greyDiffusiveRadiationMixedFvPatchScalarField.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"
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  ),
72  TName_(ptf.TName_)
73 {}
74 
75 
78 (
79  const fvPatch& p,
81  const dictionary& dict
82 )
83 :
84  mixedFvPatchScalarField(p, iF),
85  radiationCoupledBase(p, dict),
86  TName_(dict.lookupOrDefault<word>("T", "T"))
87 {
88  if (dict.found("refValue"))
89  {
91  (
92  scalarField("value", dict, p.size())
93  );
94  refValue() = scalarField("refValue", dict, p.size());
95  refGrad() = scalarField("refGradient", dict, p.size());
96  valueFraction() = scalarField("valueFraction", dict, p.size());
97  }
98  else
99  {
100  refValue() = 0.0;
101  refGrad() = 0.0;
102  valueFraction() = 1.0;
103 
104  fvPatchScalarField::operator=(refValue());
105  }
106 }
107 
108 
111 (
113 )
114 :
115  mixedFvPatchScalarField(ptf),
117  (
118  ptf.patch(),
119  ptf.emissivityMethod(),
120  ptf.emissivity_
121  ),
122  TName_(ptf.TName_)
123 {}
124 
125 
128 (
131 )
132 :
133  mixedFvPatchScalarField(ptf, iF),
135  (
136  ptf.patch(),
137  ptf.emissivityMethod(),
138  ptf.emissivity_
139  ),
140  TName_(ptf.TName_)
141 {}
142 
143 
144 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
145 
148 {
149  if (this->updated())
150  {
151  return;
152  }
153 
154  // Since we're inside initEvaluate/evaluate there might be processor
155  // comms underway. Change the tag we use.
156  int oldTag = UPstream::msgType();
157  UPstream::msgType() = oldTag+1;
158 
159  const scalarField& Tp =
160  patch().lookupPatchField<volScalarField, scalar>(TName_);
161 
162  const radiationModel& radiation =
163  db().lookupObject<radiationModel>("radiationProperties");
164 
165  const fvDOM& dom(refCast<const fvDOM>(radiation));
166 
167  label rayId = -1;
168  label lambdaId = -1;
169  dom.setRayIdLambdaId(internalField().name(), rayId, lambdaId);
170 
171  const label patchi = patch().index();
172 
173  if (dom.nLambda() != 1)
174  {
176  << " a grey boundary condition is used with a non-grey "
177  << "absorption model" << nl << exit(FatalError);
178  }
179 
180  scalarField& Iw = *this;
181  const vectorField n(patch().nf());
182 
183  radiativeIntensityRay& ray =
184  const_cast<radiativeIntensityRay&>(dom.IRay(rayId));
185 
186  const scalarField nAve(n & ray.dAve());
187 
188  ray.qr().boundaryFieldRef()[patchi] += Iw*nAve;
189 
190  const scalarField temissivity = emissivity();
191 
192  scalarField& qem = ray.qem().boundaryFieldRef()[patchi];
193  scalarField& qin = ray.qin().boundaryFieldRef()[patchi];
194 
195  const vector& myRayId = dom.IRay(rayId).d();
196 
197  // Use updated Ir while iterating over rays
198  // avoids to used lagged qin
199  scalarField Ir = dom.IRay(0).qin().boundaryField()[patchi];
200 
201  for (label rayI=1; rayI < dom.nRay(); rayI++)
202  {
203  Ir += dom.IRay(rayI).qin().boundaryField()[patchi];
204  }
205 
206  forAll(Iw, facei)
207  {
208  if ((-n[facei] & myRayId) > 0.0)
209  {
210  // direction out of the wall
211  refGrad()[facei] = 0.0;
212  valueFraction()[facei] = 1.0;
213  refValue()[facei] =
214  (
215  Ir[facei]*(scalar(1) - temissivity[facei])
216  + temissivity[facei]*physicoChemical::sigma.value()
217  * pow4(Tp[facei])
218  )/pi;
219 
220  // Emmited heat flux from this ray direction
221  qem[facei] = refValue()[facei]*nAve[facei];
222  }
223  else
224  {
225  // direction into the wall
226  valueFraction()[facei] = 0.0;
227  refGrad()[facei] = 0.0;
228  refValue()[facei] = 0.0; //not used
229 
230  // Incident heat flux on this ray direction
231  qin[facei] = Iw[facei]*nAve[facei];
232  }
233  }
234 
235  // Restore tag
236  UPstream::msgType() = oldTag;
237 
238  mixedFvPatchScalarField::updateCoeffs();
239 }
240 
241 
243 (
244  Ostream& os
245 ) const
246 {
249  writeEntryIfDifferent<word>(os, "T", "T", TName_);
250 }
251 
252 
253 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
254 
255 namespace Foam
256 {
257 namespace radiation
258 {
260  (
263  );
264 }
265 }
266 
267 
268 // ************************************************************************* //
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
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.
friend Ostream & operator(Ostream &, const fvPatchField< Type > &)
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
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.
Top level model for radiation modelling.
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
greyDiffusiveRadiationMixedFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
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.
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.
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
Namespace for OpenFOAM.