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-2019 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 )
46 :
47  mixedFvPatchScalarField(p, iF),
48  radiationCoupledBase(p, "undefined", scalarField::null()),
49  TName_("T")
50 {
51  refValue() = 0.0;
52  refGrad() = 0.0;
53  valueFraction() = 1.0;
54 }
55 
56 
59 (
61  const fvPatch& p,
63  const fvPatchFieldMapper& mapper
64 )
65 :
66  mixedFvPatchScalarField(ptf, p, iF, mapper),
68  (
69  p,
70  ptf.emissivityMethod(),
71  ptf.emissivity_,
72  mapper
73  ),
74  TName_(ptf.TName_)
75 {}
76 
77 
80 (
81  const fvPatch& p,
83  const dictionary& dict
84 )
85 :
86  mixedFvPatchScalarField(p, iF),
87  radiationCoupledBase(p, dict),
88  TName_(dict.lookupOrDefault<word>("T", "T"))
89 {
90  if (dict.found("value"))
91  {
93  (
94  scalarField("value", dict, p.size())
95  );
96  refValue() = scalarField("refValue", dict, p.size());
97  refGrad() = scalarField("refGradient", dict, p.size());
98  valueFraction() = scalarField("valueFraction", dict, p.size());
99  }
100  else
101  {
102  const scalarField& Tp =
103  patch().lookupPatchField<volScalarField, scalar>(TName_);
104 
105  refValue() =
106  4.0*physicoChemical::sigma.value()*pow4(Tp)*emissivity()/pi;
107  refGrad() = 0.0;
108 
109  fvPatchScalarField::operator=(refValue());
110  }
111 }
112 
113 
116 (
118 )
119 :
120  mixedFvPatchScalarField(ptf),
122  (
123  ptf.patch(),
124  ptf.emissivityMethod(),
125  ptf.emissivity_
126  ),
127  TName_(ptf.TName_)
128 {}
129 
130 
133 (
136 )
137 :
138  mixedFvPatchScalarField(ptf, iF),
140  (
141  ptf.patch(),
142  ptf.emissivityMethod(),
143  ptf.emissivity_
144  ),
145  TName_(ptf.TName_)
146 {}
147 
148 
149 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
150 
152 (
153  const fvPatchFieldMapper& m
154 )
155 {
156  mixedFvPatchScalarField::autoMap(m);
158 }
159 
160 
162 (
163  const fvPatchScalarField& ptf,
164  const labelList& addr
165 )
166 {
167  mixedFvPatchScalarField::rmap(ptf, addr);
168  radiationCoupledBase::rmap(ptf, addr);
169 }
170 
171 
173 {
174  if (this->updated())
175  {
176  return;
177  }
178 
179  // Since we're inside initEvaluate/evaluate there might be processor
180  // comms underway. Change the tag we use.
181  int oldTag = UPstream::msgType();
182  UPstream::msgType() = oldTag+1;
183 
184  const radiationModels::fvDOM& dom =
185  db().lookupObject<radiationModels::fvDOM>("radiationProperties");
186 
187  label rayId = -1;
188  label lambdaId = -1;
189  dom.setRayIdLambdaId(internalField().name(), rayId, lambdaId);
190 
191  const label patchi = patch().index();
192 
193  if (dom.nLambda() == 0)
194  {
196  << " a non-grey boundary condition is used with a grey "
197  << "absorption model" << nl << exit(FatalError);
198  }
199 
200  scalarField& Iw = *this;
201  const vectorField n(patch().Sf()/patch().magSf());
202 
204  const_cast<radiationModels::radiativeIntensityRay&>(dom.IRay(rayId));
205 
206  const scalarField nAve(n & ray.dAve());
207 
208  ray.qr().boundaryFieldRef()[patchi] += Iw*nAve;
209 
210  const scalarField Eb
211  (
212  dom.blackBody().bLambda(lambdaId).boundaryField()[patchi]
213  );
214 
215  scalarField temissivity = emissivity();
216 
217  scalarField& qem = ray.qem().boundaryFieldRef()[patchi];
218  scalarField& qin = ray.qin().boundaryFieldRef()[patchi];
219 
220  // Use updated Ir while iterating over rays
221  // avoids to used lagged qin
222  scalarField Ir = dom.IRay(0).qin().boundaryField()[patchi];
223 
224  for (label rayI=1; rayI < dom.nRay(); rayI++)
225  {
226  Ir += dom.IRay(rayI).qin().boundaryField()[patchi];
227  }
228 
229  forAll(Iw, facei)
230  {
231  const vector& d = dom.IRay(rayId).d();
232 
233  if ((-n[facei] & d) > 0.0)
234  {
235  // direction out of the wall
236  refGrad()[facei] = 0.0;
237  valueFraction()[facei] = 1.0;
238  refValue()[facei] =
239  (
240  Ir[facei]*(1.0 - temissivity[facei])
241  + temissivity[facei]*Eb[facei]
242  )/pi;
243 
244  // Emitted heat flux from this ray direction
245  qem[facei] = refValue()[facei]*nAve[facei];
246  }
247  else
248  {
249  // direction into the wall
250  valueFraction()[facei] = 0.0;
251  refGrad()[facei] = 0.0;
252  refValue()[facei] = 0.0; // not used
253 
254  // Incident heat flux on this ray direction
255  qin[facei] = Iw[facei]*nAve[facei];
256  }
257  }
258 
259  // Restore tag
260  UPstream::msgType() = oldTag;
261 
262  mixedFvPatchScalarField::updateCoeffs();
263 }
264 
265 
267 (
268  Ostream& os
269 ) const
270 {
273  writeEntryIfDifferent<word>(os, "T", "T", TName_);
274 }
275 
276 
277 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
278 
279 namespace Foam
280 {
282  (
285  );
286 }
287 
288 
289 // ************************************************************************* //
label nLambda() const
Number of wavelengths.
Definition: fvDOM.H:63
Collection of constants.
Finite Volume Discrete Ordinates Method. Solves the RTE equation for n directions in a participating ...
Definition: fvDOM.H:74
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:667
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
const volScalarField & bLambda(const label lambdaI) const
Black body spectrum.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
friend Ostream & operator(Ostream &, const fvPatchField< Type > &)
const vector & d() const
Return direction.
const dimensionedScalar & sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
static int & msgType()
Message tag of standard messages.
Definition: UPstream.H:476
const volScalarField & qr() const
Return const access to the boundary heat flux.
scalarField emissivity_
Emissivity.
Macros for easy insertion into run-time selection tables.
const vector & dAve() const
Return the average vector inside the solid angle.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
const radiativeIntensityRay & IRay(const label rayI) const
Ray intensity for rayI.
Definition: fvDOM.H:28
mathematical constants.
A class for handling words, derived from string.
Definition: word.H:59
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
void write(Ostream &) const
Write.
Foam::fvPatchFieldMapper.
const Type & value() const
Return const reference to value.
volScalarField & qem()
Return non-const access to the boundary emitted heat flux.
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Common functions to emissivity. It gets supplied from lookup into a dictionary or calculated by the s...
word emissivityMethod() const
Method to obtain emissivity.
virtual label size() const
Return size.
Definition: fvPatch.H:155
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
static const char nl
Definition: Ostream.H:260
void setRayIdLambdaId(const word &name, label &rayId, label &lambdaId) const
Set the rayId and lambdaId from by decomposing an intensity.
Definition: fvDOM.C:552
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
wideBandDiffusiveRadiationMixedFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
volScalarField & qin()
Return non-const access to the boundary incident heat flux.
Radiation intensity for a ray in a given direction.
label patchi
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:295
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
dimensionedScalar pow4(const dimensionedScalar &ds)
static const Field< scalar > & null()
Return a null field.
Definition: Field.H:111
label n
label nRay() const
Number of rays.
Definition: fvDOM.H:57
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
const blackBodyEmission & blackBody() const
Const access to black body.
Definition: fvDOM.H:108
Namespace for OpenFOAM.
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
This boundary condition provides a wide-band, diffusive radiation condition, where the patch temperat...