37 int nc_result=(cmd); \
38 if(nc_result!=NC_NOERR) \
39 ERRMSG("%s", nc_strerror(nc_result)); \
134 double err_formod[
ND];
137 double err_noise[
ND];
186 const char *longname,
206 gsl_matrix * s_a_inv,
207 gsl_vector * sig_eps_inv);
270 gsl_vector * sig_noise,
271 gsl_vector * sig_formod,
272 gsl_vector * sig_eps_inv);
288 static atm_t atm_apr, atm_clim, atm_i;
289 static obs_t obs_i, obs_meas;
297 double chisq, chisq_min, chisq_max, chisq_mean, z[
NP];
299 int channel[
ND], m, ntask = -1, rank, size;
306 MPI_Init(&argc, &argv);
307 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
308 MPI_Comm_size(MPI_COMM_WORLD, &size);
315 ERRMSG(
"Give parameters: <ctl> <filelist>");
322 const int nz = (int)
scan_ctl(argc, argv,
"NZ", -1,
"", NULL);
324 ERRMSG(
"Too many altitudes!");
325 for (
int iz = 0; iz < nz; iz++)
326 z[iz] =
scan_ctl(argc, argv,
"Z", iz,
"", NULL);
329 const int track0 = (int)
scan_ctl(argc, argv,
"TRACK_MIN", -1,
"0", NULL);
330 const int track1 = (int)
scan_ctl(argc, argv,
"TRACK_MAX", -1,
"134", NULL);
333 const int xtrack0 = (int)
scan_ctl(argc, argv,
"XTRACK_MIN", -1,
"0", NULL);
335 (int)
scan_ctl(argc, argv,
"XTRACK_MAX", -1,
"89", NULL);
338 int np0 = (int)
scan_ctl(argc, argv,
"NP_MIN", -1,
"0", NULL);
339 int np1 = (int)
scan_ctl(argc, argv,
"NP_MAX", -1,
"100", NULL);
340 np1 = GSL_MIN(np1, nz - 1);
343 const double sx =
scan_ctl(argc, argv,
"SX", -1,
"8", NULL);
344 const double sy =
scan_ctl(argc, argv,
"SY", -1,
"2", NULL);
347 const double sza_thresh =
scan_ctl(argc, argv,
"SZA", -1,
"96", NULL);
354 printf(
"Read filelist: %s\n", argv[2]);
355 if (!(in = fopen(argv[2],
"r")))
356 ERRMSG(
"Cannot open filelist!");
359 while (fscanf(in,
"%s", filename) != EOF) {
362 if ((++ntask) % size != rank)
366 printf(
"Retrieve file %s on rank %d of %d (with %d threads)...\n",
367 filename, rank + 1, size, omp_get_max_threads());
377 for (
int id = 0;
id < ctl.
nd;
id++) {
380 if (fabs(ctl.
nu[
id] - ncd.
l1_nu[i]) < 0.1)
383 ERRMSG(
"Cannot identify radiance channel!");
392 for (
int iz = 0; iz < nz; iz++)
393 atm_clim.
z[iz] = z[iz];
407 for (
int track = track0; track <= track1; track++) {
410 TIMER(
"retrieval", 1);
413 for (
int xtrack = xtrack0; xtrack <= xtrack1; xtrack++) {
423 for (
int id = 0;
id < ctl.
nd;
id++)
424 obs_meas.
rad[
id][0] = ncd.
l1_rad[track][xtrack][channel[
id]];
429 for (
int id = 0;
id < ctl.
nd;
id++)
430 if (ctl.
nu[
id] >= 2000)
431 obs_meas.
rad[id][0] = GSL_NAN;
434 copy_atm(&ctl, &atm_apr, &atm_clim, 0);
435 for (
int ip = 0; ip < atm_apr.
np; ip++) {
436 atm_apr.
time[ip] = obs_meas.
time[0];
437 atm_apr.
lon[ip] = obs_meas.
vplon[0];
438 atm_apr.
lat[ip] = obs_meas.
vplat[0];
442 init_l2(&ncd, track, xtrack, &ctl, &atm_apr);
446 &atm_apr, &atm_i, &chisq);
449 if (gsl_finite(chisq)) {
450 chisq_min = GSL_MIN(chisq_min, chisq);
451 chisq_max = GSL_MAX(chisq_max, chisq);
457 buffer_nc(&atm_i, chisq, &ncd, track, xtrack, np0, np1);
461 TIMER(
"retrieval", 3);
472 printf(
"chi^2: min= %g / mean= %g / max= %g / m= %d\n",
473 chisq_min, chisq_mean / m, chisq_max, m);
474 printf(
"Retrieval finished on rank %d of %d!\n", rank, size);
484 printf(
"MEMORY_ATM = %g MByte\n", 4. *
sizeof(
atm_t) / 1024. / 1024.);
485 printf(
"MEMORY_CTL = %g MByte\n", 1. *
sizeof(
ctl_t) / 1024. / 1024.);
486 printf(
"MEMORY_NCD = %g MByte\n", 1. *
sizeof(
ncd_t) / 1024. / 1024.);
487 printf(
"MEMORY_OBS = %g MByte\n", 3. *
sizeof(
atm_t) / 1024. / 1024.);
488 printf(
"MEMORY_RET = %g MByte\n", 1. *
sizeof(
ret_t) / 1024. / 1024.);
489 printf(
"MEMORY_TBL = %g MByte\n", 1. *
sizeof(
tbl_t) / 1024. / 1024.);
492 printf(
"SIZE_TASKS = %d\n", size);
493 printf(
"SIZE_THREADS = %d\n", omp_get_max_threads());
507 const char *longname,
514 if (nc_inq_varid(ncid, varname, varid) != NC_NOERR) {
517 NC(nc_def_var(ncid, varname, type, ndims, dimid, varid));
521 (ncid, *varid,
"long_name", strlen(longname), longname));
524 NC(nc_put_att_text(ncid, *varid,
"units", strlen(unit), unit));
540 ncd->
np = np1 - np0 + 1;
543 for (
int ip = np0; ip <= np1; ip++) {
544 ncd->
ret_z[ip - np0] = (float) atm->
z[ip];
547 (gsl_finite(chisq) ? (float) atm->
t[ip] : GSL_NAN);
557 gsl_vector *sig_eps_inv) {
559 double chisq_a, chisq_m = 0;
562 const size_t m = dy->size;
563 const size_t n = dx->size;
566 gsl_vector *x_aux = gsl_vector_alloc(n);
567 gsl_vector *y_aux = gsl_vector_alloc(m);
571 for (
size_t i = 0; i < m; i++)
573 gsl_pow_2(gsl_vector_get(dy, i) * gsl_vector_get(sig_eps_inv, i));
574 gsl_blas_dgemv(CblasNoTrans, 1.0, s_a_inv, dx, 0.0, x_aux);
575 gsl_blas_ddot(dx, x_aux, &chisq_a);
578 gsl_vector_free(x_aux);
579 gsl_vector_free(y_aux);
582 return (chisq_m + chisq_a) / (double) m;
595 for (
int lay = 0; lay <
L2_NLAY; lay++) {
598 for (
int track = 0; track <
L2_NTRACK; track++)
599 for (
int xtrack = 0; xtrack <
L2_NXTRACK; xtrack++) {
602 help[track][xtrack] = 0;
606 for (
int track2 = 0; track2 <
L2_NTRACK; track2++)
607 for (
int xtrack2 = 0; xtrack2 <
L2_NXTRACK; xtrack2++)
608 if (gsl_finite(x[track2][xtrack2][lay])
609 && x[track2][xtrack2][lay] > 0) {
610 const double w = exp(-gsl_pow_2((xtrack - xtrack2) / cx)
611 - gsl_pow_2((track - track2) / cy));
612 help[track][xtrack] += w * x[track2][xtrack2][lay];
618 help[track][xtrack] /= wsum;
620 help[track][xtrack] = GSL_NAN;
624 for (
int track = 0; track <
L2_NTRACK; track++)
625 for (
int xtrack = 0; xtrack <
L2_NXTRACK; xtrack++)
626 x[track][xtrack][lay] = help[track][xtrack];
639 static atm_t atm_airs;
641 double k[
NW], p, q[
NG], t, zmax = 0, zmin = 1000;
649 for (
int lay = 0; lay <
L2_NLAY; lay++)
650 if (gsl_finite(ncd->
l2_z[track][xtrack][lay])) {
651 atm_airs.
z[atm_airs.
np] = ncd->
l2_z[track][xtrack][lay];
652 atm_airs.
p[atm_airs.
np] = ncd->
l2_p[lay];
653 atm_airs.
t[atm_airs.
np] = ncd->
l2_t[track][xtrack][lay];
654 if ((++atm_airs.
np) >
NP)
655 ERRMSG(
"Too many layers!");
659 if (atm_airs.
np <= 0)
663 for (
int ip = 0; ip < atm_airs.
np; ip++) {
664 zmax = GSL_MAX(zmax, atm_airs.
z[ip]);
665 zmin = GSL_MIN(zmin, atm_airs.
z[ip]);
669 for (
int ip = 0; ip < atm->
np; ip++) {
672 intpol_atm(ctl, &atm_airs, atm->
z[ip], &p, &t, q, k);
676 if (atm->
z[ip] > zmax)
677 w = GSL_MAX(1 - (atm->
z[ip] - zmax) / 50, 0);
678 if (atm->
z[ip] < zmin)
679 w = GSL_MAX(1 - (zmin - atm->
z[ip]) / 50, 0);
682 atm->
t[ip] = w * t + (1 - w) * atm->
t[ip];
683 atm->
p[ip] = w * p + (1 - w) * atm->
p[ip];
695 const size_t n = a->size1;
698 for (
size_t i = 0; i < n && diag; i++)
699 for (
size_t j = i + 1; j < n; j++)
700 if (gsl_matrix_get(a, i, j) != 0) {
707 for (
size_t i = 0; i < n; i++)
708 gsl_matrix_set(a, i, i, 1 / gsl_matrix_get(a, i, i));
712 gsl_linalg_cholesky_decomp(a);
713 gsl_linalg_cholesky_invert(a);
726 const size_t m = a->size1;
727 const size_t n = a->size2;
730 gsl_matrix *aux = gsl_matrix_alloc(m, n);
733 if (transpose == 1) {
736 for (
size_t i = 0; i < m; i++)
737 for (
size_t j = 0; j < n; j++)
738 gsl_matrix_set(aux, i, j,
739 gsl_vector_get(b, i) * gsl_matrix_get(a, i, j));
742 gsl_blas_dgemm(CblasTrans, CblasNoTrans, 1.0, aux, aux, 0.0, c);
746 else if (transpose == 2) {
749 for (
size_t i = 0; i < m; i++)
750 for (
size_t j = 0; j < n; j++)
751 gsl_matrix_set(aux, i, j,
752 gsl_matrix_get(a, i, j) * gsl_vector_get(b, j));
755 gsl_blas_dgemm(CblasNoTrans, CblasTrans, 1.0, aux, aux, 0.0, c);
759 gsl_matrix_free(aux);
773 static int ipa[
N], iqa[
N];
775 double disq = 0, lmpar = 0.001;
782 const size_t m =
obs2y(ctl, obs_meas, NULL, NULL, NULL);
783 const size_t n =
atm2x(ctl, atm_apr, NULL, iqa, ipa);
784 if (m == 0 || n == 0) {
790 gsl_matrix *a = gsl_matrix_alloc(n, n);
791 gsl_matrix *cov = gsl_matrix_alloc(n, n);
792 gsl_matrix *k_i = gsl_matrix_alloc(m, n);
793 gsl_matrix *s_a_inv = gsl_matrix_alloc(n, n);
795 gsl_vector *b = gsl_vector_alloc(n);
796 gsl_vector *dx = gsl_vector_alloc(n);
797 gsl_vector *dy = gsl_vector_alloc(m);
798 gsl_vector *sig_eps_inv = gsl_vector_alloc(m);
799 gsl_vector *sig_formod = gsl_vector_alloc(m);
800 gsl_vector *sig_noise = gsl_vector_alloc(m);
801 gsl_vector *x_a = gsl_vector_alloc(n);
802 gsl_vector *x_i = gsl_vector_alloc(n);
803 gsl_vector *x_step = gsl_vector_alloc(n);
804 gsl_vector *y_aux = gsl_vector_alloc(m);
805 gsl_vector *y_i = gsl_vector_alloc(m);
806 gsl_vector *y_m = gsl_vector_alloc(m);
811 formod(ctl, atm_i, obs_i);
814 atm2x(ctl, atm_apr, x_a, NULL, NULL);
815 atm2x(ctl, atm_i, x_i, NULL, NULL);
816 obs2y(ctl, obs_meas, y_m, NULL, NULL);
817 obs2y(ctl, obs_i, y_i, NULL, NULL);
824 set_cov_meas(ret, ctl, obs_meas, sig_noise, sig_formod, sig_eps_inv);
827 gsl_vector_memcpy(dx, x_i);
828 gsl_vector_sub(dx, x_a);
829 gsl_vector_memcpy(dy, y_m);
830 gsl_vector_sub(dy, y_i);
836 kernel(ctl, atm_i, obs_i, k_i);
843 for (
int it = 1; it <= ret->
conv_itmax; it++) {
846 double chisq_old = *chisq;
850 kernel(ctl, atm_i, obs_i, k_i);
857 for (
size_t i = 0; i < m; i++)
858 gsl_vector_set(y_aux, i, gsl_vector_get(dy, i)
859 * gsl_pow_2(gsl_vector_get(sig_eps_inv, i)));
860 gsl_blas_dgemv(CblasTrans, 1.0, k_i, y_aux, 0.0, b);
861 gsl_blas_dgemv(CblasNoTrans, -1.0, s_a_inv, dx, 1.0, b);
864 for (
int it2 = 0; it2 < 20; it2++) {
867 gsl_matrix_memcpy(a, s_a_inv);
868 gsl_matrix_scale(a, 1 + lmpar);
869 gsl_matrix_add(a, cov);
872 gsl_linalg_cholesky_decomp(a);
873 gsl_linalg_cholesky_solve(a, b, x_step);
876 gsl_vector_add(x_i, x_step);
879 x2atm(ctl, x_i, atm_i);
882 for (
int ip = 0; ip < atm_i->
np; ip++) {
883 atm_i->
p[ip] = GSL_MIN(GSL_MAX(atm_i->
p[ip], 5e-7), 5e4);
884 atm_i->
t[ip] = GSL_MIN(GSL_MAX(atm_i->
t[ip], 100), 400);
885 for (
int ig = 0; ig < ctl->
ng; ig++)
886 atm_i->
q[ig][ip] = GSL_MIN(GSL_MAX(atm_i->
q[ig][ip], 0), 1);
887 for (
int iw = 0; iw < ctl->
nw; iw++)
888 atm_i->
k[iw][ip] = GSL_MAX(atm_i->
k[iw][ip], 0);
892 formod(ctl, atm_i, obs_i);
893 obs2y(ctl, obs_i, y_i, NULL, NULL);
896 gsl_vector_memcpy(dx, x_i);
897 gsl_vector_sub(dx, x_a);
898 gsl_vector_memcpy(dy, y_m);
899 gsl_vector_sub(dy, y_i);
905 if (*chisq > chisq_old) {
907 gsl_vector_sub(x_i, x_step);
915 gsl_blas_ddot(x_step, b, &disq);
928 gsl_matrix_free(cov);
929 gsl_matrix_free(k_i);
930 gsl_matrix_free(s_a_inv);
935 gsl_vector_free(sig_eps_inv);
936 gsl_vector_free(sig_formod);
937 gsl_vector_free(sig_noise);
938 gsl_vector_free(x_a);
939 gsl_vector_free(x_i);
940 gsl_vector_free(x_step);
941 gsl_vector_free(y_aux);
942 gsl_vector_free(y_i);
943 gsl_vector_free(y_m);
955 printf(
"Read netCDF file: %s\n", filename);
956 NC(nc_open(filename, NC_WRITE, &ncd->
ncid));
959 NC(nc_inq_varid(ncd->
ncid,
"l1_time", &varid));
961 NC(nc_inq_varid(ncd->
ncid,
"l1_lon", &varid));
962 NC(nc_get_var_double(ncd->
ncid, varid, ncd->
l1_lon[0]));
963 NC(nc_inq_varid(ncd->
ncid,
"l1_lat", &varid));
964 NC(nc_get_var_double(ncd->
ncid, varid, ncd->
l1_lat[0]));
965 NC(nc_inq_varid(ncd->
ncid,
"l1_sat_z", &varid));
967 NC(nc_inq_varid(ncd->
ncid,
"l1_sat_lon", &varid));
969 NC(nc_inq_varid(ncd->
ncid,
"l1_sat_lat", &varid));
971 NC(nc_inq_varid(ncd->
ncid,
"l1_nu", &varid));
972 NC(nc_get_var_double(ncd->
ncid, varid, ncd->
l1_nu));
973 NC(nc_inq_varid(ncd->
ncid,
"l1_rad", &varid));
974 NC(nc_get_var_float(ncd->
ncid, varid, ncd->
l1_rad[0][0]));
977 NC(nc_inq_varid(ncd->
ncid,
"l2_z", &varid));
978 NC(nc_get_var_double(ncd->
ncid, varid, ncd->
l2_z[0][0]));
979 NC(nc_inq_varid(ncd->
ncid,
"l2_press", &varid));
980 NC(nc_get_var_double(ncd->
ncid, varid, ncd->
l2_p));
981 NC(nc_inq_varid(ncd->
ncid,
"l2_temp", &varid));
982 NC(nc_get_var_double(ncd->
ncid, varid, ncd->
l2_t[0][0]));
995 (int)
scan_ctl(argc, argv,
"KERNEL_RECOMP", -1,
"3", NULL);
999 for (
int id = 0;
id < ctl->
nd;
id++)
1002 for (
int id = 0;
id < ctl->
nd;
id++)
1013 for (
int ig = 0; ig < ctl->
ng; ig++) {
1014 ret->
err_q[ig] =
scan_ctl(argc, argv,
"ERR_Q", ig,
"0", NULL);
1019 for (
int iw = 0; iw < ctl->
nw; iw++) {
1020 ret->
err_k[iw] =
scan_ctl(argc, argv,
"ERR_K", iw,
"0", NULL);
1037 const size_t n = s_a->size1;
1040 gsl_vector *x_a = gsl_vector_alloc(n);
1043 atm2x(ctl, atm, x_a, NULL, NULL);
1044 for (
size_t i = 0; i < n; i++) {
1046 gsl_vector_set(x_a, i, ret->
err_press / 100 * gsl_vector_get(x_a, i));
1048 gsl_vector_set(x_a, i, ret->
err_temp);
1049 for (
int ig = 0; ig < ctl->
ng; ig++)
1050 if (iqa[i] ==
IDXQ(ig))
1051 gsl_vector_set(x_a, i, ret->
err_q[ig] / 100 * gsl_vector_get(x_a, i));
1052 for (
int iw = 0; iw < ctl->
nw; iw++)
1053 if (iqa[i] ==
IDXK(iw))
1054 gsl_vector_set(x_a, i, ret->
err_k[iw]);
1058 for (
size_t i = 0; i < n; i++)
1059 if (gsl_pow_2(gsl_vector_get(x_a, i)) <= 0)
1060 ERRMSG(
"Check a priori data (zero standard deviation)!");
1063 gsl_matrix_set_zero(s_a);
1064 for (
size_t i = 0; i < n; i++)
1065 gsl_matrix_set(s_a, i, i, gsl_pow_2(gsl_vector_get(x_a, i)));
1068 for (
size_t i = 0; i < n; i++)
1069 for (
size_t j = 0; j < n; j++)
1070 if (i != j && iqa[i] == iqa[j]) {
1077 if (iqa[i] ==
IDXP) {
1083 if (iqa[i] ==
IDXT) {
1089 for (
int ig = 0; ig < ctl->
ng; ig++)
1090 if (iqa[i] ==
IDXQ(ig)) {
1096 for (
int iw = 0; iw < ctl->
nw; iw++)
1097 if (iqa[i] ==
IDXK(iw)) {
1103 if (cz > 0 && ch > 0) {
1106 double x0[3], x1[3];
1112 exp(-
DIST(x0, x1) / ch -
1113 fabs(atm->
z[ipa[i]] - atm->
z[ipa[j]]) / cz);
1116 gsl_matrix_set(s_a, i, j, gsl_vector_get(x_a, i)
1117 * gsl_vector_get(x_a, j) * rho);
1122 gsl_vector_free(x_a);
1131 gsl_vector *sig_noise,
1132 gsl_vector *sig_formod,
1133 gsl_vector *sig_eps_inv) {
1135 static obs_t obs_err;
1138 const size_t m = sig_eps_inv->size;
1142 for (
int ir = 0; ir < obs_err.
nr; ir++)
1143 for (
int id = 0;
id < ctl->
nd;
id++)
1145 = (gsl_finite(obs->
rad[
id][ir]) ? ret->
err_noise[
id] : GSL_NAN);
1146 obs2y(ctl, &obs_err, sig_noise, NULL, NULL);
1150 for (
int ir = 0; ir < obs_err.
nr; ir++)
1151 for (
int id = 0;
id < ctl->
nd;
id++)
1154 obs2y(ctl, &obs_err, sig_formod, NULL, NULL);
1157 for (
size_t i = 0; i < m; i++)
1158 gsl_vector_set(sig_eps_inv, i,
1159 1 / sqrt(gsl_pow_2(gsl_vector_get(sig_noise, i))
1160 + gsl_pow_2(gsl_vector_get(sig_formod, i))));
1163 for (
size_t i = 0; i < m; i++)
1164 if (gsl_vector_get(sig_eps_inv, i) <= 0)
1165 ERRMSG(
"Check measurement errors (zero standard deviation)!");
1174 int dimid[10], p_id, t_id, z_id;
1177 printf(
"Write netCDF file: %s\n", filename);
1180 NC(nc_inq_dimid(ncd->
ncid,
"L1_NTRACK", &dimid[0]));
1181 NC(nc_inq_dimid(ncd->
ncid,
"L1_NXTRACK", &dimid[1]));
1187 if (nc_inq_dimid(ncd->
ncid,
"RET_NP", &dimid[2]) != NC_NOERR)
1188 NC(nc_def_dim(ncd->
ncid,
"RET_NP", (
size_t) ncd->
np, &dimid[2]));
1191 add_var(ncd->
ncid,
"ret_z",
"km",
"altitude", NC_FLOAT, &dimid[2], &z_id,
1193 add_var(ncd->
ncid,
"ret_press",
"hPa",
"pressure", NC_FLOAT, dimid, &p_id,
1195 add_var(ncd->
ncid,
"ret_temp",
"K",
"temperature", NC_FLOAT, dimid, &t_id,
1199 NC(nc_enddef(ncd->
ncid));
1202 NC(nc_put_var_float(ncd->
ncid, z_id, ncd->
ret_z));
1203 NC(nc_put_var_float(ncd->
ncid, p_id, ncd->
ret_p));
1204 NC(nc_put_var_float(ncd->
ncid, t_id, ncd->
ret_t));
void read_ctl(int argc, char *argv[], ctl_t *ctl)
Read forward model control parameters.
void x2atm(const ctl_t *ctl, const gsl_vector *x, atm_t *atm)
Decompose parameter vector or state vector.
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.
double sza(const double sec, const double lon, const double lat)
Calculate solar zenith angle.
void formod(const ctl_t *ctl, atm_t *atm, obs_t *obs)
Determine ray paths and compute radiative transfer.
void copy_obs(const ctl_t *ctl, obs_t *obs_dest, const obs_t *obs_src, const int init)
Copy and initialize observation data.
double scan_ctl(int argc, char *argv[], const char *varname, int arridx, const char *defvalue, char *value)
Search control parameter file for variable entry.
void copy_atm(const ctl_t *ctl, atm_t *atm_dest, const atm_t *atm_src, const int init)
Copy and initialize atmospheric data.
size_t obs2y(const ctl_t *ctl, const obs_t *obs, gsl_vector *y, int *ida, int *ira)
Compose measurement vector.
size_t atm2x(const ctl_t *ctl, const atm_t *atm, gsl_vector *x, int *iqa, int *ipa)
Compose state vector or parameter vector.
void climatology(const ctl_t *ctl, atm_t *atm)
Interpolate climatological data.
void geo2cart(const double z, const double lon, const double lat, double *x)
Convert geolocation to Cartesian coordinates.
void kernel(ctl_t *ctl, atm_t *atm, obs_t *obs, gsl_matrix *k)
Compute Jacobians.
JURASSIC library declarations.
#define N
Maximum size of state vector.
#define LEN
Maximum length of ASCII data lines.
#define ERRMSG(...)
Print error message and quit program.
#define IDXK(iw)
Indices for extinction.
#define ND
Maximum number of radiance channels.
#define TIMER(name, mode)
Start or stop a timer.
#define IDXP
Index for pressure.
#define NP
Maximum number of atmospheric data points.
#define NG
Maximum number of emitters.
#define DIST(a, b)
Compute Cartesian distance between two vectors.
#define IDXQ(ig)
Indices for volume mixing ratios.
#define NW
Maximum number of spectral windows.
#define IDXT
Index for temperature.
void matrix_product(gsl_matrix *a, gsl_vector *b, int transpose, gsl_matrix *c)
Compute matrix product A^TBA or ABA^T for diagonal matrix B.
int main(int argc, char *argv[])
#define NC(cmd)
Execute netCDF library command and check result.
void set_cov_apr(ret_t *ret, ctl_t *ctl, atm_t *atm, int *iqa, int *ipa, gsl_matrix *s_a)
Set a priori covariance.
void optimal_estimation(ret_t *ret, ctl_t *ctl, obs_t *obs_meas, obs_t *obs_i, atm_t *atm_apr, atm_t *atm_i, double *chisq)
Carry out optimal estimation retrieval.
void fill_gaps(double x[L2_NTRACK][L2_NXTRACK][L2_NLAY], double cx, double cy)
Fill data gaps in L2 data.
#define L1_NXTRACK
Across-track size of AIRS radiance granule (don't change).
#define L1_NTRACK
Along-track size of AIRS radiance granule (don't change).
#define L2_NXTRACK
Across-track size of AIRS retrieval granule (don't change).
void write_nc(char *filename, ncd_t *ncd)
Write to netCDF file...
void read_nc(char *filename, ncd_t *ncd)
Read netCDF file.
void add_var(int ncid, const char *varname, const char *unit, const char *longname, int type, int dimid[], int *varid, int ndims)
Create variable in netCDF file.
#define L2_NLAY
Number of AIRS pressure layers (don't change).
void buffer_nc(atm_t *atm, double chisq, ncd_t *ncd, int track, int xtrack, int np0, int np1)
Buffer netCDF data.
void set_cov_meas(ret_t *ret, ctl_t *ctl, obs_t *obs, gsl_vector *sig_noise, gsl_vector *sig_formod, gsl_vector *sig_eps_inv)
Set measurement errors.
void read_ret_ctl(int argc, char *argv[], ctl_t *ctl, ret_t *ret)
Read retrieval control parameters.
void matrix_invert(gsl_matrix *a)
Invert symmetric matrix.
#define L1_NCHAN
Number of AIRS radiance channels (don't change).
#define L2_NTRACK
Along-track size of AIRS retrieval granule (don't change).
double cost_function(gsl_vector *dx, gsl_vector *dy, gsl_matrix *s_a_inv, gsl_vector *sig_eps_inv)
Compute cost function.
void init_l2(ncd_t *ncd, int track, int xtrack, ctl_t *ctl, atm_t *atm)
Initialize with AIRS Level-2 data.
double time[NP]
Time (seconds since 2000-01-01T00:00Z).
double k[NW][NP]
Extinction [km^-1].
double lat[NP]
Latitude [deg].
double lon[NP]
Longitude [deg].
double t[NP]
Temperature [K].
int np
Number of data points.
double z[NP]
Altitude [km].
double q[NG][NP]
Volume mixing ratio [ppv].
double p[NP]
Pressure [hPa].
Forward model control parameters.
int nw
Number of spectral windows.
double nu[ND]
Centroid wavenumber of each channel [cm^-1].
int ng
Number of emitters.
int nd
Number of radiance channels.
double l1_lat[L1_NTRACK][L1_NXTRACK]
Footprint latitude [deg].
float ret_p[L1_NTRACK *L1_NXTRACK]
Pressure [hPa].
double l1_lon[L1_NTRACK][L1_NXTRACK]
Footprint longitude [deg].
double l2_z[L2_NTRACK][L2_NXTRACK][L2_NLAY]
Altitude [km].
double l1_sat_lat[L1_NTRACK]
Satellite latitude [deg].
int np
Number of retrieval altitudes.
double l1_sat_lon[L1_NTRACK]
Satellite longitude [deg].
float ret_z[NP]
Altitude [km].
double l2_p[L2_NLAY]
Pressure [hPa].
float ret_t[L1_NTRACK *L1_NXTRACK *NP]
Temperature [K].
double l2_t[L2_NTRACK][L2_NXTRACK][L2_NLAY]
Temperature [K].
double l1_nu[L1_NCHAN]
Channel frequencies [cm^-1].
double l1_sat_z[L1_NTRACK]
Satellite altitude [km].
float l1_rad[L1_NTRACK][L1_NXTRACK][L1_NCHAN]
Radiance [W/(m^2 sr cm^-1)].
double l1_time[L1_NTRACK][L1_NXTRACK]
Time (seconds since 2000-01-01T00:00Z).
Observation geometry and radiance data.
double rad[ND][NR]
Radiance [W/(m^2 sr cm^-1)].
double vplat[NR]
View point latitude [deg].
double obslon[NR]
Observer longitude [deg].
double obslat[NR]
Observer latitude [deg].
double obsz[NR]
Observer altitude [km].
double vplon[NR]
View point longitude [deg].
double time[NR]
Time (seconds since 2000-01-01T00:00Z).
int nr
Number of ray paths.
double err_k_ch[NW]
Horizontal correlation length for extinction error [km].
double err_press_cz
Vertical correlation length for pressure error [km].
double err_k_cz[NW]
Vertical correlation length for extinction error [km].
double err_press
Pressure error [%].
double err_q_ch[NG]
Horizontal correlation length for volume mixing ratio error [km].
double err_noise[ND]
Noise error [W/(m^2 sr cm^-1)].
double err_formod[ND]
Forward model error [%].
double err_q_cz[NG]
Vertical correlation length for volume mixing ratio error [km].
double err_temp_cz
Vertical correlation length for temperature error [km].
double conv_dmin
Minimum normalized step size in state space.
double err_temp
Temperature error [K].
double err_temp_ch
Horizontal correlation length for temperature error [km].
int kernel_recomp
Recomputation of kernel matrix (number of iterations).
int conv_itmax
Maximum number of iterations.
double err_press_ch
Horizontal correlation length for pressure error [km].
double err_q[NG]
Volume mixing ratio error [%].
double err_k[NW]
Extinction error [1/km].
Emissivity look-up tables.