module_initialize_quarter_ss.F

References to this file elsewhere.
1 !IDEAL:MODEL_LAYER:INITIALIZATION
2 !
3 
4 !  This MODULE holds the routines which are used to perform various initializations
5 !  for the individual domains.  
6 
7 !  This MODULE CONTAINS the following routines:
8 
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.               
13 
14 !-----------------------------------------------------------------------
15 
16 MODULE module_initialize_ideal
17 
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
29 
30 
31 CONTAINS
32 
33 
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.
42 
43    SUBROUTINE init_domain ( grid )
44 
45    IMPLICIT NONE
46 
47    !  Input data.
48    TYPE (domain), POINTER :: grid 
49    !  Local data.
50    INTEGER                :: dyn_opt 
51    INTEGER :: idum1, idum2
52 
53    CALL nl_get_dyn_opt( 1, dyn_opt )
54    
55    CALL set_scalar_indices_from_config ( head_grid%id , idum1, idum2 )
56 
57    IF (      dyn_opt .eq. 1 &
58         .or. dyn_opt .eq. 2 &
59         .or. dyn_opt .eq. 3 &
60                                        ) THEN
61      CALL init_domain_rk( grid &
62 !
63 #include <em_actual_new_args.inc>
64 !
65                         )
66 
67    ELSE
68      WRITE(0,*)' init_domain: unknown or unimplemented dyn_opt = ',dyn_opt
69      CALL wrf_error_fatal ( ' init_domain: unknown or unimplemented dyn_opt ' )
70    ENDIF
71 
72    END SUBROUTINE init_domain
73 
74 !-------------------------------------------------------------------
75 
76    SUBROUTINE init_domain_rk ( grid &
77 !
78 # include <em_dummy_new_args.inc>
79 !
80 )
81    IMPLICIT NONE
82 
83    !  Input data.
84    TYPE (domain), POINTER :: grid
85 
86 # include <em_dummy_new_decl.inc>
87 
88    TYPE (grid_config_rec_type)              :: config_flags
89 
90    !  Local data
91    INTEGER                             ::                       &
92                                   ids, ide, jds, jde, kds, kde, &
93                                   ims, ime, jms, jme, kms, kme, &
94                                   its, ite, jts, jte, kts, kte, &
95                                   i, j, k
96 
97    ! Local data
98 
99    INTEGER, PARAMETER :: nl_max = 1000
100    REAL, DIMENSION(nl_max) :: zk, p_in, theta, rho, u, v, qv, pd_in
101    INTEGER :: nl_in
102 
103 
104    INTEGER :: icm,jcm, ii, im1, jj, jm1, loop, error, fid, nxc, nyc
105    REAL    :: u_mean,v_mean, f0, p_surf, p_level, qvf, z_at_v, z_at_u
106    REAL    :: z_scale, xrad, yrad, zrad, rad, delt, cof1, cof2
107 !   REAL, EXTERNAL :: interp_0
108    REAL    :: hm
109    REAL    :: pi
110 
111 !  stuff from original initialization that has been dropped from the Registry 
112    REAL    :: vnu, xnu, xnus, dinit0, cbh, p0_temp, t0_temp, zd, zt
113    REAL    :: qvf1, qvf2, pd_surf
114    INTEGER :: it
115    real :: thtmp, ptmp, temp(3)
116 
117    LOGICAL :: moisture_init
118    LOGICAL :: stretch_grid, dry_sounding
119 
120   INTEGER :: xs , xe , ys , ye
121   REAL :: mtn_ht
122    LOGICAL, EXTERNAL :: wrf_dm_on_monitor
123 
124 #ifdef DM_PARALLEL
125 #    include <em_data_calls.inc>
126 #endif
127 
128 
129    SELECT CASE ( model_data_order )
130          CASE ( DATA_ORDER_ZXY )
131    kds = grid%sd31 ; kde = grid%ed31 ;
132    ids = grid%sd32 ; ide = grid%ed32 ;
133    jds = grid%sd33 ; jde = grid%ed33 ;
134 
135    kms = grid%sm31 ; kme = grid%em31 ;
136    ims = grid%sm32 ; ime = grid%em32 ;
137    jms = grid%sm33 ; jme = grid%em33 ;
138 
139    kts = grid%sp31 ; kte = grid%ep31 ;   ! note that tile is entire patch
140    its = grid%sp32 ; ite = grid%ep32 ;   ! note that tile is entire patch
141    jts = grid%sp33 ; jte = grid%ep33 ;   ! note that tile is entire patch
142          CASE ( DATA_ORDER_XYZ )
143    ids = grid%sd31 ; ide = grid%ed31 ;
144    jds = grid%sd32 ; jde = grid%ed32 ;
145    kds = grid%sd33 ; kde = grid%ed33 ;
146 
147    ims = grid%sm31 ; ime = grid%em31 ;
148    jms = grid%sm32 ; jme = grid%em32 ;
149    kms = grid%sm33 ; kme = grid%em33 ;
150 
151    its = grid%sp31 ; ite = grid%ep31 ;   ! note that tile is entire patch
152    jts = grid%sp32 ; jte = grid%ep32 ;   ! note that tile is entire patch
153    kts = grid%sp33 ; kte = grid%ep33 ;   ! note that tile is entire patch
154          CASE ( DATA_ORDER_XZY )
155    ids = grid%sd31 ; ide = grid%ed31 ;
156    kds = grid%sd32 ; kde = grid%ed32 ;
157    jds = grid%sd33 ; jde = grid%ed33 ;
158 
159    ims = grid%sm31 ; ime = grid%em31 ;
160    kms = grid%sm32 ; kme = grid%em32 ;
161    jms = grid%sm33 ; jme = grid%em33 ;
162 
163    its = grid%sp31 ; ite = grid%ep31 ;   ! note that tile is entire patch
164    kts = grid%sp32 ; kte = grid%ep32 ;   ! note that tile is entire patch
165    jts = grid%sp33 ; jte = grid%ep33 ;   ! note that tile is entire patch
166 
167    END SELECT
168 
169 
170    stretch_grid = .true.
171    delt = 3.
172 !   z_scale = .50
173    z_scale = .40
174    pi = 2.*asin(1.0)
175    write(6,*) ' pi is ',pi
176    nxc = (ide-ids)/2
177    nyc = (jde-jds)/2
178 
179    CALL model_to_grid_config_rec ( grid%id , model_config_rec , config_flags )
180 
181 ! here we check to see if the boundary conditions are set properly
182 
183    CALL boundary_condition_check( config_flags, bdyzone, error, grid%id )
184 
185    moisture_init = .true.
186 
187     grid%itimestep=0
188 
189 #ifdef DM_PARALLEL
190    CALL wrf_dm_bcast_bytes( icm , IWORDSIZE )
191    CALL wrf_dm_bcast_bytes( jcm , IWORDSIZE )
192 #endif
193 
194     CALL nl_set_mminlu(1, '    ')
195     CALL nl_set_iswater(1,0)
196     CALL nl_set_cen_lat(1,40.)
197     CALL nl_set_cen_lon(1,-105.)
198     CALL nl_set_truelat1(1,0.)
199     CALL nl_set_truelat2(1,0.)
200     CALL nl_set_moad_cen_lat (1,0.)
201     CALL nl_set_stand_lon (1,0.)
202     CALL nl_set_map_proj(1,0)
203 
204 
205 !  here we initialize data we currently is not initialized 
206 !  in the input data
207 
208     DO j = jts, jte
209       DO i = its, ite
210          grid%msft(i,j)     = 1.
211          grid%msfu(i,j)     = 1.
212          grid%msfv(i,j)     = 1.
213          grid%sina(i,j)     = 0.
214          grid%cosa(i,j)     = 1.
215          grid%e(i,j)        = 0.
216          grid%f(i,j)        = 0.
217 
218       END DO
219    END DO
220 
221     DO j = jts, jte
222     DO k = kts, kte
223       DO i = its, ite
224          grid%em_ww(i,k,j)     = 0.
225       END DO
226    END DO
227    END DO
228 
229    grid%step_number = 0
230 
231 ! set up the grid
232 
233    IF (stretch_grid) THEN ! exponential stretch for eta (nearly constant dz)
234      DO k=1, kde
235       grid%em_znw(k) = (exp(-(k-1)/float(kde-1)/z_scale) - exp(-1./z_scale))/ &
236                                 (1.-exp(-1./z_scale))
237      ENDDO
238    ELSE
239      DO k=1, kde
240       grid%em_znw(k) = 1. - float(k-1)/float(kde-1)
241      ENDDO
242    ENDIF
243 
244    DO k=1, kde-1
245     grid%em_dnw(k) = grid%em_znw(k+1) - grid%em_znw(k)
246     grid%em_rdnw(k) = 1./grid%em_dnw(k)
247     grid%em_znu(k) = 0.5*(grid%em_znw(k+1)+grid%em_znw(k))
248    ENDDO
249    DO k=2, kde-1
250     grid%em_dn(k) = 0.5*(grid%em_dnw(k)+grid%em_dnw(k-1))
251     grid%em_rdn(k) = 1./grid%em_dn(k)
252     grid%em_fnp(k) = .5* grid%em_dnw(k  )/grid%em_dn(k)
253     grid%em_fnm(k) = .5* grid%em_dnw(k-1)/grid%em_dn(k)
254    ENDDO
255 
256    cof1 = (2.*grid%em_dn(2)+grid%em_dn(3))/(grid%em_dn(2)+grid%em_dn(3))*grid%em_dnw(1)/grid%em_dn(2) 
257    cof2 =     grid%em_dn(2)        /(grid%em_dn(2)+grid%em_dn(3))*grid%em_dnw(1)/grid%em_dn(3) 
258    grid%cf1  = grid%em_fnp(2) + cof1
259    grid%cf2  = grid%em_fnm(2) - cof1 - cof2
260    grid%cf3  = cof2       
261 
262    grid%cfn  = (.5*grid%em_dnw(kde-1)+grid%em_dn(kde-1))/grid%em_dn(kde-1)
263    grid%cfn1 = -.5*grid%em_dnw(kde-1)/grid%em_dn(kde-1)
264    grid%rdx = 1./config_flags%dx
265    grid%rdy = 1./config_flags%dy
266 
267 !  get the sounding from the ascii sounding file, first get dry sounding and 
268 !  calculate base state
269 
270   dry_sounding = .true.
271   IF ( wrf_dm_on_monitor() ) THEN
272   write(6,*) ' getting dry sounding for base state '
273 
274   CALL get_sounding( zk, p_in, pd_in, theta, rho, u, v, qv, dry_sounding, nl_max, nl_in )
275   ENDIF
276   CALL wrf_dm_bcast_real( zk , nl_max )
277   CALL wrf_dm_bcast_real( p_in , nl_max )
278   CALL wrf_dm_bcast_real( pd_in , nl_max )
279   CALL wrf_dm_bcast_real( theta , nl_max )
280   CALL wrf_dm_bcast_real( rho , nl_max )
281   CALL wrf_dm_bcast_real( u , nl_max )
282   CALL wrf_dm_bcast_real( v , nl_max )
283   CALL wrf_dm_bcast_real( qv , nl_max )
284   CALL wrf_dm_bcast_integer ( nl_in , 1 ) 
285 
286   write(6,*) ' returned from reading sounding, nl_in is ',nl_in
287 
288 !  find ptop for the desired ztop (ztop is input from the namelist),
289 !  and find surface pressure
290 
291   grid%p_top = interp_0( p_in, zk, config_flags%ztop, nl_in )
292 
293   DO j=jts,jte
294   DO i=its,ite
295     grid%ht(i,j) = 0.
296   ENDDO
297   ENDDO
298 
299   xs=ide/2 -3
300   xs=ids   -3
301   xe=xs + 6
302   ys=jde/2 -3
303   ye=ys + 6
304   mtn_ht = 500
305 #ifdef MTN
306   DO j=max(ys,jds),min(ye,jde-1)
307   DO i=max(xs,ids),min(xe,ide-1)
308      grid%ht(i,j) = mtn_ht * 0.25 * &
309                ( 1. + COS ( 2*pi/(xe-xs) * ( i-xs ) + pi ) ) * &
310                ( 1. + COS ( 2*pi/(ye-ys) * ( j-ys ) + pi ) )
311   ENDDO
312   ENDDO
313 #endif
314 #ifdef EW_RIDGE
315   DO j=max(ys,jds),min(ye,jde-1)
316   DO i=ids,ide
317      grid%ht(i,j) = mtn_ht * 0.50 * &
318                ( 1. + COS ( 2*pi/(ye-ys) * ( j-ys ) + pi ) )
319   ENDDO
320   ENDDO
321 #endif
322 #ifdef NS_RIDGE
323   DO j=jds,jde
324   DO i=max(xs,ids),min(xe,ide-1)
325      grid%ht(i,j) = mtn_ht * 0.50 * &
326                ( 1. + COS ( 2*pi/(xe-xs) * ( i-xs ) + pi ) )
327   ENDDO
328   ENDDO
329 #endif
330   DO j=jts,jte
331   DO i=its,ite
332     grid%em_phb(i,1,j) = g * grid%ht(i,j)
333     grid%em_ph0(i,1,j) = g * grid%ht(i,j)
334   ENDDO
335   ENDDO
336 
337   DO J = jts, jte
338   DO I = its, ite
339 
340     p_surf = interp_0( p_in, zk, grid%em_phb(i,1,j)/g, nl_in )
341     grid%em_mub(i,j) = p_surf-grid%p_top
342 
343 !  this is dry hydrostatic sounding (base state), so given grid%em_p (coordinate),
344 !  interp theta (from interp) and compute 1/rho from eqn. of state
345 
346     DO K = 1, kte-1
347       p_level = grid%em_znu(k)*(p_surf - grid%p_top) + grid%p_top
348       grid%em_pb(i,k,j) = p_level
349       grid%em_t_init(i,k,j) = interp_0( theta, p_in, p_level, nl_in ) - t0
350       grid%em_alb(i,k,j) = (r_d/p1000mb)*(grid%em_t_init(i,k,j)+t0)*(grid%em_pb(i,k,j)/p1000mb)**cvpm
351     ENDDO
352 
353 !  calc hydrostatic balance (alternatively we could interp the geopotential from the
354 !  sounding, but this assures that the base state is in exact hydrostatic balance with
355 !  respect to the model eqns.
356 
357     DO k  = 2,kte
358       grid%em_phb(i,k,j) = grid%em_phb(i,k-1,j) - grid%em_dnw(k-1)*grid%em_mub(i,j)*grid%em_alb(i,k-1,j)
359     ENDDO
360 
361   ENDDO
362   ENDDO
363 
364   IF ( wrf_dm_on_monitor() ) THEN
365     write(6,*) ' ptop is ',grid%p_top
366     write(6,*) ' base state grid%em_mub(1,1), p_surf is ',grid%em_mub(1,1),grid%em_mub(1,1)+grid%p_top
367   ENDIF
368 
369 !  calculate full state for each column - this includes moisture.
370 
371   write(6,*) ' getting moist sounding for full state '
372   dry_sounding = .false.
373   CALL get_sounding( zk, p_in, pd_in, theta, rho, u, v, qv, dry_sounding, nl_max, nl_in )
374 
375   DO J = jts, min(jde-1,jte)
376   DO I = its, min(ide-1,ite)
377 
378 !  At this point grid%p_top is already set. find the DRY mass in the column 
379 !  by interpolating the DRY pressure.  
380 
381    pd_surf = interp_0( pd_in, zk, grid%em_phb(i,1,j)/g, nl_in )
382 
383 !  compute the perturbation mass and the full mass
384 
385     grid%em_mu_1(i,j) = pd_surf-grid%p_top - grid%em_mub(i,j)
386     grid%em_mu_2(i,j) = grid%em_mu_1(i,j)
387     grid%em_mu0(i,j) = grid%em_mu_1(i,j) + grid%em_mub(i,j)
388 
389 ! given the dry pressure and coordinate system, interp the potential
390 ! temperature and qv
391 
392     do k=1,kde-1
393 
394       p_level = grid%em_znu(k)*(pd_surf - grid%p_top) + grid%p_top
395 
396       moist(i,k,j,P_QV) = interp_0( qv, pd_in, p_level, nl_in )
397       grid%em_t_1(i,k,j)          = interp_0( theta, pd_in, p_level, nl_in ) - t0
398       grid%em_t_2(i,k,j)          = grid%em_t_1(i,k,j)
399       
400 
401     enddo
402 
403 !  integrate the hydrostatic equation (from the RHS of the bigstep
404 !  vertical momentum equation) down from the top to get grid%em_p.
405 !  first from the top of the model to the top pressure
406 
407     k = kte-1  ! top level
408 
409     qvf1 = 0.5*(moist(i,k,j,P_QV)+moist(i,k,j,P_QV))
410     qvf2 = 1./(1.+qvf1)
411     qvf1 = qvf1*qvf2
412 
413 !    grid%em_p(i,k,j) = - 0.5*grid%em_mu_1(i,j)/grid%em_rdnw(k)
414     grid%em_p(i,k,j) = - 0.5*(grid%em_mu_1(i,j)+qvf1*grid%em_mub(i,j))/grid%em_rdnw(k)/qvf2
415     qvf = 1. + rvovrd*moist(i,k,j,P_QV)
416     grid%em_alt(i,k,j) = (r_d/p1000mb)*(grid%em_t_1(i,k,j)+t0)*qvf* &
417                 (((grid%em_p(i,k,j)+grid%em_pb(i,k,j))/p1000mb)**cvpm)
418     grid%em_al(i,k,j) = grid%em_alt(i,k,j) - grid%em_alb(i,k,j)
419 
420 !  down the column
421 
422     do k=kte-2,1,-1
423       qvf1 = 0.5*(moist(i,k,j,P_QV)+moist(i,k+1,j,P_QV))
424       qvf2 = 1./(1.+qvf1)
425       qvf1 = qvf1*qvf2
426       grid%em_p(i,k,j) = grid%em_p(i,k+1,j) - (grid%em_mu_1(i,j) + qvf1*grid%em_mub(i,j))/qvf2/grid%em_rdn(k+1)
427       qvf = 1. + rvovrd*moist(i,k,j,P_QV)
428       grid%em_alt(i,k,j) = (r_d/p1000mb)*(grid%em_t_1(i,k,j)+t0)*qvf* &
429                   (((grid%em_p(i,k,j)+grid%em_pb(i,k,j))/p1000mb)**cvpm)
430       grid%em_al(i,k,j) = grid%em_alt(i,k,j) - grid%em_alb(i,k,j)
431     enddo
432 
433 !  this is the hydrostatic equation used in the model after the
434 !  small timesteps.  In the model, grid%em_al (inverse density)
435 !  is computed from the geopotential.
436 
437 
438     grid%em_ph_1(i,1,j) = 0.
439     DO k  = 2,kte
440       grid%em_ph_1(i,k,j) = grid%em_ph_1(i,k-1,j) - (1./grid%em_rdnw(k-1))*(       &
441                    (grid%em_mub(i,j)+grid%em_mu_1(i,j))*grid%em_al(i,k-1,j)+ &
442                     grid%em_mu_1(i,j)*grid%em_alb(i,k-1,j)  )
443                                                    
444       grid%em_ph_2(i,k,j) = grid%em_ph_1(i,k,j) 
445       grid%em_ph0(i,k,j) = grid%em_ph_1(i,k,j) + grid%em_phb(i,k,j)
446     ENDDO
447 
448     IF ( wrf_dm_on_monitor() ) THEN
449     if((i==2) .and. (j==2)) then
450      write(6,*) ' grid%em_ph_1 calc ',grid%em_ph_1(2,1,2),grid%em_ph_1(2,2,2),&
451                               grid%em_mu_1(2,2)+grid%em_mub(2,2),grid%em_mu_1(2,2), &
452                               grid%em_alb(2,1,2),grid%em_al(1,2,1),grid%em_rdnw(1)
453     endif
454     ENDIF
455 
456   ENDDO
457   ENDDO
458 
459 !#if 0
460 
461 !  thermal perturbation to kick off convection
462 
463   write(6,*) ' nxc, nyc for perturbation ',nxc,nyc
464   write(6,*) ' delt for perturbation ',delt
465 
466   DO J = jts, min(jde-1,jte)
467     yrad = config_flags%dy*float(j-nyc)/10000.
468 !   yrad = 0.
469     DO I = its, min(ide-1,ite)
470       xrad = config_flags%dx*float(i-nxc)/10000.
471 !     xrad = 0.
472       DO K = 1, kte-1
473 
474 !  put in preturbation theta (bubble) and recalc density.  note,
475 !  the mass in the column is not changing, so when theta changes,
476 !  we recompute density and geopotential
477 
478         zrad = 0.5*(grid%em_ph_1(i,k,j)+grid%em_ph_1(i,k+1,j)  &
479                    +grid%em_phb(i,k,j)+grid%em_phb(i,k+1,j))/g
480         zrad = (zrad-1500.)/1500.
481         RAD=SQRT(xrad*xrad+yrad*yrad+zrad*zrad)
482         IF(RAD <= 1.) THEN
483            grid%em_t_1(i,k,j)=grid%em_t_1(i,k,j)+delt*COS(.5*PI*RAD)**2
484            grid%em_t_2(i,k,j)=grid%em_t_1(i,k,j)
485            qvf = 1. + rvovrd*moist(i,k,j,P_QV)
486            grid%em_alt(i,k,j) = (r_d/p1000mb)*(grid%em_t_1(i,k,j)+t0)*qvf* &
487                         (((grid%em_p(i,k,j)+grid%em_pb(i,k,j))/p1000mb)**cvpm)
488            grid%em_al(i,k,j) = grid%em_alt(i,k,j) - grid%em_alb(i,k,j)
489         ENDIF
490       ENDDO
491 
492 !  rebalance hydrostatically
493 
494       DO k  = 2,kte
495         grid%em_ph_1(i,k,j) = grid%em_ph_1(i,k-1,j) - (1./grid%em_rdnw(k-1))*(       &
496                      (grid%em_mub(i,j)+grid%em_mu_1(i,j))*grid%em_al(i,k-1,j)+ &
497                       grid%em_mu_1(i,j)*grid%em_alb(i,k-1,j)  )
498                                                    
499         grid%em_ph_2(i,k,j) = grid%em_ph_1(i,k,j) 
500         grid%em_ph0(i,k,j) = grid%em_ph_1(i,k,j) + grid%em_phb(i,k,j)
501       ENDDO
502 
503     ENDDO
504   ENDDO
505 
506 !#endif
507 
508    IF ( wrf_dm_on_monitor() ) THEN
509    write(6,*) ' grid%em_mu_1 from comp ', grid%em_mu_1(1,1)
510    write(6,*) ' full state sounding from comp, ph, grid%em_p, grid%em_al, grid%em_t_1, qv '
511    do k=1,kde-1
512      write(6,'(i3,1x,5(1x,1pe10.3))') k, grid%em_ph_1(1,k,1)+grid%em_phb(1,k,1), &
513                                       grid%em_p(1,k,1)+grid%em_pb(1,k,1), grid%em_alt(1,k,1), &
514                                       grid%em_t_1(1,k,1)+t0, moist(1,k,1,P_QV)
515    enddo
516 
517    write(6,*) ' pert state sounding from comp, grid%em_ph_1, pp, alp, grid%em_t_1, qv '
518    do k=1,kde-1
519      write(6,'(i3,1x,5(1x,1pe10.3))') k, grid%em_ph_1(1,k,1), &
520                                       grid%em_p(1,k,1), grid%em_al(1,k,1), &
521                                       grid%em_t_1(1,k,1), moist(1,k,1,P_QV)
522    enddo
523    ENDIF
524 
525 ! interp v
526 
527   DO J = jts, jte
528   DO I = its, min(ide-1,ite)
529 
530     IF (j == jds) THEN
531       z_at_v = grid%em_phb(i,1,j)/g
532     ELSE IF (j == jde) THEN
533       z_at_v = grid%em_phb(i,1,j-1)/g
534     ELSE
535       z_at_v = 0.5*(grid%em_phb(i,1,j)+grid%em_phb(i,1,j-1))/g
536     END IF
537     p_surf = interp_0( p_in, zk, z_at_v, nl_in )
538 
539     DO K = 1, kte-1
540       p_level = grid%em_znu(k)*(p_surf - grid%p_top) + grid%p_top
541       grid%em_v_1(i,k,j) = interp_0( v, p_in, p_level, nl_in )
542       grid%em_v_2(i,k,j) = grid%em_v_1(i,k,j)
543     ENDDO
544 
545   ENDDO
546   ENDDO
547 
548 ! interp u
549 
550   DO J = jts, min(jde-1,jte)
551   DO I = its, ite
552 
553     IF (i == ids) THEN
554       z_at_u = grid%em_phb(i,1,j)/g
555     ELSE IF (i == ide) THEN
556       z_at_u = grid%em_phb(i-1,1,j)/g
557     ELSE
558       z_at_u = 0.5*(grid%em_phb(i,1,j)+grid%em_phb(i-1,1,j))/g
559     END IF
560 
561     p_surf = interp_0( p_in, zk, z_at_u, nl_in )
562 
563     DO K = 1, kte-1
564       p_level = grid%em_znu(k)*(p_surf - grid%p_top) + grid%p_top
565       grid%em_u_1(i,k,j) = interp_0( u, p_in, p_level, nl_in )
566       grid%em_u_2(i,k,j) = grid%em_u_1(i,k,j)
567     ENDDO
568 
569   ENDDO
570   ENDDO
571 
572 !  set w
573 
574   DO J = jts, min(jde-1,jte)
575   DO K = kts, kte
576   DO I = its, min(ide-1,ite)
577     grid%em_w_1(i,k,j) = 0.
578     grid%em_w_2(i,k,j) = 0.
579   ENDDO
580   ENDDO
581   ENDDO
582 
583 !  set a few more things
584 
585   DO J = jts, min(jde-1,jte)
586   DO K = kts, kte-1
587   DO I = its, min(ide-1,ite)
588     grid%h_diabatic(i,k,j) = 0.
589   ENDDO
590   ENDDO
591   ENDDO
592 
593   IF ( wrf_dm_on_monitor() ) THEN
594   DO k=1,kte-1
595     grid%em_t_base(k) = grid%em_t_1(1,k,1)
596     grid%qv_base(k) = moist(1,k,1,P_QV)
597     grid%u_base(k) = grid%em_u_1(1,k,1)
598     grid%v_base(k) = grid%em_v_1(1,k,1)
599     grid%z_base(k) = 0.5*(grid%em_phb(1,k,1)+grid%em_phb(1,k+1,1)+grid%em_ph_1(1,k,1)+grid%em_ph_1(1,k+1,1))/g
600   ENDDO
601   ENDIF
602   CALL wrf_dm_bcast_real( grid%em_t_base , kte )
603   CALL wrf_dm_bcast_real( grid%qv_base , kte )
604   CALL wrf_dm_bcast_real( grid%u_base , kte )
605   CALL wrf_dm_bcast_real( grid%v_base , kte )
606   CALL wrf_dm_bcast_real( grid%z_base , kte )
607 
608   DO J = jts, min(jde-1,jte)
609   DO I = its, min(ide-1,ite)
610      thtmp   = grid%em_t_2(i,1,j)+t0
611      ptmp    = grid%em_p(i,1,j)+grid%em_pb(i,1,j)
612      temp(1) = thtmp * (ptmp/p1000mb)**rcp
613      thtmp   = grid%em_t_2(i,2,j)+t0
614      ptmp    = grid%em_p(i,2,j)+grid%em_pb(i,2,j)
615      temp(2) = thtmp * (ptmp/p1000mb)**rcp
616      thtmp   = grid%em_t_2(i,3,j)+t0
617      ptmp    = grid%em_p(i,3,j)+grid%em_pb(i,3,j)
618      temp(3) = thtmp * (ptmp/p1000mb)**rcp
619 
620      grid%tsk(I,J)=grid%cf1*temp(1)+grid%cf2*temp(2)+grid%cf3*temp(3)
621      grid%tmn(I,J)=grid%tsk(I,J)-0.5
622   ENDDO
623   ENDDO
624 
625  END SUBROUTINE init_domain_rk
626 
627    SUBROUTINE init_module_initialize
628    END SUBROUTINE init_module_initialize
629 
630 !---------------------------------------------------------------------
631 
632 !  test driver for get_sounding
633 !
634 !      implicit none
635 !      integer n
636 !      parameter(n = 1000)
637 !      real zk(n),p(n),theta(n),rho(n),u(n),v(n),qv(n),pd(n)
638 !      logical dry
639 !      integer nl,k
640 !
641 !      dry = .false.
642 !      dry = .true.
643 !      call get_sounding( zk, p, pd, theta, rho, u, v, qv, dry, n, nl )
644 !      write(6,*) ' input levels ',nl
645 !      write(6,*) ' sounding '
646 !      write(6,*) '  k  height(m)  press (Pa) pd(Pa) theta (K) den(kg/m^3)  u(m/s)     v(m/s)    qv(g/g) '
647 !      do k=1,nl
648 !        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)
649 !      enddo
650 !      end
651 !
652 !---------------------------------------------------------------------------
653 
654       subroutine get_sounding( zk, p, p_dry, theta, rho, &
655                                u, v, qv, dry, nl_max, nl_in )
656       implicit none
657 
658       integer nl_max, nl_in
659       real zk(nl_max), p(nl_max), theta(nl_max), rho(nl_max), &
660            u(nl_max), v(nl_max), qv(nl_max), p_dry(nl_max)
661       logical dry
662 
663       integer n
664       parameter(n=1000)
665       logical debug
666       parameter( debug = .true.)
667 
668 ! input sounding data
669 
670       real p_surf, th_surf, qv_surf
671       real pi_surf, pi(n)
672       real h_input(n), th_input(n), qv_input(n), u_input(n), v_input(n)
673 
674 ! diagnostics
675 
676       real rho_surf, p_input(n), rho_input(n)
677       real pm_input(n)  !  this are for full moist sounding
678 
679 ! local data
680 
681       real p1000mb,cv,cp,r,cvpm,g
682       parameter (p1000mb = 1.e+05, r = 287, cp = 1003., cv = cp-r, cvpm = -cv/cp, g=9.81 )
683       integer k, it, nl
684       real qvf, qvf1, dz
685 
686 !  first, read the sounding
687 
688       call read_sounding( p_surf, th_surf, qv_surf, &
689                           h_input, th_input, qv_input, u_input, v_input,n, nl, debug )
690 
691       if(dry) then
692        do k=1,nl
693          qv_input(k) = 0.
694        enddo
695       endif
696 
697       if(debug) write(6,*) ' number of input levels = ',nl
698 
699         nl_in = nl
700         if(nl_in .gt. nl_max ) then
701           write(6,*) ' too many levels for input arrays ',nl_in,nl_max
702           call wrf_error_fatal ( ' too many levels for input arrays ' )
703         end if
704 
705 !  compute diagnostics,
706 !  first, convert qv(g/kg) to qv(g/g)
707 
708       do k=1,nl
709         qv_input(k) = 0.001*qv_input(k)
710       enddo
711 
712       p_surf = 100.*p_surf  ! convert to pascals
713       qvf = 1. + rvovrd*qv_input(1) 
714       rho_surf = 1./((r/p1000mb)*th_surf*qvf*((p_surf/p1000mb)**cvpm))
715       pi_surf = (p_surf/p1000mb)**(r/cp)
716 
717       if(debug) then
718         write(6,*) ' surface density is ',rho_surf
719         write(6,*) ' surface pi is      ',pi_surf
720       end if
721 
722 
723 !  integrate moist sounding hydrostatically, starting from the
724 !  specified surface pressure
725 !  -> first, integrate from surface to lowest level
726 
727           qvf = 1. + rvovrd*qv_input(1) 
728           qvf1 = 1. + qv_input(1)
729           rho_input(1) = rho_surf
730           dz = h_input(1)
731           do it=1,10
732             pm_input(1) = p_surf &
733                     - 0.5*dz*(rho_surf+rho_input(1))*g*qvf1
734             rho_input(1) = 1./((r/p1000mb)*th_input(1)*qvf*((pm_input(1)/p1000mb)**cvpm))
735           enddo
736 
737 ! integrate up the column
738 
739           do k=2,nl
740             rho_input(k) = rho_input(k-1)
741             dz = h_input(k)-h_input(k-1)
742             qvf1 = 0.5*(2.+(qv_input(k-1)+qv_input(k)))
743             qvf = 1. + rvovrd*qv_input(k)   ! qv is in g/kg here
744  
745             do it=1,10
746               pm_input(k) = pm_input(k-1) &
747                       - 0.5*dz*(rho_input(k)+rho_input(k-1))*g*qvf1
748               rho_input(k) = 1./((r/p1000mb)*th_input(k)*qvf*((pm_input(k)/p1000mb)**cvpm))
749             enddo
750           enddo
751 
752 !  we have the moist sounding
753 
754 !  next, compute the dry sounding using p at the highest level from the
755 !  moist sounding and integrating down.
756 
757         p_input(nl) = pm_input(nl)
758 
759           do k=nl-1,1,-1
760             dz = h_input(k+1)-h_input(k)
761             p_input(k) = p_input(k+1) + 0.5*dz*(rho_input(k)+rho_input(k+1))*g
762           enddo
763 
764 
765         do k=1,nl
766 
767           zk(k) = h_input(k)
768           p(k) = pm_input(k)
769           p_dry(k) = p_input(k)
770           theta(k) = th_input(k)
771           rho(k) = rho_input(k)
772           u(k) = u_input(k)
773           v(k) = v_input(k)
774           qv(k) = qv_input(k)
775 
776         enddo
777 
778      if(debug) then
779       write(6,*) ' sounding '
780       write(6,*) '  k  height(m)  press (Pa) pd(Pa) theta (K) den(kg/m^3)  u(m/s)     v(m/s)    qv(g/g) '
781       do k=1,nl
782         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)
783       enddo
784 
785      end if
786 
787       end subroutine get_sounding
788 
789 !-------------------------------------------------------
790 
791       subroutine read_sounding( ps,ts,qvs,h,th,qv,u,v,n,nl,debug )
792       implicit none
793       integer n,nl
794       real ps,ts,qvs,h(n),th(n),qv(n),u(n),v(n)
795       logical end_of_file
796       logical debug
797 
798       integer k
799 
800       open(unit=10,file='input_sounding',form='formatted',status='old')
801       rewind(10)
802       read(10,*) ps, ts, qvs
803       if(debug) then
804         write(6,*) ' input sounding surface parameters '
805         write(6,*) ' surface pressure (mb) ',ps
806         write(6,*) ' surface pot. temp (K) ',ts
807         write(6,*) ' surface mixing ratio (g/kg) ',qvs
808       end if
809 
810       end_of_file = .false.
811       k = 0
812 
813       do while (.not. end_of_file)
814 
815         read(10,*,end=100) h(k+1), th(k+1), qv(k+1), u(k+1), v(k+1)
816         k = k+1
817         if(debug) write(6,'(1x,i3,5(1x,e10.3))') k, h(k), th(k), qv(k), u(k), v(k)
818         go to 110
819  100    end_of_file = .true.
820  110    continue
821       enddo
822 
823       nl = k
824 
825       close(unit=10,status = 'keep')
826 
827       end subroutine read_sounding
828 
829 END MODULE module_initialize_ideal