Loading...
Searching...
No Matches
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
11using namespace Foam;
12
13void 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 for (const label facei : cells[celli])
27 {
28 cellFld[celli] = Foam::max(cellFld[celli], faceData[facei]);
29 }
30 }
31
32 for (fvPatchScalarField& fvp : cellData.boundaryFieldRef())
33 {
34 fvp = fvp.patch().patchSlice(faceData);
35 }
37}
38
39
40void minFaceToCell
41(
42 const scalarField& faceData,
43 volScalarField& cellData
44)
45{
46 const cellList& cells = cellData.mesh().cells();
47
48 scalarField& cellFld = cellData.ref();
49
50 cellFld = GREAT;
51 forAll(cells, celli)
52 {
53 for (const label facei : cells[celli])
54 {
55 cellFld[celli] = Foam::min(cellFld[celli], faceData[facei]);
56 }
57 }
58
59 for (fvPatchScalarField& fvp : cellData.boundaryFieldRef())
60 {
61 fvp = fvp.patch().patchSlice(faceData);
62 }
64}
65
66
67void minFaceToCell
68(
69 const surfaceScalarField& faceData,
70 volScalarField& cellData,
72)
73{
74 scalarField& cellFld = cellData.primitiveFieldRef();
75
76 cellFld = GREAT;
77
78 const labelUList& own = cellData.mesh().owner();
79 const labelUList& nei = cellData.mesh().neighbour();
80
81 // Internal faces
82 forAll(own, facei)
83 {
84 cellFld[own[facei]] = Foam::min(cellFld[own[facei]], faceData[facei]);
85 cellFld[nei[facei]] = Foam::min(cellFld[nei[facei]], faceData[facei]);
86 }
87
88 // Patch faces
89 for (const fvsPatchScalarField& fvp : faceData.boundaryField())
90 {
91 label pfacei = 0;
92
93 for (const label celli : fvp.patch().faceCells())
94 {
95 cellFld[celli] = Foam::max(cellFld[celli], fvp[pfacei]);
96 ++pfacei;
97 }
98 }
99
101
102 forAll(bfld, patchi)
103 {
104 bfld[patchi] = faceData.boundaryField()[patchi];
105 }
107 {
108 cellData.correctBoundaryConditions();
109 }
110}
111
112
113void writeSurfaceField
114(
115 const fvMesh& mesh,
116 const fileName& fName,
117 const scalarField& faceData
118)
119{
120 // Write single surfaceScalarField
121
123 (
125 (
126 fName,
127 mesh.time().timeName(),
128 mesh,
132 ),
133 mesh,
135 );
136 fld.primitiveFieldRef() = faceData;
137 //fld.correctBoundaryConditions();
138 Info<< " Writing face data to " << fName << endl;
139 fld.write();
140}
141
142
144(
145 const fvMesh& mesh,
146 const wordHashSet& selectedFields,
147 const bool writeFaceFields
148)
149{
150 if (selectedFields.empty())
151 {
152 return;
153 }
154
155 Info<< "Writing fields with mesh quality parameters" << endl;
156
157 if (selectedFields.found("nonOrthoAngle"))
158 {
159 //- Face based orthogonality
160 scalarField faceOrthogonality
161 (
163 (
164 mesh,
165 mesh.faceAreas(),
166 mesh.cellCentres()
167 )
168 );
169 faceOrthogonality.clamp_range(-1, 1);
170
171 //- Face based angle
172 const scalarField nonOrthoAngle
173 (
174 radToDeg(Foam::acos(faceOrthogonality))
175 );
176
177 //- Cell field - max of either face
178 volScalarField cellNonOrthoAngle
179 (
181 (
182 "nonOrthoAngle",
183 mesh.time().timeName(),
184 mesh,
188 ),
189 mesh,
191 );
192 //- Take max
193 maxFaceToCell(nonOrthoAngle, cellNonOrthoAngle);
194 Info<< " Writing non-orthogonality (angle) to "
195 << cellNonOrthoAngle.name() << endl;
196 cellNonOrthoAngle.write();
197
198 if (writeFaceFields)
199 {
200 writeSurfaceField
201 (
202 mesh,
203 "face_nonOrthoAngle",
204 SubField<scalar>(nonOrthoAngle, mesh.nInternalFaces())
205 );
206 }
207 }
208
209 if (selectedFields.found("faceWeight"))
210 {
211 volScalarField cellWeights
212 (
214 (
215 "faceWeight",
216 mesh.time().timeName(),
217 mesh,
221 ),
222 mesh,
224 wordList // wanted bc types
225 (
226 mesh.boundary().size(),
228 ),
229 mesh.weights().boundaryField().types() // current bc types
230 );
231 //- Take min
232 minFaceToCell(mesh.weights(), cellWeights, false);
233 Info<< " Writing face interpolation weights (0..0.5) to "
234 << cellWeights.name() << endl;
235 cellWeights.write();
236
237 if (writeFaceFields)
238 {
239 writeSurfaceField(mesh, "face_faceWeight", mesh.weights());
240 }
241 }
242
243
244 // Skewness
245 // ~~~~~~~~
246
247 if (selectedFields.found("skewness"))
248 {
249 //- Face based skewness
250 const scalarField faceSkewness
251 (
253 (
254 mesh,
255 mesh.points(),
256 mesh.faceCentres(),
257 mesh.faceAreas(),
258 mesh.cellCentres()
259 )
260 );
261
262 //- Cell field - max of either face
263 volScalarField cellSkewness
264 (
266 (
267 "skewness",
268 mesh.time().timeName(),
269 mesh,
273 ),
274 mesh,
276 );
277 //- Take max
278 maxFaceToCell(faceSkewness, cellSkewness);
279 Info<< " Writing face skewness to " << cellSkewness.name() << endl;
280 cellSkewness.write();
281
282 if (writeFaceFields)
283 {
284 writeSurfaceField
285 (
286 mesh,
287 "face_skewness",
288 SubField<scalar>(faceSkewness, mesh.nInternalFaces())
289 );
290 }
291 }
292
293
294 // cellDeterminant
295 // ~~~~~~~~~~~~~~~
296
297 if (selectedFields.found("cellDeterminant"))
298 {
299 volScalarField cellDeterminant
300 (
302 (
303 "cellDeterminant",
304 mesh.time().timeName(),
305 mesh,
309 ),
310 mesh,
313 );
314 cellDeterminant.primitiveFieldRef() =
316 (
317 mesh,
318 mesh.geometricD(),
319 mesh.faceAreas(),
321 );
322 cellDeterminant.correctBoundaryConditions();
323 Info<< " Writing cell determinant to "
324 << cellDeterminant.name() << endl;
325 cellDeterminant.write();
326 }
327
328
329 // Aspect ratio
330 // ~~~~~~~~~~~~
331
332 if (selectedFields.found("aspectRatio"))
333 {
334 volScalarField aspectRatio
335 (
337 (
338 "aspectRatio",
339 mesh.time().timeName(),
340 mesh,
344 ),
345 mesh,
348 );
349
350
351 scalarField cellOpenness;
353 (
354 mesh,
355 mesh.geometricD(),
356 mesh.faceAreas(),
357 mesh.cellVolumes(),
358 cellOpenness,
359 aspectRatio.ref()
360 );
361
362 aspectRatio.correctBoundaryConditions();
363 Info<< " Writing aspect ratio to " << aspectRatio.name() << endl;
364 aspectRatio.write();
365 }
366
367 if (selectedFields.found("cellAspectRatio"))
368 {
369 volScalarField aspectRatio
370 (
372 (
373 "cellAspectRatio",
374 mesh.time().timeName(),
375 mesh,
379 ),
380 mesh,
383 );
384
385 aspectRatio.ref().field() = cellAspectRatio(mesh);
386
387 aspectRatio.correctBoundaryConditions();
388 Info<< " Writing approximate aspect ratio to "
389 << aspectRatio.name() << endl;
390 aspectRatio.write();
391 }
392
393
394 // cell type
395 // ~~~~~~~~~
396
397 if (selectedFields.found("cellShapes"))
398 {
399 volScalarField shape
400 (
402 (
403 "cellShapes",
404 mesh.time().timeName(),
405 mesh,
409 ),
410 mesh,
413 );
414 const cellShapeList& cellShapes = mesh.cellShapes();
415 forAll(cellShapes, cellI)
416 {
417 const cellModel& model = cellShapes[cellI].model();
418 shape[cellI] = model.index();
419 }
420 shape.correctBoundaryConditions();
421 Info<< " Writing cell shape (hex, tet etc.) to " << shape.name()
422 << endl;
423 shape.write();
424 }
425
426 if (selectedFields.found("cellVolume"))
427 {
429 (
431 (
432 "cellVolume",
433 mesh.time().timeName(),
434 mesh,
438 ),
439 mesh,
441 );
442 V.ref() = mesh.V();
443 Info<< " Writing cell volume to " << V.name() << endl;
444 V.write();
445 }
446
447 if (selectedFields.found("cellVolumeRatio"))
448 {
449 const scalarField faceVolumeRatio
450 (
452 (
453 mesh,
454 mesh.V()
455 )
456 );
457
458 volScalarField cellVolumeRatio
459 (
461 (
462 "cellVolumeRatio",
463 mesh.time().timeName(),
464 mesh,
468 ),
469 mesh,
471 );
472 //- Take min
473 minFaceToCell(faceVolumeRatio, cellVolumeRatio);
474 Info<< " Writing cell volume ratio to "
475 << cellVolumeRatio.name() << endl;
476 cellVolumeRatio.write();
477
478 if (writeFaceFields)
479 {
480 writeSurfaceField
481 (
482 mesh,
483 "face_cellVolumeRatio",
484 SubField<scalar>(faceVolumeRatio, mesh.nInternalFaces())
485 );
486 }
487 }
488
489 // minTetVolume
490 if (selectedFields.found("minTetVolume"))
491 {
492 volScalarField minTetVolume
493 (
495 (
496 "minTetVolume",
497 mesh.time().timeName(),
498 mesh,
502 ),
503 mesh,
504 dimensionedScalar("minTetVolume", dimless, GREAT),
506 );
507
508
509 const labelList& own = mesh.faceOwner();
510 const labelList& nei = mesh.faceNeighbour();
511 const pointField& p = mesh.points();
512 forAll(own, facei)
513 {
514 const face& f = mesh.faces()[facei];
515 const point& fc = mesh.faceCentres()[facei];
516
517 {
518 const point& ownCc = mesh.cellCentres()[own[facei]];
519 scalar& ownVol = minTetVolume[own[facei]];
520 forAll(f, fp)
521 {
522 scalar tetQual = tetPointRef
523 (
524 p[f[fp]],
525 p[f.nextLabel(fp)],
526 ownCc,
527 fc
528 ).quality();
529 ownVol = Foam::min(ownVol, tetQual);
530 }
531 }
532 if (mesh.isInternalFace(facei))
533 {
534 const point& neiCc = mesh.cellCentres()[nei[facei]];
535 scalar& neiVol = minTetVolume[nei[facei]];
536 forAll(f, fp)
537 {
538 scalar tetQual = tetPointRef
539 (
540 p[f[fp]],
541 p[f.nextLabel(fp)],
542 fc,
543 neiCc
544 ).quality();
545 neiVol = Foam::min(neiVol, tetQual);
546 }
547 }
548 }
549
550 minTetVolume.correctBoundaryConditions();
551 Info<< " Writing minTetVolume to " << minTetVolume.name() << endl;
552 minTetVolume.write();
553 }
554
555 // minPyrVolume
556 if (selectedFields.found("minPyrVolume"))
557 {
558 volScalarField minPyrVolume
559 (
561 (
562 "minPyrVolume",
563 mesh.time().timeName(),
564 mesh,
568 ),
569 mesh,
570 dimensionedScalar("minPyrVolume", dimless, GREAT),
572 );
573
574 // Get owner and neighbour pyr volumes
575 scalarField ownPyrVol(mesh.nFaces());
576 scalarField neiPyrVol(mesh.nInternalFaces());
578 (
579 mesh,
580 mesh.points(),
581 mesh.cellCentres(),
582
583 ownPyrVol,
584 neiPyrVol
585 );
586
587 // Get min pyr vol per cell
588 scalarField& cellFld = minPyrVolume.ref();
589 cellFld = GREAT;
590
591 const labelUList& own = mesh.owner();
592 const labelUList& nei = mesh.neighbour();
593
594 // Internal faces
595 forAll(own, facei)
596 {
597 cellFld[own[facei]] =
598 Foam::min(cellFld[own[facei]], ownPyrVol[facei]);
599 cellFld[nei[facei]] =
600 Foam::min(cellFld[nei[facei]], neiPyrVol[facei]);
601 }
602
603 // Patch faces
604 for (const auto& fvp : minPyrVolume.boundaryField())
605 {
606 label meshFacei = fvp.patch().start();
607
608 for (const label celli : fvp.patch().faceCells())
609 {
610 cellFld[celli] =
611 Foam::min(cellFld[celli], ownPyrVol[meshFacei]);
612
613 ++meshFacei;
614 }
615 }
616
617 minPyrVolume.correctBoundaryConditions();
618 Info<< " Writing minPyrVolume to " << minPyrVolume.name() << endl;
619 minPyrVolume.write();
620
621 if (writeFaceFields)
622 {
623 scalarField minFacePyrVol(neiPyrVol);
624 minFacePyrVol = Foam::min
625 (
626 minFacePyrVol,
627 SubField<scalar>(ownPyrVol, mesh.nInternalFaces())
628 );
629 writeSurfaceField(mesh, "face_minPyrVolume", minFacePyrVol);
630 }
631 }
632
633 if (selectedFields.found("cellRegion"))
634 {
635 volScalarField cellRegion
636 (
638 (
639 "cellRegion",
640 mesh.time().timeName(),
641 mesh,
645 ),
646 mesh,
648 );
649
650 regionSplit rs(mesh);
651 forAll(rs, celli)
652 {
653 cellRegion[celli] = rs[celli];
654 }
655 cellRegion.correctBoundaryConditions();
656 Info<< " Writing cell region to " << cellRegion.name() << endl;
657 cellRegion.write();
658 }
659
660 if (selectedFields.found("wallDistance"))
661 {
662 // See if wallDist.method entry in fvSchemes before calling factory
663 // method of wallDist. Have 'failing' version of wallDist::New instead?
664 const dictionary& schemesDict =
665 static_cast<const fvSchemes&>(mesh).schemesDict();
666 if (schemesDict.found("wallDist"))
667 {
668 if (schemesDict.subDict("wallDist").found("method"))
669 {
670 // Wall distance
671 volScalarField y("wallDistance", wallDist::New(mesh).y());
672 Info<< " Writing wall distance to " << y.name() << endl;
673 y.write();
674
675 // Wall-reflection vectors
676 //const volVectorField& n = wallDist::New(mesh).n();
677 //Info<< " Writing wall normal to " << n.name() << endl;
678 //n.write();
679 }
680 }
681 }
682
683 if (selectedFields.found("cellZone"))
684 {
686 (
688 (
689 "cellZone",
690 mesh.time().timeName(),
691 mesh,
695 ),
696 mesh,
698 );
699
700 const cellZoneMesh& czs = mesh.cellZones();
701 for (const auto& zone : czs)
702 {
704 }
705
706 cellZone.correctBoundaryConditions();
707 Info<< " Writing cell zoning to " << cellZone.name() << endl;
708 cellZone.write();
709 }
710 if (selectedFields.found("faceZone"))
711 {
712 // Determine for each face the zone index (scalar for ease of
713 // manipulation)
714 scalarField zoneID(mesh.nFaces(), -1);
715 const faceZoneMesh& czs = mesh.faceZones();
716 for (const auto& zone : czs)
717 {
719 }
720
721
722 // Split into internal and boundary values
724 (
726 (
727 "faceZone",
728 mesh.time().timeName(),
729 mesh,
733 ),
734 mesh,
736 );
737
738 faceZone.primitiveFieldRef() =
739 SubField<scalar>(zoneID, mesh.nInternalFaces());
740 surfaceScalarField::Boundary& bfld = faceZone.boundaryFieldRef();
741 for (auto& pfld : bfld)
742 {
743 const fvPatch& fvp = pfld.patch();
744 pfld == SubField<scalar>(zoneID, fvp.size(), fvp.start());
745 }
746
747 //faceZone.correctBoundaryConditions();
748 Info<< " Writing face zoning to " << faceZone.name() << endl;
749 faceZone.write();
750 }
751
752 Info<< endl;
753}
scalar y
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
const Mesh & mesh() const noexcept
Return const reference to mesh.
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef().
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
void correctBoundaryConditions()
Correct boundary field.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
bool empty() const noexcept
True if the hash table is empty.
Definition HashTable.H:353
bool found(const Key &key) const
Same as contains().
Definition HashTable.H:1370
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
static FOAM_NO_DANGLING_REFERENCE const wallDist & New(const fvMesh &mesh, Args &&... args)
SubField is a Field obtained as a section of another Field, without its own allocation....
Definition SubField.H:58
A List with indirect addressing. Like IndirectList but does not store addressing.
(Rough approximation of) cell aspect ratio
Maps a geometry to a set of cell primitives.
Definition cellModel.H:73
label index() const noexcept
Return index of model in the model list.
Definition cellModelI.H:30
A subset of mesh cells.
Definition cellZone.H:61
virtual void write(Ostream &os) const
Write (dictionary entries).
Definition cellZone.C:217
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition dictionary.C:441
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry (const access) with the given keyword.
A subset of mesh faces organised as a primitive patch.
Definition faceZone.H:63
virtual void write(Ostream &os) const
Write (dictionary entries).
Definition faceZone.C:695
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
A class for handling file names.
Definition fileName.H:75
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
static const word & calculatedType() noexcept
The type name for calculated patch fields.
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
virtual label size() const
Patch size is the number of faces, but can be overloaded.
Definition fvPatch.H:242
const polyPatch & patch() const noexcept
Return the polyPatch.
Definition fvPatch.H:202
label start() const noexcept
The patch start within the polyMesh face list.
Definition fvPatch.H:226
Selector class for finite volume differencing schemes. fvMesh is derived from fvSchemes so that all f...
Definition fvSchemes.H:54
static tmp< scalarField > faceOrthogonality(const polyMesh &mesh, const vectorField &fAreas, const vectorField &cellCtrs)
Generate orthogonality field. (1 for fully orthogonal, < 1 for non-orthogonal).
static tmp< scalarField > faceSkewness(const polyMesh &mesh, const pointField &points, const vectorField &fCtrs, const vectorField &fAreas, const vectorField &cellCtrs)
Generate skewness field.
static tmp< scalarField > volRatio(const polyMesh &mesh, const scalarField &vol)
Generate volume ratio field.
label start() const noexcept
Return start label of this patch in the polyMesh face list.
Definition polyPatch.H:446
const List< T >::subList patchSlice(const UList< T > &values) const
This patch slice from the complete list, which has size mesh::nFaces(), using the number of patch fac...
Definition polyPatch.H:502
const labelUList & faceCells() const
Return face-cell addressing.
Definition polyPatch.C:401
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.
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.
static void facePyramidVolume(const primitiveMesh &mesh, const pointField &points, const vectorField &cellCtrs, scalarField &ownPyrVol, scalarField &neiPyrVol)
Generate face pyramid volume fields.
This class separates the mesh into distinct unconnected regions, each of which is then given a label ...
static bitSet getInternalOrCoupledFaces(const polyMesh &mesh)
Get per face whether it is internal or coupled.
Definition syncTools.C:165
static const word null
An empty word.
Definition word.H:84
label index() const noexcept
The index of this zone in the zone list.
const word & name() const noexcept
The zone name.
Base class for mesh zones.
Definition zone.H:63
volScalarField & p
cellShapeList cellShapes
dynamicFvMesh & mesh
const cellShapeList & cells
const expr V(m.psi().mesh().V())
Namespace for OpenFOAM.
List< word > wordList
List of word.
Definition fileName.H:60
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
const dimensionSet dimless
Dimensionless.
HashSet< word, Hash< word > > wordHashSet
A HashSet of words, uses string hasher.
Definition HashSet.H:80
List< label > labelList
A List of labels.
Definition List.H:62
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
ZoneMesh< faceZone, polyMesh > faceZoneMesh
A ZoneMesh with faceZone content on a polyMesh.
tetrahedron< point, const point & > tetPointRef
A tetrahedron using referred points.
Definition tetrahedron.H:72
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
List< cell > cellList
List of cell.
Definition cellListFwd.H:41
ZoneMesh< cellZone, polyMesh > cellZoneMesh
A ZoneMesh with cellZone content on a polyMesh.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
vector point
Point is a vector.
Definition point.H:37
void writeFields(const fvMesh &mesh, const wordHashSet &selectedFields, const bool writeFaceFields)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
const dimensionSet dimVolume(pow3(dimLength))
vectorField pointField
pointField is a vectorField.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
UList< label > labelUList
A UList of labels.
Definition UList.H:75
fvsPatchField< scalar > fvsPatchScalarField
List< cellShape > cellShapeList
List of cellShape.
dimensionedScalar acos(const dimensionedScalar &ds)
fvPatchField< scalar > fvPatchScalarField
constexpr scalar radToDeg(const scalar rad) noexcept
Conversion from radians to degrees.
labelList f(nPoints)
cellMask correctBoundaryConditions()
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Foam::surfaceFields.