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_spline.h>
128#include <gsl/gsl_statistics.h>
129#include <math.h>
130#include <netcdf.h>
131#include <omp.h>
132#include <stdint.h>
133#include <stdio.h>
134#include <stdlib.h>
135#include <string.h>
136#include <time.h>
137#include <sys/time.h>
138
139#ifdef MPI
140#include "mpi.h"
141#endif
142
143#ifdef DD
144#include <netcdf_par.h>
145#include <stdbool.h>
146#endif
147
148#ifdef _OPENACC
149#include "openacc.h"
150#endif
151
152#ifdef CURAND
153#include "curand.h"
154#endif
155
156#ifdef ZFP
157#include "zfp.h"
158#endif
159
160#ifdef ZSTD
161#include "zstd.h"
162#endif
163
164#ifdef CMS
165#include "cmultiscale.h"
166#endif
167
168#ifdef KPP
169#include "chem_Parameters.h"
170#include "chem_Global.h"
171#include "chem_Sparse.h"
172#endif
173
174#ifdef ECCODES
175#include "eccodes.h"
176#endif
177
178/* ------------------------------------------------------------
179 Constants...
180 ------------------------------------------------------------ */
181
183#ifndef AVO
184#define AVO 6.02214076e23
185#endif
186
188#ifndef CPD
189#define CPD 1003.5
190#endif
191
193#ifndef EPS
194#define EPS (MH2O / MA)
195#endif
196
198#ifndef G0
199#define G0 9.80665
200#endif
201
203#ifndef H0
204#define H0 7.0
205#endif
206
208#ifndef LV
209#define LV 2501000.
210#endif
211
213#ifndef KARMAN
214#define KARMAN 0.40
215#endif
216
218#ifndef KB
219#define KB 1.3806504e-23
220#endif
221
223#ifndef MA
224#define MA 28.9644
225#endif
226
228#ifndef MH2O
229#define MH2O 18.01528
230#endif
231
233#ifndef MO3
234#define MO3 48.00
235#endif
236
238#ifndef P0
239#define P0 1013.25
240#endif
241
243#ifndef RA
244#define RA (1e3 * RI / MA)
245#endif
246
248#ifndef RE
249#define RE 6367.421
250#endif
251
253#ifndef RI
254#define RI 8.3144598
255#endif
256
258#ifndef T0
259#define T0 273.15
260#endif
261
263#ifndef DD_NPOLE
264#define DD_NPOLE -2
265#endif
266
267
269#ifndef DD_SPOLE
270#define DD_SPOLE -3
271#endif
272
273/* ------------------------------------------------------------
274 Dimensions...
275 ------------------------------------------------------------ */
276
278#ifndef EP
279#define EP 140
280#endif
281
283#ifndef EX
284#define EX 1444
285#endif
286
288#ifndef EY
289#define EY 724
290#endif
291
293#ifndef EP_GLOB
294#define EP_GLOB 150
295#endif
296
298#ifndef EX_GLOB
299#define EX_GLOB 1202
300#endif
301
303#ifndef EY_GLOB
304#define EY_GLOB 602
305#endif
306
308#ifndef LEN
309#define LEN 5000
310#endif
311
313#ifndef METVAR
314#define METVAR 13
315#endif
316
318#ifndef NP
319#define NP 10000000
320#endif
321
323#ifndef NQ
324#define NQ 15
325#endif
326
328#ifndef NCSI
329#define NCSI 1000000
330#endif
331
333#ifndef NENS
334#define NENS 2000
335#endif
336
338#ifndef NOBS
339#define NOBS 10000000
340#endif
341
343#ifndef NTHREADS
344#define NTHREADS 512
345#endif
346
348#ifndef CY
349#define CY 250
350#endif
351
353#ifndef CO3
354#define CO3 30
355#endif
356
358#ifndef CP
359#define CP 70
360#endif
361
363#ifndef CSZA
364#define CSZA 50
365#endif
366
368#ifndef CT
369#define CT 12
370#endif
371
373#ifndef CTS
374#define CTS 1000
375#endif
376
378#ifndef DD_NPART
379#define DD_NPART 100000
380#endif
381
383#ifndef DD_NNMAX
384#define DD_NNMAX 26
385#endif
386
387/* ------------------------------------------------------------
388 Macros...
389 ------------------------------------------------------------ */
390
410#ifdef _OPENACC
411#define ALLOC(ptr, type, n) \
412 if(acc_get_num_devices(acc_device_nvidia) <= 0) \
413 ERRMSG("Not running on a GPU device!"); \
414 if((ptr=calloc((size_t)(n), sizeof(type)))==NULL) \
415 ERRMSG("Out of memory!");
416#else
417#define ALLOC(ptr, type, n) \
418 if((ptr=calloc((size_t)(n), sizeof(type)))==NULL) \
419 ERRMSG("Out of memory!");
420#endif
421
440#define ARRAY_2D(ix, iy, ny) \
441 ((ix) * (ny) + (iy))
442
459#define ARRAY_3D(ix, iy, ny, iz, nz) \
460 (((ix)*(ny) + (iy)) * (nz) + (iz))
461
484#define ARRHENIUS(a, b, t) \
485 ((a) * exp( -(b) / (t)))
486
508#define DEG2DX(dlon, lat) \
509 (RE * DEG2RAD(dlon) * cos(DEG2RAD(lat)))
510
529#define DEG2DY(dlat) \
530 (RE * DEG2RAD(dlat))
531
546#define DEG2RAD(deg) \
547 ((deg) * (M_PI / 180.0))
548
571#define DP2DZ(dp, p) \
572 (- (dp) * H0 / (p))
573
593#define DX2DEG(dx, lat) \
594 (((lat) < -89.999 || (lat) > 89.999) ? 0 \
595 : (dx) * 180. / (M_PI * RE * cos(DEG2RAD(lat))))
596
611#define DY2DEG(dy) \
612 ((dy) * 180. / (M_PI * RE))
613
630#define DZ2DP(dz, p) \
631 (-(dz) * (p) / H0)
632
646#define DIST(a, b) \
647 sqrt(DIST2(a, b))
648
662#define DIST2(a, b) \
663 ((a[0]-b[0])*(a[0]-b[0])+(a[1]-b[1])*(a[1]-b[1])+(a[2]-b[2])*(a[2]-b[2]))
664
678#define DOTP(a, b) \
679 (a[0]*b[0]+a[1]*b[1]+a[2]*b[2])
680
692#define ECC(cmd) { \
693 int ecc_result=(cmd); \
694 if(ecc_result!=0) \
695 ERRMSG("ECCODES error: %s", codes_get_error_message(ecc_result)); \
696 }
697
711#define ECC_READ_2D(variable,target,scaling_factor,found_flag){ \
712 if( strcmp(short_name,variable)==0){ \
713 for (int ix = 0; ix < met->nx; ix++) { \
714 for (int iy = 0; iy < met->ny; iy++) { \
715 target[ix][iy] = (float)(values[iy * met->nx + ix]*scaling_factor); \
716 } \
717 } \
718 found_flag =1; \
719 } \
720 }
721
735#define ECC_READ_3D(variable,level,target,scaling_factor,found_flag){ \
736 if( strcmp(short_name,variable)==0){ \
737 for (int ix = 0; ix < met->nx; ix++) { \
738 for (int iy = 0; iy < met->ny; iy++) { \
739 target[ix][iy][level] = (float) (values[iy * met->nx + ix]*scaling_factor); \
740 } \
741 } \
742 found_flag +=1; \
743 } \
744 }
745
757#define ECC(cmd) { \
758 int ecc_result=(cmd); \
759 if(ecc_result!=0) \
760 ERRMSG("ECCODES error: %s", codes_get_error_message(ecc_result)); \
761 }
762
776#define ECC_READ_2D(variable,target,scaling_factor,found_flag){ \
777 if( strcmp(short_name,variable)==0){ \
778 for (int ix = 0; ix < met->nx; ix++) { \
779 for (int iy = 0; iy < met->ny; iy++) { \
780 target[ix][iy] = (float)(values[iy * met->nx + ix]*scaling_factor); \
781 } \
782 } \
783 found_flag =1; \
784 } \
785 }
786
800#define ECC_READ_3D(variable,level,target,scaling_factor,found_flag){ \
801 if( strcmp(short_name,variable)==0){ \
802 for (int ix = 0; ix < met->nx; ix++) { \
803 for (int iy = 0; iy < met->ny; iy++) { \
804 target[ix][iy][level] = (float) (values[iy * met->nx + ix]*scaling_factor); \
805 } \
806 } \
807 found_flag +=1; \
808 } \
809 }
810
827#define FMOD(x, y) \
828 ((x) - (int) ((x) / (y)) * (y))
829
845#define FREAD(ptr, type, size, in) { \
846 if(fread(ptr, sizeof(type), size, in)!=size) \
847 ERRMSG("Error while reading!"); \
848 }
849
865#define FWRITE(ptr, type, size, out) { \
866 if(fwrite(ptr, sizeof(type), size, out)!=size) \
867 ERRMSG("Error while writing!"); \
868 }
869
880#define INTPOL_INIT \
881 double cw[4] = {0.0, 0.0, 0.0, 0.0}; int ci[3] = {0, 0, 0};
882
894#define INTPOL_2D(var, init) \
895 intpol_met_time_2d(met0, met0->var, met1, met1->var, \
896 atm->time[ip], atm->lon[ip], atm->lat[ip], \
897 &var, ci, cw, init);
898
911#define INTPOL_3D(var, init) \
912 intpol_met_time_3d(met0, met0->var, met1, met1->var, \
913 atm->time[ip], atm->p[ip], \
914 atm->lon[ip], atm->lat[ip], \
915 &var, ci, cw, init);
916
930#define INTPOL_SPACE_ALL(p, lon, lat) { \
931 intpol_met_space_3d(met, met->z, p, lon, lat, &z, ci, cw, 1); \
932 intpol_met_space_3d(met, met->t, p, lon, lat, &t, ci, cw, 0); \
933 intpol_met_space_3d(met, met->u, p, lon, lat, &u, ci, cw, 0); \
934 intpol_met_space_3d(met, met->v, p, lon, lat, &v, ci, cw, 0); \
935 intpol_met_space_3d(met, met->w, p, lon, lat, &w, ci, cw, 0); \
936 intpol_met_space_3d(met, met->pv, p, lon, lat, &pv, ci, cw, 0); \
937 intpol_met_space_3d(met, met->h2o, p, lon, lat, &h2o, ci, cw, 0); \
938 intpol_met_space_3d(met, met->o3, p, lon, lat, &o3, ci, cw, 0); \
939 intpol_met_space_3d(met, met->lwc, p, lon, lat, &lwc, ci, cw, 0); \
940 intpol_met_space_3d(met, met->rwc, p, lon, lat, &rwc, ci, cw, 0); \
941 intpol_met_space_3d(met, met->iwc, p, lon, lat, &iwc, ci, cw, 0); \
942 intpol_met_space_3d(met, met->swc, p, lon, lat, &swc, ci, cw, 0); \
943 intpol_met_space_3d(met, met->cc, p, lon, lat, &cc, ci, cw, 0); \
944 intpol_met_space_2d(met, met->ps, lon, lat, &ps, ci, cw, 0); \
945 intpol_met_space_2d(met, met->ts, lon, lat, &ts, ci, cw, 0); \
946 intpol_met_space_2d(met, met->zs, lon, lat, &zs, ci, cw, 0); \
947 intpol_met_space_2d(met, met->us, lon, lat, &us, ci, cw, 0); \
948 intpol_met_space_2d(met, met->vs, lon, lat, &vs, ci, cw, 0); \
949 intpol_met_space_2d(met, met->ess, ess, lat, &ess, ci, cw, 0); \
950 intpol_met_space_2d(met, met->nss, nss, lat, &nss, ci, cw, 0); \
951 intpol_met_space_2d(met, met->shf, shf, lat, &shf, ci, cw, 0); \
952 intpol_met_space_2d(met, met->lsm, lon, lat, &lsm, ci, cw, 0); \
953 intpol_met_space_2d(met, met->sst, lon, lat, &sst, ci, cw, 0); \
954 intpol_met_space_2d(met, met->pbl, lon, lat, &pbl, ci, cw, 0); \
955 intpol_met_space_2d(met, met->pt, lon, lat, &pt, ci, cw, 0); \
956 intpol_met_space_2d(met, met->tt, lon, lat, &tt, ci, cw, 0); \
957 intpol_met_space_2d(met, met->zt, lon, lat, &zt, ci, cw, 0); \
958 intpol_met_space_2d(met, met->h2ot, lon, lat, &h2ot, ci, cw, 0); \
959 intpol_met_space_2d(met, met->pct, lon, lat, &pct, ci, cw, 0); \
960 intpol_met_space_2d(met, met->pcb, lon, lat, &pcb, ci, cw, 0); \
961 intpol_met_space_2d(met, met->cl, lon, lat, &cl, ci, cw, 0); \
962 intpol_met_space_2d(met, met->plcl, lon, lat, &plcl, ci, cw, 0); \
963 intpol_met_space_2d(met, met->plfc, lon, lat, &plfc, ci, cw, 0); \
964 intpol_met_space_2d(met, met->pel, lon, lat, &pel, ci, cw, 0); \
965 intpol_met_space_2d(met, met->cape, lon, lat, &cape, ci, cw, 0); \
966 intpol_met_space_2d(met, met->cin, lon, lat, &cin, ci, cw, 0); \
967 intpol_met_space_2d(met, met->o3c, lon, lat, &o3c, ci, cw, 0); \
968 }
969
984#define INTPOL_TIME_ALL(time, p, lon, lat) { \
985 intpol_met_time_3d(met0, met0->z, met1, met1->z, time, p, lon, lat, &z, ci, cw, 1); \
986 intpol_met_time_3d(met0, met0->t, met1, met1->t, time, p, lon, lat, &t, ci, cw, 0); \
987 intpol_met_time_3d(met0, met0->u, met1, met1->u, time, p, lon, lat, &u, ci, cw, 0); \
988 intpol_met_time_3d(met0, met0->v, met1, met1->v, time, p, lon, lat, &v, ci, cw, 0); \
989 intpol_met_time_3d(met0, met0->w, met1, met1->w, time, p, lon, lat, &w, ci, cw, 0); \
990 intpol_met_time_3d(met0, met0->pv, met1, met1->pv, time, p, lon, lat, &pv, ci, cw, 0); \
991 intpol_met_time_3d(met0, met0->h2o, met1, met1->h2o, time, p, lon, lat, &h2o, ci, cw, 0); \
992 intpol_met_time_3d(met0, met0->o3, met1, met1->o3, time, p, lon, lat, &o3, ci, cw, 0); \
993 intpol_met_time_3d(met0, met0->lwc, met1, met1->lwc, time, p, lon, lat, &lwc, ci, cw, 0); \
994 intpol_met_time_3d(met0, met0->rwc, met1, met1->rwc, time, p, lon, lat, &rwc, ci, cw, 0); \
995 intpol_met_time_3d(met0, met0->iwc, met1, met1->iwc, time, p, lon, lat, &iwc, ci, cw, 0); \
996 intpol_met_time_3d(met0, met0->swc, met1, met1->swc, time, p, lon, lat, &swc, ci, cw, 0); \
997 intpol_met_time_3d(met0, met0->cc, met1, met1->cc, time, p, lon, lat, &cc, ci, cw, 0); \
998 intpol_met_time_2d(met0, met0->ps, met1, met1->ps, time, lon, lat, &ps, ci, cw, 0); \
999 intpol_met_time_2d(met0, met0->ts, met1, met1->ts, time, lon, lat, &ts, ci, cw, 0); \
1000 intpol_met_time_2d(met0, met0->zs, met1, met1->zs, time, lon, lat, &zs, ci, cw, 0); \
1001 intpol_met_time_2d(met0, met0->us, met1, met1->us, time, lon, lat, &us, ci, cw, 0); \
1002 intpol_met_time_2d(met0, met0->vs, met1, met1->vs, time, lon, lat, &vs, ci, cw, 0); \
1003 intpol_met_time_2d(met0, met0->ess, met1, met1->ess, time, lon, lat, &ess, ci, cw, 0); \
1004 intpol_met_time_2d(met0, met0->nss, met1, met1->nss, time, lon, lat, &nss, ci, cw, 0); \
1005 intpol_met_time_2d(met0, met0->shf, met1, met1->shf, time, lon, lat, &shf, ci, cw, 0); \
1006 intpol_met_time_2d(met0, met0->lsm, met1, met1->lsm, time, lon, lat, &lsm, ci, cw, 0); \
1007 intpol_met_time_2d(met0, met0->sst, met1, met1->sst, time, lon, lat, &sst, ci, cw, 0); \
1008 intpol_met_time_2d(met0, met0->pbl, met1, met1->pbl, time, lon, lat, &pbl, ci, cw, 0); \
1009 intpol_met_time_2d(met0, met0->pt, met1, met1->pt, time, lon, lat, &pt, ci, cw, 0); \
1010 intpol_met_time_2d(met0, met0->tt, met1, met1->tt, time, lon, lat, &tt, ci, cw, 0); \
1011 intpol_met_time_2d(met0, met0->zt, met1, met1->zt, time, lon, lat, &zt, ci, cw, 0); \
1012 intpol_met_time_2d(met0, met0->h2ot, met1, met1->h2ot, time, lon, lat, &h2ot, ci, cw, 0); \
1013 intpol_met_time_2d(met0, met0->pct, met1, met1->pct, time, lon, lat, &pct, ci, cw, 0); \
1014 intpol_met_time_2d(met0, met0->pcb, met1, met1->pcb, time, lon, lat, &pcb, ci, cw, 0); \
1015 intpol_met_time_2d(met0, met0->cl, met1, met1->cl, time, lon, lat, &cl, ci, cw, 0); \
1016 intpol_met_time_2d(met0, met0->plcl, met1, met1->plcl, time, lon, lat, &plcl, ci, cw, 0); \
1017 intpol_met_time_2d(met0, met0->plfc, met1, met1->plfc, time, lon, lat, &plfc, ci, cw, 0); \
1018 intpol_met_time_2d(met0, met0->pel, met1, met1->pel, time, lon, lat, &pel, ci, cw, 0); \
1019 intpol_met_time_2d(met0, met0->cape, met1, met1->cape, time, lon, lat, &cape, ci, cw, 0); \
1020 intpol_met_time_2d(met0, met0->cin, met1, met1->cin, time, lon, lat, &cin, ci, cw, 0); \
1021 intpol_met_time_2d(met0, met0->o3c, met1, met1->o3c, time, lon, lat, &o3c, ci, cw, 0); \
1022 }
1023
1038#define LAPSE(p1, t1, p2, t2) \
1039 (1e3 * G0 / RA * ((t2) - (t1)) / ((t2) + (t1)) \
1040 * ((p2) + (p1)) / ((p2) - (p1)))
1041
1057#define LIN(x0, y0, x1, y1, x) \
1058 ((y0)+((y1)-(y0))/((x1)-(x0))*((x)-(x0)))
1059
1084#define MAX(a,b) \
1085 (((a)>(b))?(a):(b))
1086
1098#define MET_HEADER \
1099 fprintf(out, \
1100 "# $1 = time [s]\n" \
1101 "# $2 = altitude [km]\n" \
1102 "# $3 = longitude [deg]\n" \
1103 "# $4 = latitude [deg]\n" \
1104 "# $5 = pressure [hPa]\n" \
1105 "# $6 = temperature [K]\n" \
1106 "# $7 = zonal wind [m/s]\n" \
1107 "# $8 = meridional wind [m/s]\n" \
1108 "# $9 = vertical velocity [hPa/s]\n" \
1109 "# $10 = H2O volume mixing ratio [ppv]\n"); \
1110 fprintf(out, \
1111 "# $11 = O3 volume mixing ratio [ppv]\n" \
1112 "# $12 = geopotential height [km]\n" \
1113 "# $13 = potential vorticity [PVU]\n" \
1114 "# $14 = surface pressure [hPa]\n" \
1115 "# $15 = surface temperature [K]\n" \
1116 "# $16 = surface geopotential height [km]\n" \
1117 "# $17 = surface zonal wind [m/s]\n" \
1118 "# $18 = surface meridional wind [m/s]\n" \
1119 "# $19 = eastward turbulent surface stress [N/m^2]\n" \
1120 "# $20 = northward turbulent surface stress [N/m^2]\n"); \
1121 fprintf(out, \
1122 "# $21 = surface sensible heat flux [W/m^2]\n" \
1123 "# $22 = land-sea mask [1]\n" \
1124 "# $23 = sea surface temperature [K]\n" \
1125 "# $24 = tropopause pressure [hPa]\n" \
1126 "# $25 = tropopause geopotential height [km]\n" \
1127 "# $26 = tropopause temperature [K]\n" \
1128 "# $27 = tropopause water vapor [ppv]\n" \
1129 "# $28 = cloud liquid water content [kg/kg]\n" \
1130 "# $29 = cloud rain water content [kg/kg]\n" \
1131 "# $30 = cloud ice water content [kg/kg]\n"); \
1132 fprintf(out, \
1133 "# $31 = cloud snow water content [kg/kg]\n" \
1134 "# $32 = cloud cover [1]\n" \
1135 "# $33 = total column cloud water [kg/m^2]\n" \
1136 "# $34 = cloud top pressure [hPa]\n" \
1137 "# $35 = cloud bottom pressure [hPa]\n" \
1138 "# $36 = pressure at lifted condensation level (LCL) [hPa]\n" \
1139 "# $37 = pressure at level of free convection (LFC) [hPa]\n" \
1140 "# $38 = pressure at equilibrium level (EL) [hPa]\n" \
1141 "# $39 = convective available potential energy (CAPE) [J/kg]\n" \
1142 "# $40 = convective inhibition (CIN) [J/kg]\n"); \
1143 fprintf(out, \
1144 "# $41 = relative humidity over water [%%]\n" \
1145 "# $42 = relative humidity over ice [%%]\n" \
1146 "# $43 = dew point temperature [K]\n" \
1147 "# $44 = frost point temperature [K]\n" \
1148 "# $45 = NAT temperature [K]\n" \
1149 "# $46 = HNO3 volume mixing ratio [ppv]\n" \
1150 "# $47 = OH volume mixing ratio [ppv]\n" \
1151 "# $48 = H2O2 volume mixing ratio [ppv]\n" \
1152 "# $49 = HO2 volume mixing ratio [ppv]\n" \
1153 "# $50 = O(1D) volume mixing ratio [ppv]\n"); \
1154 fprintf(out, \
1155 "# $51 = boundary layer pressure [hPa]\n" \
1156 "# $52 = total column ozone [DU]\n" \
1157 "# $53 = number of data points\n" \
1158 "# $54 = number of tropopause data points\n" \
1159 "# $55 = number of CAPE data points\n");
1160
1185#define MIN(a,b) \
1186 (((a)<(b))?(a):(b))
1187
1200#define MOLEC_DENS(p,t) \
1201 (AVO * 1e-6 * ((p) * 100) / (RI * (t)))
1202
1214#define NC(cmd) { \
1215 int nc_result=(cmd); \
1216 if(nc_result!=NC_NOERR) \
1217 ERRMSG("%s", nc_strerror(nc_result)); \
1218 }
1219
1243#define NC_DEF_VAR(varname, type, ndims, dims, long_name, units, level, quant) { \
1244 NC(nc_def_var(ncid, varname, type, ndims, dims, &varid)); \
1245 NC(nc_put_att_text(ncid, varid, "long_name", strnlen(long_name, LEN), long_name)); \
1246 NC(nc_put_att_text(ncid, varid, "units", strnlen(units, LEN), units)); \
1247 if((quant) > 0) \
1248 NC(nc_def_var_quantize(ncid, varid, NC_QUANTIZE_GRANULARBR, quant)); \
1249 if((level) != 0) { \
1250 NC(nc_def_var_deflate(ncid, varid, 1, 1, level)); \
1251 /* unsigned int ulevel = (unsigned int)level; */ \
1252 /* NC(nc_def_var_filter(ncid, varid, 32015, 1, (unsigned int[]){ulevel})); */ \
1253 } \
1254 }
1255
1273#define NC_GET_DOUBLE(varname, ptr, force) { \
1274 if(force) { \
1275 NC(nc_inq_varid(ncid, varname, &varid)); \
1276 NC(nc_get_var_double(ncid, varid, ptr)); \
1277 } else { \
1278 if(nc_inq_varid(ncid, varname, &varid) == NC_NOERR) { \
1279 NC(nc_get_var_double(ncid, varid, ptr)); \
1280 } else \
1281 WARN("netCDF variable %s is missing!", varname); \
1282 } \
1283 }
1284
1301#define NC_INQ_DIM(dimname, ptr, min, max) { \
1302 int dimid; size_t naux; \
1303 NC(nc_inq_dimid(ncid, dimname, &dimid)); \
1304 NC(nc_inq_dimlen(ncid, dimid, &naux)); \
1305 *ptr = (int)naux; \
1306 if ((*ptr) < (min) || (*ptr) > (max)) \
1307 ERRMSG("Dimension %s is out of range!", dimname); \
1308 }
1309
1324#define NC_PUT_DOUBLE(varname, ptr, hyperslab) { \
1325 NC(nc_inq_varid(ncid, varname, &varid)); \
1326 if(hyperslab) { \
1327 NC(nc_put_vara_double(ncid, varid, start, count, ptr)); \
1328 } else { \
1329 NC(nc_put_var_double(ncid, varid, ptr)); \
1330 } \
1331 }
1332
1348#define NC_PUT_FLOAT(varname, ptr, hyperslab) { \
1349 NC(nc_inq_varid(ncid, varname, &varid)); \
1350 if(hyperslab) { \
1351 NC(nc_put_vara_float(ncid, varid, start, count, ptr)); \
1352 } else { \
1353 NC(nc_put_var_float(ncid, varid, ptr)); \
1354 } \
1355 }
1356
1371#define NC_PUT_INT(varname, ptr, hyperslab) { \
1372 NC(nc_inq_varid(ncid, varname, &varid)); \
1373 if(hyperslab) { \
1374 NC(nc_put_vara_int(ncid, varid, start, count, ptr)); \
1375 } else { \
1376 NC(nc_put_var_int(ncid, varid, ptr)); \
1377 } \
1378 }
1379
1393#define NC_PUT_ATT(varname, attname, text) { \
1394 NC(nc_inq_varid(ncid, varname, &varid)); \
1395 NC(nc_put_att_text(ncid, varid, attname, strnlen(text, LEN), text)); \
1396 }
1397
1410#define NC_PUT_ATT_GLOBAL(attname, text) \
1411 NC(nc_put_att_text(ncid, NC_GLOBAL, attname, strnlen(text, LEN), text));
1412
1430#define NN(x0, y0, x1, y1, x) \
1431 (fabs((x) - (x0)) <= fabs((x) - (x1)) ? (y0) : (y1))
1432
1445#define NORM(a) \
1446 sqrt(DOTP(a, a))
1447
1463#ifdef _OPENACC
1464#define PARTICLE_LOOP(ip0, ip1, check_dt, ...) \
1465 const int ip0_const = ip0; \
1466 const int ip1_const = ip1; \
1467 _Pragma(__VA_ARGS__) \
1468 _Pragma("acc parallel loop independent gang vector") \
1469 for (int ip = ip0_const; ip < ip1_const; ip++) \
1470 if (!check_dt || cache->dt[ip] != 0)
1471#else
1472#define PARTICLE_LOOP(ip0, ip1, check_dt, ...) \
1473 const int ip0_const = ip0; \
1474 const int ip1_const = ip1; \
1475 _Pragma("omp parallel for default(shared)") \
1476 for (int ip = ip0_const; ip < ip1_const; ip++) \
1477 if (!check_dt || cache->dt[ip] != 0)
1478#endif
1479
1502#define P(z) \
1503 (P0 * exp(-(z) / H0))
1504
1526#define PSAT(t) \
1527 (6.112 * exp(17.62 * ((t) - T0) / (243.12 + (t) - T0)))
1528
1550#define PSICE(t) \
1551 (6.112 * exp(22.46 * ((t) - T0) / (272.62 + (t) - T0)))
1552
1577#define PW(p, h2o) \
1578 ((p) * MAX((h2o), 0.1e-6) / (1. + (1. - EPS) * MAX((h2o), 0.1e-6)))
1579
1594#define RAD2DEG(rad) \
1595 ((rad) * (180.0 / M_PI))
1596
1624#define RH(p, t, h2o) \
1625 (PW(p, h2o) / PSAT(t) * 100.)
1626
1654#define RHICE(p, t, h2o) \
1655 (PW(p, h2o) / PSICE(t) * 100.)
1656
1679#define RHO(p, t) \
1680 (100. * (p) / (RA * (t)))
1681
1698#define SET_ATM(qnt, val) \
1699 if (ctl->qnt >= 0) \
1700 atm->q[ctl->qnt][ip] = val;
1701
1721#define SET_QNT(qnt, name, longname, unit) \
1722 if (strcasecmp(ctl->qnt_name[iq], name) == 0) { \
1723 ctl->qnt = iq; \
1724 sprintf(ctl->qnt_longname[iq], longname); \
1725 sprintf(ctl->qnt_unit[iq], unit); \
1726 } else
1727
1742#define SH(h2o) \
1743 (EPS * MAX((h2o), 0.1e-6))
1744
1755#define SQR(x) \
1756 ((x)*(x))
1757
1769#define SWAP(x, y, type) \
1770 do {type tmp = x; x = y; y = tmp;} while(0);
1771
1793#define TDEW(p, h2o) \
1794 (T0 + 243.12 * log(PW((p), (h2o)) / 6.112) \
1795 / (17.62 - log(PW((p), (h2o)) / 6.112)))
1796
1818#define TICE(p, h2o) \
1819 (T0 + 272.62 * log(PW((p), (h2o)) / 6.112) \
1820 / (22.46 - log(PW((p), (h2o)) / 6.112)))
1821
1842#define THETA(p, t) \
1843 ((t) * pow(1000. / (p), 0.286))
1844
1871#define THETAVIRT(p, t, h2o) \
1872 (TVIRT(THETA((p), (t)), MAX((h2o), 0.1e-6)))
1873
1892#define TOK(line, tok, format, var) { \
1893 if(((tok)=strtok((line), " \t"))) { \
1894 if(sscanf(tok, format, &(var))!=1) continue; \
1895 } else ERRMSG("Error while reading!"); \
1896 }
1897
1917#define TVIRT(t, h2o) \
1918 ((t) * (1. + (1. - EPS) * MAX((h2o), 0.1e-6)))
1919
1939#define Z(p) \
1940 (H0 * log(P0 / (p)))
1941
1970#define ZDIFF(lnp0, t0, h2o0, lnp1, t1, h2o1) \
1971 (RI / MA / G0 * 0.5 * (TVIRT((t0), (h2o0)) + TVIRT((t1), (h2o1))) \
1972 * ((lnp0) - (lnp1)))
1973
1989#define ZETA(ps, p, t) \
1990 (((p) / (ps) <= 0.3 ? 1. : \
1991 sin(M_PI / 2. * (1. - (p) / (ps)) / (1. - 0.3))) \
1992 * THETA((p), (t)))
1993
1994/* ------------------------------------------------------------
1995 Log messages...
1996 ------------------------------------------------------------ */
1997
1999#ifndef LOGLEV
2000#define LOGLEV 2
2001#endif
2002
2032#define LOG(level, ...) { \
2033 if(level >= 2) \
2034 printf(" "); \
2035 if(level <= LOGLEV) { \
2036 printf(__VA_ARGS__); \
2037 printf("\n"); \
2038 } \
2039 }
2040
2069#define WARN(...) { \
2070 printf("\nWarning (%s, %s, l%d): ", __FILE__, __func__, __LINE__); \
2071 LOG(0, __VA_ARGS__); \
2072 }
2073
2102#define ERRMSG(...) { \
2103 printf("\nError (%s, %s, l%d): ", __FILE__, __func__, __LINE__); \
2104 LOG(0, __VA_ARGS__); \
2105 exit(EXIT_FAILURE); \
2106 }
2107
2137#define PRINT(format, var) \
2138 printf("Print (%s, %s, l%d): %s= "format"\n", \
2139 __FILE__, __func__, __LINE__, #var, var);
2140
2141/* ------------------------------------------------------------
2142 Timers...
2143 ------------------------------------------------------------ */
2144
2146#define NTIMER 100
2147
2161#define PRINT_TIMERS \
2162 timer("END", "END", 1);
2163
2182#define SELECT_TIMER(id, group, color) { \
2183 NVTX_POP; \
2184 NVTX_PUSH(id, color); \
2185 timer(id, group, 0); \
2186 }
2187
2201#define START_TIMERS \
2202 NVTX_PUSH("START", NVTX_CPU);
2203
2216#define STOP_TIMERS \
2217 NVTX_POP;
2218
2219/* ------------------------------------------------------------
2220 NVIDIA Tools Extension (NVTX)...
2221 ------------------------------------------------------------ */
2222
2223#ifdef NVTX
2224#include "nvToolsExt.h"
2225
2227#define NVTX_CPU 0xFFADD8E6
2228
2230#define NVTX_GPU 0xFF00008B
2231
2233#define NVTX_H2D 0xFFFFFF00
2234
2236#define NVTX_D2H 0xFFFF8800
2237
2239#define NVTX_READ 0xFFFFCCCB
2240
2242#define NVTX_WRITE 0xFF8B0000
2243
2245#define NVTX_RECV 0xFFCCFFCB
2246
2248#define NVTX_SEND 0xFF008B00
2249
2279#define NVTX_PUSH(range_title, range_color) { \
2280 nvtxEventAttributes_t eventAttrib = {0}; \
2281 eventAttrib.version = NVTX_VERSION; \
2282 eventAttrib.size = NVTX_EVENT_ATTRIB_STRUCT_SIZE; \
2283 eventAttrib.messageType = NVTX_MESSAGE_TYPE_ASCII; \
2284 eventAttrib.colorType = NVTX_COLOR_ARGB; \
2285 eventAttrib.color = range_color; \
2286 eventAttrib.message.ascii = range_title; \
2287 nvtxRangePushEx(&eventAttrib); \
2288 }
2289
2302#define NVTX_POP { \
2303 nvtxRangePop(); \
2304 }
2305#else
2306
2307/* Empty definitions of NVTX_PUSH and NVTX_POP... */
2308#define NVTX_PUSH(range_title, range_color) {}
2309#define NVTX_POP {}
2310#endif
2311
2312/* ------------------------------------------------------------
2313 Thrust...
2314 ------------------------------------------------------------ */
2315
2341 double *__restrict__ c,
2342 int n,
2343 int *__restrict__ index);
2344
2345/* ------------------------------------------------------------
2346 Structs...
2347 ------------------------------------------------------------ */
2348
2356typedef struct {
2357
2358 /* ------------------------------------------------------------
2359 Quantity parameters...
2360 ------------------------------------------------------------ */
2361
2363 int nq;
2364
2366 char qnt_name[NQ][LEN];
2367
2369 char qnt_longname[NQ][LEN];
2370
2372 char qnt_unit[NQ][LEN];
2373
2375 char qnt_format[NQ][LEN];
2376
2379
2382
2385
2388
2391
2394
2397
2400
2403
2406
2409
2412
2415
2418
2421
2424
2427
2430
2433
2436
2439
2442
2445
2448
2451
2454
2457
2460
2463
2466
2469
2472
2475
2478
2481
2484
2487
2490
2493
2496
2499
2502
2505
2508
2511
2514
2517
2520
2523
2526
2529
2532
2535
2538
2541
2544
2547
2550
2553
2556
2559
2562
2565
2568
2571
2574
2577
2580
2583
2586
2589
2592
2595
2598
2601
2604
2607
2610
2613
2616
2619
2622
2625
2628
2631
2634
2637
2640
2643
2646
2649
2652
2655
2658
2661
2663 double t_start;
2664
2666 double t_stop;
2667
2669 double dt_mod;
2670
2671 /* ------------------------------------------------------------
2672 Meteo data parameters...
2673 ------------------------------------------------------------ */
2674
2676 char metbase[LEN];
2677
2679 double dt_met;
2680
2683
2687
2691
2694
2697
2700
2703
2706
2708 int met_zfp_prec[METVAR];
2709
2711 double met_zfp_tol[METVAR];
2712
2715
2718
2722
2725
2728
2731
2734
2737
2740
2743
2746
2749
2752
2755
2758
2761
2764
2767
2770
2773
2776
2779
2782
2785
2787 double met_p[EP];
2788
2791
2794
2796 double met_lev_hyam[EP];
2797
2799 double met_lev_hybm[EP];
2800
2803
2806
2809
2812
2815
2818
2821
2825
2828
2831
2834
2837
2840
2843
2844 /* ------------------------------------------------------------
2845 Geophysical module parameters...
2846 ------------------------------------------------------------ */
2847
2849 double sort_dt;
2850
2854
2856 char balloon[LEN];
2857
2860
2864
2867
2870
2873
2876
2879
2882
2885
2888
2891
2894
2897
2900
2903
2905 double conv_cin;
2906
2908 double conv_dt;
2909
2912
2915
2918
2921
2924
2927
2929 double bound_p0;
2930
2932 double bound_p1;
2933
2936
2939
2942
2945
2947 char species[LEN];
2948
2950 double molmass;
2951
2954
2957
2960
2962 char clim_hno3_filename[LEN];
2963
2965 char clim_oh_filename[LEN];
2966
2968 char clim_h2o2_filename[LEN];
2969
2971 char clim_ho2_filename[LEN];
2972
2974 char clim_o1d_filename[LEN];
2975
2977 char clim_o3_filename[LEN];
2978
2980 char clim_ccl4_timeseries[LEN];
2981
2983 char clim_ccl3f_timeseries[LEN];
2984
2986 char clim_ccl2f2_timeseries[LEN];
2987
2989 char clim_n2o_timeseries[LEN];
2990
2992 char clim_sf6_timeseries[LEN];
2993
2996
2999
3002
3005
3008
3011
3014
3017
3020
3023
3026
3029
3032
3035
3038
3041
3044
3047
3050
3053
3056
3059
3061 double oh_chem[4];
3062
3065
3068
3071
3073 double dt_kpp;
3074
3077
3079 double wet_depo_pre[2];
3080
3083
3086
3089
3092
3094 double wet_depo_ic_h[2];
3095
3097 double wet_depo_bc_h[2];
3098
3101
3104
3107
3110
3113
3115 double psc_h2o;
3116
3118 double psc_hno3;
3119
3120 /* ------------------------------------------------------------
3121 Output parameters...
3122 ------------------------------------------------------------ */
3123
3125 char atm_basename[LEN];
3126
3128 char atm_gpfile[LEN];
3129
3132
3135
3138
3142
3147
3150
3152 int atm_nc_quant[NQ];
3153
3156
3158 char csi_basename[LEN];
3159
3161 char csi_kernel[LEN];
3162
3165
3167 char csi_obsfile[LEN];
3168
3171
3174
3177
3179 double csi_z0;
3180
3182 double csi_z1;
3183
3186
3188 double csi_lon0;
3189
3191 double csi_lon1;
3192
3195
3197 double csi_lat0;
3198
3200 double csi_lat1;
3201
3203 int nens;
3204
3206 char ens_basename[LEN];
3207
3210
3212 char grid_basename[LEN];
3213
3215 char grid_kernel[LEN];
3216
3218 char grid_gpfile[LEN];
3219
3222
3225
3228
3230 int grid_nc_quant[NQ];
3231
3234
3237
3239 double grid_z0;
3240
3242 double grid_z1;
3243
3246
3249
3252
3255
3258
3261
3264
3266 char prof_basename[LEN];
3267
3269 char prof_obsfile[LEN];
3270
3273
3275 double prof_z0;
3276
3278 double prof_z1;
3279
3282
3285
3288
3291
3294
3297
3299 char sample_basename[LEN];
3300
3302 char sample_kernel[LEN];
3303
3305 char sample_obsfile[LEN];
3306
3309
3312
3314 char stat_basename[LEN];
3315
3317 double stat_lon;
3318
3320 double stat_lat;
3321
3323 double stat_r;
3324
3326 double stat_t0;
3327
3329 double stat_t1;
3330
3332 char vtk_basename[LEN];
3333
3336
3339
3342
3345
3348
3351
3354
3357
3360
3361} ctl_t;
3362
3371typedef struct {
3372
3374 int np;
3375
3377 double time[NP];
3378
3380 double p[NP];
3381
3383 double lon[NP];
3384
3386 double lat[NP];
3387
3389 double q[NQ][NP];
3390
3391} atm_t;
3392
3400typedef struct {
3401
3403 double time;
3404
3406 double p;
3407
3409 double lon;
3410
3412 double lat;
3413
3415 double q[NQ];
3416
3417} particle_t;
3418
3426typedef struct {
3428 int rank;
3429
3431 int size;
3432
3433#ifdef DD
3435 int neighbours[DD_NNMAX];
3436
3438 MPI_Datatype MPI_Particle;
3439#endif
3440
3441} mpi_info_t;
3442
3449typedef struct {
3450
3452 double iso_var[NP];
3453
3455 double iso_ps[NP];
3456
3458 double iso_ts[NP];
3459
3462
3464 float uvwp[NP][3];
3465
3467 double rs[3 * NP + 1];
3468
3470 double dt[NP];
3471
3472} cache_t;
3473
3481typedef struct {
3482
3484 int np;
3485
3487 int nsza;
3488
3490 int no3c;
3491
3493 double p[CP];
3494
3496 double sza[CSZA];
3497
3499 double o3c[CO3];
3500
3502 double n2o[CP][CSZA][CO3];
3503
3505 double ccl4[CP][CSZA][CO3];
3506
3508 double ccl3f[CP][CSZA][CO3];
3509
3511 double ccl2f2[CP][CSZA][CO3];
3512
3514 double o2[CP][CSZA][CO3];
3515
3517 double o3_1[CP][CSZA][CO3];
3518
3520 double o3_2[CP][CSZA][CO3];
3521
3523 double h2o2[CP][CSZA][CO3];
3524
3526 double h2o[CP][CSZA][CO3];
3527
3528} clim_photo_t;
3529
3537typedef struct {
3538
3541
3543 double time[CTS];
3544
3546 double vmr[CTS];
3547
3548} clim_ts_t;
3549
3557typedef struct {
3558
3561
3563 int nlat;
3564
3566 int np;
3567
3569 double time[CT];
3570
3572 double lat[CY];
3573
3575 double p[CP];
3576
3578 double vmr[CT][CP][CY];
3579
3580} clim_zm_t;
3581
3589typedef struct {
3590
3593
3596
3598 double tropo_time[12];
3599
3601 double tropo_lat[73];
3602
3604 double tropo[12][73];
3605
3608
3611
3614
3617
3620
3623
3626
3629
3632
3635
3638
3639} clim_t;
3640
3648typedef struct {
3649
3651 double time;
3652
3654 int nx;
3655
3657 int ny;
3658
3660 int np;
3661
3663 int npl;
3664
3665 // TODO:
3666 // They need global sizes now, maybe in the future just keep EX, EY, EP and
3667 // Introduce help data structure in read_met_grid etc.
3668
3670#ifdef DD
3671 double lon[EX_GLOB];
3672#else
3673 double lon[EX];
3674#endif
3675
3677#ifdef DD
3678 double lat[EY_GLOB];
3679#else
3680 double lat[EY];
3681#endif
3682
3684#ifdef DD
3685 double p[EP_GLOB];
3686#else
3687 double p[EP];
3688#endif
3689
3691 double hybrid[EP];
3692
3694 float ps[EX][EY];
3695
3697 float ts[EX][EY];
3698
3700 float zs[EX][EY];
3701
3703 float us[EX][EY];
3704
3706 float vs[EX][EY];
3707
3709 float ess[EX][EY];
3710
3712 float nss[EX][EY];
3713
3715 float shf[EX][EY];
3716
3718 float lsm[EX][EY];
3719
3721 float sst[EX][EY];
3722
3724 float pbl[EX][EY];
3725
3727 float pt[EX][EY];
3728
3730 float tt[EX][EY];
3731
3733 float zt[EX][EY];
3734
3736 float h2ot[EX][EY];
3737
3739 float pct[EX][EY];
3740
3742 float pcb[EX][EY];
3743
3745 float cl[EX][EY];
3746
3748 float plcl[EX][EY];
3749
3751 float plfc[EX][EY];
3752
3754 float pel[EX][EY];
3755
3757 float cape[EX][EY];
3758
3760 float cin[EX][EY];
3761
3763 float o3c[EX][EY];
3764
3766 float z[EX][EY][EP];
3767
3769 float t[EX][EY][EP];
3770
3772 float u[EX][EY][EP];
3773
3775 float v[EX][EY][EP];
3776
3778 float w[EX][EY][EP];
3779
3781 float pv[EX][EY][EP];
3782
3784 float h2o[EX][EY][EP];
3785
3787 float o3[EX][EY][EP];
3788
3790 float lwc[EX][EY][EP];
3791
3793 float rwc[EX][EY][EP];
3794
3796 float iwc[EX][EY][EP];
3797
3799 float swc[EX][EY][EP];
3800
3802 float cc[EX][EY][EP];
3803
3805 float pl[EX][EY][EP];
3806
3808 float ul[EX][EY][EP];
3809
3811 float vl[EX][EY][EP];
3812
3814 float wl[EX][EY][EP];
3815
3817 float zetal[EX][EY][EP];
3818
3820 float zeta_dotl[EX][EY][EP];
3821
3822 // TODO: Integrate this into a grid_t ?
3823
3826
3829
3832
3835
3837 size_t subdomain_start[4];
3838
3840 size_t subdomain_count[4];
3841
3842 /* Hyperslab of boundary halos start. */
3843 size_t halo_bnd_start[4];
3844
3845 /* Hyperslab of boundary halos count. */
3846 size_t halo_bnd_count[4];
3847
3848 /* Hyperslab of boundary halos count. */
3850
3851 /* Hyperslab of boundary halos count. */
3853
3856
3859
3862
3863} met_t;
3864
3865/* ------------------------------------------------------------
3866 OpenACC routines...
3867 ------------------------------------------------------------ */
3868
3869#ifdef _OPENACC
3870#pragma acc routine (clim_oh)
3871#pragma acc routine (clim_photo)
3872#pragma acc routine (clim_tropo)
3873#pragma acc routine (clim_ts)
3874#pragma acc routine (clim_zm)
3875#pragma acc routine (intpol_check_lon_lat)
3876#pragma acc routine (intpol_met_4d_coord)
3877#pragma acc routine (intpol_met_space_3d)
3878#pragma acc routine (intpol_met_space_3d_ml)
3879#pragma acc routine (intpol_met_space_2d)
3880#pragma acc routine (intpol_met_time_3d)
3881#pragma acc routine (intpol_met_time_3d_ml)
3882#pragma acc routine (intpol_met_time_2d)
3883#pragma acc routine (kernel_weight)
3884#pragma acc routine (lapse_rate)
3885#pragma acc routine (locate_irr)
3886#pragma acc routine (locate_irr_float)
3887#pragma acc routine (locate_reg)
3888#pragma acc routine (locate_vert)
3889#pragma acc routine (nat_temperature)
3890#pragma acc routine (pbl_weight)
3891#pragma acc routine (sedi)
3892#pragma acc routine (stddev)
3893#pragma acc routine (sza_calc)
3894#pragma acc routine (tropo_weight)
3895#endif
3896
3897/* ------------------------------------------------------------
3898 Functions...
3899 ------------------------------------------------------------ */
3900
3923 void *data,
3924 size_t N);
3925
3940void cart2geo(
3941 const double *x,
3942 double *z,
3943 double *lon,
3944 double *lat);
3945
3968double clim_oh(
3969 const ctl_t * ctl,
3970 const clim_t * clim,
3971 const double t,
3972 const double lon,
3973 const double lat,
3974 const double p);
3975
3995 const ctl_t * ctl,
3996 clim_t * clim);
3997
4027double clim_photo(
4028 const double rate[CP][CSZA][CO3],
4029 const clim_photo_t * photo,
4030 const double p,
4031 const double sza,
4032 const double o3c);
4033
4059double clim_tropo(
4060 const clim_t * clim,
4061 const double t,
4062 const double lat);
4063
4082void clim_tropo_init(
4083 clim_t * clim);
4084
4101double clim_ts(
4102 const clim_ts_t * ts,
4103 const double t);
4104
4126double clim_zm(
4127 const clim_zm_t * zm,
4128 const double t,
4129 const double lat,
4130 const double p);
4131
4173 const ctl_t * ctl,
4174 const char *varname,
4175 float *array,
4176 const size_t nx,
4177 const size_t ny,
4178 const size_t np,
4179 const int decompress,
4180 FILE * inout);
4181
4214void compress_pck(
4215 const char *varname,
4216 float *array,
4217 const size_t nxy,
4218 const size_t nz,
4219 const int decompress,
4220 FILE * inout);
4221
4261 const char *varname,
4262 float *array,
4263 const int nx,
4264 const int ny,
4265 const int nz,
4266 const int precision,
4267 const double tolerance,
4268 const int decompress,
4269 FILE * inout);
4270
4293 const char *varname,
4294 float *array,
4295 const size_t n,
4296 const int decompress,
4297 const int level,
4298 FILE * inout);
4299
4322void day2doy(
4323 const int year,
4324 const int mon,
4325 const int day,
4326 int *doy);
4327
4349void doy2day(
4350 const int year,
4351 const int doy,
4352 int *mon,
4353 int *day);
4354
4381void fft_help(
4382 double *fcReal,
4383 double *fcImag,
4384 const int n);
4385
4412void geo2cart(
4413 const double z,
4414 const double lon,
4415 const double lat,
4416 double *x);
4417
4442void get_met_help(
4443 const ctl_t * ctl,
4444 const double t,
4445 const int direct,
4446 const char *metbase,
4447 const double dt_met,
4448 char *filename);
4449
4473void get_met_replace(
4474 char *orig,
4475 char *search,
4476 char *repl);
4477
4514void get_tropo(
4515 const int met_tropo,
4516 ctl_t * ctl,
4517 clim_t * clim,
4518 met_t * met,
4519 const double *lons,
4520 const int nx,
4521 const double *lats,
4522 const int ny,
4523 double *pt,
4524 double *zt,
4525 double *tt,
4526 double *qt,
4527 double *o3t,
4528 double *ps,
4529 double *zs);
4530
4552 const double *lons,
4553 const int nlon,
4554 const double *lats,
4555 const int nlat,
4556 const double lon,
4557 const double lat,
4558 double *lon2,
4559 double *lat2);
4560
4603 const met_t * met0,
4604 float height0[EX][EY][EP],
4605 float array0[EX][EY][EP],
4606 const met_t * met1,
4607 float height1[EX][EY][EP],
4608 float array1[EX][EY][EP],
4609 const double ts,
4610 const double height,
4611 const double lon,
4612 const double lat,
4613 double *var,
4614 int *ci,
4615 double *cw,
4616 const int init);
4617
4653 const met_t * met,
4654 float array[EX][EY][EP],
4655 const double p,
4656 const double lon,
4657 const double lat,
4658 double *var,
4659 int *ci,
4660 double *cw,
4661 const int init);
4662
4682 const met_t * met,
4683 float zs[EX][EY][EP],
4684 float vals[EX][EY][EP],
4685 const double z,
4686 const double lon,
4687 const double lat,
4688 double *val);
4689
4725 const met_t * met,
4726 float array[EX][EY],
4727 const double lon,
4728 const double lat,
4729 double *var,
4730 int *ci,
4731 double *cw,
4732 const int init);
4733
4768 const met_t * met0,
4769 float array0[EX][EY][EP],
4770 const met_t * met1,
4771 float array1[EX][EY][EP],
4772 const double ts,
4773 const double p,
4774 const double lon,
4775 const double lat,
4776 double *var,
4777 int *ci,
4778 double *cw,
4779 const int init);
4780
4804 const met_t * met0,
4805 float zs0[EX][EY][EP],
4806 float array0[EX][EY][EP],
4807 const met_t * met1,
4808 float zs1[EX][EY][EP],
4809 float array1[EX][EY][EP],
4810 const double ts,
4811 const double p,
4812 const double lon,
4813 const double lat,
4814 double *var);
4815
4851 const met_t * met0,
4852 float array0[EX][EY],
4853 const met_t * met1,
4854 float array1[EX][EY],
4855 const double ts,
4856 const double lon,
4857 const double lat,
4858 double *var,
4859 int *ci,
4860 double *cw,
4861 const int init);
4862
4900void intpol_tropo_3d(
4901 const double time0,
4902 float array0[EX][EY],
4903 const double time1,
4904 float array1[EX][EY],
4905 const double lons[EX],
4906 const double lats[EY],
4907 const int nlon,
4908 const int nlat,
4909 const double time,
4910 const double lon,
4911 const double lat,
4912 const int method,
4913 double *var,
4914 double *sigma);
4915
4942void jsec2time(
4943 const double jsec,
4944 int *year,
4945 int *mon,
4946 int *day,
4947 int *hour,
4948 int *min,
4949 int *sec,
4950 double *remain);
4951
4978double kernel_weight(
4979 const double kz[EP],
4980 const double kw[EP],
4981 const int nk,
4982 const double p);
4983
5022double lapse_rate(
5023 const double t,
5024 const double h2o);
5025
5055 ctl_t * ctl);
5056
5076int locate_irr(
5077 const double *xx,
5078 const int n,
5079 const double x);
5080
5107 const float *xx,
5108 const int n,
5109 const double x,
5110 const int ig);
5111
5132int locate_reg(
5133 const double *xx,
5134 const int n,
5135 const double x);
5136
5158void locate_vert(
5159 float profiles[EX][EY][EP],
5160 const int np,
5161 const int lon_ap_ind,
5162 const int lat_ap_ind,
5163 const double alt_ap,
5164 int *ind);
5165
5191void module_advect(
5192 const ctl_t * ctl,
5193 const cache_t * cache,
5194 met_t * met0,
5195 met_t * met1,
5196 atm_t * atm);
5197
5221 const ctl_t * ctl,
5222 const cache_t * cache,
5223 met_t * met0,
5224 met_t * met1,
5225 atm_t * atm);
5226
5264 const ctl_t * ctl,
5265 const cache_t * cache,
5266 const clim_t * clim,
5267 met_t * met0,
5268 met_t * met1,
5269 atm_t * atm);
5270
5312void module_chem_grid(
5313 const ctl_t * ctl,
5314 met_t * met0,
5315 met_t * met1,
5316 atm_t * atm,
5317 const double t);
5318
5345void module_chem_init(
5346 const ctl_t * ctl,
5347 const cache_t * cache,
5348 const clim_t * clim,
5349 met_t * met0,
5350 met_t * met1,
5351 atm_t * atm);
5352
5377 const ctl_t * ctl,
5378 cache_t * cache,
5379 met_t * met0,
5380 met_t * met1,
5381 atm_t * atm);
5382
5409void module_decay(
5410 const ctl_t * ctl,
5411 const cache_t * cache,
5412 const clim_t * clim,
5413 atm_t * atm);
5414
5451void module_diff_meso(
5452 const ctl_t * ctl,
5453 cache_t * cache,
5454 met_t * met0,
5455 met_t * met1,
5456 atm_t * atm);
5457
5491void module_diff_pbl(
5492 const ctl_t * ctl,
5493 cache_t * cache,
5494 met_t * met0,
5495 met_t * met1,
5496 atm_t * atm);
5497
5552void module_diff_turb(
5553 const ctl_t * ctl,
5554 cache_t * cache,
5555 const clim_t * clim,
5556 met_t * met0,
5557 met_t * met1,
5558 atm_t * atm);
5559
5579void module_dry_depo(
5580 const ctl_t * ctl,
5581 const cache_t * cache,
5582 met_t * met0,
5583 met_t * met1,
5584 atm_t * atm);
5585
5618void module_h2o2_chem(
5619 const ctl_t * ctl,
5620 const cache_t * cache,
5621 const clim_t * clim,
5622 met_t * met0,
5623 met_t * met1,
5624 atm_t * atm);
5625
5646 const ctl_t * ctl,
5647 cache_t * cache,
5648 met_t * met0,
5649 met_t * met1,
5650 atm_t * atm);
5651
5669void module_isosurf(
5670 const ctl_t * ctl,
5671 const cache_t * cache,
5672 met_t * met0,
5673 met_t * met1,
5674 atm_t * atm);
5675
5708 ctl_t * ctl,
5709 cache_t * cache,
5710 clim_t * clim,
5711 met_t * met0,
5712 met_t * met1,
5713 atm_t * atm);
5714
5733void module_meteo(
5734 const ctl_t * ctl,
5735 const cache_t * cache,
5736 const clim_t * clim,
5737 met_t * met0,
5738 met_t * met1,
5739 atm_t * atm);
5740
5758void module_mixing(
5759 const ctl_t * ctl,
5760 const clim_t * clim,
5761 atm_t * atm,
5762 const double t);
5763
5792 const ctl_t * ctl,
5793 const clim_t * clim,
5794 atm_t * atm,
5795 const int *ixs,
5796 const int *iys,
5797 const int *izs,
5798 const int qnt_idx,
5799 const int use_ensemble);
5800
5833void module_oh_chem(
5834 const ctl_t * ctl,
5835 const cache_t * cache,
5836 const clim_t * clim,
5837 met_t * met0,
5838 met_t * met1,
5839 atm_t * atm);
5840
5868void module_position(
5869 const cache_t * cache,
5870 met_t * met0,
5871 met_t * met1,
5872 atm_t * atm);
5873
5898void module_rng_init(
5899 const int ntask);
5900
5926void module_rng(
5927 const ctl_t * ctl,
5928 double *rs,
5929 const size_t n,
5930 const int method);
5931
5954void module_sedi(
5955 const ctl_t * ctl,
5956 const cache_t * cache,
5957 met_t * met0,
5958 met_t * met1,
5959 atm_t * atm);
5960
5984void module_sort(
5985 const ctl_t * ctl,
5986 met_t * met0,
5987 atm_t * atm);
5988
6008void module_sort_help(
6009 double *a,
6010 const int *p,
6011 const int np);
6012
6036void module_timesteps(
6037 const ctl_t * ctl,
6038 cache_t * cache,
6039 met_t * met0,
6040 atm_t * atm,
6041 const double t);
6042
6064 ctl_t * ctl,
6065 const atm_t * atm);
6066
6100 const ctl_t * ctl,
6101 const cache_t * cache,
6102 const clim_t * clim,
6103 met_t * met0,
6104 met_t * met1,
6105 atm_t * atm);
6106
6136void module_wet_depo(
6137 const ctl_t * ctl,
6138 const cache_t * cache,
6139 met_t * met0,
6140 met_t * met1,
6141 atm_t * atm);
6142
6173void mptrac_alloc(
6174 ctl_t ** ctl,
6175 cache_t ** cache,
6176 clim_t ** clim,
6177 met_t ** met0,
6178 met_t ** met1,
6179 atm_t ** atm);
6180
6210#ifdef DD
6211void mptrac_free(
6212 ctl_t * ctl,
6213 cache_t * cache,
6214 clim_t * clim,
6215 met_t * met0,
6216 met_t * met1,
6217 atm_t * atm,
6218 mpi_info_t * mpi_info);
6219#else
6220void mptrac_free(
6221 ctl_t * ctl,
6222 cache_t * cache,
6223 clim_t * clim,
6224 met_t * met0,
6225 met_t * met1,
6226 atm_t * atm);
6227#endif
6228
6263void mptrac_get_met(
6264 ctl_t * ctl,
6265 clim_t * clim,
6266 const double t,
6267 met_t ** met0,
6268 met_t ** met1);
6269
6289void mptrac_init(
6290 ctl_t * ctl,
6291 cache_t * cache,
6292 clim_t * clim,
6293 atm_t * atm,
6294 const int ntask);
6295
6331int mptrac_read_atm(
6332 const char *filename,
6333 const ctl_t * ctl,
6334 atm_t * atm);
6335
6367void mptrac_read_clim(
6368 const ctl_t * ctl,
6369 clim_t * clim);
6370
6400void mptrac_read_ctl(
6401 const char *filename,
6402 int argc,
6403 char *argv[],
6404 ctl_t * ctl);
6405
6435int mptrac_read_met(
6436 const char *filename,
6437 const ctl_t * ctl,
6438 const clim_t * clim,
6439 met_t * met);
6440
6461#ifdef DD
6463 ctl_t * ctl,
6464 cache_t * cache,
6465 clim_t * clim,
6466 met_t ** met0,
6467 met_t ** met1,
6468 atm_t * atm,
6469 double t,
6470 mpi_info_t * mpi_info);
6471#else
6473 ctl_t * ctl,
6474 cache_t * cache,
6475 clim_t * clim,
6476 met_t ** met0,
6477 met_t ** met1,
6478 atm_t * atm,
6479 double t);
6480#endif
6481
6511void mptrac_write_atm(
6512 const char *filename,
6513 const ctl_t * ctl,
6514 const atm_t * atm,
6515 const double t);
6516
6552void mptrac_write_met(
6553 const char *filename,
6554 const ctl_t * ctl,
6555 met_t * met);
6556
6591 const char *dirname,
6592 const ctl_t * ctl,
6593 met_t * met0,
6594 met_t * met1,
6595 atm_t * atm,
6596 const double t);
6597
6629 const ctl_t * ctl,
6630 const cache_t * cache,
6631 const clim_t * clim,
6632 met_t ** met0,
6633 met_t ** met1,
6634 const atm_t * atm);
6635
6666 const ctl_t * ctl,
6667 const cache_t * cache,
6668 const clim_t * clim,
6669 met_t ** met0,
6670 met_t ** met1,
6671 const atm_t * atm);
6672
6700double nat_temperature(
6701 const double p,
6702 const double h2o,
6703 const double hno3);
6704
6725double pbl_weight(
6726 const ctl_t * ctl,
6727 const atm_t * atm,
6728 const int ip,
6729 const double pbl,
6730 const double ps);
6731
6764int read_atm_asc(
6765 const char *filename,
6766 const ctl_t * ctl,
6767 atm_t * atm);
6768
6799int read_atm_bin(
6800 const char *filename,
6801 const ctl_t * ctl,
6802 atm_t * atm);
6803
6828int read_atm_clams(
6829 const char *filename,
6830 const ctl_t * ctl,
6831 atm_t * atm);
6832
6862int read_atm_nc(
6863 const char *filename,
6864 const ctl_t * ctl,
6865 atm_t * atm);
6866
6895void read_clim_photo(
6896 const char *filename,
6897 clim_photo_t * photo);
6898
6916 const int ncid,
6917 const char *varname,
6918 const clim_photo_t * photo,
6919 double var[CP][CSZA][CO3]);
6920
6944int read_clim_ts(
6945 const char *filename,
6946 clim_ts_t * ts);
6947
6974void read_clim_zm(
6975 const char *filename,
6976 const char *varname,
6977 clim_zm_t * zm);
6978
7006void read_kernel(
7007 const char *filename,
7008 double kz[EP],
7009 double kw[EP],
7010 int *nk);
7011
7043int read_met_bin(
7044 const char *filename,
7045 const ctl_t * ctl,
7046 met_t * met);
7047
7073void read_met_bin_2d(
7074 FILE * in,
7075 const met_t * met,
7076 float var[EX][EY],
7077 const char *varname);
7078
7116void read_met_bin_3d(
7117 FILE * in,
7118 const ctl_t * ctl,
7119 const met_t * met,
7120 float var[EX][EY][EP],
7121 const char *varname,
7122 const float bound_min,
7123 const float bound_max);
7124
7151void read_met_cape(
7152 const ctl_t * ctl,
7153 const clim_t * clim,
7154 met_t * met);
7155
7178void read_met_cloud(
7179 met_t * met);
7180
7206void read_met_detrend(
7207 const ctl_t * ctl,
7208 met_t * met);
7209
7233 met_t * met);
7234
7261void read_met_geopot(
7262 const ctl_t * ctl,
7263 met_t * met);
7264
7292 const char *filename,
7293 const ctl_t * ctl,
7294 met_t * met);
7295
7317#ifdef ECCODES
7318void read_met_grib_grid(
7319 codes_handle ** handles,
7320 int count_handles,
7321 met_t * met);
7322#endif
7323
7343#ifdef ECCODES
7344void read_met_grib_levels(
7345 codes_handle ** handles,
7346 const int num_messages,
7347 const ctl_t * ctl,
7348 met_t * met);
7349#endif
7350
7380#ifdef ECCODES
7381void read_met_grib_surface(
7382 codes_handle ** handles,
7383 const int num_messages,
7384 const ctl_t * ctl,
7385 met_t * met);
7386#endif
7387
7416void read_met_ml2pl(
7417 const ctl_t * ctl,
7418 const met_t * met,
7419 float var[EX][EY][EP],
7420 const char *varname);
7421
7444 const ctl_t * ctl,
7445 met_t * met);
7446
7476int read_met_nc(
7477 const char *filename,
7478 const ctl_t * ctl,
7479 met_t * met);
7480
7512 const int ncid,
7513 const ctl_t * ctl,
7514 met_t * met);
7515
7557 const int ncid,
7558 const ctl_t * ctl,
7559 met_t * met);
7560
7592void read_met_nc_grid(
7593 const char *filename,
7594 const int ncid,
7595 const ctl_t * ctl,
7596 met_t * met);
7597
7628 const char *filename,
7629 const ctl_t * ctl,
7630 met_t * met);
7631
7665 const char *filename,
7666 const int ncid,
7667 const ctl_t * ctl,
7668 met_t * met);
7669
7701 const int ncid,
7702 const ctl_t * ctl,
7703 met_t * met);
7704
7734 const int ncid,
7735 const ctl_t * ctl,
7736 met_t * met);
7737
7767int read_met_nc_2d(
7768 const int ncid,
7769 const char *varname,
7770 const char *varname2,
7771 const char *varname3,
7772 const char *varname4,
7773 const char *varname5,
7774 const char *varname6,
7775 const ctl_t * ctl,
7776 const met_t * met,
7777 float dest[EX][EY],
7778 const float scl,
7779 const int init);
7780
7811 const int ncid,
7812 const char *varname,
7813 const char *varname2,
7814 const char *varname3,
7815 const char *varname4,
7816 const char *varname5,
7817 const char *varname6,
7818 const ctl_t * ctl,
7819 const met_t * met,
7820 float dest[EX][EY],
7821 const float scl,
7822 const int init);
7823
7854int read_met_nc_3d(
7855 const int ncid,
7856 const char *varname,
7857 const char *varname2,
7858 const char *varname3,
7859 const char *varname4,
7860 const ctl_t * ctl,
7861 const met_t * met,
7862 float dest[EX][EY][EP],
7863 const float scl);
7864
7896 const int ncid,
7897 const char *varname,
7898 const char *varname2,
7899 const char *varname3,
7900 const char *varname4,
7901 const ctl_t * ctl,
7902 const met_t * met,
7903 float dest[EX][EY][EP],
7904 const float scl);
7905
7951void read_met_pbl(
7952 const ctl_t * ctl,
7953 met_t * met);
7954
7987 met_t * met);
7988
8019 met_t * met);
8020
8051void read_met_pv(
8052 met_t * met);
8053
8076void read_met_ozone(
8077 met_t * met);
8078
8107void read_met_sample(
8108 const ctl_t * ctl,
8109 met_t * met);
8110
8139void read_met_tropo(
8140 const ctl_t * ctl,
8141 const clim_t * clim,
8142 met_t * met);
8143
8175void read_obs(
8176 const char *filename,
8177 const ctl_t * ctl,
8178 double *rt,
8179 double *rz,
8180 double *rlon,
8181 double *rlat,
8182 double *robs,
8183 int *nobs);
8184
8212void read_obs_asc(
8213 const char *filename,
8214 double *rt,
8215 double *rz,
8216 double *rlon,
8217 double *rlat,
8218 double *robs,
8219 int *nobs);
8220
8247void read_obs_nc(
8248 const char *filename,
8249 double *rt,
8250 double *rz,
8251 double *rlon,
8252 double *rlat,
8253 double *robs,
8254 int *nobs);
8255
8289double scan_ctl(
8290 const char *filename,
8291 int argc,
8292 char *argv[],
8293 const char *varname,
8294 const int arridx,
8295 const char *defvalue,
8296 char *value);
8297
8324double sedi(
8325 const double p,
8326 const double T,
8327 const double rp,
8328 const double rhop);
8329
8359void spline(
8360 const double *x,
8361 const double *y,
8362 const int n,
8363 const double *x2,
8364 double *y2,
8365 const int n2,
8366 const int method);
8367
8390float stddev(
8391 const float *data,
8392 const int n);
8393
8413double sza_calc(
8414 const double sec,
8415 const double lon,
8416 const double lat);
8417
8442void time2jsec(
8443 const int year,
8444 const int mon,
8445 const int day,
8446 const int hour,
8447 const int min,
8448 const int sec,
8449 const double remain,
8450 double *jsec);
8451
8480void timer(
8481 const char *name,
8482 const char *group,
8483 const int output);
8484
8510double time_from_filename(
8511 const char *filename,
8512 const int offset);
8513
8531double tropo_weight(
8532 const clim_t * clim,
8533 const atm_t * atm,
8534 const int ip);
8535
8558void write_atm_asc(
8559 const char *filename,
8560 const ctl_t * ctl,
8561 const atm_t * atm,
8562 const double t);
8563
8587void write_atm_bin(
8588 const char *filename,
8589 const ctl_t * ctl,
8590 const atm_t * atm);
8591
8615void write_atm_clams(
8616 const char *filename,
8617 const ctl_t * ctl,
8618 const atm_t * atm);
8619
8645 const char *dirname,
8646 const ctl_t * ctl,
8647 const atm_t * atm,
8648 const double t);
8649
8673void write_atm_nc(
8674 const char *filename,
8675 const ctl_t * ctl,
8676 const atm_t * atm);
8677
8706void write_csi(
8707 const char *filename,
8708 const ctl_t * ctl,
8709 const atm_t * atm,
8710 const double t);
8711
8753 const char *filename,
8754 const ctl_t * ctl,
8755 const atm_t * atm,
8756 const double t);
8757
8785void write_ens(
8786 const char *filename,
8787 const ctl_t * ctl,
8788 const atm_t * atm,
8789 const double t);
8790
8829void write_grid(
8830 const char *filename,
8831 const ctl_t * ctl,
8832 met_t * met0,
8833 met_t * met1,
8834 const atm_t * atm,
8835 const double t);
8836
8882void write_grid_asc(
8883 const char *filename,
8884 const ctl_t * ctl,
8885 const double *cd,
8886 double *mean[NQ],
8887 double *sigma[NQ],
8888 const double *vmr_impl,
8889 const double t,
8890 const double *z,
8891 const double *lon,
8892 const double *lat,
8893 const double *area,
8894 const double dz,
8895 const int *np);
8896
8939void write_grid_nc(
8940 const char *filename,
8941 const ctl_t * ctl,
8942 const double *cd,
8943 double *mean[NQ],
8944 double *sigma[NQ],
8945 const double *vmr_impl,
8946 const double t,
8947 const double *z,
8948 const double *lon,
8949 const double *lat,
8950 const double *area,
8951 const double dz,
8952 const int *np);
8953
8983void write_met_bin(
8984 const char *filename,
8985 const ctl_t * ctl,
8986 met_t * met);
8987
9015void write_met_bin_2d(
9016 FILE * out,
9017 met_t * met,
9018 float var[EX][EY],
9019 const char *varname);
9020
9058void write_met_bin_3d(
9059 FILE * out,
9060 const ctl_t * ctl,
9061 met_t * met,
9062 float var[EX][EY][EP],
9063 const char *varname,
9064 const int precision,
9065 const double tolerance);
9066
9094void write_met_nc(
9095 const char *filename,
9096 const ctl_t * ctl,
9097 met_t * met);
9098
9123void write_met_nc_2d(
9124 const int ncid,
9125 const char *varname,
9126 met_t * met,
9127 float var[EX][EY],
9128 const float scl);
9129
9155void write_met_nc_3d(
9156 const int ncid,
9157 const char *varname,
9158 met_t * met,
9159 float var[EX][EY][EP],
9160 const float scl);
9161
9192void write_prof(
9193 const char *filename,
9194 const ctl_t * ctl,
9195 met_t * met0,
9196 met_t * met1,
9197 const atm_t * atm,
9198 const double t);
9199
9231void write_sample(
9232 const char *filename,
9233 const ctl_t * ctl,
9234 met_t * met0,
9235 met_t * met1,
9236 const atm_t * atm,
9237 const double t);
9238
9265void write_station(
9266 const char *filename,
9267 const ctl_t * ctl,
9268 atm_t * atm,
9269 const double t);
9270
9299void write_vtk(
9300 const char *filename,
9301 const ctl_t * ctl,
9302 const atm_t * atm,
9303 const double t);
9304
9328 atm_t * atm,
9329 particle_t * particles,
9330 ctl_t * ctl,
9331 int *nparticles,
9332 cache_t * cache,
9333 int rank);
9334
9360 atm_t * atm,
9361 particle_t * particles,
9362 ctl_t * ctl,
9363 int *nparticles,
9364 cache_t * cache);
9365
9387#ifdef DD
9388void dd_register_MPI_type_particle(
9389 MPI_Datatype * MPI_Particle);
9390#endif
9391
9414#ifdef DD
9415void dd_get_rect_neighbour(
9416 const ctl_t ctl,
9417 mpi_info_t * mpi_info);
9418#endif
9419
9449#ifdef DD
9450void dd_communicate_particles(
9451 particle_t * particles,
9452 int *nparticles,
9453 MPI_Datatype MPI_Particle,
9454 int *neighbours,
9455 int nneighbours,
9456 ctl_t ctl);
9457#endif
9458
9486#ifdef DD
9487void dd_assign_rect_subdomains_atm(
9488 atm_t * atm,
9489 met_t * met,
9490 ctl_t * ctl,
9491 mpi_info_t * mpi_info,
9492 int init);
9493#endif
9494
9522#ifdef DD
9523void dd_init(
9524 ctl_t * ctl,
9525 mpi_info_t * mpi_info,
9526 atm_t * atm,
9527 met_t ** met,
9528 double t,
9529 int *dd_init);
9530#endif
9531
9559#ifdef DD
9560void module_dd(
9561 ctl_t * ctl,
9562 atm_t * atm,
9563 cache_t * cache,
9564 mpi_info_t * mpi_info,
9565 met_t ** met);
9566#endif
9567
9601#ifdef DD
9602void dd_sort(
9603 const ctl_t * ctl,
9604 met_t * met0,
9605 atm_t * atm,
9606 int *nparticles,
9607 int *rank);
9608#endif
9609
9610#endif /* LIBTRAC_H */
void read_met_geopot(const ctl_t *ctl, met_t *met)
Calculates geopotential heights from meteorological data.
Definition: mptrac.c:7148
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:9522
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:309
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:5828
void day2doy(const int year, const int mon, const int day, int *doy)
Get day of year from date.
Definition: mptrac.c:906
void read_met_extrapolate(met_t *met)
Extrapolates meteorological data.
Definition: mptrac.c:7108
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:11051
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:12437
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:4212
void read_met_sample(const ctl_t *ctl, met_t *met)
Downsamples meteorological data based on specified parameters.
Definition: mptrac.c:10050
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:12197
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:10394
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:2094
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:3934
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:5445
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:3260
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:6256
void read_met_cloud(met_t *met)
Calculates cloud-related variables for each grid point.
Definition: mptrac.c:6947
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:2647
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:10567
#define METVAR
Number of 3-D meteorological variables.
Definition: mptrac.h:314
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:1468
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:1401
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:6223
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:9658
void read_met_detrend(const ctl_t *ctl, met_t *met)
Detrends meteorological data.
Definition: mptrac.c:7004
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:8858
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:10222
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:10438
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:2532
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:2055
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:1101
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:2074
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:4301
void read_met_monotonize(const ctl_t *ctl, met_t *met)
Makes zeta and pressure profiles monotone.
Definition: mptrac.c:8742
#define EP_GLOB
Maximum number of pressure levels for meteo data.
Definition: mptrac.h:294
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:6375
void read_met_periodic(met_t *met)
Applies periodic boundary conditions to meteorological data along longitudinal axis.
Definition: mptrac.c:9795
void module_timesteps_init(ctl_t *ctl, const atm_t *atm)
Initialize start time and time interval for time-stepping.
Definition: mptrac.c:3976
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:11534
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:3365
void compress_zfp(const char *varname, float *array, const int nx, const int ny, const int nz, const int precision, const double tolerance, const int decompress, FILE *inout)
Compresses or decompresses a 3D array of floats using the ZFP library.
double clim_zm(const clim_zm_t *zm, const double t, const double lat, const double p)
Interpolates monthly mean zonal mean climatological variables.
Definition: mptrac.c:401
#define EY
Maximum number of latitudes for meteo data.
Definition: mptrac.h:289
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:3436
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:6347
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:8700
double clim_tropo(const clim_t *clim, const double t, const double lat)
Calculates the tropopause pressure based on climatological data.
Definition: mptrac.c:200
void read_obs_nc(const char *filename, double *rt, double *rz, double *rlon, double *rlat, double *robs, int *nobs)
Reads observation data from a NetCDF file.
Definition: mptrac.c:10466
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:6717
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:1991
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:3085
void level_definitions(ctl_t *ctl)
Defines pressure levels for meteorological data.
Definition: mptrac.c:1781
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:11822
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:5716
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:10709
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:4422
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:1526
void fft_help(double *fcReal, double *fcImag, const int n)
Computes the Fast Fourier Transform (FFT) of a complex sequence.
Definition: mptrac.c:955
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:4077
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:664
#define EY_GLOB
Maximum number of global latitudes for meteo data.
Definition: mptrac.h:304
double nat_temperature(const double p, const double h2o, const double hno3)
Calculates the nitric acid trihydrate (NAT) temperature.
Definition: mptrac.c:6019
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:10600
#define CP
Maximum number of pressure levels for climatological data.
Definition: mptrac.h:359
#define NQ
Maximum number of quantities per data point.
Definition: mptrac.h:324
double clim_photo(const double rate[CP][CSZA][CO3], const clim_photo_t *photo, const double p, const double sza, const double o3c)
Calculates the photolysis rate for a given set of atmospheric conditions.
Definition: mptrac.c:149
void read_clim_zm(const char *filename, const char *varname, clim_zm_t *zm)
Reads zonally averaged climatological data from a netCDF file and populates the given structure.
Definition: mptrac.c:6429
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:9174
#define EX
Maximum number of longitudes for meteo data.
Definition: mptrac.h:284
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:3809
void timer(const char *name, const char *group, const int output)
Measures and reports elapsed time for named and grouped timers.
Definition: mptrac.c:10740
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:10866
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:1343
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:1555
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:2574
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:6528
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:2274
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:10495
#define CT
Maximum number of time steps for climatological data.
Definition: mptrac.h:369
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:2247
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:3898
float stddev(const float *data, const int n)
Calculates the standard deviation of a set of data.
Definition: mptrac.c:10647
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:1613
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:9022
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:6746
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:2021
double time_from_filename(const char *filename, const int offset)
Extracts and converts a timestamp from a filename to Julian seconds.
Definition: mptrac.c:10808
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:12496
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:4512
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:2370
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:10843
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:12273
void module_rng_init(const int ntask)
Initialize random number generators for parallel tasks.
Definition: mptrac.c:3673
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:4441
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:1171
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:5772
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:5928
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:374
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:10021
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:4269
void clim_tropo_init(clim_t *clim)
Initializes the tropopause data in the climatology structure.
Definition: mptrac.c:228
void module_rng(const ctl_t *ctl, double *rs, const size_t n, const int method)
Generate random numbers using various methods and distributions.
Definition: mptrac.c:3704
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:1012
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:12885
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:319
double sza_calc(const double sec, const double lon, const double lat)
Calculates the solar zenith angle.
Definition: mptrac.c:10668
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:925
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:1580
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:3622
void clim_oh_diurnal_correction(const ctl_t *ctl, clim_t *clim)
Applies a diurnal correction to the hydroxyl radical (OH) concentration in climatology data.
Definition: mptrac.c:116
void write_met_bin_2d(FILE *out, met_t *met, float var[EX][EY], const char *varname)
Writes a 2-dimensional meteorological variable to a binary file.
Definition: mptrac.c:12168
void read_met_pv(met_t *met)
Calculates potential vorticity (PV) from meteorological data.
Definition: mptrac.c:9915
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:6107
void get_met_replace(char *orig, char *search, char *repl)
Replaces occurrences of a substring in a string with another substring.
Definition: mptrac.c:1077
#define CY
Maximum number of latitudes for climatological data.
Definition: mptrac.h:349
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:2686
void module_sort(const ctl_t *ctl, met_t *met0, atm_t *atm)
Sort particles according to box index.
Definition: mptrac.c:3838
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:8827
double clim_ts(const clim_ts_t *ts, const double t)
Interpolates a time series of climatological variables.
Definition: mptrac.c:383
void thrustSortWrapper(double *__restrict__ c, int n, int *__restrict__ index)
Wrapper to Thrust sorting function.
void jsec2time(const double jsec, int *year, int *mon, int *day, int *hour, int *min, int *sec, double *remain)
Converts Julian seconds to calendar date and time components.
Definition: mptrac.c:1704
int read_met_bin(const char *filename, const ctl_t *ctl, met_t *met)
Reads meteorological data from a binary file.
Definition: mptrac.c:6569
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:5571
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:10998
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:2888
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:6163
#define EX_GLOB
Maximum number of global longitudes for meteo data.
Definition: mptrac.h:299
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:12971
#define EP
Maximum number of pressure levels for meteo data.
Definition: mptrac.h:279
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:4007
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:4572
void read_met_polar_winds(met_t *met)
Applies a fix for polar winds in meteorological data.
Definition: mptrac.c:9856
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:3003
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:11926
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:6043
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:2763
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:12466
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:3155
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:354
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:3538
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:994
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:1737
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:9308
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:6065
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:1144
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:12723
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:11631
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:2940
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:12056
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:10948
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:6832
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:5888
#define CSZA
Maximum number of solar zenith angles for climatological data.
Definition: mptrac.h:364
double lapse_rate(const double t, const double h2o)
Calculates the moist adiabatic lapse rate in Kelvin per kilometer.
Definition: mptrac.c:1763
#define DD_NNMAX
Maximum number of neighbours to communicate with.
Definition: mptrac.h:384
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:11258
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:11209
Air parcel data.
Definition: mptrac.h:3371
int np
Number of air parcels.
Definition: mptrac.h:3374
Cache data structure.
Definition: mptrac.h:3449
int iso_n
Isosurface balloon number of data points.
Definition: mptrac.h:3461
Climatological data in the form of photolysis rates.
Definition: mptrac.h:3481
int nsza
Number of solar zenith angles.
Definition: mptrac.h:3487
int np
Number of pressure levels.
Definition: mptrac.h:3484
int no3c
Number of total ozone columns.
Definition: mptrac.h:3490
Climatological data.
Definition: mptrac.h:3589
clim_ts_t ccl2f2
CFC-12 time series.
Definition: mptrac.h:3631
clim_photo_t photo
Photolysis rates.
Definition: mptrac.h:3607
clim_zm_t ho2
HO2 zonal means.
Definition: mptrac.h:3619
clim_zm_t hno3
HNO3 zonal means.
Definition: mptrac.h:3610
int tropo_ntime
Number of tropopause timesteps.
Definition: mptrac.h:3592
clim_ts_t sf6
SF6 time series.
Definition: mptrac.h:3637
clim_ts_t ccl4
CFC-10 time series.
Definition: mptrac.h:3625
clim_ts_t ccl3f
CFC-11 time series.
Definition: mptrac.h:3628
clim_zm_t o1d
O(1D) zonal means.
Definition: mptrac.h:3622
clim_zm_t h2o2
H2O2 zonal means.
Definition: mptrac.h:3616
int tropo_nlat
Number of tropopause latitudes.
Definition: mptrac.h:3595
clim_zm_t oh
OH zonal means.
Definition: mptrac.h:3613
clim_ts_t n2o
N2O time series.
Definition: mptrac.h:3634
Climatological data in the form of time series.
Definition: mptrac.h:3537
int ntime
Number of timesteps.
Definition: mptrac.h:3540
Climatological data in the form of zonal means.
Definition: mptrac.h:3557
int np
Number of pressure levels.
Definition: mptrac.h:3566
int ntime
Number of timesteps.
Definition: mptrac.h:3560
int nlat
Number of latitudes.
Definition: mptrac.h:3563
Control parameters.
Definition: mptrac.h:2356
double grid_z0
Lower altitude of gridded data [km].
Definition: mptrac.h:3239
int qnt_o3
Quantity array index for ozone volume mixing ratio.
Definition: mptrac.h:2468
double csi_lat1
Upper latitude of gridded CSI data [deg].
Definition: mptrac.h:3200
int qnt_Coh
Quantity array index for OH volume mixing ratio (chemistry code).
Definition: mptrac.h:2618
double wet_depo_ic_a
Coefficient A for wet deposition in cloud (exponential form).
Definition: mptrac.h:3088
int met_nc_scale
Check netCDF scaling factors (0=no, 1=yes).
Definition: mptrac.h:2696
int qnt_pel
Quantity array index for pressure at equilibrium level (EL).
Definition: mptrac.h:2501
int csi_nz
Number of altitudes of gridded CSI data.
Definition: mptrac.h:3176
double molmass
Molar mass [g/mol].
Definition: mptrac.h:2950
int qnt_p
Quantity array index for pressure.
Definition: mptrac.h:2447
int qnt_Cccl2f2
Quantity array index for CFC-12 volume mixing ratio (chemistry code).
Definition: mptrac.h:2642
int dd_halos_size
Size of halos given in grid-points.
Definition: mptrac.h:3359
int mixing_nx
Number of longitudes of mixing grid.
Definition: mptrac.h:3013
double chemgrid_z1
Upper altitude of chemistry grid [km].
Definition: mptrac.h:3037
int qnt_m
Quantity array index for mass.
Definition: mptrac.h:2387
int qnt_aoa
Quantity array index for age of air.
Definition: mptrac.h:2651
int qnt_rhop
Quantity array index for particle density.
Definition: mptrac.h:2396
int qnt_swc
Quantity array index for cloud snow water content.
Definition: mptrac.h:2480
double csi_obsmin
Minimum observation index to trigger detection.
Definition: mptrac.h:3170
int qnt_pcb
Quantity array index for cloud bottom pressure.
Definition: mptrac.h:2489
double bound_dzs
Boundary conditions surface layer depth [km].
Definition: mptrac.h:2938
double csi_lon1
Upper longitude of gridded CSI data [deg].
Definition: mptrac.h:3191
int qnt_u
Quantity array index for zonal wind.
Definition: mptrac.h:2456
double stat_lon
Longitude of station [deg].
Definition: mptrac.h:3317
double mixing_trop
Interparcel exchange parameter for mixing in the troposphere.
Definition: mptrac.h:2998
double sort_dt
Time step for sorting of particle data [s].
Definition: mptrac.h:2849
double mixing_z1
Upper altitude of mixing grid [km].
Definition: mptrac.h:3010
double stat_r
Search radius around station [km].
Definition: mptrac.h:3323
double wet_depo_bc_a
Coefficient A for wet deposition below cloud (exponential form).
Definition: mptrac.h:3082
int met_zstd_level
ZSTD compression level (from -5 to 22).
Definition: mptrac.h:2705
int csi_ny
Number of latitudes of gridded CSI data.
Definition: mptrac.h:3194
int vtk_sphere
Spherical projection for VTK data (0=no, 1=yes).
Definition: mptrac.h:3347
double chemgrid_z0
Lower altitude of chemistry grid [km].
Definition: mptrac.h:3034
double met_pbl_min
Minimum depth of planetary boundary layer [km].
Definition: mptrac.h:2817
int qnt_iwc
Quantity array index for cloud ice water content.
Definition: mptrac.h:2477
double chemgrid_lat0
Lower latitude of chemistry grid [deg].
Definition: mptrac.h:3052
double conv_cape
CAPE threshold for convection module [J/kg].
Definition: mptrac.h:2902
int qnt_Co1d
Quantity array index for O(1D) volume mixing ratio (chemistry code).
Definition: mptrac.h:2630
double met_cms_eps_pv
cmultiscale compression epsilon for potential vorticity.
Definition: mptrac.h:2739
int qnt_pw
Quantity array index for partial water vapor pressure.
Definition: mptrac.h:2555
double grid_z1
Upper altitude of gridded data [km].
Definition: mptrac.h:3242
int direction
Direction flag (1=forward calculation, -1=backward calculation).
Definition: mptrac.h:2660
int qnt_Cccl4
Quantity array index for CFC-10 volume mixing ratio (chemistry code).
Definition: mptrac.h:2636
int met_dp
Stride for pressure levels.
Definition: mptrac.h:2769
double met_dt_out
Time step for sampling of meteo data along trajectories [s].
Definition: mptrac.h:2836
int qnt_h2o2
Quantity array index for H2O2 volume mixing ratio (climatology).
Definition: mptrac.h:2519
int qnt_vh
Quantity array index for horizontal wind.
Definition: mptrac.h:2585
int csi_nx
Number of longitudes of gridded CSI data.
Definition: mptrac.h:3185
double csi_lat0
Lower latitude of gridded CSI data [deg].
Definition: mptrac.h:3197
double turb_dz_trop
Vertical turbulent diffusion coefficient (troposphere) [m^2/s].
Definition: mptrac.h:2884
int met_pbl
Planetary boundary layer data (0=file, 1=z2p, 2=Richardson, 3=theta).
Definition: mptrac.h:2814
int qnt_lwc
Quantity array index for cloud liquid water content.
Definition: mptrac.h:2471
double turb_mesoz
Vertical scaling factor for mesoscale wind fluctuations.
Definition: mptrac.h:2893
int grid_nc_level
zlib compression level of netCDF grid data files (0=off).
Definition: mptrac.h:3227
int grid_nx
Number of longitudes of gridded data.
Definition: mptrac.h:3245
int atm_type
Type of atmospheric data files (0=ASCII, 1=binary, 2=netCDF, 3=CLaMS_traj, 4=CLaMS_pos).
Definition: mptrac.h:3141
double bound_mass
Boundary conditions mass per particle [kg].
Definition: mptrac.h:2911
double grid_lat0
Lower latitude of gridded data [deg].
Definition: mptrac.h:3257
int qnt_ts
Quantity array index for surface temperature.
Definition: mptrac.h:2402
int qnt_loss_rate
Quantity array index for total loss rate.
Definition: mptrac.h:2546
double met_cms_eps_h2o
cmultiscale compression epsilon for water vapor.
Definition: mptrac.h:2742
int qnt_plfc
Quantity array index for pressure at level of free convection (LCF).
Definition: mptrac.h:2498
double grid_lon0
Lower longitude of gridded data [deg].
Definition: mptrac.h:3248
int qnt_o1d
Quantity array index for O(1D) volume mixing ratio (climatology).
Definition: mptrac.h:2525
int met_tropo_spline
Tropopause interpolation method (0=linear, 1=spline).
Definition: mptrac.h:2833
int qnt_tvirt
Quantity array index for virtual temperature.
Definition: mptrac.h:2579
double dt_met
Time step of meteo data [s].
Definition: mptrac.h:2679
double chemgrid_lat1
Upper latitude of chemistry grid [deg].
Definition: mptrac.h:3055
int met_geopot_sy
Latitudinal smoothing of geopotential heights.
Definition: mptrac.h:2805
double met_cms_eps_u
cmultiscale compression epsilon for zonal wind.
Definition: mptrac.h:2730
double turb_dx_strat
Horizontal turbulent diffusion coefficient (stratosphere) [m^2/s].
Definition: mptrac.h:2878
int qnt_vmr
Quantity array index for volume mixing ratio.
Definition: mptrac.h:2390
int qnt_lsm
Quantity array index for land-sea mask.
Definition: mptrac.h:2423
int qnt_theta
Quantity array index for potential temperature.
Definition: mptrac.h:2567
double bound_lat1
Boundary conditions maximum longitude [deg].
Definition: mptrac.h:2926
double stat_t1
Stop time for station output [s].
Definition: mptrac.h:3329
double turb_dx_trop
Horizontal turbulent diffusion coefficient (troposphere) [m^2/s].
Definition: mptrac.h:2875
int grid_type
Type of grid data files (0=ASCII, 1=netCDF).
Definition: mptrac.h:3263
double csi_lon0
Lower longitude of gridded CSI data [deg].
Definition: mptrac.h:3188
int qnt_pbl
Quantity array index for boundary layer pressure.
Definition: mptrac.h:2429
int grid_stddev
Include standard deviations in grid output (0=no, 1=yes).
Definition: mptrac.h:3233
int qnt_psice
Quantity array index for saturation pressure over ice.
Definition: mptrac.h:2552
double chemgrid_lon0
Lower longitude of chemistry grid [deg].
Definition: mptrac.h:3043
int bound_pbl
Boundary conditions planetary boundary layer (0=no, 1=yes).
Definition: mptrac.h:2944
int qnt_mloss_wet
Quantity array index for total mass loss due to wet deposition.
Definition: mptrac.h:2537
int met_geopot_sx
Longitudinal smoothing of geopotential heights.
Definition: mptrac.h:2802
int met_sy
Smoothing for latitudes.
Definition: mptrac.h:2775
int qnt_ps
Quantity array index for surface pressure.
Definition: mptrac.h:2399
int rng_type
Random number generator (0=GSL, 1=Squares, 2=cuRAND).
Definition: mptrac.h:2866
int isosurf
Isosurface parameter (0=none, 1=pressure, 2=density, 3=theta, 4=balloon).
Definition: mptrac.h:2853
double bound_p1
Boundary conditions top pressure [hPa].
Definition: mptrac.h:2932
int qnt_zs
Quantity array index for surface geopotential height.
Definition: mptrac.h:2405
int prof_nz
Number of altitudes of gridded profile data.
Definition: mptrac.h:3272
double csi_dt_out
Time step for CSI output [s].
Definition: mptrac.h:3164
int met_cape
Convective available potential energy data (0=file, 1=calculate).
Definition: mptrac.h:2811
double csi_modmin
Minimum column density to trigger detection [kg/m^2].
Definition: mptrac.h:3173
int met_sx
Smoothing for longitudes.
Definition: mptrac.h:2772
double chemgrid_lon1
Upper longitude of chemistry grid [deg].
Definition: mptrac.h:3046
double turb_mesox
Horizontal scaling factor for mesoscale wind fluctuations.
Definition: mptrac.h:2890
double met_cms_eps_iwc
cmultiscale compression epsilon for cloud ice water content.
Definition: mptrac.h:2754
double met_cms_eps_swc
cmultiscale compression epsilon for cloud snow water content.
Definition: mptrac.h:2757
double met_cms_eps_v
cmultiscale compression epsilon for meridional wind.
Definition: mptrac.h:2733
double prof_z0
Lower altitude of gridded profile data [km].
Definition: mptrac.h:3275
int qnt_w
Quantity array index for vertical velocity.
Definition: mptrac.h:2462
double bound_vmr
Boundary conditions volume mixing ratio [ppv].
Definition: mptrac.h:2917
double met_tropo_pv
Dynamical tropopause potential vorticity threshold [PVU].
Definition: mptrac.h:2827
int prof_nx
Number of longitudes of gridded profile data.
Definition: mptrac.h:3281
int qnt_stat
Quantity array index for station flag.
Definition: mptrac.h:2384
int met_tropo
Tropopause definition (0=none, 1=clim, 2=cold point, 3=WMO_1st, 4=WMO_2nd, 5=dynamical).
Definition: mptrac.h:2824
int qnt_rp
Quantity array index for particle radius.
Definition: mptrac.h:2393
int met_mpi_share
Use MPI to share meteo (0=no, 1=yes).
Definition: mptrac.h:2842
double mixing_strat
Interparcel exchange parameter for mixing in the stratosphere.
Definition: mptrac.h:3001
int qnt_vz
Quantity array index for vertical velocity.
Definition: mptrac.h:2588
int qnt_ho2
Quantity array index for HO2 volume mixing ratio (climatology).
Definition: mptrac.h:2522
double csi_z1
Upper altitude of gridded CSI data [km].
Definition: mptrac.h:3182
double stat_t0
Start time for station output [s].
Definition: mptrac.h:3326
double oh_chem_beta
Beta parameter for diurnal variablity of OH.
Definition: mptrac.h:3064
double wet_depo_so2_ph
pH value used to calculate effective Henry constant of SO2.
Definition: mptrac.h:3100
double mixing_z0
Lower altitude of mixing grid [km].
Definition: mptrac.h:3007
int qnt_mloss_decay
Quantity array index for total mass loss due to exponential decay.
Definition: mptrac.h:2543
int atm_type_out
Type of atmospheric data files for output (-1=same as ATM_TYPE, 0=ASCII, 1=binary,...
Definition: mptrac.h:3146
int met_nlev
Number of meteo data model levels.
Definition: mptrac.h:2793
double dt_kpp
Time step for KPP chemistry [s].
Definition: mptrac.h:3073
double dry_depo_dp
Dry deposition surface layer [hPa].
Definition: mptrac.h:3109
int qnt_shf
Quantity array index for surface sensible heat flux.
Definition: mptrac.h:2420
int qnt_vs
Quantity array index for surface meridional wind.
Definition: mptrac.h:2411
int qnt_Cco
Quantity array index for CO volume mixing ratio (chemistry code).
Definition: mptrac.h:2615
double vtk_dt_out
Time step for VTK data output [s].
Definition: mptrac.h:3335
double t_stop
Stop time of simulation [s].
Definition: mptrac.h:2666
double conv_dt
Time interval for convection module [s].
Definition: mptrac.h:2908
int qnt_hno3
Quantity array index for HNO3 volume mixing ratio (climatology).
Definition: mptrac.h:2513
int met_clams
Read MPTRAC or CLaMS meteo data (0=MPTRAC, 1=CLaMS).
Definition: mptrac.h:2693
int qnt_h2ot
Quantity array index for tropopause water vapor volume mixing ratio.
Definition: mptrac.h:2441
int qnt_rh
Quantity array index for relative humidity over water.
Definition: mptrac.h:2561
double met_cms_eps_cc
cmultiscale compression epsilon for cloud cover.
Definition: mptrac.h:2760
double bound_lat0
Boundary conditions minimum longitude [deg].
Definition: mptrac.h:2923
double met_pbl_max
Maximum depth of planetary boundary layer [km].
Definition: mptrac.h:2820
int met_dx
Stride for longitudes.
Definition: mptrac.h:2763
int qnt_destination
Quantity array index for destination subdomain in domain decomposition.
Definition: mptrac.h:2657
int mixing_ny
Number of latitudes of mixing grid.
Definition: mptrac.h:3022
int met_convention
Meteo data layout (0=[lev, lat, lon], 1=[lon, lat, lev]).
Definition: mptrac.h:2682
int qnt_zeta_d
Quantity array index for diagnosed zeta vertical coordinate.
Definition: mptrac.h:2573
int tracer_chem
Switch for first order tracer chemistry module (0=off, 1=on).
Definition: mptrac.h:3076
double dt_mod
Time step of simulation [s].
Definition: mptrac.h:2669
int diffusion
Diffusion scheme (0=off, 1=fixed-K, 2=PBL).
Definition: mptrac.h:2869
int qnt_tnat
Quantity array index for T_NAT.
Definition: mptrac.h:2603
int qnt_tice
Quantity array index for T_ice.
Definition: mptrac.h:2597
int qnt_zg
Quantity array index for geopotential height.
Definition: mptrac.h:2444
double vtk_offset
Vertical offset for VTK data [km].
Definition: mptrac.h:3344
int qnt_v
Quantity array index for meridional wind.
Definition: mptrac.h:2459
int qnt_mloss_dry
Quantity array index for total mass loss due to dry deposition.
Definition: mptrac.h:2540
double bound_vmr_trend
Boundary conditions volume mixing ratio trend [ppv/s].
Definition: mptrac.h:2920
int met_cache
Preload meteo data into disk cache (0=no, 1=yes).
Definition: mptrac.h:2839
int qnt_oh
Quantity array index for OH volume mixing ratio (climatology).
Definition: mptrac.h:2516
int qnt_Ch
Quantity array index for H volume mixing ratio (chemistry code).
Definition: mptrac.h:2621
int met_press_level_def
Use predefined pressure levels or not.
Definition: mptrac.h:2790
int oh_chem_reaction
Reaction type for OH chemistry (0=none, 2=bimolecular, 3=termolecular).
Definition: mptrac.h:3058
int qnt_h2o
Quantity array index for water vapor volume mixing ratio.
Definition: mptrac.h:2465
int prof_ny
Number of latitudes of gridded profile data.
Definition: mptrac.h:3290
int qnt_rhice
Quantity array index for relative humidity over ice.
Definition: mptrac.h:2564
int qnt_rho
Quantity array index for density of air.
Definition: mptrac.h:2453
double sample_dz
Layer depth for sample output [km].
Definition: mptrac.h:3311
double tdec_strat
Life time of particles in the stratosphere [s].
Definition: mptrac.h:2956
int obs_type
Type of observation data files (0=ASCII, 1=netCDF).
Definition: mptrac.h:3155
double met_cms_eps_lwc
cmultiscale compression epsilon for cloud liquid water content.
Definition: mptrac.h:2748
int qnt_us
Quantity array index for surface zonal wind.
Definition: mptrac.h:2408
double met_cms_eps_z
cmultiscale compression epsilon for geopotential height.
Definition: mptrac.h:2724
double grid_lon1
Upper longitude of gridded data [deg].
Definition: mptrac.h:3251
int qnt_Cn2o
Quantity array index for N2O volume mixing ratio (chemistry code).
Definition: mptrac.h:2645
int qnt_Cccl3f
Quantity array index for CFC-11 volume mixing ratio (chemistry code).
Definition: mptrac.h:2639
double mixing_lat0
Lower latitude of mixing grid [deg].
Definition: mptrac.h:3025
int nens
Number of ensembles.
Definition: mptrac.h:3203
int qnt_pt
Quantity array index for tropopause pressure.
Definition: mptrac.h:2432
int qnt_cl
Quantity array index for total column cloud water.
Definition: mptrac.h:2492
int advect
Advection scheme (0=off, 1=Euler, 2=midpoint, 4=Runge-Kutta).
Definition: mptrac.h:2859
double prof_z1
Upper altitude of gridded profile data [km].
Definition: mptrac.h:3278
int qnt_t
Quantity array index for temperature.
Definition: mptrac.h:2450
int atm_filter
Time filter for atmospheric data output (0=none, 1=missval, 2=remove).
Definition: mptrac.h:3134
int kpp_chem
Switch for KPP chemistry module (0=off, 1=on).
Definition: mptrac.h:3070
int qnt_zeta
Quantity array index for zeta vertical coordinate.
Definition: mptrac.h:2570
double conv_pbl_trans
Depth of PBL transition layer (fraction of PBL depth).
Definition: mptrac.h:2899
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:2686
double csi_z0
Lower altitude of gridded CSI data [km].
Definition: mptrac.h:3179
int qnt_lapse
Quantity array index for lapse rate.
Definition: mptrac.h:2582
double stat_lat
Latitude of station [deg].
Definition: mptrac.h:3320
int qnt_Cho2
Quantity array index for HO2 volume mixing ratio (chemistry code).
Definition: mptrac.h:2624
int grid_ny
Number of latitudes of gridded data.
Definition: mptrac.h:3254
int qnt_Csf6
Quantity array index for SF6 volume mixing ratio (chemistry code).
Definition: mptrac.h:2648
int qnt_Ch2o
Quantity array index for H2O volume mixing ratio (chemistry code).
Definition: mptrac.h:2609
double met_detrend
FWHM of horizontal Gaussian used for detrending [km].
Definition: mptrac.h:2781
int conv_mix_pbl
Vertical mixing in the PBL (0=off, 1=on).
Definition: mptrac.h:2896
double bound_dps
Boundary conditions surface layer depth [hPa].
Definition: mptrac.h:2935
int dd_nbr_neighbours
Number of neighbours to communicate with.
Definition: mptrac.h:3356
double met_cms_eps_t
cmultiscale compression epsilon for temperature.
Definition: mptrac.h:2727
int chemgrid_nz
Number of altitudes of chemistry grid.
Definition: mptrac.h:3031
int qnt_cape
Quantity array index for convective available potential energy (CAPE).
Definition: mptrac.h:2504
int qnt_zeta_dot
Quantity array index forvelocity of zeta vertical coordinate.
Definition: mptrac.h:2576
double bound_mass_trend
Boundary conditions mass per particle trend [kg/s].
Definition: mptrac.h:2914
int mixing_nz
Number of altitudes of mixing grid.
Definition: mptrac.h:3004
int qnt_o3c
Quantity array index for total column ozone.
Definition: mptrac.h:2510
double bound_p0
Boundary conditions bottom pressure [hPa].
Definition: mptrac.h:2929
double mixing_lon0
Lower longitude of mixing grid [deg].
Definition: mptrac.h:3016
int qnt_Co3
Quantity array index for O3 volume mixing ratio (chemistry code).
Definition: mptrac.h:2612
int qnt_tsts
Quantity array index for T_STS.
Definition: mptrac.h:2600
int grid_nz
Number of altitudes of gridded data.
Definition: mptrac.h:3236
int qnt_nss
Quantity array index for northward turbulent surface stress.
Definition: mptrac.h:2417
double ens_dt_out
Time step for ensemble output [s].
Definition: mptrac.h:3209
int atm_stride
Particle index stride for atmospheric data files.
Definition: mptrac.h:3137
int met_relhum
Try to read relative humidity (0=no, 1=yes).
Definition: mptrac.h:2808
double mixing_lat1
Upper latitude of mixing grid [deg].
Definition: mptrac.h:3028
double atm_dt_out
Time step for atmospheric data output [s].
Definition: mptrac.h:3131
double prof_lat1
Upper latitude of gridded profile data [deg].
Definition: mptrac.h:3296
int met_cms_batch
cmultiscale batch size.
Definition: mptrac.h:2714
double psc_h2o
H2O volume mixing ratio for PSC analysis.
Definition: mptrac.h:3115
int met_sp
Smoothing for pressure levels.
Definition: mptrac.h:2778
double prof_lon0
Lower longitude of gridded profile data [deg].
Definition: mptrac.h:3284
int chemgrid_nx
Number of longitudes of chemistry grid.
Definition: mptrac.h:3040
int qnt_pct
Quantity array index for cloud top pressure.
Definition: mptrac.h:2486
int qnt_mloss_kpp
Quantity array index for total mass loss due to KPP chemistry.
Definition: mptrac.h:2534
int qnt_psat
Quantity array index for saturation pressure over water.
Definition: mptrac.h:2549
int qnt_subdomain
Quantity array index for current subdomain in domain decomposition.
Definition: mptrac.h:2654
double prof_lat0
Lower latitude of gridded profile data [deg].
Definition: mptrac.h:3293
int qnt_cin
Quantity array index for convective inhibition (CIN).
Definition: mptrac.h:2507
double psc_hno3
HNO3 volume mixing ratio for PSC analysis.
Definition: mptrac.h:3118
double prof_lon1
Upper longitude of gridded profile data [deg].
Definition: mptrac.h:3287
double met_cms_eps_rwc
cmultiscale compression epsilon for cloud rain water content.
Definition: mptrac.h:2751
int met_nc_quant
Number of digits for quantization of netCDF meteo files (0=off).
Definition: mptrac.h:2702
int h2o2_chem_reaction
Reaction type for H2O2 chemistry (0=none, 1=SO2).
Definition: mptrac.h:3067
int qnt_Co3p
Quantity array index for O(3P) volume mixing ratio (chemistry code).
Definition: mptrac.h:2633
double wet_depo_bc_ret_ratio
Coefficients for wet deposition below cloud: retention ratio.
Definition: mptrac.h:3106
int chemgrid_ny
Number of latitudes of chemistry grid.
Definition: mptrac.h:3049
double met_cms_eps_o3
cmultiscale compression epsilon for ozone.
Definition: mptrac.h:2745
int met_cms_zstd
cmultiscale zstd compression (0=off, 1=on).
Definition: mptrac.h:2717
int grid_sparse
Sparse output in grid data files (0=no, 1=yes).
Definition: mptrac.h:3224
double dry_depo_vdep
Dry deposition velocity [m/s].
Definition: mptrac.h:3112
int qnt_tt
Quantity array index for tropopause temperature.
Definition: mptrac.h:2435
int met_np
Number of target pressure levels.
Definition: mptrac.h:2784
int qnt_ens
Quantity array index for ensemble IDs.
Definition: mptrac.h:2381
int met_nc_level
zlib compression level of netCDF meteo files (0=off).
Definition: mptrac.h:2699
double mixing_dt
Time interval for mixing [s].
Definition: mptrac.h:2995
int qnt_mloss_h2o2
Quantity array index for total mass loss due to H2O2 chemistry.
Definition: mptrac.h:2531
double vtk_scale
Vertical scaling factor for VTK data.
Definition: mptrac.h:3341
double met_cms_eps_w
cmultiscale compression epsilon for vertical velocity.
Definition: mptrac.h:2736
double turb_dx_pbl
Horizontal turbulent diffusion coefficient (PBL) [m^2/s].
Definition: mptrac.h:2872
double conv_cin
CIN threshold for convection module [J/kg].
Definition: mptrac.h:2905
int qnt_pv
Quantity array index for potential vorticity.
Definition: mptrac.h:2591
int advect_vert_coord
Vertical velocity of air parcels (0=omega_on_plev, 1=zetadot_on_mlev, 2=omega_on_mlev).
Definition: mptrac.h:2863
int qnt_mloss_oh
Quantity array index for total mass loss due to OH chemistry.
Definition: mptrac.h:2528
int qnt_Ch2o2
Quantity array index for H2O2 volume mixing ratio (chemistry code).
Definition: mptrac.h:2627
int qnt_sst
Quantity array index for sea surface temperature.
Definition: mptrac.h:2426
double mixing_lon1
Upper longitude of mixing grid [deg].
Definition: mptrac.h:3019
int atm_nc_level
zlib compression level of netCDF atmospheric data files (0=off).
Definition: mptrac.h:3149
int met_cms_heur
cmultiscale coarsening heuristics (0=default, 1=mean diff, 2=median diff, 3=max diff).
Definition: mptrac.h:2721
double wet_depo_ic_ret_ratio
Coefficients for wet deposition in cloud: retention ratio.
Definition: mptrac.h:3103
int qnt_sh
Quantity array index for specific humidity.
Definition: mptrac.h:2558
int qnt_ess
Quantity array index for eastward turbulent surface stress.
Definition: mptrac.h:2414
double wet_depo_ic_b
Coefficient B for wet deposition in cloud (exponential form).
Definition: mptrac.h:3091
double wet_depo_bc_b
Coefficient B for wet deposition below cloud (exponential form).
Definition: mptrac.h:3085
int met_dy
Stride for latitudes.
Definition: mptrac.h:2766
int qnt_Cx
Quantity array index for trace species x volume mixing ratio (chemistry code).
Definition: mptrac.h:2606
double turb_dz_strat
Vertical turbulent diffusion coefficient (stratosphere) [m^2/s].
Definition: mptrac.h:2887
double bound_zetas
Boundary conditions surface layer zeta [K].
Definition: mptrac.h:2941
int dd_subdomains_zonal
Zonal subdomain number.
Definition: mptrac.h:3350
int qnt_idx
Quantity array index for air parcel IDs.
Definition: mptrac.h:2378
double met_tropo_theta
Dynamical tropopause potential temperature threshold [K].
Definition: mptrac.h:2830
int qnt_rwc
Quantity array index for cloud rain water content.
Definition: mptrac.h:2474
double t_start
Start time of simulation [s].
Definition: mptrac.h:2663
int nq
Number of quantities.
Definition: mptrac.h:2363
double tdec_trop
Life time of particles in the troposphere [s].
Definition: mptrac.h:2953
double sample_dx
Horizontal radius for sample output [km].
Definition: mptrac.h:3308
int vtk_stride
Particle index stride for VTK data.
Definition: mptrac.h:3338
double turb_dz_pbl
Vertical turbulent diffusion coefficient (PBL) [m^2/s].
Definition: mptrac.h:2881
double grid_lat1
Upper latitude of gridded data [deg].
Definition: mptrac.h:3260
int dd_subdomains_meridional
Meridional subdomain number.
Definition: mptrac.h:3353
int qnt_zt
Quantity array index for tropopause geopotential height.
Definition: mptrac.h:2438
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:2690
int qnt_cc
Quantity array index for cloud cover.
Definition: mptrac.h:2483
int qnt_plcl
Quantity array index for pressure at lifted condensation level (LCL).
Definition: mptrac.h:2495
double grid_dt_out
Time step for gridded data output [s].
Definition: mptrac.h:3221
int qnt_tdew
Quantity array index for dew point temperature.
Definition: mptrac.h:2594
Meteo data structure.
Definition: mptrac.h:3648
double subdomain_lon_max
Rectangular grid limit of subdomain.
Definition: mptrac.h:3825
double subdomain_lat_min
Rectangular grid limit of subdomain.
Definition: mptrac.h:3834
int nx
Number of longitudes.
Definition: mptrac.h:3654
double subdomain_lon_min
Rectangular grid limit of subdomain.
Definition: mptrac.h:3828
int ny_glob
Global sizes of meteo data.
Definition: mptrac.h:3858
int ny
Number of latitudes.
Definition: mptrac.h:3657
int halo_offset_start
Definition: mptrac.h:3849
int np
Number of pressure levels.
Definition: mptrac.h:3660
int halo_offset_end
Definition: mptrac.h:3852
int nx_glob
Global sizes of meteo data.
Definition: mptrac.h:3855
int npl
Number of model levels.
Definition: mptrac.h:3663
double time
Time [s].
Definition: mptrac.h:3651
double subdomain_lat_max
Rectangular grid limit of subdomain.
Definition: mptrac.h:3831
int np_glob
Global sizes of meteo data.
Definition: mptrac.h:3861
MPI information data.
Definition: mptrac.h:3426
int rank
Rank of node.
Definition: mptrac.h:3428
int size
Size of node.
Definition: mptrac.h:3431
Particle data.
Definition: mptrac.h:3400
double p
Pressure [hPa].
Definition: mptrac.h:3406
double lat
Latitude [deg].
Definition: mptrac.h:3412
double time
Time [s].
Definition: mptrac.h:3403
double lon
Longitude [deg].
Definition: mptrac.h:3409