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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2016-2019 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
31 #include "fvPatchFieldMapper.H"
32 #include "volFields.H"
33 
34 #include "fvDOM.H"
36 #include "constants.H"
38 
39 using namespace Foam::constant;
40 using namespace Foam::constant::mathematical;
41 
42 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43 
46 (
47  const fvPatch& p,
49 )
50 :
51  mixedFvPatchScalarField(p, iF)
52 {
53  refValue() = 0.0;
54  refGrad() = 0.0;
55  valueFraction() = 1.0;
56 }
57 
58 
61 (
63  const fvPatch& p,
65  const fvPatchFieldMapper& mapper
66 )
67 :
68  mixedFvPatchScalarField(ptf, p, iF, mapper)
69 {}
70 
71 
74 (
75  const fvPatch& p,
77  const dictionary& dict
78 )
79 :
80  mixedFvPatchScalarField(p, iF)
81 {
82  if (dict.found("refValue"))
83  {
85  (
86  scalarField("value", dict, p.size())
87  );
88  refValue() = scalarField("refValue", dict, p.size());
89  refGrad() = scalarField("refGradient", dict, p.size());
90  valueFraction() = scalarField("valueFraction", dict, p.size());
91  }
92  else
93  {
94  refValue() = 0.0;
95  refGrad() = 0.0;
96  valueFraction() = 1.0;
97 
99  }
100 }
101 
102 
105 (
107 )
108 :
109  mixedFvPatchScalarField(ptf)
110 {}
111 
112 
115 (
118 )
119 :
120  mixedFvPatchScalarField(ptf, iF)
121 {}
122 
123 
124 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
125 
128 {
129  if (this->updated())
130  {
131  return;
132  }
133 
134  // Since we're inside initEvaluate/evaluate there might be processor
135  // comms underway. Change the tag we use.
136  int oldTag = UPstream::msgType();
137  UPstream::msgType() = oldTag+1;
138 
139  const radiationModel& radiation =
140  db().lookupObject<radiationModel>("radiationProperties");
141 
142  const fvDOM& dom(refCast<const fvDOM>(radiation));
143 
144  label rayId = -1;
145  label lambdaId = -1;
146  dom.setRayIdLambdaId(internalField().name(), rayId, lambdaId);
147 
148  const label patchi = patch().index();
149 
150  if (dom.nLambda() == 0)
151  {
153  << " a non-grey boundary condition is used with a grey "
154  << "absorption model" << nl << exit(FatalError);
155  }
156 
157  scalarField& Iw = *this;
158  const vectorField n(patch().Sf()/patch().magSf());
159 
160  radiativeIntensityRay& ray =
161  const_cast<radiativeIntensityRay&>(dom.IRay(rayId));
162 
163  const scalarField nAve(n & ray.dAve());
164 
165  ray.qr().boundaryFieldRef()[patchi] += Iw*nAve;
166 
167  const scalarField Eb
168  (
169  dom.blackBody().bLambda(lambdaId).boundaryField()[patchi]
170  );
171 
172  const boundaryRadiationProperties& boundaryRadiation =
173  boundaryRadiationProperties::New(internalField().mesh());
174 
175 
176  const tmp<scalarField> temissivity
177  (
178  boundaryRadiation.emissivity(patch().index(), lambdaId)
179  );
180  const scalarField& emissivity = temissivity();
181 
182  const tmp<scalarField> ttransmissivity
183  (
184  boundaryRadiation.transmissivity(patch().index(), lambdaId)
185  );
186  const scalarField& transmissivity = ttransmissivity();
187 
188  scalarField& qem = ray.qem().boundaryFieldRef()[patchi];
189  scalarField& qin = ray.qin().boundaryFieldRef()[patchi];
190 
191  // Calculate Ir into the wall on the same lambdaId
192  scalarField Ir(patch().size(), Zero);
193  forAll(Iw, facei)
194  {
195  for (label rayi=0; rayi < dom.nRay(); rayi++)
196  {
197  const vector& d = dom.IRay(rayi).d();
198 
199  if ((-n[facei] & d) < 0.0)
200  {
201  // q into the wall
202  const scalarField& IFace =
203  dom.IRay(rayi).ILambda(lambdaId).boundaryField()[patchi];
204 
205  const vector& rayDave = dom.IRay(rayi).dAve();
206  Ir[facei] += IFace[facei]*(n[facei] & rayDave);
207  }
208  }
209  }
210 
211  if (dom.useSolarLoad())
212  {
213  // Looking for primary heat flux single band
214  Ir += patch().lookupPatchField<volScalarField,scalar>
215  (
216  dom.primaryFluxName_ + "_" + name(lambdaId)
217  );
218 
219  word qSecName = dom.relfectedFluxName_ + "_" + name(lambdaId);
220 
221  if (this->db().foundObject<volScalarField>(qSecName))
222  {
223  const volScalarField& qSec =
224  this->db().lookupObject<volScalarField>(qSecName);
225 
226  Ir += qSec.boundaryField()[patch().index()];
227  }
228  }
229 
230  scalarField Iexternal(this->size(), 0.0);
231 
232  if (dom.useExternalBeam())
233  {
234  const vector sunDir = dom.solarCalc().direction();
235  const scalar directSolarRad =
236  dom.solarCalc().directSolarRad()
237  *dom.spectralDistribution()[lambdaId];
238 
239  //label nRaysBeam = dom.nRaysBeam();
240  label SunRayId(-1);
241  scalar maxSunRay = -GREAT;
242 
243  // Looking for the ray closest to the Sun direction
244  for (label rayI=0; rayI < dom.nRay(); rayI++)
245  {
246  const vector& iD = dom.IRay(rayI).d();
247  scalar dir = sunDir & iD;
248  if (dir > maxSunRay)
249  {
250  maxSunRay = dir;
251  SunRayId = rayI;
252  }
253  }
254 
255  if (rayId == SunRayId)
256  {
257  const scalarField nAve(n & dom.IRay(rayId).dAve());
258  forAll(Iexternal, faceI)
259  {
260  Iexternal[faceI] = directSolarRad/mag(dom.IRay(rayId).dAve());
261  }
262  }
263  }
264 
265  forAll(Iw, facei)
266  {
267  const vector& d = dom.IRay(rayId).d();
268 
269  if ((-n[facei] & d) > 0.0)
270  {
271  // direction out of the wall
272  refGrad()[facei] = 0.0;
273  valueFraction()[facei] = 1.0;
274  refValue()[facei] =
275  Iexternal[facei]*transmissivity[facei]
276  + (
277  Ir[facei]*(1.0 - emissivity[facei])
278  + emissivity[facei]*Eb[facei]
279  )/pi;
280 
281  // Emitted heat flux from this ray direction (sum over lambdaId)
282  qem[facei] += refValue()[facei]*nAve[facei];
283  }
284  else
285  {
286  // direction into the wall
287  valueFraction()[facei] = 0.0;
288  refGrad()[facei] = 0.0;
289  refValue()[facei] = 0.0; //not used
290 
291  // Incident heat flux on this ray direction (sum over lambdaId)
292  qin[facei] += Iw[facei]*nAve[facei];
293  }
294  }
295 
296  // Restore tag
297  UPstream::msgType() = oldTag;
298 
299  mixedFvPatchScalarField::updateCoeffs();
300 }
301 
302 
304 (
305  Ostream& os
306 ) const
307 {
309 }
310 
311 
312 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
313 
314 namespace Foam
315 {
316 namespace radiation
317 {
319  (
322  );
323 }
324 }
325 
326 
327 // ************************************************************************* //
Foam::fvPatchField< scalar >
volFields.H
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::constant
Different types of constants.
Definition: atomicConstants.C:38
Foam::radiation::radiativeIntensityRay::qem
volScalarField & qem()
Return non-const access to the boundary emitted heat flux.
Definition: radiativeIntensityRayI.H:67
Foam::MeshObject< fvMesh, Foam::GeometricMeshObject, boundaryRadiationProperties >::New
static const boundaryRadiationProperties & New(const fvMesh &mesh, Args &&... args)
Get existing or create a new MeshObject.
Definition: MeshObject.C:48
Foam::radiation::radiativeIntensityRay
Radiation intensity for a ray in a given direction.
Definition: radiativeIntensityRay.H:57
Foam::fvPatchField< scalar >::operator
friend Ostream & operator(Ostream &, const fvPatchField< scalar > &)
fvPatchFieldMapper.H
Foam::radiation::radiativeIntensityRay::qr
const volScalarField & qr() const
Return const access to the boundary heat flux.
Definition: radiativeIntensityRayI.H:36
fvDOM.H
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::Field< scalar >
Foam::radiation::wideBandDiffusiveRadiationMixedFvPatchScalarField
This boundary condition provides a wide-band, diffusive radiation condition, where the patch temperat...
Definition: wideBandDiffusiveRadiationMixedFvPatchScalarField.H:85
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::radiation::radiativeIntensityRay::qin
volScalarField & qin()
Return non-const access to the boundary incident heat flux.
Definition: radiativeIntensityRayI.H:54
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
radiation
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
Foam::radiation::boundaryRadiationProperties
Boundary radiation properties holder.
Definition: boundaryRadiationProperties.H:56
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
os
OBJstream os(runTime.globalPath()/outputName)
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
wideBandAbsorptionEmission.H
constants.H
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::radiation::radiativeIntensityRay::dAve
const vector & dAve() const
Return the average vector inside the solid angle.
Definition: radiativeIntensityRayI.H:79
Foam::radiation::wideBandDiffusiveRadiationMixedFvPatchScalarField::wideBandDiffusiveRadiationMixedFvPatchScalarField
wideBandDiffusiveRadiationMixedFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: wideBandDiffusiveRadiationMixedFvPatchScalarField.C:46
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::constant::mathematical
Mathematical constants.
Foam::UPstream::msgType
static int & msgType() noexcept
Message tag of standard messages.
Definition: UPstream.H:540
Foam::constant::mathematical::pi
constexpr scalar pi(M_PI)
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::GeometricField::boundaryFieldRef
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Definition: GeometricField.C:783
Foam::Vector< scalar >
wideBandDiffusiveRadiationMixedFvPatchScalarField.H
Foam::radiation::radiationModel
Top level model for radiation modelling.
Definition: radiationModel.H:75
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::vtk::write
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
Definition: foamVtkOutputTemplates.C:36
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::radiation::wideBandDiffusiveRadiationMixedFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: wideBandDiffusiveRadiationMixedFvPatchScalarField.C:127
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::makePatchTypeField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Foam::radiation::wideBandDiffusiveRadiationMixedFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: wideBandDiffusiveRadiationMixedFvPatchScalarField.C:304
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::fvPatchField< scalar >::operator=
virtual void operator=(const UList< scalar > &)
Definition: fvPatchField.C:404
boundaryRadiationProperties.H
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
Foam::radiation::fvDOM
Finite Volume Discrete Ordinates Method. Solves the RTE equation for n directions in a participating ...
Definition: fvDOM.H:118