hexRef4Data.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 | Copyright (C) 2015 OpenFOAM Foundation
6  \\/ M anipulation | Copyright (C) 2015 OpenCFD Ltd.
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 "mapPolyMesh.H"
31 #include "mapDistributePolyMesh.H"
32 #include "polyMesh.H"
33 #include "syncTools.H"
34 #include "refinementHistory.H"
35 #include "fvMesh.H"
36 
37 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38 
40 {
41  {
42  IOobject rio(io);
43  rio.rename("cellLevel");
44  bool haveFile = returnReduce(rio.headerOk(), orOp<bool>());
45  if (haveFile)
46  {
47  Info<< "Reading hexRef8 data : " << rio.name() << endl;
48  cellLevelPtr_.reset(new labelIOList(rio));
49  }
50  }
51  {
52  IOobject rio(io);
53  rio.rename("pointLevel");
54  bool haveFile = returnReduce(rio.headerOk(), orOp<bool>());
55  if (haveFile)
56  {
57  Info<< "Reading hexRef8 data : " << rio.name() << endl;
58  pointLevelPtr_.reset(new labelIOList(rio));
59  }
60  }
61  {
62  IOobject rio(io);
63  rio.rename("level0Edge");
64 
65  // MEJ: temporarily (until global reading of UniformedFields)
66  // do not read level0Edge on processors that do not have it.
67  bool haveFile = rio.headerOk();
69  {
70  reduce(haveFile, orOp<bool>());
71  }
72 
73  if (haveFile)
74  {
75  Info<< "Reading hexRef8 data : " << rio.name() << endl;
77  }
78  }
79  {
80  IOobject rio(io);
81  rio.rename("refinementHistory");
82  bool haveFile = returnReduce(rio.headerOk(), orOp<bool>());
83  if (haveFile)
84  {
85  Info<< "Reading hexRef8 data : " << rio.name() << endl;
86  refHistoryPtr_.reset(new refinementHistory(rio));
87  }
88  }
89 }
90 
91 
93 (
94  const IOobject& io,
95  const hexRef8Data& data,
96  const labelList& cellMap,
97  const labelList& pointMap
98 )
99 {
100  if (data.cellLevelPtr_.valid())
101  {
102  IOobject rio(io);
103  rio.rename(data.cellLevelPtr_().name());
104 
105  cellLevelPtr_.reset
106  (
107  new labelIOList
108  (
109  rio,
110  UIndirectList<label>(data.cellLevelPtr_(), cellMap)()
111  )
112  );
113  }
114  if (data.pointLevelPtr_.valid())
115  {
116  IOobject rio(io);
117  rio.rename(data.pointLevelPtr_().name());
118 
119  pointLevelPtr_.reset
120  (
121  new labelIOList
122  (
123  rio,
124  UIndirectList<label>(data.pointLevelPtr_(), pointMap)()
125  )
126  );
127  }
128  if (data.level0EdgePtr_.valid())
129  {
130  IOobject rio(io);
131  rio.rename(data.level0EdgePtr_().name());
132 
133  level0EdgePtr_.reset
134  (
135  new uniformDimensionedScalarField(rio, data.level0EdgePtr_())
136  );
137  }
138  if (data.refHistoryPtr_.valid())
139  {
140  IOobject rio(io);
141  rio.rename(data.refHistoryPtr_().name());
142 
143  refHistoryPtr_ = data.refHistoryPtr_().clone(rio, cellMap);
144  }
145 }
146 
147 
149 (
150  const IOobject& io,
151  const UPtrList<const labelList>& cellMaps,
152  const UPtrList<const labelList>& pointMaps,
153  const UPtrList<const hexRef8Data>& procDatas
154 )
155 {
156  const polyMesh& mesh = dynamic_cast<const polyMesh&>(io.db());
157 
158  // cellLevel
159 
160  if (procDatas[0].cellLevelPtr_.valid())
161  {
162  IOobject rio(io);
163  rio.rename(procDatas[0].cellLevelPtr_().name());
164 
165  cellLevelPtr_.reset(new labelIOList(rio, mesh.nCells()));
166  labelList& cellLevel = cellLevelPtr_();
167 
168  forAll(procDatas, procI)
169  {
170  const labelList& procCellLevel = procDatas[procI].cellLevelPtr_();
171  UIndirectList<label>(cellLevel, cellMaps[procI]) = procCellLevel;
172  }
173  }
174 
175 
176  // pointLevel
177 
178  if (procDatas[0].pointLevelPtr_.valid())
179  {
180  IOobject rio(io);
181  rio.rename(procDatas[0].pointLevelPtr_().name());
182 
183  pointLevelPtr_.reset(new labelIOList(rio, mesh.nPoints()));
184  labelList& pointLevel = pointLevelPtr_();
185 
186  forAll(procDatas, procI)
187  {
188  const labelList& procPointLevel = procDatas[procI].pointLevelPtr_();
189  UIndirectList<label>(pointLevel, pointMaps[procI]) = procPointLevel;
190  }
191  }
192 
193 
194  // level0Edge
195 
196  if (procDatas[0].level0EdgePtr_.valid())
197  {
198  IOobject rio(io);
199  rio.rename(procDatas[0].level0EdgePtr_().name());
200 
201  level0EdgePtr_.reset
202  (
204  (
205  rio,
206  procDatas[0].level0EdgePtr_()
207  )
208  );
209  }
210 
211 
212  // refinementHistory
213 
214  if (procDatas[0].refHistoryPtr_.valid())
215  {
216  IOobject rio(io);
217  rio.rename(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 * * * * * * * * * * * * * * * //
239 
241 {}
242 
243 
244 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
245 
247 {
248  const polyMesh& mesh = dynamic_cast<const polyMesh&>(io.db());
249 
250  bool hasCellLevel = returnReduce(cellLevelPtr_.valid(), orOp<bool>());
251  if (hasCellLevel && !cellLevelPtr_.valid())
252  {
253  IOobject rio(io);
254  rio.rename("cellLevel");
255  rio.readOpt() = IOobject::NO_READ;
256  cellLevelPtr_.reset(new labelIOList(rio, labelList(mesh.nCells(), 0)));
257  }
258 
259  bool hasPointLevel = returnReduce(pointLevelPtr_.valid(), orOp<bool>());
260  if (hasPointLevel && !pointLevelPtr_.valid())
261  {
262  IOobject rio(io);
263  rio.rename("pointLevel");
264  rio.readOpt() = IOobject::NO_READ;
265  pointLevelPtr_.reset
266  (
267  new labelIOList(rio, labelList(mesh.nPoints(), 0))
268  );
269  }
270 
271  bool hasLevel0Edge = returnReduce(level0EdgePtr_.valid(), orOp<bool>());
272  if (hasLevel0Edge)
273  {
274  // Get master length
275  scalar masterLen = (Pstream::master() ? level0EdgePtr_().value() : 0);
276  Pstream::scatter(masterLen);
277  if (!level0EdgePtr_.valid())
278  {
279  IOobject rio(io);
280  rio.rename("level0Edge");
281  rio.readOpt() = IOobject::NO_READ;
282  level0EdgePtr_.reset
283  (
285  (
286  rio,
287  dimensionedScalar("zero", dimLength, masterLen)
288  )
289  );
290  }
291  }
292 
293  bool hasHistory = returnReduce(refHistoryPtr_.valid(), orOp<bool>());
294  if (hasHistory && !refHistoryPtr_.valid())
295  {
296  IOobject rio(io);
297  rio.rename("refinementHistory");
298  rio.readOpt() = IOobject::NO_READ;
299  refHistoryPtr_.reset(new refinementHistory(rio, mesh.nCells(), true));
300  }
301 }
302 
303 
305 {
306  if (cellLevelPtr_.valid())
307  {
308  map.cellMap().distribute(cellLevelPtr_());
309  }
310  if (pointLevelPtr_.valid())
311  {
312  map.pointMap().distribute(pointLevelPtr_());
313  }
314 
315  // No need to distribute the level0Edge
316 
317  if (refHistoryPtr_.valid() && refHistoryPtr_().active())
318  {
319  refHistoryPtr_().distribute(map);
320  }
321 }
322 
323 
325 {
326  bool ok = true;
327  if (cellLevelPtr_.valid())
328  {
329  ok = ok && cellLevelPtr_().write();
330  }
331  if (pointLevelPtr_.valid())
332  {
333  ok = ok && pointLevelPtr_().write();
334  }
335  if (level0EdgePtr_.valid())
336  {
337  ok = ok && level0EdgePtr_().write();
338  }
339  if (refHistoryPtr_.valid())
340  {
341  ok = ok && refHistoryPtr_().write();
342  }
343  return ok;
344 }
345 
346 
347 // ************************************************************************* //
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::hexRef8Data::refHistoryPtr_
autoPtr< refinementHistory > refHistoryPtr_
Definition: hexRef4Data.H:70
Foam::labelList
List< label > labelList
A List of labels.
Definition: labelList.H:56
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:86
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
mapPolyMesh.H
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:147
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
polyMesh.H
Foam::labelIOList
IOList< label > labelIOList
Label container classes.
Definition: labelIOList.H:42
syncTools.H
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Foam::mapDistributePolyMesh::cellMap
const mapDistribute & cellMap() const
Cell distribute map.
Definition: mapDistributePolyMesh.H:214
Foam::IOobject::db
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:239
Foam::primitiveMesh::nCells
label nCells() const
Definition: primitiveMeshI.H:64
Foam::IOobject::headerOk
bool headerOk()
Read and check header info.
Definition: IOobject.C:439
Foam::refinementHistory
All refinement history. Used in unrefinement.
Definition: refinementHistory.H:95
mapDistributePolyMesh.H
Foam::hexRef8Data::level0EdgePtr_
autoPtr< uniformDimensionedScalarField > level0EdgePtr_
Definition: hexRef4Data.H:68
Foam::primitiveMesh::nPoints
label nPoints() const
Definition: primitiveMeshI.H:35
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:43
Foam::UniformDimensionedField< scalar >
Foam::IOobject::readOpt
readOption readOpt() const
Definition: IOobject.H:317
Foam::orOp
Definition: ops.H:178
Foam::hexRef8Data::pointLevelPtr_
autoPtr< labelIOList > pointLevelPtr_
Definition: hexRef4Data.H:66
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::Info
messageStream Info
Foam::uniformDimensionedScalarField
UniformDimensionedField< scalar > uniformDimensionedScalarField
Definition: uniformDimensionedFields.H:47
Foam::UPtrList
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: UPtrList.H:53
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
Foam::IOobject::name
const word & name() const
Return name.
Definition: IOobject.H:273
Foam::hexRef8Data::distribute
void distribute(const mapDistributePolyMesh &)
In-place distribute.
Definition: hexRef4Data.C:304
Foam::hexRef8Data
Various for reading/decomposing/reconstructing/distributing refinement data.
Definition: hexRef4Data.H:57
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:155
refinementHistory.H
IOobject.H
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
hexRef8Data.H
fvMesh.H
Foam::IOobject::rename
virtual void rename(const word &newName)
Rename.
Definition: IOobject.H:297
Foam::hexRef8Data::sync
void sync(const IOobject &io)
Parallel synchronise. This enforces valid objects on all processors.
Definition: hexRef4Data.C:246
Foam::UPstream::master
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:399
Foam::hexRef8Data::write
bool write() const
Write.
Definition: hexRef4Data.C:324
Foam::hexRef8Data::hexRef8Data
hexRef8Data(const hexRef8Data &)
Disallow default bitwise copy construct.
UList.H
Foam::IOdictionary::name
const word & name() const
Name function is needed to disambiguate those inherited.
Definition: IOdictionary.C:181
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::UPtrList::set
bool set(const label) const
Is element set.
Foam::dictionary::clone
autoPtr< dictionary > clone() const
Construct and return clone.
Definition: dictionary.C:215
Foam::IOList< label >
Foam::hexRef8Data::~hexRef8Data
~hexRef8Data()
Destructor.
Definition: hexRef4Data.C:240
Foam::UIndirectList
A List with indirect addressing.
Definition: fvMatrix.H:106
Foam::mapDistributePolyMesh::pointMap
const mapDistribute & pointMap() const
Point distribute map.
Definition: mapDistributePolyMesh.H:202
Foam::mapDistributePolyMesh
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Definition: mapDistributePolyMesh.H:57
Foam::IOobject::READ_IF_PRESENT
@ READ_IF_PRESENT
Definition: IOobject.H:110
Foam::hexRef8Data::cellLevelPtr_
autoPtr< labelIOList > cellLevelPtr_
Definition: hexRef4Data.H:64
Foam::data
Database for solution data, solver performance and other reduced data.
Definition: data.H:52
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Foam::UPtrList::size
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:31