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
35 !-------------------------------------------------------------------
36 ! this is a wrapper for the solver-specific init_domain routines.
37 ! Also dereferences the grid variables and passes them down as arguments.
38 ! This is crucial, since the lower level routines may do message passing
39 ! and this will get fouled up on machines that insist on passing down
40 ! copies of assumed-shape arrays (by passing down as arguments, the
41 ! data are treated as assumed-size -- ie. f77 -- arrays and the copying
42 ! business is avoided). Fie on the F90 designers. Fie and a pox.
43 ! NOTE: Modified to remove all but arrays of rank 4 or more from the
44 ! argument list. Arrays with rank>3 are still problematic due to the
45 ! above-noted fie- and pox-ities. TBH 20061129.
47 SUBROUTINE init_domain ( grid )
52 TYPE (domain), POINTER :: grid
54 INTEGER :: idum1, idum2
56 CALL set_scalar_indices_from_config ( head_grid%id , idum1, idum2 )
58 CALL init_domain_rk( grid &
60 #include <actual_new_args.inc>
63 END SUBROUTINE init_domain
65 !-------------------------------------------------------------------
67 SUBROUTINE init_domain_rk ( grid &
69 # include <dummy_new_args.inc>
73 USE module_optional_input
77 TYPE (domain), POINTER :: grid
79 # include <dummy_new_decl.inc>
81 TYPE (grid_config_rec_type) :: config_flags
85 ids, ide, jds, jde, kds, kde, &
86 ims, ime, jms, jme, kms, kme, &
87 its, ite, jts, jte, kts, kte, &
90 ! JPH should add a read to a config file with:
91 ! ----- check to make sure everything is initialized from the LU index, etc.
92 ! ----- need to make a dummy category?
98 INTEGER, PARAMETER :: nl_max = 1000
99 REAL, DIMENSION(nl_max) :: zk, p_in, theta, rho, u, v, qv, pd_in
102 INTEGER :: ii, im1, jj, jm1, loop, error, fid, lm
103 REAL :: u_mean,v_mean, f0, p_surf, p_level, qvf, z_at_v, z_at_u
104 REAL :: xrad, yrad, zrad, rad, cof1, cof2
107 ! stuff from original initialization that has been dropped from the Registry
108 REAL :: vnu, xnu, xnus, dinit0, cbh, p0_temp, t0_temp, zd
109 REAL :: qvf1, qvf2, pd_surf
111 real :: thtmp, ptmp, temp(3)
114 LOGICAL :: moisture_init
115 LOGICAL :: dry_sounding
119 REAL :: tmn_input, tsk_input
120 REAL :: zs_input(100),tslb_input(100),smois_input(100)
121 LOGICAL :: real_soil = .true.
123 REAL :: zrwa(200), zwa(200)
127 # include <data_calls.inc>
131 SELECT CASE ( model_data_order )
132 CASE ( DATA_ORDER_ZXY )
133 kds = grid%sd31 ; kde = grid%ed31 ;
134 ids = grid%sd32 ; ide = grid%ed32 ;
135 jds = grid%sd33 ; jde = grid%ed33 ;
137 kms = grid%sm31 ; kme = grid%em31 ;
138 ims = grid%sm32 ; ime = grid%em32 ;
139 jms = grid%sm33 ; jme = grid%em33 ;
141 kts = grid%sp31 ; kte = grid%ep31 ; ! note that tile is entire patch
142 its = grid%sp32 ; ite = grid%ep32 ; ! note that tile is entire patch
143 jts = grid%sp33 ; jte = grid%ep33 ; ! note that tile is entire patch
144 CASE ( DATA_ORDER_XYZ )
145 ids = grid%sd31 ; ide = grid%ed31 ;
146 jds = grid%sd32 ; jde = grid%ed32 ;
147 kds = grid%sd33 ; kde = grid%ed33 ;
149 ims = grid%sm31 ; ime = grid%em31 ;
150 jms = grid%sm32 ; jme = grid%em32 ;
151 kms = grid%sm33 ; kme = grid%em33 ;
153 its = grid%sp31 ; ite = grid%ep31 ; ! note that tile is entire patch
154 jts = grid%sp32 ; jte = grid%ep32 ; ! note that tile is entire patch
155 kts = grid%sp33 ; kte = grid%ep33 ; ! note that tile is entire patch
156 CASE ( DATA_ORDER_XZY )
157 ids = grid%sd31 ; ide = grid%ed31 ;
158 kds = grid%sd32 ; kde = grid%ed32 ;
159 jds = grid%sd33 ; jde = grid%ed33 ;
161 ims = grid%sm31 ; ime = grid%em31 ;
162 kms = grid%sm32 ; kme = grid%em32 ;
163 jms = grid%sm33 ; jme = grid%em33 ;
165 its = grid%sp31 ; ite = grid%ep31 ; ! note that tile is entire patch
166 kts = grid%sp32 ; kte = grid%ep32 ; ! note that tile is entire patch
167 jts = grid%sp33 ; jte = grid%ep33 ; ! note that tile is entire patch
173 write(6,*) ' pi is ',pi
175 CALL model_to_grid_config_rec ( grid%id , model_config_rec , config_flags )
177 ! here we check to see if the boundary conditions are set properly
179 CALL boundary_condition_check( config_flags, bdyzone, error, grid%id )
181 moisture_init = .true.
185 CALL nl_set_mminlu(1, 'USGS')
186 CALL nl_set_iswater(1,16)
187 CALL nl_set_isice(1,3)
188 CALL nl_set_truelat1(1,0.)
189 CALL nl_set_truelat2(1,0.)
190 CALL nl_set_moad_cen_lat (1,0.)
191 CALL nl_set_stand_lon(1,0.)
192 CALL nl_set_map_proj(1,0)
193 ! CALL model_to_grid_config_rec(1,model_config_rec,config_flags)
194 CALL nl_get_iswater(1,grid%iswater)
196 ! here we initialize data that currently is not initialized
209 grid%msfvx_inv(i,j) = 1.
213 grid%e(i,j) = 2.0*EOMEG*cos(config_flags%scm_lat*DEGRAD)
214 grid%f(i,j) = 2.0*EOMEG*sin(config_flags%scm_lat*DEGRAD)
215 grid%xlat(i,j) = config_flags%scm_lat
216 grid%xlong(i,j) = config_flags%scm_lon
218 grid%landmask(i,j) = 1.
219 grid%lu_index(i,j) = config_flags%scm_lu_index
223 ! for Noah LSM, additional variables need to be initialized
225 other_masked_fields : SELECT CASE ( model_config_rec%sf_surface_physics(grid%id) )
231 ! JPH free of snow and ice, and only valid over land
232 DO j = jts , MIN(jde-1,jte)
233 DO i = its , MIN(ide-1,ite)
234 grid%vegfra(i,j) = config_flags%scm_vegfra
235 grid%canwat(i,j) = config_flags%scm_canwat
236 grid%isltyp(i,j) = config_flags%scm_isltyp
237 grid%ivgtyp(i,j) = config_flags%scm_lu_index
245 END SELECT other_masked_fields
249 IF ( real_soil ) THEN ! from input file
251 CALL read_soil(100,ns_input,tmn_input,tsk_input,zs_input,tslb_input,smois_input)
253 CALL init_module_optional_input(grid,config_flags)
254 num_st_levels_input = ns_input
255 num_sm_levels_input = ns_input
256 num_sw_levels_input = ns_input
258 st_levels_input(k) = zs_input(k)*100.0 ! to cm
259 sm_levels_input(k) = zs_input(k)*100.0 ! to cm
260 sw_levels_input(k) = zs_input(k)*100.0 ! to cm
261 st_input(:,k+1,:) = tslb_input(k)
262 sm_input(:,k+1,:) = smois_input(k)
263 sw_input(:,k+1,:) = smois_input(k)
270 flag_soil_layers = 0 ! go ahead and put skin temp in
271 flag_soil_levels = 0 ! go ahead and put skin moisture in
272 flag_sst = 0 ! don't modify for ocean
276 CALL process_soil_real ( grid%tsk , grid%tmn , grid%tavgsfc, &
277 grid%landmask , grid%sst , grid%ht, grid%toposoil, &
278 st_input , sm_input , sw_input , &
279 st_levels_input , sm_levels_input , sw_levels_input , &
280 grid%zs , grid%dzs , grid%tslb , grid%smois , grid%sh2o , &
281 flag_sst , flag_tavgsfc, flag_soilhgt, flag_soil_layers, flag_soil_levels, &
282 ids , ide , jds , jde , kds , kde , &
283 ims , ime , jms , jme , kms , kme , &
284 its , ite , jts , jte , kts , kte , &
285 model_config_rec%sf_surface_physics(grid%id) , &
286 model_config_rec%num_soil_layers , &
287 model_config_rec%real_data_init_type , &
288 num_st_levels_input , num_sm_levels_input , num_sw_levels_input , &
289 num_st_levels_alloc , num_sm_levels_alloc , num_sw_levels_alloc )
292 ! Process the soil; note that there are some things hard-wired into share/module_soil_pre.F
293 CALL process_soil_ideal(grid%xland,grid%xice,grid%vegfra,grid%snow,grid%canwat, &
294 grid%ivgtyp,grid%isltyp,grid%tslb,grid%smois, &
295 grid%tsk,grid%tmn,grid%zs,grid%dzs,model_config_rec%num_soil_layers, &
296 model_config_rec%sf_surface_physics(grid%id), &
297 ids,ide, jds,jde, kds,kde,&
298 ims,ime, jms,jme, kms,kme,&
299 its,ite, jts,jte, kts,kte )
311 ! this is adopted from Wayne Angevine's GABLS case
313 zrwa(kde) = exp((kde-1)/40.)
316 zrwa(k) = exp((k-1)/40.)
317 zwa(k) = (zrwa(k)-1.) * grid%ztop/(zrwa(kde)-1.)
318 grid%znw(k) = 1. - (zwa(k) / zwa(kde))
323 grid%dnw(k) = grid%znw(k+1) - grid%znw(k)
324 grid%rdnw(k) = 1./grid%dnw(k)
325 grid%znu(k) = 0.5*(grid%znw(k+1)+grid%znw(k))
328 grid%dn(k) = 0.5*(grid%dnw(k)+grid%dnw(k-1))
329 grid%rdn(k) = 1./grid%dn(k)
330 grid%fnp(k) = .5* grid%dnw(k )/grid%dn(k)
331 grid%fnm(k) = .5* grid%dnw(k-1)/grid%dn(k)
334 cof1 = (2.*grid%dn(2)+grid%dn(3))/(grid%dn(2)+grid%dn(3))*grid%dnw(1)/grid%dn(2)
335 cof2 = grid%dn(2) /(grid%dn(2)+grid%dn(3))*grid%dnw(1)/grid%dn(3)
336 grid%cf1 = grid%fnp(2) + cof1
337 grid%cf2 = grid%fnm(2) - cof1 - cof2
340 grid%cfn = (.5*grid%dnw(kde-1)+grid%dn(kde-1))/grid%dn(kde-1)
341 grid%cfn1 = -.5*grid%dnw(kde-1)/grid%dn(kde-1)
342 grid%rdx = 1./config_flags%dx
343 grid%rdy = 1./config_flags%dy
345 ! get the sounding from the ascii sounding file, first get dry sounding and
346 ! calculate base state
348 write(6,*) ' getting dry sounding for base state '
349 dry_sounding = .true.
350 CALL get_sounding( zsfc, zk, p_in, pd_in, theta, rho, u, v, qv, dry_sounding, nl_max, nl_in )
352 write(6,*) ' returned from reading sounding, nl_in is ',nl_in
354 ! find ptop for the desired ztop (ztop is input from the namelist),
355 ! and find surface pressure
357 grid%p_top = interp_0( p_in, zk, config_flags%ztop, nl_in )
360 DO i=its,ite ! flat surface
362 grid%phb(i,1,j) = grid%ht(i,j) * g
363 grid%ph0(i,1,j) = grid%ht(i,j) * g
371 p_surf = interp_0( p_in, zk, grid%phb(i,1,j)/g, nl_in )
372 grid%mub(i,j) = p_surf-grid%p_top
374 ! this is dry hydrostatic sounding (base state), so given grid%p (coordinate),
375 ! interp theta (from interp) and compute 1/rho from eqn. of state
378 p_level = grid%znu(k)*(p_surf - grid%p_top) + grid%p_top
379 grid%pb(i,k,j) = p_level
380 grid%t_init(i,k,j) = interp_0( theta, p_in, p_level, nl_in ) - t0
381 grid%alb(i,k,j) = (r_d/p1000mb)*(grid%t_init(i,k,j)+t0)*(grid%pb(i,k,j)/p1000mb)**cvpm
384 ! calc hydrostatic balance (alternatively we could interp the geopotential from the
385 ! sounding, but this assures that the base state is in exact hydrostatic balance with
386 ! respect to the model eqns.
389 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)
395 write(6,*) ' ptop is ',grid%p_top
396 write(6,*) ' base state grid%mub(1,1), p_surf is ',grid%mub(1,1),grid%mub(1,1)+grid%p_top
398 ! calculate full state for each column - this includes moisture.
400 write(6,*) ' getting moist sounding for full state '
401 dry_sounding = .false.
402 CALL get_sounding( zsfc, zk, p_in, pd_in, theta, rho, u, v, qv, dry_sounding, nl_max, nl_in )
404 DO J = jts, min(jde-1,jte)
405 DO I = its, min(ide-1,ite)
407 ! At this point grid%p_top is already set. find the DRY mass in the column
408 ! by interpolating the DRY pressure.
410 pd_surf = interp_0( pd_in, zk, grid%phb(i,1,j)/g, nl_in )
412 ! compute the perturbation mass and the full mass
414 grid%mu_1(i,j) = pd_surf-grid%p_top - grid%mub(i,j)
415 grid%mu_2(i,j) = grid%mu_1(i,j)
416 grid%mu0(i,j) = grid%mu_1(i,j) + grid%mub(i,j)
418 ! given the dry pressure and coordinate system, interp the potential
423 p_level = grid%znu(k)*(pd_surf - grid%p_top) + grid%p_top
425 moist(i,k,j,P_QV) = interp_0( qv, pd_in, p_level, nl_in )
426 grid%t_1(i,k,j) = interp_0( theta, pd_in, p_level, nl_in ) - t0
427 grid%t_2(i,k,j) = grid%t_1(i,k,j)
432 ! integrate the hydrostatic equation (from the RHS of the bigstep
433 ! vertical momentum equation) down from the top to get grid%p.
434 ! first from the top of the model to the top pressure
436 k = kte-1 ! top level
438 qvf1 = 0.5*(moist(i,k,j,P_QV)+moist(i,k,j,P_QV))
442 ! grid%p(i,k,j) = - 0.5*grid%mu_1(i,j)/grid%rdnw(k)
443 grid%p(i,k,j) = - 0.5*(grid%mu_1(i,j)+qvf1*grid%mub(i,j))/grid%rdnw(k)/qvf2
444 qvf = 1. + rvovrd*moist(i,k,j,P_QV)
445 grid%alt(i,k,j) = (r_d/p1000mb)*(grid%t_1(i,k,j)+t0)*qvf* &
446 (((grid%p(i,k,j)+grid%pb(i,k,j))/p1000mb)**cvpm)
447 grid%al(i,k,j) = grid%alt(i,k,j) - grid%alb(i,k,j)
452 qvf1 = 0.5*(moist(i,k,j,P_QV)+moist(i,k+1,j,P_QV))
455 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)
456 qvf = 1. + rvovrd*moist(i,k,j,P_QV)
457 grid%alt(i,k,j) = (r_d/p1000mb)*(grid%t_1(i,k,j)+t0)*qvf* &
458 (((grid%p(i,k,j)+grid%pb(i,k,j))/p1000mb)**cvpm)
459 grid%al(i,k,j) = grid%alt(i,k,j) - grid%alb(i,k,j)
462 ! this is the hydrostatic equation used in the model after the
463 ! small timesteps. In the model, grid%al (inverse density)
464 ! is computed from the geopotential.
467 grid%ph_1(i,1,j) = 0.
469 grid%ph_1(i,k,j) = grid%ph_1(i,k-1,j) - (1./grid%rdnw(k-1))*( &
470 (grid%mub(i,j)+grid%mu_1(i,j))*grid%al(i,k-1,j)+ &
471 grid%mu_1(i,j)*grid%alb(i,k-1,j) )
473 grid%ph_2(i,k,j) = grid%ph_1(i,k,j)
474 grid%ph0(i,k,j) = grid%ph_1(i,k,j) + grid%phb(i,k,j)
477 if((i==2) .and. (j==2)) then
478 write(6,*) ' grid%ph_1 calc ',grid%ph_1(2,1,2),grid%ph_1(2,2,2),&
479 grid%mu_1(2,2)+grid%mub(2,2),grid%mu_1(2,2), &
480 grid%alb(2,1,2),grid%al(1,2,1),grid%rdnw(1)
486 write(6,*) ' grid%mu_1 from comp ', grid%mu_1(1,1)
487 write(6,*) ' full state sounding from comp, ph, grid%p, grid%al, grid%t_1, qv '
490 write(6,'(i3,1x,5(1x,1pe10.3))') k, grid%ph_1(1,k,1)+grid%phb(1,k,1), &
491 grid%p(1,k,1)+grid%pb(1,k,1), grid%alt(1,k,1), &
492 grid%t_1(1,k,1)+t0, moist(1,k,1,P_QV)
495 write(6,*) ' pert state sounding from comp, grid%ph_1, pp, alp, grid%t_1, qv '
497 write(6,'(i3,1x,5(1x,1pe10.3))') k, grid%ph_1(1,k,1), &
498 grid%p(1,k,1), grid%al(1,k,1), &
499 grid%t_1(1,k,1), moist(1,k,1,P_QV)
505 DO I = its, min(ide-1,ite)
508 z_at_v = grid%phb(i,1,j)/g
509 ELSE IF (j == jde) THEN
510 z_at_v = grid%phb(i,1,j-1)/g
512 z_at_v = 0.5*(grid%phb(i,1,j)+grid%phb(i,1,j-1))/g
515 p_surf = interp_0( p_in, zk, z_at_v, nl_in )
518 p_level = grid%znu(k)*(p_surf - grid%p_top) + grid%p_top
519 grid%v_1(i,k,j) = interp_0( v, p_in, p_level, nl_in )
520 grid%v_2(i,k,j) = grid%v_1(i,k,j)
528 DO J = jts, min(jde-1,jte)
532 z_at_u = grid%phb(i,1,j)/g
533 ELSE IF (i == ide) THEN
534 z_at_u = grid%phb(i-1,1,j)/g
536 z_at_u = 0.5*(grid%phb(i,1,j)+grid%phb(i-1,1,j))/g
539 p_surf = interp_0( p_in, zk, z_at_u, nl_in )
542 p_level = grid%znu(k)*(p_surf - grid%p_top) + grid%p_top
543 grid%u_1(i,k,j) = interp_0( u, p_in, p_level, nl_in )
544 grid%u_2(i,k,j) = grid%u_1(i,k,j)
552 DO J = jts, min(jde-1,jte)
554 DO I = its, min(ide-1,ite)
561 ! set a few more things
563 DO J = jts, min(jde-1,jte)
565 DO I = its, min(ide-1,ite)
566 grid%h_diabatic(i,k,j) = 0.
571 ! Go ahead and initialize these from the sounding. This will allow a run
572 ! to actually succeed even if scm_force = 0
574 grid%t_base(k) = grid%t_1(1,k,1)
575 grid%qv_base(k) = moist(1,k,1,P_QV)
576 grid%u_base(k) = grid%u_1(1,k,1)
577 grid%v_base(k) = grid%v_1(1,k,1)
578 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
583 END SUBROUTINE init_domain_rk
585 SUBROUTINE init_module_initialize
586 END SUBROUTINE init_module_initialize
588 !---------------------------------------------------------------------
590 ! test driver for get_sounding
594 ! parameter(n = 1000)
595 ! real zk(n),p(n),theta(n),rho(n),u(n),v(n),qv(n),pd(n)
601 ! call get_sounding( zk, p, pd, theta, rho, u, v, qv, dry, n, nl )
602 ! write(6,*) ' input levels ',nl
603 ! write(6,*) ' sounding '
604 ! write(6,*) ' k height(m) press (Pa) pd(Pa) theta (K) den(kg/m^3) u(m/s) v(m/s) qv(g/g) '
606 ! 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)
610 !---------------------------------------------------------------------------
612 subroutine get_sounding( zsfc, zk, p, p_dry, theta, rho, &
613 u, v, qv, dry, nl_max, nl_in )
616 integer nl_max, nl_in
618 real zk(nl_max), p(nl_max), theta(nl_max), rho(nl_max), &
619 u(nl_max), v(nl_max), qv(nl_max), p_dry(nl_max)
625 parameter( debug = .true.)
627 ! input sounding data
629 real p_surf, th_surf, qv_surf
631 real h_input(n), th_input(n), qv_input(n), u_input(n), v_input(n)
635 real rho_surf, p_input(n), rho_input(n)
636 real pm_input(n) ! this are for full moist sounding
640 real p1000mb,cv,cp,r,cvpm,g
641 parameter (p1000mb = 1.e+05, r = 287, cp = 1003., cv = cp-r, cvpm = -cv/cp, g=9.81 )
645 ! first, read the sounding
647 call read_sounding( zsfc, p_surf, th_surf, qv_surf, &
648 h_input, th_input, qv_input, u_input, v_input,n, nl, debug )
656 if(debug) write(6,*) ' number of input levels = ',nl
659 if(nl_in .gt. nl_max ) then
660 write(6,*) ' too many levels for input arrays ',nl_in,nl_max
661 call wrf_error_fatal ( ' too many levels for input arrays ' )
664 ! compute diagnostics,
665 ! first, convert qv(g/kg) to qv(g/g)
668 qv_input(k) = 0.001*qv_input(k)
671 p_surf = 100.*p_surf ! convert to pascals
672 qvf = 1. + rvovrd*qv_input(1)
673 rho_surf = 1./((r/p1000mb)*th_surf*qvf*((p_surf/p1000mb)**cvpm))
674 pi_surf = (p_surf/p1000mb)**(r/cp)
677 write(6,*) ' surface density is ',rho_surf
678 write(6,*) ' surface pi is ',pi_surf
682 ! integrate moist sounding hydrostatically, starting from the
683 ! specified surface pressure
684 ! -> first, integrate from surface to lowest level
686 qvf = 1. + rvovrd*qv_input(1)
687 qvf1 = 1. + qv_input(1)
688 rho_input(1) = rho_surf
693 write(6,*) "Your first input sounding level is below the WRF terrain elevation, aborting"
694 stop "module_initialize_scm_xy:get_sounding"
697 pm_input(1) = p_surf &
698 - 0.5*dz*(rho_surf+rho_input(1))*g*qvf1
699 rho_input(1) = 1./((r/p1000mb)*th_input(1)*qvf*((pm_input(1)/p1000mb)**cvpm))
702 ! integrate up the column
705 rho_input(k) = rho_input(k-1)
706 dz = h_input(k)-h_input(k-1)
707 qvf1 = 0.5*(2.+(qv_input(k-1)+qv_input(k)))
708 qvf = 1. + rvovrd*qv_input(k) ! qv is in g/kg here
711 pm_input(k) = pm_input(k-1) &
712 - 0.5*dz*(rho_input(k)+rho_input(k-1))*g*qvf1
713 rho_input(k) = 1./((r/p1000mb)*th_input(k)*qvf*((pm_input(k)/p1000mb)**cvpm))
716 ! we have the moist sounding
718 ! next, compute the dry sounding using p at the highest level from the
719 ! moist sounding and integrating down.
721 p_input(nl) = pm_input(nl)
724 dz = h_input(k+1)-h_input(k)
725 p_input(k) = p_input(k+1) + 0.5*dz*(rho_input(k)+rho_input(k+1))*g
732 p_dry(k) = p_input(k)
733 theta(k) = th_input(k)
734 rho(k) = rho_input(k)
742 write(6,*) ' sounding '
743 write(6,*) ' k height(m) press (Pa) pd(Pa) theta (K) den(kg/m^3) u(m/s) v(m/s) qv(g/g) '
745 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)
750 end subroutine get_sounding
752 !-------------------------------------------------------
754 subroutine read_sounding( zsfc,ps,ts,qvs,h,th,qv,u,v,n,nl,debug )
757 real zsfc,ps,ts,qvs,h(n),th(n),qv(n),u(n),v(n)
764 open(unit=10,file='input_sounding',form='formatted',status='old')
767 read(10,*) zsfc, u10, v10, t2, q2, ps
773 write(6,*) ' input sounding surface parameters '
774 write(6,*) ' surface pressure (mb) ',ps
775 write(6,*) ' surface pot. temp (K) ',ts
776 write(6,*) ' surface mixing ratio (g/kg) ',qvs
779 end_of_file = .false.
782 do while (.not. end_of_file)
784 read(10,*,end=100) h(k+1), u(k+1), v(k+1), th(k+1), qv(k+1)
786 qv(k+1) = qv(k+1)*1000.0
788 if(debug) write(6,'(1x,i3,5(1x,e10.3))') k, h(k), th(k), qv(k), u(k), v(k)
790 100 end_of_file = .true.
796 close(unit=10,status = 'keep')
798 end subroutine read_sounding
800 !-------------------------------------------------------
802 subroutine read_soil( n,nl,tmn,tsk,zs,tslb,smois )
806 real zs(n),tslb(n),smois(n)
814 open(unit=11,file='input_soil',form='formatted',status='old')
817 read(11,*) zs(1),tmn,tsk
820 write(6,*) ' input deep soil temperature (K) ',tmn
821 write(6,*) ' input skin temperature (K) ',tsk
824 end_of_file = .false.
827 do while (.not. end_of_file)
829 read(11,*,end=100) zs(k+1), tslb(k+1), smois(k+1)
831 if(debug) write(6,'(1x,i3,3(1x,f16.7))') k, zs(k), tslb(k), smois(k)
833 100 end_of_file = .true.
839 close(unit=11,status = 'keep')
841 end subroutine read_soil
843 END MODULE module_initialize_ideal