42 displacementComponentLaplacianFvMotionSolver,
64 "cellDisplacement" + cmptName_,
72 cellMotionBoundaryTypes<scalar>(pointDisplacement_.boundaryField())
74 pointLocation_(
nullptr),
81 coeffDict().
found(
"frozenPointsZone")
82 ? fvMesh_.pointZones().findZoneID(coeffDict().
lookup(
"frozenPointsZone"))
89 fvMesh_.time().timeName(),
97 Info<<
"displacementComponentLaplacianFvMotionSolver:" <<
nl 98 <<
" diffusivity : " << diffusivityPtr_().type() <<
nl 99 <<
" frozenPoints zone : " << frozenPointsZone_ <<
endl;
111 fvMesh_.time().timeName(),
122 Info<<
"displacementComponentLaplacianFvMotionSolver :" 123 <<
" Read pointVectorField " 124 << pointLocation_().name()
125 <<
" to be used for boundary conditions on points." 127 <<
"Boundary conditions:" 128 << pointLocation_().boundaryField().types() <<
endl;
152 if (pointLocation_.valid())
156 Info<<
"displacementComponentLaplacianFvMotionSolver : applying " 157 <<
" boundary conditions on " << pointLocation_().name()
158 <<
" to new point location." 164 pointLocation_().primitiveFieldRef() = fvMesh_.points();
166 pointLocation_().primitiveFieldRef().replace
169 points0_ + pointDisplacement_.primitiveField()
172 pointLocation_().correctBoundaryConditions();
175 if (frozenPointsZone_ != -1)
177 const pointZone& pz = fvMesh_.pointZones()[frozenPointsZone_];
181 label pointi = pz[i];
183 pointLocation_()[pointi][cmpt_] = points0_[pointi];
199 points0_ + pointDisplacement_.primitiveField()
203 if (frozenPointsZone_ != -1)
205 const pointZone& pz = fvMesh_.pointZones()[frozenPointsZone_];
209 label pointi = pz[i];
211 curPoints[pointi][cmpt_] = points0_[pointi];
215 twoDCorrectPoints(curPoints);
226 movePoints(fvMesh_.points());
228 diffusivityPtr_->correct();
229 pointDisplacement_.boundaryFieldRef().updateCoeffs();
235 diffusivityPtr_->operator()(),
237 "laplacian(diffusivity,cellDisplacement)" 252 diffusivityPtr_.reset(
nullptr);
256 coeffDict().
lookup(
"diffusivity")
270 diffusivityPtr_.reset(
nullptr);
274 coeffDict().
lookup(
"diffusivity")
virtual void topoChange(const polyTopoChangeMap &)
Update local data for topology changes.
tmp< GeometricField< Type, pointPatchField, pointMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &) const
Interpolate volField using inverse distance weighting.
#define forAll(list, i)
Loop across all elements in list.
A list of keyword definitions, which are a keyword followed by any number of values (e...
T & ref() const
Return non-const reference or generate a fatal error.
Templated form of IOobject providing type information for file reading and header type checking...
~displacementComponentLaplacianFvMotionSolver()
Destructor.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Calculate the matrix for the laplacian of the field.
virtual void topoChange(const polyTopoChangeMap &)
Update corresponding to the given map.
void replace(const direction, const UList< cmptType > &)
Replace a component field of the field.
conserve primitiveFieldRef()+
Macros for easy insertion into run-time selection tables.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Base class for fvMesh based motionSolvers.
vectorField pointField
pointField is a vectorField.
stressControl lookup("compactNormalStress") >> compactNormalStress
static word timeName(const scalar, const int precision=curPrecision_)
Return time name of given scalar time.
A class for handling words, derived from string.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
const Time & time() const
Return time.
defineTypeNameAndDebug(combustionModel, 0)
virtual void solve()
Solve for motion.
static pointMesh & New(polyMesh &mesh)
A subset of mesh points. The labels of points in the zone can be obtained from the addressing() list...
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.
Mesh consisting of general polyhedral cells.
A class for managing temporary objects.
static autoPtr< motionDiffusivity > New(const fvMesh &mesh, Istream &mdData)
Select null constructed.
Virtual base class for displacement motion solver.
Class containing mesh-to-mesh mapping information.
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
displacementComponentLaplacianFvMotionSolver(const word &name, const polyMesh &, const dictionary &)
Construct from polyMesh and dictionary.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.