source: trunk/WRF.COMMON/WRFV3/dyn_em/module_initialize_les.F @ 2759

Last change on this file since 2759 was 2759, checked in by aslmd, 2 years ago

adding unmodified code from WRFV3.0.1.1, expurged from useless data +1M size

File size: 24.1 KB
RevLine 
[2759]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!  For LES, add randx
112   real :: randx
113
114#ifdef DM_PARALLEL
115#    include <data_calls.inc>
116#endif
117
118
119   SELECT CASE ( model_data_order )
120         CASE ( DATA_ORDER_ZXY )
121   kds = grid%sd31 ; kde = grid%ed31 ;
122   ids = grid%sd32 ; ide = grid%ed32 ;
123   jds = grid%sd33 ; jde = grid%ed33 ;
124
125   kms = grid%sm31 ; kme = grid%em31 ;
126   ims = grid%sm32 ; ime = grid%em32 ;
127   jms = grid%sm33 ; jme = grid%em33 ;
128
129   kts = grid%sp31 ; kte = grid%ep31 ;   ! note that tile is entire patch
130   its = grid%sp32 ; ite = grid%ep32 ;   ! note that tile is entire patch
131   jts = grid%sp33 ; jte = grid%ep33 ;   ! note that tile is entire patch
132         CASE ( DATA_ORDER_XYZ )
133   ids = grid%sd31 ; ide = grid%ed31 ;
134   jds = grid%sd32 ; jde = grid%ed32 ;
135   kds = grid%sd33 ; kde = grid%ed33 ;
136
137   ims = grid%sm31 ; ime = grid%em31 ;
138   jms = grid%sm32 ; jme = grid%em32 ;
139   kms = grid%sm33 ; kme = grid%em33 ;
140
141   its = grid%sp31 ; ite = grid%ep31 ;   ! note that tile is entire patch
142   jts = grid%sp32 ; jte = grid%ep32 ;   ! note that tile is entire patch
143   kts = grid%sp33 ; kte = grid%ep33 ;   ! note that tile is entire patch
144         CASE ( DATA_ORDER_XZY )
145   ids = grid%sd31 ; ide = grid%ed31 ;
146   kds = grid%sd32 ; kde = grid%ed32 ;
147   jds = grid%sd33 ; jde = grid%ed33 ;
148
149   ims = grid%sm31 ; ime = grid%em31 ;
150   kms = grid%sm32 ; kme = grid%em32 ;
151   jms = grid%sm33 ; jme = grid%em33 ;
152
153   its = grid%sp31 ; ite = grid%ep31 ;   ! note that tile is entire patch
154   kts = grid%sp32 ; kte = grid%ep32 ;   ! note that tile is entire patch
155   jts = grid%sp33 ; jte = grid%ep33 ;   ! note that tile is entire patch
156
157   END SELECT
158
159
160!  stretch_grid = .true.
161!  FOR LES, set stretch to false
162   stretch_grid = .false.
163   delt = 3.
164!   z_scale = .50
165   z_scale = .40
166   pi = 2.*asin(1.0)
167   write(6,*) ' pi is ',pi
168   nxc = (ide-ids)/2
169   nyc = (jde-jds)/2
170
171   CALL model_to_grid_config_rec ( grid%id , model_config_rec , config_flags )
172
173! here we check to see if the boundary conditions are set properly
174
175   CALL boundary_condition_check( config_flags, bdyzone, error, grid%id )
176
177   moisture_init = .true.
178
179    grid%itimestep=0
180
181#ifdef DM_PARALLEL
182   CALL wrf_dm_bcast_bytes( icm , IWORDSIZE )
183   CALL wrf_dm_bcast_bytes( jcm , IWORDSIZE )
184#endif
185
186    CALL nl_set_mminlu(1, '    ')
187    CALL nl_set_iswater(1,0)
188    CALL nl_set_cen_lat(1,40.)
189    CALL nl_set_cen_lon(1,-105.)
190    CALL nl_set_truelat1(1,0.)
191    CALL nl_set_truelat2(1,0.)
192    CALL nl_set_moad_cen_lat (1,0.)
193    CALL nl_set_stand_lon (1,0.)
194    CALL nl_set_map_proj(1,0)
195
196
197!  here we initialize data we currently is not initialized
198!  in the input data
199
200    DO j = jts, jte
201      DO i = its, ite
202         grid%msftx(i,j)    = 1.
203         grid%msfty(i,j)    = 1.
204         grid%msfux(i,j)    = 1.
205         grid%msfuy(i,j)    = 1.
206         grid%msfvx(i,j)    = 1.
207         grid%msfvx_inv(i,j)= 1.
208         grid%msfvy(i,j)    = 1.
209         grid%sina(i,j)     = 0.
210         grid%cosa(i,j)     = 1.
211         grid%e(i,j)        = 0.
212!  for LES, include Coriolis force
213         grid%f(i,j)        = 1.e-4
214
215      END DO
216   END DO
217
218    DO j = jts, jte
219    DO k = kts, kte
220      DO i = its, ite
221         grid%ww(i,k,j)     = 0.
222      END DO
223   END DO
224   END DO
225
226   grid%step_number = 0
227
228! set up the grid
229
230   IF (stretch_grid) THEN ! exponential stretch for eta (nearly constant dz)
231     DO k=1, kde
232      grid%znw(k) = (exp(-(k-1)/float(kde-1)/z_scale) - exp(-1./z_scale))/ &
233                                (1.-exp(-1./z_scale))
234     ENDDO
235   ELSE
236     DO k=1, kde
237      grid%znw(k) = 1. - float(k-1)/float(kde-1)
238     ENDDO
239   ENDIF
240
241   DO k=1, kde-1
242    grid%dnw(k) = grid%znw(k+1) - grid%znw(k)
243    grid%rdnw(k) = 1./grid%dnw(k)
244    grid%znu(k) = 0.5*(grid%znw(k+1)+grid%znw(k))
245   ENDDO
246   DO k=2, kde-1
247    grid%dn(k) = 0.5*(grid%dnw(k)+grid%dnw(k-1))
248    grid%rdn(k) = 1./grid%dn(k)
249    grid%fnp(k) = .5* grid%dnw(k  )/grid%dn(k)
250    grid%fnm(k) = .5* grid%dnw(k-1)/grid%dn(k)
251   ENDDO
252
253   cof1 = (2.*grid%dn(2)+grid%dn(3))/(grid%dn(2)+grid%dn(3))*grid%dnw(1)/grid%dn(2)
254   cof2 =     grid%dn(2)        /(grid%dn(2)+grid%dn(3))*grid%dnw(1)/grid%dn(3)
255   grid%cf1  = grid%fnp(2) + cof1
256   grid%cf2  = grid%fnm(2) - cof1 - cof2
257   grid%cf3  = cof2       
258
259   grid%cfn  = (.5*grid%dnw(kde-1)+grid%dn(kde-1))/grid%dn(kde-1)
260   grid%cfn1 = -.5*grid%dnw(kde-1)/grid%dn(kde-1)
261   grid%rdx = 1./config_flags%dx
262   grid%rdy = 1./config_flags%dy
263
264!  get the sounding from the ascii sounding file, first get dry sounding and
265!  calculate base state
266
267  dry_sounding = .true.
268  IF ( wrf_dm_on_monitor() ) THEN
269  write(6,*) ' getting dry sounding for base state '
270
271  CALL get_sounding( zk, p_in, pd_in, theta, rho, u, v, qv, dry_sounding, nl_max, nl_in )
272  ENDIF
273  CALL wrf_dm_bcast_real( zk , nl_max )
274  CALL wrf_dm_bcast_real( p_in , nl_max )
275  CALL wrf_dm_bcast_real( pd_in , nl_max )
276  CALL wrf_dm_bcast_real( theta , nl_max )
277  CALL wrf_dm_bcast_real( rho , nl_max )
278  CALL wrf_dm_bcast_real( u , nl_max )
279  CALL wrf_dm_bcast_real( v , nl_max )
280  CALL wrf_dm_bcast_real( qv , nl_max )
281  CALL wrf_dm_bcast_integer ( nl_in , 1 )
282
283  write(6,*) ' returned from reading sounding, nl_in is ',nl_in
284
285!  find ptop for the desired ztop (ztop is input from the namelist),
286!  and find surface pressure
287
288  grid%p_top = interp_0( p_in, zk, config_flags%ztop, nl_in )
289
290  DO j=jts,jte
291  DO i=its,ite
292    grid%ht(i,j) = 0.
293  ENDDO
294  ENDDO
295
296  xs=ide/2 -3
297  xs=ids   -3
298  xe=xs + 6
299  ys=jde/2 -3
300  ye=ys + 6
301  mtn_ht = 500
302#ifdef MTN
303  DO j=max(ys,jds),min(ye,jde-1)
304  DO i=max(xs,ids),min(xe,ide-1)
305     grid%ht(i,j) = mtn_ht * 0.25 * &
306               ( 1. + COS ( 2*pi/(xe-xs) * ( i-xs ) + pi ) ) * &
307               ( 1. + COS ( 2*pi/(ye-ys) * ( j-ys ) + pi ) )
308  ENDDO
309  ENDDO
310#endif
311#ifdef EW_RIDGE
312  DO j=max(ys,jds),min(ye,jde-1)
313  DO i=ids,ide
314     grid%ht(i,j) = mtn_ht * 0.50 * &
315               ( 1. + COS ( 2*pi/(ye-ys) * ( j-ys ) + pi ) )
316  ENDDO
317  ENDDO
318#endif
319#ifdef NS_RIDGE
320  DO j=jds,jde
321  DO i=max(xs,ids),min(xe,ide-1)
322     grid%ht(i,j) = mtn_ht * 0.50 * &
323               ( 1. + COS ( 2*pi/(xe-xs) * ( i-xs ) + pi ) )
324  ENDDO
325  ENDDO
326#endif
327  DO j=jts,jte
328  DO i=its,ite
329    grid%phb(i,1,j) = g * grid%ht(i,j)
330    grid%ph0(i,1,j) = g * grid%ht(i,j)
331  ENDDO
332  ENDDO
333
334  DO J = jts, jte
335  DO I = its, ite
336
337    p_surf = interp_0( p_in, zk, grid%phb(i,1,j)/g, nl_in )
338    grid%mub(i,j) = p_surf-grid%p_top
339
340!  this is dry hydrostatic sounding (base state), so given grid%p (coordinate),
341!  interp theta (from interp) and compute 1/rho from eqn. of state
342
343    DO K = 1, kte-1
344      p_level = grid%znu(k)*(p_surf - grid%p_top) + grid%p_top
345      grid%pb(i,k,j) = p_level
346      grid%t_init(i,k,j) = interp_0( theta, p_in, p_level, nl_in ) - t0
347      grid%alb(i,k,j) = (r_d/p1000mb)*(grid%t_init(i,k,j)+t0)*(grid%pb(i,k,j)/p1000mb)**cvpm
348    ENDDO
349
350!  calc hydrostatic balance (alternatively we could interp the geopotential from the
351!  sounding, but this assures that the base state is in exact hydrostatic balance with
352!  respect to the model eqns.
353
354    DO k  = 2,kte
355      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)
356    ENDDO
357
358  ENDDO
359  ENDDO
360
361  IF ( wrf_dm_on_monitor() ) THEN
362    write(6,*) ' ptop is ',grid%p_top
363    write(6,*) ' base state grid%mub(1,1), p_surf is ',grid%mub(1,1),grid%mub(1,1)+grid%p_top
364  ENDIF
365
366!  calculate full state for each column - this includes moisture.
367
368  write(6,*) ' getting moist sounding for full state '
369  dry_sounding = .false.
370  CALL get_sounding( zk, p_in, pd_in, theta, rho, u, v, qv, dry_sounding, nl_max, nl_in )
371
372  DO J = jts, min(jde-1,jte)
373  DO I = its, min(ide-1,ite)
374
375!  At this point grid%p_top is already set. find the DRY mass in the column
376!  by interpolating the DRY pressure. 
377
378   pd_surf = interp_0( pd_in, zk, grid%phb(i,1,j)/g, nl_in )
379
380!  compute the perturbation mass and the full mass
381
382    grid%mu_1(i,j) = pd_surf-grid%p_top - grid%mub(i,j)
383    grid%mu_2(i,j) = grid%mu_1(i,j)
384    grid%mu0(i,j) = grid%mu_1(i,j) + grid%mub(i,j)
385
386! given the dry pressure and coordinate system, interp the potential
387! temperature and qv
388
389    do k=1,kde-1
390
391      p_level = grid%znu(k)*(pd_surf - grid%p_top) + grid%p_top
392
393      moist(i,k,j,P_QV) = interp_0( qv, pd_in, p_level, nl_in )
394      grid%t_1(i,k,j)          = interp_0( theta, pd_in, p_level, nl_in ) - t0
395      grid%t_2(i,k,j)          = grid%t_1(i,k,j)
396     
397
398    enddo
399
400!  integrate the hydrostatic equation (from the RHS of the bigstep
401!  vertical momentum equation) down from the top to get grid%p.
402!  first from the top of the model to the top pressure
403
404    k = kte-1  ! top level
405
406    qvf1 = 0.5*(moist(i,k,j,P_QV)+moist(i,k,j,P_QV))
407    qvf2 = 1./(1.+qvf1)
408    qvf1 = qvf1*qvf2
409
410!    grid%p(i,k,j) = - 0.5*grid%mu_1(i,j)/grid%rdnw(k)
411    grid%p(i,k,j) = - 0.5*(grid%mu_1(i,j)+qvf1*grid%mub(i,j))/grid%rdnw(k)/qvf2
412    qvf = 1. + rvovrd*moist(i,k,j,P_QV)
413    grid%alt(i,k,j) = (r_d/p1000mb)*(grid%t_1(i,k,j)+t0)*qvf* &
414                (((grid%p(i,k,j)+grid%pb(i,k,j))/p1000mb)**cvpm)
415    grid%al(i,k,j) = grid%alt(i,k,j) - grid%alb(i,k,j)
416
417!  down the column
418
419    do k=kte-2,1,-1
420      qvf1 = 0.5*(moist(i,k,j,P_QV)+moist(i,k+1,j,P_QV))
421      qvf2 = 1./(1.+qvf1)
422      qvf1 = qvf1*qvf2
423      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)
424      qvf = 1. + rvovrd*moist(i,k,j,P_QV)
425      grid%alt(i,k,j) = (r_d/p1000mb)*(grid%t_1(i,k,j)+t0)*qvf* &
426                  (((grid%p(i,k,j)+grid%pb(i,k,j))/p1000mb)**cvpm)
427      grid%al(i,k,j) = grid%alt(i,k,j) - grid%alb(i,k,j)
428    enddo
429
430!  this is the hydrostatic equation used in the model after the
431!  small timesteps.  In the model, grid%al (inverse density)
432!  is computed from the geopotential.
433
434
435    grid%ph_1(i,1,j) = 0.
436    DO k  = 2,kte
437      grid%ph_1(i,k,j) = grid%ph_1(i,k-1,j) - (1./grid%rdnw(k-1))*(       &
438                   (grid%mub(i,j)+grid%mu_1(i,j))*grid%al(i,k-1,j)+ &
439                    grid%mu_1(i,j)*grid%alb(i,k-1,j)  )
440                                                   
441      grid%ph_2(i,k,j) = grid%ph_1(i,k,j)
442      grid%ph0(i,k,j) = grid%ph_1(i,k,j) + grid%phb(i,k,j)
443    ENDDO
444
445    IF ( wrf_dm_on_monitor() ) THEN
446    if((i==2) .and. (j==2)) then
447     write(6,*) ' grid%ph_1 calc ',grid%ph_1(2,1,2),grid%ph_1(2,2,2),&
448                              grid%mu_1(2,2)+grid%mub(2,2),grid%mu_1(2,2), &
449                              grid%alb(2,1,2),grid%al(1,2,1),grid%rdnw(1)
450    endif
451    ENDIF
452
453  ENDDO
454  ENDDO
455
456!#if 0
457
458!  thermal perturbation to kick off convection
459
460  write(6,*) ' nxc, nyc for perturbation ',nxc,nyc
461  write(6,*) ' delt for perturbation ',delt
462
463! For LES, change the initial random perturbations
464! For 2D test, call randx outside I-loop
465! For 3D runs, call randx inside both I-J loops
466
467  DO J = jts, min(jde-1,jte)
468!   yrad = config_flags%dy*float(j-nyc)/10000.
469    yrad = 0.
470    DO I = its, min(ide-1,ite)
471!     xrad = config_flags%dx*float(i-nxc)/10000.
472      xrad = 0.
473      call random_number (randx)
474      randx = randx - 0.5
475!     DO K = 1, kte-1
476      DO K = 1, 4
477
478!  No bubbles for LES!
479!  put in preturbation theta (bubble) and recalc density.  note,
480!  the mass in the column is not changing, so when theta changes,
481!  we recompute density and geopotential
482
483!       zrad = 0.5*(grid%ph_1(i,k,j)+grid%ph_1(i,k+1,j)  &
484!                  +grid%phb(i,k,j)+grid%phb(i,k+1,j))/g
485!       zrad = (zrad-1500.)/1500.
486        zrad = 0.
487        RAD=SQRT(xrad*xrad+yrad*yrad+zrad*zrad)
488        IF(RAD <= 1.) THEN
489!          grid%t_1(i,k,j)=grid%t_1(i,k,j)+delt*COS(.5*PI*RAD)**2
490           grid%t_1(i,k,j)=grid%t_1(i,k,j)+ 0.1 *randx
491           grid%t_2(i,k,j)=grid%t_1(i,k,j)
492           qvf = 1. + rvovrd*moist(i,k,j,P_QV)
493           grid%alt(i,k,j) = (r_d/p1000mb)*(grid%t_1(i,k,j)+t0)*qvf* &
494                        (((grid%p(i,k,j)+grid%pb(i,k,j))/p1000mb)**cvpm)
495           grid%al(i,k,j) = grid%alt(i,k,j) - grid%alb(i,k,j)
496        ENDIF
497      ENDDO
498
499!  rebalance hydrostatically
500
501      DO k  = 2,kte
502        grid%ph_1(i,k,j) = grid%ph_1(i,k-1,j) - (1./grid%rdnw(k-1))*(       &
503                     (grid%mub(i,j)+grid%mu_1(i,j))*grid%al(i,k-1,j)+ &
504                      grid%mu_1(i,j)*grid%alb(i,k-1,j)  )
505                                                   
506        grid%ph_2(i,k,j) = grid%ph_1(i,k,j)
507        grid%ph0(i,k,j) = grid%ph_1(i,k,j) + grid%phb(i,k,j)
508      ENDDO
509
510    ENDDO
511  ENDDO
512
513!#endif
514
515   IF ( wrf_dm_on_monitor() ) THEN
516   write(6,*) ' grid%mu_1 from comp ', grid%mu_1(1,1)
517   write(6,*) ' full state sounding from comp, ph, grid%p, grid%al, grid%t_1, qv '
518   do k=1,kde-1
519     write(6,'(i3,1x,5(1x,1pe10.3))') k, grid%ph_1(1,k,1)+grid%phb(1,k,1), &
520                                      grid%p(1,k,1)+grid%pb(1,k,1), grid%alt(1,k,1), &
521                                      grid%t_1(1,k,1)+t0, moist(1,k,1,P_QV)
522   enddo
523
524   write(6,*) ' pert state sounding from comp, grid%ph_1, pp, alp, grid%t_1, qv '
525   do k=1,kde-1
526     write(6,'(i3,1x,5(1x,1pe10.3))') k, grid%ph_1(1,k,1), &
527                                      grid%p(1,k,1), grid%al(1,k,1), &
528                                      grid%t_1(1,k,1), moist(1,k,1,P_QV)
529   enddo
530   ENDIF
531
532! interp v
533
534  DO J = jts, jte
535  DO I = its, min(ide-1,ite)
536
537    IF (j == jds) THEN
538      z_at_v = grid%phb(i,1,j)/g
539    ELSE IF (j == jde) THEN
540      z_at_v = grid%phb(i,1,j-1)/g
541    ELSE
542      z_at_v = 0.5*(grid%phb(i,1,j)+grid%phb(i,1,j-1))/g
543    END IF
544    p_surf = interp_0( p_in, zk, z_at_v, nl_in )
545
546    DO K = 1, kte-1
547      p_level = grid%znu(k)*(p_surf - grid%p_top) + grid%p_top
548      grid%v_1(i,k,j) = interp_0( v, p_in, p_level, nl_in )
549      grid%v_2(i,k,j) = grid%v_1(i,k,j)
550    ENDDO
551
552  ENDDO
553  ENDDO
554
555! interp u
556
557  DO J = jts, min(jde-1,jte)
558  DO I = its, ite
559
560    IF (i == ids) THEN
561      z_at_u = grid%phb(i,1,j)/g
562    ELSE IF (i == ide) THEN
563      z_at_u = grid%phb(i-1,1,j)/g
564    ELSE
565      z_at_u = 0.5*(grid%phb(i,1,j)+grid%phb(i-1,1,j))/g
566    END IF
567
568    p_surf = interp_0( p_in, zk, z_at_u, nl_in )
569
570    DO K = 1, kte-1
571      p_level = grid%znu(k)*(p_surf - grid%p_top) + grid%p_top
572      grid%u_1(i,k,j) = interp_0( u, p_in, p_level, nl_in )
573      grid%u_2(i,k,j) = grid%u_1(i,k,j)
574    ENDDO
575
576  ENDDO
577  ENDDO
578
579!  set w
580
581  DO J = jts, min(jde-1,jte)
582  DO K = kts, kte
583  DO I = its, min(ide-1,ite)
584    grid%w_1(i,k,j) = 0.
585    grid%w_2(i,k,j) = 0.
586  ENDDO
587  ENDDO
588  ENDDO
589
590!  set a few more things
591
592  DO J = jts, min(jde-1,jte)
593  DO K = kts, kte-1
594  DO I = its, min(ide-1,ite)
595    grid%h_diabatic(i,k,j) = 0.
596  ENDDO
597  ENDDO
598  ENDDO
599
600  IF ( wrf_dm_on_monitor() ) THEN
601  DO k=1,kte-1
602    grid%t_base(k) = grid%t_1(1,k,1)
603    grid%qv_base(k) = moist(1,k,1,P_QV)
604    grid%u_base(k) = grid%u_1(1,k,1)
605    grid%v_base(k) = grid%v_1(1,k,1)
606    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
607  ENDDO
608  ENDIF
609  CALL wrf_dm_bcast_real( grid%t_base , kte )
610  CALL wrf_dm_bcast_real( grid%qv_base , kte )
611  CALL wrf_dm_bcast_real( grid%u_base , kte )
612  CALL wrf_dm_bcast_real( grid%v_base , kte )
613  CALL wrf_dm_bcast_real( grid%z_base , kte )
614
615  DO J = jts, min(jde-1,jte)
616  DO I = its, min(ide-1,ite)
617     thtmp   = grid%t_2(i,1,j)+t0
618     ptmp    = grid%p(i,1,j)+grid%pb(i,1,j)
619     temp(1) = thtmp * (ptmp/p1000mb)**rcp
620     thtmp   = grid%t_2(i,2,j)+t0
621     ptmp    = grid%p(i,2,j)+grid%pb(i,2,j)
622     temp(2) = thtmp * (ptmp/p1000mb)**rcp
623     thtmp   = grid%t_2(i,3,j)+t0
624     ptmp    = grid%p(i,3,j)+grid%pb(i,3,j)
625     temp(3) = thtmp * (ptmp/p1000mb)**rcp
626
627!    For LES-CBL, add 5 degrees to the surface temperature!
628!
629!    grid%tsk(I,J)=grid%cf1*temp(1)+grid%cf2*temp(2)+grid%cf3*temp(3)
630     grid%tsk(I,J)=grid%cf1*temp(1)+grid%cf2*temp(2)+grid%cf3*temp(3)+5.
631     grid%tmn(I,J)=grid%tsk(I,J)-0.5
632  ENDDO
633  ENDDO
634
635 END SUBROUTINE init_domain_rk
636
637   SUBROUTINE init_module_initialize
638   END SUBROUTINE init_module_initialize
639
640!---------------------------------------------------------------------
641
642!  test driver for get_sounding
643!
644!      implicit none
645!      integer n
646!      parameter(n = 1000)
647!      real zk(n),p(n),theta(n),rho(n),u(n),v(n),qv(n),pd(n)
648!      logical dry
649!      integer nl,k
650!
651!      dry = .false.
652!      dry = .true.
653!      call get_sounding( zk, p, pd, theta, rho, u, v, qv, dry, n, nl )
654!      write(6,*) ' input levels ',nl
655!      write(6,*) ' sounding '
656!      write(6,*) '  k  height(m)  press (Pa) pd(Pa) theta (K) den(kg/m^3)  u(m/s)     v(m/s)    qv(g/g) '
657!      do k=1,nl
658!        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)
659!      enddo
660!      end
661!
662!---------------------------------------------------------------------------
663
664      subroutine get_sounding( zk, p, p_dry, theta, rho, &
665                               u, v, qv, dry, nl_max, nl_in )
666      implicit none
667
668      integer nl_max, nl_in
669      real zk(nl_max), p(nl_max), theta(nl_max), rho(nl_max), &
670           u(nl_max), v(nl_max), qv(nl_max), p_dry(nl_max)
671      logical dry
672
673      integer n
674      parameter(n=1000)
675      logical debug
676      parameter( debug = .true.)
677
678! input sounding data
679
680      real p_surf, th_surf, qv_surf
681      real pi_surf, pi(n)
682      real h_input(n), th_input(n), qv_input(n), u_input(n), v_input(n)
683
684! diagnostics
685
686      real rho_surf, p_input(n), rho_input(n)
687      real pm_input(n)  !  this are for full moist sounding
688
689! local data
690
691      real p1000mb,cv,cp,r,cvpm,g
692      parameter (p1000mb = 1.e+05, r = 287, cp = 1003., cv = cp-r, cvpm = -cv/cp, g=9.81 )
693      integer k, it, nl
694      real qvf, qvf1, dz
695
696!  first, read the sounding
697
698      call read_sounding( p_surf, th_surf, qv_surf, &
699                          h_input, th_input, qv_input, u_input, v_input,n, nl, debug )
700
701      if(dry) then
702       do k=1,nl
703         qv_input(k) = 0.
704       enddo
705      endif
706
707      if(debug) write(6,*) ' number of input levels = ',nl
708
709        nl_in = nl
710        if(nl_in .gt. nl_max ) then
711          write(6,*) ' too many levels for input arrays ',nl_in,nl_max
712          call wrf_error_fatal ( ' too many levels for input arrays ' )
713        end if
714
715!  compute diagnostics,
716!  first, convert qv(g/kg) to qv(g/g)
717
718      do k=1,nl
719        qv_input(k) = 0.001*qv_input(k)
720      enddo
721
722      p_surf = 100.*p_surf  ! convert to pascals
723      qvf = 1. + rvovrd*qv_input(1)
724      rho_surf = 1./((r/p1000mb)*th_surf*qvf*((p_surf/p1000mb)**cvpm))
725      pi_surf = (p_surf/p1000mb)**(r/cp)
726
727      if(debug) then
728        write(6,*) ' surface density is ',rho_surf
729        write(6,*) ' surface pi is      ',pi_surf
730      end if
731
732
733!  integrate moist sounding hydrostatically, starting from the
734!  specified surface pressure
735!  -> first, integrate from surface to lowest level
736
737          qvf = 1. + rvovrd*qv_input(1)
738          qvf1 = 1. + qv_input(1)
739          rho_input(1) = rho_surf
740          dz = h_input(1)
741          do it=1,10
742            pm_input(1) = p_surf &
743                    - 0.5*dz*(rho_surf+rho_input(1))*g*qvf1
744            rho_input(1) = 1./((r/p1000mb)*th_input(1)*qvf*((pm_input(1)/p1000mb)**cvpm))
745          enddo
746
747! integrate up the column
748
749          do k=2,nl
750            rho_input(k) = rho_input(k-1)
751            dz = h_input(k)-h_input(k-1)
752            qvf1 = 0.5*(2.+(qv_input(k-1)+qv_input(k)))
753            qvf = 1. + rvovrd*qv_input(k)   ! qv is in g/kg here
754 
755            do it=1,10
756              pm_input(k) = pm_input(k-1) &
757                      - 0.5*dz*(rho_input(k)+rho_input(k-1))*g*qvf1
758              rho_input(k) = 1./((r/p1000mb)*th_input(k)*qvf*((pm_input(k)/p1000mb)**cvpm))
759            enddo
760          enddo
761
762!  we have the moist sounding
763
764!  next, compute the dry sounding using p at the highest level from the
765!  moist sounding and integrating down.
766
767        p_input(nl) = pm_input(nl)
768
769          do k=nl-1,1,-1
770            dz = h_input(k+1)-h_input(k)
771            p_input(k) = p_input(k+1) + 0.5*dz*(rho_input(k)+rho_input(k+1))*g
772          enddo
773
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          v(k) = v_input(k)
784          qv(k) = qv_input(k)
785
786        enddo
787
788     if(debug) then
789      write(6,*) ' sounding '
790      write(6,*) '  k  height(m)  press (Pa) pd(Pa) theta (K) den(kg/m^3)  u(m/s)     v(m/s)    qv(g/g) '
791      do k=1,nl
792        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)
793      enddo
794
795     end if
796
797      end subroutine get_sounding
798
799!-------------------------------------------------------
800
801      subroutine read_sounding( ps,ts,qvs,h,th,qv,u,v,n,nl,debug )
802      implicit none
803      integer n,nl
804      real ps,ts,qvs,h(n),th(n),qv(n),u(n),v(n)
805      logical end_of_file
806      logical debug
807
808      integer k
809
810      open(unit=10,file='input_sounding',form='formatted',status='old')
811      rewind(10)
812      read(10,*) ps, ts, qvs
813      if(debug) then
814        write(6,*) ' input sounding surface parameters '
815        write(6,*) ' surface pressure (mb) ',ps
816        write(6,*) ' surface pot. temp (K) ',ts
817        write(6,*) ' surface mixing ratio (g/kg) ',qvs
818      end if
819
820      end_of_file = .false.
821      k = 0
822
823      do while (.not. end_of_file)
824
825        read(10,*,end=100) h(k+1), th(k+1), qv(k+1), u(k+1), v(k+1)
826        k = k+1
827        if(debug) write(6,'(1x,i3,5(1x,e10.3))') k, h(k), th(k), qv(k), u(k), v(k)
828        go to 110
829 100    end_of_file = .true.
830 110    continue
831      enddo
832
833      nl = k
834
835      close(unit=10,status = 'keep')
836
837      end subroutine read_sounding
838
839END MODULE module_initialize_ideal
Note: See TracBrowser for help on using the repository browser.