miniWeather_serial.cpp 26 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636
  1. //////////////////////////////////////////////////////////////////////////////////////////
  2. // miniWeather
  3. // Author: Matt Norman <normanmr@ornl.gov> , Oak Ridge National Laboratory
  4. // This code simulates dry, stratified, compressible, non-hydrostatic fluid flows
  5. // For documentation, please see the attached documentation in the "documentation" folder
  6. //////////////////////////////////////////////////////////////////////////////////////////
  7. #include <stdlib.h>
  8. #include <math.h>
  9. #include <stdio.h>
  10. #include <nvtx3/nvToolsExt.h>
  11. const double pi = 3.14159265358979323846264338327; //Pi
  12. const double grav = 9.8; //Gravitational acceleration (m / s^2)
  13. const double cp = 1004.; //Specific heat of dry air at constant pressure
  14. const double rd = 287.; //Dry air constant for equation of state (P=rho*rd*T)
  15. const double p0 = 1.e5; //Standard pressure at the surface in Pascals
  16. const double C0 = 27.5629410929725921310572974482; //Constant to translate potential temperature into pressure (P=C0*(rho*theta)**gamma)
  17. const double gamm = 1.40027894002789400278940027894; //gamma=cp/Rd , have to call this gamm because "gamma" is taken (I hate C so much)
  18. //Define domain and stability-related constants
  19. const double xlen = 2.e4; //Length of the domain in the x-direction (meters)
  20. const double zlen = 1.e4; //Length of the domain in the z-direction (meters)
  21. const double hv_beta = 0.25; //How strong to diffuse the solution: hv_beta \in [0:1]
  22. const double cfl = 1.50; //"Courant, Friedrichs, Lewy" number (for numerical stability)
  23. const double max_speed = 450; //Assumed maximum wave speed during the simulation (speed of sound + speed of wind) (meter / sec)
  24. const int hs = 2; //"Halo" size: number of cells needed for a full "stencil" of information for reconstruction
  25. const int sten_size = 4; //Size of the stencil used for interpolation
  26. //Parameters for indexing and flags
  27. const int NUM_VARS = 4; //Number of fluid state variables
  28. const int ID_DENS = 0; //index for density ("rho")
  29. const int ID_UMOM = 1; //index for momentum in the x-direction ("rho * u")
  30. const int ID_WMOM = 2; //index for momentum in the z-direction ("rho * w")
  31. const int ID_RHOT = 3; //index for density * potential temperature ("rho * theta")
  32. const int DIR_X = 1; //Integer constant to express that this operation is in the x-direction
  33. const int DIR_Z = 2; //Integer constant to express that this operation is in the z-direction
  34. const int nqpoints = 3;
  35. double qpoints[] = {0.112701665379258311482073460022E0, 0.500000000000000000000000000000E0, 0.887298334620741688517926539980E0};
  36. double qweights[] = {0.277777777777777777777777777779E0, 0.444444444444444444444444444444E0, 0.277777777777777777777777777779E0};
  37. ///////////////////////////////////////////////////////////////////////////////////////
  38. // Variables that are initialized but remain static over the course of the simulation
  39. ///////////////////////////////////////////////////////////////////////////////////////
  40. double sim_time; //total simulation time in seconds
  41. double output_freq; //frequency to perform output in seconds
  42. double dt; //Model time step (seconds)
  43. int nx, nz; //Number of local grid cells in the x- and z- dimensions
  44. double dx, dz; //Grid space length in x- and z-dimension (meters)
  45. int nx_glob, nz_glob; //Number of total grid cells in the x- and z- dimensions
  46. int i_beg, k_beg; //beginning index in the x- and z-directions
  47. int nranks, myrank; //my rank id
  48. int left_rank, right_rank; //Rank IDs that exist to my left and right in the global domain
  49. double *hy_dens_cell; //hydrostatic density (vert cell avgs). Dimensions: (1-hs:nz+hs)
  50. double *hy_dens_theta_cell; //hydrostatic rho*t (vert cell avgs). Dimensions: (1-hs:nz+hs)
  51. double *hy_dens_int; //hydrostatic density (vert cell interf). Dimensions: (1:nz+1)
  52. double *hy_dens_theta_int; //hydrostatic rho*t (vert cell interf). Dimensions: (1:nz+1)
  53. double *hy_pressure_int; //hydrostatic press (vert cell interf). Dimensions: (1:nz+1)
  54. ///////////////////////////////////////////////////////////////////////////////////////
  55. // Variables that are dynamics over the course of the simulation
  56. ///////////////////////////////////////////////////////////////////////////////////////
  57. double etime; //Elapsed model time
  58. double output_counter; //Helps determine when it's time to do output
  59. //Runtime variable arrays
  60. double *state; //Fluid state. Dimensions: (1-hs:nx+hs,1-hs:nz+hs,NUM_VARS)
  61. double *state_tmp; //Fluid state. Dimensions: (1-hs:nx+hs,1-hs:nz+hs,NUM_VARS)
  62. double *flux; //Cell interface fluxes. Dimensions: (nx+1,nz+1,NUM_VARS)
  63. double *tend; //Fluid state tendencies. Dimensions: (nx,nz,NUM_VARS)
  64. int num_out = 0; //The number of outputs performed so far
  65. int direction_switch = 1;
  66. //How is this not in the standard?!
  67. double dmin(double a, double b)
  68. {
  69. if (a < b)
  70. {
  71. return a;
  72. }
  73. else
  74. {
  75. return b;
  76. }
  77. };
  78. //Declaring the functions defined after "main"
  79. void init();
  80. void finalize();
  81. void injection(double x, double z, double &r, double &u, double &w, double &t, double &hr, double &ht);
  82. void hydro_const_theta(double z, double &r, double &t);
  83. void output(double *state, double etime);
  84. void ncwrap(int ierr, int line);
  85. void perform_timestep(double *state, double *state_tmp, double *flux, double *tend, double dt);
  86. void semi_discrete_step(double *state_init, double *state_forcing, double *state_out, double dt, int dir, double *flux, double *tend);
  87. void compute_tendencies_x(double *state, double *flux, double *tend);
  88. void compute_tendencies_z(double *state, double *flux, double *tend);
  89. void set_halo_values_x(double *state);
  90. void set_halo_values_z(double *state);
  91. ///////////////////////////////////////////////////////////////////////////////////////
  92. // THE MAIN PROGRAM STARTS HERE
  93. ///////////////////////////////////////////////////////////////////////////////////////
  94. int main(int argc, char **argv)
  95. {
  96. ///////////////////////////////////////////////////////////////////////////////////////
  97. // BEGIN USER-CONFIGURABLE PARAMETERS
  98. ///////////////////////////////////////////////////////////////////////////////////////
  99. //The x-direction length is twice as long as the z-direction length
  100. //So, you'll want to have nx_glob be twice as large as nz_glob
  101. nx_glob = 40; //Number of total cells in the x-dirction
  102. nz_glob = 20; //Number of total cells in the z-dirction
  103. sim_time = 1000; //How many seconds to run the simulation
  104. output_freq = 100; //How frequently to output data to file (in seconds)
  105. ///////////////////////////////////////////////////////////////////////////////////////
  106. // END USER-CONFIGURABLE PARAMETERS
  107. ///////////////////////////////////////////////////////////////////////////////////////
  108. if (argc == 4)
  109. {
  110. printf("The arguments supplied are %s %s %s\n", argv[1], argv[2], argv[3]);
  111. nx_glob = atoi(argv[1]);
  112. nz_glob = atoi(argv[2]);
  113. sim_time = atoi(argv[3]);
  114. }
  115. else
  116. {
  117. printf("Using default values ...\n");
  118. }
  119. nvtxRangePushA("Total");
  120. init();
  121. //Output the initial state
  122. //output(state, etime);
  123. ////////////////////////////////////////////////////
  124. // MAIN TIME STEP LOOP
  125. ////////////////////////////////////////////////////
  126. nvtxRangePushA("while");
  127. while (etime < sim_time)
  128. {
  129. //If the time step leads to exceeding the simulation time, shorten it for the last step
  130. if (etime + dt > sim_time)
  131. {
  132. dt = sim_time - etime;
  133. }
  134. //Perform a single time step
  135. nvtxRangePushA("perform_timestep");
  136. perform_timestep(state, state_tmp, flux, tend, dt);
  137. nvtxRangePop();
  138. //Inform the user
  139. printf("Elapsed Time: %lf / %lf\n", etime, sim_time);
  140. //Update the elapsed time and output counter
  141. etime = etime + dt;
  142. output_counter = output_counter + dt;
  143. //If it's time for output, reset the counter, and do output
  144. if (output_counter >= output_freq)
  145. {
  146. output_counter = output_counter - output_freq;
  147. //output(state, etime);
  148. }
  149. }
  150. nvtxRangePop();
  151. finalize();
  152. nvtxRangePop();
  153. }
  154. //Performs a single dimensionally split time step using a simple low-storate three-stage Runge-Kutta time integrator
  155. //The dimensional splitting is a second-order-accurate alternating Strang splitting in which the
  156. //order of directions is alternated each time step.
  157. //The Runge-Kutta method used here is defined as follows:
  158. // q* = q[n] + dt/3 * rhs(q[n])
  159. // q** = q[n] + dt/2 * rhs(q* )
  160. // q[n+1] = q[n] + dt/1 * rhs(q** )
  161. void perform_timestep(double *state, double *state_tmp, double *flux, double *tend, double dt)
  162. {
  163. if (direction_switch)
  164. {
  165. //x-direction first
  166. semi_discrete_step(state, state, state_tmp, dt / 3, DIR_X, flux, tend);
  167. semi_discrete_step(state, state_tmp, state_tmp, dt / 2, DIR_X, flux, tend);
  168. semi_discrete_step(state, state_tmp, state, dt / 1, DIR_X, flux, tend);
  169. //z-direction second
  170. semi_discrete_step(state, state, state_tmp, dt / 3, DIR_Z, flux, tend);
  171. semi_discrete_step(state, state_tmp, state_tmp, dt / 2, DIR_Z, flux, tend);
  172. semi_discrete_step(state, state_tmp, state, dt / 1, DIR_Z, flux, tend);
  173. }
  174. else
  175. {
  176. //z-direction second
  177. semi_discrete_step(state, state, state_tmp, dt / 3, DIR_Z, flux, tend);
  178. semi_discrete_step(state, state_tmp, state_tmp, dt / 2, DIR_Z, flux, tend);
  179. semi_discrete_step(state, state_tmp, state, dt / 1, DIR_Z, flux, tend);
  180. //x-direction first
  181. semi_discrete_step(state, state, state_tmp, dt / 3, DIR_X, flux, tend);
  182. semi_discrete_step(state, state_tmp, state_tmp, dt / 2, DIR_X, flux, tend);
  183. semi_discrete_step(state, state_tmp, state, dt / 1, DIR_X, flux, tend);
  184. }
  185. if (direction_switch)
  186. {
  187. direction_switch = 0;
  188. }
  189. else
  190. {
  191. direction_switch = 1;
  192. }
  193. }
  194. //Perform a single semi-discretized step in time with the form:
  195. //state_out = state_init + dt * rhs(state_forcing)
  196. //Meaning the step starts from state_init, computes the rhs using state_forcing, and stores the result in state_out
  197. void semi_discrete_step(double *state_init, double *state_forcing, double *state_out, double dt, int dir, double *flux, double *tend)
  198. {
  199. int i, k, ll, inds, indt;
  200. if (dir == DIR_X)
  201. {
  202. //Set the halo values in the x-direction
  203. set_halo_values_x(state_forcing);
  204. //Compute the time tendencies for the fluid state in the x-direction
  205. compute_tendencies_x(state_forcing, flux, tend);
  206. }
  207. else if (dir == DIR_Z)
  208. {
  209. //Set the halo values in the z-direction
  210. set_halo_values_z(state_forcing);
  211. //Compute the time tendencies for the fluid state in the z-direction
  212. compute_tendencies_z(state_forcing, flux, tend);
  213. }
  214. /////////////////////////////////////////////////
  215. // TODO: THREAD ME
  216. /////////////////////////////////////////////////
  217. //Apply the tendencies to the fluid state
  218. for (ll = 0; ll < NUM_VARS; ll++)
  219. {
  220. for (k = 0; k < nz; k++)
  221. {
  222. for (i = 0; i < nx; i++)
  223. {
  224. inds = ll * (nz + 2 * hs) * (nx + 2 * hs) + (k + hs) * (nx + 2 * hs) + i + hs;
  225. indt = ll * nz * nx + k * nx + i;
  226. state_out[inds] = state_init[inds] + dt * tend[indt];
  227. }
  228. }
  229. }
  230. }
  231. //Compute the time tendencies of the fluid state using forcing in the x-direction
  232. //First, compute the flux vector at each cell interface in the x-direction (including hyperviscosity)
  233. //Then, compute the tendencies using those fluxes
  234. void compute_tendencies_x(double *state, double *flux, double *tend)
  235. {
  236. int i, k, ll, s, inds, indf1, indf2, indt;
  237. double r, u, w, t, p, stencil[4], d3_vals[NUM_VARS], vals[NUM_VARS], hv_coef;
  238. //Compute the hyperviscosity coeficient
  239. hv_coef = -hv_beta * dx / (16 * dt);
  240. /////////////////////////////////////////////////
  241. // TODO: THREAD ME
  242. /////////////////////////////////////////////////
  243. //Compute fluxes in the x-direction for each cell
  244. for (k = 0; k < nz; k++)
  245. {
  246. for (i = 0; i < nx + 1; i++)
  247. {
  248. //Use fourth-order interpolation from four cell averages to compute the value at the interface in question
  249. for (ll = 0; ll < NUM_VARS; ll++)
  250. {
  251. for (s = 0; s < sten_size; s++)
  252. {
  253. inds = ll * (nz + 2 * hs) * (nx + 2 * hs) + (k + hs) * (nx + 2 * hs) + i + s;
  254. stencil[s] = state[inds];
  255. }
  256. //Fourth-order-accurate interpolation of the state
  257. vals[ll] = -stencil[0] / 12 + 7 * stencil[1] / 12 + 7 * stencil[2] / 12 - stencil[3] / 12;
  258. //First-order-accurate interpolation of the third spatial derivative of the state (for artificial viscosity)
  259. d3_vals[ll] = -stencil[0] + 3 * stencil[1] - 3 * stencil[2] + stencil[3];
  260. }
  261. //Compute density, u-wind, w-wind, potential temperature, and pressure (r,u,w,t,p respectively)
  262. r = vals[ID_DENS] + hy_dens_cell[k + hs];
  263. u = vals[ID_UMOM] / r;
  264. w = vals[ID_WMOM] / r;
  265. t = (vals[ID_RHOT] + hy_dens_theta_cell[k + hs]) / r;
  266. p = C0 * pow((r * t), gamm);
  267. //Compute the flux vector
  268. flux[ID_DENS * (nz + 1) * (nx + 1) + k * (nx + 1) + i] = r * u - hv_coef * d3_vals[ID_DENS];
  269. flux[ID_UMOM * (nz + 1) * (nx + 1) + k * (nx + 1) + i] = r * u * u + p - hv_coef * d3_vals[ID_UMOM];
  270. flux[ID_WMOM * (nz + 1) * (nx + 1) + k * (nx + 1) + i] = r * u * w - hv_coef * d3_vals[ID_WMOM];
  271. flux[ID_RHOT * (nz + 1) * (nx + 1) + k * (nx + 1) + i] = r * u * t - hv_coef * d3_vals[ID_RHOT];
  272. }
  273. }
  274. /////////////////////////////////////////////////
  275. // TODO: THREAD ME
  276. /////////////////////////////////////////////////
  277. //Use the fluxes to compute tendencies for each cell
  278. for (ll = 0; ll < NUM_VARS; ll++)
  279. {
  280. for (k = 0; k < nz; k++)
  281. {
  282. for (i = 0; i < nx; i++)
  283. {
  284. indt = ll * nz * nx + k * nx + i;
  285. indf1 = ll * (nz + 1) * (nx + 1) + k * (nx + 1) + i;
  286. indf2 = ll * (nz + 1) * (nx + 1) + k * (nx + 1) + i + 1;
  287. tend[indt] = -(flux[indf2] - flux[indf1]) / dx;
  288. }
  289. }
  290. }
  291. }
  292. //Compute the time tendencies of the fluid state using forcing in the z-direction
  293. //First, compute the flux vector at each cell interface in the z-direction (including hyperviscosity)
  294. //Then, compute the tendencies using those fluxes
  295. void compute_tendencies_z(double *state, double *flux, double *tend)
  296. {
  297. int i, k, ll, s, inds, indf1, indf2, indt;
  298. double r, u, w, t, p, stencil[4], d3_vals[NUM_VARS], vals[NUM_VARS], hv_coef;
  299. //Compute the hyperviscosity coeficient
  300. hv_coef = -hv_beta * dx / (16 * dt);
  301. /////////////////////////////////////////////////
  302. // TODO: THREAD ME
  303. /////////////////////////////////////////////////
  304. //Compute fluxes in the x-direction for each cell
  305. for (k = 0; k < nz + 1; k++)
  306. {
  307. for (i = 0; i < nx; i++)
  308. {
  309. //Use fourth-order interpolation from four cell averages to compute the value at the interface in question
  310. for (ll = 0; ll < NUM_VARS; ll++)
  311. {
  312. for (s = 0; s < sten_size; s++)
  313. {
  314. inds = ll * (nz + 2 * hs) * (nx + 2 * hs) + (k + s) * (nx + 2 * hs) + i + hs;
  315. stencil[s] = state[inds];
  316. }
  317. //Fourth-order-accurate interpolation of the state
  318. vals[ll] = -stencil[0] / 12 + 7 * stencil[1] / 12 + 7 * stencil[2] / 12 - stencil[3] / 12;
  319. //First-order-accurate interpolation of the third spatial derivative of the state
  320. d3_vals[ll] = -stencil[0] + 3 * stencil[1] - 3 * stencil[2] + stencil[3];
  321. }
  322. //Compute density, u-wind, w-wind, potential temperature, and pressure (r,u,w,t,p respectively)
  323. r = vals[ID_DENS] + hy_dens_int[k];
  324. u = vals[ID_UMOM] / r;
  325. w = vals[ID_WMOM] / r;
  326. t = (vals[ID_RHOT] + hy_dens_theta_int[k]) / r;
  327. p = C0 * pow((r * t), gamm) - hy_pressure_int[k];
  328. //Compute the flux vector with hyperviscosity
  329. flux[ID_DENS * (nz + 1) * (nx + 1) + k * (nx + 1) + i] = r * w - hv_coef * d3_vals[ID_DENS];
  330. flux[ID_UMOM * (nz + 1) * (nx + 1) + k * (nx + 1) + i] = r * w * u - hv_coef * d3_vals[ID_UMOM];
  331. flux[ID_WMOM * (nz + 1) * (nx + 1) + k * (nx + 1) + i] = r * w * w + p - hv_coef * d3_vals[ID_WMOM];
  332. flux[ID_RHOT * (nz + 1) * (nx + 1) + k * (nx + 1) + i] = r * w * t - hv_coef * d3_vals[ID_RHOT];
  333. }
  334. }
  335. /////////////////////////////////////////////////
  336. // TODO: THREAD ME
  337. /////////////////////////////////////////////////
  338. //Use the fluxes to compute tendencies for each cell
  339. for (ll = 0; ll < NUM_VARS; ll++)
  340. {
  341. for (k = 0; k < nz; k++)
  342. {
  343. for (i = 0; i < nx; i++)
  344. {
  345. indt = ll * nz * nx + k * nx + i;
  346. indf1 = ll * (nz + 1) * (nx + 1) + (k) * (nx + 1) + i;
  347. indf2 = ll * (nz + 1) * (nx + 1) + (k + 1) * (nx + 1) + i;
  348. tend[indt] = -(flux[indf2] - flux[indf1]) / dz;
  349. if (ll == ID_WMOM)
  350. {
  351. inds = ID_DENS * (nz + 2 * hs) * (nx + 2 * hs) + (k + hs) * (nx + 2 * hs) + i + hs;
  352. tend[indt] = tend[indt] - state[inds] * grav;
  353. }
  354. }
  355. }
  356. }
  357. }
  358. void set_halo_values_x(double *state)
  359. {
  360. int k, ll, ind_r, ind_u, ind_t, i;
  361. double z;
  362. for (ll = 0; ll < NUM_VARS; ll++)
  363. {
  364. for (k = 0; k < nz; k++)
  365. {
  366. state[ll * (nz + 2 * hs) * (nx + 2 * hs) + (k + hs) * (nx + 2 * hs) + 0] = state[ll * (nz + 2 * hs) * (nx + 2 * hs) + (k + hs) * (nx + 2 * hs) + nx + hs - 2];
  367. state[ll * (nz + 2 * hs) * (nx + 2 * hs) + (k + hs) * (nx + 2 * hs) + 1] = state[ll * (nz + 2 * hs) * (nx + 2 * hs) + (k + hs) * (nx + 2 * hs) + nx + hs - 1];
  368. state[ll * (nz + 2 * hs) * (nx + 2 * hs) + (k + hs) * (nx + 2 * hs) + nx + hs] = state[ll * (nz + 2 * hs) * (nx + 2 * hs) + (k + hs) * (nx + 2 * hs) + hs];
  369. state[ll * (nz + 2 * hs) * (nx + 2 * hs) + (k + hs) * (nx + 2 * hs) + nx + hs + 1] = state[ll * (nz + 2 * hs) * (nx + 2 * hs) + (k + hs) * (nx + 2 * hs) + hs + 1];
  370. }
  371. }
  372. ////////////////////////////////////////////////////
  373. if (myrank == 0)
  374. {
  375. for (k = 0; k < nz; k++)
  376. {
  377. for (i = 0; i < hs; i++)
  378. {
  379. z = (k_beg + k + 0.5) * dz;
  380. if (abs(z - 3 * zlen / 4) <= zlen / 16)
  381. {
  382. ind_r = ID_DENS * (nz + 2 * hs) * (nx + 2 * hs) + (k + hs) * (nx + 2 * hs) + i;
  383. ind_u = ID_UMOM * (nz + 2 * hs) * (nx + 2 * hs) + (k + hs) * (nx + 2 * hs) + i;
  384. ind_t = ID_RHOT * (nz + 2 * hs) * (nx + 2 * hs) + (k + hs) * (nx + 2 * hs) + i;
  385. state[ind_u] = (state[ind_r] + hy_dens_cell[k + hs]) * 50.;
  386. state[ind_t] = (state[ind_r] + hy_dens_cell[k + hs]) * 298. - hy_dens_theta_cell[k + hs];
  387. }
  388. }
  389. }
  390. }
  391. }
  392. //Set this task's halo values in the z-direction.
  393. //decomposition in the vertical direction.
  394. void set_halo_values_z(double *state)
  395. {
  396. int i, ll;
  397. const double mnt_width = xlen / 8;
  398. double x, xloc, mnt_deriv;
  399. /////////////////////////////////////////////////
  400. // TODO: THREAD ME
  401. /////////////////////////////////////////////////
  402. for (ll = 0; ll < NUM_VARS; ll++)
  403. {
  404. for (i = 0; i < nx + 2 * hs; i++)
  405. {
  406. if (ll == ID_WMOM)
  407. {
  408. state[ll * (nz + 2 * hs) * (nx + 2 * hs) + (0) * (nx + 2 * hs) + i] = 0.;
  409. state[ll * (nz + 2 * hs) * (nx + 2 * hs) + (1) * (nx + 2 * hs) + i] = 0.;
  410. state[ll * (nz + 2 * hs) * (nx + 2 * hs) + (nz + hs) * (nx + 2 * hs) + i] = 0.;
  411. state[ll * (nz + 2 * hs) * (nx + 2 * hs) + (nz + hs + 1) * (nx + 2 * hs) + i] = 0.;
  412. }
  413. else
  414. {
  415. state[ll * (nz + 2 * hs) * (nx + 2 * hs) + (0) * (nx + 2 * hs) + i] = state[ll * (nz + 2 * hs) * (nx + 2 * hs) + (hs) * (nx + 2 * hs) + i];
  416. state[ll * (nz + 2 * hs) * (nx + 2 * hs) + (1) * (nx + 2 * hs) + i] = state[ll * (nz + 2 * hs) * (nx + 2 * hs) + (hs) * (nx + 2 * hs) + i];
  417. state[ll * (nz + 2 * hs) * (nx + 2 * hs) + (nz + hs) * (nx + 2 * hs) + i] = state[ll * (nz + 2 * hs) * (nx + 2 * hs) + (nz + hs - 1) * (nx + 2 * hs) + i];
  418. state[ll * (nz + 2 * hs) * (nx + 2 * hs) + (nz + hs + 1) * (nx + 2 * hs) + i] = state[ll * (nz + 2 * hs) * (nx + 2 * hs) + (nz + hs - 1) * (nx + 2 * hs) + i];
  419. }
  420. }
  421. }
  422. }
  423. void init()
  424. {
  425. int i, k, ii, kk, ll, inds, i_end;
  426. double x, z, r, u, w, t, hr, ht, nper;
  427. //Set the cell grid size
  428. dx = xlen / nx_glob;
  429. dz = zlen / nz_glob;
  430. nranks = 1;
  431. myrank = 0;
  432. // For simpler version, replace i_beg = 0, nx = nx_glob, left_rank = 0, right_rank = 0;
  433. nper = ((double)nx_glob) / nranks;
  434. i_beg = round(nper * (myrank));
  435. i_end = round(nper * ((myrank) + 1)) - 1;
  436. nx = i_end - i_beg + 1;
  437. left_rank = myrank - 1;
  438. if (left_rank == -1)
  439. left_rank = nranks - 1;
  440. right_rank = myrank + 1;
  441. if (right_rank == nranks)
  442. right_rank = 0;
  443. ////////////////////////////////////////////////////////////////////////////////
  444. ////////////////////////////////////////////////////////////////////////////////
  445. // YOU DON'T NEED TO ALTER ANYTHING BELOW THIS POINT IN THE CODE
  446. ////////////////////////////////////////////////////////////////////////////////
  447. ////////////////////////////////////////////////////////////////////////////////
  448. k_beg = 0;
  449. nz = nz_glob;
  450. //Allocate the model data
  451. state = (double *)malloc((nx + 2 * hs) * (nz + 2 * hs) * NUM_VARS * sizeof(double));
  452. state_tmp = (double *)malloc((nx + 2 * hs) * (nz + 2 * hs) * NUM_VARS * sizeof(double));
  453. flux = (double *)malloc((nx + 1) * (nz + 1) * NUM_VARS * sizeof(double));
  454. tend = (double *)malloc(nx * nz * NUM_VARS * sizeof(double));
  455. hy_dens_cell = (double *)malloc((nz + 2 * hs) * sizeof(double));
  456. hy_dens_theta_cell = (double *)malloc((nz + 2 * hs) * sizeof(double));
  457. hy_dens_int = (double *)malloc((nz + 1) * sizeof(double));
  458. hy_dens_theta_int = (double *)malloc((nz + 1) * sizeof(double));
  459. hy_pressure_int = (double *)malloc((nz + 1) * sizeof(double));
  460. //Define the maximum stable time step based on an assumed maximum wind speed
  461. dt = dmin(dx, dz) / max_speed * cfl;
  462. //Set initial elapsed model time and output_counter to zero
  463. etime = 0.;
  464. output_counter = 0.;
  465. // Display grid information
  466. printf("nx_glob, nz_glob: %d %d\n", nx_glob, nz_glob);
  467. printf("dx,dz: %lf %lf\n", dx, dz);
  468. printf("dt: %lf\n", dt);
  469. //////////////////////////////////////////////////////////////////////////
  470. // Initialize the cell-averaged fluid state via Gauss-Legendre quadrature
  471. //////////////////////////////////////////////////////////////////////////
  472. for (k = 0; k < nz + 2 * hs; k++)
  473. {
  474. for (i = 0; i < nx + 2 * hs; i++)
  475. {
  476. //Initialize the state to zero
  477. for (ll = 0; ll < NUM_VARS; ll++)
  478. {
  479. inds = ll * (nz + 2 * hs) * (nx + 2 * hs) + k * (nx + 2 * hs) + i;
  480. state[inds] = 0.;
  481. }
  482. //Use Gauss-Legendre quadrature to initialize a hydrostatic balance + temperature perturbation
  483. for (kk = 0; kk < nqpoints; kk++)
  484. {
  485. for (ii = 0; ii < nqpoints; ii++)
  486. {
  487. //Compute the x,z location within the global domain based on cell and quadrature index
  488. x = (i_beg + i - hs + 0.5) * dx + (qpoints[ii] - 0.5) * dx;
  489. z = (k_beg + k - hs + 0.5) * dz + (qpoints[kk] - 0.5) * dz;
  490. //Set the fluid state based on the user's specification (default is injection in this example)
  491. injection(x, z, r, u, w, t, hr, ht);
  492. //Store into the fluid state array
  493. inds = ID_DENS * (nz + 2 * hs) * (nx + 2 * hs) + k * (nx + 2 * hs) + i;
  494. state[inds] = state[inds] + r * qweights[ii] * qweights[kk];
  495. inds = ID_UMOM * (nz + 2 * hs) * (nx + 2 * hs) + k * (nx + 2 * hs) + i;
  496. state[inds] = state[inds] + (r + hr) * u * qweights[ii] * qweights[kk];
  497. inds = ID_WMOM * (nz + 2 * hs) * (nx + 2 * hs) + k * (nx + 2 * hs) + i;
  498. state[inds] = state[inds] + (r + hr) * w * qweights[ii] * qweights[kk];
  499. inds = ID_RHOT * (nz + 2 * hs) * (nx + 2 * hs) + k * (nx + 2 * hs) + i;
  500. state[inds] = state[inds] + ((r + hr) * (t + ht) - hr * ht) * qweights[ii] * qweights[kk];
  501. }
  502. }
  503. for (ll = 0; ll < NUM_VARS; ll++)
  504. {
  505. inds = ll * (nz + 2 * hs) * (nx + 2 * hs) + k * (nx + 2 * hs) + i;
  506. state_tmp[inds] = state[inds];
  507. }
  508. }
  509. }
  510. //Compute the hydrostatic background state over vertical cell averages
  511. for (k = 0; k < nz + 2 * hs; k++)
  512. {
  513. hy_dens_cell[k] = 0.;
  514. hy_dens_theta_cell[k] = 0.;
  515. for (kk = 0; kk < nqpoints; kk++)
  516. {
  517. z = (k_beg + k - hs + 0.5) * dz;
  518. //Set the fluid state based on the user's specification (default is injection in this example)
  519. injection(0., z, r, u, w, t, hr, ht);
  520. hy_dens_cell[k] = hy_dens_cell[k] + hr * qweights[kk];
  521. hy_dens_theta_cell[k] = hy_dens_theta_cell[k] + hr * ht * qweights[kk];
  522. }
  523. }
  524. //Compute the hydrostatic background state at vertical cell interfaces
  525. for (k = 0; k < nz + 1; k++)
  526. {
  527. z = (k_beg + k) * dz;
  528. //Set the fluid state based on the user's specification (default is injection in this example)
  529. injection(0., z, r, u, w, t, hr, ht);
  530. hy_dens_int[k] = hr;
  531. hy_dens_theta_int[k] = hr * ht;
  532. hy_pressure_int[k] = C0 * pow((hr * ht), gamm);
  533. }
  534. }
  535. //This test case is initially balanced but injects fast, cold air from the left boundary near the model top
  536. //x and z are input coordinates at which to sample
  537. //r,u,w,t are output density, u-wind, w-wind, and potential temperature at that location
  538. //hr and ht are output background hydrostatic density and potential temperature at that location
  539. void injection(double x, double z, double &r, double &u, double &w, double &t, double &hr, double &ht)
  540. {
  541. hydro_const_theta(z, hr, ht);
  542. r = 0.;
  543. t = 0.;
  544. u = 0.;
  545. w = 0.;
  546. }
  547. //Establish hydrstatic balance using constant potential temperature (thermally neutral atmosphere)
  548. //z is the input coordinate
  549. //r and t are the output background hydrostatic density and potential temperature
  550. void hydro_const_theta(double z, double &r, double &t)
  551. {
  552. const double theta0 = 300.; //Background potential temperature
  553. const double exner0 = 1.; //Surface-level Exner pressure
  554. double p, exner, rt;
  555. //Establish hydrostatic balance first using Exner pressure
  556. t = theta0; //Potential Temperature at z
  557. exner = exner0 - grav * z / (cp * theta0); //Exner pressure at z
  558. p = p0 * pow(exner, (cp / rd)); //Pressure at z
  559. rt = pow((p / C0), (1. / gamm)); //rho*theta at z
  560. r = rt / t; //Density at z
  561. }
  562. void finalize()
  563. {
  564. free(state);
  565. free(state_tmp);
  566. free(flux);
  567. free(tend);
  568. free(hy_dens_cell);
  569. free(hy_dens_theta_cell);
  570. free(hy_dens_int);
  571. free(hy_dens_theta_int);
  572. free(hy_pressure_int);
  573. }