wideBandDiffusiveRadiationMixedFvPatchScalarField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2017 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  ),
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("value"))
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  const scalarField& Tp =
102  patch().lookupPatchField<volScalarField, scalar>(TName_);
103 
104  refValue() =
105  4.0*physicoChemical::sigma.value()*pow4(Tp)*emissivity()/pi;
106  refGrad() = 0.0;
107 
108  fvPatchScalarField::operator=(refValue());
109  }
110 }
111 
112 
115 (
117 )
118 :
119  mixedFvPatchScalarField(ptf),
121  (
122  ptf.patch(),
123  ptf.emissivityMethod(),
124  ptf.emissivity_
125  ),
126  TName_(ptf.TName_)
127 {}
128 
129 
132 (
135 )
136 :
137  mixedFvPatchScalarField(ptf, iF),
139  (
140  ptf.patch(),
141  ptf.emissivityMethod(),
142  ptf.emissivity_
143  ),
144  TName_(ptf.TName_)
145 {}
146 
147 
148 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
149 
152 {
153  if (this->updated())
154  {
155  return;
156  }
157 
158  // Since we're inside initEvaluate/evaluate there might be processor
159  // comms underway. Change the tag we use.
160  int oldTag = UPstream::msgType();
161  UPstream::msgType() = oldTag+1;
162 
163  const radiationModel& radiation =
164  db().lookupObject<radiationModel>("radiationProperties");
165 
166  const fvDOM& dom(refCast<const fvDOM>(radiation));
167 
168  label rayId = -1;
169  label lambdaId = -1;
170  dom.setRayIdLambdaId(internalField().name(), rayId, lambdaId);
171 
172  const label patchi = patch().index();
173 
174  if (dom.nLambda() == 0)
175  {
177  << " a non-grey boundary condition is used with a grey "
178  << "absorption model" << nl << exit(FatalError);
179  }
180 
181  scalarField& Iw = *this;
182  const vectorField n(patch().Sf()/patch().magSf());
183 
184  radiativeIntensityRay& ray =
185  const_cast<radiativeIntensityRay&>(dom.IRay(rayId));
186 
187  const scalarField nAve(n & ray.dAve());
188 
189  ray.qr().boundaryFieldRef()[patchi] += Iw*nAve;
190 
191  const scalarField Eb
192  (
193  dom.blackBody().bLambda(lambdaId).boundaryField()[patchi]
194  );
195 
196  scalarField temissivity = emissivity();
197 
198  scalarField& qem = ray.qem().boundaryFieldRef()[patchi];
199  scalarField& qin = ray.qin().boundaryFieldRef()[patchi];
200 
201  // Use updated Ir while iterating over rays
202  // avoids to used lagged qin
203  scalarField Ir = dom.IRay(0).qin().boundaryField()[patchi];
204 
205  for (label rayI=1; rayI < dom.nRay(); rayI++)
206  {
207  Ir += dom.IRay(rayI).qin().boundaryField()[patchi];
208  }
209 
210  forAll(Iw, facei)
211  {
212  const vector& d = dom.IRay(rayId).d();
213 
214  if ((-n[facei] & d) > 0.0)
215  {
216  // direction out of the wall
217  refGrad()[facei] = 0.0;
218  valueFraction()[facei] = 1.0;
219  refValue()[facei] =
220  (
221  Ir[facei]*(1.0 - temissivity[facei])
222  + temissivity[facei]*Eb[facei]
223  )/pi;
224 
225  // Emmited heat flux from this ray direction
226  qem[facei] = refValue()[facei]*nAve[facei];
227  }
228  else
229  {
230  // direction into the wall
231  valueFraction()[facei] = 0.0;
232  refGrad()[facei] = 0.0;
233  refValue()[facei] = 0.0; //not used
234 
235  // Incident heat flux on this ray direction
236  qin[facei] = Iw[facei]*nAve[facei];
237  }
238  }
239 
240  // Restore tag
241  UPstream::msgType() = oldTag;
242 
243  mixedFvPatchScalarField::updateCoeffs();
244 }
245 
246 
248 (
249  Ostream& os
250 ) const
251 {
254  writeEntryIfDifferent<word>(os, "T", "T", TName_);
255 }
256 
257 
258 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
259 
260 namespace Foam
261 {
262 namespace radiation
263 {
265  (
268  );
269 }
270 }
271 
272 
273 // ************************************************************************* //
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:465
scalarField emissivity_
Emissivity.
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
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:262
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
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.
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 emmited heat flux.
const vector & dAve() const
Return the average vector inside the solid angle.
Finite Volume Discrete Ordinates Method. Solves the RTE equation for n directions in a participating ...
Definition: fvDOM.H:80
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.