miniWeather_openacc.f90 26 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576
  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. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  8. !! Copyright (c) 2018, National Center for Computational Sciences, Oak Ridge National Laboratory. All rights reserved.
  9. !!
  10. !! Portions Copyright (c) 2020, NVIDIA Corporation. All rights reserved.
  11. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  12. program miniweather
  13. use nvtx
  14. !implicit none
  15. !Declare the precision for the real constants (at least 15 digits of accuracy = double precision)
  16. integer , parameter :: rp = selected_real_kind(15)
  17. !Define some physical constants to use throughout the simulation
  18. real(rp), parameter :: pi = 3.14159265358979323846264338327_rp !Pi
  19. real(rp), parameter :: grav = 9.8_rp !Gravitational acceleration (m / s^2)
  20. real(rp), parameter :: cp = 1004._rp !Specific heat of dry air at constant pressure
  21. real(rp), parameter :: rd = 287._rp !Dry air constant for equation of state (P=rho*rd*T)
  22. real(rp), parameter :: p0 = 1.e5_rp !Standard pressure at the surface in Pascals
  23. real(rp), parameter :: C0 = 27.5629410929725921310572974482_rp !Constant to translate potential temperature into pressure (P=C0*(rho*theta)**gamma)
  24. real(rp), parameter :: gamma = 1.40027894002789400278940027894_rp !gamma=cp/Rd
  25. !Define domain and stability-related constants
  26. real(rp), parameter :: xlen = 2.e4_rp !Length of the domain in the x-direction (meters)
  27. real(rp), parameter :: zlen = 1.e4_rp !Length of the domain in the z-direction (meters)
  28. real(rp), parameter :: hv_beta = 0.25_rp !How strong to diffuse the solution: hv_beta \in [0:1]
  29. real(rp), parameter :: cfl = 1.50_rp !"Courant, Friedrichs, Lewy" number (for numerical stability)
  30. real(rp), parameter :: max_speed = 450 !Assumed maximum wave speed during the simulation (speed of sound + speed of wind) (meter / sec)
  31. integer , parameter :: hs = 2 !"Halo" size: number of cells needed for a full "stencil" of information for reconstruction
  32. integer , parameter :: sten_size = 4 !Size of the stencil used for interpolation
  33. !Parameters for indexing and flags
  34. integer , parameter :: NUM_VARS = 4 !Number of fluid state variables
  35. integer , parameter :: ID_DENS = 1 !index for density ("rho")
  36. integer , parameter :: ID_UMOM = 2 !index for momentum in the x-direction ("rho * u")
  37. integer , parameter :: ID_WMOM = 3 !index for momentum in the z-direction ("rho * w")
  38. integer , parameter :: ID_RHOT = 4 !index for density * potential temperature ("rho * theta")
  39. integer , parameter :: DIR_X = 1 !Integer constant to express that this operation is in the x-direction
  40. integer , parameter :: DIR_Z = 2 !Integer constant to express that this operation is in the z-direction
  41. !Gauss-Legendre quadrature points and weights on the domain [0:1]
  42. integer , parameter :: nqpoints = 3
  43. real(rp) :: qpoints (nqpoints) = (/ 0.112701665379258311482073460022E0_rp , 0.500000000000000000000000000000E0_rp , 0.887298334620741688517926539980E0_rp /)
  44. real(rp) :: qweights(nqpoints) = (/ 0.277777777777777777777777777779E0_rp , 0.444444444444444444444444444444E0_rp , 0.277777777777777777777777777779E0_rp /)
  45. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  46. !! Variables that are initialized but remain static over the course of the simulation
  47. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  48. real(rp) :: sim_time !total simulation time in seconds
  49. real(rp) :: output_freq !frequency to perform output in seconds
  50. real(rp) :: dt !Model time step (seconds)
  51. integer :: nx, nz !Number of local grid cells in the x- and z- dimensions
  52. real(rp) :: dx, dz !Grid space length in x- and z-dimension (meters)
  53. integer :: nx_glob, nz_glob !Number of total grid cells in the x- and z- dimensions
  54. integer :: i_beg, k_beg !beginning index in the x- and z-directions
  55. integer :: nranks, myrank !my rank id
  56. integer :: left_rank, right_rank
  57. real(rp), allocatable :: hy_dens_cell (:) !hydrostatic density (vert cell avgs). Dimensions: (1-hs:nz+hs)
  58. real(rp), allocatable :: hy_dens_theta_cell(:) !hydrostatic rho*t (vert cell avgs). Dimensions: (1-hs:nz+hs)
  59. real(rp), allocatable :: hy_dens_int (:) !hydrostatic density (vert cell interf). Dimensions: (1:nz+1)
  60. real(rp), allocatable :: hy_dens_theta_int (:) !hydrostatic rho*t (vert cell interf). Dimensions: (1:nz+1)
  61. real(rp), allocatable :: hy_pressure_int (:) !hydrostatic press (vert cell interf). Dimensions: (1:nz+1)
  62. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  63. !! Variables that are dynamics over the course of the simulation
  64. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  65. real(rp) :: etime !Elapsed model time
  66. real(rp) :: output_counter !Helps determine when it's time to do output
  67. !Runtime variable arrays
  68. real(rp), allocatable :: state (:,:,:) !Fluid state. Dimensions: (1-hs:nx+hs,1-hs:nz+hs,NUM_VARS)
  69. real(rp), allocatable :: state_tmp (:,:,:) !Fluid state. Dimensions: (1-hs:nx+hs,1-hs:nz+hs,NUM_VARS)
  70. real(rp), allocatable :: flux (:,:,:) !Cell interface fluxes. Dimensions: (nx+1,nz+1,NUM_VARS)
  71. real(rp), allocatable :: tend (:,:,:) !Fluid state tendencies. Dimensions: (nx,nz,NUM_VARS)
  72. integer(8) :: t1, t2, rate !For CPU Timings
  73. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  74. !! THE MAIN PROGRAM STARTS HERE
  75. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  76. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  77. !! BEGIN USER-CONFIGURABLE PARAMETERS
  78. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  79. !The x-direction length is twice as long as the z-direction length
  80. !So, you'll want to have nx_glob be twice as large as nz_glob
  81. nx_glob = 40 !Number of total cells in the x-dirction
  82. nz_glob = 20 !Number of total cells in the z-dirction
  83. sim_time = 10 !How many seconds to run the simulation
  84. output_freq = 100 !How frequently to output data to file (in seconds)
  85. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  86. !! END USER-CONFIGURABLE PARAMETERS
  87. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  88. !Allocate arrays, initialize the grid and the data
  89. call nvtxStartRange("Total")
  90. call init()
  91. !$acc data copyin(state_tmp,hy_dens_cell,hy_dens_theta_cell,hy_dens_int,hy_dens_theta_int,hy_pressure_int) create(flux,tend) copy(state)
  92. !Output the initial state
  93. !call output(state,etime)
  94. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  95. !! MAIN TIME STEP LOOP
  96. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  97. call system_clock(t1)
  98. call nvtxStartRange("while")
  99. do while (etime < sim_time)
  100. !If the time step leads to exceeding the simulation time, shorten it for the last step
  101. if (etime + dt > sim_time) dt = sim_time - etime
  102. !Perform a single time step
  103. call nvtxStartRange("perform_timestep")
  104. call perform_timestep(state,state_tmp,flux,tend,dt)
  105. call nvtxEndRange
  106. !Inform the user
  107. write(*,*) 'Elapsed Time: ', etime , ' / ' , sim_time
  108. !Update the elapsed time and output counter
  109. etime = etime + dt
  110. output_counter = output_counter + dt
  111. !If it's time for output, reset the counter, and do output
  112. if (output_counter >= output_freq) then
  113. output_counter = output_counter - output_freq
  114. !call output(state,etime)
  115. endif
  116. enddo
  117. call nvtxEndRange
  118. call system_clock(t2,rate)
  119. write(*,*) "CPU Time: ",dble(t2-t1)/dble(rate)
  120. !$acc end data
  121. !Deallocate
  122. call finalize()
  123. call nvtxEndRange
  124. contains
  125. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  126. !! HELPER SUBROUTINES AND FUNCTIONS GO HERE
  127. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  128. !Performs a single dimensionally split time step using a simple low-storate three-stage Runge-Kutta time integrator
  129. !The dimensional splitting is a second-order-accurate alternating Strang splitting in which the
  130. !order of directions is alternated each time step.
  131. !The Runge-Kutta method used here is defined as follows:
  132. ! q* = q[n] + dt/3 * rhs(q[n])
  133. ! q** = q[n] + dt/2 * rhs(q* )
  134. ! q[n+1] = q[n] + dt/1 * rhs(q** )
  135. subroutine perform_timestep(state,state_tmp,flux,tend,dt)
  136. implicit none
  137. real(rp), intent(inout) :: state (1-hs:nx+hs,1-hs:nz+hs,NUM_VARS)
  138. real(rp), intent( out) :: state_tmp(1-hs:nx+hs,1-hs:nz+hs,NUM_VARS)
  139. real(rp), intent( out) :: flux (nx+1,nz+1,NUM_VARS)
  140. real(rp), intent( out) :: tend (nx,nz,NUM_VARS)
  141. real(rp), intent(in ) :: dt
  142. logical, save :: direction_switch = .true.
  143. if (direction_switch) then
  144. !x-direction first
  145. call semi_discrete_step( state , state , state_tmp , dt / 3 , DIR_X , flux , tend )
  146. call semi_discrete_step( state , state_tmp , state_tmp , dt / 2 , DIR_X , flux , tend )
  147. call semi_discrete_step( state , state_tmp , state , dt / 1 , DIR_X , flux , tend )
  148. !z-direction second
  149. call semi_discrete_step( state , state , state_tmp , dt / 3 , DIR_Z , flux , tend )
  150. call semi_discrete_step( state , state_tmp , state_tmp , dt / 2 , DIR_Z , flux , tend )
  151. call semi_discrete_step( state , state_tmp , state , dt / 1 , DIR_Z , flux , tend )
  152. else
  153. !z-direction second
  154. call semi_discrete_step( state , state , state_tmp , dt / 3 , DIR_Z , flux , tend )
  155. call semi_discrete_step( state , state_tmp , state_tmp , dt / 2 , DIR_Z , flux , tend )
  156. call semi_discrete_step( state , state_tmp , state , dt / 1 , DIR_Z , flux , tend )
  157. !x-direction first
  158. call semi_discrete_step( state , state , state_tmp , dt / 3 , DIR_X , flux , tend )
  159. call semi_discrete_step( state , state_tmp , state_tmp , dt / 2 , DIR_X , flux , tend )
  160. call semi_discrete_step( state , state_tmp , state , dt / 1 , DIR_X , flux , tend )
  161. endif
  162. direction_switch = .not. direction_switch
  163. end subroutine perform_timestep
  164. !Perform a single semi-discretized step in time with the form:
  165. !state_out = state_init + dt * rhs(state_forcing)
  166. !Meaning the step starts from state_init, computes the rhs using state_forcing, and stores the result in state_out
  167. subroutine semi_discrete_step( state_init , state_forcing , state_out , dt , dir , flux , tend )
  168. implicit none
  169. real(rp), intent(in ) :: state_init (1-hs:nx+hs,1-hs:nz+hs,NUM_VARS)
  170. real(rp), intent(inout) :: state_forcing(1-hs:nx+hs,1-hs:nz+hs,NUM_VARS)
  171. real(rp), intent( out) :: state_out (1-hs:nx+hs,1-hs:nz+hs,NUM_VARS)
  172. real(rp), intent( out) :: flux (nx+1,nz+1,NUM_VARS)
  173. real(rp), intent( out) :: tend (nx,nz,NUM_VARS)
  174. real(rp), intent(in ) :: dt
  175. integer , intent(in ) :: dir
  176. integer :: i,k,ll
  177. if (dir == DIR_X) then
  178. !Set the halo values in the x-direction
  179. call set_halo_values_x(state_forcing)
  180. !Compute the time tendencies for the fluid state in the x-direction
  181. call compute_tendencies_x(state_forcing,flux,tend)
  182. elseif (dir == DIR_Z) then
  183. !Set the halo values in the z-direction
  184. call set_halo_values_z(state_forcing)
  185. !Compute the time tendencies for the fluid state in the z-direction
  186. call compute_tendencies_z(state_forcing,flux,tend)
  187. endif
  188. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  189. !! TODO: THREAD ME
  190. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  191. !Apply the tendencies to the fluid state
  192. !$acc parallel loop collapse(3)
  193. do ll = 1 , NUM_VARS
  194. do k = 1 , nz
  195. do i = 1 , nx
  196. state_out(i,k,ll) = state_init(i,k,ll) + dt * tend(i,k,ll)
  197. enddo
  198. enddo
  199. enddo
  200. end subroutine semi_discrete_step
  201. !Compute the time tendencies of the fluid state using forcing in the x-direction
  202. !First, compute the flux vector at each cell interface in the x-direction (including hyperviscosity)
  203. !Then, compute the tendencies using those fluxes
  204. subroutine compute_tendencies_x(state,flux,tend)
  205. implicit none
  206. real(rp), intent(in ) :: state(1-hs:nx+hs,1-hs:nz+hs,NUM_VARS)
  207. real(rp), intent( out) :: flux (nx+1,nz+1,NUM_VARS)
  208. real(rp), intent( out) :: tend (nx,nz,NUM_VARS)
  209. integer :: i,k,ll,s
  210. real(rp) :: r,u,w,t,p, stencil(4), d3_vals(NUM_VARS), vals(NUM_VARS), hv_coef
  211. !Compute the hyperviscosity coeficient
  212. hv_coef = -hv_beta * dx / (16*dt)
  213. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  214. !! TODO: THREAD ME
  215. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  216. !Compute fluxes in the x-direction for each cell
  217. !$acc parallel loop
  218. do k = 1 , nz
  219. !$acc loop private(stencil,vals,d3_vals)
  220. do i = 1 , nx+1
  221. !Use fourth-order interpolation from four cell averages to compute the value at the interface in question
  222. do ll = 1 , NUM_VARS
  223. do s = 1 , sten_size
  224. stencil(s) = state(i-hs-1+s,k,ll)
  225. enddo
  226. !Fourth-order-accurate interpolation of the state
  227. vals(ll) = -stencil(1)/12 + 7*stencil(2)/12 + 7*stencil(3)/12 - stencil(4)/12
  228. !First-order-accurate interpolation of the third spatial derivative of the state (for artificial viscosity)
  229. d3_vals(ll) = -stencil(1) + 3*stencil(2) - 3*stencil(3) + stencil(4)
  230. enddo
  231. !Compute density, u-wind, w-wind, potential temperature, and pressure (r,u,w,t,p respectively)
  232. r = vals(ID_DENS) + hy_dens_cell(k)
  233. u = vals(ID_UMOM) / r
  234. w = vals(ID_WMOM) / r
  235. t = ( vals(ID_RHOT) + hy_dens_theta_cell(k) ) / r
  236. p = C0*(r*t)**gamma
  237. !Compute the flux vector
  238. flux(i,k,ID_DENS) = r*u - hv_coef*d3_vals(ID_DENS)
  239. flux(i,k,ID_UMOM) = r*u*u+p - hv_coef*d3_vals(ID_UMOM)
  240. flux(i,k,ID_WMOM) = r*u*w - hv_coef*d3_vals(ID_WMOM)
  241. flux(i,k,ID_RHOT) = r*u*t - hv_coef*d3_vals(ID_RHOT)
  242. enddo
  243. enddo
  244. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  245. !! TODO: THREAD ME
  246. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  247. !Use the fluxes to compute tendencies for each cell
  248. !$acc parallel loop collapse(3)
  249. do ll = 1 , NUM_VARS
  250. do k = 1 , nz
  251. do i = 1 , nx
  252. tend(i,k,ll) = -( flux(i+1,k,ll) - flux(i,k,ll) ) / dx
  253. enddo
  254. enddo
  255. enddo
  256. end subroutine compute_tendencies_x
  257. !Compute the time tendencies of the fluid state using forcing in the z-direction
  258. !First, compute the flux vector at each cell interface in the z-direction (including hyperviscosity)
  259. !Then, compute the tendencies using those fluxes
  260. subroutine compute_tendencies_z(state,flux,tend)
  261. implicit none
  262. real(rp), intent(in ) :: state(1-hs:nx+hs,1-hs:nz+hs,NUM_VARS)
  263. real(rp), intent( out) :: flux (nx+1,nz+1,NUM_VARS)
  264. real(rp), intent( out) :: tend (nx,nz,NUM_VARS)
  265. integer :: i,k,ll,s
  266. real(rp) :: r,u,w,t,p, stencil(4), d3_vals(NUM_VARS), vals(NUM_VARS), hv_coef
  267. !Compute the hyperviscosity coeficient
  268. hv_coef = -hv_beta * dz / (16*dt)
  269. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  270. !! TODO: THREAD ME
  271. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  272. !Compute fluxes in the x-direction for each cell
  273. !$acc parallel loop
  274. do k = 1 , nz+1
  275. !$acc loop private(stencil,vals,d3_vals)
  276. do i = 1 , nx
  277. !Use fourth-order interpolation from four cell averages to compute the value at the interface in question
  278. do ll = 1 , NUM_VARS
  279. do s = 1 , sten_size
  280. stencil(s) = state(i,k-hs-1+s,ll)
  281. enddo
  282. !Fourth-order-accurate interpolation of the state
  283. vals(ll) = -stencil(1)/12 + 7*stencil(2)/12 + 7*stencil(3)/12 - stencil(4)/12
  284. !First-order-accurate interpolation of the third spatial derivative of the state
  285. d3_vals(ll) = -stencil(1) + 3*stencil(2) - 3*stencil(3) + stencil(4)
  286. enddo
  287. !Compute density, u-wind, w-wind, potential temperature, and pressure (r,u,w,t,p respectively)
  288. r = vals(ID_DENS) + hy_dens_int(k)
  289. u = vals(ID_UMOM) / r
  290. w = vals(ID_WMOM) / r
  291. t = ( vals(ID_RHOT) + hy_dens_theta_int(k) ) / r
  292. p = C0*(r*t)**gamma - hy_pressure_int(k)
  293. !Enforce vertical boundary condition and exact mass conservation
  294. if (k == 1 .or. k == nz+1) then
  295. w = 0
  296. d3_vals(ID_DENS) = 0
  297. endif
  298. !Compute the flux vector with hyperviscosity
  299. flux(i,k,ID_DENS) = r*w - hv_coef*d3_vals(ID_DENS)
  300. flux(i,k,ID_UMOM) = r*w*u - hv_coef*d3_vals(ID_UMOM)
  301. flux(i,k,ID_WMOM) = r*w*w+p - hv_coef*d3_vals(ID_WMOM)
  302. flux(i,k,ID_RHOT) = r*w*t - hv_coef*d3_vals(ID_RHOT)
  303. enddo
  304. enddo
  305. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  306. !! TODO: THREAD ME
  307. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  308. !Use the fluxes to compute tendencies for each cell
  309. !$acc parallel loop collapse(3)
  310. do ll = 1 , NUM_VARS
  311. do k = 1 , nz
  312. do i = 1 , nx
  313. tend(i,k,ll) = -( flux(i,k+1,ll) - flux(i,k,ll) ) / dz
  314. if (ll == ID_WMOM) then
  315. tend(i,k,ID_WMOM) = tend(i,k,ID_WMOM) - state(i,k,ID_DENS)*grav
  316. endif
  317. enddo
  318. enddo
  319. enddo
  320. end subroutine compute_tendencies_z
  321. !Set halo values in the x-direction.
  322. subroutine set_halo_values_x(state)
  323. implicit none
  324. real(rp), intent(inout) :: state(1-hs:nx+hs,1-hs:nz+hs,NUM_VARS)
  325. integer :: k, ll
  326. real(rp) :: z
  327. !$acc parallel loop collapse(2)
  328. do ll = 1 , NUM_VARS
  329. do k = 1 , nz
  330. state(-1 ,k,ll) = state(nx-1,k,ll)
  331. state(0 ,k,ll) = state(nx ,k,ll)
  332. state(nx+1,k,ll) = state(1 ,k,ll)
  333. state(nx+2,k,ll) = state(2 ,k,ll)
  334. enddo
  335. enddo
  336. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  337. if (myrank == 0) then
  338. !$acc parallel loop
  339. do k = 1 , nz
  340. z = (k_beg-1 + k-0.5_rp)*dz
  341. if (abs(z-3*zlen/4) <= zlen/16) then
  342. state(-1:0,k,ID_UMOM) = (state(-1:0,k,ID_DENS)+hy_dens_cell(k)) * 50._rp
  343. state(-1:0,k,ID_RHOT) = (state(-1:0,k,ID_DENS)+hy_dens_cell(k)) * 298._rp - hy_dens_theta_cell(k)
  344. endif
  345. enddo
  346. endif
  347. end subroutine set_halo_values_x
  348. !Set halo values in the z-direction.
  349. !decomposition in the vertical direction
  350. subroutine set_halo_values_z(state)
  351. implicit none
  352. real(rp), intent(inout) :: state(1-hs:nx+hs,1-hs:nz+hs,NUM_VARS)
  353. integer :: i, ll
  354. real(rp), parameter :: mnt_width = xlen/8
  355. real(rp) :: x, xloc, mnt_deriv
  356. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  357. !! TODO: THREAD ME
  358. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  359. !$acc parallel loop collapse(2) private(xloc,mnt_deriv,x)
  360. do ll = 1 , NUM_VARS
  361. do i = 1-hs,nx+hs
  362. if (ll == ID_WMOM) then
  363. state(i,-1 ,ll) = 0
  364. state(i,0 ,ll) = 0
  365. state(i,nz+1,ll) = 0
  366. state(i,nz+2,ll) = 0
  367. else
  368. state(i,-1 ,ll) = state(i,1 ,ll)
  369. state(i,0 ,ll) = state(i,1 ,ll)
  370. state(i,nz+1,ll) = state(i,nz,ll)
  371. state(i,nz+2,ll) = state(i,nz,ll)
  372. endif
  373. enddo
  374. enddo
  375. end subroutine set_halo_values_z
  376. !Initialize some grid settings, allocate data, and initialize the full grid and data
  377. subroutine init()
  378. implicit none
  379. integer :: i, k, ii, kk, ll, ierr
  380. real(rp) :: x, z, r, u, w, t, hr, ht
  381. !Set the cell grid size
  382. dx = xlen / nx_glob
  383. dz = zlen / nz_glob
  384. nranks = 1
  385. myrank = 0
  386. i_beg = 1
  387. nx = nx_glob
  388. left_rank = 0
  389. right_rank = 0
  390. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  391. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  392. !! YOU DON'T NEED TO ALTER ANYTHING BELOW THIS POINT IN THE CODE
  393. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  394. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  395. k_beg = 1
  396. nz = nz_glob
  397. !Allocate the model data
  398. allocate(state (1-hs:nx+hs,1-hs:nz+hs,NUM_VARS))
  399. allocate(state_tmp (1-hs:nx+hs,1-hs:nz+hs,NUM_VARS))
  400. allocate(flux (nx+1,nz+1,NUM_VARS))
  401. allocate(tend (nx,nz,NUM_VARS))
  402. allocate(hy_dens_cell (1-hs:nz+hs))
  403. allocate(hy_dens_theta_cell(1-hs:nz+hs))
  404. allocate(hy_dens_int (nz+1))
  405. allocate(hy_dens_theta_int (nz+1))
  406. allocate(hy_pressure_int (nz+1))
  407. !Define the maximum stable time step based on an assumed maximum wind speed
  408. dt = min(dx,dz) / max_speed * cfl
  409. !Set initial elapsed model time and output_counter to zero
  410. etime = 0
  411. output_counter = 0
  412. !Display some grid information
  413. write(*,*) 'nx_glob, nz_glob:', nx_glob, nz_glob
  414. write(*,*) 'dx,dz: ',dx,dz
  415. write(*,*) 'dt: ',dt
  416. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  417. !! Initialize the cell-averaged fluid state via Gauss-Legendre quadrature
  418. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  419. state = 0
  420. do k = 1-hs , nz+hs
  421. do i = 1-hs , nx+hs
  422. !Initialize the state to zero
  423. !Use Gauss-Legendre quadrature to initialize a hydrostatic balance + temperature perturbation
  424. do kk = 1 , nqpoints
  425. do ii = 1 , nqpoints
  426. !Compute the x,z location within the global domain based on cell and quadrature index
  427. x = (i_beg-1 + i-0.5_rp) * dx + (qpoints(ii)-0.5_rp)*dx
  428. z = (k_beg-1 + k-0.5_rp) * dz + (qpoints(kk)-0.5_rp)*dz
  429. !Set the fluid state based on the user's specification
  430. call injection (x,z,r,u,w,t,hr,ht)
  431. !Store into the fluid state array
  432. state(i,k,ID_DENS) = state(i,k,ID_DENS) + r * qweights(ii)*qweights(kk)
  433. state(i,k,ID_UMOM) = state(i,k,ID_UMOM) + (r+hr)*u * qweights(ii)*qweights(kk)
  434. state(i,k,ID_WMOM) = state(i,k,ID_WMOM) + (r+hr)*w * qweights(ii)*qweights(kk)
  435. state(i,k,ID_RHOT) = state(i,k,ID_RHOT) + ( (r+hr)*(t+ht) - hr*ht ) * qweights(ii)*qweights(kk)
  436. enddo
  437. enddo
  438. do ll = 1 , NUM_VARS
  439. state_tmp(i,k,ll) = state(i,k,ll)
  440. enddo
  441. enddo
  442. enddo
  443. !Compute the hydrostatic background state over vertical cell averages
  444. hy_dens_cell = 0
  445. hy_dens_theta_cell = 0
  446. do k = 1-hs , nz+hs
  447. do kk = 1 , nqpoints
  448. z = (k_beg-1 + k-0.5_rp) * dz + (qpoints(kk)-0.5_rp)*dz
  449. !Set the fluid state based on the user's specification
  450. call injection (0._rp,z,r,u,w,t,hr,ht)
  451. hy_dens_cell(k) = hy_dens_cell(k) + hr * qweights(kk)
  452. hy_dens_theta_cell(k) = hy_dens_theta_cell(k) + hr*ht * qweights(kk)
  453. enddo
  454. enddo
  455. !Compute the hydrostatic background state at vertical cell interfaces
  456. do k = 1 , nz+1
  457. z = (k_beg-1 + k-1) * dz
  458. call injection (0._rp,z,r,u,w,t,hr,ht)
  459. hy_dens_int (k) = hr
  460. hy_dens_theta_int(k) = hr*ht
  461. hy_pressure_int (k) = C0*(hr*ht)**gamma
  462. enddo
  463. end subroutine init
  464. !This test case is initially balanced but injects fast, cold air from the left boundary near the model top
  465. subroutine injection(x,z,r,u,w,t,hr,ht)
  466. implicit none
  467. real(rp), intent(in ) :: x, z !x- and z- location of the point being sampled
  468. real(rp), intent( out) :: r, u, w, t !Density, uwind, wwind, and potential temperature
  469. real(rp), intent( out) :: hr, ht !Hydrostatic density and potential temperature
  470. call hydro_const_theta(z,hr,ht)
  471. r = 0
  472. t = 0
  473. u = 0
  474. w = 0
  475. end subroutine injection
  476. !Establish hydrstatic balance using constant potential temperature (thermally neutral atmosphere)
  477. subroutine hydro_const_theta(z,r,t)
  478. implicit none
  479. real(rp), intent(in ) :: z !x- and z- location of the point being sampled
  480. real(rp), intent( out) :: r, t !Density and potential temperature at this point
  481. real(rp), parameter :: theta0 = 300._rp !Background potential temperature
  482. real(rp), parameter :: exner0 = 1._rp !Surface-level Exner pressure
  483. real(rp) :: p,exner,rt
  484. !Establish hydrostatic balance first using Exner pressure
  485. t = theta0 !Potential Temperature at z
  486. exner = exner0 - grav * z / (cp * theta0) !Exner pressure at z
  487. p = p0 * exner**(cp/rd) !Pressure at z
  488. rt = (p / c0)**(1._rp / gamma) !rho*theta at z
  489. r = rt / t !Density at z
  490. end subroutine hydro_const_theta
  491. !Deallocate and call Finalize
  492. subroutine finalize()
  493. implicit none
  494. integer :: ierr
  495. deallocate(state )
  496. deallocate(state_tmp )
  497. deallocate(flux )
  498. deallocate(tend )
  499. deallocate(hy_dens_cell )
  500. deallocate(hy_dens_theta_cell)
  501. deallocate(hy_dens_int )
  502. deallocate(hy_dens_theta_int )
  503. deallocate(hy_pressure_int )
  504. end subroutine finalize
  505. end program miniweather