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() = 0.0;
57  refGrad() = 0.0;
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 (dict.found("refValue"))
92  {
94  (
95  scalarField("value", dict, p.size())
96  );
97  refValue() = scalarField("refValue", dict, p.size());
98  refGrad() = scalarField("refGradient", dict, p.size());
99  valueFraction() = scalarField("valueFraction", dict, p.size());
100  }
101  else
102  {
103  refValue() = 0.0;
104  refGrad() = 0.0;
105  valueFraction() = 1.0;
106 
108  }
109 }
110 
111 
114 (
116 )
117 :
118  mixedFvPatchScalarField(ptf),
119  TName_(ptf.TName_),
120  qRadExt_(ptf.qRadExt_),
121  qRadExtDir_(ptf.qRadExtDir_)
122 {}
123 
124 
127 (
130 )
131 :
132  mixedFvPatchScalarField(ptf, iF),
133  TName_(ptf.TName_),
134  qRadExt_(ptf.qRadExt_),
135  qRadExtDir_(ptf.qRadExtDir_)
136 {}
137 
138 
139 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
140 
142 updateCoeffs()
143 {
144  if (this->updated())
145  {
146  return;
147  }
148 
149  // Since we're inside initEvaluate/evaluate there might be processor
150  // comms underway. Change the tag we use.
151  int oldTag = UPstream::msgType();
152  UPstream::msgType() = oldTag+1;
153 
154  const scalarField& Tp =
155  patch().lookupPatchField<volScalarField, scalar>(TName_);
156 
157  const fvDOM& dom = db().lookupObject<fvDOM>("radiationProperties");
158 
159  label rayId = -1;
160  label lambdaId = -1;
161  dom.setRayIdLambdaId(internalField().name(), rayId, lambdaId);
162 
163  const label patchi = patch().index();
164 
165  if (dom.nLambda() != 1)
166  {
168  << " a grey boundary condition is used with a non-grey "
169  << "absorption model" << nl << exit(FatalError);
170  }
171 
172  scalarField& Iw = *this;
173 
174  const vectorField n(patch().nf());
175 
176  radiativeIntensityRay& ray =
177  const_cast<radiativeIntensityRay&>(dom.IRay(rayId));
178 
179  const scalarField nAve(n & ray.dAve());
180 
181  ray.qr().boundaryFieldRef()[patchi] += Iw*nAve;
182 
183  const boundaryRadiationProperties& boundaryRadiation =
184  boundaryRadiationProperties::New(internalField().mesh());
185 
186  const tmp<scalarField> temissivity
187  (
188  boundaryRadiation.emissivity(patch().index())
189  );
190 
191  const scalarField& emissivity = temissivity();
192 
193  const tmp<scalarField> ttransmissivity
194  (
195  boundaryRadiation.transmissivity(patch().index())
196  );
197 
198  const scalarField& transmissivity = ttransmissivity();
199 
200  scalarField& qem = ray.qem().boundaryFieldRef()[patchi];
201  scalarField& qin = ray.qin().boundaryFieldRef()[patchi];
202 
203  const vector& myRayId = dom.IRay(rayId).d();
204 
205  scalarField Ir(patch().size(), Zero);
206  forAll(Iw, facei)
207  {
208  for (label rayi=0; rayi < dom.nRay(); rayi++)
209  {
210  const vector& d = dom.IRay(rayi).d();
211 
212  if ((-n[facei] & d) < 0.0)
213  {
214  // q into the wall
215  const scalarField& IFace =
216  dom.IRay(rayi).ILambda(lambdaId).boundaryField()[patchi];
217 
218  const vector& rayDave = dom.IRay(rayi).dAve();
219  Ir[facei] += IFace[facei]*(n[facei] & rayDave);
220  }
221  }
222  }
223 
224  if (dom.useSolarLoad())
225  {
226  // Looking for primary heat flux single band
227  Ir += patch().lookupPatchField<volScalarField,scalar>
228  (
229  dom.primaryFluxName_ + "_0"
230  );
231 
232  word qSecName = dom.relfectedFluxName_ + "_0";
233 
234  if (this->db().foundObject<volScalarField>(qSecName))
235  {
236  const volScalarField& qSec =
237  this->db().lookupObject<volScalarField>(qSecName);
238 
239  Ir += qSec.boundaryField()[patch().index()];
240  }
241  }
242 
243  scalarField Iexternal(this->size(), 0.0);
244 
245  if (dom.useExternalBeam())
246  {
247  const vector sunDir = dom.solarCalc().direction();
248  const scalar directSolarRad = dom.solarCalc().directSolarRad();
249 
250  //label nRaysBeam = dom.nRaysBeam();
251  label SunRayId(-1);
252  scalar maxSunRay = -GREAT;
253 
254  // Looking for the ray closest to the Sun direction
255  for (label rayI=0; rayI < dom.nRay(); rayI++)
256  {
257  const vector& iD = dom.IRay(rayI).d();
258  scalar dir = sunDir & iD;
259  if (dir > maxSunRay)
260  {
261  maxSunRay = dir;
262  SunRayId = rayI;
263  }
264  }
265 
266  if (rayId == SunRayId)
267  {
268  const scalarField nAve(n & dom.IRay(rayId).dAve());
269  forAll(Iexternal, faceI)
270  {
271  Iexternal[faceI] = directSolarRad/mag(dom.IRay(rayId).dAve());
272  }
273  }
274  }
275 
276  scalarField Isource(this->size(), 0.0);
277 
278  if (qRadExt_ > 0)
279  {
280  if (mag(qRadExtDir_) > 0)
281  {
282  label rayqoId = -1;
283  scalar maxRay = -GREAT;
284 
285  // Looking for the ray closest to the Sun direction
286  for (label rayI = 0; rayI < dom.nRay(); ++rayI)
287  {
288  const vector& iD = dom.IRay(rayI).d();
289  const scalar dir = qRadExtDir_ & iD;
290 
291  if (dir > maxRay)
292  {
293  maxRay = dir;
294  rayqoId = rayI;
295  }
296  }
297 
298  if (rayId == rayqoId)
299  {
300  forAll(Isource, faceI)
301  {
302  Isource[faceI] += qRadExt_/mag(dom.IRay(rayId).dAve());
303  }
304  }
305  }
306  else
307  {
308  forAll(Iw, faceI)
309  {
310  label rayqoId = -1;
311  scalar maxRay = -GREAT;
312 
313  // Looking for the ray closest to the Sun direction
314  for (label rayI = 0; rayI < dom.nRay(); ++rayI)
315  {
316  const vector& iD = dom.IRay(rayI).d();
317  const scalar dir = -n[faceI] & iD;
318 
319  if (dir > maxRay)
320  {
321  maxRay = dir;
322  rayqoId = rayI;
323  }
324  }
325 
326  if (rayId == rayqoId)
327  {
328  Isource[faceI] += qRadExt_/mag(dom.IRay(rayId).dAve());
329  }
330  }
331  }
332  }
333 
334  forAll(Iw, faceI)
335  {
336  if ((-n[faceI] & myRayId) > 0.0)
337  {
338  // direction out of the wall
339  refGrad()[faceI] = 0.0;
340  valueFraction()[faceI] = 1.0;
341  refValue()[faceI] =
342  Isource[faceI]
343  + Iexternal[faceI]*transmissivity[faceI]
344  + (
345  Ir[faceI]*(scalar(1) - emissivity[faceI])
346  + emissivity[faceI]*physicoChemical::sigma.value()
347  * pow4(Tp[faceI])
348  )/pi;
349 
350  // Emitted heat flux from this ray direction
351  qem[faceI] = refValue()[faceI]*nAve[faceI];
352  }
353  else
354  {
355  // direction into the wall
356  valueFraction()[faceI] = 0.0;
357  refGrad()[faceI] = 0.0;
358  refValue()[faceI] = 0.0; //not used
359 
360  // Incident heat flux on this ray direction
361  qin[faceI] = Iw[faceI]*nAve[faceI];
362  }
363  }
364 
365  // Restore tag
366  UPstream::msgType() = oldTag;
367 
368  mixedFvPatchScalarField::updateCoeffs();
369 }
370 
371 
373 (
374  Ostream& os
375 ) const
376 {
378  os.writeEntryIfDifferent<word>("T", "T", TName_);
379  os.writeEntryIfDifferent<scalar>("qRadExt", Zero, qRadExt_);
380  os.writeEntryIfDifferent<vector>("qRadExtDir", Zero, qRadExtDir_);
381 }
382 
383 
384 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
385 
386 namespace Foam
387 {
388 namespace radiation
389 {
391  (
394  );
395 }
396 }
397 
398 
399 // ************************************************************************* //
Different types of constants.
makePatchTypeField(fvPatchScalarField, greyDiffusiveRadiationMixedFvPatchScalarField)
dictionary dict
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
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:120
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:578
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.
friend Ostream & operator(Ostream &, const fvPatchField< scalar > &)
static const boundaryRadiationProperties & New(const fvMesh &mesh, Args &&... args)
Get existing or create a new MeshObject.
Definition: MeshObject.C:41
Unit conversion functions.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:68
static int & msgType() noexcept
Message tag of standard messages.
Definition: UPstream.H:806
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:413
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:84
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
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 INVALID.
Definition: exprTraits.C:52
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:55
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:352
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.
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:157