IASI Code Collection
retrieval.c
Go to the documentation of this file.
1/*
2 This file is part of the IASI Code Collection.
3
4 the IASI Code Collections is free software: you can redistribute it
5 and/or modify it under the terms of the GNU General Public License
6 as published by the Free Software Foundation, either version 3 of
7 the License, or (at your option) any later version.
8
9 The IASI Code Collection is distributed in the hope that it will be
10 useful, but WITHOUT ANY WARRANTY; without even the implied warranty
11 of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 General Public License for more details.
13
14 You should have received a copy of the GNU General Public License
15 along with the IASI Code Collection. If not, see
16 <http://www.gnu.org/licenses/>.
17
18 Copyright (C) 2019-2025 Forschungszentrum Juelich GmbH
19*/
20
26#include <mpi.h>
27#include <omp.h>
28#include <netcdf.h>
29#include "jurassic.h"
30
31/* ------------------------------------------------------------
32 Macros...
33 ------------------------------------------------------------ */
34
36#define NC(cmd) { \
37 int nc_result=(cmd); \
38 if(nc_result!=NC_NOERR) \
39 ERRMSG("%s", nc_strerror(nc_result)); \
40}
41
42/* ------------------------------------------------------------
43 Dimensions...
44 ------------------------------------------------------------ */
45
47#define L1_NCHAN 33
48
50#define L1_NTRACK 1800
51
53#define L1_NXTRACK 60
54
56#define L2_NLAY 27
57
59#define L2_NTRACK 1800
60
62#define L2_NXTRACK 60
63
64/* ------------------------------------------------------------
65 Structs...
66 ------------------------------------------------------------ */
67
69typedef struct {
70
72 int ncid;
73
75 int np;
76
78 int ntrack;
79
81 double l1_time[L1_NTRACK][L1_NXTRACK];
82
84 double l1_lon[L1_NTRACK][L1_NXTRACK];
85
87 double l1_lat[L1_NTRACK][L1_NXTRACK];
88
90 double l1_sat_z[L1_NTRACK];
91
93 double l1_sat_lon[L1_NTRACK];
94
96 double l1_sat_lat[L1_NTRACK];
97
99 double l1_nu[L1_NCHAN];
100
103
106
108 double l2_p[L2_NLAY];
109
112
114 float ret_z[NP];
115
117 float ret_p[L1_NTRACK * L1_NXTRACK];
118
120 float ret_t[L1_NTRACK * L1_NXTRACK * NP];
121
123 float ret_chisq[L1_NTRACK * L1_NXTRACK];
124
125} ncd_t;
126
127/* ------------------------------------------------------------
128 Functions...
129 ------------------------------------------------------------ */
130
132void add_var(
133 int ncid,
134 const char *varname,
135 const char *unit,
136 const char *longname,
137 int type,
138 int dimid[],
139 int *varid,
140 int ndims);
141
143void buffer_nc(
144 atm_t * atm,
145 double chisq,
146 ncd_t * ncd,
147 int track,
148 int xtrack,
149 int np0,
150 int np1);
151
153void init_l2(
154 ncd_t * ncd,
155 int track,
156 int xtrack,
157 ctl_t * ctl,
158 atm_t * atm);
159
161void read_nc(
162 char *filename,
163 ncd_t * ncd);
164
166void write_nc(
167 char *filename,
168 ncd_t * ncd);
169
170/* ------------------------------------------------------------
171 Main...
172 ------------------------------------------------------------ */
173
175 int argc,
176 char *argv[]) {
177
178 static ctl_t ctl;
179 static atm_t atm_apr, atm_clim, atm_i;
180 static obs_t obs_i, obs_meas;
181 static ncd_t ncd;
182 static ret_t ret;
183
184 FILE *in;
185
186 char filename[LEN], filename2[LEN];
187
188 double chisq, sza_thresh, z[NP], t0;
189
190 int channel[ND], i, id, ip, iz, nz, ntask = -1, rank, size,
191 np0, np1, track, track0, track1, xtrack, xtrack0, xtrack1,
192 task0, task1, debug;
193
194 /* ------------------------------------------------------------
195 Init...
196 ------------------------------------------------------------ */
197
198 /* MPI... */
199 MPI_Init(&argc, &argv);
200 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
201 MPI_Comm_size(MPI_COMM_WORLD, &size);
202
203 /* Measure CPU time... */
204 TIMER("total", 1);
205
206 /* Check arguments... */
207 if (argc < 3)
208 ERRMSG("Give parameters: <ctl> <filelist>");
209
210 /* Read control parameters... */
211 read_ctl(argc, argv, &ctl);
212 read_ret(argc, argv, &ctl, &ret);
213 debug = (int) scan_ctl(argc, argv, "DEBUG", -1, "1", NULL);
214
215 /* Initialize look-up tables... */
216 tbl_t *tbl = read_tbl(&ctl);
217
218 /* Read retrieval grid... */
219 nz = (int) scan_ctl(argc, argv, "NZ", -1, "", NULL);
220 if (nz > NP)
221 ERRMSG("Too many altitudes!");
222 for (iz = 0; iz < nz; iz++)
223 z[iz] = scan_ctl(argc, argv, "Z", iz, "", NULL);
224
225 /* Read task range... */
226 task0 = (int) scan_ctl(argc, argv, "TASK_MIN", -1, "0", NULL);
227 task1 = (int) scan_ctl(argc, argv, "TASK_MAX", -1, "99999", NULL);
228
229 /* Read track range... */
230 track0 = (int) scan_ctl(argc, argv, "TRACK_MIN", -1, "0", NULL);
231 track1 = (int) scan_ctl(argc, argv, "TRACK_MAX", -1, "99999", NULL);
232
233 /* Read xtrack range... */
234 xtrack0 = (int) scan_ctl(argc, argv, "XTRACK_MIN", -1, "0", NULL);
235 xtrack1 = (int) scan_ctl(argc, argv, "XTRACK_MAX", -1, "59", NULL);
236
237 /* Read height range... */
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);
241
242 /* SZA threshold... */
243 sza_thresh = scan_ctl(argc, argv, "SZA", -1, "96", NULL);
244
245 /* ------------------------------------------------------------
246 Distribute granules...
247 ------------------------------------------------------------ */
248
249 /* Open filelist... */
250 printf("Read filelist: %s\n", argv[2]);
251 if (!(in = fopen(argv[2], "r")))
252 ERRMSG("Cannot open filelist!");
253
254 /* Loop over netCDF files... */
255 while (fscanf(in, "%s", filename) != EOF) {
256
257 /* Distribute files with MPI... */
258 if ((++ntask) % size != rank)
259 continue;
260
261 /* Check task range... */
262 if (ntask < task0 || ntask > task1)
263 continue;
264
265 /* Write info... */
266 printf("Retrieve file %s on rank %d of %d (with %d threads)...\n",
267 filename, rank + 1, size, omp_get_max_threads());
268
269 /* ------------------------------------------------------------
270 Initialize retrieval...
271 ------------------------------------------------------------ */
272
273 /* Read netCDF file... */
274 read_nc(filename, &ncd);
275
276 /* Adjust number of tracks... */
277 if (track1 >= ncd.ntrack)
278 track1 = ncd.ntrack - 1;
279
280 /* Identify radiance channels... */
281 for (id = 0; id < ctl.nd; id++) {
282 channel[id] = -999;
283 for (i = 0; i < L1_NCHAN; i++)
284 if (fabs(ctl.nu[id] - ncd.l1_nu[i]) < 0.1)
285 channel[id] = i;
286 if (channel[id] < 0)
287 ERRMSG("Cannot identify radiance channel!");
288 }
289
290 /* Set climatological data for center of granule... */
291 atm_clim.np = nz;
292 for (iz = 0; iz < nz; iz++)
293 atm_clim.z[iz] = z[iz];
294 climatology(&ctl, &atm_clim);
295
296 /* ------------------------------------------------------------
297 Retrieval...
298 ------------------------------------------------------------ */
299
300 /* Loop over swaths... */
301 for (track = track0; track <= track1; track++) {
302
303 /* Loop over scan... */
304 for (xtrack = xtrack0; xtrack <= xtrack1; xtrack++) {
305
306 /* Init timer... */
307 t0 = omp_get_wtime();
308
309 /* Store observation data... */
310 obs_meas.nr = 1;
311 obs_meas.time[0] = ncd.l1_time[track][xtrack];
312 obs_meas.obsz[0] = ncd.l1_sat_z[track];
313 obs_meas.obslon[0] = ncd.l1_sat_lon[track];
314 obs_meas.obslat[0] = ncd.l1_sat_lat[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]];
319
320 /* Flag out 4 micron channels for daytime measurements... */
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;
326
327 /* Prepare atmospheric data... */
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];
333 }
334
335 /* Merge Level-2 data... */
336 init_l2(&ncd, track, xtrack, &ctl, &atm_apr);
337
338 /* Retrieval... */
339 optimal_estimation(&ret, &ctl, tbl, &obs_meas, &obs_i,
340 &atm_apr, &atm_i, &chisq);
341
342 /* Buffer results... */
343 buffer_nc(&atm_i, chisq, &ncd, track, xtrack, np0, np1);
344
345 /* Write debug information... */
346 if (debug >= 1)
347 printf
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);
350 if (debug >= 2) {
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);
359 }
360 }
361 }
362
363 /* ------------------------------------------------------------
364 Finalize...
365 ------------------------------------------------------------ */
366
367 /* Write netCDF file... */
368 write_nc(filename, &ncd);
369
370 /* Write info... */
371 printf("Retrieval finished on rank %d of %d!\n", rank, size);
372 }
373
374 /* Close file list... */
375 fclose(in);
376
377 /* Measure CPU time... */
378 TIMER("total", 3);
379
380 /* Report memory usage... */
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.);
387
388 /* Report problem size... */
389 printf("SIZE_TASKS = %d\n", size);
390 printf("SIZE_THREADS = %d\n", omp_get_max_threads());
391
392 /* MPI... */
393 MPI_Finalize();
394
395 return EXIT_SUCCESS;
396}
397
398/*****************************************************************************/
399
401 int ncid,
402 const char *varname,
403 const char *unit,
404 const char *longname,
405 int type,
406 int dimid[],
407 int *varid,
408 int ndims) {
409
410 /* Check if variable exists... */
411 if (nc_inq_varid(ncid, varname, varid) != NC_NOERR) {
412
413 /* Define variable... */
414 NC(nc_def_var(ncid, varname, type, ndims, dimid, varid));
415
416 /* Set long name... */
417 NC(nc_put_att_text
418 (ncid, *varid, "long_name", strlen(longname), longname));
419
420 /* Set units... */
421 NC(nc_put_att_text(ncid, *varid, "units", strlen(unit), unit));
422 }
423}
424
425/*****************************************************************************/
426
428 atm_t *atm,
429 double chisq,
430 ncd_t *ncd,
431 int track,
432 int xtrack,
433 int np0,
434 int np1) {
435
436 int ip;
437
438 /* Set number of data points... */
439 ncd->np = np1 - np0 + 1;
440
441 /* Save retrieval data... */
442 ncd->ret_chisq[track * L1_NXTRACK + xtrack] = (float) chisq;
443 ncd->ret_p[track * L1_NXTRACK + xtrack] = (float) atm->p[np0];
444 for (ip = np0; ip <= np1; ip++) {
445 ncd->ret_z[ip - np0] = (float) atm->z[ip];
446 ncd->ret_t[(track * L1_NXTRACK + xtrack) * ncd->np + ip - np0] =
447 (gsl_finite(chisq) ? (float) atm->t[ip] : GSL_NAN);
448 }
449}
450
451/************************************************************************/
452
454 ncd_t *ncd,
455 int track,
456 int xtrack,
457 ctl_t *ctl,
458 atm_t *atm) {
459
460 static atm_t atm_iasi;
461
462 double k[NW], p, q[NG], t, w, zmax = 0, zmin = 1000;
463
464 int ip, lay;
465
466 /* Store IASI data in atmospheric data struct... */
467 atm_iasi.np = 0;
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!");
476 }
477
478 /* Check number of levels... */
479 if (atm_iasi.np < 2)
480 return;
481
482 /* Get height range of IASI data... */
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]);
486 }
487
488 /* Merge IASI data... */
489 for (ip = 0; ip < atm->np; ip++) {
490
491 /* Interpolate IASI data... */
492 intpol_atm(ctl, &atm_iasi, atm->z[ip], &p, &t, q, k);
493
494 /* Weighting factor... */
495 w = 1;
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);
500
501 /* Merge... */
502 atm->t[ip] = w * t + (1 - w) * atm->t[ip];
503 atm->p[ip] = w * p + (1 - w) * atm->p[ip];
504 }
505}
506
507/*****************************************************************************/
508
510 char *filename,
511 ncd_t *ncd) {
512
513 int dimid, varid;
514
515 size_t len;
516
517 /* Open netCDF file... */
518 printf("Read netCDF file: %s\n", filename);
519 NC(nc_open(filename, NC_WRITE, &ncd->ncid));
520
521 /* Read number of tracks... */
522 NC(nc_inq_dimid(ncd->ncid, "L1_NTRACK", &dimid));
523 NC(nc_inq_dimlen(ncd->ncid, dimid, &len));
524 ncd->ntrack = (int) len;
525
526 /* Read Level-1 data... */
527 NC(nc_inq_varid(ncd->ncid, "l1_time", &varid));
528 NC(nc_get_var_double(ncd->ncid, varid, ncd->l1_time[0]));
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));
534 NC(nc_get_var_double(ncd->ncid, varid, ncd->l1_sat_z));
535 NC(nc_inq_varid(ncd->ncid, "l1_sat_lon", &varid));
536 NC(nc_get_var_double(ncd->ncid, varid, ncd->l1_sat_lon));
537 NC(nc_inq_varid(ncd->ncid, "l1_sat_lat", &varid));
538 NC(nc_get_var_double(ncd->ncid, varid, ncd->l1_sat_lat));
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]));
543
544 /* Read Level-2 data... */
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]));
551}
552
553/*****************************************************************************/
554
556 char *filename,
557 ncd_t *ncd) {
558
559 int dimid[10], c_id, p_id, t_id, z_id;
560
561 /* Create netCDF file... */
562 printf("Write netCDF file: %s\n", filename);
563
564 /* Read existing dimensions... */
565 NC(nc_inq_dimid(ncd->ncid, "L1_NTRACK", &dimid[0]));
566 NC(nc_inq_dimid(ncd->ncid, "L1_NXTRACK", &dimid[1]));
567
568 /* Set define mode... */
569 NC(nc_redef(ncd->ncid));
570
571 /* Set new dimensions... */
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]));
574
575 /* Set new variables... */
576 add_var(ncd->ncid, "ret_z", "km", "altitude", NC_FLOAT, &dimid[2], &z_id,
577 1);
578 add_var(ncd->ncid, "ret_press", "hPa", "pressure", NC_FLOAT, dimid, &p_id,
579 2);
580 add_var(ncd->ncid, "ret_temp", "K", "temperature", NC_FLOAT, dimid, &t_id,
581 3);
582 add_var(ncd->ncid, "ret_chisq", "1", "chi^2 value of fit", NC_FLOAT, dimid,
583 &c_id, 2);
584
585 /* Leave define mode... */
586 NC(nc_enddef(ncd->ncid));
587
588 /* Write data... */
589 NC(nc_put_var_float(ncd->ncid, z_id, ncd->ret_z));
590 NC(nc_put_var_float(ncd->ncid, p_id, ncd->ret_p));
591 NC(nc_put_var_float(ncd->ncid, t_id, ncd->ret_t));
592 NC(nc_put_var_float(ncd->ncid, c_id, ncd->ret_chisq));
593
594 /* Close netCDF file... */
595 NC(nc_close(ncd->ncid));
596}
int main(int argc, char *argv[])
Definition: retrieval.c:174
#define NC(cmd)
Execute netCDF library command and check result.
Definition: retrieval.c:36
#define L1_NXTRACK
Across-track size of IASI radiance granule (don't change).
Definition: retrieval.c:53
#define L1_NTRACK
Maximum along-track size of IASI radiance granule (don't change).
Definition: retrieval.c:50
#define L2_NXTRACK
Across-track size of IASI retrieval granule (don't change).
Definition: retrieval.c:62
void write_nc(char *filename, ncd_t *ncd)
Write to netCDF file...
Definition: retrieval.c:555
void read_nc(char *filename, ncd_t *ncd)
Read netCDF file.
Definition: retrieval.c:509
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.
Definition: retrieval.c:400
#define L2_NLAY
Number of IASI pressure layers (don't change).
Definition: retrieval.c:56
void buffer_nc(atm_t *atm, double chisq, ncd_t *ncd, int track, int xtrack, int np0, int np1)
Buffer netCDF data.
Definition: retrieval.c:427
#define L1_NCHAN
Number of IASI radiance channels (don't change).
Definition: retrieval.c:47
#define L2_NTRACK
Maximum along-track size of IASI retrieval granule (don't change).
Definition: retrieval.c:59
void init_l2(ncd_t *ncd, int track, int xtrack, ctl_t *ctl, atm_t *atm)
Initialize with IASI Level-2 data.
Definition: retrieval.c:453
Buffer for netCDF data.
Definition: retrieval.c:69
float ret_z[NP]
Altitude [km].
Definition: retrieval.c:114
double l2_z[L2_NTRACK][L2_NXTRACK][L2_NLAY]
Altitude [km].
Definition: retrieval.c:105
int ntrack
Number of tacks.
Definition: retrieval.c:78
float ret_chisq[L1_NTRACK *L1_NXTRACK]
chi^2 value of fit.
Definition: retrieval.c:123
double l1_sat_lon[L1_NTRACK]
Satellite longitude [deg].
Definition: retrieval.c:93
double l1_nu[L1_NCHAN]
Channel frequencies [cm^-1].
Definition: retrieval.c:99
float ret_p[L1_NTRACK *L1_NXTRACK]
Pressure [hPa].
Definition: retrieval.c:117
int ncid
NetCDF file ID.
Definition: retrieval.c:72
int np
Number of retrieval altitudes.
Definition: retrieval.c:75
float ret_t[L1_NTRACK *L1_NXTRACK *NP]
Temperature [K].
Definition: retrieval.c:120
double l1_sat_lat[L1_NTRACK]
Satellite latitude [deg].
Definition: retrieval.c:96
double l1_sat_z[L1_NTRACK]
Satellite altitude [km].
Definition: retrieval.c:90
double l1_time[L1_NTRACK][L1_NXTRACK]
Time (seconds since 2000-01-01T00:00Z).
Definition: retrieval.c:81
double l1_lat[L1_NTRACK][L1_NXTRACK]
Footprint latitude [deg].
Definition: retrieval.c:87
double l2_p[L2_NLAY]
Pressure [hPa].
Definition: retrieval.c:108
double l2_t[L2_NTRACK][L2_NXTRACK][L2_NLAY]
Temperature [K].
Definition: retrieval.c:111
double l1_lon[L1_NTRACK][L1_NXTRACK]
Footprint longitude [deg].
Definition: retrieval.c:84
float l1_rad[L1_NTRACK][L1_NXTRACK][L1_NCHAN]
Radiance [W/(m^2 sr cm^-1)].
Definition: retrieval.c:102