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-2022 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 "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 
39 namespace Foam
40 {
41  defineTypeNameAndDebug(extrudePatchMesh, 0);
42 }
43 
44 
45 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
46 
47 Foam::extrudePatchMesh::extrudePatchMesh
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  (
59  regionName,
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 
74 Foam::extrudePatchMesh::extrudePatchMesh
75 (
76  const fvMesh& mesh,
77  const fvPatch& p,
78  const dictionary& dict,
79  const word& regionName,
80  polyPatchList& regionPatches
81 )
82 :
84 {
85  extrudeMesh(regionPatches);
86 }
87 
88 
89 Foam::extrudePatchMesh::extrudePatchMesh
90 (
91  const fvMesh& mesh,
92  const fvPatch& p,
93  const dictionary& dict,
94  const word& regionName,
95  const List<polyPatch*>& regionPatches
96 )
97 :
99 {
100  // Acquire ownership of the pointers
101  polyPatchList plist(const_cast<List<polyPatch*>&>(regionPatches));
102 
103  extrudeMesh(plist);
104 }
105 
106 
107 Foam::extrudePatchMesh::extrudePatchMesh
108 (
109  const fvMesh& mesh,
110  const fvPatch& p,
111  const dictionary& dict,
112  const word& regionName
113 )
114 :
116 {
117  polyPatchList regionPatches(3);
118  List<dictionary> dicts(regionPatches.size());
119  List<word> patchNames(regionPatches.size());
120  List<word> patchTypes(regionPatches.size());
121 
122  dicts[bottomPatchID] = dict_.subDict("bottomCoeffs");
123  dicts[sidePatchID] = dict_.subDict("sideCoeffs");
124  dicts[topPatchID] = dict_.subDict("topCoeffs");
125 
126  forAll(dicts, patchi)
127  {
128  dicts[patchi].readEntry("name", patchNames[patchi]);
129  dicts[patchi].readEntry("type", patchTypes[patchi]);
130  }
131 
132  forAll(regionPatches, patchi)
133  {
134  dictionary& patchDict = dicts[patchi];
135  patchDict.set("nFaces", 0);
136  patchDict.set("startFace", 0);
137 
138  regionPatches.set
139  (
140  patchi,
142  (
143  patchNames[patchi],
144  patchDict,
145  patchi,
146  mesh.boundaryMesh()
147  )
148  );
149  }
150 
151  extrudeMesh(regionPatches);
152 }
153 
154 
155 void Foam::extrudePatchMesh::extrudeMesh(polyPatchList& regionPatches)
156 {
157  if (this->boundaryMesh().empty())
158  {
159  const bool columnCells = dict_.get<bool>("columnCells");
160 
161  bitSet nonManifoldEdge(extrudedPatch_.nEdges());
162  for (label edgeI = 0; edgeI < extrudedPatch_.nInternalEdges(); edgeI++)
163  {
164  if (columnCells)
165  {
166  nonManifoldEdge.set(edgeI);
167  }
168  }
169 
170  autoPtr<extrudeModel> model_(extrudeModel::New(dict_));
171 
172  faceList pointGlobalRegions;
173  faceList pointLocalRegions;
174  labelList localToGlobalRegion;
175 
176  const primitiveFacePatch pp
177  (
178  extrudedPatch_, extrudedPatch_.points()
179  );
180 
182  (
183  this->globalData(),
184  pp,
185  nonManifoldEdge,
186  false,
187 
188  pointGlobalRegions,
189  pointLocalRegions,
190  localToGlobalRegion
191  );
192 
193 
194  // Per local region an originating point
195  labelList localRegionPoints(localToGlobalRegion.size());
196  forAll(pointLocalRegions, facei)
197  {
198  const face& f = extrudedPatch_.localFaces()[facei];
199  const face& pRegions = pointLocalRegions[facei];
200  forAll(pRegions, fp)
201  {
202  localRegionPoints[pRegions[fp]] = f[fp];
203  }
204  }
205 
206  // Calculate region normals by reducing local region normals
207  pointField localRegionNormals(localToGlobalRegion.size());
208  {
209  pointField localSum(localToGlobalRegion.size(), Zero);
210 
211  forAll(pointLocalRegions, facei)
212  {
213  const face& pRegions = pointLocalRegions[facei];
214  forAll(pRegions, fp)
215  {
216  label localRegionI = pRegions[fp];
217  localSum[localRegionI] +=
218  extrudedPatch_.faceNormals()[facei];
219  }
220  }
221 
222  Map<point> globalSum(2*localToGlobalRegion.size());
223 
224  forAll(localSum, localRegionI)
225  {
226  label globalRegionI = localToGlobalRegion[localRegionI];
227  globalSum.insert(globalRegionI, localSum[localRegionI]);
228  }
229 
230  // Reduce
231  Pstream::mapCombineReduce(globalSum, plusEqOp<point>());
232 
233  forAll(localToGlobalRegion, localRegionI)
234  {
235  label globalRegionI = localToGlobalRegion[localRegionI];
236  localRegionNormals[localRegionI] = globalSum[globalRegionI];
237  }
238  localRegionNormals /= mag(localRegionNormals);
239  }
240 
241 
242  // Per local region an extrusion direction
243  vectorField firstDisp(localToGlobalRegion.size());
244  forAll(firstDisp, regionI)
245  {
246  //const point& regionPt = regionCentres[regionI];
247  const point& regionPt = extrudedPatch_.points()
248  [
249  extrudedPatch_.meshPoints()
250  [
251  localRegionPoints[regionI]
252  ]
253  ];
254  const vector& n = localRegionNormals[regionI];
255  firstDisp[regionI] = model_()(regionPt, n, 1) - regionPt;
256  }
257 
258  const label nNewPatches = regionPatches.size();
259 
260  // Extrude engine
261  createShellMesh extruder
262  (
263  pp,
264  pointLocalRegions,
265  localRegionPoints
266  );
267  this->clearOut();
268  this->removeFvBoundary();
269  this->addFvPatches(regionPatches, true);
270 
271 
272  // At this point we have a valid mesh with 3 patches and zero cells.
273  // Determine:
274  // - per face the top and bottom patch (topPatchID, bottomPatchID)
275  // - per edge, per face on edge the side patch (edgePatches)
276  labelListList edgePatches(extrudedPatch_.nEdges());
277  forAll(edgePatches, edgeI)
278  {
279  const labelList& eFaces = extrudedPatch_.edgeFaces()[edgeI];
280 
281  if (eFaces.size() != 2 || nonManifoldEdge.test(edgeI))
282  {
283  edgePatches[edgeI].setSize(eFaces.size(), sidePatchID);
284  }
285  }
286 
287  polyTopoChange meshMod(nNewPatches);
288 
289  extruder.setRefinement
290  (
291  firstDisp, // first displacement
292  model_().expansionRatio(),
293  model_().nLayers(), // nLayers
294  labelList(extrudedPatch_.size(), topPatchID),
295  labelList(extrudedPatch_.size(), bottomPatchID),
296  edgePatches,
297  meshMod
298  );
299 
300  autoPtr<mapPolyMesh> map = meshMod.changeMesh
301  (
302  *this, // mesh to change
303  false // inflate
304  );
305 
306  // Update numbering on extruder.
307  extruder.updateMesh(map());
308 
309  this->setInstance(this->thisDb().time().constant());
310  this->write();
311  }
312 }
313 
314 
315 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:51
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:118
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:120
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:56
std::enable_if< std::is_same< bool, TypeT >::value, bool >::type set(const label i, bool val=true)
A bitSet::set() method for a list of bool.
Definition: List.H:463
wordList patchTypes(nPatches)
List< face > faceList
A List of faces.
Definition: faceListFwd.H:41
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:68
Mesh at a patch created on the fly. The following entry should be used on the field boundary dictiona...
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:453
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:413
dictionary()
Default construct, a top-level empty dictionary.
Definition: dictionary.C:68
PtrList< polyPatch > polyPatchList
Store lists of polyPatch as a PtrList.
Definition: polyPatch.H:56
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
Foam::word regionName(Foam::polyMesh::defaultRegion)
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
void setSize(const label n)
Alias for resize()
Definition: List.H:289
dynamicFvMesh & mesh
A class for handling words, derived from Foam::string.
Definition: word.H:63
wordList patchNames(nPatches)
Vector< scalar > vector
Definition: vector.H:57
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.
defineTypeNameAndDebug(combustionModel, 0)
labelList f(nPoints)
static autoPtr< extrudeModel > New(const dictionary &dict)
Select null constructed.
PrimitivePatch< List< face >, const pointField & > primitiveFacePatch
A PrimitivePatch with List storage for the faces, const reference for the point field.
vector point
Point is a vector.
Definition: point.H:37
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:79
const std::string patch
OpenFOAM patch number as a std::string.
label n
Field< vector > vectorField
Specialisation of Field<T> for vector.
static autoPtr< polyPatch > New(const word &patchType, const word &name, const label size, const label start, const label index, const polyBoundaryMesh &bm)
Return pointer to a new patch created on freestore from components.
Definition: polyPatchNew.C:28
List< label > labelList
A List of labels.
Definition: List.H:62
volScalarField & p
static void mapCombineReduce(Container &values, const CombineOp &cop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce) applying cop to inplace combine map values from different processo...
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:157