53 {
54
56
58
60
62
63
65
66
67 if (argc < 4)
68 ERRMSG(
"Missing or invalid command-line arguments.\n\n"
69 "Usage: met_check_dt <ctl> <dt_file> <met> [KEY VALUE ...]\n\n"
70 "Use -h for full help.");
71
72
76
77
80 const double dx = 1e3 *
scan_ctl(argv[1], argc, argv,
"DX", -1,
"", NULL);
81 const double c_max =
scan_ctl(argv[1], argc, argv,
"CMAX", -1,
"0.5", NULL);
82 const double n_max =
scan_ctl(argv[1], argc, argv,
"NMAX", -1,
"0.3", NULL);
83
84
86
87
89 ERRMSG(
"Cannot open file!");
90
91
92 FILE *out;
93 LOG(1,
"Write time step data file: %s", argv[2]);
94 if (!(out = fopen(argv[2], "w")))
95 ERRMSG(
"Cannot create file!");
96
97
98 fprintf(out,
99 "# $1 = height [km]\n"
100 "# $2 = time step for horizontal advection [s]\n"
101 "# $3 = time step for vertical advection [s]\n"
102 "# $4 = time step for horizontal diffusion [s]\n"
103 "# $5 = time step for vertical diffusion [s]\n"
104 "# $6 = time step for PBL transition diffusion [s]\n"
105 "# $7 = time step for PBL depth diffusion [s]\n\n");
106
107
108 for (
int ip = 1; ip < met->
np - 1; ip++) {
109
110
111 double dt_x_min = 1e100;
112 double dt_p_min = 1e100;
113 double dt_dx_min = 1e100;
114 double dt_dp_min = 1e100;
115 double dt_pbl_min = 1e100;
116 double dt_pbl_depth_min = 1e100;
117
118
119#pragma omp parallel for default(shared) collapse(2) reduction(min:dt_x_min,dt_p_min,dt_dx_min,dt_dp_min,dt_pbl_min,dt_pbl_depth_min)
120 for (
int ix = 0; ix < met->
nx; ix++)
121 for (
int iy = 1; iy < met->
ny - 1; iy++) {
122
123
124 const double vh =
125 sqrt(
SQR(met->
u[ix][iy][ip]) +
SQR(met->
v[ix][iy][ip]));
126 const double dt_x = fabs(c_max * dx / vh);
127 if (vh != 0)
128 dt_x_min =
MIN(dt_x, dt_x_min);
129
130 const double dp = 0.5 * fabs(met->
p[ip + 1] - met->
p[ip - 1]);
131 const double dt_p = fabs(c_max * dp / met->
w[ix][iy][ip]);
132 if (met->
w[ix][iy][ip] != 0)
133 dt_p_min =
MIN(dt_p, dt_p_min);
134
135
139 const double wpbl =
140 pbl_weight_dt(met->
p[ip], met->
pbl[ix][iy], met->
ps[ix][iy],
141 pbl_trans);
142 const double wtrop = tropo_weight_dt(met->
p[ip], pt) * (1.0 - wpbl);
143 const double wstrat = 1.0 - wpbl - wtrop;
144
145
148 if (dx_loc > 0) {
149 const double dt_dx = 0.5 *
SQR(n_max * dx) / dx_loc;
150 dt_dx_min =
MIN(dt_dx, dt_dx_min);
151 }
152
155 if (dz_loc > 0) {
156 const double dt_dp = 0.5 *
SQR(n_max * dp)
157 / (
SQR(met->
p[ip] / (100. *
H0)) * dz_loc);
158 dt_dp_min =
MIN(dt_dp, dt_dp_min);
159 }
160
161
162 if (pbl_trans > 0 && met->
ps[ix][iy] > met->
pbl[ix][iy]
163 && met->
p[ip] <= met->
pbl[ix][iy]
164 && met->
p[ip] >= met->
pbl[ix][iy]
165 - pbl_trans * (met->
ps[ix][iy] - met->
pbl[ix][iy])) {
166 const double p0 = met->
pbl[ix][iy];
167 const double p1 = p0 - pbl_trans * (met->
ps[ix][iy] - p0);
168 const double dz_trans = 1e3 * fabs(
Z(p1) -
Z(p0));
170 if (dz_trans > 0 && dz_max > 0) {
171 const double dt_trans = 0.5 *
SQR(n_max * dz_trans) / dz_max;
172 dt_pbl_min =
MIN(dt_trans, dt_pbl_min);
173 }
174 }
175
176
177 if (met->
ps[ix][iy] > met->
pbl[ix][iy]
178 && met->
p[ip] >= met->
pbl[ix][iy]) {
179 const double dz_pbl =
180 1e3 * fabs(
Z(met->
pbl[ix][iy]) -
Z(met->
ps[ix][iy]));
182 const double dt_pbl_depth =
184 dt_pbl_depth_min =
MIN(dt_pbl_depth, dt_pbl_depth_min);
185 }
186 }
187 }
188
189
190 const double out_dt_dx = dt_dx_min < 1e99 ? dt_dx_min : NAN;
191 const double out_dt_dp = dt_dp_min < 1e99 ? dt_dp_min : NAN;
192 const double out_dt_pbl = dt_pbl_min < 1e99 ? dt_pbl_min : NAN;
193 const double out_dt_pbl_depth =
194 dt_pbl_depth_min < 1e99 ? dt_pbl_depth_min : NAN;
195 fprintf(out,
"%g %g %g %g %g %g %g\n",
Z(met->
p[ip]), dt_x_min, dt_p_min,
196 out_dt_dx, out_dt_dp, out_dt_pbl, out_dt_pbl_depth);
197 }
198
199
200 fclose(out);
201
202
203 free(clim);
204 free(met);
205 free(dd);
206
207 return EXIT_SUCCESS;
208}
double clim_tropo(const clim_t *clim, const double t, const double lat)
Calculates the tropopause pressure based on climatological data.
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.
void mptrac_read_clim(const ctl_t *ctl, clim_t *clim)
Reads various climatological data and populates the given climatology structure.
int mptrac_read_met(const char *filename, const ctl_t *ctl, const clim_t *clim, met_t *met, dd_t *dd)
Reads meteorological data from a file, supporting multiple formats and MPI broadcasting.
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.
#define H0
Scale height [km].
#define MIN(a, b)
Macro to determine the minimum of two values.
#define ERRMSG(...)
Print an error message with contextual information and terminate the program.
#define USAGE
Print usage information on -h or --help.
#define Z(p)
Convert pressure to altitude.
#define ALLOC(ptr, type, n)
Allocate memory for a pointer with error handling.
#define SQR(x)
Compute the square of a value.
#define LOG(level,...)
Print a log message with a specified logging level.
#define MAX(a, b)
Macro to determine the maximum of two values.
double met_utm_ref_lat
Reference latitude [deg] for UTM grid.
double turb_dz_trop
Vertical turbulent diffusion coefficient (troposphere) [m^2/s].
double turb_dx_strat
Horizontal turbulent diffusion coefficient (stratosphere) [m^2/s].
double turb_dx_trop
Horizontal turbulent diffusion coefficient (troposphere) [m^2/s].
double turb_pbl_trans
Depth of turbulent PBL transition layer (fraction of PBL depth).
double turb_dx_pbl
Horizontal turbulent diffusion coefficient (PBL) [m^2/s].
int met_coord_type
Type of coordinates for meteo data (-1=detect, 0=lat/lon [deg], 1=UTM [m]).
double turb_dz_strat
Vertical turbulent diffusion coefficient (stratosphere) [m^2/s].
double turb_dz_pbl
Vertical turbulent diffusion coefficient (PBL) [m^2/s].
Domain decomposition data structure.
float w[EX][EY][EP]
Vertical velocity [hPa/s].
int nx
Number of longitudes.
int ny
Number of latitudes.
float ps[EX][EY]
Surface pressure [hPa].
int np
Number of pressure levels.
float u[EX][EY][EP]
Zonal wind [m/s].
float pbl[EX][EY]
Boundary layer pressure [hPa].
float v[EX][EY][EP]
Meridional wind [m/s].
double lat[EY]
Latitudes [deg].
double p[EP]
Pressure levels [hPa].