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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2015 OpenFOAM Foundation
9 Copyright (C) 2018-2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27Class
28 Foam::radiation::solarLoad
29
30Group
31 grpRadiationModels
32
33Description
34 The \c solarLoad radiation model includes Sun primary hits, their
35 reflective fluxes and diffusive sky radiative fluxes.
36
37 The primary hit rays are calculated using a face shading algorithm.
38 The first reflected fluxes can be optionally included. A view factors
39 method is needed in order to include diffusive surface to surface fluxes.
40
41 The energy is included on "visible" walls by default. The sky diffusive
42 radiation for horizontal and vertical walls is calculated following the
43 Fair Weather Conditions Method from the ASHRAE Handbook.
44
45 By default the energy is included in cells adjacent to the patches into
46 the energy equation (\c wallCoupled=false). On coupled patches the flux is
47 by default added to the wall and considered into the solid
48 (\c solidCoupled=true).
49
50 The \c solarLoad model can be used in conjuntion with \c fvDOM and
51 \c viewFactor radiation models. The flag \c useSolarLoad must be
52 \c true on the \c radiationProperties dictionary.
53
54Usage
55 Minimal examples by using \c constant/radiationProperties:
56
57 \verbatim
58 solarLoadCoeffs
59 {
60 // Mandatory entries
61 useReflectedRays true;
62 spectralDistribution (1 5 1 2);
63
64 // Optional entries
65 solidCoupled true;
66 wallCoupled false;
67 updateAbsorptivity true;
68
69 // Mandatory/Optional (inherited) entries
70 ...
71 }
72 \endverbatim
73
74 where the entries mean:
75 \table
76 Property | Description | Type | Reqd | Deflt
77 useReflectedRays | Flag to use reflected rays | bool | yes | -
78 spectralDistribution | Spectral distribution for the integrated <!--
79 --> solar heat flux | Function1<scalarField> | yes | -
80 solidCoupled | Flag to couple solids through mapped <!--
81 --> boundary patch using qr | bool | no | true
82 wallCoupled | Flag to couple wall patches using qr <!--
83 --> | bool | no | false
84 updateAbsorptivity | Flag to enable absorptivity updates <!--
85 --> | bool | no | false
86 \endtable
87
88 The inherited entries are elaborated in:
89 - \link radiationModel.H \endlink
90 - \link solarCalculator.H \endlink
91 - \link Function1.H \endlink
92
93SourceFiles
94 solarLoad.C
95
96\*---------------------------------------------------------------------------*/
97
98#ifndef radiation_solarLoad_H
99#define radiation_solarLoad_H
100
101#include "radiationModel.H"
102#include "volFields.H"
103#include "faceShading.H"
104#include "faceReflecting.H"
105#include "solarCalculator.H"
106#include "Function1.H"
107
108// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
109
110namespace Foam
111{
112namespace radiation
113{
114
115/*---------------------------------------------------------------------------*\
116 Class solarLoad Declaration
117\*---------------------------------------------------------------------------*/
118
119class solarLoad
120:
121 public radiationModel
122{
123 // Private Data
124
125 //- Solar calculator
126 solarCalculator solarCalc_;
127
128 //- Dictionary
129 dictionary dict_;
130
131 //- Net radiative heat flux [W/m2]
132 volScalarField qr_;
133
134 //- Direct hit faces Ids
135 autoPtr<faceShading> hitFaces_;
136
137 //- Reflected faces
138 autoPtr<faceReflecting> reflectedFaces_;
139
140 //- Source term for cells next to patches with flags solidCoupled
141 //- and wallCoupled false
142 DimensionedField<scalar, volMesh> Ru_;
143
144 //- Absorptivity list
145 List<List<tmp<scalarField>>> absorptivity_;
146
147 //- Spectral distribution for the integrated solar heat flux
148 scalarList spectralDistribution_;
149
150 //- Time-dependent spectral distributions
151 autoPtr<Function1<scalarField>> spectralDistributions_;
152
153 //- Primary solar radiative heat flux per band [W/m2]
154 PtrList<volScalarField> qprimaryRad_;
155
156 //- Vertical direction (default is g vector)
157 vector verticalDir_;
158
159 //- Number of bands
160 label nBands_;
161
162 //- Update Sun position index
163 label updateTimeIndex_;
164
165 //- Couple solids through mapped boundary patch using qr
166 bool solidCoupled_;
167
168 //- Couple wall patches using qr
169 bool wallCoupled_;
170
171 //- Update absorptivity
172 bool updateAbsorptivity_;
173
174 //- Include reflected rays from specular surfaces
175 bool useReflectedRays_;
176
177 //- First iteration
178 bool firstIter_;
179
180
181 // Private Member Functions
182
183 //- Initialise model parameters
184 void initialise(const dictionary&);
185
186 //- Update direct hit faces radiation
187 void updateDirectHitRadiation(const labelList&, const labelHashSet&);
188
189 //- Update reflected heat flux
190 void updateReflectedRays(const labelHashSet&);
191
192 //- Calculate diffusive heat flux
193 //void calculateQdiff(const labelHashSet&, const labelHashSet&);
194
195 //- Update Sky diffusive radiation
196 void updateSkyDiffusiveRadiation
197 (
198 const labelHashSet&,
199 const labelHashSet&
200 );
201
202 //- Update hit faces
203 bool updateHitFaces();
204
205 //- Update absorptivity
206 void updateAbsorptivity(const labelHashSet& includePatches);
207
208
209 //- No copy construct
210 solarLoad(const solarLoad&) = delete;
211
212 //- No copy assignment
213 void operator=(const solarLoad&) = delete;
214
215
216public:
217
218 //- Runtime type information
219 TypeName("solarLoad");
220
221
222 // Constructors
223
224 //- Construct from volScalarField
225 solarLoad(const volScalarField& T);
226
227 //- Construct from dictionary and volScalarField
228 solarLoad(const dictionary& dict, const volScalarField& T);
229
230
231 //- Destructor
232 virtual ~solarLoad() = default;
233
234
235 // Member Functions
236
237 // Evaluation
238
239 //- Read radiationProperties dictionary
240 bool read();
241
242 //- Solve radiation equations
243 void calculate();
244
245
246 // Access
247
248 //- Source term component (for power of T^4)
249 virtual tmp<volScalarField> Rp() const;
250
251 //- Source term component (constant)
253
254 //- Return const access to the number of bands
255 label nBands() const noexcept
256 {
257 return nBands_;
258 }
259
260 //- Return const access to the primary solar heat flux
261 const volScalarField& qprimaryRad(const label bandI) const
263 return qprimaryRad_[bandI];
264 }
265};
266
267
268// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
269
270} // End namespace radiation
271} // End namespace Foam
272
273// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
274
275#endif
276
277// ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:77
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
dictionary()
Default construct, a top-level empty dictionary.
Definition: dictionary.C:75
Top level model for radiation modelling.
The solarLoad radiation model includes Sun primary hits, their reflective fluxes and diffusive sky ra...
Definition: solarLoad.H:165
TypeName("solarLoad")
Runtime type information.
virtual tmp< DimensionedField< scalar, volMesh > > Ru() const
Source term component (constant)
Definition: solarLoad.C:950
virtual tmp< volScalarField > Rp() const
Source term component (for power of T^4)
Definition: solarLoad.C:926
const volScalarField & qprimaryRad(const label bandI) const
Return const access to the primary solar heat flux.
Definition: solarLoad.H:304
label nBands() const noexcept
Return const access to the number of bands.
Definition: solarLoad.H:298
virtual ~solarLoad()=default
Destructor.
bool read()
Read radiationProperties dictionary.
Definition: solarLoad.C:834
void calculate()
Solve radiation equations.
Definition: solarLoad.C:845
A solar calculator model providing models for the solar direction and solar loads.
A class for managing temporary objects.
Definition: tmp.H:65
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
const volScalarField & T
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
Namespace for OpenFOAM.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:64
const direction noexcept
Definition: Scalar.H:223
dictionary dict
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73