52 {
53
55
56 static atm_t atm, atm2;
57
59
60 static FILE *in, *out;
61
62 static char line[
LEN];
63
65 obs_sim, scl = 1.0, scl_err, c0, c1, cov00, cov01, cov11, sumsq,
68
71
72 static int i, ig, n, nl, ndata[
NMAX], nprof;
73
74
76
77
78 if (argc < 6)
79 ERRMSG(
"Missing or invalid command-line arguments.\n\n"
80 "Usage: invert <ctl> <prof> <inv> <atm> <kernel> [KEY VALUE ...]\n\n"
81 "Use -h for full help.");
82
83
85 const double dt =
scan_ctl(argc, argv,
"INVERT_DT", -1,
"86400", NULL);
86 const double obs_err =
87 scan_ctl(argc, argv,
"INVERT_OBS_ERR", -1,
"1.0", NULL);
88 const int data = (int)
scan_ctl(argc, argv,
"INVERT_DATA", -1,
"2", NULL);
89 const int fit = (int)
scan_ctl(argc, argv,
"INVERT_FIT", -1,
"3", NULL);
90 const int itmax =
91 (int)
scan_ctl(argc, argv,
"INVERT_ITMAX", -1,
"10", NULL);
92 const double tol =
scan_ctl(argc, argv,
"INVERT_TOL", -1,
"1e-4", NULL);
93
94
97 if (strcmp(ctl.
emitter[0],
"SO2") != 0)
98 ERRMSG(
"Set EMITTER[0] = SO2!");
99 if (strcmp(ctl.
emitter[1],
"H2O") != 0)
100 ERRMSG(
"Set EMITTER[1] = H2O!");
101 if (strcmp(ctl.
emitter[2],
"O3") != 0)
102 ERRMSG(
"Set EMITTER[2] = O3!");
103 if (strcmp(ctl.
emitter[3],
"CO2") != 0)
104 ERRMSG(
"Set EMITTER[3] = CO2!");
107
108
111
112
115
116
118
119
120
121
122
123
124 LOG(1,
"Read profile data: %s", argv[2]);
125
126
127 if (!(in = fopen(argv[2], "r")))
128 ERRMSG(
"Cannot open file!");
129
130
131 while (fgets(line,
LEN, in)) {
132
133
134 if (sscanf(line, "%lg %lg %lg %lg %g %g %g %g %g %g", &rtime[nl],
135 &rz[nl], &rlon[nl], &rlat[nl], &rp[nl], &rt[nl], &rso2[nl],
136 &rh2o[nl], &ro3[nl], &robs[nl]) != 10)
137 continue;
139 ERRMSG(
"Too many profile data points!");
140 }
141
142
143 fclose(in);
144
145
146
147
148
149
150 for (int it = 0; it < itmax; it++) {
151
152
154 for (i = 0; i <
NMAX; i++) {
155 ndata[i] = 0;
156 x[i] = y[i] = NAN;
157 }
158
159
160 for (int il = 0; il < nl; il++) {
161
162
163 if ((rtime[il] != atm.
time[0]
164 || rlon[il] != atm.
lon[0]
165 || rlat[il] != atm.
lat[0])
167
168
169 formod(&ctl, tbl, &atm, &obs);
170 obs_sim = obs.
rad[0][0] - obs.
rad[1][0];
171
172
173 i = (int) ((atm.
time[0] - rtime[0]) / dt);
174 if (i < 0 && i >=
NMAX)
175 ERRMSG(
"Time index out of range!");
176
177
178 if (data == 1) {
179 x[i] = (isfinite(x[i]) ?
MAX(x[i], obs_sim) : obs_sim);
180 y[i] = (isfinite(y[i]) ?
MAX(y[i], obs_meas) : obs_meas);
181 y_err[i] = obs_err;
182 if (isfinite(x[i]) && isfinite(y[i]))
183 ndata[i] = 1;
184 }
185
186
187 else if (data == 2) {
188 if (ndata[i] == 0) {
189 x[i] = obs_sim;
190 y[i] = obs_meas;
191 y_err[i] =
POW2(obs_meas);
192 } else {
193 x[i] += obs_sim;
194 y[i] += obs_meas;
195 y_err[i] +=
POW2(obs_meas);
196 }
197 ndata[i]++;
198 }
199
200
201 nprof++;
203 for (
int ip = 0; ip < atm.
np; ip++) {
205 atm2.
z[ip] += atm.
z[ip];
206 atm2.
lon[ip] += atm.
lon[ip];
207 atm2.
lat[ip] += atm.
lat[ip];
208 atm2.
p[ip] += atm.
p[ip];
209 atm2.
t[ip] += atm.
t[ip];
210 for (ig = 1; ig < ctl.
ng; ig++)
211 atm2.
q[ig][ip] += atm.
q[ig][ip];
212 }
213
214
216 }
217
218
219 obs_meas = robs[il];
220 atm.
time[atm.
np] = rtime[il];
221 atm.
z[atm.
np] = rz[il];
222 atm.
lon[atm.
np] = rlon[il];
223 atm.
lat[atm.
np] = rlat[il];
224 atm.
p[atm.
np] = rp[il];
225 atm.
t[atm.
np] = rt[il];
226 atm.
q[0][atm.
np] = rso2[il] * scl;
227 atm.
q[1][atm.
np] = rh2o[il];
228 atm.
q[2][atm.
np] = ro3[il];
231 ERRMSG(
"Too many data points!");
232 }
233
234
235 if (data == 2)
236 for (i = 0; i <
NMAX; i++)
237 if (ndata[i] > 0) {
238 x[i] /= ndata[i];
239 y[i] /= ndata[i];
240 y_err[i] = sqrt(
MAX(y_err[i] / ndata[i] -
POW2(y[i]), 0.0))
241 / sqrt(ndata[i]);
242 }
243
244
245 n = 0;
246 for (i = 0; i <
NMAX; i++)
247 if (ndata[i] > 0 && isfinite(x[i]) && isfinite(y[i])
248 && isfinite(y_err[i])) {
249 x2[n] = x[i];
250 y2[n] = y[i];
251 y2_err[n] = y_err[i];
252 w2[n] = 1. /
POW2(y_err[i]);
253 n++;
254 }
255
256
257 if (fit == 1)
258 gsl_fit_mul(x2, 1, y2, 1, (size_t) n, &c1, &cov11, &sumsq);
259 else if (fit == 2)
260 gsl_fit_wmul(x2, 1, w2, 1, y2, 1, (size_t) n, &c1, &cov11, &sumsq);
261 else if (fit == 3)
262 gsl_fit_linear(x2, 1, y2, 1, (size_t) n, &c0, &c1, &cov00, &cov01,
263 &cov11, &sumsq);
264 else if (fit == 4)
265 gsl_fit_wlinear(x2, 1, w2, 1, y2, 1, (size_t) n, &c0, &c1, &cov00,
266 &cov01, &cov11, &sumsq);
267 else
268 ERRMSG(
"Check INVERT_FIT!");
269
270
271 double scl_old = scl;
272 scl_err = scl * sqrt(cov11);
273 scl *= c1;
274
275
276 LOG(1,
" it= %d | scl= %g +/- %g | RMSE= %g",
277 it, scl, scl_err, sqrt(sumsq / n));
278
279
280 if (fabs(2.0 * (scl - scl_old) / (scl + scl_old)) < tol)
281 break;
282 }
283
284
285
286
287
288
289 LOG(1,
"Write inversion data: %s", argv[3]);
290
291
292 if (!(out = fopen(argv[3], "w")))
293 ERRMSG(
"Cannot create file!");
294
295
296 fprintf(out,
297 "# $1 = time (seconds since 2000-01-01T00:00Z)\n"
298 "# $2 = simulated SO2 index [K]\n"
299 "# $3 = scaled simulated SO2 index [K]\n"
300 "# $4 = error of scaled simulated SO2 index [K]\n"
301 "# $5 = observed SO2 index [K]\n"
302 "# $6 = error of observed SO2 index [K]\n\n");
303
304
305 for (i = 0; i < n; i++) {
306
307
308 if (fit == 1 || fit == 2)
309 gsl_fit_mul_est(x2[i], c1, cov11, &y2_sim[i], &y2_sim_err[i]);
310 else if (fit == 3 || fit == 4)
311 gsl_fit_linear_est(x2[i], c0, c1, cov00, cov01, cov11, &y2_sim[i],
312 &y2_sim_err[i]);
313
314
315 fprintf(out, "%.2f %g %g %g %g %g\n", rtime[0] + (i + 0.5) * dt,
316 x2[i], y2_sim[i], y2_sim_err[i], y2[i], y2_err[i]);
317 }
318
319
320 fprintf(out, "\n");
321 fprintf(out, "# scl= %g +/- %g\n", scl, scl_err);
322 fprintf(out, "# c1= %g +/- %g\n", c1, sqrt(cov11));
323 if (fit == 3 || fit == 4) {
324 fprintf(out, "# c0= %g +/- %g\n", c0, sqrt(cov00));
325 fprintf(out, "# corr= %g\n", cov01 / (sqrt(cov00) * sqrt(cov11)));
326 }
327 fprintf(out, "# RMSE= %g\n", sqrt(sumsq / n));
328 fprintf(out, "# n= %d\n", n);
329
330
331 fclose(out);
332
333
334
335
336
337
338 for (
int ip = 0; ip < atm2.
np; ip++) {
339 atm2.
time[ip] /= nprof;
341 atm2.
lon[ip] /= nprof;
342 atm2.
lat[ip] /= nprof;
345 for (ig = 0; ig < ctl.
ng; ig++)
346 atm2.
q[ig][ip] /= nprof;
347 }
348
349
350 const size_t nk =
atm2x(&ctl, &atm2, NULL, NULL, NULL);
351 const size_t mk =
obs2y(&ctl, &obs, NULL, NULL, NULL);
352
353
354 gsl_matrix *k = gsl_matrix_alloc(mk, nk);
355
356
357 kernel(&ctl, tbl, &atm2, &obs, k);
358
359
361
362
363 write_matrix(NULL, argv[5], &ctl, k, &atm2, &obs,
"y",
"x",
"r", 0);
364
365
366 gsl_matrix_free(k);
368
369 return EXIT_SUCCESS;
370}
#define NMAX
Maximum number of data points...
#define NLMAX
Maximum number of data lines...
void write_atm(const char *dirname, const char *filename, const ctl_t *ctl, const atm_t *atm, int profile)
Write atmospheric data to a file.
void tbl_free(const ctl_t *ctl, tbl_t *tbl)
Free lookup table and all internally allocated memory.
void read_ctl(int argc, char *argv[], ctl_t *ctl)
Read model control parameters from command-line and configuration input.
void write_matrix(const char *dirname, const char *filename, const ctl_t *ctl, const gsl_matrix *matrix, const atm_t *atm, const obs_t *obs, const char *rowspace, const char *colspace, const char *sort, int dataset)
Write a fully annotated matrix (e.g., Jacobian or gain matrix) to file.
void kernel(const ctl_t *ctl, const tbl_t *tbl, atm_t *atm, obs_t *obs, gsl_matrix *k)
Compute the Jacobian (kernel) matrix by finite differences.
double scan_ctl(int argc, char *argv[], const char *varname, const int arridx, const char *defvalue, char *value)
Scan control file or command-line arguments for a configuration variable.
void formod(const ctl_t *ctl, const tbl_t *tbl, atm_t *atm, obs_t *obs)
Execute the selected forward model.
tbl_t * read_tbl(const ctl_t *ctl)
Read emissivity lookup tables from disk.
size_t obs2y(const ctl_t *ctl, const obs_t *obs, gsl_vector *y, int *ida, int *ira)
Convert observation radiances into a measurement vector.
size_t atm2x(const ctl_t *ctl, const atm_t *atm, gsl_vector *x, int *iqa, int *ipa)
Convert atmospheric data to state vector elements.
#define LEN
Maximum length of ASCII data lines.
#define POW2(x)
Compute the square of a value.
#define ERRMSG(...)
Print an error message with contextual information and terminate the program.
#define USAGE
Print usage information on -h or --help.
#define CO2_FIT(TIME)
Climatological CO2 volume mixing ratio as a quadratic function of time.
#define NP
Maximum number of atmospheric data points.
#define LOG(level,...)
Print a log message with a specified logging level.
#define MAX(a, b)
Determine the maximum of two values.
Atmospheric profile data.
double time[NP]
Time (seconds since 2000-01-01T00:00Z).
double lat[NP]
Latitude [deg].
double lon[NP]
Longitude [deg].
double t[NP]
Temperature [K].
int np
Number of data points.
double z[NP]
Altitude [km].
double q[NG][NP]
Volume mixing ratio [ppv].
double p[NP]
Pressure [hPa].
int write_matrix
Write matrix file (0=no, 1=yes).
int ng
Number of emitters.
int nd
Number of radiance channels.
char emitter[NG][LEN]
Name of each emitter.
int write_bbt
Use brightness temperature instead of radiance (0=no, 1=yes).
Observation geometry and radiance data.
double rad[ND][NR]
Radiance [W/(m^2 sr cm^-1)].
double obsz[NR]
Observer altitude [km].
int nr
Number of ray paths.
Emissivity look-up tables.