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
166#ifndef LIBTRAC_H
167#define LIBTRAC_H
168
169/* ------------------------------------------------------------
170 Includes...
171 ------------------------------------------------------------ */
172
173#include <ctype.h>
174#include <gsl/gsl_fft_complex.h>
175#include <gsl/gsl_math.h>
176#include <gsl/gsl_randist.h>
177#include <gsl/gsl_rng.h>
178#include <gsl/gsl_sort.h>
179#include <gsl/gsl_spline.h>
180#include <gsl/gsl_statistics.h>
181#include <math.h>
182#include <netcdf.h>
183#include <omp.h>
184#include <stdint.h>
185#include <stdio.h>
186#include <stdlib.h>
187#include <string.h>
188#include <time.h>
189#include <sys/time.h>
190
191#ifdef MPI
192#include "mpi.h"
193#endif
194
195#ifdef DD
196#include <netcdf_par.h>
197#endif
198
199#ifdef _OPENACC
200#include "openacc.h"
201#endif
202
203#ifdef CURAND
204#include "curand.h"
205#endif
206
207#ifdef THRUST
208#include "thrustsort.h"
209#endif
210
211#ifdef ZFP
212#include "zfp.h"
213#endif
214
215#ifdef ZSTD
216#include "zstd.h"
217#endif
218
219#ifdef SZ3
220#include "SZ3c/sz3c.h"
221#endif
222
223#ifdef CMS
224#include "cmultiscale.h"
225#endif
226
227#ifdef KPP
228#include "chem_Parameters.h"
229#include "chem_Global.h"
230#include "chem_Sparse.h"
231#endif
232
233#ifdef ECCODES
234#include "eccodes.h"
235#else
237#define codes_handle void*
238#endif
239
240/* ------------------------------------------------------------
241 Constants...
242 ------------------------------------------------------------ */
243
245#ifndef AVO
246#define AVO 6.02214076e23
247#endif
248
250#ifndef CPD
251#define CPD 1003.5
252#endif
253
255#ifndef EPS
256#define EPS (MH2O / MA)
257#endif
258
260#ifndef G0
261#define G0 9.80665
262#endif
263
265#ifndef H0
266#define H0 7.0
267#endif
268
270#ifndef LV
271#define LV 2501000.
272#endif
273
275#ifndef KARMAN
276#define KARMAN 0.40
277#endif
278
280#ifndef KB
281#define KB 1.3806504e-23
282#endif
283
285#ifndef MA
286#define MA 28.9644
287#endif
288
290#ifndef MH2O
291#define MH2O 18.01528
292#endif
293
295#ifndef MO3
296#define MO3 48.00
297#endif
298
300#ifndef P0
301#define P0 1013.25
302#endif
303
305#ifndef RA
306#define RA (1e3 * RI / MA)
307#endif
308
310#ifndef RE
311#define RE 6367.421
312#endif
313
315#ifndef RI
316#define RI 8.3144598
317#endif
318
320#ifndef T0
321#define T0 273.15
322#endif
323
324/* ------------------------------------------------------------
325 Dimensions...
326 ------------------------------------------------------------ */
327
329#ifndef EP
330#define EP 140
331#endif
332
334#ifndef EX
335#define EX 1444
336#endif
337
339#ifndef EY
340#define EY 724
341#endif
342
344#ifndef LEN
345#define LEN 5000
346#endif
347
349#ifndef METVAR
350#define METVAR 13
351#endif
352
354#ifndef NP
355#define NP 10000000
356#endif
357
359#ifndef NQ
360#define NQ 15
361#endif
362
364#ifndef NCSI
365#define NCSI 1000000
366#endif
367
369#ifndef NENS
370#define NENS 2000
371#endif
372
374#ifndef NOBS
375#define NOBS 10000000
376#endif
377
379#ifndef NTHREADS
380#define NTHREADS 512
381#endif
382
384#ifndef CY
385#define CY 250
386#endif
387
389#ifndef CO3
390#define CO3 30
391#endif
392
394#ifndef CP
395#define CP 70
396#endif
397
399#ifndef CSZA
400#define CSZA 50
401#endif
402
404#ifndef CT
405#define CT 12
406#endif
407
409#ifndef CTS
410#define CTS 1000
411#endif
412
414#ifndef DD_EX_GLOB
415#define DD_EX_GLOB (EX * 16)
416#endif
417
419#ifndef DD_EY_GLOB
420#define DD_EY_GLOB (EY * 16)
421#endif
422
423/* ------------------------------------------------------------
424 Macros...
425 ------------------------------------------------------------ */
426
446#ifdef _OPENACC
447#define ALLOC(ptr, type, n) \
448 if(acc_get_num_devices(acc_device_nvidia) <= 0) \
449 ERRMSG("Not running on a GPU device!"); \
450 if((ptr=calloc((size_t)(n), sizeof(type)))==NULL) \
451 ERRMSG("Out of memory!");
452#else
453#define ALLOC(ptr, type, n) \
454 if((ptr=calloc((size_t)(n), sizeof(type)))==NULL) \
455 ERRMSG("Out of memory!");
456#endif
457
476#define ARRAY_2D(ix, iy, ny) \
477 ((ix) * (ny) + (iy))
478
495#define ARRAY_3D(ix, iy, ny, iz, nz) \
496 (((ix)*(ny) + (iy)) * (nz) + (iz))
497
520#define ARRHENIUS(a, b, t) \
521 ((a) * exp( -(b) / (t)))
522
542#define CLAMP(v, lo, hi) \
543 (((v) < (lo)) ? (lo) : (((v) > (hi)) ? (hi) : (v)))
544
566#define DEG2DX(dlon, lat) \
567 (RE * DEG2RAD(dlon) * cos(DEG2RAD(lat)))
568
587#define DEG2DY(dlat) \
588 (RE * DEG2RAD(dlat))
589
604#define DEG2RAD(deg) \
605 ((deg) * (M_PI / 180.0))
606
629#define DP2DZ(dp, p) \
630 (- (dp) * H0 / (p))
631
651#define DX2DEG(dx, lat) \
652 (((lat) < -89.999 || (lat) > 89.999) ? 0 \
653 : (dx) * 180. / (M_PI * RE * cos(DEG2RAD(lat))))
654
669#define DY2DEG(dy) \
670 ((dy) * 180. / (M_PI * RE))
671
688#define DZ2DP(dz, p) \
689 (-(dz) * (p) / H0)
690
704#define DIST(a, b) \
705 sqrt(DIST2(a, b))
706
720#define DIST2(a, b) \
721 ((a[0]-b[0])*(a[0]-b[0])+(a[1]-b[1])*(a[1]-b[1])+(a[2]-b[2])*(a[2]-b[2]))
722
736#define DOTP(a, b) \
737 (a[0]*b[0]+a[1]*b[1]+a[2]*b[2])
738
750#define ECC(cmd) { \
751 int ecc_result=(cmd); \
752 if(ecc_result!=0) \
753 ERRMSG("ECCODES error: %s", codes_get_error_message(ecc_result)); \
754 }
755
769#define ECC_READ_2D(variable, target, scaling_factor, found_flag) { \
770 if(strcmp(short_name, variable) == 0) { \
771 for (int ix = 0; ix < met->nx; ix++) \
772 for (int iy = 0; iy < met->ny; iy++) \
773 target[ix][iy] = (float)(values[iy * met->nx + ix] * scaling_factor); \
774 found_flag = 1; \
775 } \
776 }
777
792#define ECC_READ_3D(variable, level, target, scaling_factor, found_flag) { \
793 if(strcmp(short_name, variable) == 0) { \
794 for (int ix = 0; ix < met->nx; ix++) \
795 for (int iy = 0; iy < met->ny; iy++) \
796 target[ix][iy][level] = (float) (values[iy * met->nx + ix] * scaling_factor); \
797 found_flag += 1; \
798 } \
799 }
800
817#define FMOD(x, y) \
818 ((x) - (int) ((x) / (y)) * (y))
819
835#define FREAD(ptr, type, size, in) { \
836 if(fread(ptr, sizeof(type), size, in)!=size) \
837 ERRMSG("Error while reading!"); \
838 }
839
855#define FWRITE(ptr, type, size, out) { \
856 if(fwrite(ptr, sizeof(type), size, out)!=size) \
857 ERRMSG("Error while writing!"); \
858 }
859
870#define INTPOL_INIT \
871 double cw[4] = {0.0, 0.0, 0.0, 0.0}; int ci[3] = {0, 0, 0};
872
884#define INTPOL_2D(var, init) \
885 intpol_met_time_2d(met0, met0->var, met1, met1->var, \
886 atm->time[ip], atm->lon[ip], atm->lat[ip], \
887 &var, ci, cw, init);
888
901#define INTPOL_3D(var, init) \
902 intpol_met_time_3d(met0, met0->var, met1, met1->var, \
903 atm->time[ip], atm->p[ip], \
904 atm->lon[ip], atm->lat[ip], \
905 &var, ci, cw, init);
906
920#define INTPOL_SPACE_ALL(p, lon, lat) { \
921 intpol_met_space_3d(met, met->z, p, lon, lat, &z, ci, cw, 1); \
922 intpol_met_space_3d(met, met->t, p, lon, lat, &t, ci, cw, 0); \
923 intpol_met_space_3d(met, met->u, p, lon, lat, &u, ci, cw, 0); \
924 intpol_met_space_3d(met, met->v, p, lon, lat, &v, ci, cw, 0); \
925 intpol_met_space_3d(met, met->w, p, lon, lat, &w, ci, cw, 0); \
926 intpol_met_space_3d(met, met->pv, p, lon, lat, &pv, ci, cw, 0); \
927 intpol_met_space_3d(met, met->h2o, p, lon, lat, &h2o, ci, cw, 0); \
928 intpol_met_space_3d(met, met->o3, p, lon, lat, &o3, ci, cw, 0); \
929 intpol_met_space_3d(met, met->lwc, p, lon, lat, &lwc, ci, cw, 0); \
930 intpol_met_space_3d(met, met->rwc, p, lon, lat, &rwc, ci, cw, 0); \
931 intpol_met_space_3d(met, met->iwc, p, lon, lat, &iwc, ci, cw, 0); \
932 intpol_met_space_3d(met, met->swc, p, lon, lat, &swc, ci, cw, 0); \
933 intpol_met_space_3d(met, met->cc, p, lon, lat, &cc, ci, cw, 0); \
934 intpol_met_space_2d(met, met->ps, lon, lat, &ps, ci, cw, 0); \
935 intpol_met_space_2d(met, met->ts, lon, lat, &ts, ci, cw, 0); \
936 intpol_met_space_2d(met, met->zs, lon, lat, &zs, ci, cw, 0); \
937 intpol_met_space_2d(met, met->us, lon, lat, &us, ci, cw, 0); \
938 intpol_met_space_2d(met, met->vs, lon, lat, &vs, ci, cw, 0); \
939 intpol_met_space_2d(met, met->ess, ess, lat, &ess, ci, cw, 0); \
940 intpol_met_space_2d(met, met->nss, nss, lat, &nss, ci, cw, 0); \
941 intpol_met_space_2d(met, met->shf, shf, lat, &shf, ci, cw, 0); \
942 intpol_met_space_2d(met, met->lsm, lon, lat, &lsm, ci, cw, 0); \
943 intpol_met_space_2d(met, met->sst, lon, lat, &sst, ci, cw, 0); \
944 intpol_met_space_2d(met, met->pbl, lon, lat, &pbl, ci, cw, 0); \
945 intpol_met_space_2d(met, met->pt, lon, lat, &pt, ci, cw, 0); \
946 intpol_met_space_2d(met, met->tt, lon, lat, &tt, ci, cw, 0); \
947 intpol_met_space_2d(met, met->zt, lon, lat, &zt, ci, cw, 0); \
948 intpol_met_space_2d(met, met->h2ot, lon, lat, &h2ot, ci, cw, 0); \
949 intpol_met_space_2d(met, met->pct, lon, lat, &pct, ci, cw, 0); \
950 intpol_met_space_2d(met, met->pcb, lon, lat, &pcb, ci, cw, 0); \
951 intpol_met_space_2d(met, met->cl, lon, lat, &cl, ci, cw, 0); \
952 intpol_met_space_2d(met, met->plcl, lon, lat, &plcl, ci, cw, 0); \
953 intpol_met_space_2d(met, met->plfc, lon, lat, &plfc, ci, cw, 0); \
954 intpol_met_space_2d(met, met->pel, lon, lat, &pel, ci, cw, 0); \
955 intpol_met_space_2d(met, met->cape, lon, lat, &cape, ci, cw, 0); \
956 intpol_met_space_2d(met, met->cin, lon, lat, &cin, ci, cw, 0); \
957 intpol_met_space_2d(met, met->o3c, lon, lat, &o3c, ci, cw, 0); \
958 }
959
974#define INTPOL_TIME_ALL(time, p, lon, lat) { \
975 intpol_met_time_3d(met0, met0->z, met1, met1->z, time, p, lon, lat, &z, ci, cw, 1); \
976 intpol_met_time_3d(met0, met0->t, met1, met1->t, time, p, lon, lat, &t, ci, cw, 0); \
977 intpol_met_time_3d(met0, met0->u, met1, met1->u, time, p, lon, lat, &u, ci, cw, 0); \
978 intpol_met_time_3d(met0, met0->v, met1, met1->v, time, p, lon, lat, &v, ci, cw, 0); \
979 intpol_met_time_3d(met0, met0->w, met1, met1->w, time, p, lon, lat, &w, ci, cw, 0); \
980 intpol_met_time_3d(met0, met0->pv, met1, met1->pv, time, p, lon, lat, &pv, ci, cw, 0); \
981 intpol_met_time_3d(met0, met0->h2o, met1, met1->h2o, time, p, lon, lat, &h2o, ci, cw, 0); \
982 intpol_met_time_3d(met0, met0->o3, met1, met1->o3, time, p, lon, lat, &o3, ci, cw, 0); \
983 intpol_met_time_3d(met0, met0->lwc, met1, met1->lwc, time, p, lon, lat, &lwc, ci, cw, 0); \
984 intpol_met_time_3d(met0, met0->rwc, met1, met1->rwc, time, p, lon, lat, &rwc, ci, cw, 0); \
985 intpol_met_time_3d(met0, met0->iwc, met1, met1->iwc, time, p, lon, lat, &iwc, ci, cw, 0); \
986 intpol_met_time_3d(met0, met0->swc, met1, met1->swc, time, p, lon, lat, &swc, ci, cw, 0); \
987 intpol_met_time_3d(met0, met0->cc, met1, met1->cc, time, p, lon, lat, &cc, ci, cw, 0); \
988 intpol_met_time_2d(met0, met0->ps, met1, met1->ps, time, lon, lat, &ps, ci, cw, 0); \
989 intpol_met_time_2d(met0, met0->ts, met1, met1->ts, time, lon, lat, &ts, ci, cw, 0); \
990 intpol_met_time_2d(met0, met0->zs, met1, met1->zs, time, lon, lat, &zs, ci, cw, 0); \
991 intpol_met_time_2d(met0, met0->us, met1, met1->us, time, lon, lat, &us, ci, cw, 0); \
992 intpol_met_time_2d(met0, met0->vs, met1, met1->vs, time, lon, lat, &vs, ci, cw, 0); \
993 intpol_met_time_2d(met0, met0->ess, met1, met1->ess, time, lon, lat, &ess, ci, cw, 0); \
994 intpol_met_time_2d(met0, met0->nss, met1, met1->nss, time, lon, lat, &nss, ci, cw, 0); \
995 intpol_met_time_2d(met0, met0->shf, met1, met1->shf, time, lon, lat, &shf, ci, cw, 0); \
996 intpol_met_time_2d(met0, met0->lsm, met1, met1->lsm, time, lon, lat, &lsm, ci, cw, 0); \
997 intpol_met_time_2d(met0, met0->sst, met1, met1->sst, time, lon, lat, &sst, ci, cw, 0); \
998 intpol_met_time_2d(met0, met0->pbl, met1, met1->pbl, time, lon, lat, &pbl, ci, cw, 0); \
999 intpol_met_time_2d(met0, met0->pt, met1, met1->pt, time, lon, lat, &pt, ci, cw, 0); \
1000 intpol_met_time_2d(met0, met0->tt, met1, met1->tt, time, lon, lat, &tt, ci, cw, 0); \
1001 intpol_met_time_2d(met0, met0->zt, met1, met1->zt, time, lon, lat, &zt, ci, cw, 0); \
1002 intpol_met_time_2d(met0, met0->h2ot, met1, met1->h2ot, time, lon, lat, &h2ot, ci, cw, 0); \
1003 intpol_met_time_2d(met0, met0->pct, met1, met1->pct, time, lon, lat, &pct, ci, cw, 0); \
1004 intpol_met_time_2d(met0, met0->pcb, met1, met1->pcb, time, lon, lat, &pcb, ci, cw, 0); \
1005 intpol_met_time_2d(met0, met0->cl, met1, met1->cl, time, lon, lat, &cl, ci, cw, 0); \
1006 intpol_met_time_2d(met0, met0->plcl, met1, met1->plcl, time, lon, lat, &plcl, ci, cw, 0); \
1007 intpol_met_time_2d(met0, met0->plfc, met1, met1->plfc, time, lon, lat, &plfc, ci, cw, 0); \
1008 intpol_met_time_2d(met0, met0->pel, met1, met1->pel, time, lon, lat, &pel, ci, cw, 0); \
1009 intpol_met_time_2d(met0, met0->cape, met1, met1->cape, time, lon, lat, &cape, ci, cw, 0); \
1010 intpol_met_time_2d(met0, met0->cin, met1, met1->cin, time, lon, lat, &cin, ci, cw, 0); \
1011 intpol_met_time_2d(met0, met0->o3c, met1, met1->o3c, time, lon, lat, &o3c, ci, cw, 0); \
1012 }
1013
1028#define LAPSE(p1, t1, p2, t2) \
1029 (1e3 * G0 / RA * ((t2) - (t1)) / ((t2) + (t1)) \
1030 * ((p2) + (p1)) / ((p2) - (p1)))
1031
1047#define LIN(x0, y0, x1, y1, x) \
1048 ((y0)+((y1)-(y0))/((x1)-(x0))*((x)-(x0)))
1049
1074#define MAX(a,b) \
1075 (((a)>(b))?(a):(b))
1076
1088#define MET_HEADER \
1089 fprintf(out, \
1090 "# $1 = time [s]\n" \
1091 "# $2 = altitude [km]\n" \
1092 "# $3 = longitude [deg]\n" \
1093 "# $4 = latitude [deg]\n" \
1094 "# $5 = pressure [hPa]\n" \
1095 "# $6 = temperature [K]\n" \
1096 "# $7 = zonal wind [m/s]\n" \
1097 "# $8 = meridional wind [m/s]\n" \
1098 "# $9 = vertical velocity [hPa/s]\n" \
1099 "# $10 = H2O volume mixing ratio [ppv]\n"); \
1100 fprintf(out, \
1101 "# $11 = O3 volume mixing ratio [ppv]\n" \
1102 "# $12 = geopotential height [km]\n" \
1103 "# $13 = potential vorticity [PVU]\n" \
1104 "# $14 = surface pressure [hPa]\n" \
1105 "# $15 = surface temperature [K]\n" \
1106 "# $16 = surface geopotential height [km]\n" \
1107 "# $17 = surface zonal wind [m/s]\n" \
1108 "# $18 = surface meridional wind [m/s]\n" \
1109 "# $19 = eastward turbulent surface stress [N/m^2]\n" \
1110 "# $20 = northward turbulent surface stress [N/m^2]\n"); \
1111 fprintf(out, \
1112 "# $21 = surface sensible heat flux [W/m^2]\n" \
1113 "# $22 = land-sea mask [1]\n" \
1114 "# $23 = sea surface temperature [K]\n" \
1115 "# $24 = tropopause pressure [hPa]\n" \
1116 "# $25 = tropopause geopotential height [km]\n" \
1117 "# $26 = tropopause temperature [K]\n" \
1118 "# $27 = tropopause water vapor [ppv]\n" \
1119 "# $28 = cloud liquid water content [kg/kg]\n" \
1120 "# $29 = cloud rain water content [kg/kg]\n" \
1121 "# $30 = cloud ice water content [kg/kg]\n"); \
1122 fprintf(out, \
1123 "# $31 = cloud snow water content [kg/kg]\n" \
1124 "# $32 = cloud cover [1]\n" \
1125 "# $33 = total column cloud water [kg/m^2]\n" \
1126 "# $34 = cloud top pressure [hPa]\n" \
1127 "# $35 = cloud bottom pressure [hPa]\n" \
1128 "# $36 = pressure at lifted condensation level (LCL) [hPa]\n" \
1129 "# $37 = pressure at level of free convection (LFC) [hPa]\n" \
1130 "# $38 = pressure at equilibrium level (EL) [hPa]\n" \
1131 "# $39 = convective available potential energy (CAPE) [J/kg]\n" \
1132 "# $40 = convective inhibition (CIN) [J/kg]\n"); \
1133 fprintf(out, \
1134 "# $41 = relative humidity over water [%%]\n" \
1135 "# $42 = relative humidity over ice [%%]\n" \
1136 "# $43 = dew point temperature [K]\n" \
1137 "# $44 = frost point temperature [K]\n" \
1138 "# $45 = NAT temperature [K]\n" \
1139 "# $46 = HNO3 volume mixing ratio [ppv]\n" \
1140 "# $47 = OH volume mixing ratio [ppv]\n" \
1141 "# $48 = H2O2 volume mixing ratio [ppv]\n" \
1142 "# $49 = HO2 volume mixing ratio [ppv]\n" \
1143 "# $50 = O(1D) volume mixing ratio [ppv]\n"); \
1144 fprintf(out, \
1145 "# $51 = boundary layer pressure [hPa]\n" \
1146 "# $52 = total column ozone [DU]\n" \
1147 "# $53 = number of data points\n" \
1148 "# $54 = number of tropopause data points\n" \
1149 "# $55 = number of CAPE data points\n");
1150
1175#define MIN(a,b) \
1176 (((a)<(b))?(a):(b))
1177
1190#define MOLEC_DENS(p,t) \
1191 (AVO * 1e-6 * ((p) * 100) / (RI * (t)))
1192
1204#define NC(cmd) { \
1205 int nc_result=(cmd); \
1206 if(nc_result!=NC_NOERR) \
1207 ERRMSG("%s", nc_strerror(nc_result)); \
1208 }
1209
1233#define NC_DEF_VAR(varname, type, ndims, dims, long_name, units, level, quant) { \
1234 NC(nc_def_var(ncid, varname, type, ndims, dims, &varid)); \
1235 NC(nc_put_att_text(ncid, varid, "long_name", strnlen(long_name, LEN), long_name)); \
1236 NC(nc_put_att_text(ncid, varid, "units", strnlen(units, LEN), units)); \
1237 if((quant) > 0) \
1238 NC(nc_def_var_quantize(ncid, varid, NC_QUANTIZE_GRANULARBR, quant)); \
1239 if((level) != 0) { \
1240 NC(nc_def_var_deflate(ncid, varid, 1, 1, level)); \
1241 /* unsigned int ulevel = (unsigned int)level; */ \
1242 /* NC(nc_def_var_filter(ncid, varid, 32015, 1, (unsigned int[]){ulevel})); */ \
1243 } \
1244 }
1245
1263#define NC_GET_DOUBLE(varname, ptr, force) { \
1264 if(force) { \
1265 NC(nc_inq_varid(ncid, varname, &varid)); \
1266 NC(nc_get_var_double(ncid, varid, ptr)); \
1267 } else { \
1268 if(nc_inq_varid(ncid, varname, &varid) == NC_NOERR) { \
1269 NC(nc_get_var_double(ncid, varid, ptr)); \
1270 } else \
1271 WARN("netCDF variable %s is missing!", varname); \
1272 } \
1273 }
1274
1293#define NC_INQ_DIM(dimname, ptr, min, max, check) { \
1294 int dimid; size_t naux; \
1295 NC(nc_inq_dimid(ncid, dimname, &dimid)); \
1296 NC(nc_inq_dimlen(ncid, dimid, &naux)); \
1297 *ptr = (int)naux; \
1298 if (check) \
1299 if ((*ptr) < (min) || (*ptr) > (max)) \
1300 ERRMSG("Dimension %s is out of range!", dimname); \
1301 }
1302
1317#define NC_PUT_DOUBLE(varname, ptr, hyperslab) { \
1318 NC(nc_inq_varid(ncid, varname, &varid)); \
1319 if(hyperslab) { \
1320 NC(nc_put_vara_double(ncid, varid, start, count, ptr)); \
1321 } else { \
1322 NC(nc_put_var_double(ncid, varid, ptr)); \
1323 } \
1324 }
1325
1341#define NC_PUT_FLOAT(varname, ptr, hyperslab) { \
1342 NC(nc_inq_varid(ncid, varname, &varid)); \
1343 if(hyperslab) { \
1344 NC(nc_put_vara_float(ncid, varid, start, count, ptr)); \
1345 } else { \
1346 NC(nc_put_var_float(ncid, varid, ptr)); \
1347 } \
1348 }
1349
1364#define NC_PUT_INT(varname, ptr, hyperslab) { \
1365 NC(nc_inq_varid(ncid, varname, &varid)); \
1366 if(hyperslab) { \
1367 NC(nc_put_vara_int(ncid, varid, start, count, ptr)); \
1368 } else { \
1369 NC(nc_put_var_int(ncid, varid, ptr)); \
1370 } \
1371 }
1372
1386#define NC_PUT_ATT(varname, attname, text) { \
1387 NC(nc_inq_varid(ncid, varname, &varid)); \
1388 NC(nc_put_att_text(ncid, varid, attname, strnlen(text, LEN), text)); \
1389 }
1390
1403#define NC_PUT_ATT_GLOBAL(attname, text) \
1404 NC(nc_put_att_text(ncid, NC_GLOBAL, attname, strnlen(text, LEN), text));
1405
1423#define NN(x0, y0, x1, y1, x) \
1424 (fabs((x) - (x0)) <= fabs((x) - (x1)) ? (y0) : (y1))
1425
1441#ifdef _OPENACC
1442#define PARTICLE_LOOP(ip0, ip1, check_dt, ...) \
1443 const int ip0_const = ip0; \
1444 const int ip1_const = ip1; \
1445 _Pragma(__VA_ARGS__) \
1446 _Pragma("acc parallel loop independent gang vector") \
1447 for (int ip = ip0_const; ip < ip1_const; ip++) \
1448 if (!check_dt || cache->dt[ip] != 0)
1449#else
1450#define PARTICLE_LOOP(ip0, ip1, check_dt, ...) \
1451 const int ip0_const = ip0; \
1452 const int ip1_const = ip1; \
1453 _Pragma("omp parallel for default(shared)") \
1454 for (int ip = ip0_const; ip < ip1_const; ip++) \
1455 if (!check_dt || cache->dt[ip] != 0)
1456#endif
1457
1480#define P(z) \
1481 (P0 * exp(-(z) / H0))
1482
1504#define PSAT(t) \
1505 (6.112 * exp(17.62 * ((t) - T0) / (243.12 + (t) - T0)))
1506
1528#define PSICE(t) \
1529 (6.112 * exp(22.46 * ((t) - T0) / (272.62 + (t) - T0)))
1530
1555#define PW(p, h2o) \
1556 ((p) * MAX((h2o), 0.1e-6) / (1. + (1. - EPS) * MAX((h2o), 0.1e-6)))
1557
1572#define RAD2DEG(rad) \
1573 ((rad) * (180.0 / M_PI))
1574
1602#define RH(p, t, h2o) \
1603 (PW(p, h2o) / PSAT(t) * 100.)
1604
1632#define RHICE(p, t, h2o) \
1633 (PW(p, h2o) / PSICE(t) * 100.)
1634
1657#define RHO(p, t) \
1658 (100. * (p) / (RA * (t)))
1659
1676#define SET_ATM(qnt, val) \
1677 if (ctl->qnt >= 0) \
1678 atm->q[ctl->qnt][ip] = val;
1679
1699#define SET_QNT(qnt, name, longname, unit) \
1700 if (strcasecmp(ctl->qnt_name[iq], name) == 0) { \
1701 ctl->qnt = iq; \
1702 sprintf(ctl->qnt_longname[iq], longname); \
1703 sprintf(ctl->qnt_unit[iq], unit); \
1704 } else
1705
1720#define SH(h2o) \
1721 (EPS * MAX((h2o), 0.1e-6))
1722
1733#define SQR(x) \
1734 ((x)*(x))
1735
1747#define SWAP(x, y, type) \
1748 do {type tmp = x; x = y; y = tmp;} while(0);
1749
1771#define TDEW(p, h2o) \
1772 (T0 + 243.12 * log(PW((p), (h2o)) / 6.112) \
1773 / (17.62 - log(PW((p), (h2o)) / 6.112)))
1774
1796#define TICE(p, h2o) \
1797 (T0 + 272.62 * log(PW((p), (h2o)) / 6.112) \
1798 / (22.46 - log(PW((p), (h2o)) / 6.112)))
1799
1820#define THETA(p, t) \
1821 ((t) * pow(1000. / (p), 0.286))
1822
1849#define THETAVIRT(p, t, h2o) \
1850 (TVIRT(THETA((p), (t)), MAX((h2o), 0.1e-6)))
1851
1870#define TOK(line, tok, format, var) { \
1871 if(((tok)=strtok((line), " \t"))) { \
1872 if(sscanf(tok, format, &(var))!=1) continue; \
1873 } else ERRMSG("Error while reading!"); \
1874 }
1875
1895#define TVIRT(t, h2o) \
1896 ((t) * (1. + (1. - EPS) * MAX((h2o), 0.1e-6)))
1897
1909#define USAGE \
1910 do { \
1911 int iusage; \
1912 for (iusage = 1; iusage < argc; iusage++) \
1913 if (!strcmp(argv[iusage], "-h") \
1914 || !strcmp(argv[iusage], "--help")) { \
1915 usage(); \
1916 return EXIT_SUCCESS; \
1917 } \
1918 } while (0)
1919
1939#define Z(p) \
1940 (H0 * log(P0 / (p)))
1941
1970#define ZDIFF(lnp0, t0, h2o0, lnp1, t1, h2o1) \
1971 (RI / MA / G0 * 0.5 * (TVIRT((t0), (h2o0)) + TVIRT((t1), (h2o1))) \
1972 * ((lnp0) - (lnp1)))
1973
1989#define ZETA(ps, p, t) \
1990 (((p) / (ps) <= 0.3 ? 1. : \
1991 sin(M_PI / 2. * (1. - (p) / (ps)) / (1. - 0.3))) \
1992 * THETA((p), (t)))
1993
1994/* ------------------------------------------------------------
1995 Log messages...
1996 ------------------------------------------------------------ */
1997
1999#ifndef LOGLEV
2000#define LOGLEV 2
2001#endif
2002
2032#define LOG(level, ...) { \
2033 if(level >= 2) \
2034 printf(" "); \
2035 if(level <= LOGLEV) { \
2036 printf(__VA_ARGS__); \
2037 printf("\n"); \
2038 } \
2039 }
2040
2069#define WARN(...) { \
2070 printf("\nWarning (%s, %s, l%d): ", __FILE__, __func__, __LINE__); \
2071 LOG(0, __VA_ARGS__); \
2072 }
2073
2102#define ERRMSG(...) { \
2103 printf("\nError (%s, %s, l%d): ", __FILE__, __func__, __LINE__); \
2104 LOG(0, __VA_ARGS__); \
2105 exit(EXIT_FAILURE); \
2106 }
2107
2137#define PRINT(format, var) \
2138 printf("Print (%s, %s, l%d): %s= "format"\n", \
2139 __FILE__, __func__, __LINE__, #var, var);
2140
2141/* ------------------------------------------------------------
2142 Timers...
2143 ------------------------------------------------------------ */
2144
2146#define NTIMER 100
2147
2161#define PRINT_TIMERS \
2162 timer("END", "END", 1);
2163
2176#define SELECT_TIMER(id, group) \
2177 timer(id, group, 0);
2178
2179/* ------------------------------------------------------------
2180 Structs...
2181 ------------------------------------------------------------ */
2182
2190typedef struct {
2191
2192 /* ------------------------------------------------------------
2193 Quantity parameters...
2194 ------------------------------------------------------------ */
2195
2197 int nq;
2198
2200 char qnt_name[NQ][LEN];
2201
2203 char qnt_longname[NQ][LEN];
2204
2206 char qnt_unit[NQ][LEN];
2207
2209 char qnt_format[NQ][LEN];
2210
2213
2216
2219
2222
2225
2228
2231
2234
2237
2240
2243
2246
2249
2252
2255
2258
2261
2264
2267
2270
2273
2276
2279
2282
2285
2288
2291
2294
2297
2300
2303
2306
2309
2312
2315
2318
2321
2324
2327
2330
2333
2336
2339
2342
2345
2348
2351
2354
2357
2360
2363
2366
2369
2372
2375
2378
2381
2384
2387
2390
2393
2396
2399
2402
2405
2408
2411
2414
2417
2420
2423
2426
2429
2432
2435
2438
2441
2444
2447
2450
2453
2456
2459
2462
2465
2468
2471
2474
2477
2480
2483
2486
2489
2492
2495
2498
2501
2504
2507
2510
2513
2516
2519
2521 double t_start;
2522
2524 double t_stop;
2525
2527 double dt_mod;
2528
2529 /* ------------------------------------------------------------
2530 Meteo data parameters...
2531 ------------------------------------------------------------ */
2532
2534 char metbase[LEN];
2535
2537 double dt_met;
2538
2541
2545
2549
2552
2555
2558
2561
2564
2566 int met_comp_prec[METVAR];
2567
2569 double met_comp_tol[METVAR];
2570
2573
2576
2579
2582
2585
2588
2591
2594
2597
2600
2603
2606
2609
2612
2615
2618
2621
2624
2627
2630
2633
2636
2639
2642
2645
2648
2650 double met_p[EP];
2651
2654
2657
2659 double met_lev_hyam[EP];
2660
2662 double met_lev_hybm[EP];
2663
2666
2669
2672
2675
2678
2681
2684
2688
2691
2694
2697
2700
2703
2706
2707 /* ------------------------------------------------------------
2708 Geophysical module parameters...
2709 ------------------------------------------------------------ */
2710
2712 double sort_dt;
2713
2717
2719 char balloon[LEN];
2720
2723
2727
2730
2733
2736
2739
2742
2745
2748
2751
2754
2757
2760
2763
2766
2768 double conv_cin;
2769
2771 double conv_dt;
2772
2775
2778
2781
2784
2787
2790
2792 double bound_p0;
2793
2795 double bound_p1;
2796
2799
2802
2805
2808
2810 char species[LEN];
2811
2813 double molmass;
2814
2817
2820
2823
2825 char clim_hno3_filename[LEN];
2826
2828 char clim_oh_filename[LEN];
2829
2831 char clim_h2o2_filename[LEN];
2832
2834 char clim_ho2_filename[LEN];
2835
2837 char clim_o1d_filename[LEN];
2838
2840 char clim_o3_filename[LEN];
2841
2843 char clim_ccl4_timeseries[LEN];
2844
2846 char clim_ccl3f_timeseries[LEN];
2847
2849 char clim_ccl2f2_timeseries[LEN];
2850
2852 char clim_n2o_timeseries[LEN];
2853
2855 char clim_sf6_timeseries[LEN];
2856
2859
2862
2865
2868
2871
2874
2877
2880
2883
2886
2889
2892
2895
2898
2901
2904
2907
2910
2913
2916
2919
2922
2924 double oh_chem[4];
2925
2928
2931
2934
2936 double dt_kpp;
2937
2940
2943
2945 double wet_depo_pre[2];
2946
2949
2952
2955
2958
2960 double wet_depo_ic_h[2];
2961
2963 double wet_depo_bc_h[2];
2964
2967
2970
2973
2976
2979
2981 double psc_h2o;
2982
2984 double psc_hno3;
2985
2986 /* ------------------------------------------------------------
2987 Output parameters...
2988 ------------------------------------------------------------ */
2989
2991 char atm_basename[LEN];
2992
2994 char atm_gpfile[LEN];
2995
2998
3001
3004
3008
3013
3016
3018 int atm_nc_quant[NQ];
3019
3022
3024 char csi_basename[LEN];
3025
3027 char csi_kernel[LEN];
3028
3031
3033 char csi_obsfile[LEN];
3034
3037
3040
3043
3045 double csi_z0;
3046
3048 double csi_z1;
3049
3052
3054 double csi_lon0;
3055
3057 double csi_lon1;
3058
3061
3063 double csi_lat0;
3064
3066 double csi_lat1;
3067
3069 int nens;
3070
3072 char ens_basename[LEN];
3073
3076
3078 char grid_basename[LEN];
3079
3081 char grid_kernel[LEN];
3082
3084 char grid_gpfile[LEN];
3085
3088
3091
3094
3096 int grid_nc_quant[NQ];
3097
3100
3103
3105 double grid_z0;
3106
3108 double grid_z1;
3109
3112
3115
3118
3121
3124
3127
3130
3132 char prof_basename[LEN];
3133
3135 char prof_obsfile[LEN];
3136
3139
3141 double prof_z0;
3142
3144 double prof_z1;
3145
3148
3151
3154
3157
3160
3163
3165 char sample_basename[LEN];
3166
3168 char sample_kernel[LEN];
3169
3171 char sample_obsfile[LEN];
3172
3175
3178
3180 char stat_basename[LEN];
3181
3183 double stat_lon;
3184
3186 double stat_lat;
3187
3189 double stat_r;
3190
3192 double stat_t0;
3193
3195 double stat_t1;
3196
3198 char vtk_basename[LEN];
3199
3202
3205
3208
3211
3214
3215 /* ------------------------------------------------------------
3216 Domain decomposition...
3217 ------------------------------------------------------------ */
3218
3220 int dd;
3221
3224
3227
3230
3231} ctl_t;
3232
3241typedef struct {
3242
3244 int np;
3245
3247 double time[NP];
3248
3250 double p[NP];
3251
3253 double lon[NP];
3254
3256 double lat[NP];
3257
3259 double q[NQ][NP];
3260
3261} atm_t;
3262
3270typedef struct {
3271
3273 double time;
3274
3276 double p;
3277
3279 double lon;
3280
3282 double lat;
3283
3285 double q[NQ];
3286
3287} particle_t;
3288
3289
3296typedef struct {
3297
3299 double iso_var[NP];
3300
3302 double iso_ps[NP];
3303
3305 double iso_ts[NP];
3306
3309
3311 float uvwp[NP][3];
3312
3314 double rs[3 * NP + 1];
3315
3317 double dt[NP];
3318
3319} cache_t;
3320
3328typedef struct {
3329
3331 int np;
3332
3334 int nsza;
3335
3337 int no3c;
3338
3340 double p[CP];
3341
3343 double sza[CSZA];
3344
3346 double o3c[CO3];
3347
3349 double n2o[CP][CSZA][CO3];
3350
3352 double ccl4[CP][CSZA][CO3];
3353
3355 double ccl3f[CP][CSZA][CO3];
3356
3358 double ccl2f2[CP][CSZA][CO3];
3359
3361 double o2[CP][CSZA][CO3];
3362
3364 double o3_1[CP][CSZA][CO3];
3365
3367 double o3_2[CP][CSZA][CO3];
3368
3370 double h2o2[CP][CSZA][CO3];
3371
3373 double h2o[CP][CSZA][CO3];
3374
3375} clim_photo_t;
3376
3384typedef struct {
3385
3388
3390 double time[CTS];
3391
3393 double vmr[CTS];
3394
3395} clim_ts_t;
3396
3404typedef struct {
3405
3408
3410 int nlat;
3411
3413 int np;
3414
3416 double time[CT];
3417
3419 double lat[CY];
3420
3422 double p[CP];
3423
3425 double vmr[CT][CP][CY];
3426
3427} clim_zm_t;
3428
3436typedef struct {
3437
3440
3443
3445 double tropo_time[12];
3446
3448 double tropo_lat[73];
3449
3451 double tropo[12][73];
3452
3455
3458
3461
3464
3467
3470
3473
3476
3479
3482
3485
3486} clim_t;
3487
3495typedef struct {
3496
3498 double time;
3499
3501 int nx;
3502
3504 int ny;
3505
3507 int np;
3508
3510 int npl;
3511
3513 double lon[EX];
3514
3516 double lat[EY];
3517
3519 double p[EP];
3520
3522 double hybrid[EP];
3523
3525 double hyam[EP];
3526
3528 double hybm[EP];
3529
3531 double eta[EP];
3532
3534 float ps[EX][EY];
3535
3537 float ts[EX][EY];
3538
3540 float zs[EX][EY];
3541
3543 float us[EX][EY];
3544
3546 float vs[EX][EY];
3547
3549 float ess[EX][EY];
3550
3552 float nss[EX][EY];
3553
3555 float shf[EX][EY];
3556
3558 float lsm[EX][EY];
3559
3561 float sst[EX][EY];
3562
3564 float pbl[EX][EY];
3565
3567 float pt[EX][EY];
3568
3570 float tt[EX][EY];
3571
3573 float zt[EX][EY];
3574
3576 float h2ot[EX][EY];
3577
3579 float pct[EX][EY];
3580
3582 float pcb[EX][EY];
3583
3585 float cl[EX][EY];
3586
3588 float plcl[EX][EY];
3589
3591 float plfc[EX][EY];
3592
3594 float pel[EX][EY];
3595
3597 float cape[EX][EY];
3598
3600 float cin[EX][EY];
3601
3603 float o3c[EX][EY];
3604
3606 float z[EX][EY][EP];
3607
3609 float t[EX][EY][EP];
3610
3612 float u[EX][EY][EP];
3613
3615 float v[EX][EY][EP];
3616
3618 float w[EX][EY][EP];
3619
3621 float pv[EX][EY][EP];
3622
3624 float h2o[EX][EY][EP];
3625
3627 float o3[EX][EY][EP];
3628
3630 float lwc[EX][EY][EP];
3631
3633 float rwc[EX][EY][EP];
3634
3636 float iwc[EX][EY][EP];
3637
3639 float swc[EX][EY][EP];
3640
3642 float cc[EX][EY][EP];
3643
3645 float pl[EX][EY][EP];
3646
3648 float ul[EX][EY][EP];
3649
3651 float vl[EX][EY][EP];
3652
3654 float wl[EX][EY][EP];
3655
3657 float zetal[EX][EY][EP];
3658
3660 float zeta_dotl[EX][EY][EP];
3661
3662} met_t;
3663
3669typedef struct {
3670
3671 /* ------------------------------------------------------------
3672 Global grid...
3673 ------------------------------------------------------------ */
3674
3677
3680
3682 double lon_glob[DD_EX_GLOB];
3683
3685 double lat_glob[DD_EY_GLOB];
3686
3687 /* ------------------------------------------------------------
3688 Subdomains...
3689 ------------------------------------------------------------ */
3690
3692 size_t subdomain_start[4];
3693
3695 size_t subdomain_count[4];
3696
3698 size_t halo_bnd_start[4];
3699
3701 size_t halo_bnd_count[4];
3702
3705
3708
3709 /* ------------------------------------------------------------
3710 Helpers...
3711 ------------------------------------------------------------ */
3712
3713#ifdef DD
3714
3716 MPI_Datatype MPI_Particle;
3717
3719 double sort_key[NP];
3720
3722 int perm[NP];
3723
3725 double tmp[NP];
3726
3727#endif
3728
3729} dd_t;
3730
3731/* ------------------------------------------------------------
3732 Functions...
3733 ------------------------------------------------------------ */
3734
3757 void *data,
3758 size_t N);
3759
3774void cart2geo(
3775 const double *x,
3776 double *z,
3777 double *lon,
3778 double *lat);
3779
3802double clim_oh(
3803 const ctl_t * ctl,
3804 const clim_t * clim,
3805 const double t,
3806 const double lon,
3807 const double lat,
3808 const double p);
3809
3829 const ctl_t * ctl,
3830 clim_t * clim);
3831
3861double clim_photo(
3862 const double rate[CP][CSZA][CO3],
3863 const clim_photo_t * photo,
3864 const double p,
3865 const double sza,
3866 const double o3c);
3867
3893double clim_tropo(
3894 const clim_t * clim,
3895 const double t,
3896 const double lat);
3897
3916void clim_tropo_init(
3917 clim_t * clim);
3918
3935double clim_ts(
3936 const clim_ts_t * ts,
3937 const double t);
3938
3960double clim_zm(
3961 const clim_zm_t * zm,
3962 const double t,
3963 const double lat,
3964 const double p);
3965
4010 const ctl_t * ctl,
4011 const char *varname,
4012 float *array,
4013 const size_t nx,
4014 const size_t ny,
4015 const size_t np,
4016 const double *plev,
4017 const int decompress,
4018 FILE * inout);
4019
4052void compress_pck(
4053 const char *varname,
4054 float *array,
4055 const size_t nxy,
4056 const size_t nz,
4057 const int decompress,
4058 FILE * inout);
4059
4090 const char *varname,
4091 float *array,
4092 const int nx,
4093 const int ny,
4094 const int nz,
4095 const int precision,
4096 const double tolerance,
4097 const int decompress,
4098 FILE * inout);
4099
4139 const char *varname,
4140 float *array,
4141 const int nx,
4142 const int ny,
4143 const int nz,
4144 const int precision,
4145 const double tolerance,
4146 const int decompress,
4147 FILE * inout);
4148
4171 const char *varname,
4172 float *array,
4173 const size_t n,
4174 const int decompress,
4175 const int level,
4176 FILE * inout);
4177
4202double cos_sza(
4203 const double sec,
4204 const double lon,
4205 const double lat);
4206
4229void day2doy(
4230 const int year,
4231 const int mon,
4232 const int day,
4233 int *doy);
4234
4288 const ctl_t * ctl,
4289 const dd_t * dd,
4290 atm_t * atm,
4291 const int init);
4292
4345 const ctl_t * ctl,
4346 cache_t * cache,
4347 atm_t * atm,
4348 particle_t * particles,
4349 const int npart);
4350
4380 const ctl_t * ctl,
4381 const dd_t * dd,
4382 const double lon,
4383 const double lat);
4384
4441 const ctl_t * ctl,
4442 const dd_t * dd,
4443 particle_t ** particles,
4444 int *npart,
4445 int *capacity);
4446
4486 const ctl_t * ctl,
4487 dd_t * dd,
4488 atm_t * atm);
4489
4524 const dd_t * dd,
4525 double *lon,
4526 double *lat);
4527
4576 const ctl_t * ctl,
4577 cache_t * cache,
4578 const particle_t * particles,
4579 const int npart,
4580 atm_t * atm);
4581
4640 const ctl_t * ctl,
4641 const met_t * met0,
4642 atm_t * atm,
4643 dd_t * dd,
4644 int *npart);
4645
4675 double *a,
4676 dd_t * dd,
4677 const int np);
4678
4700void doy2day(
4701 const int year,
4702 const int doy,
4703 int *mon,
4704 int *day);
4705
4732void fft_help(
4733 double *fcReal,
4734 double *fcImag,
4735 const int n);
4736
4763void geo2cart(
4764 const double z,
4765 const double lon,
4766 const double lat,
4767 double *x);
4768
4793void get_met_filename(
4794 const ctl_t * ctl,
4795 const double t,
4796 const int direct,
4797 const char *metbase,
4798 const double dt_met,
4799 char *filename);
4800
4824void get_met_replace(
4825 char *orig,
4826 const char *search,
4827 const char *repl);
4828
4865void get_tropo(
4866 const int met_tropo,
4867 ctl_t * ctl,
4868 const clim_t * clim,
4869 met_t * met,
4870 const double *lons,
4871 const int nx,
4872 const double *lats,
4873 const int ny,
4874 double *pt,
4875 double *zt,
4876 double *tt,
4877 double *qt,
4878 double *o3t,
4879 double *ps,
4880 double *zs);
4881
4903 const double *lons,
4904 const int nlon,
4905 const double *lats,
4906 const int nlat,
4907 const double lon,
4908 const double lat,
4909 double *lon2,
4910 double *lat2);
4911
4954 const met_t * met0,
4955 float height0[EX][EY][EP],
4956 float array0[EX][EY][EP],
4957 const met_t * met1,
4958 float height1[EX][EY][EP],
4959 float array1[EX][EY][EP],
4960 const double ts,
4961 const double height,
4962 const double lon,
4963 const double lat,
4964 double *var,
4965 int *ci,
4966 double *cw,
4967 const int init);
4968
5004 const met_t * met,
5005 float array[EX][EY][EP],
5006 const double p,
5007 const double lon,
5008 const double lat,
5009 double *var,
5010 int *ci,
5011 double *cw,
5012 const int init);
5013
5049 const met_t * met,
5050 float array[EX][EY],
5051 const double lon,
5052 const double lat,
5053 double *var,
5054 int *ci,
5055 double *cw,
5056 const int init);
5057
5092 const met_t * met0,
5093 float array0[EX][EY][EP],
5094 const met_t * met1,
5095 float array1[EX][EY][EP],
5096 const double ts,
5097 const double p,
5098 const double lon,
5099 const double lat,
5100 double *var,
5101 int *ci,
5102 double *cw,
5103 const int init);
5104
5140 const met_t * met0,
5141 float array0[EX][EY],
5142 const met_t * met1,
5143 float array1[EX][EY],
5144 const double ts,
5145 const double lon,
5146 const double lat,
5147 double *var,
5148 int *ci,
5149 double *cw,
5150 const int init);
5151
5189void intpol_tropo_3d(
5190 const double time0,
5191 float array0[EX][EY],
5192 const double time1,
5193 float array1[EX][EY],
5194 const double lons[EX],
5195 const double lats[EY],
5196 const int nlon,
5197 const int nlat,
5198 const double time,
5199 const double lon,
5200 const double lat,
5201 const int method,
5202 double *var,
5203 double *sigma);
5204
5231void jsec2time(
5232 const double jsec,
5233 int *year,
5234 int *mon,
5235 int *day,
5236 int *hour,
5237 int *min,
5238 int *sec,
5239 double *remain);
5240
5267double kernel_weight(
5268 const double kz[EP],
5269 const double kw[EP],
5270 const int nk,
5271 const double p);
5272
5311double lapse_rate(
5312 const double t,
5313 const double h2o);
5314
5342 ctl_t * ctl);
5343
5363int locate_irr(
5364 const double *xx,
5365 const int n,
5366 const double x);
5367
5394 const float *xx,
5395 const int n,
5396 const double x,
5397 const int ig);
5398
5419int locate_reg(
5420 const double *xx,
5421 const int n,
5422 const double x);
5423
5445void locate_vert(
5446 float profiles[EX][EY][EP],
5447 const int np,
5448 const int lon_ap_ind,
5449 const int lat_ap_ind,
5450 const double alt_ap,
5451 int *ind);
5452
5478void module_advect(
5479 const ctl_t * ctl,
5480 const cache_t * cache,
5481 met_t * met0,
5482 met_t * met1,
5483 atm_t * atm);
5484
5508 const ctl_t * ctl,
5509 const cache_t * cache,
5510 met_t * met0,
5511 met_t * met1,
5512 atm_t * atm);
5513
5551 const ctl_t * ctl,
5552 const cache_t * cache,
5553 const clim_t * clim,
5554 met_t * met0,
5555 met_t * met1,
5556 atm_t * atm);
5557
5599void module_chem_grid(
5600 const ctl_t * ctl,
5601 met_t * met0,
5602 met_t * met1,
5603 atm_t * atm,
5604 const double t);
5605
5633void module_chem_init(
5634 const ctl_t * ctl,
5635 const cache_t * cache,
5636 const clim_t * clim,
5637 met_t * met0,
5638 met_t * met1,
5639 atm_t * atm);
5640
5665 const ctl_t * ctl,
5666 cache_t * cache,
5667 met_t * met0,
5668 met_t * met1,
5669 atm_t * atm);
5670
5715 const ctl_t * ctl,
5716 cache_t * cache,
5717 dd_t * dd,
5718 atm_t * atm,
5719 met_t ** met);
5720
5747void module_decay(
5748 const ctl_t * ctl,
5749 const cache_t * cache,
5750 const clim_t * clim,
5751 atm_t * atm);
5752
5789void module_diff_meso(
5790 const ctl_t * ctl,
5791 cache_t * cache,
5792 met_t * met0,
5793 met_t * met1,
5794 atm_t * atm);
5795
5829void module_diff_pbl(
5830 const ctl_t * ctl,
5831 cache_t * cache,
5832 met_t * met0,
5833 met_t * met1,
5834 atm_t * atm);
5835
5890void module_diff_turb(
5891 const ctl_t * ctl,
5892 cache_t * cache,
5893 const clim_t * clim,
5894 met_t * met0,
5895 met_t * met1,
5896 atm_t * atm);
5897
5917void module_dry_depo(
5918 const ctl_t * ctl,
5919 const cache_t * cache,
5920 met_t * met0,
5921 met_t * met1,
5922 atm_t * atm);
5923
5956void module_h2o2_chem(
5957 const ctl_t * ctl,
5958 const cache_t * cache,
5959 const clim_t * clim,
5960 met_t * met0,
5961 met_t * met1,
5962 atm_t * atm);
5963
5984 const ctl_t * ctl,
5985 cache_t * cache,
5986 met_t * met0,
5987 met_t * met1,
5988 atm_t * atm);
5989
6007void module_isosurf(
6008 const ctl_t * ctl,
6009 const cache_t * cache,
6010 met_t * met0,
6011 met_t * met1,
6012 atm_t * atm);
6013
6046 ctl_t * ctl,
6047 cache_t * cache,
6048 clim_t * clim,
6049 met_t * met0,
6050 met_t * met1,
6051 atm_t * atm);
6052
6071void module_meteo(
6072 const ctl_t * ctl,
6073 const cache_t * cache,
6074 const clim_t * clim,
6075 met_t * met0,
6076 met_t * met1,
6077 atm_t * atm);
6078
6096void module_mixing(
6097 const ctl_t * ctl,
6098 const clim_t * clim,
6099 atm_t * atm,
6100 const double t);
6101
6130 const ctl_t * ctl,
6131 const clim_t * clim,
6132 atm_t * atm,
6133 const int *ixs,
6134 const int *iys,
6135 const int *izs,
6136 const int qnt_idx,
6137 const int use_ensemble);
6138
6171void module_oh_chem(
6172 const ctl_t * ctl,
6173 const cache_t * cache,
6174 const clim_t * clim,
6175 met_t * met0,
6176 met_t * met1,
6177 atm_t * atm);
6178
6206void module_position(
6207 const cache_t * cache,
6208 met_t * met0,
6209 met_t * met1,
6210 atm_t * atm);
6211
6236void module_rng_init(
6237 const int ntask);
6238
6264void module_rng(
6265 const ctl_t * ctl,
6266 double *rs,
6267 const size_t n,
6268 const int method);
6269
6305 const ctl_t * ctl,
6306 const cache_t * cache,
6307 atm_t * atm);
6308
6331void module_sedi(
6332 const ctl_t * ctl,
6333 const cache_t * cache,
6334 met_t * met0,
6335 met_t * met1,
6336 atm_t * atm);
6337
6361void module_sort(
6362 const ctl_t * ctl,
6363 const met_t * met0,
6364 atm_t * atm);
6365
6385void module_sort_help(
6386 double *a,
6387 const int *p,
6388 const int np);
6389
6413void module_timesteps(
6414 const ctl_t * ctl,
6415 cache_t * cache,
6416 met_t * met0,
6417 atm_t * atm,
6418 const double t);
6419
6441 ctl_t * ctl,
6442 const atm_t * atm);
6443
6477 const ctl_t * ctl,
6478 const cache_t * cache,
6479 const clim_t * clim,
6480 met_t * met0,
6481 met_t * met1,
6482 atm_t * atm);
6483
6513void module_wet_depo(
6514 const ctl_t * ctl,
6515 const cache_t * cache,
6516 met_t * met0,
6517 met_t * met1,
6518 atm_t * atm);
6519
6551void mptrac_alloc(
6552 ctl_t ** ctl,
6553 cache_t ** cache,
6554 clim_t ** clim,
6555 met_t ** met0,
6556 met_t ** met1,
6557 atm_t ** atm,
6558 dd_t ** dd);
6559
6590void mptrac_free(
6591 ctl_t * ctl,
6592 cache_t * cache,
6593 clim_t * clim,
6594 met_t * met0,
6595 met_t * met1,
6596 atm_t * atm,
6597 dd_t * dd);
6598
6634void mptrac_get_met(
6635 ctl_t * ctl,
6636 clim_t * clim,
6637 const double t,
6638 met_t ** met0,
6639 met_t ** met1,
6640 dd_t * dd);
6641
6661void mptrac_init(
6662 ctl_t * ctl,
6663 cache_t * cache,
6664 clim_t * clim,
6665 atm_t * atm,
6666 const int ntask);
6667
6703int mptrac_read_atm(
6704 const char *filename,
6705 const ctl_t * ctl,
6706 atm_t * atm);
6707
6739void mptrac_read_clim(
6740 const ctl_t * ctl,
6741 clim_t * clim);
6742
6772void mptrac_read_ctl(
6773 const char *filename,
6774 int argc,
6775 char *argv[],
6776 ctl_t * ctl);
6777
6808int mptrac_read_met(
6809 const char *filename,
6810 const ctl_t * ctl,
6811 const clim_t * clim,
6812 met_t * met,
6813 dd_t * dd);
6814
6836 ctl_t * ctl,
6837 cache_t * cache,
6838 clim_t * clim,
6839 met_t ** met0,
6840 met_t ** met1,
6841 atm_t * atm,
6842 double t,
6843 dd_t * dd);
6844
6874void mptrac_write_atm(
6875 const char *filename,
6876 const ctl_t * ctl,
6877 const atm_t * atm,
6878 const double t);
6879
6917void mptrac_write_met(
6918 const char *filename,
6919 const ctl_t * ctl,
6920 met_t * met);
6921
6956 const char *dirname,
6957 const ctl_t * ctl,
6958 met_t * met0,
6959 met_t * met1,
6960 atm_t * atm,
6961 const double t);
6962
6994 const ctl_t * ctl,
6995 const cache_t * cache,
6996 const clim_t * clim,
6997 met_t ** met0,
6998 met_t ** met1,
6999 const atm_t * atm);
7000
7031 const ctl_t * ctl,
7032 const cache_t * cache,
7033 const clim_t * clim,
7034 met_t ** met0,
7035 met_t ** met1,
7036 const atm_t * atm);
7037
7065double nat_temperature(
7066 const double p,
7067 const double h2o,
7068 const double hno3);
7069
7090double pbl_weight(
7091 const ctl_t * ctl,
7092 const atm_t * atm,
7093 const int ip,
7094 const double pbl,
7095 const double ps);
7096
7129int read_atm_asc(
7130 const char *filename,
7131 const ctl_t * ctl,
7132 atm_t * atm);
7133
7164int read_atm_bin(
7165 const char *filename,
7166 const ctl_t * ctl,
7167 atm_t * atm);
7168
7193int read_atm_clams(
7194 const char *filename,
7195 const ctl_t * ctl,
7196 atm_t * atm);
7197
7227int read_atm_nc(
7228 const char *filename,
7229 const ctl_t * ctl,
7230 atm_t * atm);
7231
7260void read_clim_photo(
7261 const char *filename,
7262 clim_photo_t * photo);
7263
7281 const int ncid,
7282 const char *varname,
7283 const clim_photo_t * photo,
7284 double var[CP][CSZA][CO3]);
7285
7309int read_clim_ts(
7310 const char *filename,
7311 clim_ts_t * ts);
7312
7339void read_clim_zm(
7340 const char *filename,
7341 const char *varname,
7342 clim_zm_t * zm);
7343
7371void read_kernel(
7372 const char *filename,
7373 double kz[EP],
7374 double kw[EP],
7375 int *nk);
7376
7408int read_met_bin(
7409 const char *filename,
7410 const ctl_t * ctl,
7411 met_t * met);
7412
7438void read_met_bin_2d(
7439 FILE * in,
7440 const met_t * met,
7441 float var[EX][EY],
7442 const char *varname);
7443
7481void read_met_bin_3d(
7482 FILE * in,
7483 const ctl_t * ctl,
7484 const met_t * met,
7485 float var[EX][EY][EP],
7486 const char *varname,
7487 const float bound_min,
7488 const float bound_max);
7489
7517void read_met_cape(
7518 const ctl_t * ctl,
7519 const clim_t * clim,
7520 met_t * met);
7521
7544void read_met_cloud(
7545 met_t * met);
7546
7572void read_met_detrend(
7573 const ctl_t * ctl,
7574 met_t * met);
7575
7599 met_t * met);
7600
7627void read_met_geopot(
7628 const ctl_t * ctl,
7629 met_t * met);
7630
7658 const char *filename,
7659 const ctl_t * ctl,
7660 met_t * met);
7661
7683 codes_handle ** handles,
7684 int count_handles,
7685 met_t * met);
7686
7707 codes_handle ** handles,
7708 const int num_messages,
7709 const ctl_t * ctl,
7710 met_t * met);
7711
7742 codes_handle ** handles,
7743 const int num_messages,
7744 const ctl_t * ctl,
7745 met_t * met);
7746
7775void read_met_ml2pl(
7776 const ctl_t * ctl,
7777 const met_t * met,
7778 float var[EX][EY][EP],
7779 const char *varname);
7780
7803 const ctl_t * ctl,
7804 met_t * met);
7805
7836int read_met_nc(
7837 const char *filename,
7838 const ctl_t * ctl,
7839 met_t * met,
7840 dd_t * dd);
7841
7876void read_met_nc_grid(
7877 const char *filename,
7878 const int ncid,
7879 const ctl_t * ctl,
7880 met_t * met,
7881 dd_t * dd);
7882
7928 dd_t * dd,
7929 const ctl_t * ctl,
7930 met_t * met,
7931 const int ncid);
7932
7965 const int ncid,
7966 const ctl_t * ctl,
7967 met_t * met,
7968 dd_t * dd);
7969
8000 const int ncid,
8001 const ctl_t * ctl,
8002 met_t * met,
8003 dd_t * dd);
8004
8037int read_met_nc_2d(
8038 const int ncid,
8039 const char *varname,
8040 const char *varname2,
8041 const char *varname3,
8042 const char *varname4,
8043 const char *varname5,
8044 const char *varname6,
8045 const ctl_t * ctl,
8046 const met_t * met,
8047 dd_t * dd,
8048 float dest[EX][EY],
8049 const float scl,
8050 const int init);
8051
8081int read_met_nc_3d(
8082 const int ncid,
8083 const char *varname,
8084 const char *varname2,
8085 const char *varname3,
8086 const char *varname4,
8087 const ctl_t * ctl,
8088 const met_t * met,
8089 dd_t * dd,
8090 float dest[EX][EY][EP],
8091 const float scl);
8092
8138void read_met_pbl(
8139 const ctl_t * ctl,
8140 met_t * met);
8141
8174 met_t * met);
8175
8206 met_t * met);
8207
8238void read_met_pv(
8239 met_t * met);
8240
8263void read_met_ozone(
8264 met_t * met);
8265
8294void read_met_sample(
8295 const ctl_t * ctl,
8296 met_t * met);
8297
8326void read_met_tropo(
8327 const ctl_t * ctl,
8328 const clim_t * clim,
8329 met_t * met);
8330
8362void read_obs(
8363 const char *filename,
8364 const ctl_t * ctl,
8365 double *rt,
8366 double *rz,
8367 double *rlon,
8368 double *rlat,
8369 double *robs,
8370 int *nobs);
8371
8399void read_obs_asc(
8400 const char *filename,
8401 double *rt,
8402 double *rz,
8403 double *rlon,
8404 double *rlat,
8405 double *robs,
8406 int *nobs);
8407
8434void read_obs_nc(
8435 const char *filename,
8436 double *rt,
8437 double *rz,
8438 double *rlon,
8439 double *rlat,
8440 double *robs,
8441 int *nobs);
8442
8476double scan_ctl(
8477 const char *filename,
8478 int argc,
8479 char *argv[],
8480 const char *varname,
8481 const int arridx,
8482 const char *defvalue,
8483 char *value);
8484
8511double sedi(
8512 const double p,
8513 const double T,
8514 const double rp,
8515 const double rhop);
8516
8546void spline(
8547 const double *x,
8548 const double *y,
8549 const int n,
8550 const double *x2,
8551 double *y2,
8552 const int n2,
8553 const int method);
8554
8577float stddev(
8578 const float *data,
8579 const int n);
8580
8605void time2jsec(
8606 const int year,
8607 const int mon,
8608 const int day,
8609 const int hour,
8610 const int min,
8611 const int sec,
8612 const double remain,
8613 double *jsec);
8614
8643void timer(
8644 const char *name,
8645 const char *group,
8646 const int output);
8647
8673double time_from_filename(
8674 const char *filename,
8675 const int offset);
8676
8694double tropo_weight(
8695 const clim_t * clim,
8696 const atm_t * atm,
8697 const int ip);
8698
8721void write_atm_asc(
8722 const char *filename,
8723 const ctl_t * ctl,
8724 const atm_t * atm,
8725 const double t);
8726
8750void write_atm_bin(
8751 const char *filename,
8752 const ctl_t * ctl,
8753 const atm_t * atm);
8754
8778void write_atm_clams(
8779 const char *filename,
8780 const ctl_t * ctl,
8781 const atm_t * atm);
8782
8808 const char *dirname,
8809 const ctl_t * ctl,
8810 const atm_t * atm,
8811 const double t);
8812
8836void write_atm_nc(
8837 const char *filename,
8838 const ctl_t * ctl,
8839 const atm_t * atm);
8840
8869void write_csi(
8870 const char *filename,
8871 const ctl_t * ctl,
8872 const atm_t * atm,
8873 const double t);
8874
8916 const char *filename,
8917 const ctl_t * ctl,
8918 const atm_t * atm,
8919 const double t);
8920
8948void write_ens(
8949 const char *filename,
8950 const ctl_t * ctl,
8951 const atm_t * atm,
8952 const double t);
8953
8992void write_grid(
8993 const char *filename,
8994 const ctl_t * ctl,
8995 met_t * met0,
8996 met_t * met1,
8997 const atm_t * atm,
8998 const double t);
8999
9045void write_grid_asc(
9046 const char *filename,
9047 const ctl_t * ctl,
9048 const double *cd,
9049 double *mean[NQ],
9050 double *sigma[NQ],
9051 const double *vmr_impl,
9052 const double t,
9053 const double *z,
9054 const double *lon,
9055 const double *lat,
9056 const double *area,
9057 const double dz,
9058 const int *np);
9059
9102void write_grid_nc(
9103 const char *filename,
9104 const ctl_t * ctl,
9105 const double *cd,
9106 double *mean[NQ],
9107 double *sigma[NQ],
9108 const double *vmr_impl,
9109 const double t,
9110 const double *z,
9111 const double *lon,
9112 const double *lat,
9113 const double *area,
9114 const double dz,
9115 const int *np);
9116
9146void write_met_bin(
9147 const char *filename,
9148 const ctl_t * ctl,
9149 met_t * met);
9150
9178void write_met_bin_2d(
9179 FILE * out,
9180 met_t * met,
9181 float var[EX][EY],
9182 const char *varname);
9183
9221void write_met_bin_3d(
9222 FILE * out,
9223 const ctl_t * ctl,
9224 met_t * met,
9225 float var[EX][EY][EP],
9226 const char *varname,
9227 const int precision,
9228 const double tolerance);
9229
9257void write_met_nc(
9258 const char *filename,
9259 const ctl_t * ctl,
9260 met_t * met);
9261
9284void write_met_nc_2d(
9285 const int ncid,
9286 const char *varname,
9287 met_t * met,
9288 float var[EX][EY],
9289 const float scl);
9290
9314void write_met_nc_3d(
9315 const int ncid,
9316 const char *varname,
9317 met_t * met,
9318 float var[EX][EY][EP],
9319 const float scl);
9320
9351void write_prof(
9352 const char *filename,
9353 const ctl_t * ctl,
9354 met_t * met0,
9355 met_t * met1,
9356 const atm_t * atm,
9357 const double t);
9358
9390void write_sample(
9391 const char *filename,
9392 const ctl_t * ctl,
9393 met_t * met0,
9394 met_t * met1,
9395 const atm_t * atm,
9396 const double t);
9397
9424void write_station(
9425 const char *filename,
9426 const ctl_t * ctl,
9427 atm_t * atm,
9428 const double t);
9429
9458void write_vtk(
9459 const char *filename,
9460 const ctl_t * ctl,
9461 const atm_t * atm,
9462 const double t);
9463
9464/* ------------------------------------------------------------
9465 OpenACC routines...
9466 ------------------------------------------------------------ */
9467
9468#ifdef _OPENACC
9469#pragma acc routine (clim_oh)
9470#pragma acc routine (clim_photo)
9471#pragma acc routine (clim_tropo)
9472#pragma acc routine (clim_ts)
9473#pragma acc routine (clim_zm)
9474#pragma acc routine (cos_sza)
9475#pragma acc routine (dd_calc_subdomain_from_coords)
9476#pragma acc routine (dd_normalize_lon_lat)
9477#pragma acc routine (intpol_check_lon_lat)
9478#pragma acc routine (intpol_met_4d_zeta)
9479#pragma acc routine (intpol_met_space_3d)
9480#pragma acc routine (intpol_met_space_2d)
9481#pragma acc routine (intpol_met_time_3d)
9482#pragma acc routine (intpol_met_time_2d)
9483#pragma acc routine (kernel_weight)
9484#pragma acc routine (lapse_rate)
9485#pragma acc routine (locate_irr)
9486#pragma acc routine (locate_irr_float)
9487#pragma acc routine (locate_reg)
9488#pragma acc routine (locate_vert)
9489#pragma acc routine (nat_temperature)
9490#pragma acc routine (pbl_weight)
9491#pragma acc routine (sedi)
9492#pragma acc routine (stddev)
9493#pragma acc routine (tropo_weight)
9494#endif
9495
9496#endif /* LIBTRAC_H */
void read_met_geopot(const ctl_t *ctl, met_t *met)
Calculates geopotential heights from meteorological data.
Definition: mptrac.c:7870
void dd_init(const ctl_t *ctl, dd_t *dd, atm_t *atm)
Initialize the domain decomposition infrastructure.
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:345
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:4887
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:6526
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:7830
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:11260
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:8500
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:12661
void read_met_sample(const ctl_t *ctl, met_t *met)
Downsamples meteorological data based on specified parameters.
Definition: mptrac.c:10299
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:12405
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:10644
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:2627
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:4603
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:3852
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:6958
void read_met_cloud(met_t *met)
Calculates cloud-related variables for each grid point.
Definition: mptrac.c:7669
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:3238
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:10817
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:350
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:2083
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:8822
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:6925
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:9907
void read_met_detrend(const ctl_t *ctl, met_t *met)
Detrends meteorological data.
Definition: mptrac.c:7726
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:10472
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:10688
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:3075
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:2588
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:2607
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:8303
void read_met_monotonize(const ctl_t *ctl, met_t *met)
Makes zeta and pressure profiles monotone.
Definition: mptrac.c:9606
int dd_calc_subdomain_from_coords(const ctl_t *ctl, const dd_t *dd, const double lon, const double lat)
Determine MPI subdomain from particle coordinates.
#define DD_EY_GLOB
Maximum number of latitudes of global meteo data.
Definition: mptrac.h:420
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:7077
void read_met_periodic(met_t *met)
Applies periodic boundary conditions to meteorological data along longitudinal axis.
Definition: mptrac.c:10044
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:6260
void module_timesteps_init(ctl_t *ctl, const atm_t *atm)
Initialize start time and time interval for time-stepping.
Definition: mptrac.c:4650
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:11742
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:3959
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:340
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:4031
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:7049
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:9564
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:10716
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:7419
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:2524
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:3676
void level_definitions(ctl_t *ctl)
Defines pressure levels for meteorological data.
Definition: mptrac.c:2371
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:1853
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:12030
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:6414
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:10918
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:5098
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:2141
#define codes_handle
Placeholder when ECCODES is not available.
Definition: mptrac.h:237
void get_met_filename(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:1692
void fft_help(double *fcReal, double *fcImag, const int n)
Computes the Fast Fourier Transform (FFT) of a complex sequence.
Definition: mptrac.c:1635
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:4752
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:6721
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:10850
#define CP
Maximum number of pressure levels for climatological data.
Definition: mptrac.h:395
#define NQ
Maximum number of quantities per data point.
Definition: mptrac.h:360
void dd_assign_subdomains(const ctl_t *ctl, const dd_t *dd, atm_t *atm, const int init)
Assign or update particle subdomain ownership.
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:7131
#define EX
Maximum number of longitudes for meteo data.
Definition: mptrac.h:335
void read_met_nc_grid_dd_naive(dd_t *dd, const ctl_t *ctl, met_t *met, const int ncid)
Read meteorological grid information and construct the domain-decomposed grid with halo regions.
Definition: mptrac.c:9731
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:4462
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:9691
void timer(const char *name, const char *group, const int output)
Measures and reports elapsed time for named and grouped timers.
Definition: mptrac.c:10949
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:11075
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:2025
void module_sort(const ctl_t *ctl, const met_t *met0, atm_t *atm)
Sort particles according to box index.
Definition: mptrac.c:4491
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:3117
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:7230
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:2817
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:10745
void dd_normalize_lon_lat(const dd_t *dd, double *lon, double *lat)
Normalize geographic coordinates to the global grid convention.
#define CT
Maximum number of time steps for climatological data.
Definition: mptrac.h:405
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:2790
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:4268
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 dd_sort(const ctl_t *ctl, const met_t *met0, atm_t *atm, dd_t *dd, int *npart)
Sort local atmospheric particles and determine export counts for domain decomposition.
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:4976
void module_sort_help(double *a, const int *p, const int np)
Reorder an array based on a given permutation.
Definition: mptrac.c:4565
void dd_particles2atm(const ctl_t *ctl, cache_t *cache, const particle_t *particles, const int npart, atm_t *atm)
Copy received particles from the communication buffer into the atmospheric state.
float stddev(const float *data, const int n)
Calculates the standard deviation of a set of data.
Definition: mptrac.c:10897
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:2203
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:7448
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:2554
double time_from_filename(const char *filename, const int offset)
Extracts and converts a timestamp from a filename to Julian seconds.
Definition: mptrac.c:11017
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:12720
void mptrac_read_clim(const ctl_t *ctl, clim_t *clim)
Reads various climatological data and populates the given climatology structure.
Definition: mptrac.c:5188
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:2913
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:11052
void write_met_nc(const char *filename, const ctl_t *ctl, met_t *met)
Writes meteorological data to a NetCDF file.
Definition: mptrac.c:12497
void module_rng_init(const int ntask)
Initialize random number generators for parallel tasks.
Definition: mptrac.c:4326
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:5117
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:6470
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:6630
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:410
void read_met_ozone(met_t *met)
Calculates the total column ozone from meteorological ozone data.
Definition: mptrac.c:10270
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:8165
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:4357
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:13109
void cart2geo(const double *x, double *z, double *lon, double *lat)
Converts Cartesian coordinates to geographic coordinates.
Definition: mptrac.c:74
void dd_atm2particles(const ctl_t *ctl, cache_t *cache, atm_t *atm, particle_t *particles, const int npart)
Copy migratable atmospheric particles from the ATM state into a particle buffer.
#define NP
Maximum number of atmospheric data points.
Definition: mptrac.h:355
void dd_communicate_particles(const ctl_t *ctl, const dd_t *dd, particle_t **particles, int *npart, int *capacity)
Exchange particles between MPI ranks according to their destination rank.
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:1605
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:2170
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:4217
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 module_dd(const ctl_t *ctl, cache_t *cache, dd_t *dd, atm_t *atm, met_t **met)
Perform domain decomposition and exchange particles between MPI ranks.
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:12376
void read_met_pv(met_t *met)
Calculates potential vorticity (PV) from meteorological data.
Definition: mptrac.c:10164
#define DD_EX_GLOB
Maximum number of longitudes of global meteo data.
Definition: mptrac.h:415
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:6809
#define CY
Maximum number of latitudes for climatological data.
Definition: mptrac.h:385
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:3277
void dd_sort_help(double *a, dd_t *dd, const int np)
Apply the sorting permutation to a particle data array.
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:2294
int read_met_bin(const char *filename, const ctl_t *ctl, met_t *met)
Reads meteorological data from a binary file.
Definition: mptrac.c:7271
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:11207
void get_met_replace(char *orig, const char *search, const char *repl)
Replaces occurrences of a substring in a string with another substring.
Definition: mptrac.c:1759
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:3479
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:6865
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:6151
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:13195
#define EP
Maximum number of pressure levels for meteo data.
Definition: mptrac.h:330
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:4681
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:5248
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:4940
void read_met_polar_winds(met_t *met)
Applies a fix for polar winds in meteorological data.
Definition: mptrac.c:10105
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:3594
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:12134
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:6745
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:3354
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:12690
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:3746
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:390
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:4133
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:1674
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:7998
void get_tropo(const int met_tropo, ctl_t *ctl, const 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:1783
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:2327
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.
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:6767
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:1826
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:12947
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:11839
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:3531
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:12264
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:11157
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:7554
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:6586
#define CSZA
Maximum number of solar zenith angles for climatological data.
Definition: mptrac.h:400
double lapse_rate(const double t, const double h2o)
Calculates the moist adiabatic lapse rate in Kelvin per kilometer.
Definition: mptrac.c:2353
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:11467
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:11418
Air parcel data.
Definition: mptrac.h:3241
int np
Number of air parcels.
Definition: mptrac.h:3244
Cache data structure.
Definition: mptrac.h:3296
int iso_n
Isosurface balloon number of data points.
Definition: mptrac.h:3308
Climatological data in the form of photolysis rates.
Definition: mptrac.h:3328
int nsza
Number of solar zenith angles.
Definition: mptrac.h:3334
int np
Number of pressure levels.
Definition: mptrac.h:3331
int no3c
Number of total ozone columns.
Definition: mptrac.h:3337
Climatological data.
Definition: mptrac.h:3436
clim_ts_t ccl2f2
CFC-12 time series.
Definition: mptrac.h:3478
clim_photo_t photo
Photolysis rates.
Definition: mptrac.h:3454
clim_zm_t ho2
HO2 zonal means.
Definition: mptrac.h:3466
clim_zm_t hno3
HNO3 zonal means.
Definition: mptrac.h:3457
int tropo_ntime
Number of tropopause timesteps.
Definition: mptrac.h:3439
clim_ts_t sf6
SF6 time series.
Definition: mptrac.h:3484
clim_ts_t ccl4
CFC-10 time series.
Definition: mptrac.h:3472
clim_ts_t ccl3f
CFC-11 time series.
Definition: mptrac.h:3475
clim_zm_t o1d
O(1D) zonal means.
Definition: mptrac.h:3469
clim_zm_t h2o2
H2O2 zonal means.
Definition: mptrac.h:3463
int tropo_nlat
Number of tropopause latitudes.
Definition: mptrac.h:3442
clim_zm_t oh
OH zonal means.
Definition: mptrac.h:3460
clim_ts_t n2o
N2O time series.
Definition: mptrac.h:3481
Climatological data in the form of time series.
Definition: mptrac.h:3384
int ntime
Number of timesteps.
Definition: mptrac.h:3387
Climatological data in the form of zonal means.
Definition: mptrac.h:3404
int np
Number of pressure levels.
Definition: mptrac.h:3413
int ntime
Number of timesteps.
Definition: mptrac.h:3407
int nlat
Number of latitudes.
Definition: mptrac.h:3410
Control parameters.
Definition: mptrac.h:2190
double grid_z0
Lower altitude of gridded data [km].
Definition: mptrac.h:3105
int qnt_o3
Quantity array index for ozone volume mixing ratio.
Definition: mptrac.h:2302
double csi_lat1
Upper latitude of gridded CSI data [deg].
Definition: mptrac.h:3066
int qnt_Coh
Quantity array index for OH volume mixing ratio (chemistry code).
Definition: mptrac.h:2458
double wet_depo_ic_a
Coefficient A for wet deposition in cloud (exponential form).
Definition: mptrac.h:2954
int met_nc_scale
Check netCDF scaling factors (0=no, 1=yes).
Definition: mptrac.h:2554
int qnt_pel
Quantity array index for pressure at equilibrium level (EL).
Definition: mptrac.h:2335
int csi_nz
Number of altitudes of gridded CSI data.
Definition: mptrac.h:3042
double molmass
Molar mass [g/mol].
Definition: mptrac.h:2813
int qnt_p
Quantity array index for pressure.
Definition: mptrac.h:2281
int qnt_Cccl2f2
Quantity array index for CFC-12 volume mixing ratio (chemistry code).
Definition: mptrac.h:2482
int dd_halos_size
Domain decomposition size of halos given in grid-points.
Definition: mptrac.h:3229
int mixing_nx
Number of longitudes of mixing grid.
Definition: mptrac.h:2876
double chemgrid_z1
Upper altitude of chemistry grid [km].
Definition: mptrac.h:2900
int qnt_m
Quantity array index for mass.
Definition: mptrac.h:2221
int qnt_aoa
Quantity array index for age of air.
Definition: mptrac.h:2491
int qnt_rhop
Quantity array index for particle density.
Definition: mptrac.h:2230
int qnt_swc
Quantity array index for cloud snow water content.
Definition: mptrac.h:2314
double csi_obsmin
Minimum observation index to trigger detection.
Definition: mptrac.h:3036
int qnt_pcb
Quantity array index for cloud bottom pressure.
Definition: mptrac.h:2323
double bound_dzs
Boundary conditions surface layer depth [km].
Definition: mptrac.h:2801
double csi_lon1
Upper longitude of gridded CSI data [deg].
Definition: mptrac.h:3057
int qnt_u
Quantity array index for zonal wind.
Definition: mptrac.h:2290
double stat_lon
Longitude of station [deg].
Definition: mptrac.h:3183
double mixing_trop
Interparcel exchange parameter for mixing in the troposphere.
Definition: mptrac.h:2861
double sort_dt
Time step for sorting of particle data [s].
Definition: mptrac.h:2712
double mixing_z1
Upper altitude of mixing grid [km].
Definition: mptrac.h:2873
double stat_r
Search radius around station [km].
Definition: mptrac.h:3189
double wet_depo_bc_a
Coefficient A for wet deposition below cloud (exponential form).
Definition: mptrac.h:2948
int met_zstd_level
ZSTD compression level (from -5 to 22).
Definition: mptrac.h:2563
int csi_ny
Number of latitudes of gridded CSI data.
Definition: mptrac.h:3060
int vtk_sphere
Spherical projection for VTK data (0=no, 1=yes).
Definition: mptrac.h:3213
double chemgrid_z0
Lower altitude of chemistry grid [km].
Definition: mptrac.h:2897
double met_pbl_min
Minimum depth of planetary boundary layer [km].
Definition: mptrac.h:2680
int qnt_iwc
Quantity array index for cloud ice water content.
Definition: mptrac.h:2311
double chemgrid_lat0
Lower latitude of chemistry grid [deg].
Definition: mptrac.h:2915
double conv_cape
CAPE threshold for convection module [J/kg].
Definition: mptrac.h:2765
int qnt_Co1d
Quantity array index for O(1D) volume mixing ratio (chemistry code).
Definition: mptrac.h:2470
double met_cms_eps_pv
cmultiscale compression epsilon for potential vorticity.
Definition: mptrac.h:2602
int qnt_pw
Quantity array index for partial water vapor pressure.
Definition: mptrac.h:2389
double grid_z1
Upper altitude of gridded data [km].
Definition: mptrac.h:3108
int direction
Direction flag (1=forward calculation, -1=backward calculation).
Definition: mptrac.h:2518
int qnt_Cccl4
Quantity array index for CFC-10 volume mixing ratio (chemistry code).
Definition: mptrac.h:2476
int met_dp
Stride for pressure levels.
Definition: mptrac.h:2632
double met_dt_out
Time step for sampling of meteo data along trajectories [s].
Definition: mptrac.h:2699
int qnt_h2o2
Quantity array index for H2O2 volume mixing ratio (climatology).
Definition: mptrac.h:2353
int qnt_vh
Quantity array index for horizontal wind.
Definition: mptrac.h:2425
int csi_nx
Number of longitudes of gridded CSI data.
Definition: mptrac.h:3051
double csi_lat0
Lower latitude of gridded CSI data [deg].
Definition: mptrac.h:3063
double turb_dz_trop
Vertical turbulent diffusion coefficient (troposphere) [m^2/s].
Definition: mptrac.h:2747
int met_pbl
Planetary boundary layer data (0=file, 1=z2p, 2=Richardson, 3=theta).
Definition: mptrac.h:2677
int qnt_lwc
Quantity array index for cloud liquid water content.
Definition: mptrac.h:2305
double turb_mesoz
Vertical scaling factor for mesoscale wind fluctuations.
Definition: mptrac.h:2756
int grid_nc_level
zlib compression level of netCDF grid data files (0=off).
Definition: mptrac.h:3093
int grid_nx
Number of longitudes of gridded data.
Definition: mptrac.h:3111
int atm_type
Type of atmospheric data files (0=ASCII, 1=binary, 2=netCDF, 3=CLaMS_traj, 4=CLaMS_pos).
Definition: mptrac.h:3007
double bound_mass
Boundary conditions mass per particle [kg].
Definition: mptrac.h:2774
double grid_lat0
Lower latitude of gridded data [deg].
Definition: mptrac.h:3123
int qnt_ts
Quantity array index for surface temperature.
Definition: mptrac.h:2236
int qnt_loss_rate
Quantity array index for total loss rate.
Definition: mptrac.h:2380
double met_cms_eps_h2o
cmultiscale compression epsilon for water vapor.
Definition: mptrac.h:2605
int qnt_plfc
Quantity array index for pressure at level of free convection (LCF).
Definition: mptrac.h:2332
int qnt_Acs137
Quantity array index for radioactive activity of Cs-137.
Definition: mptrac.h:2503
double grid_lon0
Lower longitude of gridded data [deg].
Definition: mptrac.h:3114
int qnt_o1d
Quantity array index for O(1D) volume mixing ratio (climatology).
Definition: mptrac.h:2359
int met_tropo_spline
Tropopause interpolation method (0=linear, 1=spline).
Definition: mptrac.h:2696
int qnt_tvirt
Quantity array index for virtual temperature.
Definition: mptrac.h:2419
double dt_met
Time step of meteo data [s].
Definition: mptrac.h:2537
double chemgrid_lat1
Upper latitude of chemistry grid [deg].
Definition: mptrac.h:2918
int met_geopot_sy
Latitudinal smoothing of geopotential heights.
Definition: mptrac.h:2668
double met_cms_eps_u
cmultiscale compression epsilon for zonal wind.
Definition: mptrac.h:2593
double turb_dx_strat
Horizontal turbulent diffusion coefficient (stratosphere) [m^2/s].
Definition: mptrac.h:2741
int qnt_vmr
Quantity array index for volume mixing ratio.
Definition: mptrac.h:2224
int qnt_lsm
Quantity array index for land-sea mask.
Definition: mptrac.h:2257
int qnt_theta
Quantity array index for potential temperature.
Definition: mptrac.h:2401
double bound_lat1
Boundary conditions maximum longitude [deg].
Definition: mptrac.h:2789
double stat_t1
Stop time for station output [s].
Definition: mptrac.h:3195
double turb_dx_trop
Horizontal turbulent diffusion coefficient (troposphere) [m^2/s].
Definition: mptrac.h:2738
int grid_type
Type of grid data files (0=ASCII, 1=netCDF).
Definition: mptrac.h:3129
double csi_lon0
Lower longitude of gridded CSI data [deg].
Definition: mptrac.h:3054
int qnt_pbl
Quantity array index for boundary layer pressure.
Definition: mptrac.h:2263
int grid_stddev
Include standard deviations in grid output (0=no, 1=yes).
Definition: mptrac.h:3099
int qnt_psice
Quantity array index for saturation pressure over ice.
Definition: mptrac.h:2386
double chemgrid_lon0
Lower longitude of chemistry grid [deg].
Definition: mptrac.h:2906
int bound_pbl
Boundary conditions planetary boundary layer (0=no, 1=yes).
Definition: mptrac.h:2807
int qnt_mloss_wet
Quantity array index for total mass loss due to wet deposition.
Definition: mptrac.h:2371
int radio_decay
Switch for radioactive decay module (0=off, 1=on).
Definition: mptrac.h:2942
int met_geopot_sx
Longitudinal smoothing of geopotential heights.
Definition: mptrac.h:2665
int met_sy
Smoothing for latitudes.
Definition: mptrac.h:2638
int qnt_ps
Quantity array index for surface pressure.
Definition: mptrac.h:2233
int rng_type
Random number generator (0=GSL, 1=Squares, 2=cuRAND).
Definition: mptrac.h:2729
int isosurf
Isosurface parameter (0=none, 1=pressure, 2=density, 3=theta, 4=balloon).
Definition: mptrac.h:2716
double bound_p1
Boundary conditions top pressure [hPa].
Definition: mptrac.h:2795
int qnt_zs
Quantity array index for surface geopotential height.
Definition: mptrac.h:2239
int prof_nz
Number of altitudes of gridded profile data.
Definition: mptrac.h:3138
double csi_dt_out
Time step for CSI output [s].
Definition: mptrac.h:3030
int met_cape
Convective available potential energy data (0=file, 1=calculate).
Definition: mptrac.h:2674
double csi_modmin
Minimum column density to trigger detection [kg/m^2].
Definition: mptrac.h:3039
int met_sx
Smoothing for longitudes.
Definition: mptrac.h:2635
double chemgrid_lon1
Upper longitude of chemistry grid [deg].
Definition: mptrac.h:2909
double turb_mesox
Horizontal scaling factor for mesoscale wind fluctuations.
Definition: mptrac.h:2753
double met_cms_eps_iwc
cmultiscale compression epsilon for cloud ice water content.
Definition: mptrac.h:2617
double met_cms_eps_swc
cmultiscale compression epsilon for cloud snow water content.
Definition: mptrac.h:2620
double met_cms_eps_v
cmultiscale compression epsilon for meridional wind.
Definition: mptrac.h:2596
double prof_z0
Lower altitude of gridded profile data [km].
Definition: mptrac.h:3141
int qnt_w
Quantity array index for vertical velocity.
Definition: mptrac.h:2296
double bound_vmr
Boundary conditions volume mixing ratio [ppv].
Definition: mptrac.h:2780
double met_tropo_pv
Dynamical tropopause potential vorticity threshold [PVU].
Definition: mptrac.h:2690
int prof_nx
Number of longitudes of gridded profile data.
Definition: mptrac.h:3147
int qnt_stat
Quantity array index for station flag.
Definition: mptrac.h:2218
int met_tropo
Tropopause definition (0=none, 1=clim, 2=cold point, 3=WMO_1st, 4=WMO_2nd, 5=dynamical).
Definition: mptrac.h:2687
int qnt_rp
Quantity array index for particle radius.
Definition: mptrac.h:2227
int met_mpi_share
Use MPI to share meteo (0=no, 1=yes).
Definition: mptrac.h:2705
double mixing_strat
Interparcel exchange parameter for mixing in the stratosphere.
Definition: mptrac.h:2864
int qnt_vz
Quantity array index for vertical velocity.
Definition: mptrac.h:2428
int qnt_ho2
Quantity array index for HO2 volume mixing ratio (climatology).
Definition: mptrac.h:2356
double csi_z1
Upper altitude of gridded CSI data [km].
Definition: mptrac.h:3048
double stat_t0
Start time for station output [s].
Definition: mptrac.h:3192
double oh_chem_beta
Beta parameter for diurnal variablity of OH.
Definition: mptrac.h:2927
int dd
Domain decomposition (0=no, 1=yes, with 2x2 if not specified).
Definition: mptrac.h:3220
int qnt_eta
Quantity array index for eta vertical coordinate.
Definition: mptrac.h:2413
double wet_depo_so2_ph
pH value used to calculate effective Henry constant of SO2.
Definition: mptrac.h:2966
double mixing_z0
Lower altitude of mixing grid [km].
Definition: mptrac.h:2870
int qnt_mloss_decay
Quantity array index for total mass loss due to exponential decay.
Definition: mptrac.h:2377
int atm_type_out
Type of atmospheric data files for output (-1=same as ATM_TYPE, 0=ASCII, 1=binary,...
Definition: mptrac.h:3012
int met_cms_nd0x
cmultiscale number of cells of coarsest grid in x-direction.
Definition: mptrac.h:2578
int met_nlev
Number of meteo data model levels.
Definition: mptrac.h:2656
double dt_kpp
Time step for KPP chemistry [s].
Definition: mptrac.h:2936
double dry_depo_dp
Dry deposition surface layer [hPa].
Definition: mptrac.h:2975
int qnt_shf
Quantity array index for surface sensible heat flux.
Definition: mptrac.h:2254
int qnt_vs
Quantity array index for surface meridional wind.
Definition: mptrac.h:2245
int qnt_Cco
Quantity array index for CO volume mixing ratio (chemistry code).
Definition: mptrac.h:2455
double vtk_dt_out
Time step for VTK data output [s].
Definition: mptrac.h:3201
double t_stop
Stop time of simulation [s].
Definition: mptrac.h:2524
double conv_dt
Time interval for convection module [s].
Definition: mptrac.h:2771
int qnt_hno3
Quantity array index for HNO3 volume mixing ratio (climatology).
Definition: mptrac.h:2347
int met_clams
Read MPTRAC or CLaMS meteo data (0=MPTRAC, 1=CLaMS).
Definition: mptrac.h:2551
int qnt_h2ot
Quantity array index for tropopause water vapor volume mixing ratio.
Definition: mptrac.h:2275
int qnt_rh
Quantity array index for relative humidity over water.
Definition: mptrac.h:2395
double met_cms_eps_cc
cmultiscale compression epsilon for cloud cover.
Definition: mptrac.h:2623
double bound_lat0
Boundary conditions minimum longitude [deg].
Definition: mptrac.h:2786
double met_pbl_max
Maximum depth of planetary boundary layer [km].
Definition: mptrac.h:2683
int met_dx
Stride for longitudes.
Definition: mptrac.h:2626
int qnt_destination
Quantity array index for destination subdomain in domain decomposition.
Definition: mptrac.h:2515
int mixing_ny
Number of latitudes of mixing grid.
Definition: mptrac.h:2885
int met_convention
Meteo data layout (0=[lev, lat, lon], 1=[lon, lat, lev]).
Definition: mptrac.h:2540
int qnt_zeta_d
Quantity array index for diagnosed zeta vertical coordinate.
Definition: mptrac.h:2407
int tracer_chem
Switch for first order tracer chemistry module (0=off, 1=on).
Definition: mptrac.h:2939
double dt_mod
Time step of simulation [s].
Definition: mptrac.h:2527
int diffusion
Diffusion scheme (0=off, 1=fixed-K, 2=PBL).
Definition: mptrac.h:2732
int qnt_tnat
Quantity array index for T_NAT.
Definition: mptrac.h:2443
int qnt_eta_dot
Quantity array index for velocity of eta vertical coordinate.
Definition: mptrac.h:2416
int qnt_tice
Quantity array index for T_ice.
Definition: mptrac.h:2437
int qnt_zg
Quantity array index for geopotential height.
Definition: mptrac.h:2278
double vtk_offset
Vertical offset for VTK data [km].
Definition: mptrac.h:3210
int qnt_v
Quantity array index for meridional wind.
Definition: mptrac.h:2293
int qnt_mloss_dry
Quantity array index for total mass loss due to dry deposition.
Definition: mptrac.h:2374
double bound_vmr_trend
Boundary conditions volume mixing ratio trend [ppv/s].
Definition: mptrac.h:2783
int met_cache
Preload meteo data into disk cache (0=no, 1=yes).
Definition: mptrac.h:2702
int qnt_oh
Quantity array index for OH volume mixing ratio (climatology).
Definition: mptrac.h:2350
int qnt_Ch
Quantity array index for H volume mixing ratio (chemistry code).
Definition: mptrac.h:2461
int met_press_level_def
Use predefined pressure levels or not.
Definition: mptrac.h:2653
int oh_chem_reaction
Reaction type for OH chemistry (0=none, 2=bimolecular, 3=termolecular).
Definition: mptrac.h:2921
int qnt_h2o
Quantity array index for water vapor volume mixing ratio.
Definition: mptrac.h:2299
int prof_ny
Number of latitudes of gridded profile data.
Definition: mptrac.h:3156
int qnt_rhice
Quantity array index for relative humidity over ice.
Definition: mptrac.h:2398
int qnt_rho
Quantity array index for density of air.
Definition: mptrac.h:2287
double sample_dz
Layer depth for sample output [km].
Definition: mptrac.h:3177
double tdec_strat
Life time of particles in the stratosphere [s].
Definition: mptrac.h:2819
int obs_type
Type of observation data files (0=ASCII, 1=netCDF).
Definition: mptrac.h:3021
double met_cms_eps_lwc
cmultiscale compression epsilon for cloud liquid water content.
Definition: mptrac.h:2611
int qnt_us
Quantity array index for surface zonal wind.
Definition: mptrac.h:2242
double met_cms_eps_z
cmultiscale compression epsilon for geopotential height.
Definition: mptrac.h:2587
double grid_lon1
Upper longitude of gridded data [deg].
Definition: mptrac.h:3117
int qnt_Cn2o
Quantity array index for N2O volume mixing ratio (chemistry code).
Definition: mptrac.h:2485
int qnt_Cccl3f
Quantity array index for CFC-11 volume mixing ratio (chemistry code).
Definition: mptrac.h:2479
double mixing_lat0
Lower latitude of mixing grid [deg].
Definition: mptrac.h:2888
int nens
Number of ensembles.
Definition: mptrac.h:3069
int qnt_pt
Quantity array index for tropopause pressure.
Definition: mptrac.h:2266
int qnt_cl
Quantity array index for total column cloud water.
Definition: mptrac.h:2326
int advect
Advection scheme (0=off, 1=Euler, 2=midpoint, 4=Runge-Kutta).
Definition: mptrac.h:2722
double prof_z1
Upper altitude of gridded profile data [km].
Definition: mptrac.h:3144
int qnt_t
Quantity array index for temperature.
Definition: mptrac.h:2284
int atm_filter
Time filter for atmospheric data output (0=none, 1=missval, 2=remove).
Definition: mptrac.h:3000
int kpp_chem
Switch for KPP chemistry module (0=off, 1=on).
Definition: mptrac.h:2933
int qnt_zeta
Quantity array index for zeta vertical coordinate.
Definition: mptrac.h:2404
double conv_pbl_trans
Depth of PBL transition layer (fraction of PBL depth).
Definition: mptrac.h:2762
int qnt_Ai131
Quantity array index for radioactive activity of I-131.
Definition: mptrac.h:2506
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:2544
double csi_z0
Lower altitude of gridded CSI data [km].
Definition: mptrac.h:3045
int qnt_lapse
Quantity array index for lapse rate.
Definition: mptrac.h:2422
int qnt_Apb210
Quantity array index for radioactive activity of Pb-210.
Definition: mptrac.h:2497
double stat_lat
Latitude of station [deg].
Definition: mptrac.h:3186
int qnt_Cho2
Quantity array index for HO2 volume mixing ratio (chemistry code).
Definition: mptrac.h:2464
int grid_ny
Number of latitudes of gridded data.
Definition: mptrac.h:3120
int qnt_Csf6
Quantity array index for SF6 volume mixing ratio (chemistry code).
Definition: mptrac.h:2488
int qnt_Ch2o
Quantity array index for H2O volume mixing ratio (chemistry code).
Definition: mptrac.h:2449
double met_detrend
FWHM of horizontal Gaussian used for detrending [km].
Definition: mptrac.h:2644
int conv_mix_pbl
Vertical mixing in the PBL (0=off, 1=on).
Definition: mptrac.h:2759
double bound_dps
Boundary conditions surface layer depth [hPa].
Definition: mptrac.h:2798
double met_cms_eps_t
cmultiscale compression epsilon for temperature.
Definition: mptrac.h:2590
int chemgrid_nz
Number of altitudes of chemistry grid.
Definition: mptrac.h:2894
int qnt_cape
Quantity array index for convective available potential energy (CAPE).
Definition: mptrac.h:2338
int qnt_zeta_dot
Quantity array index for velocity of zeta vertical coordinate.
Definition: mptrac.h:2410
double bound_mass_trend
Boundary conditions mass per particle trend [kg/s].
Definition: mptrac.h:2777
int met_cms_nd0y
cmultiscale number of cells of coarsest grid in y-direction.
Definition: mptrac.h:2581
int mixing_nz
Number of altitudes of mixing grid.
Definition: mptrac.h:2867
int qnt_o3c
Quantity array index for total column ozone.
Definition: mptrac.h:2344
double bound_p0
Boundary conditions bottom pressure [hPa].
Definition: mptrac.h:2792
double mixing_lon0
Lower longitude of mixing grid [deg].
Definition: mptrac.h:2879
int qnt_Co3
Quantity array index for O3 volume mixing ratio (chemistry code).
Definition: mptrac.h:2452
int qnt_tsts
Quantity array index for T_STS.
Definition: mptrac.h:2440
int grid_nz
Number of altitudes of gridded data.
Definition: mptrac.h:3102
int qnt_nss
Quantity array index for northward turbulent surface stress.
Definition: mptrac.h:2251
double ens_dt_out
Time step for ensemble output [s].
Definition: mptrac.h:3075
int atm_stride
Particle index stride for atmospheric data files.
Definition: mptrac.h:3003
int met_relhum
Try to read relative humidity (0=no, 1=yes).
Definition: mptrac.h:2671
double mixing_lat1
Upper latitude of mixing grid [deg].
Definition: mptrac.h:2891
double atm_dt_out
Time step for atmospheric data output [s].
Definition: mptrac.h:2997
double prof_lat1
Upper latitude of gridded profile data [deg].
Definition: mptrac.h:3162
int met_cms_batch
cmultiscale batch size.
Definition: mptrac.h:2572
double psc_h2o
H2O volume mixing ratio for PSC analysis.
Definition: mptrac.h:2981
int met_sp
Smoothing for pressure levels.
Definition: mptrac.h:2641
double prof_lon0
Lower longitude of gridded profile data [deg].
Definition: mptrac.h:3150
int qnt_Axe133
Quantity array index for radioactive activity of Xe-133.
Definition: mptrac.h:2509
int chemgrid_nx
Number of longitudes of chemistry grid.
Definition: mptrac.h:2903
int qnt_pct
Quantity array index for cloud top pressure.
Definition: mptrac.h:2320
int qnt_mloss_kpp
Quantity array index for total mass loss due to KPP chemistry.
Definition: mptrac.h:2368
int qnt_psat
Quantity array index for saturation pressure over water.
Definition: mptrac.h:2383
int qnt_subdomain
Quantity array index for current subdomain in domain decomposition.
Definition: mptrac.h:2512
double prof_lat0
Lower latitude of gridded profile data [deg].
Definition: mptrac.h:3159
int qnt_cin
Quantity array index for convective inhibition (CIN).
Definition: mptrac.h:2341
double psc_hno3
HNO3 volume mixing ratio for PSC analysis.
Definition: mptrac.h:2984
double prof_lon1
Upper longitude of gridded profile data [deg].
Definition: mptrac.h:3153
double met_cms_eps_rwc
cmultiscale compression epsilon for cloud rain water content.
Definition: mptrac.h:2614
int met_nc_quant
Number of digits for quantization of netCDF meteo files (0=off).
Definition: mptrac.h:2560
int h2o2_chem_reaction
Reaction type for H2O2 chemistry (0=none, 1=SO2).
Definition: mptrac.h:2930
int qnt_Co3p
Quantity array index for O(3P) volume mixing ratio (chemistry code).
Definition: mptrac.h:2473
double wet_depo_bc_ret_ratio
Coefficients for wet deposition below cloud: retention ratio.
Definition: mptrac.h:2972
int chemgrid_ny
Number of latitudes of chemistry grid.
Definition: mptrac.h:2912
int qnt_Abe7
Quantity array index for radioactive activity of Be-7.
Definition: mptrac.h:2500
double met_cms_eps_o3
cmultiscale compression epsilon for ozone.
Definition: mptrac.h:2608
int met_cms_zstd
cmultiscale ZSTD compression (0=off, 1=on).
Definition: mptrac.h:2575
int met_cms_maxlev
cmultiscale maximum refinement level.
Definition: mptrac.h:2584
int grid_sparse
Sparse output in grid data files (0=no, 1=yes).
Definition: mptrac.h:3090
double dry_depo_vdep
Dry deposition velocity [m/s].
Definition: mptrac.h:2978
int qnt_tt
Quantity array index for tropopause temperature.
Definition: mptrac.h:2269
int met_np
Number of target pressure levels.
Definition: mptrac.h:2647
int qnt_ens
Quantity array index for ensemble IDs.
Definition: mptrac.h:2215
int met_nc_level
zlib compression level of netCDF meteo files (0=off).
Definition: mptrac.h:2557
double mixing_dt
Time interval for mixing [s].
Definition: mptrac.h:2858
int qnt_Arn222
Quantity array index for radioactive activity of Rn-222.
Definition: mptrac.h:2494
int qnt_mloss_h2o2
Quantity array index for total mass loss due to H2O2 chemistry.
Definition: mptrac.h:2365
double vtk_scale
Vertical scaling factor for VTK data.
Definition: mptrac.h:3207
double met_cms_eps_w
cmultiscale compression epsilon for vertical velocity.
Definition: mptrac.h:2599
double turb_dx_pbl
Horizontal turbulent diffusion coefficient (PBL) [m^2/s].
Definition: mptrac.h:2735
double conv_cin
CIN threshold for convection module [J/kg].
Definition: mptrac.h:2768
int qnt_pv
Quantity array index for potential vorticity.
Definition: mptrac.h:2431
int advect_vert_coord
Vertical velocity of air parcels (0=omega_on_plev, 1=zetadot_on_mlev, 2=omega_on_mlev,...
Definition: mptrac.h:2726
int qnt_mloss_oh
Quantity array index for total mass loss due to OH chemistry.
Definition: mptrac.h:2362
int qnt_Ch2o2
Quantity array index for H2O2 volume mixing ratio (chemistry code).
Definition: mptrac.h:2467
int qnt_sst
Quantity array index for sea surface temperature.
Definition: mptrac.h:2260
double mixing_lon1
Upper longitude of mixing grid [deg].
Definition: mptrac.h:2882
int atm_nc_level
zlib compression level of netCDF atmospheric data files (0=off).
Definition: mptrac.h:3015
double wet_depo_ic_ret_ratio
Coefficients for wet deposition in cloud: retention ratio.
Definition: mptrac.h:2969
int qnt_sh
Quantity array index for specific humidity.
Definition: mptrac.h:2392
int qnt_ess
Quantity array index for eastward turbulent surface stress.
Definition: mptrac.h:2248
double wet_depo_ic_b
Coefficient B for wet deposition in cloud (exponential form).
Definition: mptrac.h:2957
double wet_depo_bc_b
Coefficient B for wet deposition below cloud (exponential form).
Definition: mptrac.h:2951
int met_dy
Stride for latitudes.
Definition: mptrac.h:2629
int qnt_Cx
Quantity array index for trace species x volume mixing ratio (chemistry code).
Definition: mptrac.h:2446
double turb_dz_strat
Vertical turbulent diffusion coefficient (stratosphere) [m^2/s].
Definition: mptrac.h:2750
double bound_zetas
Boundary conditions surface layer zeta [K].
Definition: mptrac.h:2804
int dd_subdomains_zonal
Domain decomposition zonal subdomain number.
Definition: mptrac.h:3223
int qnt_idx
Quantity array index for air parcel IDs.
Definition: mptrac.h:2212
double met_tropo_theta
Dynamical tropopause potential temperature threshold [K].
Definition: mptrac.h:2693
int qnt_rwc
Quantity array index for cloud rain water content.
Definition: mptrac.h:2308
double t_start
Start time of simulation [s].
Definition: mptrac.h:2521
int nq
Number of quantities.
Definition: mptrac.h:2197
double tdec_trop
Life time of particles in the troposphere [s].
Definition: mptrac.h:2816
double sample_dx
Horizontal radius for sample output [km].
Definition: mptrac.h:3174
int vtk_stride
Particle index stride for VTK data.
Definition: mptrac.h:3204
double turb_dz_pbl
Vertical turbulent diffusion coefficient (PBL) [m^2/s].
Definition: mptrac.h:2744
double grid_lat1
Upper latitude of gridded data [deg].
Definition: mptrac.h:3126
int dd_subdomains_meridional
Domain decomposition meridional subdomain number.
Definition: mptrac.h:3226
int qnt_zt
Quantity array index for tropopause geopotential height.
Definition: mptrac.h:2272
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:2548
int qnt_cc
Quantity array index for cloud cover.
Definition: mptrac.h:2317
int qnt_plcl
Quantity array index for pressure at lifted condensation level (LCL).
Definition: mptrac.h:2329
double grid_dt_out
Time step for gridded data output [s].
Definition: mptrac.h:3087
int qnt_tdew
Quantity array index for dew point temperature.
Definition: mptrac.h:2434
Domain decomposition data structure.
Definition: mptrac.h:3669
int halo_offset_end
Offset of the periodic halo block at the end of the local x-array.
Definition: mptrac.h:3707
int nx_glob
Number of global longitudes.
Definition: mptrac.h:3676
int halo_offset_start
Offset of the periodic halo block at the beginning of the local x-array.
Definition: mptrac.h:3704
int ny_glob
Number of global latitudes.
Definition: mptrac.h:3679
Meteo data structure.
Definition: mptrac.h:3495
int nx
Number of longitudes.
Definition: mptrac.h:3501
int ny
Number of latitudes.
Definition: mptrac.h:3504
int np
Number of pressure levels.
Definition: mptrac.h:3507
int npl
Number of model levels.
Definition: mptrac.h:3510
double time
Time [s].
Definition: mptrac.h:3498
Particle data.
Definition: mptrac.h:3270
double p
Pressure [hPa].
Definition: mptrac.h:3276
double lat
Latitude [deg].
Definition: mptrac.h:3282
double time
Time [s].
Definition: mptrac.h:3273
double lon
Longitude [deg].
Definition: mptrac.h:3279