37 int nc_result=(cmd); \
38 if(nc_result!=NC_NOERR) \
39 ERRMSG("%s", nc_strerror(nc_result)); \
130 const char *longname,
175 static atm_t atm_apr, atm_clim, atm_i;
176 static obs_t obs_i, obs_meas;
180 static FILE *in, *out;
182 static char filename[LEN], filename2[2 * LEN];
187 static int channel[ND], n, ntask = -1, rank, size;
194 MPI_Init(&argc, &argv);
195 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
196 MPI_Comm_size(MPI_COMM_WORLD, &size);
203 ERRMSG(
"Give parameters: <ctl> <filelist>");
206 read_ctl(argc, argv, &ctl);
207 read_ret(argc, argv, &ctl, &ret);
210 tbl_t *tbl = read_tbl(&ctl);
213 const int nz = (int) scan_ctl(argc, argv,
"NZ", -1,
"", NULL);
215 ERRMSG(
"Too many altitudes!");
216 for (
int iz = 0; iz < nz; iz++)
217 z[iz] = scan_ctl(argc, argv,
"Z", iz,
"", NULL);
220 const int track0 = (int) scan_ctl(argc, argv,
"TRACK_MIN", -1,
"0", NULL);
221 const int track1 = (int) scan_ctl(argc, argv,
"TRACK_MAX", -1,
"134", NULL);
224 const int xtrack0 = (int) scan_ctl(argc, argv,
"XTRACK_MIN", -1,
"0", NULL);
226 (int) scan_ctl(argc, argv,
"XTRACK_MAX", -1,
"89", NULL);
229 const double sx = scan_ctl(argc, argv,
"SX", -1,
"8", NULL);
230 const double sy = scan_ctl(argc, argv,
"SY", -1,
"2", NULL);
237 printf(
"Read filelist: %s\n", argv[2]);
238 if (!(in = fopen(argv[2],
"r")))
239 ERRMSG(
"Cannot open filelist!");
242 while (fscanf(in,
"%s", filename) != EOF) {
245 if ((++ntask) % size != rank)
249 printf(
"Retrieve file %s on rank %d of %d (with %d threads)...\n",
250 filename, rank + 1, size, omp_get_max_threads());
260 for (
int id = 0;
id < ctl.nd;
id++) {
263 if (fabs(ctl.nu[
id] - ncd.
l1_nu[i]) < 0.1)
266 ERRMSG(
"Cannot identify radiance channel!");
275 for (
int iz = 0; iz < nz; iz++)
276 atm_clim.z[iz] = z[iz];
277 climatology(&ctl, &atm_clim);
284 for (
int track = track0; track <= track1; track++) {
287 TIMER(
"retrieval", 1);
290 for (
int xtrack = xtrack0; xtrack <= xtrack1; xtrack++) {
294 obs_meas.time[0] = ncd.
l1_time[track][xtrack];
295 obs_meas.obsz[0] = ncd.
l1_sat_z[track];
298 obs_meas.vplon[0] = ncd.
l1_lon[track][xtrack];
299 obs_meas.vplat[0] = ncd.
l1_lat[track][xtrack];
300 for (
int id = 0;
id < ctl.nd;
id++)
301 obs_meas.rad[
id][0] = ncd.
l1_rad[track][xtrack][channel[
id]];
304 for (
int id = 0;
id < ctl.nd;
id++)
305 if (ctl.nu[
id] >= 2000)
306 obs_meas.rad[id][0] = GSL_NAN;
309 copy_atm(&ctl, &atm_apr, &atm_clim, 0);
310 for (
int ip = 0; ip < atm_apr.np; ip++) {
311 atm_apr.time[ip] = obs_meas.time[0];
312 atm_apr.lon[ip] = obs_meas.vplon[0];
313 atm_apr.lat[ip] = obs_meas.vplat[0];
317 init_l2(&ncd, track, xtrack, &ctl, &atm_apr);
321 &atm_apr, &atm_i, &chisq[track][xtrack]);
324 for (
int id = 0;
id < ctl.nd;
id++)
325 obs_meas.rad[
id][0] = obs_i.rad[
id][0]
326 = ncd.
l1_rad[track][xtrack][channel[
id]];
327 formod(&ctl, tbl, &atm_i, &obs_i);
331 ni[track][xtrack] = 0;
332 for (
int id = 0;
id < ctl.nd;
id++)
333 if (ctl.nu[
id] >= 2000 && gsl_finite(obs_meas.rad[
id][0])) {
335 (BRIGHT(obs_meas.rad[
id][0], ncd.
l1_nu[channel[
id]])
336 - BRIGHT(obs_i.rad[
id][0], ncd.
l1_nu[channel[
id]]));
339 ni[track][xtrack] /= n;
343 TIMER(
"retrieval", 3);
351 sprintf(filename2,
"%s.nlte.tab", filename);
354 printf(
"Write non-LTE data: %s\n", filename2);
355 if (!(out = fopen(filename2,
"w")))
356 ERRMSG(
"Cannot create file!");
360 "# $1 = time (seconds since 2000-01-01, 00:00 UTC)\n"
361 "# $2 = longitude [deg]\n"
362 "# $3 = latitude [deg]\n"
363 "# $4 = solar zenith angle [deg]\n"
364 "# $5 = non-LTE index [K]\n" "# $6 = chi^2 of retrieval fit\n");
367 for (
int track = track0; track <= track1; track++) {
369 for (
int xtrack = xtrack0; xtrack <= xtrack1; xtrack++)
370 fprintf(out,
"%.2f %g %g %g %g %g\n",
372 ncd.
l1_lon[track][xtrack],
373 ncd.
l1_lat[track][xtrack],
375 ncd.
l1_lat[track][xtrack]), ni[track][xtrack],
376 chisq[track][xtrack]);
383 printf(
"Retrieval finished on rank %d of %d!\n", rank, size);
396 printf(
"MEMORY_ATM = %g MByte\n", 4. *
sizeof(atm_t) / 1024. / 1024.);
397 printf(
"MEMORY_CTL = %g MByte\n", 1. *
sizeof(ctl_t) / 1024. / 1024.);
398 printf(
"MEMORY_NCD = %g MByte\n", 1. *
sizeof(
ncd_t) / 1024. / 1024.);
399 printf(
"MEMORY_OBS = %g MByte\n", 3. *
sizeof(atm_t) / 1024. / 1024.);
400 printf(
"MEMORY_RET = %g MByte\n", 1. *
sizeof(ret_t) / 1024. / 1024.);
401 printf(
"MEMORY_TBL = %g MByte\n", 1. *
sizeof(tbl_t) / 1024. / 1024.);
404 printf(
"SIZE_TASKS = %d\n", size);
405 printf(
"SIZE_THREADS = %d\n", omp_get_max_threads());
419 const char *longname,
426 if (nc_inq_varid(ncid, varname, varid) != NC_NOERR) {
429 NC(nc_def_var(ncid, varname, type, ndims, dimid, varid));
433 (ncid, *varid,
"long_name", strlen(longname), longname));
436 NC(nc_put_att_text(ncid, *varid,
"units", strlen(unit), unit));
450 for (
int lay = 0; lay <
L2_NLAY; lay++) {
453 for (
int track = 0; track <
L2_NTRACK; track++)
454 for (
int xtrack = 0; xtrack <
L2_NXTRACK; xtrack++) {
457 help[track][xtrack] = 0;
461 for (
int track2 = 0; track2 <
L2_NTRACK; track2++)
462 for (
int xtrack2 = 0; xtrack2 <
L2_NXTRACK; xtrack2++)
463 if (gsl_finite(x[track2][xtrack2][lay])
464 && x[track2][xtrack2][lay] > 0) {
465 const double w = exp(-gsl_pow_2((xtrack - xtrack2) / cx)
466 - gsl_pow_2((track - track2) / cy));
467 help[track][xtrack] += w * x[track2][xtrack2][lay];
473 help[track][xtrack] /= wsum;
475 help[track][xtrack] = GSL_NAN;
479 for (
int track = 0; track <
L2_NTRACK; track++)
480 for (
int xtrack = 0; xtrack <
L2_NXTRACK; xtrack++)
481 x[track][xtrack][lay] = help[track][xtrack];
494 static atm_t atm_airs;
496 double k[NW], p, q[NG], t, zmax = 0, zmin = 1000;
504 for (
int lay = 0; lay <
L2_NLAY; lay++)
505 if (gsl_finite(ncd->
l2_z[track][xtrack][lay])) {
506 atm_airs.z[atm_airs.np] = ncd->
l2_z[track][xtrack][lay];
507 atm_airs.p[atm_airs.np] = ncd->
l2_p[lay];
508 atm_airs.t[atm_airs.np] = ncd->
l2_t[track][xtrack][lay];
509 if ((++atm_airs.np) > NP)
510 ERRMSG(
"Too many layers!");
514 if (atm_airs.np <= 0)
518 for (
int ip = 0; ip < atm_airs.np; ip++) {
519 zmax = GSL_MAX(zmax, atm_airs.z[ip]);
520 zmin = GSL_MIN(zmin, atm_airs.z[ip]);
524 for (
int ip = 0; ip < atm->np; ip++) {
527 intpol_atm(ctl, &atm_airs, atm->z[ip], &p, &t, q, k);
531 if (atm->z[ip] > zmax)
532 w = GSL_MAX(1 - (atm->z[ip] - zmax) / 50, 0);
533 if (atm->z[ip] < zmin)
534 w = GSL_MAX(1 - (zmin - atm->z[ip]) / 50, 0);
537 atm->t[ip] = w * t + (1 - w) * atm->t[ip];
538 atm->p[ip] = w * p + (1 - w) * atm->p[ip];
554 static int ipa[N], iqa[N];
556 double disq = 0, lmpar = 0.001;
563 const size_t m = obs2y(ctl, obs_meas, NULL, NULL, NULL);
564 const size_t n = atm2x(ctl, atm_apr, NULL, iqa, ipa);
565 if (m == 0 || n == 0) {
571 gsl_matrix *a = gsl_matrix_alloc(n, n);
572 gsl_matrix *cov = gsl_matrix_alloc(n, n);
573 gsl_matrix *k_i = gsl_matrix_alloc(m, n);
574 gsl_matrix *s_a_inv = gsl_matrix_alloc(n, n);
576 gsl_vector *b = gsl_vector_alloc(n);
577 gsl_vector *dx = gsl_vector_alloc(n);
578 gsl_vector *dy = gsl_vector_alloc(m);
579 gsl_vector *sig_eps_inv = gsl_vector_alloc(m);
580 gsl_vector *sig_formod = gsl_vector_alloc(m);
581 gsl_vector *sig_noise = gsl_vector_alloc(m);
582 gsl_vector *x_a = gsl_vector_alloc(n);
583 gsl_vector *x_i = gsl_vector_alloc(n);
584 gsl_vector *x_step = gsl_vector_alloc(n);
585 gsl_vector *y_aux = gsl_vector_alloc(m);
586 gsl_vector *y_i = gsl_vector_alloc(m);
587 gsl_vector *y_m = gsl_vector_alloc(m);
590 copy_atm(ctl, atm_i, atm_apr, 0);
591 copy_obs(ctl, obs_i, obs_meas, 0);
592 formod(ctl, tbl, atm_i, obs_i);
595 atm2x(ctl, atm_apr, x_a, NULL, NULL);
596 atm2x(ctl, atm_i, x_i, NULL, NULL);
597 obs2y(ctl, obs_meas, y_m, NULL, NULL);
598 obs2y(ctl, obs_i, y_i, NULL, NULL);
601 set_cov_apr(ret, ctl, atm_apr, iqa, ipa, s_a_inv);
602 matrix_invert(s_a_inv);
605 set_cov_meas(ret, ctl, obs_meas, sig_noise, sig_formod, sig_eps_inv);
608 gsl_vector_memcpy(dx, x_i);
609 gsl_vector_sub(dx, x_a);
610 gsl_vector_memcpy(dy, y_m);
611 gsl_vector_sub(dy, y_i);
614 *chisq = cost_function(dx, dy, s_a_inv, sig_eps_inv);
617 kernel(ctl, tbl, atm_i, obs_i, k_i);
624 for (
int it = 1; it <= ret->conv_itmax; it++) {
627 double chisq_old = *chisq;
630 if (it > 1 && it % ret->kernel_recomp == 0)
631 kernel(ctl, tbl, atm_i, obs_i, k_i);
634 if (it == 1 || it % ret->kernel_recomp == 0)
635 matrix_product(k_i, sig_eps_inv, 1, cov);
638 for (
size_t i = 0; i < m; i++)
639 gsl_vector_set(y_aux, i, gsl_vector_get(dy, i)
640 * gsl_pow_2(gsl_vector_get(sig_eps_inv, i)));
641 gsl_blas_dgemv(CblasTrans, 1.0, k_i, y_aux, 0.0, b);
642 gsl_blas_dgemv(CblasNoTrans, -1.0, s_a_inv, dx, 1.0, b);
645 for (
int it2 = 0; it2 < 20; it2++) {
648 gsl_matrix_memcpy(a, s_a_inv);
649 gsl_matrix_scale(a, 1 + lmpar);
650 gsl_matrix_add(a, cov);
653 gsl_linalg_cholesky_decomp(a);
654 gsl_linalg_cholesky_solve(a, b, x_step);
657 gsl_vector_add(x_i, x_step);
658 copy_atm(ctl, atm_i, atm_apr, 0);
659 copy_obs(ctl, obs_i, obs_meas, 0);
660 x2atm(ctl, x_i, atm_i);
663 for (
int ip = 0; ip < atm_i->np; ip++) {
664 atm_i->p[ip] = GSL_MIN(GSL_MAX(atm_i->p[ip], 5e-7), 5e4);
665 atm_i->t[ip] = GSL_MIN(GSL_MAX(atm_i->t[ip], 100), 400);
666 for (
int ig = 0; ig < ctl->ng; ig++)
667 atm_i->q[ig][ip] = GSL_MIN(GSL_MAX(atm_i->q[ig][ip], 0), 1);
668 for (
int iw = 0; iw < ctl->nw; iw++)
669 atm_i->k[iw][ip] = GSL_MAX(atm_i->k[iw][ip], 0);
673 formod(ctl, tbl, atm_i, obs_i);
674 obs2y(ctl, obs_i, y_i, NULL, NULL);
677 gsl_vector_memcpy(dx, x_i);
678 gsl_vector_sub(dx, x_a);
679 gsl_vector_memcpy(dy, y_m);
680 gsl_vector_sub(dy, y_i);
683 *chisq = cost_function(dx, dy, s_a_inv, sig_eps_inv);
686 if (*chisq > chisq_old) {
688 gsl_vector_sub(x_i, x_step);
696 gsl_blas_ddot(x_step, b, &disq);
700 if ((it == 1 || it % ret->kernel_recomp == 0) && disq < ret->conv_dmin)
709 gsl_matrix_free(cov);
710 gsl_matrix_free(k_i);
711 gsl_matrix_free(s_a_inv);
716 gsl_vector_free(sig_eps_inv);
717 gsl_vector_free(sig_formod);
718 gsl_vector_free(sig_noise);
719 gsl_vector_free(x_a);
720 gsl_vector_free(x_i);
721 gsl_vector_free(x_step);
722 gsl_vector_free(y_aux);
723 gsl_vector_free(y_i);
724 gsl_vector_free(y_m);
736 printf(
"Read netCDF file: %s\n", filename);
737 NC(nc_open(filename, NC_WRITE, &ncd->
ncid));
740 NC(nc_inq_varid(ncd->
ncid,
"l1_time", &varid));
742 NC(nc_inq_varid(ncd->
ncid,
"l1_lon", &varid));
743 NC(nc_get_var_double(ncd->
ncid, varid, ncd->
l1_lon[0]));
744 NC(nc_inq_varid(ncd->
ncid,
"l1_lat", &varid));
745 NC(nc_get_var_double(ncd->
ncid, varid, ncd->
l1_lat[0]));
746 NC(nc_inq_varid(ncd->
ncid,
"l1_sat_z", &varid));
748 NC(nc_inq_varid(ncd->
ncid,
"l1_sat_lon", &varid));
750 NC(nc_inq_varid(ncd->
ncid,
"l1_sat_lat", &varid));
752 NC(nc_inq_varid(ncd->
ncid,
"l1_nu", &varid));
753 NC(nc_get_var_double(ncd->
ncid, varid, ncd->
l1_nu));
754 NC(nc_inq_varid(ncd->
ncid,
"l1_rad", &varid));
755 NC(nc_get_var_float(ncd->
ncid, varid, ncd->
l1_rad[0][0]));
758 NC(nc_inq_varid(ncd->
ncid,
"l2_z", &varid));
759 NC(nc_get_var_double(ncd->
ncid, varid, ncd->
l2_z[0][0]));
760 NC(nc_inq_varid(ncd->
ncid,
"l2_press", &varid));
761 NC(nc_get_var_double(ncd->
ncid, varid, ncd->
l2_p));
762 NC(nc_inq_varid(ncd->
ncid,
"l2_temp", &varid));
763 NC(nc_get_var_double(ncd->
ncid, varid, ncd->
l2_t[0][0]));
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 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.
#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].
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].
double l1_sat_lon[L1_NTRACK]
Satellite longitude [deg].
double l2_p[L2_NLAY]
Pressure [hPa].
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).