JURASSIC
Data Structures | Functions
retrieval.c File Reference

JURASSIC retrieval processor. More...

#include "jurassic.h"

Go to the source code of this file.

Data Structures

struct  ret_t
 Retrieval control parameters. More...
 

Functions

void analyze_avk (ret_t *ret, ctl_t *ctl, atm_t *atm, int *iqa, int *ipa, gsl_matrix *avk)
 Compute information content and resolution. More...
 
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. More...
 
double cost_function (gsl_vector *dx, gsl_vector *dy, gsl_matrix *s_a_inv, gsl_vector *sig_eps_inv)
 Compute cost function. More...
 
void matrix_invert (gsl_matrix *a)
 Invert symmetric matrix. More...
 
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. More...
 
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. More...
 
void read_ret (int argc, char *argv[], ctl_t *ctl, ret_t *ret)
 Read retrieval control parameters. More...
 
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. More...
 
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. More...
 
void write_stddev (const char *quantity, ret_t *ret, ctl_t *ctl, atm_t *atm, gsl_matrix *s)
 Write retrieval error to file. More...
 
int main (int argc, char *argv[])
 

Detailed Description

JURASSIC retrieval processor.

Definition in file retrieval.c.

Function Documentation

◆ analyze_avk()

void analyze_avk ( ret_t ret,
ctl_t ctl,
atm_t atm,
int *  iqa,
int *  ipa,
gsl_matrix *  avk 
)

Compute information content and resolution.

Definition at line 257 of file retrieval.c.

263 {
264
265 static atm_t atm_cont, atm_res;
266
267 int icl, ig, iq, isf, iw;
268
269 size_t i, n, n0[NQ], n1[NQ];
270
271 /* Get sizes... */
272 n = avk->size1;
273
274 /* Find sub-matrices for different quantities... */
275 for (iq = 0; iq < NQ; iq++) {
276 n0[iq] = N;
277 for (i = 0; i < n; i++) {
278 if (iqa[i] == iq && n0[iq] == N)
279 n0[iq] = i;
280 if (iqa[i] == iq)
281 n1[iq] = i - n0[iq] + 1;
282 }
283 }
284
285 /* Initialize... */
286 copy_atm(ctl, &atm_cont, atm, 1);
287 copy_atm(ctl, &atm_res, atm, 1);
288
289 /* Analyze quantities... */
290 analyze_avk_quantity(avk, IDXP, ipa, n0, n1, atm_cont.p, atm_res.p);
291 analyze_avk_quantity(avk, IDXT, ipa, n0, n1, atm_cont.t, atm_res.t);
292 for (ig = 0; ig < ctl->ng; ig++)
293 analyze_avk_quantity(avk, IDXQ(ig), ipa, n0, n1,
294 atm_cont.q[ig], atm_res.q[ig]);
295 for (iw = 0; iw < ctl->nw; iw++)
296 analyze_avk_quantity(avk, IDXK(iw), ipa, n0, n1,
297 atm_cont.k[iw], atm_res.k[iw]);
298 analyze_avk_quantity(avk, IDXCLZ, ipa, n0, n1, &atm_cont.clz, &atm_res.clz);
299 analyze_avk_quantity(avk, IDXCLDZ, ipa, n0, n1, &atm_cont.cldz,
300 &atm_res.cldz);
301 for (icl = 0; icl < ctl->ncl; icl++)
302 analyze_avk_quantity(avk, IDXCLK(icl), ipa, n0, n1,
303 &atm_cont.clk[icl], &atm_res.clk[icl]);
304 analyze_avk_quantity(avk, IDXSFZ, ipa, n0, n1, &atm_cont.sfz, &atm_res.sfz);
305 analyze_avk_quantity(avk, IDXSFP, ipa, n0, n1, &atm_cont.sfp, &atm_res.sfp);
306 analyze_avk_quantity(avk, IDXSFT, ipa, n0, n1, &atm_cont.sft, &atm_res.sft);
307 for (isf = 0; isf < ctl->nsf; isf++)
308 analyze_avk_quantity(avk, IDXSFEPS(isf), ipa, n0, n1,
309 &atm_cont.sfeps[isf], &atm_res.sfeps[isf]);
310
311 /* Write results to disk... */
312 write_atm(ret->dir, "atm_cont.tab", ctl, &atm_cont);
313 write_atm(ret->dir, "atm_res.tab", ctl, &atm_res);
314}
void write_atm(const char *dirname, const char *filename, const ctl_t *ctl, const atm_t *atm)
Write atmospheric data.
Definition: jurassic.c:5339
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
#define N
Maximum size of state vector.
Definition: jurassic.h:393
#define IDXCLZ
Index for cloud layer height.
Definition: jurassic.h:463
#define IDXCLDZ
Index for cloud layer depth.
Definition: jurassic.h:466
#define IDXK(iw)
Indices for extinction.
Definition: jurassic.h:460
#define NQ
Maximum number of quantities.
Definition: jurassic.h:398
#define IDXSFT
Index for surface layer temperature.
Definition: jurassic.h:478
#define IDXSFEPS(isf)
Indices for surface layer emissivity.
Definition: jurassic.h:481
#define IDXCLK(icl)
Indices for cloud layer extinction.
Definition: jurassic.h:469
#define IDXP
Index for pressure.
Definition: jurassic.h:451
#define IDXSFP
Index for surface layer pressure.
Definition: jurassic.h:475
#define IDXQ(ig)
Indices for volume mixing ratios.
Definition: jurassic.h:457
#define IDXT
Index for temperature.
Definition: jurassic.h:454
#define IDXSFZ
Index for surface layer height.
Definition: jurassic.h:472
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.
Definition: retrieval.c:318
Atmospheric data.
Definition: jurassic.h:488
double sfeps[NSF]
Surface emissivity.
Definition: jurassic.h:536
double k[NW][NP]
Extinction [km^-1].
Definition: jurassic.h:515
double sfz
Surface height [km].
Definition: jurassic.h:527
double t[NP]
Temperature [K].
Definition: jurassic.h:509
double sfp
Surface pressure [hPa].
Definition: jurassic.h:530
double clz
Cloud layer height [km].
Definition: jurassic.h:518
double cldz
Cloud layer depth [km].
Definition: jurassic.h:521
double sft
Surface temperature [K].
Definition: jurassic.h:533
double clk[NCL]
Cloud layer extinction [km^-1].
Definition: jurassic.h:524
double q[NG][NP]
Volume mixing ratio [ppv].
Definition: jurassic.h:512
double p[NP]
Pressure [hPa].
Definition: jurassic.h:506
int nw
Number of spectral windows.
Definition: jurassic.h:556
int ng
Number of emitters.
Definition: jurassic.h:544
int ncl
Number of cloud layer spectral grid points.
Definition: jurassic.h:562
int nsf
Number of surface layer spectral grid points.
Definition: jurassic.h:568
char dir[LEN]
Working directory.
Definition: retrieval.c:35
Here is the call graph for this function:

◆ analyze_avk_quantity()

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.

Definition at line 318 of file retrieval.c.

325 {
326
327 /* Loop over state vector elements... */
328 if (n0[iq] < N)
329 for (size_t i = 0; i < n1[iq]; i++) {
330
331 /* Get area of averaging kernel... */
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);
334
335 /* Get information density... */
336 res[ipa[n0[iq] + i]] = 1 / gsl_matrix_get(avk, n0[iq] + i, n0[iq] + i);
337 }
338}

◆ cost_function()

double cost_function ( gsl_vector *  dx,
gsl_vector *  dy,
gsl_matrix *  s_a_inv,
gsl_vector *  sig_eps_inv 
)

Compute cost function.

Definition at line 342 of file retrieval.c.

346 {
347
348 gsl_vector *x_aux, *y_aux;
349
350 double chisq_a, chisq_m = 0;
351
352 /* Get sizes... */
353 size_t m = dy->size;
354 size_t n = dx->size;
355
356 /* Allocate... */
357 x_aux = gsl_vector_alloc(n);
358 y_aux = gsl_vector_alloc(m);
359
360 /* Determine normalized cost function...
361 (chi^2 = 1/m * [dy^T * S_eps^{-1} * dy + dx^T * S_a^{-1} * dx]) */
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);
366
367 /* Free... */
368 gsl_vector_free(x_aux);
369 gsl_vector_free(y_aux);
370
371 /* Return cost function value... */
372 return (chisq_m + chisq_a) / (double) m;
373}
#define POW2(x)
Compute x^2.
Definition: jurassic.h:187

◆ matrix_invert()

void matrix_invert ( gsl_matrix *  a)

Invert symmetric matrix.

Definition at line 377 of file retrieval.c.

378 {
379
380 size_t diag = 1;
381
382 /* Get size... */
383 size_t n = a->size1;
384
385 /* Check if matrix is diagonal... */
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) {
389 diag = 0;
390 break;
391 }
392
393 /* Quick inversion of diagonal matrix... */
394 if (diag)
395 for (size_t i = 0; i < n; i++)
396 gsl_matrix_set(a, i, i, 1 / gsl_matrix_get(a, i, i));
397
398 /* Matrix inversion by means of Cholesky decomposition... */
399 else {
400 gsl_linalg_cholesky_decomp(a);
401 gsl_linalg_cholesky_invert(a);
402 }
403}

◆ matrix_product()

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.

Definition at line 407 of file retrieval.c.

411 {
412
413 gsl_matrix *aux;
414
415 /* Set sizes... */
416 size_t m = a->size1;
417 size_t n = a->size2;
418
419 /* Allocate... */
420 aux = gsl_matrix_alloc(m, n);
421
422 /* Compute A^T B A... */
423 if (transpose == 1) {
424
425 /* Compute B^1/2 A... */
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));
430
431 /* Compute A^T B A = (B^1/2 A)^T (B^1/2 A)... */
432 gsl_blas_dgemm(CblasTrans, CblasNoTrans, 1.0, aux, aux, 0.0, c);
433 }
434
435 /* Compute A B A^T... */
436 else if (transpose == 2) {
437
438 /* Compute A B^1/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));
443
444 /* Compute A B A^T = (A B^1/2) (A B^1/2)^T... */
445 gsl_blas_dgemm(CblasNoTrans, CblasTrans, 1.0, aux, aux, 0.0, c);
446 }
447
448 /* Free... */
449 gsl_matrix_free(aux);
450}

◆ optimal_estimation()

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.

Definition at line 454 of file retrieval.c.

460 {
461
462 static int ipa[N], iqa[N];
463
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;
467
468 FILE *out;
469
470 char filename[2 * LEN];
471
472 double chisq, disq = 0, lmpar = 0.001;
473
474 int icl, ig, ip, isf, it = 0, it2, iw;
475
476 size_t i, m, n;
477
478 /* ------------------------------------------------------------
479 Initialize...
480 ------------------------------------------------------------ */
481
482 /* Get sizes... */
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!");
487
488 /* Write info... */
489 LOG(1, "Problem size: m= %d / n= %d "
490 "(alloc= %.4g MB / stat= %.4g MB)",
491 (int) m, (int) n,
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.);
496
497 /* Allocate... */
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);
502
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);
515
516 /* Set initial state... */
517 copy_atm(ctl, atm_i, atm_apr, 0);
518 copy_obs(ctl, obs_i, obs_meas, 0);
519 formod(ctl, atm_i, obs_i);
520
521 /* Set state vectors and observation vectors... */
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);
526
527 /* Set inverse a priori covariance S_a^-1... */
528 set_cov_apr(ret, ctl, atm_apr, iqa, ipa, s_a_inv);
529 write_matrix(ret->dir, "matrix_cov_apr.tab", ctl, s_a_inv,
530 atm_i, obs_i, "x", "x", "r");
531 matrix_invert(s_a_inv);
532
533 /* Get measurement errors... */
534 set_cov_meas(ret, ctl, obs_meas, sig_noise, sig_formod, sig_eps_inv);
535
536 /* Create cost function file... */
537 sprintf(filename, "%s/costs.tab", ret->dir);
538 if (!(out = fopen(filename, "w")))
539 ERRMSG("Cannot create cost function file!");
540
541 /* Write header... */
542 fprintf(out,
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");
547
548 /* Determine dx = x_i - x_a and dy = y - F(x_i) ... */
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);
553
554 /* Compute cost function... */
555 chisq = cost_function(dx, dy, s_a_inv, sig_eps_inv);
556
557 /* Write info... */
558 LOG(1, "it= %d / chi^2/m= %g", it, chisq);
559
560 /* Write to cost function file... */
561 fprintf(out, "%d %g %d %d\n", it, chisq, (int) m, (int) n);
562
563 /* Compute initial kernel... */
564 kernel(ctl, atm_i, obs_i, k_i);
565
566 /* ------------------------------------------------------------
567 Levenberg-Marquardt minimization...
568 ------------------------------------------------------------ */
569
570 /* Outer loop... */
571 for (it = 1; it <= ret->conv_itmax; it++) {
572
573 /* Store current cost function value... */
574 double chisq_old = chisq;
575
576 /* Compute kernel matrix K_i... */
577 if (it > 1 && it % ret->kernel_recomp == 0)
578 kernel(ctl, atm_i, obs_i, k_i);
579
580 /* Compute K_i^T * S_eps^{-1} * K_i ... */
581 if (it == 1 || it % ret->kernel_recomp == 0)
582 matrix_product(k_i, sig_eps_inv, 1, cov);
583
584 /* Determine b = K_i^T * S_eps^{-1} * dy - S_a^{-1} * dx ... */
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);
590
591 /* Inner loop... */
592 for (it2 = 0; it2 < 20; it2++) {
593
594 /* Compute A = (1 + lmpar) * S_a^{-1} + K_i^T * S_eps^{-1} * K_i ... */
595 gsl_matrix_memcpy(a, s_a_inv);
596 gsl_matrix_scale(a, 1 + lmpar);
597 gsl_matrix_add(a, cov);
598
599 /* Solve A * x_step = b by means of Cholesky decomposition... */
600 gsl_linalg_cholesky_decomp(a);
601 gsl_linalg_cholesky_solve(a, b, x_step);
602
603 /* Update atmospheric state... */
604 gsl_vector_add(x_i, x_step);
605 copy_atm(ctl, atm_i, atm_apr, 0);
606 copy_obs(ctl, obs_i, obs_meas, 0);
607 x2atm(ctl, x_i, atm_i);
608
609 /* Check atmospheric state... */
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);
617 }
618 atm_i->clz = MAX(atm_i->clz, 0);
619 atm_i->cldz = MAX(atm_i->cldz, 0.1);
620 for (icl = 0; icl < ctl->ncl; icl++)
621 atm_i->clk[icl] = MAX(atm_i->clk[icl], 0);
622 atm_i->sfz = MAX(atm_i->sfz, 0);
623 atm_i->sfp = MAX(atm_i->sfp, 0);
624 atm_i->sft = MIN(MAX(atm_i->sft, 100), 400);
625 for (isf = 0; isf < ctl->nsf; isf++)
626 atm_i->sfeps[isf] = MIN(MAX(atm_i->sfeps[isf], 0), 1);
627
628 /* Forward calculation... */
629 formod(ctl, atm_i, obs_i);
630 obs2y(ctl, obs_i, y_i, NULL, NULL);
631
632 /* Determine dx = x_i - x_a and dy = y - F(x_i) ... */
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);
637
638 /* Compute cost function... */
639 chisq = cost_function(dx, dy, s_a_inv, sig_eps_inv);
640
641 /* Modify Levenberg-Marquardt parameter... */
642 if (chisq > chisq_old) {
643 lmpar *= 10;
644 gsl_vector_sub(x_i, x_step);
645 } else {
646 lmpar /= 10;
647 break;
648 }
649 }
650
651 /* Write info... */
652 LOG(1, "it= %d / chi^2/m= %g", it, chisq);
653
654 /* Write to cost function file... */
655 fprintf(out, "%d %g %d %d\n", it, chisq, (int) m, (int) n);
656
657 /* Get normalized step size in state space... */
658 gsl_blas_ddot(x_step, b, &disq);
659 disq /= (double) n;
660
661 /* Convergence test... */
662 if ((it == 1 || it % ret->kernel_recomp == 0) && disq < ret->conv_dmin)
663 break;
664 }
665
666 /* Close cost function file... */
667 fclose(out);
668
669 /* Store results... */
670 write_atm(ret->dir, "atm_final.tab", ctl, atm_i);
671 write_obs(ret->dir, "obs_final.tab", ctl, obs_i);
672 write_matrix(ret->dir, "matrix_kernel.tab", ctl, k_i,
673 atm_i, obs_i, "y", "x", "r");
674
675 /* ------------------------------------------------------------
676 Analysis of retrieval results...
677 ------------------------------------------------------------ */
678
679 /* Check if error analysis is requested... */
680 if (ret->err_ana) {
681
682 /* Allocate... */
683 auxnm = gsl_matrix_alloc(n, m);
684 corr = gsl_matrix_alloc(n, n);
685 gain = gsl_matrix_alloc(n, m);
686
687 /* Compute inverse retrieval covariance...
688 cov^{-1} = S_a^{-1} + K_i^T * S_eps^{-1} * K_i */
689 matrix_product(k_i, sig_eps_inv, 1, cov);
690 gsl_matrix_add(cov, s_a_inv);
691
692 /* Compute retrieval covariance... */
693 matrix_invert(cov);
694 write_matrix(ret->dir, "matrix_cov_ret.tab", ctl, cov,
695 atm_i, obs_i, "x", "x", "r");
696 write_stddev("total", ret, ctl, atm_i, cov);
697
698 /* Compute correlation matrix... */
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)));
704 write_matrix(ret->dir, "matrix_corr.tab", ctl, corr,
705 atm_i, obs_i, "x", "x", "r");
706
707 /* Compute gain matrix...
708 G = cov * K^T * S_eps^{-1} */
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);
714 write_matrix(ret->dir, "matrix_gain.tab", ctl, gain,
715 atm_i, obs_i, "x", "y", "c");
716
717 /* Compute retrieval error due to noise... */
718 matrix_product(gain, sig_noise, 2, a);
719 write_stddev("noise", ret, ctl, atm_i, a);
720
721 /* Compute retrieval error due to forward model errors... */
722 matrix_product(gain, sig_formod, 2, a);
723 write_stddev("formod", ret, ctl, atm_i, a);
724
725 /* Compute averaging kernel matrix
726 A = G * K ... */
727 gsl_blas_dgemm(CblasNoTrans, CblasNoTrans, 1.0, gain, k_i, 0.0, a);
728 write_matrix(ret->dir, "matrix_avk.tab", ctl, a,
729 atm_i, obs_i, "x", "x", "r");
730
731 /* Analyze averaging kernel matrix... */
732 analyze_avk(ret, ctl, atm_i, iqa, ipa, a);
733
734 /* Free... */
735 gsl_matrix_free(auxnm);
736 gsl_matrix_free(corr);
737 gsl_matrix_free(gain);
738 }
739
740 /* ------------------------------------------------------------
741 Finalize...
742 ------------------------------------------------------------ */
743
744 gsl_matrix_free(a);
745 gsl_matrix_free(cov);
746 gsl_matrix_free(k_i);
747 gsl_matrix_free(s_a_inv);
748
749 gsl_vector_free(b);
750 gsl_vector_free(dx);
751 gsl_vector_free(dy);
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);
761}
void x2atm(const ctl_t *ctl, const gsl_vector *x, atm_t *atm)
Decompose parameter vector or state vector.
Definition: jurassic.c:5899
void formod(const ctl_t *ctl, atm_t *atm, obs_t *obs)
Determine ray paths and compute radiative transfer.
Definition: jurassic.c:3026
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 write_obs(const char *dirname, const char *filename, const ctl_t *ctl, const obs_t *obs)
Write observation data.
Definition: jurassic.c:5673
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
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
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
void kernel(ctl_t *ctl, atm_t *atm, obs_t *obs, gsl_matrix *k)
Compute Jacobians.
Definition: jurassic.c:4006
#define LEN
Maximum length of ASCII data lines.
Definition: jurassic.h:383
#define MIN(a, b)
Macro to determine the minimum of two values.
Definition: jurassic.h:160
#define ERRMSG(...)
Print error message and quit program.
Definition: jurassic.h:237
#define LOG(level,...)
Print log message.
Definition: jurassic.h:221
#define MAX(a, b)
Macro to determine the maximum of two values.
Definition: jurassic.h:156
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.
Definition: retrieval.c:407
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.
Definition: retrieval.c:820
void analyze_avk(ret_t *ret, ctl_t *ctl, atm_t *atm, int *iqa, int *ipa, gsl_matrix *avk)
Compute information content and resolution.
Definition: retrieval.c:257
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.
Definition: retrieval.c:938
void write_stddev(const char *quantity, ret_t *ret, ctl_t *ctl, atm_t *atm, gsl_matrix *s)
Write retrieval error to file.
Definition: retrieval.c:982
void matrix_invert(gsl_matrix *a)
Invert symmetric matrix.
Definition: retrieval.c:377
double cost_function(gsl_vector *dx, gsl_vector *dy, gsl_matrix *s_a_inv, gsl_vector *sig_eps_inv)
Compute cost function.
Definition: retrieval.c:342
int np
Number of data points.
Definition: jurassic.h:491
Observation geometry and radiance data.
Definition: jurassic.h:734
int err_ana
Carry out error analysis (0=no, 1=yes).
Definition: retrieval.c:47
double conv_dmin
Minimum normalized step size in state space.
Definition: retrieval.c:44
int kernel_recomp
Re-computation of kernel matrix (number of iterations).
Definition: retrieval.c:38
int conv_itmax
Maximum number of iterations.
Definition: retrieval.c:41
Here is the call graph for this function:

◆ read_ret()

void read_ret ( int  argc,
char *  argv[],
ctl_t ctl,
ret_t ret 
)

Read retrieval control parameters.

Definition at line 765 of file retrieval.c.

769 {
770
771 /* Iteration control... */
772 ret->kernel_recomp =
773 (int) scan_ctl(argc, argv, "KERNEL_RECOMP", -1, "3", NULL);
774 ret->conv_itmax = (int) scan_ctl(argc, argv, "CONV_ITMAX", -1, "30", NULL);
775 ret->conv_dmin = scan_ctl(argc, argv, "CONV_DMIN", -1, "0.1", NULL);
776
777 /* Error analysis... */
778 ret->err_ana = (int) scan_ctl(argc, argv, "ERR_ANA", -1, "1", NULL);
779
780 for (int id = 0; id < ctl->nd; id++)
781 ret->err_formod[id] = scan_ctl(argc, argv, "ERR_FORMOD", id, "0", NULL);
782
783 for (int id = 0; id < ctl->nd; id++)
784 ret->err_noise[id] = scan_ctl(argc, argv, "ERR_NOISE", id, "0", NULL);
785
786 ret->err_press = scan_ctl(argc, argv, "ERR_PRESS", -1, "0", NULL);
787 ret->err_press_cz = scan_ctl(argc, argv, "ERR_PRESS_CZ", -1, "-999", NULL);
788 ret->err_press_ch = scan_ctl(argc, argv, "ERR_PRESS_CH", -1, "-999", NULL);
789
790 ret->err_temp = scan_ctl(argc, argv, "ERR_TEMP", -1, "0", NULL);
791 ret->err_temp_cz = scan_ctl(argc, argv, "ERR_TEMP_CZ", -1, "-999", NULL);
792 ret->err_temp_ch = scan_ctl(argc, argv, "ERR_TEMP_CH", -1, "-999", NULL);
793
794 for (int ig = 0; ig < ctl->ng; ig++) {
795 ret->err_q[ig] = scan_ctl(argc, argv, "ERR_Q", ig, "0", NULL);
796 ret->err_q_cz[ig] = scan_ctl(argc, argv, "ERR_Q_CZ", ig, "-999", NULL);
797 ret->err_q_ch[ig] = scan_ctl(argc, argv, "ERR_Q_CH", ig, "-999", NULL);
798 }
799
800 for (int iw = 0; iw < ctl->nw; iw++) {
801 ret->err_k[iw] = scan_ctl(argc, argv, "ERR_K", iw, "0", NULL);
802 ret->err_k_cz[iw] = scan_ctl(argc, argv, "ERR_K_CZ", iw, "-999", NULL);
803 ret->err_k_ch[iw] = scan_ctl(argc, argv, "ERR_K_CH", iw, "-999", NULL);
804 }
805
806 ret->err_clz = scan_ctl(argc, argv, "ERR_CLZ", -1, "0", NULL);
807 ret->err_cldz = scan_ctl(argc, argv, "ERR_CLDZ", -1, "0", NULL);
808 for (int icl = 0; icl < ctl->ncl; icl++)
809 ret->err_clk[icl] = scan_ctl(argc, argv, "ERR_CLK", icl, "0", NULL);
810
811 ret->err_sfz = scan_ctl(argc, argv, "ERR_SFZ", -1, "0", NULL);
812 ret->err_sfp = scan_ctl(argc, argv, "ERR_SFP", -1, "0", NULL);
813 ret->err_sft = scan_ctl(argc, argv, "ERR_SFT", -1, "0", NULL);
814 for (int isf = 0; isf < ctl->nsf; isf++)
815 ret->err_sfeps[isf] = scan_ctl(argc, argv, "ERR_SFEPS", isf, "0", NULL);
816}
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
int nd
Number of radiance channels.
Definition: jurassic.h:550
double err_sfz
Surface height error [km].
Definition: retrieval.c:101
double err_press_cz
Vertical correlation length for pressure error [km].
Definition: retrieval.c:59
double err_sfp
Surface pressure error [hPa].
Definition: retrieval.c:104
double err_press
Pressure error [%].
Definition: retrieval.c:56
double err_k_cz[NW]
Vertical correlation length for extinction error [km].
Definition: retrieval.c:86
double err_k_ch[NW]
Horizontal correlation length for extinction error [km].
Definition: retrieval.c:89
double err_temp_cz
Vertical correlation length for temperature error [km].
Definition: retrieval.c:68
double err_formod[ND]
Forward model error [%].
Definition: retrieval.c:50
double err_temp
Temperature error [K].
Definition: retrieval.c:65
double err_q_cz[NG]
Vertical correlation length for volume mixing ratio error [km].
Definition: retrieval.c:77
double err_noise[ND]
Noise error [W/(m^2 sr cm^-1)].
Definition: retrieval.c:53
double err_clz
Cloud height error [km].
Definition: retrieval.c:92
double err_sft
Surface temperature error [K].
Definition: retrieval.c:107
double err_clk[NCL]
Cloud extinction error [km^-1].
Definition: retrieval.c:98
double err_temp_ch
Horizontal correlation length for temperature error [km].
Definition: retrieval.c:71
double err_cldz
Cloud depth error [km].
Definition: retrieval.c:95
double err_press_ch
Horizontal correlation length for pressure error [km].
Definition: retrieval.c:62
double err_q_ch[NG]
Horizontal correlation length for volume mixing ratio error [km].
Definition: retrieval.c:80
double err_q[NG]
Volume mixing ratio error [%].
Definition: retrieval.c:74
double err_sfeps[NSF]
Surface emissivity error.
Definition: retrieval.c:110
double err_k[NW]
Extinction error [km^-1].
Definition: retrieval.c:83
Here is the call graph for this function:

◆ set_cov_apr()

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.

Definition at line 820 of file retrieval.c.

826 {
827
828 gsl_vector *x_a;
829
830 double x0[3], x1[3];
831
832 /* Get sizes... */
833 size_t n = s_a->size1;
834
835 /* Allocate... */
836 x_a = gsl_vector_alloc(n);
837
838 /* Get sigma vector... */
839 atm2x(ctl, atm, x_a, NULL, NULL);
840 for (size_t i = 0; i < n; i++) {
841 if (iqa[i] == IDXP)
842 gsl_vector_set(x_a, i, ret->err_press / 100 * gsl_vector_get(x_a, i));
843 if (iqa[i] == IDXT)
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]);
851 if (iqa[i] == IDXCLZ)
852 gsl_vector_set(x_a, i, ret->err_clz);
853 if (iqa[i] == IDXCLDZ)
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]);
858 if (iqa[i] == IDXSFZ)
859 gsl_vector_set(x_a, i, ret->err_sfz);
860 if (iqa[i] == IDXSFP)
861 gsl_vector_set(x_a, i, ret->err_sfp);
862 if (iqa[i] == IDXSFT)
863 gsl_vector_set(x_a, i, ret->err_sft);
864 for (int isf = 0; isf < ctl->nsf; isf++)
865 if (iqa[i] == IDXSFEPS(isf))
866 gsl_vector_set(x_a, i, ret->err_sfeps[isf]);
867 }
868
869 /* Check standard deviations... */
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)!");
873
874 /* Initialize diagonal covariance... */
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)));
878
879 /* Loop over matrix elements... */
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]) {
883
884 /* Initialize... */
885 double cz = 0;
886 double ch = 0;
887
888 /* Set correlation lengths for pressure... */
889 if (iqa[i] == IDXP) {
890 cz = ret->err_press_cz;
891 ch = ret->err_press_ch;
892 }
893
894 /* Set correlation lengths for temperature... */
895 if (iqa[i] == IDXT) {
896 cz = ret->err_temp_cz;
897 ch = ret->err_temp_ch;
898 }
899
900 /* Set correlation lengths for volume mixing ratios... */
901 for (int ig = 0; ig < ctl->ng; ig++)
902 if (iqa[i] == IDXQ(ig)) {
903 cz = ret->err_q_cz[ig];
904 ch = ret->err_q_ch[ig];
905 }
906
907 /* Set correlation lengths for extinction... */
908 for (int iw = 0; iw < ctl->nw; iw++)
909 if (iqa[i] == IDXK(iw)) {
910 cz = ret->err_k_cz[iw];
911 ch = ret->err_k_ch[iw];
912 }
913
914 /* Compute correlations... */
915 if (cz > 0 && ch > 0) {
916
917 /* Get Cartesian coordinates... */
918 geo2cart(0, atm->lon[ipa[i]], atm->lat[ipa[i]], x0);
919 geo2cart(0, atm->lon[ipa[j]], atm->lat[ipa[j]], x1);
920
921 /* Compute correlations... */
922 double rho =
923 exp(-DIST(x0, x1) / ch -
924 fabs(atm->z[ipa[i]] - atm->z[ipa[j]]) / cz);
925
926 /* Set covariance... */
927 gsl_matrix_set(s_a, i, j, gsl_vector_get(x_a, i)
928 * gsl_vector_get(x_a, j) * rho);
929 }
930 }
931
932 /* Free... */
933 gsl_vector_free(x_a);
934}
void geo2cart(const double z, const double lon, const double lat, double *x)
Convert geolocation to Cartesian coordinates.
Definition: jurassic.c:3500
#define DIST(a, b)
Compute Cartesian distance between two vectors.
Definition: jurassic.h:134
double lat[NP]
Latitude [deg].
Definition: jurassic.h:503
double lon[NP]
Longitude [deg].
Definition: jurassic.h:500
double z[NP]
Altitude [km].
Definition: jurassic.h:497
Here is the call graph for this function:

◆ set_cov_meas()

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.

Definition at line 938 of file retrieval.c.

944 {
945
946 static obs_t obs_err;
947
948 /* Get size... */
949 size_t m = sig_eps_inv->size;
950
951 /* Noise error (always considered in retrieval fit)... */
952 copy_obs(ctl, &obs_err, obs, 1);
953 for (int ir = 0; ir < obs_err.nr; ir++)
954 for (int id = 0; id < ctl->nd; id++)
955 obs_err.rad[id][ir]
956 = (isfinite(obs->rad[id][ir]) ? ret->err_noise[id] : NAN);
957 obs2y(ctl, &obs_err, sig_noise, NULL, NULL);
958
959 /* Forward model error (always considered in retrieval fit)... */
960 copy_obs(ctl, &obs_err, obs, 1);
961 for (int ir = 0; ir < obs_err.nr; ir++)
962 for (int id = 0; id < ctl->nd; id++)
963 obs_err.rad[id][ir]
964 = fabs(ret->err_formod[id] / 100 * obs->rad[id][ir]);
965 obs2y(ctl, &obs_err, sig_formod, NULL, NULL);
966
967 /* Total error... */
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))
970 +
971 POW2(gsl_vector_get
972 (sig_formod, i))));
973
974 /* Check standard deviations... */
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)!");
978}
double rad[ND][NR]
Radiance [W/(m^2 sr cm^-1)].
Definition: jurassic.h:773
int nr
Number of ray paths.
Definition: jurassic.h:737
Here is the call graph for this function:

◆ write_stddev()

void write_stddev ( const char *  quantity,
ret_t ret,
ctl_t ctl,
atm_t atm,
gsl_matrix *  s 
)

Write retrieval error to file.

Definition at line 982 of file retrieval.c.

987 {
988
989 static atm_t atm_aux;
990
991 gsl_vector *x_aux;
992
993 char filename[LEN];
994
995 /* Get sizes... */
996 size_t n = s->size1;
997
998 /* Allocate... */
999 x_aux = gsl_vector_alloc(n);
1000
1001 /* Compute standard deviation... */
1002 for (size_t i = 0; i < n; i++)
1003 gsl_vector_set(x_aux, i, sqrt(gsl_matrix_get(s, i, i)));
1004
1005 /* Write to disk... */
1006 copy_atm(ctl, &atm_aux, atm, 1);
1007 x2atm(ctl, x_aux, &atm_aux);
1008 sprintf(filename, "atm_err_%s.tab", quantity);
1009 write_atm(ret->dir, filename, ctl, &atm_aux);
1010
1011 /* Free... */
1012 gsl_vector_free(x_aux);
1013}
Here is the call graph for this function:

◆ main()

int main ( int  argc,
char *  argv[] 
)

Definition at line 201 of file retrieval.c.

203 {
204
205 static atm_t atm_i, atm_apr;
206 static ctl_t ctl;
207 static obs_t obs_i, obs_meas;
208 static ret_t ret;
209
210 FILE *dirlist;
211
212 /* Check arguments... */
213 if (argc < 3)
214 ERRMSG("Give parameters: <ctl> <dirlist>");
215
216 /* Measure CPU-time... */
217 TIMER("total", 1);
218
219 /* Read control parameters... */
220 read_ctl(argc, argv, &ctl);
221 read_ret(argc, argv, &ctl, &ret);
222
223 /* Open directory list... */
224 if (!(dirlist = fopen(argv[2], "r")))
225 ERRMSG("Cannot open directory list!");
226
227 /* Loop over directories... */
228 while (fscanf(dirlist, "%s", ret.dir) != EOF) {
229
230 /* Write info... */
231 LOG(1, "\nRetrieve in directory %s...\n", ret.dir);
232
233 /* Read atmospheric data... */
234 read_atm(ret.dir, "atm_apr.tab", &ctl, &atm_apr);
235
236 /* Read observation data... */
237 read_obs(ret.dir, "obs_meas.tab", &ctl, &obs_meas);
238
239 /* Run retrieval... */
240 optimal_estimation(&ret, &ctl, &obs_meas, &obs_i, &atm_apr, &atm_i);
241
242 /* Measure CPU-time... */
243 TIMER("total", 2);
244 }
245
246 /* Write info... */
247 LOG(1, "\nRetrieval done...");
248
249 /* Measure CPU-time... */
250 TIMER("total", 3);
251
252 return EXIT_SUCCESS;
253}
void read_ctl(int argc, char *argv[], ctl_t *ctl)
Read forward model control parameters.
Definition: jurassic.c:4550
void read_atm(const char *dirname, const char *filename, const ctl_t *ctl, atm_t *atm)
Read atmospheric data.
Definition: jurassic.c:4445
void read_obs(const char *dirname, const char *filename, const ctl_t *ctl, obs_t *obs)
Read observation data.
Definition: jurassic.c:4700
#define TIMER(name, mode)
Start or stop a timer.
Definition: jurassic.h:201
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.
Definition: retrieval.c:454
void read_ret(int argc, char *argv[], ctl_t *ctl, ret_t *ret)
Read retrieval control parameters.
Definition: retrieval.c:765
Forward model control parameters.
Definition: jurassic.h:541
Retrieval control parameters.
Definition: retrieval.c:32
Here is the call graph for this function: