50 double err_formod[
ND];
141 gsl_matrix * s_a_inv,
142 gsl_vector * sig_eps_inv);
185 gsl_vector * sig_noise,
186 gsl_vector * sig_formod,
187 gsl_vector * sig_eps_inv);
191 const char *quantity,
205 static atm_t atm_i, atm_apr;
207 static obs_t obs_i, obs_meas;
214 ERRMSG(
"Give parameters: <ctl> <dirlist>");
224 if (!(dirlist = fopen(argv[2],
"r")))
225 ERRMSG(
"Cannot open directory list!");
228 while (fscanf(dirlist,
"%s", ret.
dir) != EOF) {
231 LOG(1,
"\nRetrieve in directory %s...\n", ret.
dir);
237 read_obs(ret.
dir,
"obs_meas.tab", &ctl, &obs_meas);
247 LOG(1,
"\nRetrieval done...");
265 static atm_t atm_cont, atm_res;
267 int icl, ig, iq, isf, iw;
269 size_t i, n, n0[
NQ], n1[
NQ];
275 for (iq = 0; iq <
NQ; iq++) {
277 for (i = 0; i < n; i++) {
278 if (iqa[i] == iq && n0[iq] ==
N)
281 n1[iq] = i - n0[iq] + 1;
292 for (ig = 0; ig < ctl->
ng; ig++)
294 atm_cont.
q[ig], atm_res.
q[ig]);
295 for (iw = 0; iw < ctl->
nw; iw++)
297 atm_cont.
k[iw], atm_res.
k[iw]);
301 for (icl = 0; icl < ctl->
ncl; icl++)
303 &atm_cont.
clk[icl], &atm_res.
clk[icl]);
307 for (isf = 0; isf < ctl->
nsf; isf++)
329 for (
size_t i = 0; i < n1[iq]; i++) {
332 for (
size_t j = 0; j < n1[iq]; j++)
333 cont[ipa[n0[iq] + i]] += gsl_matrix_get(avk, n0[iq] + i, n0[iq] + j);
336 res[ipa[n0[iq] + i]] = 1 / gsl_matrix_get(avk, n0[iq] + i, n0[iq] + i);
345 gsl_matrix * s_a_inv,
346 gsl_vector * sig_eps_inv) {
348 gsl_vector *x_aux, *y_aux;
350 double chisq_a, chisq_m = 0;
357 x_aux = gsl_vector_alloc(n);
358 y_aux = gsl_vector_alloc(m);
362 for (
size_t i = 0; i < m; i++)
363 chisq_m +=
POW2(gsl_vector_get(dy, i) * gsl_vector_get(sig_eps_inv, i));
364 gsl_blas_dgemv(CblasNoTrans, 1.0, s_a_inv, dx, 0.0, x_aux);
365 gsl_blas_ddot(dx, x_aux, &chisq_a);
368 gsl_vector_free(x_aux);
369 gsl_vector_free(y_aux);
372 return (chisq_m + chisq_a) / (double) m;
386 for (
size_t i = 0; i < n && diag; i++)
387 for (
size_t j = i + 1; j < n; j++)
388 if (gsl_matrix_get(a, i, j) != 0) {
395 for (
size_t i = 0; i < n; i++)
396 gsl_matrix_set(a, i, i, 1 / gsl_matrix_get(a, i, i));
400 gsl_linalg_cholesky_decomp(a);
401 gsl_linalg_cholesky_invert(a);
420 aux = gsl_matrix_alloc(m, n);
423 if (transpose == 1) {
426 for (
size_t i = 0; i < m; i++)
427 for (
size_t j = 0; j < n; j++)
428 gsl_matrix_set(aux, i, j,
429 gsl_vector_get(b, i) * gsl_matrix_get(a, i, j));
432 gsl_blas_dgemm(CblasTrans, CblasNoTrans, 1.0, aux, aux, 0.0, c);
436 else if (transpose == 2) {
439 for (
size_t i = 0; i < m; i++)
440 for (
size_t j = 0; j < n; j++)
441 gsl_matrix_set(aux, i, j,
442 gsl_matrix_get(a, i, j) * gsl_vector_get(b, j));
445 gsl_blas_dgemm(CblasNoTrans, CblasTrans, 1.0, aux, aux, 0.0, c);
449 gsl_matrix_free(aux);
462 static int ipa[
N], iqa[
N];
464 gsl_matrix *a, *auxnm, *corr, *cov, *gain, *k_i, *s_a_inv;
465 gsl_vector *b, *dx, *dy, *sig_eps_inv, *sig_formod, *sig_noise,
466 *x_a, *x_i, *x_step, *y_aux, *y_i, *y_m;
470 char filename[2 *
LEN];
472 double chisq, disq = 0, lmpar = 0.001;
474 int icl, ig, ip, isf, it = 0, it2, iw;
483 m =
obs2y(ctl, obs_meas, NULL, NULL, NULL);
484 n =
atm2x(ctl, atm_apr, NULL, iqa, ipa);
485 if (m == 0 || n == 0)
486 ERRMSG(
"Check problem definition!");
489 LOG(1,
"Problem size: m= %d / n= %d "
490 "(alloc= %.4g MB / stat= %.4g MB)",
492 (
double) (3 * m * n + 4 * n * n + 8 * m +
493 8 * n) *
sizeof(
double) / 1024. / 1024.,
494 (
double) (5 *
sizeof(
atm_t) + 3 *
sizeof(
obs_t)
495 + 2 *
N *
sizeof(
int)) / 1024. / 1024.);
498 a = gsl_matrix_alloc(n, n);
499 cov = gsl_matrix_alloc(n, n);
500 k_i = gsl_matrix_alloc(m, n);
501 s_a_inv = gsl_matrix_alloc(n, n);
503 b = gsl_vector_alloc(n);
504 dx = gsl_vector_alloc(n);
505 dy = gsl_vector_alloc(m);
506 sig_eps_inv = gsl_vector_alloc(m);
507 sig_formod = gsl_vector_alloc(m);
508 sig_noise = gsl_vector_alloc(m);
509 x_a = gsl_vector_alloc(n);
510 x_i = gsl_vector_alloc(n);
511 x_step = gsl_vector_alloc(n);
512 y_aux = gsl_vector_alloc(m);
513 y_i = gsl_vector_alloc(m);
514 y_m = gsl_vector_alloc(m);
519 formod(ctl, atm_i, obs_i);
522 atm2x(ctl, atm_apr, x_a, NULL, NULL);
523 atm2x(ctl, atm_i, x_i, NULL, NULL);
524 obs2y(ctl, obs_meas, y_m, NULL, NULL);
525 obs2y(ctl, obs_i, y_i, NULL, NULL);
530 atm_i, obs_i,
"x",
"x",
"r");
534 set_cov_meas(ret, ctl, obs_meas, sig_noise, sig_formod, sig_eps_inv);
537 sprintf(filename,
"%s/costs.tab", ret->
dir);
538 if (!(out = fopen(filename,
"w")))
539 ERRMSG(
"Cannot create cost function file!");
543 "# $1 = iteration number\n"
544 "# $2 = normalized cost function\n"
545 "# $3 = number of measurements\n"
546 "# $4 = number of state vector elements\n\n");
549 gsl_vector_memcpy(dx, x_i);
550 gsl_vector_sub(dx, x_a);
551 gsl_vector_memcpy(dy, y_m);
552 gsl_vector_sub(dy, y_i);
558 LOG(1,
"it= %d / chi^2/m= %g", it, chisq);
561 fprintf(out,
"%d %g %d %d\n", it, chisq, (
int) m, (
int) n);
564 kernel(ctl, atm_i, obs_i, k_i);
574 double chisq_old = chisq;
578 kernel(ctl, atm_i, obs_i, k_i);
585 for (i = 0; i < m; i++)
586 gsl_vector_set(y_aux, i, gsl_vector_get(dy, i)
587 *
POW2(gsl_vector_get(sig_eps_inv, i)));
588 gsl_blas_dgemv(CblasTrans, 1.0, k_i, y_aux, 0.0, b);
589 gsl_blas_dgemv(CblasNoTrans, -1.0, s_a_inv, dx, 1.0, b);
592 for (it2 = 0; it2 < 20; it2++) {
595 gsl_matrix_memcpy(a, s_a_inv);
596 gsl_matrix_scale(a, 1 + lmpar);
597 gsl_matrix_add(a, cov);
600 gsl_linalg_cholesky_decomp(a);
601 gsl_linalg_cholesky_solve(a, b, x_step);
604 gsl_vector_add(x_i, x_step);
607 x2atm(ctl, x_i, atm_i);
610 for (ip = 0; ip < atm_i->
np; ip++) {
611 atm_i->
p[ip] =
MIN(
MAX(atm_i->
p[ip], 5e-7), 5e4);
612 atm_i->
t[ip] =
MIN(
MAX(atm_i->
t[ip], 100), 400);
613 for (ig = 0; ig < ctl->
ng; ig++)
614 atm_i->
q[ig][ip] =
MIN(
MAX(atm_i->
q[ig][ip], 0), 1);
615 for (iw = 0; iw < ctl->
nw; iw++)
616 atm_i->
k[iw][ip] =
MAX(atm_i->
k[iw][ip], 0);
620 for (icl = 0; icl < ctl->
ncl; icl++)
621 atm_i->
clk[icl] =
MAX(atm_i->
clk[icl], 0);
625 for (isf = 0; isf < ctl->
nsf; isf++)
629 formod(ctl, atm_i, obs_i);
630 obs2y(ctl, obs_i, y_i, NULL, NULL);
633 gsl_vector_memcpy(dx, x_i);
634 gsl_vector_sub(dx, x_a);
635 gsl_vector_memcpy(dy, y_m);
636 gsl_vector_sub(dy, y_i);
642 if (chisq > chisq_old) {
644 gsl_vector_sub(x_i, x_step);
652 LOG(1,
"it= %d / chi^2/m= %g", it, chisq);
655 fprintf(out,
"%d %g %d %d\n", it, chisq, (
int) m, (
int) n);
658 gsl_blas_ddot(x_step, b, &disq);
673 atm_i, obs_i,
"y",
"x",
"r");
683 auxnm = gsl_matrix_alloc(n, m);
684 corr = gsl_matrix_alloc(n, n);
685 gain = gsl_matrix_alloc(n, m);
690 gsl_matrix_add(cov, s_a_inv);
695 atm_i, obs_i,
"x",
"x",
"r");
699 for (i = 0; i < n; i++)
700 for (
size_t j = 0; j < n; j++)
701 gsl_matrix_set(corr, i, j, gsl_matrix_get(cov, i, j)
702 / sqrt(gsl_matrix_get(cov, i, i))
703 / sqrt(gsl_matrix_get(cov, j, j)));
705 atm_i, obs_i,
"x",
"x",
"r");
709 for (i = 0; i < n; i++)
710 for (
size_t j = 0; j < m; j++)
711 gsl_matrix_set(auxnm, i, j, gsl_matrix_get(k_i, j, i)
712 *
POW2(gsl_vector_get(sig_eps_inv, j)));
713 gsl_blas_dgemm(CblasNoTrans, CblasNoTrans, 1.0, cov, auxnm, 0.0, gain);
715 atm_i, obs_i,
"x",
"y",
"c");
727 gsl_blas_dgemm(CblasNoTrans, CblasNoTrans, 1.0, gain, k_i, 0.0, a);
729 atm_i, obs_i,
"x",
"x",
"r");
735 gsl_matrix_free(auxnm);
736 gsl_matrix_free(corr);
737 gsl_matrix_free(gain);
745 gsl_matrix_free(cov);
746 gsl_matrix_free(k_i);
747 gsl_matrix_free(s_a_inv);
752 gsl_vector_free(sig_eps_inv);
753 gsl_vector_free(sig_formod);
754 gsl_vector_free(sig_noise);
755 gsl_vector_free(x_a);
756 gsl_vector_free(x_i);
757 gsl_vector_free(x_step);
758 gsl_vector_free(y_aux);
759 gsl_vector_free(y_i);
760 gsl_vector_free(y_m);
773 (int)
scan_ctl(argc, argv,
"KERNEL_RECOMP", -1,
"3", NULL);
780 for (
int id = 0;
id < ctl->
nd;
id++)
783 for (
int id = 0;
id < ctl->
nd;
id++)
794 for (
int ig = 0; ig < ctl->
ng; ig++) {
795 ret->
err_q[ig] =
scan_ctl(argc, argv,
"ERR_Q", ig,
"0", NULL);
800 for (
int iw = 0; iw < ctl->
nw; iw++) {
801 ret->
err_k[iw] =
scan_ctl(argc, argv,
"ERR_K", iw,
"0", NULL);
808 for (
int icl = 0; icl < ctl->
ncl; icl++)
814 for (
int isf = 0; isf < ctl->
nsf; isf++)
833 size_t n = s_a->size1;
836 x_a = gsl_vector_alloc(n);
839 atm2x(ctl, atm, x_a, NULL, NULL);
840 for (
size_t i = 0; i < n; i++) {
842 gsl_vector_set(x_a, i, ret->
err_press / 100 * gsl_vector_get(x_a, i));
844 gsl_vector_set(x_a, i, ret->
err_temp);
845 for (
int ig = 0; ig < ctl->
ng; ig++)
846 if (iqa[i] ==
IDXQ(ig))
847 gsl_vector_set(x_a, i, ret->
err_q[ig] / 100 * gsl_vector_get(x_a, i));
848 for (
int iw = 0; iw < ctl->
nw; iw++)
849 if (iqa[i] ==
IDXK(iw))
850 gsl_vector_set(x_a, i, ret->
err_k[iw]);
852 gsl_vector_set(x_a, i, ret->
err_clz);
854 gsl_vector_set(x_a, i, ret->
err_cldz);
855 for (
int icl = 0; icl < ctl->
ncl; icl++)
856 if (iqa[i] ==
IDXCLK(icl))
857 gsl_vector_set(x_a, i, ret->
err_clk[icl]);
859 gsl_vector_set(x_a, i, ret->
err_sfz);
861 gsl_vector_set(x_a, i, ret->
err_sfp);
863 gsl_vector_set(x_a, i, ret->
err_sft);
864 for (
int isf = 0; isf < ctl->
nsf; isf++)
866 gsl_vector_set(x_a, i, ret->
err_sfeps[isf]);
870 for (
size_t i = 0; i < n; i++)
871 if (
POW2(gsl_vector_get(x_a, i)) <= 0)
872 ERRMSG(
"Check a priori data (zero standard deviation)!");
875 gsl_matrix_set_zero(s_a);
876 for (
size_t i = 0; i < n; i++)
877 gsl_matrix_set(s_a, i, i,
POW2(gsl_vector_get(x_a, i)));
880 for (
size_t i = 0; i < n; i++)
881 for (
size_t j = 0; j < n; j++)
882 if (i != j && iqa[i] == iqa[j]) {
889 if (iqa[i] ==
IDXP) {
895 if (iqa[i] ==
IDXT) {
901 for (
int ig = 0; ig < ctl->
ng; ig++)
902 if (iqa[i] ==
IDXQ(ig)) {
908 for (
int iw = 0; iw < ctl->
nw; iw++)
909 if (iqa[i] ==
IDXK(iw)) {
915 if (cz > 0 && ch > 0) {
923 exp(-
DIST(x0, x1) / ch -
924 fabs(atm->
z[ipa[i]] - atm->
z[ipa[j]]) / cz);
927 gsl_matrix_set(s_a, i, j, gsl_vector_get(x_a, i)
928 * gsl_vector_get(x_a, j) * rho);
933 gsl_vector_free(x_a);
942 gsl_vector * sig_noise,
943 gsl_vector * sig_formod,
944 gsl_vector * sig_eps_inv) {
946 static obs_t obs_err;
949 size_t m = sig_eps_inv->size;
953 for (
int ir = 0; ir < obs_err.
nr; ir++)
954 for (
int id = 0;
id < ctl->
nd;
id++)
956 = (isfinite(obs->
rad[
id][ir]) ? ret->
err_noise[
id] : NAN);
957 obs2y(ctl, &obs_err, sig_noise, NULL, NULL);
961 for (
int ir = 0; ir < obs_err.
nr; ir++)
962 for (
int id = 0;
id < ctl->
nd;
id++)
965 obs2y(ctl, &obs_err, sig_formod, NULL, NULL);
968 for (
size_t i = 0; i < m; i++)
969 gsl_vector_set(sig_eps_inv, i, 1 / sqrt(
POW2(gsl_vector_get(sig_noise, i))
975 for (
size_t i = 0; i < m; i++)
976 if (gsl_vector_get(sig_eps_inv, i) <= 0)
977 ERRMSG(
"Check measurement errors (zero standard deviation)!");
983 const char *quantity,
989 static atm_t atm_aux;
999 x_aux = gsl_vector_alloc(n);
1002 for (
size_t i = 0; i < n; i++)
1003 gsl_vector_set(x_aux, i, sqrt(gsl_matrix_get(s, i, i)));
1007 x2atm(ctl, x_aux, &atm_aux);
1008 sprintf(filename,
"atm_err_%s.tab", quantity);
1012 gsl_vector_free(x_aux);
void write_atm(const char *dirname, const char *filename, const ctl_t *ctl, const atm_t *atm)
Write atmospheric data.
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 formod(const ctl_t *ctl, atm_t *atm, obs_t *obs)
Determine ray paths and compute radiative transfer.
void read_atm(const char *dirname, const char *filename, const ctl_t *ctl, atm_t *atm)
Read atmospheric data.
void copy_obs(const ctl_t *ctl, obs_t *obs_dest, const obs_t *obs_src, const int init)
Copy and initialize observation data.
void read_obs(const char *dirname, const char *filename, const ctl_t *ctl, obs_t *obs)
Read observation data.
void write_obs(const char *dirname, const char *filename, const ctl_t *ctl, const obs_t *obs)
Write 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 geo2cart(const double z, const double lon, const double lat, double *x)
Convert geolocation to Cartesian coordinates.
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.
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 IDXCLZ
Index for cloud layer height.
#define LEN
Maximum length of ASCII data lines.
#define POW2(x)
Compute x^2.
#define MIN(a, b)
Macro to determine the minimum of two values.
#define IDXCLDZ
Index for cloud layer depth.
#define ERRMSG(...)
Print error message and quit program.
#define IDXK(iw)
Indices for extinction.
#define NQ
Maximum number of quantities.
#define ND
Maximum number of radiance channels.
#define IDXSFT
Index for surface layer temperature.
#define IDXSFEPS(isf)
Indices for surface layer emissivity.
#define IDXCLK(icl)
Indices for cloud layer extinction.
#define TIMER(name, mode)
Start or stop a timer.
#define IDXP
Index for pressure.
#define IDXSFP
Index for surface layer pressure.
#define NG
Maximum number of emitters.
#define LOG(level,...)
Print log message.
#define NSF
Maximum number of surface layer spectral grid points.
#define NCL
Maximum number of cloud layer spectral grid points.
#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 MAX(a, b)
Macro to determine the maximum of two values.
#define IDXT
Index for temperature.
#define IDXSFZ
Index for surface layer height.
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[])
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)
Carry out optimal estimation retrieval.
void analyze_avk(ret_t *ret, ctl_t *ctl, atm_t *atm, int *iqa, int *ipa, gsl_matrix *avk)
Compute information content and resolution.
void read_ret(int argc, char *argv[], ctl_t *ctl, ret_t *ret)
Read retrieval control parameters.
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 write_stddev(const char *quantity, ret_t *ret, ctl_t *ctl, atm_t *atm, gsl_matrix *s)
Write retrieval error to file.
void matrix_invert(gsl_matrix *a)
Invert symmetric matrix.
double cost_function(gsl_vector *dx, gsl_vector *dy, gsl_matrix *s_a_inv, gsl_vector *sig_eps_inv)
Compute cost function.
void analyze_avk_quantity(gsl_matrix *avk, int iq, int *ipa, size_t *n0, size_t *n1, double *cont, double *res)
Analyze averaging kernels for individual retrieval target.
double sfeps[NSF]
Surface emissivity.
double k[NW][NP]
Extinction [km^-1].
double sfz
Surface height [km].
double lat[NP]
Latitude [deg].
double lon[NP]
Longitude [deg].
double t[NP]
Temperature [K].
double sfp
Surface pressure [hPa].
double clz
Cloud layer height [km].
int np
Number of data points.
double cldz
Cloud layer depth [km].
double sft
Surface temperature [K].
double z[NP]
Altitude [km].
double clk[NCL]
Cloud layer extinction [km^-1].
double q[NG][NP]
Volume mixing ratio [ppv].
double p[NP]
Pressure [hPa].
Forward model control parameters.
int nw
Number of spectral windows.
int ng
Number of emitters.
int nd
Number of radiance channels.
int ncl
Number of cloud layer spectral grid points.
int nsf
Number of surface layer spectral grid points.
Observation geometry and radiance data.
double rad[ND][NR]
Radiance [W/(m^2 sr cm^-1)].
int nr
Number of ray paths.
Retrieval control parameters.
double err_sfz
Surface height error [km].
double err_press_cz
Vertical correlation length for pressure error [km].
double err_sfp
Surface pressure error [hPa].
double err_press
Pressure error [%].
double err_k_cz[NW]
Vertical correlation length for extinction error [km].
int err_ana
Carry out error analysis (0=no, 1=yes).
double err_k_ch[NW]
Horizontal correlation length for extinction error [km].
double err_temp_cz
Vertical correlation length for temperature error [km].
double err_formod[ND]
Forward model error [%].
double conv_dmin
Minimum normalized step size in state space.
double err_temp
Temperature error [K].
double err_q_cz[NG]
Vertical correlation length for volume mixing ratio error [km].
double err_noise[ND]
Noise error [W/(m^2 sr cm^-1)].
double err_clz
Cloud height error [km].
double err_sft
Surface temperature error [K].
double err_clk[NCL]
Cloud extinction error [km^-1].
double err_temp_ch
Horizontal correlation length for temperature error [km].
int kernel_recomp
Re-computation of kernel matrix (number of iterations).
int conv_itmax
Maximum number of iterations.
double err_cldz
Cloud depth error [km].
double err_press_ch
Horizontal correlation length for pressure error [km].
double err_q_ch[NG]
Horizontal correlation length for volume mixing ratio error [km].
double err_q[NG]
Volume mixing ratio error [%].
char dir[LEN]
Working directory.
double err_sfeps[NSF]
Surface emissivity error.
double err_k[NW]
Extinction error [km^-1].