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-2021 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 #include "hexRef8Data.H"
32 #include "mapPolyMesh.H"
33 #include "mapDistributePolyMesh.H"
34 #include "polyMesh.H"
35 #include "syncTools.H"
36 #include "refinementHistory.H"
37 #include "fvMesh.H"
38 
39 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40 
41 Foam::hexRef8Data::hexRef8Data(const IOobject& io)
42 {
43  {
44  typedef labelIOList Type;
45  IOobject rio(io, "cellLevel");
46 
47  // haveFile
48  if (returnReduce(rio.typeHeaderOk<Type>(true), orOp<bool>()))
49  {
50  Info<< "Reading hexRef8 data : " << rio.name() << endl;
51  cellLevelPtr_.reset(new Type(rio));
52  }
53  }
54  {
55  typedef labelIOList Type;
56  IOobject rio(io, "pointLevel");
57 
58  // haveFile
59  if (returnReduce(rio.typeHeaderOk<Type>(true), orOp<bool>()))
60  {
61  Info<< "Reading hexRef8 data : " << rio.name() << endl;
62  pointLevelPtr_.reset(new Type(rio));
63  }
64  }
65  {
66  typedef uniformDimensionedScalarField Type;
67  IOobject rio(io, "level0Edge");
68 
69  // haveFile
70  if (returnReduce(rio.typeHeaderOk<Type>(true), orOp<bool>()))
71  {
72  Info<< "Reading hexRef8 data : " << rio.name() << endl;
73  level0EdgePtr_.reset(new Type(rio));
74  }
75  }
76  {
77  typedef refinementHistory Type;
78  IOobject rio(io, "refinementHistory");
79 
80  // haveFile
81  if (returnReduce(rio.typeHeaderOk<Type>(true), orOp<bool>()))
82  {
83  Info<< "Reading hexRef8 data : " << rio.name() << endl;
84  refHistoryPtr_.reset(new Type(rio));
85  }
86  }
87 }
88 
89 
90 Foam::hexRef8Data::hexRef8Data
91 (
92  const IOobject& io,
93  const hexRef8Data& data,
94  const labelList& cellMap,
95  const labelList& pointMap
96 )
97 {
98  if (data.cellLevelPtr_)
99  {
100  IOobject rio(io, data.cellLevelPtr_().name());
101 
102  cellLevelPtr_.reset
103  (
104  new labelIOList
105  (
106  rio,
107  labelUIndList(data.cellLevelPtr_(), cellMap)()
108  )
109  );
110  }
111  if (data.pointLevelPtr_)
112  {
113  IOobject rio(io, data.pointLevelPtr_().name());
114 
115  pointLevelPtr_.reset
116  (
117  new labelIOList
118  (
119  rio,
120  labelUIndList(data.pointLevelPtr_(), pointMap)()
121  )
122  );
123  }
124  if (data.level0EdgePtr_)
125  {
126  IOobject rio(io, data.level0EdgePtr_().name());
127 
128  level0EdgePtr_.reset
129  (
130  new uniformDimensionedScalarField(rio, data.level0EdgePtr_())
131  );
132  }
133  if (data.refHistoryPtr_)
134  {
135  IOobject rio(io, data.refHistoryPtr_().name());
136 
137  refHistoryPtr_ = data.refHistoryPtr_().clone(rio, cellMap);
138  }
139 }
140 
141 
142 Foam::hexRef8Data::hexRef8Data
143 (
144  const IOobject& io,
145  const UPtrList<const labelList>& cellMaps,
146  const UPtrList<const labelList>& pointMaps,
147  const UPtrList<const hexRef8Data>& procDatas
148 )
149 {
150  const polyMesh& mesh = dynamic_cast<const polyMesh&>(io.db());
151 
152  // cellLevel
153 
154  if (procDatas[0].cellLevelPtr_)
155  {
156  IOobject rio(io, procDatas[0].cellLevelPtr_().name());
157 
158  cellLevelPtr_.reset(new labelIOList(rio, mesh.nCells()));
159  auto& cellLevel = *cellLevelPtr_;
160 
161  forAll(procDatas, procI)
162  {
163  const labelList& procCellLevel = procDatas[procI].cellLevelPtr_();
164  labelUIndList(cellLevel, cellMaps[procI]) = procCellLevel;
165  }
166  }
167 
168 
169  // pointLevel
170 
171  if (procDatas[0].pointLevelPtr_)
172  {
173  IOobject rio(io, procDatas[0].pointLevelPtr_().name());
174 
175  pointLevelPtr_.reset(new labelIOList(rio, mesh.nPoints()));
176  auto& pointLevel = *pointLevelPtr_;
177 
178  forAll(procDatas, procI)
179  {
180  const labelList& procPointLevel = procDatas[procI].pointLevelPtr_();
181  labelUIndList(pointLevel, pointMaps[procI]) = procPointLevel;
182  }
183  }
184 
185 
186  // level0Edge
187 
188  if (procDatas[0].level0EdgePtr_)
189  {
190  IOobject rio(io, procDatas[0].level0EdgePtr_().name());
191 
192  level0EdgePtr_.reset
193  (
195  (
196  rio,
197  procDatas[0].level0EdgePtr_()
198  )
199  );
200  }
201 
202 
203  // refinementHistory
204 
205  if (procDatas[0].refHistoryPtr_)
206  {
207  IOobject rio(io, procDatas[0].refHistoryPtr_().name());
208 
209  UPtrList<const refinementHistory> procRefs(procDatas.size());
210  forAll(procDatas, i)
211  {
212  procRefs.set(i, &procDatas[i].refHistoryPtr_());
213  }
214 
215  refHistoryPtr_.reset
216  (
218  (
219  rio,
220  cellMaps,
221  procRefs
222  )
223  );
224  }
225 }
226 
227 
228 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
229 
231 {} // refinementHistory forward declared
232 
233 
234 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
235 
237 {
238  const polyMesh& mesh = dynamic_cast<const polyMesh&>(io.db());
239 
240  bool hasCellLevel = returnReduce(bool(cellLevelPtr_), orOp<bool>());
241  if (hasCellLevel && !cellLevelPtr_)
242  {
243  IOobject rio(io, "cellLevel");
245  cellLevelPtr_.reset
246  (
247  new labelIOList(rio, labelList(mesh.nCells(), Zero))
248  );
249  }
250 
251  bool hasPointLevel = returnReduce(bool(pointLevelPtr_), orOp<bool>());
252  if (hasPointLevel && !pointLevelPtr_)
253  {
254  IOobject rio(io, "pointLevel");
256  pointLevelPtr_.reset
257  (
258  new labelIOList(rio, labelList(mesh.nPoints(), Zero))
259  );
260  }
261 
262  bool hasLevel0Edge = returnReduce(bool(level0EdgePtr_), orOp<bool>());
263  if (hasLevel0Edge)
264  {
265  // Get master length
266  scalar masterLen = (Pstream::master() ? level0EdgePtr_().value() : 0);
267  Pstream::scatter(masterLen);
268  if (!level0EdgePtr_)
269  {
270  IOobject rio(io, "level0Edge");
272  level0EdgePtr_.reset
273  (
275  (
276  rio,
277  dimensionedScalar(rio.name(), dimLength, masterLen)
278  )
279  );
280  }
281  }
282 
283  bool hasHistory = returnReduce(bool(refHistoryPtr_), orOp<bool>());
284  if (hasHistory && !refHistoryPtr_)
285  {
286  IOobject rio(io, "refinementHistory");
288  refHistoryPtr_.reset(new refinementHistory(rio, mesh.nCells(), true));
289  }
290 }
291 
292 
294 {
295  // Sanity check
296  if
297  (
298  (cellLevelPtr_ && cellLevelPtr_().size() != map.nOldCells())
299  || (pointLevelPtr_ && pointLevelPtr_().size() != map.nOldPoints())
300  )
301  {
302  cellLevelPtr_.clear();
303  pointLevelPtr_.clear();
304  level0EdgePtr_.clear();
305  refHistoryPtr_.clear();
306  return;
307  }
308 
309 
310  if (cellLevelPtr_)
311  {
312  const labelList& cellMap = map.cellMap();
313  labelList& cellLevel = cellLevelPtr_();
314 
315  labelList newCellLevel(cellMap.size());
316  forAll(cellMap, newCelli)
317  {
318  label oldCelli = cellMap[newCelli];
319 
320  if (oldCelli == -1)
321  {
322  newCellLevel[newCelli] = 0;
323  }
324  else
325  {
326  newCellLevel[newCelli] = cellLevel[oldCelli];
327  }
328  }
329  cellLevel.transfer(newCellLevel);
330  cellLevelPtr_().instance() = map.mesh().facesInstance();
331  }
332  if (pointLevelPtr_)
333  {
334  const labelList& pointMap = map.pointMap();
335  labelList& pointLevel = pointLevelPtr_();
336 
337  labelList newPointLevel(pointMap.size());
338  forAll(pointMap, newPointi)
339  {
340  label oldPointi = pointMap[newPointi];
341 
342  if (oldPointi == -1)
343  {
344  newPointLevel[newPointi] = 0;
345  }
346  else
347  {
348  newPointLevel[newPointi] = pointLevel[oldPointi];
349  }
350  }
351  pointLevel.transfer(newPointLevel);
352  pointLevelPtr_().instance() = map.mesh().facesInstance();
353  }
354 
355 
356  if (refHistoryPtr_ && refHistoryPtr_().active())
357  {
358  refHistoryPtr_().updateMesh(map);
359  refHistoryPtr_().instance() = map.mesh().facesInstance();
360  }
361 }
362 
363 
365 {
366  if (cellLevelPtr_)
367  {
368  map.cellMap().distribute(*cellLevelPtr_);
369  }
370  if (pointLevelPtr_)
371  {
372  map.pointMap().distribute(*pointLevelPtr_);
373  }
374 
375  // No need to distribute the level0Edge
376 
377  if (refHistoryPtr_ && refHistoryPtr_().active())
378  {
379  refHistoryPtr_().distribute(map);
380  }
381 }
382 
383 
385 {
386  bool ok = true;
387  if (cellLevelPtr_)
388  {
389  ok = ok && cellLevelPtr_().write();
390  }
391  if (pointLevelPtr_)
392  {
393  ok = ok && pointLevelPtr_().write();
394  }
395  if (level0EdgePtr_)
396  {
397  ok = ok && level0EdgePtr_().write();
398  }
399  if (refHistoryPtr_)
400  {
401  ok = ok && refHistoryPtr_().write();
402  }
403  return ok;
404 }
405 
406 
407 // ************************************************************************* //
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::UPtrList::size
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:106
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:52
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
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:852
Foam::UPstream::master
static bool master(const label communicator=worldComm)
Am I the master process.
Definition: UPstream.H:457
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:369
Foam::baseIOdictionary::name
const word & name() const
Definition: baseIOdictionary.C:85
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
Foam::primitiveMesh::nPoints
label nPoints() const noexcept
Number of mesh points.
Definition: primitiveMeshI.H:37
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::mapDistributePolyMesh::cellMap
const mapDistribute & cellMap() const
Cell distribute map.
Definition: mapDistributePolyMesh.H:223
Foam::refinementHistory
All refinement history. Used in unrefinement.
Definition: refinementHistory.H:102
mapDistributePolyMesh.H
Foam::primitiveMesh::nCells
label nCells() const noexcept
Number of mesh cells.
Definition: primitiveMeshI.H:96
Foam::UniformDimensionedField< scalar >
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::mapPolyMesh::nOldPoints
label nOldPoints() const
Number of old points.
Definition: mapPolyMesh.H:369
Foam::uniformDimensionedScalarField
UniformDimensionedField< scalar > uniformDimensionedScalarField
Definition: uniformDimensionedFields.H:49
Foam::List::transfer
void transfer(List< T > &list)
Definition: List.C:456
Foam::hexRef8Data::updateMesh
void updateMesh(const mapPolyMesh &)
Update local numbering for changed mesh.
Definition: hexRef8Data.C:293
Foam::UPtrList
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition: UPtrList.H:62
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
Foam::hexRef8Data::distribute
void distribute(const mapDistributePolyMesh &)
In-place distribute.
Definition: hexRef8Data.C:364
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::hexRef8Data::sync
void sync(const IOobject &io)
Parallel synchronise. This enforces valid objects on all processors.
Definition: hexRef8Data.C:236
Foam::IOobject::readOpt
readOption readOpt() const noexcept
The read option.
Definition: IOobjectI.H:164
Foam::IOobject::name
const word & name() const noexcept
Return name.
Definition: IOobjectI.H:65
Foam::hexRef8Data::write
bool write() const
Write.
Definition: hexRef8Data.C:384
UList.H
Foam::List< label >
Foam::mapPolyMesh::pointMap
const labelList & pointMap() const
Old point map.
Definition: mapPolyMesh.H:396
Foam::mapPolyMesh::nOldCells
label nOldCells() const
Number of old cells.
Definition: mapPolyMesh.H:387
Foam::dictionary::clone
autoPtr< dictionary > clone() const
Construct and return clone.
Definition: dictionary.C:172
Foam::IOList< label >
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:161
Foam::hexRef8Data::~hexRef8Data
~hexRef8Data()
Destructor.
Definition: hexRef8Data.C:230
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
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:363
Foam::IOobject::NO_READ
Definition: IOobject.H:188
Foam::orOp
Definition: ops.H:234
Foam::data
Database for solution data, solver performance and other reduced data.
Definition: data.H:55
Foam::mapPolyMesh::cellMap
const labelList & cellMap() const
Old cell map.
Definition: mapPolyMesh.H:435
Foam::labelUIndList
UIndirectList< label > labelUIndList
UIndirectList of labels.
Definition: UIndirectList.H:58
Foam::IOobject::db
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition: IOobject.C:487