source: lmdz_wrf/WRFV3/dyn_em/module_initialize_scm_xy.F @ 1

Last change on this file since 1 was 1, checked in by lfita, 10 years ago
  • -- --- Opening of the WRF+LMDZ coupling repository --- -- -

WRF: version v3.3
LMDZ: version v1818

More details in:

File size: 25.6 KB
Line 
1!IDEAL:MODEL_LAYER:INITIALIZATION
2!
3
4!  This MODULE holds the routines which are used to perform various initializations
5!  for the individual domains. 
6
7!  This MODULE CONTAINS the following routines:
8
9!  initialize_field_test - 1. Set different fields to different constant
10!                             values.  This is only a test.  If the correct
11!                             domain is not found (based upon the "id")
12!                             then a fatal error is issued.               
13
14!-----------------------------------------------------------------------
15
16MODULE module_initialize_ideal
17
18   USE module_domain
19   USE module_io_domain
20   USE module_state_description
21   USE module_model_constants
22   USE module_bc
23   USE module_timing
24   USE module_configure
25   USE module_init_utilities
26   USE module_soil_pre
27#ifdef DM_PARALLEL
28   USE module_dm
29#endif
30
31
32CONTAINS
33
34
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. 
46
47   SUBROUTINE init_domain ( grid )
48
49   IMPLICIT NONE
50
51   !  Input data.
52   TYPE (domain), POINTER :: grid
53   !  Local data.
54   INTEGER :: idum1, idum2
55
56   CALL set_scalar_indices_from_config ( head_grid%id , idum1, idum2 )
57
58     CALL init_domain_rk( grid &
59!
60#include <actual_new_args.inc>
61!
62                        )
63   END SUBROUTINE init_domain
64
65!-------------------------------------------------------------------
66
67   SUBROUTINE init_domain_rk ( grid &
68!
69# include <dummy_new_args.inc>
70!
71)
72
73   USE module_optional_input
74   IMPLICIT NONE
75
76   !  Input data.
77   TYPE (domain), POINTER :: grid
78
79# include <dummy_new_decl.inc>
80
81   TYPE (grid_config_rec_type)              :: config_flags
82
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
89
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
96
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
101
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
106
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
113
114   LOGICAL :: moisture_init
115   LOGICAL :: dry_sounding
116   character (len=256) :: mminlu2
117
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.
123
124   REAL    :: zrwa(200), zwa(200)
125
126   
127#ifdef DM_PARALLEL
128#    include <data_calls.inc>
129#endif
130
131
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 ;
137
138   kms = grid%sm31 ; kme = grid%em31 ;
139   ims = grid%sm32 ; ime = grid%em32 ;
140   jms = grid%sm33 ; jme = grid%em33 ;
141
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 ;
149
150   ims = grid%sm31 ; ime = grid%em31 ;
151   jms = grid%sm32 ; jme = grid%em32 ;
152   kms = grid%sm33 ; kme = grid%em33 ;
153
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 ;
161
162   ims = grid%sm31 ; ime = grid%em31 ;
163   kms = grid%sm32 ; kme = grid%em32 ;
164   jms = grid%sm33 ; jme = grid%em33 ;
165
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
169
170   END SELECT
171
172
173   pi = 2.*asin(1.0)
174   write(6,*) ' pi is ',pi
175
176   CALL model_to_grid_config_rec ( grid%id , model_config_rec , config_flags )
177
178! here we check to see if the boundary conditions are set properly
179
180   CALL boundary_condition_check( config_flags, bdyzone, error, grid%id )
181
182   moisture_init = .true.
183
184    grid%itimestep=0
185
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)
201
202!  here we initialize data that currently is not initialized
203!  in the input data
204
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
228
229! for LSM, additional variables need to be initialized
230
231!  other_masked_fields : SELECT CASE ( model_config_rec%sf_surface_physics(grid%id) )
232
233!     CASE (SLABSCHEME)
234
235!     CASE (LSMSCHEME)
236
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
248
249!     CASE (RUCLSMSCHEME)
250
251!  END SELECT other_masked_fields
252
253   grid%step_number = 0
254
255   IF ( real_soil ) THEN ! from input file
256
257      IF (config_flags%sf_surface_physics .NE. 2) WRITE (6, *)   &
258         'If using LSM option other than Noah, must edit input_soil file in test/em_scm_xy/ directory'
259   
260      CALL read_soil(100,ns_input,tmn_input,tsk_input,zs_input,tslb_input,smois_input)
261
262      CALL init_module_optional_input(grid,config_flags)
263      num_st_levels_input = ns_input
264      num_sm_levels_input = ns_input
265      num_sw_levels_input = ns_input
266      DO k = 1,ns_input
267         st_levels_input(k) = zs_input(k)*100.0 ! to cm
268         sm_levels_input(k) = zs_input(k)*100.0 ! to cm
269         sw_levels_input(k) = zs_input(k)*100.0 ! to cm
270         st_input(:,k+1,:) = tslb_input(k)
271         sm_input(:,k+1,:) = smois_input(k)
272         sw_input(:,k+1,:) = smois_input(k)
273      ENDDO
274 
275      grid%tsk = tsk_input
276      grid%sst = tsk_input
277      grid%tmn = tmn_input
278
279      flag_soil_layers  = 0 ! go ahead and put skin temp in
280      flag_soil_levels  = 0 ! go ahead and put skin moisture in
281      flag_sst          = 0 ! don't modify for ocean
282      flag_tavgsfc      = 0
283      flag_soilhgt      = 0
284
285      CALL process_soil_real ( grid%tsk , grid%tmn , grid%tavgsfc, &
286                   grid%landmask , grid%sst , grid%ht, grid%toposoil, &
287                   st_input , sm_input , sw_input , &
288                   st_levels_input , sm_levels_input , sw_levels_input , &
289                   grid%zs , grid%dzs , grid%tslb , grid%smois , grid%sh2o , &
290                   flag_sst , flag_tavgsfc, flag_soilhgt, flag_soil_layers, flag_soil_levels,  &
291                   ids , ide , jds , jde , kds , kde , &
292                   ims , ime , jms , jme , kms , kme , &
293                   its , ite , jts , jte , kts , kte , &
294                   model_config_rec%sf_surface_physics(grid%id) , &
295                   model_config_rec%num_soil_layers , &
296                   model_config_rec%real_data_init_type , &
297                   num_st_levels_input , num_sm_levels_input , num_sw_levels_input , &
298                   num_st_levels_alloc , num_sm_levels_alloc , num_sw_levels_alloc )
299
300   ELSE ! ideal soil
301! Process the soil; note that there are some things hard-wired into share/module_soil_pre.F
302      CALL process_soil_ideal(grid%xland,grid%xice,grid%vegfra,grid%snow,grid%canwat, &
303                     grid%ivgtyp,grid%isltyp,grid%tslb,grid%smois, &
304                     grid%tsk,grid%tmn,grid%zs,grid%dzs,model_config_rec%num_soil_layers, &
305                     model_config_rec%sf_surface_physics(grid%id), &
306                                   ids,ide, jds,jde, kds,kde,&
307                                   ims,ime, jms,jme, kms,kme,&
308                                   its,ite, jts,jte, kts,kte )
309
310    ENDIF
311
312    DO j = jts, jte
313     DO k = kts, kte
314       DO i = its, ite
315          grid%ww(i,k,j)     = 0.
316       END DO
317     END DO
318    END DO
319
320! this is adopted from Wayne Angevine's GABLS case
321    grid%znw(1) = 1.0
322    zrwa(kde) = exp((kde-1)/40.)
323    zwa(kde) = grid%ztop
324    DO k=2, kde-1
325       zrwa(k) = exp((k-1)/40.)
326       zwa(k) = (zrwa(k)-1.) * grid%ztop/(zrwa(kde)-1.)
327       grid%znw(k) = 1. - (zwa(k) / zwa(kde))
328    ENDDO
329    grid%znw(kde) = 0.
330
331   DO k=1, kde-1
332    grid%dnw(k) = grid%znw(k+1) - grid%znw(k)
333    grid%rdnw(k) = 1./grid%dnw(k)
334    grid%znu(k) = 0.5*(grid%znw(k+1)+grid%znw(k))
335   ENDDO
336   DO k=2, kde-1
337    grid%dn(k) = 0.5*(grid%dnw(k)+grid%dnw(k-1))
338    grid%rdn(k) = 1./grid%dn(k)
339    grid%fnp(k) = .5* grid%dnw(k  )/grid%dn(k)
340    grid%fnm(k) = .5* grid%dnw(k-1)/grid%dn(k)
341   ENDDO
342
343   cof1 = (2.*grid%dn(2)+grid%dn(3))/(grid%dn(2)+grid%dn(3))*grid%dnw(1)/grid%dn(2)
344   cof2 =     grid%dn(2)        /(grid%dn(2)+grid%dn(3))*grid%dnw(1)/grid%dn(3)
345   grid%cf1  = grid%fnp(2) + cof1
346   grid%cf2  = grid%fnm(2) - cof1 - cof2
347   grid%cf3  = cof2       
348
349   grid%cfn  = (.5*grid%dnw(kde-1)+grid%dn(kde-1))/grid%dn(kde-1)
350   grid%cfn1 = -.5*grid%dnw(kde-1)/grid%dn(kde-1)
351   grid%rdx = 1./config_flags%dx
352   grid%rdy = 1./config_flags%dy
353
354!  get the sounding from the ascii sounding file, first get dry sounding and
355!  calculate base state
356
357  write(6,*) ' getting dry sounding for base state '
358  dry_sounding = .true.
359  CALL get_sounding( zsfc, zk, p_in, pd_in, theta, rho, u, v, qv, dry_sounding, nl_max, nl_in )
360
361  write(6,*) ' returned from reading sounding, nl_in is ',nl_in
362
363!  find ptop for the desired ztop (ztop is input from the namelist),
364!  and find surface pressure
365
366  grid%p_top = interp_0( p_in, zk, config_flags%ztop, nl_in )
367
368  DO j=jts,jte
369  DO i=its,ite  ! flat surface
370    grid%ht(i,j) = zsfc
371    grid%phb(i,1,j) = grid%ht(i,j) * g
372    grid%ph0(i,1,j) = grid%ht(i,j) * g
373    grid%php(i,1,j) = 0.
374  ENDDO
375  ENDDO
376
377  DO J = jts, jte
378  DO I = its, ite
379
380    p_surf = interp_0( p_in, zk, grid%phb(i,1,j)/g, nl_in )
381    grid%mub(i,j) = p_surf-grid%p_top
382
383!  this is dry hydrostatic sounding (base state), so given grid%p (coordinate),
384!  interp theta (from interp) and compute 1/rho from eqn. of state
385
386    DO K = 1, kte-1
387      p_level = grid%znu(k)*(p_surf - grid%p_top) + grid%p_top
388      grid%pb(i,k,j) = p_level
389      grid%t_init(i,k,j) = interp_0( theta, p_in, p_level, nl_in ) - t0
390      grid%alb(i,k,j) = (r_d/p1000mb)*(grid%t_init(i,k,j)+t0)*(grid%pb(i,k,j)/p1000mb)**cvpm
391    ENDDO
392
393!  calc hydrostatic balance (alternatively we could interp the geopotential from the
394!  sounding, but this assures that the base state is in exact hydrostatic balance with
395!  respect to the model eqns.
396
397    DO k  = 2,kte
398      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)
399    ENDDO
400
401  ENDDO
402  ENDDO
403
404  write(6,*) ' ptop is ',grid%p_top
405  write(6,*) ' base state grid%mub(1,1), p_surf is ',grid%mub(1,1),grid%mub(1,1)+grid%p_top
406
407!  calculate full state for each column - this includes moisture.
408
409  write(6,*) ' getting moist sounding for full state '
410  dry_sounding = .false.
411  CALL get_sounding( zsfc, zk, p_in, pd_in, theta, rho, u, v, qv, dry_sounding, nl_max, nl_in )
412
413  DO J = jts, min(jde-1,jte)
414  DO I = its, min(ide-1,ite)
415
416!  At this point grid%p_top is already set. find the DRY mass in the column
417!  by interpolating the DRY pressure. 
418
419   pd_surf = interp_0( pd_in, zk, grid%phb(i,1,j)/g, nl_in )
420
421!  compute the perturbation mass and the full mass
422
423    grid%mu_1(i,j) = pd_surf-grid%p_top - grid%mub(i,j)
424    grid%mu_2(i,j) = grid%mu_1(i,j)
425    grid%mu0(i,j) = grid%mu_1(i,j) + grid%mub(i,j)
426
427! given the dry pressure and coordinate system, interp the potential
428! temperature and qv
429
430    do k=1,kde-1
431
432      p_level = grid%znu(k)*(pd_surf - grid%p_top) + grid%p_top
433
434      moist(i,k,j,P_QV) = interp_0( qv, pd_in, p_level, nl_in )
435      grid%t_1(i,k,j)          = interp_0( theta, pd_in, p_level, nl_in ) - t0
436      grid%t_2(i,k,j)          = grid%t_1(i,k,j)
437     
438
439    enddo
440
441!  integrate the hydrostatic equation (from the RHS of the bigstep
442!  vertical momentum equation) down from the top to get grid%p.
443!  first from the top of the model to the top pressure
444
445    k = kte-1  ! top level
446
447    qvf1 = 0.5*(moist(i,k,j,P_QV)+moist(i,k,j,P_QV))
448    qvf2 = 1./(1.+qvf1)
449    qvf1 = qvf1*qvf2
450
451!    grid%p(i,k,j) = - 0.5*grid%mu_1(i,j)/grid%rdnw(k)
452    grid%p(i,k,j) = - 0.5*(grid%mu_1(i,j)+qvf1*grid%mub(i,j))/grid%rdnw(k)/qvf2
453    qvf = 1. + rvovrd*moist(i,k,j,P_QV)
454    grid%alt(i,k,j) = (r_d/p1000mb)*(grid%t_1(i,k,j)+t0)*qvf* &
455                (((grid%p(i,k,j)+grid%pb(i,k,j))/p1000mb)**cvpm)
456    grid%al(i,k,j) = grid%alt(i,k,j) - grid%alb(i,k,j)
457
458!  down the column
459
460    do k=kte-2,1,-1
461      qvf1 = 0.5*(moist(i,k,j,P_QV)+moist(i,k+1,j,P_QV))
462      qvf2 = 1./(1.+qvf1)
463      qvf1 = qvf1*qvf2
464      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)
465      qvf = 1. + rvovrd*moist(i,k,j,P_QV)
466      grid%alt(i,k,j) = (r_d/p1000mb)*(grid%t_1(i,k,j)+t0)*qvf* &
467                  (((grid%p(i,k,j)+grid%pb(i,k,j))/p1000mb)**cvpm)
468      grid%al(i,k,j) = grid%alt(i,k,j) - grid%alb(i,k,j)
469    enddo
470
471!  this is the hydrostatic equation used in the model after the
472!  small timesteps.  In the model, grid%al (inverse density)
473!  is computed from the geopotential.
474
475
476    grid%ph_1(i,1,j) = 0.
477    DO k  = 2,kte
478      grid%ph_1(i,k,j) = grid%ph_1(i,k-1,j) - (1./grid%rdnw(k-1))*(       &
479                   (grid%mub(i,j)+grid%mu_1(i,j))*grid%al(i,k-1,j)+ &
480                    grid%mu_1(i,j)*grid%alb(i,k-1,j)  )
481                                                   
482      grid%ph_2(i,k,j) = grid%ph_1(i,k,j)
483      grid%ph0(i,k,j) = grid%ph_1(i,k,j) + grid%phb(i,k,j)
484    ENDDO
485
486    if((i==2) .and. (j==2)) then
487     write(6,*) ' grid%ph_1 calc ',grid%ph_1(2,1,2),grid%ph_1(2,2,2),&
488                              grid%mu_1(2,2)+grid%mub(2,2),grid%mu_1(2,2), &
489                              grid%alb(2,1,2),grid%al(1,2,1),grid%rdnw(1)
490    endif
491
492  ENDDO
493  ENDDO
494
495   write(6,*) ' grid%mu_1 from comp ', grid%mu_1(1,1)
496   write(6,*) ' full state sounding from comp, ph, grid%p, grid%al, grid%t_1, qv '
497
498   do k=1,kde-1
499     write(6,'(i3,1x,5(1x,1pe10.3))') k, grid%ph_1(1,k,1)+grid%phb(1,k,1), &
500                                      grid%p(1,k,1)+grid%pb(1,k,1), grid%alt(1,k,1), &
501                                      grid%t_1(1,k,1)+t0, moist(1,k,1,P_QV)
502   enddo
503
504   write(6,*) ' pert state sounding from comp, grid%ph_1, pp, alp, grid%t_1, qv '
505   do k=1,kde-1
506     write(6,'(i3,1x,5(1x,1pe10.3))') k, grid%ph_1(1,k,1), &
507                                      grid%p(1,k,1), grid%al(1,k,1), &
508                                      grid%t_1(1,k,1), moist(1,k,1,P_QV)
509   enddo
510
511! interp v
512
513  DO J = jts, jte
514  DO I = its, min(ide-1,ite)
515
516    IF (j == jds) THEN
517      z_at_v = grid%phb(i,1,j)/g
518    ELSE IF (j == jde) THEN
519      z_at_v = grid%phb(i,1,j-1)/g
520    ELSE
521      z_at_v = 0.5*(grid%phb(i,1,j)+grid%phb(i,1,j-1))/g
522    END IF
523
524    p_surf = interp_0( p_in, zk, z_at_v, nl_in )
525
526    DO K = 1, kte
527      p_level = grid%znu(k)*(p_surf - grid%p_top) + grid%p_top
528      grid%v_1(i,k,j) = interp_0( v, p_in, p_level, nl_in )
529      grid%v_2(i,k,j) = grid%v_1(i,k,j)
530    ENDDO
531
532  ENDDO
533  ENDDO
534
535! interp u
536
537  DO J = jts, min(jde-1,jte)
538  DO I = its, ite
539
540    IF (i == ids) THEN
541      z_at_u = grid%phb(i,1,j)/g
542    ELSE IF (i == ide) THEN
543      z_at_u = grid%phb(i-1,1,j)/g
544    ELSE
545      z_at_u = 0.5*(grid%phb(i,1,j)+grid%phb(i-1,1,j))/g
546    END IF
547
548    p_surf = interp_0( p_in, zk, z_at_u, nl_in )
549
550    DO K = 1, kte
551      p_level = grid%znu(k)*(p_surf - grid%p_top) + grid%p_top
552      grid%u_1(i,k,j) = interp_0( u, p_in, p_level, nl_in )
553      grid%u_2(i,k,j) = grid%u_1(i,k,j)
554    ENDDO
555
556  ENDDO
557  ENDDO
558
559!  set w
560
561  DO J = jts, min(jde-1,jte)
562  DO K = kts, kte
563  DO I = its, min(ide-1,ite)
564    grid%w_1(i,k,j) = 0.
565    grid%w_2(i,k,j) = 0.
566  ENDDO
567  ENDDO
568  ENDDO
569
570!  set a few more things
571
572  DO J = jts, min(jde-1,jte)
573  DO K = kts, kte-1
574  DO I = its, min(ide-1,ite)
575    grid%h_diabatic(i,k,j) = 0.
576  ENDDO
577  ENDDO
578  ENDDO
579
580! Go ahead and initialize these from the sounding.  This will allow a run
581! to actually succeed even if scm_force = 0
582  DO k=1,kte-1
583    grid%t_base(k) = grid%t_1(1,k,1)
584    grid%qv_base(k) = moist(1,k,1,P_QV)
585    grid%u_base(k) = grid%u_1(1,k,1)
586    grid%v_base(k) = grid%v_1(1,k,1)
587    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
588  ENDDO
589
590  RETURN
591
592 END SUBROUTINE init_domain_rk
593
594   SUBROUTINE init_module_initialize
595   END SUBROUTINE init_module_initialize
596
597!---------------------------------------------------------------------
598
599!  test driver for get_sounding
600!
601!      implicit none
602!      integer n
603!      parameter(n = 1000)
604!      real zk(n),p(n),theta(n),rho(n),u(n),v(n),qv(n),pd(n)
605!      logical dry
606!      integer nl,k
607!
608!      dry = .false.
609!      dry = .true.
610!      call get_sounding( zk, p, pd, theta, rho, u, v, qv, dry, n, nl )
611!      write(6,*) ' input levels ',nl
612!      write(6,*) ' sounding '
613!      write(6,*) '  k  height(m)  press (Pa) pd(Pa) theta (K) den(kg/m^3)  u(m/s)     v(m/s)    qv(g/g) '
614!      do k=1,nl
615!        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)
616!      enddo
617!      end
618!
619!---------------------------------------------------------------------------
620
621      subroutine get_sounding( zsfc, zk, p, p_dry, theta, rho, &
622                               u, v, qv, dry, nl_max, nl_in )
623      implicit none
624
625      integer nl_max, nl_in
626      real zsfc
627      real zk(nl_max), p(nl_max), theta(nl_max), rho(nl_max), &
628           u(nl_max), v(nl_max), qv(nl_max), p_dry(nl_max)
629      logical dry
630
631      integer n
632      parameter(n=3000)
633      logical debug
634      parameter( debug = .true.)
635
636! input sounding data
637
638      real p_surf, th_surf, qv_surf
639      real pi_surf, pi(n)
640      real h_input(n), th_input(n), qv_input(n), u_input(n), v_input(n)
641
642! diagnostics
643
644      real rho_surf, p_input(n), rho_input(n)
645      real pm_input(n)  !  this are for full moist sounding
646
647! local data
648
649      real r
650      parameter (r = r_d)
651      integer k, it, nl
652      real qvf, qvf1, dz
653
654!  first, read the sounding
655
656      call read_sounding( zsfc, p_surf, th_surf, qv_surf, &
657                          h_input, th_input, qv_input, u_input, v_input,n, nl, debug )
658
659      if(dry) then
660       do k=1,nl
661         qv_input(k) = 0.
662       enddo
663      endif
664
665      if(debug) write(6,*) ' number of input levels = ',nl
666
667        nl_in = nl
668        if(nl_in .gt. nl_max ) then
669          write(6,*) ' too many levels for input arrays ',nl_in,nl_max
670          call wrf_error_fatal ( ' too many levels for input arrays ' )
671        end if
672
673!  compute diagnostics,
674!  first, convert qv(g/kg) to qv(g/g)
675
676      do k=1,nl
677        qv_input(k) = 0.001*qv_input(k)
678      enddo
679
680      p_surf = 100.*p_surf  ! convert to pascals
681      qvf = 1. + rvovrd*qv_input(1)
682      rho_surf = 1./((r/p1000mb)*th_surf*qvf*((p_surf/p1000mb)**cvpm))
683      pi_surf = (p_surf/p1000mb)**(r/cp)
684
685      if(debug) then
686        write(6,*) ' surface density is ',rho_surf
687        write(6,*) ' surface pi is      ',pi_surf
688      end if
689
690
691!  integrate moist sounding hydrostatically, starting from the
692!  specified surface pressure
693!  -> first, integrate from surface to lowest level
694
695          qvf = 1. + rvovrd*qv_input(1)
696          qvf1 = 1. + qv_input(1)
697          rho_input(1) = rho_surf
698          dz = h_input(1)-zsfc
699
700! error check here
701          if ( dz < 0.0 ) then
702            write(6,*) "Your first input sounding level is below the WRF terrain elevation, aborting"
703            stop "module_initialize_scm_xy:get_sounding"
704          endif
705          do it=1,10
706            pm_input(1) = p_surf &
707                    - 0.5*dz*(rho_surf+rho_input(1))*g*qvf1
708            rho_input(1) = 1./((r/p1000mb)*th_input(1)*qvf*((pm_input(1)/p1000mb)**cvpm))
709          enddo
710
711! integrate up the column
712
713          do k=2,nl
714            rho_input(k) = rho_input(k-1)
715            dz = h_input(k)-h_input(k-1)
716            qvf1 = 0.5*(2.+(qv_input(k-1)+qv_input(k)))
717            qvf = 1. + rvovrd*qv_input(k)   ! qv is in g/kg here
718 
719            do it=1,10
720              pm_input(k) = pm_input(k-1) &
721                      - 0.5*dz*(rho_input(k)+rho_input(k-1))*g*qvf1
722              rho_input(k) = 1./((r/p1000mb)*th_input(k)*qvf*((pm_input(k)/p1000mb)**cvpm))
723            enddo
724          enddo
725!  we have the moist sounding
726
727!  next, compute the dry sounding using p at the highest level from the
728!  moist sounding and integrating down.
729
730        p_input(nl) = pm_input(nl)
731
732          do k=nl-1,1,-1
733            dz = h_input(k+1)-h_input(k)
734            p_input(k) = p_input(k+1) + 0.5*dz*(rho_input(k)+rho_input(k+1))*g
735          enddo
736
737
738        do k=1,nl
739          zk(k) = h_input(k)
740          p(k) = pm_input(k)
741          p_dry(k) = p_input(k)
742          theta(k) = th_input(k)
743          rho(k) = rho_input(k)
744          u(k) = u_input(k)
745          v(k) = v_input(k)
746          qv(k) = qv_input(k)
747
748        enddo
749
750     if(debug) then
751      write(6,*) ' sounding '
752      write(6,*) '  k  height(m)  press (Pa) pd(Pa) theta (K) den(kg/m^3)  u(m/s)     v(m/s)    qv(g/g) '
753      do k=1,nl
754        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)
755      enddo
756
757     end if
758
759      end subroutine get_sounding
760
761!-------------------------------------------------------
762
763      subroutine read_sounding( zsfc,ps,ts,qvs,h,th,qv,u,v,n,nl,debug )
764      implicit none
765      integer n,nl
766      real zsfc,ps,ts,qvs,h(n),th(n),qv(n),u(n),v(n)
767      real u10,v10,t2,q2
768      logical end_of_file
769      logical debug
770
771      integer k
772
773      open(unit=10,file='input_sounding',form='formatted',status='old')
774      rewind(10)
775
776      read(10,*) zsfc, u10, v10, t2, q2, ps
777      ps = ps/100.0
778      ts = t2
779      qvs = q2*1000
780
781      if(debug) then
782        write(6,*) ' input sounding surface parameters '
783        write(6,*) ' surface pressure (mb) ',ps
784        write(6,*) ' surface pot. temp (K) ',ts
785        write(6,*) ' surface mixing ratio (g/kg) ',qvs
786      end if
787
788      end_of_file = .false.
789      k = 0
790
791      do while (.not. end_of_file)
792
793        read(10,*,end=100) h(k+1), u(k+1), v(k+1), th(k+1), qv(k+1)
794
795        qv(k+1) = qv(k+1)*1000.0
796        k = k+1
797        if(debug) write(6,'(1x,i3,5(1x,e10.3))') k, h(k), th(k), qv(k), u(k), v(k)
798        go to 110
799 100    end_of_file = .true.
800 110    continue
801      enddo
802
803      nl = k
804
805      close(unit=10,status = 'keep')
806
807      end subroutine read_sounding
808
809!-------------------------------------------------------
810
811      subroutine read_soil( n,nl,tmn,tsk,zs,tslb,smois )
812      implicit none
813      integer n,nl
814      real tmn,tsk
815      real zs(n),tslb(n),smois(n)
816      logical end_of_file
817      logical debug
818
819      integer k
820   
821      debug = .true.
822
823      open(unit=11,file='input_soil',form='formatted',status='old')
824      rewind(11)
825
826      read(11,*) zs(1),tmn,tsk
827
828      if(debug) then
829        write(6,*) ' input deep soil temperature (K) ',tmn
830        write(6,*) ' input skin temperature (K) ',tsk
831      end if
832
833      end_of_file = .false.
834      k = 0
835
836      do while (.not. end_of_file)
837
838        read(11,*,end=100) zs(k+1), tslb(k+1), smois(k+1)
839        k = k+1
840        if(debug) write(6,'(1x,i3,3(1x,f16.7))') k, zs(k), tslb(k), smois(k)
841        go to 110
842 100    end_of_file = .true.
843 110    continue
844      enddo
845
846      nl = k
847
848      close(unit=11,status = 'keep')
849
850      end subroutine read_soil
851
852END MODULE module_initialize_ideal
Note: See TracBrowser for help on using the repository browser.