blockMesh.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-2017 OpenFOAM Foundation
9  Copyright (C) 2018-2021 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 "blockMesh.H"
30 #include "transform.H"
31 #include "Time.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37  defineDebugSwitch(blockMesh, 0);
38 }
39 
41 
42 
44 Foam::blockMesh::strategyNames_
45 ({
46  { mergeStrategy::MERGE_TOPOLOGY, "topology" },
47  { mergeStrategy::MERGE_POINTS, "points" },
48 });
49 
50 
51 // * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
52 
53 namespace Foam
54 {
55 
56 // Define/refine the verbosity level
57 // Command-line options have precedence over dictionary setting
58 static int getVerbosity(const dictionary& dict, int verbosity)
59 {
60  if (verbosity < 0)
61  {
62  // Forced as 'off'
63  verbosity = 0;
64  }
65  else if (!verbosity)
66  {
67  // Not specified: use dictionary value or static default
68  verbosity = dict.getOrDefault("verbose", blockMesh::verboseOutput);
69  }
70 
71  return verbosity;
72 }
73 
74 
75 // Process dictionary entry with single scalar or vector quantity
76 // - return 0 if scaling is not needed. Eg, Not found or is unity
77 // - return 1 for uniform scaling
78 // - return 3 for non-uniform scaling
79 
80 int readScaling(const entry *eptr, vector& scale)
81 {
82  int nCmpt = 0;
83 
84  if (!eptr)
85  {
86  return nCmpt;
87  }
88 
89  ITstream& is = eptr->stream();
90 
91  if (is.peek().isNumber())
92  {
93  // scalar value
94  scalar val;
95  is >> val;
96 
97  if ((val > 0) && !equal(val, 1))
98  {
99  // Uniform scaling
100  nCmpt = 1;
101  scale = vector::uniform(val);
102  }
103  }
104  else
105  {
106  // vector value
107  is >> scale;
108 
109  bool nonUnity = false;
110  for (direction cmpt=0; cmpt < vector::nComponents; ++cmpt)
111  {
112  if (scale[cmpt] <= 0)
113  {
114  scale[cmpt] = 1;
115  }
116  else if (!equal(scale[cmpt], 1))
117  {
118  nonUnity = true;
119  }
120  }
121 
122  if (nonUnity)
123  {
124  if (equal(scale.x(), scale.y()) && equal(scale.x(), scale.z()))
125  {
126  // Uniform scaling
127  nCmpt = 1;
128  }
129  else
130  {
131  nCmpt = 3;
132  }
133  }
134  }
135 
136  eptr->checkITstream(is);
137 
138  return nCmpt;
139 }
140 
141 } // End namespace Foam
142 
143 
144 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
145 
146 bool Foam::blockMesh::readPointTransforms(const dictionary& dict)
147 {
148  transformType_ = transformTypes::NO_TRANSFORM;
149 
150  // Optional cartesian coordinate system transform, since JUL-2021
151  if (const dictionary* dictptr = dict.findDict("transform"))
152  {
153  transform_ = coordSystem::cartesian(*dictptr);
154 
155  constexpr scalar tol = VSMALL;
156 
157  // Check for non-zero origin
158  const vector& o = transform_.origin();
159  if
160  (
161  (mag(o.x()) > tol)
162  || (mag(o.y()) > tol)
163  || (mag(o.z()) > tol)
164  )
165  {
166  transformType_ |= transformTypes::TRANSLATION;
167  }
168 
169  // Check for non-identity rotation
170  const tensor& r = transform_.R();
171  if
172  (
173  (mag(r.xx() - 1) > tol)
174  || (mag(r.xy()) > tol)
175  || (mag(r.xz()) > tol)
176 
177  || (mag(r.yx()) > tol)
178  || (mag(r.yy() - 1) > tol)
179  || (mag(r.yz()) > tol)
180 
181  || (mag(r.zx()) > tol)
182  || (mag(r.zy()) > tol)
183  || (mag(r.zz() - 1) > tol)
184  )
185  {
186  transformType_ |= transformTypes::ROTATION;
187  }
188  }
189  else
190  {
191  transform_.clear();
192  }
193 
194 
195  // Optional 'prescale' factor.
196  {
197  prescaling_ = vector::uniform(1);
198 
199  const int scaleType = readScaling
200  (
201  dict.findEntry("prescale", keyType::LITERAL),
202  prescaling_
203  );
204 
205  if (scaleType == 1)
206  {
207  transformType_ |= transformTypes::PRESCALING;
208  }
209  else if (scaleType == 3)
210  {
211  transformType_ |= transformTypes::PRESCALING3;
212  }
213  }
214 
215 
216  // Optional 'scale' factor. Was 'convertToMeters' until OCT-2008
217  {
218  scaling_ = vector::uniform(1);
219 
220  const int scaleType = readScaling
221  (
222  // Mark as changed from 2010 onwards
223  dict.findCompat
224  (
225  "scale",
226  {{"convertToMeters", 1012}},
228  ),
229  scaling_
230  );
231 
232  if (scaleType == 1)
233  {
234  transformType_ |= transformTypes::SCALING;
235  }
236  else if (scaleType == 3)
237  {
238  transformType_ |= transformTypes::SCALING3;
239  }
240  }
241 
242  return bool(transformType_);
243 }
244 
245 
246 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
247 
248 Foam::blockMesh::blockMesh
249 (
250  const IOdictionary& dict,
251  const word& regionName,
252  mergeStrategy strategy,
253  int verbosity
254 )
255 :
256  meshDict_(dict),
257  verbose_(getVerbosity(dict, verbosity)),
258  checkFaceCorrespondence_
259  (
260  meshDict_.getOrDefault("checkFaceCorrespondence", true)
261  ),
262  mergeStrategy_(strategy),
263  transformType_(transformTypes::NO_TRANSFORM),
264  geometry_
265  (
266  IOobject
267  (
268  "geometry", // dummy name
269  meshDict_.time().constant(), // instance
270  "geometry", // local
271  meshDict_.time(), // registry
274  ),
275  meshDict_.found("geometry")
276  ? meshDict_.subDict("geometry")
277  : dictionary(),
278  true
279  ),
280  blockVertices_
281  (
282  meshDict_.lookup("vertices"),
283  blockVertex::iNew(meshDict_, geometry_)
284  ),
285  vertices_(Foam::vertices(blockVertices_)),
286  prescaling_(vector::uniform(1)),
287  scaling_(vector::uniform(1)),
288  transform_(),
289  topologyPtr_(createTopology(meshDict_, regionName))
290 {
291  if (mergeStrategy_ == mergeStrategy::DEFAULT_MERGE)
292  {
293  strategyNames_.readIfPresent("mergeType", meshDict_, mergeStrategy_);
294 
295  // Warn about fairly obscure old "fastMerge" option?
296 
297  if
298  (
299  mergeStrategy_ == mergeStrategy::DEFAULT_MERGE
300  && checkDegenerate()
301  )
302  {
303  Info<< nl
304  << "Detected collapsed blocks "
305  << "- using merge points instead of merge topology" << nl
306  << endl;
307 
308  mergeStrategy_ = mergeStrategy::MERGE_POINTS;
309  }
310  }
311 
312  if (mergeStrategy_ == mergeStrategy::MERGE_POINTS)
313  {
314  // MERGE_POINTS
315  calcGeometricalMerge();
316  }
317  else
318  {
319  // MERGE_TOPOLOGY
320  calcTopologicalMerge();
321  }
322 }
323 
324 
325 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
326 
327 bool Foam::blockMesh::valid() const noexcept
328 {
329  return bool(topologyPtr_);
330 }
331 
332 
333 int Foam::blockMesh::verbose() const noexcept
334 {
335  return verbose_;
336 }
337 
338 
339 int Foam::blockMesh::verbose(const int level) noexcept
340 {
341  int old(verbose_);
342  verbose_ = level;
343  return old;
344 }
345 
346 
348 {
349  return vertices_;
350 }
351 
352 
354 Foam::blockMesh::vertices(bool applyTransform) const
355 {
356  if (applyTransform && hasPointTransforms())
357  {
358  auto tpts = tmp<pointField>::New(vertices_);
359 
360  inplacePointTransforms(tpts.ref());
361 
362  return tpts;
363  }
364 
365  return vertices_;
366 }
367 
368 
370 {
371  const polyPatchList& topoPatches = topology().boundaryMesh();
372 
373  PtrList<dictionary> patchDicts(topoPatches.size());
374 
375  forAll(topoPatches, patchi)
376  {
378  topoPatches[patchi].write(os);
379  IStringStream is(os.str());
380  patchDicts.set(patchi, new dictionary(is));
381  }
382  return patchDicts;
383 }
384 
385 
387 {
388  if (points_.empty())
389  {
390  createPoints();
391  }
392 
393  return points_;
394 }
395 
396 
398 {
399  if (cells_.empty())
400  {
401  createCells();
402  }
403 
404  return cells_;
405 }
406 
407 
409 {
410  if (patches_.empty())
411  {
412  createPatches();
413  }
414 
415  return patches_;
416 }
417 
418 
420 {
421  return topology().boundaryMesh().names();
422 }
423 
424 
425 //Foam::wordList Foam::blockMesh::patchTypes() const
426 //{
427 // return topology().boundaryMesh().types();
428 //}
429 //
430 //
431 //Foam::wordList Foam::blockMesh::patchPhysicalTypes() const
432 //{
433 // return topology().boundaryMesh().physicalTypes();
434 //}
435 
436 
438 {
439  const blockList& blocks = *this;
440 
441  label count = 0;
442 
443  for (const block& blk : blocks)
444  {
445  if (blk.zoneName().size())
446  {
447  ++count;
448  }
449  }
450 
451  return count;
452 }
453 
454 
456 {
457  return bool(transformType_);
458 }
459 
460 
462 {
463  if (!transformType_)
464  {
465  return false;
466  }
467 
468  if (transformType_ & transformTypes::PRESCALING)
469  {
470  for (point& p : pts)
471  {
472  p *= prescaling_.x();
473  }
474  }
475  else if (transformType_ & transformTypes::PRESCALING3)
476  {
477  for (point& p : pts)
478  {
479  p = cmptMultiply(p, prescaling_);
480  }
481  }
482 
483  if (transformType_ & transformTypes::ROTATION)
484  {
485  const tensor rot(transform_.R());
486 
487  if (transformType_ & transformTypes::TRANSLATION)
488  {
489  const point origin(transform_.origin());
490 
491  for (point& p : pts)
492  {
493  p = Foam::transform(rot, p) + origin;
494  }
495  }
496  else
497  {
498  for (point& p : pts)
499  {
500  p = Foam::transform(rot, p);
501  }
502  }
503  }
504  else if (transformType_ & transformTypes::TRANSLATION)
505  {
506  const point origin(transform_.origin());
507 
508  for (point& p : pts)
509  {
510  p += origin;
511  }
512  }
513 
514  if (transformType_ & transformTypes::SCALING)
515  {
516  for (point& p : pts)
517  {
518  p *= scaling_.x();
519  }
520  }
521  else if (transformType_ & transformTypes::SCALING3)
522  {
523  for (point& p : pts)
524  {
525  p = cmptMultiply(p, scaling_);
526  }
527  }
528 
529  return true;
530 }
531 
532 
535 {
536  if (hasPointTransforms())
537  {
538  auto tpts = tmp<pointField>::New(localPoints);
539 
540  inplacePointTransforms(tpts.ref());
541 
542  return tpts;
543  }
544  else
545  {
546  return localPoints;
547  }
548 }
549 
550 
551 // ************************************************************************* //
Foam::entry
A keyword and a list of tokens is an 'entry'.
Definition: entry.H:67
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:54
Foam::Tensor< scalar >
Foam::blockMesh::patchDicts
PtrList< dictionary > patchDicts() const
Patch information from the topology mesh.
Definition: blockMesh.C:369
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::block
Creates a single block of cells from point coordinates, numbers of cells in each direction and an exp...
Definition: block.H:58
Foam::blockMesh::valid
bool valid() const noexcept
True if the blockMesh topology exists.
Definition: blockMesh.C:327
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::Enum
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: IOstreamOption.H:57
Foam::blockMesh::patches
const faceListList & patches() const
Return the patch face lists.
Definition: blockMesh.C:408
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::VectorSpace< Vector< scalar >, scalar, 3 >::uniform
static Vector< scalar > uniform(const scalar &s)
Return a VectorSpace with all elements = s.
Definition: VectorSpaceI.H:164
Foam::entry::stream
virtual ITstream & stream() const =0
Return token stream, if entry is a primitive entry.
Foam::cmptMultiply
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::coordinateSystem::clear
virtual void clear()
Reset origin and rotation to an identity coordinateSystem.
Definition: coordinateSystem.C:286
Foam::blockMesh::verbose
int verbose() const noexcept
Output verbosity level.
Definition: blockMesh.C:333
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::blockMesh::mergeStrategy
mergeStrategy
The block merging strategy.
Definition: blockMesh.H:180
Foam::blockMesh::vertices
const pointField & vertices() const noexcept
Definition: blockMesh.C:347
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
blockMesh.H
Foam::transform
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:521
Foam::readScaling
int readScaling(const entry *eptr, vector &scale)
Definition: blockMesh.C:80
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
patchDicts
PtrList< dictionary > patchDicts
Definition: readKivaGrid.H:532
Foam::token::isNumber
bool isNumber() const noexcept
Token is LABEL, FLOAT or DOUBLE.
Definition: tokenI.H:587
Foam::coordinateSystem::origin
virtual const point & origin() const
Return origin.
Definition: coordinateSystem.H:469
regionName
Foam::word regionName
Definition: createNamedDynamicFvMesh.H:1
Foam::Field< vector >
Foam::ITstream
An input stream of tokens.
Definition: ITstream.H:52
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::defineDebugSwitch
defineDebugSwitch(blockMesh, 0)
Foam::PtrList< Foam::dictionary >
Foam::blockMesh::inplacePointTransforms
bool inplacePointTransforms(pointField &pts) const
Apply coordinate transforms and scaling.
Definition: blockMesh.C:461
Foam::blockMesh::globalPosition
tmp< pointField > globalPosition(const pointField &localPoints) const
Apply coordinate transforms and scaling.
Definition: blockMesh.C:534
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::vertices
pointField vertices(const blockVertexList &bvl)
Definition: blockVertexList.H:49
os
OBJstream os(runTime.globalPath()/outputName)
Foam::IStringStream
Input from string buffer, using a ISstream. Always UNCOMPRESSED.
Definition: StringStream.H:108
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
Foam::blockMesh::hasPointTransforms
bool hasPointTransforms() const noexcept
True if scaling and/or transformations are needed.
Definition: blockMesh.C:455
Time.H
Foam::blockMesh::numZonedBlocks
label numZonedBlocks() const
Number of blocks with specified zones.
Definition: blockMesh.C:437
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::blockMesh::points
const pointField & points() const
Definition: blockMesh.C:386
Foam::BitOps::count
unsigned int count(const UList< bool > &bools, const bool val=true)
Count number of 'true' entries.
Definition: BitOps.H:77
Foam::Vector< scalar >
Foam::List< cellShape >
Foam::blockMesh::patchNames
wordList patchNames() const
Return the patch names.
Definition: blockMesh.C:419
Foam::OStringStream
Output to string buffer, using a OSstream. Always UNCOMPRESSED.
Definition: StringStream.H:227
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::tmp::New
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
bool
bool
Definition: EEqn.H:20
Foam::direction
uint8_t direction
Definition: direction.H:52
Foam::blockVertex::iNew
Class used for the read-construction of.
Definition: blockVertex.H:95
Foam::keyType::LITERAL
String literal.
Definition: keyType.H:81
transform.H
3D tensor transformation operations.
Foam::entry::checkITstream
void checkITstream(const ITstream &is) const
Definition: entry.C:110
Foam::dictionary::getOrDefault
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:148
Foam::blockMesh::cells
const cellShapeList & cells() const
Return cell shapes list.
Definition: blockMesh.C:397
Foam::getVerbosity
static int getVerbosity(const dictionary &dict, int verbosity)
Definition: blockMesh.C:58
Foam::ITstream::peek
const token & peek() const
Failsafe peek at what the next read would return,.
Definition: ITstream.C:352
Foam::equal
bool equal(const T &s1, const T &s2)
Compare two values for equality.
Definition: doubleFloat.H:46
Foam::tensor
Tensor< scalar > tensor
Tensor of scalars, i.e. Tensor<scalar>.
Definition: symmTensor.H:61
Foam::blockMesh::verboseOutput
static bool verboseOutput
The default verbosity (true)
Definition: blockMesh.H:339
Foam::VectorSpace< Vector< scalar >, scalar, 3 >::nComponents
static constexpr direction nComponents
Number of components in this vector space.
Definition: VectorSpace.H:101
Foam::IOobject::MUST_READ
Definition: IOobject.H:185
Foam::coordinateSystem::R
virtual const tensor & R() const
Return const reference to the rotation tensor.
Definition: coordinateSystem.H:475