Loading...
Searching...
No Matches
displacementSmartPointSmoothingMotionSolver.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) 2024-2025 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
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
30#include "syncTools.H"
31#include "pointConstraints.H"
32#include "motionSmootherAlgo.H"
33
34//#include "fvMesh.H"
35//#include "fvGeometryScheme.H"
36#include "OBJstream.H"
39// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40
41namespace Foam
42{
44
46 (
50 );
51
53 (
56 displacement
57 );
58}
59
60
61// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
62
64(
65 const labelHashSet& changedFaces,
66 labelHashSet& affectedFaces
67)
68{
69 bitSet affectedPoints(mesh().nPoints());
70
71 for (const label facei : changedFaces)
72 {
73 affectedPoints.set(mesh().faces()[facei]);
74 }
75
77 (
78 mesh(),
79 affectedPoints,
80 orEqOp<unsigned int>(),
81 0U
82 );
83
84 for (const label pointi : affectedPoints)
85 {
86 for (const label celli : mesh().pointCells()[pointi])
87 {
88 affectedFaces.insert(mesh().cells()[celli]);
89 }
90 }
91}
92
93
95{
96 if
97 (
98 relaxationFactors_.empty()
99 || (relaxationFactors_.size() == 1 && relaxationFactors_[0] == 1.0)
100 )
101 {
102 relaxedPoints_ = points0() + pointDisplacement().internalField();
103 return true;
104 }
105
106
107 const pointField oldRelaxedPoints(relaxedPoints_);
108
109 labelHashSet affectedFaces(facesToMove_);
110
111 // Create a list of relaxation levels
112 // -1 indicates a point which is not to be moved
113 // 0 is the starting value for a moving point
114 labelList relaxationLevel(mesh().nPoints(), -1);
115 for (const label facei : affectedFaces)
116 {
117 for (const label pointi : mesh().faces()[facei])
118 {
119 relaxationLevel[pointi] = 0;
120 }
121 }
122
124 (
125 mesh(),
126 relaxationLevel,
128 label(-1)
129 );
130
131 // Loop whilst relaxation levels are being incremented
132 bool complete(false);
133 while (!complete)
134 {
135 //scalar nAffectedFaces(affectedFaces.size());
136 //reduce(nAffectedFaces, sumOp<scalar>());
137 //Info << " Moving " << nAffectedFaces << " faces" << endl;
138
139 // Move the points
140 forAll(relaxationLevel, pointI)
141 {
142 if (relaxationLevel[pointI] >= 0)
143 {
144 const scalar x
145 (
146 relaxationFactors_[relaxationLevel[pointI]]
147 );
148
149 relaxedPoints_[pointI] =
150 (1 - x)*oldRelaxedPoints[pointI]
151 + x*(points0()[pointI] + pointDisplacement()[pointI]);
152 }
153 }
154
155 // Get a list of changed faces
156 labelHashSet markedFaces;
157 markAffectedFaces(affectedFaces, markedFaces);
158 labelList markedFacesList(markedFaces.toc());
159
160 // Update the geometry
161 meshGeometry_.correct(relaxedPoints_, markedFacesList);
162
163 // Check the modified face quality
164 if (false)
165 {
166 // Use snappyHexMesh compatible checks
167 markedFaces.clear();
169 (
170 false,
171 meshQualityDict_,
172 meshGeometry_,
173 relaxedPoints_,
174 markedFacesList,
175 markedFaces
176 );
177
178 // Mark the affected faces
179 affectedFaces.clear();
180 markAffectedFaces(markedFaces, affectedFaces);
181 }
182 else
183 {
184 // Use pointSmoother specific
185 tmp<scalarField> tfaceQ
186 (
187 pointUntangler_->faceQuality
188 (
189 relaxedPoints_,
190 meshGeometry_.faceCentres(),
191 meshGeometry_.faceAreas(),
192 meshGeometry_.cellCentres(),
193 meshGeometry_.cellVolumes()
194 )
195 );
196
197 if (debug)
198 {
199 MinMax<scalar> range(gMinMax(tfaceQ()));
200 Pout<< " min:" << range.min() << nl
201 << " max:" << range.max() << endl;
202 }
203
204 labelList order;
205 Foam::sortedOrder(tfaceQ(), order);
206
207 label nUntangle = 0;
208 forAll(order, i)
209 {
210 if (tfaceQ()[order[i]] > untangleQ_)
211 {
212 nUntangle = i;
213 break;
214 }
215 }
216
217 affectedFaces = labelList(SubList<label>(order, nUntangle));
218 }
219
220 // Increase relaxation and check convergence
221 bitSet pointsToRelax(mesh().nPoints());
222 complete = true;
223 for (const label facei : affectedFaces)
224 {
225 pointsToRelax.set(mesh().faces()[facei]);
226 }
227
228 for (const label pointI : pointsToRelax)
229 {
230 if (relaxationLevel[pointI] < relaxationFactors_.size() - 1)
231 {
232 ++ relaxationLevel[pointI];
233
234 complete = false;
235 }
236 }
237
238 // Synchronise relaxation levels
240 (
241 mesh(),
242 relaxationLevel,
244 label(0)
245 );
246
247 // Synchronise completion
248 UPstream::reduceAnd(complete);
249 }
250
251 // Check for convergence
252 bool converged(true);
253 forAll(mesh().faces(), facei)
254 {
255 for (const label pointi : mesh().faces()[facei])
256 {
257 if (relaxationLevel[pointi] > 0)
258 {
259 facesToMove_.insert(facei);
260
261 converged = false;
262
263 break;
264 }
265 }
266 }
267
268 // Syncronise convergence
269 UPstream::reduceAnd(converged);
270
271 //if (converged)
272 //{
273 // Info<< "... Converged" << endl << endl;
274 //}
275 //else
276 //{
277 // Info<< "... Not converged" << endl << endl;
278 //}
279
280 return converged;
281}
282
283
285(
286 const dictionary& dict
287)
288{
289 if (dict.getOrDefault<bool>("moveInternalFaces", true))
290 {
291 facesToMove_.resize(2*mesh().nFaces());
292 forAll(mesh().faces(), faceI)
293 {
294 facesToMove_.insert(faceI);
295 }
296 }
297 else
298 {
299 facesToMove_.resize(2*(mesh().nBoundaryFaces()));
300 for
301 (
302 label faceI = mesh().nInternalFaces();
303 faceI < mesh().nFaces();
304 ++ faceI
305 )
307 facesToMove_.insert(faceI);
308 }
309 }
310}
311
312
314(
315 pointVectorField& pointDisplacement
316)
317{
318 // Assume empty point patches are already in correct location
319 // so knock out any off-plane displacement.
320 auto& fld = pointDisplacement.primitiveFieldRef();
321 for (const auto& ppf : pointDisplacement.boundaryField())
322 {
324 {
325 const auto& mp = ppf.patch().meshPoints();
326 forAll(mp, i)
327 {
329 ppf.patch().applyConstraint(i, pc);
330 fld[mp[i]] = pc.constrainDisplacement(fld[mp[i]]);
331 }
332 }
333 }
334
335 pointField wantedPoints(points0() + fld);
336 twoDCorrectPoints(wantedPoints);
337 fld = wantedPoints-points0();
338}
339
340
341// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
342
345(
346 const polyMesh& mesh,
347 const IOdictionary& dict
348)
349:
351 meshGeometry_(mesh),
352 pointUntangler_
353 (
355 (
356 coeffDict().get<word>("untangler"),
357 mesh,
358 coeffDict()
359 )
360 ),
361 untangleQ_(coeffDict().get<scalar>("untangleQ")),
362 minQ_(coeffDict().get<scalar>("minQ")),
363 pointSmoother_(pointSmoother::New(mesh, coeffDict())),
364 nPointSmootherIter_
365 (
366 coeffDict().get<label>("nPointSmootherIter")
367 ),
368 relaxedPoints_(mesh.points())
369{
370 if (coeffDict().readIfPresent("relaxationFactors", relaxationFactors_))
372 meshQualityDict_ = coeffDict().subDict("meshQuality");
373 }
375}
376
377
380(
381 const polyMesh& mesh,
382 const IOdictionary& dict,
383 const pointVectorField& pointDisplacement,
384 const pointIOField& points0
385)
386:
387 displacementMotionSolver(mesh, dict, pointDisplacement, points0, typeName),
388 meshGeometry_(mesh),
389 pointUntangler_
390 (
391 pointSmoother::New
392 (
393 coeffDict().get<word>("untangler"),
394 mesh,
395 coeffDict()
396 )
397 ),
398 untangleQ_(coeffDict().get<scalar>("untangleQ")),
399 minQ_(coeffDict().get<scalar>("minQ")),
400 pointSmoother_
401 (
402 pointSmoother::New
403 (
404 mesh,
405 coeffDict()
406 )
407 ),
408 nPointSmootherIter_
409 (
410 coeffDict().get<label>("nPointSmootherIter")
411 ),
412 relaxedPoints_(mesh.points())
413{
414 if (coeffDict().readIfPresent("relaxationFactors", relaxationFactors_))
415 {
416 meshQualityDict_ = coeffDict().subDict("meshQuality");
417 }
419}
420
421
422// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
423
427 //Note: twoDCorrect already done by ::solve
428
429 return relaxedPoints_;
430}
431
432
434{
435 movePoints(curPoints());
436
437 // Update values on pointDisplacement. Note: should also evaluate? Since
438 // e.g. cellMotionBC uses pointDisplacement value.
439 pointDisplacement().boundaryFieldRef().updateCoeffs();
440
441
442 fileName debugDir;
443 if (debug & 2)
444 {
445 debugDir = mesh().time().timePath();
446 mkDir(debugDir);
447 OBJstream os(debugDir/"bc.obj");
448
449 const pointField wantedPoints
450 (
451 points0()
452 + pointDisplacement().internalField()
453 );
454 const auto& pbm = pointDisplacement().mesh().boundary();
455 for (const auto& ppp : pbm)
456 {
457 if (!isA<emptyPointPatch>(ppp))
458 {
459 const auto& mp = ppp.meshPoints();
460 for (const label pointi : mp)
461 {
462 os.write
463 (
465 (
466 points0()[pointi],
467 wantedPoints[pointi]
468 )
469 );
470 }
471 }
472 }
473 Pout<< "Written " << os.nVertices() << " initial displacements to "
474 << os.name() << endl;
475 }
476
477
478 // Extend: face-to-point-to-cell-to-faces
479 labelHashSet affectedFaces;
480 markAffectedFaces(facesToMove_, affectedFaces);
481
482
483 for (label nIter = 0; nIter < nPointSmootherIter_; ++nIter)
484 {
485 const pointField wantedPoints
486 (
487 points0()
488 + pointDisplacement().internalField()
489 );
490
491 meshGeometry_.correct
492 (
493 wantedPoints,
494 affectedFaces.toc()
495 );
496
497 //{
498 // // Debugging: check meshGeometry consistent with fvGeometryScheme
499 // const auto& geom =
500 // reinterpret_cast<const fvMesh&>(mesh()).geometry();
501 // pointField faceCentres(mesh().nFaces());
502 // vectorField faceAreas(mesh().nFaces());
503 // pointField cellCentres(mesh().nCells());
504 // scalarField cellVolumes(mesh().nCells());
505 // geom.updateGeom
506 // (
507 // wantedPoints,
508 // mesh().points(), // old points
509 // faceCentres,
510 // faceAreas,
511 // cellCentres,
512 // cellVolumes
513 // );
514 // forAll(faceCentres, facei)
515 // {
516 // const point& meshFc = mesh().faceCentres()[facei];
517 // const point& meshGeomFc = meshGeometry_.faceCentres()[facei];
518 // const point& updatedFc = faceCentres[facei];
519 //
520 // if (updatedFc != meshGeomFc)
521 // {
522 // const face& f = mesh().faces()[facei];
523 //
524 // Pout<< "At face:" << facei << nl
525 // << " old :" << meshFc << nl
526 // << " new :" << updatedFc << nl
527 // << " polyMeshGeom:" << meshGeomFc << nl
528 // << " oldPoints :"
529 // << UIndirectList<point>(mesh().points(), f) << nl
530 // << " wantedPoints:"
531 // << UIndirectList<point>(wantedPoints, f) << nl
532 // << endl;
533 // }
534 // }
535 //}
536
537 // Get measure of face quality
538 tmp<scalarField> tfaceQ
539 (
540 pointUntangler_->faceQuality
541 (
542 wantedPoints,
543 meshGeometry_.faceCentres(),
544 meshGeometry_.faceAreas(),
545 meshGeometry_.cellCentres(),
546 meshGeometry_.cellVolumes()
547 )
548 );
549
550
551 if (debug)
552 {
553 MinMax<scalar> range(gMinMax(tfaceQ()));
554 Pout<< " min:" << range.min() << nl
555 << " max:" << range.max() << endl;
556 }
557
558 labelList order;
559 Foam::sortedOrder(tfaceQ(), order);
560
561 label nUntangle = 0;
562 forAll(order, i)
563 {
564 if (tfaceQ()[order[i]] > untangleQ_)
565 {
566 nUntangle = i;
567 break;
568 }
569 }
570 label nLow = 0;
571 forAll(order, i)
572 {
573 if (tfaceQ()[order[i]] > minQ_)
574 {
575 nLow = i;
576 break;
577 }
578 }
579
580
581 if (debug)
582 {
583 Pout<< " nUntangle:" << returnReduce(nUntangle, sumOp<label>())
584 << nl
585 << " nLow :" << returnReduce(nLow, sumOp<label>())
586 << nl;
587 }
588
589
590 if (returnReduce(nUntangle, sumOp<label>()))
591 {
592 // Start untangling
593 labelList lowQFaces(SubList<label>(order, nUntangle));
594 //{
595 // // Grow set (non parallel)
596 // bitSet isMarkedFace(mesh().nFaces());
597 // for (const label facei : lowQFaces)
598 // {
599 // for (const label pointi : mesh().faces()[facei])
600 // {
601 // isMarkedFace.set(mesh().pointFaces()[pointi]);
602 // }
603 // }
604 // lowQFaces = isMarkedFace.sortedToc();
605 //}
606
607 //Pout<< " untangling "
608 // << returnReduce(lowQFaces.size(), sumOp<label>())
609 // << " faces" << endl;
610 pointUntangler_->update
611 (
612 lowQFaces,
613 points0(),
614 wantedPoints,
615 meshGeometry_,
616 pointDisplacement()
617 //false // ! do NOT apply bcs, constraints
618 );
619
620 // Keep points on empty patches. Note: since pointConstraints
621 // does not implement constraints on emptyPointPatches and
622 // emptyPointPatchField does not either.
623 emptyCorrectPoints(pointDisplacement());
624
625 if (debug & 2)
626 {
627 OBJstream os(debugDir/"untangle_" + Foam::name(nIter) + ".obj");
628
629 const pointField wantedPoints
630 (
631 points0()
632 + pointDisplacement().internalField()
633 );
634 forAll(wantedPoints, pointi)
635 {
636 os.write
637 (
639 (
640 points0()[pointi],
641 wantedPoints[pointi]
642 )
643 );
644 }
645 Pout<< "Written " << os.nVertices() << " wanted untangle to "
646 << os.name() << endl;
647 }
648 }
649 else if (returnReduce(nLow, sumOp<label>()))
650 {
651 labelList lowQFaces(SubList<label>(order, nLow));
652 //{
653 // // Grow set (non parallel)
654 // bitSet isMarkedFace(mesh().nFaces());
655 // for (const label facei : lowQFaces)
656 // {
657 // for (const label pointi : mesh().faces()[facei])
658 // {
659 // isMarkedFace.set(mesh().pointFaces()[pointi]);
660 // }
661 // }
662 // lowQFaces = isMarkedFace.sortedToc();
663 //}
664
665 //Pout<< " smoothing "
666 // << returnReduce(lowQFaces.size(), sumOp<label>())
667 // << " faces" << endl;
668
669 pointSmoother_->update
670 (
671 lowQFaces,
672 points0(),
673 wantedPoints,
674 meshGeometry_,
675 pointDisplacement()
676 );
677 // Keep points on empty patches
678 emptyCorrectPoints(pointDisplacement());
679 }
680 else
681 {
682 //Pout<< "** converged" << endl;
683 break;
684 }
685 }
686
687
688 relax();
689 //relaxedPoints_ = points0() + pointDisplacement().internalField();
690
691 twoDCorrectPoints(relaxedPoints_);
692
693 // Update pointDisplacement for actual relaxedPoints. Keep fixed-value
694 // bcs.
695 pointDisplacement().primitiveFieldRef() = relaxedPoints_-points0();
696
697 // Adhere to multi-point constraints. Does correctBoundaryConditions +
698 // multi-patch issues.
699 const pointConstraints& pcs =
700 pointConstraints::New(pointDisplacement().mesh());
701 pcs.constrainDisplacement(pointDisplacement(), false);
702}
703
704
705// ************************************************************************* //
scalar range
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
propsDict readIfPresent("fields", acceptFields)
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
const polyBoundaryMesh & pbm
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition HashSet.H:229
List< Key > toc() const
The table of contents (the keys) in unsorted order.
Definition HashTable.C:141
void clear()
Remove all entries from table.
Definition HashTable.C:742
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
static FOAM_NO_DANGLING_REFERENCE const pointConstraints & New(const pointMesh &mesh, Args &&... args)
A min/max value pair with additional methods. In addition to conveniently storing values,...
Definition MinMax.H:106
An OFstream that keeps track of vertices and provides convenience output methods for OBJ files.
Definition OBJstream.H:58
label nVertices() const noexcept
Return the number of vertices written.
Definition OBJstream.H:120
virtual Ostream & write(const char c) override
Write character.
Definition OBJstream.C:69
virtual const fileName & name() const override
Read/write access to the name of the stream.
Definition OSstream.H:134
A non-owning sub-view of a List (allocated or unallocated storage).
Definition SubList.H:61
fileName timePath() const
Return current time path = path/timeName.
Definition Time.H:512
static void reduceAnd(bool &value, const int communicator=worldComm)
Logical (and) reduction (MPI_AllReduce).
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition bitSet.H:61
void set(const bitSet &bitset)
Set specified bits from another bitset.
Definition bitSetI.H:502
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition dictionary.C:441
Virtual base class for displacement motion solver.
pointVectorField & pointDisplacement() noexcept
Return reference to the point motion displacement field.
displacementMotionSolver(const displacementMotionSolver &)=delete
No copy construct.
static autoPtr< displacementMotionSolver > New(const word &solverTypeName, const polyMesh &, const IOdictionary &, const pointVectorField &pointDisplacement, const pointIOField &points0)
Select constructed from polyMesh, dictionary and components.
Quality-based under-relaxation for run-time selectable point smoothing. WIP.
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
labelHashSet facesToMove_
Set of the faces which are to be moved.
scalarList relaxationFactors_
Relaxation factors to use in each iteration.
void markAffectedFaces(const labelHashSet &changedFaces, labelHashSet &affectedFaces)
Mark affected faces.
const label nPointSmootherIter_
Number of point smoother iterations per timestep.
void emptyCorrectPoints(pointVectorField &pointDisplacement)
Handle 2D & empty bcs. Assume in both cases the starting mesh.
virtual void setFacesToMove(const dictionary &)
Set all the faces to be moved.
displacementSmartPointSmoothingMotionSolver(const polyMesh &, const IOdictionary &)
Construct from a polyMesh and an IOdictionary.
A class for handling file names.
Definition fileName.H:75
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
static bool checkMesh(const bool report, const polyMesh &mesh, const dictionary &dict, labelHashSet &wrongFaces, const bool dryRun=false)
Check mesh with mesh settings in dict. Collects incorrect faces.
Virtual base class for mesh motion solver.
const polyMesh & mesh() const
Return reference to mesh.
virtual void twoDCorrectPoints(pointField &) const
const dictionary & coeffDict() const
Const access to the coefficients dictionary.
Smooth ATC in cells having a point to a set of patches supplied by type.
Definition pointCells.H:55
Accumulates point constraints through successive applications of the applyConstraint function.
void applyConstraint(const vector &cd)
Apply and accumulate the effect of the given constraint direction.
vector constrainDisplacement(const vector &disp) const
Constrain a displacement.
Application of (multi-)patch point constraints.
void constrainDisplacement(pointVectorField &displacement, const bool overrideValue=false) const
Apply boundary conditions (single-patch constraints),.
Abstract base class for point smoothing methods. Handles parallel communication via reset and average...
pointField & points0() noexcept
Return reference to the reference ('0') pointField.
virtual void movePoints(const pointField &)
Update local data for geometry changes.
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
label nFaces() const noexcept
Number of mesh faces.
static void syncPointList(const polyMesh &mesh, List< T > &pointValues, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh points.
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
U
Definition pEqn.H:72
UEqn relax()
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
const pointField & points
label nPoints
const cellShapeList & cells
const dimensionedScalar mp
Proton mass.
Namespace for handling debugging switches.
Definition debug.C:45
Namespace for OpenFOAM.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
List< label > labelList
A List of labels.
Definition List.H:62
bool mkDir(const fileName &pathName, mode_t mode=0777)
Make a directory and return an error if it could not be created.
Definition POSIX.C:616
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
GeometricField< vector, pointPatchField, pointMesh > pointVectorField
vectorIOField pointIOField
pointIOField is a vectorIOField.
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
labelList sortedOrder(const UList< T > &input)
Return the (stable) sort order for the list.
line< point, const point & > linePointRef
A line using referred points.
Definition line.H:66
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
pointField points0(pointIOField(IOobject("points", mesh.time().constant(), polyMesh::meshSubDir, mesh, IOobject::MUST_READ, IOobject::NO_WRITE, IOobject::NO_REGISTER)))