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-2022 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 "blockMesh.H"
30#include "transform.H"
31#include "Time.H"
32
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34
35namespace Foam
36{
38}
39
41
42
44Foam::blockMesh::strategyNames_
45({
46 { mergeStrategy::MERGE_TOPOLOGY, "topology" },
47 { mergeStrategy::MERGE_POINTS, "points" },
48});
49
50
51// * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
52
53namespace Foam
54{
55
56// Define/refine the verbosity level
57// Command-line options have precedence over dictionary setting
58static 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
80int 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
146bool Foam::blockMesh::readPointTransforms(const dictionary& dict)
147{
148 transformType_ = transformTypes::NO_TRANSFORM;
149
150 const dictionary* dictptr;
151
152 // Optional cartesian coordinate system transform, since JUL-2021
153 if ((dictptr = dict.findDict("transform", keyType::LITERAL)) != nullptr)
154 {
155 transform_ = coordSystem::cartesian(*dictptr);
156
157 // Non-zero origin?
158 if (magSqr(transform_.origin()) > ROOTVSMALL)
159 {
160 transformType_ |= transformTypes::TRANSLATION;
161 }
162
163 // Non-identity rotation?
164 if (!transform_.R().is_identity(ROOTVSMALL))
165 {
166 transformType_ |= transformTypes::ROTATION;
167 }
168 }
169 else
170 {
171 transform_.clear();
172 }
173
174
175 // Optional 'prescale' factor.
176 {
177 prescaling_ = vector::uniform(1);
178
179 const int scaleType = readScaling
180 (
181 dict.findEntry("prescale", keyType::LITERAL),
182 prescaling_
183 );
184
185 if (scaleType == 1)
186 {
187 transformType_ |= transformTypes::PRESCALING;
188 }
189 else if (scaleType == 3)
190 {
191 transformType_ |= transformTypes::PRESCALING3;
192 }
193 }
194
195
196 // Optional 'scale' factor. Was 'convertToMeters' until OCT-2008
197 {
198 scaling_ = vector::uniform(1);
199
200 const int scaleType = readScaling
201 (
202 // Mark as changed from 2010 onwards
204 (
205 "scale",
206 {{"convertToMeters", 1012}},
208 ),
209 scaling_
210 );
211
212 if (scaleType == 1)
213 {
214 transformType_ |= transformTypes::SCALING;
215 }
216 else if (scaleType == 3)
217 {
218 transformType_ |= transformTypes::SCALING3;
219 }
220 }
221
222 return bool(transformType_);
223}
224
225
226// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
227
229(
230 const IOdictionary& dict,
231 const word& regionName,
232 mergeStrategy strategy,
233 int verbosity
234)
235:
236 meshDict_(dict),
237 verbose_(getVerbosity(dict, verbosity)),
238 checkFaceCorrespondence_
239 (
240 meshDict_.getOrDefault("checkFaceCorrespondence", true)
241 ),
242 mergeStrategy_(strategy),
243 transformType_(transformTypes::NO_TRANSFORM),
244 geometry_
245 (
247 (
248 "geometry", // dummy name
249 meshDict_.time().constant(), // instance
250 "geometry", // local
251 meshDict_.time(), // registry
252 IOobject::MUST_READ,
253 IOobject::NO_WRITE
254 ),
255 meshDict_.found("geometry")
256 ? meshDict_.subDict("geometry")
257 : dictionary(),
258 true
259 ),
260 blockVertices_
261 (
262 meshDict_.lookup("vertices"),
263 blockVertex::iNew(meshDict_, geometry_)
264 ),
265 vertices_(Foam::vertices(blockVertices_)),
266 prescaling_(vector::uniform(1)),
267 scaling_(vector::uniform(1)),
268 transform_(),
269 topologyPtr_(createTopology(meshDict_, regionName))
270{
271 if (mergeStrategy_ == mergeStrategy::DEFAULT_MERGE)
272 {
273 strategyNames_.readIfPresent("mergeType", meshDict_, mergeStrategy_);
274
275 // Warn about fairly obscure old "fastMerge" option?
276
277 if
278 (
279 mergeStrategy_ == mergeStrategy::DEFAULT_MERGE
280 && checkDegenerate()
281 )
282 {
283 Info<< nl
284 << "Detected collapsed blocks "
285 << "- using merge points instead of merge topology" << nl
286 << endl;
287
288 mergeStrategy_ = mergeStrategy::MERGE_POINTS;
289 }
290 }
291
292 if (mergeStrategy_ == mergeStrategy::MERGE_POINTS)
293 {
294 // MERGE_POINTS
295 calcGeometricalMerge();
296 }
297 else
298 {
299 // MERGE_TOPOLOGY
300 calcTopologicalMerge();
301 }
302}
303
304
305// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
306
308{
309 return bool(topologyPtr_);
310}
311
312
314{
315 return verbose_;
316}
317
318
319int Foam::blockMesh::verbose(const int level) noexcept
320{
321 int old(verbose_);
322 verbose_ = level;
323 return old;
324}
325
326
328{
329 return vertices_;
330}
331
332
334Foam::blockMesh::vertices(bool applyTransform) const
335{
336 if (applyTransform && hasPointTransforms())
337 {
338 auto tpts = tmp<pointField>::New(vertices_);
339
340 inplacePointTransforms(tpts.ref());
341
342 return tpts;
343 }
344
345 return vertices_;
346}
347
348
350{
351 const polyPatchList& topoPatches = topology().boundaryMesh();
352
353 PtrList<dictionary> patchDicts(topoPatches.size());
354
355 forAll(topoPatches, patchi)
356 {
358 topoPatches[patchi].write(os);
359 IStringStream is(os.str());
360 patchDicts.set(patchi, new dictionary(is));
361 }
362 return patchDicts;
363}
364
365
367{
368 if (points_.empty())
369 {
370 createPoints();
371 }
372
373 return points_;
374}
375
376
378{
379 if (cells_.empty())
380 {
381 createCells();
382 }
383
384 return cells_;
385}
386
387
389{
390 if (patches_.empty())
391 {
392 createPatches();
393 }
394
395 return patches_;
396}
397
398
400{
401 return topology().boundaryMesh().names();
402}
403
404
405//Foam::wordList Foam::blockMesh::patchTypes() const
406//{
407// return topology().boundaryMesh().types();
408//}
409//
410//
411//Foam::wordList Foam::blockMesh::patchPhysicalTypes() const
412//{
413// return topology().boundaryMesh().physicalTypes();
414//}
415
416
418{
419 const blockList& blocks = *this;
420
421 label count = 0;
422
423 for (const block& blk : blocks)
424 {
425 if (blk.zoneName().size())
426 {
427 ++count;
428 }
429 }
430
431 return count;
432}
433
434
436{
437 return bool(transformType_);
438}
439
440
442{
443 if (!transformType_)
444 {
445 return false;
446 }
447
448 if (transformType_ & transformTypes::PRESCALING)
449 {
450 for (point& p : pts)
451 {
452 p *= prescaling_.x();
453 }
454 }
455 else if (transformType_ & transformTypes::PRESCALING3)
456 {
457 for (point& p : pts)
458 {
459 p = cmptMultiply(p, prescaling_);
460 }
461 }
462
463 if (transformType_ & transformTypes::ROTATION)
464 {
465 const tensor rot(transform_.R());
466
467 if (transformType_ & transformTypes::TRANSLATION)
468 {
469 const point origin(transform_.origin());
470
471 for (point& p : pts)
472 {
473 p = Foam::transform(rot, p) + origin;
474 }
475 }
476 else
477 {
478 for (point& p : pts)
479 {
480 p = Foam::transform(rot, p);
481 }
482 }
483 }
484 else if (transformType_ & transformTypes::TRANSLATION)
485 {
486 const point origin(transform_.origin());
487
488 for (point& p : pts)
489 {
490 p += origin;
491 }
492 }
493
494 if (transformType_ & transformTypes::SCALING)
495 {
496 for (point& p : pts)
497 {
498 p *= scaling_.x();
499 }
500 }
501 else if (transformType_ & transformTypes::SCALING3)
502 {
503 for (point& p : pts)
504 {
505 p = cmptMultiply(p, scaling_);
506 }
507 }
508
509 return true;
510}
511
512
515{
516 if (hasPointTransforms())
517 {
518 auto tpts = tmp<pointField>::New(localPoints);
519
520 inplacePointTransforms(tpts.ref());
521
522 return tpts;
523 }
524 else
525 {
526 return localPoints;
527 }
528}
529
530
531// ************************************************************************* //
bool found
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: Enum.H:61
bool readIfPresent(const word &key, const dictionary &dict, EnumType &val) const
Find an entry if present, and assign to T val.
Definition: EnumI.H:132
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
Input from string buffer, using a ISstream. Always UNCOMPRESSED.
Definition: StringStream.H:112
An input stream of tokens.
Definition: ITstream.H:56
const token & peek() const
Failsafe peek at what the next read would return,.
Definition: ITstream.C:352
Output to string buffer, using a OSstream. Always UNCOMPRESSED.
Definition: StringStream.H:231
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:106
A multi-block mesh generator.
Definition: blockMesh.H:168
const faceListList & patches() const
Return the patch face lists.
Definition: blockMesh.C:388
bool hasPointTransforms() const noexcept
True if scaling and/or transformations are needed.
Definition: blockMesh.C:435
const pointField & vertices() const noexcept
Definition: blockMesh.C:327
static bool verboseOutput
The default verbosity (true)
Definition: blockMesh.H:339
wordList patchNames() const
Return the patch names.
Definition: blockMesh.C:399
int verbose() const noexcept
Output verbosity level.
Definition: blockMesh.C:313
tmp< pointField > globalPosition(const pointField &localPoints) const
Apply coordinate transforms and scaling.
Definition: blockMesh.C:514
mergeStrategy
The block merging strategy.
Definition: blockMesh.H:181
@ DEFAULT_MERGE
Default (TOPOLOGY), not selectable.
Definition: blockMesh.H:182
@ MERGE_POINTS
"points" merge by point geometry
Definition: blockMesh.H:184
label numZonedBlocks() const
Number of blocks with specified zones.
Definition: blockMesh.C:417
bool valid() const noexcept
True if the blockMesh topology exists.
Definition: blockMesh.C:307
bool inplacePointTransforms(pointField &pts) const
Apply coordinate transforms and scaling.
Definition: blockMesh.C:441
const cellShapeList & cells() const
Return cell shapes list.
Definition: blockMesh.C:377
const pointField & points() const
Definition: blockMesh.C:366
PtrList< dictionary > patchDicts() const
Patch information from the topology mesh.
Definition: blockMesh.C:349
Define a block vertex.
Definition: blockVertex.H:52
Creates a single block of cells from point coordinates, numbers of cells in each direction and an exp...
Definition: block.H:61
virtual const point & origin() const
Return origin.
virtual const tensor & R() const
Return const reference to the rotation tensor.
virtual void clear()
Reset origin and rotation to an identity coordinateSystem.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
const entry * findCompat(const word &keyword, std::initializer_list< std::pair< const char *, int > > compat, enum keyType::option matchOpt) const
dictionary * findDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX)
Find and return a sub-dictionary pointer if present.
Definition: dictionaryI.H:127
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
entry * findEntry(const word &keyword, enum keyType::option matchOpt=keyType::REGEX)
Find for an entry (non-const access) with the given keyword.
Definition: dictionaryI.H:97
Particle-size distribution model wherein random samples are drawn from the doubly-truncated uniform p...
Definition: uniform.H:164
A keyword and a list of tokens is an 'entry'.
Definition: entry.H:70
virtual ITstream & stream() const =0
Return token stream, if entry is a primitive entry.
void checkITstream(const ITstream &is) const
Definition: entry.C:110
@ LITERAL
String literal.
Definition: keyType.H:81
static constexpr direction nComponents
Number of components in bool is 1.
Definition: bool.H:98
constant condensation/saturation model.
Lookup type of boundary radiation properties.
Definition: lookup.H:66
A class for managing temporary objects.
Definition: tmp.H:65
bool isNumber() const noexcept
Token is LABEL, FLOAT or DOUBLE.
Definition: tokenI.H:587
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
bool
Definition: EEqn.H:20
Foam::word regionName(Foam::polyMesh::defaultRegion)
#define defineDebugSwitch(Type, Value)
Define the debug information.
OBJstream os(runTime.globalPath()/outputName)
Namespace for OpenFOAM.
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:536
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
bool equal(const T &s1, const T &s2)
Compare two values for equality.
Definition: doubleFloat.H:46
int readScaling(const entry *eptr, vector &scale)
Definition: blockMesh.C:80
pointField vertices(const blockVertexList &bvl)
uint8_t direction
Definition: direction.H:56
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
const direction noexcept
Definition: Scalar.H:223
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
static int getVerbosity(const dictionary &dict, int verbosity)
Definition: blockMesh.C:58
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
PtrList< dictionary > patchDicts
Definition: readKivaGrid.H:532
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
3D tensor transformation operations.