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 | Copyright (C) 2015-2016 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 "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 
39 Foam::hexRef8Data::hexRef8Data(const IOobject& io)
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  bool haveFile = returnReduce(rio.headerOk(), orOp<bool>());
65  if (haveFile)
66  {
67  Info<< "Reading hexRef8 data : " << rio.name() << endl;
68  level0EdgePtr_.reset(new uniformDimensionedScalarField(rio));
69  }
70  }
71  {
72  IOobject rio(io);
73  rio.rename("refinementHistory");
74  bool haveFile = returnReduce(rio.headerOk(), orOp<bool>());
75  if (haveFile)
76  {
77  Info<< "Reading hexRef8 data : " << rio.name() << endl;
78  refHistoryPtr_.reset(new refinementHistory(rio));
79  }
80  }
81 }
82 
83 
84 Foam::hexRef8Data::hexRef8Data
85 (
86  const IOobject& io,
87  const hexRef8Data& data,
88  const labelList& cellMap,
89  const labelList& pointMap
90 )
91 {
92  if (data.cellLevelPtr_.valid())
93  {
94  IOobject rio(io);
95  rio.rename(data.cellLevelPtr_().name());
96 
97  cellLevelPtr_.reset
98  (
99  new labelIOList
100  (
101  rio,
102  UIndirectList<label>(data.cellLevelPtr_(), cellMap)()
103  )
104  );
105  }
106  if (data.pointLevelPtr_.valid())
107  {
108  IOobject rio(io);
109  rio.rename(data.pointLevelPtr_().name());
110 
111  pointLevelPtr_.reset
112  (
113  new labelIOList
114  (
115  rio,
116  UIndirectList<label>(data.pointLevelPtr_(), pointMap)()
117  )
118  );
119  }
120  if (data.level0EdgePtr_.valid())
121  {
122  IOobject rio(io);
123  rio.rename(data.level0EdgePtr_().name());
124 
125  level0EdgePtr_.reset
126  (
127  new uniformDimensionedScalarField(rio, data.level0EdgePtr_())
128  );
129  }
130  if (data.refHistoryPtr_.valid())
131  {
132  IOobject rio(io);
133  rio.rename(data.refHistoryPtr_().name());
134 
135  refHistoryPtr_ = data.refHistoryPtr_().clone(rio, cellMap);
136  }
137 }
138 
139 
140 Foam::hexRef8Data::hexRef8Data
141 (
142  const IOobject& io,
143  const UPtrList<const labelList>& cellMaps,
144  const UPtrList<const labelList>& pointMaps,
145  const UPtrList<const hexRef8Data>& procDatas
146 )
147 {
148  const polyMesh& mesh = dynamic_cast<const polyMesh&>(io.db());
149 
150  // cellLevel
151 
152  if (procDatas[0].cellLevelPtr_.valid())
153  {
154  IOobject rio(io);
155  rio.rename(procDatas[0].cellLevelPtr_().name());
156 
157  cellLevelPtr_.reset(new labelIOList(rio, mesh.nCells()));
158  labelList& cellLevel = cellLevelPtr_();
159 
160  forAll(procDatas, procI)
161  {
162  const labelList& procCellLevel = procDatas[procI].cellLevelPtr_();
163  UIndirectList<label>(cellLevel, cellMaps[procI]) = procCellLevel;
164  }
165  }
166 
167 
168  // pointLevel
169 
170  if (procDatas[0].pointLevelPtr_.valid())
171  {
172  IOobject rio(io);
173  rio.rename(procDatas[0].pointLevelPtr_().name());
174 
175  pointLevelPtr_.reset(new labelIOList(rio, mesh.nPoints()));
176  labelList& pointLevel = pointLevelPtr_();
177 
178  forAll(procDatas, procI)
179  {
180  const labelList& procPointLevel = procDatas[procI].pointLevelPtr_();
181  UIndirectList<label>(pointLevel, pointMaps[procI]) = procPointLevel;
182  }
183  }
184 
185 
186  // level0Edge
187 
188  if (procDatas[0].level0EdgePtr_.valid())
189  {
190  IOobject rio(io);
191  rio.rename(procDatas[0].level0EdgePtr_().name());
192 
193  level0EdgePtr_.reset
194  (
196  (
197  rio,
198  procDatas[0].level0EdgePtr_()
199  )
200  );
201  }
202 
203 
204  // refinementHistory
205 
206  if (procDatas[0].refHistoryPtr_.valid())
207  {
208  IOobject rio(io);
209  rio.rename(procDatas[0].refHistoryPtr_().name());
210 
211  UPtrList<const refinementHistory> procRefs(procDatas.size());
212  forAll(procDatas, i)
213  {
214  procRefs.set(i, &procDatas[i].refHistoryPtr_());
215  }
216 
217  refHistoryPtr_.reset
218  (
220  (
221  rio,
222  cellMaps,
223  procRefs
224  )
225  );
226  }
227 }
228 
229 
230 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
231 
233 {}
234 
235 
236 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
237 
239 {
240  const polyMesh& mesh = dynamic_cast<const polyMesh&>(io.db());
241 
242  bool hasCellLevel = returnReduce(cellLevelPtr_.valid(), orOp<bool>());
243  if (hasCellLevel && !cellLevelPtr_.valid())
244  {
245  IOobject rio(io);
246  rio.rename("cellLevel");
247  rio.readOpt() = IOobject::NO_READ;
248  cellLevelPtr_.reset(new labelIOList(rio, labelList(mesh.nCells(), 0)));
249  }
250 
251  bool hasPointLevel = returnReduce(pointLevelPtr_.valid(), orOp<bool>());
252  if (hasPointLevel && !pointLevelPtr_.valid())
253  {
254  IOobject rio(io);
255  rio.rename("pointLevel");
256  rio.readOpt() = IOobject::NO_READ;
257  pointLevelPtr_.reset
258  (
259  new labelIOList(rio, labelList(mesh.nPoints(), 0))
260  );
261  }
262 
263  bool hasLevel0Edge = returnReduce(level0EdgePtr_.valid(), orOp<bool>());
264  if (hasLevel0Edge)
265  {
266  // Get master length
267  scalar masterLen = level0EdgePtr_().value();
268  Pstream::scatter(masterLen);
269  if (!level0EdgePtr_.valid())
270  {
271  IOobject rio(io);
272  rio.rename("level0Edge");
273  rio.readOpt() = IOobject::NO_READ;
274  level0EdgePtr_.reset
275  (
277  (
278  rio,
279  dimensionedScalar("zero", dimLength, masterLen)
280  )
281  );
282  }
283  }
284 
285  bool hasHistory = returnReduce(refHistoryPtr_.valid(), orOp<bool>());
286  if (hasHistory && !refHistoryPtr_.valid())
287  {
288  IOobject rio(io);
289  rio.rename("refinementHistory");
290  rio.readOpt() = IOobject::NO_READ;
291  refHistoryPtr_.reset(new refinementHistory(rio, mesh.nCells(), true));
292  }
293 }
294 
295 
297 {
298  if (cellLevelPtr_.valid())
299  {
300  map.cellMap().distribute(cellLevelPtr_());
301  }
302  if (pointLevelPtr_.valid())
303  {
304  map.pointMap().distribute(pointLevelPtr_());
305  }
306 
307  // No need to distribute the level0Edge
308 
309  if (refHistoryPtr_.valid() && refHistoryPtr_().active())
310  {
311  refHistoryPtr_().distribute(map);
312  }
313 }
314 
315 
317 {
318  bool ok = true;
319  if (cellLevelPtr_.valid())
320  {
321  ok = ok && cellLevelPtr_().write();
322  }
323  if (pointLevelPtr_.valid())
324  {
325  ok = ok && pointLevelPtr_().write();
326  }
327  if (level0EdgePtr_.valid())
328  {
329  ok = ok && level0EdgePtr_().write();
330  }
331  if (refHistoryPtr_.valid())
332  {
333  ok = ok && refHistoryPtr_().write();
334  }
335  return ok;
336 }
337 
338 
339 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:221
void sync(const IOobject &io)
Parallel synchronise. This enforces valid objects on all processors.
Definition: hexRef8Data.C:238
bool write() const
Write.
Definition: hexRef8Data.C:316
virtual void rename(const word &newName)
Rename.
Definition: IOobject.H:284
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
UniformDimensionedField< scalar > uniformDimensionedScalarField
label nCells() const
dynamicFvMesh & mesh
List< label > labelList
A List of labels.
Definition: labelList.H:56
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: UPtrList.H:54
static void scatter(const List< commsStruct > &comms, T &Value, const int tag, const label comm)
Scatter data. Distribute without modification. Reverse of gather.
const mapDistribute & pointMap() const
Point distribute map.
void distribute(List< T > &fld, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Distribute data using default commsType.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
const mapDistribute & cellMap() const
Cell distribute map.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A List with indirect addressing.
Definition: fvMatrix.H:106
readOption readOpt() const
Definition: IOobject.H:304
label nPoints() const
messageStream Info
bool headerOk()
Read and check header info.
Definition: IOobject.C:400
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
void distribute(const mapDistributePolyMesh &)
In-place distribute.
Definition: hexRef8Data.C:296
Various for reading/decomposing/reconstructing/distributing refinement data.
Definition: hexRef8Data.H:57
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
All refinement history. Used in unrefinement.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
~hexRef8Data()
Destructor.
Definition: hexRef8Data.C:232
const word & name() const
Return name.
Definition: IOobject.H:260
IOList< label > labelIOList
Label container classes.
Definition: labelIOList.H:42
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29