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