module_initialize_b_wave.F

References to this file elsewhere.
1 !IDEAL:MODEL_LAYER:INITIALIZATION
2 
3 !  This MODULE holds the routines which are used to perform various initializations
4 !  for the individual domains.  
5 
6 !-----------------------------------------------------------------------
7 
8 MODULE module_initialize_ideal
9 
10    USE module_domain
11    USE module_io_domain
12    USE module_state_description
13    USE module_model_constants
14    USE module_bc
15    USE module_timing
16    USE module_configure
17    USE module_init_utilities
18 #ifdef DM_PARALLEL
19    USE module_dm
20 #endif
21 
22 
23 CONTAINS
24 
25 
26 !-------------------------------------------------------------------
27 ! this is a wrapper for the solver-specific init_domain routines.
28 ! Also dereferences the grid variables and passes them down as arguments.
29 ! This is crucial, since the lower level routines may do message passing
30 ! and this will get fouled up on machines that insist on passing down
31 ! copies of assumed-shape arrays (by passing down as arguments, the 
32 ! data are treated as assumed-size -- ie. f77 -- arrays and the copying
33 ! business is avoided).  Fie on the F90 designers.  Fie and a pox.
34 
35    SUBROUTINE init_domain ( grid )
36 
37    IMPLICIT NONE
38 
39    !  Input data.
40    TYPE (domain), POINTER :: grid 
41    !  Local data.
42    INTEGER                :: dyn_opt 
43    INTEGER :: idum1, idum2
44 
45    CALL nl_get_dyn_opt( 1, dyn_opt )
46    
47    CALL set_scalar_indices_from_config ( head_grid%id , idum1, idum2 )
48 
49    IF (      dyn_opt .eq. 1 &
50         .or. dyn_opt .eq. 2 &
51         .or. dyn_opt .eq. 3 &
52                                        ) THEN
53      CALL init_domain_rk( grid &
54 !
55 #include <em_actual_new_args.inc>
56 !
57                         )
58 
59    ELSE
60      WRITE(0,*)' init_domain: unknown or unimplemented dyn_opt = ',dyn_opt
61      call wrf_error_fatal ( ' init_domain: unknown or unimplemented dyn_opt ' )
62    ENDIF
63 
64    END SUBROUTINE init_domain
65 
66 !-------------------------------------------------------------------
67 
68    SUBROUTINE init_domain_rk ( grid &
69 !
70 # include <em_dummy_new_args.inc>
71 !
72 )
73    IMPLICIT NONE
74 
75    !  Input data.
76    TYPE (domain), POINTER :: grid
77 
78 # include <em_dummy_decl.inc>
79 
80    TYPE (grid_config_rec_type)              :: config_flags
81 
82    !  Local data
83    INTEGER                             ::                       &
84                                   ids, ide, jds, jde, kds, kde, &
85                                   ims, ime, jms, jme, kms, kme, &
86                                   its, ite, jts, jte, kts, kte, &
87                                   i, j, k
88 
89    ! Local data
90 
91    INTEGER, PARAMETER :: nl_max = 1000
92    REAL, DIMENSION(nl_max) :: zk, p_in, theta, rho, u, v, qv, pd_in
93    INTEGER :: nl_in
94 
95    INTEGER :: icm,jcm, ii, im1, jj, jm1, loop, error, fid, nxc, nyc
96    REAL    :: u_mean,v_mean, f0, p_surf, p_level, qvf, z_at_v, z_at_u
97    REAL    :: z_scale, xrad, yrad, zrad, rad, delt, cof1, cof2
98 !   REAL, EXTERNAL :: interp_0
99    REAL    :: hm
100    REAL    :: pi
101 
102 !  stuff from original initialization that has been dropped from the Registry 
103    REAL    :: vnu, xnu, xnus, dinit0, cbh, p0_temp, t0_temp, zd, zt
104    REAL    :: qvf1, qvf2, pd_surf
105    INTEGER :: it
106 
107    LOGICAL :: moisture_init
108    LOGICAL :: stretch_grid, dry_sounding, debug
109 
110 !  kludge space for initial jet
111 
112    INTEGER, parameter :: nz_jet=64, ny_jet=80
113    REAL, DIMENSION(nz_jet, ny_jet) :: u_jet, rho_jet, th_jet, z_jet
114 
115 !  perturbation parameters
116 
117    REAL, PARAMETER :: htbub=8000., radbub=2000000., radz=8000., tpbub=1.0
118    REAL :: piov2, tp
119    INTEGER :: icen, jcen
120    real :: thtmp, ptmp, temp(3)
121 
122 #ifdef DM_PARALLEL
123 #    include <em_data_calls.inc>
124 #endif
125 
126    SELECT CASE ( model_data_order )
127          CASE ( DATA_ORDER_ZXY )
128    kds = grid%sd31 ; kde = grid%ed31 ;
129    ids = grid%sd32 ; ide = grid%ed32 ;
130    jds = grid%sd33 ; jde = grid%ed33 ;
131 
132    kms = grid%sm31 ; kme = grid%em31 ;
133    ims = grid%sm32 ; ime = grid%em32 ;
134    jms = grid%sm33 ; jme = grid%em33 ;
135 
136    kts = grid%sp31 ; kte = grid%ep31 ;   ! note that tile is entire patch
137    its = grid%sp32 ; ite = grid%ep32 ;   ! note that tile is entire patch
138    jts = grid%sp33 ; jte = grid%ep33 ;   ! note that tile is entire patch
139          CASE ( DATA_ORDER_XYZ )
140    ids = grid%sd31 ; ide = grid%ed31 ;
141    jds = grid%sd32 ; jde = grid%ed32 ;
142    kds = grid%sd33 ; kde = grid%ed33 ;
143 
144    ims = grid%sm31 ; ime = grid%em31 ;
145    jms = grid%sm32 ; jme = grid%em32 ;
146    kms = grid%sm33 ; kme = grid%em33 ;
147 
148    its = grid%sp31 ; ite = grid%ep31 ;   ! note that tile is entire patch
149    jts = grid%sp32 ; jte = grid%ep32 ;   ! note that tile is entire patch
150    kts = grid%sp33 ; kte = grid%ep33 ;   ! note that tile is entire patch
151          CASE ( DATA_ORDER_XZY )
152    ids = grid%sd31 ; ide = grid%ed31 ;
153    kds = grid%sd32 ; kde = grid%ed32 ;
154    jds = grid%sd33 ; jde = grid%ed33 ;
155 
156    ims = grid%sm31 ; ime = grid%em31 ;
157    kms = grid%sm32 ; kme = grid%em32 ;
158    jms = grid%sm33 ; jme = grid%em33 ;
159 
160    its = grid%sp31 ; ite = grid%ep31 ;   ! note that tile is entire patch
161    kts = grid%sp32 ; kte = grid%ep32 ;   ! note that tile is entire patch
162    jts = grid%sp33 ; jte = grid%ep33 ;   ! note that tile is entire patch
163 
164    END SELECT
165 
166    piov2 = 2.*atan(1.0)
167    icen = ide/4
168    jcen = jde/2
169 
170    stretch_grid = .true.
171    delt = 0.
172    z_scale = .50
173    pi = 2.*asin(1.0)
174    write(6,*) ' pi is ',pi
175    nxc = (ide-ids)/4
176    nyc = (jde-jds)/2
177 
178    CALL model_to_grid_config_rec ( grid%id , model_config_rec , config_flags )
179 
180 ! here we check to see if the boundary conditions are set properly
181 
182    CALL boundary_condition_check( config_flags, bdyzone, error, grid%id )
183 
184    moisture_init = .true.
185 
186     grid%itimestep=0
187 
188 #ifdef DM_PARALLEL
189    CALL wrf_dm_bcast_bytes( icm , IWORDSIZE )
190    CALL wrf_dm_bcast_bytes( jcm , IWORDSIZE )
191 #endif
192 
193     CALL nl_set_mminlu(1,'    ')
194     CALL nl_set_iswater(1,0)
195     CALL nl_set_cen_lat(1,40.)
196     CALL nl_set_cen_lon(1,-105.)
197     CALL nl_set_truelat1(1,0.)
198     CALL nl_set_truelat2(1,0.)
199     CALL nl_set_moad_cen_lat (1,0.)
200     CALL nl_set_stand_lon (1,0.)
201     CALL nl_set_map_proj(1,0)
202 
203 
204 !  here we initialize data we currently is not initialized 
205 !  in the input data
206 
207     DO j = jts, jte
208       DO i = its, ite
209 
210          grid%ht(i,j)       = 0.
211          grid%msftx(i,j)    = 1.
212          grid%msfty(i,j)    = 1.
213          grid%msfux(i,j)    = 1.
214          grid%msfuy(i,j)    = 1.
215          grid%msfvx(i,j)    = 1.
216          grid%msfvx_inv(i,j)= 1.
217          grid%msfvy(i,j)    = 1.
218          grid%sina(i,j)     = 0.
219          grid%cosa(i,j)     = 1.
220          grid%e(i,j)        = 0.
221          grid%f(i,j)        = 1.e-04
222 
223       END DO
224    END DO
225 
226     DO j = jts, jte
227     DO k = kts, kte
228       DO i = its, ite
229          grid%em_ww(i,k,j)     = 0.
230       END DO
231    END DO
232    END DO
233 
234    grid%step_number = 0
235 
236 ! set up the grid
237 
238    IF (stretch_grid) THEN ! exponential stretch for eta (nearly constant dz)
239      DO k=1, kde
240       grid%em_znw(k) = (exp(-(k-1)/float(kde-1)/z_scale) - exp(-1./z_scale))/ &
241                                 (1.-exp(-1./z_scale))
242      ENDDO
243    ELSE
244      DO k=1, kde
245       grid%em_znw(k) = 1. - float(k-1)/float(kde-1)
246      ENDDO
247    ENDIF
248 
249    DO k=1, kde-1
250     grid%em_dnw(k) = grid%em_znw(k+1) - grid%em_znw(k)
251     grid%em_rdnw(k) = 1./grid%em_dnw(k)
252     grid%em_znu(k) = 0.5*(grid%em_znw(k+1)+grid%em_znw(k))
253    ENDDO
254    DO k=2, kde-1
255     grid%em_dn(k) = 0.5*(grid%em_dnw(k)+grid%em_dnw(k-1))
256     grid%em_rdn(k) = 1./grid%em_dn(k)
257     grid%em_fnp(k) = .5* grid%em_dnw(k  )/grid%em_dn(k)
258     grid%em_fnm(k) = .5* grid%em_dnw(k-1)/grid%em_dn(k)
259    ENDDO
260 
261    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) 
262    cof2 =     grid%em_dn(2)        /(grid%em_dn(2)+grid%em_dn(3))*grid%em_dnw(1)/grid%em_dn(3) 
263    grid%cf1  = grid%em_fnp(2) + cof1
264    grid%cf2  = grid%em_fnm(2) - cof1 - cof2
265    grid%cf3  = cof2       
266 
267    grid%cfn  = (.5*grid%em_dnw(kde-1)+grid%em_dn(kde-1))/grid%em_dn(kde-1)
268    grid%cfn1 = -.5*grid%em_dnw(kde-1)/grid%em_dn(kde-1)
269    grid%rdx = 1./config_flags%dx
270    grid%rdy = 1./config_flags%dy
271 
272 !  get the sounding from the ascii sounding file, first get dry sounding and 
273 !  calculate base state
274 
275   write(6,*) ' reading input jet sounding '
276   call read_input_jet( u_jet, rho_jet, th_jet, z_jet, nz_jet, ny_jet )
277 
278   write(6,*) ' getting dry sounding for base state '
279   write(6,*) ' using middle column in jet sounding, j = ',ny_jet/2
280   dry_sounding = .true.
281 
282   dry_sounding   = .true.
283   debug = .true.  !  this will produce print of the sounding
284   CALL get_sounding( zk, p_in, pd_in, theta, rho, u, v, qv, dry_sounding, &
285                       nl_max, nl_in, u_jet, rho_jet, th_jet, z_jet,      &
286                       nz_jet, ny_jet, ny_jet/2, debug                   )
287 
288   write(6,*) ' returned from reading sounding, nl_in is ',nl_in
289 
290 !  find ptop for the desired ztop (ztop is input from the namelist),
291 !  and find surface pressure
292 
293 !  For the jet, using the middle column for the base state means that
294 !  we will be extrapolating above the highest height data to the south
295 !  of the centerline.
296 
297   grid%p_top = interp_0( p_in, zk, config_flags%ztop, nl_in )
298 
299   DO j=jts,jte
300   DO i=its,ite  ! flat surface
301     grid%em_phb(i,1,j) = 0.
302     grid%em_php(i,1,j) = 0.
303     grid%em_ph0(i,1,j) = 0.
304     grid%ht(i,j) = 0.
305   ENDDO
306   ENDDO
307 
308   DO J = jts, jte
309   DO I = its, ite
310 
311     p_surf = interp_0( p_in, zk, grid%em_phb(i,1,j)/g, nl_in )
312     grid%em_mub(i,j) = p_surf-grid%p_top
313 
314 !  this is dry hydrostatic sounding (base state), so given grid%em_p (coordinate),
315 !  interp theta (from interp) and compute 1/rho from eqn. of state
316 
317     DO K = 1, kte-1
318       p_level = grid%em_znu(k)*(p_surf - grid%p_top) + grid%p_top
319       grid%em_pb(i,k,j) = p_level
320       grid%em_t_init(i,k,j) = interp_0( theta, p_in, p_level, nl_in ) - t0
321       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
322     ENDDO
323 
324 !  calc hydrostatic balance (alternatively we could interp the geopotential from the
325 !  sounding, but this assures that the base state is in exact hydrostatic balance with
326 !  respect to the model eqns.
327 
328     DO k  = 2,kte
329       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)
330     ENDDO
331 
332   ENDDO
333   ENDDO
334 
335   write(6,*) ' ptop is ',grid%p_top
336   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
337 
338 !  calculate full state for each column - this includes moisture.
339 
340   write(6,*) ' getting grid%moist sounding for full state '
341 
342   dry_sounding = .true.
343   IF (config_flags%mp_physics /= 0)  dry_sounding = .false.
344 
345   DO J = jts, min(jde-1,jte)
346 
347 !  get sounding for this point
348 
349   debug = .false.  !  this will turn off print of the sounding
350   CALL get_sounding( zk, p_in, pd_in, theta, rho, u, v, qv, dry_sounding, &
351                       nl_max, nl_in, u_jet, rho_jet, th_jet, z_jet,      &
352                       nz_jet, ny_jet, j, debug                          )
353 
354   DO I = its, min(ide-1,ite)
355 
356 !   we could just do the first point in "i" and copy from there, but we'll
357 !   be lazy and do all the points as if they are all, independent
358 
359 !   At this point grid%p_top is already set. find the DRY mass in the column 
360 !   by interpolating the DRY pressure.  
361 
362     pd_surf = interp_0( pd_in, zk, grid%em_phb(i,1,j)/g, nl_in )
363 
364 !   compute the perturbation mass and the full mass
365 
366     grid%em_mu_1(i,j) = pd_surf-grid%p_top - grid%em_mub(i,j)
367     grid%em_mu_2(i,j) = grid%em_mu_1(i,j)
368     grid%em_mu0(i,j) = grid%em_mu_1(i,j) + grid%em_mub(i,j)
369 
370 !   given the dry pressure and coordinate system, interp the potential
371 !   temperature and qv
372 
373     do k=1,kde-1
374 
375       p_level = grid%em_znu(k)*(pd_surf - grid%p_top) + grid%p_top
376 
377       grid%moist(i,k,j,P_QV) = interp_0( qv, pd_in, p_level, nl_in )
378       grid%em_t_1(i,k,j)          = interp_0( theta, pd_in, p_level, nl_in ) - t0
379       grid%em_t_2(i,k,j)          = grid%em_t_1(i,k,j)
380       
381 
382     enddo
383 
384 !   integrate the hydrostatic equation (from the RHS of the bigstep
385 !   vertical momentum equation) down from the top to get grid%em_p.
386 !   first from the top of the model to the top pressure
387 
388     k = kte-1  ! top level
389 
390     qvf1 = 0.5*(grid%moist(i,k,j,P_QV)+grid%moist(i,k,j,P_QV))
391     qvf2 = 1./(1.+qvf1)
392     qvf1 = qvf1*qvf2
393 
394 !    grid%em_p(i,k,j) = - 0.5*grid%em_mu_1(i,j)/grid%em_rdnw(k)
395     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
396     qvf = 1. + rvovrd*grid%moist(i,k,j,P_QV)
397     grid%em_alt(i,k,j) = (r_d/p1000mb)*(grid%em_t_1(i,k,j)+t0)*qvf* &
398                 (((grid%em_p(i,k,j)+grid%em_pb(i,k,j))/p1000mb)**cvpm)
399     grid%em_al(i,k,j) = grid%em_alt(i,k,j) - grid%em_alb(i,k,j)
400 
401 !  down the column
402 
403     do k=kte-2,1,-1
404       qvf1 = 0.5*(grid%moist(i,k,j,P_QV)+grid%moist(i,k+1,j,P_QV))
405       qvf2 = 1./(1.+qvf1)
406       qvf1 = qvf1*qvf2
407       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)
408       qvf = 1. + rvovrd*grid%moist(i,k,j,P_QV)
409       grid%em_alt(i,k,j) = (r_d/p1000mb)*(grid%em_t_1(i,k,j)+t0)*qvf* &
410                   (((grid%em_p(i,k,j)+grid%em_pb(i,k,j))/p1000mb)**cvpm)
411       grid%em_al(i,k,j) = grid%em_alt(i,k,j) - grid%em_alb(i,k,j)
412     enddo
413 
414 !  this is the hydrostatic equation used in the model after the
415 !  small timesteps.  In the model, grid%em_al (inverse density)
416 !  is computed from the geopotential.
417 
418 
419     grid%em_ph_1(i,1,j) = 0.
420     DO k  = 2,kte
421       grid%em_ph_1(i,k,j) = grid%em_ph_1(i,k-1,j) - (1./grid%em_rdnw(k-1))*(       &
422                    (grid%em_mub(i,j)+grid%em_mu_1(i,j))*grid%em_al(i,k-1,j)+ &
423                     grid%em_mu_1(i,j)*grid%em_alb(i,k-1,j)  )
424                                                    
425       grid%em_ph_2(i,k,j) = grid%em_ph_1(i,k,j) 
426       grid%em_ph0(i,k,j) = grid%em_ph_1(i,k,j) + grid%em_phb(i,k,j)
427     ENDDO
428 
429 ! interp u
430 
431     DO K = 1, kte
432       p_level = grid%em_znu(k)*(p_surf - grid%p_top) + grid%p_top
433       grid%em_u_1(i,k,j) = interp_0( u, p_in, p_level, nl_in )
434       grid%em_u_2(i,k,j) = grid%em_u_1(i,k,j)
435     ENDDO
436 
437   ENDDO
438   ENDDO
439 
440 !  thermal perturbation to kick off convection
441 
442   write(6,*) ' nxc, nyc for perturbation ',nxc,nyc
443   write(6,*) ' delt for perturbation ',tpbub
444 
445   DO J = jts, min(jde-1,jte)
446     yrad = config_flags%dy*float(j-jde/2-1)/radbub
447     DO I = its, min(ide-1,ite)
448       xrad = float(i-1)/float(ide-ids)
449 
450       DO K = 1, kte-1
451 
452 !  put in preturbation theta (bubble) and recalc density.  note,
453 !  the mass in the column is not changing, so when theta changes,
454 !  we recompute density and geopotential
455 
456         zrad = 0.5*(grid%em_ph_1(i,k,j)+grid%em_ph_1(i,k+1,j)  &
457                    +grid%em_phb(i,k,j)+grid%em_phb(i,k+1,j))/g
458         zrad = (zrad-htbub)/radz
459         RAD=SQRT(yrad*yrad+zrad*zrad)
460         IF(RAD <= 1.) THEN
461            tp = tpbub*cos(rad*piov2)*cos(rad*piov2)*cos(xrad*2*pi+pi)
462            grid%em_t_1(i,k,j)=grid%em_t_1(i,k,j)+tp
463            grid%em_t_2(i,k,j)=grid%em_t_1(i,k,j)
464            qvf = 1. + rvovrd*grid%moist(i,k,j,P_QV)
465            grid%em_alt(i,k,j) = (r_d/p1000mb)*(grid%em_t_1(i,k,j)+t0)*qvf* &
466                         (((grid%em_p(i,k,j)+grid%em_pb(i,k,j))/p1000mb)**cvpm)
467            grid%em_al(i,k,j) = grid%em_alt(i,k,j) - grid%em_alb(i,k,j)
468         ENDIF
469       ENDDO
470 
471 !  rebalance hydrostatically
472 
473       DO k  = 2,kte
474         grid%em_ph_1(i,k,j) = grid%em_ph_1(i,k-1,j) - (1./grid%em_rdnw(k-1))*(       &
475                      (grid%em_mub(i,j)+grid%em_mu_1(i,j))*grid%em_al(i,k-1,j)+ &
476                       grid%em_mu_1(i,j)*grid%em_alb(i,k-1,j)  )
477                                                    
478         grid%em_ph_2(i,k,j) = grid%em_ph_1(i,k,j) 
479         grid%em_ph0(i,k,j) = grid%em_ph_1(i,k,j) + grid%em_phb(i,k,j)
480       ENDDO
481 
482     ENDDO
483   ENDDO
484 
485 !#endif
486 
487    write(6,*) ' grid%em_mu_1 from comp ', grid%em_mu_1(1,1)
488    write(6,*) ' pert state sounding from comp, grid%em_ph_1, pp, grid%em_al, grid%em_t_1, qv '
489    do k=1,kde-1
490      write(6,'(i3,1x,5(1x,1pe10.3))') k, grid%em_ph_1(1,k,1),grid%em_p(1,k,1), grid%em_al(1,k,1), &
491                      grid%em_t_1(1,k,1), grid%moist(1,k,1,P_QV)
492    enddo
493 
494    write(6,*) ' grid%em_mu_1 from comp ', grid%em_mu_1(1,1)
495    write(6,*) ' full state sounding from comp, ph, grid%em_p, grid%em_al, grid%em_t_1, qv '
496    write(6,*) ' at j = 1 '
497    do k=1,kde-1
498      write(6,'(i3,1x,5(1x,1pe10.3))') k, grid%em_ph_1(1,k,1)+grid%em_phb(1,k,1), &
499                      grid%em_p(1,k,1)+grid%em_pb(1,k,1), grid%em_al(1,k,1)+grid%em_alb(1,k,1), &
500                      grid%em_t_1(1,k,1)+t0, grid%moist(1,k,1,P_QV)
501    enddo
502 
503 
504    write(6,*) ' full state sounding from comp, ph, grid%em_p, grid%em_al, grid%em_t_1, qv '
505    write(6,*) ' at j = jde/2 '
506    do k=1,kde-1
507      write(6,'(i3,1x,5(1x,1pe10.3))') k, grid%em_ph_1(1,k,jde/2)+grid%em_phb(1,k,jde/2), &
508                      grid%em_p(1,k,jde/2)+grid%em_pb(1,k,jde/2), grid%em_al(1,k,jde/2)+grid%em_alb(1,k,jde/2), &
509                      grid%em_t_1(1,k,jde/2)+t0, grid%moist(1,k,jde/2,P_QV)
510    enddo
511 
512    write(6,*) ' full state sounding from comp, ph, grid%em_p, grid%em_al, grid%em_t_1, qv '
513    write(6,*) ' at j = jde-1 '
514    do k=1,kde-1
515      write(6,'(i3,1x,5(1x,1pe10.3))') k, grid%em_ph_1(1,k,jde-1)+grid%em_phb(1,k,jde-1), &
516                      grid%em_p(1,k,jde-1)+grid%em_pb(1,k,jde-1), grid%em_al(1,k,jde-1)+grid%em_alb(1,k,jde-1), &
517                      grid%em_t_1(1,k,jde-1)+t0, grid%moist(1,k,jde-1,P_QV)
518    enddo
519 
520 ! set v
521 
522   DO J = jts, jte
523   DO I = its, min(ide-1,ite)
524 
525     DO K = 1, kte
526       grid%em_v_1(i,k,j) = 0.
527       grid%em_v_2(i,k,j) = grid%em_v_1(i,k,j)
528     ENDDO
529 
530   ENDDO
531   ENDDO
532 
533 !  fill out last i row for u
534 
535   DO J = jts, min(jde-1,jte)
536   DO I = ite, ite
537 
538     DO K = 1, kte
539       grid%em_u_1(i,k,j) = grid%em_u_1(its,k,j)
540       grid%em_u_2(i,k,j) = grid%em_u_2(its,k,j)
541     ENDDO
542 
543   ENDDO
544   ENDDO
545 
546 !  set w
547 
548   DO J = jts, min(jde-1,jte)
549   DO K = kts, kte
550   DO I = its, min(ide-1,ite)
551     grid%em_w_1(i,k,j) = 0.
552     grid%em_w_2(i,k,j) = 0.
553   ENDDO
554   ENDDO
555   ENDDO
556 
557 !  set a few more things
558 
559   DO J = jts, min(jde-1,jte)
560   DO K = kts, kte-1
561   DO I = its, min(ide-1,ite)
562     grid%h_diabatic(i,k,j) = 0.
563   ENDDO
564   ENDDO
565   ENDDO
566 
567   DO k=1,kte-1
568     grid%em_t_base(k) = grid%em_t_1(1,k,1)
569     grid%qv_base(k) = grid%moist(1,k,1,P_QV)
570     grid%u_base(k) = grid%em_u_1(1,k,1)
571     grid%v_base(k) = grid%em_v_1(1,k,1)
572   ENDDO
573 
574   DO J = jts, min(jde-1,jte)
575   DO I = its, min(ide-1,ite)
576      thtmp   = grid%em_t_2(i,1,j)+t0
577      ptmp    = grid%em_p(i,1,j)+grid%em_pb(i,1,j)
578      temp(1) = thtmp * (ptmp/p1000mb)**rcp
579      thtmp   = grid%em_t_2(i,2,j)+t0
580      ptmp    = grid%em_p(i,2,j)+grid%em_pb(i,2,j)
581      temp(2) = thtmp * (ptmp/p1000mb)**rcp
582      thtmp   = grid%em_t_2(i,3,j)+t0
583      ptmp    = grid%em_p(i,3,j)+grid%em_pb(i,3,j)
584      temp(3) = thtmp * (ptmp/p1000mb)**rcp
585 
586      grid%tsk(I,J)=grid%cf1*temp(1)+grid%cf2*temp(2)+grid%cf3*temp(3)
587      if (i .eq. 1) print*,'sfctem',j,temp(1),temp(2),temp(3),grid%tsk(I,J)
588      grid%tmn(I,J)=grid%tsk(I,J)-0.5
589   ENDDO
590   ENDDO
591 
592   RETURN
593 
594  END SUBROUTINE init_domain_rk
595 
596 !---------------------------------------------------------------------
597 
598  SUBROUTINE init_module_initialize
599  END SUBROUTINE init_module_initialize
600 
601 !---------------------------------------------------------------------
602 #if 0
603 ! TEST DRIVER FOR "read_input_jet" and "get_sounding"
604   implicit none 
605   integer, parameter :: nz_jet=64, ny_jet=80
606   real, dimension(nz_jet,ny_jet) :: u_jet, rho_jet, &
607                                     th_jet, z_jet
608 
609   real, dimension(nz_jet,ny_jet) :: zk,p,p_dry,theta,rho,u,v,qv
610   logical :: dry, debug
611   integer :: j, nl
612 
613   call read_input_jet( u_jet, rho_jet, th_jet, z_jet, nz_jet, ny_jet )
614 
615   call opngks
616   call parray( u_jet, nz_jet, ny_jet)
617   call parray( rho_jet, nz_jet, ny_jet)
618   call parray( th_jet, nz_jet, ny_jet)
619 !  call clsgks
620 
621 !  set up initial jet
622 
623   debug = .true.
624   dry = .true.
625   do j=1,ny_jet
626 
627     call get_sounding( zk(:,j),p(:,j),p_dry(:,j),theta(:,j),      &
628                        rho(:,j),u(:,j), v(:,j),  qv(:,j),        &
629                        dry, nz_jet, nl, u_jet, rho_jet, th_jet,  &
630                        z_jet, nz_jet, ny_jet, j, debug          )
631     debug = .false.
632 
633   enddo
634 
635   write(6,*) ' lowest level p, th, and rho, highest level p '
636 
637   do j=1,ny_jet
638     write(6,*) j, p(1,j),theta(1,j),rho(1,j), p(nz_jet,j)
639 !    write(6,*) j, p(1,j),theta(1,j)-th_jet(1,j),rho(1,j)-rho_jet(1,j)
640   enddo
641 
642   call parray( p, nz_jet, ny_jet)
643   call parray( p_dry, nz_jet, ny_jet)
644   call parray( theta, nz_jet, ny_jet)
645 
646   call clsgks
647 
648   end
649 
650 !---------------------------------
651 
652       subroutine parray(a,m,n)
653       dimension a(m,n)
654       dimension b(n,m)
655 
656     do i=1,m
657     do j=1,n
658       b(j,i) = a(i,j)
659     enddo
660     enddo
661       
662       write(6,'(''  dimensions m,n  '',2i6)')m,n
663         call set(.05,.95,.05,.95,0.,1.,0.,1.,1)
664         call perim(4,5,4,5)
665         call setusv('LW',2000)
666 !        CALL CONREC(a,m,m,n,cmax,cmin,cinc,-1,-638,-922)
667         CALL CONREC(b,n,n,m,0.,0.,0.,-1,-638,-922)
668         call frame
669       return
670       end
671 
672 ! END TEST DRIVER FOR "read_input_jet" and "get_sounding"
673 #endif
674 
675 !------------------------------------------------------------------
676 
677     subroutine get_sounding( zk, p, p_dry, theta, rho,       &
678                              u, v, qv, dry, nl_max, nl_in,  &
679                              u_jet, rho_jet, th_jet, z_jet, &
680                              nz_jet, ny_jet, j_point, debug )
681     implicit none
682 
683     integer nl_max, nl_in
684     real zk(nl_max), p(nl_max), theta(nl_max), rho(nl_max), &
685          u(nl_max), v(nl_max), qv(nl_max), p_dry(nl_max)
686     logical dry
687 
688     integer nz_jet, ny_jet, j_point
689     real, dimension(nz_jet, ny_jet) :: u_jet, rho_jet, th_jet, z_jet
690 
691     integer n
692     parameter(n=1000)
693     logical debug
694 
695 ! input sounding data
696 
697     real p_surf, th_surf, qv_surf
698     real pi_surf, pi(n)
699     real h_input(n), th_input(n), qv_input(n), u_input(n), v_input(n)
700 
701 ! diagnostics
702 
703     real rho_surf, p_input(n), rho_input(n)
704     real pm_input(n)  !  this are for full moist sounding
705 
706 ! local data
707 
708     real p1000mb,cv,cp,r,cvpm,g
709     parameter (p1000mb = 1.e+05, r = 287, cp = 1003., cv = cp-r, cvpm = -cv/cp, g=9.81 )
710     integer k, it, nl
711     real qvf, qvf1, dz
712 
713 !  first, read the sounding
714 
715 !    call read_sounding( p_surf, th_surf, qv_surf, &
716 !                          h_input, th_input, qv_input, u_input, v_input,n, nl, debug )
717 
718    call calc_jet_sounding( p_surf, th_surf, qv_surf,                             &
719                            h_input, th_input, qv_input, u_input, v_input,        &
720                            n, nl, debug, u_jet, rho_jet, th_jet, z_jet, j_point, &
721                            nz_jet, ny_jet, dry                                  )
722 
723    nl = nz_jet
724 
725     if(dry) then
726      do k=1,nl
727        qv_input(k) = 0.
728      enddo
729     endif
730 
731     if(debug) write(6,*) ' number of input levels = ',nl
732 
733       nl_in = nl
734       if(nl_in .gt. nl_max ) then
735         write(6,*) ' too many levels for input arrays ',nl_in,nl_max
736         call wrf_error_fatal ( ' too many levels for input arrays ' )
737       end if
738 
739 !  compute diagnostics,
740 !  first, convert qv(g/kg) to qv(g/g)
741 !
742 !      do k=1,nl
743 !        qv_input(k) = 0.001*qv_input(k)
744 !      enddo
745 !      p_surf = 100.*p_surf  ! convert to pascals
746 
747     qvf = 1. + rvovrd*qv_input(1) 
748     rho_surf = 1./((r/p1000mb)*th_surf*qvf*((p_surf/p1000mb)**cvpm))
749     pi_surf = (p_surf/p1000mb)**(r/cp)
750 
751     if(debug) then
752       write(6,*) ' surface density is ',rho_surf
753       write(6,*) ' surface pi is    ',pi_surf
754     end if
755 
756 
757 !  integrate moist sounding hydrostatically, starting from the
758 !  specified surface pressure
759 !  -> first, integrate from surface to lowest level
760 
761         qvf = 1. + rvovrd*qv_input(1) 
762         qvf1 = 1. + qv_input(1)
763         rho_input(1) = rho_surf
764         dz = h_input(1)
765         do it=1,10
766           pm_input(1) = p_surf &
767                   - 0.5*dz*(rho_surf+rho_input(1))*g*qvf1
768           rho_input(1) = 1./((r/p1000mb)*th_input(1)*qvf*((pm_input(1)/p1000mb)**cvpm))
769         enddo
770 
771 ! integrate up the column
772 
773         do k=2,nl
774           rho_input(k) = rho_input(k-1)
775           dz = h_input(k)-h_input(k-1)
776           qvf1 = 0.5*(2.+(qv_input(k-1)+qv_input(k)))
777           qvf = 1. + rvovrd*qv_input(k)   ! qv is in g/kg here
778  
779           do it=1,10
780             pm_input(k) = pm_input(k-1) &
781                     - 0.5*dz*(rho_input(k)+rho_input(k-1))*g*qvf1
782             rho_input(k) = 1./((r/p1000mb)*th_input(k)*qvf*((pm_input(k)/p1000mb)**cvpm))
783           enddo
784         enddo
785 
786 !  we have the moist sounding
787 
788 !  next, compute the dry sounding using p at the highest level from the
789 !  moist sounding and integrating down.
790 
791         p_input(nl) = pm_input(nl)
792 
793           do k=nl-1,1,-1
794             dz = h_input(k+1)-h_input(k)
795             p_input(k) = p_input(k+1) + 0.5*dz*(rho_input(k)+rho_input(k+1))*g
796           enddo
797 
798 
799         do k=1,nl
800 
801           zk(k) = h_input(k)
802           p(k) = pm_input(k)
803           p_dry(k) = p_input(k)
804           theta(k) = th_input(k)
805           rho(k) = rho_input(k)
806           u(k) = u_input(k)
807           v(k) = v_input(k)
808           qv(k) = qv_input(k)
809 
810         enddo
811 
812      if(debug) then
813       write(6,*) ' sounding '
814       write(6,*) '  k  height(m)  press (Pa)   pd(Pa)   theta (K)  den(kg/m^3)  u(m/s)     v(m/s)    qv(g/g) '
815       do k=1,nl
816         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)
817       enddo
818 
819      end if
820 
821      end subroutine get_sounding
822 
823 !------------------------------------------------------------------
824 
825   subroutine calc_jet_sounding( p_surf, th_surf, qv_surf,      &
826                                 h, th, qv, u, v, n, nl, debug, &
827                                 u_jet, rho_jet, th_jet, z_jet, &
828                                 jp, nz_jet, ny_jet, dry       )
829   implicit none
830   integer :: n, nl, jp, nz_jet, ny_jet
831 
832   real, dimension(nz_jet, ny_jet) :: u_jet, rho_jet, th_jet, z_jet
833   real, dimension(n) :: h,th,qv,u,v
834   real :: p_surf, th_surf, qv_surf
835   logical :: debug, dry
836 
837   real, dimension(1:nz_jet) :: rho, rel_hum, p
838   integer :: k
839 
840 !  some local stuff
841 
842   real :: tmppi, es, qvs, temperature
843   real, parameter :: p1000mb=1.e+05, rcp=287./1004.5, svpt0=273.15, &
844                      svp3 = 29.65, ep_2=287./461.6, r_d = 287., &
845                      cpovcv = 1004./(1004.-287.),               &
846                      svp1 = 0.6112, svp2 = 17.67
847 
848 !  get sounding from column jp
849 
850    do k=1,nz_jet
851      h(k)  = z_jet(k,jp)
852      th(k) = th_jet(k,jp)
853      qv(k) = 0.
854      rho(k) = rho_jet(k,jp)
855      u(k) = u_jet(k,jp)
856      v(k) = 0.
857    enddo
858 
859    if (.not.dry) then
860      DO k=1,nz_jet
861        if(h(k) .gt. 8000.) then
862          rel_hum(k)=0.1
863        else
864          rel_hum(k)=(1.-0.90*(h(k)/8000.)**1.25)
865        end if
866        rel_hum(k) = min(0.7,rel_hum(k))
867      ENDDO
868    else
869      do k=1,nz_jet
870        rel_hum(k) = 0.
871      enddo
872    endif
873 
874 !  next, compute pressure
875 
876    do k=1,nz_jet
877      p(k) = p1000mb*(R_d*rho(k)*th(k)/p1000mb)**cpovcv
878    enddo
879 
880 !  here we adjust for fixed moisture profile
881 
882      IF (.not.dry)  THEN
883 
884 !  here we assume the input theta is th_v, so we reset theta accordingly
885 
886        DO k=1,nz_jet
887          tmppi=(p(k)/p1000mb)**rcp
888          temperature = tmppi*th(k)
889          if (temperature .gt. svpt0) then
890             es  = 1000.*svp1*exp(svp2*(temperature-svpt0)/(temperature-svp3))
891             qvs = ep_2*es/(p(k)-es)
892          else
893             es  = 1000.*svp1*exp( 21.8745584*(temperature-273.16)/(temperature-7.66) )
894             qvs = ep_2*es/(p(k)-es)
895          endif
896          qv(k) = rel_hum(k)*qvs
897          th(k) = th(k)/(1.+.61*qv(k))
898        ENDDO
899  
900      ENDIF
901 
902 !  finally, set the surface data. We'll just do a simple extrapolation
903 
904    p_surf = 1.5*p(1) - 0.5*p(2)
905    th_surf = 1.5*th(1) - 0.5*th(2)
906    qv_surf = 1.5*qv(1) - 0.5*qv(2)
907 
908    end subroutine calc_jet_sounding
909 
910 !---------------------------------------------------------------------
911 
912  SUBROUTINE read_input_jet( u, r, t, zk, nz, ny )
913  implicit none
914 
915  integer, intent(in) :: nz,ny
916  real, dimension(nz,ny), intent(out) :: u,r,t,zk
917  integer :: ny_in, nz_in, j,k
918  real, dimension(ny,nz) :: field_in
919 
920 ! this code assumes it is called on processor 0 only
921 
922    OPEN(unit=10, file='input_jet', form='unformatted', status='old' )
923    REWIND(10) 
924    read(10) ny_in,nz_in
925    if((ny_in /= ny ) .or. (nz_in /= nz)) then
926      write(0,*) ' error in input jet dimensions '
927      write(0,*) ' ny, ny_input, nz, nz_input ', ny, ny_in, nz,nz_in
928      write(0,*) ' error exit '
929      call wrf_error_fatal ( ' error in input jet dimensions ' )
930    end if
931    read(10) field_in
932    do j=1,ny
933    do k=1,nz
934      u(k,j) = field_in(j,k)
935    enddo
936    enddo
937    read(10) field_in
938    do j=1,ny
939    do k=1,nz
940      t(k,j) = field_in(j,k)
941    enddo
942    enddo
943 
944    read(10) field_in
945    do j=1,ny
946    do k=1,nz
947      r(k,j) = field_in(j,k)
948    enddo
949    enddo
950 
951    do j=1,ny
952    do k=1,nz
953      zk(k,j) = 125. + 250.*float(k-1)
954    enddo
955    enddo
956 
957  end subroutine read_input_jet
958 
959 END MODULE module_initialize_ideal