Go to the documentation of this file.
42 void Foam::wallDist::constructn()
const
44 n_ = tmp<volVectorField>
56 patchDistMethod::patchTypes<vector>(
mesh(), patchIDs_)
62 volVectorField::Boundary& nbf = n_.ref().boundaryFieldRef();
64 for (
const label patchi : patchIDs_)
66 nbf[patchi] ==
patches[patchi].nf();
73 Foam::wallDist::wallDist
77 const word& patchTypeName
82 patchTypeName_(patchTypeName),
87 patchTypeName_ &
"Dist"
103 "y" & patchTypeName_,
109 patchDistMethod::patchTypes<scalar>(
mesh, patchIDs_)
111 nRequired_(dict_.getOrDefault(
"nRequired",
false)),
113 updateInterval_(dict_.getOrDefault<label>(
"updateInterval", 1)),
125 Foam::wallDist::wallDist
128 const word& defaultPatchDistMethod,
130 const word& patchTypeName
135 patchTypeName_(patchTypeName),
140 patchTypeName_ &
"Dist"
150 defaultPatchDistMethod
157 "y" & patchTypeName_,
163 patchDistMethod::patchTypes<scalar>(
mesh, patchIDs_)
165 nRequired_(dict_.getOrDefault(
"nRequired",
false)),
167 updateInterval_(dict_.getOrDefault<label>(
"updateInterval", 1)),
203 <<
"n requested but 'nRequired' not specified in the "
204 << (patchTypeName_ &
"Dist") <<
" dictionary" <<
nl
205 <<
" Recalculating y and n fields." <<
endl;
209 pdm_->correct(y_, n_.ref());
220 (updateInterval_ != 0)
221 && ((mesh_.time().timeIndex() % updateInterval_) == 0)
224 requireUpdate_ =
true;
227 if (requireUpdate_ && pdm_->movePoints())
231 requireUpdate_ =
false;
235 return pdm_->correct(y_, n_.ref());
239 return pdm_->correct(y_);
249 pdm_->updateMesh(mpm);
255 requireUpdate_ =
true;
Defines the attributes of an object for which implicit objectRegistry management is supported,...
virtual void updateMesh(const mapPolyMesh &)
Update the y-field when the mesh changes.
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Dimensionless.
A class for handling words, derived from Foam::string.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
static constexpr const zero Zero
Global zero (0)
Selector class for finite volume differencing schemes. fvMesh is derived from fvShemes so that all fi...
static word timeName(const scalar t, const int precision=precision_)
IOobject(const IOobject &)=default
Copy construct.
Ostream & endl(Ostream &os)
Add newline and flush stream.
const Time & time() const
Return time.
static autoPtr< patchDistMethod > New(const dictionary &dict, const fvMesh &mesh, const labelHashSet &patchIDs)
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
PtrList< fvPatch > fvPatchList
container classes for fvPatch
Interface to run-time selectable methods to calculate the distance-to-wall and normal-to-wall fields.
const fvMesh & mesh() const
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual bool movePoints()
Update the y-field when the mesh moves.
const volVectorField & n() const
Return reference to cached normal-to-wall field.
Mesh data needed to do the Finite Volume discretisation.
GeometricField< vector, fvPatchField, volMesh > volVectorField
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
#define DebugInfo
Report an information message using Foam::Info.
const polyBoundaryMesh & patches
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
const Time & time() const
Return the top-level database.
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Templated abstract base-class for optional mesh objects used to automate their allocation to the mesh...
bool isNull(const T *ptr)
True if ptr is a pointer (of type T) to the nullObject.
static const GeometricField< vector, fvPatchField, volMesh > & null()
Return a null geometric field.
defineTypeNameAndDebug(combustionModel, 0)
#define WarningInFunction
Report a warning using Foam::Warning.
virtual ~wallDist()
Destructor.