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