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-2025 Forschungszentrum Juelich GmbH
18*/
19
97#ifndef LIBTRAC_H
98#define LIBTRAC_H
99
100/* ------------------------------------------------------------
101 Includes...
102 ------------------------------------------------------------ */
103
104#include <ctype.h>
105#include <gsl/gsl_fft_complex.h>
106#include <gsl/gsl_math.h>
107#include <gsl/gsl_randist.h>
108#include <gsl/gsl_rng.h>
109#include <gsl/gsl_spline.h>
110#include <gsl/gsl_statistics.h>
111#include <math.h>
112#include <netcdf.h>
113#include <omp.h>
114#include <stdint.h>
115#include <stdio.h>
116#include <stdlib.h>
117#include <string.h>
118#include <time.h>
119#include <sys/time.h>
120
121#ifdef MPI
122#include "mpi.h"
123#endif
124
125#ifdef _OPENACC
126#include "openacc.h"
127#endif
128
129#ifdef CURAND
130#include "curand.h"
131#endif
132
133#ifdef ZFP
134#include "zfp.h"
135#endif
136
137#ifdef ZSTD
138#include "zstd.h"
139#endif
140
141#ifdef CMS
142#include "cmultiscale.h"
143#endif
144
145#ifdef KPP
146#include "chem_Parameters.h"
147#include "chem_Global.h"
148#include "chem_Sparse.h"
149#endif
150
151/* ------------------------------------------------------------
152 Constants...
153 ------------------------------------------------------------ */
154
156#ifndef AVO
157#define AVO 6.02214076e23
158#endif
159
161#ifndef CPD
162#define CPD 1003.5
163#endif
164
166#ifndef EPS
167#define EPS (MH2O / MA)
168#endif
169
171#ifndef G0
172#define G0 9.80665
173#endif
174
176#ifndef H0
177#define H0 7.0
178#endif
179
181#ifndef LV
182#define LV 2501000.
183#endif
184
186#ifndef KARMAN
187#define KARMAN 0.40
188#endif
189
191#ifndef KB
192#define KB 1.3806504e-23
193#endif
194
196#ifndef MA
197#define MA 28.9644
198#endif
199
201#ifndef MH2O
202#define MH2O 18.01528
203#endif
204
206#ifndef MO3
207#define MO3 48.00
208#endif
209
211#ifndef P0
212#define P0 1013.25
213#endif
214
216#ifndef RA
217#define RA (1e3 * RI / MA)
218#endif
219
221#ifndef RE
222#define RE 6367.421
223#endif
224
226#ifndef RI
227#define RI 8.3144598
228#endif
229
231#ifndef T0
232#define T0 273.15
233#endif
234
235/* ------------------------------------------------------------
236 Dimensions...
237 ------------------------------------------------------------ */
238
240#ifndef LEN
241#define LEN 5000
242#endif
243
245#ifndef NP
246#define NP 10000000
247#endif
248
250#ifndef NQ
251#define NQ 15
252#endif
253
255#ifndef NCSI
256#define NCSI 1000000
257#endif
258
260#ifndef EP
261#define EP 140
262#endif
263
265#ifndef EX
266#define EX 1202
267#endif
268
270#ifndef EY
271#define EY 602
272#endif
273
275#ifndef NENS
276#define NENS 2000
277#endif
278
280#ifndef NOBS
281#define NOBS 10000000
282#endif
283
285#ifndef NTHREADS
286#define NTHREADS 512
287#endif
288
290#ifndef CY
291#define CY 250
292#endif
293
295#ifndef CO3
296#define CO3 30
297#endif
298
300#ifndef CP
301#define CP 70
302#endif
303
305#ifndef CSZA
306#define CSZA 50
307#endif
308
310#ifndef CT
311#define CT 12
312#endif
313
315#ifndef CTS
316#define CTS 1000
317#endif
318
319/* ------------------------------------------------------------
320 Macros...
321 ------------------------------------------------------------ */
322
342#ifdef _OPENACC
343#define ALLOC(ptr, type, n) \
344 if(acc_get_num_devices(acc_device_nvidia) <= 0) \
345 ERRMSG("Not running on a GPU device!"); \
346 if((ptr=calloc((size_t)(n), sizeof(type)))==NULL) \
347 ERRMSG("Out of memory!");
348#else
349#define ALLOC(ptr, type, n) \
350 if((ptr=calloc((size_t)(n), sizeof(type)))==NULL) \
351 ERRMSG("Out of memory!");
352#endif
353
372#define ARRAY_2D(ix, iy, ny) \
373 ((ix) * (ny) + (iy))
374
391#define ARRAY_3D(ix, iy, ny, iz, nz) \
392 (((ix)*(ny) + (iy)) * (nz) + (iz))
393
416#define ARRHENIUS(a, b, t) \
417 ((a) * exp( -(b) / (t)))
418
440#define DEG2DX(dlon, lat) \
441 (RE * DEG2RAD(dlon) * cos(DEG2RAD(lat)))
442
461#define DEG2DY(dlat) \
462 (RE * DEG2RAD(dlat))
463
478#define DEG2RAD(deg) \
479 ((deg) * (M_PI / 180.0))
480
503#define DP2DZ(dp, p) \
504 (- (dp) * H0 / (p))
505
525#define DX2DEG(dx, lat) \
526 (((lat) < -89.999 || (lat) > 89.999) ? 0 \
527 : (dx) * 180. / (M_PI * RE * cos(DEG2RAD(lat))))
528
543#define DY2DEG(dy) \
544 ((dy) * 180. / (M_PI * RE))
545
562#define DZ2DP(dz, p) \
563 (-(dz) * (p) / H0)
564
578#define DIST(a, b) \
579 sqrt(DIST2(a, b))
580
594#define DIST2(a, b) \
595 ((a[0]-b[0])*(a[0]-b[0])+(a[1]-b[1])*(a[1]-b[1])+(a[2]-b[2])*(a[2]-b[2]))
596
610#define DOTP(a, b) \
611 (a[0]*b[0]+a[1]*b[1]+a[2]*b[2])
612
629#define FMOD(x, y) \
630 ((x) - (int) ((x) / (y)) * (y))
631
647#define FREAD(ptr, type, size, in) { \
648 if(fread(ptr, sizeof(type), size, in)!=size) \
649 ERRMSG("Error while reading!"); \
650 }
651
667#define FWRITE(ptr, type, size, out) { \
668 if(fwrite(ptr, sizeof(type), size, out)!=size) \
669 ERRMSG("Error while writing!"); \
670 }
671
682#define INTPOL_INIT \
683 double cw[4] = {0.0, 0.0, 0.0, 0.0}; int ci[3] = {0, 0, 0};
684
696#define INTPOL_2D(var, init) \
697 intpol_met_time_2d(met0, met0->var, met1, met1->var, \
698 atm->time[ip], atm->lon[ip], atm->lat[ip], \
699 &var, ci, cw, init);
700
713#define INTPOL_3D(var, init) \
714 intpol_met_time_3d(met0, met0->var, met1, met1->var, \
715 atm->time[ip], atm->p[ip], \
716 atm->lon[ip], atm->lat[ip], \
717 &var, ci, cw, init);
718
732#define INTPOL_SPACE_ALL(p, lon, lat) { \
733 intpol_met_space_3d(met, met->z, p, lon, lat, &z, ci, cw, 1); \
734 intpol_met_space_3d(met, met->t, p, lon, lat, &t, ci, cw, 0); \
735 intpol_met_space_3d(met, met->u, p, lon, lat, &u, ci, cw, 0); \
736 intpol_met_space_3d(met, met->v, p, lon, lat, &v, ci, cw, 0); \
737 intpol_met_space_3d(met, met->w, p, lon, lat, &w, ci, cw, 0); \
738 intpol_met_space_3d(met, met->pv, p, lon, lat, &pv, ci, cw, 0); \
739 intpol_met_space_3d(met, met->h2o, p, lon, lat, &h2o, ci, cw, 0); \
740 intpol_met_space_3d(met, met->o3, p, lon, lat, &o3, ci, cw, 0); \
741 intpol_met_space_3d(met, met->lwc, p, lon, lat, &lwc, ci, cw, 0); \
742 intpol_met_space_3d(met, met->rwc, p, lon, lat, &rwc, ci, cw, 0); \
743 intpol_met_space_3d(met, met->iwc, p, lon, lat, &iwc, ci, cw, 0); \
744 intpol_met_space_3d(met, met->swc, p, lon, lat, &swc, ci, cw, 0); \
745 intpol_met_space_3d(met, met->cc, p, lon, lat, &cc, ci, cw, 0); \
746 intpol_met_space_2d(met, met->ps, lon, lat, &ps, ci, cw, 0); \
747 intpol_met_space_2d(met, met->ts, lon, lat, &ts, ci, cw, 0); \
748 intpol_met_space_2d(met, met->zs, lon, lat, &zs, ci, cw, 0); \
749 intpol_met_space_2d(met, met->us, lon, lat, &us, ci, cw, 0); \
750 intpol_met_space_2d(met, met->vs, lon, lat, &vs, ci, cw, 0); \
751 intpol_met_space_2d(met, met->ess, ess, lat, &ess, ci, cw, 0); \
752 intpol_met_space_2d(met, met->nss, nss, lat, &nss, ci, cw, 0); \
753 intpol_met_space_2d(met, met->shf, shf, lat, &shf, ci, cw, 0); \
754 intpol_met_space_2d(met, met->lsm, lon, lat, &lsm, ci, cw, 0); \
755 intpol_met_space_2d(met, met->sst, lon, lat, &sst, ci, cw, 0); \
756 intpol_met_space_2d(met, met->pbl, lon, lat, &pbl, ci, cw, 0); \
757 intpol_met_space_2d(met, met->pt, lon, lat, &pt, ci, cw, 0); \
758 intpol_met_space_2d(met, met->tt, lon, lat, &tt, ci, cw, 0); \
759 intpol_met_space_2d(met, met->zt, lon, lat, &zt, ci, cw, 0); \
760 intpol_met_space_2d(met, met->h2ot, lon, lat, &h2ot, ci, cw, 0); \
761 intpol_met_space_2d(met, met->pct, lon, lat, &pct, ci, cw, 0); \
762 intpol_met_space_2d(met, met->pcb, lon, lat, &pcb, ci, cw, 0); \
763 intpol_met_space_2d(met, met->cl, lon, lat, &cl, ci, cw, 0); \
764 intpol_met_space_2d(met, met->plcl, lon, lat, &plcl, ci, cw, 0); \
765 intpol_met_space_2d(met, met->plfc, lon, lat, &plfc, ci, cw, 0); \
766 intpol_met_space_2d(met, met->pel, lon, lat, &pel, ci, cw, 0); \
767 intpol_met_space_2d(met, met->cape, lon, lat, &cape, ci, cw, 0); \
768 intpol_met_space_2d(met, met->cin, lon, lat, &cin, ci, cw, 0); \
769 intpol_met_space_2d(met, met->o3c, lon, lat, &o3c, ci, cw, 0); \
770 }
771
786#define INTPOL_TIME_ALL(time, p, lon, lat) { \
787 intpol_met_time_3d(met0, met0->z, met1, met1->z, time, p, lon, lat, &z, ci, cw, 1); \
788 intpol_met_time_3d(met0, met0->t, met1, met1->t, time, p, lon, lat, &t, ci, cw, 0); \
789 intpol_met_time_3d(met0, met0->u, met1, met1->u, time, p, lon, lat, &u, ci, cw, 0); \
790 intpol_met_time_3d(met0, met0->v, met1, met1->v, time, p, lon, lat, &v, ci, cw, 0); \
791 intpol_met_time_3d(met0, met0->w, met1, met1->w, time, p, lon, lat, &w, ci, cw, 0); \
792 intpol_met_time_3d(met0, met0->pv, met1, met1->pv, time, p, lon, lat, &pv, ci, cw, 0); \
793 intpol_met_time_3d(met0, met0->h2o, met1, met1->h2o, time, p, lon, lat, &h2o, ci, cw, 0); \
794 intpol_met_time_3d(met0, met0->o3, met1, met1->o3, time, p, lon, lat, &o3, ci, cw, 0); \
795 intpol_met_time_3d(met0, met0->lwc, met1, met1->lwc, time, p, lon, lat, &lwc, ci, cw, 0); \
796 intpol_met_time_3d(met0, met0->rwc, met1, met1->rwc, time, p, lon, lat, &rwc, ci, cw, 0); \
797 intpol_met_time_3d(met0, met0->iwc, met1, met1->iwc, time, p, lon, lat, &iwc, ci, cw, 0); \
798 intpol_met_time_3d(met0, met0->swc, met1, met1->swc, time, p, lon, lat, &swc, ci, cw, 0); \
799 intpol_met_time_3d(met0, met0->cc, met1, met1->cc, time, p, lon, lat, &cc, ci, cw, 0); \
800 intpol_met_time_2d(met0, met0->ps, met1, met1->ps, time, lon, lat, &ps, ci, cw, 0); \
801 intpol_met_time_2d(met0, met0->ts, met1, met1->ts, time, lon, lat, &ts, ci, cw, 0); \
802 intpol_met_time_2d(met0, met0->zs, met1, met1->zs, time, lon, lat, &zs, ci, cw, 0); \
803 intpol_met_time_2d(met0, met0->us, met1, met1->us, time, lon, lat, &us, ci, cw, 0); \
804 intpol_met_time_2d(met0, met0->vs, met1, met1->vs, time, lon, lat, &vs, ci, cw, 0); \
805 intpol_met_time_2d(met0, met0->ess, met1, met1->ess, time, lon, lat, &ess, ci, cw, 0); \
806 intpol_met_time_2d(met0, met0->nss, met1, met1->nss, time, lon, lat, &nss, ci, cw, 0); \
807 intpol_met_time_2d(met0, met0->shf, met1, met1->shf, time, lon, lat, &shf, ci, cw, 0); \
808 intpol_met_time_2d(met0, met0->lsm, met1, met1->lsm, time, lon, lat, &lsm, ci, cw, 0); \
809 intpol_met_time_2d(met0, met0->sst, met1, met1->sst, time, lon, lat, &sst, ci, cw, 0); \
810 intpol_met_time_2d(met0, met0->pbl, met1, met1->pbl, time, lon, lat, &pbl, ci, cw, 0); \
811 intpol_met_time_2d(met0, met0->pt, met1, met1->pt, time, lon, lat, &pt, ci, cw, 0); \
812 intpol_met_time_2d(met0, met0->tt, met1, met1->tt, time, lon, lat, &tt, ci, cw, 0); \
813 intpol_met_time_2d(met0, met0->zt, met1, met1->zt, time, lon, lat, &zt, ci, cw, 0); \
814 intpol_met_time_2d(met0, met0->h2ot, met1, met1->h2ot, time, lon, lat, &h2ot, ci, cw, 0); \
815 intpol_met_time_2d(met0, met0->pct, met1, met1->pct, time, lon, lat, &pct, ci, cw, 0); \
816 intpol_met_time_2d(met0, met0->pcb, met1, met1->pcb, time, lon, lat, &pcb, ci, cw, 0); \
817 intpol_met_time_2d(met0, met0->cl, met1, met1->cl, time, lon, lat, &cl, ci, cw, 0); \
818 intpol_met_time_2d(met0, met0->plcl, met1, met1->plcl, time, lon, lat, &plcl, ci, cw, 0); \
819 intpol_met_time_2d(met0, met0->plfc, met1, met1->plfc, time, lon, lat, &plfc, ci, cw, 0); \
820 intpol_met_time_2d(met0, met0->pel, met1, met1->pel, time, lon, lat, &pel, ci, cw, 0); \
821 intpol_met_time_2d(met0, met0->cape, met1, met1->cape, time, lon, lat, &cape, ci, cw, 0); \
822 intpol_met_time_2d(met0, met0->cin, met1, met1->cin, time, lon, lat, &cin, ci, cw, 0); \
823 intpol_met_time_2d(met0, met0->o3c, met1, met1->o3c, time, lon, lat, &o3c, ci, cw, 0); \
824 }
825
840#define LAPSE(p1, t1, p2, t2) \
841 (1e3 * G0 / RA * ((t2) - (t1)) / ((t2) + (t1)) \
842 * ((p2) + (p1)) / ((p2) - (p1)))
843
859#define LIN(x0, y0, x1, y1, x) \
860 ((y0)+((y1)-(y0))/((x1)-(x0))*((x)-(x0)))
861
886#define MAX(a,b) \
887 (((a)>(b))?(a):(b))
888
900#define MET_HEADER \
901 fprintf(out, \
902 "# $1 = time [s]\n" \
903 "# $2 = altitude [km]\n" \
904 "# $3 = longitude [deg]\n" \
905 "# $4 = latitude [deg]\n" \
906 "# $5 = pressure [hPa]\n" \
907 "# $6 = temperature [K]\n" \
908 "# $7 = zonal wind [m/s]\n" \
909 "# $8 = meridional wind [m/s]\n" \
910 "# $9 = vertical velocity [hPa/s]\n" \
911 "# $10 = H2O volume mixing ratio [ppv]\n"); \
912 fprintf(out, \
913 "# $11 = O3 volume mixing ratio [ppv]\n" \
914 "# $12 = geopotential height [km]\n" \
915 "# $13 = potential vorticity [PVU]\n" \
916 "# $14 = surface pressure [hPa]\n" \
917 "# $15 = surface temperature [K]\n" \
918 "# $16 = surface geopotential height [km]\n" \
919 "# $17 = surface zonal wind [m/s]\n" \
920 "# $18 = surface meridional wind [m/s]\n" \
921 "# $19 = eastward turbulent surface stress [N/m^2]\n" \
922 "# $20 = northward turbulent surface stress [N/m^2]\n"); \
923 fprintf(out, \
924 "# $21 = surface sensible heat flux [W/m^2]\n" \
925 "# $22 = land-sea mask [1]\n" \
926 "# $23 = sea surface temperature [K]\n" \
927 "# $24 = tropopause pressure [hPa]\n" \
928 "# $25 = tropopause geopotential height [km]\n" \
929 "# $26 = tropopause temperature [K]\n" \
930 "# $27 = tropopause water vapor [ppv]\n" \
931 "# $28 = cloud liquid water content [kg/kg]\n" \
932 "# $29 = cloud rain water content [kg/kg]\n" \
933 "# $30 = cloud ice water content [kg/kg]\n"); \
934 fprintf(out, \
935 "# $31 = cloud snow water content [kg/kg]\n" \
936 "# $32 = cloud cover [1]\n" \
937 "# $33 = total column cloud water [kg/m^2]\n" \
938 "# $34 = cloud top pressure [hPa]\n" \
939 "# $35 = cloud bottom pressure [hPa]\n" \
940 "# $36 = pressure at lifted condensation level (LCL) [hPa]\n" \
941 "# $37 = pressure at level of free convection (LFC) [hPa]\n" \
942 "# $38 = pressure at equilibrium level (EL) [hPa]\n" \
943 "# $39 = convective available potential energy (CAPE) [J/kg]\n" \
944 "# $40 = convective inhibition (CIN) [J/kg]\n"); \
945 fprintf(out, \
946 "# $41 = relative humidity over water [%%]\n" \
947 "# $42 = relative humidity over ice [%%]\n" \
948 "# $43 = dew point temperature [K]\n" \
949 "# $44 = frost point temperature [K]\n" \
950 "# $45 = NAT temperature [K]\n" \
951 "# $46 = HNO3 volume mixing ratio [ppv]\n" \
952 "# $47 = OH volume mixing ratio [ppv]\n" \
953 "# $48 = H2O2 volume mixing ratio [ppv]\n" \
954 "# $49 = HO2 volume mixing ratio [ppv]\n" \
955 "# $50 = O(1D) volume mixing ratio [ppv]\n"); \
956 fprintf(out, \
957 "# $51 = boundary layer pressure [hPa]\n" \
958 "# $52 = total column ozone [DU]\n" \
959 "# $53 = number of data points\n" \
960 "# $54 = number of tropopause data points\n" \
961 "# $55 = number of CAPE data points\n");
962
987#define MIN(a,b) \
988 (((a)<(b))?(a):(b))
989
1002#define MOLEC_DENS(p,t) \
1003 (AVO * 1e-6 * ((p) * 100) / (RI * (t)))
1004
1016#define NC(cmd) { \
1017 int nc_result=(cmd); \
1018 if(nc_result!=NC_NOERR) \
1019 ERRMSG("%s", nc_strerror(nc_result)); \
1020 }
1021
1045#define NC_DEF_VAR(varname, type, ndims, dims, long_name, units, level, quant) { \
1046 NC(nc_def_var(ncid, varname, type, ndims, dims, &varid)); \
1047 NC(nc_put_att_text(ncid, varid, "long_name", strnlen(long_name, LEN), long_name)); \
1048 NC(nc_put_att_text(ncid, varid, "units", strnlen(units, LEN), units)); \
1049 if((quant) > 0) \
1050 NC(nc_def_var_quantize(ncid, varid, NC_QUANTIZE_GRANULARBR, quant)); \
1051 if((level) != 0) { \
1052 NC(nc_def_var_deflate(ncid, varid, 1, 1, level)); \
1053 /* unsigned int ulevel = (unsigned int)level; */ \
1054 /* NC(nc_def_var_filter(ncid, varid, 32015, 1, (unsigned int[]){ulevel})); */ \
1055 } \
1056 }
1057
1075#define NC_GET_DOUBLE(varname, ptr, force) { \
1076 if(force) { \
1077 NC(nc_inq_varid(ncid, varname, &varid)); \
1078 NC(nc_get_var_double(ncid, varid, ptr)); \
1079 } else { \
1080 if(nc_inq_varid(ncid, varname, &varid) == NC_NOERR) { \
1081 NC(nc_get_var_double(ncid, varid, ptr)); \
1082 } else \
1083 WARN("netCDF variable %s is missing!", varname); \
1084 } \
1085 }
1086
1103#define NC_INQ_DIM(dimname, ptr, min, max) { \
1104 int dimid; size_t naux; \
1105 NC(nc_inq_dimid(ncid, dimname, &dimid)); \
1106 NC(nc_inq_dimlen(ncid, dimid, &naux)); \
1107 *ptr = (int)naux; \
1108 if ((*ptr) < (min) || (*ptr) > (max)) \
1109 ERRMSG("Dimension %s is out of range!", dimname); \
1110 }
1111
1126#define NC_PUT_DOUBLE(varname, ptr, hyperslab) { \
1127 NC(nc_inq_varid(ncid, varname, &varid)); \
1128 if(hyperslab) { \
1129 NC(nc_put_vara_double(ncid, varid, start, count, ptr)); \
1130 } else { \
1131 NC(nc_put_var_double(ncid, varid, ptr)); \
1132 } \
1133 }
1134
1150#define NC_PUT_FLOAT(varname, ptr, hyperslab) { \
1151 NC(nc_inq_varid(ncid, varname, &varid)); \
1152 if(hyperslab) { \
1153 NC(nc_put_vara_float(ncid, varid, start, count, ptr)); \
1154 } else { \
1155 NC(nc_put_var_float(ncid, varid, ptr)); \
1156 } \
1157 }
1158
1173#define NC_PUT_INT(varname, ptr, hyperslab) { \
1174 NC(nc_inq_varid(ncid, varname, &varid)); \
1175 if(hyperslab) { \
1176 NC(nc_put_vara_int(ncid, varid, start, count, ptr)); \
1177 } else { \
1178 NC(nc_put_var_int(ncid, varid, ptr)); \
1179 } \
1180 }
1181
1195#define NC_PUT_ATT(varname, attname, text) { \
1196 NC(nc_inq_varid(ncid, varname, &varid)); \
1197 NC(nc_put_att_text(ncid, varid, attname, strnlen(text, LEN), text)); \
1198 }
1199
1212#define NC_PUT_ATT_GLOBAL(attname, text) \
1213 NC(nc_put_att_text(ncid, NC_GLOBAL, attname, strnlen(text, LEN), text));
1214
1232#define NN(x0, y0, x1, y1, x) \
1233 (fabs((x) - (x0)) <= fabs((x) - (x1)) ? (y0) : (y1))
1234
1247#define NORM(a) \
1248 sqrt(DOTP(a, a))
1249
1265#ifdef _OPENACC
1266#define PARTICLE_LOOP(ip0, ip1, check_dt, ...) \
1267 const int ip0_const = ip0; \
1268 const int ip1_const = ip1; \
1269 _Pragma(__VA_ARGS__) \
1270 _Pragma("acc parallel loop independent gang vector") \
1271 for (int ip = ip0_const; ip < ip1_const; ip++) \
1272 if (!check_dt || cache->dt[ip] != 0)
1273#else
1274#define PARTICLE_LOOP(ip0, ip1, check_dt, ...) \
1275 const int ip0_const = ip0; \
1276 const int ip1_const = ip1; \
1277 _Pragma("omp parallel for default(shared)") \
1278 for (int ip = ip0_const; ip < ip1_const; ip++) \
1279 if (!check_dt || cache->dt[ip] != 0)
1280#endif
1281
1304#define P(z) \
1305 (P0 * exp(-(z) / H0))
1306
1328#define PSAT(t) \
1329 (6.112 * exp(17.62 * ((t) - T0) / (243.12 + (t) - T0)))
1330
1352#define PSICE(t) \
1353 (6.112 * exp(22.46 * ((t) - T0) / (272.62 + (t) - T0)))
1354
1379#define PW(p, h2o) \
1380 ((p) * MAX((h2o), 0.1e-6) / (1. + (1. - EPS) * MAX((h2o), 0.1e-6)))
1381
1396#define RAD2DEG(rad) \
1397 ((rad) * (180.0 / M_PI))
1398
1426#define RH(p, t, h2o) \
1427 (PW(p, h2o) / PSAT(t) * 100.)
1428
1456#define RHICE(p, t, h2o) \
1457 (PW(p, h2o) / PSICE(t) * 100.)
1458
1481#define RHO(p, t) \
1482 (100. * (p) / (RA * (t)))
1483
1500#define SET_ATM(qnt, val) \
1501 if (ctl->qnt >= 0) \
1502 atm->q[ctl->qnt][ip] = val;
1503
1523#define SET_QNT(qnt, name, longname, unit) \
1524 if (strcasecmp(ctl->qnt_name[iq], name) == 0) { \
1525 ctl->qnt = iq; \
1526 sprintf(ctl->qnt_longname[iq], longname); \
1527 sprintf(ctl->qnt_unit[iq], unit); \
1528 } else
1529
1544#define SH(h2o) \
1545 (EPS * MAX((h2o), 0.1e-6))
1546
1557#define SQR(x) \
1558 ((x)*(x))
1559
1571#define SWAP(x, y, type) \
1572 do {type tmp = x; x = y; y = tmp;} while(0);
1573
1595#define TDEW(p, h2o) \
1596 (T0 + 243.12 * log(PW((p), (h2o)) / 6.112) \
1597 / (17.62 - log(PW((p), (h2o)) / 6.112)))
1598
1620#define TICE(p, h2o) \
1621 (T0 + 272.62 * log(PW((p), (h2o)) / 6.112) \
1622 / (22.46 - log(PW((p), (h2o)) / 6.112)))
1623
1644#define THETA(p, t) \
1645 ((t) * pow(1000. / (p), 0.286))
1646
1673#define THETAVIRT(p, t, h2o) \
1674 (TVIRT(THETA((p), (t)), MAX((h2o), 0.1e-6)))
1675
1694#define TOK(line, tok, format, var) { \
1695 if(((tok)=strtok((line), " \t"))) { \
1696 if(sscanf(tok, format, &(var))!=1) continue; \
1697 } else ERRMSG("Error while reading!"); \
1698 }
1699
1719#define TVIRT(t, h2o) \
1720 ((t) * (1. + (1. - EPS) * MAX((h2o), 0.1e-6)))
1721
1741#define Z(p) \
1742 (H0 * log(P0 / (p)))
1743
1772#define ZDIFF(lnp0, t0, h2o0, lnp1, t1, h2o1) \
1773 (RI / MA / G0 * 0.5 * (TVIRT((t0), (h2o0)) + TVIRT((t1), (h2o1))) \
1774 * ((lnp0) - (lnp1)))
1775
1791#define ZETA(ps, p, t) \
1792 (((p) / (ps) <= 0.3 ? 1. : \
1793 sin(M_PI / 2. * (1. - (p) / (ps)) / (1. - 0.3))) \
1794 * THETA((p), (t)))
1795
1796/* ------------------------------------------------------------
1797 Log messages...
1798 ------------------------------------------------------------ */
1799
1801#ifndef LOGLEV
1802#define LOGLEV 2
1803#endif
1804
1834#define LOG(level, ...) { \
1835 if(level >= 2) \
1836 printf(" "); \
1837 if(level <= LOGLEV) { \
1838 printf(__VA_ARGS__); \
1839 printf("\n"); \
1840 } \
1841 }
1842
1871#define WARN(...) { \
1872 printf("\nWarning (%s, %s, l%d): ", __FILE__, __func__, __LINE__); \
1873 LOG(0, __VA_ARGS__); \
1874 }
1875
1904#define ERRMSG(...) { \
1905 printf("\nError (%s, %s, l%d): ", __FILE__, __func__, __LINE__); \
1906 LOG(0, __VA_ARGS__); \
1907 exit(EXIT_FAILURE); \
1908 }
1909
1939#define PRINT(format, var) \
1940 printf("Print (%s, %s, l%d): %s= "format"\n", \
1941 __FILE__, __func__, __LINE__, #var, var);
1942
1943/* ------------------------------------------------------------
1944 Timers...
1945 ------------------------------------------------------------ */
1946
1948#define NTIMER 100
1949
1963#define PRINT_TIMERS \
1964 timer("END", "END", 1);
1965
1984#define SELECT_TIMER(id, group, color) { \
1985 NVTX_POP; \
1986 NVTX_PUSH(id, color); \
1987 timer(id, group, 0); \
1988 }
1989
2003#define START_TIMERS \
2004 NVTX_PUSH("START", NVTX_CPU);
2005
2018#define STOP_TIMERS \
2019 NVTX_POP;
2020
2021/* ------------------------------------------------------------
2022 NVIDIA Tools Extension (NVTX)...
2023 ------------------------------------------------------------ */
2024
2025#ifdef NVTX
2026#include "nvToolsExt.h"
2027
2029#define NVTX_CPU 0xFFADD8E6
2030
2032#define NVTX_GPU 0xFF00008B
2033
2035#define NVTX_H2D 0xFFFFFF00
2036
2038#define NVTX_D2H 0xFFFF8800
2039
2041#define NVTX_READ 0xFFFFCCCB
2042
2044#define NVTX_WRITE 0xFF8B0000
2045
2047#define NVTX_RECV 0xFFCCFFCB
2048
2050#define NVTX_SEND 0xFF008B00
2051
2081#define NVTX_PUSH(range_title, range_color) { \
2082 nvtxEventAttributes_t eventAttrib = {0}; \
2083 eventAttrib.version = NVTX_VERSION; \
2084 eventAttrib.size = NVTX_EVENT_ATTRIB_STRUCT_SIZE; \
2085 eventAttrib.messageType = NVTX_MESSAGE_TYPE_ASCII; \
2086 eventAttrib.colorType = NVTX_COLOR_ARGB; \
2087 eventAttrib.color = range_color; \
2088 eventAttrib.message.ascii = range_title; \
2089 nvtxRangePushEx(&eventAttrib); \
2090 }
2091
2104#define NVTX_POP { \
2105 nvtxRangePop(); \
2106 }
2107#else
2108
2109/* Empty definitions of NVTX_PUSH and NVTX_POP... */
2110#define NVTX_PUSH(range_title, range_color) {}
2111#define NVTX_POP {}
2112#endif
2113
2114/* ------------------------------------------------------------
2115 Thrust...
2116 ------------------------------------------------------------ */
2117
2143 double *__restrict__ c,
2144 int n,
2145 int *__restrict__ index);
2146
2147/* ------------------------------------------------------------
2148 Structs...
2149 ------------------------------------------------------------ */
2150
2158typedef struct {
2159
2160 /* ------------------------------------------------------------
2161 Quantity parameters...
2162 ------------------------------------------------------------ */
2163
2165 int nq;
2166
2168 char qnt_name[NQ][LEN];
2169
2171 char qnt_longname[NQ][LEN];
2172
2174 char qnt_unit[NQ][LEN];
2175
2177 char qnt_format[NQ][LEN];
2178
2181
2184
2187
2190
2193
2196
2199
2202
2205
2208
2211
2214
2217
2220
2223
2226
2229
2232
2235
2238
2241
2244
2247
2250
2253
2256
2259
2262
2265
2268
2271
2274
2277
2280
2283
2286
2289
2292
2295
2298
2301
2304
2307
2310
2313
2316
2319
2322
2325
2328
2331
2334
2337
2340
2343
2346
2349
2352
2355
2358
2361
2364
2367
2370
2373
2376
2379
2382
2385
2388
2391
2394
2397
2400
2403
2406
2409
2412
2415
2418
2421
2424
2427
2430
2433
2436
2439
2442
2445
2448
2451
2454
2456 double t_start;
2457
2459 double t_stop;
2460
2462 double dt_mod;
2463
2464 /* ------------------------------------------------------------
2465 Meteo data parameters...
2466 ------------------------------------------------------------ */
2467
2469 char metbase[LEN];
2470
2472 double dt_met;
2473
2476
2480
2484
2487
2490
2493
2496
2499
2502
2505
2508
2511
2515
2518
2521
2524
2527
2530
2533
2536
2539
2542
2545
2548
2551
2554
2557
2560
2563
2566
2569
2572
2575
2578
2580 double met_p[EP];
2581
2584
2587
2589 double met_lev_hyam[EP];
2590
2592 double met_lev_hybm[EP];
2593
2596
2599
2602
2605
2608
2611
2614
2618
2621
2624
2627
2630
2633
2636
2637 /* ------------------------------------------------------------
2638 Geophysical module parameters...
2639 ------------------------------------------------------------ */
2640
2642 double sort_dt;
2643
2647
2649 char balloon[LEN];
2650
2653
2657
2660
2663
2666
2669
2672
2675
2678
2681
2684
2687
2690
2693
2696
2698 double conv_cin;
2699
2701 double conv_dt;
2702
2705
2708
2711
2714
2717
2720
2722 double bound_p0;
2723
2725 double bound_p1;
2726
2729
2732
2735
2738
2740 char species[LEN];
2741
2743 double molmass;
2744
2747
2750
2753
2755 char clim_hno3_filename[LEN];
2756
2758 char clim_oh_filename[LEN];
2759
2761 char clim_h2o2_filename[LEN];
2762
2764 char clim_ho2_filename[LEN];
2765
2767 char clim_o1d_filename[LEN];
2768
2770 char clim_o3_filename[LEN];
2771
2773 char clim_ccl4_timeseries[LEN];
2774
2776 char clim_ccl3f_timeseries[LEN];
2777
2779 char clim_ccl2f2_timeseries[LEN];
2780
2782 char clim_n2o_timeseries[LEN];
2783
2785 char clim_sf6_timeseries[LEN];
2786
2789
2792
2795
2798
2801
2804
2807
2810
2813
2816
2819
2822
2825
2828
2831
2834
2837
2840
2843
2846
2849
2852
2854 double oh_chem[4];
2855
2858
2861
2864
2866 double dt_kpp;
2867
2870
2872 double wet_depo_pre[2];
2873
2876
2879
2882
2885
2887 double wet_depo_ic_h[2];
2888
2890 double wet_depo_bc_h[2];
2891
2894
2897
2900
2903
2906
2908 double psc_h2o;
2909
2911 double psc_hno3;
2912
2913 /* ------------------------------------------------------------
2914 Output parameters...
2915 ------------------------------------------------------------ */
2916
2918 char atm_basename[LEN];
2919
2921 char atm_gpfile[LEN];
2922
2925
2928
2931
2935
2940
2943
2945 int atm_nc_quant[NQ];
2946
2949
2951 char csi_basename[LEN];
2952
2954 char csi_kernel[LEN];
2955
2958
2960 char csi_obsfile[LEN];
2961
2964
2967
2970
2972 double csi_z0;
2973
2975 double csi_z1;
2976
2979
2981 double csi_lon0;
2982
2984 double csi_lon1;
2985
2988
2990 double csi_lat0;
2991
2993 double csi_lat1;
2994
2996 char ens_basename[LEN];
2997
3000
3002 char grid_basename[LEN];
3003
3005 char grid_kernel[LEN];
3006
3008 char grid_gpfile[LEN];
3009
3012
3015
3018
3020 int grid_nc_quant[NQ];
3021
3024
3027
3029 double grid_z0;
3030
3032 double grid_z1;
3033
3036
3039
3042
3045
3048
3051
3054
3056 char prof_basename[LEN];
3057
3059 char prof_obsfile[LEN];
3060
3063
3065 double prof_z0;
3066
3068 double prof_z1;
3069
3072
3075
3078
3081
3084
3087
3089 char sample_basename[LEN];
3090
3092 char sample_kernel[LEN];
3093
3095 char sample_obsfile[LEN];
3096
3099
3102
3104 char stat_basename[LEN];
3105
3107 double stat_lon;
3108
3110 double stat_lat;
3111
3113 double stat_r;
3114
3116 double stat_t0;
3117
3119 double stat_t1;
3120
3122 char vtk_basename[LEN];
3123
3126
3129
3132
3135
3138
3139} ctl_t;
3140
3149typedef struct {
3150
3152 int np;
3153
3155 double time[NP];
3156
3158 double p[NP];
3159
3161 double lon[NP];
3162
3164 double lat[NP];
3165
3167 double q[NQ][NP];
3168
3169} atm_t;
3170
3177typedef struct {
3178
3180 double iso_var[NP];
3181
3183 double iso_ps[NP];
3184
3186 double iso_ts[NP];
3187
3190
3192 float uvwp[NP][3];
3193
3195 double rs[3 * NP + 1];
3196
3198 double dt[NP];
3199
3200} cache_t;
3201
3209typedef struct {
3210
3212 int np;
3213
3215 int nsza;
3216
3218 int no3c;
3219
3221 double p[CP];
3222
3224 double sza[CSZA];
3225
3227 double o3c[CO3];
3228
3230 double n2o[CP][CSZA][CO3];
3231
3233 double ccl4[CP][CSZA][CO3];
3234
3236 double ccl3f[CP][CSZA][CO3];
3237
3239 double ccl2f2[CP][CSZA][CO3];
3240
3242 double o2[CP][CSZA][CO3];
3243
3245 double o3_1[CP][CSZA][CO3];
3246
3248 double o3_2[CP][CSZA][CO3];
3249
3251 double h2o2[CP][CSZA][CO3];
3252
3254 double h2o[CP][CSZA][CO3];
3255
3256} clim_photo_t;
3257
3265typedef struct {
3266
3269
3271 double time[CTS];
3272
3274 double vmr[CTS];
3275
3276} clim_ts_t;
3277
3285typedef struct {
3286
3289
3291 int nlat;
3292
3294 int np;
3295
3297 double time[CT];
3298
3300 double lat[CY];
3301
3303 double p[CP];
3304
3306 double vmr[CT][CP][CY];
3307
3308} clim_zm_t;
3309
3317typedef struct {
3318
3321
3324
3326 double tropo_time[12];
3327
3329 double tropo_lat[73];
3330
3332 double tropo[12][73];
3333
3336
3339
3342
3345
3348
3351
3354
3357
3360
3363
3366
3367} clim_t;
3368
3376typedef struct {
3377
3379 double time;
3380
3382 int nx;
3383
3385 int ny;
3386
3388 int np;
3389
3391 int npl;
3392
3394 double lon[EX];
3395
3397 double lat[EY];
3398
3400 double p[EP];
3401
3403 double hybrid[EP];
3404
3406 float ps[EX][EY];
3407
3409 float ts[EX][EY];
3410
3412 float zs[EX][EY];
3413
3415 float us[EX][EY];
3416
3418 float vs[EX][EY];
3419
3421 float ess[EX][EY];
3422
3424 float nss[EX][EY];
3425
3427 float shf[EX][EY];
3428
3430 float lsm[EX][EY];
3431
3433 float sst[EX][EY];
3434
3436 float pbl[EX][EY];
3437
3439 float pt[EX][EY];
3440
3442 float tt[EX][EY];
3443
3445 float zt[EX][EY];
3446
3448 float h2ot[EX][EY];
3449
3451 float pct[EX][EY];
3452
3454 float pcb[EX][EY];
3455
3457 float cl[EX][EY];
3458
3460 float plcl[EX][EY];
3461
3463 float plfc[EX][EY];
3464
3466 float pel[EX][EY];
3467
3469 float cape[EX][EY];
3470
3472 float cin[EX][EY];
3473
3475 float o3c[EX][EY];
3476
3478 float z[EX][EY][EP];
3479
3481 float t[EX][EY][EP];
3482
3484 float u[EX][EY][EP];
3485
3487 float v[EX][EY][EP];
3488
3490 float w[EX][EY][EP];
3491
3493 float pv[EX][EY][EP];
3494
3496 float h2o[EX][EY][EP];
3497
3499 float o3[EX][EY][EP];
3500
3502 float lwc[EX][EY][EP];
3503
3505 float rwc[EX][EY][EP];
3506
3508 float iwc[EX][EY][EP];
3509
3511 float swc[EX][EY][EP];
3512
3514 float cc[EX][EY][EP];
3515
3517 float pl[EX][EY][EP];
3518
3520 float ul[EX][EY][EP];
3521
3523 float vl[EX][EY][EP];
3524
3526 float wl[EX][EY][EP];
3527
3529 float zetal[EX][EY][EP];
3530
3532 float zeta_dotl[EX][EY][EP];
3533
3534} met_t;
3535
3536/* ------------------------------------------------------------
3537 OpenACC routines...
3538 ------------------------------------------------------------ */
3539
3540#ifdef _OPENACC
3541#pragma acc routine (clim_oh)
3542#pragma acc routine (clim_photo)
3543#pragma acc routine (clim_tropo)
3544#pragma acc routine (clim_ts)
3545#pragma acc routine (clim_zm)
3546#pragma acc routine (intpol_check_lon_lat)
3547#pragma acc routine (intpol_met_4d_coord)
3548#pragma acc routine (intpol_met_space_3d)
3549#pragma acc routine (intpol_met_space_3d_ml)
3550#pragma acc routine (intpol_met_space_2d)
3551#pragma acc routine (intpol_met_time_3d)
3552#pragma acc routine (intpol_met_time_3d_ml)
3553#pragma acc routine (intpol_met_time_2d)
3554#pragma acc routine (kernel_weight)
3555#pragma acc routine (lapse_rate)
3556#pragma acc routine (locate_irr)
3557#pragma acc routine (locate_irr_float)
3558#pragma acc routine (locate_reg)
3559#pragma acc routine (locate_vert)
3560#pragma acc routine (nat_temperature)
3561#pragma acc routine (pbl_weight)
3562#pragma acc routine (sedi)
3563#pragma acc routine (stddev)
3564#pragma acc routine (sza_calc)
3565#pragma acc routine (tropo_weight)
3566#endif
3567
3568/* ------------------------------------------------------------
3569 Functions...
3570 ------------------------------------------------------------ */
3571
3594 void *data,
3595 size_t N);
3596
3611void cart2geo(
3612 const double *x,
3613 double *z,
3614 double *lon,
3615 double *lat);
3616
3639double clim_oh(
3640 const ctl_t * ctl,
3641 const clim_t * clim,
3642 const double t,
3643 const double lon,
3644 const double lat,
3645 const double p);
3646
3666 const ctl_t * ctl,
3667 clim_t * clim);
3668
3698double clim_photo(
3699 const double rate[CP][CSZA][CO3],
3700 const clim_photo_t * photo,
3701 const double p,
3702 const double sza,
3703 const double o3c);
3704
3730double clim_tropo(
3731 const clim_t * clim,
3732 const double t,
3733 const double lat);
3734
3753void clim_tropo_init(
3754 clim_t * clim);
3755
3772double clim_ts(
3773 const clim_ts_t * ts,
3774 const double t);
3775
3797double clim_zm(
3798 const clim_zm_t * zm,
3799 const double t,
3800 const double lat,
3801 const double p);
3802
3844 const ctl_t * ctl,
3845 const char *varname,
3846 float *array,
3847 const size_t nx,
3848 const size_t ny,
3849 const size_t np,
3850 const int decompress,
3851 FILE * inout);
3852
3885void compress_pck(
3886 const char *varname,
3887 float *array,
3888 const size_t nxy,
3889 const size_t nz,
3890 const int decompress,
3891 FILE * inout);
3892
3932 const char *varname,
3933 float *array,
3934 const int nx,
3935 const int ny,
3936 const int nz,
3937 const int precision,
3938 const double tolerance,
3939 const int decompress,
3940 FILE * inout);
3941
3978 const char *varname,
3979 float *array,
3980 const size_t n,
3981 const int decompress,
3982 FILE * inout);
3983
4006void day2doy(
4007 const int year,
4008 const int mon,
4009 const int day,
4010 int *doy);
4011
4033void doy2day(
4034 const int year,
4035 const int doy,
4036 int *mon,
4037 int *day);
4038
4065void fft_help(
4066 double *fcReal,
4067 double *fcImag,
4068 const int n);
4069
4096void geo2cart(
4097 const double z,
4098 const double lon,
4099 const double lat,
4100 double *x);
4101
4126void get_met_help(
4127 const ctl_t * ctl,
4128 const double t,
4129 const int direct,
4130 const char *metbase,
4131 const double dt_met,
4132 char *filename);
4133
4157void get_met_replace(
4158 char *orig,
4159 char *search,
4160 char *repl);
4161
4198void get_tropo(
4199 const int met_tropo,
4200 ctl_t * ctl,
4201 clim_t * clim,
4202 met_t * met,
4203 const double *lons,
4204 const int nx,
4205 const double *lats,
4206 const int ny,
4207 double *pt,
4208 double *zt,
4209 double *tt,
4210 double *qt,
4211 double *o3t,
4212 double *ps,
4213 double *zs);
4214
4236 const double *lons,
4237 const int nlon,
4238 const double *lats,
4239 const int nlat,
4240 const double lon,
4241 const double lat,
4242 double *lon2,
4243 double *lat2);
4244
4287 const met_t * met0,
4288 float height0[EX][EY][EP],
4289 float array0[EX][EY][EP],
4290 const met_t * met1,
4291 float height1[EX][EY][EP],
4292 float array1[EX][EY][EP],
4293 const double ts,
4294 const double height,
4295 const double lon,
4296 const double lat,
4297 double *var,
4298 int *ci,
4299 double *cw,
4300 const int init);
4301
4337 const met_t * met,
4338 float array[EX][EY][EP],
4339 const double p,
4340 const double lon,
4341 const double lat,
4342 double *var,
4343 int *ci,
4344 double *cw,
4345 const int init);
4346
4366 const met_t * met,
4367 float zs[EX][EY][EP],
4368 float vals[EX][EY][EP],
4369 const double z,
4370 const double lon,
4371 const double lat,
4372 double *val);
4373
4409 const met_t * met,
4410 float array[EX][EY],
4411 const double lon,
4412 const double lat,
4413 double *var,
4414 int *ci,
4415 double *cw,
4416 const int init);
4417
4452 const met_t * met0,
4453 float array0[EX][EY][EP],
4454 const met_t * met1,
4455 float array1[EX][EY][EP],
4456 const double ts,
4457 const double p,
4458 const double lon,
4459 const double lat,
4460 double *var,
4461 int *ci,
4462 double *cw,
4463 const int init);
4464
4488 const met_t * met0,
4489 float zs0[EX][EY][EP],
4490 float array0[EX][EY][EP],
4491 const met_t * met1,
4492 float zs1[EX][EY][EP],
4493 float array1[EX][EY][EP],
4494 const double ts,
4495 const double p,
4496 const double lon,
4497 const double lat,
4498 double *var);
4499
4535 const met_t * met0,
4536 float array0[EX][EY],
4537 const met_t * met1,
4538 float array1[EX][EY],
4539 const double ts,
4540 const double lon,
4541 const double lat,
4542 double *var,
4543 int *ci,
4544 double *cw,
4545 const int init);
4546
4584void intpol_tropo_3d(
4585 const double time0,
4586 float array0[EX][EY],
4587 const double time1,
4588 float array1[EX][EY],
4589 const double lons[EX],
4590 const double lats[EY],
4591 const int nlon,
4592 const int nlat,
4593 const double time,
4594 const double lon,
4595 const double lat,
4596 const int method,
4597 double *var,
4598 double *sigma);
4599
4626void jsec2time(
4627 const double jsec,
4628 int *year,
4629 int *mon,
4630 int *day,
4631 int *hour,
4632 int *min,
4633 int *sec,
4634 double *remain);
4635
4662double kernel_weight(
4663 const double kz[EP],
4664 const double kw[EP],
4665 const int nk,
4666 const double p);
4667
4706double lapse_rate(
4707 const double t,
4708 const double h2o);
4709
4739 ctl_t * ctl);
4740
4760int locate_irr(
4761 const double *xx,
4762 const int n,
4763 const double x);
4764
4791 const float *xx,
4792 const int n,
4793 const double x,
4794 const int ig);
4795
4816int locate_reg(
4817 const double *xx,
4818 const int n,
4819 const double x);
4820
4842void locate_vert(
4843 float profiles[EX][EY][EP],
4844 const int np,
4845 const int lon_ap_ind,
4846 const int lat_ap_ind,
4847 const double alt_ap,
4848 int *ind);
4849
4875void module_advect(
4876 const ctl_t * ctl,
4877 const cache_t * cache,
4878 met_t * met0,
4879 met_t * met1,
4880 atm_t * atm);
4881
4905 const ctl_t * ctl,
4906 const cache_t * cache,
4907 met_t * met0,
4908 met_t * met1,
4909 atm_t * atm);
4910
4948 const ctl_t * ctl,
4949 const cache_t * cache,
4950 const clim_t * clim,
4951 met_t * met0,
4952 met_t * met1,
4953 atm_t * atm);
4954
4972void module_chem_grid(
4973 const ctl_t * ctl,
4974 met_t * met0,
4975 met_t * met1,
4976 atm_t * atm,
4977 const double t);
4978
5005void module_chem_init(
5006 const ctl_t * ctl,
5007 const cache_t * cache,
5008 const clim_t * clim,
5009 met_t * met0,
5010 met_t * met1,
5011 atm_t * atm);
5012
5037 const ctl_t * ctl,
5038 cache_t * cache,
5039 met_t * met0,
5040 met_t * met1,
5041 atm_t * atm);
5042
5069void module_decay(
5070 const ctl_t * ctl,
5071 const cache_t * cache,
5072 const clim_t * clim,
5073 atm_t * atm);
5074
5111void module_diff_meso(
5112 const ctl_t * ctl,
5113 cache_t * cache,
5114 met_t * met0,
5115 met_t * met1,
5116 atm_t * atm);
5117
5151void module_diff_pbl(
5152 const ctl_t * ctl,
5153 cache_t * cache,
5154 met_t * met0,
5155 met_t * met1,
5156 atm_t * atm);
5157
5212void module_diff_turb(
5213 const ctl_t * ctl,
5214 cache_t * cache,
5215 const clim_t * clim,
5216 met_t * met0,
5217 met_t * met1,
5218 atm_t * atm);
5219
5239void module_dry_depo(
5240 const ctl_t * ctl,
5241 const cache_t * cache,
5242 met_t * met0,
5243 met_t * met1,
5244 atm_t * atm);
5245
5278void module_h2o2_chem(
5279 const ctl_t * ctl,
5280 const cache_t * cache,
5281 const clim_t * clim,
5282 met_t * met0,
5283 met_t * met1,
5284 atm_t * atm);
5285
5306 const ctl_t * ctl,
5307 cache_t * cache,
5308 met_t * met0,
5309 met_t * met1,
5310 atm_t * atm);
5311
5329void module_isosurf(
5330 const ctl_t * ctl,
5331 const cache_t * cache,
5332 met_t * met0,
5333 met_t * met1,
5334 atm_t * atm);
5335
5368 ctl_t * ctl,
5369 cache_t * cache,
5370 clim_t * clim,
5371 met_t * met0,
5372 met_t * met1,
5373 atm_t * atm);
5374
5393void module_meteo(
5394 const ctl_t * ctl,
5395 const cache_t * cache,
5396 const clim_t * clim,
5397 met_t * met0,
5398 met_t * met1,
5399 atm_t * atm);
5400
5418void module_mixing(
5419 const ctl_t * ctl,
5420 const clim_t * clim,
5421 atm_t * atm,
5422 const double t);
5423
5444 const ctl_t * ctl,
5445 const clim_t * clim,
5446 atm_t * atm,
5447 const int *ixs,
5448 const int *iys,
5449 const int *izs,
5450 const int qnt_idx);
5451
5484void module_oh_chem(
5485 const ctl_t * ctl,
5486 const cache_t * cache,
5487 const clim_t * clim,
5488 met_t * met0,
5489 met_t * met1,
5490 atm_t * atm);
5491
5519void module_position(
5520 const cache_t * cache,
5521 met_t * met0,
5522 met_t * met1,
5523 atm_t * atm);
5524
5549void module_rng_init(
5550 const int ntask);
5551
5577void module_rng(
5578 const ctl_t * ctl,
5579 double *rs,
5580 const size_t n,
5581 const int method);
5582
5605void module_sedi(
5606 const ctl_t * ctl,
5607 const cache_t * cache,
5608 met_t * met0,
5609 met_t * met1,
5610 atm_t * atm);
5611
5635void module_sort(
5636 const ctl_t * ctl,
5637 met_t * met0,
5638 atm_t * atm);
5639
5659void module_sort_help(
5660 double *a,
5661 const int *p,
5662 const int np);
5663
5687void module_timesteps(
5688 const ctl_t * ctl,
5689 cache_t * cache,
5690 met_t * met0,
5691 atm_t * atm,
5692 const double t);
5693
5715 ctl_t * ctl,
5716 const atm_t * atm);
5717
5751 const ctl_t * ctl,
5752 const cache_t * cache,
5753 const clim_t * clim,
5754 met_t * met0,
5755 met_t * met1,
5756 atm_t * atm);
5757
5787void module_wet_depo(
5788 const ctl_t * ctl,
5789 const cache_t * cache,
5790 met_t * met0,
5791 met_t * met1,
5792 atm_t * atm);
5793
5824void mptrac_alloc(
5825 ctl_t ** ctl,
5826 cache_t ** cache,
5827 clim_t ** clim,
5828 met_t ** met0,
5829 met_t ** met1,
5830 atm_t ** atm);
5831
5861void mptrac_free(
5862 ctl_t * ctl,
5863 cache_t * cache,
5864 clim_t * clim,
5865 met_t * met0,
5866 met_t * met1,
5867 atm_t * atm);
5868
5903void mptrac_get_met(
5904 ctl_t * ctl,
5905 clim_t * clim,
5906 const double t,
5907 met_t ** met0,
5908 met_t ** met1);
5909
5929void mptrac_init(
5930 ctl_t * ctl,
5931 cache_t * cache,
5932 clim_t * clim,
5933 atm_t * atm,
5934 const int ntask);
5935
5971int mptrac_read_atm(
5972 const char *filename,
5973 const ctl_t * ctl,
5974 atm_t * atm);
5975
6007void mptrac_read_clim(
6008 const ctl_t * ctl,
6009 clim_t * clim);
6010
6040void mptrac_read_ctl(
6041 const char *filename,
6042 int argc,
6043 char *argv[],
6044 ctl_t * ctl);
6045
6074int mptrac_read_met(
6075 const char *filename,
6076 const ctl_t * ctl,
6077 const clim_t * clim,
6078 met_t * met);
6079
6100 ctl_t * ctl,
6101 cache_t * cache,
6102 clim_t * clim,
6103 met_t ** met0,
6104 met_t ** met1,
6105 atm_t * atm,
6106 double t);
6107
6137void mptrac_write_atm(
6138 const char *filename,
6139 const ctl_t * ctl,
6140 const atm_t * atm,
6141 const double t);
6142
6178void mptrac_write_met(
6179 const char *filename,
6180 const ctl_t * ctl,
6181 met_t * met);
6182
6217 const char *dirname,
6218 const ctl_t * ctl,
6219 met_t * met0,
6220 met_t * met1,
6221 atm_t * atm,
6222 const double t);
6223
6255 const ctl_t * ctl,
6256 const cache_t * cache,
6257 const clim_t * clim,
6258 met_t ** met0,
6259 met_t ** met1,
6260 const atm_t * atm);
6261
6292 const ctl_t * ctl,
6293 const cache_t * cache,
6294 const clim_t * clim,
6295 met_t ** met0,
6296 met_t ** met1,
6297 const atm_t * atm);
6298
6326double nat_temperature(
6327 const double p,
6328 const double h2o,
6329 const double hno3);
6330
6351double pbl_weight(
6352 const ctl_t * ctl,
6353 const atm_t * atm,
6354 const int ip,
6355 const double pbl,
6356 const double ps);
6357
6390int read_atm_asc(
6391 const char *filename,
6392 const ctl_t * ctl,
6393 atm_t * atm);
6394
6425int read_atm_bin(
6426 const char *filename,
6427 const ctl_t * ctl,
6428 atm_t * atm);
6429
6454int read_atm_clams(
6455 const char *filename,
6456 const ctl_t * ctl,
6457 atm_t * atm);
6458
6488int read_atm_nc(
6489 const char *filename,
6490 const ctl_t * ctl,
6491 atm_t * atm);
6492
6521void read_clim_photo(
6522 const char *filename,
6523 clim_photo_t * photo);
6524
6542 const int ncid,
6543 const char *varname,
6544 const clim_photo_t * photo,
6545 double var[CP][CSZA][CO3]);
6546
6570int read_clim_ts(
6571 const char *filename,
6572 clim_ts_t * ts);
6573
6600void read_clim_zm(
6601 const char *filename,
6602 const char *varname,
6603 clim_zm_t * zm);
6604
6632void read_kernel(
6633 const char *filename,
6634 double kz[EP],
6635 double kw[EP],
6636 int *nk);
6637
6669int read_met_bin(
6670 const char *filename,
6671 const ctl_t * ctl,
6672 met_t * met);
6673
6699void read_met_bin_2d(
6700 FILE * in,
6701 const met_t * met,
6702 float var[EX][EY],
6703 const char *varname);
6704
6742void read_met_bin_3d(
6743 FILE * in,
6744 const ctl_t * ctl,
6745 const met_t * met,
6746 float var[EX][EY][EP],
6747 const char *varname,
6748 const float bound_min,
6749 const float bound_max);
6750
6777void read_met_cape(
6778 const ctl_t * ctl,
6779 const clim_t * clim,
6780 met_t * met);
6781
6804void read_met_cloud(
6805 met_t * met);
6806
6832void read_met_detrend(
6833 const ctl_t * ctl,
6834 met_t * met);
6835
6859 met_t * met);
6860
6887void read_met_geopot(
6888 const ctl_t * ctl,
6889 met_t * met);
6890
6922void read_met_grid(
6923 const char *filename,
6924 const int ncid,
6925 const ctl_t * ctl,
6926 met_t * met);
6927
6958void read_met_levels(
6959 const int ncid,
6960 const ctl_t * ctl,
6961 met_t * met);
6962
6991void read_met_ml2pl(
6992 const ctl_t * ctl,
6993 const met_t * met,
6994 float var[EX][EY][EP],
6995 const char *varname);
6996
7019 const ctl_t * ctl,
7020 met_t * met);
7021
7053int read_met_nc(
7054 const char *filename,
7055 const ctl_t * ctl,
7056 const clim_t * clim,
7057 met_t * met);
7058
7088int read_met_nc_2d(
7089 const int ncid,
7090 const char *varname,
7091 const char *varname2,
7092 const char *varname3,
7093 const char *varname4,
7094 const char *varname5,
7095 const char *varname6,
7096 const ctl_t * ctl,
7097 const met_t * met,
7098 float dest[EX][EY],
7099 const float scl,
7100 const int init);
7101
7132int read_met_nc_3d(
7133 const int ncid,
7134 const char *varname,
7135 const char *varname2,
7136 const char *varname3,
7137 const char *varname4,
7138 const ctl_t * ctl,
7139 const met_t * met,
7140 float dest[EX][EY][EP],
7141 const float scl);
7142
7188void read_met_pbl(
7189 const ctl_t * ctl,
7190 met_t * met);
7191
7224 met_t * met);
7225
7256 met_t * met);
7257
7288void read_met_pv(
7289 met_t * met);
7290
7313void read_met_ozone(
7314 met_t * met);
7315
7344void read_met_sample(
7345 const ctl_t * ctl,
7346 met_t * met);
7347
7388void read_met_surface(
7389 const int ncid,
7390 const ctl_t * ctl,
7391 met_t * met);
7392
7421void read_met_tropo(
7422 const ctl_t * ctl,
7423 const clim_t * clim,
7424 met_t * met);
7425
7457void read_obs(
7458 const char *filename,
7459 const ctl_t * ctl,
7460 double *rt,
7461 double *rz,
7462 double *rlon,
7463 double *rlat,
7464 double *robs,
7465 int *nobs);
7466
7494void read_obs_asc(
7495 const char *filename,
7496 double *rt,
7497 double *rz,
7498 double *rlon,
7499 double *rlat,
7500 double *robs,
7501 int *nobs);
7502
7529void read_obs_nc(
7530 const char *filename,
7531 double *rt,
7532 double *rz,
7533 double *rlon,
7534 double *rlat,
7535 double *robs,
7536 int *nobs);
7537
7571double scan_ctl(
7572 const char *filename,
7573 int argc,
7574 char *argv[],
7575 const char *varname,
7576 const int arridx,
7577 const char *defvalue,
7578 char *value);
7579
7606double sedi(
7607 const double p,
7608 const double T,
7609 const double rp,
7610 const double rhop);
7611
7641void spline(
7642 const double *x,
7643 const double *y,
7644 const int n,
7645 const double *x2,
7646 double *y2,
7647 const int n2,
7648 const int method);
7649
7672float stddev(
7673 const float *data,
7674 const int n);
7675
7695double sza_calc(
7696 const double sec,
7697 const double lon,
7698 const double lat);
7699
7724void time2jsec(
7725 const int year,
7726 const int mon,
7727 const int day,
7728 const int hour,
7729 const int min,
7730 const int sec,
7731 const double remain,
7732 double *jsec);
7733
7762void timer(
7763 const char *name,
7764 const char *group,
7765 const int output);
7766
7792double time_from_filename(
7793 const char *filename,
7794 const int offset);
7795
7813double tropo_weight(
7814 const clim_t * clim,
7815 const atm_t * atm,
7816 const int ip);
7817
7840void write_atm_asc(
7841 const char *filename,
7842 const ctl_t * ctl,
7843 const atm_t * atm,
7844 const double t);
7845
7869void write_atm_bin(
7870 const char *filename,
7871 const ctl_t * ctl,
7872 const atm_t * atm);
7873
7897void write_atm_clams(
7898 const char *filename,
7899 const ctl_t * ctl,
7900 const atm_t * atm);
7901
7927 const char *dirname,
7928 const ctl_t * ctl,
7929 const atm_t * atm,
7930 const double t);
7931
7955void write_atm_nc(
7956 const char *filename,
7957 const ctl_t * ctl,
7958 const atm_t * atm);
7959
7988void write_csi(
7989 const char *filename,
7990 const ctl_t * ctl,
7991 const atm_t * atm,
7992 const double t);
7993
8021void write_ens(
8022 const char *filename,
8023 const ctl_t * ctl,
8024 const atm_t * atm,
8025 const double t);
8026
8065void write_grid(
8066 const char *filename,
8067 const ctl_t * ctl,
8068 met_t * met0,
8069 met_t * met1,
8070 const atm_t * atm,
8071 const double t);
8072
8118void write_grid_asc(
8119 const char *filename,
8120 const ctl_t * ctl,
8121 const double *cd,
8122 double *mean[NQ],
8123 double *sigma[NQ],
8124 const double *vmr_impl,
8125 const double t,
8126 const double *z,
8127 const double *lon,
8128 const double *lat,
8129 const double *area,
8130 const double dz,
8131 const int *np);
8132
8175void write_grid_nc(
8176 const char *filename,
8177 const ctl_t * ctl,
8178 const double *cd,
8179 double *mean[NQ],
8180 double *sigma[NQ],
8181 const double *vmr_impl,
8182 const double t,
8183 const double *z,
8184 const double *lon,
8185 const double *lat,
8186 const double *area,
8187 const double dz,
8188 const int *np);
8189
8219void write_met_bin(
8220 const char *filename,
8221 const ctl_t * ctl,
8222 met_t * met);
8223
8251void write_met_bin_2d(
8252 FILE * out,
8253 met_t * met,
8254 float var[EX][EY],
8255 const char *varname);
8256
8294void write_met_bin_3d(
8295 FILE * out,
8296 const ctl_t * ctl,
8297 met_t * met,
8298 float var[EX][EY][EP],
8299 const char *varname,
8300 const int precision,
8301 const double tolerance);
8302
8330void write_met_nc(
8331 const char *filename,
8332 const ctl_t * ctl,
8333 met_t * met);
8334
8359void write_met_nc_2d(
8360 const int ncid,
8361 const char *varname,
8362 met_t * met,
8363 float var[EX][EY],
8364 const float scl);
8365
8391void write_met_nc_3d(
8392 const int ncid,
8393 const char *varname,
8394 met_t * met,
8395 float var[EX][EY][EP],
8396 const float scl);
8397
8428void write_prof(
8429 const char *filename,
8430 const ctl_t * ctl,
8431 met_t * met0,
8432 met_t * met1,
8433 const atm_t * atm,
8434 const double t);
8435
8467void write_sample(
8468 const char *filename,
8469 const ctl_t * ctl,
8470 met_t * met0,
8471 met_t * met1,
8472 const atm_t * atm,
8473 const double t);
8474
8501void write_station(
8502 const char *filename,
8503 const ctl_t * ctl,
8504 atm_t * atm,
8505 const double t);
8506
8535void write_vtk(
8536 const char *filename,
8537 const ctl_t * ctl,
8538 const atm_t * atm,
8539 const double t);
8540
8541#endif /* LIBTRAC_H */
void read_met_geopot(const ctl_t *ctl, met_t *met)
Calculates geopotential heights from meteorological data.
Definition: mptrac.c:6993
#define LEN
Maximum length of ASCII data lines.
Definition: mptrac.h:241
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:5677
void day2doy(const int year, const int mon, const int day, int *doy)
Get day of year from date.
Definition: mptrac.c:900
void read_met_extrapolate(met_t *met)
Extrapolates meteorological data.
Definition: mptrac.c:6953
void read_met_levels(const int ncid, const ctl_t *ctl, met_t *met)
Reads meteorological variables at different vertical levels from a NetCDF file.
Definition: mptrac.c:7253
int read_met_nc(const char *filename, const ctl_t *ctl, const clim_t *clim, met_t *met)
Reads meteorological data from a NetCDF file and processes it.
Definition: mptrac.c:7581
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:9485
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:10845
void mptrac_alloc(ctl_t **ctl, cache_t **cache, clim_t **clim, met_t **met0, met_t **met1, atm_t **atm)
Allocates and initializes memory resources for MPTRAC.
Definition: mptrac.c:4195
void read_met_sample(const ctl_t *ctl, met_t *met)
Downsamples meteorological data based on specified parameters.
Definition: mptrac.c:8348
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:10606
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:8828
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:2084
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:3922
int mptrac_read_met(const char *filename, const ctl_t *ctl, const clim_t *clim, met_t *met)
Reads meteorological data from a file, supporting multiple formats and MPI broadcasting.
Definition: mptrac.c:5373
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:3238
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:6101
void read_met_cloud(met_t *met)
Calculates cloud-related variables for each grid point.
Definition: mptrac.c:6792
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:2624
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:9001
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:1462
void intpol_met_space_3d_ml(const met_t *met, float zs[EX][EY][EP], float vals[EX][EY][EP], const double z, const double lon, const double lat, double *val)
Interpolates meteorological data in 3D space.
Definition: mptrac.c:1395
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:6068
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:7956
void read_met_detrend(const ctl_t *ctl, met_t *met)
Detrends meteorological data.
Definition: mptrac.c:6849
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, float dest[EX][EY], const float scl, const int init)
Reads a 2-dimensional meteorological variable from a NetCDF file.
Definition: mptrac.c:7652
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:8656
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:8872
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:2509
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:2045
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:1095
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:2064
void mptrac_get_met(ctl_t *ctl, clim_t *clim, const double t, met_t **met0, met_t **met1)
Retrieves meteorological data for the specified time.
Definition: mptrac.c:4268
void read_met_monotonize(const ctl_t *ctl, met_t *met)
Makes zeta and pressure profiles monotone.
Definition: mptrac.c:7496
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:6220
void read_met_periodic(met_t *met)
Applies periodic boundary conditions to meteorological data along longitudinal axis.
Definition: mptrac.c:8093
void module_timesteps_init(ctl_t *ctl, const atm_t *atm)
Initialize start time and time interval for time-stepping.
Definition: mptrac.c:3959
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:9954
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 int decompress, FILE *inout)
Compresses or decompresses a 3D array of floats using a custom multiscale compression algorithm.
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:3342
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:401
#define EY
Maximum number of latitudes for meteo data.
Definition: mptrac.h:271
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:6192
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:7454
double clim_tropo(const clim_t *clim, const double t, const double lat)
Calculates the tropopause pressure based on climatological data.
Definition: mptrac.c:200
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:8900
void read_met_grid(const char *filename, const int ncid, const ctl_t *ctl, met_t *met)
Reads meteorological grid information from a NetCDF file.
Definition: mptrac.c:7121
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:6562
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:1981
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:3063
void level_definitions(ctl_t *ctl)
Defines pressure levels for meteorological data.
Definition: mptrac.c:1775
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:10242
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:5565
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:9143
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:4389
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:1520
void fft_help(double *fcReal, double *fcImag, const int n)
Computes the Fast Fourier Transform (FFT) of a complex sequence.
Definition: mptrac.c:949
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:4060
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:665
double nat_temperature(const double p, const double h2o, const double hno3)
Calculates the nitric acid trihydrate (NAT) temperature.
Definition: mptrac.c:5868
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:9034
#define CP
Maximum number of pressure levels for climatological data.
Definition: mptrac.h:301
#define NQ
Maximum number of quantities per data point.
Definition: mptrac.h:251
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:149
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:6274
#define EX
Maximum number of longitudes for meteo data.
Definition: mptrac.h:266
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:3797
void timer(const char *name, const char *group, const int output)
Measures and reports elapsed time for named and grouped timers.
Definition: mptrac.c:9174
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:9300
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:1337
void read_met_surface(const int ncid, const ctl_t *ctl, met_t *met)
Reads surface meteorological data from a netCDF file and stores it in the meteorological data structu...
Definition: mptrac.c:8520
void intpol_met_time_3d_ml(const met_t *met0, float zs0[EX][EY][EP], float array0[EX][EY][EP], const met_t *met1, float zs1[EX][EY][EP], float array1[EX][EY][EP], const double ts, const double p, const double lon, const double lat, double *var)
Interpolates meteorological data in both space and time.
Definition: mptrac.c:1549
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:2551
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:6373
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:2266
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:8929
#define CT
Maximum number of time steps for climatological data.
Definition: mptrac.h:311
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:2239
void module_sort_help(double *a, const int *p, const int np)
Reorder an array based on a given permutation.
Definition: mptrac.c:3886
float stddev(const float *data, const int n)
Calculates the standard deviation of a set of data.
Definition: mptrac.c:9081
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:1607
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, float dest[EX][EY][EP], const float scl)
Reads a 3-dimensional meteorological variable from a NetCDF file.
Definition: mptrac.c:7810
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:6591
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:2011
double time_from_filename(const char *filename, const int offset)
Extracts and converts a timestamp from a filename to Julian seconds.
Definition: mptrac.c:9242
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:10904
void mptrac_read_clim(const ctl_t *ctl, clim_t *clim)
Reads various climatological data and populates the given climatology structure.
Definition: mptrac.c:4479
void module_chem_grid(const ctl_t *ctl, met_t *met0, met_t *met1, atm_t *atm, const double t)
Calculate grid data for chemistry modules.
Definition: mptrac.c:2362
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:9277
void write_met_nc(const char *filename, const ctl_t *ctl, met_t *met)
Writes meteorological data to a NetCDF file.
Definition: mptrac.c:10682
void module_rng_init(const int ntask)
Initialize random number generators for parallel tasks.
Definition: mptrac.c:3661
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:4408
void intpol_met_4d_coord(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:1165
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:5621
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:5777
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
#define CTS
Maximum number of data points of climatological time series.
Definition: mptrac.h:316
void read_met_ozone(met_t *met)
Calculates the total column ozone from meteorological ozone data.
Definition: mptrac.c:8319
void broadcast_large_data(void *data, size_t N)
Broadcasts large data across all processes in an MPI communicator.
void mptrac_free(ctl_t *ctl, cache_t *cache, clim_t *clim, met_t *met0, met_t *met1, atm_t *atm)
Frees memory resources allocated for MPTRAC.
Definition: mptrac.c:4242
void clim_tropo_init(clim_t *clim)
Initializes the tropopause data in the climatology structure.
Definition: mptrac.c:228
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:3692
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:1006
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:11293
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:246
double sza_calc(const double sec, const double lon, const double lat)
Calculates the solar zenith angle.
Definition: mptrac.c:9102
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)
Perform interparcel mixing for a specific quantity.
Definition: mptrac.c:3432
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:919
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:1574
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:3610
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:116
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:10577
void read_met_pv(met_t *met)
Calculates potential vorticity (PV) from meteorological data.
Definition: mptrac.c:8213
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:5956
void get_met_replace(char *orig, char *search, char *repl)
Replaces occurrences of a substring in a string with another substring.
Definition: mptrac.c:1071
#define CY
Maximum number of latitudes for climatological data.
Definition: mptrac.h:291
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:2663
void module_sort(const ctl_t *ctl, met_t *met0, atm_t *atm)
Sort particles according to box index.
Definition: mptrac.c:3826
double clim_ts(const clim_ts_t *ts, const double t)
Interpolates a time series of climatological variables.
Definition: mptrac.c:383
void thrustSortWrapper(double *__restrict__ c, int n, int *__restrict__ index)
Wrapper to Thrust sorting function.
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:1698
int read_met_bin(const char *filename, const ctl_t *ctl, met_t *met)
Reads meteorological data from a binary file.
Definition: mptrac.c:6414
void mptrac_run_timestep(ctl_t *ctl, cache_t *cache, clim_t *clim, met_t **met0, met_t **met1, atm_t *atm, double t)
Executes a single timestep of the MPTRAC model simulation.
Definition: mptrac.c:5428
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:9432
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:2865
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:6012
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:11379
#define EP
Maximum number of pressure levels for meteo data.
Definition: mptrac.h:261
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:3990
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:4539
void read_met_polar_winds(met_t *met)
Applies a fix for polar winds in meteorological data.
Definition: mptrac.c:8154
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:2980
void compress_zstd(const char *varname, float *array, const size_t n, const int decompress, FILE *inout)
Compresses or decompresses an array of floats using the Zstandard (ZSTD) library.
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:10346
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:5892
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:2740
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:10874
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:3133
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:296
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:3526
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:988
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:1731
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:5914
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:1138
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:11131
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:10051
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:2917
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:10476
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:9382
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:6677
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:5737
#define CSZA
Maximum number of solar zenith angles for climatological data.
Definition: mptrac.h:306
double lapse_rate(const double t, const double h2o)
Calculates the moist adiabatic lapse rate in Kelvin per kilometer.
Definition: mptrac.c:1757
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:9692
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:9643
Air parcel data.
Definition: mptrac.h:3149
int np
Number of air parcels.
Definition: mptrac.h:3152
Cache data structure.
Definition: mptrac.h:3177
int iso_n
Isosurface balloon number of data points.
Definition: mptrac.h:3189
Climatological data in the form of photolysis rates.
Definition: mptrac.h:3209
int nsza
Number of solar zenith angles.
Definition: mptrac.h:3215
int np
Number of pressure levels.
Definition: mptrac.h:3212
int no3c
Number of total ozone columns.
Definition: mptrac.h:3218
Climatological data.
Definition: mptrac.h:3317
clim_ts_t ccl2f2
CFC-12 time series.
Definition: mptrac.h:3359
clim_photo_t photo
Photolysis rates.
Definition: mptrac.h:3335
clim_zm_t ho2
HO2 zonal means.
Definition: mptrac.h:3347
clim_zm_t hno3
HNO3 zonal means.
Definition: mptrac.h:3338
int tropo_ntime
Number of tropopause timesteps.
Definition: mptrac.h:3320
clim_ts_t sf6
SF6 time series.
Definition: mptrac.h:3365
clim_ts_t ccl4
CFC-10 time series.
Definition: mptrac.h:3353
clim_ts_t ccl3f
CFC-11 time series.
Definition: mptrac.h:3356
clim_zm_t o1d
O(1D) zonal means.
Definition: mptrac.h:3350
clim_zm_t h2o2
H2O2 zonal means.
Definition: mptrac.h:3344
int tropo_nlat
Number of tropopause latitudes.
Definition: mptrac.h:3323
clim_zm_t oh
OH zonal means.
Definition: mptrac.h:3341
clim_ts_t n2o
N2O time series.
Definition: mptrac.h:3362
Climatological data in the form of time series.
Definition: mptrac.h:3265
int ntime
Number of timesteps.
Definition: mptrac.h:3268
Climatological data in the form of zonal means.
Definition: mptrac.h:3285
int np
Number of pressure levels.
Definition: mptrac.h:3294
int ntime
Number of timesteps.
Definition: mptrac.h:3288
int nlat
Number of latitudes.
Definition: mptrac.h:3291
Control parameters.
Definition: mptrac.h:2158
double grid_z0
Lower altitude of gridded data [km].
Definition: mptrac.h:3029
int qnt_o3
Quantity array index for ozone volume mixing ratio.
Definition: mptrac.h:2270
double csi_lat1
Upper latitude of gridded CSI data [deg].
Definition: mptrac.h:2993
int qnt_Coh
Quantity array index for OH volume mixing ratio (chemistry code).
Definition: mptrac.h:2417
double wet_depo_ic_a
Coefficient A for wet deposition in cloud (exponential form).
Definition: mptrac.h:2881
int met_nc_scale
Check netCDF scaling factors (0=no, 1=yes).
Definition: mptrac.h:2489
int qnt_pel
Quantity array index for pressure at equilibrium level (EL).
Definition: mptrac.h:2303
int csi_nz
Number of altitudes of gridded CSI data.
Definition: mptrac.h:2969
double molmass
Molar mass [g/mol].
Definition: mptrac.h:2743
int qnt_p
Quantity array index for pressure.
Definition: mptrac.h:2249
int qnt_Cccl2f2
Quantity array index for CFC-12 volume mixing ratio (chemistry code).
Definition: mptrac.h:2441
int mixing_nx
Number of longitudes of mixing grid.
Definition: mptrac.h:2806
double chemgrid_z1
Upper altitude of chemistry grid [km].
Definition: mptrac.h:2830
int qnt_m
Quantity array index for mass.
Definition: mptrac.h:2189
int qnt_aoa
Quantity array index for age of air.
Definition: mptrac.h:2450
int qnt_rhop
Quantity array index for particle density.
Definition: mptrac.h:2198
int qnt_swc
Quantity array index for cloud snow water content.
Definition: mptrac.h:2282
double csi_obsmin
Minimum observation index to trigger detection.
Definition: mptrac.h:2963
int qnt_pcb
Quantity array index for cloud bottom pressure.
Definition: mptrac.h:2291
double bound_dzs
Boundary conditions surface layer depth [km].
Definition: mptrac.h:2731
double csi_lon1
Upper longitude of gridded CSI data [deg].
Definition: mptrac.h:2984
int qnt_u
Quantity array index for zonal wind.
Definition: mptrac.h:2258
double stat_lon
Longitude of station [deg].
Definition: mptrac.h:3107
double mixing_trop
Interparcel exchange parameter for mixing in the troposphere.
Definition: mptrac.h:2791
double sort_dt
Time step for sorting of particle data [s].
Definition: mptrac.h:2642
double mixing_z1
Upper altitude of mixing grid [km].
Definition: mptrac.h:2803
double stat_r
Search radius around station [km].
Definition: mptrac.h:3113
double wet_depo_bc_a
Coefficient A for wet deposition below cloud (exponential form).
Definition: mptrac.h:2875
int csi_ny
Number of latitudes of gridded CSI data.
Definition: mptrac.h:2987
int vtk_sphere
Spherical projection for VTK data (0=no, 1=yes).
Definition: mptrac.h:3137
double chemgrid_z0
Lower altitude of chemistry grid [km].
Definition: mptrac.h:2827
double met_pbl_min
Minimum depth of planetary boundary layer [km].
Definition: mptrac.h:2610
int qnt_iwc
Quantity array index for cloud ice water content.
Definition: mptrac.h:2279
double chemgrid_lat0
Lower latitude of chemistry grid [deg].
Definition: mptrac.h:2845
double conv_cape
CAPE threshold for convection module [J/kg].
Definition: mptrac.h:2695
int qnt_Co1d
Quantity array index for O(1D) volume mixing ratio (chemistry code).
Definition: mptrac.h:2429
double met_cms_eps_pv
cmultiscale compression epsilon for potential vorticity.
Definition: mptrac.h:2532
int qnt_pw
Quantity array index for partial water vapor pressure.
Definition: mptrac.h:2357
double grid_z1
Upper altitude of gridded data [km].
Definition: mptrac.h:3032
int direction
Direction flag (1=forward calculation, -1=backward calculation).
Definition: mptrac.h:2453
int qnt_Cccl4
Quantity array index for CFC-10 volume mixing ratio (chemistry code).
Definition: mptrac.h:2435
int met_dp
Stride for pressure levels.
Definition: mptrac.h:2562
double met_dt_out
Time step for sampling of meteo data along trajectories [s].
Definition: mptrac.h:2629
int qnt_h2o2
Quantity array index for H2O2 volume mixing ratio (climatology).
Definition: mptrac.h:2321
int qnt_vh
Quantity array index for horizontal wind.
Definition: mptrac.h:2384
int csi_nx
Number of longitudes of gridded CSI data.
Definition: mptrac.h:2978
double csi_lat0
Lower latitude of gridded CSI data [deg].
Definition: mptrac.h:2990
double turb_dz_trop
Vertical turbulent diffusion coefficient (troposphere) [m^2/s].
Definition: mptrac.h:2677
int met_pbl
Planetary boundary layer data (0=file, 1=z2p, 2=Richardson, 3=theta).
Definition: mptrac.h:2607
int qnt_lwc
Quantity array index for cloud liquid water content.
Definition: mptrac.h:2273
double turb_mesoz
Vertical scaling factor for mesoscale wind fluctuations.
Definition: mptrac.h:2686
int grid_nc_level
zlib compression level of netCDF grid data files (0=off).
Definition: mptrac.h:3017
int grid_nx
Number of longitudes of gridded data.
Definition: mptrac.h:3035
int atm_type
Type of atmospheric data files (0=ASCII, 1=binary, 2=netCDF, 3=CLaMS_traj, 4=CLaMS_pos).
Definition: mptrac.h:2934
double bound_mass
Boundary conditions mass per particle [kg].
Definition: mptrac.h:2704
double grid_lat0
Lower latitude of gridded data [deg].
Definition: mptrac.h:3047
int qnt_ts
Quantity array index for surface temperature.
Definition: mptrac.h:2204
int qnt_loss_rate
Quantity array index for total loss rate.
Definition: mptrac.h:2348
double met_cms_eps_h2o
cmultiscale compression epsilon for water vapor.
Definition: mptrac.h:2535
int qnt_plfc
Quantity array index for pressure at level of free convection (LCF).
Definition: mptrac.h:2300
double grid_lon0
Lower longitude of gridded data [deg].
Definition: mptrac.h:3038
int qnt_o1d
Quantity array index for O(1D) volume mixing ratio (climatology).
Definition: mptrac.h:2327
int met_tropo_spline
Tropopause interpolation method (0=linear, 1=spline).
Definition: mptrac.h:2626
int qnt_tvirt
Quantity array index for virtual temperature.
Definition: mptrac.h:2378
double dt_met
Time step of meteo data [s].
Definition: mptrac.h:2472
double chemgrid_lat1
Upper latitude of chemistry grid [deg].
Definition: mptrac.h:2848
int met_geopot_sy
Latitudinal smoothing of geopotential heights.
Definition: mptrac.h:2598
double met_cms_eps_u
cmultiscale compression epsilon for zonal wind.
Definition: mptrac.h:2523
double turb_dx_strat
Horizontal turbulent diffusion coefficient (stratosphere) [m^2/s].
Definition: mptrac.h:2671
int qnt_vmr
Quantity array index for volume mixing ratio.
Definition: mptrac.h:2192
int qnt_lsm
Quantity array index for land-sea mask.
Definition: mptrac.h:2225
int qnt_theta
Quantity array index for potential temperature.
Definition: mptrac.h:2369
double bound_lat1
Boundary conditions maximum longitude [deg].
Definition: mptrac.h:2719
double stat_t1
Stop time for station output [s].
Definition: mptrac.h:3119
double turb_dx_trop
Horizontal turbulent diffusion coefficient (troposphere) [m^2/s].
Definition: mptrac.h:2668
int grid_type
Type of grid data files (0=ASCII, 1=netCDF).
Definition: mptrac.h:3053
double csi_lon0
Lower longitude of gridded CSI data [deg].
Definition: mptrac.h:2981
int qnt_pbl
Quantity array index for boundary layer pressure.
Definition: mptrac.h:2231
int grid_stddev
Include standard deviations in grid output (0=no, 1=yes).
Definition: mptrac.h:3023
int qnt_psice
Quantity array index for saturation pressure over ice.
Definition: mptrac.h:2354
double chemgrid_lon0
Lower longitude of chemistry grid [deg].
Definition: mptrac.h:2836
int bound_pbl
Boundary conditions planetary boundary layer (0=no, 1=yes).
Definition: mptrac.h:2737
int qnt_mloss_wet
Quantity array index for total mass loss due to wet deposition.
Definition: mptrac.h:2339
int met_geopot_sx
Longitudinal smoothing of geopotential heights.
Definition: mptrac.h:2595
int met_sy
Smoothing for latitudes.
Definition: mptrac.h:2568
int qnt_ps
Quantity array index for surface pressure.
Definition: mptrac.h:2201
int rng_type
Random number generator (0=GSL, 1=Squares, 2=cuRAND).
Definition: mptrac.h:2659
int isosurf
Isosurface parameter (0=none, 1=pressure, 2=density, 3=theta, 4=balloon).
Definition: mptrac.h:2646
double bound_p1
Boundary conditions top pressure [hPa].
Definition: mptrac.h:2725
int qnt_zs
Quantity array index for surface geopotential height.
Definition: mptrac.h:2207
int prof_nz
Number of altitudes of gridded profile data.
Definition: mptrac.h:3062
double csi_dt_out
Time step for CSI output [s].
Definition: mptrac.h:2957
int met_cape
Convective available potential energy data (0=file, 1=calculate).
Definition: mptrac.h:2604
double csi_modmin
Minimum column density to trigger detection [kg/m^2].
Definition: mptrac.h:2966
int met_sx
Smoothing for longitudes.
Definition: mptrac.h:2565
double chemgrid_lon1
Upper longitude of chemistry grid [deg].
Definition: mptrac.h:2839
double turb_mesox
Horizontal scaling factor for mesoscale wind fluctuations.
Definition: mptrac.h:2683
double met_cms_eps_iwc
cmultiscale compression epsilon for cloud ice water content.
Definition: mptrac.h:2547
double met_cms_eps_swc
cmultiscale compression epsilon for cloud snow water content.
Definition: mptrac.h:2550
int met_zfp_prec
ZFP compression precision for all variables, except z and T.
Definition: mptrac.h:2498
double met_cms_eps_v
cmultiscale compression epsilon for meridional wind.
Definition: mptrac.h:2526
double prof_z0
Lower altitude of gridded profile data [km].
Definition: mptrac.h:3065
int qnt_w
Quantity array index for vertical velocity.
Definition: mptrac.h:2264
double bound_vmr
Boundary conditions volume mixing ratio [ppv].
Definition: mptrac.h:2710
double met_tropo_pv
Dynamical tropopause potential vorticity threshold [PVU].
Definition: mptrac.h:2620
int prof_nx
Number of longitudes of gridded profile data.
Definition: mptrac.h:3071
int qnt_stat
Quantity array index for station flag.
Definition: mptrac.h:2186
int met_tropo
Tropopause definition (0=none, 1=clim, 2=cold point, 3=WMO_1st, 4=WMO_2nd, 5=dynamical).
Definition: mptrac.h:2617
int qnt_rp
Quantity array index for particle radius.
Definition: mptrac.h:2195
int met_mpi_share
Use MPI to share meteo (0=no, 1=yes).
Definition: mptrac.h:2635
double mixing_strat
Interparcel exchange parameter for mixing in the stratosphere.
Definition: mptrac.h:2794
int qnt_vz
Quantity array index for vertical velocity.
Definition: mptrac.h:2387
int qnt_ho2
Quantity array index for HO2 volume mixing ratio (climatology).
Definition: mptrac.h:2324
double csi_z1
Upper altitude of gridded CSI data [km].
Definition: mptrac.h:2975
double stat_t0
Start time for station output [s].
Definition: mptrac.h:3116
double oh_chem_beta
Beta parameter for diurnal variablity of OH.
Definition: mptrac.h:2857
double wet_depo_so2_ph
pH value used to calculate effective Henry constant of SO2.
Definition: mptrac.h:2893
double mixing_z0
Lower altitude of mixing grid [km].
Definition: mptrac.h:2800
int qnt_mloss_decay
Quantity array index for total mass loss due to exponential decay.
Definition: mptrac.h:2345
int atm_type_out
Type of atmospheric data files for output (-1=same as ATM_TYPE, 0=ASCII, 1=binary,...
Definition: mptrac.h:2939
int met_nlev
Number of meteo data model levels.
Definition: mptrac.h:2586
double dt_kpp
Time step for KPP chemistry [s].
Definition: mptrac.h:2866
double dry_depo_dp
Dry deposition surface layer [hPa].
Definition: mptrac.h:2902
int qnt_shf
Quantity array index for surface sensible heat flux.
Definition: mptrac.h:2222
int qnt_vs
Quantity array index for surface meridional wind.
Definition: mptrac.h:2213
int qnt_Cco
Quantity array index for CO volume mixing ratio (chemistry code).
Definition: mptrac.h:2414
double vtk_dt_out
Time step for VTK data output [s].
Definition: mptrac.h:3125
double t_stop
Stop time of simulation [s].
Definition: mptrac.h:2459
double conv_dt
Time interval for convection module [s].
Definition: mptrac.h:2701
int qnt_hno3
Quantity array index for HNO3 volume mixing ratio (climatology).
Definition: mptrac.h:2315
int met_clams
Read MPTRAC or CLaMS meteo data (0=MPTRAC, 1=CLaMS).
Definition: mptrac.h:2486
int qnt_h2ot
Quantity array index for tropopause water vapor volume mixing ratio.
Definition: mptrac.h:2243
int qnt_rh
Quantity array index for relative humidity over water.
Definition: mptrac.h:2363
double met_cms_eps_cc
cmultiscale compression epsilon for cloud cover.
Definition: mptrac.h:2553
double bound_lat0
Boundary conditions minimum longitude [deg].
Definition: mptrac.h:2716
double met_pbl_max
Maximum depth of planetary boundary layer [km].
Definition: mptrac.h:2613
int met_dx
Stride for longitudes.
Definition: mptrac.h:2556
int mixing_ny
Number of latitudes of mixing grid.
Definition: mptrac.h:2815
int met_convention
Meteo data layout (0=[lev, lat, lon], 1=[lon, lat, lev]).
Definition: mptrac.h:2475
int qnt_zeta_d
Quantity array index for diagnosed zeta vertical coordinate.
Definition: mptrac.h:2375
int tracer_chem
Switch for first order tracer chemistry module (0=off, 1=on).
Definition: mptrac.h:2869
double dt_mod
Time step of simulation [s].
Definition: mptrac.h:2462
int diffusion
Diffusion scheme (0=off, 1=fixed-K, 2=PBL).
Definition: mptrac.h:2662
int qnt_tnat
Quantity array index for T_NAT.
Definition: mptrac.h:2402
int qnt_tice
Quantity array index for T_ice.
Definition: mptrac.h:2396
int qnt_zg
Quantity array index for geopotential height.
Definition: mptrac.h:2246
double vtk_offset
Vertical offset for VTK data [km].
Definition: mptrac.h:3134
int qnt_v
Quantity array index for meridional wind.
Definition: mptrac.h:2261
int qnt_mloss_dry
Quantity array index for total mass loss due to dry deposition.
Definition: mptrac.h:2342
double bound_vmr_trend
Boundary conditions volume mixing ratio trend [ppv/s].
Definition: mptrac.h:2713
int met_cache
Preload meteo data into disk cache (0=no, 1=yes).
Definition: mptrac.h:2632
int qnt_oh
Quantity array index for OH volume mixing ratio (climatology).
Definition: mptrac.h:2318
int qnt_Ch
Quantity array index for H volume mixing ratio (chemistry code).
Definition: mptrac.h:2420
int met_press_level_def
Use predefined pressure levels or not.
Definition: mptrac.h:2583
int oh_chem_reaction
Reaction type for OH chemistry (0=none, 2=bimolecular, 3=termolecular).
Definition: mptrac.h:2851
int qnt_h2o
Quantity array index for water vapor volume mixing ratio.
Definition: mptrac.h:2267
int prof_ny
Number of latitudes of gridded profile data.
Definition: mptrac.h:3080
int qnt_rhice
Quantity array index for relative humidity over ice.
Definition: mptrac.h:2366
int qnt_rho
Quantity array index for density of air.
Definition: mptrac.h:2255
double sample_dz
Layer depth for sample output [km].
Definition: mptrac.h:3101
double tdec_strat
Life time of particles in the stratosphere [s].
Definition: mptrac.h:2749
int obs_type
Type of observation data files (0=ASCII, 1=netCDF).
Definition: mptrac.h:2948
double met_cms_eps_lwc
cmultiscale compression epsilon for cloud liquid water content.
Definition: mptrac.h:2541
int qnt_us
Quantity array index for surface zonal wind.
Definition: mptrac.h:2210
double met_cms_eps_z
cmultiscale compression epsilon for geopotential height.
Definition: mptrac.h:2517
double grid_lon1
Upper longitude of gridded data [deg].
Definition: mptrac.h:3041
int qnt_Cn2o
Quantity array index for N2O volume mixing ratio (chemistry code).
Definition: mptrac.h:2444
int qnt_Cccl3f
Quantity array index for CFC-11 volume mixing ratio (chemistry code).
Definition: mptrac.h:2438
double mixing_lat0
Lower latitude of mixing grid [deg].
Definition: mptrac.h:2818
int qnt_pt
Quantity array index for tropopause pressure.
Definition: mptrac.h:2234
int qnt_cl
Quantity array index for total column cloud water.
Definition: mptrac.h:2294
int advect
Advection scheme (0=off, 1=Euler, 2=midpoint, 4=Runge-Kutta).
Definition: mptrac.h:2652
double prof_z1
Upper altitude of gridded profile data [km].
Definition: mptrac.h:3068
int qnt_t
Quantity array index for temperature.
Definition: mptrac.h:2252
int atm_filter
Time filter for atmospheric data output (0=none, 1=missval, 2=remove).
Definition: mptrac.h:2927
int kpp_chem
Switch for KPP chemistry module (0=off, 1=on).
Definition: mptrac.h:2863
int qnt_zeta
Quantity array index for zeta vertical coordinate.
Definition: mptrac.h:2372
double conv_pbl_trans
Depth of PBL transition layer (fraction of PBL depth).
Definition: mptrac.h:2692
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:2479
double csi_z0
Lower altitude of gridded CSI data [km].
Definition: mptrac.h:2972
int qnt_lapse
Quantity array index for lapse rate.
Definition: mptrac.h:2381
double stat_lat
Latitude of station [deg].
Definition: mptrac.h:3110
int qnt_Cho2
Quantity array index for HO2 volume mixing ratio (chemistry code).
Definition: mptrac.h:2423
int grid_ny
Number of latitudes of gridded data.
Definition: mptrac.h:3044
int qnt_Csf6
Quantity array index for SF6 volume mixing ratio (chemistry code).
Definition: mptrac.h:2447
int qnt_Ch2o
Quantity array index for H2O volume mixing ratio (chemistry code).
Definition: mptrac.h:2408
double met_detrend
FWHM of horizontal Gaussian used for detrending [km].
Definition: mptrac.h:2574
int conv_mix_pbl
Vertical mixing in the PBL (0=off, 1=on).
Definition: mptrac.h:2689
double bound_dps
Boundary conditions surface layer depth [hPa].
Definition: mptrac.h:2728
double met_cms_eps_t
cmultiscale compression epsilon for temperature.
Definition: mptrac.h:2520
int chemgrid_nz
Number of altitudes of chemistry grid.
Definition: mptrac.h:2824
int qnt_cape
Quantity array index for convective available potential energy (CAPE).
Definition: mptrac.h:2306
double bound_mass_trend
Boundary conditions mass per particle trend [kg/s].
Definition: mptrac.h:2707
int mixing_nz
Number of altitudes of mixing grid.
Definition: mptrac.h:2797
int qnt_o3c
Quantity array index for total column ozone.
Definition: mptrac.h:2312
double bound_p0
Boundary conditions bottom pressure [hPa].
Definition: mptrac.h:2722
double mixing_lon0
Lower longitude of mixing grid [deg].
Definition: mptrac.h:2809
int qnt_Co3
Quantity array index for O3 volume mixing ratio (chemistry code).
Definition: mptrac.h:2411
int qnt_tsts
Quantity array index for T_STS.
Definition: mptrac.h:2399
int grid_nz
Number of altitudes of gridded data.
Definition: mptrac.h:3026
int qnt_nss
Quantity array index for northward turbulent surface stress.
Definition: mptrac.h:2219
double ens_dt_out
Time step for ensemble output [s].
Definition: mptrac.h:2999
int atm_stride
Particle index stride for atmospheric data files.
Definition: mptrac.h:2930
int met_relhum
Try to read relative humidity (0=no, 1=yes).
Definition: mptrac.h:2601
double mixing_lat1
Upper latitude of mixing grid [deg].
Definition: mptrac.h:2821
double atm_dt_out
Time step for atmospheric data output [s].
Definition: mptrac.h:2924
double prof_lat1
Upper latitude of gridded profile data [deg].
Definition: mptrac.h:3086
int met_cms_batch
cmultiscale batch size.
Definition: mptrac.h:2507
double psc_h2o
H2O volume mixing ratio for PSC analysis.
Definition: mptrac.h:2908
int met_sp
Smoothing for pressure levels.
Definition: mptrac.h:2571
double prof_lon0
Lower longitude of gridded profile data [deg].
Definition: mptrac.h:3074
int chemgrid_nx
Number of longitudes of chemistry grid.
Definition: mptrac.h:2833
int qnt_pct
Quantity array index for cloud top pressure.
Definition: mptrac.h:2288
int qnt_mloss_kpp
Quantity array index for total mass loss due to KPP chemistry.
Definition: mptrac.h:2336
int qnt_psat
Quantity array index for saturation pressure over water.
Definition: mptrac.h:2351
double prof_lat0
Lower latitude of gridded profile data [deg].
Definition: mptrac.h:3083
int qnt_cin
Quantity array index for convective inhibition (CIN).
Definition: mptrac.h:2309
double psc_hno3
HNO3 volume mixing ratio for PSC analysis.
Definition: mptrac.h:2911
double prof_lon1
Upper longitude of gridded profile data [deg].
Definition: mptrac.h:3077
double met_cms_eps_rwc
cmultiscale compression epsilon for cloud rain water content.
Definition: mptrac.h:2544
int met_nc_quant
Number of digits for quantization of netCDF meteo files (0=off).
Definition: mptrac.h:2495
int h2o2_chem_reaction
Reaction type for H2O2 chemistry (0=none, 1=SO2).
Definition: mptrac.h:2860
int qnt_Co3p
Quantity array index for O(3P) volume mixing ratio (chemistry code).
Definition: mptrac.h:2432
double wet_depo_bc_ret_ratio
Coefficients for wet deposition below cloud: retention ratio.
Definition: mptrac.h:2899
int chemgrid_ny
Number of latitudes of chemistry grid.
Definition: mptrac.h:2842
double met_cms_eps_o3
cmultiscale compression epsilon for ozone.
Definition: mptrac.h:2538
int met_cms_zstd
cmultiscale zstd compression (0=off, 1=on).
Definition: mptrac.h:2510
int grid_sparse
Sparse output in grid data files (0=no, 1=yes).
Definition: mptrac.h:3014
double dry_depo_vdep
Dry deposition velocity [m/s].
Definition: mptrac.h:2905
int qnt_tt
Quantity array index for tropopause temperature.
Definition: mptrac.h:2237
int met_np
Number of target pressure levels.
Definition: mptrac.h:2577
int qnt_ens
Quantity array index for ensemble IDs.
Definition: mptrac.h:2183
int met_nc_level
zlib compression level of netCDF meteo files (0=off).
Definition: mptrac.h:2492
double met_zfp_tol_t
ZFP compression tolerance for temperature.
Definition: mptrac.h:2501
double mixing_dt
Time interval for mixing [s].
Definition: mptrac.h:2788
int qnt_mloss_h2o2
Quantity array index for total mass loss due to H2O2 chemistry.
Definition: mptrac.h:2333
double met_zfp_tol_z
ZFP compression tolerance for geopotential height.
Definition: mptrac.h:2504
double vtk_scale
Vertical scaling factor for VTK data.
Definition: mptrac.h:3131
double met_cms_eps_w
cmultiscale compression epsilon for vertical velocity.
Definition: mptrac.h:2529
double turb_dx_pbl
Horizontal turbulent diffusion coefficient (PBL) [m^2/s].
Definition: mptrac.h:2665
double conv_cin
CIN threshold for convection module [J/kg].
Definition: mptrac.h:2698
int qnt_pv
Quantity array index for potential vorticity.
Definition: mptrac.h:2390
int advect_vert_coord
Vertical velocity of air parcels (0=omega_on_plev, 1=zetadot_on_mlev, 2=omega_on_mlev).
Definition: mptrac.h:2656
int qnt_mloss_oh
Quantity array index for total mass loss due to OH chemistry.
Definition: mptrac.h:2330
int qnt_Ch2o2
Quantity array index for H2O2 volume mixing ratio (chemistry code).
Definition: mptrac.h:2426
int qnt_sst
Quantity array index for sea surface temperature.
Definition: mptrac.h:2228
double mixing_lon1
Upper longitude of mixing grid [deg].
Definition: mptrac.h:2812
int atm_nc_level
zlib compression level of netCDF atmospheric data files (0=off).
Definition: mptrac.h:2942
int met_cms_heur
cmultiscale coarsening heuristics (0=default, 1=mean diff, 2=median diff, 3=max diff).
Definition: mptrac.h:2514
double wet_depo_ic_ret_ratio
Coefficients for wet deposition in cloud: retention ratio.
Definition: mptrac.h:2896
int qnt_sh
Quantity array index for specific humidity.
Definition: mptrac.h:2360
int qnt_ess
Quantity array index for eastward turbulent surface stress.
Definition: mptrac.h:2216
double wet_depo_ic_b
Coefficient B for wet deposition in cloud (exponential form).
Definition: mptrac.h:2884
double wet_depo_bc_b
Coefficient B for wet deposition below cloud (exponential form).
Definition: mptrac.h:2878
int met_dy
Stride for latitudes.
Definition: mptrac.h:2559
int qnt_Cx
Quantity array index for trace species x volume mixing ratio (chemistry code).
Definition: mptrac.h:2405
double turb_dz_strat
Vertical turbulent diffusion coefficient (stratosphere) [m^2/s].
Definition: mptrac.h:2680
double bound_zetas
Boundary conditions surface layer zeta [K].
Definition: mptrac.h:2734
int qnt_idx
Quantity array index for air parcel IDs.
Definition: mptrac.h:2180
double met_tropo_theta
Dynamical tropopause potential temperature threshold [K].
Definition: mptrac.h:2623
int qnt_rwc
Quantity array index for cloud rain water content.
Definition: mptrac.h:2276
double t_start
Start time of simulation [s].
Definition: mptrac.h:2456
int nq
Number of quantities.
Definition: mptrac.h:2165
double tdec_trop
Life time of particles in the troposphere [s].
Definition: mptrac.h:2746
double sample_dx
Horizontal radius for sample output [km].
Definition: mptrac.h:3098
int vtk_stride
Particle index stride for VTK data.
Definition: mptrac.h:3128
double turb_dz_pbl
Vertical turbulent diffusion coefficient (PBL) [m^2/s].
Definition: mptrac.h:2674
double grid_lat1
Upper latitude of gridded data [deg].
Definition: mptrac.h:3050
int qnt_zt
Quantity array index for tropopause geopotential height.
Definition: mptrac.h:2240
int met_type
Type of meteo data files (0=netCDF, 1=binary, 2=pck, 3=zfp, 4=zstd, 5=cms).
Definition: mptrac.h:2483
int qnt_cc
Quantity array index for cloud cover.
Definition: mptrac.h:2285
int qnt_plcl
Quantity array index for pressure at lifted condensation level (LCL).
Definition: mptrac.h:2297
double grid_dt_out
Time step for gridded data output [s].
Definition: mptrac.h:3011
int qnt_tdew
Quantity array index for dew point temperature.
Definition: mptrac.h:2393
Meteo data structure.
Definition: mptrac.h:3376
int nx
Number of longitudes.
Definition: mptrac.h:3382
int ny
Number of latitudes.
Definition: mptrac.h:3385
int np
Number of pressure levels.
Definition: mptrac.h:3388
int npl
Number of model levels.
Definition: mptrac.h:3391
double time
Time [s].
Definition: mptrac.h:3379