GPS Code Collection
libgps.h
Go to the documentation of this file.
1/*
2 This file is part of the GPS Code Collection.
3
4 the GPS Code Collections is free software: you can redistribute it
5 and/or modify it under the terms of the GNU General Public License
6 as published by the Free Software Foundation, either version 3 of
7 the License, or (at your option) any later version.
8
9 The GPS Code Collection is distributed in the hope that it will be
10 useful, but WITHOUT ANY WARRANTY; without even the implied warranty
11 of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 General Public License for more details.
13
14 You should have received a copy of the GNU General Public License
15 along with the GPS Code Collection. If not, see
16 <http://www.gnu.org/licenses/>.
17
18 Copyright (C) 2019-2025 Forschungszentrum Juelich GmbH
19*/
20
80#include <netcdf.h>
81#include <gsl/gsl_multifit.h>
82#include <gsl/gsl_poly.h>
83#include <gsl/gsl_spline.h>
84#include "jurassic.h"
85
86/* ------------------------------------------------------------
87 Dimensions...
88 ------------------------------------------------------------ */
89
91#define EP 73
92
94#define EX 721
95
97#define EY 361
98
100#define NDS 10000
101
103#define NZ 5000
104
105/* ------------------------------------------------------------
106 Constants...
107 ------------------------------------------------------------ */
108
110#ifndef MA
111#define MA 28.9644
112#endif
113
115#ifndef RA
116#define RA (1e3 * RI / MA)
117#endif
118
120#ifndef RI
121#define RI 8.3144598
122#endif
123
124/* ------------------------------------------------------------
125 Macros...
126 ------------------------------------------------------------ */
127
129#define LAPSE(p1, t1, p2, t2) \
130 (1e3 * G0 / RA * ((t2) - (t1)) / ((t2) + (t1)) \
131 * ((p2) + (p1)) / ((p2) - (p1)))
132
134#define NC(cmd) { \
135 int nc_result=(cmd); \
136 if(nc_result!=NC_NOERR) \
137 ERRMSG("%s", nc_strerror(nc_result)); \
138 }
139
141#define P(z) \
142 (P0 * exp(-(z) / H0))
143
145#define Z(p) \
146 (H0 * log(P0 / (p)))
147
148/* ------------------------------------------------------------
149 Structs...
150 ------------------------------------------------------------ */
151
153typedef struct {
154
156 int nds;
157
159 int nz[NDS];
160
162 double time[NDS];
163
165 double z[NDS][NZ];
166
168 double lon[NDS][NZ];
169
171 double lat[NDS][NZ];
172
174 double p[NDS][NZ];
175
177 double t[NDS][NZ];
178
180 double wv[NDS][NZ];
181
183 double pt[NDS][NZ];
184
186 double th[NDS];
187
189 double tp[NDS];
190
192 double tt[NDS];
193
195 double tq[NDS];
196
198 double tlon[NDS];
199
201 double tlat[NDS];
202
203} gps_t;
204
206typedef struct {
207
209 double time;
210
212 int nx;
213
215 int ny;
216
218 int np;
219
221 double lon[EX];
222
224 double lat[EY];
225
227 double p[EP];
228
230 float t[EX][EY][EP];
231
232} met_t;
233
234/* ------------------------------------------------------------
235 Functions...
236 ------------------------------------------------------------ */
237
239void add_var(
240 int ncid,
241 const char *varname,
242 const char *unit,
243 const char *longname,
244 int type,
245 int dimid[],
246 int *varid,
247 int ndims);
248
250void detrend_met(
251 gps_t * gps,
252 char *metbase,
253 double dt_met);
254
256void gauss(
257 gps_t * gps,
258 double dx,
259 double dy);
260
262void grid_gps(
263 gps_t * gps,
264 double zmin,
265 double zmax,
266 int nz);
267
269void get_met(
270 char *metbase,
271 double dt_met,
272 double t,
273 met_t * met0,
274 met_t * met1);
275
277void get_met_help(
278 double t,
279 int direct,
280 char *metbase,
281 double dt_met,
282 char *filename);
283
285void intpol_met_3d(
286 float array[EX][EY][EP],
287 int ip,
288 int ix,
289 int iy,
290 double wp,
291 double wx,
292 double wy,
293 double *var);
294
297 met_t * met,
298 double p,
299 double lon,
300 double lat,
301 double *t);
302
304void intpol_met_time(
305 met_t * met0,
306 met_t * met1,
307 double ts,
308 double p,
309 double lon,
310 double lat,
311 double *t);
312
315 gps_t * gps,
316 double dz);
317
320 gps_t * gps,
321 double dz);
322
324void poly(
325 gps_t * gps,
326 int dim,
327 double zmin,
328 double zmax);
329
331void poly_help(
332 double *xx,
333 double *yy,
334 int n,
335 int dim,
336 double xmin,
337 double xmax);
338
340void read_gps_prof(
341 char *filename,
342 gps_t * gps);
343
345void read_gps(
346 char *filename,
347 gps_t * gps);
348
350void read_met(
351 char *filename,
352 met_t * met);
353
356 met_t * met);
357
359void read_met_help(
360 int ncid,
361 char *varname,
362 char *varname2,
363 met_t * met,
364 float dest[EX][EY][EP],
365 float scl);
366
369 met_t * met);
370
372void spline(
373 const double *x,
374 const double *y,
375 const int n,
376 const double *x2,
377 double *y2,
378 const int n2,
379 const int method);
380
382void tropopause(
383 gps_t * gps);
384
387 gps_t * gps,
388 int met_tropo);
389
391void write_gps(
392 char *filename,
393 gps_t * gps);
void read_met_extrapolate(met_t *met)
Extrapolate meteorological data at lower boundary.
Definition: libgps.c:762
void intpol_met_space(met_t *met, double p, double lon, double lat, double *t)
Spatial interpolation of meteorological data.
Definition: libgps.c:274
void write_gps(char *filename, gps_t *gps)
Write GPS-RO data file.
Definition: libgps.c:1060
void read_met_periodic(met_t *met)
Create meteorological data with periodic boundary conditions.
Definition: libgps.c:816
#define EY
Maximum number of latitudes for meteorological data.
Definition: libgps.h:97
void gauss(gps_t *gps, double dx, double dy)
Calculate horizontal Gaussian mean to extract perturbations.
Definition: libgps.c:95
void read_met(char *filename, met_t *met)
Read meteorological data file.
Definition: libgps.c:682
void intpol_met_time(met_t *met0, met_t *met1, double ts, double p, double lon, double lat, double *t)
Temporal interpolation of meteorological data.
Definition: libgps.c:301
void spline(const double *x, const double *y, const int n, const double *x2, double *y2, const int n2, const int method)
Performs spline interpolation or linear interpolation.
Definition: libgps.c:839
void poly(gps_t *gps, int dim, double zmin, double zmax)
Remove polynomial fit from perturbation profile.
Definition: libgps.c:435
#define NZ
Maximum number of altitudes per GPS-RO profile.
Definition: libgps.h:103
void grid_gps(gps_t *gps, double zmin, double zmax, int nz)
Interpolate GPS data to regular altitude grid.
Definition: libgps.c:133
#define EX
Maximum number of longitudes for meteorological data.
Definition: libgps.h:94
void intpol_met_3d(float array[EX][EY][EP], int ip, int ix, int iy, double wp, double wx, double wy, double *var)
Linear interpolation of 3-D meteorological data.
Definition: libgps.c:246
#define NDS
Maximum number of GPS-RO profiles.
Definition: libgps.h:100
void add_var(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: libgps.c:30
void get_met(char *metbase, double dt_met, double t, met_t *met0, met_t *met1)
Get meteorological data for given timestep.
Definition: libgps.c:188
void read_gps(char *filename, gps_t *gps)
Read GPS-RO data file.
Definition: libgps.c:611
void hamming_low_pass(gps_t *gps, double dz)
Apply vertical Hamming filter to extract perturbations.
Definition: libgps.c:325
void detrend_met(gps_t *gps, char *metbase, double dt_met)
Detrending by means of meteo data.
Definition: libgps.c:57
void read_gps_prof(char *filename, gps_t *gps)
Read GPS-RO profile.
Definition: libgps.c:520
void get_met_help(double t, int direct, char *metbase, double dt_met, char *filename)
Get meteorological data for timestep.
Definition: libgps.c:220
void hamming_high_pass(gps_t *gps, double dz)
Apply vertical Hamming filter to reduce noise.
Definition: libgps.c:382
void tropopause_spline(gps_t *gps, int met_tropo)
Find tropopause height using cubic spline interpolation.
Definition: libgps.c:922
void poly_help(double *xx, double *yy, int n, int dim, double xmin, double xmax)
Auxiliary function for polynomial interpolation.
Definition: libgps.c:461
void read_met_help(int ncid, char *varname, char *varname2, met_t *met, float dest[EX][EY][EP], float scl)
Read and convert variable from meteorological data file.
Definition: libgps.c:784
void tropopause(gps_t *gps)
Find tropopause height.
Definition: libgps.c:886
#define EP
Maximum number of pressure levels for meteorological data.
Definition: libgps.h:91
GPS-RO profile data.
Definition: libgps.h:153
int nds
Number of profiles.
Definition: libgps.h:156
Meteorological data.
Definition: libgps.h:206
int nx
Number of longitudes.
Definition: libgps.h:212
int ny
Number of latitudes.
Definition: libgps.h:215
int np
Number of pressure levels.
Definition: libgps.h:218
double time
Time [s].
Definition: libgps.h:209