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