Loading...
Searching...
No Matches
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-2025 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 faMeshI.H
50
51Author
52 Zeljko Tukovic, FMENA
53 Hrvoje Jasak, Wikki Ltd.
54
55\*---------------------------------------------------------------------------*/
56
57#ifndef Foam_faMesh_H
58#define Foam_faMesh_H
59
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
77// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
78
79namespace Foam
80{
81
82// Forward Declarations
85class faMeshMapper;
86class faPatchData;
87
88/*---------------------------------------------------------------------------*\
89 Class faMeshRegistry Declaration
90\*---------------------------------------------------------------------------*/
91
92//- The objectRegistry for faMesh.
93// This is a separate class to ensure it will be fully constructed before
94// other data members use it, which ensures that its virtual methods are
95// callable during construction (gcc).
97:
98 public objectRegistry
99{
100public:
101
102 // Constructors
103
104 //- Construct an objectRegistry for given area region name.
105 //- Treat an empty name like polyMesh::defaultRegion
106 faMeshRegistry(const word& areaName, const polyMesh& mesh);
107
108
109 // Member Functions
110
111 //- True - thisDb() is a valid registry
112 virtual bool hasDb() const { return true; }
113
114 //- Reference to the mesh database
115 virtual const objectRegistry& thisDb() const { return *this; }
116
117 //- Local directory path of the objectRegistry relative to Time
118 //- with override for the single-region case
119 virtual const fileName& dbDir() const;
120};
121
122
123/*---------------------------------------------------------------------------*\
124 Class faMesh Declaration
125\*---------------------------------------------------------------------------*/
126
127class faMesh
128:
129 public faMeshRegistry,
130 public lduMesh,
131 public faSchemes,
132 public faSolution,
133 public edgeInterpolation // may need input from faSchemes
134{
135 // Private (internal) classes/structures
136
137 //- A (proc, patchi, patchEdgei, meshFacei) tuple used internally
138 //- for managing patch/patch bookkeeping during construction.
139 // Finite-area patches are stored with negated indices (offset -2),
140 // which makes them readily identifiable and always sort before normal
141 // patches.
142 struct patchTuple
143 :
144 public FixedList<label, 4>
145 {
146
147 // Constructors
148
149 // Inherit constructors
150 using FixedList<label, 4>::FixedList;
151
152 //- Default construct as 'invalid'
153 patchTuple()
154 {
155 clear();
156 }
157
158
159 // Static Member Functions
160
161 // Globally consistent ordering
162 //
163 // 1. sort left/right as lower/higher processor connection
164 // 2. sort by proc/patch/patch index
165 static void sort(UList<Pair<patchTuple>>& list)
166 {
167 for (auto& tuples : list)
168 {
169 tuples.sort();
170 }
171 Foam::stableSort(list);
172 }
173
174
175 // Member Functions
176
177 //- Reset to 'invalid'
178 void clear()
179 {
180 procNo(-1);
181 patchi(labelMax);
182 patchEdgei(-1);
183 meshFacei(-1);
184 }
185
186 //- Valid if proc and patch (or patch edge) are non-negative
187 bool valid() const
188 {
189 return (procNo() >= 0 && patchi() != -1);
190 }
191
192 // Processor is the first sort index
193 label procNo() const { return (*this)[0]; }
194 void procNo(label val) { (*this)[0] = val; }
195
196 // PatchId is the second sort index (finiteArea patches are < -1)
197 label patchi() const { return (*this)[1]; }
198 void patchi(label val) { (*this)[1] = val; }
199
200 // The patch edge index (on the finiteArea patch)
201 // is the third sort index
202 label patchEdgei() const { return (*this)[2]; }
203 void patchEdgei(label val) { (*this)[2] = val; }
204
205 // The processor-local mesh face is the fourth sort index
206 label meshFacei() const { return (*this)[3]; }
207 void meshFacei(label val) { (*this)[3] = val; }
208
209 //- Return the real patchId
210 label realPatchi() const
211 {
212 const label id = patchi();
213 return (id < -1 ? -(id + 2) : id);
214 }
215
216 //- Set patchId as finiteArea
217 void faPatchi(label val)
218 {
219 patchi(-(val + 2));
220 }
221
222 //- Considered to be finiteArea if (patchi < -1)
223 bool is_finiteArea() const noexcept
224 {
225 return (patchi() < -1);
226 }
227
228 //- Considered to be processor local
229 bool is_localProc() const
230 {
231 return (procNo() == UPstream::myProcNo());
232 }
233 };
234
235
236 // Private Data
237
238 //- Face labels
239 labelIOList faceLabels_;
240
241 //- Boundary mesh
242 faBoundaryMesh boundary_;
243
244
245 // Primitive mesh data
246
247 //- Edges, addressing into local point list
248 edgeList edges_;
249
250 //- Edge owner
251 labelList edgeOwner_;
252
253 //- Edge neighbour
254 labelList edgeNeighbour_;
255
256
257 // Primitive size data
258
259 //- Total number of points
260 mutable label nPoints_;
261
262 //- Total number of edges
263 mutable label nEdges_;
264
265 //- Number of internal edges
266 mutable label nInternalEdges_;
267
268 //- Number of faces
269 mutable label nFaces_;
270
271
272 // Communication support, updating
273
274 //- Communicator used for parallel communication
275 label comm_;
276
277 //- Current time index for motion
278 // Note. The whole mechanism will be replaced once the
279 // dimensionedField is created and the dimensionedField
280 // will take care of the old-time levels.
281 mutable label curTimeIndex_;
282
283
284 // Demand-driven data
285
286 //- Primitive patch
287 mutable std::unique_ptr<uindirectPrimitivePatch> patchPtr_;
288
289 //- List of poly patch/patch-face for faceLabels
290 mutable std::unique_ptr<List<labelPair>> polyPatchFacesPtr_;
291
292 //- List of polyPatch ids used by the area mesh
293 mutable std::unique_ptr<labelList> polyPatchIdsPtr_;
294
295 //- List of proc/mesh-face for boundary edge neighbours
296 mutable std::unique_ptr<List<labelPair>> bndConnectPtr_;
297
298 //- Ldu addressing data
299 mutable std::unique_ptr<faMeshLduAddressing> lduPtr_;
300
301
302 // Geometric Data
303
304 //- Face areas
305 mutable std::unique_ptr<DimensionedField<scalar, areaMesh>> SPtr_;
306
307 //- Face areas old time level
308 mutable std::unique_ptr<DimensionedField<scalar, areaMesh>> S0Ptr_;
309
310 //- Face areas old-old time level
311 mutable std::unique_ptr<DimensionedField<scalar, areaMesh>> S00Ptr_;
312
313 //- Patch starts in the edge list
314 mutable std::unique_ptr<labelList> patchStartsPtr_;
315
316 //- Edge length vectors
317 mutable std::unique_ptr<edgeVectorField> LePtr_;
318
319 //- Mag edge length vectors
320 mutable std::unique_ptr<edgeScalarField> magLePtr_;
321
322 //- Face centres
323 mutable std::unique_ptr<areaVectorField> faceCentresPtr_;
324
325 //- Edge centres
326 mutable std::unique_ptr<edgeVectorField> edgeCentresPtr_;
327
328 //- Face area normals
329 mutable std::unique_ptr<areaVectorField> faceAreaNormalsPtr_;
330
331 //- Edge area normals
332 mutable std::unique_ptr<edgeVectorField> edgeAreaNormalsPtr_;
333
334 //- Point area normals
335 mutable std::unique_ptr<vectorField> pointAreaNormalsPtr_;
336
337 //- Face curvatures
338 mutable std::unique_ptr<areaScalarField> faceCurvaturesPtr_;
339
340 //- Edge transformation tensors
341 mutable std::unique_ptr<FieldField<Field, tensor>>
342 edgeTransformTensorsPtr_;
343
344 //- Whether point normals must be corrected for a patch
345 mutable std::unique_ptr<boolList> correctPatchPointNormalsPtr_;
346
347
348 // Other mesh-related data
349
350 //- Parallel info
351 mutable autoPtr<faGlobalMeshData> globalMeshDataPtr_;
352
353 //- Mapping/swapping for boundary edge halo neighbours
354 mutable std::unique_ptr<faMeshBoundaryHalo> haloMapPtr_;
355
356 //- Face centres for boundary edge halo neighbours
357 mutable std::unique_ptr<pointField> haloFaceCentresPtr_;
358
359 //- Face normals for boundary edge halo neighbours
360 mutable std::unique_ptr<vectorField> haloFaceNormalsPtr_;
361
362
363 // Static Private Data
364
365 //- The prefix to local: %finite-area
366 static const word prefix_;
367
368 //- Quadrics fit for pointAreaNormals (experimental)
369 static const int quadricsFit_;
370
371
372 // Private Member Functions
373
374 //- Set indirect patch, removing any old one.
375 // No communication
376 void initPatch() const;
377
378 //- Set primitive mesh data.
379 // No communication
380 void setPrimitiveMeshData();
381
382 //- Get list of (proc/patchi/patchEdgei/meshFacei) tuple pairs in an
383 //- globally consistent ordering
384 List<Pair<patchTuple>> getBoundaryEdgeConnections() const;
385
386 //- Determine the boundary edge neighbour connections
387 void calcBoundaryConnections() const;
388
389 //- Define boundary edge neighbours (proc/face) based on
390 //- gathered topology information
391 void setBoundaryConnections
392 (
393 const List<Pair<patchTuple>>& bndEdgeConnections
394 ) const;
395
396
397 // Private Member Functions to calculate demand driven data
398
399 //- Calculate ldu addressing
400 void calcLduAddressing() const;
401
402 //- Calculate patch starts in the edge list
403 void calcPatchStarts() const;
404
405 //- Calculate which polyPatches, polyPatch/local-face
406 //- are related to the areaMesh
407 void calcWhichPatchFaces() const;
408
409 //- Calculate the edge normals (direct calculation)
410 //- with flat boundary addressing
411 // Possible communication
412 tmp<vectorField> calcRawEdgeNormals(int calcType) const;
413
414 //- Calculate edge lengths
415 // Triggers communication via calcEdgeAreaNormals
416 void calcLe() const;
417
418 //- Calculate mag edge lengths
419 // No communication
420 void calcMagLe() const;
421
422 //- Calculate face centres
423 // No communication
424 void calcFaceCentres() const;
425
426 //- Calculate edge centres
427 // No communication
428 void calcEdgeCentres() const;
429
430 //- Calculate face areas
431 // No communication
432 void calcS() const;
433
434 //- Calculate face area normals
435 // Triggers communication via calcEdgeAreaNormals
436 void calcFaceAreaNormals() const;
437
438 //- Calculate edge area normals
439 // Triggers communication via pointAreaNormals
440 void calcEdgeAreaNormals() const;
441
442 //- Calculate point area normals
443 // Triggers communication for processor patches
444 void calcPointAreaNormals(vectorField& result) const;
445
446 //- Calculate point area normals by quadrics fit
447 void calcPointAreaNormalsByQuadricsFit(vectorField& result) const;
448
449 //- Calculate face curvatures
450 // Triggers communication via edgeLengthCorrection
451 void calcFaceCurvatures() const;
452
453 //- Calculate edge transformation tensors
454 void calcEdgeTransformTensors() const;
455
456 //- Clear geometry but not the face areas
457 void clearGeomNotAreas() const;
458
459 //- Has halo face centers/normals
460 bool hasHaloFaceGeometry() const noexcept;
461
462 //- Clear boundary halo information
463 void clearHalo() const;
464
465 //- Clear geometry
466 void clearGeom() const;
467
468 //- Clear addressing
469 void clearAddressing() const;
470
471 //- Clear demand-driven data
472 void clearOut() const;
473
474
475 // Halo handling
476
477 //- Calculate halo centres/normals
478 void calcHaloFaceGeometry() const;
479
480
481 // Helpers
482
483 //- Create a single patch
484 faPatchList createOnePatch
485 (
486 const word& patchName,
487 const word& patchType = ""
488 ) const;
489
490 //- Create list of patches from boundary definition
491 faPatchList createPatchList
492 (
493 const dictionary& bndDict,
494 const word& emptyPatchName = "",
495 const dictionary* defaultPatchDefinition = nullptr
496 ) const;
497
498
499 //- Fatal error if edge labels are out of range
500 void checkBoundaryEdgeLabelRange(const labelUList& edgeLabels) const;
501
502 //- Extract list from contiguous (unordered) boundary data
503 //- to the locally sorted order.
504 template<class T>
505 List<T> boundarySubset
506 (
507 const UList<T>& bndField,
508 const labelUList& edgeLabels
509 ) const
510 {
511 #ifdef FULLDEBUG
512 checkBoundaryEdgeLabelRange(edgeLabels);
513 #endif
514 // Like UIndirectList but with an offset
515 List<T> result(edgeLabels.size());
516 forAll(edgeLabels, i)
517 {
518 result[i] = bndField[edgeLabels[i] - nInternalEdges_];
519 }
520 return result;
521 }
522
523
524 // Static Functions
525
526 //- Test if faSchemes/faSolution files are available
527 static bool hasSystemFiles
528 (
529 const word& areaName,
530 const polyMesh& pMesh
531 );
532
533 //- Test if mesh files needed for read construction are available
534 static bool hasMeshFiles
535 (
536 const word& areaName,
537 const polyMesh& pMesh
538 );
539
540public:
541
542 // Public Typedefs
543
544 //- The mesh type
545 typedef faMesh Mesh;
546
547 //- The boundary type associated with the mesh
549
550
551 // Tuning switches
552
553 //- Geometry treatment
554 static int geometryOrder_;
555
556
557 //- Runtime type information
558 TypeName("faMesh");
559
560 //- The prefix to the parent registry name: %finite-area
561 static const word& prefix() noexcept;
562
563 //- The mesh sub-directory name (usually "faMesh")
564 static word meshSubDir;
565
566
567 // Generated Methods
568
569 //- No copy construct
570 faMesh(const faMesh&) = delete;
571
572 //- No copy assignment
573 void operator=(const faMesh&) = delete;
574
575
576 // Constructors
577
578 //- Construct zero-sized from polyMesh
579 // Boundary is added using addFaPatches() member function
580 faMesh(const word& areaName, const polyMesh&, Foam::zero);
581
582 //- Construct zero-sized from polyMesh
583 // Boundary is added using addFaPatches() member function
584 faMesh(const polyMesh&, Foam::zero);
585
586 //- Read construct from polyMesh, using its IOobject properties
587 faMesh(const word& areaName, const polyMesh&, const bool doInit = true);
588
589 //- Read construct from polyMesh, using its IOobject properties
590 explicit faMesh(const polyMesh&, const bool doInit = true);
591
592 //- Construct from components (face labels) without boundary,
593 //- using specified read properties.
594 // Boundary is added using addFaPatches() member function.
595 faMesh
596 (
597 const word& areaName,
598 const polyMesh& pMesh,
600 IOobjectOption ioOpt
601 );
602
603 //- Construct from components (face labels) without boundary,
604 //- using IOobject properties from polyMesh.
605 // Boundary is added using addFaPatches() member function.
606 faMesh(const word& areaName, const polyMesh&, labelList&& faceLabels);
607
608 //- Construct from components (face labels) without boundary,
609 //- using specified read properties.
610 // Boundary is added using addFaPatches() member function.
612
613 //- Construct from components (face labels) without boundary,
614 //- using IOobject properties from polyMesh.
615 // Boundary is added using addFaPatches() member function.
617
618 //- Construct as copy (for dictionaries) and zero-sized
619 //- without boundary.
620 // Boundary is added using addFaPatches() member function
621 faMesh(const word& areaName, const faMesh& baseMesh, Foam::zero);
622
623 //- Construct as copy (for dictionaries) and zero-sized
624 //- without boundary, using IOobject properties from polyMesh.
625 // Boundary is added using addFaPatches() member function
626 faMesh(const faMesh& baseMesh, Foam::zero);
627
628 //- Construct as copy (for dictionaries) and faceLabels
629 //- without boundary, using specified read properties.
630 // Boundary is added using addFaPatches() member function
631 faMesh
632 (
633 const word& areaName,
634 const faMesh& baseMesh,
636 IOobjectOption ioOpt
637 );
638
639 //- Construct as copy (for dictionaries) and faceLabels
640 //- without boundary, using specified read properties.
641 // Boundary is added using addFaPatches() member function
642 faMesh
643 (
644 const faMesh& baseMesh,
646 IOobjectOption ioOpt
647 );
648
649 //- Construct as copy (for dictionaries) and faceLabels
650 //- without boundary, using IOobject properties from polyMesh.
651 // Boundary is added using addFaPatches() member function
652 faMesh
653 (
654 const word& areaName,
655 const faMesh& baseMesh,
657 );
658
659 //- Construct as copy (for dictionaries) and faceLabels
660 //- without boundary, using read properties from baseMesh.
661 // Boundary is added using addFaPatches() member function
662 faMesh(const faMesh& baseMesh, labelList&& faceLabels);
663
664 //- Construct from single polyPatch
665 faMesh
666 (
667 const word& areaName,
668 const polyPatch& pp,
669 const bool doInit = true
670 );
671
672 //- Construct from single polyPatch
673 explicit faMesh(const polyPatch& pp, const bool doInit = true);
674
675 //- Construct from definition
676 faMesh
677 (
678 const word& areaName,
679 const polyMesh& pMesh,
680 const dictionary& faMeshDefinition,
681 const bool doInit = true
682 );
683
684 //- Construct from definition
685 faMesh
686 (
687 const polyMesh& pMesh,
688 const dictionary& faMeshDefinition,
689 const bool doInit = true
690 );
691
692
693 //- Destructor
694 virtual ~faMesh();
695
696
697 // Static Functions
698
699 //- Return the current geometry treatment
700 // A zero level or negative is with restricted neighbour information
701 static int geometryOrder() noexcept
702 {
703 return geometryOrder_;
704 }
705
706 //- Set the preferred geometry treatment
707 // \return the previous value
708 static int geometryOrder(int order) noexcept
709 {
710 int old(geometryOrder_);
711 geometryOrder_ = order;
712 return old;
713 }
714
715 //- Read construction from polyMesh if all files are available
717 (
718 const word& areaName,
719 const polyMesh& pMesh
720 );
722 //- Read construction from polyMesh if all files are available
723 static autoPtr<faMesh> TryNew(const polyMesh& pMesh);
724
725
726 // Member Functions
727
728 // Topological Change
729
730 //- Add boundary patches. Constructor helper
731 void addFaPatches
732 (
733 faPatchList& plist,
734 const bool validBoundary = true
735 );
736
737 //- Add boundary patches. Constructor helper
738 void addFaPatches
739 (
741 const bool validBoundary = true
742 );
743
744 //- Initialise non-demand-driven data etc
745 bool init(const bool doInit);
746
747 //- Processor/processor synchronisation for geometry fields.
748 // Largely internal use only (slightly hacky).
749 void syncGeom();
751
752 // Database
753
754 //- Find the singleton parent registry (on the polyMesh)
755 //- that contains all objects related to finite-area.
756 static const objectRegistry* registry(const polyMesh& pMesh);
757
758 //- Return the singleton parent registry (on the polyMesh)
759 //- that contains all objects related to finite-area.
760 static const objectRegistry& Registry(const polyMesh& pMesh);
761
762 //- The single-region finite-area region on the polyMesh.
763 //- Uses lookupObject semantics - Fatal if non-existent
764 static const faMesh& mesh(const polyMesh& pMesh);
765
766 //- The single-region or specified finite-area region on the polyMesh.
767 //- Uses lookupObject semantics - Fatal if non-existent
768 static const faMesh& mesh(const polyMesh&, const word& areaRegion);
769
770 //- Return access to polyMesh
771 const polyMesh& mesh() const;
772
773 //- Return the local mesh directory (dbDir()/meshSubDir)
774 fileName meshDir() const;
775
776 //- Return reference to time
777 const Time& time() const;
778
779 //- Return the current instance directory for points
780 // Used in the construction of geometric mesh data dependent
781 // on points
782 const fileName& pointsInstance() const;
783
784 //- Return the current instance directory for faces
785 const fileName& facesInstance() const;
786
787 //- Const reference to the mesh and solver state data
788 const meshState& data() const { return mesh().data(); }
789
790 //- Reference to the mesh and solver state data
791 meshState& data() { return const_cast<polyMesh&>(mesh()).data(); }
792
793
794 // Parallel
795
796 //- Return communicator used for parallel communication
797 label comm() const noexcept { return comm_; }
798
799 //- Return communicator used for parallel communication
800 label& comm() noexcept { return comm_; }
801
802 //- Is demand-driven parallel info available?
803 bool hasGlobalData() const noexcept;
804
805 //- Return parallel info (demand-driven)
806 const faGlobalMeshData& globalData() const;
807
808
809 // Solution Control
810
811 //- Non-null if faSchemes exists (can test as bool).
812 const faSchemes* hasSchemes() const;
813
814 //- Non-null if faSolution exists (can test as bool).
815 const faSolution* hasSolution() const;
816
817 //- Read-access to the faSchemes controls
818 const faSchemes& schemes() const;
819
820 //- Read/write-access to the faSchemes controls
822
823 //- Read-access to the faSolution controls
824 const faSolution& solution() const;
825
826 //- Read/write-access to the faSolution controls
828
829
830 // Access: Mesh size parameters
831
832 //- Number of local mesh points
833 inline label nPoints() const noexcept;
834
835 //- Number of local mesh edges
836 inline label nEdges() const noexcept;
837
838 //- Number of internal faces
839 inline label nInternalEdges() const noexcept;
840
841 //- Number of boundary edges (== nEdges - nInternalEdges)
842 inline label nBoundaryEdges() const noexcept;
843
844 //- Number of patch faces
845 inline label nFaces() const noexcept;
846
847
848 // Access: Primitive mesh data
849
850 //- Return local points
851 inline const pointField& points() const;
852
853 //- Return local edges with reordered boundary
854 inline const edgeList& edges() const noexcept;
855
856 //- Sub-list of local internal edges
857 inline const edgeList::subList internalEdges() const;
858
859 //- Return local faces
860 inline const faceList& faces() const;
861
862 //- Edge owner addressing
863 inline const labelList& edgeOwner() const noexcept;
864
865 //- Edge neighbour addressing
866 inline const labelList& edgeNeighbour() const noexcept;
867
868 //- True if the internalEdges use an ordering that does not
869 //- correspond 1-to-1 with the patch internalEdges.
870 inline bool hasInternalEdgeLabels() const noexcept;
871
872
873 // Database
874
875 //- True - thisDb() is a valid registry
876 virtual bool hasDb() const { return true; }
877
878 //- Reference to the mesh database
879 virtual const objectRegistry& thisDb() const
880 {
881 return faMeshRegistry::thisDb();
882 }
883
884 //- Local directory path of the objectRegistry relative to Time
885 //- with override for the single-region case
886 virtual const fileName& dbDir() const
887 {
888 return faMeshRegistry::dbDir();
889 }
890
891 //- The name of the area mesh
892 const word& name() const
893 {
894 return faMeshRegistry::thisDb().name();
895 }
896
897
898 // Regions
899
900 //- Local registry directory path (relative to Time) for specified
901 //- area mesh (of a single-region volume mesh)
902 static fileName dbDir(const word& areaRegion);
903
904 //- Local registry directory path (relative to Time) for specified
905 //- volume mesh and area mesh combination
906 static fileName dbDir(const word& volRegion, const word& areaRegion);
907
908 //- Local registry directory path (relative to Time) for specified
909 //- volume mesh and area mesh combination
910 static fileName dbDir
911 (
912 const polyMesh& pMesh,
913 const word& areaRegion = word::null
914 );
915
916 //- The local mesh directory name (eg, "faMesh") for specified
917 //- area mesh (of a single-region volume mesh)
918 static fileName meshDir(const word& areaRegion);
919
920 //- The local mesh directory name (eg, "faMesh") for specified
921 //- volume mesh and area mesh combination
922 static fileName meshDir(const word& volRegion, const word& areaRegion);
923
924 //- The local mesh directory name (eg, "faMesh") for specified
925 //- volume mesh and area mesh combination
926 static fileName meshDir
927 (
928 const polyMesh& pMesh,
929 const word& areaRegion = word::null
930 );
931
932 //- The region name or word::null if polyMesh::defaultRegion
933 const word& regionName() const;
934
935
936 // Access
937
938 //- Return constant reference to boundary mesh
939 inline const faBoundaryMesh& boundary() const noexcept;
940
941 //- Return the underlying polyMesh face labels
942 inline const labelList& faceLabels() const noexcept;
943
944 //- The area-face corresponding to the mesh-face, -1 if not found
945 inline label whichFace(const label meshFacei) const;
946
947 //- The polyPatches related to the areaMesh, in sorted order
948 inline const labelList& whichPolyPatches() const;
949
950 //- The polyPatch/local-face for each faceLabels()
951 inline const List<labelPair>& whichPatchFaces() const;
952
953 //- Return ldu addressing
954 virtual const lduAddressing& lduAddr() const;
955
956 //- Return a list of pointers for each patch
957 // with only those pointing to interfaces being set
958 virtual lduInterfacePtrsList interfaces() const
959 {
960 return boundary().interfaces();
961 }
962
963 //- Internal face owner
964 const labelUList& owner() const
965 {
966 return lduAddr().lowerAddr();
967 }
968
969 //- Internal face neighbour
970 const labelUList& neighbour() const
971 {
972 return lduAddr().upperAddr();
973 }
974
975 //- True if given edge label is internal to the mesh
976 bool isInternalEdge(const label edgeIndex) const noexcept
977 {
978 return (edgeIndex < nInternalEdges_);
979 }
980
981 //- List of proc/face for the boundary edge neighbours
982 //- using primitive patch edge numbering.
983 inline const List<labelPair>& boundaryConnections() const;
984
985 //- Boundary edge neighbour processors
986 //- (does not include own proc)
987 labelList boundaryProcs() const;
988
989 //- List of proc/size for the boundary edge neighbour processors
990 //- (does not include own proc)
991 List<labelPair> boundaryProcSizes() const;
992
993 //- Mapping/swapping for boundary halo neighbours
994 const faMeshBoundaryHalo& boundaryHaloMap() const;
995
996 //- Face centres of boundary halo neighbours
997 const pointField& haloFaceCentres() const;
998
999 //- Face unit-normals of boundary halo neighbours
1000 const vectorField& haloFaceNormals() const;
1001
1002 //- Face centres of boundary halo neighbours for specified patch
1003 tmp<pointField> haloFaceCentres(const label patchi) const;
1004
1005 //- Face unit-normals of boundary halo neighbours for specified patch
1006 tmp<vectorField> haloFaceNormals(const label patchi) const;
1007
1008
1009 // Interfacing
1010
1011 //- The volume (owner) cells associated with the area-mesh
1012 labelList faceCells() const;
1013
1014
1015 // Storage Management
1016
1017 //- Remove all files from mesh instance
1018 void removeFiles(const fileName& instanceDir) const;
1019
1020 //- Remove all files from mesh instance()
1021 void removeFiles() const;
1022
1023
1024 //- Has face areas: S()
1025 bool hasFaceAreas() const noexcept { return bool(SPtr_); }
1026
1027 //- Has face centres: areaCentres()
1028 bool hasAreaCentres() const noexcept { return bool(faceCentresPtr_); }
1029
1030 //- Has edge centres: edgeCentres()
1031 bool hasEdgeCentres() const noexcept { return bool(edgeCentresPtr_); }
1032
1033 //- Has edge length vectors: Le()
1034 bool hasLe() const noexcept { return bool(LePtr_); }
1035
1036 //- Has edge length magnitudes: magLe()
1037 bool hasMagLe() const noexcept { return bool(magLePtr_); }
1038
1039 //- Has face area normals: faceAreaNormals()
1040 bool hasFaceAreaNormals() const noexcept
1041 {
1042 return bool(faceAreaNormalsPtr_);
1043 }
1044
1045 //- Has edge area normals: edgeAreaNormals()
1046 bool hasEdgeAreaNormals() const noexcept
1047 {
1048 return bool(edgeAreaNormalsPtr_);
1049 }
1050
1051 //- Has point area normals: pointAreaNormals()
1052 bool hasPointAreaNormals() const noexcept
1053 {
1054 return bool(pointAreaNormalsPtr_);
1055 }
1056
1057 //- Has face curvatures: faceCurvatures()
1058 bool hasFaceCurvatures() const noexcept
1059 {
1060 return bool(faceCurvaturesPtr_);
1061 }
1062
1063 //- Has patch point normals corrections
1065 {
1066 return bool(correctPatchPointNormalsPtr_);
1067 }
1069
1070 // Mesh motion and morphing
1071
1072 //- Is mesh moving
1073 bool moving() const
1074 {
1075 return mesh().moving();
1076 }
1077
1078 //- Update after mesh motion
1079 virtual bool movePoints();
1080
1081 //- Update after topo change
1082 virtual void updateMesh(const mapPolyMesh&);
1083
1084
1085 // Mapping
1087 //- Map all fields in time using given map.
1088 virtual void mapFields(const faMeshMapper& mapper) const;
1089
1090 //- Map face areas in time using given map.
1091 virtual void mapOldAreas(const faMeshMapper& mapper) const;
1092
1093
1094 // Demand-driven data
1095
1096 //- Return constant reference to primitive patch
1097 inline const uindirectPrimitivePatch& patch() const;
1098
1099 //- Return reference to primitive patch
1101
1102 //- Return patch starts
1103 const labelList& patchStarts() const;
1104
1105 //- Return edge length vectors
1106 const edgeVectorField& Le() const;
1107
1108 //- Return edge length magnitudes
1109 const edgeScalarField& magLe() const;
1110
1111 //- Return normalised edge length vectors
1113
1114 //- Return face centres as areaVectorField
1115 const areaVectorField& areaCentres() const;
1116
1117 //- Return edge centres as edgeVectorField
1118 const edgeVectorField& edgeCentres() const;
1119
1120 //- Return face areas
1122
1123 //- Return old-time face areas
1125
1126 //- Return old-old-time face areas
1128
1129 //- Return face area normals
1130 const areaVectorField& faceAreaNormals() const;
1131
1132 //- Return edge area normals
1133 const edgeVectorField& edgeAreaNormals() const;
1134
1135 //- Return point area normals
1136 const vectorField& pointAreaNormals() const;
1137
1138 //- Return face curvatures
1139 const areaScalarField& faceCurvatures() const;
1140
1141 //- Return edge transformation tensors
1143
1144 //- Return internal point labels
1145 labelList internalPoints() const;
1146
1147 //- Return boundary point labels
1148 labelList boundaryPoints() const;
1149
1150 //- Return edge length correction
1152
1153 //- Whether point normals should be corrected for a patch
1154 bool correctPatchPointNormals(const label patchID) const;
1155
1156 //- Set whether point normals should be corrected for a patch
1158
1159
1160 // Storage management
1161
1162
1163 //- Write mesh
1164 virtual bool write(const bool writeOnProc = true) const;
1165
1166
1167 // Member Operators
1168
1169 bool operator!=(const faMesh& m) const;
1170
1171 bool operator==(const faMesh& m) const;
1172
1173
1174 // Housekeeping
1175
1176 //- No call operator. Prior to 2312 was used to obtain polyMesh
1177 void operator()() const = delete;
1178};
1179
1180
1181// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1182
1183} // End namespace Foam
1184
1185
1186// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1187
1188#include "faMeshI.H"
1189
1190#ifdef NoRepository
1191 #include "faPatchFaMeshTemplates.C"
1192#endif
1193
1194// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1195
1196#endif
1197
1198// ************************************************************************* //
Forwards and collection of common area field types.
labelList faceLabels(nFaceLabels)
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
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:77
A 1D vector of objects of type <T> with a fixed length <N>.
Definition FixedList.H:73
friend Ostream & operator(Ostream &, const HashTable< regIOobject *, word, Foam::Hash< word > > &tbl)
A simple container of IOobject preferences. Can also be used for general handling of read/no-read/rea...
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
An ordered pair of two objects of type <T> with first() and second() elements.
Definition Pair.H:66
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
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:89
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Negative if the process is not a...
Definition UPstream.H:1706
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
edgeInterpolation(const edgeInterpolation &)=delete
No copy construct.
Finite area boundary mesh, which is a faPatch list with registered IO, a reference to the associated ...
lduInterfacePtrsList interfaces() const
Return a list of pointers for each patch with only those pointing to interfaces being set.
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.
The objectRegistry for faMesh.
Definition faMesh.H:97
virtual const objectRegistry & thisDb() const
Reference to the mesh database.
Definition faMesh.H:119
virtual const fileName & dbDir() const
Local directory path of the objectRegistry relative to Time with override for the single-region case.
faMeshRegistry(const word &areaName, const polyMesh &mesh)
Construct an objectRegistry for given area region name. Treat an empty name like polyMesh::defaultReg...
virtual bool hasDb() const
True - thisDb() is a valid registry.
Definition faMesh.H:114
Finite area mesh (used for 2-D non-Euclidian finite area method) defined using a patch of faces on a ...
Definition faMesh.H:140
virtual bool movePoints()
Update after mesh motion.
Definition faMesh.C:1280
boolList & correctPatchPointNormals() const
Set whether point normals should be corrected for a patch.
Definition faMesh.C:1357
labelList internalPoints() const
Return internal point labels.
virtual const objectRegistry & thisDb() const
Reference to the mesh database.
Definition faMesh.H:1209
bool operator==(const faMesh &m) const
Definition faMesh.C:1386
const List< labelPair > & whichPatchFaces() const
The polyPatch/local-face for each faceLabels().
Definition faMeshI.H:138
static int geometryOrder() noexcept
Return the current geometry treatment.
Definition faMesh.H:944
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition faMesh.C:1040
bool hasAreaCentres() const noexcept
Has face centres: areaCentres().
Definition faMesh.H:1423
const vectorField & haloFaceNormals() const
Face unit-normals of boundary halo neighbours.
const labelList & patchStarts() const
Return patch starts.
Definition faMesh.C:1083
label & comm() noexcept
Return communicator used for parallel communication.
Definition faMesh.H:1086
bool hasGlobalData() const noexcept
Is demand-driven parallel info available?
Definition faMesh.C:1252
bool hasPatchPointNormalsCorrection() const noexcept
Has patch point normals corrections.
Definition faMesh.H:1475
virtual void mapOldAreas(const faMeshMapper &mapper) const
Map face areas in time using given map.
faBoundaryMesh BoundaryMesh
The boundary type associated with the mesh.
Definition faMesh.H:726
const labelList & edgeNeighbour() const noexcept
Edge neighbour addressing.
Definition faMeshI.H:90
bool hasInternalEdgeLabels() const noexcept
True if the internalEdges use an ordering that does not correspond 1-to-1 with the patch internalEdge...
Definition faMeshI.H:148
const Time & time() const
Return reference to time.
Definition faMesh.C:1028
tmp< edgeVectorField > unitLe() const
Return normalised edge length vectors.
const faBoundaryMesh & boundary() const noexcept
Return constant reference to boundary mesh.
Definition faMeshI.H:24
bool hasPointAreaNormals() const noexcept
Has point area normals: pointAreaNormals().
Definition faMesh.H:1459
faMesh(const faMesh &)=delete
No copy construct.
virtual const lduAddressing & lduAddr() const
Return ldu addressing.
Definition faMesh.C:1269
fileName meshDir() const
Return the local mesh directory (dbDir()/meshSubDir).
Definition faMesh.C:1022
const pointField & haloFaceCentres() const
Face centres of boundary halo neighbours.
virtual const fileName & dbDir() const
Local directory path of the objectRegistry relative to Time with override for the single-region case.
Definition faMesh.H:1218
bool hasFaceAreas() const noexcept
Has face areas: S().
Definition faMesh.H:1418
label nEdges() const noexcept
Number of local mesh edges.
Definition faMeshI.H:36
bool hasLe() const noexcept
Has edge length vectors: Le().
Definition faMesh.H:1433
const edgeVectorField & edgeCentres() const
Return edge centres as edgeVectorField.
Definition faMesh.C:1127
const DimensionedField< scalar, areaMesh > & S00() const
Return old-old-time face areas.
Definition faMesh.C:1165
label nBoundaryEdges() const noexcept
Number of boundary edges (== nEdges - nInternalEdges).
Definition faMeshI.H:48
const faceList & faces() const
Return local faces.
Definition faMeshI.H:78
bool isInternalEdge(const label edgeIndex) const noexcept
True if given edge label is internal to the mesh.
Definition faMesh.H:1345
const polyMesh & mesh() const
Return access to polyMesh.
Definition faMesh.C:1016
bool init(const bool doInit)
Initialise non-demand-driven data etc.
Definition faMesh.C:423
const edgeScalarField & magLe() const
Return edge length magnitudes.
Definition faMesh.C:1105
const DimensionedField< scalar, areaMesh > & S() const
Return face areas.
Definition faMesh.C:1139
void operator=(const faMesh &)=delete
No copy assignment.
const edgeList::subList internalEdges() const
Sub-list of local internal edges.
Definition faMeshI.H:72
const word & name() const
The name of the area mesh.
Definition faMesh.H:1226
label whichFace(const label meshFacei) const
The area-face corresponding to the mesh-face, -1 if not found.
Definition faMeshI.H:122
const faSchemes * hasSchemes() const
Non-null if faSchemes exists (can test as bool).
Definition faMesh.C:980
meshState & data()
Reference to the mesh and solver state data.
Definition faMesh.H:1073
static const objectRegistry * registry(const polyMesh &pMesh)
Find the singleton parent registry (on the polyMesh) that contains all objects related to finite-area...
Definition faMesh.C:145
virtual void mapFields(const faMeshMapper &mapper) const
Map all fields in time using given map.
const faSolution & solution() const
Read-access to the faSolution controls.
Definition faMesh.C:1004
const edgeList & edges() const noexcept
Return local edges with reordered boundary.
Definition faMeshI.H:66
bool hasEdgeCentres() const noexcept
Has edge centres: edgeCentres().
Definition faMesh.H:1428
const meshState & data() const
Const reference to the mesh and solver state data.
Definition faMesh.H:1068
const uindirectPrimitivePatch & patch() const
Return constant reference to primitive patch.
Definition faMeshI.H:102
const faSchemes & schemes() const
Read-access to the faSchemes controls.
Definition faMesh.C:992
const edgeVectorField & Le() const
Return edge length vectors.
Definition faMesh.C:1094
const FieldField< Field, tensor > & edgeTransformTensors() const
Return edge transformation tensors.
Definition faMesh.C:1241
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition faMesh.C:1034
const labelUList & owner() const
Internal face owner.
Definition faMesh.H:1329
static int geometryOrder_
Geometry treatment.
Definition faMesh.H:734
bool operator!=(const faMesh &m) const
Definition faMesh.C:1380
label nInternalEdges() const noexcept
Number of internal faces.
Definition faMeshI.H:42
const faGlobalMeshData & globalData() const
Return parallel info (demand-driven).
Definition faMesh.C:1258
const DimensionedField< scalar, areaMesh > & S0() const
Return old-time face areas.
Definition faMesh.C:1151
static int geometryOrder(int order) noexcept
Set the preferred geometry treatment.
Definition faMesh.H:954
const areaScalarField & faceCurvatures() const
Return face curvatures.
Definition faMesh.C:1229
tmp< edgeScalarField > edgeLengthCorrection() const
Return edge length correction.
label nPoints() const noexcept
Number of local mesh points.
Definition faMeshI.H:30
virtual lduInterfacePtrsList interfaces() const
Return a list of pointers for each patch.
Definition faMesh.H:1321
const labelUList & neighbour() const
Internal face neighbour.
Definition faMesh.H:1337
labelList faceCells() const
The volume (owner) cells associated with the area-mesh.
Definition faMesh.C:1052
void syncGeom()
Processor/processor synchronisation for geometry fields.
Definition faMesh.C:402
virtual bool hasDb() const
True - thisDb() is a valid registry.
Definition faMesh.H:1204
label nFaces() const noexcept
Number of patch faces.
Definition faMeshI.H:54
const List< labelPair > & boundaryConnections() const
List of proc/face for the boundary edge neighbours using primitive patch edge numbering.
Definition faMeshI.H:155
void operator()() const =delete
No call operator. Prior to 2312 was used to obtain polyMesh.
static autoPtr< faMesh > TryNew(const word &areaName, const polyMesh &pMesh)
Read construction from polyMesh if all files are available.
Definition faMeshNew.C:188
faMesh Mesh
The mesh type.
Definition faMesh.H:721
static const word & prefix() noexcept
The prefix to the parent registry name: finite-area.
Definition faMesh.C:66
labelList boundaryProcs() const
Boundary edge neighbour processors (does not include own proc).
const edgeVectorField & edgeAreaNormals() const
Return edge area normals.
Definition faMesh.C:1200
label comm() const noexcept
Return communicator used for parallel communication.
Definition faMesh.H:1081
const labelList & whichPolyPatches() const
The polyPatches related to the areaMesh, in sorted order.
Definition faMeshI.H:128
static const objectRegistry & Registry(const polyMesh &pMesh)
Return the singleton parent registry (on the polyMesh) that contains all objects related to finite-ar...
Definition faMesh.C:157
const pointField & points() const
Return local points.
Definition faMeshI.H:60
void addFaPatches(faPatchList &plist, const bool validBoundary=true)
Add boundary patches. Constructor helper.
void removeFiles() const
Remove all files from mesh instance().
Definition faMesh.C:1077
bool hasFaceCurvatures() const noexcept
Has face curvatures: faceCurvatures().
Definition faMesh.H:1467
const labelList & faceLabels() const noexcept
Return the underlying polyMesh face labels.
Definition faMeshI.H:96
const labelList & edgeOwner() const noexcept
Edge owner addressing.
Definition faMeshI.H:84
List< labelPair > boundaryProcSizes() const
List of proc/size for the boundary edge neighbour processors (does not include own proc).
bool hasEdgeAreaNormals() const noexcept
Has edge area normals: edgeAreaNormals().
Definition faMesh.H:1451
virtual void updateMesh(const mapPolyMesh &)
Update after topo change.
const vectorField & pointAreaNormals() const
Return point area normals.
Definition faMesh.C:1211
virtual ~faMesh()
Destructor.
Definition faMesh.C:972
static word meshSubDir
The mesh sub-directory name (usually "faMesh").
Definition faMesh.H:750
bool hasMagLe() const noexcept
Has edge length magnitudes: magLe().
Definition faMesh.H:1438
TypeName("faMesh")
Runtime type information.
bool hasFaceAreaNormals() const noexcept
Has face area normals: faceAreaNormals().
Definition faMesh.H:1443
const word & regionName() const
The region name or word::null if polyMesh::defaultRegion.
Definition faMesh.C:1046
const areaVectorField & faceAreaNormals() const
Return face area normals.
Definition faMesh.C:1189
const faMeshBoundaryHalo & boundaryHaloMap() const
Mapping/swapping for boundary halo neighbours.
const faSolution * hasSolution() const
Non-null if faSolution exists (can test as bool).
Definition faMesh.C:986
labelList boundaryPoints() const
Return boundary point labels.
bool moving() const
Is mesh moving.
Definition faMesh.H:1486
const areaVectorField & areaCentres() const
Return face centres as areaVectorField.
Definition faMesh.C:1116
Helper class for holding data during faPatch construction. Most data members are exposed at the momen...
Definition faPatchData.H:49
Selector class for finite area differencing schemes. faMesh is derived from faSchemes so that all fie...
Definition faSchemes.H:54
faSchemes(const faSchemes &)=delete
No copy construct.
Selector class for finite area solution. faMesh is derived from faSolution so that all fields have ac...
Definition faSolution.H:54
faSolution(const faSolution &)=delete
No copy construct.
A class for handling file names.
Definition fileName.H:75
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:54
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Database for mesh data, solution data, solver performance and other reduced data.
Definition meshState.H:56
Registry of regIOobjects.
void clear()
Clear all entries from the registry.
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
bool moving() const noexcept
Is mesh moving.
Definition polyMesh.H:732
virtual const meshState & data() const noexcept
Const reference to the mesh and solver state data.
Definition polyMesh.H:559
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
Selector class for relaxation factors, solver type and solution.
Definition solution.H:95
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
static const word null
An empty word.
Definition word.H:84
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition zero.H:58
volScalarField & p
dynamicFvMesh & mesh
Forwards for edge field types.
const pointField & points
label nPoints
Namespace for OpenFOAM.
List< edge > edgeList
List of edge.
Definition edgeList.H:32
Pair< label > labelPair
A pair of labels.
Definition Pair.H:54
GeometricField< vector, faePatchField, edgeMesh > edgeVectorField
List< label > labelList
A List of labels.
Definition List.H:62
IOList< label > labelIOList
IO for a List of label.
Definition labelIOList.H:32
UPtrList< const lduInterface > lduInterfacePtrsList
Store lists of lduInterface as a UPtrList.
List< face > faceList
List of faces.
Definition faceListFwd.H:41
constexpr label labelMax
Definition label.H:55
GeometricField< vector, faPatchField, areaMesh > areaVectorField
GeometricField< scalar, faePatchField, edgeMesh > edgeScalarField
GeometricField< scalar, faPatchField, areaMesh > areaScalarField
void sort(UList< T > &list)
Sort the list.
Definition UList.C:283
Field< vector > vectorField
Specialisation of Field<T> for vector.
List< bool > boolList
A List of bools.
Definition List.H:60
const direction noexcept
Definition scalarImpl.H:265
PrimitivePatch< UIndirectList< face >, const pointField & > uindirectPrimitivePatch
A PrimitivePatch with UIndirectList for the faces, const reference for the point field.
vectorField pointField
pointField is a vectorField.
PtrList< faPatch > faPatchList
Store lists of faPatch as a PtrList.
Definition faPatch.H:64
UList< label > labelUList
A UList of labels.
Definition UList.H:75
void stableSort(UList< T > &list)
Stable sort the list.
Definition UList.C:299
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:299
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition typeInfo.H:68