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 -------------------------------------------------------------------------------
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  const List<polyPatch*>& 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 )
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 
141 void 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::mapCombineGather(globalSum, plusEqOp<point>());
218  Pstream::mapCombineScatter(globalSum);
219 
220  forAll(localToGlobalRegion, localRegionI)
221  {
222  label globalRegionI = localToGlobalRegion[localRegionI];
223  localRegionNormals[localRegionI] = globalSum[globalRegionI];
224  }
225  localRegionNormals /= mag(localRegionNormals);
226  }
227 
228 
229  // Per local region an extrusion direction
230  vectorField firstDisp(localToGlobalRegion.size());
231  forAll(firstDisp, regionI)
232  {
233  //const point& regionPt = regionCentres[regionI];
234  const point& regionPt = extrudedPatch_.points()
235  [
236  extrudedPatch_.meshPoints()
237  [
238  localRegionPoints[regionI]
239  ]
240  ];
241  const vector& n = localRegionNormals[regionI];
242  firstDisp[regionI] = model_()(regionPt, n, 1) - regionPt;
243  }
244 
245  const label nNewPatches = regionPatches.size();
246 
247  // Extrude engine
248  createShellMesh extruder
249  (
250  pp,
251  pointLocalRegions,
252  localRegionPoints
253  );
254  this->clearOut();
255  this->removeFvBoundary();
256  this->addFvPatches(regionPatches, true);
257 
258 
259  // At this point we have a valid mesh with 3 patches and zero cells.
260  // Determine:
261  // - per face the top and bottom patch (topPatchID, bottomPatchID)
262  // - per edge, per face on edge the side patch (edgePatches)
263  labelListList edgePatches(extrudedPatch_.nEdges());
264  forAll(edgePatches, edgeI)
265  {
266  const labelList& eFaces = extrudedPatch_.edgeFaces()[edgeI];
267 
268  if (eFaces.size() != 2 || nonManifoldEdge.test(edgeI))
269  {
270  edgePatches[edgeI].setSize(eFaces.size(), sidePatchID);
271  }
272  }
273 
274  polyTopoChange meshMod(nNewPatches);
275 
276  extruder.setRefinement
277  (
278  firstDisp, // first displacement
279  model_().expansionRatio(),
280  model_().nLayers(), // nLayers
281  labelList(extrudedPatch_.size(), topPatchID),
282  labelList(extrudedPatch_.size(), bottomPatchID),
283  edgePatches,
284  meshMod
285  );
286 
287  autoPtr<mapPolyMesh> map = meshMod.changeMesh
288  (
289  *this, // mesh to change
290  false // inflate
291  );
292 
293  // Update numbering on extruder.
294  extruder.updateMesh(map());
295 
296  this->setInstance(this->thisDb().time().constant());
297  this->write();
298  }
299 }
300 
301 
302 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
Foam::Pstream::mapCombineGather
static void mapCombineGather(const List< commsStruct > &comms, Container &Values, const CombineOp &cop, const int tag, const label comm)
Definition: combineGatherScatter.C:551
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::bitSet
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:63
Foam::primitiveFacePatch
PrimitivePatch< List< face >, const pointField & > primitiveFacePatch
A PrimitivePatch with List storage for the faces, const reference for the point field.
Definition: primitiveFacePatch.H:51
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
wedgePolyPatch.H
polyTopoChange.H
wallPolyPatch.H
Foam::bitSet::set
void set(const bitSet &bitset)
Set specified bits from another bitset.
Definition: bitSetI.H:574
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:444
Foam::dictionary::set
entry * set(entry *entryPtr)
Assign a new entry, overwriting any existing entry.
Definition: dictionary.C:780
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
patchTypes
wordList patchTypes(nPatches)
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
n
label n
Definition: TABSMDCalcMethod2.H:31
regionName
Foam::word regionName
Definition: createNamedDynamicFvMesh.H:1
createShellMesh.H
Foam::List::setSize
void setSize(const label n)
Alias for resize()
Definition: List.H:222
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::extrudeModel::New
static autoPtr< extrudeModel > New(const dictionary &dict)
Select null constructed.
Definition: extrudeModelNew.C:34
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:59
patchNames
wordList patchNames(nPatches)
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::polyPatch::New
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:35
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
emptyPolyPatch.H
Foam::labelListList
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
Foam::Pstream::mapCombineScatter
static void mapCombineScatter(const List< commsStruct > &comms, Container &Values, const int tag, const label comm)
Scatter data. Reverse of combineGather.
Definition: combineGatherScatter.C:666
f
labelList f(nPoints)
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::faceList
List< face > faceList
A List of faces.
Definition: faceListFwd.H:47
Foam::createShellMesh::calcPointRegions
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.
Definition: createShellMesh.C:157
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:63
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::extrudePatchMesh
Mesh at a patch created on the fly. The following entry should be used on the field boundary dictiona...
Definition: extrudePatchMesh.H:93
Foam::vtk::write
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
Definition: foamVtkOutputTemplates.C:36
Foam::boundaryMesh
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Definition: boundaryMesh.H:62
Foam::point
vector point
Point is a vector.
Definition: point.H:43
constant
constant condensation/saturation model.
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
extrudePatchMesh.H