PageRenderTime 49ms CodeModel.GetById 13ms RepoModel.GetById 1ms app.codeStats 0ms

/wrfv2_fire/dyn_em/module_initialize_grav2d_x.F

http://github.com/jbeezley/wrf-fire
FORTRAN Legacy | 834 lines | 474 code | 201 blank | 159 comment | 13 complexity | 0b25372741f95a550725c629e0ba230c MD5 | raw file
Possible License(s): AGPL-1.0
  1. !IDEAL:MODEL_LAYER:INITIALIZATION
  2. !
  3. ! This MODULE holds the routines which are used to perform various initializations
  4. ! for the individual domains.
  5. ! This MODULE CONTAINS the following routines:
  6. ! initialize_field_test - 1. Set different fields to different constant
  7. ! values. This is only a test. If the correct
  8. ! domain is not found (based upon the "id")
  9. ! then a fatal error is issued.
  10. !-----------------------------------------------------------------------
  11. MODULE module_initialize_ideal
  12. USE module_domain
  13. USE module_io_domain
  14. USE module_state_description
  15. USE module_model_constants
  16. USE module_bc
  17. USE module_timing
  18. USE module_configure
  19. USE module_init_utilities
  20. #ifdef DM_PARALLEL
  21. USE module_dm
  22. #endif
  23. CONTAINS
  24. !-------------------------------------------------------------------
  25. ! this is a wrapper for the solver-specific init_domain routines.
  26. ! Also dereferences the grid variables and passes them down as arguments.
  27. ! This is crucial, since the lower level routines may do message passing
  28. ! and this will get fouled up on machines that insist on passing down
  29. ! copies of assumed-shape arrays (by passing down as arguments, the
  30. ! data are treated as assumed-size -- ie. f77 -- arrays and the copying
  31. ! business is avoided). Fie on the F90 designers. Fie and a pox.
  32. SUBROUTINE init_domain ( grid )
  33. IMPLICIT NONE
  34. ! Input data.
  35. TYPE (domain), POINTER :: grid
  36. ! Local data.
  37. INTEGER :: idum1, idum2
  38. CALL set_scalar_indices_from_config ( head_grid%id , idum1, idum2 )
  39. CALL init_domain_rk( grid &
  40. !
  41. #include <actual_new_args.inc>
  42. !
  43. )
  44. END SUBROUTINE init_domain
  45. !-------------------------------------------------------------------
  46. SUBROUTINE init_domain_rk ( grid &
  47. !
  48. # include <dummy_new_args.inc>
  49. !
  50. )
  51. IMPLICIT NONE
  52. ! Input data.
  53. TYPE (domain), POINTER :: grid
  54. # include <dummy_new_decl.inc>
  55. TYPE (grid_config_rec_type) :: config_flags
  56. ! Local data
  57. INTEGER :: &
  58. ids, ide, jds, jde, kds, kde, &
  59. ims, ime, jms, jme, kms, kme, &
  60. its, ite, jts, jte, kts, kte, &
  61. i, j, k
  62. ! Local data
  63. INTEGER, PARAMETER :: nl_max = 1000
  64. REAL, DIMENSION(nl_max) :: zk, p_in, theta, rho, u, v, qv, pd_in
  65. INTEGER :: nl_in
  66. INTEGER :: icm,jcm, ii, im1, jj, jm1, loop, error, fid, nxc, nyc
  67. REAL :: u_mean,v_mean, f0, p_surf, p_level, qvf, z_at_v, z_at_u
  68. REAL :: z_scale, xrad, yrad, zrad, rad, delt, cof1, cof2, t_min, t_max
  69. ! REAL, EXTERNAL :: interp_0
  70. REAL :: hm, xa, xpos, xposml, xpospl
  71. REAL :: pi
  72. ! stuff from original initialization that has been dropped from the Registry
  73. REAL :: vnu, xnu, xnus, dinit0, cbh, p0_temp, t0_temp, zd, zt
  74. REAL :: qvf1, qvf2, pd_surf
  75. INTEGER :: it
  76. real :: thtmp, ptmp, temp(3)
  77. LOGICAL :: moisture_init
  78. LOGICAL :: stretch_grid, dry_sounding
  79. REAL :: xa1, xal1,pii,hm1 ! data for intercomparison setup from dale
  80. #ifdef DM_PARALLEL
  81. # include <data_calls.inc>
  82. #endif
  83. SELECT CASE ( model_data_order )
  84. CASE ( DATA_ORDER_ZXY )
  85. kds = grid%sd31 ; kde = grid%ed31 ;
  86. ids = grid%sd32 ; ide = grid%ed32 ;
  87. jds = grid%sd33 ; jde = grid%ed33 ;
  88. kms = grid%sm31 ; kme = grid%em31 ;
  89. ims = grid%sm32 ; ime = grid%em32 ;
  90. jms = grid%sm33 ; jme = grid%em33 ;
  91. kts = grid%sp31 ; kte = grid%ep31 ; ! note that tile is entire patch
  92. its = grid%sp32 ; ite = grid%ep32 ; ! note that tile is entire patch
  93. jts = grid%sp33 ; jte = grid%ep33 ; ! note that tile is entire patch
  94. CASE ( DATA_ORDER_XYZ )
  95. ids = grid%sd31 ; ide = grid%ed31 ;
  96. jds = grid%sd32 ; jde = grid%ed32 ;
  97. kds = grid%sd33 ; kde = grid%ed33 ;
  98. ims = grid%sm31 ; ime = grid%em31 ;
  99. jms = grid%sm32 ; jme = grid%em32 ;
  100. kms = grid%sm33 ; kme = grid%em33 ;
  101. its = grid%sp31 ; ite = grid%ep31 ; ! note that tile is entire patch
  102. jts = grid%sp32 ; jte = grid%ep32 ; ! note that tile is entire patch
  103. kts = grid%sp33 ; kte = grid%ep33 ; ! note that tile is entire patch
  104. CASE ( DATA_ORDER_XZY )
  105. ids = grid%sd31 ; ide = grid%ed31 ;
  106. kds = grid%sd32 ; kde = grid%ed32 ;
  107. jds = grid%sd33 ; jde = grid%ed33 ;
  108. ims = grid%sm31 ; ime = grid%em31 ;
  109. kms = grid%sm32 ; kme = grid%em32 ;
  110. jms = grid%sm33 ; jme = grid%em33 ;
  111. its = grid%sp31 ; ite = grid%ep31 ; ! note that tile is entire patch
  112. kts = grid%sp32 ; kte = grid%ep32 ; ! note that tile is entire patch
  113. jts = grid%sp33 ; jte = grid%ep33 ; ! note that tile is entire patch
  114. END SELECT
  115. hm = 000.
  116. xa = 5.0
  117. icm = ide/2
  118. xa1 = 5000./500.
  119. xal1 = 4000./500.
  120. pii = 2.*asin(1.0)
  121. hm1 = 250.
  122. ! hm1 = 1000.
  123. stretch_grid = .true.
  124. ! z_scale = .50
  125. z_scale = 1.675
  126. pi = 2.*asin(1.0)
  127. write(6,*) ' pi is ',pi
  128. nxc = (ide-ids)/2
  129. nyc = (jde-jds)/2
  130. CALL model_to_grid_config_rec ( grid%id , model_config_rec , config_flags )
  131. ! here we check to see if the boundary conditions are set properly
  132. CALL boundary_condition_check( config_flags, bdyzone, error, grid%id )
  133. moisture_init = .true.
  134. grid%itimestep=0
  135. #ifdef DM_PARALLEL
  136. CALL wrf_dm_bcast_bytes( icm , IWORDSIZE )
  137. CALL wrf_dm_bcast_bytes( jcm , IWORDSIZE )
  138. #endif
  139. CALL nl_set_mminlu(1, ' ')
  140. CALL nl_set_iswater(1,0)
  141. CALL nl_set_cen_lat(1,40.)
  142. CALL nl_set_cen_lon(1,-105.)
  143. CALL nl_set_truelat1(1,0.)
  144. CALL nl_set_truelat2(1,0.)
  145. CALL nl_set_moad_cen_lat (1,0.)
  146. CALL nl_set_stand_lon (1,0.)
  147. CALL nl_set_pole_lon (1,0.)
  148. CALL nl_set_pole_lat (1,90.)
  149. CALL nl_set_map_proj(1,0)
  150. ! here we initialize data we currently is not initialized
  151. ! in the input data
  152. DO j = jts, jte
  153. DO i = its, ite
  154. grid%msftx(i,j) = 1.
  155. grid%msfty(i,j) = 1.
  156. grid%msfux(i,j) = 1.
  157. grid%msfuy(i,j) = 1.
  158. grid%msfvx(i,j) = 1.
  159. grid%msfvx_inv(i,j)= 1.
  160. grid%msfvy(i,j) = 1.
  161. grid%sina(i,j) = 0.
  162. grid%cosa(i,j) = 1.
  163. grid%e(i,j) = 0.
  164. grid%f(i,j) = 0.
  165. END DO
  166. END DO
  167. DO j = jts, jte
  168. DO k = kts, kte
  169. DO i = its, ite
  170. grid%ww(i,k,j) = 0.
  171. END DO
  172. END DO
  173. END DO
  174. grid%step_number = 0
  175. ! set up the grid
  176. IF (stretch_grid) THEN ! exponential stretch for eta (nearly constant dz)
  177. DO k=1, kde
  178. grid%znw(k) = (exp(-(k-1)/float(kde-1)/z_scale) - exp(-1./z_scale))/ &
  179. (1.-exp(-1./z_scale))
  180. ENDDO
  181. ELSE
  182. DO k=1, kde
  183. grid%znw(k) = 1. - float(k-1)/float(kde-1)
  184. ENDDO
  185. ENDIF
  186. DO k=1, kde-1
  187. grid%dnw(k) = grid%znw(k+1) - grid%znw(k)
  188. grid%rdnw(k) = 1./grid%dnw(k)
  189. grid%znu(k) = 0.5*(grid%znw(k+1)+grid%znw(k))
  190. ENDDO
  191. DO k=2, kde-1
  192. grid%dn(k) = 0.5*(grid%dnw(k)+grid%dnw(k-1))
  193. grid%rdn(k) = 1./grid%dn(k)
  194. grid%fnp(k) = .5* grid%dnw(k )/grid%dn(k)
  195. grid%fnm(k) = .5* grid%dnw(k-1)/grid%dn(k)
  196. ENDDO
  197. cof1 = (2.*grid%dn(2)+grid%dn(3))/(grid%dn(2)+grid%dn(3))*grid%dnw(1)/grid%dn(2)
  198. cof2 = grid%dn(2) /(grid%dn(2)+grid%dn(3))*grid%dnw(1)/grid%dn(3)
  199. grid%cf1 = grid%fnp(2) + cof1
  200. grid%cf2 = grid%fnm(2) - cof1 - cof2
  201. grid%cf3 = cof2
  202. grid%cfn = (.5*grid%dnw(kde-1)+grid%dn(kde-1))/grid%dn(kde-1)
  203. grid%cfn1 = -.5*grid%dnw(kde-1)/grid%dn(kde-1)
  204. grid%rdx = 1./config_flags%dx
  205. grid%rdy = 1./config_flags%dy
  206. ! get the sounding from the ascii sounding file, first get dry sounding and
  207. ! calculate base state
  208. write(6,*) ' getting dry sounding for base state '
  209. dry_sounding = .true.
  210. CALL get_sounding( zk, p_in, pd_in, theta, rho, u, v, qv, dry_sounding, &
  211. nl_max, nl_in, .true.)
  212. write(6,*) ' returned from reading sounding, nl_in is ',nl_in
  213. ! find ptop for the desired ztop (ztop is input from the namelist),
  214. ! and find surface pressure
  215. grid%p_top = interp_0( p_in, zk, config_flags%ztop, nl_in )
  216. DO j=jts,jte
  217. DO i=its,ite ! flat surface
  218. !! grid%ht(i,j) = 0.
  219. grid%ht(i,j) = hm/(1.+(float(i-icm)/xa)**2)
  220. ! grid%ht(i,j) = hm1*exp(-(( float(i-icm)/xa1)**2)) &
  221. ! *( (cos(pii*float(i-icm)/xal1))**2 )
  222. grid%phb(i,1,j) = g*grid%ht(i,j)
  223. grid%php(i,1,j) = 0.
  224. grid%ph0(i,1,j) = grid%phb(i,1,j)
  225. ENDDO
  226. ENDDO
  227. DO J = jts, jte
  228. DO I = its, ite
  229. p_surf = interp_0( p_in, zk, grid%phb(i,1,j)/g, nl_in )
  230. grid%mub(i,j) = p_surf-grid%p_top
  231. ! this is dry hydrostatic sounding (base state), so given p (coordinate),
  232. ! interp theta (from interp) and compute 1/rho from eqn. of state
  233. DO K = 1, kte-1
  234. p_level = grid%znu(k)*(p_surf - grid%p_top) + grid%p_top
  235. grid%pb(i,k,j) = p_level
  236. grid%t_init(i,k,j) = interp_0( theta, p_in, p_level, nl_in ) - t0
  237. grid%alb(i,k,j) = (r_d/p1000mb)*(grid%t_init(i,k,j)+t0)*(grid%pb(i,k,j)/p1000mb)**cvpm
  238. ENDDO
  239. ! calc hydrostatic balance (alternatively we could interp the geopotential from the
  240. ! sounding, but this assures that the base state is in exact hydrostatic balance with
  241. ! respect to the model eqns.
  242. DO k = 2,kte
  243. grid%phb(i,k,j) = grid%phb(i,k-1,j) - grid%dnw(k-1)*grid%mub(i,j)*grid%alb(i,k-1,j)
  244. ENDDO
  245. ENDDO
  246. ENDDO
  247. write(6,*) ' ptop is ',grid%p_top
  248. write(6,*) ' base state mub(1,1), p_surf is ',grid%mub(1,1),grid%mub(1,1)+grid%p_top
  249. ! calculate full state for each column - this includes moisture.
  250. write(6,*) ' getting moist sounding for full state '
  251. dry_sounding = .false.
  252. CALL get_sounding( zk, p_in, pd_in, theta, rho, u, v, qv, dry_sounding, &
  253. nl_max, nl_in, .false. )
  254. DO J = jts, min(jde-1,jte)
  255. DO I = its, min(ide-1,ite)
  256. ! At this point grid%p_top is already set. find the DRY mass in the column
  257. ! by interpolating the DRY pressure.
  258. pd_surf = interp_0( pd_in, zk, grid%phb(i,1,j)/g, nl_in )
  259. ! compute the perturbation mass and the full mass
  260. grid%mu_1(i,j) = pd_surf-grid%p_top - grid%mub(i,j)
  261. grid%mu_2(i,j) = grid%mu_1(i,j)
  262. grid%mu0(i,j) = grid%mu_1(i,j) + grid%mub(i,j)
  263. ! given the dry pressure and coordinate system, interp the potential
  264. ! temperature and qv
  265. do k=1,kde-1
  266. p_level = grid%znu(k)*(pd_surf - grid%p_top) + grid%p_top
  267. grid%moist(i,k,j,P_QV) = interp_0( qv, pd_in, p_level, nl_in )
  268. grid%t_1(i,k,j) = interp_0( theta, pd_in, p_level, nl_in ) - t0
  269. grid%t_2(i,k,j) = grid%t_1(i,k,j)
  270. enddo
  271. ! integrate the hydrostatic equation (from the RHS of the bigstep
  272. ! vertical momentum equation) down from the top to get p.
  273. ! first from the top of the model to the top pressure
  274. k = kte-1 ! top level
  275. qvf1 = 0.5*(grid%moist(i,k,j,P_QV)+grid%moist(i,k,j,P_QV))
  276. qvf2 = 1./(1.+qvf1)
  277. qvf1 = qvf1*qvf2
  278. ! grid%p(i,k,j) = - 0.5*grid%mu_1(i,j)/grid%rdnw(k)
  279. grid%p(i,k,j) = - 0.5*(grid%mu_1(i,j)+qvf1*grid%mub(i,j))/grid%rdnw(k)/qvf2
  280. qvf = 1. + rvovrd*moist(i,k,j,P_QV)
  281. grid%alt(i,k,j) = (r_d/p1000mb)*(grid%t_1(i,k,j)+t0)*qvf* &
  282. (((grid%p(i,k,j)+grid%pb(i,k,j))/p1000mb)**cvpm)
  283. grid%al(i,k,j) = grid%alt(i,k,j) - grid%alb(i,k,j)
  284. ! down the column
  285. do k=kte-2,1,-1
  286. qvf1 = 0.5*(moist(i,k,j,P_QV)+moist(i,k+1,j,P_QV))
  287. qvf2 = 1./(1.+qvf1)
  288. qvf1 = qvf1*qvf2
  289. grid%p(i,k,j) = grid%p(i,k+1,j) - (grid%mu_1(i,j) + qvf1*grid%mub(i,j))/qvf2/grid%rdn(k+1)
  290. qvf = 1. + rvovrd*moist(i,k,j,P_QV)
  291. grid%alt(i,k,j) = (r_d/p1000mb)*(grid%t_1(i,k,j)+t0)*qvf* &
  292. (((grid%p(i,k,j)+grid%pb(i,k,j))/p1000mb)**cvpm)
  293. grid%al(i,k,j) = grid%alt(i,k,j) - grid%alb(i,k,j)
  294. enddo
  295. ! this is the hydrostatic equation used in the model after the
  296. ! small timesteps. In the model, al (inverse density)
  297. ! is computed from the geopotential.
  298. grid%ph_1(i,1,j) = 0.
  299. DO k = 2,kte
  300. grid%ph_1(i,k,j) = grid%ph_1(i,k-1,j) - (1./grid%rdnw(k-1))*( &
  301. (grid%mub(i,j)+grid%mu_1(i,j))*grid%al(i,k-1,j)+ &
  302. grid%mu_1(i,j)*grid%alb(i,k-1,j) )
  303. grid%ph_2(i,k,j) = grid%ph_1(i,k,j)
  304. grid%ph0(i,k,j) = grid%ph_1(i,k,j) + grid%phb(i,k,j)
  305. ENDDO
  306. if((i==2) .and. (j==2)) then
  307. write(6,*) ' ph_1 calc ',grid%ph_1(2,1,2),grid%ph_1(2,2,2),&
  308. grid%mu_1(2,2)+grid%mub(2,2),grid%mu_1(2,2), &
  309. grid%alb(2,1,2),grid%al(1,2,1),grid%rdnw(1)
  310. endif
  311. ENDDO
  312. ENDDO
  313. ! cold bubble input (from straka et al, IJNMF, vol 17, 1993 pp 1-22)
  314. t_min = grid%t_1(its,kts,jts)
  315. t_max = t_min
  316. u_mean = 00.
  317. xpos = config_flags%dx*nxc - u_mean*900.
  318. xposml = xpos - config_flags%dx*(ide-1)
  319. xpospl = xpos + config_flags%dx*(ide-1)
  320. DO J = jts, min(jde-1,jte)
  321. DO I = its, min(ide-1,ite)
  322. ! xrad = config_flags%dx*float(i-nxc)/4000. ! 4000 meter horizontal radius
  323. ! ! centered in the domain
  324. xrad = min( abs(config_flags%dx*float(i)-xpos), &
  325. abs(config_flags%dx*float(i)-xposml), &
  326. abs(config_flags%dx*float(i)-xpospl))/4000.
  327. DO K = 1, kte-1
  328. ! put in preturbation theta (bubble) and recalc density. note,
  329. ! the mass in the column is not changing, so when theta changes,
  330. ! we recompute density and geopotential
  331. zrad = 0.5*(grid%ph_1(i,k,j)+grid%ph_1(i,k+1,j) &
  332. +grid%phb(i,k,j)+grid%phb(i,k+1,j))/g
  333. zrad = (zrad-3000.)/2000. ! 2000 meter vertical radius,
  334. ! centered at z=3000,
  335. RAD=SQRT(xrad*xrad+zrad*zrad)
  336. IF(RAD <= 1.) THEN
  337. ! perturbation temperature is 15 C, convert to potential temperature
  338. delt = -15.0 / ((grid%p(i,k,j)+grid%pb(i,k,j))/p1000mb)**rcp
  339. grid%T_1(i,k,j)=grid%T_1(i,k,j)+delt*(COS(PI*RAD)+1.0)/2.
  340. grid%T_2(i,k,j)=grid%T_1(i,k,j)
  341. qvf = 1. + rvovrd*moist(i,k,j,P_QV)
  342. grid%alt(i,k,j) = (r_d/p1000mb)*(grid%t_1(i,k,j)+t0)*qvf* &
  343. (((grid%p(i,k,j)+grid%pb(i,k,j))/p1000mb)**cvpm)
  344. grid%al(i,k,j) = grid%alt(i,k,j) - grid%alb(i,k,j)
  345. ENDIF
  346. t_min = min(t_min, grid%t_1(i,k,j))
  347. t_max = max(t_max, grid%t_1(i,k,j))
  348. ENDDO
  349. ! rebalance hydrostatically
  350. DO k = 2,kte
  351. grid%ph_1(i,k,j) = grid%ph_1(i,k-1,j) - (1./grid%rdnw(k-1))*( &
  352. (grid%mub(i,j)+grid%mu_1(i,j))*grid%al(i,k-1,j)+ &
  353. grid%mu_1(i,j)*grid%alb(i,k-1,j) )
  354. grid%ph_2(i,k,j) = grid%ph_1(i,k,j)
  355. grid%ph0(i,k,j) = grid%ph_1(i,k,j) + grid%phb(i,k,j)
  356. ENDDO
  357. ENDDO
  358. ENDDO
  359. write(6,*) ' min and max theta perturbation ',t_min,t_max
  360. ! -- end bubble insert
  361. write(6,*) ' mu_1 from comp ', grid%mu_1(1,1)
  362. write(6,*) ' full state sounding from comp, ph, p, al, t_1, qv '
  363. do k=1,kde-1
  364. write(6,'(i3,1x,5(1x,1pe10.3))') k, grid%ph_1(1,k,1)+grid%phb(1,k,1), &
  365. grid%p(1,k,1)+grid%pb(1,k,1), grid%alt(1,k,1), &
  366. grid%t_1(1,k,1)+t0, moist(1,k,1,P_QV)
  367. enddo
  368. write(6,*) ' pert state sounding from comp, ph_1, pp, alp, t_1, qv '
  369. do k=1,kde-1
  370. write(6,'(i3,1x,5(1x,1pe10.3))') k, grid%ph_1(1,k,1), &
  371. grid%p(1,k,1), grid%al(1,k,1), &
  372. grid%t_1(1,k,1), moist(1,k,1,P_QV)
  373. enddo
  374. write(6,*) ' '
  375. write(6,*) ' k, model level, dz '
  376. do k=1,kde-1
  377. write(6,'(i3,1x,e12.5,1x,f10.2)') k, &
  378. .5*(grid%ph_1(1,k,1)+grid%phb(1,k,1)+grid%ph_1(1,k+1,1)+grid%phb(1,k+1,1))/g, &
  379. (grid%ph_1(1,k+1,1)+grid%phb(1,k+1,1)-grid%ph_1(1,k,1)-grid%phb(1,k,1))/g
  380. enddo
  381. write(6,*) ' model top (m) is ', (grid%ph_1(1,kde,1)+grid%phb(1,kde,1))/g
  382. ! interp v
  383. DO J = jts, jte
  384. DO I = its, min(ide-1,ite)
  385. IF (j == jds) THEN
  386. z_at_v = grid%phb(i,1,j)/g
  387. ELSE IF (j == jde) THEN
  388. z_at_v = grid%phb(i,1,j-1)/g
  389. ELSE
  390. z_at_v = 0.5*(grid%phb(i,1,j)+grid%phb(i,1,j-1))/g
  391. END IF
  392. p_surf = interp_0( p_in, zk, z_at_v, nl_in )
  393. DO K = 1, kte
  394. p_level = grid%znu(k)*(p_surf - grid%p_top) + grid%p_top
  395. grid%v_1(i,k,j) = interp_0( v, p_in, p_level, nl_in )
  396. grid%v_2(i,k,j) = grid%v_1(i,k,j)
  397. ENDDO
  398. ENDDO
  399. ENDDO
  400. ! interp u
  401. DO J = jts, min(jde-1,jte)
  402. DO I = its, ite
  403. IF (i == ids) THEN
  404. z_at_u = grid%phb(i,1,j)/g
  405. ELSE IF (i == ide) THEN
  406. z_at_u = grid%phb(i-1,1,j)/g
  407. ELSE
  408. z_at_u = 0.5*(grid%phb(i,1,j)+grid%phb(i-1,1,j))/g
  409. END IF
  410. p_surf = interp_0( p_in, zk, z_at_u, nl_in )
  411. DO K = 1, kte
  412. p_level = grid%znu(k)*(p_surf - grid%p_top) + grid%p_top
  413. grid%u_1(i,k,j) = interp_0( u, p_in, p_level, nl_in )
  414. grid%u_2(i,k,j) = grid%u_1(i,k,j)
  415. ENDDO
  416. ENDDO
  417. ENDDO
  418. ! set w
  419. DO J = jts, min(jde-1,jte)
  420. DO K = kts, kte
  421. DO I = its, min(ide-1,ite)
  422. grid%w_1(i,k,j) = 0.
  423. grid%w_2(i,k,j) = 0.
  424. ENDDO
  425. ENDDO
  426. ENDDO
  427. ! set a few more things
  428. DO J = jts, min(jde-1,jte)
  429. DO K = kts, kte-1
  430. DO I = its, min(ide-1,ite)
  431. grid%h_diabatic(i,k,j) = 0.
  432. ENDDO
  433. ENDDO
  434. ENDDO
  435. DO k=1,kte-1
  436. grid%t_base(k) = grid%t_1(1,k,1)
  437. grid%qv_base(k) = moist(1,k,1,P_QV)
  438. grid%u_base(k) = grid%u_1(1,k,1)
  439. grid%v_base(k) = grid%v_1(1,k,1)
  440. grid%z_base(k) = 0.5*(grid%phb(1,k,1)+grid%phb(1,k+1,1)+grid%ph_1(1,k,1)+grid%ph_1(1,k+1,1))/g
  441. ENDDO
  442. DO J = jts, min(jde-1,jte)
  443. DO I = its, min(ide-1,ite)
  444. thtmp = grid%t_2(i,1,j)+t0
  445. ptmp = grid%p(i,1,j)+grid%pb(i,1,j)
  446. temp(1) = thtmp * (ptmp/p1000mb)**rcp
  447. thtmp = grid%t_2(i,2,j)+t0
  448. ptmp = grid%p(i,2,j)+grid%pb(i,2,j)
  449. temp(2) = thtmp * (ptmp/p1000mb)**rcp
  450. thtmp = grid%t_2(i,3,j)+t0
  451. ptmp = grid%p(i,3,j)+grid%pb(i,3,j)
  452. temp(3) = thtmp * (ptmp/p1000mb)**rcp
  453. grid%TSK(I,J)=grid%cf1*temp(1)+grid%cf2*temp(2)+grid%cf3*temp(3)
  454. grid%TMN(I,J)=grid%TSK(I,J)-0.5
  455. ENDDO
  456. ENDDO
  457. RETURN
  458. END SUBROUTINE init_domain_rk
  459. SUBROUTINE init_module_initialize
  460. END SUBROUTINE init_module_initialize
  461. !---------------------------------------------------------------------
  462. ! test driver for get_sounding
  463. !
  464. ! implicit none
  465. ! integer n
  466. ! parameter(n = 1000)
  467. ! real zk(n),p(n),theta(n),rho(n),u(n),v(n),qv(n),pd(n)
  468. ! logical dry
  469. ! integer nl,k
  470. !
  471. ! dry = .false.
  472. ! dry = .true.
  473. ! call get_sounding( zk, p, pd, theta, rho, u, v, qv, dry, n, nl )
  474. ! write(6,*) ' input levels ',nl
  475. ! write(6,*) ' sounding '
  476. ! write(6,*) ' k height(m) press (Pa) pd(Pa) theta (K) den(kg/m^3) u(m/s) v(m/s) qv(g/g) '
  477. ! do k=1,nl
  478. ! write(6,'(1x,i3,8(1x,1pe10.3))') k, zk(k), p(k), pd(k), theta(k), rho(k), u(k), v(k), qv(k)
  479. ! enddo
  480. ! end
  481. !
  482. !---------------------------------------------------------------------------
  483. subroutine get_sounding( zk, p, p_dry, theta, rho, &
  484. u, v, qv, dry, nl_max, nl_in, base_state )
  485. implicit none
  486. integer nl_max, nl_in
  487. real zk(nl_max), p(nl_max), theta(nl_max), rho(nl_max), &
  488. u(nl_max), v(nl_max), qv(nl_max), p_dry(nl_max)
  489. logical dry
  490. logical base_state
  491. integer n, iz
  492. parameter(n=1000)
  493. logical debug
  494. parameter( debug = .false.)
  495. ! input sounding data
  496. real p_surf, th_surf, qv_surf
  497. real pi_surf, pi(n)
  498. real h_input(n), th_input(n), qv_input(n), u_input(n), v_input(n)
  499. ! diagnostics
  500. real rho_surf, p_input(n), rho_input(n)
  501. real pm_input(n) ! this are for full moist sounding
  502. ! local data
  503. real r
  504. parameter (r = r_d)
  505. integer k, it, nl
  506. real qvf, qvf1, dz
  507. ! first, read the sounding
  508. call read_sounding( p_surf, th_surf, qv_surf, &
  509. h_input, th_input, qv_input, u_input, v_input,n, nl, debug )
  510. ! iz = 1
  511. ! do k=2,nl
  512. ! if(h_input(k) .lt. 12000.) iz = k
  513. ! enddo
  514. ! write(6,*) " tropopause ",iz,h_input(iz)
  515. ! if(dry) then
  516. ! write(6,*) ' nl is ',nl
  517. ! do k=1,nl
  518. ! th_input(k) = th_input(k)+10.+10*float(k)/nl
  519. ! enddo
  520. ! write(6,*) ' finished adjusting theta '
  521. ! endif
  522. ! do k=1,nl
  523. ! u_input(k) = 2*u_input(k)
  524. ! enddo
  525. !
  526. ! end if
  527. if(dry) then
  528. do k=1,nl
  529. qv_input(k) = 0.
  530. enddo
  531. endif
  532. if(debug) write(6,*) ' number of input levels = ',nl
  533. nl_in = nl
  534. if(nl_in .gt. nl_max ) then
  535. write(6,*) ' too many levels for input arrays ',nl_in,nl_max
  536. call wrf_error_fatal ( ' too many levels for input arrays ' )
  537. end if
  538. ! compute diagnostics,
  539. ! first, convert qv(g/kg) to qv(g/g)
  540. do k=1,nl
  541. qv_input(k) = 0.001*qv_input(k)
  542. enddo
  543. p_surf = 100.*p_surf ! convert to pascals
  544. qvf = 1. + rvovrd*qv_input(1)
  545. rho_surf = 1./((r/p1000mb)*th_surf*qvf*((p_surf/p1000mb)**cvpm))
  546. pi_surf = (p_surf/p1000mb)**(r/cp)
  547. if(debug) then
  548. write(6,*) ' surface density is ',rho_surf
  549. write(6,*) ' surface pi is ',pi_surf
  550. end if
  551. ! integrate moist sounding hydrostatically, starting from the
  552. ! specified surface pressure
  553. ! -> first, integrate from surface to lowest level
  554. qvf = 1. + rvovrd*qv_input(1)
  555. qvf1 = 1. + qv_input(1)
  556. rho_input(1) = rho_surf
  557. dz = h_input(1)
  558. do it=1,10
  559. pm_input(1) = p_surf &
  560. - 0.5*dz*(rho_surf+rho_input(1))*g*qvf1
  561. rho_input(1) = 1./((r/p1000mb)*th_input(1)*qvf*((pm_input(1)/p1000mb)**cvpm))
  562. enddo
  563. ! integrate up the column
  564. do k=2,nl
  565. rho_input(k) = rho_input(k-1)
  566. dz = h_input(k)-h_input(k-1)
  567. qvf1 = 0.5*(2.+(qv_input(k-1)+qv_input(k)))
  568. qvf = 1. + rvovrd*qv_input(k) ! qv is in g/kg here
  569. do it=1,20
  570. pm_input(k) = pm_input(k-1) &
  571. - 0.5*dz*(rho_input(k)+rho_input(k-1))*g*qvf1
  572. rho_input(k) = 1./((r/p1000mb)*th_input(k)*qvf*((pm_input(k)/p1000mb)**cvpm))
  573. enddo
  574. enddo
  575. ! we have the moist sounding
  576. ! next, compute the dry sounding using p at the highest level from the
  577. ! moist sounding and integrating down.
  578. p_input(nl) = pm_input(nl)
  579. do k=nl-1,1,-1
  580. dz = h_input(k+1)-h_input(k)
  581. p_input(k) = p_input(k+1) + 0.5*dz*(rho_input(k)+rho_input(k+1))*g
  582. enddo
  583. ! write(6,*) ' zeroing u input '
  584. do k=1,nl
  585. zk(k) = h_input(k)
  586. p(k) = pm_input(k)
  587. p_dry(k) = p_input(k)
  588. theta(k) = th_input(k)
  589. rho(k) = rho_input(k)
  590. u(k) = u_input(k)
  591. ! u(k) = 0.
  592. v(k) = v_input(k)
  593. qv(k) = qv_input(k)
  594. enddo
  595. if(debug) then
  596. write(6,*) ' sounding '
  597. write(6,*) ' k height(m) press (Pa) pd(Pa) theta (K) den(kg/m^3) u(m/s) v(m/s) qv(g/g) '
  598. do k=1,nl
  599. write(6,'(1x,i3,8(1x,1pe10.3))') k, zk(k), p(k), p_dry(k), theta(k), rho(k), u(k), v(k), qv(k)
  600. enddo
  601. end if
  602. end subroutine get_sounding
  603. !-------------------------------------------------------
  604. subroutine read_sounding( ps,ts,qvs,h,th,qv,u,v,n,nl,debug )
  605. implicit none
  606. integer n,nl
  607. real ps,ts,qvs,h(n),th(n),qv(n),u(n),v(n)
  608. logical end_of_file
  609. logical debug
  610. integer k
  611. open(unit=10,file='input_sounding',form='formatted',status='old')
  612. rewind(10)
  613. read(10,*) ps, ts, qvs
  614. if(debug) then
  615. write(6,*) ' input sounding surface parameters '
  616. write(6,*) ' surface pressure (mb) ',ps
  617. write(6,*) ' surface pot. temp (K) ',ts
  618. write(6,*) ' surface mixing ratio (g/kg) ',qvs
  619. end if
  620. end_of_file = .false.
  621. k = 0
  622. do while (.not. end_of_file)
  623. read(10,*,end=100) h(k+1), th(k+1), qv(k+1), u(k+1), v(k+1)
  624. k = k+1
  625. if(debug) write(6,'(1x,i3,5(1x,e10.3))') k, h(k), th(k), qv(k), u(k), v(k)
  626. go to 110
  627. 100 end_of_file = .true.
  628. 110 continue
  629. enddo
  630. nl = k
  631. close(unit=10,status = 'keep')
  632. end subroutine read_sounding
  633. END MODULE module_initialize_ideal