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-2015 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(dimensionedInternalField().name(), rayId, lambdaId);
170 
171  const label patchI = patch().index();
172 
173  if (dom.nLambda() != 1)
174  {
176  (
177  "Foam::radiation::"
178  "greyDiffusiveRadiationMixedFvPatchScalarField::updateCoeffs"
179  ) << " a grey boundary condition is used with a non-grey "
180  << "absorption model" << nl << exit(FatalError);
181  }
182 
183  scalarField& Iw = *this;
184  const vectorField n(patch().nf());
185 
186  radiativeIntensityRay& ray =
187  const_cast<radiativeIntensityRay&>(dom.IRay(rayId));
188 
189  const scalarField nAve(n & ray.dAve());
190 
191  ray.Qr().boundaryField()[patchI] += Iw*nAve;
192 
193  const scalarField temissivity = emissivity();
194 
195  scalarField& Qem = ray.Qem().boundaryField()[patchI];
196  scalarField& Qin = ray.Qin().boundaryField()[patchI];
197 
198  const vector& myRayId = dom.IRay(rayId).d();
199 
200  // Use updated Ir while iterating over rays
201  // avoids to used lagged Qin
202  scalarField Ir = dom.IRay(0).Qin().boundaryField()[patchI];
203 
204  for (label rayI=1; rayI < dom.nRay(); rayI++)
205  {
206  Ir += dom.IRay(rayI).Qin().boundaryField()[patchI];
207  }
208 
209  forAll(Iw, faceI)
210  {
211  if ((-n[faceI] & myRayId) > 0.0)
212  {
213  // direction out of the wall
214  refGrad()[faceI] = 0.0;
215  valueFraction()[faceI] = 1.0;
216  refValue()[faceI] =
217  (
218  Ir[faceI]*(scalar(1.0) - temissivity[faceI])
219  + temissivity[faceI]*physicoChemical::sigma.value()
220  * pow4(Tp[faceI])
221  )/pi;
222 
223  // Emmited heat flux from this ray direction
224  Qem[faceI] = refValue()[faceI]*nAve[faceI];
225  }
226  else
227  {
228  // direction into the wall
229  valueFraction()[faceI] = 0.0;
230  refGrad()[faceI] = 0.0;
231  refValue()[faceI] = 0.0; //not used
232 
233  // Incident heat flux on this ray direction
234  Qin[faceI] = Iw[faceI]*nAve[faceI];
235  }
236  }
237 
238  // Restore tag
239  UPstream::msgType() = oldTag;
240 
241  mixedFvPatchScalarField::updateCoeffs();
242 }
243 
244 
246 (
247  Ostream& os
248 ) const
249 {
252  writeEntryIfDifferent<word>(os, "T", "T", TName_);
253 }
254 
255 
256 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
257 
258 namespace Foam
259 {
260 namespace radiation
261 {
263  (
266  );
267 }
268 }
269 
270 
271 // ************************************************************************* //
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:306
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:387
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Top level model for radiation modelling.
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Collection of constants.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
volScalarField & Qem()
Return non-const access to the boundary emmited heat flux.
Finite Volume Discrete Ordinates Method. Solves the RTE equation for n directions in a participating ...
Definition: fvDOM.H:83
A class for handling words, derived from string.
Definition: word.H:59
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
Foam::fvPatchFieldMapper.
void write(Ostream &) const
Write.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
word emissivityMethod() const
Method to obtain emissivity.
Namespace for OpenFOAM.
runTime write()
greyDiffusiveRadiationMixedFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static const Field< scalar > & null()
Return a null field.
Definition: Field.H:100
label n
static const char nl
Definition: Ostream.H:260
const vector & dAve() const
Return the average vector inside the solid angle.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
#define forAll(list, i)
Definition: UList.H:421
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.
static int & msgType()
Message tag of standard messages.
Definition: UPstream.H:451
Macros for easy insertion into run-time selection tables.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:314
scalarField emissivity_
Emissivity.
rDeltaT dimensionedInternalField()
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
error FatalError
volScalarField & Qin()
Return non-const access to the boundary incident heat flux.
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
Radiation intensity for a ray in a given direction.
const volScalarField & Qr() const
Return const access to the boundary heat flux.
dimensionedScalar pow4(const dimensionedScalar &ds)
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
mathematical constants.
virtual label size() const
Return size.
Definition: fvPatch.H:161
friend Ostream & operator(Ostream &, const fvPatchField< Type > &)