Loading...
Searching...
No Matches
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-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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 "List.H"
31#include "hexRef8Data.H"
32#include "mapPolyMesh.H"
34#include "polyMesh.H"
35#include "syncTools.H"
36#include "refinementHistory.H"
37#include "fvMesh.H"
38
39// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40
42{
43 if
44 (
45 // Or: !io.anyRead() ?
46 io.readOpt() != IOobject::MUST_READ
47 && io.readOpt() != IOobject::LAZY_READ
48 )
49 {
50 return;
51 }
52
53 {
54 typedef labelIOList Type;
55 IOobject rio(io, "cellLevel");
56
57 // haveFile
58 if (returnReduceOr(rio.typeHeaderOk<Type>(true)))
59 {
60 Info<< "Reading hexRef8 data : " << rio.name() << endl;
61 cellLevelPtr_.reset(new Type(rio));
62 }
63 }
64 {
65 typedef labelIOList Type;
66 IOobject rio(io, "pointLevel");
67
68 // haveFile
69 if (returnReduceOr(rio.typeHeaderOk<Type>(true)))
70 {
71 Info<< "Reading hexRef8 data : " << rio.name() << endl;
72 pointLevelPtr_.reset(new Type(rio));
73 }
74 }
75 {
77 IOobject rio(io, "level0Edge");
78
79 // haveFile
80 if (returnReduceOr(rio.typeHeaderOk<Type>(true)))
81 {
82 Info<< "Reading hexRef8 data : " << rio.name() << endl;
83 level0EdgePtr_.reset(new Type(rio));
84 }
85 }
86 {
87 typedef refinementHistory Type;
88 IOobject rio(io, "refinementHistory");
89
90 // haveFile
91 if (returnReduceOr(rio.typeHeaderOk<Type>(true)))
92 {
93 Info<< "Reading hexRef8 data : " << rio.name() << endl;
94 refHistoryPtr_.reset(new Type(rio));
95 }
96 }
97}
98
99
101(
102 const IOobject& io,
103 const hexRef8Data& data,
104 const labelUList& cellMap,
105 const labelUList& pointMap
106)
107{
108 if (data.cellLevelPtr_)
109 {
110 IOobject rio(io, data.cellLevelPtr_().name());
111
112 cellLevelPtr_.reset
113 (
114 new labelIOList
115 (
116 rio,
117 labelUIndList(data.cellLevelPtr_(), cellMap)()
118 )
119 );
120 }
121 if (data.pointLevelPtr_)
122 {
123 IOobject rio(io, data.pointLevelPtr_().name());
124
125 pointLevelPtr_.reset
126 (
127 new labelIOList
128 (
129 rio,
130 labelUIndList(data.pointLevelPtr_(), pointMap)()
131 )
132 );
133 }
134 if (data.level0EdgePtr_)
135 {
136 IOobject rio(io, data.level0EdgePtr_().name());
137
138 level0EdgePtr_.reset
139 (
140 new uniformDimensionedScalarField(rio, data.level0EdgePtr_())
141 );
142 }
143 if (data.refHistoryPtr_)
144 {
145 IOobject rio(io, data.refHistoryPtr_().name());
146
147 refHistoryPtr_ = data.refHistoryPtr_().clone(rio, cellMap);
148 }
149}
150
151
153(
154 const IOobject& io,
155 const UPtrList<const labelList>& cellMaps,
156 const UPtrList<const labelList>& pointMaps,
157 const UPtrList<const hexRef8Data>& procDatas
158)
159{
160 const polyMesh& mesh = dynamic_cast<const polyMesh&>(io.db());
161
162 // cellLevel
163
164 if (procDatas[0].cellLevelPtr_)
165 {
166 IOobject rio(io, procDatas[0].cellLevelPtr_().name());
167
168 cellLevelPtr_.reset(new labelIOList(rio, mesh.nCells()));
169 auto& cellLevel = *cellLevelPtr_;
170
171 forAll(procDatas, procI)
172 {
173 const labelList& procCellLevel = procDatas[procI].cellLevelPtr_();
174 labelUIndList(cellLevel, cellMaps[procI]) = procCellLevel;
175 }
176 }
177
178
179 // pointLevel
180
181 if (procDatas[0].pointLevelPtr_)
182 {
183 IOobject rio(io, procDatas[0].pointLevelPtr_().name());
184
185 pointLevelPtr_.reset(new labelIOList(rio, mesh.nPoints()));
186 auto& pointLevel = *pointLevelPtr_;
187
188 forAll(procDatas, procI)
189 {
190 const labelList& procPointLevel = procDatas[procI].pointLevelPtr_();
191 labelUIndList(pointLevel, pointMaps[procI]) = procPointLevel;
192 }
193 }
194
195
196 // level0Edge
197
198 if (procDatas[0].level0EdgePtr_)
199 {
200 IOobject rio(io, procDatas[0].level0EdgePtr_().name());
201
202 level0EdgePtr_.reset
203 (
205 (
206 rio,
207 procDatas[0].level0EdgePtr_()
208 )
209 );
210 }
211
212
213 // refinementHistory
214
215 if (procDatas[0].refHistoryPtr_)
216 {
217 IOobject rio(io, procDatas[0].refHistoryPtr_().name());
218
219 UPtrList<const refinementHistory> procRefs(procDatas.size());
220 forAll(procDatas, i)
221 {
222 procRefs.set(i, &procDatas[i].refHistoryPtr_());
223 }
224
225 refHistoryPtr_.reset
226 (
228 (
229 rio,
230 cellMaps,
231 procRefs
232 )
233 );
234 }
235}
236
237
238// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
241{} // refinementHistory forward declared
242
243
244// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
245
247{
248 const polyMesh& mesh = dynamic_cast<const polyMesh&>(io.db());
249
250 if (returnReduceOr(cellLevelPtr_) && !cellLevelPtr_)
251 {
252 IOobject rio(io, "cellLevel");
254 cellLevelPtr_.reset
255 (
256 new labelIOList(rio, labelList(mesh.nCells(), Zero))
257 );
258 }
259 // Take over IOobject settings
260 if (cellLevelPtr_)
261 {
262 cellLevelPtr_->IOobjectOption::operator=(io);
263 }
264
265 if (returnReduceOr(pointLevelPtr_) && !pointLevelPtr_)
266 {
267 IOobject rio(io, "pointLevel");
269 pointLevelPtr_.reset
270 (
272 );
273 }
274 // Take over IOobject settings
275 if (pointLevelPtr_)
276 {
277 pointLevelPtr_->IOobjectOption::operator=(io);
278 }
279
280 if (returnReduceOr(level0EdgePtr_))
281 {
282 // Get master length
283 scalar masterLen = (UPstream::master() ? level0EdgePtr_().value() : 0);
284 Pstream::broadcast(masterLen);
285 if (!level0EdgePtr_)
286 {
287 IOobject rio(io, "level0Edge");
289 level0EdgePtr_.reset
290 (
291 new uniformDimensionedScalarField(rio, dimLength, masterLen)
292 );
293 }
294 }
295 // Take over IOobject settings
296 if (level0EdgePtr_)
297 {
298 level0EdgePtr_->IOobjectOption::operator=(io);
299 }
300
301 if (returnReduceOr(refHistoryPtr_) && !refHistoryPtr_)
302 {
303 IOobject rio(io, "refinementHistory");
305 refHistoryPtr_.reset(new refinementHistory(rio, mesh.nCells(), true));
306 }
307 // Take over IOobject settings
308 if (refHistoryPtr_)
309 {
310 refHistoryPtr_->IOobjectOption::operator=(io);
311 }
312}
313
314
315void Foam::hexRef8Data::updateMesh(const mapPolyMesh& map)
316{
317 // Sanity check
318 if
319 (
320 (cellLevelPtr_ && cellLevelPtr_().size() != map.nOldCells())
321 || (pointLevelPtr_ && pointLevelPtr_().size() != map.nOldPoints())
322 )
323 {
324 cellLevelPtr_.clear();
325 pointLevelPtr_.clear();
326 level0EdgePtr_.clear();
327 refHistoryPtr_.clear();
328 return;
329 }
330
331
332 if (cellLevelPtr_)
333 {
334 const labelList& cellMap = map.cellMap();
335 labelList& cellLevel = cellLevelPtr_();
336
337 labelList newCellLevel(cellMap.size());
338 forAll(cellMap, newCelli)
339 {
340 label oldCelli = cellMap[newCelli];
341
342 if (oldCelli == -1)
343 {
344 newCellLevel[newCelli] = 0;
345 }
346 else
347 {
348 newCellLevel[newCelli] = cellLevel[oldCelli];
349 }
350 }
351 cellLevel.transfer(newCellLevel);
352 cellLevelPtr_().instance() = map.mesh().facesInstance();
353 }
354 if (pointLevelPtr_)
355 {
356 const labelList& pointMap = map.pointMap();
357 labelList& pointLevel = pointLevelPtr_();
358
359 labelList newPointLevel(pointMap.size());
360 forAll(pointMap, newPointi)
361 {
362 label oldPointi = pointMap[newPointi];
363
364 if (oldPointi == -1)
365 {
366 newPointLevel[newPointi] = 0;
367 }
368 else
369 {
370 newPointLevel[newPointi] = pointLevel[oldPointi];
371 }
372 }
373 pointLevel.transfer(newPointLevel);
374 pointLevelPtr_().instance() = map.mesh().facesInstance();
375 }
376
377
378 if (refHistoryPtr_ && refHistoryPtr_().active())
380 refHistoryPtr_().updateMesh(map);
381 refHistoryPtr_().instance() = map.mesh().facesInstance();
382 }
383}
384
385
386void Foam::hexRef8Data::distribute(const mapDistributePolyMesh& map)
387{
388 if (cellLevelPtr_)
389 {
390 map.cellMap().distribute(*cellLevelPtr_);
391 }
392 if (pointLevelPtr_)
393 {
394 map.pointMap().distribute(*pointLevelPtr_);
395 }
396
397 // No need to distribute the level0Edge
398
399 if (refHistoryPtr_ && refHistoryPtr_().active())
400 {
401 refHistoryPtr_().distribute(map);
402 }
403}
404
405
406bool Foam::hexRef8Data::write() const
407{
408 bool ok = true;
409 if (cellLevelPtr_)
410 {
411 ok = ok && cellLevelPtr_().write();
412 }
413 if (pointLevelPtr_)
414 {
415 ok = ok && pointLevelPtr_().write();
416 }
417 if (level0EdgePtr_)
418 {
419 ok = ok && level0EdgePtr_().write();
420 }
421 if (refHistoryPtr_)
422 {
423 ok = ok && refHistoryPtr_().write();
424 }
425 return ok;
426}
427
428
429// ************************************************************************* //
void transfer(HashTable< T, Key, Hash > &rhs)
Transfer contents into this table.
Definition HashTable.C:794
readOption readOpt() const noexcept
Get the read option.
@ NO_READ
Nothing to be read.
@ MUST_READ
Reading required.
@ LAZY_READ
Reading is optional [identical to READ_IF_PRESENT].
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
bool typeHeaderOk(const bool checkType=true, const bool search=true, const bool verbose=true)
Read header (respects is_globalIOobject trait) and check its info. A void type suppresses trait and t...
void transfer(List< T > &list)
Transfer the contents of the argument List into this list and annul the argument list.
Definition List.C:347
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition UPstream.H:1714
@ broadcast
broadcast [MPI]
Definition UPstream.H:189
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition UPtrList.H:101
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie,...
Definition UPtrList.H:366
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
Various for reading/decomposing/reconstructing/distributing refinement data.
Definition hexRef8Data.H:57
hexRef8Data(const hexRef8Data &)=delete
No copy construct.
void distribute(const mapDistributePolyMesh &)
In-place distribute.
void sync(const IOobject &io)
Parallel synchronise. This enforces valid objects on all processors (even if they don't have a mesh)....
~hexRef8Data()
Destructor.
bool write() const
Write.
void updateMesh(const mapPolyMesh &)
Update local numbering for changed mesh.
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
const mapDistribute & cellMap() const noexcept
Cell distribute map.
const mapDistribute & pointMap() const noexcept
Point distribute map.
void distribute(List< T > &fld, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Distribute List data using default commsType, default flip/negate operator.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
label nOldCells() const noexcept
Number of old cells.
const polyMesh & mesh() const noexcept
Return polyMesh.
const labelList & cellMap() const noexcept
Old cell map.
const labelList & pointMap() const noexcept
Old point map.
label nOldPoints() const noexcept
Number of old points.
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition polyMesh.C:844
label nPoints() const noexcept
Number of mesh points.
label nCells() const noexcept
Number of mesh cells.
void reset(const label nPoints, const label nInternalFaces, const label nFaces, const label nCells)
Reset this primitiveMesh given the primitive array sizes.
All refinement history. Used in unrefinement.
dynamicFvMesh & mesh
const auto & io
auto & name
bool returnReduceOr(const bool value, const int communicator=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
List< label > labelList
A List of labels.
Definition List.H:62
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
UIndirectList< label > labelUIndList
UIndirectList of labels.
IOList< label > labelIOList
IO for a List of label.
Definition labelIOList.H:32
UniformDimensionedField< scalar > uniformDimensionedScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
UList< label > labelUList
A UList of labels.
Definition UList.H:75
label newPointi
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299