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-2023 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  const dictionary& dict
45 )
46 :
47  mixedFvPatchScalarField(p, iF, dict, false),
49  TName_(dict.lookupOrDefault<word>("T", "T"))
50 {
51  if (dict.found("refValue"))
52  {
54  (
55  scalarField("value", dict, p.size())
56  );
57  refValue() = scalarField("refValue", dict, p.size());
58  refGrad() = scalarField("refGradient", dict, p.size());
59  valueFraction() = scalarField("valueFraction", dict, p.size());
60  }
61  else
62  {
63  refValue() = 0.0;
64  refGrad() = 0.0;
65  valueFraction() = 1.0;
66 
68  }
69 }
70 
71 
74 (
76  const fvPatch& p,
78  const fvPatchFieldMapper& mapper
79 )
80 :
81  mixedFvPatchScalarField(ptf, p, iF, mapper),
83  (
84  p,
85  ptf.emissivityMethod(),
86  ptf.emissivity_,
87  mapper
88  ),
89  TName_(ptf.TName_)
90 {}
91 
92 
95 (
98 )
99 :
100  mixedFvPatchScalarField(ptf, iF),
102  (
103  ptf.patch(),
104  ptf.emissivityMethod(),
105  ptf.emissivity_
106  ),
107  TName_(ptf.TName_)
108 {}
109 
110 
111 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
112 
114 (
115  const fvPatchScalarField& ptf,
116  const fvPatchFieldMapper& mapper
117 )
118 {
119  mixedFvPatchScalarField::map(ptf, mapper);
120  radiationCoupledBase::map(ptf, mapper);
121 }
122 
123 
125 (
126  const fvPatchScalarField& ptf
127 )
128 {
129  mixedFvPatchScalarField::reset(ptf);
131 }
132 
133 
135 {
136  if (this->updated())
137  {
138  return;
139  }
140 
141  // Since we're inside initEvaluate/evaluate there might be processor
142  // comms underway. Change the tag we use.
143  int oldTag = UPstream::msgType();
144  UPstream::msgType() = oldTag+1;
145 
146  const scalarField& Tp =
147  patch().lookupPatchField<volScalarField, scalar>(TName_);
148 
149  const radiationModels::fvDOM& dom =
150  db().lookupObject<radiationModels::fvDOM>("radiationProperties");
151 
152  label rayId = -1;
153  label lambdaId = -1;
154  dom.setRayIdLambdaId(internalField().name(), rayId, lambdaId);
155 
156  const label patchi = patch().index();
157 
158  if (dom.nLambda() != 1)
159  {
161  << " a grey boundary condition is used with a non-grey "
162  << "absorption model" << nl << exit(FatalError);
163  }
164 
165  scalarField& Iw = *this;
166  const vectorField n(patch().nf());
167 
169  const_cast<radiationModels::radiativeIntensityRay&>(dom.IRay(rayId));
170 
171  const scalarField nAve(n & ray.dAve());
172 
173  ray.qr().boundaryFieldRef()[patchi] += Iw*nAve;
174 
175  const scalarField emissivity(this->emissivity());
176 
177  scalarField& qem = ray.qem().boundaryFieldRef()[patchi];
178  scalarField& qin = ray.qin().boundaryFieldRef()[patchi];
179 
180  const vector& myRayId = dom.IRay(rayId).dAve();
181 
182  // Use updated Ir while iterating over rays
183  // avoids to used lagged qin
184  scalarField Ir = dom.IRay(0).qin().boundaryField()[patchi];
185 
186  for (label rayI=1; rayI < dom.nRay(); rayI++)
187  {
188  Ir += dom.IRay(rayI).qin().boundaryField()[patchi];
189  }
190 
191  forAll(Iw, facei)
192  {
193  if ((-n[facei] & myRayId) > 0.0)
194  {
195  // direction out of the wall
196  refGrad()[facei] = 0.0;
197  valueFraction()[facei] = 1.0;
198  refValue()[facei] =
199  (
200  Ir[facei]*(scalar(1) - emissivity[facei])
201  + emissivity[facei]*physicoChemical::sigma.value()
202  * pow4(Tp[facei])
203  )/pi;
204 
205  // Emitted heat flux from this ray direction
206  qem[facei] = refValue()[facei]*nAve[facei];
207  }
208  else
209  {
210  // direction into the wall
211  valueFraction()[facei] = 0.0;
212  refGrad()[facei] = 0.0;
213  refValue()[facei] = 0.0; // not used
214 
215  // Incident heat flux on this ray direction
216  qin[facei] = Iw[facei]*nAve[facei];
217  }
218  }
219 
220  // Restore tag
221  UPstream::msgType() = oldTag;
222 
223  mixedFvPatchScalarField::updateCoeffs();
224 }
225 
226 
228 (
229  Ostream& os
230 ) const
231 {
234  writeEntryIfDifferent<word>(os, "T", "T", TName_);
235 }
236 
237 
238 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
239 
240 namespace Foam
241 {
243  (
246  );
247 }
248 
249 
250 // ************************************************************************* //
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
static int & msgType()
Message tag of standard messages.
Definition: UPstream.H:476
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Foam::fvPatchFieldMapper.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:87
friend Ostream & operator(Ostream &, const fvPatchField< Type > &)
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:251
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
This boundary condition provides a grey-diffuse condition for radiation intensity,...
greyDiffusiveRadiationMixedFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
virtual void reset(const fvPatchScalarField &)
Reset the fvPatchField to the given fvPatchField.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void map(const fvPatchScalarField &, const fvPatchFieldMapper &)
Map the given fvPatchField onto this fvPatchField.
Common functions to emissivity. It gets supplied from lookup into a dictionary or calculated by the s...
void write(Ostream &) const
Write.
virtual void reset(const fvPatchScalarField &)
Reset the fvPatchField to the given fvPatchField.
virtual void map(const fvPatchScalarField &, const fvPatchFieldMapper &)
Map the given fvPatchField onto this fvPatchField.
Finite Volume Discrete Ordinates Method. Solves the RTE equation for n directions in a participating ...
Definition: fvDOM.H:77
void setRayIdLambdaId(const word &name, label &rayId, label &lambdaId) const
Set the rayId and lambdaId from by decomposing an intensity.
Definition: fvDOM.C:576
Radiation intensity for a ray in a given direction.
const vector & dAve() const
Return the average vector inside the solid angle.
const volScalarField & qr() const
Return const access to the boundary heat flux.
volScalarField & qin()
Return non-const access to the boundary incident heat flux.
volScalarField & qem()
Return non-const access to the boundary emitted heat flux.
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
label patchi
mathematical constants.
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
Collection of constants.
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
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
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61
dimensionedScalar pow4(const dimensionedScalar &ds)
error FatalError
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
static const char nl
Definition: Ostream.H:260
dictionary dict
volScalarField & p