Loading...
Searching...
No Matches
PrimitivePatch.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) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2016-2023 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::PrimitivePatch
29
30Description
31 A list of faces which address into the list of points.
32
33 The class is templated on the face type (e.g. triangle, polygon etc.)
34 and on the list type of faces and points so that it can refer to
35 existing lists using UList and const pointField& or hold the storage
36 using List and pointField.
37
38SourceFiles
39 PrimitivePatch.C
40 PrimitivePatchAddressing.C
41 PrimitivePatchBdryFaces.C
42 PrimitivePatchBdryPoints.C
43 PrimitivePatchCheck.C
44 PrimitivePatchClear.C
45 PrimitivePatchEdgeLoops.C
46 PrimitivePatchLocalPointOrder.C
47 PrimitivePatchMeshData.C
48 PrimitivePatchMeshEdges.C
49 PrimitivePatchPointAddressing.C
50 PrimitivePatchProjectPoints.C
51
52\*---------------------------------------------------------------------------*/
53
54#ifndef Foam_PrimitivePatch_H
55#define Foam_PrimitivePatch_H
56
57#include "boolList.H"
58#include "labelList.H"
59#include "edgeList.H"
60#include "face.H"
61#include "pointField.H"
62#include "intersection.H"
63#include "HashSet.H"
64#include "SubList.H"
65#include "objectHit.H"
66#include "PrimitivePatchBase.H"
67
68// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
69
70namespace Foam
71{
72
73// Forward Declarations
74template<class T> class Map;
75template<class FaceList, class PointField> class PrimitivePatch;
76
77
78/*---------------------------------------------------------------------------*\
79 Class PrimitivePatch Declaration
80\*---------------------------------------------------------------------------*/
81
82template<class FaceList, class PointField>
84:
85 public PrimitivePatchBase,
86 public FaceList
87{
88public:
89
90 // Public Typedefs
92 //- The face type
93 typedef typename
94 std::remove_reference<FaceList>::type::value_type face_type;
95
96 //- The point type
97 typedef typename
98 std::remove_reference<PointField>::type::value_type point_type;
99
100 //- The face list type
101 typedef FaceList FaceListType;
103 //- The point field type
105
106 //- Deprecated(2020-03) prefer face_type typedef
107 // \deprecated(2020-03) prefer face_type typedef
108 typedef face_type FaceType;
109
110
111 // Public Data Types
112
113 //- Enumeration defining the surface type. Used in check routines.
114 // Ordered from 'good' to 'bad'
115 enum surfaceTopo
116 {
117 MANIFOLD,
118 OPEN,
119 ILLEGAL
120 };
121
122private:
123
124 // Private Data
125
126 //- Reference to global list of points
127 PointField points_;
129
130 // Demand-driven Private Data
131
132 //- Edges of the patch; address into local point list;
133 // sorted with internal edges first in upper-triangular order
134 // and external edges last.
135 mutable std::unique_ptr<edgeList> edgesPtr_;
136
137 //- Which part of edgesPtr_ is internal edges.
138 mutable label nInternalEdges_;
139
140 //- Boundary point labels, addressing into local point list
141 mutable std::unique_ptr<labelList> boundaryPointsPtr_;
142
143 //- Face-face addressing
144 mutable std::unique_ptr<labelListList> faceFacesPtr_;
145
146 //- Edge-face addressing
147 mutable std::unique_ptr<labelListList> edgeFacesPtr_;
148
149 //- Face-edge addressing
150 mutable std::unique_ptr<labelListList> faceEdgesPtr_;
151
152 //- Point-edge addressing
153 mutable std::unique_ptr<labelListList> pointEdgesPtr_;
154
155 //- Point-face addressing
156 mutable std::unique_ptr<labelListList> pointFacesPtr_;
157
158 //- Faces addressing into local point list
159 mutable std::unique_ptr<List<face_type>> localFacesPtr_;
160
161 //- Labels of mesh points
162 mutable std::unique_ptr<labelList> meshPointsPtr_;
163
164 //- Mesh point map. Given the global point index find its
165 //- location in the patch
166 mutable std::unique_ptr<Map<label>> meshPointMapPtr_;
167
168 //- Outside edge loops
169 mutable std::unique_ptr<labelListList> edgeLoopsPtr_;
170
171 //- Points local to patch
172 mutable std::unique_ptr<Field<point_type>> localPointsPtr_;
173
174 //- Local point order for most efficient search
175 mutable std::unique_ptr<labelList> localPointOrderPtr_;
176
177 //- Face centres
178 mutable std::unique_ptr<Field<point_type>> faceCentresPtr_;
179
180 //- Face area vectors
181 mutable std::unique_ptr<Field<point_type>> faceAreasPtr_;
182
183 //- Mag face area
184 mutable std::unique_ptr<Field<scalar>> magFaceAreasPtr_;
185
186 //- Face unit normals
187 mutable std::unique_ptr<Field<point_type>> faceNormalsPtr_;
188
189 //- Point unit normals
190 mutable std::unique_ptr<Field<point_type>> pointNormalsPtr_;
191
192
193 // Private Member Functions
194
195 //- Calculate internal points on a patch
196 void calcInternPoints() const;
197
198 //- Calculate boundary points on a patch
199 void calcBdryPoints() const;
200
201 //- Calculate addressing
202 void calcAddressing() const;
203
204 //- Calculate point-edge addressing
205 void calcPointEdges() const;
206
207 //- Calculate point-face addressing
208 void calcPointFaces() const;
209
210 //- Calculate mesh addressing
211 void calcMeshData() const;
212
213 //- Calculate mesh point map
214 void calcMeshPointMap() const;
215
216 //- Calculate outside edge loops
217 void calcEdgeLoops() const;
218
219 //- Calculate local points
220 void calcLocalPoints() const;
221
222 //- Calculate local point order
223 void calcLocalPointOrder() const;
224
225 //- Calculate face centres
226 void calcFaceCentres() const;
227
228 //- Calculate face area vectors
229 void calcFaceAreas() const;
230
231 //- Calculate face area magnitudes
232 void calcMagFaceAreas() const;
233
234 //- Calculate unit face normals
235 void calcFaceNormals() const;
236
237 //- Calculate unit point normals
238 void calcPointNormals() const;
239
240
241 //- Face-edge-face walk while remaining on a patch point.
242 // Used to determine if surface multiply connected through point.
243 void visitPointRegion
244 (
245 const label pointi,
246 const labelUList& pFaces,
247 const label startFacei,
248 const label startEdgei,
249 UList<bool>& pFacesVisited
250 ) const;
251
252
253public:
254
255 // Constructors
256
257 //- Construct from components
259 (
260 const FaceList& faces,
261 const PointField& points
262 );
263
264 //- Construct from components, transferring faces
266 (
267 FaceList&& faces,
268 const PointField& points
269 );
270
271 //- Construct from components, reuse storage
273 (
274 FaceList& faces,
276 const bool reuse
277 );
278
279 //- Copy construct
281
282
283 //- Destructor
284 virtual ~PrimitivePatch();
285
286 void clearOut();
287
288 void clearGeom();
289
290 void clearTopology();
291
292 void clearPatchMeshAddr();
293
294
295 // Member Functions
296
297 //- Suppress direct swapping, since storage containers may be const
298 void swap(PrimitivePatch&) = delete;
299
300
301 // Access
302
303 //- Return reference to global points
304 const Field<point_type>& points() const noexcept
305 {
306 return points_;
307 }
308
309 //- Number of faces in the patch
310 label nFaces() const noexcept
311 {
312 return FaceList::size();
313 }
314
315
316 // Access functions for demand-driven data
317
318 // Topological data; no mesh required.
319
320 //- Number of points supporting patch faces
321 label nPoints() const
322 {
323 return meshPoints().size();
324 }
325
326 //- Number of edges in patch
327 label nEdges() const
328 {
329 return edges().size();
330 }
331
332 //- Return list of edges, address into LOCAL point list
333 const edgeList& edges() const;
334
335 //- Return sub-list of internal edges, address into LOCAL point list
336 const edgeList::subList internalEdges() const;
337
338 //- Return sub-list of boundary edges, address into LOCAL point list
339 const edgeList::subList boundaryEdges() const;
340
341 //- Number of internal edges
342 label nInternalEdges() const;
344 //- Number of boundary edges == (nEdges() - nInternalEdges())
345 label nBoundaryEdges() const;
346
347 //- Is internal edge?
348 bool isInternalEdge(const label edgei) const
349 {
350 return edgei < nInternalEdges();
351 }
353 //- Return list of boundary points, address into LOCAL point list
354 // Uses edge addressing (if it exists) or calculates directly
355 // from localFaces()
356 const labelList& boundaryPoints() const;
357
358 //- Return face-face addressing
359 const labelListList& faceFaces() const;
360
361 //- Return edge-face addressing
362 const labelListList& edgeFaces() const;
363
364 //- Return face-edge addressing
365 const labelListList& faceEdges() const;
366
367 //- Return point-edge addressing
368 const labelListList& pointEdges() const;
369
370 //- Return point-face addressing
371 const labelListList& pointFaces() const;
372
373 //- Return patch faces addressing into local point list
374 const List<face_type>& localFaces() const;
375
376 //- Extract list of local faces corresponding to
377 //- the boundary edges.
378 labelList boundaryFaces() const;
380 //- Extract sorted list of unique local faces associated with
381 //- the boundary edges.
384
385 // Addressing into mesh
386
387 //- Return labelList of mesh points in patch.
388 // They are constructed by walking through the faces in
389 // incremental order and not sorted anymore.
390 const labelList& meshPoints() const;
391
392 //- Mesh point map.
393 // Given the global point index find its location in the patch
394 const Map<label>& meshPointMap() const;
395
396 //- Return pointField of points in patch
397 const Field<point_type>& localPoints() const;
398
399 //- Return orders the local points for most efficient search
400 const labelList& localPointOrder() const;
402 //- Given a global point index, return the local point index.
403 // If the point is not found, return -1
404 label whichPoint(const label gp) const;
405
406 //- From patch edge to global edge using meshPoints.
407 edge meshEdge(const label edgei) const;
408
409 //- From patch edge to global edge using meshPoints.
410 edge meshEdge(const edge& e) const;
411
412 //- Search for edge (local point labels) and return its
413 //- index in the edge list or -1 if not found.
414 // Ignores invalid or out-of-range edges
415 label findEdge(const edge& e) const;
416
417 //- Return labels of patch edges in the global edge list using
418 //- cell addressing
420 (
421 const edgeList& allEdges,
422 const labelListList& cellEdges,
423 const labelList& faceCells
424 ) const;
425
426 //- Return labels of patch edges into the global edge list using
427 //- basic edge addressing.
429 (
430 const edgeList& allEdges,
432 ) const;
433
434 //- Return label of the local patch edge
435 //- into the global edge list using basic edge addressing.
436 label meshEdge
437 (
438 const label edgei,
439 const edgeList& allEdges,
441 ) const;
442
443 //- Return labels of specified patch edges
444 //- into the global edge list using basic edge addressing.
446 (
447 const labelUList& edgeLabels,
448 const edgeList& allEdges,
450 ) const;
451
452
453 //- Return face centres for patch
454 const Field<point_type>& faceCentres() const;
455
456 //- Return face area vectors for patch
457 const Field<point_type>& faceAreas() const;
459 //- Return face area magnitudes for patch
460 const Field<scalar>& magFaceAreas() const;
461
462 //- Return face unit normals for patch
464
465 //- Return point normals for patch
466 const Field<point_type>& pointNormals() const;
467
468 //- The enclosing (bounding) box for the patch points
469 Pair<point_type> box() const;
470
471 //- The enclosing (bounding) sphere radius^2 for specified face
472 scalar sphere(const label facei) const;
473
475 // Storage Management
476
477 bool hasFaceAreas() const { return bool(faceAreasPtr_); }
478 bool hasFaceCentres() const { return bool(faceCentresPtr_); }
479 bool hasFaceNormals() const { return bool(faceNormalsPtr_); }
480 bool hasPointNormals() const { return bool(pointNormalsPtr_); }
481
482 bool hasBoundaryPoints() const { return bool(boundaryPointsPtr_); }
483
484 // These ones are currently all calculated together:
485 // - edges(), faceFaces(), edgeFaces(), faceEdges()
486
487 bool hasEdges() const { return bool(edgesPtr_); }
488 bool hasFaceFaces() const { return bool(faceFacesPtr_); }
489 bool hasEdgeFaces() const { return bool(edgeFacesPtr_); }
490 bool hasFaceEdges() const { return bool(faceEdgesPtr_); }
491
492 bool hasPointEdges() const { return bool(pointEdgesPtr_); }
493 bool hasPointFaces() const { return bool(pointFacesPtr_); }
495 bool hasMeshPoints() const { return bool(meshPointsPtr_); }
496 bool hasMeshPointMap() const { return bool(meshPointMapPtr_); }
497
498
499 // Other patch operations
500
501 //- Project vertices of patch onto another patch
502 template<class ToPatch>
503 List<objectHit> projectPoints
505 const ToPatch& targetPatch,
506 const Field<point_type>& projectionDirection,
509 ) const;
511 //- Project vertices of patch onto another patch
512 template<class ToPatch>
513 List<objectHit> projectFaceCentres
514 (
515 const ToPatch& targetPatch,
516 const Field<point_type>& projectionDirection,
519 ) const;
520
521 //- Return list of closed loops of boundary vertices.
522 // Edge loops are given as ordered lists of vertices
523 // in local addressing
524 const labelListList& edgeLoops() const;
525
526
527 // Check
528
529 //- Calculate surface type formed by patch,
530 //- optionally recording the indices of illegal edges.
531 //
532 // Return types:
533 // - all edges have two neighbours (manifold)
534 // - some edges have more than two neighbours (illegal)
535 // - other (open)
537
538 //- Check surface formed by patch for manifoldness (see above).
539 // Return true if any incorrect edges are found.
540 // Insert vertices of incorrect edges into set.
541 bool checkTopology
542 (
543 const bool report = false,
545 ) const;
546
547 //- Checks primitivePatch for faces sharing point but not edge.
548 // This denotes a surface that is pinched at a single point
549 // (test for pinched at single edge is already in PrimitivePatch)
550 // Returns true if this situation found and puts conflicting
551 // (mesh)point in set. Based on all the checking routines in
552 // primitiveMesh.
554 (
555 const bool report = false,
557 ) const;
558
559
560 // Edit
562 //- Correct patch after moving points
563 virtual void movePoints(const Field<point_type>&);
564
565
566 // Member Operators
567
568 //- Copy assign faces. Leave points alone (could be a reference).
570
571 //- Move assign faces. Leave points alone (could be a reference).
573
574
575 // Housekeeping
576
577 //- Identical to findEdge
578 label whichEdge(const edge& e) const { return this->findEdge(e); }
579};
580
581
582// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
583
584} // End namespace Foam
585
586// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
587
588#ifdef NoRepository
589 #include "PrimitivePatch.C"
590#endif
591
592// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
593
594#endif
595
596// ************************************************************************* //
labelHashSet * pointSetPtr
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
labelHashSet * badEdgesPtr
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
SubList< edge > subList
Definition List.H:129
A HashTable to objects of type <T> with a label key.
Definition Map.H:54
An ordered pair of two objects of type <T> with first() and second() elements.
Definition Pair.H:66
PrimitivePatchBase()=default
Default construct.
A list of faces which address into the list of points.
label nEdges() const
Number of edges in patch.
bool hasPointFaces() const
label nPoints() const
Number of points supporting patch faces.
label nFaces() const noexcept
Number of faces in the patch.
label nBoundaryEdges() const
Number of boundary edges == (nEdges() - nInternalEdges()).
const labelListList & pointEdges() const
Return point-edge addressing.
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
label nInternalEdges() const
Number of internal edges.
bool checkTopology(const bool report=false, labelHashSet *pointSetPtr=nullptr) const
Check surface formed by patch for manifoldness (see above).
bool hasBoundaryPoints() const
edge meshEdge(const edge &e) const
From patch edge to global edge using meshPoints.
PrimitivePatch(const PrimitivePatch< FaceList, PointField > &pp)
Copy construct.
labelList meshEdges(const edgeList &allEdges, const labelListList &cellEdges, const labelList &faceCells) const
Return labels of patch edges in the global edge list using cell addressing.
bool hasFaceNormals() const
const Map< label > & meshPointMap() const
Mesh point map.
const edgeList::subList internalEdges() const
Return sub-list of internal edges, address into LOCAL point list.
List< objectHit > projectFaceCentres(const ToPatch &targetPatch, const Field< point_type > &projectionDirection, const intersection::algorithm=intersection::FULL_RAY, const intersection::direction=intersection::VECTOR) const
Project vertices of patch onto another patch.
const labelList & meshPoints() const
Return labelList of mesh points in patch.
bool hasFaceAreas() const
edge meshEdge(const label edgei) const
From patch edge to global edge using meshPoints.
void operator=(const PrimitivePatch< FaceList, PointField > &rhs)
Copy assign faces. Leave points alone (could be a reference).
bool hasFaceCentres() const
const Field< point_type > & localPoints() const
Return pointField of points in patch.
label meshEdge(const label edgei, const edgeList &allEdges, const labelListList &pointEdges) const
Return label of the local patch edge into the global edge list using basic edge addressing.
const Field< point_type > & faceNormals() const
Return face unit normals for patch.
const labelList & localPointOrder() const
Return orders the local points for most efficient search.
PrimitivePatch(FaceList &faces, PointField &points, const bool reuse)
Construct from components, reuse storage.
PrimitivePatch(const FaceList &faces, const PointField &points)
Construct from components.
const Field< point_type > & pointNormals() const
Return point normals for patch.
const edgeList::subList boundaryEdges() const
Return sub-list of boundary edges, address into LOCAL point list.
std::remove_reference< const pointField >::type::value_type point_type
bool hasMeshPointMap() const
labelList meshEdges(const edgeList &allEdges, const labelListList &pointEdges) const
Return labels of patch edges into the global edge list using basic edge addressing.
labelList boundaryFaces() const
Extract list of local faces corresponding to the boundary edges.
const Field< point_type > & points() const noexcept
bool checkPointManifold(const bool report=false, labelHashSet *pointSetPtr=nullptr) const
Checks primitivePatch for faces sharing point but not edge.
Pair< point_type > box() const
The enclosing (bounding) box for the patch points.
const Field< point_type > & faceAreas() const
Return face area vectors for patch.
label whichPoint(const label gp) const
Given a global point index, return the local point index.
const labelListList & edgeLoops() const
Return list of closed loops of boundary vertices.
scalar sphere(const label facei) const
The enclosing (bounding) sphere radius^2 for specified face.
virtual void movePoints(const Field< point_type > &)
Correct patch after moving points.
bool isInternalEdge(const label edgei) const
Is internal edge?
const labelList & boundaryPoints() const
Return list of boundary points, address into LOCAL point list.
const Field< point_type > & faceCentres() const
Return face centres for patch.
const labelListList & faceFaces() const
Return face-face addressing.
virtual ~PrimitivePatch()
Destructor.
bool hasPointNormals() const
void swap(PrimitivePatch &)=delete
Suppress direct swapping, since storage containers may be const.
const labelListList & pointFaces() const
Return point-face addressing.
void operator=(PrimitivePatch< FaceList, PointField > &&rhs)
Move assign faces. Leave points alone (could be a reference).
surfaceTopo surfaceType(labelHashSet *badEdgesPtr=nullptr) const
Calculate surface type formed by patch, optionally recording the indices of illegal edges.
std::remove_reference< List< face > >::type::value_type face_type
bool hasFaceEdges() const
const labelListList & edgeFaces() const
Return edge-face addressing.
const labelListList & faceEdges() const
Return face-edge addressing.
const Field< scalar > & magFaceAreas() const
Return face area magnitudes for patch.
List< objectHit > projectPoints(const ToPatch &targetPatch, const Field< point_type > &projectionDirection, const intersection::algorithm=intersection::FULL_RAY, const intersection::direction=intersection::VECTOR) const
Project vertices of patch onto another patch.
label whichEdge(const edge &e) const
Identical to findEdge.
labelList meshEdges(const labelUList &edgeLabels, const edgeList &allEdges, const labelListList &pointEdges) const
Return labels of specified patch edges into the global edge list using basic edge addressing.
bool hasEdgeFaces() const
label findEdge(const edge &e) const
Search for edge (local point labels) and return its index in the edge list or -1 if not found.
const List< face_type > & localFaces() const
Return patch faces addressing into local point list.
bool hasMeshPoints() const
bool hasPointEdges() const
PrimitivePatch(FaceList &&faces, const PointField &points)
Construct from components, transferring faces.
surfaceTopo
Enumeration defining the surface type. Used in check routines.
bool hasFaceFaces() const
labelList uniqBoundaryFaces() const
Extract sorted list of unique local faces associated with the boundary edges.
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
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition edge.H:62
Smooth ATC in cells next to a set of patches supplied by type.
Definition faceCells.H:55
Namespace for OpenFOAM.
List< edge > edgeList
List of edge.
Definition edgeList.H:32
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
void rhs(fvMatrix< typename Expr::value_type > &m, const Expr &expression)
const direction noexcept
Definition scalarImpl.H:265
GeometricField< Type, pointPatchField, pointMesh > PointField
A point field for a given type.
UList< label > labelUList
A UList of labels.
Definition UList.H:75
Info<< "Finished reading KIVA file"<< endl;cellShapeList cellShapes(nPoints);labelList cellZoning(nPoints, -1);const cellModel &hex=cellModel::ref(cellModel::HEX);labelList hexLabels(8);label activeCells=0;labelList pointMap(nPoints);forAll(pointMap, i){ pointMap[i]=i;}for(label i=0;i< nPoints;i++){ if(f[i] > 0.0) { hexLabels[0]=i;hexLabels[1]=i1tab[i];hexLabels[2]=i3tab[i1tab[i]];hexLabels[3]=i3tab[i];hexLabels[4]=i8tab[i];hexLabels[5]=i1tab[i8tab[i]];hexLabels[6]=i3tab[i1tab[i8tab[i]]];hexLabels[7]=i3tab[i8tab[i]];cellShapes[activeCells].reset(hex, hexLabels);edgeList edges=cellShapes[activeCells].edges();forAll(edges, ei) { if(edges[ei].mag(points)< SMALL) { label start=pointMap[edges[ei].start()];while(start !=pointMap[start]) { start=pointMap[start];} label end=pointMap[edges[ei].end()];while(end !=pointMap[end]) { end=pointMap[end];} pointMap[start]=pointMap[end]=Foam::min(start, end);} } cellZoning[activeCells]=idreg[i];activeCells++;}}cellShapes.setSize(activeCells);cellZoning.setSize(activeCells);forAll(cellShapes, celli){ cellShape &cs=cellShapes[celli];forAll(cs, i) { cs[i]=pointMap[cs[i]];} cs.collapse();}label bcIDs[11]={-1, 0, 2, 4, -1, 5, -1, 6, 7, 8, 9};constexpr label nBCs=12;const word *kivaPatchTypes[nBCs]={ &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &symmetryPolyPatch::typeName, &wedgePolyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &symmetryPolyPatch::typeName, &oldCyclicPolyPatch::typeName};enum patchTypeNames{ PISTON, VALVE, LINER, CYLINDERHEAD, AXIS, WEDGE, INFLOW, OUTFLOW, PRESIN, PRESOUT, SYMMETRYPLANE, CYCLIC};const char *kivaPatchNames[nBCs]={ "piston", "valve", "liner", "cylinderHead", "axis", "wedge", "inflow", "outflow", "presin", "presout", "symmetryPlane", "cyclic"};List< DynamicList< face > > pFaces[nBCs]
volScalarField & e