Loading...
Searching...
No Matches
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 "solarLoadBase.H"
103#include "volFields.H"
104#include "faceShading.H"
105#include "faceReflecting.H"
106#include "solarCalculator.H"
107#include "Function1.H"
108
109// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
110
111namespace Foam
112{
113namespace radiation
114{
115
116/*---------------------------------------------------------------------------*\
117 Class solarLoad Declaration
118\*---------------------------------------------------------------------------*/
119
120class solarLoad
121:
122 public radiationModel,
123 public solarLoadBase
124{
125 // Private Data
126
127 //- Solar calculator
128 solarCalculator solarCalc_;
129
130 //- Dictionary
131 dictionary dict_;
132
133 //- Net radiative heat flux [W/m2]
134 volScalarField qr_;
135
136 //- Direct hit faces Ids
137 autoPtr<faceShading> hitFaces_;
138
139 //- Reflected faces
140 autoPtr<faceReflecting> reflectedFaces_;
141
142 //- Source term for cells next to patches with flags solidCoupled
143 //- and wallCoupled false
144 DimensionedField<scalar, volMesh> Ru_;
145
146 //- Absorptivity list
147 List<List<tmp<scalarField>>> absorptivity_;
148
149 //- Spectral distribution for the integrated solar heat flux
150 scalarList spectralDistribution_;
151
152 //- Time-dependent spectral distributions
153 autoPtr<Function1<scalarField>> spectralDistributions_;
154
155 //- Primary solar radiative heat flux per band [W/m2]
156 PtrList<volScalarField> qprimaryRad_;
157
158 //- Number of bands
159 label nBands_;
160
161 //- Update Sun position index
162 label updateTimeIndex_;
163
164 //- Couple solids through mapped boundary patch using qr
165 bool solidCoupled_;
166
167 //- Couple wall patches using qr
168 bool wallCoupled_;
169
170 //- Update absorptivity
171 bool updateAbsorptivity_;
172
173 //- Include reflected rays from specular surfaces
174 bool useReflectedRays_;
175
176 //- First iteration
177 bool firstIter_;
178
179
180 // Private Member Functions
181
182 //- Initialise model parameters
183 void initialise(const dictionary&);
184
185 //- Update direct hit faces radiation
186 void updateDirectHitRadiation(const labelList&, const labelHashSet&);
187
188 //- Update reflected heat flux
189 void updateReflectedRays(const labelHashSet&);
190
191 //- Calculate diffusive heat flux
192 //void calculateQdiff(const labelHashSet&, const labelHashSet&);
193
194 //- Update Sky diffusive radiation
195 void updateSkyDiffusiveRadiation
196 (
197 const labelHashSet&,
198 const labelHashSet&
199 );
200
201 //- Update hit faces
202 bool updateHitFaces();
203
204 //- Update absorptivity
205 void updateAbsorptivity(const labelHashSet& includePatches);
206
207
208 //- No copy construct
209 solarLoad(const solarLoad&) = delete;
210
211 //- No copy assignment
212 void operator=(const solarLoad&) = delete;
213
214
215public:
216
217 //- Runtime type information
218 TypeName("solarLoad");
219
220
221 // Constructors
222
223 //- Construct from volScalarField
224 solarLoad(const volScalarField& T);
225
226 //- Construct from dictionary and volScalarField
227 solarLoad(const dictionary& dict, const volScalarField& T);
228
229
230 //- Destructor
231 virtual ~solarLoad() = default;
232
233
234 // Member Functions
235
236 // Evaluation
237
238 //- Read radiationProperties dictionary
239 bool read();
240
241 //- Solve radiation equations
242 void calculate();
243
244
245 // Access
246
247 //- Source term component (for power of T^4)
248 virtual tmp<volScalarField> Rp() const;
249
250 //- Source term component (constant)
252
253 //- Return const access to the number of bands
254 label nBands() const noexcept
255 {
256 return nBands_;
257 }
258
259 //- Return const access to the primary solar heat flux
260 const volScalarField& qprimaryRad(const label bandI) const
261 {
262 return qprimaryRad_[bandI];
263 }
264
265 //- Return const reference to the solar calculator
266 virtual const solarCalculator& solarCalculatorRef() const;
267
268 //- Return const reference to the face shading calculator
269 virtual const faceShading& faceShadingRef() const;
270};
271
272
273// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
274
275} // End namespace radiation
276} // End namespace Foam
277
278// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
279
280#endif
281
282// ************************************************************************* //
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:72
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
dictionary()
Default construct, a top-level empty dictionary.
Definition dictionary.C:68
Helper class to calculate visible faces for global, sun-like illumination.
Definition faceShading.H:58
Top level model for radiation modelling.
const volScalarField & T() const noexcept
Return access to the temperature field.
Base class for solarLoad models.
solarLoadBase(const fvMesh &mesh)
Construct.
The solarLoad radiation model includes Sun primary hits, their reflective fluxes and diffusive sky ra...
Definition solarLoad.H:163
TypeName("solarLoad")
Runtime type information.
virtual tmp< DimensionedField< scalar, volMesh > > Ru() const
Source term component (constant).
Definition solarLoad.C:922
virtual tmp< volScalarField > Rp() const
Source term component (for power of T^4).
Definition solarLoad.C:905
virtual const solarCalculator & solarCalculatorRef() const
Return const reference to the solar calculator.
Definition solarLoad.C:929
const volScalarField & qprimaryRad(const label bandI) const
Return const access to the primary solar heat flux.
Definition solarLoad.H:372
virtual const faceShading & faceShadingRef() const
Return const reference to the face shading calculator.
Definition solarLoad.C:936
label nBands() const noexcept
Return const access to the number of bands.
Definition solarLoad.H:364
virtual ~solarLoad()=default
Destructor.
bool read()
Read radiationProperties dictionary.
Definition solarLoad.C:813
void calculate()
Solve radiation equations.
Definition solarLoad.C:824
A solar calculator model providing models for the solar direction and solar loads.
A class for managing temporary objects.
Definition tmp.H:75
Namespace for radiation modelling.
Namespace for OpenFOAM.
List< label > labelList
A List of labels.
Definition List.H:62
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const direction noexcept
Definition scalarImpl.H:265
List< scalar > scalarList
List of scalar.
Definition scalarList.H:32
dictionary dict
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition typeInfo.H:68