greyDiffusiveRadiationMixedFvPatchScalarField.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 "constants.H"
33 
34 using namespace Foam::constant;
35 using namespace Foam::constant::mathematical;
36 
37 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38 
41 (
42  const fvPatch& p,
44 )
45 :
46  mixedFvPatchScalarField(p, iF),
47  radiationCoupledBase(p, "undefined", scalarField::null()),
48  TName_("T")
49 {
50  refValue() = 0.0;
51  refGrad() = 0.0;
52  valueFraction() = 1.0;
53 }
54 
55 
58 (
60  const fvPatch& p,
62  const fvPatchFieldMapper& mapper
63 )
64 :
65  mixedFvPatchScalarField(ptf, p, iF, mapper),
67  (
68  p,
69  ptf.emissivityMethod(),
70  ptf.emissivity_,
71  mapper
72  ),
73  TName_(ptf.TName_)
74 {}
75 
76 
79 (
80  const fvPatch& p,
82  const dictionary& dict
83 )
84 :
85  mixedFvPatchScalarField(p, iF),
86  radiationCoupledBase(p, dict),
87  TName_(dict.lookupOrDefault<word>("T", "T"))
88 {
89  if (dict.found("refValue"))
90  {
92  (
93  scalarField("value", dict, p.size())
94  );
95  refValue() = scalarField("refValue", dict, p.size());
96  refGrad() = scalarField("refGradient", dict, p.size());
97  valueFraction() = scalarField("valueFraction", dict, p.size());
98  }
99  else
100  {
101  refValue() = 0.0;
102  refGrad() = 0.0;
103  valueFraction() = 1.0;
104 
105  fvPatchScalarField::operator=(refValue());
106  }
107 }
108 
109 
112 (
114 )
115 :
116  mixedFvPatchScalarField(ptf),
118  (
119  ptf.patch(),
120  ptf.emissivityMethod(),
121  ptf.emissivity_
122  ),
123  TName_(ptf.TName_)
124 {}
125 
126 
129 (
132 )
133 :
134  mixedFvPatchScalarField(ptf, iF),
136  (
137  ptf.patch(),
138  ptf.emissivityMethod(),
139  ptf.emissivity_
140  ),
141  TName_(ptf.TName_)
142 {}
143 
144 
145 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
146 
148 (
149  const fvPatchFieldMapper& m
150 )
151 {
152  mixedFvPatchScalarField::autoMap(m);
154 }
155 
156 
158 (
159  const fvPatchScalarField& ptf,
160  const labelList& addr
161 )
162 {
163  mixedFvPatchScalarField::rmap(ptf, addr);
164  radiationCoupledBase::rmap(ptf, addr);
165 }
166 
167 
169 {
170  if (this->updated())
171  {
172  return;
173  }
174 
175  // Since we're inside initEvaluate/evaluate there might be processor
176  // comms underway. Change the tag we use.
177  int oldTag = UPstream::msgType();
178  UPstream::msgType() = oldTag+1;
179 
180  const scalarField& Tp =
181  patch().lookupPatchField<volScalarField, scalar>(TName_);
182 
183  const radiationModels::fvDOM& dom =
184  db().lookupObject<radiationModels::fvDOM>("radiationProperties");
185 
186  label rayId = -1;
187  label lambdaId = -1;
188  dom.setRayIdLambdaId(internalField().name(), rayId, lambdaId);
189 
190  const label patchi = patch().index();
191 
192  if (dom.nLambda() != 1)
193  {
195  << " a grey boundary condition is used with a non-grey "
196  << "absorption model" << nl << exit(FatalError);
197  }
198 
199  scalarField& Iw = *this;
200  const vectorField n(patch().nf());
201 
203  const_cast<radiationModels::radiativeIntensityRay&>(dom.IRay(rayId));
204 
205  const scalarField nAve(n & ray.dAve());
206 
207  ray.qr().boundaryFieldRef()[patchi] += Iw*nAve;
208 
209  const scalarField temissivity = emissivity();
210 
211  scalarField& qem = ray.qem().boundaryFieldRef()[patchi];
212  scalarField& qin = ray.qin().boundaryFieldRef()[patchi];
213 
214  const vector& myRayId = dom.IRay(rayId).dAve();
215 
216  // Use updated Ir while iterating over rays
217  // avoids to used lagged qin
218  scalarField Ir = dom.IRay(0).qin().boundaryField()[patchi];
219 
220  for (label rayI=1; rayI < dom.nRay(); rayI++)
221  {
222  Ir += dom.IRay(rayI).qin().boundaryField()[patchi];
223  }
224 
225  forAll(Iw, facei)
226  {
227  if ((-n[facei] & myRayId) > 0.0)
228  {
229  // direction out of the wall
230  refGrad()[facei] = 0.0;
231  valueFraction()[facei] = 1.0;
232  refValue()[facei] =
233  (
234  Ir[facei]*(scalar(1) - temissivity[facei])
235  + temissivity[facei]*physicoChemical::sigma.value()
236  * pow4(Tp[facei])
237  )/pi;
238 
239  // Emitted heat flux from this ray direction
240  qem[facei] = refValue()[facei]*nAve[facei];
241  }
242  else
243  {
244  // direction into the wall
245  valueFraction()[facei] = 0.0;
246  refGrad()[facei] = 0.0;
247  refValue()[facei] = 0.0; // not used
248 
249  // Incident heat flux on this ray direction
250  qin[facei] = Iw[facei]*nAve[facei];
251  }
252  }
253 
254  // Restore tag
255  UPstream::msgType() = oldTag;
256 
257  mixedFvPatchScalarField::updateCoeffs();
258 }
259 
260 
262 (
263  Ostream& os
264 ) const
265 {
268  writeEntryIfDifferent<word>(os, "T", "T", TName_);
269 }
270 
271 
272 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
273 
274 namespace Foam
275 {
277  (
280  );
281 }
282 
283 
284 // ************************************************************************* //
Collection of constants.
Finite Volume Discrete Ordinates Method. Solves the RTE equation for n directions in a participating ...
Definition: fvDOM.H:83
This boundary condition provides a grey-diffuse condition for radiation intensity, I, for use with the finite-volume discrete-ordinates model (fvDOM), in which the radiation temperature is retrieved from the temperature field boundary condition.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:438
#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
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/m^2/K^4].
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
greyDiffusiveRadiationMixedFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
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.
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
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.
volScalarField & qem()
Return non-const access to the boundary emitted heat flux.
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:53
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
static const char nl
Definition: Ostream.H:265
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.
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
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Namespace for OpenFOAM.