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
425#define DEG2RAD(deg) ((deg) * (M_PI / 180.0))
426
443#define DIST(a, b) sqrt(DIST2(a, b))
444
460#define DIST2(a, b) \
461 ((a[0]-b[0])*(a[0]-b[0])+(a[1]-b[1])*(a[1]-b[1])+(a[2]-b[2])*(a[2]-b[2]))
462
478#define DOTP(a, b) (a[0]*b[0]+a[1]*b[1]+a[2]*b[2])
479
498#define FREAD(ptr, type, size, out) { \
499 if(fread(ptr, sizeof(type), size, out)!=size) \
500 ERRMSG("Error while reading!"); \
501 }
502
521#define FWRITE(ptr, type, size, out) { \
522 if(fwrite(ptr, sizeof(type), size, out)!=size) \
523 ERRMSG("Error while writing!"); \
524 }
525
544#define LIN(x0, y0, x1, y1, x) \
545 ((y0)+((y1)-(y0))/((x1)-(x0))*((x)-(x0)))
546
565#define LOGX(x0, y0, x1, y1, x) \
566 (((x)/(x0)>0 && (x1)/(x0)>0) \
567 ? ((y0)+((y1)-(y0))*log((x)/(x0))/log((x1)/(x0))) \
568 : LIN(x0, y0, x1, y1, x))
569
588#define LOGY(x0, y0, x1, y1, x) \
589 (((y1)/(y0)>0) \
590 ? ((y0)*exp(log((y1)/(y0))/((x1)-(x0))*((x)-(x0)))) \
591 : LIN(x0, y0, x1, y1, x))
592
609#define MAX(a,b) (((a)>(b))?(a):(b))
610
627#define MIN(a,b) (((a)<(b))?(a):(b))
628
640#define NC(cmd) { \
641 int nc_result=(cmd); \
642 if(nc_result!=NC_NOERR) \
643 ERRMSG("%s", nc_strerror(nc_result)); \
644 }
645
669#define NC_DEF_VAR(varname, type, ndims, dims, long_name, units, level, quant) { \
670 NC(nc_def_var(ncid, varname, type, ndims, dims, &varid)); \
671 NC(nc_put_att_text(ncid, varid, "long_name", strnlen(long_name, LEN), long_name)); \
672 NC(nc_put_att_text(ncid, varid, "units", strnlen(units, LEN), units)); \
673 if((quant) > 0) \
674 NC(nc_def_var_quantize(ncid, varid, NC_QUANTIZE_GRANULARBR, quant)); \
675 if((level) != 0) { \
676 NC(nc_def_var_deflate(ncid, varid, 1, 1, level)); \
677 /* unsigned int ulevel = (unsigned int)level; */ \
678 /* NC(nc_def_var_filter(ncid, varid, 32015, 1, (unsigned int[]){ulevel})); */ \
679 } \
680 }
681
699#define NC_GET_DOUBLE(varname, ptr, force) { \
700 if(force) { \
701 NC(nc_inq_varid(ncid, varname, &varid)); \
702 NC(nc_get_var_double(ncid, varid, ptr)); \
703 } else { \
704 if(nc_inq_varid(ncid, varname, &varid) == NC_NOERR) { \
705 NC(nc_get_var_double(ncid, varid, ptr)); \
706 } else \
707 WARN("netCDF variable %s is missing!", varname); \
708 } \
709 }
710
729#define NC_INQ_DIM(dimname, ptr, min, max, check) { \
730 int dimid; size_t naux; \
731 NC(nc_inq_dimid(ncid, dimname, &dimid)); \
732 NC(nc_inq_dimlen(ncid, dimid, &naux)); \
733 *ptr = (int)naux; \
734 if (check) \
735 if ((*ptr) < (min) || (*ptr) > (max)) \
736 ERRMSG("Dimension %s is out of range!", dimname); \
737 }
738
753#define NC_PUT_DOUBLE(varname, ptr, hyperslab) { \
754 NC(nc_inq_varid(ncid, varname, &varid)); \
755 if(hyperslab) { \
756 NC(nc_put_vara_double(ncid, varid, start, count, ptr)); \
757 } else { \
758 NC(nc_put_var_double(ncid, varid, ptr)); \
759 } \
760 }
761
777#define NC_PUT_FLOAT(varname, ptr, hyperslab) { \
778 NC(nc_inq_varid(ncid, varname, &varid)); \
779 if(hyperslab) { \
780 NC(nc_put_vara_float(ncid, varid, start, count, ptr)); \
781 } else { \
782 NC(nc_put_var_float(ncid, varid, ptr)); \
783 } \
784 }
785
800#define NC_PUT_INT(varname, ptr, hyperslab) { \
801 NC(nc_inq_varid(ncid, varname, &varid)); \
802 if(hyperslab) { \
803 NC(nc_put_vara_int(ncid, varid, start, count, ptr)); \
804 } else { \
805 NC(nc_put_var_int(ncid, varid, ptr)); \
806 } \
807 }
808
822#define NC_PUT_ATT(varname, attname, text) { \
823 NC(nc_inq_varid(ncid, varname, &varid)); \
824 NC(nc_put_att_text(ncid, varid, attname, strnlen(text, LEN), text)); \
825 }
826
839#define NC_PUT_ATT_GLOBAL(attname, text) \
840 NC(nc_put_att_text(ncid, NC_GLOBAL, attname, strnlen(text, LEN), text));
841
871#define NEDT(t_bg, nesr, nu) \
872 (BRIGHT(PLANCK((t_bg), (nu)) + (nesr), (nu)) - (t_bg))
873
900#define NESR(t_bg, nedt, nu) \
901 (PLANCK((t_bg) + (nedt), (nu)) - PLANCK((t_bg), (nu)))
902
917#define NORM(a) sqrt(DOTP(a, a))
918
948#define PLANCK(T, nu) \
949 (C1 * POW3(nu) / gsl_expm1(C2 * (nu) / (T)))
950
962#define POW2(x) ((x)*(x))
963
975#define POW3(x) ((x)*(x)*(x))
976
991#define RAD2DEG(rad) ((rad) * (180.0 / M_PI))
992
1007#define REFRAC(p, T) (7.753e-05 * (p) / (T))
1008
1023#define TBL_LOG(tbl, id, ig) \
1024 do { \
1025 for (int ip = 0; ip < (tbl)->np[(id)][(ig)]; ip++) \
1026 LOG(2, \
1027 "p[%2d]= %.5e hPa | T[0:%2d]= %.2f ... %.2f K | " \
1028 "u[0:%3d]= %.5e ... %.5e molec/cm^2 | " \
1029 "eps[0:%3d]= %.5e ... %.5e", \
1030 (ip), \
1031 (tbl)->p[(id)][(ig)][(ip)], \
1032 (tbl)->nt[(id)][(ig)][(ip)] - 1, \
1033 (tbl)->t[(id)][(ig)][(ip)][0], \
1034 (tbl)->t[(id)][(ig)][(ip)] \
1035 [(tbl)->nt[(id)][(ig)][(ip)] - 1], \
1036 (tbl)->nu[(id)][(ig)][(ip)][0] - 1, \
1037 exp((tbl)->logu[(id)][(ig)][(ip)][0][0]), \
1038 exp((tbl)->logu[(id)][(ig)][(ip)][0] \
1039 [(tbl)->nu[(id)][(ig)][(ip)][0] - 1]), \
1040 (tbl)->nu[(id)][(ig)][(ip)][0] - 1, \
1041 exp((tbl)->logeps[(id)][(ig)][(ip)][0][0]), \
1042 exp((tbl)->logeps[(id)][(ig)][(ip)][0] \
1043 [(tbl)->nu[(id)][(ig)][(ip)][0] - 1])); \
1044 } while (0)
1045
1059#define TIMER(name, mode) \
1060 {timer(name, __FILE__, __func__, __LINE__, mode);}
1061
1080#define TOK(line, tok, format, var) { \
1081 if(((tok)=strtok((line), " \t"))) { \
1082 if(sscanf(tok, format, &(var))!=1) continue; \
1083 } else ERRMSG("Error while reading!"); \
1084 }
1085
1086/* ------------------------------------------------------------
1087 Log messages...
1088 ------------------------------------------------------------ */
1089
1091#ifndef LOGLEV
1092#define LOGLEV 2
1093#endif
1094
1124#define LOG(level, ...) { \
1125 if(level >= 2) \
1126 printf(" "); \
1127 if(level <= LOGLEV) { \
1128 printf(__VA_ARGS__); \
1129 printf("\n"); \
1130 } \
1131 }
1132
1161#define WARN(...) { \
1162 printf("\nWarning (%s, %s, l%d): ", __FILE__, __func__, __LINE__); \
1163 LOG(0, __VA_ARGS__); \
1164 }
1165
1194#define ERRMSG(...) { \
1195 printf("\nError (%s, %s, l%d): ", __FILE__, __func__, __LINE__); \
1196 LOG(0, __VA_ARGS__); \
1197 exit(EXIT_FAILURE); \
1198 }
1199
1229#define PRINT(format, var) \
1230 printf("Print (%s, %s, l%d): %s= "format"\n", \
1231 __FILE__, __func__, __LINE__, #var, var);
1232
1233/* ------------------------------------------------------------
1234 Structs...
1235 ------------------------------------------------------------ */
1236
1244typedef struct {
1245
1247 int np;
1248
1250 double time[NP];
1251
1253 double z[NP];
1254
1256 double lon[NP];
1257
1259 double lat[NP];
1260
1262 double p[NP];
1263
1265 double t[NP];
1266
1268 double q[NG][NP];
1269
1271 double k[NW][NP];
1272
1274 double clz;
1275
1277 double cldz;
1278
1280 double clk[NCL];
1281
1283 double sft;
1284
1286 double sfeps[NSF];
1287
1288} atm_t;
1289
1297typedef struct {
1298
1300 int ng;
1301
1303 char emitter[NG][LEN];
1304
1307
1310
1313
1316
1318 int nd;
1319
1321 double nu[ND];
1322
1324 int nw;
1325
1327 int window[ND];
1328
1330 int ncl;
1331
1333 double clnu[NCL];
1334
1336 int nsf;
1337
1339 double sfnu[NSF];
1340
1343
1345 double sfsza;
1346
1348 char tblbase[LEN];
1349
1352
1355
1358
1360 double hydz;
1361
1364
1367
1370
1373
1376
1378 double rayds;
1379
1381 double raydz;
1382
1384 char fov[LEN];
1385
1387 double fov_dz[NSHAPE];
1388
1390 double fov_w[NSHAPE];
1391
1394
1397
1400
1403
1406
1408 double retq_zmin[NG];
1409
1411 double retq_zmax[NG];
1412
1414 double retk_zmin[NW];
1415
1417 double retk_zmax[NW];
1418
1421
1424
1427
1430
1433
1436
1439
1442
1444 char rfmbin[LEN];
1445
1447 char rfmhit[LEN];
1448
1450 char rfmxsc[NG][LEN];
1451
1452} ctl_t;
1453
1461typedef struct {
1462
1464 int np;
1465
1467 double z[NLOS];
1468
1470 double lon[NLOS];
1471
1473 double lat[NLOS];
1474
1476 double p[NLOS];
1477
1479 double t[NLOS];
1480
1482 double q[NLOS][NG];
1483
1485 double k[NLOS][ND];
1486
1488 double sft;
1489
1491 double sfeps[ND];
1492
1494 double ds[NLOS];
1495
1497 double u[NLOS][NG];
1498
1500 double cgp[NLOS][NG];
1501
1503 double cgt[NLOS][NG];
1504
1506 double cgu[NLOS][NG];
1507
1509 double eps[NLOS][ND];
1510
1512 double src[NLOS][ND];
1513
1514} los_t;
1515
1523typedef struct {
1524
1526 int nr;
1527
1529 double time[NR];
1530
1532 double obsz[NR];
1533
1535 double obslon[NR];
1536
1538 double obslat[NR];
1539
1541 double vpz[NR];
1542
1544 double vplon[NR];
1545
1547 double vplat[NR];
1548
1550 double tpz[NR];
1551
1553 double tplon[NR];
1554
1556 double tplat[NR];
1557
1559 double tau[ND][NR];
1560
1562 double rad[ND][NR];
1563
1564} obs_t;
1565
1575typedef struct {
1576
1578 char dir[LEN];
1579
1582
1585
1588
1591
1593 double err_formod[ND];
1594
1596 double err_noise[ND];
1597
1600
1603
1606
1608 double err_temp;
1609
1612
1615
1617 double err_q[NG];
1618
1620 double err_q_cz[NG];
1621
1623 double err_q_ch[NG];
1624
1626 double err_k[NW];
1627
1629 double err_k_cz[NW];
1630
1632 double err_k_ch[NW];
1633
1635 double err_clz;
1636
1638 double err_cldz;
1639
1641 double err_clk[NCL];
1642
1644 double err_sft;
1645
1647 double err_sfeps[NSF];
1648
1649} ret_t;
1650
1657typedef struct {
1658
1660 int np[ND][NG];
1661
1663 int nt[ND][NG][TBLNP];
1664
1666 int nu[ND][NG][TBLNP][TBLNT];
1667
1669 double p[ND][NG][TBLNP];
1670
1672 double t[ND][NG][TBLNP][TBLNT];
1673
1675 float *logu[ND][NG][TBLNP][TBLNT];
1676
1678 float *logeps[ND][NG][TBLNP][TBLNT];
1679
1681 int filt_n[ND];
1682
1684 double filt_nu[ND][NSHAPE];
1685
1687 double filt_f[ND][NSHAPE];
1688
1690 double st[TBLNS];
1691
1693 double sr[TBLNS][ND];
1694
1695} tbl_t;
1696
1697/* ------------------------------------------------------------
1698 Functions...
1699 ------------------------------------------------------------ */
1700
1752void analyze_avk(
1753 const ret_t * ret,
1754 const ctl_t * ctl,
1755 const atm_t * atm,
1756 const int *iqa,
1757 const int *ipa,
1758 const gsl_matrix * avk);
1759
1816 const gsl_matrix * avk,
1817 const int iq,
1818 const int *ipa,
1819 const size_t *n0,
1820 const size_t *n1,
1821 double *cont,
1822 double *res);
1823
1850size_t atm2x(
1851 const ctl_t * ctl,
1852 const atm_t * atm,
1853 gsl_vector * x,
1854 int *iqa,
1855 int *ipa);
1856
1878void atm2x_help(
1879 const double value,
1880 const int value_iqa,
1881 const int value_ip,
1882 gsl_vector * x,
1883 int *iqa,
1884 int *ipa,
1885 size_t *n);
1886
1901void cart2geo(
1902 const double *x,
1903 double *z,
1904 double *lon,
1905 double *lat);
1906
1922void climatology(
1923 const ctl_t * ctl,
1924 atm_t * atm);
1925
1950double cos_sza(
1951 const double sec,
1952 const double lon,
1953 const double lat);
1954
1994double cost_function(
1995 const gsl_vector * dx,
1996 const gsl_vector * dy,
1997 const gsl_matrix * s_a_inv,
1998 const gsl_vector * sig_eps_inv);
1999
2005double ctmco2(
2006 const double nu,
2007 const double p,
2008 const double t,
2009 const double u);
2010
2016double ctmh2o(
2017 const double nu,
2018 const double p,
2019 const double t,
2020 const double q,
2021 const double u);
2022
2052double ctmn2(
2053 const double nu,
2054 const double p,
2055 const double t);
2056
2086double ctmo2(
2087 const double nu,
2088 const double p,
2089 const double t);
2090
2115void copy_atm(
2116 const ctl_t * ctl,
2117 atm_t * atm_dest,
2118 const atm_t * atm_src,
2119 const int init);
2120
2145void copy_obs(
2146 const ctl_t * ctl,
2147 obs_t * obs_dest,
2148 const obs_t * obs_src,
2149 const int init);
2150
2171void day2doy(
2172 int year,
2173 int mon,
2174 int day,
2175 int *doy);
2176
2199void doy2day(
2200 int year,
2201 int doy,
2202 int *mon,
2203 int *day);
2204
2224int find_emitter(
2225 const ctl_t * ctl,
2226 const char *emitter);
2227
2253void formod(
2254 const ctl_t * ctl,
2255 const tbl_t * tbl,
2256 atm_t * atm,
2257 obs_t * obs);
2258
2284void formod_continua(
2285 const ctl_t * ctl,
2286 const los_t * los,
2287 const int ip,
2288 double *beta);
2289
2311void formod_fov(
2312 const ctl_t * ctl,
2313 obs_t * obs);
2314
2344void formod_pencil(
2345 const ctl_t * ctl,
2346 const tbl_t * tbl,
2347 const atm_t * atm,
2348 obs_t * obs,
2349 const int ir);
2350
2389void formod_rfm(
2390 const ctl_t * ctl,
2391 const atm_t * atm,
2392 obs_t * obs);
2393
2418void formod_srcfunc(
2419 const ctl_t * ctl,
2420 const tbl_t * tbl,
2421 const double t,
2422 double *src);
2423
2450void geo2cart(
2451 const double z,
2452 const double lon,
2453 const double lat,
2454 double *x);
2455
2483void hydrostatic(
2484 const ctl_t * ctl,
2485 atm_t * atm);
2486
2510void idx2name(
2511 const ctl_t * ctl,
2512 const int idx,
2513 char *quantity);
2514
2544void init_srcfunc(
2545 const ctl_t * ctl,
2546 tbl_t * tbl);
2547
2572void intpol_atm(
2573 const ctl_t * ctl,
2574 const atm_t * atm,
2575 const double z,
2576 double *p,
2577 double *t,
2578 double *q,
2579 double *k);
2580
2614void intpol_tbl_cga(
2615 const ctl_t * ctl,
2616 const tbl_t * tbl,
2617 const los_t * los,
2618 const int ip,
2619 double tau_path[ND][NG],
2620 double tau_seg[ND]);
2621
2660void intpol_tbl_ega(
2661 const ctl_t * ctl,
2662 const tbl_t * tbl,
2663 const los_t * los,
2664 const int ip,
2665 double tau_path[ND][NG],
2666 double tau_seg[ND]);
2667
2694double intpol_tbl_eps(
2695 const tbl_t * tbl,
2696 const int ig,
2697 const int id,
2698 const int ip,
2699 const int it,
2700 const double logu);
2701
2732double intpol_tbl_u(
2733 const tbl_t * tbl,
2734 const int ig,
2735 const int id,
2736 const int ip,
2737 const int it,
2738 const double logeps);
2739
2766void jsec2time(
2767 const double jsec,
2768 int *year,
2769 int *mon,
2770 int *day,
2771 int *hour,
2772 int *min,
2773 int *sec,
2774 double *remain);
2775
2817void kernel(
2818 const ctl_t * ctl,
2819 const tbl_t * tbl,
2820 atm_t * atm,
2821 obs_t * obs,
2822 gsl_matrix * k);
2823
2850int locate_irr(
2851 const double *xx,
2852 const int n,
2853 const double x);
2854
2881int locate_reg(
2882 const double *xx,
2883 const int n,
2884 const double x);
2885
2911int locate_tbl(
2912 const float *xx,
2913 const int n,
2914 const double x);
2915
2954void matrix_invert(
2955 gsl_matrix * a);
2956
3001void matrix_product(
3002 const gsl_matrix * a,
3003 const gsl_vector * b,
3004 const int transpose,
3005 gsl_matrix * c);
3006
3034size_t obs2y(
3035 const ctl_t * ctl,
3036 const obs_t * obs,
3037 gsl_vector * y,
3038 int *ida,
3039 int *ira);
3040
3092 ret_t * ret,
3093 ctl_t * ctl,
3094 tbl_t * tbl,
3095 obs_t * obs_meas,
3096 obs_t * obs_i,
3097 atm_t * atm_apr,
3098 atm_t * atm_i,
3099 double *chisq);
3100
3139void raytrace(
3140 const ctl_t * ctl,
3141 const atm_t * atm,
3142 obs_t * obs,
3143 los_t * los,
3144 const int ir);
3145
3191void read_atm(
3192 const char *dirname,
3193 const char *filename,
3194 const ctl_t * ctl,
3195 atm_t * atm);
3196
3243void read_atm_asc(
3244 const char *filename,
3245 const ctl_t * ctl,
3246 atm_t * atm);
3247
3297void read_atm_bin(
3298 const char *filename,
3299 const ctl_t * ctl,
3300 atm_t * atm);
3301
3341void read_atm_nc(
3342 const char *filename,
3343 const ctl_t * ctl,
3344 atm_t * atm,
3345 int profile);
3346
3387void read_ctl(
3388 int argc,
3389 char *argv[],
3390 ctl_t * ctl);
3391
3423void read_matrix(
3424 const char *dirname,
3425 const char *filename,
3426 gsl_matrix * matrix);
3427
3462void read_obs(
3463 const char *dirname,
3464 const char *filename,
3465 const ctl_t * ctl,
3466 obs_t * obs);
3467
3504void read_obs_asc(
3505 const char *filename,
3506 const ctl_t * ctl,
3507 obs_t * obs);
3508
3549void read_obs_bin(
3550 const char *filename,
3551 const ctl_t * ctl,
3552 obs_t * obs);
3553
3585void read_obs_nc(
3586 const char *filename,
3587 const ctl_t * ctl,
3588 obs_t * obs,
3589 const int profile);
3590
3630double read_obs_rfm(
3631 const char *basename,
3632 const double z,
3633 const double *nu,
3634 const double *f,
3635 const int n);
3636
3701void read_ret(
3702 int argc,
3703 char *argv[],
3704 const ctl_t * ctl,
3705 ret_t * ret);
3706
3747void read_rfm_spec(
3748 const char *filename,
3749 double *nu,
3750 double *rad,
3751 int *npts);
3752
3786void read_shape(
3787 const char *filename,
3788 double *x,
3789 double *y,
3790 int *n);
3791
3828 const ctl_t * ctl);
3829
3865void read_tbl_asc(
3866 const ctl_t * ctl,
3867 tbl_t * tbl,
3868 const int id,
3869 const int ig);
3870
3902void read_tbl_bin(
3903 const ctl_t * ctl,
3904 tbl_t * tbl,
3905 const int id,
3906 const int ig);
3907
3951 const ctl_t * ctl,
3952 tbl_t * tbl,
3953 int id,
3954 int ig,
3955 int ncid);
3956
4011double scan_ctl(
4012 int argc,
4013 char *argv[],
4014 const char *varname,
4015 const int arridx,
4016 const char *defvalue,
4017 char *value);
4018
4072void set_cov_apr(
4073 const ret_t * ret,
4074 const ctl_t * ctl,
4075 const atm_t * atm,
4076 const int *iqa,
4077 const int *ipa,
4078 gsl_matrix * s_a);
4079
4133void set_cov_meas(
4134 const ret_t * ret,
4135 const ctl_t * ctl,
4136 const obs_t * obs,
4137 gsl_vector * sig_noise,
4138 gsl_vector * sig_formod,
4139 gsl_vector * sig_eps_inv);
4140
4176double sza(
4177 double sec,
4178 double lon,
4179 double lat);
4180
4220void tangent_point(
4221 const los_t * los,
4222 double *tpz,
4223 double *tplon,
4224 double *tplat);
4225
4247void tbl_free(
4248 const ctl_t * ctl,
4249 tbl_t * tbl);
4250
4290void tbl_pack(
4291 const tbl_t * tbl,
4292 int id,
4293 int ig,
4294 uint8_t * buf,
4295 size_t *bytes_used);
4296
4307size_t tbl_packed_size(
4308 const tbl_t * tbl,
4309 int id,
4310 int ig);
4311
4350size_t tbl_unpack(
4351 tbl_t * tbl,
4352 int id,
4353 int ig,
4354 const uint8_t * buf);
4355
4356
4381void time2jsec(
4382 const int year,
4383 const int mon,
4384 const int day,
4385 const int hour,
4386 const int min,
4387 const int sec,
4388 const double remain,
4389 double *jsec);
4390
4432void timer(
4433 const char *name,
4434 const char *file,
4435 const char *func,
4436 int line,
4437 int mode);
4438
4483void write_atm(
4484 const char *dirname,
4485 const char *filename,
4486 const ctl_t * ctl,
4487 const atm_t * atm);
4488
4540void write_atm_asc(
4541 const char *filename,
4542 const ctl_t * ctl,
4543 const atm_t * atm);
4544
4598void write_atm_bin(
4599 const char *filename,
4600 const ctl_t * ctl,
4601 const atm_t * atm);
4602
4646void write_atm_nc(
4647 const char *filename,
4648 const ctl_t * ctl,
4649 const atm_t * atm,
4650 int profile);
4651
4703void write_atm_rfm(
4704 const char *filename,
4705 const ctl_t * ctl,
4706 const atm_t * atm);
4707
4755void write_matrix(
4756 const char *dirname,
4757 const char *filename,
4758 const ctl_t * ctl,
4759 const gsl_matrix * matrix,
4760 const atm_t * atm,
4761 const obs_t * obs,
4762 const char *rowspace,
4763 const char *colspace,
4764 const char *sort);
4765
4798void write_obs(
4799 const char *dirname,
4800 const char *filename,
4801 const ctl_t * ctl,
4802 const obs_t * obs);
4803
4842void write_obs_asc(
4843 const char *filename,
4844 const ctl_t * ctl,
4845 const obs_t * obs);
4846
4883void write_obs_bin(
4884 const char *filename,
4885 const ctl_t * ctl,
4886 const obs_t * obs);
4887
4923void write_obs_nc(
4924 const char *filename,
4925 const ctl_t * ctl,
4926 const obs_t * obs,
4927 const int profile);
4928
4963void write_shape(
4964 const char *filename,
4965 const double *x,
4966 const double *y,
4967 const int n);
4968
5014void write_stddev(
5015 const char *quantity,
5016 const ret_t * ret,
5017 const ctl_t * ctl,
5018 const atm_t * atm,
5019 const gsl_matrix * s);
5020
5056void write_tbl(
5057 const ctl_t * ctl,
5058 const tbl_t * tbl);
5059
5085void write_tbl_asc(
5086 const ctl_t * ctl,
5087 const tbl_t * tbl,
5088 const int id,
5089 const int ig);
5090
5118void write_tbl_bin(
5119 const ctl_t * ctl,
5120 const tbl_t * tbl,
5121 const int id,
5122 const int ig);
5123
5156void write_tbl_nc(
5157 const ctl_t * ctl,
5158 const tbl_t * tbl,
5159 const int id,
5160 const int ig);
5161
5208void x2atm(
5209 const ctl_t * ctl,
5210 const gsl_vector * x,
5211 atm_t * atm);
5212
5241void x2atm_help(
5242 double *value,
5243 const gsl_vector * x,
5244 size_t *n);
5245
5278void y2obs(
5279 const ctl_t * ctl,
5280 const gsl_vector * y,
5281 obs_t * obs);
5282
5283#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:4225
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:6299
#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:5892
void read_matrix(const char *dirname, const char *filename, gsl_matrix *matrix)
Read a numerical matrix from an ASCII file.
Definition: jurassic.c:5574
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:6882
void read_rfm_spec(const char *filename, double *nu, double *rad, int *npts)
Read a Reference Forward Model (RFM) ASCII spectrum.
Definition: jurassic.c:6003
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:6687
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:5252
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:7614
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:6920
int locate_reg(const double *xx, const int n, const double x)
Locate index for interpolation on a regular (uniform) grid.
Definition: jurassic.c:4474
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:4128
void formod_rfm(const ctl_t *ctl, const atm_t *atm, obs_t *obs)
Interface routine for the Reference Forward Model (RFM).
Definition: jurassic.c:3693
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:7738
void tbl_free(const ctl_t *ctl, tbl_t *tbl)
Free lookup table and all internally allocated memory.
Definition: jurassic.c:6648
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:8019
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:8103
void idx2name(const ctl_t *ctl, const int idx, char *quantity)
Convert a quantity index to a descriptive name string.
Definition: jurassic.c:3931
void read_ctl(int argc, char *argv[], ctl_t *ctl)
Read model control parameters from command-line and configuration input.
Definition: jurassic.c:5450
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:5736
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:4889
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:4545
int locate_irr(const double *xx, const int n, const double x)
Locate index for interpolation on an irregular grid.
Definition: jurassic.c:4444
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:7142
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:6851
void x2atm(const ctl_t *ctl, const gsl_vector *x, atm_t *atm)
Map retrieval state vector back to atmospheric structure.
Definition: jurassic.c:8168
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:8061
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:4042
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:7325
#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:6560
#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:7976
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:6775
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:4017
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:4617
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:4270
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:6604
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:6991
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:5126
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:5614
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:7541
#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:7670
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:8218
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:4349
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:5814
void init_srcfunc(const ctl_t *ctl, tbl_t *tbl)
Initialize source function lookup tables from emissivity data.
Definition: jurassic.c:3970
void matrix_invert(gsl_matrix *a)
Invert a square matrix, optimized for diagonal or symmetric positive-definite matrices.
Definition: jurassic.c:4515
void read_ret(int argc, char *argv[], const ctl_t *ctl, ret_t *ret)
Read retrieval configuration and error parameters.
Definition: jurassic.c:5950
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:6380
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:7913
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:5200
#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:6176
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:3834
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:4316
#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:6106
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:6343
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:5346
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:4590
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:8230
void hydrostatic(const ctl_t *ctl, atm_t *atm)
Adjust pressure profile using the hydrostatic equation.
Definition: jurassic.c:3871
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:6064
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:3851
#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:7363
#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:6449
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:7943
#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:6741
#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:4493
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:7058
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:5691
#define NW
Maximum number of spectral windows.
Definition: jurassic.h:263
Atmospheric profile data.
Definition: jurassic.h:1244
double clz
Cloud layer height [km].
Definition: jurassic.h:1274
int np
Number of data points.
Definition: jurassic.h:1247
double cldz
Cloud layer depth [km].
Definition: jurassic.h:1277
double sft
Surface temperature [K].
Definition: jurassic.h:1283
Control parameters.
Definition: jurassic.h:1297
int write_matrix
Write matrix file (0=no, 1=yes).
Definition: jurassic.h:1438
int nw
Number of spectral windows.
Definition: jurassic.h:1324
int atmfmt
Atmospheric data file format (1=ASCII, 2=binary, 3=netCDF).
Definition: jurassic.h:1354
double retp_zmin
Minimum altitude for pressure retrieval [km].
Definition: jurassic.h:1396
int ig_co2
Emitter index of CO2.
Definition: jurassic.h:1306
double hydz
Reference height for hydrostatic pressure profile (-999 to skip) [km].
Definition: jurassic.h:1360
int ctm_co2
Compute CO2 continuum (0=no, 1=yes).
Definition: jurassic.h:1363
double rett_zmax
Maximum altitude for temperature retrieval [km].
Definition: jurassic.h:1405
int ret_sfeps
Retrieve surface layer emissivity (0=no, 1=yes).
Definition: jurassic.h:1432
int ret_sft
Retrieve surface layer temperature (0=no, 1=yes).
Definition: jurassic.h:1429
int ret_clz
Retrieve cloud layer height (0=no, 1=yes).
Definition: jurassic.h:1420
int ctm_n2
Compute N2 continuum (0=no, 1=yes).
Definition: jurassic.h:1369
int ctm_h2o
Compute H2O continuum (0=no, 1=yes).
Definition: jurassic.h:1366
int formod
Forward model (0=CGA, 1=EGA, 2=RFM).
Definition: jurassic.h:1441
int ng
Number of emitters.
Definition: jurassic.h:1300
int refrac
Take into account refractivity (0=no, 1=yes).
Definition: jurassic.h:1375
int ig_o2
Emitter index of O2.
Definition: jurassic.h:1315
double rett_zmin
Minimum altitude for temperature retrieval [km].
Definition: jurassic.h:1402
double sfsza
Solar zenith angle at the surface [deg] (-999=auto).
Definition: jurassic.h:1345
int nd
Number of radiance channels.
Definition: jurassic.h:1318
int fov_n
Field-of-view number of data points.
Definition: jurassic.h:1393
int sftype
Surface treatment (0=none, 1=emissions, 2=downward, 3=solar).
Definition: jurassic.h:1342
int ncl
Number of cloud layer spectral grid points.
Definition: jurassic.h:1330
int ig_n2
Emitter index of N2.
Definition: jurassic.h:1312
int obsfmt
Observation data file format (1=ASCII, 2=binary, 3=netCDF).
Definition: jurassic.h:1357
int ctm_o2
Compute O2 continuum (0=no, 1=yes).
Definition: jurassic.h:1372
int ret_clk
Retrieve cloud layer extinction (0=no, 1=yes).
Definition: jurassic.h:1426
int nsf
Number of surface layer spectral grid points.
Definition: jurassic.h:1336
double rayds
Maximum step length for raytracing [km].
Definition: jurassic.h:1378
int ret_cldz
Retrieve cloud layer depth (0=no, 1=yes).
Definition: jurassic.h:1423
int ig_h2o
Emitter index of H2O.
Definition: jurassic.h:1309
int tblfmt
Look-up table file format (1=ASCII, 2=binary, 3=netCDF).
Definition: jurassic.h:1351
double raydz
Vertical step length for raytracing [km].
Definition: jurassic.h:1381
int write_bbt
Use brightness temperature instead of radiance (0=no, 1=yes).
Definition: jurassic.h:1435
double retp_zmax
Maximum altitude for pressure retrieval [km].
Definition: jurassic.h:1399
Line-of-sight data.
Definition: jurassic.h:1461
double sft
Surface temperature [K].
Definition: jurassic.h:1488
int np
Number of LOS points.
Definition: jurassic.h:1464
Observation geometry and radiance data.
Definition: jurassic.h:1523
int nr
Number of ray paths.
Definition: jurassic.h:1526
Retrieval control parameters.
Definition: jurassic.h:1575
double err_press_cz
Vertical correlation length for pressure error [km].
Definition: jurassic.h:1602
double err_press
Pressure error [%].
Definition: jurassic.h:1599
int err_ana
Carry out error analysis (0=no, 1=yes).
Definition: jurassic.h:1590
double err_temp_cz
Vertical correlation length for temperature error [km].
Definition: jurassic.h:1611
double conv_dmin
Minimum normalized step size in state space.
Definition: jurassic.h:1587
double err_temp
Temperature error [K].
Definition: jurassic.h:1608
double err_clz
Cloud height error [km].
Definition: jurassic.h:1635
double err_sft
Surface temperature error [K].
Definition: jurassic.h:1644
double err_temp_ch
Horizontal correlation length for temperature error [km].
Definition: jurassic.h:1614
int kernel_recomp
Re-computation of kernel matrix (number of iterations).
Definition: jurassic.h:1581
int conv_itmax
Maximum number of iterations.
Definition: jurassic.h:1584
double err_cldz
Cloud depth error [km].
Definition: jurassic.h:1638
double err_press_ch
Horizontal correlation length for pressure error [km].
Definition: jurassic.h:1605
Emissivity look-up tables.
Definition: jurassic.h:1657