Carry out optimal estimation retrieval.
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
480
481
482
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
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
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
519 formod(ctl, atm_i, obs_i);
520
521
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
530 atm_i, obs_i, "x", "x", "r");
532
533
534 set_cov_meas(ret, ctl, obs_meas, sig_noise, sig_formod, sig_eps_inv);
535
536
537 sprintf(filename,
"%s/costs.tab", ret->
dir);
538 if (!(out = fopen(filename, "w")))
539 ERRMSG(
"Cannot create cost function file!");
540
541
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
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
556
557
558 LOG(1,
"it= %d / chi^2/m= %g", it, chisq);
559
560
561 fprintf(out, "%d %g %d %d\n", it, chisq, (int) m, (int) n);
562
563
564 kernel(ctl, atm_i, obs_i, k_i);
565
566
567
568
569
570
572
573
574 double chisq_old = chisq;
575
576
578 kernel(ctl, atm_i, obs_i, k_i);
579
580
583
584
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
592 for (it2 = 0; it2 < 20; it2++) {
593
594
595 gsl_matrix_memcpy(a, s_a_inv);
596 gsl_matrix_scale(a, 1 + lmpar);
597 gsl_matrix_add(a, cov);
598
599
600 gsl_linalg_cholesky_decomp(a);
601 gsl_linalg_cholesky_solve(a, b, x_step);
602
603
604 gsl_vector_add(x_i, x_step);
607 x2atm(ctl, x_i, atm_i);
608
609
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 }
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++)
627
628
629 formod(ctl, atm_i, obs_i);
630 obs2y(ctl, obs_i, y_i, NULL, NULL);
631
632
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
640
641
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
652 LOG(1,
"it= %d / chi^2/m= %g", it, chisq);
653
654
655 fprintf(out, "%d %g %d %d\n", it, chisq, (int) m, (int) n);
656
657
658 gsl_blas_ddot(x_step, b, &disq);
659 disq /= (double) n;
660
661
663 break;
664 }
665
666
667 fclose(out);
668
669
673 atm_i, obs_i, "y", "x", "r");
674
675
676
677
678
679
681
682
683 auxnm = gsl_matrix_alloc(n, m);
684 corr = gsl_matrix_alloc(n, n);
685 gain = gsl_matrix_alloc(n, m);
686
687
688
690 gsl_matrix_add(cov, s_a_inv);
691
692
695 atm_i, obs_i, "x", "x", "r");
697
698
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");
706
707
708
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");
716
717
720
721
724
725
726
727 gsl_blas_dgemm(CblasNoTrans, CblasNoTrans, 1.0, gain, k_i, 0.0, a);
729 atm_i, obs_i, "x", "x", "r");
730
731
733
734
735 gsl_matrix_free(auxnm);
736 gsl_matrix_free(corr);
737 gsl_matrix_free(gain);
738 }
739
740
741
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.
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.
void write_obs(const char *dirname, const char *filename, const ctl_t *ctl, const obs_t *obs)
Write observation 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 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.
#define LEN
Maximum length of ASCII data lines.
#define MIN(a, b)
Macro to determine the minimum of two values.
#define ERRMSG(...)
Print error message and quit program.
#define LOG(level,...)
Print log message.
#define MAX(a, b)
Macro to determine the maximum of two values.
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.
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 analyze_avk(ret_t *ret, ctl_t *ctl, atm_t *atm, int *iqa, int *ipa, gsl_matrix *avk)
Compute information content and resolution.
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.
int np
Number of data points.
Observation geometry and radiance data.
int err_ana
Carry out error analysis (0=no, 1=yes).
double conv_dmin
Minimum normalized step size in state space.
int kernel_recomp
Re-computation of kernel matrix (number of iterations).
int conv_itmax
Maximum number of iterations.