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