faMesh.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2016-2017 Wikki Ltd
9 Copyright (C) 2020-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "faMesh.H"
30#include "faMeshBoundaryHalo.H"
31#include "faGlobalMeshData.H"
32#include "Time.H"
33#include "polyMesh.H"
34#include "primitiveMesh.H"
35#include "demandDrivenData.H"
36#include "IndirectList.H"
37#include "areaFields.H"
38#include "edgeFields.H"
39#include "faMeshLduAddressing.H"
40#include "processorFaPatch.H"
41#include "wedgeFaPatch.H"
42#include "faPatchData.H"
43
44// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
45
46namespace Foam
47{
49}
50
51
52const Foam::word Foam::faMesh::prefix("finite-area");
53
55
56int Foam::faMesh::geometryOrder_ = 1; // 1: Standard treatment
57
58const int Foam::faMesh::quadricsFit_ = 0; // Tuning (experimental)
59
60
61// * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
62
63namespace Foam
64{
65
66// Convert patch names to face labels. Preserve patch order
68(
69 const polyBoundaryMesh& pbm,
70 const wordRes& polyPatchNames
71)
72{
73 const labelList patchIDs
74 (
75 pbm.patchSet
76 (
77 polyPatchNames,
78 false, // warnNotFound
79 true // useGroups
80 ).sortedToc()
81 );
82
83 if (patchIDs.empty())
84 {
86 << "No matching patches: " << polyPatchNames << nl
87 << exit(FatalError);
88 }
89
90 label nFaceLabels = 0;
91 for (const label patchi : patchIDs)
92 {
93 nFaceLabels += pbm[patchi].size();
94 }
95
96 labelList faceLabels(nFaceLabels);
97
98 nFaceLabels = 0;
99 for (const label patchi : patchIDs)
100 {
101 for (const label facei : pbm[patchi].range())
102 {
103 faceLabels[nFaceLabels] = facei;
104 ++nFaceLabels;
105 }
106 }
107
108 return faceLabels;
109}
110
111} // End namespace Foam
112
113
114// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
115
116void Foam::faMesh::checkBoundaryEdgeLabelRange
117(
118 const labelUList& edgeLabels
119) const
120{
121 label nErrors = 0;
122
123 for (const label edgei : edgeLabels)
124 {
125 if (edgei < nInternalEdges_ || edgei >= nEdges_)
126 {
127 if (!nErrors++)
128 {
130 << "Boundary edge label out of range "
131 << nInternalEdges_ << ".." << (nEdges_-1) << nl
132 << " ";
133 }
134
135 FatalError<< ' ' << edgei;
136 }
137 }
138
139 if (nErrors)
140 {
142 }
143}
144
145
146void Foam::faMesh::initPatch() const
147{
148 patchPtr_.reset
149 (
151 (
152 UIndirectList<face>(mesh().faces(), faceLabels_),
153 mesh().points()
154 )
155 );
156 bndConnectPtr_.reset(nullptr);
157 haloMapPtr_.reset(nullptr);
158 haloFaceCentresPtr_.reset(nullptr);
159 haloFaceNormalsPtr_.reset(nullptr);
160}
161
162
163void Foam::faMesh::setPrimitiveMeshData()
164{
165 DebugInFunction << "Setting primitive data" << endl;
166
167 const uindirectPrimitivePatch& bp = patch();
168 const labelListList& edgeFaces = bp.edgeFaces();
169
170 // Dimensions
171
172 nEdges_ = bp.nEdges();
173 nInternalEdges_ = bp.nInternalEdges();
174 nFaces_ = bp.size();
175 nPoints_ = bp.nPoints();
176
177 edges_.resize(nEdges_);
178 edgeOwner_.resize(nEdges_);
179 edgeNeighbour_.resize(nInternalEdges_);
180
181 // Internal edges
182 for (label edgei = 0; edgei < nInternalEdges_; ++edgei)
183 {
184 edges_[edgei] = bp.edges()[edgei];
185
186 edgeOwner_[edgei] = edgeFaces[edgei][0];
187
188 edgeNeighbour_[edgei] = edgeFaces[edgei][1];
189 }
190
191 // Continue with boundary edges
192 label edgei = nInternalEdges_;
193
194 for (const faPatch& p : boundary())
195 {
196 for (const label patchEdgei : p.edgeLabels())
197 {
198 edges_[edgei] = bp.edges()[patchEdgei];
199
200 edgeOwner_[edgei] = edgeFaces[patchEdgei][0];
201
202 ++edgei;
203 }
204 }
205}
206
207
208void Foam::faMesh::clearHalo() const
209{
210 DebugInFunction << "Clearing halo information" << endl;
211
212 haloMapPtr_.reset(nullptr);
213 haloFaceCentresPtr_.reset(nullptr);
214 haloFaceNormalsPtr_.reset(nullptr);
215}
216
217
218void Foam::faMesh::clearGeomNotAreas() const
219{
220 DebugInFunction << "Clearing geometry" << endl;
221
222 clearHalo();
223 patchPtr_.reset(nullptr);
224 bndConnectPtr_.reset(nullptr);
226 deleteDemandDrivenData(patchStartsPtr_);
228 deleteDemandDrivenData(magLePtr_);
229 deleteDemandDrivenData(centresPtr_);
230 deleteDemandDrivenData(edgeCentresPtr_);
231 deleteDemandDrivenData(faceAreaNormalsPtr_);
232 deleteDemandDrivenData(edgeAreaNormalsPtr_);
233 pointAreaNormalsPtr_.reset(nullptr);
234 deleteDemandDrivenData(faceCurvaturesPtr_);
235 deleteDemandDrivenData(edgeTransformTensorsPtr_);
236}
237
238
239void Foam::faMesh::clearGeom() const
240{
241 DebugInFunction << "Clearing geometry" << endl;
242
243 clearGeomNotAreas();
245 deleteDemandDrivenData(S00Ptr_);
246 deleteDemandDrivenData(correctPatchPointNormalsPtr_);
247}
248
249
250void Foam::faMesh::clearAddressing() const
251{
252 DebugInFunction << "Clearing addressing" << endl;
253
254 deleteDemandDrivenData(lduPtr_);
255}
256
257
258void Foam::faMesh::clearOut() const
259{
260 clearGeom();
261 clearAddressing();
262 globalMeshDataPtr_.reset(nullptr);
263}
264
265
266// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
267
268bool Foam::faMesh::init(const bool doInit)
269{
270 if (doInit)
271 {
272 setPrimitiveMeshData();
273 }
274
275 // Create global mesh data
276 if (Pstream::parRun())
277 {
278 (void)globalData();
279 }
280
281 // Calculate topology for the patches (processor-processor comms etc.)
282 boundary_.updateMesh();
283
284 // Calculate the geometry for the patches (transformation tensors etc.)
285 boundary_.calcGeometry();
286
287 // Ensure processor/processor information is properly synchronised
288 if (Pstream::parRun())
289 {
290 const_cast<areaVectorField&>(areaCentres()).boundaryFieldRef()
291 .evaluateCoupled<processorFaPatch>();
292
293 // This roughly corresponds to what OpenFOAM-v2112 (and earlier) had,
294 // but should nominally be unnecessary.
295 //
298 }
299
300 return false;
301}
302
303
305:
306 faMesh(pMesh, labelList(), static_cast<const IOobject&>(pMesh))
307{}
308
309
311:
312 faMesh(baseMesh, labelList())
313{}
314
315
317(
318 const polyMesh& pMesh,
319 const bool doInit
320)
321:
323 faSchemes(mesh()),
324 edgeInterpolation(*this),
325 faSolution(mesh()),
326 data(mesh()), // Always NO_READ, NO_WRITE
327 faceLabels_
328 (
330 (
331 "faceLabels",
332 time().findInstance(meshDir(), "faceLabels"),
333 faMesh::meshSubDir,
334 mesh(),
335 IOobject::MUST_READ,
336 IOobject::NO_WRITE
337 )
338 ),
339 boundary_
340 (
342 (
343 "faBoundary",
344 // Allow boundary file that is newer than faceLabels
345 time().findInstance
346 (
347 meshDir(),
348 "faBoundary",
349 IOobject::MUST_READ,
350 faceLabels_.instance()
351 ),
352 faMesh::meshSubDir,
353 mesh(),
354 IOobject::MUST_READ,
355 IOobject::NO_WRITE
356 ),
357 *this
358 ),
359 comm_(Pstream::worldComm),
360 curTimeIndex_(time().timeIndex()),
361
362 patchPtr_(nullptr),
363 bndConnectPtr_(nullptr),
364 lduPtr_(nullptr),
365
366 SPtr_(nullptr),
367 S0Ptr_(nullptr),
368 S00Ptr_(nullptr),
369 patchStartsPtr_(nullptr),
370 LePtr_(nullptr),
371 magLePtr_(nullptr),
372 centresPtr_(nullptr),
373 edgeCentresPtr_(nullptr),
374 faceAreaNormalsPtr_(nullptr),
375 edgeAreaNormalsPtr_(nullptr),
376 pointAreaNormalsPtr_(nullptr),
377 faceCurvaturesPtr_(nullptr),
378 edgeTransformTensorsPtr_(nullptr),
379 correctPatchPointNormalsPtr_(nullptr),
380 globalMeshDataPtr_(nullptr),
381
382 haloMapPtr_(nullptr),
383 haloFaceCentresPtr_(nullptr),
384 haloFaceNormalsPtr_(nullptr)
385{
386 DebugInFunction << "Creating from IOobject" << endl;
387
388 setPrimitiveMeshData();
389
390 if (doInit)
391 {
392 faMesh::init(false); // do not init lower levels
393 }
394
395 if (doInit && fileHandler().isFile(pMesh.time().timePath()/"S0"))
396 {
398 (
400 (
401 "S0",
402 time().timeName(),
404 mesh(),
407 ),
408 *this
409 );
410 }
411}
412
413
415(
416 const polyMesh& pMesh,
417 labelList&& faceLabels
418)
419:
420 faMesh(pMesh, std::move(faceLabels), static_cast<const IOobject&>(pMesh))
421{}
422
423
425(
426 const polyMesh& pMesh,
427 labelList&& faceLabels,
428 const IOobject& io
429)
430:
432 faSchemes(mesh(), io.readOpt()),
433 edgeInterpolation(*this),
434 faSolution(mesh(), io.readOpt()),
435 data(mesh()), // Always NO_READ, NO_WRITE
436 faceLabels_
437 (
439 (
440 "faceLabels",
441 mesh().facesInstance(),
442 faMesh::meshSubDir,
443 mesh(),
444 IOobject::NO_READ,
445 IOobject::NO_WRITE
446 ),
447 std::move(faceLabels)
448 ),
449 boundary_
450 (
452 (
453 "faBoundary",
454 mesh().facesInstance(),
455 faMesh::meshSubDir,
456 mesh(),
457 IOobject::NO_READ,
458 IOobject::NO_WRITE
459 ),
460 *this,
461 label(0)
462 ),
463 comm_(Pstream::worldComm),
464 curTimeIndex_(time().timeIndex()),
465
466 patchPtr_(nullptr),
467 bndConnectPtr_(nullptr),
468 lduPtr_(nullptr),
469
470 SPtr_(nullptr),
471 S0Ptr_(nullptr),
472 S00Ptr_(nullptr),
473 patchStartsPtr_(nullptr),
474 LePtr_(nullptr),
475 magLePtr_(nullptr),
476 centresPtr_(nullptr),
477 edgeCentresPtr_(nullptr),
478 faceAreaNormalsPtr_(nullptr),
479 edgeAreaNormalsPtr_(nullptr),
480 pointAreaNormalsPtr_(nullptr),
481 faceCurvaturesPtr_(nullptr),
482 edgeTransformTensorsPtr_(nullptr),
483 correctPatchPointNormalsPtr_(nullptr),
484 globalMeshDataPtr_(nullptr),
485
486 haloMapPtr_(nullptr),
487 haloFaceCentresPtr_(nullptr),
488 haloFaceNormalsPtr_(nullptr)
489{}
490
491
493(
494 const faMesh& baseMesh,
495 labelList&& faceLabels
496)
497:
500 (
501 mesh(),
502 static_cast<const faSchemes&>(baseMesh)
503 ),
504 edgeInterpolation(*this),
506 (
507 mesh(),
508 static_cast<const faSolution&>(baseMesh)
509 ),
510 data
511 (
512 mesh(),
513 static_cast<const data&>(baseMesh)
514 ),
515 faceLabels_
516 (
518 (
519 "faceLabels",
520 mesh().facesInstance(),
521 faMesh::meshSubDir,
522 mesh(),
523 IOobject::NO_READ,
524 IOobject::NO_WRITE
525 ),
526 std::move(faceLabels)
527 ),
528 boundary_
529 (
531 (
532 "faBoundary",
533 mesh().facesInstance(),
534 faMesh::meshSubDir,
535 mesh(),
536 IOobject::NO_READ,
537 IOobject::NO_WRITE
538 ),
539 *this,
540 label(0)
541 ),
542 comm_(Pstream::worldComm),
543 curTimeIndex_(time().timeIndex()),
544
545 patchPtr_(nullptr),
546 bndConnectPtr_(nullptr),
547 lduPtr_(nullptr),
548
549 SPtr_(nullptr),
550 S0Ptr_(nullptr),
551 S00Ptr_(nullptr),
552 patchStartsPtr_(nullptr),
553 LePtr_(nullptr),
554 magLePtr_(nullptr),
555 centresPtr_(nullptr),
556 edgeCentresPtr_(nullptr),
557 faceAreaNormalsPtr_(nullptr),
558 edgeAreaNormalsPtr_(nullptr),
559 pointAreaNormalsPtr_(nullptr),
560 faceCurvaturesPtr_(nullptr),
561 edgeTransformTensorsPtr_(nullptr),
562 correctPatchPointNormalsPtr_(nullptr),
563 globalMeshDataPtr_(nullptr),
564
565 haloMapPtr_(nullptr),
566 haloFaceCentresPtr_(nullptr),
567 haloFaceNormalsPtr_(nullptr)
568{}
569
570
571Foam::faMesh::faMesh(const polyPatch& pp, const bool doInit)
572:
573 faMesh
574 (
575 pp.boundaryMesh().mesh(),
576 identity(pp.range())
577 )
578{
579 DebugInFunction << "Creating from polyPatch:" << pp.name() << endl;
580
581 // Add single faPatch "default", but with processor connections
582 faPatchList newPatches
583 (
584 createOnePatch("default")
585 );
586
587 addFaPatches(newPatches);
588
589 setPrimitiveMeshData();
590
591 if (doInit)
592 {
593 faMesh::init(false); // do not init lower levels
594 }
595}
596
597
599(
600 const polyMesh& pMesh,
601 const dictionary& faMeshDefinition,
602 const bool doInit
603)
604:
605 faMesh
606 (
607 pMesh,
609 (
610 pMesh.boundaryMesh(),
611 faMeshDefinition.get<wordRes>("polyMeshPatches")
612 )
613 )
614{
615 DebugInFunction << "Creating from definition (dictionary)" << endl;
616
617 faPatchList newPatches
618 (
619 createPatchList
620 (
621 faMeshDefinition.subDict("boundary"),
622
623 // Optional 'empty' patch
624 faMeshDefinition.getOrDefault<word>("emptyPatch", word::null),
625
626 // Optional specification for default patch
627 faMeshDefinition.findDict("defaultPatch")
628 )
629 );
630
631 addFaPatches(newPatches);
632
633 if (doInit)
634 {
635 faMesh::init(false); // do not init lower levels
636 }
637
638 if (doInit && fileHandler().isFile(pMesh.time().timePath()/"S0"))
639 {
641 (
643 (
644 "S0",
645 time().timeName(),
647 mesh(),
650 ),
651 *this
652 );
653 }
654}
655
656
657// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
658
660{
661 clearOut();
662}
663
664
665// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
666
668{
669 return mesh().dbDir()/faMesh::meshSubDir;
670}
671
672
674{
675 return mesh().time();
676}
677
678
680{
681 return mesh().pointsInstance();
682}
683
684
686{
687 return mesh().facesInstance();
688}
689
690
692{
693 return true;
694}
695
696
698{
699 return mesh().thisDb();
700}
701
702
704{
705 return polyMesh::regionName(thisDb().name());
706}
707
708
709void Foam::faMesh::removeFiles(const fileName& instanceDir) const
710{
711 fileName meshFilesPath = thisDb().time().path()/instanceDir/meshDir();
712
713 Foam::rm(meshFilesPath/"faceLabels");
714 Foam::rm(meshFilesPath/"faBoundary");
715}
716
717
719{
720 removeFiles(mesh().instance());
721}
722
723
725{
726 if (!patchStartsPtr_)
727 {
728 calcPatchStarts();
729 }
730
731 return *patchStartsPtr_;
732}
733
734
736{
737 if (!LePtr_)
738 {
739 calcLe();
740 }
741
742 return *LePtr_;
743}
744
745
747{
748 if (!magLePtr_)
749 {
750 calcMagLe();
751 }
752
753 return *magLePtr_;
754}
755
756
758{
759 if (!centresPtr_)
760 {
761 calcAreaCentres();
762 }
763
764 return *centresPtr_;
765}
766
767
769{
770 if (!edgeCentresPtr_)
771 {
772 calcEdgeCentres();
773 }
774
775 return *edgeCentresPtr_;
776}
777
778
781{
782 if (!SPtr_)
783 {
784 calcS();
785 }
786
787 return *SPtr_;
788}
789
790
793{
794 if (!S0Ptr_)
795 {
797 << "S0 is not available"
798 << abort(FatalError);
799 }
800
801 return *S0Ptr_;
802}
803
804
807{
808 if (!S00Ptr_)
809 {
811 (
813 (
814 "S00",
815 time().timeName(),
816 mesh(),
819 ),
820 S0()
821 );
822
823 S0Ptr_->writeOpt(IOobject::AUTO_WRITE);
824 }
825
826 return *S00Ptr_;
827}
828
829
831{
832 if (!faceAreaNormalsPtr_)
833 {
834 calcFaceAreaNormals();
835 }
836
837 return *faceAreaNormalsPtr_;
838}
839
840
842{
843 if (!edgeAreaNormalsPtr_)
844 {
845 calcEdgeAreaNormals();
846 }
847
848 return *edgeAreaNormalsPtr_;
849}
850
851
853{
854 if (!pointAreaNormalsPtr_)
855 {
856 pointAreaNormalsPtr_.reset(new vectorField(nPoints()));
857
858 calcPointAreaNormals(*pointAreaNormalsPtr_);
859
860 if (quadricsFit_ > 0)
861 {
862 calcPointAreaNormalsByQuadricsFit(*pointAreaNormalsPtr_);
863 }
864 }
865
866 return *pointAreaNormalsPtr_;
867}
868
869
871{
872 if (!faceCurvaturesPtr_)
873 {
874 calcFaceCurvatures();
875 }
876
877 return *faceCurvaturesPtr_;
878}
879
880
883{
884 if (!edgeTransformTensorsPtr_)
885 {
886 calcEdgeTransformTensors();
887 }
888
889 return *edgeTransformTensorsPtr_;
890}
891
892
894{
895 if (!globalMeshDataPtr_)
896 {
897 globalMeshDataPtr_.reset(new faGlobalMeshData(*this));
898 }
899
900 return *globalMeshDataPtr_;
901}
902
903
905{
906 if (!lduPtr_)
907 {
908 calcLduAddressing();
909 }
910
911 return *lduPtr_;
912}
913
914
916{
917 // Grab point motion from polyMesh
918 const vectorField& newPoints = mesh().points();
919
920 // Grab old time areas if the time has been incremented
921 if (curTimeIndex_ < time().timeIndex())
922 {
923 if (S00Ptr_ && S0Ptr_)
924 {
925 DebugInfo<< "Copy old-old S" << endl;
926 *S00Ptr_ = *S0Ptr_;
927 }
928
929 if (S0Ptr_)
930 {
931 DebugInfo<< "Copy old S" << endl;
932 *S0Ptr_ = S();
933 }
934 else
935 {
936 DebugInfo<< "Creating old cell volumes." << endl;
937
939 (
941 (
942 "S0",
943 time().timeName(),
944 mesh(),
947 false
948 ),
949 S()
950 );
951 }
952
953 curTimeIndex_ = time().timeIndex();
954 }
955
956 clearGeomNotAreas();
957
958 // To satisfy the motion interface for MeshObject, const cast is needed
959 if (patchPtr_)
960 {
961 patchPtr_->movePoints(newPoints);
962 }
963
964 // Move boundary points
965 boundary_.movePoints(newPoints);
966
967 // Move interpolation
969
970 // Note: Fluxes were dummy?
971
972 return true;
973}
974
975
977{
978 if ((patchID < 0) || (patchID >= boundary().size()))
979 {
981 << "patchID is not in the valid range"
982 << abort(FatalError);
983 }
984
985 if (correctPatchPointNormalsPtr_)
986 {
987 return (*correctPatchPointNormalsPtr_)[patchID];
988 }
989
990 return false;
991}
992
993
995{
996 if (!correctPatchPointNormalsPtr_)
997 {
998 correctPatchPointNormalsPtr_ = new boolList(boundary().size(), false);
999 }
1000
1001 return *correctPatchPointNormalsPtr_;
1002}
1003
1004
1005bool Foam::faMesh::write(const bool valid) const
1006{
1007 faceLabels_.write();
1008 boundary_.write();
1009
1010 return false;
1011}
1012
1013
1014// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
1015
1017{
1018 return &m != this;
1019}
1020
1021
1023{
1024 return &m == this;
1025}
1026
1027
1028// ************************************************************************* //
scalar range
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
A field of fields is a PtrList of fields with reference counting.
Definition: FieldField.H:80
List< Key > sortedToc() const
The table of contents (the keys) in sorted order.
Definition: HashTable.C:137
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
Templated abstract base-class for optional mesh objects used to automate their allocation to the mesh...
Definition: MeshObject.H:91
Inter-processor communications stream.
Definition: Pstream.H:63
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
fileName timePath() const
Return current time path.
Definition: Time.H:375
bool empty() const noexcept
True if the UList is empty (ie, size() is zero)
Definition: UListI.H:427
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:433
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:106
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Definition: boundaryMesh.H:63
Database for solution data, solver performance and other reduced data.
Definition: data.H:58
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
dictionary * findDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX)
Find and return a sub-dictionary pointer if present.
Definition: dictionaryI.H:127
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:460
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Face to edge interpolation scheme. Included in faMesh.
bool movePoints() const
Do what is necessary if the mesh has moved.
Various mesh related information for a parallel run.
Finite area mesh (used for 2-D non-Euclidian finite area method) defined using a patch of faces on a ...
Definition: faMesh.H:100
virtual bool movePoints()
Update after mesh motion.
Definition: faMesh.C:915
boolList & correctPatchPointNormals() const
Set whether point normals should be corrected for a patch.
Definition: faMesh.C:994
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: faMesh.C:685
const labelList & patchStarts() const
Return patch starts.
Definition: faMesh.C:724
const Time & time() const
Return reference to time.
Definition: faMesh.C:673
virtual const lduAddressing & lduAddr() const
Return ldu addressing.
Definition: faMesh.C:904
fileName meshDir() const
Return the local mesh directory (dbDir()/meshSubDir)
Definition: faMesh.C:667
static const word prefix
The prefix to local: finite-area.
Definition: faMesh.H:501
const edgeVectorField & edgeCentres() const
Return edge centres as edgeVectorField.
Definition: faMesh.C:768
const DimensionedField< scalar, areaMesh > & S00() const
Return old-old-time face areas.
Definition: faMesh.C:806
const polyMesh & mesh() const
Return access to polyMesh.
Definition: faMeshI.H:31
const edgeScalarField & magLe() const
Return edge length magnitudes.
Definition: faMesh.C:746
const DimensionedField< scalar, areaMesh > & S() const
Return face areas.
Definition: faMesh.C:780
const edgeVectorField & Le() const
Return edge length vectors.
Definition: faMesh.C:735
const FieldField< Field, tensor > & edgeTransformTensors() const
Return edge transformation tensors.
Definition: faMesh.C:882
virtual bool hasDb() const
Return true if thisDb() is a valid DB.
Definition: faMesh.C:691
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: faMesh.C:679
virtual const objectRegistry & thisDb() const
Return reference to the mesh database.
Definition: faMesh.C:697
const faGlobalMeshData & globalData() const
Return parallel info.
Definition: faMesh.C:893
const DimensionedField< scalar, areaMesh > & S0() const
Return old-time face areas.
Definition: faMesh.C:792
const areaScalarField & faceCurvatures() const
Return face curvatures.
Definition: faMesh.C:870
const edgeVectorField & edgeAreaNormals() const
Return edge area normals.
Definition: faMesh.C:841
void addFaPatches(faPatchList &plist, const bool validBoundary=true)
Add boundary patches. Constructor helper.
Definition: faMeshPatches.C:38
void removeFiles() const
Remove all files from mesh instance()
Definition: faMesh.C:718
const vectorField & pointAreaNormals() const
Return point area normals.
Definition: faMesh.C:852
virtual ~faMesh()
Destructor.
Definition: faMesh.C:659
static word meshSubDir
The mesh sub-directory name (usually "faMesh")
Definition: faMesh.H:504
const word & regionName() const
The mesh region name or word::null if polyMesh::defaultRegion.
Definition: faMesh.C:703
const areaVectorField & faceAreaNormals() const
Return face area normals.
Definition: faMesh.C:830
const areaVectorField & areaCentres() const
Return face centres as areaVectorField.
Definition: faMesh.C:757
Selector class for finite area differencing schemes. faMesh is derived from faSchemes so that all fie...
Definition: faSchemes.H:60
Selector class for finite area solution. faMesh is derived from faSolution so that all fields have ac...
Definition: faSolution.H:60
A class for handling file names.
Definition: fileName.H:76
static std::string path(const std::string &str)
Return directory path name (part before last /)
Definition: fileNameI.H:176
virtual bool write()
Write the output fields.
virtual const objectRegistry & thisDb() const
Return the object registry - resolve conflict polyMesh/lduMesh.
Definition: fvMesh.H:302
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
The class contains the addressing required by the lduMatrix: upper, lower and losort.
Registry of regIOobjects.
const Time & time() const noexcept
Return time registry.
const word & name() const noexcept
The patch name.
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
labelHashSet patchSet(const UList< wordRe > &patchNames, const bool warnNotFound=true, const bool useGroups=true) const
Return the set of patch IDs corresponding to the given names.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:866
virtual const fileName & dbDir() const
Override the objectRegistry dbDir for a single-region case.
Definition: polyMesh.C:837
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:860
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1083
const word & regionName() const
The mesh region name or word::null if polyMesh::defaultRegion.
Definition: polyMesh.C:854
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
friend bool operator!=(const refineCell &rc1, const refineCell &rc2)
Definition: refineCell.H:106
friend bool operator==(const refineCell &rc1, const refineCell &rc2)
Definition: refineCell.H:97
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:54
A class for handling words, derived from Foam::string.
Definition: word.H:68
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition: zero.H:63
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
volScalarField & p
faceListList boundary
dynamicFvMesh & mesh
Template functions to aid in the implementation of demand driven data.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
const pointField & points
label nPoints
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, false)
word timeName
Definition: getTimeIndex.H:3
#define DebugInfo
Report an information message using Foam::Info.
#define DebugInFunction
Report an information message using Foam::Info.
const std::string patch
OpenFOAM patch number as a std::string.
Namespace for OpenFOAM.
bool rm(const fileName &file)
Remove a file (or its gz equivalent), returning true if successful.
Definition: MSwindows.C:1012
const fileOperation & fileHandler()
Get current file handler.
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i)
Definition: labelList.C:38
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
static labelList selectPatchFaces(const polyBoundaryMesh &pbm, const wordRes &polyPatchNames)
Definition: faMesh.C:68
PrimitivePatch< UIndirectList< face >, const pointField & > uindirectPrimitivePatch
A PrimitivePatch with UIndirectList for the faces, const reference for the point field.
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Field< vector > vectorField
Specialisation of Field<T> for vector.
List< bool > boolList
A List of bools.
Definition: List.H:64
bool isFile(const fileName &name, const bool checkGzip=true, const bool followLink=true)
Does the name exist as a FILE in the file system?
Definition: MSwindows.C:666
error FatalError
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
void deleteDemandDrivenData(DataPtr &dataPtr)
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
label timeIndex
Definition: getTimeIndex.H:30