JURASSIC
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
1108 double hydz;
1109
1112
1115
1118
1121
1124
1126 double rayds;
1127
1129 double raydz;
1130
1132 char fov[LEN];
1133
1135 double fov_dz[NSHAPE];
1136
1138 double fov_w[NSHAPE];
1139
1142
1145
1148
1151
1154
1156 double retq_zmin[NG];
1157
1159 double retq_zmax[NG];
1160
1162 double retk_zmin[NW];
1163
1165 double retk_zmax[NW];
1166
1169
1172
1175
1178
1181
1184
1187
1190
1192 char rfmbin[LEN];
1193
1195 char rfmhit[LEN];
1196
1198 char rfmxsc[NG][LEN];
1199
1200} ctl_t;
1201
1209typedef struct {
1210
1212 int np;
1213
1215 double z[NLOS];
1216
1218 double lon[NLOS];
1219
1221 double lat[NLOS];
1222
1224 double p[NLOS];
1225
1227 double t[NLOS];
1228
1230 double q[NLOS][NG];
1231
1233 double k[NLOS][ND];
1234
1236 double sft;
1237
1239 double sfeps[ND];
1240
1242 double ds[NLOS];
1243
1245 double u[NLOS][NG];
1246
1248 double cgp[NLOS][NG];
1249
1251 double cgt[NLOS][NG];
1252
1254 double cgu[NLOS][NG];
1255
1257 double eps[NLOS][ND];
1258
1260 double src[NLOS][ND];
1261
1262} los_t;
1263
1271typedef struct {
1272
1274 int nr;
1275
1277 double time[NR];
1278
1280 double obsz[NR];
1281
1283 double obslon[NR];
1284
1286 double obslat[NR];
1287
1289 double vpz[NR];
1290
1292 double vplon[NR];
1293
1295 double vplat[NR];
1296
1298 double tpz[NR];
1299
1301 double tplon[NR];
1302
1304 double tplat[NR];
1305
1307 double tau[ND][NR];
1308
1310 double rad[ND][NR];
1311
1312} obs_t;
1313
1323typedef struct {
1324
1326 char dir[LEN];
1327
1330
1333
1336
1339
1341 double err_formod[ND];
1342
1344 double err_noise[ND];
1345
1348
1351
1354
1356 double err_temp;
1357
1360
1363
1365 double err_q[NG];
1366
1368 double err_q_cz[NG];
1369
1371 double err_q_ch[NG];
1372
1374 double err_k[NW];
1375
1377 double err_k_cz[NW];
1378
1380 double err_k_ch[NW];
1381
1383 double err_clz;
1384
1386 double err_cldz;
1387
1389 double err_clk[NCL];
1390
1392 double err_sft;
1393
1395 double err_sfeps[NSF];
1396
1397} ret_t;
1398
1405typedef struct {
1406
1408 int np[ND][NG];
1409
1411 int nt[ND][NG][TBLNP];
1412
1414 int nu[ND][NG][TBLNP][TBLNT];
1415
1417 double p[ND][NG][TBLNP];
1418
1420 double t[ND][NG][TBLNP][TBLNT];
1421
1423 float u[ND][NG][TBLNP][TBLNT][TBLNU];
1424
1426 float eps[ND][NG][TBLNP][TBLNT][TBLNU];
1427
1429 double st[TBLNS];
1430
1432 double sr[TBLNS][ND];
1433
1434} tbl_t;
1435
1442typedef struct {
1443
1445 double freq;
1446
1448 int64_t offset;
1449
1451 int64_t size;
1452
1454
1462typedef struct {
1463
1465 FILE *fp;
1466
1468 int32_t ntables;
1469
1472
1475
1476} tbl_gas_t;
1477
1478/* ------------------------------------------------------------
1479 Functions...
1480 ------------------------------------------------------------ */
1481
1533void analyze_avk(
1534 const ret_t * ret,
1535 const ctl_t * ctl,
1536 const atm_t * atm,
1537 const int *iqa,
1538 const int *ipa,
1539 const gsl_matrix * avk);
1540
1597 const gsl_matrix * avk,
1598 const int iq,
1599 const int *ipa,
1600 const size_t *n0,
1601 const size_t *n1,
1602 double *cont,
1603 double *res);
1604
1631size_t atm2x(
1632 const ctl_t * ctl,
1633 const atm_t * atm,
1634 gsl_vector * x,
1635 int *iqa,
1636 int *ipa);
1637
1659void atm2x_help(
1660 const double value,
1661 const int value_iqa,
1662 const int value_ip,
1663 gsl_vector * x,
1664 int *iqa,
1665 int *ipa,
1666 size_t *n);
1667
1682void cart2geo(
1683 const double *x,
1684 double *z,
1685 double *lon,
1686 double *lat);
1687
1703void climatology(
1704 const ctl_t * ctl,
1705 atm_t * atm);
1706
1731double cos_sza(
1732 const double sec,
1733 const double lon,
1734 const double lat);
1735
1775double cost_function(
1776 const gsl_vector * dx,
1777 const gsl_vector * dy,
1778 const gsl_matrix * s_a_inv,
1779 const gsl_vector * sig_eps_inv);
1780
1786double ctmco2(
1787 const double nu,
1788 const double p,
1789 const double t,
1790 const double u);
1791
1797double ctmh2o(
1798 const double nu,
1799 const double p,
1800 const double t,
1801 const double q,
1802 const double u);
1803
1833double ctmn2(
1834 const double nu,
1835 const double p,
1836 const double t);
1837
1867double ctmo2(
1868 const double nu,
1869 const double p,
1870 const double t);
1871
1896void copy_atm(
1897 const ctl_t * ctl,
1898 atm_t * atm_dest,
1899 const atm_t * atm_src,
1900 const int init);
1901
1926void copy_obs(
1927 const ctl_t * ctl,
1928 obs_t * obs_dest,
1929 const obs_t * obs_src,
1930 const int init);
1931
1951int find_emitter(
1952 const ctl_t * ctl,
1953 const char *emitter);
1954
1980void formod(
1981 const ctl_t * ctl,
1982 const tbl_t * tbl,
1983 atm_t * atm,
1984 obs_t * obs);
1985
2011void formod_continua(
2012 const ctl_t * ctl,
2013 const los_t * los,
2014 const int ip,
2015 double *beta);
2016
2038void formod_fov(
2039 const ctl_t * ctl,
2040 obs_t * obs);
2041
2071void formod_pencil(
2072 const ctl_t * ctl,
2073 const tbl_t * tbl,
2074 const atm_t * atm,
2075 obs_t * obs,
2076 const int ir);
2077
2116void formod_rfm(
2117 const ctl_t * ctl,
2118 const atm_t * atm,
2119 obs_t * obs);
2120
2145void formod_srcfunc(
2146 const ctl_t * ctl,
2147 const tbl_t * tbl,
2148 const double t,
2149 double *src);
2150
2177void geo2cart(
2178 const double z,
2179 const double lon,
2180 const double lat,
2181 double *x);
2182
2210void hydrostatic(
2211 const ctl_t * ctl,
2212 atm_t * atm);
2213
2237void idx2name(
2238 const ctl_t * ctl,
2239 const int idx,
2240 char *quantity);
2241
2270void init_srcfunc(
2271 const ctl_t * ctl,
2272 tbl_t * tbl);
2273
2298void intpol_atm(
2299 const ctl_t * ctl,
2300 const atm_t * atm,
2301 const double z,
2302 double *p,
2303 double *t,
2304 double *q,
2305 double *k);
2306
2340void intpol_tbl_cga(
2341 const ctl_t * ctl,
2342 const tbl_t * tbl,
2343 const los_t * los,
2344 const int ip,
2345 double tau_path[ND][NG],
2346 double tau_seg[ND]);
2347
2386void intpol_tbl_ega(
2387 const ctl_t * ctl,
2388 const tbl_t * tbl,
2389 const los_t * los,
2390 const int ip,
2391 double tau_path[ND][NG],
2392 double tau_seg[ND]);
2393
2426double intpol_tbl_eps(
2427 const tbl_t * tbl,
2428 const int ig,
2429 const int id,
2430 const int ip,
2431 const int it,
2432 const double u);
2433
2466double intpol_tbl_u(
2467 const tbl_t * tbl,
2468 const int ig,
2469 const int id,
2470 const int ip,
2471 const int it,
2472 const double eps);
2473
2500void jsec2time(
2501 const double jsec,
2502 int *year,
2503 int *mon,
2504 int *day,
2505 int *hour,
2506 int *min,
2507 int *sec,
2508 double *remain);
2509
2551void kernel(
2552 const ctl_t * ctl,
2553 const tbl_t * tbl,
2554 atm_t * atm,
2555 obs_t * obs,
2556 gsl_matrix * k);
2557
2584int locate_irr(
2585 const double *xx,
2586 const int n,
2587 const double x);
2588
2615int locate_reg(
2616 const double *xx,
2617 const int n,
2618 const double x);
2619
2645int locate_tbl(
2646 const float *xx,
2647 const int n,
2648 const double x);
2649
2688void matrix_invert(
2689 gsl_matrix * a);
2690
2735void matrix_product(
2736 const gsl_matrix * a,
2737 const gsl_vector * b,
2738 const int transpose,
2739 gsl_matrix * c);
2740
2768size_t obs2y(
2769 const ctl_t * ctl,
2770 const obs_t * obs,
2771 gsl_vector * y,
2772 int *ida,
2773 int *ira);
2774
2813void raytrace(
2814 const ctl_t * ctl,
2815 const atm_t * atm,
2816 obs_t * obs,
2817 los_t * los,
2818 const int ir);
2819
2857void read_atm(
2858 const char *dirname,
2859 const char *filename,
2860 const ctl_t * ctl,
2861 atm_t * atm);
2862
2903void read_ctl(
2904 int argc,
2905 char *argv[],
2906 ctl_t * ctl);
2907
2939void read_matrix(
2940 const char *dirname,
2941 const char *filename,
2942 gsl_matrix * matrix);
2943
2980void read_obs(
2981 const char *dirname,
2982 const char *filename,
2983 const ctl_t * ctl,
2984 obs_t * obs);
2985
3025double read_obs_rfm(
3026 const char *basename,
3027 const double z,
3028 const double *nu,
3029 const double *f,
3030 const int n);
3031
3096void read_ret(
3097 int argc,
3098 char *argv[],
3099 const ctl_t * ctl,
3100 ret_t * ret);
3101
3142void read_rfm_spec(
3143 const char *filename,
3144 double *nu,
3145 double *rad,
3146 int *npts);
3147
3181void read_shape(
3182 const char *filename,
3183 double *x,
3184 double *y,
3185 int *n);
3186
3211 const ctl_t * ctl);
3212
3236void read_tbl_asc(
3237 const ctl_t * ctl,
3238 tbl_t * tbl,
3239 const int id,
3240 const int ig);
3241
3268void read_tbl_bin(
3269 const ctl_t * ctl,
3270 tbl_t * tbl,
3271 const int id,
3272 const int ig);
3273
3291void read_tbl_gas(
3292 const ctl_t * ctl,
3293 tbl_t * tbl,
3294 const int id,
3295 const int ig);
3296
3312 tbl_gas_t * g);
3313
3332 const char *path,
3333 tbl_gas_t * g);
3334
3367 const tbl_gas_t * g,
3368 const double freq,
3369 tbl_t * tbl,
3370 const int id,
3371 const int ig);
3372
3428double scan_ctl(
3429 int argc,
3430 char *argv[],
3431 const char *varname,
3432 const int arridx,
3433 const char *defvalue,
3434 char *value);
3435
3490void set_cov_apr(
3491 const ret_t * ret,
3492 const ctl_t * ctl,
3493 const atm_t * atm,
3494 const int *iqa,
3495 const int *ipa,
3496 gsl_matrix * s_a);
3497
3552void set_cov_meas(
3553 const ret_t * ret,
3554 const ctl_t * ctl,
3555 const obs_t * obs,
3556 gsl_vector * sig_noise,
3557 gsl_vector * sig_formod,
3558 gsl_vector * sig_eps_inv);
3559
3596double sza(
3597 double sec,
3598 double lon,
3599 double lat);
3600
3641void tangent_point(
3642 const los_t * los,
3643 double *tpz,
3644 double *tplon,
3645 double *tplat);
3646
3671void time2jsec(
3672 const int year,
3673 const int mon,
3674 const int day,
3675 const int hour,
3676 const int min,
3677 const int sec,
3678 const double remain,
3679 double *jsec);
3680
3723void timer(
3724 const char *name,
3725 const char *file,
3726 const char *func,
3727 int line,
3728 int mode);
3729
3772void write_atm(
3773 const char *dirname,
3774 const char *filename,
3775 const ctl_t * ctl,
3776 const atm_t * atm);
3777
3830void write_atm_rfm(
3831 const char *filename,
3832 const ctl_t * ctl,
3833 const atm_t * atm);
3834
3883void write_matrix(
3884 const char *dirname,
3885 const char *filename,
3886 const ctl_t * ctl,
3887 const gsl_matrix * matrix,
3888 const atm_t * atm,
3889 const obs_t * obs,
3890 const char *rowspace,
3891 const char *colspace,
3892 const char *sort);
3893
3951void write_obs(
3952 const char *dirname,
3953 const char *filename,
3954 const ctl_t * ctl,
3955 const obs_t * obs);
3956
3992void write_shape(
3993 const char *filename,
3994 const double *x,
3995 const double *y,
3996 const int n);
3997
4044void write_stddev(
4045 const char *quantity,
4046 const ret_t * ret,
4047 const ctl_t * ctl,
4048 const atm_t * atm,
4049 const gsl_matrix * s);
4050
4070void write_tbl(
4071 const ctl_t * ctl,
4072 const tbl_t * tbl);
4073
4098void write_tbl_asc(
4099 const ctl_t * ctl,
4100 const tbl_t * tbl);
4101
4127void write_tbl_bin(
4128 const ctl_t * ctl,
4129 const tbl_t * tbl);
4130
4158void write_tbl_gas(
4159 const ctl_t * ctl,
4160 const tbl_t * tbl);
4161
4182 const char *path);
4183
4222 tbl_gas_t * g,
4223 const double freq,
4224 const tbl_t * tbl,
4225 const int id,
4226 const int ig);
4227
4275void x2atm(
4276 const ctl_t * ctl,
4277 const gsl_vector * x,
4278 atm_t * atm);
4279
4309void x2atm_help(
4310 double *value,
4311 const gsl_vector * x,
4312 size_t *n);
4313
4347void y2obs(
4348 const ctl_t * ctl,
4349 const gsl_vector * y,
4350 obs_t * obs);
4351
4352#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:5508
#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
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:5144
void read_matrix(const char *dirname, const char *filename, gsl_matrix *matrix)
Read a numerical matrix from an ASCII file.
Definition: jurassic.c:5003
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:6037
void read_rfm_spec(const char *filename, double *nu, double *rad, int *npts)
Read a Reference Forward Model (RFM) ASCII spectrum.
Definition: jurassic.c:5255
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 profile data to a text file.
Definition: jurassic.c:6075
int locate_reg(const double *xx, const int n, const double x)
Locate index for interpolation on a regular (uniform) grid.
Definition: jurassic.c:4408
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:4090
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:3645
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:6775
void idx2name(const ctl_t *ctl, const int idx, char *quantity)
Convert a quantity index to a descriptive name string.
Definition: jurassic.c:3883
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:5660
void read_ctl(int argc, char *argv[], ctl_t *ctl)
Read model control parameters from command-line and configuration input.
Definition: jurassic.c:4889
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:3409
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:4551
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:4479
int locate_irr(const double *xx, const int n, const double x)
Locate index for interpolation on an irregular grid.
Definition: jurassic.c:4378
void write_tbl_asc(const ctl_t *ctl, const tbl_t *tbl)
Write all lookup tables in human-readable ASCII format.
Definition: jurassic.c:6604
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:6006
void x2atm(const ctl_t *ctl, const gsl_vector *x, atm_t *atm)
Map retrieval state vector back to atmospheric structure.
Definition: jurassic.c:6872
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:4001
void write_tbl_bin(const ctl_t *ctl, const tbl_t *tbl)
Write all lookup tables in compact binary format.
Definition: jurassic.c:6649
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:6192
#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:5918
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:4218
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:4186
#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:6581
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:6704
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:3976
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:5568
int find_emitter(const ctl_t *ctl, const char *emitter)
Find gas species index by name.
Definition: jurassic.c:3346
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:5962
#define TBLNU
Maximum number of column densities in emissivity tables.
Definition: jurassic.h:312
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:3510
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_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)
Write observation geometry and radiance data to a text file.
Definition: jurassic.c:6408
#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:6922
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:3445
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:4283
void init_srcfunc(const ctl_t *ctl, tbl_t *tbl)
Initialize the source-function (Planck radiance) lookup table.
Definition: jurassic.c:3922
void matrix_invert(gsl_matrix *a)
Invert a square matrix, optimized for diagonal or symmetric positive-definite matrices.
Definition: jurassic.c:4449
void read_ret(int argc, char *argv[], const ctl_t *ctl, ret_t *ret)
Read retrieval configuration and error parameters.
Definition: jurassic.c:5202
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:5738
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:6518
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
#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:5410
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:3786
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:4250
#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:5358
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:4524
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:6934
void hydrostatic(const ctl_t *ctl, atm_t *atm)
Adjust pressure profile using the hydrostatic equation.
Definition: jurassic.c:3823
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:5316
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:5598
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:5631
#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:3803
#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:6230
#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:5807
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:6548
#define NLOS
Maximum number of LOS points.
Definition: jurassic.h:287
#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:4427
#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:6744
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:1186
int nw
Number of spectral windows.
Definition: jurassic.h:1078
double retp_zmin
Minimum altitude for pressure retrieval [km].
Definition: jurassic.h:1144
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:1108
int ctm_co2
Compute CO2 continuum (0=no, 1=yes).
Definition: jurassic.h:1111
double rett_zmax
Maximum altitude for temperature retrieval [km].
Definition: jurassic.h:1153
int ret_sfeps
Retrieve surface layer emissivity (0=no, 1=yes).
Definition: jurassic.h:1180
int ret_sft
Retrieve surface layer temperature (0=no, 1=yes).
Definition: jurassic.h:1177
int ret_clz
Retrieve cloud layer height (0=no, 1=yes).
Definition: jurassic.h:1168
int ctm_n2
Compute N2 continuum (0=no, 1=yes).
Definition: jurassic.h:1117
int ctm_h2o
Compute H2O continuum (0=no, 1=yes).
Definition: jurassic.h:1114
int formod
Forward model (0=CGA, 1=EGA, 2=RFM).
Definition: jurassic.h:1189
int ng
Number of emitters.
Definition: jurassic.h:1054
int refrac
Take into account refractivity (0=no, 1=yes).
Definition: jurassic.h:1123
int ig_o2
Emitter index of O2.
Definition: jurassic.h:1069
double rett_zmin
Minimum altitude for temperature retrieval [km].
Definition: jurassic.h:1150
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:1141
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 ctm_o2
Compute O2 continuum (0=no, 1=yes).
Definition: jurassic.h:1120
int ret_clk
Retrieve cloud layer extinction (0=no, 1=yes).
Definition: jurassic.h:1174
int nsf
Number of surface layer spectral grid points.
Definition: jurassic.h:1090
double rayds
Maximum step length for raytracing [km].
Definition: jurassic.h:1126
int ret_cldz
Retrieve cloud layer depth (0=no, 1=yes).
Definition: jurassic.h:1171
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:1129
int write_bbt
Use brightness temperature instead of radiance (0=no, 1=yes).
Definition: jurassic.h:1183
double retp_zmax
Maximum altitude for pressure retrieval [km].
Definition: jurassic.h:1147
Line-of-sight data.
Definition: jurassic.h:1209
double sft
Surface temperature [K].
Definition: jurassic.h:1236
int np
Number of LOS points.
Definition: jurassic.h:1212
Observation geometry and radiance data.
Definition: jurassic.h:1271
int nr
Number of ray paths.
Definition: jurassic.h:1274
Retrieval control parameters.
Definition: jurassic.h:1323
double err_press_cz
Vertical correlation length for pressure error [km].
Definition: jurassic.h:1350
double err_press
Pressure error [%].
Definition: jurassic.h:1347
int err_ana
Carry out error analysis (0=no, 1=yes).
Definition: jurassic.h:1338
double err_temp_cz
Vertical correlation length for temperature error [km].
Definition: jurassic.h:1359
double conv_dmin
Minimum normalized step size in state space.
Definition: jurassic.h:1335
double err_temp
Temperature error [K].
Definition: jurassic.h:1356
double err_clz
Cloud height error [km].
Definition: jurassic.h:1383
double err_sft
Surface temperature error [K].
Definition: jurassic.h:1392
double err_temp_ch
Horizontal correlation length for temperature error [km].
Definition: jurassic.h:1362
int kernel_recomp
Re-computation of kernel matrix (number of iterations).
Definition: jurassic.h:1329
int conv_itmax
Maximum number of iterations.
Definition: jurassic.h:1332
double err_cldz
Cloud depth error [km].
Definition: jurassic.h:1386
double err_press_ch
Horizontal correlation length for pressure error [km].
Definition: jurassic.h:1353
On-disk index entry describing one frequency table block in a gas file.
Definition: jurassic.h:1442
int64_t offset
Byte offset in file where the serialized block begins.
Definition: jurassic.h:1448
int64_t size
Size of the serialized block (in bytes).
Definition: jurassic.h:1451
double freq
Frequency identifier ν_j for this table block.
Definition: jurassic.h:1445
In-memory representation of an open per-gas lookup-table file.
Definition: jurassic.h:1462
FILE * fp
Open file handle ("rb+"), NULL if not open.
Definition: jurassic.h:1465
int dirty
Definition: jurassic.h:1474
int32_t ntables
Number of index entries currently in use.
Definition: jurassic.h:1468
tbl_gas_index_t * index
In-memory index table of length MAX_TABLES.
Definition: jurassic.h:1471
Emissivity look-up tables.
Definition: jurassic.h:1405