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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "extrudePatchMesh.H"
29 
30 #include "createShellMesh.H"
31 #include "polyTopoChange.H"
32 #include "wallPolyPatch.H"
33 #include "emptyPolyPatch.H"
34 #include "wedgePolyPatch.H"
35 
36 
37 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41 
42 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43 
44 defineTypeNameAndDebug(extrudePatchMesh, 0);
45 
46 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47 
49 (
50  const fvMesh& mesh,
51  const fvPatch& patch,
52  const dictionary& dict,
53  const word regionName,
54  const List<polyPatch*>& regionPatches
55 )
56 :
57  fvMesh
58  (
59  IOobject
60  (
61  regionName,
63  mesh,
66  true
67  ),
68  Zero,
69  false
70  ),
71  extrudedPatch_(patch.patch()),
72  dict_(dict)
73 {
74  extrudeMesh(regionPatches);
75 }
76 
77 
79 (
80  const fvMesh& mesh,
81  const fvPatch& patch,
82  const dictionary& dict,
83  const word regionName
84 )
85 :
86  fvMesh
87  (
88  IOobject
89  (
90  regionName,
92  mesh,
95  true
96  ),
97  Zero,
98  false
99  ),
100  extrudedPatch_(patch.patch()),
101  dict_(dict)
102 {
103  List<polyPatch*> regionPatches(3);
104  List<word> patchNames(regionPatches.size());
105  List<word> patchTypes(regionPatches.size());
106  PtrList<dictionary> dicts(regionPatches.size());
107 
108  forAll(dicts, patchi)
109  {
110  if (!dicts.set(patchi))
111  {
112  dicts.set(patchi, new dictionary());
113  }
114  }
115 
116  dicts[bottomPatchID] = dict_.subDict("bottomCoeffs");
117  dicts[sidePatchID] = dict_.subDict("sideCoeffs");
118  dicts[topPatchID] = dict_.subDict("topCoeffs");
119 
120  forAll(dicts, patchi)
121  {
122  dicts[patchi].readEntry("name", patchNames[patchi]);
123  dicts[patchi].readEntry("type", patchTypes[patchi]);
124  }
125 
126  forAll(regionPatches, patchi)
127  {
128  dictionary& patchDict = dicts[patchi];
129  patchDict.set("nFaces", 0);
130  patchDict.set("startFace", 0);
131 
132  regionPatches[patchi] = polyPatch::New
133  (
134  patchNames[patchi],
135  patchDict,
136  patchi,
138  ).ptr();
139 
140  }
141 
142  extrudeMesh(regionPatches);
143 }
144 
145 
146 void extrudePatchMesh::extrudeMesh(const List<polyPatch*>& regionPatches)
147 {
148  if (this->boundaryMesh().size() == 0)
149  {
150  const bool columnCells = dict_.get<bool>("columnCells");
151 
152  bitSet nonManifoldEdge(extrudedPatch_.nEdges());
153  for (label edgeI = 0; edgeI < extrudedPatch_.nInternalEdges(); edgeI++)
154  {
155  if (columnCells)
156  {
157  nonManifoldEdge.set(edgeI);
158  }
159  }
160 
161  autoPtr<extrudeModel> model_(extrudeModel::New(dict_));
162 
163  faceList pointGlobalRegions;
164  faceList pointLocalRegions;
165  labelList localToGlobalRegion;
166 
167  const primitiveFacePatch pp
168  (
169  extrudedPatch_, extrudedPatch_.points()
170  );
171 
173  (
174  this->globalData(),
175  pp,
176  nonManifoldEdge,
177  false,
178 
179  pointGlobalRegions,
180  pointLocalRegions,
181  localToGlobalRegion
182  );
183 
184 
185  // Per local region an originating point
186  labelList localRegionPoints(localToGlobalRegion.size());
187  forAll(pointLocalRegions, facei)
188  {
189  const face& f = extrudedPatch_.localFaces()[facei];
190  const face& pRegions = pointLocalRegions[facei];
191  forAll(pRegions, fp)
192  {
193  localRegionPoints[pRegions[fp]] = f[fp];
194  }
195  }
196 
197  // Calculate region normals by reducing local region normals
198  pointField localRegionNormals(localToGlobalRegion.size());
199  {
200  pointField localSum(localToGlobalRegion.size(), Zero);
201 
202  forAll(pointLocalRegions, facei)
203  {
204  const face& pRegions = pointLocalRegions[facei];
205  forAll(pRegions, fp)
206  {
207  label localRegionI = pRegions[fp];
208  localSum[localRegionI] +=
209  extrudedPatch_.faceNormals()[facei];
210  }
211  }
212 
213  Map<point> globalSum(2*localToGlobalRegion.size());
214 
215  forAll(localSum, localRegionI)
216  {
217  label globalRegionI = localToGlobalRegion[localRegionI];
218  globalSum.insert(globalRegionI, localSum[localRegionI]);
219  }
220 
221  // Reduce
222  Pstream::mapCombineGather(globalSum, plusEqOp<point>());
223  Pstream::mapCombineScatter(globalSum);
224 
225  forAll(localToGlobalRegion, localRegionI)
226  {
227  label globalRegionI = localToGlobalRegion[localRegionI];
228  localRegionNormals[localRegionI] = globalSum[globalRegionI];
229  }
230  localRegionNormals /= mag(localRegionNormals);
231  }
232 
233 
234  // Per local region an extrusion direction
235  vectorField firstDisp(localToGlobalRegion.size());
236  forAll(firstDisp, regionI)
237  {
238  //const point& regionPt = regionCentres[regionI];
239  const point& regionPt = extrudedPatch_.points()
240  [
241  extrudedPatch_.meshPoints()
242  [
243  localRegionPoints[regionI]
244  ]
245  ];
246  const vector& n = localRegionNormals[regionI];
247  firstDisp[regionI] = model_()(regionPt, n, 1) - regionPt;
248  }
249 
250 
251  // Extrude engine
252  createShellMesh extruder
253  (
254  pp,
255  pointLocalRegions,
256  localRegionPoints
257  );
258 /*
259  List<polyPatch*> regionPatches(3);
260  List<word> patchNames(regionPatches.size());
261  List<word> patchTypes(regionPatches.size());
262  PtrList<dictionary> dicts(regionPatches.size());
263 
264  forAll(dicts, patchi)
265  {
266  if (!dicts.set(patchi))
267  {
268  dicts.set(patchi, new dictionary());
269  }
270  }
271 
272  dicts[bottomPatchID] = dict_.subDict("bottomCoeffs");
273  dicts[sidePatchID] = dict_.subDict("sideCoeffs");
274  dicts[topPatchID] = dict_.subDict("topCoeffs");
275 
276  forAll(dicts, patchi)
277  {
278  dicts[patchi].readEntry("name", patchNames[patchi]);
279  dicts[patchi].readEntry("type", patchTypes[patchi]);
280  }
281 
282  forAll(regionPatches, patchi)
283  {
284  dictionary& patchDict = dicts[patchi];
285  patchDict.set("nFaces", 0);
286  patchDict.set("startFace", 0);
287 
288  regionPatches[patchi] = polyPatch::New
289  (
290  patchNames[patchi],
291  patchDict,
292  patchi,
293  mesh.boundaryMesh()
294  ).ptr();
295 
296  }
297 */
298  this->clearOut();
299  this->removeFvBoundary();
300  this->addFvPatches(regionPatches, true);
301 
302 
303  // At this point we have a valid mesh with 3 patches and zero cells.
304  // Determine:
305  // - per face the top and bottom patch (topPatchID, bottomPatchID)
306  // - per edge, per face on edge the side patch (edgePatches)
307  labelListList edgePatches(extrudedPatch_.nEdges());
308  forAll(edgePatches, edgeI)
309  {
310  const labelList& eFaces = extrudedPatch_.edgeFaces()[edgeI];
311 
312  if (eFaces.size() != 2 || nonManifoldEdge.test(edgeI))
313  {
314  edgePatches[edgeI].setSize(eFaces.size(), sidePatchID);
315  }
316  }
317 
318  polyTopoChange meshMod(regionPatches.size());
319 
320  extruder.setRefinement
321  (
322  firstDisp, // first displacement
323  model_().expansionRatio(),
324  model_().nLayers(), // nLayers
325  labelList(extrudedPatch_.size(), topPatchID),
326  labelList(extrudedPatch_.size(), bottomPatchID),
327  edgePatches,
328  meshMod
329  );
330 
331  autoPtr<mapPolyMesh> map = meshMod.changeMesh
332  (
333  *this, // mesh to change
334  false // inflate
335  );
336 
337  // Update numbering on extruder.
338  extruder.updateMesh(map());
339 
340  this->setInstance(this->thisDb().time().constant());
341  this->write();
342  }
343 }
344 
345 
346 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
347 
348 } // End namespace Foam
349 
350 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
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::IOobject::NO_WRITE
Definition: IOobject.H:130
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:74
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Foam::HashTable< regIOobject * >::size
label size() const noexcept
The number of elements in table.
Definition: HashTableI.H:52
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::PrimitivePatch::points
const Field< PointType > & points() const
Return reference to global points.
Definition: PrimitivePatch.H:300
Foam::PrimitivePatch::edgeFaces
const labelListList & edgeFaces() const
Return edge-face addressing.
Definition: PrimitivePatch.C:318
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::fvMesh::write
virtual bool write(const bool valid=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:939
Foam::bitSet
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:64
Foam::primitiveFacePatch
PrimitivePatch< face, List, const pointField & > primitiveFacePatch
Foam::primitiveFacePatch.
Definition: primitiveFacePatch.H:47
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::fvMesh::thisDb
virtual const objectRegistry & thisDb() const
Return the object registry - resolve conflict polyMesh/lduMesh.
Definition: fvMesh.H:258
wedgePolyPatch.H
Foam::PrimitivePatch::nEdges
label nEdges() const
Return number of edges in patch.
Definition: PrimitivePatch.H:317
polyTopoChange.H
wallPolyPatch.H
Foam::polyMesh::facesInstance
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:821
Foam::bitSet::set
void set(const bitSet &bitset)
Set specified bits from another bitset.
Definition: bitSetI.H:563
Foam::fvMesh::removeFvBoundary
void removeFvBoundary()
Remove boundary patches. Warning: fvPatchFields hold ref to.
Definition: fvMesh.C:522
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:435
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:81
Foam::dictionary::set
entry * set(entry *entryPtr)
Assign a new entry, overwriting any existing entry.
Definition: dictionary.C:842
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
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::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:62
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
Foam::IOobject::READ_IF_PRESENT
Definition: IOobject.H:122
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:65
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:121
Foam::fvMesh::addFvPatches
void addFvPatches(PtrList< polyPatch > &plist, const bool validBoundary=true)
Add boundary patches. Constructor helper.
Definition: fvMesh.C:492
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 a pointer to a new patch created on freestore from.
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:84
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::PrimitivePatch::nInternalEdges
label nInternalEdges() const
Number of internal edges.
Definition: PrimitivePatch.C:258
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
Foam::PrimitivePatch::faceNormals
const Field< PointType > & faceNormals() const
Return face unit normals for patch.
Definition: PrimitivePatch.C:577
Foam::extrudePatchMesh::extrudePatchMesh
extrudePatchMesh(const fvMesh &, const fvPatch &, const dictionary &, const word)
Construct from mesh, patch and dictionary.
Definition: extrudePatchMesh.C:79
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: HashTable.H:102
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:246
Foam::PrimitivePatch::localFaces
const List< Face > & localFaces() const
Return patch faces addressing into local point list.
Definition: PrimitivePatch.C:398
Foam::polyMesh::globalData
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1241
Foam::point
vector point
Point is a vector.
Definition: point.H:43
Foam::PrimitivePatch::meshPoints
const labelList & meshPoints() const
Return labelList of mesh points in patch.
Definition: PrimitivePatch.C:418
Foam::List::setSize
void setSize(const label newSize)
Alias for resize(const label)
Definition: ListI.H:146
constant
constant condensation/saturation model.
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
extrudePatchMesh.H
Foam::polyMesh::setInstance
void setInstance(const fileName &instance, const IOobject::writeOption wOpt=IOobject::AUTO_WRITE)
Set the instance for mesh files.
Definition: polyMeshIO.C:35
Foam::fvMesh::clearOut
void clearOut()
Clear all geometry and addressing.
Definition: fvMesh.C:234