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-2024 Forschungszentrum Juelich GmbH
18*/
19
100#ifndef JURASSIC_H
101#define JURASSIC_H
102
103#include <gsl/gsl_math.h>
104#include <gsl/gsl_blas.h>
105#include <gsl/gsl_linalg.h>
106#include <gsl/gsl_randist.h>
107#include <gsl/gsl_rng.h>
108#include <gsl/gsl_statistics.h>
109#include <math.h>
110#include <omp.h>
111#include <stdio.h>
112#include <stdlib.h>
113#include <string.h>
114#include <time.h>
115
116/* ------------------------------------------------------------
117 Macros...
118 ------------------------------------------------------------ */
119
121#define ALLOC(ptr, type, n) \
122 if((ptr=malloc((size_t)(n)*sizeof(type)))==NULL) \
123 ERRMSG("Out of memory!");
124
126#define BRIGHT(rad, nu) \
127 (C2 * (nu) / gsl_log1p(C1 * POW3(nu) / (rad)))
128
130#define DEG2RAD(deg) \
131 ((deg) * (M_PI / 180.0))
132
134#define DIST(a, b) sqrt(DIST2(a, b))
135
137#define DIST2(a, b) \
138 ((a[0]-b[0])*(a[0]-b[0])+(a[1]-b[1])*(a[1]-b[1])+(a[2]-b[2])*(a[2]-b[2]))
139
141#define DOTP(a, b) (a[0]*b[0]+a[1]*b[1]+a[2]*b[2])
142
144#define FREAD(ptr, type, size, out) { \
145 if(fread(ptr, sizeof(type), size, out)!=size) \
146 ERRMSG("Error while reading!"); \
147 }
148
150#define FWRITE(ptr, type, size, out) { \
151 if(fwrite(ptr, sizeof(type), size, out)!=size) \
152 ERRMSG("Error while writing!"); \
153 }
154
156#define MAX(a,b) \
157 (((a)>(b))?(a):(b))
158
160#define MIN(a,b) \
161 (((a)<(b))?(a):(b))
162
164#define LIN(x0, y0, x1, y1, x) \
165 ((y0)+((y1)-(y0))/((x1)-(x0))*((x)-(x0)))
166
168#define LOGX(x0, y0, x1, y1, x) \
169 (((x)/(x0)>0 && (x1)/(x0)>0) \
170 ? ((y0)+((y1)-(y0))*log((x)/(x0))/log((x1)/(x0))) \
171 : LIN(x0, y0, x1, y1, x))
172
174#define LOGY(x0, y0, x1, y1, x) \
175 (((y1)/(y0)>0) \
176 ? ((y0)*exp(log((y1)/(y0))/((x1)-(x0))*((x)-(x0)))) \
177 : LIN(x0, y0, x1, y1, x))
178
180#define NORM(a) sqrt(DOTP(a, a))
181
183#define PLANCK(T, nu) \
184 (C1 * POW3(nu) / gsl_expm1(C2 * (nu) / (T)))
185
187#define POW2(x) ((x)*(x))
188
190#define POW3(x) ((x)*(x)*(x))
191
193#define RAD2DEG(rad) \
194 ((rad) * (180.0 / M_PI))
195
197#define REFRAC(p, T) \
198 (7.753e-05 * (p) / (T))
199
201#define TIMER(name, mode) \
202 {timer(name, __FILE__, __func__, __LINE__, mode);}
203
205#define TOK(line, tok, format, var) { \
206 if(((tok)=strtok((line), " \t"))) { \
207 if(sscanf(tok, format, &(var))!=1) continue; \
208 } else ERRMSG("Error while reading!"); \
209 }
210
211/* ------------------------------------------------------------
212 Log messages...
213 ------------------------------------------------------------ */
214
216#ifndef LOGLEV
217#define LOGLEV 2
218#endif
219
221#define LOG(level, ...) { \
222 if(level >= 2) \
223 printf(" "); \
224 if(level <= LOGLEV) { \
225 printf(__VA_ARGS__); \
226 printf("\n"); \
227 } \
228 }
229
231#define WARN(...) { \
232 printf("\nWarning (%s, %s, l%d): ", __FILE__, __func__, __LINE__); \
233 LOG(0, __VA_ARGS__); \
234 }
235
237#define ERRMSG(...) { \
238 printf("\nError (%s, %s, l%d): ", __FILE__, __func__, __LINE__); \
239 LOG(0, __VA_ARGS__); \
240 exit(EXIT_FAILURE); \
241 }
242
244#define PRINT(format, var) \
245 printf("Print (%s, %s, l%d): %s= "format"\n", \
246 __FILE__, __func__, __LINE__, #var, var);
247
248/* ------------------------------------------------------------
249 Constants...
250 ------------------------------------------------------------ */
251
253#ifndef C1
254#define C1 1.19104259e-8
255#endif
256
258#ifndef C2
259#define C2 1.43877506
260#endif
261
263#ifndef EPSMIN
264#define EPSMIN 0
265#endif
266
268#ifndef EPSMAX
269#define EPSMAX 1
270#endif
271
273#ifndef G0
274#define G0 9.80665
275#endif
276
278#ifndef H0
279#define H0 7.0
280#endif
281
283#ifndef KB
284#define KB 1.3806504e-23
285#endif
286
288#ifndef ME
289#define ME 5.976e24
290#endif
291
293#ifndef NA
294#define NA 6.02214199e23
295#endif
296
298#ifndef RE
299#define RE 6367.421
300#endif
301
303#ifndef RI
304#define RI 8.3144598
305#endif
306
308#ifndef P0
309#define P0 1013.25
310#endif
311
313#ifndef T0
314#define T0 273.15
315#endif
316
318#ifndef TMIN
319#define TMIN 100.
320#endif
321
323#ifndef TMAX
324#define TMAX 400.
325#endif
326
328#ifndef TSUN
329#define TSUN 5780.
330#endif
331
333#ifndef UMIN
334#define UMIN 0
335#endif
336
338#ifndef UMAX
339#define UMAX 1e30
340#endif
341
342/* ------------------------------------------------------------
343 Dimensions...
344 ------------------------------------------------------------ */
345
347#ifndef NCL
348#define NCL 8
349#endif
350
352#ifndef ND
353#define ND 128
354#endif
355
357#ifndef NG
358#define NG 8
359#endif
360
362#ifndef NP
363#define NP 256
364#endif
365
367#ifndef NR
368#define NR 256
369#endif
370
372#ifndef NSF
373#define NSF 8
374#endif
375
377#ifndef NW
378#define NW 4
379#endif
380
382#ifndef LEN
383#define LEN 10000
384#endif
385
387#ifndef M
388#define M (NR*ND)
389#endif
390
392#ifndef N
393#define N ((2+NG+NW)*NP+NCL+NSF+5)
394#endif
395
397#ifndef NQ
398#define NQ (7+NG+NW+NCL+NSF)
399#endif
400
402#ifndef NLOS
403#define NLOS 4096
404#endif
405
407#ifndef NSHAPE
408#define NSHAPE 20000
409#endif
410
412#ifndef NFOV
413#define NFOV 5
414#endif
415
417#ifndef TBLNP
418#define TBLNP 41
419#endif
420
422#ifndef TBLNT
423#define TBLNT 30
424#endif
425
427#ifndef TBLNU
428#define TBLNU 320
429#endif
430
432#ifndef TBLNS
433#define TBLNS 1200
434#endif
435
437#ifndef RFMNPTS
438#define RFMNPTS 10000000
439#endif
440
442#ifndef RFMLINE
443#define RFMLINE 100000
444#endif
445
446/* ------------------------------------------------------------
447 Quantity indices...
448 ------------------------------------------------------------ */
449
451#define IDXP 0
452
454#define IDXT 1
455
457#define IDXQ(ig) (2+ig)
458
460#define IDXK(iw) (2+ctl->ng+iw)
461
463#define IDXCLZ (2+ctl->ng+ctl->nw)
464
466#define IDXCLDZ (3+ctl->ng+ctl->nw)
467
469#define IDXCLK(icl) (4+ctl->ng+ctl->nw+icl)
470
472#define IDXSFZ (4+ctl->ng+ctl->nw+ctl->ncl)
473
475#define IDXSFP (5+ctl->ng+ctl->nw+ctl->ncl)
476
478#define IDXSFT (6+ctl->ng+ctl->nw+ctl->ncl)
479
481#define IDXSFEPS(isf) (7+ctl->ng+ctl->nw+ctl->ncl+isf)
482
483/* ------------------------------------------------------------
484 Structs...
485 ------------------------------------------------------------ */
486
488typedef struct {
489
491 int np;
492
494 double time[NP];
495
497 double z[NP];
498
500 double lon[NP];
501
503 double lat[NP];
504
506 double p[NP];
507
509 double t[NP];
510
512 double q[NG][NP];
513
515 double k[NW][NP];
516
518 double clz;
519
521 double cldz;
522
524 double clk[NCL];
525
527 double sfz;
528
530 double sfp;
531
533 double sft;
534
536 double sfeps[NSF];
537
538} atm_t;
539
541typedef struct {
542
544 int ng;
545
547 char emitter[NG][LEN];
548
550 int nd;
551
553 double nu[ND];
554
556 int nw;
557
559 int window[ND];
560
562 int ncl;
563
565 double clnu[NCL];
566
568 int nsf;
569
571 double sfnu[NSF];
572
575
577 double sfsza;
578
580 char tblbase[LEN];
581
584
586 double hydz;
587
590
593
596
599
602
604 double rayds;
605
607 double raydz;
608
610 char fov[LEN];
611
613 double retp_zmin;
614
616 double retp_zmax;
617
619 double rett_zmin;
620
622 double rett_zmax;
623
625 double retq_zmin[NG];
626
628 double retq_zmax[NG];
629
631 double retk_zmin[NW];
632
634 double retk_zmax[NW];
635
638
641
644
647
650
653
656
659
662
665
667 char rfmbin[LEN];
668
670 char rfmhit[LEN];
671
673 char rfmxsc[NG][LEN];
674
675} ctl_t;
676
678typedef struct {
679
681 int np;
682
684 double z[NLOS];
685
687 double lon[NLOS];
688
690 double lat[NLOS];
691
693 double p[NLOS];
694
696 double t[NLOS];
697
699 double q[NLOS][NG];
700
702 double k[NLOS][ND];
703
705 double sft;
706
708 double sfeps[ND];
709
711 double ds[NLOS];
712
714 double u[NLOS][NG];
715
717 double cgp[NLOS][NG];
718
720 double cgt[NLOS][NG];
721
723 double cgu[NLOS][NG];
724
726 double eps[NLOS][ND];
727
729 double src[NLOS][ND];
730
731} los_t;
732
734typedef struct {
735
737 int nr;
738
740 double time[NR];
741
743 double obsz[NR];
744
746 double obslon[NR];
747
749 double obslat[NR];
750
752 double vpz[NR];
753
755 double vplon[NR];
756
758 double vplat[NR];
759
761 double tpz[NR];
762
764 double tplon[NR];
765
767 double tplat[NR];
768
770 double tau[ND][NR];
771
773 double rad[ND][NR];
774
775} obs_t;
776
778typedef struct {
779
781 int np[ND][NG];
782
784 int nt[ND][NG][TBLNP];
785
787 int nu[ND][NG][TBLNP][TBLNT];
788
790 double p[ND][NG][TBLNP];
791
793 double t[ND][NG][TBLNP][TBLNT];
794
796 float u[ND][NG][TBLNP][TBLNT][TBLNU];
797
799 float eps[ND][NG][TBLNP][TBLNT][TBLNU];
800
802 double st[TBLNS];
803
805 double sr[TBLNS][ND];
806
807} tbl_t;
808
809/* ------------------------------------------------------------
810 Functions...
811 ------------------------------------------------------------ */
812
814size_t atm2x(
815 const ctl_t * ctl,
816 const atm_t * atm,
817 gsl_vector * x,
818 int *iqa,
819 int *ipa);
820
822void atm2x_help(
823 const double value,
824 const int value_iqa,
825 const int value_ip,
826 gsl_vector * x,
827 int *iqa,
828 int *ipa,
829 size_t *n);
830
832void cart2geo(
833 const double *x,
834 double *z,
835 double *lon,
836 double *lat);
837
839void climatology(
840 const ctl_t * ctl,
841 atm_t * atm_mean);
842
844double ctmco2(
845 const double nu,
846 const double p,
847 const double t,
848 const double u);
849
851double ctmh2o(
852 const double nu,
853 const double p,
854 const double t,
855 const double q,
856 const double u);
857
859double ctmn2(
860 const double nu,
861 const double p,
862 const double t);
863
865double ctmo2(
866 const double nu,
867 const double p,
868 const double t);
869
871void copy_atm(
872 const ctl_t * ctl,
873 atm_t * atm_dest,
874 const atm_t * atm_src,
875 const int init);
876
878void copy_obs(
879 const ctl_t * ctl,
880 obs_t * obs_dest,
881 const obs_t * obs_src,
882 const int init);
883
885int find_emitter(
886 const ctl_t * ctl,
887 const char *emitter);
888
890void formod(
891 const ctl_t * ctl,
892 atm_t * atm,
893 obs_t * obs);
894
896void formod_continua(
897 const ctl_t * ctl,
898 const los_t * los,
899 const int ip,
900 double *beta);
901
903void formod_fov(
904 const ctl_t * ctl,
905 obs_t * obs);
906
908void formod_pencil(
909 const ctl_t * ctl,
910 const atm_t * atm,
911 obs_t * obs,
912 const int ir);
913
915void formod_rfm(
916 const ctl_t * ctl,
917 const atm_t * atm,
918 obs_t * obs);
919
921void formod_srcfunc(
922 const ctl_t * ctl,
923 const tbl_t * tbl,
924 const double t,
925 double *src);
926
928void geo2cart(
929 const double z,
930 const double lon,
931 const double lat,
932 double *x);
933
935void hydrostatic(
936 const ctl_t * ctl,
937 atm_t * atm);
938
940void idx2name(
941 const ctl_t * ctl,
942 const int idx,
943 char *quantity);
944
946void init_srcfunc(
947 const ctl_t * ctl,
948 tbl_t * tbl);
949
951void intpol_atm(
952 const ctl_t * ctl,
953 const atm_t * atm,
954 const double z,
955 double *p,
956 double *t,
957 double *q,
958 double *k);
959
961void intpol_tbl_cga(
962 const ctl_t * ctl,
963 const tbl_t * tbl,
964 const los_t * los,
965 const int ip,
966 double tau_path[ND][NG],
967 double tau_seg[ND]);
968
970void intpol_tbl_ega(
971 const ctl_t * ctl,
972 const tbl_t * tbl,
973 const los_t * los,
974 const int ip,
975 double tau_path[ND][NG],
976 double tau_seg[ND]);
977
979double intpol_tbl_eps(
980 const tbl_t * tbl,
981 const int ig,
982 const int id,
983 const int ip,
984 const int it,
985 const double u);
986
988double intpol_tbl_u(
989 const tbl_t * tbl,
990 const int ig,
991 const int id,
992 const int ip,
993 const int it,
994 const double eps);
995
997void jsec2time(
998 const double jsec,
999 int *year,
1000 int *mon,
1001 int *day,
1002 int *hour,
1003 int *min,
1004 int *sec,
1005 double *remain);
1006
1008void kernel(
1009 ctl_t * ctl,
1010 atm_t * atm,
1011 obs_t * obs,
1012 gsl_matrix * k);
1013
1015int locate_irr(
1016 const double *xx,
1017 const int n,
1018 const double x);
1019
1021int locate_reg(
1022 const double *xx,
1023 const int n,
1024 const double x);
1025
1027int locate_tbl(
1028 const float *xx,
1029 const int n,
1030 const double x);
1031
1033size_t obs2y(
1034 const ctl_t * ctl,
1035 const obs_t * obs,
1036 gsl_vector * y,
1037 int *ida,
1038 int *ira);
1039
1041void raytrace(
1042 const ctl_t * ctl,
1043 const atm_t * atm,
1044 obs_t * obs,
1045 los_t * los,
1046 const int ir);
1047
1049void read_atm(
1050 const char *dirname,
1051 const char *filename,
1052 const ctl_t * ctl,
1053 atm_t * atm);
1054
1056void read_ctl(
1057 int argc,
1058 char *argv[],
1059 ctl_t * ctl);
1060
1062void read_matrix(
1063 const char *dirname,
1064 const char *filename,
1065 gsl_matrix * matrix);
1066
1068void read_obs(
1069 const char *dirname,
1070 const char *filename,
1071 const ctl_t * ctl,
1072 obs_t * obs);
1073
1075double read_obs_rfm(
1076 const char *basename,
1077 const double z,
1078 double *nu,
1079 double *f,
1080 int n);
1081
1083void read_rfm_spec(
1084 const char *filename,
1085 double *nu,
1086 double *rad,
1087 int *npts);
1088
1090void read_shape(
1091 const char *filename,
1092 double *x,
1093 double *y,
1094 int *n);
1095
1097void read_tbl(
1098 const ctl_t * ctl,
1099 tbl_t * tbl);
1100
1102double scan_ctl(
1103 int argc,
1104 char *argv[],
1105 const char *varname,
1106 int arridx,
1107 const char *defvalue,
1108 char *value);
1109
1111double sza(
1112 double sec,
1113 double lon,
1114 double lat);
1115
1117void tangent_point(
1118 const los_t * los,
1119 double *tpz,
1120 double *tplon,
1121 double *tplat);
1122
1124void time2jsec(
1125 const int year,
1126 const int mon,
1127 const int day,
1128 const int hour,
1129 const int min,
1130 const int sec,
1131 const double remain,
1132 double *jsec);
1133
1135void timer(
1136 const char *name,
1137 const char *file,
1138 const char *func,
1139 int line,
1140 int mode);
1141
1143void write_atm(
1144 const char *dirname,
1145 const char *filename,
1146 const ctl_t * ctl,
1147 const atm_t * atm);
1148
1150void write_atm_rfm(
1151 const char *filename,
1152 const ctl_t * ctl,
1153 const atm_t * atm);
1154
1156void write_matrix(
1157 const char *dirname,
1158 const char *filename,
1159 const ctl_t * ctl,
1160 const gsl_matrix * matrix,
1161 const atm_t * atm,
1162 const obs_t * obs,
1163 const char *rowspace,
1164 const char *colspace,
1165 const char *sort);
1166
1168void write_obs(
1169 const char *dirname,
1170 const char *filename,
1171 const ctl_t * ctl,
1172 const obs_t * obs);
1173
1175void write_shape(
1176 const char *filename,
1177 const double *x,
1178 const double *y,
1179 const int n);
1180
1182void write_tbl(
1183 const ctl_t * ctl,
1184 const tbl_t * tbl);
1185
1187void x2atm(
1188 const ctl_t * ctl,
1189 const gsl_vector * x,
1190 atm_t * atm);
1191
1193void x2atm_help(
1194 double *value,
1195 const gsl_vector * x,
1196 size_t *n);
1197
1199void y2obs(
1200 const ctl_t * ctl,
1201 const gsl_vector * y,
1202 obs_t * obs);
1203
1204#endif
#define LEN
Maximum length of ASCII data lines.
Definition: jurassic.h:383
void read_matrix(const char *dirname, const char *filename, gsl_matrix *matrix)
Read matrix.
Definition: jurassic.c:4660
void timer(const char *name, const char *file, const char *func, int line, int mode)
Measure wall-clock time.
Definition: jurassic.c:5301
void read_rfm_spec(const char *filename, double *nu, double *rad, int *npts)
Read RFM spectrum.
Definition: jurassic.c:4859
void read_tbl(const ctl_t *ctl, tbl_t *tbl)
Read look-up table data.
Definition: jurassic.c:4952
void write_atm(const char *dirname, const char *filename, const ctl_t *ctl, const atm_t *atm)
Write atmospheric data.
Definition: jurassic.c:5339
int locate_reg(const double *xx, const int n, const double x)
Find array index for regular grid.
Definition: jurassic.c:4137
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:3799
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:4801
void formod_rfm(const ctl_t *ctl, const atm_t *atm, obs_t *obs)
Apply RFM for radiative transfer calculations.
Definition: jurassic.c:3342
double ctmo2(const double nu, const double p, const double t)
Compute oxygen continuum (absorption coefficient).
Definition: jurassic.c:2861
void idx2name(const ctl_t *ctl, const int idx, char *quantity)
Determine name of state vector quantity for given index.
Definition: jurassic.c:3586
void read_ctl(int argc, char *argv[], ctl_t *ctl)
Read forward model control parameters.
Definition: jurassic.c:4550
void formod_continua(const ctl_t *ctl, const los_t *los, const int ip, double *beta)
Compute absorption coefficient of continua.
Definition: jurassic.c:3075
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:4205
int locate_irr(const double *xx, const int n, const double x)
Find array index for irregular grid.
Definition: jurassic.c:4107
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:5270
void x2atm(const ctl_t *ctl, const gsl_vector *x, atm_t *atm)
Decompose parameter vector or state vector.
Definition: jurassic.c:5899
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:1746
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:3710
void write_atm_rfm(const char *filename, const ctl_t *ctl, const atm_t *atm)
Write atmospheric data in RFM format.
Definition: jurassic.c:5457
#define ND
Maximum number of radiance channels.
Definition: jurassic.h:353
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:3933
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:3894
void write_tbl(const ctl_t *ctl, const tbl_t *tbl)
Write look-up table data.
Definition: jurassic.c:5813
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:3685
int find_emitter(const ctl_t *ctl, const char *emitter)
Find index of an emitter.
Definition: jurassic.c:3013
void formod_pencil(const ctl_t *ctl, const atm_t *atm, obs_t *obs, const int ir)
Compute radiative transfer for a pencil beam.
Definition: jurassic.c:3196
void tangent_point(const los_t *los, double *tpz, double *tplon, double *tplat)
Find tangent point of a given LOS.
Definition: jurassic.c:5226
#define TBLNU
Maximum number of column densities in emissivity tables.
Definition: jurassic.h:428
void formod(const ctl_t *ctl, atm_t *atm, obs_t *obs)
Determine ray paths and compute radiative transfer.
Definition: jurassic.c:3026
void climatology(const ctl_t *ctl, atm_t *atm_mean)
Interpolate climatological data.
Definition: jurassic.c:123
double ctmco2(const double nu, const double p, const double t, const double u)
Compute carbon dioxide continuum (optical depth).
Definition: jurassic.c:886
void read_atm(const char *dirname, const char *filename, const ctl_t *ctl, atm_t *atm)
Read atmospheric data.
Definition: jurassic.c:4445
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:2975
void read_obs(const char *dirname, const char *filename, const ctl_t *ctl, obs_t *obs)
Read observation data.
Definition: jurassic.c:4700
void write_obs(const char *dirname, const char *filename, const ctl_t *ctl, const obs_t *obs)
Write observation data.
Definition: jurassic.c:5673
#define TBLNT
Maximum number of temperatures in emissivity tables.
Definition: jurassic.h:423
void x2atm_help(double *value, const gsl_vector *x, size_t *n)
Get element from state vector.
Definition: jurassic.c:5953
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:363
void formod_fov(const ctl_t *ctl, obs_t *obs)
Apply field of view convolution.
Definition: jurassic.c:3121
void init_srcfunc(const ctl_t *ctl, tbl_t *tbl)
Initialize source function table.
Definition: jurassic.c:3631
void write_shape(const char *filename, const double *x, const double *y, const int n)
Write shape function.
Definition: jurassic.c:5783
#define NG
Maximum number of emitters.
Definition: jurassic.h:358
void formod_srcfunc(const ctl_t *ctl, const tbl_t *tbl, const double t, double *src)
Compute Planck source function.
Definition: jurassic.c:3483
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:3973
#define TBLNS
Maximum number of source function temperature levels.
Definition: jurassic.h:433
double scan_ctl(int argc, char *argv[], const char *varname, int arridx, const char *defvalue, char *value)
Search control parameter file for variable entry.
Definition: jurassic.c:5117
double sza(double sec, double lon, double lat)
Calculate solar zenith angle.
Definition: jurassic.c:5186
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:2921
size_t obs2y(const ctl_t *ctl, const obs_t *obs, gsl_vector *y, int *ida, int *ira)
Compose measurement vector.
Definition: jurassic.c:4178
void y2obs(const ctl_t *ctl, const gsl_vector *y, obs_t *obs)
Decompose measurement vector.
Definition: jurassic.c:5965
void hydrostatic(const ctl_t *ctl, atm_t *atm)
Set hydrostatic equilibrium.
Definition: jurassic.c:3520
void read_shape(const char *filename, double *x, double *y, int *n)
Read shape function.
Definition: jurassic.c:4910
double ctmn2(const double nu, const double p, const double t)
Compute nitrogen continuum (absorption coefficient).
Definition: jurassic.c:2796
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:418
void geo2cart(const double z, const double lon, const double lat, double *x)
Convert geolocation to Cartesian coordinates.
Definition: jurassic.c:3500
#define NSF
Maximum number of surface layer spectral grid points.
Definition: jurassic.h:373
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:5495
#define NCL
Maximum number of cloud layer spectral grid points.
Definition: jurassic.h:348
#define NLOS
Maximum number of LOS points.
Definition: jurassic.h:403
void kernel(ctl_t *ctl, atm_t *atm, obs_t *obs, gsl_matrix *k)
Compute Jacobians.
Definition: jurassic.c:4006
#define NR
Maximum number of ray paths.
Definition: jurassic.h:368
int locate_tbl(const float *xx, const int n, const double x)
Find array index in float array.
Definition: jurassic.c:4156
#define NW
Maximum number of spectral windows.
Definition: jurassic.h:378
Atmospheric data.
Definition: jurassic.h:488
double sfz
Surface height [km].
Definition: jurassic.h:527
double sfp
Surface pressure [hPa].
Definition: jurassic.h:530
double clz
Cloud layer height [km].
Definition: jurassic.h:518
int np
Number of data points.
Definition: jurassic.h:491
double cldz
Cloud layer depth [km].
Definition: jurassic.h:521
double sft
Surface temperature [K].
Definition: jurassic.h:533
Forward model control parameters.
Definition: jurassic.h:541
int write_matrix
Write matrix file (0=no, 1=yes).
Definition: jurassic.h:661
int nw
Number of spectral windows.
Definition: jurassic.h:556
double retp_zmin
Minimum altitude for pressure retrieval [km].
Definition: jurassic.h:613
double hydz
Reference height for hydrostatic pressure profile (-999 to skip) [km].
Definition: jurassic.h:586
int ret_sfz
Retrieve surface layer height (0=no, 1=yes).
Definition: jurassic.h:646
int ctm_co2
Compute CO2 continuum (0=no, 1=yes).
Definition: jurassic.h:589
double rett_zmax
Maximum altitude for temperature retrieval [km].
Definition: jurassic.h:622
int ret_sfeps
Retrieve surface layer emissivity (0=no, 1=yes).
Definition: jurassic.h:655
int ret_sft
Retrieve surface layer temperature (0=no, 1=yes).
Definition: jurassic.h:652
int ret_clz
Retrieve cloud layer height (0=no, 1=yes).
Definition: jurassic.h:637
int ctm_n2
Compute N2 continuum (0=no, 1=yes).
Definition: jurassic.h:595
int ctm_h2o
Compute H2O continuum (0=no, 1=yes).
Definition: jurassic.h:592
int formod
Forward model (0=CGA, 1=EGA, 2=RFM).
Definition: jurassic.h:664
int ng
Number of emitters.
Definition: jurassic.h:544
int refrac
Take into account refractivity (0=no, 1=yes).
Definition: jurassic.h:601
double rett_zmin
Minimum altitude for temperature retrieval [km].
Definition: jurassic.h:619
double sfsza
Solar zenith angle at the surface [deg] (-999=auto).
Definition: jurassic.h:577
int nd
Number of radiance channels.
Definition: jurassic.h:550
int sftype
Surface treatment (0=none, 1=emissions, 2=downward, 3=solar).
Definition: jurassic.h:574
int ret_sfp
Retrieve surface layer pressure (0=no, 1=yes).
Definition: jurassic.h:649
int ncl
Number of cloud layer spectral grid points.
Definition: jurassic.h:562
int ctm_o2
Compute O2 continuum (0=no, 1=yes).
Definition: jurassic.h:598
int ret_clk
Retrieve cloud layer extinction (0=no, 1=yes).
Definition: jurassic.h:643
int nsf
Number of surface layer spectral grid points.
Definition: jurassic.h:568
double rayds
Maximum step length for raytracing [km].
Definition: jurassic.h:604
int ret_cldz
Retrieve cloud layer depth (0=no, 1=yes).
Definition: jurassic.h:640
int tblfmt
Look-up table file format (1=ASCII, 2=binary).
Definition: jurassic.h:583
double raydz
Vertical step length for raytracing [km].
Definition: jurassic.h:607
int write_bbt
Use brightness temperature instead of radiance (0=no, 1=yes).
Definition: jurassic.h:658
double retp_zmax
Maximum altitude for pressure retrieval [km].
Definition: jurassic.h:616
Line-of-sight data.
Definition: jurassic.h:678
double sft
Surface temperature [K].
Definition: jurassic.h:705
int np
Number of LOS points.
Definition: jurassic.h:681
Observation geometry and radiance data.
Definition: jurassic.h:734
int nr
Number of ray paths.
Definition: jurassic.h:737
Emissivity look-up tables.
Definition: jurassic.h:778