Test-hexRef8.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) 2013-2015 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 Application
25  Test-hexRef8
26 
27 Description
28  Test app for refinement and unrefinement. Runs a few iterations refining
29  and unrefining.
30 
31 \*---------------------------------------------------------------------------*/
32 
33 #include "argList.H"
34 #include "Time.H"
35 #include "volFields.H"
36 #include "surfaceFields.H"
37 #include "pointFields.H"
38 #include "hexRef8.H"
39 #include "mapPolyMesh.H"
40 #include "polyTopoChange.H"
41 #include "Random.H"
44 #include "pointConstraints.H"
45 #include "fvcDiv.H"
46 
47 using namespace Foam;
48 
49 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50 
51 bool notEqual(const scalar s1, const scalar s2, const scalar tol)
52 {
53  return mag(s1-s2) > tol;
54 }
55 
56 
57 // Main program:
58 int main(int argc, char *argv[])
59 {
60  #include "addTimeOptions.H"
61  argList::validArgs.append("inflate (true|false)");
62  #include "setRootCase.H"
63  #include "createTime.H"
64  #include "createMesh.H"
65 
66 
68 
69  const Switch inflate(args.args()[1]);
70 
71  if (inflate)
72  {
73  Info<< "Splitting/deleting cells using inflation/deflation" << nl
74  << endl;
75  }
76  else
77  {
78  Info<< "Splitting/deleting cells, introducing points at new position"
79  << nl << endl;
80  }
81 
82 
83  Random rndGen(0);
84 
85 
86  // Force generation of V()
87  (void)mesh.V();
88 
89 
90  // Test mapping
91  // ------------
92 
93  // 1. uniform field stays uniform
95  (
96  IOobject
97  (
98  "one",
99  runTime.timeName(),
100  mesh,
103  ),
104  mesh,
105  dimensionedScalar("one", dimless, 1.0),
106  zeroGradientFvPatchScalarField::typeName
107  );
108  Info<< "Writing one field "
109  << one.name() << " in " << runTime.timeName() << endl;
110  one.write();
111 
112 
113  // 2. linear profile gets preserved
114  volScalarField ccX
115  (
116  IOobject
117  (
118  "ccX",
119  runTime.timeName(),
120  mesh,
123  ),
124  mesh.C().component(0)
125  );
126  Info<< "Writing x component of cell centres to "
127  << ccX.name()
128  << " in " << runTime.timeName() << endl;
129  ccX.write();
130 
131 
132  // Uniform surface field
133  surfaceScalarField surfaceOne
134  (
135  IOobject
136  (
137  "surfaceOne",
138  runTime.timeName(),
139  mesh,
142  ),
143  mesh,
144  dimensionedScalar("one", dimless, 1.0),
145  calculatedFvsPatchScalarField::typeName
146  );
147  Info<< "Writing surface one field "
148  << surfaceOne.name() << " in " << runTime.timeName() << endl;
149  surfaceOne.write();
150 
151 
152  // Uniform point field
153  pointScalarField pointX
154  (
155  IOobject
156  (
157  "pointX",
158  runTime.timeName(),
159  mesh,
162  ),
164  dimensionedScalar("one", dimless, 1.0),
165  calculatedPointPatchScalarField::typeName
166  );
167  pointX.internalField() = mesh.points().component(0);
168  pointX.correctBoundaryConditions();
169  Info<< "Writing x-component field "
170  << pointX.name() << " in " << runTime.timeName() << endl;
171  pointX.write();
172 
173 
174 
175  // Force allocation of V. Important for any mesh changes since otherwise
176  // old time volumes are not stored
177  const scalar totalVol = gSum(mesh.V());
178 
179 
180  // Construct refiner. Read initial cell and point levels.
182 
183 
184  while (runTime.loop())
185  {
186  Info<< "Time = " << runTime.timeName() << nl << endl;
187 
188  if (mesh.globalData().nTotalCells() == 0)
189  {
190  break;
191  }
192 
193 
194  mesh.moving(false);
195  mesh.topoChanging(false);
196 
197 
198  label action = rndGen.integer(0, 5);
199 
200 
201  if (action == 0)
202  {
203  Info<< nl << "-- moving only" << endl;
205  }
206  else if (action == 1 || action == 2)
207  {
208  // Mesh changing engine.
209  polyTopoChange meshMod(mesh);
210 
211  if (action == 1)
212  {
213  // Refine
214  label nRefine = mesh.nCells()/20;
215  DynamicList<label> refineCandidates(nRefine);
216 
217  for (label i=0; i<nRefine; i++)
218  {
219  refineCandidates.append(rndGen.integer(0, mesh.nCells()-1));
220  }
221 
222  labelList cellsToRefine
223  (
224  meshCutter.consistentRefinement
225  (
226  refineCandidates,
227  true // buffer layer
228  )
229  );
230  Info<< nl << "-- selected "
231  << returnReduce(cellsToRefine.size(), sumOp<label>())
232  << " cells out of " << mesh.globalData().nTotalCells()
233  << " for refinement" << endl;
234 
235  // Play refinement commands into mesh changer.
236  meshCutter.setRefinement(cellsToRefine, meshMod);
237  }
238  else
239  {
240  // Unrefine
241  labelList allSplitPoints(meshCutter.getSplitPoints());
242 
243  label nUnrefine = allSplitPoints.size()/20;
244  labelHashSet candidates(2*nUnrefine);
245 
246  for (label i=0; i<nUnrefine; i++)
247  {
248  label index = rndGen.integer(0, allSplitPoints.size()-1);
249  candidates.insert(allSplitPoints[index]);
250  }
251 
252  labelList splitPoints = meshCutter.consistentUnrefinement
253  (
254  candidates.toc(),
255  false
256  );
257  Info<< nl << "-- selected "
258  << returnReduce(splitPoints.size(), sumOp<label>())
259  << " points out of "
260  << returnReduce(allSplitPoints.size(), sumOp<label>())
261  << " for unrefinement" << endl;
262 
263  // Play refinement commands into mesh changer.
264  meshCutter.setUnrefinement(splitPoints, meshMod);
265  }
266 
267 
268 
269 
270  // Create mesh, return map from old to new mesh.
271  Info<< nl << "-- actually changing mesh" << endl;
272  autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh, inflate);
273 
274  // Update fields
275  Info<< nl << "-- mapping mesh data" << endl;
276  mesh.updateMesh(map);
277 
278  // Inflate mesh
279  if (map().hasMotionPoints())
280  {
281  Info<< nl << "-- moving mesh" << endl;
282  mesh.movePoints(map().preMotionPoints());
283  }
284 
285  // Update numbering of cells/vertices.
286  Info<< nl << "-- mapping hexRef8 data" << endl;
287  meshCutter.updateMesh(map);
288  }
289 
290 
291  Info<< nl<< "-- Mesh : moving:" << mesh.moving()
292  << " topoChanging:" << mesh.topoChanging()
293  << " changing:" << mesh.changing()
294  << endl;
295 
296 
297 
298  Info<< "Writing fields" << nl << endl;
299  runTime.write();
300 
301 
302 
303  // Check mesh volume conservation
304  if (mesh.moving())
305  {
306  #include "volContinuity.H"
307  }
308  else
309  {
310  if (mesh.V().size() != mesh.nCells())
311  {
313  << "Volume not mapped. V:" << mesh.V().size()
314  << " nCells:" << mesh.nCells()
315  << exit(FatalError);
316  }
317 
318  const scalar newVol = gSum(mesh.V());
319  Info<< "Initial volume = " << totalVol
320  << " New volume = " << newVol
321  << endl;
322 
323  if (mag(newVol-totalVol)/totalVol > 1e-10)
324  {
326  << "Volume loss: old volume:" << totalVol
327  << " new volume:" << newVol
328  << exit(FatalError);
329  }
330  else
331  {
332  Info<< "Volume check OK" << nl << endl;
333  }
334  }
335 
336 
337  // Check constant profile
338  {
339  const scalar max = gMax(one);
340  const scalar min = gMin(one);
341 
342  Info<< "Uniform one field min = " << min
343  << " max = " << max << endl;
344 
345  if (notEqual(max, 1.0, 1e-10) || notEqual(min, 1.0, 1e-10))
346  {
348  << "Uniform volVectorField not preserved."
349  << " Min and max should both be 1.0. min:" << min
350  << " max:" << max
351  << exit(FatalError);
352  }
353  else
354  {
355  Info<< "Uniform field mapping check OK" << nl << endl;
356  }
357  }
358 
359  // Check linear profile
360  {
361  const scalarField diff = ccX-mesh.C().component(0);
362 
363  const scalar max = gMax(diff);
364  const scalar min = gMin(diff);
365 
366  Info<< "Linear profile field min = " << min
367  << " max = " << max << endl;
368 
369  if (notEqual(max, 0.0, 1e-10) || notEqual(min, 0.0, 1e-10))
370  {
371  Info<< "Linear profile not preserved."
372  << " Min and max should both be 0.0. min:" << min
373  << " max:" << max << nl << endl;
374  }
375  else
376  {
377  Info<< "Linear profile mapping check OK" << nl << endl;
378  }
379  }
380 
381  // Check face field mapping
382  if (surfaceOne.size())
383  {
384  const scalar max = gMax(surfaceOne.internalField());
385  const scalar min = gMin(surfaceOne.internalField());
386 
387  Info<< "Uniform surface field min = " << min
388  << " max = " << max << endl;
389 
390  if (notEqual(max, 1.0, 1e-10) || notEqual(min, 1.0, 1e-10))
391  {
393  << "Uniform surfaceScalarField not preserved."
394  << " Min and max should both be 1.0. min:" << min
395  << " max:" << max
396  << exit(FatalError);
397  }
398  else
399  {
400  Info<< "Uniform surfaceScalarField mapping check OK" << nl
401  << endl;
402  }
403  }
404 
405 
406  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
407  << " ClockTime = " << runTime.elapsedClockTime() << " s"
408  << nl << endl;
409  }
410 
411  Info<< "pc:" << pc.patchPatchPointConstraintPoints().size() << endl;
412 
413 
414  Info<< "End\n" << endl;
415 
416  return 0;
417 }
418 
419 
420 // ************************************************************************* //
Foam::argList::validArgs
static SLList< string > validArgs
A list of valid (mandatory) arguments.
Definition: argList.H:143
volFields.H
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
Foam::polyMesh::changing
bool changing() const
Is mesh changing (topology changing and/or moving)
Definition: polyMesh.H:521
Foam::Random
Simple random number generator.
Definition: Random.H:49
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:979
Foam::polyMesh::topoChanging
bool topoChanging() const
Is mesh topology changing.
Definition: polyMesh.H:507
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::GeometricField::component
tmp< GeometricField< cmptType, PatchField, GeoMesh > > component(const direction) const
Return a component of the field.
Foam::Switch
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:60
main
int main(int argc, char *argv[])
Definition: Test-hexRef8.C:55
Foam::meshCutter
Cuts (splits) cells.
Definition: meshCutter.H:134
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
Foam::argList::args
const stringList & args() const
Return arguments.
Definition: argListI.H:66
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::HashTable::toc
List< Key > toc() const
Return the table of contents.
Definition: HashTable.C:201
Foam::IOobject::AUTO_WRITE
@ AUTO_WRITE
Definition: IOobject.H:117
Foam::polyTopoChange::changeMesh
autoPtr< mapPolyMesh > changeMesh(polyMesh &mesh, const bool inflate, const bool syncParallel=true, const bool orderCells=false, const bool orderPoints=false)
Inplace changes mesh without change of patches.
Definition: polyTopoChange.C:3039
Foam::DynamicList< label >
calculatedPointPatchFields.H
Foam::polyMesh::moving
bool moving() const
Is mesh moving.
Definition: polyMesh.H:493
mapPolyMesh.H
polyTopoChange.H
Foam::polyTopoChange
Direct mesh changes based on v1.3 polyTopoChange syntax.
Definition: polyTopoChange.H:97
Foam::globalMeshData::nTotalCells
label nTotalCells() const
Return total number of cells in decomposed mesh.
Definition: globalMeshData.H:394
fvcDiv.H
Calculate the divergence of the given field.
Foam::one
A class representing the concept of 1 (scalar(1.0)) used to avoid unnecessary manipulations for objec...
Definition: one.H:45
Foam::fvMesh::movePoints
virtual tmp< scalarField > movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
Definition: fvMesh.C:726
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
surfaceFields.H
Foam::surfaceFields.
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:564
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::meshCutter::updateMesh
void updateMesh(const mapPolyMesh &)
Force recalculation of locally stored data on topological change.
Definition: meshCutter.C:972
Foam::HashSet< label, Hash< label > >
addTimeOptions.H
Foam::MeshObject< pointMesh, UpdateableMeshObject, pointConstraints >::New
static const pointConstraints & New(const pointMesh &mesh)
Definition: MeshObject.C:44
Foam::primitiveMesh::nCells
label nCells() const
Definition: primitiveMeshI.H:64
Foam::diff
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
Definition: triad.C:407
Foam::fvMesh::V
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
Definition: fvMeshGeometry.C:199
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::Ostream::write
virtual Ostream & write(const token &)=0
Write next token to stream.
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
hexRef8.H
pointConstraints.H
argList.H
Foam::fvMesh::updateMesh
virtual void updateMesh(const mapPolyMesh &mpm)
Update mesh corresponding to the given map.
Definition: fvMesh.C:802
Foam::GeometricField::internalField
InternalField & internalField()
Return internal field.
Definition: GeometricField.C:724
Foam::FatalError
error FatalError
Foam::pointConstraints::patchPatchPointConstraintPoints
const labelList & patchPatchPointConstraintPoints() const
Mesh points on which to apply special constraints.
Definition: pointConstraints.H:111
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::Field::component
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Definition: Field.C:644
Foam::e
const double e
Elementary charge.
Definition: doubleFloat.H:94
Random.H
Foam::DynamicList::append
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::GeometricField::correctBoundaryConditions
void correctBoundaryConditions()
Correct boundary field.
Definition: GeometricField.C:885
Foam::pointConstraints
Application of (multi-)patch point contraints.
Definition: pointConstraints.H:62
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::autoPtr
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
setRootCase.H
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::fvMesh::C
const volVectorField & C() const
Return cell centres as volVectorField.
Definition: fvMeshGeometry.C:369
Foam::sumOp
Definition: ops.H:162
Foam::hexRef8
Refinement of (split) hexes using polyTopoChange.
Definition: hexRef4.H:64
volContinuity.H
Foam::meshCutter::setRefinement
void setRefinement(const cellCuts &cuts, polyTopoChange &meshMod)
Do actual cutting with cut description. Inserts mesh changes.
Definition: meshCutter.C:524
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
createMesh.H
Foam::HashSet::insert
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:116
createTime.H
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::gMin
Type gMin(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:563
Foam::polyMesh::globalData
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1143
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
args
Foam::argList args(argc, argv)
Foam::notEqual
bool notEqual(const Scalar s1, const Scalar s2)
Definition: Scalar.H:155
zeroGradientFvPatchFields.H
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
rndGen
cachedRandom rndGen(label(0), -1)
pointFields.H
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:562