JURASSIC
jurassic.h
Go to the documentation of this file.
1/*
2 This file is part of JURASSIC.
3
4 JURASSIC is free software: you can redistribute it and/or modify
5 it under the terms of the GNU General Public License as published by
6 the Free Software Foundation, either version 3 of the License, or
7 (at your option) any later version.
8
9 JURASSIC is distributed in the hope that it will be useful,
10 but WITHOUT ANY WARRANTY; without even the implied warranty of
11 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12 GNU General Public License for more details.
13
14 You should have received a copy of the GNU General Public License
15 along with JURASSIC. If not, see <http://www.gnu.org/licenses/>.
16
17 Copyright (C) 2003-2025 Forschungszentrum Juelich GmbH
18*/
19
100#ifndef JURASSIC_H
101#define JURASSIC_H
102
103#include <errno.h>
104#include <gsl/gsl_math.h>
105#include <gsl/gsl_blas.h>
106#include <gsl/gsl_linalg.h>
107#include <gsl/gsl_randist.h>
108#include <gsl/gsl_rng.h>
109#include <gsl/gsl_statistics.h>
110#include <math.h>
111#include <omp.h>
112#include <stdio.h>
113#include <stdlib.h>
114#include <string.h>
115#include <time.h>
116
117/* ------------------------------------------------------------
118 Macros...
119 ------------------------------------------------------------ */
120
122#define ALLOC(ptr, type, n) \
123 if((ptr=malloc((size_t)(n)*sizeof(type)))==NULL) \
124 ERRMSG("Out of memory!");
125
127#define BRIGHT(rad, nu) \
128 (C2 * (nu) / gsl_log1p(C1 * POW3(nu) / (rad)))
129
131#define DEG2RAD(deg) \
132 ((deg) * (M_PI / 180.0))
133
135#define DIST(a, b) sqrt(DIST2(a, b))
136
138#define DIST2(a, b) \
139 ((a[0]-b[0])*(a[0]-b[0])+(a[1]-b[1])*(a[1]-b[1])+(a[2]-b[2])*(a[2]-b[2]))
140
142#define DOTP(a, b) (a[0]*b[0]+a[1]*b[1]+a[2]*b[2])
143
145#define FREAD(ptr, type, size, out) { \
146 if(fread(ptr, sizeof(type), size, out)!=size) \
147 ERRMSG("Error while reading!"); \
148 }
149
151#define FWRITE(ptr, type, size, out) { \
152 if(fwrite(ptr, sizeof(type), size, out)!=size) \
153 ERRMSG("Error while writing!"); \
154 }
155
157#define MAX(a,b) \
158 (((a)>(b))?(a):(b))
159
161#define MIN(a,b) \
162 (((a)<(b))?(a):(b))
163
165#define LIN(x0, y0, x1, y1, x) \
166 ((y0)+((y1)-(y0))/((x1)-(x0))*((x)-(x0)))
167
169#define LOGX(x0, y0, x1, y1, x) \
170 (((x)/(x0)>0 && (x1)/(x0)>0) \
171 ? ((y0)+((y1)-(y0))*log((x)/(x0))/log((x1)/(x0))) \
172 : LIN(x0, y0, x1, y1, x))
173
175#define LOGY(x0, y0, x1, y1, x) \
176 (((y1)/(y0)>0) \
177 ? ((y0)*exp(log((y1)/(y0))/((x1)-(x0))*((x)-(x0)))) \
178 : LIN(x0, y0, x1, y1, x))
179
181#define NORM(a) sqrt(DOTP(a, a))
182
184#define PLANCK(T, nu) \
185 (C1 * POW3(nu) / gsl_expm1(C2 * (nu) / (T)))
186
188#define POW2(x) ((x)*(x))
189
191#define POW3(x) ((x)*(x)*(x))
192
194#define RAD2DEG(rad) \
195 ((rad) * (180.0 / M_PI))
196
198#define REFRAC(p, T) \
199 (7.753e-05 * (p) / (T))
200
202#define TIMER(name, mode) \
203 {timer(name, __FILE__, __func__, __LINE__, mode);}
204
206#define TOK(line, tok, format, var) { \
207 if(((tok)=strtok((line), " \t"))) { \
208 if(sscanf(tok, format, &(var))!=1) continue; \
209 } else ERRMSG("Error while reading!"); \
210 }
211
212/* ------------------------------------------------------------
213 Log messages...
214 ------------------------------------------------------------ */
215
217#ifndef LOGLEV
218#define LOGLEV 2
219#endif
220
222#define LOG(level, ...) { \
223 if(level >= 2) \
224 printf(" "); \
225 if(level <= LOGLEV) { \
226 printf(__VA_ARGS__); \
227 printf("\n"); \
228 } \
229 }
230
232#define WARN(...) { \
233 printf("\nWarning (%s, %s, l%d): ", __FILE__, __func__, __LINE__); \
234 LOG(0, __VA_ARGS__); \
235 }
236
238#define ERRMSG(...) { \
239 printf("\nError (%s, %s, l%d): ", __FILE__, __func__, __LINE__); \
240 LOG(0, __VA_ARGS__); \
241 exit(EXIT_FAILURE); \
242 }
243
245#define PRINT(format, var) \
246 printf("Print (%s, %s, l%d): %s= "format"\n", \
247 __FILE__, __func__, __LINE__, #var, var);
248
249/* ------------------------------------------------------------
250 Constants...
251 ------------------------------------------------------------ */
252
254#ifndef C1
255#define C1 1.19104259e-8
256#endif
257
259#ifndef C2
260#define C2 1.43877506
261#endif
262
264#ifndef EPSMIN
265#define EPSMIN 0
266#endif
267
269#ifndef EPSMAX
270#define EPSMAX 1
271#endif
272
274#ifndef G0
275#define G0 9.80665
276#endif
277
279#ifndef H0
280#define H0 7.0
281#endif
282
284#ifndef KB
285#define KB 1.3806504e-23
286#endif
287
289#ifndef ME
290#define ME 5.976e24
291#endif
292
294#ifndef NA
295#define NA 6.02214199e23
296#endif
297
299#ifndef N2
300#define N2 0.78084
301#endif
302
304#ifndef O2
305#define O2 0.20946
306#endif
307
309#ifndef P0
310#define P0 1013.25
311#endif
312
314#ifndef RE
315#define RE 6367.421
316#endif
317
319#ifndef RI
320#define RI 8.3144598
321#endif
322
324#ifndef T0
325#define T0 273.15
326#endif
327
329#ifndef TMIN
330#define TMIN 100.
331#endif
332
334#ifndef TMAX
335#define TMAX 400.
336#endif
337
339#ifndef TSUN
340#define TSUN 5780.
341#endif
342
344#ifndef UMIN
345#define UMIN 0
346#endif
347
349#ifndef UMAX
350#define UMAX 1e30
351#endif
352
353/* ------------------------------------------------------------
354 Dimensions...
355 ------------------------------------------------------------ */
356
358#ifndef NCL
359#define NCL 8
360#endif
361
363#ifndef ND
364#define ND 128
365#endif
366
368#ifndef NG
369#define NG 8
370#endif
371
373#ifndef NP
374#define NP 256
375#endif
376
378#ifndef NR
379#define NR 256
380#endif
381
383#ifndef NSF
384#define NSF 8
385#endif
386
388#ifndef NW
389#define NW 4
390#endif
391
393#ifndef LEN
394#define LEN 10000
395#endif
396
398#ifndef M
399#define M (NR*ND)
400#endif
401
403#ifndef N
404#define N ((2+NG+NW)*NP+NCL+NSF+5)
405#endif
406
408#ifndef NQ
409#define NQ (7+NG+NW+NCL+NSF)
410#endif
411
413#ifndef NLOS
414#define NLOS 4096
415#endif
416
418#ifndef NSHAPE
419#define NSHAPE 20000
420#endif
421
423#ifndef NFOV
424#define NFOV 5
425#endif
426
428#ifndef TBLNP
429#define TBLNP 41
430#endif
431
433#ifndef TBLNT
434#define TBLNT 30
435#endif
436
438#ifndef TBLNU
439#define TBLNU 320
440#endif
441
443#ifndef TBLNS
444#define TBLNS 1200
445#endif
446
448#ifndef RFMNPTS
449#define RFMNPTS 10000000
450#endif
451
452/* ------------------------------------------------------------
453 Quantity indices...
454 ------------------------------------------------------------ */
455
457#define IDXP 0
458
460#define IDXT 1
461
463#define IDXQ(ig) (2+ig)
464
466#define IDXK(iw) (2+ctl->ng+iw)
467
469#define IDXCLZ (2+ctl->ng+ctl->nw)
470
472#define IDXCLDZ (3+ctl->ng+ctl->nw)
473
475#define IDXCLK(icl) (4+ctl->ng+ctl->nw+icl)
476
478#define IDXSFZ (4+ctl->ng+ctl->nw+ctl->ncl)
479
481#define IDXSFP (5+ctl->ng+ctl->nw+ctl->ncl)
482
484#define IDXSFT (6+ctl->ng+ctl->nw+ctl->ncl)
485
487#define IDXSFEPS(isf) (7+ctl->ng+ctl->nw+ctl->ncl+isf)
488
489/* ------------------------------------------------------------
490 Structs...
491 ------------------------------------------------------------ */
492
494typedef struct {
495
497 int np;
498
500 double time[NP];
501
503 double z[NP];
504
506 double lon[NP];
507
509 double lat[NP];
510
512 double p[NP];
513
515 double t[NP];
516
518 double q[NG][NP];
519
521 double k[NW][NP];
522
524 double clz;
525
527 double cldz;
528
530 double clk[NCL];
531
533 double sfz;
534
536 double sfp;
537
539 double sft;
540
542 double sfeps[NSF];
543
544} atm_t;
545
547typedef struct {
548
550 int ng;
551
553 char emitter[NG][LEN];
554
557
560
562 int ig_n2;
563
565 int ig_o2;
566
568 int nd;
569
571 double nu[ND];
572
574 int nw;
575
577 int window[ND];
578
580 int ncl;
581
583 double clnu[NCL];
584
586 int nsf;
587
589 double sfnu[NSF];
590
593
595 double sfsza;
596
598 char tblbase[LEN];
599
602
604 double hydz;
605
608
611
614
617
620
622 double rayds;
623
625 double raydz;
626
628 char fov[LEN];
629
631 double fov_dz[NSHAPE];
632
634 double fov_w[NSHAPE];
635
637 int fov_n;
638
640 double retp_zmin;
641
643 double retp_zmax;
644
646 double rett_zmin;
647
649 double rett_zmax;
650
652 double retq_zmin[NG];
653
655 double retq_zmax[NG];
656
658 double retk_zmin[NW];
659
661 double retk_zmax[NW];
662
665
668
671
674
677
680
683
686
689
692
694 char rfmbin[LEN];
695
697 char rfmhit[LEN];
698
700 char rfmxsc[NG][LEN];
701
702} ctl_t;
703
705typedef struct {
706
708 int np;
709
711 double z[NLOS];
712
714 double lon[NLOS];
715
717 double lat[NLOS];
718
720 double p[NLOS];
721
723 double t[NLOS];
724
726 double q[NLOS][NG];
727
729 double k[NLOS][ND];
730
732 double sft;
733
735 double sfeps[ND];
736
738 double ds[NLOS];
739
741 double u[NLOS][NG];
742
744 double cgp[NLOS][NG];
745
747 double cgt[NLOS][NG];
748
750 double cgu[NLOS][NG];
751
753 double eps[NLOS][ND];
754
756 double src[NLOS][ND];
757
758} los_t;
759
761typedef struct {
762
764 int nr;
765
767 double time[NR];
768
770 double obsz[NR];
771
773 double obslon[NR];
774
776 double obslat[NR];
777
779 double vpz[NR];
780
782 double vplon[NR];
783
785 double vplat[NR];
786
788 double tpz[NR];
789
791 double tplon[NR];
792
794 double tplat[NR];
795
797 double tau[ND][NR];
798
800 double rad[ND][NR];
801
802} obs_t;
803
805typedef struct {
806
808 int np[ND][NG];
809
811 int nt[ND][NG][TBLNP];
812
814 int nu[ND][NG][TBLNP][TBLNT];
815
817 double p[ND][NG][TBLNP];
818
820 double t[ND][NG][TBLNP][TBLNT];
821
823 float u[ND][NG][TBLNP][TBLNT][TBLNU];
824
826 float eps[ND][NG][TBLNP][TBLNT][TBLNU];
827
829 double st[TBLNS];
830
832 double sr[TBLNS][ND];
833
834} tbl_t;
835
836/* ------------------------------------------------------------
837 Functions...
838 ------------------------------------------------------------ */
839
841size_t atm2x(
842 const ctl_t * ctl,
843 const atm_t * atm,
844 gsl_vector * x,
845 int *iqa,
846 int *ipa);
847
849void atm2x_help(
850 const double value,
851 const int value_iqa,
852 const int value_ip,
853 gsl_vector * x,
854 int *iqa,
855 int *ipa,
856 size_t *n);
857
859void cart2geo(
860 const double *x,
861 double *z,
862 double *lon,
863 double *lat);
864
866void climatology(
867 const ctl_t * ctl,
868 atm_t * atm_mean);
869
871double ctmco2(
872 const double nu,
873 const double p,
874 const double t,
875 const double u);
876
878double ctmh2o(
879 const double nu,
880 const double p,
881 const double t,
882 const double q,
883 const double u);
884
886double ctmn2(
887 const double nu,
888 const double p,
889 const double t);
890
892double ctmo2(
893 const double nu,
894 const double p,
895 const double t);
896
898void copy_atm(
899 const ctl_t * ctl,
900 atm_t * atm_dest,
901 const atm_t * atm_src,
902 const int init);
903
905void copy_obs(
906 const ctl_t * ctl,
907 obs_t * obs_dest,
908 const obs_t * obs_src,
909 const int init);
910
912int find_emitter(
913 const ctl_t * ctl,
914 const char *emitter);
915
917void formod(
918 const ctl_t * ctl,
919 const tbl_t * tbl,
920 atm_t * atm,
921 obs_t * obs);
922
924void formod_continua(
925 const ctl_t * ctl,
926 const los_t * los,
927 const int ip,
928 double *beta);
929
931void formod_fov(
932 const ctl_t * ctl,
933 obs_t * obs);
934
936void formod_pencil(
937 const ctl_t * ctl,
938 const tbl_t * tbl,
939 const atm_t * atm,
940 obs_t * obs,
941 const int ir);
942
944void formod_rfm(
945 const ctl_t * ctl,
946 const atm_t * atm,
947 obs_t * obs);
948
950void formod_srcfunc(
951 const ctl_t * ctl,
952 const tbl_t * tbl,
953 const double t,
954 double *src);
955
957void geo2cart(
958 const double z,
959 const double lon,
960 const double lat,
961 double *x);
962
964void hydrostatic(
965 const ctl_t * ctl,
966 atm_t * atm);
967
969void idx2name(
970 const ctl_t * ctl,
971 const int idx,
972 char *quantity);
973
975void init_srcfunc(
976 const ctl_t * ctl,
977 tbl_t * tbl);
978
980void intpol_atm(
981 const ctl_t * ctl,
982 const atm_t * atm,
983 const double z,
984 double *p,
985 double *t,
986 double *q,
987 double *k);
988
990void intpol_tbl_cga(
991 const ctl_t * ctl,
992 const tbl_t * tbl,
993 const los_t * los,
994 const int ip,
995 double tau_path[ND][NG],
996 double tau_seg[ND]);
997
999void intpol_tbl_ega(
1000 const ctl_t * ctl,
1001 const tbl_t * tbl,
1002 const los_t * los,
1003 const int ip,
1004 double tau_path[ND][NG],
1005 double tau_seg[ND]);
1006
1008double intpol_tbl_eps(
1009 const tbl_t * tbl,
1010 const int ig,
1011 const int id,
1012 const int ip,
1013 const int it,
1014 const double u);
1015
1017double intpol_tbl_u(
1018 const tbl_t * tbl,
1019 const int ig,
1020 const int id,
1021 const int ip,
1022 const int it,
1023 const double eps);
1024
1026void jsec2time(
1027 const double jsec,
1028 int *year,
1029 int *mon,
1030 int *day,
1031 int *hour,
1032 int *min,
1033 int *sec,
1034 double *remain);
1035
1037void kernel(
1038 const ctl_t * ctl,
1039 const tbl_t * tbl,
1040 atm_t * atm,
1041 obs_t * obs,
1042 gsl_matrix * k);
1043
1045int locate_irr(
1046 const double *xx,
1047 const int n,
1048 const double x);
1049
1051int locate_reg(
1052 const double *xx,
1053 const int n,
1054 const double x);
1055
1057int locate_tbl(
1058 const float *xx,
1059 const int n,
1060 const double x);
1061
1063size_t obs2y(
1064 const ctl_t * ctl,
1065 const obs_t * obs,
1066 gsl_vector * y,
1067 int *ida,
1068 int *ira);
1069
1071void raytrace(
1072 const ctl_t * ctl,
1073 const atm_t * atm,
1074 obs_t * obs,
1075 los_t * los,
1076 const int ir);
1077
1079void read_atm(
1080 const char *dirname,
1081 const char *filename,
1082 const ctl_t * ctl,
1083 atm_t * atm);
1084
1086void read_ctl(
1087 int argc,
1088 char *argv[],
1089 ctl_t * ctl);
1090
1092void read_matrix(
1093 const char *dirname,
1094 const char *filename,
1095 gsl_matrix * matrix);
1096
1098void read_obs(
1099 const char *dirname,
1100 const char *filename,
1101 const ctl_t * ctl,
1102 obs_t * obs);
1103
1105double read_obs_rfm(
1106 const char *basename,
1107 const double z,
1108 double *nu,
1109 double *f,
1110 int n);
1111
1113void read_rfm_spec(
1114 const char *filename,
1115 double *nu,
1116 double *rad,
1117 int *npts);
1118
1120void read_shape(
1121 const char *filename,
1122 double *x,
1123 double *y,
1124 int *n);
1125
1128 const ctl_t * ctl);
1129
1131double scan_ctl(
1132 int argc,
1133 char *argv[],
1134 const char *varname,
1135 const int arridx,
1136 const char *defvalue,
1137 char *value);
1138
1140double sza(
1141 double sec,
1142 double lon,
1143 double lat);
1144
1146void tangent_point(
1147 const los_t * los,
1148 double *tpz,
1149 double *tplon,
1150 double *tplat);
1151
1153void time2jsec(
1154 const int year,
1155 const int mon,
1156 const int day,
1157 const int hour,
1158 const int min,
1159 const int sec,
1160 const double remain,
1161 double *jsec);
1162
1164void timer(
1165 const char *name,
1166 const char *file,
1167 const char *func,
1168 int line,
1169 int mode);
1170
1172void write_atm(
1173 const char *dirname,
1174 const char *filename,
1175 const ctl_t * ctl,
1176 const atm_t * atm);
1177
1179void write_atm_rfm(
1180 const char *filename,
1181 const ctl_t * ctl,
1182 const atm_t * atm);
1183
1185void write_matrix(
1186 const char *dirname,
1187 const char *filename,
1188 const ctl_t * ctl,
1189 const gsl_matrix * matrix,
1190 const atm_t * atm,
1191 const obs_t * obs,
1192 const char *rowspace,
1193 const char *colspace,
1194 const char *sort);
1195
1197void write_obs(
1198 const char *dirname,
1199 const char *filename,
1200 const ctl_t * ctl,
1201 const obs_t * obs);
1202
1204void write_shape(
1205 const char *filename,
1206 const double *x,
1207 const double *y,
1208 const int n);
1209
1211void write_tbl(
1212 const ctl_t * ctl,
1213 const tbl_t * tbl);
1214
1216void x2atm(
1217 const ctl_t * ctl,
1218 const gsl_vector * x,
1219 atm_t * atm);
1220
1222void x2atm_help(
1223 double *value,
1224 const gsl_vector * x,
1225 size_t *n);
1226
1228void y2obs(
1229 const ctl_t * ctl,
1230 const gsl_vector * y,
1231 obs_t * obs);
1232
1233#endif
#define LEN
Maximum length of ASCII data lines.
Definition: jurassic.h:394
void read_matrix(const char *dirname, const char *filename, gsl_matrix *matrix)
Read matrix.
Definition: jurassic.c:4802
void timer(const char *name, const char *file, const char *func, int line, int mode)
Measure wall-clock time.
Definition: jurassic.c:5462
void read_rfm_spec(const char *filename, double *nu, double *rad, int *npts)
Read RFM spectrum.
Definition: jurassic.c:5001
void write_atm(const char *dirname, const char *filename, const ctl_t *ctl, const atm_t *atm)
Write atmospheric data.
Definition: jurassic.c:5500
int locate_reg(const double *xx, const int n, const double x)
Find array index for regular grid.
Definition: jurassic.c:4271
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])
Get transmittance from look-up tables (EGA method).
Definition: jurassic.c:3949
double read_obs_rfm(const char *basename, const double z, double *nu, double *f, int n)
Read observation data in RFM format.
Definition: jurassic.c:4943
void formod_rfm(const ctl_t *ctl, const atm_t *atm, obs_t *obs)
Apply RFM for radiative transfer calculations.
Definition: jurassic.c:3498
double ctmo2(const double nu, const double p, const double t)
Compute oxygen continuum (absorption coefficient).
Definition: jurassic.c:3045
void idx2name(const ctl_t *ctl, const int idx, char *quantity)
Determine name of state vector quantity for given index.
Definition: jurassic.c:3736
void read_ctl(int argc, char *argv[], ctl_t *ctl)
Read forward model control parameters.
Definition: jurassic.c:4686
void formod_continua(const ctl_t *ctl, const los_t *los, const int ip, double *beta)
Compute absorption coefficient of continua.
Definition: jurassic.c:3262
void raytrace(const ctl_t *ctl, const atm_t *atm, obs_t *obs, los_t *los, const int ir)
Do ray-tracing to determine LOS.
Definition: jurassic.c:4339
int locate_irr(const double *xx, const int n, const double x)
Find array index for irregular grid.
Definition: jurassic.c:4241
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)
Convert date to seconds.
Definition: jurassic.c:5431
void x2atm(const ctl_t *ctl, const gsl_vector *x, atm_t *atm)
Decompose parameter vector or state vector.
Definition: jurassic.c:6060
void atm2x_help(const double value, const int value_iqa, const int value_ip, gsl_vector *x, int *iqa, int *ipa, size_t *n)
Add element to state vector.
Definition: jurassic.c:87
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:1924
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])
Get transmittance from look-up tables (CGA method).
Definition: jurassic.c:3860
void write_atm_rfm(const char *filename, const ctl_t *ctl, const atm_t *atm)
Write atmospheric data in RFM format.
Definition: jurassic.c:5618
#define ND
Maximum number of radiance channels.
Definition: jurassic.h:364
double intpol_tbl_u(const tbl_t *tbl, const int ig, const int id, const int ip, const int it, const double eps)
Interpolate column density from look-up tables.
Definition: jurassic.c:4077
double intpol_tbl_eps(const tbl_t *tbl, const int ig, const int id, const int ip, const int it, const double u)
Interpolate emissivity from look-up tables.
Definition: jurassic.c:4045
#define NSHAPE
Maximum number of shape function grid points.
Definition: jurassic.h:419
void write_tbl(const ctl_t *ctl, const tbl_t *tbl)
Write look-up table data.
Definition: jurassic.c:5974
void intpol_atm(const ctl_t *ctl, const atm_t *atm, const double z, double *p, double *t, double *q, double *k)
Interpolate atmospheric data.
Definition: jurassic.c:3835
int find_emitter(const ctl_t *ctl, const char *emitter)
Find index of an emitter.
Definition: jurassic.c:3199
void tangent_point(const los_t *los, double *tpz, double *tplon, double *tplat)
Find tangent point of a given LOS.
Definition: jurassic.c:5387
#define TBLNU
Maximum number of column densities in emissivity tables.
Definition: jurassic.h:439
void climatology(const ctl_t *ctl, atm_t *atm_mean)
Interpolate climatological data.
Definition: jurassic.c:123
void formod_pencil(const ctl_t *ctl, const tbl_t *tbl, const atm_t *atm, obs_t *obs, const int ir)
Compute radiative transfer for a pencil beam.
Definition: jurassic.c:3363
double ctmco2(const double nu, const double p, const double t, const double u)
Compute carbon dioxide continuum (optical depth).
Definition: jurassic.c:1061
void read_atm(const char *dirname, const char *filename, const ctl_t *ctl, atm_t *atm)
Read atmospheric data.
Definition: jurassic.c:4581
void copy_obs(const ctl_t *ctl, obs_t *obs_dest, const obs_t *obs_src, const int init)
Copy and initialize observation data.
Definition: jurassic.c:3161
void read_obs(const char *dirname, const char *filename, const ctl_t *ctl, obs_t *obs)
Read observation data.
Definition: jurassic.c:4842
void write_obs(const char *dirname, const char *filename, const ctl_t *ctl, const obs_t *obs)
Write observation data.
Definition: jurassic.c:5834
#define TBLNT
Maximum number of temperatures in emissivity tables.
Definition: jurassic.h:434
void x2atm_help(double *value, const gsl_vector *x, size_t *n)
Get element from state vector.
Definition: jurassic.c:6114
void cart2geo(const double *x, double *z, double *lon, double *lat)
Convert Cartesian coordinates to geolocation.
Definition: jurassic.c:108
#define NP
Maximum number of atmospheric data points.
Definition: jurassic.h:374
void formod_fov(const ctl_t *ctl, obs_t *obs)
Apply field of view convolution.
Definition: jurassic.c:3298
void kernel(const ctl_t *ctl, const tbl_t *tbl, atm_t *atm, obs_t *obs, gsl_matrix *k)
Compute Jacobians.
Definition: jurassic.c:4142
void init_srcfunc(const ctl_t *ctl, tbl_t *tbl)
Initialize source function table.
Definition: jurassic.c:3781
double scan_ctl(int argc, char *argv[], const char *varname, const int arridx, const char *defvalue, char *value)
Search control parameter file for variable entry.
Definition: jurassic.c:5278
void write_shape(const char *filename, const double *x, const double *y, const int n)
Write shape function.
Definition: jurassic.c:5944
void formod(const ctl_t *ctl, const tbl_t *tbl, atm_t *atm, obs_t *obs)
Determine ray paths and compute radiative transfer.
Definition: jurassic.c:3212
#define NG
Maximum number of emitters.
Definition: jurassic.h:369
void formod_srcfunc(const ctl_t *ctl, const tbl_t *tbl, const double t, double *src)
Compute Planck source function.
Definition: jurassic.c:3639
void jsec2time(const double jsec, int *year, int *mon, int *day, int *hour, int *min, int *sec, double *remain)
Convert seconds to date.
Definition: jurassic.c:4109
#define TBLNS
Maximum number of source function temperature levels.
Definition: jurassic.h:444
tbl_t * read_tbl(const ctl_t *ctl)
Read look-up table data.
Definition: jurassic.c:5104
double sza(double sec, double lon, double lat)
Calculate solar zenith angle.
Definition: jurassic.c:5347
void copy_atm(const ctl_t *ctl, atm_t *atm_dest, const atm_t *atm_src, const int init)
Copy and initialize atmospheric data.
Definition: jurassic.c:3107
size_t obs2y(const ctl_t *ctl, const obs_t *obs, gsl_vector *y, int *ida, int *ira)
Compose measurement vector.
Definition: jurassic.c:4312
void y2obs(const ctl_t *ctl, const gsl_vector *y, obs_t *obs)
Decompose measurement vector.
Definition: jurassic.c:6126
void hydrostatic(const ctl_t *ctl, atm_t *atm)
Set hydrostatic equilibrium.
Definition: jurassic.c:3676
void read_shape(const char *filename, double *x, double *y, int *n)
Read shape function.
Definition: jurassic.c:5062
double ctmn2(const double nu, const double p, const double t)
Compute nitrogen continuum (absorption coefficient).
Definition: jurassic.c:2976
size_t atm2x(const ctl_t *ctl, const atm_t *atm, gsl_vector *x, int *iqa, int *ipa)
Compose state vector or parameter vector.
Definition: jurassic.c:29
#define TBLNP
Maximum number of pressure levels in emissivity tables.
Definition: jurassic.h:429
void geo2cart(const double z, const double lon, const double lat, double *x)
Convert geolocation to Cartesian coordinates.
Definition: jurassic.c:3656
#define NSF
Maximum number of surface layer spectral grid points.
Definition: jurassic.h:384
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 matrix.
Definition: jurassic.c:5656
#define NCL
Maximum number of cloud layer spectral grid points.
Definition: jurassic.h:359
#define NLOS
Maximum number of LOS points.
Definition: jurassic.h:414
#define NR
Maximum number of ray paths.
Definition: jurassic.h:379
int locate_tbl(const float *xx, const int n, const double x)
Find array index in float array.
Definition: jurassic.c:4290
#define NW
Maximum number of spectral windows.
Definition: jurassic.h:389
Atmospheric data.
Definition: jurassic.h:494
double sfz
Surface height [km].
Definition: jurassic.h:533
double sfp
Surface pressure [hPa].
Definition: jurassic.h:536
double clz
Cloud layer height [km].
Definition: jurassic.h:524
int np
Number of data points.
Definition: jurassic.h:497
double cldz
Cloud layer depth [km].
Definition: jurassic.h:527
double sft
Surface temperature [K].
Definition: jurassic.h:539
Forward model control parameters.
Definition: jurassic.h:547
int write_matrix
Write matrix file (0=no, 1=yes).
Definition: jurassic.h:688
int nw
Number of spectral windows.
Definition: jurassic.h:574
double retp_zmin
Minimum altitude for pressure retrieval [km].
Definition: jurassic.h:640
int ig_co2
Emitter index of CO2.
Definition: jurassic.h:556
double hydz
Reference height for hydrostatic pressure profile (-999 to skip) [km].
Definition: jurassic.h:604
int ret_sfz
Retrieve surface layer height (0=no, 1=yes).
Definition: jurassic.h:673
int ctm_co2
Compute CO2 continuum (0=no, 1=yes).
Definition: jurassic.h:607
double rett_zmax
Maximum altitude for temperature retrieval [km].
Definition: jurassic.h:649
int ret_sfeps
Retrieve surface layer emissivity (0=no, 1=yes).
Definition: jurassic.h:682
int ret_sft
Retrieve surface layer temperature (0=no, 1=yes).
Definition: jurassic.h:679
int ret_clz
Retrieve cloud layer height (0=no, 1=yes).
Definition: jurassic.h:664
int ctm_n2
Compute N2 continuum (0=no, 1=yes).
Definition: jurassic.h:613
int ctm_h2o
Compute H2O continuum (0=no, 1=yes).
Definition: jurassic.h:610
int formod
Forward model (0=CGA, 1=EGA, 2=RFM).
Definition: jurassic.h:691
int ng
Number of emitters.
Definition: jurassic.h:550
int refrac
Take into account refractivity (0=no, 1=yes).
Definition: jurassic.h:619
int ig_o2
Emitter index of O2.
Definition: jurassic.h:565
double rett_zmin
Minimum altitude for temperature retrieval [km].
Definition: jurassic.h:646
double sfsza
Solar zenith angle at the surface [deg] (-999=auto).
Definition: jurassic.h:595
int nd
Number of radiance channels.
Definition: jurassic.h:568
int fov_n
Field-of-view number of data points.
Definition: jurassic.h:637
int sftype
Surface treatment (0=none, 1=emissions, 2=downward, 3=solar).
Definition: jurassic.h:592
int ret_sfp
Retrieve surface layer pressure (0=no, 1=yes).
Definition: jurassic.h:676
int ncl
Number of cloud layer spectral grid points.
Definition: jurassic.h:580
int ig_n2
Emitter index of N2.
Definition: jurassic.h:562
int ctm_o2
Compute O2 continuum (0=no, 1=yes).
Definition: jurassic.h:616
int ret_clk
Retrieve cloud layer extinction (0=no, 1=yes).
Definition: jurassic.h:670
int nsf
Number of surface layer spectral grid points.
Definition: jurassic.h:586
double rayds
Maximum step length for raytracing [km].
Definition: jurassic.h:622
int ret_cldz
Retrieve cloud layer depth (0=no, 1=yes).
Definition: jurassic.h:667
int ig_h2o
Emitter index of H2O.
Definition: jurassic.h:559
int tblfmt
Look-up table file format (1=ASCII, 2=binary).
Definition: jurassic.h:601
double raydz
Vertical step length for raytracing [km].
Definition: jurassic.h:625
int write_bbt
Use brightness temperature instead of radiance (0=no, 1=yes).
Definition: jurassic.h:685
double retp_zmax
Maximum altitude for pressure retrieval [km].
Definition: jurassic.h:643
Line-of-sight data.
Definition: jurassic.h:705
double sft
Surface temperature [K].
Definition: jurassic.h:732
int np
Number of LOS points.
Definition: jurassic.h:708
Observation geometry and radiance data.
Definition: jurassic.h:761
int nr
Number of ray paths.
Definition: jurassic.h:764
Emissivity look-up tables.
Definition: jurassic.h:805