Go to the documentation of this file.
43 velocityDisplacementMotionSolver,
52 Foam::velocityDisplacementMotionSolver::pointDisplacementBoundaryTypes()
const
54 const pointVectorField::Boundary& pmUbf(
pointMotionU().boundaryField());
60 if (
isA<fixedValuePointPatchField<vector>>(pmUbf[patchI]))
62 cmUbf[patchI] = fixedValuePointPatchField<vector>::typeName;
72 Foam::velocityDisplacementMotionSolver::velocityDisplacementMotionSolver
79 displacementMotionSolverPtr_()
87 "pointVelocityDisplacement",
91 pointMotionU().mesh(),
93 pointDisplacementBoundaryTypes()
98 displacementMotionSolverPtr_.reset
104 coeffDict().get<word>(
"solver"),
135 return displacementMotionSolverPtr_->curPoints();
143 const scalar deltaT(
mesh().time().deltaTValue());
148 displacementMotionSolverPtr_->pointDisplacement()
152 mesh().
points() - displacementMotionSolverPtr_->points0()
156 pointMotionU().correctBoundaryConditions();
158 pointVectorField::Boundary& dispBf = displacement.boundaryFieldRef();
161 forAll(pointMotionU().boundaryField(), patchI)
165 pointMotionU().boundaryField()[patchI]
174 displacementMotionSolverPtr_->solve();
177 pointMotionU().primitiveFieldRef() =
178 (displacement.primitiveField() - displacementOld)/deltaT;
186 displacementMotionSolverPtr_->movePoints(
p);
197 displacementMotionSolverPtr_->updateMesh(mpm);
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
virtual const pointField & points() const
Return raw points.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
virtual void movePoints(const pointField &)
Update geometry.
static autoPtr< displacementMotionSolver > New(const word &solverTypeName, const polyMesh &, const IOdictionary &, const pointVectorField &pointDisplacement, const pointIOField &points0)
Select constructed from polyMesh, dictionary and components.
A primitive field of type <T> with automated input and output.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
A class for managing temporary objects.
virtual void movePoints(const pointField &)
Update local data for geometry changes.
static word timeName(const scalar t, const int precision=precision_)
Virtual base class for velocity motion solver.
~velocityDisplacementMotionSolver()
Destructor.
Mesh consisting of general polyhedral cells.
#define forAll(list, i)
Loop across all elements in list.
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
static IOobject points0IO(const polyMesh &mesh)
Return IO object for points0.
List< word > wordList
A List of words.
const fileName & name() const
The dictionary name.
virtual void updateMesh(const mapPolyMesh &)
Update local data for topology changes.
Virtual base class for displacement motion solver.
Macros for easy insertion into run-time selection tables.
virtual void solve()
Solve for motion.
virtual void updateMesh(const mapPolyMesh &)
Update topology.
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
pointField points0(pointIOField(IOobject("points", mesh.time().constant(), polyMesh::meshSubDir, mesh, IOobject::MUST_READ, IOobject::NO_WRITE, false)))
const TargetType * isA(const Type &t)
Check if dynamic_cast to TargetType is possible.
pointVectorField & pointMotionU()
Return reference to the point motion velocity field.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
tmp< Field< Type > > patchInternalField() const
Return field created from appropriate internal field values.
const Time & time() const
Return the top-level database.
const word & constant() const
Return constant name.
defineTypeNameAndDebug(combustionModel, 0)