50 void Foam::displacementLayeredMotionSolver::walkLayers
52 const polyPatch& startPatch,
58 const labelList& startMeshPoints = startPatch.meshPoints();
68 List<pointEdgeStructuredWalk> startInfo(startMeshPoints.size());
71 startInfo[i] = pointEdgeStructuredWalk
81 List<pointEdgeStructuredWalk> allPointInfo(
mesh().
nPoints());
82 forAll(allPointInfo, pointi)
84 allPointInfo[pointi] = pointEdgeStructuredWalk
94 List<pointEdgeStructuredWalk> allEdgeInfo(
mesh().nEdges());
97 allEdgeInfo[edgei] = pointEdgeStructuredWalk
107 PointEdgeWave<pointEdgeStructuredWalk> walk
114 mesh().globalData().nTotalPoints()
118 forAll(allPointInfo, pointi)
120 distance[pointi] = allPointInfo[pointi].dist();
121 displacement[pointi] = allPointInfo[pointi].data();
136 oppositePatchNames_(
dict.lookup(
"oppositePatches")),
158 points0() + pointDisplacement_.primitiveField()
171 pointDisplacement_.boundaryFieldRef().updateCoeffs();
177 walkLayers(patch0, patchDist0, patchDisp0);
183 walkLayers(patch1, patchDist1, patchDisp1);
187 const scalarField w(patchDist0/(patchDist0 + patchDist1 + small));
190 pointDisplacement_.primitiveFieldRef() = (1 - w)*patchDisp0 + w*patchDisp1;
#define forAll(list, i)
Loop across all elements in list.
Macros for easy insertion into run-time selection tables.
static pointConstraints & New(const word &name, const pointMesh &mesh)
Construct and return the named DemandDrivenMeshObject.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Interpolating motion solver for extruded/layered meshes.
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
virtual void topoChange(const polyTopoChangeMap &)
Update topology (not implemented)
~displacementLayeredMotionSolver()
Destructor.
displacementLayeredMotionSolver(const word &name, const polyMesh &, const dictionary &)
Construct from polyMesh and dictionary.
virtual void solve()
Solve for motion.
Virtual base class for displacement motion solver.
pointVectorField pointDisplacement_
Point motion field.
Virtual base class for mesh motion solver.
const polyMesh & mesh() const
Return reference to mesh.
void constrainDisplacement(pointVectorField &displacement, const bool overrideValue=false) const
Apply boundary conditions (single-patch constraints),.
pointField & points0()
Return reference to the reference field.
Mesh consisting of general polyhedral cells.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
A patch is a list of labels that address the faces in the global face list.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
A class for managing temporary objects.
A class for handling words, derived from string.
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
volScalarField scalarField(fieldObject, mesh)
volVectorField vectorField(fieldObject, mesh)
List< label > labelList
A List of labels.
addToRunTimeSelectionTable(polyPatch, mergedCyclicPolyPatch, word)
labelList second(const UList< labelPair > &p)
labelList first(const UList< labelPair > &p)
defineTypeNameAndDebug(combustionModel, 0)
Field< vector > vectorField
Specialisation of Field<T> for vector.
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.