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-2026 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 <netcdf.h>
116#include <omp.h>
117#include <stdint.h>
118#include <stdio.h>
119#include <stdlib.h>
120#include <string.h>
121#include <time.h>
122
123/* ------------------------------------------------------------
124 Constants...
125 ------------------------------------------------------------ */
126
128#ifndef C1
129#define C1 1.19104259e-8
130#endif
131
133#ifndef C2
134#define C2 1.43877506
135#endif
136
138#ifndef EPSMIN
139#define EPSMIN 0
140#endif
141
143#ifndef EPSMAX
144#define EPSMAX 1
145#endif
146
148#ifndef G0
149#define G0 9.80665
150#endif
151
153#ifndef H0
154#define H0 7.0
155#endif
156
158#ifndef KB
159#define KB 1.3806504e-23
160#endif
161
163#ifndef ME
164#define ME 5.976e24
165#endif
166
168#ifndef NA
169#define NA 6.02214199e23
170#endif
171
173#ifndef N2
174#define N2 0.78084
175#endif
176
178#ifndef O2
179#define O2 0.20946
180#endif
181
183#ifndef P0
184#define P0 1013.25
185#endif
186
188#ifndef RE
189#define RE 6367.421
190#endif
191
193#ifndef RI
194#define RI 8.3144598
195#endif
196
198#ifndef T0
199#define T0 273.15
200#endif
201
203#ifndef TMIN
204#define TMIN 100.
205#endif
206
208#ifndef TMAX
209#define TMAX 400.
210#endif
211
213#ifndef TSUN
214#define TSUN 5780.
215#endif
216
218#ifndef UMIN
219#define UMIN 0
220#endif
221
223#ifndef UMAX
224#define UMAX 1e30
225#endif
226
227/* ------------------------------------------------------------
228 Dimensions...
229 ------------------------------------------------------------ */
230
232#ifndef NCL
233#define NCL 8
234#endif
235
237#ifndef ND
238#define ND 128
239#endif
240
242#ifndef NG
243#define NG 8
244#endif
245
247#ifndef NP
248#define NP 256
249#endif
250
252#ifndef NR
253#define NR 256
254#endif
255
257#ifndef NSF
258#define NSF 8
259#endif
260
262#ifndef NW
263#define NW 4
264#endif
265
267#ifndef LEN
268#define LEN 10000
269#endif
270
272#ifndef M
273#define M (NR*ND)
274#endif
275
277#ifndef N
278#define N ((2 + NG + NW) * NP + NCL + NSF + 3)
279#endif
280
282#ifndef NQ
283#define NQ (5 + NG + NW + NCL + NSF)
284#endif
285
287#ifndef NLOS
288#define NLOS 4096
289#endif
290
292#ifndef NSHAPE
293#define NSHAPE 20000
294#endif
295
297#ifndef NFOV
298#define NFOV 5
299#endif
300
302#ifndef TBLNP
303#define TBLNP 41
304#endif
305
307#ifndef TBLNT
308#define TBLNT 30
309#endif
310
312#ifndef TBLNU
313#define TBLNU 512
314#endif
315
317#ifndef TBLNS
318#define TBLNS 1200
319#endif
320
322#ifndef RFMNPTS
323#define RFMNPTS 10000000
324#endif
325
326/* ------------------------------------------------------------
327 Quantity indices...
328 ------------------------------------------------------------ */
329
331#define IDXP 0
332
334#define IDXT 1
335
337#define IDXQ(ig) (2 + (ig))
338
340#define IDXK(iw) (2 + (ctl->ng) + (iw))
341
343#define IDXCLZ (2 + (ctl->ng) + (ctl->nw))
344
346#define IDXCLDZ (3 + (ctl->ng) + (ctl->nw))
347
349#define IDXCLK(icl) (4 + (ctl->ng) + (ctl->nw) + (icl))
350
352#define IDXSFT (4 + (ctl->ng) + (ctl->nw) + (ctl->ncl))
353
355#define IDXSFEPS(isf) (5 + (ctl->ng) + (ctl->nw) + (ctl->ncl) + (isf))
356
357/* ------------------------------------------------------------
358 Macros...
359 ------------------------------------------------------------ */
360
378#define ALLOC(ptr, type, n) \
379 if((ptr=calloc((size_t)(n), sizeof(type)))==NULL) \
380 ERRMSG("Out of memory!");
381
408#define BRIGHT(rad, nu) \
409 (C2 * (nu) / gsl_log1p(C1 * POW3(nu) / (rad)))
410
430#define CLAMP(v, lo, hi) \
431 (((v) < (lo)) ? (lo) : (((v) > (hi)) ? (hi) : (v)))
432
447#define DEG2RAD(deg) ((deg) * (M_PI / 180.0))
448
465#define DIST(a, b) sqrt(DIST2(a, b))
466
482#define DIST2(a, b) \
483 ((a[0]-b[0])*(a[0]-b[0])+(a[1]-b[1])*(a[1]-b[1])+(a[2]-b[2])*(a[2]-b[2]))
484
500#define DOTP(a, b) (a[0]*b[0]+a[1]*b[1]+a[2]*b[2])
501
520#define FREAD(ptr, type, size, out) { \
521 if(fread(ptr, sizeof(type), size, out)!=size) \
522 ERRMSG("Error while reading!"); \
523 }
524
543#define FWRITE(ptr, type, size, out) { \
544 if(fwrite(ptr, sizeof(type), size, out)!=size) \
545 ERRMSG("Error while writing!"); \
546 }
547
566#define LIN(x0, y0, x1, y1, x) \
567 ((y0)+((y1)-(y0))/((x1)-(x0))*((x)-(x0)))
568
587#define LOGY(x0, y0, x1, y1, x) \
588 (((y1)/(y0)>0) \
589 ? ((y0)*exp(log((y1)/(y0))/((x1)-(x0))*((x)-(x0)))) \
590 : LIN(x0, y0, x1, y1, x))
591
608#define MAX(a,b) (((a)>(b))?(a):(b))
609
626#define MIN(a,b) (((a)<(b))?(a):(b))
627
639#define NC(cmd) { \
640 int nc_result=(cmd); \
641 if(nc_result!=NC_NOERR) \
642 ERRMSG("%s", nc_strerror(nc_result)); \
643 }
644
668#define NC_DEF_VAR(varname, type, ndims, dims, long_name, units, level, quant) { \
669 NC(nc_def_var(ncid, varname, type, ndims, dims, &varid)); \
670 NC(nc_put_att_text(ncid, varid, "long_name", strnlen(long_name, LEN), long_name)); \
671 NC(nc_put_att_text(ncid, varid, "units", strnlen(units, LEN), units)); \
672 if((quant) > 0) \
673 NC(nc_def_var_quantize(ncid, varid, NC_QUANTIZE_GRANULARBR, quant)); \
674 if((level) != 0) { \
675 NC(nc_def_var_deflate(ncid, varid, 1, 1, level)); \
676 /* unsigned int ulevel = (unsigned int)level; */ \
677 /* NC(nc_def_var_filter(ncid, varid, 32015, 1, (unsigned int[]){ulevel})); */ \
678 } \
679 }
680
698#define NC_GET_DOUBLE(varname, ptr, force) { \
699 if(force) { \
700 NC(nc_inq_varid(ncid, varname, &varid)); \
701 NC(nc_get_var_double(ncid, varid, ptr)); \
702 } else { \
703 if(nc_inq_varid(ncid, varname, &varid) == NC_NOERR) { \
704 NC(nc_get_var_double(ncid, varid, ptr)); \
705 } else \
706 WARN("netCDF variable %s is missing!", varname); \
707 } \
708 }
709
728#define NC_INQ_DIM(dimname, ptr, min, max, check) { \
729 int dimid; size_t naux; \
730 NC(nc_inq_dimid(ncid, dimname, &dimid)); \
731 NC(nc_inq_dimlen(ncid, dimid, &naux)); \
732 *ptr = (int)naux; \
733 if (check) \
734 if ((*ptr) < (min) || (*ptr) > (max)) \
735 ERRMSG("Dimension %s is out of range!", dimname); \
736 }
737
752#define NC_PUT_DOUBLE(varname, ptr, hyperslab) { \
753 NC(nc_inq_varid(ncid, varname, &varid)); \
754 if(hyperslab) { \
755 NC(nc_put_vara_double(ncid, varid, start, count, ptr)); \
756 } else { \
757 NC(nc_put_var_double(ncid, varid, ptr)); \
758 } \
759 }
760
776#define NC_PUT_FLOAT(varname, ptr, hyperslab) { \
777 NC(nc_inq_varid(ncid, varname, &varid)); \
778 if(hyperslab) { \
779 NC(nc_put_vara_float(ncid, varid, start, count, ptr)); \
780 } else { \
781 NC(nc_put_var_float(ncid, varid, ptr)); \
782 } \
783 }
784
799#define NC_PUT_INT(varname, ptr, hyperslab) { \
800 NC(nc_inq_varid(ncid, varname, &varid)); \
801 if(hyperslab) { \
802 NC(nc_put_vara_int(ncid, varid, start, count, ptr)); \
803 } else { \
804 NC(nc_put_var_int(ncid, varid, ptr)); \
805 } \
806 }
807
821#define NC_PUT_ATT(varname, attname, text) { \
822 NC(nc_inq_varid(ncid, varname, &varid)); \
823 NC(nc_put_att_text(ncid, varid, attname, strnlen(text, LEN), text)); \
824 }
825
838#define NC_PUT_ATT_GLOBAL(attname, text) \
839 NC(nc_put_att_text(ncid, NC_GLOBAL, attname, strnlen(text, LEN), text));
840
870#define NEDT(t_bg, nesr, nu) \
871 (BRIGHT(PLANCK((t_bg), (nu)) + (nesr), (nu)) - (t_bg))
872
899#define NESR(t_bg, nedt, nu) \
900 (PLANCK((t_bg) + (nedt), (nu)) - PLANCK((t_bg), (nu)))
901
916#define NORM(a) sqrt(DOTP(a, a))
917
947#define PLANCK(T, nu) \
948 (C1 * POW3(nu) / gsl_expm1(C2 * (nu) / (T)))
949
961#define POW2(x) ((x)*(x))
962
974#define POW3(x) ((x)*(x)*(x))
975
990#define RAD2DEG(rad) ((rad) * (180.0 / M_PI))
991
1006#define REFRAC(p, T) (7.753e-05 * (p) / (T))
1007
1022#define TBL_LOG(tbl, id, ig) \
1023 do { \
1024 for (int ip = 0; ip < (tbl)->np[(id)][(ig)]; ip++) \
1025 LOG(2, \
1026 "p[%2d]= %.5e hPa | T[0:%2d]= %.2f ... %.2f K | " \
1027 "u[0:%3d]= %.5e ... %.5e molec/cm^2 | " \
1028 "eps[0:%3d]= %.5e ... %.5e", \
1029 (ip), \
1030 (tbl)->p[(id)][(ig)][(ip)], \
1031 (tbl)->nt[(id)][(ig)][(ip)] - 1, \
1032 (tbl)->t[(id)][(ig)][(ip)][0], \
1033 (tbl)->t[(id)][(ig)][(ip)] \
1034 [(tbl)->nt[(id)][(ig)][(ip)] - 1], \
1035 (tbl)->nu[(id)][(ig)][(ip)][0] - 1, \
1036 exp((tbl)->logu[(id)][(ig)][(ip)][0][0]), \
1037 exp((tbl)->logu[(id)][(ig)][(ip)][0] \
1038 [(tbl)->nu[(id)][(ig)][(ip)][0] - 1]), \
1039 (tbl)->nu[(id)][(ig)][(ip)][0] - 1, \
1040 exp((tbl)->logeps[(id)][(ig)][(ip)][0][0]), \
1041 exp((tbl)->logeps[(id)][(ig)][(ip)][0] \
1042 [(tbl)->nu[(id)][(ig)][(ip)][0] - 1])); \
1043 } while (0)
1044
1058#define TIMER(name, mode) \
1059 {timer(name, __FILE__, __func__, __LINE__, mode);}
1060
1079#define TOK(line, tok, format, var) { \
1080 if(((tok)=strtok((line), " \t"))) { \
1081 if(sscanf(tok, format, &(var))!=1) continue; \
1082 } else ERRMSG("Error while reading!"); \
1083 }
1084
1085/* ------------------------------------------------------------
1086 Log messages...
1087 ------------------------------------------------------------ */
1088
1090#ifndef LOGLEV
1091#define LOGLEV 2
1092#endif
1093
1123#define LOG(level, ...) { \
1124 if(level >= 2) \
1125 printf(" "); \
1126 if(level <= LOGLEV) { \
1127 printf(__VA_ARGS__); \
1128 printf("\n"); \
1129 } \
1130 }
1131
1160#define WARN(...) { \
1161 printf("\nWarning (%s, %s, l%d): ", __FILE__, __func__, __LINE__); \
1162 LOG(0, __VA_ARGS__); \
1163 }
1164
1193#define ERRMSG(...) { \
1194 printf("\nError (%s, %s, l%d): ", __FILE__, __func__, __LINE__); \
1195 LOG(0, __VA_ARGS__); \
1196 exit(EXIT_FAILURE); \
1197 }
1198
1228#define PRINT(format, var) \
1229 printf("Print (%s, %s, l%d): %s= "format"\n", \
1230 __FILE__, __func__, __LINE__, #var, var);
1231
1232/* ------------------------------------------------------------
1233 Structs...
1234 ------------------------------------------------------------ */
1235
1243typedef struct {
1244
1246 int np;
1247
1249 double time[NP];
1250
1252 double z[NP];
1253
1255 double lon[NP];
1256
1258 double lat[NP];
1259
1261 double p[NP];
1262
1264 double t[NP];
1265
1267 double q[NG][NP];
1268
1270 double k[NW][NP];
1271
1273 double clz;
1274
1276 double cldz;
1277
1279 double clk[NCL];
1280
1282 double sft;
1283
1285 double sfeps[NSF];
1286
1287} atm_t;
1288
1296typedef struct {
1297
1299 int ng;
1300
1302 char emitter[NG][LEN];
1303
1306
1309
1312
1315
1317 int nd;
1318
1320 double nu[ND];
1321
1323 int nw;
1324
1326 int window[ND];
1327
1329 int ncl;
1330
1332 double clnu[NCL];
1333
1335 int nsf;
1336
1338 double sfnu[NSF];
1339
1342
1344 double sfsza;
1345
1347 char tblbase[LEN];
1348
1351
1354
1357
1359 double hydz;
1360
1363
1366
1369
1372
1375
1377 double rayds;
1378
1380 double raydz;
1381
1383 char fov[LEN];
1384
1386 double fov_dz[NSHAPE];
1387
1389 double fov_w[NSHAPE];
1390
1393
1396
1399
1402
1405
1407 double retq_zmin[NG];
1408
1410 double retq_zmax[NG];
1411
1413 double retk_zmin[NW];
1414
1416 double retk_zmax[NW];
1417
1420
1423
1426
1429
1432
1435
1438
1441
1443 char rfmbin[LEN];
1444
1446 char rfmhit[LEN];
1447
1449 char rfmxsc[NG][LEN];
1450
1451} ctl_t;
1452
1460typedef struct {
1461
1463 int np;
1464
1466 double z[NLOS];
1467
1469 double lon[NLOS];
1470
1472 double lat[NLOS];
1473
1475 double p[NLOS];
1476
1478 double t[NLOS];
1479
1481 double q[NLOS][NG];
1482
1484 double k[NLOS][ND];
1485
1487 double sft;
1488
1490 double sfeps[ND];
1491
1493 double ds[NLOS];
1494
1496 double u[NLOS][NG];
1497
1499 double cgp[NLOS][NG];
1500
1502 double cgt[NLOS][NG];
1503
1505 double cgu[NLOS][NG];
1506
1508 double eps[NLOS][ND];
1509
1511 double src[NLOS][ND];
1512
1513} los_t;
1514
1522typedef struct {
1523
1525 int nr;
1526
1528 double time[NR];
1529
1531 double obsz[NR];
1532
1534 double obslon[NR];
1535
1537 double obslat[NR];
1538
1540 double vpz[NR];
1541
1543 double vplon[NR];
1544
1546 double vplat[NR];
1547
1549 double tpz[NR];
1550
1552 double tplon[NR];
1553
1555 double tplat[NR];
1556
1558 double tau[ND][NR];
1559
1561 double rad[ND][NR];
1562
1563} obs_t;
1564
1574typedef struct {
1575
1577 char dir[LEN];
1578
1581
1584
1587
1590
1592 double err_formod[ND];
1593
1595 double err_noise[ND];
1596
1599
1602
1605
1607 double err_temp;
1608
1611
1614
1616 double err_q[NG];
1617
1619 double err_q_cz[NG];
1620
1622 double err_q_ch[NG];
1623
1625 double err_k[NW];
1626
1628 double err_k_cz[NW];
1629
1631 double err_k_ch[NW];
1632
1634 double err_clz;
1635
1637 double err_cldz;
1638
1640 double err_clk[NCL];
1641
1643 double err_sft;
1644
1646 double err_sfeps[NSF];
1647
1648} ret_t;
1649
1656typedef struct {
1657
1659 int np[ND][NG];
1660
1662 int nt[ND][NG][TBLNP];
1663
1665 int nu[ND][NG][TBLNP][TBLNT];
1666
1668 double p[ND][NG][TBLNP];
1669
1671 double lnp[ND][NG][TBLNP];
1672
1674 double t[ND][NG][TBLNP][TBLNT];
1675
1677 float *logu[ND][NG][TBLNP][TBLNT];
1678
1680 float *logeps[ND][NG][TBLNP][TBLNT];
1681
1683 int filt_n[ND];
1684
1686 double filt_nu[ND][NSHAPE];
1687
1689 double filt_f[ND][NSHAPE];
1690
1692 double st[TBLNS];
1693
1695 double sr[TBLNS][ND];
1696
1697} tbl_t;
1698
1699/* ------------------------------------------------------------
1700 Functions...
1701 ------------------------------------------------------------ */
1702
1754void analyze_avk(
1755 const ret_t * ret,
1756 const ctl_t * ctl,
1757 const atm_t * atm,
1758 const int *iqa,
1759 const int *ipa,
1760 const gsl_matrix * avk);
1761
1818 const gsl_matrix * avk,
1819 const int iq,
1820 const int *ipa,
1821 const size_t *n0,
1822 const size_t *n1,
1823 double *cont,
1824 double *res);
1825
1852size_t atm2x(
1853 const ctl_t * ctl,
1854 const atm_t * atm,
1855 gsl_vector * x,
1856 int *iqa,
1857 int *ipa);
1858
1880void atm2x_help(
1881 const double value,
1882 const int value_iqa,
1883 const int value_ip,
1884 gsl_vector * x,
1885 int *iqa,
1886 int *ipa,
1887 size_t *n);
1888
1903void cart2geo(
1904 const double *x,
1905 double *z,
1906 double *lon,
1907 double *lat);
1908
1924void climatology(
1925 const ctl_t * ctl,
1926 atm_t * atm);
1927
1952double cos_sza(
1953 const double sec,
1954 const double lon,
1955 const double lat);
1956
1996double cost_function(
1997 const gsl_vector * dx,
1998 const gsl_vector * dy,
1999 const gsl_matrix * s_a_inv,
2000 const gsl_vector * sig_eps_inv);
2001
2007double ctmco2(
2008 const double nu,
2009 const double p,
2010 const double t,
2011 const double u);
2012
2018double ctmh2o(
2019 const double nu,
2020 const double p,
2021 const double t,
2022 const double q,
2023 const double u);
2024
2054double ctmn2(
2055 const double nu,
2056 const double p,
2057 const double t);
2058
2088double ctmo2(
2089 const double nu,
2090 const double p,
2091 const double t);
2092
2117void copy_atm(
2118 const ctl_t * ctl,
2119 atm_t * atm_dest,
2120 const atm_t * atm_src,
2121 const int init);
2122
2147void copy_obs(
2148 const ctl_t * ctl,
2149 obs_t * obs_dest,
2150 const obs_t * obs_src,
2151 const int init);
2152
2173void day2doy(
2174 int year,
2175 int mon,
2176 int day,
2177 int *doy);
2178
2201void doy2day(
2202 int year,
2203 int doy,
2204 int *mon,
2205 int *day);
2206
2226int find_emitter(
2227 const ctl_t * ctl,
2228 const char *emitter);
2229
2255void formod(
2256 const ctl_t * ctl,
2257 const tbl_t * tbl,
2258 atm_t * atm,
2259 obs_t * obs);
2260
2286void formod_continua(
2287 const ctl_t * ctl,
2288 const los_t * los,
2289 const int ip,
2290 double *beta);
2291
2313void formod_fov(
2314 const ctl_t * ctl,
2315 obs_t * obs);
2316
2346void formod_pencil(
2347 const ctl_t * ctl,
2348 const tbl_t * tbl,
2349 const atm_t * atm,
2350 obs_t * obs,
2351 const int ir);
2352
2404void formod_rfm(
2405 const ctl_t * ctl,
2406 const tbl_t * tbl,
2407 const atm_t * atm,
2408 obs_t * obs);
2409
2434void formod_srcfunc(
2435 const ctl_t * ctl,
2436 const tbl_t * tbl,
2437 const double t,
2438 double *src);
2439
2466void geo2cart(
2467 const double z,
2468 const double lon,
2469 const double lat,
2470 double *x);
2471
2499void hydrostatic(
2500 const ctl_t * ctl,
2501 atm_t * atm);
2502
2526void idx2name(
2527 const ctl_t * ctl,
2528 const int idx,
2529 char *quantity);
2530
2560void init_srcfunc(
2561 const ctl_t * ctl,
2562 tbl_t * tbl);
2563
2588void intpol_atm(
2589 const ctl_t * ctl,
2590 const atm_t * atm,
2591 const double z,
2592 double *p,
2593 double *t,
2594 double *q,
2595 double *k);
2596
2630void intpol_tbl_cga(
2631 const ctl_t * ctl,
2632 const tbl_t * tbl,
2633 const los_t * los,
2634 const int ip,
2635 double tau_path[ND][NG],
2636 double tau_seg[ND]);
2637
2676void intpol_tbl_ega(
2677 const ctl_t * ctl,
2678 const tbl_t * tbl,
2679 const los_t * los,
2680 const int ip,
2681 double tau_path[ND][NG],
2682 double tau_seg[ND]);
2683
2710double intpol_tbl_eps(
2711 const tbl_t * tbl,
2712 const int ig,
2713 const int id,
2714 const int ip,
2715 const int it,
2716 const double logu);
2717
2748double intpol_tbl_u(
2749 const tbl_t * tbl,
2750 const int ig,
2751 const int id,
2752 const int ip,
2753 const int it,
2754 const double logeps);
2755
2782void jsec2time(
2783 const double jsec,
2784 int *year,
2785 int *mon,
2786 int *day,
2787 int *hour,
2788 int *min,
2789 int *sec,
2790 double *remain);
2791
2833void kernel(
2834 const ctl_t * ctl,
2835 const tbl_t * tbl,
2836 atm_t * atm,
2837 obs_t * obs,
2838 gsl_matrix * k);
2839
2866int locate_irr(
2867 const double *xx,
2868 const int n,
2869 const double x);
2870
2897int locate_reg(
2898 const double *xx,
2899 const int n,
2900 const double x);
2901
2927int locate_tbl(
2928 const float *xx,
2929 const int n,
2930 const double x);
2931
2970void matrix_invert(
2971 gsl_matrix * a);
2972
3017void matrix_product(
3018 const gsl_matrix * a,
3019 const gsl_vector * b,
3020 const int transpose,
3021 gsl_matrix * c);
3022
3050size_t obs2y(
3051 const ctl_t * ctl,
3052 const obs_t * obs,
3053 gsl_vector * y,
3054 int *ida,
3055 int *ira);
3056
3108 ret_t * ret,
3109 ctl_t * ctl,
3110 tbl_t * tbl,
3111 obs_t * obs_meas,
3112 obs_t * obs_i,
3113 atm_t * atm_apr,
3114 atm_t * atm_i,
3115 double *chisq);
3116
3155void raytrace(
3156 const ctl_t * ctl,
3157 const atm_t * atm,
3158 obs_t * obs,
3159 los_t * los,
3160 const int ir);
3161
3207void read_atm(
3208 const char *dirname,
3209 const char *filename,
3210 const ctl_t * ctl,
3211 atm_t * atm);
3212
3259void read_atm_asc(
3260 const char *filename,
3261 const ctl_t * ctl,
3262 atm_t * atm);
3263
3313void read_atm_bin(
3314 const char *filename,
3315 const ctl_t * ctl,
3316 atm_t * atm);
3317
3357void read_atm_nc(
3358 const char *filename,
3359 const ctl_t * ctl,
3360 atm_t * atm,
3361 int profile);
3362
3403void read_ctl(
3404 int argc,
3405 char *argv[],
3406 ctl_t * ctl);
3407
3439void read_matrix(
3440 const char *dirname,
3441 const char *filename,
3442 gsl_matrix * matrix);
3443
3478void read_obs(
3479 const char *dirname,
3480 const char *filename,
3481 const ctl_t * ctl,
3482 obs_t * obs);
3483
3520void read_obs_asc(
3521 const char *filename,
3522 const ctl_t * ctl,
3523 obs_t * obs);
3524
3565void read_obs_bin(
3566 const char *filename,
3567 const ctl_t * ctl,
3568 obs_t * obs);
3569
3601void read_obs_nc(
3602 const char *filename,
3603 const ctl_t * ctl,
3604 obs_t * obs,
3605 const int profile);
3606
3646double read_obs_rfm(
3647 const char *basename,
3648 const double z,
3649 const double *nu,
3650 const double *f,
3651 const int n);
3652
3717void read_ret(
3718 int argc,
3719 char *argv[],
3720 const ctl_t * ctl,
3721 ret_t * ret);
3722
3763void read_rfm_spec(
3764 const char *filename,
3765 double *nu,
3766 double *rad,
3767 int *npts);
3768
3802void read_shape(
3803 const char *filename,
3804 double *x,
3805 double *y,
3806 int *n);
3807
3844 const ctl_t * ctl);
3845
3881void read_tbl_asc(
3882 const ctl_t * ctl,
3883 tbl_t * tbl,
3884 const int id,
3885 const int ig);
3886
3918void read_tbl_bin(
3919 const ctl_t * ctl,
3920 tbl_t * tbl,
3921 const int id,
3922 const int ig);
3923
3967 const ctl_t * ctl,
3968 tbl_t * tbl,
3969 int id,
3970 int ig,
3971 int ncid);
3972
4027double scan_ctl(
4028 int argc,
4029 char *argv[],
4030 const char *varname,
4031 const int arridx,
4032 const char *defvalue,
4033 char *value);
4034
4088void set_cov_apr(
4089 const ret_t * ret,
4090 const ctl_t * ctl,
4091 const atm_t * atm,
4092 const int *iqa,
4093 const int *ipa,
4094 gsl_matrix * s_a);
4095
4149void set_cov_meas(
4150 const ret_t * ret,
4151 const ctl_t * ctl,
4152 const obs_t * obs,
4153 gsl_vector * sig_noise,
4154 gsl_vector * sig_formod,
4155 gsl_vector * sig_eps_inv);
4156
4192double sza(
4193 double sec,
4194 double lon,
4195 double lat);
4196
4236void tangent_point(
4237 const los_t * los,
4238 double *tpz,
4239 double *tplon,
4240 double *tplat);
4241
4263void tbl_free(
4264 const ctl_t * ctl,
4265 tbl_t * tbl);
4266
4306void tbl_pack(
4307 const tbl_t * tbl,
4308 int id,
4309 int ig,
4310 uint8_t * buf,
4311 size_t *bytes_used);
4312
4323size_t tbl_packed_size(
4324 const tbl_t * tbl,
4325 int id,
4326 int ig);
4327
4366size_t tbl_unpack(
4367 tbl_t * tbl,
4368 int id,
4369 int ig,
4370 const uint8_t * buf);
4371
4372
4397void time2jsec(
4398 const int year,
4399 const int mon,
4400 const int day,
4401 const int hour,
4402 const int min,
4403 const int sec,
4404 const double remain,
4405 double *jsec);
4406
4448void timer(
4449 const char *name,
4450 const char *file,
4451 const char *func,
4452 int line,
4453 int mode);
4454
4499void write_atm(
4500 const char *dirname,
4501 const char *filename,
4502 const ctl_t * ctl,
4503 const atm_t * atm);
4504
4556void write_atm_asc(
4557 const char *filename,
4558 const ctl_t * ctl,
4559 const atm_t * atm);
4560
4614void write_atm_bin(
4615 const char *filename,
4616 const ctl_t * ctl,
4617 const atm_t * atm);
4618
4662void write_atm_nc(
4663 const char *filename,
4664 const ctl_t * ctl,
4665 const atm_t * atm,
4666 int profile);
4667
4719void write_atm_rfm(
4720 const char *filename,
4721 const ctl_t * ctl,
4722 const atm_t * atm);
4723
4771void write_matrix(
4772 const char *dirname,
4773 const char *filename,
4774 const ctl_t * ctl,
4775 const gsl_matrix * matrix,
4776 const atm_t * atm,
4777 const obs_t * obs,
4778 const char *rowspace,
4779 const char *colspace,
4780 const char *sort);
4781
4814void write_obs(
4815 const char *dirname,
4816 const char *filename,
4817 const ctl_t * ctl,
4818 const obs_t * obs);
4819
4858void write_obs_asc(
4859 const char *filename,
4860 const ctl_t * ctl,
4861 const obs_t * obs);
4862
4899void write_obs_bin(
4900 const char *filename,
4901 const ctl_t * ctl,
4902 const obs_t * obs);
4903
4939void write_obs_nc(
4940 const char *filename,
4941 const ctl_t * ctl,
4942 const obs_t * obs,
4943 const int profile);
4944
4979void write_shape(
4980 const char *filename,
4981 const double *x,
4982 const double *y,
4983 const int n);
4984
5030void write_stddev(
5031 const char *quantity,
5032 const ret_t * ret,
5033 const ctl_t * ctl,
5034 const atm_t * atm,
5035 const gsl_matrix * s);
5036
5072void write_tbl(
5073 const ctl_t * ctl,
5074 const tbl_t * tbl);
5075
5101void write_tbl_asc(
5102 const ctl_t * ctl,
5103 const tbl_t * tbl,
5104 const int id,
5105 const int ig);
5106
5134void write_tbl_bin(
5135 const ctl_t * ctl,
5136 const tbl_t * tbl,
5137 const int id,
5138 const int ig);
5139
5172void write_tbl_nc(
5173 const ctl_t * ctl,
5174 const tbl_t * tbl,
5175 const int id,
5176 const int ig);
5177
5224void x2atm(
5225 const ctl_t * ctl,
5226 const gsl_vector * x,
5227 atm_t * atm);
5228
5257void x2atm_help(
5258 double *value,
5259 const gsl_vector * x,
5260 size_t *n);
5261
5294void y2obs(
5295 const ctl_t * ctl,
5296 const gsl_vector * y,
5297 obs_t * obs);
5298
5299#endif
double intpol_tbl_eps(const tbl_t *tbl, const int ig, const int id, const int ip, const int it, const double logu)
Interpolate gas emissivity as a function of column amount.
Definition: jurassic.c:4252
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 compact binary emissivity lookup table.
Definition: jurassic.c:6332
#define LEN
Maximum length of ASCII data lines.
Definition: jurassic.h:268
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:5919
void read_matrix(const char *dirname, const char *filename, gsl_matrix *matrix)
Read a numerical matrix from an ASCII file.
Definition: jurassic.c:5601
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:6915
void read_rfm_spec(const char *filename, double *nu, double *rad, int *npts)
Read a Reference Forward Model (RFM) ASCII spectrum.
Definition: jurassic.c:6030
void tbl_pack(const tbl_t *tbl, int id, int ig, uint8_t *buf, size_t *bytes_used)
Pack a lookup table into a contiguous binary buffer.
Definition: jurassic.c:6720
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 read_atm_bin(const char *filename, const ctl_t *ctl, atm_t *atm)
Read atmospheric data in binary format.
Definition: jurassic.c:5279
void write_obs_asc(const char *filename, const ctl_t *ctl, const obs_t *obs)
Write observation data to an ASCII text file.
Definition: jurassic.c:7647
void formod_rfm(const ctl_t *ctl, const tbl_t *tbl, const atm_t *atm, obs_t *obs)
Forward-model radiance and transmittance with the Reference Forward Model (RFM).
Definition: jurassic.c:3693
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:6953
int locate_reg(const double *xx, const int n, const double x)
Locate index for interpolation on a regular (uniform) grid.
Definition: jurassic.c:4501
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:4153
double ctmo2(const double nu, const double p, const double t)
Compute O₂ collision-induced absorption coefficient.
Definition: jurassic.c:3194
void write_obs_nc(const char *filename, const ctl_t *ctl, const obs_t *obs, const int profile)
Write one observation profile to a NetCDF file.
Definition: jurassic.c:7771
void tbl_free(const ctl_t *ctl, tbl_t *tbl)
Free lookup table and all internally allocated memory.
Definition: jurassic.c:6681
void write_tbl_asc(const ctl_t *ctl, const tbl_t *tbl, const int id, const int ig)
Write one emissivity lookup table in human-readable ASCII format.
Definition: jurassic.c:8052
void write_tbl_nc(const ctl_t *ctl, const tbl_t *tbl, const int id, const int ig)
Write one packed lookup table to a NetCDF file.
Definition: jurassic.c:8136
void idx2name(const ctl_t *ctl, const int idx, char *quantity)
Convert a quantity index to a descriptive name string.
Definition: jurassic.c:3952
void read_ctl(int argc, char *argv[], ctl_t *ctl)
Read model control parameters from command-line and configuration input.
Definition: jurassic.c:5477
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:3457
void read_obs_bin(const char *filename, const ctl_t *ctl, obs_t *obs)
Read binary-formatted observation data from a file.
Definition: jurassic.c:5763
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:4916
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:4572
int locate_irr(const double *xx, const int n, const double x)
Locate index for interpolation on an irregular grid.
Definition: jurassic.c:4471
void day2doy(int year, int mon, int day, int *doy)
Convert a calendar date to day-of-year.
Definition: jurassic.c:3344
void write_atm_nc(const char *filename, const ctl_t *ctl, const atm_t *atm, int profile)
Write one atmospheric profile to a netCDF file.
Definition: jurassic.c:7175
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:6884
void x2atm(const ctl_t *ctl, const gsl_vector *x, atm_t *atm)
Map retrieval state vector back to atmospheric structure.
Definition: jurassic.c:8201
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:2073
void write_tbl_bin(const ctl_t *ctl, const tbl_t *tbl, const int id, const int ig)
Write one emissivity lookup table in compact binary format.
Definition: jurassic.c:8094
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:4063
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:7358
#define ND
Maximum number of radiance channels.
Definition: jurassic.h:238
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:6593
#define NSHAPE
Maximum number of shape function grid points.
Definition: jurassic.h:293
void write_tbl(const ctl_t *ctl, const tbl_t *tbl)
Write emissivity lookup tables to disk.
Definition: jurassic.c:8009
size_t tbl_unpack(tbl_t *tbl, int id, int ig, const uint8_t *buf)
Unpack a lookup table from a contiguous binary buffer.
Definition: jurassic.c:6808
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:4038
int find_emitter(const ctl_t *ctl, const char *emitter)
Find gas species index by name.
Definition: jurassic.c:3394
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:4644
double intpol_tbl_u(const tbl_t *tbl, const int ig, const int id, const int ip, const int it, const double logeps)
Interpolate column amount as a function of emissivity.
Definition: jurassic.c:4297
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:6637
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:3558
double ctmco2(const double nu, const double p, const double t, const double u)
Compute carbon dioxide continuum (optical depth).
Definition: jurassic.c:1210
void write_atm_asc(const char *filename, const ctl_t *ctl, const atm_t *atm)
Write atmospheric data to an ASCII file.
Definition: jurassic.c:7024
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:5153
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:3306
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:5641
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:7574
#define TBLNT
Maximum number of temperatures in emissivity tables.
Definition: jurassic.h:308
void write_obs_bin(const char *filename, const ctl_t *ctl, const obs_t *obs)
Write observation data in binary format to a file.
Definition: jurassic.c:7703
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:8251
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:248
void formod_fov(const ctl_t *ctl, obs_t *obs)
Apply field-of-view (FOV) convolution to modeled radiances.
Definition: jurassic.c:3493
void doy2day(int year, int doy, int *mon, int *day)
Convert a day-of-year value to a calendar date.
Definition: jurassic.c:3364
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:4376
void read_obs_nc(const char *filename, const ctl_t *ctl, obs_t *obs, const int profile)
Read one observation profile from a NetCDF file.
Definition: jurassic.c:5841
void init_srcfunc(const ctl_t *ctl, tbl_t *tbl)
Initialize source function lookup tables from emissivity data.
Definition: jurassic.c:3991
void matrix_invert(gsl_matrix *a)
Invert a square matrix, optimized for diagonal or symmetric positive-definite matrices.
Definition: jurassic.c:4542
void read_ret(int argc, char *argv[], const ctl_t *ctl, ret_t *ret)
Read retrieval configuration and error parameters.
Definition: jurassic.c:5977
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:6413
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:7946
void formod(const ctl_t *ctl, const tbl_t *tbl, atm_t *atm, obs_t *obs)
Execute the selected forward model.
Definition: jurassic.c:3407
void read_atm_asc(const char *filename, const ctl_t *ctl, atm_t *atm)
Read atmospheric data in ASCII format.
Definition: jurassic.c:5227
#define NG
Maximum number of emitters.
Definition: jurassic.h:243
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:6209
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:3855
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:4343
#define TBLNS
Maximum number of source function temperature levels.
Definition: jurassic.h:318
tbl_t * read_tbl(const ctl_t *ctl)
Read emissivity lookup tables from disk.
Definition: jurassic.c:6133
double sza(double sec, double lon, double lat)
Compute the solar zenith angle for a given time and location.
void read_tbl_nc_channel(const ctl_t *ctl, tbl_t *tbl, int id, int ig, int ncid)
Read one packed emissivity lookup table from an open NetCDF file.
Definition: jurassic.c:6376
void read_atm_nc(const char *filename, const ctl_t *ctl, atm_t *atm, int profile)
Read one atmospheric profile from a netCDF file.
Definition: jurassic.c:5373
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:3256
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:4617
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:8263
void hydrostatic(const ctl_t *ctl, atm_t *atm)
Adjust pressure profile using the hydrostatic equation.
Definition: jurassic.c:3892
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:6091
double ctmn2(const double nu, const double p, const double t)
Compute N₂ collision-induced absorption coefficient.
Definition: jurassic.c:3125
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
#define TBLNP
Maximum number of pressure levels in emissivity tables.
Definition: jurassic.h:303
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:3872
#define NSF
Maximum number of surface layer spectral grid points.
Definition: jurassic.h:258
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:7396
#define NCL
Maximum number of cloud layer spectral grid points.
Definition: jurassic.h:233
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:6482
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:7976
#define NLOS
Maximum number of LOS points.
Definition: jurassic.h:288
size_t tbl_packed_size(const tbl_t *tbl, int id, int ig)
Compute required buffer size (in bytes) for tbl_pack().
Definition: jurassic.c:6774
#define NR
Maximum number of ray paths.
Definition: jurassic.h:253
int locate_tbl(const float *xx, const int n, const double x)
Locate index for interpolation within emissivity table grids.
Definition: jurassic.c:4520
void write_atm_bin(const char *filename, const ctl_t *ctl, const atm_t *atm)
Write atmospheric data to a binary file.
Definition: jurassic.c:7091
void read_obs_asc(const char *filename, const ctl_t *ctl, obs_t *obs)
Read ASCII-formatted observation data from a file.
Definition: jurassic.c:5718
#define NW
Maximum number of spectral windows.
Definition: jurassic.h:263
Atmospheric profile data.
Definition: jurassic.h:1243
double clz
Cloud layer height [km].
Definition: jurassic.h:1273
int np
Number of data points.
Definition: jurassic.h:1246
double cldz
Cloud layer depth [km].
Definition: jurassic.h:1276
double sft
Surface temperature [K].
Definition: jurassic.h:1282
Control parameters.
Definition: jurassic.h:1296
int write_matrix
Write matrix file (0=no, 1=yes).
Definition: jurassic.h:1437
int nw
Number of spectral windows.
Definition: jurassic.h:1323
int atmfmt
Atmospheric data file format (1=ASCII, 2=binary, 3=netCDF).
Definition: jurassic.h:1353
double retp_zmin
Minimum altitude for pressure retrieval [km].
Definition: jurassic.h:1395
int ig_co2
Emitter index of CO2.
Definition: jurassic.h:1305
double hydz
Reference height for hydrostatic pressure profile (-999 to skip) [km].
Definition: jurassic.h:1359
int ctm_co2
Compute CO2 continuum (0=no, 1=yes).
Definition: jurassic.h:1362
double rett_zmax
Maximum altitude for temperature retrieval [km].
Definition: jurassic.h:1404
int ret_sfeps
Retrieve surface layer emissivity (0=no, 1=yes).
Definition: jurassic.h:1431
int ret_sft
Retrieve surface layer temperature (0=no, 1=yes).
Definition: jurassic.h:1428
int ret_clz
Retrieve cloud layer height (0=no, 1=yes).
Definition: jurassic.h:1419
int ctm_n2
Compute N2 continuum (0=no, 1=yes).
Definition: jurassic.h:1368
int ctm_h2o
Compute H2O continuum (0=no, 1=yes).
Definition: jurassic.h:1365
int formod
Forward model (0=CGA, 1=EGA, 2=RFM).
Definition: jurassic.h:1440
int ng
Number of emitters.
Definition: jurassic.h:1299
int refrac
Take into account refractivity (0=no, 1=yes).
Definition: jurassic.h:1374
int ig_o2
Emitter index of O2.
Definition: jurassic.h:1314
double rett_zmin
Minimum altitude for temperature retrieval [km].
Definition: jurassic.h:1401
double sfsza
Solar zenith angle at the surface [deg] (-999=auto).
Definition: jurassic.h:1344
int nd
Number of radiance channels.
Definition: jurassic.h:1317
int fov_n
Field-of-view number of data points.
Definition: jurassic.h:1392
int sftype
Surface treatment (0=none, 1=emissions, 2=downward, 3=solar).
Definition: jurassic.h:1341
int ncl
Number of cloud layer spectral grid points.
Definition: jurassic.h:1329
int ig_n2
Emitter index of N2.
Definition: jurassic.h:1311
int obsfmt
Observation data file format (1=ASCII, 2=binary, 3=netCDF).
Definition: jurassic.h:1356
int ctm_o2
Compute O2 continuum (0=no, 1=yes).
Definition: jurassic.h:1371
int ret_clk
Retrieve cloud layer extinction (0=no, 1=yes).
Definition: jurassic.h:1425
int nsf
Number of surface layer spectral grid points.
Definition: jurassic.h:1335
double rayds
Maximum step length for raytracing [km].
Definition: jurassic.h:1377
int ret_cldz
Retrieve cloud layer depth (0=no, 1=yes).
Definition: jurassic.h:1422
int ig_h2o
Emitter index of H2O.
Definition: jurassic.h:1308
int tblfmt
Look-up table file format (1=ASCII, 2=binary, 3=netCDF).
Definition: jurassic.h:1350
double raydz
Vertical step length for raytracing [km].
Definition: jurassic.h:1380
int write_bbt
Use brightness temperature instead of radiance (0=no, 1=yes).
Definition: jurassic.h:1434
double retp_zmax
Maximum altitude for pressure retrieval [km].
Definition: jurassic.h:1398
Line-of-sight data.
Definition: jurassic.h:1460
double sft
Surface temperature [K].
Definition: jurassic.h:1487
int np
Number of LOS points.
Definition: jurassic.h:1463
Observation geometry and radiance data.
Definition: jurassic.h:1522
int nr
Number of ray paths.
Definition: jurassic.h:1525
Retrieval control parameters.
Definition: jurassic.h:1574
double err_press_cz
Vertical correlation length for pressure error [km].
Definition: jurassic.h:1601
double err_press
Pressure error [%].
Definition: jurassic.h:1598
int err_ana
Carry out error analysis (0=no, 1=yes).
Definition: jurassic.h:1589
double err_temp_cz
Vertical correlation length for temperature error [km].
Definition: jurassic.h:1610
double conv_dmin
Minimum normalized step size in state space.
Definition: jurassic.h:1586
double err_temp
Temperature error [K].
Definition: jurassic.h:1607
double err_clz
Cloud height error [km].
Definition: jurassic.h:1634
double err_sft
Surface temperature error [K].
Definition: jurassic.h:1643
double err_temp_ch
Horizontal correlation length for temperature error [km].
Definition: jurassic.h:1613
int kernel_recomp
Re-computation of kernel matrix (number of iterations).
Definition: jurassic.h:1580
int conv_itmax
Maximum number of iterations.
Definition: jurassic.h:1583
double err_cldz
Cloud depth error [km].
Definition: jurassic.h:1637
double err_press_ch
Horizontal correlation length for pressure error [km].
Definition: jurassic.h:1604
Emissivity look-up tables.
Definition: jurassic.h:1656