Loading...
Searching...
No Matches
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}
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 // Optional (cartesian) coordinate system transform
151 auto csysPtr = coordinateSystem::NewIfPresent(dict, "transform");
152
153 if (csysPtr)
154 {
155 transform_ = csysPtr();
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
203 dict.findCompat
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 }
222 return bool(transformType_);
223}
224
225
226// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
227
228Foam::blockMesh::blockMesh
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 * * * * * * * * * * * * * //
308{
309 return bool(topologyPtr_);
310}
311
314{
315 return verbose_;
316}
317
318
319int Foam::blockMesh::verbose(const int level) noexcept
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 OCharStream os;
356 ISpanStream is;
357
358 forAll(topoPatches, patchi)
359 {
360 os.rewind();
361 topoPatches[patchi].write(os);
362
363 is.reset(os.view());
364 patchDicts.emplace_set(patchi, is);
365 }
366 return patchDicts;
367}
368
369
371{
372 if (points_.empty())
373 {
374 createPoints();
375 }
376
377 return points_;
378}
379
380
382{
383 if (cells_.empty())
384 {
385 createCells();
386 }
387
388 return cells_;
389}
390
391
393{
394 if (patches_.empty())
395 {
396 createPatches();
397 }
398
399 return patches_;
400}
401
402
404{
405 return topology().boundaryMesh().names();
406}
407
408
409//Foam::wordList Foam::blockMesh::patchTypes() const
410//{
411// return topology().boundaryMesh().types();
412//}
413//
415//Foam::wordList Foam::blockMesh::patchPhysicalTypes() const
416//{
417// return topology().boundaryMesh().physicalTypes();
418//}
419
420
421Foam::label Foam::blockMesh::numZonedBlocks() const
422{
423 const blockList& blocks = *this;
424
425 label count = 0;
426
427 for (const block& blk : blocks)
428 {
429 if (blk.zoneName().size())
430 {
431 ++count;
433 }
434
435 return count;
436}
437
440{
441 return bool(transformType_);
442}
443
444
446{
447 if (!transformType_)
448 {
449 return false;
450 }
451
452 if (transformType_ & transformTypes::PRESCALING)
453 {
454 for (point& p : pts)
455 {
456 p *= prescaling_.x();
457 }
458 }
459 else if (transformType_ & transformTypes::PRESCALING3)
460 {
461 for (point& p : pts)
462 {
463 p = cmptMultiply(p, prescaling_);
464 }
465 }
466
467 if (transformType_ & transformTypes::ROTATION)
468 {
469 const tensor rot(transform_.R());
470
471 if (transformType_ & transformTypes::TRANSLATION)
472 {
473 const point origin(transform_.origin());
474
475 for (point& p : pts)
476 {
477 p = Foam::transform(rot, p) + origin;
478 }
479 }
480 else
481 {
482 for (point& p : pts)
483 {
484 p = Foam::transform(rot, p);
485 }
486 }
487 }
488 else if (transformType_ & transformTypes::TRANSLATION)
489 {
490 const point origin(transform_.origin());
491
492 for (point& p : pts)
493 {
494 p += origin;
495 }
496 }
497
498 if (transformType_ & transformTypes::SCALING)
499 {
500 for (point& p : pts)
501 {
502 p *= scaling_.x();
503 }
504 }
505 else if (transformType_ & transformTypes::SCALING3)
506 {
507 for (point& p : pts)
508 {
509 p = cmptMultiply(p, scaling_);
510 }
512
513 return true;
514}
515
516
518Foam::blockMesh::globalPosition(const pointField& localPoints) const
519{
520 if (hasPointTransforms())
521 {
522 auto tpts = tmp<pointField>::New(localPoints);
523
524 inplacePointTransforms(tpts.ref());
525
526 return tpts;
527 }
528 else
529 {
530 return localPoints;
531 }
532}
533
534
535// ************************************************************************* //
bool found
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition Enum.H:57
bool readIfPresent(const word &key, const dictionary &dict, EnumType &val, const bool warnOnly=false) const
Find an entry if present, and assign to T val.
Definition EnumI.H:111
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
Similar to IStringStream but using an externally managed buffer for its input. This allows the input ...
void reset(const char *buffer, size_t nbytes)
Reset input area, position to buffer start and clear errors.
An input stream of tokens.
Definition ITstream.H:56
const token & peek() const noexcept
Failsafe peek at what the next read would return, including handling of any putback.
Definition ITstream.C:327
An OSstream with internal List storage.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
constexpr PtrList() noexcept
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
label count() const noexcept
The number of non-nullptr entries in the list.
Definition UPtrList.H:890
static constexpr direction nComponents
Number of components in this vector space.
static Form uniform(const Cmpt &s)
Return a VectorSpace with all elements = s.
A multi-block mesh generator.
Definition blockMesh.H:164
const faceListList & patches() const
Return the patch face lists.
Definition blockMesh.C:385
bool hasPointTransforms() const noexcept
True if scaling and/or transformations are needed.
Definition blockMesh.C:432
const pointField & vertices() const noexcept
Reference to point field defining the blocks, these points are unscaled and non-transformed.
Definition blockMesh.C:320
static bool verboseOutput
The default verbosity (true).
Definition blockMesh.H:396
wordList patchNames() const
Return the patch names.
Definition blockMesh.C:396
int verbose() const noexcept
Output verbosity level.
Definition blockMesh.C:306
bool good() const noexcept
True if the blockMesh topology exists.
Definition blockMesh.C:300
tmp< pointField > globalPosition(const pointField &localPoints) const
Apply coordinate transforms and scaling.
Definition blockMesh.C:511
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:414
bool inplacePointTransforms(pointField &pts) const
Apply coordinate transforms and scaling.
Definition blockMesh.C:438
const cellShapeList & cells() const
Return cell shapes list.
Definition blockMesh.C:374
const pointField & points() const
The points for the entire mesh. These points are scaled and transformed.
Definition blockMesh.C:363
PtrList< dictionary > patchDicts() const
Patch information from the topology mesh.
Definition blockMesh.C:342
const polyMesh & topology() const
The blockMesh topology as a polyMesh unscaled and non-transformed.
PtrList< block > blockList
The list of blocks is stored as a PtrList.
Definition blockMesh.H:172
Define a block vertex.
Definition blockVertex.H:48
Creates a single block of cells from point coordinates, numbers of cells in each direction and an exp...
Definition block.H:57
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A keyword and a list of tokens is an 'entry'.
Definition entry.H:66
virtual ITstream & stream() const =0
Return token stream, if entry is a primitive entry.
void checkITstream(const ITstream &is) const
Check after reading if the input token stream has unconsumed tokens remaining or if there were no tok...
Definition entry.C:103
@ LITERAL
String literal.
Definition keyType.H:82
Lookup type of boundary radiation properties.
Definition lookup.H:60
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
bool isNumber() const noexcept
Token is (signed/unsigned) integer type, FLOAT or DOUBLE.
Definition tokenI.H:992
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
Foam::word regionName(args.getOrDefault< word >("region", Foam::polyMesh::defaultRegion))
#define defineDebugSwitch(Type, Value)
Define the debug information.
OBJstream os(runTime.globalPath()/outputName)
Different types of constants.
Namespace for OpenFOAM.
PtrList< polyPatch > polyPatchList
Store lists of polyPatch as a PtrList.
Definition polyPatch.H:61
List< word > wordList
List of word.
Definition fileName.H:60
bool equal(const T &a, const T &b)
Compare two values for equality.
Definition label.H:180
refinementData transform(const tensor &, const refinementData val)
No-op rotational transform for base types.
messageStream Info
Information stream (stdout output on master, null elsewhere).
Tensor< scalar > tensor
Definition symmTensor.H:57
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
int readScaling(const entry *eptr, vector &scale)
Definition blockMesh.C:73
List< faceList > faceListList
List of faceList.
Definition faceListFwd.H:44
pointField vertices(const blockVertexList &bvl)
uint8_t direction
Definition direction.H:49
vector point
Point is a vector.
Definition point.H:37
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
const direction noexcept
Definition scalarImpl.H:265
vectorField pointField
pointField is a vectorField.
Vector< scalar > vector
Definition vector.H:57
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
List< cellShape > cellShapeList
List of cellShape.
static int getVerbosity(const dictionary &dict, int verbosity)
Definition blockMesh.C:51
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
PtrList< dictionary > patchDicts
dictionary dict
const pointField & pts
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
3D tensor transformation operations.