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 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 "faGlobalMeshData.H"
31 #include "Time.H"
32 #include "polyMesh.H"
33 #include "primitiveMesh.H"
34 #include "demandDrivenData.H"
35 #include "IndirectList.H"
36 #include "areaFields.H"
37 #include "edgeFields.H"
38 #include "faMeshLduAddressing.H"
39 #include "wedgeFaPatch.H"
40 #include "faPatchData.H"
41 #include "SortableList.H"
42 
43 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
44 
45 namespace Foam
46 {
47  defineTypeNameAndDebug(faMesh, 0);
48 }
49 
51 
52 const int Foam::faMesh::quadricsFit_ = 0;
53 
54 
55 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
56 
57 void Foam::faMesh::setPrimitiveMeshData()
58 {
59  DebugInFunction << "Setting primitive data" << endl;
60 
61  const indirectPrimitivePatch& bp = patch();
62 
63  // Set faMesh edges
64  edges_.setSize(bp.nEdges());
65 
66  label edgeI = -1;
67 
68  label nIntEdges = bp.nInternalEdges();
69 
70  for (label curEdge = 0; curEdge < nIntEdges; ++curEdge)
71  {
72  edges_[++edgeI] = bp.edges()[curEdge];
73  }
74 
75  forAll(boundary(), patchI)
76  {
77  const labelList& curP = boundary()[patchI];
78 
79  forAll(curP, eI)
80  {
81  edges_[++edgeI] = bp.edges()[curP[eI]];
82  }
83  }
84 
85  nEdges_ = edges_.size();
86  nInternalEdges_ = nIntEdges;
87 
88 
89  // Set edge owner and neighbour
90  edgeOwner_.setSize(nEdges());
91  edgeNeighbour_.setSize(nInternalEdges());
92 
93  edgeI = -1;
94 
95  const labelListList& bpef = bp.edgeFaces();
96 
97  for (label curEdge = 0; curEdge < nIntEdges; ++curEdge)
98  {
99  edgeOwner_[++edgeI] = bpef[curEdge][0];
100 
101  edgeNeighbour_[edgeI] = bpef[curEdge][1];
102  }
103 
104  forAll(boundary(), patchI)
105  {
106  const labelList& curP = boundary()[patchI];
107 
108  forAll(curP, eI)
109  {
110  edgeOwner_[++edgeI] = bpef[curP[eI]][0];
111  }
112  }
113 
114  // Set number of faces
115  nFaces_ = bp.size();
116 
117  // Set number of points
118  nPoints_ = bp.nPoints();
119 }
120 
121 
122 void Foam::faMesh::clearGeomNotAreas() const
123 {
124  DebugInFunction << "Clearing geometry" << endl;
125 
126  deleteDemandDrivenData(SPtr_);
127  deleteDemandDrivenData(patchPtr_);
128  deleteDemandDrivenData(patchStartsPtr_);
129  deleteDemandDrivenData(LePtr_);
130  deleteDemandDrivenData(magLePtr_);
131  deleteDemandDrivenData(centresPtr_);
132  deleteDemandDrivenData(edgeCentresPtr_);
133  deleteDemandDrivenData(faceAreaNormalsPtr_);
134  deleteDemandDrivenData(edgeAreaNormalsPtr_);
135  deleteDemandDrivenData(pointAreaNormalsPtr_);
136  deleteDemandDrivenData(faceCurvaturesPtr_);
137  deleteDemandDrivenData(edgeTransformTensorsPtr_);
138 }
139 
140 
141 void Foam::faMesh::clearGeom() const
142 {
143  DebugInFunction << "Clearing geometry" << endl;
144 
145  clearGeomNotAreas();
146  deleteDemandDrivenData(S0Ptr_);
147  deleteDemandDrivenData(S00Ptr_);
148  deleteDemandDrivenData(correctPatchPointNormalsPtr_);
149 }
150 
151 
152 void Foam::faMesh::clearAddressing() const
153 {
154  DebugInFunction << "Clearing addressing" << endl;
155 
156  deleteDemandDrivenData(lduPtr_);
157 }
158 
159 
160 void Foam::faMesh::clearOut() const
161 {
162  clearGeom();
163  clearAddressing();
164  deleteDemandDrivenData(globalMeshDataPtr_);
165 }
166 
167 
168 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
169 
170 Foam::faMesh::faMesh(const polyMesh& pMesh)
171 :
172  GeoMesh<polyMesh>(pMesh),
174  edgeInterpolation(*this),
175  faSchemes(mesh()),
176  faSolution(mesh()),
177  data(mesh()),
178  faceLabels_
179  (
180  IOobject
181  (
182  "faceLabels",
183  time().findInstance(meshDir(), "faceLabels"),
184  meshSubDir,
185  mesh(),
186  IOobject::MUST_READ,
187  IOobject::NO_WRITE
188  )
189  ),
190  boundary_
191  (
192  IOobject
193  (
194  "faBoundary",
195  time().findInstance(meshDir(), "faBoundary"),
196  meshSubDir,
197  mesh(),
198  IOobject::MUST_READ,
199  IOobject::NO_WRITE
200  ),
201  *this
202  ),
203  comm_(Pstream::worldComm),
204  patchPtr_(nullptr),
205  lduPtr_(nullptr),
206  curTimeIndex_(time().timeIndex()),
207  SPtr_(nullptr),
208  S0Ptr_(nullptr),
209  S00Ptr_(nullptr),
210  patchStartsPtr_(nullptr),
211  LePtr_(nullptr),
212  magLePtr_(nullptr),
213  centresPtr_(nullptr),
214  edgeCentresPtr_(nullptr),
215  faceAreaNormalsPtr_(nullptr),
216  edgeAreaNormalsPtr_(nullptr),
217  pointAreaNormalsPtr_(nullptr),
218  faceCurvaturesPtr_(nullptr),
219  edgeTransformTensorsPtr_(nullptr),
220  correctPatchPointNormalsPtr_(nullptr),
221  globalMeshDataPtr_(nullptr)
222 {
223  DebugInFunction << "Creating faMesh from IOobject" << endl;
224 
225  setPrimitiveMeshData();
226 
227  // Create global mesh data
228  if (Pstream::parRun())
229  {
230  globalData();
231  }
232 
233  // Calculate topology for the patches (processor-processor comms etc.)
234  boundary_.updateMesh();
235 
236  // Calculate the geometry for the patches (transformation tensors etc.)
237  boundary_.calcGeometry();
238 
239  if (isFile(pMesh.time().timePath()/mesh().dbDir()/"S0"))
240  {
242  (
243  IOobject
244  (
245  "S0",
246  time().timeName(),
247  meshSubDir,
248  mesh(),
251  ),
252  *this
253  );
254  }
255 }
256 
257 
258 Foam::faMesh::faMesh
259 (
260  const polyMesh& pMesh,
261  const labelList& faceLabels
262 )
263 :
264  GeoMesh<polyMesh>(pMesh),
266  edgeInterpolation(*this),
267  faSchemes(mesh()),
268  faSolution(mesh()),
269  data(mesh()),
270  faceLabels_
271  (
272  IOobject
273  (
274  "faceLabels",
275  mesh().facesInstance(),
276  meshSubDir,
277  mesh(),
280  ),
281  faceLabels
282  ),
283  boundary_
284  (
285  IOobject
286  (
287  "faBoundary",
288  mesh().facesInstance(),
289  meshSubDir,
290  mesh(),
293  ),
294  *this,
295  0
296  ),
297  comm_(Pstream::worldComm),
298  patchPtr_(nullptr),
299  lduPtr_(nullptr),
300  curTimeIndex_(time().timeIndex()),
301  SPtr_(nullptr),
302  S0Ptr_(nullptr),
303  S00Ptr_(nullptr),
304  patchStartsPtr_(nullptr),
305  LePtr_(nullptr),
306  magLePtr_(nullptr),
307  centresPtr_(nullptr),
308  edgeCentresPtr_(nullptr),
309  faceAreaNormalsPtr_(nullptr),
310  edgeAreaNormalsPtr_(nullptr),
311  pointAreaNormalsPtr_(nullptr),
312  faceCurvaturesPtr_(nullptr),
313  edgeTransformTensorsPtr_(nullptr),
314  correctPatchPointNormalsPtr_(nullptr),
315  globalMeshDataPtr_(nullptr)
316 {
317  DebugInFunction << "Creating faMesh from components" << endl;
318 }
319 
320 
321 Foam::faMesh::faMesh
322 (
323  const polyMesh& pMesh,
324  const fileName& defFile
325 )
326 :
327  GeoMesh<polyMesh>(pMesh),
329  edgeInterpolation(*this),
330  faSchemes(mesh()),
331  faSolution(mesh()),
332  data(mesh()),
333  faceLabels_
334  (
335  IOobject
336  (
337  "faceLabels",
338  mesh().facesInstance(),
339  meshSubDir,
340  mesh(),
343  ),
344  List<label>(0)
345  ),
346  boundary_
347  (
348  IOobject
349  (
350  "faBoundary",
351  mesh().facesInstance(),
352  meshSubDir,
353  mesh(),
356  ),
357  *this,
358  0
359  ),
360  comm_(Pstream::worldComm),
361  patchPtr_(nullptr),
362  lduPtr_(nullptr),
363  curTimeIndex_(time().timeIndex()),
364  SPtr_(nullptr),
365  S0Ptr_(nullptr),
366  S00Ptr_(nullptr),
367  patchStartsPtr_(nullptr),
368  LePtr_(nullptr),
369  magLePtr_(nullptr),
370  centresPtr_(nullptr),
371  edgeCentresPtr_(nullptr),
372  faceAreaNormalsPtr_(nullptr),
373  edgeAreaNormalsPtr_(nullptr),
374  pointAreaNormalsPtr_(nullptr),
375  faceCurvaturesPtr_(nullptr),
376  edgeTransformTensorsPtr_(nullptr),
377  correctPatchPointNormalsPtr_(nullptr),
378  globalMeshDataPtr_(nullptr)
379 {
380  DebugInFunction << "Creating faMesh from definition file" << endl;
381 
382  // Reading faMeshDefinition dictionary
383  IOdictionary faMeshDefinition
384  (
385  IOobject
386  (
387  defFile,
388  mesh().time().constant(),
389  meshSubDir,
390  mesh(),
393  )
394  );
395 
396  const wordList polyMeshPatches
397  (
398  faMeshDefinition.get<wordList>("polyMeshPatches")
399  );
400 
401  const dictionary& bndDict = faMeshDefinition.subDict("boundary");
402 
403  const wordList faPatchNames(bndDict.toc());
404 
405  List<faPatchData> faPatches(faPatchNames.size() + 1);
406 
407  const polyBoundaryMesh& pbm = pMesh.boundaryMesh();
408 
409  forAll(faPatchNames, patchI)
410  {
411  dictionary curPatchDict = bndDict.subDict(faPatchNames[patchI]);
412 
413  faPatches[patchI].name_ = faPatchNames[patchI];
414 
415  faPatches[patchI].type_ = curPatchDict.get<word>("type");
416 
417  faPatches[patchI].ownPolyPatchID_ =
418  pbm.findPatchID(curPatchDict.get<word>("ownerPolyPatch"));
419 
420  faPatches[patchI].ngbPolyPatchID_ =
421  pbm.findPatchID(curPatchDict.get<word>("neighbourPolyPatch"));
422  }
423 
424 
425  // Setting faceLabels list size
426  label size = 0;
427 
428  labelList patchIDs(polyMeshPatches.size(), -1);
429 
430  forAll(polyMeshPatches, patchI)
431  {
432  patchIDs[patchI] = pbm.findPatchID(polyMeshPatches[patchI]);
433 
434  size += pbm[patchIDs[patchI]].size();
435  }
436 
437  faceLabels_ = labelList(size, -1);
438 
439 
440  // Filling of faceLabels list
441  label faceI = -1;
442 
443  sort(patchIDs);
444 
445  forAll(polyMeshPatches, patchI)
446  {
447  label start = pbm[patchIDs[patchI]].start();
448  label size = pbm[patchIDs[patchI]].size();
449 
450  for (label i = 0; i < size; ++i)
451  {
452  faceLabels_[++faceI] = start + i;
453  }
454  }
455 
456 
457  // Determination of faPatch ID for each boundary edge.
458  // Result is in the bndEdgeFaPatchIDs list
459  labelList faceCells(faceLabels_.size(), -1);
460 
461  forAll(faceCells, faceI)
462  {
463  label faceID = faceLabels_[faceI];
464 
465  faceCells[faceI] = mesh().faceOwner()[faceID];
466  }
467 
468  labelList meshEdges =
469  patch().meshEdges
470  (
471  mesh().edges(),
472  mesh().cellEdges(),
473  faceCells
474  );
475 
476  const labelListList& edgeFaces = mesh().edgeFaces();
477 
478  const label nTotalEdges = patch().nEdges();
479  const label nInternalEdges = patch().nInternalEdges();
480 
481  labelList bndEdgeFaPatchIDs(nTotalEdges - nInternalEdges, -1);
482 
483  for (label edgeI = nInternalEdges; edgeI < nTotalEdges; ++edgeI)
484  {
485  label curMeshEdge = meshEdges[edgeI];
486 
487  labelList curEdgePatchIDs(2, label(-1));
488 
489  label patchI = -1;
490 
491  forAll(edgeFaces[curMeshEdge], faceI)
492  {
493  label curFace = edgeFaces[curMeshEdge][faceI];
494 
495  label curPatchID = pbm.whichPatch(curFace);
496 
497  if (curPatchID != -1)
498  {
499  curEdgePatchIDs[++patchI] = curPatchID;
500  }
501  }
502 
503  for (label pI = 0; pI < faPatches.size() - 1; ++pI)
504  {
505  if
506  (
507  (
508  curEdgePatchIDs[0] == faPatches[pI].ownPolyPatchID_
509  && curEdgePatchIDs[1] == faPatches[pI].ngbPolyPatchID_
510  )
511  ||
512  (
513  curEdgePatchIDs[1] == faPatches[pI].ownPolyPatchID_
514  && curEdgePatchIDs[0] == faPatches[pI].ngbPolyPatchID_
515  )
516  )
517  {
518  bndEdgeFaPatchIDs[edgeI - nInternalEdges] = pI;
519  break;
520  }
521  }
522  }
523 
524 
525  // Set edgeLabels for each faPatch
526  for (label pI = 0; pI < (faPatches.size() - 1); ++pI)
527  {
528  SLList<label> tmpList;
529 
530  forAll(bndEdgeFaPatchIDs, eI)
531  {
532  if (bndEdgeFaPatchIDs[eI] == pI)
533  {
534  tmpList.append(nInternalEdges + eI);
535  }
536  }
537 
538  faPatches[pI].edgeLabels_ = tmpList;
539  }
540 
541  // Check for undefined edges
542  SLList<label> tmpList;
543 
544  forAll(bndEdgeFaPatchIDs, eI)
545  {
546  if (bndEdgeFaPatchIDs[eI] == -1)
547  {
548  tmpList.append(nInternalEdges + eI);
549  }
550  }
551 
552  if (tmpList.size() > 0)
553  {
554  // Check for processor edges
555  labelList allUndefEdges(tmpList);
556  labelList ngbPolyPatch(allUndefEdges.size(), -1);
557  forAll(ngbPolyPatch, edgeI)
558  {
559  label curEdge = allUndefEdges[edgeI];
560 
561  label curPMeshEdge = meshEdges[curEdge];
562 
563  forAll(edgeFaces[curPMeshEdge], faceI)
564  {
565  label curFace = edgeFaces[curPMeshEdge][faceI];
566 
567  if (faceLabels_.found(curFace))
568  {
569  label polyPatchID = pbm.whichPatch(curFace);
570 
571  if (polyPatchID != -1)
572  {
573  ngbPolyPatch[edgeI] = polyPatchID;
574  }
575  }
576  }
577  }
578 
579  // Count ngb processorPolyPatch-es
580  labelHashSet processorPatchSet;
581  forAll(ngbPolyPatch, edgeI)
582  {
583  if (ngbPolyPatch[edgeI] != -1)
584  {
585  if (isA<processorPolyPatch>(pbm[ngbPolyPatch[edgeI]]))
586  {
587  processorPatchSet.insert(ngbPolyPatch[edgeI]);
588  }
589  }
590  }
591  labelList processorPatches(processorPatchSet.toc());
592  faPatches.setSize(faPatches.size() + processorPatches.size());
593 
594  for (label i = 0; i < processorPatches.size(); ++i)
595  {
596  SLList<label> tmpLst;
597 
598  forAll(ngbPolyPatch, eI)
599  {
600  if (ngbPolyPatch[eI] == processorPatches[i])
601  {
602  tmpLst.append(allUndefEdges[eI]);
603  }
604  }
605 
606  faPatches[faPatchNames.size() + i].edgeLabels_ = tmpLst;
607 
608  faPatches[faPatchNames.size() + i].name_ =
609  pbm[processorPatches[i]].name();
610 
611  faPatches[faPatchNames.size() + i].type_ =
612  processorFaPatch::typeName;
613 
614  faPatches[faPatchNames.size() + i].ngbPolyPatchID_ =
615  processorPatches[i];
616  }
617 
618  // Remaining undefined edges
619  SLList<label> undefEdges;
620  forAll(ngbPolyPatch, eI)
621  {
622  if (ngbPolyPatch[eI] == -1)
623  {
624  undefEdges.append(allUndefEdges[eI]);
625  }
626  else if (!isA<processorPolyPatch>(pbm[ngbPolyPatch[eI]]))
627  {
628  undefEdges.append(allUndefEdges[eI]);
629  }
630  }
631 
632  if (undefEdges.size() > 0)
633  {
634  label pI = faPatches.size()-1;
635 
636  faPatches[pI].name_ = "undefined";
637  faPatches[pI].type_ = "patch";
638  faPatches[pI].edgeLabels_ = undefEdges;
639  }
640  else
641  {
642  faPatches.setSize(faPatches.size() - 1);
643  }
644  }
645  else
646  {
647  faPatches.setSize(faPatches.size() - 1);
648  }
649 
650 
651  // Reorder processorFaPatch using
652  // ordering of ngb processorPolyPatch
653  forAll(faPatches, patchI)
654  {
655  if (faPatches[patchI].type_ == processorFaPatch::typeName)
656  {
657  labelList ngbFaces(faPatches[patchI].edgeLabels_.size(), -1);
658 
659  forAll(ngbFaces, edgeI)
660  {
661  label curEdge = faPatches[patchI].edgeLabels_[edgeI];
662 
663  label curPMeshEdge = meshEdges[curEdge];
664 
665  forAll(edgeFaces[curPMeshEdge], faceI)
666  {
667  label curFace = edgeFaces[curPMeshEdge][faceI];
668 
669  label curPatchID = pbm.whichPatch(curFace);
670 
671  if (curPatchID == faPatches[patchI].ngbPolyPatchID_)
672  {
673  ngbFaces[edgeI] = curFace;
674  }
675  }
676  }
677 
678  SortableList<label> sortedNgbFaces(ngbFaces);
679  labelList reorderedEdgeLabels(ngbFaces.size(), -1);
680  for (label i = 0; i < reorderedEdgeLabels.size(); ++i)
681  {
682  reorderedEdgeLabels[i] =
683  faPatches[patchI].edgeLabels_
684  [
685  sortedNgbFaces.indices()[i]
686  ];
687  }
688 
689  faPatches[patchI].edgeLabels_ = reorderedEdgeLabels;
690  }
691  }
692 
693 
694  // Add good patches to faMesh
695  SLList<faPatch*> faPatchLst;
696 
697  for (label pI = 0; pI < faPatches.size(); ++pI)
698  {
699  faPatches[pI].dict_.add("type", faPatches[pI].type_);
700  faPatches[pI].dict_.add("edgeLabels", faPatches[pI].edgeLabels_);
701  faPatches[pI].dict_.add
702  (
703  "ngbPolyPatchIndex",
704  faPatches[pI].ngbPolyPatchID_
705  );
706 
707  if (faPatches[pI].type_ == processorFaPatch::typeName)
708  {
709  if (faPatches[pI].ngbPolyPatchID_ == -1)
710  {
712  << "ngbPolyPatch is not defined for processorFaPatch: "
713  << faPatches[pI].name_
714  << abort(FatalError);
715  }
716 
717  const processorPolyPatch& procPolyPatch =
718  refCast<const processorPolyPatch>
719  (
720  pbm[faPatches[pI].ngbPolyPatchID_]
721  );
722 
723  faPatches[pI].dict_.add("myProcNo", procPolyPatch.myProcNo());
724  faPatches[pI].dict_.add
725  (
726  "neighbProcNo",
727  procPolyPatch.neighbProcNo()
728  );
729  }
730 
731  faPatchLst.append
732  (
734  (
735  faPatches[pI].name_,
736  faPatches[pI].dict_,
737  pI,
738  boundary()
739  ).ptr()
740  );
741  }
742 
743  addFaPatches(List<faPatch*>(faPatchLst));
744 
745  // Create global mesh data
746  if (Pstream::parRun())
747  {
748  globalData();
749  }
750 
751  // Calculate topology for the patches (processor-processor comms etc.)
752  boundary_.updateMesh();
753 
754  // Calculate the geometry for the patches (transformation tensors etc.)
755  boundary_.calcGeometry();
756 
757  if (isFile(mesh().time().timePath()/"S0"))
758  {
760  (
761  IOobject
762  (
763  "S0",
764  time().timeName(),
765  meshSubDir,
766  mesh(),
769  ),
770  *this
771  );
772  }
773 }
774 
775 
776 Foam::faMesh::faMesh
777 (
778  const polyMesh& pMesh,
779  const label polyPatchID
780 )
781 :
782  GeoMesh<polyMesh>(pMesh),
784  edgeInterpolation(*this),
785  faSchemes(mesh()),
786  faSolution(mesh()),
787  data(mesh()),
788  faceLabels_
789  (
790  IOobject
791  (
792  "faceLabels",
793  mesh().facesInstance(),
794  meshSubDir,
795  mesh(),
798  ),
799  labelList(pMesh.boundaryMesh()[polyPatchID].size(), -1)
800  ),
801  boundary_
802  (
803  IOobject
804  (
805  "faBoundary",
806  mesh().facesInstance(),
807  meshSubDir,
808  mesh(),
811  ),
812  *this,
813  0
814  ),
815  comm_(Pstream::worldComm),
816  patchPtr_(nullptr),
817  lduPtr_(nullptr),
818  curTimeIndex_(time().timeIndex()),
819  SPtr_(nullptr),
820  S0Ptr_(nullptr),
821  S00Ptr_(nullptr),
822  patchStartsPtr_(nullptr),
823  LePtr_(nullptr),
824  magLePtr_(nullptr),
825  centresPtr_(nullptr),
826  edgeCentresPtr_(nullptr),
827  faceAreaNormalsPtr_(nullptr),
828  edgeAreaNormalsPtr_(nullptr),
829  pointAreaNormalsPtr_(nullptr),
830  faceCurvaturesPtr_(nullptr),
831  edgeTransformTensorsPtr_(nullptr),
832  correctPatchPointNormalsPtr_(nullptr),
833  globalMeshDataPtr_(nullptr)
834 {
835  DebugInFunction << "Creating faMesh from polyPatch" << endl;
836 
837  const polyBoundaryMesh& pbm = pMesh.boundaryMesh();
838 
839  // Set faceLabels
840  forAll(faceLabels_, faceI)
841  {
842  faceLabels_[faceI] = pbm[polyPatchID].start() + faceI;
843  }
844 
845  // Add one faPatch
846  const indirectPrimitivePatch& bp = patch();
847 
848  const label nTotalEdges = bp.nEdges();
849 
850  const label nInternalEdges = bp.nInternalEdges();
851 
852  labelList edgeLabels(nTotalEdges - nInternalEdges, -1);
853 
854  forAll(edgeLabels, edgeI)
855  {
856  edgeLabels[edgeI] = nInternalEdges + edgeI;
857  }
858 
859  dictionary patchDict;
860 
861  patchDict.add("type", "patch");
862  patchDict.add("edgeLabels", edgeLabels);
863  patchDict.add("ngbPolyPatchIndex", -1);
864 
865  List<faPatch*> faPatchLst(1);
866 
867  faPatchLst[0] = faPatch::New("default", patchDict, 0, boundary()).ptr();
868 
869  addFaPatches(faPatchLst);
870 
871  setPrimitiveMeshData();
872 
873  // Calculate topology for the patches (processor-processor comms etc.)
874  boundary_.updateMesh();
875 
876  // Calculate the geometry for the patches (transformation tensors etc.)
877  boundary_.calcGeometry();
878 }
879 
880 
881 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
882 
884 {
885  clearOut();
886 }
887 
888 
889 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
890 
892 {
893  return mesh().dbDir()/meshSubDir;
894 }
895 
896 
898 {
899  return mesh().time();
900 }
901 
902 
904 {
905  return mesh().pointsInstance();
906 }
907 
908 
910 {
911  return mesh().facesInstance();
912 }
913 
914 
916 {
917  if (!patchPtr_)
918  {
919  patchPtr_ = new indirectPrimitivePatch
920  (
922  (
923  mesh().faces(),
924  faceLabels_
925  ),
926  mesh().points()
927  );
928  }
929 
930  return *patchPtr_;
931 }
932 
933 
935 {
936  if (!patchPtr_)
937  {
938  patchPtr_ = new indirectPrimitivePatch
939  (
941  (
942  mesh().faces(),
943  faceLabels_
944  ),
945  mesh().points()
946  );
947  }
948 
949  return *patchPtr_;
950 }
951 
952 
954 {
955  return patch().localPoints();
956 }
957 
958 
960 {
961  return edges_;
962 }
963 
964 
966 {
967  return patch().localFaces();
968 }
969 
970 
972 {
973  DebugInFunction << "Adding patches to faMesh" << endl;
974 
975  if (boundary().size() > 0)
976  {
978  << "boundary already exists"
979  << abort(FatalError);
980  }
981 
982  boundary_.setSize(p.size());
983 
984  forAll(p, patchI)
985  {
986  boundary_.set(patchI, p[patchI]);
987  }
988 
989  setPrimitiveMeshData();
990 
991  boundary_.checkDefinition();
992 }
993 
994 
995 Foam::label Foam::faMesh::comm() const
996 {
997  return comm_;
998 }
999 
1000 
1001 Foam::label& Foam::faMesh::comm()
1002 {
1003  return comm_;
1004 }
1005 
1006 
1008 {
1009  return true;
1010 }
1011 
1012 
1014 {
1015  return mesh().thisDb();
1016 }
1017 
1018 
1020 {
1021  return boundary_;
1022 }
1023 
1024 
1026 {
1027  if (!patchStartsPtr_)
1028  {
1029  calcPatchStarts();
1030  }
1031 
1032  return *patchStartsPtr_;
1033 }
1034 
1035 
1037 {
1038  if (!LePtr_)
1039  {
1040  calcLe();
1041  }
1042 
1043  return *LePtr_;
1044 }
1045 
1046 
1048 {
1049  if (!magLePtr_)
1050  {
1051  calcMagLe();
1052  }
1053 
1054  return *magLePtr_;
1055 }
1056 
1057 
1059 {
1060  if (!centresPtr_)
1061  {
1062  calcAreaCentres();
1063  }
1064 
1065  return *centresPtr_;
1066 }
1067 
1068 
1070 {
1071  if (!edgeCentresPtr_)
1072  {
1073  calcEdgeCentres();
1074  }
1075 
1076  return *edgeCentresPtr_;
1077 }
1078 
1079 
1082 {
1083  if (!SPtr_)
1084  {
1085  calcS();
1086  }
1087 
1088  return *SPtr_;
1089 }
1090 
1091 
1094 {
1095  if (!S0Ptr_)
1096  {
1098  << "S0 is not available"
1099  << abort(FatalError);
1100  }
1101 
1102  return *S0Ptr_;
1103 }
1104 
1105 
1108 {
1109  if (!S00Ptr_)
1110  {
1112  (
1113  IOobject
1114  (
1115  "S00",
1116  time().timeName(),
1117  mesh(),
1120  ),
1121  S0()
1122  );
1123 
1124  S0Ptr_->writeOpt() = IOobject::AUTO_WRITE;
1125  }
1126 
1127  return *S00Ptr_;
1128 }
1129 
1130 
1132 {
1133  if (!faceAreaNormalsPtr_)
1134  {
1135  calcFaceAreaNormals();
1136  }
1137 
1138  return *faceAreaNormalsPtr_;
1139 }
1140 
1141 
1143 {
1144  if (!edgeAreaNormalsPtr_)
1145  {
1146  calcEdgeAreaNormals();
1147  }
1148 
1149  return *edgeAreaNormalsPtr_;
1150 }
1151 
1152 
1154 {
1155  if (!pointAreaNormalsPtr_)
1156  {
1157  calcPointAreaNormals();
1158 
1159  if (quadricsFit_ > 0)
1160  {
1161  calcPointAreaNormalsByQuadricsFit();
1162  }
1163  }
1164 
1165  return *pointAreaNormalsPtr_;
1166 }
1167 
1168 
1170 {
1171  if (!faceCurvaturesPtr_)
1172  {
1173  calcFaceCurvatures();
1174  }
1175 
1176  return *faceCurvaturesPtr_;
1177 }
1178 
1179 
1182 {
1183  if (!edgeTransformTensorsPtr_)
1184  {
1185  calcEdgeTransformTensors();
1186  }
1187 
1188  return *edgeTransformTensorsPtr_;
1189 }
1190 
1191 
1193 {
1194  if (!globalMeshDataPtr_)
1195  {
1196  globalMeshDataPtr_ = new faGlobalMeshData(*this);
1197  }
1198 
1199  return *globalMeshDataPtr_;
1200 }
1201 
1202 
1204 {
1205  if (!lduPtr_)
1206  {
1207  calcLduAddressing();
1208  }
1209 
1210  return *lduPtr_;
1211 }
1212 
1213 
1215 {
1216  // Grab point motion from polyMesh
1217  const vectorField& newPoints = mesh().points();
1218 
1219  // Grab old time areas if the time has been incremented
1220  if (curTimeIndex_ < time().timeIndex())
1221  {
1222  if (S00Ptr_ && S0Ptr_)
1223  {
1224  DebugInfo<< "Copy old-old S" << endl;
1225  *S00Ptr_ = *S0Ptr_;
1226  }
1227 
1228  if (S0Ptr_)
1229  {
1230  DebugInfo<< "Copy old S" << endl;
1231  *S0Ptr_ = S();
1232  }
1233  else
1234  {
1235  DebugInfo<< "Creating old cell volumes." << endl;
1236 
1238  (
1239  IOobject
1240  (
1241  "S0",
1242  time().timeName(),
1243  mesh(),
1246  false
1247  ),
1248  S()
1249  );
1250  }
1251 
1252  curTimeIndex_ = time().timeIndex();
1253  }
1254 
1255  clearGeomNotAreas();
1256 
1257  // To satisfy the motion interface for MeshObject, const cast is needed
1258  if (patchPtr_)
1259  {
1260  patchPtr_->movePoints(newPoints);
1261  }
1262 
1263  // Move boundary points
1264  const_cast<faBoundaryMesh&>(boundary_).movePoints(newPoints);
1265 
1266  // Move interpolation
1267  const edgeInterpolation& cei = *this;
1269 
1270  // Note: Fluxes were dummy?
1271 
1272  return true;
1273 }
1274 
1275 
1277 {
1278  if ((patchID < 0) || (patchID >= boundary().size()))
1279  {
1281  << "patchID is not in the valid range"
1282  << abort(FatalError);
1283  }
1284 
1285  if (correctPatchPointNormalsPtr_)
1286  {
1287  return (*correctPatchPointNormalsPtr_)[patchID];
1288  }
1289 
1290  return false;
1291 }
1292 
1293 
1295 {
1296  if (!correctPatchPointNormalsPtr_)
1297  {
1298  correctPatchPointNormalsPtr_ = new boolList(boundary().size(), false);
1299  }
1300 
1301  return *correctPatchPointNormalsPtr_;
1302 }
1303 
1304 
1305 bool Foam::faMesh::write(const bool valid) const
1306 {
1307  faceLabels_.write();
1308  boundary_.write();
1309 
1310  return false;
1311 }
1312 
1313 
1314 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
1315 
1316 bool Foam::faMesh::operator!=(const faMesh& m) const
1317 {
1318  return &m != this;
1319 }
1320 
1321 
1322 bool Foam::faMesh::operator==(const faMesh& m) const
1323 {
1324  return &m == this;
1325 }
1326 
1327 
1328 // ************************************************************************* //
Foam::lduAddressing
The class contains the addressing required by the lduMatrix: upper, lower and losort.
Definition: lduAddressing.H:114
Foam::LList::append
void append(const T &item)
Add copy at tail of list.
Definition: LList.H:237
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:71
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:54
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1069
Foam::faMesh::meshDir
fileName meshDir() const
Return the local mesh directory (dbDir()/meshSubDir)
Definition: faMesh.C:891
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::faMesh::Le
const edgeVectorField & Le() const
Return edge length vectors.
Definition: faMesh.C:1036
Foam::faBoundaryMesh::updateMesh
void updateMesh()
Correct faBoundaryMesh after topology update.
Definition: faBoundaryMesh.C:412
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:129
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:73
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::faMesh::faces
const faceList & faces() const
Return faces.
Definition: faMesh.C:965
Foam::fileName
A class for handling file names.
Definition: fileName.H:69
Foam::FieldField
A field of fields is a PtrList of fields with reference counting.
Definition: FieldField.H:53
Foam::faMesh::write
virtual bool write(const bool valid=true) const
Write mesh.
Definition: faMesh.C:1305
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:62
Foam::HashTable::toc
List< Key > toc() const
The table of contents (the keys) in unsorted order.
Definition: HashTable.C:121
Foam::faSchemes
Selector class for finite area differencing schemes. faMesh is derived from faShemes so that all fiel...
Definition: faSchemes.H:52
Foam::faMesh::patch
const indirectPrimitivePatch & patch() const
Return constant reference to primitive patch.
Definition: faMesh.C:915
Foam::fvMesh::thisDb
virtual const objectRegistry & thisDb() const
Return the object registry - resolve conflict polyMesh/lduMesh.
Definition: fvMesh.H:287
demandDrivenData.H
Template functions to aid in the implementation of demand driven data.
Foam::primitiveMesh::edgeFaces
const labelListList & edgeFaces() const
Definition: primitiveMeshEdgeFaces.C:33
Foam::polyMesh::dbDir
virtual const fileName & dbDir() const
Override the objectRegistry dbDir for a single-region case.
Definition: polyMesh.C:829
Foam::PrimitivePatch::nEdges
label nEdges() const
Return number of edges in patch.
Definition: PrimitivePatch.H:322
Foam::objectRegistry::time
const Time & time() const
Return time.
Definition: objectRegistry.H:186
Foam::UPstream::parRun
static bool & parRun()
Test if this a parallel run, or allow modify access.
Definition: UPstream.H:434
Foam::faBoundaryMesh::calcGeometry
void calcGeometry()
Calculate the geometry for the patches.
Definition: faBoundaryMesh.C:175
Foam::polyMesh::facesInstance
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:852
Foam::processorPolyPatch::neighbProcNo
int neighbProcNo() const
Return neighbour processor number.
Definition: processorPolyPatch.H:274
Foam::faMesh::edges
const edgeList & edges() const
Return edges.
Definition: faMesh.C:959
Foam::boolList
List< bool > boolList
A List of bools.
Definition: List.H:69
Foam::isFile
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:658
Foam::faGlobalMeshData
Various mesh related information for a parallel run.
Definition: faGlobalMeshData.H:54
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:444
primitiveMesh.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
faMesh.H
Foam::faPatch::New
static autoPtr< faPatch > New(const word &name, const dictionary &dict, const label index, const faBoundaryMesh &bm)
Definition: faPatchNew.C:35
Foam::faMesh::pointAreaNormals
const vectorField & pointAreaNormals() const
Return point area normals.
Definition: faMesh.C:1153
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:81
Foam::faMesh::~faMesh
virtual ~faMesh()
Destructor.
Definition: faMesh.C:883
Foam::faMesh::S00
const DimensionedField< scalar, areaMesh > & S00() const
Return old-old-time face areas.
Definition: faMesh.C:1107
Foam::DynamicID< polyBoundaryMesh >
polyMesh.H
Foam::HashSet< label, Hash< label > >
Foam::UpdateableMeshObject
Definition: MeshObject.H:241
Foam::polyBoundaryMesh::start
label start() const
The start label of the boundary faces in the polyMesh face list.
Definition: polyBoundaryMesh.C:638
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
Foam::faMesh::patchStarts
const labelList & patchStarts() const
Return patch starts.
Definition: faMesh.C:1025
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::faMesh::globalData
const faGlobalMeshData & globalData() const
Return parallel info.
Definition: faMesh.C:1192
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:60
Foam::LList
Template class for non-intrusive linked lists.
Definition: LList.H:54
Foam::deleteDemandDrivenData
void deleteDemandDrivenData(DataPtr &dataPtr)
Definition: demandDrivenData.H:42
Foam::polyMesh::pointsInstance
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:846
SortableList.H
Foam::faBoundaryMesh
Finite area boundary mesh.
Definition: faBoundaryMesh.H:66
Foam::faMesh::edgeTransformTensors
const FieldField< Field, tensor > & edgeTransformTensors() const
Return edge transformation tensors.
Definition: faMesh.C:1181
Foam::faMesh::S0
const DimensionedField< scalar, areaMesh > & S0() const
Return old-time face areas.
Definition: faMesh.C:1093
Foam::Field< vector >
Foam::faMesh::S
const DimensionedField< scalar, areaMesh > & S() const
Return face areas.
Definition: faMesh.C:1081
Foam::faMesh::operator==
bool operator==(const faMesh &m) const
Definition: faMesh.C:1322
Foam::faMesh::lduAddr
virtual const lduAddressing & lduAddr() const
Return ldu addressing.
Definition: faMesh.C:1203
DebugInFunction
#define DebugInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:365
edgeFields.H
Foam::faMesh::addFaPatches
void addFaPatches(const List< faPatch * > &)
Add boundary patches. Constructor helper.
Definition: faMesh.C:971
Foam::Time::timePath
fileName timePath() const
Return current time path.
Definition: Time.H:375
Foam::faMesh::faceAreaNormals
const areaVectorField & faceAreaNormals() const
Return face area normals.
Definition: faMesh.C:1131
Foam::faMesh::time
const Time & time() const
Return reference to time.
Definition: faMesh.C:897
patchID
label patchID
Definition: boundaryProcessorFaPatchPoints.H:5
Foam::polyMesh::faceOwner
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1107
Foam::dictionary::subDict
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:528
Foam::sort
void sort(UList< T > &a)
Definition: UList.C:254
Foam::IndirectList
A List with indirect addressing.
Definition: IndirectList.H:56
IndirectList.H
Foam::processorPolyPatch
Neighbour processor patch.
Definition: processorPolyPatch.H:58
Foam::faSolution
Selector class for finite area solution. faMesh is derived from faSolution so that all fields have ac...
Definition: faSolution.H:53
Foam::faMesh::magLe
const edgeScalarField & magLe() const
Return edge length magnitudes.
Definition: faMesh.C:1047
Foam::faMesh::facesInstance
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: faMesh.C:909
areaFields.H
timeName
word timeName
Definition: getTimeIndex.H:3
Foam::faMesh::nEdges
label nEdges() const
Definition: faMesh.H:356
Foam::faMesh::operator!=
bool operator!=(const faMesh &m) const
Definition: faMesh.C:1316
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
faPatchData.H
Foam::SortableList
A list that is sorted upon construction or when explicitly requested with the sort() method.
Definition: List.H:63
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::faMesh::points
const pointField & points() const
Return mesh points.
Definition: faMesh.C:953
Foam::GeoMesh
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:48
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::faMesh::nInternalEdges
label nInternalEdges() const
Definition: faMesh.H:361
Foam::edgeInterpolation::movePoints
bool movePoints() const
Do what is necessary if the mesh has moved.
Definition: edgeInterpolation.C:161
Foam::faMesh::edgeAreaNormals
const edgeVectorField & edgeAreaNormals() const
Return edge area normals.
Definition: faMesh.C:1142
Foam::faMesh::mesh
const polyMesh & mesh() const
Return access to polyMesh.
Definition: faMesh.H:323
Foam::faMesh::correctPatchPointNormals
boolList & correctPatchPointNormals() const
Set whether point normals should be corrected for a patch.
Definition: faMesh.C:1294
Foam::indirectPrimitivePatch
PrimitivePatch< IndirectList< face >, const pointField & > indirectPrimitivePatch
A PrimitivePatch with an IndirectList for the faces, const reference for the point field.
Definition: indirectPrimitivePatch.H:49
Foam::PrimitivePatch::nInternalEdges
label nInternalEdges() const
Number of internal edges.
Definition: PrimitivePatch.C:203
Foam::faMesh::thisDb
virtual const objectRegistry & thisDb() const
Return reference to the mesh database.
Definition: faMesh.C:1013
Time.H
Foam::labelListList
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
Foam::faMesh::hasDb
virtual bool hasDb() const
Return true if thisDb() is a valid DB.
Definition: faMesh.C:1007
DebugInfo
#define DebugInfo
Report an information message using Foam::Info.
Definition: messageStream.H:359
Foam::faMesh::faceCurvatures
const areaScalarField & faceCurvatures() const
Return face curvatures.
Definition: faMesh.C:1169
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::polyPatchID
DynamicID< polyBoundaryMesh > polyPatchID
Foam::polyPatchID.
Definition: polyPatchID.H:39
Foam::Pstream
Inter-processor communications stream.
Definition: Pstream.H:56
Foam::UPstream::worldComm
static label worldComm
Default communicator (all processors)
Definition: UPstream.H:295
Foam::List< label >
faMeshLduAddressing.H
Foam::faMesh::meshSubDir
static word meshSubDir
Return the mesh sub-directory name (usually "faMesh")
Definition: faMesh.H:277
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::faMesh::boundary
const faBoundaryMesh & boundary() const
Return constant reference to boundary mesh.
Definition: faMesh.C:1019
Foam::faMesh::edgeCentres
const edgeVectorField & edgeCentres() const
Return edge centres as edgeVectorField.
Definition: faMesh.C:1069
Foam::HashSet::insert
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:181
Foam::SortableList::indices
const labelList & indices() const
Return the list of sorted indices. Updated every sort.
Definition: SortableList.H:113
Foam::faMesh
Finite area mesh. Used for 2-D non-Euclidian finite area method.
Definition: faMesh.H:77
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:275
Foam::MeshObject
Templated abstract base-class for optional mesh objects used to automate their allocation to the mesh...
Definition: MeshObject.H:88
Foam::dictionary::add
entry * add(entry *entryPtr, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:708
timeIndex
label timeIndex
Definition: getTimeIndex.H:4
Foam::processorPolyPatch::myProcNo
int myProcNo() const
Return processor number.
Definition: processorPolyPatch.H:268
Foam::dictionary::toc
wordList toc() const
Return the table of contents.
Definition: dictionary.C:670
Foam::faMesh::pointsInstance
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: faMesh.C:903
Foam::List::setSize
void setSize(const label newSize)
Alias for resize(const label)
Definition: ListI.H:146
Foam::GeometricField
Generic GeometricField class.
Definition: areaFieldsFwd.H:53
Foam::IOobject::NO_READ
Definition: IOobject.H:123
constant
constant condensation/saturation model.
wedgeFaPatch.H
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::faceCells
Smooth ATC in cells next to a set of patches supplied by type.
Definition: faceCells.H:56
Foam::edgeInterpolation
Face to edge interpolation scheme. Included in faMesh.
Definition: edgeInterpolation.H:60
Foam::data
Database for solution data, solver performance and other reduced data.
Definition: data.H:55
boundary
faceListList boundary
Definition: createBlockMesh.H:4
Foam::faMesh::movePoints
virtual bool movePoints()
Update after mesh motion.
Definition: faMesh.C:1214
Foam::faMesh::comm
label comm() const
Return communicator used for parallel communication.
Definition: faMesh.C:995
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatch.H:85
Foam::faMesh::areaCentres
const areaVectorField & areaCentres() const
Return face centres as areaVectorField.
Definition: faMesh.C:1058
Foam::IOobject::MUST_READ
Definition: IOobject.H:120
faGlobalMeshData.H