source: trunk/WRF.COMMON/WRFV2/dyn_em/module_initialize_squall2d_y.F @ 3026

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

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

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