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"
31 
32 #include "fvDOM.H"
33 #include "constants.H"
34 
35 using namespace Foam::constant;
36 using namespace Foam::constant::mathematical;
37 
38 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
39 
42 (
43  const fvPatch& p,
45 )
46 :
47  mixedFvPatchScalarField(p, iF),
48  TName_("T"),
49  solarLoad_(false)
50 {
51  refValue() = 0.0;
52  refGrad() = 0.0;
53  valueFraction() = 1.0;
54 }
55 
56 
59 (
61  const fvPatch& p,
63  const fvPatchFieldMapper& mapper
64 )
65 :
66  mixedFvPatchScalarField(ptf, p, iF, mapper),
67  TName_(ptf.TName_),
68  solarLoad_(ptf.solarLoad_)
69 {}
70 
71 
74 (
75  const fvPatch& p,
77  const dictionary& dict
78 )
79 :
80  mixedFvPatchScalarField(p, iF),
81  TName_(dict.lookupOrDefault<word>("T", "T")),
82  solarLoad_(dict.lookupOrDefault<bool>("solarLoad", false))
83 {
84  if (dict.found("refValue"))
85  {
87  (
88  scalarField("value", dict, p.size())
89  );
90  refValue() = scalarField("refValue", dict, p.size());
91  refGrad() = scalarField("refGradient", dict, p.size());
92  valueFraction() = scalarField("valueFraction", dict, p.size());
93  }
94  else
95  {
96  refValue() = 0.0;
97  refGrad() = 0.0;
98  valueFraction() = 1.0;
99 
100  fvPatchScalarField::operator=(refValue());
101  }
102 }
103 
104 
107 (
109 )
110 :
111  mixedFvPatchScalarField(ptf),
112  TName_(ptf.TName_),
113  solarLoad_(ptf.solarLoad_)
114 {}
115 
116 
119 (
122 )
123 :
124  mixedFvPatchScalarField(ptf, iF),
125  TName_(ptf.TName_),
126  solarLoad_(ptf.solarLoad_)
127 {}
128 
129 
130 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
131 
134 {
135  if (this->updated())
136  {
137  return;
138  }
139 
140  // Since we're inside initEvaluate/evaluate there might be processor
141  // comms underway. Change the tag we use.
142  int oldTag = UPstream::msgType();
143  UPstream::msgType() = oldTag+1;
144 
145  const scalarField& Tp =
146  patch().lookupPatchField<volScalarField, scalar>(TName_);
147 
148  const radiationModel& radiation =
149  db().lookupObject<radiationModel>("radiationProperties");
150 
151  const fvDOM& dom(refCast<const fvDOM>(radiation));
152 
153  label rayId = -1;
154  label lambdaId = -1;
155  dom.setRayIdLambdaId(dimensionedInternalField().name(), rayId, lambdaId);
156 
157  const label patchI = patch().index();
158 
159  if (dom.nLambda() != 1)
160  {
162  << " a grey boundary condition is used with a non-grey "
163  << "absorption model" << nl << exit(FatalError);
164  }
165 
166  scalarField& Iw = *this;
167 
168  const vectorField n(patch().nf());
169 
170  radiativeIntensityRay& ray =
171  const_cast<radiativeIntensityRay&>(dom.IRay(rayId));
172 
173  const scalarField nAve(n & ray.dAve());
174 
175  ray.Qr().boundaryField()[patchI] += Iw*nAve;
176 
177  const boundaryRadiationProperties& boundaryRadiation =
179 
180  const tmp<scalarField> temissivity
181  (
182  boundaryRadiation.emissivity(patch().index())
183  );
184 
185  const scalarField& emissivity = temissivity();
186 
187  scalarField& Qem = ray.Qem().boundaryField()[patchI];
188  scalarField& Qin = ray.Qin().boundaryField()[patchI];
189 
190  const vector& myRayId = dom.IRay(rayId).d();
191 
192  // Use updated Ir while iterating over rays
193  // avoids to used lagged Qin
194  scalarField Ir = dom.IRay(0).Qin().boundaryField()[patchI];
195 
196  for (label rayI=1; rayI < dom.nRay(); rayI++)
197  {
198  Ir += dom.IRay(rayI).Qin().boundaryField()[patchI];
199  }
200 
201  if (solarLoad_)
202  {
203  Ir += patch().lookupPatchField<volScalarField,scalar>
204  (
205  radiation.externalRadHeatFieldName_
206  );
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) - emissivity[faceI])
219  + emissivity[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 {
251  writeEntryIfDifferent<word>(os, "T", "T", TName_);
252  os.writeKeyword("solarLoad") << solarLoad_ << token::END_STATEMENT << nl;
253 }
254 
255 
256 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
257 
258 namespace Foam
259 {
260 namespace radiation
261 {
263  (
266  );
267 }
268 }
269 
270 
271 // ************************************************************************* //
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
Foam::token::END_STATEMENT
@ END_STATEMENT
Definition: token.H:99
p
p
Definition: pEqn.H:62
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
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::greyDiffusiveRadiationMixedFvPatchScalarField::greyDiffusiveRadiationMixedFvPatchScalarField
greyDiffusiveRadiationMixedFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: greyDiffusiveRadiationMixedFvPatchScalarField.C:42
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
Foam::radiation::greyDiffusiveRadiationMixedFvPatchScalarField
This boundary condition provides a grey-diffuse condition for radiation intensity,...
Definition: greyDiffusiveRadiationMixedFvPatchScalarField.H:88
fvDOM.H
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::pow4
dimensionedScalar pow4(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:98
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::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
Foam::radiation::radiativeIntensityRay::Qem
volScalarField & Qem()
Return non-const access to the boundary emmited heat flux.
Definition: radiativeIntensityRayI.H:59
Foam::radiation::greyDiffusiveRadiationMixedFvPatchScalarField::TName_
word TName_
Name of temperature field.
Definition: greyDiffusiveRadiationMixedFvPatchScalarField.H:95
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)
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::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::radiation::greyDiffusiveRadiationMixedFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: greyDiffusiveRadiationMixedFvPatchScalarField.C:246
Foam::Vector< scalar >
Foam::radiation::radiationModel
Top level model for radiation modelling.
Definition: radiationModel.H:73
Foam::Ostream::writeKeyword
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:59
Foam::radiation::greyDiffusiveRadiationMixedFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: greyDiffusiveRadiationMixedFvPatchScalarField.C:133
Foam::constant::physicoChemical::sigma
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
Foam::constant::mathematical::pi
const scalar pi(M_PI)
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:45
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
greyDiffusiveRadiationMixedFvPatchScalarField.H
Foam::radiation::greyDiffusiveRadiationMixedFvPatchScalarField::solarLoad_
bool solarLoad_
Activate solar load.
Definition: greyDiffusiveRadiationMixedFvPatchScalarField.H:98
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