displacementLayeredMotionMotionSolver.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  Copyright (C) 2015-2019 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 
32 #include "pointFields.H"
33 #include "PointEdgeWave.H"
34 #include "syncTools.H"
35 #include "interpolationTable.H"
36 #include "mapPolyMesh.H"
37 #include "pointConstraints.H"
38 
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 
41 namespace Foam
42 {
43  defineTypeNameAndDebug(displacementLayeredMotionMotionSolver, 0);
44 
46  (
47  motionSolver,
48  displacementLayeredMotionMotionSolver,
49  dictionary
50  );
51 
53  (
54  displacementMotionSolver,
55  displacementLayeredMotionMotionSolver,
56  displacement
57  );
58 }
59 
60 
61 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
62 
64 Foam::displacementLayeredMotionMotionSolver::getCylindrical
65 (
66  const label cellZoneI,
67  const dictionary& zoneDict
68 )
69 {
70  auto iter = cylSystems_.cfind(cellZoneI);
71 
72  if (iter.found())
73  {
74  return *(iter.val());
75  }
76 
77  cylSystems_.set(cellZoneI, new coordSystem::cylindrical(zoneDict));
78 
79  return *cylSystems_[cellZoneI];
80 }
81 
82 
83 void Foam::displacementLayeredMotionMotionSolver::calcZoneMask
84 (
85  const label cellZoneI,
86  bitSet& isZonePoint,
87  bitSet& isZoneEdge
88 ) const
89 {
90  isZonePoint.resize(mesh().nPoints());
91  isZoneEdge.resize(mesh().nEdges());
92 
93  if (cellZoneI == -1)
94  {
95  isZonePoint = true;
96  isZoneEdge = true;
97  return;
98  }
99 
100 
101  isZonePoint.reset();
102  isZoneEdge.reset();
103 
104  const cellZone& cz = mesh().cellZones()[cellZoneI];
105 
106  // Mark points, edges inside cellZone
107  for (const label celli : cz)
108  {
109  isZonePoint.set(mesh().cellPoints(celli));
110  isZoneEdge.set(mesh().cellEdges(celli));
111  }
112 
114  (
115  mesh(),
116  isZonePoint,
117  orEqOp<unsigned int>(),
118  0
119  );
120 
122  (
123  mesh(),
124  isZoneEdge,
125  orEqOp<unsigned int>(),
126  0
127  );
128 
129  DebugInfo
130  << "On cellZone " << cz.name() << " marked "
131  << returnReduce(isZonePoint.count(), sumOp<label>()) << " points and "
132  << returnReduce(isZoneEdge.count(), sumOp<label>()) << " edges" << nl;
133 }
134 
135 
136 // Find distance to starting point
137 void Foam::displacementLayeredMotionMotionSolver::walkStructured
138 (
139  const label cellZoneI,
140  const bitSet& isZonePoint,
141  const bitSet& isZoneEdge,
142  const labelList& seedPoints,
143  const vectorField& seedData,
145  vectorField& data,
146  labelField& patchPoints
147 ) const
148 {
149  List<pointEdgeStructuredWalk> seedInfo(seedPoints.size());
150 
151  forAll(seedPoints, i)
152  {
153  const label pointi = seedPoints[i];
154 
155  seedInfo[i] = pointEdgeStructuredWalk
156  (
157  points0()[pointi], // location of data
158  points0()[pointi], // previous location
159  0.0,
160  seedData[i],
161  pointi // pass thru seed point id
162  );
163  }
164 
165  // Current info on points
166  List<pointEdgeStructuredWalk> allPointInfo(mesh().nPoints());
167 
168  // Mark points inside cellZone.
169  // Note that we use points0, not mesh.points()
170  // so as not to accumulate errors.
171  for (const label pointi : isZonePoint)
172  {
173  allPointInfo[pointi] = pointEdgeStructuredWalk
174  (
175  points0()[pointi], // location of data
176  point::max, // not valid
177  0.0,
178  Zero // passive data
179  );
180  }
181 
182 
183  // Current info on edges
184  List<pointEdgeStructuredWalk> allEdgeInfo(mesh().nEdges());
185 
186  // Mark edges inside cellZone
187  for (const label edgei : isZoneEdge)
188  {
189  allEdgeInfo[edgei] = pointEdgeStructuredWalk
190  (
191  mesh().edges()[edgei].centre(points0()), // location of data
192  point::max, // not valid
193  0.0,
194  Zero
195  );
196  }
197 
198  // Walk
199  PointEdgeWave<pointEdgeStructuredWalk> wallCalc
200  (
201  mesh(),
202  seedPoints,
203  seedInfo,
204 
205  allPointInfo,
206  allEdgeInfo,
207  mesh().globalData().nTotalPoints() // max iterations
208  );
209 
210  // Extract distance and passive data
211  for (const label pointi : isZonePoint)
212  {
213  const auto& pointInfo = allPointInfo[pointi];
214 
215  distance[pointi] = pointInfo.dist();
216  data[pointi] = pointInfo.data();
217 
218  // Optional information
219  if (patchPoints.size())
220  {
221  patchPoints[pointi] = pointInfo.index();
222  }
223  }
224 }
225 
226 
227 // Evaluate faceZone patch
229 Foam::displacementLayeredMotionMotionSolver::faceZoneEvaluate
230 (
231  const faceZone& fz,
232  const labelList& meshPoints,
233  const dictionary& dict,
234  const PtrList<pointVectorField>& patchDisp,
235  const label patchi
236 ) const
237 {
238  auto tfld = tmp<vectorField>::New(meshPoints.size());
239  auto& fld = tfld.ref();
240 
241  const word type(dict.get<word>("type"));
242 
243  if (type == "fixedValue")
244  {
245  fld = vectorField("value", dict, meshPoints.size());
246  }
247  else if (type == "timeVaryingUniformFixedValue")
248  {
249  interpolationTable<vector> timeSeries(dict);
250 
251  fld = timeSeries(mesh().time().timeOutputValue());
252  }
253  else if (type == "slip")
254  {
255  if ((patchi % 2) != 1)
256  {
258  << "FaceZone:" << fz.name()
259  << exit(FatalIOError);
260  }
261  // Use field set by previous bc
262  fld = vectorField(patchDisp[patchi - 1], meshPoints);
263  }
264  else if (type == "follow")
265  {
266  // Only on boundary faces - follow boundary conditions
267  fld = vectorField(pointDisplacement_, meshPoints);
268  }
269  else if (type == "uniformFollow")
270  {
271  // Reads name of name of patch. Then get average point displacement on
272  // patch. That becomes the value of fld.
273  const word patchName(dict.get<word>("patch"));
274  const label patchID = mesh().boundaryMesh().findPatchID(patchName);
275  pointField pdf
276  (
277  pointDisplacement_.boundaryField()[patchID].patchInternalField()
278  );
279  fld = gAverage(pdf);
280  }
281  else
282  {
284  << "Unknown faceZonePatch type " << type
285  << " for faceZone " << fz.name() << nl
286  << exit(FatalIOError);
287  }
288 
289  return tfld;
290 }
291 
292 
293 void Foam::displacementLayeredMotionMotionSolver::cellZoneSolve
294 (
295  const label cellZoneI,
296  const dictionary& zoneDict
297 )
298 {
299  bitSet isZonePoint, isZoneEdge;
300  calcZoneMask(cellZoneI, isZonePoint, isZoneEdge);
301 
302  const dictionary& patchesDict = zoneDict.subDict("boundaryField");
303 
304  if (patchesDict.size() != 2)
305  {
307  << "Two faceZones (patches) must be specified per cellZone. "
308  << " cellZone:" << cellZoneI
309  << " patches:" << patchesDict.toc()
310  << exit(FatalIOError);
311  }
312 
313  PtrList<labelField> patchPoints(patchesDict.size());
314  PtrList<scalarField> patchDist(patchesDict.size());
315  PtrList<pointVectorField> patchDisp(patchesDict.size());
316 
317  // Allocate the fields
318  label patchi = 0;
319  for (const entry& dEntry : patchesDict)
320  {
321  const word& faceZoneName = dEntry.keyword();
322  const label zoneI = mesh().faceZones().findZoneID(faceZoneName);
323  if (zoneI == -1)
324  {
326  << "Cannot find faceZone " << faceZoneName
327  << endl << "Valid zones are " << mesh().faceZones().names()
328  << exit(FatalIOError);
329  }
330 
331  // Determine the points of the faceZone within the cellZone
332  const faceZone& fz = mesh().faceZones()[zoneI];
333 
334  patchPoints.set(patchi, new labelField(mesh().nPoints(), label(-1)));
335  patchDist.set(patchi, new scalarField(mesh().nPoints()));
336  patchDisp.set
337  (
338  patchi,
339  new pointVectorField
340  (
341  IOobject
342  (
343  mesh().cellZones()[cellZoneI].name() + "_" + fz.name(),
344  mesh().time().timeName(),
345  mesh(),
348  false
349  ),
350  pointDisplacement_ // to inherit the boundary conditions
351  )
352  );
353 
354  ++patchi;
355  }
356 
357 
358 
359  // 'correctBoundaryConditions'
360  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~
361  // Loops over all the faceZones and walks their boundary values
362 
363  // Make sure we can pick up bc values from field
364  pointDisplacement_.correctBoundaryConditions();
365 
366  patchi = 0;
367  for (const entry& dEntry : patchesDict)
368  {
369  const word& faceZoneName = dEntry.keyword();
370  const dictionary& faceZoneDict = dEntry.dict();
371 
372  // Determine the points of the faceZone within the cellZone
373  const faceZone& fz = mesh().faceZones()[faceZoneName];
374  const labelList& fzMeshPoints = fz().meshPoints();
375  DynamicList<label> meshPoints(fzMeshPoints.size());
376  forAll(fzMeshPoints, i)
377  {
378  if (isZonePoint[fzMeshPoints[i]])
379  {
380  meshPoints.append(fzMeshPoints[i]);
381  }
382  }
383 
384  // Get initial value for all the faceZone points
385  tmp<vectorField> tseed = faceZoneEvaluate
386  (
387  fz,
388  meshPoints,
389  faceZoneDict,
390  patchDisp,
391  patchi
392  );
393 
394  DebugInfo
395  << "For cellZone:" << cellZoneI
396  << " for faceZone:" << fz.name()
397  << " nPoints:" << tseed().size()
398  << " have patchField:"
399  << " max:" << gMax(tseed())
400  << " min:" << gMin(tseed())
401  << " avg:" << gAverage(tseed())
402  << endl;
403 
404 
405  // Set distance and transported value
406  walkStructured
407  (
408  cellZoneI,
409  isZonePoint,
410  isZoneEdge,
411 
412  meshPoints,
413  tseed,
414  patchDist[patchi],
415  patchDisp[patchi],
416  patchPoints[patchi]
417  );
418 
419  // Implement real bc.
420  patchDisp[patchi].correctBoundaryConditions();
421 
422  ++patchi;
423  }
424 
425 
426  // Solve
427  // ~~~~~
428 
429  if (debug)
430  {
431  // Normalised distance
433  (
434  IOobject
435  (
436  mesh().cellZones()[cellZoneI].name() + ":distance",
437  mesh().time().timeName(),
438  mesh(),
441  false
442  ),
443  pointMesh::New(mesh()),
445  );
446 
447  for (const label pointi : isZonePoint)
448  {
449  const scalar d1 = patchDist[0][pointi];
450  const scalar d2 = patchDist[1][pointi];
451  if (d1 + d2 > SMALL)
452  {
453  distance[pointi] = d1/(d1 + d2);
454  }
455  }
456 
457  Info<< "Writing " << pointScalarField::typeName
458  << distance.name() << " to "
459  << mesh().time().timeName() << endl;
460  distance.write();
461  }
462 
463 
464  const word interpolationScheme(zoneDict.get<word>("interpolationScheme"));
465 
466  if (interpolationScheme == "oneSided")
467  {
468  for (const label pointi : isZonePoint)
469  {
470  pointDisplacement_[pointi] = patchDisp[0][pointi];
471  }
472  }
473  else if (interpolationScheme == "linear")
474  {
475  for (const label pointi : isZonePoint)
476  {
477  const scalar d1 = patchDist[0][pointi];
478  const scalar d2 = patchDist[1][pointi];
479  const scalar s = d1/(d1 + d2 + VSMALL);
480 
481  const vector& pd1 = patchDisp[0][pointi];
482  const vector& pd2 = patchDisp[1][pointi];
483 
484  pointDisplacement_[pointi] = (1 - s)*pd1 + s*pd2;
485  }
486  }
487  else if (interpolationScheme == "cylindrical")
488  {
489  const coordSystem::cylindrical& cs =
490  this->getCylindrical(cellZoneI, zoneDict);
491 
492  // May wish to implement alternative distance calculation here
493 
494 
496  << "cylindrical interpolation not yet available" << nl
497  << "coordinate system " << cs << nl
498  << exit(FatalError);
499  }
500  else
501  {
503  << "Invalid interpolationScheme: " << interpolationScheme
504  << ". Valid schemes: (oneSided linear cylindrical)" << nl
505  << exit(FatalError);
506  }
507 }
508 
509 
510 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
511 
512 Foam::displacementLayeredMotionMotionSolver::
513 displacementLayeredMotionMotionSolver
514 (
515  const polyMesh& mesh,
516  const IOdictionary& dict
517 )
518 :
520 {}
521 
522 
523 Foam::displacementLayeredMotionMotionSolver::
524 displacementLayeredMotionMotionSolver
525 (
526  const polyMesh& mesh,
527  const IOdictionary& dict,
528  const pointVectorField& pointDisplacement,
529  const pointIOField& points0
530 )
531 :
532  displacementMotionSolver(mesh, dict, pointDisplacement, points0, typeName)
533 {}
534 
535 
536 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
537 
540 {
541  tmp<pointField> tcurPoints
542  (
543  points0() + pointDisplacement_.primitiveField()
544  );
545 
546  return tcurPoints;
547 }
548 
549 
551 {
552  // The points have moved so before interpolation update the motionSolver
553  movePoints(mesh().points());
554 
555  // Apply boundary conditions
556  pointDisplacement_.boundaryFieldRef().updateCoeffs();
557 
558  // Solve motion on all regions (=cellZones)
559  for (const entry& dEntry : coeffDict().subDict("regions"))
560  {
561  const word& cellZoneName = dEntry.keyword();
562  const dictionary& regionDict = dEntry.dict();
563 
564  const label zoneI = mesh().cellZones().findZoneID(cellZoneName);
565 
566  Info<< "solving for zone: " << cellZoneName << endl;
567 
568  if (zoneI == -1)
569  {
571  << "Cannot find cellZone " << cellZoneName
572  << endl << "Valid zones are " << mesh().cellZones().names()
573  << exit(FatalIOError);
574  }
575 
576  cellZoneSolve(zoneI, regionDict);
577  }
578 
579  // Update pointDisplacement for solved values
580  const pointConstraints& pcs =
581  pointConstraints::New(pointDisplacement_.mesh());
582  pcs.constrainDisplacement(pointDisplacement_, false);
583 }
584 
585 
587 (
588  const mapPolyMesh& mpm
589 )
590 {
592 
593  const vectorField displacement(this->newPoints() - points0_);
594 
595  forAll(points0_, pointi)
596  {
597  const label oldPointi = mpm.pointMap()[pointi];
598 
599  if (oldPointi >= 0)
600  {
601  const label masterPointi = mpm.reversePointMap()[oldPointi];
602 
603  if ((masterPointi != pointi))
604  {
605  // newly inserted point in this cellZone
606 
607  // need to set point0 so that it represents the position that
608  // it would have had if it had existed for all time
609  points0_[pointi] -= displacement[pointi];
610  }
611  }
612  }
613 }
614 
615 
616 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
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::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
Foam::pointVectorField
GeometricField< vector, pointPatchField, pointMesh > pointVectorField
Definition: pointFieldsFwd.H:61
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:54
Foam::syncTools::syncEdgeList
static void syncEdgeList(const polyMesh &mesh, List< T > &edgeValues, const CombineOp &cop, const T &nullValue, const TransformOp &top, const FlipOp &fop)
Synchronize values on all mesh edges.
Definition: syncToolsTemplates.C:800
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::displacementLayeredMotionMotionSolver::curPoints
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
Definition: displacementLayeredMotionMotionSolver.C:539
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::IOField
A primitive field of type <T> with automated input and output.
Definition: foamVtkLagrangianWriter.H:61
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:94
s
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Definition: gmvOutputSpray.H:25
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::gAverage
Type gAverage(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:604
mapPolyMesh.H
Foam::displacementLayeredMotionMotionSolver::updateMesh
virtual void updateMesh(const mapPolyMesh &mpm)
Update topology.
Definition: displacementLayeredMotionMotionSolver.C:587
Foam::MeshObject< polyMesh, UpdateableMeshObject, pointMesh >::New
static const pointMesh & New(const polyMesh &mesh, Args &&... args)
Get existing or create a new MeshObject.
Definition: MeshObject.C:48
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
interpolationTable.H
Foam::FatalIOError
IOerror FatalIOError
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:444
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:107
Foam::pointConstraints::constrainDisplacement
void constrainDisplacement(pointVectorField &displacement, const bool overrideValue=false) const
Apply boundary conditions (single-patch constraints),.
Definition: pointConstraints.C:380
syncTools.H
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
Foam::syncTools::syncPointList
static void syncPointList(const polyMesh &mesh, List< T > &pointValues, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh points.
Definition: syncToolsTemplates.C:721
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
nPoints
label nPoints
Definition: gmvOutputHeader.H:2
displacementLayeredMotionMotionSolver.H
Foam::Field< vector >
Foam::labelField
Field< label > labelField
Specialisation of Field<T> for label.
Definition: labelField.H:52
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
pointConstraints.H
PointEdgeWave.H
Foam::polyMesh::cellZones
const cellZoneMesh & cellZones() const noexcept
Return cell zone mesh.
Definition: polyMesh.H:492
patchID
label patchID
Definition: boundaryProcessorFaPatchPoints.H:5
Foam::polyMesh::faceZones
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition: polyMesh.H:486
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
Foam::polyBoundaryMesh::findPatchID
label findPatchID(const word &patchName, const bool allowNotFound=true) const
Find patch index given a name, return -1 if not found.
Definition: polyBoundaryMesh.C:765
fld
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Definition: gmvOutputLagrangian.H:23
Foam::coordSystem::cylindrical
A cylindrical coordinate system (r-theta-z). The coordinate system angle theta is always in radians.
Definition: cylindricalCS.H:71
timeName
word timeName
Definition: getTimeIndex.H:3
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::displacementLayeredMotionMotionSolver::solve
virtual void solve()
Solve for motion.
Definition: displacementLayeredMotionMotionSolver.C:550
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::displacementMotionSolver
Virtual base class for displacement motion solver.
Definition: displacementMotionSolver.H:53
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::ZoneMesh::findZoneID
label findZoneID(const word &zoneName) const
Find zone index by name, return -1 if not found.
Definition: ZoneMesh.C:519
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
Foam::distance
scalar distance(const vector &p1, const vector &p2)
Definition: curveTools.C:12
Foam::ZoneMesh::names
wordList names() const
A list of the zone names.
Definition: ZoneMesh.C:305
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::pointConstraints
Application of (multi-)patch point constraints.
Definition: pointConstraints.H:64
points0
pointField points0(pointIOField(IOobject("points", mesh.time().constant(), polyMesh::meshSubDir, mesh, IOobject::MUST_READ, IOobject::NO_WRITE, false)))
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
DebugInfo
#define DebugInfo
Report an information message using Foam::Info.
Definition: messageStream.H:382
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::VectorSpace< Vector< Cmpt >, Cmpt, 3 >::max
static const Vector< Cmpt > max
Definition: VectorSpace.H:117
Foam::points0MotionSolver::updateMesh
virtual void updateMesh(const mapPolyMesh &)
Update local data for topology changes.
Definition: points0MotionSolver.C:157
Foam::mapPolyMesh::reversePointMap
const labelList & reversePointMap() const
Reverse point map.
Definition: mapPolyMesh.H:469
pointEdgeStructuredWalk.H
Foam::mapPolyMesh::pointMap
const labelList & pointMap() const
Old point map.
Definition: mapPolyMesh.H:396
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::tmp::New
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:161
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:280
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:473
Foam::gMin
Type gMin(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:593
Foam::GeometricField< vector, pointPatchField, pointMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:188
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
pointFields.H
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:592
Foam::pointScalarField
GeometricField< scalar, pointPatchField, pointMesh > pointScalarField
Definition: pointFieldsFwd.H:56