writeFields.C
Go to the documentation of this file.
1 #include "writeFields.H"
2 #include "volFields.H"
3 #include "surfaceFields.H"
4 #include "polyMeshTools.H"
5 #include "syncTools.H"
6 #include "tetrahedron.H"
7 #include "regionSplit.H"
8 #include "wallDist.H"
9 #include "cellAspectRatio.H"
10 
11 using namespace Foam;
12 
13 void maxFaceToCell
14 (
15  const scalarField& faceData,
16  volScalarField& cellData
17 )
18 {
19  const cellList& cells = cellData.mesh().cells();
20 
21  scalarField& cellFld = cellData.ref();
22 
23  cellFld = -GREAT;
24  forAll(cells, cellI)
25  {
26  const cell& cFaces = cells[cellI];
27  forAll(cFaces, i)
28  {
29  cellFld[cellI] = max(cellFld[cellI], faceData[cFaces[i]]);
30  }
31  }
32 
33  forAll(cellData.boundaryField(), patchI)
34  {
35  fvPatchScalarField& fvp = cellData.boundaryFieldRef()[patchI];
36 
37  fvp = fvp.patch().patchSlice(faceData);
38  }
39  cellData.correctBoundaryConditions();
40 }
41 
42 
43 void minFaceToCell
44 (
45  const scalarField& faceData,
46  volScalarField& cellData
47 )
48 {
49  const cellList& cells = cellData.mesh().cells();
50 
51  scalarField& cellFld = cellData.ref();
52 
53  cellFld = GREAT;
54  forAll(cells, cellI)
55  {
56  const cell& cFaces = cells[cellI];
57  forAll(cFaces, i)
58  {
59  cellFld[cellI] = min(cellFld[cellI], faceData[cFaces[i]]);
60  }
61  }
62 
63  forAll(cellData.boundaryField(), patchI)
64  {
65  fvPatchScalarField& fvp = cellData.boundaryFieldRef()[patchI];
66 
67  fvp = fvp.patch().patchSlice(faceData);
68  }
69  cellData.correctBoundaryConditions();
70 }
71 
72 
73 void minFaceToCell
74 (
75  const surfaceScalarField& faceData,
76  volScalarField& cellData,
77  const bool correctBoundaryConditions
78 )
79 {
80  scalarField& cellFld = cellData.ref();
81 
82  cellFld = GREAT;
83 
84  const labelUList& own = cellData.mesh().owner();
85  const labelUList& nei = cellData.mesh().neighbour();
86 
87  // Internal faces
88  forAll(own, facei)
89  {
90  cellFld[own[facei]] = min(cellFld[own[facei]], faceData[facei]);
91  cellFld[nei[facei]] = min(cellFld[nei[facei]], faceData[facei]);
92  }
93 
94  // Patch faces
95  forAll(faceData.boundaryField(), patchi)
96  {
97  const fvsPatchScalarField& fvp = faceData.boundaryField()[patchi];
98  const labelUList& fc = fvp.patch().faceCells();
99 
100  forAll(fc, i)
101  {
102  cellFld[fc[i]] = min(cellFld[fc[i]], fvp[i]);
103  }
104  }
105 
106  volScalarField::Boundary& bfld = cellData.boundaryFieldRef();
107 
108  forAll(bfld, patchi)
109  {
110  bfld[patchi] = faceData.boundaryField()[patchi];
111  }
113  {
114  cellData.correctBoundaryConditions();
115  }
116 }
117 
118 
119 void writeSurfaceField
120 (
121  const fvMesh& mesh,
122  const fileName& fName,
123  const scalarField& faceData
124 )
125 {
126  // Write single surfaceScalarField
127 
129  (
130  IOobject
131  (
132  fName,
133  mesh.time().timeName(),
134  mesh,
138  ),
139  mesh,
141  );
142  fld.primitiveFieldRef() = faceData;
143  //fld.correctBoundaryConditions();
144  Info<< " Writing face data to " << fName << endl;
145  fld.write();
146 }
147 
148 
150 (
151  const fvMesh& mesh,
152  const wordHashSet& selectedFields,
153  const bool writeFaceFields
154 )
155 {
156  if (selectedFields.empty())
157  {
158  return;
159  }
160 
161  Info<< "Writing fields with mesh quality parameters" << endl;
162 
163  if (selectedFields.found("nonOrthoAngle"))
164  {
165  //- Face based orthogonality
166  const scalarField faceOrthogonality
167  (
169  (
170  mesh,
171  mesh.faceAreas(),
172  mesh.cellCentres()
173  )
174  );
175 
176  //- Face based angle
177  const scalarField nonOrthoAngle
178  (
179  radToDeg
180  (
181  Foam::acos(min(scalar(1), max(scalar(-1), faceOrthogonality)))
182  )
183  );
184 
185  //- Cell field - max of either face
186  volScalarField cellNonOrthoAngle
187  (
188  IOobject
189  (
190  "nonOrthoAngle",
191  mesh.time().timeName(),
192  mesh,
196  ),
197  mesh,
199  );
200  //- Take max
201  maxFaceToCell(nonOrthoAngle, cellNonOrthoAngle);
202  Info<< " Writing non-orthogonality (angle) to "
203  << cellNonOrthoAngle.name() << endl;
204  cellNonOrthoAngle.write();
205 
206  if (writeFaceFields)
207  {
208  writeSurfaceField
209  (
210  mesh,
211  "face_nonOrthoAngle",
212  SubField<scalar>(nonOrthoAngle, mesh.nInternalFaces())
213  );
214  }
215  }
216 
217  if (selectedFields.found("faceWeight"))
218  {
219  volScalarField cellWeights
220  (
221  IOobject
222  (
223  "faceWeight",
224  mesh.time().timeName(),
225  mesh,
229  ),
230  mesh,
232  wordList // wanted bc types
233  (
234  mesh.boundary().size(),
236  ),
237  mesh.weights().boundaryField().types() // current bc types
238  );
239  //- Take min
240  minFaceToCell(mesh.weights(), cellWeights, false);
241  Info<< " Writing face interpolation weights (0..0.5) to "
242  << cellWeights.name() << endl;
243  cellWeights.write();
244 
245  if (writeFaceFields)
246  {
247  writeSurfaceField(mesh, "face_faceWeight", mesh.weights());
248  }
249  }
250 
251 
252  // Skewness
253  // ~~~~~~~~
254 
255  if (selectedFields.found("skewness"))
256  {
257  //- Face based skewness
258  const scalarField faceSkewness
259  (
261  (
262  mesh,
263  mesh.points(),
264  mesh.faceCentres(),
265  mesh.faceAreas(),
266  mesh.cellCentres()
267  )
268  );
269 
270  //- Cell field - max of either face
271  volScalarField cellSkewness
272  (
273  IOobject
274  (
275  "skewness",
276  mesh.time().timeName(),
277  mesh,
281  ),
282  mesh,
284  );
285  //- Take max
286  maxFaceToCell(faceSkewness, cellSkewness);
287  Info<< " Writing face skewness to " << cellSkewness.name() << endl;
288  cellSkewness.write();
289 
290  if (writeFaceFields)
291  {
292  writeSurfaceField
293  (
294  mesh,
295  "face_skewness",
296  SubField<scalar>(faceSkewness, mesh.nInternalFaces())
297  );
298  }
299  }
300 
301 
302  // cellDeterminant
303  // ~~~~~~~~~~~~~~~
304 
305  if (selectedFields.found("cellDeterminant"))
306  {
307  volScalarField cellDeterminant
308  (
309  IOobject
310  (
311  "cellDeterminant",
312  mesh.time().timeName(),
313  mesh,
317  ),
318  mesh,
321  );
322  cellDeterminant.primitiveFieldRef() =
324  (
325  mesh,
326  mesh.geometricD(),
327  mesh.faceAreas(),
329  );
330  cellDeterminant.correctBoundaryConditions();
331  Info<< " Writing cell determinant to "
332  << cellDeterminant.name() << endl;
333  cellDeterminant.write();
334  }
335 
336 
337  // Aspect ratio
338  // ~~~~~~~~~~~~
339 
340  if (selectedFields.found("aspectRatio"))
341  {
342  volScalarField aspectRatio
343  (
344  IOobject
345  (
346  "aspectRatio",
347  mesh.time().timeName(),
348  mesh,
352  ),
353  mesh,
356  );
357 
358 
359  scalarField cellOpenness;
361  (
362  mesh,
363  mesh.geometricD(),
364  mesh.faceAreas(),
365  mesh.cellVolumes(),
366  cellOpenness,
367  aspectRatio.ref()
368  );
369 
370  aspectRatio.correctBoundaryConditions();
371  Info<< " Writing aspect ratio to " << aspectRatio.name() << endl;
372  aspectRatio.write();
373  }
374 
375  if (selectedFields.found("cellAspectRatio"))
376  {
377  volScalarField aspectRatio
378  (
379  IOobject
380  (
381  "cellAspectRatio",
382  mesh.time().timeName(),
383  mesh,
387  ),
388  mesh,
391  );
392 
393  aspectRatio.ref().field() = cellAspectRatio(mesh);
394 
395  aspectRatio.correctBoundaryConditions();
396  Info<< " Writing approximate aspect ratio to "
397  << aspectRatio.name() << endl;
398  aspectRatio.write();
399  }
400 
401 
402  // cell type
403  // ~~~~~~~~~
404 
405  if (selectedFields.found("cellShapes"))
406  {
407  volScalarField shape
408  (
409  IOobject
410  (
411  "cellShapes",
412  mesh.time().timeName(),
413  mesh,
417  ),
418  mesh,
421  );
423  forAll(cellShapes, cellI)
424  {
425  const cellModel& model = cellShapes[cellI].model();
426  shape[cellI] = model.index();
427  }
428  shape.correctBoundaryConditions();
429  Info<< " Writing cell shape (hex, tet etc.) to " << shape.name()
430  << endl;
431  shape.write();
432  }
433 
434  if (selectedFields.found("cellVolume"))
435  {
437  (
438  IOobject
439  (
440  "cellVolume",
441  mesh.time().timeName(),
442  mesh,
446  ),
447  mesh,
449  );
450  V.ref() = mesh.V();
451  Info<< " Writing cell volume to " << V.name() << endl;
452  V.write();
453  }
454 
455  if (selectedFields.found("cellVolumeRatio"))
456  {
457  const scalarField faceVolumeRatio
458  (
460  (
461  mesh,
462  mesh.V()
463  )
464  );
465 
466  volScalarField cellVolumeRatio
467  (
468  IOobject
469  (
470  "cellVolumeRatio",
471  mesh.time().timeName(),
472  mesh,
476  ),
477  mesh,
479  );
480  //- Take min
481  minFaceToCell(faceVolumeRatio, cellVolumeRatio);
482  Info<< " Writing cell volume ratio to "
483  << cellVolumeRatio.name() << endl;
484  cellVolumeRatio.write();
485 
486  if (writeFaceFields)
487  {
488  writeSurfaceField
489  (
490  mesh,
491  "face_cellVolumeRatio",
492  SubField<scalar>(faceVolumeRatio, mesh.nInternalFaces())
493  );
494  }
495  }
496 
497  // minTetVolume
498  if (selectedFields.found("minTetVolume"))
499  {
500  volScalarField minTetVolume
501  (
502  IOobject
503  (
504  "minTetVolume",
505  mesh.time().timeName(),
506  mesh,
510  ),
511  mesh,
512  dimensionedScalar("minTetVolume", dimless, GREAT),
514  );
515 
516 
517  const labelList& own = mesh.faceOwner();
518  const labelList& nei = mesh.faceNeighbour();
519  const pointField& p = mesh.points();
520  forAll(own, facei)
521  {
522  const face& f = mesh.faces()[facei];
523  const point& fc = mesh.faceCentres()[facei];
524 
525  {
526  const point& ownCc = mesh.cellCentres()[own[facei]];
527  scalar& ownVol = minTetVolume[own[facei]];
528  forAll(f, fp)
529  {
530  scalar tetQual = tetPointRef
531  (
532  p[f[fp]],
533  p[f.nextLabel(fp)],
534  ownCc,
535  fc
536  ).quality();
537  ownVol = min(ownVol, tetQual);
538  }
539  }
540  if (mesh.isInternalFace(facei))
541  {
542  const point& neiCc = mesh.cellCentres()[nei[facei]];
543  scalar& neiVol = minTetVolume[nei[facei]];
544  forAll(f, fp)
545  {
546  scalar tetQual = tetPointRef
547  (
548  p[f[fp]],
549  p[f.nextLabel(fp)],
550  fc,
551  neiCc
552  ).quality();
553  neiVol = min(neiVol, tetQual);
554  }
555  }
556  }
557 
558  minTetVolume.correctBoundaryConditions();
559  Info<< " Writing minTetVolume to " << minTetVolume.name() << endl;
560  minTetVolume.write();
561  }
562 
563  // minPyrVolume
564  if (selectedFields.found("minPyrVolume"))
565  {
566  volScalarField minPyrVolume
567  (
568  IOobject
569  (
570  "minPyrVolume",
571  mesh.time().timeName(),
572  mesh,
576  ),
577  mesh,
578  dimensionedScalar("minPyrVolume", dimless, GREAT),
580  );
581 
582  // Get owner and neighbour pyr volumes
583  scalarField ownPyrVol(mesh.nFaces());
584  scalarField neiPyrVol(mesh.nInternalFaces());
586  (
587  mesh,
588  mesh.points(),
589  mesh.cellCentres(),
590 
591  ownPyrVol,
592  neiPyrVol
593  );
594 
595  // Get min pyr vol per cell
596  scalarField& cellFld = minPyrVolume.ref();
597  cellFld = GREAT;
598 
599  const labelUList& own = mesh.owner();
600  const labelUList& nei = mesh.neighbour();
601 
602  // Internal faces
603  forAll(own, facei)
604  {
605  cellFld[own[facei]] = min(cellFld[own[facei]], ownPyrVol[facei]);
606  cellFld[nei[facei]] = min(cellFld[nei[facei]], neiPyrVol[facei]);
607  }
608 
609  // Patch faces
610  for (const auto& fvp : minPyrVolume.boundaryField())
611  {
612  const labelUList& fc = fvp.patch().faceCells();
613 
614  forAll(fc, i)
615  {
616  const label meshFacei = fvp.patch().start();
617  cellFld[fc[i]] = min(cellFld[fc[i]], ownPyrVol[meshFacei]);
618  }
619  }
620 
621  minPyrVolume.correctBoundaryConditions();
622  Info<< " Writing minPyrVolume to " << minPyrVolume.name() << endl;
623  minPyrVolume.write();
624 
625  if (writeFaceFields)
626  {
627  scalarField minFacePyrVol(neiPyrVol);
628  minFacePyrVol = min
629  (
630  minFacePyrVol,
631  SubField<scalar>(ownPyrVol, mesh.nInternalFaces())
632  );
633  writeSurfaceField(mesh, "face_minPyrVolume", minFacePyrVol);
634  }
635  }
636 
637  if (selectedFields.found("cellRegion"))
638  {
639  volScalarField cellRegion
640  (
641  IOobject
642  (
643  "cellRegion",
644  mesh.time().timeName(),
645  mesh,
649  ),
650  mesh,
652  );
653 
654  regionSplit rs(mesh);
655  forAll(rs, celli)
656  {
657  cellRegion[celli] = rs[celli];
658  }
659  cellRegion.correctBoundaryConditions();
660  Info<< " Writing cell region to " << cellRegion.name() << endl;
661  cellRegion.write();
662  }
663 
664  if (selectedFields.found("wallDistance"))
665  {
666  // See if wallDist.method entry in fvSchemes before calling factory
667  // method of wallDist. Have 'failing' version of wallDist::New instead?
668  const dictionary& schemesDict =
669  static_cast<const fvSchemes&>(mesh).schemesDict();
670  if (schemesDict.found("wallDist"))
671  {
672  if (schemesDict.subDict("wallDist").found("method"))
673  {
674  // Wall distance
675  volScalarField y("wallDistance", wallDist::New(mesh).y());
676  Info<< " Writing wall distance to " << y.name() << endl;
677  y.write();
678 
679  // Wall-reflection vectors
680  //const volVectorField& n = wallDist::New(mesh).n();
681  //Info<< " Writing wall normal to " << n.name() << endl;
682  //n.write();
683  }
684  }
685  }
686 
687  if (selectedFields.found("cellZone"))
688  {
690  (
691  IOobject
692  (
693  "cellZone",
694  mesh.time().timeName(),
695  mesh,
699  ),
700  mesh,
702  );
703 
704  const cellZoneMesh& czs = mesh.cellZones();
705  for (const auto& zone : czs)
706  {
708  }
709 
710  cellZone.correctBoundaryConditions();
711  Info<< " Writing cell zoning to " << cellZone.name() << endl;
712  cellZone.write();
713  }
714  if (selectedFields.found("faceZone"))
715  {
716  // Determine for each face the zone index (scalar for ease of
717  // manipulation)
718  scalarField zoneID(mesh.nFaces(), -1);
719  const faceZoneMesh& czs = mesh.faceZones();
720  for (const auto& zone : czs)
721  {
722  UIndirectList<scalar>(zoneID, zone) = zone.index();
723  }
724 
725 
726  // Split into internal and boundary values
728  (
729  IOobject
730  (
731  "faceZone",
732  mesh.time().timeName(),
733  mesh,
737  ),
738  mesh,
740  );
741 
742  faceZone.primitiveFieldRef() =
744  surfaceScalarField::Boundary& bfld = faceZone.boundaryFieldRef();
745  for (auto& pfld : bfld)
746  {
747  const fvPatch& fvp = pfld.patch();
748  pfld == SubField<scalar>(zoneID, fvp.size(), fvp.start());
749  }
750 
751  //faceZone.correctBoundaryConditions();
752  Info<< " Writing face zoning to " << faceZone.name() << endl;
753  faceZone.write();
754  }
755 
756  Info<< endl;
757 }
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
Definition: fvPatchField.H:218
Foam::surfaceFields.
label start() const noexcept
The patch start within the polyMesh face list.
Definition: fvPatch.H:226
This class separates the mesh into distinct unconnected regions, each of which is then given a label ...
Definition: regionSplit.H:136
const List< T >::subList patchSlice(const List< T > &values) const
This patch slice from the complete list, which has size mesh::nFaces(), using the virtual patch size...
Definition: fvPatch.H:271
dimensionedScalar acos(const dimensionedScalar &ds)
A class for handling file names.
Definition: fileName.H:71
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:68
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:120
const fvPatch & patch() const noexcept
Return the patch.
Definition: fvPatchField.H:266
bool found(const Key &key) const
Same as contains()
Definition: HashTable.H:1333
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1117
tetrahedron< point, const point & > tetPointRef
A tetrahedron using referred points.
Definition: tetrahedron.H:72
static const wallDist & New(const fvMesh &mesh, Args &&... args)
Get existing or create a new MeshObject. Registered with typeName.
Definition: MeshObject.C:53
(Rough approximation of) cell aspect ratio
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
const cellShapeList & cellShapes() const
Return cell shapes.
wordList types() const
Return a list of the patch types.
cellMask correctBoundaryConditions()
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
const fvPatch & patch() const noexcept
Return the patch.
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:362
label nFaces() const noexcept
Number of mesh faces.
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
SubField is a Field obtained as a section of another Field, without its own allocation. SubField is derived from a SubList rather than a List.
Definition: Field.H:63
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:453
static bitSet getInternalOrCoupledFaces(const polyMesh &mesh)
Get per face whether it is internal or coupled.
Definition: syncTools.C:169
bool isInternalFace(const label faceIndex) const noexcept
Return true if given face label is internal to the mesh.
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry (const access) with the given keyword.
Definition: dictionaryI.H:100
const labelUList & neighbour() const
Internal face neighbour.
Definition: fvMesh.H:549
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1073
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:414
const Vector< label > & geometricD() const
Return the vector of geometric directions in mesh.
Definition: polyMesh.C:860
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
scalar y
virtual const surfaceScalarField & weights() const
Return reference to linear difference weighting factors.
dynamicFvMesh & mesh
const cellShapeList & cells
Base class for mesh zones.
Definition: zone.H:59
virtual const labelUList & faceCells() const
Return faceCells.
Definition: fvPatch.C:107
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:113
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1111
static const word null
An empty word.
Definition: word.H:84
cellShapeList cellShapes
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
static tmp< scalarField > faceSkewness(const polyMesh &mesh, const pointField &points, const vectorField &fCtrs, const vectorField &fAreas, const vectorField &cellCtrs)
Generate skewness field.
Definition: polyMeshTools.C:84
label nInternalFaces() const noexcept
Number of internal faces.
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1098
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
const vectorField & cellCentres() const
static tmp< scalarField > volRatio(const polyMesh &mesh, const scalarField &vol)
Generate volume ratio field.
void writeFields(const fvMesh &mesh, const wordHashSet &selectedFields, const bool writeFaceFields)
virtual label size() const
Patch size is the number of faces, but can be overloaded.
Definition: fvPatch.H:234
virtual void write(Ostream &os) const
Write.
Definition: faceZone.C:609
const Mesh & mesh() const noexcept
Return mesh.
bool empty() const noexcept
True if the hash table is empty.
Definition: HashTable.H:326
label index() const noexcept
The index of this zone in the zone list.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
static word timeName(const scalar t, const int precision=precision_)
Return time name of given scalar time formatted with the given precision.
Definition: Time.C:770
static void facePyramidVolume(const primitiveMesh &mesh, const pointField &points, const vectorField &cellCtrs, scalarField &ownPyrVol, scalarField &neiPyrVol)
Generate face pyramid volume fields.
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition: polyMesh.H:646
static tmp< scalarField > cellDeterminant(const primitiveMesh &mesh, const Vector< label > &directions, const vectorField &faceAreas, const bitSet &internalOrCoupledFace)
Generate cell determinant field. Normalised to 1 for an internal cube.
labelList f(nPoints)
A subset of mesh cells.
Definition: cellZone.H:58
virtual void write(Ostream &os) const
Write.
Definition: zone.C:203
const labelUList & owner() const
Internal face owner. Note bypassing virtual mechanism so.
Definition: fvMesh.H:541
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
const Field< Type > & field() const noexcept
Return const-reference to the field values.
const vectorField & faceCentres() const
static tmp< scalarField > faceOrthogonality(const polyMesh &mesh, const vectorField &fAreas, const vectorField &cellCtrs)
Generate orthogonality field. (1 for fully orthogonal, < 1 for non-orthogonal)
Definition: polyMeshTools.C:30
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
label index() const noexcept
Return index of model in the model list.
Definition: cellModelI.H:30
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
constexpr scalar radToDeg(const scalar rad) noexcept
Conversion from radians to degrees.
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:53
Selector class for finite volume differencing schemes. fvMesh is derived from fvSchemes so that all f...
Definition: fvSchemes.H:51
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:79
A List with indirect addressing. Like IndirectList but does not store addressing. ...
Definition: faMatrix.H:56
const vectorField & faceAreas() const
const word & name() const noexcept
The zone name.
Nothing to be read.
static const word & calculatedType() noexcept
The type name for calculated patch fields.
Definition: fvPatchField.H:201
Maps a geometry to a set of cell primitives.
Definition: cellModel.H:72
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
const cellZoneMesh & cellZones() const noexcept
Return cell zone mesh.
Definition: polyMesh.H:654
void correctBoundaryConditions()
Correct boundary field.
messageStream Info
Information stream (stdout output on master, null elsewhere)
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef()
const fvBoundaryMesh & boundary() const noexcept
Return reference to boundary mesh.
Definition: fvMesh.H:397
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:60
volScalarField & p
static void cellClosedness(const primitiveMesh &mesh, const Vector< label > &meshD, const vectorField &areas, const scalarField &vols, scalarField &openness, scalarField &aratio)
Generate cell openness and cell aspect ratio field.
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:171
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
Do not request registration (bool: false)
Namespace for OpenFOAM.
const scalarField & cellVolumes() const
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:133