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
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];
187
188 double chisq, chisq_min, chisq_max, chisq_mean, z[NP];
189
190 int channel[ND], m, ntask = -1, rank, size;
191
192 /* ------------------------------------------------------------
193 Init...
194 ------------------------------------------------------------ */
195
196 /* MPI... */
197 MPI_Init(&argc, &argv);
198 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
199 MPI_Comm_size(MPI_COMM_WORLD, &size);
200
201 /* Measure CPU time... */
202 TIMER("total", 1);
203
204 /* Check arguments... */
205 if (argc < 3)
206 ERRMSG("Give parameters: <ctl> <filelist>");
207
208 /* Read control parameters... */
209 read_ctl(argc, argv, &ctl);
210 read_ret(argc, argv, &ctl, &ret);
211
212 /* Initialize look-up tables... */
213 tbl_t *tbl = read_tbl(&ctl);
214
215 /* Read retrieval grid... */
216 const int nz = (int) scan_ctl(argc, argv, "NZ", -1, "", NULL);
217 if (nz > NP)
218 ERRMSG("Too many altitudes!");
219 for (int iz = 0; iz < nz; iz++)
220 z[iz] = scan_ctl(argc, argv, "Z", iz, "", NULL);
221
222 /* Read track range... */
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);
225
226 /* Read xtrack range... */
227 const int xtrack0 = (int) scan_ctl(argc, argv, "XTRACK_MIN", -1, "0", NULL);
228 const int xtrack1 =
229 (int) scan_ctl(argc, argv, "XTRACK_MAX", -1, "89", NULL);
230
231 /* Read height range... */
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);
235
236 /* Background smoothing... */
237 const double sx = scan_ctl(argc, argv, "SX", -1, "8", NULL);
238 const double sy = scan_ctl(argc, argv, "SY", -1, "2", NULL);
239
240 /* SZA threshold... */
241 const double sza_thresh = scan_ctl(argc, argv, "SZA", -1, "96", NULL);
242
243 /* ------------------------------------------------------------
244 Distribute granules...
245 ------------------------------------------------------------ */
246
247 /* Open filelist... */
248 printf("Read filelist: %s\n", argv[2]);
249 if (!(in = fopen(argv[2], "r")))
250 ERRMSG("Cannot open filelist!");
251
252 /* Loop over netCDF files... */
253 while (fscanf(in, "%s", filename) != EOF) {
254
255 /* Distribute files with MPI... */
256 if ((++ntask) % size != rank)
257 continue;
258
259 /* Write info... */
260 printf("Retrieve file %s on rank %d of %d (with %d threads)...\n",
261 filename, rank + 1, size, omp_get_max_threads());
262
263 /* ------------------------------------------------------------
264 Initialize retrieval...
265 ------------------------------------------------------------ */
266
267 /* Read netCDF file... */
268 read_nc(filename, &ncd);
269
270 /* Identify radiance channels... */
271 for (int id = 0; id < ctl.nd; id++) {
272 channel[id] = -999;
273 for (int i = 0; i < L1_NCHAN; i++)
274 if (fabs(ctl.nu[id] - ncd.l1_nu[i]) < 0.1)
275 channel[id] = i;
276 if (channel[id] < 0)
277 ERRMSG("Cannot identify radiance channel!");
278 }
279
280 /* Fill data gaps... */
281 fill_gaps(ncd.l2_t, sx, sy);
282 fill_gaps(ncd.l2_z, sx, sy);
283
284 /* Set climatological data for center of granule... */
285 atm_clim.np = nz;
286 for (int iz = 0; iz < nz; iz++)
287 atm_clim.z[iz] = z[iz];
288 climatology(&ctl, &atm_clim);
289
290 /* ------------------------------------------------------------
291 Retrieval...
292 ------------------------------------------------------------ */
293
294 /* Get chi^2 statistics... */
295 chisq_min = 1e100;
296 chisq_max = -1e100;
297 chisq_mean = 0;
298 m = 0;
299
300 /* Loop over swaths... */
301 for (int track = track0; track <= track1; track++) {
302
303 /* Measure CPU time... */
304 TIMER("retrieval", 1);
305
306 /* Loop over scan... */
307 for (int xtrack = xtrack0; xtrack <= xtrack1; xtrack++) {
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 (int 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 (int 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 (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];
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 /* Get chi^2 statistics... */
343 if (gsl_finite(chisq)) {
344 chisq_min = GSL_MIN(chisq_min, chisq);
345 chisq_max = GSL_MAX(chisq_max, chisq);
346 chisq_mean += chisq;
347 m++;
348 }
349
350 /* Buffer results... */
351 buffer_nc(&atm_i, chisq, &ncd, track, xtrack, np0, np1);
352 }
353
354 /* Measure CPU time... */
355 TIMER("retrieval", 3);
356 }
357
358 /* ------------------------------------------------------------
359 Finalize...
360 ------------------------------------------------------------ */
361
362 /* Write netCDF file... */
363 write_nc(filename, &ncd);
364
365 /* Write info... */
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);
369 }
370
371 /* Close file list... */
372 fclose(in);
373
374 /* Measure CPU time... */
375 TIMER("total", 3);
376
377 /* Report memory usage... */
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.);
384
385 /* Report problem size... */
386 printf("SIZE_TASKS = %d\n", size);
387 printf("SIZE_THREADS = %d\n", omp_get_max_threads());
388
389 /* MPI... */
390 MPI_Finalize();
391
392 return EXIT_SUCCESS;
393}
394
395/*****************************************************************************/
396
398 int ncid,
399 const char *varname,
400 const char *unit,
401 const char *longname,
402 int type,
403 int dimid[],
404 int *varid,
405 int ndims) {
406
407 /* Check if variable exists... */
408 if (nc_inq_varid(ncid, varname, varid) != NC_NOERR) {
409
410 /* Define variable... */
411 NC(nc_def_var(ncid, varname, type, ndims, dimid, varid));
412
413 /* Set long name... */
414 NC(nc_put_att_text
415 (ncid, *varid, "long_name", strlen(longname), longname));
416
417 /* Set units... */
418 NC(nc_put_att_text(ncid, *varid, "units", strlen(unit), unit));
419 }
420}
421
422/*****************************************************************************/
423
425 atm_t *atm,
426 double chisq,
427 ncd_t *ncd,
428 int track,
429 int xtrack,
430 int np0,
431 int np1) {
432
433 /* Set number of data points... */
434 ncd->np = np1 - np0 + 1;
435
436 /* Save retrieval data... */
437 for (int ip = np0; ip <= np1; ip++) {
438 ncd->ret_z[ip - np0] = (float) atm->z[ip];
439 ncd->ret_p[track * L1_NXTRACK + xtrack] = (float) atm->p[np0];
440 ncd->ret_t[(track * L1_NXTRACK + xtrack) * ncd->np + ip - np0] =
441 (gsl_finite(chisq) ? (float) atm->t[ip] : GSL_NAN);
442 }
443}
444
445/************************************************************************/
446
448 double x[L2_NTRACK][L2_NXTRACK][L2_NLAY],
449 double cx,
450 double cy) {
451
452 double help[L2_NTRACK][L2_NXTRACK];
453
454 /* Loop over layers... */
455 for (int lay = 0; lay < L2_NLAY; lay++) {
456
457 /* Loop over grid points... */
458 for (int track = 0; track < L2_NTRACK; track++)
459 for (int xtrack = 0; xtrack < L2_NXTRACK; xtrack++) {
460
461 /* Init... */
462 help[track][xtrack] = 0;
463 double wsum = 0;
464
465 /* Averrage data points... */
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];
473 wsum += w;
474 }
475
476 /* Normalize... */
477 if (wsum > 0)
478 help[track][xtrack] /= wsum;
479 else
480 help[track][xtrack] = GSL_NAN;
481 }
482
483 /* Copy grid points... */
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];
487 }
488}
489
490/************************************************************************/
491
493 ncd_t *ncd,
494 int track,
495 int xtrack,
496 ctl_t *ctl,
497 atm_t *atm) {
498
499 static atm_t atm_airs;
500
501 double k[NW], p, q[NG], t, zmax = 0, zmin = 1000;
502
503 /* Reset track- and xtrack-index to match Level-2 data... */
504 track /= 3;
505 xtrack /= 3;
506
507 /* Store AIRS data in atmospheric data struct... */
508 atm_airs.np = 0;
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!");
516 }
517
518 /* Check number of levels... */
519 if (atm_airs.np <= 0)
520 return;
521
522 /* Get height range of AIRS data... */
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]);
526 }
527
528 /* Merge AIRS data... */
529 for (int ip = 0; ip < atm->np; ip++) {
530
531 /* Interpolate AIRS data... */
532 intpol_atm(ctl, &atm_airs, atm->z[ip], &p, &t, q, k);
533
534 /* Weighting factor... */
535 double w = 1;
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);
540
541 /* Merge... */
542 atm->t[ip] = w * t + (1 - w) * atm->t[ip];
543 atm->p[ip] = w * p + (1 - w) * atm->p[ip];
544 }
545}
546
547/*****************************************************************************/
548
550 char *filename,
551 ncd_t *ncd) {
552
553 int varid;
554
555 /* Open netCDF file... */
556 printf("Read netCDF file: %s\n", filename);
557 NC(nc_open(filename, NC_WRITE, &ncd->ncid));
558
559 /* Read Level-1 data... */
560 NC(nc_inq_varid(ncd->ncid, "l1_time", &varid));
561 NC(nc_get_var_double(ncd->ncid, varid, ncd->l1_time[0]));
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));
567 NC(nc_get_var_double(ncd->ncid, varid, ncd->l1_sat_z));
568 NC(nc_inq_varid(ncd->ncid, "l1_sat_lon", &varid));
569 NC(nc_get_var_double(ncd->ncid, varid, ncd->l1_sat_lon));
570 NC(nc_inq_varid(ncd->ncid, "l1_sat_lat", &varid));
571 NC(nc_get_var_double(ncd->ncid, varid, ncd->l1_sat_lat));
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]));
576
577 /* Read Level-2 data... */
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]));
584}
585
586/*****************************************************************************/
587
589 char *filename,
590 ncd_t *ncd) {
591
592 int dimid[10], p_id, t_id, z_id;
593
594 /* Create netCDF file... */
595 printf("Write netCDF file: %s\n", filename);
596
597 /* Read existing dimensions... */
598 NC(nc_inq_dimid(ncd->ncid, "L1_NTRACK", &dimid[0]));
599 NC(nc_inq_dimid(ncd->ncid, "L1_NXTRACK", &dimid[1]));
600
601 /* Set define mode... */
602 NC(nc_redef(ncd->ncid));
603
604 /* Set new dimensions... */
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]));
607
608 /* Set new variables... */
609 add_var(ncd->ncid, "ret_z", "km", "altitude", NC_FLOAT, &dimid[2], &z_id,
610 1);
611 add_var(ncd->ncid, "ret_press", "hPa", "pressure", NC_FLOAT, dimid, &p_id,
612 2);
613 add_var(ncd->ncid, "ret_temp", "K", "temperature", NC_FLOAT, dimid, &t_id,
614 3);
615
616 /* Leave define mode... */
617 NC(nc_enddef(ncd->ncid));
618
619 /* Write data... */
620 NC(nc_put_var_float(ncd->ncid, z_id, ncd->ret_z));
621 NC(nc_put_var_float(ncd->ncid, p_id, ncd->ret_p));
622 NC(nc_put_var_float(ncd->ncid, t_id, ncd->ret_t));
623
624 /* Close netCDF file... */
625 NC(nc_close(ncd->ncid));
626}
int main(int argc, char *argv[])
Definition: retrieval.c:174
#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:447
#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:588
void read_nc(char *filename, ncd_t *ncd)
Read netCDF file.
Definition: retrieval.c:549
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:397
#define L2_NLAY
Number of AIRS 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:424
#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:492
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