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-2014 OpenFOAM Foundation
6  \\/ M anipulation |
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 }
49 
50 
51 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
52 
53 void Foam::displacementLayeredMotionMotionSolver::calcZoneMask
54 (
55  const label cellZoneI,
56  PackedBoolList& isZonePoint,
57  PackedBoolList& isZoneEdge
58 ) const
59 {
60  if (cellZoneI == -1)
61  {
62  isZonePoint.setSize(mesh().nPoints());
63  isZonePoint = 1;
64 
65  isZoneEdge.setSize(mesh().nEdges());
66  isZoneEdge = 1;
67  }
68  else
69  {
70  const cellZone& cz = mesh().cellZones()[cellZoneI];
71 
72  label nPoints = 0;
73  forAll(cz, i)
74  {
75  const labelList& cPoints = mesh().cellPoints(cz[i]);
76  forAll(cPoints, cPointI)
77  {
78  if (!isZonePoint[cPoints[cPointI]])
79  {
80  isZonePoint[cPoints[cPointI]] = 1;
81  nPoints++;
82  }
83  }
84  }
86  (
87  mesh(),
88  isZonePoint,
89  orEqOp<unsigned int>(),
90  0
91  );
92 
93 
94  // Mark edge inside cellZone
95  label nEdges = 0;
96  forAll(cz, i)
97  {
98  const labelList& cEdges = mesh().cellEdges(cz[i]);
99  forAll(cEdges, cEdgeI)
100  {
101  if (!isZoneEdge[cEdges[cEdgeI]])
102  {
103  isZoneEdge[cEdges[cEdgeI]] = 1;
104  nEdges++;
105  }
106  }
107  }
109  (
110  mesh(),
111  isZoneEdge,
112  orEqOp<unsigned int>(),
113  0
114  );
115 
116  if (debug)
117  {
118  Info<< "On cellZone " << cz.name()
119  << " marked " << returnReduce(nPoints, sumOp<label>())
120  << " points and " << returnReduce(nEdges, sumOp<label>())
121  << " edges." << endl;
122  }
123  }
124 }
125 
126 
127 // Find distance to starting point
128 void Foam::displacementLayeredMotionMotionSolver::walkStructured
129 (
130  const label cellZoneI,
131  const PackedBoolList& isZonePoint,
132  const PackedBoolList& isZoneEdge,
133  const labelList& seedPoints,
134  const vectorField& seedData,
135  scalarField& distance,
136  vectorField& data
137 ) const
138 {
139  List<pointEdgeStructuredWalk> seedInfo(seedPoints.size());
140 
141  forAll(seedPoints, i)
142  {
143  seedInfo[i] = pointEdgeStructuredWalk
144  (
145  points0()[seedPoints[i]], // location of data
146  points0()[seedPoints[i]], // previous location
147  0.0,
148  seedData[i]
149  );
150  }
151 
152  // Current info on points
153  List<pointEdgeStructuredWalk> allPointInfo(mesh().nPoints());
154 
155  // Mark points inside cellZone.
156  // Note that we use points0, not mesh.points()
157  // so as not to accumulate errors.
158  forAll(isZonePoint, pointI)
159  {
160  if (isZonePoint[pointI])
161  {
162  allPointInfo[pointI] = pointEdgeStructuredWalk
163  (
164  points0()[pointI], // location of data
165  vector::max, // not valid
166  0.0,
167  vector::zero // passive data
168  );
169  }
170  }
171 
172  // Current info on edges
173  List<pointEdgeStructuredWalk> allEdgeInfo(mesh().nEdges());
174 
175  // Mark edges inside cellZone
176  forAll(isZoneEdge, edgeI)
177  {
178  if (isZoneEdge[edgeI])
179  {
180  allEdgeInfo[edgeI] = pointEdgeStructuredWalk
181  (
182  mesh().edges()[edgeI].centre(points0()), // location of data
183  vector::max, // not valid
184  0.0,
186  );
187  }
188  }
189 
190  // Walk
191  PointEdgeWave<pointEdgeStructuredWalk> wallCalc
192  (
193  mesh(),
194  seedPoints,
195  seedInfo,
196 
197  allPointInfo,
198  allEdgeInfo,
199  mesh().globalData().nTotalPoints() // max iterations
200  );
201 
202  // Extract distance and passive data
203  forAll(allPointInfo, pointI)
204  {
205  if (isZonePoint[pointI])
206  {
207  distance[pointI] = allPointInfo[pointI].dist();
208  data[pointI] = allPointInfo[pointI].data();
209  }
210  }
211 }
212 
213 
214 // Evaluate faceZone patch
216 Foam::displacementLayeredMotionMotionSolver::faceZoneEvaluate
217 (
218  const faceZone& fz,
219  const labelList& meshPoints,
220  const dictionary& dict,
221  const PtrList<pointVectorField>& patchDisp,
222  const label patchI
223 ) const
224 {
225  tmp<vectorField> tfld(new vectorField(meshPoints.size()));
226  vectorField& fld = tfld();
227 
228  const word type(dict.lookup("type"));
229 
230  if (type == "fixedValue")
231  {
232  fld = vectorField("value", dict, meshPoints.size());
233  }
234  else if (type == "timeVaryingUniformFixedValue")
235  {
236  interpolationTable<vector> timeSeries(dict);
237 
238  fld = timeSeries(mesh().time().timeOutputValue());
239  }
240  else if (type == "slip")
241  {
242  if ((patchI % 2) != 1)
243  {
245  (
246  "displacementLayeredMotionMotionSolver::faceZoneEvaluate"
247  "("
248  "const faceZone&, "
249  "const labelList&, "
250  "const dictionary&, "
251  "const PtrList<pointVectorField>&, "
252  "const label"
253  ") const",
254  *this
255  ) << "slip can only be used on second faceZone patch of pair. "
256  << "FaceZone:" << fz.name()
257  << exit(FatalIOError);
258  }
259  // Use field set by previous bc
260  fld = vectorField(patchDisp[patchI - 1], meshPoints);
261  }
262  else if (type == "follow")
263  {
264  // Only on boundary faces - follow boundary conditions
265  fld = vectorField(pointDisplacement_, meshPoints);
266  }
267  else if (type == "uniformFollow")
268  {
269  // Reads name of name of patch. Then get average point dislacement on
270  // patch. That becomes the value of fld.
271  const word patchName(dict.lookup("patch"));
272  label patchID = mesh().boundaryMesh().findPatchID(patchName);
273  pointField pdf
274  (
275  pointDisplacement_.boundaryField()[patchID].patchInternalField()
276  );
277  fld = gAverage(pdf);
278  }
279  else
280  {
282  (
283  "displacementLayeredMotionMotionSolver::faceZoneEvaluate"
284  "("
285  "const faceZone&, "
286  "const labelList&, "
287  "const dictionary&, "
288  "const PtrList<pointVectorField>&, "
289  "const label"
290  ") const",
291  *this
292  ) << "Unknown faceZonePatch type " << type << " for faceZone "
293  << fz.name() << exit(FatalIOError);
294  }
295  return tfld;
296 }
297 
298 
299 void Foam::displacementLayeredMotionMotionSolver::cellZoneSolve
300 (
301  const label cellZoneI,
302  const dictionary& zoneDict
303 )
304 {
305  PackedBoolList isZonePoint(mesh().nPoints());
306  PackedBoolList isZoneEdge(mesh().nEdges());
307  calcZoneMask(cellZoneI, isZonePoint, isZoneEdge);
308 
309  const dictionary& patchesDict = zoneDict.subDict("boundaryField");
310 
311  if (patchesDict.size() != 2)
312  {
314  (
315  "displacementLayeredMotionMotionSolver::"
316  "cellZoneSolve(const label, const dictionary&)",
317  *this
318  ) << "Two faceZones (patches) must be specifed per cellZone. "
319  << " cellZone:" << cellZoneI
320  << " patches:" << patchesDict.toc()
321  << exit(FatalIOError);
322  }
323 
324  PtrList<scalarField> patchDist(patchesDict.size());
325  PtrList<pointVectorField> patchDisp(patchesDict.size());
326 
327  // Allocate the fields
328  label patchI = 0;
329  forAllConstIter(dictionary, patchesDict, patchIter)
330  {
331  const word& faceZoneName = patchIter().keyword();
332  label zoneI = mesh().faceZones().findZoneID(faceZoneName);
333  if (zoneI == -1)
334  {
336  (
337  "displacementLayeredMotionMotionSolver::"
338  "cellZoneSolve(const label, const dictionary&)",
339  *this
340  ) << "Cannot find faceZone " << faceZoneName
341  << endl << "Valid zones are " << mesh().faceZones().names()
342  << exit(FatalIOError);
343  }
344 
345  // Determine the points of the faceZone within the cellZone
346  const faceZone& fz = mesh().faceZones()[zoneI];
347 
348  patchDist.set(patchI, new scalarField(mesh().nPoints()));
349  patchDisp.set
350  (
351  patchI,
352  new pointVectorField
353  (
354  IOobject
355  (
356  mesh().cellZones()[cellZoneI].name() + "_" + fz.name(),
357  mesh().time().timeName(),
358  mesh(),
361  false
362  ),
363  pointDisplacement_ // to inherit the boundary conditions
364  )
365  );
366 
367  patchI++;
368  }
369 
370 
371 
372  // 'correctBoundaryConditions'
373  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~
374  // Loops over all the faceZones and walks their boundary values
375 
376  // Make sure we can pick up bc values from field
377  pointDisplacement_.correctBoundaryConditions();
378 
379  patchI = 0;
380  forAllConstIter(dictionary, patchesDict, patchIter)
381  {
382  const word& faceZoneName = patchIter().keyword();
383  const dictionary& faceZoneDict = patchIter().dict();
384 
385  // Determine the points of the faceZone within the cellZone
386  const faceZone& fz = mesh().faceZones()[faceZoneName];
387  const labelList& fzMeshPoints = fz().meshPoints();
388  DynamicList<label> meshPoints(fzMeshPoints.size());
389  forAll(fzMeshPoints, i)
390  {
391  if (isZonePoint[fzMeshPoints[i]])
392  {
393  meshPoints.append(fzMeshPoints[i]);
394  }
395  }
396 
397  // Get initial value for all the faceZone points
398  tmp<vectorField> tseed = faceZoneEvaluate
399  (
400  fz,
401  meshPoints,
402  faceZoneDict,
403  patchDisp,
404  patchI
405  );
406 
407  if (debug)
408  {
409  Info<< "For cellZone:" << cellZoneI
410  << " for faceZone:" << fz.name()
411  << " nPoints:" << tseed().size()
412  << " have patchField:"
413  << " max:" << gMax(tseed())
414  << " min:" << gMin(tseed())
415  << " avg:" << gAverage(tseed())
416  << endl;
417  }
418 
419  // Set distance and transported value
420  walkStructured
421  (
422  cellZoneI,
423  isZonePoint,
424  isZoneEdge,
425 
426  meshPoints,
427  tseed,
428  patchDist[patchI],
429  patchDisp[patchI]
430  );
431 
432  // Implement real bc.
433  patchDisp[patchI].correctBoundaryConditions();
434 
435  patchI++;
436  }
437 
438 
439  // Solve
440  // ~~~~~
441 
442  if (debug)
443  {
444  // Normalised distance
446  (
447  IOobject
448  (
449  mesh().cellZones()[cellZoneI].name() + ":distance",
450  mesh().time().timeName(),
451  mesh(),
454  false
455  ),
456  pointMesh::New(mesh()),
457  dimensionedScalar("zero", dimLength, 0.0)
458  );
459 
460  forAll(distance, pointI)
461  {
462  scalar d1 = patchDist[0][pointI];
463  scalar d2 = patchDist[1][pointI];
464  if (d1 + d2 > SMALL)
465  {
466  scalar s = d1/(d1 + d2);
467  distance[pointI] = s;
468  }
469  }
470 
471  Info<< "Writing " << pointScalarField::typeName
472  << distance.name() << " to "
473  << mesh().time().timeName() << endl;
474  distance.write();
475  }
476 
477 
478  const word interpolationScheme = zoneDict.lookup("interpolationScheme");
479 
480  if (interpolationScheme == "oneSided")
481  {
482  forAll(pointDisplacement_, pointI)
483  {
484  if (isZonePoint[pointI])
485  {
486  pointDisplacement_[pointI] = patchDisp[0][pointI];
487  }
488  }
489  }
490  else if (interpolationScheme == "linear")
491  {
492  forAll(pointDisplacement_, pointI)
493  {
494  if (isZonePoint[pointI])
495  {
496  scalar d1 = patchDist[0][pointI];
497  scalar d2 = patchDist[1][pointI];
498  scalar s = d1/(d1 + d2 + VSMALL);
499 
500  const vector& pd1 = patchDisp[0][pointI];
501  const vector& pd2 = patchDisp[1][pointI];
502 
503  pointDisplacement_[pointI] = (1 - s)*pd1 + s*pd2;
504  }
505  }
506  }
507  else
508  {
510  (
511  "displacementLayeredMotionMotionSolver::"
512  "cellZoneSolve(const label, const dictionary&)"
513  )
514  << "Invalid interpolationScheme: " << interpolationScheme
515  << ". Valid schemes are 'oneSided' and 'linear'"
516  << exit(FatalError);
517  }
518 }
519 
520 
521 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
522 
523 Foam::displacementLayeredMotionMotionSolver::
524 displacementLayeredMotionMotionSolver
525 (
526  const polyMesh& mesh,
527  const IOdictionary& dict
528 )
529 :
530  displacementMotionSolver(mesh, dict, typeName)
531 {}
532 
533 
534 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
535 
538 {}
539 
540 
541 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
542 
545 {
546  tmp<pointField> tcurPoints
547  (
548  points0() + pointDisplacement_.internalField()
549  );
550 
551  return tcurPoints;
552 }
553 
554 
556 {
557  // The points have moved so before interpolation update the motionSolver
558  movePoints(mesh().points());
559 
560  // Apply boundary conditions
561  pointDisplacement_.boundaryField().updateCoeffs();
562 
563  // Solve motion on all regions (=cellZones)
564  const dictionary& regionDicts = coeffDict().subDict("regions");
565  forAllConstIter(dictionary, regionDicts, regionIter)
566  {
567  const word& cellZoneName = regionIter().keyword();
568  const dictionary& regionDict = regionIter().dict();
569 
570  label zoneI = mesh().cellZones().findZoneID(cellZoneName);
571 
572  Info<< "solving for zone: " << cellZoneName << endl;
573 
574  if (zoneI == -1)
575  {
577  (
578  "displacementLayeredMotionMotionSolver::solve()",
579  *this
580  ) << "Cannot find cellZone " << cellZoneName
581  << endl << "Valid zones are " << mesh().cellZones().names()
582  << exit(FatalIOError);
583  }
584 
585  cellZoneSolve(zoneI, regionDict);
586  }
587 
588  // Update pointDisplacement for solved values
589  const pointConstraints& pcs =
590  pointConstraints::New(pointDisplacement_.mesh());
591  pcs.constrainDisplacement(pointDisplacement_, false);
592 }
593 
594 
596 (
597  const mapPolyMesh& mpm
598 )
599 {
601 
602  const vectorField displacement(this->newPoints() - points0_);
603 
604  forAll(points0_, pointI)
605  {
606  label oldPointI = mpm.pointMap()[pointI];
607 
608  if (oldPointI >= 0)
609  {
610  label masterPointI = mpm.reversePointMap()[oldPointI];
611 
612  if ((masterPointI != pointI))
613  {
614  // newly inserted point in this cellZone
615 
616  // need to set point0 so that it represents the position that
617  // it would have had if it had existed for all time
618  points0_[pointI] -= displacement[pointI];
619  }
620  }
621  }
622 }
623 
624 
625 // ************************************************************************* //
const cellZoneMesh & cellZones() const
Return cell zone mesh.
Definition: polyMesh.H:469
const labelListList & cellPoints() const
const pointField & points
const faceZoneMesh & faceZones() const
Return face zone mesh.
Definition: polyMesh.H:463
virtual void updateMesh(const mapPolyMesh &)
Update local data for topology changes.
static const char *const typeName
Definition: Field.H:94
bool set(const label) const
Is element set.
Definition: PtrListI.H:107
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
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 ))
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
const labelList & reversePointMap() const
Reverse point map.
Definition: mapPolyMesh.H:465
static const pointMesh & New(const polyMesh &mesh)
virtual void updateMesh(const mapPolyMesh &)
Update topology.
A class for handling words, derived from string.
Definition: word.H:59
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
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:741
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const Vector max
Definition: Vector.H:82
static void syncEdgeList(const polyMesh &, List< T > &, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh edges.
messageStream Info
dynamicFvMesh & mesh
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:638
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
const labelListList & cellEdges() const
Namespace for OpenFOAM.
scalar distance(const vector &p1, const vector &p2)
Definition: curveTools.C:12
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
IOerror FatalIOError
void constrainDisplacement(pointVectorField &displacement, const bool overrideValue=false) const
Apply boundary conditions (single-patch constraints),.
Field< vector > vectorField
Specialisation of Field<T> for vector.
#define forAll(list, i)
Definition: UList.H:421
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
virtual Ostream & write(const token &)=0
Write next token to stream.
Macros for easy insertion into run-time selection tables.
Type gMin(const FieldField< Field, Type > &f)
static void syncPointList(const polyMesh &, List< T > &, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh points.
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:589
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:314
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
error FatalError
Type gAverage(const FieldField< Field, Type > &f)
List< label > labelList
A List of labels.
Definition: labelList.H:56
static const Vector zero
Definition: Vector.H:80
const labelList & pointMap() const
Old point map.
Definition: mapPolyMesh.H:392
GeometricField< vector, pointPatchField, pointMesh > pointVectorField
Definition: pointFields.H:50
label findZoneID(const word &zoneName) const
Find zone index given a name.
Definition: ZoneMesh.C:354
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
wordList names() const
Return a list of zone names.
Definition: ZoneMesh.C:269
GeometricField< scalar, pointPatchField, pointMesh > pointScalarField
Definition: pointFields.H:49
label nPoints
label findPatchID(const word &patchName) const
Find patch index given a name.
#define FatalIOErrorIn(functionName, ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:325
Application of (multi-)patch point contraints.
Type gMax(const FieldField< Field, Type > &f)
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A class for managing temporary objects.
Definition: PtrList.H:118
defineTypeNameAndDebug(combustionModel, 0)
word timeName
Definition: getTimeIndex.H:3
Virtual base class for displacement motion solver.