CrIS Code Collection
jurassic.h
Go to the documentation of this file.
1/*
2 This file is part of JURASSIC.
3
4 JURASSIC is free software: you can redistribute it and/or modify
5 it under the terms of the GNU General Public License as published by
6 the Free Software Foundation, either version 3 of the License, or
7 (at your option) any later version.
8
9 JURASSIC is distributed in the hope that it will be useful,
10 but WITHOUT ANY WARRANTY; without even the implied warranty of
11 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12 GNU General Public License for more details.
13
14 You should have received a copy of the GNU General Public License
15 along with JURASSIC. If not, see <http://www.gnu.org/licenses/>.
16
17 Copyright (C) 2003-2025 Forschungszentrum Juelich GmbH
18*/
19
100#ifndef JURASSIC_H
101#define JURASSIC_H
102
103/* ------------------------------------------------------------
104 Includes...
105 ------------------------------------------------------------ */
106
107#include <errno.h>
108#include <gsl/gsl_math.h>
109#include <gsl/gsl_blas.h>
110#include <gsl/gsl_linalg.h>
111#include <gsl/gsl_randist.h>
112#include <gsl/gsl_rng.h>
113#include <gsl/gsl_statistics.h>
114#include <math.h>
115#include <omp.h>
116#include <stdint.h>
117#include <stdio.h>
118#include <stdlib.h>
119#include <string.h>
120#include <time.h>
121
122/* ------------------------------------------------------------
123 Constants...
124 ------------------------------------------------------------ */
125
127#ifndef C1
128#define C1 1.19104259e-8
129#endif
130
132#ifndef C2
133#define C2 1.43877506
134#endif
135
137#ifndef EPSMIN
138#define EPSMIN 0
139#endif
140
142#ifndef EPSMAX
143#define EPSMAX 1
144#endif
145
147#ifndef G0
148#define G0 9.80665
149#endif
150
152#ifndef H0
153#define H0 7.0
154#endif
155
157#ifndef KB
158#define KB 1.3806504e-23
159#endif
160
162#ifndef ME
163#define ME 5.976e24
164#endif
165
167#ifndef NA
168#define NA 6.02214199e23
169#endif
170
172#ifndef N2
173#define N2 0.78084
174#endif
175
177#ifndef O2
178#define O2 0.20946
179#endif
180
182#ifndef P0
183#define P0 1013.25
184#endif
185
187#ifndef RE
188#define RE 6367.421
189#endif
190
192#ifndef RI
193#define RI 8.3144598
194#endif
195
197#ifndef T0
198#define T0 273.15
199#endif
200
202#ifndef TMIN
203#define TMIN 100.
204#endif
205
207#ifndef TMAX
208#define TMAX 400.
209#endif
210
212#ifndef TSUN
213#define TSUN 5780.
214#endif
215
217#ifndef UMIN
218#define UMIN 0
219#endif
220
222#ifndef UMAX
223#define UMAX 1e30
224#endif
225
226/* ------------------------------------------------------------
227 Dimensions...
228 ------------------------------------------------------------ */
229
231#ifndef NCL
232#define NCL 8
233#endif
234
236#ifndef ND
237#define ND 128
238#endif
239
241#ifndef NG
242#define NG 8
243#endif
244
246#ifndef NP
247#define NP 256
248#endif
249
251#ifndef NR
252#define NR 256
253#endif
254
256#ifndef NSF
257#define NSF 8
258#endif
259
261#ifndef NW
262#define NW 4
263#endif
264
266#ifndef LEN
267#define LEN 10000
268#endif
269
271#ifndef M
272#define M (NR*ND)
273#endif
274
276#ifndef N
277#define N ((2 + NG + NW) * NP + NCL + NSF + 3)
278#endif
279
281#ifndef NQ
282#define NQ (5 + NG + NW + NCL + NSF)
283#endif
284
286#ifndef NLOS
287#define NLOS 4096
288#endif
289
291#ifndef NSHAPE
292#define NSHAPE 20000
293#endif
294
296#ifndef NFOV
297#define NFOV 5
298#endif
299
301#ifndef TBLNP
302#define TBLNP 41
303#endif
304
306#ifndef TBLNT
307#define TBLNT 30
308#endif
309
311#ifndef TBLNU
312#define TBLNU 320
313#endif
314
316#ifndef TBLNS
317#define TBLNS 1200
318#endif
319
321#ifndef MAX_TABLES
322#define MAX_TABLES 10000
323#endif
324
326#ifndef RFMNPTS
327#define RFMNPTS 10000000
328#endif
329
330/* ------------------------------------------------------------
331 Quantity indices...
332 ------------------------------------------------------------ */
333
335#define IDXP 0
336
338#define IDXT 1
339
341#define IDXQ(ig) (2 + (ig))
342
344#define IDXK(iw) (2 + (ctl->ng) + (iw))
345
347#define IDXCLZ (2 + (ctl->ng) + (ctl->nw))
348
350#define IDXCLDZ (3 + (ctl->ng) + (ctl->nw))
351
353#define IDXCLK(icl) (4 + (ctl->ng) + (ctl->nw) + (icl))
354
356#define IDXSFT (4 + (ctl->ng) + (ctl->nw) + (ctl->ncl))
357
359#define IDXSFEPS(isf) (5 + (ctl->ng) + (ctl->nw) + (ctl->ncl) + (isf))
360
361/* ------------------------------------------------------------
362 Macros...
363 ------------------------------------------------------------ */
364
382#define ALLOC(ptr, type, n) \
383 if((ptr=calloc((size_t)(n), sizeof(type)))==NULL) \
384 ERRMSG("Out of memory!");
385
412#define BRIGHT(rad, nu) \
413 (C2 * (nu) / gsl_log1p(C1 * POW3(nu) / (rad)))
414
429#define DEG2RAD(deg) ((deg) * (M_PI / 180.0))
430
447#define DIST(a, b) sqrt(DIST2(a, b))
448
464#define DIST2(a, b) \
465 ((a[0]-b[0])*(a[0]-b[0])+(a[1]-b[1])*(a[1]-b[1])+(a[2]-b[2])*(a[2]-b[2]))
466
482#define DOTP(a, b) (a[0]*b[0]+a[1]*b[1]+a[2]*b[2])
483
502#define FREAD(ptr, type, size, out) { \
503 if(fread(ptr, sizeof(type), size, out)!=size) \
504 ERRMSG("Error while reading!"); \
505 }
506
525#define FWRITE(ptr, type, size, out) { \
526 if(fwrite(ptr, sizeof(type), size, out)!=size) \
527 ERRMSG("Error while writing!"); \
528 }
529
548#define LIN(x0, y0, x1, y1, x) \
549 ((y0)+((y1)-(y0))/((x1)-(x0))*((x)-(x0)))
550
569#define LOGX(x0, y0, x1, y1, x) \
570 (((x)/(x0)>0 && (x1)/(x0)>0) \
571 ? ((y0)+((y1)-(y0))*log((x)/(x0))/log((x1)/(x0))) \
572 : LIN(x0, y0, x1, y1, x))
573
592#define LOGY(x0, y0, x1, y1, x) \
593 (((y1)/(y0)>0) \
594 ? ((y0)*exp(log((y1)/(y0))/((x1)-(x0))*((x)-(x0)))) \
595 : LIN(x0, y0, x1, y1, x))
596
613#define MAX(a,b) (((a)>(b))?(a):(b))
614
631#define MIN(a,b) (((a)<(b))?(a):(b))
632
662#define NEDT(t_bg, nesr, nu) \
663 (BRIGHT(PLANCK((t_bg), (nu)) + (nesr), (nu)) - (t_bg))
664
691#define NESR(t_bg, nedt, nu) \
692 (PLANCK((t_bg) + (nedt), (nu)) - PLANCK((t_bg), (nu)))
693
708#define NORM(a) sqrt(DOTP(a, a))
709
739#define PLANCK(T, nu) \
740 (C1 * POW3(nu) / gsl_expm1(C2 * (nu) / (T)))
741
753#define POW2(x) ((x)*(x))
754
766#define POW3(x) ((x)*(x)*(x))
767
782#define RAD2DEG(rad) ((rad) * (180.0 / M_PI))
783
798#define REFRAC(p, T) (7.753e-05 * (p) / (T))
799
813#define TIMER(name, mode) \
814 {timer(name, __FILE__, __func__, __LINE__, mode);}
815
834#define TOK(line, tok, format, var) { \
835 if(((tok)=strtok((line), " \t"))) { \
836 if(sscanf(tok, format, &(var))!=1) continue; \
837 } else ERRMSG("Error while reading!"); \
838 }
839
840/* ------------------------------------------------------------
841 Log messages...
842 ------------------------------------------------------------ */
843
845#ifndef LOGLEV
846#define LOGLEV 2
847#endif
848
878#define LOG(level, ...) { \
879 if(level >= 2) \
880 printf(" "); \
881 if(level <= LOGLEV) { \
882 printf(__VA_ARGS__); \
883 printf("\n"); \
884 } \
885 }
886
915#define WARN(...) { \
916 printf("\nWarning (%s, %s, l%d): ", __FILE__, __func__, __LINE__); \
917 LOG(0, __VA_ARGS__); \
918 }
919
948#define ERRMSG(...) { \
949 printf("\nError (%s, %s, l%d): ", __FILE__, __func__, __LINE__); \
950 LOG(0, __VA_ARGS__); \
951 exit(EXIT_FAILURE); \
952 }
953
983#define PRINT(format, var) \
984 printf("Print (%s, %s, l%d): %s= "format"\n", \
985 __FILE__, __func__, __LINE__, #var, var);
986
987/* ------------------------------------------------------------
988 Structs...
989 ------------------------------------------------------------ */
990
998typedef struct {
999
1001 int np;
1002
1004 double time[NP];
1005
1007 double z[NP];
1008
1010 double lon[NP];
1011
1013 double lat[NP];
1014
1016 double p[NP];
1017
1019 double t[NP];
1020
1022 double q[NG][NP];
1023
1025 double k[NW][NP];
1026
1028 double clz;
1029
1031 double cldz;
1032
1034 double clk[NCL];
1035
1037 double sft;
1038
1040 double sfeps[NSF];
1041
1042} atm_t;
1043
1051typedef struct {
1052
1054 int ng;
1055
1057 char emitter[NG][LEN];
1058
1061
1064
1067
1070
1072 int nd;
1073
1075 double nu[ND];
1076
1078 int nw;
1079
1081 int window[ND];
1082
1084 int ncl;
1085
1087 double clnu[NCL];
1088
1090 int nsf;
1091
1093 double sfnu[NSF];
1094
1097
1099 double sfsza;
1100
1102 char tblbase[LEN];
1103
1106
1109
1112
1114 double hydz;
1115
1118
1121
1124
1127
1130
1132 double rayds;
1133
1135 double raydz;
1136
1138 char fov[LEN];
1139
1141 double fov_dz[NSHAPE];
1142
1144 double fov_w[NSHAPE];
1145
1148
1151
1154
1157
1160
1162 double retq_zmin[NG];
1163
1165 double retq_zmax[NG];
1166
1168 double retk_zmin[NW];
1169
1171 double retk_zmax[NW];
1172
1175
1178
1181
1184
1187
1190
1193
1196
1198 char rfmbin[LEN];
1199
1201 char rfmhit[LEN];
1202
1204 char rfmxsc[NG][LEN];
1205
1206} ctl_t;
1207
1215typedef struct {
1216
1218 int np;
1219
1221 double z[NLOS];
1222
1224 double lon[NLOS];
1225
1227 double lat[NLOS];
1228
1230 double p[NLOS];
1231
1233 double t[NLOS];
1234
1236 double q[NLOS][NG];
1237
1239 double k[NLOS][ND];
1240
1242 double sft;
1243
1245 double sfeps[ND];
1246
1248 double ds[NLOS];
1249
1251 double u[NLOS][NG];
1252
1254 double cgp[NLOS][NG];
1255
1257 double cgt[NLOS][NG];
1258
1260 double cgu[NLOS][NG];
1261
1263 double eps[NLOS][ND];
1264
1266 double src[NLOS][ND];
1267
1268} los_t;
1269
1277typedef struct {
1278
1280 int nr;
1281
1283 double time[NR];
1284
1286 double obsz[NR];
1287
1289 double obslon[NR];
1290
1292 double obslat[NR];
1293
1295 double vpz[NR];
1296
1298 double vplon[NR];
1299
1301 double vplat[NR];
1302
1304 double tpz[NR];
1305
1307 double tplon[NR];
1308
1310 double tplat[NR];
1311
1313 double tau[ND][NR];
1314
1316 double rad[ND][NR];
1317
1318} obs_t;
1319
1329typedef struct {
1330
1332 char dir[LEN];
1333
1336
1339
1342
1345
1347 double err_formod[ND];
1348
1350 double err_noise[ND];
1351
1354
1357
1360
1362 double err_temp;
1363
1366
1369
1371 double err_q[NG];
1372
1374 double err_q_cz[NG];
1375
1377 double err_q_ch[NG];
1378
1380 double err_k[NW];
1381
1383 double err_k_cz[NW];
1384
1386 double err_k_ch[NW];
1387
1389 double err_clz;
1390
1392 double err_cldz;
1393
1395 double err_clk[NCL];
1396
1398 double err_sft;
1399
1401 double err_sfeps[NSF];
1402
1403} ret_t;
1404
1411typedef struct {
1412
1414 int np[ND][NG];
1415
1417 int nt[ND][NG][TBLNP];
1418
1420 int nu[ND][NG][TBLNP][TBLNT];
1421
1423 double p[ND][NG][TBLNP];
1424
1426 double t[ND][NG][TBLNP][TBLNT];
1427
1429 float u[ND][NG][TBLNP][TBLNT][TBLNU];
1430
1432 float eps[ND][NG][TBLNP][TBLNT][TBLNU];
1433
1435 double st[TBLNS];
1436
1438 double sr[TBLNS][ND];
1439
1440} tbl_t;
1441
1448typedef struct {
1449
1451 double freq;
1452
1454 int64_t offset;
1455
1457 int64_t size;
1458
1460
1468typedef struct {
1469
1471 FILE *fp;
1472
1474 int32_t ntables;
1475
1478
1481
1482} tbl_gas_t;
1483
1484/* ------------------------------------------------------------
1485 Functions...
1486 ------------------------------------------------------------ */
1487
1539void analyze_avk(
1540 const ret_t * ret,
1541 const ctl_t * ctl,
1542 const atm_t * atm,
1543 const int *iqa,
1544 const int *ipa,
1545 const gsl_matrix * avk);
1546
1603 const gsl_matrix * avk,
1604 const int iq,
1605 const int *ipa,
1606 const size_t *n0,
1607 const size_t *n1,
1608 double *cont,
1609 double *res);
1610
1637size_t atm2x(
1638 const ctl_t * ctl,
1639 const atm_t * atm,
1640 gsl_vector * x,
1641 int *iqa,
1642 int *ipa);
1643
1665void atm2x_help(
1666 const double value,
1667 const int value_iqa,
1668 const int value_ip,
1669 gsl_vector * x,
1670 int *iqa,
1671 int *ipa,
1672 size_t *n);
1673
1688void cart2geo(
1689 const double *x,
1690 double *z,
1691 double *lon,
1692 double *lat);
1693
1709void climatology(
1710 const ctl_t * ctl,
1711 atm_t * atm);
1712
1737double cos_sza(
1738 const double sec,
1739 const double lon,
1740 const double lat);
1741
1781double cost_function(
1782 const gsl_vector * dx,
1783 const gsl_vector * dy,
1784 const gsl_matrix * s_a_inv,
1785 const gsl_vector * sig_eps_inv);
1786
1792double ctmco2(
1793 const double nu,
1794 const double p,
1795 const double t,
1796 const double u);
1797
1803double ctmh2o(
1804 const double nu,
1805 const double p,
1806 const double t,
1807 const double q,
1808 const double u);
1809
1839double ctmn2(
1840 const double nu,
1841 const double p,
1842 const double t);
1843
1873double ctmo2(
1874 const double nu,
1875 const double p,
1876 const double t);
1877
1902void copy_atm(
1903 const ctl_t * ctl,
1904 atm_t * atm_dest,
1905 const atm_t * atm_src,
1906 const int init);
1907
1932void copy_obs(
1933 const ctl_t * ctl,
1934 obs_t * obs_dest,
1935 const obs_t * obs_src,
1936 const int init);
1937
1958void day2doy(
1959 int year,
1960 int mon,
1961 int day,
1962 int *doy);
1963
1986void doy2day(
1987 int year,
1988 int doy,
1989 int *mon,
1990 int *day);
1991
2011int find_emitter(
2012 const ctl_t * ctl,
2013 const char *emitter);
2014
2040void formod(
2041 const ctl_t * ctl,
2042 const tbl_t * tbl,
2043 atm_t * atm,
2044 obs_t * obs);
2045
2071void formod_continua(
2072 const ctl_t * ctl,
2073 const los_t * los,
2074 const int ip,
2075 double *beta);
2076
2098void formod_fov(
2099 const ctl_t * ctl,
2100 obs_t * obs);
2101
2131void formod_pencil(
2132 const ctl_t * ctl,
2133 const tbl_t * tbl,
2134 const atm_t * atm,
2135 obs_t * obs,
2136 const int ir);
2137
2176void formod_rfm(
2177 const ctl_t * ctl,
2178 const atm_t * atm,
2179 obs_t * obs);
2180
2205void formod_srcfunc(
2206 const ctl_t * ctl,
2207 const tbl_t * tbl,
2208 const double t,
2209 double *src);
2210
2237void geo2cart(
2238 const double z,
2239 const double lon,
2240 const double lat,
2241 double *x);
2242
2270void hydrostatic(
2271 const ctl_t * ctl,
2272 atm_t * atm);
2273
2297void idx2name(
2298 const ctl_t * ctl,
2299 const int idx,
2300 char *quantity);
2301
2330void init_srcfunc(
2331 const ctl_t * ctl,
2332 tbl_t * tbl);
2333
2358void intpol_atm(
2359 const ctl_t * ctl,
2360 const atm_t * atm,
2361 const double z,
2362 double *p,
2363 double *t,
2364 double *q,
2365 double *k);
2366
2400void intpol_tbl_cga(
2401 const ctl_t * ctl,
2402 const tbl_t * tbl,
2403 const los_t * los,
2404 const int ip,
2405 double tau_path[ND][NG],
2406 double tau_seg[ND]);
2407
2446void intpol_tbl_ega(
2447 const ctl_t * ctl,
2448 const tbl_t * tbl,
2449 const los_t * los,
2450 const int ip,
2451 double tau_path[ND][NG],
2452 double tau_seg[ND]);
2453
2486double intpol_tbl_eps(
2487 const tbl_t * tbl,
2488 const int ig,
2489 const int id,
2490 const int ip,
2491 const int it,
2492 const double u);
2493
2526double intpol_tbl_u(
2527 const tbl_t * tbl,
2528 const int ig,
2529 const int id,
2530 const int ip,
2531 const int it,
2532 const double eps);
2533
2560void jsec2time(
2561 const double jsec,
2562 int *year,
2563 int *mon,
2564 int *day,
2565 int *hour,
2566 int *min,
2567 int *sec,
2568 double *remain);
2569
2611void kernel(
2612 const ctl_t * ctl,
2613 const tbl_t * tbl,
2614 atm_t * atm,
2615 obs_t * obs,
2616 gsl_matrix * k);
2617
2644int locate_irr(
2645 const double *xx,
2646 const int n,
2647 const double x);
2648
2675int locate_reg(
2676 const double *xx,
2677 const int n,
2678 const double x);
2679
2705int locate_tbl(
2706 const float *xx,
2707 const int n,
2708 const double x);
2709
2748void matrix_invert(
2749 gsl_matrix * a);
2750
2795void matrix_product(
2796 const gsl_matrix * a,
2797 const gsl_vector * b,
2798 const int transpose,
2799 gsl_matrix * c);
2800
2828size_t obs2y(
2829 const ctl_t * ctl,
2830 const obs_t * obs,
2831 gsl_vector * y,
2832 int *ida,
2833 int *ira);
2834
2886 ret_t * ret,
2887 ctl_t * ctl,
2888 tbl_t * tbl,
2889 obs_t * obs_meas,
2890 obs_t * obs_i,
2891 atm_t * atm_apr,
2892 atm_t * atm_i,
2893 double *chisq);
2894
2933void raytrace(
2934 const ctl_t * ctl,
2935 const atm_t * atm,
2936 obs_t * obs,
2937 los_t * los,
2938 const int ir);
2939
2985void read_atm(
2986 const char *dirname,
2987 const char *filename,
2988 const ctl_t * ctl,
2989 atm_t * atm);
2990
3037void read_atm_asc(
3038 FILE * in,
3039 const ctl_t * ctl,
3040 atm_t * atm);
3041
3091void read_atm_bin(
3092 FILE * in,
3093 const ctl_t * ctl,
3094 atm_t * atm);
3095
3136void read_ctl(
3137 int argc,
3138 char *argv[],
3139 ctl_t * ctl);
3140
3172void read_matrix(
3173 const char *dirname,
3174 const char *filename,
3175 gsl_matrix * matrix);
3176
3211void read_obs(
3212 const char *dirname,
3213 const char *filename,
3214 const ctl_t * ctl,
3215 obs_t * obs);
3216
3249void read_obs_asc(
3250 FILE * in,
3251 const ctl_t * ctl,
3252 obs_t * obs);
3253
3290void read_obs_bin(
3291 FILE * in,
3292 const ctl_t * ctl,
3293 obs_t * obs);
3294
3334double read_obs_rfm(
3335 const char *basename,
3336 const double z,
3337 const double *nu,
3338 const double *f,
3339 const int n);
3340
3405void read_ret(
3406 int argc,
3407 char *argv[],
3408 const ctl_t * ctl,
3409 ret_t * ret);
3410
3451void read_rfm_spec(
3452 const char *filename,
3453 double *nu,
3454 double *rad,
3455 int *npts);
3456
3490void read_shape(
3491 const char *filename,
3492 double *x,
3493 double *y,
3494 int *n);
3495
3520 const ctl_t * ctl);
3521
3544void read_tbl_asc(
3545 const ctl_t * ctl,
3546 tbl_t * tbl,
3547 const int id,
3548 const int ig);
3549
3575void read_tbl_bin(
3576 const ctl_t * ctl,
3577 tbl_t * tbl,
3578 const int id,
3579 const int ig);
3580
3597void read_tbl_gas(
3598 const ctl_t * ctl,
3599 tbl_t * tbl,
3600 const int id,
3601 const int ig);
3602
3617 tbl_gas_t * g);
3618
3636 const char *path,
3637 tbl_gas_t * g);
3638
3670 const tbl_gas_t * g,
3671 const double freq,
3672 tbl_t * tbl,
3673 const int id,
3674 const int ig);
3675
3730double scan_ctl(
3731 int argc,
3732 char *argv[],
3733 const char *varname,
3734 const int arridx,
3735 const char *defvalue,
3736 char *value);
3737
3791void set_cov_apr(
3792 const ret_t * ret,
3793 const ctl_t * ctl,
3794 const atm_t * atm,
3795 const int *iqa,
3796 const int *ipa,
3797 gsl_matrix * s_a);
3798
3852void set_cov_meas(
3853 const ret_t * ret,
3854 const ctl_t * ctl,
3855 const obs_t * obs,
3856 gsl_vector * sig_noise,
3857 gsl_vector * sig_formod,
3858 gsl_vector * sig_eps_inv);
3859
3895double sza(
3896 double sec,
3897 double lon,
3898 double lat);
3899
3939void tangent_point(
3940 const los_t * los,
3941 double *tpz,
3942 double *tplon,
3943 double *tplat);
3944
3969void time2jsec(
3970 const int year,
3971 const int mon,
3972 const int day,
3973 const int hour,
3974 const int min,
3975 const int sec,
3976 const double remain,
3977 double *jsec);
3978
4020void timer(
4021 const char *name,
4022 const char *file,
4023 const char *func,
4024 int line,
4025 int mode);
4026
4071void write_atm(
4072 const char *dirname,
4073 const char *filename,
4074 const ctl_t * ctl,
4075 const atm_t * atm);
4076
4128void write_atm_asc(
4129 FILE * out,
4130 const ctl_t * ctl,
4131 const atm_t * atm);
4132
4187void write_atm_bin(
4188 FILE * out,
4189 const ctl_t * ctl,
4190 const atm_t * atm);
4191
4243void write_atm_rfm(
4244 const char *filename,
4245 const ctl_t * ctl,
4246 const atm_t * atm);
4247
4295void write_matrix(
4296 const char *dirname,
4297 const char *filename,
4298 const ctl_t * ctl,
4299 const gsl_matrix * matrix,
4300 const atm_t * atm,
4301 const obs_t * obs,
4302 const char *rowspace,
4303 const char *colspace,
4304 const char *sort);
4305
4338void write_obs(
4339 const char *dirname,
4340 const char *filename,
4341 const ctl_t * ctl,
4342 const obs_t * obs);
4343
4378void write_obs_asc(
4379 FILE * out,
4380 const ctl_t * ctl,
4381 const obs_t * obs);
4382
4414void write_obs_bin(
4415 FILE * out,
4416 const ctl_t * ctl,
4417 const obs_t * obs);
4418
4453void write_shape(
4454 const char *filename,
4455 const double *x,
4456 const double *y,
4457 const int n);
4458
4504void write_stddev(
4505 const char *quantity,
4506 const ret_t * ret,
4507 const ctl_t * ctl,
4508 const atm_t * atm,
4509 const gsl_matrix * s);
4510
4529void write_tbl(
4530 const ctl_t * ctl,
4531 const tbl_t * tbl);
4532
4556void write_tbl_asc(
4557 const ctl_t * ctl,
4558 const tbl_t * tbl);
4559
4584void write_tbl_bin(
4585 const ctl_t * ctl,
4586 const tbl_t * tbl);
4587
4614void write_tbl_gas(
4615 const ctl_t * ctl,
4616 const tbl_t * tbl);
4617
4637 const char *path);
4638
4676 tbl_gas_t * g,
4677 const double freq,
4678 const tbl_t * tbl,
4679 const int id,
4680 const int ig);
4681
4728void x2atm(
4729 const ctl_t * ctl,
4730 const gsl_vector * x,
4731 atm_t * atm);
4732
4761void x2atm_help(
4762 double *value,
4763 const gsl_vector * x,
4764 size_t *n);
4765
4798void y2obs(
4799 const ctl_t * ctl,
4800 const gsl_vector * y,
4801 obs_t * obs);
4802
4803#endif
void analyze_avk(const ret_t *ret, const ctl_t *ctl, const atm_t *atm, const int *iqa, const int *ipa, const gsl_matrix *avk)
Analyze averaging kernel (AVK) matrix for retrieval diagnostics.
Definition: jurassic.c:29
void read_tbl_bin(const ctl_t *ctl, tbl_t *tbl, const int id, const int ig)
Read a single binary emissivity lookup table.
Definition: jurassic.c:6037
#define LEN
Maximum length of ASCII data lines.
Definition: jurassic.h:267
double cost_function(const gsl_vector *dx, const gsl_vector *dy, const gsl_matrix *s_a_inv, const gsl_vector *sig_eps_inv)
Compute the normalized quadratic cost function for optimal estimation.
Definition: jurassic.c:1179
void read_atm_bin(FILE *in, const ctl_t *ctl, atm_t *atm)
Read atmospheric data in binary format.
Definition: jurassic.c:5236
double read_obs_rfm(const char *basename, const double z, const double *nu, const double *f, const int n)
Read and spectrally convolve an RFM output spectrum.
Definition: jurassic.c:5673
void read_matrix(const char *dirname, const char *filename, gsl_matrix *matrix)
Read a numerical matrix from an ASCII file.
Definition: jurassic.c:5440
void timer(const char *name, const char *file, const char *func, int line, int mode)
Simple wall-clock timer for runtime diagnostics.
Definition: jurassic.c:6566
void read_rfm_spec(const char *filename, double *nu, double *rad, int *npts)
Read a Reference Forward Model (RFM) ASCII spectrum.
Definition: jurassic.c:5784
double cos_sza(const double sec, const double lon, const double lat)
Calculates the cosine of the solar zenith angle.
Definition: jurassic.c:1138
void write_atm(const char *dirname, const char *filename, const ctl_t *ctl, const atm_t *atm)
Write atmospheric data to a file.
Definition: jurassic.c:6604
int locate_reg(const double *xx, const int n, const double x)
Locate index for interpolation on a regular (uniform) grid.
Definition: jurassic.c:4458
void intpol_tbl_ega(const ctl_t *ctl, const tbl_t *tbl, const los_t *los, const int ip, double tau_path[ND][NG], double tau_seg[ND])
Interpolate emissivities and transmittances using the Emissivity Growth Approximation (EGA).
Definition: jurassic.c:4140
void formod_rfm(const ctl_t *ctl, const atm_t *atm, obs_t *obs)
Interface routine for the Reference Forward Model (RFM).
Definition: jurassic.c:3695
double ctmo2(const double nu, const double p, const double t)
Compute O₂ collision-induced absorption coefficient.
Definition: jurassic.c:3196
int write_tbl_gas_single(tbl_gas_t *g, const double freq, const tbl_t *tbl, const int id, const int ig)
Append or overwrite a single frequency-table block in a per-gas file.
Definition: jurassic.c:7480
void idx2name(const ctl_t *ctl, const int idx, char *quantity)
Convert a quantity index to a descriptive name string.
Definition: jurassic.c:3933
int read_tbl_gas_single(const tbl_gas_t *g, const double freq, tbl_t *tbl, const int id, const int ig)
Read one emissivity table block from a per-gas table file.
Definition: jurassic.c:6189
void read_ctl(int argc, char *argv[], ctl_t *ctl)
Read model control parameters from command-line and configuration input.
Definition: jurassic.c:5322
void formod_continua(const ctl_t *ctl, const los_t *los, const int ip, double *beta)
Compute total extinction including gaseous continua.
Definition: jurassic.c:3459
void write_obs_bin(FILE *out, const ctl_t *ctl, const obs_t *obs)
Write observation data in binary format to an output file stream.
Definition: jurassic.c:7163
void raytrace(const ctl_t *ctl, const atm_t *atm, obs_t *obs, los_t *los, const int ir)
Perform line-of-sight (LOS) ray tracing through the atmosphere.
Definition: jurassic.c:4873
void matrix_product(const gsl_matrix *a, const gsl_vector *b, const int transpose, gsl_matrix *c)
Compute structured matrix products of the form or .
Definition: jurassic.c:4529
int locate_irr(const double *xx, const int n, const double x)
Locate index for interpolation on an irregular grid.
Definition: jurassic.c:4428
void write_tbl_asc(const ctl_t *ctl, const tbl_t *tbl)
Write all lookup tables in human-readable ASCII format.
Definition: jurassic.c:7309
void day2doy(int year, int mon, int day, int *doy)
Convert a calendar date to day-of-year.
Definition: jurassic.c:3346
void time2jsec(const int year, const int mon, const int day, const int hour, const int min, const int sec, const double remain, double *jsec)
Converts time components to seconds since January 1, 2000, 12:00:00 UTC.
Definition: jurassic.c:6535
void x2atm(const ctl_t *ctl, const gsl_vector *x, atm_t *atm)
Map retrieval state vector back to atmospheric structure.
Definition: jurassic.c:7577
void atm2x_help(const double value, const int value_iqa, const int value_ip, gsl_vector *x, int *iqa, int *ipa, size_t *n)
Append a single atmospheric value to the state vector.
Definition: jurassic.c:164
double ctmh2o(const double nu, const double p, const double t, const double q, const double u)
Compute water vapor continuum (optical depth).
Definition: jurassic.c:2075
void intpol_tbl_cga(const ctl_t *ctl, const tbl_t *tbl, const los_t *los, const int ip, double tau_path[ND][NG], double tau_seg[ND])
Interpolate emissivities and transmittances using the Curtis–Godson approximation (CGA).
Definition: jurassic.c:4051
void write_tbl_bin(const ctl_t *ctl, const tbl_t *tbl)
Write all lookup tables in compact binary format.
Definition: jurassic.c:7354
void write_atm_rfm(const char *filename, const ctl_t *ctl, const atm_t *atm)
Write atmospheric profile in RFM-compatible format.
Definition: jurassic.c:6817
#define ND
Maximum number of radiance channels.
Definition: jurassic.h:237
void set_cov_meas(const ret_t *ret, const ctl_t *ctl, const obs_t *obs, gsl_vector *sig_noise, gsl_vector *sig_formod, gsl_vector *sig_eps_inv)
Construct measurement error standard deviations and their inverse.
Definition: jurassic.c:6447
double intpol_tbl_u(const tbl_t *tbl, const int ig, const int id, const int ip, const int it, const double eps)
Interpolate column density from lookup tables as a function of emissivity.
Definition: jurassic.c:4268
double intpol_tbl_eps(const tbl_t *tbl, const int ig, const int id, const int ip, const int it, const double u)
Interpolate emissivity from lookup tables as a function of column density.
Definition: jurassic.c:4236
#define NSHAPE
Maximum number of shape function grid points.
Definition: jurassic.h:292
void write_tbl(const ctl_t *ctl, const tbl_t *tbl)
Write all emissivity lookup tables in the format specified by the control structure.
Definition: jurassic.c:7286
void write_obs_asc(FILE *out, const ctl_t *ctl, const obs_t *obs)
Write observation data to an ASCII text file.
Definition: jurassic.c:7115
void write_tbl_gas(const ctl_t *ctl, const tbl_t *tbl)
Write lookup tables into per-gas binary table files with indexed blocks.
Definition: jurassic.c:7409
void intpol_atm(const ctl_t *ctl, const atm_t *atm, const double z, double *p, double *t, double *q, double *k)
Interpolate atmospheric state variables at a given altitude.
Definition: jurassic.c:4026
void read_tbl_gas(const ctl_t *ctl, tbl_t *tbl, const int id, const int ig)
Read one frequency block from a per-gas binary table file.
Definition: jurassic.c:6097
int find_emitter(const ctl_t *ctl, const char *emitter)
Find gas species index by name.
Definition: jurassic.c:3396
void optimal_estimation(ret_t *ret, ctl_t *ctl, tbl_t *tbl, obs_t *obs_meas, obs_t *obs_i, atm_t *atm_apr, atm_t *atm_i, double *chisq)
Perform optimal estimation retrieval using Levenberg–Marquardt minimization.
Definition: jurassic.c:4601
void tangent_point(const los_t *los, double *tpz, double *tplon, double *tplat)
Determine the tangent point along a line of sight (LOS).
Definition: jurassic.c:6491
#define TBLNU
Maximum number of column densities in emissivity tables.
Definition: jurassic.h:312
void read_obs_asc(FILE *in, const ctl_t *ctl, obs_t *obs)
Read ASCII-formatted observation data from an open file stream.
Definition: jurassic.c:5566
void formod_pencil(const ctl_t *ctl, const tbl_t *tbl, const atm_t *atm, obs_t *obs, const int ir)
Compute line-of-sight radiances using the pencil-beam forward model.
Definition: jurassic.c:3560
double ctmco2(const double nu, const double p, const double t, const double u)
Compute carbon dioxide continuum (optical depth).
Definition: jurassic.c:1212
void read_obs_bin(FILE *in, const ctl_t *ctl, obs_t *obs)
Read binary-formatted observation data from an open file stream.
Definition: jurassic.c:5603
void read_atm(const char *dirname, const char *filename, const ctl_t *ctl, atm_t *atm)
Read atmospheric input data from a file.
Definition: jurassic.c:5109
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 data from an input file.
Definition: jurassic.c:5480
void write_obs(const char *dirname, const char *filename, const ctl_t *ctl, const obs_t *obs)
Write observation data to an output file in ASCII or binary format.
Definition: jurassic.c:7033
#define TBLNT
Maximum number of temperatures in emissivity tables.
Definition: jurassic.h:307
void x2atm_help(double *value, const gsl_vector *x, size_t *n)
Helper function to extract a single value from the retrieval state vector.
Definition: jurassic.c:7627
void cart2geo(const double *x, double *z, double *lon, double *lat)
Converts Cartesian coordinates to geographic coordinates.
Definition: jurassic.c:185
#define NP
Maximum number of atmospheric data points.
Definition: jurassic.h:247
void formod_fov(const ctl_t *ctl, obs_t *obs)
Apply field-of-view (FOV) convolution to modeled radiances.
Definition: jurassic.c:3495
void doy2day(int year, int doy, int *mon, int *day)
Convert a day-of-year value to a calendar date.
Definition: jurassic.c:3366
void kernel(const ctl_t *ctl, const tbl_t *tbl, atm_t *atm, obs_t *obs, gsl_matrix *k)
Compute the Jacobian (kernel) matrix by finite differences.
Definition: jurassic.c:4333
void init_srcfunc(const ctl_t *ctl, tbl_t *tbl)
Initialize the source-function (Planck radiance) lookup table.
Definition: jurassic.c:3972
void matrix_invert(gsl_matrix *a)
Invert a square matrix, optimized for diagonal or symmetric positive-definite matrices.
Definition: jurassic.c:4499
void read_ret(int argc, char *argv[], const ctl_t *ctl, ret_t *ret)
Read retrieval configuration and error parameters.
Definition: jurassic.c:5731
void write_atm_bin(FILE *out, const ctl_t *ctl, const atm_t *atm)
Write atmospheric data to a binary file.
Definition: jurassic.c:6741
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:6267
void write_shape(const char *filename, const double *x, const double *y, const int n)
Write tabulated shape function data to a text file.
Definition: jurassic.c:7223
void formod(const ctl_t *ctl, const tbl_t *tbl, atm_t *atm, obs_t *obs)
Execute the selected forward model.
Definition: jurassic.c:3409
#define NG
Maximum number of emitters.
Definition: jurassic.h:242
void read_tbl_asc(const ctl_t *ctl, tbl_t *tbl, const int id, const int ig)
Read a single ASCII emissivity lookup table.
Definition: jurassic.c:5939
void formod_srcfunc(const ctl_t *ctl, const tbl_t *tbl, const double t, double *src)
Interpolate the source function (Planck radiance) at a given temperature.
Definition: jurassic.c:3836
void analyze_avk_quantity(const gsl_matrix *avk, const int iq, const int *ipa, const size_t *n0, const size_t *n1, double *cont, double *res)
Analyze averaging kernel submatrix for a specific retrieved quantity.
Definition: jurassic.c:86
void jsec2time(const double jsec, int *year, int *mon, int *day, int *hour, int *min, int *sec, double *remain)
Converts Julian seconds to calendar date and time components.
Definition: jurassic.c:4300
#define TBLNS
Maximum number of source function temperature levels.
Definition: jurassic.h:317
tbl_t * read_tbl(const ctl_t *ctl)
Read all emissivity lookup tables for all gases and frequencies.
Definition: jurassic.c:5887
double sza(double sec, double lon, double lat)
Compute the solar zenith angle for a given time and location.
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
size_t obs2y(const ctl_t *ctl, const obs_t *obs, gsl_vector *y, int *ida, int *ira)
Convert observation radiances into a measurement vector.
Definition: jurassic.c:4574
void y2obs(const ctl_t *ctl, const gsl_vector *y, obs_t *obs)
Copy elements from the measurement vector y into the observation structure.
Definition: jurassic.c:7639
void hydrostatic(const ctl_t *ctl, atm_t *atm)
Adjust pressure profile using the hydrostatic equation.
Definition: jurassic.c:3873
void read_shape(const char *filename, double *x, double *y, int *n)
Read a two-column shape function from an ASCII file.
Definition: jurassic.c:5845
void read_atm_asc(FILE *in, const ctl_t *ctl, atm_t *atm)
Read atmospheric data in ASCII format.
Definition: jurassic.c:5192
double ctmn2(const double nu, const double p, const double t)
Compute N₂ collision-induced absorption coefficient.
Definition: jurassic.c:3127
int read_tbl_gas_close(tbl_gas_t *g)
Close a per-gas binary table file and optionally rewrite metadata.
Definition: jurassic.c:6127
size_t atm2x(const ctl_t *ctl, const atm_t *atm, gsl_vector *x, int *iqa, int *ipa)
Convert atmospheric data to state vector elements.
Definition: jurassic.c:110
int read_tbl_gas_open(const char *path, tbl_gas_t *g)
Open a per-gas binary table file for reading and writing.
Definition: jurassic.c:6160
#define TBLNP
Maximum number of pressure levels in emissivity tables.
Definition: jurassic.h:302
void climatology(const ctl_t *ctl, atm_t *atm)
Initializes atmospheric climatology profiles.
Definition: jurassic.c:200
void geo2cart(const double z, const double lon, const double lat, double *x)
Converts geographic coordinates (longitude, latitude, altitude) to Cartesian coordinates.
Definition: jurassic.c:3853
#define NSF
Maximum number of surface layer spectral grid points.
Definition: jurassic.h:257
void write_matrix(const char *dirname, const char *filename, const ctl_t *ctl, const gsl_matrix *matrix, const atm_t *atm, const obs_t *obs, const char *rowspace, const char *colspace, const char *sort)
Write a fully annotated matrix (e.g., Jacobian or gain matrix) to file.
Definition: jurassic.c:6855
#define NCL
Maximum number of cloud layer spectral grid points.
Definition: jurassic.h:232
void set_cov_apr(const ret_t *ret, const ctl_t *ctl, const atm_t *atm, const int *iqa, const int *ipa, gsl_matrix *s_a)
Construct the a priori covariance matrix for retrieval parameters.
Definition: jurassic.c:6336
void write_stddev(const char *quantity, const ret_t *ret, const ctl_t *ctl, const atm_t *atm, const gsl_matrix *s)
Write retrieval standard deviation profiles to disk.
Definition: jurassic.c:7253
#define NLOS
Maximum number of LOS points.
Definition: jurassic.h:287
void write_atm_asc(FILE *out, const ctl_t *ctl, const atm_t *atm)
Write atmospheric data to an ASCII file.
Definition: jurassic.c:6680
#define NR
Maximum number of ray paths.
Definition: jurassic.h:252
int locate_tbl(const float *xx, const int n, const double x)
Locate index for interpolation within emissivity table grids.
Definition: jurassic.c:4477
#define NW
Maximum number of spectral windows.
Definition: jurassic.h:262
int write_tbl_gas_create(const char *path)
Create a new per-gas table file with an empty index.
Definition: jurassic.c:7449
Atmospheric profile data.
Definition: jurassic.h:998
double clz
Cloud layer height [km].
Definition: jurassic.h:1028
int np
Number of data points.
Definition: jurassic.h:1001
double cldz
Cloud layer depth [km].
Definition: jurassic.h:1031
double sft
Surface temperature [K].
Definition: jurassic.h:1037
Control parameters.
Definition: jurassic.h:1051
int write_matrix
Write matrix file (0=no, 1=yes).
Definition: jurassic.h:1192
int nw
Number of spectral windows.
Definition: jurassic.h:1078
int atmfmt
Atmospheric data file format (1=ASCII, 2=binary).
Definition: jurassic.h:1108
double retp_zmin
Minimum altitude for pressure retrieval [km].
Definition: jurassic.h:1150
int ig_co2
Emitter index of CO2.
Definition: jurassic.h:1060
double hydz
Reference height for hydrostatic pressure profile (-999 to skip) [km].
Definition: jurassic.h:1114
int ctm_co2
Compute CO2 continuum (0=no, 1=yes).
Definition: jurassic.h:1117
double rett_zmax
Maximum altitude for temperature retrieval [km].
Definition: jurassic.h:1159
int ret_sfeps
Retrieve surface layer emissivity (0=no, 1=yes).
Definition: jurassic.h:1186
int ret_sft
Retrieve surface layer temperature (0=no, 1=yes).
Definition: jurassic.h:1183
int ret_clz
Retrieve cloud layer height (0=no, 1=yes).
Definition: jurassic.h:1174
int ctm_n2
Compute N2 continuum (0=no, 1=yes).
Definition: jurassic.h:1123
int ctm_h2o
Compute H2O continuum (0=no, 1=yes).
Definition: jurassic.h:1120
int formod
Forward model (0=CGA, 1=EGA, 2=RFM).
Definition: jurassic.h:1195
int ng
Number of emitters.
Definition: jurassic.h:1054
int refrac
Take into account refractivity (0=no, 1=yes).
Definition: jurassic.h:1129
int ig_o2
Emitter index of O2.
Definition: jurassic.h:1069
double rett_zmin
Minimum altitude for temperature retrieval [km].
Definition: jurassic.h:1156
double sfsza
Solar zenith angle at the surface [deg] (-999=auto).
Definition: jurassic.h:1099
int nd
Number of radiance channels.
Definition: jurassic.h:1072
int fov_n
Field-of-view number of data points.
Definition: jurassic.h:1147
int sftype
Surface treatment (0=none, 1=emissions, 2=downward, 3=solar).
Definition: jurassic.h:1096
int ncl
Number of cloud layer spectral grid points.
Definition: jurassic.h:1084
int ig_n2
Emitter index of N2.
Definition: jurassic.h:1066
int obsfmt
Observation data file format (1=ASCII, 2=binary).
Definition: jurassic.h:1111
int ctm_o2
Compute O2 continuum (0=no, 1=yes).
Definition: jurassic.h:1126
int ret_clk
Retrieve cloud layer extinction (0=no, 1=yes).
Definition: jurassic.h:1180
int nsf
Number of surface layer spectral grid points.
Definition: jurassic.h:1090
double rayds
Maximum step length for raytracing [km].
Definition: jurassic.h:1132
int ret_cldz
Retrieve cloud layer depth (0=no, 1=yes).
Definition: jurassic.h:1177
int ig_h2o
Emitter index of H2O.
Definition: jurassic.h:1063
int tblfmt
Look-up table file format (1=ASCII, 2=binary).
Definition: jurassic.h:1105
double raydz
Vertical step length for raytracing [km].
Definition: jurassic.h:1135
int write_bbt
Use brightness temperature instead of radiance (0=no, 1=yes).
Definition: jurassic.h:1189
double retp_zmax
Maximum altitude for pressure retrieval [km].
Definition: jurassic.h:1153
Line-of-sight data.
Definition: jurassic.h:1215
double sft
Surface temperature [K].
Definition: jurassic.h:1242
int np
Number of LOS points.
Definition: jurassic.h:1218
Observation geometry and radiance data.
Definition: jurassic.h:1277
int nr
Number of ray paths.
Definition: jurassic.h:1280
Retrieval control parameters.
Definition: jurassic.h:1329
double err_press_cz
Vertical correlation length for pressure error [km].
Definition: jurassic.h:1356
double err_press
Pressure error [%].
Definition: jurassic.h:1353
int err_ana
Carry out error analysis (0=no, 1=yes).
Definition: jurassic.h:1344
double err_temp_cz
Vertical correlation length for temperature error [km].
Definition: jurassic.h:1365
double conv_dmin
Minimum normalized step size in state space.
Definition: jurassic.h:1341
double err_temp
Temperature error [K].
Definition: jurassic.h:1362
double err_clz
Cloud height error [km].
Definition: jurassic.h:1389
double err_sft
Surface temperature error [K].
Definition: jurassic.h:1398
double err_temp_ch
Horizontal correlation length for temperature error [km].
Definition: jurassic.h:1368
int kernel_recomp
Re-computation of kernel matrix (number of iterations).
Definition: jurassic.h:1335
int conv_itmax
Maximum number of iterations.
Definition: jurassic.h:1338
double err_cldz
Cloud depth error [km].
Definition: jurassic.h:1392
double err_press_ch
Horizontal correlation length for pressure error [km].
Definition: jurassic.h:1359
On-disk index entry describing one frequency table block in a gas file.
Definition: jurassic.h:1448
int64_t offset
Byte offset in file where the serialized block begins.
Definition: jurassic.h:1454
int64_t size
Size of the serialized block (in bytes).
Definition: jurassic.h:1457
double freq
Frequency identifier ν_j for this table block.
Definition: jurassic.h:1451
In-memory representation of an open per-gas lookup-table file.
Definition: jurassic.h:1468
FILE * fp
Open file handle ("rb+"), NULL if not open.
Definition: jurassic.h:1471
int dirty
Definition: jurassic.h:1480
int32_t ntables
Number of index entries currently in use.
Definition: jurassic.h:1474
tbl_gas_index_t * index
In-memory index table of length MAX_TABLES.
Definition: jurassic.h:1477
Emissivity look-up tables.
Definition: jurassic.h:1411