37 int nc_result=(cmd); \
38 if(nc_result!=NC_NOERR) \
39 ERRMSG("%s", nc_strerror(nc_result)); \
136 const char *longname,
179 static atm_t atm_apr, atm_clim, atm_i;
180 static obs_t obs_i, obs_meas;
186 char filename[LEN], filename2[LEN];
188 double chisq, sza_thresh, z[NP], t0;
190 int channel[ND], i, id, ip, iz, nz, ntask = -1, rank, size,
191 np0, np1, track, track0, track1, xtrack, xtrack0, xtrack1,
199 MPI_Init(&argc, &argv);
200 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
201 MPI_Comm_size(MPI_COMM_WORLD, &size);
208 ERRMSG(
"Give parameters: <ctl> <filelist>");
211 read_ctl(argc, argv, &ctl);
212 read_ret(argc, argv, &ctl, &ret);
213 debug = (int) scan_ctl(argc, argv,
"DEBUG", -1,
"1", NULL);
216 tbl_t *tbl = read_tbl(&ctl);
219 nz = (int) scan_ctl(argc, argv,
"NZ", -1,
"", NULL);
221 ERRMSG(
"Too many altitudes!");
222 for (iz = 0; iz < nz; iz++)
223 z[iz] = scan_ctl(argc, argv,
"Z", iz,
"", NULL);
226 task0 = (int) scan_ctl(argc, argv,
"TASK_MIN", -1,
"0", NULL);
227 task1 = (int) scan_ctl(argc, argv,
"TASK_MAX", -1,
"99999", NULL);
230 track0 = (int) scan_ctl(argc, argv,
"TRACK_MIN", -1,
"0", NULL);
231 track1 = (int) scan_ctl(argc, argv,
"TRACK_MAX", -1,
"99999", NULL);
234 xtrack0 = (int) scan_ctl(argc, argv,
"XTRACK_MIN", -1,
"0", NULL);
235 xtrack1 = (int) scan_ctl(argc, argv,
"XTRACK_MAX", -1,
"59", NULL);
238 np0 = (int) scan_ctl(argc, argv,
"NP_MIN", -1,
"0", NULL);
239 np1 = (int) scan_ctl(argc, argv,
"NP_MAX", -1,
"100", NULL);
240 np1 = GSL_MIN(np1, nz - 1);
243 sza_thresh = scan_ctl(argc, argv,
"SZA", -1,
"96", NULL);
250 printf(
"Read filelist: %s\n", argv[2]);
251 if (!(in = fopen(argv[2],
"r")))
252 ERRMSG(
"Cannot open filelist!");
255 while (fscanf(in,
"%s", filename) != EOF) {
258 if ((++ntask) % size != rank)
262 if (ntask < task0 || ntask > task1)
266 printf(
"Retrieve file %s on rank %d of %d (with %d threads)...\n",
267 filename, rank + 1, size, omp_get_max_threads());
281 for (
id = 0;
id < ctl.nd;
id++) {
284 if (fabs(ctl.nu[
id] - ncd.
l1_nu[i]) < 0.1)
287 ERRMSG(
"Cannot identify radiance channel!");
292 for (iz = 0; iz < nz; iz++)
293 atm_clim.z[iz] = z[iz];
294 climatology(&ctl, &atm_clim);
301 for (track = track0; track <= track1; track++) {
304 for (xtrack = xtrack0; xtrack <= xtrack1; xtrack++) {
307 t0 = omp_get_wtime();
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 (
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 (
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 (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 buffer_nc(&atm_i, chisq, &ncd, track, xtrack, np0, np1);
348 (
" task= %4d | track= %5d | xtrack= %3d | chi^2= %8.3f | time= %8.3f s\n",
349 ntask, track, xtrack, chisq, omp_get_wtime() - t0);
351 sprintf(filename2,
"atm_apr_%d_%d_%d.tab", ntask, track, xtrack);
352 write_atm(NULL, filename2, &ctl, &atm_apr);
353 sprintf(filename2,
"atm_i_%d_%d_%d.tab", ntask, track, xtrack);
354 write_atm(NULL, filename2, &ctl, &atm_i);
355 sprintf(filename2,
"obs_meas_%d_%d_%d.tab", ntask, track, xtrack);
356 write_obs(NULL, filename2, &ctl, &obs_meas);
357 sprintf(filename2,
"obs_i_%d_%d_%d.tab", ntask, track, xtrack);
358 write_obs(NULL, filename2, &ctl, &obs_i);
371 printf(
"Retrieval finished on rank %d of %d!\n", rank, size);
381 printf(
"MEMORY_ATM = %g MByte\n", 4. *
sizeof(atm_t) / 1024. / 1024.);
382 printf(
"MEMORY_CTL = %g MByte\n", 1. *
sizeof(ctl_t) / 1024. / 1024.);
383 printf(
"MEMORY_NCD = %g MByte\n", 1. *
sizeof(
ncd_t) / 1024. / 1024.);
384 printf(
"MEMORY_OBS = %g MByte\n", 3. *
sizeof(atm_t) / 1024. / 1024.);
385 printf(
"MEMORY_RET = %g MByte\n", 1. *
sizeof(ret_t) / 1024. / 1024.);
386 printf(
"MEMORY_TBL = %g MByte\n", 1. *
sizeof(tbl_t) / 1024. / 1024.);
389 printf(
"SIZE_TASKS = %d\n", size);
390 printf(
"SIZE_THREADS = %d\n", omp_get_max_threads());
404 const char *longname,
411 if (nc_inq_varid(ncid, varname, varid) != NC_NOERR) {
414 NC(nc_def_var(ncid, varname, type, ndims, dimid, varid));
418 (ncid, *varid,
"long_name", strlen(longname), longname));
421 NC(nc_put_att_text(ncid, *varid,
"units", strlen(unit), unit));
439 ncd->
np = np1 - np0 + 1;
444 for (ip = np0; ip <= np1; ip++) {
445 ncd->
ret_z[ip - np0] = (float) atm->z[ip];
447 (gsl_finite(chisq) ? (float) atm->t[ip] : GSL_NAN);
460 static atm_t atm_iasi;
462 double k[NW], p, q[NG], t, w, zmax = 0, zmin = 1000;
468 for (lay = 0; lay <
L2_NLAY; lay++)
469 if (gsl_finite(ncd->
l2_z[track][xtrack][lay])
470 && ncd->
l2_z[track][xtrack][lay] <= 60.) {
471 atm_iasi.z[atm_iasi.np] = ncd->
l2_z[track][xtrack][lay];
472 atm_iasi.p[atm_iasi.np] = ncd->
l2_p[lay];
473 atm_iasi.t[atm_iasi.np] = ncd->
l2_t[track][xtrack][lay];
474 if ((++atm_iasi.np) > NP)
475 ERRMSG(
"Too many layers!");
483 for (ip = 0; ip < atm_iasi.np; ip++) {
484 zmax = GSL_MAX(zmax, atm_iasi.z[ip]);
485 zmin = GSL_MIN(zmin, atm_iasi.z[ip]);
489 for (ip = 0; ip < atm->np; ip++) {
492 intpol_atm(ctl, &atm_iasi, atm->z[ip], &p, &t, q, k);
496 if (atm->z[ip] > zmax)
497 w = GSL_MAX(1 - (atm->z[ip] - zmax) / 50, 0);
498 if (atm->z[ip] < zmin)
499 w = GSL_MAX(1 - (zmin - atm->z[ip]) / 50, 0);
502 atm->t[ip] = w * t + (1 - w) * atm->t[ip];
503 atm->p[ip] = w * p + (1 - w) * atm->p[ip];
518 printf(
"Read netCDF file: %s\n", filename);
519 NC(nc_open(filename, NC_WRITE, &ncd->
ncid));
522 NC(nc_inq_dimid(ncd->
ncid,
"L1_NTRACK", &dimid));
523 NC(nc_inq_dimlen(ncd->
ncid, dimid, &len));
527 NC(nc_inq_varid(ncd->
ncid,
"l1_time", &varid));
529 NC(nc_inq_varid(ncd->
ncid,
"l1_lon", &varid));
530 NC(nc_get_var_double(ncd->
ncid, varid, ncd->
l1_lon[0]));
531 NC(nc_inq_varid(ncd->
ncid,
"l1_lat", &varid));
532 NC(nc_get_var_double(ncd->
ncid, varid, ncd->
l1_lat[0]));
533 NC(nc_inq_varid(ncd->
ncid,
"l1_sat_z", &varid));
535 NC(nc_inq_varid(ncd->
ncid,
"l1_sat_lon", &varid));
537 NC(nc_inq_varid(ncd->
ncid,
"l1_sat_lat", &varid));
539 NC(nc_inq_varid(ncd->
ncid,
"l1_nu", &varid));
540 NC(nc_get_var_double(ncd->
ncid, varid, ncd->
l1_nu));
541 NC(nc_inq_varid(ncd->
ncid,
"l1_rad", &varid));
542 NC(nc_get_var_float(ncd->
ncid, varid, ncd->
l1_rad[0][0]));
545 NC(nc_inq_varid(ncd->
ncid,
"l2_z", &varid));
546 NC(nc_get_var_double(ncd->
ncid, varid, ncd->
l2_z[0][0]));
547 NC(nc_inq_varid(ncd->
ncid,
"l2_press", &varid));
548 NC(nc_get_var_double(ncd->
ncid, varid, ncd->
l2_p));
549 NC(nc_inq_varid(ncd->
ncid,
"l2_temp", &varid));
550 NC(nc_get_var_double(ncd->
ncid, varid, ncd->
l2_t[0][0]));
559 int dimid[10], c_id, p_id, t_id, z_id;
562 printf(
"Write netCDF file: %s\n", filename);
565 NC(nc_inq_dimid(ncd->
ncid,
"L1_NTRACK", &dimid[0]));
566 NC(nc_inq_dimid(ncd->
ncid,
"L1_NXTRACK", &dimid[1]));
572 if (nc_inq_dimid(ncd->
ncid,
"RET_NP", &dimid[2]) != NC_NOERR)
573 NC(nc_def_dim(ncd->
ncid,
"RET_NP", (
size_t) ncd->
np, &dimid[2]));
576 add_var(ncd->
ncid,
"ret_z",
"km",
"altitude", NC_FLOAT, &dimid[2], &z_id,
578 add_var(ncd->
ncid,
"ret_press",
"hPa",
"pressure", NC_FLOAT, dimid, &p_id,
580 add_var(ncd->
ncid,
"ret_temp",
"K",
"temperature", NC_FLOAT, dimid, &t_id,
582 add_var(ncd->
ncid,
"ret_chisq",
"1",
"chi^2 value of fit", NC_FLOAT, dimid,
int main(int argc, char *argv[])
#define NC(cmd)
Execute netCDF library command and check result.
#define L1_NXTRACK
Across-track size of IASI radiance granule (don't change).
#define L1_NTRACK
Maximum along-track size of IASI radiance granule (don't change).
#define L2_NXTRACK
Across-track size of IASI 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 IASI 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 IASI radiance channels (don't change).
#define L2_NTRACK
Maximum along-track size of IASI retrieval granule (don't change).
void init_l2(ncd_t *ncd, int track, int xtrack, ctl_t *ctl, atm_t *atm)
Initialize with IASI Level-2 data.
float ret_z[NP]
Altitude [km].
double l2_z[L2_NTRACK][L2_NXTRACK][L2_NLAY]
Altitude [km].
int ntrack
Number of tacks.
float ret_chisq[L1_NTRACK *L1_NXTRACK]
chi^2 value of fit.
double l1_sat_lon[L1_NTRACK]
Satellite longitude [deg].
double l1_nu[L1_NCHAN]
Channel frequencies [cm^-1].
float ret_p[L1_NTRACK *L1_NXTRACK]
Pressure [hPa].
int np
Number of retrieval altitudes.
float ret_t[L1_NTRACK *L1_NXTRACK *NP]
Temperature [K].
double l1_sat_lat[L1_NTRACK]
Satellite latitude [deg].
double l1_sat_z[L1_NTRACK]
Satellite altitude [km].
double l1_time[L1_NTRACK][L1_NXTRACK]
Time (seconds since 2000-01-01T00:00Z).
double l1_lat[L1_NTRACK][L1_NXTRACK]
Footprint latitude [deg].
double l2_p[L2_NLAY]
Pressure [hPa].
double l2_t[L2_NTRACK][L2_NXTRACK][L2_NLAY]
Temperature [K].
double l1_lon[L1_NTRACK][L1_NXTRACK]
Footprint longitude [deg].
float l1_rad[L1_NTRACK][L1_NXTRACK][L1_NCHAN]
Radiance [W/(m^2 sr cm^-1)].