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-2024 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 "fieldMapper.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  {
54  refValue() = scalarField("refValue", iF.dimensions(), dict, p.size());
55  refGrad() =
57  (
58  "refGradient",
59  iF.dimensions()/dimLength,
60  dict,
61  p.size()
62  );
63  valueFraction() =
64  scalarField("valueFraction", unitFraction, dict, p.size());
65 
67  (
68  scalarField("value", iF.dimensions(), dict, p.size())
69  );
70  }
71  else
72  {
73  const scalarField& Tp =
74  patch().lookupPatchField<volScalarField, scalar>(TName_);
75 
76  refValue() =
78  refGrad() = 0.0;
79 
81  }
82 }
83 
84 
87 (
89  const fvPatch& p,
91  const fieldMapper& mapper
92 )
93 :
94  mixedFvPatchScalarField(ptf, p, iF, mapper),
96  (
97  p,
98  ptf.emissivityMethod(),
99  ptf.emissivity_,
100  mapper
101  ),
102  TName_(ptf.TName_)
103 {}
104 
105 
108 (
111 )
112 :
113  mixedFvPatchScalarField(ptf, iF),
115  (
116  ptf.patch(),
117  ptf.emissivityMethod(),
118  ptf.emissivity_
119  ),
120  TName_(ptf.TName_)
121 {}
122 
123 
124 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
125 
127 (
128  const fvPatchScalarField& ptf,
129  const fieldMapper& mapper
130 )
131 {
132  mixedFvPatchScalarField::map(ptf, mapper);
133  radiationCoupledBase::map(ptf, mapper);
134 }
135 
136 
138 (
139  const fvPatchScalarField& ptf
140 )
141 {
142  mixedFvPatchScalarField::reset(ptf);
144 }
145 
146 
148 {
149  if (this->updated())
150  {
151  return;
152  }
153 
154  // Since we're inside initEvaluate/evaluate there might be processor
155  // comms underway. Change the tag we use.
156  int oldTag = UPstream::msgType();
157  UPstream::msgType() = oldTag+1;
158 
159  const radiationModels::fvDOM& dom =
160  db().lookupObject<radiationModels::fvDOM>("radiationProperties");
161 
162  label rayId = -1;
163  label lambdaId = -1;
164  dom.setRayIdLambdaId(internalField().name(), rayId, lambdaId);
165 
166  const label patchi = patch().index();
167 
168  if (dom.nLambda() == 0)
169  {
171  << " a non-grey boundary condition is used with a grey "
172  << "absorption model" << nl << exit(FatalError);
173  }
174 
175  scalarField& Iw = *this;
176  const vectorField n(patch().Sf()/patch().magSf());
177 
179  const_cast<radiationModels::radiativeIntensityRay&>(dom.IRay(rayId));
180 
181  const scalarField nAve(n & ray.dAve());
182 
183  ray.qr().boundaryFieldRef()[patchi] += Iw*nAve;
184 
185  const scalarField Eb
186  (
187  dom.blackBody().bLambda(lambdaId).boundaryField()[patchi]
188  );
189 
190  scalarField emissivity(this->emissivity());
191 
192  scalarField& qem = ray.qem().boundaryFieldRef()[patchi];
193  scalarField& qin = ray.qin().boundaryFieldRef()[patchi];
194 
195  // Use updated Ir while iterating over rays
196  // avoids to used lagged qin
197  scalarField Ir = dom.IRay(0).qin().boundaryField()[patchi];
198 
199  for (label rayI=1; rayI < dom.nRay(); rayI++)
200  {
201  Ir += dom.IRay(rayI).qin().boundaryField()[patchi];
202  }
203 
204  forAll(Iw, facei)
205  {
206  const vector& d = dom.IRay(rayId).d();
207 
208  if ((-n[facei] & d) > 0.0)
209  {
210  // direction out of the wall
211  refGrad()[facei] = 0.0;
212  valueFraction()[facei] = 1.0;
213  refValue()[facei] =
214  (
215  Ir[facei]*(1.0 - emissivity[facei])
216  + emissivity[facei]*Eb[facei]
217  )/pi;
218 
219  // Emitted heat flux from this ray direction
220  qem[facei] = refValue()[facei]*nAve[facei];
221  }
222  else
223  {
224  // direction into the wall
225  valueFraction()[facei] = 0.0;
226  refGrad()[facei] = 0.0;
227  refValue()[facei] = 0.0; // not used
228 
229  // Incident heat flux on this ray direction
230  qin[facei] = Iw[facei]*nAve[facei];
231  }
232  }
233 
234  // Restore tag
235  UPstream::msgType() = oldTag;
236 
237  mixedFvPatchScalarField::updateCoeffs();
238 }
239 
240 
242 (
243  Ostream& os
244 ) const
245 {
248  writeEntryIfDifferent<word>(os, "T", "T", TName_);
249 }
250 
251 
252 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
253 
254 namespace Foam
255 {
257  (
260  );
261 }
262 
263 
264 // ************************************************************************* //
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 dimensionSet & dimensions() const
Return dimensions.
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:162
const Type & value() const
Return const reference to value.
Abstract base class for field mapping.
Definition: fieldMapper.H:48
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:88
friend Ostream & operator(Ostream &, const fvPatchField< Type > &)
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:249
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 fieldMapper &)
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:577
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 fieldMapper &)
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:334
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
const dimensionSet dimLength
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:64
dimensionedScalar pow4(const dimensionedScalar &ds)
error FatalError
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
static const char nl
Definition: Ostream.H:266
const unitConversion unitFraction
dictionary dict
volScalarField & p