Loading...
Searching...
No Matches
wideBandDiffusiveRadiationMixedFvPatchScalarField.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2016-2025 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
31#include "fvPatchFieldMapper.H"
32#include "volFields.H"
33
34#include "fvDOM.H"
36#include "constants.H"
38
39using namespace Foam::constant;
40using namespace Foam::constant::mathematical;
41
42// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43
46(
47 const fvPatch& p,
49)
50:
51 mixedFvPatchScalarField(p, iF)
63 const fvPatch& p,
65 const fvPatchFieldMapper& mapper
66)
67:
68 mixedFvPatchScalarField(ptf, p, iF, mapper)
69{}
70
71
74(
75 const fvPatch& p,
77 const dictionary& dict
78)
79:
80 mixedFvPatchScalarField(p, iF, dict, IOobjectOption::NO_READ)
81{
82 if (this->readMixedEntries(dict))
83 {
84 // Full restart
85 this->readValueEntry(dict, IOobjectOption::MUST_READ);
86 }
87 else
88 {
89 refValue() = Zero;
90 refGrad() = Zero;
100(
110(
113)
115 mixedFvPatchScalarField(ptf, iF)
116{}
117
118
119// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
120
123{
124 if (this->updated())
125 {
126 return;
127 }
128
129 // Since we're inside initEvaluate/evaluate there might be processor
130 // comms underway. Change the tag we use.
131 const int oldTag = UPstream::incrMsgType();
132
134 db().lookupObject<radiationModel>("radiationProperties");
135
137
138 label rayId = -1;
139 label lambdaId = -1;
140 dom.setRayIdLambdaId(internalField().name(), rayId, lambdaId);
141
142 const label patchi = patch().index();
143
144 if (dom.nLambda() == 0)
145 {
147 << " a non-grey boundary condition is used with a grey "
148 << "absorption model" << nl << exit(FatalError);
149 }
150
151 scalarField& Iw = *this;
152 const vectorField n(patch().Sf()/patch().magSf());
153
155 const_cast<radiativeIntensityRay&>(dom.IRay(rayId));
156
157 const scalarField nAve(n & ray.dAve());
158
159 ray.qr().boundaryFieldRef()[patchi] += Iw*nAve;
160
161 const scalarField Eb
162 (
163 dom.blackBody().bLambda(lambdaId).boundaryField()[patchi]
164 );
165
166 const boundaryRadiationProperties& boundaryRadiation =
167 boundaryRadiationProperties::New(internalField().mesh());
168
169
170 const auto& Tp = radiation.T().boundaryField()[patchi];
171
172 const tmp<scalarField> temissivity
173 (
174 boundaryRadiation.emissivity(patchi, lambdaId, nullptr, &Tp)
175 );
176 const scalarField& emissivity = temissivity();
177
178 const tmp<scalarField> ttransmissivity
179 (
180 boundaryRadiation.transmissivity(patchi, lambdaId, nullptr, &Tp)
181 );
182 const scalarField& transmissivity = ttransmissivity();
183
184 scalarField& qem = ray.qem().boundaryFieldRef()[patchi];
185 scalarField& qin = ray.qin().boundaryFieldRef()[patchi];
186
187 // Calculate Ir into the wall on the same lambdaId
188 scalarField Ir(patch().size(), Zero);
189 forAll(Iw, facei)
190 {
191 for (label rayi=0; rayi < dom.nRay(); rayi++)
192 {
193 const vector& d = dom.IRay(rayi).d();
194
195 if ((-n[facei] & d) < 0.0)
196 {
197 // q into the wall
198 const scalarField& IFace =
199 dom.IRay(rayi).ILambda(lambdaId).boundaryField()[patchi];
200
201 const vector& rayDave = dom.IRay(rayi).dAve();
202 Ir[facei] += IFace[facei]*(n[facei] & rayDave);
203 }
204 }
205 }
206
207 if (dom.useSolarLoad())
208 {
209 // Looking for primary heat flux single band
210 Ir += patch().lookupPatchField<volScalarField>
211 (
212 dom.primaryFluxName_ + "_" + name(lambdaId)
213 );
214
215 if
216 (
217 const auto* qSec
218 = patch().cfindPatchField<volScalarField>
219 (
220 dom.relfectedFluxName_ + "_" + name(lambdaId)
221 )
222 )
223 {
224 Ir += *qSec;
225 }
226 }
227
228 scalarField Iexternal(this->size(), Zero);
229
230 if (dom.useExternalBeam())
231 {
232 const vector sunDir = dom.solarCalc().direction();
233 const scalar directSolarRad =
234 dom.solarCalc().directSolarRad()
235 *dom.spectralDistribution()[lambdaId];
236
237 //label nRaysBeam = dom.nRaysBeam();
238 label SunRayId(-1);
239 scalar maxSunRay = -GREAT;
240
241 // Looking for the ray closest to the Sun direction
242 for (label rayI=0; rayI < dom.nRay(); rayI++)
243 {
244 const vector& iD = dom.IRay(rayI).d();
245 scalar dir = sunDir & iD;
246 if (dir > maxSunRay)
247 {
248 maxSunRay = dir;
249 SunRayId = rayI;
250 }
251 }
252
253 if (rayId == SunRayId)
254 {
255 const scalarField nAve(n & dom.IRay(rayId).dAve());
256 forAll(Iexternal, faceI)
257 {
258 Iexternal[faceI] = directSolarRad/mag(dom.IRay(rayId).dAve());
259 }
260 }
261 }
262
263 forAll(Iw, facei)
264 {
265 const vector& d = dom.IRay(rayId).d();
266
267 if ((-n[facei] & d) > 0.0)
268 {
269 // direction out of the wall
270 refGrad()[facei] = 0.0;
271 valueFraction()[facei] = 1.0;
272 refValue()[facei] =
273 Iexternal[facei]*transmissivity[facei]
274 + (
275 Ir[facei]*(1.0 - emissivity[facei])
276 + emissivity[facei]*Eb[facei]
277 )/pi;
278
279 // Emitted heat flux from this ray direction (sum over lambdaId)
280 qem[facei] += refValue()[facei]*nAve[facei];
281 }
282 else
283 {
284 // direction into the wall
285 valueFraction()[facei] = 0.0;
286 refGrad()[facei] = 0.0;
287 refValue()[facei] = 0.0; //not used
288
289 // Incident heat flux on this ray direction (sum over lambdaId)
290 qin[facei] += Iw[facei]*nAve[facei];
291 }
292 }
294 UPstream::msgType(oldTag); // Restore tag
295
296 mixedFvPatchScalarField::updateCoeffs();
297}
298
299
301(
302 Ostream& os
303) const
304{
306}
307
309// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
310
311namespace Foam
312{
313namespace radiation
314{
316 (
319 );
320}
321}
322
323
324// ************************************************************************* //
label n
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.
A simple container of IOobject preferences. Can also be used for general handling of read/no-read/rea...
@ MUST_READ
Reading required.
static FOAM_NO_DANGLING_REFERENCE const boundaryRadiationProperties & New(const fvMesh &mesh, Args &&... args)
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
static int & msgType() noexcept
Message tag of standard messages.
Definition UPstream.H:1926
static int incrMsgType(int val=1) noexcept
Increment the message tag for standard messages.
Definition UPstream.H:1948
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A FieldMapper for finite-volume patch fields.
virtual void operator=(const UList< scalar > &)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
virtual void write(Ostream &) const
Write.
tmp< scalarField > transmissivity(const label patchI, const label bandI=0, const vectorField *incomingDirection=nullptr, const scalarField *T=nullptr) const
Access boundary transmissivity on patch.
tmp< scalarField > emissivity(const label patchI, const label bandI=0, const vectorField *incomingDirection=nullptr, const scalarField *T=nullptr) const
Access boundary emissivity on patch.
Finite Volume Discrete Ordinates Method. Solves the RTE equation for n directions in a participating ...
Definition fvDOM.H:117
Top level model for radiation modelling.
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.
This boundary condition provides a wide-band, diffusive radiation condition, where the patch temperat...
wideBandDiffusiveRadiationMixedFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A class for managing temporary objects.
Definition tmp.H:75
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
volScalarField & p
dynamicFvMesh & mesh
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete fvPatchField type and add to run-time tables Example, (fvPatchScalarField,...
auto & name
Mathematical constants.
constexpr scalar pi(M_PI)
Different types of constants.
const std::string patch
OpenFOAM patch number as a std::string.
Namespace for radiation modelling.
Namespace for OpenFOAM.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Field< vector > vectorField
Specialisation of Field<T> for vector.
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
Vector< scalar > vector
Definition vector.H:57
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
fvPatchField< scalar > fvPatchScalarField
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299