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-2021 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 (
60  const fvPatch& p,
62  const dictionary& dict
63 )
64 :
65  mixedFvPatchScalarField(p, iF),
66  radiationCoupledBase(p, dict),
67  TName_(dict.lookupOrDefault<word>("T", "T"))
68 {
69  if (dict.found("value"))
70  {
72  (
73  scalarField("value", dict, p.size())
74  );
75  refValue() = scalarField("refValue", dict, p.size());
76  refGrad() = scalarField("refGradient", dict, p.size());
77  valueFraction() = scalarField("valueFraction", dict, p.size());
78  }
79  else
80  {
81  const scalarField& Tp =
82  patch().lookupPatchField<volScalarField, scalar>(TName_);
83 
84  refValue() =
85  4.0*physicoChemical::sigma.value()*pow4(Tp)*emissivity()/pi;
86  refGrad() = 0.0;
87 
89  }
90 }
91 
92 
95 (
97  const fvPatch& p,
99  const fvPatchFieldMapper& mapper
100 )
101 :
102  mixedFvPatchScalarField(ptf, p, iF, mapper),
104  (
105  p,
106  ptf.emissivityMethod(),
107  ptf.emissivity_,
108  mapper
109  ),
110  TName_(ptf.TName_)
111 {}
112 
113 
116 (
119 )
120 :
121  mixedFvPatchScalarField(ptf, iF),
123  (
124  ptf.patch(),
125  ptf.emissivityMethod(),
126  ptf.emissivity_
127  ),
128  TName_(ptf.TName_)
129 {}
130 
131 
132 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
133 
135 (
136  const fvPatchFieldMapper& m
137 )
138 {
139  mixedFvPatchScalarField::autoMap(m);
141 }
142 
143 
145 (
146  const fvPatchScalarField& ptf,
147  const labelList& addr
148 )
149 {
150  mixedFvPatchScalarField::rmap(ptf, addr);
151  radiationCoupledBase::rmap(ptf, addr);
152 }
153 
154 
156 {
157  if (this->updated())
158  {
159  return;
160  }
161 
162  // Since we're inside initEvaluate/evaluate there might be processor
163  // comms underway. Change the tag we use.
164  int oldTag = UPstream::msgType();
165  UPstream::msgType() = oldTag+1;
166 
167  const radiationModels::fvDOM& dom =
168  db().lookupObject<radiationModels::fvDOM>("radiationProperties");
169 
170  label rayId = -1;
171  label lambdaId = -1;
172  dom.setRayIdLambdaId(internalField().name(), rayId, lambdaId);
173 
174  const label patchi = patch().index();
175 
176  if (dom.nLambda() == 0)
177  {
179  << " a non-grey boundary condition is used with a grey "
180  << "absorption model" << nl << exit(FatalError);
181  }
182 
183  scalarField& Iw = *this;
184  const vectorField n(patch().Sf()/patch().magSf());
185 
187  const_cast<radiationModels::radiativeIntensityRay&>(dom.IRay(rayId));
188 
189  const scalarField nAve(n & ray.dAve());
190 
191  ray.qr().boundaryFieldRef()[patchi] += Iw*nAve;
192 
193  const scalarField Eb
194  (
195  dom.blackBody().bLambda(lambdaId).boundaryField()[patchi]
196  );
197 
198  scalarField temissivity = emissivity();
199 
200  scalarField& qem = ray.qem().boundaryFieldRef()[patchi];
201  scalarField& qin = ray.qin().boundaryFieldRef()[patchi];
202 
203  // Use updated Ir while iterating over rays
204  // avoids to used lagged qin
205  scalarField Ir = dom.IRay(0).qin().boundaryField()[patchi];
206 
207  for (label rayI=1; rayI < dom.nRay(); rayI++)
208  {
209  Ir += dom.IRay(rayI).qin().boundaryField()[patchi];
210  }
211 
212  forAll(Iw, facei)
213  {
214  const vector& d = dom.IRay(rayId).d();
215 
216  if ((-n[facei] & d) > 0.0)
217  {
218  // direction out of the wall
219  refGrad()[facei] = 0.0;
220  valueFraction()[facei] = 1.0;
221  refValue()[facei] =
222  (
223  Ir[facei]*(1.0 - temissivity[facei])
224  + temissivity[facei]*Eb[facei]
225  )/pi;
226 
227  // Emitted heat flux from this ray direction
228  qem[facei] = refValue()[facei]*nAve[facei];
229  }
230  else
231  {
232  // direction into the wall
233  valueFraction()[facei] = 0.0;
234  refGrad()[facei] = 0.0;
235  refValue()[facei] = 0.0; // not used
236 
237  // Incident heat flux on this ray direction
238  qin[facei] = Iw[facei]*nAve[facei];
239  }
240  }
241 
242  // Restore tag
243  UPstream::msgType() = oldTag;
244 
245  mixedFvPatchScalarField::updateCoeffs();
246 }
247 
248 
250 (
251  Ostream& os
252 ) const
253 {
256  writeEntryIfDifferent<word>(os, "T", "T", TName_);
257 }
258 
259 
260 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
261 
262 namespace Foam
263 {
265  (
268  );
269 }
270 
271 
272 // ************************************************************************* //
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:643
#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:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:323
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.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:62
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
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
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:156
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:576
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:275
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...