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 | Copyright (C) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation | Copyright (C) 2015 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
29 #include "pointFields.H"
30 #include "PointEdgeWave.H"
31 #include "syncTools.H"
32 #include "interpolationTable.H"
33 #include "mapPolyMesh.H"
34 #include "pointConstraints.H"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40  defineTypeNameAndDebug(displacementLayeredMotionMotionSolver, 0);
41 
43  (
44  motionSolver,
45  displacementLayeredMotionMotionSolver,
46  dictionary
47  );
48 
50  (
51  displacementMotionSolver,
52  displacementLayeredMotionMotionSolver,
53  displacement
54  );
55 }
56 
57 
58 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
59 
61 (
62  const label cellZoneI,
63  PackedBoolList& isZonePoint,
64  PackedBoolList& isZoneEdge
65 ) const
66 {
67  if (cellZoneI == -1)
68  {
69  isZonePoint.setSize(mesh().nPoints());
70  isZonePoint = 1;
71 
72  isZoneEdge.setSize(mesh().nEdges());
73  isZoneEdge = 1;
74  }
75  else
76  {
77  const cellZone& cz = mesh().cellZones()[cellZoneI];
78 
79  label nPoints = 0;
80  forAll(cz, i)
81  {
82  const labelList& cPoints = mesh().cellPoints(cz[i]);
83  forAll(cPoints, cPointI)
84  {
85  if (!isZonePoint[cPoints[cPointI]])
86  {
87  isZonePoint[cPoints[cPointI]] = 1;
88  nPoints++;
89  }
90  }
91  }
92  syncTools::syncPointList
93  (
94  mesh(),
95  isZonePoint,
97  0
98  );
99 
100 
101  // Mark edge inside cellZone
102  label nEdges = 0;
103  forAll(cz, i)
104  {
105  const labelList& cEdges = mesh().cellEdges(cz[i]);
106  forAll(cEdges, cEdgeI)
107  {
108  if (!isZoneEdge[cEdges[cEdgeI]])
109  {
110  isZoneEdge[cEdges[cEdgeI]] = 1;
111  nEdges++;
112  }
113  }
114  }
115  syncTools::syncEdgeList
116  (
117  mesh(),
118  isZoneEdge,
120  0
121  );
122 
123  if (debug)
124  {
125  Info<< "On cellZone " << cz.name()
126  << " marked " << returnReduce(nPoints, sumOp<label>())
127  << " points and " << returnReduce(nEdges, sumOp<label>())
128  << " edges." << endl;
129  }
130  }
131 }
132 
133 
134 // Find distance to starting point
136 (
137  const label cellZoneI,
138  const PackedBoolList& isZonePoint,
139  const PackedBoolList& isZoneEdge,
140  const labelList& seedPoints,
141  const vectorField& seedData,
144 ) const
145 {
146  List<pointEdgeStructuredWalk> seedInfo(seedPoints.size());
147 
148  forAll(seedPoints, i)
149  {
150  seedInfo[i] = pointEdgeStructuredWalk
151  (
152  points0()[seedPoints[i]], // location of data
153  points0()[seedPoints[i]], // previous location
154  0.0,
155  seedData[i]
156  );
157  }
158 
159  // Current info on points
160  List<pointEdgeStructuredWalk> allPointInfo(mesh().nPoints());
161 
162  // Mark points inside cellZone.
163  // Note that we use points0, not mesh.points()
164  // so as not to accumulate errors.
165  forAll(isZonePoint, pointI)
166  {
167  if (isZonePoint[pointI])
168  {
169  allPointInfo[pointI] = pointEdgeStructuredWalk
170  (
171  points0()[pointI], // location of data
172  vector::max, // not valid
173  0.0,
174  vector::zero // passive data
175  );
176  }
177  }
178 
179  // Current info on edges
180  List<pointEdgeStructuredWalk> allEdgeInfo(mesh().nEdges());
181 
182  // Mark edges inside cellZone
183  forAll(isZoneEdge, edgeI)
184  {
185  if (isZoneEdge[edgeI])
186  {
187  allEdgeInfo[edgeI] = pointEdgeStructuredWalk
188  (
189  mesh().edges()[edgeI].centre(points0()), // location of data
190  vector::max, // not valid
191  0.0,
192  vector::zero
193  );
194  }
195  }
196 
197  // Walk
199  (
200  mesh(),
201  seedPoints,
202  seedInfo,
203 
204  allPointInfo,
205  allEdgeInfo,
206  mesh().globalData().nTotalPoints() // max iterations
207  );
208 
209  // Extract distance and passive data
210  forAll(allPointInfo, pointI)
211  {
212  if (isZonePoint[pointI])
213  {
214  distance[pointI] = allPointInfo[pointI].dist();
215  data[pointI] = allPointInfo[pointI].data();
216  }
217  }
218 }
219 
220 
221 // Evaluate faceZone patch
224 (
225  const faceZone& fz,
226  const labelList& meshPoints,
227  const dictionary& dict,
228  const PtrList<pointVectorField>& patchDisp,
229  const label patchI
230 ) const
231 {
232  tmp<vectorField> tfld(new vectorField(meshPoints.size()));
233  vectorField& fld = tfld();
234 
235  const word type(dict.lookup("type"));
236 
237  if (type == "fixedValue")
238  {
239  fld = vectorField("value", dict, meshPoints.size());
240  }
241  else if (type == "timeVaryingUniformFixedValue")
242  {
243  interpolationTable<vector> timeSeries(dict);
244 
245  fld = timeSeries(mesh().time().timeOutputValue());
246  }
247  else if (type == "slip")
248  {
249  if ((patchI % 2) != 1)
250  {
252  << "FaceZone:" << fz.name()
253  << exit(FatalIOError);
254  }
255  // Use field set by previous bc
256  fld = vectorField(patchDisp[patchI - 1], meshPoints);
257  }
258  else if (type == "follow")
259  {
260  // Only on boundary faces - follow boundary conditions
261  fld = vectorField(pointDisplacement_, meshPoints);
262  }
263  else if (type == "uniformFollow")
264  {
265  // Reads name of name of patch. Then get average point dislacement on
266  // patch. That becomes the value of fld.
267  const word patchName(dict.lookup("patch"));
268  label patchID = mesh().boundaryMesh().findPatchID(patchName);
269  pointField pdf
270  (
271  pointDisplacement_.boundaryField()[patchID].patchInternalField()
272  );
273  fld = gAverage(pdf);
274  }
275  else
276  {
278  << "Unknown faceZonePatch type " << type << " for faceZone "
279  << fz.name() << exit(FatalIOError);
280  }
281  return tfld;
282 }
283 
284 
286 (
287  const label cellZoneI,
288  const dictionary& zoneDict
289 )
290 {
291  PackedBoolList isZonePoint(mesh().nPoints());
292  PackedBoolList isZoneEdge(mesh().nEdges());
293  calcZoneMask(cellZoneI, isZonePoint, isZoneEdge);
294 
295  const dictionary& patchesDict = zoneDict.subDict("boundaryField");
296 
297  if (patchesDict.size() != 2)
298  {
300  << "Two faceZones (patches) must be specifed per cellZone. "
301  << " cellZone:" << cellZoneI
302  << " patches:" << patchesDict.toc()
303  << exit(FatalIOError);
304  }
305 
306  PtrList<scalarField> patchDist(patchesDict.size());
307  PtrList<pointVectorField> patchDisp(patchesDict.size());
308 
309  // Allocate the fields
310  label patchI = 0;
311  forAllConstIter(dictionary, patchesDict, patchIter)
312  {
313  const word& faceZoneName = patchIter().keyword();
314  label zoneI = mesh().faceZones().findZoneID(faceZoneName);
315  if (zoneI == -1)
316  {
318  << "Cannot find faceZone " << faceZoneName
319  << endl << "Valid zones are " << mesh().faceZones().names()
320  << exit(FatalIOError);
321  }
322 
323  // Determine the points of the faceZone within the cellZone
324  const faceZone& fz = mesh().faceZones()[zoneI];
325 
326  patchDist.set(patchI, new scalarField(mesh().nPoints()));
327  patchDisp.set
328  (
329  patchI,
330  new pointVectorField
331  (
332  IOobject
333  (
334  mesh().cellZones()[cellZoneI].name() + "_" + fz.name(),
335  mesh().time().timeName(),
336  mesh(),
337  IOobject::NO_READ,
338  IOobject::NO_WRITE,
339  false
340  ),
341  pointDisplacement_ // to inherit the boundary conditions
342  )
343  );
344 
345  patchI++;
346  }
347 
348 
349 
350  // 'correctBoundaryConditions'
351  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~
352  // Loops over all the faceZones and walks their boundary values
353 
354  // Make sure we can pick up bc values from field
355  pointDisplacement_.correctBoundaryConditions();
356 
357  patchI = 0;
358  forAllConstIter(dictionary, patchesDict, patchIter)
359  {
360  const word& faceZoneName = patchIter().keyword();
361  const dictionary& faceZoneDict = patchIter().dict();
362 
363  // Determine the points of the faceZone within the cellZone
364  const faceZone& fz = mesh().faceZones()[faceZoneName];
365  const labelList& fzMeshPoints = fz().meshPoints();
366  DynamicList<label> meshPoints(fzMeshPoints.size());
367  forAll(fzMeshPoints, i)
368  {
369  if (isZonePoint[fzMeshPoints[i]])
370  {
371  meshPoints.append(fzMeshPoints[i]);
372  }
373  }
374 
375  // Get initial value for all the faceZone points
376  tmp<vectorField> tseed = faceZoneEvaluate
377  (
378  fz,
379  meshPoints,
380  faceZoneDict,
381  patchDisp,
382  patchI
383  );
384 
385  if (debug)
386  {
387  Info<< "For cellZone:" << cellZoneI
388  << " for faceZone:" << fz.name()
389  << " nPoints:" << tseed().size()
390  << " have patchField:"
391  << " max:" << gMax(tseed())
392  << " min:" << gMin(tseed())
393  << " avg:" << gAverage(tseed())
394  << endl;
395  }
396 
397  // Set distance and transported value
398  walkStructured
399  (
400  cellZoneI,
401  isZonePoint,
402  isZoneEdge,
403 
404  meshPoints,
405  tseed,
406  patchDist[patchI],
407  patchDisp[patchI]
408  );
409 
410  // Implement real bc.
411  patchDisp[patchI].correctBoundaryConditions();
412 
413  patchI++;
414  }
415 
416 
417  // Solve
418  // ~~~~~
419 
420  if (debug)
421  {
422  // Normalised distance
424  (
425  IOobject
426  (
427  mesh().cellZones()[cellZoneI].name() + ":distance",
428  mesh().time().timeName(),
429  mesh(),
430  IOobject::NO_READ,
431  IOobject::NO_WRITE,
432  false
433  ),
434  pointMesh::New(mesh()),
435  dimensionedScalar("zero", dimLength, 0.0)
436  );
437 
438  forAll(distance, pointI)
439  {
440  scalar d1 = patchDist[0][pointI];
441  scalar d2 = patchDist[1][pointI];
442  if (d1 + d2 > SMALL)
443  {
444  scalar s = d1/(d1 + d2);
445  distance[pointI] = s;
446  }
447  }
448 
449  Info<< "Writing " << pointScalarField::typeName
450  << distance.name() << " to "
451  << mesh().time().timeName() << endl;
452  distance.write();
453  }
454 
455 
456  const word interpolationScheme = zoneDict.lookup("interpolationScheme");
457 
458  if (interpolationScheme == "oneSided")
459  {
460  forAll(pointDisplacement_, pointI)
461  {
462  if (isZonePoint[pointI])
463  {
464  pointDisplacement_[pointI] = patchDisp[0][pointI];
465  }
466  }
467  }
468  else if (interpolationScheme == "linear")
469  {
470  forAll(pointDisplacement_, pointI)
471  {
472  if (isZonePoint[pointI])
473  {
474  scalar d1 = patchDist[0][pointI];
475  scalar d2 = patchDist[1][pointI];
476  scalar s = d1/(d1 + d2 + VSMALL);
477 
478  const vector& pd1 = patchDisp[0][pointI];
479  const vector& pd2 = patchDisp[1][pointI];
480 
481  pointDisplacement_[pointI] = (1 - s)*pd1 + s*pd2;
482  }
483  }
484  }
485  else
486  {
488  << "Invalid interpolationScheme: " << interpolationScheme
489  << ". Valid schemes are 'oneSided' and 'linear'"
490  << exit(FatalError);
491  }
492 }
493 
494 
495 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
496 
499 (
500  const polyMesh& mesh,
501  const IOdictionary& dict
502 )
503 :
505 {}
506 
507 
510 (
511  const polyMesh& mesh,
512  const IOdictionary& dict,
513  const pointVectorField& pointDisplacement,
514  const pointIOField& points0
515 )
516 :
517  displacementMotionSolver(mesh, dict, pointDisplacement, points0, typeName)
518 {}
519 
520 
521 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
522 
525 {}
526 
527 
528 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
529 
532 {
533  tmp<pointField> tcurPoints
534  (
535  points0() + pointDisplacement_.internalField()
536  );
537 
538  return tcurPoints;
539 }
540 
541 
543 {
544  // The points have moved so before interpolation update the motionSolver
545  movePoints(mesh().points());
546 
547  // Apply boundary conditions
548  pointDisplacement_.boundaryField().updateCoeffs();
549 
550  // Solve motion on all regions (=cellZones)
551  const dictionary& regionDicts = coeffDict().subDict("regions");
552  forAllConstIter(dictionary, regionDicts, regionIter)
553  {
554  const word& cellZoneName = regionIter().keyword();
555  const dictionary& regionDict = regionIter().dict();
556 
557  label zoneI = mesh().cellZones().findZoneID(cellZoneName);
558 
559  Info<< "solving for zone: " << cellZoneName << endl;
560 
561  if (zoneI == -1)
562  {
564  << "Cannot find cellZone " << cellZoneName
565  << endl << "Valid zones are " << mesh().cellZones().names()
566  << exit(FatalIOError);
567  }
568 
569  cellZoneSolve(zoneI, regionDict);
570  }
571 
572  // Update pointDisplacement for solved values
573  const pointConstraints& pcs =
574  pointConstraints::New(pointDisplacement_.mesh());
575  pcs.constrainDisplacement(pointDisplacement_, false);
576 }
577 
578 
580 (
581  const mapPolyMesh& mpm
582 )
583 {
585 
586  const vectorField displacement(this->newPoints() - points0_);
587 
588  forAll(points0_, pointI)
589  {
590  label oldPointI = mpm.pointMap()[pointI];
591 
592  if (oldPointI >= 0)
593  {
594  label masterPointI = mpm.reversePointMap()[oldPointI];
595 
596  if ((masterPointI != pointI))
597  {
598  // newly inserted point in this cellZone
599 
600  // need to set point0 so that it represents the position that
601  // it would have had if it had existed for all time
602  points0_[pointI] -= displacement[pointI];
603  }
604  }
605  }
606 }
607 
608 
609 // ************************************************************************* //
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::PackedBoolList
A bit-packed bool list.
Definition: PackedBoolList.H:63
Foam::displacementLayeredMotionMotionSolver::curPoints
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
Definition: displacementLayeredMotionMotionSolver.C:531
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::pointEdgeStructuredWalk
Determines length of string of edges walked to point.
Definition: pointEdgeStructuredWalk.H:54
Foam::compressible::New
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Definition: turbulentFluidThermoModel.C:36
Foam::IOField< vector >
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:86
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::displacementLayeredMotionMotionSolver::displacementLayeredMotionMotionSolver
displacementLayeredMotionMotionSolver(const displacementLayeredMotionMotionSolver &)
Disallow default bitwise copy construct.
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::DynamicList< label >
Foam::gAverage
Type gAverage(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:571
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
Foam::PointEdgeWave
Wave propagation of information through grid. Every iteration information goes through one layer of e...
Definition: PointEdgeWave.H:85
mapPolyMesh.H
Foam::dictionary::lookup
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:449
Foam::polyMesh::cellZones
const cellZoneMesh & cellZones() const
Return cell zone mesh.
Definition: polyMesh.H:469
Foam::displacementMotionSolver::updateMesh
virtual void updateMesh(const mapPolyMesh &)
Update local data for topology changes.
Definition: displacementMotionSolver.C:248
interpolationTable.H
Foam::FatalIOError
IOerror FatalIOError
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::pointConstraints::constrainDisplacement
void constrainDisplacement(pointVectorField &displacement, const bool overrideValue=false) const
Apply boundary conditions (single-patch constraints),.
Definition: pointConstraints.C:369
syncTools.H
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Foam::MeshObject< pointMesh, UpdateableMeshObject, pointConstraints >::New
static const pointConstraints & New(const pointMesh &mesh)
Definition: MeshObject.C:44
Foam::PackedList::setSize
void setSize(const label, const unsigned int &val=0u)
Alias for resize()
Definition: PackedListI.H:823
Foam::cellZone
A subset of mesh cells.
Definition: cellZone.H:61
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
nPoints
label nPoints
Definition: gmvOutputHeader.H:2
displacementLayeredMotionMotionSolver.H
Foam::PtrList::set
bool set(const label) const
Is element set.
Foam::displacementLayeredMotionMotionSolver::walkStructured
void walkStructured(const label cellZoneI, const PackedBoolList &isZonePoint, const PackedBoolList &isZoneEdge, const labelList &seedPoints, const vectorField &seedData, scalarField &distance, vectorField &data) const
Definition: displacementLayeredMotionMotionSolver.C:136
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::displacementLayeredMotionMotionSolver::faceZoneEvaluate
tmp< vectorField > faceZoneEvaluate(const faceZone &fz, const labelList &meshPoints, const dictionary &dict, const PtrList< pointVectorField > &patchDisp, const label patchI) const
Definition: displacementLayeredMotionMotionSolver.C:224
Foam::faceZone
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:64
Foam::Info
messageStream Info
pointConstraints.H
PointEdgeWave.H
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
vectorField
volVectorField vectorField(fieldObject, mesh)
Foam::PtrList
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:61
Foam::orEqOp
Definition: ops.H:82
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::displacementLayeredMotionMotionSolver::solve
virtual void solve()
Solve for motion.
Definition: displacementLayeredMotionMotionSolver.C:542
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::displacementMotionSolver
Virtual base class for displacement motion solver.
Definition: displacementMotionSolver.H:55
Foam::displacementLayeredMotionMotionSolver::calcZoneMask
void calcZoneMask(const label cellZoneI, PackedBoolList &isZonePoint, PackedBoolList &isZoneEdge) const
Definition: displacementLayeredMotionMotionSolver.C:61
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
fld
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< ' ';}gmvFile<< nl;forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::zone::name
const word & name() const
Return name.
Definition: zone.H:150
Foam::ZoneMesh::findZoneID
label findZoneID(const word &zoneName) const
Find zone index given a name.
Definition: ZoneMesh.C:348
s
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().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))
Foam::distance
scalar distance(const vector &p1, const vector &p2)
Definition: curveTools.C:12
Foam::ZoneMesh::names
wordList names() const
Return a list of zone names.
Definition: ZoneMesh.C:263
Foam::DynamicList::append
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::pointConstraints
Application of (multi-)patch point contraints.
Definition: pointConstraints.H:62
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::sumOp
Definition: ops.H:162
Foam::displacementLayeredMotionMotionSolver::updateMesh
virtual void updateMesh(const mapPolyMesh &)
Update topology.
Definition: displacementLayeredMotionMotionSolver.C:580
Foam::mapPolyMesh::reversePointMap
const labelList & reversePointMap() const
Reverse point map.
Definition: mapPolyMesh.H:465
Foam::Vector< scalar >
pointEdgeStructuredWalk.H
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
scalarField
volScalarField scalarField(fieldObject, mesh)
Foam::mapPolyMesh::pointMap
const labelList & pointMap() const
Old point map.
Definition: mapPolyMesh.H:392
points
const pointField & points
Definition: gmvOutputHeader.H:1
timeName
word timeName
Definition: getTimeIndex.H:3
Foam::displacementLayeredMotionMotionSolver::cellZoneSolve
void cellZoneSolve(const label cellZoneI, const dictionary &zoneDict)
Definition: displacementLayeredMotionMotionSolver.C:286
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:330
Foam::displacementLayeredMotionMotionSolver::~displacementLayeredMotionMotionSolver
~displacementLayeredMotionMotionSolver()
Destructor.
Definition: displacementLayeredMotionMotionSolver.C:524
Foam::interpolationTable
An interpolation/look-up table of scalar vs <Type> values. The reference scalar values must be monoto...
Definition: interpolationTable.H:77
Foam::dictionary::toc
wordList toc() const
Return the table of contents.
Definition: dictionary.C:697
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::gMin
Type gMin(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:563
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::dictionary::subDict
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:631
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::data
Database for solution data, solver performance and other reduced data.
Definition: data.H:52
pointFields.H
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:562
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47