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"
6 #include "syncTools.H"
7 #include "tetPointRef.H"
8 #include "regionSplit.H"
9 #include "wallDist.H"
10 #include "cellAspectRatio.H"
11 
12 using namespace Foam;
13 
14 void maxFaceToCell
15 (
16  const scalarField& faceData,
17  volScalarField& cellData
18 )
19 {
20  const cellList& cells = cellData.mesh().cells();
21 
22  scalarField& cellFld = cellData.ref();
23 
24  cellFld = -GREAT;
25  forAll(cells, cellI)
26  {
27  const cell& cFaces = cells[cellI];
28  forAll(cFaces, i)
29  {
30  cellFld[cellI] = max(cellFld[cellI], faceData[cFaces[i]]);
31  }
32  }
33 
34  forAll(cellData.boundaryField(), patchI)
35  {
36  fvPatchScalarField& fvp = cellData.boundaryFieldRef()[patchI];
37 
38  fvp = fvp.patch().patchSlice(faceData);
39  }
40  cellData.correctBoundaryConditions();
41 }
42 
43 
44 void minFaceToCell
45 (
46  const scalarField& faceData,
47  volScalarField& cellData
48 )
49 {
50  const cellList& cells = cellData.mesh().cells();
51 
52  scalarField& cellFld = cellData.ref();
53 
54  cellFld = GREAT;
55  forAll(cells, cellI)
56  {
57  const cell& cFaces = cells[cellI];
58  forAll(cFaces, i)
59  {
60  cellFld[cellI] = min(cellFld[cellI], faceData[cFaces[i]]);
61  }
62  }
63 
64  forAll(cellData.boundaryField(), patchI)
65  {
66  fvPatchScalarField& fvp = cellData.boundaryFieldRef()[patchI];
67 
68  fvp = fvp.patch().patchSlice(faceData);
69  }
70  cellData.correctBoundaryConditions();
71 }
72 
73 
74 void minFaceToCell
75 (
76  const surfaceScalarField& faceData,
77  volScalarField& cellData,
78  const bool correctBoundaryConditions
79 )
80 {
81  scalarField& cellFld = cellData.ref();
82 
83  cellFld = GREAT;
84 
85  const labelUList& own = cellData.mesh().owner();
86  const labelUList& nei = cellData.mesh().neighbour();
87 
88  // Internal faces
89  forAll(own, facei)
90  {
91  cellFld[own[facei]] = min(cellFld[own[facei]], faceData[facei]);
92  cellFld[nei[facei]] = min(cellFld[nei[facei]], faceData[facei]);
93  }
94 
95  // Patch faces
96  forAll(faceData.boundaryField(), patchi)
97  {
98  const fvsPatchScalarField& fvp = faceData.boundaryField()[patchi];
99  const labelUList& fc = fvp.patch().faceCells();
100 
101  forAll(fc, i)
102  {
103  cellFld[fc[i]] = min(cellFld[fc[i]], fvp[i]);
104  }
105  }
106 
107  volScalarField::Boundary& bfld = cellData.boundaryFieldRef();
108 
109  forAll(bfld, patchi)
110  {
111  bfld[patchi] = faceData.boundaryField()[patchi];
112  }
114  {
115  cellData.correctBoundaryConditions();
116  }
117 }
118 
119 
120 void writeSurfaceField
121 (
122  const fvMesh& mesh,
123  const fileName& fName,
124  const scalarField& faceData
125 )
126 {
127  // Write single surfaceScalarField
128 
130  (
131  IOobject
132  (
133  fName,
134  mesh.time().timeName(),
135  mesh,
138  false
139  ),
140  mesh,
142  calculatedFvsPatchScalarField::typeName
143  );
144  fld.primitiveFieldRef() = faceData;
145  //fld.correctBoundaryConditions();
146  Info<< " Writing face data to " << fName << endl;
147  fld.write();
148 }
149 
150 
152 (
153  const fvMesh& mesh,
154  const wordHashSet& selectedFields,
155  const bool writeFaceFields
156 )
157 {
158  if (selectedFields.empty())
159  {
160  return;
161  }
162 
163  Info<< "Writing fields with mesh quality parameters" << endl;
164 
165  if (selectedFields.found("nonOrthoAngle"))
166  {
167  //- Face based orthogonality
168  const scalarField faceOrthogonality
169  (
171  (
172  mesh,
173  mesh.faceAreas(),
174  mesh.cellCentres()
175  )
176  );
177 
178  //- Face based angle
179  const scalarField nonOrthoAngle
180  (
181  radToDeg
182  (
183  Foam::acos(min(scalar(1), max(scalar(-1), faceOrthogonality)))
184  )
185  );
186 
187  //- Cell field - max of either face
188  volScalarField cellNonOrthoAngle
189  (
190  IOobject
191  (
192  "nonOrthoAngle",
193  mesh.time().timeName(),
194  mesh,
197  false
198  ),
199  mesh,
201  calculatedFvPatchScalarField::typeName
202  );
203  //- Take max
204  maxFaceToCell(nonOrthoAngle, cellNonOrthoAngle);
205  Info<< " Writing non-orthogonality (angle) to "
206  << cellNonOrthoAngle.name() << endl;
207  cellNonOrthoAngle.write();
208 
209  if (writeFaceFields)
210  {
211  writeSurfaceField
212  (
213  mesh,
214  "face_nonOrthoAngle",
215  SubField<scalar>(nonOrthoAngle, mesh.nInternalFaces())
216  );
217  }
218  }
219 
220  if (selectedFields.found("faceWeight"))
221  {
222  volScalarField cellWeights
223  (
224  IOobject
225  (
226  "faceWeight",
227  mesh.time().timeName(),
228  mesh,
231  false
232  ),
233  mesh,
235  wordList // wanted bc types
236  (
237  mesh.boundary().size(),
238  calculatedFvPatchScalarField::typeName
239  ),
240  mesh.weights().boundaryField().types() // current bc types
241  );
242  //- Take min
243  minFaceToCell(mesh.weights(), cellWeights, false);
244  Info<< " Writing face interpolation weights (0..0.5) to "
245  << cellWeights.name() << endl;
246  cellWeights.write();
247 
248  if (writeFaceFields)
249  {
250  writeSurfaceField(mesh, "face_faceWeight", mesh.weights());
251  }
252  }
253 
254 
255  // Skewness
256  // ~~~~~~~~
257 
258  if (selectedFields.found("skewness"))
259  {
260  //- Face based skewness
261  const scalarField faceSkewness
262  (
264  (
265  mesh,
266  mesh.points(),
267  mesh.faceCentres(),
268  mesh.faceAreas(),
269  mesh.cellCentres()
270  )
271  );
272 
273  //- Cell field - max of either face
274  volScalarField cellSkewness
275  (
276  IOobject
277  (
278  "skewness",
279  mesh.time().timeName(),
280  mesh,
283  false
284  ),
285  mesh,
287  calculatedFvPatchScalarField::typeName
288  );
289  //- Take max
290  maxFaceToCell(faceSkewness, cellSkewness);
291  Info<< " Writing face skewness to " << cellSkewness.name() << endl;
292  cellSkewness.write();
293 
294  if (writeFaceFields)
295  {
296  writeSurfaceField
297  (
298  mesh,
299  "face_skewness",
300  SubField<scalar>(faceSkewness, mesh.nInternalFaces())
301  );
302  }
303  }
304 
305 
306  // cellDeterminant
307  // ~~~~~~~~~~~~~~~
308 
309  if (selectedFields.found("cellDeterminant"))
310  {
311  volScalarField cellDeterminant
312  (
313  IOobject
314  (
315  "cellDeterminant",
316  mesh.time().timeName(),
317  mesh,
320  false
321  ),
322  mesh,
324  zeroGradientFvPatchScalarField::typeName
325  );
326  cellDeterminant.primitiveFieldRef() =
328  (
329  mesh,
330  mesh.geometricD(),
331  mesh.faceAreas(),
333  );
334  cellDeterminant.correctBoundaryConditions();
335  Info<< " Writing cell determinant to "
336  << cellDeterminant.name() << endl;
337  cellDeterminant.write();
338  }
339 
340 
341  // Aspect ratio
342  // ~~~~~~~~~~~~
343 
344  if (selectedFields.found("aspectRatio"))
345  {
346  volScalarField aspectRatio
347  (
348  IOobject
349  (
350  "aspectRatio",
351  mesh.time().timeName(),
352  mesh,
355  false
356  ),
357  mesh,
359  zeroGradientFvPatchScalarField::typeName
360  );
361 
362 
363  scalarField cellOpenness;
365  (
366  mesh,
367  mesh.geometricD(),
368  mesh.faceAreas(),
369  mesh.cellVolumes(),
370  cellOpenness,
371  aspectRatio.ref()
372  );
373 
374  aspectRatio.correctBoundaryConditions();
375  Info<< " Writing aspect ratio to " << aspectRatio.name() << endl;
376  aspectRatio.write();
377  }
378 
379  if (selectedFields.found("cellAspectRatio"))
380  {
381  volScalarField aspectRatio
382  (
383  IOobject
384  (
385  "cellAspectRatio",
386  mesh.time().timeName(),
387  mesh,
390  false
391  ),
392  mesh,
394  zeroGradientFvPatchScalarField::typeName
395  );
396 
397  aspectRatio.ref().field() = cellAspectRatio(mesh);
398 
399  aspectRatio.correctBoundaryConditions();
400  Info<< " Writing approximate aspect ratio to "
401  << aspectRatio.name() << endl;
402  aspectRatio.write();
403  }
404 
405 
406  // cell type
407  // ~~~~~~~~~
408 
409  if (selectedFields.found("cellShapes"))
410  {
411  volScalarField shape
412  (
413  IOobject
414  (
415  "cellShapes",
416  mesh.time().timeName(),
417  mesh,
420  false
421  ),
422  mesh,
424  zeroGradientFvPatchScalarField::typeName
425  );
427  forAll(cellShapes, cellI)
428  {
429  const cellModel& model = cellShapes[cellI].model();
430  shape[cellI] = model.index();
431  }
432  shape.correctBoundaryConditions();
433  Info<< " Writing cell shape (hex, tet etc.) to " << shape.name()
434  << endl;
435  shape.write();
436  }
437 
438  if (selectedFields.found("cellVolume"))
439  {
441  (
442  IOobject
443  (
444  "cellVolume",
445  mesh.time().timeName(),
446  mesh,
449  false
450  ),
451  mesh,
453  calculatedFvPatchScalarField::typeName
454  );
455  V.ref() = mesh.V();
456  Info<< " Writing cell volume to " << V.name() << endl;
457  V.write();
458  }
459 
460  if (selectedFields.found("cellVolumeRatio"))
461  {
462  const scalarField faceVolumeRatio
463  (
465  (
466  mesh,
467  mesh.V()
468  )
469  );
470 
471  volScalarField cellVolumeRatio
472  (
473  IOobject
474  (
475  "cellVolumeRatio",
476  mesh.time().timeName(),
477  mesh,
480  false
481  ),
482  mesh,
484  calculatedFvPatchScalarField::typeName
485  );
486  //- Take min
487  minFaceToCell(faceVolumeRatio, cellVolumeRatio);
488  Info<< " Writing cell volume ratio to "
489  << cellVolumeRatio.name() << endl;
490  cellVolumeRatio.write();
491 
492  if (writeFaceFields)
493  {
494  writeSurfaceField
495  (
496  mesh,
497  "face_cellVolumeRatio",
498  SubField<scalar>(faceVolumeRatio, mesh.nInternalFaces())
499  );
500  }
501  }
502 
503  // minTetVolume
504  if (selectedFields.found("minTetVolume"))
505  {
506  volScalarField minTetVolume
507  (
508  IOobject
509  (
510  "minTetVolume",
511  mesh.time().timeName(),
512  mesh,
515  false
516  ),
517  mesh,
518  dimensionedScalar("minTetVolume", dimless, GREAT),
519  zeroGradientFvPatchScalarField::typeName
520  );
521 
522 
523  const labelList& own = mesh.faceOwner();
524  const labelList& nei = mesh.faceNeighbour();
525  const pointField& p = mesh.points();
526  forAll(own, facei)
527  {
528  const face& f = mesh.faces()[facei];
529  const point& fc = mesh.faceCentres()[facei];
530 
531  {
532  const point& ownCc = mesh.cellCentres()[own[facei]];
533  scalar& ownVol = minTetVolume[own[facei]];
534  forAll(f, fp)
535  {
536  scalar tetQual = tetPointRef
537  (
538  p[f[fp]],
539  p[f.nextLabel(fp)],
540  ownCc,
541  fc
542  ).quality();
543  ownVol = min(ownVol, tetQual);
544  }
545  }
546  if (mesh.isInternalFace(facei))
547  {
548  const point& neiCc = mesh.cellCentres()[nei[facei]];
549  scalar& neiVol = minTetVolume[nei[facei]];
550  forAll(f, fp)
551  {
552  scalar tetQual = tetPointRef
553  (
554  p[f[fp]],
555  p[f.nextLabel(fp)],
556  fc,
557  neiCc
558  ).quality();
559  neiVol = min(neiVol, tetQual);
560  }
561  }
562  }
563 
564  minTetVolume.correctBoundaryConditions();
565  Info<< " Writing minTetVolume to " << minTetVolume.name() << endl;
566  minTetVolume.write();
567  }
568 
569  // minPyrVolume
570  if (selectedFields.found("minPyrVolume"))
571  {
572  volScalarField minPyrVolume
573  (
574  IOobject
575  (
576  "minPyrVolume",
577  mesh.time().timeName(),
578  mesh,
581  false
582  ),
583  mesh,
584  dimensionedScalar("minPyrVolume", dimless, GREAT),
585  zeroGradientFvPatchScalarField::typeName
586  );
587 
588  // Get owner and neighbour pyr volumes
589  scalarField ownPyrVol(mesh.nFaces());
590  scalarField neiPyrVol(mesh.nInternalFaces());
592  (
593  mesh,
594  mesh.points(),
595  mesh.cellCentres(),
596 
597  ownPyrVol,
598  neiPyrVol
599  );
600 
601  // Get min pyr vol per cell
602  scalarField& cellFld = minPyrVolume.ref();
603  cellFld = GREAT;
604 
605  const labelUList& own = mesh.owner();
606  const labelUList& nei = mesh.neighbour();
607 
608  // Internal faces
609  forAll(own, facei)
610  {
611  cellFld[own[facei]] = min(cellFld[own[facei]], ownPyrVol[facei]);
612  cellFld[nei[facei]] = min(cellFld[nei[facei]], neiPyrVol[facei]);
613  }
614 
615  // Patch faces
616  for (const auto& fvp : minPyrVolume.boundaryField())
617  {
618  const labelUList& fc = fvp.patch().faceCells();
619 
620  forAll(fc, i)
621  {
622  const label meshFacei = fvp.patch().start();
623  cellFld[fc[i]] = min(cellFld[fc[i]], ownPyrVol[meshFacei]);
624  }
625  }
626 
627  minPyrVolume.correctBoundaryConditions();
628  Info<< " Writing minPyrVolume to " << minPyrVolume.name() << endl;
629  minPyrVolume.write();
630 
631  if (writeFaceFields)
632  {
633  scalarField minFacePyrVol(neiPyrVol);
634  minFacePyrVol = min
635  (
636  minFacePyrVol,
637  SubField<scalar>(ownPyrVol, mesh.nInternalFaces())
638  );
639  writeSurfaceField(mesh, "face_minPyrVolume", minFacePyrVol);
640  }
641  }
642 
643  if (selectedFields.found("cellRegion"))
644  {
645  volScalarField cellRegion
646  (
647  IOobject
648  (
649  "cellRegion",
650  mesh.time().timeName(),
651  mesh,
654  false
655  ),
656  mesh,
658  calculatedFvPatchScalarField::typeName
659  );
660 
661  regionSplit rs(mesh);
662  forAll(rs, celli)
663  {
664  cellRegion[celli] = rs[celli];
665  }
666  cellRegion.correctBoundaryConditions();
667  Info<< " Writing cell region to " << cellRegion.name() << endl;
668  cellRegion.write();
669  }
670 
671  if (selectedFields.found("wallDistance"))
672  {
673  // See if wallDist.method entry in fvSchemes before calling factory
674  // method of wallDist. Have 'failing' version of wallDist::New instead?
675  const dictionary& schemesDict =
676  static_cast<const fvSchemes&>(mesh).schemesDict();
677  if (schemesDict.found("wallDist"))
678  {
679  if (schemesDict.subDict("wallDist").found("method"))
680  {
681  // Wall distance
682  volScalarField y("wallDistance", wallDist::New(mesh).y());
683  Info<< " Writing wall distance to " << y.name() << endl;
684  y.write();
685 
686  // Wall-reflection vectors
687  //const volVectorField& n = wallDist::New(mesh).n();
688  //Info<< " Writing wall normal to " << n.name() << endl;
689  //n.write();
690  }
691  }
692  }
693 
694  if (selectedFields.found("cellZone"))
695  {
697  (
698  IOobject
699  (
700  "cellZone",
701  mesh.time().timeName(),
702  mesh,
705  false
706  ),
707  mesh,
708  dimensionedScalar(scalar(-1)),
709  calculatedFvPatchScalarField::typeName
710  );
711 
712  const cellZoneMesh& czs = mesh.cellZones();
713  for (const auto& zone : czs)
714  {
716  }
717 
718  cellZone.correctBoundaryConditions();
719  Info<< " Writing cell zoning to " << cellZone.name() << endl;
720  cellZone.write();
721  }
722  if (selectedFields.found("faceZone"))
723  {
724  // Determine for each face the zone index (scalar for ease of
725  // manipulation)
726  scalarField zoneID(mesh.nFaces(), -1);
727  const faceZoneMesh& czs = mesh.faceZones();
728  for (const auto& zone : czs)
729  {
731  }
732 
733 
734  // Split into internal and boundary values
736  (
737  IOobject
738  (
739  "faceZone",
740  mesh.time().timeName(),
741  mesh,
744  false
745  ),
746  mesh,
747  dimensionedScalar(scalar(-1)),
748  calculatedFvsPatchScalarField::typeName
749  );
750 
751  faceZone.primitiveFieldRef() =
753  surfaceScalarField::Boundary& bfld = faceZone.boundaryFieldRef();
754  for (auto& pfld : bfld)
755  {
756  const fvPatch& fvp = pfld.patch();
757  pfld == SubField<scalar>(zoneID, fvp.size(), fvp.start());
758  }
759 
760  //faceZone.correctBoundaryConditions();
761  Info<< " Writing face zoning to " << faceZone.name() << endl;
762  faceZone.write();
763  }
764 
765  Info<< endl;
766 }
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: volSurfaceMapping.H:47
wallDist.H
Foam::cellModel::index
label index() const noexcept
Definition: cellModelI.H:30
volFields.H
Foam::polyMesh::points
virtual const pointField & points() const
Definition: polyMesh.C:1062
Foam::primitiveMeshTools::cellDeterminant
static tmp< scalarField > cellDeterminant(const primitiveMesh &mesh, const Vector< label > &directions, const vectorField &faceAreas, const bitSet &internalOrCoupledFace)
Definition: primitiveMeshTools.C:621
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:165
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::IOobject::AUTO_WRITE
@ AUTO_WRITE
Definition: IOobject.H:190
Foam::syncTools::getInternalOrCoupledFaces
static bitSet getInternalOrCoupledFaces(const polyMesh &mesh)
Definition: syncTools.C:169
Foam::fileName
A class for handling file names.
Definition: fileName.H:71
Foam::Zero
static constexpr const zero Zero
Definition: zero.H:131
cellAspectRatio.H
Foam::dictionary::found
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryI.H:80
Foam::fvPatch::start
virtual label start() const
Definition: fvPatch.H:169
Foam::fvSchemes
Selector class for finite volume differencing schemes. fvMesh is derived from fvSchemes so that all f...
Definition: fvSchemes.H:48
Foam::polyMeshTools::volRatio
static tmp< scalarField > volRatio(const polyMesh &mesh, const scalarField &vol)
Definition: polyMeshTools.C:228
Foam::MeshObject< fvMesh, UpdateableMeshObject, wallDist >::New
static const wallDist & New(const fvMesh &mesh, Args &&... args)
Definition: MeshObject.C:41
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:773
Foam::fvsPatchField
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:64
Foam::zone
Base class for mesh zones.
Definition: zone.H:59
tetPointRef.H
Foam::endl
Ostream & endl(Ostream &os)
Definition: Ostream.H:381
surfaceFields.H
Foam::surfaceFields.
Foam::polyMeshTools::faceSkewness
static tmp< scalarField > faceSkewness(const polyMesh &mesh, const pointField &points, const vectorField &fCtrs, const vectorField &fAreas, const vectorField &cellCtrs)
Definition: polyMeshTools.C:85
Foam::HashSet
A HashTable with keys but without contents that is similar to std::unordered_set.
Definition: HashSet.H:73
Foam::polyMesh::geometricD
const Vector< label > & geometricD() const
Definition: polyMesh.C:851
syncTools.H
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Definition: hashSets.C:26
Foam::tetPointRef
tetrahedron< point, const point & > tetPointRef
A tetrahedron using referred points.
Definition: tetPointRef.H:41
Foam::fvPatch::patchSlice
const List< T >::subList patchSlice(const List< T > &l) const
Definition: fvPatch.H:206
forAll
#define forAll(list, i)
Definition: stdFoam.H:349
Foam::cellZone
A subset of mesh cells.
Definition: cellZone.H:58
Foam::SubField
SubField is a Field obtained as a section of another Field, without its own allocation....
Definition: Field.H:60
Foam::fvPatch::size
virtual label size() const
Definition: fvPatch.H:175
correctBoundaryConditions
cellMask correctBoundaryConditions()
Foam::zone::write
virtual void write(Ostream &os) const
Definition: zone.C:201
Foam::Field
Generic templated field type.
Definition: Field.H:59
Foam::faceZone
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:60
Foam::Info
messageStream Info
Foam::polyMesh::cellZones
const cellZoneMesh & cellZones() const noexcept
Definition: polyMesh.H:488
Foam::polyMesh::faceOwner
virtual const labelList & faceOwner() const
Definition: polyMesh.C:1100
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
regionSplit.H
Foam::dictionary::subDict
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionary.C:453
Foam::ZoneMesh< cellZone, polyMesh >
Foam::polyMesh::faceZones
const faceZoneMesh & faceZones() const noexcept
Definition: polyMesh.H:482
Foam::GeometricField::Boundary::types
wordList types() const
Definition: GeometricBoundaryField.C:463
Foam::radToDeg
constexpr scalar radToDeg(const scalar rad) noexcept
Definition: unitConversion.H:52
Foam::regionSplit
This class separates the mesh into distinct unconnected regions, each of which is then given a label ...
Definition: regionSplit.H:136
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Definition: hashSets.C:40
fld
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))
Definition: gmvOutputLagrangian.H:23
writeFields.H
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:119
Foam::fvMesh::neighbour
const labelUList & neighbour() const
Definition: fvMesh.H:409
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:36
zoneID
const labelIOList & zoneID
Definition: interpolatedFaces.H:22
Foam::Ostream::write
virtual bool write(const token &tok)=0
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:81
Foam
Definition: atmBoundaryLayer.C:26
Foam::primitiveMesh::cellVolumes
const scalarField & cellVolumes() const
Definition: primitiveMeshCellCentresAndVols.C:90
Foam::cellAspectRatio
(Rough approximation of) cell aspect ratio
Definition: cellAspectRatio.H:47
Foam::primitiveMesh::cellShapes
const cellShapeList & cellShapes() const
Definition: primitiveMesh.C:346
Foam::fvMesh::owner
const labelUList & owner() const
Definition: fvMesh.H:403
Foam::GeometricField::correctBoundaryConditions
void correctBoundaryConditions()
Definition: GeometricField.C:933
Foam::zoneIdentifier::index
label index() const noexcept
Definition: zoneIdentifier.H:131
Foam::fvMesh::boundary
const fvBoundaryMesh & boundary() const
Definition: fvMesh.C:678
polyMeshTools.H
Foam::primitiveMesh::cellCentres
const vectorField & cellCentres() const
Definition: primitiveMeshCellCentresAndVols.C:78
Foam::primitiveMesh::isInternalFace
bool isInternalFace(const label faceIndex) const noexcept
Definition: primitiveMeshI.H:96
Foam::polyMesh::faces
virtual const faceList & faces() const
Definition: polyMesh.C:1087
Foam::GeometricField::Boundary
Definition: GeometricField.H:111
Foam::fvsPatchField::patch
const fvPatch & patch() const
Definition: fvsPatchField.H:277
Foam::zoneIdentifier::name
const word & name() const noexcept
Definition: zoneIdentifier.H:119
f
labelList f(nPoints)
Foam::GeometricField::ref
Internal & ref(const bool updateAccessTime=true)
Definition: GeometricField.C:742
Foam::GeometricField::boundaryFieldRef
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Definition: GeometricField.C:776
Foam::Vector< scalar >
Foam::primitiveMesh::nInternalFaces
label nInternalFaces() const noexcept
Definition: primitiveMeshI.H:71
Foam::List< cell >
Foam::primitiveMeshTools::cellClosedness
static void cellClosedness(const primitiveMesh &mesh, const Vector< label > &meshD, const vectorField &areas, const scalarField &vols, scalarField &openness, scalarField &aratio)
Definition: primitiveMeshTools.C:407
Foam::acos
dimensionedScalar acos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:261
Foam::surfaceInterpolation::weights
virtual const surfaceScalarField & weights() const
Definition: surfaceInterpolation.C:96
Foam::primitiveMesh::faceCentres
const vectorField & faceCentres() const
Definition: primitiveMeshFaceCentresAndAreas.C:71
Foam::UList
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:99
Foam::fvPatch::faceCells
virtual const labelUList & faceCells() const
Definition: fvPatch.C:106
Foam::fvMesh::time
const Time & time() const
Definition: fvMesh.H:276
Foam::UIndirectList
A List with indirect addressing.
Definition: faMatrix.H:56
Foam::primitiveMesh::nFaces
label nFaces() const noexcept
Definition: primitiveMeshI.H:83
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:68
cellShapes
cellShapeList cellShapes
Definition: createBlockMesh.H:3
Foam::faceZone::write
virtual void write(Ostream &os) const
Definition: faceZone.C:609
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
Foam::cellModel
Maps a geometry to a set of cell primitives.
Definition: cellModel.H:68
Foam::fvPatchField::patch
const fvPatch & patch() const
Definition: fvPatchField.H:353
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
Foam::GeometricField
Generic GeometricField class.
Definition: areaFieldsFwd.H:49
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:184
zeroGradientFvPatchFields.H
Foam::cell
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:50
Foam::polyMesh::faceNeighbour
virtual const labelList & faceNeighbour() const
Definition: polyMesh.C:1106
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Definition: GeometricFieldI.H:55
Foam::writeFields
void writeFields(const fvMesh &mesh, const wordHashSet &selectedFields, const bool writeFaceFields)
Foam::dimless
const dimensionSet dimless
Definition: dimensionSets.C:182
y
scalar y
Definition: LISASMDCalcMethod1.H:14
Foam::fvMesh::V
const DimensionedField< scalar, volMesh > & V() const
Definition: fvMeshGeometry.C:172
Foam::primitiveMeshTools::facePyramidVolume
static void facePyramidVolume(const primitiveMesh &mesh, const pointField &points, const vectorField &cellCtrs, scalarField &ownPyrVol, scalarField &neiPyrVol)
Definition: primitiveMeshTools.C:368
Foam::polyMeshTools::faceOrthogonality
static tmp< scalarField > faceOrthogonality(const polyMesh &mesh, const vectorField &fAreas, const vectorField &cellCtrs)
Definition: polyMeshTools.C:30
Foam::primitiveMesh::faceAreas
const vectorField & faceAreas() const
Definition: primitiveMeshFaceCentresAndAreas.C:83