solidAbsorption.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) 2015-2018 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "solidAbsorption.H"
30 #include "mappedPatchBase.H"
31 #include "radiationModel.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38  namespace radiation
39  {
40  defineTypeNameAndDebug(solidAbsorption, 0);
41 
43  (
44  wallAbsorptionEmissionModel,
45  solidAbsorption,
46  dictionary
47  );
48  }
49 }
50 
51 
52 // * * * * * * * * * * * * * * * * Private members * * * * * * * * * * * * * //
53 
54 const Foam::fvMesh& Foam::radiation::solidAbsorption::nbrRegion() const
55 {
56  const mappedPatchBase& mpp = refCast<const mappedPatchBase>(pp_);
57  return (refCast<const fvMesh>(mpp.sampleMesh()));
58 }
59 
60 
61 Foam::label Foam::radiation::solidAbsorption::nbrPatchIndex() const
62 {
63  const mappedPatchBase& mpp = refCast<const mappedPatchBase>(pp_);
64  return (mpp.samplePolyPatch().index());
65 }
66 
67 
68 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
69 
71 (
72  const dictionary& dict,
73  const polyPatch& pp
74 )
75 :
77 {
78  if (!isA<mappedPatchBase>(pp))
79  {
81  << "\n patch type '" << pp.type()
82  << "' not type '" << mappedPatchBase::typeName << "'"
83  << "\n for patch " << pp.name()
84  << abort(FatalIOError);
85  }
86 }
87 
88 
89 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
90 
92 {}
93 
94 
95 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
96 
98 (
99  const label bandI,
100  vectorField* incomingDirection,
101  scalarField* T
102 ) const
103 {
104  // Since we're inside initEvaluate/evaluate there might be processor
105  // comms underway. Change the tag we use.
106  int oldTag = UPstream::msgType();
107  UPstream::msgType() = oldTag+1;
108 
109  const fvMesh& nbrMesh = nbrRegion();
110 
113  (
114  "radiationProperties"
115  );
116 
117  scalarField absorptivity
118  (
119  radiation.absorptionEmission().a
120  (
121  bandI
122  )().boundaryField()
123  [
124  nbrPatchIndex()
125  ]
126  );
127 
128  const mappedPatchBase& mpp = refCast<const mappedPatchBase>(pp_);
129 
130  mpp.distribute(absorptivity);
131 
132  // Restore tag
133  UPstream::msgType() = oldTag;
134 
135  return tmp<scalarField>::New(std::move(absorptivity));
136 }
137 
138 
140 (
141  const label faceI,
142  const label bandI,
143  const vector dir,
144  const scalar T
145 ) const
146 {
147  return a(bandI, nullptr, nullptr)()[faceI];
148 }
149 
151 (
152  const label bandI,
153  vectorField* incomingDirection,
154  scalarField* T
155 ) const
156 {
157 
158  // Since we're inside initEvaluate/evaluate there might be processor
159  // comms underway. Change the tag we use.
160  int oldTag = UPstream::msgType();
161  UPstream::msgType() = oldTag+1;
162 
163  const fvMesh& nbrMesh = nbrRegion();
164 
167  (
168  "radiationProperties"
169  );
170 
171  scalarField emissivity
172  (
173  radiation.absorptionEmission().e
174  (
175  bandI
176  )().boundaryField()
177  [
178  nbrPatchIndex()
179  ]
180  );
181 
182  const mappedPatchBase& mpp = refCast<const mappedPatchBase>(pp_);
183 
184  mpp.distribute(emissivity);
185 
186  // Restore tag
187  UPstream::msgType() = oldTag;
188 
189  return tmp<scalarField>::New(std::move(emissivity));
190 }
191 
192 
194 (
195  const label faceI,
196  const label bandI,
197  const vector dir,
198  const scalar T
199 ) const
200 {
201  return e(bandI, nullptr, nullptr)()[faceI];
202 }
203 
204 
206 {
207  const fvMesh& nbrMesh = nbrRegion();
208 
211  (
212  "radiationProperties"
213  );
214 
215  return (radiation.absorptionEmission().nBands());
216 }
217 // ************************************************************************* //
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::radiation::solidAbsorption::solidAbsorption
solidAbsorption(const dictionary &dict, const polyPatch &pp)
Construct from components.
Definition: solidAbsorption.C:71
Foam::FatalIOError
IOerror FatalIOError
Foam::radiation::solidAbsorption::~solidAbsorption
virtual ~solidAbsorption()
Destructor.
Definition: solidAbsorption.C:91
Foam::radiation::solidAbsorption::e
tmp< scalarField > e(const label bandI=0, vectorField *incomingDirection=nullptr, scalarField *T=nullptr) const
Return emission coefficient.
Definition: solidAbsorption.C:151
Foam::radiation::wallAbsorptionEmissionModel::pp_
const polyPatch & pp_
Reference to the polyPatch.
Definition: wallAbsorptionEmissionModel.H:60
Foam::mappedPatchBase
Determines a mapping between patch face centres and mesh cell or face centres and processors they're ...
Definition: mappedPatchBase.H:112
Foam::Field< vector >
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:68
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
solidAbsorption.H
Foam::objectRegistry::lookupObject
const Type & lookupObject(const word &name, const bool recursive=false) const
Definition: objectRegistryTemplates.C:434
radiation
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::radiation::wallAbsorptionEmissionModel
Based class for wall absorption emission models.
Definition: wallAbsorptionEmissionModel.H:52
Foam::radiation::defineTypeNameAndDebug
defineTypeNameAndDebug(cloudAbsorptionEmission, 0)
Foam::radiation::solidAbsorption::a
tmp< scalarField > a(const label bandI=0, vectorField *incomingDirection=nullptr, scalarField *T=nullptr) const
absorptivity coefficient
Definition: solidAbsorption.C:98
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::radiation::solidAbsorption::nBands
label nBands() const
Number of bands.
Definition: solidAbsorption.C:205
Foam::UPstream::msgType
static int & msgType() noexcept
Message tag of standard messages.
Definition: UPstream.H:540
absorptionEmissionModel.H
Foam::Vector< scalar >
Foam::radiation::radiationModel
Top level model for radiation modelling.
Definition: radiationModel.H:75
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::tmp::New
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
Foam::patchIdentifier::name
const word & name() const noexcept
The patch name.
Definition: patchIdentifier.H:135
radiationModel.H
Foam::radiation::addToRunTimeSelectionTable
addToRunTimeSelectionTable(absorptionEmissionModel, cloudAbsorptionEmission, dictionary)
mappedPatchBase.H