Go to the documentation of this file.
35 namespace functionObjects
63 mesh_.time().timeName(),
73 mesh_.objectRegistry::store(procFieldPtr);
86 doCells_ =
dict.getOrDefault(
"calculateCells",
true);
96 mesh_.time().constant(),
102 dict.subDict(
"geometry"),
118 volScalarField::Boundary& bfld =
distance.boundaryFieldRef();
138 dist[i] =
mag(nearestInfo[i].hitPoint()-fc[i]);
140 bfld[patchi] == dist;
158 forAll(nearestInfo, celli)
160 distance[celli] =
mag(nearestInfo[celli].hitPoint()-cc[celli]);
163 distance.correctBoundaryConditions();
171 Log <<
" functionObjects::" <<
type() <<
" " <<
name()
172 <<
" writing distance-to-surface field" <<
endl;
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
virtual bool write()
Write the interpolated fields.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
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)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
#define forAll(list, i)
Loop across all elements in list.
Switch doCells_
Switch to calculate distance-to-cells.
word name(const complex &c)
Return string representation of complex.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual bool execute()
Calculate the interpolated fields.
static bool constraintType(const word &pt)
Return true if the given type is a constraint type.
const volVectorField & C() const
Return cell centres as volVectorField.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Macros for easy insertion into run-time selection tables.
surfaceDistance(const word &name, const Time &runTime, const dictionary &dict)
Construct for given objectRegistry and dictionary.
scalar distance(const vector &p1, const vector &p2)
Type & lookupObjectRef(const word &name, const bool recursive=false) const
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
const std::string patch
OpenFOAM patch number as a std::string.
virtual const word & type() const =0
Runtime type information.
Container for searchableSurfaces. The collection is specified as a dictionary. For example,...
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
const fvMesh & mesh_
Reference to the fvMesh.
autoPtr< searchableSurfaces > geomPtr_
Geometry.
virtual bool read(const dictionary &)
Read the controls.
defineTypeNameAndDebug(ObukhovLength, 0)
const Boundary & boundaryField() const
Return const-reference to the boundary field.