Go to the documentation of this file.
42 displacementComponentLaplacianFvMotionSolver,
63 "cellDisplacement" + cmptName_,
64 mesh.time().timeName(),
66 IOobject::READ_IF_PRESENT,
73 pointDisplacement_.dimensions(),
76 cellMotionBoundaryTypes<scalar>(pointDisplacement_.boundaryField())
81 coeffDict().
found(
"interpolation")
91 coeffDict().
found(
"frozenPointsZone")
92 ? fvMesh_.pointZones().findZoneID(coeffDict().
lookup(
"frozenPointsZone"))
98 coeffDict().lookupOrDefault
100 "applyPointLocation",
105 if (applyPointLocation)
114 fvMesh_.time().timeName(),
125 Info<<
"displacementComponentLaplacianFvMotionSolver :"
126 <<
" Read pointVectorField "
127 << pointLocation_().name()
128 <<
" to be used for boundary conditions on points."
130 <<
"Boundary conditions:"
131 << pointLocation_().boundaryField().types() <<
endl;
149 interpolationPtr_->interpolate
155 if (pointLocation_.valid())
159 Info<<
"displacementComponentLaplacianFvMotionSolver : applying "
160 <<
" boundary conditions on " << pointLocation_().name()
161 <<
" to new point location."
167 pointLocation_().internalField() = fvMesh_.points();
169 pointLocation_().internalField().replace
172 points0_ + pointDisplacement_.internalField()
175 pointLocation_().correctBoundaryConditions();
178 if (frozenPointsZone_ != -1)
180 const pointZone& pz = fvMesh_.pointZones()[frozenPointsZone_];
184 label pointI = pz[i];
186 pointLocation_()[pointI][cmpt_] = points0_[pointI];
201 points0_ + pointDisplacement_.internalField()
205 if (frozenPointsZone_ != -1)
207 const pointZone& pz = fvMesh_.pointZones()[frozenPointsZone_];
211 label pointI = pz[i];
213 tcurPoints()[pointI][cmpt_] = points0_[pointI];
217 twoDCorrectPoints(tcurPoints());
228 movePoints(fvMesh_.points());
230 diffusivityPtr_->correct();
231 pointDisplacement_.boundaryField().updateCoeffs();
237 diffusivityPtr_->operator()(),
239 "laplacian(diffusivity,cellDisplacement)"
254 diffusivityPtr_.reset(NULL);
258 coeffDict().
lookup(
"diffusivity")
vectorField pointField
pointField is a vectorField.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
displacementComponentLaplacianFvMotionSolver(const displacementComponentLaplacianFvMotionSolver &)
Disallow default bitwise copy construct.
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
#define forAll(list, i)
Loop across all elements in list.
A class for managing temporary objects.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
A subset of mesh points. The labels of points in the zone can be obtained from the addressing() list.
Virtual base class for displacement motion solver.
virtual void updateMesh(const mapPolyMesh &)
Update topology.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Mesh consisting of general polyhedral cells.
~displacementComponentLaplacianFvMotionSolver()
Destructor.
virtual void updateMesh(const mapPolyMesh &)
Update local data for topology changes.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Base class for fvMesh based motionSolvers.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Calculate the matrix for the laplacian of the field.
static autoPtr< motionDiffusivity > New(const fvMesh &mesh, Istream &mdData)
Select null constructed.
Macros for easy insertion into run-time selection tables.
Generic dimensioned Type class.
conserve internalField()+
SolverPerformance< Type > solve(fvMatrix< Type > &, const dictionary &)
Solve returning the solution statistics given convergence tolerance.
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
virtual void solve()
Solve for motion.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Generic GeometricField class.
defineTypeNameAndDebug(combustionModel, 0)
stressControl lookup("compactNormalStress") >> compactNormalStress