miniWeather_serial.cpp 30 KB

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