merge standard release WRF/WPS V3.0.1.1 into wrffire
[wrffire.git] / wrfv2_fire / dyn_em / module_initialize_quarter_ss.F
blobfb9213f69bd05826b8a30b4bde3db8ebb9d6b47b
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                         )
60    END SUBROUTINE init_domain
62 !-------------------------------------------------------------------
64    SUBROUTINE init_domain_rk ( grid &
66 # include <dummy_new_args.inc>
69    IMPLICIT NONE
71    !  Input data.
72    TYPE (domain), POINTER :: grid
74 # include <dummy_new_decl.inc>
76    TYPE (grid_config_rec_type)              :: config_flags
78    !  Local data
79    INTEGER                             ::                       &
80                                   ids, ide, jds, jde, kds, kde, &
81                                   ims, ime, jms, jme, kms, kme, &
82                                   its, ite, jts, jte, kts, kte, &
83                                   i, j, k
85    ! Local data
87    INTEGER, PARAMETER :: nl_max = 1000
88    REAL, DIMENSION(nl_max) :: zk, p_in, theta, rho, u, v, qv, pd_in
89    INTEGER :: nl_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
96    REAL    :: hm
97    REAL    :: pi
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
102    INTEGER :: it
103    real :: thtmp, ptmp, temp(3)
105    LOGICAL :: moisture_init
106    LOGICAL :: stretch_grid, dry_sounding
108   INTEGER :: xs , xe , ys , ye
109   REAL :: mtn_ht
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
150    END SELECT
153    stretch_grid = .true.
154    delt = 3.
155 !   z_scale = .50
156    z_scale = .40
157    pi = 2.*asin(1.0)
158    write(6,*) ' pi is ',pi
159    nxc = (ide-ids)/2
160    nyc = (jde-jds)/2
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.
170     grid%itimestep=0
172 #ifdef DM_PARALLEL
173    CALL wrf_dm_bcast_bytes( icm , IWORDSIZE )
174    CALL wrf_dm_bcast_bytes( jcm , IWORDSIZE )
175 #endif
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_map_proj(1,0)
188 !  here we initialize data we currently is not initialized 
189 !  in the input data
191     DO j = jts, jte
192       DO i = its, ite
193          grid%msftx(i,j)    = 1.
194          grid%msfty(i,j)    = 1.
195          grid%msfux(i,j)    = 1.
196          grid%msfuy(i,j)    = 1.
197          grid%msfvx(i,j)    = 1.
198          grid%msfvx_inv(i,j)= 1.
199          grid%msfvy(i,j)    = 1.
200          grid%sina(i,j)     = 0.
201          grid%cosa(i,j)     = 1.
202          grid%e(i,j)        = 0.
203          grid%f(i,j)        = 0.
205       END DO
206    END DO
208     DO j = jts, jte
209     DO k = kts, kte
210       DO i = its, ite
211          grid%ww(i,k,j)     = 0.
212       END DO
213    END DO
214    END DO
216    grid%step_number = 0
218 ! set up the grid
220    IF (stretch_grid) THEN ! exponential stretch for eta (nearly constant dz)
221      DO k=1, kde
222       grid%znw(k) = (exp(-(k-1)/float(kde-1)/z_scale) - exp(-1./z_scale))/ &
223                                 (1.-exp(-1./z_scale))
224      ENDDO
225    ELSE
226      DO k=1, kde
227       grid%znw(k) = 1. - float(k-1)/float(kde-1)
228      ENDDO
229    ENDIF
231    DO k=1, kde-1
232     grid%dnw(k) = grid%znw(k+1) - grid%znw(k)
233     grid%rdnw(k) = 1./grid%dnw(k)
234     grid%znu(k) = 0.5*(grid%znw(k+1)+grid%znw(k))
235    ENDDO
236    DO k=2, kde-1
237     grid%dn(k) = 0.5*(grid%dnw(k)+grid%dnw(k-1))
238     grid%rdn(k) = 1./grid%dn(k)
239     grid%fnp(k) = .5* grid%dnw(k  )/grid%dn(k)
240     grid%fnm(k) = .5* grid%dnw(k-1)/grid%dn(k)
241    ENDDO
243    cof1 = (2.*grid%dn(2)+grid%dn(3))/(grid%dn(2)+grid%dn(3))*grid%dnw(1)/grid%dn(2) 
244    cof2 =     grid%dn(2)        /(grid%dn(2)+grid%dn(3))*grid%dnw(1)/grid%dn(3) 
245    grid%cf1  = grid%fnp(2) + cof1
246    grid%cf2  = grid%fnm(2) - cof1 - cof2
247    grid%cf3  = cof2       
249    grid%cfn  = (.5*grid%dnw(kde-1)+grid%dn(kde-1))/grid%dn(kde-1)
250    grid%cfn1 = -.5*grid%dnw(kde-1)/grid%dn(kde-1)
251    grid%rdx = 1./config_flags%dx
252    grid%rdy = 1./config_flags%dy
254 !  get the sounding from the ascii sounding file, first get dry sounding and 
255 !  calculate base state
257   dry_sounding = .true.
258   IF ( wrf_dm_on_monitor() ) THEN
259   write(6,*) ' getting dry sounding for base state '
261   CALL get_sounding( zk, p_in, pd_in, theta, rho, u, v, qv, dry_sounding, nl_max, nl_in )
262   ENDIF
263   CALL wrf_dm_bcast_real( zk , nl_max )
264   CALL wrf_dm_bcast_real( p_in , nl_max )
265   CALL wrf_dm_bcast_real( pd_in , nl_max )
266   CALL wrf_dm_bcast_real( theta , nl_max )
267   CALL wrf_dm_bcast_real( rho , nl_max )
268   CALL wrf_dm_bcast_real( u , nl_max )
269   CALL wrf_dm_bcast_real( v , nl_max )
270   CALL wrf_dm_bcast_real( qv , nl_max )
271   CALL wrf_dm_bcast_integer ( nl_in , 1 ) 
273   write(6,*) ' returned from reading sounding, nl_in is ',nl_in
275 !  find ptop for the desired ztop (ztop is input from the namelist),
276 !  and find surface pressure
278   grid%p_top = interp_0( p_in, zk, config_flags%ztop, nl_in )
280   DO j=jts,jte
281   DO i=its,ite
282     grid%ht(i,j) = 0.
283   ENDDO
284   ENDDO
286   xs=ide/2 -3
287   xs=ids   -3
288   xe=xs + 6
289   ys=jde/2 -3
290   ye=ys + 6
291   mtn_ht = 500
292 #ifdef MTN
293   DO j=max(ys,jds),min(ye,jde-1)
294   DO i=max(xs,ids),min(xe,ide-1)
295      grid%ht(i,j) = mtn_ht * 0.25 * &
296                ( 1. + COS ( 2*pi/(xe-xs) * ( i-xs ) + pi ) ) * &
297                ( 1. + COS ( 2*pi/(ye-ys) * ( j-ys ) + pi ) )
298   ENDDO
299   ENDDO
300 #endif
301 #ifdef EW_RIDGE
302   DO j=max(ys,jds),min(ye,jde-1)
303   DO i=ids,ide
304      grid%ht(i,j) = mtn_ht * 0.50 * &
305                ( 1. + COS ( 2*pi/(ye-ys) * ( j-ys ) + pi ) )
306   ENDDO
307   ENDDO
308 #endif
309 #ifdef NS_RIDGE
310   DO j=jds,jde
311   DO i=max(xs,ids),min(xe,ide-1)
312      grid%ht(i,j) = mtn_ht * 0.50 * &
313                ( 1. + COS ( 2*pi/(xe-xs) * ( i-xs ) + pi ) )
314   ENDDO
315   ENDDO
316 #endif
317   DO j=jts,jte
318   DO i=its,ite
319     grid%phb(i,1,j) = g * grid%ht(i,j)
320     grid%ph0(i,1,j) = g * grid%ht(i,j)
321   ENDDO
322   ENDDO
324   DO J = jts, jte
325   DO I = its, ite
327     p_surf = interp_0( p_in, zk, grid%phb(i,1,j)/g, nl_in )
328     grid%mub(i,j) = p_surf-grid%p_top
330 !  this is dry hydrostatic sounding (base state), so given grid%p (coordinate),
331 !  interp theta (from interp) and compute 1/rho from eqn. of state
333     DO K = 1, kte-1
334       p_level = grid%znu(k)*(p_surf - grid%p_top) + grid%p_top
335       grid%pb(i,k,j) = p_level
336       grid%t_init(i,k,j) = interp_0( theta, p_in, p_level, nl_in ) - t0
337       grid%alb(i,k,j) = (r_d/p1000mb)*(grid%t_init(i,k,j)+t0)*(grid%pb(i,k,j)/p1000mb)**cvpm
338     ENDDO
340 !  calc hydrostatic balance (alternatively we could interp the geopotential from the
341 !  sounding, but this assures that the base state is in exact hydrostatic balance with
342 !  respect to the model eqns.
344     DO k  = 2,kte
345       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)
346     ENDDO
348   ENDDO
349   ENDDO
351   IF ( wrf_dm_on_monitor() ) THEN
352     write(6,*) ' ptop is ',grid%p_top
353     write(6,*) ' base state grid%mub(1,1), p_surf is ',grid%mub(1,1),grid%mub(1,1)+grid%p_top
354   ENDIF
356 !  calculate full state for each column - this includes moisture.
358   write(6,*) ' getting moist sounding for full state '
359   dry_sounding = .false.
360   CALL get_sounding( zk, p_in, pd_in, theta, rho, u, v, qv, dry_sounding, nl_max, nl_in )
362   DO J = jts, min(jde-1,jte)
363   DO I = its, min(ide-1,ite)
365 !  At this point grid%p_top is already set. find the DRY mass in the column 
366 !  by interpolating the DRY pressure.  
368    pd_surf = interp_0( pd_in, zk, grid%phb(i,1,j)/g, nl_in )
370 !  compute the perturbation mass and the full mass
372     grid%mu_1(i,j) = pd_surf-grid%p_top - grid%mub(i,j)
373     grid%mu_2(i,j) = grid%mu_1(i,j)
374     grid%mu0(i,j) = grid%mu_1(i,j) + grid%mub(i,j)
376 ! given the dry pressure and coordinate system, interp the potential
377 ! temperature and qv
379     do k=1,kde-1
381       p_level = grid%znu(k)*(pd_surf - grid%p_top) + grid%p_top
383       moist(i,k,j,P_QV) = interp_0( qv, pd_in, p_level, nl_in )
384       grid%t_1(i,k,j)          = interp_0( theta, pd_in, p_level, nl_in ) - t0
385       grid%t_2(i,k,j)          = grid%t_1(i,k,j)
386       
388     enddo
390 !  integrate the hydrostatic equation (from the RHS of the bigstep
391 !  vertical momentum equation) down from the top to get grid%p.
392 !  first from the top of the model to the top pressure
394     k = kte-1  ! top level
396     qvf1 = 0.5*(moist(i,k,j,P_QV)+moist(i,k,j,P_QV))
397     qvf2 = 1./(1.+qvf1)
398     qvf1 = qvf1*qvf2
400 !    grid%p(i,k,j) = - 0.5*grid%mu_1(i,j)/grid%rdnw(k)
401     grid%p(i,k,j) = - 0.5*(grid%mu_1(i,j)+qvf1*grid%mub(i,j))/grid%rdnw(k)/qvf2
402     qvf = 1. + rvovrd*moist(i,k,j,P_QV)
403     grid%alt(i,k,j) = (r_d/p1000mb)*(grid%t_1(i,k,j)+t0)*qvf* &
404                 (((grid%p(i,k,j)+grid%pb(i,k,j))/p1000mb)**cvpm)
405     grid%al(i,k,j) = grid%alt(i,k,j) - grid%alb(i,k,j)
407 !  down the column
409     do k=kte-2,1,-1
410       qvf1 = 0.5*(moist(i,k,j,P_QV)+moist(i,k+1,j,P_QV))
411       qvf2 = 1./(1.+qvf1)
412       qvf1 = qvf1*qvf2
413       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)
414       qvf = 1. + rvovrd*moist(i,k,j,P_QV)
415       grid%alt(i,k,j) = (r_d/p1000mb)*(grid%t_1(i,k,j)+t0)*qvf* &
416                   (((grid%p(i,k,j)+grid%pb(i,k,j))/p1000mb)**cvpm)
417       grid%al(i,k,j) = grid%alt(i,k,j) - grid%alb(i,k,j)
418     enddo
420 !  this is the hydrostatic equation used in the model after the
421 !  small timesteps.  In the model, grid%al (inverse density)
422 !  is computed from the geopotential.
425     grid%ph_1(i,1,j) = 0.
426     DO k  = 2,kte
427       grid%ph_1(i,k,j) = grid%ph_1(i,k-1,j) - (1./grid%rdnw(k-1))*(       &
428                    (grid%mub(i,j)+grid%mu_1(i,j))*grid%al(i,k-1,j)+ &
429                     grid%mu_1(i,j)*grid%alb(i,k-1,j)  )
430                                                    
431       grid%ph_2(i,k,j) = grid%ph_1(i,k,j) 
432       grid%ph0(i,k,j) = grid%ph_1(i,k,j) + grid%phb(i,k,j)
433     ENDDO
435     IF ( wrf_dm_on_monitor() ) THEN
436     if((i==2) .and. (j==2)) then
437      write(6,*) ' grid%ph_1 calc ',grid%ph_1(2,1,2),grid%ph_1(2,2,2),&
438                               grid%mu_1(2,2)+grid%mub(2,2),grid%mu_1(2,2), &
439                               grid%alb(2,1,2),grid%al(1,2,1),grid%rdnw(1)
440     endif
441     ENDIF
443   ENDDO
444   ENDDO
446 !#if 0
448 !  thermal perturbation to kick off convection
450   write(6,*) ' nxc, nyc for perturbation ',nxc,nyc
451   write(6,*) ' delt for perturbation ',delt
453   DO J = jts, min(jde-1,jte)
454     yrad = config_flags%dy*float(j-nyc)/10000.
455 !   yrad = 0.
456     DO I = its, min(ide-1,ite)
457       xrad = config_flags%dx*float(i-nxc)/10000.
458 !     xrad = 0.
459       DO K = 1, kte-1
461 !  put in preturbation theta (bubble) and recalc density.  note,
462 !  the mass in the column is not changing, so when theta changes,
463 !  we recompute density and geopotential
465         zrad = 0.5*(grid%ph_1(i,k,j)+grid%ph_1(i,k+1,j)  &
466                    +grid%phb(i,k,j)+grid%phb(i,k+1,j))/g
467         zrad = (zrad-1500.)/1500.
468         RAD=SQRT(xrad*xrad+yrad*yrad+zrad*zrad)
469         IF(RAD <= 1.) THEN
470            grid%t_1(i,k,j)=grid%t_1(i,k,j)+delt*COS(.5*PI*RAD)**2
471            grid%t_2(i,k,j)=grid%t_1(i,k,j)
472            qvf = 1. + rvovrd*moist(i,k,j,P_QV)
473            grid%alt(i,k,j) = (r_d/p1000mb)*(grid%t_1(i,k,j)+t0)*qvf* &
474                         (((grid%p(i,k,j)+grid%pb(i,k,j))/p1000mb)**cvpm)
475            grid%al(i,k,j) = grid%alt(i,k,j) - grid%alb(i,k,j)
476         ENDIF
477       ENDDO
479 !  rebalance hydrostatically
481       DO k  = 2,kte
482         grid%ph_1(i,k,j) = grid%ph_1(i,k-1,j) - (1./grid%rdnw(k-1))*(       &
483                      (grid%mub(i,j)+grid%mu_1(i,j))*grid%al(i,k-1,j)+ &
484                       grid%mu_1(i,j)*grid%alb(i,k-1,j)  )
485                                                    
486         grid%ph_2(i,k,j) = grid%ph_1(i,k,j) 
487         grid%ph0(i,k,j) = grid%ph_1(i,k,j) + grid%phb(i,k,j)
488       ENDDO
490     ENDDO
491   ENDDO
493 !#endif
495    IF ( wrf_dm_on_monitor() ) THEN
496    write(6,*) ' grid%mu_1 from comp ', grid%mu_1(1,1)
497    write(6,*) ' full state sounding from comp, ph, grid%p, grid%al, grid%t_1, qv '
498    do k=1,kde-1
499      write(6,'(i3,1x,5(1x,1pe10.3))') k, grid%ph_1(1,k,1)+grid%phb(1,k,1), &
500                                       grid%p(1,k,1)+grid%pb(1,k,1), grid%alt(1,k,1), &
501                                       grid%t_1(1,k,1)+t0, moist(1,k,1,P_QV)
502    enddo
504    write(6,*) ' pert state sounding from comp, grid%ph_1, pp, alp, grid%t_1, qv '
505    do k=1,kde-1
506      write(6,'(i3,1x,5(1x,1pe10.3))') k, grid%ph_1(1,k,1), &
507                                       grid%p(1,k,1), grid%al(1,k,1), &
508                                       grid%t_1(1,k,1), moist(1,k,1,P_QV)
509    enddo
510    ENDIF
512 ! interp v
514   DO J = jts, jte
515   DO I = its, min(ide-1,ite)
517     IF (j == jds) THEN
518       z_at_v = grid%phb(i,1,j)/g
519     ELSE IF (j == jde) THEN
520       z_at_v = grid%phb(i,1,j-1)/g
521     ELSE
522       z_at_v = 0.5*(grid%phb(i,1,j)+grid%phb(i,1,j-1))/g
523     END IF
524     p_surf = interp_0( p_in, zk, z_at_v, nl_in )
526     DO K = 1, kte-1
527       p_level = grid%znu(k)*(p_surf - grid%p_top) + grid%p_top
528       grid%v_1(i,k,j) = interp_0( v, p_in, p_level, nl_in )
529       grid%v_2(i,k,j) = grid%v_1(i,k,j)
530     ENDDO
532   ENDDO
533   ENDDO
535 ! interp u
537   DO J = jts, min(jde-1,jte)
538   DO I = its, ite
540     IF (i == ids) THEN
541       z_at_u = grid%phb(i,1,j)/g
542     ELSE IF (i == ide) THEN
543       z_at_u = grid%phb(i-1,1,j)/g
544     ELSE
545       z_at_u = 0.5*(grid%phb(i,1,j)+grid%phb(i-1,1,j))/g
546     END IF
548     p_surf = interp_0( p_in, zk, z_at_u, nl_in )
550     DO K = 1, kte-1
551       p_level = grid%znu(k)*(p_surf - grid%p_top) + grid%p_top
552       grid%u_1(i,k,j) = interp_0( u, p_in, p_level, nl_in )
553       grid%u_2(i,k,j) = grid%u_1(i,k,j)
554     ENDDO
556   ENDDO
557   ENDDO
559 !  set w
561   DO J = jts, min(jde-1,jte)
562   DO K = kts, kte
563   DO I = its, min(ide-1,ite)
564     grid%w_1(i,k,j) = 0.
565     grid%w_2(i,k,j) = 0.
566   ENDDO
567   ENDDO
568   ENDDO
570 !  set a few more things
572   DO J = jts, min(jde-1,jte)
573   DO K = kts, kte-1
574   DO I = its, min(ide-1,ite)
575     grid%h_diabatic(i,k,j) = 0.
576   ENDDO
577   ENDDO
578   ENDDO
580   IF ( wrf_dm_on_monitor() ) THEN
581   DO k=1,kte-1
582     grid%t_base(k) = grid%t_1(1,k,1)
583     grid%qv_base(k) = moist(1,k,1,P_QV)
584     grid%u_base(k) = grid%u_1(1,k,1)
585     grid%v_base(k) = grid%v_1(1,k,1)
586     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
587   ENDDO
588   ENDIF
589   CALL wrf_dm_bcast_real( grid%t_base , kte )
590   CALL wrf_dm_bcast_real( grid%qv_base , kte )
591   CALL wrf_dm_bcast_real( grid%u_base , kte )
592   CALL wrf_dm_bcast_real( grid%v_base , kte )
593   CALL wrf_dm_bcast_real( grid%z_base , kte )
595   DO J = jts, min(jde-1,jte)
596   DO I = its, min(ide-1,ite)
597      thtmp   = grid%t_2(i,1,j)+t0
598      ptmp    = grid%p(i,1,j)+grid%pb(i,1,j)
599      temp(1) = thtmp * (ptmp/p1000mb)**rcp
600      thtmp   = grid%t_2(i,2,j)+t0
601      ptmp    = grid%p(i,2,j)+grid%pb(i,2,j)
602      temp(2) = thtmp * (ptmp/p1000mb)**rcp
603      thtmp   = grid%t_2(i,3,j)+t0
604      ptmp    = grid%p(i,3,j)+grid%pb(i,3,j)
605      temp(3) = thtmp * (ptmp/p1000mb)**rcp
607      grid%tsk(I,J)=grid%cf1*temp(1)+grid%cf2*temp(2)+grid%cf3*temp(3)
608      grid%tmn(I,J)=grid%tsk(I,J)-0.5
609   ENDDO
610   ENDDO
612  END SUBROUTINE init_domain_rk
614    SUBROUTINE init_module_initialize
615    END SUBROUTINE init_module_initialize
617 !---------------------------------------------------------------------
619 !  test driver for get_sounding
621 !      implicit none
622 !      integer n
623 !      parameter(n = 1000)
624 !      real zk(n),p(n),theta(n),rho(n),u(n),v(n),qv(n),pd(n)
625 !      logical dry
626 !      integer nl,k
628 !      dry = .false.
629 !      dry = .true.
630 !      call get_sounding( zk, p, pd, theta, rho, u, v, qv, dry, n, nl )
631 !      write(6,*) ' input levels ',nl
632 !      write(6,*) ' sounding '
633 !      write(6,*) '  k  height(m)  press (Pa) pd(Pa) theta (K) den(kg/m^3)  u(m/s)     v(m/s)    qv(g/g) '
634 !      do k=1,nl
635 !        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)
636 !      enddo
637 !      end
639 !---------------------------------------------------------------------------
641       subroutine get_sounding( zk, p, p_dry, theta, rho, &
642                                u, v, qv, dry, nl_max, nl_in )
643       implicit none
645       integer nl_max, nl_in
646       real zk(nl_max), p(nl_max), theta(nl_max), rho(nl_max), &
647            u(nl_max), v(nl_max), qv(nl_max), p_dry(nl_max)
648       logical dry
650       integer n
651       parameter(n=1000)
652       logical debug
653       parameter( debug = .true.)
655 ! input sounding data
657       real p_surf, th_surf, qv_surf
658       real pi_surf, pi(n)
659       real h_input(n), th_input(n), qv_input(n), u_input(n), v_input(n)
661 ! diagnostics
663       real rho_surf, p_input(n), rho_input(n)
664       real pm_input(n)  !  this are for full moist sounding
666 ! local data
668       real p1000mb,cv,cp,r,cvpm,g
669       parameter (p1000mb = 1.e+05, r = 287, cp = 1003., cv = cp-r, cvpm = -cv/cp, g=9.81 )
670       integer k, it, nl
671       real qvf, qvf1, dz
673 !  first, read the sounding
675       call read_sounding( p_surf, th_surf, qv_surf, &
676                           h_input, th_input, qv_input, u_input, v_input,n, nl, debug )
678       if(dry) then
679        do k=1,nl
680          qv_input(k) = 0.
681        enddo
682       endif
684       if(debug) write(6,*) ' number of input levels = ',nl
686         nl_in = nl
687         if(nl_in .gt. nl_max ) then
688           write(6,*) ' too many levels for input arrays ',nl_in,nl_max
689           call wrf_error_fatal ( ' too many levels for input arrays ' )
690         end if
692 !  compute diagnostics,
693 !  first, convert qv(g/kg) to qv(g/g)
695       do k=1,nl
696         qv_input(k) = 0.001*qv_input(k)
697       enddo
699       p_surf = 100.*p_surf  ! convert to pascals
700       qvf = 1. + rvovrd*qv_input(1) 
701       rho_surf = 1./((r/p1000mb)*th_surf*qvf*((p_surf/p1000mb)**cvpm))
702       pi_surf = (p_surf/p1000mb)**(r/cp)
704       if(debug) then
705         write(6,*) ' surface density is ',rho_surf
706         write(6,*) ' surface pi is      ',pi_surf
707       end if
710 !  integrate moist sounding hydrostatically, starting from the
711 !  specified surface pressure
712 !  -> first, integrate from surface to lowest level
714           qvf = 1. + rvovrd*qv_input(1) 
715           qvf1 = 1. + qv_input(1)
716           rho_input(1) = rho_surf
717           dz = h_input(1)
718           do it=1,10
719             pm_input(1) = p_surf &
720                     - 0.5*dz*(rho_surf+rho_input(1))*g*qvf1
721             rho_input(1) = 1./((r/p1000mb)*th_input(1)*qvf*((pm_input(1)/p1000mb)**cvpm))
722           enddo
724 ! integrate up the column
726           do k=2,nl
727             rho_input(k) = rho_input(k-1)
728             dz = h_input(k)-h_input(k-1)
729             qvf1 = 0.5*(2.+(qv_input(k-1)+qv_input(k)))
730             qvf = 1. + rvovrd*qv_input(k)   ! qv is in g/kg here
732             do it=1,10
733               pm_input(k) = pm_input(k-1) &
734                       - 0.5*dz*(rho_input(k)+rho_input(k-1))*g*qvf1
735               rho_input(k) = 1./((r/p1000mb)*th_input(k)*qvf*((pm_input(k)/p1000mb)**cvpm))
736             enddo
737           enddo
739 !  we have the moist sounding
741 !  next, compute the dry sounding using p at the highest level from the
742 !  moist sounding and integrating down.
744         p_input(nl) = pm_input(nl)
746           do k=nl-1,1,-1
747             dz = h_input(k+1)-h_input(k)
748             p_input(k) = p_input(k+1) + 0.5*dz*(rho_input(k)+rho_input(k+1))*g
749           enddo
752         do k=1,nl
754           zk(k) = h_input(k)
755           p(k) = pm_input(k)
756           p_dry(k) = p_input(k)
757           theta(k) = th_input(k)
758           rho(k) = rho_input(k)
759           u(k) = u_input(k)
760           v(k) = v_input(k)
761           qv(k) = qv_input(k)
763         enddo
765      if(debug) then
766       write(6,*) ' sounding '
767       write(6,*) '  k  height(m)  press (Pa) pd(Pa) theta (K) den(kg/m^3)  u(m/s)     v(m/s)    qv(g/g) '
768       do k=1,nl
769         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)
770       enddo
772      end if
774       end subroutine get_sounding
776 !-------------------------------------------------------
778       subroutine read_sounding( ps,ts,qvs,h,th,qv,u,v,n,nl,debug )
779       implicit none
780       integer n,nl
781       real ps,ts,qvs,h(n),th(n),qv(n),u(n),v(n)
782       logical end_of_file
783       logical debug
785       integer k
787       open(unit=10,file='input_sounding',form='formatted',status='old')
788       rewind(10)
789       read(10,*) ps, ts, qvs
790       if(debug) then
791         write(6,*) ' input sounding surface parameters '
792         write(6,*) ' surface pressure (mb) ',ps
793         write(6,*) ' surface pot. temp (K) ',ts
794         write(6,*) ' surface mixing ratio (g/kg) ',qvs
795       end if
797       end_of_file = .false.
798       k = 0
800       do while (.not. end_of_file)
802         read(10,*,end=100) h(k+1), th(k+1), qv(k+1), u(k+1), v(k+1)
803         k = k+1
804         if(debug) write(6,'(1x,i3,5(1x,e10.3))') k, h(k), th(k), qv(k), u(k), v(k)
805         go to 110
806  100    end_of_file = .true.
807  110    continue
808       enddo
810       nl = k
812       close(unit=10,status = 'keep')
814       end subroutine read_sounding
816 END MODULE module_initialize_ideal