CrIS Code Collection
libcris.h
Go to the documentation of this file.
1#include <netcdf.h>
2#include <gsl/gsl_randist.h>
3#include <gsl/gsl_fft_complex.h>
4#include <gsl/gsl_multifit.h>
5#include <gsl/gsl_poly.h>
6#include <gsl/gsl_sort.h>
7#include <gsl/gsl_spline.h>
8#include "jurassic.h"
9
10/* ------------------------------------------------------------
11 Dimensions...
12 ------------------------------------------------------------ */
13
15#define NDS 200000
16
18#define NPG 30
19
21#define L1_NTRACK 45
22
24#define L1_NXTRACK 30
25
27#define L1_NFOV 9
28
30#define L1_NCHAN_LW 717
31
33#define L1_NCHAN_MW 869
34
36#define L1_NCHAN_SW 637
37
39#define PERT_NTRACK 44000
40
42#define PERT_NXTRACK 120
43
45#define PERT_NFOV 9
46
48#define WX 300
49
51#define WY 33000
52
54#define PMAX 512
55
56/* ------------------------------------------------------------
57 Macros...
58 ------------------------------------------------------------ */
59
61#define NC(cmd) { \
62 int nc_result=(cmd); \
63 if(nc_result!=NC_NOERR) \
64 ERRMSG("%s", nc_strerror(nc_result)); \
65}
66
67/* ------------------------------------------------------------
68 Structs...
69 ------------------------------------------------------------ */
70
72typedef struct {
73
75 double time[L1_NTRACK][L1_NXTRACK];
76
79
82
84 double sat_z[L1_NTRACK];
85
87 double sat_lon[L1_NTRACK];
88
90 double sat_lat[L1_NTRACK];
91
93 double nu_lw[L1_NCHAN_LW];
94
97
99 float nedn_lw[L1_NFOV][L1_NCHAN_LW];
100
102 short qual_lw[L1_NTRACK][L1_NXTRACK][L1_NFOV];
103
105 double nu_mw[L1_NCHAN_MW];
106
109
111 float nedn_mw[L1_NFOV][L1_NCHAN_MW];
112
114 short qual_mw[L1_NTRACK][L1_NXTRACK][L1_NFOV];
115
117 double nu_sw[L1_NCHAN_SW];
118
121
123 float nedn_sw[L1_NFOV][L1_NCHAN_SW];
124
126 short qual_sw[L1_NTRACK][L1_NXTRACK][L1_NFOV];
127
128} cris_l1_t;
129
131typedef struct {
132
135
138
140 int nfov;
141
144
147
150
153
156
159
162
163} pert_t;
164
166typedef struct {
167
169 int nds;
170
172 int np;
173
175 double time[NDS][NPG];
176
178 double z[NDS][NPG];
179
181 double lon[NDS][NPG];
182
184 double lat[NDS][NPG];
185
187 double p[NDS][NPG];
188
190 double t[NDS][NPG];
191
193 double t_apr[NDS][NPG];
194
196 double t_tot[NDS][NPG];
197
199 double t_noise[NDS][NPG];
200
202 double t_fm[NDS][NPG];
203
205 double t_cont[NDS][NPG];
206
208 double t_res[NDS][NPG];
209
211 double chisq[NDS];
212
213} retr_t;
214
216typedef struct {
217
219 int nx;
220
222 int ny;
223
225 double time;
226
228 double z;
229
231 double lon[WX][WY];
232
234 double lat[WX][WY];
235
237 double x[WX];
238
240 double y[WY];
241
243 double temp[WX][WY];
244
246 double bg[WX][WY];
247
249 double pt[WX][WY];
250
252 double fit[WX][WY];
253
255 double var[WX][WY];
256
257} wave_t;
258
259/* ------------------------------------------------------------
260 Functions...
261 ------------------------------------------------------------ */
262
264void add_att(
265 const int ncid,
266 const int varid,
267 const char *unit,
268 const char *long_name);
269
271void add_var(
272 const int ncid,
273 const char *varname,
274 const char *unit,
275 const char *longname,
276 int type,
277 int dimid[],
278 int *varid,
279 int ndims);
280
282void background_poly(
283 wave_t * wave,
284 int dim_x,
285 int dim_y);
286
289 const double *xx,
290 double *yy,
291 const int n,
292 const int dim);
293
296 wave_t * wave,
297 int npts_x,
298 int npts_y);
299
302 wave_t * wave);
303
305void create_noise(
306 wave_t * wave,
307 double nedt);
308
310void create_wave(
311 wave_t * wave,
312 double amp,
313 double lx,
314 double ly,
315 double phi,
316 double fwhm);
317
319void fit_wave(
320 wave_t * wave,
321 double amp,
322 double phi,
323 double kx,
324 double ky,
325 double *chisq);
326
328void fft_help(
329 double *fcReal,
330 double *fcImag,
331 int n);
332
334void fft(
335 wave_t * wave,
336 double *Amax,
337 double *phimax,
338 double *lhmax,
339 double *kxmax,
340 double *kymax,
341 double *alphamax,
342 double *betamax,
343 char *filename);
344
346void gauss(
347 wave_t * wave,
348 double fwhm);
349
351void hamming(
352 wave_t * wave,
353 int nit);
354
356void intpol_x(
357 wave_t * wave,
358 int n);
359
361void median(
362 wave_t * wave,
363 int dx);
364
366void merge_y(
367 wave_t * wave1,
368 wave_t * wave2);
369
371void noise(
372 wave_t * wave,
373 double *mu,
374 double *sig);
375
377void noise_pert(
378 pert_t * pert,
379 int track0,
380 int track1,
381 double *mu,
382 double *sig);
383
385void period(
386 wave_t * wave,
387 double lxymax,
388 double dlxy,
389 double *Amax,
390 double *phimax,
391 double *lhmax,
392 double *kxmax,
393 double *kymax,
394 double *alphamax,
395 double *betamax,
396 char *filename);
397
400 pert_t * pert,
401 wave_t * wave,
402 int track0,
403 int track1,
404 int xtrack0,
405 int xtrack1);
406
408int read_cris_l1(
409 char *filename,
410 cris_l1_t * l1,
411 int apo);
412
414void read_pert(
415 char *filename,
416 char *pertname,
417 int dc,
418 pert_t * pert);
419
421void read_retr(
422 char *filename,
423 retr_t * ret);
424
426void read_retr_help(
427 double *help,
428 int nds,
429 int np,
430 double mat[NDS][NPG]);
431
433void read_wave(
434 char *filename,
435 wave_t * wave);
436
438void ret2wave(
439 retr_t * ret,
440 wave_t * wave,
441 int dataset,
442 int ip);
443
445void variance(
446 wave_t * wave,
447 double dh);
448
450void write_wave(
451 char *filename,
452 wave_t * wave);
JURASSIC library declarations.
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:42
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:30
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:24
#define L1_NTRACK
Along-track size of CrIS radiance granule.
Definition: libcris.h:21
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:27
#define WX
Across-track size of wave analysis data.
Definition: libcris.h:48
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:39
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:15
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:33
#define L1_NCHAN_SW
Number of CrIS shortwave radiance channels.
Definition: libcris.h:36
#define NPG
Maximum number of data points per granule.
Definition: libcris.h:18
#define PERT_NFOV
Number of field of views of perturbation data.
Definition: libcris.h:45
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:51
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:72
Perturbation data.
Definition: libcris.h:131
int nfov
Number of field of views.
Definition: libcris.h:140
int ntrack
Number of along-track values.
Definition: libcris.h:134
int nxtrack
Number of across-track values.
Definition: libcris.h:137
Retrieval results.
Definition: libcris.h:166
int np
Number of data points.
Definition: libcris.h:172
int nds
Number of data sets.
Definition: libcris.h:169
Wave analysis data.
Definition: libcris.h:216
int nx
Number of across-track values.
Definition: libcris.h:219
int ny
Number of along-track values.
Definition: libcris.h:222
double z
Altitude [km].
Definition: libcris.h:228
double time
Time (seconds since 2000-01-01T00:00Z).
Definition: libcris.h:225