AIRS Code Collection
Functions
var1d.c File Reference

Estimate horizontal wavelength sensitivity. More...

#include "libairs.h"

Go to the source code of this file.

Functions

int main (int argc, char *argv[])
 

Detailed Description

Estimate horizontal wavelength sensitivity.

Definition in file var1d.c.

Function Documentation

◆ main()

int main ( int  argc,
char *  argv[] 
)

Definition at line 28 of file var1d.c.

30 {
31
32 static double chisq;
33
34 /* Check arguments... */
35 if (argc != 8)
36 ERRMSG("Give parameters: <width> <n> <lxmin> <lxmax> <dlx> <fwhm> <dim>");
37
38 /* Get arguments... */
39 const double width = atof(argv[1]);
40 const int n = atoi(argv[2]);
41 const double lxmin = atof(argv[3]);
42 const double lxmax = atof(argv[4]);
43 const double dlx = atoi(argv[5]);
44 const double fwhm = atof(argv[6]);
45 const int dim = atoi(argv[7]);
46
47 /* Initialize... */
48 gsl_vector *c = gsl_vector_alloc((size_t) dim);
49 gsl_matrix *cov = gsl_matrix_alloc((size_t) dim, (size_t) dim);
50 gsl_multifit_linear_workspace *work
51 = gsl_multifit_linear_alloc((size_t) n, (size_t) dim);
52 gsl_matrix *X = gsl_matrix_alloc((size_t) n, (size_t) dim);
53 gsl_vector *xvec = gsl_vector_alloc((size_t) n);
54 gsl_vector *yvec = gsl_vector_alloc((size_t) n);
55 gsl_vector *yfit = gsl_vector_alloc((size_t) n);
56
57 /* Loop over wavelengths... */
58 for (double lx = lxmin; lx <= lxmax; lx += dlx) {
59
60 /* Initialize... */
61 double vmean = 0;
62 double vmean2 = 0;
63
64 /* Loop over phases... */
65 for (double phi = 0; phi < 2 * M_PI; phi += M_PI / 180) {
66
67 /* Initialize... */
68 double var = 0, var2 = 0, wsum = 0;
69
70 /* Set wave... */
71 for (int i = 0; i < n; i++) {
72 gsl_vector_set(xvec, (size_t) i, width / (n - 1.0) * i - width / 2.);
73 gsl_vector_set(yvec, (size_t) i,
74 sin(2 * M_PI / lx * gsl_vector_get(xvec, (size_t) i) +
75 phi));
76 if (fwhm > 0) {
77 double w = gsl_ran_gaussian_pdf(gsl_vector_get(xvec, (size_t) i),
78 fwhm * lx / 2.3548);
79 gsl_vector_set(yvec, (size_t) i,
80 w * gsl_vector_get(yvec, (size_t) i));
81 wsum += w;
82 }
83 }
84 if (wsum > 0)
85 gsl_vector_scale(yvec, 1 / wsum);
86
87 /* Detrending... */
88 for (int i = 0; i < n; i++)
89 for (int i2 = 0; i2 < dim; i2++)
90 gsl_matrix_set(X, (size_t) i, (size_t) i2,
91 pow(gsl_vector_get(xvec, (size_t) i), 1. * i2));
92 gsl_multifit_linear(X, yvec, c, cov, &chisq, work);
93 for (int i = 0; i < n; i++)
94 gsl_vector_set(yfit, (size_t) i, gsl_vector_get(yvec, (size_t) i)
95 - gsl_poly_eval(c->data, (int) dim,
96 gsl_vector_get(xvec, (size_t) i)));
97
98 /* Compute variances... */
99 for (int i = 0; i < n; i++) {
100 var += gsl_pow_2(gsl_vector_get(yfit, (size_t) i)) / (double) n;
101 var2 += gsl_pow_2(gsl_vector_get(yvec, (size_t) i)) / (double) n;
102 }
103 vmean += var;
104 vmean2 += var2;
105 }
106
107 /* Write output... */
108 printf("%g %g\n", lx, 100 * vmean / vmean2);
109 }
110
111 return EXIT_SUCCESS;
112}
#define ERRMSG(...)
Print error message and quit program.
Definition: jurassic.h:237