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-2021 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 (
59  const fvPatch& p,
61  const dictionary& dict
62 )
63 :
64  mixedFvPatchScalarField(p, iF),
65  radiationCoupledBase(p, dict),
66  TName_(dict.lookupOrDefault<word>("T", "T"))
67 {
68  if (dict.found("refValue"))
69  {
71  (
72  scalarField("value", dict, p.size())
73  );
74  refValue() = scalarField("refValue", dict, p.size());
75  refGrad() = scalarField("refGradient", dict, p.size());
76  valueFraction() = scalarField("valueFraction", dict, p.size());
77  }
78  else
79  {
80  refValue() = 0.0;
81  refGrad() = 0.0;
82  valueFraction() = 1.0;
83 
85  }
86 }
87 
88 
91 (
93  const fvPatch& p,
95  const fvPatchFieldMapper& mapper
96 )
97 :
98  mixedFvPatchScalarField(ptf, p, iF, mapper),
100  (
101  p,
102  ptf.emissivityMethod(),
103  ptf.emissivity_,
104  mapper
105  ),
106  TName_(ptf.TName_)
107 {}
108 
109 
112 (
115 )
116 :
117  mixedFvPatchScalarField(ptf, iF),
119  (
120  ptf.patch(),
121  ptf.emissivityMethod(),
122  ptf.emissivity_
123  ),
124  TName_(ptf.TName_)
125 {}
126 
127 
128 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
129 
131 (
132  const fvPatchFieldMapper& m
133 )
134 {
135  mixedFvPatchScalarField::autoMap(m);
137 }
138 
139 
141 (
142  const fvPatchScalarField& ptf,
143  const labelList& addr
144 )
145 {
146  mixedFvPatchScalarField::rmap(ptf, addr);
147  radiationCoupledBase::rmap(ptf, addr);
148 }
149 
150 
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 scalarField& Tp =
164  patch().lookupPatchField<volScalarField, scalar>(TName_);
165 
166  const radiationModels::fvDOM& dom =
167  db().lookupObject<radiationModels::fvDOM>("radiationProperties");
168 
169  label rayId = -1;
170  label lambdaId = -1;
171  dom.setRayIdLambdaId(internalField().name(), rayId, lambdaId);
172 
173  const label patchi = patch().index();
174 
175  if (dom.nLambda() != 1)
176  {
178  << " a grey boundary condition is used with a non-grey "
179  << "absorption model" << nl << exit(FatalError);
180  }
181 
182  scalarField& Iw = *this;
183  const vectorField n(patch().nf());
184 
186  const_cast<radiationModels::radiativeIntensityRay&>(dom.IRay(rayId));
187 
188  const scalarField nAve(n & ray.dAve());
189 
190  ray.qr().boundaryFieldRef()[patchi] += Iw*nAve;
191 
192  const scalarField temissivity = emissivity();
193 
194  scalarField& qem = ray.qem().boundaryFieldRef()[patchi];
195  scalarField& qin = ray.qin().boundaryFieldRef()[patchi];
196 
197  const vector& myRayId = dom.IRay(rayId).dAve();
198 
199  // Use updated Ir while iterating over rays
200  // avoids to used lagged qin
201  scalarField Ir = dom.IRay(0).qin().boundaryField()[patchi];
202 
203  for (label rayI=1; rayI < dom.nRay(); rayI++)
204  {
205  Ir += dom.IRay(rayI).qin().boundaryField()[patchi];
206  }
207 
208  forAll(Iw, facei)
209  {
210  if ((-n[facei] & myRayId) > 0.0)
211  {
212  // direction out of the wall
213  refGrad()[facei] = 0.0;
214  valueFraction()[facei] = 1.0;
215  refValue()[facei] =
216  (
217  Ir[facei]*(scalar(1) - temissivity[facei])
218  + temissivity[facei]*physicoChemical::sigma.value()
219  * pow4(Tp[facei])
220  )/pi;
221 
222  // Emitted heat flux from this ray direction
223  qem[facei] = refValue()[facei]*nAve[facei];
224  }
225  else
226  {
227  // direction into the wall
228  valueFraction()[facei] = 0.0;
229  refGrad()[facei] = 0.0;
230  refValue()[facei] = 0.0; // not used
231 
232  // Incident heat flux on this ray direction
233  qin[facei] = Iw[facei]*nAve[facei];
234  }
235  }
236 
237  // Restore tag
238  UPstream::msgType() = oldTag;
239 
240  mixedFvPatchScalarField::updateCoeffs();
241 }
242 
243 
245 (
246  Ostream& os
247 ) const
248 {
251  writeEntryIfDifferent<word>(os, "T", "T", TName_);
252 }
253 
254 
255 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
256 
257 namespace Foam
258 {
260  (
263  );
264 }
265 
266 
267 // ************************************************************************* //
Collection of constants.
Finite Volume Discrete Ordinates Method. Solves the RTE equation for n directions in a participating ...
Definition: fvDOM.H:74
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:643
#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:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:323
friend Ostream & operator(Ostream &, const fvPatchField< Type > &)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:62
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:57
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
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.
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
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:156
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
static const char nl
Definition: Ostream.H:260
void setRayIdLambdaId(const word &name, label &rayId, label &lambdaId) const
Set the rayId and lambdaId from by decomposing an intensity.
Definition: fvDOM.C:576
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:275
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.