JURASSIC
Macros | Functions
invert.c File Reference

Inversion tool for MPTRAC. More...

#include "jurassic.h"
#include <gsl/gsl_fit.h>

Go to the source code of this file.

Macros

#define NLMAX   30000000
 Maximum number of data lines... More...
 
#define NMAX   1000
 Maximum number of data points... More...
 

Functions

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

Detailed Description

Inversion tool for MPTRAC.

Definition in file invert.c.

Macro Definition Documentation

◆ NLMAX

#define NLMAX   30000000

Maximum number of data lines...

Definition at line 33 of file invert.c.

◆ NMAX

#define NMAX   1000

Maximum number of data points...

Definition at line 36 of file invert.c.

Function Documentation

◆ main()

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

Definition at line 42 of file invert.c.

44 {
45
46 static ctl_t ctl;
47
48 static atm_t atm, atm2;
49
50 static obs_t obs;
51
52 static gsl_matrix *k;
53
54 static FILE *in, *out;
55
56 static char line[LEN];
57
58 static double rtime[NLMAX], rz[NLMAX], rlon[NLMAX], rlat[NLMAX], obs_meas,
59 obs_sim, scl = 1.0, scl_err, c0, c1, cov00, cov01, cov11, sumsq,
60 x[NMAX], x2[NMAX], y[NMAX], y_err[NMAX], y2[NMAX], y2_err[NMAX],
61 y2_sim[NMAX], y2_sim_err[NMAX], w2[NMAX], dt, tol, obs_err;
62
63 static float rp[NLMAX], rt[NLMAX], rso2[NLMAX], rh2o[NLMAX],
64 ro3[NLMAX], robs[NLMAX];
65
66 static int data, fit, i, ig, il, ip, it, itmax, n, nl, ndata[NMAX], nprof;
67
68 static size_t mk, nk;
69
70 /* Check arguments... */
71 if (argc < 6)
72 ERRMSG("Give parameters: <ctl> <prof> <inv> <atm> <kernel>");
73
74 /* Read control parameters... */
75 read_ctl(argc, argv, &ctl);
76 dt = scan_ctl(argc, argv, "INVERT_DT", -1, "86400", NULL);
77 obs_err = scan_ctl(argc, argv, "INVERT_OBS_ERR", -1, "1.0", NULL);
78 data = (int) scan_ctl(argc, argv, "INVERT_DATA", -1, "2", NULL);
79 fit = (int) scan_ctl(argc, argv, "INVERT_FIT", -1, "3", NULL);
80 itmax = (int) scan_ctl(argc, argv, "INVERT_ITMAX", -1, "10", NULL);
81 tol = scan_ctl(argc, argv, "INVERT_TOL", -1, "1e-4", NULL);
82
83 /* Check control parameters... */
84 if (ctl.ng != 4)
85 ERRMSG("Set NG = 4!");
86 if (strcmp(ctl.emitter[0], "SO2") != 0)
87 ERRMSG("Set EMITTER[0] = SO2!");
88 if (strcmp(ctl.emitter[1], "H2O") != 0)
89 ERRMSG("Set EMITTER[1] = H2O!");
90 if (strcmp(ctl.emitter[2], "O3") != 0)
91 ERRMSG("Set EMITTER[2] = O3!");
92 if (strcmp(ctl.emitter[3], "CO2") != 0)
93 ERRMSG("Set EMITTER[3] = CO2!");
94 if (ctl.nd != 2)
95 ERRMSG("Set ND = 2!");
96
97 /* Set control parameters... */
98 ctl.write_bbt = 1;
99 ctl.write_matrix = 1;
100
101 /* Set observation data... */
102 obs.nr = 1;
103 obs.obsz[0] = 705;
104
105 /* ------------------------------------------------------------
106 Read profiles...
107 ------------------------------------------------------------ */
108
109 /* Read profile data... */
110 LOG(1, "Read profile data: %s", argv[2]);
111
112 /* Open file... */
113 if (!(in = fopen(argv[2], "r")))
114 ERRMSG("Cannot open file!");
115
116 /* Read file... */
117 while (fgets(line, LEN, in)) {
118
119 /* Read data... */
120 if (sscanf(line, "%lg %lg %lg %lg %g %g %g %g %g %g", &rtime[nl],
121 &rz[nl], &rlon[nl], &rlat[nl], &rp[nl], &rt[nl], &rso2[nl],
122 &rh2o[nl], &ro3[nl], &robs[nl]) != 10)
123 continue;
124 if ((++nl) > NLMAX)
125 ERRMSG("Too many profile data points!");
126 }
127
128 /* Close files... */
129 fclose(in);
130
131 /* ------------------------------------------------------------
132 Fit scaling factor for total mass...
133 ------------------------------------------------------------ */
134
135 /* Iterations... */
136 for (it = 0; it < itmax; it++) {
137
138 /* Init... */
139 atm.np = n = 0;
140 for (i = 0; i < NMAX; i++) {
141 ndata[i] = 0;
142 x[i] = y[i] = NAN;
143 }
144
145 /* Loop over lines... */
146 for (il = 0; il < nl; il++) {
147
148 /* Check for new profile... */
149 if ((rtime[il] != atm.time[0]
150 || rlon[il] != atm.lon[0]
151 || rlat[il] != atm.lat[0])
152 && atm.np > 0) {
153
154 /* Call forward model... */
155 formod(&ctl, &atm, &obs);
156 obs_sim = obs.rad[0][0] - obs.rad[1][0];
157
158 /* Get time index... */
159 i = (int) ((atm.time[0] - rtime[0]) / dt);
160 if (i < 0 && i >= NMAX)
161 ERRMSG("Time index out of range!");
162
163 /* Get maxima... */
164 if (data == 1) {
165 x[i] = (isfinite(x[i]) ? MAX(x[i], obs_sim) : obs_sim);
166 y[i] = (isfinite(y[i]) ? MAX(y[i], obs_meas) : obs_meas);
167 y_err[i] = obs_err;
168 if (isfinite(x[i]) && isfinite(y[i]))
169 ndata[i] = 1;
170 }
171
172 /* Get means... */
173 else if (data == 2) {
174 if (ndata[i] == 0) {
175 x[i] = obs_sim;
176 y[i] = obs_meas;
177 y_err[i] = POW2(obs_meas);
178 } else {
179 x[i] += obs_sim;
180 y[i] += obs_meas;
181 y_err[i] += POW2(obs_meas);
182 }
183 ndata[i]++;
184 }
185
186 /* Calculate mean atmospheric profile... */
187 nprof++;
188 atm2.np = atm.np;
189 for (ip = 0; ip < atm.np; ip++) {
190 atm2.time[ip] += atm.time[ip];
191 atm2.z[ip] += atm.z[ip];
192 atm2.lon[ip] += atm.lon[ip];
193 atm2.lat[ip] += atm.lat[ip];
194 atm2.p[ip] += atm.p[ip];
195 atm2.t[ip] += atm.t[ip];
196 for (ig = 1; ig < ctl.ng; ig++)
197 atm2.q[ig][ip] += atm.q[ig][ip];
198 }
199
200 /* Reset counter... */
201 atm.np = 0;
202 }
203
204 /* Save data... */
205 obs_meas = robs[il];
206 atm.time[atm.np] = rtime[il];
207 atm.z[atm.np] = rz[il];
208 atm.lon[atm.np] = rlon[il];
209 atm.lat[atm.np] = rlat[il];
210 atm.p[atm.np] = rp[il];
211 atm.t[atm.np] = rt[il];
212 atm.q[0][atm.np] = rso2[il] * scl;
213 atm.q[1][atm.np] = rh2o[il];
214 atm.q[2][atm.np] = ro3[il];
215 atm.q[3][atm.np] = 371.789948e-6 + 2.026214e-6
216 * (atm.time[atm.np] - 63158400.) / 31557600.;
217 if ((++atm.np) > NP)
218 ERRMSG("Too many data points!");
219 }
220
221 /* Calculate means... */
222 if (data == 2)
223 for (i = 0; i < NMAX; i++)
224 if (ndata[i] > 0) {
225 x[i] /= ndata[i];
226 y[i] /= ndata[i];
227 y_err[i] = sqrt(MAX(y_err[i] / ndata[i] - POW2(y[i]), 0.0))
228 / sqrt(ndata[i]); /* standard error! */
229 }
230
231 /* Filter data... */
232 n = 0;
233 for (i = 0; i < NMAX; i++)
234 if (ndata[i] > 0 && isfinite(x[i]) && isfinite(y[i])
235 && isfinite(y_err[i])) {
236 x2[n] = x[i];
237 y2[n] = y[i];
238 y2_err[n] = y_err[i];
239 w2[n] = 1. / POW2(y_err[i]);
240 n++;
241 }
242
243 /* Fit radiance data... */
244 if (fit == 1)
245 gsl_fit_mul(x2, 1, y2, 1, (size_t) n, &c1, &cov11, &sumsq);
246 else if (fit == 2)
247 gsl_fit_wmul(x2, 1, w2, 1, y2, 1, (size_t) n, &c1, &cov11, &sumsq);
248 else if (fit == 3)
249 gsl_fit_linear(x2, 1, y2, 1, (size_t) n, &c0, &c1, &cov00, &cov01,
250 &cov11, &sumsq);
251 else if (fit == 4)
252 gsl_fit_wlinear(x2, 1, w2, 1, y2, 1, (size_t) n, &c0, &c1, &cov00,
253 &cov01, &cov11, &sumsq);
254 else
255 ERRMSG("Check INVERT_FIT!");
256
257 /* Get new scaling factor... */
258 double scl_old = scl;
259 scl_err = scl * sqrt(cov11);
260 scl *= c1;
261
262 /* Write info... */
263 LOG(1, " it= %d | scl= %g +/- %g | RMSE= %g",
264 it, scl, scl_err, sqrt(sumsq / n));
265
266 /* Convergence test... */
267 if (fabs(2.0 * (scl - scl_old) / (scl + scl_old)) < tol)
268 break;
269 }
270
271 /* ------------------------------------------------------------
272 Write inversion data...
273 ------------------------------------------------------------ */
274
275 /* Write info... */
276 LOG(1, "Write inversion data: %s", argv[3]);
277
278 /* Create file... */
279 if (!(out = fopen(argv[3], "w")))
280 ERRMSG("Cannot create file!");
281
282 /* Write header... */
283 fprintf(out,
284 "# $1 = time (seconds since 2000-01-01T00:00Z)\n"
285 "# $2 = simulated SO2 index [K]\n"
286 "# $3 = scaled simulated SO2 index [K]\n"
287 "# $4 = error of scaled simulated SO2 index [K]\n"
288 "# $5 = observed SO2 index [K]\n"
289 "# $6 = error of observed SO2 index [K]\n\n");
290
291 /* Write data... */
292 for (i = 0; i < n; i++) {
293
294 /* Calculate scaled SO2 index... */
295 if (fit == 1 || fit == 2)
296 gsl_fit_mul_est(x2[i], c1, cov11, &y2_sim[i], &y2_sim_err[i]);
297 else if (fit == 3 || fit == 4)
298 gsl_fit_linear_est(x2[i], c0, c1, cov00, cov01, cov11, &y2_sim[i],
299 &y2_sim_err[i]);
300
301 /* Write output... */
302 fprintf(out, "%.2f %g %g %g %g %g\n", rtime[0] + (i + 0.5) * dt,
303 x2[i], y2_sim[i], y2_sim_err[i], y2[i], y2_err[i]);
304 }
305
306 /* Report scaling factor for total mass... */
307 fprintf(out, "\n");
308 fprintf(out, "# scl= %g +/- %g\n", scl, scl_err);
309 fprintf(out, "# c1= %g +/- %g\n", c1, sqrt(cov11));
310 if (fit == 3 || fit == 4) {
311 fprintf(out, "# c0= %g +/- %g\n", c0, sqrt(cov00));
312 fprintf(out, "# corr= %g\n", cov01 / (sqrt(cov00) * sqrt(cov11)));
313 }
314 fprintf(out, "# RMSE= %g\n", sqrt(sumsq / n));
315 fprintf(out, "# n= %d\n", n);
316
317 /* Close files... */
318 fclose(out);
319
320 /* ------------------------------------------------------------
321 Calculate kernel...
322 ------------------------------------------------------------ */
323
324 /* Set atmospheric data... */
325 for (ip = 0; ip < atm2.np; ip++) {
326 atm2.time[ip] /= nprof;
327 atm2.z[ip] /= nprof;
328 atm2.lon[ip] /= nprof;
329 atm2.lat[ip] /= nprof;
330 atm2.p[ip] /= nprof;
331 atm2.t[ip] /= nprof;
332 for (ig = 0; ig < ctl.ng; ig++)
333 atm2.q[ig][ip] /= nprof;
334 }
335
336 /* Get sizes... */
337 nk = atm2x(&ctl, &atm2, NULL, NULL, NULL);
338 mk = obs2y(&ctl, &obs, NULL, NULL, NULL);
339
340 /* Allocate... */
341 k = gsl_matrix_alloc(mk, nk);
342
343 /* Compute kernel matrix... */
344 kernel(&ctl, &atm2, &obs, k);
345
346 /* Write atmospheric data... */
347 write_atm(NULL, argv[4], &ctl, &atm);
348
349 /* Write matrix to file... */
350 write_matrix(NULL, argv[5], &ctl, k, &atm2, &obs, "y", "x", "r");
351
352 /* Free... */
353 gsl_matrix_free(k);
354
355 return EXIT_SUCCESS;
356}
#define NMAX
Maximum number of data points...
Definition: invert.c:36
#define NLMAX
Maximum number of data lines...
Definition: invert.c:33
void read_ctl(int argc, char *argv[], ctl_t *ctl)
Read forward model control parameters.
Definition: jurassic.c:4561
void write_atm(const char *dirname, const char *filename, ctl_t *ctl, atm_t *atm)
Write atmospheric data.
Definition: jurassic.c:5361
size_t obs2y(ctl_t *ctl, obs_t *obs, gsl_vector *y, int *ida, int *ira)
Compose measurement vector.
Definition: jurassic.c:4180
void write_matrix(const char *dirname, const char *filename, ctl_t *ctl, gsl_matrix *matrix, atm_t *atm, obs_t *obs, const char *rowspace, const char *colspace, const char *sort)
Write matrix.
Definition: jurassic.c:5519
double scan_ctl(int argc, char *argv[], const char *varname, int arridx, const char *defvalue, char *value)
Search control parameter file for variable entry.
Definition: jurassic.c:5138
void formod(ctl_t *ctl, atm_t *atm, obs_t *obs)
Determine ray paths and compute radiative transfer.
Definition: jurassic.c:3033
void kernel(ctl_t *ctl, atm_t *atm, obs_t *obs, gsl_matrix *k)
Compute Jacobians.
Definition: jurassic.c:4008
size_t atm2x(ctl_t *ctl, atm_t *atm, gsl_vector *x, int *iqa, int *ipa)
Compose state vector or parameter vector.
Definition: jurassic.c:29
#define LEN
Maximum length of ASCII data lines.
Definition: jurassic.h:363
#define POW2(x)
Compute x^2.
Definition: jurassic.h:175
#define ERRMSG(...)
Print error message and quit program.
Definition: jurassic.h:217
#define NP
Maximum number of atmospheric data points.
Definition: jurassic.h:343
#define LOG(level,...)
Print log message.
Definition: jurassic.h:201
#define MAX(a, b)
Macro to determine the maximum of two values.
Definition: jurassic.h:148
Atmospheric data.
Definition: jurassic.h:468
double time[NP]
Time (seconds since 2000-01-01T00:00Z).
Definition: jurassic.h:474
double lat[NP]
Latitude [deg].
Definition: jurassic.h:483
double lon[NP]
Longitude [deg].
Definition: jurassic.h:480
double t[NP]
Temperature [K].
Definition: jurassic.h:489
int np
Number of data points.
Definition: jurassic.h:471
double z[NP]
Altitude [km].
Definition: jurassic.h:477
double q[NG][NP]
Volume mixing ratio [ppv].
Definition: jurassic.h:492
double p[NP]
Pressure [hPa].
Definition: jurassic.h:486
Forward model control parameters.
Definition: jurassic.h:521
int write_matrix
Write matrix file (0=no, 1=yes).
Definition: jurassic.h:641
int ng
Number of emitters.
Definition: jurassic.h:524
int nd
Number of radiance channels.
Definition: jurassic.h:530
char emitter[NG][LEN]
Name of each emitter.
Definition: jurassic.h:527
int write_bbt
Use brightness temperature instead of radiance (0=no, 1=yes).
Definition: jurassic.h:638
Observation geometry and radiance data.
Definition: jurassic.h:714
double rad[ND][NR]
Radiance [W/(m^2 sr cm^-1)].
Definition: jurassic.h:753
double obsz[NR]
Observer altitude [km].
Definition: jurassic.h:723
int nr
Number of ray paths.
Definition: jurassic.h:717
Here is the call graph for this function: