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
133#ifndef JURASSIC_H
134#define JURASSIC_H
135
136/* ------------------------------------------------------------
137 Includes...
138 ------------------------------------------------------------ */
139
140#include <errno.h>
141#include <fcntl.h>
142#include <gsl/gsl_math.h>
143#include <gsl/gsl_blas.h>
144#include <gsl/gsl_linalg.h>
145#include <gsl/gsl_randist.h>
146#include <gsl/gsl_rng.h>
147#include <gsl/gsl_statistics.h>
148#include <math.h>
149#include <netcdf.h>
150#include <omp.h>
151#include <stdint.h>
152#include <stdio.h>
153#include <stdlib.h>
154#include <string.h>
155#include <time.h>
156#include <unistd.h>
157
158/* ------------------------------------------------------------
159 Constants...
160 ------------------------------------------------------------ */
161
163#ifndef C1
164#define C1 1.19104259e-8
165#endif
166
168#ifndef C2
169#define C2 1.43877506
170#endif
171
173#ifndef EPSMIN
174#define EPSMIN 0
175#endif
176
178#ifndef EPSMAX
179#define EPSMAX 1
180#endif
181
183#ifndef G0
184#define G0 9.80665
185#endif
186
188#ifndef H0
189#define H0 7.0
190#endif
191
193#ifndef KB
194#define KB 1.3806504e-23
195#endif
196
198#ifndef ME
199#define ME 5.976e24
200#endif
201
203#ifndef NA
204#define NA 6.02214199e23
205#endif
206
208#ifndef N2
209#define N2 0.78084
210#endif
211
213#ifndef O2
214#define O2 0.20946
215#endif
216
218#ifndef P0
219#define P0 1013.25
220#endif
221
223#ifndef RE
224#define RE 6367.421
225#endif
226
228#ifndef RI
229#define RI 8.3144598
230#endif
231
233#ifndef T0
234#define T0 273.15
235#endif
236
238#ifndef TMAX
239#define TMAX 400.
240#endif
241
243#ifndef TMIN
244#define TMIN 100.
245#endif
246
248#ifndef TSUN
249#define TSUN 5780.
250#endif
251
253#ifndef UMIN
254#define UMIN 0
255#endif
256
258#ifndef UMAX
259#define UMAX 1e30
260#endif
261
262/* ------------------------------------------------------------
263 Dimensions...
264 ------------------------------------------------------------ */
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 NCL
283#define NCL 8
284#endif
285
287#ifndef ND
288#define ND 128
289#endif
290
292#ifndef NFOV
293#define NFOV 5
294#endif
295
297#ifndef NG
298#define NG 8
299#endif
300
302#ifndef NLOS
303#define NLOS 4096
304#endif
305
307#ifndef NP
308#define NP 256
309#endif
310
312#ifndef NQ
313#define NQ (5 + NG + NW + NCL + NSF)
314#endif
315
317#ifndef NR
318#define NR 256
319#endif
320
322#ifndef NSF
323#define NSF 8
324#endif
325
327#ifndef NSHAPE
328#define NSHAPE 20000
329#endif
330
332#ifndef NTIMER
333#define NTIMER 100
334#endif
335
337#ifndef NW
338#define NW 4
339#endif
340
342#ifndef RFMNPTS
343#define RFMNPTS 10000000
344#endif
345
347#ifndef TBLNP
348#define TBLNP 41
349#endif
350
352#ifndef TBLNS
353#define TBLNS 1200
354#endif
355
357#ifndef TBLNT
358#define TBLNT 30
359#endif
360
362#ifndef TBLNU
363#define TBLNU 512
364#endif
365
366/* ------------------------------------------------------------
367 Quantity indices...
368 ------------------------------------------------------------ */
369
371#define IDXP 0
372
374#define IDXT 1
375
377#define IDXQ(ig) (2 + (ig))
378
380#define IDXK(iw) (2 + (ctl->ng) + (iw))
381
383#define IDXCLZ (2 + (ctl->ng) + (ctl->nw))
384
386#define IDXCLDZ (3 + (ctl->ng) + (ctl->nw))
387
389#define IDXCLK(icl) (4 + (ctl->ng) + (ctl->nw) + (icl))
390
392#define IDXSFT (4 + (ctl->ng) + (ctl->nw) + (ctl->ncl))
393
395#define IDXSFEPS(isf) (5 + (ctl->ng) + (ctl->nw) + (ctl->ncl) + (isf))
396
397/* ------------------------------------------------------------
398 Macros...
399 ------------------------------------------------------------ */
400
418#define ALLOC(ptr, type, n) \
419 if((ptr=calloc((size_t)(n), sizeof(type)))==NULL) \
420 ERRMSG("Out of memory!");
421
448#define BRIGHT(rad, nu) \
449 (C2 * (nu) / gsl_log1p(C1 * POW3(nu) / (rad)))
450
470#define CLAMP(v, lo, hi) \
471 (((v) < (lo)) ? (lo) : (((v) > (hi)) ? (hi) : (v)))
472
485#define CO2_FIT(TIME) \
486 (373.87884731310215e-6 \
487 + 1.9021498612395775e-6 * (((TIME) - 94694400.) / 31557600.) \
488 + 0.016490044598136553e-6 \
489 * POW2((((TIME) - 94694400.) / 31557600.)))
490
505#define DEG2RAD(deg) ((deg) * (M_PI / 180.0))
506
523#define DIST(a, b) sqrt(DIST2(a, b))
524
540#define DIST2(a, b) \
541 ((a[0]-b[0])*(a[0]-b[0])+(a[1]-b[1])*(a[1]-b[1])+(a[2]-b[2])*(a[2]-b[2]))
542
558#define DOTP(a, b) (a[0]*b[0]+a[1]*b[1]+a[2]*b[2])
559
578#define FREAD(ptr, type, size, out) { \
579 if(fread(ptr, sizeof(type), size, out)!=size) \
580 ERRMSG("Error while reading!"); \
581 }
582
601#define FWRITE(ptr, type, size, out) { \
602 if(fwrite(ptr, sizeof(type), size, out)!=size) \
603 ERRMSG("Error while writing!"); \
604 }
605
624#define LIN(x0, y0, x1, y1, x) \
625 ((y0)+((y1)-(y0))/((x1)-(x0))*((x)-(x0)))
626
645#define LOGY(x0, y0, x1, y1, x) \
646 (((y1)/(y0)>0) \
647 ? ((y0)*exp(log((y1)/(y0))/((x1)-(x0))*((x)-(x0)))) \
648 : LIN(x0, y0, x1, y1, x))
649
666#define MAX(a,b) (((a)>(b))?(a):(b))
667
684#define MIN(a,b) (((a)<(b))?(a):(b))
685
697#define NC(cmd) { \
698 int nc_result=(cmd); \
699 if(nc_result!=NC_NOERR) \
700 ERRMSG("%s", nc_strerror(nc_result)); \
701 }
702
726#define NC_DEF_VAR(varname, type, ndims, dims, long_name, units, level, quant) { \
727 NC(nc_def_var(ncid, varname, type, ndims, dims, &varid)); \
728 NC(nc_put_att_text(ncid, varid, "long_name", strnlen(long_name, LEN), long_name)); \
729 NC(nc_put_att_text(ncid, varid, "units", strnlen(units, LEN), units)); \
730 if((quant) > 0) \
731 NC(nc_def_var_quantize(ncid, varid, NC_QUANTIZE_GRANULARBR, quant)); \
732 if((level) != 0) { \
733 NC(nc_def_var_deflate(ncid, varid, 1, 1, level)); \
734 /* unsigned int ulevel = (unsigned int)level; */ \
735 /* NC(nc_def_var_filter(ncid, varid, 32015, 1, (unsigned int[]){ulevel})); */ \
736 } \
737 }
738
756#define NC_GET_DOUBLE(varname, ptr, force) { \
757 if(force) { \
758 NC(nc_inq_varid(ncid, varname, &varid)); \
759 NC(nc_get_var_double(ncid, varid, ptr)); \
760 } else { \
761 if(nc_inq_varid(ncid, varname, &varid) == NC_NOERR) { \
762 NC(nc_get_var_double(ncid, varid, ptr)); \
763 } else \
764 WARN("netCDF variable %s is missing!", varname); \
765 } \
766 }
767
785#define NC_INQ_DIM(dimname, ptr, min, max, check) { \
786 int dimid; size_t naux; \
787 NC(nc_inq_dimid(ncid, dimname, &dimid)); \
788 NC(nc_inq_dimlen(ncid, dimid, &naux)); \
789 *ptr = (int)naux; \
790 if (check) \
791 if ((*ptr) < (min) || (*ptr) > (max)) \
792 ERRMSG("Dimension %s is out of range!", dimname); \
793 }
794
809#define NC_PUT_DOUBLE(varname, ptr, hyperslab) { \
810 NC(nc_inq_varid(ncid, varname, &varid)); \
811 if(hyperslab) { \
812 NC(nc_put_vara_double(ncid, varid, start, count, ptr)); \
813 } else { \
814 NC(nc_put_var_double(ncid, varid, ptr)); \
815 } \
816 }
817
833#define NC_PUT_FLOAT(varname, ptr, hyperslab) { \
834 NC(nc_inq_varid(ncid, varname, &varid)); \
835 if(hyperslab) { \
836 NC(nc_put_vara_float(ncid, varid, start, count, ptr)); \
837 } else { \
838 NC(nc_put_var_float(ncid, varid, ptr)); \
839 } \
840 }
841
856#define NC_PUT_INT(varname, ptr, hyperslab) { \
857 NC(nc_inq_varid(ncid, varname, &varid)); \
858 if(hyperslab) { \
859 NC(nc_put_vara_int(ncid, varid, start, count, ptr)); \
860 } else { \
861 NC(nc_put_var_int(ncid, varid, ptr)); \
862 } \
863 }
864
878#define NC_PUT_ATT(varname, attname, text) { \
879 NC(nc_inq_varid(ncid, varname, &varid)); \
880 NC(nc_put_att_text(ncid, varid, attname, strnlen(text, LEN), text)); \
881 }
882
895#define NC_PUT_ATT_GLOBAL(attname, text) \
896 NC(nc_put_att_text(ncid, NC_GLOBAL, attname, strnlen(text, LEN), text));
897
927#define NEDT(t_bg, nesr, nu) \
928 (BRIGHT(PLANCK((t_bg), (nu)) + (nesr), (nu)) - (t_bg))
929
956#define NESR(t_bg, nedt, nu) \
957 (PLANCK((t_bg) + (nedt), (nu)) - PLANCK((t_bg), (nu)))
958
973#define NORM(a) sqrt(DOTP(a, a))
974
1004#define PLANCK(T, nu) \
1005 (C1 * POW3(nu) / gsl_expm1(C2 * (nu) / (T)))
1006
1018#define POW2(x) ((x)*(x))
1019
1031#define POW3(x) ((x)*(x)*(x))
1032
1047#define RAD2DEG(rad) ((rad) * (180.0 / M_PI))
1048
1077#define REFRAC_SW53(p, T, h2o) \
1078 (1e-6 * (77.6 * (p) / (T) + 3.73e5 * (h2o) * (p) / ((T) * (T))))
1079
1104#define REFRAC_EDLEN(p, T) \
1105 ((1e-8 \
1106 * (8342.54 + 2406147.0 / (130.0 - 1.0 / (10.0 * 10.0)) \
1107 + 15998.0 / (38.9 - 1.0 / (10.0 * 10.0)))) \
1108 * ((p) / 1013.25) * (288.15 / (T)))
1109
1124#define TBL_LOG(tbl, id, ig) \
1125 do { \
1126 for (int ip = 0; ip < (tbl)->np[(id)][(ig)]; ip++) \
1127 LOG(2, \
1128 "p[%2d]= %.5e hPa | T[0:%2d]= %.2f ... %.2f K | " \
1129 "u[0:%3d]= %.5e ... %.5e molec/cm^2 | " \
1130 "eps[0:%3d]= %.5e ... %.5e", \
1131 (ip), \
1132 (tbl)->p[(id)][(ig)][(ip)], \
1133 (tbl)->nt[(id)][(ig)][(ip)] - 1, \
1134 (tbl)->t[(id)][(ig)][(ip)][0], \
1135 (tbl)->t[(id)][(ig)][(ip)] \
1136 [(tbl)->nt[(id)][(ig)][(ip)] - 1], \
1137 (tbl)->nu[(id)][(ig)][(ip)][0] - 1, \
1138 exp((tbl)->logu[(id)][(ig)][(ip)][0][0]), \
1139 exp((tbl)->logu[(id)][(ig)][(ip)][0] \
1140 [(tbl)->nu[(id)][(ig)][(ip)][0] - 1]), \
1141 (tbl)->nu[(id)][(ig)][(ip)][0] - 1, \
1142 exp((tbl)->logeps[(id)][(ig)][(ip)][0][0]), \
1143 exp((tbl)->logeps[(id)][(ig)][(ip)][0] \
1144 [(tbl)->nu[(id)][(ig)][(ip)][0] - 1])); \
1145 } while (0)
1146
1158#define SELECT_TIMER(id, group) \
1159 {timer_group(id, group, 0);}
1160
1168#define PRINT_TIMERS \
1169 {timer_group("END", "END", 1);}
1170
1189#define TOK(line, tok, format, var) { \
1190 if(((tok)=strtok((line), " \t"))) { \
1191 if(sscanf(tok, format, &(var))!=1) continue; \
1192 } else ERRMSG("Error while reading!"); \
1193 }
1194
1206#define USAGE \
1207 do { \
1208 int iusage; \
1209 for (iusage = 1; iusage < argc; iusage++) \
1210 if (!strcmp(argv[iusage], "-h") \
1211 || !strcmp(argv[iusage], "--help")) { \
1212 usage(); \
1213 return EXIT_SUCCESS; \
1214 } \
1215 } while (0)
1216
1217/* ------------------------------------------------------------
1218 Log messages...
1219 ------------------------------------------------------------ */
1220
1222#ifndef LOGLEV
1223#define LOGLEV 2
1224#endif
1225
1255#define LOG(level, ...) { \
1256 if(level >= 2) \
1257 printf(" "); \
1258 if(level <= LOGLEV) { \
1259 printf(__VA_ARGS__); \
1260 printf("\n"); \
1261 } \
1262 }
1263
1292#define WARN(...) { \
1293 printf("\nWarning (%s, %s, l%d): ", __FILE__, __func__, __LINE__); \
1294 LOG(0, __VA_ARGS__); \
1295 }
1296
1325#define ERRMSG(...) { \
1326 printf("\nError (%s, %s, l%d): ", __FILE__, __func__, __LINE__); \
1327 LOG(0, __VA_ARGS__); \
1328 exit(EXIT_FAILURE); \
1329 }
1330
1360#define PRINT(format, var) \
1361 printf("Print (%s, %s, l%d): %s= "format"\n", \
1362 __FILE__, __func__, __LINE__, #var, var);
1363
1364/* ------------------------------------------------------------
1365 Structs...
1366 ------------------------------------------------------------ */
1367
1375typedef struct {
1376
1378 int np;
1379
1381 double time[NP];
1382
1384 double z[NP];
1385
1387 double lon[NP];
1388
1390 double lat[NP];
1391
1393 double p[NP];
1394
1396 double t[NP];
1397
1399 double q[NG][NP];
1400
1402 double k[NW][NP];
1403
1405 double clz;
1406
1408 double cldz;
1409
1411 double clk[NCL];
1412
1414 double sft;
1415
1417 double sfeps[NSF];
1418
1419} atm_t;
1420
1428typedef struct {
1429
1431 int ng;
1432
1434 char emitter[NG][LEN];
1435
1438
1441
1444
1447
1449 int nd;
1450
1452 double nu[ND];
1453
1455 int nw;
1456
1458 int window[ND];
1459
1461 int ncl;
1462
1464 double clnu[NCL];
1465
1467 int nsf;
1468
1470 double sfnu[NSF];
1471
1474
1476 double sfsza;
1477
1479 char tblbase[LEN];
1480
1483
1486
1489
1492
1494 double hydz;
1495
1498
1501
1504
1507
1510
1512 double rayds;
1513
1515 double raydz;
1516
1518 char fov[LEN];
1519
1521 double fov_dz[NSHAPE];
1522
1524 double fov_w[NSHAPE];
1525
1528
1531
1534
1537
1540
1542 double retq_zmin[NG];
1543
1545 double retq_zmax[NG];
1546
1548 double retk_zmin[NW];
1549
1551 double retk_zmax[NW];
1552
1555
1558
1561
1564
1567
1570
1573
1576
1578 char rfmbin[LEN];
1579
1581 char rfmhit[LEN];
1582
1584 char rfmxsc[NG][LEN];
1585
1586} ctl_t;
1587
1595typedef struct {
1596
1598 int np;
1599
1601 double z[NLOS];
1602
1604 double lon[NLOS];
1605
1607 double lat[NLOS];
1608
1610 double p[NLOS];
1611
1613 double t[NLOS];
1614
1616 double q[NLOS][NG];
1617
1619 double k[NLOS][ND];
1620
1622 double sft;
1623
1625 double sfeps[ND];
1626
1628 double ds[NLOS];
1629
1631 double u[NLOS][NG];
1632
1634 double cgp[NLOS][NG];
1635
1637 double cgt[NLOS][NG];
1638
1640 double cgu[NLOS][NG];
1641
1643 double eps[NLOS][ND];
1644
1646 double src[NLOS][ND];
1647
1648} los_t;
1649
1657typedef struct {
1658
1660 int nr;
1661
1663 double time[NR];
1664
1666 double obsz[NR];
1667
1669 double obslon[NR];
1670
1672 double obslat[NR];
1673
1675 double vpz[NR];
1676
1678 double vplon[NR];
1679
1681 double vplat[NR];
1682
1684 double tpz[NR];
1685
1687 double tplon[NR];
1688
1690 double tplat[NR];
1691
1693 double tau[ND][NR];
1694
1696 double rad[ND][NR];
1697
1698} obs_t;
1699
1709typedef struct {
1710
1712 char dir[LEN];
1713
1716
1719
1722
1725
1727 double err_formod[ND];
1728
1730 double err_noise[ND];
1731
1734
1737
1740
1742 double err_temp;
1743
1746
1749
1751 double err_q[NG];
1752
1754 double err_q_cz[NG];
1755
1757 double err_q_ch[NG];
1758
1760 double err_k[NW];
1761
1763 double err_k_cz[NW];
1764
1766 double err_k_ch[NW];
1767
1769 double err_clz;
1770
1772 double err_cldz;
1773
1775 double err_clk[NCL];
1776
1778 double err_sft;
1779
1781 double err_sfeps[NSF];
1782
1785
1787 char shared_io_proflist[LEN];
1788
1790 char shared_io_atm_apr_file[LEN];
1791
1793 char shared_io_obs_meas_file[LEN];
1794
1796 char shared_io_atm_final_file[LEN];
1797
1799 char shared_io_obs_final_file[LEN];
1800
1802 char shared_io_matrix_cov_apr_file[LEN];
1803
1805 char shared_io_matrix_kernel_file[LEN];
1806
1808 char shared_io_matrix_cov_ret_file[LEN];
1809
1811 char shared_io_matrix_corr_file[LEN];
1812
1814 char shared_io_matrix_gain_file[LEN];
1815
1817 char shared_io_matrix_avk_file[LEN];
1818
1820 char shared_io_atm_err_total_file[LEN];
1821
1823 char shared_io_atm_err_noise_file[LEN];
1824
1826 char shared_io_atm_err_formod_file[LEN];
1827
1829 char shared_io_atm_cont_file[LEN];
1830
1832 char shared_io_atm_res_file[LEN];
1833
1834} ret_t;
1835
1842typedef struct {
1843
1845 int np[ND][NG];
1846
1848 int nt[ND][NG][TBLNP];
1849
1851 int nu[ND][NG][TBLNP][TBLNT];
1852
1854 double p[ND][NG][TBLNP];
1855
1857 double lnp[ND][NG][TBLNP];
1858
1860 double t[ND][NG][TBLNP][TBLNT];
1861
1863 float *logu[ND][NG][TBLNP][TBLNT];
1864
1866 float *logeps[ND][NG][TBLNP][TBLNT];
1867
1869 int filt_n[ND];
1870
1872 double filt_nu[ND][NSHAPE];
1873
1875 double filt_f[ND][NSHAPE];
1876
1878 double st[TBLNS];
1879
1881 double sr[TBLNS][ND];
1882
1883} tbl_t;
1884
1885/* ------------------------------------------------------------
1886 Functions...
1887 ------------------------------------------------------------ */
1888
1940void analyze_avk(
1941 const ret_t * ret,
1942 const ctl_t * ctl,
1943 const atm_t * atm,
1944 const int *iqa,
1945 const int *ipa,
1946 const gsl_matrix * avk);
1947
2004 const gsl_matrix * avk,
2005 const int iq,
2006 const int *ipa,
2007 const size_t *n0,
2008 const size_t *n1,
2009 double *cont,
2010 double *res);
2011
2038size_t atm2x(
2039 const ctl_t * ctl,
2040 const atm_t * atm,
2041 gsl_vector * x,
2042 int *iqa,
2043 int *ipa);
2044
2066void atm2x_help(
2067 const double value,
2068 const int value_iqa,
2069 const int value_ip,
2070 gsl_vector * x,
2071 int *iqa,
2072 int *ipa,
2073 size_t *n);
2074
2089void cart2geo(
2090 const double *x,
2091 double *z,
2092 double *lon,
2093 double *lat);
2094
2110void climatology(
2111 const ctl_t * ctl,
2112 atm_t * atm);
2113
2138double cos_sza(
2139 const double sec,
2140 const double lon,
2141 const double lat);
2142
2182double cost_function(
2183 const gsl_vector * dx,
2184 const gsl_vector * dy,
2185 const gsl_matrix * s_a_inv,
2186 const gsl_vector * sig_eps_inv);
2187
2193double ctmco2(
2194 const double nu,
2195 const double p,
2196 const double t,
2197 const double u);
2198
2204double ctmh2o(
2205 const double nu,
2206 const double p,
2207 const double t,
2208 const double q,
2209 const double u);
2210
2240double ctmn2(
2241 const double nu,
2242 const double p,
2243 const double t);
2244
2274double ctmo2(
2275 const double nu,
2276 const double p,
2277 const double t);
2278
2303void copy_atm(
2304 const ctl_t * ctl,
2305 atm_t * atm_dest,
2306 const atm_t * atm_src,
2307 const int init);
2308
2333void copy_obs(
2334 const ctl_t * ctl,
2335 obs_t * obs_dest,
2336 const obs_t * obs_src,
2337 const int init);
2338
2359void day2doy(
2360 int year,
2361 int mon,
2362 int day,
2363 int *doy);
2364
2387void doy2day(
2388 int year,
2389 int doy,
2390 int *mon,
2391 int *day);
2392
2412int find_emitter(
2413 const ctl_t * ctl,
2414 const char *emitter);
2415
2441void formod(
2442 const ctl_t * ctl,
2443 const tbl_t * tbl,
2444 atm_t * atm,
2445 obs_t * obs);
2446
2472void formod_continua(
2473 const ctl_t * ctl,
2474 const los_t * los,
2475 const int ip,
2476 double *beta);
2477
2499void formod_fov(
2500 const ctl_t * ctl,
2501 obs_t * obs);
2502
2532void formod_pencil(
2533 const ctl_t * ctl,
2534 const tbl_t * tbl,
2535 const atm_t * atm,
2536 obs_t * obs,
2537 const int ir);
2538
2590void formod_rfm(
2591 const ctl_t * ctl,
2592 const tbl_t * tbl,
2593 const atm_t * atm,
2594 obs_t * obs);
2595
2620void formod_srcfunc(
2621 const ctl_t * ctl,
2622 const tbl_t * tbl,
2623 const double t,
2624 double *src);
2625
2652void geo2cart(
2653 const double z,
2654 const double lon,
2655 const double lat,
2656 double *x);
2657
2685void hydrostatic(
2686 const ctl_t * ctl,
2687 atm_t * atm);
2688
2712void idx2name(
2713 const ctl_t * ctl,
2714 const int idx,
2715 char *quantity);
2716
2746void init_srcfunc(
2747 const ctl_t * ctl,
2748 tbl_t * tbl);
2749
2774void intpol_atm(
2775 const ctl_t * ctl,
2776 const atm_t * atm,
2777 const double z,
2778 double *p,
2779 double *t,
2780 double *q,
2781 double *k);
2782
2816void intpol_tbl_cga(
2817 const ctl_t * ctl,
2818 const tbl_t * tbl,
2819 const los_t * los,
2820 const int ip,
2821 double tau_path[ND][NG],
2822 double tau_seg[ND]);
2823
2862void intpol_tbl_ega(
2863 const ctl_t * ctl,
2864 const tbl_t * tbl,
2865 const los_t * los,
2866 const int ip,
2867 double tau_path[ND][NG],
2868 double tau_seg[ND]);
2869
2896double intpol_tbl_eps(
2897 const tbl_t * tbl,
2898 const int ig,
2899 const int id,
2900 const int ip,
2901 const int it,
2902 const double logu);
2903
2934double intpol_tbl_u(
2935 const tbl_t * tbl,
2936 const int ig,
2937 const int id,
2938 const int ip,
2939 const int it,
2940 const double logeps);
2941
2968void jsec2time(
2969 const double jsec,
2970 int *year,
2971 int *mon,
2972 int *day,
2973 int *hour,
2974 int *min,
2975 int *sec,
2976 double *remain);
2977
3019void kernel(
3020 const ctl_t * ctl,
3021 const tbl_t * tbl,
3022 atm_t * atm,
3023 obs_t * obs,
3024 gsl_matrix * k);
3025
3052int locate_irr(
3053 const double *xx,
3054 const int n,
3055 const double x);
3056
3083int locate_reg(
3084 const double *xx,
3085 const int n,
3086 const double x);
3087
3113int locate_tbl(
3114 const float *xx,
3115 const int n,
3116 const double x);
3117
3156void matrix_invert(
3157 gsl_matrix * a);
3158
3203void matrix_product(
3204 const gsl_matrix * a,
3205 const gsl_vector * b,
3206 const int transpose,
3207 gsl_matrix * c);
3208
3236size_t obs2y(
3237 const ctl_t * ctl,
3238 const obs_t * obs,
3239 gsl_vector * y,
3240 int *ida,
3241 int *ira);
3242
3294 ret_t * ret,
3295 ctl_t * ctl,
3296 tbl_t * tbl,
3297 obs_t * obs_meas,
3298 obs_t * obs_i,
3299 atm_t * atm_apr,
3300 atm_t * atm_i,
3301 double *chisq);
3302
3344void raytrace(
3345 const ctl_t * ctl,
3346 const atm_t * atm,
3347 obs_t * obs,
3348 los_t * los,
3349 const int ir);
3350
3400void read_atm(
3401 const char *dirname,
3402 const char *filename,
3403 const ctl_t * ctl,
3404 atm_t * atm,
3405 int profile);
3406
3453void read_atm_asc(
3454 const char *filename,
3455 const ctl_t * ctl,
3456 atm_t * atm);
3457
3507void read_atm_bin(
3508 const char *filename,
3509 const ctl_t * ctl,
3510 atm_t * atm);
3511
3553void read_atm_nc(
3554 const char *filename,
3555 const ctl_t * ctl,
3556 atm_t * atm,
3557 int profile);
3558
3599void read_ctl(
3600 int argc,
3601 char *argv[],
3602 ctl_t * ctl);
3603
3622void read_matrix(
3623 const char *dirname,
3624 const char *filename,
3625 const ctl_t * ctl,
3626 gsl_matrix * matrix,
3627 int dataset);
3628
3660void read_matrix_asc(
3661 const char *dirname,
3662 const char *filename,
3663 gsl_matrix * matrix);
3664
3682void read_matrix_bin(
3683 const char *dirname,
3684 const char *filename,
3685 gsl_matrix * matrix);
3686
3708void read_matrix_nc(
3709 const char *dirname,
3710 const char *filename,
3711 gsl_matrix * matrix,
3712 int dataset);
3713
3750void read_obs(
3751 const char *dirname,
3752 const char *filename,
3753 const ctl_t * ctl,
3754 obs_t * obs,
3755 int profile);
3756
3794void read_obs_asc(
3795 const char *filename,
3796 const ctl_t * ctl,
3797 obs_t * obs);
3798
3840void read_obs_bin(
3841 const char *filename,
3842 const ctl_t * ctl,
3843 obs_t * obs);
3844
3877void read_obs_nc(
3878 const char *filename,
3879 const ctl_t * ctl,
3880 obs_t * obs,
3881 const int profile);
3882
3922double read_obs_rfm(
3923 const char *basename,
3924 const double z,
3925 const double *nu,
3926 const double *f,
3927 const int n);
3928
4003void read_ret(
4004 int argc,
4005 char *argv[],
4006 const ctl_t * ctl,
4007 ret_t * ret);
4008
4049void read_rfm_spec(
4050 const char *filename,
4051 double *nu,
4052 double *rad,
4053 int *npts);
4054
4088void read_shape(
4089 const char *filename,
4090 double *x,
4091 double *y,
4092 int *n);
4093
4130 const ctl_t * ctl);
4131
4167void read_tbl_asc(
4168 const ctl_t * ctl,
4169 tbl_t * tbl,
4170 const int id,
4171 const int ig);
4172
4204void read_tbl_bin(
4205 const ctl_t * ctl,
4206 tbl_t * tbl,
4207 const int id,
4208 const int ig);
4209
4253 const ctl_t * ctl,
4254 tbl_t * tbl,
4255 int id,
4256 int ig,
4257 int ncid);
4258
4313double scan_ctl(
4314 int argc,
4315 char *argv[],
4316 const char *varname,
4317 const int arridx,
4318 const char *defvalue,
4319 char *value);
4320
4374void set_cov_apr(
4375 const ret_t * ret,
4376 const ctl_t * ctl,
4377 const atm_t * atm,
4378 const int *iqa,
4379 const int *ipa,
4380 gsl_matrix * s_a);
4381
4435void set_cov_meas(
4436 const ret_t * ret,
4437 const ctl_t * ctl,
4438 const obs_t * obs,
4439 gsl_vector * sig_noise,
4440 gsl_vector * sig_formod,
4441 gsl_vector * sig_eps_inv);
4442
4469const char *shared_io_input_target(
4470 const ret_t * ret,
4471 const char *shared_file,
4472 const char *legacy_file,
4473 const char **dirname,
4474 int *profile);
4475
4502const char *shared_io_output_target(
4503 const ret_t * ret,
4504 const char *shared_file,
4505 const char *legacy_file,
4506 const char **dirname,
4507 int *profile);
4508
4529int shared_io_lock(
4530 const ret_t * ret);
4531
4546void shared_io_unlock(
4547 int fd);
4548
4569const char *shared_io_output_file(
4570 const ret_t * ret);
4571
4611void tangent_point(
4612 const los_t * los,
4613 double *tpz,
4614 double *tplon,
4615 double *tplat);
4616
4638void tbl_free(
4639 const ctl_t * ctl,
4640 tbl_t * tbl);
4641
4681void tbl_pack(
4682 const tbl_t * tbl,
4683 int id,
4684 int ig,
4685 uint8_t * buf,
4686 size_t *bytes_used);
4687
4698size_t tbl_packed_size(
4699 const tbl_t * tbl,
4700 int id,
4701 int ig);
4702
4741size_t tbl_unpack(
4742 tbl_t * tbl,
4743 int id,
4744 int ig,
4745 const uint8_t * buf);
4746
4747
4772void time2jsec(
4773 const int year,
4774 const int mon,
4775 const int day,
4776 const int hour,
4777 const int min,
4778 const int sec,
4779 const double remain,
4780 double *jsec);
4781
4796void timer_group(
4797 const char *name,
4798 const char *group,
4799 int output);
4800
4849void write_atm(
4850 const char *dirname,
4851 const char *filename,
4852 const ctl_t * ctl,
4853 const atm_t * atm,
4854 int profile);
4855
4907void write_atm_asc(
4908 const char *filename,
4909 const ctl_t * ctl,
4910 const atm_t * atm);
4911
4965void write_atm_bin(
4966 const char *filename,
4967 const ctl_t * ctl,
4968 const atm_t * atm);
4969
5017void write_atm_nc(
5018 const char *filename,
5019 const ctl_t * ctl,
5020 const atm_t * atm,
5021 int profile);
5022
5074void write_atm_rfm(
5075 const char *filename,
5076 const ctl_t * ctl,
5077 const atm_t * atm);
5078
5102void write_matrix(
5103 const char *dirname,
5104 const char *filename,
5105 const ctl_t * ctl,
5106 const gsl_matrix * matrix,
5107 const atm_t * atm,
5108 const obs_t * obs,
5109 const char *rowspace,
5110 const char *colspace,
5111 const char *sort,
5112 int dataset);
5113
5161void write_matrix_asc(
5162 const char *dirname,
5163 const char *filename,
5164 const ctl_t * ctl,
5165 const gsl_matrix * matrix,
5166 const atm_t * atm,
5167 const obs_t * obs,
5168 const char *rowspace,
5169 const char *colspace,
5170 const char *sort);
5171
5186void write_matrix_bin(
5187 const char *dirname,
5188 const char *filename,
5189 const gsl_matrix * matrix);
5190
5217void write_matrix_nc(
5218 const char *dirname,
5219 const char *filename,
5220 const ctl_t * ctl,
5221 const gsl_matrix * matrix,
5222 const atm_t * atm,
5223 const obs_t * obs,
5224 const char *rowspace,
5225 const char *colspace,
5226 const char *sort,
5227 int dataset);
5228
5264void write_obs(
5265 const char *dirname,
5266 const char *filename,
5267 const ctl_t * ctl,
5268 const obs_t * obs,
5269 int profile);
5270
5309void write_obs_asc(
5310 const char *filename,
5311 const ctl_t * ctl,
5312 const obs_t * obs);
5313
5351void write_obs_bin(
5352 const char *filename,
5353 const ctl_t * ctl,
5354 const obs_t * obs);
5355
5393void write_obs_nc(
5394 const char *filename,
5395 const ctl_t * ctl,
5396 const obs_t * obs,
5397 const int profile);
5398
5433void write_shape(
5434 const char *filename,
5435 const double *x,
5436 const double *y,
5437 const int n);
5438
5484void write_stddev(
5485 const char *quantity,
5486 const ret_t * ret,
5487 const ctl_t * ctl,
5488 const atm_t * atm,
5489 const gsl_matrix * s);
5490
5526void write_tbl(
5527 const ctl_t * ctl,
5528 const tbl_t * tbl);
5529
5555void write_tbl_asc(
5556 const ctl_t * ctl,
5557 const tbl_t * tbl,
5558 const int id,
5559 const int ig);
5560
5588void write_tbl_bin(
5589 const ctl_t * ctl,
5590 const tbl_t * tbl,
5591 const int id,
5592 const int ig);
5593
5626void write_tbl_nc(
5627 const ctl_t * ctl,
5628 const tbl_t * tbl,
5629 const int id,
5630 const int ig);
5631
5678void x2atm(
5679 const ctl_t * ctl,
5680 const gsl_vector * x,
5681 atm_t * atm);
5682
5711void x2atm_help(
5712 double *value,
5713 const gsl_vector * x,
5714 size_t *n);
5715
5748void y2obs(
5749 const ctl_t * ctl,
5750 const gsl_vector * y,
5751 obs_t * obs);
5752
5753#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:4263
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:6531
#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:1193
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:6086
void write_atm(const char *dirname, const char *filename, const ctl_t *ctl, const atm_t *atm, int profile)
Write atmospheric data to a file.
Definition: jurassic.c:7307
void read_rfm_spec(const char *filename, double *nu, double *rad, int *npts)
Read a Reference Forward Model (RFM) ASCII spectrum.
Definition: jurassic.c:6231
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:7041
void read_matrix_nc(const char *dirname, const char *filename, gsl_matrix *matrix, int dataset)
Read a numerical matrix from a netCDF file.
Definition: jurassic.c:5761
double cos_sza(const double sec, const double lon, const double lat)
Calculates the cosine of the solar zenith angle.
Definition: jurassic.c:1152
void read_atm_bin(const char *filename, const ctl_t *ctl, atm_t *atm)
Read atmospheric data in binary format.
Definition: jurassic.c:5318
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:8346
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:3707
int locate_reg(const double *xx, const int n, const double x)
Locate index for interpolation on a regular (uniform) grid.
Definition: jurassic.c:4511
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:4164
double ctmo2(const double nu, const double p, const double t)
Compute O₂ collision-induced absorption coefficient.
Definition: jurassic.c:3208
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:8470
void tbl_free(const ctl_t *ctl, tbl_t *tbl)
Free lookup table and all internally allocated memory.
Definition: jurassic.c:7002
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:8765
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:8849
void idx2name(const ctl_t *ctl, const int idx, char *quantity)
Convert a quantity index to a descriptive name string.
Definition: jurassic.c:3963
void read_ctl(int argc, char *argv[], ctl_t *ctl)
Read model control parameters from command-line and configuration input.
Definition: jurassic.c:5516
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:3470
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:5929
int shared_io_lock(const ret_t *ret)
Acquire an exclusive lock for shared retrieval output files.
Definition: jurassic.c:6880
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:4965
void read_obs(const char *dirname, const char *filename, const ctl_t *ctl, obs_t *obs, int profile)
Read observation data from an input file.
Definition: jurassic.c:5810
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:4582
int locate_irr(const double *xx, const int n, const double x)
Locate index for interpolation on an irregular grid.
Definition: jurassic.c:4481
const char * shared_io_input_target(const ret_t *ret, const char *shared_file, const char *legacy_file, const char **dirname, int *profile)
Resolve retrieval input target for legacy-directory and shared-file modes.
Definition: jurassic.c:6840
void write_matrix_asc(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 an ASCII file.
Definition: jurassic.c:7790
void day2doy(int year, int mon, int day, int *doy)
Convert a calendar date to day-of-year.
Definition: jurassic.c:3358
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:7529
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, int dataset)
Write a fully annotated matrix (e.g., Jacobian or gain matrix) to file.
Definition: jurassic.c:7753
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:7205
void x2atm(const ctl_t *ctl, const gsl_vector *x, atm_t *atm)
Map retrieval state vector back to atmospheric structure.
Definition: jurassic.c:8914
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:179
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:2087
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:8807
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:4074
void shared_io_unlock(int fd)
Release a shared retrieval output lock.
Definition: jurassic.c:6907
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:7716
#define ND
Maximum number of radiance channels.
Definition: jurassic.h:288
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:6791
void write_obs(const char *dirname, const char *filename, const ctl_t *ctl, const obs_t *obs, int profile)
Write observation data to an output file in ASCII or binary format.
Definition: jurassic.c:8275
#define NSHAPE
Maximum number of shape function grid points.
Definition: jurassic.h:328
void write_tbl(const ctl_t *ctl, const tbl_t *tbl)
Write emissivity lookup tables to disk.
Definition: jurassic.c:8722
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:7129
void read_matrix(const char *dirname, const char *filename, const ctl_t *ctl, gsl_matrix *matrix, int dataset)
Read a numerical matrix from a file.
Definition: jurassic.c:5643
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:4049
int find_emitter(const ctl_t *ctl, const char *emitter)
Find gas species index by name.
Definition: jurassic.c:3408
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:4654
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:4308
void timer_group(const char *name, const char *group, int output)
Aggregate wall-clock timings for named code phases.
Definition: jurassic.c:7236
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:6959
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:3570
double ctmco2(const double nu, const double p, const double t, const double u)
Compute carbon dioxide continuum (optical depth).
Definition: jurassic.c:1224
void read_atm(const char *dirname, const char *filename, const ctl_t *ctl, atm_t *atm, int profile)
Read atmospheric input data from a file.
Definition: jurassic.c:5193
void write_matrix_nc(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, int dataset)
Write a numerical matrix to a netCDF file.
Definition: jurassic.c:8012
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:7378
void write_matrix_bin(const char *dirname, const char *filename, const gsl_matrix *matrix)
Write a numerical matrix to a binary file.
Definition: jurassic.c:7963
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:3320
#define TBLNT
Maximum number of temperatures in emissivity tables.
Definition: jurassic.h:358
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:8402
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:8964
void cart2geo(const double *x, double *z, double *lon, double *lat)
Converts Cartesian coordinates to geographic coordinates.
Definition: jurassic.c:200
#define NP
Maximum number of atmospheric data points.
Definition: jurassic.h:308
void formod_fov(const ctl_t *ctl, obs_t *obs)
Apply field-of-view (FOV) convolution to modeled radiances.
Definition: jurassic.c:3506
void doy2day(int year, int doy, int *mon, int *day)
Convert a day-of-year value to a calendar date.
Definition: jurassic.c:3378
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:4387
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:6007
void init_srcfunc(const ctl_t *ctl, tbl_t *tbl)
Initialize source function lookup tables from emissivity data.
Definition: jurassic.c:4002
void matrix_invert(gsl_matrix *a)
Invert a square matrix, optimized for diagonal or symmetric positive-definite matrices.
Definition: jurassic.c:4552
void read_ret(int argc, char *argv[], const ctl_t *ctl, ret_t *ret)
Read retrieval configuration and error parameters.
Definition: jurassic.c:6143
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:6612
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:8650
void formod(const ctl_t *ctl, const tbl_t *tbl, atm_t *atm, obs_t *obs)
Execute the selected forward model.
Definition: jurassic.c:3421
void read_atm_asc(const char *filename, const ctl_t *ctl, atm_t *atm)
Read atmospheric data in ASCII format.
Definition: jurassic.c:5267
#define NG
Maximum number of emitters.
Definition: jurassic.h:298
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:6408
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:3868
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:101
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:4354
#define TBLNS
Maximum number of source function temperature levels.
Definition: jurassic.h:353
void read_matrix_asc(const char *dirname, const char *filename, gsl_matrix *matrix)
Read a numerical matrix from an ASCII file.
Definition: jurassic.c:5669
tbl_t * read_tbl(const ctl_t *ctl)
Read emissivity lookup tables from disk.
Definition: jurassic.c:6332
void read_matrix_bin(const char *dirname, const char *filename, gsl_matrix *matrix)
Read a numerical matrix from a binary file.
Definition: jurassic.c:5708
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:6575
const char * shared_io_output_target(const ret_t *ret, const char *shared_file, const char *legacy_file, const char **dirname, int *profile)
Resolve retrieval output target for legacy-directory and shared-file modes.
Definition: jurassic.c:6860
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:5412
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:3270
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:4627
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:8976
void hydrostatic(const ctl_t *ctl, atm_t *atm)
Adjust pressure profile using the hydrostatic equation.
Definition: jurassic.c:3903
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:6291
double ctmn2(const double nu, const double p, const double t)
Compute N₂ collision-induced absorption coefficient.
Definition: jurassic.c:3139
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:125
#define TBLNP
Maximum number of pressure levels in emissivity tables.
Definition: jurassic.h:348
void climatology(const ctl_t *ctl, atm_t *atm)
Initializes atmospheric climatology profiles.
Definition: jurassic.c:215
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:3885
#define NSF
Maximum number of surface layer spectral grid points.
Definition: jurassic.h:323
#define NCL
Maximum number of cloud layer spectral grid points.
Definition: jurassic.h:283
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:6680
const char * shared_io_output_file(const ret_t *ret)
Return the primary shared retrieval output file configured for locking.
Definition: jurassic.c:6925
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:8679
#define NLOS
Maximum number of LOS points.
Definition: jurassic.h:303
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:7095
#define NR
Maximum number of ray paths.
Definition: jurassic.h:318
int locate_tbl(const float *xx, const int n, const double x)
Locate index for interpolation within emissivity table grids.
Definition: jurassic.c:4530
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:7445
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:5885
#define NW
Maximum number of spectral windows.
Definition: jurassic.h:338
Atmospheric profile data.
Definition: jurassic.h:1375
double clz
Cloud layer height [km].
Definition: jurassic.h:1405
int np
Number of data points.
Definition: jurassic.h:1378
double cldz
Cloud layer depth [km].
Definition: jurassic.h:1408
double sft
Surface temperature [K].
Definition: jurassic.h:1414
Control parameters.
Definition: jurassic.h:1428
int write_matrix
Write matrix file (0=no, 1=yes).
Definition: jurassic.h:1572
int nw
Number of spectral windows.
Definition: jurassic.h:1455
int atmfmt
Atmospheric data file format (1=ASCII, 2=binary, 3=netCDF).
Definition: jurassic.h:1485
double retp_zmin
Minimum altitude for pressure retrieval [km].
Definition: jurassic.h:1530
int ig_co2
Emitter index of CO2.
Definition: jurassic.h:1437
double hydz
Reference height for hydrostatic pressure profile (-999 to skip) [km].
Definition: jurassic.h:1494
int ctm_co2
Compute CO2 continuum (0=no, 1=yes).
Definition: jurassic.h:1497
double rett_zmax
Maximum altitude for temperature retrieval [km].
Definition: jurassic.h:1539
int ret_sfeps
Retrieve surface layer emissivity (0=no, 1=yes).
Definition: jurassic.h:1566
int ret_sft
Retrieve surface layer temperature (0=no, 1=yes).
Definition: jurassic.h:1563
int ret_clz
Retrieve cloud layer height (0=no, 1=yes).
Definition: jurassic.h:1554
int ctm_n2
Compute N2 continuum (0=no, 1=yes).
Definition: jurassic.h:1503
int ctm_h2o
Compute H2O continuum (0=no, 1=yes).
Definition: jurassic.h:1500
int formod
Forward model (0=CGA, 1=EGA, 2=RFM).
Definition: jurassic.h:1575
int ng
Number of emitters.
Definition: jurassic.h:1431
int refrac
Take into account refractivity (0=no, 1=yes).
Definition: jurassic.h:1509
int ig_o2
Emitter index of O2.
Definition: jurassic.h:1446
double rett_zmin
Minimum altitude for temperature retrieval [km].
Definition: jurassic.h:1536
double sfsza
Solar zenith angle at the surface [deg] (-999=auto).
Definition: jurassic.h:1476
int nd
Number of radiance channels.
Definition: jurassic.h:1449
int fov_n
Field-of-view number of data points.
Definition: jurassic.h:1527
int sftype
Surface treatment (0=none, 1=emissions, 2=downward, 3=solar).
Definition: jurassic.h:1473
int matrixfmt
Matrix data file format (1=ASCII, 2=binary, 3=netCDF).
Definition: jurassic.h:1491
int ncl
Number of cloud layer spectral grid points.
Definition: jurassic.h:1461
int ig_n2
Emitter index of N2.
Definition: jurassic.h:1443
int obsfmt
Observation data file format (1=ASCII, 2=binary, 3=netCDF).
Definition: jurassic.h:1488
int ctm_o2
Compute O2 continuum (0=no, 1=yes).
Definition: jurassic.h:1506
int ret_clk
Retrieve cloud layer extinction (0=no, 1=yes).
Definition: jurassic.h:1560
int nsf
Number of surface layer spectral grid points.
Definition: jurassic.h:1467
double rayds
Maximum step length for raytracing [km].
Definition: jurassic.h:1512
int ret_cldz
Retrieve cloud layer depth (0=no, 1=yes).
Definition: jurassic.h:1557
int ig_h2o
Emitter index of H2O.
Definition: jurassic.h:1440
int tblfmt
Look-up table file format (1=ASCII, 2=binary, 3=netCDF).
Definition: jurassic.h:1482
double raydz
Vertical step length for raytracing [km].
Definition: jurassic.h:1515
int write_bbt
Use brightness temperature instead of radiance (0=no, 1=yes).
Definition: jurassic.h:1569
double retp_zmax
Maximum altitude for pressure retrieval [km].
Definition: jurassic.h:1533
Line-of-sight data.
Definition: jurassic.h:1595
double sft
Surface temperature [K].
Definition: jurassic.h:1622
int np
Number of LOS points.
Definition: jurassic.h:1598
Observation geometry and radiance data.
Definition: jurassic.h:1657
int nr
Number of ray paths.
Definition: jurassic.h:1660
Retrieval control parameters.
Definition: jurassic.h:1709
double err_press_cz
Vertical correlation length for pressure error [km].
Definition: jurassic.h:1736
double err_press
Pressure error [%].
Definition: jurassic.h:1733
int err_ana
Carry out error analysis (0=no, 1=yes).
Definition: jurassic.h:1724
double err_temp_cz
Vertical correlation length for temperature error [km].
Definition: jurassic.h:1745
double conv_dmin
Minimum normalized step size in state space.
Definition: jurassic.h:1721
double err_temp
Temperature error [K].
Definition: jurassic.h:1742
double err_clz
Cloud height error [km].
Definition: jurassic.h:1769
double err_sft
Surface temperature error [K].
Definition: jurassic.h:1778
double err_temp_ch
Horizontal correlation length for temperature error [km].
Definition: jurassic.h:1748
int kernel_recomp
Re-computation of kernel matrix (number of iterations).
Definition: jurassic.h:1715
int conv_itmax
Maximum number of iterations.
Definition: jurassic.h:1718
double err_cldz
Cloud depth error [km].
Definition: jurassic.h:1772
double err_press_ch
Horizontal correlation length for pressure error [km].
Definition: jurassic.h:1739
int shared_io_profile
Profile index used for shared netCDF retrieval inputs and outputs.
Definition: jurassic.h:1784
Emissivity look-up tables.
Definition: jurassic.h:1842