NASA Logo
Ocean Color Science Software

ocssw V2022
dtranbrdf.cpp
Go to the documentation of this file.
1 #include <cmath>
2 #include <string>
3 #include <iostream>
4 #include <fstream>
5 #include <sstream>
6 
7 #define NBIG 10U
8 #define NRAD 31U
9 #define NPHI 4U
10 #define NPHASE 16U
11 #define NWAVE 6U
12 #define NGAUS (2U * NRAD - 1U)
13 #define NNG 50U
14 #define NG (2U * NNG)
15 #define NUM 75U
16 #define NCHL 6U
17 #define NSUN 6U
18 #ifndef M_PI
19 #define M_PI 3.141592653589793238462643383279502884L
20 #endif
21 
22 static float PHSA[NPHASE][NWAVE][NBIG][NGAUS][NPHI], PHSR[NBIG][NGAUS][NPHI];
23 static float APHSRADA[NPHASE][NWAVE][NBIG][NGAUS][NPHI];
24 static float MU[NGAUS], PDIV[NG], PWT[NG], THETA[NGAUS];
25 static float PHSRADA[NGAUS][NGAUS][NPHI], PHSRADR[NGAUS][NGAUS][NPHI];
26 static float TAUA_RAT[NPHASE][NWAVE];
27 static float RAD1[NRAD][NPHI];
28 static float tst[2], tdf[2];
29 static float TSTARR[NWAVE][NRAD], TSTARA[NPHASE][NWAVE][NRAD];
30 namespace morel_vars {
31  float BRDF[NWAVE][NSUN][NCHL][NRAD][NPHI];
32  float SmBRDF[NWAVE][NSUN][NCHL][NBIG][NPHI];
34  float LChl[NCHL];
35 } // namespace morel_vars
36 bool initialized = false;
37 
38 const float aindex = 1.334f;
39 
40 float fresref(float muair, float index = aindex) {
41  float theta1 = std::acos(muair);
42  float stheta1 = std::sin(theta1);
43  float stheta2 = stheta1 / index;
44  float theta2 = std::asin(stheta2);
45  float fresref = ((index - 1) / (index + 1)) * ((index - 1) / (index + 1));
46  if (theta1 > 0.0) {
47  float tan_ratio = std::tan(theta1 - theta2) / std::tan(theta1 + theta2);
48  float sin_ratio = std::sin(theta1 - theta2) / std::sin(theta1 + theta2);
49  tan_ratio *= tan_ratio;
50  sin_ratio *= sin_ratio;
51  fresref = 0.5f * (tan_ratio + sin_ratio);
52  }
53  return fresref;
54 }
55 
56 void read_morel_data(const std::string &data_path) {
57  using namespace morel_vars;
58  std::string INFL_FOURIER31 =data_path + "/NEW_Morel_NRAD31-1-EDITED";
59  std::string INFL_FOURIER10 = data_path + "/NEW_Morel_SMALL-1-EDITED";
61 
62  {
63  std::ifstream lut_file(INFL_FOURIER31.c_str(), std::ios::in);
64  if (!lut_file.good()) {
65  std::cerr << "-Error opening the Morel LUT txt file: " << INFL_FOURIER31 << std::endl;
66  exit(EXIT_FAILURE);
67  }
68  for (size_t iwave = 0; iwave < NWAVE; iwave++) {
69  for (size_t isun = 0; isun < NSUN; isun++) {
70  for (size_t ichl = 0; ichl < NCHL; ichl++) {
71  std::getline(lut_file, line);
72  // std::cout << line << std::endl;
73  {
74  std::istringstream iss(line);
75  iss >> Wave[iwave] >> Theta0[isun] >> Chl[ichl];
76  }
77  for (size_t iphi = 0; iphi < NPHI; iphi++) {
78  std::getline(lut_file, line);
79  // std::cout << line << std::endl;
80  std::string number_line;
81  for (size_t il = 0; il < 7; il++) {
82  std::getline(lut_file, line);
83  number_line += line;
84  }
85  // std::cout << number_line << std::endl;
86  std::istringstream iss(number_line);
87  for (size_t irad = 0; irad < NRAD; irad++) {
88  iss >> BRDF[iwave][isun][ichl][irad][iphi];
89  // std::cout << BRDF[iwave][isun][ichl][irad][iphi]
90  // << std::endl;
91  }
92  }
93  }
94  }
95  }
96  lut_file.close();
97  }
98 
99  {
100  std::ifstream lut_file(INFL_FOURIER10.c_str(), std::ios::in);
101  if (!lut_file.good()) {
102  std::cerr << "-Error opening the Morel LUT txt file: " << INFL_FOURIER10 << std::endl;
103  exit(EXIT_FAILURE);
104  }
105  for (size_t iwave = 0; iwave < NWAVE; iwave++) {
106  for (size_t isun = 0; isun < NSUN; isun++) {
107  for (size_t ichl = 0; ichl < NCHL; ichl++) {
108  std::getline(lut_file, line);
109  // std::cout << line << std::endl;
110  {
111  std::istringstream iss(line);
112  iss >> Wave[iwave] >> Theta0[isun] >> Chl[ichl];
113  }
114  for (size_t iphi = 0; iphi < NPHI; iphi++) {
115  std::getline(lut_file, line);
116  // std::cout << line << std::endl;
117  std::string number_line;
118  for (size_t il = 0; il < 2; il++) {
119  std::getline(lut_file, line);
120  number_line += line;
121  }
122  // std::cout << number_line << std::endl;
123  std::istringstream iss(number_line);
124  for (size_t ibig = 0; ibig < NBIG; ibig++) {
125  iss >> SmBRDF[iwave][isun][ichl][ibig][iphi];
126  // std::cout <<
127  // SmBRDF[iwave][isun][ichl][ibig][iphi]
128  // << std::endl;
129  }
130  }
131  }
132  }
133  }
134  lut_file.close();
135  }
136  for (size_t i = 0; i < NCHL; i++) {
137  LChl[i] = std::log10(Chl[i]);
138  }
139 }
140 
142  std::string ocdataroot = std::getenv("OCDATAROOT");
143  if (ocdataroot.empty()) {
144  std::cerr << "-Error- : OCDATAROOT is not defined. Exiting ..." << std::endl;
145  exit(EXIT_FAILURE);
146  }
147  std::string common_path = ocdataroot + "/eval" + "/common" + "/dtran_brdf";
148  // file pathes
149  std::string infl_aer = common_path +"/Aerosols_Partial_Inegr.dat";
150  std::string infl_ray = common_path + "/Rayleigh_Partial_Inegr.dat";
151  std::string infl_rat = common_path + "/spec_var_EDITED.dat";
152  //
153  std::ifstream aer(infl_aer.c_str(), std::ios::in);
154  if (!aer.good()) {
155  std::cerr << "-Error-: opening the Aerosol LUT txt file: " << infl_aer << std::endl;
156  exit(EXIT_FAILURE);
157  }
159  size_t line_counter = 0;
160  // read the text Aerosol LUT file line by line.
161  // let's convert it to netcdf
162  for (auto &iphase: PHSA) {
163  for (auto &iwave: iphase) {
164  for (size_t m = 0; m < NPHI; m++) {
165  for (size_t jup = 0; jup < NGAUS; jup++) {
166  line_counter++;
167  std::getline(aer, line); // read the header
168  {
169  line_counter++;
170  std::getline(aer, line);
171  // std::cout << line << std::endl;
172  std::istringstream iss(line);
173  if (!(iss >> iwave[0][jup][m] >> iwave[1][jup][m] >> iwave[2][jup][m] >>
174  iwave[3][jup][m] >> iwave[4][jup][m])) {
175  std::cerr << "-Error-: Couldn't parse the line # " << line_counter << " : "
176  << line << std::endl;
177  exit(EXIT_FAILURE);
178  };
179  }
180  {
181  line_counter++;
182  std::getline(aer, line);
183  // std::cout << line << std::endl;
184  std::istringstream iss(line);
185  if (!(iss >> iwave[5][jup][m] >> iwave[6][jup][m] >> iwave[7][jup][m] >>
186  iwave[8][jup][m] >> iwave[9][jup][m])) {
187  std::cerr << "-Error-: Couldn't parse the line # " << line_counter << " : "
188  << line << std::endl;
189  exit(EXIT_FAILURE);
190  };
191  }
192  }
193  }
194  }
195  }
196  aer.close();
197 
198  // read the text Rayleigh LUT file line by line.
199  line_counter = 0;
200  std::ifstream rayleigh(infl_ray.c_str(), std::ios::in);
201  if (!rayleigh.good()) {
202  std::cerr << "-Error-: opening the Rayleigh LUT txt file: " << infl_ray << std::endl;
203  exit(EXIT_FAILURE);
204  }
205 
206  for (size_t m = 0; m < NPHI; m++) {
207  for (size_t jup = 0; jup < NGAUS; jup++) {
208  line_counter++;
209  std::getline(rayleigh, line); // read the header
210  {
211  line_counter++;
212  std::getline(rayleigh, line);
213  // std::cout << line << std::endl;
214  std::istringstream iss(line);
215  if (!(iss >> PHSR[0][jup][m] >> PHSR[1][jup][m] >> PHSR[2][jup][m] >> PHSR[3][jup][m] >>
216  PHSR[4][jup][m])) {
217  std::cerr << "-Error-: Couldn't parse the line # " << line_counter << " : " << line
218  << std::endl;
219  exit(EXIT_FAILURE);
220  };
221  }
222  {
223  line_counter++;
224  std::getline(rayleigh, line);
225  // std::cout << line << std::endl;
226  std::istringstream iss(line);
227  if (!(iss >> PHSR[5][jup][m] >> PHSR[6][jup][m] >> PHSR[7][jup][m] >> PHSR[8][jup][m] >>
228  PHSR[9][jup][m])) {
229  std::cerr << "-Error-: Couldn't parse the line # " << line_counter << " : " << line
230  << std::endl;
231  exit(EXIT_FAILURE);
232  };
233  }
234  }
235  }
236  rayleigh.close();
237  // read the text Tau LUT file line by line.
238  std::ifstream rat(infl_rat.c_str(), std::ios::in);
239  if (!rat.good()) {
240  std::cerr << "-Error-: opening the Tau LUT txt file: " << infl_rat << std::endl;
241  exit(EXIT_FAILURE);
242  }
243 
244  for (auto &iphase: TAUA_RAT) {
245  std::getline(rat, line); // read the header
246  for (float &iwave: iphase) {
247  std::getline(rat, line);
248  // std::cout << line << std::endl;
249  std::istringstream iss(line);
250  size_t wv;
251  iss >> wv >> iwave;
252  // std::cout << wv << " = tau is " << TAUA_RAT[iphase][iwave]
253  // << std::endl;
254  }
255  }
256  rat.close();
257 
258  // read transmittance LUTs
259 
260  std::string t_aer_path = common_path + "/tstar_aerosol.dat";
261  std::string t_ray_path = common_path + "/tstar_rayleigh.dat";
262 
263  std::ifstream t_aer_file(t_aer_path.c_str(), std::ios::in);
264  if (!t_aer_file.good()) {
265  std::cerr << "-Error-: opening the transmittance LUT txt file: " << t_aer_path << std::endl;
266  exit(EXIT_FAILURE);
267  }
268 
269  for (auto &iphase: TSTARA) {
270  for (auto &iwave: iphase) {
271  std::getline(t_aer_file, line);
272  // std::cout << line << std::endl;
273  std::string number_line;
274  for (size_t il = 0; il < 6; il++) {
275  std::getline(t_aer_file, line);
276  number_line += line;
277  }
278  // std::cout << number_line << std::endl;
279  std::istringstream iss(number_line);
280  for (size_t irad = 0; irad < NRAD - 3; irad++) {
281  iss >> iwave[irad];
282  std::cout << iwave[irad] << std::endl;
283  }
284  }
285  }
286  t_aer_file.close();
287 
288  std::ifstream t_ray_file(t_ray_path.c_str(), std::ios::in);
289  if (!t_ray_file.good()) {
290  std::cerr << "-Error-: opening the transmittance LUT txt file: " << t_ray_path << std::endl;
291  exit(EXIT_FAILURE);
292  }
293  for (auto &iwave: TSTARR) {
294  std::getline(t_ray_file, line);
295  // std::cout << line << std::endl;
296  std::string number_line;
297  for (size_t il = 0; il < 6; il++) {
298  std::getline(t_ray_file, line);
299  number_line += line;
300  }
301  // std::cout << number_line << std::endl;
302  std::istringstream iss(number_line);
303  for (size_t irad = 0; irad < NRAD - 3; irad++) {
304  iss >> iwave[irad];
305  // std::cout << TSTARR[iwave][irad] << std::endl;
306  }
307  }
308  t_ray_file.close();
309  read_morel_data(common_path.c_str());
310 }
311 
312 void Morel_BRDF(float Sun, float Chlor, float MBRDF[NWAVE][NRAD][NPHI], float SmMBRDF[NWAVE][NBIG][NPHI]) {
313  using namespace morel_vars;
314  int32_t isun = -1;
315  int32_t ichl = -1;
316  for (size_t i = 0; i < NSUN; i++) {
317  if (Sun < Theta0[i]) {
318  isun = i;
319  break;
320  }
321  }
322  for (size_t i = 0; i < NCHL; i++) {
323  if (Chlor < Chl[i]) {
324  ichl = i;
325  break;
326  }
327  }
328  float LChlor = std::log10(Chlor);
329  // C Interpolate in the tables to get the water radiance
330 
331  float sun_interp = (Sun - Theta0[isun - 1]) / (Theta0[isun] - Theta0[isun - 1]);
332  float chl_interp = (LChlor - LChl[ichl - 1]) / (LChl[ichl] - LChl[ichl - 1]);
333 
334  for (size_t iwave = 0; iwave < NWAVE; iwave++) {
335  for (size_t i = 0; i < NRAD; i++) {
336  for (size_t m = 0; m < NPHI; m++) {
337  float interp_term =
338  (1. - sun_interp) * (1. - chl_interp) * BRDF[iwave][isun - 1][ichl - 1][i][m];
339  interp_term =
340  interp_term + sun_interp * (1. - chl_interp) * BRDF[iwave][isun][ichl - 1][i][m];
341  interp_term =
342  interp_term + (1. - sun_interp) * chl_interp * BRDF[iwave][isun - 1][ichl][i][m];
343  interp_term = interp_term + sun_interp * chl_interp * BRDF[iwave][isun][ichl][i][m];
344  MBRDF[iwave][i][m] = interp_term;
345  }
346  }
347  }
348  for (size_t iwave = 0; iwave < NWAVE; iwave++) {
349  for (size_t i = 0; i < NBIG; i++) {
350  for (size_t m = 0; m < NPHI; m++) {
351  float interp_term =
352  (1. - sun_interp) * (1. - chl_interp) * SmBRDF[iwave][isun - 1][ichl - 1][i][m];
353  interp_term =
354  interp_term + sun_interp * (1. - chl_interp) * SmBRDF[iwave][isun][ichl - 1][i][m];
355  interp_term =
356  interp_term + (1. - sun_interp) * chl_interp * SmBRDF[iwave][isun - 1][ichl][i][m];
357  interp_term = interp_term + sun_interp * chl_interp * SmBRDF[iwave][isun][ichl][i][m];
358  SmMBRDF[iwave][i][m] = interp_term;
359  }
360  }
361  }
362 }
363 
364 extern "C" void diff_tran_corr_(int *iphase, float *solz, float *senz, float *phi, float *chl, float *taua,
365  float *correct) {
366  std::cout
367  << "\n\n\nDifftran BRDF exits; Uncomment to continue. Exiting ... Source file - dtranbrdf.cpp. Now we know the usage!\n\n\n\n"
368  << std::endl;
369  exit(EXIT_FAILURE);
370  float MBRDF[NWAVE][NRAD][NPHI];
371  float SmMBRDF[NWAVE][NBIG][NPHI];
372  static float tr[] = {0.3132, 0.2336, 0.1547, 0.1330, 0.0957, 0.0446};
373  int32_t iview = -1;
374  if (!initialized) {
375  read_luts_txt();
376  for (size_t i = 0; i < NGAUS; i++) {
377  THETA[i] = static_cast<float>(i) * 90.0f / static_cast<float>(NRAD - 1);
378  MU[i] = std::cos(THETA[i] * M_PI / 180.0);
379  }
380  initialized = true;
381  }
382  Morel_BRDF(*solz, *chl, MBRDF, SmMBRDF);
383  for (size_t i = 0; i < NRAD; i++) {
384  if (*senz < THETA[i]) {
385  iview = i;
386  break;
387  }
388  }
389  float Ta = *taua;
390  size_t mgaus = NGAUS;
391  size_t maxphi = NPHI;
392  for (size_t iwave = 0; iwave < NWAVE; iwave++) {
393  float tstar1, tstar2, ans1, ans2, tdiff1, tdiff2;
394  for (size_t jup = iview - 1; jup <= iview; jup++) {
395  // float adelphi = (180.0 / M_PI) * (*phi);
396  size_t jdn = mgaus + 1 - jup;
397  float ya1 = 0.0f;
398  float yr1 = 0.0f;
399  for (size_t m = 0; m < maxphi; m++) {
400  float order = static_cast<float>(m);
401  float x1 = 0;
402  float v1 = 0;
403  for (size_t i = 0; i < NBIG; i++) {
404  float x2 = PHSA[*iphase][iwave][i][jup][m] * SmMBRDF[iwave][i][m];
405  float v2 = PHSR[i][jup][m] * SmMBRDF[iwave][i][m];
406  x1 += x2;
407  v1 += v2;
408  }
409  if (m > 0) {
410  x1 *= 2;
411  v1 *= 2;
412  }
413  ya1 += std::cos(order * *phi) * x1;
414  yr1 += std::cos(order * *phi) * v1;
415  }
416  float za1 = 0.0f;
417  float zr1 = 0.0f;
418  for (size_t m = 0; m < maxphi; m++) {
419  float order = static_cast<float>(m);
420  float x1 = 0;
421  float v1 = 0;
422  for (size_t i = 0; i < NBIG; i++) {
423  float x2 = PHSA[*iphase][iwave][i][jdn][m] * SmMBRDF[iwave][i][m];
424  float v2 = PHSR[i][jdn][m] * SmMBRDF[iwave][i][m];
425  x1 += x2;
426  v1 += v2;
427  }
428  if (m > 0) {
429  x1 *= 2;
430  v1 *= 2;
431  }
432  za1 += std::cos(order * *phi) * x1;
433  zr1 += std::cos(order * *phi) * v1;
434  }
435  float yl = 0.0f;
436  for (size_t m = 0; m < maxphi; m++) {
437  float order = static_cast<float>(m);
438  float fac = m > 0 ? 2.0f : 1.0f;
439  yl += MBRDF[iwave][jup][m] * fac * std::cos(order * *phi);
440  }
441  float rfres = fresref(MU[jup]);
442  float aint_pl_a = (ya1 + rfres * za1) / yl / (1.0f - rfres);
443  float aint_pl_r = (yr1 + rfres * zr1) / yl / (1.0f - rfres);
444  float plterm_a = (1.0f - aint_pl_a / 2.0f) / MU[jup];
445  float plterm_r = (1.0f - aint_pl_r / 2.0f) / MU[jup];
446  float t_diff_a = std::exp(-plterm_a * Ta * TAUA_RAT[*iphase][iwave]);
447  float t_diff_r = std::exp(-plterm_r * tr[iwave]);
448  if (jup == iview - 1) {
449  tstar1 =
450  TSTARR[iwave][jup] * std::pow(TSTARA[*iphase][iwave][jup], Ta * TAUA_RAT[*iphase][iwave]);
451  tdiff1 = t_diff_a * t_diff_r;
452  ans1 = (tdiff1 - tstar1) / (tstar1);
453  } else {
454  tstar2 =
455  TSTARR[iwave][jup] * std::pow(TSTARA[*iphase][iwave][jup], Ta * TAUA_RAT[*iphase][iwave]);
456  tdiff2 = t_diff_a * t_diff_r;
457  ans2 = (tdiff2 - tstar2) / (tstar2);
458  }
459  }
460  float slope_ans = (ans2 - ans1) / (THETA[iview] - THETA[iview - 1]);
461 // float slope_tst = (tstar2 - tstar1) / (THETA[iview] - THETA[iview - 1]);
462 // float slope_tdf = (tdiff2 - tdiff1) / (THETA[iview] - THETA[iview - 1]);
463  correct[iwave] = ans1 + slope_ans * (*senz - THETA[iview - 1]);
464  }
465 };
an array had not been initialized Several spelling and grammar corrections were which is read from the appropriate MCF the above metadata values were hard coded A problem calculating the average background DN for SWIR bands when the moon is in the space view port was corrected The new algorithm used to calculate the average background DN for all reflective bands when the moon is in the space view port is now the same as the algorithm employed by the thermal bands For non SWIR changes in the averages are typically less than Also for non SWIR the black body DNs remain a backup in case the SV DNs are not available For SWIR the changes in computed averages were larger because the old which used the black body suffered from contamination by the micron leak As a consequence of the if SV DNs are not available for the SWIR the EV pixels will not be the granule time is used to identify the appropriate tables within the set given for one LUT the first two or last two tables respectively will be used for the interpolation If there is only one LUT in the set of it will be treated as a constant LUT The manner in which Earth View data is checked for saturation was changed Previously the raw Earth View DNs and Space View DNs were checked against the lookup table values contained in the table dn_sat The change made is to check the raw Earth and Space View DNs to be sure they are less than the maximum saturation value and to check the Space View subtracted Earth View dns against a set of values contained in the new lookup table dn_sat_ev The metadata configuration and ASSOCIATEDINSTRUMENTSHORTNAME from the MOD02HKM product The same metatdata with extensions and were removed from the MOD021KM and MOD02OBC products ASSOCIATEDSENSORSHORTNAME was set to MODIS in all products These changes are reflected in new File Specification which users may consult for exact the pow functions were eliminated in Emissive_Cal and Emissive bands replaced by more efficient code Other calculations throughout the code were also made more efficient Aside from a few round off there was no difference to the product The CPU time decreased by about for a day case and for a night case A minor bug in calculating the uncertainty index for emissive bands was corrected The frame index(0-based) was previously being used the frame number(1-based) should have been used. There were only a few minor changes to the uncertainty index(maximum of 1 digit). 3. Some inefficient arrays(Sigma_RVS_norm_sq) were eliminated and some code lines in Preprocess_L1A_Data were moved into Process_OBCEng_Emiss. There were no changes to the product. Required RAM was reduced by 20 MB. Now
#define NCHL
Definition: dtranbrdf.cpp:16
#define NPHASE
Definition: dtranbrdf.cpp:10
void Morel_BRDF(float Sun, float Chlor, float MBRDF[NWAVE][NRAD][NPHI], float SmMBRDF[NWAVE][NBIG][NPHI])
Definition: dtranbrdf.cpp:312
float LChl[NCHL]
Definition: dtranbrdf.cpp:34
void rayleigh(l1str *l1rec, int32_t ip)
Definition: rayleigh.c:169
@ string
float Theta0[NSUN]
Definition: dtranbrdf.cpp:33
#define M_PI
Definition: dtranbrdf.cpp:19
float fresref(float muair, float index=aindex)
Definition: dtranbrdf.cpp:40
#define fac
#define NG
Definition: dtranbrdf.cpp:14
float BRDF[NWAVE][NSUN][NCHL][NRAD][NPHI]
Definition: dtranbrdf.cpp:31
#define NRAD
Definition: dtranbrdf.cpp:8
double precision function f(R1)
Definition: tmd.lp.f:1454
float SmBRDF[NWAVE][NSUN][NCHL][NBIG][NPHI]
Definition: dtranbrdf.cpp:32
const char * ocdataroot
#define NPHI
Definition: dtranbrdf.cpp:9
void diff_tran_corr_(int *iphase, float *solz, float *senz, float *phi, float *chl, float *taua, float *correct)
Definition: dtranbrdf.cpp:364
#define NGAUS
Definition: dtranbrdf.cpp:12
float Thetav[NRAD]
Definition: dtranbrdf.cpp:33
float Wave[NWAVE]
Definition: dtranbrdf.cpp:33
float Chl[NCHL]
Definition: dtranbrdf.cpp:33
an array had not been initialized Several spelling and grammar corrections were which is read from the appropriate MCF the above metadata values were hard coded A problem calculating the average background DN for SWIR bands when the moon is in the space view port was corrected The new algorithm used to calculate the average background DN for all reflective bands when the moon is in the space view port is now the same as the algorithm employed by the thermal bands For non SWIR changes in the averages are typically less than Also for non SWIR the black body DNs remain a backup in case the SV DNs are not available For SWIR the changes in computed averages were larger because the old which used the black body suffered from contamination by the micron leak As a consequence of the if SV DNs are not available for the SWIR the EV pixels will not be the granule time is used to identify the appropriate tables within the set given for one LUT the first two or last two tables respectively will be used for the interpolation If there is only one LUT in the set of it will be treated as a constant LUT The manner in which Earth View data is checked for saturation was changed Previously the raw Earth View DNs and Space View DNs were checked against the lookup table values contained in the table dn_sat The change made is to check the raw Earth and Space View DNs to be sure they are less than the maximum saturation value and to check the Space View subtracted Earth View dns against a set of values contained in the new lookup table dn_sat_ev The metadata configuration and ASSOCIATEDINSTRUMENTSHORTNAME from the MOD02HKM product The same metatdata with extensions and were removed from the MOD021KM and MOD02OBC products ASSOCIATEDSENSORSHORTNAME was set to MODIS in all products These changes are reflected in new File Specification which users may consult for exact the pow functions were eliminated in Emissive_Cal and Emissive bands replaced by more efficient code Other calculations throughout the code were also made more efficient Aside from a few round off there was no difference to the product The CPU time decreased by about for a day case and for a night case A minor bug in calculating the uncertainty index for emissive bands was corrected The frame the required RAM for each execution is MB on the DEC ALPHA and MB on the SGI Octane v2
Definition: HISTORY.txt:728
bool initialized
Definition: dtranbrdf.cpp:36
void read_luts_txt()
Definition: dtranbrdf.cpp:141
const float aindex
Definition: dtranbrdf.cpp:38
#define NWAVE
Definition: dtranbrdf.cpp:11
int i
Definition: decode_rs.h:71
#define NBIG
Definition: dtranbrdf.cpp:7
#define NSUN
Definition: dtranbrdf.cpp:17
void read_morel_data(const std::string &data_path)
Definition: dtranbrdf.cpp:56