solarHeatLoad.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::solarHeatLoad
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 solarHeatLoad 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 "solarHeatLoad" must be set to
51  true.
52 
53 
54 SourceFiles
55  solarHeatLoad.C
56 
57 \*---------------------------------------------------------------------------*/
58 
59 #ifndef radiationModelsolarHeatLoad_H
60 #define radiationModelsolarHeatLoad_H
61 
62 #include "radiationModel.H"
63 #include "singleCellFvMesh.H"
64 #include "scalarListIOList.H"
65 #include "volFields.H"
66 #include "faceHeatShading.H"
67 #include "solarCalculator.H"
68 #include "IOmapDistribute.H"
69 
70 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
71 
72 namespace Foam
73 {
74 namespace radiation
75 {
76 
77 /*---------------------------------------------------------------------------*\
78  Class solarHeatLoad Declaration
79 \*---------------------------------------------------------------------------*/
80 
81 class solarHeatLoad
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  // Private data
95 
96 
97  //- Agglomeration List
99 
100  //- Coarse mesh
102 
103  //- direct solar radiative heat flux [W/m2] 净辐射量,直射辐射
105 
106  //- reflect solar radiative heat flux [W/m2] 面面之间的辐射量,中间量,不用保存
108 
109  //- diffuse solar radiative heat flux [W/m2] 散射辐射量
111 
112  //- 总辐射量total solar radiative heat flux [W/m2] 总的辐射量
114 
115  //- Direct hit faces Ids
117 
118  //- Constant source term
120 
121  //- Solar calculator
123 
124  //- Vertical direction (Default is g vector)
126 
127  //- Include diffuse reflected heat fluxes from direct heat flux
128  bool useVFbeamToDiffuse_;
129 
130  //- Selected patches to apply solar load
132 
133  //- Chached coarse to fine mapping for coarse mesh
135 
136  //-Number of bands
137  label nBands_;
138 
139  //- Spectral distribution for the integrated solar heat flux
141 
142  //- Map distribute
144 
145  //- Face-compact map
147 
148  //- Couple solids through mapped boundary patch using Qr (default:true)
149  bool solidCoupled_;
150 
151  //- Couple wall patches using Qr (default:false)
152  bool wallCoupled_;
153 
154  //- Absorptivity list
156 
157  //- Update absorptivity
158  bool updateAbsorptivity_;
159 
160  //- First iteration
161  bool firstIter_;
162 
163  //- Update Sun position index
165 
166 
167  // Private Member Functions
168 
169 
170  //- Initialise
171  void initialise(const dictionary&);
172 
173  //- Update direct hit faces radiation
174  void updateDirectHitRadiation(const labelList&, const labelHashSet&);
175 
176  //- Calculate diffusive heat flux
177  void calculateQdiff(const labelHashSet&, const labelHashSet&);
178 
179  //- Update Sky diffusive radiation
181  (
182  const labelHashSet&,
183  const labelHashSet&
184  );
185 
186  //- Update hit faces
187  bool updateHitFaces();
188 
189  //- Update absorptivity
190  void updateAbsorptivity(const labelHashSet& includePatches);
191 
192  //- Disallow default bitwise copy construct
194 
195  //- Disallow default bitwise assignment
196  void operator=(const solarHeatLoad&);
197 
198 
199 public:
200 
201  //- Runtime type information
202  TypeName("solarHeatLoad");
203 
204 
205  // Constructors
206 
207  //- Construct from volScalarField
209 
210  //- Construct from dictionary and volScalarField
211  solarHeatLoad(const dictionary& dict, const volScalarField& T);
212 
213  //- Constructor from local components. Does not create a radiationModel.
214  // radWallFieldName is the solar heat field name
216  (
217  const dictionary& dict,
218  const volScalarField& T,
219  const word radWallFieldName
220  );
221 
222 
223  //- Destructor
224  virtual ~solarHeatLoad();
225 
226 
227  // Member functions
228 
229  // Edit
230 
231  //- Solve
232  void calculate();
233 
234  //- Read radiation properties dictionary
235  bool read();
236 
237  //- Source term component (for power of T^4)
238  virtual tmp<volScalarField> Rp() const;
239 
240  //- Source term component (constant)
241  virtual tmp<DimensionedField<scalar, volMesh> > Ru() const;
242 
243 
244  // Access
245 
246  //- Number of bands
247  label nBands() const;
248 };
249 
250 
251 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
252 
253 } // End namespace radiation
254 } // End namespace Foam
255 
256 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
257 
258 #endif
259 
260 // ************************************************************************* //
Foam::radiation::solarHeatLoad::calculate
void calculate()
Solve.
Definition: solarHeatLoad.C:1054
Foam::radiation::solarHeatLoad::spectralDistribution_
scalarList spectralDistribution_
Spectral distribution for the integrated solar heat flux.
Definition: solarHeatLoad.H:139
volFields.H
Foam::radiation::solarHeatLoad::absorptivity_
List< List< tmp< scalarField > > > absorptivity_
Absorptivity list.
Definition: solarHeatLoad.H:154
Foam::radiation::solarHeatLoad::updateDirectHitRadiation
void updateDirectHitRadiation(const labelList &, const labelHashSet &)
Update direct hit faces radiation.
Definition: solarHeatLoad.C:121
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::solarHeatLoad::finalAgglom_
labelListIOList finalAgglom_
Agglomeration List.
Definition: solarHeatLoad.H:97
Foam::radiation::solarHeatLoad::solarCalc_
solarCalculator solarCalc_
Solar calculator.
Definition: solarHeatLoad.H:121
Foam::radiation::solarHeatLoad::firstIter_
bool firstIter_
First iteration.
Definition: solarHeatLoad.H:160
Foam::radiation::solarHeatLoad::Rp
virtual tmp< volScalarField > Rp() const
Source term component (for power of T^4)
Definition: solarHeatLoad.C:1131
Foam::radiation::solarHeatLoad::coarseToFine_
List< labelListList > coarseToFine_
Chached coarse to fine mapping for coarse mesh.
Definition: solarHeatLoad.H:133
faceHeatShading.H
Foam::radiation::solarHeatLoad::visibleFaceFaces_
labelListIOList visibleFaceFaces_
Face-compact map.
Definition: solarHeatLoad.H:145
Foam::radiation::solarHeatLoad::useVFbeamToDiffuse_
bool useVFbeamToDiffuse_
Include diffuse reflected heat fluxes from direct heat flux.
Definition: solarHeatLoad.H:127
Foam::HashSet< label, Hash< label > >
Foam::radiation::solarHeatLoad::TypeName
TypeName("solarHeatLoad")
Runtime type information.
Foam::radiation::solarHeatLoad::updateSkyDiffusiveRadiation
void updateSkyDiffusiveRadiation(const labelHashSet &, const labelHashSet &)
Update Sky diffusive radiation.
Definition: solarHeatLoad.C:183
Foam::solarCalculator
The solar calculator model provides information about the Sun direction and Sun load model....
Definition: solarCalculator.H:101
Foam::radiation::solarHeatLoad::verticalDir_
vector verticalDir_
Vertical direction (Default is g vector)
Definition: solarHeatLoad.H:124
Foam::radiation::solarHeatLoad::~solarHeatLoad
virtual ~solarHeatLoad()
Destructor.
Definition: solarHeatLoad.C:1029
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::solarHeatLoad::Ru
virtual tmp< DimensionedField< scalar, volMesh > > Ru() const
Source term component (constant)
Definition: solarHeatLoad.C:1159
Foam::radiation::solarHeatLoad::QdiffRad_
volScalarField QdiffRad_
diffuse solar radiative heat flux [W/m2] 散射辐射量
Definition: solarHeatLoad.H:109
Foam::radiation::solarHeatLoad::nBands_
label nBands_
Number of bands.
Definition: solarHeatLoad.H:136
Foam::radiation::solarHeatLoad::nBands
label nBands() const
Number of bands.
Definition: solarHeatLoad.C:1035
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:55
singleCellFvMesh.H
radiation
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::radiation::solarHeatLoad::updateAbsorptivity
void updateAbsorptivity(const labelHashSet &includePatches)
Update absorptivity.
Definition: solarHeatLoad.C:100
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::radiation::solarHeatLoad::operator=
void operator=(const solarHeatLoad &)
Disallow default bitwise assignment.
Foam::radiation::solarHeatLoad::viewFactorWalls
static const word viewFactorWalls
Static name for view factor walls.
Definition: solarHeatLoad.H:89
Foam::radiation::solarHeatLoad::map_
autoPtr< IOmapDistribute > map_
Map distribute.
Definition: solarHeatLoad.H:142
Foam::radiation::solarHeatLoad::updateHitFaces
bool updateHitFaces()
Update hit faces.
Definition: solarHeatLoad.C:52
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::radiation::solarHeatLoad
The solar load radiation model includes Sun primary hits, their reflective fluxes and diffusive sky r...
Definition: solarHeatLoad.H:80
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
IOmapDistribute.H
scalarListIOList.H
Foam::radiation::solarHeatLoad::solarHeatLoad
solarHeatLoad(const solarHeatLoad &)
Disallow default bitwise copy construct.
Foam::Vector< scalar >
Foam::radiation::solarHeatLoad::Qr_
volScalarField Qr_
direct solar radiative heat flux [W/m2] 净辐射量,直射辐射
Definition: solarHeatLoad.H:103
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::radiationModel
Top level model for radiation modelling.
Definition: radiationModel.H:73
Foam::radiation::solarHeatLoad::updateAbsorptivity_
bool updateAbsorptivity_
Update absorptivity.
Definition: solarHeatLoad.H:157
Foam::radiation::radiationModel::initialise
void initialise()
Initialise.
Definition: radiationModel.C:77
solarCalculator.H
Foam::IOList< labelList >
Foam::radiation::solarHeatLoad::solidCoupled_
bool solidCoupled_
Couple solids through mapped boundary patch using Qr (default:true)
Definition: solarHeatLoad.H:148
Foam::radiation::solarHeatLoad::Ru_
DimensionedField< scalar, volMesh > Ru_
Constant source term.
Definition: solarHeatLoad.H:118
Foam::radiation::solarHeatLoad::coarseMesh_
autoPtr< singleCellFvMesh > coarseMesh_
Coarse mesh.
Definition: solarHeatLoad.H:100
Foam::radiation::solarHeatLoad::hitFaces_
autoPtr< faceHeatShading > hitFaces_
Direct hit faces Ids.
Definition: solarHeatLoad.H:115
Foam::radiation::solarHeatLoad::QtotalRad_
volScalarField QtotalRad_
总辐射量total solar radiative heat flux [W/m2] 总的辐射量
Definition: solarHeatLoad.H:112
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::radiation::solarHeatLoad::QrefecFormOtherWallRad_
volScalarField QrefecFormOtherWallRad_
reflect solar radiative heat flux [W/m2] 面面之间的辐射量,中间量,不用保存
Definition: solarHeatLoad.H:106
Foam::radiation::solarHeatLoad::read
bool read()
Read radiation properties dictionary.
Definition: solarHeatLoad.C:1041
Foam::radiation::solarHeatLoad::includePatches_
labelList includePatches_
Selected patches to apply solar load.
Definition: solarHeatLoad.H:130
radiationModel.H
Foam::radiation::solarHeatLoad::updateTimeIndex_
label updateTimeIndex_
Update Sun position index.
Definition: solarHeatLoad.H:163
Foam::radiation::solarHeatLoad::calculateQdiff
void calculateQdiff(const labelHashSet &, const labelHashSet &)
Calculate diffusive heat flux.
Definition: solarHeatLoad.C:390
Foam::radiation::solarHeatLoad::wallCoupled_
bool wallCoupled_
Couple wall patches using Qr (default:false)
Definition: solarHeatLoad.H:151
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:51