foamyHexMeshBackgroundMesh.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) 2012-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  foamyHexMeshBackGroundMesh
26 
27 Description
28  Writes out background mesh as constructed by foamyHexMesh and constructs
29  distanceSurface.
30 
31 \*---------------------------------------------------------------------------*/
32 
33 #include "PatchTools.H"
34 #include "argList.H"
35 #include "Time.H"
36 #include "triSurface.H"
37 #include "searchableSurfaces.H"
38 #include "conformationSurfaces.H"
39 #include "cellShapeControl.H"
41 #include "cellShape.H"
42 #include "cellModeller.H"
43 #include "DynamicField.H"
44 #include "isoSurfaceCell.H"
45 #include "vtkSurfaceWriter.H"
46 #include "syncTools.H"
47 #include "decompositionModel.H"
48 
49 using namespace Foam;
50 
51 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52 
53 // Tolerance (as fraction of the bounding box). Needs to be fairly lax since
54 // usually meshes get written with limited precision (6 digits)
55 static const scalar defaultMergeTol = 1e-6;
56 
57 // Get merging distance when matching face centres
58 scalar getMergeDistance
59 (
60  const argList& args,
61  const Time& runTime,
62  const boundBox& bb
63 )
64 {
65  scalar mergeTol = defaultMergeTol;
66  args.optionReadIfPresent("mergeTol", mergeTol);
67 
68  scalar writeTol =
69  Foam::pow(scalar(10.0), -scalar(IOstream::defaultPrecision()));
70 
71  Info<< "Merge tolerance : " << mergeTol << nl
72  << "Write tolerance : " << writeTol << endl;
73 
74  if (runTime.writeFormat() == IOstream::ASCII && mergeTol < writeTol)
75  {
77  << "Your current settings specify ASCII writing with "
78  << IOstream::defaultPrecision() << " digits precision." << endl
79  << "Your merging tolerance (" << mergeTol << ") is finer than this."
80  << endl
81  << "Please change your writeFormat to binary"
82  << " or increase the writePrecision" << endl
83  << "or adjust the merge tolerance (-mergeTol)."
84  << exit(FatalError);
85  }
86 
87  scalar mergeDist = mergeTol * bb.mag();
88 
89  Info<< "Overall meshes bounding box : " << bb << nl
90  << "Relative tolerance : " << mergeTol << nl
91  << "Absolute matching distance : " << mergeDist << nl
92  << endl;
93 
94  return mergeDist;
95 }
96 
97 
98 void printMeshData(const polyMesh& mesh)
99 {
100  // Collect all data on master
101 
102  globalIndex globalCells(mesh.nCells());
103  labelListList patchNeiProcNo(Pstream::nProcs());
104  labelListList patchSize(Pstream::nProcs());
105  const labelList& pPatches = mesh.globalData().processorPatches();
106  patchNeiProcNo[Pstream::myProcNo()].setSize(pPatches.size());
107  patchSize[Pstream::myProcNo()].setSize(pPatches.size());
108  forAll(pPatches, i)
109  {
110  const processorPolyPatch& ppp = refCast<const processorPolyPatch>
111  (
112  mesh.boundaryMesh()[pPatches[i]]
113  );
114  patchNeiProcNo[Pstream::myProcNo()][i] = ppp.neighbProcNo();
115  patchSize[Pstream::myProcNo()][i] = ppp.size();
116  }
117  Pstream::gatherList(patchNeiProcNo);
118  Pstream::gatherList(patchSize);
119 
120 
121  // Print stats
122 
123  globalIndex globalBoundaryFaces(mesh.nFaces()-mesh.nInternalFaces());
124 
125  label maxProcCells = 0;
126  label totProcFaces = 0;
127  label maxProcPatches = 0;
128  label totProcPatches = 0;
129  label maxProcFaces = 0;
130 
131  for (label procI = 0; procI < Pstream::nProcs(); procI++)
132  {
133  Info<< endl
134  << "Processor " << procI << nl
135  << " Number of cells = " << globalCells.localSize(procI)
136  << endl;
137 
138  label nProcFaces = 0;
139 
140  const labelList& nei = patchNeiProcNo[procI];
141 
142  forAll(patchNeiProcNo[procI], i)
143  {
144  Info<< " Number of faces shared with processor "
145  << patchNeiProcNo[procI][i] << " = " << patchSize[procI][i]
146  << endl;
147 
148  nProcFaces += patchSize[procI][i];
149  }
150 
151  Info<< " Number of processor patches = " << nei.size() << nl
152  << " Number of processor faces = " << nProcFaces << nl
153  << " Number of boundary faces = "
154  << globalBoundaryFaces.localSize(procI) << endl;
155 
156  maxProcCells = max(maxProcCells, globalCells.localSize(procI));
157  totProcFaces += nProcFaces;
158  totProcPatches += nei.size();
159  maxProcPatches = max(maxProcPatches, nei.size());
160  maxProcFaces = max(maxProcFaces, nProcFaces);
161  }
162 
163  // Stats
164 
165  scalar avgProcCells = scalar(globalCells.size())/Pstream::nProcs();
166  scalar avgProcPatches = scalar(totProcPatches)/Pstream::nProcs();
167  scalar avgProcFaces = scalar(totProcFaces)/Pstream::nProcs();
168 
169  // In case of all faces on one processor. Just to avoid division by 0.
170  if (totProcPatches == 0)
171  {
172  avgProcPatches = 1;
173  }
174  if (totProcFaces == 0)
175  {
176  avgProcFaces = 1;
177  }
178 
179  Info<< nl
180  << "Number of processor faces = " << totProcFaces/2 << nl
181  << "Max number of cells = " << maxProcCells
182  << " (" << 100.0*(maxProcCells-avgProcCells)/avgProcCells
183  << "% above average " << avgProcCells << ")" << nl
184  << "Max number of processor patches = " << maxProcPatches
185  << " (" << 100.0*(maxProcPatches-avgProcPatches)/avgProcPatches
186  << "% above average " << avgProcPatches << ")" << nl
187  << "Max number of faces between processors = " << maxProcFaces
188  << " (" << 100.0*(maxProcFaces-avgProcFaces)/avgProcFaces
189  << "% above average " << avgProcFaces << ")" << nl
190  << endl;
191 }
192 
193 
194 // Return cellID
195 label cellLabel
196 (
197  const Vector<label>& nCells,
198  const label i,
199  const label j,
200  const label k
201 )
202 {
203  return i*nCells[1]*nCells[2]+j*nCells[2]+k;
204 }
205 
206 label vtxLabel
207 (
208  const Vector<label>& nCells,
209  const label i,
210  const label j,
211  const label k
212 )
213 {
215  (
216  nCells[0]+1,
217  nCells[1]+1,
218  nCells[2]+1
219  );
220  return i*nPoints[1]*nPoints[2]+j*nPoints[2]+k;
221 }
222 
223 
224 autoPtr<polyMesh> generateHexMesh
225 (
226  const IOobject& io,
227  const point& origin,
228  const vector& cellSize,
229  const Vector<label>& nCells
230 )
231 {
233  if (nCells[0]+nCells[1]+nCells[2] > 0)
234  {
235  points.setSize((nCells[0]+1)*(nCells[1]+1)*(nCells[2]+1));
236 
237  // Generate points
238  for (label i = 0; i <= nCells[0]; i++)
239  {
240  for (label j = 0; j <= nCells[1]; j++)
241  {
242  for (label k = 0; k <= nCells[2]; k++)
243  {
244  point pt = origin;
245  pt.x() += i*cellSize[0];
246  pt.y() += j*cellSize[1];
247  pt.z() += k*cellSize[2];
248  points[vtxLabel(nCells, i, j, k)] = pt;
249  }
250  }
251  }
252  }
253 
254 
255  const cellModel& hex = *(cellModeller::lookup("hex"));
256  cellShapeList cellShapes(nCells[0]*nCells[1]*nCells[2]);
257 
258  labelList hexPoints(8);
259  label cellI = 0;
260  for (label i = 0; i < nCells[0]; i++)
261  {
262  for (label j = 0; j < nCells[1]; j++)
263  {
264  for (label k = 0; k < nCells[2]; k++)
265  {
266  hexPoints[0] = vtxLabel(nCells, i, j, k);
267  hexPoints[1] = vtxLabel(nCells, i+1, j, k);
268  hexPoints[2] = vtxLabel(nCells, i+1, j+1, k);
269  hexPoints[3] = vtxLabel(nCells, i, j+1, k);
270  hexPoints[4] = vtxLabel(nCells, i, j, k+1);
271  hexPoints[5] = vtxLabel(nCells, i+1, j, k+1);
272  hexPoints[6] = vtxLabel(nCells, i+1, j+1, k+1);
273  hexPoints[7] = vtxLabel(nCells, i, j+1, k+1);
274  cellShapes[cellI++] = cellShape(hex, hexPoints);
275  }
276  }
277  }
278 
280  wordList patchNames(0);
281  wordList patchTypes(0);
282  word defaultFacesName = "defaultFaces";
283  word defaultFacesType = polyPatch::typeName;
284  wordList patchPhysicalTypes(0);
285 
286  return autoPtr<polyMesh>
287  (
288  new polyMesh
289  (
290  io,
291  xferMoveTo<pointField>(points),
292  cellShapes,
293  boundary,
294  patchNames,
295  patchTypes,
298  patchPhysicalTypes
299  )
300  );
301 }
302 
303 
304 // Determine for every point a signed distance to the nearest surface
305 // (outside is positive)
306 tmp<scalarField> signedDistance
307 (
308  const scalarField& distSqr,
309  const pointField& points,
310  const searchableSurfaces& geometry,
311  const labelList& surfaces
312 )
313 {
314  tmp<scalarField> tfld(new scalarField(points.size(), Foam::sqr(GREAT)));
315  scalarField& fld = tfld();
316 
317  // Find nearest
318  List<pointIndexHit> nearest;
319  labelList nearestSurfaces;
321  (
322  geometry,
323  surfaces,
324  points,
325  scalarField(points.size(), Foam::sqr(GREAT)),//distSqr
326  nearestSurfaces,
327  nearest
328  );
329 
330  // Determine sign of nearest. Sort by surface to do this.
331  DynamicField<point> surfPoints(points.size());
332  DynamicList<label> surfIndices(points.size());
333 
334  forAll(surfaces, surfI)
335  {
336  // Extract points on this surface
337  surfPoints.clear();
338  surfIndices.clear();
339  forAll(nearestSurfaces, i)
340  {
341  if (nearestSurfaces[i] == surfI)
342  {
343  surfPoints.append(points[i]);
344  surfIndices.append(i);
345  }
346  }
347 
348  // Calculate sideness of these surface points
349  label geomI = surfaces[surfI];
350  List<volumeType> volType;
351  geometry[geomI].getVolumeType(surfPoints, volType);
352 
353  // Push back to original
354  forAll(volType, i)
355  {
356  label pointI = surfIndices[i];
357  scalar dist = mag(points[pointI] - nearest[pointI].hitPoint());
358 
359  volumeType vT = volType[i];
360 
361  if (vT == volumeType::OUTSIDE)
362  {
363  fld[pointI] = dist;
364  }
365  else if (vT == volumeType::INSIDE)
366  {
367  fld[i] = -dist;
368  }
369  else
370  {
372  << "getVolumeType failure, neither INSIDE or OUTSIDE"
373  << exit(FatalError);
374  }
375  }
376  }
377 
378  return tfld;
379 }
380 
381 
382 
383 // Main program:
384 
385 int main(int argc, char *argv[])
386 {
388  (
389  "Generate foamyHexMesh-consistent representation of surfaces"
390  );
392  (
393  "writeMesh",
394  "write the resulting mesh and distance fields"
395  );
397  (
398  "mergeTol",
399  "scalar",
400  "specify the merge distance relative to the bounding box size "
401  "(default 1e-6)"
402  );
403 
404  #include "setRootCase.H"
405  #include "createTime.H"
406  runTime.functionObjects().off();
407 
408  const bool writeMesh = args.optionFound("writeMesh");
409 
410  if (writeMesh)
411  {
412  Info<< "Writing resulting mesh and cellDistance, pointDistance fields."
413  << nl << endl;
414  }
415 
416 
417  IOdictionary foamyHexMeshDict
418  (
419  IOobject
420  (
421  "foamyHexMeshDict",
422  runTime.system(),
423  runTime,
426  )
427  );
428 
429  // Define/load all geometry
430  searchableSurfaces allGeometry
431  (
432  IOobject
433  (
434  "cvSearchableSurfaces",
435  runTime.constant(),
436  "triSurface",
437  runTime,
440  ),
441  foamyHexMeshDict.subDict("geometry"),
442  foamyHexMeshDict.lookupOrDefault("singleRegionName", true)
443  );
444 
446 
447  conformationSurfaces geometryToConformTo
448  (
449  runTime,
450  rndGen,
451  allGeometry,
452  foamyHexMeshDict.subDict("surfaceConformation")
453  );
454 
455  cellShapeControl cellShapeControls
456  (
457  runTime,
458  foamyHexMeshDict.subDict("motionControl"),
459  allGeometry,
460  geometryToConformTo
461  );
462 
463 
464  // Generate starting block mesh
465  vector cellSize;
466  {
467  const treeBoundBox& bb = geometryToConformTo.globalBounds();
468 
469  // Determine the number of cells in each direction.
470  const vector span = bb.span();
471  vector nScalarCells = span/cellShapeControls.defaultCellSize();
472 
473  // Calculate initial cell size to be a little bit smaller than the
474  // defaultCellSize to avoid initial refinement triggering.
476  (
477  label(nScalarCells.x())+2,
478  label(nScalarCells.y())+2,
479  label(nScalarCells.z())+2
480  );
481  cellSize = vector
482  (
483  span[0]/nCells[0],
484  span[1]/nCells[1],
485  span[2]/nCells[2]
486  );
487 
488  Info<< "Generating initial hex mesh with" << nl
489  << " bounding box : " << bb << nl
490  << " nCells : " << nCells << nl
491  << " cellSize : " << cellSize << nl
492  << endl;
493 
495  (
496  generateHexMesh
497  (
498  IOobject
499  (
501  runTime.constant(),
502  runTime
503  ),
504  bb.min(),
505  cellSize,
506  (
508  ? nCells
509  : Vector<label>(0, 0, 0)
510  )
511  )
512  );
513  Info<< "Writing initial hex mesh to " << meshPtr().instance() << nl
514  << endl;
515  meshPtr().write();
516  }
517 
518  // Distribute the initial mesh
519  if (Pstream::parRun())
520  {
521  #include "createMesh.H"
522  Info<< "Loaded mesh:" << endl;
524 
525  // Allow override of decomposeParDict location
526  fileName decompDictFile;
527  if (args.optionReadIfPresent("decomposeParDict", decompDictFile))
528  {
529  if (isDir(decompDictFile))
530  {
531  decompDictFile = decompDictFile / "decomposeParDict";
532  }
533  }
534 
536  (
537  mesh,
538  decompDictFile
540 
541  // Global matching tolerance
542  const scalar tolDim = getMergeDistance
543  (
544  args,
545  runTime,
546  mesh.bounds()
547  );
548 
549  // Mesh distribution engine
550  fvMeshDistribute distributor(mesh, tolDim);
551 
552  Info<< "Wanted distribution:"
553  << distributor.countCells(decomp) << nl << endl;
554 
555  // Do actual sending/receiving of mesh
556  autoPtr<mapDistributePolyMesh> map = distributor.distribute(decomp);
557 
558  // Print some statistics
559  //Info<< "After distribution:" << endl;
560  //printMeshData(mesh);
561 
562  mesh.setInstance(runTime.constant());
563  Info<< "Writing redistributed mesh" << nl << endl;
564  mesh.write();
565  }
566 
567 
568  Info<< "Refining backgroud mesh according to cell size specification" << nl
569  << endl;
570 
571  const dictionary& backgroundMeshDict =
572  foamyHexMeshDict.subDict("backgroundMeshDecomposition");
573 
574  backgroundMeshDecomposition backgroundMesh
575  (
576  runTime,
577  rndGen,
578  geometryToConformTo,
579  backgroundMeshDict
580  );
581 
582  if (writeMesh)
583  {
584  runTime++;
585  Info<< "Writing mesh to " << runTime.timeName() << endl;
586  backgroundMesh.mesh().write();
587  }
588 
589  const scalar tolDim = getMergeDistance
590  (
591  args,
592  runTime,
593  backgroundMesh.mesh().bounds()
594  );
595 
596 
597  faceList isoFaces;
598  pointField isoPoints;
599 
600  {
601  // Apply a distanceSurface to it.
602  const fvMesh& fvm = backgroundMesh.mesh();
603 
604  volScalarField cellDistance
605  (
606  IOobject
607  (
608  "cellDistance",
609  fvm.time().timeName(),
610  fvm.time(),
613  false
614  ),
615  fvm,
616  dimensionedScalar("zero", dimLength, 0)
617  );
618 
619  const searchableSurfaces& geometry = geometryToConformTo.geometry();
620  const labelList& surfaces = geometryToConformTo.surfaces();
621 
622 
623  // Get maximum search size per cell
624  scalarField distSqr(cellDistance.size());
625 
626  const labelList& cellLevel = backgroundMesh.cellLevel();
627  forAll(cellLevel, cellI)
628  {
629  // The largest edge of the cell will always be less than the
630  // span of the bounding box of the cell.
631  distSqr[cellI] = magSqr(cellSize)/pow(2, cellLevel[cellI]);
632  }
633 
634  {
635  // Internal field
636  cellDistance.internalField() = signedDistance
637  (
638  distSqr,
639  fvm.C(),
640  geometry,
641  surfaces
642  );
643  // Patch fields
644  forAll(fvm.C().boundaryField(), patchI)
645  {
646  const pointField& cc = fvm.C().boundaryField()[patchI];
647  fvPatchScalarField& fld = cellDistance.boundaryField()[patchI];
648  scalarField patchDistSqr
649  (
650  fld.patch().patchInternalField(distSqr)
651  );
652  fld = signedDistance(patchDistSqr, cc, geometry, surfaces);
653  }
654 
655  // On processor patches the fvm.C() will already be the cell centre
656  // on the opposite side so no need to swap cellDistance.
657 
658  if (writeMesh)
659  {
660  cellDistance.write();
661  }
662  }
663 
664 
665  // Distance to points
666  pointScalarField pointDistance
667  (
668  IOobject
669  (
670  "pointDistance",
671  fvm.time().timeName(),
672  fvm.time(),
675  false
676  ),
677  pointMesh::New(fvm),
678  dimensionedScalar("zero", dimLength, 0)
679  );
680  {
681  scalarField pointDistSqr(fvm.nPoints(), -sqr(GREAT));
682  for (label faceI = 0; faceI < fvm.nInternalFaces(); faceI++)
683  {
684  label own = fvm.faceOwner()[faceI];
685  label ownDistSqr = distSqr[own];
686 
687  const face& f = fvm.faces()[faceI];
688  forAll(f, fp)
689  {
690  pointDistSqr[f[fp]] = max(pointDistSqr[f[fp]], ownDistSqr);
691  }
692  }
694  (
695  fvm,
696  pointDistSqr,
697  maxEqOp<scalar>(),
698  -sqr(GREAT) // null value
699  );
700 
701  pointDistance.internalField() = signedDistance
702  (
703  pointDistSqr,
704  fvm.points(),
705  geometry,
706  surfaces
707  );
708 
709  if (writeMesh)
710  {
711  pointDistance.write();
712  }
713  }
714 
715  isoSurfaceCell iso
716  (
717  fvm,
718  cellDistance,
719  pointDistance,
720  0, //distance,
721  false //regularise
722  );
723 
724  isoFaces.setSize(iso.size());
725  forAll(isoFaces, i)
726  {
727  isoFaces[i] = iso[i].triFaceFace();
728  }
729  isoPoints = iso.points();
730  }
731 
732 
733  pointField mergedPoints;
734  faceList mergedFaces;
735  labelList pointMergeMap;
737  (
738  tolDim,
740  (
741  SubList<face>(isoFaces, isoFaces.size()),
742  isoPoints
743  ),
744  mergedPoints,
745  mergedFaces,
746  pointMergeMap
747  );
748 
750  writer.write
751  (
752  runTime.path(),
753  "iso",
754  mergedPoints,
755  mergedFaces
756  );
757 
758  Info<< "End\n" << endl;
759 
760  return 0;
761 }
762 
763 
764 // ************************************************************************* //
cellShapes
const cellShapeList & cellShapes
Definition: getFieldScalar.H:35
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
Foam::Random
Simple random number generator.
Definition: Random.H:49
Foam::volumeType::OUTSIDE
@ OUTSIDE
Definition: volumeType.H:64
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:979
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:48
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::volumeType::INSIDE
@ INSIDE
Definition: volumeType.H:63
Foam::boundBox::mag
scalar mag() const
The magnitude of the bounding box span.
Definition: boundBoxI.H:90
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
cellShapeControl.H
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::fileName
A class for handling file names.
Definition: fileName.H:69
Foam::argList::addOption
static void addOption(const word &opt, const string &param="", const string &usage="")
Add to an option to validOptions with usage information.
Definition: argList.C:108
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
Foam::argList::addNote
static void addNote(const string &)
Add extra notes for the usage information.
Definition: argList.C:139
Foam::polyMesh::defaultRegion
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:306
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::DynamicList< label >
Foam::globalMeshData::processorPatches
const labelList & processorPatches() const
Return list of processor patch labels.
Definition: globalMeshData.H:404
isoSurfaceCell.H
Foam::Time::writeFormat
IOstream::streamFormat writeFormat() const
Default write format.
Definition: Time.H:303
Foam::SubList
A List obtained as a section of another List.
Definition: SubList.H:53
Foam::UPstream::nProcs
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:387
Foam::treeBoundBox
Standard boundBox + extra functionality for use in octree.
Definition: treeBoundBox.H:75
boundary
faceListList boundary(nPatches)
Foam::UPstream::parRun
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:377
printMeshData
void printMeshData(const polyMesh &mesh)
Definition: redistributePar.C:128
Foam::maxEqOp
Definition: ops.H:77
Foam::writer::write
virtual void write(const coordSet &, const wordList &, const List< const Field< Type > * > &, Ostream &) const =0
General entry point for writing.
Foam::argList::addBoolOption
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:98
Foam::Time::functionObjects
const functionObjectList & functionObjects() const
Return the list of function objects.
Definition: Time.H:435
Foam::GeometricField::boundaryField
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Definition: GeometricField.C:735
Foam::processorPolyPatch::neighbProcNo
int neighbProcNo() const
Return neigbour processor number.
Definition: processorPolyPatch.H:255
Foam::argList
Extract command arguments and options from the supplied argc and argv parameters.
Definition: argList.H:97
Foam::IOobject::MUST_READ
@ MUST_READ
Definition: IOobject.H:108
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::IOobject::instance
const fileName & instance() const
Definition: IOobject.H:350
defaultFacesType
word defaultFacesType
Definition: readKivaGrid.H:461
syncTools.H
Foam::conformationSurfaces
Definition: conformationSurfaces.H:53
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
Foam::MeshObject< polyMesh, UpdateableMeshObject, pointMesh >::New
static const pointMesh & New(const polyMesh &mesh)
Definition: MeshObject.C:44
Foam::isoSurfaceCell
A surface formed by the iso value. After "Polygonising A Scalar Field Using Tetrahedrons",...
Definition: isoSurfaceCell.H:64
patchTypes
wordList patchTypes(nPatches)
Foam::DynamicField
Dynamically sized Field.
Definition: DynamicField.H:49
nPoints
label nPoints
Definition: gmvOutputHeader.H:2
Foam::primitiveMesh::nCells
label nCells() const
Definition: primitiveMeshI.H:64
Foam::boundBox::span
vector span() const
The bounding box span (from minimum to maximum)
Definition: boundBoxI.H:84
Foam::fvMesh::write
virtual bool write() const
Write mesh using IO settings from time.
Definition: fvMesh.C:873
Foam::IOobject::MUST_READ_IF_MODIFIED
@ MUST_READ_IF_MODIFIED
Definition: IOobject.H:109
getMergeDistance
scalar getMergeDistance(const argList &args, const Time &runTime, const boundBox &bb)
Definition: redistributePar.C:89
searchableSurfaces.H
Foam::primitiveMesh::nPoints
label nPoints() const
Definition: primitiveMeshI.H:35
Foam::IOstream::ASCII
@ ASCII
Definition: IOstream.H:88
backgroundMeshDecomposition.H
cellModeller.H
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::volumeType
Definition: volumeType.H:54
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::TimePaths::system
const word & system() const
Return system name.
Definition: TimePaths.H:120
Foam::functionObjectList::off
virtual void off()
Switch the function objects off.
Definition: functionObjectList.C:178
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
Foam::cellModeller::lookup
static const cellModel * lookup(const word &)
Look up a model by name and return a pointer to the model or NULL.
Definition: cellModeller.C:91
argList.H
Foam::decompositionModel::decomposer
decompositionMethod & decomposer() const
Definition: decompositionModel.H:110
Foam::polyMesh::faceOwner
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1017
conformationSurfaces.H
Foam::processorPolyPatch
Neighbour processor patch.
Definition: processorPolyPatch.H:55
main
int main(int argc, char *argv[])
Definition: postCalc.C:54
patchNames
wordList patchNames(nPatches)
Foam::boundBox::min
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:54
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:73
meshPtr
static fvMesh * meshPtr
Definition: globalFoam.H:52
Foam::Vector::x
const Cmpt & x() const
Definition: VectorI.H:65
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::primitiveMesh::nInternalFaces
label nInternalFaces() const
Definition: primitiveMeshI.H:52
Foam::cellShapeControl
Definition: cellShapeControl.H:62
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
fld
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< ' ';}gmvFile<< nl;forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam::writer
Base class for graphics format writing. Entry points are.
Definition: writer.H:78
Foam::cellShape
An analytical geometric cellShape.
Definition: cellShape.H:69
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::Vector::z
const Cmpt & z() const
Definition: VectorI.H:77
Foam::isDir
bool isDir(const fileName &)
Does the name exist as a DIRECTORY in the file system?
Definition: POSIX.C:615
Foam::globalIndex
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:63
Foam::e
const double e
Elementary charge.
Definition: doubleFloat.H:94
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::hex
IOstream & hex(IOstream &io)
Definition: IOstream.H:564
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:405
Foam::UPstream::master
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:399
Foam::polyMesh::setInstance
void setInstance(const fileName &)
Set the instance for mesh files.
Definition: polyMeshIO.C:32
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::polyMesh::bounds
const boundBox & bounds() const
Return mesh bounding box.
Definition: polyMesh.H:427
Foam::decompositionModel::New
static const decompositionModel & New(const polyMesh &mesh, const fileName &decompDictFile="")
Read (optionallly from absolute path) & register on mesh.
Definition: decompositionModel.C:108
Foam::List::setSize
void setSize(const label)
Reset size of List.
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
Foam::decompositionMethod::decompose
virtual labelList decompose(const pointField &points, const scalarField &pointWeights)
Return for every coordinate the wanted processor number.
Definition: decompositionMethod.H:126
Foam::primitiveMesh::cellCentres
const vectorField & cellCentres() const
Definition: primitiveMeshCellCentresAndVols.C:211
Foam::polyMesh::faces
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1004
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::Pstream::gatherList
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
Definition: gatherScatterList.C:49
DynamicField.H
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
Foam::fvMesh::C
const volVectorField & C() const
Return cell centres as volVectorField.
Definition: fvMeshGeometry.C:369
Foam::primitiveMesh::nFaces
label nFaces() const
Definition: primitiveMeshI.H:58
Foam::Time::path
fileName path() const
Return path.
Definition: Time.H:281
defaultFacesName
word defaultFacesName
Definition: readKivaGrid.H:460
Foam::IOstream::defaultPrecision
static unsigned int defaultPrecision()
Return the default precision.
Definition: IOstream.H:461
f
labelList f(nPoints)
Foam::Vector< label >
Foam::Time::timeName
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:741
vtkSurfaceWriter.H
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::searchableSurfaces
Container for searchableSurfaces.
Definition: searchableSurfaces.H:53
Foam::TimePaths::constant
const word & constant() const
Return constant name.
Definition: TimePaths.H:130
Foam::backgroundMeshDecomposition
Store a background polyMesh to use for the decomposition of space and queries for parallel conformalV...
Definition: backgroundMeshDecomposition.H:92
createMesh.H
points
const pointField & points
Definition: gmvOutputHeader.H:1
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
createTime.H
Foam::argList::optionFound
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:108
Foam::boundBox
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:55
cellShape.H
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
Foam::Vector::y
const Cmpt & y() const
Definition: VectorI.H:71
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
decompositionModel.H
Foam::vtkSurfaceWriter
A surfaceWriter for VTK legacy format.
Definition: vtkSurfaceWriter.H:48
Foam::cellModel
Maps a geometry to a set of cell primitives, which enables geometric cell data to be calculated witho...
Definition: cellModel.H:64
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::polyMesh::globalData
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1143
Foam::PatchTools::gatherAndMerge
static void gatherAndMerge(const scalar mergeDist, const PrimitivePatch< Face, FaceList, PointField, PointType > &p, Field< PointType > &mergedPoints, List< Face > &mergedFaces, labelList &pointMergeMap)
Gather points and faces onto master and merge into single patch.
Definition: PatchToolsGatherAndMerge.C:43
defaultMergeTol
static const scalar defaultMergeTol
Definition: reconstructParMesh.C:57
Foam::dictionary::subDict
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:631
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
args
Foam::argList args(argc, argv)
Foam::argList::optionReadIfPresent
bool optionReadIfPresent(const word &opt, T &) const
Read a value from the named option if present.
Definition: argListI.H:198
Foam::fvMeshDistribute
Sends/receives parts of mesh+fvfields to neighbouring processors. Used in load balancing.
Definition: fvMeshDistribute.H:70
Foam::searchableSurfacesQueries::findNearest
static void findNearest(const PtrList< searchableSurface > &, const labelList &surfacesToTest, const pointField &, const scalarField &nearestDistSqr, labelList &surfaces, List< pointIndexHit > &)
Find nearest. Return -1 (and a miss()) or surface and nearest.
Definition: searchableSurfacesQueries.C:625
Foam::magSqr
dimensioned< scalar > magSqr(const dimensioned< Type > &)
rndGen
cachedRandom rndGen(label(0), -1)
Foam::syncTools::syncPointList
static void syncPointList(const polyMesh &, List< T > &, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh points.
Definition: syncToolsTemplates.C:984
writeMesh
void writeMesh(const string &msg, const meshRefinement &meshRefiner, const meshRefinement::debugType debugLevel, const meshRefinement::writeType writeLevel)
Definition: snappyHexMesh.C:580
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatchTemplate.H:88