Go to the documentation of this file.
64 this->
db().foundObject<compressible::turbulenceModel>
115 return turbModel.
kappaEff(patch().index());
127 return thermo.kappa(patch().index());
134 <<
" on mesh " << this->db().name() <<
" patch "
136 <<
" could not find a method in. Methods are: "
137 << methodTypeNames_.toc()
138 <<
" Not turbulenceModel or thermophysicalProperties"
160 const fvPatch& nbrPatch = regionCoupledPatch_.neighbFvPatch();
176 regionCoupledPatch_.regionCoupledPatch().interpolate
184 regionCoupledPatch_.regionCoupledPatch().interpolate
186 nbrPatch.
nf() & nbrPatch.
delta()
190 const scalarField nbrAlphaDelta(nbrAlpha/nbrDeltas);
197 scalar di = alphaDelta[faceI];
198 scalar dni = nbrAlphaDelta[faceI];
200 w[faceI] = di/(di + dni);
217 regionCoupledPatch_(refCast<const regionCoupledBaseFvPatch>(
p)),
234 regionCoupledPatch_(refCast<const regionCoupledBaseFvPatch>(
p)),
250 regionCoupledPatch_(refCast<const regionCoupledBaseFvPatch>(
p)),
256 if (!isA<regionCoupledBase>(this->patch().patch()))
259 <<
"' not type '" << regionCoupledBase::typeName <<
"'"
260 <<
"\n for patch " <<
p.name()
304 regionCoupledPatch_.patch().deltaCoeffs()
305 *(*
this - patchInternalField());
337 lWeights*patchInternalTemperatureField()
338 + (1.0 - lWeights)*patchNeighbourTemperatureField(),
351 const fvPatch& nbrPatch = regionCoupledPatch_.neighbFvPatch();
359 nbrThermoPtr_->T().internalField(), nbrFaceCells
364 regionCoupledPatch_.regionCoupledPatch().interpolate(nbrIntT)
378 const fvPatch& nbrPatch = regionCoupledPatch_.neighbFvPatch();
384 nbrThermoPtr_->T().internalField(), nbrFaceCells
388 regionCoupledPatch_.regionCoupledPatch().interpolate(nbrIntT);
397 const labelUList& faceCells = regionCoupledPatch_.faceCells();
401 new scalarField(thermoPtr_->T().internalField(), faceCells)
427 myHE = thermoPtr_->he(pp, Tp,
patchi);
435 myHE[facei] = psiInternal[regionCoupledPatch_.faceCells()[facei]];
440 const labelUList& faceCells = regionCoupledPatch_.faceCells();
444 result[faceCells[elemI]] -= coeffs[elemI]*myHE[elemI];
464 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.
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::blocking)
Evaluate the patch field, sets Updated to false.
volScalarField w(IOobject("w", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE), mesh, dimensionedScalar("w", dimensionSet(0, 0, 0, 0, 0, 0, 0), 0.0))
virtual const polyMesh & nbrMesh() const
Returns neighbour polyMesh.
static const NamedEnum< kappaMethodType, 3 > methodTypeNames_
Methof to extract kappa.
energyRegionCoupledFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
#define forAll(list, i)
Loop across all elements in list.
A class for managing temporary objects.
kappaMethodType method_
How to get K.
static const word dictName
virtual tmp< vectorField > delta() const
Return cell-centre to face-centre vector.
rDeltaT dimensionedInternalField()
tmp< scalarField > patchNeighbourTemperatureField() const
Return nbr temperature internal field.
friend Ostream & operator(Ostream &, const Field< Type > &)
static const word propertiesName
Default name of the turbulence properties dictionary.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
virtual void evaluate(const Pstream::commsTypes commsType)
Evaluate the patch field.
Abstract base-class for fluid and solid thermodynamic properties.
const basicThermo * thermoPtr_
AutoPtr to my thermo.
const basicThermo * nbrThermoPtr_
AutoPtr to nbr thermo.
tmp< scalarField > patchInternalTemperatureField() const
Return local temperature internal field.
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
#define Debug(var)
Report a variable name and value.
tmp< vectorField > nf() const
Return face normals.
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Pre-declare SubField and related Field type.
tmp< Field< Type > > T() const
Return the field transpose (only defined for second rank tensors)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
tmp< scalarField > weights() const
Local weight for this coupled field.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
virtual void updateInterfaceMatrix(scalarField &, const scalarField &, const scalarField &, const direction, const Pstream::commsTypes commsType) const =0
Inherit updateInterfaceMatrix from lduInterfaceField.
A list of keyword definitions, which are a keyword followed by any number of values (e....
const regionCoupledBaseFvPatch & regionCoupledPatch_
Local reference to region couple patch.
Macros for easy insertion into run-time selection tables.
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
To & refCast(From &r)
Reference type cast template function.
errorManipArg< error, int > exit(error &err, const int errNo=1)
conserve internalField()+
virtual tmp< volScalarField > kappaEff() const
Return the effective turbulent thermal diffusivity for temperature.
commsTypes
Types of communications.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
tmp< scalarField > kappa() const
Return kappa.
const objectRegistry & db() const
Return local objectRegistry.
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
void setMethod() const
Set method.
Energy region coupled implicit boundary condition. The fvPatch is treated as uncoupled from the delta...
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
const GeometricField::PatchFieldType & lookupPatchField(const word &name, const GeometricField *=NULL, const Type *=NULL) const
Lookup and return the patchField of the named field from the.
virtual const labelUList & faceCells() const
Return faceCells.
const polyPatch & patch() const
Return the polyPatch.
Foam::fvPatchFieldMapper.
virtual void write(Ostream &) const
Write.
virtual tmp< scalarField > patchNeighbourField() const
Return neighbour coupled internal cell data.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Initialise the NamedEnum HashTable from the static list of names.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
virtual tmp< scalarField > snGrad() const
Return patch-normal gradient.