merge standard release WRF/WPS V3.0.1.1 into wrffire
[wrffire.git] / wrfv2_fire / dyn_em / module_initialize_squall2d_y.F
blob7d38105d32c9c7b7027be69fae16af52c4078c01
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
18    USE module_domain
19    USE module_io_domain
20    USE module_state_description
21    USE module_model_constants
22    USE module_bc
23    USE module_timing
24    USE module_configure
25    USE module_init_utilities
26 #ifdef DM_PARALLEL
27    USE module_dm
28 #endif
31 CONTAINS
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 )
45    IMPLICIT NONE
47    !  Input data.
48    TYPE (domain), POINTER :: grid 
49    !  Local data.
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>
58                         )
59    END SUBROUTINE init_domain
61 !-------------------------------------------------------------------
63    SUBROUTINE init_domain_rk ( grid &
65 # include <dummy_new_args.inc>
68    IMPLICIT NONE
70    !  Input data.
71    TYPE (domain), POINTER :: grid
73 # include <dummy_new_decl.inc>
75    TYPE (grid_config_rec_type)              :: config_flags
77    !  Local data
78    INTEGER                             ::                       &
79                                   ids, ide, jds, jde, kds, kde, &
80                                   ims, ime, jms, jme, kms, kme, &
81                                   its, ite, jts, jte, kts, kte, &
82                                   i, j, k
84    ! Local data
86    INTEGER, PARAMETER :: nl_max = 1000
87    REAL, DIMENSION(nl_max) :: zk, p_in, theta, rho, u, v, qv, pd_in
88    INTEGER :: nl_in
91    INTEGER :: icm,jcm, ii, im1, jj, jm1, loop, error, fid, nxc, nyc
92    REAL    :: u_mean,v_mean, f0, p_surf, p_level, qvf, z_at_v, z_at_u
93    REAL    :: z_scale, xrad, yrad, zrad, rad, delt, cof1, cof2
94 !   REAL, EXTERNAL :: interp_0
95    REAL    :: hm
96    REAL    :: pi
98 !  stuff from original initialization that has been dropped from the Registry 
99    REAL    :: vnu, xnu, xnus, dinit0, cbh, p0_temp, t0_temp, zd, zt
100    REAL    :: qvf1, qvf2, pd_surf
101    INTEGER :: it
102    real :: thtmp, ptmp, temp(3)
104    LOGICAL :: moisture_init
105    LOGICAL :: stretch_grid, dry_sounding
107    SELECT CASE ( model_data_order )
108          CASE ( DATA_ORDER_ZXY )
109    kds = grid%sd31 ; kde = grid%ed31 ;
110    ids = grid%sd32 ; ide = grid%ed32 ;
111    jds = grid%sd33 ; jde = grid%ed33 ;
113    kms = grid%sm31 ; kme = grid%em31 ;
114    ims = grid%sm32 ; ime = grid%em32 ;
115    jms = grid%sm33 ; jme = grid%em33 ;
117    kts = grid%sp31 ; kte = grid%ep31 ;   ! note that tile is entire patch
118    its = grid%sp32 ; ite = grid%ep32 ;   ! note that tile is entire patch
119    jts = grid%sp33 ; jte = grid%ep33 ;   ! note that tile is entire patch
120          CASE ( DATA_ORDER_XYZ )
121    ids = grid%sd31 ; ide = grid%ed31 ;
122    jds = grid%sd32 ; jde = grid%ed32 ;
123    kds = grid%sd33 ; kde = grid%ed33 ;
125    ims = grid%sm31 ; ime = grid%em31 ;
126    jms = grid%sm32 ; jme = grid%em32 ;
127    kms = grid%sm33 ; kme = grid%em33 ;
129    its = grid%sp31 ; ite = grid%ep31 ;   ! note that tile is entire patch
130    jts = grid%sp32 ; jte = grid%ep32 ;   ! note that tile is entire patch
131    kts = grid%sp33 ; kte = grid%ep33 ;   ! note that tile is entire patch
132          CASE ( DATA_ORDER_XZY )
133    ids = grid%sd31 ; ide = grid%ed31 ;
134    kds = grid%sd32 ; kde = grid%ed32 ;
135    jds = grid%sd33 ; jde = grid%ed33 ;
137    ims = grid%sm31 ; ime = grid%em31 ;
138    kms = grid%sm32 ; kme = grid%em32 ;
139    jms = grid%sm33 ; jme = grid%em33 ;
141    its = grid%sp31 ; ite = grid%ep31 ;   ! note that tile is entire patch
142    kts = grid%sp32 ; kte = grid%ep32 ;   ! note that tile is entire patch
143    jts = grid%sp33 ; jte = grid%ep33 ;   ! note that tile is entire patch
145    END SELECT
148    stretch_grid = .true.
149    delt = 3.
150 !   z_scale = .50
151    z_scale = .40
152    pi = 2.*asin(1.0)
153    write(6,*) ' pi is ',pi
154    nxc = (ide-ids)/2
155    nyc = (jde-jds)/2
157    CALL model_to_grid_config_rec ( grid%id , model_config_rec , config_flags )
159 ! here we check to see if the boundary conditions are set properly
161    CALL boundary_condition_check( config_flags, bdyzone, error, grid%id )
163    moisture_init = .true.
165     grid%itimestep=0
167 #ifdef DM_PARALLEL
168    CALL wrf_dm_bcast_bytes( icm , IWORDSIZE )
169    CALL wrf_dm_bcast_bytes( jcm , IWORDSIZE )
170 #endif
172     CALL nl_set_mminlu(1, '    ')
173     CALL nl_set_iswater(1,0)
174     CALL nl_set_cen_lat(1,40.)
175     CALL nl_set_cen_lon(1,-105.)
176     CALL nl_set_truelat1(1,0.)
177     CALL nl_set_truelat2(1,0.)
178     CALL nl_set_moad_cen_lat (1,0.)
179     CALL nl_set_stand_lon (1,0.)
180     CALL nl_set_map_proj(1,0)
183 !  here we initialize data we currently is not initialized 
184 !  in the input data
186     DO j = jts, jte
187       DO i = its, ite
188          grid%msftx(i,j)    = 1.
189          grid%msfty(i,j)    = 1.
190          grid%msfux(i,j)    = 1.
191          grid%msfuy(i,j)    = 1.
192          grid%msfvx(i,j)    = 1.
193          grid%msfvx_inv(i,j)= 1.
194          grid%msfvy(i,j)    = 1.
195          grid%sina(i,j)     = 0.
196          grid%cosa(i,j)     = 1.
197          grid%e(i,j)        = 0.
198          grid%f(i,j)        = 0.
200       END DO
201    END DO
203     DO j = jts, jte
204     DO k = kts, kte
205       DO i = its, ite
206          grid%ww(i,k,j)     = 0.
207       END DO
208    END DO
209    END DO
211    grid%step_number = 0
213 ! set up the grid
215    IF (stretch_grid) THEN ! exponential stretch for eta (nearly constant dz)
216      DO k=1, kde
217       grid%znw(k) = (exp(-(k-1)/float(kde-1)/z_scale) - exp(-1./z_scale))/ &
218                                 (1.-exp(-1./z_scale))
219      ENDDO
220    ELSE
221      DO k=1, kde
222       grid%znw(k) = 1. - float(k-1)/float(kde-1)
223      ENDDO
224    ENDIF
226    DO k=1, kde-1
227     grid%dnw(k) = grid%znw(k+1) - grid%znw(k)
228     grid%rdnw(k) = 1./grid%dnw(k)
229     grid%znu(k) = 0.5*(grid%znw(k+1)+grid%znw(k))
230    ENDDO
231    DO k=2, kde-1
232     grid%dn(k) = 0.5*(grid%dnw(k)+grid%dnw(k-1))
233     grid%rdn(k) = 1./grid%dn(k)
234     grid%fnp(k) = .5* grid%dnw(k  )/grid%dn(k)
235     grid%fnm(k) = .5* grid%dnw(k-1)/grid%dn(k)
236    ENDDO
238    cof1 = (2.*grid%dn(2)+grid%dn(3))/(grid%dn(2)+grid%dn(3))*grid%dnw(1)/grid%dn(2) 
239    cof2 =     grid%dn(2)        /(grid%dn(2)+grid%dn(3))*grid%dnw(1)/grid%dn(3) 
240    grid%cf1  = grid%fnp(2) + cof1
241    grid%cf2  = grid%fnm(2) - cof1 - cof2
242    grid%cf3  = cof2       
244    grid%cfn  = (.5*grid%dnw(kde-1)+grid%dn(kde-1))/grid%dn(kde-1)
245    grid%cfn1 = -.5*grid%dnw(kde-1)/grid%dn(kde-1)
246    grid%rdx = 1./config_flags%dx
247    grid%rdy = 1./config_flags%dy
249 !  get the sounding from the ascii sounding file, first get dry sounding and 
250 !  calculate base state
252   write(6,*) ' getting dry sounding for base state '
253   dry_sounding = .true.
254   CALL get_sounding( zk, p_in, pd_in, theta, rho, u, v, qv, dry_sounding, nl_max, nl_in )
256   write(6,*) ' returned from reading sounding, nl_in is ',nl_in
258 !  find ptop for the desired ztop (ztop is input from the namelist),
259 !  and find surface pressure
261   grid%p_top = interp_0( p_in, zk, config_flags%ztop, nl_in )
263   DO j=jts,jte
264   DO i=its,ite  ! flat surface
265     grid%phb(i,1,j) = 0.
266     grid%php(i,1,j) = 0.
267     grid%ph0(i,1,j) = 0.
268     grid%ht(i,j) = 0.
269   ENDDO
270   ENDDO
272   DO J = jts, jte
273   DO I = its, ite
275     p_surf = interp_0( p_in, zk, grid%phb(i,1,j)/g, nl_in )
276     grid%mub(i,j) = p_surf-grid%p_top
278 !  this is dry hydrostatic sounding (base state), so given grid%p (coordinate),
279 !  interp theta (from interp) and compute 1/rho from eqn. of state
281     DO K = 1, kte-1
282       p_level = grid%znu(k)*(p_surf - grid%p_top) + grid%p_top
283       grid%pb(i,k,j) = p_level
284       grid%t_init(i,k,j) = interp_0( theta, p_in, p_level, nl_in ) - t0
285       grid%alb(i,k,j) = (r_d/p1000mb)*(grid%t_init(i,k,j)+t0)*(grid%pb(i,k,j)/p1000mb)**cvpm
286     ENDDO
288 !  calc hydrostatic balance (alternatively we could interp the geopotential from the
289 !  sounding, but this assures that the base state is in exact hydrostatic balance with
290 !  respect to the model eqns.
292     DO k  = 2,kte
293       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)
294     ENDDO
296   ENDDO
297   ENDDO
299   write(6,*) ' ptop is ',grid%p_top
300   write(6,*) ' base state grid%mub(1,1), p_surf is ',grid%mub(1,1),grid%mub(1,1)+grid%p_top
302 !  calculate full state for each column - this includes moisture.
304   write(6,*) ' getting moist sounding for full state '
305   dry_sounding = .false.
306   CALL get_sounding( zk, p_in, pd_in, theta, rho, u, v, qv, dry_sounding, nl_max, nl_in )
308   DO J = jts, min(jde-1,jte)
309   DO I = its, min(ide-1,ite)
311 !  At this point grid%p_top is already set. find the DRY mass in the column 
312 !  by interpolating the DRY pressure.  
314    pd_surf = interp_0( pd_in, zk, grid%phb(i,1,j)/g, nl_in )
316 !  compute the perturbation mass and the full mass
318     grid%mu_1(i,j) = pd_surf-grid%p_top - grid%mub(i,j)
319     grid%mu_2(i,j) = grid%mu_1(i,j)
320     grid%mu0(i,j) = grid%mu_1(i,j) + grid%mub(i,j)
322 ! given the dry pressure and coordinate system, interp the potential
323 ! temperature and qv
325     do k=1,kde-1
327       p_level = grid%znu(k)*(pd_surf - grid%p_top) + grid%p_top
329       moist(i,k,j,P_QV) = interp_0( qv, pd_in, p_level, nl_in )
330       grid%t_1(i,k,j)          = interp_0( theta, pd_in, p_level, nl_in ) - t0
331       grid%t_2(i,k,j)          = grid%t_1(i,k,j)
332       
334     enddo
336 !  integrate the hydrostatic equation (from the RHS of the bigstep
337 !  vertical momentum equation) down from the top to get grid%p.
338 !  first from the top of the model to the top pressure
340     k = kte-1  ! top level
342     qvf1 = 0.5*(moist(i,k,j,P_QV)+moist(i,k,j,P_QV))
343     qvf2 = 1./(1.+qvf1)
344     qvf1 = qvf1*qvf2
346 !    grid%p(i,k,j) = - 0.5*grid%mu_1(i,j)/grid%rdnw(k)
347     grid%p(i,k,j) = - 0.5*(grid%mu_1(i,j)+qvf1*grid%mub(i,j))/grid%rdnw(k)/qvf2
348     qvf = 1. + rvovrd*moist(i,k,j,P_QV)
349     grid%alt(i,k,j) = (r_d/p1000mb)*(grid%t_1(i,k,j)+t0)*qvf* &
350                 (((grid%p(i,k,j)+grid%pb(i,k,j))/p1000mb)**cvpm)
351     grid%al(i,k,j) = grid%alt(i,k,j) - grid%alb(i,k,j)
353 !  down the column
355     do k=kte-2,1,-1
356       qvf1 = 0.5*(moist(i,k,j,P_QV)+moist(i,k+1,j,P_QV))
357       qvf2 = 1./(1.+qvf1)
358       qvf1 = qvf1*qvf2
359       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)
360       qvf = 1. + rvovrd*moist(i,k,j,P_QV)
361       grid%alt(i,k,j) = (r_d/p1000mb)*(grid%t_1(i,k,j)+t0)*qvf* &
362                   (((grid%p(i,k,j)+grid%pb(i,k,j))/p1000mb)**cvpm)
363       grid%al(i,k,j) = grid%alt(i,k,j) - grid%alb(i,k,j)
364     enddo
366 !  this is the hydrostatic equation used in the model after the
367 !  small timesteps.  In the model, grid%al (inverse density)
368 !  is computed from the geopotential.
371     grid%ph_1(i,1,j) = 0.
372     DO k  = 2,kte
373       grid%ph_1(i,k,j) = grid%ph_1(i,k-1,j) - (1./grid%rdnw(k-1))*(       &
374                    (grid%mub(i,j)+grid%mu_1(i,j))*grid%al(i,k-1,j)+ &
375                     grid%mu_1(i,j)*grid%alb(i,k-1,j)  )
376                                                    
377       grid%ph_2(i,k,j) = grid%ph_1(i,k,j) 
378       grid%ph0(i,k,j) = grid%ph_1(i,k,j) + grid%phb(i,k,j)
379     ENDDO
381     if((i==2) .and. (j==2)) then
382      write(6,*) ' grid%ph_1 calc ',grid%ph_1(2,1,2),grid%ph_1(2,2,2),&
383                               grid%mu_1(2,2)+grid%mub(2,2),grid%mu_1(2,2), &
384                               grid%alb(2,1,2),grid%al(1,2,1),grid%rdnw(1)
385     endif
387   ENDDO
388   ENDDO
390 !#if 0
392 !  thermal perturbation to kick off convection
394   write(6,*) ' nxc, nyc for perturbation ',nxc,nyc
395   write(6,*) ' delt for perturbation ',delt
397   DO J = jts, min(jde-1,jte)
398     yrad = config_flags%dy*float(j-nyc)/4000.
399 !     yrad = 0.
400     DO I = its, min(ide-1,ite)
401 !      xrad = config_flags%dx*float(i-nxc)/4000.
402      xrad = 0.
403       DO K = 1, kte-1
405 !  put in preturbation theta (bubble) and recalc density.  note,
406 !  the mass in the column is not changing, so when theta changes,
407 !  we recompute density and geopotential
409         zrad = 0.5*(grid%ph_1(i,k,j)+grid%ph_1(i,k+1,j)  &
410                    +grid%phb(i,k,j)+grid%phb(i,k+1,j))/g
411         zrad = (zrad-1500.)/1500.
412         RAD=SQRT(xrad*xrad+yrad*yrad+zrad*zrad)
413         IF(RAD <= 1.) THEN
414            grid%t_1(i,k,j)=grid%t_1(i,k,j)+delt*COS(.5*PI*RAD)**2
415            grid%t_2(i,k,j)=grid%t_1(i,k,j)
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)
420         ENDIF
421       ENDDO
423 !  rebalance hydrostatically
425       DO k  = 2,kte
426         grid%ph_1(i,k,j) = grid%ph_1(i,k-1,j) - (1./grid%rdnw(k-1))*(       &
427                      (grid%mub(i,j)+grid%mu_1(i,j))*grid%al(i,k-1,j)+ &
428                       grid%mu_1(i,j)*grid%alb(i,k-1,j)  )
429                                                    
430         grid%ph_2(i,k,j) = grid%ph_1(i,k,j) 
431         grid%ph0(i,k,j) = grid%ph_1(i,k,j) + grid%phb(i,k,j)
432       ENDDO
434     ENDDO
435   ENDDO
437 !#endif
439    write(6,*) ' grid%mu_1 from comp ', grid%mu_1(1,1)
440    write(6,*) ' full state sounding from comp, ph, grid%p, grid%al, grid%t_1, qv '
441    do k=1,kde-1
442      write(6,'(i3,1x,5(1x,1pe10.3))') k, grid%ph_1(1,k,1)+grid%phb(1,k,1), &
443                                       grid%p(1,k,1)+grid%pb(1,k,1), grid%alt(1,k,1), &
444                                       grid%t_1(1,k,1)+t0, moist(1,k,1,P_QV)
445    enddo
447    write(6,*) ' pert state sounding from comp, grid%ph_1, pp, alp, grid%t_1, qv '
448    do k=1,kde-1
449      write(6,'(i3,1x,5(1x,1pe10.3))') k, grid%ph_1(1,k,1), &
450                                       grid%p(1,k,1), grid%al(1,k,1), &
451                                       grid%t_1(1,k,1), moist(1,k,1,P_QV)
452    enddo
454 ! interp v
456   DO J = jts, jte
457   DO I = its, min(ide-1,ite)
459     IF (j == jds) THEN
460       z_at_v = grid%phb(i,1,j)/g
461     ELSE IF (j == jde) THEN
462       z_at_v = grid%phb(i,1,j-1)/g
463     ELSE
464       z_at_v = 0.5*(grid%phb(i,1,j)+grid%phb(i,1,j-1))/g
465     END IF
467     p_surf = interp_0( p_in, zk, z_at_v, nl_in )
469     DO K = 1, kte
470       p_level = grid%znu(k)*(p_surf - grid%p_top) + grid%p_top
471       grid%v_1(i,k,j) = interp_0( v, p_in, p_level, nl_in )
472       grid%v_2(i,k,j) = grid%v_1(i,k,j)
473     ENDDO
475   ENDDO
476   ENDDO
478 ! interp u
480   DO J = jts, min(jde-1,jte)
481   DO I = its, ite
483     IF (i == ids) THEN
484       z_at_u = grid%phb(i,1,j)/g
485     ELSE IF (i == ide) THEN
486       z_at_u = grid%phb(i-1,1,j)/g
487     ELSE
488       z_at_u = 0.5*(grid%phb(i,1,j)+grid%phb(i-1,1,j))/g
489     END IF
491     p_surf = interp_0( p_in, zk, z_at_u, nl_in )
493     DO K = 1, kte
494       p_level = grid%znu(k)*(p_surf - grid%p_top) + grid%p_top
495       grid%u_1(i,k,j) = interp_0( u, p_in, p_level, nl_in )
496       grid%u_2(i,k,j) = grid%u_1(i,k,j)
497     ENDDO
499   ENDDO
500   ENDDO
502 !  set w
504   DO J = jts, min(jde-1,jte)
505   DO K = kts, kte
506   DO I = its, min(ide-1,ite)
507     grid%w_1(i,k,j) = 0.
508     grid%w_2(i,k,j) = 0.
509   ENDDO
510   ENDDO
511   ENDDO
513 !  set a few more things
515   DO J = jts, min(jde-1,jte)
516   DO K = kts, kte-1
517   DO I = its, min(ide-1,ite)
518     grid%h_diabatic(i,k,j) = 0.
519   ENDDO
520   ENDDO
521   ENDDO
523   DO k=1,kte-1
524     grid%t_base(k) = grid%t_1(1,k,1)
525     grid%qv_base(k) = moist(1,k,1,P_QV)
526     grid%u_base(k) = grid%u_1(1,k,1)
527     grid%v_base(k) = grid%v_1(1,k,1)
528     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
529   ENDDO
531   DO J = jts, min(jde-1,jte)
532   DO I = its, min(ide-1,ite)
533      thtmp   = grid%t_2(i,1,j)+t0
534      ptmp    = grid%p(i,1,j)+grid%pb(i,1,j)
535      temp(1) = thtmp * (ptmp/p1000mb)**rcp
536      thtmp   = grid%t_2(i,2,j)+t0
537      ptmp    = grid%p(i,2,j)+grid%pb(i,2,j)
538      temp(2) = thtmp * (ptmp/p1000mb)**rcp
539      thtmp   = grid%t_2(i,3,j)+t0
540      ptmp    = grid%p(i,3,j)+grid%pb(i,3,j)
541      temp(3) = thtmp * (ptmp/p1000mb)**rcp
543      grid%tsk(I,J)=grid%cf1*temp(1)+grid%cf2*temp(2)+grid%cf3*temp(3)
544      grid%tmn(I,J)=grid%tsk(I,J)-0.5
545   ENDDO
546   ENDDO
548   RETURN
550  END SUBROUTINE init_domain_rk
552    SUBROUTINE init_module_initialize
553    END SUBROUTINE init_module_initialize
555 !---------------------------------------------------------------------
557 !  test driver for get_sounding
559 !      implicit none
560 !      integer n
561 !      parameter(n = 1000)
562 !      real zk(n),p(n),theta(n),rho(n),u(n),v(n),qv(n),pd(n)
563 !      logical dry
564 !      integer nl,k
566 !      dry = .false.
567 !      dry = .true.
568 !      call get_sounding( zk, p, pd, theta, rho, u, v, qv, dry, n, nl )
569 !      write(6,*) ' input levels ',nl
570 !      write(6,*) ' sounding '
571 !      write(6,*) '  k  height(m)  press (Pa) pd(Pa) theta (K) den(kg/m^3)  u(m/s)     v(m/s)    qv(g/g) '
572 !      do k=1,nl
573 !        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)
574 !      enddo
575 !      end
577 !---------------------------------------------------------------------------
579       subroutine get_sounding( zk, p, p_dry, theta, rho, &
580                                u, v, qv, dry, nl_max, nl_in )
581       implicit none
583       integer nl_max, nl_in
584       real zk(nl_max), p(nl_max), theta(nl_max), rho(nl_max), &
585            u(nl_max), v(nl_max), qv(nl_max), p_dry(nl_max)
586       logical dry
588       integer n
589       parameter(n=1000)
590       logical debug
591       parameter( debug = .true.)
593 ! input sounding data
595       real p_surf, th_surf, qv_surf
596       real pi_surf, pi(n)
597       real h_input(n), th_input(n), qv_input(n), u_input(n), v_input(n)
599 ! diagnostics
601       real rho_surf, p_input(n), rho_input(n)
602       real pm_input(n)  !  this are for full moist sounding
604 ! local data
606       real p1000mb,cv,cp,r,cvpm,g
607       parameter (p1000mb = 1.e+05, r = 287, cp = 1003., cv = cp-r, cvpm = -cv/cp, g=9.81 )
608       integer k, it, nl
609       real qvf, qvf1, dz
611 !  first, read the sounding
613       call read_sounding( p_surf, th_surf, qv_surf, &
614                           h_input, th_input, qv_input, u_input, v_input,n, nl, debug )
616       if(dry) then
617        do k=1,nl
618          qv_input(k) = 0.
619        enddo
620       endif
622       if(debug) write(6,*) ' number of input levels = ',nl
624         nl_in = nl
625         if(nl_in .gt. nl_max ) then
626           write(6,*) ' too many levels for input arrays ',nl_in,nl_max
627           call wrf_error_fatal ( ' too many levels for input arrays ' )
628         end if
630 !  compute diagnostics,
631 !  first, convert qv(g/kg) to qv(g/g)
633       do k=1,nl
634         qv_input(k) = 0.001*qv_input(k)
635       enddo
637       p_surf = 100.*p_surf  ! convert to pascals
638       qvf = 1. + rvovrd*qv_input(1) 
639       rho_surf = 1./((r/p1000mb)*th_surf*qvf*((p_surf/p1000mb)**cvpm))
640       pi_surf = (p_surf/p1000mb)**(r/cp)
642       if(debug) then
643         write(6,*) ' surface density is ',rho_surf
644         write(6,*) ' surface pi is      ',pi_surf
645       end if
648 !  integrate moist sounding hydrostatically, starting from the
649 !  specified surface pressure
650 !  -> first, integrate from surface to lowest level
652           qvf = 1. + rvovrd*qv_input(1) 
653           qvf1 = 1. + qv_input(1)
654           rho_input(1) = rho_surf
655           dz = h_input(1)
656           do it=1,10
657             pm_input(1) = p_surf &
658                     - 0.5*dz*(rho_surf+rho_input(1))*g*qvf1
659             rho_input(1) = 1./((r/p1000mb)*th_input(1)*qvf*((pm_input(1)/p1000mb)**cvpm))
660           enddo
662 ! integrate up the column
664           do k=2,nl
665             rho_input(k) = rho_input(k-1)
666             dz = h_input(k)-h_input(k-1)
667             qvf1 = 0.5*(2.+(qv_input(k-1)+qv_input(k)))
668             qvf = 1. + rvovrd*qv_input(k)   ! qv is in g/kg here
670             do it=1,10
671               pm_input(k) = pm_input(k-1) &
672                       - 0.5*dz*(rho_input(k)+rho_input(k-1))*g*qvf1
673               rho_input(k) = 1./((r/p1000mb)*th_input(k)*qvf*((pm_input(k)/p1000mb)**cvpm))
674             enddo
675           enddo
677 !  we have the moist sounding
679 !  next, compute the dry sounding using p at the highest level from the
680 !  moist sounding and integrating down.
682         p_input(nl) = pm_input(nl)
684           do k=nl-1,1,-1
685             dz = h_input(k+1)-h_input(k)
686             p_input(k) = p_input(k+1) + 0.5*dz*(rho_input(k)+rho_input(k+1))*g
687           enddo
690         do k=1,nl
692           zk(k) = h_input(k)
693           p(k) = pm_input(k)
694           p_dry(k) = p_input(k)
695           theta(k) = th_input(k)
696           rho(k) = rho_input(k)
697           u(k) = u_input(k)
698           v(k) = v_input(k)
699           qv(k) = qv_input(k)
701         enddo
703      if(debug) then
704       write(6,*) ' sounding '
705       write(6,*) '  k  height(m)  press (Pa) pd(Pa) theta (K) den(kg/m^3)  u(m/s)     v(m/s)    qv(g/g) '
706       do k=1,nl
707         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)
708       enddo
710      end if
712       end subroutine get_sounding
714 !-------------------------------------------------------
716       subroutine read_sounding( ps,ts,qvs,h,th,qv,u,v,n,nl,debug )
717       implicit none
718       integer n,nl
719       real ps,ts,qvs,h(n),th(n),qv(n),u(n),v(n)
720       logical end_of_file
721       logical debug
723       integer k
725       open(unit=10,file='input_sounding',form='formatted',status='old')
726       rewind(10)
727       read(10,*) ps, ts, qvs
728       if(debug) then
729         write(6,*) ' input sounding surface parameters '
730         write(6,*) ' surface pressure (mb) ',ps
731         write(6,*) ' surface pot. temp (K) ',ts
732         write(6,*) ' surface mixing ratio (g/kg) ',qvs
733       end if
735       end_of_file = .false.
736       k = 0
738       do while (.not. end_of_file)
740         read(10,*,end=100) h(k+1), th(k+1), qv(k+1), u(k+1), v(k+1)
741         k = k+1
742         if(debug) write(6,'(1x,i3,5(1x,e10.3))') k, h(k), th(k), qv(k), u(k), v(k)
743         go to 110
744  100    end_of_file = .true.
745  110    continue
746       enddo
748       nl = k
750       close(unit=10,status = 'keep')
752       end subroutine read_sounding
754 END MODULE module_initialize_ideal