solarLoad.H
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 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 Class
25  Foam::radiation::solarLoad
26 
27 Description
28 
29  The solar load radiation model includes Sun primary hits, their
30  reflective fluxes and diffusive sky radiative fluxes.
31 
32  The primary hit rays are calculated using a face shading algorithm.
33  The reflected fluxes are considered diffusive and use a view factors method
34  to deposit the energy on "visible" walls. The sky diffusive radiation for
35  horizontal and vertical walls is calculated following the Fair Weather
36  Conditions Method from the ASHRAE Handbook.
37 
38  By default the energy is included in cells adjacent to the patches into
39  the energy Equation (wallCoupled = false). On coupled patches the flux is
40  by default added to the wall and considered into the solid
41  (solidCoupled = true).
42 
43  The reflected fluxes uses a grey absoprtion/emission model wich is weighted
44  by the spectral distribution. The flag useVFbeamToDiffuse should be
45  switched on and the view factors should be calculated using the
46  'viewFactorsGen' application.
47 
48  The solarLoad model can be used in conjuntion with fvDOM and viewFactor
49  radiation models but only using a single band spectrum. On the
50  corresponding BC's for these models the flag "solarLoad" must be set to
51  true.
52 
53 
54 SourceFiles
55  solarLoad.C
56 
57 \*---------------------------------------------------------------------------*/
58 
59 #ifndef radiationModelsolarLoad_H
60 #define radiationModelsolarLoad_H
61 
62 #include "radiationModel.H"
63 #include "singleCellFvMesh.H"
64 #include "scalarListIOList.H"
65 #include "volFields.H"
66 #include "faceShading.H"
67 #include "solarCalculator.H"
68 #include "IOmapDistribute.H"
69 
70 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
71 
72 namespace Foam
73 {
74 namespace radiation
75 {
76 
77 /*---------------------------------------------------------------------------*\
78  Class solarLoad Declaration
79 \*---------------------------------------------------------------------------*/
80 
81 class solarLoad
82 :
83  public radiationModel
84 {
85 protected:
86 
87  // Static data
88 
89  //- Static name for view factor walls
90  static const word viewFactorWalls;
91 
92 
93 private:
94 
95  // Private data
96 
97 
98  //- Agglomeration List
100 
101  //- Coarse mesh
103 
104  //- Net radiative heat flux [W/m2]
106 
107  //- Secondary solar radiative heat flux [W/m2]
109 
110  //- Direct hit faces Ids
112 
113  //- Constant source term
115 
116  //- Solar calculator
118 
119  //- Vertical direction (Default is g vector)
121 
122  //- Include diffuse reflected heat fluxes from direct heat flux
123  bool useVFbeamToDiffuse_;
124 
125  //- Selected patches to apply solar load
127 
128  //- Chached coarse to fine mapping for coarse mesh
130 
131  //-Number of bands
132  label nBands_;
133 
134  //- Spectral distribution for the integrated solar heat flux
136 
137  //- Map distribute
139 
140  //- Face-compact map
142 
143  //- Couple solids through mapped boundary patch using Qr (default:true)
144  bool solidCoupled_;
145 
146  //- Couple wall patches using Qr (default:false)
147  bool wallCoupled_;
148 
149  //- Absorptivity list
151 
152  //- Update absorptivity
153  bool updateAbsorptivity_;
154 
155  //- First iteration
156  bool firstIter_;
157 
158  //- Update Sun position index
160 
161 
162  // Private Member Functions
163 
164 
165  //- Initialise
166  void initialise(const dictionary&);
167 
168  //- Update direct hit faces radiation
169  void updateDirectHitRadiation(const labelList&, const labelHashSet&);
170 
171  //- Calculate diffusive heat flux
172  void calculateQdiff(const labelHashSet&, const labelHashSet&);
173 
174  //- Update Sky diffusive radiation
176  (
177  const labelHashSet&,
178  const labelHashSet&
179  );
180 
181  //- Update hit faces
182  bool updateHitFaces();
183 
184  //- Update absorptivity
185  void updateAbsorptivity(const labelHashSet& includePatches);
186 
187  //- Disallow default bitwise copy construct
188  solarLoad(const solarLoad&);
189 
190  //- Disallow default bitwise assignment
191  void operator=(const solarLoad&);
192 
193 
194 public:
195 
196  //- Runtime type information
197  TypeName("solarLoad");
198 
199 
200  // Constructors
201 
202  //- Construct from volScalarField
203  solarLoad(const volScalarField& T);
204 
205  //- Construct from dictionary and volScalarField
206  solarLoad(const dictionary& dict, const volScalarField& T);
207 
208  //- Constructor from local components. Does not create a radiationModel.
209  // radWallFieldName is the solar heat field name
210  solarLoad
211  (
212  const dictionary& dict,
213  const volScalarField& T,
214  const word radWallFieldName
215  );
216 
217 
218  //- Destructor
219  virtual ~solarLoad();
220 
221 
222  // Member functions
223 
224  // Edit
225 
226  //- Solve
227  void calculate();
228 
229  //- Read radiation properties dictionary
230  bool read();
231 
232  //- Source term component (for power of T^4)
233  virtual tmp<volScalarField> Rp() const;
234 
235  //- Source term component (constant)
236  virtual tmp<DimensionedField<scalar, volMesh> > Ru() const;
237 
238 
239  // Access
240 
241  //- Number of bands
242  label nBands() const;
243 };
244 
245 
246 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
247 
248 } // End namespace radiation
249 } // End namespace Foam
250 
251 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
252 
253 #endif
254 
255 // ************************************************************************* //
Foam::radiation::solarLoad::nBands_
label nBands_
Number of bands.
Definition: solarLoad.H:131
volFields.H
Foam::radiation::solarLoad::nBands
label nBands() const
Number of bands.
Definition: solarLoad.C:921
Foam::radiation::solarLoad::updateAbsorptivity
void updateAbsorptivity(const labelHashSet &includePatches)
Update absorptivity.
Definition: solarLoad.C:94
Foam::radiation::solarLoad::coarseToFine_
List< labelListList > coarseToFine_
Chached coarse to fine mapping for coarse mesh.
Definition: solarLoad.H:128
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::radiation::solarLoad::firstIter_
bool firstIter_
First iteration.
Definition: solarLoad.H:155
Foam::radiation::solarLoad::updateHitFaces
bool updateHitFaces()
Update hit faces.
Definition: solarLoad.C:52
Foam::radiation::solarLoad::solidCoupled_
bool solidCoupled_
Couple solids through mapped boundary patch using Qr (default:true)
Definition: solarLoad.H:143
Foam::radiation::solarLoad::updateAbsorptivity_
bool updateAbsorptivity_
Update absorptivity.
Definition: solarLoad.H:152
Foam::radiation::solarLoad::solarCalc_
solarCalculator solarCalc_
Solar calculator.
Definition: solarLoad.H:116
Foam::HashSet< label, Hash< label > >
Foam::radiation::solarLoad::Ru
virtual tmp< DimensionedField< scalar, volMesh > > Ru() const
Source term component (constant)
Definition: solarLoad.C:1038
Foam::radiation::solarLoad::calculateQdiff
void calculateQdiff(const labelHashSet &, const labelHashSet &)
Calculate diffusive heat flux.
Definition: solarLoad.C:358
Foam::radiation::solarLoad::TypeName
TypeName("solarLoad")
Runtime type information.
Foam::radiation::solarLoad::absorptivity_
List< List< tmp< scalarField > > > absorptivity_
Absorptivity list.
Definition: solarLoad.H:149
Foam::solarCalculator
The solar calculator model provides information about the Sun direction and Sun load model....
Definition: solarCalculator.H:101
Foam::radiation::solarLoad::Ru_
DimensionedField< scalar, volMesh > Ru_
Constant source term.
Definition: solarLoad.H:113
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::radiation::solarLoad::Qr_
volScalarField Qr_
Net radiative heat flux [W/m2].
Definition: solarLoad.H:104
Foam::radiation::solarLoad::QsecondRad_
volScalarField QsecondRad_
Secondary solar radiative heat flux [W/m2].
Definition: solarLoad.H:107
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:55
Foam::radiation::solarLoad::updateSkyDiffusiveRadiation
void updateSkyDiffusiveRadiation(const labelHashSet &, const labelHashSet &)
Update Sky diffusive radiation.
Definition: solarLoad.C:163
faceShading.H
singleCellFvMesh.H
Foam::radiation::solarLoad::verticalDir_
vector verticalDir_
Vertical direction (Default is g vector)
Definition: solarLoad.H:119
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 any number of values (e....
Definition: dictionary.H:137
Foam::radiation::solarLoad::includePatches_
labelList includePatches_
Selected patches to apply solar load.
Definition: solarLoad.H:125
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::radiation::solarLoad::calculate
void calculate()
Solve.
Definition: solarLoad.C:940
Foam::radiation::solarLoad::operator=
void operator=(const solarLoad &)
Disallow default bitwise assignment.
Foam::radiation::solarLoad::solarLoad
solarLoad(const solarLoad &)
Disallow default bitwise copy construct.
Foam::radiation::solarLoad::wallCoupled_
bool wallCoupled_
Couple wall patches using Qr (default:false)
Definition: solarLoad.H:146
Foam::autoPtr
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
Foam::radiation::solarLoad::map_
autoPtr< IOmapDistribute > map_
Map distribute.
Definition: solarLoad.H:137
IOmapDistribute.H
Foam::radiation::solarLoad
The solar load radiation model includes Sun primary hits, their reflective fluxes and diffusive sky r...
Definition: solarLoad.H:80
scalarListIOList.H
Foam::Vector< scalar >
Foam::radiation::solarLoad::~solarLoad
virtual ~solarLoad()
Destructor.
Definition: solarLoad.C:915
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::radiation::solarLoad::updateTimeIndex_
label updateTimeIndex_
Update Sun position index.
Definition: solarLoad.H:158
Foam::radiation::radiationModel
Top level model for radiation modelling.
Definition: radiationModel.H:73
Foam::radiation::radiationModel::initialise
void initialise()
Initialise.
Definition: radiationModel.C:77
solarCalculator.H
Foam::IOList< labelList >
Foam::radiation::solarLoad::hitFaces_
autoPtr< faceShading > hitFaces_
Direct hit faces Ids.
Definition: solarLoad.H:110
Foam::radiation::solarLoad::visibleFaceFaces_
labelListIOList visibleFaceFaces_
Face-compact map.
Definition: solarLoad.H:140
Foam::radiation::solarLoad::viewFactorWalls
static const word viewFactorWalls
Static name for view factor walls.
Definition: solarLoad.H:89
Foam::radiation::solarLoad::useVFbeamToDiffuse_
bool useVFbeamToDiffuse_
Include diffuse reflected heat fluxes from direct heat flux.
Definition: solarLoad.H:122
Foam::radiation::solarLoad::finalAgglom_
labelListIOList finalAgglom_
Agglomeration List.
Definition: solarLoad.H:98
Foam::radiation::solarLoad::read
bool read()
Read radiation properties dictionary.
Definition: solarLoad.C:927
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::radiation::solarLoad::updateDirectHitRadiation
void updateDirectHitRadiation(const labelList &, const labelHashSet &)
Update direct hit faces radiation.
Definition: solarLoad.C:115
Foam::radiation::solarLoad::coarseMesh_
autoPtr< singleCellFvMesh > coarseMesh_
Coarse mesh.
Definition: solarLoad.H:101
Foam::radiation::solarLoad::spectralDistribution_
scalarList spectralDistribution_
Spectral distribution for the integrated solar heat flux.
Definition: solarLoad.H:134
radiationModel.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::solarLoad::Rp
virtual tmp< volScalarField > Rp() const
Source term component (for power of T^4)
Definition: solarLoad.C:1010