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-2019 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"
33 
34 #include "fvDOM.H"
36 #include "constants.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 {
53  refValue() = Zero;
54  refGrad() = Zero;
55  valueFraction() = 1.0;
56 }
57 
58 
61 (
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;
91  valueFraction() = 1.0;
92 
94  }
95 }
96 
97 
100 (
102 )
103 :
104  mixedFvPatchScalarField(ptf)
105 {}
106 
107 
110 (
113 )
114 :
115  mixedFvPatchScalarField(ptf, iF)
116 {}
117 
118 
119 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
120 
122 updateCoeffs()
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 
133  const radiationModel& radiation =
134  db().lookupObject<radiationModel>("radiationProperties");
135 
136  const fvDOM& dom(refCast<const fvDOM>(radiation));
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 
154  radiativeIntensityRay& ray =
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  word qSecName = dom.relfectedFluxName_ + "_" + name(lambdaId);
216 
217  if (this->db().foundObject<volScalarField>(qSecName))
218  {
219  const volScalarField& qSec =
220  this->db().lookupObject<volScalarField>(qSecName);
221 
222  Ir += qSec.boundaryField()[patchi];
223  }
224  }
225 
226  scalarField Iexternal(this->size(), 0.0);
227 
228  if (dom.useExternalBeam())
229  {
230  const vector sunDir = dom.solarCalc().direction();
231  const scalar directSolarRad =
232  dom.solarCalc().directSolarRad()
233  *dom.spectralDistribution()[lambdaId];
234 
235  //label nRaysBeam = dom.nRaysBeam();
236  label SunRayId(-1);
237  scalar maxSunRay = -GREAT;
238 
239  // Looking for the ray closest to the Sun direction
240  for (label rayI=0; rayI < dom.nRay(); rayI++)
241  {
242  const vector& iD = dom.IRay(rayI).d();
243  scalar dir = sunDir & iD;
244  if (dir > maxSunRay)
245  {
246  maxSunRay = dir;
247  SunRayId = rayI;
248  }
249  }
250 
251  if (rayId == SunRayId)
252  {
253  const scalarField nAve(n & dom.IRay(rayId).dAve());
254  forAll(Iexternal, faceI)
255  {
256  Iexternal[faceI] = directSolarRad/mag(dom.IRay(rayId).dAve());
257  }
258  }
259  }
260 
261  forAll(Iw, facei)
262  {
263  const vector& d = dom.IRay(rayId).d();
264 
265  if ((-n[facei] & d) > 0.0)
266  {
267  // direction out of the wall
268  refGrad()[facei] = 0.0;
269  valueFraction()[facei] = 1.0;
270  refValue()[facei] =
271  Iexternal[facei]*transmissivity[facei]
272  + (
273  Ir[facei]*(1.0 - emissivity[facei])
274  + emissivity[facei]*Eb[facei]
275  )/pi;
276 
277  // Emitted heat flux from this ray direction (sum over lambdaId)
278  qem[facei] += refValue()[facei]*nAve[facei];
279  }
280  else
281  {
282  // direction into the wall
283  valueFraction()[facei] = 0.0;
284  refGrad()[facei] = 0.0;
285  refValue()[facei] = 0.0; //not used
286 
287  // Incident heat flux on this ray direction (sum over lambdaId)
288  qin[facei] += Iw[facei]*nAve[facei];
289  }
290  }
291 
292  UPstream::msgType(oldTag); // Restore tag
293 
294  mixedFvPatchScalarField::updateCoeffs();
295 }
296 
297 
299 (
300  Ostream& os
301 ) const
302 {
304 }
305 
306 
307 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
308 
309 namespace Foam
310 {
311 namespace radiation
312 {
314  (
317  );
318 }
319 }
320 
321 
322 // ************************************************************************* //
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:1274
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:608
static const boundaryRadiationProperties & New(const fvMesh &mesh, Args &&... args)
Get existing or create MeshObject registered with typeName.
Definition: MeshObject.C:53
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
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:1252
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:72
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.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
A FieldMapper for finite-volume patch fields.
constexpr scalar pi(M_PI)
Top level model for radiation modelling.
Vector< scalar > vector
Definition: vector.H:57
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...
Definition: Ostream.H:56
OBJstream os(runTime.globalPath()/outputName)
This boundary condition provides a wide-band, diffusive radiation condition, where the patch temperat...
virtual void operator=(const UList< scalar > &)
Definition: fvPatchField.C:391
A simple container of IOobject preferences. Can also be used for general handling of read/no-read/rea...
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
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.
wideBandDiffusiveRadiationMixedFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127