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 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 
27 Class
28  Foam::radiation::solarLoad
29 
30 Group
31  grpRadiationModels
32 
33 Description
34  The solar load 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 (wallCoupled = false). On coupled patches the flux is
47  by default added to the wall and considered into the solid
48  (solidCoupled = true).
49 
50  The solarLoad model can be used in conjuntion with fvDOM and viewFactor
51  radiation models. The flag useSolarLoad must be true on the rediation
52  dictionary.
53 
54 
55 SourceFiles
56  solarLoad.C
57 
58 \*---------------------------------------------------------------------------*/
59 
60 #ifndef radiation_solarLoad_H
61 #define radiation_solarLoad_H
62 
63 #include "radiationModel.H"
64 #include "volFields.H"
65 #include "faceShading.H"
66 #include "faceReflecting.H"
67 #include "solarCalculator.H"
68 
69 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
70 
71 namespace Foam
72 {
73 namespace radiation
74 {
75 
76 /*---------------------------------------------------------------------------*\
77  Class solarLoad Declaration
78 \*---------------------------------------------------------------------------*/
79 
80 class solarLoad
81 :
82  public radiationModel
83 {
84 
85 private:
86 
87  // Private data
88 
89  //- Dictionary
90  dictionary dict_;
91 
92  //- Net radiative heat flux [W/m2]
93  volScalarField qr_;
94 
95  //- Direct hit faces Ids
96  autoPtr<faceShading> hitFaces_;
97 
98  //- Reflected faces
99  autoPtr<faceReflecting> reflectedFaces_;
100 
101  //- Source term for cells next to patches with flags solidCoupled
102  //- and wallCoupled false
104 
105  //- Solar calculator
106  solarCalculator solarCalc_;
107 
108  //- Vertical direction (Default is g vector)
109  vector verticalDir_;
110 
111  //- Include reflected rays from specular surfaces
112  bool useReflectedRays_;
113 
114  //- Spectral distribution for the integrated solar heat flux
115  scalarList spectralDistribution_;
116 
117  //-Number of bands
118  label nBands_;
119 
120  //- Primary solar radiative heat flux per band [W/m2]
121  PtrList<volScalarField> qprimaryRad_;
122 
123  //- Couple solids through mapped boundary patch using qr (default:true)
124  bool solidCoupled_;
125 
126  //- Couple wall patches using qr (default:false)
127  bool wallCoupled_;
128 
129  //- Absorptivity list
130  List<List<tmp<scalarField>>> absorptivity_;
131 
132  //- Update absorptivity
133  bool updateAbsorptivity_;
134 
135  //- First iteration
136  bool firstIter_;
137 
138  //- Update Sun position index
139  label updateTimeIndex_;
140 
141 
142  // Private Member Functions
143 
144  //- Initialise
145  void initialise(const dictionary&);
146 
147  //- Update direct hit faces radiation
148  void updateDirectHitRadiation(const labelList&, const labelHashSet&);
149 
150  //- Update reflected heat flux
151  void updateReflectedRays(const labelHashSet&);
152 
153  //- Calculate diffusive heat flux
154  //void calculateQdiff(const labelHashSet&, const labelHashSet&);
155 
156  //- Update Sky diffusive radiation
157  void updateSkyDiffusiveRadiation
158  (
159  const labelHashSet&,
160  const labelHashSet&
161  );
162 
163  //- Update hit faces
164  bool updateHitFaces();
165 
166  //- Update absorptivity
167  void updateAbsorptivity(const labelHashSet& includePatches);
168 
169  //- No copy construct
170  solarLoad(const solarLoad&) = delete;
171 
172  //- No copy assignment
173  void operator=(const solarLoad&) = delete;
174 
175 
176 public:
177 
178  //- Runtime type information
179  TypeName("solarLoad");
180 
181 
182  // Constructors
183 
184  //- Construct from volScalarField
185  solarLoad(const volScalarField& T);
186 
187  //- Construct from dictionary and volScalarField
188  solarLoad(const dictionary& dict, const volScalarField& T);
189 
190 
191  //- Destructor
192  virtual ~solarLoad() = default;
193 
194 
195  // Member functions
196 
197  // Edit
198 
199  //- Solve
200  void calculate();
201 
202  //- Read radiation properties dictionary
203  bool read();
204 
205  //- Source term component (for power of T^4)
206  virtual tmp<volScalarField> Rp() const;
207 
208  //- Source term component (constant)
209  virtual tmp<DimensionedField<scalar, volMesh>> Ru() const;
210 
211 
212  // Access
213 
214  //- Number of bands
215  label nBands() const;
216 
217  //- Primary solar heat flux
218  const volScalarField& qprimaryRad(const label bandI) const
219  {
220  return qprimaryRad_[bandI];
221  }
222 };
223 
224 
225 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
226 
227 } // End namespace radiation
228 } // End namespace Foam
229 
230 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
231 
232 #endif
233 
234 // ************************************************************************* //
volFields.H
Foam::radiation::solarLoad::nBands
label nBands() const
Number of bands.
Definition: solarLoad.C:824
faceReflecting.H
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::HashSet< label, Hash< label > >
Foam::radiation::solarLoad::Ru
virtual tmp< DimensionedField< scalar, volMesh > > Ru() const
Source term component (constant)
Definition: solarLoad.C:934
Foam::radiation::solarLoad::TypeName
TypeName("solarLoad")
Runtime type information.
Foam::solarCalculator
The solar calculator model provides information about the Sun direction and Sun load model....
Definition: solarCalculator.H:105
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:62
faceShading.H
radiation
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
Foam::radiation::solarLoad::qprimaryRad
const volScalarField & qprimaryRad(const label bandI) const
Primary solar heat flux.
Definition: solarLoad.H:217
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:121
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::radiation::solarLoad::calculate
void calculate()
Solve.
Definition: solarLoad.C:841
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::radiation::solarLoad
The solar load radiation model includes Sun primary hits, their reflective fluxes and diffusive sky r...
Definition: solarLoad.H:79
Foam::Vector< scalar >
Foam::List< scalar >
Foam::radiation::radiationModel
Top level model for radiation modelling.
Definition: radiationModel.H:75
Foam::radiation::solarLoad::~solarLoad
virtual ~solarLoad()=default
Destructor.
solarCalculator.H
Foam::radiation::solarLoad::read
bool read()
Read radiation properties dictionary.
Definition: solarLoad.C:830
Foam::GeometricField< scalar, fvPatchField, volMesh >
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:54
Foam::radiation::solarLoad::Rp
virtual tmp< volScalarField > Rp() const
Source term component (for power of T^4)
Definition: solarLoad.C:910