Go to the documentation of this file.
45 fixedGradientFvPatchScalarField(
p, iF),
46 heatFlux_(
p.size(), 0.0)
59 fixedGradientFvPatchScalarField(ptf,
p, iF, mapper),
60 heatFlux_(ptf.heatFlux_, mapper)
71 fixedGradientFvPatchScalarField(
p, iF),
72 heatFlux_(
"heatFlux",
dict,
p.size())
77 Info <<
"initialeze the heatFlux"<<
endl;
80 if (
dict.found(
"value"))
96 Info <<
"value can not found,"<<
endl;
104 if (
dict.found(
"gradient"))
112 fixedGradientFvPatchScalarField::updateCoeffs();
113 fixedGradientFvPatchScalarField::evaluate();
120 Info <<
"gardient can not found,"<<
endl;
121 Info <<
"the value is equal to the internal value near cell" <<
endl;
122 Info <<
"gradient default is 0" <<
endl;
135 fixedGradientFvPatchScalarField(wbppsf),
136 heatFlux_(wbppsf.heatFlux_)
146 fixedGradientFvPatchScalarField(wbppsf, iF),
147 heatFlux_(wbppsf.heatFlux_)
159 Info <<
"update the boundary patch field Coeffs now ... " <<
endl;
168 = db().lookupObject<
IOdictionary>(
"transportProperties");
185 const scalar Ap =
gSum(patch().magSf());
192 gradient() =
heatFlux_ / (Ap*alphaEffW * Cp0.value() *
rho0.value() );
196 Info <<
"input Ap is : " << Ap <<
endl;
197 Info <<
"input alphaEffw is : " << alphaEffW <<
endl;
198 Info <<
"input Cp0 is : " << Cp0 <<
endl;
200 Info <<
"the gradient of temperature is : " << gradient() <<
endl;
203 fixedGradientFvPatchScalarField::updateCoeffs();
212 gradient().writeEntry(
"gradient", os);
213 this->writeEntry(
"value", os);
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
virtual void write(Ostream &) const
Write.
virtual void operator=(const UList< Type > &)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
friend Ostream & operator(Ostream &, const fvPatchField< Type > &)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Type gSum(const FieldField< Field, Type > &f)
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
virtual void write(Ostream &) const
Write.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
A list of keyword definitions, which are a keyword followed by any number of values (e....
Macros for easy insertion into run-time selection tables.
Generic dimensioned Type class.
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
void writeEntry(const word &keyword, Ostream &os) const
Write the field as a dictionary entry.
wallHeatFluxFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Foam::fvPatchFieldMapper.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...