merge standard release WRF/WPS V3.0.1.1 into wrffire
[wrffire.git] / wrfv2_fire / dyn_em / module_initialize_b_wave.F
blob70946b0c7be2f3f9cf4a1ad46b69bf3477dc239b
1 !IDEAL:MODEL_LAYER:INITIALIZATION
3 !  This MODULE holds the routines which are used to perform various initializations
4 !  for the individual domains.  
6 !-----------------------------------------------------------------------
8 MODULE module_initialize_ideal
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
23 CONTAINS
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.
35    SUBROUTINE init_domain ( grid )
37    IMPLICIT NONE
39    !  Input data.
40    TYPE (domain), POINTER :: grid 
41    !  Local data.
42    INTEGER :: idum1, idum2
44    CALL set_scalar_indices_from_config ( head_grid%id , idum1, idum2 )
46      CALL init_domain_rk( grid &
48 #include <actual_new_args.inc>
50                         )
52    END SUBROUTINE init_domain
54 !-------------------------------------------------------------------
56    SUBROUTINE init_domain_rk ( grid &
58 # include <dummy_new_args.inc>
61    IMPLICIT NONE
63    !  Input data.
64    TYPE (domain), POINTER :: grid
66 # include <dummy_decl.inc>
68    TYPE (grid_config_rec_type)              :: config_flags
70    !  Local data
71    INTEGER                             ::                       &
72                                   ids, ide, jds, jde, kds, kde, &
73                                   ims, ime, jms, jme, kms, kme, &
74                                   its, ite, jts, jte, kts, kte, &
75                                   i, j, k
77    ! Local data
79    INTEGER, PARAMETER :: nl_max = 1000
80    REAL, DIMENSION(nl_max) :: zk, p_in, theta, rho, u, v, qv, pd_in
81    INTEGER :: nl_in
83    INTEGER :: icm,jcm, ii, im1, jj, jm1, loop, error, fid, nxc, nyc
84    REAL    :: u_mean,v_mean, f0, p_surf, p_level, qvf, z_at_v, z_at_u
85    REAL    :: z_scale, xrad, yrad, zrad, rad, delt, cof1, cof2
86 !   REAL, EXTERNAL :: interp_0
87    REAL    :: hm
88    REAL    :: pi
90 !  stuff from original initialization that has been dropped from the Registry 
91    REAL    :: vnu, xnu, xnus, dinit0, cbh, p0_temp, t0_temp, zd, zt
92    REAL    :: qvf1, qvf2, pd_surf
93    INTEGER :: it
95    LOGICAL :: moisture_init
96    LOGICAL :: stretch_grid, dry_sounding, debug
98 !  kludge space for initial jet
100    INTEGER, parameter :: nz_jet=64, ny_jet=80
101    REAL, DIMENSION(nz_jet, ny_jet) :: u_jet, rho_jet, th_jet, z_jet
103 !  perturbation parameters
105    REAL, PARAMETER :: htbub=8000., radbub=2000000., radz=8000., tpbub=1.0
106    REAL :: piov2, tp
107    INTEGER :: icen, jcen
108    real :: thtmp, ptmp, temp(3)
110    SELECT CASE ( model_data_order )
111          CASE ( DATA_ORDER_ZXY )
112    kds = grid%sd31 ; kde = grid%ed31 ;
113    ids = grid%sd32 ; ide = grid%ed32 ;
114    jds = grid%sd33 ; jde = grid%ed33 ;
116    kms = grid%sm31 ; kme = grid%em31 ;
117    ims = grid%sm32 ; ime = grid%em32 ;
118    jms = grid%sm33 ; jme = grid%em33 ;
120    kts = grid%sp31 ; kte = grid%ep31 ;   ! note that tile is entire patch
121    its = grid%sp32 ; ite = grid%ep32 ;   ! note that tile is entire patch
122    jts = grid%sp33 ; jte = grid%ep33 ;   ! note that tile is entire patch
123          CASE ( DATA_ORDER_XYZ )
124    ids = grid%sd31 ; ide = grid%ed31 ;
125    jds = grid%sd32 ; jde = grid%ed32 ;
126    kds = grid%sd33 ; kde = grid%ed33 ;
128    ims = grid%sm31 ; ime = grid%em31 ;
129    jms = grid%sm32 ; jme = grid%em32 ;
130    kms = grid%sm33 ; kme = grid%em33 ;
132    its = grid%sp31 ; ite = grid%ep31 ;   ! note that tile is entire patch
133    jts = grid%sp32 ; jte = grid%ep32 ;   ! note that tile is entire patch
134    kts = grid%sp33 ; kte = grid%ep33 ;   ! note that tile is entire patch
135          CASE ( DATA_ORDER_XZY )
136    ids = grid%sd31 ; ide = grid%ed31 ;
137    kds = grid%sd32 ; kde = grid%ed32 ;
138    jds = grid%sd33 ; jde = grid%ed33 ;
140    ims = grid%sm31 ; ime = grid%em31 ;
141    kms = grid%sm32 ; kme = grid%em32 ;
142    jms = grid%sm33 ; jme = grid%em33 ;
144    its = grid%sp31 ; ite = grid%ep31 ;   ! note that tile is entire patch
145    kts = grid%sp32 ; kte = grid%ep32 ;   ! note that tile is entire patch
146    jts = grid%sp33 ; jte = grid%ep33 ;   ! note that tile is entire patch
148    END SELECT
150    piov2 = 2.*atan(1.0)
151    icen = ide/4
152    jcen = jde/2
154    stretch_grid = .true.
155    delt = 0.
156    z_scale = .50
157    pi = 2.*asin(1.0)
158    write(6,*) ' pi is ',pi
159    nxc = (ide-ids)/4
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
194          grid%ht(i,j)       = 0.
195          grid%msftx(i,j)    = 1.
196          grid%msfty(i,j)    = 1.
197          grid%msfux(i,j)    = 1.
198          grid%msfuy(i,j)    = 1.
199          grid%msfvx(i,j)    = 1.
200          grid%msfvx_inv(i,j)= 1.
201          grid%msfvy(i,j)    = 1.
202          grid%sina(i,j)     = 0.
203          grid%cosa(i,j)     = 1.
204          grid%e(i,j)        = 0.
205          grid%f(i,j)        = 1.e-04
207       END DO
208    END DO
210     DO j = jts, jte
211     DO k = kts, kte
212       DO i = its, ite
213          grid%ww(i,k,j)     = 0.
214       END DO
215    END DO
216    END DO
218    grid%step_number = 0
220 ! set up the grid
222    IF (stretch_grid) THEN ! exponential stretch for eta (nearly constant dz)
223      DO k=1, kde
224       grid%znw(k) = (exp(-(k-1)/float(kde-1)/z_scale) - exp(-1./z_scale))/ &
225                                 (1.-exp(-1./z_scale))
226      ENDDO
227    ELSE
228      DO k=1, kde
229       grid%znw(k) = 1. - float(k-1)/float(kde-1)
230      ENDDO
231    ENDIF
233    DO k=1, kde-1
234     grid%dnw(k) = grid%znw(k+1) - grid%znw(k)
235     grid%rdnw(k) = 1./grid%dnw(k)
236     grid%znu(k) = 0.5*(grid%znw(k+1)+grid%znw(k))
237    ENDDO
238    DO k=2, kde-1
239     grid%dn(k) = 0.5*(grid%dnw(k)+grid%dnw(k-1))
240     grid%rdn(k) = 1./grid%dn(k)
241     grid%fnp(k) = .5* grid%dnw(k  )/grid%dn(k)
242     grid%fnm(k) = .5* grid%dnw(k-1)/grid%dn(k)
243    ENDDO
245    cof1 = (2.*grid%dn(2)+grid%dn(3))/(grid%dn(2)+grid%dn(3))*grid%dnw(1)/grid%dn(2) 
246    cof2 =     grid%dn(2)        /(grid%dn(2)+grid%dn(3))*grid%dnw(1)/grid%dn(3) 
247    grid%cf1  = grid%fnp(2) + cof1
248    grid%cf2  = grid%fnm(2) - cof1 - cof2
249    grid%cf3  = cof2       
251    grid%cfn  = (.5*grid%dnw(kde-1)+grid%dn(kde-1))/grid%dn(kde-1)
252    grid%cfn1 = -.5*grid%dnw(kde-1)/grid%dn(kde-1)
253    grid%rdx = 1./config_flags%dx
254    grid%rdy = 1./config_flags%dy
256 !  get the sounding from the ascii sounding file, first get dry sounding and 
257 !  calculate base state
259   write(6,*) ' reading input jet sounding '
260   call read_input_jet( u_jet, rho_jet, th_jet, z_jet, nz_jet, ny_jet )
262   write(6,*) ' getting dry sounding for base state '
263   write(6,*) ' using middle column in jet sounding, j = ',ny_jet/2
264   dry_sounding = .true.
266   dry_sounding   = .true.
267   debug = .true.  !  this will produce print of the sounding
268   CALL get_sounding( zk, p_in, pd_in, theta, rho, u, v, qv, dry_sounding, &
269                       nl_max, nl_in, u_jet, rho_jet, th_jet, z_jet,      &
270                       nz_jet, ny_jet, ny_jet/2, debug                   )
272   write(6,*) ' returned from reading sounding, nl_in is ',nl_in
274 !  find ptop for the desired ztop (ztop is input from the namelist),
275 !  and find surface pressure
277 !  For the jet, using the middle column for the base state means that
278 !  we will be extrapolating above the highest height data to the south
279 !  of the centerline.
281   grid%p_top = interp_0( p_in, zk, config_flags%ztop, nl_in )
283   DO j=jts,jte
284   DO i=its,ite  ! flat surface
285     grid%phb(i,1,j) = 0.
286     grid%php(i,1,j) = 0.
287     grid%ph0(i,1,j) = 0.
288     grid%ht(i,j) = 0.
289   ENDDO
290   ENDDO
292   DO J = jts, jte
293   DO I = its, ite
295     p_surf = interp_0( p_in, zk, grid%phb(i,1,j)/g, nl_in )
296     grid%mub(i,j) = p_surf-grid%p_top
298 !  this is dry hydrostatic sounding (base state), so given grid%p (coordinate),
299 !  interp theta (from interp) and compute 1/rho from eqn. of state
301     DO K = 1, kte-1
302       p_level = grid%znu(k)*(p_surf - grid%p_top) + grid%p_top
303       grid%pb(i,k,j) = p_level
304       grid%t_init(i,k,j) = interp_0( theta, p_in, p_level, nl_in ) - t0
305       grid%alb(i,k,j) = (r_d/p1000mb)*(grid%t_init(i,k,j)+t0)*(grid%pb(i,k,j)/p1000mb)**cvpm
306     ENDDO
308 !  calc hydrostatic balance (alternatively we could interp the geopotential from the
309 !  sounding, but this assures that the base state is in exact hydrostatic balance with
310 !  respect to the model eqns.
312     DO k  = 2,kte
313       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)
314     ENDDO
316   ENDDO
317   ENDDO
319   write(6,*) ' ptop is ',grid%p_top
320   write(6,*) ' base state grid%mub(1,1), p_surf is ',grid%mub(1,1),grid%mub(1,1)+grid%p_top
322 !  calculate full state for each column - this includes moisture.
324   write(6,*) ' getting grid%moist sounding for full state '
326   dry_sounding = .true.
327   IF (config_flags%mp_physics /= 0)  dry_sounding = .false.
329   DO J = jts, min(jde-1,jte)
331 !  get sounding for this point
333   debug = .false.  !  this will turn off print of the sounding
334   CALL get_sounding( zk, p_in, pd_in, theta, rho, u, v, qv, dry_sounding, &
335                       nl_max, nl_in, u_jet, rho_jet, th_jet, z_jet,      &
336                       nz_jet, ny_jet, j, debug                          )
338   DO I = its, min(ide-1,ite)
340 !   we could just do the first point in "i" and copy from there, but we'll
341 !   be lazy and do all the points as if they are all, independent
343 !   At this point grid%p_top is already set. find the DRY mass in the column 
344 !   by interpolating the DRY pressure.  
346     pd_surf = interp_0( pd_in, zk, grid%phb(i,1,j)/g, nl_in )
348 !   compute the perturbation mass and the full mass
350     grid%mu_1(i,j) = pd_surf-grid%p_top - grid%mub(i,j)
351     grid%mu_2(i,j) = grid%mu_1(i,j)
352     grid%mu0(i,j) = grid%mu_1(i,j) + grid%mub(i,j)
354 !   given the dry pressure and coordinate system, interp the potential
355 !   temperature and qv
357     do k=1,kde-1
359       p_level = grid%znu(k)*(pd_surf - grid%p_top) + grid%p_top
361       grid%moist(i,k,j,P_QV) = interp_0( qv, pd_in, p_level, nl_in )
362       grid%t_1(i,k,j)          = interp_0( theta, pd_in, p_level, nl_in ) - t0
363       grid%t_2(i,k,j)          = grid%t_1(i,k,j)
364       
366     enddo
368 !   integrate the hydrostatic equation (from the RHS of the bigstep
369 !   vertical momentum equation) down from the top to get grid%p.
370 !   first from the top of the model to the top pressure
372     k = kte-1  ! top level
374     qvf1 = 0.5*(grid%moist(i,k,j,P_QV)+grid%moist(i,k,j,P_QV))
375     qvf2 = 1./(1.+qvf1)
376     qvf1 = qvf1*qvf2
378 !    grid%p(i,k,j) = - 0.5*grid%mu_1(i,j)/grid%rdnw(k)
379     grid%p(i,k,j) = - 0.5*(grid%mu_1(i,j)+qvf1*grid%mub(i,j))/grid%rdnw(k)/qvf2
380     qvf = 1. + rvovrd*grid%moist(i,k,j,P_QV)
381     grid%alt(i,k,j) = (r_d/p1000mb)*(grid%t_1(i,k,j)+t0)*qvf* &
382                 (((grid%p(i,k,j)+grid%pb(i,k,j))/p1000mb)**cvpm)
383     grid%al(i,k,j) = grid%alt(i,k,j) - grid%alb(i,k,j)
385 !  down the column
387     do k=kte-2,1,-1
388       qvf1 = 0.5*(grid%moist(i,k,j,P_QV)+grid%moist(i,k+1,j,P_QV))
389       qvf2 = 1./(1.+qvf1)
390       qvf1 = qvf1*qvf2
391       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)
392       qvf = 1. + rvovrd*grid%moist(i,k,j,P_QV)
393       grid%alt(i,k,j) = (r_d/p1000mb)*(grid%t_1(i,k,j)+t0)*qvf* &
394                   (((grid%p(i,k,j)+grid%pb(i,k,j))/p1000mb)**cvpm)
395       grid%al(i,k,j) = grid%alt(i,k,j) - grid%alb(i,k,j)
396     enddo
398 !  this is the hydrostatic equation used in the model after the
399 !  small timesteps.  In the model, grid%al (inverse density)
400 !  is computed from the geopotential.
403     grid%ph_1(i,1,j) = 0.
404     DO k  = 2,kte
405       grid%ph_1(i,k,j) = grid%ph_1(i,k-1,j) - (1./grid%rdnw(k-1))*(       &
406                    (grid%mub(i,j)+grid%mu_1(i,j))*grid%al(i,k-1,j)+ &
407                     grid%mu_1(i,j)*grid%alb(i,k-1,j)  )
408                                                    
409       grid%ph_2(i,k,j) = grid%ph_1(i,k,j) 
410       grid%ph0(i,k,j) = grid%ph_1(i,k,j) + grid%phb(i,k,j)
411     ENDDO
413 ! interp u
415     DO K = 1, kte
416       p_level = grid%znu(k)*(p_surf - grid%p_top) + grid%p_top
417       grid%u_1(i,k,j) = interp_0( u, p_in, p_level, nl_in )
418       grid%u_2(i,k,j) = grid%u_1(i,k,j)
419     ENDDO
421   ENDDO
422   ENDDO
424 !  thermal perturbation to kick off convection
426   write(6,*) ' nxc, nyc for perturbation ',nxc,nyc
427   write(6,*) ' delt for perturbation ',tpbub
429   DO J = jts, min(jde-1,jte)
430     yrad = config_flags%dy*float(j-jde/2-1)/radbub
431     DO I = its, min(ide-1,ite)
432       xrad = float(i-1)/float(ide-ids)
434       DO K = 1, kte-1
436 !  put in preturbation theta (bubble) and recalc density.  note,
437 !  the mass in the column is not changing, so when theta changes,
438 !  we recompute density and geopotential
440         zrad = 0.5*(grid%ph_1(i,k,j)+grid%ph_1(i,k+1,j)  &
441                    +grid%phb(i,k,j)+grid%phb(i,k+1,j))/g
442         zrad = (zrad-htbub)/radz
443         RAD=SQRT(yrad*yrad+zrad*zrad)
444         IF(RAD <= 1.) THEN
445            tp = tpbub*cos(rad*piov2)*cos(rad*piov2)*cos(xrad*2*pi+pi)
446            grid%t_1(i,k,j)=grid%t_1(i,k,j)+tp
447            grid%t_2(i,k,j)=grid%t_1(i,k,j)
448            qvf = 1. + rvovrd*grid%moist(i,k,j,P_QV)
449            grid%alt(i,k,j) = (r_d/p1000mb)*(grid%t_1(i,k,j)+t0)*qvf* &
450                         (((grid%p(i,k,j)+grid%pb(i,k,j))/p1000mb)**cvpm)
451            grid%al(i,k,j) = grid%alt(i,k,j) - grid%alb(i,k,j)
452         ENDIF
453       ENDDO
455 !  rebalance hydrostatically
457       DO k  = 2,kte
458         grid%ph_1(i,k,j) = grid%ph_1(i,k-1,j) - (1./grid%rdnw(k-1))*(       &
459                      (grid%mub(i,j)+grid%mu_1(i,j))*grid%al(i,k-1,j)+ &
460                       grid%mu_1(i,j)*grid%alb(i,k-1,j)  )
461                                                    
462         grid%ph_2(i,k,j) = grid%ph_1(i,k,j) 
463         grid%ph0(i,k,j) = grid%ph_1(i,k,j) + grid%phb(i,k,j)
464       ENDDO
466     ENDDO
467   ENDDO
469 !#endif
471    write(6,*) ' grid%mu_1 from comp ', grid%mu_1(1,1)
472    write(6,*) ' pert state sounding from comp, grid%ph_1, pp, grid%al, grid%t_1, qv '
473    do k=1,kde-1
474      write(6,'(i3,1x,5(1x,1pe10.3))') k, grid%ph_1(1,k,1),grid%p(1,k,1), grid%al(1,k,1), &
475                      grid%t_1(1,k,1), grid%moist(1,k,1,P_QV)
476    enddo
478    write(6,*) ' grid%mu_1 from comp ', grid%mu_1(1,1)
479    write(6,*) ' full state sounding from comp, ph, grid%p, grid%al, grid%t_1, qv '
480    write(6,*) ' at j = 1 '
481    do k=1,kde-1
482      write(6,'(i3,1x,5(1x,1pe10.3))') k, grid%ph_1(1,k,1)+grid%phb(1,k,1), &
483                      grid%p(1,k,1)+grid%pb(1,k,1), grid%al(1,k,1)+grid%alb(1,k,1), &
484                      grid%t_1(1,k,1)+t0, grid%moist(1,k,1,P_QV)
485    enddo
488    write(6,*) ' full state sounding from comp, ph, grid%p, grid%al, grid%t_1, qv '
489    write(6,*) ' at j = jde/2 '
490    do k=1,kde-1
491      write(6,'(i3,1x,5(1x,1pe10.3))') k, grid%ph_1(1,k,jde/2)+grid%phb(1,k,jde/2), &
492                      grid%p(1,k,jde/2)+grid%pb(1,k,jde/2), grid%al(1,k,jde/2)+grid%alb(1,k,jde/2), &
493                      grid%t_1(1,k,jde/2)+t0, grid%moist(1,k,jde/2,P_QV)
494    enddo
496    write(6,*) ' full state sounding from comp, ph, grid%p, grid%al, grid%t_1, qv '
497    write(6,*) ' at j = jde-1 '
498    do k=1,kde-1
499      write(6,'(i3,1x,5(1x,1pe10.3))') k, grid%ph_1(1,k,jde-1)+grid%phb(1,k,jde-1), &
500                      grid%p(1,k,jde-1)+grid%pb(1,k,jde-1), grid%al(1,k,jde-1)+grid%alb(1,k,jde-1), &
501                      grid%t_1(1,k,jde-1)+t0, grid%moist(1,k,jde-1,P_QV)
502    enddo
504 ! set v
506   DO J = jts, jte
507   DO I = its, min(ide-1,ite)
509     DO K = 1, kte
510       grid%v_1(i,k,j) = 0.
511       grid%v_2(i,k,j) = grid%v_1(i,k,j)
512     ENDDO
514   ENDDO
515   ENDDO
517 !  fill out last i row for u
519   DO J = jts, min(jde-1,jte)
520   DO I = ite, ite
522     DO K = 1, kte
523       grid%u_1(i,k,j) = grid%u_1(its,k,j)
524       grid%u_2(i,k,j) = grid%u_2(its,k,j)
525     ENDDO
527   ENDDO
528   ENDDO
530 !  set w
532   DO J = jts, min(jde-1,jte)
533   DO K = kts, kte
534   DO I = its, min(ide-1,ite)
535     grid%w_1(i,k,j) = 0.
536     grid%w_2(i,k,j) = 0.
537   ENDDO
538   ENDDO
539   ENDDO
541 !  set a few more things
543   DO J = jts, min(jde-1,jte)
544   DO K = kts, kte-1
545   DO I = its, min(ide-1,ite)
546     grid%h_diabatic(i,k,j) = 0.
547   ENDDO
548   ENDDO
549   ENDDO
551   DO k=1,kte-1
552     grid%t_base(k) = grid%t_1(1,k,1)
553     grid%qv_base(k) = grid%moist(1,k,1,P_QV)
554     grid%u_base(k) = grid%u_1(1,k,1)
555     grid%v_base(k) = grid%v_1(1,k,1)
556   ENDDO
558   DO J = jts, min(jde-1,jte)
559   DO I = its, min(ide-1,ite)
560      thtmp   = grid%t_2(i,1,j)+t0
561      ptmp    = grid%p(i,1,j)+grid%pb(i,1,j)
562      temp(1) = thtmp * (ptmp/p1000mb)**rcp
563      thtmp   = grid%t_2(i,2,j)+t0
564      ptmp    = grid%p(i,2,j)+grid%pb(i,2,j)
565      temp(2) = thtmp * (ptmp/p1000mb)**rcp
566      thtmp   = grid%t_2(i,3,j)+t0
567      ptmp    = grid%p(i,3,j)+grid%pb(i,3,j)
568      temp(3) = thtmp * (ptmp/p1000mb)**rcp
570      grid%tsk(I,J)=grid%cf1*temp(1)+grid%cf2*temp(2)+grid%cf3*temp(3)
571      if (i .eq. 1) print*,'sfctem',j,temp(1),temp(2),temp(3),grid%tsk(I,J)
572      grid%tmn(I,J)=grid%tsk(I,J)-0.5
573   ENDDO
574   ENDDO
576   RETURN
578  END SUBROUTINE init_domain_rk
580 !---------------------------------------------------------------------
582  SUBROUTINE init_module_initialize
583  END SUBROUTINE init_module_initialize
585 !---------------------------------------------------------------------
586 #if 0
587 ! TEST DRIVER FOR "read_input_jet" and "get_sounding"
588   implicit none 
589   integer, parameter :: nz_jet=64, ny_jet=80
590   real, dimension(nz_jet,ny_jet) :: u_jet, rho_jet, &
591                                     th_jet, z_jet
593   real, dimension(nz_jet,ny_jet) :: zk,p,p_dry,theta,rho,u,v,qv
594   logical :: dry, debug
595   integer :: j, nl
597   call read_input_jet( u_jet, rho_jet, th_jet, z_jet, nz_jet, ny_jet )
599   call opngks
600   call parray( u_jet, nz_jet, ny_jet)
601   call parray( rho_jet, nz_jet, ny_jet)
602   call parray( th_jet, nz_jet, ny_jet)
603 !  call clsgks
605 !  set up initial jet
607   debug = .true.
608   dry = .true.
609   do j=1,ny_jet
611     call get_sounding( zk(:,j),p(:,j),p_dry(:,j),theta(:,j),      &
612                        rho(:,j),u(:,j), v(:,j),  qv(:,j),        &
613                        dry, nz_jet, nl, u_jet, rho_jet, th_jet,  &
614                        z_jet, nz_jet, ny_jet, j, debug          )
615     debug = .false.
617   enddo
619   write(6,*) ' lowest level p, th, and rho, highest level p '
621   do j=1,ny_jet
622     write(6,*) j, p(1,j),theta(1,j),rho(1,j), p(nz_jet,j)
623 !    write(6,*) j, p(1,j),theta(1,j)-th_jet(1,j),rho(1,j)-rho_jet(1,j)
624   enddo
626   call parray( p, nz_jet, ny_jet)
627   call parray( p_dry, nz_jet, ny_jet)
628   call parray( theta, nz_jet, ny_jet)
630   call clsgks
632   end
634 !---------------------------------
636       subroutine parray(a,m,n)
637       dimension a(m,n)
638       dimension b(n,m)
640     do i=1,m
641     do j=1,n
642       b(j,i) = a(i,j)
643     enddo
644     enddo
645       
646       write(6,'(''  dimensions m,n  '',2i6)')m,n
647         call set(.05,.95,.05,.95,0.,1.,0.,1.,1)
648         call perim(4,5,4,5)
649         call setusv('LW',2000)
650 !        CALL CONREC(a,m,m,n,cmax,cmin,cinc,-1,-638,-922)
651         CALL CONREC(b,n,n,m,0.,0.,0.,-1,-638,-922)
652         call frame
653       return
654       end
656 ! END TEST DRIVER FOR "read_input_jet" and "get_sounding"
657 #endif
659 !------------------------------------------------------------------
661     subroutine get_sounding( zk, p, p_dry, theta, rho,       &
662                              u, v, qv, dry, nl_max, nl_in,  &
663                              u_jet, rho_jet, th_jet, z_jet, &
664                              nz_jet, ny_jet, j_point, debug )
665     implicit none
667     integer nl_max, nl_in
668     real zk(nl_max), p(nl_max), theta(nl_max), rho(nl_max), &
669          u(nl_max), v(nl_max), qv(nl_max), p_dry(nl_max)
670     logical dry
672     integer nz_jet, ny_jet, j_point
673     real, dimension(nz_jet, ny_jet) :: u_jet, rho_jet, th_jet, z_jet
675     integer n
676     parameter(n=1000)
677     logical debug
679 ! input sounding data
681     real p_surf, th_surf, qv_surf
682     real pi_surf, pi(n)
683     real h_input(n), th_input(n), qv_input(n), u_input(n), v_input(n)
685 ! diagnostics
687     real rho_surf, p_input(n), rho_input(n)
688     real pm_input(n)  !  this are for full moist sounding
690 ! local data
692     real p1000mb,cv,cp,r,cvpm,g
693     parameter (p1000mb = 1.e+05, r = 287, cp = 1003., cv = cp-r, cvpm = -cv/cp, g=9.81 )
694     integer k, it, nl
695     real qvf, qvf1, dz
697 !  first, read the sounding
699 !    call read_sounding( p_surf, th_surf, qv_surf, &
700 !                          h_input, th_input, qv_input, u_input, v_input,n, nl, debug )
702    call calc_jet_sounding( p_surf, th_surf, qv_surf,                             &
703                            h_input, th_input, qv_input, u_input, v_input,        &
704                            n, nl, debug, u_jet, rho_jet, th_jet, z_jet, j_point, &
705                            nz_jet, ny_jet, dry                                  )
707    nl = nz_jet
709     if(dry) then
710      do k=1,nl
711        qv_input(k) = 0.
712      enddo
713     endif
715     if(debug) write(6,*) ' number of input levels = ',nl
717       nl_in = nl
718       if(nl_in .gt. nl_max ) then
719         write(6,*) ' too many levels for input arrays ',nl_in,nl_max
720         call wrf_error_fatal ( ' too many levels for input arrays ' )
721       end if
723 !  compute diagnostics,
724 !  first, convert qv(g/kg) to qv(g/g)
726 !      do k=1,nl
727 !        qv_input(k) = 0.001*qv_input(k)
728 !      enddo
729 !      p_surf = 100.*p_surf  ! convert to pascals
731     qvf = 1. + rvovrd*qv_input(1) 
732     rho_surf = 1./((r/p1000mb)*th_surf*qvf*((p_surf/p1000mb)**cvpm))
733     pi_surf = (p_surf/p1000mb)**(r/cp)
735     if(debug) then
736       write(6,*) ' surface density is ',rho_surf
737       write(6,*) ' surface pi is    ',pi_surf
738     end if
741 !  integrate moist sounding hydrostatically, starting from the
742 !  specified surface pressure
743 !  -> first, integrate from surface to lowest level
745         qvf = 1. + rvovrd*qv_input(1) 
746         qvf1 = 1. + qv_input(1)
747         rho_input(1) = rho_surf
748         dz = h_input(1)
749         do it=1,10
750           pm_input(1) = p_surf &
751                   - 0.5*dz*(rho_surf+rho_input(1))*g*qvf1
752           rho_input(1) = 1./((r/p1000mb)*th_input(1)*qvf*((pm_input(1)/p1000mb)**cvpm))
753         enddo
755 ! integrate up the column
757         do k=2,nl
758           rho_input(k) = rho_input(k-1)
759           dz = h_input(k)-h_input(k-1)
760           qvf1 = 0.5*(2.+(qv_input(k-1)+qv_input(k)))
761           qvf = 1. + rvovrd*qv_input(k)   ! qv is in g/kg here
763           do it=1,10
764             pm_input(k) = pm_input(k-1) &
765                     - 0.5*dz*(rho_input(k)+rho_input(k-1))*g*qvf1
766             rho_input(k) = 1./((r/p1000mb)*th_input(k)*qvf*((pm_input(k)/p1000mb)**cvpm))
767           enddo
768         enddo
770 !  we have the moist sounding
772 !  next, compute the dry sounding using p at the highest level from the
773 !  moist sounding and integrating down.
775         p_input(nl) = pm_input(nl)
777           do k=nl-1,1,-1
778             dz = h_input(k+1)-h_input(k)
779             p_input(k) = p_input(k+1) + 0.5*dz*(rho_input(k)+rho_input(k+1))*g
780           enddo
783         do k=1,nl
785           zk(k) = h_input(k)
786           p(k) = pm_input(k)
787           p_dry(k) = p_input(k)
788           theta(k) = th_input(k)
789           rho(k) = rho_input(k)
790           u(k) = u_input(k)
791           v(k) = v_input(k)
792           qv(k) = qv_input(k)
794         enddo
796      if(debug) then
797       write(6,*) ' sounding '
798       write(6,*) '  k  height(m)  press (Pa)   pd(Pa)   theta (K)  den(kg/m^3)  u(m/s)     v(m/s)    qv(g/g) '
799       do k=1,nl
800         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)
801       enddo
803      end if
805      end subroutine get_sounding
807 !------------------------------------------------------------------
809   subroutine calc_jet_sounding( p_surf, th_surf, qv_surf,      &
810                                 h, th, qv, u, v, n, nl, debug, &
811                                 u_jet, rho_jet, th_jet, z_jet, &
812                                 jp, nz_jet, ny_jet, dry       )
813   implicit none
814   integer :: n, nl, jp, nz_jet, ny_jet
816   real, dimension(nz_jet, ny_jet) :: u_jet, rho_jet, th_jet, z_jet
817   real, dimension(n) :: h,th,qv,u,v
818   real :: p_surf, th_surf, qv_surf
819   logical :: debug, dry
821   real, dimension(1:nz_jet) :: rho, rel_hum, p
822   integer :: k
824 !  some local stuff
826   real :: tmppi, es, qvs, temperature
827   real, parameter :: p1000mb=1.e+05, rcp=287./1004.5, svpt0=273.15, &
828                      svp3 = 29.65, ep_2=287./461.6, r_d = 287., &
829                      cpovcv = 1004./(1004.-287.),               &
830                      svp1 = 0.6112, svp2 = 17.67
832 !  get sounding from column jp
834    do k=1,nz_jet
835      h(k)  = z_jet(k,jp)
836      th(k) = th_jet(k,jp)
837      qv(k) = 0.
838      rho(k) = rho_jet(k,jp)
839      u(k) = u_jet(k,jp)
840      v(k) = 0.
841    enddo
843    if (.not.dry) then
844      DO k=1,nz_jet
845        if(h(k) .gt. 8000.) then
846          rel_hum(k)=0.1
847        else
848          rel_hum(k)=(1.-0.90*(h(k)/8000.)**1.25)
849        end if
850        rel_hum(k) = min(0.7,rel_hum(k))
851      ENDDO
852    else
853      do k=1,nz_jet
854        rel_hum(k) = 0.
855      enddo
856    endif
858 !  next, compute pressure
860    do k=1,nz_jet
861      p(k) = p1000mb*(R_d*rho(k)*th(k)/p1000mb)**cpovcv
862    enddo
864 !  here we adjust for fixed moisture profile
866      IF (.not.dry)  THEN
868 !  here we assume the input theta is th_v, so we reset theta accordingly
870        DO k=1,nz_jet
871          tmppi=(p(k)/p1000mb)**rcp
872          temperature = tmppi*th(k)
873          if (temperature .gt. svpt0) then
874             es  = 1000.*svp1*exp(svp2*(temperature-svpt0)/(temperature-svp3))
875             qvs = ep_2*es/(p(k)-es)
876          else
877             es  = 1000.*svp1*exp( 21.8745584*(temperature-273.16)/(temperature-7.66) )
878             qvs = ep_2*es/(p(k)-es)
879          endif
880          qv(k) = rel_hum(k)*qvs
881          th(k) = th(k)/(1.+.61*qv(k))
882        ENDDO
884      ENDIF
886 !  finally, set the surface data. We'll just do a simple extrapolation
888    p_surf = 1.5*p(1) - 0.5*p(2)
889    th_surf = 1.5*th(1) - 0.5*th(2)
890    qv_surf = 1.5*qv(1) - 0.5*qv(2)
892    end subroutine calc_jet_sounding
894 !---------------------------------------------------------------------
896  SUBROUTINE read_input_jet( u, r, t, zk, nz, ny )
897  implicit none
899  integer, intent(in) :: nz,ny
900  real, dimension(nz,ny), intent(out) :: u,r,t,zk
901  integer :: ny_in, nz_in, j,k
902  real, dimension(ny,nz) :: field_in
904 ! this code assumes it is called on processor 0 only
906    OPEN(unit=10, file='input_jet', form='unformatted', status='old' )
907    REWIND(10) 
908    read(10) ny_in,nz_in
909    if((ny_in /= ny ) .or. (nz_in /= nz)) then
910      write(0,*) ' error in input jet dimensions '
911      write(0,*) ' ny, ny_input, nz, nz_input ', ny, ny_in, nz,nz_in
912      write(0,*) ' error exit '
913      call wrf_error_fatal ( ' error in input jet dimensions ' )
914    end if
915    read(10) field_in
916    do j=1,ny
917    do k=1,nz
918      u(k,j) = field_in(j,k)
919    enddo
920    enddo
921    read(10) field_in
922    do j=1,ny
923    do k=1,nz
924      t(k,j) = field_in(j,k)
925    enddo
926    enddo
928    read(10) field_in
929    do j=1,ny
930    do k=1,nz
931      r(k,j) = field_in(j,k)
932    enddo
933    enddo
935    do j=1,ny
936    do k=1,nz
937      zk(k,j) = 125. + 250.*float(k-1)
938    enddo
939    enddo
941  end subroutine read_input_jet
943 END MODULE module_initialize_ideal