source: lmdz_wrf/WRFV3/dyn_em/module_initialize_quarter_ss.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: 23.5 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#ifdef DM_PARALLEL
27   USE module_dm
28#endif
29
30
31CONTAINS
32
33
34!-------------------------------------------------------------------
35! this is a wrapper for the solver-specific init_domain routines.
36! Also dereferences the grid variables and passes them down as arguments.
37! This is crucial, since the lower level routines may do message passing
38! and this will get fouled up on machines that insist on passing down
39! copies of assumed-shape arrays (by passing down as arguments, the
40! data are treated as assumed-size -- ie. f77 -- arrays and the copying
41! business is avoided).  Fie on the F90 designers.  Fie and a pox.
42
43   SUBROUTINE init_domain ( grid )
44
45   IMPLICIT NONE
46
47   !  Input data.
48   TYPE (domain), POINTER :: grid
49   !  Local data.
50   INTEGER :: idum1, idum2
51
52   CALL set_scalar_indices_from_config ( head_grid%id , idum1, idum2 )
53
54     CALL init_domain_rk( grid &
55!
56#include <actual_new_args.inc>
57!
58                        )
59
60   END SUBROUTINE init_domain
61
62!-------------------------------------------------------------------
63
64   SUBROUTINE init_domain_rk ( grid &
65!
66# include <dummy_new_args.inc>
67!
68)
69   IMPLICIT NONE
70
71   !  Input data.
72   TYPE (domain), POINTER :: grid
73
74# include <dummy_new_decl.inc>
75
76   TYPE (grid_config_rec_type)              :: config_flags
77
78   !  Local data
79   INTEGER                             ::                       &
80                                  ids, ide, jds, jde, kds, kde, &
81                                  ims, ime, jms, jme, kms, kme, &
82                                  its, ite, jts, jte, kts, kte, &
83                                  i, j, k
84
85   ! Local data
86
87   INTEGER, PARAMETER :: nl_max = 1000
88   REAL, DIMENSION(nl_max) :: zk, p_in, theta, rho, u, v, qv, pd_in
89   INTEGER :: nl_in
90
91
92   INTEGER :: icm,jcm, ii, im1, jj, jm1, loop, error, fid, nxc, nyc
93   REAL    :: u_mean,v_mean, f0, p_surf, p_level, qvf, z_at_v, z_at_u
94   REAL    :: z_scale, xrad, yrad, zrad, rad, delt, cof1, cof2
95!   REAL, EXTERNAL :: interp_0
96   REAL    :: hm
97   REAL    :: pi
98
99!  stuff from original initialization that has been dropped from the Registry
100   REAL    :: vnu, xnu, xnus, dinit0, cbh, p0_temp, t0_temp, zd, zt
101   REAL    :: qvf1, qvf2, pd_surf
102   INTEGER :: it
103   real :: thtmp, ptmp, temp(3)
104
105   LOGICAL :: moisture_init
106   LOGICAL :: stretch_grid, dry_sounding
107
108  INTEGER :: xs , xe , ys , ye
109  REAL :: mtn_ht
110   LOGICAL, EXTERNAL :: wrf_dm_on_monitor
111
112   SELECT CASE ( model_data_order )
113         CASE ( DATA_ORDER_ZXY )
114   kds = grid%sd31 ; kde = grid%ed31 ;
115   ids = grid%sd32 ; ide = grid%ed32 ;
116   jds = grid%sd33 ; jde = grid%ed33 ;
117
118   kms = grid%sm31 ; kme = grid%em31 ;
119   ims = grid%sm32 ; ime = grid%em32 ;
120   jms = grid%sm33 ; jme = grid%em33 ;
121
122   kts = grid%sp31 ; kte = grid%ep31 ;   ! note that tile is entire patch
123   its = grid%sp32 ; ite = grid%ep32 ;   ! note that tile is entire patch
124   jts = grid%sp33 ; jte = grid%ep33 ;   ! note that tile is entire patch
125         CASE ( DATA_ORDER_XYZ )
126   ids = grid%sd31 ; ide = grid%ed31 ;
127   jds = grid%sd32 ; jde = grid%ed32 ;
128   kds = grid%sd33 ; kde = grid%ed33 ;
129
130   ims = grid%sm31 ; ime = grid%em31 ;
131   jms = grid%sm32 ; jme = grid%em32 ;
132   kms = grid%sm33 ; kme = grid%em33 ;
133
134   its = grid%sp31 ; ite = grid%ep31 ;   ! note that tile is entire patch
135   jts = grid%sp32 ; jte = grid%ep32 ;   ! note that tile is entire patch
136   kts = grid%sp33 ; kte = grid%ep33 ;   ! note that tile is entire patch
137         CASE ( DATA_ORDER_XZY )
138   ids = grid%sd31 ; ide = grid%ed31 ;
139   kds = grid%sd32 ; kde = grid%ed32 ;
140   jds = grid%sd33 ; jde = grid%ed33 ;
141
142   ims = grid%sm31 ; ime = grid%em31 ;
143   kms = grid%sm32 ; kme = grid%em32 ;
144   jms = grid%sm33 ; jme = grid%em33 ;
145
146   its = grid%sp31 ; ite = grid%ep31 ;   ! note that tile is entire patch
147   kts = grid%sp32 ; kte = grid%ep32 ;   ! note that tile is entire patch
148   jts = grid%sp33 ; jte = grid%ep33 ;   ! note that tile is entire patch
149
150   END SELECT
151
152
153   stretch_grid = .true.
154   delt = 3.
155!   z_scale = .50
156   z_scale = .40
157   pi = 2.*asin(1.0)
158   write(6,*) ' pi is ',pi
159   nxc = (ide-ids)/2
160   nyc = (jde-jds)/2
161
162   CALL model_to_grid_config_rec ( grid%id , model_config_rec , config_flags )
163
164! here we check to see if the boundary conditions are set properly
165
166   CALL boundary_condition_check( config_flags, bdyzone, error, grid%id )
167
168   moisture_init = .true.
169
170    grid%itimestep=0
171
172#ifdef DM_PARALLEL
173   CALL wrf_dm_bcast_bytes( icm , IWORDSIZE )
174   CALL wrf_dm_bcast_bytes( jcm , IWORDSIZE )
175#endif
176
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)
188
189
190!  here we initialize data we currently is not initialized
191!  in the input data
192
193    DO j = jts, jte
194      DO i = its, ite
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)        = 0.
206
207      END DO
208   END DO
209
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
217
218   grid%step_number = 0
219
220! set up the grid
221
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
232
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
244
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       
250
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
255
256!  get the sounding from the ascii sounding file, first get dry sounding and
257!  calculate base state
258
259  dry_sounding = .true.
260  IF ( wrf_dm_on_monitor() ) THEN
261  write(6,*) ' getting dry sounding for base state '
262
263  CALL get_sounding( zk, p_in, pd_in, theta, rho, u, v, qv, dry_sounding, nl_max, nl_in )
264  ENDIF
265  CALL wrf_dm_bcast_real( zk , nl_max )
266  CALL wrf_dm_bcast_real( p_in , nl_max )
267  CALL wrf_dm_bcast_real( pd_in , nl_max )
268  CALL wrf_dm_bcast_real( theta , nl_max )
269  CALL wrf_dm_bcast_real( rho , nl_max )
270  CALL wrf_dm_bcast_real( u , nl_max )
271  CALL wrf_dm_bcast_real( v , nl_max )
272  CALL wrf_dm_bcast_real( qv , nl_max )
273  CALL wrf_dm_bcast_integer ( nl_in , 1 )
274
275  write(6,*) ' returned from reading sounding, nl_in is ',nl_in
276
277!  find ptop for the desired ztop (ztop is input from the namelist),
278!  and find surface pressure
279
280  grid%p_top = interp_0( p_in, zk, config_flags%ztop, nl_in )
281
282  DO j=jts,jte
283  DO i=its,ite
284    grid%ht(i,j) = 0.
285  ENDDO
286  ENDDO
287
288  xs=ide/2 -3
289  xs=ids   -3
290  xe=xs + 6
291  ys=jde/2 -3
292  ye=ys + 6
293  mtn_ht = 500
294#ifdef MTN
295  DO j=max(ys,jds),min(ye,jde-1)
296  DO i=max(xs,ids),min(xe,ide-1)
297     grid%ht(i,j) = mtn_ht * 0.25 * &
298               ( 1. + COS ( 2*pi/(xe-xs) * ( i-xs ) + pi ) ) * &
299               ( 1. + COS ( 2*pi/(ye-ys) * ( j-ys ) + pi ) )
300  ENDDO
301  ENDDO
302#endif
303#ifdef EW_RIDGE
304  DO j=max(ys,jds),min(ye,jde-1)
305  DO i=ids,ide
306     grid%ht(i,j) = mtn_ht * 0.50 * &
307               ( 1. + COS ( 2*pi/(ye-ys) * ( j-ys ) + pi ) )
308  ENDDO
309  ENDDO
310#endif
311#ifdef NS_RIDGE
312  DO j=jds,jde
313  DO i=max(xs,ids),min(xe,ide-1)
314     grid%ht(i,j) = mtn_ht * 0.50 * &
315               ( 1. + COS ( 2*pi/(xe-xs) * ( i-xs ) + pi ) )
316  ENDDO
317  ENDDO
318#endif
319  DO j=jts,jte
320  DO i=its,ite
321    grid%phb(i,1,j) = g * grid%ht(i,j)
322    grid%ph0(i,1,j) = g * grid%ht(i,j)
323  ENDDO
324  ENDDO
325
326  DO J = jts, jte
327  DO I = its, ite
328
329    p_surf = interp_0( p_in, zk, grid%phb(i,1,j)/g, nl_in )
330    grid%mub(i,j) = p_surf-grid%p_top
331
332!  this is dry hydrostatic sounding (base state), so given grid%p (coordinate),
333!  interp theta (from interp) and compute 1/rho from eqn. of state
334
335    DO K = 1, kte-1
336      p_level = grid%znu(k)*(p_surf - grid%p_top) + grid%p_top
337      grid%pb(i,k,j) = p_level
338      grid%t_init(i,k,j) = interp_0( theta, p_in, p_level, nl_in ) - t0
339      grid%alb(i,k,j) = (r_d/p1000mb)*(grid%t_init(i,k,j)+t0)*(grid%pb(i,k,j)/p1000mb)**cvpm
340    ENDDO
341
342!  calc hydrostatic balance (alternatively we could interp the geopotential from the
343!  sounding, but this assures that the base state is in exact hydrostatic balance with
344!  respect to the model eqns.
345
346    DO k  = 2,kte
347      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)
348    ENDDO
349
350  ENDDO
351  ENDDO
352
353  IF ( wrf_dm_on_monitor() ) THEN
354    write(6,*) ' ptop is ',grid%p_top
355    write(6,*) ' base state grid%mub(1,1), p_surf is ',grid%mub(1,1),grid%mub(1,1)+grid%p_top
356  ENDIF
357
358!  calculate full state for each column - this includes moisture.
359
360  write(6,*) ' getting moist sounding for full state '
361  dry_sounding = .false.
362  CALL get_sounding( zk, p_in, pd_in, theta, rho, u, v, qv, dry_sounding, nl_max, nl_in )
363
364  DO J = jts, min(jde-1,jte)
365  DO I = its, min(ide-1,ite)
366
367!  At this point grid%p_top is already set. find the DRY mass in the column
368!  by interpolating the DRY pressure. 
369
370   pd_surf = interp_0( pd_in, zk, grid%phb(i,1,j)/g, nl_in )
371
372!  compute the perturbation mass and the full mass
373
374    grid%mu_1(i,j) = pd_surf-grid%p_top - grid%mub(i,j)
375    grid%mu_2(i,j) = grid%mu_1(i,j)
376    grid%mu0(i,j) = grid%mu_1(i,j) + grid%mub(i,j)
377
378! given the dry pressure and coordinate system, interp the potential
379! temperature and qv
380
381    do k=1,kde-1
382
383      p_level = grid%znu(k)*(pd_surf - grid%p_top) + grid%p_top
384
385      moist(i,k,j,P_QV) = interp_0( qv, pd_in, p_level, nl_in )
386      grid%t_1(i,k,j)          = interp_0( theta, pd_in, p_level, nl_in ) - t0
387      grid%t_2(i,k,j)          = grid%t_1(i,k,j)
388     
389
390    enddo
391
392!  integrate the hydrostatic equation (from the RHS of the bigstep
393!  vertical momentum equation) down from the top to get grid%p.
394!  first from the top of the model to the top pressure
395
396    k = kte-1  ! top level
397
398    qvf1 = 0.5*(moist(i,k,j,P_QV)+moist(i,k,j,P_QV))
399    qvf2 = 1./(1.+qvf1)
400    qvf1 = qvf1*qvf2
401
402!    grid%p(i,k,j) = - 0.5*grid%mu_1(i,j)/grid%rdnw(k)
403    grid%p(i,k,j) = - 0.5*(grid%mu_1(i,j)+qvf1*grid%mub(i,j))/grid%rdnw(k)/qvf2
404    qvf = 1. + rvovrd*moist(i,k,j,P_QV)
405    grid%alt(i,k,j) = (r_d/p1000mb)*(grid%t_1(i,k,j)+t0)*qvf* &
406                (((grid%p(i,k,j)+grid%pb(i,k,j))/p1000mb)**cvpm)
407    grid%al(i,k,j) = grid%alt(i,k,j) - grid%alb(i,k,j)
408
409!  down the column
410
411    do k=kte-2,1,-1
412      qvf1 = 0.5*(moist(i,k,j,P_QV)+moist(i,k+1,j,P_QV))
413      qvf2 = 1./(1.+qvf1)
414      qvf1 = qvf1*qvf2
415      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)
416      qvf = 1. + rvovrd*moist(i,k,j,P_QV)
417      grid%alt(i,k,j) = (r_d/p1000mb)*(grid%t_1(i,k,j)+t0)*qvf* &
418                  (((grid%p(i,k,j)+grid%pb(i,k,j))/p1000mb)**cvpm)
419      grid%al(i,k,j) = grid%alt(i,k,j) - grid%alb(i,k,j)
420    enddo
421
422!  this is the hydrostatic equation used in the model after the
423!  small timesteps.  In the model, grid%al (inverse density)
424!  is computed from the geopotential.
425
426
427    grid%ph_1(i,1,j) = 0.
428    DO k  = 2,kte
429      grid%ph_1(i,k,j) = grid%ph_1(i,k-1,j) - (1./grid%rdnw(k-1))*(       &
430                   (grid%mub(i,j)+grid%mu_1(i,j))*grid%al(i,k-1,j)+ &
431                    grid%mu_1(i,j)*grid%alb(i,k-1,j)  )
432                                                   
433      grid%ph_2(i,k,j) = grid%ph_1(i,k,j)
434      grid%ph0(i,k,j) = grid%ph_1(i,k,j) + grid%phb(i,k,j)
435    ENDDO
436
437    IF ( wrf_dm_on_monitor() ) THEN
438    if((i==2) .and. (j==2)) then
439     write(6,*) ' grid%ph_1 calc ',grid%ph_1(2,1,2),grid%ph_1(2,2,2),&
440                              grid%mu_1(2,2)+grid%mub(2,2),grid%mu_1(2,2), &
441                              grid%alb(2,1,2),grid%al(1,2,1),grid%rdnw(1)
442    endif
443    ENDIF
444
445  ENDDO
446  ENDDO
447
448!#if 0
449
450!  thermal perturbation to kick off convection
451
452  write(6,*) ' nxc, nyc for perturbation ',nxc,nyc
453  write(6,*) ' delt for perturbation ',delt
454
455  DO J = jts, min(jde-1,jte)
456    yrad = config_flags%dy*float(j-nyc)/10000.
457!   yrad = 0.
458    DO I = its, min(ide-1,ite)
459      xrad = config_flags%dx*float(i-nxc)/10000.
460!     xrad = 0.
461      DO K = 1, kte-1
462
463!  put in preturbation theta (bubble) and recalc density.  note,
464!  the mass in the column is not changing, so when theta changes,
465!  we recompute density and geopotential
466
467        zrad = 0.5*(grid%ph_1(i,k,j)+grid%ph_1(i,k+1,j)  &
468                   +grid%phb(i,k,j)+grid%phb(i,k+1,j))/g
469        zrad = (zrad-1500.)/1500.
470        RAD=SQRT(xrad*xrad+yrad*yrad+zrad*zrad)
471        IF(RAD <= 1.) THEN
472           grid%t_1(i,k,j)=grid%t_1(i,k,j)+delt*COS(.5*PI*RAD)**2
473           grid%t_2(i,k,j)=grid%t_1(i,k,j)
474           qvf = 1. + rvovrd*moist(i,k,j,P_QV)
475           grid%alt(i,k,j) = (r_d/p1000mb)*(grid%t_1(i,k,j)+t0)*qvf* &
476                        (((grid%p(i,k,j)+grid%pb(i,k,j))/p1000mb)**cvpm)
477           grid%al(i,k,j) = grid%alt(i,k,j) - grid%alb(i,k,j)
478        ENDIF
479      ENDDO
480
481!  rebalance hydrostatically
482
483      DO k  = 2,kte
484        grid%ph_1(i,k,j) = grid%ph_1(i,k-1,j) - (1./grid%rdnw(k-1))*(       &
485                     (grid%mub(i,j)+grid%mu_1(i,j))*grid%al(i,k-1,j)+ &
486                      grid%mu_1(i,j)*grid%alb(i,k-1,j)  )
487                                                   
488        grid%ph_2(i,k,j) = grid%ph_1(i,k,j)
489        grid%ph0(i,k,j) = grid%ph_1(i,k,j) + grid%phb(i,k,j)
490      ENDDO
491
492    ENDDO
493  ENDDO
494
495!#endif
496
497   IF ( wrf_dm_on_monitor() ) THEN
498   write(6,*) ' grid%mu_1 from comp ', grid%mu_1(1,1)
499   write(6,*) ' full state sounding from comp, ph, grid%p, grid%al, grid%t_1, qv '
500   do k=1,kde-1
501     write(6,'(i3,1x,5(1x,1pe10.3))') k, grid%ph_1(1,k,1)+grid%phb(1,k,1), &
502                                      grid%p(1,k,1)+grid%pb(1,k,1), grid%alt(1,k,1), &
503                                      grid%t_1(1,k,1)+t0, moist(1,k,1,P_QV)
504   enddo
505
506   write(6,*) ' pert state sounding from comp, grid%ph_1, pp, alp, grid%t_1, qv '
507   do k=1,kde-1
508     write(6,'(i3,1x,5(1x,1pe10.3))') k, grid%ph_1(1,k,1), &
509                                      grid%p(1,k,1), grid%al(1,k,1), &
510                                      grid%t_1(1,k,1), moist(1,k,1,P_QV)
511   enddo
512   ENDIF
513
514! interp v
515
516  DO J = jts, jte
517  DO I = its, min(ide-1,ite)
518
519    IF (j == jds) THEN
520      z_at_v = grid%phb(i,1,j)/g
521    ELSE IF (j == jde) THEN
522      z_at_v = grid%phb(i,1,j-1)/g
523    ELSE
524      z_at_v = 0.5*(grid%phb(i,1,j)+grid%phb(i,1,j-1))/g
525    END IF
526    p_surf = interp_0( p_in, zk, z_at_v, nl_in )
527
528    DO K = 1, kte-1
529      p_level = grid%znu(k)*(p_surf - grid%p_top) + grid%p_top
530      grid%v_1(i,k,j) = interp_0( v, p_in, p_level, nl_in )
531      grid%v_2(i,k,j) = grid%v_1(i,k,j)
532    ENDDO
533
534  ENDDO
535  ENDDO
536
537! interp u
538
539  DO J = jts, min(jde-1,jte)
540  DO I = its, ite
541
542    IF (i == ids) THEN
543      z_at_u = grid%phb(i,1,j)/g
544    ELSE IF (i == ide) THEN
545      z_at_u = grid%phb(i-1,1,j)/g
546    ELSE
547      z_at_u = 0.5*(grid%phb(i,1,j)+grid%phb(i-1,1,j))/g
548    END IF
549
550    p_surf = interp_0( p_in, zk, z_at_u, nl_in )
551
552    DO K = 1, kte-1
553      p_level = grid%znu(k)*(p_surf - grid%p_top) + grid%p_top
554      grid%u_1(i,k,j) = interp_0( u, p_in, p_level, nl_in )
555      grid%u_2(i,k,j) = grid%u_1(i,k,j)
556    ENDDO
557
558  ENDDO
559  ENDDO
560
561!  set w
562
563  DO J = jts, min(jde-1,jte)
564  DO K = kts, kte
565  DO I = its, min(ide-1,ite)
566    grid%w_1(i,k,j) = 0.
567    grid%w_2(i,k,j) = 0.
568  ENDDO
569  ENDDO
570  ENDDO
571
572!  set a few more things
573
574  DO J = jts, min(jde-1,jte)
575  DO K = kts, kte-1
576  DO I = its, min(ide-1,ite)
577    grid%h_diabatic(i,k,j) = 0.
578  ENDDO
579  ENDDO
580  ENDDO
581
582  IF ( wrf_dm_on_monitor() ) THEN
583  DO k=1,kte-1
584    grid%t_base(k) = grid%t_1(1,k,1)
585    grid%qv_base(k) = moist(1,k,1,P_QV)
586    grid%u_base(k) = grid%u_1(1,k,1)
587    grid%v_base(k) = grid%v_1(1,k,1)
588    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
589  ENDDO
590  ENDIF
591  CALL wrf_dm_bcast_real( grid%t_base , kte )
592  CALL wrf_dm_bcast_real( grid%qv_base , kte )
593  CALL wrf_dm_bcast_real( grid%u_base , kte )
594  CALL wrf_dm_bcast_real( grid%v_base , kte )
595  CALL wrf_dm_bcast_real( grid%z_base , kte )
596
597  DO J = jts, min(jde-1,jte)
598  DO I = its, min(ide-1,ite)
599     thtmp   = grid%t_2(i,1,j)+t0
600     ptmp    = grid%p(i,1,j)+grid%pb(i,1,j)
601     temp(1) = thtmp * (ptmp/p1000mb)**rcp
602     thtmp   = grid%t_2(i,2,j)+t0
603     ptmp    = grid%p(i,2,j)+grid%pb(i,2,j)
604     temp(2) = thtmp * (ptmp/p1000mb)**rcp
605     thtmp   = grid%t_2(i,3,j)+t0
606     ptmp    = grid%p(i,3,j)+grid%pb(i,3,j)
607     temp(3) = thtmp * (ptmp/p1000mb)**rcp
608
609     grid%tsk(I,J)=grid%cf1*temp(1)+grid%cf2*temp(2)+grid%cf3*temp(3)
610     grid%tmn(I,J)=grid%tsk(I,J)-0.5
611  ENDDO
612  ENDDO
613
614 END SUBROUTINE init_domain_rk
615
616   SUBROUTINE init_module_initialize
617   END SUBROUTINE init_module_initialize
618
619!---------------------------------------------------------------------
620
621!  test driver for get_sounding
622!
623!      implicit none
624!      integer n
625!      parameter(n = 1000)
626!      real zk(n),p(n),theta(n),rho(n),u(n),v(n),qv(n),pd(n)
627!      logical dry
628!      integer nl,k
629!
630!      dry = .false.
631!      dry = .true.
632!      call get_sounding( zk, p, pd, theta, rho, u, v, qv, dry, n, nl )
633!      write(6,*) ' input levels ',nl
634!      write(6,*) ' sounding '
635!      write(6,*) '  k  height(m)  press (Pa) pd(Pa) theta (K) den(kg/m^3)  u(m/s)     v(m/s)    qv(g/g) '
636!      do k=1,nl
637!        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)
638!      enddo
639!      end
640!
641!---------------------------------------------------------------------------
642
643      subroutine get_sounding( zk, p, p_dry, theta, rho, &
644                               u, v, qv, dry, nl_max, nl_in )
645      implicit none
646
647      integer nl_max, nl_in
648      real zk(nl_max), p(nl_max), theta(nl_max), rho(nl_max), &
649           u(nl_max), v(nl_max), qv(nl_max), p_dry(nl_max)
650      logical dry
651
652      integer n
653      parameter(n=1000)
654      logical debug
655      parameter( debug = .true.)
656
657! input sounding data
658
659      real p_surf, th_surf, qv_surf
660      real pi_surf, pi(n)
661      real h_input(n), th_input(n), qv_input(n), u_input(n), v_input(n)
662
663! diagnostics
664
665      real rho_surf, p_input(n), rho_input(n)
666      real pm_input(n)  !  this are for full moist sounding
667
668! local data
669
670      real r
671      parameter (r = r_d)
672      integer k, it, nl
673      real qvf, qvf1, dz
674
675!  first, read the sounding
676
677      call read_sounding( p_surf, th_surf, qv_surf, &
678                          h_input, th_input, qv_input, u_input, v_input,n, nl, debug )
679
680      if(dry) then
681       do k=1,nl
682         qv_input(k) = 0.
683       enddo
684      endif
685
686      if(debug) write(6,*) ' number of input levels = ',nl
687
688        nl_in = nl
689        if(nl_in .gt. nl_max ) then
690          write(6,*) ' too many levels for input arrays ',nl_in,nl_max
691          call wrf_error_fatal ( ' too many levels for input arrays ' )
692        end if
693
694!  compute diagnostics,
695!  first, convert qv(g/kg) to qv(g/g)
696
697      do k=1,nl
698        qv_input(k) = 0.001*qv_input(k)
699      enddo
700
701      p_surf = 100.*p_surf  ! convert to pascals
702      qvf = 1. + rvovrd*qv_input(1)
703      rho_surf = 1./((r/p1000mb)*th_surf*qvf*((p_surf/p1000mb)**cvpm))
704      pi_surf = (p_surf/p1000mb)**(r/cp)
705
706      if(debug) then
707        write(6,*) ' surface density is ',rho_surf
708        write(6,*) ' surface pi is      ',pi_surf
709      end if
710
711
712!  integrate moist sounding hydrostatically, starting from the
713!  specified surface pressure
714!  -> first, integrate from surface to lowest level
715
716          qvf = 1. + rvovrd*qv_input(1)
717          qvf1 = 1. + qv_input(1)
718          rho_input(1) = rho_surf
719          dz = h_input(1)
720          do it=1,10
721            pm_input(1) = p_surf &
722                    - 0.5*dz*(rho_surf+rho_input(1))*g*qvf1
723            rho_input(1) = 1./((r/p1000mb)*th_input(1)*qvf*((pm_input(1)/p1000mb)**cvpm))
724          enddo
725
726! integrate up the column
727
728          do k=2,nl
729            rho_input(k) = rho_input(k-1)
730            dz = h_input(k)-h_input(k-1)
731            qvf1 = 0.5*(2.+(qv_input(k-1)+qv_input(k)))
732            qvf = 1. + rvovrd*qv_input(k)   ! qv is in g/kg here
733 
734            do it=1,10
735              pm_input(k) = pm_input(k-1) &
736                      - 0.5*dz*(rho_input(k)+rho_input(k-1))*g*qvf1
737              rho_input(k) = 1./((r/p1000mb)*th_input(k)*qvf*((pm_input(k)/p1000mb)**cvpm))
738            enddo
739          enddo
740
741!  we have the moist sounding
742
743!  next, compute the dry sounding using p at the highest level from the
744!  moist sounding and integrating down.
745
746        p_input(nl) = pm_input(nl)
747
748          do k=nl-1,1,-1
749            dz = h_input(k+1)-h_input(k)
750            p_input(k) = p_input(k+1) + 0.5*dz*(rho_input(k)+rho_input(k+1))*g
751          enddo
752
753
754        do k=1,nl
755
756          zk(k) = h_input(k)
757          p(k) = pm_input(k)
758          p_dry(k) = p_input(k)
759          theta(k) = th_input(k)
760          rho(k) = rho_input(k)
761          u(k) = u_input(k)
762          v(k) = v_input(k)
763          qv(k) = qv_input(k)
764
765        enddo
766
767     if(debug) then
768      write(6,*) ' sounding '
769      write(6,*) '  k  height(m)  press (Pa) pd(Pa) theta (K) den(kg/m^3)  u(m/s)     v(m/s)    qv(g/g) '
770      do k=1,nl
771        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)
772      enddo
773
774     end if
775
776      end subroutine get_sounding
777
778!-------------------------------------------------------
779
780      subroutine read_sounding( ps,ts,qvs,h,th,qv,u,v,n,nl,debug )
781      implicit none
782      integer n,nl
783      real ps,ts,qvs,h(n),th(n),qv(n),u(n),v(n)
784      logical end_of_file
785      logical debug
786
787      integer k
788
789      open(unit=10,file='input_sounding',form='formatted',status='old')
790      rewind(10)
791      read(10,*) ps, ts, qvs
792      if(debug) then
793        write(6,*) ' input sounding surface parameters '
794        write(6,*) ' surface pressure (mb) ',ps
795        write(6,*) ' surface pot. temp (K) ',ts
796        write(6,*) ' surface mixing ratio (g/kg) ',qvs
797      end if
798
799      end_of_file = .false.
800      k = 0
801
802      do while (.not. end_of_file)
803
804        read(10,*,end=100) h(k+1), th(k+1), qv(k+1), u(k+1), v(k+1)
805        k = k+1
806        if(debug) write(6,'(1x,i3,5(1x,e10.3))') k, h(k), th(k), qv(k), u(k), v(k)
807        go to 110
808 100    end_of_file = .true.
809 110    continue
810      enddo
811
812      nl = k
813
814      close(unit=10,status = 'keep')
815
816      end subroutine read_sounding
817
818END MODULE module_initialize_ideal
Note: See TracBrowser for help on using the repository browser.