Go to the documentation of this file.
51 mixedFvPatchScalarField(
p, iF)
55 valueFraction() = 1.0;
68 mixedFvPatchScalarField(ptf,
p, iF, mapper)
80 mixedFvPatchScalarField(
p, iF)
82 if (
dict.found(
"refValue"))
96 valueFraction() = 1.0;
109 mixedFvPatchScalarField(ptf)
120 mixedFvPatchScalarField(ptf, iF)
146 dom.setRayIdLambdaId(internalField().
name(), rayId, lambdaId);
148 const label patchi =
patch().index();
150 if (dom.nLambda() == 0)
153 <<
" a non-grey boundary condition is used with a grey "
169 dom.blackBody().bLambda(lambdaId).boundaryField()[patchi]
178 boundaryRadiation.emissivity(
patch().index(), lambdaId)
184 boundaryRadiation.transmissivity(
patch().index(), lambdaId)
186 const scalarField& transmissivity = ttransmissivity();
195 for (label rayi=0; rayi < dom.nRay(); rayi++)
197 const vector& d = dom.IRay(rayi).d();
199 if ((-
n[facei] & d) < 0.0)
203 dom.IRay(rayi).ILambda(lambdaId).boundaryField()[patchi];
205 const vector& rayDave = dom.IRay(rayi).dAve();
206 Ir[facei] += IFace[facei]*(
n[facei] & rayDave);
211 if (dom.useSolarLoad())
216 dom.primaryFluxName_ +
"_" +
name(lambdaId)
219 word qSecName = dom.relfectedFluxName_ +
"_" +
name(lambdaId);
221 if (this->db().foundObject<volScalarField>(qSecName))
232 if (dom.useExternalBeam())
234 const vector sunDir = dom.solarCalc().direction();
235 const scalar directSolarRad =
236 dom.solarCalc().directSolarRad()
237 *dom.spectralDistribution()[lambdaId];
241 scalar maxSunRay = -GREAT;
244 for (label rayI=0; rayI < dom.nRay(); rayI++)
246 const vector& iD = dom.IRay(rayI).d();
247 scalar dir = sunDir & iD;
255 if (rayId == SunRayId)
260 Iexternal[faceI] = directSolarRad/
mag(dom.IRay(rayId).dAve());
267 const vector& d = dom.IRay(rayId).d();
269 if ((-
n[facei] & d) > 0.0)
272 refGrad()[facei] = 0.0;
273 valueFraction()[facei] = 1.0;
275 Iexternal[facei]*transmissivity[facei]
277 Ir[facei]*(1.0 - emissivity[facei])
278 + emissivity[facei]*Eb[facei]
282 qem[facei] += refValue()[facei]*nAve[facei];
287 valueFraction()[facei] = 0.0;
288 refGrad()[facei] = 0.0;
289 refValue()[facei] = 0.0;
292 qin[facei] += Iw[facei]*nAve[facei];
299 mixedFvPatchScalarField::updateCoeffs();
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
A class for handling words, derived from Foam::string.
A class for managing temporary objects.
static constexpr const zero Zero
Global zero (0)
Different types of constants.
volScalarField & qem()
Return non-const access to the boundary emitted heat flux.
static const boundaryRadiationProperties & New(const fvMesh &mesh, Args &&... args)
Get existing or create a new MeshObject.
Radiation intensity for a ray in a given direction.
friend Ostream & operator(Ostream &, const fvPatchField< scalar > &)
const volScalarField & qr() const
Return const access to the boundary heat flux.
#define forAll(list, i)
Loop across all elements in list.
word name(const complex &c)
Return string representation of complex.
This boundary condition provides a wide-band, diffusive radiation condition, where the patch temperat...
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
volScalarField & qin()
Return non-const access to the boundary incident heat flux.
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,...
Macros for easy insertion into run-time selection tables.
errorManipArg< error, int > exit(error &err, const int errNo=1)
const vector & dAve() const
Return the average vector inside the solid angle.
wideBandDiffusiveRadiationMixedFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
static int & msgType()
Message tag of standard messages.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
constexpr scalar pi(M_PI)
const std::string patch
OpenFOAM patch number as a std::string.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Top level model for radiation modelling.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
Foam::fvPatchFieldMapper.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
virtual void write(Ostream &) const
Write.
virtual void operator=(const UList< scalar > &)
const Boundary & boundaryField() const
Return const-reference to the boundary field.
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 ...