r4627 | gill | 2010-12-29 16:29:58 -0700 (Wed, 29 Dec 2010) | 5 lines
[wrffire.git] / wrfv2_fire / dyn_em / module_initialize_scm_xy.F
blob5e66a349b9e28635a1a1a39b1ad7f27b6538dc95
1 !IDEAL:MODEL_LAYER:INITIALIZATION
4 !  This MODULE holds the routines which are used to perform various initializations
5 !  for the individual domains.  
7 !  This MODULE CONTAINS the following routines:
9 !  initialize_field_test - 1. Set different fields to different constant
10 !                             values.  This is only a test.  If the correct
11 !                             domain is not found (based upon the "id")
12 !                             then a fatal error is issued.               
14 !-----------------------------------------------------------------------
16 MODULE module_initialize_ideal
18    USE module_domain
19    USE module_io_domain
20    USE module_state_description
21    USE module_model_constants
22    USE module_bc
23    USE module_timing
24    USE module_configure
25    USE module_init_utilities
26    USE module_soil_pre
27 #ifdef DM_PARALLEL
28    USE module_dm
29 #endif
32 CONTAINS
35 !-------------------------------------------------------------------
36 ! this is a wrapper for the solver-specific init_domain routines.
37 ! Also dereferences the grid variables and passes them down as arguments.
38 ! This is crucial, since the lower level routines may do message passing
39 ! and this will get fouled up on machines that insist on passing down
40 ! copies of assumed-shape arrays (by passing down as arguments, the 
41 ! data are treated as assumed-size -- ie. f77 -- arrays and the copying
42 ! business is avoided).  Fie on the F90 designers.  Fie and a pox.
43 ! NOTE:  Modified to remove all but arrays of rank 4 or more from the 
44 !        argument list.  Arrays with rank>3 are still problematic due to the 
45 !        above-noted fie- and pox-ities.  TBH 20061129.  
47    SUBROUTINE init_domain ( grid )
49    IMPLICIT NONE
51    !  Input data.
52    TYPE (domain), POINTER :: grid 
53    !  Local data.
54    INTEGER :: idum1, idum2
56    CALL set_scalar_indices_from_config ( head_grid%id , idum1, idum2 )
58      CALL init_domain_rk( grid &
60 #include <actual_new_args.inc>
62                         )
63    END SUBROUTINE init_domain
65 !-------------------------------------------------------------------
67    SUBROUTINE init_domain_rk ( grid &
69 # include <dummy_new_args.inc>
73    USE module_optional_input
74    IMPLICIT NONE
76    !  Input data.
77    TYPE (domain), POINTER :: grid
79 # include <dummy_new_decl.inc>
81    TYPE (grid_config_rec_type)              :: config_flags
83    !  Local data
84    INTEGER                             ::                       &
85                                   ids, ide, jds, jde, kds, kde, &
86                                   ims, ime, jms, jme, kms, kme, &
87                                   its, ite, jts, jte, kts, kte, &
88                                   i, j, k
90 ! JPH should add a read to a config file with:
91 ! ----- check to make sure everything is initialized from the LU index, etc.
92 ! ----- need to make a dummy category?
93 ! cen_lat, cen_lon
94 ! land-use category
95 ! soil category
97    ! Local data
98    INTEGER, PARAMETER :: nl_max = 1000
99    REAL, DIMENSION(nl_max) :: zk, p_in, theta, rho, u, v, qv, pd_in
100    INTEGER :: nl_in
102    INTEGER :: ii, im1, jj, jm1, loop, error, fid, lm
103    REAL    :: u_mean,v_mean, f0, p_surf, p_level, qvf, z_at_v, z_at_u
104    REAL    :: xrad, yrad, zrad, rad, cof1, cof2
105    REAL    :: pi, rnd
107 !  stuff from original initialization that has been dropped from the Registry 
108    REAL    :: vnu, xnu, xnus, dinit0, cbh, p0_temp, t0_temp, zd
109    REAL    :: qvf1, qvf2, pd_surf
110    INTEGER :: it
111    real :: thtmp, ptmp, temp(3)
112    real :: zsfc
114    LOGICAL :: moisture_init
115    LOGICAL :: dry_sounding
116    character (len=256) :: mminlu2
118 ! soil input
119    INTEGER :: ns_input
120    REAL    :: tmn_input, tsk_input
121    REAL    :: zs_input(100),tslb_input(100),smois_input(100)
122    LOGICAL :: real_soil = .true.
124    REAL    :: zrwa(200), zwa(200)
126    
127 #ifdef DM_PARALLEL
128 #    include <data_calls.inc>
129 #endif
132    SELECT CASE ( model_data_order )
133          CASE ( DATA_ORDER_ZXY )
134    kds = grid%sd31 ; kde = grid%ed31 ;
135    ids = grid%sd32 ; ide = grid%ed32 ;
136    jds = grid%sd33 ; jde = grid%ed33 ;
138    kms = grid%sm31 ; kme = grid%em31 ;
139    ims = grid%sm32 ; ime = grid%em32 ;
140    jms = grid%sm33 ; jme = grid%em33 ;
142    kts = grid%sp31 ; kte = grid%ep31 ;   ! note that tile is entire patch
143    its = grid%sp32 ; ite = grid%ep32 ;   ! note that tile is entire patch
144    jts = grid%sp33 ; jte = grid%ep33 ;   ! note that tile is entire patch
145          CASE ( DATA_ORDER_XYZ )
146    ids = grid%sd31 ; ide = grid%ed31 ;
147    jds = grid%sd32 ; jde = grid%ed32 ;
148    kds = grid%sd33 ; kde = grid%ed33 ;
150    ims = grid%sm31 ; ime = grid%em31 ;
151    jms = grid%sm32 ; jme = grid%em32 ;
152    kms = grid%sm33 ; kme = grid%em33 ;
154    its = grid%sp31 ; ite = grid%ep31 ;   ! note that tile is entire patch
155    jts = grid%sp32 ; jte = grid%ep32 ;   ! note that tile is entire patch
156    kts = grid%sp33 ; kte = grid%ep33 ;   ! note that tile is entire patch
157          CASE ( DATA_ORDER_XZY )
158    ids = grid%sd31 ; ide = grid%ed31 ;
159    kds = grid%sd32 ; kde = grid%ed32 ;
160    jds = grid%sd33 ; jde = grid%ed33 ;
162    ims = grid%sm31 ; ime = grid%em31 ;
163    kms = grid%sm32 ; kme = grid%em32 ;
164    jms = grid%sm33 ; jme = grid%em33 ;
166    its = grid%sp31 ; ite = grid%ep31 ;   ! note that tile is entire patch
167    kts = grid%sp32 ; kte = grid%ep32 ;   ! note that tile is entire patch
168    jts = grid%sp33 ; jte = grid%ep33 ;   ! note that tile is entire patch
170    END SELECT
173    pi = 2.*asin(1.0)
174    write(6,*) ' pi is ',pi
176    CALL model_to_grid_config_rec ( grid%id , model_config_rec , config_flags )
178 ! here we check to see if the boundary conditions are set properly
180    CALL boundary_condition_check( config_flags, bdyzone, error, grid%id )
182    moisture_init = .true.
184     grid%itimestep=0
186     mminlu2 = ' '
187     mminlu2(1:4) = 'USGS'
188     CALL nl_set_mminlu(1, mminlu2)
189 !   CALL nl_set_mminlu(1, 'USGS')
190     CALL nl_set_iswater(1,16)
191     CALL nl_set_isice(1,3)
192     CALL nl_set_truelat1(1,0.)
193     CALL nl_set_truelat2(1,0.)
194     CALL nl_set_moad_cen_lat (1,0.)
195     CALL nl_set_stand_lon(1,0.)
196     CALL nl_set_pole_lon (1,0.)
197     CALL nl_set_pole_lat (1,90.)
198     CALL nl_set_map_proj(1,0)
199 !   CALL model_to_grid_config_rec(1,model_config_rec,config_flags)
200     CALL nl_get_iswater(1,grid%iswater)
202 !  here we initialize data that currently is not initialized 
203 !  in the input data
205     DO j = jts, jte
206       DO i = its, ite
207          grid%msft(i,j)     = 1.
208          grid%msfu(i,j)     = 1.
209          grid%msfv(i,j)     = 1.
210          grid%msftx(i,j)    = 1.
211          grid%msfty(i,j)    = 1.
212          grid%msfux(i,j)    = 1.
213          grid%msfuy(i,j)    = 1.
214          grid%msfvx(i,j)    = 1.
215          grid%msfvx_inv(i,j)    = 1.
216          grid%msfvy(i,j)    = 1.
217          grid%sina(i,j)     = 0.
218          grid%cosa(i,j)     = 1.
219          grid%e(i,j)        = 2.0*EOMEG*cos(config_flags%scm_lat*DEGRAD)
220          grid%f(i,j)        = 2.0*EOMEG*sin(config_flags%scm_lat*DEGRAD)
221          grid%xlat(i,j)     = config_flags%scm_lat
222          grid%xlong(i,j)    = config_flags%scm_lon
223          grid%xland(i,j)     = 1.
224          grid%landmask(i,j)  = 1.
225          grid%lu_index(i,j)  = config_flags%scm_lu_index
226       END DO
227    END DO
229 ! for Noah LSM, additional variables need to be initialized
231    other_masked_fields : SELECT CASE ( model_config_rec%sf_surface_physics(grid%id) )
233       CASE (SLABSCHEME)
235       CASE (LSMSCHEME)
237 ! JPH free of snow and ice, and only valid over land
238         DO j = jts , MIN(jde-1,jte)
239            DO i = its , MIN(ide-1,ite)
240               grid%vegfra(i,j) = config_flags%scm_vegfra
241               grid%canwat(i,j) = config_flags%scm_canwat
242               grid%isltyp(i,j)  = config_flags%scm_isltyp
243               grid%ivgtyp(i,j)  = config_flags%scm_lu_index
244               grid%xice(i,j) = 0.
245               grid%snow(i,j) = 0.
246            END DO
247         END DO
249       CASE (RUCLSMSCHEME)
251    END SELECT other_masked_fields
253    grid%step_number = 0
255    IF ( real_soil ) THEN ! from input file
256     
257       CALL read_soil(100,ns_input,tmn_input,tsk_input,zs_input,tslb_input,smois_input)
259       CALL init_module_optional_input(grid,config_flags)
260       num_st_levels_input = ns_input
261       num_sm_levels_input = ns_input
262       num_sw_levels_input = ns_input
263       DO k = 1,ns_input
264          st_levels_input(k) = zs_input(k)*100.0 ! to cm
265          sm_levels_input(k) = zs_input(k)*100.0 ! to cm
266          sw_levels_input(k) = zs_input(k)*100.0 ! to cm
267          st_input(:,k+1,:) = tslb_input(k)
268          sm_input(:,k+1,:) = smois_input(k)
269          sw_input(:,k+1,:) = smois_input(k)
270       ENDDO
271   
272       grid%tsk = tsk_input
273       grid%sst = tsk_input
274       grid%tmn = tmn_input
276       flag_soil_layers  = 0 ! go ahead and put skin temp in
277       flag_soil_levels  = 0 ! go ahead and put skin moisture in
278       flag_sst          = 0 ! don't modify for ocean
279       flag_tavgsfc      = 0 
280       flag_soilhgt      = 0 
282       CALL process_soil_real ( grid%tsk , grid%tmn , grid%tavgsfc, &
283                    grid%landmask , grid%sst , grid%ht, grid%toposoil, &
284                    st_input , sm_input , sw_input , &
285                    st_levels_input , sm_levels_input , sw_levels_input , &
286                    grid%zs , grid%dzs , grid%tslb , grid%smois , grid%sh2o , &
287                    flag_sst , flag_tavgsfc, flag_soilhgt, flag_soil_layers, flag_soil_levels,  &
288                    ids , ide , jds , jde , kds , kde , &
289                    ims , ime , jms , jme , kms , kme , &
290                    its , ite , jts , jte , kts , kte , &
291                    model_config_rec%sf_surface_physics(grid%id) , &
292                    model_config_rec%num_soil_layers , &
293                    model_config_rec%real_data_init_type , &
294                    num_st_levels_input , num_sm_levels_input , num_sw_levels_input , &
295                    num_st_levels_alloc , num_sm_levels_alloc , num_sw_levels_alloc )
297    ELSE ! ideal soil
298 ! Process the soil; note that there are some things hard-wired into share/module_soil_pre.F
299       CALL process_soil_ideal(grid%xland,grid%xice,grid%vegfra,grid%snow,grid%canwat, &
300                      grid%ivgtyp,grid%isltyp,grid%tslb,grid%smois, &
301                      grid%tsk,grid%tmn,grid%zs,grid%dzs,model_config_rec%num_soil_layers, &
302                      model_config_rec%sf_surface_physics(grid%id), &
303                                    ids,ide, jds,jde, kds,kde,&
304                                    ims,ime, jms,jme, kms,kme,&
305                                    its,ite, jts,jte, kts,kte )
307     ENDIF
309     DO j = jts, jte
310      DO k = kts, kte
311        DO i = its, ite
312           grid%ww(i,k,j)     = 0.
313        END DO
314      END DO
315     END DO
317 ! this is adopted from Wayne Angevine's GABLS case
318     grid%znw(1) = 1.0
319     zrwa(kde) = exp((kde-1)/40.)
320     zwa(kde) = grid%ztop
321     DO k=2, kde-1
322        zrwa(k) = exp((k-1)/40.)
323        zwa(k) = (zrwa(k)-1.) * grid%ztop/(zrwa(kde)-1.)
324        grid%znw(k) = 1. - (zwa(k) / zwa(kde))
325     ENDDO
326     grid%znw(kde) = 0.
328    DO k=1, kde-1
329     grid%dnw(k) = grid%znw(k+1) - grid%znw(k)
330     grid%rdnw(k) = 1./grid%dnw(k)
331     grid%znu(k) = 0.5*(grid%znw(k+1)+grid%znw(k))
332    ENDDO
333    DO k=2, kde-1
334     grid%dn(k) = 0.5*(grid%dnw(k)+grid%dnw(k-1))
335     grid%rdn(k) = 1./grid%dn(k)
336     grid%fnp(k) = .5* grid%dnw(k  )/grid%dn(k)
337     grid%fnm(k) = .5* grid%dnw(k-1)/grid%dn(k)
338    ENDDO
340    cof1 = (2.*grid%dn(2)+grid%dn(3))/(grid%dn(2)+grid%dn(3))*grid%dnw(1)/grid%dn(2) 
341    cof2 =     grid%dn(2)        /(grid%dn(2)+grid%dn(3))*grid%dnw(1)/grid%dn(3) 
342    grid%cf1  = grid%fnp(2) + cof1
343    grid%cf2  = grid%fnm(2) - cof1 - cof2
344    grid%cf3  = cof2       
346    grid%cfn  = (.5*grid%dnw(kde-1)+grid%dn(kde-1))/grid%dn(kde-1)
347    grid%cfn1 = -.5*grid%dnw(kde-1)/grid%dn(kde-1)
348    grid%rdx = 1./config_flags%dx
349    grid%rdy = 1./config_flags%dy
351 !  get the sounding from the ascii sounding file, first get dry sounding and 
352 !  calculate base state
354   write(6,*) ' getting dry sounding for base state '
355   dry_sounding = .true.
356   CALL get_sounding( zsfc, zk, p_in, pd_in, theta, rho, u, v, qv, dry_sounding, nl_max, nl_in )
358   write(6,*) ' returned from reading sounding, nl_in is ',nl_in
360 !  find ptop for the desired ztop (ztop is input from the namelist),
361 !  and find surface pressure
363   grid%p_top = interp_0( p_in, zk, config_flags%ztop, nl_in )
365   DO j=jts,jte
366   DO i=its,ite  ! flat surface
367     grid%ht(i,j) = zsfc
368     grid%phb(i,1,j) = grid%ht(i,j) * g
369     grid%ph0(i,1,j) = grid%ht(i,j) * g
370     grid%php(i,1,j) = 0.
371   ENDDO
372   ENDDO
374   DO J = jts, jte
375   DO I = its, ite
377     p_surf = interp_0( p_in, zk, grid%phb(i,1,j)/g, nl_in )
378     grid%mub(i,j) = p_surf-grid%p_top
380 !  this is dry hydrostatic sounding (base state), so given grid%p (coordinate),
381 !  interp theta (from interp) and compute 1/rho from eqn. of state
383     DO K = 1, kte-1
384       p_level = grid%znu(k)*(p_surf - grid%p_top) + grid%p_top
385       grid%pb(i,k,j) = p_level
386       grid%t_init(i,k,j) = interp_0( theta, p_in, p_level, nl_in ) - t0
387       grid%alb(i,k,j) = (r_d/p1000mb)*(grid%t_init(i,k,j)+t0)*(grid%pb(i,k,j)/p1000mb)**cvpm
388     ENDDO
390 !  calc hydrostatic balance (alternatively we could interp the geopotential from the
391 !  sounding, but this assures that the base state is in exact hydrostatic balance with
392 !  respect to the model eqns.
394     DO k  = 2,kte
395       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)
396     ENDDO
398   ENDDO
399   ENDDO
401   write(6,*) ' ptop is ',grid%p_top
402   write(6,*) ' base state grid%mub(1,1), p_surf is ',grid%mub(1,1),grid%mub(1,1)+grid%p_top
404 !  calculate full state for each column - this includes moisture.
406   write(6,*) ' getting moist sounding for full state '
407   dry_sounding = .false.
408   CALL get_sounding( zsfc, zk, p_in, pd_in, theta, rho, u, v, qv, dry_sounding, nl_max, nl_in )
410   DO J = jts, min(jde-1,jte)
411   DO I = its, min(ide-1,ite)
413 !  At this point grid%p_top is already set. find the DRY mass in the column 
414 !  by interpolating the DRY pressure.  
416    pd_surf = interp_0( pd_in, zk, grid%phb(i,1,j)/g, nl_in )
418 !  compute the perturbation mass and the full mass
420     grid%mu_1(i,j) = pd_surf-grid%p_top - grid%mub(i,j)
421     grid%mu_2(i,j) = grid%mu_1(i,j)
422     grid%mu0(i,j) = grid%mu_1(i,j) + grid%mub(i,j)
424 ! given the dry pressure and coordinate system, interp the potential
425 ! temperature and qv
427     do k=1,kde-1
429       p_level = grid%znu(k)*(pd_surf - grid%p_top) + grid%p_top
431       moist(i,k,j,P_QV) = interp_0( qv, pd_in, p_level, nl_in )
432       grid%t_1(i,k,j)          = interp_0( theta, pd_in, p_level, nl_in ) - t0
433       grid%t_2(i,k,j)          = grid%t_1(i,k,j)
434       
436     enddo
438 !  integrate the hydrostatic equation (from the RHS of the bigstep
439 !  vertical momentum equation) down from the top to get grid%p.
440 !  first from the top of the model to the top pressure
442     k = kte-1  ! top level
444     qvf1 = 0.5*(moist(i,k,j,P_QV)+moist(i,k,j,P_QV))
445     qvf2 = 1./(1.+qvf1)
446     qvf1 = qvf1*qvf2
448 !    grid%p(i,k,j) = - 0.5*grid%mu_1(i,j)/grid%rdnw(k)
449     grid%p(i,k,j) = - 0.5*(grid%mu_1(i,j)+qvf1*grid%mub(i,j))/grid%rdnw(k)/qvf2
450     qvf = 1. + rvovrd*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)
455 !  down the column
457     do k=kte-2,1,-1
458       qvf1 = 0.5*(moist(i,k,j,P_QV)+moist(i,k+1,j,P_QV))
459       qvf2 = 1./(1.+qvf1)
460       qvf1 = qvf1*qvf2
461       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)
462       qvf = 1. + rvovrd*moist(i,k,j,P_QV)
463       grid%alt(i,k,j) = (r_d/p1000mb)*(grid%t_1(i,k,j)+t0)*qvf* &
464                   (((grid%p(i,k,j)+grid%pb(i,k,j))/p1000mb)**cvpm)
465       grid%al(i,k,j) = grid%alt(i,k,j) - grid%alb(i,k,j)
466     enddo
468 !  this is the hydrostatic equation used in the model after the
469 !  small timesteps.  In the model, grid%al (inverse density)
470 !  is computed from the geopotential.
473     grid%ph_1(i,1,j) = 0.
474     DO k  = 2,kte
475       grid%ph_1(i,k,j) = grid%ph_1(i,k-1,j) - (1./grid%rdnw(k-1))*(       &
476                    (grid%mub(i,j)+grid%mu_1(i,j))*grid%al(i,k-1,j)+ &
477                     grid%mu_1(i,j)*grid%alb(i,k-1,j)  )
478                                                    
479       grid%ph_2(i,k,j) = grid%ph_1(i,k,j) 
480       grid%ph0(i,k,j) = grid%ph_1(i,k,j) + grid%phb(i,k,j)
481     ENDDO
483     if((i==2) .and. (j==2)) then
484      write(6,*) ' grid%ph_1 calc ',grid%ph_1(2,1,2),grid%ph_1(2,2,2),&
485                               grid%mu_1(2,2)+grid%mub(2,2),grid%mu_1(2,2), &
486                               grid%alb(2,1,2),grid%al(1,2,1),grid%rdnw(1)
487     endif
489   ENDDO
490   ENDDO
492    write(6,*) ' grid%mu_1 from comp ', grid%mu_1(1,1)
493    write(6,*) ' full state sounding from comp, ph, grid%p, grid%al, grid%t_1, qv '
495    do k=1,kde-1
496      write(6,'(i3,1x,5(1x,1pe10.3))') k, grid%ph_1(1,k,1)+grid%phb(1,k,1), &
497                                       grid%p(1,k,1)+grid%pb(1,k,1), grid%alt(1,k,1), &
498                                       grid%t_1(1,k,1)+t0, moist(1,k,1,P_QV)
499    enddo
501    write(6,*) ' pert state sounding from comp, grid%ph_1, pp, alp, grid%t_1, qv '
502    do k=1,kde-1
503      write(6,'(i3,1x,5(1x,1pe10.3))') k, grid%ph_1(1,k,1), &
504                                       grid%p(1,k,1), grid%al(1,k,1), &
505                                       grid%t_1(1,k,1), moist(1,k,1,P_QV)
506    enddo
508 ! interp v
510   DO J = jts, jte
511   DO I = its, min(ide-1,ite)
513     IF (j == jds) THEN
514       z_at_v = grid%phb(i,1,j)/g
515     ELSE IF (j == jde) THEN
516       z_at_v = grid%phb(i,1,j-1)/g
517     ELSE
518       z_at_v = 0.5*(grid%phb(i,1,j)+grid%phb(i,1,j-1))/g
519     END IF
521     p_surf = interp_0( p_in, zk, z_at_v, nl_in )
523     DO K = 1, kte
524       p_level = grid%znu(k)*(p_surf - grid%p_top) + grid%p_top
525       grid%v_1(i,k,j) = interp_0( v, p_in, p_level, nl_in )
526       grid%v_2(i,k,j) = grid%v_1(i,k,j)
527     ENDDO
529   ENDDO
530   ENDDO
532 ! interp u
534   DO J = jts, min(jde-1,jte)
535   DO I = its, ite
537     IF (i == ids) THEN
538       z_at_u = grid%phb(i,1,j)/g
539     ELSE IF (i == ide) THEN
540       z_at_u = grid%phb(i-1,1,j)/g
541     ELSE
542       z_at_u = 0.5*(grid%phb(i,1,j)+grid%phb(i-1,1,j))/g
543     END IF
545     p_surf = interp_0( p_in, zk, z_at_u, nl_in )
547     DO K = 1, kte
548       p_level = grid%znu(k)*(p_surf - grid%p_top) + grid%p_top
549       grid%u_1(i,k,j) = interp_0( u, p_in, p_level, nl_in )
550       grid%u_2(i,k,j) = grid%u_1(i,k,j)
551     ENDDO
553   ENDDO
554   ENDDO
556 !  set w
558   DO J = jts, min(jde-1,jte)
559   DO K = kts, kte
560   DO I = its, min(ide-1,ite)
561     grid%w_1(i,k,j) = 0.
562     grid%w_2(i,k,j) = 0.
563   ENDDO
564   ENDDO
565   ENDDO
567 !  set a few more things
569   DO J = jts, min(jde-1,jte)
570   DO K = kts, kte-1
571   DO I = its, min(ide-1,ite)
572     grid%h_diabatic(i,k,j) = 0.
573   ENDDO
574   ENDDO
575   ENDDO
577 ! Go ahead and initialize these from the sounding.  This will allow a run
578 ! to actually succeed even if scm_force = 0
579   DO k=1,kte-1
580     grid%t_base(k) = grid%t_1(1,k,1)
581     grid%qv_base(k) = moist(1,k,1,P_QV)
582     grid%u_base(k) = grid%u_1(1,k,1)
583     grid%v_base(k) = grid%v_1(1,k,1)
584     grid%z_base(k) = 0.5*(grid%phb(1,k,1)+grid%phb(1,k+1,1)+grid%ph_1(1,k,1)+grid%ph_1(1,k+1,1))/g
585   ENDDO
587   RETURN
589  END SUBROUTINE init_domain_rk
591    SUBROUTINE init_module_initialize
592    END SUBROUTINE init_module_initialize
594 !---------------------------------------------------------------------
596 !  test driver for get_sounding
598 !      implicit none
599 !      integer n
600 !      parameter(n = 1000)
601 !      real zk(n),p(n),theta(n),rho(n),u(n),v(n),qv(n),pd(n)
602 !      logical dry
603 !      integer nl,k
605 !      dry = .false.
606 !      dry = .true.
607 !      call get_sounding( zk, p, pd, theta, rho, u, v, qv, dry, n, nl )
608 !      write(6,*) ' input levels ',nl
609 !      write(6,*) ' sounding '
610 !      write(6,*) '  k  height(m)  press (Pa) pd(Pa) theta (K) den(kg/m^3)  u(m/s)     v(m/s)    qv(g/g) '
611 !      do k=1,nl
612 !        write(6,'(1x,i3,8(1x,1pe10.3))') k, zk(k), p(k), pd(k), theta(k), rho(k), u(k), v(k), qv(k)
613 !      enddo
614 !      end
616 !---------------------------------------------------------------------------
618       subroutine get_sounding( zsfc, zk, p, p_dry, theta, rho, &
619                                u, v, qv, dry, nl_max, nl_in )
620       implicit none
622       integer nl_max, nl_in
623       real zsfc
624       real zk(nl_max), p(nl_max), theta(nl_max), rho(nl_max), &
625            u(nl_max), v(nl_max), qv(nl_max), p_dry(nl_max)
626       logical dry
628       integer n
629       parameter(n=3000)
630       logical debug
631       parameter( debug = .true.)
633 ! input sounding data
635       real p_surf, th_surf, qv_surf
636       real pi_surf, pi(n)
637       real h_input(n), th_input(n), qv_input(n), u_input(n), v_input(n)
639 ! diagnostics
641       real rho_surf, p_input(n), rho_input(n)
642       real pm_input(n)  !  this are for full moist sounding
644 ! local data
646       real r
647       parameter (r = r_d)
648       integer k, it, nl
649       real qvf, qvf1, dz
651 !  first, read the sounding
653       call read_sounding( zsfc, p_surf, th_surf, qv_surf, &
654                           h_input, th_input, qv_input, u_input, v_input,n, nl, debug )
656       if(dry) then
657        do k=1,nl
658          qv_input(k) = 0.
659        enddo
660       endif
662       if(debug) write(6,*) ' number of input levels = ',nl
664         nl_in = nl
665         if(nl_in .gt. nl_max ) then
666           write(6,*) ' too many levels for input arrays ',nl_in,nl_max
667           call wrf_error_fatal ( ' too many levels for input arrays ' )
668         end if
670 !  compute diagnostics,
671 !  first, convert qv(g/kg) to qv(g/g)
673       do k=1,nl
674         qv_input(k) = 0.001*qv_input(k)
675       enddo
677       p_surf = 100.*p_surf  ! convert to pascals
678       qvf = 1. + rvovrd*qv_input(1) 
679       rho_surf = 1./((r/p1000mb)*th_surf*qvf*((p_surf/p1000mb)**cvpm))
680       pi_surf = (p_surf/p1000mb)**(r/cp)
682       if(debug) then
683         write(6,*) ' surface density is ',rho_surf
684         write(6,*) ' surface pi is      ',pi_surf
685       end if
688 !  integrate moist sounding hydrostatically, starting from the
689 !  specified surface pressure
690 !  -> first, integrate from surface to lowest level
692           qvf = 1. + rvovrd*qv_input(1) 
693           qvf1 = 1. + qv_input(1)
694           rho_input(1) = rho_surf
695           dz = h_input(1)-zsfc
697 ! error check here
698           if ( dz < 0.0 ) then
699             write(6,*) "Your first input sounding level is below the WRF terrain elevation, aborting"
700             stop "module_initialize_scm_xy:get_sounding"
701           endif
702           do it=1,10
703             pm_input(1) = p_surf &
704                     - 0.5*dz*(rho_surf+rho_input(1))*g*qvf1
705             rho_input(1) = 1./((r/p1000mb)*th_input(1)*qvf*((pm_input(1)/p1000mb)**cvpm))
706           enddo
708 ! integrate up the column
710           do k=2,nl
711             rho_input(k) = rho_input(k-1)
712             dz = h_input(k)-h_input(k-1)
713             qvf1 = 0.5*(2.+(qv_input(k-1)+qv_input(k)))
714             qvf = 1. + rvovrd*qv_input(k)   ! qv is in g/kg here
716             do it=1,10
717               pm_input(k) = pm_input(k-1) &
718                       - 0.5*dz*(rho_input(k)+rho_input(k-1))*g*qvf1
719               rho_input(k) = 1./((r/p1000mb)*th_input(k)*qvf*((pm_input(k)/p1000mb)**cvpm))
720             enddo
721           enddo
722 !  we have the moist sounding
724 !  next, compute the dry sounding using p at the highest level from the
725 !  moist sounding and integrating down.
727         p_input(nl) = pm_input(nl)
729           do k=nl-1,1,-1
730             dz = h_input(k+1)-h_input(k)
731             p_input(k) = p_input(k+1) + 0.5*dz*(rho_input(k)+rho_input(k+1))*g
732           enddo
735         do k=1,nl
736           zk(k) = h_input(k)
737           p(k) = pm_input(k)
738           p_dry(k) = p_input(k)
739           theta(k) = th_input(k)
740           rho(k) = rho_input(k)
741           u(k) = u_input(k)
742           v(k) = v_input(k)
743           qv(k) = qv_input(k) 
745         enddo
747      if(debug) then
748       write(6,*) ' sounding '
749       write(6,*) '  k  height(m)  press (Pa) pd(Pa) theta (K) den(kg/m^3)  u(m/s)     v(m/s)    qv(g/g) '
750       do k=1,nl
751         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)
752       enddo
754      end if
756       end subroutine get_sounding
758 !-------------------------------------------------------
760       subroutine read_sounding( zsfc,ps,ts,qvs,h,th,qv,u,v,n,nl,debug )
761       implicit none
762       integer n,nl
763       real zsfc,ps,ts,qvs,h(n),th(n),qv(n),u(n),v(n)
764       real u10,v10,t2,q2
765       logical end_of_file
766       logical debug
768       integer k
770       open(unit=10,file='input_sounding',form='formatted',status='old')
771       rewind(10)
773       read(10,*) zsfc, u10, v10, t2, q2, ps
774       ps = ps/100.0
775       ts = t2
776       qvs = q2*1000
778       if(debug) then
779         write(6,*) ' input sounding surface parameters '
780         write(6,*) ' surface pressure (mb) ',ps
781         write(6,*) ' surface pot. temp (K) ',ts
782         write(6,*) ' surface mixing ratio (g/kg) ',qvs
783       end if
785       end_of_file = .false.
786       k = 0
788       do while (.not. end_of_file)
790         read(10,*,end=100) h(k+1), u(k+1), v(k+1), th(k+1), qv(k+1)
792         qv(k+1) = qv(k+1)*1000.0
793         k = k+1
794         if(debug) write(6,'(1x,i3,5(1x,e10.3))') k, h(k), th(k), qv(k), u(k), v(k)
795         go to 110
796  100    end_of_file = .true.
797  110    continue
798       enddo
800       nl = k
802       close(unit=10,status = 'keep')
804       end subroutine read_sounding
806 !-------------------------------------------------------
808       subroutine read_soil( n,nl,tmn,tsk,zs,tslb,smois )
809       implicit none
810       integer n,nl
811       real tmn,tsk
812       real zs(n),tslb(n),smois(n)
813       logical end_of_file
814       logical debug 
816       integer k
817     
818       debug = .true.
820       open(unit=11,file='input_soil',form='formatted',status='old')
821       rewind(11)
823       read(11,*) zs(1),tmn,tsk
825       if(debug) then
826         write(6,*) ' input deep soil temperature (K) ',tmn
827         write(6,*) ' input skin temperature (K) ',tsk
828       end if
830       end_of_file = .false.
831       k = 0
833       do while (.not. end_of_file)
835         read(11,*,end=100) zs(k+1), tslb(k+1), smois(k+1)
836         k = k+1
837         if(debug) write(6,'(1x,i3,3(1x,f16.7))') k, zs(k), tslb(k), smois(k)
838         go to 110
839  100    end_of_file = .true.
840  110    continue
841       enddo
843       nl = k
845       close(unit=11,status = 'keep')
847       end subroutine read_soil
849 END MODULE module_initialize_ideal