Go to the documentation of this file.
54 this->refValue() =
Zero;
55 this->refGrad() =
Zero;
56 this->valueFraction() = 0.0;
70 phiName_(ptf.phiName_),
71 rhoName_(ptf.rhoName_),
72 fieldInf_(ptf.fieldInf_),
103 this->refValue() = *
this;
104 this->refGrad() =
Zero;
105 this->valueFraction() = 0.0;
114 <<
"unphysical lInf specified (lInf < 0)" <<
nl
115 <<
" on patch " << this->
patch().name()
116 <<
" of field " << this->internalField().name()
117 <<
" in file " << this->internalField().objectPath()
160 this->db().objectRegistry::template lookupObject<surfaceScalarField>
164 this->
patch().template lookupPatchField<surfaceScalarField, scalar>
172 this->
patch().template lookupPatchField<volScalarField, scalar>
177 return phip/(rhop*this->
patch().magSf());
181 return phip/this->
patch().magSf();
194 const fvMesh&
mesh = this->internalField().mesh();
200 scalar deltaT = this->db().time().deltaTValue();
203 this->db().objectRegistry::template
204 lookupObject<GeometricField<Type, fvPatchField, volMesh>>
206 this->internalField().name()
216 label patchi = this->
patch().index();
233 field.oldTime().boundaryField()[patchi] +
k*fieldInf_
236 this->valueFraction() = (1.0 +
k)/(1.0 +
alpha +
k);
242 2.0*
field.oldTime().boundaryField()[patchi]
243 - 0.5*
field.oldTime().oldTime().boundaryField()[patchi]
247 this->valueFraction() = (1.5 +
k)/(1.5 +
alpha +
k);
268 field.oldTime().boundaryField()[patchi] +
k*fieldInf_
271 this->valueFraction() = (1.0 +
k)/(1.0 +
alpha +
k);
277 <<
" on patch " << this->
patch().name()
278 <<
" of field " << this->internalField().name()
279 <<
" in file " << this->internalField().objectPath()
291 this->refValue() =
field.oldTime().boundaryField()[patchi];
293 this->valueFraction() = 1.0/(1.0 +
alpha);
299 2.0*
field.oldTime().boundaryField()[patchi]
300 - 0.5*
field.oldTime().oldTime().boundaryField()[patchi]
303 this->valueFraction() = 1.5/(1.5 +
alpha);
319 this->refValue() =
field.oldTime().boundaryField()[patchi];
321 this->valueFraction() = 1.0/(1.0 +
alpha);
327 <<
"\n on patch " << this->
patch().name()
328 <<
" of field " << this->internalField().name()
329 <<
" in file " << this->internalField().objectPath()
352 this->writeEntry(
"value", os);
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
virtual void write(Ostream &) const
Write.
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
Type fieldInf_
Field value of the far-field.
static const volScalarField & localRDeltaT(const fvMesh &mesh)
Return the reciprocal of the local time-step.
A class for handling words, derived from Foam::string.
A class for managing temporary objects.
static constexpr const zero Zero
Global zero (0)
const dimensionSet dimVelocity
ITstream & ddtScheme(const word &name) const
const dimensionSet dimDensity
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Local time-step first-order Euler implicit/explicit ddt.
const dimensionSet dimArea(sqr(dimLength))
Generic templated field type.
Second-oder Crank-Nicolson implicit ddt using the current and previous time-step fields as well as th...
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
virtual void write(Ostream &) const
Write.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
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.
Mesh data needed to do the Finite Volume discretisation.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
This boundary condition provides a base class for 'mixed' type boundary conditions,...
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
advectiveFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
errorManipArg< error, int > exit(error &err, const int errNo=1)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Second-order backward-differencing ddt using the current and two previous time-step values.
const std::string patch
OpenFOAM patch number as a std::string.
label k
Boltzmann constant.
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Foam::fvPatchFieldMapper.
virtual tmp< scalarField > advectionSpeed() const
Calculate and return the advection speed at the boundary.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
scalar lInf_
Relaxation length-scale.
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
virtual void operator=(const UList< Type > &)
Basic first-order Euler implicit/explicit ddt using only the current and previous time-step values.
This boundary condition provides an advective outflow condition, based on solving DDt(W,...
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
word phiName_
Name of the flux transporting the field.