miniWeather_openacc.f90 26 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574
  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 = 1000 !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. !$acc wait
  119. call system_clock(t2,rate)
  120. write(*,*) "CPU Time: ",dble(t2-t1)/dble(rate)
  121. !$acc end data
  122. !Deallocate
  123. call finalize()
  124. call nvtxEndRange
  125. contains
  126. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  127. !! HELPER SUBROUTINES AND FUNCTIONS GO HERE
  128. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  129. !Performs a single dimensionally split time step using a simple low-storate three-stage Runge-Kutta time integrator
  130. !The dimensional splitting is a second-order-accurate alternating Strang splitting in which the
  131. !order of directions is alternated each time step.
  132. !The Runge-Kutta method used here is defined as follows:
  133. ! q* = q[n] + dt/3 * rhs(q[n])
  134. ! q** = q[n] + dt/2 * rhs(q* )
  135. ! q[n+1] = q[n] + dt/1 * rhs(q** )
  136. subroutine perform_timestep(state,state_tmp,flux,tend,dt)
  137. implicit none
  138. real(rp), intent(inout) :: state (1-hs:nx+hs,1-hs:nz+hs,NUM_VARS)
  139. real(rp), intent( out) :: state_tmp(1-hs:nx+hs,1-hs:nz+hs,NUM_VARS)
  140. real(rp), intent( out) :: flux (nx+1,nz+1,NUM_VARS)
  141. real(rp), intent( out) :: tend (nx,nz,NUM_VARS)
  142. real(rp), intent(in ) :: dt
  143. logical, save :: direction_switch = .true.
  144. if (direction_switch) then
  145. !x-direction first
  146. call semi_discrete_step( state , state , state_tmp , dt / 3 , DIR_X , flux , tend )
  147. call semi_discrete_step( state , state_tmp , state_tmp , dt / 2 , DIR_X , flux , tend )
  148. call semi_discrete_step( state , state_tmp , state , dt / 1 , DIR_X , flux , tend )
  149. !z-direction second
  150. call semi_discrete_step( state , state , state_tmp , dt / 3 , DIR_Z , flux , tend )
  151. call semi_discrete_step( state , state_tmp , state_tmp , dt / 2 , DIR_Z , flux , tend )
  152. call semi_discrete_step( state , state_tmp , state , dt / 1 , DIR_Z , flux , tend )
  153. else
  154. !z-direction second
  155. call semi_discrete_step( state , state , state_tmp , dt / 3 , DIR_Z , flux , tend )
  156. call semi_discrete_step( state , state_tmp , state_tmp , dt / 2 , DIR_Z , flux , tend )
  157. call semi_discrete_step( state , state_tmp , state , dt / 1 , DIR_Z , flux , tend )
  158. !x-direction first
  159. call semi_discrete_step( state , state , state_tmp , dt / 3 , DIR_X , flux , tend )
  160. call semi_discrete_step( state , state_tmp , state_tmp , dt / 2 , DIR_X , flux , tend )
  161. call semi_discrete_step( state , state_tmp , state , dt / 1 , DIR_X , flux , tend )
  162. endif
  163. direction_switch = .not. direction_switch
  164. end subroutine perform_timestep
  165. !Perform a single semi-discretized step in time with the form:
  166. !state_out = state_init + dt * rhs(state_forcing)
  167. !Meaning the step starts from state_init, computes the rhs using state_forcing, and stores the result in state_out
  168. subroutine semi_discrete_step( state_init , state_forcing , state_out , dt , dir , flux , tend )
  169. implicit none
  170. real(rp), intent(in ) :: state_init (1-hs:nx+hs,1-hs:nz+hs,NUM_VARS)
  171. real(rp), intent(inout) :: state_forcing(1-hs:nx+hs,1-hs:nz+hs,NUM_VARS)
  172. real(rp), intent( out) :: state_out (1-hs:nx+hs,1-hs:nz+hs,NUM_VARS)
  173. real(rp), intent( out) :: flux (nx+1,nz+1,NUM_VARS)
  174. real(rp), intent( out) :: tend (nx,nz,NUM_VARS)
  175. real(rp), intent(in ) :: dt
  176. integer , intent(in ) :: dir
  177. integer :: i,k,ll
  178. if (dir == DIR_X) then
  179. !Set the halo values in the x-direction
  180. call set_halo_values_x(state_forcing)
  181. !Compute the time tendencies for the fluid state in the x-direction
  182. call compute_tendencies_x(state_forcing,flux,tend)
  183. elseif (dir == DIR_Z) then
  184. !Set the halo values in the z-direction
  185. call set_halo_values_z(state_forcing)
  186. !Compute the time tendencies for the fluid state in the z-direction
  187. call compute_tendencies_z(state_forcing,flux,tend)
  188. endif
  189. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  190. !! TODO: THREAD ME
  191. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  192. !Apply the tendencies to the fluid state
  193. !$acc parallel loop collapse(3)
  194. do ll = 1 , NUM_VARS
  195. do k = 1 , nz
  196. do i = 1 , nx
  197. state_out(i,k,ll) = state_init(i,k,ll) + dt * tend(i,k,ll)
  198. enddo
  199. enddo
  200. enddo
  201. end subroutine semi_discrete_step
  202. !Compute the time tendencies of the fluid state using forcing in the x-direction
  203. !First, compute the flux vector at each cell interface in the x-direction (including hyperviscosity)
  204. !Then, compute the tendencies using those fluxes
  205. subroutine compute_tendencies_x(state,flux,tend)
  206. implicit none
  207. real(rp), intent(in ) :: state(1-hs:nx+hs,1-hs:nz+hs,NUM_VARS)
  208. real(rp), intent( out) :: flux (nx+1,nz+1,NUM_VARS)
  209. real(rp), intent( out) :: tend (nx,nz,NUM_VARS)
  210. integer :: i,k,ll,s
  211. real(rp) :: r,u,w,t,p, stencil(4), d3_vals(NUM_VARS), vals(NUM_VARS), hv_coef
  212. !Compute the hyperviscosity coeficient
  213. hv_coef = -hv_beta * dx / (16*dt)
  214. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  215. !! TODO: THREAD ME
  216. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  217. !Compute fluxes in the x-direction for each cell
  218. !$acc parallel loop gang
  219. do k = 1 , nz
  220. !$acc loop vector private(stencil,vals,d3_vals)
  221. do i = 1 , nx+1
  222. !Use fourth-order interpolation from four cell averages to compute the value at the interface in question
  223. do ll = 1 , NUM_VARS
  224. do s = 1 , sten_size
  225. stencil(s) = state(i-hs-1+s,k,ll)
  226. enddo
  227. !Fourth-order-accurate interpolation of the state
  228. vals(ll) = -stencil(1)/12 + 7*stencil(2)/12 + 7*stencil(3)/12 - stencil(4)/12
  229. !First-order-accurate interpolation of the third spatial derivative of the state (for artificial viscosity)
  230. d3_vals(ll) = -stencil(1) + 3*stencil(2) - 3*stencil(3) + stencil(4)
  231. enddo
  232. !Compute density, u-wind, w-wind, potential temperature, and pressure (r,u,w,t,p respectively)
  233. r = vals(ID_DENS) + hy_dens_cell(k)
  234. u = vals(ID_UMOM) / r
  235. w = vals(ID_WMOM) / r
  236. t = ( vals(ID_RHOT) + hy_dens_theta_cell(k) ) / r
  237. p = C0*(r*t)**gamma
  238. !Compute the flux vector
  239. flux(i,k,ID_DENS) = r*u - hv_coef*d3_vals(ID_DENS)
  240. flux(i,k,ID_UMOM) = r*u*u+p - hv_coef*d3_vals(ID_UMOM)
  241. flux(i,k,ID_WMOM) = r*u*w - hv_coef*d3_vals(ID_WMOM)
  242. flux(i,k,ID_RHOT) = r*u*t - hv_coef*d3_vals(ID_RHOT)
  243. enddo
  244. enddo
  245. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  246. !! TODO: THREAD ME
  247. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  248. !Use the fluxes to compute tendencies for each cell
  249. !$acc parallel loop collapse(3)
  250. do ll = 1 , NUM_VARS
  251. do k = 1 , nz
  252. do i = 1 , nx
  253. tend(i,k,ll) = -( flux(i+1,k,ll) - flux(i,k,ll) ) / dx
  254. enddo
  255. enddo
  256. enddo
  257. end subroutine compute_tendencies_x
  258. !Compute the time tendencies of the fluid state using forcing in the z-direction
  259. !First, compute the flux vector at each cell interface in the z-direction (including hyperviscosity)
  260. !Then, compute the tendencies using those fluxes
  261. subroutine compute_tendencies_z(state,flux,tend)
  262. implicit none
  263. real(rp), intent(in ) :: state(1-hs:nx+hs,1-hs:nz+hs,NUM_VARS)
  264. real(rp), intent( out) :: flux (nx+1,nz+1,NUM_VARS)
  265. real(rp), intent( out) :: tend (nx,nz,NUM_VARS)
  266. integer :: i,k,ll,s
  267. real(rp) :: r,u,w,t,p, stencil(4), d3_vals(NUM_VARS), vals(NUM_VARS), hv_coef
  268. !Compute the hyperviscosity coeficient
  269. hv_coef = -hv_beta * dz / (16*dt)
  270. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  271. !! TODO: THREAD ME
  272. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  273. !Compute fluxes in the x-direction for each cell
  274. !$acc parallel loop gang
  275. do k = 1 , nz+1
  276. !$acc loop vector private(stencil,vals,d3_vals)
  277. do i = 1 , nx
  278. !Use fourth-order interpolation from four cell averages to compute the value at the interface in question
  279. do ll = 1 , NUM_VARS
  280. do s = 1 , sten_size
  281. stencil(s) = state(i,k-hs-1+s,ll)
  282. enddo
  283. !Fourth-order-accurate interpolation of the state
  284. vals(ll) = -stencil(1)/12 + 7*stencil(2)/12 + 7*stencil(3)/12 - stencil(4)/12
  285. !First-order-accurate interpolation of the third spatial derivative of the state
  286. d3_vals(ll) = -stencil(1) + 3*stencil(2) - 3*stencil(3) + stencil(4)
  287. enddo
  288. !Compute density, u-wind, w-wind, potential temperature, and pressure (r,u,w,t,p respectively)
  289. r = vals(ID_DENS) + hy_dens_int(k)
  290. u = vals(ID_UMOM) / r
  291. w = vals(ID_WMOM) / r
  292. t = ( vals(ID_RHOT) + hy_dens_theta_int(k) ) / r
  293. p = C0*(r*t)**gamma - hy_pressure_int(k)
  294. !Enforce vertical boundary condition and exact mass conservation
  295. if (k == 1 .or. k == nz+1) then
  296. w = 0
  297. d3_vals(ID_DENS) = 0
  298. endif
  299. !Compute the flux vector with hyperviscosity
  300. flux(i,k,ID_DENS) = r*w - hv_coef*d3_vals(ID_DENS)
  301. flux(i,k,ID_UMOM) = r*w*u - hv_coef*d3_vals(ID_UMOM)
  302. flux(i,k,ID_WMOM) = r*w*w+p - hv_coef*d3_vals(ID_WMOM)
  303. flux(i,k,ID_RHOT) = r*w*t - hv_coef*d3_vals(ID_RHOT)
  304. enddo
  305. enddo
  306. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  307. !! TODO: THREAD ME
  308. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  309. !Use the fluxes to compute tendencies for each cell
  310. !$acc parallel loop collapse(3)
  311. do ll = 1 , NUM_VARS
  312. do k = 1 , nz
  313. do i = 1 , nx
  314. tend(i,k,ll) = -( flux(i,k+1,ll) - flux(i,k,ll) ) / dz
  315. if (ll == ID_WMOM) then
  316. tend(i,k,ID_WMOM) = tend(i,k,ID_WMOM) - state(i,k,ID_DENS)*grav
  317. endif
  318. enddo
  319. enddo
  320. enddo
  321. end subroutine compute_tendencies_z
  322. !Set halo values in the x-direction.
  323. subroutine set_halo_values_x(state)
  324. implicit none
  325. real(rp), intent(inout) :: state(1-hs:nx+hs,1-hs:nz+hs,NUM_VARS)
  326. integer :: k, ll
  327. real(rp) :: z
  328. !$acc parallel loop collapse(2)
  329. do ll = 1 , NUM_VARS
  330. do k = 1 , nz
  331. state(-1 ,k,ll) = state(nx-1,k,ll)
  332. state(0 ,k,ll) = state(nx ,k,ll)
  333. state(nx+1,k,ll) = state(1 ,k,ll)
  334. state(nx+2,k,ll) = state(2 ,k,ll)
  335. enddo
  336. enddo
  337. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  338. if (myrank == 0) then
  339. !$acc parallel loop
  340. do k = 1 , nz
  341. z = (k_beg-1 + k-0.5_rp)*dz
  342. if (abs(z-3*zlen/4) <= zlen/16) then
  343. state(-1:0,k,ID_UMOM) = (state(-1:0,k,ID_DENS)+hy_dens_cell(k)) * 50._rp
  344. state(-1:0,k,ID_RHOT) = (state(-1:0,k,ID_DENS)+hy_dens_cell(k)) * 298._rp - hy_dens_theta_cell(k)
  345. endif
  346. enddo
  347. endif
  348. end subroutine set_halo_values_x
  349. !Set halo values in the z-direction.
  350. !decomposition in the vertical direction
  351. subroutine set_halo_values_z(state)
  352. implicit none
  353. real(rp), intent(inout) :: state(1-hs:nx+hs,1-hs:nz+hs,NUM_VARS)
  354. integer :: i, ll
  355. real(rp), parameter :: mnt_width = xlen/8
  356. real(rp) :: x, xloc, mnt_deriv
  357. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  358. !! TODO: THREAD ME
  359. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  360. !$acc parallel loop collapse(2) private(xloc,mnt_deriv,x)
  361. do ll = 1 , NUM_VARS
  362. do i = 1-hs,nx+hs
  363. if (ll == ID_WMOM) then
  364. state(i,-1 ,ll) = 0
  365. state(i,0 ,ll) = 0
  366. state(i,nz+1,ll) = 0
  367. state(i,nz+2,ll) = 0
  368. else
  369. state(i,-1 ,ll) = state(i,1 ,ll)
  370. state(i,0 ,ll) = state(i,1 ,ll)
  371. state(i,nz+1,ll) = state(i,nz,ll)
  372. state(i,nz+2,ll) = state(i,nz,ll)
  373. endif
  374. enddo
  375. enddo
  376. end subroutine set_halo_values_z
  377. !Initialize some grid settings, allocate data, and initialize the full grid and data
  378. subroutine init()
  379. implicit none
  380. integer :: i, k, ii, kk, ll, ierr
  381. real(rp) :: x, z, r, u, w, t, hr, ht
  382. !Set the cell grid size
  383. dx = xlen / nx_glob
  384. dz = zlen / nz_glob
  385. nranks = 1
  386. myrank = 0
  387. i_beg = 1
  388. nx = nx_glob
  389. left_rank = 0
  390. right_rank = 0
  391. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  392. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  393. !! YOU DON'T NEED TO ALTER ANYTHING BELOW THIS POINT IN THE CODE
  394. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  395. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  396. k_beg = 1
  397. nz = nz_glob
  398. !Allocate the model data
  399. allocate(state (1-hs:nx+hs,1-hs:nz+hs,NUM_VARS))
  400. allocate(state_tmp (1-hs:nx+hs,1-hs:nz+hs,NUM_VARS))
  401. allocate(flux (nx+1,nz+1,NUM_VARS))
  402. allocate(tend (nx,nz,NUM_VARS))
  403. allocate(hy_dens_cell (1-hs:nz+hs))
  404. allocate(hy_dens_theta_cell(1-hs:nz+hs))
  405. allocate(hy_dens_int (nz+1))
  406. allocate(hy_dens_theta_int (nz+1))
  407. allocate(hy_pressure_int (nz+1))
  408. !Define the maximum stable time step based on an assumed maximum wind speed
  409. dt = min(dx,dz) / max_speed * cfl
  410. !Set initial elapsed model time and output_counter to zero
  411. etime = 0
  412. output_counter = 0
  413. !Display some grid information
  414. write(*,*) 'nx_glob, nz_glob:', nx_glob, nz_glob
  415. write(*,*) 'dx,dz: ',dx,dz
  416. write(*,*) 'dt: ',dt
  417. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  418. !! Initialize the cell-averaged fluid state via Gauss-Legendre quadrature
  419. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  420. state = 0
  421. do k = 1-hs , nz+hs
  422. do i = 1-hs , nx+hs
  423. !Initialize the state to zero
  424. !Use Gauss-Legendre quadrature to initialize a hydrostatic balance + temperature perturbation
  425. do kk = 1 , nqpoints
  426. do ii = 1 , nqpoints
  427. !Compute the x,z location within the global domain based on cell and quadrature index
  428. x = (i_beg-1 + i-0.5_rp) * dx + (qpoints(ii)-0.5_rp)*dx
  429. z = (k_beg-1 + k-0.5_rp) * dz + (qpoints(kk)-0.5_rp)*dz
  430. !Set the fluid state based on the user's specification
  431. call injection (x,z,r,u,w,t,hr,ht)
  432. !Store into the fluid state array
  433. state(i,k,ID_DENS) = state(i,k,ID_DENS) + r * qweights(ii)*qweights(kk)
  434. state(i,k,ID_UMOM) = state(i,k,ID_UMOM) + (r+hr)*u * qweights(ii)*qweights(kk)
  435. state(i,k,ID_WMOM) = state(i,k,ID_WMOM) + (r+hr)*w * qweights(ii)*qweights(kk)
  436. state(i,k,ID_RHOT) = state(i,k,ID_RHOT) + ( (r+hr)*(t+ht) - hr*ht ) * qweights(ii)*qweights(kk)
  437. enddo
  438. enddo
  439. do ll = 1 , NUM_VARS
  440. state_tmp(i,k,ll) = state(i,k,ll)
  441. enddo
  442. enddo
  443. enddo
  444. !Compute the hydrostatic background state over vertical cell averages
  445. hy_dens_cell = 0
  446. hy_dens_theta_cell = 0
  447. do k = 1-hs , nz+hs
  448. do kk = 1 , nqpoints
  449. z = (k_beg-1 + k-0.5_rp) * dz + (qpoints(kk)-0.5_rp)*dz
  450. !Set the fluid state based on the user's specification
  451. call injection (0._rp,z,r,u,w,t,hr,ht)
  452. hy_dens_cell(k) = hy_dens_cell(k) + hr * qweights(kk)
  453. hy_dens_theta_cell(k) = hy_dens_theta_cell(k) + hr*ht * qweights(kk)
  454. enddo
  455. enddo
  456. !Compute the hydrostatic background state at vertical cell interfaces
  457. do k = 1 , nz+1
  458. z = (k_beg-1 + k-1) * dz
  459. call injection (0._rp,z,r,u,w,t,hr,ht)
  460. hy_dens_int (k) = hr
  461. hy_dens_theta_int(k) = hr*ht
  462. hy_pressure_int (k) = C0*(hr*ht)**gamma
  463. enddo
  464. end subroutine init
  465. !This test case is initially balanced but injects fast, cold air from the left boundary near the model top
  466. subroutine injection(x,z,r,u,w,t,hr,ht)
  467. implicit none
  468. real(rp), intent(in ) :: x, z !x- and z- location of the point being sampled
  469. real(rp), intent( out) :: r, u, w, t !Density, uwind, wwind, and potential temperature
  470. real(rp), intent( out) :: hr, ht !Hydrostatic density and potential temperature
  471. call hydro_const_theta(z,hr,ht)
  472. r = 0
  473. t = 0
  474. u = 0
  475. w = 0
  476. end subroutine injection
  477. !Establish hydrstatic balance using constant potential temperature (thermally neutral atmosphere)
  478. subroutine hydro_const_theta(z,r,t)
  479. implicit none
  480. real(rp), intent(in ) :: z !x- and z- location of the point being sampled
  481. real(rp), intent( out) :: r, t !Density and potential temperature at this point
  482. real(rp), parameter :: theta0 = 300._rp !Background potential temperature
  483. real(rp), parameter :: exner0 = 1._rp !Surface-level Exner pressure
  484. real(rp) :: p,exner,rt
  485. !Establish hydrostatic balance first using Exner pressure
  486. t = theta0 !Potential Temperature at z
  487. exner = exner0 - grav * z / (cp * theta0) !Exner pressure at z
  488. p = p0 * exner**(cp/rd) !Pressure at z
  489. rt = (p / c0)**(1._rp / gamma) !rho*theta at z
  490. r = rt / t !Density at z
  491. end subroutine hydro_const_theta
  492. !Deallocate and call Finalize
  493. subroutine finalize()
  494. implicit none
  495. integer :: ierr
  496. deallocate(state )
  497. deallocate(state_tmp )
  498. deallocate(flux )
  499. deallocate(tend )
  500. deallocate(hy_dens_cell )
  501. deallocate(hy_dens_theta_cell)
  502. deallocate(hy_dens_int )
  503. deallocate(hy_dens_theta_int )
  504. deallocate(hy_pressure_int )
  505. end subroutine finalize
  506. end program miniweather