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-2018 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"
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 autoMap
153 (
154  const fvPatchFieldMapper& m
155 )
156 {
157  mixedFvPatchScalarField::autoMap(m);
159 }
160 
161 
163 (
164  const fvPatchScalarField& ptf,
165  const labelList& addr
166 )
167 {
168  mixedFvPatchScalarField::rmap(ptf, addr);
169  radiationCoupledBase::rmap(ptf, addr);
170 }
171 
174 {
175  if (this->updated())
176  {
177  return;
178  }
179 
180  // Since we're inside initEvaluate/evaluate there might be processor
181  // comms underway. Change the tag we use.
182  int oldTag = UPstream::msgType();
183  UPstream::msgType() = oldTag+1;
184 
185  const radiationModel& radiation =
186  db().lookupObject<radiationModel>("radiationProperties");
187 
188  const fvDOM& dom(refCast<const fvDOM>(radiation));
189 
190  label rayId = -1;
191  label lambdaId = -1;
192  dom.setRayIdLambdaId(internalField().name(), rayId, lambdaId);
193 
194  const label patchi = patch().index();
195 
196  if (dom.nLambda() == 0)
197  {
199  << " a non-grey boundary condition is used with a grey "
200  << "absorption model" << nl << exit(FatalError);
201  }
202 
203  scalarField& Iw = *this;
204  const vectorField n(patch().Sf()/patch().magSf());
205 
206  radiativeIntensityRay& ray =
207  const_cast<radiativeIntensityRay&>(dom.IRay(rayId));
208 
209  const scalarField nAve(n & ray.dAve());
210 
211  ray.qr().boundaryFieldRef()[patchi] += Iw*nAve;
212 
213  const scalarField Eb
214  (
215  dom.blackBody().bLambda(lambdaId).boundaryField()[patchi]
216  );
217 
218  scalarField temissivity = emissivity();
219 
220  scalarField& qem = ray.qem().boundaryFieldRef()[patchi];
221  scalarField& qin = ray.qin().boundaryFieldRef()[patchi];
222 
223  // Use updated Ir while iterating over rays
224  // avoids to used lagged qin
225  scalarField Ir = dom.IRay(0).qin().boundaryField()[patchi];
226 
227  for (label rayI=1; rayI < dom.nRay(); rayI++)
228  {
229  Ir += dom.IRay(rayI).qin().boundaryField()[patchi];
230  }
231 
232  forAll(Iw, facei)
233  {
234  const vector& d = dom.IRay(rayId).d();
235 
236  if ((-n[facei] & d) > 0.0)
237  {
238  // direction out of the wall
239  refGrad()[facei] = 0.0;
240  valueFraction()[facei] = 1.0;
241  refValue()[facei] =
242  (
243  Ir[facei]*(1.0 - temissivity[facei])
244  + temissivity[facei]*Eb[facei]
245  )/pi;
246 
247  // Emitted heat flux from this ray direction
248  qem[facei] = refValue()[facei]*nAve[facei];
249  }
250  else
251  {
252  // direction into the wall
253  valueFraction()[facei] = 0.0;
254  refGrad()[facei] = 0.0;
255  refValue()[facei] = 0.0; // not used
256 
257  // Incident heat flux on this ray direction
258  qin[facei] = Iw[facei]*nAve[facei];
259  }
260  }
261 
262  // Restore tag
263  UPstream::msgType() = oldTag;
264 
265  mixedFvPatchScalarField::updateCoeffs();
266 }
267 
268 
270 (
271  Ostream& os
272 ) const
273 {
276  writeEntryIfDifferent<word>(os, "T", "T", TName_);
277 }
278 
279 
280 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
281 
282 namespace Foam
283 {
284 namespace radiation
285 {
287  (
290  );
291 }
292 }
293 
294 
295 // ************************************************************************* //
Collection of constants.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:431
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
volScalarField & qin()
Return non-const access to the boundary incident heat flux.
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
friend Ostream & operator(Ostream &, const fvPatchField< Type > &)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
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:474
scalarField emissivity_
Emissivity.
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
void write(Ostream &, const label, const dictionary &)
Write with dictionary lookup.
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.
Top level model for radiation modelling.
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:161
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
static const char nl
Definition: Ostream.H:265
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Radiation intensity for a ray in a given direction.
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
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.
This boundary condition provides a wide-band, diffusive radiation condition, where the patch temperat...
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:395
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:100
label n
volScalarField & qem()
Return non-const access to the boundary emitted heat flux.
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
const vector & dAve() const
Return the average vector inside the solid angle.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Finite Volume Discrete Ordinates Method. Solves the RTE equation for n directions in a participating ...
Definition: fvDOM.H:83
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
const volScalarField & qr() const
Return const access to the boundary heat flux.
wideBandDiffusiveRadiationMixedFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
Namespace for OpenFOAM.