extrudePatchMesh.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) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2020 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
27\*---------------------------------------------------------------------------*/
28
29#include "extrudePatchMesh.H"
30
31#include "createShellMesh.H"
32#include "polyTopoChange.H"
33#include "wallPolyPatch.H"
34#include "emptyPolyPatch.H"
35#include "wedgePolyPatch.H"
36
37// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38
39namespace Foam
40{
42}
43
44
45// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
46
48(
49 const word& regionName,
50 const fvMesh& mesh,
51 const fvPatch& p,
52 const dictionary& dict
53)
54:
55 fvMesh
56 (
57 IOobject
58 (
60 mesh.facesInstance(),
61 mesh,
62 IOobject::READ_IF_PRESENT,
63 IOobject::NO_WRITE,
64 true
65 ),
66 Zero,
67 false
68 ),
69 extrudedPatch_(p.patch()),
70 dict_(dict)
71{}
72
73
75(
76 const fvMesh& mesh,
77 const fvPatch& p,
78 const dictionary& dict,
79 const word& regionName,
80 const List<polyPatch*>& regionPatches
81)
82:
84{
85 extrudeMesh(regionPatches);
86}
87
88
90(
91 const fvMesh& mesh,
92 const fvPatch& p,
93 const dictionary& dict,
94 const word& regionName
95)
96:
98{
99 List<polyPatch*> regionPatches(3);
100 List<word> patchNames(regionPatches.size());
101 List<word> patchTypes(regionPatches.size());
102 PtrList<dictionary> dicts(regionPatches.size());
103
104 forAll(dicts, patchi)
105 {
106 if (!dicts.set(patchi))
107 {
108 dicts.set(patchi, new dictionary());
109 }
110 }
111
112 dicts[bottomPatchID] = dict_.subDict("bottomCoeffs");
113 dicts[sidePatchID] = dict_.subDict("sideCoeffs");
114 dicts[topPatchID] = dict_.subDict("topCoeffs");
115
116 forAll(dicts, patchi)
117 {
118 dicts[patchi].readEntry("name", patchNames[patchi]);
119 dicts[patchi].readEntry("type", patchTypes[patchi]);
120 }
121
122 forAll(regionPatches, patchi)
123 {
124 dictionary& patchDict = dicts[patchi];
125 patchDict.set("nFaces", 0);
126 patchDict.set("startFace", 0);
127
128 regionPatches[patchi] = polyPatch::New
129 (
130 patchNames[patchi],
131 patchDict,
132 patchi,
134 ).ptr();
135 }
136
137 extrudeMesh(regionPatches);
138}
139
140
141void Foam::extrudePatchMesh::extrudeMesh(const List<polyPatch*>& regionPatches)
142{
143 if (this->boundaryMesh().size() == 0)
144 {
145 const bool columnCells = dict_.get<bool>("columnCells");
146
147 bitSet nonManifoldEdge(extrudedPatch_.nEdges());
148 for (label edgeI = 0; edgeI < extrudedPatch_.nInternalEdges(); edgeI++)
149 {
150 if (columnCells)
151 {
152 nonManifoldEdge.set(edgeI);
153 }
154 }
155
156 autoPtr<extrudeModel> model_(extrudeModel::New(dict_));
157
158 faceList pointGlobalRegions;
159 faceList pointLocalRegions;
160 labelList localToGlobalRegion;
161
162 const primitiveFacePatch pp
163 (
164 extrudedPatch_, extrudedPatch_.points()
165 );
166
168 (
169 this->globalData(),
170 pp,
171 nonManifoldEdge,
172 false,
173
174 pointGlobalRegions,
175 pointLocalRegions,
176 localToGlobalRegion
177 );
178
179
180 // Per local region an originating point
181 labelList localRegionPoints(localToGlobalRegion.size());
182 forAll(pointLocalRegions, facei)
183 {
184 const face& f = extrudedPatch_.localFaces()[facei];
185 const face& pRegions = pointLocalRegions[facei];
186 forAll(pRegions, fp)
187 {
188 localRegionPoints[pRegions[fp]] = f[fp];
189 }
190 }
191
192 // Calculate region normals by reducing local region normals
193 pointField localRegionNormals(localToGlobalRegion.size());
194 {
195 pointField localSum(localToGlobalRegion.size(), Zero);
196
197 forAll(pointLocalRegions, facei)
198 {
199 const face& pRegions = pointLocalRegions[facei];
200 forAll(pRegions, fp)
201 {
202 label localRegionI = pRegions[fp];
203 localSum[localRegionI] +=
204 extrudedPatch_.faceNormals()[facei];
205 }
206 }
207
208 Map<point> globalSum(2*localToGlobalRegion.size());
209
210 forAll(localSum, localRegionI)
211 {
212 label globalRegionI = localToGlobalRegion[localRegionI];
213 globalSum.insert(globalRegionI, localSum[localRegionI]);
214 }
215
216 // Reduce
217 Pstream::mapCombineAllGather(globalSum, plusEqOp<point>());
218
219 forAll(localToGlobalRegion, localRegionI)
220 {
221 label globalRegionI = localToGlobalRegion[localRegionI];
222 localRegionNormals[localRegionI] = globalSum[globalRegionI];
223 }
224 localRegionNormals /= mag(localRegionNormals);
225 }
226
227
228 // Per local region an extrusion direction
229 vectorField firstDisp(localToGlobalRegion.size());
230 forAll(firstDisp, regionI)
231 {
232 //const point& regionPt = regionCentres[regionI];
233 const point& regionPt = extrudedPatch_.points()
234 [
235 extrudedPatch_.meshPoints()
236 [
237 localRegionPoints[regionI]
238 ]
239 ];
240 const vector& n = localRegionNormals[regionI];
241 firstDisp[regionI] = model_()(regionPt, n, 1) - regionPt;
242 }
243
244 const label nNewPatches = regionPatches.size();
245
246 // Extrude engine
247 createShellMesh extruder
248 (
249 pp,
250 pointLocalRegions,
251 localRegionPoints
252 );
253 this->clearOut();
254 this->removeFvBoundary();
255 this->addFvPatches(regionPatches, true);
256
257
258 // At this point we have a valid mesh with 3 patches and zero cells.
259 // Determine:
260 // - per face the top and bottom patch (topPatchID, bottomPatchID)
261 // - per edge, per face on edge the side patch (edgePatches)
262 labelListList edgePatches(extrudedPatch_.nEdges());
263 forAll(edgePatches, edgeI)
264 {
265 const labelList& eFaces = extrudedPatch_.edgeFaces()[edgeI];
266
267 if (eFaces.size() != 2 || nonManifoldEdge.test(edgeI))
268 {
269 edgePatches[edgeI].setSize(eFaces.size(), sidePatchID);
270 }
271 }
272
273 polyTopoChange meshMod(nNewPatches);
274
275 extruder.setRefinement
276 (
277 firstDisp, // first displacement
278 model_().expansionRatio(),
279 model_().nLayers(), // nLayers
280 labelList(extrudedPatch_.size(), topPatchID),
281 labelList(extrudedPatch_.size(), bottomPatchID),
282 edgePatches,
283 meshMod
284 );
285
286 autoPtr<mapPolyMesh> map = meshMod.changeMesh
287 (
288 *this, // mesh to change
289 false // inflate
290 );
291
292 // Update numbering on extruder.
293 extruder.updateMesh(map());
294
295 this->setInstance(this->thisDb().time().constant());
296 this->write();
297 }
298}
299
300
301// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
label n
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:77
void setSize(const label n)
Alias for resize()
Definition: List.H:218
static void mapCombineAllGather(const List< commsStruct > &comms, Container &values, const CombineOp &cop, const int tag, const label comm)
After completion all processors have the same data.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
const T * set(const label i) const
Definition: PtrList.H:138
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:66
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Definition: boundaryMesh.H:63
static void calcPointRegions(const globalMeshData &globalData, const primitiveFacePatch &patch, const bitSet &nonManifoldEdge, const bool syncNonCollocated, faceList &pointGlobalRegions, faceList &pointLocalRegions, labelList &localToGlobalRegion)
Helper: calculate point regions. The point region is the.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:460
entry * set(entry *entryPtr)
Assign a new entry, overwriting any existing entry.
Definition: dictionary.C:780
Mesh at a patch created on the fly. The following entry should be used on the field boundary dictiona...
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
constant condensation/saturation model.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:456
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
volScalarField & p
dynamicFvMesh & mesh
Foam::word regionName(Foam::polyMesh::defaultRegion)
Namespace for OpenFOAM.
List< label > labelList
A List of labels.
Definition: List.H:66
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
vector point
Point is a vector.
Definition: point.H:43
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
Field< vector > vectorField
Specialisation of Field<T> for vector.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
PrimitivePatch< List< face >, const pointField & > primitiveFacePatch
A PrimitivePatch with List storage for the faces, const reference for the point field.
List< face > faceList
A List of faces.
Definition: faceListFwd.H:47
runTime write()
wordList patchTypes(nPatches)
wordList patchNames(nPatches)
labelList f(nPoints)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333