wideBandDiffusiveRadiationMixedFvPatchScalarField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2023 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
28 #include "fvPatchFieldMapper.H"
29 #include "volFields.H"
30 
31 #include "fvDOM.H"
32 #include "wideBand.H"
33 #include "constants.H"
34 
35 using namespace Foam::constant;
36 using namespace Foam::constant::mathematical;
37 
38 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
39 
42 (
43  const fvPatch& p,
45  const dictionary& dict
46 )
47 :
48  mixedFvPatchScalarField(p, iF, dict, false),
50  TName_(dict.lookupOrDefault<word>("T", "T"))
51 {
52  if (dict.found("value"))
53  {
55  (
56  scalarField("value", dict, p.size())
57  );
58  refValue() = scalarField("refValue", dict, p.size());
59  refGrad() = scalarField("refGradient", dict, p.size());
60  valueFraction() = scalarField("valueFraction", dict, p.size());
61  }
62  else
63  {
64  const scalarField& Tp =
65  patch().lookupPatchField<volScalarField, scalar>(TName_);
66 
67  refValue() =
69  refGrad() = 0.0;
70 
72  }
73 }
74 
75 
78 (
80  const fvPatch& p,
82  const fvPatchFieldMapper& mapper
83 )
84 :
85  mixedFvPatchScalarField(ptf, p, iF, mapper),
87  (
88  p,
89  ptf.emissivityMethod(),
90  ptf.emissivity_,
91  mapper
92  ),
93  TName_(ptf.TName_)
94 {}
95 
96 
99 (
102 )
103 :
104  mixedFvPatchScalarField(ptf, iF),
106  (
107  ptf.patch(),
108  ptf.emissivityMethod(),
109  ptf.emissivity_
110  ),
111  TName_(ptf.TName_)
112 {}
113 
114 
115 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
116 
118 (
119  const fvPatchScalarField& ptf,
120  const fvPatchFieldMapper& mapper
121 )
122 {
123  mixedFvPatchScalarField::map(ptf, mapper);
124  radiationCoupledBase::map(ptf, mapper);
125 }
126 
127 
129 (
130  const fvPatchScalarField& ptf
131 )
132 {
133  mixedFvPatchScalarField::reset(ptf);
135 }
136 
137 
139 {
140  if (this->updated())
141  {
142  return;
143  }
144 
145  // Since we're inside initEvaluate/evaluate there might be processor
146  // comms underway. Change the tag we use.
147  int oldTag = UPstream::msgType();
148  UPstream::msgType() = oldTag+1;
149 
150  const radiationModels::fvDOM& dom =
151  db().lookupObject<radiationModels::fvDOM>("radiationProperties");
152 
153  label rayId = -1;
154  label lambdaId = -1;
155  dom.setRayIdLambdaId(internalField().name(), rayId, lambdaId);
156 
157  const label patchi = patch().index();
158 
159  if (dom.nLambda() == 0)
160  {
162  << " a non-grey boundary condition is used with a grey "
163  << "absorption model" << nl << exit(FatalError);
164  }
165 
166  scalarField& Iw = *this;
167  const vectorField n(patch().Sf()/patch().magSf());
168 
170  const_cast<radiationModels::radiativeIntensityRay&>(dom.IRay(rayId));
171 
172  const scalarField nAve(n & ray.dAve());
173 
174  ray.qr().boundaryFieldRef()[patchi] += Iw*nAve;
175 
176  const scalarField Eb
177  (
178  dom.blackBody().bLambda(lambdaId).boundaryField()[patchi]
179  );
180 
181  scalarField emissivity(this->emissivity());
182 
183  scalarField& qem = ray.qem().boundaryFieldRef()[patchi];
184  scalarField& qin = ray.qin().boundaryFieldRef()[patchi];
185 
186  // Use updated Ir while iterating over rays
187  // avoids to used lagged qin
188  scalarField Ir = dom.IRay(0).qin().boundaryField()[patchi];
189 
190  for (label rayI=1; rayI < dom.nRay(); rayI++)
191  {
192  Ir += dom.IRay(rayI).qin().boundaryField()[patchi];
193  }
194 
195  forAll(Iw, facei)
196  {
197  const vector& d = dom.IRay(rayId).d();
198 
199  if ((-n[facei] & d) > 0.0)
200  {
201  // direction out of the wall
202  refGrad()[facei] = 0.0;
203  valueFraction()[facei] = 1.0;
204  refValue()[facei] =
205  (
206  Ir[facei]*(1.0 - emissivity[facei])
207  + emissivity[facei]*Eb[facei]
208  )/pi;
209 
210  // Emitted heat flux from this ray direction
211  qem[facei] = refValue()[facei]*nAve[facei];
212  }
213  else
214  {
215  // direction into the wall
216  valueFraction()[facei] = 0.0;
217  refGrad()[facei] = 0.0;
218  refValue()[facei] = 0.0; // not used
219 
220  // Incident heat flux on this ray direction
221  qin[facei] = Iw[facei]*nAve[facei];
222  }
223  }
224 
225  // Restore tag
226  UPstream::msgType() = oldTag;
227 
228  mixedFvPatchScalarField::updateCoeffs();
229 }
230 
231 
233 (
234  Ostream& os
235 ) const
236 {
239  writeEntryIfDifferent<word>(os, "T", "T", TName_);
240 }
241 
242 
243 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
244 
245 namespace Foam
246 {
248  (
251  );
252 }
253 
254 
255 // ************************************************************************* //
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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...
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
static int & msgType()
Message tag of standard messages.
Definition: UPstream.H:476
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
const Type & value() const
Return const reference to value.
Foam::fvPatchFieldMapper.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:87
friend Ostream & operator(Ostream &, const fvPatchField< Type > &)
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:251
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
Common functions to emissivity. It gets supplied from lookup into a dictionary or calculated by the s...
void write(Ostream &) const
Write.
virtual void reset(const fvPatchScalarField &)
Reset the fvPatchField to the given fvPatchField.
tmp< scalarField > emissivity() const
Calculate corresponding emissivity field.
virtual void map(const fvPatchScalarField &, const fvPatchFieldMapper &)
Map the given fvPatchField onto this fvPatchField.
const volScalarField & bLambda(const label lambdaI) const
Black body spectrum.
Finite Volume Discrete Ordinates Method. Solves the RTE equation for n directions in a participating ...
Definition: fvDOM.H:77
const blackBodyEmission & blackBody() const
Const access to black body.
Definition: fvDOMI.H:107
label nRay() const
Number of rays.
Definition: fvDOMI.H:56
label nLambda() const
Number of wavelengths.
Definition: fvDOMI.H:62
void setRayIdLambdaId(const word &name, label &rayId, label &lambdaId) const
Set the rayId and lambdaId from by decomposing an intensity.
Definition: fvDOM.C:576
const radiativeIntensityRay & IRay(const label rayI) const
Ray intensity for rayI.
Definition: fvDOMI.H:27
Radiation intensity for a ray in a given direction.
const vector & dAve() const
Return the average vector inside the solid angle.
const vector & d() const
Return direction.
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...
virtual void reset(const fvPatchScalarField &)
Reset the fvPatchField to the given fvPatchField.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void map(const fvPatchScalarField &, const fvPatchFieldMapper &)
Map the given fvPatchField onto this fvPatchField.
wideBandDiffusiveRadiationMixedFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
label patchi
mathematical constants.
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
Collection of constants.
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61
dimensionedScalar pow4(const dimensionedScalar &ds)
error FatalError
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
static const char nl
Definition: Ostream.H:260
dictionary dict
volScalarField & p