hexRef8Data.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2015-2017 OpenFOAM Foundation
9  Copyright (C) 2017 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "IOobject.H"
30 #include "UList.H"
31 
32 #include "hexRef8Data.H"
33 #include "mapPolyMesh.H"
34 #include "mapDistributePolyMesh.H"
35 #include "polyMesh.H"
36 #include "syncTools.H"
37 #include "refinementHistory.H"
38 #include "fvMesh.H"
39 
40 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
41 
42 Foam::hexRef8Data::hexRef8Data(const IOobject& io)
43 {
44  {
45  IOobject rio(io);
46  rio.rename("cellLevel");
47  bool haveFile = returnReduce
48  (
49  rio.typeHeaderOk<labelIOList>(true),
50  orOp<bool>()
51  );
52  if (haveFile)
53  {
54  Info<< "Reading hexRef8 data : " << rio.name() << endl;
55  cellLevelPtr_.reset(new labelIOList(rio));
56  }
57  }
58  {
59  IOobject rio(io);
60  rio.rename("pointLevel");
61  bool haveFile = returnReduce
62  (
63  rio.typeHeaderOk<labelIOList>(true),
64  orOp<bool>()
65  );
66  if (haveFile)
67  {
68  Info<< "Reading hexRef8 data : " << rio.name() << endl;
69  pointLevelPtr_.reset(new labelIOList(rio));
70  }
71  }
72  {
73  IOobject rio(io);
74  rio.rename("level0Edge");
75  bool haveFile = returnReduce
76  (
78  orOp<bool>()
79  );
80  if (haveFile)
81  {
82  Info<< "Reading hexRef8 data : " << rio.name() << endl;
83  level0EdgePtr_.reset(new uniformDimensionedScalarField(rio));
84  }
85  }
86  {
87  IOobject rio(io);
88  rio.rename("refinementHistory");
89  bool haveFile = returnReduce
90  (
92  orOp<bool>()
93  );
94  if (haveFile)
95  {
96  Info<< "Reading hexRef8 data : " << rio.name() << endl;
97  refHistoryPtr_.reset(new refinementHistory(rio));
98  }
99  }
100 }
101 
102 
103 Foam::hexRef8Data::hexRef8Data
104 (
105  const IOobject& io,
106  const hexRef8Data& data,
107  const labelList& cellMap,
108  const labelList& pointMap
109 )
110 {
111  if (data.cellLevelPtr_.valid())
112  {
113  IOobject rio(io);
114  rio.rename(data.cellLevelPtr_().name());
115 
116  cellLevelPtr_.reset
117  (
118  new labelIOList
119  (
120  rio,
121  labelUIndList(data.cellLevelPtr_(), cellMap)()
122  )
123  );
124  }
125  if (data.pointLevelPtr_.valid())
126  {
127  IOobject rio(io);
128  rio.rename(data.pointLevelPtr_().name());
129 
130  pointLevelPtr_.reset
131  (
132  new labelIOList
133  (
134  rio,
135  labelUIndList(data.pointLevelPtr_(), pointMap)()
136  )
137  );
138  }
139  if (data.level0EdgePtr_.valid())
140  {
141  IOobject rio(io);
142  rio.rename(data.level0EdgePtr_().name());
143 
144  level0EdgePtr_.reset
145  (
146  new uniformDimensionedScalarField(rio, data.level0EdgePtr_())
147  );
148  }
149  if (data.refHistoryPtr_.valid())
150  {
151  IOobject rio(io);
152  rio.rename(data.refHistoryPtr_().name());
153 
154  refHistoryPtr_ = data.refHistoryPtr_().clone(rio, cellMap);
155  }
156 }
157 
158 
159 Foam::hexRef8Data::hexRef8Data
160 (
161  const IOobject& io,
162  const UPtrList<const labelList>& cellMaps,
163  const UPtrList<const labelList>& pointMaps,
164  const UPtrList<const hexRef8Data>& procDatas
165 )
166 {
167  const polyMesh& mesh = dynamic_cast<const polyMesh&>(io.db());
168 
169  // cellLevel
170 
171  if (procDatas[0].cellLevelPtr_.valid())
172  {
173  IOobject rio(io);
174  rio.rename(procDatas[0].cellLevelPtr_().name());
175 
176  cellLevelPtr_.reset(new labelIOList(rio, mesh.nCells()));
177  labelList& cellLevel = cellLevelPtr_();
178 
179  forAll(procDatas, procI)
180  {
181  const labelList& procCellLevel = procDatas[procI].cellLevelPtr_();
182  labelUIndList(cellLevel, cellMaps[procI]) = procCellLevel;
183  }
184  }
185 
186 
187  // pointLevel
188 
189  if (procDatas[0].pointLevelPtr_.valid())
190  {
191  IOobject rio(io);
192  rio.rename(procDatas[0].pointLevelPtr_().name());
193 
194  pointLevelPtr_.reset(new labelIOList(rio, mesh.nPoints()));
195  labelList& pointLevel = pointLevelPtr_();
196 
197  forAll(procDatas, procI)
198  {
199  const labelList& procPointLevel = procDatas[procI].pointLevelPtr_();
200  labelUIndList(pointLevel, pointMaps[procI]) = procPointLevel;
201  }
202  }
203 
204 
205  // level0Edge
206 
207  if (procDatas[0].level0EdgePtr_.valid())
208  {
209  IOobject rio(io);
210  rio.rename(procDatas[0].level0EdgePtr_().name());
211 
212  level0EdgePtr_.reset
213  (
215  (
216  rio,
217  procDatas[0].level0EdgePtr_()
218  )
219  );
220  }
221 
222 
223  // refinementHistory
224 
225  if (procDatas[0].refHistoryPtr_.valid())
226  {
227  IOobject rio(io);
228  rio.rename(procDatas[0].refHistoryPtr_().name());
229 
230  UPtrList<const refinementHistory> procRefs(procDatas.size());
231  forAll(procDatas, i)
232  {
233  procRefs.set(i, &procDatas[i].refHistoryPtr_());
234  }
235 
236  refHistoryPtr_.reset
237  (
239  (
240  rio,
241  cellMaps,
242  procRefs
243  )
244  );
245  }
246 }
247 
248 
249 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
250 
252 {}
253 
254 
255 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
256 
258 {
259  const polyMesh& mesh = dynamic_cast<const polyMesh&>(io.db());
260 
261  bool hasCellLevel = returnReduce(cellLevelPtr_.valid(), orOp<bool>());
262  if (hasCellLevel && !cellLevelPtr_.valid())
263  {
264  IOobject rio(io);
265  rio.rename("cellLevel");
266  rio.readOpt() = IOobject::NO_READ;
267  cellLevelPtr_.reset
268  (
269  new labelIOList(rio, labelList(mesh.nCells(), Zero))
270  );
271  }
272 
273  bool hasPointLevel = returnReduce(pointLevelPtr_.valid(), orOp<bool>());
274  if (hasPointLevel && !pointLevelPtr_.valid())
275  {
276  IOobject rio(io);
277  rio.rename("pointLevel");
278  rio.readOpt() = IOobject::NO_READ;
279  pointLevelPtr_.reset
280  (
281  new labelIOList(rio, labelList(mesh.nPoints(), Zero))
282  );
283  }
284 
285  bool hasLevel0Edge = returnReduce(level0EdgePtr_.valid(), orOp<bool>());
286  if (hasLevel0Edge)
287  {
288  // Get master length
289  scalar masterLen = (Pstream::master() ? level0EdgePtr_().value() : 0);
290  Pstream::scatter(masterLen);
291  if (!level0EdgePtr_.valid())
292  {
293  IOobject rio(io);
294  rio.rename("level0Edge");
295  rio.readOpt() = IOobject::NO_READ;
296  level0EdgePtr_.reset
297  (
299  (
300  rio,
301  dimensionedScalar(rio.name(), dimLength, masterLen)
302  )
303  );
304  }
305  }
306 
307  bool hasHistory = returnReduce(refHistoryPtr_.valid(), orOp<bool>());
308  if (hasHistory && !refHistoryPtr_.valid())
309  {
310  IOobject rio(io);
311  rio.rename("refinementHistory");
312  rio.readOpt() = IOobject::NO_READ;
313  refHistoryPtr_.reset(new refinementHistory(rio, mesh.nCells(), true));
314  }
315 }
316 
317 
319 {
320  // Sanity check
321  if
322  (
323  (cellLevelPtr_.valid() && cellLevelPtr_().size() != map.nOldCells())
324  || (pointLevelPtr_.valid() && pointLevelPtr_().size() != map.nOldPoints())
325  )
326  {
327  cellLevelPtr_.clear();
328  pointLevelPtr_.clear();
329  level0EdgePtr_.clear();
330  refHistoryPtr_.clear();
331  return;
332  }
333 
334 
335  if (cellLevelPtr_.valid())
336  {
337  const labelList& cellMap = map.cellMap();
338  labelList& cellLevel = cellLevelPtr_();
339 
340  labelList newCellLevel(cellMap.size());
341  forAll(cellMap, newCelli)
342  {
343  label oldCelli = cellMap[newCelli];
344 
345  if (oldCelli == -1)
346  {
347  newCellLevel[newCelli] = 0;
348  }
349  else
350  {
351  newCellLevel[newCelli] = cellLevel[oldCelli];
352  }
353  }
354  cellLevel.transfer(newCellLevel);
355  cellLevelPtr_().instance() = map.mesh().facesInstance();
356  }
357  if (pointLevelPtr_.valid())
358  {
359  const labelList& pointMap = map.pointMap();
360  labelList& pointLevel = pointLevelPtr_();
361 
362  labelList newPointLevel(pointMap.size());
363  forAll(pointMap, newPointi)
364  {
365  label oldPointi = pointMap[newPointi];
366 
367  if (oldPointi == -1)
368  {
369  newPointLevel[newPointi] = 0;
370  }
371  else
372  {
373  newPointLevel[newPointi] = pointLevel[oldPointi];
374  }
375  }
376  pointLevel.transfer(newPointLevel);
377  pointLevelPtr_().instance() = map.mesh().facesInstance();
378  }
379 
380 
381  if (refHistoryPtr_.valid() && refHistoryPtr_().active())
382  {
383  refHistoryPtr_().updateMesh(map);
384  refHistoryPtr_().instance() = map.mesh().facesInstance();
385  }
386 }
387 
388 
390 {
391  if (cellLevelPtr_.valid())
392  {
393  map.cellMap().distribute(cellLevelPtr_());
394  }
395  if (pointLevelPtr_.valid())
396  {
397  map.pointMap().distribute(pointLevelPtr_());
398  }
399 
400  // No need to distribute the level0Edge
401 
402  if (refHistoryPtr_.valid() && refHistoryPtr_().active())
403  {
404  refHistoryPtr_().distribute(map);
405  }
406 }
407 
408 
410 {
411  bool ok = true;
412  if (cellLevelPtr_.valid())
413  {
414  ok = ok && cellLevelPtr_().write();
415  }
416  if (pointLevelPtr_.valid())
417  {
418  ok = ok && pointLevelPtr_().write();
419  }
420  if (level0EdgePtr_.valid())
421  {
422  ok = ok && level0EdgePtr_().write();
423  }
424  if (refHistoryPtr_.valid())
425  {
426  ok = ok && refHistoryPtr_().write();
427  }
428  return ok;
429 }
430 
431 
432 // ************************************************************************* //
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:74
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::UPtrList::size
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:90
Foam::IOobject::name
const word & name() const
Return name.
Definition: IOobjectI.H:46
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:94
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::IOobject::typeHeaderOk
bool typeHeaderOk(const bool checkType=true, const bool search=true, const bool verbose=true)
Read header (uses typeFilePath to find file) and check its info.
Definition: IOobjectTemplates.C:39
mapPolyMesh.H
Foam::polyMesh::facesInstance
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:821
Foam::Pstream::scatter
static void scatter(const List< commsStruct > &comms, T &Value, const int tag, const label comm)
Scatter data. Distribute without modification. Reverse of gather.
Definition: gatherScatter.C:150
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::baseIOdictionary::name
const word & name() const
Name function is needed to disambiguate those inherited.
Definition: baseIOdictionary.C:88
polyMesh.H
Foam::labelIOList
IOList< label > labelIOList
Label container classes.
Definition: labelIOList.H:44
syncTools.H
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::mapDistributePolyMesh::cellMap
const mapDistribute & cellMap() const
Cell distribute map.
Definition: mapDistributePolyMesh.H:223
Foam::IOobject::db
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:432
Foam::refinementHistory
All refinement history. Used in unrefinement.
Definition: refinementHistory.H:106
mapDistributePolyMesh.H
Foam::primitiveMesh::nCells
label nCells() const
Number of mesh cells.
Definition: primitiveMeshI.H:96
Foam::UniformDimensionedField< scalar >
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:62
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::mapPolyMesh::nOldPoints
label nOldPoints() const
Number of old points.
Definition: mapPolyMesh.H:368
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::uniformDimensionedScalarField
UniformDimensionedField< scalar > uniformDimensionedScalarField
Definition: uniformDimensionedFields.H:49
Foam::List::transfer
void transfer(List< T > &list)
Definition: List.C:436
Foam::hexRef8Data::updateMesh
void updateMesh(const mapPolyMesh &)
Update local numbering for changed mesh.
Definition: hexRef8Data.C:318
Foam::UPtrList
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition: UPtrList.H:63
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
Foam::hexRef8Data::distribute
void distribute(const mapDistributePolyMesh &)
In-place distribute.
Definition: hexRef8Data.C:389
Foam::hexRef8Data
Various for reading/decomposing/reconstructing/distributing refinement data.
Definition: hexRef8Data.H:60
Foam::mapDistribute::distribute
void distribute(List< T > &fld, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Distribute data using default commsType.
Definition: mapDistributeTemplates.C:152
newPointi
label newPointi
Definition: readKivaGrid.H:496
IOobject.H
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
hexRef8Data.H
fvMesh.H
Foam::IOobject::rename
virtual void rename(const word &newName)
Rename.
Definition: IOobject.H:368
Foam::hexRef8Data::sync
void sync(const IOobject &io)
Parallel synchronise. This enforces valid objects on all processors.
Definition: hexRef8Data.C:257
Foam::UPstream::master
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:438
Foam::hexRef8Data::write
bool write() const
Write.
Definition: hexRef8Data.C:409
UList.H
Foam::List< label >
Foam::mapPolyMesh::pointMap
const labelList & pointMap() const
Old point map.
Definition: mapPolyMesh.H:395
Foam::IOobject::readOpt
readOption readOpt() const
The read option.
Definition: IOobjectI.H:141
Foam::mapPolyMesh::nOldCells
label nOldCells() const
Number of old cells.
Definition: mapPolyMesh.H:386
Foam::dictionary::clone
autoPtr< dictionary > clone() const
Construct and return clone.
Definition: dictionary.C:178
Foam::primitiveMesh::nPoints
label nPoints() const
Number of mesh points.
Definition: primitiveMeshI.H:37
Foam::IOList< label >
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:160
Foam::hexRef8Data::~hexRef8Data
~hexRef8Data()
Destructor.
Definition: hexRef8Data.C:251
Foam::mapDistributePolyMesh::pointMap
const mapDistribute & pointMap() const
Point distribute map.
Definition: mapDistributePolyMesh.H:211
refinementHistory.H
Foam::mapDistributePolyMesh
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Definition: mapDistributePolyMesh.H:66
Foam::mapPolyMesh::mesh
const polyMesh & mesh() const
Return polyMesh.
Definition: mapPolyMesh.H:362
Foam::IOobject::NO_READ
Definition: IOobject.H:123
Foam::orOp
Definition: ops.H:234
Foam::data
Database for solution data, solver performance and other reduced data.
Definition: data.H:54
Foam::mapPolyMesh::cellMap
const labelList & cellMap() const
Old cell map.
Definition: mapPolyMesh.H:434
Foam::labelUIndList
UIndirectList< label > labelUIndList
UIndirectList of labels.
Definition: UIndirectList.H:59