CrIS Code Collection
libcris.h
Go to the documentation of this file.
1
68#include <netcdf.h>
69#include <gsl/gsl_randist.h>
70#include <gsl/gsl_fft_complex.h>
71#include <gsl/gsl_multifit.h>
72#include <gsl/gsl_poly.h>
73#include <gsl/gsl_sort.h>
74#include <gsl/gsl_spline.h>
75#include "jurassic.h"
76
77/* ------------------------------------------------------------
78 Dimensions...
79 ------------------------------------------------------------ */
80
82#define NDS 200000
83
85#define NPG 30
86
88#define L1_NTRACK 45
89
91#define L1_NXTRACK 30
92
94#define L1_NFOV 9
95
97#define L1_NCHAN_LW 717
98
100#define L1_NCHAN_MW 869
101
103#define L1_NCHAN_SW 637
104
106#define PERT_NTRACK 44000
107
109#define PERT_NXTRACK 120
110
112#define PERT_NFOV 9
113
115#define WX 300
116
118#define WY 33000
119
121#define PMAX 512
122
123/* ------------------------------------------------------------
124 Macros...
125 ------------------------------------------------------------ */
126
128#define NC(cmd) { \
129 int nc_result=(cmd); \
130 if(nc_result!=NC_NOERR) \
131 ERRMSG("%s", nc_strerror(nc_result)); \
132}
133
134/* ------------------------------------------------------------
135 Structs...
136 ------------------------------------------------------------ */
137
139typedef struct {
140
142 double time[L1_NTRACK][L1_NXTRACK];
143
146
149
151 double sat_z[L1_NTRACK];
152
154 double sat_lon[L1_NTRACK];
155
157 double sat_lat[L1_NTRACK];
158
160 double nu_lw[L1_NCHAN_LW];
161
164
166 float nedn_lw[L1_NFOV][L1_NCHAN_LW];
167
169 short qual_lw[L1_NTRACK][L1_NXTRACK][L1_NFOV];
170
172 double nu_mw[L1_NCHAN_MW];
173
176
178 float nedn_mw[L1_NFOV][L1_NCHAN_MW];
179
181 short qual_mw[L1_NTRACK][L1_NXTRACK][L1_NFOV];
182
184 double nu_sw[L1_NCHAN_SW];
185
188
190 float nedn_sw[L1_NFOV][L1_NCHAN_SW];
191
193 short qual_sw[L1_NTRACK][L1_NXTRACK][L1_NFOV];
194
195} cris_l1_t;
196
198typedef struct {
199
202
205
207 int nfov;
208
211
214
217
220
223
226
229
230} pert_t;
231
233typedef struct {
234
236 int nds;
237
239 int np;
240
242 double time[NDS][NPG];
243
245 double z[NDS][NPG];
246
248 double lon[NDS][NPG];
249
251 double lat[NDS][NPG];
252
254 double p[NDS][NPG];
255
257 double t[NDS][NPG];
258
260 double t_apr[NDS][NPG];
261
263 double t_tot[NDS][NPG];
264
266 double t_noise[NDS][NPG];
267
269 double t_fm[NDS][NPG];
270
272 double t_cont[NDS][NPG];
273
275 double t_res[NDS][NPG];
276
278 double chisq[NDS];
279
280} retr_t;
281
283typedef struct {
284
286 int nx;
287
289 int ny;
290
292 double time;
293
295 double z;
296
298 double lon[WX][WY];
299
301 double lat[WX][WY];
302
304 double x[WX];
305
307 double y[WY];
308
310 double temp[WX][WY];
311
313 double bg[WX][WY];
314
316 double pt[WX][WY];
317
319 double fit[WX][WY];
320
322 double var[WX][WY];
323
324} wave_t;
325
326/* ------------------------------------------------------------
327 Functions...
328 ------------------------------------------------------------ */
329
331void add_att(
332 const int ncid,
333 const int varid,
334 const char *unit,
335 const char *long_name);
336
338void add_var(
339 const int ncid,
340 const char *varname,
341 const char *unit,
342 const char *longname,
343 int type,
344 int dimid[],
345 int *varid,
346 int ndims);
347
349void background_poly(
350 wave_t * wave,
351 int dim_x,
352 int dim_y);
353
356 const double *xx,
357 double *yy,
358 const int n,
359 const int dim);
360
363 wave_t * wave,
364 int npts_x,
365 int npts_y);
366
369 wave_t * wave);
370
372void create_noise(
373 wave_t * wave,
374 double nedt);
375
377void create_wave(
378 wave_t * wave,
379 double amp,
380 double lx,
381 double ly,
382 double phi,
383 double fwhm);
384
386void fit_wave(
387 wave_t * wave,
388 double amp,
389 double phi,
390 double kx,
391 double ky,
392 double *chisq);
393
395void fft_help(
396 double *fcReal,
397 double *fcImag,
398 int n);
399
401void fft(
402 wave_t * wave,
403 double *Amax,
404 double *phimax,
405 double *lhmax,
406 double *kxmax,
407 double *kymax,
408 double *alphamax,
409 double *betamax,
410 char *filename);
411
413void gauss(
414 wave_t * wave,
415 double fwhm);
416
418void hamming(
419 wave_t * wave,
420 int nit);
421
423void intpol_x(
424 wave_t * wave,
425 int n);
426
428void median(
429 wave_t * wave,
430 int dx);
431
433void merge_y(
434 wave_t * wave1,
435 wave_t * wave2);
436
438void noise(
439 wave_t * wave,
440 double *mu,
441 double *sig);
442
444void noise_pert(
445 pert_t * pert,
446 int track0,
447 int track1,
448 double *mu,
449 double *sig);
450
452void period(
453 wave_t * wave,
454 double lxymax,
455 double dlxy,
456 double *Amax,
457 double *phimax,
458 double *lhmax,
459 double *kxmax,
460 double *kymax,
461 double *alphamax,
462 double *betamax,
463 char *filename);
464
467 pert_t * pert,
468 wave_t * wave,
469 int track0,
470 int track1,
471 int xtrack0,
472 int xtrack1);
473
475int read_cris_l1(
476 char *filename,
477 cris_l1_t * l1,
478 int apo);
479
481void read_pert(
482 char *filename,
483 char *pertname,
484 int dc,
485 pert_t * pert);
486
488void read_retr(
489 char *filename,
490 retr_t * ret);
491
493void read_retr_help(
494 double *help,
495 int nds,
496 int np,
497 double mat[NDS][NPG]);
498
500void read_wave(
501 char *filename,
502 wave_t * wave);
503
505void ret2wave(
506 retr_t * ret,
507 wave_t * wave,
508 int dataset,
509 int ip);
510
512void variance(
513 wave_t * wave,
514 double dh);
515
517void write_wave(
518 char *filename,
519 wave_t * wave);
void hamming(wave_t *wave, int nit)
Apply Hamming filter to perturbations...
Definition: libcris.c:546
void merge_y(wave_t *wave1, wave_t *wave2)
Merge wave structs in y-direction.
Definition: libcris.c:694
void fit_wave(wave_t *wave, double amp, double phi, double kx, double ky, double *chisq)
Evaluate wave fit...
Definition: libcris.c:279
void create_wave(wave_t *wave, double amp, double lx, double ly, double phi, double fwhm)
Add linear wave pattern...
Definition: libcris.c:250
void noise_pert(pert_t *pert, int track0, int track1, double *mu, double *sig)
Estimate noise from perurbations.
Definition: libcris.c:772
void background_poly_help(const double *xx, double *yy, const int n, const int dim)
Get background based on polynomial fits.
Definition: libcris.c:47
#define PERT_NXTRACK
Across-track size of perturbation data.
Definition: libcris.h:109
void fft_help(double *fcReal, double *fcImag, int n)
Calculate 1-D FFT...
Definition: libcris.c:305
#define L1_NCHAN_LW
Number of CrIS longwave radiance channels.
Definition: libcris.h:97
void background_smooth(wave_t *wave, int npts_x, int npts_y)
Smooth background.
Definition: libcris.c:151
#define L1_NXTRACK
Across-track size of CrIS radiance granule.
Definition: libcris.h:91
#define L1_NTRACK
Along-track size of CrIS radiance granule.
Definition: libcris.h:88
void intpol_x(wave_t *wave, int n)
Interpolate to regular grid in x-direction.
Definition: libcris.c:575
void fft(wave_t *wave, double *Amax, double *phimax, double *lhmax, double *kxmax, double *kymax, double *alphamax, double *betamax, char *filename)
Calculate 2-D FFT...
Definition: libcris.c:344
#define L1_NFOV
Number of field of views of CrIS radiance granule.
Definition: libcris.h:94
#define WX
Across-track size of wave analysis data.
Definition: libcris.h:115
void noise(wave_t *wave, double *mu, double *sig)
Estimate noise.
Definition: libcris.c:727
void read_retr_help(double *help, int nds, int np, double mat[NDS][NPG])
Convert array.
Definition: libcris.c:1518
#define PERT_NTRACK
Along-track size of perturbation data.
Definition: libcris.h:106
void ret2wave(retr_t *ret, wave_t *wave, int dataset, int ip)
Convert CrIS retrieval results to wave analysis struct.
Definition: libcris.c:1594
void background_poly(wave_t *wave, int dim_x, int dim_y)
Get background based on polynomial fits.
Definition: libcris.c:101
void create_noise(wave_t *wave, double nedt)
Add noise to perturbations and temperatures...
Definition: libcris.c:229
int read_cris_l1(char *filename, cris_l1_t *l1, int apo)
Read CrIS Level-1 data.
Definition: libcris.c:1051
void variance(wave_t *wave, double dh)
Compute local variance.
Definition: libcris.c:1648
#define NDS
Maximum number of data sets per granule.
Definition: libcris.h:82
void read_retr(char *filename, retr_t *ret)
Read CrIS retrieval data.
Definition: libcris.c:1352
void pert2wave(pert_t *pert, wave_t *wave, int track0, int track1, int xtrack0, int xtrack1)
Convert radiance perturbation data to wave analysis struct.
#define L1_NCHAN_MW
Number of CrIS midwave radiance channels.
Definition: libcris.h:100
#define L1_NCHAN_SW
Number of CrIS shortwave radiance channels.
Definition: libcris.h:103
#define NPG
Maximum number of data points per granule.
Definition: libcris.h:85
#define PERT_NFOV
Number of field of views of perturbation data.
Definition: libcris.h:112
void add_att(const int ncid, const int varid, const char *unit, const char *long_name)
Add variable attributes to netCDF file.
Definition: libcris.c:5
void create_background(wave_t *wave)
Set background...
Definition: libcris.c:204
void read_pert(char *filename, char *pertname, int dc, pert_t *pert)
Read radiance perturbation data.
Definition: libcris.c:1221
#define WY
Along-track size of wave analysis data.
Definition: libcris.h:118
void read_wave(char *filename, wave_t *wave)
Read wave analysis data.
Definition: libcris.c:1533
void gauss(wave_t *wave, double fwhm)
Apply Gaussian filter to perturbations...
Definition: libcris.c:506
void add_var(const int ncid, const char *varname, const char *unit, const char *longname, int type, int dimid[], int *varid, int ndims)
Add variable to netCDF file.
Definition: libcris.c:20
void median(wave_t *wave, int dx)
Apply median filter to perturbations...
Definition: libcris.c:653
void write_wave(char *filename, wave_t *wave)
Write wave analysis data.
Definition: libcris.c:1698
void period(wave_t *wave, double lxymax, double dlxy, double *Amax, double *phimax, double *lhmax, double *kxmax, double *kymax, double *alphamax, double *betamax, char *filename)
Compute periodogram.
Definition: libcris.c:819
CrIS Level-1 data.
Definition: libcris.h:139
Perturbation data.
Definition: libcris.h:198
int nfov
Number of field of views.
Definition: libcris.h:207
int ntrack
Number of along-track values.
Definition: libcris.h:201
int nxtrack
Number of across-track values.
Definition: libcris.h:204
Retrieval results.
Definition: libcris.h:233
int np
Number of data points.
Definition: libcris.h:239
int nds
Number of data sets.
Definition: libcris.h:236
Wave analysis data.
Definition: libcris.h:283
int nx
Number of across-track values.
Definition: libcris.h:286
int ny
Number of along-track values.
Definition: libcris.h:289
double z
Altitude [km].
Definition: libcris.h:295
double time
Time (seconds since 2000-01-01T00:00Z).
Definition: libcris.h:292