boundaryRadiationPropertiesFvPatchField.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) 2015 OpenCFD Ltd
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 
27 #include "volFields.H"
28 #include "mappedPatchBase.H"
29 #include "fvPatchFieldMapper.H"
30 #include "radiationModel.H"
33 
34 // * * * * * * * * * * * * * Static Member Data * * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38  template<>
39  const char* Foam::NamedEnum
40  <
42  3
43  >::names[] =
44  {
45  "solidRadiation",
46  "lookup",
47  "model"
48  };
49 }
50 
51 const Foam::NamedEnum
52 <
54  3
56 
57 
58 // * * * * * * * * * * * * * * * * Private functions * * * * * * * * * * * * //
59 
62 {
63  // Get the coupling information from the mappedPatchBase
64  const mappedPatchBase& mpp =
65  refCast<const mappedPatchBase>(patch().patch());
66 
67  return (mpp.samplePolyPatch().index());
68 }
69 
70 
71 const Foam::fvMesh&
73 {
74  const mappedPatchBase& mpp =
75  refCast<const mappedPatchBase>(patch().patch());
76 
77  return (refCast<const fvMesh>(mpp.sampleMesh()));
78 }
79 
80 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
81 
82 
85 (
86  const fvPatch& p,
88 )
89 :
90  calculatedFvPatchScalarField(p, iF),
91  method_(LOOKUP),
92  dict_(),
93  absorptionEmission_(NULL),
94  transmissivity_(NULL)
95 {}
96 
97 
100 (
101  const fvPatch& p,
103  const dictionary& dict
104 )
105 :
106  calculatedFvPatchScalarField(p, iF),
107  method_(methodTypeNames_.read(dict.lookup("mode"))),
108  dict_(dict),
109  absorptionEmission_(NULL),
110  transmissivity_(NULL)
111 {
112  switch (method_)
113  {
114  case SOLIDRADIATION:
115  {
116  if (!isA<mappedPatchBase>(p.patch()))
117  {
119  << "\n patch type '" << p.type()
120  << "' not type '" << mappedPatchBase::typeName << "'"
121  << "\n for patch " << p.name()
122  << abort(FatalIOError);
123  }
124  }
125  break;
126 
127  case MODEL:
128  {
129  const fvMesh& mesh = this->dimensionedInternalField().mesh();
130 
131  //if (dict.found("absorptionEmissionModel"))
132  {
133  absorptionEmission_.reset
134  (
136  );
137  }
138 
139  // if (dict.found("transmissivityModel"))
140  {
141  transmissivity_.reset
142  (
144  );
145  }
146  }
147  case LOOKUP:
148  {
149  //Do nothing
150  }
151  break;
152  }
153 
154  if (dict.found("value"))
155  {
157  (
158  scalarField("value", dict, p.size())
159  );
160 
161  }
162  else
163  {
165  }
166 }
167 
168 
171 (
173  const fvPatch& p,
175  const fvPatchFieldMapper& mapper
176 )
177 :
178  calculatedFvPatchScalarField(ptf, p, iF, mapper),
179  method_(ptf.method_),
180  dict_(ptf.dict_),
181  absorptionEmission_(NULL),
182  transmissivity_(NULL)
183 {}
184 
185 
188 (
190 )
191 :
192  calculatedFvPatchScalarField(ptf),
193  method_(ptf.method_),
194  dict_(ptf.dict_),
195  absorptionEmission_(NULL),
196  transmissivity_(NULL)
197 {}
198 
199 
202 (
205 )
206 :
207  calculatedFvPatchScalarField(ptf, iF),
208  method_(ptf.method_),
209  dict_(ptf.dict_),
210  absorptionEmission_(NULL),
211  transmissivity_(NULL)
212 {}
213 
214 
215 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
216 
220 {
221  return absorptionEmission_();
222 }
223 
224 
228 {
229  return transmissivity_();
230 }
231 
232 
235 emissivity(const label bandI) const
236 {
237  switch (method_)
238  {
239  case SOLIDRADIATION:
240  {
241  const fvMesh& nbrMesh = nbrRegion();
242 
245  (
246  "radiationProperties"
247  );
248 
249  scalarField emissivity
250  (
251  radiation.absorptionEmission().e(bandI)().boundaryField()
252  [
253  nbrPatchIndex()
254  ]
255  );
256 
257  const mappedPatchBase& mpp =
258  refCast<const mappedPatchBase>(patch().patch());
259 
260  mpp.distribute(emissivity);
261 
262  const tmp<scalarField> te(new scalarField(emissivity));
263 
264  return te;
265 
266  }
267  break;
268 
269  case LOOKUP:
270  {
272  (
273  new scalarField("emissivity", dict_, patch().size())
274  );
275 
276  return e;
277  }
278 
279  case MODEL:
280  {
281  const label index = patch().index();
282 
284  (
285  new scalarField
286  (
287  absorptionEmission_->e(bandI)().boundaryField()[index]
288  )
289  );
290 
291  return e;
292  }
293 
294  default:
295  {
297  << "Please set 'mode' to one of " << methodTypeNames_.toc()
298  << exit(FatalError);
299  }
300  break;
301  }
302 
303  return scalarField(0);
304 }
305 
306 
309 (
310  const label bandI
311 ) const
312 {
313  switch (method_)
314  {
315  case SOLIDRADIATION:
316  {
317  const fvMesh& nbrMesh = nbrRegion();
318 
321  (
322  "radiationProperties"
323  );
324 
325  scalarField absorp
326  (
327  radiation.absorptionEmission().a(bandI)().boundaryField()
328  [
329  nbrPatchIndex()
330  ]
331  );
332 
333  const mappedPatchBase& mpp =
334  refCast<const mappedPatchBase>(patch().patch());
335 
336  mpp.distribute(absorp);
337 
338  const tmp<scalarField> ta(new scalarField(absorp));
339 
340  return ta;
341 
342  }
343  break;
344 
345  case MODEL:
346  {
347  const label index = patch().index();
349  (
350  new scalarField
351  (
352  absorptionEmission_->a(bandI)().boundaryField()[index]
353  )
354  );
355  return a;
356  }
357 
358  case LOOKUP:
359  {
361  (
362  new scalarField("absorptivity", dict_, patch().size())
363  );
364 
365  return a;
366  }
367 
368  default:
369  {
371  << "Unimplemented method " << method_ << endl
372  << "Please set 'mode' to one of "
373  << methodTypeNames_.toc()
374  << exit(FatalError);
375  }
376  break;
377  }
378 
379  return scalarField(0);
380 }
381 
382 
385 transmissivity(const label bandI) const
386 {
387  switch (method_)
388  {
389  case SOLIDRADIATION:
390  {
391  const fvMesh& nbrMesh = nbrRegion();
392 
395  (
396  "radiationProperties"
397  );
398 
400  (
401  radiation.transmissivity().tauEff(bandI)().boundaryField()
402  [
403  nbrPatchIndex()
404  ]
405  );
406 
407  const mappedPatchBase& mpp =
408  refCast<const mappedPatchBase>(patch().patch());
409 
410  mpp.distribute(trans);
411 
412  const tmp<scalarField> tt(new scalarField(trans));
413 
414  return tt;
415 
416  }
417  break;
418 
419  case MODEL:
420  {
421  const label index = patch().index();
422  tmp<scalarField> tau
423  (
424  new scalarField
425  (
426  transmissivity_->tauEff(bandI)().boundaryField()[index]
427  )
428  );
429  return tau;
430  }
431 
432  case LOOKUP:
433  {
434  tmp<scalarField> tau
435  (
436  new scalarField
437  (
438  "transmissivity", dict_, patch().size()
439  )
440  );
441  return tau;
442  }
443 
444  default:
445  {
447  << "Unimplemented method " << method_ << endl
448  << "Please set 'mode' to one of "
449  << methodTypeNames_.toc()
450  << exit(FatalError);
451  }
452  break;
453  }
454 
455  return scalarField(0);
456 }
457 
458 
459 
462 reflectivity(const label bandI) const
463 {
464  const tmp<scalarField> tt = transmissivity(bandI);
465  const tmp<scalarField> ta = absorptivity(bandI);
466 
467  return (1.0 - tt - ta);
468 }
469 
470 
472 write(Ostream& os) const
473 {
475 
476  os.writeKeyword("mode") << methodTypeNames_[method_]
477  << token::END_STATEMENT << nl;
478 
479  switch (method_)
480  {
481  case MODEL:
482  {
483  word modelType
484  (
485  word(dict_.lookup("absorptionEmissionModel"))
486  );
487 
488  os.writeKeyword("absorptionEmissionModel") << modelType
489  << token::END_STATEMENT << nl;
490 
491  word modelCoeffs(modelType + word("Coeffs"));
492  os.writeKeyword(modelCoeffs);
493 
494  dict_.subDict(modelCoeffs).write(os);
495 
496  modelType = word(dict_.lookup("transmissivityModel"));
497 
498  os.writeKeyword("transmissivityModel") << modelType
499  << token::END_STATEMENT << nl;
500 
501  modelCoeffs = modelType + word("Coeffs");
502 
503  os.writeKeyword(modelCoeffs);
504 
505  dict_.subDict(modelCoeffs).write(os);
506 
507  break;
508  }
509 
510  case LOOKUP:
511  {
512  const scalarField emissivity("emissivity", dict_, patch().size());
513  emissivity.writeEntry("emissivity", os);
514 
515  const scalarField absorptivity
516  (
517  "absorptivity", dict_, patch().size()
518  );
519  absorptivity.writeEntry("absorptivity", os);
520 
521  const scalarField transmissivity
522  (
523  "transmissivity", dict_, patch().size()
524  );
525  transmissivity.writeEntry("transmissivity", os);
526 
527  break;
528  }
529 
530  case SOLIDRADIATION:
531  {
532  }
533  }
534 }
535 
536 
537 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
538 
539 namespace Foam
540 {
541 namespace radiation
542 {
544  (
547  );
548 }
549 }
550 
551 // ************************************************************************* //
boundaryRadiationPropertiesFvPatchField.H
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
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::radiation::absorptionEmissionModel::New
static autoPtr< absorptionEmissionModel > New(const dictionary &dict, const fvMesh &mesh)
Selector.
Definition: absorptionEmissionModelNew.C:33
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::radiation::boundaryRadiationPropertiesFvPatchField::methodType
methodType
Definition: boundaryRadiationPropertiesFvPatchField.H:66
Foam::radiation::boundaryRadiationPropertiesFvPatchField::boundaryRadiationPropertiesFvPatchField
boundaryRadiationPropertiesFvPatchField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: boundaryRadiationPropertiesFvPatchField.C:85
dimensionedInternalField
rDeltaT dimensionedInternalField()
Foam::mappedPatchBase::samplePolyPatch
const polyPatch & samplePolyPatch() const
Get the patch on the region.
Definition: mappedPatchBase.C:1246
Foam::radiation::boundaryRadiationPropertiesFvPatchField::absorptivity
tmp< scalarField > absorptivity(const label bandI=0) const
Calculate corresponding absorptivity field for bandI.
Definition: boundaryRadiationPropertiesFvPatchField.C:309
Foam::radiation::makePatchTypeField
makePatchTypeField(fvPatchScalarField, boundaryRadiationPropertiesFvPatchField)
Foam::radiation::boundaryRadiationPropertiesFvPatchField::nbrRegion
const fvMesh & nbrRegion() const
Return nbr mesh.
Definition: boundaryRadiationPropertiesFvPatchField.C:72
Foam::fvPatchField::operator
friend Ostream & operator(Ostream &, const fvPatchField< Type > &)
Foam::FatalIOError
IOerror FatalIOError
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
fvPatchFieldMapper.H
Foam::radiation::boundaryRadiationPropertiesFvPatchField::reflectivity
tmp< scalarField > reflectivity(const label bandI=0) const
Calculate corresponding reflectivity field.
Definition: boundaryRadiationPropertiesFvPatchField.C:462
Foam::mappedPatchBase
Determines a mapping between patch face centres and mesh cell or face centres and processors they're ...
Definition: mappedPatchBase.H:101
Foam::radiation::boundaryRadiationPropertiesFvPatchField::absorptionEmission
const absorptionEmissionModel & absorptionEmission() const
Return absorptionEmissionModel.
Definition: boundaryRadiationPropertiesFvPatchField.C:219
Foam::radiation::transmissivityModel::New
static autoPtr< transmissivityModel > New(const dictionary &dict, const fvMesh &mesh)
Definition: transmissivityModelNew.C:33
Foam::radiation::transmissivityModel
Base class for radiation scattering.
Definition: transmissivityModel.H:50
Foam::radiation::boundaryRadiationPropertiesFvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: boundaryRadiationPropertiesFvPatchField.C:472
Foam::mappedPatchBase::sampleMesh
const polyMesh & sampleMesh() const
Get the region mesh.
Definition: mappedPatchBase.C:1237
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::Ostream::write
virtual Ostream & write(const token &)=0
Write next token to stream.
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::mappedPatchBase::distribute
void distribute(List< Type > &lst) const
Wrapper around map/interpolate data distribution.
Definition: mappedPatchBaseTemplates.C:27
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
radiation
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
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::radiation::boundaryRadiationPropertiesFvPatchField::transmissivity
tmp< scalarField > transmissivity(const label bandI=0) const
Calculate corresponding transmissivity field for bandI.
Definition: boundaryRadiationPropertiesFvPatchField.C:385
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Foam::e
const double e
Elementary charge.
Definition: doubleFloat.H:94
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::Field::writeEntry
void writeEntry(const word &keyword, Ostream &os) const
Write the field as a dictionary entry.
Definition: Field.C:700
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
absorptionEmissionModel.H
Foam::radiation::boundaryRadiationPropertiesFvPatchField::transmissiveModel
const transmissivityModel & transmissiveModel() const
Return transmissivityModel.
Definition: boundaryRadiationPropertiesFvPatchField.C:227
Foam::radiation::radiationModel
Top level model for radiation modelling.
Definition: radiationModel.H:73
Foam::radiation::boundaryRadiationPropertiesFvPatchField::dict_
const dictionary dict_
Dictionary.
Definition: boundaryRadiationPropertiesFvPatchField.H:84
Foam::Ostream::writeKeyword
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:59
Foam::radiation::boundaryRadiationPropertiesFvPatchField::methodTypeNames_
static const NamedEnum< methodType, 3 > methodTypeNames_
Definition: boundaryRadiationPropertiesFvPatchField.H:78
Foam::radiation::boundaryRadiationPropertiesFvPatchField
Definition: boundaryRadiationPropertiesFvPatchField.H:59
Foam::radiation::boundaryRadiationPropertiesFvPatchField::emissivity
tmp< scalarField > emissivity(const label bandI=0) const
Calculate corresponding emissivity field for bandI.
Definition: boundaryRadiationPropertiesFvPatchField.C:235
Foam::radiation::absorptionEmissionModel
Model to supply absorption and emission coefficients for radiation modelling.
Definition: absorptionEmissionModel.H:52
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:45
Foam::objectRegistry::lookupObject
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
Definition: objectRegistryTemplates.C:165
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::trans
dimensionSet trans(const dimensionSet &)
Definition: dimensionSet.C:438
Foam::patchIdentifier::index
label index() const
Return the index of this patch in the boundaryMesh.
Definition: patchIdentifier.H:133
Foam::radiation::boundaryRadiationPropertiesFvPatchField::method_
const methodType method_
How to get property.
Definition: boundaryRadiationPropertiesFvPatchField.H:81
radiationModel.H
write
Tcoeff write()
Foam::radiation::boundaryRadiationPropertiesFvPatchField::nbrPatchIndex
label nbrPatchIndex() const
Return nbr patch index.
Definition: boundaryRadiationPropertiesFvPatchField.C:61
Foam::NamedEnum
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:52
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:51
mappedPatchBase.H