solarCalculator.C
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-2017 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "solarCalculator.H"
29 #include "Time.H"
30 #include "unitConversion.H"
31 #include "constants.H"
32 
33 using namespace Foam::constant;
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39  defineTypeNameAndDebug(solarCalculator, 0);
40 }
41 
42 
43 const Foam::Enum
44 <
46 >
48 ({
49  { sunDirModel::mSunDirConstant, "sunDirConstant" },
50  { sunDirModel::mSunDirTracking, "sunDirTracking" },
51 });
52 
53 
54 const Foam::Enum
55 <
57 >
59 ({
60  { sunLModel::mSunLoadConstant, "sunLoadConstant" },
61  {
62  sunLModel::mSunLoadFairWeatherConditions,
63  "sunLoadFairWeatherConditions"
64  },
65  { sunLModel::mSunLoadTheoreticalMaximum, "sunLoadTheoreticalMaximum" },
66 });
67 
68 
69 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
70 
71 void Foam::solarCalculator::calculateBetaTheta()
72 {
73  scalar runTime = 0.0;
74  switch (sunDirectionModel_)
75  {
76  case mSunDirTracking:
77  {
78  runTime = mesh_.time().value();
79  break;
80  }
81  case mSunDirConstant:
82  {
83  break;
84  }
85  }
86 
87  scalar LSM = 15.0*(dict_.get<scalar>("localStandardMeridian"));
88 
89  scalar D = dict_.get<scalar>("startDay") + runTime/86400.0;
90  scalar M = 6.24004 + 0.0172*D;
91  scalar EOT = -7.659*sin(M) + 9.863*sin(2*M + 3.5932);
92 
93  dict_.readEntry("startTime", startTime_);
94 
95  scalar LST = startTime_ + runTime/3600.0;
96 
97  scalar LON = dict_.get<scalar>("longitude");
98 
99  scalar AST = LST + EOT/60.0 + (LON - LSM)/15;
100 
101  scalar delta = 23.45*sin(degToRad((360*(284 + D))/365));
102 
103  scalar H = degToRad(15*(AST - 12));
104 
105  scalar L = degToRad(dict_.get<scalar>("latitude"));
106 
107  scalar deltaRad = degToRad(delta);
108  beta_ = max(asin(cos(L)*cos(deltaRad)*cos(H) + sin(L)*sin(deltaRad)), 1e-3);
109  theta_ = acos((sin(beta_)*sin(L) - sin(deltaRad))/(cos(beta_)*cos(L)));
110 
111  // theta is the angle between the SOUTH axis and the Sun
112  // If the hour angle is lower than zero (morning) the Sun is positioned
113  // on the East side.
114  if (H < 0)
115  {
116  theta_ += 2*(constant::mathematical::pi - theta_);
117  }
118 
119  if (debug)
120  {
121  Info << tab << "altitude : " << radToDeg(beta_) << endl;
122  Info << tab << "azimuth : " << radToDeg(theta_) << endl;
123  }
124 }
125 
126 
127 void Foam::solarCalculator::calculateSunDirection()
128 {
129  gridUp_ = normalised(dict_.get<vector>("gridUp"));
130  eastDir_ = normalised(dict_.get<vector>("gridEast"));
131 
132  coord_.reset
133  (
134  new coordinateSystem("grid", Zero, gridUp_, eastDir_)
135  );
136 
137  // Assuming 'z' vertical, 'y' North and 'x' East
138  direction_.z() = -sin(beta_);
139  direction_.y() = cos(beta_)*cos(theta_); // South axis
140  direction_.x() = cos(beta_)*sin(theta_); // West axis
141 
142  direction_.normalise();
143 
144  if (debug)
145  {
146  Info<< "Sun direction in absolute coordinates : " << direction_ <<endl;
147  }
148 
149  // Transform to actual coordinate system
150  direction_ = coord_->transform(direction_);
151 
152  if (debug)
153  {
154  Info<< "Sun direction in the Grid coordinates : " << direction_ <<endl;
155  }
156 }
157 
158 
159 void Foam::solarCalculator::init()
160 {
161  switch (sunDirectionModel_)
162  {
163  case mSunDirConstant:
164  {
165  if (dict_.readIfPresent("sunDirection", direction_))
166  {
167  direction_.normalise();
168  }
169  else
170  {
171  calculateBetaTheta();
172  calculateSunDirection();
173  }
174 
175  break;
176  }
177  case mSunDirTracking:
178  {
179  if (word(mesh_.ddtScheme("default")) == "steadyState")
180  {
182  << " Sun direction model can not be sunDirtracking if the "
183  << " case is steady " << nl << exit(FatalError);
184  }
185 
186  dict_.readEntry
187  (
188  "sunTrackingUpdateInterval",
189  sunTrackingUpdateInterval_
190  );
191 
192  calculateBetaTheta();
193  calculateSunDirection();
194  break;
195  }
196  }
197 
198  switch (sunLoadModel_)
199  {
200  case mSunLoadConstant:
201  {
202  dict_.readEntry("directSolarRad", directSolarRad_);
203  dict_.readEntry("diffuseSolarRad", diffuseSolarRad_);
204  break;
205  }
206  case mSunLoadFairWeatherConditions:
207  {
208  dict_.readIfPresent
209  (
210  "skyCloudCoverFraction",
211  skyCloudCoverFraction_
212  );
213 
214  dict_.readEntry("A", A_);
215  dict_.readEntry("B", B_);
216 
217  if (!dict_.readIfPresent("beta", beta_))
218  {
219  calculateBetaTheta();
220  }
221 
222  directSolarRad_ =
223  (1.0 - 0.75*pow(skyCloudCoverFraction_, 3.0))
224  * A_/exp(B_/sin(beta_));
225 
226  dict_.readEntry("groundReflectivity", groundReflectivity_);
227  break;
228  }
229  case mSunLoadTheoreticalMaximum:
230  {
231  dict_.readEntry("Setrn", Setrn_);
232  dict_.readEntry("SunPrime", SunPrime_);
233  directSolarRad_ = Setrn_*SunPrime_;
234 
235  dict_.readEntry("groundReflectivity", groundReflectivity_);
236  break;
237  }
238  }
239 }
240 
241 
242 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
243 
244 Foam::solarCalculator::solarCalculator
245 (
246  const dictionary& dict,
247  const fvMesh& mesh
248 )
249 :
250  mesh_(mesh),
251  dict_(dict),
252  direction_(Zero),
253  directSolarRad_(0.0),
254  diffuseSolarRad_(0.0),
255  groundReflectivity_(0.0),
256  A_(0.0),
257  B_(0.0),
258  beta_(0.0),
259  theta_(0.0),
260  skyCloudCoverFraction_(0.0),
261  Setrn_(0.0),
262  SunPrime_(0.0),
263  C_(dict.get<scalar>("C")),
264  sunDirectionModel_
265  (
266  sunDirectionModelTypeNames_.get("sunDirectionModel", dict)
267  ),
268  sunLoadModel_(sunLoadModelTypeNames_.get("sunLoadModel", dict)),
269  coord_()
270 {
271  init();
272 }
273 
274 
275 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
276 
278 {
279  switch (sunDirectionModel_)
280  {
281  case mSunDirConstant:
282  {
283  break;
284  }
285  case mSunDirTracking:
286  {
287  calculateBetaTheta();
288  calculateSunDirection();
289  directSolarRad_ = A_/exp(B_/sin(max(beta_, ROOTVSMALL)));
290  break;
291  }
292  }
293 }
294 
295 
296 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
runTime
engineTime & runTime
Definition: createEngineTime.H:13
L
const vector L(dict.get< vector >("L"))
Foam::Enum
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: IOstreamOption.H:51
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::constant
Different types of constants.
Definition: atomicConstants.C:38
Foam::sin
dimensionedScalar sin(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:264
H
volScalarField H(IOobject("H", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), mesh, dimensionedScalar(dimLength, Zero))
unitConversion.H
Unit conversion functions.
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:261
Foam::solarCalculator::sunLModel
sunLModel
Direct sun load models.
Definition: solarCalculator.H:119
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
delta
scalar delta
Definition: LISASMDCalcMethod2.H:8
Foam::radToDeg
constexpr scalar radToDeg(const scalar rad) noexcept
Conversion from radians to degrees.
Definition: unitConversion.H:54
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:84
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
Foam::degToRad
constexpr scalar degToRad(const scalar deg) noexcept
Conversion from degrees to radians.
Definition: unitConversion.H:48
constants.H
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::normalised
VectorSpace< Form, Cmpt, Ncmpts > normalised(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
Definition: VectorSpaceI.H:496
Time.H
Foam::solarCalculator::correctSunDirection
void correctSunDirection()
Recalculate.
Definition: solarCalculator.C:277
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::tab
constexpr char tab
Definition: Ostream.H:371
Foam::constant::mathematical::pi
constexpr scalar pi(M_PI)
Foam::nl
constexpr char nl
Definition: Ostream.H:372
Foam::solarCalculator::sunDirectionModelTypeNames_
static const Enum< sunDirModel > sunDirectionModelTypeNames_
Sun direction models.
Definition: solarCalculator.H:130
Foam::acos
dimensionedScalar acos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:268
M
#define M(I)
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
solarCalculator.H
D
const dimensionedScalar & D
Definition: solveBulkSurfactant.H:4
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::solarCalculator::sunDirModel
sunDirModel
Sun direction models.
Definition: solarCalculator.H:112
Foam::solarCalculator::sunLoadModelTypeNames_
static const Enum< sunLModel > sunLoadModelTypeNames_
Sun load models.
Definition: solarCalculator.H:133
Foam::asin
dimensionedScalar asin(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:267
Foam::cos
dimensionedScalar cos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:265