AIRS Code Collection
Data Structures | Macros | Functions
nlte.c File Reference

Retrieval of non-LTE index. More...

#include <mpi.h>
#include <omp.h>
#include <netcdf.h>
#include "jurassic.h"

Go to the source code of this file.

Data Structures

struct  ncd_t
 Buffer for netCDF data. More...
 
struct  ret_t
 Retrieval results. More...
 

Macros

#define NC(cmd)
 Execute netCDF library command and check result. More...
 
#define L1_NCHAN   34
 Number of AIRS radiance channels (don't change). More...
 
#define L1_NTRACK   135
 Along-track size of AIRS radiance granule (don't change). More...
 
#define L1_NXTRACK   90
 Across-track size of AIRS radiance granule (don't change). More...
 
#define L2_NLAY   27
 Number of AIRS pressure layers (don't change). More...
 
#define L2_NTRACK   45
 Along-track size of AIRS retrieval granule (don't change). More...
 
#define L2_NXTRACK   30
 Across-track size of AIRS retrieval granule (don't change). More...
 

Functions

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. 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 fill_gaps (double x[L2_NTRACK][L2_NXTRACK][L2_NLAY], double cx, double cy)
 Fill data gaps in L2 data. More...
 
void init_l2 (ncd_t *ncd, int track, int xtrack, ctl_t *ctl, atm_t *atm)
 Initialize with AIRS Level-2 data. 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, tbl_t *tbl, obs_t *obs_meas, obs_t *obs_i, atm_t *atm_apr, atm_t *atm_i, double *chisq)
 Carry out optimal estimation retrieval. More...
 
void read_nc (char *filename, ncd_t *ncd)
 Read netCDF file. More...
 
void read_ret_ctl (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...
 
int main (int argc, char *argv[])
 

Detailed Description

Retrieval of non-LTE index.

Definition in file nlte.c.

Macro Definition Documentation

◆ NC

#define NC (   cmd)
Value:
{ \
int nc_result=(cmd); \
if(nc_result!=NC_NOERR) \
ERRMSG("%s", nc_strerror(nc_result)); \
}

Execute netCDF library command and check result.

Definition at line 36 of file nlte.c.

◆ L1_NCHAN

#define L1_NCHAN   34

Number of AIRS radiance channels (don't change).

Definition at line 47 of file nlte.c.

◆ L1_NTRACK

#define L1_NTRACK   135

Along-track size of AIRS radiance granule (don't change).

Definition at line 50 of file nlte.c.

◆ L1_NXTRACK

#define L1_NXTRACK   90

Across-track size of AIRS radiance granule (don't change).

Definition at line 53 of file nlte.c.

◆ L2_NLAY

#define L2_NLAY   27

Number of AIRS pressure layers (don't change).

Definition at line 56 of file nlte.c.

◆ L2_NTRACK

#define L2_NTRACK   45

Along-track size of AIRS retrieval granule (don't change).

Definition at line 59 of file nlte.c.

◆ L2_NXTRACK

#define L2_NXTRACK   30

Across-track size of AIRS retrieval granule (don't change).

Definition at line 62 of file nlte.c.

Function Documentation

◆ add_var()

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.

Add variable to netCDF file.

Definition at line 514 of file nlte.c.

522 {
523
524 /* Check if variable exists... */
525 if (nc_inq_varid(ncid, varname, varid) != NC_NOERR) {
526
527 /* Define variable... */
528 NC(nc_def_var(ncid, varname, type, ndims, dimid, varid));
529
530 /* Set long name... */
531 NC(nc_put_att_text
532 (ncid, *varid, "long_name", strlen(longname), longname));
533
534 /* Set units... */
535 NC(nc_put_att_text(ncid, *varid, "units", strlen(unit), unit));
536 }
537}
#define NC(cmd)
Execute netCDF library command and check result.
Definition: nlte.c:36

◆ 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 541 of file nlte.c.

545 {
546
547 double chisq_a, chisq_m = 0;
548
549 /* Get sizes... */
550 const size_t m = dy->size;
551 const size_t n = dx->size;
552
553 /* Allocate... */
554 gsl_vector *x_aux = gsl_vector_alloc(n);
555 gsl_vector *y_aux = gsl_vector_alloc(m);
556
557 /* Determine normalized cost function...
558 (chi^2 = 1/m * [dy^T * S_eps^{-1} * dy + dx^T * S_a^{-1} * dx]) */
559 for (size_t i = 0; i < m; i++)
560 chisq_m +=
561 gsl_pow_2(gsl_vector_get(dy, i) * gsl_vector_get(sig_eps_inv, i));
562 gsl_blas_dgemv(CblasNoTrans, 1.0, s_a_inv, dx, 0.0, x_aux);
563 gsl_blas_ddot(dx, x_aux, &chisq_a);
564
565 /* Free... */
566 gsl_vector_free(x_aux);
567 gsl_vector_free(y_aux);
568
569 /* Return cost function value... */
570 return (chisq_m + chisq_a) / (double) m;
571}

◆ fill_gaps()

void fill_gaps ( double  x[L2_NTRACK][L2_NXTRACK][L2_NLAY],
double  cx,
double  cy 
)

Fill data gaps in L2 data.

Definition at line 575 of file nlte.c.

578 {
579
580 double help[L2_NTRACK][L2_NXTRACK];
581
582 /* Loop over layers... */
583 for (int lay = 0; lay < L2_NLAY; lay++) {
584
585 /* Loop over grid points... */
586 for (int track = 0; track < L2_NTRACK; track++)
587 for (int xtrack = 0; xtrack < L2_NXTRACK; xtrack++) {
588
589 /* Init... */
590 help[track][xtrack] = 0;
591 double wsum = 0;
592
593 /* Averrage data points... */
594 for (int track2 = 0; track2 < L2_NTRACK; track2++)
595 for (int xtrack2 = 0; xtrack2 < L2_NXTRACK; xtrack2++)
596 if (gsl_finite(x[track2][xtrack2][lay])
597 && x[track2][xtrack2][lay] > 0) {
598 const double w = exp(-gsl_pow_2((xtrack - xtrack2) / cx)
599 - gsl_pow_2((track - track2) / cy));
600 help[track][xtrack] += w * x[track2][xtrack2][lay];
601 wsum += w;
602 }
603
604 /* Normalize... */
605 if (wsum > 0)
606 help[track][xtrack] /= wsum;
607 else
608 help[track][xtrack] = GSL_NAN;
609 }
610
611 /* Copy grid points... */
612 for (int track = 0; track < L2_NTRACK; track++)
613 for (int xtrack = 0; xtrack < L2_NXTRACK; xtrack++)
614 x[track][xtrack][lay] = help[track][xtrack];
615 }
616}
#define L2_NXTRACK
Across-track size of AIRS retrieval granule (don't change).
Definition: nlte.c:62
#define L2_NLAY
Number of AIRS pressure layers (don't change).
Definition: nlte.c:56
#define L2_NTRACK
Along-track size of AIRS retrieval granule (don't change).
Definition: nlte.c:59

◆ init_l2()

void init_l2 ( ncd_t ncd,
int  track,
int  xtrack,
ctl_t *  ctl,
atm_t *  atm 
)

Initialize with AIRS Level-2 data.

Definition at line 620 of file nlte.c.

625 {
626
627 static atm_t atm_airs;
628
629 double k[NW], p, q[NG], t, zmax = 0, zmin = 1000;
630
631 /* Reset track- and xtrack-index to match Level-2 data... */
632 track /= 3;
633 xtrack /= 3;
634
635 /* Store AIRS data in atmospheric data struct... */
636 atm_airs.np = 0;
637 for (int lay = 0; lay < L2_NLAY; lay++)
638 if (gsl_finite(ncd->l2_z[track][xtrack][lay])) {
639 atm_airs.z[atm_airs.np] = ncd->l2_z[track][xtrack][lay];
640 atm_airs.p[atm_airs.np] = ncd->l2_p[lay];
641 atm_airs.t[atm_airs.np] = ncd->l2_t[track][xtrack][lay];
642 if ((++atm_airs.np) > NP)
643 ERRMSG("Too many layers!");
644 }
645
646 /* Check number of levels... */
647 if (atm_airs.np <= 0)
648 return;
649
650 /* Get height range of AIRS data... */
651 for (int ip = 0; ip < atm_airs.np; ip++) {
652 zmax = GSL_MAX(zmax, atm_airs.z[ip]);
653 zmin = GSL_MIN(zmin, atm_airs.z[ip]);
654 }
655
656 /* Merge AIRS data... */
657 for (int ip = 0; ip < atm->np; ip++) {
658
659 /* Interpolate AIRS data... */
660 intpol_atm(ctl, &atm_airs, atm->z[ip], &p, &t, q, k);
661
662 /* Weighting factor... */
663 double w = 1;
664 if (atm->z[ip] > zmax)
665 w = GSL_MAX(1 - (atm->z[ip] - zmax) / 50, 0);
666 if (atm->z[ip] < zmin)
667 w = GSL_MAX(1 - (zmin - atm->z[ip]) / 50, 0);
668
669 /* Merge... */
670 atm->t[ip] = w * t + (1 - w) * atm->t[ip];
671 atm->p[ip] = w * p + (1 - w) * atm->p[ip];
672 }
673}
double l2_z[L2_NTRACK][L2_NXTRACK][L2_NLAY]
Altitude [km].
Definition: diff_apr.c:101
double l2_p[L2_NLAY]
Pressure [hPa].
Definition: diff_apr.c:104
double l2_t[L2_NTRACK][L2_NXTRACK][L2_NLAY]
Temperature [K].
Definition: diff_apr.c:107

◆ matrix_invert()

void matrix_invert ( gsl_matrix *  a)

Invert symmetric matrix.

Definition at line 677 of file nlte.c.

678 {
679
680 size_t diag = 1;
681
682 /* Get size... */
683 const size_t n = a->size1;
684
685 /* Check if matrix is diagonal... */
686 for (size_t i = 0; i < n && diag; i++)
687 for (size_t j = i + 1; j < n; j++)
688 if (gsl_matrix_get(a, i, j) != 0) {
689 diag = 0;
690 break;
691 }
692
693 /* Quick inversion of diagonal matrix... */
694 if (diag)
695 for (size_t i = 0; i < n; i++)
696 gsl_matrix_set(a, i, i, 1 / gsl_matrix_get(a, i, i));
697
698 /* Matrix inversion by means of Cholesky decomposition... */
699 else {
700 gsl_linalg_cholesky_decomp(a);
701 gsl_linalg_cholesky_invert(a);
702 }
703}

◆ 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 707 of file nlte.c.

711 {
712
713 /* Set sizes... */
714 const size_t m = a->size1;
715 const size_t n = a->size2;
716
717 /* Allocate... */
718 gsl_matrix *aux = gsl_matrix_alloc(m, n);
719
720 /* Compute A^T B A... */
721 if (transpose == 1) {
722
723 /* Compute B^1/2 A... */
724 for (size_t i = 0; i < m; i++)
725 for (size_t j = 0; j < n; j++)
726 gsl_matrix_set(aux, i, j,
727 gsl_vector_get(b, i) * gsl_matrix_get(a, i, j));
728
729 /* Compute A^T B A = (B^1/2 A)^T (B^1/2 A)... */
730 gsl_blas_dgemm(CblasTrans, CblasNoTrans, 1.0, aux, aux, 0.0, c);
731 }
732
733 /* Compute A B A^T... */
734 else if (transpose == 2) {
735
736 /* Compute A B^1/2... */
737 for (size_t i = 0; i < m; i++)
738 for (size_t j = 0; j < n; j++)
739 gsl_matrix_set(aux, i, j,
740 gsl_matrix_get(a, i, j) * gsl_vector_get(b, j));
741
742 /* Compute A B A^T = (A B^1/2) (A B^1/2)^T... */
743 gsl_blas_dgemm(CblasNoTrans, CblasTrans, 1.0, aux, aux, 0.0, c);
744 }
745
746 /* Free... */
747 gsl_matrix_free(aux);
748}

◆ optimal_estimation()

void optimal_estimation ( ret_t ret,
ctl_t *  ctl,
tbl_t *  tbl,
obs_t *  obs_meas,
obs_t *  obs_i,
atm_t *  atm_apr,
atm_t *  atm_i,
double *  chisq 
)

Carry out optimal estimation retrieval.

Definition at line 752 of file nlte.c.

760 {
761
762 static int ipa[N], iqa[N];
763
764 double disq = 0, lmpar = 0.001;
765
766 /* ------------------------------------------------------------
767 Initialize...
768 ------------------------------------------------------------ */
769
770 /* Get sizes... */
771 const size_t m = obs2y(ctl, obs_meas, NULL, NULL, NULL);
772 const size_t n = atm2x(ctl, atm_apr, NULL, iqa, ipa);
773 if (m == 0 || n == 0) {
774 *chisq = GSL_NAN;
775 return;
776 }
777
778 /* Allocate... */
779 gsl_matrix *a = gsl_matrix_alloc(n, n);
780 gsl_matrix *cov = gsl_matrix_alloc(n, n);
781 gsl_matrix *k_i = gsl_matrix_alloc(m, n);
782 gsl_matrix *s_a_inv = gsl_matrix_alloc(n, n);
783
784 gsl_vector *b = gsl_vector_alloc(n);
785 gsl_vector *dx = gsl_vector_alloc(n);
786 gsl_vector *dy = gsl_vector_alloc(m);
787 gsl_vector *sig_eps_inv = gsl_vector_alloc(m);
788 gsl_vector *sig_formod = gsl_vector_alloc(m);
789 gsl_vector *sig_noise = gsl_vector_alloc(m);
790 gsl_vector *x_a = gsl_vector_alloc(n);
791 gsl_vector *x_i = gsl_vector_alloc(n);
792 gsl_vector *x_step = gsl_vector_alloc(n);
793 gsl_vector *y_aux = gsl_vector_alloc(m);
794 gsl_vector *y_i = gsl_vector_alloc(m);
795 gsl_vector *y_m = gsl_vector_alloc(m);
796
797 /* Set initial state... */
798 copy_atm(ctl, atm_i, atm_apr, 0);
799 copy_obs(ctl, obs_i, obs_meas, 0);
800 formod(ctl, tbl, atm_i, obs_i);
801
802 /* Set state vectors and observation vectors... */
803 atm2x(ctl, atm_apr, x_a, NULL, NULL);
804 atm2x(ctl, atm_i, x_i, NULL, NULL);
805 obs2y(ctl, obs_meas, y_m, NULL, NULL);
806 obs2y(ctl, obs_i, y_i, NULL, NULL);
807
808 /* Set inverse a priori covariance S_a^-1... */
809 set_cov_apr(ret, ctl, atm_apr, iqa, ipa, s_a_inv);
810 matrix_invert(s_a_inv);
811
812 /* Get measurement errors... */
813 set_cov_meas(ret, ctl, obs_meas, sig_noise, sig_formod, sig_eps_inv);
814
815 /* Determine dx = x_i - x_a and dy = y - F(x_i) ... */
816 gsl_vector_memcpy(dx, x_i);
817 gsl_vector_sub(dx, x_a);
818 gsl_vector_memcpy(dy, y_m);
819 gsl_vector_sub(dy, y_i);
820
821 /* Compute cost function... */
822 *chisq = cost_function(dx, dy, s_a_inv, sig_eps_inv);
823
824 /* Compute initial kernel... */
825 kernel(ctl, tbl, atm_i, obs_i, k_i);
826
827 /* ------------------------------------------------------------
828 Levenberg-Marquardt minimization...
829 ------------------------------------------------------------ */
830
831 /* Outer loop... */
832 for (int it = 1; it <= ret->conv_itmax; it++) {
833
834 /* Store current cost function value... */
835 double chisq_old = *chisq;
836
837 /* Compute kernel matrix K_i... */
838 if (it > 1 && it % ret->kernel_recomp == 0)
839 kernel(ctl, tbl, atm_i, obs_i, k_i);
840
841 /* Compute K_i^T * S_eps^{-1} * K_i ... */
842 if (it == 1 || it % ret->kernel_recomp == 0)
843 matrix_product(k_i, sig_eps_inv, 1, cov);
844
845 /* Determine b = K_i^T * S_eps^{-1} * dy - S_a^{-1} * dx ... */
846 for (size_t i = 0; i < m; i++)
847 gsl_vector_set(y_aux, i, gsl_vector_get(dy, i)
848 * gsl_pow_2(gsl_vector_get(sig_eps_inv, i)));
849 gsl_blas_dgemv(CblasTrans, 1.0, k_i, y_aux, 0.0, b);
850 gsl_blas_dgemv(CblasNoTrans, -1.0, s_a_inv, dx, 1.0, b);
851
852 /* Inner loop... */
853 for (int it2 = 0; it2 < 20; it2++) {
854
855 /* Compute A = (1 + lmpar) * S_a^{-1} + K_i^T * S_eps^{-1} * K_i ... */
856 gsl_matrix_memcpy(a, s_a_inv);
857 gsl_matrix_scale(a, 1 + lmpar);
858 gsl_matrix_add(a, cov);
859
860 /* Solve A * x_step = b by means of Cholesky decomposition... */
861 gsl_linalg_cholesky_decomp(a);
862 gsl_linalg_cholesky_solve(a, b, x_step);
863
864 /* Update atmospheric state... */
865 gsl_vector_add(x_i, x_step);
866 copy_atm(ctl, atm_i, atm_apr, 0);
867 copy_obs(ctl, obs_i, obs_meas, 0);
868 x2atm(ctl, x_i, atm_i);
869
870 /* Check atmospheric state... */
871 for (int ip = 0; ip < atm_i->np; ip++) {
872 atm_i->p[ip] = GSL_MIN(GSL_MAX(atm_i->p[ip], 5e-7), 5e4);
873 atm_i->t[ip] = GSL_MIN(GSL_MAX(atm_i->t[ip], 100), 400);
874 for (int ig = 0; ig < ctl->ng; ig++)
875 atm_i->q[ig][ip] = GSL_MIN(GSL_MAX(atm_i->q[ig][ip], 0), 1);
876 for (int iw = 0; iw < ctl->nw; iw++)
877 atm_i->k[iw][ip] = GSL_MAX(atm_i->k[iw][ip], 0);
878 }
879
880 /* Forward calculation... */
881 formod(ctl, tbl, atm_i, obs_i);
882 obs2y(ctl, obs_i, y_i, NULL, NULL);
883
884 /* Determine dx = x_i - x_a and dy = y - F(x_i) ... */
885 gsl_vector_memcpy(dx, x_i);
886 gsl_vector_sub(dx, x_a);
887 gsl_vector_memcpy(dy, y_m);
888 gsl_vector_sub(dy, y_i);
889
890 /* Compute cost function... */
891 *chisq = cost_function(dx, dy, s_a_inv, sig_eps_inv);
892
893 /* Modify Levenberg-Marquardt parameter... */
894 if (*chisq > chisq_old) {
895 lmpar *= 10;
896 gsl_vector_sub(x_i, x_step);
897 } else {
898 lmpar /= 10;
899 break;
900 }
901 }
902
903 /* Get normalized step size in state space... */
904 gsl_blas_ddot(x_step, b, &disq);
905 disq /= (double) n;
906
907 /* Convergence test... */
908 if ((it == 1 || it % ret->kernel_recomp == 0) && disq < ret->conv_dmin)
909 break;
910 }
911
912 /* ------------------------------------------------------------
913 Finalize...
914 ------------------------------------------------------------ */
915
916 gsl_matrix_free(a);
917 gsl_matrix_free(cov);
918 gsl_matrix_free(k_i);
919 gsl_matrix_free(s_a_inv);
920
921 gsl_vector_free(b);
922 gsl_vector_free(dx);
923 gsl_vector_free(dy);
924 gsl_vector_free(sig_eps_inv);
925 gsl_vector_free(sig_formod);
926 gsl_vector_free(sig_noise);
927 gsl_vector_free(x_a);
928 gsl_vector_free(x_i);
929 gsl_vector_free(x_step);
930 gsl_vector_free(y_aux);
931 gsl_vector_free(y_i);
932 gsl_vector_free(y_m);
933}
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: nlte.c:707
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: nlte.c:1017
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: nlte.c:1115
void matrix_invert(gsl_matrix *a)
Invert symmetric matrix.
Definition: nlte.c:677
double cost_function(gsl_vector *dx, gsl_vector *dy, gsl_matrix *s_a_inv, gsl_vector *sig_eps_inv)
Compute cost function.
Definition: nlte.c:541
double conv_dmin
Minimum normalized step size in state space.
Definition: nlte.c:131
int kernel_recomp
Recomputation of kernel matrix (number of iterations).
Definition: nlte.c:125
int conv_itmax
Maximum number of iterations.
Definition: nlte.c:128
Here is the call graph for this function:

◆ read_nc()

void read_nc ( char *  filename,
ncd_t ncd 
)

Read netCDF file.

Definition at line 937 of file nlte.c.

939 {
940
941 int varid;
942
943 /* Open netCDF file... */
944 printf("Read netCDF file: %s\n", filename);
945 NC(nc_open(filename, NC_WRITE, &ncd->ncid));
946
947 /* Read Level-1 data... */
948 NC(nc_inq_varid(ncd->ncid, "l1_time", &varid));
949 NC(nc_get_var_double(ncd->ncid, varid, ncd->l1_time[0]));
950 NC(nc_inq_varid(ncd->ncid, "l1_lon", &varid));
951 NC(nc_get_var_double(ncd->ncid, varid, ncd->l1_lon[0]));
952 NC(nc_inq_varid(ncd->ncid, "l1_lat", &varid));
953 NC(nc_get_var_double(ncd->ncid, varid, ncd->l1_lat[0]));
954 NC(nc_inq_varid(ncd->ncid, "l1_sat_z", &varid));
955 NC(nc_get_var_double(ncd->ncid, varid, ncd->l1_sat_z));
956 NC(nc_inq_varid(ncd->ncid, "l1_sat_lon", &varid));
957 NC(nc_get_var_double(ncd->ncid, varid, ncd->l1_sat_lon));
958 NC(nc_inq_varid(ncd->ncid, "l1_sat_lat", &varid));
959 NC(nc_get_var_double(ncd->ncid, varid, ncd->l1_sat_lat));
960 NC(nc_inq_varid(ncd->ncid, "l1_nu", &varid));
961 NC(nc_get_var_double(ncd->ncid, varid, ncd->l1_nu));
962 NC(nc_inq_varid(ncd->ncid, "l1_rad", &varid));
963 NC(nc_get_var_float(ncd->ncid, varid, ncd->l1_rad[0][0]));
964
965 /* Read Level-2 data... */
966 NC(nc_inq_varid(ncd->ncid, "l2_z", &varid));
967 NC(nc_get_var_double(ncd->ncid, varid, ncd->l2_z[0][0]));
968 NC(nc_inq_varid(ncd->ncid, "l2_press", &varid));
969 NC(nc_get_var_double(ncd->ncid, varid, ncd->l2_p));
970 NC(nc_inq_varid(ncd->ncid, "l2_temp", &varid));
971 NC(nc_get_var_double(ncd->ncid, varid, ncd->l2_t[0][0]));
972}
double l1_lat[L1_NTRACK][L1_NXTRACK]
Footprint latitude [deg].
Definition: diff_apr.c:83
double l1_lon[L1_NTRACK][L1_NXTRACK]
Footprint longitude [deg].
Definition: diff_apr.c:80
double l1_sat_lat[L1_NTRACK]
Satellite latitude [deg].
Definition: diff_apr.c:92
int ncid
NetCDF file ID.
Definition: diff_apr.c:71
double l1_sat_lon[L1_NTRACK]
Satellite longitude [deg].
Definition: diff_apr.c:89
double l1_nu[L1_NCHAN]
Channel frequencies [cm^-1].
Definition: diff_apr.c:95
double l1_sat_z[L1_NTRACK]
Satellite altitude [km].
Definition: diff_apr.c:86
float l1_rad[L1_NTRACK][L1_NXTRACK][L1_NCHAN]
Radiance [W/(m^2 sr cm^-1)].
Definition: diff_apr.c:98
double l1_time[L1_NTRACK][L1_NXTRACK]
Time (seconds since 2000-01-01T00:00Z).
Definition: diff_apr.c:77

◆ read_ret_ctl()

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

Read retrieval control parameters.

Definition at line 976 of file nlte.c.

980 {
981
982 /* Iteration control... */
983 ret->kernel_recomp =
984 (int) scan_ctl(argc, argv, "KERNEL_RECOMP", -1, "3", NULL);
985 ret->conv_itmax = (int) scan_ctl(argc, argv, "CONV_ITMAX", -1, "30", NULL);
986 ret->conv_dmin = scan_ctl(argc, argv, "CONV_DMIN", -1, "0.1", NULL);
987
988 for (int id = 0; id < ctl->nd; id++)
989 ret->err_formod[id] = scan_ctl(argc, argv, "ERR_FORMOD", id, "0", NULL);
990
991 for (int id = 0; id < ctl->nd; id++)
992 ret->err_noise[id] = scan_ctl(argc, argv, "ERR_NOISE", id, "0", NULL);
993
994 ret->err_press = scan_ctl(argc, argv, "ERR_PRESS", -1, "0", NULL);
995 ret->err_press_cz = scan_ctl(argc, argv, "ERR_PRESS_CZ", -1, "-999", NULL);
996 ret->err_press_ch = scan_ctl(argc, argv, "ERR_PRESS_CH", -1, "-999", NULL);
997
998 ret->err_temp = scan_ctl(argc, argv, "ERR_TEMP", -1, "0", NULL);
999 ret->err_temp_cz = scan_ctl(argc, argv, "ERR_TEMP_CZ", -1, "-999", NULL);
1000 ret->err_temp_ch = scan_ctl(argc, argv, "ERR_TEMP_CH", -1, "-999", NULL);
1001
1002 for (int ig = 0; ig < ctl->ng; ig++) {
1003 ret->err_q[ig] = scan_ctl(argc, argv, "ERR_Q", ig, "0", NULL);
1004 ret->err_q_cz[ig] = scan_ctl(argc, argv, "ERR_Q_CZ", ig, "-999", NULL);
1005 ret->err_q_ch[ig] = scan_ctl(argc, argv, "ERR_Q_CH", ig, "-999", NULL);
1006 }
1007
1008 for (int iw = 0; iw < ctl->nw; iw++) {
1009 ret->err_k[iw] = scan_ctl(argc, argv, "ERR_K", iw, "0", NULL);
1010 ret->err_k_cz[iw] = scan_ctl(argc, argv, "ERR_K_CZ", iw, "-999", NULL);
1011 ret->err_k_ch[iw] = scan_ctl(argc, argv, "ERR_K_CH", iw, "-999", NULL);
1012 }
1013}
double err_k_ch[NW]
Horizontal correlation length for extinction error [km].
Definition: nlte.c:173
double err_press_cz
Vertical correlation length for pressure error [km].
Definition: nlte.c:143
double err_k_cz[NW]
Vertical correlation length for extinction error [km].
Definition: nlte.c:170
double err_press
Pressure error [%].
Definition: nlte.c:140
double err_q_ch[NG]
Horizontal correlation length for volume mixing ratio error [km].
Definition: nlte.c:164
double err_noise[ND]
Noise error [W/(m^2 sr cm^-1)].
Definition: nlte.c:137
double err_formod[ND]
Forward model error [%].
Definition: nlte.c:134
double err_q_cz[NG]
Vertical correlation length for volume mixing ratio error [km].
Definition: nlte.c:161
double err_temp_cz
Vertical correlation length for temperature error [km].
Definition: nlte.c:152
double err_temp
Temperature error [K].
Definition: nlte.c:149
double err_temp_ch
Horizontal correlation length for temperature error [km].
Definition: nlte.c:155
double err_press_ch
Horizontal correlation length for pressure error [km].
Definition: nlte.c:146
double err_q[NG]
Volume mixing ratio error [%].
Definition: nlte.c:158
double err_k[NW]
Extinction error [1/km].
Definition: nlte.c:167

◆ 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 1017 of file nlte.c.

1023 {
1024
1025 /* Get sizes... */
1026 const size_t n = s_a->size1;
1027
1028 /* Allocate... */
1029 gsl_vector *x_a = gsl_vector_alloc(n);
1030
1031 /* Get sigma vector... */
1032 atm2x(ctl, atm, x_a, NULL, NULL);
1033 for (size_t i = 0; i < n; i++) {
1034 if (iqa[i] == IDXP)
1035 gsl_vector_set(x_a, i, ret->err_press / 100 * gsl_vector_get(x_a, i));
1036 if (iqa[i] == IDXT)
1037 gsl_vector_set(x_a, i, ret->err_temp);
1038 for (int ig = 0; ig < ctl->ng; ig++)
1039 if (iqa[i] == IDXQ(ig))
1040 gsl_vector_set(x_a, i, ret->err_q[ig] / 100 * gsl_vector_get(x_a, i));
1041 for (int iw = 0; iw < ctl->nw; iw++)
1042 if (iqa[i] == IDXK(iw))
1043 gsl_vector_set(x_a, i, ret->err_k[iw]);
1044 }
1045
1046 /* Check standard deviations... */
1047 for (size_t i = 0; i < n; i++)
1048 if (gsl_pow_2(gsl_vector_get(x_a, i)) <= 0)
1049 ERRMSG("Check a priori data (zero standard deviation)!");
1050
1051 /* Initialize diagonal covariance... */
1052 gsl_matrix_set_zero(s_a);
1053 for (size_t i = 0; i < n; i++)
1054 gsl_matrix_set(s_a, i, i, gsl_pow_2(gsl_vector_get(x_a, i)));
1055
1056 /* Loop over matrix elements... */
1057 for (size_t i = 0; i < n; i++)
1058 for (size_t j = 0; j < n; j++)
1059 if (i != j && iqa[i] == iqa[j]) {
1060
1061 /* Initialize... */
1062 double cz = 0, ch = 0;
1063
1064 /* Set correlation lengths for pressure... */
1065 if (iqa[i] == IDXP) {
1066 cz = ret->err_press_cz;
1067 ch = ret->err_press_ch;
1068 }
1069
1070 /* Set correlation lengths for temperature... */
1071 if (iqa[i] == IDXT) {
1072 cz = ret->err_temp_cz;
1073 ch = ret->err_temp_ch;
1074 }
1075
1076 /* Set correlation lengths for volume mixing ratios... */
1077 for (int ig = 0; ig < ctl->ng; ig++)
1078 if (iqa[i] == IDXQ(ig)) {
1079 cz = ret->err_q_cz[ig];
1080 ch = ret->err_q_ch[ig];
1081 }
1082
1083 /* Set correlation lengths for extinction... */
1084 for (int iw = 0; iw < ctl->nw; iw++)
1085 if (iqa[i] == IDXK(iw)) {
1086 cz = ret->err_k_cz[iw];
1087 ch = ret->err_k_ch[iw];
1088 }
1089
1090 /* Compute correlations... */
1091 if (cz > 0 && ch > 0) {
1092
1093 /* Get Cartesian coordinates... */
1094 double x0[3], x1[3];
1095 geo2cart(0, atm->lon[ipa[i]], atm->lat[ipa[i]], x0);
1096 geo2cart(0, atm->lon[ipa[j]], atm->lat[ipa[j]], x1);
1097
1098 /* Compute correlations... */
1099 const double rho =
1100 exp(-DIST(x0, x1) / ch -
1101 fabs(atm->z[ipa[i]] - atm->z[ipa[j]]) / cz);
1102
1103 /* Set covariance... */
1104 gsl_matrix_set(s_a, i, j, gsl_vector_get(x_a, i)
1105 * gsl_vector_get(x_a, j) * rho);
1106 }
1107 }
1108
1109 /* Free... */
1110 gsl_vector_free(x_a);
1111}

◆ 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 1115 of file nlte.c.

1121 {
1122
1123 static obs_t obs_err;
1124
1125 /* Get size... */
1126 const size_t m = sig_eps_inv->size;
1127
1128 /* Noise error (always considered in retrieval fit)... */
1129 copy_obs(ctl, &obs_err, obs, 1);
1130 for (int ir = 0; ir < obs_err.nr; ir++)
1131 for (int id = 0; id < ctl->nd; id++)
1132 obs_err.rad[id][ir]
1133 = (gsl_finite(obs->rad[id][ir]) ? ret->err_noise[id] : GSL_NAN);
1134 obs2y(ctl, &obs_err, sig_noise, NULL, NULL);
1135
1136 /* Forward model error (always considered in retrieval fit)... */
1137 copy_obs(ctl, &obs_err, obs, 1);
1138 for (int ir = 0; ir < obs_err.nr; ir++)
1139 for (int id = 0; id < ctl->nd; id++)
1140 obs_err.rad[id][ir]
1141 = fabs(ret->err_formod[id] / 100 * obs->rad[id][ir]);
1142 obs2y(ctl, &obs_err, sig_formod, NULL, NULL);
1143
1144 /* Total error... */
1145 for (size_t i = 0; i < m; i++)
1146 gsl_vector_set(sig_eps_inv, i,
1147 1 / sqrt(gsl_pow_2(gsl_vector_get(sig_noise, i))
1148 + gsl_pow_2(gsl_vector_get(sig_formod, i))));
1149
1150 /* Check standard deviations... */
1151 for (size_t i = 0; i < m; i++)
1152 if (gsl_vector_get(sig_eps_inv, i) <= 0)
1153 ERRMSG("Check measurement errors (zero standard deviation)!");
1154}

◆ main()

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

Definition at line 269 of file nlte.c.

271 {
272
273 static ctl_t ctl;
274 static atm_t atm_apr, atm_clim, atm_i;
275 static obs_t obs_i, obs_meas;
276 static ncd_t ncd;
277 static ret_t ret;
278
279 static FILE *in, *out;
280
281 static char filename[LEN], filename2[2 * LEN];
282
283 static double chisq[L1_NTRACK][L1_NXTRACK], ni[L1_NTRACK][L1_NXTRACK],
284 z[NP];
285
286 static int channel[ND], n, ntask = -1, rank, size;
287
288 /* ------------------------------------------------------------
289 Init...
290 ------------------------------------------------------------ */
291
292 /* MPI... */
293 MPI_Init(&argc, &argv);
294 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
295 MPI_Comm_size(MPI_COMM_WORLD, &size);
296
297 /* Measure CPU time... */
298 TIMER("total", 1);
299
300 /* Check arguments... */
301 if (argc < 3)
302 ERRMSG("Give parameters: <ctl> <filelist>");
303
304 /* Read control parameters... */
305 read_ctl(argc, argv, &ctl);
306 read_ret_ctl(argc, argv, &ctl, &ret);
307
308 /* Initialize look-up tables... */
309 tbl_t *tbl = read_tbl(&ctl);
310
311 /* Read retrieval grid... */
312 const int nz = (int) scan_ctl(argc, argv, "NZ", -1, "", NULL);
313 if (nz > NP)
314 ERRMSG("Too many altitudes!");
315 for (int iz = 0; iz < nz; iz++)
316 z[iz] = scan_ctl(argc, argv, "Z", iz, "", NULL);
317
318 /* Read track range... */
319 const int track0 = (int) scan_ctl(argc, argv, "TRACK_MIN", -1, "0", NULL);
320 const int track1 = (int) scan_ctl(argc, argv, "TRACK_MAX", -1, "134", NULL);
321
322 /* Read xtrack range... */
323 const int xtrack0 = (int) scan_ctl(argc, argv, "XTRACK_MIN", -1, "0", NULL);
324 const int xtrack1 =
325 (int) scan_ctl(argc, argv, "XTRACK_MAX", -1, "89", NULL);
326
327 /* Background smoothing... */
328 const double sx = scan_ctl(argc, argv, "SX", -1, "8", NULL);
329 const double sy = scan_ctl(argc, argv, "SY", -1, "2", NULL);
330
331 /* ------------------------------------------------------------
332 Distribute granules...
333 ------------------------------------------------------------ */
334
335 /* Open filelist... */
336 printf("Read filelist: %s\n", argv[2]);
337 if (!(in = fopen(argv[2], "r")))
338 ERRMSG("Cannot open filelist!");
339
340 /* Loop over netCDF files... */
341 while (fscanf(in, "%s", filename) != EOF) {
342
343 /* Distribute files with MPI... */
344 if ((++ntask) % size != rank)
345 continue;
346
347 /* Write info... */
348 printf("Retrieve file %s on rank %d of %d (with %d threads)...\n",
349 filename, rank + 1, size, omp_get_max_threads());
350
351 /* ------------------------------------------------------------
352 Initialize retrieval...
353 ------------------------------------------------------------ */
354
355 /* Read netCDF file... */
356 read_nc(filename, &ncd);
357
358 /* Identify radiance channels... */
359 for (int id = 0; id < ctl.nd; id++) {
360 channel[id] = -999;
361 for (int i = 0; i < L1_NCHAN; i++)
362 if (fabs(ctl.nu[id] - ncd.l1_nu[i]) < 0.1)
363 channel[id] = i;
364 if (channel[id] < 0)
365 ERRMSG("Cannot identify radiance channel!");
366 }
367
368 /* Fill data gaps... */
369 fill_gaps(ncd.l2_t, sx, sy);
370 fill_gaps(ncd.l2_z, sx, sy);
371
372 /* Set climatological data for center of granule... */
373 atm_clim.np = nz;
374 for (int iz = 0; iz < nz; iz++)
375 atm_clim.z[iz] = z[iz];
376 climatology(&ctl, &atm_clim);
377
378 /* ------------------------------------------------------------
379 Retrieval...
380 ------------------------------------------------------------ */
381
382 /* Loop over swaths... */
383 for (int track = track0; track <= track1; track++) {
384
385 /* Measure CPU time... */
386 TIMER("retrieval", 1);
387
388 /* Loop over scan... */
389 for (int xtrack = xtrack0; xtrack <= xtrack1; xtrack++) {
390
391 /* Store observation data... */
392 obs_meas.nr = 1;
393 obs_meas.time[0] = ncd.l1_time[track][xtrack];
394 obs_meas.obsz[0] = ncd.l1_sat_z[track];
395 obs_meas.obslon[0] = ncd.l1_sat_lon[track];
396 obs_meas.obslat[0] = ncd.l1_sat_lat[track];
397 obs_meas.vplon[0] = ncd.l1_lon[track][xtrack];
398 obs_meas.vplat[0] = ncd.l1_lat[track][xtrack];
399 for (int id = 0; id < ctl.nd; id++)
400 obs_meas.rad[id][0] = ncd.l1_rad[track][xtrack][channel[id]];
401
402 /* Flag out 4 micron channels... */
403 for (int id = 0; id < ctl.nd; id++)
404 if (ctl.nu[id] >= 2000)
405 obs_meas.rad[id][0] = GSL_NAN;
406
407 /* Prepare atmospheric data... */
408 copy_atm(&ctl, &atm_apr, &atm_clim, 0);
409 for (int ip = 0; ip < atm_apr.np; ip++) {
410 atm_apr.time[ip] = obs_meas.time[0];
411 atm_apr.lon[ip] = obs_meas.vplon[0];
412 atm_apr.lat[ip] = obs_meas.vplat[0];
413 }
414
415 /* Merge Level-2 data... */
416 init_l2(&ncd, track, xtrack, &ctl, &atm_apr);
417
418 /* Retrieval... */
419 optimal_estimation(&ret, &ctl, tbl, &obs_meas, &obs_i,
420 &atm_apr, &atm_i, &chisq[track][xtrack]);
421
422 /* Run forward model including 4 micron channels... */
423 for (int id = 0; id < ctl.nd; id++)
424 obs_meas.rad[id][0] = obs_i.rad[id][0]
425 = ncd.l1_rad[track][xtrack][channel[id]];
426 formod(&ctl, tbl, &atm_i, &obs_i);
427
428 /* Calculate non-LTE index... */
429 n = 0;
430 ni[track][xtrack] = 0;
431 for (int id = 0; id < ctl.nd; id++)
432 if (ctl.nu[id] >= 2000 && gsl_finite(obs_meas.rad[id][0])) {
433 ni[track][xtrack] +=
434 (BRIGHT(obs_meas.rad[id][0], ncd.l1_nu[channel[id]])
435 - BRIGHT(obs_i.rad[id][0], ncd.l1_nu[channel[id]]));
436 n++;
437 }
438 ni[track][xtrack] /= n;
439 }
440
441 /* Measure CPU time... */
442 TIMER("retrieval", 3);
443 }
444
445 /* ------------------------------------------------------------
446 Write output...
447 ------------------------------------------------------------ */
448
449 /* Set filename... */
450 sprintf(filename2, "%s.nlte.tab", filename);
451
452 /* Create output file... */
453 printf("Write non-LTE data: %s\n", filename2);
454 if (!(out = fopen(filename2, "w")))
455 ERRMSG("Cannot create file!");
456
457 /* Write header... */
458 fprintf(out,
459 "# $1 = time (seconds since 2000-01-01, 00:00 UTC)\n"
460 "# $2 = longitude [deg]\n"
461 "# $3 = latitude [deg]\n"
462 "# $4 = solar zenith angle [deg]\n"
463 "# $5 = non-LTE index [K]\n" "# $6 = chi^2 of retrieval fit\n");
464
465 /* Write data... */
466 for (int track = track0; track <= track1; track++) {
467 fprintf(out, "\n");
468 for (int xtrack = xtrack0; xtrack <= xtrack1; xtrack++)
469 fprintf(out, "%.2f %g %g %g %g %g\n",
470 ncd.l1_time[track][xtrack],
471 ncd.l1_lon[track][xtrack],
472 ncd.l1_lat[track][xtrack],
473 sza(ncd.l1_time[track][xtrack], ncd.l1_lon[track][xtrack],
474 ncd.l1_lat[track][xtrack]), ni[track][xtrack],
475 chisq[track][xtrack]);
476 }
477
478 /* Close output file... */
479 fclose(out);
480
481 /* Write info... */
482 printf("Retrieval finished on rank %d of %d!\n", rank, size);
483 }
484
485 /* Close file list... */
486 fclose(in);
487
488 /* Free... */
489 free(tbl);
490
491 /* Measure CPU time... */
492 TIMER("total", 3);
493
494 /* Report memory usage... */
495 printf("MEMORY_ATM = %g MByte\n", 4. * sizeof(atm_t) / 1024. / 1024.);
496 printf("MEMORY_CTL = %g MByte\n", 1. * sizeof(ctl_t) / 1024. / 1024.);
497 printf("MEMORY_NCD = %g MByte\n", 1. * sizeof(ncd_t) / 1024. / 1024.);
498 printf("MEMORY_OBS = %g MByte\n", 3. * sizeof(atm_t) / 1024. / 1024.);
499 printf("MEMORY_RET = %g MByte\n", 1. * sizeof(ret_t) / 1024. / 1024.);
500 printf("MEMORY_TBL = %g MByte\n", 1. * sizeof(tbl_t) / 1024. / 1024.);
501
502 /* Report problem size... */
503 printf("SIZE_TASKS = %d\n", size);
504 printf("SIZE_THREADS = %d\n", omp_get_max_threads());
505
506 /* MPI... */
507 MPI_Finalize();
508
509 return EXIT_SUCCESS;
510}
void fill_gaps(double x[L2_NTRACK][L2_NXTRACK][L2_NLAY], double cx, double cy)
Fill data gaps in L2 data.
Definition: nlte.c:575
#define L1_NXTRACK
Across-track size of AIRS radiance granule (don't change).
Definition: nlte.c:53
#define L1_NTRACK
Along-track size of AIRS radiance granule (don't change).
Definition: nlte.c:50
void read_nc(char *filename, ncd_t *ncd)
Read netCDF file.
Definition: nlte.c:937
void optimal_estimation(ret_t *ret, ctl_t *ctl, tbl_t *tbl, obs_t *obs_meas, obs_t *obs_i, atm_t *atm_apr, atm_t *atm_i, double *chisq)
Carry out optimal estimation retrieval.
Definition: nlte.c:752
void read_ret_ctl(int argc, char *argv[], ctl_t *ctl, ret_t *ret)
Read retrieval control parameters.
Definition: nlte.c:976
#define L1_NCHAN
Number of AIRS radiance channels (don't change).
Definition: nlte.c:47
void init_l2(ncd_t *ncd, int track, int xtrack, ctl_t *ctl, atm_t *atm)
Initialize with AIRS Level-2 data.
Definition: nlte.c:620
Buffer for netCDF data.
Definition: diff_apr.c:68
Retrieval results.
Definition: libairs.h:237
Here is the call graph for this function: