Loading...
Searching...
No Matches
snappyLayerDriver.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-2015 OpenFOAM Foundation
9 Copyright (C) 2017 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::snappyLayerDriver
29
30Description
31 All to do with adding layers
32
33SourceFiles
34 snappyLayerDriver.C
35 snappyLayerDriverOneByOne.C
36
37\*---------------------------------------------------------------------------*/
38
39#ifndef snappyLayerDriver_H
40#define snappyLayerDriver_H
41
42#include "meshRefinement.H"
43#include "scalarIOField.H"
44
45// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46
47namespace Foam
48{
49
50// Forward declaration of classes
51class removePoints;
52class pointSet;
53class motionSmoother;
55class faceSet;
56class layerParameters;
59/*---------------------------------------------------------------------------*\
60 Class snappyLayerDriver Declaration
61\*---------------------------------------------------------------------------*/
62
63class snappyLayerDriver
64{
65public:
66
67 // Public data types
68
69 //- Extrusion controls
71 {
72 NOEXTRUDE,
73 EXTRUDE,
76 };
77
78private:
79
80 // Private classes
81
82 //- Combine operator class to combine normal with other normal.
83 class nomalsCombine
84 {
85 public:
86
87 void operator()(vector& x, const vector& y) const
88 {
89 if (y != point::max)
90 {
91 if (x == point::max)
92 {
93 x = y;
94 }
95 else
96 {
97 x *= (x&y);
98 }
99 }
100 }
101 };
102
103
104 // Private data
105
106 //- Mesh+surface
107 meshRefinement& meshRefiner_;
108
109 //- From surface region to patch
110 const labelList globalToMasterPatch_;
111
112 //- From surface region to patch
113 const labelList globalToSlavePatch_;
114
115 //- Are we operating in test mode?
116 const bool dryRun_;
117
118
119 // Private Member Functions
120
121 // Layers
122
123 //- For debugging: Dump displacement to .obj files
124 static void dumpDisplacement
125 (
126 const fileName&,
128 const vectorField&,
129 const List<extrudeMode>&
130 );
131
132 //- Average point wise data to face wise
133 static tmp<scalarField> avgPointData
134 (
136 const scalarField& pointFld
137 );
138
139 //- Check that primitivePatch is not multiply connected.
140 // Collect non-manifold points in pointSet.
141 static void checkManifold
142 (
144 pointSet& nonManifoldPoints
145 );
146
147 //- Check that mesh outside is not multiply connected.
148 void checkMeshManifold() const;
149
150
151 // Static extrusion setup
152
153 //- Unset extrusion on point. Returns true if anything unset.
154 static bool unmarkExtrusion
155 (
156 const label patchPointi,
157 pointField& patchDisp,
158 labelList& patchNLayers,
159 List<extrudeMode>& extrudeStatus
160 );
161
162 //- Unset extrusion on face. Returns true if anything unset.
163 static bool unmarkExtrusion
164 (
165 const face& localFace,
166 pointField& patchDisp,
167 labelList& patchNLayers,
168 List<extrudeMode>& extrudeStatus
169 );
170
171 //- Truncate index in face
172 static label constrainFp(const label sz, const label fp);
173
174 //- Count common points between face and its neighbours
175 void countCommonPoints
176 (
178 const label facei,
179 Map<label>&
180 ) const;
181
182 //- Check if any common points form single string. Return
183 // false if not.
184 bool checkCommonOrder
185 (
186 const label nCommon,
187 const face& curFace,
188 const face& nbFace
189 ) const;
190
191 //- Check if any common points form single string; unmark
192 // points on face if not
193 void checkCommonOrder
194 (
196 const label facei,
197 const Map<label>& nCommonPoints,
198 pointField& patchDisp,
199 labelList& patchNLayers,
200 List<extrudeMode>& extrudeStatus
201 ) const;
202
203 //- Check if any common points form single string; unmark
204 // points on face if not
205 void handleNonStringConnected
206 (
208 pointField& patchDisp,
209 labelList& patchNLayers,
210 List<extrudeMode>& extrudeStatus
211 ) const;
212
213 //- No extrusion at non-manifold points.
214 void handleNonManifolds
215 (
217 const labelList& meshEdges,
218 const labelListList& edgeGlobalFaces,
219 pointField& patchDisp,
220 labelList& patchNLayers,
221 List<extrudeMode>& extrudeStatus
222 ) const;
223
224 //- No extrusion on feature edges. Assumes non-manifold
225 // edges already handled.
226 void handleFeatureAngle
227 (
229 const labelList& meshEdges,
230 const scalar minAngle,
231 pointField& patchDisp,
232 labelList& patchNLayers,
233 List<extrudeMode>& extrudeStatus
234 ) const;
235
236 //- No extrusion on warped faces
237 void handleWarpedFaces
238 (
240 const scalar faceRatio,
241 const boolList& relativeSizes,
242 const scalar edge0Len,
243 const labelList& cellLevel,
244 pointField& patchDisp,
245 labelList& patchNLayers,
246 List<extrudeMode>& extrudeStatus
247 ) const;
248
249 //- Determine the number of layers per point from the number of
250 // layers per surface.
251 static void setNumLayers
252 (
253 meshRefinement& meshRefiner,
254 const labelList& patchToNLayers,
255 const labelList& patchIDs,
257
258 labelList& patchNLayers,
259 List<extrudeMode>& extrudeStatus,
260 label& nIdealAddedCells
261 );
262
263 //- Determine number of layers per point; include static checks
264 //- on invalid extrusion (e.g. non-manifold)
265 label setPointNumLayers
266 (
267 const layerParameters& layerParams,
268
269 const labelList& numLayers,
270 const labelList& patchIDs,
272 const labelListList& edgeGlobalFaces,
273
274 vectorField& patchDisp,
275 labelList& patchNLayers,
276 List<extrudeMode>&
277 ) const;
278 autoPtr<externalDisplacementMeshMover> makeMeshMover
279 (
280 const layerParameters& layerParams,
281 const dictionary& motionDict,
282 const labelList& internalFaceZones,
283 const scalarIOField& minThickness,
284 pointVectorField& displacement
285 ) const;
286 void updatePatch
287 (
288 const labelList& patchIDs,
289 const mapPolyMesh& map,
290 autoPtr<indirectPrimitivePatch>& pp,
291 labelList& newToOldPatchPoints
292 ) const;
293
294
295 //- Helper function to make a pointVectorField with correct
296 // bcs for layer addition:
297 // - numLayers > 0 : fixedValue
298 // - numLayers == 0 : fixedValue (always zero)
299 // - processor : calculated (so free to move)
300 // - cyclic/wedge/symmetry : slip
301 // - other : slip
302 static tmp<pointVectorField> makeLayerDisplacementField
303 (
304 const pointMesh& pMesh,
305 const labelList& numLayers
306 );
307
308 //- Grow no-extrusion layer.
309 void growNoExtrusion
310 (
312 pointField& patchDisp,
313 labelList& patchNLayers,
314 List<extrudeMode>& extrudeStatus
315 ) const;
316
317 //- Re-merge points/faces on faceZones. Opposite of
318 //- dupFaceZonePoints above
319 void mergeFaceZonePoints
320 (
321 const labelList& pointToMaster,
322 labelList& cellNLayers,
323 scalarField& faceRealThickness,
324 scalarField& faceWantedThickness
325 );
326
327 //- Calculate pointwise wanted and minimum thickness.
328 // thickness: wanted thickness
329 // minthickness: when to give up and not extrude
330 // Gets per patch parameters and determine pp pointwise
331 // parameters.
332 void calculateLayerThickness
333 (
335 const labelList& patchIDs,
336 const layerParameters& layerParams,
337 const labelList& cellLevel,
338 const labelList& patchNLayers,
339 const scalar edge0Len,
340
341 scalarField& thickness,
342 scalarField& minThickness,
343 scalarField& expansionRatio
344 ) const;
345
346
347 // Extrusion execution
348
349 //- Synchronize displacement among coupled patches.
350 static void syncPatchDisplacement
351 (
352 const fvMesh& mesh,
354 const scalarField& minThickness,
355 pointField& patchDisp,
356 labelList& patchNLayers,
357 List<extrudeMode>& extrudeStatus
358 );
359
360 //- Get nearest point on surface to snap to
361 void getPatchDisplacement
362 (
364 const scalarField& thickness,
365 const scalarField& minThickness,
366 const scalarField& expansionRatio,
367
368 pointField& patchDisp,
369 labelList& patchNLayers,
370 List<extrudeMode>& extrudeStatus
371 ) const;
372
373 //- For truncateDisplacement: find strings of edges
374 bool sameEdgeNeighbour
375 (
376 const labelListList& globalEdgeFaces,
377 const label myGlobalFacei,
378 const label nbrGlobFacei,
379 const label edgeI
380 ) const;
381
382 //- For truncateDisplacement: find strings of edges
383 void getVertexString
384 (
386 const labelListList& globalEdgeFaces,
387 const label facei,
388 const label edgeI,
389 const label myGlobFacei,
390 const label nbrGlobFacei,
391 DynamicList<label>& vertices
392 ) const;
393
394 //- Truncates displacement
395 // - for all patchFaces in the faceset displacement gets set
396 // to zero
397 // - all displacement < minThickness gets set to zero
398 // - all non-consecutive extrusions get set to 0
399 label truncateDisplacement
400 (
401 const globalIndex& globalFaces,
402 const labelListList& edgeGlobalFaces,
404 const scalarField& minThickness,
405 const faceSet& illegalPatchFaces,
406 pointField& patchDisp,
407 labelList& patchNLayers,
408 List<extrudeMode>& extrudeStatus
409 ) const;
410
411 //- Setup layer information (at points and faces) to
412 // modify mesh topology in
413 // regions where layer mesh terminates. Guarantees an
414 // optional slow decreasing of the number of layers.
415 // Returns the number of layers per face and per point
416 // to go into the actual layer addition engine.
417 void setupLayerInfoTruncation
418 (
420 const labelList& patchNLayers,
421 const List<extrudeMode>& extrudeStatus,
422 const label nBufferCellsNoExtrude,
423 labelList& nPatchPointLayers,
424 labelList& nPatchFaceLayers
425 ) const;
426
427 //- Does any of the cells use a face from faces?
428 static bool cellsUseFace
429 (
430 const polyMesh& mesh,
431 const labelList& cellLabels,
432 const labelHashSet& faces
433 );
434
435 //- Checks the newly added cells and locally unmarks points
436 // so they will not get extruded next time round. Returns
437 // global number of unmarked points (0 if all was fine)
438 static label checkAndUnmark
439 (
440 const addPatchCellLayer& addLayer,
441 const dictionary& motionDict,
442 const bool additionalReporting,
443 const List<labelPair>& baffles,
445 const fvMesh&,
446
447 pointField& patchDisp,
448 labelList& patchNLayers,
449 List<extrudeMode>& extrudeStatus
450 );
451
452 //- Count global number of extruded faces
453 static label countExtrusion
454 (
456 const List<extrudeMode>& extrudeStatus
457 );
458
459 //- After adding to mesh get the new baffles
460 static List<labelPair> getBafflesOnAddedMesh
461 (
462 const polyMesh& mesh,
463 const labelList& newToOldFaces,
464 const List<labelPair>& baffles
465 );
466
467 //- Collect layer faces and layer cells into bools
468 // for ease of handling
469 static void getLayerCellsFaces
470 (
471 const polyMesh&,
472 const addPatchCellLayer&,
473 const scalarField& oldRealThickness,
474
475 labelList& cellStatus,
476 scalarField& faceRealThickness
477 );
478
479 //- Print layer coverage table
480 void printLayerData
481 (
482 const fvMesh& mesh,
483 const labelList& patchIDs,
484 const labelList& cellNLayers,
485 const scalarField& faceWantedThickness,
486 const scalarField& faceRealThickness,
487 const layerParameters& layerParams
488 ) const;
489
490 //- Write cellSet,faceSet for layers
491 bool writeLayerSets
492 (
493 const fvMesh& mesh,
494 const labelList& cellNLayers,
495 const scalarField& faceRealThickness
496 ) const;
497
498 //- Write volFields,cellSet,faceSet for layers depending
499 // on write level
500 bool writeLayerData
501 (
502 const fvMesh& mesh,
503 const labelList& patchIDs,
504 const labelList& cellNLayers,
505 const scalarField& faceWantedThickness,
506 const scalarField& faceRealThickness
507 ) const;
508
509
510 // Mesh shrinking (to create space for layers)
511
512 //- Average field (over all subset of mesh points) by
513 // summing contribution from edges. Global parallel since only
514 // does master edges for coupled edges.
515 template<class Type>
516 static void averageNeighbours
517 (
518 const polyMesh& mesh,
519 const bitSet& isMasterEdge,
520 const labelList& meshEdges,
521 const labelList& meshPoints,
522 const edgeList& edges,
523 const scalarField& invSumWeight,
524 const Field<Type>& data,
525 Field<Type>& average
526 );
527
528 //- Calculate inverse sum of edge weights (currently always 1.0)
529 void sumWeights
530 (
531 const bitSet& isMasterEdge,
532 const labelList& meshEdges,
533 const labelList& meshPoints,
534 const edgeList& edges,
535 scalarField& invSumWeight
536 ) const;
537
538 //- Smooth scalar field on patch
539 void smoothField
540 (
541 const motionSmoother& meshMover,
542 const bitSet& isMasterPoint,
543 const bitSet& isMasterEdge,
544 const labelList& meshEdges,
545 const scalarField& fieldMin,
546 const label nSmoothDisp,
548 ) const;
549
550 //- Smooth normals on patch.
551 void smoothPatchNormals
552 (
553 const motionSmoother& meshMover,
554 const bitSet& isMasterPoint,
555 const bitSet& isMasterEdge,
556 const labelList& meshEdges,
557 const label nSmoothDisp,
558 pointField& normals
559 ) const;
560
561 //- Smooth normals in interior.
562 void smoothNormals
563 (
564 const label nSmoothDisp,
565 const bitSet& isMasterPoint,
566 const bitSet& isMasterEdge,
567 const labelList& fixedPoints,
568 pointVectorField& normals
569 ) const;
570
571 //- Stop layer growth where mesh wraps around edge with a
572 // large feature angle
573 void handleFeatureAngleLayerTerminations
574 (
575 const scalar minCos,
576 const bitSet& isMasterPoint,
578 const labelList& meshEdges,
579
580 List<extrudeMode>& extrudeStatus,
581 pointField& patchDisp,
582 labelList& patchNLayers,
583 label& nPointCounter
584 ) const;
585
586 //- Find isolated islands (points, edges and faces and
587 // layer terminations)
588 // in the layer mesh and stop any layer growth at these points.
589 void findIsolatedRegions
590 (
591 const scalar minCosLayerTermination,
592 const bitSet& isMasterPoint,
593 const bitSet& isMasterEdge,
595 const labelList& meshEdges,
596 const scalarField& minThickness,
597
598 List<extrudeMode>& extrudeStatus,
599 pointField& patchDisp,
600 labelList& patchNLayers
601 ) const;
602
603
604 //- No copy construct
605 snappyLayerDriver(const snappyLayerDriver&) = delete;
606
607 //- No copy assignment
608 void operator=(const snappyLayerDriver&) = delete;
609
610
611public:
612
613 //- Runtime type information
614 ClassName("snappyLayerDriver");
615
616 // Constructors
617
618 //- Construct from components
619 snappyLayerDriver
620 (
621 meshRefinement& meshRefiner,
622 const labelList& globalToMasterPatch,
623 const labelList& globalToSlavePatch,
624 const bool dryRun = false
625 );
626
627
628 // Member Functions
629
630 //- Merge patch faces on same cell. Return total number of
631 //- faces/edges changed.
633 (
634 const layerParameters& layerParams,
635 const dictionary& motionDict,
636 const meshRefinement::FaceMergeType mergeType
637 );
638
639 void addLayers
640 (
641 const layerParameters& layerParams,
642 const label nLayerIter,
643
644 const dictionary& motionDict,
645 const label nRelaxedIter,
646 const label nAllowableErrors,
647
648 const labelList& patchIDs,
649 const labelList& internalFaceZones,
650 const List<labelPair>& baffles,
651 const labelList& numLayers,
652 const label nIdealTotAddedCells,
653
654 const globalIndex& globalFaces,
656 const labelListList& edgeGlobalFaces,
657 const labelList& edgePatchID,
658 const labelList& edgeZoneID,
659 const boolList& edgeFlip,
660 const labelList& inflateFaceID,
661 const scalarField& thickness,
662 const scalarIOField& minThickness,
663 const scalarField& expansionRatio,
664 vectorField& patchDisp,
665 labelList& patchNLayers,
666 List<extrudeMode>& extrudeStatus,
667 polyTopoChange& savedMeshMod,
668 labelList& cellNLayers,
669 scalarField& faceRealThickness
670 );
671
672 //- Add cell layers
673 void addLayers
674 (
675 const layerParameters& layerParams,
676 const dictionary& motionDict,
677 const labelList& patchIDs,
678 const label nAllowableErrors,
679 decompositionMethod& decomposer,
680 fvMeshDistribute& distributor
681 );
682
683 //- For debugging. Can be removed.
685 (
686 const layerParameters& layerParams,
687 const dictionary& motionDict,
688 const labelList& patchIDs,
689 const label nAllowableErrors,
690 decompositionMethod& decomposer,
691 fvMeshDistribute& distributor
692 );
693
694 //- Add layers according to the dictionary settings
695 void doLayers
696 (
697 const dictionary& shrinkDict,
698 const dictionary& motionDict,
699 const layerParameters& layerParams,
700 const meshRefinement::FaceMergeType mergeType,
701 const bool preBalance, // balance before adding?
702 decompositionMethod& decomposer,
703 fvMeshDistribute& distributor
704 );
705
706 //- Helper: see what zones and patches edges should be extruded into
707 static void determineSidePatches
708 (
709 meshRefinement& meshRefiner,
710 const globalIndex& globalFaces,
711 const labelListList& edgeGlobalFaces,
713
714 labelList& edgePatchID,
715 labelList& edgeZoneID,
716 boolList& edgeFlip,
717 labelList& inflateFaceID
718 );
719
720 //- Duplicate points on faceZones with layers. Re-used when adding
721 //- buffer layers. Can be made private again once multi-side
722 //- layer addition working.
723 static autoPtr<mapPolyMesh> dupFaceZonePoints
724 (
725 meshRefinement& meshRefiner,
726 const labelList& patchIDs, // patch indices
727 const labelList& numLayers, // number of layers per patch
728 List<labelPair> baffles,
729 labelList& pointToMaster
730 );
731
732 //- Map numbering after adding cell layers
733 static void mapFaceZonePoints
734 (
735 meshRefinement& meshRefiner,
736 const mapPolyMesh& map,
737 labelPairList& baffles,
738 labelList& pointToMaster
739 );
740};
741
742
743// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
744
745} // End namespace Foam
746
747// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
748
749#ifdef NoRepository
751#endif
752
753// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
754
755#endif
756
757// ************************************************************************* //
scalar y
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
labelList patchIDs
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
static const Form max
Adds layers of cells to outside (or inside) of polyMesh. Can optionally create stand-alone extruded m...
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
Abstract base class for domain decomposition.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Virtual base class for mesh movers with externally provided displacement field giving the boundary co...
A list of face labels.
Definition faceSet.H:50
Sends/receives parts of mesh+fvfields to neighbouring processors. Used in load balancing.
Calculates a non-overlapping list of offsets based on an input size (eg, number of cells) from differ...
Definition globalIndex.H:77
Simple container to keep together layer specific information.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Helper class which maintains intersections of (changing) mesh with (static) surfaces.
FaceMergeType
Enumeration for what to do with co-planar patch faces on a single.
Given a displacement moves the mesh by scaling the displacement back until there are no more mesh err...
A set of point labels.
Definition pointSet.H:50
Direct mesh changes based on v1.3 polyTopoChange syntax.
Removes selected points from mesh and updates faces using these points.
void addLayersSinglePass(const layerParameters &layerParams, const dictionary &motionDict, const labelList &patchIDs, const label nAllowableErrors, decompositionMethod &decomposer, fvMeshDistribute &distributor)
For debugging. Can be removed.
static autoPtr< mapPolyMesh > dupFaceZonePoints(meshRefinement &meshRefiner, const labelList &patchIDs, const labelList &numLayers, List< labelPair > baffles, labelList &pointToMaster)
Duplicate points on faceZones with layers. Re-used when adding buffer layers. Can be made private aga...
void addLayers(const layerParameters &layerParams, const label nLayerIter, const dictionary &motionDict, const label nRelaxedIter, const label nAllowableErrors, const labelList &patchIDs, const labelList &internalFaceZones, const List< labelPair > &baffles, const labelList &numLayers, const label nIdealTotAddedCells, const globalIndex &globalFaces, indirectPrimitivePatch &pp, const labelListList &edgeGlobalFaces, const labelList &edgePatchID, const labelList &edgeZoneID, const boolList &edgeFlip, const labelList &inflateFaceID, const scalarField &thickness, const scalarIOField &minThickness, const scalarField &expansionRatio, vectorField &patchDisp, labelList &patchNLayers, List< extrudeMode > &extrudeStatus, polyTopoChange &savedMeshMod, labelList &cellNLayers, scalarField &faceRealThickness)
static void determineSidePatches(meshRefinement &meshRefiner, const globalIndex &globalFaces, const labelListList &edgeGlobalFaces, const indirectPrimitivePatch &pp, labelList &edgePatchID, labelList &edgeZoneID, boolList &edgeFlip, labelList &inflateFaceID)
Helper: see what zones and patches edges should be extruded into.
static void mapFaceZonePoints(meshRefinement &meshRefiner, const mapPolyMesh &map, labelPairList &baffles, labelList &pointToMaster)
Map numbering after adding cell layers.
void doLayers(const dictionary &shrinkDict, const dictionary &motionDict, const layerParameters &layerParams, const meshRefinement::FaceMergeType mergeType, const bool preBalance, decompositionMethod &decomposer, fvMeshDistribute &distributor)
Add layers according to the dictionary settings.
label mergePatchFacesUndo(const layerParameters &layerParams, const dictionary &motionDict, const meshRefinement::FaceMergeType mergeType)
Merge patch faces on same cell. Return total number of faces/edges changed.
ClassName("snappyLayerDriver")
Runtime type information.
extrudeMode
Extrusion controls.
@ NOEXTRUDE
Do not extrude. No layers added.
#define ClassName(TypeNameString)
Add typeName information from argument TypeNameString to a class.
Definition className.H:74
rDeltaTY field()
dynamicFvMesh & mesh
Namespace for OpenFOAM.
List< edge > edgeList
List of edge.
Definition edgeList.H:32
List< labelPair > labelPairList
List of labelPair.
Definition labelPair.H:33
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
GeometricField< vector, pointPatchField, pointMesh > pointVectorField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
PrimitivePatch< IndirectList< face >, const pointField & > indirectPrimitivePatch
A PrimitivePatch with an IndirectList for the faces, const reference for the point field.
pointField vertices(const blockVertexList &bvl)
Field< vector > vectorField
Specialisation of Field<T> for vector.
List< bool > boolList
A List of bools.
Definition List.H:60
vectorField pointField
pointField is a vectorField.
Vector< scalar > vector
Definition vector.H:57
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &f1, const label comm)
IOField< scalar > scalarIOField
IO for a Field of scalar.