| 1 | 
greg | 
2.1 | 
#ifndef lint | 
| 2 | 
greg | 
2.18 | 
static const char RCSid[] = "$Id: gendaymtx.c,v 2.17 2014/10/26 17:35:53 greg Exp $"; | 
| 3 | 
greg | 
2.1 | 
#endif | 
| 4 | 
  | 
  | 
/* | 
| 5 | 
  | 
  | 
 *  gendaymtx.c | 
| 6 | 
  | 
  | 
 *   | 
| 7 | 
  | 
  | 
 *  Generate a daylight matrix based on Perez Sky Model. | 
| 8 | 
  | 
  | 
 * | 
| 9 | 
  | 
  | 
 *  Most of this code is borrowed (see copyright below) from Ian Ashdown's | 
| 10 | 
  | 
  | 
 *  excellent re-implementation of Jean-Jacques Delaunay's gendaylit.c | 
| 11 | 
  | 
  | 
 * | 
| 12 | 
  | 
  | 
 *  Created by Greg Ward on 1/16/13. | 
| 13 | 
  | 
  | 
 */ | 
| 14 | 
  | 
  | 
 | 
| 15 | 
  | 
  | 
/********************************************************************* | 
| 16 | 
  | 
  | 
 * | 
| 17 | 
  | 
  | 
 *  H32_gendaylit.CPP - Perez Sky Model Calculation | 
| 18 | 
  | 
  | 
 * | 
| 19 | 
  | 
  | 
 *  Version:    1.00A | 
| 20 | 
  | 
  | 
 * | 
| 21 | 
  | 
  | 
 *  History:    09/10/01 - Created. | 
| 22 | 
  | 
  | 
 *                              11/10/08 - Modified for Unix compilation. | 
| 23 | 
  | 
  | 
 *                              11/10/12 - Fixed conditional __max directive. | 
| 24 | 
  | 
  | 
 *                              1/11/13 - Tweaks and optimizations (G.Ward) | 
| 25 | 
  | 
  | 
 * | 
| 26 | 
  | 
  | 
 *  Compilers:  Microsoft Visual C/C++ Professional V10.0     | 
| 27 | 
  | 
  | 
 * | 
| 28 | 
  | 
  | 
 *  Author:     Ian Ashdown, P.Eng. | 
| 29 | 
  | 
  | 
 *              byHeart Consultants Limited | 
| 30 | 
  | 
  | 
 *              620 Ballantree Road | 
| 31 | 
  | 
  | 
 *              West Vancouver, B.C. | 
| 32 | 
  | 
  | 
 *              Canada V7S 1W3 | 
| 33 | 
  | 
  | 
 *              e-mail: [email protected] | 
| 34 | 
  | 
  | 
 * | 
| 35 | 
  | 
  | 
 *  References: Perez, R., P. Ineichen, R. Seals, J. Michalsky, and R. | 
| 36 | 
  | 
  | 
 *                              Stewart. 1990. ìModeling Daylight Availability and | 
| 37 | 
  | 
  | 
 *                              Irradiance Components from Direct and Global | 
| 38 | 
  | 
  | 
 *                              Irradiance,î Solar Energy 44(5):271-289. | 
| 39 | 
  | 
  | 
 * | 
| 40 | 
  | 
  | 
 *                              Perez, R., R. Seals, and J. Michalsky. 1993. | 
| 41 | 
  | 
  | 
 *                              ìAll-Weather Model for Sky Luminance Distribution - | 
| 42 | 
  | 
  | 
 *                              Preliminary Configuration and Validation,î Solar Energy | 
| 43 | 
  | 
  | 
 *                              50(3):235-245. | 
| 44 | 
  | 
  | 
 * | 
| 45 | 
  | 
  | 
 *                              Perez, R., R. Seals, and J. Michalsky. 1993. "ERRATUM to | 
| 46 | 
  | 
  | 
 *                              All-Weather Model for Sky Luminance Distribution - | 
| 47 | 
  | 
  | 
 *                              Preliminary Configuration and Validation,î Solar Energy | 
| 48 | 
  | 
  | 
 *                              51(5):423. | 
| 49 | 
  | 
  | 
 * | 
| 50 | 
  | 
  | 
 *  NOTE:               This program is a completely rewritten version of | 
| 51 | 
  | 
  | 
 *                              gendaylit.c written by Jean-Jacques Delaunay (1994). | 
| 52 | 
  | 
  | 
 * | 
| 53 | 
  | 
  | 
 *  Copyright 2009-2012 byHeart Consultants Limited. All rights | 
| 54 | 
  | 
  | 
 *  reserved. | 
| 55 | 
  | 
  | 
 * | 
| 56 | 
  | 
  | 
 *  Redistribution and use in source and binary forms, with or without | 
| 57 | 
  | 
  | 
 *  modification, are permitted for personal and commercial purposes | 
| 58 | 
  | 
  | 
 *  provided that redistribution of source code must retain the above | 
| 59 | 
  | 
  | 
 *  copyright notice, this list of conditions and the following | 
| 60 | 
  | 
  | 
 *  disclaimer: | 
| 61 | 
  | 
  | 
 * | 
| 62 | 
  | 
  | 
 *    THIS SOFTWARE IS PROVIDED "AS IS" AND ANY EXPRESSED OR IMPLIED | 
| 63 | 
  | 
  | 
 *    WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED WARRANTIES | 
| 64 | 
  | 
  | 
 *    OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE | 
| 65 | 
  | 
  | 
 *    DISCLAIMED. IN NO EVENT SHALL byHeart Consultants Limited OR | 
| 66 | 
  | 
  | 
 *    ITS EMPLOYEES BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, | 
| 67 | 
  | 
  | 
 *    SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT | 
| 68 | 
  | 
  | 
 *    LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF | 
| 69 | 
  | 
  | 
 *    USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED | 
| 70 | 
  | 
  | 
 *    AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT | 
| 71 | 
  | 
  | 
 *    LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN | 
| 72 | 
  | 
  | 
 *    ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE | 
| 73 | 
  | 
  | 
 *    POSSIBILITY OF SUCH DAMAGE. | 
| 74 | 
  | 
  | 
 * | 
| 75 | 
  | 
  | 
 *********************************************************************/ | 
| 76 | 
  | 
  | 
 | 
| 77 | 
  | 
  | 
/* Zenith is along the Z-axis */ | 
| 78 | 
  | 
  | 
/* X-axis points east */ | 
| 79 | 
  | 
  | 
/* Y-axis points north */ | 
| 80 | 
  | 
  | 
/* azimuth is measured as degrees or radians east of North */ | 
| 81 | 
  | 
  | 
 | 
| 82 | 
  | 
  | 
/* Include files */ | 
| 83 | 
  | 
  | 
#define _USE_MATH_DEFINES | 
| 84 | 
  | 
  | 
#include <stdio.h> | 
| 85 | 
  | 
  | 
#include <stdlib.h> | 
| 86 | 
  | 
  | 
#include <string.h> | 
| 87 | 
  | 
  | 
#include <ctype.h> | 
| 88 | 
  | 
  | 
#include "rtmath.h" | 
| 89 | 
greg | 
2.14 | 
#include "platform.h" | 
| 90 | 
greg | 
2.1 | 
#include "color.h" | 
| 91 | 
greg | 
2.17 | 
#include "resolu.h" | 
| 92 | 
greg | 
2.1 | 
 | 
| 93 | 
  | 
  | 
char *progname;                                                         /* Program name */ | 
| 94 | 
  | 
  | 
char errmsg[128];                                                       /* Error message buffer */ | 
| 95 | 
  | 
  | 
const double DC_SolarConstantE = 1367.0;        /* Solar constant W/m^2 */ | 
| 96 | 
  | 
  | 
const double DC_SolarConstantL = 127.5;         /* Solar constant klux */ | 
| 97 | 
  | 
  | 
 | 
| 98 | 
  | 
  | 
double altitude;                        /* Solar altitude (radians) */ | 
| 99 | 
  | 
  | 
double azimuth;                         /* Solar azimuth (radians) */ | 
| 100 | 
  | 
  | 
double apwc;                            /* Atmospheric precipitable water content */ | 
| 101 | 
  | 
  | 
double dew_point = 11.0;                /* Surface dew point temperature (deg. C) */ | 
| 102 | 
  | 
  | 
double diff_illum;                      /* Diffuse illuminance */ | 
| 103 | 
  | 
  | 
double diff_irrad;                      /* Diffuse irradiance */ | 
| 104 | 
  | 
  | 
double dir_illum;                       /* Direct illuminance */ | 
| 105 | 
  | 
  | 
double dir_irrad;                       /* Direct irradiance */ | 
| 106 | 
  | 
  | 
int julian_date;                        /* Julian date */ | 
| 107 | 
  | 
  | 
double perez_param[5];                  /* Perez sky model parameters */ | 
| 108 | 
  | 
  | 
double sky_brightness;                  /* Sky brightness */ | 
| 109 | 
  | 
  | 
double sky_clearness;                   /* Sky clearness */ | 
| 110 | 
  | 
  | 
double solar_rad;                       /* Solar radiance */ | 
| 111 | 
  | 
  | 
double sun_zenith;                      /* Sun zenith angle (radians) */ | 
| 112 | 
  | 
  | 
int     input = 0;                              /* Input type */ | 
| 113 | 
greg | 
2.11 | 
int     output = 0;                             /* Output type */ | 
| 114 | 
greg | 
2.1 | 
 | 
| 115 | 
  | 
  | 
extern double dmax( double, double ); | 
| 116 | 
  | 
  | 
extern double CalcAirMass(); | 
| 117 | 
  | 
  | 
extern double CalcDiffuseIllumRatio( int ); | 
| 118 | 
  | 
  | 
extern double CalcDiffuseIrradiance(); | 
| 119 | 
  | 
  | 
extern double CalcDirectIllumRatio( int ); | 
| 120 | 
  | 
  | 
extern double CalcDirectIrradiance(); | 
| 121 | 
  | 
  | 
extern double CalcEccentricity(); | 
| 122 | 
  | 
  | 
extern double CalcPrecipWater( double ); | 
| 123 | 
  | 
  | 
extern double CalcRelHorzIllum( float *parr ); | 
| 124 | 
  | 
  | 
extern double CalcRelLuminance( double, double ); | 
| 125 | 
  | 
  | 
extern double CalcSkyBrightness(); | 
| 126 | 
  | 
  | 
extern double CalcSkyClearness(); | 
| 127 | 
  | 
  | 
extern int CalcSkyParamFromIllum(); | 
| 128 | 
  | 
  | 
extern int GetCategoryIndex(); | 
| 129 | 
  | 
  | 
extern void CalcPerezParam( double, double, double, int ); | 
| 130 | 
  | 
  | 
extern void CalcSkyPatchLumin( float *parr ); | 
| 131 | 
  | 
  | 
extern void ComputeSky( float *parr ); | 
| 132 | 
  | 
  | 
 | 
| 133 | 
  | 
  | 
/* Degrees into radians */ | 
| 134 | 
  | 
  | 
#define DegToRad(deg)   ((deg)*(PI/180.)) | 
| 135 | 
  | 
  | 
 | 
| 136 | 
  | 
  | 
/* Radiuans into degrees */ | 
| 137 | 
  | 
  | 
#define RadToDeg(rad)   ((rad)*(180./PI)) | 
| 138 | 
  | 
  | 
 | 
| 139 | 
  | 
  | 
 | 
| 140 | 
  | 
  | 
/* Perez sky model coefficients */ | 
| 141 | 
  | 
  | 
 | 
| 142 | 
  | 
  | 
/* Reference:   Perez, R., R. Seals, and J. Michalsky, 1993. "All- */ | 
| 143 | 
  | 
  | 
/*                              Weather Model for Sky Luminance Distribution - */ | 
| 144 | 
  | 
  | 
/*                              Preliminary Configuration and Validation," Solar */ | 
| 145 | 
  | 
  | 
/*                              Energy 50(3):235-245, Table 1. */ | 
| 146 | 
  | 
  | 
 | 
| 147 | 
  | 
  | 
static const double PerezCoeff[8][20] = | 
| 148 | 
  | 
  | 
{ | 
| 149 | 
  | 
  | 
        /* Sky clearness (epsilon): 1.000 to 1.065 */ | 
| 150 | 
  | 
  | 
        {   1.3525,  -0.2576,  -0.2690,  -1.4366,   -0.7670, | 
| 151 | 
  | 
  | 
            0.0007,   1.2734,  -0.1233,   2.8000,    0.6004, | 
| 152 | 
  | 
  | 
            1.2375,   1.0000,   1.8734,   0.6297,    0.9738, | 
| 153 | 
  | 
  | 
            0.2809,   0.0356,  -0.1246,  -0.5718,    0.9938 }, | 
| 154 | 
  | 
  | 
    /* Sky clearness (epsilon): 1.065 to 1.230 */ | 
| 155 | 
  | 
  | 
        {  -1.2219,  -0.7730,   1.4148,   1.1016,   -0.2054, | 
| 156 | 
  | 
  | 
            0.0367,  -3.9128,   0.9156,   6.9750,    0.1774, | 
| 157 | 
  | 
  | 
                6.4477,  -0.1239,  -1.5798,  -0.5081,   -1.7812, | 
| 158 | 
  | 
  | 
                0.1080,   0.2624,   0.0672,  -0.2190,   -0.4285 }, | 
| 159 | 
  | 
  | 
    /* Sky clearness (epsilon): 1.230 to 1.500 */ | 
| 160 | 
  | 
  | 
        {  -1.1000,  -0.2515,   0.8952,   0.0156,    0.2782, | 
| 161 | 
  | 
  | 
           -0.1812, - 4.5000,   1.1766,  24.7219,  -13.0812, | 
| 162 | 
  | 
  | 
          -37.7000,  34.8438,  -5.0000,   1.5218,    3.9229, | 
| 163 | 
  | 
  | 
           -2.6204,  -0.0156,   0.1597,   0.4199,   -0.5562 }, | 
| 164 | 
  | 
  | 
    /* Sky clearness (epsilon): 1.500 to 1.950 */ | 
| 165 | 
  | 
  | 
        {  -0.5484,  -0.6654,  -0.2672,   0.7117,   0.7234, | 
| 166 | 
  | 
  | 
           -0.6219,  -5.6812,   2.6297,  33.3389, -18.3000, | 
| 167 | 
  | 
  | 
          -62.2500,  52.0781,  -3.5000,   0.0016,   1.1477, | 
| 168 | 
  | 
  | 
            0.1062,   0.4659,  -0.3296,  -0.0876,  -0.0329 }, | 
| 169 | 
  | 
  | 
    /* Sky clearness (epsilon): 1.950 to 2.800 */ | 
| 170 | 
  | 
  | 
        {  -0.6000,  -0.3566,  -2.5000,   2.3250,   0.2937, | 
| 171 | 
  | 
  | 
            0.0496,  -5.6812,   1.8415,  21.0000,  -4.7656 , | 
| 172 | 
  | 
  | 
          -21.5906,   7.2492,  -3.5000,  -0.1554,   1.4062, | 
| 173 | 
  | 
  | 
            0.3988,   0.0032,   0.0766,  -0.0656,  -0.1294 }, | 
| 174 | 
  | 
  | 
    /* Sky clearness (epsilon): 2.800 to 4.500 */ | 
| 175 | 
  | 
  | 
        {  -1.0156,  -0.3670,   1.0078,   1.4051,   0.2875, | 
| 176 | 
  | 
  | 
           -0.5328,  -3.8500,   3.3750,  14.0000,  -0.9999, | 
| 177 | 
  | 
  | 
           -7.1406,   7.5469,  -3.4000,  -0.1078,  -1.0750, | 
| 178 | 
  | 
  | 
            1.5702,  -0.0672,   0.4016,   0.3017,  -0.4844 }, | 
| 179 | 
  | 
  | 
    /* Sky clearness (epsilon): 4.500 to 6.200 */ | 
| 180 | 
  | 
  | 
        {  -1.0000,   0.0211,   0.5025,  -0.5119,  -0.3000, | 
| 181 | 
  | 
  | 
            0.1922,   0.7023,  -1.6317,  19.0000,  -5.0000, | 
| 182 | 
  | 
  | 
                1.2438,  -1.9094,  -4.0000,   0.0250,   0.3844, | 
| 183 | 
  | 
  | 
                0.2656,   1.0468,  -0.3788,  -2.4517,   1.4656 }, | 
| 184 | 
  | 
  | 
    /* Sky clearness (epsilon): 6.200 to ... */ | 
| 185 | 
  | 
  | 
        {  -1.0500,   0.0289,   0.4260,   0.3590,  -0.3250, | 
| 186 | 
  | 
  | 
            0.1156,   0.7781,   0.0025,  31.0625, -14.5000, | 
| 187 | 
  | 
  | 
          -46.1148,  55.3750,  -7.2312,   0.4050,  13.3500, | 
| 188 | 
  | 
  | 
            0.6234,   1.5000,  -0.6426,   1.8564,   0.5636 } | 
| 189 | 
  | 
  | 
}; | 
| 190 | 
  | 
  | 
 | 
| 191 | 
  | 
  | 
/* Perez irradiance component model coefficients */ | 
| 192 | 
  | 
  | 
 | 
| 193 | 
  | 
  | 
/* Reference:   Perez, R., P. Ineichen, R. Seals, J. Michalsky, and R. */ | 
| 194 | 
  | 
  | 
/*                              Stewart. 1990. ìModeling Daylight Availability and */ | 
| 195 | 
  | 
  | 
/*                              Irradiance Components from Direct and Global */ | 
| 196 | 
  | 
  | 
/*                              Irradiance,î Solar Energy 44(5):271-289. */ | 
| 197 | 
  | 
  | 
 | 
| 198 | 
  | 
  | 
typedef struct | 
| 199 | 
  | 
  | 
{ | 
| 200 | 
  | 
  | 
        double lower;   /* Lower bound */ | 
| 201 | 
  | 
  | 
        double upper;   /* Upper bound */ | 
| 202 | 
  | 
  | 
} CategoryBounds; | 
| 203 | 
  | 
  | 
 | 
| 204 | 
  | 
  | 
/* Perez sky clearness (epsilon) categories (Table 1) */ | 
| 205 | 
  | 
  | 
static const CategoryBounds SkyClearCat[8] = | 
| 206 | 
  | 
  | 
{ | 
| 207 | 
  | 
  | 
        { 1.000, 1.065 },       /* Overcast */ | 
| 208 | 
  | 
  | 
        { 1.065, 1.230 }, | 
| 209 | 
  | 
  | 
        { 1.230, 1.500 }, | 
| 210 | 
  | 
  | 
        { 1.500, 1.950 }, | 
| 211 | 
  | 
  | 
        { 1.950, 2.800 }, | 
| 212 | 
  | 
  | 
        { 2.800, 4.500 }, | 
| 213 | 
  | 
  | 
        { 4.500, 6.200 }, | 
| 214 | 
greg | 
2.12 | 
        { 6.200, 12.01 }        /* Clear */ | 
| 215 | 
greg | 
2.1 | 
}; | 
| 216 | 
  | 
  | 
 | 
| 217 | 
  | 
  | 
/* Luminous efficacy model coefficients */ | 
| 218 | 
  | 
  | 
typedef struct | 
| 219 | 
  | 
  | 
{ | 
| 220 | 
  | 
  | 
        double a; | 
| 221 | 
  | 
  | 
        double b; | 
| 222 | 
  | 
  | 
        double c; | 
| 223 | 
  | 
  | 
        double d; | 
| 224 | 
  | 
  | 
} ModelCoeff; | 
| 225 | 
  | 
  | 
 | 
| 226 | 
  | 
  | 
/* Diffuse luminous efficacy model coefficients (Table 4, Eqn. 7) */ | 
| 227 | 
  | 
  | 
static const ModelCoeff DiffuseLumEff[8] = | 
| 228 | 
  | 
  | 
{ | 
| 229 | 
  | 
  | 
        {  97.24, -0.46,  12.00,  -8.91 }, | 
| 230 | 
  | 
  | 
        { 107.22,  1.15,   0.59,  -3.95 }, | 
| 231 | 
  | 
  | 
        { 104.97,  2.96,  -5.53,  -8.77 }, | 
| 232 | 
  | 
  | 
        { 102.39,  5.59, -13.95, -13.90 }, | 
| 233 | 
  | 
  | 
        { 100.71,  5.94, -22.75, -23.74 }, | 
| 234 | 
  | 
  | 
        { 106.42,  3.83, -36.15, -28.83 }, | 
| 235 | 
  | 
  | 
        { 141.88,  1.90, -53.24, -14.03 }, | 
| 236 | 
  | 
  | 
        { 152.23,  0.35, -45.27,  -7.98 } | 
| 237 | 
  | 
  | 
}; | 
| 238 | 
  | 
  | 
 | 
| 239 | 
  | 
  | 
/* Direct luminous efficacy model coefficients (Table 4, Eqn. 8) */ | 
| 240 | 
  | 
  | 
static const ModelCoeff DirectLumEff[8] = | 
| 241 | 
  | 
  | 
{ | 
| 242 | 
  | 
  | 
        {  57.20, -4.55, -2.98, 117.12 }, | 
| 243 | 
  | 
  | 
        {  98.99, -3.46, -1.21,  12.38 }, | 
| 244 | 
  | 
  | 
        { 109.83, -4.90, -1.71,  -8.81 }, | 
| 245 | 
  | 
  | 
        { 110.34, -5.84, -1.99,  -4.56 }, | 
| 246 | 
  | 
  | 
        { 106.36, -3.97, -1.75,  -6.16 }, | 
| 247 | 
  | 
  | 
        { 107.19, -1.25, -1.51, -26.73 }, | 
| 248 | 
  | 
  | 
        { 105.75,  0.77, -1.26, -34.44 }, | 
| 249 | 
  | 
  | 
        { 101.18,  1.58, -1.10,  -8.29 } | 
| 250 | 
  | 
  | 
}; | 
| 251 | 
  | 
  | 
 | 
| 252 | 
greg | 
2.3 | 
#ifndef NSUNPATCH | 
| 253 | 
greg | 
2.10 | 
#define NSUNPATCH       4               /* max. # patches to spread sun into */ | 
| 254 | 
greg | 
2.3 | 
#endif | 
| 255 | 
  | 
  | 
 | 
| 256 | 
greg | 
2.1 | 
extern int jdate(int month, int day); | 
| 257 | 
  | 
  | 
extern double stadj(int  jd); | 
| 258 | 
  | 
  | 
extern double sdec(int  jd); | 
| 259 | 
  | 
  | 
extern double salt(double sd, double st); | 
| 260 | 
  | 
  | 
extern double sazi(double sd, double st); | 
| 261 | 
  | 
  | 
                                        /* sun calculation constants */ | 
| 262 | 
  | 
  | 
extern double  s_latitude; | 
| 263 | 
  | 
  | 
extern double  s_longitude; | 
| 264 | 
  | 
  | 
extern double  s_meridian; | 
| 265 | 
  | 
  | 
 | 
| 266 | 
greg | 
2.10 | 
int             nsuns = NSUNPATCH;      /* number of sun patches to use */ | 
| 267 | 
  | 
  | 
double          fixed_sun_sa = -1;      /* fixed solid angle per sun? */ | 
| 268 | 
  | 
  | 
 | 
| 269 | 
greg | 
2.1 | 
int             verbose = 0;            /* progress reports to stderr? */ | 
| 270 | 
  | 
  | 
 | 
| 271 | 
  | 
  | 
int             outfmt = 'a';           /* output format */ | 
| 272 | 
  | 
  | 
 | 
| 273 | 
  | 
  | 
int             rhsubdiv = 1;           /* Reinhart sky subdivisions */ | 
| 274 | 
  | 
  | 
 | 
| 275 | 
greg | 
2.4 | 
COLOR           skycolor = {.96, 1.004, 1.118}; /* sky coloration */ | 
| 276 | 
  | 
  | 
COLOR           suncolor = {1., 1., 1.};        /* sun color */ | 
| 277 | 
  | 
  | 
COLOR           grefl = {.2, .2, .2};           /* ground reflectance */ | 
| 278 | 
greg | 
2.1 | 
 | 
| 279 | 
  | 
  | 
int             nskypatch;              /* number of Reinhart patches */ | 
| 280 | 
  | 
  | 
float           *rh_palt;               /* sky patch altitudes (radians) */ | 
| 281 | 
  | 
  | 
float           *rh_pazi;               /* sky patch azimuths (radians) */ | 
| 282 | 
  | 
  | 
float           *rh_dom;                /* sky patch solid angle (sr) */ | 
| 283 | 
  | 
  | 
 | 
| 284 | 
  | 
  | 
#define         vector(v,alt,azi)       (       (v)[1] = tcos(alt), \ | 
| 285 | 
  | 
  | 
                                                (v)[0] = (v)[1]*tsin(azi), \ | 
| 286 | 
  | 
  | 
                                                (v)[1] *= tcos(azi), \ | 
| 287 | 
  | 
  | 
                                                (v)[2] = tsin(alt) ) | 
| 288 | 
  | 
  | 
 | 
| 289 | 
  | 
  | 
#define         rh_vector(v,i)          vector(v,rh_palt[i],rh_pazi[i]) | 
| 290 | 
  | 
  | 
 | 
| 291 | 
  | 
  | 
#define         rh_cos(i)               tsin(rh_palt[i]) | 
| 292 | 
  | 
  | 
 | 
| 293 | 
  | 
  | 
extern int      rh_init(void); | 
| 294 | 
  | 
  | 
extern float *  resize_dmatrix(float *mtx_data, int nsteps, int npatch); | 
| 295 | 
  | 
  | 
extern void     AddDirect(float *parr); | 
| 296 | 
  | 
  | 
 | 
| 297 | 
greg | 
2.14 | 
 | 
| 298 | 
  | 
  | 
static const char * | 
| 299 | 
  | 
  | 
getfmtname(int fmt) | 
| 300 | 
  | 
  | 
{ | 
| 301 | 
  | 
  | 
        switch (fmt) { | 
| 302 | 
  | 
  | 
        case 'a': | 
| 303 | 
  | 
  | 
                return("ascii"); | 
| 304 | 
  | 
  | 
        case 'f': | 
| 305 | 
  | 
  | 
                return("float"); | 
| 306 | 
  | 
  | 
        case 'd': | 
| 307 | 
  | 
  | 
                return("double"); | 
| 308 | 
  | 
  | 
        } | 
| 309 | 
  | 
  | 
        return("unknown"); | 
| 310 | 
  | 
  | 
} | 
| 311 | 
  | 
  | 
 | 
| 312 | 
  | 
  | 
 | 
| 313 | 
greg | 
2.1 | 
int | 
| 314 | 
  | 
  | 
main(int argc, char *argv[]) | 
| 315 | 
  | 
  | 
{ | 
| 316 | 
  | 
  | 
        char    buf[256]; | 
| 317 | 
greg | 
2.14 | 
        int     doheader = 1;           /* output header? */ | 
| 318 | 
greg | 
2.8 | 
        double  rotation = 0;           /* site rotation (degrees) */ | 
| 319 | 
greg | 
2.1 | 
        double  elevation;              /* site elevation (meters) */ | 
| 320 | 
  | 
  | 
        int     dir_is_horiz;           /* direct is meas. on horizontal? */ | 
| 321 | 
  | 
  | 
        float   *mtx_data = NULL;       /* our matrix data */ | 
| 322 | 
  | 
  | 
        int     ntsteps = 0;            /* number of rows in matrix */ | 
| 323 | 
greg | 
2.16 | 
        int     step_alloc = 0; | 
| 324 | 
greg | 
2.1 | 
        int     last_monthly = 0;       /* month of last report */ | 
| 325 | 
  | 
  | 
        int     mo, da;                 /* month (1-12) and day (1-31) */ | 
| 326 | 
  | 
  | 
        double  hr;                     /* hour (local standard time) */ | 
| 327 | 
  | 
  | 
        double  dir, dif;               /* direct and diffuse values */ | 
| 328 | 
  | 
  | 
        int     mtx_offset; | 
| 329 | 
  | 
  | 
        int     i, j; | 
| 330 | 
  | 
  | 
 | 
| 331 | 
  | 
  | 
        progname = argv[0]; | 
| 332 | 
  | 
  | 
                                        /* get options */ | 
| 333 | 
  | 
  | 
        for (i = 1; i < argc && argv[i][0] == '-'; i++) | 
| 334 | 
  | 
  | 
                switch (argv[i][1]) { | 
| 335 | 
greg | 
2.4 | 
                case 'g':                       /* ground reflectance */ | 
| 336 | 
  | 
  | 
                        grefl[0] = atof(argv[++i]); | 
| 337 | 
  | 
  | 
                        grefl[1] = atof(argv[++i]); | 
| 338 | 
  | 
  | 
                        grefl[2] = atof(argv[++i]); | 
| 339 | 
greg | 
2.1 | 
                        break; | 
| 340 | 
greg | 
2.4 | 
                case 'v':                       /* verbose progress reports */ | 
| 341 | 
greg | 
2.1 | 
                        verbose++; | 
| 342 | 
  | 
  | 
                        break; | 
| 343 | 
greg | 
2.14 | 
                case 'h':                       /* turn off header */ | 
| 344 | 
  | 
  | 
                        doheader = 0; | 
| 345 | 
  | 
  | 
                        break; | 
| 346 | 
greg | 
2.4 | 
                case 'o':                       /* output format */ | 
| 347 | 
greg | 
2.1 | 
                        switch (argv[i][2]) { | 
| 348 | 
  | 
  | 
                        case 'f': | 
| 349 | 
  | 
  | 
                        case 'd': | 
| 350 | 
  | 
  | 
                        case 'a': | 
| 351 | 
  | 
  | 
                                outfmt = argv[i][2]; | 
| 352 | 
  | 
  | 
                                break; | 
| 353 | 
  | 
  | 
                        default: | 
| 354 | 
  | 
  | 
                                goto userr; | 
| 355 | 
  | 
  | 
                        } | 
| 356 | 
  | 
  | 
                        break; | 
| 357 | 
greg | 
2.11 | 
                case 'O':                       /* output type */ | 
| 358 | 
  | 
  | 
                        switch (argv[i][2]) { | 
| 359 | 
  | 
  | 
                        case '0': | 
| 360 | 
  | 
  | 
                                output = 0; | 
| 361 | 
  | 
  | 
                                break; | 
| 362 | 
  | 
  | 
                        case '1': | 
| 363 | 
  | 
  | 
                                output = 1; | 
| 364 | 
  | 
  | 
                                break; | 
| 365 | 
  | 
  | 
                        default: | 
| 366 | 
  | 
  | 
                                goto userr; | 
| 367 | 
  | 
  | 
                        } | 
| 368 | 
  | 
  | 
                        if (argv[i][3]) | 
| 369 | 
  | 
  | 
                                goto userr; | 
| 370 | 
  | 
  | 
                        break; | 
| 371 | 
greg | 
2.4 | 
                case 'm':                       /* Reinhart subdivisions */ | 
| 372 | 
greg | 
2.1 | 
                        rhsubdiv = atoi(argv[++i]); | 
| 373 | 
  | 
  | 
                        break; | 
| 374 | 
greg | 
2.4 | 
                case 'c':                       /* sky color */ | 
| 375 | 
greg | 
2.1 | 
                        skycolor[0] = atof(argv[++i]); | 
| 376 | 
  | 
  | 
                        skycolor[1] = atof(argv[++i]); | 
| 377 | 
  | 
  | 
                        skycolor[2] = atof(argv[++i]); | 
| 378 | 
  | 
  | 
                        break; | 
| 379 | 
greg | 
2.4 | 
                case 'd':                       /* solar (direct) only */ | 
| 380 | 
greg | 
2.1 | 
                        skycolor[0] = skycolor[1] = skycolor[2] = 0; | 
| 381 | 
greg | 
2.4 | 
                        if (suncolor[1] <= 1e-4) | 
| 382 | 
  | 
  | 
                                suncolor[0] = suncolor[1] = suncolor[2] = 1; | 
| 383 | 
greg | 
2.1 | 
                        break; | 
| 384 | 
greg | 
2.4 | 
                case 's':                       /* sky only (no direct) */ | 
| 385 | 
  | 
  | 
                        suncolor[0] = suncolor[1] = suncolor[2] = 0; | 
| 386 | 
greg | 
2.1 | 
                        if (skycolor[1] <= 1e-4) | 
| 387 | 
  | 
  | 
                                skycolor[0] = skycolor[1] = skycolor[2] = 1; | 
| 388 | 
  | 
  | 
                        break; | 
| 389 | 
greg | 
2.8 | 
                case 'r':                       /* rotate distribution */ | 
| 390 | 
  | 
  | 
                        if (argv[i][2] && argv[i][2] != 'z') | 
| 391 | 
  | 
  | 
                                goto userr; | 
| 392 | 
  | 
  | 
                        rotation = atof(argv[++i]); | 
| 393 | 
  | 
  | 
                        break; | 
| 394 | 
greg | 
2.10 | 
                case '5':                       /* 5-phase calculation */ | 
| 395 | 
  | 
  | 
                        nsuns = 1; | 
| 396 | 
  | 
  | 
                        fixed_sun_sa = 6.797e-05; | 
| 397 | 
  | 
  | 
                        break; | 
| 398 | 
greg | 
2.1 | 
                default: | 
| 399 | 
  | 
  | 
                        goto userr; | 
| 400 | 
  | 
  | 
                } | 
| 401 | 
  | 
  | 
        if (i < argc-1) | 
| 402 | 
  | 
  | 
                goto userr; | 
| 403 | 
  | 
  | 
        if (i == argc-1 && freopen(argv[i], "r", stdin) == NULL) { | 
| 404 | 
  | 
  | 
                fprintf(stderr, "%s: cannot open '%s' for input\n", | 
| 405 | 
  | 
  | 
                                progname, argv[i]); | 
| 406 | 
  | 
  | 
                exit(1); | 
| 407 | 
  | 
  | 
        } | 
| 408 | 
  | 
  | 
        if (verbose) { | 
| 409 | 
  | 
  | 
                if (i == argc-1) | 
| 410 | 
  | 
  | 
                        fprintf(stderr, "%s: reading weather tape '%s'\n", | 
| 411 | 
  | 
  | 
                                        progname, argv[i]); | 
| 412 | 
  | 
  | 
                else | 
| 413 | 
  | 
  | 
                        fprintf(stderr, "%s: reading weather tape from <stdin>\n", | 
| 414 | 
  | 
  | 
                                        progname); | 
| 415 | 
  | 
  | 
        } | 
| 416 | 
  | 
  | 
                                        /* read weather tape header */ | 
| 417 | 
greg | 
2.2 | 
        if (scanf("place %[^\r\n] ", buf) != 1) | 
| 418 | 
greg | 
2.1 | 
                goto fmterr; | 
| 419 | 
  | 
  | 
        if (scanf("latitude %lf\n", &s_latitude) != 1) | 
| 420 | 
  | 
  | 
                goto fmterr; | 
| 421 | 
  | 
  | 
        if (scanf("longitude %lf\n", &s_longitude) != 1) | 
| 422 | 
  | 
  | 
                goto fmterr; | 
| 423 | 
  | 
  | 
        if (scanf("time_zone %lf\n", &s_meridian) != 1) | 
| 424 | 
  | 
  | 
                goto fmterr; | 
| 425 | 
  | 
  | 
        if (scanf("site_elevation %lf\n", &elevation) != 1) | 
| 426 | 
  | 
  | 
                goto fmterr; | 
| 427 | 
  | 
  | 
        if (scanf("weather_data_file_units %d\n", &input) != 1) | 
| 428 | 
  | 
  | 
                goto fmterr; | 
| 429 | 
  | 
  | 
        switch (input) {                /* translate units */ | 
| 430 | 
  | 
  | 
        case 1: | 
| 431 | 
  | 
  | 
                input = 1;              /* radiometric quantities */ | 
| 432 | 
  | 
  | 
                dir_is_horiz = 0;       /* direct is perpendicular meas. */ | 
| 433 | 
  | 
  | 
                break; | 
| 434 | 
  | 
  | 
        case 2: | 
| 435 | 
  | 
  | 
                input = 1;              /* radiometric quantities */ | 
| 436 | 
  | 
  | 
                dir_is_horiz = 1;       /* solar measured horizontally */ | 
| 437 | 
  | 
  | 
                break; | 
| 438 | 
  | 
  | 
        case 3: | 
| 439 | 
  | 
  | 
                input = 2;              /* photometric quantities */ | 
| 440 | 
  | 
  | 
                dir_is_horiz = 0;       /* direct is perpendicular meas. */ | 
| 441 | 
  | 
  | 
                break; | 
| 442 | 
  | 
  | 
        default: | 
| 443 | 
  | 
  | 
                goto fmterr; | 
| 444 | 
  | 
  | 
        } | 
| 445 | 
  | 
  | 
        rh_init();                      /* initialize sky patches */ | 
| 446 | 
  | 
  | 
        if (verbose) { | 
| 447 | 
  | 
  | 
                fprintf(stderr, "%s: location '%s'\n", progname, buf); | 
| 448 | 
  | 
  | 
                fprintf(stderr, "%s: (lat,long)=(%.1f,%.1f) degrees north, west\n", | 
| 449 | 
  | 
  | 
                                progname, s_latitude, s_longitude); | 
| 450 | 
  | 
  | 
                fprintf(stderr, "%s: %d sky patches per time step\n", | 
| 451 | 
  | 
  | 
                                progname, nskypatch); | 
| 452 | 
greg | 
2.8 | 
                if (rotation != 0) | 
| 453 | 
  | 
  | 
                        fprintf(stderr, "%s: rotating output %.0f degrees\n", | 
| 454 | 
  | 
  | 
                                        progname, rotation); | 
| 455 | 
greg | 
2.1 | 
        } | 
| 456 | 
greg | 
2.2 | 
                                        /* convert quantities to radians */ | 
| 457 | 
  | 
  | 
        s_latitude = DegToRad(s_latitude); | 
| 458 | 
  | 
  | 
        s_longitude = DegToRad(s_longitude); | 
| 459 | 
  | 
  | 
        s_meridian = DegToRad(s_meridian); | 
| 460 | 
greg | 
2.1 | 
                                        /* process each time step in tape */ | 
| 461 | 
  | 
  | 
        while (scanf("%d %d %lf %lf %lf\n", &mo, &da, &hr, &dir, &dif) == 5) { | 
| 462 | 
  | 
  | 
                double          sda, sta; | 
| 463 | 
  | 
  | 
                                        /* make space for next time step */ | 
| 464 | 
  | 
  | 
                mtx_offset = 3*nskypatch*ntsteps++; | 
| 465 | 
greg | 
2.15 | 
                if (ntsteps > step_alloc) { | 
| 466 | 
greg | 
2.16 | 
                        step_alloc += (step_alloc>>1) + ntsteps + 7; | 
| 467 | 
greg | 
2.15 | 
                        mtx_data = resize_dmatrix(mtx_data, step_alloc, nskypatch); | 
| 468 | 
  | 
  | 
                } | 
| 469 | 
greg | 
2.1 | 
                if (dif <= 1e-4) { | 
| 470 | 
  | 
  | 
                        memset(mtx_data+mtx_offset, 0, sizeof(float)*3*nskypatch); | 
| 471 | 
  | 
  | 
                        continue; | 
| 472 | 
  | 
  | 
                } | 
| 473 | 
  | 
  | 
                if (verbose && mo != last_monthly) | 
| 474 | 
  | 
  | 
                        fprintf(stderr, "%s: stepping through month %d...\n", | 
| 475 | 
  | 
  | 
                                                progname, last_monthly=mo); | 
| 476 | 
  | 
  | 
                                        /* compute solar position */ | 
| 477 | 
  | 
  | 
                julian_date = jdate(mo, da); | 
| 478 | 
  | 
  | 
                sda = sdec(julian_date); | 
| 479 | 
  | 
  | 
                sta = stadj(julian_date); | 
| 480 | 
  | 
  | 
                altitude = salt(sda, hr+sta); | 
| 481 | 
greg | 
2.8 | 
                azimuth = sazi(sda, hr+sta) + PI - DegToRad(rotation); | 
| 482 | 
greg | 
2.1 | 
                                        /* convert measured values */ | 
| 483 | 
  | 
  | 
                if (dir_is_horiz && altitude > 0.) | 
| 484 | 
  | 
  | 
                        dir /= sin(altitude); | 
| 485 | 
  | 
  | 
                if (input == 1) { | 
| 486 | 
  | 
  | 
                        dir_irrad = dir; | 
| 487 | 
  | 
  | 
                        diff_irrad = dif; | 
| 488 | 
  | 
  | 
                } else /* input == 2 */ { | 
| 489 | 
  | 
  | 
                        dir_illum = dir; | 
| 490 | 
  | 
  | 
                        diff_illum = dif; | 
| 491 | 
  | 
  | 
                } | 
| 492 | 
  | 
  | 
                                        /* compute sky patch values */ | 
| 493 | 
  | 
  | 
                ComputeSky(mtx_data+mtx_offset); | 
| 494 | 
greg | 
2.4 | 
                AddDirect(mtx_data+mtx_offset); | 
| 495 | 
greg | 
2.1 | 
        } | 
| 496 | 
  | 
  | 
                                        /* check for junk at end */ | 
| 497 | 
  | 
  | 
        while ((i = fgetc(stdin)) != EOF) | 
| 498 | 
  | 
  | 
                if (!isspace(i)) { | 
| 499 | 
  | 
  | 
                        fprintf(stderr, "%s: warning - unexpected data past EOT: ", | 
| 500 | 
  | 
  | 
                                        progname); | 
| 501 | 
  | 
  | 
                        buf[0] = i; buf[1] = '\0'; | 
| 502 | 
  | 
  | 
                        fgets(buf+1, sizeof(buf)-1, stdin); | 
| 503 | 
  | 
  | 
                        fputs(buf, stderr); fputc('\n', stderr); | 
| 504 | 
  | 
  | 
                        break; | 
| 505 | 
  | 
  | 
                } | 
| 506 | 
  | 
  | 
                                        /* write out matrix */ | 
| 507 | 
greg | 
2.14 | 
        if (outfmt != 'a') | 
| 508 | 
  | 
  | 
                SET_FILE_BINARY(stdout); | 
| 509 | 
greg | 
2.1 | 
#ifdef getc_unlocked | 
| 510 | 
  | 
  | 
        flockfile(stdout); | 
| 511 | 
  | 
  | 
#endif | 
| 512 | 
  | 
  | 
        if (verbose) | 
| 513 | 
  | 
  | 
                fprintf(stderr, "%s: writing %smatrix with %d time steps...\n", | 
| 514 | 
  | 
  | 
                                progname, outfmt=='a' ? "" : "binary ", ntsteps); | 
| 515 | 
greg | 
2.14 | 
        if (doheader) { | 
| 516 | 
  | 
  | 
                newheader("RADIANCE", stdout); | 
| 517 | 
  | 
  | 
                printargs(argc, argv, stdout); | 
| 518 | 
  | 
  | 
                printf("LATLONG= %.8f %.8f\n", RadToDeg(s_latitude), | 
| 519 | 
  | 
  | 
                                        -RadToDeg(s_longitude)); | 
| 520 | 
  | 
  | 
                printf("NROWS=%d\n", nskypatch); | 
| 521 | 
  | 
  | 
                printf("NCOLS=%d\n", ntsteps); | 
| 522 | 
  | 
  | 
                printf("NCOMP=3\n"); | 
| 523 | 
greg | 
2.18 | 
                fputformat((char *)getfmtname(outfmt), stdout); | 
| 524 | 
greg | 
2.14 | 
                putchar('\n'); | 
| 525 | 
  | 
  | 
        } | 
| 526 | 
greg | 
2.1 | 
                                        /* patches are rows (outer sort) */ | 
| 527 | 
  | 
  | 
        for (i = 0; i < nskypatch; i++) { | 
| 528 | 
  | 
  | 
                mtx_offset = 3*i; | 
| 529 | 
  | 
  | 
                switch (outfmt) { | 
| 530 | 
  | 
  | 
                case 'a': | 
| 531 | 
  | 
  | 
                        for (j = 0; j < ntsteps; j++) { | 
| 532 | 
greg | 
2.3 | 
                                printf("%.3g %.3g %.3g\n", mtx_data[mtx_offset], | 
| 533 | 
greg | 
2.1 | 
                                                mtx_data[mtx_offset+1], | 
| 534 | 
  | 
  | 
                                                mtx_data[mtx_offset+2]); | 
| 535 | 
  | 
  | 
                                mtx_offset += 3*nskypatch; | 
| 536 | 
  | 
  | 
                        } | 
| 537 | 
greg | 
2.2 | 
                        if (ntsteps > 1) | 
| 538 | 
  | 
  | 
                                fputc('\n', stdout); | 
| 539 | 
greg | 
2.1 | 
                        break; | 
| 540 | 
  | 
  | 
                case 'f': | 
| 541 | 
  | 
  | 
                        for (j = 0; j < ntsteps; j++) { | 
| 542 | 
  | 
  | 
                                fwrite(mtx_data+mtx_offset, sizeof(float), 3, | 
| 543 | 
  | 
  | 
                                                stdout); | 
| 544 | 
  | 
  | 
                                mtx_offset += 3*nskypatch; | 
| 545 | 
  | 
  | 
                        } | 
| 546 | 
  | 
  | 
                        break; | 
| 547 | 
  | 
  | 
                case 'd': | 
| 548 | 
  | 
  | 
                        for (j = 0; j < ntsteps; j++) { | 
| 549 | 
  | 
  | 
                                double  ment[3]; | 
| 550 | 
  | 
  | 
                                ment[0] = mtx_data[mtx_offset]; | 
| 551 | 
  | 
  | 
                                ment[1] = mtx_data[mtx_offset+1]; | 
| 552 | 
  | 
  | 
                                ment[2] = mtx_data[mtx_offset+2]; | 
| 553 | 
  | 
  | 
                                fwrite(ment, sizeof(double), 3, stdout); | 
| 554 | 
  | 
  | 
                                mtx_offset += 3*nskypatch; | 
| 555 | 
  | 
  | 
                        } | 
| 556 | 
  | 
  | 
                        break; | 
| 557 | 
  | 
  | 
                } | 
| 558 | 
  | 
  | 
                if (ferror(stdout)) | 
| 559 | 
  | 
  | 
                        goto writerr; | 
| 560 | 
  | 
  | 
        } | 
| 561 | 
  | 
  | 
        if (fflush(stdout) == EOF) | 
| 562 | 
  | 
  | 
                goto writerr; | 
| 563 | 
  | 
  | 
        if (verbose) | 
| 564 | 
  | 
  | 
                fprintf(stderr, "%s: done.\n", progname); | 
| 565 | 
  | 
  | 
        exit(0); | 
| 566 | 
  | 
  | 
userr: | 
| 567 | 
greg | 
2.14 | 
        fprintf(stderr, "Usage: %s [-v][-h][-d|-s][-r deg][-m N][-g r g b][-c r g b][-o{f|d}][-O{0|1}] [tape.wea]\n", | 
| 568 | 
greg | 
2.1 | 
                        progname); | 
| 569 | 
  | 
  | 
        exit(1); | 
| 570 | 
  | 
  | 
fmterr: | 
| 571 | 
  | 
  | 
        fprintf(stderr, "%s: input weather tape format error\n", progname); | 
| 572 | 
  | 
  | 
        exit(1); | 
| 573 | 
  | 
  | 
writerr: | 
| 574 | 
  | 
  | 
        fprintf(stderr, "%s: write error on output\n", progname); | 
| 575 | 
  | 
  | 
        exit(1); | 
| 576 | 
  | 
  | 
} | 
| 577 | 
  | 
  | 
 | 
| 578 | 
  | 
  | 
/* Return maximum of two doubles */ | 
| 579 | 
  | 
  | 
double dmax( double a, double b ) | 
| 580 | 
  | 
  | 
{ return (a > b) ? a : b; } | 
| 581 | 
  | 
  | 
 | 
| 582 | 
  | 
  | 
/* Compute sky patch radiance values (modified by GW) */ | 
| 583 | 
  | 
  | 
void | 
| 584 | 
  | 
  | 
ComputeSky(float *parr) | 
| 585 | 
  | 
  | 
{ | 
| 586 | 
  | 
  | 
        int index;                      /* Category index */ | 
| 587 | 
  | 
  | 
        double norm_diff_illum;         /* Normalized diffuse illuimnance */ | 
| 588 | 
  | 
  | 
        int i; | 
| 589 | 
  | 
  | 
         | 
| 590 | 
  | 
  | 
        /* Calculate atmospheric precipitable water content */ | 
| 591 | 
  | 
  | 
        apwc = CalcPrecipWater(dew_point); | 
| 592 | 
  | 
  | 
 | 
| 593 | 
greg | 
2.6 | 
        /* Calculate sun zenith angle (don't let it dip below horizon) */ | 
| 594 | 
  | 
  | 
        /* Also limit minimum angle to keep circumsolar off zenith */ | 
| 595 | 
  | 
  | 
        if (altitude <= 0.0) | 
| 596 | 
  | 
  | 
                sun_zenith = DegToRad(90.0); | 
| 597 | 
  | 
  | 
        else if (altitude >= DegToRad(87.0)) | 
| 598 | 
  | 
  | 
                sun_zenith = DegToRad(3.0); | 
| 599 | 
  | 
  | 
        else | 
| 600 | 
  | 
  | 
                sun_zenith = DegToRad(90.0) - altitude; | 
| 601 | 
greg | 
2.1 | 
 | 
| 602 | 
  | 
  | 
        /* Compute the inputs for the calculation of the sky distribution */ | 
| 603 | 
  | 
  | 
         | 
| 604 | 
  | 
  | 
        if (input == 0)                                 /* XXX never used */ | 
| 605 | 
  | 
  | 
        { | 
| 606 | 
  | 
  | 
                /* Calculate irradiance */ | 
| 607 | 
  | 
  | 
                diff_irrad = CalcDiffuseIrradiance(); | 
| 608 | 
  | 
  | 
                dir_irrad = CalcDirectIrradiance(); | 
| 609 | 
  | 
  | 
                 | 
| 610 | 
  | 
  | 
                /* Calculate illuminance */ | 
| 611 | 
  | 
  | 
                index = GetCategoryIndex(); | 
| 612 | 
  | 
  | 
                diff_illum = diff_irrad * CalcDiffuseIllumRatio(index); | 
| 613 | 
  | 
  | 
                dir_illum = dir_irrad * CalcDirectIllumRatio(index); | 
| 614 | 
  | 
  | 
        } | 
| 615 | 
  | 
  | 
        else if (input == 1) | 
| 616 | 
  | 
  | 
        { | 
| 617 | 
  | 
  | 
                sky_brightness = CalcSkyBrightness(); | 
| 618 | 
  | 
  | 
                sky_clearness =  CalcSkyClearness(); | 
| 619 | 
  | 
  | 
 | 
| 620 | 
greg | 
2.9 | 
                /* Limit sky clearness */ | 
| 621 | 
  | 
  | 
                if (sky_clearness > 11.9) | 
| 622 | 
  | 
  | 
                        sky_clearness = 11.9; | 
| 623 | 
  | 
  | 
 | 
| 624 | 
  | 
  | 
                /* Limit sky brightness */ | 
| 625 | 
  | 
  | 
                if (sky_brightness < 0.01) | 
| 626 | 
greg | 
2.11 | 
                        sky_brightness = 0.01; | 
| 627 | 
greg | 
2.9 | 
 | 
| 628 | 
greg | 
2.1 | 
                /* Calculate illuminance */ | 
| 629 | 
  | 
  | 
                index = GetCategoryIndex(); | 
| 630 | 
  | 
  | 
                diff_illum = diff_irrad * CalcDiffuseIllumRatio(index); | 
| 631 | 
  | 
  | 
                dir_illum = dir_irrad * CalcDirectIllumRatio(index); | 
| 632 | 
  | 
  | 
        } | 
| 633 | 
  | 
  | 
        else if (input == 2) | 
| 634 | 
  | 
  | 
        { | 
| 635 | 
  | 
  | 
                /* Calculate sky brightness and clearness from illuminance values */ | 
| 636 | 
  | 
  | 
                index = CalcSkyParamFromIllum(); | 
| 637 | 
  | 
  | 
        } | 
| 638 | 
  | 
  | 
 | 
| 639 | 
greg | 
2.11 | 
        if (output == 1) {                      /* hack for solar radiance */ | 
| 640 | 
  | 
  | 
                diff_illum = diff_irrad * WHTEFFICACY; | 
| 641 | 
  | 
  | 
                dir_illum = dir_irrad * WHTEFFICACY; | 
| 642 | 
  | 
  | 
        } | 
| 643 | 
  | 
  | 
 | 
| 644 | 
greg | 
2.2 | 
        if (bright(skycolor) <= 1e-4) {                 /* 0 sky component? */ | 
| 645 | 
  | 
  | 
                memset(parr, 0, sizeof(float)*3*nskypatch); | 
| 646 | 
  | 
  | 
                return; | 
| 647 | 
  | 
  | 
        } | 
| 648 | 
greg | 
2.1 | 
        /* Compute ground radiance (include solar contribution if any) */ | 
| 649 | 
greg | 
2.3 | 
        parr[0] = diff_illum; | 
| 650 | 
greg | 
2.1 | 
        if (altitude > 0) | 
| 651 | 
greg | 
2.3 | 
                parr[0] += dir_illum * sin(altitude); | 
| 652 | 
greg | 
2.4 | 
        parr[2] = parr[1] = parr[0] *= (1./PI/WHTEFFICACY); | 
| 653 | 
  | 
  | 
        multcolor(parr, grefl); | 
| 654 | 
greg | 
2.1 | 
 | 
| 655 | 
  | 
  | 
        /* Calculate Perez sky model parameters */ | 
| 656 | 
  | 
  | 
        CalcPerezParam(sun_zenith, sky_clearness, sky_brightness, index); | 
| 657 | 
  | 
  | 
 | 
| 658 | 
  | 
  | 
        /* Calculate sky patch luminance values */ | 
| 659 | 
  | 
  | 
        CalcSkyPatchLumin(parr); | 
| 660 | 
  | 
  | 
 | 
| 661 | 
  | 
  | 
        /* Calculate relative horizontal illuminance */ | 
| 662 | 
  | 
  | 
        norm_diff_illum = CalcRelHorzIllum(parr); | 
| 663 | 
  | 
  | 
 | 
| 664 | 
greg | 
2.13 | 
        /* Check for zero sky -- make uniform in that case */ | 
| 665 | 
  | 
  | 
        if (norm_diff_illum <= FTINY) { | 
| 666 | 
  | 
  | 
                for (i = 1; i < nskypatch; i++) | 
| 667 | 
  | 
  | 
                        setcolor(parr+3*i, 1., 1., 1.); | 
| 668 | 
  | 
  | 
                norm_diff_illum = PI; | 
| 669 | 
  | 
  | 
        } | 
| 670 | 
greg | 
2.1 | 
        /* Normalization coefficient */ | 
| 671 | 
  | 
  | 
        norm_diff_illum = diff_illum / norm_diff_illum; | 
| 672 | 
  | 
  | 
 | 
| 673 | 
  | 
  | 
        /* Apply to sky patches to get absolute radiance values */ | 
| 674 | 
  | 
  | 
        for (i = 1; i < nskypatch; i++) { | 
| 675 | 
greg | 
2.7 | 
                scalecolor(parr+3*i, norm_diff_illum*(1./WHTEFFICACY)); | 
| 676 | 
greg | 
2.1 | 
                multcolor(parr+3*i, skycolor); | 
| 677 | 
  | 
  | 
        } | 
| 678 | 
  | 
  | 
} | 
| 679 | 
  | 
  | 
 | 
| 680 | 
  | 
  | 
/* Add in solar direct to nearest sky patches (GW) */ | 
| 681 | 
  | 
  | 
void | 
| 682 | 
  | 
  | 
AddDirect(float *parr) | 
| 683 | 
  | 
  | 
{ | 
| 684 | 
  | 
  | 
        FVECT   svec; | 
| 685 | 
greg | 
2.3 | 
        double  near_dprod[NSUNPATCH]; | 
| 686 | 
  | 
  | 
        int     near_patch[NSUNPATCH]; | 
| 687 | 
  | 
  | 
        double  wta[NSUNPATCH], wtot; | 
| 688 | 
greg | 
2.1 | 
        int     i, j, p; | 
| 689 | 
  | 
  | 
 | 
| 690 | 
greg | 
2.4 | 
        if (dir_illum <= 1e-4 || bright(suncolor) <= 1e-4) | 
| 691 | 
greg | 
2.1 | 
                return; | 
| 692 | 
greg | 
2.10 | 
                                        /* identify nsuns closest patches */ | 
| 693 | 
  | 
  | 
        if (nsuns > NSUNPATCH) | 
| 694 | 
  | 
  | 
                nsuns = NSUNPATCH; | 
| 695 | 
  | 
  | 
        else if (nsuns <= 0) | 
| 696 | 
  | 
  | 
                nsuns = 1; | 
| 697 | 
  | 
  | 
        for (i = nsuns; i--; ) | 
| 698 | 
greg | 
2.1 | 
                near_dprod[i] = -1.; | 
| 699 | 
  | 
  | 
        vector(svec, altitude, azimuth); | 
| 700 | 
  | 
  | 
        for (p = 1; p < nskypatch; p++) { | 
| 701 | 
  | 
  | 
                FVECT   pvec; | 
| 702 | 
  | 
  | 
                double  dprod; | 
| 703 | 
  | 
  | 
                rh_vector(pvec, p); | 
| 704 | 
  | 
  | 
                dprod = DOT(pvec, svec); | 
| 705 | 
greg | 
2.10 | 
                for (i = 0; i < nsuns; i++) | 
| 706 | 
greg | 
2.1 | 
                        if (dprod > near_dprod[i]) { | 
| 707 | 
greg | 
2.10 | 
                                for (j = nsuns; --j > i; ) { | 
| 708 | 
greg | 
2.1 | 
                                        near_dprod[j] = near_dprod[j-1]; | 
| 709 | 
  | 
  | 
                                        near_patch[j] = near_patch[j-1]; | 
| 710 | 
  | 
  | 
                                } | 
| 711 | 
  | 
  | 
                                near_dprod[i] = dprod; | 
| 712 | 
  | 
  | 
                                near_patch[i] = p; | 
| 713 | 
  | 
  | 
                                break; | 
| 714 | 
  | 
  | 
                        } | 
| 715 | 
  | 
  | 
        } | 
| 716 | 
  | 
  | 
        wtot = 0;                       /* weight by proximity */ | 
| 717 | 
greg | 
2.10 | 
        for (i = nsuns; i--; ) | 
| 718 | 
greg | 
2.1 | 
                wtot += wta[i] = 1./(1.002 - near_dprod[i]); | 
| 719 | 
  | 
  | 
                                        /* add to nearest patch radiances */ | 
| 720 | 
greg | 
2.10 | 
        for (i = nsuns; i--; ) { | 
| 721 | 
greg | 
2.2 | 
                float   *pdest = parr + 3*near_patch[i]; | 
| 722 | 
greg | 
2.10 | 
                float   val_add = wta[i] * dir_illum / (WHTEFFICACY * wtot); | 
| 723 | 
  | 
  | 
 | 
| 724 | 
  | 
  | 
                val_add /= (fixed_sun_sa > 0)   ? fixed_sun_sa  | 
| 725 | 
  | 
  | 
                                                : rh_dom[near_patch[i]] ; | 
| 726 | 
greg | 
2.4 | 
                *pdest++ += val_add*suncolor[0]; | 
| 727 | 
  | 
  | 
                *pdest++ += val_add*suncolor[1]; | 
| 728 | 
  | 
  | 
                *pdest++ += val_add*suncolor[2]; | 
| 729 | 
greg | 
2.2 | 
        } | 
| 730 | 
greg | 
2.1 | 
} | 
| 731 | 
  | 
  | 
 | 
| 732 | 
  | 
  | 
/* Initialize Reinhart sky patch positions (GW) */ | 
| 733 | 
  | 
  | 
int | 
| 734 | 
  | 
  | 
rh_init(void) | 
| 735 | 
  | 
  | 
{ | 
| 736 | 
  | 
  | 
#define NROW    7 | 
| 737 | 
  | 
  | 
        static const int        tnaz[NROW] = {30, 30, 24, 24, 18, 12, 6}; | 
| 738 | 
  | 
  | 
        const double            alpha = (PI/2.)/(NROW*rhsubdiv + .5); | 
| 739 | 
  | 
  | 
        int                     p, i, j; | 
| 740 | 
  | 
  | 
                                        /* allocate patch angle arrays */ | 
| 741 | 
  | 
  | 
        nskypatch = 0; | 
| 742 | 
  | 
  | 
        for (p = 0; p < NROW; p++) | 
| 743 | 
  | 
  | 
                nskypatch += tnaz[p]; | 
| 744 | 
  | 
  | 
        nskypatch *= rhsubdiv*rhsubdiv; | 
| 745 | 
  | 
  | 
        nskypatch += 2; | 
| 746 | 
  | 
  | 
        rh_palt = (float *)malloc(sizeof(float)*nskypatch); | 
| 747 | 
  | 
  | 
        rh_pazi = (float *)malloc(sizeof(float)*nskypatch); | 
| 748 | 
  | 
  | 
        rh_dom = (float *)malloc(sizeof(float)*nskypatch); | 
| 749 | 
  | 
  | 
        if ((rh_palt == NULL) | (rh_pazi == NULL) | (rh_dom == NULL)) { | 
| 750 | 
  | 
  | 
                fprintf(stderr, "%s: out of memory in rh_init()\n", progname); | 
| 751 | 
  | 
  | 
                exit(1); | 
| 752 | 
  | 
  | 
        } | 
| 753 | 
  | 
  | 
        rh_palt[0] = -PI/2.;            /* ground & zenith patches */ | 
| 754 | 
  | 
  | 
        rh_pazi[0] = 0.; | 
| 755 | 
  | 
  | 
        rh_dom[0] = 2.*PI; | 
| 756 | 
  | 
  | 
        rh_palt[nskypatch-1] = PI/2.; | 
| 757 | 
  | 
  | 
        rh_pazi[nskypatch-1] = 0.; | 
| 758 | 
  | 
  | 
        rh_dom[nskypatch-1] = 2.*PI*(1. - cos(alpha*.5)); | 
| 759 | 
  | 
  | 
        p = 1;                          /* "normal" patches */ | 
| 760 | 
  | 
  | 
        for (i = 0; i < NROW*rhsubdiv; i++) { | 
| 761 | 
  | 
  | 
                const float     ralt = alpha*(i + .5); | 
| 762 | 
  | 
  | 
                const int       ninrow = tnaz[i/rhsubdiv]*rhsubdiv; | 
| 763 | 
greg | 
2.3 | 
                const float     dom = 2.*PI*(sin(alpha*(i+1)) - sin(alpha*i)) / | 
| 764 | 
  | 
  | 
                                                (double)ninrow; | 
| 765 | 
greg | 
2.1 | 
                for (j = 0; j < ninrow; j++) { | 
| 766 | 
  | 
  | 
                        rh_palt[p] = ralt; | 
| 767 | 
  | 
  | 
                        rh_pazi[p] = 2.*PI * j / (double)ninrow; | 
| 768 | 
  | 
  | 
                        rh_dom[p++] = dom; | 
| 769 | 
  | 
  | 
                } | 
| 770 | 
  | 
  | 
        } | 
| 771 | 
  | 
  | 
        return nskypatch; | 
| 772 | 
  | 
  | 
#undef NROW | 
| 773 | 
  | 
  | 
} | 
| 774 | 
  | 
  | 
 | 
| 775 | 
  | 
  | 
/* Resize daylight matrix (GW) */ | 
| 776 | 
  | 
  | 
float * | 
| 777 | 
  | 
  | 
resize_dmatrix(float *mtx_data, int nsteps, int npatch) | 
| 778 | 
  | 
  | 
{ | 
| 779 | 
  | 
  | 
        if (mtx_data == NULL) | 
| 780 | 
  | 
  | 
                mtx_data = (float *)malloc(sizeof(float)*3*nsteps*npatch); | 
| 781 | 
  | 
  | 
        else | 
| 782 | 
  | 
  | 
                mtx_data = (float *)realloc(mtx_data, | 
| 783 | 
  | 
  | 
                                        sizeof(float)*3*nsteps*npatch); | 
| 784 | 
  | 
  | 
        if (mtx_data == NULL) { | 
| 785 | 
  | 
  | 
                fprintf(stderr, "%s: out of memory in resize_dmatrix(%d,%d)\n", | 
| 786 | 
  | 
  | 
                                progname, nsteps, npatch); | 
| 787 | 
  | 
  | 
                exit(1); | 
| 788 | 
  | 
  | 
        } | 
| 789 | 
  | 
  | 
        return(mtx_data); | 
| 790 | 
  | 
  | 
} | 
| 791 | 
  | 
  | 
 | 
| 792 | 
  | 
  | 
/* Determine category index */ | 
| 793 | 
  | 
  | 
int GetCategoryIndex() | 
| 794 | 
  | 
  | 
{ | 
| 795 | 
  | 
  | 
        int index;      /* Loop index */ | 
| 796 | 
  | 
  | 
 | 
| 797 | 
  | 
  | 
        for (index = 0; index < 8; index++) | 
| 798 | 
  | 
  | 
                if ((sky_clearness >= SkyClearCat[index].lower) && | 
| 799 | 
  | 
  | 
                                (sky_clearness < SkyClearCat[index].upper)) | 
| 800 | 
  | 
  | 
                        break; | 
| 801 | 
  | 
  | 
 | 
| 802 | 
  | 
  | 
        return index; | 
| 803 | 
  | 
  | 
} | 
| 804 | 
  | 
  | 
 | 
| 805 | 
  | 
  | 
/* Calculate diffuse illuminance to diffuse irradiance ratio */ | 
| 806 | 
  | 
  | 
 | 
| 807 | 
  | 
  | 
/* Reference:   Perez, R., P. Ineichen, R. Seals, J. Michalsky, and R. */ | 
| 808 | 
  | 
  | 
/*                              Stewart. 1990. ìModeling Daylight Availability and */ | 
| 809 | 
  | 
  | 
/*                              Irradiance Components from Direct and Global */ | 
| 810 | 
  | 
  | 
/*                              Irradiance,î Solar Energy 44(5):271-289, Eqn. 7. */ | 
| 811 | 
  | 
  | 
 | 
| 812 | 
  | 
  | 
double CalcDiffuseIllumRatio( int index ) | 
| 813 | 
  | 
  | 
{ | 
| 814 | 
  | 
  | 
        ModelCoeff const *pnle; /* Category coefficient pointer */ | 
| 815 | 
  | 
  | 
         | 
| 816 | 
  | 
  | 
        /* Get category coefficient pointer */ | 
| 817 | 
  | 
  | 
        pnle = &(DiffuseLumEff[index]); | 
| 818 | 
  | 
  | 
 | 
| 819 | 
  | 
  | 
        return pnle->a + pnle->b * apwc + pnle->c * cos(sun_zenith) + | 
| 820 | 
  | 
  | 
                        pnle->d * log(sky_brightness); | 
| 821 | 
  | 
  | 
} | 
| 822 | 
  | 
  | 
 | 
| 823 | 
  | 
  | 
/* Calculate direct illuminance to direct irradiance ratio */ | 
| 824 | 
  | 
  | 
 | 
| 825 | 
  | 
  | 
/* Reference:   Perez, R., P. Ineichen, R. Seals, J. Michalsky, and R. */ | 
| 826 | 
  | 
  | 
/*                              Stewart. 1990. ìModeling Daylight Availability and */ | 
| 827 | 
  | 
  | 
/*                              Irradiance Components from Direct and Global */ | 
| 828 | 
  | 
  | 
/*                              Irradiance,î Solar Energy 44(5):271-289, Eqn. 8. */ | 
| 829 | 
  | 
  | 
 | 
| 830 | 
  | 
  | 
double CalcDirectIllumRatio( int index ) | 
| 831 | 
  | 
  | 
{ | 
| 832 | 
  | 
  | 
        ModelCoeff const *pnle; /* Category coefficient pointer */ | 
| 833 | 
  | 
  | 
 | 
| 834 | 
  | 
  | 
        /* Get category coefficient pointer */ | 
| 835 | 
  | 
  | 
        pnle = &(DirectLumEff[index]); | 
| 836 | 
  | 
  | 
 | 
| 837 | 
  | 
  | 
        /* Calculate direct illuminance from direct irradiance */ | 
| 838 | 
  | 
  | 
         | 
| 839 | 
  | 
  | 
        return dmax((pnle->a + pnle->b * apwc + pnle->c * exp(5.73 * | 
| 840 | 
  | 
  | 
                        sun_zenith - 5.0) + pnle->d * sky_brightness), | 
| 841 | 
  | 
  | 
                        0.0); | 
| 842 | 
  | 
  | 
} | 
| 843 | 
  | 
  | 
 | 
| 844 | 
  | 
  | 
/* Calculate sky brightness */ | 
| 845 | 
  | 
  | 
 | 
| 846 | 
  | 
  | 
/* Reference:   Perez, R., P. Ineichen, R. Seals, J. Michalsky, and R. */ | 
| 847 | 
  | 
  | 
/*                              Stewart. 1990. ìModeling Daylight Availability and */ | 
| 848 | 
  | 
  | 
/*                              Irradiance Components from Direct and Global */ | 
| 849 | 
  | 
  | 
/*                              Irradiance,î Solar Energy 44(5):271-289, Eqn. 2. */ | 
| 850 | 
  | 
  | 
 | 
| 851 | 
  | 
  | 
double CalcSkyBrightness() | 
| 852 | 
  | 
  | 
{ | 
| 853 | 
  | 
  | 
        return diff_irrad * CalcAirMass() / (DC_SolarConstantE * | 
| 854 | 
  | 
  | 
                        CalcEccentricity()); | 
| 855 | 
  | 
  | 
} | 
| 856 | 
  | 
  | 
 | 
| 857 | 
  | 
  | 
/* Calculate sky clearness */ | 
| 858 | 
  | 
  | 
 | 
| 859 | 
  | 
  | 
/* Reference:   Perez, R., P. Ineichen, R. Seals, J. Michalsky, and R. */ | 
| 860 | 
  | 
  | 
/*                              Stewart. 1990. ìModeling Daylight Availability and */ | 
| 861 | 
  | 
  | 
/*                              Irradiance Components from Direct and Global */ | 
| 862 | 
  | 
  | 
/*                              Irradiance,î Solar Energy 44(5):271-289, Eqn. 1. */ | 
| 863 | 
  | 
  | 
 | 
| 864 | 
  | 
  | 
double CalcSkyClearness() | 
| 865 | 
  | 
  | 
{ | 
| 866 | 
  | 
  | 
        double sz_cubed;        /* Sun zenith angle cubed */ | 
| 867 | 
  | 
  | 
 | 
| 868 | 
  | 
  | 
        /* Calculate sun zenith angle cubed */ | 
| 869 | 
greg | 
2.11 | 
        sz_cubed = sun_zenith*sun_zenith*sun_zenith; | 
| 870 | 
greg | 
2.1 | 
 | 
| 871 | 
  | 
  | 
        return ((diff_irrad + dir_irrad) / diff_irrad + 1.041 * | 
| 872 | 
  | 
  | 
                        sz_cubed) / (1.0 + 1.041 * sz_cubed); | 
| 873 | 
  | 
  | 
} | 
| 874 | 
  | 
  | 
 | 
| 875 | 
  | 
  | 
/* Calculate diffuse horizontal irradiance from Perez sky brightness */ | 
| 876 | 
  | 
  | 
 | 
| 877 | 
  | 
  | 
/* Reference:   Perez, R., P. Ineichen, R. Seals, J. Michalsky, and R. */ | 
| 878 | 
  | 
  | 
/*                              Stewart. 1990. ìModeling Daylight Availability and */ | 
| 879 | 
  | 
  | 
/*                              Irradiance Components from Direct and Global */ | 
| 880 | 
  | 
  | 
/*                              Irradiance,î Solar Energy 44(5):271-289, Eqn. 2 */ | 
| 881 | 
  | 
  | 
/*                              (inverse). */ | 
| 882 | 
  | 
  | 
 | 
| 883 | 
  | 
  | 
double CalcDiffuseIrradiance() | 
| 884 | 
  | 
  | 
{ | 
| 885 | 
  | 
  | 
        return sky_brightness * DC_SolarConstantE * CalcEccentricity() / | 
| 886 | 
  | 
  | 
                        CalcAirMass(); | 
| 887 | 
  | 
  | 
} | 
| 888 | 
  | 
  | 
 | 
| 889 | 
  | 
  | 
/* Calculate direct normal irradiance from Perez sky clearness */ | 
| 890 | 
  | 
  | 
 | 
| 891 | 
  | 
  | 
/* Reference:   Perez, R., P. Ineichen, R. Seals, J. Michalsky, and R. */ | 
| 892 | 
  | 
  | 
/*                              Stewart. 1990. ìModeling Daylight Availability and */ | 
| 893 | 
  | 
  | 
/*                              Irradiance Components from Direct and Global */ | 
| 894 | 
  | 
  | 
/*                              Irradiance,î Solar Energy 44(5):271-289, Eqn. 1 */ | 
| 895 | 
  | 
  | 
/*                              (inverse). */ | 
| 896 | 
  | 
  | 
 | 
| 897 | 
  | 
  | 
double CalcDirectIrradiance() | 
| 898 | 
  | 
  | 
{ | 
| 899 | 
  | 
  | 
        return CalcDiffuseIrradiance() * ((sky_clearness - 1.0) * (1 + 1.041 | 
| 900 | 
greg | 
2.11 | 
                        * sun_zenith*sun_zenith*sun_zenith)); | 
| 901 | 
greg | 
2.1 | 
} | 
| 902 | 
  | 
  | 
 | 
| 903 | 
  | 
  | 
/* Calculate sky brightness and clearness from illuminance values */ | 
| 904 | 
  | 
  | 
int CalcSkyParamFromIllum() | 
| 905 | 
  | 
  | 
{ | 
| 906 | 
  | 
  | 
        double test1 = 0.1; | 
| 907 | 
  | 
  | 
        double test2 = 0.1; | 
| 908 | 
  | 
  | 
        int     counter = 0; | 
| 909 | 
  | 
  | 
        int index = 0;                  /* Category index */ | 
| 910 | 
  | 
  | 
 | 
| 911 | 
  | 
  | 
        /* Convert illuminance to irradiance */ | 
| 912 | 
  | 
  | 
        diff_irrad = diff_illum * DC_SolarConstantE / | 
| 913 | 
  | 
  | 
                        (DC_SolarConstantL * 1000.0); | 
| 914 | 
  | 
  | 
        dir_irrad = dir_illum * DC_SolarConstantE / | 
| 915 | 
  | 
  | 
                        (DC_SolarConstantL * 1000.0); | 
| 916 | 
  | 
  | 
 | 
| 917 | 
  | 
  | 
        /* Calculate sky brightness and clearness */ | 
| 918 | 
  | 
  | 
        sky_brightness = CalcSkyBrightness(); | 
| 919 | 
  | 
  | 
        sky_clearness =  CalcSkyClearness();  | 
| 920 | 
  | 
  | 
 | 
| 921 | 
  | 
  | 
        /* Limit sky clearness */ | 
| 922 | 
  | 
  | 
        if (sky_clearness > 12.0) | 
| 923 | 
  | 
  | 
                sky_clearness = 12.0; | 
| 924 | 
  | 
  | 
 | 
| 925 | 
  | 
  | 
        /* Limit sky brightness */ | 
| 926 | 
greg | 
2.9 | 
        if (sky_brightness < 0.01) | 
| 927 | 
greg | 
2.1 | 
                        sky_brightness = 0.01;  | 
| 928 | 
  | 
  | 
 | 
| 929 | 
  | 
  | 
        while (((fabs(diff_irrad - test1) > 10.0) || | 
| 930 | 
  | 
  | 
                        (fabs(dir_irrad - test2) > 10.0)) && !(counter == 5)) | 
| 931 | 
  | 
  | 
        { | 
| 932 | 
  | 
  | 
                test1 = diff_irrad; | 
| 933 | 
  | 
  | 
                test2 = dir_irrad;       | 
| 934 | 
  | 
  | 
                counter++; | 
| 935 | 
  | 
  | 
         | 
| 936 | 
  | 
  | 
                /* Convert illuminance to irradiance */ | 
| 937 | 
  | 
  | 
                index = GetCategoryIndex(); | 
| 938 | 
  | 
  | 
                diff_irrad = diff_illum / CalcDiffuseIllumRatio(index); | 
| 939 | 
  | 
  | 
                dir_irrad = dir_illum / CalcDirectIllumRatio(index); | 
| 940 | 
  | 
  | 
         | 
| 941 | 
  | 
  | 
                /* Calculate sky brightness and clearness */ | 
| 942 | 
  | 
  | 
                sky_brightness = CalcSkyBrightness(); | 
| 943 | 
  | 
  | 
                sky_clearness =  CalcSkyClearness(); | 
| 944 | 
  | 
  | 
 | 
| 945 | 
  | 
  | 
                /* Limit sky clearness */ | 
| 946 | 
  | 
  | 
                if (sky_clearness > 12.0) | 
| 947 | 
  | 
  | 
                        sky_clearness = 12.0; | 
| 948 | 
  | 
  | 
         | 
| 949 | 
  | 
  | 
                /* Limit sky brightness */ | 
| 950 | 
greg | 
2.9 | 
                if (sky_brightness < 0.01) | 
| 951 | 
greg | 
2.1 | 
                        sky_brightness = 0.01;  | 
| 952 | 
  | 
  | 
        } | 
| 953 | 
  | 
  | 
 | 
| 954 | 
  | 
  | 
        return GetCategoryIndex(); | 
| 955 | 
  | 
  | 
}                | 
| 956 | 
  | 
  | 
 | 
| 957 | 
  | 
  | 
/* Calculate relative luminance */ | 
| 958 | 
  | 
  | 
 | 
| 959 | 
  | 
  | 
/* Reference:   Perez, R., R. Seals, and J. Michalsky. 1993. */ | 
| 960 | 
  | 
  | 
/*                              ìAll-Weather Model for Sky Luminance Distribution - */ | 
| 961 | 
  | 
  | 
/*                              Preliminary Configuration and Validation,î Solar Energy */ | 
| 962 | 
  | 
  | 
/*                              50(3):235-245, Eqn. 1. */ | 
| 963 | 
  | 
  | 
 | 
| 964 | 
  | 
  | 
double CalcRelLuminance( double gamma, double zeta ) | 
| 965 | 
  | 
  | 
{ | 
| 966 | 
  | 
  | 
        return (1.0 + perez_param[0] * exp(perez_param[1] / cos(zeta))) * | 
| 967 | 
  | 
  | 
                    (1.0 + perez_param[2] * exp(perez_param[3] * gamma) + | 
| 968 | 
  | 
  | 
                        perez_param[4] * cos(gamma) * cos(gamma)); | 
| 969 | 
  | 
  | 
} | 
| 970 | 
  | 
  | 
 | 
| 971 | 
  | 
  | 
/* Calculate Perez sky model parameters */ | 
| 972 | 
  | 
  | 
 | 
| 973 | 
  | 
  | 
/* Reference:   Perez, R., R. Seals, and J. Michalsky. 1993. */ | 
| 974 | 
  | 
  | 
/*                              ìAll-Weather Model for Sky Luminance Distribution - */ | 
| 975 | 
  | 
  | 
/*                              Preliminary Configuration and Validation,î Solar Energy */ | 
| 976 | 
  | 
  | 
/*                              50(3):235-245, Eqns. 6 - 8. */ | 
| 977 | 
  | 
  | 
 | 
| 978 | 
  | 
  | 
void CalcPerezParam( double sz, double epsilon, double delta, | 
| 979 | 
  | 
  | 
                int index ) | 
| 980 | 
  | 
  | 
{ | 
| 981 | 
  | 
  | 
        double x[5][4];         /* Coefficents a, b, c, d, e */ | 
| 982 | 
  | 
  | 
        int i, j;                       /* Loop indices */ | 
| 983 | 
  | 
  | 
 | 
| 984 | 
  | 
  | 
        /* Limit sky brightness */ | 
| 985 | 
  | 
  | 
        if (epsilon > 1.065 && epsilon < 2.8) | 
| 986 | 
  | 
  | 
        { | 
| 987 | 
  | 
  | 
                if (delta < 0.2) | 
| 988 | 
  | 
  | 
                        delta = 0.2; | 
| 989 | 
  | 
  | 
        } | 
| 990 | 
  | 
  | 
 | 
| 991 | 
  | 
  | 
        /* Get Perez coefficients */ | 
| 992 | 
  | 
  | 
        for (i = 0; i < 5; i++) | 
| 993 | 
  | 
  | 
                for (j = 0; j < 4; j++) | 
| 994 | 
  | 
  | 
                        x[i][j] = PerezCoeff[index][4 * i + j]; | 
| 995 | 
  | 
  | 
 | 
| 996 | 
  | 
  | 
        if (index != 0) | 
| 997 | 
  | 
  | 
        { | 
| 998 | 
  | 
  | 
                /* Calculate parameter a, b, c, d and e (Eqn. 6) */ | 
| 999 | 
  | 
  | 
                for (i = 0; i < 5; i++) | 
| 1000 | 
  | 
  | 
                        perez_param[i] = x[i][0] + x[i][1] * sz + delta * (x[i][2] + | 
| 1001 | 
  | 
  | 
                                        x[i][3] * sz); | 
| 1002 | 
  | 
  | 
        } | 
| 1003 | 
  | 
  | 
        else | 
| 1004 | 
  | 
  | 
        { | 
| 1005 | 
  | 
  | 
                /* Parameters a, b and e (Eqn. 6) */ | 
| 1006 | 
  | 
  | 
                perez_param[0] = x[0][0] + x[0][1] * sz + delta * (x[0][2] + | 
| 1007 | 
  | 
  | 
                                x[0][3] * sz); | 
| 1008 | 
  | 
  | 
                perez_param[1] = x[1][0] + x[1][1] * sz + delta * (x[1][2] + | 
| 1009 | 
  | 
  | 
                                x[1][3] * sz); | 
| 1010 | 
  | 
  | 
                perez_param[4] = x[4][0] + x[4][1] * sz + delta * (x[4][2] + | 
| 1011 | 
  | 
  | 
                                x[4][3] * sz); | 
| 1012 | 
  | 
  | 
 | 
| 1013 | 
  | 
  | 
                /* Parameter c (Eqn. 7) */ | 
| 1014 | 
  | 
  | 
                perez_param[2] = exp(pow(delta * (x[2][0] + x[2][1] * sz), | 
| 1015 | 
  | 
  | 
                                x[2][2])) - x[2][3]; | 
| 1016 | 
  | 
  | 
 | 
| 1017 | 
  | 
  | 
                /* Parameter d (Eqn. 8) */ | 
| 1018 | 
  | 
  | 
                perez_param[3] = -exp(delta * (x[3][0] + x[3][1] * sz)) +  | 
| 1019 | 
  | 
  | 
                                x[3][2] + delta * x[3][3]; | 
| 1020 | 
  | 
  | 
        } | 
| 1021 | 
  | 
  | 
} | 
| 1022 | 
  | 
  | 
 | 
| 1023 | 
  | 
  | 
/* Calculate relative horizontal illuminance (modified by GW) */ | 
| 1024 | 
  | 
  | 
 | 
| 1025 | 
  | 
  | 
/* Reference:   Perez, R., R. Seals, and J. Michalsky. 1993. */ | 
| 1026 | 
  | 
  | 
/*                              ìAll-Weather Model for Sky Luminance Distribution - */ | 
| 1027 | 
  | 
  | 
/*                              Preliminary Configuration and Validation,î Solar Energy */ | 
| 1028 | 
  | 
  | 
/*                              50(3):235-245, Eqn. 3. */ | 
| 1029 | 
  | 
  | 
 | 
| 1030 | 
  | 
  | 
double CalcRelHorzIllum( float *parr ) | 
| 1031 | 
  | 
  | 
{ | 
| 1032 | 
  | 
  | 
        int i; | 
| 1033 | 
  | 
  | 
        double rh_illum = 0.0;  /* Relative horizontal illuminance */ | 
| 1034 | 
  | 
  | 
 | 
| 1035 | 
  | 
  | 
        for (i = 1; i < nskypatch; i++) | 
| 1036 | 
greg | 
2.7 | 
                rh_illum += parr[3*i+1] * rh_cos(i) * rh_dom[i]; | 
| 1037 | 
greg | 
2.1 | 
 | 
| 1038 | 
greg | 
2.7 | 
        return rh_illum; | 
| 1039 | 
greg | 
2.1 | 
} | 
| 1040 | 
  | 
  | 
 | 
| 1041 | 
  | 
  | 
/* Calculate earth orbit eccentricity correction factor */ | 
| 1042 | 
  | 
  | 
 | 
| 1043 | 
  | 
  | 
/* Reference:   Sen, Z. 2008. Solar Energy Fundamental and Modeling  */ | 
| 1044 | 
  | 
  | 
/*                              Techniques. Springer, p. 72. */ | 
| 1045 | 
  | 
  | 
 | 
| 1046 | 
  | 
  | 
double CalcEccentricity() | 
| 1047 | 
  | 
  | 
{ | 
| 1048 | 
  | 
  | 
        double day_angle;       /* Day angle (radians) */ | 
| 1049 | 
  | 
  | 
        double E0;                      /* Eccentricity */ | 
| 1050 | 
  | 
  | 
 | 
| 1051 | 
  | 
  | 
        /* Calculate day angle */ | 
| 1052 | 
  | 
  | 
        day_angle  = (julian_date - 1.0) * (2.0 * PI / 365.0); | 
| 1053 | 
  | 
  | 
 | 
| 1054 | 
  | 
  | 
        /* Calculate eccentricity */ | 
| 1055 | 
  | 
  | 
        E0 = 1.00011 + 0.034221 * cos(day_angle) + 0.00128 * sin(day_angle) | 
| 1056 | 
  | 
  | 
                        + 0.000719 * cos(2.0 * day_angle) + 0.000077 * sin(2.0 * | 
| 1057 | 
  | 
  | 
                        day_angle); | 
| 1058 | 
  | 
  | 
 | 
| 1059 | 
  | 
  | 
        return E0; | 
| 1060 | 
  | 
  | 
} | 
| 1061 | 
  | 
  | 
 | 
| 1062 | 
  | 
  | 
/* Calculate atmospheric precipitable water content */ | 
| 1063 | 
  | 
  | 
 | 
| 1064 | 
  | 
  | 
/* Reference:   Perez, R., P. Ineichen, R. Seals, J. Michalsky, and R. */ | 
| 1065 | 
  | 
  | 
/*                              Stewart. 1990. ìModeling Daylight Availability and */ | 
| 1066 | 
  | 
  | 
/*                              Irradiance Components from Direct and Global */ | 
| 1067 | 
  | 
  | 
/*                              Irradiance,î Solar Energy 44(5):271-289, Eqn. 3. */ | 
| 1068 | 
  | 
  | 
 | 
| 1069 | 
  | 
  | 
/* Note:        The default surface dew point temperature is 11 deg. C */ | 
| 1070 | 
  | 
  | 
/*                      (52 deg. F). Typical values are: */ | 
| 1071 | 
  | 
  | 
 | 
| 1072 | 
  | 
  | 
/*                      Celsius         Fahrenheit              Human Perception */ | 
| 1073 | 
  | 
  | 
/*                      > 24            > 75                    Extremely uncomfortable */ | 
| 1074 | 
  | 
  | 
/*                      21 - 24         70 - 74                 Very humid */ | 
| 1075 | 
  | 
  | 
/*                      18 - 21         65 - 69                 Somewhat uncomfortable */ | 
| 1076 | 
  | 
  | 
/*                      16 - 18         60 - 64                 OK for most people */ | 
| 1077 | 
  | 
  | 
/*                      13 - 16         55 - 59                 Comfortable */ | 
| 1078 | 
  | 
  | 
/*                      10 - 12         50 - 54                 Very comfortable */ | 
| 1079 | 
  | 
  | 
/*                      < 10            < 49                    A bit dry for some */ | 
| 1080 | 
  | 
  | 
 | 
| 1081 | 
  | 
  | 
double CalcPrecipWater( double dpt ) | 
| 1082 | 
  | 
  | 
{ return exp(0.07 * dpt - 0.075); } | 
| 1083 | 
  | 
  | 
 | 
| 1084 | 
  | 
  | 
/* Calculate relative air mass */ | 
| 1085 | 
  | 
  | 
 | 
| 1086 | 
  | 
  | 
/* Reference:   Kasten, F. 1966. "A New Table and Approximation Formula */ | 
| 1087 | 
  | 
  | 
/*                              for the Relative Optical Air Mass," Arch. Meteorol. */ | 
| 1088 | 
  | 
  | 
/*                              Geophys. Bioklimataol. Ser. B14, pp. 206-233. */ | 
| 1089 | 
  | 
  | 
 | 
| 1090 | 
  | 
  | 
/* Note:                More sophisticated relative air mass models are */ | 
| 1091 | 
  | 
  | 
/*                              available, but they differ significantly only for */ | 
| 1092 | 
  | 
  | 
/*                              sun zenith angles greater than 80 degrees. */ | 
| 1093 | 
  | 
  | 
 | 
| 1094 | 
  | 
  | 
double CalcAirMass() | 
| 1095 | 
  | 
  | 
{ | 
| 1096 | 
  | 
  | 
        return (1.0 / (cos(sun_zenith) + 0.15 * pow(93.885 - | 
| 1097 | 
  | 
  | 
                        RadToDeg(sun_zenith), -1.253))); | 
| 1098 | 
  | 
  | 
} | 
| 1099 | 
  | 
  | 
 | 
| 1100 | 
  | 
  | 
/* Calculate Perez All-Weather sky patch luminances (modified by GW) */ | 
| 1101 | 
  | 
  | 
 | 
| 1102 | 
  | 
  | 
/* NOTE: The sky patches centers are determined in accordance with the */ | 
| 1103 | 
  | 
  | 
/*       BRE-IDMP sky luminance measurement procedures. (See for example */ | 
| 1104 | 
  | 
  | 
/*       Mardaljevic, J. 2001. "The BRE-IDMP Dataset: A New Benchmark */ | 
| 1105 | 
  | 
  | 
/*       for the Validation of Illuminance Prediction Techniques," */ | 
| 1106 | 
  | 
  | 
/*       Lighting Research & Technology 33(2):117-136.) */ | 
| 1107 | 
  | 
  | 
 | 
| 1108 | 
  | 
  | 
void CalcSkyPatchLumin( float *parr ) | 
| 1109 | 
  | 
  | 
{ | 
| 1110 | 
  | 
  | 
        int i; | 
| 1111 | 
  | 
  | 
        double aas;                             /* Sun-sky point azimuthal angle */ | 
| 1112 | 
  | 
  | 
        double sspa;                    /* Sun-sky point angle */ | 
| 1113 | 
  | 
  | 
        double zsa;                             /* Zenithal sun angle */ | 
| 1114 | 
  | 
  | 
 | 
| 1115 | 
  | 
  | 
        for (i = 1; i < nskypatch; i++) | 
| 1116 | 
  | 
  | 
        { | 
| 1117 | 
  | 
  | 
                /* Calculate sun-sky point azimuthal angle */ | 
| 1118 | 
  | 
  | 
                aas = fabs(rh_pazi[i] - azimuth); | 
| 1119 | 
  | 
  | 
 | 
| 1120 | 
  | 
  | 
                /* Calculate zenithal sun angle */ | 
| 1121 | 
  | 
  | 
                zsa = PI * 0.5 - rh_palt[i]; | 
| 1122 | 
  | 
  | 
 | 
| 1123 | 
  | 
  | 
                /* Calculate sun-sky point angle (Equation 8-20) */ | 
| 1124 | 
  | 
  | 
                sspa = acos(cos(sun_zenith) * cos(zsa) + sin(sun_zenith) * | 
| 1125 | 
  | 
  | 
                                sin(zsa) * cos(aas)); | 
| 1126 | 
  | 
  | 
 | 
| 1127 | 
  | 
  | 
                /* Calculate patch luminance */ | 
| 1128 | 
  | 
  | 
                parr[3*i] = CalcRelLuminance(sspa, zsa); | 
| 1129 | 
  | 
  | 
                if (parr[3*i] < 0) parr[3*i] = 0; | 
| 1130 | 
  | 
  | 
                parr[3*i+2] = parr[3*i+1] = parr[3*i]; | 
| 1131 | 
  | 
  | 
        } | 
| 1132 | 
  | 
  | 
} |