37 int nc_result=(cmd); \
38 if(nc_result!=NC_NOERR) \
39 ERRMSG("%s", nc_strerror(nc_result)); \
130 const char *longname,
190 static atm_t atm_apr, atm_clim, atm_i;
191 static obs_t obs_i, obs_meas;
199 double chisq, chisq_min, chisq_max, chisq_mean, z[NP];
201 int channel[ND], m, ntask = -1, rank, size;
208 MPI_Init(&argc, &argv);
209 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
210 MPI_Comm_size(MPI_COMM_WORLD, &size);
217 ERRMSG(
"Give parameters: <ctl> <filelist>");
220 read_ctl(argc, argv, &ctl);
221 read_ret(argc, argv, &ctl, &ret);
224 tbl_t *tbl = read_tbl(&ctl);
227 const int nz = (int) scan_ctl(argc, argv,
"NZ", -1,
"", NULL);
229 ERRMSG(
"Too many altitudes!");
230 for (
int iz = 0; iz < nz; iz++)
231 z[iz] = scan_ctl(argc, argv,
"Z", iz,
"", NULL);
234 const int track0 = (int) scan_ctl(argc, argv,
"TRACK_MIN", -1,
"0", NULL);
235 const int track1 = (int) scan_ctl(argc, argv,
"TRACK_MAX", -1,
"134", NULL);
238 const int xtrack0 = (int) scan_ctl(argc, argv,
"XTRACK_MIN", -1,
"0", NULL);
240 (int) scan_ctl(argc, argv,
"XTRACK_MAX", -1,
"89", NULL);
243 int np0 = (int) scan_ctl(argc, argv,
"NP_MIN", -1,
"0", NULL);
244 int np1 = (int) scan_ctl(argc, argv,
"NP_MAX", -1,
"100", NULL);
245 np1 = GSL_MIN(np1, nz - 1);
248 const double sx = scan_ctl(argc, argv,
"SX", -1,
"8", NULL);
249 const double sy = scan_ctl(argc, argv,
"SY", -1,
"2", NULL);
252 const double sza_thresh = scan_ctl(argc, argv,
"SZA", -1,
"96", NULL);
259 printf(
"Read filelist: %s\n", argv[2]);
260 if (!(in = fopen(argv[2],
"r")))
261 ERRMSG(
"Cannot open filelist!");
264 while (fscanf(in,
"%s", filename) != EOF) {
267 if ((++ntask) % size != rank)
271 printf(
"Retrieve file %s on rank %d of %d (with %d threads)...\n",
272 filename, rank + 1, size, omp_get_max_threads());
282 for (
int id = 0;
id < ctl.nd;
id++) {
285 if (fabs(ctl.nu[
id] - ncd.
l1_nu[i]) < 0.1)
288 ERRMSG(
"Cannot identify radiance channel!");
297 for (
int iz = 0; iz < nz; iz++)
298 atm_clim.z[iz] = z[iz];
299 climatology(&ctl, &atm_clim);
312 for (
int track = track0; track <= track1; track++) {
315 TIMER(
"retrieval", 1);
318 for (
int xtrack = xtrack0; xtrack <= xtrack1; xtrack++) {
322 obs_meas.time[0] = ncd.
l1_time[track][xtrack];
323 obs_meas.obsz[0] = ncd.
l1_sat_z[track];
326 obs_meas.vplon[0] = ncd.
l1_lon[track][xtrack];
327 obs_meas.vplat[0] = ncd.
l1_lat[track][xtrack];
328 for (
int id = 0;
id < ctl.nd;
id++)
329 obs_meas.rad[
id][0] = ncd.
l1_rad[track][xtrack][channel[
id]];
332 if (sza(obs_meas.time[0], obs_meas.obslon[0], obs_meas.obslat[0])
334 for (
int id = 0;
id < ctl.nd;
id++)
335 if (ctl.nu[
id] >= 2000)
336 obs_meas.rad[id][0] = GSL_NAN;
339 copy_atm(&ctl, &atm_apr, &atm_clim, 0);
340 for (
int ip = 0; ip < atm_apr.np; ip++) {
341 atm_apr.time[ip] = obs_meas.time[0];
342 atm_apr.lon[ip] = obs_meas.vplon[0];
343 atm_apr.lat[ip] = obs_meas.vplat[0];
347 init_l2(&ncd, track, xtrack, &ctl, &atm_apr);
351 &atm_apr, &atm_i, &chisq);
354 if (gsl_finite(chisq)) {
355 chisq_min = GSL_MIN(chisq_min, chisq);
356 chisq_max = GSL_MAX(chisq_max, chisq);
362 buffer_nc(&atm_i, chisq, &ncd, track, xtrack, np0, np1);
366 TIMER(
"retrieval", 3);
377 printf(
"chi^2: min= %g / mean= %g / max= %g / m= %d\n",
378 chisq_min, chisq_mean / m, chisq_max, m);
379 printf(
"Retrieval finished on rank %d of %d!\n", rank, size);
389 printf(
"MEMORY_ATM = %g MByte\n", 4. *
sizeof(atm_t) / 1024. / 1024.);
390 printf(
"MEMORY_CTL = %g MByte\n", 1. *
sizeof(ctl_t) / 1024. / 1024.);
391 printf(
"MEMORY_NCD = %g MByte\n", 1. *
sizeof(
ncd_t) / 1024. / 1024.);
392 printf(
"MEMORY_OBS = %g MByte\n", 3. *
sizeof(atm_t) / 1024. / 1024.);
393 printf(
"MEMORY_RET = %g MByte\n", 1. *
sizeof(ret_t) / 1024. / 1024.);
394 printf(
"MEMORY_TBL = %g MByte\n", 1. *
sizeof(tbl_t) / 1024. / 1024.);
397 printf(
"SIZE_TASKS = %d\n", size);
398 printf(
"SIZE_THREADS = %d\n", omp_get_max_threads());
412 const char *longname,
419 if (nc_inq_varid(ncid, varname, varid) != NC_NOERR) {
422 NC(nc_def_var(ncid, varname, type, ndims, dimid, varid));
426 (ncid, *varid,
"long_name", strlen(longname), longname));
429 NC(nc_put_att_text(ncid, *varid,
"units", strlen(unit), unit));
445 ncd->
np = np1 - np0 + 1;
448 for (
int ip = np0; ip <= np1; ip++) {
449 ncd->
ret_z[ip - np0] = (float) atm->z[ip];
452 (gsl_finite(chisq) ? (float) atm->t[ip] : GSL_NAN);
466 for (
int lay = 0; lay <
L2_NLAY; lay++) {
469 for (
int track = 0; track <
L2_NTRACK; track++)
470 for (
int xtrack = 0; xtrack <
L2_NXTRACK; xtrack++) {
473 help[track][xtrack] = 0;
477 for (
int track2 = 0; track2 <
L2_NTRACK; track2++)
478 for (
int xtrack2 = 0; xtrack2 <
L2_NXTRACK; xtrack2++)
479 if (gsl_finite(x[track2][xtrack2][lay])
480 && x[track2][xtrack2][lay] > 0) {
481 const double w = exp(-gsl_pow_2((xtrack - xtrack2) / cx)
482 - gsl_pow_2((track - track2) / cy));
483 help[track][xtrack] += w * x[track2][xtrack2][lay];
489 help[track][xtrack] /= wsum;
491 help[track][xtrack] = GSL_NAN;
495 for (
int track = 0; track <
L2_NTRACK; track++)
496 for (
int xtrack = 0; xtrack <
L2_NXTRACK; xtrack++)
497 x[track][xtrack][lay] = help[track][xtrack];
510 static atm_t atm_airs;
512 double k[NW], p, q[NG], t, zmax = 0, zmin = 1000;
520 for (
int lay = 0; lay <
L2_NLAY; lay++)
521 if (gsl_finite(ncd->
l2_z[track][xtrack][lay])) {
522 atm_airs.z[atm_airs.np] = ncd->
l2_z[track][xtrack][lay];
523 atm_airs.p[atm_airs.np] = ncd->
l2_p[lay];
524 atm_airs.t[atm_airs.np] = ncd->
l2_t[track][xtrack][lay];
525 if ((++atm_airs.np) > NP)
526 ERRMSG(
"Too many layers!");
530 if (atm_airs.np <= 0)
534 for (
int ip = 0; ip < atm_airs.np; ip++) {
535 zmax = GSL_MAX(zmax, atm_airs.z[ip]);
536 zmin = GSL_MIN(zmin, atm_airs.z[ip]);
540 for (
int ip = 0; ip < atm->np; ip++) {
543 intpol_atm(ctl, &atm_airs, atm->z[ip], &p, &t, q, k);
547 if (atm->z[ip] > zmax)
548 w = GSL_MAX(1 - (atm->z[ip] - zmax) / 50, 0);
549 if (atm->z[ip] < zmin)
550 w = GSL_MAX(1 - (zmin - atm->z[ip]) / 50, 0);
553 atm->t[ip] = w * t + (1 - w) * atm->t[ip];
554 atm->p[ip] = w * p + (1 - w) * atm->p[ip];
570 static int ipa[N], iqa[N];
572 double disq = 0, lmpar = 0.001;
579 const size_t m = obs2y(ctl, obs_meas, NULL, NULL, NULL);
580 const size_t n = atm2x(ctl, atm_apr, NULL, iqa, ipa);
581 if (m == 0 || n == 0) {
587 gsl_matrix *a = gsl_matrix_alloc(n, n);
588 gsl_matrix *cov = gsl_matrix_alloc(n, n);
589 gsl_matrix *k_i = gsl_matrix_alloc(m, n);
590 gsl_matrix *s_a_inv = gsl_matrix_alloc(n, n);
592 gsl_vector *b = gsl_vector_alloc(n);
593 gsl_vector *dx = gsl_vector_alloc(n);
594 gsl_vector *dy = gsl_vector_alloc(m);
595 gsl_vector *sig_eps_inv = gsl_vector_alloc(m);
596 gsl_vector *sig_formod = gsl_vector_alloc(m);
597 gsl_vector *sig_noise = gsl_vector_alloc(m);
598 gsl_vector *x_a = gsl_vector_alloc(n);
599 gsl_vector *x_i = gsl_vector_alloc(n);
600 gsl_vector *x_step = gsl_vector_alloc(n);
601 gsl_vector *y_aux = gsl_vector_alloc(m);
602 gsl_vector *y_i = gsl_vector_alloc(m);
603 gsl_vector *y_m = gsl_vector_alloc(m);
606 copy_atm(ctl, atm_i, atm_apr, 0);
607 copy_obs(ctl, obs_i, obs_meas, 0);
608 formod(ctl, tbl, atm_i, obs_i);
611 atm2x(ctl, atm_apr, x_a, NULL, NULL);
612 atm2x(ctl, atm_i, x_i, NULL, NULL);
613 obs2y(ctl, obs_meas, y_m, NULL, NULL);
614 obs2y(ctl, obs_i, y_i, NULL, NULL);
617 set_cov_apr(ret, ctl, atm_apr, iqa, ipa, s_a_inv);
618 matrix_invert(s_a_inv);
621 set_cov_meas(ret, ctl, obs_meas, sig_noise, sig_formod, sig_eps_inv);
624 gsl_vector_memcpy(dx, x_i);
625 gsl_vector_sub(dx, x_a);
626 gsl_vector_memcpy(dy, y_m);
627 gsl_vector_sub(dy, y_i);
630 *chisq = cost_function(dx, dy, s_a_inv, sig_eps_inv);
633 kernel(ctl, tbl, atm_i, obs_i, k_i);
640 for (
int it = 1; it <= ret->conv_itmax; it++) {
643 double chisq_old = *chisq;
646 if (it > 1 && it % ret->kernel_recomp == 0)
647 kernel(ctl, tbl, atm_i, obs_i, k_i);
650 if (it == 1 || it % ret->kernel_recomp == 0)
651 matrix_product(k_i, sig_eps_inv, 1, cov);
654 for (
size_t i = 0; i < m; i++)
655 gsl_vector_set(y_aux, i, gsl_vector_get(dy, i)
656 * gsl_pow_2(gsl_vector_get(sig_eps_inv, i)));
657 gsl_blas_dgemv(CblasTrans, 1.0, k_i, y_aux, 0.0, b);
658 gsl_blas_dgemv(CblasNoTrans, -1.0, s_a_inv, dx, 1.0, b);
661 for (
int it2 = 0; it2 < 20; it2++) {
664 gsl_matrix_memcpy(a, s_a_inv);
665 gsl_matrix_scale(a, 1 + lmpar);
666 gsl_matrix_add(a, cov);
669 gsl_linalg_cholesky_decomp(a);
670 gsl_linalg_cholesky_solve(a, b, x_step);
673 gsl_vector_add(x_i, x_step);
674 copy_atm(ctl, atm_i, atm_apr, 0);
675 copy_obs(ctl, obs_i, obs_meas, 0);
676 x2atm(ctl, x_i, atm_i);
679 for (
int ip = 0; ip < atm_i->np; ip++) {
680 atm_i->p[ip] = GSL_MIN(GSL_MAX(atm_i->p[ip], 5e-7), 5e4);
681 atm_i->t[ip] = GSL_MIN(GSL_MAX(atm_i->t[ip], 100), 400);
682 for (
int ig = 0; ig < ctl->ng; ig++)
683 atm_i->q[ig][ip] = GSL_MIN(GSL_MAX(atm_i->q[ig][ip], 0), 1);
684 for (
int iw = 0; iw < ctl->nw; iw++)
685 atm_i->k[iw][ip] = GSL_MAX(atm_i->k[iw][ip], 0);
689 formod(ctl, tbl, atm_i, obs_i);
690 obs2y(ctl, obs_i, y_i, NULL, NULL);
693 gsl_vector_memcpy(dx, x_i);
694 gsl_vector_sub(dx, x_a);
695 gsl_vector_memcpy(dy, y_m);
696 gsl_vector_sub(dy, y_i);
699 *chisq = cost_function(dx, dy, s_a_inv, sig_eps_inv);
702 if (*chisq > chisq_old) {
704 gsl_vector_sub(x_i, x_step);
712 gsl_blas_ddot(x_step, b, &disq);
716 if ((it == 1 || it % ret->kernel_recomp == 0) && disq < ret->conv_dmin)
725 gsl_matrix_free(cov);
726 gsl_matrix_free(k_i);
727 gsl_matrix_free(s_a_inv);
732 gsl_vector_free(sig_eps_inv);
733 gsl_vector_free(sig_formod);
734 gsl_vector_free(sig_noise);
735 gsl_vector_free(x_a);
736 gsl_vector_free(x_i);
737 gsl_vector_free(x_step);
738 gsl_vector_free(y_aux);
739 gsl_vector_free(y_i);
740 gsl_vector_free(y_m);
752 printf(
"Read netCDF file: %s\n", filename);
753 NC(nc_open(filename, NC_WRITE, &ncd->
ncid));
756 NC(nc_inq_varid(ncd->
ncid,
"l1_time", &varid));
758 NC(nc_inq_varid(ncd->
ncid,
"l1_lon", &varid));
759 NC(nc_get_var_double(ncd->
ncid, varid, ncd->
l1_lon[0]));
760 NC(nc_inq_varid(ncd->
ncid,
"l1_lat", &varid));
761 NC(nc_get_var_double(ncd->
ncid, varid, ncd->
l1_lat[0]));
762 NC(nc_inq_varid(ncd->
ncid,
"l1_sat_z", &varid));
764 NC(nc_inq_varid(ncd->
ncid,
"l1_sat_lon", &varid));
766 NC(nc_inq_varid(ncd->
ncid,
"l1_sat_lat", &varid));
768 NC(nc_inq_varid(ncd->
ncid,
"l1_nu", &varid));
769 NC(nc_get_var_double(ncd->
ncid, varid, ncd->
l1_nu));
770 NC(nc_inq_varid(ncd->
ncid,
"l1_rad", &varid));
771 NC(nc_get_var_float(ncd->
ncid, varid, ncd->
l1_rad[0][0]));
774 NC(nc_inq_varid(ncd->
ncid,
"l2_z", &varid));
775 NC(nc_get_var_double(ncd->
ncid, varid, ncd->
l2_z[0][0]));
776 NC(nc_inq_varid(ncd->
ncid,
"l2_press", &varid));
777 NC(nc_get_var_double(ncd->
ncid, varid, ncd->
l2_p));
778 NC(nc_inq_varid(ncd->
ncid,
"l2_temp", &varid));
779 NC(nc_get_var_double(ncd->
ncid, varid, ncd->
l2_t[0][0]));
788 int dimid[10], p_id, t_id, z_id;
791 printf(
"Write netCDF file: %s\n", filename);
794 NC(nc_inq_dimid(ncd->
ncid,
"L1_NTRACK", &dimid[0]));
795 NC(nc_inq_dimid(ncd->
ncid,
"L1_NXTRACK", &dimid[1]));
801 if (nc_inq_dimid(ncd->
ncid,
"RET_NP", &dimid[2]) != NC_NOERR)
802 NC(nc_def_dim(ncd->
ncid,
"RET_NP", (
size_t) ncd->
np, &dimid[2]));
805 add_var(ncd->
ncid,
"ret_z",
"km",
"altitude", NC_FLOAT, &dimid[2], &z_id,
807 add_var(ncd->
ncid,
"ret_press",
"hPa",
"pressure", NC_FLOAT, dimid, &p_id,
809 add_var(ncd->
ncid,
"ret_temp",
"K",
"temperature", NC_FLOAT, dimid, &t_id,
int main(int argc, char *argv[])
#define NC(cmd)
Execute netCDF library command and check result.
void fill_gaps(double x[L2_NTRACK][L2_NXTRACK][L2_NLAY], double cx, double cy)
Fill data gaps in L2 data.
#define L1_NXTRACK
Across-track size of AIRS radiance granule (don't change).
#define L1_NTRACK
Along-track size of AIRS radiance granule (don't change).
#define L2_NXTRACK
Across-track size of AIRS retrieval granule (don't change).
void write_nc(char *filename, ncd_t *ncd)
Write to netCDF file...
void read_nc(char *filename, ncd_t *ncd)
Read netCDF file.
void add_var(int ncid, const char *varname, const char *unit, const char *longname, int type, int dimid[], int *varid, int ndims)
Create variable in netCDF file.
#define L2_NLAY
Number of AIRS pressure layers (don't change).
void 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.
void buffer_nc(atm_t *atm, double chisq, ncd_t *ncd, int track, int xtrack, int np0, int np1)
Buffer netCDF data.
#define L1_NCHAN
Number of AIRS radiance channels (don't change).
#define L2_NTRACK
Along-track size of AIRS retrieval granule (don't change).
void init_l2(ncd_t *ncd, int track, int xtrack, ctl_t *ctl, atm_t *atm)
Initialize with AIRS Level-2 data.
double l1_lat[L1_NTRACK][L1_NXTRACK]
Footprint latitude [deg].
float ret_p[L1_NTRACK *L1_NXTRACK]
Pressure [hPa].
double l1_lon[L1_NTRACK][L1_NXTRACK]
Footprint longitude [deg].
double l2_z[L2_NTRACK][L2_NXTRACK][L2_NLAY]
Altitude [km].
double l1_sat_lat[L1_NTRACK]
Satellite latitude [deg].
int np
Number of retrieval altitudes.
double l1_sat_lon[L1_NTRACK]
Satellite longitude [deg].
float ret_z[NP]
Altitude [km].
double l2_p[L2_NLAY]
Pressure [hPa].
float ret_t[L1_NTRACK *L1_NXTRACK *NP]
Temperature [K].
double l2_t[L2_NTRACK][L2_NXTRACK][L2_NLAY]
Temperature [K].
double l1_nu[L1_NCHAN]
Channel frequencies [cm^-1].
double l1_sat_z[L1_NTRACK]
Satellite altitude [km].
float l1_rad[L1_NTRACK][L1_NXTRACK][L1_NCHAN]
Radiance [W/(m^2 sr cm^-1)].
double l1_time[L1_NTRACK][L1_NXTRACK]
Time (seconds since 2000-01-01T00:00Z).