faMesh.H
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) 2021-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
27Class
28 Foam::faMesh
29
30Description
31 Finite area mesh (used for 2-D non-Euclidian finite area method)
32 defined using a \em patch of faces on a polyMesh
33 (ie, uindirectPrimitivePatch).
34
35 The ordering of faces and points on the faMesh corresponds to
36 the localFaces and localPoints as per Foam::PrimitivePatch but
37 the edge addressing is handled slightly differently.
38 The internal edges of the faMesh will generally correspond identically
39 to the internalEdges of the PrimitivePatch (may change in the future)
40 but the boundaryEdges will be reordered compared to the PrimitivePatch
41 to allow edge boundary slices to be obtained.
42
43SourceFiles
44 faMesh.C
45 faMeshDemandDrivenData.C
46 faMeshPatches.C
47 faMeshTopology.C
48 faMeshUpdate.C
49
50Author
51 Zeljko Tukovic, FMENA
52 Hrvoje Jasak, Wikki Ltd.
53
54\*---------------------------------------------------------------------------*/
55
56#ifndef Foam_faMesh_H
57#define Foam_faMesh_H
58
59#include "MeshObject.H"
60#include "polyMesh.H"
61#include "lduMesh.H"
62#include "faBoundaryMesh.H"
63#include "edgeList.H"
64#include "faceList.H"
65#include "primitiveFieldsFwd.H"
66#include "DimensionedField.H"
67#include "areaFieldsFwd.H"
68#include "edgeFieldsFwd.H"
70#include "edgeInterpolation.H"
71#include "labelIOList.H"
72#include "FieldFields.H"
73#include "faGlobalMeshData.H"
74#include "faSchemes.H"
75#include "faSolution.H"
76#include "data.H"
77
78// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
79
80namespace Foam
81{
82
83// Forward Declarations
84class faMeshBoundaryHalo;
85class faMeshLduAddressing;
86class faMeshMapper;
87class faPatchData;
88
89/*---------------------------------------------------------------------------*\
90 Class faMesh Declaration
91\*---------------------------------------------------------------------------*/
93class faMesh
94:
95 public MeshObject<polyMesh, Foam::UpdateableMeshObject, faMesh>,
96 public lduMesh,
97 public faSchemes,
98 public edgeInterpolation, // may need input from faSchemes
99 public faSolution,
100 public data
101{
102 // Private (internal) classes/structures
103
104 //- A (proc, patchi, patchEdgei) tuple used internally for managing
105 //- patch/patch bookkeeping during construction.
106 // Finite-area patches are stored with negated indices, which makes
107 // them readily identifiable and always sort before normal patches.
108 // Note
109 struct patchTuple
110 :
111 public FixedList<label, 4>
112 {
113
114 // Constructors
115
116 // Inherit constructors
117 using FixedList<label, 4>::FixedList;
118
119 //- Default construct as 'invalid'
120 patchTuple()
121 {
122 clear();
123 }
124
125
126 // Static Member Functions
127
128 // Globally consistent ordering
129 //
130 // 1. sort left/right as lower/higher processor connection
131 // 2. sort by proc/patch/patch index
132 static void sort(UList<Pair<patchTuple>>& list)
133 {
134 for (auto& tuples : list)
135 {
136 tuples.sort();
137 }
138 Foam::stableSort(list);
139 }
140
141
142 // Member Functions
143
144 //- Reset to 'invalid'
145 void clear()
146 {
147 procNo(-1);
148 patchi(labelMax);
149 patchEdgei(-1);
150 meshFacei(-1);
151 }
152
153 //- Valid if proc and edge are non-negative
154 bool valid() const noexcept
155 {
156 return (procNo() >= 0 && patchEdgei() >= 0);
157 }
158
159 // Processor is the first sort index
160 label procNo() const { return (*this)[0]; }
161 void procNo(label val) { (*this)[0] = val; }
162
163 // PatchId (-ve for finiteArea patches) is the second sort index
164 label patchi() const { return (*this)[1]; }
165 void patchi(label val) { (*this)[1] = val; }
166
167 // The patch edge index (on the finiteArea patch)
168 // is the third sort index
169 label patchEdgei() const { return (*this)[2]; }
170 void patchEdgei(label val) { (*this)[2] = val; }
171
172 // The processor-local mesh face is the fourth sort index
173 label meshFacei() const { return (*this)[3]; }
174 void meshFacei(label val) { (*this)[3] = val; }
175
176 //- Return the real patchId
177 label realPatchi() const
178 {
179 const label id = patchi();
180 return (id < 0 ? -(id + 1) : id);
181 }
182
183 //- Set patchId as finiteArea
184 void faPatchi(label val)
185 {
186 patchi(-(val + 1));
187 }
188
189 //- Considered to be finiteArea if patchi < 0
190 bool is_finiteArea() const noexcept
191 {
192 return (patchi() < 0);
193 }
194
195 //- Considered to be processor local
196 bool is_localProc() const noexcept
197 {
198 return (procNo() == Pstream::myProcNo());
199 }
200 };
201
202
203 // Private Data
204
205 //- Face labels
206 labelIOList faceLabels_;
207
208 //- Boundary mesh
209 faBoundaryMesh boundary_;
210
211
212 // Primitive mesh data
213
214 //- Edges, addressing into local point list
215 edgeList edges_;
216
217 //- Edge owner
218 labelList edgeOwner_;
219
220 //- Edge neighbour
221 labelList edgeNeighbour_;
222
223
224 // Primitive size data
225
226 //- Total number of points
227 mutable label nPoints_;
228
229 //- Total number of edges
230 mutable label nEdges_;
231
232 //- Number of internal edges
233 mutable label nInternalEdges_;
234
235 //- Number of faces
236 mutable label nFaces_;
237
238
239 // Communication support, updating
240
241 //- Communicator used for parallel communication
242 label comm_;
243
244 //- Current time index for motion
245 // Note. The whole mechanism will be replaced once the
246 // dimensionedField is created and the dimensionedField
247 // will take care of the old-time levels.
248 mutable label curTimeIndex_;
249
250
251 // Demand-driven data
252
253 //- Primitive patch
254 mutable std::unique_ptr<uindirectPrimitivePatch> patchPtr_;
255
256 //- List of proc/mesh-face for boundary edge neighbours
257 mutable std::unique_ptr<List<labelPair>> bndConnectPtr_;
258
259 //- Ldu addressing data
260 mutable faMeshLduAddressing* lduPtr_;
261
262
263 // Geometric Data
264
265 //- Face areas
267
268 //- Face areas old time level
270
271 //- Face areas old-old time level
272 mutable DimensionedField<scalar, areaMesh>* S00Ptr_;
273
274 //- Patch starts in the edge list
275 mutable labelList* patchStartsPtr_;
276
277 //- Edge length vectors
278 mutable edgeVectorField* LePtr_;
279
280 //- Mag edge length vectors
281 mutable edgeScalarField* magLePtr_;
282
283 //- Face centres
284 mutable areaVectorField* centresPtr_;
285
286 //- Edge centres
287 mutable edgeVectorField* edgeCentresPtr_;
288
289 //- Face area normals
290 mutable areaVectorField* faceAreaNormalsPtr_;
291
292 //- Edge area normals
293 mutable edgeVectorField* edgeAreaNormalsPtr_;
294
295 //- Point area normals
296 mutable std::unique_ptr<vectorField> pointAreaNormalsPtr_;
297
298 //- Face curvatures
299 mutable areaScalarField* faceCurvaturesPtr_;
300
301 //- Edge transformation tensors
302 mutable FieldField<Field, tensor>* edgeTransformTensorsPtr_;
303
304 //- Whether point normals must be corrected for a patch
305 mutable boolList* correctPatchPointNormalsPtr_;
306
307
308 // Other mesh-related data
309
310 //- Parallel info
311 mutable autoPtr<faGlobalMeshData> globalMeshDataPtr_;
312
313 //- Mapping/swapping for boundary edge halo neighbours
314 mutable std::unique_ptr<faMeshBoundaryHalo> haloMapPtr_;
315
316 //- Face centres for boundary edge halo neighbours
317 mutable std::unique_ptr<pointField> haloFaceCentresPtr_;
318
319 //- Face normals for boundary edge halo neighbours
320 mutable std::unique_ptr<vectorField> haloFaceNormalsPtr_;
321
322
323 // Static Private Data
324
325 //- Geometry treatment (0: primitive, 1: standard)
326 static int geometryOrder_;
327
328 //- Quadrics fit for pointAreaNormals (experimental)
329 static const int quadricsFit_;
330
331
332 // Private Member Functions
333
334 //- No copy construct
335 faMesh(const faMesh&) = delete;
336
337 //- No copy assignment
338 void operator=(const faMesh&) = delete;
339
340 //- Set indirect patch, removing any old one.
341 // No communication
342 void initPatch() const;
343
344 //- Set primitive mesh data.
345 // No communication
346 void setPrimitiveMeshData();
347
348 //- Get list of (proc/patchi/patchEdgei/meshFacei) tuple pairs in an
349 //- globally consistent ordering
350 List<Pair<patchTuple>> getBoundaryEdgeConnections() const;
351
352 //- Determine the boundary edge neighbour connections
353 void calcBoundaryConnections() const;
354
355 //- Define boundary edge neighbours (proc/face) based on
356 //- gathered topology information
357 void setBoundaryConnections
358 (
359 const List<Pair<patchTuple>>& bndEdgeConnections
360 ) const;
361
362
363 // Private member functions to calculate demand driven data
364
365 //- Calculate ldu addressing
366 void calcLduAddressing() const;
367
368 //- Calculate patch starts in the edge list
369 void calcPatchStarts() const;
370
371 //- Calculate edge lengths
372 // Triggers communication via calcEdgeAreaNormals
373 void calcLe() const;
374
375 //- Calculate mag edge lengths
376 // No communication
377 void calcMagLe() const;
378
379 //- Calculate face centres
380 // No communication
381 void calcAreaCentres() const;
382
383 //- Calculate edge centres
384 // No communication
385 void calcEdgeCentres() const;
386
387 //- Calculate face areas
388 // No communication
389 void calcS() const;
390
391 //- Calculate face area normals
392 // Triggers communication via calcEdgeAreaNormals
393 void calcFaceAreaNormals() const;
394
395 //- Calculate edge area normals
396 // Triggers communication via pointAreaNormals
397 void calcEdgeAreaNormals() const;
398
399 //- Calculate point area normals
400 // Triggers communication for processor patches
401 void calcPointAreaNormals(vectorField& result) const;
402
403 //- Calculate point area normals by quadrics fit
404 void calcPointAreaNormalsByQuadricsFit(vectorField& result) const;
405
406 //- Calculate face curvatures
407 // Triggers communication via edgeLengthCorrection
408 void calcFaceCurvatures() const;
409
410 //- Calculate edge transformation tensors
411 void calcEdgeTransformTensors() const;
412
413 //- Clear geometry but not the face areas
414 void clearGeomNotAreas() const;
415
416 //- Clear boundary halo information
417 void clearHalo() const;
418
419 //- Clear geometry
420 void clearGeom() const;
421
422 //- Clear addressing
423 void clearAddressing() const;
424
425 //- Clear demand-driven data
426 void clearOut() const;
427
428
429 // Halo handling
430
431 //- Calculate halo centres/normals
432 void calcHaloFaceGeometry() const;
433
434
435 // Helpers
436
437 //- Create a single patch
438 faPatchList createOnePatch
439 (
440 const word& patchName,
441 const word& patchType = ""
442 ) const;
443
444 //- Create list of patches from boundary definition
445 faPatchList createPatchList
446 (
447 const dictionary& bndDict,
448 const word& emptyPatchName = "",
449 const dictionary* defaultPatchDefinition = nullptr
450 ) const;
451
452
453 //- Fatal error if edge labels are out of range
454 void checkBoundaryEdgeLabelRange(const labelUList& edgeLabels) const;
455
456 //- Extract list from contiguous (unordered) boundary data
457 //- to the locally sorted order.
458 template<class T>
459 List<T> boundarySubset
460 (
461 const UList<T>& bndField,
462 const labelUList& edgeLabels
463 ) const
464 {
465 #ifdef FULLDEBUG
466 checkBoundaryEdgeLabelRange(edgeLabels);
467 #endif
468 // Like UIndirectList but with an offset
469 List<T> result(edgeLabels.size());
470 forAll(edgeLabels, i)
471 {
472 result[i] = bndField[edgeLabels[i] - nInternalEdges_];
473 }
474 return result;
475 }
476
477
478 // Static Functions
479
480 //- Test if faSchemes/faSolution files are available
481 static bool hasSystemFiles(const polyMesh& pMesh);
482
483 //- Test if all files needed for read construction are available
484 static bool hasFiles(const polyMesh& pMesh);
485
486
487public:
488
489 // Public Typedefs
490
491 //- The mesh type
492 typedef faMesh Mesh;
493
494 //- The boundary type associated with the mesh
496
497
498 //- Runtime type information
499 TypeName("faMesh");
500
501 //- The prefix to local: %finite-area
502 static const word prefix;
503
504 //- The mesh sub-directory name (usually "faMesh")
505 static word meshSubDir;
506
507
508 // Constructors
509
510 //- Read construct from polyMesh, using its IOobject properties
511 explicit faMesh(const polyMesh& pMesh, const bool doInit = true);
512
513 //- Construct zero-sized from polyMesh
514 // Boundary is added using addFaPatches() member function
515 faMesh(const polyMesh& pMesh, const Foam::zero);
516
517 //- Construct as copy (for dictionaries) and zero-sized
518 //- without boundary, using IOobject properties from polyMesh.
519 // Boundary is added using addFaPatches() member function
520 faMesh(const faMesh& baseMesh, const Foam::zero);
521
522 //- Construct as copy (for dictionaries) and faceLabels
523 //- without boundary, using IOobject properties from polyMesh.
524 // Boundary is added using addFaPatches() member function
525 faMesh(const faMesh& baseMesh, labelList&& faceLabels);
526
527 //- Construct from components (face labels) without boundary,
528 //- using IOobject properties from polyMesh.
529 // Boundary is added using addFaPatches() member function.
530 faMesh(const polyMesh& pMesh, labelList&& faceLabels);
531
532 //- Construct from components (face labels) without boundary,
533 //- using alternative IOobject properties
534 //- (primarily the readOption).
535 // Boundary is added using addFaPatches() member function.
536 faMesh
537 (
538 const polyMesh& pMesh,
540 const IOobject& io
541 );
542
543 //- Construct from single polyPatch
544 explicit faMesh(const polyPatch& pp, const bool doInit = true);
545
546 //- Construct from definition
547 faMesh
548 (
549 const polyMesh& pMesh,
550 const dictionary& faMeshDefinition,
551 const bool doInit = true
552 );
553
554
555 //- Destructor
556 virtual ~faMesh();
557
558
559 // Static Functions
560
561 //- Return the current geometry treatment (0: primitive, 1: standard)
562 // A zero level is with restricted neighbour information
563 static int geometryOrder() noexcept
564 {
565 return geometryOrder_;
566 }
567
568 //- Set the preferred geometry treatment
569 // \return the previous value
570 static int geometryOrder(const int order) noexcept
571 {
572 int old(geometryOrder_);
573 geometryOrder_ = order;
574 return old;
575 }
576
577 //- Read construction from polyMesh if all files are available
578 static autoPtr<faMesh> TryNew(const polyMesh& pMesh);
579
580
581 // Member Functions
582
583 // Topological Change
584
585 //- Add boundary patches. Constructor helper
586 void addFaPatches
587 (
588 faPatchList& plist,
589 const bool validBoundary = true
590 );
591
592 //- Add boundary patches. Constructor helper
593 void addFaPatches
594 (
595 const List<faPatch*>& p,
596 const bool validBoundary = true
597 );
598
599 //- Initialise non-demand-driven data etc
600 bool init(const bool doInit);
601
602
603 // Database
604
605 //- Return access to polyMesh
606 inline const polyMesh& mesh() const;
607
608 //- Interface to referenced polyMesh (similar to GeoMesh)
609 const polyMesh& operator()() const { return mesh(); }
610
611 //- Return the local mesh directory (dbDir()/meshSubDir)
612 fileName meshDir() const;
613
614 //- Return reference to time
615 const Time& time() const;
616
617 //- Return the current instance directory for points
618 // Used in the construction of geometric mesh data dependent
619 // on points
620 const fileName& pointsInstance() const;
621
622 //- Return the current instance directory for faces
623 const fileName& facesInstance() const;
624
625
626 // Communication support
627
628 //- Return communicator used for parallel communication
629 inline label comm() const noexcept;
630
631 //- Return communicator used for parallel communication
632 inline label& comm() noexcept;
633
634
635 // Access: Mesh size parameters
636
637 //- Number of local mesh points
638 inline label nPoints() const noexcept;
639
640 //- Number of local mesh edges
641 inline label nEdges() const noexcept;
642
643 //- Number of internal faces
644 inline label nInternalEdges() const noexcept;
645
646 //- Number of boundary edges (== nEdges - nInternalEdges)
647 inline label nBoundaryEdges() const noexcept;
648
649 //- Number of patch faces
650 inline label nFaces() const noexcept;
651
652
653 // Access: Primitive mesh data
654
655 //- Return local points
656 inline const pointField& points() const;
657
658 //- Return local edges with reordered boundary
659 inline const edgeList& edges() const noexcept;
660
661 //- Sub-list of local internal edges
662 inline const edgeList::subList internalEdges() const;
663
664 //- Return local faces
665 inline const faceList& faces() const;
666
667 //- Edge owner addressing
668 inline const labelList& edgeOwner() const noexcept;
669
670 //- Edge neighbour addressing
671 inline const labelList& edgeNeighbour() const noexcept;
672
673 //- True if the internalEdges use an ordering that does not
674 //- correspond 1-to-1 with the patch internalEdges.
675 inline bool hasInternalEdgeLabels() const noexcept;
676
677
678 // Database
679
680 //- Return true if thisDb() is a valid DB
681 virtual bool hasDb() const;
682
683 //- Return reference to the mesh database
684 virtual const objectRegistry& thisDb() const;
685
686 //- Name function is needed to disambiguate those inherited
687 //- from base classes
688 const word& name() const
689 {
690 return thisDb().name();
691 }
692
693 //- The mesh region name or word::null if polyMesh::defaultRegion
694 const word& regionName() const;
695
696
697 // Access
698
699 //- Return constant reference to boundary mesh
700 inline const faBoundaryMesh& boundary() const noexcept;
701
702 //- Return the underlying polyMesh face labels
703 inline const labelList& faceLabels() const noexcept;
704
705
706 //- Return parallel info
707 const faGlobalMeshData& globalData() const;
708
709 //- Return ldu addressing
710 virtual const lduAddressing& lduAddr() const;
711
712 //- Return a list of pointers for each patch
713 // with only those pointing to interfaces being set
714 virtual lduInterfacePtrsList interfaces() const
715 {
716 return boundary().interfaces();
717 }
718
719 //- Internal face owner
720 const labelUList& owner() const
721 {
722 return lduAddr().lowerAddr();
723 }
724
725 //- Internal face neighbour
726 const labelUList& neighbour() const
727 {
728 return lduAddr().upperAddr();
729 }
730
731 //- True if given edge label is internal to the mesh
732 bool isInternalEdge(const label edgeIndex) const noexcept
733 {
734 return (edgeIndex < nInternalEdges_);
735 }
736
737 //- List of proc/face for the boundary edge neighbours
738 //- using primitive patch edge numbering.
739 inline const List<labelPair>& boundaryConnections() const;
740
741 //- Boundary edge neighbour processors
742 //- (does not include own proc)
743 labelList boundaryProcs() const;
744
745 //- List of proc/size for the boundary edge neighbour processors
746 //- (does not include own proc)
748
749 //- Mapping/swapping for boundary halo neighbours
750 const faMeshBoundaryHalo& boundaryHaloMap() const;
751
752 //- Face centres of boundary halo neighbours
753 const pointField& haloFaceCentres() const;
754
755 //- Face normals of boundary halo neighbours
756 const vectorField& haloFaceNormals() const;
757
758 //- Face centres of boundary halo neighbours for specified patch
759 tmp<pointField> haloFaceCentres(const label patchi) const;
760
761 //- Face normals of boundary halo neighbours for specified patch
762 tmp<vectorField> haloFaceNormals(const label patchi) const;
763
764
765 // Storage management
766
767 //- Remove all files from mesh instance
768 void removeFiles(const fileName& instanceDir) const;
769
770 //- Remove all files from mesh instance()
771 void removeFiles() const;
772
773
774 // Mesh motion and morphing
775
776 //- Is mesh moving
777 bool moving() const
778 {
779 return mesh().moving();
780 }
781
782 //- Update after mesh motion
783 virtual bool movePoints();
784
785 //- Update after topo change
786 virtual void updateMesh(const mapPolyMesh&);
787
788
789 // Mapping
790
791 //- Map all fields in time using given map.
792 virtual void mapFields(const faMeshMapper& mapper) const;
793
794 //- Map face areas in time using given map.
795 virtual void mapOldAreas(const faMeshMapper& mapper) const;
796
797
798 // Demand-driven data
799
800 //- Return constant reference to primitive patch
801 inline const uindirectPrimitivePatch& patch() const;
802
803 //- Return reference to primitive patch
805
806 //- Return patch starts
807 const labelList& patchStarts() const;
808
809 //- Return edge length vectors
810 const edgeVectorField& Le() const;
811
812 //- Return edge length magnitudes
813 const edgeScalarField& magLe() const;
814
815 //- Return face centres as areaVectorField
816 const areaVectorField& areaCentres() const;
817
818 //- Return edge centres as edgeVectorField
819 const edgeVectorField& edgeCentres() const;
820
821 //- Return face areas
823
824 //- Return old-time face areas
826
827 //- Return old-old-time face areas
829
830 //- Return face area normals
831 const areaVectorField& faceAreaNormals() const;
832
833 //- Return edge area normals
834 const edgeVectorField& edgeAreaNormals() const;
835
836 //- Return point area normals
837 const vectorField& pointAreaNormals() const;
838
839 //- Return face curvatures
840 const areaScalarField& faceCurvatures() const;
841
842 //- Return edge transformation tensors
844
845 //- Return internal point labels
847
848 //- Return boundary point labels
850
851 //- Return edge length correction
853
854 //- Whether point normals should be corrected for a patch
855 bool correctPatchPointNormals(const label patchID) const;
856
857 //- Set whether point normals should be corrected for a patch
859
860
861 // Storage management
862
863
864 //- Write mesh
865 virtual bool write(const bool valid = true) const;
866
867
868 // Member Operators
869
870 bool operator!=(const faMesh& m) const;
871
872 bool operator==(const faMesh& m) const;
873};
874
875
876// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
877
878} // End namespace Foam
879
880
881// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
882
883#include "faMeshI.H"
884
885#ifdef NoRepository
886 #include "faPatchFaMeshTemplates.C"
887#endif
888
889// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
890
891#endif
892
893// ************************************************************************* //
Forwards and collection of common area field types.
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
A 1D vector of objects of type <T> with a fixed length <N>.
Definition: FixedList.H:81
FixedList()=default
Default construct.
Generic GeometricField class.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:65
Templated abstract base-class for optional mesh objects used to automate their allocation to the mesh...
Definition: MeshObject.H:91
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: Pair.H:69
A list of faces which address into the list of points.
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
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: UList.H:94
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
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
Face to edge interpolation scheme. Included in faMesh.
Finite area boundary mesh.
lduInterfacePtrsList interfaces() const
Various mesh related information for a parallel run.
Class for obtaining halo face data for the boundary edges. The ordering follows that natural edge ord...
lduAddressing wrapper for faMesh
Class holds all the necessary information for mapping fields associated with faMesh.
Definition: faMeshMapper.H:69
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
labelList internalPoints() const
Return internal point labels.
bool operator==(const faMesh &m) const
Definition: faMesh.C:1022
static int geometryOrder() noexcept
Return the current geometry treatment (0: primitive, 1: standard)
Definition: faMesh.H:562
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: faMesh.C:685
const vectorField & haloFaceNormals() const
Face normals of boundary halo neighbours.
const labelList & patchStarts() const
Return patch starts.
Definition: faMesh.C:724
virtual void mapOldAreas(const faMeshMapper &mapper) const
Map face areas in time using given map.
Definition: faMeshUpdate.C:173
faBoundaryMesh BoundaryMesh
The boundary type associated with the mesh.
Definition: faMesh.H:494
const labelList & edgeNeighbour() const noexcept
Edge neighbour addressing.
Definition: faMeshI.H:116
bool hasInternalEdgeLabels() const noexcept
Definition: faMeshI.H:148
const Time & time() const
Return reference to time.
Definition: faMesh.C:673
const faBoundaryMesh & boundary() const noexcept
Return constant reference to boundary mesh.
Definition: faMeshI.H:38
label comm() const noexcept
Return communicator used for parallel communication.
Definition: faMeshI.H:44
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
const pointField & haloFaceCentres() const
Face centres of boundary halo neighbours.
label nEdges() const noexcept
Number of local mesh edges.
Definition: faMeshI.H:62
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
label nBoundaryEdges() const noexcept
Number of boundary edges (== nEdges - nInternalEdges)
Definition: faMeshI.H:74
const faceList & faces() const
Return local faces.
Definition: faMeshI.H:104
bool isInternalEdge(const label edgeIndex) const noexcept
True if given edge label is internal to the mesh.
Definition: faMesh.H:731
const polyMesh & mesh() const
Return access to polyMesh.
Definition: faMeshI.H:31
bool init(const bool doInit)
Initialise non-demand-driven data etc.
Definition: faMesh.C:268
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 edgeList::subList internalEdges() const
Sub-list of local internal edges.
Definition: faMeshI.H:98
const word & name() const
Definition: faMesh.H:687
virtual void mapFields(const faMeshMapper &mapper) const
Map all fields in time using given map.
Definition: faMeshUpdate.C:151
const edgeList & edges() const noexcept
Return local edges with reordered boundary.
Definition: faMeshI.H:92
const uindirectPrimitivePatch & patch() const
Return constant reference to primitive patch.
Definition: faMeshI.H:128
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
const labelUList & owner() const
Internal face owner.
Definition: faMesh.H:719
virtual const objectRegistry & thisDb() const
Return reference to the mesh database.
Definition: faMesh.C:697
bool operator!=(const faMesh &m) const
Definition: faMesh.C:1016
label nInternalEdges() const noexcept
Number of internal faces.
Definition: faMeshI.H:68
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 polyMesh & operator()() const
Interface to referenced polyMesh (similar to GeoMesh)
Definition: faMesh.H:608
const areaScalarField & faceCurvatures() const
Return face curvatures.
Definition: faMesh.C:870
tmp< edgeScalarField > edgeLengthCorrection() const
Return edge length correction.
label nPoints() const noexcept
Number of local mesh points.
Definition: faMeshI.H:56
virtual lduInterfacePtrsList interfaces() const
Return a list of pointers for each patch.
Definition: faMesh.H:713
const labelUList & neighbour() const
Internal face neighbour.
Definition: faMesh.H:725
label nFaces() const noexcept
Number of patch faces.
Definition: faMeshI.H:80
const List< labelPair > & boundaryConnections() const
Definition: faMeshI.H:155
faMesh Mesh
The mesh type.
Definition: faMesh.H:491
labelList boundaryProcs() const
const edgeVectorField & edgeAreaNormals() const
Return edge area normals.
Definition: faMesh.C:841
static autoPtr< faMesh > TryNew(const polyMesh &pMesh)
Read construction from polyMesh if all files are available.
Definition: faMeshNew.C:146
const pointField & points() const
Return local points.
Definition: faMeshI.H:86
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 labelList & faceLabels() const noexcept
Return the underlying polyMesh face labels.
Definition: faMeshI.H:122
const labelList & edgeOwner() const noexcept
Edge owner addressing.
Definition: faMeshI.H:110
List< labelPair > boundaryProcSizes() const
virtual void updateMesh(const mapPolyMesh &)
Update after topo change.
Definition: faMeshUpdate.C:38
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
TypeName("faMesh")
Runtime type information.
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 faMeshBoundaryHalo & boundaryHaloMap() const
Mapping/swapping for boundary halo neighbours.
labelList boundaryPoints() const
Return boundary point labels.
bool moving() const
Is mesh moving.
Definition: faMesh.H:776
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
The class contains the addressing required by the lduMatrix: upper, lower and losort.
virtual const labelUList & upperAddr() const =0
Return upper addressing.
virtual const labelUList & lowerAddr() const =0
Return lower addressing.
Abstract base class for meshes which provide LDU addressing for the construction of lduMatrix and LDU...
Definition: lduMesh.H:63
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:162
Registry of regIOobjects.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
bool moving() const noexcept
Is mesh moving.
Definition: polyMesh.H:534
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
int myProcNo() const noexcept
Return processor number.
A class for managing temporary objects.
Definition: tmp.H:65
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
volScalarField & p
patchWriters clear()
Forwards for edge field types.
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, false)
Namespace for OpenFOAM.
constexpr label labelMax
Definition: label.H:61
void sort(UList< T > &list)
Sort the list.
Definition: UList.C:342
const direction noexcept
Definition: Scalar.H:223
PtrList< faPatch > faPatchList
Store lists of faPatch as a PtrList.
Definition: faPatch.H:66
void stableSort(UList< T > &list)
Stable sort the list.
Definition: UList.C:356
runTime write()
Forward declarations of the specialisations of Field<T> for scalar, vector and tensor.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73