Go to the documentation of this file.
51 mixedFvPatchScalarField(
p, iF),
56 valueFraction() = 1.0;
69 mixedFvPatchScalarField(ptf,
p, iF, mapper),
82 mixedFvPatchScalarField(
p, iF),
83 TName_(
dict.getOrDefault<
word>(
"T",
"T"))
85 if (
dict.found(
"refValue"))
99 valueFraction() = 1.0;
113 mixedFvPatchScalarField(ptf),
125 mixedFvPatchScalarField(ptf, iF),
148 const fvDOM& dom = db().lookupObject<
fvDOM>(
"radiationProperties");
154 const label patchi =
patch().index();
156 if (dom.nLambda() != 1)
159 <<
" a grey boundary condition is used with a non-grey "
167 radiativeIntensityRay& ray =
168 const_cast<radiativeIntensityRay&
>(dom.IRay(rayId));
172 ray.qr().boundaryFieldRef()[patchi] += Iw*nAve;
174 const boundaryRadiationProperties& boundaryRadiation =
177 const tmp<scalarField> temissivity
179 boundaryRadiation.emissivity(
patch().index())
184 const tmp<scalarField> ttransmissivity
186 boundaryRadiation.transmissivity(
patch().index())
189 const scalarField& transmissivity = ttransmissivity();
191 scalarField& qem = ray.qem().boundaryFieldRef()[patchi];
192 scalarField& qin = ray.qin().boundaryFieldRef()[patchi];
194 const vector& myRayId = dom.IRay(rayId).d();
199 for (label rayi=0; rayi < dom.nRay(); rayi++)
201 const vector& d = dom.IRay(rayi).d();
203 if ((-
n[facei] & d) < 0.0)
207 dom.IRay(rayi).ILambda(lambdaId).boundaryField()[patchi];
209 const vector& rayDave = dom.IRay(rayi).dAve();
210 Ir[facei] += IFace[facei]*(
n[facei] & rayDave);
215 if (dom.useSolarLoad())
220 dom.primaryFluxName_ +
"_0"
223 word qSecName = dom.relfectedFluxName_ +
"_0";
236 if (dom.useExternalBeam())
238 const vector sunDir = dom.solarCalc().direction();
239 const scalar directSolarRad = dom.solarCalc().directSolarRad();
243 scalar maxSunRay = -GREAT;
246 for (label rayI=0; rayI < dom.nRay(); rayI++)
248 const vector& iD = dom.IRay(rayI).d();
249 scalar dir = sunDir & iD;
257 if (rayId == SunRayId)
262 Iexternal[faceI] = directSolarRad/
mag(dom.IRay(rayId).dAve());
269 if ((-
n[faceI] & myRayId) > 0.0)
272 refGrad()[faceI] = 0.0;
273 valueFraction()[faceI] = 1.0;
275 Iexternal[faceI]*transmissivity[faceI]
277 Ir[faceI]*(scalar(1) - emissivity[faceI])
283 qem[faceI] = refValue()[faceI]*nAve[faceI];
288 valueFraction()[faceI] = 0.0;
289 refGrad()[faceI] = 0.0;
290 refValue()[faceI] = 0.0;
293 qin[faceI] = Iw[faceI]*nAve[faceI];
300 mixedFvPatchScalarField::updateCoeffs();
310 os.writeEntryIfDifferent<
word>(
"T",
"T", TName_);
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
virtual void operator=(const UList< Type > &)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
A class for handling words, derived from Foam::string.
const volScalarField & qr() const
makePatchTypeField(fvPatchScalarField, greyDiffusiveRadiationMixedFvPatchScalarField)
A class for managing temporary objects.
static constexpr const zero Zero
Different types of constants.
static const boundaryRadiationProperties & New(const fvMesh &mesh, Args &&... args)
greyDiffusiveRadiationMixedFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Radiation intensity for a ray in a given direction.
Unit conversion functions.
friend Ostream & operator(Ostream &, const fvPatchField< Type > &)
This boundary condition provides a grey-diffuse condition for radiation intensity,...
Field< vector > vectorField
Specialisation of Field<T> for vector.
tmp< scalarField > emissivity(const label patchI, const label bandI=0, vectorField *incomingDirection=nullptr, scalarField *T=nullptr) const
dimensionedScalar pow4(const dimensionedScalar &ds)
Generic templated field type.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
const vector & dAve() const
GeometricField< scalar, fvPatchField, volMesh > volScalarField
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
Boundary radiation properties holder.
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.
Vector< scalar > vector
A scalar version of the templated Vector.
errorManipArg< error, int > exit(error &err, const int errNo=1)
#define FatalErrorInFunction
static int & msgType() noexcept
constexpr scalar pi(M_PI)
virtual void write(Ostream &) const
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
virtual void updateCoeffs()
const dimensionedScalar sigma
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Foam::fvPatchFieldMapper.
tmp< scalarField > transmissivity(const label patchI, const label bandI=0, vectorField *incomingDirection=nullptr, scalarField *T=nullptr) const
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Generic GeometricField class.
void setRayIdLambdaId(const word &name, label &rayId, label &lambdaId) const
const Boundary & boundaryField() const
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Finite Volume Discrete Ordinates Method. Solves the RTE equation for n directions in a participating ...