53 void Foam::displacementLayeredMotionMotionSolver::calcZoneMask
55 const label cellZoneI,
56 PackedBoolList& isZonePoint,
57 PackedBoolList& isZoneEdge
65 isZoneEdge.setSize(
mesh().nEdges());
78 if (!isZonePoint[cPoints[cPointi]])
80 isZonePoint[cPoints[cPointi]] = 1;
89 orEqOp<unsigned int>(),
101 if (!isZoneEdge[cEdges[cEdgeI]])
103 isZoneEdge[cEdges[cEdgeI]] = 1;
112 orEqOp<unsigned int>(),
118 Info<<
"On cellZone " << cz.name()
120 <<
" points and " <<
returnReduce(nEdges, sumOp<label>())
121 <<
" edges." <<
endl;
128 void Foam::displacementLayeredMotionMotionSolver::walkStructured
130 const label cellZoneI,
131 const PackedBoolList& isZonePoint,
132 const PackedBoolList& isZoneEdge,
139 List<pointEdgeStructuredWalk> seedInfo(seedPoints.size());
143 seedInfo[i] = pointEdgeStructuredWalk
145 points0()[seedPoints[i]],
146 points0()[seedPoints[i]],
153 List<pointEdgeStructuredWalk> allPointInfo(mesh().
nPoints());
158 forAll(isZonePoint, pointi)
160 if (isZonePoint[pointi])
162 allPointInfo[pointi] = pointEdgeStructuredWalk
173 List<pointEdgeStructuredWalk> allEdgeInfo(mesh().nEdges());
178 if (isZoneEdge[edgeI])
180 allEdgeInfo[edgeI] = pointEdgeStructuredWalk
182 mesh().edges()[edgeI].centre(points0()),
191 PointEdgeWave<pointEdgeStructuredWalk> wallCalc
199 mesh().globalData().nTotalPoints()
203 forAll(allPointInfo, pointi)
205 if (isZonePoint[pointi])
207 distance[pointi] = allPointInfo[pointi].dist();
208 data[pointi] = allPointInfo[pointi].data();
216 Foam::displacementLayeredMotionMotionSolver::faceZoneEvaluate
220 const dictionary&
dict,
221 const PtrList<pointVectorField>& patchDisp,
225 tmp<vectorField> tfld(
new vectorField(meshPoints.size()));
230 if (
type ==
"fixedValue")
234 else if (
type ==
"timeVaryingUniformFixedValue")
237 Function1s::Table<vector>
240 mesh().time().userUnits(),
243 ).value(mesh().time().value());
245 else if (
type ==
"slip")
250 <<
"FaceZone:" << fz.name()
257 else if (
type ==
"follow")
262 else if (
type ==
"uniformFollow")
267 label patchID = mesh().boundaryMesh().findIndex(patchName);
270 pointDisplacement_.boundaryField()[patchID].patchInternalField()
277 <<
"Unknown faceZonePatch type " <<
type <<
" for faceZone "
284 void Foam::displacementLayeredMotionMotionSolver::cellZoneSolve
286 const label cellZoneI,
287 const dictionary& zoneDict
290 PackedBoolList isZonePoint(mesh().
nPoints());
291 PackedBoolList isZoneEdge(mesh().nEdges());
292 calcZoneMask(cellZoneI, isZonePoint, isZoneEdge);
294 const dictionary& patchesDict = zoneDict.subDict(
"boundaryField");
296 if (patchesDict.size() != 2)
299 <<
"Two faceZones (patches) must be specified per cellZone. "
300 <<
" cellZone:" << cellZoneI
301 <<
" patches:" << patchesDict.toc()
305 PtrList<scalarField> patchDist(patchesDict.size());
306 PtrList<pointVectorField> patchDisp(patchesDict.size());
312 const word& faceZoneName = patchiter().keyword();
313 label zoneI = mesh().faceZones().findIndex(faceZoneName);
317 <<
"Cannot find faceZone " << faceZoneName
318 <<
endl <<
"Valid zones are " << mesh().faceZones().toc()
323 const faceZone& fz = mesh().faceZones()[zoneI];
333 mesh().cellZones()[cellZoneI].
name() +
"_" + fz.name(),
334 mesh().time().
name(),
354 pointDisplacement_.correctBoundaryConditions();
359 const word& faceZoneName = patchiter().keyword();
360 const dictionary& faceZoneDict = patchiter().dict();
363 const faceZone& fz = mesh().faceZones()[faceZoneName];
364 const labelList& fzMeshPoints = fz.patch().meshPoints();
365 DynamicList<label> meshPoints(fzMeshPoints.size());
368 if (isZonePoint[fzMeshPoints[i]])
370 meshPoints.append(fzMeshPoints[i]);
375 tmp<vectorField> tseed = faceZoneEvaluate
386 Info<<
"For cellZone:" << cellZoneI
387 <<
" for faceZone:" << fz.name()
388 <<
" nPoints:" << tseed().size()
389 <<
" have patchField:"
390 <<
" max:" <<
gMax(tseed())
391 <<
" min:" <<
gMin(tseed())
410 patchDisp[
patchi].correctBoundaryConditions();
426 mesh().cellZones()[cellZoneI].
name() +
":distance",
427 mesh().time().
name(),
439 scalar d1 = patchDist[0][pointi];
440 scalar d2 = patchDist[1][pointi];
443 scalar
s = d1/(d1 + d2);
444 distance[pointi] =
s;
449 << distance.name() <<
" to "
450 << mesh().time().name() <<
endl;
455 const word interpolationScheme = zoneDict.lookup(
"interpolationScheme");
457 if (interpolationScheme ==
"oneSided")
459 forAll(pointDisplacement_, pointi)
461 if (isZonePoint[pointi])
463 pointDisplacement_[pointi] = patchDisp[0][pointi];
467 else if (interpolationScheme ==
"linear")
469 forAll(pointDisplacement_, pointi)
471 if (isZonePoint[pointi])
473 scalar d1 = patchDist[0][pointi];
474 scalar d2 = patchDist[1][pointi];
475 scalar
s = d1/(d1 + d2 + vSmall);
477 const vector& pd1 = patchDisp[0][pointi];
478 const vector& pd2 = patchDisp[1][pointi];
480 pointDisplacement_[pointi] = (1 -
s)*pd1 +
s*pd2;
487 <<
"Invalid interpolationScheme: " << interpolationScheme
488 <<
". Valid schemes are 'oneSided' and 'linear'"
522 points0() + pointDisplacement_.primitiveField()
532 movePoints(mesh().
points());
535 pointDisplacement_.boundaryFieldRef().updateCoeffs();
541 const word& cellZoneName = regionIter().keyword();
542 const dictionary& regionDict = regionIter().dict();
544 label zoneI = mesh().cellZones().findIndex(cellZoneName);
546 Info<<
"solving for zone: " << cellZoneName <<
endl;
551 <<
"Cannot find cellZone " << cellZoneName
552 <<
endl <<
"Valid zones are " << mesh().cellZones().toc()
556 cellZoneSolve(zoneI, regionDict);
572 <<
"Probably inconsistent with points0MotionSolver" <<
nl
573 <<
" Needs to be updated and tested."
578 const vectorField displacement(this->newPoints() - points0_);
588 if ((masterPointi != pointi))
594 points0_[pointi] -= displacement[pointi];
#define forAll(list, i)
Loop across all elements in list.
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Macros for easy insertion into run-time selection tables.
static pointMesh & New(const word &name, const polyMesh &mesh)
Construct and return the named DemandDrivenMeshObject.
static const char *const typeName
virtual Ostream & write(const char)=0
Write character.
A list of keyword definitions, which are a keyword followed by any number of values (e....
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Mesh motion solver for an (multi-block) extruded fvMesh. Gets given the structure of the mesh blocks ...
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
displacementLayeredMotionMotionSolver(const word &name, const polyMesh &, const dictionary &)
Construct from polyMesh and dictionary.
virtual void topoChange(const polyTopoChangeMap &)
Update topology.
~displacementLayeredMotionMotionSolver()
Destructor.
virtual void solve()
Solve for motion.
Virtual base class for displacement motion solver.
Virtual base class for mesh motion solver.
const polyMesh & mesh() const
Return reference to mesh.
Application of (multi-)patch point constraints.
void constrainDisplacement(pointVectorField &displacement, const bool overrideValue=false) const
Apply boundary conditions (single-patch constraints),.
virtual void topoChange(const polyTopoChangeMap &)
Update local data for topology changes.
Mesh consisting of general polyhedral cells.
const cellZoneList & cellZones() const
Return cell zones.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
const labelList & reversePointMap() const
Reverse point map.
const labelList & pointMap() const
Old point map.
const labelListList & cellEdges() const
const labelListList & cellPoints() const
A class for managing temporary objects.
A class for handling words, derived from string.
static const word null
An empty word.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
gmvFile<< "tracers "<< particles.size()<< nl;{ pointField positions(particles.size());label particlei=0;forAllConstIter(Cloud< passiveParticle >, particles, iter) { positions[particlei++]=iter().position(mesh);} for(i=0;i< pTraits< point >::nComponents;i++) { forAll(positions, particlei) { gmvFile<< component(positions[particlei], i)<< ' ';} gmvFile<< nl;}}forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
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.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
errorManipArg< error, int > exit(error &err, const int errNo=1)
List< label > labelList
A List of labels.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
addToRunTimeSelectionTable(polyPatch, mergedCyclicPolyPatch, word)
Ostream & endl(Ostream &os)
Add newline and flush stream.
word name(const bool)
Return a word representation of a bool.
const dimensionSet dimLength
vectorField pointField
pointField is a vectorField.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
PointField< scalar > pointScalarField
Vector< scalar > vector
A scalar version of the templated Vector.
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
defineTypeNameAndDebug(combustionModel, 0)
PointField< vector > pointVectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Type gAverage(const FieldField< Field, Type > &f)
Type gMin(const FieldField< Field, Type > &f)
Type gMax(const FieldField< Field, Type > &f)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.