44Foam::radiation::greyDiffusiveRadiationMixedFvPatchScalarField::
45greyDiffusiveRadiationMixedFvPatchScalarField
51 mixedFvPatchScalarField(
p, iF),
58 valueFraction() = 1.0;
62Foam::radiation::greyDiffusiveRadiationMixedFvPatchScalarField::
63greyDiffusiveRadiationMixedFvPatchScalarField
71 mixedFvPatchScalarField(ptf,
p, iF, mapper),
73 qRadExt_(ptf.qRadExt_),
74 qRadExtDir_(ptf.qRadExtDir_)
78Foam::radiation::greyDiffusiveRadiationMixedFvPatchScalarField::
79greyDiffusiveRadiationMixedFvPatchScalarField
86 mixedFvPatchScalarField(
p, iF),
87 TName_(
dict.getOrDefault<
word>(
"T",
"T")),
88 qRadExt_(
dict.getOrDefault<scalar>(
"qRadExt", 0)),
105 valueFraction() = 1.0;
112Foam::radiation::greyDiffusiveRadiationMixedFvPatchScalarField::
113greyDiffusiveRadiationMixedFvPatchScalarField
118 mixedFvPatchScalarField(ptf),
120 qRadExt_(ptf.qRadExt_),
121 qRadExtDir_(ptf.qRadExtDir_)
125Foam::radiation::greyDiffusiveRadiationMixedFvPatchScalarField::
126greyDiffusiveRadiationMixedFvPatchScalarField
132 mixedFvPatchScalarField(ptf, iF),
134 qRadExt_(ptf.qRadExt_),
135 qRadExtDir_(ptf.qRadExtDir_)
157 const fvDOM& dom = db().lookupObject<
fvDOM>(
"radiationProperties");
163 const label patchi = patch().index();
165 if (dom.nLambda() != 1)
168 <<
" a grey boundary condition is used with a non-grey "
198 const scalarField& transmissivity = ttransmissivity();
203 const vector& myRayId = dom.IRay(rayId).d();
208 for (label rayi=0; rayi < dom.nRay(); rayi++)
210 const vector& d = dom.IRay(rayi).d();
212 if ((-
n[facei] & d) < 0.0)
216 dom.IRay(rayi).ILambda(lambdaId).boundaryField()[patchi];
218 const vector& rayDave = dom.IRay(rayi).dAve();
219 Ir[facei] += IFace[facei]*(
n[facei] & rayDave);
224 if (dom.useSolarLoad())
229 dom.primaryFluxName_ +
"_0"
232 word qSecName = dom.relfectedFluxName_ +
"_0";
245 if (dom.useExternalBeam())
247 const vector sunDir = dom.solarCalc().direction();
248 const scalar directSolarRad = dom.solarCalc().directSolarRad();
252 scalar maxSunRay = -GREAT;
255 for (label rayI=0; rayI < dom.nRay(); rayI++)
257 const vector& iD = dom.IRay(rayI).d();
258 scalar dir = sunDir & iD;
266 if (rayId == SunRayId)
271 Iexternal[faceI] = directSolarRad/
mag(dom.IRay(rayId).dAve());
280 if (
mag(qRadExtDir_) > 0)
283 scalar maxRay = -GREAT;
286 for (label rayI = 0; rayI < dom.nRay(); ++rayI)
288 const vector& iD = dom.IRay(rayI).d();
289 const scalar dir = qRadExtDir_ & iD;
298 if (rayId == rayqoId)
302 Isource[faceI] += qRadExt_/
mag(dom.IRay(rayId).dAve());
311 scalar maxRay = -GREAT;
314 for (label rayI = 0; rayI < dom.nRay(); ++rayI)
316 const vector& iD = dom.IRay(rayI).d();
317 const scalar dir = -
n[faceI] & iD;
326 if (rayId == rayqoId)
328 Isource[faceI] += qRadExt_/
mag(dom.IRay(rayId).dAve());
336 if ((-
n[faceI] & myRayId) > 0.0)
339 refGrad()[faceI] = 0.0;
340 valueFraction()[faceI] = 1.0;
343 + Iexternal[faceI]*transmissivity[faceI]
345 Ir[faceI]*(scalar(1) - emissivity[faceI])
351 qem[faceI] = refValue()[faceI]*nAve[faceI];
356 valueFraction()[faceI] = 0.0;
357 refGrad()[faceI] = 0.0;
358 refValue()[faceI] = 0.0;
361 qin[faceI] = Iw[faceI]*nAve[faceI];
368 mixedFvPatchScalarField::updateCoeffs();
377 mixedFvPatchScalarField::write(
os);
378 os.writeEntryIfDifferent<
word>(
"T",
"T", TName_);
379 os.writeEntryIfDifferent<scalar>(
"qRadExt",
Zero, qRadExt_);
380 os.writeEntryIfDifferent<
vector>(
"qRadExtDir",
Zero, qRadExtDir_);
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
void size(const label n)
Older name for setAddressableSize.
static int & msgType() noexcept
Message tag of standard messages.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
virtual bool write()
Write the output fields.
A FieldMapper for finite-volume patch fields.
virtual void operator=(const UList< scalar > &)
friend Ostream & operator(Ostream &, const fvPatchField< scalar > &)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Boundary radiation properties holder.
tmp< scalarField > emissivity(const label patchI, const label bandI=0, vectorField *incomingDirection=nullptr, scalarField *T=nullptr) const
Access boundary emissivity on patch.
tmp< scalarField > transmissivity(const label patchI, const label bandI=0, vectorField *incomingDirection=nullptr, scalarField *T=nullptr) const
Access boundary transmissivity on patch.
Finite Volume Discrete Ordinates Method. Solves the RTE equation for n directions in a participating ...
void setRayIdLambdaId(const word &name, label &rayId, label &lambdaId) const
Set the rayId and lambdaId from by decomposing an intensity.
This boundary condition provides a grey-diffuse condition for radiation intensity,...
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Radiation intensity for a ray in a given direction.
const vector & dAve() const
Return the average vector inside the solid angle.
const volScalarField & qr() const
Return const access to the boundary heat flux.
volScalarField & qin()
Return non-const access to the boundary incident heat flux.
volScalarField & qem()
Return non-const access to the boundary emitted heat flux.
A class for managing temporary objects.
A class for handling words, derived from Foam::string.
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
constexpr scalar pi(M_PI)
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
Different types of constants.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
dimensionedScalar pow4(const dimensionedScalar &ds)
static constexpr const zero Zero
Global zero (0)
errorManipArg< error, int > exit(error &err, const int errNo=1)
constexpr char nl
The newline '\n' character (0x0a)
#define forAll(list, i)
Loop across all elements in list.
Unit conversion functions.