| 1 |
greg |
2.2 |
#ifndef lint |
| 2 |
greg |
2.3 |
static const char RCSid[] = "$Id: atmos.c,v 2.2 2024/07/19 23:38:28 greg Exp $"; |
| 3 |
greg |
2.2 |
#endif |
| 4 |
|
|
/* |
| 5 |
|
|
The Radiance Software License, Version 2.0 |
| 6 |
|
|
|
| 7 |
|
|
Radiance v5.4 Copyright (c) 1990 to 2022, The Regents of the University of |
| 8 |
|
|
California, through Lawrence Berkeley National Laboratory (subject to receipt |
| 9 |
|
|
of any required approvals from the U.S. Dept. of Energy). All rights reserved. |
| 10 |
|
|
|
| 11 |
|
|
Redistribution and use in source and binary forms, with or without |
| 12 |
|
|
modification, are permitted provided that the following conditions are met: |
| 13 |
|
|
|
| 14 |
|
|
(1) Redistributions of source code must retain the above copyright notice, |
| 15 |
|
|
this list of conditions and the following disclaimer. |
| 16 |
|
|
|
| 17 |
|
|
(2) Redistributions in binary form must reproduce the above copyright |
| 18 |
|
|
notice, this list of conditions and the following disclaimer in the |
| 19 |
|
|
documentation and/or other materials provided with the distribution. |
| 20 |
|
|
|
| 21 |
|
|
(3) Neither the name of the University of California, Lawrence Berkeley |
| 22 |
|
|
National Laboratory, U.S. Dept. of Energy nor the names of its contributors |
| 23 |
|
|
may be used to endorse or promote products derived from this software |
| 24 |
|
|
without specific prior written permission. |
| 25 |
|
|
|
| 26 |
|
|
THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" |
| 27 |
|
|
AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE |
| 28 |
|
|
IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE |
| 29 |
|
|
ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER OR CONTRIBUTORS BE |
| 30 |
|
|
LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR |
| 31 |
|
|
CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF |
| 32 |
|
|
SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS |
| 33 |
|
|
INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN |
| 34 |
|
|
CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) |
| 35 |
|
|
ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE |
| 36 |
|
|
POSSIBILITY OF SUCH DAMAGE. |
| 37 |
|
|
|
| 38 |
|
|
You are under no obligation whatsoever to provide any bug fixes, patches, |
| 39 |
|
|
or upgrades to the features, functionality or performance of the source |
| 40 |
|
|
code ("Enhancements") to anyone; however, if you choose to make your |
| 41 |
|
|
Enhancements available either publicly, or directly to Lawrence Berkeley |
| 42 |
|
|
National Laboratory, without imposing a separate written license agreement |
| 43 |
|
|
for such Enhancements, then you hereby grant the following license: a |
| 44 |
|
|
non-exclusive, royalty-free perpetual license to install, use, modify, |
| 45 |
|
|
prepare derivative works, incorporate into other computer software, |
| 46 |
|
|
distribute, and sublicense such enhancements or derivative works thereof, |
| 47 |
|
|
in binary and source code form. |
| 48 |
|
|
|
| 49 |
|
|
================================================================================= |
| 50 |
|
|
|
| 51 |
|
|
Copyright (c) 2017 Eric Bruneton |
| 52 |
|
|
All rights reserved. |
| 53 |
|
|
|
| 54 |
|
|
Redistribution and use in source and binary forms, with or without |
| 55 |
|
|
modification, are permitted provided that the following conditions are met: |
| 56 |
|
|
|
| 57 |
|
|
1. Redistributions of source code must retain the above copyright notice, this |
| 58 |
|
|
list of conditions and the following disclaimer. |
| 59 |
|
|
|
| 60 |
|
|
2. Redistributions in binary form must reproduce the above copyright notice, |
| 61 |
|
|
this list of conditions and the following disclaimer in the documentation |
| 62 |
|
|
and/or other materials provided with the distribution. |
| 63 |
|
|
|
| 64 |
|
|
3. Neither the name of the copyright holder nor the names of its contributors |
| 65 |
|
|
may be used to endorse or promote products derived from this software without |
| 66 |
|
|
specific prior written permission. |
| 67 |
|
|
|
| 68 |
|
|
THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" |
| 69 |
|
|
AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE |
| 70 |
|
|
IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE |
| 71 |
|
|
DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE LIABLE |
| 72 |
|
|
FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL |
| 73 |
|
|
DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR |
| 74 |
|
|
SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER |
| 75 |
|
|
CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, |
| 76 |
|
|
OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE |
| 77 |
|
|
OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE. |
| 78 |
|
|
|
| 79 |
|
|
Additional modifications by T. Wang |
| 80 |
|
|
- N-dimensional data storage and lookup using Radiance data structure |
| 81 |
|
|
- Mie scattering coefficients precomputed using libRadTrans and interpolated at runtime |
| 82 |
|
|
- Translation to C from C++ with multi-threading |
| 83 |
|
|
*/ |
| 84 |
|
|
|
| 85 |
greg |
2.1 |
#include "atmos.h" |
| 86 |
|
|
#include "data.h" |
| 87 |
|
|
|
| 88 |
|
|
#if defined(_WIN32) || defined(_WIN64) |
| 89 |
|
|
#include <windows.h> |
| 90 |
|
|
#else |
| 91 |
|
|
#include <pthread.h> |
| 92 |
|
|
#endif |
| 93 |
|
|
|
| 94 |
|
|
#if defined(_WIN32) || defined(_WIN64) |
| 95 |
|
|
#define THREAD HANDLE |
| 96 |
|
|
#define CREATE_THREAD(thread, func, context) \ |
| 97 |
|
|
((*(thread) = CreateThread(NULL, 0, (func), \ |
| 98 |
|
|
(context), 0, NULL)) != NULL) |
| 99 |
|
|
#define THREAD_RETURN DWORD WINAPI |
| 100 |
|
|
#define THREAD_JOIN(thread) WaitForSingleObject(thread, INFINITE) |
| 101 |
|
|
#define THREAD_CLOSE(thread) CloseHandle(thread) |
| 102 |
|
|
#else |
| 103 |
|
|
#define THREAD pthread_t |
| 104 |
|
|
#define CREATE_THREAD(thread, func, context) \ |
| 105 |
|
|
(pthread_create((thread), NULL, (func), (context)) == 0) |
| 106 |
|
|
#define THREAD_RETURN void * |
| 107 |
|
|
#define THREAD_JOIN(thread) pthread_join(thread, NULL) |
| 108 |
|
|
#endif |
| 109 |
|
|
|
| 110 |
|
|
typedef struct { |
| 111 |
|
|
int start_l; |
| 112 |
|
|
int end_l; |
| 113 |
|
|
const Atmosphere *atmos; |
| 114 |
|
|
DATARRAY *tau_dp; |
| 115 |
|
|
DATARRAY *delta_rayleigh_scattering_dp; |
| 116 |
|
|
DATARRAY *delta_mie_scattering_dp; |
| 117 |
|
|
DATARRAY *scattering_dp; |
| 118 |
|
|
} Scat1Tdat; |
| 119 |
|
|
|
| 120 |
|
|
typedef struct { |
| 121 |
|
|
int start_l; |
| 122 |
|
|
int end_l; |
| 123 |
|
|
DATARRAY *tau_dp; |
| 124 |
|
|
DATARRAY *delta_scattering_density_dp; |
| 125 |
|
|
DATARRAY *delta_multiple_scattering_dp; |
| 126 |
|
|
DATARRAY *scattering_dp; |
| 127 |
|
|
} ScatNTdat; |
| 128 |
|
|
|
| 129 |
|
|
typedef struct { |
| 130 |
|
|
int start_l; |
| 131 |
|
|
int end_l; |
| 132 |
|
|
const Atmosphere *atmos; |
| 133 |
|
|
DATARRAY *tau_dp; |
| 134 |
|
|
DATARRAY *delta_rayleigh_scattering_dp; |
| 135 |
|
|
DATARRAY *delta_mie_scattering_dp; |
| 136 |
|
|
DATARRAY *delta_multiple_scattering_dp; |
| 137 |
|
|
DATARRAY *delta_irradiance_dp; |
| 138 |
|
|
int scattering_order; |
| 139 |
|
|
DATARRAY *delta_scattering_density_dp; |
| 140 |
|
|
} ScatDenseTdat; |
| 141 |
|
|
|
| 142 |
greg |
2.2 |
const double ER = 6360000.0; /* Earth radius in meters */ |
| 143 |
|
|
const double AR = 6420000; /* Atmosphere radius in meters */ |
| 144 |
|
|
const double AH = 60000; /* Atmosphere thickness in meters */ |
| 145 |
greg |
2.1 |
|
| 146 |
|
|
const float START_WVL = 390.0; |
| 147 |
|
|
const float END_WVL = 770.0; |
| 148 |
|
|
const double SUNRAD = 0.004625; |
| 149 |
|
|
|
| 150 |
greg |
2.2 |
/* The cosine of the maximum Sun zenith angle */ |
| 151 |
greg |
2.1 |
const double MU_S_MIN = -0.2; |
| 152 |
|
|
|
| 153 |
greg |
2.2 |
/* Scale heights (m) */ |
| 154 |
|
|
/* Rayleigh scattering */ |
| 155 |
|
|
const double HR_MS = 8820; /* Midlatitude summer */ |
| 156 |
|
|
const double HR_MW = 8100; /* Midlatitude winter */ |
| 157 |
|
|
const double HR_SS = 8550; /* Subarctic summer */ |
| 158 |
|
|
const double HR_SW = 7600; /* Subarctic winter */ |
| 159 |
|
|
const double HR_T = 9050; /* Tropical */ |
| 160 |
greg |
2.1 |
|
| 161 |
greg |
2.2 |
const double AOD0_CA = 0.115; /* Broadband AOD for continental average */ |
| 162 |
greg |
2.1 |
|
| 163 |
greg |
2.2 |
const double SOLOMG = 6.7967e-5; /* .533 apex angle */ |
| 164 |
|
|
const int WVLSPAN = 400; /* 380nm to 780nm */ |
| 165 |
greg |
2.1 |
|
| 166 |
greg |
2.2 |
/* Aerosol optical depth at 550nm for continental clean */ |
| 167 |
greg |
2.1 |
const double AOD_CC_550 = 0.05352; |
| 168 |
|
|
|
| 169 |
|
|
/** 380nm to 780nm at 20nm intervals **/ |
| 170 |
greg |
2.2 |
/* Extraterrestrial solar W/m^2/nm */ |
| 171 |
greg |
2.1 |
const float EXTSOL[NSSAMP] = {1.11099345, 1.75363199, 1.67126921, 1.97235667, |
| 172 |
|
|
2.01863015, 1.93612482, 1.87310758, 1.88792588, |
| 173 |
|
|
1.86274213, 1.84433248, 1.80370942, 1.73020456, |
| 174 |
|
|
1.67036654, 1.57482149, 1.53646383, 1.45515319, |
| 175 |
|
|
1.38207435, 1.3265141, 1.26648111, 1.20592043}; |
| 176 |
|
|
|
| 177 |
greg |
2.2 |
/* Rayleight scattering coefficients at sea level in m^-1 */ |
| 178 |
greg |
2.1 |
const float BR0_MS[NSSAMP] = { |
| 179 |
|
|
4.63889967e-05, 3.76703293e-05, 3.09193474e-05, 2.56223081e-05, |
| 180 |
|
|
2.14198405e-05, 1.80483398e-05, 1.53177709e-05, 1.30867298e-05, |
| 181 |
|
|
1.12487269e-05, 9.72370830e-06, 8.44932642e-06, 7.37769845e-06, |
| 182 |
|
|
6.47116653e-06, 5.70005327e-06, 5.04091603e-06, 4.47430240e-06, |
| 183 |
|
|
3.98549839e-06, 3.56178512e-06, 3.19293761e-06, 2.87072599e-06}; |
| 184 |
|
|
|
| 185 |
|
|
const float BR0_MW[NSSAMP] = { |
| 186 |
|
|
5.03857564e-05, 4.09159104e-05, 3.35832809e-05, 2.78298619e-05, |
| 187 |
|
|
2.32653203e-05, 1.96033395e-05, 1.66375116e-05, 1.42142496e-05, |
| 188 |
|
|
1.22178890e-05, 1.05614786e-05, 9.17729921e-06, 8.01334246e-06, |
| 189 |
|
|
7.02870602e-06, 6.19115557e-06, 5.47522872e-06, 4.85979708e-06, |
| 190 |
|
|
4.32887894e-06, 3.86865960e-06, 3.46803311e-06, 3.11806054e-06}; |
| 191 |
|
|
|
| 192 |
|
|
const float BR0_SS[NSSAMP] = { |
| 193 |
|
|
4.73789047e-05, 3.84741872e-05, 3.15791442e-05, 2.61690698e-05, |
| 194 |
|
|
2.18769246e-05, 1.84334785e-05, 1.56446412e-05, 1.33659912e-05, |
| 195 |
|
|
1.14887667e-05, 9.93120528e-06, 8.62962901e-06, 7.53513326e-06, |
| 196 |
|
|
6.60925660e-06, 5.82168833e-06, 5.14848557e-06, 4.56978082e-06, |
| 197 |
|
|
4.07054608e-06, 3.63779107e-06, 3.26107261e-06, 2.93198523e-06}; |
| 198 |
|
|
|
| 199 |
|
|
const float BR0_SW[NSSAMP] = { |
| 200 |
|
|
5.30623659e-05, 4.30894595e-05, 3.53673035e-05, 2.93082494e-05, |
| 201 |
|
|
2.45012287e-05, 2.06447149e-05, 1.75213353e-05, 1.49693438e-05, |
| 202 |
|
|
1.28669320e-05, 1.11225291e-05, 9.66481886e-06, 8.43902999e-06, |
| 203 |
|
|
7.40208736e-06, 6.52004427e-06, 5.76608570e-06, 5.11796089e-06, |
| 204 |
|
|
4.55883914e-06, 4.07417187e-06, 3.65226316e-06, 3.28369923e-06}; |
| 205 |
|
|
const float BR0_T[NSSAMP] = { |
| 206 |
|
|
4.55376661e-05, 3.69790036e-05, 3.03519157e-05, 2.51520876e-05, |
| 207 |
|
|
2.10267438e-05, 1.77171168e-05, 1.50366593e-05, 1.28465622e-05, |
| 208 |
|
|
1.10422904e-05, 9.54525887e-06, 8.29426444e-06, 7.24230298e-06, |
| 209 |
|
|
6.35240772e-06, 5.59544593e-06, 4.94840517e-06, 4.39219003e-06, |
| 210 |
|
|
3.91235655e-06, 3.49641927e-06, 3.13434084e-06, 2.81804244e-06}; |
| 211 |
|
|
|
| 212 |
greg |
2.2 |
/* Ozone extinction coefficients at sea level */ |
| 213 |
greg |
2.1 |
const float AO0[NSSAMP] = { |
| 214 |
|
|
3.5661864e-09, 1.4848362e-08, 4.5415674e-08, 1.2446184e-07, 2.6461576e-07, |
| 215 |
|
|
4.8451984e-07, 8.609148e-07, 1.480537e-06, 1.8809e-06, 2.5107328e-06, |
| 216 |
|
|
2.5263174e-06, 2.313507e-06, 1.727741e-06, 1.2027012e-06, 7.915902e-07, |
| 217 |
|
|
5.0639202e-07, 3.5285684e-07, 2.23021e-07, 1.7395638e-07, 1.5052574e-07}; |
| 218 |
|
|
|
| 219 |
greg |
2.2 |
const double MIE_G = 0.7; /* Mie eccentricity */ |
| 220 |
greg |
2.1 |
|
| 221 |
|
|
const int TRANSMITTANCE_TEXTURE_WIDTH = 256; |
| 222 |
|
|
const int TRANSMITTANCE_TEXTURE_HEIGHT = 64; |
| 223 |
|
|
|
| 224 |
greg |
2.2 |
/* Resolution halved by original impl */ |
| 225 |
greg |
2.1 |
const int SCATTERING_TEXTURE_R_SIZE = 16; |
| 226 |
|
|
const int SCATTERING_TEXTURE_MU_SIZE = 64; |
| 227 |
|
|
const int SCATTERING_TEXTURE_MU_S_SIZE = 16; |
| 228 |
|
|
const int SCATTERING_TEXTURE_NU_SIZE = 4; |
| 229 |
|
|
|
| 230 |
|
|
const int IRRADIANCE_TEXTURE_WIDTH = 64; |
| 231 |
|
|
const int IRRADIANCE_TEXTURE_HEIGHT = 16; |
| 232 |
|
|
|
| 233 |
|
|
static inline double clamp_cosine(double mu) { |
| 234 |
|
|
return fmax(-1.0, fmin(1.0, mu)); |
| 235 |
|
|
} |
| 236 |
|
|
|
| 237 |
|
|
static inline double clamp_distance(double d) { return fmax(d, 0.0); } |
| 238 |
|
|
|
| 239 |
|
|
static inline double clamp_radius(double r) { return fmax(ER, fmin(r, AR)); } |
| 240 |
|
|
|
| 241 |
|
|
static inline double safe_sqrt(double a) { return sqrt(fmax(a, 0.0)); } |
| 242 |
|
|
|
| 243 |
|
|
static inline double distance_to_space(const double radius, const double mu) { |
| 244 |
|
|
const double descriminant = radius * radius * (mu * mu - 1.0) + AR * AR; |
| 245 |
|
|
return clamp_distance(-radius * mu + safe_sqrt(descriminant)); |
| 246 |
|
|
} |
| 247 |
|
|
|
| 248 |
|
|
static inline double distance_to_earth(const double radius, const double mu) { |
| 249 |
|
|
const double descriminant = radius * radius * (mu * mu - 1.0) + ER * ER; |
| 250 |
|
|
return clamp_distance(-radius * mu - safe_sqrt(descriminant)); |
| 251 |
|
|
} |
| 252 |
|
|
|
| 253 |
|
|
static inline int ray_hits_ground(double radius, double mu) { |
| 254 |
|
|
return mu < 0.0 && radius * radius * (mu * mu - 1.0) + ER * ER >= 0.0; |
| 255 |
|
|
} |
| 256 |
|
|
|
| 257 |
|
|
static inline double get_layer_density(const DensityProfileLayer *layer, |
| 258 |
|
|
const double altitude) { |
| 259 |
|
|
double density = layer->exp_term * exp(layer->exp_scale * altitude) + |
| 260 |
|
|
layer->linear_term * altitude + layer->constant_term; |
| 261 |
|
|
return fmax(0.0, fmin(1.0, density)); |
| 262 |
|
|
} |
| 263 |
|
|
|
| 264 |
|
|
static inline double get_profile_density(const DensityProfile *profile, |
| 265 |
|
|
const double altitude) { |
| 266 |
|
|
if (altitude < profile->layers[0].width) { |
| 267 |
|
|
return get_layer_density(&(profile->layers[0]), altitude); |
| 268 |
|
|
} else { |
| 269 |
|
|
return get_layer_density(&(profile->layers[1]), altitude); |
| 270 |
|
|
} |
| 271 |
|
|
} |
| 272 |
|
|
|
| 273 |
|
|
static int compute_optical_length_to_space_mie(DATARRAY *mdp, const double radius, |
| 274 |
|
|
const double mu, double *result) { |
| 275 |
greg |
2.2 |
/* Number of intervals for the numerical integration. */ |
| 276 |
greg |
2.1 |
const int nsamp = 500; |
| 277 |
greg |
2.2 |
/* The integration step, i.e. the length of each integration interval. */ |
| 278 |
greg |
2.1 |
const double dx = distance_to_space(radius, mu) / nsamp; |
| 279 |
greg |
2.2 |
/* Integration loop. */ |
| 280 |
greg |
2.1 |
for (int i = 0; i <= nsamp; ++i) { |
| 281 |
|
|
double d_i = i * dx; |
| 282 |
greg |
2.2 |
/* Distance between the current sample point and the planet center. */ |
| 283 |
greg |
2.1 |
double r_i = sqrt(d_i * d_i + 2.0 * radius * mu * d_i + radius * radius); |
| 284 |
|
|
double pt[1] = {r_i - ER}; |
| 285 |
|
|
DATARRAY *mie = datavector(mdp, pt); |
| 286 |
|
|
double weight_i = i == 0 || i == nsamp ? 0.5 : 1.0; |
| 287 |
|
|
for (int j = 0; j < NSSAMP; ++j) { |
| 288 |
|
|
result[j] += mie->arr.d[j] * weight_i * dx; |
| 289 |
|
|
} |
| 290 |
|
|
free(mie); |
| 291 |
|
|
} |
| 292 |
|
|
return 1; |
| 293 |
|
|
} |
| 294 |
|
|
|
| 295 |
|
|
static double compute_optical_length_to_space(const DensityProfile *profile, |
| 296 |
|
|
const double radius, const double mu) { |
| 297 |
greg |
2.2 |
/* Number of intervals for the numerical integration. */ |
| 298 |
greg |
2.1 |
const int SAMPLE_COUNT = 500; |
| 299 |
greg |
2.2 |
/* The integration step, i.e. the length of each integration interval. */ |
| 300 |
greg |
2.1 |
const double dx = distance_to_space(radius, mu) / SAMPLE_COUNT; |
| 301 |
greg |
2.2 |
/* Integration loop. */ |
| 302 |
greg |
2.1 |
double result = 0.0; |
| 303 |
|
|
for (int i = 0; i <= SAMPLE_COUNT; ++i) { |
| 304 |
|
|
double d_i = i * dx; |
| 305 |
greg |
2.2 |
/* Distance between the current sample point and the planet center. */ |
| 306 |
greg |
2.1 |
double r_i = sqrt(d_i * d_i + 2.0 * radius * mu * d_i + radius * radius); |
| 307 |
greg |
2.2 |
/* Number density at the current sample point (divided by the number density */ |
| 308 |
|
|
/* at the bottom of the atmosphere, yielding a dimensionless number). */ |
| 309 |
greg |
2.1 |
double y_i = get_profile_density(profile, r_i - ER); |
| 310 |
greg |
2.2 |
/* Sample weight (from the trapezoidal rule). */ |
| 311 |
greg |
2.1 |
double weight_i = i == 0 || i == SAMPLE_COUNT ? 0.5 : 1.0; |
| 312 |
|
|
result += y_i * weight_i * dx; |
| 313 |
|
|
} |
| 314 |
|
|
return result; |
| 315 |
|
|
} |
| 316 |
|
|
|
| 317 |
|
|
static void compute_transmittance_to_space(const Atmosphere *atmos, |
| 318 |
|
|
const double radius, const double mu, |
| 319 |
|
|
float *result) { |
| 320 |
|
|
const double taur = |
| 321 |
|
|
compute_optical_length_to_space(&atmos->rayleigh_density, radius, mu); |
| 322 |
|
|
double taum[NSSAMP] = {0}; |
| 323 |
|
|
compute_optical_length_to_space_mie(atmos->beta_m, radius, mu, taum); |
| 324 |
|
|
const double tauo = |
| 325 |
|
|
compute_optical_length_to_space(&atmos->ozone_density, radius, mu); |
| 326 |
|
|
for (int i = 0; i < NSSAMP; ++i) { |
| 327 |
|
|
result[i] = exp(-(taur * atmos->beta_r0[i] + taum[i] * atmos->beta_scale + |
| 328 |
|
|
tauo * AO0[i])); |
| 329 |
|
|
} |
| 330 |
|
|
} |
| 331 |
|
|
|
| 332 |
|
|
static inline double get_texture_coord_from_unit_range(const double x, |
| 333 |
|
|
const int texture_size) { |
| 334 |
|
|
return 0.5 / (double)texture_size + x * (1.0 - 1.0 / (double)texture_size); |
| 335 |
|
|
} |
| 336 |
|
|
|
| 337 |
|
|
static inline double get_unit_range_from_texture_coord(const double u, |
| 338 |
|
|
const int texture_size) { |
| 339 |
|
|
return (u - 0.5 / (double)texture_size) / (1.0 - 1.0 / (double)texture_size); |
| 340 |
|
|
} |
| 341 |
|
|
|
| 342 |
|
|
static void to_transmittance_uv(const double radius, const double mu, double *u, |
| 343 |
|
|
double *v) { |
| 344 |
greg |
2.2 |
/* Distance to top atmosphere boundary for a horizontal ray at ground level. */ |
| 345 |
greg |
2.1 |
double H = sqrt(AR * AR - ER * ER); |
| 346 |
greg |
2.2 |
/* Distance to the horizon. */ |
| 347 |
greg |
2.1 |
double rho = safe_sqrt(radius * radius - ER * ER); |
| 348 |
greg |
2.2 |
/* Distance to the top atmosphere boundary for the ray (r,mu), and its minimum */ |
| 349 |
|
|
/* and maximum values over all mu - obtained for (r,1) and (r,mu_horizon). */ |
| 350 |
greg |
2.1 |
double d = distance_to_space(radius, mu); |
| 351 |
|
|
double d_min = AR - radius; |
| 352 |
|
|
double d_max = rho + H; |
| 353 |
|
|
double x_mu = (d - d_min) / (d_max - d_min); |
| 354 |
|
|
double x_r = rho / H; |
| 355 |
|
|
*u = fmax(0.0, fmin(1.0, x_mu)); |
| 356 |
|
|
*v = fmax(0.0, fmin(1.0, x_r)); |
| 357 |
|
|
} |
| 358 |
|
|
|
| 359 |
|
|
static void from_transmittance_uv(const double u, const double v, double *radius, |
| 360 |
|
|
double *mu) { |
| 361 |
|
|
double x_mu = u; |
| 362 |
|
|
double x_r = v; |
| 363 |
greg |
2.2 |
/* Distance to top atmosphere boundary for a horizontal ray at ground level. */ |
| 364 |
greg |
2.1 |
double H = sqrt(AR * AR - ER * ER); |
| 365 |
greg |
2.2 |
/* Distance to the horizon, from which we can compute r: */ |
| 366 |
greg |
2.1 |
double rho = H * x_r; |
| 367 |
|
|
*radius = sqrt(rho * rho + ER * ER); |
| 368 |
greg |
2.2 |
/* Distance to the top atmosphere boundary for the ray (r,mu), and its minimum */ |
| 369 |
|
|
/* and maximum values over all mu - obtained for (r,1) and (r,mu_horizon) - */ |
| 370 |
|
|
/* from which we can recover mu: */ |
| 371 |
greg |
2.1 |
double d_min = AR - *radius; |
| 372 |
|
|
double d_max = rho + H; |
| 373 |
|
|
double d = d_min + x_mu * (d_max - d_min); |
| 374 |
|
|
*mu = d == 0.0 ? 1.0 : (H * H - rho * rho - d * d) / (2.0 * *radius * d); |
| 375 |
|
|
*mu = clamp_cosine(*mu); |
| 376 |
|
|
} |
| 377 |
|
|
|
| 378 |
|
|
static DATARRAY *get_transmittance_to_space(DATARRAY *dp, const double radius, |
| 379 |
|
|
const double mu) { |
| 380 |
|
|
double u, v; |
| 381 |
|
|
to_transmittance_uv(radius, mu, &u, &v); |
| 382 |
|
|
|
| 383 |
|
|
double pt[2] = {v, u}; |
| 384 |
|
|
return datavector(dp, pt); |
| 385 |
|
|
} |
| 386 |
|
|
|
| 387 |
|
|
static void get_transmittance(DATARRAY *tau_dp, const double radius, const double mu, const double d, |
| 388 |
|
|
const int intersects_ground, double *result) { |
| 389 |
|
|
|
| 390 |
|
|
DATARRAY *result1; |
| 391 |
|
|
DATARRAY *result2; |
| 392 |
|
|
|
| 393 |
|
|
const double r_d = clamp_radius(sqrt(d * d + 2.0 * radius * mu * d + radius * radius)); |
| 394 |
|
|
const double mu_d = clamp_cosine((radius * mu + d) / r_d); |
| 395 |
|
|
|
| 396 |
|
|
if (intersects_ground) { |
| 397 |
|
|
result1 = get_transmittance_to_space(tau_dp, r_d, -mu_d); |
| 398 |
|
|
result2 = get_transmittance_to_space(tau_dp, radius, -mu); |
| 399 |
|
|
} else { |
| 400 |
|
|
result1 = get_transmittance_to_space(tau_dp, radius, mu); |
| 401 |
|
|
result2 = get_transmittance_to_space(tau_dp, r_d, mu_d); |
| 402 |
|
|
} |
| 403 |
|
|
for (int i = 0; i < NSSAMP; ++i) { |
| 404 |
|
|
result[i] = fmin(result1->arr.d[i] / result2->arr.d[i], 1.0); |
| 405 |
|
|
} |
| 406 |
|
|
free(result1); |
| 407 |
|
|
free(result2); |
| 408 |
|
|
} |
| 409 |
|
|
|
| 410 |
|
|
void get_transmittance_to_sun(DATARRAY *tau_dp, const double r, |
| 411 |
|
|
const double mu_s, double *result) { |
| 412 |
|
|
const double sin_theta_h = ER / r; |
| 413 |
|
|
const double cos_theta_h = -sqrt(fmax(1.0 - sin_theta_h * sin_theta_h, 0.0)); |
| 414 |
|
|
DATARRAY *tau_sun = get_transmittance_to_space(tau_dp, r, mu_s); |
| 415 |
greg |
2.2 |
/* Smooth step */ |
| 416 |
greg |
2.1 |
const double edge0 = -sin_theta_h * SUNRAD; |
| 417 |
|
|
const double edge1 = edge0 * -1; |
| 418 |
|
|
double x = mu_s - cos_theta_h; |
| 419 |
|
|
x = (x - edge0) / (edge1 - edge0); |
| 420 |
|
|
x = x < 0.0 ? 0.0 : x > 1.0 ? 1.0 : x; |
| 421 |
|
|
const double st = x * x * (3.0 - 2.0 * x); |
| 422 |
|
|
for (int i = 0; i < NSSAMP; ++i) { |
| 423 |
|
|
result[i] = tau_sun->arr.d[i] * st; |
| 424 |
|
|
} |
| 425 |
|
|
free(tau_sun); |
| 426 |
|
|
} |
| 427 |
|
|
|
| 428 |
|
|
static void compute_single_scattering_integrand( |
| 429 |
|
|
const Atmosphere *atmos, DATARRAY *tau_dp, const double radius, const double mu, |
| 430 |
|
|
const double mu_s, const double nu, const double distance, |
| 431 |
|
|
const int r_mu_intersects_ground, double *rayleigh, double *mie) { |
| 432 |
|
|
|
| 433 |
|
|
const double radius_d = clamp_radius(sqrt(distance * distance + 2.0 * radius * mu * distance + radius * radius)); |
| 434 |
|
|
const double mu_s_d = clamp_cosine((radius * mu_s + distance * nu) / radius_d); |
| 435 |
greg |
2.2 |
/* double transmittance[NSSAMP]; */ |
| 436 |
greg |
2.1 |
double tau_r_mu[NSSAMP] = {0}; |
| 437 |
|
|
double tau_sun[NSSAMP] = {0}; |
| 438 |
|
|
get_transmittance(tau_dp, radius, mu, distance, r_mu_intersects_ground, tau_r_mu); |
| 439 |
|
|
get_transmittance_to_sun(tau_dp, radius_d, mu_s_d, tau_sun); |
| 440 |
|
|
const double rayleigh_profile_density = |
| 441 |
|
|
get_profile_density(&atmos->rayleigh_density, radius_d - ER); |
| 442 |
|
|
DATARRAY *mie_scat; |
| 443 |
|
|
double pt[1] = {radius_d - ER}; |
| 444 |
|
|
mie_scat = datavector(atmos->beta_m, pt); |
| 445 |
|
|
for (int i = 0; i < NSSAMP; ++i) { |
| 446 |
|
|
double transmittance = tau_r_mu[i] * tau_sun[i]; |
| 447 |
|
|
rayleigh[i] = transmittance * rayleigh_profile_density; |
| 448 |
|
|
mie[i] = transmittance * mie_scat->arr.d[i]; |
| 449 |
|
|
} |
| 450 |
|
|
free(mie_scat); |
| 451 |
|
|
} |
| 452 |
|
|
|
| 453 |
|
|
static double |
| 454 |
|
|
distance_to_nearst_atmosphere_boundary(const double radius, const double mu, |
| 455 |
|
|
const int hits_ground) { |
| 456 |
|
|
if (hits_ground) { |
| 457 |
|
|
return distance_to_earth(radius, mu); |
| 458 |
|
|
} else { |
| 459 |
|
|
return distance_to_space(radius, mu); |
| 460 |
|
|
} |
| 461 |
|
|
} |
| 462 |
|
|
|
| 463 |
|
|
static void compute_single_scattering(const Atmosphere *atmos, DATARRAY *tau_dp, |
| 464 |
|
|
const double radius, const double mu, |
| 465 |
|
|
const double mu_s, const double nu, |
| 466 |
|
|
const int hits_ground, |
| 467 |
|
|
float *rayleigh, float *mie) { |
| 468 |
|
|
|
| 469 |
|
|
const int nsamp = 50; |
| 470 |
|
|
double distance_to_boundary; |
| 471 |
|
|
|
| 472 |
|
|
if (hits_ground) { |
| 473 |
|
|
distance_to_boundary = distance_to_earth(radius, mu); |
| 474 |
|
|
} else { |
| 475 |
|
|
distance_to_boundary = distance_to_space(radius, mu); |
| 476 |
|
|
} |
| 477 |
|
|
const double dx = distance_to_boundary / (double)nsamp; |
| 478 |
|
|
|
| 479 |
|
|
double rayleigh_sum[NSSAMP] = {0}; |
| 480 |
|
|
double mie_sum[NSSAMP] = {0}; |
| 481 |
|
|
|
| 482 |
|
|
for (int i = 0; i <= nsamp; ++i) { |
| 483 |
|
|
const double d_i = i * dx; |
| 484 |
|
|
double rayleigh_i[NSSAMP] = {0}; |
| 485 |
|
|
double mie_i[NSSAMP] = {0}; |
| 486 |
|
|
compute_single_scattering_integrand(atmos, tau_dp, radius, mu, mu_s, nu, d_i, |
| 487 |
|
|
hits_ground, rayleigh_i, |
| 488 |
|
|
mie_i); |
| 489 |
|
|
double weight_i = (i == 0 || i == nsamp) ? 0.5 : 1.0; |
| 490 |
|
|
for (int j = 0; j < NSSAMP; ++j) { |
| 491 |
|
|
rayleigh_sum[j] += rayleigh_i[j] * weight_i; |
| 492 |
|
|
mie_sum[j] += mie_i[j] * weight_i; |
| 493 |
|
|
} |
| 494 |
|
|
} |
| 495 |
|
|
for (int i = 0; i < NSSAMP; ++i) { |
| 496 |
|
|
rayleigh[i] = rayleigh_sum[i] * dx * EXTSOL[i] * atmos->beta_r0[i]; |
| 497 |
|
|
mie[i] = dx * EXTSOL[i] * (mie_sum[i] * atmos->beta_scale); |
| 498 |
|
|
} |
| 499 |
|
|
} |
| 500 |
|
|
|
| 501 |
|
|
inline static double rayleigh_phase_function(const double nu) { |
| 502 |
|
|
const double k = 3.0 / (16.0 * PI); |
| 503 |
|
|
return k * (1.0 + nu * nu); |
| 504 |
|
|
} |
| 505 |
|
|
|
| 506 |
|
|
inline static double mie_phase_function(const double g, const double nu) { |
| 507 |
|
|
const double k = 3.0 / (8.0 * PI) * (1.0 - g * g) / (2.0 + g * g); |
| 508 |
|
|
return k * (1.0 + nu * nu) / pow(1.0 + g * g - 2.0 * g * nu, 1.5); |
| 509 |
|
|
} |
| 510 |
|
|
|
| 511 |
|
|
inline static double cloud_phase_function(const double g, const double nu) { |
| 512 |
|
|
const double alpha = 0.5; |
| 513 |
|
|
const double g1 = -0.5; |
| 514 |
|
|
const double hg0 = mie_phase_function(g, nu); |
| 515 |
|
|
const double hg1 = mie_phase_function(g1, nu); |
| 516 |
|
|
return alpha * hg0 + (1.0 - alpha) * hg1; |
| 517 |
|
|
} |
| 518 |
|
|
|
| 519 |
|
|
static void to_scattering_uvwz(const double radius, const double mu, const double mu_s, const double nu, |
| 520 |
|
|
const int r_mu_hits_ground, double *u, double *v, |
| 521 |
|
|
double *w, double *z) { |
| 522 |
|
|
|
| 523 |
|
|
double H = sqrt(AR * AR - ER * ER); |
| 524 |
|
|
double rho = safe_sqrt(radius * radius - ER * ER); |
| 525 |
|
|
double u_r = rho / H; |
| 526 |
|
|
|
| 527 |
|
|
double r_mu = radius * mu; |
| 528 |
|
|
double discriminant = r_mu * r_mu - radius * radius + ER * ER; |
| 529 |
|
|
|
| 530 |
|
|
double u_mu; |
| 531 |
|
|
if (r_mu_hits_ground) { |
| 532 |
|
|
double d = -r_mu - safe_sqrt(discriminant); |
| 533 |
|
|
double d_min = radius - ER; |
| 534 |
|
|
double d_max = rho; |
| 535 |
|
|
u_mu = 0.5 - 0.5 * get_texture_coord_from_unit_range( |
| 536 |
|
|
d_max == d_min ? 0.0 : (d - d_min) / (d_max - d_min), |
| 537 |
|
|
SCATTERING_TEXTURE_MU_SIZE / 2); |
| 538 |
|
|
} else { |
| 539 |
|
|
double d = -r_mu + sqrt(discriminant + H * H); |
| 540 |
|
|
double d_min = AR - radius; |
| 541 |
|
|
double d_max = rho + H; |
| 542 |
|
|
u_mu = 0.5 + 0.5 * get_texture_coord_from_unit_range( |
| 543 |
|
|
(d - d_min) / (d_max - d_min), |
| 544 |
|
|
SCATTERING_TEXTURE_MU_SIZE / 2); |
| 545 |
|
|
} |
| 546 |
|
|
|
| 547 |
|
|
double d = distance_to_space(ER, mu_s); |
| 548 |
|
|
double d_min = AH; |
| 549 |
|
|
double d_max = H; |
| 550 |
|
|
double a = (d - d_min) / (d_max - d_min); |
| 551 |
|
|
double D = distance_to_space(ER, MU_S_MIN); |
| 552 |
|
|
double A = (D - d_min) / (d_max - d_min); |
| 553 |
|
|
double u_mu_s = get_texture_coord_from_unit_range( |
| 554 |
|
|
fmax(1.0 - a / A, 0.0) / (1.0 + a), SCATTERING_TEXTURE_MU_S_SIZE); |
| 555 |
|
|
double u_nu = (nu + 1.0) / 2; |
| 556 |
|
|
*u = u_nu; |
| 557 |
|
|
*v = u_mu_s; |
| 558 |
|
|
*w = u_mu; |
| 559 |
|
|
*z = u_r; |
| 560 |
|
|
} |
| 561 |
|
|
|
| 562 |
|
|
static void from_scattering_uvwz(const double x, const double y, const double z, const double w, |
| 563 |
|
|
double *radius, double *mu, double *mu_s, |
| 564 |
|
|
double *nu, int *r_mu_hits_ground) { |
| 565 |
|
|
|
| 566 |
|
|
|
| 567 |
|
|
const double H = sqrt(AR * AR - ER * ER); |
| 568 |
|
|
const double rho = H * w; |
| 569 |
|
|
*radius = sqrt(rho * rho + ER * ER); |
| 570 |
|
|
|
| 571 |
|
|
if (z < 0.5) { |
| 572 |
greg |
2.2 |
/* Distance to the ground for the ray (r,mu), and its minimum and maximum */ |
| 573 |
|
|
/* values over all mu - obtained for (r,-1) and (r,mu_horizon) - from which */ |
| 574 |
|
|
/* we can recover mu: */ |
| 575 |
greg |
2.1 |
const double d_min = *radius - ER; |
| 576 |
|
|
const double d_max = rho; |
| 577 |
|
|
const double d = d_min + (d_max - d_min) * |
| 578 |
|
|
get_unit_range_from_texture_coord( |
| 579 |
|
|
1.0 - 2.0 * z, SCATTERING_TEXTURE_MU_SIZE / 2); |
| 580 |
|
|
*mu = d == 0.0 ? -1.0 : clamp_cosine(-(rho * rho + d * d) / (2.0 * *radius * d)); |
| 581 |
|
|
*r_mu_hits_ground = 1; |
| 582 |
|
|
} else { |
| 583 |
greg |
2.2 |
/* Distance to the top atmosphere boundary for the ray (r,mu), and its */ |
| 584 |
|
|
/* minimum and maximum values over all mu - obtained for (r,1) and */ |
| 585 |
|
|
/* (r,mu_horizon) - from which we can recover mu: */ |
| 586 |
greg |
2.1 |
const double d_min = AR - *radius; |
| 587 |
|
|
const double d_max = rho + H; |
| 588 |
|
|
const double d = d_min + (d_max - d_min) * |
| 589 |
|
|
get_unit_range_from_texture_coord( |
| 590 |
|
|
2.0 * z - 1.0, SCATTERING_TEXTURE_MU_SIZE / 2); |
| 591 |
|
|
*mu = d == 0.0 |
| 592 |
|
|
? 1.0 |
| 593 |
|
|
: clamp_cosine((H * H - rho * rho - d * d) / (2.0 * (*radius) * d)); |
| 594 |
|
|
*r_mu_hits_ground = 0; |
| 595 |
|
|
} |
| 596 |
|
|
|
| 597 |
|
|
const double x_mu_s = y; |
| 598 |
|
|
const double d_min = AH; |
| 599 |
|
|
const double d_max = H; |
| 600 |
|
|
const double D = distance_to_space(ER, MU_S_MIN); |
| 601 |
|
|
const double A = (D - d_min) / (d_max - d_min); |
| 602 |
|
|
const double a = (A - x_mu_s * A) / (1.0 + x_mu_s * A); |
| 603 |
|
|
const double d = d_min + fmin(a, A) * (d_max - d_min); |
| 604 |
|
|
*mu_s = d == 0.0 ? 1.0 : clamp_cosine((H * H - d * d) / (2.0 * ER * d)); |
| 605 |
|
|
*nu = clamp_cosine(x * 2.0 - 1.0); |
| 606 |
|
|
} |
| 607 |
|
|
|
| 608 |
|
|
static DATARRAY *interpolate_scattering(DATARRAY *dp, const double radius, |
| 609 |
|
|
const double mu, const double mu_s, |
| 610 |
|
|
const double nu, |
| 611 |
|
|
const int r_mu_hits_ground) { |
| 612 |
|
|
double u, v, w, z; |
| 613 |
|
|
to_scattering_uvwz(radius, mu, mu_s, nu, r_mu_hits_ground, &u, &v, &w, |
| 614 |
|
|
&z); |
| 615 |
|
|
double pt[4] = {z, w, v, u}; |
| 616 |
|
|
return datavector(dp, pt); |
| 617 |
|
|
} |
| 618 |
|
|
|
| 619 |
|
|
static void get_scattering(DATARRAY *scat1r, DATARRAY *scat1m, DATARRAY *scat, |
| 620 |
|
|
const double radius, const double mu, const double mu_s, |
| 621 |
|
|
const double nu, |
| 622 |
|
|
const int ray_r_mu_intersects_ground, |
| 623 |
|
|
const int scattering_order, double *result) { |
| 624 |
|
|
if (scattering_order == 1) { |
| 625 |
|
|
DATARRAY *rayleigh = interpolate_scattering(scat1r, radius, mu, mu_s, nu, |
| 626 |
|
|
ray_r_mu_intersects_ground); |
| 627 |
|
|
DATARRAY *mie = interpolate_scattering(scat1m, radius, mu, mu_s, nu, |
| 628 |
|
|
ray_r_mu_intersects_ground); |
| 629 |
|
|
const double rayleigh_phase = rayleigh_phase_function(nu); |
| 630 |
|
|
const double mie_phase = mie_phase_function(MIE_G, nu); |
| 631 |
|
|
for (int i = 0; i < NSSAMP; ++i) { |
| 632 |
|
|
result[i] = |
| 633 |
|
|
rayleigh->arr.d[i] * rayleigh_phase + mie->arr.d[i] * mie_phase; |
| 634 |
|
|
} |
| 635 |
|
|
free(rayleigh); |
| 636 |
|
|
free(mie); |
| 637 |
|
|
} else { |
| 638 |
|
|
DATARRAY *scattering = interpolate_scattering(scat, radius, mu, mu_s, nu, |
| 639 |
|
|
ray_r_mu_intersects_ground); |
| 640 |
|
|
for (int i = 0; i < NSSAMP; ++i) { |
| 641 |
|
|
result[i] = scattering->arr.d[i]; |
| 642 |
|
|
} |
| 643 |
|
|
free(scattering); |
| 644 |
|
|
} |
| 645 |
|
|
} |
| 646 |
|
|
|
| 647 |
|
|
static void to_irradiance_uv(const double radius, const double mu_s, double *u, |
| 648 |
|
|
double *v) { |
| 649 |
|
|
const double x_r = (radius - ER) / AH; |
| 650 |
|
|
const double x_mu_s = mu_s * 0.5 + 0.5; |
| 651 |
|
|
*u = x_mu_s; |
| 652 |
|
|
*v = x_r; |
| 653 |
|
|
} |
| 654 |
|
|
|
| 655 |
|
|
static void from_irradiance_uv(const double u, const double v, double *radius, |
| 656 |
|
|
double *mu_s) { |
| 657 |
|
|
*radius = v * AH + ER; |
| 658 |
|
|
*mu_s = clamp_cosine(u * 2.0 - 1.0); |
| 659 |
|
|
} |
| 660 |
|
|
|
| 661 |
|
|
DATARRAY *get_indirect_irradiance(DATARRAY *dp, const double radius, |
| 662 |
|
|
const double mu_s) { |
| 663 |
|
|
double u, v; |
| 664 |
|
|
to_irradiance_uv(radius, mu_s, &u, &v); |
| 665 |
|
|
double pt[2] = {v, u}; |
| 666 |
|
|
return datavector(dp, pt); |
| 667 |
|
|
} |
| 668 |
|
|
|
| 669 |
|
|
static void |
| 670 |
|
|
compute_scattering_density(const Atmosphere *atmos, DATARRAY *tau_dp, |
| 671 |
|
|
DATARRAY *scat1r, DATARRAY *scat1m, DATARRAY *scat, |
| 672 |
|
|
DATARRAY *irrad_dp, const double radius, const double mu, |
| 673 |
|
|
const double mu_s, const double nu, |
| 674 |
|
|
const int scattering_order, double *result) { |
| 675 |
greg |
2.2 |
/* Compute unit direction vectors for the zenith, the view direction omega and */ |
| 676 |
|
|
/* and the sun direction omega_s, such that the cosine of the view-zenith */ |
| 677 |
|
|
/* angle is mu, the cosine of the sun-zenith angle is mu_s, and the cosine of */ |
| 678 |
|
|
/* the view-sun angle is nu. The goal is to simplify computations below. */ |
| 679 |
greg |
2.1 |
const double height = radius - ER; |
| 680 |
|
|
const double epsilon = 1e-6; |
| 681 |
|
|
const FVECT zenith_direction = {0.0, 0.0, 1.0}; |
| 682 |
|
|
double omega_0 = 1.0 - mu * mu; |
| 683 |
|
|
omega_0 = fabs(omega_0) <= epsilon ? 0.0 : sqrt(omega_0); |
| 684 |
|
|
const FVECT omega = {omega_0, 0.0, mu}; |
| 685 |
|
|
const double sun_dir_x = omega[0] == 0.0 ? 0.0 : (nu - mu * mu_s) / omega[0]; |
| 686 |
|
|
const double sun_dir_y = |
| 687 |
|
|
sqrt(fmax(1.0 - sun_dir_x * sun_dir_x - mu_s * mu_s, 0.0)); |
| 688 |
|
|
const FVECT omega_s = {sun_dir_x, sun_dir_y, mu_s}; |
| 689 |
|
|
|
| 690 |
|
|
const int SAMPLE_COUNT = 16; |
| 691 |
|
|
const double dphi = PI / SAMPLE_COUNT; |
| 692 |
|
|
const double dtheta = PI / SAMPLE_COUNT; |
| 693 |
|
|
|
| 694 |
greg |
2.2 |
/* Nested loops for the integral over all the incident directions omega_i. */ |
| 695 |
greg |
2.1 |
for (int l = 0; l < SAMPLE_COUNT; ++l) { |
| 696 |
|
|
const double theta = (l + 0.5) * dtheta; |
| 697 |
|
|
const double cos_theta = cos(theta); |
| 698 |
|
|
const double sin_theta = sin(theta); |
| 699 |
|
|
const int r_theta_hits_ground = ray_hits_ground(radius, cos_theta); |
| 700 |
|
|
|
| 701 |
greg |
2.2 |
/* The distance and transmittance to the ground only depend on theta, so we */ |
| 702 |
|
|
/* can compute them in the outer loop for efficiency. */ |
| 703 |
greg |
2.1 |
double distance_to_ground = 0.0; |
| 704 |
|
|
double transmittance_to_ground[NSSAMP] = {0}; |
| 705 |
|
|
double ground_albedo = 0; |
| 706 |
|
|
if (r_theta_hits_ground) { |
| 707 |
|
|
distance_to_ground = distance_to_earth(radius, cos_theta); |
| 708 |
greg |
2.2 |
/* assert(distance_to_ground >= 0.0 && distance_to_ground <= HMAX); */ |
| 709 |
greg |
2.1 |
get_transmittance(tau_dp, radius, cos_theta, distance_to_ground, r_theta_hits_ground, |
| 710 |
|
|
transmittance_to_ground); |
| 711 |
|
|
ground_albedo = atmos->grefl; |
| 712 |
|
|
} |
| 713 |
|
|
|
| 714 |
|
|
for (int m = 0; m < 2 * SAMPLE_COUNT; ++m) { |
| 715 |
|
|
const double phi = (m + 0.5) * dphi; |
| 716 |
|
|
const FVECT omega_i = {cos(phi) * sin_theta, sin(phi) * sin_theta, |
| 717 |
|
|
cos_theta}; |
| 718 |
|
|
const double domega_i = dtheta * dphi * sin(theta); |
| 719 |
|
|
|
| 720 |
greg |
2.2 |
/* The radiance L_i arriving from direction omega_i after n-1 bounces is */ |
| 721 |
|
|
/* the sum of a term given by the precomputed scattering texture for the */ |
| 722 |
|
|
/* (n-1)-th order: */ |
| 723 |
greg |
2.1 |
double nu1 = fdot(omega_s, omega_i); |
| 724 |
|
|
if (nu1 <= -1.0) { |
| 725 |
|
|
nu1 += 0.1; |
| 726 |
|
|
} else if (nu1 >= 1.0) { |
| 727 |
|
|
nu1 -= 0.1; |
| 728 |
|
|
} |
| 729 |
|
|
double incident_radiance[NSSAMP] = {0}; |
| 730 |
|
|
get_scattering(scat1r, scat1m, scat, radius, omega_i[2], mu_s, nu1, |
| 731 |
|
|
r_theta_hits_ground, scattering_order - 1, |
| 732 |
|
|
incident_radiance); |
| 733 |
|
|
|
| 734 |
greg |
2.2 |
/* and of the contribution from the light paths with n-1 bounces and whose */ |
| 735 |
|
|
/* last bounce is on the ground. This contribution is the product of the */ |
| 736 |
|
|
/* transmittance to the ground, the ground albedo, the ground BRDF, and */ |
| 737 |
|
|
/* the irradiance received on the ground after n-2 bounces. */ |
| 738 |
greg |
2.1 |
FVECT ground_normal = { |
| 739 |
|
|
zenith_direction[0] * radius + omega_i[0] * distance_to_ground, |
| 740 |
|
|
zenith_direction[1] * radius + omega_i[1] * distance_to_ground, |
| 741 |
|
|
zenith_direction[2] * radius + omega_i[2] * distance_to_ground}; |
| 742 |
|
|
normalize(ground_normal); |
| 743 |
|
|
DATARRAY *ground_irradiance = |
| 744 |
|
|
get_indirect_irradiance(irrad_dp, ER, fdot(ground_normal, omega_s)); |
| 745 |
|
|
for (int i = 0; i < NSSAMP; ++i) { |
| 746 |
|
|
incident_radiance[i] += transmittance_to_ground[i] * ground_albedo * |
| 747 |
|
|
(1.0 / PI) * ground_irradiance->arr.d[i]; |
| 748 |
|
|
} |
| 749 |
|
|
free(ground_irradiance); |
| 750 |
|
|
|
| 751 |
greg |
2.2 |
/* The radiance finally scattered from direction omega_i towards direction */ |
| 752 |
|
|
/* -omega is the product of the incident radiance, the scattering */ |
| 753 |
|
|
/* coefficient, and the phase function for directions omega and omega_i */ |
| 754 |
|
|
/* (all this summed over all particle types, i.e. Rayleigh and Mie). */ |
| 755 |
greg |
2.1 |
const double nu2 = fdot(omega, omega_i); |
| 756 |
|
|
const double rayleigh_density = |
| 757 |
|
|
get_profile_density(&atmos->rayleigh_density, height); |
| 758 |
|
|
DATARRAY *mie_scat; |
| 759 |
|
|
double pt[1] = {height}; |
| 760 |
|
|
mie_scat = datavector(atmos->beta_m, pt); |
| 761 |
|
|
const double rayleigh_phase = rayleigh_phase_function(nu2); |
| 762 |
|
|
const double mie_phase = mie_phase_function(MIE_G, nu2); |
| 763 |
|
|
for (int j = 0; j < NSSAMP; ++j) { |
| 764 |
|
|
result[j] += incident_radiance[j] * domega_i * |
| 765 |
|
|
(atmos->beta_r0[j] * rayleigh_density * rayleigh_phase + |
| 766 |
|
|
mie_scat->arr.d[j] * mie_phase * atmos->beta_scale); |
| 767 |
|
|
} |
| 768 |
|
|
free(mie_scat); |
| 769 |
|
|
} |
| 770 |
|
|
} |
| 771 |
|
|
} |
| 772 |
|
|
|
| 773 |
|
|
static void compute_multi_scattering(DATARRAY *tau_dp, |
| 774 |
|
|
DATARRAY *scattering_density_dp, |
| 775 |
|
|
const double radius, const double mu, |
| 776 |
|
|
const double mu_s, const double nu, |
| 777 |
|
|
const int r_mu_hits_ground, |
| 778 |
|
|
double *result) { |
| 779 |
|
|
|
| 780 |
|
|
|
| 781 |
greg |
2.2 |
/* Numerical integration sample count. */ |
| 782 |
greg |
2.1 |
const int nsamp = 50; |
| 783 |
|
|
const double dx = distance_to_nearst_atmosphere_boundary( |
| 784 |
|
|
radius, mu, r_mu_hits_ground) / |
| 785 |
|
|
(double)nsamp; |
| 786 |
|
|
|
| 787 |
|
|
for (int i = 0; i <= nsamp; ++i) { |
| 788 |
|
|
const double d_i = i * dx; |
| 789 |
|
|
|
| 790 |
greg |
2.2 |
/* The r, mu and mu_s parameters at the current integration point (see the */ |
| 791 |
|
|
/* single scattering section for a detailed explanation). */ |
| 792 |
greg |
2.1 |
const double r_i = |
| 793 |
|
|
clamp_radius(sqrt(d_i * d_i + 2.0 * radius * mu * d_i + radius * radius)); |
| 794 |
|
|
const double mu_i = clamp_cosine((radius * mu + d_i) / r_i); |
| 795 |
|
|
const double mu_s_i = clamp_cosine((radius * mu_s + d_i * nu) / r_i); |
| 796 |
|
|
|
| 797 |
greg |
2.2 |
/* The Rayleigh and Mie multiple scattering at the current sample point. */ |
| 798 |
greg |
2.1 |
double transmittance[NSSAMP] = {0}; |
| 799 |
|
|
DATARRAY *rayleigh_mie_s = |
| 800 |
|
|
interpolate_scattering(scattering_density_dp, r_i, mu_i, mu_s_i, nu, |
| 801 |
|
|
r_mu_hits_ground); |
| 802 |
|
|
get_transmittance(tau_dp, radius, mu, d_i, r_mu_hits_ground, |
| 803 |
|
|
transmittance); |
| 804 |
greg |
2.2 |
/* Sample weight (from the trapezoidal rule). */ |
| 805 |
greg |
2.1 |
const double weight_i = (i == 0 || i == nsamp) ? 0.5 : 1.0; |
| 806 |
|
|
for (int j = 0; j < NSSAMP; ++j) { |
| 807 |
|
|
result[j] += rayleigh_mie_s->arr.d[j] * transmittance[j] * dx * weight_i; |
| 808 |
|
|
} |
| 809 |
|
|
free(rayleigh_mie_s); |
| 810 |
|
|
} |
| 811 |
|
|
} |
| 812 |
|
|
|
| 813 |
|
|
static void compute_direct_irradiance(DATARRAY *tau_dp, const double r, |
| 814 |
|
|
const double mu_s, float *result) { |
| 815 |
|
|
|
| 816 |
greg |
2.2 |
/* Approximate the average of the cosine factor mu_s over the visible fraction */ |
| 817 |
|
|
/* of the sun disc */ |
| 818 |
greg |
2.1 |
const double average_cosine_factor = |
| 819 |
|
|
mu_s < -SUNRAD |
| 820 |
|
|
? 0.0 |
| 821 |
|
|
: (mu_s > SUNRAD |
| 822 |
|
|
? mu_s |
| 823 |
|
|
: (mu_s + SUNRAD) * (mu_s + SUNRAD) / (4.0 * SUNRAD)); |
| 824 |
|
|
|
| 825 |
|
|
DATARRAY *transmittance = get_transmittance_to_space(tau_dp, r, mu_s); |
| 826 |
|
|
for (int i = 0; i < NSSAMP; ++i) { |
| 827 |
|
|
result[i] = EXTSOL[i] * transmittance->arr.d[i] * average_cosine_factor; |
| 828 |
|
|
} |
| 829 |
|
|
free(transmittance); |
| 830 |
|
|
} |
| 831 |
|
|
|
| 832 |
|
|
static void compute_indirect_irradiance(DATARRAY *scat1r, DATARRAY *scat1m, |
| 833 |
|
|
DATARRAY *scat, const double radius, |
| 834 |
|
|
const double mu_s, |
| 835 |
|
|
const int scattering_order, |
| 836 |
|
|
double *result) { |
| 837 |
|
|
|
| 838 |
|
|
const int SAMPLE_COUNT = 32; |
| 839 |
|
|
const double dphi = PI / SAMPLE_COUNT; |
| 840 |
|
|
const double dtheta = PI / SAMPLE_COUNT; |
| 841 |
|
|
|
| 842 |
|
|
const FVECT omega_s = {sqrt(1.0 - mu_s * mu_s), 0.0, mu_s}; |
| 843 |
|
|
for (int j = 0; j < SAMPLE_COUNT / 2; ++j) { |
| 844 |
|
|
const double theta = (j + 0.5) * dtheta; |
| 845 |
|
|
for (int i = 0; i < 2 * SAMPLE_COUNT; ++i) { |
| 846 |
|
|
const double phi = (i + 0.5) * dphi; |
| 847 |
|
|
const FVECT omega = {sin(theta) * cos(phi), sin(theta) * sin(phi), |
| 848 |
|
|
cos(theta)}; |
| 849 |
|
|
const double domega = dtheta * dphi * sin(theta); |
| 850 |
|
|
const double nu = fdot(omega, omega_s); |
| 851 |
|
|
double result1[NSSAMP] = {0}; |
| 852 |
|
|
get_scattering(scat1r, scat1m, scat, radius, omega[2], mu_s, nu, 0, |
| 853 |
|
|
scattering_order, result1); |
| 854 |
|
|
for (int k = 0; k < NSSAMP; ++k) { |
| 855 |
|
|
result[k] += result1[k] * omega[2] * domega; |
| 856 |
|
|
} |
| 857 |
|
|
} |
| 858 |
|
|
} |
| 859 |
|
|
} |
| 860 |
|
|
|
| 861 |
|
|
DATARRAY *allocate_3d_datarray(const char *name, const int ri, const int rj, |
| 862 |
|
|
const int rk) { |
| 863 |
|
|
DATARRAY *dp; |
| 864 |
|
|
if ((dp = (DATARRAY *)malloc(sizeof(DATARRAY))) == NULL) |
| 865 |
|
|
goto memerr; |
| 866 |
|
|
int asize = ri * rj * rk; |
| 867 |
|
|
dp->name = savestr(name); |
| 868 |
|
|
dp->type = DATATY; |
| 869 |
|
|
dp->nd = 3; |
| 870 |
|
|
dp->dim[0].org = 0; |
| 871 |
|
|
dp->dim[0].siz = 1; |
| 872 |
|
|
dp->dim[0].ne = ri; |
| 873 |
|
|
dp->dim[0].p = NULL; |
| 874 |
|
|
dp->dim[1].org = 0; |
| 875 |
|
|
dp->dim[1].siz = 1; |
| 876 |
|
|
dp->dim[1].ne = rj; |
| 877 |
|
|
dp->dim[1].p = NULL; |
| 878 |
|
|
dp->dim[2].org = 0; |
| 879 |
|
|
dp->dim[2].siz = rk - 1; |
| 880 |
|
|
dp->dim[2].ne = rk; |
| 881 |
|
|
dp->dim[2].p = NULL; |
| 882 |
greg |
2.2 |
/* dp->arr.p = */ |
| 883 |
|
|
/* malloc(sizeof(DATATYPE)*dp->dim[0].ne*dp->dim[1].ne*dp->dim[2].ne); */ |
| 884 |
greg |
2.1 |
if ((dp->arr.d = (DATATYPE *)malloc(asize * sizeof(DATATYPE))) == NULL) |
| 885 |
|
|
goto memerr; |
| 886 |
|
|
return (dp); |
| 887 |
|
|
|
| 888 |
|
|
memerr: |
| 889 |
|
|
fprintf(stderr, "Memory allocation error in allocate_3d_datarray\n"); |
| 890 |
|
|
return (NULL); |
| 891 |
|
|
} |
| 892 |
|
|
|
| 893 |
|
|
DATARRAY *allocate_5d_datarray(const char *name, const int ri, const int rj, |
| 894 |
|
|
const int rk, const int rl, const int rm) { |
| 895 |
|
|
DATARRAY *dp; |
| 896 |
|
|
int asize = ri * rj * rk * rl * rm; |
| 897 |
|
|
if ((dp = (DATARRAY *)malloc(sizeof(DATARRAY))) == NULL) |
| 898 |
|
|
goto memerr; |
| 899 |
|
|
dp->name = savestr(name); |
| 900 |
|
|
dp->type = DATATY; |
| 901 |
|
|
dp->nd = 5; |
| 902 |
|
|
dp->dim[0].org = 0; |
| 903 |
|
|
dp->dim[0].siz = 1; |
| 904 |
|
|
dp->dim[0].ne = ri; |
| 905 |
|
|
dp->dim[0].p = NULL; |
| 906 |
|
|
dp->dim[1].org = 0; |
| 907 |
|
|
dp->dim[1].siz = 1; |
| 908 |
|
|
dp->dim[1].ne = rj; |
| 909 |
|
|
dp->dim[1].p = NULL; |
| 910 |
|
|
dp->dim[2].org = 0; |
| 911 |
|
|
dp->dim[2].siz = 1; |
| 912 |
|
|
dp->dim[2].ne = rk; |
| 913 |
|
|
dp->dim[2].p = NULL; |
| 914 |
|
|
dp->dim[3].org = 0; |
| 915 |
|
|
dp->dim[3].siz = 1; |
| 916 |
|
|
dp->dim[3].ne = rl; |
| 917 |
|
|
dp->dim[3].p = NULL; |
| 918 |
|
|
dp->dim[4].org = 0; |
| 919 |
|
|
dp->dim[4].siz = rm - 1; |
| 920 |
|
|
dp->dim[4].ne = rm; |
| 921 |
|
|
dp->dim[4].p = NULL; |
| 922 |
|
|
if ((dp->arr.d = (DATATYPE *)malloc(asize * sizeof(DATATYPE))) == NULL) |
| 923 |
|
|
goto memerr; |
| 924 |
|
|
return (dp); |
| 925 |
|
|
memerr: |
| 926 |
|
|
fprintf(stderr, "Memory allocation error in allocate_5d_datarray\n"); |
| 927 |
|
|
return (NULL); |
| 928 |
|
|
} |
| 929 |
|
|
|
| 930 |
|
|
/* save to a .dat file */ |
| 931 |
|
|
void savedata(DATARRAY *dp) { |
| 932 |
|
|
if (dp == NULL) |
| 933 |
|
|
return; |
| 934 |
|
|
FILE *fp = fopen(dp->name, "w"); |
| 935 |
|
|
if (fp == NULL) { |
| 936 |
|
|
fprintf(stderr, "Error opening file %s\n", dp->name); |
| 937 |
|
|
return; |
| 938 |
|
|
} |
| 939 |
|
|
int i, j; |
| 940 |
|
|
int nvals = 1; |
| 941 |
|
|
for (i = 0; i < dp->nd; i++) { |
| 942 |
|
|
nvals *= dp->dim[i].ne; |
| 943 |
|
|
} |
| 944 |
|
|
fprintf(fp, "%d\n", dp->nd); |
| 945 |
|
|
for (i = 0; i < dp->nd; i++) { |
| 946 |
|
|
fprintf(fp, "%f %f %d\n", dp->dim[i].org, dp->dim[i].siz, dp->dim[i].ne); |
| 947 |
|
|
} |
| 948 |
|
|
for (i = 0; i < nvals; i++) { |
| 949 |
greg |
2.3 |
fprintf(fp, "%.15e\n", dp->arr.d[i]); |
| 950 |
greg |
2.1 |
} |
| 951 |
|
|
fclose(fp); |
| 952 |
|
|
} |
| 953 |
|
|
|
| 954 |
|
|
void increment_dp(DATARRAY *dp1, DATARRAY *dp2) { |
| 955 |
|
|
if (dp1 == NULL || dp2 == NULL) |
| 956 |
|
|
perror("null pointer in increment_dp\n"); |
| 957 |
|
|
|
| 958 |
|
|
if (dp1->nd != dp2->nd) |
| 959 |
|
|
perror("dimension mismatch in increment_dp\n"); |
| 960 |
|
|
|
| 961 |
|
|
int i; |
| 962 |
|
|
int nvals1 = 1; |
| 963 |
|
|
int nvals2 = 1; |
| 964 |
|
|
for (i = 0; i < dp1->nd; i++) { |
| 965 |
|
|
nvals1 *= dp1->dim[i].ne; |
| 966 |
|
|
nvals2 *= dp2->dim[i].ne; |
| 967 |
|
|
} |
| 968 |
|
|
if (nvals1 != nvals2) |
| 969 |
|
|
perror("size mismatch in increment_dp\n"); |
| 970 |
|
|
for (i = 0; i < nvals1; i++) { |
| 971 |
|
|
dp1->arr.d[i] += dp2->arr.d[i]; |
| 972 |
|
|
} |
| 973 |
|
|
} |
| 974 |
|
|
|
| 975 |
|
|
THREAD_RETURN compute_scattering_density_thread(void *arg) { |
| 976 |
|
|
ScatDenseTdat *tdata = (ScatDenseTdat *)arg; |
| 977 |
|
|
for (unsigned int l = tdata->start_l; l < tdata->end_l; ++l) { |
| 978 |
|
|
for (unsigned int k = 0; k < SCATTERING_TEXTURE_MU_SIZE; ++k) { |
| 979 |
|
|
for (unsigned int j = 0; j < SCATTERING_TEXTURE_MU_S_SIZE; ++j) { |
| 980 |
|
|
for (unsigned int i = 0; i < SCATTERING_TEXTURE_NU_SIZE; ++i) { |
| 981 |
|
|
double scattering_density[NSSAMP] = {0}; |
| 982 |
|
|
double r, mu, mu_s, nu; |
| 983 |
|
|
int ray_r_mu_intersects_ground; |
| 984 |
|
|
double xr = (double)i / SCATTERING_TEXTURE_NU_SIZE; |
| 985 |
|
|
double yr = (double)j / SCATTERING_TEXTURE_MU_S_SIZE; |
| 986 |
|
|
double zr = (double)k / SCATTERING_TEXTURE_MU_SIZE; |
| 987 |
|
|
double wr = (double)l / SCATTERING_TEXTURE_R_SIZE; |
| 988 |
|
|
from_scattering_uvwz(xr, yr, zr, wr, &r, &mu, &mu_s, &nu, |
| 989 |
|
|
&ray_r_mu_intersects_ground); |
| 990 |
|
|
nu = |
| 991 |
|
|
fmax(mu * mu_s - sqrt((1.0 - mu * mu) * (1.0 - mu_s * mu_s)), |
| 992 |
|
|
fmin(mu * mu_s + sqrt((1.0 - mu * mu) * (1.0 - mu_s * mu_s)), |
| 993 |
|
|
nu)); |
| 994 |
|
|
compute_scattering_density( |
| 995 |
|
|
tdata->atmos, tdata->tau_dp, tdata->delta_rayleigh_scattering_dp, |
| 996 |
|
|
tdata->delta_mie_scattering_dp, |
| 997 |
|
|
tdata->delta_multiple_scattering_dp, tdata->delta_irradiance_dp, |
| 998 |
|
|
r, mu, mu_s, nu, tdata->scattering_order, scattering_density); |
| 999 |
|
|
for (int m = 0; m < NSSAMP; ++m) { |
| 1000 |
|
|
int idx = l * SCATTERING_TEXTURE_MU_SIZE * |
| 1001 |
|
|
SCATTERING_TEXTURE_MU_S_SIZE * |
| 1002 |
|
|
SCATTERING_TEXTURE_NU_SIZE * NSSAMP + |
| 1003 |
|
|
k * SCATTERING_TEXTURE_MU_S_SIZE * |
| 1004 |
|
|
SCATTERING_TEXTURE_NU_SIZE * NSSAMP + |
| 1005 |
|
|
j * SCATTERING_TEXTURE_NU_SIZE * NSSAMP + i * NSSAMP + m; |
| 1006 |
|
|
tdata->delta_scattering_density_dp->arr.d[idx] = |
| 1007 |
|
|
scattering_density[m]; |
| 1008 |
|
|
} |
| 1009 |
|
|
} |
| 1010 |
|
|
} |
| 1011 |
|
|
} |
| 1012 |
|
|
} |
| 1013 |
|
|
return 0; |
| 1014 |
|
|
} |
| 1015 |
|
|
|
| 1016 |
|
|
THREAD_RETURN compute_single_scattering_thread(void *arg) { |
| 1017 |
|
|
Scat1Tdat *tdata = (Scat1Tdat *)arg; |
| 1018 |
|
|
for (unsigned int l = tdata->start_l; l < tdata->end_l; ++l) { |
| 1019 |
|
|
for (int k = 0; k < SCATTERING_TEXTURE_MU_SIZE; ++k) { |
| 1020 |
|
|
for (int j = 0; j < SCATTERING_TEXTURE_MU_S_SIZE; ++j) { |
| 1021 |
|
|
for (int i = 0; i < SCATTERING_TEXTURE_NU_SIZE; ++i) { |
| 1022 |
|
|
float rayleigh[NSSAMP] = {0}; |
| 1023 |
|
|
float mie[NSSAMP] = {0}; |
| 1024 |
|
|
double r, mu, mu_s, nu; |
| 1025 |
|
|
int r_mu_hits_ground; |
| 1026 |
|
|
const double xr = (double)i / SCATTERING_TEXTURE_NU_SIZE; |
| 1027 |
|
|
const double yr = (double)j / SCATTERING_TEXTURE_MU_S_SIZE; |
| 1028 |
|
|
const double zr = (double)k / SCATTERING_TEXTURE_MU_SIZE; |
| 1029 |
|
|
const double wr = (double)l / SCATTERING_TEXTURE_R_SIZE; |
| 1030 |
|
|
from_scattering_uvwz(xr, yr, zr, wr, &r, &mu, &mu_s, &nu, |
| 1031 |
|
|
&r_mu_hits_ground); |
| 1032 |
|
|
const double mu_mu_s = mu * mu_s; |
| 1033 |
|
|
const double sqrt_term = sqrt((1.0 - mu * mu) * (1.0 - mu_s * mu_s)); |
| 1034 |
|
|
nu = fmax(mu_mu_s - sqrt_term, fmin(mu_mu_s + sqrt_term, nu)); |
| 1035 |
|
|
compute_single_scattering(tdata->atmos, tdata->tau_dp, r, mu, mu_s, |
| 1036 |
|
|
nu, r_mu_hits_ground, rayleigh, |
| 1037 |
|
|
mie); |
| 1038 |
|
|
for (int m = 0; m < NSSAMP; ++m) { |
| 1039 |
|
|
int idx = l * SCATTERING_TEXTURE_MU_SIZE * |
| 1040 |
|
|
SCATTERING_TEXTURE_MU_S_SIZE * |
| 1041 |
|
|
SCATTERING_TEXTURE_NU_SIZE * NSSAMP + |
| 1042 |
|
|
k * SCATTERING_TEXTURE_MU_S_SIZE * |
| 1043 |
|
|
SCATTERING_TEXTURE_NU_SIZE * NSSAMP + |
| 1044 |
|
|
j * SCATTERING_TEXTURE_NU_SIZE * NSSAMP + i * NSSAMP + m; |
| 1045 |
|
|
tdata->delta_rayleigh_scattering_dp->arr.d[idx] = rayleigh[m]; |
| 1046 |
|
|
tdata->delta_mie_scattering_dp->arr.d[idx] = mie[m]; |
| 1047 |
|
|
tdata->scattering_dp->arr.d[idx] = rayleigh[m]; |
| 1048 |
|
|
} |
| 1049 |
|
|
} |
| 1050 |
|
|
} |
| 1051 |
|
|
} |
| 1052 |
|
|
} |
| 1053 |
|
|
return 0; |
| 1054 |
|
|
} |
| 1055 |
|
|
|
| 1056 |
|
|
THREAD_RETURN compute_multi_scattering_thread(void *arg) { |
| 1057 |
|
|
ScatNTdat *tdata = (ScatNTdat *)arg; |
| 1058 |
|
|
for (unsigned int l = tdata->start_l; l < tdata->end_l; ++l) { |
| 1059 |
|
|
for (unsigned int k = 0; k < SCATTERING_TEXTURE_MU_SIZE; ++k) { |
| 1060 |
|
|
for (unsigned int j = 0; j < SCATTERING_TEXTURE_MU_S_SIZE; ++j) { |
| 1061 |
|
|
for (unsigned int i = 0; i < SCATTERING_TEXTURE_NU_SIZE; ++i) { |
| 1062 |
|
|
double delta_multiple_scattering[NSSAMP] = {0}; |
| 1063 |
|
|
double r, mu, mu_s, nu; |
| 1064 |
|
|
int r_mu_hits_ground; |
| 1065 |
|
|
const double xr = (double)i / SCATTERING_TEXTURE_NU_SIZE; |
| 1066 |
|
|
const double yr = (double)j / SCATTERING_TEXTURE_MU_S_SIZE; |
| 1067 |
|
|
const double zr = (double)k / SCATTERING_TEXTURE_MU_SIZE; |
| 1068 |
|
|
const double wr = (double)l / SCATTERING_TEXTURE_R_SIZE; |
| 1069 |
|
|
from_scattering_uvwz(xr, yr, zr, wr, &r, &mu, &mu_s, &nu, |
| 1070 |
|
|
&r_mu_hits_ground); |
| 1071 |
|
|
const double mu_mu_s = mu * mu_s; |
| 1072 |
|
|
const double sqrt_term = sqrt((1.0 - mu * mu) * (1.0 - mu_s * mu_s)); |
| 1073 |
|
|
nu = fmax(mu_mu_s - sqrt_term, fmin(mu_mu_s + sqrt_term, nu)); |
| 1074 |
|
|
compute_multi_scattering( |
| 1075 |
|
|
tdata->tau_dp, tdata->delta_scattering_density_dp, r, mu, mu_s, |
| 1076 |
|
|
nu, r_mu_hits_ground, delta_multiple_scattering); |
| 1077 |
|
|
const double rayleigh_phase = rayleigh_phase_function(nu); |
| 1078 |
|
|
for (int m = 0; m < NSSAMP; ++m) { |
| 1079 |
|
|
int idx = l * SCATTERING_TEXTURE_MU_SIZE * |
| 1080 |
|
|
SCATTERING_TEXTURE_MU_S_SIZE * |
| 1081 |
|
|
SCATTERING_TEXTURE_NU_SIZE * NSSAMP + |
| 1082 |
|
|
k * SCATTERING_TEXTURE_MU_S_SIZE * |
| 1083 |
|
|
SCATTERING_TEXTURE_NU_SIZE * NSSAMP + |
| 1084 |
|
|
j * SCATTERING_TEXTURE_NU_SIZE * NSSAMP + i * NSSAMP + m; |
| 1085 |
|
|
tdata->delta_multiple_scattering_dp->arr.d[idx] = |
| 1086 |
|
|
delta_multiple_scattering[m]; |
| 1087 |
|
|
tdata->scattering_dp->arr.d[idx] += |
| 1088 |
|
|
delta_multiple_scattering[m] * (1.0 / rayleigh_phase); |
| 1089 |
|
|
} |
| 1090 |
|
|
} |
| 1091 |
|
|
} |
| 1092 |
|
|
} |
| 1093 |
|
|
} |
| 1094 |
|
|
return 0; |
| 1095 |
|
|
} |
| 1096 |
|
|
|
| 1097 |
|
|
int precompute(const int sorder, const DpPaths dppaths, const Atmosphere *atmos, |
| 1098 |
|
|
int num_threads) { |
| 1099 |
|
|
num_threads = (num_threads < SCATTERING_TEXTURE_R_SIZE) |
| 1100 |
|
|
? num_threads |
| 1101 |
|
|
: SCATTERING_TEXTURE_R_SIZE; |
| 1102 |
|
|
const int tchunk = SCATTERING_TEXTURE_R_SIZE / num_threads; |
| 1103 |
|
|
const int tremainder = SCATTERING_TEXTURE_R_SIZE % num_threads; |
| 1104 |
|
|
|
| 1105 |
|
|
if (sorder < 2) { |
| 1106 |
|
|
fprintf(stderr, "scattering order must be at least 2\n"); |
| 1107 |
|
|
return 0; |
| 1108 |
|
|
} |
| 1109 |
|
|
|
| 1110 |
|
|
DATARRAY *tau_dp = |
| 1111 |
|
|
allocate_3d_datarray(dppaths.tau, TRANSMITTANCE_TEXTURE_HEIGHT, |
| 1112 |
|
|
TRANSMITTANCE_TEXTURE_WIDTH, NSSAMP); |
| 1113 |
|
|
DATARRAY *delta_irradiance_dp = allocate_3d_datarray( |
| 1114 |
|
|
"ssky_delta_irradiance.dat", IRRADIANCE_TEXTURE_HEIGHT, |
| 1115 |
|
|
IRRADIANCE_TEXTURE_WIDTH, NSSAMP); |
| 1116 |
|
|
DATARRAY *irradiance_dp = |
| 1117 |
|
|
allocate_3d_datarray(dppaths.irrad, IRRADIANCE_TEXTURE_HEIGHT, |
| 1118 |
|
|
IRRADIANCE_TEXTURE_WIDTH, NSSAMP); |
| 1119 |
|
|
|
| 1120 |
|
|
DATARRAY *delta_rayleigh_scattering_dp = allocate_5d_datarray( |
| 1121 |
|
|
"ssky_delta_rayleigh_scattering_dp.dat", SCATTERING_TEXTURE_R_SIZE, |
| 1122 |
|
|
SCATTERING_TEXTURE_MU_SIZE, SCATTERING_TEXTURE_MU_S_SIZE, |
| 1123 |
|
|
SCATTERING_TEXTURE_NU_SIZE, NSSAMP); |
| 1124 |
|
|
DATARRAY *delta_mie_scattering_dp = allocate_5d_datarray( |
| 1125 |
|
|
dppaths.scat1m, SCATTERING_TEXTURE_R_SIZE, SCATTERING_TEXTURE_MU_SIZE, |
| 1126 |
|
|
SCATTERING_TEXTURE_MU_S_SIZE, SCATTERING_TEXTURE_NU_SIZE, NSSAMP); |
| 1127 |
|
|
DATARRAY *scattering_dp = allocate_5d_datarray( |
| 1128 |
|
|
dppaths.scat, SCATTERING_TEXTURE_R_SIZE, SCATTERING_TEXTURE_MU_SIZE, |
| 1129 |
|
|
SCATTERING_TEXTURE_MU_S_SIZE, SCATTERING_TEXTURE_NU_SIZE, NSSAMP); |
| 1130 |
|
|
DATARRAY *delta_multiple_scattering_dp = allocate_5d_datarray( |
| 1131 |
|
|
"ssky_delta_multiple_scattering.dat", SCATTERING_TEXTURE_R_SIZE, |
| 1132 |
|
|
SCATTERING_TEXTURE_MU_SIZE, SCATTERING_TEXTURE_MU_S_SIZE, |
| 1133 |
|
|
SCATTERING_TEXTURE_NU_SIZE, NSSAMP); |
| 1134 |
|
|
DATARRAY *delta_scattering_density_dp = allocate_5d_datarray( |
| 1135 |
|
|
"ssky_delta_scattering_density.dat", SCATTERING_TEXTURE_R_SIZE, |
| 1136 |
|
|
SCATTERING_TEXTURE_MU_SIZE, SCATTERING_TEXTURE_MU_S_SIZE, |
| 1137 |
|
|
SCATTERING_TEXTURE_NU_SIZE, NSSAMP); |
| 1138 |
|
|
|
| 1139 |
|
|
int idx = 0; |
| 1140 |
greg |
2.2 |
/* Computing transmittance... */ |
| 1141 |
greg |
2.1 |
for (int j = 0; j < TRANSMITTANCE_TEXTURE_HEIGHT; ++j) { |
| 1142 |
|
|
for (int i = 0; i < TRANSMITTANCE_TEXTURE_WIDTH; ++i) { |
| 1143 |
|
|
double r, mu; |
| 1144 |
|
|
float tau[NSSAMP] = {0}; |
| 1145 |
|
|
const float xr = (float)i / TRANSMITTANCE_TEXTURE_WIDTH; |
| 1146 |
|
|
const float yr = (float)j / TRANSMITTANCE_TEXTURE_HEIGHT; |
| 1147 |
|
|
from_transmittance_uv(xr, yr, &r, &mu); |
| 1148 |
|
|
compute_transmittance_to_space(atmos, r, mu, tau); |
| 1149 |
|
|
for (int k = 0; k < NSSAMP; ++k) { |
| 1150 |
|
|
tau_dp->arr.d[idx] = tau[k]; |
| 1151 |
|
|
idx += 1; |
| 1152 |
|
|
} |
| 1153 |
|
|
} |
| 1154 |
|
|
} |
| 1155 |
|
|
savedata(tau_dp); |
| 1156 |
|
|
|
| 1157 |
greg |
2.2 |
/* Computing direct irradiance... */ |
| 1158 |
greg |
2.1 |
idx = 0; |
| 1159 |
|
|
for (int j = 0; j < IRRADIANCE_TEXTURE_HEIGHT; ++j) { |
| 1160 |
|
|
for (int i = 0; i < IRRADIANCE_TEXTURE_WIDTH; ++i) { |
| 1161 |
|
|
double r, mu_s; |
| 1162 |
|
|
float result[NSSAMP] = {0}; |
| 1163 |
|
|
const double xr = (double)i / IRRADIANCE_TEXTURE_WIDTH; |
| 1164 |
|
|
const double yr = (double)j / IRRADIANCE_TEXTURE_HEIGHT; |
| 1165 |
|
|
from_irradiance_uv(xr, yr, &r, &mu_s); |
| 1166 |
|
|
compute_direct_irradiance(tau_dp, r, mu_s, result); |
| 1167 |
|
|
for (int k = 0; k < NSSAMP; ++k) { |
| 1168 |
|
|
delta_irradiance_dp->arr.d[idx] = result[k]; |
| 1169 |
|
|
irradiance_dp->arr.d[idx] = 0.0; |
| 1170 |
|
|
idx += 1; |
| 1171 |
|
|
} |
| 1172 |
|
|
} |
| 1173 |
|
|
} |
| 1174 |
|
|
|
| 1175 |
greg |
2.2 |
/* Computing single scattering... */ |
| 1176 |
greg |
2.1 |
THREAD *threads = (THREAD *)malloc(num_threads * sizeof(THREAD)); |
| 1177 |
|
|
Scat1Tdat *tdata = (Scat1Tdat *)malloc(num_threads * sizeof(Scat1Tdat)); |
| 1178 |
|
|
for (int i = 0; i < num_threads; i++) { |
| 1179 |
|
|
tdata[i].atmos = atmos; |
| 1180 |
|
|
tdata[i].tau_dp = tau_dp; |
| 1181 |
|
|
tdata[i].delta_rayleigh_scattering_dp = delta_rayleigh_scattering_dp; |
| 1182 |
|
|
tdata[i].delta_mie_scattering_dp = delta_mie_scattering_dp; |
| 1183 |
|
|
tdata[i].scattering_dp = scattering_dp; |
| 1184 |
|
|
tdata[i].start_l = i * tchunk; |
| 1185 |
|
|
tdata[i].end_l = (i + 1) * tchunk; |
| 1186 |
|
|
if (i == num_threads - 1) |
| 1187 |
|
|
tdata[i].end_l += tremainder; |
| 1188 |
|
|
CREATE_THREAD(&threads[i], compute_single_scattering_thread, |
| 1189 |
|
|
(void *)&tdata[i]); |
| 1190 |
|
|
} |
| 1191 |
|
|
|
| 1192 |
|
|
for (int i = 0; i < num_threads; i++) { |
| 1193 |
|
|
THREAD_JOIN(threads[i]); |
| 1194 |
|
|
#if defined(_WIN32) || defined(_WIN64) |
| 1195 |
|
|
THREAD_CLOSE(threads[i]); |
| 1196 |
|
|
#endif |
| 1197 |
|
|
} |
| 1198 |
|
|
free(tdata); |
| 1199 |
|
|
free(threads); |
| 1200 |
|
|
savedata(delta_mie_scattering_dp); |
| 1201 |
|
|
|
| 1202 |
|
|
|
| 1203 |
|
|
for (int scattering_order = 2; scattering_order <= sorder; ++scattering_order) { |
| 1204 |
|
|
THREAD *threads = (THREAD *)malloc(num_threads * sizeof(THREAD)); |
| 1205 |
|
|
ScatDenseTdat *tdata = (ScatDenseTdat *)malloc(num_threads * sizeof(ScatDenseTdat)); |
| 1206 |
|
|
for (int i = 0; i < num_threads; i++) { |
| 1207 |
|
|
tdata[i].atmos = atmos; |
| 1208 |
|
|
tdata[i].tau_dp = tau_dp; |
| 1209 |
|
|
tdata[i].delta_rayleigh_scattering_dp = delta_rayleigh_scattering_dp; |
| 1210 |
|
|
tdata[i].delta_mie_scattering_dp = delta_mie_scattering_dp; |
| 1211 |
|
|
tdata[i].delta_multiple_scattering_dp = delta_multiple_scattering_dp; |
| 1212 |
|
|
tdata[i].delta_irradiance_dp = delta_irradiance_dp; |
| 1213 |
|
|
tdata[i].scattering_order = scattering_order; |
| 1214 |
|
|
tdata[i].delta_scattering_density_dp = delta_scattering_density_dp; |
| 1215 |
|
|
tdata[i].start_l = i * tchunk; |
| 1216 |
|
|
tdata[i].end_l = (i + 1) * tchunk; |
| 1217 |
|
|
if (i == num_threads - 1) |
| 1218 |
|
|
tdata[i].end_l += tremainder; |
| 1219 |
|
|
CREATE_THREAD(&threads[i], compute_scattering_density_thread, |
| 1220 |
|
|
(void *)&tdata[i]); |
| 1221 |
|
|
} |
| 1222 |
|
|
|
| 1223 |
|
|
for (int i = 0; i < num_threads; i++) { |
| 1224 |
|
|
THREAD_JOIN(threads[i]); |
| 1225 |
|
|
#if defined(_WIN32) || defined(_WIN64) |
| 1226 |
|
|
THREAD_CLOSE(threads[i]); |
| 1227 |
|
|
#endif |
| 1228 |
|
|
} |
| 1229 |
|
|
free(tdata); |
| 1230 |
|
|
free(threads); |
| 1231 |
|
|
|
| 1232 |
|
|
idx = 0; |
| 1233 |
|
|
for (unsigned int j = 0; j < IRRADIANCE_TEXTURE_HEIGHT; ++j) { |
| 1234 |
|
|
for (unsigned int i = 0; i < IRRADIANCE_TEXTURE_WIDTH; ++i) { |
| 1235 |
|
|
double delta_irradiance[NSSAMP] = {0}; |
| 1236 |
|
|
double r, mu_s; |
| 1237 |
|
|
const double xr = (double)i / IRRADIANCE_TEXTURE_WIDTH; |
| 1238 |
|
|
const double yr = (double)j / IRRADIANCE_TEXTURE_HEIGHT; |
| 1239 |
|
|
from_irradiance_uv(xr, yr, &r, &mu_s); |
| 1240 |
|
|
compute_indirect_irradiance(delta_rayleigh_scattering_dp, |
| 1241 |
|
|
delta_mie_scattering_dp, |
| 1242 |
|
|
delta_multiple_scattering_dp, r, mu_s, |
| 1243 |
|
|
scattering_order - 1, delta_irradiance); |
| 1244 |
|
|
for (int k = 0; k < NSSAMP; ++k) { |
| 1245 |
|
|
delta_irradiance_dp->arr.d[idx] = delta_irradiance[k]; |
| 1246 |
|
|
idx += 1; |
| 1247 |
|
|
} |
| 1248 |
|
|
} |
| 1249 |
|
|
} |
| 1250 |
|
|
|
| 1251 |
|
|
increment_dp(irradiance_dp, delta_irradiance_dp); |
| 1252 |
|
|
|
| 1253 |
|
|
THREAD *threads2 = (THREAD *)malloc(num_threads * sizeof(THREAD)); |
| 1254 |
|
|
ScatNTdat *tdata2 = (ScatNTdat *)malloc(num_threads * sizeof(ScatNTdat)); |
| 1255 |
|
|
for (int i = 0; i < num_threads; i++) { |
| 1256 |
|
|
tdata2[i].tau_dp = tau_dp; |
| 1257 |
|
|
tdata2[i].delta_multiple_scattering_dp = delta_multiple_scattering_dp; |
| 1258 |
|
|
tdata2[i].delta_scattering_density_dp = delta_scattering_density_dp; |
| 1259 |
|
|
tdata2[i].scattering_dp = scattering_dp; |
| 1260 |
|
|
tdata2[i].start_l = i * tchunk; |
| 1261 |
|
|
tdata2[i].end_l = (i + 1) * tchunk; |
| 1262 |
|
|
if (i == num_threads - 1) |
| 1263 |
|
|
tdata2[i].end_l += tremainder; |
| 1264 |
|
|
CREATE_THREAD(&threads2[i], compute_multi_scattering_thread, |
| 1265 |
|
|
(void *)&tdata2[i]); |
| 1266 |
|
|
} |
| 1267 |
|
|
|
| 1268 |
|
|
for (int i = 0; i < num_threads; i++) { |
| 1269 |
|
|
THREAD_JOIN(threads2[i]); |
| 1270 |
|
|
#if defined(_WIN32) || defined(_WIN64) |
| 1271 |
|
|
THREAD_CLOSE(threads2[i]); |
| 1272 |
|
|
#endif |
| 1273 |
|
|
} |
| 1274 |
|
|
free(tdata2); |
| 1275 |
|
|
free(threads2); |
| 1276 |
|
|
} |
| 1277 |
|
|
savedata(scattering_dp); |
| 1278 |
|
|
savedata(irradiance_dp); |
| 1279 |
|
|
|
| 1280 |
|
|
freedata(tau_dp); |
| 1281 |
|
|
freedata(delta_rayleigh_scattering_dp); |
| 1282 |
|
|
freedata(delta_mie_scattering_dp); |
| 1283 |
|
|
freedata(scattering_dp); |
| 1284 |
|
|
freedata(delta_multiple_scattering_dp); |
| 1285 |
|
|
freedata(delta_scattering_density_dp); |
| 1286 |
|
|
freedata(delta_irradiance_dp); |
| 1287 |
|
|
freedata(irradiance_dp); |
| 1288 |
|
|
return 1; |
| 1289 |
|
|
} |
| 1290 |
|
|
|
| 1291 |
|
|
int compute_sundir(const int year, const int month, const int day, |
| 1292 |
|
|
const double hour, const int tsolar, double sundir[3]) { |
| 1293 |
|
|
double sd, st = hour; |
| 1294 |
|
|
|
| 1295 |
|
|
if (year) { |
| 1296 |
|
|
/* Michalsky algorithm? */ |
| 1297 |
|
|
double mjd = mjdate(year, month, day, hour); |
| 1298 |
|
|
if (tsolar) |
| 1299 |
|
|
sd = msdec(mjd, NULL); |
| 1300 |
|
|
else |
| 1301 |
|
|
sd = msdec(mjd, &st); |
| 1302 |
|
|
} else { |
| 1303 |
|
|
int jd = jdate(month, day); /* Julian date */ |
| 1304 |
|
|
sd = sdec(jd); /* solar declination */ |
| 1305 |
|
|
if (!tsolar) /* get solar time? */ |
| 1306 |
|
|
st = hour + stadj(jd); |
| 1307 |
|
|
} |
| 1308 |
|
|
const double altitude = salt(sd, st); |
| 1309 |
|
|
const double azimuth = sazi(sd, st); |
| 1310 |
|
|
sundir[0] = -sin(azimuth) * cos(altitude); |
| 1311 |
|
|
sundir[1] = -cos(azimuth) * cos(altitude); |
| 1312 |
|
|
sundir[2] = sin(altitude); |
| 1313 |
|
|
return 1; |
| 1314 |
|
|
} |
| 1315 |
|
|
|
| 1316 |
|
|
void get_sky_transmittance(DATARRAY *tau, const double radius, const double mu, float *result) { |
| 1317 |
|
|
DATARRAY *trans = get_transmittance_to_space(tau, radius, mu); |
| 1318 |
|
|
for (int i = 0; i < NSSAMP; ++i) { |
| 1319 |
|
|
result[i] = trans->arr.d[i]; |
| 1320 |
|
|
} |
| 1321 |
|
|
free(trans); |
| 1322 |
|
|
} |
| 1323 |
|
|
|
| 1324 |
|
|
void get_sky_radiance(DATARRAY *scat, DATARRAY *scat1m, const double radius, |
| 1325 |
|
|
const double mu, const double mu_s, const double nu, float *result) { |
| 1326 |
|
|
double u, v, w, z; |
| 1327 |
|
|
to_scattering_uvwz(radius, mu, mu_s, nu, 0, &u, &v, &w, &z); |
| 1328 |
|
|
double pt[4] = {z, w, v, u}; |
| 1329 |
|
|
|
| 1330 |
|
|
DATARRAY *scattering = datavector(scat, pt); |
| 1331 |
|
|
DATARRAY *single_mie_scattering = datavector(scat1m, pt); |
| 1332 |
|
|
|
| 1333 |
|
|
const double rayleigh_phase = rayleigh_phase_function(nu); |
| 1334 |
|
|
const double mie_phase = mie_phase_function(MIE_G, nu); |
| 1335 |
|
|
|
| 1336 |
|
|
for (int i = 0; i < NSSAMP; ++i) { |
| 1337 |
|
|
result[i] = scattering->arr.d[i] * rayleigh_phase + |
| 1338 |
|
|
single_mie_scattering->arr.d[i] * mie_phase; |
| 1339 |
|
|
} |
| 1340 |
|
|
free(single_mie_scattering); |
| 1341 |
|
|
free(scattering); |
| 1342 |
|
|
} |
| 1343 |
|
|
|
| 1344 |
|
|
static void get_irradiance(DATARRAY *tau, DATARRAY *irrad, const double radius, |
| 1345 |
|
|
const FVECT point, const FVECT normal, |
| 1346 |
|
|
const FVECT sundir, double *result) { |
| 1347 |
|
|
double mu_s = fdot(point, sundir) / radius; |
| 1348 |
|
|
|
| 1349 |
|
|
double point_trans_sun[NSSAMP] = {0}; |
| 1350 |
|
|
get_transmittance_to_sun(tau, radius, mu_s, point_trans_sun); |
| 1351 |
|
|
|
| 1352 |
|
|
DATARRAY *indirect_irradiance = get_indirect_irradiance(irrad, radius, mu_s); |
| 1353 |
|
|
const double sun_ct = fmax(fdot(normal, sundir), 0.0); |
| 1354 |
|
|
|
| 1355 |
|
|
for (int i = 0; i < NSSAMP; ++i) { |
| 1356 |
|
|
result[i] = indirect_irradiance->arr.d[i] + point_trans_sun[i] * EXTSOL[i] * sun_ct; |
| 1357 |
|
|
} |
| 1358 |
|
|
|
| 1359 |
|
|
free(indirect_irradiance); |
| 1360 |
|
|
} |
| 1361 |
|
|
|
| 1362 |
|
|
|
| 1363 |
|
|
void get_ground_radiance(DATARRAY *tau, DATARRAY *scat, DATARRAY *scat1m, DATARRAY *irrad, |
| 1364 |
|
|
const FVECT view_point, const FVECT view_direction, |
| 1365 |
|
|
const double radius, const double mu, const double sun_ct, |
| 1366 |
|
|
const double nu, const double grefl, const FVECT sundir, |
| 1367 |
|
|
float *ground_radiance) { |
| 1368 |
|
|
|
| 1369 |
|
|
const int intersect = ray_hits_ground(radius, mu); |
| 1370 |
|
|
|
| 1371 |
|
|
if (intersect) { |
| 1372 |
|
|
FVECT point, normal; |
| 1373 |
|
|
const double distance = distance_to_earth(radius, mu); |
| 1374 |
|
|
|
| 1375 |
greg |
2.2 |
/* direct + indirect irradiance */ |
| 1376 |
greg |
2.1 |
VSUM(point, view_point, view_direction, distance); |
| 1377 |
|
|
VCOPY(normal, point); |
| 1378 |
|
|
normalize(normal); |
| 1379 |
|
|
double irradiance[NSSAMP] = {0}; |
| 1380 |
|
|
get_irradiance(tau, irrad, ER, point, normal, sundir, irradiance); |
| 1381 |
|
|
|
| 1382 |
greg |
2.2 |
/* transmittance between view point and ground point */ |
| 1383 |
greg |
2.1 |
double trans[NSSAMP] = {0}; |
| 1384 |
|
|
get_transmittance(tau, radius, mu, distance, intersect, trans); |
| 1385 |
|
|
|
| 1386 |
greg |
2.2 |
/* inscattering */ |
| 1387 |
greg |
2.1 |
float inscatter[NSSAMP] = {0}; |
| 1388 |
|
|
get_sky_radiance(scat, scat1m, radius, mu, sun_ct, nu, inscatter); |
| 1389 |
|
|
|
| 1390 |
|
|
for (int i = 0; i < NSSAMP; ++i) { |
| 1391 |
|
|
ground_radiance[i] = inscatter[i] + irradiance[i] * trans[i] * grefl / PI; |
| 1392 |
|
|
} |
| 1393 |
|
|
} |
| 1394 |
|
|
} |
| 1395 |
|
|
|
| 1396 |
|
|
void get_solar_radiance(DATARRAY *tau, DATARRAY *scat, DATARRAY *scat1m, const FVECT sundir, const double radius, const double sun_ct, double *sun_radiance) { |
| 1397 |
|
|
double trans_sun[NSSAMP] = {0}; |
| 1398 |
|
|
float sky_radiance_sun[NSSAMP] = {0}; |
| 1399 |
|
|
double u, v, w, z; |
| 1400 |
|
|
const double nu = fdot(sundir, sundir); |
| 1401 |
|
|
get_transmittance_to_sun(tau, radius, sun_ct, trans_sun); |
| 1402 |
|
|
get_sky_radiance(scat, scat1m, radius, sun_ct, sun_ct, nu, sky_radiance_sun); |
| 1403 |
|
|
for (int i = 0; i < NSSAMP; ++i) { |
| 1404 |
|
|
sun_radiance[i] = sky_radiance_sun[i] + trans_sun[i] * EXTSOL[i] / SOLOMG; |
| 1405 |
|
|
} |
| 1406 |
|
|
} |