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#else
195#define MPI_Datatype void*
196#endif
197
198#ifdef DD
199#include <netcdf_par.h>
200#endif
201
202#ifdef _OPENACC
203#include "openacc.h"
204#endif
205
206#ifdef CURAND
207#include "curand.h"
208#endif
209
210#ifdef THRUST
211#include "thrustsort.h"
212#endif
213
214#ifdef ZFP
215#include "zfp.h"
216#endif
217
218#ifdef ZSTD
219#include "zstd.h"
220#endif
221
222#ifdef SZ3
223#include "SZ3c/sz3c.h"
224#endif
225
226#ifdef CMS
227#include "cmultiscale.h"
228#endif
229
230#ifdef KPP
231#include "chem_Parameters.h"
232#include "chem_Global.h"
233#include "chem_Sparse.h"
234#endif
235
236#ifdef ECCODES
237#include "eccodes.h"
238#else
240#define codes_handle void*
241#endif
242
243/* ------------------------------------------------------------
244 Constants...
245 ------------------------------------------------------------ */
246
248#ifndef AVO
249#define AVO 6.02214076e23
250#endif
251
253#ifndef CPD
254#define CPD 1003.5
255#endif
256
258#ifndef EPS
259#define EPS (MH2O / MA)
260#endif
261
263#ifndef G0
264#define G0 9.80665
265#endif
266
268#ifndef H0
269#define H0 7.0
270#endif
271
273#ifndef LV
274#define LV 2501000.
275#endif
276
278#ifndef KARMAN
279#define KARMAN 0.40
280#endif
281
283#ifndef KB
284#define KB 1.3806504e-23
285#endif
286
288#ifndef MA
289#define MA 28.9644
290#endif
291
293#ifndef MH2O
294#define MH2O 18.01528
295#endif
296
298#ifndef MO3
299#define MO3 48.00
300#endif
301
303#ifndef P0
304#define P0 1013.25
305#endif
306
308#ifndef RA
309#define RA (1e3 * RI / MA)
310#endif
311
313#ifndef RE
314#define RE 6367.421
315#endif
316
318#ifndef RI
319#define RI 8.3144598
320#endif
321
323#ifndef T0
324#define T0 273.15
325#endif
326
327/* ------------------------------------------------------------
328 Dimensions...
329 ------------------------------------------------------------ */
330
332#ifndef EP
333#define EP 140
334#endif
335
337#ifndef EX
338#define EX 1444
339#endif
340
342#ifndef EY
343#define EY 724
344#endif
345
347#ifndef LEN
348#define LEN 5000
349#endif
350
352#ifndef METVAR
353#define METVAR 13
354#endif
355
357#ifndef NP
358#define NP 10000000
359#endif
360
362#ifndef NQ
363#define NQ 15
364#endif
365
367#ifndef NCSI
368#define NCSI 1000000
369#endif
370
372#ifndef NENS
373#define NENS 2000
374#endif
375
377#ifndef NOBS
378#define NOBS 10000000
379#endif
380
382#ifndef NTHREADS
383#define NTHREADS 512
384#endif
385
387#ifndef CY
388#define CY 250
389#endif
390
392#ifndef CO3
393#define CO3 30
394#endif
395
397#ifndef CP
398#define CP 70
399#endif
400
402#ifndef CSZA
403#define CSZA 50
404#endif
405
407#ifndef CT
408#define CT 12
409#endif
410
412#ifndef CTS
413#define CTS 1000
414#endif
415
417#ifndef DD_NPART
418#define DD_NPART 100000
419#endif
420
422#ifndef DD_NNMAX
423#define DD_NNMAX 26
424#endif
425
427#ifndef DD_NPOLE
428#define DD_NPOLE -2
429#endif
430
432#ifndef DD_SPOLE
433#define DD_SPOLE -3
434#endif
435
436/* ------------------------------------------------------------
437 Macros...
438 ------------------------------------------------------------ */
439
459#ifdef _OPENACC
460#define ALLOC(ptr, type, n) \
461 if(acc_get_num_devices(acc_device_nvidia) <= 0) \
462 ERRMSG("Not running on a GPU device!"); \
463 if((ptr=calloc((size_t)(n), sizeof(type)))==NULL) \
464 ERRMSG("Out of memory!");
465#else
466#define ALLOC(ptr, type, n) \
467 if((ptr=calloc((size_t)(n), sizeof(type)))==NULL) \
468 ERRMSG("Out of memory!");
469#endif
470
489#define ARRAY_2D(ix, iy, ny) \
490 ((ix) * (ny) + (iy))
491
508#define ARRAY_3D(ix, iy, ny, iz, nz) \
509 (((ix)*(ny) + (iy)) * (nz) + (iz))
510
533#define ARRHENIUS(a, b, t) \
534 ((a) * exp( -(b) / (t)))
535
555#define CLAMP(v, lo, hi) \
556 (((v) < (lo)) ? (lo) : (((v) > (hi)) ? (hi) : (v)))
557
579#define DEG2DX(dlon, lat) \
580 (RE * DEG2RAD(dlon) * cos(DEG2RAD(lat)))
581
600#define DEG2DY(dlat) \
601 (RE * DEG2RAD(dlat))
602
617#define DEG2RAD(deg) \
618 ((deg) * (M_PI / 180.0))
619
642#define DP2DZ(dp, p) \
643 (- (dp) * H0 / (p))
644
664#define DX2DEG(dx, lat) \
665 (((lat) < -89.999 || (lat) > 89.999) ? 0 \
666 : (dx) * 180. / (M_PI * RE * cos(DEG2RAD(lat))))
667
682#define DY2DEG(dy) \
683 ((dy) * 180. / (M_PI * RE))
684
701#define DZ2DP(dz, p) \
702 (-(dz) * (p) / H0)
703
717#define DIST(a, b) \
718 sqrt(DIST2(a, b))
719
733#define DIST2(a, b) \
734 ((a[0]-b[0])*(a[0]-b[0])+(a[1]-b[1])*(a[1]-b[1])+(a[2]-b[2])*(a[2]-b[2]))
735
749#define DOTP(a, b) \
750 (a[0]*b[0]+a[1]*b[1]+a[2]*b[2])
751
763#define ECC(cmd) { \
764 int ecc_result=(cmd); \
765 if(ecc_result!=0) \
766 ERRMSG("ECCODES error: %s", codes_get_error_message(ecc_result)); \
767 }
768
782#define ECC_READ_2D(variable, target, scaling_factor, found_flag) { \
783 if(strcmp(short_name, variable) == 0) { \
784 for (int ix = 0; ix < met->nx; ix++) \
785 for (int iy = 0; iy < met->ny; iy++) \
786 target[ix][iy] = (float)(values[iy * met->nx + ix] * scaling_factor); \
787 found_flag = 1; \
788 } \
789 }
790
805#define ECC_READ_3D(variable, level, target, scaling_factor, found_flag) { \
806 if(strcmp(short_name, variable) == 0) { \
807 for (int ix = 0; ix < met->nx; ix++) \
808 for (int iy = 0; iy < met->ny; iy++) \
809 target[ix][iy][level] = (float) (values[iy * met->nx + ix] * scaling_factor); \
810 found_flag += 1; \
811 } \
812 }
813
830#define FMOD(x, y) \
831 ((x) - (int) ((x) / (y)) * (y))
832
848#define FREAD(ptr, type, size, in) { \
849 if(fread(ptr, sizeof(type), size, in)!=size) \
850 ERRMSG("Error while reading!"); \
851 }
852
868#define FWRITE(ptr, type, size, out) { \
869 if(fwrite(ptr, sizeof(type), size, out)!=size) \
870 ERRMSG("Error while writing!"); \
871 }
872
883#define INTPOL_INIT \
884 double cw[4] = {0.0, 0.0, 0.0, 0.0}; int ci[3] = {0, 0, 0};
885
897#define INTPOL_2D(var, init) \
898 intpol_met_time_2d(met0, met0->var, met1, met1->var, \
899 atm->time[ip], atm->lon[ip], atm->lat[ip], \
900 &var, ci, cw, init);
901
914#define INTPOL_3D(var, init) \
915 intpol_met_time_3d(met0, met0->var, met1, met1->var, \
916 atm->time[ip], atm->p[ip], \
917 atm->lon[ip], atm->lat[ip], \
918 &var, ci, cw, init);
919
933#define INTPOL_SPACE_ALL(p, lon, lat) { \
934 intpol_met_space_3d(met, met->z, p, lon, lat, &z, ci, cw, 1); \
935 intpol_met_space_3d(met, met->t, p, lon, lat, &t, ci, cw, 0); \
936 intpol_met_space_3d(met, met->u, p, lon, lat, &u, ci, cw, 0); \
937 intpol_met_space_3d(met, met->v, p, lon, lat, &v, ci, cw, 0); \
938 intpol_met_space_3d(met, met->w, p, lon, lat, &w, ci, cw, 0); \
939 intpol_met_space_3d(met, met->pv, p, lon, lat, &pv, ci, cw, 0); \
940 intpol_met_space_3d(met, met->h2o, p, lon, lat, &h2o, ci, cw, 0); \
941 intpol_met_space_3d(met, met->o3, p, lon, lat, &o3, ci, cw, 0); \
942 intpol_met_space_3d(met, met->lwc, p, lon, lat, &lwc, ci, cw, 0); \
943 intpol_met_space_3d(met, met->rwc, p, lon, lat, &rwc, ci, cw, 0); \
944 intpol_met_space_3d(met, met->iwc, p, lon, lat, &iwc, ci, cw, 0); \
945 intpol_met_space_3d(met, met->swc, p, lon, lat, &swc, ci, cw, 0); \
946 intpol_met_space_3d(met, met->cc, p, lon, lat, &cc, ci, cw, 0); \
947 intpol_met_space_2d(met, met->ps, lon, lat, &ps, ci, cw, 0); \
948 intpol_met_space_2d(met, met->ts, lon, lat, &ts, ci, cw, 0); \
949 intpol_met_space_2d(met, met->zs, lon, lat, &zs, ci, cw, 0); \
950 intpol_met_space_2d(met, met->us, lon, lat, &us, ci, cw, 0); \
951 intpol_met_space_2d(met, met->vs, lon, lat, &vs, ci, cw, 0); \
952 intpol_met_space_2d(met, met->ess, ess, lat, &ess, ci, cw, 0); \
953 intpol_met_space_2d(met, met->nss, nss, lat, &nss, ci, cw, 0); \
954 intpol_met_space_2d(met, met->shf, shf, lat, &shf, ci, cw, 0); \
955 intpol_met_space_2d(met, met->lsm, lon, lat, &lsm, ci, cw, 0); \
956 intpol_met_space_2d(met, met->sst, lon, lat, &sst, ci, cw, 0); \
957 intpol_met_space_2d(met, met->pbl, lon, lat, &pbl, ci, cw, 0); \
958 intpol_met_space_2d(met, met->pt, lon, lat, &pt, ci, cw, 0); \
959 intpol_met_space_2d(met, met->tt, lon, lat, &tt, ci, cw, 0); \
960 intpol_met_space_2d(met, met->zt, lon, lat, &zt, ci, cw, 0); \
961 intpol_met_space_2d(met, met->h2ot, lon, lat, &h2ot, ci, cw, 0); \
962 intpol_met_space_2d(met, met->pct, lon, lat, &pct, ci, cw, 0); \
963 intpol_met_space_2d(met, met->pcb, lon, lat, &pcb, ci, cw, 0); \
964 intpol_met_space_2d(met, met->cl, lon, lat, &cl, ci, cw, 0); \
965 intpol_met_space_2d(met, met->plcl, lon, lat, &plcl, ci, cw, 0); \
966 intpol_met_space_2d(met, met->plfc, lon, lat, &plfc, ci, cw, 0); \
967 intpol_met_space_2d(met, met->pel, lon, lat, &pel, ci, cw, 0); \
968 intpol_met_space_2d(met, met->cape, lon, lat, &cape, ci, cw, 0); \
969 intpol_met_space_2d(met, met->cin, lon, lat, &cin, ci, cw, 0); \
970 intpol_met_space_2d(met, met->o3c, lon, lat, &o3c, ci, cw, 0); \
971 }
972
987#define INTPOL_TIME_ALL(time, p, lon, lat) { \
988 intpol_met_time_3d(met0, met0->z, met1, met1->z, time, p, lon, lat, &z, ci, cw, 1); \
989 intpol_met_time_3d(met0, met0->t, met1, met1->t, time, p, lon, lat, &t, ci, cw, 0); \
990 intpol_met_time_3d(met0, met0->u, met1, met1->u, time, p, lon, lat, &u, ci, cw, 0); \
991 intpol_met_time_3d(met0, met0->v, met1, met1->v, time, p, lon, lat, &v, ci, cw, 0); \
992 intpol_met_time_3d(met0, met0->w, met1, met1->w, time, p, lon, lat, &w, ci, cw, 0); \
993 intpol_met_time_3d(met0, met0->pv, met1, met1->pv, time, p, lon, lat, &pv, ci, cw, 0); \
994 intpol_met_time_3d(met0, met0->h2o, met1, met1->h2o, time, p, lon, lat, &h2o, ci, cw, 0); \
995 intpol_met_time_3d(met0, met0->o3, met1, met1->o3, time, p, lon, lat, &o3, ci, cw, 0); \
996 intpol_met_time_3d(met0, met0->lwc, met1, met1->lwc, time, p, lon, lat, &lwc, ci, cw, 0); \
997 intpol_met_time_3d(met0, met0->rwc, met1, met1->rwc, time, p, lon, lat, &rwc, ci, cw, 0); \
998 intpol_met_time_3d(met0, met0->iwc, met1, met1->iwc, time, p, lon, lat, &iwc, ci, cw, 0); \
999 intpol_met_time_3d(met0, met0->swc, met1, met1->swc, time, p, lon, lat, &swc, ci, cw, 0); \
1000 intpol_met_time_3d(met0, met0->cc, met1, met1->cc, time, p, lon, lat, &cc, ci, cw, 0); \
1001 intpol_met_time_2d(met0, met0->ps, met1, met1->ps, time, lon, lat, &ps, ci, cw, 0); \
1002 intpol_met_time_2d(met0, met0->ts, met1, met1->ts, time, lon, lat, &ts, ci, cw, 0); \
1003 intpol_met_time_2d(met0, met0->zs, met1, met1->zs, time, lon, lat, &zs, ci, cw, 0); \
1004 intpol_met_time_2d(met0, met0->us, met1, met1->us, time, lon, lat, &us, ci, cw, 0); \
1005 intpol_met_time_2d(met0, met0->vs, met1, met1->vs, time, lon, lat, &vs, ci, cw, 0); \
1006 intpol_met_time_2d(met0, met0->ess, met1, met1->ess, time, lon, lat, &ess, ci, cw, 0); \
1007 intpol_met_time_2d(met0, met0->nss, met1, met1->nss, time, lon, lat, &nss, ci, cw, 0); \
1008 intpol_met_time_2d(met0, met0->shf, met1, met1->shf, time, lon, lat, &shf, ci, cw, 0); \
1009 intpol_met_time_2d(met0, met0->lsm, met1, met1->lsm, time, lon, lat, &lsm, ci, cw, 0); \
1010 intpol_met_time_2d(met0, met0->sst, met1, met1->sst, time, lon, lat, &sst, ci, cw, 0); \
1011 intpol_met_time_2d(met0, met0->pbl, met1, met1->pbl, time, lon, lat, &pbl, ci, cw, 0); \
1012 intpol_met_time_2d(met0, met0->pt, met1, met1->pt, time, lon, lat, &pt, ci, cw, 0); \
1013 intpol_met_time_2d(met0, met0->tt, met1, met1->tt, time, lon, lat, &tt, ci, cw, 0); \
1014 intpol_met_time_2d(met0, met0->zt, met1, met1->zt, time, lon, lat, &zt, ci, cw, 0); \
1015 intpol_met_time_2d(met0, met0->h2ot, met1, met1->h2ot, time, lon, lat, &h2ot, ci, cw, 0); \
1016 intpol_met_time_2d(met0, met0->pct, met1, met1->pct, time, lon, lat, &pct, ci, cw, 0); \
1017 intpol_met_time_2d(met0, met0->pcb, met1, met1->pcb, time, lon, lat, &pcb, ci, cw, 0); \
1018 intpol_met_time_2d(met0, met0->cl, met1, met1->cl, time, lon, lat, &cl, ci, cw, 0); \
1019 intpol_met_time_2d(met0, met0->plcl, met1, met1->plcl, time, lon, lat, &plcl, ci, cw, 0); \
1020 intpol_met_time_2d(met0, met0->plfc, met1, met1->plfc, time, lon, lat, &plfc, ci, cw, 0); \
1021 intpol_met_time_2d(met0, met0->pel, met1, met1->pel, time, lon, lat, &pel, ci, cw, 0); \
1022 intpol_met_time_2d(met0, met0->cape, met1, met1->cape, time, lon, lat, &cape, ci, cw, 0); \
1023 intpol_met_time_2d(met0, met0->cin, met1, met1->cin, time, lon, lat, &cin, ci, cw, 0); \
1024 intpol_met_time_2d(met0, met0->o3c, met1, met1->o3c, time, lon, lat, &o3c, ci, cw, 0); \
1025 }
1026
1041#define LAPSE(p1, t1, p2, t2) \
1042 (1e3 * G0 / RA * ((t2) - (t1)) / ((t2) + (t1)) \
1043 * ((p2) + (p1)) / ((p2) - (p1)))
1044
1060#define LIN(x0, y0, x1, y1, x) \
1061 ((y0)+((y1)-(y0))/((x1)-(x0))*((x)-(x0)))
1062
1087#define MAX(a,b) \
1088 (((a)>(b))?(a):(b))
1089
1101#define MET_HEADER \
1102 fprintf(out, \
1103 "# $1 = time [s]\n" \
1104 "# $2 = altitude [km]\n" \
1105 "# $3 = longitude [deg]\n" \
1106 "# $4 = latitude [deg]\n" \
1107 "# $5 = pressure [hPa]\n" \
1108 "# $6 = temperature [K]\n" \
1109 "# $7 = zonal wind [m/s]\n" \
1110 "# $8 = meridional wind [m/s]\n" \
1111 "# $9 = vertical velocity [hPa/s]\n" \
1112 "# $10 = H2O volume mixing ratio [ppv]\n"); \
1113 fprintf(out, \
1114 "# $11 = O3 volume mixing ratio [ppv]\n" \
1115 "# $12 = geopotential height [km]\n" \
1116 "# $13 = potential vorticity [PVU]\n" \
1117 "# $14 = surface pressure [hPa]\n" \
1118 "# $15 = surface temperature [K]\n" \
1119 "# $16 = surface geopotential height [km]\n" \
1120 "# $17 = surface zonal wind [m/s]\n" \
1121 "# $18 = surface meridional wind [m/s]\n" \
1122 "# $19 = eastward turbulent surface stress [N/m^2]\n" \
1123 "# $20 = northward turbulent surface stress [N/m^2]\n"); \
1124 fprintf(out, \
1125 "# $21 = surface sensible heat flux [W/m^2]\n" \
1126 "# $22 = land-sea mask [1]\n" \
1127 "# $23 = sea surface temperature [K]\n" \
1128 "# $24 = tropopause pressure [hPa]\n" \
1129 "# $25 = tropopause geopotential height [km]\n" \
1130 "# $26 = tropopause temperature [K]\n" \
1131 "# $27 = tropopause water vapor [ppv]\n" \
1132 "# $28 = cloud liquid water content [kg/kg]\n" \
1133 "# $29 = cloud rain water content [kg/kg]\n" \
1134 "# $30 = cloud ice water content [kg/kg]\n"); \
1135 fprintf(out, \
1136 "# $31 = cloud snow water content [kg/kg]\n" \
1137 "# $32 = cloud cover [1]\n" \
1138 "# $33 = total column cloud water [kg/m^2]\n" \
1139 "# $34 = cloud top pressure [hPa]\n" \
1140 "# $35 = cloud bottom pressure [hPa]\n" \
1141 "# $36 = pressure at lifted condensation level (LCL) [hPa]\n" \
1142 "# $37 = pressure at level of free convection (LFC) [hPa]\n" \
1143 "# $38 = pressure at equilibrium level (EL) [hPa]\n" \
1144 "# $39 = convective available potential energy (CAPE) [J/kg]\n" \
1145 "# $40 = convective inhibition (CIN) [J/kg]\n"); \
1146 fprintf(out, \
1147 "# $41 = relative humidity over water [%%]\n" \
1148 "# $42 = relative humidity over ice [%%]\n" \
1149 "# $43 = dew point temperature [K]\n" \
1150 "# $44 = frost point temperature [K]\n" \
1151 "# $45 = NAT temperature [K]\n" \
1152 "# $46 = HNO3 volume mixing ratio [ppv]\n" \
1153 "# $47 = OH volume mixing ratio [ppv]\n" \
1154 "# $48 = H2O2 volume mixing ratio [ppv]\n" \
1155 "# $49 = HO2 volume mixing ratio [ppv]\n" \
1156 "# $50 = O(1D) volume mixing ratio [ppv]\n"); \
1157 fprintf(out, \
1158 "# $51 = boundary layer pressure [hPa]\n" \
1159 "# $52 = total column ozone [DU]\n" \
1160 "# $53 = number of data points\n" \
1161 "# $54 = number of tropopause data points\n" \
1162 "# $55 = number of CAPE data points\n");
1163
1188#define MIN(a,b) \
1189 (((a)<(b))?(a):(b))
1190
1203#define MOLEC_DENS(p,t) \
1204 (AVO * 1e-6 * ((p) * 100) / (RI * (t)))
1205
1217#define NC(cmd) { \
1218 int nc_result=(cmd); \
1219 if(nc_result!=NC_NOERR) \
1220 ERRMSG("%s", nc_strerror(nc_result)); \
1221 }
1222
1246#define NC_DEF_VAR(varname, type, ndims, dims, long_name, units, level, quant) { \
1247 NC(nc_def_var(ncid, varname, type, ndims, dims, &varid)); \
1248 NC(nc_put_att_text(ncid, varid, "long_name", strnlen(long_name, LEN), long_name)); \
1249 NC(nc_put_att_text(ncid, varid, "units", strnlen(units, LEN), units)); \
1250 if((quant) > 0) \
1251 NC(nc_def_var_quantize(ncid, varid, NC_QUANTIZE_GRANULARBR, quant)); \
1252 if((level) != 0) { \
1253 NC(nc_def_var_deflate(ncid, varid, 1, 1, level)); \
1254 /* unsigned int ulevel = (unsigned int)level; */ \
1255 /* NC(nc_def_var_filter(ncid, varid, 32015, 1, (unsigned int[]){ulevel})); */ \
1256 } \
1257 }
1258
1276#define NC_GET_DOUBLE(varname, ptr, force) { \
1277 if(force) { \
1278 NC(nc_inq_varid(ncid, varname, &varid)); \
1279 NC(nc_get_var_double(ncid, varid, ptr)); \
1280 } else { \
1281 if(nc_inq_varid(ncid, varname, &varid) == NC_NOERR) { \
1282 NC(nc_get_var_double(ncid, varid, ptr)); \
1283 } else \
1284 WARN("netCDF variable %s is missing!", varname); \
1285 } \
1286 }
1287
1306#define NC_INQ_DIM(dimname, ptr, min, max, check) { \
1307 int dimid; size_t naux; \
1308 NC(nc_inq_dimid(ncid, dimname, &dimid)); \
1309 NC(nc_inq_dimlen(ncid, dimid, &naux)); \
1310 *ptr = (int)naux; \
1311 if (check) \
1312 if ((*ptr) < (min) || (*ptr) > (max)) \
1313 ERRMSG("Dimension %s is out of range!", dimname); \
1314 }
1315
1330#define NC_PUT_DOUBLE(varname, ptr, hyperslab) { \
1331 NC(nc_inq_varid(ncid, varname, &varid)); \
1332 if(hyperslab) { \
1333 NC(nc_put_vara_double(ncid, varid, start, count, ptr)); \
1334 } else { \
1335 NC(nc_put_var_double(ncid, varid, ptr)); \
1336 } \
1337 }
1338
1354#define NC_PUT_FLOAT(varname, ptr, hyperslab) { \
1355 NC(nc_inq_varid(ncid, varname, &varid)); \
1356 if(hyperslab) { \
1357 NC(nc_put_vara_float(ncid, varid, start, count, ptr)); \
1358 } else { \
1359 NC(nc_put_var_float(ncid, varid, ptr)); \
1360 } \
1361 }
1362
1377#define NC_PUT_INT(varname, ptr, hyperslab) { \
1378 NC(nc_inq_varid(ncid, varname, &varid)); \
1379 if(hyperslab) { \
1380 NC(nc_put_vara_int(ncid, varid, start, count, ptr)); \
1381 } else { \
1382 NC(nc_put_var_int(ncid, varid, ptr)); \
1383 } \
1384 }
1385
1399#define NC_PUT_ATT(varname, attname, text) { \
1400 NC(nc_inq_varid(ncid, varname, &varid)); \
1401 NC(nc_put_att_text(ncid, varid, attname, strnlen(text, LEN), text)); \
1402 }
1403
1416#define NC_PUT_ATT_GLOBAL(attname, text) \
1417 NC(nc_put_att_text(ncid, NC_GLOBAL, attname, strnlen(text, LEN), text));
1418
1436#define NN(x0, y0, x1, y1, x) \
1437 (fabs((x) - (x0)) <= fabs((x) - (x1)) ? (y0) : (y1))
1438
1454#ifdef _OPENACC
1455#define PARTICLE_LOOP(ip0, ip1, check_dt, ...) \
1456 const int ip0_const = ip0; \
1457 const int ip1_const = ip1; \
1458 _Pragma(__VA_ARGS__) \
1459 _Pragma("acc parallel loop independent gang vector") \
1460 for (int ip = ip0_const; ip < ip1_const; ip++) \
1461 if (!check_dt || cache->dt[ip] != 0)
1462#else
1463#define PARTICLE_LOOP(ip0, ip1, check_dt, ...) \
1464 const int ip0_const = ip0; \
1465 const int ip1_const = ip1; \
1466 _Pragma("omp parallel for default(shared)") \
1467 for (int ip = ip0_const; ip < ip1_const; ip++) \
1468 if (!check_dt || cache->dt[ip] != 0)
1469#endif
1470
1493#define P(z) \
1494 (P0 * exp(-(z) / H0))
1495
1517#define PSAT(t) \
1518 (6.112 * exp(17.62 * ((t) - T0) / (243.12 + (t) - T0)))
1519
1541#define PSICE(t) \
1542 (6.112 * exp(22.46 * ((t) - T0) / (272.62 + (t) - T0)))
1543
1568#define PW(p, h2o) \
1569 ((p) * MAX((h2o), 0.1e-6) / (1. + (1. - EPS) * MAX((h2o), 0.1e-6)))
1570
1585#define RAD2DEG(rad) \
1586 ((rad) * (180.0 / M_PI))
1587
1615#define RH(p, t, h2o) \
1616 (PW(p, h2o) / PSAT(t) * 100.)
1617
1645#define RHICE(p, t, h2o) \
1646 (PW(p, h2o) / PSICE(t) * 100.)
1647
1670#define RHO(p, t) \
1671 (100. * (p) / (RA * (t)))
1672
1689#define SET_ATM(qnt, val) \
1690 if (ctl->qnt >= 0) \
1691 atm->q[ctl->qnt][ip] = val;
1692
1712#define SET_QNT(qnt, name, longname, unit) \
1713 if (strcasecmp(ctl->qnt_name[iq], name) == 0) { \
1714 ctl->qnt = iq; \
1715 sprintf(ctl->qnt_longname[iq], longname); \
1716 sprintf(ctl->qnt_unit[iq], unit); \
1717 } else
1718
1733#define SH(h2o) \
1734 (EPS * MAX((h2o), 0.1e-6))
1735
1746#define SQR(x) \
1747 ((x)*(x))
1748
1760#define SWAP(x, y, type) \
1761 do {type tmp = x; x = y; y = tmp;} while(0);
1762
1784#define TDEW(p, h2o) \
1785 (T0 + 243.12 * log(PW((p), (h2o)) / 6.112) \
1786 / (17.62 - log(PW((p), (h2o)) / 6.112)))
1787
1809#define TICE(p, h2o) \
1810 (T0 + 272.62 * log(PW((p), (h2o)) / 6.112) \
1811 / (22.46 - log(PW((p), (h2o)) / 6.112)))
1812
1833#define THETA(p, t) \
1834 ((t) * pow(1000. / (p), 0.286))
1835
1862#define THETAVIRT(p, t, h2o) \
1863 (TVIRT(THETA((p), (t)), MAX((h2o), 0.1e-6)))
1864
1883#define TOK(line, tok, format, var) { \
1884 if(((tok)=strtok((line), " \t"))) { \
1885 if(sscanf(tok, format, &(var))!=1) continue; \
1886 } else ERRMSG("Error while reading!"); \
1887 }
1888
1908#define TVIRT(t, h2o) \
1909 ((t) * (1. + (1. - EPS) * MAX((h2o), 0.1e-6)))
1910
1922#define USAGE \
1923 do { \
1924 int iusage; \
1925 for (iusage = 1; iusage < argc; iusage++) \
1926 if (!strcmp(argv[iusage], "-h") \
1927 || !strcmp(argv[iusage], "--help")) { \
1928 usage(); \
1929 return EXIT_SUCCESS; \
1930 } \
1931 } while (0)
1932
1952#define Z(p) \
1953 (H0 * log(P0 / (p)))
1954
1983#define ZDIFF(lnp0, t0, h2o0, lnp1, t1, h2o1) \
1984 (RI / MA / G0 * 0.5 * (TVIRT((t0), (h2o0)) + TVIRT((t1), (h2o1))) \
1985 * ((lnp0) - (lnp1)))
1986
2002#define ZETA(ps, p, t) \
2003 (((p) / (ps) <= 0.3 ? 1. : \
2004 sin(M_PI / 2. * (1. - (p) / (ps)) / (1. - 0.3))) \
2005 * THETA((p), (t)))
2006
2007/* ------------------------------------------------------------
2008 Log messages...
2009 ------------------------------------------------------------ */
2010
2012#ifndef LOGLEV
2013#define LOGLEV 2
2014#endif
2015
2045#define LOG(level, ...) { \
2046 if(level >= 2) \
2047 printf(" "); \
2048 if(level <= LOGLEV) { \
2049 printf(__VA_ARGS__); \
2050 printf("\n"); \
2051 } \
2052 }
2053
2082#define WARN(...) { \
2083 printf("\nWarning (%s, %s, l%d): ", __FILE__, __func__, __LINE__); \
2084 LOG(0, __VA_ARGS__); \
2085 }
2086
2115#define ERRMSG(...) { \
2116 printf("\nError (%s, %s, l%d): ", __FILE__, __func__, __LINE__); \
2117 LOG(0, __VA_ARGS__); \
2118 exit(EXIT_FAILURE); \
2119 }
2120
2150#define PRINT(format, var) \
2151 printf("Print (%s, %s, l%d): %s= "format"\n", \
2152 __FILE__, __func__, __LINE__, #var, var);
2153
2154/* ------------------------------------------------------------
2155 Timers...
2156 ------------------------------------------------------------ */
2157
2159#define NTIMER 100
2160
2174#define PRINT_TIMERS \
2175 timer("END", "END", 1);
2176
2189#define SELECT_TIMER(id, group) \
2190 timer(id, group, 0);
2191
2192/* ------------------------------------------------------------
2193 Structs...
2194 ------------------------------------------------------------ */
2195
2203typedef struct {
2204
2205 /* ------------------------------------------------------------
2206 Quantity parameters...
2207 ------------------------------------------------------------ */
2208
2210 int nq;
2211
2213 char qnt_name[NQ][LEN];
2214
2216 char qnt_longname[NQ][LEN];
2217
2219 char qnt_unit[NQ][LEN];
2220
2222 char qnt_format[NQ][LEN];
2223
2226
2229
2232
2235
2238
2241
2244
2247
2250
2253
2256
2259
2262
2265
2268
2271
2274
2277
2280
2283
2286
2289
2292
2295
2298
2301
2304
2307
2310
2313
2316
2319
2322
2325
2328
2331
2334
2337
2340
2343
2346
2349
2352
2355
2358
2361
2364
2367
2370
2373
2376
2379
2382
2385
2388
2391
2394
2397
2400
2403
2406
2409
2412
2415
2418
2421
2424
2427
2430
2433
2436
2439
2442
2445
2448
2451
2454
2457
2460
2463
2466
2469
2472
2475
2478
2481
2484
2487
2490
2493
2496
2499
2502
2505
2508
2511
2514
2517
2520
2523
2526
2529
2532
2534 double t_start;
2535
2537 double t_stop;
2538
2540 double dt_mod;
2541
2542 /* ------------------------------------------------------------
2543 Meteo data parameters...
2544 ------------------------------------------------------------ */
2545
2547 char metbase[LEN];
2548
2550 double dt_met;
2551
2554
2558
2562
2565
2568
2571
2574
2577
2579 int met_comp_prec[METVAR];
2580
2582 double met_comp_tol[METVAR];
2583
2586
2589
2592
2595
2598
2601
2604
2607
2610
2613
2616
2619
2622
2625
2628
2631
2634
2637
2640
2643
2646
2649
2652
2655
2658
2661
2663 double met_p[EP];
2664
2667
2670
2672 double met_lev_hyam[EP];
2673
2675 double met_lev_hybm[EP];
2676
2679
2682
2685
2688
2691
2694
2697
2701
2704
2707
2710
2713
2716
2719
2720 /* ------------------------------------------------------------
2721 Geophysical module parameters...
2722 ------------------------------------------------------------ */
2723
2725 double sort_dt;
2726
2730
2732 char balloon[LEN];
2733
2736
2740
2743
2746
2749
2752
2755
2758
2761
2764
2767
2770
2773
2776
2779
2781 double conv_cin;
2782
2784 double conv_dt;
2785
2788
2791
2794
2797
2800
2803
2805 double bound_p0;
2806
2808 double bound_p1;
2809
2812
2815
2818
2821
2823 char species[LEN];
2824
2826 double molmass;
2827
2830
2833
2836
2838 char clim_hno3_filename[LEN];
2839
2841 char clim_oh_filename[LEN];
2842
2844 char clim_h2o2_filename[LEN];
2845
2847 char clim_ho2_filename[LEN];
2848
2850 char clim_o1d_filename[LEN];
2851
2853 char clim_o3_filename[LEN];
2854
2856 char clim_ccl4_timeseries[LEN];
2857
2859 char clim_ccl3f_timeseries[LEN];
2860
2862 char clim_ccl2f2_timeseries[LEN];
2863
2865 char clim_n2o_timeseries[LEN];
2866
2868 char clim_sf6_timeseries[LEN];
2869
2872
2875
2878
2881
2884
2887
2890
2893
2896
2899
2902
2905
2908
2911
2914
2917
2920
2923
2926
2929
2932
2935
2937 double oh_chem[4];
2938
2941
2944
2947
2949 double dt_kpp;
2950
2953
2956
2958 double wet_depo_pre[2];
2959
2962
2965
2968
2971
2973 double wet_depo_ic_h[2];
2974
2976 double wet_depo_bc_h[2];
2977
2980
2983
2986
2989
2992
2994 double psc_h2o;
2995
2997 double psc_hno3;
2998
2999 /* ------------------------------------------------------------
3000 Output parameters...
3001 ------------------------------------------------------------ */
3002
3004 char atm_basename[LEN];
3005
3007 char atm_gpfile[LEN];
3008
3011
3014
3017
3021
3026
3029
3031 int atm_nc_quant[NQ];
3032
3035
3037 char csi_basename[LEN];
3038
3040 char csi_kernel[LEN];
3041
3044
3046 char csi_obsfile[LEN];
3047
3050
3053
3056
3058 double csi_z0;
3059
3061 double csi_z1;
3062
3065
3067 double csi_lon0;
3068
3070 double csi_lon1;
3071
3074
3076 double csi_lat0;
3077
3079 double csi_lat1;
3080
3082 int nens;
3083
3085 char ens_basename[LEN];
3086
3089
3091 char grid_basename[LEN];
3092
3094 char grid_kernel[LEN];
3095
3097 char grid_gpfile[LEN];
3098
3101
3104
3107
3109 int grid_nc_quant[NQ];
3110
3113
3116
3118 double grid_z0;
3119
3121 double grid_z1;
3122
3125
3128
3131
3134
3137
3140
3143
3145 char prof_basename[LEN];
3146
3148 char prof_obsfile[LEN];
3149
3152
3154 double prof_z0;
3155
3157 double prof_z1;
3158
3161
3164
3167
3170
3173
3176
3178 char sample_basename[LEN];
3179
3181 char sample_kernel[LEN];
3182
3184 char sample_obsfile[LEN];
3185
3188
3191
3193 char stat_basename[LEN];
3194
3196 double stat_lon;
3197
3199 double stat_lat;
3200
3202 double stat_r;
3203
3205 double stat_t0;
3206
3208 double stat_t1;
3209
3211 char vtk_basename[LEN];
3212
3215
3218
3221
3224
3227
3228 /* ------------------------------------------------------------
3229 Domain decomposition...
3230 ------------------------------------------------------------ */
3231
3233 int dd;
3234
3237
3240
3243
3246
3247} ctl_t;
3248
3257typedef struct {
3258
3260 int np;
3261
3263 double time[NP];
3264
3266 double p[NP];
3267
3269 double lon[NP];
3270
3272 double lat[NP];
3273
3275 double q[NQ][NP];
3276
3277} atm_t;
3278
3286typedef struct {
3287
3289 double time;
3290
3292 double p;
3293
3295 double lon;
3296
3298 double lat;
3299
3301 double q[NQ];
3302
3303} particle_t;
3304
3305
3312typedef struct {
3313
3315 double iso_var[NP];
3316
3318 double iso_ps[NP];
3319
3321 double iso_ts[NP];
3322
3325
3327 float uvwp[NP][3];
3328
3330 double rs[3 * NP + 1];
3331
3333 double dt[NP];
3334
3335} cache_t;
3336
3344typedef struct {
3345
3347 int np;
3348
3350 int nsza;
3351
3353 int no3c;
3354
3356 double p[CP];
3357
3359 double sza[CSZA];
3360
3362 double o3c[CO3];
3363
3365 double n2o[CP][CSZA][CO3];
3366
3368 double ccl4[CP][CSZA][CO3];
3369
3371 double ccl3f[CP][CSZA][CO3];
3372
3374 double ccl2f2[CP][CSZA][CO3];
3375
3377 double o2[CP][CSZA][CO3];
3378
3380 double o3_1[CP][CSZA][CO3];
3381
3383 double o3_2[CP][CSZA][CO3];
3384
3386 double h2o2[CP][CSZA][CO3];
3387
3389 double h2o[CP][CSZA][CO3];
3390
3391} clim_photo_t;
3392
3400typedef struct {
3401
3404
3406 double time[CTS];
3407
3409 double vmr[CTS];
3410
3411} clim_ts_t;
3412
3420typedef struct {
3421
3424
3426 int nlat;
3427
3429 int np;
3430
3432 double time[CT];
3433
3435 double lat[CY];
3436
3438 double p[CP];
3439
3441 double vmr[CT][CP][CY];
3442
3443} clim_zm_t;
3444
3452typedef struct {
3453
3456
3459
3461 double tropo_time[12];
3462
3464 double tropo_lat[73];
3465
3467 double tropo[12][73];
3468
3471
3474
3477
3480
3483
3486
3489
3492
3495
3498
3501
3502} clim_t;
3503
3511typedef struct {
3512
3514 double time;
3515
3517 int nx;
3518
3520 int ny;
3521
3523 int np;
3524
3526 int npl;
3527
3529 double lon[EX];
3530
3532 double lat[EY];
3533
3535 double p[EP];
3536
3538 double hybrid[EP];
3539
3541 double hyam[EP];
3542
3544 double hybm[EP];
3545
3547 double eta[EP];
3548
3550 float ps[EX][EY];
3551
3553 float ts[EX][EY];
3554
3556 float zs[EX][EY];
3557
3559 float us[EX][EY];
3560
3562 float vs[EX][EY];
3563
3565 float ess[EX][EY];
3566
3568 float nss[EX][EY];
3569
3571 float shf[EX][EY];
3572
3574 float lsm[EX][EY];
3575
3577 float sst[EX][EY];
3578
3580 float pbl[EX][EY];
3581
3583 float pt[EX][EY];
3584
3586 float tt[EX][EY];
3587
3589 float zt[EX][EY];
3590
3592 float h2ot[EX][EY];
3593
3595 float pct[EX][EY];
3596
3598 float pcb[EX][EY];
3599
3601 float cl[EX][EY];
3602
3604 float plcl[EX][EY];
3605
3607 float plfc[EX][EY];
3608
3610 float pel[EX][EY];
3611
3613 float cape[EX][EY];
3614
3616 float cin[EX][EY];
3617
3619 float o3c[EX][EY];
3620
3622 float z[EX][EY][EP];
3623
3625 float t[EX][EY][EP];
3626
3628 float u[EX][EY][EP];
3629
3631 float v[EX][EY][EP];
3632
3634 float w[EX][EY][EP];
3635
3637 float pv[EX][EY][EP];
3638
3640 float h2o[EX][EY][EP];
3641
3643 float o3[EX][EY][EP];
3644
3646 float lwc[EX][EY][EP];
3647
3649 float rwc[EX][EY][EP];
3650
3652 float iwc[EX][EY][EP];
3653
3655 float swc[EX][EY][EP];
3656
3658 float cc[EX][EY][EP];
3659
3661 float pl[EX][EY][EP];
3662
3664 float ul[EX][EY][EP];
3665
3667 float vl[EX][EY][EP];
3668
3670 float wl[EX][EY][EP];
3671
3673 float zetal[EX][EY][EP];
3674
3676 float zeta_dotl[EX][EY][EP];
3677
3678} met_t;
3679
3685typedef struct {
3686
3687 /* ------------------------------------------------------------
3688 MPI Information
3689 ------------------------------------------------------------ */
3690
3692 int rank;
3693
3695 int size;
3696
3698 int neighbours[DD_NNMAX];
3699
3700#ifdef DD
3702 MPI_Datatype MPI_Particle;
3703#endif
3704
3705 /* ------------------------------------------------------------
3706 Properties of subdomains
3707 ------------------------------------------------------------ */
3708
3711
3714
3717
3720
3722 size_t subdomain_start[4];
3723
3725 size_t subdomain_count[4];
3726
3728 size_t halo_bnd_start[4];
3729
3731 size_t halo_bnd_count[4];
3732
3735
3738
3739 /* ------------------------------------------------------------
3740 Caches
3741 ------------------------------------------------------------ */
3742
3744 int init;
3745
3746#ifdef DD
3747
3749 double a[NP];
3750
3752 int p[NP];
3753
3755 double help[NP];
3756
3757#endif
3758
3759} dd_t;
3760
3761/* ------------------------------------------------------------
3762 OpenACC routines...
3763 ------------------------------------------------------------ */
3764
3765#ifdef _OPENACC
3766#pragma acc routine (clim_oh)
3767#pragma acc routine (clim_photo)
3768#pragma acc routine (clim_tropo)
3769#pragma acc routine (clim_ts)
3770#pragma acc routine (clim_zm)
3771#pragma acc routine (cos_sza)
3772#pragma acc routine (intpol_check_lon_lat)
3773#pragma acc routine (intpol_met_4d_zeta)
3774#pragma acc routine (intpol_met_space_3d)
3775#pragma acc routine (intpol_met_space_2d)
3776#pragma acc routine (intpol_met_time_3d)
3777#pragma acc routine (intpol_met_time_2d)
3778#pragma acc routine (kernel_weight)
3779#pragma acc routine (lapse_rate)
3780#pragma acc routine (locate_irr)
3781#pragma acc routine (locate_irr_float)
3782#pragma acc routine (locate_reg)
3783#pragma acc routine (locate_vert)
3784#pragma acc routine (nat_temperature)
3785#pragma acc routine (pbl_weight)
3786#pragma acc routine (sedi)
3787#pragma acc routine (stddev)
3788#pragma acc routine (tropo_weight)
3789#endif
3790
3791/* ------------------------------------------------------------
3792 Functions...
3793 ------------------------------------------------------------ */
3794
3817 void *data,
3818 size_t N);
3819
3834void cart2geo(
3835 const double *x,
3836 double *z,
3837 double *lon,
3838 double *lat);
3839
3862double clim_oh(
3863 const ctl_t * ctl,
3864 const clim_t * clim,
3865 const double t,
3866 const double lon,
3867 const double lat,
3868 const double p);
3869
3889 const ctl_t * ctl,
3890 clim_t * clim);
3891
3921double clim_photo(
3922 const double rate[CP][CSZA][CO3],
3923 const clim_photo_t * photo,
3924 const double p,
3925 const double sza,
3926 const double o3c);
3927
3953double clim_tropo(
3954 const clim_t * clim,
3955 const double t,
3956 const double lat);
3957
3976void clim_tropo_init(
3977 clim_t * clim);
3978
3995double clim_ts(
3996 const clim_ts_t * ts,
3997 const double t);
3998
4020double clim_zm(
4021 const clim_zm_t * zm,
4022 const double t,
4023 const double lat,
4024 const double p);
4025
4070 const ctl_t * ctl,
4071 const char *varname,
4072 float *array,
4073 const size_t nx,
4074 const size_t ny,
4075 const size_t np,
4076 const double *plev,
4077 const int decompress,
4078 FILE * inout);
4079
4112void compress_pck(
4113 const char *varname,
4114 float *array,
4115 const size_t nxy,
4116 const size_t nz,
4117 const int decompress,
4118 FILE * inout);
4119
4150 const char *varname,
4151 float *array,
4152 const int nx,
4153 const int ny,
4154 const int nz,
4155 const int precision,
4156 const double tolerance,
4157 const int decompress,
4158 FILE * inout);
4159
4199 const char *varname,
4200 float *array,
4201 const int nx,
4202 const int ny,
4203 const int nz,
4204 const int precision,
4205 const double tolerance,
4206 const int decompress,
4207 FILE * inout);
4208
4231 const char *varname,
4232 float *array,
4233 const size_t n,
4234 const int decompress,
4235 const int level,
4236 FILE * inout);
4237
4262double cos_sza(
4263 const double sec,
4264 const double lon,
4265 const double lat);
4266
4289void day2doy(
4290 const int year,
4291 const int mon,
4292 const int day,
4293 int *doy);
4294
4342 atm_t * atm,
4343 ctl_t * ctl,
4344 dd_t * dd,
4345 int init);
4346
4412 atm_t * atm,
4413 particle_t * particles,
4414 ctl_t * ctl,
4415 int *nparticles,
4416 cache_t * cache,
4417 int rank);
4418
4483 double lon,
4484 double lat,
4485 met_t * met,
4486 ctl_t * ctl,
4487 int mpi_size,
4488 int nx_glob,
4489 int ny_glob);
4490
4521 particle_t * particles,
4522 int *nparticles,
4523 MPI_Datatype MPI_Particle,
4524 int *neighbours,
4525 int nneighbours,
4526 ctl_t ctl);
4527
4551 const ctl_t ctl,
4552 dd_t * dd);
4553
4579 ctl_t * ctl,
4580 dd_t * dd,
4581 atm_t * atm);
4582
4607 met_t * met,
4608 int nx_glob);
4609
4635 atm_t * atm,
4636 particle_t * particles,
4637 ctl_t * ctl,
4638 int *nparticles,
4639 cache_t * cache);
4640
4663 MPI_Datatype * MPI_Particle);
4664
4700 const ctl_t * ctl,
4701 met_t * met0,
4702 atm_t * atm,
4703 dd_t * dd,
4704 int *nparticles,
4705 int *rank);
4706
4735 double *a,
4736 dd_t * dd,
4737 const int np);
4738
4760void doy2day(
4761 const int year,
4762 const int doy,
4763 int *mon,
4764 int *day);
4765
4792void fft_help(
4793 double *fcReal,
4794 double *fcImag,
4795 const int n);
4796
4823void geo2cart(
4824 const double z,
4825 const double lon,
4826 const double lat,
4827 double *x);
4828
4853void get_met_help(
4854 const ctl_t * ctl,
4855 const double t,
4856 const int direct,
4857 const char *metbase,
4858 const double dt_met,
4859 char *filename);
4860
4884void get_met_replace(
4885 char *orig,
4886 char *search,
4887 char *repl);
4888
4925void get_tropo(
4926 const int met_tropo,
4927 ctl_t * ctl,
4928 clim_t * clim,
4929 met_t * met,
4930 const double *lons,
4931 const int nx,
4932 const double *lats,
4933 const int ny,
4934 double *pt,
4935 double *zt,
4936 double *tt,
4937 double *qt,
4938 double *o3t,
4939 double *ps,
4940 double *zs);
4941
4963 const double *lons,
4964 const int nlon,
4965 const double *lats,
4966 const int nlat,
4967 const double lon,
4968 const double lat,
4969 double *lon2,
4970 double *lat2);
4971
5014 const met_t * met0,
5015 float height0[EX][EY][EP],
5016 float array0[EX][EY][EP],
5017 const met_t * met1,
5018 float height1[EX][EY][EP],
5019 float array1[EX][EY][EP],
5020 const double ts,
5021 const double height,
5022 const double lon,
5023 const double lat,
5024 double *var,
5025 int *ci,
5026 double *cw,
5027 const int init);
5028
5064 const met_t * met,
5065 float array[EX][EY][EP],
5066 const double p,
5067 const double lon,
5068 const double lat,
5069 double *var,
5070 int *ci,
5071 double *cw,
5072 const int init);
5073
5109 const met_t * met,
5110 float array[EX][EY],
5111 const double lon,
5112 const double lat,
5113 double *var,
5114 int *ci,
5115 double *cw,
5116 const int init);
5117
5152 const met_t * met0,
5153 float array0[EX][EY][EP],
5154 const met_t * met1,
5155 float array1[EX][EY][EP],
5156 const double ts,
5157 const double p,
5158 const double lon,
5159 const double lat,
5160 double *var,
5161 int *ci,
5162 double *cw,
5163 const int init);
5164
5200 const met_t * met0,
5201 float array0[EX][EY],
5202 const met_t * met1,
5203 float array1[EX][EY],
5204 const double ts,
5205 const double lon,
5206 const double lat,
5207 double *var,
5208 int *ci,
5209 double *cw,
5210 const int init);
5211
5249void intpol_tropo_3d(
5250 const double time0,
5251 float array0[EX][EY],
5252 const double time1,
5253 float array1[EX][EY],
5254 const double lons[EX],
5255 const double lats[EY],
5256 const int nlon,
5257 const int nlat,
5258 const double time,
5259 const double lon,
5260 const double lat,
5261 const int method,
5262 double *var,
5263 double *sigma);
5264
5291void jsec2time(
5292 const double jsec,
5293 int *year,
5294 int *mon,
5295 int *day,
5296 int *hour,
5297 int *min,
5298 int *sec,
5299 double *remain);
5300
5327double kernel_weight(
5328 const double kz[EP],
5329 const double kw[EP],
5330 const int nk,
5331 const double p);
5332
5371double lapse_rate(
5372 const double t,
5373 const double h2o);
5374
5404 ctl_t * ctl);
5405
5425int locate_irr(
5426 const double *xx,
5427 const int n,
5428 const double x);
5429
5456 const float *xx,
5457 const int n,
5458 const double x,
5459 const int ig);
5460
5481int locate_reg(
5482 const double *xx,
5483 const int n,
5484 const double x);
5485
5507void locate_vert(
5508 float profiles[EX][EY][EP],
5509 const int np,
5510 const int lon_ap_ind,
5511 const int lat_ap_ind,
5512 const double alt_ap,
5513 int *ind);
5514
5540void module_advect(
5541 const ctl_t * ctl,
5542 const cache_t * cache,
5543 met_t * met0,
5544 met_t * met1,
5545 atm_t * atm);
5546
5570 const ctl_t * ctl,
5571 const cache_t * cache,
5572 met_t * met0,
5573 met_t * met1,
5574 atm_t * atm);
5575
5613 const ctl_t * ctl,
5614 const cache_t * cache,
5615 const clim_t * clim,
5616 met_t * met0,
5617 met_t * met1,
5618 atm_t * atm);
5619
5661void module_chem_grid(
5662 const ctl_t * ctl,
5663 met_t * met0,
5664 met_t * met1,
5665 atm_t * atm,
5666 const double t);
5667
5695void module_chem_init(
5696 const ctl_t * ctl,
5697 const cache_t * cache,
5698 const clim_t * clim,
5699 met_t * met0,
5700 met_t * met1,
5701 atm_t * atm);
5702
5727 const ctl_t * ctl,
5728 cache_t * cache,
5729 met_t * met0,
5730 met_t * met1,
5731 atm_t * atm);
5732
5761 ctl_t * ctl,
5762 atm_t * atm,
5763 cache_t * cache,
5764 dd_t * dd,
5765 met_t ** met);
5766
5793void module_decay(
5794 const ctl_t * ctl,
5795 const cache_t * cache,
5796 const clim_t * clim,
5797 atm_t * atm);
5798
5835void module_diff_meso(
5836 const ctl_t * ctl,
5837 cache_t * cache,
5838 met_t * met0,
5839 met_t * met1,
5840 atm_t * atm);
5841
5875void module_diff_pbl(
5876 const ctl_t * ctl,
5877 cache_t * cache,
5878 met_t * met0,
5879 met_t * met1,
5880 atm_t * atm);
5881
5936void module_diff_turb(
5937 const ctl_t * ctl,
5938 cache_t * cache,
5939 const clim_t * clim,
5940 met_t * met0,
5941 met_t * met1,
5942 atm_t * atm);
5943
5963void module_dry_depo(
5964 const ctl_t * ctl,
5965 const cache_t * cache,
5966 met_t * met0,
5967 met_t * met1,
5968 atm_t * atm);
5969
6002void module_h2o2_chem(
6003 const ctl_t * ctl,
6004 const cache_t * cache,
6005 const clim_t * clim,
6006 met_t * met0,
6007 met_t * met1,
6008 atm_t * atm);
6009
6030 const ctl_t * ctl,
6031 cache_t * cache,
6032 met_t * met0,
6033 met_t * met1,
6034 atm_t * atm);
6035
6053void module_isosurf(
6054 const ctl_t * ctl,
6055 const cache_t * cache,
6056 met_t * met0,
6057 met_t * met1,
6058 atm_t * atm);
6059
6092 ctl_t * ctl,
6093 cache_t * cache,
6094 clim_t * clim,
6095 met_t * met0,
6096 met_t * met1,
6097 atm_t * atm);
6098
6117void module_meteo(
6118 const ctl_t * ctl,
6119 const cache_t * cache,
6120 const clim_t * clim,
6121 met_t * met0,
6122 met_t * met1,
6123 atm_t * atm);
6124
6142void module_mixing(
6143 const ctl_t * ctl,
6144 const clim_t * clim,
6145 atm_t * atm,
6146 const double t);
6147
6176 const ctl_t * ctl,
6177 const clim_t * clim,
6178 atm_t * atm,
6179 const int *ixs,
6180 const int *iys,
6181 const int *izs,
6182 const int qnt_idx,
6183 const int use_ensemble);
6184
6217void module_oh_chem(
6218 const ctl_t * ctl,
6219 const cache_t * cache,
6220 const clim_t * clim,
6221 met_t * met0,
6222 met_t * met1,
6223 atm_t * atm);
6224
6252void module_position(
6253 const cache_t * cache,
6254 met_t * met0,
6255 met_t * met1,
6256 atm_t * atm);
6257
6282void module_rng_init(
6283 const int ntask);
6284
6310void module_rng(
6311 const ctl_t * ctl,
6312 double *rs,
6313 const size_t n,
6314 const int method);
6315
6351 const ctl_t * ctl,
6352 const cache_t * cache,
6353 atm_t * atm);
6354
6377void module_sedi(
6378 const ctl_t * ctl,
6379 const cache_t * cache,
6380 met_t * met0,
6381 met_t * met1,
6382 atm_t * atm);
6383
6407void module_sort(
6408 const ctl_t * ctl,
6409 met_t * met0,
6410 atm_t * atm);
6411
6431void module_sort_help(
6432 double *a,
6433 const int *p,
6434 const int np);
6435
6459void module_timesteps(
6460 const ctl_t * ctl,
6461 cache_t * cache,
6462 met_t * met0,
6463 atm_t * atm,
6464 const double t);
6465
6487 ctl_t * ctl,
6488 const atm_t * atm);
6489
6523 const ctl_t * ctl,
6524 const cache_t * cache,
6525 const clim_t * clim,
6526 met_t * met0,
6527 met_t * met1,
6528 atm_t * atm);
6529
6559void module_wet_depo(
6560 const ctl_t * ctl,
6561 const cache_t * cache,
6562 met_t * met0,
6563 met_t * met1,
6564 atm_t * atm);
6565
6597void mptrac_alloc(
6598 ctl_t ** ctl,
6599 cache_t ** cache,
6600 clim_t ** clim,
6601 met_t ** met0,
6602 met_t ** met1,
6603 atm_t ** atm,
6604 dd_t ** dd);
6605
6636void mptrac_free(
6637 ctl_t * ctl,
6638 cache_t * cache,
6639 clim_t * clim,
6640 met_t * met0,
6641 met_t * met1,
6642 atm_t * atm,
6643 dd_t * dd);
6644
6680void mptrac_get_met(
6681 ctl_t * ctl,
6682 clim_t * clim,
6683 const double t,
6684 met_t ** met0,
6685 met_t ** met1,
6686 dd_t * dd);
6687
6707void mptrac_init(
6708 ctl_t * ctl,
6709 cache_t * cache,
6710 clim_t * clim,
6711 atm_t * atm,
6712 const int ntask);
6713
6749int mptrac_read_atm(
6750 const char *filename,
6751 const ctl_t * ctl,
6752 atm_t * atm);
6753
6785void mptrac_read_clim(
6786 const ctl_t * ctl,
6787 clim_t * clim);
6788
6818void mptrac_read_ctl(
6819 const char *filename,
6820 int argc,
6821 char *argv[],
6822 ctl_t * ctl);
6823
6854int mptrac_read_met(
6855 const char *filename,
6856 const ctl_t * ctl,
6857 const clim_t * clim,
6858 met_t * met,
6859 dd_t * dd);
6860
6882 ctl_t * ctl,
6883 cache_t * cache,
6884 clim_t * clim,
6885 met_t ** met0,
6886 met_t ** met1,
6887 atm_t * atm,
6888 double t,
6889 dd_t * dd);
6890
6920void mptrac_write_atm(
6921 const char *filename,
6922 const ctl_t * ctl,
6923 const atm_t * atm,
6924 const double t);
6925
6963void mptrac_write_met(
6964 const char *filename,
6965 const ctl_t * ctl,
6966 met_t * met);
6967
7002 const char *dirname,
7003 const ctl_t * ctl,
7004 met_t * met0,
7005 met_t * met1,
7006 atm_t * atm,
7007 const double t);
7008
7040 const ctl_t * ctl,
7041 const cache_t * cache,
7042 const clim_t * clim,
7043 met_t ** met0,
7044 met_t ** met1,
7045 const atm_t * atm);
7046
7077 const ctl_t * ctl,
7078 const cache_t * cache,
7079 const clim_t * clim,
7080 met_t ** met0,
7081 met_t ** met1,
7082 const atm_t * atm);
7083
7111double nat_temperature(
7112 const double p,
7113 const double h2o,
7114 const double hno3);
7115
7136double pbl_weight(
7137 const ctl_t * ctl,
7138 const atm_t * atm,
7139 const int ip,
7140 const double pbl,
7141 const double ps);
7142
7175int read_atm_asc(
7176 const char *filename,
7177 const ctl_t * ctl,
7178 atm_t * atm);
7179
7210int read_atm_bin(
7211 const char *filename,
7212 const ctl_t * ctl,
7213 atm_t * atm);
7214
7239int read_atm_clams(
7240 const char *filename,
7241 const ctl_t * ctl,
7242 atm_t * atm);
7243
7273int read_atm_nc(
7274 const char *filename,
7275 const ctl_t * ctl,
7276 atm_t * atm);
7277
7306void read_clim_photo(
7307 const char *filename,
7308 clim_photo_t * photo);
7309
7327 const int ncid,
7328 const char *varname,
7329 const clim_photo_t * photo,
7330 double var[CP][CSZA][CO3]);
7331
7355int read_clim_ts(
7356 const char *filename,
7357 clim_ts_t * ts);
7358
7385void read_clim_zm(
7386 const char *filename,
7387 const char *varname,
7388 clim_zm_t * zm);
7389
7417void read_kernel(
7418 const char *filename,
7419 double kz[EP],
7420 double kw[EP],
7421 int *nk);
7422
7454int read_met_bin(
7455 const char *filename,
7456 const ctl_t * ctl,
7457 met_t * met);
7458
7484void read_met_bin_2d(
7485 FILE * in,
7486 const met_t * met,
7487 float var[EX][EY],
7488 const char *varname);
7489
7527void read_met_bin_3d(
7528 FILE * in,
7529 const ctl_t * ctl,
7530 const met_t * met,
7531 float var[EX][EY][EP],
7532 const char *varname,
7533 const float bound_min,
7534 const float bound_max);
7535
7563void read_met_cape(
7564 const ctl_t * ctl,
7565 const clim_t * clim,
7566 met_t * met);
7567
7590void read_met_cloud(
7591 met_t * met);
7592
7618void read_met_detrend(
7619 const ctl_t * ctl,
7620 met_t * met);
7621
7645 met_t * met);
7646
7673void read_met_geopot(
7674 const ctl_t * ctl,
7675 met_t * met);
7676
7704 const char *filename,
7705 const ctl_t * ctl,
7706 met_t * met);
7707
7729 codes_handle ** handles,
7730 int count_handles,
7731 met_t * met);
7732
7753 codes_handle ** handles,
7754 const int num_messages,
7755 const ctl_t * ctl,
7756 met_t * met);
7757
7788 codes_handle ** handles,
7789 const int num_messages,
7790 const ctl_t * ctl,
7791 met_t * met);
7792
7821void read_met_ml2pl(
7822 const ctl_t * ctl,
7823 const met_t * met,
7824 float var[EX][EY][EP],
7825 const char *varname);
7826
7849 const ctl_t * ctl,
7850 met_t * met);
7851
7882int read_met_nc(
7883 const char *filename,
7884 const ctl_t * ctl,
7885 met_t * met,
7886 dd_t * dd);
7887
7922void read_met_nc_grid(
7923 const char *filename,
7924 const int ncid,
7925 const ctl_t * ctl,
7926 met_t * met,
7927 dd_t * dd);
7928
7974 dd_t * dd,
7975 const ctl_t * ctl,
7976 met_t * met,
7977 const int ncid);
7978
8011 const int ncid,
8012 const ctl_t * ctl,
8013 met_t * met,
8014 dd_t * dd);
8015
8046 const int ncid,
8047 const ctl_t * ctl,
8048 met_t * met,
8049 dd_t * dd);
8050
8083int read_met_nc_2d(
8084 const int ncid,
8085 const char *varname,
8086 const char *varname2,
8087 const char *varname3,
8088 const char *varname4,
8089 const char *varname5,
8090 const char *varname6,
8091 const ctl_t * ctl,
8092 const met_t * met,
8093 dd_t * dd,
8094 float dest[EX][EY],
8095 const float scl,
8096 const int init);
8097
8127int read_met_nc_3d(
8128 const int ncid,
8129 const char *varname,
8130 const char *varname2,
8131 const char *varname3,
8132 const char *varname4,
8133 const ctl_t * ctl,
8134 const met_t * met,
8135 dd_t * dd,
8136 float dest[EX][EY][EP],
8137 const float scl);
8138
8184void read_met_pbl(
8185 const ctl_t * ctl,
8186 met_t * met);
8187
8220 met_t * met);
8221
8252 met_t * met);
8253
8284void read_met_pv(
8285 met_t * met);
8286
8309void read_met_ozone(
8310 met_t * met);
8311
8340void read_met_sample(
8341 const ctl_t * ctl,
8342 met_t * met);
8343
8372void read_met_tropo(
8373 const ctl_t * ctl,
8374 const clim_t * clim,
8375 met_t * met);
8376
8408void read_obs(
8409 const char *filename,
8410 const ctl_t * ctl,
8411 double *rt,
8412 double *rz,
8413 double *rlon,
8414 double *rlat,
8415 double *robs,
8416 int *nobs);
8417
8445void read_obs_asc(
8446 const char *filename,
8447 double *rt,
8448 double *rz,
8449 double *rlon,
8450 double *rlat,
8451 double *robs,
8452 int *nobs);
8453
8480void read_obs_nc(
8481 const char *filename,
8482 double *rt,
8483 double *rz,
8484 double *rlon,
8485 double *rlat,
8486 double *robs,
8487 int *nobs);
8488
8522double scan_ctl(
8523 const char *filename,
8524 int argc,
8525 char *argv[],
8526 const char *varname,
8527 const int arridx,
8528 const char *defvalue,
8529 char *value);
8530
8557double sedi(
8558 const double p,
8559 const double T,
8560 const double rp,
8561 const double rhop);
8562
8592void spline(
8593 const double *x,
8594 const double *y,
8595 const int n,
8596 const double *x2,
8597 double *y2,
8598 const int n2,
8599 const int method);
8600
8623float stddev(
8624 const float *data,
8625 const int n);
8626
8651void time2jsec(
8652 const int year,
8653 const int mon,
8654 const int day,
8655 const int hour,
8656 const int min,
8657 const int sec,
8658 const double remain,
8659 double *jsec);
8660
8689void timer(
8690 const char *name,
8691 const char *group,
8692 const int output);
8693
8719double time_from_filename(
8720 const char *filename,
8721 const int offset);
8722
8740double tropo_weight(
8741 const clim_t * clim,
8742 const atm_t * atm,
8743 const int ip);
8744
8767void write_atm_asc(
8768 const char *filename,
8769 const ctl_t * ctl,
8770 const atm_t * atm,
8771 const double t);
8772
8796void write_atm_bin(
8797 const char *filename,
8798 const ctl_t * ctl,
8799 const atm_t * atm);
8800
8824void write_atm_clams(
8825 const char *filename,
8826 const ctl_t * ctl,
8827 const atm_t * atm);
8828
8854 const char *dirname,
8855 const ctl_t * ctl,
8856 const atm_t * atm,
8857 const double t);
8858
8882void write_atm_nc(
8883 const char *filename,
8884 const ctl_t * ctl,
8885 const atm_t * atm);
8886
8915void write_csi(
8916 const char *filename,
8917 const ctl_t * ctl,
8918 const atm_t * atm,
8919 const double t);
8920
8962 const char *filename,
8963 const ctl_t * ctl,
8964 const atm_t * atm,
8965 const double t);
8966
8994void write_ens(
8995 const char *filename,
8996 const ctl_t * ctl,
8997 const atm_t * atm,
8998 const double t);
8999
9038void write_grid(
9039 const char *filename,
9040 const ctl_t * ctl,
9041 met_t * met0,
9042 met_t * met1,
9043 const atm_t * atm,
9044 const double t);
9045
9091void write_grid_asc(
9092 const char *filename,
9093 const ctl_t * ctl,
9094 const double *cd,
9095 double *mean[NQ],
9096 double *sigma[NQ],
9097 const double *vmr_impl,
9098 const double t,
9099 const double *z,
9100 const double *lon,
9101 const double *lat,
9102 const double *area,
9103 const double dz,
9104 const int *np);
9105
9148void write_grid_nc(
9149 const char *filename,
9150 const ctl_t * ctl,
9151 const double *cd,
9152 double *mean[NQ],
9153 double *sigma[NQ],
9154 const double *vmr_impl,
9155 const double t,
9156 const double *z,
9157 const double *lon,
9158 const double *lat,
9159 const double *area,
9160 const double dz,
9161 const int *np);
9162
9192void write_met_bin(
9193 const char *filename,
9194 const ctl_t * ctl,
9195 met_t * met);
9196
9224void write_met_bin_2d(
9225 FILE * out,
9226 met_t * met,
9227 float var[EX][EY],
9228 const char *varname);
9229
9267void write_met_bin_3d(
9268 FILE * out,
9269 const ctl_t * ctl,
9270 met_t * met,
9271 float var[EX][EY][EP],
9272 const char *varname,
9273 const int precision,
9274 const double tolerance);
9275
9303void write_met_nc(
9304 const char *filename,
9305 const ctl_t * ctl,
9306 met_t * met);
9307
9330void write_met_nc_2d(
9331 const int ncid,
9332 const char *varname,
9333 met_t * met,
9334 float var[EX][EY],
9335 const float scl);
9336
9360void write_met_nc_3d(
9361 const int ncid,
9362 const char *varname,
9363 met_t * met,
9364 float var[EX][EY][EP],
9365 const float scl);
9366
9397void write_prof(
9398 const char *filename,
9399 const ctl_t * ctl,
9400 met_t * met0,
9401 met_t * met1,
9402 const atm_t * atm,
9403 const double t);
9404
9436void write_sample(
9437 const char *filename,
9438 const ctl_t * ctl,
9439 met_t * met0,
9440 met_t * met1,
9441 const atm_t * atm,
9442 const double t);
9443
9470void write_station(
9471 const char *filename,
9472 const ctl_t * ctl,
9473 atm_t * atm,
9474 const double t);
9475
9504void write_vtk(
9505 const char *filename,
9506 const ctl_t * ctl,
9507 const atm_t * atm,
9508 const double t);
9509
9510#endif /* LIBTRAC_H */
void read_met_geopot(const ctl_t *ctl, met_t *met)
Calculates geopotential heights from meteorological data.
Definition: mptrac.c:8184
int dd_init(ctl_t *ctl, dd_t *dd, atm_t *atm)
Initializes domain decomposition for parallel processing.
void compress_zstd(const char *varname, float *array, const size_t n, const int decompress, const int level, FILE *inout)
Compresses or decompresses a float array using Zstandard (ZSTD).
#define LEN
Maximum length of ASCII data lines.
Definition: mptrac.h:348
void mptrac_alloc(ctl_t **ctl, cache_t **cache, clim_t **clim, met_t **met0, met_t **met1, atm_t **atm, dd_t **dd)
Allocates and initializes memory resources for MPTRAC.
Definition: mptrac.c:5203
void mptrac_write_atm(const char *filename, const ctl_t *ctl, const atm_t *atm, const double t)
Writes air parcel data to a file in various formats.
Definition: mptrac.c:6840
void day2doy(const int year, const int mon, const int day, int *doy)
Get day of year from date.
Definition: mptrac.c:1052
void read_met_grib_surface(codes_handle **handles, const int num_messages, const ctl_t *ctl, met_t *met)
Reads surface meteorological data from a grib file and stores it in the meteorological data structure...
void read_met_extrapolate(met_t *met)
Extrapolates meteorological data.
Definition: mptrac.c:8144
void write_atm_clams_traj(const char *dirname, const ctl_t *ctl, const atm_t *atm, const double t)
Writes CLaMS trajectory data to a NetCDF file.
Definition: mptrac.c:11517
int read_met_nc_2d(const int ncid, const char *varname, const char *varname2, const char *varname3, const char *varname4, const char *varname5, const char *varname6, const ctl_t *ctl, const met_t *met, dd_t *dd, float dest[EX][EY], const float scl, const int init)
Reads a 2-dimensional meteorological variable from a NetCDF file.
Definition: mptrac.c:8814
void write_met_nc_2d(const int ncid, const char *varname, met_t *met, float var[EX][EY], const float scl)
Writes a 2D meteorological variable to a NetCDF file.
Definition: mptrac.c:12916
void read_met_sample(const ctl_t *ctl, met_t *met)
Downsamples meteorological data based on specified parameters.
Definition: mptrac.c:10556
void write_met_bin_3d(FILE *out, const ctl_t *ctl, met_t *met, float var[EX][EY][EP], const char *varname, const int precision, const double tolerance)
Writes a 3-dimensional meteorological variable to a binary file.
Definition: mptrac.c:12662
void read_obs(const char *filename, const ctl_t *ctl, double *rt, double *rz, double *rlon, double *rlat, double *robs, int *nobs)
Reads observation data from a file and stores it in arrays.
Definition: mptrac.c:10901
void module_advect(const ctl_t *ctl, const cache_t *cache, met_t *met0, met_t *met1, atm_t *atm)
Advances particle positions using different advection schemes.
Definition: mptrac.c:2960
void module_timesteps(const ctl_t *ctl, cache_t *cache, met_t *met0, atm_t *atm, const double t)
Calculate time steps for air parcels based on specified conditions.
Definition: mptrac.c:4919
void module_meteo(const ctl_t *ctl, const cache_t *cache, const clim_t *clim, met_t *met0, met_t *met1, atm_t *atm)
Update atmospheric properties using meteorological data.
Definition: mptrac.c:4177
void read_clim_photo(const char *filename, clim_photo_t *photo)
Reads photolysis rates from a NetCDF file and populates the given photolysis structure.
Definition: mptrac.c:7272
void read_met_cloud(met_t *met)
Calculates cloud-related variables for each grid point.
Definition: mptrac.c:7983
void module_decay(const ctl_t *ctl, const cache_t *cache, const clim_t *clim, atm_t *atm)
Simulate exponential decay processes for atmospheric particles.
Definition: mptrac.c:3563
double sedi(const double p, const double T, const double rp, const double rhop)
Calculates the sedimentation velocity of a particle in air.
Definition: mptrac.c:11074
double cos_sza(const double sec, const double lon, const double lat)
Calculates the cosine of the solar zenith angle.
Definition: mptrac.c:1011
#define METVAR
Number of 3-D meteorological variables.
Definition: mptrac.h:353
void intpol_met_space_2d(const met_t *met, float array[EX][EY], const double lon, const double lat, double *var, int *ci, double *cw, const int init)
Interpolates meteorological variables in 2D space.
Definition: mptrac.c:2359
int read_met_nc_3d(const int ncid, const char *varname, const char *varname2, const char *varname3, const char *varname4, const ctl_t *ctl, const met_t *met, dd_t *dd, float dest[EX][EY][EP], const float scl)
Reads a 3-dimensional meteorological variable from a NetCDF file.
Definition: mptrac.c:9076
void dd_sort(const ctl_t *ctl, met_t *met0, atm_t *atm, dd_t *dd, int *nparticles, int *rank)
Sort particles according to box index and target rank for neighbours.
int read_atm_nc(const char *filename, const ctl_t *ctl, atm_t *atm)
Reads air parcel data from a generic netCDF file and populates the given atmospheric structure.
Definition: mptrac.c:7239
void write_csi_ens(const char *filename, const ctl_t *ctl, const atm_t *atm, const double t)
Writes ensemble-based Critical Success Index (CSI) and other verification statistics to an output fil...
void read_met_pbl(const ctl_t *ctl, met_t *met)
Computes the planetary boundary layer (PBL) pressure based on meteorological data.
Definition: mptrac.c:10164
void read_met_detrend(const ctl_t *ctl, met_t *met)
Detrends meteorological data.
Definition: mptrac.c:8040
void read_met_tropo(const ctl_t *ctl, const clim_t *clim, met_t *met)
Calculates the tropopause and related meteorological variables based on various methods and stores th...
Definition: mptrac.c:10729
void read_obs_asc(const char *filename, double *rt, double *rz, double *rlon, double *rlat, double *robs, int *nobs)
Reads observation data from an ASCII file.
Definition: mptrac.c:10945
void module_chem_init(const ctl_t *ctl, const cache_t *cache, const clim_t *clim, met_t *met0, met_t *met1, atm_t *atm)
Initializes the chemistry modules by setting atmospheric composition.
Definition: mptrac.c:3408
int locate_reg(const double *xx, const int n, const double x)
Locate the index of the interval containing a given value in a regular grid.
Definition: mptrac.c:2921
void get_tropo(const int met_tropo, ctl_t *ctl, clim_t *clim, met_t *met, const double *lons, const int nx, const double *lats, const int ny, double *pt, double *zt, double *tt, double *qt, double *o3t, double *ps, double *zs)
Calculate tropopause data.
Definition: mptrac.c:2059
void locate_vert(float profiles[EX][EY][EP], const int np, const int lon_ap_ind, const int lat_ap_ind, const double alt_ap, int *ind)
Locate the four vertical indizes of a box for a given height value.
Definition: mptrac.c:2940
void compress_cms(const ctl_t *ctl, const char *varname, float *array, const size_t nx, const size_t ny, const size_t np, const double *plev, const int decompress, FILE *inout)
Compresses or decompresses a 3-D meteorological field using cmultiscale.
void read_met_nc_levels(const int ncid, const ctl_t *ctl, met_t *met, dd_t *dd)
Reads and processes meteorological level data from NetCDF files with domain decomposition.
Definition: mptrac.c:8617
void read_met_monotonize(const ctl_t *ctl, met_t *met)
Makes zeta and pressure profiles monotone.
Definition: mptrac.c:9810
int read_clim_ts(const char *filename, clim_ts_t *ts)
Reads a climatological time series from a file and populates the given time series structure.
Definition: mptrac.c:7391
void read_met_periodic(met_t *met)
Applies periodic boundary conditions to meteorological data along longitudinal axis.
Definition: mptrac.c:10301
void mptrac_run_timestep(ctl_t *ctl, cache_t *cache, clim_t *clim, met_t **met0, met_t **met1, atm_t *atm, double t, dd_t *dd)
Executes a single timestep of the MPTRAC model simulation.
Definition: mptrac.c:6576
void module_timesteps_init(ctl_t *ctl, const atm_t *atm)
Initialize start time and time interval for time-stepping.
Definition: mptrac.c:4966
void write_ens(const char *filename, const ctl_t *ctl, const atm_t *atm, const double t)
Writes ensemble data to a file.
Definition: mptrac.c:11999
void module_mixing(const ctl_t *ctl, const clim_t *clim, atm_t *atm, const double t)
Update atmospheric properties through interparcel mixing.
Definition: mptrac.c:4284
void compress_zfp(const char *varname, float *array, const int nx, const int ny, const int nz, const int precision, const double tolerance, const int decompress, FILE *inout)
Compresses or decompresses a 3D array of floats using the ZFP library.
double clim_zm(const clim_zm_t *zm, const double t, const double lat, const double p)
Interpolates monthly mean zonal mean climatological variables.
Definition: mptrac.c:405
#define EY
Maximum number of latitudes for meteo data.
Definition: mptrac.h:343
void module_mixing_help(const ctl_t *ctl, const clim_t *clim, atm_t *atm, const int *ixs, const int *iys, const int *izs, const int qnt_idx, const int use_ensemble)
Perform subgrid-scale interparcel mixing of a given quantity.
Definition: mptrac.c:4356
void read_clim_photo_help(const int ncid, const char *varname, const clim_photo_t *photo, double var[CP][CSZA][CO3])
Reads a 3D climatological photochemistry variable from a NetCDF file.
Definition: mptrac.c:7363
void read_met_ml2pl(const ctl_t *ctl, const met_t *met, float var[EX][EY][EP], const char *varname)
Interpolates meteorological data to specified pressure levels.
Definition: mptrac.c:9768
double clim_tropo(const clim_t *clim, const double t, const double lat)
Calculates the tropopause pressure based on climatological data.
Definition: mptrac.c:204
void read_obs_nc(const char *filename, double *rt, double *rz, double *rlon, double *rlat, double *robs, int *nobs)
Reads observation data from a NetCDF file.
Definition: mptrac.c:10973
void dd_particles2atm(atm_t *atm, particle_t *particles, ctl_t *ctl, int *nparticles, cache_t *cache)
Converts particle data to atmospheric data.
void read_met_bin_2d(FILE *in, const met_t *met, float var[EX][EY], const char *varname)
Reads a 2-dimensional meteorological variable from a binary file and stores it in the provided array.
Definition: mptrac.c:7733
int locate_irr(const double *xx, const int n, const double x)
Locate the index of the interval containing a given value in a sorted array.
Definition: mptrac.c:2857
void module_isosurf_init(const ctl_t *ctl, cache_t *cache, met_t *met0, met_t *met1, atm_t *atm)
Initialize the isosurface module based on atmospheric data.
Definition: mptrac.c:4001
void level_definitions(ctl_t *ctl)
Defines pressure levels for meteorological data.
Definition: mptrac.c:2647
void intpol_met_4d_zeta(const met_t *met0, float height0[EX][EY][EP], float array0[EX][EY][EP], const met_t *met1, float height1[EX][EY][EP], float array1[EX][EY][EP], const double ts, const double height, const double lon, const double lat, double *var, int *ci, double *cw, const int init)
Interpolates meteorological variables to a given position and time.
Definition: mptrac.c:2129
void write_grid_asc(const char *filename, const ctl_t *ctl, const double *cd, double *mean[NQ], double *sigma[NQ], const double *vmr_impl, const double t, const double *z, const double *lon, const double *lat, const double *area, const double dz, const int *np)
Writes grid data to an ASCII file.
Definition: mptrac.c:12287
void mptrac_update_device(const ctl_t *ctl, const cache_t *cache, const clim_t *clim, met_t **met0, met_t **met1, const atm_t *atm)
Updates device memory for specified data structures.
Definition: mptrac.c:6728
void time2jsec(const int year, const int mon, const int day, const int hour, const int min, const int sec, const double remain, double *jsec)
Converts time components to seconds since January 1, 2000, 12:00:00 UTC.
Definition: mptrac.c:11175
void mptrac_init(ctl_t *ctl, cache_t *cache, clim_t *clim, atm_t *atm, const int ntask)
Initializes the MPTRAC model and its associated components.
Definition: mptrac.c:5411
void intpol_met_time_3d(const met_t *met0, float array0[EX][EY][EP], const met_t *met1, float array1[EX][EY][EP], const double ts, const double p, const double lon, const double lat, double *var, int *ci, double *cw, const int init)
Interpolates meteorological data in 3D space and time.
Definition: mptrac.c:2417
#define codes_handle
Placeholder when ECCODES is not available.
Definition: mptrac.h:240
void fft_help(double *fcReal, double *fcImag, const int n)
Computes the Fast Fourier Transform (FFT) of a complex sequence.
Definition: mptrac.c:1911
void module_wet_depo(const ctl_t *ctl, const cache_t *cache, met_t *met0, met_t *met1, atm_t *atm)
Perform wet deposition calculations for air parcels.
Definition: mptrac.c:5068
void compress_pck(const char *varname, float *array, const size_t nxy, const size_t nz, const int decompress, FILE *inout)
Compresses or decompresses a 3D array of floats.
Definition: mptrac.c:680
double nat_temperature(const double p, const double h2o, const double hno3)
Calculates the nitric acid trihydrate (NAT) temperature.
Definition: mptrac.c:7035
void spline(const double *x, const double *y, const int n, const double *x2, double *y2, const int n2, const int method)
Performs spline interpolation or linear interpolation.
Definition: mptrac.c:11107
#define CP
Maximum number of pressure levels for climatological data.
Definition: mptrac.h:398
#define NQ
Maximum number of quantities per data point.
Definition: mptrac.h:363
double clim_photo(const double rate[CP][CSZA][CO3], const clim_photo_t *photo, const double p, const double sza, const double o3c)
Calculates the photolysis rate for a given set of atmospheric conditions.
Definition: mptrac.c:147
void read_clim_zm(const char *filename, const char *varname, clim_zm_t *zm)
Reads zonally averaged climatological data from a netCDF file and populates the given structure.
Definition: mptrac.c:7445
#define EX
Maximum number of longitudes for meteo data.
Definition: mptrac.h:338
void read_met_nc_grid_dd_naive(dd_t *dd, const ctl_t *ctl, met_t *met, const int ncid)
Read meteorological grid data from a NetCDF file and set up subdomain decomposition with halos.
Definition: mptrac.c:9935
void module_sedi(const ctl_t *ctl, const cache_t *cache, met_t *met0, met_t *met1, atm_t *atm)
Simulate sedimentation of particles in the atmosphere.
Definition: mptrac.c:4787
int read_met_nc(const char *filename, const ctl_t *ctl, met_t *met, dd_t *dd)
Reads meteorological data from a NetCDF file and processes it.
Definition: mptrac.c:9895
void timer(const char *name, const char *group, const int output)
Measures and reports elapsed time for named and grouped timers.
Definition: mptrac.c:11206
void write_atm_asc(const char *filename, const ctl_t *ctl, const atm_t *atm, const double t)
Writes air parcel data to an ASCII file or gnuplot.
Definition: mptrac.c:11332
void intpol_met_space_3d(const met_t *met, float array[EX][EY][EP], const double p, const double lon, const double lat, double *var, int *ci, double *cw, const int init)
Interpolates meteorological variables in 3D space.
Definition: mptrac.c:2301
void module_convection(const ctl_t *ctl, cache_t *cache, met_t *met0, met_t *met1, atm_t *atm)
Performs convective mixing of atmospheric particles.
Definition: mptrac.c:3450
void read_kernel(const char *filename, double kz[EP], double kw[EP], int *nk)
Reads kernel function data from a file and populates the provided arrays.
Definition: mptrac.c:7544
void module_bound_cond(const ctl_t *ctl, const cache_t *cache, const clim_t *clim, met_t *met0, met_t *met1, atm_t *atm)
Apply boundary conditions to particles based on meteorological and climatological data.
Definition: mptrac.c:3150
double scan_ctl(const char *filename, int argc, char *argv[], const char *varname, const int arridx, const char *defvalue, char *value)
Scans a control file or command-line arguments for a specified variable.
Definition: mptrac.c:11002
#define CT
Maximum number of time steps for climatological data.
Definition: mptrac.h:408
void module_advect_init(const ctl_t *ctl, const cache_t *cache, met_t *met0, met_t *met1, atm_t *atm)
Initializes the advection module by setting up pressure fields.
Definition: mptrac.c:3123
int dd_is_periodic_longitude(met_t *met, int nx_glob)
Check whether the longitude grid is periodic (global coverage).
void module_radio_decay(const ctl_t *ctl, const cache_t *cache, atm_t *atm)
Apply radioactive decay to atmospheric tracer species.
Definition: mptrac.c:4593
int read_met_grib(const char *filename, const ctl_t *ctl, met_t *met)
Reads meteorological data from a grib file and processes it.
void mptrac_get_met(ctl_t *ctl, clim_t *clim, const double t, met_t **met0, met_t **met1, dd_t *dd)
Retrieves meteorological data for the specified time.
Definition: mptrac.c:5289
void module_sort_help(double *a, const int *p, const int np)
Reorder an array based on a given permutation.
Definition: mptrac.c:4881
float stddev(const float *data, const int n)
Calculates the standard deviation of a set of data.
Definition: mptrac.c:11154
void intpol_tropo_3d(const double time0, float array0[EX][EY], const double time1, float array1[EX][EY], const double lons[EX], const double lats[EY], const int nlon, const int nlat, const double time, const double lon, const double lat, const int method, double *var, double *sigma)
Interpolates tropopause data in 3D (latitude, longitude, and time).
Definition: mptrac.c:2479
void read_met_bin_3d(FILE *in, const ctl_t *ctl, const met_t *met, float var[EX][EY][EP], const char *varname, const float bound_min, const float bound_max)
Reads 3D meteorological data from a binary file, potentially using different compression methods.
Definition: mptrac.c:7762
int locate_irr_float(const float *xx, const int n, const double x, const int ig)
Locate the index of the interval containing a given value in an irregularly spaced array.
Definition: mptrac.c:2887
void dd_get_rect_neighbour(const ctl_t ctl, dd_t *dd)
Determines rectangular neighbouring ranks for MPI processes.
double time_from_filename(const char *filename, const int offset)
Extracts and converts a timestamp from a filename to Julian seconds.
Definition: mptrac.c:11274
void write_prof(const char *filename, const ctl_t *ctl, met_t *met0, met_t *met1, const atm_t *atm, const double t)
Writes profile data to a specified file.
Definition: mptrac.c:12975
void dd_atm2particles(atm_t *atm, particle_t *particles, ctl_t *ctl, int *nparticles, cache_t *cache, int rank)
Extracts particles from an atmospheric state and prepares them for inter-domain transfer.
void mptrac_read_clim(const ctl_t *ctl, clim_t *clim)
Reads various climatological data and populates the given climatology structure.
Definition: mptrac.c:5501
void module_chem_grid(const ctl_t *ctl, met_t *met0, met_t *met1, atm_t *atm, const double t)
Computes gridded chemical tracer concentrations (volume mixing ratio) from individual air parcel mass...
Definition: mptrac.c:3246
double tropo_weight(const clim_t *clim, const atm_t *atm, const int ip)
Computes a weighting factor based on tropopause pressure.
Definition: mptrac.c:11309
void write_met_nc(const char *filename, const ctl_t *ctl, met_t *met)
Writes meteorological data to a NetCDF file.
Definition: mptrac.c:12752
void module_rng_init(const int ntask)
Initialize random number generators for parallel tasks.
Definition: mptrac.c:4651
int mptrac_read_atm(const char *filename, const ctl_t *ctl, atm_t *atm)
Reads air parcel data from a specified file into the given atmospheric structure.
Definition: mptrac.c:5430
void mptrac_update_host(const ctl_t *ctl, const cache_t *cache, const clim_t *clim, met_t **met0, met_t **met1, const atm_t *atm)
Updates host memory for specified data structures.
Definition: mptrac.c:6784
void mptrac_write_output(const char *dirname, const ctl_t *ctl, met_t *met0, met_t *met1, atm_t *atm, const double t)
Writes various types of output data to files in a specified directory.
Definition: mptrac.c:6944
double clim_oh(const ctl_t *ctl, const clim_t *clim, const double t, const double lon, const double lat, const double p)
Calculates the hydroxyl radical (OH) concentration from climatology data, with an optional diurnal co...
Definition: mptrac.c:89
void read_met_grib_levels(codes_handle **handles, const int num_messages, const ctl_t *ctl, met_t *met)
Reads meteorological variables at different vertical levels from a grib file.
#define CTS
Maximum number of data points of climatological time series.
Definition: mptrac.h:413
void read_met_ozone(met_t *met)
Calculates the total column ozone from meteorological ozone data.
Definition: mptrac.c:10527
int dd_calc_subdomain_from_coords(double lon, double lat, met_t *met, ctl_t *ctl, int mpi_size, int nx_glob, int ny_glob)
Computes the destination subdomain (MPI rank) for a particle based on its geographic coordinates.
void broadcast_large_data(void *data, size_t N)
Broadcasts large data across all processes in an MPI communicator.
void read_met_nc_surface(const int ncid, const ctl_t *ctl, met_t *met, dd_t *dd)
Reads and processes surface meteorological data from NetCDF files with domain decomposition.
Definition: mptrac.c:8479
void clim_tropo_init(clim_t *clim)
Initializes the tropopause data in the climatology structure.
Definition: mptrac.c:232
void module_rng(const ctl_t *ctl, double *rs, const size_t n, const int method)
Generate random numbers using various methods and distributions.
Definition: mptrac.c:4682
void get_met_help(const ctl_t *ctl, const double t, const int direct, const char *metbase, const double dt_met, char *filename)
Generates a formatted filename for meteorological data files based on the input parameters.
Definition: mptrac.c:1968
void write_station(const char *filename, const ctl_t *ctl, atm_t *atm, const double t)
Writes station data to a specified file.
Definition: mptrac.c:13364
void cart2geo(const double *x, double *z, double *lon, double *lat)
Converts Cartesian coordinates to geographic coordinates.
Definition: mptrac.c:74
#define NP
Maximum number of atmospheric data points.
Definition: mptrac.h:358
void doy2day(const int year, const int doy, int *mon, int *day)
Converts a given day of the year (DOY) to a date (month and day).
Definition: mptrac.c:1881
void dd_communicate_particles(particle_t *particles, int *nparticles, MPI_Datatype MPI_Particle, int *neighbours, int nneighbours, ctl_t ctl)
Communicates particles between MPI processes.
void intpol_met_time_2d(const met_t *met0, float array0[EX][EY], const met_t *met1, float array1[EX][EY], const double ts, const double lon, const double lat, double *var, int *ci, double *cw, const int init)
Interpolates meteorological data in 2D space and time.
Definition: mptrac.c:2446
void module_position(const cache_t *cache, met_t *met0, met_t *met1, atm_t *atm)
Update the positions and pressure levels of atmospheric particles.
Definition: mptrac.c:4542
void clim_oh_diurnal_correction(const ctl_t *ctl, clim_t *clim)
Applies a diurnal correction to the hydroxyl radical (OH) concentration in climatology data.
Definition: mptrac.c:115
void write_met_bin_2d(FILE *out, met_t *met, float var[EX][EY], const char *varname)
Writes a 2-dimensional meteorological variable to a binary file.
Definition: mptrac.c:12633
void read_met_pv(met_t *met)
Calculates potential vorticity (PV) from meteorological data.
Definition: mptrac.c:10421
int read_atm_bin(const char *filename, const ctl_t *ctl, atm_t *atm)
Reads air parcel data from a binary file and populates the given atmospheric structure.
Definition: mptrac.c:7123
void get_met_replace(char *orig, char *search, char *repl)
Replaces occurrences of a substring in a string with another substring.
Definition: mptrac.c:2035
#define CY
Maximum number of latitudes for climatological data.
Definition: mptrac.h:388
void module_diff_meso(const ctl_t *ctl, cache_t *cache, met_t *met0, met_t *met1, atm_t *atm)
Simulate mesoscale diffusion for atmospheric particles.
Definition: mptrac.c:3602
void module_dd(ctl_t *ctl, atm_t *atm, cache_t *cache, dd_t *dd, met_t **met)
Manages domain decomposition and particle communication in parallel processing.
void dd_sort_help(double *a, dd_t *dd, const int np)
Reorder an array according to a permutation vector.
void module_sort(const ctl_t *ctl, met_t *met0, atm_t *atm)
Sort particles according to box index.
Definition: mptrac.c:4816
double clim_ts(const clim_ts_t *ts, const double t)
Interpolates a time series of climatological variables.
Definition: mptrac.c:387
void jsec2time(const double jsec, int *year, int *mon, int *day, int *hour, int *min, int *sec, double *remain)
Converts Julian seconds to calendar date and time components.
Definition: mptrac.c:2570
int read_met_bin(const char *filename, const ctl_t *ctl, met_t *met)
Reads meteorological data from a binary file.
Definition: mptrac.c:7585
void write_atm_clams(const char *filename, const ctl_t *ctl, const atm_t *atm)
Writes air parcel data to a NetCDF file in the CLaMS format.
Definition: mptrac.c:11464
void module_diff_turb(const ctl_t *ctl, cache_t *cache, const clim_t *clim, met_t *met0, met_t *met1, atm_t *atm)
Applies turbulent diffusion processes to atmospheric particles.
Definition: mptrac.c:3804
int read_atm_clams(const char *filename, const ctl_t *ctl, atm_t *atm)
Reads atmospheric data from a CLAMS NetCDF file.
Definition: mptrac.c:7179
void dd_assign_rect_subdomains_atm(atm_t *atm, ctl_t *ctl, dd_t *dd, int init)
Assign atmospheric particles to rectangular subdomains.
int mptrac_read_met(const char *filename, const ctl_t *ctl, const clim_t *clim, met_t *met, dd_t *dd)
Reads meteorological data from a file, supporting multiple formats and MPI broadcasting.
Definition: mptrac.c:6467
void write_vtk(const char *filename, const ctl_t *ctl, const atm_t *atm, const double t)
Writes VTK (Visualization Toolkit) data to a specified file.
Definition: mptrac.c:13450
#define MPI_Datatype
Placeholder when MPI is not available.
Definition: mptrac.h:195
#define EP
Maximum number of pressure levels for meteo data.
Definition: mptrac.h:333
void module_tracer_chem(const ctl_t *ctl, const cache_t *cache, const clim_t *clim, met_t *met0, met_t *met1, atm_t *atm)
Simulate chemical reactions involving long-lived atmospheric tracers.
Definition: mptrac.c:4997
void mptrac_read_ctl(const char *filename, int argc, char *argv[], ctl_t *ctl)
Reads control parameters from a configuration file and populates the given structure.
Definition: mptrac.c:5561
void mptrac_free(ctl_t *ctl, cache_t *cache, clim_t *clim, met_t *met0, met_t *met1, atm_t *atm, dd_t *dd)
Frees memory resources allocated for MPTRAC.
Definition: mptrac.c:5256
void read_met_polar_winds(met_t *met)
Applies a fix for polar winds in meteorological data.
Definition: mptrac.c:10362
void module_h2o2_chem(const ctl_t *ctl, const cache_t *cache, const clim_t *clim, met_t *met0, met_t *met1, atm_t *atm)
Perform chemical reactions involving H2O2 within cloud particles.
Definition: mptrac.c:3919
void write_grid_nc(const char *filename, const ctl_t *ctl, const double *cd, double *mean[NQ], double *sigma[NQ], const double *vmr_impl, const double t, const double *z, const double *lon, const double *lat, const double *area, const double dz, const int *np)
Writes grid data to a NetCDF file.
Definition: mptrac.c:12391
double pbl_weight(const ctl_t *ctl, const atm_t *atm, const int ip, const double pbl, const double ps)
Computes a weighting factor based on planetary boundary layer pressure.
Definition: mptrac.c:7059
void module_diff_pbl(const ctl_t *ctl, cache_t *cache, met_t *met0, met_t *met1, atm_t *atm)
Computes particle diffusion within the planetary boundary layer (PBL).
Definition: mptrac.c:3679
void write_met_nc_3d(const int ncid, const char *varname, met_t *met, float var[EX][EY][EP], const float scl)
Writes a 3D meteorological variable to a NetCDF file.
Definition: mptrac.c:12945
void module_isosurf(const ctl_t *ctl, const cache_t *cache, met_t *met0, met_t *met1, atm_t *atm)
Apply the isosurface module to adjust atmospheric properties.
Definition: mptrac.c:4071
void module_kpp_chem(ctl_t *ctl, cache_t *cache, clim_t *clim, met_t *met0, met_t *met1, atm_t *atm)
KPP chemistry module.
#define CO3
Maximum number of total column ozone data for climatological data.
Definition: mptrac.h:393
void module_oh_chem(const ctl_t *ctl, const cache_t *cache, const clim_t *clim, met_t *met0, met_t *met1, atm_t *atm)
Perform hydroxyl chemistry calculations for atmospheric particles.
Definition: mptrac.c:4458
void read_met_grib_grid(codes_handle **handles, int count_handles, met_t *met)
Reads global meteorological information from a grib file.
void geo2cart(const double z, const double lon, const double lat, double *x)
Converts geographic coordinates (longitude, latitude, altitude) to Cartesian coordinates.
Definition: mptrac.c:1950
void read_met_nc_grid(const char *filename, const int ncid, const ctl_t *ctl, met_t *met, dd_t *dd)
Reads meteorological grid data from NetCDF files with domain decomposition.
Definition: mptrac.c:8312
double kernel_weight(const double kz[EP], const double kw[EP], const int nk, const double p)
Calculates the kernel weight based on altitude and given kernel data.
Definition: mptrac.c:2603
void compress_sz3(const char *varname, float *array, const int nx, const int ny, const int nz, const int precision, const double tolerance, const int decompress, FILE *inout)
Compresses or decompresses a 3-D float array using the SZ3 library.
void dd_register_MPI_type_particle(MPI_Datatype *MPI_Particle)
Registers a custom MPI datatype for particle structures.
int read_atm_asc(const char *filename, const ctl_t *ctl, atm_t *atm)
Reads air parcel data from an ASCII file and populates the given atmospheric structure.
Definition: mptrac.c:7081
void intpol_check_lon_lat(const double *lons, const int nlon, const double *lats, const int nlat, const double lon, const double lat, double *lon2, double *lat2)
Adjusts longitude and latitude to ensure they fall within valid bounds.
Definition: mptrac.c:2102
void write_sample(const char *filename, const ctl_t *ctl, met_t *met0, met_t *met1, const atm_t *atm, const double t)
Writes sample data to a specified file.
Definition: mptrac.c:13202
void write_grid(const char *filename, const ctl_t *ctl, met_t *met0, met_t *met1, const atm_t *atm, const double t)
Writes grid data to a file in ASCII or netCDF format.
Definition: mptrac.c:12096
void module_dry_depo(const ctl_t *ctl, const cache_t *cache, met_t *met0, met_t *met1, atm_t *atm)
Simulate dry deposition of atmospheric particles.
Definition: mptrac.c:3856
void write_met_bin(const char *filename, const ctl_t *ctl, met_t *met)
Writes meteorological data in binary format to a specified file.
Definition: mptrac.c:12521
void write_atm_bin(const char *filename, const ctl_t *ctl, const atm_t *atm)
Writes air parcel data to a binary file.
Definition: mptrac.c:11414
void read_met_cape(const ctl_t *ctl, const clim_t *clim, met_t *met)
Calculates Convective Available Potential Energy (CAPE) for each grid point.
Definition: mptrac.c:7868
void mptrac_write_met(const char *filename, const ctl_t *ctl, met_t *met)
Writes meteorological data to a file, supporting multiple formats and compression options.
Definition: mptrac.c:6900
#define CSZA
Maximum number of solar zenith angles for climatological data.
Definition: mptrac.h:403
double lapse_rate(const double t, const double h2o)
Calculates the moist adiabatic lapse rate in Kelvin per kilometer.
Definition: mptrac.c:2629
#define DD_NNMAX
Maximum number of neighbours to communicate with in domain decomposition.
Definition: mptrac.h:423
void write_csi(const char *filename, const ctl_t *ctl, const atm_t *atm, const double t)
Writes Critical Success Index (CSI) data to a file.
Definition: mptrac.c:11724
void write_atm_nc(const char *filename, const ctl_t *ctl, const atm_t *atm)
Writes air parcel data to a NetCDF file.
Definition: mptrac.c:11675
Air parcel data.
Definition: mptrac.h:3257
int np
Number of air parcels.
Definition: mptrac.h:3260
Cache data structure.
Definition: mptrac.h:3312
int iso_n
Isosurface balloon number of data points.
Definition: mptrac.h:3324
Climatological data in the form of photolysis rates.
Definition: mptrac.h:3344
int nsza
Number of solar zenith angles.
Definition: mptrac.h:3350
int np
Number of pressure levels.
Definition: mptrac.h:3347
int no3c
Number of total ozone columns.
Definition: mptrac.h:3353
Climatological data.
Definition: mptrac.h:3452
clim_ts_t ccl2f2
CFC-12 time series.
Definition: mptrac.h:3494
clim_photo_t photo
Photolysis rates.
Definition: mptrac.h:3470
clim_zm_t ho2
HO2 zonal means.
Definition: mptrac.h:3482
clim_zm_t hno3
HNO3 zonal means.
Definition: mptrac.h:3473
int tropo_ntime
Number of tropopause timesteps.
Definition: mptrac.h:3455
clim_ts_t sf6
SF6 time series.
Definition: mptrac.h:3500
clim_ts_t ccl4
CFC-10 time series.
Definition: mptrac.h:3488
clim_ts_t ccl3f
CFC-11 time series.
Definition: mptrac.h:3491
clim_zm_t o1d
O(1D) zonal means.
Definition: mptrac.h:3485
clim_zm_t h2o2
H2O2 zonal means.
Definition: mptrac.h:3479
int tropo_nlat
Number of tropopause latitudes.
Definition: mptrac.h:3458
clim_zm_t oh
OH zonal means.
Definition: mptrac.h:3476
clim_ts_t n2o
N2O time series.
Definition: mptrac.h:3497
Climatological data in the form of time series.
Definition: mptrac.h:3400
int ntime
Number of timesteps.
Definition: mptrac.h:3403
Climatological data in the form of zonal means.
Definition: mptrac.h:3420
int np
Number of pressure levels.
Definition: mptrac.h:3429
int ntime
Number of timesteps.
Definition: mptrac.h:3423
int nlat
Number of latitudes.
Definition: mptrac.h:3426
Control parameters.
Definition: mptrac.h:2203
double grid_z0
Lower altitude of gridded data [km].
Definition: mptrac.h:3118
int qnt_o3
Quantity array index for ozone volume mixing ratio.
Definition: mptrac.h:2315
double csi_lat1
Upper latitude of gridded CSI data [deg].
Definition: mptrac.h:3079
int qnt_Coh
Quantity array index for OH volume mixing ratio (chemistry code).
Definition: mptrac.h:2471
double wet_depo_ic_a
Coefficient A for wet deposition in cloud (exponential form).
Definition: mptrac.h:2967
int met_nc_scale
Check netCDF scaling factors (0=no, 1=yes).
Definition: mptrac.h:2567
int qnt_pel
Quantity array index for pressure at equilibrium level (EL).
Definition: mptrac.h:2348
int csi_nz
Number of altitudes of gridded CSI data.
Definition: mptrac.h:3055
double molmass
Molar mass [g/mol].
Definition: mptrac.h:2826
int qnt_p
Quantity array index for pressure.
Definition: mptrac.h:2294
int qnt_Cccl2f2
Quantity array index for CFC-12 volume mixing ratio (chemistry code).
Definition: mptrac.h:2495
int dd_halos_size
Domain decomposition size of halos given in grid-points.
Definition: mptrac.h:3245
int mixing_nx
Number of longitudes of mixing grid.
Definition: mptrac.h:2889
double chemgrid_z1
Upper altitude of chemistry grid [km].
Definition: mptrac.h:2913
int qnt_m
Quantity array index for mass.
Definition: mptrac.h:2234
int qnt_aoa
Quantity array index for age of air.
Definition: mptrac.h:2504
int qnt_rhop
Quantity array index for particle density.
Definition: mptrac.h:2243
int qnt_swc
Quantity array index for cloud snow water content.
Definition: mptrac.h:2327
double csi_obsmin
Minimum observation index to trigger detection.
Definition: mptrac.h:3049
int qnt_pcb
Quantity array index for cloud bottom pressure.
Definition: mptrac.h:2336
double bound_dzs
Boundary conditions surface layer depth [km].
Definition: mptrac.h:2814
double csi_lon1
Upper longitude of gridded CSI data [deg].
Definition: mptrac.h:3070
int qnt_u
Quantity array index for zonal wind.
Definition: mptrac.h:2303
double stat_lon
Longitude of station [deg].
Definition: mptrac.h:3196
double mixing_trop
Interparcel exchange parameter for mixing in the troposphere.
Definition: mptrac.h:2874
double sort_dt
Time step for sorting of particle data [s].
Definition: mptrac.h:2725
double mixing_z1
Upper altitude of mixing grid [km].
Definition: mptrac.h:2886
double stat_r
Search radius around station [km].
Definition: mptrac.h:3202
double wet_depo_bc_a
Coefficient A for wet deposition below cloud (exponential form).
Definition: mptrac.h:2961
int met_zstd_level
ZSTD compression level (from -5 to 22).
Definition: mptrac.h:2576
int csi_ny
Number of latitudes of gridded CSI data.
Definition: mptrac.h:3073
int vtk_sphere
Spherical projection for VTK data (0=no, 1=yes).
Definition: mptrac.h:3226
double chemgrid_z0
Lower altitude of chemistry grid [km].
Definition: mptrac.h:2910
double met_pbl_min
Minimum depth of planetary boundary layer [km].
Definition: mptrac.h:2693
int qnt_iwc
Quantity array index for cloud ice water content.
Definition: mptrac.h:2324
double chemgrid_lat0
Lower latitude of chemistry grid [deg].
Definition: mptrac.h:2928
double conv_cape
CAPE threshold for convection module [J/kg].
Definition: mptrac.h:2778
int qnt_Co1d
Quantity array index for O(1D) volume mixing ratio (chemistry code).
Definition: mptrac.h:2483
double met_cms_eps_pv
cmultiscale compression epsilon for potential vorticity.
Definition: mptrac.h:2615
int qnt_pw
Quantity array index for partial water vapor pressure.
Definition: mptrac.h:2402
double grid_z1
Upper altitude of gridded data [km].
Definition: mptrac.h:3121
int direction
Direction flag (1=forward calculation, -1=backward calculation).
Definition: mptrac.h:2531
int qnt_Cccl4
Quantity array index for CFC-10 volume mixing ratio (chemistry code).
Definition: mptrac.h:2489
int met_dp
Stride for pressure levels.
Definition: mptrac.h:2645
double met_dt_out
Time step for sampling of meteo data along trajectories [s].
Definition: mptrac.h:2712
int qnt_h2o2
Quantity array index for H2O2 volume mixing ratio (climatology).
Definition: mptrac.h:2366
int qnt_vh
Quantity array index for horizontal wind.
Definition: mptrac.h:2438
int csi_nx
Number of longitudes of gridded CSI data.
Definition: mptrac.h:3064
double csi_lat0
Lower latitude of gridded CSI data [deg].
Definition: mptrac.h:3076
double turb_dz_trop
Vertical turbulent diffusion coefficient (troposphere) [m^2/s].
Definition: mptrac.h:2760
int met_pbl
Planetary boundary layer data (0=file, 1=z2p, 2=Richardson, 3=theta).
Definition: mptrac.h:2690
int qnt_lwc
Quantity array index for cloud liquid water content.
Definition: mptrac.h:2318
double turb_mesoz
Vertical scaling factor for mesoscale wind fluctuations.
Definition: mptrac.h:2769
int grid_nc_level
zlib compression level of netCDF grid data files (0=off).
Definition: mptrac.h:3106
int grid_nx
Number of longitudes of gridded data.
Definition: mptrac.h:3124
int atm_type
Type of atmospheric data files (0=ASCII, 1=binary, 2=netCDF, 3=CLaMS_traj, 4=CLaMS_pos).
Definition: mptrac.h:3020
double bound_mass
Boundary conditions mass per particle [kg].
Definition: mptrac.h:2787
double grid_lat0
Lower latitude of gridded data [deg].
Definition: mptrac.h:3136
int qnt_ts
Quantity array index for surface temperature.
Definition: mptrac.h:2249
int qnt_loss_rate
Quantity array index for total loss rate.
Definition: mptrac.h:2393
double met_cms_eps_h2o
cmultiscale compression epsilon for water vapor.
Definition: mptrac.h:2618
int qnt_plfc
Quantity array index for pressure at level of free convection (LCF).
Definition: mptrac.h:2345
int qnt_Acs137
Quantity array index for radioactive activity of Cs-137.
Definition: mptrac.h:2516
double grid_lon0
Lower longitude of gridded data [deg].
Definition: mptrac.h:3127
int qnt_o1d
Quantity array index for O(1D) volume mixing ratio (climatology).
Definition: mptrac.h:2372
int met_tropo_spline
Tropopause interpolation method (0=linear, 1=spline).
Definition: mptrac.h:2709
int qnt_tvirt
Quantity array index for virtual temperature.
Definition: mptrac.h:2432
double dt_met
Time step of meteo data [s].
Definition: mptrac.h:2550
double chemgrid_lat1
Upper latitude of chemistry grid [deg].
Definition: mptrac.h:2931
int met_geopot_sy
Latitudinal smoothing of geopotential heights.
Definition: mptrac.h:2681
double met_cms_eps_u
cmultiscale compression epsilon for zonal wind.
Definition: mptrac.h:2606
double turb_dx_strat
Horizontal turbulent diffusion coefficient (stratosphere) [m^2/s].
Definition: mptrac.h:2754
int qnt_vmr
Quantity array index for volume mixing ratio.
Definition: mptrac.h:2237
int qnt_lsm
Quantity array index for land-sea mask.
Definition: mptrac.h:2270
int qnt_theta
Quantity array index for potential temperature.
Definition: mptrac.h:2414
double bound_lat1
Boundary conditions maximum longitude [deg].
Definition: mptrac.h:2802
double stat_t1
Stop time for station output [s].
Definition: mptrac.h:3208
double turb_dx_trop
Horizontal turbulent diffusion coefficient (troposphere) [m^2/s].
Definition: mptrac.h:2751
int grid_type
Type of grid data files (0=ASCII, 1=netCDF).
Definition: mptrac.h:3142
double csi_lon0
Lower longitude of gridded CSI data [deg].
Definition: mptrac.h:3067
int qnt_pbl
Quantity array index for boundary layer pressure.
Definition: mptrac.h:2276
int grid_stddev
Include standard deviations in grid output (0=no, 1=yes).
Definition: mptrac.h:3112
int qnt_psice
Quantity array index for saturation pressure over ice.
Definition: mptrac.h:2399
double chemgrid_lon0
Lower longitude of chemistry grid [deg].
Definition: mptrac.h:2919
int bound_pbl
Boundary conditions planetary boundary layer (0=no, 1=yes).
Definition: mptrac.h:2820
int qnt_mloss_wet
Quantity array index for total mass loss due to wet deposition.
Definition: mptrac.h:2384
int radio_decay
Switch for radioactive decay module (0=off, 1=on).
Definition: mptrac.h:2955
int met_geopot_sx
Longitudinal smoothing of geopotential heights.
Definition: mptrac.h:2678
int met_sy
Smoothing for latitudes.
Definition: mptrac.h:2651
int qnt_ps
Quantity array index for surface pressure.
Definition: mptrac.h:2246
int rng_type
Random number generator (0=GSL, 1=Squares, 2=cuRAND).
Definition: mptrac.h:2742
int isosurf
Isosurface parameter (0=none, 1=pressure, 2=density, 3=theta, 4=balloon).
Definition: mptrac.h:2729
double bound_p1
Boundary conditions top pressure [hPa].
Definition: mptrac.h:2808
int qnt_zs
Quantity array index for surface geopotential height.
Definition: mptrac.h:2252
int prof_nz
Number of altitudes of gridded profile data.
Definition: mptrac.h:3151
double csi_dt_out
Time step for CSI output [s].
Definition: mptrac.h:3043
int met_cape
Convective available potential energy data (0=file, 1=calculate).
Definition: mptrac.h:2687
double csi_modmin
Minimum column density to trigger detection [kg/m^2].
Definition: mptrac.h:3052
int met_sx
Smoothing for longitudes.
Definition: mptrac.h:2648
double chemgrid_lon1
Upper longitude of chemistry grid [deg].
Definition: mptrac.h:2922
double turb_mesox
Horizontal scaling factor for mesoscale wind fluctuations.
Definition: mptrac.h:2766
double met_cms_eps_iwc
cmultiscale compression epsilon for cloud ice water content.
Definition: mptrac.h:2630
double met_cms_eps_swc
cmultiscale compression epsilon for cloud snow water content.
Definition: mptrac.h:2633
double met_cms_eps_v
cmultiscale compression epsilon for meridional wind.
Definition: mptrac.h:2609
double prof_z0
Lower altitude of gridded profile data [km].
Definition: mptrac.h:3154
int qnt_w
Quantity array index for vertical velocity.
Definition: mptrac.h:2309
double bound_vmr
Boundary conditions volume mixing ratio [ppv].
Definition: mptrac.h:2793
double met_tropo_pv
Dynamical tropopause potential vorticity threshold [PVU].
Definition: mptrac.h:2703
int prof_nx
Number of longitudes of gridded profile data.
Definition: mptrac.h:3160
int qnt_stat
Quantity array index for station flag.
Definition: mptrac.h:2231
int met_tropo
Tropopause definition (0=none, 1=clim, 2=cold point, 3=WMO_1st, 4=WMO_2nd, 5=dynamical).
Definition: mptrac.h:2700
int qnt_rp
Quantity array index for particle radius.
Definition: mptrac.h:2240
int met_mpi_share
Use MPI to share meteo (0=no, 1=yes).
Definition: mptrac.h:2718
double mixing_strat
Interparcel exchange parameter for mixing in the stratosphere.
Definition: mptrac.h:2877
int qnt_vz
Quantity array index for vertical velocity.
Definition: mptrac.h:2441
int qnt_ho2
Quantity array index for HO2 volume mixing ratio (climatology).
Definition: mptrac.h:2369
double csi_z1
Upper altitude of gridded CSI data [km].
Definition: mptrac.h:3061
double stat_t0
Start time for station output [s].
Definition: mptrac.h:3205
double oh_chem_beta
Beta parameter for diurnal variablity of OH.
Definition: mptrac.h:2940
int dd
Domain decomposition (0=no, 1=yes, with 2x2 if not specified).
Definition: mptrac.h:3233
int qnt_eta
Quantity array index for eta vertical coordinate.
Definition: mptrac.h:2426
double wet_depo_so2_ph
pH value used to calculate effective Henry constant of SO2.
Definition: mptrac.h:2979
double mixing_z0
Lower altitude of mixing grid [km].
Definition: mptrac.h:2883
int qnt_mloss_decay
Quantity array index for total mass loss due to exponential decay.
Definition: mptrac.h:2390
int atm_type_out
Type of atmospheric data files for output (-1=same as ATM_TYPE, 0=ASCII, 1=binary,...
Definition: mptrac.h:3025
int met_cms_nd0x
cmultiscale number of cells of coarsest grid in x-direction.
Definition: mptrac.h:2591
int met_nlev
Number of meteo data model levels.
Definition: mptrac.h:2669
double dt_kpp
Time step for KPP chemistry [s].
Definition: mptrac.h:2949
double dry_depo_dp
Dry deposition surface layer [hPa].
Definition: mptrac.h:2988
int qnt_shf
Quantity array index for surface sensible heat flux.
Definition: mptrac.h:2267
int qnt_vs
Quantity array index for surface meridional wind.
Definition: mptrac.h:2258
int qnt_Cco
Quantity array index for CO volume mixing ratio (chemistry code).
Definition: mptrac.h:2468
double vtk_dt_out
Time step for VTK data output [s].
Definition: mptrac.h:3214
double t_stop
Stop time of simulation [s].
Definition: mptrac.h:2537
double conv_dt
Time interval for convection module [s].
Definition: mptrac.h:2784
int qnt_hno3
Quantity array index for HNO3 volume mixing ratio (climatology).
Definition: mptrac.h:2360
int met_clams
Read MPTRAC or CLaMS meteo data (0=MPTRAC, 1=CLaMS).
Definition: mptrac.h:2564
int qnt_h2ot
Quantity array index for tropopause water vapor volume mixing ratio.
Definition: mptrac.h:2288
int qnt_rh
Quantity array index for relative humidity over water.
Definition: mptrac.h:2408
double met_cms_eps_cc
cmultiscale compression epsilon for cloud cover.
Definition: mptrac.h:2636
double bound_lat0
Boundary conditions minimum longitude [deg].
Definition: mptrac.h:2799
double met_pbl_max
Maximum depth of planetary boundary layer [km].
Definition: mptrac.h:2696
int met_dx
Stride for longitudes.
Definition: mptrac.h:2639
int qnt_destination
Quantity array index for destination subdomain in domain decomposition.
Definition: mptrac.h:2528
int mixing_ny
Number of latitudes of mixing grid.
Definition: mptrac.h:2898
int met_convention
Meteo data layout (0=[lev, lat, lon], 1=[lon, lat, lev]).
Definition: mptrac.h:2553
int qnt_zeta_d
Quantity array index for diagnosed zeta vertical coordinate.
Definition: mptrac.h:2420
int tracer_chem
Switch for first order tracer chemistry module (0=off, 1=on).
Definition: mptrac.h:2952
double dt_mod
Time step of simulation [s].
Definition: mptrac.h:2540
int diffusion
Diffusion scheme (0=off, 1=fixed-K, 2=PBL).
Definition: mptrac.h:2745
int qnt_tnat
Quantity array index for T_NAT.
Definition: mptrac.h:2456
int qnt_eta_dot
Quantity array index for velocity of eta vertical coordinate.
Definition: mptrac.h:2429
int qnt_tice
Quantity array index for T_ice.
Definition: mptrac.h:2450
int qnt_zg
Quantity array index for geopotential height.
Definition: mptrac.h:2291
double vtk_offset
Vertical offset for VTK data [km].
Definition: mptrac.h:3223
int qnt_v
Quantity array index for meridional wind.
Definition: mptrac.h:2306
int qnt_mloss_dry
Quantity array index for total mass loss due to dry deposition.
Definition: mptrac.h:2387
double bound_vmr_trend
Boundary conditions volume mixing ratio trend [ppv/s].
Definition: mptrac.h:2796
int met_cache
Preload meteo data into disk cache (0=no, 1=yes).
Definition: mptrac.h:2715
int qnt_oh
Quantity array index for OH volume mixing ratio (climatology).
Definition: mptrac.h:2363
int qnt_Ch
Quantity array index for H volume mixing ratio (chemistry code).
Definition: mptrac.h:2474
int met_press_level_def
Use predefined pressure levels or not.
Definition: mptrac.h:2666
int oh_chem_reaction
Reaction type for OH chemistry (0=none, 2=bimolecular, 3=termolecular).
Definition: mptrac.h:2934
int qnt_h2o
Quantity array index for water vapor volume mixing ratio.
Definition: mptrac.h:2312
int prof_ny
Number of latitudes of gridded profile data.
Definition: mptrac.h:3169
int qnt_rhice
Quantity array index for relative humidity over ice.
Definition: mptrac.h:2411
int qnt_rho
Quantity array index for density of air.
Definition: mptrac.h:2300
double sample_dz
Layer depth for sample output [km].
Definition: mptrac.h:3190
double tdec_strat
Life time of particles in the stratosphere [s].
Definition: mptrac.h:2832
int obs_type
Type of observation data files (0=ASCII, 1=netCDF).
Definition: mptrac.h:3034
double met_cms_eps_lwc
cmultiscale compression epsilon for cloud liquid water content.
Definition: mptrac.h:2624
int qnt_us
Quantity array index for surface zonal wind.
Definition: mptrac.h:2255
double met_cms_eps_z
cmultiscale compression epsilon for geopotential height.
Definition: mptrac.h:2600
double grid_lon1
Upper longitude of gridded data [deg].
Definition: mptrac.h:3130
int qnt_Cn2o
Quantity array index for N2O volume mixing ratio (chemistry code).
Definition: mptrac.h:2498
int qnt_Cccl3f
Quantity array index for CFC-11 volume mixing ratio (chemistry code).
Definition: mptrac.h:2492
double mixing_lat0
Lower latitude of mixing grid [deg].
Definition: mptrac.h:2901
int nens
Number of ensembles.
Definition: mptrac.h:3082
int qnt_pt
Quantity array index for tropopause pressure.
Definition: mptrac.h:2279
int qnt_cl
Quantity array index for total column cloud water.
Definition: mptrac.h:2339
int advect
Advection scheme (0=off, 1=Euler, 2=midpoint, 4=Runge-Kutta).
Definition: mptrac.h:2735
double prof_z1
Upper altitude of gridded profile data [km].
Definition: mptrac.h:3157
int qnt_t
Quantity array index for temperature.
Definition: mptrac.h:2297
int atm_filter
Time filter for atmospheric data output (0=none, 1=missval, 2=remove).
Definition: mptrac.h:3013
int kpp_chem
Switch for KPP chemistry module (0=off, 1=on).
Definition: mptrac.h:2946
int qnt_zeta
Quantity array index for zeta vertical coordinate.
Definition: mptrac.h:2417
double conv_pbl_trans
Depth of PBL transition layer (fraction of PBL depth).
Definition: mptrac.h:2775
int qnt_Ai131
Quantity array index for radioactive activity of I-131.
Definition: mptrac.h:2519
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:2557
double csi_z0
Lower altitude of gridded CSI data [km].
Definition: mptrac.h:3058
int qnt_lapse
Quantity array index for lapse rate.
Definition: mptrac.h:2435
int qnt_Apb210
Quantity array index for radioactive activity of Pb-210.
Definition: mptrac.h:2510
double stat_lat
Latitude of station [deg].
Definition: mptrac.h:3199
int qnt_Cho2
Quantity array index for HO2 volume mixing ratio (chemistry code).
Definition: mptrac.h:2477
int grid_ny
Number of latitudes of gridded data.
Definition: mptrac.h:3133
int qnt_Csf6
Quantity array index for SF6 volume mixing ratio (chemistry code).
Definition: mptrac.h:2501
int qnt_Ch2o
Quantity array index for H2O volume mixing ratio (chemistry code).
Definition: mptrac.h:2462
double met_detrend
FWHM of horizontal Gaussian used for detrending [km].
Definition: mptrac.h:2657
int conv_mix_pbl
Vertical mixing in the PBL (0=off, 1=on).
Definition: mptrac.h:2772
double bound_dps
Boundary conditions surface layer depth [hPa].
Definition: mptrac.h:2811
int dd_nbr_neighbours
Domain decomposition number of neighbours to communicate with.
Definition: mptrac.h:3242
double met_cms_eps_t
cmultiscale compression epsilon for temperature.
Definition: mptrac.h:2603
int chemgrid_nz
Number of altitudes of chemistry grid.
Definition: mptrac.h:2907
int qnt_cape
Quantity array index for convective available potential energy (CAPE).
Definition: mptrac.h:2351
int qnt_zeta_dot
Quantity array index for velocity of zeta vertical coordinate.
Definition: mptrac.h:2423
double bound_mass_trend
Boundary conditions mass per particle trend [kg/s].
Definition: mptrac.h:2790
int met_cms_nd0y
cmultiscale number of cells of coarsest grid in y-direction.
Definition: mptrac.h:2594
int mixing_nz
Number of altitudes of mixing grid.
Definition: mptrac.h:2880
int qnt_o3c
Quantity array index for total column ozone.
Definition: mptrac.h:2357
double bound_p0
Boundary conditions bottom pressure [hPa].
Definition: mptrac.h:2805
double mixing_lon0
Lower longitude of mixing grid [deg].
Definition: mptrac.h:2892
int qnt_Co3
Quantity array index for O3 volume mixing ratio (chemistry code).
Definition: mptrac.h:2465
int qnt_tsts
Quantity array index for T_STS.
Definition: mptrac.h:2453
int grid_nz
Number of altitudes of gridded data.
Definition: mptrac.h:3115
int qnt_nss
Quantity array index for northward turbulent surface stress.
Definition: mptrac.h:2264
double ens_dt_out
Time step for ensemble output [s].
Definition: mptrac.h:3088
int atm_stride
Particle index stride for atmospheric data files.
Definition: mptrac.h:3016
int met_relhum
Try to read relative humidity (0=no, 1=yes).
Definition: mptrac.h:2684
double mixing_lat1
Upper latitude of mixing grid [deg].
Definition: mptrac.h:2904
double atm_dt_out
Time step for atmospheric data output [s].
Definition: mptrac.h:3010
double prof_lat1
Upper latitude of gridded profile data [deg].
Definition: mptrac.h:3175
int met_cms_batch
cmultiscale batch size.
Definition: mptrac.h:2585
double psc_h2o
H2O volume mixing ratio for PSC analysis.
Definition: mptrac.h:2994
int met_sp
Smoothing for pressure levels.
Definition: mptrac.h:2654
double prof_lon0
Lower longitude of gridded profile data [deg].
Definition: mptrac.h:3163
int qnt_Axe133
Quantity array index for radioactive activity of Xe-133.
Definition: mptrac.h:2522
int chemgrid_nx
Number of longitudes of chemistry grid.
Definition: mptrac.h:2916
int qnt_pct
Quantity array index for cloud top pressure.
Definition: mptrac.h:2333
int qnt_mloss_kpp
Quantity array index for total mass loss due to KPP chemistry.
Definition: mptrac.h:2381
int qnt_psat
Quantity array index for saturation pressure over water.
Definition: mptrac.h:2396
int qnt_subdomain
Quantity array index for current subdomain in domain decomposition.
Definition: mptrac.h:2525
double prof_lat0
Lower latitude of gridded profile data [deg].
Definition: mptrac.h:3172
int qnt_cin
Quantity array index for convective inhibition (CIN).
Definition: mptrac.h:2354
double psc_hno3
HNO3 volume mixing ratio for PSC analysis.
Definition: mptrac.h:2997
double prof_lon1
Upper longitude of gridded profile data [deg].
Definition: mptrac.h:3166
double met_cms_eps_rwc
cmultiscale compression epsilon for cloud rain water content.
Definition: mptrac.h:2627
int met_nc_quant
Number of digits for quantization of netCDF meteo files (0=off).
Definition: mptrac.h:2573
int h2o2_chem_reaction
Reaction type for H2O2 chemistry (0=none, 1=SO2).
Definition: mptrac.h:2943
int qnt_Co3p
Quantity array index for O(3P) volume mixing ratio (chemistry code).
Definition: mptrac.h:2486
double wet_depo_bc_ret_ratio
Coefficients for wet deposition below cloud: retention ratio.
Definition: mptrac.h:2985
int chemgrid_ny
Number of latitudes of chemistry grid.
Definition: mptrac.h:2925
int qnt_Abe7
Quantity array index for radioactive activity of Be-7.
Definition: mptrac.h:2513
double met_cms_eps_o3
cmultiscale compression epsilon for ozone.
Definition: mptrac.h:2621
int met_cms_zstd
cmultiscale ZSTD compression (0=off, 1=on).
Definition: mptrac.h:2588
int met_cms_maxlev
cmultiscale maximum refinement level.
Definition: mptrac.h:2597
int grid_sparse
Sparse output in grid data files (0=no, 1=yes).
Definition: mptrac.h:3103
double dry_depo_vdep
Dry deposition velocity [m/s].
Definition: mptrac.h:2991
int qnt_tt
Quantity array index for tropopause temperature.
Definition: mptrac.h:2282
int met_np
Number of target pressure levels.
Definition: mptrac.h:2660
int qnt_ens
Quantity array index for ensemble IDs.
Definition: mptrac.h:2228
int met_nc_level
zlib compression level of netCDF meteo files (0=off).
Definition: mptrac.h:2570
double mixing_dt
Time interval for mixing [s].
Definition: mptrac.h:2871
int qnt_Arn222
Quantity array index for radioactive activity of Rn-222.
Definition: mptrac.h:2507
int qnt_mloss_h2o2
Quantity array index for total mass loss due to H2O2 chemistry.
Definition: mptrac.h:2378
double vtk_scale
Vertical scaling factor for VTK data.
Definition: mptrac.h:3220
double met_cms_eps_w
cmultiscale compression epsilon for vertical velocity.
Definition: mptrac.h:2612
double turb_dx_pbl
Horizontal turbulent diffusion coefficient (PBL) [m^2/s].
Definition: mptrac.h:2748
double conv_cin
CIN threshold for convection module [J/kg].
Definition: mptrac.h:2781
int qnt_pv
Quantity array index for potential vorticity.
Definition: mptrac.h:2444
int advect_vert_coord
Vertical velocity of air parcels (0=omega_on_plev, 1=zetadot_on_mlev, 2=omega_on_mlev,...
Definition: mptrac.h:2739
int qnt_mloss_oh
Quantity array index for total mass loss due to OH chemistry.
Definition: mptrac.h:2375
int qnt_Ch2o2
Quantity array index for H2O2 volume mixing ratio (chemistry code).
Definition: mptrac.h:2480
int qnt_sst
Quantity array index for sea surface temperature.
Definition: mptrac.h:2273
double mixing_lon1
Upper longitude of mixing grid [deg].
Definition: mptrac.h:2895
int atm_nc_level
zlib compression level of netCDF atmospheric data files (0=off).
Definition: mptrac.h:3028
double wet_depo_ic_ret_ratio
Coefficients for wet deposition in cloud: retention ratio.
Definition: mptrac.h:2982
int qnt_sh
Quantity array index for specific humidity.
Definition: mptrac.h:2405
int qnt_ess
Quantity array index for eastward turbulent surface stress.
Definition: mptrac.h:2261
double wet_depo_ic_b
Coefficient B for wet deposition in cloud (exponential form).
Definition: mptrac.h:2970
double wet_depo_bc_b
Coefficient B for wet deposition below cloud (exponential form).
Definition: mptrac.h:2964
int met_dy
Stride for latitudes.
Definition: mptrac.h:2642
int qnt_Cx
Quantity array index for trace species x volume mixing ratio (chemistry code).
Definition: mptrac.h:2459
double turb_dz_strat
Vertical turbulent diffusion coefficient (stratosphere) [m^2/s].
Definition: mptrac.h:2763
double bound_zetas
Boundary conditions surface layer zeta [K].
Definition: mptrac.h:2817
int dd_subdomains_zonal
Domain decomposition zonal subdomain number.
Definition: mptrac.h:3236
int qnt_idx
Quantity array index for air parcel IDs.
Definition: mptrac.h:2225
double met_tropo_theta
Dynamical tropopause potential temperature threshold [K].
Definition: mptrac.h:2706
int qnt_rwc
Quantity array index for cloud rain water content.
Definition: mptrac.h:2321
double t_start
Start time of simulation [s].
Definition: mptrac.h:2534
int nq
Number of quantities.
Definition: mptrac.h:2210
double tdec_trop
Life time of particles in the troposphere [s].
Definition: mptrac.h:2829
double sample_dx
Horizontal radius for sample output [km].
Definition: mptrac.h:3187
int vtk_stride
Particle index stride for VTK data.
Definition: mptrac.h:3217
double turb_dz_pbl
Vertical turbulent diffusion coefficient (PBL) [m^2/s].
Definition: mptrac.h:2757
double grid_lat1
Upper latitude of gridded data [deg].
Definition: mptrac.h:3139
int dd_subdomains_meridional
Domain decomposition meridional subdomain number.
Definition: mptrac.h:3239
int qnt_zt
Quantity array index for tropopause geopotential height.
Definition: mptrac.h:2285
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:2561
int qnt_cc
Quantity array index for cloud cover.
Definition: mptrac.h:2330
int qnt_plcl
Quantity array index for pressure at lifted condensation level (LCL).
Definition: mptrac.h:2342
double grid_dt_out
Time step for gridded data output [s].
Definition: mptrac.h:3100
int qnt_tdew
Quantity array index for dew point temperature.
Definition: mptrac.h:2447
Domain decomposition data structure.
Definition: mptrac.h:3685
int halo_offset_end
Hyperslab of boundary halos count.
Definition: mptrac.h:3737
int rank
Rank of node.
Definition: mptrac.h:3692
double subdomain_lon_min
Rectangular grid limit of subdomain.
Definition: mptrac.h:3713
double subdomain_lat_max
Rectangular grid limit of subdomain.
Definition: mptrac.h:3716
int init
Shows if domain decomposition was initialized.
Definition: mptrac.h:3744
double subdomain_lon_max
Rectangular grid limit of subdomain.
Definition: mptrac.h:3710
int halo_offset_start
Hyperslab of boundary halos count.
Definition: mptrac.h:3734
int size
Size of node.
Definition: mptrac.h:3695
double subdomain_lat_min
Rectangular grid limit of subdomain.
Definition: mptrac.h:3719
Meteo data structure.
Definition: mptrac.h:3511
int nx
Number of longitudes.
Definition: mptrac.h:3517
int ny
Number of latitudes.
Definition: mptrac.h:3520
int np
Number of pressure levels.
Definition: mptrac.h:3523
int npl
Number of model levels.
Definition: mptrac.h:3526
double time
Time [s].
Definition: mptrac.h:3514
Particle data.
Definition: mptrac.h:3286
double p
Pressure [hPa].
Definition: mptrac.h:3292
double lat
Latitude [deg].
Definition: mptrac.h:3298
double time
Time [s].
Definition: mptrac.h:3289
double lon
Longitude [deg].
Definition: mptrac.h:3295