1 !IDEAL:MODEL_LAYER:INITIALIZATION
4 ! This MODULE holds the routines which are used to perform various initializations
5 ! for the individual domains.
7 ! This MODULE CONTAINS the following routines:
9 ! initialize_field_test - 1. Set different fields to different constant
10 ! values. This is only a test. If the correct
11 ! domain is not found (based upon the "id")
12 ! then a fatal error is issued.
14 !-----------------------------------------------------------------------
16 MODULE module_initialize_ideal
20 USE module_state_description
21 USE module_model_constants
25 USE module_init_utilities
34 !-------------------------------------------------------------------
35 ! this is a wrapper for the solver-specific init_domain routines.
36 ! Also dereferences the grid variables and passes them down as arguments.
37 ! This is crucial, since the lower level routines may do message passing
38 ! and this will get fouled up on machines that insist on passing down
39 ! copies of assumed-shape arrays (by passing down as arguments, the
40 ! data are treated as assumed-size -- ie. f77 -- arrays and the copying
41 ! business is avoided). Fie on the F90 designers. Fie and a pox.
43 SUBROUTINE init_domain ( grid )
48 TYPE (domain), POINTER :: grid
50 INTEGER :: idum1, idum2
52 CALL set_scalar_indices_from_config ( head_grid%id , idum1, idum2 )
54 CALL init_domain_rk( grid &
56 #include <actual_new_args.inc>
60 END SUBROUTINE init_domain
62 !-------------------------------------------------------------------
64 SUBROUTINE init_domain_rk ( grid &
66 # include <dummy_new_args.inc>
72 TYPE (domain), POINTER :: grid
74 # include <dummy_new_decl.inc>
76 TYPE (grid_config_rec_type) :: config_flags
80 ids, ide, jds, jde, kds, kde, &
81 ims, ime, jms, jme, kms, kme, &
82 its, ite, jts, jte, kts, kte, &
87 INTEGER, PARAMETER :: nl_max = 1000
88 REAL, DIMENSION(nl_max) :: zk, p_in, theta, rho, u, v, qv, pd_in
92 INTEGER :: icm,jcm, ii, im1, jj, jm1, loop, error, fid, nxc, nyc
93 REAL :: u_mean,v_mean, f0, p_surf, p_level, qvf, z_at_v, z_at_u
94 REAL :: z_scale, xrad, yrad, zrad, rad, delt, cof1, cof2
95 ! REAL, EXTERNAL :: interp_0
99 ! stuff from original initialization that has been dropped from the Registry
100 REAL :: vnu, xnu, xnus, dinit0, cbh, p0_temp, t0_temp, zd, zt
101 REAL :: qvf1, qvf2, pd_surf
103 real :: thtmp, ptmp, temp(3)
105 LOGICAL :: moisture_init
106 LOGICAL :: stretch_grid, dry_sounding
108 INTEGER :: xs , xe , ys , ye
110 LOGICAL, EXTERNAL :: wrf_dm_on_monitor
112 SELECT CASE ( model_data_order )
113 CASE ( DATA_ORDER_ZXY )
114 kds = grid%sd31 ; kde = grid%ed31 ;
115 ids = grid%sd32 ; ide = grid%ed32 ;
116 jds = grid%sd33 ; jde = grid%ed33 ;
118 kms = grid%sm31 ; kme = grid%em31 ;
119 ims = grid%sm32 ; ime = grid%em32 ;
120 jms = grid%sm33 ; jme = grid%em33 ;
122 kts = grid%sp31 ; kte = grid%ep31 ; ! note that tile is entire patch
123 its = grid%sp32 ; ite = grid%ep32 ; ! note that tile is entire patch
124 jts = grid%sp33 ; jte = grid%ep33 ; ! note that tile is entire patch
125 CASE ( DATA_ORDER_XYZ )
126 ids = grid%sd31 ; ide = grid%ed31 ;
127 jds = grid%sd32 ; jde = grid%ed32 ;
128 kds = grid%sd33 ; kde = grid%ed33 ;
130 ims = grid%sm31 ; ime = grid%em31 ;
131 jms = grid%sm32 ; jme = grid%em32 ;
132 kms = grid%sm33 ; kme = grid%em33 ;
134 its = grid%sp31 ; ite = grid%ep31 ; ! note that tile is entire patch
135 jts = grid%sp32 ; jte = grid%ep32 ; ! note that tile is entire patch
136 kts = grid%sp33 ; kte = grid%ep33 ; ! note that tile is entire patch
137 CASE ( DATA_ORDER_XZY )
138 ids = grid%sd31 ; ide = grid%ed31 ;
139 kds = grid%sd32 ; kde = grid%ed32 ;
140 jds = grid%sd33 ; jde = grid%ed33 ;
142 ims = grid%sm31 ; ime = grid%em31 ;
143 kms = grid%sm32 ; kme = grid%em32 ;
144 jms = grid%sm33 ; jme = grid%em33 ;
146 its = grid%sp31 ; ite = grid%ep31 ; ! note that tile is entire patch
147 kts = grid%sp32 ; kte = grid%ep32 ; ! note that tile is entire patch
148 jts = grid%sp33 ; jte = grid%ep33 ; ! note that tile is entire patch
153 stretch_grid = .true.
158 write(6,*) ' pi is ',pi
162 CALL model_to_grid_config_rec ( grid%id , model_config_rec , config_flags )
164 ! here we check to see if the boundary conditions are set properly
166 CALL boundary_condition_check( config_flags, bdyzone, error, grid%id )
168 moisture_init = .true.
173 CALL wrf_dm_bcast_bytes( icm , IWORDSIZE )
174 CALL wrf_dm_bcast_bytes( jcm , IWORDSIZE )
177 CALL nl_set_mminlu(1, ' ')
178 CALL nl_set_iswater(1,0)
179 CALL nl_set_cen_lat(1,40.)
180 CALL nl_set_cen_lon(1,-105.)
181 CALL nl_set_truelat1(1,0.)
182 CALL nl_set_truelat2(1,0.)
183 CALL nl_set_moad_cen_lat (1,0.)
184 CALL nl_set_stand_lon (1,0.)
185 CALL nl_set_pole_lon (1,0.)
186 CALL nl_set_pole_lat (1,90.)
187 CALL nl_set_map_proj(1,0)
190 ! here we initialize data we currently is not initialized
200 grid%msfvx_inv(i,j)= 1.
222 IF (stretch_grid) THEN ! exponential stretch for eta (nearly constant dz)
224 grid%znw(k) = (exp(-(k-1)/float(kde-1)/z_scale) - exp(-1./z_scale))/ &
225 (1.-exp(-1./z_scale))
229 grid%znw(k) = 1. - float(k-1)/float(kde-1)
234 grid%dnw(k) = grid%znw(k+1) - grid%znw(k)
235 grid%rdnw(k) = 1./grid%dnw(k)
236 grid%znu(k) = 0.5*(grid%znw(k+1)+grid%znw(k))
239 grid%dn(k) = 0.5*(grid%dnw(k)+grid%dnw(k-1))
240 grid%rdn(k) = 1./grid%dn(k)
241 grid%fnp(k) = .5* grid%dnw(k )/grid%dn(k)
242 grid%fnm(k) = .5* grid%dnw(k-1)/grid%dn(k)
245 cof1 = (2.*grid%dn(2)+grid%dn(3))/(grid%dn(2)+grid%dn(3))*grid%dnw(1)/grid%dn(2)
246 cof2 = grid%dn(2) /(grid%dn(2)+grid%dn(3))*grid%dnw(1)/grid%dn(3)
247 grid%cf1 = grid%fnp(2) + cof1
248 grid%cf2 = grid%fnm(2) - cof1 - cof2
251 grid%cfn = (.5*grid%dnw(kde-1)+grid%dn(kde-1))/grid%dn(kde-1)
252 grid%cfn1 = -.5*grid%dnw(kde-1)/grid%dn(kde-1)
253 grid%rdx = 1./config_flags%dx
254 grid%rdy = 1./config_flags%dy
256 ! get the sounding from the ascii sounding file, first get dry sounding and
257 ! calculate base state
259 dry_sounding = .true.
260 IF ( wrf_dm_on_monitor() ) THEN
261 write(6,*) ' getting dry sounding for base state '
263 CALL get_sounding( zk, p_in, pd_in, theta, rho, u, v, qv, dry_sounding, nl_max, nl_in )
265 CALL wrf_dm_bcast_real( zk , nl_max )
266 CALL wrf_dm_bcast_real( p_in , nl_max )
267 CALL wrf_dm_bcast_real( pd_in , nl_max )
268 CALL wrf_dm_bcast_real( theta , nl_max )
269 CALL wrf_dm_bcast_real( rho , nl_max )
270 CALL wrf_dm_bcast_real( u , nl_max )
271 CALL wrf_dm_bcast_real( v , nl_max )
272 CALL wrf_dm_bcast_real( qv , nl_max )
273 CALL wrf_dm_bcast_integer ( nl_in , 1 )
275 write(6,*) ' returned from reading sounding, nl_in is ',nl_in
277 ! find ptop for the desired ztop (ztop is input from the namelist),
278 ! and find surface pressure
280 grid%p_top = interp_0( p_in, zk, config_flags%ztop, nl_in )
295 DO j=max(ys,jds),min(ye,jde-1)
296 DO i=max(xs,ids),min(xe,ide-1)
297 grid%ht(i,j) = mtn_ht * 0.25 * &
298 ( 1. + COS ( 2*pi/(xe-xs) * ( i-xs ) + pi ) ) * &
299 ( 1. + COS ( 2*pi/(ye-ys) * ( j-ys ) + pi ) )
304 DO j=max(ys,jds),min(ye,jde-1)
306 grid%ht(i,j) = mtn_ht * 0.50 * &
307 ( 1. + COS ( 2*pi/(ye-ys) * ( j-ys ) + pi ) )
313 DO i=max(xs,ids),min(xe,ide-1)
314 grid%ht(i,j) = mtn_ht * 0.50 * &
315 ( 1. + COS ( 2*pi/(xe-xs) * ( i-xs ) + pi ) )
321 grid%phb(i,1,j) = g * grid%ht(i,j)
322 grid%ph0(i,1,j) = g * grid%ht(i,j)
329 p_surf = interp_0( p_in, zk, grid%phb(i,1,j)/g, nl_in )
330 grid%mub(i,j) = p_surf-grid%p_top
332 ! this is dry hydrostatic sounding (base state), so given grid%p (coordinate),
333 ! interp theta (from interp) and compute 1/rho from eqn. of state
336 p_level = grid%znu(k)*(p_surf - grid%p_top) + grid%p_top
337 grid%pb(i,k,j) = p_level
338 grid%t_init(i,k,j) = interp_0( theta, p_in, p_level, nl_in ) - t0
339 grid%alb(i,k,j) = (r_d/p1000mb)*(grid%t_init(i,k,j)+t0)*(grid%pb(i,k,j)/p1000mb)**cvpm
342 ! calc hydrostatic balance (alternatively we could interp the geopotential from the
343 ! sounding, but this assures that the base state is in exact hydrostatic balance with
344 ! respect to the model eqns.
347 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)
353 IF ( wrf_dm_on_monitor() ) THEN
354 write(6,*) ' ptop is ',grid%p_top
355 write(6,*) ' base state grid%mub(1,1), p_surf is ',grid%mub(1,1),grid%mub(1,1)+grid%p_top
358 ! calculate full state for each column - this includes moisture.
360 write(6,*) ' getting moist sounding for full state '
361 dry_sounding = .false.
362 CALL get_sounding( zk, p_in, pd_in, theta, rho, u, v, qv, dry_sounding, nl_max, nl_in )
364 DO J = jts, min(jde-1,jte)
365 DO I = its, min(ide-1,ite)
367 ! At this point grid%p_top is already set. find the DRY mass in the column
368 ! by interpolating the DRY pressure.
370 pd_surf = interp_0( pd_in, zk, grid%phb(i,1,j)/g, nl_in )
372 ! compute the perturbation mass and the full mass
374 grid%mu_1(i,j) = pd_surf-grid%p_top - grid%mub(i,j)
375 grid%mu_2(i,j) = grid%mu_1(i,j)
376 grid%mu0(i,j) = grid%mu_1(i,j) + grid%mub(i,j)
378 ! given the dry pressure and coordinate system, interp the potential
383 p_level = grid%znu(k)*(pd_surf - grid%p_top) + grid%p_top
385 moist(i,k,j,P_QV) = interp_0( qv, pd_in, p_level, nl_in )
386 grid%t_1(i,k,j) = interp_0( theta, pd_in, p_level, nl_in ) - t0
387 grid%t_2(i,k,j) = grid%t_1(i,k,j)
392 ! integrate the hydrostatic equation (from the RHS of the bigstep
393 ! vertical momentum equation) down from the top to get grid%p.
394 ! first from the top of the model to the top pressure
396 k = kte-1 ! top level
398 qvf1 = 0.5*(moist(i,k,j,P_QV)+moist(i,k,j,P_QV))
402 ! grid%p(i,k,j) = - 0.5*grid%mu_1(i,j)/grid%rdnw(k)
403 grid%p(i,k,j) = - 0.5*(grid%mu_1(i,j)+qvf1*grid%mub(i,j))/grid%rdnw(k)/qvf2
404 qvf = 1. + rvovrd*moist(i,k,j,P_QV)
405 grid%alt(i,k,j) = (r_d/p1000mb)*(grid%t_1(i,k,j)+t0)*qvf* &
406 (((grid%p(i,k,j)+grid%pb(i,k,j))/p1000mb)**cvpm)
407 grid%al(i,k,j) = grid%alt(i,k,j) - grid%alb(i,k,j)
412 qvf1 = 0.5*(moist(i,k,j,P_QV)+moist(i,k+1,j,P_QV))
415 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)
416 qvf = 1. + rvovrd*moist(i,k,j,P_QV)
417 grid%alt(i,k,j) = (r_d/p1000mb)*(grid%t_1(i,k,j)+t0)*qvf* &
418 (((grid%p(i,k,j)+grid%pb(i,k,j))/p1000mb)**cvpm)
419 grid%al(i,k,j) = grid%alt(i,k,j) - grid%alb(i,k,j)
422 ! this is the hydrostatic equation used in the model after the
423 ! small timesteps. In the model, grid%al (inverse density)
424 ! is computed from the geopotential.
427 grid%ph_1(i,1,j) = 0.
429 grid%ph_1(i,k,j) = grid%ph_1(i,k-1,j) - (1./grid%rdnw(k-1))*( &
430 (grid%mub(i,j)+grid%mu_1(i,j))*grid%al(i,k-1,j)+ &
431 grid%mu_1(i,j)*grid%alb(i,k-1,j) )
433 grid%ph_2(i,k,j) = grid%ph_1(i,k,j)
434 grid%ph0(i,k,j) = grid%ph_1(i,k,j) + grid%phb(i,k,j)
437 IF ( wrf_dm_on_monitor() ) THEN
438 if((i==2) .and. (j==2)) then
439 write(6,*) ' grid%ph_1 calc ',grid%ph_1(2,1,2),grid%ph_1(2,2,2),&
440 grid%mu_1(2,2)+grid%mub(2,2),grid%mu_1(2,2), &
441 grid%alb(2,1,2),grid%al(1,2,1),grid%rdnw(1)
450 ! thermal perturbation to kick off convection
452 write(6,*) ' nxc, nyc for perturbation ',nxc,nyc
453 write(6,*) ' delt for perturbation ',delt
455 DO J = jts, min(jde-1,jte)
456 yrad = config_flags%dy*float(j-nyc)/10000.
458 DO I = its, min(ide-1,ite)
459 xrad = config_flags%dx*float(i-nxc)/10000.
463 ! put in preturbation theta (bubble) and recalc density. note,
464 ! the mass in the column is not changing, so when theta changes,
465 ! we recompute density and geopotential
467 zrad = 0.5*(grid%ph_1(i,k,j)+grid%ph_1(i,k+1,j) &
468 +grid%phb(i,k,j)+grid%phb(i,k+1,j))/g
469 zrad = (zrad-1500.)/1500.
470 RAD=SQRT(xrad*xrad+yrad*yrad+zrad*zrad)
472 grid%t_1(i,k,j)=grid%t_1(i,k,j)+delt*COS(.5*PI*RAD)**2
473 grid%t_2(i,k,j)=grid%t_1(i,k,j)
474 qvf = 1. + rvovrd*moist(i,k,j,P_QV)
475 grid%alt(i,k,j) = (r_d/p1000mb)*(grid%t_1(i,k,j)+t0)*qvf* &
476 (((grid%p(i,k,j)+grid%pb(i,k,j))/p1000mb)**cvpm)
477 grid%al(i,k,j) = grid%alt(i,k,j) - grid%alb(i,k,j)
481 ! rebalance hydrostatically
484 grid%ph_1(i,k,j) = grid%ph_1(i,k-1,j) - (1./grid%rdnw(k-1))*( &
485 (grid%mub(i,j)+grid%mu_1(i,j))*grid%al(i,k-1,j)+ &
486 grid%mu_1(i,j)*grid%alb(i,k-1,j) )
488 grid%ph_2(i,k,j) = grid%ph_1(i,k,j)
489 grid%ph0(i,k,j) = grid%ph_1(i,k,j) + grid%phb(i,k,j)
497 IF ( wrf_dm_on_monitor() ) THEN
498 write(6,*) ' grid%mu_1 from comp ', grid%mu_1(1,1)
499 write(6,*) ' full state sounding from comp, ph, grid%p, grid%al, grid%t_1, qv '
501 write(6,'(i3,1x,5(1x,1pe10.3))') k, grid%ph_1(1,k,1)+grid%phb(1,k,1), &
502 grid%p(1,k,1)+grid%pb(1,k,1), grid%alt(1,k,1), &
503 grid%t_1(1,k,1)+t0, moist(1,k,1,P_QV)
506 write(6,*) ' pert state sounding from comp, grid%ph_1, pp, alp, grid%t_1, qv '
508 write(6,'(i3,1x,5(1x,1pe10.3))') k, grid%ph_1(1,k,1), &
509 grid%p(1,k,1), grid%al(1,k,1), &
510 grid%t_1(1,k,1), moist(1,k,1,P_QV)
517 DO I = its, min(ide-1,ite)
520 z_at_v = grid%phb(i,1,j)/g
521 ELSE IF (j == jde) THEN
522 z_at_v = grid%phb(i,1,j-1)/g
524 z_at_v = 0.5*(grid%phb(i,1,j)+grid%phb(i,1,j-1))/g
526 p_surf = interp_0( p_in, zk, z_at_v, nl_in )
529 p_level = grid%znu(k)*(p_surf - grid%p_top) + grid%p_top
530 grid%v_1(i,k,j) = interp_0( v, p_in, p_level, nl_in )
531 grid%v_2(i,k,j) = grid%v_1(i,k,j)
539 DO J = jts, min(jde-1,jte)
543 z_at_u = grid%phb(i,1,j)/g
544 ELSE IF (i == ide) THEN
545 z_at_u = grid%phb(i-1,1,j)/g
547 z_at_u = 0.5*(grid%phb(i,1,j)+grid%phb(i-1,1,j))/g
550 p_surf = interp_0( p_in, zk, z_at_u, nl_in )
553 p_level = grid%znu(k)*(p_surf - grid%p_top) + grid%p_top
554 grid%u_1(i,k,j) = interp_0( u, p_in, p_level, nl_in )
555 grid%u_2(i,k,j) = grid%u_1(i,k,j)
563 DO J = jts, min(jde-1,jte)
565 DO I = its, min(ide-1,ite)
572 ! set a few more things
574 DO J = jts, min(jde-1,jte)
576 DO I = its, min(ide-1,ite)
577 grid%h_diabatic(i,k,j) = 0.
582 IF ( wrf_dm_on_monitor() ) THEN
584 grid%t_base(k) = grid%t_1(1,k,1)
585 grid%qv_base(k) = moist(1,k,1,P_QV)
586 grid%u_base(k) = grid%u_1(1,k,1)
587 grid%v_base(k) = grid%v_1(1,k,1)
588 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
591 CALL wrf_dm_bcast_real( grid%t_base , kte )
592 CALL wrf_dm_bcast_real( grid%qv_base , kte )
593 CALL wrf_dm_bcast_real( grid%u_base , kte )
594 CALL wrf_dm_bcast_real( grid%v_base , kte )
595 CALL wrf_dm_bcast_real( grid%z_base , kte )
597 DO J = jts, min(jde-1,jte)
598 DO I = its, min(ide-1,ite)
599 thtmp = grid%t_2(i,1,j)+t0
600 ptmp = grid%p(i,1,j)+grid%pb(i,1,j)
601 temp(1) = thtmp * (ptmp/p1000mb)**rcp
602 thtmp = grid%t_2(i,2,j)+t0
603 ptmp = grid%p(i,2,j)+grid%pb(i,2,j)
604 temp(2) = thtmp * (ptmp/p1000mb)**rcp
605 thtmp = grid%t_2(i,3,j)+t0
606 ptmp = grid%p(i,3,j)+grid%pb(i,3,j)
607 temp(3) = thtmp * (ptmp/p1000mb)**rcp
609 grid%tsk(I,J)=grid%cf1*temp(1)+grid%cf2*temp(2)+grid%cf3*temp(3)
610 grid%tmn(I,J)=grid%tsk(I,J)-0.5
614 END SUBROUTINE init_domain_rk
616 SUBROUTINE init_module_initialize
617 END SUBROUTINE init_module_initialize
619 !---------------------------------------------------------------------
621 ! test driver for get_sounding
625 ! parameter(n = 1000)
626 ! real zk(n),p(n),theta(n),rho(n),u(n),v(n),qv(n),pd(n)
632 ! call get_sounding( zk, p, pd, theta, rho, u, v, qv, dry, n, nl )
633 ! write(6,*) ' input levels ',nl
634 ! write(6,*) ' sounding '
635 ! write(6,*) ' k height(m) press (Pa) pd(Pa) theta (K) den(kg/m^3) u(m/s) v(m/s) qv(g/g) '
637 ! 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)
641 !---------------------------------------------------------------------------
643 subroutine get_sounding( zk, p, p_dry, theta, rho, &
644 u, v, qv, dry, nl_max, nl_in )
647 integer nl_max, nl_in
648 real zk(nl_max), p(nl_max), theta(nl_max), rho(nl_max), &
649 u(nl_max), v(nl_max), qv(nl_max), p_dry(nl_max)
655 parameter( debug = .true.)
657 ! input sounding data
659 real p_surf, th_surf, qv_surf
661 real h_input(n), th_input(n), qv_input(n), u_input(n), v_input(n)
665 real rho_surf, p_input(n), rho_input(n)
666 real pm_input(n) ! this are for full moist sounding
675 ! first, read the sounding
677 call read_sounding( p_surf, th_surf, qv_surf, &
678 h_input, th_input, qv_input, u_input, v_input,n, nl, debug )
686 if(debug) write(6,*) ' number of input levels = ',nl
689 if(nl_in .gt. nl_max ) then
690 write(6,*) ' too many levels for input arrays ',nl_in,nl_max
691 call wrf_error_fatal ( ' too many levels for input arrays ' )
694 ! compute diagnostics,
695 ! first, convert qv(g/kg) to qv(g/g)
698 qv_input(k) = 0.001*qv_input(k)
701 p_surf = 100.*p_surf ! convert to pascals
702 qvf = 1. + rvovrd*qv_input(1)
703 rho_surf = 1./((r/p1000mb)*th_surf*qvf*((p_surf/p1000mb)**cvpm))
704 pi_surf = (p_surf/p1000mb)**(r/cp)
707 write(6,*) ' surface density is ',rho_surf
708 write(6,*) ' surface pi is ',pi_surf
712 ! integrate moist sounding hydrostatically, starting from the
713 ! specified surface pressure
714 ! -> first, integrate from surface to lowest level
716 qvf = 1. + rvovrd*qv_input(1)
717 qvf1 = 1. + qv_input(1)
718 rho_input(1) = rho_surf
721 pm_input(1) = p_surf &
722 - 0.5*dz*(rho_surf+rho_input(1))*g*qvf1
723 rho_input(1) = 1./((r/p1000mb)*th_input(1)*qvf*((pm_input(1)/p1000mb)**cvpm))
726 ! integrate up the column
729 rho_input(k) = rho_input(k-1)
730 dz = h_input(k)-h_input(k-1)
731 qvf1 = 0.5*(2.+(qv_input(k-1)+qv_input(k)))
732 qvf = 1. + rvovrd*qv_input(k) ! qv is in g/kg here
735 pm_input(k) = pm_input(k-1) &
736 - 0.5*dz*(rho_input(k)+rho_input(k-1))*g*qvf1
737 rho_input(k) = 1./((r/p1000mb)*th_input(k)*qvf*((pm_input(k)/p1000mb)**cvpm))
741 ! we have the moist sounding
743 ! next, compute the dry sounding using p at the highest level from the
744 ! moist sounding and integrating down.
746 p_input(nl) = pm_input(nl)
749 dz = h_input(k+1)-h_input(k)
750 p_input(k) = p_input(k+1) + 0.5*dz*(rho_input(k)+rho_input(k+1))*g
758 p_dry(k) = p_input(k)
759 theta(k) = th_input(k)
760 rho(k) = rho_input(k)
768 write(6,*) ' sounding '
769 write(6,*) ' k height(m) press (Pa) pd(Pa) theta (K) den(kg/m^3) u(m/s) v(m/s) qv(g/g) '
771 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)
776 end subroutine get_sounding
778 !-------------------------------------------------------
780 subroutine read_sounding( ps,ts,qvs,h,th,qv,u,v,n,nl,debug )
783 real ps,ts,qvs,h(n),th(n),qv(n),u(n),v(n)
789 open(unit=10,file='input_sounding',form='formatted',status='old')
791 read(10,*) ps, ts, qvs
793 write(6,*) ' input sounding surface parameters '
794 write(6,*) ' surface pressure (mb) ',ps
795 write(6,*) ' surface pot. temp (K) ',ts
796 write(6,*) ' surface mixing ratio (g/kg) ',qvs
799 end_of_file = .false.
802 do while (.not. end_of_file)
804 read(10,*,end=100) h(k+1), th(k+1), qv(k+1), u(k+1), v(k+1)
806 if(debug) write(6,'(1x,i3,5(1x,e10.3))') k, h(k), th(k), qv(k), u(k), v(k)
808 100 end_of_file = .true.
814 close(unit=10,status = 'keep')
816 end subroutine read_sounding
818 END MODULE module_initialize_ideal