source: trunk/WRF.COMMON/WRFV2/dyn_em/module_initialize_grav2d_x.F @ 3555

Last change on this file since 3555 was 11, checked in by aslmd, 14 years ago

spiga@svn-planeto:ajoute le modele meso-echelle martien

File size: 24.9 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
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                :: dyn_opt
51   INTEGER :: idum1, idum2
52
53   CALL nl_get_dyn_opt( 1, dyn_opt )
54   
55   CALL set_scalar_indices_from_config ( head_grid%id , idum1, idum2 )
56
57   IF (      dyn_opt .eq. 1 &
58        .or. dyn_opt .eq. 2 &
59        .or. dyn_opt .eq. 3 &
60                                       ) THEN
61     CALL init_domain_rk( grid &
62!
63#include <em_actual_new_args.inc>
64!
65                        )
66
67   ELSE
68     WRITE(0,*)' init_domain: unknown or unimplemented dyn_opt = ',dyn_opt
69     CALL wrf_error_fatal ( ' init_domain: unknown or unimplemented dyn_opt ' )
70   ENDIF
71
72   END SUBROUTINE init_domain
73
74!-------------------------------------------------------------------
75
76   SUBROUTINE init_domain_rk ( grid &
77!
78# include <em_dummy_new_args.inc>
79!
80)
81   IMPLICIT NONE
82
83   !  Input data.
84   TYPE (domain), POINTER :: grid
85
86# include <em_dummy_new_decl.inc>
87
88   TYPE (grid_config_rec_type)              :: config_flags
89
90   !  Local data
91   INTEGER                             ::                       &
92                                  ids, ide, jds, jde, kds, kde, &
93                                  ims, ime, jms, jme, kms, kme, &
94                                  its, ite, jts, jte, kts, kte, &
95                                  i, j, k
96
97   ! Local data
98
99   INTEGER, PARAMETER :: nl_max = 1000
100   REAL, DIMENSION(nl_max) :: zk, p_in, theta, rho, u, v, qv, pd_in
101   INTEGER :: nl_in
102
103
104   INTEGER :: icm,jcm, ii, im1, jj, jm1, loop, error, fid, nxc, nyc
105   REAL    :: u_mean,v_mean, f0, p_surf, p_level, qvf, z_at_v, z_at_u
106   REAL    :: z_scale, xrad, yrad, zrad, rad, delt, cof1, cof2, t_min, t_max
107!   REAL, EXTERNAL :: interp_0
108   REAL    :: hm, xa, xpos, xposml, xpospl
109   REAL    :: pi
110
111!  stuff from original initialization that has been dropped from the Registry
112   REAL    :: vnu, xnu, xnus, dinit0, cbh, p0_temp, t0_temp, zd, zt
113   REAL    :: qvf1, qvf2, pd_surf
114   INTEGER :: it
115   real :: thtmp, ptmp, temp(3)
116
117   LOGICAL :: moisture_init
118   LOGICAL :: stretch_grid, dry_sounding
119
120   REAL    :: xa1, xal1,pii,hm1  !  data for intercomparison setup from dale
121
122#ifdef DM_PARALLEL
123#    include <em_data_calls.inc>
124#endif
125
126
127   SELECT CASE ( model_data_order )
128         CASE ( DATA_ORDER_ZXY )
129   kds = grid%sd31 ; kde = grid%ed31 ;
130   ids = grid%sd32 ; ide = grid%ed32 ;
131   jds = grid%sd33 ; jde = grid%ed33 ;
132
133   kms = grid%sm31 ; kme = grid%em31 ;
134   ims = grid%sm32 ; ime = grid%em32 ;
135   jms = grid%sm33 ; jme = grid%em33 ;
136
137   kts = grid%sp31 ; kte = grid%ep31 ;   ! note that tile is entire patch
138   its = grid%sp32 ; ite = grid%ep32 ;   ! note that tile is entire patch
139   jts = grid%sp33 ; jte = grid%ep33 ;   ! note that tile is entire patch
140         CASE ( DATA_ORDER_XYZ )
141   ids = grid%sd31 ; ide = grid%ed31 ;
142   jds = grid%sd32 ; jde = grid%ed32 ;
143   kds = grid%sd33 ; kde = grid%ed33 ;
144
145   ims = grid%sm31 ; ime = grid%em31 ;
146   jms = grid%sm32 ; jme = grid%em32 ;
147   kms = grid%sm33 ; kme = grid%em33 ;
148
149   its = grid%sp31 ; ite = grid%ep31 ;   ! note that tile is entire patch
150   jts = grid%sp32 ; jte = grid%ep32 ;   ! note that tile is entire patch
151   kts = grid%sp33 ; kte = grid%ep33 ;   ! note that tile is entire patch
152         CASE ( DATA_ORDER_XZY )
153   ids = grid%sd31 ; ide = grid%ed31 ;
154   kds = grid%sd32 ; kde = grid%ed32 ;
155   jds = grid%sd33 ; jde = grid%ed33 ;
156
157   ims = grid%sm31 ; ime = grid%em31 ;
158   kms = grid%sm32 ; kme = grid%em32 ;
159   jms = grid%sm33 ; jme = grid%em33 ;
160
161   its = grid%sp31 ; ite = grid%ep31 ;   ! note that tile is entire patch
162   kts = grid%sp32 ; kte = grid%ep32 ;   ! note that tile is entire patch
163   jts = grid%sp33 ; jte = grid%ep33 ;   ! note that tile is entire patch
164
165   END SELECT
166
167
168   hm = 000.
169   xa = 5.0
170
171   icm = ide/2
172
173
174   xa1  = 5000./500.
175   xal1 = 4000./500.
176   pii  = 2.*asin(1.0)
177   hm1  = 250.
178!   hm1  = 1000.
179
180
181   stretch_grid = .true.
182!   z_scale = .50
183   z_scale = 1.675
184   pi = 2.*asin(1.0)
185   write(6,*) ' pi is ',pi
186   nxc = (ide-ids)/2
187   nyc = (jde-jds)/2
188
189   CALL model_to_grid_config_rec ( grid%id , model_config_rec , config_flags )
190
191! here we check to see if the boundary conditions are set properly
192
193   CALL boundary_condition_check( config_flags, bdyzone, error, grid%id )
194
195   moisture_init = .true.
196
197    grid%itimestep=0
198
199#ifdef DM_PARALLEL
200   CALL wrf_dm_bcast_bytes( icm , IWORDSIZE )
201   CALL wrf_dm_bcast_bytes( jcm , IWORDSIZE )
202#endif
203
204    CALL nl_set_mminlu(1, '    ')
205    CALL nl_set_iswater(1,0)
206    CALL nl_set_cen_lat(1,40.)
207    CALL nl_set_cen_lon(1,-105.)
208    CALL nl_set_truelat1(1,0.)
209    CALL nl_set_truelat2(1,0.)
210    CALL nl_set_moad_cen_lat (1,0.)
211    CALL nl_set_stand_lon (1,0.)
212    CALL nl_set_map_proj(1,0)
213
214
215!  here we initialize data we currently is not initialized
216!  in the input data
217
218    DO j = jts, jte
219      DO i = its, ite
220         grid%msft(i,j)     = 1.
221         grid%msfu(i,j)     = 1.
222         grid%msfv(i,j)     = 1.
223         grid%sina(i,j)     = 0.
224         grid%cosa(i,j)     = 1.
225         grid%e(i,j)        = 0.
226         grid%f(i,j)        = 0.
227
228      END DO
229   END DO
230
231    DO j = jts, jte
232    DO k = kts, kte
233      DO i = its, ite
234         grid%em_ww(i,k,j)     = 0.
235      END DO
236   END DO
237   END DO
238
239   grid%step_number = 0
240
241! set up the grid
242
243   IF (stretch_grid) THEN ! exponential stretch for eta (nearly constant dz)
244     DO k=1, kde
245      grid%em_znw(k) = (exp(-(k-1)/float(kde-1)/z_scale) - exp(-1./z_scale))/ &
246                                (1.-exp(-1./z_scale))
247     ENDDO
248   ELSE
249     DO k=1, kde
250      grid%em_znw(k) = 1. - float(k-1)/float(kde-1)
251     ENDDO
252   ENDIF
253
254   DO k=1, kde-1
255    grid%em_dnw(k) = grid%em_znw(k+1) - grid%em_znw(k)
256    grid%em_rdnw(k) = 1./grid%em_dnw(k)
257    grid%em_znu(k) = 0.5*(grid%em_znw(k+1)+grid%em_znw(k))
258   ENDDO
259   DO k=2, kde-1
260    grid%em_dn(k) = 0.5*(grid%em_dnw(k)+grid%em_dnw(k-1))
261    grid%em_rdn(k) = 1./grid%em_dn(k)
262    grid%em_fnp(k) = .5* grid%em_dnw(k  )/grid%em_dn(k)
263    grid%em_fnm(k) = .5* grid%em_dnw(k-1)/grid%em_dn(k)
264   ENDDO
265
266   cof1 = (2.*grid%em_dn(2)+grid%em_dn(3))/(grid%em_dn(2)+grid%em_dn(3))*grid%em_dnw(1)/grid%em_dn(2)
267   cof2 =     grid%em_dn(2)        /(grid%em_dn(2)+grid%em_dn(3))*grid%em_dnw(1)/grid%em_dn(3)
268   grid%cf1  = grid%em_fnp(2) + cof1
269   grid%cf2  = grid%em_fnm(2) - cof1 - cof2
270   grid%cf3  = cof2       
271
272   grid%cfn  = (.5*grid%em_dnw(kde-1)+grid%em_dn(kde-1))/grid%em_dn(kde-1)
273   grid%cfn1 = -.5*grid%em_dnw(kde-1)/grid%em_dn(kde-1)
274   grid%rdx = 1./config_flags%dx
275   grid%rdy = 1./config_flags%dy
276
277!  get the sounding from the ascii sounding file, first get dry sounding and
278!  calculate base state
279
280  write(6,*) ' getting dry sounding for base state '
281  dry_sounding = .true.
282  CALL get_sounding( zk, p_in, pd_in, theta, rho, u, v, qv, dry_sounding, &
283                     nl_max, nl_in, .true.)
284
285  write(6,*) ' returned from reading sounding, nl_in is ',nl_in
286
287
288!  find ptop for the desired ztop (ztop is input from the namelist),
289!  and find surface pressure
290
291  grid%p_top = interp_0( p_in, zk, config_flags%ztop, nl_in )
292
293  DO j=jts,jte
294  DO i=its,ite  ! flat surface
295!!    grid%ht(i,j) = 0.
296    grid%ht(i,j) = hm/(1.+(float(i-icm)/xa)**2)
297!    grid%ht(i,j) = hm1*exp(-(( float(i-icm)/xa1)**2))   &
298!               *( (cos(pii*float(i-icm)/xal1))**2 )
299    grid%em_phb(i,1,j) = g*grid%ht(i,j)
300    grid%em_php(i,1,j) = 0.
301    grid%em_ph0(i,1,j) = grid%em_phb(i,1,j)
302  ENDDO
303  ENDDO
304
305  DO J = jts, jte
306  DO I = its, ite
307
308    p_surf = interp_0( p_in, zk, grid%em_phb(i,1,j)/g, nl_in )
309    grid%em_mub(i,j) = p_surf-grid%p_top
310
311!  this is dry hydrostatic sounding (base state), so given p (coordinate),
312!  interp theta (from interp) and compute 1/rho from eqn. of state
313
314    DO K = 1, kte-1
315      p_level = grid%em_znu(k)*(p_surf - grid%p_top) + grid%p_top
316      grid%em_pb(i,k,j) = p_level
317      grid%em_t_init(i,k,j) = interp_0( theta, p_in, p_level, nl_in ) - t0
318      grid%em_alb(i,k,j) = (r_d/p1000mb)*(grid%em_t_init(i,k,j)+t0)*(grid%em_pb(i,k,j)/p1000mb)**cvpm
319    ENDDO
320
321!  calc hydrostatic balance (alternatively we could interp the geopotential from the
322!  sounding, but this assures that the base state is in exact hydrostatic balance with
323!  respect to the model eqns.
324
325    DO k  = 2,kte
326      grid%em_phb(i,k,j) = grid%em_phb(i,k-1,j) - grid%em_dnw(k-1)*grid%em_mub(i,j)*grid%em_alb(i,k-1,j)
327    ENDDO
328
329  ENDDO
330  ENDDO
331
332  write(6,*) ' ptop is ',grid%p_top
333  write(6,*) ' base state mub(1,1), p_surf is ',grid%em_mub(1,1),grid%em_mub(1,1)+grid%p_top
334
335!  calculate full state for each column - this includes moisture.
336
337  write(6,*) ' getting moist sounding for full state '
338  dry_sounding = .false.
339  CALL get_sounding( zk, p_in, pd_in, theta, rho, u, v, qv, dry_sounding, &
340                     nl_max, nl_in, .false. )
341
342  DO J = jts, min(jde-1,jte)
343  DO I = its, min(ide-1,ite)
344
345!  At this point grid%p_top is already set. find the DRY mass in the column
346!  by interpolating the DRY pressure. 
347
348   pd_surf = interp_0( pd_in, zk, grid%em_phb(i,1,j)/g, nl_in )
349
350!  compute the perturbation mass and the full mass
351
352    grid%em_mu_1(i,j) = pd_surf-grid%p_top - grid%em_mub(i,j)
353    grid%em_mu_2(i,j) = grid%em_mu_1(i,j)
354    grid%em_mu0(i,j) = grid%em_mu_1(i,j) + grid%em_mub(i,j)
355
356! given the dry pressure and coordinate system, interp the potential
357! temperature and qv
358
359    do k=1,kde-1
360
361      p_level = grid%em_znu(k)*(pd_surf - grid%p_top) + grid%p_top
362
363      grid%moist(i,k,j,P_QV) = interp_0( qv, pd_in, p_level, nl_in )
364      grid%em_t_1(i,k,j)          = interp_0( theta, pd_in, p_level, nl_in ) - t0
365      grid%em_t_2(i,k,j)          = grid%em_t_1(i,k,j)
366     
367
368    enddo
369
370!  integrate the hydrostatic equation (from the RHS of the bigstep
371!  vertical momentum equation) down from the top to get p.
372!  first from the top of the model to the top pressure
373
374    k = kte-1  ! top level
375
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
379
380!    grid%em_p(i,k,j) = - 0.5*grid%em_mu_1(i,j)/grid%em_rdnw(k)
381    grid%em_p(i,k,j) = - 0.5*(grid%em_mu_1(i,j)+qvf1*grid%em_mub(i,j))/grid%em_rdnw(k)/qvf2
382    qvf = 1. + rvovrd*moist(i,k,j,P_QV)
383    grid%em_alt(i,k,j) = (r_d/p1000mb)*(grid%em_t_1(i,k,j)+t0)*qvf* &
384                (((grid%em_p(i,k,j)+grid%em_pb(i,k,j))/p1000mb)**cvpm)
385    grid%em_al(i,k,j) = grid%em_alt(i,k,j) - grid%em_alb(i,k,j)
386
387!  down the column
388
389    do k=kte-2,1,-1
390      qvf1 = 0.5*(moist(i,k,j,P_QV)+moist(i,k+1,j,P_QV))
391      qvf2 = 1./(1.+qvf1)
392      qvf1 = qvf1*qvf2
393      grid%em_p(i,k,j) = grid%em_p(i,k+1,j) - (grid%em_mu_1(i,j) + qvf1*grid%em_mub(i,j))/qvf2/grid%em_rdn(k+1)
394      qvf = 1. + rvovrd*moist(i,k,j,P_QV)
395      grid%em_alt(i,k,j) = (r_d/p1000mb)*(grid%em_t_1(i,k,j)+t0)*qvf* &
396                  (((grid%em_p(i,k,j)+grid%em_pb(i,k,j))/p1000mb)**cvpm)
397      grid%em_al(i,k,j) = grid%em_alt(i,k,j) - grid%em_alb(i,k,j)
398    enddo
399
400!  this is the hydrostatic equation used in the model after the
401!  small timesteps.  In the model, al (inverse density)
402!  is computed from the geopotential.
403
404
405    grid%em_ph_1(i,1,j) = 0.
406    DO k  = 2,kte
407      grid%em_ph_1(i,k,j) = grid%em_ph_1(i,k-1,j) - (1./grid%em_rdnw(k-1))*(       &
408                   (grid%em_mub(i,j)+grid%em_mu_1(i,j))*grid%em_al(i,k-1,j)+ &
409                    grid%em_mu_1(i,j)*grid%em_alb(i,k-1,j)  )
410                                                   
411      grid%em_ph_2(i,k,j) = grid%em_ph_1(i,k,j)
412      grid%em_ph0(i,k,j) = grid%em_ph_1(i,k,j) + grid%em_phb(i,k,j)
413    ENDDO
414
415    if((i==2) .and. (j==2)) then
416     write(6,*) ' ph_1 calc ',grid%em_ph_1(2,1,2),grid%em_ph_1(2,2,2),&
417                              grid%em_mu_1(2,2)+grid%em_mub(2,2),grid%em_mu_1(2,2), &
418                              grid%em_alb(2,1,2),grid%em_al(1,2,1),grid%em_rdnw(1)
419    endif
420
421  ENDDO
422  ENDDO
423
424!  cold bubble input  (from straka et al, IJNMF, vol 17, 1993 pp 1-22)
425
426  t_min = grid%em_t_1(its,kts,jts)
427  t_max = t_min
428  u_mean = 00.
429
430  xpos = config_flags%dx*nxc - u_mean*900.
431  xposml = xpos - config_flags%dx*(ide-1)
432  xpospl = xpos + config_flags%dx*(ide-1)
433
434  DO J = jts, min(jde-1,jte)
435    DO I = its, min(ide-1,ite)
436!      xrad = config_flags%dx*float(i-nxc)/4000.  !  4000 meter horizontal radius
437!                                    !  centered in the domain
438
439       xrad = min( abs(config_flags%dx*float(i)-xpos),   &
440                   abs(config_flags%dx*float(i)-xposml), &
441                   abs(config_flags%dx*float(i)-xpospl))/4000.
442
443      DO K = 1, kte-1
444
445!  put in preturbation theta (bubble) and recalc density.  note,
446!  the mass in the column is not changing, so when theta changes,
447!  we recompute density and geopotential
448
449        zrad = 0.5*(grid%em_ph_1(i,k,j)+grid%em_ph_1(i,k+1,j)  &
450                   +grid%em_phb(i,k,j)+grid%em_phb(i,k+1,j))/g
451        zrad = (zrad-3000.)/2000. !  2000 meter vertical radius,
452                                  !  centered at z=3000,
453        RAD=SQRT(xrad*xrad+zrad*zrad)
454        IF(RAD <= 1.) THEN
455
456           !  perturbation temperature is 15 C, convert to potential temperature
457
458           delt = -15.0 / ((grid%em_p(i,k,j)+grid%em_pb(i,k,j))/p1000mb)**rcp
459
460           grid%em_T_1(i,k,j)=grid%em_T_1(i,k,j)+delt*(COS(PI*RAD)+1.0)/2.
461           grid%em_T_2(i,k,j)=grid%em_T_1(i,k,j)
462           qvf = 1. + rvovrd*moist(i,k,j,P_QV)
463           grid%em_alt(i,k,j) = (r_d/p1000mb)*(grid%em_t_1(i,k,j)+t0)*qvf* &
464                        (((grid%em_p(i,k,j)+grid%em_pb(i,k,j))/p1000mb)**cvpm)
465           grid%em_al(i,k,j) = grid%em_alt(i,k,j) - grid%em_alb(i,k,j)
466        ENDIF
467
468        t_min = min(t_min, grid%em_t_1(i,k,j))
469        t_max = max(t_max, grid%em_t_1(i,k,j))
470      ENDDO
471
472!  rebalance hydrostatically
473
474      DO k  = 2,kte
475        grid%em_ph_1(i,k,j) = grid%em_ph_1(i,k-1,j) - (1./grid%em_rdnw(k-1))*(       &
476                     (grid%em_mub(i,j)+grid%em_mu_1(i,j))*grid%em_al(i,k-1,j)+ &
477                      grid%em_mu_1(i,j)*grid%em_alb(i,k-1,j)  )
478                                                   
479        grid%em_ph_2(i,k,j) = grid%em_ph_1(i,k,j)
480        grid%em_ph0(i,k,j) = grid%em_ph_1(i,k,j) + grid%em_phb(i,k,j)
481      ENDDO
482
483    ENDDO
484  ENDDO
485
486  write(6,*) ' min and max theta perturbation ',t_min,t_max
487
488
489
490
491! -- end bubble insert
492
493   write(6,*) ' mu_1 from comp ', grid%em_mu_1(1,1)
494   write(6,*) ' full state sounding from comp, ph, p, al, t_1, qv '
495   do k=1,kde-1
496     write(6,'(i3,1x,5(1x,1pe10.3))') k, grid%em_ph_1(1,k,1)+grid%em_phb(1,k,1), &
497                                      grid%em_p(1,k,1)+grid%em_pb(1,k,1), grid%em_alt(1,k,1), &
498                                      grid%em_t_1(1,k,1)+t0, moist(1,k,1,P_QV)
499   enddo
500
501   write(6,*) ' pert state sounding from comp, ph_1, pp, alp, t_1, qv '
502   do k=1,kde-1
503     write(6,'(i3,1x,5(1x,1pe10.3))') k, grid%em_ph_1(1,k,1), &
504                                      grid%em_p(1,k,1), grid%em_al(1,k,1), &
505                                      grid%em_t_1(1,k,1), moist(1,k,1,P_QV)
506   enddo
507
508   write(6,*) ' '
509   write(6,*) ' k, model level, dz '
510   do k=1,kde-1
511     write(6,'(i3,1x,e12.5,1x,f10.2)') k,  &
512      .5*(grid%em_ph_1(1,k,1)+grid%em_phb(1,k,1)+grid%em_ph_1(1,k+1,1)+grid%em_phb(1,k+1,1))/g, &
513      (grid%em_ph_1(1,k+1,1)+grid%em_phb(1,k+1,1)-grid%em_ph_1(1,k,1)-grid%em_phb(1,k,1))/g
514   enddo
515   write(6,*) ' model top (m) is ', (grid%em_ph_1(1,kde,1)+grid%em_phb(1,kde,1))/g
516
517
518! interp v
519
520  DO J = jts, jte
521  DO I = its, min(ide-1,ite)
522
523    IF (j == jds) THEN
524      z_at_v = grid%em_phb(i,1,j)/g
525    ELSE IF (j == jde) THEN
526      z_at_v = grid%em_phb(i,1,j-1)/g
527    ELSE
528      z_at_v = 0.5*(grid%em_phb(i,1,j)+grid%em_phb(i,1,j-1))/g
529    END IF
530
531    p_surf = interp_0( p_in, zk, z_at_v, nl_in )
532
533    DO K = 1, kte
534      p_level = grid%em_znu(k)*(p_surf - grid%p_top) + grid%p_top
535      grid%em_v_1(i,k,j) = interp_0( v, p_in, p_level, nl_in )
536      grid%em_v_2(i,k,j) = grid%em_v_1(i,k,j)
537    ENDDO
538
539  ENDDO
540  ENDDO
541
542! interp u
543
544  DO J = jts, min(jde-1,jte)
545  DO I = its, ite
546
547    IF (i == ids) THEN
548      z_at_u = grid%em_phb(i,1,j)/g
549    ELSE IF (i == ide) THEN
550      z_at_u = grid%em_phb(i-1,1,j)/g
551    ELSE
552      z_at_u = 0.5*(grid%em_phb(i,1,j)+grid%em_phb(i-1,1,j))/g
553    END IF
554
555    p_surf = interp_0( p_in, zk, z_at_u, nl_in )
556
557    DO K = 1, kte
558      p_level = grid%em_znu(k)*(p_surf - grid%p_top) + grid%p_top
559      grid%em_u_1(i,k,j) = interp_0( u, p_in, p_level, nl_in )
560      grid%em_u_2(i,k,j) = grid%em_u_1(i,k,j)
561    ENDDO
562
563  ENDDO
564  ENDDO
565
566!  set w
567
568  DO J = jts, min(jde-1,jte)
569  DO K = kts, kte
570  DO I = its, min(ide-1,ite)
571    grid%em_w_1(i,k,j) = 0.
572    grid%em_w_2(i,k,j) = 0.
573  ENDDO
574  ENDDO
575  ENDDO
576
577!  set a few more things
578
579  DO J = jts, min(jde-1,jte)
580  DO K = kts, kte-1
581  DO I = its, min(ide-1,ite)
582    grid%h_diabatic(i,k,j) = 0.
583  ENDDO
584  ENDDO
585  ENDDO
586
587  DO k=1,kte-1
588    grid%em_t_base(k) = grid%em_t_1(1,k,1)
589    grid%qv_base(k) = moist(1,k,1,P_QV)
590    grid%u_base(k) = grid%em_u_1(1,k,1)
591    grid%v_base(k) = grid%em_v_1(1,k,1)
592    grid%z_base(k) = 0.5*(grid%em_phb(1,k,1)+grid%em_phb(1,k+1,1)+grid%em_ph_1(1,k,1)+grid%em_ph_1(1,k+1,1))/g
593  ENDDO
594
595  DO J = jts, min(jde-1,jte)
596  DO I = its, min(ide-1,ite)
597     thtmp   = grid%em_t_2(i,1,j)+t0
598     ptmp    = grid%em_p(i,1,j)+grid%em_pb(i,1,j)
599     temp(1) = thtmp * (ptmp/p1000mb)**rcp
600     thtmp   = grid%em_t_2(i,2,j)+t0
601     ptmp    = grid%em_p(i,2,j)+grid%em_pb(i,2,j)
602     temp(2) = thtmp * (ptmp/p1000mb)**rcp
603     thtmp   = grid%em_t_2(i,3,j)+t0
604     ptmp    = grid%em_p(i,3,j)+grid%em_pb(i,3,j)
605     temp(3) = thtmp * (ptmp/p1000mb)**rcp
606
607     grid%TSK(I,J)=grid%cf1*temp(1)+grid%cf2*temp(2)+grid%cf3*temp(3)
608     grid%TMN(I,J)=grid%TSK(I,J)-0.5
609  ENDDO
610  ENDDO
611
612  RETURN
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, base_state )
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      logical base_state
652
653      integer n, iz
654      parameter(n=1000)
655      logical debug
656      parameter( debug = .false.)
657
658! input sounding data
659
660      real p_surf, th_surf, qv_surf
661      real pi_surf, pi(n)
662      real h_input(n), th_input(n), qv_input(n), u_input(n), v_input(n)
663
664! diagnostics
665
666      real rho_surf, p_input(n), rho_input(n)
667      real pm_input(n)  !  this are for full moist sounding
668
669! local data
670
671      real p1000mb,cv,cp,r,cvpm,g
672      parameter (p1000mb = 1.e+05, r = 287, cp = 1003., cv = cp-r, cvpm = -cv/cp, g=9.81 )
673      integer k, it, nl
674      real qvf, qvf1, dz
675
676!  first, read the sounding
677
678      call read_sounding( p_surf, th_surf, qv_surf, &
679                          h_input, th_input, qv_input, u_input, v_input,n, nl, debug )
680
681!        iz = 1
682!        do k=2,nl
683!          if(h_input(k) .lt. 12000.) iz = k
684!        enddo
685!        write(6,*) " tropopause ",iz,h_input(iz)
686!        if(dry) then
687!        write(6,*) ' nl is ',nl
688!        do k=1,nl
689!          th_input(k) = th_input(k)+10.+10*float(k)/nl
690!        enddo
691!        write(6,*) ' finished adjusting theta '
692!        endif
693
694!        do k=1,nl
695!          u_input(k) = 2*u_input(k)
696!        enddo
697!
698!      end if
699
700      if(dry) then
701       do k=1,nl
702         qv_input(k) = 0.
703       enddo
704      endif
705
706      if(debug) write(6,*) ' number of input levels = ',nl
707
708        nl_in = nl
709        if(nl_in .gt. nl_max ) then
710          write(6,*) ' too many levels for input arrays ',nl_in,nl_max
711          call wrf_error_fatal ( ' too many levels for input arrays ' )
712        end if
713
714!  compute diagnostics,
715!  first, convert qv(g/kg) to qv(g/g)
716
717      do k=1,nl
718        qv_input(k) = 0.001*qv_input(k)
719      enddo
720
721      p_surf = 100.*p_surf  ! convert to pascals
722      qvf = 1. + rvovrd*qv_input(1)
723      rho_surf = 1./((r/p1000mb)*th_surf*qvf*((p_surf/p1000mb)**cvpm))
724      pi_surf = (p_surf/p1000mb)**(r/cp)
725
726      if(debug) then
727        write(6,*) ' surface density is ',rho_surf
728        write(6,*) ' surface pi is      ',pi_surf
729      end if
730
731
732!  integrate moist sounding hydrostatically, starting from the
733!  specified surface pressure
734!  -> first, integrate from surface to lowest level
735
736          qvf = 1. + rvovrd*qv_input(1)
737          qvf1 = 1. + qv_input(1)
738          rho_input(1) = rho_surf
739          dz = h_input(1)
740          do it=1,10
741            pm_input(1) = p_surf &
742                    - 0.5*dz*(rho_surf+rho_input(1))*g*qvf1
743            rho_input(1) = 1./((r/p1000mb)*th_input(1)*qvf*((pm_input(1)/p1000mb)**cvpm))
744          enddo
745
746! integrate up the column
747
748          do k=2,nl
749            rho_input(k) = rho_input(k-1)
750            dz = h_input(k)-h_input(k-1)
751            qvf1 = 0.5*(2.+(qv_input(k-1)+qv_input(k)))
752            qvf = 1. + rvovrd*qv_input(k)   ! qv is in g/kg here
753 
754            do it=1,20
755              pm_input(k) = pm_input(k-1) &
756                      - 0.5*dz*(rho_input(k)+rho_input(k-1))*g*qvf1
757              rho_input(k) = 1./((r/p1000mb)*th_input(k)*qvf*((pm_input(k)/p1000mb)**cvpm))
758            enddo
759          enddo
760
761!  we have the moist sounding
762
763!  next, compute the dry sounding using p at the highest level from the
764!  moist sounding and integrating down.
765
766        p_input(nl) = pm_input(nl)
767
768          do k=nl-1,1,-1
769            dz = h_input(k+1)-h_input(k)
770            p_input(k) = p_input(k+1) + 0.5*dz*(rho_input(k)+rho_input(k+1))*g
771          enddo
772
773!      write(6,*) ' zeroing u input '
774
775        do k=1,nl
776
777          zk(k) = h_input(k)
778          p(k) = pm_input(k)
779          p_dry(k) = p_input(k)
780          theta(k) = th_input(k)
781          rho(k) = rho_input(k)
782          u(k) = u_input(k)
783!          u(k) = 0.
784          v(k) = v_input(k)
785          qv(k) = qv_input(k)
786
787        enddo
788
789     if(debug) then
790      write(6,*) ' sounding '
791      write(6,*) '  k  height(m)  press (Pa) pd(Pa) theta (K) den(kg/m^3)  u(m/s)     v(m/s)    qv(g/g) '
792      do k=1,nl
793        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)
794      enddo
795
796     end if
797
798      end subroutine get_sounding
799
800!-------------------------------------------------------
801
802      subroutine read_sounding( ps,ts,qvs,h,th,qv,u,v,n,nl,debug )
803      implicit none
804      integer n,nl
805      real ps,ts,qvs,h(n),th(n),qv(n),u(n),v(n)
806      logical end_of_file
807      logical debug
808
809      integer k
810
811      open(unit=10,file='input_sounding',form='formatted',status='old')
812      rewind(10)
813      read(10,*) ps, ts, qvs
814      if(debug) then
815        write(6,*) ' input sounding surface parameters '
816        write(6,*) ' surface pressure (mb) ',ps
817        write(6,*) ' surface pot. temp (K) ',ts
818        write(6,*) ' surface mixing ratio (g/kg) ',qvs
819      end if
820
821      end_of_file = .false.
822      k = 0
823
824      do while (.not. end_of_file)
825
826        read(10,*,end=100) h(k+1), th(k+1), qv(k+1), u(k+1), v(k+1)
827        k = k+1
828        if(debug) write(6,'(1x,i3,5(1x,e10.3))') k, h(k), th(k), qv(k), u(k), v(k)
829        go to 110
830 100    end_of_file = .true.
831 110    continue
832      enddo
833
834      nl = k
835
836      close(unit=10,status = 'keep')
837
838      end subroutine read_sounding
839
840END MODULE module_initialize
Note: See TracBrowser for help on using the repository browser.