AIRS Code Collection
retrieval.c
Go to the documentation of this file.
1/*
2 This file is part of the AIRS Code Collection.
3
4 the AIRS 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 AIRS 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 AIRS 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 34
48
50#define L1_NTRACK 135
51
53#define L1_NXTRACK 90
54
56#define L2_NLAY 27
57
59#define L2_NTRACK 45
60
62#define L2_NXTRACK 30
63
64/* ------------------------------------------------------------
65 Structs...
66 ------------------------------------------------------------ */
67
69typedef struct {
70
72 int ncid;
73
75 int np;
76
78 double l1_time[L1_NTRACK][L1_NXTRACK];
79
81 double l1_lon[L1_NTRACK][L1_NXTRACK];
82
84 double l1_lat[L1_NTRACK][L1_NXTRACK];
85
87 double l1_sat_z[L1_NTRACK];
88
90 double l1_sat_lon[L1_NTRACK];
91
93 double l1_sat_lat[L1_NTRACK];
94
96 double l1_nu[L1_NCHAN];
97
99 float l1_rad[L1_NTRACK][L1_NXTRACK][L1_NCHAN];
100
102 double l2_z[L2_NTRACK][L2_NXTRACK][L2_NLAY];
103
105 double l2_p[L2_NLAY];
106
108 double l2_t[L2_NTRACK][L2_NXTRACK][L2_NLAY];
109
111 float ret_z[NP];
112
114 float ret_p[L1_NTRACK * L1_NXTRACK];
115
117 float ret_t[L1_NTRACK * L1_NXTRACK * NP];
118
119} ncd_t;
120
121/* ------------------------------------------------------------
122 Functions...
123 ------------------------------------------------------------ */
124
126void add_var(
127 int ncid,
128 const char *varname,
129 const char *unit,
130 const char *longname,
131 int type,
132 int dimid[],
133 int *varid,
134 int ndims);
135
137void buffer_nc(
138 atm_t * atm,
139 double chisq,
140 ncd_t * ncd,
141 int track,
142 int xtrack,
143 int np0,
144 int np1);
145
147void fill_gaps(
148 double x[L2_NTRACK][L2_NXTRACK][L2_NLAY],
149 double cx,
150 double cy);
151
153void init_l2(
154 ncd_t * ncd,
155 int track,
156 int xtrack,
157 ctl_t * ctl,
158 atm_t * atm);
159
162 ret_t * ret,
163 ctl_t * ctl,
164 tbl_t * tbl,
165 obs_t * obs_meas,
166 obs_t * obs_i,
167 atm_t * atm_apr,
168 atm_t * atm_i,
169 double *chisq);
170
172void read_nc(
173 char *filename,
174 ncd_t * ncd);
175
177void write_nc(
178 char *filename,
179 ncd_t * ncd);
180
181/* ------------------------------------------------------------
182 Main...
183 ------------------------------------------------------------ */
184
186 int argc,
187 char *argv[]) {
188
189 static ctl_t ctl;
190 static atm_t atm_apr, atm_clim, atm_i;
191 static obs_t obs_i, obs_meas;
192 static ncd_t ncd;
193 static ret_t ret;
194
195 FILE *in;
196
197 char filename[LEN];
198
199 double chisq, chisq_min, chisq_max, chisq_mean, z[NP];
200
201 int channel[ND], m, ntask = -1, rank, size;
202
203 /* ------------------------------------------------------------
204 Init...
205 ------------------------------------------------------------ */
206
207 /* MPI... */
208 MPI_Init(&argc, &argv);
209 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
210 MPI_Comm_size(MPI_COMM_WORLD, &size);
211
212 /* Measure CPU time... */
213 TIMER("total", 1);
214
215 /* Check arguments... */
216 if (argc < 3)
217 ERRMSG("Give parameters: <ctl> <filelist>");
218
219 /* Read control parameters... */
220 read_ctl(argc, argv, &ctl);
221 read_ret(argc, argv, &ctl, &ret);
222
223 /* Initialize look-up tables... */
224 tbl_t *tbl = read_tbl(&ctl);
225
226 /* Read retrieval grid... */
227 const int nz = (int) scan_ctl(argc, argv, "NZ", -1, "", NULL);
228 if (nz > NP)
229 ERRMSG("Too many altitudes!");
230 for (int iz = 0; iz < nz; iz++)
231 z[iz] = scan_ctl(argc, argv, "Z", iz, "", NULL);
232
233 /* Read track range... */
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);
236
237 /* Read xtrack range... */
238 const int xtrack0 = (int) scan_ctl(argc, argv, "XTRACK_MIN", -1, "0", NULL);
239 const int xtrack1 =
240 (int) scan_ctl(argc, argv, "XTRACK_MAX", -1, "89", NULL);
241
242 /* Read height range... */
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);
246
247 /* Background smoothing... */
248 const double sx = scan_ctl(argc, argv, "SX", -1, "8", NULL);
249 const double sy = scan_ctl(argc, argv, "SY", -1, "2", NULL);
250
251 /* SZA threshold... */
252 const double sza_thresh = scan_ctl(argc, argv, "SZA", -1, "96", NULL);
253
254 /* ------------------------------------------------------------
255 Distribute granules...
256 ------------------------------------------------------------ */
257
258 /* Open filelist... */
259 printf("Read filelist: %s\n", argv[2]);
260 if (!(in = fopen(argv[2], "r")))
261 ERRMSG("Cannot open filelist!");
262
263 /* Loop over netCDF files... */
264 while (fscanf(in, "%s", filename) != EOF) {
265
266 /* Distribute files with MPI... */
267 if ((++ntask) % size != rank)
268 continue;
269
270 /* Write info... */
271 printf("Retrieve file %s on rank %d of %d (with %d threads)...\n",
272 filename, rank + 1, size, omp_get_max_threads());
273
274 /* ------------------------------------------------------------
275 Initialize retrieval...
276 ------------------------------------------------------------ */
277
278 /* Read netCDF file... */
279 read_nc(filename, &ncd);
280
281 /* Identify radiance channels... */
282 for (int id = 0; id < ctl.nd; id++) {
283 channel[id] = -999;
284 for (int i = 0; i < L1_NCHAN; i++)
285 if (fabs(ctl.nu[id] - ncd.l1_nu[i]) < 0.1)
286 channel[id] = i;
287 if (channel[id] < 0)
288 ERRMSG("Cannot identify radiance channel!");
289 }
290
291 /* Fill data gaps... */
292 fill_gaps(ncd.l2_t, sx, sy);
293 fill_gaps(ncd.l2_z, sx, sy);
294
295 /* Set climatological data for center of granule... */
296 atm_clim.np = nz;
297 for (int iz = 0; iz < nz; iz++)
298 atm_clim.z[iz] = z[iz];
299 climatology(&ctl, &atm_clim);
300
301 /* ------------------------------------------------------------
302 Retrieval...
303 ------------------------------------------------------------ */
304
305 /* Get chi^2 statistics... */
306 chisq_min = 1e100;
307 chisq_max = -1e100;
308 chisq_mean = 0;
309 m = 0;
310
311 /* Loop over swaths... */
312 for (int track = track0; track <= track1; track++) {
313
314 /* Measure CPU time... */
315 TIMER("retrieval", 1);
316
317 /* Loop over scan... */
318 for (int xtrack = xtrack0; xtrack <= xtrack1; xtrack++) {
319
320 /* Store observation data... */
321 obs_meas.nr = 1;
322 obs_meas.time[0] = ncd.l1_time[track][xtrack];
323 obs_meas.obsz[0] = ncd.l1_sat_z[track];
324 obs_meas.obslon[0] = ncd.l1_sat_lon[track];
325 obs_meas.obslat[0] = ncd.l1_sat_lat[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]];
330
331 /* Flag out 4 micron channels for daytime measurements... */
332 if (sza(obs_meas.time[0], obs_meas.obslon[0], obs_meas.obslat[0])
333 < sza_thresh)
334 for (int id = 0; id < ctl.nd; id++)
335 if (ctl.nu[id] >= 2000)
336 obs_meas.rad[id][0] = GSL_NAN;
337
338 /* Prepare atmospheric data... */
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];
344 }
345
346 /* Merge Level-2 data... */
347 init_l2(&ncd, track, xtrack, &ctl, &atm_apr);
348
349 /* Retrieval... */
350 optimal_estimation(&ret, &ctl, tbl, &obs_meas, &obs_i,
351 &atm_apr, &atm_i, &chisq);
352
353 /* Get chi^2 statistics... */
354 if (gsl_finite(chisq)) {
355 chisq_min = GSL_MIN(chisq_min, chisq);
356 chisq_max = GSL_MAX(chisq_max, chisq);
357 chisq_mean += chisq;
358 m++;
359 }
360
361 /* Buffer results... */
362 buffer_nc(&atm_i, chisq, &ncd, track, xtrack, np0, np1);
363 }
364
365 /* Measure CPU time... */
366 TIMER("retrieval", 3);
367 }
368
369 /* ------------------------------------------------------------
370 Finalize...
371 ------------------------------------------------------------ */
372
373 /* Write netCDF file... */
374 write_nc(filename, &ncd);
375
376 /* Write info... */
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);
380 }
381
382 /* Close file list... */
383 fclose(in);
384
385 /* Measure CPU time... */
386 TIMER("total", 3);
387
388 /* Report memory usage... */
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.);
395
396 /* Report problem size... */
397 printf("SIZE_TASKS = %d\n", size);
398 printf("SIZE_THREADS = %d\n", omp_get_max_threads());
399
400 /* MPI... */
401 MPI_Finalize();
402
403 return EXIT_SUCCESS;
404}
405
406/*****************************************************************************/
407
409 int ncid,
410 const char *varname,
411 const char *unit,
412 const char *longname,
413 int type,
414 int dimid[],
415 int *varid,
416 int ndims) {
417
418 /* Check if variable exists... */
419 if (nc_inq_varid(ncid, varname, varid) != NC_NOERR) {
420
421 /* Define variable... */
422 NC(nc_def_var(ncid, varname, type, ndims, dimid, varid));
423
424 /* Set long name... */
425 NC(nc_put_att_text
426 (ncid, *varid, "long_name", strlen(longname), longname));
427
428 /* Set units... */
429 NC(nc_put_att_text(ncid, *varid, "units", strlen(unit), unit));
430 }
431}
432
433/*****************************************************************************/
434
436 atm_t *atm,
437 double chisq,
438 ncd_t *ncd,
439 int track,
440 int xtrack,
441 int np0,
442 int np1) {
443
444 /* Set number of data points... */
445 ncd->np = np1 - np0 + 1;
446
447 /* Save retrieval data... */
448 for (int ip = np0; ip <= np1; ip++) {
449 ncd->ret_z[ip - np0] = (float) atm->z[ip];
450 ncd->ret_p[track * L1_NXTRACK + xtrack] = (float) atm->p[np0];
451 ncd->ret_t[(track * L1_NXTRACK + xtrack) * ncd->np + ip - np0] =
452 (gsl_finite(chisq) ? (float) atm->t[ip] : GSL_NAN);
453 }
454}
455
456/************************************************************************/
457
459 double x[L2_NTRACK][L2_NXTRACK][L2_NLAY],
460 double cx,
461 double cy) {
462
463 double help[L2_NTRACK][L2_NXTRACK];
464
465 /* Loop over layers... */
466 for (int lay = 0; lay < L2_NLAY; lay++) {
467
468 /* Loop over grid points... */
469 for (int track = 0; track < L2_NTRACK; track++)
470 for (int xtrack = 0; xtrack < L2_NXTRACK; xtrack++) {
471
472 /* Init... */
473 help[track][xtrack] = 0;
474 double wsum = 0;
475
476 /* Averrage data points... */
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];
484 wsum += w;
485 }
486
487 /* Normalize... */
488 if (wsum > 0)
489 help[track][xtrack] /= wsum;
490 else
491 help[track][xtrack] = GSL_NAN;
492 }
493
494 /* Copy grid points... */
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];
498 }
499}
500
501/************************************************************************/
502
504 ncd_t *ncd,
505 int track,
506 int xtrack,
507 ctl_t *ctl,
508 atm_t *atm) {
509
510 static atm_t atm_airs;
511
512 double k[NW], p, q[NG], t, zmax = 0, zmin = 1000;
513
514 /* Reset track- and xtrack-index to match Level-2 data... */
515 track /= 3;
516 xtrack /= 3;
517
518 /* Store AIRS data in atmospheric data struct... */
519 atm_airs.np = 0;
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!");
527 }
528
529 /* Check number of levels... */
530 if (atm_airs.np <= 0)
531 return;
532
533 /* Get height range of AIRS data... */
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]);
537 }
538
539 /* Merge AIRS data... */
540 for (int ip = 0; ip < atm->np; ip++) {
541
542 /* Interpolate AIRS data... */
543 intpol_atm(ctl, &atm_airs, atm->z[ip], &p, &t, q, k);
544
545 /* Weighting factor... */
546 double w = 1;
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);
551
552 /* Merge... */
553 atm->t[ip] = w * t + (1 - w) * atm->t[ip];
554 atm->p[ip] = w * p + (1 - w) * atm->p[ip];
555 }
556}
557
558/*****************************************************************************/
559
561 ret_t *ret,
562 ctl_t *ctl,
563 tbl_t *tbl,
564 obs_t *obs_meas,
565 obs_t *obs_i,
566 atm_t *atm_apr,
567 atm_t *atm_i,
568 double *chisq) {
569
570 static int ipa[N], iqa[N];
571
572 double disq = 0, lmpar = 0.001;
573
574 /* ------------------------------------------------------------
575 Initialize...
576 ------------------------------------------------------------ */
577
578 /* Get sizes... */
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) {
582 *chisq = GSL_NAN;
583 return;
584 }
585
586 /* Allocate... */
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);
591
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);
604
605 /* Set initial state... */
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);
609
610 /* Set state vectors and observation vectors... */
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);
615
616 /* Set inverse a priori covariance S_a^-1... */
617 set_cov_apr(ret, ctl, atm_apr, iqa, ipa, s_a_inv);
618 matrix_invert(s_a_inv);
619
620 /* Get measurement errors... */
621 set_cov_meas(ret, ctl, obs_meas, sig_noise, sig_formod, sig_eps_inv);
622
623 /* Determine dx = x_i - x_a and dy = y - F(x_i) ... */
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);
628
629 /* Compute cost function... */
630 *chisq = cost_function(dx, dy, s_a_inv, sig_eps_inv);
631
632 /* Compute initial kernel... */
633 kernel(ctl, tbl, atm_i, obs_i, k_i);
634
635 /* ------------------------------------------------------------
636 Levenberg-Marquardt minimization...
637 ------------------------------------------------------------ */
638
639 /* Outer loop... */
640 for (int it = 1; it <= ret->conv_itmax; it++) {
641
642 /* Store current cost function value... */
643 double chisq_old = *chisq;
644
645 /* Compute kernel matrix K_i... */
646 if (it > 1 && it % ret->kernel_recomp == 0)
647 kernel(ctl, tbl, atm_i, obs_i, k_i);
648
649 /* Compute K_i^T * S_eps^{-1} * K_i ... */
650 if (it == 1 || it % ret->kernel_recomp == 0)
651 matrix_product(k_i, sig_eps_inv, 1, cov);
652
653 /* Determine b = K_i^T * S_eps^{-1} * dy - S_a^{-1} * dx ... */
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);
659
660 /* Inner loop... */
661 for (int it2 = 0; it2 < 20; it2++) {
662
663 /* Compute A = (1 + lmpar) * S_a^{-1} + K_i^T * S_eps^{-1} * K_i ... */
664 gsl_matrix_memcpy(a, s_a_inv);
665 gsl_matrix_scale(a, 1 + lmpar);
666 gsl_matrix_add(a, cov);
667
668 /* Solve A * x_step = b by means of Cholesky decomposition... */
669 gsl_linalg_cholesky_decomp(a);
670 gsl_linalg_cholesky_solve(a, b, x_step);
671
672 /* Update atmospheric state... */
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);
677
678 /* Check atmospheric state... */
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);
686 }
687
688 /* Forward calculation... */
689 formod(ctl, tbl, atm_i, obs_i);
690 obs2y(ctl, obs_i, y_i, NULL, NULL);
691
692 /* Determine dx = x_i - x_a and dy = y - F(x_i) ... */
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);
697
698 /* Compute cost function... */
699 *chisq = cost_function(dx, dy, s_a_inv, sig_eps_inv);
700
701 /* Modify Levenberg-Marquardt parameter... */
702 if (*chisq > chisq_old) {
703 lmpar *= 10;
704 gsl_vector_sub(x_i, x_step);
705 } else {
706 lmpar /= 10;
707 break;
708 }
709 }
710
711 /* Get normalized step size in state space... */
712 gsl_blas_ddot(x_step, b, &disq);
713 disq /= (double) n;
714
715 /* Convergence test... */
716 if ((it == 1 || it % ret->kernel_recomp == 0) && disq < ret->conv_dmin)
717 break;
718 }
719
720 /* ------------------------------------------------------------
721 Finalize...
722 ------------------------------------------------------------ */
723
724 gsl_matrix_free(a);
725 gsl_matrix_free(cov);
726 gsl_matrix_free(k_i);
727 gsl_matrix_free(s_a_inv);
728
729 gsl_vector_free(b);
730 gsl_vector_free(dx);
731 gsl_vector_free(dy);
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);
741}
742
743/*****************************************************************************/
744
746 char *filename,
747 ncd_t *ncd) {
748
749 int varid;
750
751 /* Open netCDF file... */
752 printf("Read netCDF file: %s\n", filename);
753 NC(nc_open(filename, NC_WRITE, &ncd->ncid));
754
755 /* Read Level-1 data... */
756 NC(nc_inq_varid(ncd->ncid, "l1_time", &varid));
757 NC(nc_get_var_double(ncd->ncid, varid, ncd->l1_time[0]));
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));
763 NC(nc_get_var_double(ncd->ncid, varid, ncd->l1_sat_z));
764 NC(nc_inq_varid(ncd->ncid, "l1_sat_lon", &varid));
765 NC(nc_get_var_double(ncd->ncid, varid, ncd->l1_sat_lon));
766 NC(nc_inq_varid(ncd->ncid, "l1_sat_lat", &varid));
767 NC(nc_get_var_double(ncd->ncid, varid, ncd->l1_sat_lat));
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]));
772
773 /* Read Level-2 data... */
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]));
780}
781
782/*****************************************************************************/
783
785 char *filename,
786 ncd_t *ncd) {
787
788 int dimid[10], p_id, t_id, z_id;
789
790 /* Create netCDF file... */
791 printf("Write netCDF file: %s\n", filename);
792
793 /* Read existing dimensions... */
794 NC(nc_inq_dimid(ncd->ncid, "L1_NTRACK", &dimid[0]));
795 NC(nc_inq_dimid(ncd->ncid, "L1_NXTRACK", &dimid[1]));
796
797 /* Set define mode... */
798 NC(nc_redef(ncd->ncid));
799
800 /* Set new dimensions... */
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]));
803
804 /* Set new variables... */
805 add_var(ncd->ncid, "ret_z", "km", "altitude", NC_FLOAT, &dimid[2], &z_id,
806 1);
807 add_var(ncd->ncid, "ret_press", "hPa", "pressure", NC_FLOAT, dimid, &p_id,
808 2);
809 add_var(ncd->ncid, "ret_temp", "K", "temperature", NC_FLOAT, dimid, &t_id,
810 3);
811
812 /* Leave define mode... */
813 NC(nc_enddef(ncd->ncid));
814
815 /* Write data... */
816 NC(nc_put_var_float(ncd->ncid, z_id, ncd->ret_z));
817 NC(nc_put_var_float(ncd->ncid, p_id, ncd->ret_p));
818 NC(nc_put_var_float(ncd->ncid, t_id, ncd->ret_t));
819
820 /* Close netCDF file... */
821 NC(nc_close(ncd->ncid));
822}
int main(int argc, char *argv[])
Definition: retrieval.c:185
#define NC(cmd)
Execute netCDF library command and check result.
Definition: retrieval.c:36
void fill_gaps(double x[L2_NTRACK][L2_NXTRACK][L2_NLAY], double cx, double cy)
Fill data gaps in L2 data.
Definition: retrieval.c:458
#define L1_NXTRACK
Across-track size of AIRS radiance granule (don't change).
Definition: retrieval.c:53
#define L1_NTRACK
Along-track size of AIRS radiance granule (don't change).
Definition: retrieval.c:50
#define L2_NXTRACK
Across-track size of AIRS retrieval granule (don't change).
Definition: retrieval.c:62
void write_nc(char *filename, ncd_t *ncd)
Write to netCDF file...
Definition: retrieval.c:784
void read_nc(char *filename, ncd_t *ncd)
Read netCDF file.
Definition: retrieval.c:745
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:408
#define L2_NLAY
Number of AIRS pressure layers (don't change).
Definition: retrieval.c:56
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.
Definition: retrieval.c:560
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:435
#define L1_NCHAN
Number of AIRS radiance channels (don't change).
Definition: retrieval.c:47
#define L2_NTRACK
Along-track size of AIRS 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 AIRS Level-2 data.
Definition: retrieval.c:503
Buffer for netCDF data.
Definition: diff_apr.c:68
double l1_lat[L1_NTRACK][L1_NXTRACK]
Footprint latitude [deg].
Definition: diff_apr.c:83
float ret_p[L1_NTRACK *L1_NXTRACK]
Pressure [hPa].
Definition: diff_apr.c:113
double l1_lon[L1_NTRACK][L1_NXTRACK]
Footprint longitude [deg].
Definition: diff_apr.c:80
double l2_z[L2_NTRACK][L2_NXTRACK][L2_NLAY]
Altitude [km].
Definition: diff_apr.c:101
double l1_sat_lat[L1_NTRACK]
Satellite latitude [deg].
Definition: diff_apr.c:92
int ncid
NetCDF file ID.
Definition: diff_apr.c:71
int np
Number of retrieval altitudes.
Definition: diff_apr.c:74
double l1_sat_lon[L1_NTRACK]
Satellite longitude [deg].
Definition: diff_apr.c:89
float ret_z[NP]
Altitude [km].
Definition: diff_apr.c:110
double l2_p[L2_NLAY]
Pressure [hPa].
Definition: diff_apr.c:104
float ret_t[L1_NTRACK *L1_NXTRACK *NP]
Temperature [K].
Definition: diff_apr.c:116
double l2_t[L2_NTRACK][L2_NXTRACK][L2_NLAY]
Temperature [K].
Definition: diff_apr.c:107
double l1_nu[L1_NCHAN]
Channel frequencies [cm^-1].
Definition: diff_apr.c:95
double l1_sat_z[L1_NTRACK]
Satellite altitude [km].
Definition: diff_apr.c:86
float l1_rad[L1_NTRACK][L1_NXTRACK][L1_NCHAN]
Radiance [W/(m^2 sr cm^-1)].
Definition: diff_apr.c:98
double l1_time[L1_NTRACK][L1_NXTRACK]
Time (seconds since 2000-01-01T00:00Z).
Definition: diff_apr.c:77