Go to the documentation of this file.
60 phiName_(
dict.getOrDefault<
word>(
"phi",
"phi")),
61 rhoName_(
dict.getOrDefault<
word>(
"rho",
"rho")),
64 length_(
dict.
get<scalar>(
"length")),
65 uniformJump_(
dict.getOrDefault(
"uniformJump", false))
83 phiName_(ptf.phiName_),
84 rhoName_(ptf.rhoName_),
88 uniformJump_(ptf.uniformJump_)
99 phiName_(ptf.phiName_),
100 rhoName_(ptf.rhoName_),
103 length_(ptf.length_),
104 uniformJump_(ptf.uniformJump_)
115 phiName_(ptf.phiName_),
116 rhoName_(ptf.rhoName_),
119 length_(ptf.length_),
120 uniformJump_(ptf.uniformJump_)
135 const fvsPatchField<scalar>& phip =
156 internalField().
group()
160 const scalar t = db().time().timeOutputValue();
161 const scalar
D = D_->
value(t);
162 const scalar
I = I_->value(t);
168 D*turbModel.nu(
patch().index())
177 jump()*
patch().lookupPatchField<volScalarField, scalar>(rhoName_)
183 scalar avePressureJump =
gAverage(jump());
186 Info<<
patch().boundaryMesh().mesh().name() <<
':'
187 <<
patch().name() <<
':'
188 <<
" Average pressure drop :" << avePressureJump
189 <<
" Average velocity :" << aveVelocity
216 porousBafflePressureFvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
fvPatchField< scalar > fvPatchScalarField
scalar timeOutputValue() const
const dimensionSet dimPressure
const DimensionedField< Type, volMesh > & internalField() const
virtual void write(Ostream &) const
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual tmp< volScalarField > nu() const =0
virtual void setJump(const Field< scalar > &jump)
A class for handling words, derived from Foam::string.
constexpr const char *const group
Type gAverage(const FieldField< Field, Type > &f)
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
static const word propertiesName
Ostream & endl(Ostream &os)
const Type & value() const
virtual const word & name() const
const fvMesh & mesh() const noexcept
dimensionedScalar sign(const dimensionedScalar &ds)
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
const fvBoundaryMesh & boundaryMesh() const
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
virtual void write(Ostream &) const
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
This boundary condition provides a jump condition, using the cyclic condition as a base.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const Type & lookupObject(const word &name, const bool recursive=false) const
const GeometricField::Patch & lookupPatchField(const word &name, const GeometricField *=nullptr, const Type *=nullptr) const
porousBafflePressureFvPatchField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Abstract base class for turbulence models (RAS, LES and laminar).
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
OBJstream os(runTime.globalPath()/outputName)
Macros for easy insertion into run-time selection tables.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
virtual void updateCoeffs()
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
virtual void updateCoeffs()
Abstract base class for cyclic coupled interfaces.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
const objectRegistry & db() const
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
const dimensionedScalar & D
const GeometricField::Patch & patchField(const GeometricField &) const
Ostream & writeEntry(const keyType &key, const T &value)
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
Foam::fvPatchFieldMapper.
const fvPatch & patch() const
static word groupName(StringType base, const word &group)
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Generic GeometricField class.
const Time & time() const noexcept
static const Identity< scalar > I
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
virtual tmp< Field< scalar > > jump() const
const word & name() const