JURASSIC
Functions
formod.c File Reference

JURASSIC forward model. More...

#include "jurassic.h"

Go to the source code of this file.

Functions

void call_formod (ctl_t *ctl, const tbl_t *tbl, const char *wrkdir, const char *obsfile, const char *atmfile, const char *radfile, const char *task, const char *obsref)
 Perform forward model calculations in a single directory. More...
 
void compute_rel_errors (const ctl_t *ctl, const obs_t *obs_test, const obs_t *obs_ref, double *mre, double *sdre, double *minre, double *maxre)
 Calculate relative errors. More...
 
int main (int argc, char *argv[])
 

Detailed Description

JURASSIC forward model.

Definition in file formod.c.

Function Documentation

◆ call_formod()

void call_formod ( ctl_t ctl,
const tbl_t tbl,
const char *  wrkdir,
const char *  obsfile,
const char *  atmfile,
const char *  radfile,
const char *  task,
const char *  obsref 
)

Perform forward model calculations in a single directory.

Definition at line 143 of file formod.c.

151 {
152
153 static atm_t atm, atm2;
154 static obs_t obs, obs2;
155
156 /* Read atmospheric data... */
157 read_atm(wrkdir, atmfile, ctl, &atm);
158
159 /* Read observation geometry... */
160 read_obs(wrkdir, obsfile, ctl, &obs);
161
162 /* Compute multiple profiles... */
163 if (task[0] == 'p' || task[0] == 'P') {
164
165 /* Loop over ray paths... */
166 for (int ir = 0; ir < obs.nr; ir++) {
167
168 /* Get atmospheric data... */
169 atm2.np = 0;
170 for (int ip = 0; ip < atm.np; ip++)
171 if (atm.time[ip] == obs.time[ir]) {
172 atm2.time[atm2.np] = atm.time[ip];
173 atm2.z[atm2.np] = atm.z[ip];
174 atm2.lon[atm2.np] = atm.lon[ip];
175 atm2.lat[atm2.np] = atm.lat[ip];
176 atm2.p[atm2.np] = atm.p[ip];
177 atm2.t[atm2.np] = atm.t[ip];
178 for (int ig = 0; ig < ctl->ng; ig++)
179 atm2.q[ig][atm2.np] = atm.q[ig][ip];
180 for (int iw = 0; iw < ctl->nw; iw++)
181 atm2.k[iw][atm2.np] = atm.k[iw][ip];
182 atm2.np++;
183 }
184
185 /* Get observation data... */
186 obs2.nr = 1;
187 obs2.time[0] = obs.time[ir];
188 obs2.vpz[0] = obs.vpz[ir];
189 obs2.vplon[0] = obs.vplon[ir];
190 obs2.vplat[0] = obs.vplat[ir];
191 obs2.obsz[0] = obs.obsz[ir];
192 obs2.obslon[0] = obs.obslon[ir];
193 obs2.obslat[0] = obs.obslat[ir];
194
195 /* Check number of data points... */
196 if (atm2.np > 0) {
197
198 /* Call forward model... */
199 formod(ctl, tbl, &atm2, &obs2);
200
201 /* Save radiance data... */
202 for (int id = 0; id < ctl->nd; id++) {
203 obs.rad[id][ir] = obs2.rad[id][0];
204 obs.tau[id][ir] = obs2.tau[id][0];
205 }
206 }
207 }
208
209 /* Write radiance data... */
210 write_obs(wrkdir, radfile, ctl, &obs);
211 }
212
213 /* Compute single profile... */
214 else {
215
216 /* Call forward model... */
217 formod(ctl, tbl, &atm, &obs);
218
219 /* Save radiance data... */
220 write_obs(wrkdir, radfile, ctl, &obs);
221
222 /* Evaluate results... */
223 if (obsref[0] != '-') {
224
225 /* Read reference data... */
226 read_obs(wrkdir, obsref, ctl, &obs2);
227
228 /* Calculate relative errors... */
229 double mre[ND], sdre[ND], minre[ND], maxre[ND];
230 compute_rel_errors(ctl, &obs, &obs2, mre, sdre, minre, maxre);
231
232 /* Write results... */
233 for (int id = 0; id < ctl->nd; id++)
234 printf
235 ("EVAL: nu= %.4f cm^-1 | MRE= %g %% | SDRE= %g %% | MinRE= %g %% | MaxRE= %g %%\n",
236 ctl->nu[id], mre[id], sdre[id], minre[id], maxre[id]);
237 }
238
239 /* Compute contributions of emitters... */
240 if (task[0] == 'c' || task[0] == 'C') {
241
242 char filename[LEN];
243
244 /* Switch off continua... */
245 ctl->ctm_co2 = 0;
246 ctl->ctm_h2o = 0;
247 ctl->ctm_n2 = 0;
248 ctl->ctm_o2 = 0;
249
250 /* Loop over emitters... */
251 for (int ig = 0; ig < ctl->ng; ig++) {
252
253 /* Copy atmospheric data... */
254 copy_atm(ctl, &atm2, &atm, 0);
255
256 /* Set extinction to zero... */
257 for (int iw = 0; iw < ctl->nw; iw++)
258 for (int ip = 0; ip < atm2.np; ip++)
259 atm2.k[iw][ip] = 0;
260
261 /* Select emitter... */
262 for (int ig2 = 0; ig2 < ctl->ng; ig2++)
263 if (ig2 != ig)
264 for (int ip = 0; ip < atm2.np; ip++)
265 atm2.q[ig2][ip] = 0;
266
267 /* Call forward model... */
268 formod(ctl, tbl, &atm2, &obs);
269
270 /* Save radiance data... */
271 sprintf(filename, "%s.%s", radfile, ctl->emitter[ig]);
272 write_obs(wrkdir, filename, ctl, &obs);
273 }
274
275 /* Copy atmospheric data... */
276 copy_atm(ctl, &atm2, &atm, 0);
277
278 /* Set volume mixing ratios to zero, keep extinction... */
279 for (int ig = 0; ig < ctl->ng; ig++)
280 for (int ip = 0; ip < atm2.np; ip++)
281 atm2.q[ig][ip] = 0;
282
283 /* Call forward model... */
284 formod(ctl, tbl, &atm2, &obs);
285
286 /* Save radiance data... */
287 sprintf(filename, "%s.EXTINCT", radfile);
288 write_obs(wrkdir, filename, ctl, &obs);
289 }
290
291 /* Measure CPU-time... */
292 if (task[0] == 't' || task[0] == 'T') {
293
294 /* Init... */
295 double t_min, t_max, t_mean = 0, t_sd = 0;
296 int n = 0;
297
298 /* Initialize random number generator... */
299 gsl_rng_env_setup();
300 gsl_rng *rng = gsl_rng_alloc(gsl_rng_default);
301
302 /* Loop over profiles... */
303 do {
304
305 /* Create random atmosphere... */
306 copy_atm(ctl, &atm2, &atm, 0);
307 double dtemp = 40. * (gsl_rng_uniform(rng) - 0.5);
308 double dpress = 1. - 0.1 * gsl_rng_uniform(rng);
309 double dq[NG];
310 for (int ig = 0; ig < ctl->ng; ig++)
311 dq[ig] = 0.8 + 0.4 * gsl_rng_uniform(rng);
312 for (int ip = 0; ip < atm2.np; ip++) {
313 atm2.t[ip] += dtemp;
314 atm2.p[ip] *= dpress;
315 for (int ig = 0; ig < ctl->ng; ig++)
316 atm2.q[ig][ip] *= dq[ig];
317 }
318
319 /* Measure runtime... */
320 double t0 = omp_get_wtime();
321 formod(ctl, tbl, &atm2, &obs);
322 double dt = omp_get_wtime() - t0;
323
324 /* Get runtime statistics... */
325 t_mean += dt;
326 t_sd += POW2(dt);
327 if (n == 0 || dt < t_min)
328 t_min = dt;
329 if (n == 0 || dt > t_max)
330 t_max = dt;
331 n++;
332
333 } while (t_mean < 10.0);
334
335 /* Write results... */
336 t_mean /= (double) n;
337 t_sd = sqrt(t_sd / (double) n - POW2(t_mean));
338 printf("RUNTIME: mean= %g s | stddev= %g s | min= %g s | max= %g s\n",
339 t_mean, t_sd, t_min, t_max);
340
341 /* Free... */
342 gsl_rng_free(rng);
343 }
344
345 /* Analyze impact of step size... */
346 if (task[0] == 's' || task[0] == 'S') {
347
348 /* Reference run... */
349 ctl->rayds = 0.1;
350 ctl->raydz = 0.01;
351 formod(ctl, tbl, &atm, &obs);
352 copy_obs(ctl, &obs2, &obs, 0);
353
354 /* Loop over step size... */
355 for (double dz = 0.01; dz <= 2; dz *= 1.1)
356 for (double ds = 0.1; ds <= 50; ds *= 1.1) {
357
358 /* Set step size... */
359 ctl->rayds = ds;
360 ctl->raydz = dz;
361
362 /* Measure runtime... */
363 double t0 = omp_get_wtime();
364 formod(ctl, tbl, &atm, &obs);
365 double dt = omp_get_wtime() - t0;
366
367 /* Calculate relative errors... */
368 double mre[ND], sdre[ND], minre[ND], maxre[ND];
369 compute_rel_errors(ctl, &obs, &obs2, mre, sdre, minre, maxre);
370
371 /* Write results... */
372 for (int id = 0; id < ctl->nd; id++)
373 printf
374 ("STEPSIZE: ds= %.4f km | dz= %g km | t= %g s | nu= %.4f cm^-1"
375 " | MRE= %g %% | SDRE= %g %% | MinRE= %g %% | MaxRE= %g %%\n",
376 ds, dz, dt, ctl->nu[id], mre[id], sdre[id], minre[id], maxre[id]);
377 }
378 }
379 }
380}
void compute_rel_errors(const ctl_t *ctl, const obs_t *obs_test, const obs_t *obs_ref, double *mre, double *sdre, double *minre, double *maxre)
Calculate relative errors.
Definition: formod.c:384
void read_atm(const char *dirname, const char *filename, const ctl_t *ctl, atm_t *atm)
Read atmospheric profile data from an ASCII file.
Definition: jurassic.c:4787
void copy_obs(const ctl_t *ctl, obs_t *obs_dest, const obs_t *obs_src, const int init)
Copy or initialize observation geometry and radiance data.
Definition: jurassic.c:3308
void read_obs(const char *dirname, const char *filename, const ctl_t *ctl, obs_t *obs)
Read observation geometry and radiance data from an ASCII file.
Definition: jurassic.c:5043
void write_obs(const char *dirname, const char *filename, const ctl_t *ctl, const obs_t *obs)
Definition: jurassic.c:6202
void formod(const ctl_t *ctl, const tbl_t *tbl, atm_t *atm, obs_t *obs)
Execute the selected forward model.
Definition: jurassic.c:3359
void copy_atm(const ctl_t *ctl, atm_t *atm_dest, const atm_t *atm_src, const int init)
Copy or initialize atmospheric profile data.
Definition: jurassic.c:3258
#define LEN
Maximum length of ASCII data lines.
Definition: jurassic.h:266
#define POW2(x)
Compute the square of a value.
Definition: jurassic.h:669
#define ND
Maximum number of radiance channels.
Definition: jurassic.h:236
#define NG
Maximum number of emitters.
Definition: jurassic.h:241
Atmospheric profile data.
Definition: jurassic.h:914
double time[NP]
Time (seconds since 2000-01-01T00:00Z).
Definition: jurassic.h:920
double k[NW][NP]
Extinction [km^-1].
Definition: jurassic.h:941
double lat[NP]
Latitude [deg].
Definition: jurassic.h:929
double lon[NP]
Longitude [deg].
Definition: jurassic.h:926
double t[NP]
Temperature [K].
Definition: jurassic.h:935
int np
Number of data points.
Definition: jurassic.h:917
double z[NP]
Altitude [km].
Definition: jurassic.h:923
double q[NG][NP]
Volume mixing ratio [ppv].
Definition: jurassic.h:938
double p[NP]
Pressure [hPa].
Definition: jurassic.h:932
int nw
Number of spectral windows.
Definition: jurassic.h:994
double nu[ND]
Centroid wavenumber of each channel [cm^-1].
Definition: jurassic.h:991
int ctm_co2
Compute CO2 continuum (0=no, 1=yes).
Definition: jurassic.h:1027
int ctm_n2
Compute N2 continuum (0=no, 1=yes).
Definition: jurassic.h:1033
int ctm_h2o
Compute H2O continuum (0=no, 1=yes).
Definition: jurassic.h:1030
int ng
Number of emitters.
Definition: jurassic.h:970
int nd
Number of radiance channels.
Definition: jurassic.h:988
int ctm_o2
Compute O2 continuum (0=no, 1=yes).
Definition: jurassic.h:1036
char emitter[NG][LEN]
Name of each emitter.
Definition: jurassic.h:973
double rayds
Maximum step length for raytracing [km].
Definition: jurassic.h:1042
double raydz
Vertical step length for raytracing [km].
Definition: jurassic.h:1045
Observation geometry and radiance data.
Definition: jurassic.h:1187
double tau[ND][NR]
Transmittance of ray path.
Definition: jurassic.h:1223
double rad[ND][NR]
Radiance [W/(m^2 sr cm^-1)].
Definition: jurassic.h:1226
double vpz[NR]
View point altitude [km].
Definition: jurassic.h:1205
double vplat[NR]
View point latitude [deg].
Definition: jurassic.h:1211
double obslon[NR]
Observer longitude [deg].
Definition: jurassic.h:1199
double obslat[NR]
Observer latitude [deg].
Definition: jurassic.h:1202
double obsz[NR]
Observer altitude [km].
Definition: jurassic.h:1196
double vplon[NR]
View point longitude [deg].
Definition: jurassic.h:1208
double time[NR]
Time (seconds since 2000-01-01T00:00Z).
Definition: jurassic.h:1193
int nr
Number of ray paths.
Definition: jurassic.h:1190
Here is the call graph for this function:

◆ compute_rel_errors()

void compute_rel_errors ( const ctl_t ctl,
const obs_t obs_test,
const obs_t obs_ref,
double *  mre,
double *  sdre,
double *  minre,
double *  maxre 
)

Calculate relative errors.

Definition at line 384 of file formod.c.

391 {
392
393 /* Loop over channels... */
394 for (int id = 0; id < ctl->nd; id++) {
395
396 double sum = 0, sum2 = 0;
397 minre[id] = +1e9;
398 maxre[id] = -1e9;
399 int n = 0;
400
401 /* Loop over ray paths... */
402 for (int ir = 0; ir < obs_test->nr; ir++) {
403
404 /* Check for zero... */
405 if (obs_ref->rad[id][ir] == 0)
406 continue;
407
408 /* Calculate relative error... */
409 double err = 100.0 * (obs_test->rad[id][ir] - obs_ref->rad[id][ir])
410 / obs_ref->rad[id][ir];
411
412 /* Get statistics... */
413 sum += err;
414 sum2 += err * err;
415 if (err > maxre[id])
416 maxre[id] = err;
417 if (err < minre[id])
418 minre[id] = err;
419 n++;
420 }
421
422 /* Get mean and standard deviaton... */
423 mre[id] = sum / n;
424 sdre[id] = sqrt(sum2 / n - mre[id] * mre[id]);
425 }
426}

◆ main()

int main ( int  argc,
char *  argv[] 
)

Definition at line 59 of file formod.c.

61 {
62
63 static ctl_t ctl;
64
65 /* Check arguments... */
66 if (argc < 5)
67 ERRMSG("Give parameters: <ctl> <obs> <atm> <rad>");
68
69 /* Read control parameters... */
70 read_ctl(argc, argv, &ctl);
71
72#ifdef UNIFIED
73
74 static atm_t atm;
75 static obs_t obs;
76
77 /* Read observation geometry... */
78 read_obs(NULL, argv[2], &ctl, &obs);
79
80 /* Read atmospheric data... */
81 read_atm(NULL, argv[3], &ctl, &atm);
82
83 /* Call forward model... */
84 jur_unified_init(argc, argv);
85 jur_unified_formod_multiple_packages(&atm, &obs, 1, NULL);
86
87 /* Save radiance data... */
88 write_obs(NULL, argv[4], &ctl, &obs);
89
90#else
91
92 char dirlist[LEN], obsref[LEN], task[LEN];
93
94 /* Initialize look-up tables... */
95 tbl_t *tbl = read_tbl(&ctl);
96
97 /* Get task... */
98 scan_ctl(argc, argv, "TASK", -1, "-", task);
99
100 /* Get dirlist... */
101 scan_ctl(argc, argv, "DIRLIST", -1, "-", dirlist);
102
103 /* Get reference data... */
104 scan_ctl(argc, argv, "OBSREF", -1, "-", obsref);
105
106 /* Single forward calculation... */
107 if (dirlist[0] == '-')
108 call_formod(&ctl, tbl, NULL, argv[2], argv[3], argv[4], task, obsref);
109
110 /* Work on directory list... */
111 else {
112
113 /* Open directory list... */
114 FILE *in;
115 if (!(in = fopen(dirlist, "r")))
116 ERRMSG("Cannot open directory list!");
117
118 /* Loop over directories... */
119 char wrkdir[LEN];
120 while (fscanf(in, "%4999s", wrkdir) != EOF) {
121
122 /* Write info... */
123 LOG(1, "\nWorking directory: %s", wrkdir);
124
125 /* Call forward model... */
126 call_formod(&ctl, tbl, wrkdir, argv[2], argv[3], argv[4], task, obsref);
127 }
128
129 /* Close dirlist... */
130 fclose(in);
131 }
132
133#endif
134
135 /* Free... */
136 free(tbl);
137
138 return EXIT_SUCCESS;
139}
void call_formod(ctl_t *ctl, const tbl_t *tbl, const char *wrkdir, const char *obsfile, const char *atmfile, const char *radfile, const char *task, const char *obsref)
Perform forward model calculations in a single directory.
Definition: formod.c:143
void read_ctl(int argc, char *argv[], ctl_t *ctl)
Read model control parameters from command-line and configuration input.
Definition: jurassic.c:4889
double scan_ctl(int argc, char *argv[], const char *varname, const int arridx, const char *defvalue, char *value)
Scan control file or command-line arguments for a configuration variable.
Definition: jurassic.c:5532
tbl_t * read_tbl(const ctl_t *ctl)
Read gas emissivity look-up tables for all channels and emitters.
Definition: jurassic.c:5358
#define ERRMSG(...)
Print an error message with contextual information and terminate the program.
Definition: jurassic.h:864
#define LOG(level,...)
Print a log message with a specified logging level.
Definition: jurassic.h:794
Control parameters.
Definition: jurassic.h:967
Emissivity look-up tables.
Definition: jurassic.h:1321
Here is the call graph for this function: