MPTRAC
mptrac.h
Go to the documentation of this file.
1/*
2 This file is part of MPTRAC.
3
4 MPTRAC is free software: you can redistribute it and/or modify it
5 under the terms of the GNU General Public License as published by
6 the Free Software Foundation, either version 3 of the License, or
7 (at your option) any later version.
8
9 MPTRAC is distributed in the hope that it will be useful,
10 but WITHOUT ANY WARRANTY; without even the implied warranty of
11 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12 GNU General Public License for more details.
13
14 You should have received a copy of the GNU General Public License
15 along with MPTRAC. If not, see <http://www.gnu.org/licenses/>.
16
17 Copyright (C) 2013-2026 Forschungszentrum Juelich GmbH
18*/
19
116#ifndef LIBTRAC_H
117#define LIBTRAC_H
118
119/* ------------------------------------------------------------
120 Includes...
121 ------------------------------------------------------------ */
122
123#include <ctype.h>
124#include <gsl/gsl_fft_complex.h>
125#include <gsl/gsl_math.h>
126#include <gsl/gsl_randist.h>
127#include <gsl/gsl_rng.h>
128#include <gsl/gsl_sort.h>
129#include <gsl/gsl_spline.h>
130#include <gsl/gsl_statistics.h>
131#include <math.h>
132#include <netcdf.h>
133#include <omp.h>
134#include <stdint.h>
135#include <stdio.h>
136#include <stdlib.h>
137#include <string.h>
138#include <time.h>
139#include <sys/time.h>
140
141#ifdef MPI
142#include "mpi.h"
143#else
145#define MPI_Datatype void*
146#endif
147
148#ifdef DD
149#include <netcdf_par.h>
150#endif
151
152#ifdef _OPENACC
153#include "openacc.h"
154#endif
155
156#ifdef CURAND
157#include "curand.h"
158#endif
159
160#ifdef THRUST
161#include "thrustsort.h"
162#endif
163
164#ifdef ZFP
165#include "zfp.h"
166#endif
167
168#ifdef ZSTD
169#include "zstd.h"
170#endif
171
172#ifdef SZ3
173#include "SZ3c/sz3c.h"
174#endif
175
176#ifdef CMS
177#include "cmultiscale.h"
178#endif
179
180#ifdef KPP
181#include "chem_Parameters.h"
182#include "chem_Global.h"
183#include "chem_Sparse.h"
184#endif
185
186#ifdef ECCODES
187#include "eccodes.h"
188#else
190#define codes_handle void*
191#endif
192
193/* ------------------------------------------------------------
194 Constants...
195 ------------------------------------------------------------ */
196
198#ifndef AVO
199#define AVO 6.02214076e23
200#endif
201
203#ifndef CPD
204#define CPD 1003.5
205#endif
206
208#ifndef EPS
209#define EPS (MH2O / MA)
210#endif
211
213#ifndef G0
214#define G0 9.80665
215#endif
216
218#ifndef H0
219#define H0 7.0
220#endif
221
223#ifndef LV
224#define LV 2501000.
225#endif
226
228#ifndef KARMAN
229#define KARMAN 0.40
230#endif
231
233#ifndef KB
234#define KB 1.3806504e-23
235#endif
236
238#ifndef MA
239#define MA 28.9644
240#endif
241
243#ifndef MH2O
244#define MH2O 18.01528
245#endif
246
248#ifndef MO3
249#define MO3 48.00
250#endif
251
253#ifndef P0
254#define P0 1013.25
255#endif
256
258#ifndef RA
259#define RA (1e3 * RI / MA)
260#endif
261
263#ifndef RE
264#define RE 6367.421
265#endif
266
268#ifndef RI
269#define RI 8.3144598
270#endif
271
273#ifndef T0
274#define T0 273.15
275#endif
276
277/* ------------------------------------------------------------
278 Dimensions...
279 ------------------------------------------------------------ */
280
282#ifndef EP
283#define EP 140
284#endif
285
287#ifndef EX
288#define EX 1444
289#endif
290
292#ifndef EY
293#define EY 724
294#endif
295
297#ifndef LEN
298#define LEN 5000
299#endif
300
302#ifndef METVAR
303#define METVAR 13
304#endif
305
307#ifndef NP
308#define NP 10000000
309#endif
310
312#ifndef NQ
313#define NQ 15
314#endif
315
317#ifndef NCSI
318#define NCSI 1000000
319#endif
320
322#ifndef NENS
323#define NENS 2000
324#endif
325
327#ifndef NOBS
328#define NOBS 10000000
329#endif
330
332#ifndef NTHREADS
333#define NTHREADS 512
334#endif
335
337#ifndef CY
338#define CY 250
339#endif
340
342#ifndef CO3
343#define CO3 30
344#endif
345
347#ifndef CP
348#define CP 70
349#endif
350
352#ifndef CSZA
353#define CSZA 50
354#endif
355
357#ifndef CT
358#define CT 12
359#endif
360
362#ifndef CTS
363#define CTS 1000
364#endif
365
367#ifndef DD_NPART
368#define DD_NPART 100000
369#endif
370
372#ifndef DD_NNMAX
373#define DD_NNMAX 26
374#endif
375
377#ifndef DD_NPOLE
378#define DD_NPOLE -2
379#endif
380
382#ifndef DD_SPOLE
383#define DD_SPOLE -3
384#endif
385
386/* ------------------------------------------------------------
387 Macros...
388 ------------------------------------------------------------ */
389
409#ifdef _OPENACC
410#define ALLOC(ptr, type, n) \
411 if(acc_get_num_devices(acc_device_nvidia) <= 0) \
412 ERRMSG("Not running on a GPU device!"); \
413 if((ptr=calloc((size_t)(n), sizeof(type)))==NULL) \
414 ERRMSG("Out of memory!");
415#else
416#define ALLOC(ptr, type, n) \
417 if((ptr=calloc((size_t)(n), sizeof(type)))==NULL) \
418 ERRMSG("Out of memory!");
419#endif
420
439#define ARRAY_2D(ix, iy, ny) \
440 ((ix) * (ny) + (iy))
441
458#define ARRAY_3D(ix, iy, ny, iz, nz) \
459 (((ix)*(ny) + (iy)) * (nz) + (iz))
460
483#define ARRHENIUS(a, b, t) \
484 ((a) * exp( -(b) / (t)))
485
505#define CLAMP(v, lo, hi) \
506 (((v) < (lo)) ? (lo) : (((v) > (hi)) ? (hi) : (v)))
507
529#define DEG2DX(dlon, lat) \
530 (RE * DEG2RAD(dlon) * cos(DEG2RAD(lat)))
531
550#define DEG2DY(dlat) \
551 (RE * DEG2RAD(dlat))
552
567#define DEG2RAD(deg) \
568 ((deg) * (M_PI / 180.0))
569
592#define DP2DZ(dp, p) \
593 (- (dp) * H0 / (p))
594
614#define DX2DEG(dx, lat) \
615 (((lat) < -89.999 || (lat) > 89.999) ? 0 \
616 : (dx) * 180. / (M_PI * RE * cos(DEG2RAD(lat))))
617
632#define DY2DEG(dy) \
633 ((dy) * 180. / (M_PI * RE))
634
651#define DZ2DP(dz, p) \
652 (-(dz) * (p) / H0)
653
667#define DIST(a, b) \
668 sqrt(DIST2(a, b))
669
683#define DIST2(a, b) \
684 ((a[0]-b[0])*(a[0]-b[0])+(a[1]-b[1])*(a[1]-b[1])+(a[2]-b[2])*(a[2]-b[2]))
685
699#define DOTP(a, b) \
700 (a[0]*b[0]+a[1]*b[1]+a[2]*b[2])
701
713#define ECC(cmd) { \
714 int ecc_result=(cmd); \
715 if(ecc_result!=0) \
716 ERRMSG("ECCODES error: %s", codes_get_error_message(ecc_result)); \
717 }
718
732#define ECC_READ_2D(variable, target, scaling_factor, found_flag) { \
733 if(strcmp(short_name, variable) == 0) { \
734 for (int ix = 0; ix < met->nx; ix++) \
735 for (int iy = 0; iy < met->ny; iy++) \
736 target[ix][iy] = (float)(values[iy * met->nx + ix] * scaling_factor); \
737 found_flag = 1; \
738 } \
739 }
740
755#define ECC_READ_3D(variable, level, target, scaling_factor, found_flag) { \
756 if(strcmp(short_name, variable) == 0) { \
757 for (int ix = 0; ix < met->nx; ix++) \
758 for (int iy = 0; iy < met->ny; iy++) \
759 target[ix][iy][level] = (float) (values[iy * met->nx + ix] * scaling_factor); \
760 found_flag += 1; \
761 } \
762 }
763
780#define FMOD(x, y) \
781 ((x) - (int) ((x) / (y)) * (y))
782
798#define FREAD(ptr, type, size, in) { \
799 if(fread(ptr, sizeof(type), size, in)!=size) \
800 ERRMSG("Error while reading!"); \
801 }
802
818#define FWRITE(ptr, type, size, out) { \
819 if(fwrite(ptr, sizeof(type), size, out)!=size) \
820 ERRMSG("Error while writing!"); \
821 }
822
833#define INTPOL_INIT \
834 double cw[4] = {0.0, 0.0, 0.0, 0.0}; int ci[3] = {0, 0, 0};
835
847#define INTPOL_2D(var, init) \
848 intpol_met_time_2d(met0, met0->var, met1, met1->var, \
849 atm->time[ip], atm->lon[ip], atm->lat[ip], \
850 &var, ci, cw, init);
851
864#define INTPOL_3D(var, init) \
865 intpol_met_time_3d(met0, met0->var, met1, met1->var, \
866 atm->time[ip], atm->p[ip], \
867 atm->lon[ip], atm->lat[ip], \
868 &var, ci, cw, init);
869
883#define INTPOL_SPACE_ALL(p, lon, lat) { \
884 intpol_met_space_3d(met, met->z, p, lon, lat, &z, ci, cw, 1); \
885 intpol_met_space_3d(met, met->t, p, lon, lat, &t, ci, cw, 0); \
886 intpol_met_space_3d(met, met->u, p, lon, lat, &u, ci, cw, 0); \
887 intpol_met_space_3d(met, met->v, p, lon, lat, &v, ci, cw, 0); \
888 intpol_met_space_3d(met, met->w, p, lon, lat, &w, ci, cw, 0); \
889 intpol_met_space_3d(met, met->pv, p, lon, lat, &pv, ci, cw, 0); \
890 intpol_met_space_3d(met, met->h2o, p, lon, lat, &h2o, ci, cw, 0); \
891 intpol_met_space_3d(met, met->o3, p, lon, lat, &o3, ci, cw, 0); \
892 intpol_met_space_3d(met, met->lwc, p, lon, lat, &lwc, ci, cw, 0); \
893 intpol_met_space_3d(met, met->rwc, p, lon, lat, &rwc, ci, cw, 0); \
894 intpol_met_space_3d(met, met->iwc, p, lon, lat, &iwc, ci, cw, 0); \
895 intpol_met_space_3d(met, met->swc, p, lon, lat, &swc, ci, cw, 0); \
896 intpol_met_space_3d(met, met->cc, p, lon, lat, &cc, ci, cw, 0); \
897 intpol_met_space_2d(met, met->ps, lon, lat, &ps, ci, cw, 0); \
898 intpol_met_space_2d(met, met->ts, lon, lat, &ts, ci, cw, 0); \
899 intpol_met_space_2d(met, met->zs, lon, lat, &zs, ci, cw, 0); \
900 intpol_met_space_2d(met, met->us, lon, lat, &us, ci, cw, 0); \
901 intpol_met_space_2d(met, met->vs, lon, lat, &vs, ci, cw, 0); \
902 intpol_met_space_2d(met, met->ess, ess, lat, &ess, ci, cw, 0); \
903 intpol_met_space_2d(met, met->nss, nss, lat, &nss, ci, cw, 0); \
904 intpol_met_space_2d(met, met->shf, shf, lat, &shf, ci, cw, 0); \
905 intpol_met_space_2d(met, met->lsm, lon, lat, &lsm, ci, cw, 0); \
906 intpol_met_space_2d(met, met->sst, lon, lat, &sst, ci, cw, 0); \
907 intpol_met_space_2d(met, met->pbl, lon, lat, &pbl, ci, cw, 0); \
908 intpol_met_space_2d(met, met->pt, lon, lat, &pt, ci, cw, 0); \
909 intpol_met_space_2d(met, met->tt, lon, lat, &tt, ci, cw, 0); \
910 intpol_met_space_2d(met, met->zt, lon, lat, &zt, ci, cw, 0); \
911 intpol_met_space_2d(met, met->h2ot, lon, lat, &h2ot, ci, cw, 0); \
912 intpol_met_space_2d(met, met->pct, lon, lat, &pct, ci, cw, 0); \
913 intpol_met_space_2d(met, met->pcb, lon, lat, &pcb, ci, cw, 0); \
914 intpol_met_space_2d(met, met->cl, lon, lat, &cl, ci, cw, 0); \
915 intpol_met_space_2d(met, met->plcl, lon, lat, &plcl, ci, cw, 0); \
916 intpol_met_space_2d(met, met->plfc, lon, lat, &plfc, ci, cw, 0); \
917 intpol_met_space_2d(met, met->pel, lon, lat, &pel, ci, cw, 0); \
918 intpol_met_space_2d(met, met->cape, lon, lat, &cape, ci, cw, 0); \
919 intpol_met_space_2d(met, met->cin, lon, lat, &cin, ci, cw, 0); \
920 intpol_met_space_2d(met, met->o3c, lon, lat, &o3c, ci, cw, 0); \
921 }
922
937#define INTPOL_TIME_ALL(time, p, lon, lat) { \
938 intpol_met_time_3d(met0, met0->z, met1, met1->z, time, p, lon, lat, &z, ci, cw, 1); \
939 intpol_met_time_3d(met0, met0->t, met1, met1->t, time, p, lon, lat, &t, ci, cw, 0); \
940 intpol_met_time_3d(met0, met0->u, met1, met1->u, time, p, lon, lat, &u, ci, cw, 0); \
941 intpol_met_time_3d(met0, met0->v, met1, met1->v, time, p, lon, lat, &v, ci, cw, 0); \
942 intpol_met_time_3d(met0, met0->w, met1, met1->w, time, p, lon, lat, &w, ci, cw, 0); \
943 intpol_met_time_3d(met0, met0->pv, met1, met1->pv, time, p, lon, lat, &pv, ci, cw, 0); \
944 intpol_met_time_3d(met0, met0->h2o, met1, met1->h2o, time, p, lon, lat, &h2o, ci, cw, 0); \
945 intpol_met_time_3d(met0, met0->o3, met1, met1->o3, time, p, lon, lat, &o3, ci, cw, 0); \
946 intpol_met_time_3d(met0, met0->lwc, met1, met1->lwc, time, p, lon, lat, &lwc, ci, cw, 0); \
947 intpol_met_time_3d(met0, met0->rwc, met1, met1->rwc, time, p, lon, lat, &rwc, ci, cw, 0); \
948 intpol_met_time_3d(met0, met0->iwc, met1, met1->iwc, time, p, lon, lat, &iwc, ci, cw, 0); \
949 intpol_met_time_3d(met0, met0->swc, met1, met1->swc, time, p, lon, lat, &swc, ci, cw, 0); \
950 intpol_met_time_3d(met0, met0->cc, met1, met1->cc, time, p, lon, lat, &cc, ci, cw, 0); \
951 intpol_met_time_2d(met0, met0->ps, met1, met1->ps, time, lon, lat, &ps, ci, cw, 0); \
952 intpol_met_time_2d(met0, met0->ts, met1, met1->ts, time, lon, lat, &ts, ci, cw, 0); \
953 intpol_met_time_2d(met0, met0->zs, met1, met1->zs, time, lon, lat, &zs, ci, cw, 0); \
954 intpol_met_time_2d(met0, met0->us, met1, met1->us, time, lon, lat, &us, ci, cw, 0); \
955 intpol_met_time_2d(met0, met0->vs, met1, met1->vs, time, lon, lat, &vs, ci, cw, 0); \
956 intpol_met_time_2d(met0, met0->ess, met1, met1->ess, time, lon, lat, &ess, ci, cw, 0); \
957 intpol_met_time_2d(met0, met0->nss, met1, met1->nss, time, lon, lat, &nss, ci, cw, 0); \
958 intpol_met_time_2d(met0, met0->shf, met1, met1->shf, time, lon, lat, &shf, ci, cw, 0); \
959 intpol_met_time_2d(met0, met0->lsm, met1, met1->lsm, time, lon, lat, &lsm, ci, cw, 0); \
960 intpol_met_time_2d(met0, met0->sst, met1, met1->sst, time, lon, lat, &sst, ci, cw, 0); \
961 intpol_met_time_2d(met0, met0->pbl, met1, met1->pbl, time, lon, lat, &pbl, ci, cw, 0); \
962 intpol_met_time_2d(met0, met0->pt, met1, met1->pt, time, lon, lat, &pt, ci, cw, 0); \
963 intpol_met_time_2d(met0, met0->tt, met1, met1->tt, time, lon, lat, &tt, ci, cw, 0); \
964 intpol_met_time_2d(met0, met0->zt, met1, met1->zt, time, lon, lat, &zt, ci, cw, 0); \
965 intpol_met_time_2d(met0, met0->h2ot, met1, met1->h2ot, time, lon, lat, &h2ot, ci, cw, 0); \
966 intpol_met_time_2d(met0, met0->pct, met1, met1->pct, time, lon, lat, &pct, ci, cw, 0); \
967 intpol_met_time_2d(met0, met0->pcb, met1, met1->pcb, time, lon, lat, &pcb, ci, cw, 0); \
968 intpol_met_time_2d(met0, met0->cl, met1, met1->cl, time, lon, lat, &cl, ci, cw, 0); \
969 intpol_met_time_2d(met0, met0->plcl, met1, met1->plcl, time, lon, lat, &plcl, ci, cw, 0); \
970 intpol_met_time_2d(met0, met0->plfc, met1, met1->plfc, time, lon, lat, &plfc, ci, cw, 0); \
971 intpol_met_time_2d(met0, met0->pel, met1, met1->pel, time, lon, lat, &pel, ci, cw, 0); \
972 intpol_met_time_2d(met0, met0->cape, met1, met1->cape, time, lon, lat, &cape, ci, cw, 0); \
973 intpol_met_time_2d(met0, met0->cin, met1, met1->cin, time, lon, lat, &cin, ci, cw, 0); \
974 intpol_met_time_2d(met0, met0->o3c, met1, met1->o3c, time, lon, lat, &o3c, ci, cw, 0); \
975 }
976
991#define LAPSE(p1, t1, p2, t2) \
992 (1e3 * G0 / RA * ((t2) - (t1)) / ((t2) + (t1)) \
993 * ((p2) + (p1)) / ((p2) - (p1)))
994
1010#define LIN(x0, y0, x1, y1, x) \
1011 ((y0)+((y1)-(y0))/((x1)-(x0))*((x)-(x0)))
1012
1037#define MAX(a,b) \
1038 (((a)>(b))?(a):(b))
1039
1051#define MET_HEADER \
1052 fprintf(out, \
1053 "# $1 = time [s]\n" \
1054 "# $2 = altitude [km]\n" \
1055 "# $3 = longitude [deg]\n" \
1056 "# $4 = latitude [deg]\n" \
1057 "# $5 = pressure [hPa]\n" \
1058 "# $6 = temperature [K]\n" \
1059 "# $7 = zonal wind [m/s]\n" \
1060 "# $8 = meridional wind [m/s]\n" \
1061 "# $9 = vertical velocity [hPa/s]\n" \
1062 "# $10 = H2O volume mixing ratio [ppv]\n"); \
1063 fprintf(out, \
1064 "# $11 = O3 volume mixing ratio [ppv]\n" \
1065 "# $12 = geopotential height [km]\n" \
1066 "# $13 = potential vorticity [PVU]\n" \
1067 "# $14 = surface pressure [hPa]\n" \
1068 "# $15 = surface temperature [K]\n" \
1069 "# $16 = surface geopotential height [km]\n" \
1070 "# $17 = surface zonal wind [m/s]\n" \
1071 "# $18 = surface meridional wind [m/s]\n" \
1072 "# $19 = eastward turbulent surface stress [N/m^2]\n" \
1073 "# $20 = northward turbulent surface stress [N/m^2]\n"); \
1074 fprintf(out, \
1075 "# $21 = surface sensible heat flux [W/m^2]\n" \
1076 "# $22 = land-sea mask [1]\n" \
1077 "# $23 = sea surface temperature [K]\n" \
1078 "# $24 = tropopause pressure [hPa]\n" \
1079 "# $25 = tropopause geopotential height [km]\n" \
1080 "# $26 = tropopause temperature [K]\n" \
1081 "# $27 = tropopause water vapor [ppv]\n" \
1082 "# $28 = cloud liquid water content [kg/kg]\n" \
1083 "# $29 = cloud rain water content [kg/kg]\n" \
1084 "# $30 = cloud ice water content [kg/kg]\n"); \
1085 fprintf(out, \
1086 "# $31 = cloud snow water content [kg/kg]\n" \
1087 "# $32 = cloud cover [1]\n" \
1088 "# $33 = total column cloud water [kg/m^2]\n" \
1089 "# $34 = cloud top pressure [hPa]\n" \
1090 "# $35 = cloud bottom pressure [hPa]\n" \
1091 "# $36 = pressure at lifted condensation level (LCL) [hPa]\n" \
1092 "# $37 = pressure at level of free convection (LFC) [hPa]\n" \
1093 "# $38 = pressure at equilibrium level (EL) [hPa]\n" \
1094 "# $39 = convective available potential energy (CAPE) [J/kg]\n" \
1095 "# $40 = convective inhibition (CIN) [J/kg]\n"); \
1096 fprintf(out, \
1097 "# $41 = relative humidity over water [%%]\n" \
1098 "# $42 = relative humidity over ice [%%]\n" \
1099 "# $43 = dew point temperature [K]\n" \
1100 "# $44 = frost point temperature [K]\n" \
1101 "# $45 = NAT temperature [K]\n" \
1102 "# $46 = HNO3 volume mixing ratio [ppv]\n" \
1103 "# $47 = OH volume mixing ratio [ppv]\n" \
1104 "# $48 = H2O2 volume mixing ratio [ppv]\n" \
1105 "# $49 = HO2 volume mixing ratio [ppv]\n" \
1106 "# $50 = O(1D) volume mixing ratio [ppv]\n"); \
1107 fprintf(out, \
1108 "# $51 = boundary layer pressure [hPa]\n" \
1109 "# $52 = total column ozone [DU]\n" \
1110 "# $53 = number of data points\n" \
1111 "# $54 = number of tropopause data points\n" \
1112 "# $55 = number of CAPE data points\n");
1113
1138#define MIN(a,b) \
1139 (((a)<(b))?(a):(b))
1140
1153#define MOLEC_DENS(p,t) \
1154 (AVO * 1e-6 * ((p) * 100) / (RI * (t)))
1155
1167#define NC(cmd) { \
1168 int nc_result=(cmd); \
1169 if(nc_result!=NC_NOERR) \
1170 ERRMSG("%s", nc_strerror(nc_result)); \
1171 }
1172
1196#define NC_DEF_VAR(varname, type, ndims, dims, long_name, units, level, quant) { \
1197 NC(nc_def_var(ncid, varname, type, ndims, dims, &varid)); \
1198 NC(nc_put_att_text(ncid, varid, "long_name", strnlen(long_name, LEN), long_name)); \
1199 NC(nc_put_att_text(ncid, varid, "units", strnlen(units, LEN), units)); \
1200 if((quant) > 0) \
1201 NC(nc_def_var_quantize(ncid, varid, NC_QUANTIZE_GRANULARBR, quant)); \
1202 if((level) != 0) { \
1203 NC(nc_def_var_deflate(ncid, varid, 1, 1, level)); \
1204 /* unsigned int ulevel = (unsigned int)level; */ \
1205 /* NC(nc_def_var_filter(ncid, varid, 32015, 1, (unsigned int[]){ulevel})); */ \
1206 } \
1207 }
1208
1226#define NC_GET_DOUBLE(varname, ptr, force) { \
1227 if(force) { \
1228 NC(nc_inq_varid(ncid, varname, &varid)); \
1229 NC(nc_get_var_double(ncid, varid, ptr)); \
1230 } else { \
1231 if(nc_inq_varid(ncid, varname, &varid) == NC_NOERR) { \
1232 NC(nc_get_var_double(ncid, varid, ptr)); \
1233 } else \
1234 WARN("netCDF variable %s is missing!", varname); \
1235 } \
1236 }
1237
1256#define NC_INQ_DIM(dimname, ptr, min, max, check) { \
1257 int dimid; size_t naux; \
1258 NC(nc_inq_dimid(ncid, dimname, &dimid)); \
1259 NC(nc_inq_dimlen(ncid, dimid, &naux)); \
1260 *ptr = (int)naux; \
1261 if (check) \
1262 if ((*ptr) < (min) || (*ptr) > (max)) \
1263 ERRMSG("Dimension %s is out of range!", dimname); \
1264 }
1265
1280#define NC_PUT_DOUBLE(varname, ptr, hyperslab) { \
1281 NC(nc_inq_varid(ncid, varname, &varid)); \
1282 if(hyperslab) { \
1283 NC(nc_put_vara_double(ncid, varid, start, count, ptr)); \
1284 } else { \
1285 NC(nc_put_var_double(ncid, varid, ptr)); \
1286 } \
1287 }
1288
1304#define NC_PUT_FLOAT(varname, ptr, hyperslab) { \
1305 NC(nc_inq_varid(ncid, varname, &varid)); \
1306 if(hyperslab) { \
1307 NC(nc_put_vara_float(ncid, varid, start, count, ptr)); \
1308 } else { \
1309 NC(nc_put_var_float(ncid, varid, ptr)); \
1310 } \
1311 }
1312
1327#define NC_PUT_INT(varname, ptr, hyperslab) { \
1328 NC(nc_inq_varid(ncid, varname, &varid)); \
1329 if(hyperslab) { \
1330 NC(nc_put_vara_int(ncid, varid, start, count, ptr)); \
1331 } else { \
1332 NC(nc_put_var_int(ncid, varid, ptr)); \
1333 } \
1334 }
1335
1349#define NC_PUT_ATT(varname, attname, text) { \
1350 NC(nc_inq_varid(ncid, varname, &varid)); \
1351 NC(nc_put_att_text(ncid, varid, attname, strnlen(text, LEN), text)); \
1352 }
1353
1366#define NC_PUT_ATT_GLOBAL(attname, text) \
1367 NC(nc_put_att_text(ncid, NC_GLOBAL, attname, strnlen(text, LEN), text));
1368
1386#define NN(x0, y0, x1, y1, x) \
1387 (fabs((x) - (x0)) <= fabs((x) - (x1)) ? (y0) : (y1))
1388
1404#ifdef _OPENACC
1405#define PARTICLE_LOOP(ip0, ip1, check_dt, ...) \
1406 const int ip0_const = ip0; \
1407 const int ip1_const = ip1; \
1408 _Pragma(__VA_ARGS__) \
1409 _Pragma("acc parallel loop independent gang vector") \
1410 for (int ip = ip0_const; ip < ip1_const; ip++) \
1411 if (!check_dt || cache->dt[ip] != 0)
1412#else
1413#define PARTICLE_LOOP(ip0, ip1, check_dt, ...) \
1414 const int ip0_const = ip0; \
1415 const int ip1_const = ip1; \
1416 _Pragma("omp parallel for default(shared)") \
1417 for (int ip = ip0_const; ip < ip1_const; ip++) \
1418 if (!check_dt || cache->dt[ip] != 0)
1419#endif
1420
1443#define P(z) \
1444 (P0 * exp(-(z) / H0))
1445
1467#define PSAT(t) \
1468 (6.112 * exp(17.62 * ((t) - T0) / (243.12 + (t) - T0)))
1469
1491#define PSICE(t) \
1492 (6.112 * exp(22.46 * ((t) - T0) / (272.62 + (t) - T0)))
1493
1518#define PW(p, h2o) \
1519 ((p) * MAX((h2o), 0.1e-6) / (1. + (1. - EPS) * MAX((h2o), 0.1e-6)))
1520
1535#define RAD2DEG(rad) \
1536 ((rad) * (180.0 / M_PI))
1537
1565#define RH(p, t, h2o) \
1566 (PW(p, h2o) / PSAT(t) * 100.)
1567
1595#define RHICE(p, t, h2o) \
1596 (PW(p, h2o) / PSICE(t) * 100.)
1597
1620#define RHO(p, t) \
1621 (100. * (p) / (RA * (t)))
1622
1639#define SET_ATM(qnt, val) \
1640 if (ctl->qnt >= 0) \
1641 atm->q[ctl->qnt][ip] = val;
1642
1662#define SET_QNT(qnt, name, longname, unit) \
1663 if (strcasecmp(ctl->qnt_name[iq], name) == 0) { \
1664 ctl->qnt = iq; \
1665 sprintf(ctl->qnt_longname[iq], longname); \
1666 sprintf(ctl->qnt_unit[iq], unit); \
1667 } else
1668
1683#define SH(h2o) \
1684 (EPS * MAX((h2o), 0.1e-6))
1685
1696#define SQR(x) \
1697 ((x)*(x))
1698
1710#define SWAP(x, y, type) \
1711 do {type tmp = x; x = y; y = tmp;} while(0);
1712
1734#define TDEW(p, h2o) \
1735 (T0 + 243.12 * log(PW((p), (h2o)) / 6.112) \
1736 / (17.62 - log(PW((p), (h2o)) / 6.112)))
1737
1759#define TICE(p, h2o) \
1760 (T0 + 272.62 * log(PW((p), (h2o)) / 6.112) \
1761 / (22.46 - log(PW((p), (h2o)) / 6.112)))
1762
1783#define THETA(p, t) \
1784 ((t) * pow(1000. / (p), 0.286))
1785
1812#define THETAVIRT(p, t, h2o) \
1813 (TVIRT(THETA((p), (t)), MAX((h2o), 0.1e-6)))
1814
1833#define TOK(line, tok, format, var) { \
1834 if(((tok)=strtok((line), " \t"))) { \
1835 if(sscanf(tok, format, &(var))!=1) continue; \
1836 } else ERRMSG("Error while reading!"); \
1837 }
1838
1858#define TVIRT(t, h2o) \
1859 ((t) * (1. + (1. - EPS) * MAX((h2o), 0.1e-6)))
1860
1880#define Z(p) \
1881 (H0 * log(P0 / (p)))
1882
1911#define ZDIFF(lnp0, t0, h2o0, lnp1, t1, h2o1) \
1912 (RI / MA / G0 * 0.5 * (TVIRT((t0), (h2o0)) + TVIRT((t1), (h2o1))) \
1913 * ((lnp0) - (lnp1)))
1914
1930#define ZETA(ps, p, t) \
1931 (((p) / (ps) <= 0.3 ? 1. : \
1932 sin(M_PI / 2. * (1. - (p) / (ps)) / (1. - 0.3))) \
1933 * THETA((p), (t)))
1934
1935/* ------------------------------------------------------------
1936 Log messages...
1937 ------------------------------------------------------------ */
1938
1940#ifndef LOGLEV
1941#define LOGLEV 2
1942#endif
1943
1973#define LOG(level, ...) { \
1974 if(level >= 2) \
1975 printf(" "); \
1976 if(level <= LOGLEV) { \
1977 printf(__VA_ARGS__); \
1978 printf("\n"); \
1979 } \
1980 }
1981
2010#define WARN(...) { \
2011 printf("\nWarning (%s, %s, l%d): ", __FILE__, __func__, __LINE__); \
2012 LOG(0, __VA_ARGS__); \
2013 }
2014
2043#define ERRMSG(...) { \
2044 printf("\nError (%s, %s, l%d): ", __FILE__, __func__, __LINE__); \
2045 LOG(0, __VA_ARGS__); \
2046 exit(EXIT_FAILURE); \
2047 }
2048
2078#define PRINT(format, var) \
2079 printf("Print (%s, %s, l%d): %s= "format"\n", \
2080 __FILE__, __func__, __LINE__, #var, var);
2081
2082/* ------------------------------------------------------------
2083 Timers...
2084 ------------------------------------------------------------ */
2085
2087#define NTIMER 100
2088
2102#define PRINT_TIMERS \
2103 timer("END", "END", 1);
2104
2117#define SELECT_TIMER(id, group) \
2118 timer(id, group, 0);
2119
2120/* ------------------------------------------------------------
2121 Structs...
2122 ------------------------------------------------------------ */
2123
2131typedef struct {
2132
2133 /* ------------------------------------------------------------
2134 Quantity parameters...
2135 ------------------------------------------------------------ */
2136
2138 int nq;
2139
2141 char qnt_name[NQ][LEN];
2142
2144 char qnt_longname[NQ][LEN];
2145
2147 char qnt_unit[NQ][LEN];
2148
2150 char qnt_format[NQ][LEN];
2151
2154
2157
2160
2163
2166
2169
2172
2175
2178
2181
2184
2187
2190
2193
2196
2199
2202
2205
2208
2211
2214
2217
2220
2223
2226
2229
2232
2235
2238
2241
2244
2247
2250
2253
2256
2259
2262
2265
2268
2271
2274
2277
2280
2283
2286
2289
2292
2295
2298
2301
2304
2307
2310
2313
2316
2319
2322
2325
2328
2331
2334
2337
2340
2343
2346
2349
2352
2355
2358
2361
2364
2367
2370
2373
2376
2379
2382
2385
2388
2391
2394
2397
2400
2403
2406
2409
2412
2415
2418
2421
2424
2427
2430
2433
2436
2439
2442
2445
2448
2451
2454
2457
2460
2462 double t_start;
2463
2465 double t_stop;
2466
2468 double dt_mod;
2469
2470 /* ------------------------------------------------------------
2471 Meteo data parameters...
2472 ------------------------------------------------------------ */
2473
2475 char metbase[LEN];
2476
2478 double dt_met;
2479
2482
2486
2490
2493
2496
2499
2502
2505
2507 int met_comp_prec[METVAR];
2508
2510 double met_comp_tol[METVAR];
2511
2514
2517
2520
2523
2526
2529
2532
2535
2538
2541
2544
2547
2550
2553
2556
2559
2562
2565
2568
2571
2574
2577
2580
2583
2586
2589
2591 double met_p[EP];
2592
2595
2598
2600 double met_lev_hyam[EP];
2601
2603 double met_lev_hybm[EP];
2604
2607
2610
2613
2616
2619
2622
2625
2629
2632
2635
2638
2641
2644
2647
2648 /* ------------------------------------------------------------
2649 Geophysical module parameters...
2650 ------------------------------------------------------------ */
2651
2653 double sort_dt;
2654
2658
2660 char balloon[LEN];
2661
2664
2668
2671
2674
2677
2680
2683
2686
2689
2692
2695
2698
2701
2704
2707
2709 double conv_cin;
2710
2712 double conv_dt;
2713
2716
2719
2722
2725
2728
2731
2733 double bound_p0;
2734
2736 double bound_p1;
2737
2740
2743
2746
2749
2751 char species[LEN];
2752
2754 double molmass;
2755
2758
2761
2764
2766 char clim_hno3_filename[LEN];
2767
2769 char clim_oh_filename[LEN];
2770
2772 char clim_h2o2_filename[LEN];
2773
2775 char clim_ho2_filename[LEN];
2776
2778 char clim_o1d_filename[LEN];
2779
2781 char clim_o3_filename[LEN];
2782
2784 char clim_ccl4_timeseries[LEN];
2785
2787 char clim_ccl3f_timeseries[LEN];
2788
2790 char clim_ccl2f2_timeseries[LEN];
2791
2793 char clim_n2o_timeseries[LEN];
2794
2796 char clim_sf6_timeseries[LEN];
2797
2800
2803
2806
2809
2812
2815
2818
2821
2824
2827
2830
2833
2836
2839
2842
2845
2848
2851
2854
2857
2860
2863
2865 double oh_chem[4];
2866
2869
2872
2875
2877 double dt_kpp;
2878
2881
2884
2886 double wet_depo_pre[2];
2887
2890
2893
2896
2899
2901 double wet_depo_ic_h[2];
2902
2904 double wet_depo_bc_h[2];
2905
2908
2911
2914
2917
2920
2922 double psc_h2o;
2923
2925 double psc_hno3;
2926
2927 /* ------------------------------------------------------------
2928 Output parameters...
2929 ------------------------------------------------------------ */
2930
2932 char atm_basename[LEN];
2933
2935 char atm_gpfile[LEN];
2936
2939
2942
2945
2949
2954
2957
2959 int atm_nc_quant[NQ];
2960
2963
2965 char csi_basename[LEN];
2966
2968 char csi_kernel[LEN];
2969
2972
2974 char csi_obsfile[LEN];
2975
2978
2981
2984
2986 double csi_z0;
2987
2989 double csi_z1;
2990
2993
2995 double csi_lon0;
2996
2998 double csi_lon1;
2999
3002
3004 double csi_lat0;
3005
3007 double csi_lat1;
3008
3010 int nens;
3011
3013 char ens_basename[LEN];
3014
3017
3019 char grid_basename[LEN];
3020
3022 char grid_kernel[LEN];
3023
3025 char grid_gpfile[LEN];
3026
3029
3032
3035
3037 int grid_nc_quant[NQ];
3038
3041
3044
3046 double grid_z0;
3047
3049 double grid_z1;
3050
3053
3056
3059
3062
3065
3068
3071
3073 char prof_basename[LEN];
3074
3076 char prof_obsfile[LEN];
3077
3080
3082 double prof_z0;
3083
3085 double prof_z1;
3086
3089
3092
3095
3098
3101
3104
3106 char sample_basename[LEN];
3107
3109 char sample_kernel[LEN];
3110
3112 char sample_obsfile[LEN];
3113
3116
3119
3121 char stat_basename[LEN];
3122
3124 double stat_lon;
3125
3127 double stat_lat;
3128
3130 double stat_r;
3131
3133 double stat_t0;
3134
3136 double stat_t1;
3137
3139 char vtk_basename[LEN];
3140
3143
3146
3149
3152
3155
3156 /* ------------------------------------------------------------
3157 Domain decomposition...
3158 ------------------------------------------------------------ */
3159
3161 int dd;
3162
3165
3168
3171
3174
3175} ctl_t;
3176
3185typedef struct {
3186
3188 int np;
3189
3191 double time[NP];
3192
3194 double p[NP];
3195
3197 double lon[NP];
3198
3200 double lat[NP];
3201
3203 double q[NQ][NP];
3204
3205} atm_t;
3206
3214typedef struct {
3215
3217 double time;
3218
3220 double p;
3221
3223 double lon;
3224
3226 double lat;
3227
3229 double q[NQ];
3230
3231} particle_t;
3232
3233
3240typedef struct {
3241
3243 double iso_var[NP];
3244
3246 double iso_ps[NP];
3247
3249 double iso_ts[NP];
3250
3253
3255 float uvwp[NP][3];
3256
3258 double rs[3 * NP + 1];
3259
3261 double dt[NP];
3262
3263} cache_t;
3264
3272typedef struct {
3273
3275 int np;
3276
3278 int nsza;
3279
3281 int no3c;
3282
3284 double p[CP];
3285
3287 double sza[CSZA];
3288
3290 double o3c[CO3];
3291
3293 double n2o[CP][CSZA][CO3];
3294
3296 double ccl4[CP][CSZA][CO3];
3297
3299 double ccl3f[CP][CSZA][CO3];
3300
3302 double ccl2f2[CP][CSZA][CO3];
3303
3305 double o2[CP][CSZA][CO3];
3306
3308 double o3_1[CP][CSZA][CO3];
3309
3311 double o3_2[CP][CSZA][CO3];
3312
3314 double h2o2[CP][CSZA][CO3];
3315
3317 double h2o[CP][CSZA][CO3];
3318
3319} clim_photo_t;
3320
3328typedef struct {
3329
3332
3334 double time[CTS];
3335
3337 double vmr[CTS];
3338
3339} clim_ts_t;
3340
3348typedef struct {
3349
3352
3354 int nlat;
3355
3357 int np;
3358
3360 double time[CT];
3361
3363 double lat[CY];
3364
3366 double p[CP];
3367
3369 double vmr[CT][CP][CY];
3370
3371} clim_zm_t;
3372
3380typedef struct {
3381
3384
3387
3389 double tropo_time[12];
3390
3392 double tropo_lat[73];
3393
3395 double tropo[12][73];
3396
3399
3402
3405
3408
3411
3414
3417
3420
3423
3426
3429
3430} clim_t;
3431
3439typedef struct {
3440
3442 double time;
3443
3445 int nx;
3446
3448 int ny;
3449
3451 int np;
3452
3454 int npl;
3455
3457 double lon[EX];
3458
3460 double lat[EY];
3461
3463 double p[EP];
3464
3466 double hybrid[EP];
3467
3469 double hyam[EP];
3470
3472 double hybm[EP];
3473
3475 double eta[EP];
3476
3478 float ps[EX][EY];
3479
3481 float ts[EX][EY];
3482
3484 float zs[EX][EY];
3485
3487 float us[EX][EY];
3488
3490 float vs[EX][EY];
3491
3493 float ess[EX][EY];
3494
3496 float nss[EX][EY];
3497
3499 float shf[EX][EY];
3500
3502 float lsm[EX][EY];
3503
3505 float sst[EX][EY];
3506
3508 float pbl[EX][EY];
3509
3511 float pt[EX][EY];
3512
3514 float tt[EX][EY];
3515
3517 float zt[EX][EY];
3518
3520 float h2ot[EX][EY];
3521
3523 float pct[EX][EY];
3524
3526 float pcb[EX][EY];
3527
3529 float cl[EX][EY];
3530
3532 float plcl[EX][EY];
3533
3535 float plfc[EX][EY];
3536
3538 float pel[EX][EY];
3539
3541 float cape[EX][EY];
3542
3544 float cin[EX][EY];
3545
3547 float o3c[EX][EY];
3548
3550 float z[EX][EY][EP];
3551
3553 float t[EX][EY][EP];
3554
3556 float u[EX][EY][EP];
3557
3559 float v[EX][EY][EP];
3560
3562 float w[EX][EY][EP];
3563
3565 float pv[EX][EY][EP];
3566
3568 float h2o[EX][EY][EP];
3569
3571 float o3[EX][EY][EP];
3572
3574 float lwc[EX][EY][EP];
3575
3577 float rwc[EX][EY][EP];
3578
3580 float iwc[EX][EY][EP];
3581
3583 float swc[EX][EY][EP];
3584
3586 float cc[EX][EY][EP];
3587
3589 float pl[EX][EY][EP];
3590
3592 float ul[EX][EY][EP];
3593
3595 float vl[EX][EY][EP];
3596
3598 float wl[EX][EY][EP];
3599
3601 float zetal[EX][EY][EP];
3602
3604 float zeta_dotl[EX][EY][EP];
3605
3606} met_t;
3607
3613typedef struct {
3614
3615 /* ------------------------------------------------------------
3616 MPI Information
3617 ------------------------------------------------------------ */
3618
3620 int rank;
3621
3623 int size;
3624
3626 int neighbours[DD_NNMAX];
3627
3628#ifdef DD
3630 MPI_Datatype MPI_Particle;
3631#endif
3632
3633 /* ------------------------------------------------------------
3634 Properties of subdomains
3635 ------------------------------------------------------------ */
3636
3639
3642
3645
3648
3650 size_t subdomain_start[4];
3651
3653 size_t subdomain_count[4];
3654
3656 size_t halo_bnd_start[4];
3657
3659 size_t halo_bnd_count[4];
3660
3663
3666
3667 /* ------------------------------------------------------------
3668 Caches
3669 ------------------------------------------------------------ */
3670
3672 int init;
3673
3674#ifdef DD
3675
3677 double a[NP];
3678
3680 int p[NP];
3681
3683 double help[NP];
3684
3685#endif
3686
3687} dd_t;
3688
3689/* ------------------------------------------------------------
3690 OpenACC routines...
3691 ------------------------------------------------------------ */
3692
3693#ifdef _OPENACC
3694#pragma acc routine (clim_oh)
3695#pragma acc routine (clim_photo)
3696#pragma acc routine (clim_tropo)
3697#pragma acc routine (clim_ts)
3698#pragma acc routine (clim_zm)
3699#pragma acc routine (cos_sza)
3700#pragma acc routine (intpol_check_lon_lat)
3701#pragma acc routine (intpol_met_4d_zeta)
3702#pragma acc routine (intpol_met_space_3d)
3703#pragma acc routine (intpol_met_space_2d)
3704#pragma acc routine (intpol_met_time_3d)
3705#pragma acc routine (intpol_met_time_2d)
3706#pragma acc routine (kernel_weight)
3707#pragma acc routine (lapse_rate)
3708#pragma acc routine (locate_irr)
3709#pragma acc routine (locate_irr_float)
3710#pragma acc routine (locate_reg)
3711#pragma acc routine (locate_vert)
3712#pragma acc routine (nat_temperature)
3713#pragma acc routine (pbl_weight)
3714#pragma acc routine (sedi)
3715#pragma acc routine (stddev)
3716#pragma acc routine (tropo_weight)
3717#endif
3718
3719/* ------------------------------------------------------------
3720 Functions...
3721 ------------------------------------------------------------ */
3722
3745 void *data,
3746 size_t N);
3747
3762void cart2geo(
3763 const double *x,
3764 double *z,
3765 double *lon,
3766 double *lat);
3767
3790double clim_oh(
3791 const ctl_t * ctl,
3792 const clim_t * clim,
3793 const double t,
3794 const double lon,
3795 const double lat,
3796 const double p);
3797
3817 const ctl_t * ctl,
3818 clim_t * clim);
3819
3849double clim_photo(
3850 const double rate[CP][CSZA][CO3],
3851 const clim_photo_t * photo,
3852 const double p,
3853 const double sza,
3854 const double o3c);
3855
3881double clim_tropo(
3882 const clim_t * clim,
3883 const double t,
3884 const double lat);
3885
3904void clim_tropo_init(
3905 clim_t * clim);
3906
3923double clim_ts(
3924 const clim_ts_t * ts,
3925 const double t);
3926
3948double clim_zm(
3949 const clim_zm_t * zm,
3950 const double t,
3951 const double lat,
3952 const double p);
3953
3998 const ctl_t * ctl,
3999 const char *varname,
4000 float *array,
4001 const size_t nx,
4002 const size_t ny,
4003 const size_t np,
4004 const double *plev,
4005 const int decompress,
4006 FILE * inout);
4007
4040void compress_pck(
4041 const char *varname,
4042 float *array,
4043 const size_t nxy,
4044 const size_t nz,
4045 const int decompress,
4046 FILE * inout);
4047
4078 const char *varname,
4079 float *array,
4080 const int nx,
4081 const int ny,
4082 const int nz,
4083 const int precision,
4084 const double tolerance,
4085 const int decompress,
4086 FILE * inout);
4087
4127 const char *varname,
4128 float *array,
4129 const int nx,
4130 const int ny,
4131 const int nz,
4132 const int precision,
4133 const double tolerance,
4134 const int decompress,
4135 FILE * inout);
4136
4159 const char *varname,
4160 float *array,
4161 const size_t n,
4162 const int decompress,
4163 const int level,
4164 FILE * inout);
4165
4190double cos_sza(
4191 const double sec,
4192 const double lon,
4193 const double lat);
4194
4217void day2doy(
4218 const int year,
4219 const int mon,
4220 const int day,
4221 int *doy);
4222
4270 atm_t * atm,
4271 ctl_t * ctl,
4272 dd_t * dd,
4273 int init);
4274
4340 atm_t * atm,
4341 particle_t * particles,
4342 ctl_t * ctl,
4343 int *nparticles,
4344 cache_t * cache,
4345 int rank);
4346
4411 double lon,
4412 double lat,
4413 met_t * met,
4414 ctl_t * ctl,
4415 int mpi_size,
4416 int nx_glob,
4417 int ny_glob);
4418
4449 particle_t * particles,
4450 int *nparticles,
4451 MPI_Datatype MPI_Particle,
4452 int *neighbours,
4453 int nneighbours,
4454 ctl_t ctl);
4455
4479 const ctl_t ctl,
4480 dd_t * dd);
4481
4507 ctl_t * ctl,
4508 dd_t * dd,
4509 atm_t * atm);
4510
4535 met_t * met,
4536 int nx_glob);
4537
4563 atm_t * atm,
4564 particle_t * particles,
4565 ctl_t * ctl,
4566 int *nparticles,
4567 cache_t * cache);
4568
4591 MPI_Datatype * MPI_Particle);
4592
4628 const ctl_t * ctl,
4629 met_t * met0,
4630 atm_t * atm,
4631 dd_t * dd,
4632 int *nparticles,
4633 int *rank);
4634
4663 double *a,
4664 dd_t * dd,
4665 const int np);
4666
4688void doy2day(
4689 const int year,
4690 const int doy,
4691 int *mon,
4692 int *day);
4693
4720void fft_help(
4721 double *fcReal,
4722 double *fcImag,
4723 const int n);
4724
4751void geo2cart(
4752 const double z,
4753 const double lon,
4754 const double lat,
4755 double *x);
4756
4781void get_met_help(
4782 const ctl_t * ctl,
4783 const double t,
4784 const int direct,
4785 const char *metbase,
4786 const double dt_met,
4787 char *filename);
4788
4812void get_met_replace(
4813 char *orig,
4814 char *search,
4815 char *repl);
4816
4853void get_tropo(
4854 const int met_tropo,
4855 ctl_t * ctl,
4856 clim_t * clim,
4857 met_t * met,
4858 const double *lons,
4859 const int nx,
4860 const double *lats,
4861 const int ny,
4862 double *pt,
4863 double *zt,
4864 double *tt,
4865 double *qt,
4866 double *o3t,
4867 double *ps,
4868 double *zs);
4869
4891 const double *lons,
4892 const int nlon,
4893 const double *lats,
4894 const int nlat,
4895 const double lon,
4896 const double lat,
4897 double *lon2,
4898 double *lat2);
4899
4942 const met_t * met0,
4943 float height0[EX][EY][EP],
4944 float array0[EX][EY][EP],
4945 const met_t * met1,
4946 float height1[EX][EY][EP],
4947 float array1[EX][EY][EP],
4948 const double ts,
4949 const double height,
4950 const double lon,
4951 const double lat,
4952 double *var,
4953 int *ci,
4954 double *cw,
4955 const int init);
4956
4992 const met_t * met,
4993 float array[EX][EY][EP],
4994 const double p,
4995 const double lon,
4996 const double lat,
4997 double *var,
4998 int *ci,
4999 double *cw,
5000 const int init);
5001
5037 const met_t * met,
5038 float array[EX][EY],
5039 const double lon,
5040 const double lat,
5041 double *var,
5042 int *ci,
5043 double *cw,
5044 const int init);
5045
5080 const met_t * met0,
5081 float array0[EX][EY][EP],
5082 const met_t * met1,
5083 float array1[EX][EY][EP],
5084 const double ts,
5085 const double p,
5086 const double lon,
5087 const double lat,
5088 double *var,
5089 int *ci,
5090 double *cw,
5091 const int init);
5092
5128 const met_t * met0,
5129 float array0[EX][EY],
5130 const met_t * met1,
5131 float array1[EX][EY],
5132 const double ts,
5133 const double lon,
5134 const double lat,
5135 double *var,
5136 int *ci,
5137 double *cw,
5138 const int init);
5139
5177void intpol_tropo_3d(
5178 const double time0,
5179 float array0[EX][EY],
5180 const double time1,
5181 float array1[EX][EY],
5182 const double lons[EX],
5183 const double lats[EY],
5184 const int nlon,
5185 const int nlat,
5186 const double time,
5187 const double lon,
5188 const double lat,
5189 const int method,
5190 double *var,
5191 double *sigma);
5192
5219void jsec2time(
5220 const double jsec,
5221 int *year,
5222 int *mon,
5223 int *day,
5224 int *hour,
5225 int *min,
5226 int *sec,
5227 double *remain);
5228
5255double kernel_weight(
5256 const double kz[EP],
5257 const double kw[EP],
5258 const int nk,
5259 const double p);
5260
5299double lapse_rate(
5300 const double t,
5301 const double h2o);
5302
5332 ctl_t * ctl);
5333
5353int locate_irr(
5354 const double *xx,
5355 const int n,
5356 const double x);
5357
5384 const float *xx,
5385 const int n,
5386 const double x,
5387 const int ig);
5388
5409int locate_reg(
5410 const double *xx,
5411 const int n,
5412 const double x);
5413
5435void locate_vert(
5436 float profiles[EX][EY][EP],
5437 const int np,
5438 const int lon_ap_ind,
5439 const int lat_ap_ind,
5440 const double alt_ap,
5441 int *ind);
5442
5468void module_advect(
5469 const ctl_t * ctl,
5470 const cache_t * cache,
5471 met_t * met0,
5472 met_t * met1,
5473 atm_t * atm);
5474
5498 const ctl_t * ctl,
5499 const cache_t * cache,
5500 met_t * met0,
5501 met_t * met1,
5502 atm_t * atm);
5503
5541 const ctl_t * ctl,
5542 const cache_t * cache,
5543 const clim_t * clim,
5544 met_t * met0,
5545 met_t * met1,
5546 atm_t * atm);
5547
5589void module_chem_grid(
5590 const ctl_t * ctl,
5591 met_t * met0,
5592 met_t * met1,
5593 atm_t * atm,
5594 const double t);
5595
5623void module_chem_init(
5624 const ctl_t * ctl,
5625 const cache_t * cache,
5626 const clim_t * clim,
5627 met_t * met0,
5628 met_t * met1,
5629 atm_t * atm);
5630
5655 const ctl_t * ctl,
5656 cache_t * cache,
5657 met_t * met0,
5658 met_t * met1,
5659 atm_t * atm);
5660
5689 ctl_t * ctl,
5690 atm_t * atm,
5691 cache_t * cache,
5692 dd_t * dd,
5693 met_t ** met);
5694
5721void module_decay(
5722 const ctl_t * ctl,
5723 const cache_t * cache,
5724 const clim_t * clim,
5725 atm_t * atm);
5726
5763void module_diff_meso(
5764 const ctl_t * ctl,
5765 cache_t * cache,
5766 met_t * met0,
5767 met_t * met1,
5768 atm_t * atm);
5769
5803void module_diff_pbl(
5804 const ctl_t * ctl,
5805 cache_t * cache,
5806 met_t * met0,
5807 met_t * met1,
5808 atm_t * atm);
5809
5864void module_diff_turb(
5865 const ctl_t * ctl,
5866 cache_t * cache,
5867 const clim_t * clim,
5868 met_t * met0,
5869 met_t * met1,
5870 atm_t * atm);
5871
5891void module_dry_depo(
5892 const ctl_t * ctl,
5893 const cache_t * cache,
5894 met_t * met0,
5895 met_t * met1,
5896 atm_t * atm);
5897
5930void module_h2o2_chem(
5931 const ctl_t * ctl,
5932 const cache_t * cache,
5933 const clim_t * clim,
5934 met_t * met0,
5935 met_t * met1,
5936 atm_t * atm);
5937
5958 const ctl_t * ctl,
5959 cache_t * cache,
5960 met_t * met0,
5961 met_t * met1,
5962 atm_t * atm);
5963
5981void module_isosurf(
5982 const ctl_t * ctl,
5983 const cache_t * cache,
5984 met_t * met0,
5985 met_t * met1,
5986 atm_t * atm);
5987
6020 ctl_t * ctl,
6021 cache_t * cache,
6022 clim_t * clim,
6023 met_t * met0,
6024 met_t * met1,
6025 atm_t * atm);
6026
6045void module_meteo(
6046 const ctl_t * ctl,
6047 const cache_t * cache,
6048 const clim_t * clim,
6049 met_t * met0,
6050 met_t * met1,
6051 atm_t * atm);
6052
6070void module_mixing(
6071 const ctl_t * ctl,
6072 const clim_t * clim,
6073 atm_t * atm,
6074 const double t);
6075
6104 const ctl_t * ctl,
6105 const clim_t * clim,
6106 atm_t * atm,
6107 const int *ixs,
6108 const int *iys,
6109 const int *izs,
6110 const int qnt_idx,
6111 const int use_ensemble);
6112
6145void module_oh_chem(
6146 const ctl_t * ctl,
6147 const cache_t * cache,
6148 const clim_t * clim,
6149 met_t * met0,
6150 met_t * met1,
6151 atm_t * atm);
6152
6180void module_position(
6181 const cache_t * cache,
6182 met_t * met0,
6183 met_t * met1,
6184 atm_t * atm);
6185
6210void module_rng_init(
6211 const int ntask);
6212
6238void module_rng(
6239 const ctl_t * ctl,
6240 double *rs,
6241 const size_t n,
6242 const int method);
6243
6279 const ctl_t * ctl,
6280 const cache_t * cache,
6281 atm_t * atm);
6282
6305void module_sedi(
6306 const ctl_t * ctl,
6307 const cache_t * cache,
6308 met_t * met0,
6309 met_t * met1,
6310 atm_t * atm);
6311
6335void module_sort(
6336 const ctl_t * ctl,
6337 met_t * met0,
6338 atm_t * atm);
6339
6359void module_sort_help(
6360 double *a,
6361 const int *p,
6362 const int np);
6363
6387void module_timesteps(
6388 const ctl_t * ctl,
6389 cache_t * cache,
6390 met_t * met0,
6391 atm_t * atm,
6392 const double t);
6393
6415 ctl_t * ctl,
6416 const atm_t * atm);
6417
6451 const ctl_t * ctl,
6452 const cache_t * cache,
6453 const clim_t * clim,
6454 met_t * met0,
6455 met_t * met1,
6456 atm_t * atm);
6457
6487void module_wet_depo(
6488 const ctl_t * ctl,
6489 const cache_t * cache,
6490 met_t * met0,
6491 met_t * met1,
6492 atm_t * atm);
6493
6525void mptrac_alloc(
6526 ctl_t ** ctl,
6527 cache_t ** cache,
6528 clim_t ** clim,
6529 met_t ** met0,
6530 met_t ** met1,
6531 atm_t ** atm,
6532 dd_t ** dd);
6533
6564void mptrac_free(
6565 ctl_t * ctl,
6566 cache_t * cache,
6567 clim_t * clim,
6568 met_t * met0,
6569 met_t * met1,
6570 atm_t * atm,
6571 dd_t * dd);
6572
6608void mptrac_get_met(
6609 ctl_t * ctl,
6610 clim_t * clim,
6611 const double t,
6612 met_t ** met0,
6613 met_t ** met1,
6614 dd_t * dd);
6615
6635void mptrac_init(
6636 ctl_t * ctl,
6637 cache_t * cache,
6638 clim_t * clim,
6639 atm_t * atm,
6640 const int ntask);
6641
6677int mptrac_read_atm(
6678 const char *filename,
6679 const ctl_t * ctl,
6680 atm_t * atm);
6681
6713void mptrac_read_clim(
6714 const ctl_t * ctl,
6715 clim_t * clim);
6716
6746void mptrac_read_ctl(
6747 const char *filename,
6748 int argc,
6749 char *argv[],
6750 ctl_t * ctl);
6751
6782int mptrac_read_met(
6783 const char *filename,
6784 const ctl_t * ctl,
6785 const clim_t * clim,
6786 met_t * met,
6787 dd_t * dd);
6788
6810 ctl_t * ctl,
6811 cache_t * cache,
6812 clim_t * clim,
6813 met_t ** met0,
6814 met_t ** met1,
6815 atm_t * atm,
6816 double t,
6817 dd_t * dd);
6818
6848void mptrac_write_atm(
6849 const char *filename,
6850 const ctl_t * ctl,
6851 const atm_t * atm,
6852 const double t);
6853
6891void mptrac_write_met(
6892 const char *filename,
6893 const ctl_t * ctl,
6894 met_t * met);
6895
6930 const char *dirname,
6931 const ctl_t * ctl,
6932 met_t * met0,
6933 met_t * met1,
6934 atm_t * atm,
6935 const double t);
6936
6968 const ctl_t * ctl,
6969 const cache_t * cache,
6970 const clim_t * clim,
6971 met_t ** met0,
6972 met_t ** met1,
6973 const atm_t * atm);
6974
7005 const ctl_t * ctl,
7006 const cache_t * cache,
7007 const clim_t * clim,
7008 met_t ** met0,
7009 met_t ** met1,
7010 const atm_t * atm);
7011
7039double nat_temperature(
7040 const double p,
7041 const double h2o,
7042 const double hno3);
7043
7064double pbl_weight(
7065 const ctl_t * ctl,
7066 const atm_t * atm,
7067 const int ip,
7068 const double pbl,
7069 const double ps);
7070
7103int read_atm_asc(
7104 const char *filename,
7105 const ctl_t * ctl,
7106 atm_t * atm);
7107
7138int read_atm_bin(
7139 const char *filename,
7140 const ctl_t * ctl,
7141 atm_t * atm);
7142
7167int read_atm_clams(
7168 const char *filename,
7169 const ctl_t * ctl,
7170 atm_t * atm);
7171
7201int read_atm_nc(
7202 const char *filename,
7203 const ctl_t * ctl,
7204 atm_t * atm);
7205
7234void read_clim_photo(
7235 const char *filename,
7236 clim_photo_t * photo);
7237
7255 const int ncid,
7256 const char *varname,
7257 const clim_photo_t * photo,
7258 double var[CP][CSZA][CO3]);
7259
7283int read_clim_ts(
7284 const char *filename,
7285 clim_ts_t * ts);
7286
7313void read_clim_zm(
7314 const char *filename,
7315 const char *varname,
7316 clim_zm_t * zm);
7317
7345void read_kernel(
7346 const char *filename,
7347 double kz[EP],
7348 double kw[EP],
7349 int *nk);
7350
7382int read_met_bin(
7383 const char *filename,
7384 const ctl_t * ctl,
7385 met_t * met);
7386
7412void read_met_bin_2d(
7413 FILE * in,
7414 const met_t * met,
7415 float var[EX][EY],
7416 const char *varname);
7417
7455void read_met_bin_3d(
7456 FILE * in,
7457 const ctl_t * ctl,
7458 const met_t * met,
7459 float var[EX][EY][EP],
7460 const char *varname,
7461 const float bound_min,
7462 const float bound_max);
7463
7491void read_met_cape(
7492 const ctl_t * ctl,
7493 const clim_t * clim,
7494 met_t * met);
7495
7518void read_met_cloud(
7519 met_t * met);
7520
7546void read_met_detrend(
7547 const ctl_t * ctl,
7548 met_t * met);
7549
7573 met_t * met);
7574
7601void read_met_geopot(
7602 const ctl_t * ctl,
7603 met_t * met);
7604
7632 const char *filename,
7633 const ctl_t * ctl,
7634 met_t * met);
7635
7657 codes_handle ** handles,
7658 int count_handles,
7659 met_t * met);
7660
7681 codes_handle ** handles,
7682 const int num_messages,
7683 const ctl_t * ctl,
7684 met_t * met);
7685
7716 codes_handle ** handles,
7717 const int num_messages,
7718 const ctl_t * ctl,
7719 met_t * met);
7720
7749void read_met_ml2pl(
7750 const ctl_t * ctl,
7751 const met_t * met,
7752 float var[EX][EY][EP],
7753 const char *varname);
7754
7777 const ctl_t * ctl,
7778 met_t * met);
7779
7810int read_met_nc(
7811 const char *filename,
7812 const ctl_t * ctl,
7813 met_t * met,
7814 dd_t * dd);
7815
7850void read_met_nc_grid(
7851 const char *filename,
7852 const int ncid,
7853 const ctl_t * ctl,
7854 met_t * met,
7855 dd_t * dd);
7856
7902 dd_t * dd,
7903 const ctl_t * ctl,
7904 met_t * met,
7905 const int ncid);
7906
7939 const int ncid,
7940 const ctl_t * ctl,
7941 met_t * met,
7942 dd_t * dd);
7943
7974 const int ncid,
7975 const ctl_t * ctl,
7976 met_t * met,
7977 dd_t * dd);
7978
8011int read_met_nc_2d(
8012 const int ncid,
8013 const char *varname,
8014 const char *varname2,
8015 const char *varname3,
8016 const char *varname4,
8017 const char *varname5,
8018 const char *varname6,
8019 const ctl_t * ctl,
8020 const met_t * met,
8021 dd_t * dd,
8022 float dest[EX][EY],
8023 const float scl,
8024 const int init);
8025
8055int read_met_nc_3d(
8056 const int ncid,
8057 const char *varname,
8058 const char *varname2,
8059 const char *varname3,
8060 const char *varname4,
8061 const ctl_t * ctl,
8062 const met_t * met,
8063 dd_t * dd,
8064 float dest[EX][EY][EP],
8065 const float scl);
8066
8112void read_met_pbl(
8113 const ctl_t * ctl,
8114 met_t * met);
8115
8148 met_t * met);
8149
8180 met_t * met);
8181
8212void read_met_pv(
8213 met_t * met);
8214
8237void read_met_ozone(
8238 met_t * met);
8239
8268void read_met_sample(
8269 const ctl_t * ctl,
8270 met_t * met);
8271
8300void read_met_tropo(
8301 const ctl_t * ctl,
8302 const clim_t * clim,
8303 met_t * met);
8304
8336void read_obs(
8337 const char *filename,
8338 const ctl_t * ctl,
8339 double *rt,
8340 double *rz,
8341 double *rlon,
8342 double *rlat,
8343 double *robs,
8344 int *nobs);
8345
8373void read_obs_asc(
8374 const char *filename,
8375 double *rt,
8376 double *rz,
8377 double *rlon,
8378 double *rlat,
8379 double *robs,
8380 int *nobs);
8381
8408void read_obs_nc(
8409 const char *filename,
8410 double *rt,
8411 double *rz,
8412 double *rlon,
8413 double *rlat,
8414 double *robs,
8415 int *nobs);
8416
8450double scan_ctl(
8451 const char *filename,
8452 int argc,
8453 char *argv[],
8454 const char *varname,
8455 const int arridx,
8456 const char *defvalue,
8457 char *value);
8458
8485double sedi(
8486 const double p,
8487 const double T,
8488 const double rp,
8489 const double rhop);
8490
8520void spline(
8521 const double *x,
8522 const double *y,
8523 const int n,
8524 const double *x2,
8525 double *y2,
8526 const int n2,
8527 const int method);
8528
8551float stddev(
8552 const float *data,
8553 const int n);
8554
8579void time2jsec(
8580 const int year,
8581 const int mon,
8582 const int day,
8583 const int hour,
8584 const int min,
8585 const int sec,
8586 const double remain,
8587 double *jsec);
8588
8617void timer(
8618 const char *name,
8619 const char *group,
8620 const int output);
8621
8647double time_from_filename(
8648 const char *filename,
8649 const int offset);
8650
8668double tropo_weight(
8669 const clim_t * clim,
8670 const atm_t * atm,
8671 const int ip);
8672
8695void write_atm_asc(
8696 const char *filename,
8697 const ctl_t * ctl,
8698 const atm_t * atm,
8699 const double t);
8700
8724void write_atm_bin(
8725 const char *filename,
8726 const ctl_t * ctl,
8727 const atm_t * atm);
8728
8752void write_atm_clams(
8753 const char *filename,
8754 const ctl_t * ctl,
8755 const atm_t * atm);
8756
8782 const char *dirname,
8783 const ctl_t * ctl,
8784 const atm_t * atm,
8785 const double t);
8786
8810void write_atm_nc(
8811 const char *filename,
8812 const ctl_t * ctl,
8813 const atm_t * atm);
8814
8843void write_csi(
8844 const char *filename,
8845 const ctl_t * ctl,
8846 const atm_t * atm,
8847 const double t);
8848
8890 const char *filename,
8891 const ctl_t * ctl,
8892 const atm_t * atm,
8893 const double t);
8894
8922void write_ens(
8923 const char *filename,
8924 const ctl_t * ctl,
8925 const atm_t * atm,
8926 const double t);
8927
8966void write_grid(
8967 const char *filename,
8968 const ctl_t * ctl,
8969 met_t * met0,
8970 met_t * met1,
8971 const atm_t * atm,
8972 const double t);
8973
9019void write_grid_asc(
9020 const char *filename,
9021 const ctl_t * ctl,
9022 const double *cd,
9023 double *mean[NQ],
9024 double *sigma[NQ],
9025 const double *vmr_impl,
9026 const double t,
9027 const double *z,
9028 const double *lon,
9029 const double *lat,
9030 const double *area,
9031 const double dz,
9032 const int *np);
9033
9076void write_grid_nc(
9077 const char *filename,
9078 const ctl_t * ctl,
9079 const double *cd,
9080 double *mean[NQ],
9081 double *sigma[NQ],
9082 const double *vmr_impl,
9083 const double t,
9084 const double *z,
9085 const double *lon,
9086 const double *lat,
9087 const double *area,
9088 const double dz,
9089 const int *np);
9090
9120void write_met_bin(
9121 const char *filename,
9122 const ctl_t * ctl,
9123 met_t * met);
9124
9152void write_met_bin_2d(
9153 FILE * out,
9154 met_t * met,
9155 float var[EX][EY],
9156 const char *varname);
9157
9195void write_met_bin_3d(
9196 FILE * out,
9197 const ctl_t * ctl,
9198 met_t * met,
9199 float var[EX][EY][EP],
9200 const char *varname,
9201 const int precision,
9202 const double tolerance);
9203
9231void write_met_nc(
9232 const char *filename,
9233 const ctl_t * ctl,
9234 met_t * met);
9235
9258void write_met_nc_2d(
9259 const int ncid,
9260 const char *varname,
9261 met_t * met,
9262 float var[EX][EY],
9263 const float scl);
9264
9288void write_met_nc_3d(
9289 const int ncid,
9290 const char *varname,
9291 met_t * met,
9292 float var[EX][EY][EP],
9293 const float scl);
9294
9325void write_prof(
9326 const char *filename,
9327 const ctl_t * ctl,
9328 met_t * met0,
9329 met_t * met1,
9330 const atm_t * atm,
9331 const double t);
9332
9364void write_sample(
9365 const char *filename,
9366 const ctl_t * ctl,
9367 met_t * met0,
9368 met_t * met1,
9369 const atm_t * atm,
9370 const double t);
9371
9398void write_station(
9399 const char *filename,
9400 const ctl_t * ctl,
9401 atm_t * atm,
9402 const double t);
9403
9432void write_vtk(
9433 const char *filename,
9434 const ctl_t * ctl,
9435 const atm_t * atm,
9436 const double t);
9437
9438#endif /* LIBTRAC_H */
void read_met_geopot(const ctl_t *ctl, met_t *met)
Calculates geopotential heights from meteorological data.
Definition: mptrac.c:8184
int dd_init(ctl_t *ctl, dd_t *dd, atm_t *atm)
Initializes domain decomposition for parallel processing.
void compress_zstd(const char *varname, float *array, const size_t n, const int decompress, const int level, FILE *inout)
Compresses or decompresses a float array using Zstandard (ZSTD).
#define LEN
Maximum length of ASCII data lines.
Definition: mptrac.h:298
void mptrac_alloc(ctl_t **ctl, cache_t **cache, clim_t **clim, met_t **met0, met_t **met1, atm_t **atm, dd_t **dd)
Allocates and initializes memory resources for MPTRAC.
Definition: mptrac.c:5203
void mptrac_write_atm(const char *filename, const ctl_t *ctl, const atm_t *atm, const double t)
Writes air parcel data to a file in various formats.
Definition: mptrac.c:6840
void day2doy(const int year, const int mon, const int day, int *doy)
Get day of year from date.
Definition: mptrac.c:1052
void read_met_grib_surface(codes_handle **handles, const int num_messages, const ctl_t *ctl, met_t *met)
Reads surface meteorological data from a grib file and stores it in the meteorological data structure...
void read_met_extrapolate(met_t *met)
Extrapolates meteorological data.
Definition: mptrac.c:8144
void write_atm_clams_traj(const char *dirname, const ctl_t *ctl, const atm_t *atm, const double t)
Writes CLaMS trajectory data to a NetCDF file.
Definition: mptrac.c:11517
int read_met_nc_2d(const int ncid, const char *varname, const char *varname2, const char *varname3, const char *varname4, const char *varname5, const char *varname6, const ctl_t *ctl, const met_t *met, dd_t *dd, float dest[EX][EY], const float scl, const int init)
Reads a 2-dimensional meteorological variable from a NetCDF file.
Definition: mptrac.c:8814
void write_met_nc_2d(const int ncid, const char *varname, met_t *met, float var[EX][EY], const float scl)
Writes a 2D meteorological variable to a NetCDF file.
Definition: mptrac.c:12916
void read_met_sample(const ctl_t *ctl, met_t *met)
Downsamples meteorological data based on specified parameters.
Definition: mptrac.c:10556
void write_met_bin_3d(FILE *out, const ctl_t *ctl, met_t *met, float var[EX][EY][EP], const char *varname, const int precision, const double tolerance)
Writes a 3-dimensional meteorological variable to a binary file.
Definition: mptrac.c:12662
void read_obs(const char *filename, const ctl_t *ctl, double *rt, double *rz, double *rlon, double *rlat, double *robs, int *nobs)
Reads observation data from a file and stores it in arrays.
Definition: mptrac.c:10901
void module_advect(const ctl_t *ctl, const cache_t *cache, met_t *met0, met_t *met1, atm_t *atm)
Advances particle positions using different advection schemes.
Definition: mptrac.c:2960
void module_timesteps(const ctl_t *ctl, cache_t *cache, met_t *met0, atm_t *atm, const double t)
Calculate time steps for air parcels based on specified conditions.
Definition: mptrac.c:4919
void module_meteo(const ctl_t *ctl, const cache_t *cache, const clim_t *clim, met_t *met0, met_t *met1, atm_t *atm)
Update atmospheric properties using meteorological data.
Definition: mptrac.c:4177
void read_clim_photo(const char *filename, clim_photo_t *photo)
Reads photolysis rates from a NetCDF file and populates the given photolysis structure.
Definition: mptrac.c:7272
void read_met_cloud(met_t *met)
Calculates cloud-related variables for each grid point.
Definition: mptrac.c:7983
void module_decay(const ctl_t *ctl, const cache_t *cache, const clim_t *clim, atm_t *atm)
Simulate exponential decay processes for atmospheric particles.
Definition: mptrac.c:3563
double sedi(const double p, const double T, const double rp, const double rhop)
Calculates the sedimentation velocity of a particle in air.
Definition: mptrac.c:11074
double cos_sza(const double sec, const double lon, const double lat)
Calculates the cosine of the solar zenith angle.
Definition: mptrac.c:1011
#define METVAR
Number of 3-D meteorological variables.
Definition: mptrac.h:303
void intpol_met_space_2d(const met_t *met, float array[EX][EY], const double lon, const double lat, double *var, int *ci, double *cw, const int init)
Interpolates meteorological variables in 2D space.
Definition: mptrac.c:2359
int read_met_nc_3d(const int ncid, const char *varname, const char *varname2, const char *varname3, const char *varname4, const ctl_t *ctl, const met_t *met, dd_t *dd, float dest[EX][EY][EP], const float scl)
Reads a 3-dimensional meteorological variable from a NetCDF file.
Definition: mptrac.c:9076
void dd_sort(const ctl_t *ctl, met_t *met0, atm_t *atm, dd_t *dd, int *nparticles, int *rank)
Sort particles according to box index and target rank for neighbours.
int read_atm_nc(const char *filename, const ctl_t *ctl, atm_t *atm)
Reads air parcel data from a generic netCDF file and populates the given atmospheric structure.
Definition: mptrac.c:7239
void write_csi_ens(const char *filename, const ctl_t *ctl, const atm_t *atm, const double t)
Writes ensemble-based Critical Success Index (CSI) and other verification statistics to an output fil...
void read_met_pbl(const ctl_t *ctl, met_t *met)
Computes the planetary boundary layer (PBL) pressure based on meteorological data.
Definition: mptrac.c:10164
void read_met_detrend(const ctl_t *ctl, met_t *met)
Detrends meteorological data.
Definition: mptrac.c:8040
void read_met_tropo(const ctl_t *ctl, const clim_t *clim, met_t *met)
Calculates the tropopause and related meteorological variables based on various methods and stores th...
Definition: mptrac.c:10729
void read_obs_asc(const char *filename, double *rt, double *rz, double *rlon, double *rlat, double *robs, int *nobs)
Reads observation data from an ASCII file.
Definition: mptrac.c:10945
void module_chem_init(const ctl_t *ctl, const cache_t *cache, const clim_t *clim, met_t *met0, met_t *met1, atm_t *atm)
Initializes the chemistry modules by setting atmospheric composition.
Definition: mptrac.c:3408
int locate_reg(const double *xx, const int n, const double x)
Locate the index of the interval containing a given value in a regular grid.
Definition: mptrac.c:2921
void get_tropo(const int met_tropo, ctl_t *ctl, clim_t *clim, met_t *met, const double *lons, const int nx, const double *lats, const int ny, double *pt, double *zt, double *tt, double *qt, double *o3t, double *ps, double *zs)
Calculate tropopause data.
Definition: mptrac.c:2059
void locate_vert(float profiles[EX][EY][EP], const int np, const int lon_ap_ind, const int lat_ap_ind, const double alt_ap, int *ind)
Locate the four vertical indizes of a box for a given height value.
Definition: mptrac.c:2940
void compress_cms(const ctl_t *ctl, const char *varname, float *array, const size_t nx, const size_t ny, const size_t np, const double *plev, const int decompress, FILE *inout)
Compresses or decompresses a 3-D meteorological field using cmultiscale.
void read_met_nc_levels(const int ncid, const ctl_t *ctl, met_t *met, dd_t *dd)
Reads and processes meteorological level data from NetCDF files with domain decomposition.
Definition: mptrac.c:8617
void read_met_monotonize(const ctl_t *ctl, met_t *met)
Makes zeta and pressure profiles monotone.
Definition: mptrac.c:9810
int read_clim_ts(const char *filename, clim_ts_t *ts)
Reads a climatological time series from a file and populates the given time series structure.
Definition: mptrac.c:7391
void read_met_periodic(met_t *met)
Applies periodic boundary conditions to meteorological data along longitudinal axis.
Definition: mptrac.c:10301
void mptrac_run_timestep(ctl_t *ctl, cache_t *cache, clim_t *clim, met_t **met0, met_t **met1, atm_t *atm, double t, dd_t *dd)
Executes a single timestep of the MPTRAC model simulation.
Definition: mptrac.c:6576
void module_timesteps_init(ctl_t *ctl, const atm_t *atm)
Initialize start time and time interval for time-stepping.
Definition: mptrac.c:4966
void write_ens(const char *filename, const ctl_t *ctl, const atm_t *atm, const double t)
Writes ensemble data to a file.
Definition: mptrac.c:11999
void module_mixing(const ctl_t *ctl, const clim_t *clim, atm_t *atm, const double t)
Update atmospheric properties through interparcel mixing.
Definition: mptrac.c:4284
void compress_zfp(const char *varname, float *array, const int nx, const int ny, const int nz, const int precision, const double tolerance, const int decompress, FILE *inout)
Compresses or decompresses a 3D array of floats using the ZFP library.
double clim_zm(const clim_zm_t *zm, const double t, const double lat, const double p)
Interpolates monthly mean zonal mean climatological variables.
Definition: mptrac.c:405
#define EY
Maximum number of latitudes for meteo data.
Definition: mptrac.h:293
void module_mixing_help(const ctl_t *ctl, const clim_t *clim, atm_t *atm, const int *ixs, const int *iys, const int *izs, const int qnt_idx, const int use_ensemble)
Perform subgrid-scale interparcel mixing of a given quantity.
Definition: mptrac.c:4356
void read_clim_photo_help(const int ncid, const char *varname, const clim_photo_t *photo, double var[CP][CSZA][CO3])
Reads a 3D climatological photochemistry variable from a NetCDF file.
Definition: mptrac.c:7363
void read_met_ml2pl(const ctl_t *ctl, const met_t *met, float var[EX][EY][EP], const char *varname)
Interpolates meteorological data to specified pressure levels.
Definition: mptrac.c:9768
double clim_tropo(const clim_t *clim, const double t, const double lat)
Calculates the tropopause pressure based on climatological data.
Definition: mptrac.c:204
void read_obs_nc(const char *filename, double *rt, double *rz, double *rlon, double *rlat, double *robs, int *nobs)
Reads observation data from a NetCDF file.
Definition: mptrac.c:10973
void dd_particles2atm(atm_t *atm, particle_t *particles, ctl_t *ctl, int *nparticles, cache_t *cache)
Converts particle data to atmospheric data.
void read_met_bin_2d(FILE *in, const met_t *met, float var[EX][EY], const char *varname)
Reads a 2-dimensional meteorological variable from a binary file and stores it in the provided array.
Definition: mptrac.c:7733
int locate_irr(const double *xx, const int n, const double x)
Locate the index of the interval containing a given value in a sorted array.
Definition: mptrac.c:2857
void module_isosurf_init(const ctl_t *ctl, cache_t *cache, met_t *met0, met_t *met1, atm_t *atm)
Initialize the isosurface module based on atmospheric data.
Definition: mptrac.c:4001
void level_definitions(ctl_t *ctl)
Defines pressure levels for meteorological data.
Definition: mptrac.c:2647
void intpol_met_4d_zeta(const met_t *met0, float height0[EX][EY][EP], float array0[EX][EY][EP], const met_t *met1, float height1[EX][EY][EP], float array1[EX][EY][EP], const double ts, const double height, const double lon, const double lat, double *var, int *ci, double *cw, const int init)
Interpolates meteorological variables to a given position and time.
Definition: mptrac.c:2129
void write_grid_asc(const char *filename, const ctl_t *ctl, const double *cd, double *mean[NQ], double *sigma[NQ], const double *vmr_impl, const double t, const double *z, const double *lon, const double *lat, const double *area, const double dz, const int *np)
Writes grid data to an ASCII file.
Definition: mptrac.c:12287
void mptrac_update_device(const ctl_t *ctl, const cache_t *cache, const clim_t *clim, met_t **met0, met_t **met1, const atm_t *atm)
Updates device memory for specified data structures.
Definition: mptrac.c:6728
void time2jsec(const int year, const int mon, const int day, const int hour, const int min, const int sec, const double remain, double *jsec)
Converts time components to seconds since January 1, 2000, 12:00:00 UTC.
Definition: mptrac.c:11175
void mptrac_init(ctl_t *ctl, cache_t *cache, clim_t *clim, atm_t *atm, const int ntask)
Initializes the MPTRAC model and its associated components.
Definition: mptrac.c:5411
void intpol_met_time_3d(const met_t *met0, float array0[EX][EY][EP], const met_t *met1, float array1[EX][EY][EP], const double ts, const double p, const double lon, const double lat, double *var, int *ci, double *cw, const int init)
Interpolates meteorological data in 3D space and time.
Definition: mptrac.c:2417
#define codes_handle
Placeholder when ECCODES is not available.
Definition: mptrac.h:190
void fft_help(double *fcReal, double *fcImag, const int n)
Computes the Fast Fourier Transform (FFT) of a complex sequence.
Definition: mptrac.c:1911
void module_wet_depo(const ctl_t *ctl, const cache_t *cache, met_t *met0, met_t *met1, atm_t *atm)
Perform wet deposition calculations for air parcels.
Definition: mptrac.c:5068
void compress_pck(const char *varname, float *array, const size_t nxy, const size_t nz, const int decompress, FILE *inout)
Compresses or decompresses a 3D array of floats.
Definition: mptrac.c:680
double nat_temperature(const double p, const double h2o, const double hno3)
Calculates the nitric acid trihydrate (NAT) temperature.
Definition: mptrac.c:7035
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: mptrac.c:11107
#define CP
Maximum number of pressure levels for climatological data.
Definition: mptrac.h:348
#define NQ
Maximum number of quantities per data point.
Definition: mptrac.h:313
double clim_photo(const double rate[CP][CSZA][CO3], const clim_photo_t *photo, const double p, const double sza, const double o3c)
Calculates the photolysis rate for a given set of atmospheric conditions.
Definition: mptrac.c:147
void read_clim_zm(const char *filename, const char *varname, clim_zm_t *zm)
Reads zonally averaged climatological data from a netCDF file and populates the given structure.
Definition: mptrac.c:7445
#define EX
Maximum number of longitudes for meteo data.
Definition: mptrac.h:288
void read_met_nc_grid_dd_naive(dd_t *dd, const ctl_t *ctl, met_t *met, const int ncid)
Read meteorological grid data from a NetCDF file and set up subdomain decomposition with halos.
Definition: mptrac.c:9935
void module_sedi(const ctl_t *ctl, const cache_t *cache, met_t *met0, met_t *met1, atm_t *atm)
Simulate sedimentation of particles in the atmosphere.
Definition: mptrac.c:4787
int read_met_nc(const char *filename, const ctl_t *ctl, met_t *met, dd_t *dd)
Reads meteorological data from a NetCDF file and processes it.
Definition: mptrac.c:9895
void timer(const char *name, const char *group, const int output)
Measures and reports elapsed time for named and grouped timers.
Definition: mptrac.c:11206
void write_atm_asc(const char *filename, const ctl_t *ctl, const atm_t *atm, const double t)
Writes air parcel data to an ASCII file or gnuplot.
Definition: mptrac.c:11332
void intpol_met_space_3d(const met_t *met, float array[EX][EY][EP], const double p, const double lon, const double lat, double *var, int *ci, double *cw, const int init)
Interpolates meteorological variables in 3D space.
Definition: mptrac.c:2301
void module_convection(const ctl_t *ctl, cache_t *cache, met_t *met0, met_t *met1, atm_t *atm)
Performs convective mixing of atmospheric particles.
Definition: mptrac.c:3450
void read_kernel(const char *filename, double kz[EP], double kw[EP], int *nk)
Reads kernel function data from a file and populates the provided arrays.
Definition: mptrac.c:7544
void module_bound_cond(const ctl_t *ctl, const cache_t *cache, const clim_t *clim, met_t *met0, met_t *met1, atm_t *atm)
Apply boundary conditions to particles based on meteorological and climatological data.
Definition: mptrac.c:3150
double scan_ctl(const char *filename, int argc, char *argv[], const char *varname, const int arridx, const char *defvalue, char *value)
Scans a control file or command-line arguments for a specified variable.
Definition: mptrac.c:11002
#define CT
Maximum number of time steps for climatological data.
Definition: mptrac.h:358
void module_advect_init(const ctl_t *ctl, const cache_t *cache, met_t *met0, met_t *met1, atm_t *atm)
Initializes the advection module by setting up pressure fields.
Definition: mptrac.c:3123
int dd_is_periodic_longitude(met_t *met, int nx_glob)
Check whether the longitude grid is periodic (global coverage).
void module_radio_decay(const ctl_t *ctl, const cache_t *cache, atm_t *atm)
Apply radioactive decay to atmospheric tracer species.
Definition: mptrac.c:4593
int read_met_grib(const char *filename, const ctl_t *ctl, met_t *met)
Reads meteorological data from a grib file and processes it.
void mptrac_get_met(ctl_t *ctl, clim_t *clim, const double t, met_t **met0, met_t **met1, dd_t *dd)
Retrieves meteorological data for the specified time.
Definition: mptrac.c:5289
void module_sort_help(double *a, const int *p, const int np)
Reorder an array based on a given permutation.
Definition: mptrac.c:4881
float stddev(const float *data, const int n)
Calculates the standard deviation of a set of data.
Definition: mptrac.c:11154
void intpol_tropo_3d(const double time0, float array0[EX][EY], const double time1, float array1[EX][EY], const double lons[EX], const double lats[EY], const int nlon, const int nlat, const double time, const double lon, const double lat, const int method, double *var, double *sigma)
Interpolates tropopause data in 3D (latitude, longitude, and time).
Definition: mptrac.c:2479
void read_met_bin_3d(FILE *in, const ctl_t *ctl, const met_t *met, float var[EX][EY][EP], const char *varname, const float bound_min, const float bound_max)
Reads 3D meteorological data from a binary file, potentially using different compression methods.
Definition: mptrac.c:7762
int locate_irr_float(const float *xx, const int n, const double x, const int ig)
Locate the index of the interval containing a given value in an irregularly spaced array.
Definition: mptrac.c:2887
void dd_get_rect_neighbour(const ctl_t ctl, dd_t *dd)
Determines rectangular neighbouring ranks for MPI processes.
double time_from_filename(const char *filename, const int offset)
Extracts and converts a timestamp from a filename to Julian seconds.
Definition: mptrac.c:11274
void write_prof(const char *filename, const ctl_t *ctl, met_t *met0, met_t *met1, const atm_t *atm, const double t)
Writes profile data to a specified file.
Definition: mptrac.c:12975
void dd_atm2particles(atm_t *atm, particle_t *particles, ctl_t *ctl, int *nparticles, cache_t *cache, int rank)
Extracts particles from an atmospheric state and prepares them for inter-domain transfer.
void mptrac_read_clim(const ctl_t *ctl, clim_t *clim)
Reads various climatological data and populates the given climatology structure.
Definition: mptrac.c:5501
void module_chem_grid(const ctl_t *ctl, met_t *met0, met_t *met1, atm_t *atm, const double t)
Computes gridded chemical tracer concentrations (volume mixing ratio) from individual air parcel mass...
Definition: mptrac.c:3246
double tropo_weight(const clim_t *clim, const atm_t *atm, const int ip)
Computes a weighting factor based on tropopause pressure.
Definition: mptrac.c:11309
void write_met_nc(const char *filename, const ctl_t *ctl, met_t *met)
Writes meteorological data to a NetCDF file.
Definition: mptrac.c:12752
void module_rng_init(const int ntask)
Initialize random number generators for parallel tasks.
Definition: mptrac.c:4651
int mptrac_read_atm(const char *filename, const ctl_t *ctl, atm_t *atm)
Reads air parcel data from a specified file into the given atmospheric structure.
Definition: mptrac.c:5430
void mptrac_update_host(const ctl_t *ctl, const cache_t *cache, const clim_t *clim, met_t **met0, met_t **met1, const atm_t *atm)
Updates host memory for specified data structures.
Definition: mptrac.c:6784
void mptrac_write_output(const char *dirname, const ctl_t *ctl, met_t *met0, met_t *met1, atm_t *atm, const double t)
Writes various types of output data to files in a specified directory.
Definition: mptrac.c:6944
double clim_oh(const ctl_t *ctl, const clim_t *clim, const double t, const double lon, const double lat, const double p)
Calculates the hydroxyl radical (OH) concentration from climatology data, with an optional diurnal co...
Definition: mptrac.c:89
void read_met_grib_levels(codes_handle **handles, const int num_messages, const ctl_t *ctl, met_t *met)
Reads meteorological variables at different vertical levels from a grib file.
#define CTS
Maximum number of data points of climatological time series.
Definition: mptrac.h:363
void read_met_ozone(met_t *met)
Calculates the total column ozone from meteorological ozone data.
Definition: mptrac.c:10527
int dd_calc_subdomain_from_coords(double lon, double lat, met_t *met, ctl_t *ctl, int mpi_size, int nx_glob, int ny_glob)
Computes the destination subdomain (MPI rank) for a particle based on its geographic coordinates.
void broadcast_large_data(void *data, size_t N)
Broadcasts large data across all processes in an MPI communicator.
void read_met_nc_surface(const int ncid, const ctl_t *ctl, met_t *met, dd_t *dd)
Reads and processes surface meteorological data from NetCDF files with domain decomposition.
Definition: mptrac.c:8479
void clim_tropo_init(clim_t *clim)
Initializes the tropopause data in the climatology structure.
Definition: mptrac.c:232
void module_rng(const ctl_t *ctl, double *rs, const size_t n, const int method)
Generate random numbers using various methods and distributions.
Definition: mptrac.c:4682
void get_met_help(const ctl_t *ctl, const double t, const int direct, const char *metbase, const double dt_met, char *filename)
Generates a formatted filename for meteorological data files based on the input parameters.
Definition: mptrac.c:1968
void write_station(const char *filename, const ctl_t *ctl, atm_t *atm, const double t)
Writes station data to a specified file.
Definition: mptrac.c:13364
void cart2geo(const double *x, double *z, double *lon, double *lat)
Converts Cartesian coordinates to geographic coordinates.
Definition: mptrac.c:74
#define NP
Maximum number of atmospheric data points.
Definition: mptrac.h:308
void doy2day(const int year, const int doy, int *mon, int *day)
Converts a given day of the year (DOY) to a date (month and day).
Definition: mptrac.c:1881
void dd_communicate_particles(particle_t *particles, int *nparticles, MPI_Datatype MPI_Particle, int *neighbours, int nneighbours, ctl_t ctl)
Communicates particles between MPI processes.
void intpol_met_time_2d(const met_t *met0, float array0[EX][EY], const met_t *met1, float array1[EX][EY], const double ts, const double lon, const double lat, double *var, int *ci, double *cw, const int init)
Interpolates meteorological data in 2D space and time.
Definition: mptrac.c:2446
void module_position(const cache_t *cache, met_t *met0, met_t *met1, atm_t *atm)
Update the positions and pressure levels of atmospheric particles.
Definition: mptrac.c:4542
void clim_oh_diurnal_correction(const ctl_t *ctl, clim_t *clim)
Applies a diurnal correction to the hydroxyl radical (OH) concentration in climatology data.
Definition: mptrac.c:115
void write_met_bin_2d(FILE *out, met_t *met, float var[EX][EY], const char *varname)
Writes a 2-dimensional meteorological variable to a binary file.
Definition: mptrac.c:12633
void read_met_pv(met_t *met)
Calculates potential vorticity (PV) from meteorological data.
Definition: mptrac.c:10421
int read_atm_bin(const char *filename, const ctl_t *ctl, atm_t *atm)
Reads air parcel data from a binary file and populates the given atmospheric structure.
Definition: mptrac.c:7123
void get_met_replace(char *orig, char *search, char *repl)
Replaces occurrences of a substring in a string with another substring.
Definition: mptrac.c:2035
#define CY
Maximum number of latitudes for climatological data.
Definition: mptrac.h:338
void module_diff_meso(const ctl_t *ctl, cache_t *cache, met_t *met0, met_t *met1, atm_t *atm)
Simulate mesoscale diffusion for atmospheric particles.
Definition: mptrac.c:3602
void module_dd(ctl_t *ctl, atm_t *atm, cache_t *cache, dd_t *dd, met_t **met)
Manages domain decomposition and particle communication in parallel processing.
void dd_sort_help(double *a, dd_t *dd, const int np)
Reorder an array according to a permutation vector.
void module_sort(const ctl_t *ctl, met_t *met0, atm_t *atm)
Sort particles according to box index.
Definition: mptrac.c:4816
double clim_ts(const clim_ts_t *ts, const double t)
Interpolates a time series of climatological variables.
Definition: mptrac.c:387
void jsec2time(const double jsec, int *year, int *mon, int *day, int *hour, int *min, int *sec, double *remain)
Converts Julian seconds to calendar date and time components.
Definition: mptrac.c:2570
int read_met_bin(const char *filename, const ctl_t *ctl, met_t *met)
Reads meteorological data from a binary file.
Definition: mptrac.c:7585
void write_atm_clams(const char *filename, const ctl_t *ctl, const atm_t *atm)
Writes air parcel data to a NetCDF file in the CLaMS format.
Definition: mptrac.c:11464
void module_diff_turb(const ctl_t *ctl, cache_t *cache, const clim_t *clim, met_t *met0, met_t *met1, atm_t *atm)
Applies turbulent diffusion processes to atmospheric particles.
Definition: mptrac.c:3804
int read_atm_clams(const char *filename, const ctl_t *ctl, atm_t *atm)
Reads atmospheric data from a CLAMS NetCDF file.
Definition: mptrac.c:7179
void dd_assign_rect_subdomains_atm(atm_t *atm, ctl_t *ctl, dd_t *dd, int init)
Assign atmospheric particles to rectangular subdomains.
int mptrac_read_met(const char *filename, const ctl_t *ctl, const clim_t *clim, met_t *met, dd_t *dd)
Reads meteorological data from a file, supporting multiple formats and MPI broadcasting.
Definition: mptrac.c:6467
void write_vtk(const char *filename, const ctl_t *ctl, const atm_t *atm, const double t)
Writes VTK (Visualization Toolkit) data to a specified file.
Definition: mptrac.c:13450
#define MPI_Datatype
Placeholder when MPI is not available.
Definition: mptrac.h:145
#define EP
Maximum number of pressure levels for meteo data.
Definition: mptrac.h:283
void module_tracer_chem(const ctl_t *ctl, const cache_t *cache, const clim_t *clim, met_t *met0, met_t *met1, atm_t *atm)
Simulate chemical reactions involving long-lived atmospheric tracers.
Definition: mptrac.c:4997
void mptrac_read_ctl(const char *filename, int argc, char *argv[], ctl_t *ctl)
Reads control parameters from a configuration file and populates the given structure.
Definition: mptrac.c:5561
void mptrac_free(ctl_t *ctl, cache_t *cache, clim_t *clim, met_t *met0, met_t *met1, atm_t *atm, dd_t *dd)
Frees memory resources allocated for MPTRAC.
Definition: mptrac.c:5256
void read_met_polar_winds(met_t *met)
Applies a fix for polar winds in meteorological data.
Definition: mptrac.c:10362
void module_h2o2_chem(const ctl_t *ctl, const cache_t *cache, const clim_t *clim, met_t *met0, met_t *met1, atm_t *atm)
Perform chemical reactions involving H2O2 within cloud particles.
Definition: mptrac.c:3919
void write_grid_nc(const char *filename, const ctl_t *ctl, const double *cd, double *mean[NQ], double *sigma[NQ], const double *vmr_impl, const double t, const double *z, const double *lon, const double *lat, const double *area, const double dz, const int *np)
Writes grid data to a NetCDF file.
Definition: mptrac.c:12391
double pbl_weight(const ctl_t *ctl, const atm_t *atm, const int ip, const double pbl, const double ps)
Computes a weighting factor based on planetary boundary layer pressure.
Definition: mptrac.c:7059
void module_diff_pbl(const ctl_t *ctl, cache_t *cache, met_t *met0, met_t *met1, atm_t *atm)
Computes particle diffusion within the planetary boundary layer (PBL).
Definition: mptrac.c:3679
void write_met_nc_3d(const int ncid, const char *varname, met_t *met, float var[EX][EY][EP], const float scl)
Writes a 3D meteorological variable to a NetCDF file.
Definition: mptrac.c:12945
void module_isosurf(const ctl_t *ctl, const cache_t *cache, met_t *met0, met_t *met1, atm_t *atm)
Apply the isosurface module to adjust atmospheric properties.
Definition: mptrac.c:4071
void module_kpp_chem(ctl_t *ctl, cache_t *cache, clim_t *clim, met_t *met0, met_t *met1, atm_t *atm)
KPP chemistry module.
#define CO3
Maximum number of total column ozone data for climatological data.
Definition: mptrac.h:343
void module_oh_chem(const ctl_t *ctl, const cache_t *cache, const clim_t *clim, met_t *met0, met_t *met1, atm_t *atm)
Perform hydroxyl chemistry calculations for atmospheric particles.
Definition: mptrac.c:4458
void read_met_grib_grid(codes_handle **handles, int count_handles, met_t *met)
Reads global meteorological information from a grib file.
void geo2cart(const double z, const double lon, const double lat, double *x)
Converts geographic coordinates (longitude, latitude, altitude) to Cartesian coordinates.
Definition: mptrac.c:1950
void read_met_nc_grid(const char *filename, const int ncid, const ctl_t *ctl, met_t *met, dd_t *dd)
Reads meteorological grid data from NetCDF files with domain decomposition.
Definition: mptrac.c:8312
double kernel_weight(const double kz[EP], const double kw[EP], const int nk, const double p)
Calculates the kernel weight based on altitude and given kernel data.
Definition: mptrac.c:2603
void compress_sz3(const char *varname, float *array, const int nx, const int ny, const int nz, const int precision, const double tolerance, const int decompress, FILE *inout)
Compresses or decompresses a 3-D float array using the SZ3 library.
void dd_register_MPI_type_particle(MPI_Datatype *MPI_Particle)
Registers a custom MPI datatype for particle structures.
int read_atm_asc(const char *filename, const ctl_t *ctl, atm_t *atm)
Reads air parcel data from an ASCII file and populates the given atmospheric structure.
Definition: mptrac.c:7081
void intpol_check_lon_lat(const double *lons, const int nlon, const double *lats, const int nlat, const double lon, const double lat, double *lon2, double *lat2)
Adjusts longitude and latitude to ensure they fall within valid bounds.
Definition: mptrac.c:2102
void write_sample(const char *filename, const ctl_t *ctl, met_t *met0, met_t *met1, const atm_t *atm, const double t)
Writes sample data to a specified file.
Definition: mptrac.c:13202
void write_grid(const char *filename, const ctl_t *ctl, met_t *met0, met_t *met1, const atm_t *atm, const double t)
Writes grid data to a file in ASCII or netCDF format.
Definition: mptrac.c:12096
void module_dry_depo(const ctl_t *ctl, const cache_t *cache, met_t *met0, met_t *met1, atm_t *atm)
Simulate dry deposition of atmospheric particles.
Definition: mptrac.c:3856
void write_met_bin(const char *filename, const ctl_t *ctl, met_t *met)
Writes meteorological data in binary format to a specified file.
Definition: mptrac.c:12521
void write_atm_bin(const char *filename, const ctl_t *ctl, const atm_t *atm)
Writes air parcel data to a binary file.
Definition: mptrac.c:11414
void read_met_cape(const ctl_t *ctl, const clim_t *clim, met_t *met)
Calculates Convective Available Potential Energy (CAPE) for each grid point.
Definition: mptrac.c:7868
void mptrac_write_met(const char *filename, const ctl_t *ctl, met_t *met)
Writes meteorological data to a file, supporting multiple formats and compression options.
Definition: mptrac.c:6900
#define CSZA
Maximum number of solar zenith angles for climatological data.
Definition: mptrac.h:353
double lapse_rate(const double t, const double h2o)
Calculates the moist adiabatic lapse rate in Kelvin per kilometer.
Definition: mptrac.c:2629
#define DD_NNMAX
Maximum number of neighbours to communicate with in domain decomposition.
Definition: mptrac.h:373
void write_csi(const char *filename, const ctl_t *ctl, const atm_t *atm, const double t)
Writes Critical Success Index (CSI) data to a file.
Definition: mptrac.c:11724
void write_atm_nc(const char *filename, const ctl_t *ctl, const atm_t *atm)
Writes air parcel data to a NetCDF file.
Definition: mptrac.c:11675
Air parcel data.
Definition: mptrac.h:3185
int np
Number of air parcels.
Definition: mptrac.h:3188
Cache data structure.
Definition: mptrac.h:3240
int iso_n
Isosurface balloon number of data points.
Definition: mptrac.h:3252
Climatological data in the form of photolysis rates.
Definition: mptrac.h:3272
int nsza
Number of solar zenith angles.
Definition: mptrac.h:3278
int np
Number of pressure levels.
Definition: mptrac.h:3275
int no3c
Number of total ozone columns.
Definition: mptrac.h:3281
Climatological data.
Definition: mptrac.h:3380
clim_ts_t ccl2f2
CFC-12 time series.
Definition: mptrac.h:3422
clim_photo_t photo
Photolysis rates.
Definition: mptrac.h:3398
clim_zm_t ho2
HO2 zonal means.
Definition: mptrac.h:3410
clim_zm_t hno3
HNO3 zonal means.
Definition: mptrac.h:3401
int tropo_ntime
Number of tropopause timesteps.
Definition: mptrac.h:3383
clim_ts_t sf6
SF6 time series.
Definition: mptrac.h:3428
clim_ts_t ccl4
CFC-10 time series.
Definition: mptrac.h:3416
clim_ts_t ccl3f
CFC-11 time series.
Definition: mptrac.h:3419
clim_zm_t o1d
O(1D) zonal means.
Definition: mptrac.h:3413
clim_zm_t h2o2
H2O2 zonal means.
Definition: mptrac.h:3407
int tropo_nlat
Number of tropopause latitudes.
Definition: mptrac.h:3386
clim_zm_t oh
OH zonal means.
Definition: mptrac.h:3404
clim_ts_t n2o
N2O time series.
Definition: mptrac.h:3425
Climatological data in the form of time series.
Definition: mptrac.h:3328
int ntime
Number of timesteps.
Definition: mptrac.h:3331
Climatological data in the form of zonal means.
Definition: mptrac.h:3348
int np
Number of pressure levels.
Definition: mptrac.h:3357
int ntime
Number of timesteps.
Definition: mptrac.h:3351
int nlat
Number of latitudes.
Definition: mptrac.h:3354
Control parameters.
Definition: mptrac.h:2131
double grid_z0
Lower altitude of gridded data [km].
Definition: mptrac.h:3046
int qnt_o3
Quantity array index for ozone volume mixing ratio.
Definition: mptrac.h:2243
double csi_lat1
Upper latitude of gridded CSI data [deg].
Definition: mptrac.h:3007
int qnt_Coh
Quantity array index for OH volume mixing ratio (chemistry code).
Definition: mptrac.h:2399
double wet_depo_ic_a
Coefficient A for wet deposition in cloud (exponential form).
Definition: mptrac.h:2895
int met_nc_scale
Check netCDF scaling factors (0=no, 1=yes).
Definition: mptrac.h:2495
int qnt_pel
Quantity array index for pressure at equilibrium level (EL).
Definition: mptrac.h:2276
int csi_nz
Number of altitudes of gridded CSI data.
Definition: mptrac.h:2983
double molmass
Molar mass [g/mol].
Definition: mptrac.h:2754
int qnt_p
Quantity array index for pressure.
Definition: mptrac.h:2222
int qnt_Cccl2f2
Quantity array index for CFC-12 volume mixing ratio (chemistry code).
Definition: mptrac.h:2423
int dd_halos_size
Domain decomposition size of halos given in grid-points.
Definition: mptrac.h:3173
int mixing_nx
Number of longitudes of mixing grid.
Definition: mptrac.h:2817
double chemgrid_z1
Upper altitude of chemistry grid [km].
Definition: mptrac.h:2841
int qnt_m
Quantity array index for mass.
Definition: mptrac.h:2162
int qnt_aoa
Quantity array index for age of air.
Definition: mptrac.h:2432
int qnt_rhop
Quantity array index for particle density.
Definition: mptrac.h:2171
int qnt_swc
Quantity array index for cloud snow water content.
Definition: mptrac.h:2255
double csi_obsmin
Minimum observation index to trigger detection.
Definition: mptrac.h:2977
int qnt_pcb
Quantity array index for cloud bottom pressure.
Definition: mptrac.h:2264
double bound_dzs
Boundary conditions surface layer depth [km].
Definition: mptrac.h:2742
double csi_lon1
Upper longitude of gridded CSI data [deg].
Definition: mptrac.h:2998
int qnt_u
Quantity array index for zonal wind.
Definition: mptrac.h:2231
double stat_lon
Longitude of station [deg].
Definition: mptrac.h:3124
double mixing_trop
Interparcel exchange parameter for mixing in the troposphere.
Definition: mptrac.h:2802
double sort_dt
Time step for sorting of particle data [s].
Definition: mptrac.h:2653
double mixing_z1
Upper altitude of mixing grid [km].
Definition: mptrac.h:2814
double stat_r
Search radius around station [km].
Definition: mptrac.h:3130
double wet_depo_bc_a
Coefficient A for wet deposition below cloud (exponential form).
Definition: mptrac.h:2889
int met_zstd_level
ZSTD compression level (from -5 to 22).
Definition: mptrac.h:2504
int csi_ny
Number of latitudes of gridded CSI data.
Definition: mptrac.h:3001
int vtk_sphere
Spherical projection for VTK data (0=no, 1=yes).
Definition: mptrac.h:3154
double chemgrid_z0
Lower altitude of chemistry grid [km].
Definition: mptrac.h:2838
double met_pbl_min
Minimum depth of planetary boundary layer [km].
Definition: mptrac.h:2621
int qnt_iwc
Quantity array index for cloud ice water content.
Definition: mptrac.h:2252
double chemgrid_lat0
Lower latitude of chemistry grid [deg].
Definition: mptrac.h:2856
double conv_cape
CAPE threshold for convection module [J/kg].
Definition: mptrac.h:2706
int qnt_Co1d
Quantity array index for O(1D) volume mixing ratio (chemistry code).
Definition: mptrac.h:2411
double met_cms_eps_pv
cmultiscale compression epsilon for potential vorticity.
Definition: mptrac.h:2543
int qnt_pw
Quantity array index for partial water vapor pressure.
Definition: mptrac.h:2330
double grid_z1
Upper altitude of gridded data [km].
Definition: mptrac.h:3049
int direction
Direction flag (1=forward calculation, -1=backward calculation).
Definition: mptrac.h:2459
int qnt_Cccl4
Quantity array index for CFC-10 volume mixing ratio (chemistry code).
Definition: mptrac.h:2417
int met_dp
Stride for pressure levels.
Definition: mptrac.h:2573
double met_dt_out
Time step for sampling of meteo data along trajectories [s].
Definition: mptrac.h:2640
int qnt_h2o2
Quantity array index for H2O2 volume mixing ratio (climatology).
Definition: mptrac.h:2294
int qnt_vh
Quantity array index for horizontal wind.
Definition: mptrac.h:2366
int csi_nx
Number of longitudes of gridded CSI data.
Definition: mptrac.h:2992
double csi_lat0
Lower latitude of gridded CSI data [deg].
Definition: mptrac.h:3004
double turb_dz_trop
Vertical turbulent diffusion coefficient (troposphere) [m^2/s].
Definition: mptrac.h:2688
int met_pbl
Planetary boundary layer data (0=file, 1=z2p, 2=Richardson, 3=theta).
Definition: mptrac.h:2618
int qnt_lwc
Quantity array index for cloud liquid water content.
Definition: mptrac.h:2246
double turb_mesoz
Vertical scaling factor for mesoscale wind fluctuations.
Definition: mptrac.h:2697
int grid_nc_level
zlib compression level of netCDF grid data files (0=off).
Definition: mptrac.h:3034
int grid_nx
Number of longitudes of gridded data.
Definition: mptrac.h:3052
int atm_type
Type of atmospheric data files (0=ASCII, 1=binary, 2=netCDF, 3=CLaMS_traj, 4=CLaMS_pos).
Definition: mptrac.h:2948
double bound_mass
Boundary conditions mass per particle [kg].
Definition: mptrac.h:2715
double grid_lat0
Lower latitude of gridded data [deg].
Definition: mptrac.h:3064
int qnt_ts
Quantity array index for surface temperature.
Definition: mptrac.h:2177
int qnt_loss_rate
Quantity array index for total loss rate.
Definition: mptrac.h:2321
double met_cms_eps_h2o
cmultiscale compression epsilon for water vapor.
Definition: mptrac.h:2546
int qnt_plfc
Quantity array index for pressure at level of free convection (LCF).
Definition: mptrac.h:2273
int qnt_Acs137
Quantity array index for radioactive activity of Cs-137.
Definition: mptrac.h:2444
double grid_lon0
Lower longitude of gridded data [deg].
Definition: mptrac.h:3055
int qnt_o1d
Quantity array index for O(1D) volume mixing ratio (climatology).
Definition: mptrac.h:2300
int met_tropo_spline
Tropopause interpolation method (0=linear, 1=spline).
Definition: mptrac.h:2637
int qnt_tvirt
Quantity array index for virtual temperature.
Definition: mptrac.h:2360
double dt_met
Time step of meteo data [s].
Definition: mptrac.h:2478
double chemgrid_lat1
Upper latitude of chemistry grid [deg].
Definition: mptrac.h:2859
int met_geopot_sy
Latitudinal smoothing of geopotential heights.
Definition: mptrac.h:2609
double met_cms_eps_u
cmultiscale compression epsilon for zonal wind.
Definition: mptrac.h:2534
double turb_dx_strat
Horizontal turbulent diffusion coefficient (stratosphere) [m^2/s].
Definition: mptrac.h:2682
int qnt_vmr
Quantity array index for volume mixing ratio.
Definition: mptrac.h:2165
int qnt_lsm
Quantity array index for land-sea mask.
Definition: mptrac.h:2198
int qnt_theta
Quantity array index for potential temperature.
Definition: mptrac.h:2342
double bound_lat1
Boundary conditions maximum longitude [deg].
Definition: mptrac.h:2730
double stat_t1
Stop time for station output [s].
Definition: mptrac.h:3136
double turb_dx_trop
Horizontal turbulent diffusion coefficient (troposphere) [m^2/s].
Definition: mptrac.h:2679
int grid_type
Type of grid data files (0=ASCII, 1=netCDF).
Definition: mptrac.h:3070
double csi_lon0
Lower longitude of gridded CSI data [deg].
Definition: mptrac.h:2995
int qnt_pbl
Quantity array index for boundary layer pressure.
Definition: mptrac.h:2204
int grid_stddev
Include standard deviations in grid output (0=no, 1=yes).
Definition: mptrac.h:3040
int qnt_psice
Quantity array index for saturation pressure over ice.
Definition: mptrac.h:2327
double chemgrid_lon0
Lower longitude of chemistry grid [deg].
Definition: mptrac.h:2847
int bound_pbl
Boundary conditions planetary boundary layer (0=no, 1=yes).
Definition: mptrac.h:2748
int qnt_mloss_wet
Quantity array index for total mass loss due to wet deposition.
Definition: mptrac.h:2312
int radio_decay
Switch for radioactive decay module (0=off, 1=on).
Definition: mptrac.h:2883
int met_geopot_sx
Longitudinal smoothing of geopotential heights.
Definition: mptrac.h:2606
int met_sy
Smoothing for latitudes.
Definition: mptrac.h:2579
int qnt_ps
Quantity array index for surface pressure.
Definition: mptrac.h:2174
int rng_type
Random number generator (0=GSL, 1=Squares, 2=cuRAND).
Definition: mptrac.h:2670
int isosurf
Isosurface parameter (0=none, 1=pressure, 2=density, 3=theta, 4=balloon).
Definition: mptrac.h:2657
double bound_p1
Boundary conditions top pressure [hPa].
Definition: mptrac.h:2736
int qnt_zs
Quantity array index for surface geopotential height.
Definition: mptrac.h:2180
int prof_nz
Number of altitudes of gridded profile data.
Definition: mptrac.h:3079
double csi_dt_out
Time step for CSI output [s].
Definition: mptrac.h:2971
int met_cape
Convective available potential energy data (0=file, 1=calculate).
Definition: mptrac.h:2615
double csi_modmin
Minimum column density to trigger detection [kg/m^2].
Definition: mptrac.h:2980
int met_sx
Smoothing for longitudes.
Definition: mptrac.h:2576
double chemgrid_lon1
Upper longitude of chemistry grid [deg].
Definition: mptrac.h:2850
double turb_mesox
Horizontal scaling factor for mesoscale wind fluctuations.
Definition: mptrac.h:2694
double met_cms_eps_iwc
cmultiscale compression epsilon for cloud ice water content.
Definition: mptrac.h:2558
double met_cms_eps_swc
cmultiscale compression epsilon for cloud snow water content.
Definition: mptrac.h:2561
double met_cms_eps_v
cmultiscale compression epsilon for meridional wind.
Definition: mptrac.h:2537
double prof_z0
Lower altitude of gridded profile data [km].
Definition: mptrac.h:3082
int qnt_w
Quantity array index for vertical velocity.
Definition: mptrac.h:2237
double bound_vmr
Boundary conditions volume mixing ratio [ppv].
Definition: mptrac.h:2721
double met_tropo_pv
Dynamical tropopause potential vorticity threshold [PVU].
Definition: mptrac.h:2631
int prof_nx
Number of longitudes of gridded profile data.
Definition: mptrac.h:3088
int qnt_stat
Quantity array index for station flag.
Definition: mptrac.h:2159
int met_tropo
Tropopause definition (0=none, 1=clim, 2=cold point, 3=WMO_1st, 4=WMO_2nd, 5=dynamical).
Definition: mptrac.h:2628
int qnt_rp
Quantity array index for particle radius.
Definition: mptrac.h:2168
int met_mpi_share
Use MPI to share meteo (0=no, 1=yes).
Definition: mptrac.h:2646
double mixing_strat
Interparcel exchange parameter for mixing in the stratosphere.
Definition: mptrac.h:2805
int qnt_vz
Quantity array index for vertical velocity.
Definition: mptrac.h:2369
int qnt_ho2
Quantity array index for HO2 volume mixing ratio (climatology).
Definition: mptrac.h:2297
double csi_z1
Upper altitude of gridded CSI data [km].
Definition: mptrac.h:2989
double stat_t0
Start time for station output [s].
Definition: mptrac.h:3133
double oh_chem_beta
Beta parameter for diurnal variablity of OH.
Definition: mptrac.h:2868
int dd
Domain decomposition (0=no, 1=yes, with 2x2 if not specified).
Definition: mptrac.h:3161
int qnt_eta
Quantity array index for eta vertical coordinate.
Definition: mptrac.h:2354
double wet_depo_so2_ph
pH value used to calculate effective Henry constant of SO2.
Definition: mptrac.h:2907
double mixing_z0
Lower altitude of mixing grid [km].
Definition: mptrac.h:2811
int qnt_mloss_decay
Quantity array index for total mass loss due to exponential decay.
Definition: mptrac.h:2318
int atm_type_out
Type of atmospheric data files for output (-1=same as ATM_TYPE, 0=ASCII, 1=binary,...
Definition: mptrac.h:2953
int met_cms_nd0x
cmultiscale number of cells of coarsest grid in x-direction.
Definition: mptrac.h:2519
int met_nlev
Number of meteo data model levels.
Definition: mptrac.h:2597
double dt_kpp
Time step for KPP chemistry [s].
Definition: mptrac.h:2877
double dry_depo_dp
Dry deposition surface layer [hPa].
Definition: mptrac.h:2916
int qnt_shf
Quantity array index for surface sensible heat flux.
Definition: mptrac.h:2195
int qnt_vs
Quantity array index for surface meridional wind.
Definition: mptrac.h:2186
int qnt_Cco
Quantity array index for CO volume mixing ratio (chemistry code).
Definition: mptrac.h:2396
double vtk_dt_out
Time step for VTK data output [s].
Definition: mptrac.h:3142
double t_stop
Stop time of simulation [s].
Definition: mptrac.h:2465
double conv_dt
Time interval for convection module [s].
Definition: mptrac.h:2712
int qnt_hno3
Quantity array index for HNO3 volume mixing ratio (climatology).
Definition: mptrac.h:2288
int met_clams
Read MPTRAC or CLaMS meteo data (0=MPTRAC, 1=CLaMS).
Definition: mptrac.h:2492
int qnt_h2ot
Quantity array index for tropopause water vapor volume mixing ratio.
Definition: mptrac.h:2216
int qnt_rh
Quantity array index for relative humidity over water.
Definition: mptrac.h:2336
double met_cms_eps_cc
cmultiscale compression epsilon for cloud cover.
Definition: mptrac.h:2564
double bound_lat0
Boundary conditions minimum longitude [deg].
Definition: mptrac.h:2727
double met_pbl_max
Maximum depth of planetary boundary layer [km].
Definition: mptrac.h:2624
int met_dx
Stride for longitudes.
Definition: mptrac.h:2567
int qnt_destination
Quantity array index for destination subdomain in domain decomposition.
Definition: mptrac.h:2456
int mixing_ny
Number of latitudes of mixing grid.
Definition: mptrac.h:2826
int met_convention
Meteo data layout (0=[lev, lat, lon], 1=[lon, lat, lev]).
Definition: mptrac.h:2481
int qnt_zeta_d
Quantity array index for diagnosed zeta vertical coordinate.
Definition: mptrac.h:2348
int tracer_chem
Switch for first order tracer chemistry module (0=off, 1=on).
Definition: mptrac.h:2880
double dt_mod
Time step of simulation [s].
Definition: mptrac.h:2468
int diffusion
Diffusion scheme (0=off, 1=fixed-K, 2=PBL).
Definition: mptrac.h:2673
int qnt_tnat
Quantity array index for T_NAT.
Definition: mptrac.h:2384
int qnt_eta_dot
Quantity array index for velocity of eta vertical coordinate.
Definition: mptrac.h:2357
int qnt_tice
Quantity array index for T_ice.
Definition: mptrac.h:2378
int qnt_zg
Quantity array index for geopotential height.
Definition: mptrac.h:2219
double vtk_offset
Vertical offset for VTK data [km].
Definition: mptrac.h:3151
int qnt_v
Quantity array index for meridional wind.
Definition: mptrac.h:2234
int qnt_mloss_dry
Quantity array index for total mass loss due to dry deposition.
Definition: mptrac.h:2315
double bound_vmr_trend
Boundary conditions volume mixing ratio trend [ppv/s].
Definition: mptrac.h:2724
int met_cache
Preload meteo data into disk cache (0=no, 1=yes).
Definition: mptrac.h:2643
int qnt_oh
Quantity array index for OH volume mixing ratio (climatology).
Definition: mptrac.h:2291
int qnt_Ch
Quantity array index for H volume mixing ratio (chemistry code).
Definition: mptrac.h:2402
int met_press_level_def
Use predefined pressure levels or not.
Definition: mptrac.h:2594
int oh_chem_reaction
Reaction type for OH chemistry (0=none, 2=bimolecular, 3=termolecular).
Definition: mptrac.h:2862
int qnt_h2o
Quantity array index for water vapor volume mixing ratio.
Definition: mptrac.h:2240
int prof_ny
Number of latitudes of gridded profile data.
Definition: mptrac.h:3097
int qnt_rhice
Quantity array index for relative humidity over ice.
Definition: mptrac.h:2339
int qnt_rho
Quantity array index for density of air.
Definition: mptrac.h:2228
double sample_dz
Layer depth for sample output [km].
Definition: mptrac.h:3118
double tdec_strat
Life time of particles in the stratosphere [s].
Definition: mptrac.h:2760
int obs_type
Type of observation data files (0=ASCII, 1=netCDF).
Definition: mptrac.h:2962
double met_cms_eps_lwc
cmultiscale compression epsilon for cloud liquid water content.
Definition: mptrac.h:2552
int qnt_us
Quantity array index for surface zonal wind.
Definition: mptrac.h:2183
double met_cms_eps_z
cmultiscale compression epsilon for geopotential height.
Definition: mptrac.h:2528
double grid_lon1
Upper longitude of gridded data [deg].
Definition: mptrac.h:3058
int qnt_Cn2o
Quantity array index for N2O volume mixing ratio (chemistry code).
Definition: mptrac.h:2426
int qnt_Cccl3f
Quantity array index for CFC-11 volume mixing ratio (chemistry code).
Definition: mptrac.h:2420
double mixing_lat0
Lower latitude of mixing grid [deg].
Definition: mptrac.h:2829
int nens
Number of ensembles.
Definition: mptrac.h:3010
int qnt_pt
Quantity array index for tropopause pressure.
Definition: mptrac.h:2207
int qnt_cl
Quantity array index for total column cloud water.
Definition: mptrac.h:2267
int advect
Advection scheme (0=off, 1=Euler, 2=midpoint, 4=Runge-Kutta).
Definition: mptrac.h:2663
double prof_z1
Upper altitude of gridded profile data [km].
Definition: mptrac.h:3085
int qnt_t
Quantity array index for temperature.
Definition: mptrac.h:2225
int atm_filter
Time filter for atmospheric data output (0=none, 1=missval, 2=remove).
Definition: mptrac.h:2941
int kpp_chem
Switch for KPP chemistry module (0=off, 1=on).
Definition: mptrac.h:2874
int qnt_zeta
Quantity array index for zeta vertical coordinate.
Definition: mptrac.h:2345
double conv_pbl_trans
Depth of PBL transition layer (fraction of PBL depth).
Definition: mptrac.h:2703
int qnt_Ai131
Quantity array index for radioactive activity of I-131.
Definition: mptrac.h:2447
int met_vert_coord
Vertical coordinate of input meteo data (0=plev, 1=mlev_p_file, 2=mlev_ab_file, 3=mlev_ab_full,...
Definition: mptrac.h:2485
double csi_z0
Lower altitude of gridded CSI data [km].
Definition: mptrac.h:2986
int qnt_lapse
Quantity array index for lapse rate.
Definition: mptrac.h:2363
int qnt_Apb210
Quantity array index for radioactive activity of Pb-210.
Definition: mptrac.h:2438
double stat_lat
Latitude of station [deg].
Definition: mptrac.h:3127
int qnt_Cho2
Quantity array index for HO2 volume mixing ratio (chemistry code).
Definition: mptrac.h:2405
int grid_ny
Number of latitudes of gridded data.
Definition: mptrac.h:3061
int qnt_Csf6
Quantity array index for SF6 volume mixing ratio (chemistry code).
Definition: mptrac.h:2429
int qnt_Ch2o
Quantity array index for H2O volume mixing ratio (chemistry code).
Definition: mptrac.h:2390
double met_detrend
FWHM of horizontal Gaussian used for detrending [km].
Definition: mptrac.h:2585
int conv_mix_pbl
Vertical mixing in the PBL (0=off, 1=on).
Definition: mptrac.h:2700
double bound_dps
Boundary conditions surface layer depth [hPa].
Definition: mptrac.h:2739
int dd_nbr_neighbours
Domain decomposition number of neighbours to communicate with.
Definition: mptrac.h:3170
double met_cms_eps_t
cmultiscale compression epsilon for temperature.
Definition: mptrac.h:2531
int chemgrid_nz
Number of altitudes of chemistry grid.
Definition: mptrac.h:2835
int qnt_cape
Quantity array index for convective available potential energy (CAPE).
Definition: mptrac.h:2279
int qnt_zeta_dot
Quantity array index for velocity of zeta vertical coordinate.
Definition: mptrac.h:2351
double bound_mass_trend
Boundary conditions mass per particle trend [kg/s].
Definition: mptrac.h:2718
int met_cms_nd0y
cmultiscale number of cells of coarsest grid in y-direction.
Definition: mptrac.h:2522
int mixing_nz
Number of altitudes of mixing grid.
Definition: mptrac.h:2808
int qnt_o3c
Quantity array index for total column ozone.
Definition: mptrac.h:2285
double bound_p0
Boundary conditions bottom pressure [hPa].
Definition: mptrac.h:2733
double mixing_lon0
Lower longitude of mixing grid [deg].
Definition: mptrac.h:2820
int qnt_Co3
Quantity array index for O3 volume mixing ratio (chemistry code).
Definition: mptrac.h:2393
int qnt_tsts
Quantity array index for T_STS.
Definition: mptrac.h:2381
int grid_nz
Number of altitudes of gridded data.
Definition: mptrac.h:3043
int qnt_nss
Quantity array index for northward turbulent surface stress.
Definition: mptrac.h:2192
double ens_dt_out
Time step for ensemble output [s].
Definition: mptrac.h:3016
int atm_stride
Particle index stride for atmospheric data files.
Definition: mptrac.h:2944
int met_relhum
Try to read relative humidity (0=no, 1=yes).
Definition: mptrac.h:2612
double mixing_lat1
Upper latitude of mixing grid [deg].
Definition: mptrac.h:2832
double atm_dt_out
Time step for atmospheric data output [s].
Definition: mptrac.h:2938
double prof_lat1
Upper latitude of gridded profile data [deg].
Definition: mptrac.h:3103
int met_cms_batch
cmultiscale batch size.
Definition: mptrac.h:2513
double psc_h2o
H2O volume mixing ratio for PSC analysis.
Definition: mptrac.h:2922
int met_sp
Smoothing for pressure levels.
Definition: mptrac.h:2582
double prof_lon0
Lower longitude of gridded profile data [deg].
Definition: mptrac.h:3091
int qnt_Axe133
Quantity array index for radioactive activity of Xe-133.
Definition: mptrac.h:2450
int chemgrid_nx
Number of longitudes of chemistry grid.
Definition: mptrac.h:2844
int qnt_pct
Quantity array index for cloud top pressure.
Definition: mptrac.h:2261
int qnt_mloss_kpp
Quantity array index for total mass loss due to KPP chemistry.
Definition: mptrac.h:2309
int qnt_psat
Quantity array index for saturation pressure over water.
Definition: mptrac.h:2324
int qnt_subdomain
Quantity array index for current subdomain in domain decomposition.
Definition: mptrac.h:2453
double prof_lat0
Lower latitude of gridded profile data [deg].
Definition: mptrac.h:3100
int qnt_cin
Quantity array index for convective inhibition (CIN).
Definition: mptrac.h:2282
double psc_hno3
HNO3 volume mixing ratio for PSC analysis.
Definition: mptrac.h:2925
double prof_lon1
Upper longitude of gridded profile data [deg].
Definition: mptrac.h:3094
double met_cms_eps_rwc
cmultiscale compression epsilon for cloud rain water content.
Definition: mptrac.h:2555
int met_nc_quant
Number of digits for quantization of netCDF meteo files (0=off).
Definition: mptrac.h:2501
int h2o2_chem_reaction
Reaction type for H2O2 chemistry (0=none, 1=SO2).
Definition: mptrac.h:2871
int qnt_Co3p
Quantity array index for O(3P) volume mixing ratio (chemistry code).
Definition: mptrac.h:2414
double wet_depo_bc_ret_ratio
Coefficients for wet deposition below cloud: retention ratio.
Definition: mptrac.h:2913
int chemgrid_ny
Number of latitudes of chemistry grid.
Definition: mptrac.h:2853
int qnt_Abe7
Quantity array index for radioactive activity of Be-7.
Definition: mptrac.h:2441
double met_cms_eps_o3
cmultiscale compression epsilon for ozone.
Definition: mptrac.h:2549
int met_cms_zstd
cmultiscale ZSTD compression (0=off, 1=on).
Definition: mptrac.h:2516
int met_cms_maxlev
cmultiscale maximum refinement level.
Definition: mptrac.h:2525
int grid_sparse
Sparse output in grid data files (0=no, 1=yes).
Definition: mptrac.h:3031
double dry_depo_vdep
Dry deposition velocity [m/s].
Definition: mptrac.h:2919
int qnt_tt
Quantity array index for tropopause temperature.
Definition: mptrac.h:2210
int met_np
Number of target pressure levels.
Definition: mptrac.h:2588
int qnt_ens
Quantity array index for ensemble IDs.
Definition: mptrac.h:2156
int met_nc_level
zlib compression level of netCDF meteo files (0=off).
Definition: mptrac.h:2498
double mixing_dt
Time interval for mixing [s].
Definition: mptrac.h:2799
int qnt_Arn222
Quantity array index for radioactive activity of Rn-222.
Definition: mptrac.h:2435
int qnt_mloss_h2o2
Quantity array index for total mass loss due to H2O2 chemistry.
Definition: mptrac.h:2306
double vtk_scale
Vertical scaling factor for VTK data.
Definition: mptrac.h:3148
double met_cms_eps_w
cmultiscale compression epsilon for vertical velocity.
Definition: mptrac.h:2540
double turb_dx_pbl
Horizontal turbulent diffusion coefficient (PBL) [m^2/s].
Definition: mptrac.h:2676
double conv_cin
CIN threshold for convection module [J/kg].
Definition: mptrac.h:2709
int qnt_pv
Quantity array index for potential vorticity.
Definition: mptrac.h:2372
int advect_vert_coord
Vertical velocity of air parcels (0=omega_on_plev, 1=zetadot_on_mlev, 2=omega_on_mlev,...
Definition: mptrac.h:2667
int qnt_mloss_oh
Quantity array index for total mass loss due to OH chemistry.
Definition: mptrac.h:2303
int qnt_Ch2o2
Quantity array index for H2O2 volume mixing ratio (chemistry code).
Definition: mptrac.h:2408
int qnt_sst
Quantity array index for sea surface temperature.
Definition: mptrac.h:2201
double mixing_lon1
Upper longitude of mixing grid [deg].
Definition: mptrac.h:2823
int atm_nc_level
zlib compression level of netCDF atmospheric data files (0=off).
Definition: mptrac.h:2956
double wet_depo_ic_ret_ratio
Coefficients for wet deposition in cloud: retention ratio.
Definition: mptrac.h:2910
int qnt_sh
Quantity array index for specific humidity.
Definition: mptrac.h:2333
int qnt_ess
Quantity array index for eastward turbulent surface stress.
Definition: mptrac.h:2189
double wet_depo_ic_b
Coefficient B for wet deposition in cloud (exponential form).
Definition: mptrac.h:2898
double wet_depo_bc_b
Coefficient B for wet deposition below cloud (exponential form).
Definition: mptrac.h:2892
int met_dy
Stride for latitudes.
Definition: mptrac.h:2570
int qnt_Cx
Quantity array index for trace species x volume mixing ratio (chemistry code).
Definition: mptrac.h:2387
double turb_dz_strat
Vertical turbulent diffusion coefficient (stratosphere) [m^2/s].
Definition: mptrac.h:2691
double bound_zetas
Boundary conditions surface layer zeta [K].
Definition: mptrac.h:2745
int dd_subdomains_zonal
Domain decomposition zonal subdomain number.
Definition: mptrac.h:3164
int qnt_idx
Quantity array index for air parcel IDs.
Definition: mptrac.h:2153
double met_tropo_theta
Dynamical tropopause potential temperature threshold [K].
Definition: mptrac.h:2634
int qnt_rwc
Quantity array index for cloud rain water content.
Definition: mptrac.h:2249
double t_start
Start time of simulation [s].
Definition: mptrac.h:2462
int nq
Number of quantities.
Definition: mptrac.h:2138
double tdec_trop
Life time of particles in the troposphere [s].
Definition: mptrac.h:2757
double sample_dx
Horizontal radius for sample output [km].
Definition: mptrac.h:3115
int vtk_stride
Particle index stride for VTK data.
Definition: mptrac.h:3145
double turb_dz_pbl
Vertical turbulent diffusion coefficient (PBL) [m^2/s].
Definition: mptrac.h:2685
double grid_lat1
Upper latitude of gridded data [deg].
Definition: mptrac.h:3067
int dd_subdomains_meridional
Domain decomposition meridional subdomain number.
Definition: mptrac.h:3167
int qnt_zt
Quantity array index for tropopause geopotential height.
Definition: mptrac.h:2213
int met_type
Type of meteo data files (0=netCDF, 1=binary, 2=pck, 3=ZFP, 4=ZSTD, 5=cms, 6=grib,...
Definition: mptrac.h:2489
int qnt_cc
Quantity array index for cloud cover.
Definition: mptrac.h:2258
int qnt_plcl
Quantity array index for pressure at lifted condensation level (LCL).
Definition: mptrac.h:2270
double grid_dt_out
Time step for gridded data output [s].
Definition: mptrac.h:3028
int qnt_tdew
Quantity array index for dew point temperature.
Definition: mptrac.h:2375
Domain decomposition data structure.
Definition: mptrac.h:3613
int halo_offset_end
Hyperslab of boundary halos count.
Definition: mptrac.h:3665
int rank
Rank of node.
Definition: mptrac.h:3620
double subdomain_lon_min
Rectangular grid limit of subdomain.
Definition: mptrac.h:3641
double subdomain_lat_max
Rectangular grid limit of subdomain.
Definition: mptrac.h:3644
int init
Shows if domain decomposition was initialized.
Definition: mptrac.h:3672
double subdomain_lon_max
Rectangular grid limit of subdomain.
Definition: mptrac.h:3638
int halo_offset_start
Hyperslab of boundary halos count.
Definition: mptrac.h:3662
int size
Size of node.
Definition: mptrac.h:3623
double subdomain_lat_min
Rectangular grid limit of subdomain.
Definition: mptrac.h:3647
Meteo data structure.
Definition: mptrac.h:3439
int nx
Number of longitudes.
Definition: mptrac.h:3445
int ny
Number of latitudes.
Definition: mptrac.h:3448
int np
Number of pressure levels.
Definition: mptrac.h:3451
int npl
Number of model levels.
Definition: mptrac.h:3454
double time
Time [s].
Definition: mptrac.h:3442
Particle data.
Definition: mptrac.h:3214
double p
Pressure [hPa].
Definition: mptrac.h:3220
double lat
Latitude [deg].
Definition: mptrac.h:3226
double time
Time [s].
Definition: mptrac.h:3217
double lon
Longitude [deg].
Definition: mptrac.h:3223