r4627 | gill | 2010-12-29 16:29:58 -0700 (Wed, 29 Dec 2010) | 5 lines
[wrffire.git] / wrfv2_fire / dyn_em / module_initialize_b_wave.F
blob44465a22099b54cbfb568381819f18852543da58
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_pole_lon (1,0.)
186     CALL nl_set_pole_lat (1,90.)
187     CALL nl_set_map_proj(1,0)
190 !  here we initialize data we currently is not initialized 
191 !  in the input data
193     DO j = jts, jte
194       DO i = its, ite
196          grid%ht(i,j)       = 0.
197          grid%msftx(i,j)    = 1.
198          grid%msfty(i,j)    = 1.
199          grid%msfux(i,j)    = 1.
200          grid%msfuy(i,j)    = 1.
201          grid%msfvx(i,j)    = 1.
202          grid%msfvx_inv(i,j)= 1.
203          grid%msfvy(i,j)    = 1.
204          grid%sina(i,j)     = 0.
205          grid%cosa(i,j)     = 1.
206          grid%e(i,j)        = 0.
207          grid%f(i,j)        = 1.e-04
209       END DO
210    END DO
212     DO j = jts, jte
213     DO k = kts, kte
214       DO i = its, ite
215          grid%ww(i,k,j)     = 0.
216       END DO
217    END DO
218    END DO
220    grid%step_number = 0
222 ! set up the grid
224    IF (stretch_grid) THEN ! exponential stretch for eta (nearly constant dz)
225      DO k=1, kde
226       grid%znw(k) = (exp(-(k-1)/float(kde-1)/z_scale) - exp(-1./z_scale))/ &
227                                 (1.-exp(-1./z_scale))
228      ENDDO
229    ELSE
230      DO k=1, kde
231       grid%znw(k) = 1. - float(k-1)/float(kde-1)
232      ENDDO
233    ENDIF
235    DO k=1, kde-1
236     grid%dnw(k) = grid%znw(k+1) - grid%znw(k)
237     grid%rdnw(k) = 1./grid%dnw(k)
238     grid%znu(k) = 0.5*(grid%znw(k+1)+grid%znw(k))
239    ENDDO
240    DO k=2, kde-1
241     grid%dn(k) = 0.5*(grid%dnw(k)+grid%dnw(k-1))
242     grid%rdn(k) = 1./grid%dn(k)
243     grid%fnp(k) = .5* grid%dnw(k  )/grid%dn(k)
244     grid%fnm(k) = .5* grid%dnw(k-1)/grid%dn(k)
245    ENDDO
247    cof1 = (2.*grid%dn(2)+grid%dn(3))/(grid%dn(2)+grid%dn(3))*grid%dnw(1)/grid%dn(2) 
248    cof2 =     grid%dn(2)        /(grid%dn(2)+grid%dn(3))*grid%dnw(1)/grid%dn(3) 
249    grid%cf1  = grid%fnp(2) + cof1
250    grid%cf2  = grid%fnm(2) - cof1 - cof2
251    grid%cf3  = cof2       
253    grid%cfn  = (.5*grid%dnw(kde-1)+grid%dn(kde-1))/grid%dn(kde-1)
254    grid%cfn1 = -.5*grid%dnw(kde-1)/grid%dn(kde-1)
255    grid%rdx = 1./config_flags%dx
256    grid%rdy = 1./config_flags%dy
258 !  get the sounding from the ascii sounding file, first get dry sounding and 
259 !  calculate base state
261   write(6,*) ' reading input jet sounding '
262   call read_input_jet( u_jet, rho_jet, th_jet, z_jet, nz_jet, ny_jet )
264   write(6,*) ' getting dry sounding for base state '
265   write(6,*) ' using middle column in jet sounding, j = ',ny_jet/2
266   dry_sounding = .true.
268   dry_sounding   = .true.
269   debug = .true.  !  this will produce print of the sounding
270   CALL get_sounding( zk, p_in, pd_in, theta, rho, u, v, qv, dry_sounding, &
271                       nl_max, nl_in, u_jet, rho_jet, th_jet, z_jet,      &
272                       nz_jet, ny_jet, ny_jet/2, debug                   )
274   write(6,*) ' returned from reading sounding, nl_in is ',nl_in
276 !  find ptop for the desired ztop (ztop is input from the namelist),
277 !  and find surface pressure
279 !  For the jet, using the middle column for the base state means that
280 !  we will be extrapolating above the highest height data to the south
281 !  of the centerline.
283   grid%p_top = interp_0( p_in, zk, config_flags%ztop, nl_in )
285   DO j=jts,jte
286   DO i=its,ite  ! flat surface
287     grid%phb(i,1,j) = 0.
288     grid%php(i,1,j) = 0.
289     grid%ph0(i,1,j) = 0.
290     grid%ht(i,j) = 0.
291   ENDDO
292   ENDDO
294   DO J = jts, jte
295   DO I = its, ite
297     p_surf = interp_0( p_in, zk, grid%phb(i,1,j)/g, nl_in )
298     grid%mub(i,j) = p_surf-grid%p_top
300 !  this is dry hydrostatic sounding (base state), so given grid%p (coordinate),
301 !  interp theta (from interp) and compute 1/rho from eqn. of state
303     DO K = 1, kte-1
304       p_level = grid%znu(k)*(p_surf - grid%p_top) + grid%p_top
305       grid%pb(i,k,j) = p_level
306       grid%t_init(i,k,j) = interp_0( theta, p_in, p_level, nl_in ) - t0
307       grid%alb(i,k,j) = (r_d/p1000mb)*(grid%t_init(i,k,j)+t0)*(grid%pb(i,k,j)/p1000mb)**cvpm
308     ENDDO
310 !  calc hydrostatic balance (alternatively we could interp the geopotential from the
311 !  sounding, but this assures that the base state is in exact hydrostatic balance with
312 !  respect to the model eqns.
314     DO k  = 2,kte
315       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)
316     ENDDO
318   ENDDO
319   ENDDO
321   write(6,*) ' ptop is ',grid%p_top
322   write(6,*) ' base state grid%mub(1,1), p_surf is ',grid%mub(1,1),grid%mub(1,1)+grid%p_top
324 !  calculate full state for each column - this includes moisture.
326   write(6,*) ' getting grid%moist sounding for full state '
328   dry_sounding = .true.
329   IF (config_flags%mp_physics /= 0)  dry_sounding = .false.
331   DO J = jts, min(jde-1,jte)
333 !  get sounding for this point
335   debug = .false.  !  this will turn off print of the sounding
336   CALL get_sounding( zk, p_in, pd_in, theta, rho, u, v, qv, dry_sounding, &
337                       nl_max, nl_in, u_jet, rho_jet, th_jet, z_jet,      &
338                       nz_jet, ny_jet, j, debug                          )
340   DO I = its, min(ide-1,ite)
342 !   we could just do the first point in "i" and copy from there, but we'll
343 !   be lazy and do all the points as if they are all, independent
345 !   At this point grid%p_top is already set. find the DRY mass in the column 
346 !   by interpolating the DRY pressure.  
348     pd_surf = interp_0( pd_in, zk, grid%phb(i,1,j)/g, nl_in )
350 !   compute the perturbation mass and the full mass
352     grid%mu_1(i,j) = pd_surf-grid%p_top - grid%mub(i,j)
353     grid%mu_2(i,j) = grid%mu_1(i,j)
354     grid%mu0(i,j) = grid%mu_1(i,j) + grid%mub(i,j)
356 !   given the dry pressure and coordinate system, interp the potential
357 !   temperature and qv
359     do k=1,kde-1
361       p_level = grid%znu(k)*(pd_surf - grid%p_top) + grid%p_top
363       grid%moist(i,k,j,P_QV) = interp_0( qv, pd_in, p_level, nl_in )
364       grid%t_1(i,k,j)          = interp_0( theta, pd_in, p_level, nl_in ) - t0
365       grid%t_2(i,k,j)          = grid%t_1(i,k,j)
366       
368     enddo
370 !   integrate the hydrostatic equation (from the RHS of the bigstep
371 !   vertical momentum equation) down from the top to get grid%p.
372 !   first from the top of the model to the top pressure
374     k = kte-1  ! top level
376     qvf1 = 0.5*(grid%moist(i,k,j,P_QV)+grid%moist(i,k,j,P_QV))
377     qvf2 = 1./(1.+qvf1)
378     qvf1 = qvf1*qvf2
380 !    grid%p(i,k,j) = - 0.5*grid%mu_1(i,j)/grid%rdnw(k)
381     grid%p(i,k,j) = - 0.5*(grid%mu_1(i,j)+qvf1*grid%mub(i,j))/grid%rdnw(k)/qvf2
382     qvf = 1. + rvovrd*grid%moist(i,k,j,P_QV)
383     grid%alt(i,k,j) = (r_d/p1000mb)*(grid%t_1(i,k,j)+t0)*qvf* &
384                 (((grid%p(i,k,j)+grid%pb(i,k,j))/p1000mb)**cvpm)
385     grid%al(i,k,j) = grid%alt(i,k,j) - grid%alb(i,k,j)
387 !  down the column
389     do k=kte-2,1,-1
390       qvf1 = 0.5*(grid%moist(i,k,j,P_QV)+grid%moist(i,k+1,j,P_QV))
391       qvf2 = 1./(1.+qvf1)
392       qvf1 = qvf1*qvf2
393       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)
394       qvf = 1. + rvovrd*grid%moist(i,k,j,P_QV)
395       grid%alt(i,k,j) = (r_d/p1000mb)*(grid%t_1(i,k,j)+t0)*qvf* &
396                   (((grid%p(i,k,j)+grid%pb(i,k,j))/p1000mb)**cvpm)
397       grid%al(i,k,j) = grid%alt(i,k,j) - grid%alb(i,k,j)
398     enddo
400 !  this is the hydrostatic equation used in the model after the
401 !  small timesteps.  In the model, grid%al (inverse density)
402 !  is computed from the geopotential.
405     grid%ph_1(i,1,j) = 0.
406     DO k  = 2,kte
407       grid%ph_1(i,k,j) = grid%ph_1(i,k-1,j) - (1./grid%rdnw(k-1))*(       &
408                    (grid%mub(i,j)+grid%mu_1(i,j))*grid%al(i,k-1,j)+ &
409                     grid%mu_1(i,j)*grid%alb(i,k-1,j)  )
410                                                    
411       grid%ph_2(i,k,j) = grid%ph_1(i,k,j) 
412       grid%ph0(i,k,j) = grid%ph_1(i,k,j) + grid%phb(i,k,j)
413     ENDDO
415 ! interp u
417     DO K = 1, kte
418       p_level = grid%znu(k)*(p_surf - grid%p_top) + grid%p_top
419       grid%u_1(i,k,j) = interp_0( u, p_in, p_level, nl_in )
420       grid%u_2(i,k,j) = grid%u_1(i,k,j)
421     ENDDO
423   ENDDO
424   ENDDO
426 !  thermal perturbation to kick off convection
428   write(6,*) ' nxc, nyc for perturbation ',nxc,nyc
429   write(6,*) ' delt for perturbation ',tpbub
431   DO J = jts, min(jde-1,jte)
432     yrad = config_flags%dy*float(j-jde/2-1)/radbub
433     DO I = its, min(ide-1,ite)
434       xrad = float(i-1)/float(ide-ids)
436       DO K = 1, kte-1
438 !  put in preturbation theta (bubble) and recalc density.  note,
439 !  the mass in the column is not changing, so when theta changes,
440 !  we recompute density and geopotential
442         zrad = 0.5*(grid%ph_1(i,k,j)+grid%ph_1(i,k+1,j)  &
443                    +grid%phb(i,k,j)+grid%phb(i,k+1,j))/g
444         zrad = (zrad-htbub)/radz
445         RAD=SQRT(yrad*yrad+zrad*zrad)
446         IF(RAD <= 1.) THEN
447            tp = tpbub*cos(rad*piov2)*cos(rad*piov2)*cos(xrad*2*pi+pi)
448            grid%t_1(i,k,j)=grid%t_1(i,k,j)+tp
449            grid%t_2(i,k,j)=grid%t_1(i,k,j)
450            qvf = 1. + rvovrd*grid%moist(i,k,j,P_QV)
451            grid%alt(i,k,j) = (r_d/p1000mb)*(grid%t_1(i,k,j)+t0)*qvf* &
452                         (((grid%p(i,k,j)+grid%pb(i,k,j))/p1000mb)**cvpm)
453            grid%al(i,k,j) = grid%alt(i,k,j) - grid%alb(i,k,j)
454         ENDIF
455       ENDDO
457 !  rebalance hydrostatically
459       DO k  = 2,kte
460         grid%ph_1(i,k,j) = grid%ph_1(i,k-1,j) - (1./grid%rdnw(k-1))*(       &
461                      (grid%mub(i,j)+grid%mu_1(i,j))*grid%al(i,k-1,j)+ &
462                       grid%mu_1(i,j)*grid%alb(i,k-1,j)  )
463                                                    
464         grid%ph_2(i,k,j) = grid%ph_1(i,k,j) 
465         grid%ph0(i,k,j) = grid%ph_1(i,k,j) + grid%phb(i,k,j)
466       ENDDO
468     ENDDO
469   ENDDO
471 !#endif
473    write(6,*) ' grid%mu_1 from comp ', grid%mu_1(1,1)
474    write(6,*) ' pert state sounding from comp, grid%ph_1, pp, grid%al, grid%t_1, qv '
475    do k=1,kde-1
476      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), &
477                      grid%t_1(1,k,1), grid%moist(1,k,1,P_QV)
478    enddo
480    write(6,*) ' grid%mu_1 from comp ', grid%mu_1(1,1)
481    write(6,*) ' full state sounding from comp, ph, grid%p, grid%al, grid%t_1, qv '
482    write(6,*) ' at j = 1 '
483    do k=1,kde-1
484      write(6,'(i3,1x,5(1x,1pe10.3))') k, grid%ph_1(1,k,1)+grid%phb(1,k,1), &
485                      grid%p(1,k,1)+grid%pb(1,k,1), grid%al(1,k,1)+grid%alb(1,k,1), &
486                      grid%t_1(1,k,1)+t0, grid%moist(1,k,1,P_QV)
487    enddo
490    write(6,*) ' full state sounding from comp, ph, grid%p, grid%al, grid%t_1, qv '
491    write(6,*) ' at j = jde/2 '
492    do k=1,kde-1
493      write(6,'(i3,1x,5(1x,1pe10.3))') k, grid%ph_1(1,k,jde/2)+grid%phb(1,k,jde/2), &
494                      grid%p(1,k,jde/2)+grid%pb(1,k,jde/2), grid%al(1,k,jde/2)+grid%alb(1,k,jde/2), &
495                      grid%t_1(1,k,jde/2)+t0, grid%moist(1,k,jde/2,P_QV)
496    enddo
498    write(6,*) ' full state sounding from comp, ph, grid%p, grid%al, grid%t_1, qv '
499    write(6,*) ' at j = jde-1 '
500    do k=1,kde-1
501      write(6,'(i3,1x,5(1x,1pe10.3))') k, grid%ph_1(1,k,jde-1)+grid%phb(1,k,jde-1), &
502                      grid%p(1,k,jde-1)+grid%pb(1,k,jde-1), grid%al(1,k,jde-1)+grid%alb(1,k,jde-1), &
503                      grid%t_1(1,k,jde-1)+t0, grid%moist(1,k,jde-1,P_QV)
504    enddo
506 ! set v
508   DO J = jts, jte
509   DO I = its, min(ide-1,ite)
511     DO K = 1, kte
512       grid%v_1(i,k,j) = 0.
513       grid%v_2(i,k,j) = grid%v_1(i,k,j)
514     ENDDO
516   ENDDO
517   ENDDO
519 !  fill out last i row for u
521   DO J = jts, min(jde-1,jte)
522   DO I = ite, ite
524     DO K = 1, kte
525       grid%u_1(i,k,j) = grid%u_1(its,k,j)
526       grid%u_2(i,k,j) = grid%u_2(its,k,j)
527     ENDDO
529   ENDDO
530   ENDDO
532 !  set w
534   DO J = jts, min(jde-1,jte)
535   DO K = kts, kte
536   DO I = its, min(ide-1,ite)
537     grid%w_1(i,k,j) = 0.
538     grid%w_2(i,k,j) = 0.
539   ENDDO
540   ENDDO
541   ENDDO
543 !  set a few more things
545   DO J = jts, min(jde-1,jte)
546   DO K = kts, kte-1
547   DO I = its, min(ide-1,ite)
548     grid%h_diabatic(i,k,j) = 0.
549   ENDDO
550   ENDDO
551   ENDDO
553   DO k=1,kte-1
554     grid%t_base(k) = grid%t_1(1,k,1)
555     grid%qv_base(k) = grid%moist(1,k,1,P_QV)
556     grid%u_base(k) = grid%u_1(1,k,1)
557     grid%v_base(k) = grid%v_1(1,k,1)
558   ENDDO
560   DO J = jts, min(jde-1,jte)
561   DO I = its, min(ide-1,ite)
562      thtmp   = grid%t_2(i,1,j)+t0
563      ptmp    = grid%p(i,1,j)+grid%pb(i,1,j)
564      temp(1) = thtmp * (ptmp/p1000mb)**rcp
565      thtmp   = grid%t_2(i,2,j)+t0
566      ptmp    = grid%p(i,2,j)+grid%pb(i,2,j)
567      temp(2) = thtmp * (ptmp/p1000mb)**rcp
568      thtmp   = grid%t_2(i,3,j)+t0
569      ptmp    = grid%p(i,3,j)+grid%pb(i,3,j)
570      temp(3) = thtmp * (ptmp/p1000mb)**rcp
572      grid%tsk(I,J)=grid%cf1*temp(1)+grid%cf2*temp(2)+grid%cf3*temp(3)
573      if (i .eq. 1) print*,'sfctem',j,temp(1),temp(2),temp(3),grid%tsk(I,J)
574      grid%tmn(I,J)=grid%tsk(I,J)-0.5
575   ENDDO
576   ENDDO
578   RETURN
580  END SUBROUTINE init_domain_rk
582 !---------------------------------------------------------------------
584  SUBROUTINE init_module_initialize
585  END SUBROUTINE init_module_initialize
587 !---------------------------------------------------------------------
588 #if 0
589 ! TEST DRIVER FOR "read_input_jet" and "get_sounding"
590   implicit none 
591   integer, parameter :: nz_jet=64, ny_jet=80
592   real, dimension(nz_jet,ny_jet) :: u_jet, rho_jet, &
593                                     th_jet, z_jet
595   real, dimension(nz_jet,ny_jet) :: zk,p,p_dry,theta,rho,u,v,qv
596   logical :: dry, debug
597   integer :: j, nl
599   call read_input_jet( u_jet, rho_jet, th_jet, z_jet, nz_jet, ny_jet )
601   call opngks
602   call parray( u_jet, nz_jet, ny_jet)
603   call parray( rho_jet, nz_jet, ny_jet)
604   call parray( th_jet, nz_jet, ny_jet)
605 !  call clsgks
607 !  set up initial jet
609   debug = .true.
610   dry = .true.
611   do j=1,ny_jet
613     call get_sounding( zk(:,j),p(:,j),p_dry(:,j),theta(:,j),      &
614                        rho(:,j),u(:,j), v(:,j),  qv(:,j),        &
615                        dry, nz_jet, nl, u_jet, rho_jet, th_jet,  &
616                        z_jet, nz_jet, ny_jet, j, debug          )
617     debug = .false.
619   enddo
621   write(6,*) ' lowest level p, th, and rho, highest level p '
623   do j=1,ny_jet
624     write(6,*) j, p(1,j),theta(1,j),rho(1,j), p(nz_jet,j)
625 !    write(6,*) j, p(1,j),theta(1,j)-th_jet(1,j),rho(1,j)-rho_jet(1,j)
626   enddo
628   call parray( p, nz_jet, ny_jet)
629   call parray( p_dry, nz_jet, ny_jet)
630   call parray( theta, nz_jet, ny_jet)
632   call clsgks
634   end
636 !---------------------------------
638       subroutine parray(a,m,n)
639       dimension a(m,n)
640       dimension b(n,m)
642     do i=1,m
643     do j=1,n
644       b(j,i) = a(i,j)
645     enddo
646     enddo
647       
648       write(6,'(''  dimensions m,n  '',2i6)')m,n
649         call set(.05,.95,.05,.95,0.,1.,0.,1.,1)
650         call perim(4,5,4,5)
651         call setusv('LW',2000)
652 !        CALL CONREC(a,m,m,n,cmax,cmin,cinc,-1,-638,-922)
653         CALL CONREC(b,n,n,m,0.,0.,0.,-1,-638,-922)
654         call frame
655       return
656       end
658 ! END TEST DRIVER FOR "read_input_jet" and "get_sounding"
659 #endif
661 !------------------------------------------------------------------
663     subroutine get_sounding( zk, p, p_dry, theta, rho,       &
664                              u, v, qv, dry, nl_max, nl_in,  &
665                              u_jet, rho_jet, th_jet, z_jet, &
666                              nz_jet, ny_jet, j_point, debug )
667     implicit none
669     integer nl_max, nl_in
670     real zk(nl_max), p(nl_max), theta(nl_max), rho(nl_max), &
671          u(nl_max), v(nl_max), qv(nl_max), p_dry(nl_max)
672     logical dry
674     integer nz_jet, ny_jet, j_point
675     real, dimension(nz_jet, ny_jet) :: u_jet, rho_jet, th_jet, z_jet
677     integer n
678     parameter(n=1000)
679     logical debug
681 ! input sounding data
683     real p_surf, th_surf, qv_surf
684     real pi_surf, pi(n)
685     real h_input(n), th_input(n), qv_input(n), u_input(n), v_input(n)
687 ! diagnostics
689     real rho_surf, p_input(n), rho_input(n)
690     real pm_input(n)  !  this are for full moist sounding
692 ! local data
694     real r
695     parameter (r = r_d)
696     integer k, it, nl
697     real qvf, qvf1, dz
699 !  first, read the sounding
701 !    call read_sounding( p_surf, th_surf, qv_surf, &
702 !                          h_input, th_input, qv_input, u_input, v_input,n, nl, debug )
704    call calc_jet_sounding( p_surf, th_surf, qv_surf,                             &
705                            h_input, th_input, qv_input, u_input, v_input,        &
706                            n, nl, debug, u_jet, rho_jet, th_jet, z_jet, j_point, &
707                            nz_jet, ny_jet, dry                                  )
709    nl = nz_jet
711     if(dry) then
712      do k=1,nl
713        qv_input(k) = 0.
714      enddo
715     endif
717     if(debug) write(6,*) ' number of input levels = ',nl
719       nl_in = nl
720       if(nl_in .gt. nl_max ) then
721         write(6,*) ' too many levels for input arrays ',nl_in,nl_max
722         call wrf_error_fatal ( ' too many levels for input arrays ' )
723       end if
725 !  compute diagnostics,
726 !  first, convert qv(g/kg) to qv(g/g)
728 !      do k=1,nl
729 !        qv_input(k) = 0.001*qv_input(k)
730 !      enddo
731 !      p_surf = 100.*p_surf  ! convert to pascals
733     qvf = 1. + rvovrd*qv_input(1) 
734     rho_surf = 1./((r/p1000mb)*th_surf*qvf*((p_surf/p1000mb)**cvpm))
735     pi_surf = (p_surf/p1000mb)**(r/cp)
737     if(debug) then
738       write(6,*) ' surface density is ',rho_surf
739       write(6,*) ' surface pi is    ',pi_surf
740     end if
743 !  integrate moist sounding hydrostatically, starting from the
744 !  specified surface pressure
745 !  -> first, integrate from surface to lowest level
747         qvf = 1. + rvovrd*qv_input(1) 
748         qvf1 = 1. + qv_input(1)
749         rho_input(1) = rho_surf
750         dz = h_input(1)
751         do it=1,10
752           pm_input(1) = p_surf &
753                   - 0.5*dz*(rho_surf+rho_input(1))*g*qvf1
754           rho_input(1) = 1./((r/p1000mb)*th_input(1)*qvf*((pm_input(1)/p1000mb)**cvpm))
755         enddo
757 ! integrate up the column
759         do k=2,nl
760           rho_input(k) = rho_input(k-1)
761           dz = h_input(k)-h_input(k-1)
762           qvf1 = 0.5*(2.+(qv_input(k-1)+qv_input(k)))
763           qvf = 1. + rvovrd*qv_input(k)   ! qv is in g/kg here
765           do it=1,10
766             pm_input(k) = pm_input(k-1) &
767                     - 0.5*dz*(rho_input(k)+rho_input(k-1))*g*qvf1
768             rho_input(k) = 1./((r/p1000mb)*th_input(k)*qvf*((pm_input(k)/p1000mb)**cvpm))
769           enddo
770         enddo
772 !  we have the moist sounding
774 !  next, compute the dry sounding using p at the highest level from the
775 !  moist sounding and integrating down.
777         p_input(nl) = pm_input(nl)
779           do k=nl-1,1,-1
780             dz = h_input(k+1)-h_input(k)
781             p_input(k) = p_input(k+1) + 0.5*dz*(rho_input(k)+rho_input(k+1))*g
782           enddo
785         do k=1,nl
787           zk(k) = h_input(k)
788           p(k) = pm_input(k)
789           p_dry(k) = p_input(k)
790           theta(k) = th_input(k)
791           rho(k) = rho_input(k)
792           u(k) = u_input(k)
793           v(k) = v_input(k)
794           qv(k) = qv_input(k)
796         enddo
798      if(debug) then
799       write(6,*) ' sounding '
800       write(6,*) '  k  height(m)  press (Pa)   pd(Pa)   theta (K)  den(kg/m^3)  u(m/s)     v(m/s)    qv(g/g) '
801       do k=1,nl
802         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)
803       enddo
805      end if
807      end subroutine get_sounding
809 !------------------------------------------------------------------
811   subroutine calc_jet_sounding( p_surf, th_surf, qv_surf,      &
812                                 h, th, qv, u, v, n, nl, debug, &
813                                 u_jet, rho_jet, th_jet, z_jet, &
814                                 jp, nz_jet, ny_jet, dry       )
815   implicit none
816   integer :: n, nl, jp, nz_jet, ny_jet
818   real, dimension(nz_jet, ny_jet) :: u_jet, rho_jet, th_jet, z_jet
819   real, dimension(n) :: h,th,qv,u,v
820   real :: p_surf, th_surf, qv_surf
821   logical :: debug, dry
823   real, dimension(1:nz_jet) :: rho, rel_hum, p
824   integer :: k
826 !  some local stuff
828   real :: tmppi, es, qvs, temperature
830 !  get sounding from column jp
832    do k=1,nz_jet
833      h(k)  = z_jet(k,jp)
834      th(k) = th_jet(k,jp)
835      qv(k) = 0.
836      rho(k) = rho_jet(k,jp)
837      u(k) = u_jet(k,jp)
838      v(k) = 0.
839    enddo
841    if (.not.dry) then
842      DO k=1,nz_jet
843        if(h(k) .gt. 8000.) then
844          rel_hum(k)=0.1
845        else
846          rel_hum(k)=(1.-0.90*(h(k)/8000.)**1.25)
847        end if
848        rel_hum(k) = min(0.7,rel_hum(k))
849      ENDDO
850    else
851      do k=1,nz_jet
852        rel_hum(k) = 0.
853      enddo
854    endif
856 !  next, compute pressure
858    do k=1,nz_jet
859      p(k) = p1000mb*(R_d*rho(k)*th(k)/p1000mb)**cpovcv
860    enddo
862 !  here we adjust for fixed moisture profile
864      IF (.not.dry)  THEN
866 !  here we assume the input theta is th_v, so we reset theta accordingly
868        DO k=1,nz_jet
869          tmppi=(p(k)/p1000mb)**rcp
870          temperature = tmppi*th(k)
871          if (temperature .gt. svpt0) then
872             es  = 1000.*svp1*exp(svp2*(temperature-svpt0)/(temperature-svp3))
873             qvs = ep_2*es/(p(k)-es)
874          else
875             es  = 1000.*svp1*exp( 21.8745584*(temperature-273.16)/(temperature-7.66) )
876             qvs = ep_2*es/(p(k)-es)
877          endif
878          qv(k) = rel_hum(k)*qvs
879          th(k) = th(k)/(1.+.61*qv(k))
880        ENDDO
882      ENDIF
884 !  finally, set the surface data. We'll just do a simple extrapolation
886    p_surf = 1.5*p(1) - 0.5*p(2)
887    th_surf = 1.5*th(1) - 0.5*th(2)
888    qv_surf = 1.5*qv(1) - 0.5*qv(2)
890    end subroutine calc_jet_sounding
892 !---------------------------------------------------------------------
894  SUBROUTINE read_input_jet( u, r, t, zk, nz, ny )
895  implicit none
897  integer, intent(in) :: nz,ny
898  real, dimension(nz,ny), intent(out) :: u,r,t,zk
899  integer :: ny_in, nz_in, j,k
900  real, dimension(ny,nz) :: field_in
902 ! this code assumes it is called on processor 0 only
904    OPEN(unit=10, file='input_jet', form='unformatted', status='old' )
905    REWIND(10) 
906    read(10) ny_in,nz_in
907    if((ny_in /= ny ) .or. (nz_in /= nz)) then
908      write(0,*) ' error in input jet dimensions '
909      write(0,*) ' ny, ny_input, nz, nz_input ', ny, ny_in, nz,nz_in
910      write(0,*) ' error exit '
911      call wrf_error_fatal ( ' error in input jet dimensions ' )
912    end if
913    read(10) field_in
914    do j=1,ny
915    do k=1,nz
916      u(k,j) = field_in(j,k)
917    enddo
918    enddo
919    read(10) field_in
920    do j=1,ny
921    do k=1,nz
922      t(k,j) = field_in(j,k)
923    enddo
924    enddo
926    read(10) field_in
927    do j=1,ny
928    do k=1,nz
929      r(k,j) = field_in(j,k)
930    enddo
931    enddo
933    do j=1,ny
934    do k=1,nz
935      zk(k,j) = 125. + 250.*float(k-1)
936    enddo
937    enddo
939  end subroutine read_input_jet
941 END MODULE module_initialize_ideal