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 | 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"
33 #include "constants.H"
35 
36 using namespace Foam::constant;
37 using namespace Foam::constant::mathematical;
38 
39 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40 
43 (
44  const fvPatch& p,
46 )
47 :
48  mixedFvPatchScalarField(p, iF)
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)
66 {}
67 
68 
71 (
72  const fvPatch& p,
74  const dictionary& dict
75 )
76 :
77  mixedFvPatchScalarField(p, iF)
78 {
79  if (dict.found("value"))
80  {
82  (
83  scalarField("value", dict, p.size())
84  );
85  refValue() = scalarField("refValue", dict, p.size());
86  refGrad() = scalarField("refGradient", dict, p.size());
87  valueFraction() = scalarField("valueFraction", dict, p.size());
88  }
89  else
90  {
91  refValue() = 0.0;
92  refGrad() = 0.0;
93  valueFraction() = 1.0;
94 
96  }
97 }
98 
99 
102 (
104 )
105 :
106  mixedFvPatchScalarField(ptf)
107 {}
108 
109 
112 (
115 )
116 :
117  mixedFvPatchScalarField(ptf, iF)
118 {}
119 
120 
121 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
122 
125 {
126  if (this->updated())
127  {
128  return;
129  }
130 
131  // Since we're inside initEvaluate/evaluate there might be processor
132  // comms underway. Change the tag we use.
133  int oldTag = UPstream::msgType();
134  UPstream::msgType() = oldTag+1;
135 
136  const radiationModel& radiation =
137  db().lookupObject<radiationModel>("radiationProperties");
138 
139  const fvDOM& dom(refCast<const fvDOM>(radiation));
140 
141  label rayId = -1;
142  label lambdaId = -1;
143  dom.setRayIdLambdaId(dimensionedInternalField().name(), rayId, lambdaId);
144 
145  const label patchI = patch().index();
146 
147  if (dom.nLambda() == 0)
148  {
150  << " a non-grey boundary condition is used with a grey "
151  << "absorption model" << nl << exit(FatalError);
152  }
153 
154  scalarField& Iw = *this;
155  const vectorField n(patch().Sf()/patch().magSf());
156 
157  radiativeIntensityRay& ray =
158  const_cast<radiativeIntensityRay&>(dom.IRay(rayId));
159 
160  const scalarField nAve(n & ray.dAve());
161 
162  ray.Qr().boundaryField()[patchI] += Iw*nAve;
163 
164  const scalarField Eb
165  (
166  dom.blackBody().bLambda(lambdaId).boundaryField()[patchI]
167  );
168 
169  const boundaryRadiationProperties& boundaryRadiation =
171 
172 
173  const tmp<scalarField> temissivity
174  (
175  boundaryRadiation.emissivity(patch().index(), lambdaId)
176  );
177 
178  const scalarField& emissivity = temissivity();
179 
180  scalarField& Qem = ray.Qem().boundaryField()[patchI];
181  scalarField& Qin = ray.Qin().boundaryField()[patchI];
182 
183  // Use updated Ir while iterating over rays
184  // avoids to used lagged Qin
185  scalarField Ir = dom.IRay(0).Qin().boundaryField()[patchI];
186 
187  for (label rayI=1; rayI < dom.nRay(); rayI++)
188  {
189  Ir += dom.IRay(rayI).Qin().boundaryField()[patchI];
190  }
191 
192  forAll(Iw, faceI)
193  {
194  const vector& d = dom.IRay(rayId).d();
195 
196  if ((-n[faceI] & d) > 0.0)
197  {
198  // direction out of the wall
199  refGrad()[faceI] = 0.0;
200  valueFraction()[faceI] = 1.0;
201  refValue()[faceI] =
202  (
203  Ir[faceI]*(1.0 - emissivity[faceI])
204  + emissivity[faceI]*Eb[faceI]
205  )/pi;
206 
207  // Emmited heat flux from this ray direction
208  Qem[faceI] = refValue()[faceI]*nAve[faceI];
209  }
210  else
211  {
212  // direction into the wall
213  valueFraction()[faceI] = 0.0;
214  refGrad()[faceI] = 0.0;
215  refValue()[faceI] = 0.0; //not used
216 
217  // Incident heat flux on this ray direction
218  Qin[faceI] = Iw[faceI]*nAve[faceI];
219  }
220  }
221 
222  // Restore tag
223  UPstream::msgType() = oldTag;
224 
225  mixedFvPatchScalarField::updateCoeffs();
226 }
227 
228 
230 (
231  Ostream& os
232 ) const
233 {
235 }
236 
237 
238 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
239 
240 namespace Foam
241 {
242 namespace radiation
243 {
245  (
248  );
249 }
250 }
251 
252 
253 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
Foam::radiation::boundaryRadiationProperties::emissivity
tmp< scalarField > emissivity(const label patchId, const label bandI=0) const
Access boundary emissivity on patch.
Definition: boundaryRadiationProperties.C:110
volFields.H
Foam::fvPatchField::operator=
virtual void operator=(const UList< Type > &)
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:48
p
p
Definition: pEqn.H:62
Foam::compressible::New
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Definition: turbulentFluidThermoModel.C:36
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::constant
Collection of constants.
Definition: atomicConstants.C:37
dimensionedInternalField
rDeltaT dimensionedInternalField()
Foam::radiation::radiativeIntensityRay
Radiation intensity for a ray in a given direction.
Definition: radiativeIntensityRay.H:55
Foam::GeometricField::boundaryField
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Definition: GeometricField.C:735
Foam::fvPatchField::operator
friend Ostream & operator(Ostream &, const fvPatchField< Type > &)
fvPatchFieldMapper.H
fvDOM.H
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::label
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
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::radiation::wideBandDiffusiveRadiationMixedFvPatchScalarField
This boundary condition provides a wide-band, diffusive radiation condition, where the patch temperat...
Definition: wideBandDiffusiveRadiationMixedFvPatchScalarField.H:84
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Foam::radiation::radiativeIntensityRay::Qem
volScalarField & Qem()
Return non-const access to the boundary emmited heat flux.
Definition: radiativeIntensityRayI.H:59
radiation
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
Foam::radiation::boundaryRadiationProperties
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 any number of values (e....
Definition: dictionary.H:137
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::makePatchTypeField
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
wideBandAbsorptionEmission.H
constants.H
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::radiation::radiativeIntensityRay::dAve
const vector & dAve() const
Return the average vector inside the solid angle.
Definition: radiativeIntensityRayI.H:77
Foam::radiation::wideBandDiffusiveRadiationMixedFvPatchScalarField::wideBandDiffusiveRadiationMixedFvPatchScalarField
wideBandDiffusiveRadiationMixedFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: wideBandDiffusiveRadiationMixedFvPatchScalarField.C:43
Foam::UPstream::msgType
static int & msgType()
Message tag of standard messages.
Definition: UPstream.H:452
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::constant::mathematical
mathematical constants.
Foam::Vector< scalar >
wideBandDiffusiveRadiationMixedFvPatchScalarField.H
Foam::radiation::radiationModel
Top level model for radiation modelling.
Definition: radiationModel.H:73
Foam::constant::mathematical::pi
const scalar pi(M_PI)
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:45
Foam::radiation::wideBandDiffusiveRadiationMixedFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: wideBandDiffusiveRadiationMixedFvPatchScalarField.C:124
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::radiation::radiativeIntensityRay::Qr
const volScalarField & Qr() const
Return const access to the boundary heat flux.
Definition: radiativeIntensityRayI.H:34
Foam::radiation::wideBandDiffusiveRadiationMixedFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: wideBandDiffusiveRadiationMixedFvPatchScalarField.C:230
write
Tcoeff write()
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
boundaryRadiationProperties.H
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:51
Foam::radiation::radiativeIntensityRay::Qin
volScalarField & Qin()
Return non-const access to the boundary incident heat flux.
Definition: radiativeIntensityRayI.H:46
Foam::radiation::fvDOM
Finite Volume Discrete Ordinates Method. Solves the RTE equation for n directions in a participating ...
Definition: fvDOM.H:91