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