miniWeather_serial.f90 25 KB

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