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 -------------------------------------------------------------------------------
11 License
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"
37 #include "unitConversion.H"
38 
39 using namespace Foam::constant;
40 using 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)
55 {
56  refValue() = Zero;
57  refGrad() = Zero;
58  valueFraction() = 1.0;
59 }
60 
61 
64 (
66  const fvPatch& p,
68  const fvPatchFieldMapper& mapper
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;
100  valueFraction() = 1.0;
101 
103  }
104 }
105 
106 
109 (
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 
137 updateCoeffs()
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 
169  radiativeIntensityRay& ray =
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())
182  );
183 
184  const scalarField& emissivity = temissivity();
185 
186  const tmp<scalarField> ttransmissivity
187  (
188  boundaryRadiation.transmissivity(patch().index())
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  word qSecName = dom.relfectedFluxName_ + "_0";
226 
227  if (this->db().foundObject<volScalarField>(qSecName))
228  {
229  const volScalarField& qSec =
230  this->db().lookupObject<volScalarField>(qSecName);
231 
232  Ir += qSec.boundaryField()[patch().index()];
233  }
234  }
235 
236  scalarField Iexternal(this->size(), 0.0);
237 
238  if (dom.useExternalBeam())
239  {
240  const vector sunDir = dom.solarCalc().direction();
241  const scalar directSolarRad = dom.solarCalc().directSolarRad();
242 
243  //label nRaysBeam = dom.nRaysBeam();
244  label SunRayId(-1);
245  scalar maxSunRay = -GREAT;
246 
247  // Looking for the ray closest to the Sun direction
248  for (label rayI=0; rayI < dom.nRay(); rayI++)
249  {
250  const vector& iD = dom.IRay(rayI).d();
251  scalar dir = sunDir & iD;
252  if (dir > maxSunRay)
253  {
254  maxSunRay = dir;
255  SunRayId = rayI;
256  }
257  }
258 
259  if (rayId == SunRayId)
260  {
261  const scalarField nAve(n & dom.IRay(rayId).dAve());
262  forAll(Iexternal, faceI)
263  {
264  Iexternal[faceI] = directSolarRad/mag(dom.IRay(rayId).dAve());
265  }
266  }
267  }
268 
269  scalarField Isource(this->size(), 0.0);
270 
271  if (qRadExt_ > 0)
272  {
273  if (mag(qRadExtDir_) > 0)
274  {
275  label rayqoId = -1;
276  scalar maxRay = -GREAT;
277 
278  // Looking for the ray closest to the Sun direction
279  for (label rayI = 0; rayI < dom.nRay(); ++rayI)
280  {
281  const vector& iD = dom.IRay(rayI).d();
282  const scalar dir = qRadExtDir_ & iD;
283 
284  if (dir > maxRay)
285  {
286  maxRay = dir;
287  rayqoId = rayI;
288  }
289  }
290 
291  if (rayId == rayqoId)
292  {
293  forAll(Isource, faceI)
294  {
295  Isource[faceI] += qRadExt_/mag(dom.IRay(rayId).dAve());
296  }
297  }
298  }
299  else
300  {
301  forAll(Iw, faceI)
302  {
303  label rayqoId = -1;
304  scalar maxRay = -GREAT;
305 
306  // Looking for the ray closest to the Sun direction
307  for (label rayI = 0; rayI < dom.nRay(); ++rayI)
308  {
309  const vector& iD = dom.IRay(rayI).d();
310  const scalar dir = -n[faceI] & iD;
311 
312  if (dir > maxRay)
313  {
314  maxRay = dir;
315  rayqoId = rayI;
316  }
317  }
318 
319  if (rayId == rayqoId)
320  {
321  Isource[faceI] += qRadExt_/mag(dom.IRay(rayId).dAve());
322  }
323  }
324  }
325  }
326 
327  forAll(Iw, faceI)
328  {
329  if ((-n[faceI] & myRayId) > 0.0)
330  {
331  // direction out of the wall
332  refGrad()[faceI] = 0.0;
333  valueFraction()[faceI] = 1.0;
334  refValue()[faceI] =
335  Isource[faceI]
336  + Iexternal[faceI]*transmissivity[faceI]
337  + (
338  Ir[faceI]*(scalar(1) - emissivity[faceI])
339  + emissivity[faceI]*physicoChemical::sigma.value()
340  * pow4(Tp[faceI])
341  )/pi;
342 
343  // Emitted heat flux from this ray direction
344  qem[faceI] = refValue()[faceI]*nAve[faceI];
345  }
346  else
347  {
348  // direction into the wall
349  valueFraction()[faceI] = 0.0;
350  refGrad()[faceI] = 0.0;
351  refValue()[faceI] = 0.0; //not used
352 
353  // Incident heat flux on this ray direction
354  qin[faceI] = Iw[faceI]*nAve[faceI];
355  }
356  }
357 
358  UPstream::msgType(oldTag); // Restore tag
359 
360  mixedFvPatchScalarField::updateCoeffs();
361 }
362 
363 
365 (
366  Ostream& os
367 ) const
368 {
370  os.writeEntryIfDifferent<word>("T", "T", TName_);
371  os.writeEntryIfDifferent<scalar>("qRadExt", Zero, qRadExt_);
372  os.writeEntryIfDifferent<vector>("qRadExtDir", Zero, qRadExtDir_);
373 }
374 
375 
376 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
377 
378 namespace Foam
379 {
380 namespace radiation
381 {
383  (
386  );
387 }
388 }
389 
390 
391 // ************************************************************************* //
Different types of constants.
makePatchTypeField(fvPatchScalarField, greyDiffusiveRadiationMixedFvPatchScalarField)
dictionary dict
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
static int incrMsgType(int val=1) noexcept
Increment the message tag for standard messages.
Definition: UPstream.H:1251
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
This boundary condition provides a grey-diffuse condition for radiation intensity, I, for use with the finite-volume discrete-ordinates model (fvDOM), in which the radiation temperature is retrieved from the temperature field boundary condition.
static const boundaryRadiationProperties & New(const fvMesh &mesh, Args &&... args)
Get existing or create a new MeshObject. Registered with typeName.
Definition: MeshObject.C:53
Unit conversion functions.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
static int & msgType() noexcept
Message tag of standard messages.
Definition: UPstream.H:1229
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
void setRayIdLambdaId(const word &name, label &rayId, label &lambdaId) const
Set the rayId and lambdaId from by decomposing an intensity.
Definition: fvDOM.C:754
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
Mathematical constants.
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
A FieldMapper for finite-volume patch fields.
constexpr scalar pi(M_PI)
Vector< scalar > vector
Definition: vector.H:57
greyDiffusiveRadiationMixedFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
OBJstream os(runTime.globalPath()/outputName)
virtual void operator=(const UList< scalar > &)
Definition: fvPatchField.C:391
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
dimensionedScalar pow4(const dimensionedScalar &ds)
const std::string patch
OpenFOAM patch number as a std::string.
label n
Field< vector > vectorField
Specialisation of Field<T> for vector.
virtual void write(Ostream &) const
Write.
volScalarField & p
Finite Volume Discrete Ordinates Method. Solves the RTE equation for n directions in a participating ...
Definition: fvDOM.H:114
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127