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
115#ifndef LIBTRAC_H
116#define LIBTRAC_H
117
118/* ------------------------------------------------------------
119 Includes...
120 ------------------------------------------------------------ */
121
122#include <ctype.h>
123#include <gsl/gsl_fft_complex.h>
124#include <gsl/gsl_math.h>
125#include <gsl/gsl_randist.h>
126#include <gsl/gsl_rng.h>
127#include <gsl/gsl_sort.h>
128#include <gsl/gsl_spline.h>
129#include <gsl/gsl_statistics.h>
130#include <math.h>
131#include <netcdf.h>
132#include <omp.h>
133#include <stdint.h>
134#include <stdio.h>
135#include <stdlib.h>
136#include <string.h>
137#include <time.h>
138#include <sys/time.h>
139
140#ifdef MPI
141#include "mpi.h"
142#endif
143
144#ifdef DD
145#include <netcdf_par.h>
146#include <stdbool.h>
147#endif
148
149#ifdef _OPENACC
150#include "openacc.h"
151#endif
152
153#ifdef CURAND
154#include "curand.h"
155#endif
156
157#ifdef THRUST
158#include "thrustsort.h"
159#endif
160
161#ifdef ZFP
162#include "zfp.h"
163#endif
164
165#ifdef ZSTD
166#include "zstd.h"
167#endif
168
169#ifdef SZ3
170#include "SZ3c/sz3c.h"
171#endif
172
173#ifdef CMS
174#include "cmultiscale.h"
175#endif
176
177#ifdef KPP
178#include "chem_Parameters.h"
179#include "chem_Global.h"
180#include "chem_Sparse.h"
181#endif
182
183#ifdef ECCODES
184#include "eccodes.h"
185#endif
186
187/* ------------------------------------------------------------
188 Constants...
189 ------------------------------------------------------------ */
190
192#ifndef AVO
193#define AVO 6.02214076e23
194#endif
195
197#ifndef CPD
198#define CPD 1003.5
199#endif
200
202#ifndef EPS
203#define EPS (MH2O / MA)
204#endif
205
207#ifndef G0
208#define G0 9.80665
209#endif
210
212#ifndef H0
213#define H0 7.0
214#endif
215
217#ifndef LV
218#define LV 2501000.
219#endif
220
222#ifndef KARMAN
223#define KARMAN 0.40
224#endif
225
227#ifndef KB
228#define KB 1.3806504e-23
229#endif
230
232#ifndef MA
233#define MA 28.9644
234#endif
235
237#ifndef MH2O
238#define MH2O 18.01528
239#endif
240
242#ifndef MO3
243#define MO3 48.00
244#endif
245
247#ifndef P0
248#define P0 1013.25
249#endif
250
252#ifndef RA
253#define RA (1e3 * RI / MA)
254#endif
255
257#ifndef RE
258#define RE 6367.421
259#endif
260
262#ifndef RI
263#define RI 8.3144598
264#endif
265
267#ifndef T0
268#define T0 273.15
269#endif
270
272#ifndef DD_NPOLE
273#define DD_NPOLE -2
274#endif
275
276
278#ifndef DD_SPOLE
279#define DD_SPOLE -3
280#endif
281
282/* ------------------------------------------------------------
283 Dimensions...
284 ------------------------------------------------------------ */
285
287#ifndef EP
288#define EP 140
289#endif
290
292#ifndef EX
293#define EX 1444
294#endif
295
297#ifndef EY
298#define EY 724
299#endif
300
302#ifndef EP_GLOB
303#define EP_GLOB 150
304#endif
305
307#ifndef EX_GLOB
308#define EX_GLOB 1202
309#endif
310
312#ifndef EY_GLOB
313#define EY_GLOB 602
314#endif
315
317#ifndef LEN
318#define LEN 5000
319#endif
320
322#ifndef METVAR
323#define METVAR 13
324#endif
325
327#ifndef NP
328#define NP 10000000
329#endif
330
332#ifndef NQ
333#define NQ 15
334#endif
335
337#ifndef NCSI
338#define NCSI 1000000
339#endif
340
342#ifndef NENS
343#define NENS 2000
344#endif
345
347#ifndef NOBS
348#define NOBS 10000000
349#endif
350
352#ifndef NTHREADS
353#define NTHREADS 512
354#endif
355
357#ifndef CY
358#define CY 250
359#endif
360
362#ifndef CO3
363#define CO3 30
364#endif
365
367#ifndef CP
368#define CP 70
369#endif
370
372#ifndef CSZA
373#define CSZA 50
374#endif
375
377#ifndef CT
378#define CT 12
379#endif
380
382#ifndef CTS
383#define CTS 1000
384#endif
385
387#ifndef DD_NPART
388#define DD_NPART 100000
389#endif
390
392#ifndef DD_NNMAX
393#define DD_NNMAX 26
394#endif
395
396/* ------------------------------------------------------------
397 Macros...
398 ------------------------------------------------------------ */
399
419#ifdef _OPENACC
420#define ALLOC(ptr, type, n) \
421 if(acc_get_num_devices(acc_device_nvidia) <= 0) \
422 ERRMSG("Not running on a GPU device!"); \
423 if((ptr=calloc((size_t)(n), sizeof(type)))==NULL) \
424 ERRMSG("Out of memory!");
425#else
426#define ALLOC(ptr, type, n) \
427 if((ptr=calloc((size_t)(n), sizeof(type)))==NULL) \
428 ERRMSG("Out of memory!");
429#endif
430
449#define ARRAY_2D(ix, iy, ny) \
450 ((ix) * (ny) + (iy))
451
468#define ARRAY_3D(ix, iy, ny, iz, nz) \
469 (((ix)*(ny) + (iy)) * (nz) + (iz))
470
493#define ARRHENIUS(a, b, t) \
494 ((a) * exp( -(b) / (t)))
495
517#define DEG2DX(dlon, lat) \
518 (RE * DEG2RAD(dlon) * cos(DEG2RAD(lat)))
519
538#define DEG2DY(dlat) \
539 (RE * DEG2RAD(dlat))
540
555#define DEG2RAD(deg) \
556 ((deg) * (M_PI / 180.0))
557
580#define DP2DZ(dp, p) \
581 (- (dp) * H0 / (p))
582
602#define DX2DEG(dx, lat) \
603 (((lat) < -89.999 || (lat) > 89.999) ? 0 \
604 : (dx) * 180. / (M_PI * RE * cos(DEG2RAD(lat))))
605
620#define DY2DEG(dy) \
621 ((dy) * 180. / (M_PI * RE))
622
639#define DZ2DP(dz, p) \
640 (-(dz) * (p) / H0)
641
655#define DIST(a, b) \
656 sqrt(DIST2(a, b))
657
671#define DIST2(a, b) \
672 ((a[0]-b[0])*(a[0]-b[0])+(a[1]-b[1])*(a[1]-b[1])+(a[2]-b[2])*(a[2]-b[2]))
673
687#define DOTP(a, b) \
688 (a[0]*b[0]+a[1]*b[1]+a[2]*b[2])
689
701#define ECC(cmd) { \
702 int ecc_result=(cmd); \
703 if(ecc_result!=0) \
704 ERRMSG("ECCODES error: %s", codes_get_error_message(ecc_result)); \
705 }
706
720#define ECC_READ_2D(variable, target, scaling_factor, found_flag) { \
721 if(strcmp(short_name, variable) == 0) { \
722 for (int ix = 0; ix < met->nx; ix++) \
723 for (int iy = 0; iy < met->ny; iy++) \
724 target[ix][iy] = (float)(values[iy * met->nx + ix] * scaling_factor); \
725 found_flag = 1; \
726 } \
727 }
728
743#define ECC_READ_3D(variable, level, target, scaling_factor, found_flag) { \
744 if(strcmp(short_name, variable) == 0) { \
745 for (int ix = 0; ix < met->nx; ix++) \
746 for (int iy = 0; iy < met->ny; iy++) \
747 target[ix][iy][level] = (float) (values[iy * met->nx + ix] * scaling_factor); \
748 found_flag += 1; \
749 } \
750 }
751
768#define FMOD(x, y) \
769 ((x) - (int) ((x) / (y)) * (y))
770
786#define FREAD(ptr, type, size, in) { \
787 if(fread(ptr, sizeof(type), size, in)!=size) \
788 ERRMSG("Error while reading!"); \
789 }
790
806#define FWRITE(ptr, type, size, out) { \
807 if(fwrite(ptr, sizeof(type), size, out)!=size) \
808 ERRMSG("Error while writing!"); \
809 }
810
821#define INTPOL_INIT \
822 double cw[4] = {0.0, 0.0, 0.0, 0.0}; int ci[3] = {0, 0, 0};
823
835#define INTPOL_2D(var, init) \
836 intpol_met_time_2d(met0, met0->var, met1, met1->var, \
837 atm->time[ip], atm->lon[ip], atm->lat[ip], \
838 &var, ci, cw, init);
839
852#define INTPOL_3D(var, init) \
853 intpol_met_time_3d(met0, met0->var, met1, met1->var, \
854 atm->time[ip], atm->p[ip], \
855 atm->lon[ip], atm->lat[ip], \
856 &var, ci, cw, init);
857
871#define INTPOL_SPACE_ALL(p, lon, lat) { \
872 intpol_met_space_3d(met, met->z, p, lon, lat, &z, ci, cw, 1); \
873 intpol_met_space_3d(met, met->t, p, lon, lat, &t, ci, cw, 0); \
874 intpol_met_space_3d(met, met->u, p, lon, lat, &u, ci, cw, 0); \
875 intpol_met_space_3d(met, met->v, p, lon, lat, &v, ci, cw, 0); \
876 intpol_met_space_3d(met, met->w, p, lon, lat, &w, ci, cw, 0); \
877 intpol_met_space_3d(met, met->pv, p, lon, lat, &pv, ci, cw, 0); \
878 intpol_met_space_3d(met, met->h2o, p, lon, lat, &h2o, ci, cw, 0); \
879 intpol_met_space_3d(met, met->o3, p, lon, lat, &o3, ci, cw, 0); \
880 intpol_met_space_3d(met, met->lwc, p, lon, lat, &lwc, ci, cw, 0); \
881 intpol_met_space_3d(met, met->rwc, p, lon, lat, &rwc, ci, cw, 0); \
882 intpol_met_space_3d(met, met->iwc, p, lon, lat, &iwc, ci, cw, 0); \
883 intpol_met_space_3d(met, met->swc, p, lon, lat, &swc, ci, cw, 0); \
884 intpol_met_space_3d(met, met->cc, p, lon, lat, &cc, ci, cw, 0); \
885 intpol_met_space_2d(met, met->ps, lon, lat, &ps, ci, cw, 0); \
886 intpol_met_space_2d(met, met->ts, lon, lat, &ts, ci, cw, 0); \
887 intpol_met_space_2d(met, met->zs, lon, lat, &zs, ci, cw, 0); \
888 intpol_met_space_2d(met, met->us, lon, lat, &us, ci, cw, 0); \
889 intpol_met_space_2d(met, met->vs, lon, lat, &vs, ci, cw, 0); \
890 intpol_met_space_2d(met, met->ess, ess, lat, &ess, ci, cw, 0); \
891 intpol_met_space_2d(met, met->nss, nss, lat, &nss, ci, cw, 0); \
892 intpol_met_space_2d(met, met->shf, shf, lat, &shf, ci, cw, 0); \
893 intpol_met_space_2d(met, met->lsm, lon, lat, &lsm, ci, cw, 0); \
894 intpol_met_space_2d(met, met->sst, lon, lat, &sst, ci, cw, 0); \
895 intpol_met_space_2d(met, met->pbl, lon, lat, &pbl, ci, cw, 0); \
896 intpol_met_space_2d(met, met->pt, lon, lat, &pt, ci, cw, 0); \
897 intpol_met_space_2d(met, met->tt, lon, lat, &tt, ci, cw, 0); \
898 intpol_met_space_2d(met, met->zt, lon, lat, &zt, ci, cw, 0); \
899 intpol_met_space_2d(met, met->h2ot, lon, lat, &h2ot, ci, cw, 0); \
900 intpol_met_space_2d(met, met->pct, lon, lat, &pct, ci, cw, 0); \
901 intpol_met_space_2d(met, met->pcb, lon, lat, &pcb, ci, cw, 0); \
902 intpol_met_space_2d(met, met->cl, lon, lat, &cl, ci, cw, 0); \
903 intpol_met_space_2d(met, met->plcl, lon, lat, &plcl, ci, cw, 0); \
904 intpol_met_space_2d(met, met->plfc, lon, lat, &plfc, ci, cw, 0); \
905 intpol_met_space_2d(met, met->pel, lon, lat, &pel, ci, cw, 0); \
906 intpol_met_space_2d(met, met->cape, lon, lat, &cape, ci, cw, 0); \
907 intpol_met_space_2d(met, met->cin, lon, lat, &cin, ci, cw, 0); \
908 intpol_met_space_2d(met, met->o3c, lon, lat, &o3c, ci, cw, 0); \
909 }
910
925#define INTPOL_TIME_ALL(time, p, lon, lat) { \
926 intpol_met_time_3d(met0, met0->z, met1, met1->z, time, p, lon, lat, &z, ci, cw, 1); \
927 intpol_met_time_3d(met0, met0->t, met1, met1->t, time, p, lon, lat, &t, ci, cw, 0); \
928 intpol_met_time_3d(met0, met0->u, met1, met1->u, time, p, lon, lat, &u, ci, cw, 0); \
929 intpol_met_time_3d(met0, met0->v, met1, met1->v, time, p, lon, lat, &v, ci, cw, 0); \
930 intpol_met_time_3d(met0, met0->w, met1, met1->w, time, p, lon, lat, &w, ci, cw, 0); \
931 intpol_met_time_3d(met0, met0->pv, met1, met1->pv, time, p, lon, lat, &pv, ci, cw, 0); \
932 intpol_met_time_3d(met0, met0->h2o, met1, met1->h2o, time, p, lon, lat, &h2o, ci, cw, 0); \
933 intpol_met_time_3d(met0, met0->o3, met1, met1->o3, time, p, lon, lat, &o3, ci, cw, 0); \
934 intpol_met_time_3d(met0, met0->lwc, met1, met1->lwc, time, p, lon, lat, &lwc, ci, cw, 0); \
935 intpol_met_time_3d(met0, met0->rwc, met1, met1->rwc, time, p, lon, lat, &rwc, ci, cw, 0); \
936 intpol_met_time_3d(met0, met0->iwc, met1, met1->iwc, time, p, lon, lat, &iwc, ci, cw, 0); \
937 intpol_met_time_3d(met0, met0->swc, met1, met1->swc, time, p, lon, lat, &swc, ci, cw, 0); \
938 intpol_met_time_3d(met0, met0->cc, met1, met1->cc, time, p, lon, lat, &cc, ci, cw, 0); \
939 intpol_met_time_2d(met0, met0->ps, met1, met1->ps, time, lon, lat, &ps, ci, cw, 0); \
940 intpol_met_time_2d(met0, met0->ts, met1, met1->ts, time, lon, lat, &ts, ci, cw, 0); \
941 intpol_met_time_2d(met0, met0->zs, met1, met1->zs, time, lon, lat, &zs, ci, cw, 0); \
942 intpol_met_time_2d(met0, met0->us, met1, met1->us, time, lon, lat, &us, ci, cw, 0); \
943 intpol_met_time_2d(met0, met0->vs, met1, met1->vs, time, lon, lat, &vs, ci, cw, 0); \
944 intpol_met_time_2d(met0, met0->ess, met1, met1->ess, time, lon, lat, &ess, ci, cw, 0); \
945 intpol_met_time_2d(met0, met0->nss, met1, met1->nss, time, lon, lat, &nss, ci, cw, 0); \
946 intpol_met_time_2d(met0, met0->shf, met1, met1->shf, time, lon, lat, &shf, ci, cw, 0); \
947 intpol_met_time_2d(met0, met0->lsm, met1, met1->lsm, time, lon, lat, &lsm, ci, cw, 0); \
948 intpol_met_time_2d(met0, met0->sst, met1, met1->sst, time, lon, lat, &sst, ci, cw, 0); \
949 intpol_met_time_2d(met0, met0->pbl, met1, met1->pbl, time, lon, lat, &pbl, ci, cw, 0); \
950 intpol_met_time_2d(met0, met0->pt, met1, met1->pt, time, lon, lat, &pt, ci, cw, 0); \
951 intpol_met_time_2d(met0, met0->tt, met1, met1->tt, time, lon, lat, &tt, ci, cw, 0); \
952 intpol_met_time_2d(met0, met0->zt, met1, met1->zt, time, lon, lat, &zt, ci, cw, 0); \
953 intpol_met_time_2d(met0, met0->h2ot, met1, met1->h2ot, time, lon, lat, &h2ot, ci, cw, 0); \
954 intpol_met_time_2d(met0, met0->pct, met1, met1->pct, time, lon, lat, &pct, ci, cw, 0); \
955 intpol_met_time_2d(met0, met0->pcb, met1, met1->pcb, time, lon, lat, &pcb, ci, cw, 0); \
956 intpol_met_time_2d(met0, met0->cl, met1, met1->cl, time, lon, lat, &cl, ci, cw, 0); \
957 intpol_met_time_2d(met0, met0->plcl, met1, met1->plcl, time, lon, lat, &plcl, ci, cw, 0); \
958 intpol_met_time_2d(met0, met0->plfc, met1, met1->plfc, time, lon, lat, &plfc, ci, cw, 0); \
959 intpol_met_time_2d(met0, met0->pel, met1, met1->pel, time, lon, lat, &pel, ci, cw, 0); \
960 intpol_met_time_2d(met0, met0->cape, met1, met1->cape, time, lon, lat, &cape, ci, cw, 0); \
961 intpol_met_time_2d(met0, met0->cin, met1, met1->cin, time, lon, lat, &cin, ci, cw, 0); \
962 intpol_met_time_2d(met0, met0->o3c, met1, met1->o3c, time, lon, lat, &o3c, ci, cw, 0); \
963 }
964
979#define LAPSE(p1, t1, p2, t2) \
980 (1e3 * G0 / RA * ((t2) - (t1)) / ((t2) + (t1)) \
981 * ((p2) + (p1)) / ((p2) - (p1)))
982
998#define LIN(x0, y0, x1, y1, x) \
999 ((y0)+((y1)-(y0))/((x1)-(x0))*((x)-(x0)))
1000
1025#define MAX(a,b) \
1026 (((a)>(b))?(a):(b))
1027
1039#define MET_HEADER \
1040 fprintf(out, \
1041 "# $1 = time [s]\n" \
1042 "# $2 = altitude [km]\n" \
1043 "# $3 = longitude [deg]\n" \
1044 "# $4 = latitude [deg]\n" \
1045 "# $5 = pressure [hPa]\n" \
1046 "# $6 = temperature [K]\n" \
1047 "# $7 = zonal wind [m/s]\n" \
1048 "# $8 = meridional wind [m/s]\n" \
1049 "# $9 = vertical velocity [hPa/s]\n" \
1050 "# $10 = H2O volume mixing ratio [ppv]\n"); \
1051 fprintf(out, \
1052 "# $11 = O3 volume mixing ratio [ppv]\n" \
1053 "# $12 = geopotential height [km]\n" \
1054 "# $13 = potential vorticity [PVU]\n" \
1055 "# $14 = surface pressure [hPa]\n" \
1056 "# $15 = surface temperature [K]\n" \
1057 "# $16 = surface geopotential height [km]\n" \
1058 "# $17 = surface zonal wind [m/s]\n" \
1059 "# $18 = surface meridional wind [m/s]\n" \
1060 "# $19 = eastward turbulent surface stress [N/m^2]\n" \
1061 "# $20 = northward turbulent surface stress [N/m^2]\n"); \
1062 fprintf(out, \
1063 "# $21 = surface sensible heat flux [W/m^2]\n" \
1064 "# $22 = land-sea mask [1]\n" \
1065 "# $23 = sea surface temperature [K]\n" \
1066 "# $24 = tropopause pressure [hPa]\n" \
1067 "# $25 = tropopause geopotential height [km]\n" \
1068 "# $26 = tropopause temperature [K]\n" \
1069 "# $27 = tropopause water vapor [ppv]\n" \
1070 "# $28 = cloud liquid water content [kg/kg]\n" \
1071 "# $29 = cloud rain water content [kg/kg]\n" \
1072 "# $30 = cloud ice water content [kg/kg]\n"); \
1073 fprintf(out, \
1074 "# $31 = cloud snow water content [kg/kg]\n" \
1075 "# $32 = cloud cover [1]\n" \
1076 "# $33 = total column cloud water [kg/m^2]\n" \
1077 "# $34 = cloud top pressure [hPa]\n" \
1078 "# $35 = cloud bottom pressure [hPa]\n" \
1079 "# $36 = pressure at lifted condensation level (LCL) [hPa]\n" \
1080 "# $37 = pressure at level of free convection (LFC) [hPa]\n" \
1081 "# $38 = pressure at equilibrium level (EL) [hPa]\n" \
1082 "# $39 = convective available potential energy (CAPE) [J/kg]\n" \
1083 "# $40 = convective inhibition (CIN) [J/kg]\n"); \
1084 fprintf(out, \
1085 "# $41 = relative humidity over water [%%]\n" \
1086 "# $42 = relative humidity over ice [%%]\n" \
1087 "# $43 = dew point temperature [K]\n" \
1088 "# $44 = frost point temperature [K]\n" \
1089 "# $45 = NAT temperature [K]\n" \
1090 "# $46 = HNO3 volume mixing ratio [ppv]\n" \
1091 "# $47 = OH volume mixing ratio [ppv]\n" \
1092 "# $48 = H2O2 volume mixing ratio [ppv]\n" \
1093 "# $49 = HO2 volume mixing ratio [ppv]\n" \
1094 "# $50 = O(1D) volume mixing ratio [ppv]\n"); \
1095 fprintf(out, \
1096 "# $51 = boundary layer pressure [hPa]\n" \
1097 "# $52 = total column ozone [DU]\n" \
1098 "# $53 = number of data points\n" \
1099 "# $54 = number of tropopause data points\n" \
1100 "# $55 = number of CAPE data points\n");
1101
1126#define MIN(a,b) \
1127 (((a)<(b))?(a):(b))
1128
1141#define MOLEC_DENS(p,t) \
1142 (AVO * 1e-6 * ((p) * 100) / (RI * (t)))
1143
1155#define NC(cmd) { \
1156 int nc_result=(cmd); \
1157 if(nc_result!=NC_NOERR) \
1158 ERRMSG("%s", nc_strerror(nc_result)); \
1159 }
1160
1184#define NC_DEF_VAR(varname, type, ndims, dims, long_name, units, level, quant) { \
1185 NC(nc_def_var(ncid, varname, type, ndims, dims, &varid)); \
1186 NC(nc_put_att_text(ncid, varid, "long_name", strnlen(long_name, LEN), long_name)); \
1187 NC(nc_put_att_text(ncid, varid, "units", strnlen(units, LEN), units)); \
1188 if((quant) > 0) \
1189 NC(nc_def_var_quantize(ncid, varid, NC_QUANTIZE_GRANULARBR, quant)); \
1190 if((level) != 0) { \
1191 NC(nc_def_var_deflate(ncid, varid, 1, 1, level)); \
1192 /* unsigned int ulevel = (unsigned int)level; */ \
1193 /* NC(nc_def_var_filter(ncid, varid, 32015, 1, (unsigned int[]){ulevel})); */ \
1194 } \
1195 }
1196
1214#define NC_GET_DOUBLE(varname, ptr, force) { \
1215 if(force) { \
1216 NC(nc_inq_varid(ncid, varname, &varid)); \
1217 NC(nc_get_var_double(ncid, varid, ptr)); \
1218 } else { \
1219 if(nc_inq_varid(ncid, varname, &varid) == NC_NOERR) { \
1220 NC(nc_get_var_double(ncid, varid, ptr)); \
1221 } else \
1222 WARN("netCDF variable %s is missing!", varname); \
1223 } \
1224 }
1225
1242#define NC_INQ_DIM(dimname, ptr, min, max) { \
1243 int dimid; size_t naux; \
1244 NC(nc_inq_dimid(ncid, dimname, &dimid)); \
1245 NC(nc_inq_dimlen(ncid, dimid, &naux)); \
1246 *ptr = (int)naux; \
1247 if ((*ptr) < (min) || (*ptr) > (max)) \
1248 ERRMSG("Dimension %s is out of range!", dimname); \
1249 }
1250
1265#define NC_PUT_DOUBLE(varname, ptr, hyperslab) { \
1266 NC(nc_inq_varid(ncid, varname, &varid)); \
1267 if(hyperslab) { \
1268 NC(nc_put_vara_double(ncid, varid, start, count, ptr)); \
1269 } else { \
1270 NC(nc_put_var_double(ncid, varid, ptr)); \
1271 } \
1272 }
1273
1289#define NC_PUT_FLOAT(varname, ptr, hyperslab) { \
1290 NC(nc_inq_varid(ncid, varname, &varid)); \
1291 if(hyperslab) { \
1292 NC(nc_put_vara_float(ncid, varid, start, count, ptr)); \
1293 } else { \
1294 NC(nc_put_var_float(ncid, varid, ptr)); \
1295 } \
1296 }
1297
1312#define NC_PUT_INT(varname, ptr, hyperslab) { \
1313 NC(nc_inq_varid(ncid, varname, &varid)); \
1314 if(hyperslab) { \
1315 NC(nc_put_vara_int(ncid, varid, start, count, ptr)); \
1316 } else { \
1317 NC(nc_put_var_int(ncid, varid, ptr)); \
1318 } \
1319 }
1320
1334#define NC_PUT_ATT(varname, attname, text) { \
1335 NC(nc_inq_varid(ncid, varname, &varid)); \
1336 NC(nc_put_att_text(ncid, varid, attname, strnlen(text, LEN), text)); \
1337 }
1338
1351#define NC_PUT_ATT_GLOBAL(attname, text) \
1352 NC(nc_put_att_text(ncid, NC_GLOBAL, attname, strnlen(text, LEN), text));
1353
1371#define NN(x0, y0, x1, y1, x) \
1372 (fabs((x) - (x0)) <= fabs((x) - (x1)) ? (y0) : (y1))
1373
1389#ifdef _OPENACC
1390#define PARTICLE_LOOP(ip0, ip1, check_dt, ...) \
1391 const int ip0_const = ip0; \
1392 const int ip1_const = ip1; \
1393 _Pragma(__VA_ARGS__) \
1394 _Pragma("acc parallel loop independent gang vector") \
1395 for (int ip = ip0_const; ip < ip1_const; ip++) \
1396 if (!check_dt || cache->dt[ip] != 0)
1397#else
1398#define PARTICLE_LOOP(ip0, ip1, check_dt, ...) \
1399 const int ip0_const = ip0; \
1400 const int ip1_const = ip1; \
1401 _Pragma("omp parallel for default(shared)") \
1402 for (int ip = ip0_const; ip < ip1_const; ip++) \
1403 if (!check_dt || cache->dt[ip] != 0)
1404#endif
1405
1428#define P(z) \
1429 (P0 * exp(-(z) / H0))
1430
1452#define PSAT(t) \
1453 (6.112 * exp(17.62 * ((t) - T0) / (243.12 + (t) - T0)))
1454
1476#define PSICE(t) \
1477 (6.112 * exp(22.46 * ((t) - T0) / (272.62 + (t) - T0)))
1478
1503#define PW(p, h2o) \
1504 ((p) * MAX((h2o), 0.1e-6) / (1. + (1. - EPS) * MAX((h2o), 0.1e-6)))
1505
1520#define RAD2DEG(rad) \
1521 ((rad) * (180.0 / M_PI))
1522
1550#define RH(p, t, h2o) \
1551 (PW(p, h2o) / PSAT(t) * 100.)
1552
1580#define RHICE(p, t, h2o) \
1581 (PW(p, h2o) / PSICE(t) * 100.)
1582
1605#define RHO(p, t) \
1606 (100. * (p) / (RA * (t)))
1607
1624#define SET_ATM(qnt, val) \
1625 if (ctl->qnt >= 0) \
1626 atm->q[ctl->qnt][ip] = val;
1627
1647#define SET_QNT(qnt, name, longname, unit) \
1648 if (strcasecmp(ctl->qnt_name[iq], name) == 0) { \
1649 ctl->qnt = iq; \
1650 sprintf(ctl->qnt_longname[iq], longname); \
1651 sprintf(ctl->qnt_unit[iq], unit); \
1652 } else
1653
1668#define SH(h2o) \
1669 (EPS * MAX((h2o), 0.1e-6))
1670
1681#define SQR(x) \
1682 ((x)*(x))
1683
1695#define SWAP(x, y, type) \
1696 do {type tmp = x; x = y; y = tmp;} while(0);
1697
1719#define TDEW(p, h2o) \
1720 (T0 + 243.12 * log(PW((p), (h2o)) / 6.112) \
1721 / (17.62 - log(PW((p), (h2o)) / 6.112)))
1722
1744#define TICE(p, h2o) \
1745 (T0 + 272.62 * log(PW((p), (h2o)) / 6.112) \
1746 / (22.46 - log(PW((p), (h2o)) / 6.112)))
1747
1768#define THETA(p, t) \
1769 ((t) * pow(1000. / (p), 0.286))
1770
1797#define THETAVIRT(p, t, h2o) \
1798 (TVIRT(THETA((p), (t)), MAX((h2o), 0.1e-6)))
1799
1818#define TOK(line, tok, format, var) { \
1819 if(((tok)=strtok((line), " \t"))) { \
1820 if(sscanf(tok, format, &(var))!=1) continue; \
1821 } else ERRMSG("Error while reading!"); \
1822 }
1823
1843#define TVIRT(t, h2o) \
1844 ((t) * (1. + (1. - EPS) * MAX((h2o), 0.1e-6)))
1845
1865#define Z(p) \
1866 (H0 * log(P0 / (p)))
1867
1896#define ZDIFF(lnp0, t0, h2o0, lnp1, t1, h2o1) \
1897 (RI / MA / G0 * 0.5 * (TVIRT((t0), (h2o0)) + TVIRT((t1), (h2o1))) \
1898 * ((lnp0) - (lnp1)))
1899
1915#define ZETA(ps, p, t) \
1916 (((p) / (ps) <= 0.3 ? 1. : \
1917 sin(M_PI / 2. * (1. - (p) / (ps)) / (1. - 0.3))) \
1918 * THETA((p), (t)))
1919
1920/* ------------------------------------------------------------
1921 Log messages...
1922 ------------------------------------------------------------ */
1923
1925#ifndef LOGLEV
1926#define LOGLEV 2
1927#endif
1928
1958#define LOG(level, ...) { \
1959 if(level >= 2) \
1960 printf(" "); \
1961 if(level <= LOGLEV) { \
1962 printf(__VA_ARGS__); \
1963 printf("\n"); \
1964 } \
1965 }
1966
1995#define WARN(...) { \
1996 printf("\nWarning (%s, %s, l%d): ", __FILE__, __func__, __LINE__); \
1997 LOG(0, __VA_ARGS__); \
1998 }
1999
2028#define ERRMSG(...) { \
2029 printf("\nError (%s, %s, l%d): ", __FILE__, __func__, __LINE__); \
2030 LOG(0, __VA_ARGS__); \
2031 exit(EXIT_FAILURE); \
2032 }
2033
2063#define PRINT(format, var) \
2064 printf("Print (%s, %s, l%d): %s= "format"\n", \
2065 __FILE__, __func__, __LINE__, #var, var);
2066
2067/* ------------------------------------------------------------
2068 Timers...
2069 ------------------------------------------------------------ */
2070
2072#define NTIMER 100
2073
2087#define PRINT_TIMERS \
2088 timer("END", "END", 1);
2089
2108#define SELECT_TIMER(id, group, color) { \
2109 NVTX_POP; \
2110 NVTX_PUSH(id, color); \
2111 timer(id, group, 0); \
2112 }
2113
2127#define START_TIMERS \
2128 NVTX_PUSH("START", NVTX_CPU);
2129
2142#define STOP_TIMERS \
2143 NVTX_POP;
2144
2145/* ------------------------------------------------------------
2146 NVIDIA Tools Extension (NVTX)...
2147 ------------------------------------------------------------ */
2148
2149#ifdef NVTX
2150#include "nvToolsExt.h"
2151
2153#define NVTX_CPU 0xFFADD8E6
2154
2156#define NVTX_GPU 0xFF00008B
2157
2159#define NVTX_H2D 0xFFFFFF00
2160
2162#define NVTX_D2H 0xFFFF8800
2163
2165#define NVTX_READ 0xFFFFCCCB
2166
2168#define NVTX_WRITE 0xFF8B0000
2169
2171#define NVTX_RECV 0xFFCCFFCB
2172
2174#define NVTX_SEND 0xFF008B00
2175
2205#define NVTX_PUSH(range_title, range_color) { \
2206 nvtxEventAttributes_t eventAttrib = {0}; \
2207 eventAttrib.version = NVTX_VERSION; \
2208 eventAttrib.size = NVTX_EVENT_ATTRIB_STRUCT_SIZE; \
2209 eventAttrib.messageType = NVTX_MESSAGE_TYPE_ASCII; \
2210 eventAttrib.colorType = NVTX_COLOR_ARGB; \
2211 eventAttrib.color = range_color; \
2212 eventAttrib.message.ascii = range_title; \
2213 nvtxRangePushEx(&eventAttrib); \
2214 }
2215
2228#define NVTX_POP { \
2229 nvtxRangePop(); \
2230 }
2231#else
2232
2233/* Empty definitions of NVTX_PUSH and NVTX_POP... */
2234#define NVTX_PUSH(range_title, range_color) {}
2235#define NVTX_POP {}
2236#endif
2237
2238/* ------------------------------------------------------------
2239 Structs...
2240 ------------------------------------------------------------ */
2241
2249typedef struct {
2250
2251 /* ------------------------------------------------------------
2252 Quantity parameters...
2253 ------------------------------------------------------------ */
2254
2256 int nq;
2257
2259 char qnt_name[NQ][LEN];
2260
2262 char qnt_longname[NQ][LEN];
2263
2265 char qnt_unit[NQ][LEN];
2266
2268 char qnt_format[NQ][LEN];
2269
2272
2275
2278
2281
2284
2287
2290
2293
2296
2299
2302
2305
2308
2311
2314
2317
2320
2323
2326
2329
2332
2335
2338
2341
2344
2347
2350
2353
2356
2359
2362
2365
2368
2371
2374
2377
2380
2383
2386
2389
2392
2395
2398
2401
2404
2407
2410
2413
2416
2419
2422
2425
2428
2431
2434
2437
2440
2443
2446
2449
2452
2455
2458
2461
2464
2467
2470
2473
2476
2479
2482
2485
2488
2491
2494
2497
2500
2503
2506
2509
2512
2515
2518
2521
2524
2527
2530
2533
2536
2539
2542
2545
2548
2551
2554
2556 double t_start;
2557
2559 double t_stop;
2560
2562 double dt_mod;
2563
2564 /* ------------------------------------------------------------
2565 Meteo data parameters...
2566 ------------------------------------------------------------ */
2567
2569 char metbase[LEN];
2570
2572 double dt_met;
2573
2576
2580
2584
2587
2590
2593
2596
2599
2601 int met_comp_prec[METVAR];
2602
2604 double met_comp_tol[METVAR];
2605
2608
2611
2615
2618
2621
2624
2627
2630
2633
2636
2639
2642
2645
2648
2651
2654
2657
2660
2663
2666
2669
2672
2675
2678
2680 double met_p[EP];
2681
2684
2687
2689 double met_lev_hyam[EP];
2690
2692 double met_lev_hybm[EP];
2693
2696
2699
2702
2705
2708
2711
2714
2718
2721
2724
2727
2730
2733
2736
2737 /* ------------------------------------------------------------
2738 Geophysical module parameters...
2739 ------------------------------------------------------------ */
2740
2742 double sort_dt;
2743
2747
2749 char balloon[LEN];
2750
2753
2757
2760
2763
2766
2769
2772
2775
2778
2781
2784
2787
2790
2793
2796
2798 double conv_cin;
2799
2801 double conv_dt;
2802
2805
2808
2811
2814
2817
2820
2822 double bound_p0;
2823
2825 double bound_p1;
2826
2829
2832
2835
2838
2840 char species[LEN];
2841
2843 double molmass;
2844
2847
2850
2853
2855 char clim_hno3_filename[LEN];
2856
2858 char clim_oh_filename[LEN];
2859
2861 char clim_h2o2_filename[LEN];
2862
2864 char clim_ho2_filename[LEN];
2865
2867 char clim_o1d_filename[LEN];
2868
2870 char clim_o3_filename[LEN];
2871
2873 char clim_ccl4_timeseries[LEN];
2874
2876 char clim_ccl3f_timeseries[LEN];
2877
2879 char clim_ccl2f2_timeseries[LEN];
2880
2882 char clim_n2o_timeseries[LEN];
2883
2885 char clim_sf6_timeseries[LEN];
2886
2889
2892
2895
2898
2901
2904
2907
2910
2913
2916
2919
2922
2925
2928
2931
2934
2937
2940
2943
2946
2949
2952
2954 double oh_chem[4];
2955
2958
2961
2964
2966 double dt_kpp;
2967
2970
2972 double wet_depo_pre[2];
2973
2976
2979
2982
2985
2987 double wet_depo_ic_h[2];
2988
2990 double wet_depo_bc_h[2];
2991
2994
2997
3000
3003
3006
3008 double psc_h2o;
3009
3011 double psc_hno3;
3012
3013 /* ------------------------------------------------------------
3014 Output parameters...
3015 ------------------------------------------------------------ */
3016
3018 char atm_basename[LEN];
3019
3021 char atm_gpfile[LEN];
3022
3025
3028
3031
3035
3040
3043
3045 int atm_nc_quant[NQ];
3046
3049
3051 char csi_basename[LEN];
3052
3054 char csi_kernel[LEN];
3055
3058
3060 char csi_obsfile[LEN];
3061
3064
3067
3070
3072 double csi_z0;
3073
3075 double csi_z1;
3076
3079
3081 double csi_lon0;
3082
3084 double csi_lon1;
3085
3088
3090 double csi_lat0;
3091
3093 double csi_lat1;
3094
3096 int nens;
3097
3099 char ens_basename[LEN];
3100
3103
3105 char grid_basename[LEN];
3106
3108 char grid_kernel[LEN];
3109
3111 char grid_gpfile[LEN];
3112
3115
3118
3121
3123 int grid_nc_quant[NQ];
3124
3127
3130
3132 double grid_z0;
3133
3135 double grid_z1;
3136
3139
3142
3145
3148
3151
3154
3157
3159 char prof_basename[LEN];
3160
3162 char prof_obsfile[LEN];
3163
3166
3168 double prof_z0;
3169
3171 double prof_z1;
3172
3175
3178
3181
3184
3187
3190
3192 char sample_basename[LEN];
3193
3195 char sample_kernel[LEN];
3196
3198 char sample_obsfile[LEN];
3199
3202
3205
3207 char stat_basename[LEN];
3208
3210 double stat_lon;
3211
3213 double stat_lat;
3214
3216 double stat_r;
3217
3219 double stat_t0;
3220
3222 double stat_t1;
3223
3225 char vtk_basename[LEN];
3226
3229
3232
3235
3238
3241
3244
3247
3250
3253
3254} ctl_t;
3255
3264typedef struct {
3265
3267 int np;
3268
3270 double time[NP];
3271
3273 double p[NP];
3274
3276 double lon[NP];
3277
3279 double lat[NP];
3280
3282 double q[NQ][NP];
3283
3284} atm_t;
3285
3293typedef struct {
3294
3296 double time;
3297
3299 double p;
3300
3302 double lon;
3303
3305 double lat;
3306
3308 double q[NQ];
3309
3310} particle_t;
3311
3319typedef struct {
3320
3322 int rank;
3323
3325 int size;
3326
3327#ifdef DD
3329 int neighbours[DD_NNMAX];
3330
3332 MPI_Datatype MPI_Particle;
3333#endif
3334
3335} mpi_info_t;
3336
3343typedef struct {
3344
3346 double iso_var[NP];
3347
3349 double iso_ps[NP];
3350
3352 double iso_ts[NP];
3353
3356
3358 float uvwp[NP][3];
3359
3361 double rs[3 * NP + 1];
3362
3364 double dt[NP];
3365
3366} cache_t;
3367
3375typedef struct {
3376
3378 int np;
3379
3381 int nsza;
3382
3384 int no3c;
3385
3387 double p[CP];
3388
3390 double sza[CSZA];
3391
3393 double o3c[CO3];
3394
3396 double n2o[CP][CSZA][CO3];
3397
3399 double ccl4[CP][CSZA][CO3];
3400
3402 double ccl3f[CP][CSZA][CO3];
3403
3405 double ccl2f2[CP][CSZA][CO3];
3406
3408 double o2[CP][CSZA][CO3];
3409
3411 double o3_1[CP][CSZA][CO3];
3412
3414 double o3_2[CP][CSZA][CO3];
3415
3417 double h2o2[CP][CSZA][CO3];
3418
3420 double h2o[CP][CSZA][CO3];
3421
3422} clim_photo_t;
3423
3431typedef struct {
3432
3435
3437 double time[CTS];
3438
3440 double vmr[CTS];
3441
3442} clim_ts_t;
3443
3451typedef struct {
3452
3455
3457 int nlat;
3458
3460 int np;
3461
3463 double time[CT];
3464
3466 double lat[CY];
3467
3469 double p[CP];
3470
3472 double vmr[CT][CP][CY];
3473
3474} clim_zm_t;
3475
3483typedef struct {
3484
3487
3490
3492 double tropo_time[12];
3493
3495 double tropo_lat[73];
3496
3498 double tropo[12][73];
3499
3502
3505
3508
3511
3514
3517
3520
3523
3526
3529
3532
3533} clim_t;
3534
3542typedef struct {
3543
3545 double time;
3546
3548 int nx;
3549
3551 int ny;
3552
3554 int np;
3555
3557 int npl;
3558
3559 // TODO:
3560 // They need global sizes now, maybe in the future just keep EX, EY, EP and
3561 // Introduce help data structure in read_met_grid etc.
3562
3564#ifdef DD
3565 double lon[EX_GLOB];
3566#else
3567 double lon[EX];
3568#endif
3569
3571#ifdef DD
3572 double lat[EY_GLOB];
3573#else
3574 double lat[EY];
3575#endif
3576
3578#ifdef DD
3579 double p[EP_GLOB];
3580#else
3581 double p[EP];
3582#endif
3583
3585 double hybrid[EP];
3586
3588 float ps[EX][EY];
3589
3591 float ts[EX][EY];
3592
3594 float zs[EX][EY];
3595
3597 float us[EX][EY];
3598
3600 float vs[EX][EY];
3601
3603 float ess[EX][EY];
3604
3606 float nss[EX][EY];
3607
3609 float shf[EX][EY];
3610
3612 float lsm[EX][EY];
3613
3615 float sst[EX][EY];
3616
3618 float pbl[EX][EY];
3619
3621 float pt[EX][EY];
3622
3624 float tt[EX][EY];
3625
3627 float zt[EX][EY];
3628
3630 float h2ot[EX][EY];
3631
3633 float pct[EX][EY];
3634
3636 float pcb[EX][EY];
3637
3639 float cl[EX][EY];
3640
3642 float plcl[EX][EY];
3643
3645 float plfc[EX][EY];
3646
3648 float pel[EX][EY];
3649
3651 float cape[EX][EY];
3652
3654 float cin[EX][EY];
3655
3657 float o3c[EX][EY];
3658
3660 float z[EX][EY][EP];
3661
3663 float t[EX][EY][EP];
3664
3666 float u[EX][EY][EP];
3667
3669 float v[EX][EY][EP];
3670
3672 float w[EX][EY][EP];
3673
3675 float pv[EX][EY][EP];
3676
3678 float h2o[EX][EY][EP];
3679
3681 float o3[EX][EY][EP];
3682
3684 float lwc[EX][EY][EP];
3685
3687 float rwc[EX][EY][EP];
3688
3690 float iwc[EX][EY][EP];
3691
3693 float swc[EX][EY][EP];
3694
3696 float cc[EX][EY][EP];
3697
3699 float pl[EX][EY][EP];
3700
3702 float ul[EX][EY][EP];
3703
3705 float vl[EX][EY][EP];
3706
3708 float wl[EX][EY][EP];
3709
3711 float zetal[EX][EY][EP];
3712
3714 float zeta_dotl[EX][EY][EP];
3715
3716 // TODO: Integrate this into a grid_t ?
3717
3720
3723
3726
3729
3731 size_t subdomain_start[4];
3732
3734 size_t subdomain_count[4];
3735
3736 /* Hyperslab of boundary halos start. */
3737 size_t halo_bnd_start[4];
3738
3739 /* Hyperslab of boundary halos count. */
3740 size_t halo_bnd_count[4];
3741
3742 /* Hyperslab of boundary halos count. */
3744
3745 /* Hyperslab of boundary halos count. */
3747
3750
3753
3756
3757} met_t;
3758
3759/* ------------------------------------------------------------
3760 OpenACC routines...
3761 ------------------------------------------------------------ */
3762
3763#ifdef _OPENACC
3764#pragma acc routine (clim_oh)
3765#pragma acc routine (clim_photo)
3766#pragma acc routine (clim_tropo)
3767#pragma acc routine (clim_ts)
3768#pragma acc routine (clim_zm)
3769#pragma acc routine (intpol_check_lon_lat)
3770#pragma acc routine (intpol_met_4d_coord)
3771#pragma acc routine (intpol_met_space_3d)
3772#pragma acc routine (intpol_met_space_3d_ml)
3773#pragma acc routine (intpol_met_space_2d)
3774#pragma acc routine (intpol_met_time_3d)
3775#pragma acc routine (intpol_met_time_3d_ml)
3776#pragma acc routine (intpol_met_time_2d)
3777#pragma acc routine (kernel_weight)
3778#pragma acc routine (lapse_rate)
3779#pragma acc routine (locate_irr)
3780#pragma acc routine (locate_irr_float)
3781#pragma acc routine (locate_reg)
3782#pragma acc routine (locate_vert)
3783#pragma acc routine (nat_temperature)
3784#pragma acc routine (pbl_weight)
3785#pragma acc routine (sedi)
3786#pragma acc routine (stddev)
3787#pragma acc routine (sza_calc)
3788#pragma acc routine (tropo_weight)
3789#endif
3790
3791/* ------------------------------------------------------------
3792 Functions...
3793 ------------------------------------------------------------ */
3794
3817 void *data,
3818 size_t N);
3819
3834void cart2geo(
3835 const double *x,
3836 double *z,
3837 double *lon,
3838 double *lat);
3839
3862double clim_oh(
3863 const ctl_t * ctl,
3864 const clim_t * clim,
3865 const double t,
3866 const double lon,
3867 const double lat,
3868 const double p);
3869
3889 const ctl_t * ctl,
3890 clim_t * clim);
3891
3921double clim_photo(
3922 const double rate[CP][CSZA][CO3],
3923 const clim_photo_t * photo,
3924 const double p,
3925 const double sza,
3926 const double o3c);
3927
3953double clim_tropo(
3954 const clim_t * clim,
3955 const double t,
3956 const double lat);
3957
3976void clim_tropo_init(
3977 clim_t * clim);
3978
3995double clim_ts(
3996 const clim_ts_t * ts,
3997 const double t);
3998
4020double clim_zm(
4021 const clim_zm_t * zm,
4022 const double t,
4023 const double lat,
4024 const double p);
4025
4067 const ctl_t * ctl,
4068 const char *varname,
4069 float *array,
4070 const size_t nx,
4071 const size_t ny,
4072 const size_t np,
4073 const int decompress,
4074 FILE * inout);
4075
4108void compress_pck(
4109 const char *varname,
4110 float *array,
4111 const size_t nxy,
4112 const size_t nz,
4113 const int decompress,
4114 FILE * inout);
4115
4146 const char *varname,
4147 float *array,
4148 const int nx,
4149 const int ny,
4150 const int nz,
4151 const int precision,
4152 const double tolerance,
4153 const int decompress,
4154 FILE * inout);
4155
4195 const char *varname,
4196 float *array,
4197 const int nx,
4198 const int ny,
4199 const int nz,
4200 const int precision,
4201 const double tolerance,
4202 const int decompress,
4203 FILE * inout);
4204
4227 const char *varname,
4228 float *array,
4229 const size_t n,
4230 const int decompress,
4231 const int level,
4232 FILE * inout);
4233
4256void day2doy(
4257 const int year,
4258 const int mon,
4259 const int day,
4260 int *doy);
4261
4283void doy2day(
4284 const int year,
4285 const int doy,
4286 int *mon,
4287 int *day);
4288
4315void fft_help(
4316 double *fcReal,
4317 double *fcImag,
4318 const int n);
4319
4346void geo2cart(
4347 const double z,
4348 const double lon,
4349 const double lat,
4350 double *x);
4351
4376void get_met_help(
4377 const ctl_t * ctl,
4378 const double t,
4379 const int direct,
4380 const char *metbase,
4381 const double dt_met,
4382 char *filename);
4383
4407void get_met_replace(
4408 char *orig,
4409 char *search,
4410 char *repl);
4411
4448void get_tropo(
4449 const int met_tropo,
4450 ctl_t * ctl,
4451 clim_t * clim,
4452 met_t * met,
4453 const double *lons,
4454 const int nx,
4455 const double *lats,
4456 const int ny,
4457 double *pt,
4458 double *zt,
4459 double *tt,
4460 double *qt,
4461 double *o3t,
4462 double *ps,
4463 double *zs);
4464
4486 const double *lons,
4487 const int nlon,
4488 const double *lats,
4489 const int nlat,
4490 const double lon,
4491 const double lat,
4492 double *lon2,
4493 double *lat2);
4494
4537 const met_t * met0,
4538 float height0[EX][EY][EP],
4539 float array0[EX][EY][EP],
4540 const met_t * met1,
4541 float height1[EX][EY][EP],
4542 float array1[EX][EY][EP],
4543 const double ts,
4544 const double height,
4545 const double lon,
4546 const double lat,
4547 double *var,
4548 int *ci,
4549 double *cw,
4550 const int init);
4551
4587 const met_t * met,
4588 float array[EX][EY][EP],
4589 const double p,
4590 const double lon,
4591 const double lat,
4592 double *var,
4593 int *ci,
4594 double *cw,
4595 const int init);
4596
4616 const met_t * met,
4617 float zs[EX][EY][EP],
4618 float vals[EX][EY][EP],
4619 const double z,
4620 const double lon,
4621 const double lat,
4622 double *val);
4623
4659 const met_t * met,
4660 float array[EX][EY],
4661 const double lon,
4662 const double lat,
4663 double *var,
4664 int *ci,
4665 double *cw,
4666 const int init);
4667
4702 const met_t * met0,
4703 float array0[EX][EY][EP],
4704 const met_t * met1,
4705 float array1[EX][EY][EP],
4706 const double ts,
4707 const double p,
4708 const double lon,
4709 const double lat,
4710 double *var,
4711 int *ci,
4712 double *cw,
4713 const int init);
4714
4738 const met_t * met0,
4739 float zs0[EX][EY][EP],
4740 float array0[EX][EY][EP],
4741 const met_t * met1,
4742 float zs1[EX][EY][EP],
4743 float array1[EX][EY][EP],
4744 const double ts,
4745 const double p,
4746 const double lon,
4747 const double lat,
4748 double *var);
4749
4785 const met_t * met0,
4786 float array0[EX][EY],
4787 const met_t * met1,
4788 float array1[EX][EY],
4789 const double ts,
4790 const double lon,
4791 const double lat,
4792 double *var,
4793 int *ci,
4794 double *cw,
4795 const int init);
4796
4834void intpol_tropo_3d(
4835 const double time0,
4836 float array0[EX][EY],
4837 const double time1,
4838 float array1[EX][EY],
4839 const double lons[EX],
4840 const double lats[EY],
4841 const int nlon,
4842 const int nlat,
4843 const double time,
4844 const double lon,
4845 const double lat,
4846 const int method,
4847 double *var,
4848 double *sigma);
4849
4876void jsec2time(
4877 const double jsec,
4878 int *year,
4879 int *mon,
4880 int *day,
4881 int *hour,
4882 int *min,
4883 int *sec,
4884 double *remain);
4885
4912double kernel_weight(
4913 const double kz[EP],
4914 const double kw[EP],
4915 const int nk,
4916 const double p);
4917
4956double lapse_rate(
4957 const double t,
4958 const double h2o);
4959
4989 ctl_t * ctl);
4990
5010int locate_irr(
5011 const double *xx,
5012 const int n,
5013 const double x);
5014
5041 const float *xx,
5042 const int n,
5043 const double x,
5044 const int ig);
5045
5066int locate_reg(
5067 const double *xx,
5068 const int n,
5069 const double x);
5070
5092void locate_vert(
5093 float profiles[EX][EY][EP],
5094 const int np,
5095 const int lon_ap_ind,
5096 const int lat_ap_ind,
5097 const double alt_ap,
5098 int *ind);
5099
5125void module_advect(
5126 const ctl_t * ctl,
5127 const cache_t * cache,
5128 met_t * met0,
5129 met_t * met1,
5130 atm_t * atm);
5131
5155 const ctl_t * ctl,
5156 const cache_t * cache,
5157 met_t * met0,
5158 met_t * met1,
5159 atm_t * atm);
5160
5198 const ctl_t * ctl,
5199 const cache_t * cache,
5200 const clim_t * clim,
5201 met_t * met0,
5202 met_t * met1,
5203 atm_t * atm);
5204
5246void module_chem_grid(
5247 const ctl_t * ctl,
5248 met_t * met0,
5249 met_t * met1,
5250 atm_t * atm,
5251 const double t);
5252
5279void module_chem_init(
5280 const ctl_t * ctl,
5281 const cache_t * cache,
5282 const clim_t * clim,
5283 met_t * met0,
5284 met_t * met1,
5285 atm_t * atm);
5286
5311 const ctl_t * ctl,
5312 cache_t * cache,
5313 met_t * met0,
5314 met_t * met1,
5315 atm_t * atm);
5316
5343void module_decay(
5344 const ctl_t * ctl,
5345 const cache_t * cache,
5346 const clim_t * clim,
5347 atm_t * atm);
5348
5385void module_diff_meso(
5386 const ctl_t * ctl,
5387 cache_t * cache,
5388 met_t * met0,
5389 met_t * met1,
5390 atm_t * atm);
5391
5425void module_diff_pbl(
5426 const ctl_t * ctl,
5427 cache_t * cache,
5428 met_t * met0,
5429 met_t * met1,
5430 atm_t * atm);
5431
5486void module_diff_turb(
5487 const ctl_t * ctl,
5488 cache_t * cache,
5489 const clim_t * clim,
5490 met_t * met0,
5491 met_t * met1,
5492 atm_t * atm);
5493
5513void module_dry_depo(
5514 const ctl_t * ctl,
5515 const cache_t * cache,
5516 met_t * met0,
5517 met_t * met1,
5518 atm_t * atm);
5519
5552void module_h2o2_chem(
5553 const ctl_t * ctl,
5554 const cache_t * cache,
5555 const clim_t * clim,
5556 met_t * met0,
5557 met_t * met1,
5558 atm_t * atm);
5559
5580 const ctl_t * ctl,
5581 cache_t * cache,
5582 met_t * met0,
5583 met_t * met1,
5584 atm_t * atm);
5585
5603void module_isosurf(
5604 const ctl_t * ctl,
5605 const cache_t * cache,
5606 met_t * met0,
5607 met_t * met1,
5608 atm_t * atm);
5609
5642 ctl_t * ctl,
5643 cache_t * cache,
5644 clim_t * clim,
5645 met_t * met0,
5646 met_t * met1,
5647 atm_t * atm);
5648
5667void module_meteo(
5668 const ctl_t * ctl,
5669 const cache_t * cache,
5670 const clim_t * clim,
5671 met_t * met0,
5672 met_t * met1,
5673 atm_t * atm);
5674
5692void module_mixing(
5693 const ctl_t * ctl,
5694 const clim_t * clim,
5695 atm_t * atm,
5696 const double t);
5697
5726 const ctl_t * ctl,
5727 const clim_t * clim,
5728 atm_t * atm,
5729 const int *ixs,
5730 const int *iys,
5731 const int *izs,
5732 const int qnt_idx,
5733 const int use_ensemble);
5734
5767void module_oh_chem(
5768 const ctl_t * ctl,
5769 const cache_t * cache,
5770 const clim_t * clim,
5771 met_t * met0,
5772 met_t * met1,
5773 atm_t * atm);
5774
5802void module_position(
5803 const cache_t * cache,
5804 met_t * met0,
5805 met_t * met1,
5806 atm_t * atm);
5807
5832void module_rng_init(
5833 const int ntask);
5834
5860void module_rng(
5861 const ctl_t * ctl,
5862 double *rs,
5863 const size_t n,
5864 const int method);
5865
5888void module_sedi(
5889 const ctl_t * ctl,
5890 const cache_t * cache,
5891 met_t * met0,
5892 met_t * met1,
5893 atm_t * atm);
5894
5918void module_sort(
5919 const ctl_t * ctl,
5920 met_t * met0,
5921 atm_t * atm);
5922
5942void module_sort_help(
5943 double *a,
5944 const int *p,
5945 const int np);
5946
5970void module_timesteps(
5971 const ctl_t * ctl,
5972 cache_t * cache,
5973 met_t * met0,
5974 atm_t * atm,
5975 const double t);
5976
5998 ctl_t * ctl,
5999 const atm_t * atm);
6000
6034 const ctl_t * ctl,
6035 const cache_t * cache,
6036 const clim_t * clim,
6037 met_t * met0,
6038 met_t * met1,
6039 atm_t * atm);
6040
6070void module_wet_depo(
6071 const ctl_t * ctl,
6072 const cache_t * cache,
6073 met_t * met0,
6074 met_t * met1,
6075 atm_t * atm);
6076
6107void mptrac_alloc(
6108 ctl_t ** ctl,
6109 cache_t ** cache,
6110 clim_t ** clim,
6111 met_t ** met0,
6112 met_t ** met1,
6113 atm_t ** atm);
6114
6144#ifdef DD
6145void mptrac_free(
6146 ctl_t * ctl,
6147 cache_t * cache,
6148 clim_t * clim,
6149 met_t * met0,
6150 met_t * met1,
6151 atm_t * atm,
6152 mpi_info_t * mpi_info);
6153#else
6154void mptrac_free(
6155 ctl_t * ctl,
6156 cache_t * cache,
6157 clim_t * clim,
6158 met_t * met0,
6159 met_t * met1,
6160 atm_t * atm);
6161#endif
6162
6197void mptrac_get_met(
6198 ctl_t * ctl,
6199 clim_t * clim,
6200 const double t,
6201 met_t ** met0,
6202 met_t ** met1);
6203
6223void mptrac_init(
6224 ctl_t * ctl,
6225 cache_t * cache,
6226 clim_t * clim,
6227 atm_t * atm,
6228 const int ntask);
6229
6265int mptrac_read_atm(
6266 const char *filename,
6267 const ctl_t * ctl,
6268 atm_t * atm);
6269
6301void mptrac_read_clim(
6302 const ctl_t * ctl,
6303 clim_t * clim);
6304
6334void mptrac_read_ctl(
6335 const char *filename,
6336 int argc,
6337 char *argv[],
6338 ctl_t * ctl);
6339
6369int mptrac_read_met(
6370 const char *filename,
6371 const ctl_t * ctl,
6372 const clim_t * clim,
6373 met_t * met);
6374
6395#ifdef DD
6397 ctl_t * ctl,
6398 cache_t * cache,
6399 clim_t * clim,
6400 met_t ** met0,
6401 met_t ** met1,
6402 atm_t * atm,
6403 double t,
6404 mpi_info_t * mpi_info);
6405#else
6407 ctl_t * ctl,
6408 cache_t * cache,
6409 clim_t * clim,
6410 met_t ** met0,
6411 met_t ** met1,
6412 atm_t * atm,
6413 double t);
6414#endif
6415
6445void mptrac_write_atm(
6446 const char *filename,
6447 const ctl_t * ctl,
6448 const atm_t * atm,
6449 const double t);
6450
6488void mptrac_write_met(
6489 const char *filename,
6490 const ctl_t * ctl,
6491 met_t * met);
6492
6527 const char *dirname,
6528 const ctl_t * ctl,
6529 met_t * met0,
6530 met_t * met1,
6531 atm_t * atm,
6532 const double t);
6533
6565 const ctl_t * ctl,
6566 const cache_t * cache,
6567 const clim_t * clim,
6568 met_t ** met0,
6569 met_t ** met1,
6570 const atm_t * atm);
6571
6602 const ctl_t * ctl,
6603 const cache_t * cache,
6604 const clim_t * clim,
6605 met_t ** met0,
6606 met_t ** met1,
6607 const atm_t * atm);
6608
6636double nat_temperature(
6637 const double p,
6638 const double h2o,
6639 const double hno3);
6640
6661double pbl_weight(
6662 const ctl_t * ctl,
6663 const atm_t * atm,
6664 const int ip,
6665 const double pbl,
6666 const double ps);
6667
6700int read_atm_asc(
6701 const char *filename,
6702 const ctl_t * ctl,
6703 atm_t * atm);
6704
6735int read_atm_bin(
6736 const char *filename,
6737 const ctl_t * ctl,
6738 atm_t * atm);
6739
6764int read_atm_clams(
6765 const char *filename,
6766 const ctl_t * ctl,
6767 atm_t * atm);
6768
6798int read_atm_nc(
6799 const char *filename,
6800 const ctl_t * ctl,
6801 atm_t * atm);
6802
6831void read_clim_photo(
6832 const char *filename,
6833 clim_photo_t * photo);
6834
6852 const int ncid,
6853 const char *varname,
6854 const clim_photo_t * photo,
6855 double var[CP][CSZA][CO3]);
6856
6880int read_clim_ts(
6881 const char *filename,
6882 clim_ts_t * ts);
6883
6910void read_clim_zm(
6911 const char *filename,
6912 const char *varname,
6913 clim_zm_t * zm);
6914
6942void read_kernel(
6943 const char *filename,
6944 double kz[EP],
6945 double kw[EP],
6946 int *nk);
6947
6979int read_met_bin(
6980 const char *filename,
6981 const ctl_t * ctl,
6982 met_t * met);
6983
7009void read_met_bin_2d(
7010 FILE * in,
7011 const met_t * met,
7012 float var[EX][EY],
7013 const char *varname);
7014
7052void read_met_bin_3d(
7053 FILE * in,
7054 const ctl_t * ctl,
7055 const met_t * met,
7056 float var[EX][EY][EP],
7057 const char *varname,
7058 const float bound_min,
7059 const float bound_max);
7060
7087void read_met_cape(
7088 const ctl_t * ctl,
7089 const clim_t * clim,
7090 met_t * met);
7091
7114void read_met_cloud(
7115 met_t * met);
7116
7142void read_met_detrend(
7143 const ctl_t * ctl,
7144 met_t * met);
7145
7169 met_t * met);
7170
7197void read_met_geopot(
7198 const ctl_t * ctl,
7199 met_t * met);
7200
7228 const char *filename,
7229 const ctl_t * ctl,
7230 met_t * met);
7231
7252#ifdef ECCODES
7253void read_met_grib_grid(
7254 codes_handle ** handles,
7255 int count_handles,
7256 met_t * met);
7257#endif
7258
7278#ifdef ECCODES
7279void read_met_grib_levels(
7280 codes_handle ** handles,
7281 const int num_messages,
7282 const ctl_t * ctl,
7283 met_t * met);
7284#endif
7285
7315#ifdef ECCODES
7316void read_met_grib_surface(
7317 codes_handle ** handles,
7318 const int num_messages,
7319 const ctl_t * ctl,
7320 met_t * met);
7321#endif
7322
7351void read_met_ml2pl(
7352 const ctl_t * ctl,
7353 const met_t * met,
7354 float var[EX][EY][EP],
7355 const char *varname);
7356
7379 const ctl_t * ctl,
7380 met_t * met);
7381
7411int read_met_nc(
7412 const char *filename,
7413 const ctl_t * ctl,
7414 met_t * met);
7415
7447 const int ncid,
7448 const ctl_t * ctl,
7449 met_t * met);
7450
7492 const int ncid,
7493 const ctl_t * ctl,
7494 met_t * met);
7495
7527void read_met_nc_grid(
7528 const char *filename,
7529 const int ncid,
7530 const ctl_t * ctl,
7531 met_t * met);
7532
7563 const char *filename,
7564 const ctl_t * ctl,
7565 met_t * met);
7566
7600 const char *filename,
7601 const int ncid,
7602 const ctl_t * ctl,
7603 met_t * met);
7604
7636 const int ncid,
7637 const ctl_t * ctl,
7638 met_t * met);
7639
7669 const int ncid,
7670 const ctl_t * ctl,
7671 met_t * met);
7672
7702int read_met_nc_2d(
7703 const int ncid,
7704 const char *varname,
7705 const char *varname2,
7706 const char *varname3,
7707 const char *varname4,
7708 const char *varname5,
7709 const char *varname6,
7710 const ctl_t * ctl,
7711 const met_t * met,
7712 float dest[EX][EY],
7713 const float scl,
7714 const int init);
7715
7746 const int ncid,
7747 const char *varname,
7748 const char *varname2,
7749 const char *varname3,
7750 const char *varname4,
7751 const char *varname5,
7752 const char *varname6,
7753 const ctl_t * ctl,
7754 const met_t * met,
7755 float dest[EX][EY],
7756 const float scl,
7757 const int init);
7758
7789int read_met_nc_3d(
7790 const int ncid,
7791 const char *varname,
7792 const char *varname2,
7793 const char *varname3,
7794 const char *varname4,
7795 const ctl_t * ctl,
7796 const met_t * met,
7797 float dest[EX][EY][EP],
7798 const float scl);
7799
7831 const int ncid,
7832 const char *varname,
7833 const char *varname2,
7834 const char *varname3,
7835 const char *varname4,
7836 const ctl_t * ctl,
7837 const met_t * met,
7838 float dest[EX][EY][EP],
7839 const float scl);
7840
7886void read_met_pbl(
7887 const ctl_t * ctl,
7888 met_t * met);
7889
7922 met_t * met);
7923
7954 met_t * met);
7955
7986void read_met_pv(
7987 met_t * met);
7988
8011void read_met_ozone(
8012 met_t * met);
8013
8042void read_met_sample(
8043 const ctl_t * ctl,
8044 met_t * met);
8045
8074void read_met_tropo(
8075 const ctl_t * ctl,
8076 const clim_t * clim,
8077 met_t * met);
8078
8110void read_obs(
8111 const char *filename,
8112 const ctl_t * ctl,
8113 double *rt,
8114 double *rz,
8115 double *rlon,
8116 double *rlat,
8117 double *robs,
8118 int *nobs);
8119
8147void read_obs_asc(
8148 const char *filename,
8149 double *rt,
8150 double *rz,
8151 double *rlon,
8152 double *rlat,
8153 double *robs,
8154 int *nobs);
8155
8182void read_obs_nc(
8183 const char *filename,
8184 double *rt,
8185 double *rz,
8186 double *rlon,
8187 double *rlat,
8188 double *robs,
8189 int *nobs);
8190
8224double scan_ctl(
8225 const char *filename,
8226 int argc,
8227 char *argv[],
8228 const char *varname,
8229 const int arridx,
8230 const char *defvalue,
8231 char *value);
8232
8259double sedi(
8260 const double p,
8261 const double T,
8262 const double rp,
8263 const double rhop);
8264
8294void spline(
8295 const double *x,
8296 const double *y,
8297 const int n,
8298 const double *x2,
8299 double *y2,
8300 const int n2,
8301 const int method);
8302
8325float stddev(
8326 const float *data,
8327 const int n);
8328
8348double sza_calc(
8349 const double sec,
8350 const double lon,
8351 const double lat);
8352
8377void time2jsec(
8378 const int year,
8379 const int mon,
8380 const int day,
8381 const int hour,
8382 const int min,
8383 const int sec,
8384 const double remain,
8385 double *jsec);
8386
8415void timer(
8416 const char *name,
8417 const char *group,
8418 const int output);
8419
8445double time_from_filename(
8446 const char *filename,
8447 const int offset);
8448
8466double tropo_weight(
8467 const clim_t * clim,
8468 const atm_t * atm,
8469 const int ip);
8470
8493void write_atm_asc(
8494 const char *filename,
8495 const ctl_t * ctl,
8496 const atm_t * atm,
8497 const double t);
8498
8522void write_atm_bin(
8523 const char *filename,
8524 const ctl_t * ctl,
8525 const atm_t * atm);
8526
8550void write_atm_clams(
8551 const char *filename,
8552 const ctl_t * ctl,
8553 const atm_t * atm);
8554
8580 const char *dirname,
8581 const ctl_t * ctl,
8582 const atm_t * atm,
8583 const double t);
8584
8608void write_atm_nc(
8609 const char *filename,
8610 const ctl_t * ctl,
8611 const atm_t * atm);
8612
8641void write_csi(
8642 const char *filename,
8643 const ctl_t * ctl,
8644 const atm_t * atm,
8645 const double t);
8646
8688 const char *filename,
8689 const ctl_t * ctl,
8690 const atm_t * atm,
8691 const double t);
8692
8720void write_ens(
8721 const char *filename,
8722 const ctl_t * ctl,
8723 const atm_t * atm,
8724 const double t);
8725
8764void write_grid(
8765 const char *filename,
8766 const ctl_t * ctl,
8767 met_t * met0,
8768 met_t * met1,
8769 const atm_t * atm,
8770 const double t);
8771
8817void write_grid_asc(
8818 const char *filename,
8819 const ctl_t * ctl,
8820 const double *cd,
8821 double *mean[NQ],
8822 double *sigma[NQ],
8823 const double *vmr_impl,
8824 const double t,
8825 const double *z,
8826 const double *lon,
8827 const double *lat,
8828 const double *area,
8829 const double dz,
8830 const int *np);
8831
8874void write_grid_nc(
8875 const char *filename,
8876 const ctl_t * ctl,
8877 const double *cd,
8878 double *mean[NQ],
8879 double *sigma[NQ],
8880 const double *vmr_impl,
8881 const double t,
8882 const double *z,
8883 const double *lon,
8884 const double *lat,
8885 const double *area,
8886 const double dz,
8887 const int *np);
8888
8918void write_met_bin(
8919 const char *filename,
8920 const ctl_t * ctl,
8921 met_t * met);
8922
8950void write_met_bin_2d(
8951 FILE * out,
8952 met_t * met,
8953 float var[EX][EY],
8954 const char *varname);
8955
8993void write_met_bin_3d(
8994 FILE * out,
8995 const ctl_t * ctl,
8996 met_t * met,
8997 float var[EX][EY][EP],
8998 const char *varname,
8999 const int precision,
9000 const double tolerance);
9001
9029void write_met_nc(
9030 const char *filename,
9031 const ctl_t * ctl,
9032 met_t * met);
9033
9058void write_met_nc_2d(
9059 const int ncid,
9060 const char *varname,
9061 met_t * met,
9062 float var[EX][EY],
9063 const float scl);
9064
9090void write_met_nc_3d(
9091 const int ncid,
9092 const char *varname,
9093 met_t * met,
9094 float var[EX][EY][EP],
9095 const float scl);
9096
9127void write_prof(
9128 const char *filename,
9129 const ctl_t * ctl,
9130 met_t * met0,
9131 met_t * met1,
9132 const atm_t * atm,
9133 const double t);
9134
9166void write_sample(
9167 const char *filename,
9168 const ctl_t * ctl,
9169 met_t * met0,
9170 met_t * met1,
9171 const atm_t * atm,
9172 const double t);
9173
9200void write_station(
9201 const char *filename,
9202 const ctl_t * ctl,
9203 atm_t * atm,
9204 const double t);
9205
9234void write_vtk(
9235 const char *filename,
9236 const ctl_t * ctl,
9237 const atm_t * atm,
9238 const double t);
9239
9263 atm_t * atm,
9264 particle_t * particles,
9265 ctl_t * ctl,
9266 int *nparticles,
9267 cache_t * cache,
9268 int rank);
9269
9295 atm_t * atm,
9296 particle_t * particles,
9297 ctl_t * ctl,
9298 int *nparticles,
9299 cache_t * cache);
9300
9322#ifdef DD
9323void dd_register_MPI_type_particle(
9324 MPI_Datatype * MPI_Particle);
9325#endif
9326
9349#ifdef DD
9350void dd_get_rect_neighbour(
9351 const ctl_t ctl,
9352 mpi_info_t * mpi_info);
9353#endif
9354
9384#ifdef DD
9385void dd_communicate_particles(
9386 particle_t * particles,
9387 int *nparticles,
9388 MPI_Datatype MPI_Particle,
9389 int *neighbours,
9390 int nneighbours,
9391 ctl_t ctl);
9392#endif
9393
9421#ifdef DD
9422void dd_assign_rect_subdomains_atm(
9423 atm_t * atm,
9424 met_t * met,
9425 ctl_t * ctl,
9426 mpi_info_t * mpi_info,
9427 int init);
9428#endif
9429
9456#ifdef DD
9457void dd_init(
9458 ctl_t * ctl,
9459 mpi_info_t * mpi_info,
9460 atm_t * atm,
9461 met_t ** met,
9462 int *dd_init);
9463#endif
9464
9492#ifdef DD
9493void module_dd(
9494 ctl_t * ctl,
9495 atm_t * atm,
9496 cache_t * cache,
9497 mpi_info_t * mpi_info,
9498 met_t ** met);
9499#endif
9500
9534#ifdef DD
9535void dd_sort(
9536 const ctl_t * ctl,
9537 met_t * met0,
9538 atm_t * atm,
9539 int *nparticles,
9540 int *rank);
9541#endif
9542
9543#endif /* LIBTRAC_H */
void read_met_geopot(const ctl_t *ctl, met_t *met)
Calculates geopotential heights from meteorological data.
Definition: mptrac.c:7283
void read_met_nc_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:9679
void compress_zstd(const char *varname, float *array, const size_t n, const int decompress, const int level, FILE *inout)
Compresses or decompresses a float array using Zstandard (ZSTD).
#define LEN
Maximum length of ASCII data lines.
Definition: mptrac.h:318
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:5939
void day2doy(const int year, const int mon, const int day, int *doy)
Get day of year from date.
Definition: mptrac.c:1004
void read_met_extrapolate(met_t *met)
Extrapolates meteorological data.
Definition: mptrac.c:7243
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:11209
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:12612
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:4325
void read_met_sample(const ctl_t *ctl, met_t *met)
Downsamples meteorological data based on specified parameters.
Definition: mptrac.c:10207
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:12358
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:10552
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:2194
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:4042
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:5558
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:3361
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:6371
void read_met_cloud(met_t *met)
Calculates cloud-related variables for each grid point.
Definition: mptrac.c:7082
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:2747
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:10725
#define METVAR
Number of 3-D meteorological variables.
Definition: mptrac.h:323
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:1568
int read_met_nc_dd(const char *filename, const ctl_t *ctl, met_t *met)
Reads meteorological data from a NetCDF file and processes it.
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:1501
int read_met_nc_3d_dd(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.
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:6338
void write_csi_ens(const char *filename, const ctl_t *ctl, const atm_t *atm, const double t)
Writes ensemble-based Critical Success Index (CSI) and other verification statistics to an output fil...
void read_met_pbl(const ctl_t *ctl, met_t *met)
Computes the planetary boundary layer (PBL) pressure based on meteorological data.
Definition: mptrac.c:9815
void read_met_detrend(const ctl_t *ctl, met_t *met)
Detrends meteorological data.
Definition: mptrac.c:7139
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:9015
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:10380
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:10596
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:2632
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:2155
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:1201
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:2174
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:4414
void read_met_monotonize(const ctl_t *ctl, met_t *met)
Makes zeta and pressure profiles monotone.
Definition: mptrac.c:8899
#define EP_GLOB
Maximum number of pressure levels for meteo data.
Definition: mptrac.h:303
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:6490
void read_met_periodic(met_t *met)
Applies periodic boundary conditions to meteorological data along longitudinal axis.
Definition: mptrac.c:9952
void module_timesteps_init(ctl_t *ctl, const atm_t *atm)
Initialize start time and time interval for time-stepping.
Definition: mptrac.c:4089
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:11695
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:3466
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:406
#define EY
Maximum number of latitudes for meteo data.
Definition: mptrac.h:298
void module_mixing_help(const ctl_t *ctl, const clim_t *clim, atm_t *atm, const int *ixs, const int *iys, const int *izs, const int qnt_idx, const int use_ensemble)
Perform subgrid-scale interparcel mixing of a given quantity.
Definition: mptrac.c:3537
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:6462
void read_met_ml2pl(const ctl_t *ctl, const met_t *met, float var[EX][EY][EP], const char *varname)
Reads global meteorological information from a grib file.
Definition: mptrac.c:8857
double clim_tropo(const clim_t *clim, const double t, const double lat)
Calculates the tropopause pressure based on climatological data.
Definition: mptrac.c:205
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:10624
void dd_particles2atm(atm_t *atm, particle_t *particles, ctl_t *ctl, int *nparticles, cache_t *cache)
Converts particle data to atmospheric data.
void read_met_bin_2d(FILE *in, const met_t *met, float var[EX][EY], const char *varname)
Reads a 2-dimensional meteorological variable from a binary file and stores it in the provided array.
Definition: mptrac.c:6832
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:2091
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:3185
void level_definitions(ctl_t *ctl)
Defines pressure levels for meteorological data.
Definition: mptrac.c:1881
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:11983
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:5827
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:10867
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:4535
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:1626
void fft_help(double *fcReal, double *fcImag, const int n)
Computes the Fast Fourier Transform (FFT) of a complex sequence.
Definition: mptrac.c:1053
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:4190
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:673
#define EY_GLOB
Maximum number of global latitudes for meteo data.
Definition: mptrac.h:313
double nat_temperature(const double p, const double h2o, const double hno3)
Calculates the nitric acid trihydrate (NAT) temperature.
Definition: mptrac.c:6134
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:10758
#define CP
Maximum number of pressure levels for climatological data.
Definition: mptrac.h:368
#define NQ
Maximum number of quantities per data point.
Definition: mptrac.h:333
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:148
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:6544
void read_met_nc_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:9331
#define EX
Maximum number of longitudes for meteo data.
Definition: mptrac.h:293
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:3910
void timer(const char *name, const char *group, const int output)
Measures and reports elapsed time for named and grouped timers.
Definition: mptrac.c:10898
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:11024
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:1443
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:1655
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:2674
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:6643
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:2374
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:10653
#define CT
Maximum number of time steps for climatological data.
Definition: mptrac.h:378
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:2347
int read_met_grib(const char *filename, const ctl_t *ctl, met_t *met)
Reads meteorological data from a grib file and processes it.
void module_sort_help(double *a, const int *p, const int np)
Reorder an array based on a given permutation.
Definition: mptrac.c:4004
float stddev(const float *data, const int n)
Calculates the standard deviation of a set of data.
Definition: mptrac.c:10805
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:1713
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:9179
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:6861
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:2121
double time_from_filename(const char *filename, const int offset)
Extracts and converts a timestamp from a filename to Julian seconds.
Definition: mptrac.c:10966
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:12671
void dd_atm2particles(atm_t *atm, particle_t *particles, ctl_t *ctl, int *nparticles, cache_t *cache, int rank)
Converts atmospheric data to particle data.
void mptrac_read_clim(const ctl_t *ctl, clim_t *clim)
Reads various climatological data and populates the given climatology structure.
Definition: mptrac.c:4625
void module_chem_grid(const ctl_t *ctl, met_t *met0, met_t *met1, atm_t *atm, const double t)
Computes gridded chemical tracer concentrations (volume mixing ratio) from individual air parcel mass...
Definition: mptrac.c:2470
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:11001
void read_met_nc_grid_dd(const char *filename, const int ncid, const ctl_t *ctl, met_t *met)
Reads meteorological grid data from NetCDF files with domain decomposition.
void write_met_nc(const char *filename, const ctl_t *ctl, met_t *met)
Writes meteorological data to a NetCDF file.
Definition: mptrac.c:12448
void module_rng_init(const int ntask)
Initialize random number generators for parallel tasks.
Definition: mptrac.c:3774
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:4554
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:1271
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:5883
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:6043
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:383
void read_met_nc_surface_dd(const int ncid, const ctl_t *ctl, met_t *met)
Reads and processes surface meteorological data from NetCDF files with domain decomposition.
void read_met_ozone(met_t *met)
Calculates the total column ozone from meteorological ozone data.
Definition: mptrac.c:10178
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:4382
void clim_tropo_init(clim_t *clim)
Initializes the tropopause data in the climatology structure.
Definition: mptrac.c:233
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:3805
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:1110
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:13060
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:328
double sza_calc(const double sec, const double lon, const double lat)
Calculates the solar zenith angle.
Definition: mptrac.c:10826
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:1023
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:1680
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:3723
void clim_oh_diurnal_correction(const ctl_t *ctl, clim_t *clim)
Applies a diurnal correction to the hydroxyl radical (OH) concentration in climatology data.
Definition: mptrac.c:115
void write_met_bin_2d(FILE *out, met_t *met, float var[EX][EY], const char *varname)
Writes a 2-dimensional meteorological variable to a binary file.
Definition: mptrac.c:12329
void read_met_pv(met_t *met)
Calculates potential vorticity (PV) from meteorological data.
Definition: mptrac.c:10072
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:6222
void get_met_replace(char *orig, char *search, char *repl)
Replaces occurrences of a substring in a string with another substring.
Definition: mptrac.c:1177
#define CY
Maximum number of latitudes for climatological data.
Definition: mptrac.h:358
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:2786
void module_sort(const ctl_t *ctl, met_t *met0, atm_t *atm)
Sort particles according to box index.
Definition: mptrac.c:3939
int read_met_nc(const char *filename, const ctl_t *ctl, met_t *met)
Reads meteorological data from a NetCDF file and processes it.
Definition: mptrac.c:8984
double clim_ts(const clim_ts_t *ts, const double t)
Interpolates a time series of climatological variables.
Definition: mptrac.c:388
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:1804
int read_met_bin(const char *filename, const ctl_t *ctl, met_t *met)
Reads meteorological data from a binary file.
Definition: mptrac.c:6684
int read_met_nc_2d_dd(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.
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:5682
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:11156
void read_met_nc_levels_dd(const int ncid, const ctl_t *ctl, met_t *met)
Reads and processes meteorological level data from NetCDF files with domain decomposition.
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:2988
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:6278
#define EX_GLOB
Maximum number of global longitudes for meteo data.
Definition: mptrac.h:308
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:13146
#define EP
Maximum number of pressure levels for meteo data.
Definition: mptrac.h:288
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:4120
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:4685
void read_met_polar_winds(met_t *met)
Applies a fix for polar winds in meteorological data.
Definition: mptrac.c:10013
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:3103
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:12087
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:6158
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:2863
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:12641
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:3255
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:363
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:3639
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:1092
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:1837
void compress_sz3(const char *varname, float *array, const int nx, const int ny, const int nz, const int precision, const double tolerance, const int decompress, FILE *inout)
Compresses or decompresses a 3-D float array using the SZ3 library.
void read_met_nc_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:9465
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:6180
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:1244
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:12898
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:11792
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:3040
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:12217
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:11106
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:6967
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:5999
#define CSZA
Maximum number of solar zenith angles for climatological data.
Definition: mptrac.h:373
double lapse_rate(const double t, const double h2o)
Calculates the moist adiabatic lapse rate in Kelvin per kilometer.
Definition: mptrac.c:1863
#define DD_NNMAX
Maximum number of neighbours to communicate with.
Definition: mptrac.h:393
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:11416
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:11367
Air parcel data.
Definition: mptrac.h:3264
int np
Number of air parcels.
Definition: mptrac.h:3267
Cache data structure.
Definition: mptrac.h:3343
int iso_n
Isosurface balloon number of data points.
Definition: mptrac.h:3355
Climatological data in the form of photolysis rates.
Definition: mptrac.h:3375
int nsza
Number of solar zenith angles.
Definition: mptrac.h:3381
int np
Number of pressure levels.
Definition: mptrac.h:3378
int no3c
Number of total ozone columns.
Definition: mptrac.h:3384
Climatological data.
Definition: mptrac.h:3483
clim_ts_t ccl2f2
CFC-12 time series.
Definition: mptrac.h:3525
clim_photo_t photo
Photolysis rates.
Definition: mptrac.h:3501
clim_zm_t ho2
HO2 zonal means.
Definition: mptrac.h:3513
clim_zm_t hno3
HNO3 zonal means.
Definition: mptrac.h:3504
int tropo_ntime
Number of tropopause timesteps.
Definition: mptrac.h:3486
clim_ts_t sf6
SF6 time series.
Definition: mptrac.h:3531
clim_ts_t ccl4
CFC-10 time series.
Definition: mptrac.h:3519
clim_ts_t ccl3f
CFC-11 time series.
Definition: mptrac.h:3522
clim_zm_t o1d
O(1D) zonal means.
Definition: mptrac.h:3516
clim_zm_t h2o2
H2O2 zonal means.
Definition: mptrac.h:3510
int tropo_nlat
Number of tropopause latitudes.
Definition: mptrac.h:3489
clim_zm_t oh
OH zonal means.
Definition: mptrac.h:3507
clim_ts_t n2o
N2O time series.
Definition: mptrac.h:3528
Climatological data in the form of time series.
Definition: mptrac.h:3431
int ntime
Number of timesteps.
Definition: mptrac.h:3434
Climatological data in the form of zonal means.
Definition: mptrac.h:3451
int np
Number of pressure levels.
Definition: mptrac.h:3460
int ntime
Number of timesteps.
Definition: mptrac.h:3454
int nlat
Number of latitudes.
Definition: mptrac.h:3457
Control parameters.
Definition: mptrac.h:2249
double grid_z0
Lower altitude of gridded data [km].
Definition: mptrac.h:3132
int qnt_o3
Quantity array index for ozone volume mixing ratio.
Definition: mptrac.h:2361
double csi_lat1
Upper latitude of gridded CSI data [deg].
Definition: mptrac.h:3093
int qnt_Coh
Quantity array index for OH volume mixing ratio (chemistry code).
Definition: mptrac.h:2511
double wet_depo_ic_a
Coefficient A for wet deposition in cloud (exponential form).
Definition: mptrac.h:2981
int met_nc_scale
Check netCDF scaling factors (0=no, 1=yes).
Definition: mptrac.h:2589
int qnt_pel
Quantity array index for pressure at equilibrium level (EL).
Definition: mptrac.h:2394
int csi_nz
Number of altitudes of gridded CSI data.
Definition: mptrac.h:3069
double molmass
Molar mass [g/mol].
Definition: mptrac.h:2843
int qnt_p
Quantity array index for pressure.
Definition: mptrac.h:2340
int qnt_Cccl2f2
Quantity array index for CFC-12 volume mixing ratio (chemistry code).
Definition: mptrac.h:2535
int dd_halos_size
Size of halos given in grid-points.
Definition: mptrac.h:3252
int mixing_nx
Number of longitudes of mixing grid.
Definition: mptrac.h:2906
double chemgrid_z1
Upper altitude of chemistry grid [km].
Definition: mptrac.h:2930
int qnt_m
Quantity array index for mass.
Definition: mptrac.h:2280
int qnt_aoa
Quantity array index for age of air.
Definition: mptrac.h:2544
int qnt_rhop
Quantity array index for particle density.
Definition: mptrac.h:2289
int qnt_swc
Quantity array index for cloud snow water content.
Definition: mptrac.h:2373
double csi_obsmin
Minimum observation index to trigger detection.
Definition: mptrac.h:3063
int qnt_pcb
Quantity array index for cloud bottom pressure.
Definition: mptrac.h:2382
double bound_dzs
Boundary conditions surface layer depth [km].
Definition: mptrac.h:2831
double csi_lon1
Upper longitude of gridded CSI data [deg].
Definition: mptrac.h:3084
int qnt_u
Quantity array index for zonal wind.
Definition: mptrac.h:2349
double stat_lon
Longitude of station [deg].
Definition: mptrac.h:3210
double mixing_trop
Interparcel exchange parameter for mixing in the troposphere.
Definition: mptrac.h:2891
double sort_dt
Time step for sorting of particle data [s].
Definition: mptrac.h:2742
double mixing_z1
Upper altitude of mixing grid [km].
Definition: mptrac.h:2903
double stat_r
Search radius around station [km].
Definition: mptrac.h:3216
double wet_depo_bc_a
Coefficient A for wet deposition below cloud (exponential form).
Definition: mptrac.h:2975
int met_zstd_level
ZSTD compression level (from -5 to 22).
Definition: mptrac.h:2598
int csi_ny
Number of latitudes of gridded CSI data.
Definition: mptrac.h:3087
int vtk_sphere
Spherical projection for VTK data (0=no, 1=yes).
Definition: mptrac.h:3240
double chemgrid_z0
Lower altitude of chemistry grid [km].
Definition: mptrac.h:2927
double met_pbl_min
Minimum depth of planetary boundary layer [km].
Definition: mptrac.h:2710
int qnt_iwc
Quantity array index for cloud ice water content.
Definition: mptrac.h:2370
double chemgrid_lat0
Lower latitude of chemistry grid [deg].
Definition: mptrac.h:2945
double conv_cape
CAPE threshold for convection module [J/kg].
Definition: mptrac.h:2795
int qnt_Co1d
Quantity array index for O(1D) volume mixing ratio (chemistry code).
Definition: mptrac.h:2523
double met_cms_eps_pv
cmultiscale compression epsilon for potential vorticity.
Definition: mptrac.h:2632
int qnt_pw
Quantity array index for partial water vapor pressure.
Definition: mptrac.h:2448
double grid_z1
Upper altitude of gridded data [km].
Definition: mptrac.h:3135
int direction
Direction flag (1=forward calculation, -1=backward calculation).
Definition: mptrac.h:2553
int qnt_Cccl4
Quantity array index for CFC-10 volume mixing ratio (chemistry code).
Definition: mptrac.h:2529
int met_dp
Stride for pressure levels.
Definition: mptrac.h:2662
double met_dt_out
Time step for sampling of meteo data along trajectories [s].
Definition: mptrac.h:2729
int qnt_h2o2
Quantity array index for H2O2 volume mixing ratio (climatology).
Definition: mptrac.h:2412
int qnt_vh
Quantity array index for horizontal wind.
Definition: mptrac.h:2478
int csi_nx
Number of longitudes of gridded CSI data.
Definition: mptrac.h:3078
double csi_lat0
Lower latitude of gridded CSI data [deg].
Definition: mptrac.h:3090
double turb_dz_trop
Vertical turbulent diffusion coefficient (troposphere) [m^2/s].
Definition: mptrac.h:2777
int met_pbl
Planetary boundary layer data (0=file, 1=z2p, 2=Richardson, 3=theta).
Definition: mptrac.h:2707
int qnt_lwc
Quantity array index for cloud liquid water content.
Definition: mptrac.h:2364
double turb_mesoz
Vertical scaling factor for mesoscale wind fluctuations.
Definition: mptrac.h:2786
int grid_nc_level
zlib compression level of netCDF grid data files (0=off).
Definition: mptrac.h:3120
int grid_nx
Number of longitudes of gridded data.
Definition: mptrac.h:3138
int atm_type
Type of atmospheric data files (0=ASCII, 1=binary, 2=netCDF, 3=CLaMS_traj, 4=CLaMS_pos).
Definition: mptrac.h:3034
double bound_mass
Boundary conditions mass per particle [kg].
Definition: mptrac.h:2804
double grid_lat0
Lower latitude of gridded data [deg].
Definition: mptrac.h:3150
int qnt_ts
Quantity array index for surface temperature.
Definition: mptrac.h:2295
int qnt_loss_rate
Quantity array index for total loss rate.
Definition: mptrac.h:2439
double met_cms_eps_h2o
cmultiscale compression epsilon for water vapor.
Definition: mptrac.h:2635
int qnt_plfc
Quantity array index for pressure at level of free convection (LCF).
Definition: mptrac.h:2391
double grid_lon0
Lower longitude of gridded data [deg].
Definition: mptrac.h:3141
int qnt_o1d
Quantity array index for O(1D) volume mixing ratio (climatology).
Definition: mptrac.h:2418
int met_tropo_spline
Tropopause interpolation method (0=linear, 1=spline).
Definition: mptrac.h:2726
int qnt_tvirt
Quantity array index for virtual temperature.
Definition: mptrac.h:2472
double dt_met
Time step of meteo data [s].
Definition: mptrac.h:2572
double chemgrid_lat1
Upper latitude of chemistry grid [deg].
Definition: mptrac.h:2948
int met_geopot_sy
Latitudinal smoothing of geopotential heights.
Definition: mptrac.h:2698
double met_cms_eps_u
cmultiscale compression epsilon for zonal wind.
Definition: mptrac.h:2623
double turb_dx_strat
Horizontal turbulent diffusion coefficient (stratosphere) [m^2/s].
Definition: mptrac.h:2771
int qnt_vmr
Quantity array index for volume mixing ratio.
Definition: mptrac.h:2283
int qnt_lsm
Quantity array index for land-sea mask.
Definition: mptrac.h:2316
int qnt_theta
Quantity array index for potential temperature.
Definition: mptrac.h:2460
double bound_lat1
Boundary conditions maximum longitude [deg].
Definition: mptrac.h:2819
double stat_t1
Stop time for station output [s].
Definition: mptrac.h:3222
double turb_dx_trop
Horizontal turbulent diffusion coefficient (troposphere) [m^2/s].
Definition: mptrac.h:2768
int grid_type
Type of grid data files (0=ASCII, 1=netCDF).
Definition: mptrac.h:3156
double csi_lon0
Lower longitude of gridded CSI data [deg].
Definition: mptrac.h:3081
int qnt_pbl
Quantity array index for boundary layer pressure.
Definition: mptrac.h:2322
int grid_stddev
Include standard deviations in grid output (0=no, 1=yes).
Definition: mptrac.h:3126
int qnt_psice
Quantity array index for saturation pressure over ice.
Definition: mptrac.h:2445
double chemgrid_lon0
Lower longitude of chemistry grid [deg].
Definition: mptrac.h:2936
int bound_pbl
Boundary conditions planetary boundary layer (0=no, 1=yes).
Definition: mptrac.h:2837
int qnt_mloss_wet
Quantity array index for total mass loss due to wet deposition.
Definition: mptrac.h:2430
int met_geopot_sx
Longitudinal smoothing of geopotential heights.
Definition: mptrac.h:2695
int met_sy
Smoothing for latitudes.
Definition: mptrac.h:2668
int qnt_ps
Quantity array index for surface pressure.
Definition: mptrac.h:2292
int rng_type
Random number generator (0=GSL, 1=Squares, 2=cuRAND).
Definition: mptrac.h:2759
int isosurf
Isosurface parameter (0=none, 1=pressure, 2=density, 3=theta, 4=balloon).
Definition: mptrac.h:2746
double bound_p1
Boundary conditions top pressure [hPa].
Definition: mptrac.h:2825
int qnt_zs
Quantity array index for surface geopotential height.
Definition: mptrac.h:2298
int prof_nz
Number of altitudes of gridded profile data.
Definition: mptrac.h:3165
double csi_dt_out
Time step for CSI output [s].
Definition: mptrac.h:3057
int met_cape
Convective available potential energy data (0=file, 1=calculate).
Definition: mptrac.h:2704
double csi_modmin
Minimum column density to trigger detection [kg/m^2].
Definition: mptrac.h:3066
int met_sx
Smoothing for longitudes.
Definition: mptrac.h:2665
double chemgrid_lon1
Upper longitude of chemistry grid [deg].
Definition: mptrac.h:2939
double turb_mesox
Horizontal scaling factor for mesoscale wind fluctuations.
Definition: mptrac.h:2783
double met_cms_eps_iwc
cmultiscale compression epsilon for cloud ice water content.
Definition: mptrac.h:2647
double met_cms_eps_swc
cmultiscale compression epsilon for cloud snow water content.
Definition: mptrac.h:2650
double met_cms_eps_v
cmultiscale compression epsilon for meridional wind.
Definition: mptrac.h:2626
double prof_z0
Lower altitude of gridded profile data [km].
Definition: mptrac.h:3168
int qnt_w
Quantity array index for vertical velocity.
Definition: mptrac.h:2355
double bound_vmr
Boundary conditions volume mixing ratio [ppv].
Definition: mptrac.h:2810
double met_tropo_pv
Dynamical tropopause potential vorticity threshold [PVU].
Definition: mptrac.h:2720
int prof_nx
Number of longitudes of gridded profile data.
Definition: mptrac.h:3174
int qnt_stat
Quantity array index for station flag.
Definition: mptrac.h:2277
int met_tropo
Tropopause definition (0=none, 1=clim, 2=cold point, 3=WMO_1st, 4=WMO_2nd, 5=dynamical).
Definition: mptrac.h:2717
int qnt_rp
Quantity array index for particle radius.
Definition: mptrac.h:2286
int met_mpi_share
Use MPI to share meteo (0=no, 1=yes).
Definition: mptrac.h:2735
double mixing_strat
Interparcel exchange parameter for mixing in the stratosphere.
Definition: mptrac.h:2894
int qnt_vz
Quantity array index for vertical velocity.
Definition: mptrac.h:2481
int qnt_ho2
Quantity array index for HO2 volume mixing ratio (climatology).
Definition: mptrac.h:2415
double csi_z1
Upper altitude of gridded CSI data [km].
Definition: mptrac.h:3075
double stat_t0
Start time for station output [s].
Definition: mptrac.h:3219
double oh_chem_beta
Beta parameter for diurnal variablity of OH.
Definition: mptrac.h:2957
double wet_depo_so2_ph
pH value used to calculate effective Henry constant of SO2.
Definition: mptrac.h:2993
double mixing_z0
Lower altitude of mixing grid [km].
Definition: mptrac.h:2900
int qnt_mloss_decay
Quantity array index for total mass loss due to exponential decay.
Definition: mptrac.h:2436
int atm_type_out
Type of atmospheric data files for output (-1=same as ATM_TYPE, 0=ASCII, 1=binary,...
Definition: mptrac.h:3039
int met_nlev
Number of meteo data model levels.
Definition: mptrac.h:2686
double dt_kpp
Time step for KPP chemistry [s].
Definition: mptrac.h:2966
double dry_depo_dp
Dry deposition surface layer [hPa].
Definition: mptrac.h:3002
int qnt_shf
Quantity array index for surface sensible heat flux.
Definition: mptrac.h:2313
int qnt_vs
Quantity array index for surface meridional wind.
Definition: mptrac.h:2304
int qnt_Cco
Quantity array index for CO volume mixing ratio (chemistry code).
Definition: mptrac.h:2508
double vtk_dt_out
Time step for VTK data output [s].
Definition: mptrac.h:3228
double t_stop
Stop time of simulation [s].
Definition: mptrac.h:2559
double conv_dt
Time interval for convection module [s].
Definition: mptrac.h:2801
int qnt_hno3
Quantity array index for HNO3 volume mixing ratio (climatology).
Definition: mptrac.h:2406
int met_clams
Read MPTRAC or CLaMS meteo data (0=MPTRAC, 1=CLaMS).
Definition: mptrac.h:2586
int qnt_h2ot
Quantity array index for tropopause water vapor volume mixing ratio.
Definition: mptrac.h:2334
int qnt_rh
Quantity array index for relative humidity over water.
Definition: mptrac.h:2454
double met_cms_eps_cc
cmultiscale compression epsilon for cloud cover.
Definition: mptrac.h:2653
double bound_lat0
Boundary conditions minimum longitude [deg].
Definition: mptrac.h:2816
double met_pbl_max
Maximum depth of planetary boundary layer [km].
Definition: mptrac.h:2713
int met_dx
Stride for longitudes.
Definition: mptrac.h:2656
int qnt_destination
Quantity array index for destination subdomain in domain decomposition.
Definition: mptrac.h:2550
int mixing_ny
Number of latitudes of mixing grid.
Definition: mptrac.h:2915
int met_convention
Meteo data layout (0=[lev, lat, lon], 1=[lon, lat, lev]).
Definition: mptrac.h:2575
int qnt_zeta_d
Quantity array index for diagnosed zeta vertical coordinate.
Definition: mptrac.h:2466
int tracer_chem
Switch for first order tracer chemistry module (0=off, 1=on).
Definition: mptrac.h:2969
double dt_mod
Time step of simulation [s].
Definition: mptrac.h:2562
int diffusion
Diffusion scheme (0=off, 1=fixed-K, 2=PBL).
Definition: mptrac.h:2762
int qnt_tnat
Quantity array index for T_NAT.
Definition: mptrac.h:2496
int qnt_tice
Quantity array index for T_ice.
Definition: mptrac.h:2490
int qnt_zg
Quantity array index for geopotential height.
Definition: mptrac.h:2337
double vtk_offset
Vertical offset for VTK data [km].
Definition: mptrac.h:3237
int qnt_v
Quantity array index for meridional wind.
Definition: mptrac.h:2352
int qnt_mloss_dry
Quantity array index for total mass loss due to dry deposition.
Definition: mptrac.h:2433
double bound_vmr_trend
Boundary conditions volume mixing ratio trend [ppv/s].
Definition: mptrac.h:2813
int met_cache
Preload meteo data into disk cache (0=no, 1=yes).
Definition: mptrac.h:2732
int qnt_oh
Quantity array index for OH volume mixing ratio (climatology).
Definition: mptrac.h:2409
int qnt_Ch
Quantity array index for H volume mixing ratio (chemistry code).
Definition: mptrac.h:2514
int met_press_level_def
Use predefined pressure levels or not.
Definition: mptrac.h:2683
int oh_chem_reaction
Reaction type for OH chemistry (0=none, 2=bimolecular, 3=termolecular).
Definition: mptrac.h:2951
int qnt_h2o
Quantity array index for water vapor volume mixing ratio.
Definition: mptrac.h:2358
int prof_ny
Number of latitudes of gridded profile data.
Definition: mptrac.h:3183
int qnt_rhice
Quantity array index for relative humidity over ice.
Definition: mptrac.h:2457
int qnt_rho
Quantity array index for density of air.
Definition: mptrac.h:2346
double sample_dz
Layer depth for sample output [km].
Definition: mptrac.h:3204
double tdec_strat
Life time of particles in the stratosphere [s].
Definition: mptrac.h:2849
int obs_type
Type of observation data files (0=ASCII, 1=netCDF).
Definition: mptrac.h:3048
double met_cms_eps_lwc
cmultiscale compression epsilon for cloud liquid water content.
Definition: mptrac.h:2641
int qnt_us
Quantity array index for surface zonal wind.
Definition: mptrac.h:2301
double met_cms_eps_z
cmultiscale compression epsilon for geopotential height.
Definition: mptrac.h:2617
double grid_lon1
Upper longitude of gridded data [deg].
Definition: mptrac.h:3144
int qnt_Cn2o
Quantity array index for N2O volume mixing ratio (chemistry code).
Definition: mptrac.h:2538
int qnt_Cccl3f
Quantity array index for CFC-11 volume mixing ratio (chemistry code).
Definition: mptrac.h:2532
double mixing_lat0
Lower latitude of mixing grid [deg].
Definition: mptrac.h:2918
int nens
Number of ensembles.
Definition: mptrac.h:3096
int qnt_pt
Quantity array index for tropopause pressure.
Definition: mptrac.h:2325
int qnt_cl
Quantity array index for total column cloud water.
Definition: mptrac.h:2385
int advect
Advection scheme (0=off, 1=Euler, 2=midpoint, 4=Runge-Kutta).
Definition: mptrac.h:2752
double prof_z1
Upper altitude of gridded profile data [km].
Definition: mptrac.h:3171
int qnt_t
Quantity array index for temperature.
Definition: mptrac.h:2343
int atm_filter
Time filter for atmospheric data output (0=none, 1=missval, 2=remove).
Definition: mptrac.h:3027
int kpp_chem
Switch for KPP chemistry module (0=off, 1=on).
Definition: mptrac.h:2963
int qnt_zeta
Quantity array index for zeta vertical coordinate.
Definition: mptrac.h:2463
double conv_pbl_trans
Depth of PBL transition layer (fraction of PBL depth).
Definition: mptrac.h:2792
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:2579
double csi_z0
Lower altitude of gridded CSI data [km].
Definition: mptrac.h:3072
int qnt_lapse
Quantity array index for lapse rate.
Definition: mptrac.h:2475
double stat_lat
Latitude of station [deg].
Definition: mptrac.h:3213
int qnt_Cho2
Quantity array index for HO2 volume mixing ratio (chemistry code).
Definition: mptrac.h:2517
int grid_ny
Number of latitudes of gridded data.
Definition: mptrac.h:3147
int qnt_Csf6
Quantity array index for SF6 volume mixing ratio (chemistry code).
Definition: mptrac.h:2541
int qnt_Ch2o
Quantity array index for H2O volume mixing ratio (chemistry code).
Definition: mptrac.h:2502
double met_detrend
FWHM of horizontal Gaussian used for detrending [km].
Definition: mptrac.h:2674
int conv_mix_pbl
Vertical mixing in the PBL (0=off, 1=on).
Definition: mptrac.h:2789
double bound_dps
Boundary conditions surface layer depth [hPa].
Definition: mptrac.h:2828
int dd_nbr_neighbours
Number of neighbours to communicate with.
Definition: mptrac.h:3249
double met_cms_eps_t
cmultiscale compression epsilon for temperature.
Definition: mptrac.h:2620
int chemgrid_nz
Number of altitudes of chemistry grid.
Definition: mptrac.h:2924
int qnt_cape
Quantity array index for convective available potential energy (CAPE).
Definition: mptrac.h:2397
int qnt_zeta_dot
Quantity array index forvelocity of zeta vertical coordinate.
Definition: mptrac.h:2469
double bound_mass_trend
Boundary conditions mass per particle trend [kg/s].
Definition: mptrac.h:2807
int mixing_nz
Number of altitudes of mixing grid.
Definition: mptrac.h:2897
int qnt_o3c
Quantity array index for total column ozone.
Definition: mptrac.h:2403
double bound_p0
Boundary conditions bottom pressure [hPa].
Definition: mptrac.h:2822
double mixing_lon0
Lower longitude of mixing grid [deg].
Definition: mptrac.h:2909
int qnt_Co3
Quantity array index for O3 volume mixing ratio (chemistry code).
Definition: mptrac.h:2505
int qnt_tsts
Quantity array index for T_STS.
Definition: mptrac.h:2493
int grid_nz
Number of altitudes of gridded data.
Definition: mptrac.h:3129
int qnt_nss
Quantity array index for northward turbulent surface stress.
Definition: mptrac.h:2310
double ens_dt_out
Time step for ensemble output [s].
Definition: mptrac.h:3102
int atm_stride
Particle index stride for atmospheric data files.
Definition: mptrac.h:3030
int met_relhum
Try to read relative humidity (0=no, 1=yes).
Definition: mptrac.h:2701
double mixing_lat1
Upper latitude of mixing grid [deg].
Definition: mptrac.h:2921
double atm_dt_out
Time step for atmospheric data output [s].
Definition: mptrac.h:3024
double prof_lat1
Upper latitude of gridded profile data [deg].
Definition: mptrac.h:3189
int met_cms_batch
cmultiscale batch size.
Definition: mptrac.h:2607
double psc_h2o
H2O volume mixing ratio for PSC analysis.
Definition: mptrac.h:3008
int met_sp
Smoothing for pressure levels.
Definition: mptrac.h:2671
double prof_lon0
Lower longitude of gridded profile data [deg].
Definition: mptrac.h:3177
int chemgrid_nx
Number of longitudes of chemistry grid.
Definition: mptrac.h:2933
int qnt_pct
Quantity array index for cloud top pressure.
Definition: mptrac.h:2379
int qnt_mloss_kpp
Quantity array index for total mass loss due to KPP chemistry.
Definition: mptrac.h:2427
int qnt_psat
Quantity array index for saturation pressure over water.
Definition: mptrac.h:2442
int qnt_subdomain
Quantity array index for current subdomain in domain decomposition.
Definition: mptrac.h:2547
double prof_lat0
Lower latitude of gridded profile data [deg].
Definition: mptrac.h:3186
int qnt_cin
Quantity array index for convective inhibition (CIN).
Definition: mptrac.h:2400
double psc_hno3
HNO3 volume mixing ratio for PSC analysis.
Definition: mptrac.h:3011
double prof_lon1
Upper longitude of gridded profile data [deg].
Definition: mptrac.h:3180
double met_cms_eps_rwc
cmultiscale compression epsilon for cloud rain water content.
Definition: mptrac.h:2644
int met_nc_quant
Number of digits for quantization of netCDF meteo files (0=off).
Definition: mptrac.h:2595
int h2o2_chem_reaction
Reaction type for H2O2 chemistry (0=none, 1=SO2).
Definition: mptrac.h:2960
int qnt_Co3p
Quantity array index for O(3P) volume mixing ratio (chemistry code).
Definition: mptrac.h:2526
double wet_depo_bc_ret_ratio
Coefficients for wet deposition below cloud: retention ratio.
Definition: mptrac.h:2999
int chemgrid_ny
Number of latitudes of chemistry grid.
Definition: mptrac.h:2942
double met_cms_eps_o3
cmultiscale compression epsilon for ozone.
Definition: mptrac.h:2638
int met_cms_zstd
cmultiscale ZSTD compression (0=off, 1=on).
Definition: mptrac.h:2610
int grid_sparse
Sparse output in grid data files (0=no, 1=yes).
Definition: mptrac.h:3117
double dry_depo_vdep
Dry deposition velocity [m/s].
Definition: mptrac.h:3005
int qnt_tt
Quantity array index for tropopause temperature.
Definition: mptrac.h:2328
int met_np
Number of target pressure levels.
Definition: mptrac.h:2677
int qnt_ens
Quantity array index for ensemble IDs.
Definition: mptrac.h:2274
int met_nc_level
zlib compression level of netCDF meteo files (0=off).
Definition: mptrac.h:2592
double mixing_dt
Time interval for mixing [s].
Definition: mptrac.h:2888
int qnt_mloss_h2o2
Quantity array index for total mass loss due to H2O2 chemistry.
Definition: mptrac.h:2424
double vtk_scale
Vertical scaling factor for VTK data.
Definition: mptrac.h:3234
double met_cms_eps_w
cmultiscale compression epsilon for vertical velocity.
Definition: mptrac.h:2629
double turb_dx_pbl
Horizontal turbulent diffusion coefficient (PBL) [m^2/s].
Definition: mptrac.h:2765
double conv_cin
CIN threshold for convection module [J/kg].
Definition: mptrac.h:2798
int qnt_pv
Quantity array index for potential vorticity.
Definition: mptrac.h:2484
int advect_vert_coord
Vertical velocity of air parcels (0=omega_on_plev, 1=zetadot_on_mlev, 2=omega_on_mlev).
Definition: mptrac.h:2756
int qnt_mloss_oh
Quantity array index for total mass loss due to OH chemistry.
Definition: mptrac.h:2421
int qnt_Ch2o2
Quantity array index for H2O2 volume mixing ratio (chemistry code).
Definition: mptrac.h:2520
int qnt_sst
Quantity array index for sea surface temperature.
Definition: mptrac.h:2319
double mixing_lon1
Upper longitude of mixing grid [deg].
Definition: mptrac.h:2912
int atm_nc_level
zlib compression level of netCDF atmospheric data files (0=off).
Definition: mptrac.h:3042
int met_cms_heur
cmultiscale coarsening heuristics (0=default, 1=mean diff, 2=median diff, 3=max diff).
Definition: mptrac.h:2614
double wet_depo_ic_ret_ratio
Coefficients for wet deposition in cloud: retention ratio.
Definition: mptrac.h:2996
int qnt_sh
Quantity array index for specific humidity.
Definition: mptrac.h:2451
int qnt_ess
Quantity array index for eastward turbulent surface stress.
Definition: mptrac.h:2307
double wet_depo_ic_b
Coefficient B for wet deposition in cloud (exponential form).
Definition: mptrac.h:2984
double wet_depo_bc_b
Coefficient B for wet deposition below cloud (exponential form).
Definition: mptrac.h:2978
int met_dy
Stride for latitudes.
Definition: mptrac.h:2659
int qnt_Cx
Quantity array index for trace species x volume mixing ratio (chemistry code).
Definition: mptrac.h:2499
double turb_dz_strat
Vertical turbulent diffusion coefficient (stratosphere) [m^2/s].
Definition: mptrac.h:2780
double bound_zetas
Boundary conditions surface layer zeta [K].
Definition: mptrac.h:2834
int dd_subdomains_zonal
Zonal subdomain number.
Definition: mptrac.h:3243
int qnt_idx
Quantity array index for air parcel IDs.
Definition: mptrac.h:2271
double met_tropo_theta
Dynamical tropopause potential temperature threshold [K].
Definition: mptrac.h:2723
int qnt_rwc
Quantity array index for cloud rain water content.
Definition: mptrac.h:2367
double t_start
Start time of simulation [s].
Definition: mptrac.h:2556
int nq
Number of quantities.
Definition: mptrac.h:2256
double tdec_trop
Life time of particles in the troposphere [s].
Definition: mptrac.h:2846
double sample_dx
Horizontal radius for sample output [km].
Definition: mptrac.h:3201
int vtk_stride
Particle index stride for VTK data.
Definition: mptrac.h:3231
double turb_dz_pbl
Vertical turbulent diffusion coefficient (PBL) [m^2/s].
Definition: mptrac.h:2774
double grid_lat1
Upper latitude of gridded data [deg].
Definition: mptrac.h:3153
int dd_subdomains_meridional
Meridional subdomain number.
Definition: mptrac.h:3246
int qnt_zt
Quantity array index for tropopause geopotential height.
Definition: mptrac.h:2331
int met_type
Type of meteo data files (0=netCDF, 1=binary, 2=pck, 3=ZFP, 4=ZSTD, 5=cms, 6=grib,...
Definition: mptrac.h:2583
int qnt_cc
Quantity array index for cloud cover.
Definition: mptrac.h:2376
int qnt_plcl
Quantity array index for pressure at lifted condensation level (LCL).
Definition: mptrac.h:2388
double grid_dt_out
Time step for gridded data output [s].
Definition: mptrac.h:3114
int qnt_tdew
Quantity array index for dew point temperature.
Definition: mptrac.h:2487
Meteo data structure.
Definition: mptrac.h:3542
double subdomain_lon_max
Rectangular grid limit of subdomain.
Definition: mptrac.h:3719
double subdomain_lat_min
Rectangular grid limit of subdomain.
Definition: mptrac.h:3728
int nx
Number of longitudes.
Definition: mptrac.h:3548
double subdomain_lon_min
Rectangular grid limit of subdomain.
Definition: mptrac.h:3722
int ny_glob
Global sizes of meteo data.
Definition: mptrac.h:3752
int ny
Number of latitudes.
Definition: mptrac.h:3551
int halo_offset_start
Definition: mptrac.h:3743
int np
Number of pressure levels.
Definition: mptrac.h:3554
int halo_offset_end
Definition: mptrac.h:3746
int nx_glob
Global sizes of meteo data.
Definition: mptrac.h:3749
int npl
Number of model levels.
Definition: mptrac.h:3557
double time
Time [s].
Definition: mptrac.h:3545
double subdomain_lat_max
Rectangular grid limit of subdomain.
Definition: mptrac.h:3725
int np_glob
Global sizes of meteo data.
Definition: mptrac.h:3755
MPI information data.
Definition: mptrac.h:3319
int rank
Rank of node.
Definition: mptrac.h:3322
int size
Size of node.
Definition: mptrac.h:3325
Particle data.
Definition: mptrac.h:3293
double p
Pressure [hPa].
Definition: mptrac.h:3299
double lat
Latitude [deg].
Definition: mptrac.h:3305
double time
Time [s].
Definition: mptrac.h:3296
double lon
Longitude [deg].
Definition: mptrac.h:3302