Loading...
Searching...
No Matches
greyDiffusiveRadiationMixedFvPatchScalarField.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-2018 OpenFOAM Foundation
9 Copyright (C) 2016-2022 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"
34
35#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),
52 TName_("T"),
53 qRadExt_(0),
54 qRadExtDir_(Zero)
69)
70:
71 mixedFvPatchScalarField(ptf, p, iF, mapper),
72 TName_(ptf.TName_),
73 qRadExt_(ptf.qRadExt_),
74 qRadExtDir_(ptf.qRadExtDir_)
75{}
76
77
80(
81 const fvPatch& p,
83 const dictionary& dict
84)
85:
86 mixedFvPatchScalarField(p, iF),
87 TName_(dict.getOrDefault<word>("T", "T")),
88 qRadExt_(dict.getOrDefault<scalar>("qRadExt", 0)),
89 qRadExtDir_(dict.getOrDefault<vector>("qRadExtDir", Zero))
90{
91 if (this->readMixedEntries(dict))
92 {
93 // Full restart
94 this->readValueEntry(dict, IOobjectOption::MUST_READ);
95 }
96 else
97 {
98 refValue() = Zero;
99 refGrad() = Zero;
111)
112:
113 mixedFvPatchScalarField(ptf),
114 TName_(ptf.TName_),
115 qRadExt_(ptf.qRadExt_),
116 qRadExtDir_(ptf.qRadExtDir_)
117{}
118
119
122(
125)
126:
127 mixedFvPatchScalarField(ptf, iF),
128 TName_(ptf.TName_),
129 qRadExt_(ptf.qRadExt_),
130 qRadExtDir_(ptf.qRadExtDir_)
131{}
132
133
134// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
135
138{
139 if (this->updated())
140 {
141 return;
142 }
143
144 // Since we're inside initEvaluate/evaluate there might be processor
145 // comms underway. Change the tag we use.
146 const int oldTag = UPstream::incrMsgType();
147
148 const auto& Tp = patch().lookupPatchField<volScalarField>(TName_);
149
150 const fvDOM& dom = db().lookupObject<fvDOM>("radiationProperties");
151
152 label rayId = -1;
153 label lambdaId = -1;
154 dom.setRayIdLambdaId(internalField().name(), rayId, lambdaId);
155
156 const label patchi = patch().index();
157
158 if (dom.nLambda() != 1)
159 {
161 << " a grey boundary condition is used with a non-grey "
162 << "absorption model" << nl << exit(FatalError);
163 }
164
165 scalarField& Iw = *this;
166
167 const vectorField n(patch().nf());
168
170 const_cast<radiativeIntensityRay&>(dom.IRay(rayId));
171
172 const scalarField nAve(n & ray.dAve());
173
174 ray.qr().boundaryFieldRef()[patchi] += Iw*nAve;
175
176 const boundaryRadiationProperties& boundaryRadiation =
177 boundaryRadiationProperties::New(internalField().mesh());
178
179 const tmp<scalarField> temissivity
180 (
181 boundaryRadiation.emissivity(patch().index(), 0, nullptr, &Tp)
182 );
183
184 const scalarField& emissivity = temissivity();
185
186 const tmp<scalarField> ttransmissivity
187 (
188 boundaryRadiation.transmissivity(patch().index(), 0, nullptr, &Tp)
189 );
190
191 const scalarField& transmissivity = ttransmissivity();
192
193 scalarField& qem = ray.qem().boundaryFieldRef()[patchi];
194 scalarField& qin = ray.qin().boundaryFieldRef()[patchi];
195
196 const vector& myRayId = dom.IRay(rayId).d();
197
198 scalarField Ir(patch().size(), Zero);
199 forAll(Iw, facei)
200 {
201 for (label rayi=0; rayi < dom.nRay(); rayi++)
202 {
203 const vector& d = dom.IRay(rayi).d();
204
205 if ((-n[facei] & d) < 0.0)
206 {
207 // q into the wall
208 const scalarField& IFace =
209 dom.IRay(rayi).ILambda(lambdaId).boundaryField()[patchi];
210
211 const vector& rayDave = dom.IRay(rayi).dAve();
212 Ir[facei] += IFace[facei]*(n[facei] & rayDave);
213 }
214 }
215 }
216
217 if (dom.useSolarLoad())
218 {
219 // Looking for primary heat flux single band
220 Ir += patch().lookupPatchField<volScalarField>
221 (
222 dom.primaryFluxName_ + "_0"
223 );
224
225 if
226 (
227 const auto* qSec
228 = patch().cfindPatchField<volScalarField>
229 (
230 dom.relfectedFluxName_ + "_0"
231 )
232 )
233 {
234 Ir += *qSec;
235 }
236 }
237
238 scalarField Iexternal(this->size(), 0.0);
239
240 if (dom.useExternalBeam())
241 {
242 const vector sunDir = dom.solarCalc().direction();
243 const scalar directSolarRad = dom.solarCalc().directSolarRad();
244
245 //label nRaysBeam = dom.nRaysBeam();
246 label SunRayId(-1);
247 scalar maxSunRay = -GREAT;
248
249 // Looking for the ray closest to the Sun direction
250 for (label rayI=0; rayI < dom.nRay(); rayI++)
251 {
252 const vector& iD = dom.IRay(rayI).d();
253 scalar dir = sunDir & iD;
254 if (dir > maxSunRay)
255 {
256 maxSunRay = dir;
257 SunRayId = rayI;
258 }
259 }
260
261 if (rayId == SunRayId)
262 {
263 const scalarField nAve(n & dom.IRay(rayId).dAve());
264 forAll(Iexternal, faceI)
265 {
266 Iexternal[faceI] = directSolarRad/mag(dom.IRay(rayId).dAve());
267 }
268 }
269 }
270
271 scalarField Isource(this->size(), 0.0);
272
273 if (qRadExt_ > 0)
274 {
275 if (mag(qRadExtDir_) > 0)
276 {
277 label rayqoId = -1;
278 scalar maxRay = -GREAT;
279
280 // Looking for the ray closest to the Sun direction
281 for (label rayI = 0; rayI < dom.nRay(); ++rayI)
282 {
283 const vector& iD = dom.IRay(rayI).d();
284 const scalar dir = qRadExtDir_ & iD;
285
286 if (dir > maxRay)
287 {
288 maxRay = dir;
289 rayqoId = rayI;
290 }
291 }
292
293 if (rayId == rayqoId)
294 {
295 forAll(Isource, faceI)
296 {
297 Isource[faceI] += qRadExt_/mag(dom.IRay(rayId).dAve());
298 }
299 }
300 }
301 else
302 {
303 forAll(Iw, faceI)
304 {
305 label rayqoId = -1;
306 scalar maxRay = -GREAT;
307
308 // Looking for the ray closest to the Sun direction
309 for (label rayI = 0; rayI < dom.nRay(); ++rayI)
310 {
311 const vector& iD = dom.IRay(rayI).d();
312 const scalar dir = -n[faceI] & iD;
313
314 if (dir > maxRay)
315 {
316 maxRay = dir;
317 rayqoId = rayI;
318 }
319 }
320
321 if (rayId == rayqoId)
322 {
323 Isource[faceI] += qRadExt_/mag(dom.IRay(rayId).dAve());
324 }
325 }
326 }
327 }
328
329 forAll(Iw, faceI)
330 {
331 if ((-n[faceI] & myRayId) > 0.0)
332 {
333 // direction out of the wall
334 refGrad()[faceI] = 0.0;
335 valueFraction()[faceI] = 1.0;
336 refValue()[faceI] =
337 Isource[faceI]
338 + Iexternal[faceI]*transmissivity[faceI]
339 + (
340 Ir[faceI]*(scalar(1) - emissivity[faceI])
341 + emissivity[faceI]*physicoChemical::sigma.value()
342 * pow4(Tp[faceI])
343 )/pi;
344
345 // Emitted heat flux from this ray direction
346 qem[faceI] = refValue()[faceI]*nAve[faceI];
347 }
348 else
349 {
350 // direction into the wall
351 valueFraction()[faceI] = 0.0;
352 refGrad()[faceI] = 0.0;
353 refValue()[faceI] = 0.0; //not used
354
355 // Incident heat flux on this ray direction
356 qin[faceI] = Iw[faceI]*nAve[faceI];
357 }
358 }
360 UPstream::msgType(oldTag); // Restore tag
361
362 mixedFvPatchScalarField::updateCoeffs();
363}
364
365
367(
368 Ostream& os
369) const
370{
372 os.writeEntryIfDifferent<word>("T", "T", TName_);
373 os.writeEntryIfDifferent<scalar>("qRadExt", Zero, qRadExt_);
374 os.writeEntryIfDifferent<vector>("qRadExtDir", Zero, qRadExtDir_);
375}
376
378// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
379
380namespace Foam
381{
382namespace radiation
383{
385 (
388 );
389}
390}
391
392
393// ************************************************************************* //
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.
@ 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
const Type & value() const noexcept
Return const reference to value.
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
void setRayIdLambdaId(const word &name, label &rayId, label &lambdaId) const
Set the rayId and lambdaId from by decomposing an intensity.
Definition fvDOM.C:738
This boundary condition provides a grey-diffuse condition for radiation intensity,...
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
greyDiffusiveRadiationMixedFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal 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.
Definition tmp.H:75
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
dynamicFvMesh & mesh
#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)
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
Different types of constants.
const std::string patch
OpenFOAM patch number as a std::string.
Namespace for radiation modelling.
Namespace for OpenFOAM.
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)
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...
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
Unit conversion functions.