37 int nc_result=(cmd); \
38 if(nc_result!=NC_NOERR) \
39 ERRMSG("%s", nc_strerror(nc_result)); \
130 const char *longname,
179 static atm_t atm_apr, atm_clim, atm_i;
180 static obs_t obs_i, obs_meas;
188 double chisq, chisq_min, chisq_max, chisq_mean, z[NP];
190 int channel[ND], m, ntask = -1, rank, size;
197 MPI_Init(&argc, &argv);
198 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
199 MPI_Comm_size(MPI_COMM_WORLD, &size);
206 ERRMSG(
"Give parameters: <ctl> <filelist>");
209 read_ctl(argc, argv, &ctl);
210 read_ret(argc, argv, &ctl, &ret);
213 tbl_t *tbl = read_tbl(&ctl);
216 const int nz = (int) scan_ctl(argc, argv,
"NZ", -1,
"", NULL);
218 ERRMSG(
"Too many altitudes!");
219 for (
int iz = 0; iz < nz; iz++)
220 z[iz] = scan_ctl(argc, argv,
"Z", iz,
"", NULL);
223 const int track0 = (int) scan_ctl(argc, argv,
"TRACK_MIN", -1,
"0", NULL);
224 const int track1 = (int) scan_ctl(argc, argv,
"TRACK_MAX", -1,
"134", NULL);
227 const int xtrack0 = (int) scan_ctl(argc, argv,
"XTRACK_MIN", -1,
"0", NULL);
229 (int) scan_ctl(argc, argv,
"XTRACK_MAX", -1,
"89", NULL);
232 int np0 = (int) scan_ctl(argc, argv,
"NP_MIN", -1,
"0", NULL);
233 int np1 = (int) scan_ctl(argc, argv,
"NP_MAX", -1,
"100", NULL);
234 np1 = GSL_MIN(np1, nz - 1);
237 const double sx = scan_ctl(argc, argv,
"SX", -1,
"8", NULL);
238 const double sy = scan_ctl(argc, argv,
"SY", -1,
"2", NULL);
241 const double sza_thresh = scan_ctl(argc, argv,
"SZA", -1,
"96", NULL);
248 printf(
"Read filelist: %s\n", argv[2]);
249 if (!(in = fopen(argv[2],
"r")))
250 ERRMSG(
"Cannot open filelist!");
253 while (fscanf(in,
"%s", filename) != EOF) {
256 if ((++ntask) % size != rank)
260 printf(
"Retrieve file %s on rank %d of %d (with %d threads)...\n",
261 filename, rank + 1, size, omp_get_max_threads());
271 for (
int id = 0;
id < ctl.nd;
id++) {
274 if (fabs(ctl.nu[
id] - ncd.
l1_nu[i]) < 0.1)
277 ERRMSG(
"Cannot identify radiance channel!");
286 for (
int iz = 0; iz < nz; iz++)
287 atm_clim.z[iz] = z[iz];
288 climatology(&ctl, &atm_clim);
301 for (
int track = track0; track <= track1; track++) {
304 TIMER(
"retrieval", 1);
307 for (
int xtrack = xtrack0; xtrack <= xtrack1; xtrack++) {
311 obs_meas.time[0] = ncd.
l1_time[track][xtrack];
312 obs_meas.obsz[0] = ncd.
l1_sat_z[track];
315 obs_meas.vplon[0] = ncd.
l1_lon[track][xtrack];
316 obs_meas.vplat[0] = ncd.
l1_lat[track][xtrack];
317 for (
int id = 0;
id < ctl.nd;
id++)
318 obs_meas.rad[
id][0] = ncd.
l1_rad[track][xtrack][channel[
id]];
321 if (RAD2DEG(acos(cos_sza(obs_meas.time[0], obs_meas.obslon[0],
322 obs_meas.obslat[0]))) < sza_thresh)
323 for (
int id = 0;
id < ctl.nd;
id++)
324 if (ctl.nu[
id] >= 2000)
325 obs_meas.rad[id][0] = GSL_NAN;
328 copy_atm(&ctl, &atm_apr, &atm_clim, 0);
329 for (
int ip = 0; ip < atm_apr.np; ip++) {
330 atm_apr.time[ip] = obs_meas.time[0];
331 atm_apr.lon[ip] = obs_meas.vplon[0];
332 atm_apr.lat[ip] = obs_meas.vplat[0];
336 init_l2(&ncd, track, xtrack, &ctl, &atm_apr);
339 optimal_estimation(&ret, &ctl, tbl, &obs_meas, &obs_i,
340 &atm_apr, &atm_i, &chisq);
343 if (gsl_finite(chisq)) {
344 chisq_min = GSL_MIN(chisq_min, chisq);
345 chisq_max = GSL_MAX(chisq_max, chisq);
351 buffer_nc(&atm_i, chisq, &ncd, track, xtrack, np0, np1);
355 TIMER(
"retrieval", 3);
366 printf(
"chi^2: min= %g / mean= %g / max= %g / m= %d\n",
367 chisq_min, chisq_mean / m, chisq_max, m);
368 printf(
"Retrieval finished on rank %d of %d!\n", rank, size);
378 printf(
"MEMORY_ATM = %g MByte\n", 4. *
sizeof(atm_t) / 1024. / 1024.);
379 printf(
"MEMORY_CTL = %g MByte\n", 1. *
sizeof(ctl_t) / 1024. / 1024.);
380 printf(
"MEMORY_NCD = %g MByte\n", 1. *
sizeof(
ncd_t) / 1024. / 1024.);
381 printf(
"MEMORY_OBS = %g MByte\n", 3. *
sizeof(atm_t) / 1024. / 1024.);
382 printf(
"MEMORY_RET = %g MByte\n", 1. *
sizeof(ret_t) / 1024. / 1024.);
383 printf(
"MEMORY_TBL = %g MByte\n", 1. *
sizeof(tbl_t) / 1024. / 1024.);
386 printf(
"SIZE_TASKS = %d\n", size);
387 printf(
"SIZE_THREADS = %d\n", omp_get_max_threads());
401 const char *longname,
408 if (nc_inq_varid(ncid, varname, varid) != NC_NOERR) {
411 NC(nc_def_var(ncid, varname, type, ndims, dimid, varid));
415 (ncid, *varid,
"long_name", strlen(longname), longname));
418 NC(nc_put_att_text(ncid, *varid,
"units", strlen(unit), unit));
434 ncd->
np = np1 - np0 + 1;
437 for (
int ip = np0; ip <= np1; ip++) {
438 ncd->
ret_z[ip - np0] = (float) atm->z[ip];
441 (gsl_finite(chisq) ? (float) atm->t[ip] : GSL_NAN);
455 for (
int lay = 0; lay <
L2_NLAY; lay++) {
458 for (
int track = 0; track <
L2_NTRACK; track++)
459 for (
int xtrack = 0; xtrack <
L2_NXTRACK; xtrack++) {
462 help[track][xtrack] = 0;
466 for (
int track2 = 0; track2 <
L2_NTRACK; track2++)
467 for (
int xtrack2 = 0; xtrack2 <
L2_NXTRACK; xtrack2++)
468 if (gsl_finite(x[track2][xtrack2][lay])
469 && x[track2][xtrack2][lay] > 0) {
470 const double w = exp(-gsl_pow_2((xtrack - xtrack2) / cx)
471 - gsl_pow_2((track - track2) / cy));
472 help[track][xtrack] += w * x[track2][xtrack2][lay];
478 help[track][xtrack] /= wsum;
480 help[track][xtrack] = GSL_NAN;
484 for (
int track = 0; track <
L2_NTRACK; track++)
485 for (
int xtrack = 0; xtrack <
L2_NXTRACK; xtrack++)
486 x[track][xtrack][lay] = help[track][xtrack];
499 static atm_t atm_airs;
501 double k[NW], p, q[NG], t, zmax = 0, zmin = 1000;
509 for (
int lay = 0; lay <
L2_NLAY; lay++)
510 if (gsl_finite(ncd->
l2_z[track][xtrack][lay])) {
511 atm_airs.z[atm_airs.np] = ncd->
l2_z[track][xtrack][lay];
512 atm_airs.p[atm_airs.np] = ncd->
l2_p[lay];
513 atm_airs.t[atm_airs.np] = ncd->
l2_t[track][xtrack][lay];
514 if ((++atm_airs.np) > NP)
515 ERRMSG(
"Too many layers!");
519 if (atm_airs.np <= 0)
523 for (
int ip = 0; ip < atm_airs.np; ip++) {
524 zmax = GSL_MAX(zmax, atm_airs.z[ip]);
525 zmin = GSL_MIN(zmin, atm_airs.z[ip]);
529 for (
int ip = 0; ip < atm->np; ip++) {
532 intpol_atm(ctl, &atm_airs, atm->z[ip], &p, &t, q, k);
536 if (atm->z[ip] > zmax)
537 w = GSL_MAX(1 - (atm->z[ip] - zmax) / 50, 0);
538 if (atm->z[ip] < zmin)
539 w = GSL_MAX(1 - (zmin - atm->z[ip]) / 50, 0);
542 atm->t[ip] = w * t + (1 - w) * atm->t[ip];
543 atm->p[ip] = w * p + (1 - w) * atm->p[ip];
556 printf(
"Read netCDF file: %s\n", filename);
557 NC(nc_open(filename, NC_WRITE, &ncd->
ncid));
560 NC(nc_inq_varid(ncd->
ncid,
"l1_time", &varid));
562 NC(nc_inq_varid(ncd->
ncid,
"l1_lon", &varid));
563 NC(nc_get_var_double(ncd->
ncid, varid, ncd->
l1_lon[0]));
564 NC(nc_inq_varid(ncd->
ncid,
"l1_lat", &varid));
565 NC(nc_get_var_double(ncd->
ncid, varid, ncd->
l1_lat[0]));
566 NC(nc_inq_varid(ncd->
ncid,
"l1_sat_z", &varid));
568 NC(nc_inq_varid(ncd->
ncid,
"l1_sat_lon", &varid));
570 NC(nc_inq_varid(ncd->
ncid,
"l1_sat_lat", &varid));
572 NC(nc_inq_varid(ncd->
ncid,
"l1_nu", &varid));
573 NC(nc_get_var_double(ncd->
ncid, varid, ncd->
l1_nu));
574 NC(nc_inq_varid(ncd->
ncid,
"l1_rad", &varid));
575 NC(nc_get_var_float(ncd->
ncid, varid, ncd->
l1_rad[0][0]));
578 NC(nc_inq_varid(ncd->
ncid,
"l2_z", &varid));
579 NC(nc_get_var_double(ncd->
ncid, varid, ncd->
l2_z[0][0]));
580 NC(nc_inq_varid(ncd->
ncid,
"l2_press", &varid));
581 NC(nc_get_var_double(ncd->
ncid, varid, ncd->
l2_p));
582 NC(nc_inq_varid(ncd->
ncid,
"l2_temp", &varid));
583 NC(nc_get_var_double(ncd->
ncid, varid, ncd->
l2_t[0][0]));
592 int dimid[10], p_id, t_id, z_id;
595 printf(
"Write netCDF file: %s\n", filename);
598 NC(nc_inq_dimid(ncd->
ncid,
"L1_NTRACK", &dimid[0]));
599 NC(nc_inq_dimid(ncd->
ncid,
"L1_NXTRACK", &dimid[1]));
605 if (nc_inq_dimid(ncd->
ncid,
"RET_NP", &dimid[2]) != NC_NOERR)
606 NC(nc_def_dim(ncd->
ncid,
"RET_NP", (
size_t) ncd->
np, &dimid[2]));
609 add_var(ncd->
ncid,
"ret_z",
"km",
"altitude", NC_FLOAT, &dimid[2], &z_id,
611 add_var(ncd->
ncid,
"ret_press",
"hPa",
"pressure", NC_FLOAT, dimid, &p_id,
613 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 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).