53Foam::velocityLaplacianFvMotionSolver::velocityLaplacianFvMotionSolver
73 cellMotionBoundaryTypes<
vector>(pointMotionU_.boundaryField())
77 coeffDict().
found(
"interpolation")
99 interpolationPtr_->interpolate
108 + fvMesh_.time().deltaTValue()*pointMotionU_.primitiveField()
121 movePoints(fvMesh_.points());
123 diffusivityPtr_->correct();
124 pointMotionU_.boundaryFieldRef().updateCoeffs();
130 getOrDefault<label>(
"nNonOrthogonalCorrectors", 1)
140 * diffusivityPtr_->operator()(),
142 "laplacian(diffusivity,cellMotionU)"
149 UEqn.solveSegregatedOrCoupled();
157 const mapPolyMesh& mpm
164 diffusivityPtr_.reset(
nullptr);
168 coeffDict().lookup(
"diffusivity")
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
const Time & time() const noexcept
Return Time associated with the objectRegistry.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
Base class for fvMesh based motionSolvers.
const fvMesh & fvMesh_
The fvMesh to be moved.
wordList cellMotionBoundaryTypes(const typename GeometricField< Type, pointPatchField, pointMesh >::Boundary &pmUbf) const
Create the corresponding patch types for cellMotion from those.
fvMotionSolver(const polyMesh &)
Construct from polyMesh.
const fvMesh & mesh() const
Return reference to the fvMesh to be moved.
Finite-volume options, which is an IOdictionary of values and a fv::optionList.
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present otherwise lookup and return.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Abstract base class for cell-centre mesh motion diffusivity.
static autoPtr< motionDiffusivity > New(const fvMesh &mesh, Istream &mdData)
Select null constructed.
Base class for interpolation of cell displacement fields, generated by fvMotionSolvers,...
Virtual base class for mesh motion solver.
static autoPtr< motionSolver > New(const polyMesh &)
Select constructed from polyMesh.
virtual void twoDCorrectPoints(pointField &) const
const dictionary & coeffDict() const
Const access to the coefficients dictionary.
Mesh consisting of general polyhedral cells.
Lookup type of boundary radiation properties.
A class for managing temporary objects.
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Mesh motion solver for an fvMesh. Based on solving the cell-centre Laplacian for the motion velocity.
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
virtual void updateMesh(const mapPolyMesh &)
Update topology.
~velocityLaplacianFvMotionSolver()
Destructor.
virtual void solve()
Solve for motion.
Virtual base class for velocity motion solver.
velocityMotionSolver(const velocityMotionSolver &)=delete
No copy construct.
pointVectorField pointMotionU_
Point motion field.
virtual void movePoints(const pointField &)
Update local data for geometry changes.
virtual void updateMesh(const mapPolyMesh &)
Update local data for topology changes.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Calculate the matrix for the laplacian of the field.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
const dimensionSet dimViscosity
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
fvMatrix< vector > fvVectorMatrix
static constexpr const zero Zero
Global zero (0).
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.