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