source: trunk/WRF.COMMON/WRFV2/dyn_em/module_initialize_hill2d_x.F @ 3555

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

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

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