Go to the documentation of this file.
41 velocityComponentLaplacianFvMotionSolver,
62 "cellMotionU" + cmptName_,
63 mesh.time().timeName(),
65 IOobject::READ_IF_PRESENT,
72 pointMotionU_.dimensions(),
75 cellMotionBoundaryTypes<scalar>(pointMotionU_.boundaryField())
79 coeffDict().
found(
"interpolation")
102 interpolationPtr_->interpolate
114 + fvMesh_.time().deltaTValue()*pointMotionU_.internalField()
117 twoDCorrectPoints(tcurPoints());
127 movePoints(fvMesh_.points());
129 diffusivityPtr_->correct();
130 pointMotionU_.boundaryField().updateCoeffs();
136 diffusivityPtr_->operator()(),
138 "laplacian(diffusivity,cellMotionU)"
153 diffusivityPtr_.reset(NULL);
157 coeffDict().
lookup(
"diffusivity")
virtual void updateMesh(const mapPolyMesh &)
Update local data for topology changes.
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...
velocityComponentLaplacianFvMotionSolver(const velocityComponentLaplacianFvMotionSolver &)
Disallow default bitwise copy construct.
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
A class for managing temporary objects.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
~velocityComponentLaplacianFvMotionSolver()
Destructor.
Mesh consisting of general polyhedral cells.
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.
Virtual base class for velocity motion solver.
Macros for easy insertion into run-time selection tables.
Generic dimensioned Type class.
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
SolverPerformance< Type > solve(fvMatrix< Type > &, const dictionary &)
Solve returning the solution statistics given convergence tolerance.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
virtual void solve()
Solve for motion.
defineTypeNameAndDebug(combustionModel, 0)
stressControl lookup("compactNormalStress") >> compactNormalStress
virtual void updateMesh(const mapPolyMesh &)
Update topology.