source: trunk/mesoscale/LMD_LES_MARS/modif_mars/module_initialize_les.F @ 58

Last change on this file since 58 was 18, checked in by aslmd, 14 years ago

spiga:correction, v17 majeure (ajout LES), celle-ci v18 mineure

File size: 27.2 KB
Line 
1!IDEAL:MODEL_LAYER:INITIALIZATION
2!
3
4!  This MODULE holds the routines which are used to perform various initializations
5!  for the individual domains. 
6
7!  This MODULE CONTAINS the following routines:
8
9!  initialize_field_test - 1. Set different fields to different constant
10!                             values.  This is only a test.  If the correct
11!                             domain is not found (based upon the "id")
12!                             then a fatal error is issued.               
13
14!-----------------------------------------------------------------------
15
16MODULE module_initialize_ideal
17
18   USE module_domain
19   USE module_io_domain
20   USE module_state_description
21   USE module_model_constants
22   USE module_bc
23   USE module_timing
24   USE module_configure
25   USE module_init_utilities
26#ifdef DM_PARALLEL
27   USE module_dm
28#endif
29
30
31CONTAINS
32
33
34!-------------------------------------------------------------------
35! this is a wrapper for the solver-specific init_domain routines.
36! Also dereferences the grid variables and passes them down as arguments.
37! This is crucial, since the lower level routines may do message passing
38! and this will get fouled up on machines that insist on passing down
39! copies of assumed-shape arrays (by passing down as arguments, the
40! data are treated as assumed-size -- ie. f77 -- arrays and the copying
41! business is avoided).  Fie on the F90 designers.  Fie and a pox.
42
43   SUBROUTINE init_domain ( grid )
44
45   IMPLICIT NONE
46
47   !  Input data.
48   TYPE (domain), POINTER :: grid
49   !  Local data.
50   INTEGER :: idum1, idum2
51
52   CALL set_scalar_indices_from_config ( head_grid%id , idum1, idum2 )
53
54     CALL init_domain_rk( grid &
55!
56#include <actual_new_args.inc>
57!
58                        )
59
60   END SUBROUTINE init_domain
61
62!-------------------------------------------------------------------
63
64   SUBROUTINE init_domain_rk ( grid &
65!
66# include <dummy_new_args.inc>
67!
68)
69   IMPLICIT NONE
70
71   !  Input data.
72   TYPE (domain), POINTER :: grid
73
74# include <dummy_new_decl.inc>
75
76   TYPE (grid_config_rec_type)              :: config_flags
77
78   !  Local data
79   INTEGER                             ::                       &
80                                  ids, ide, jds, jde, kds, kde, &
81                                  ims, ime, jms, jme, kms, kme, &
82                                  its, ite, jts, jte, kts, kte, &
83                                  i, j, k
84
85   ! Local data
86
87   INTEGER, PARAMETER :: nl_max = 1000
88   REAL, DIMENSION(nl_max) :: zk, p_in, theta, rho, u, v, qv, pd_in
89   INTEGER :: nl_in
90
91
92   INTEGER :: icm,jcm, ii, im1, jj, jm1, loop, error, fid, nxc, nyc
93   REAL    :: u_mean,v_mean, f0, p_surf, p_level, qvf, z_at_v, z_at_u
94   REAL    :: z_scale, xrad, yrad, zrad, rad, delt, cof1, cof2
95!   REAL, EXTERNAL :: interp_0
96   REAL    :: hm
97   REAL    :: pi
98
99!  stuff from original initialization that has been dropped from the Registry
100   REAL    :: vnu, xnu, xnus, dinit0, cbh, p0_temp, t0_temp, zd, zt
101   REAL    :: qvf1, qvf2, pd_surf
102   INTEGER :: it
103   real :: thtmp, ptmp, temp(3)
104
105   LOGICAL :: moisture_init
106   LOGICAL :: stretch_grid, dry_sounding
107
108  INTEGER :: xs , xe , ys , ye
109  REAL :: mtn_ht
110   LOGICAL, EXTERNAL :: wrf_dm_on_monitor
111!  For LES, add randx
112   real :: randx
113
114!!MARS
115 REAL :: lon_input, lat_input, alt_input, tsurf_input
116!!MARS
117
118#ifdef DM_PARALLEL
119#    include <data_calls.inc>
120#endif
121
122
123   SELECT CASE ( model_data_order )
124         CASE ( DATA_ORDER_ZXY )
125   kds = grid%sd31 ; kde = grid%ed31 ;
126   ids = grid%sd32 ; ide = grid%ed32 ;
127   jds = grid%sd33 ; jde = grid%ed33 ;
128
129   kms = grid%sm31 ; kme = grid%em31 ;
130   ims = grid%sm32 ; ime = grid%em32 ;
131   jms = grid%sm33 ; jme = grid%em33 ;
132
133   kts = grid%sp31 ; kte = grid%ep31 ;   ! note that tile is entire patch
134   its = grid%sp32 ; ite = grid%ep32 ;   ! note that tile is entire patch
135   jts = grid%sp33 ; jte = grid%ep33 ;   ! note that tile is entire patch
136         CASE ( DATA_ORDER_XYZ )
137   ids = grid%sd31 ; ide = grid%ed31 ;
138   jds = grid%sd32 ; jde = grid%ed32 ;
139   kds = grid%sd33 ; kde = grid%ed33 ;
140
141   ims = grid%sm31 ; ime = grid%em31 ;
142   jms = grid%sm32 ; jme = grid%em32 ;
143   kms = grid%sm33 ; kme = grid%em33 ;
144
145   its = grid%sp31 ; ite = grid%ep31 ;   ! note that tile is entire patch
146   jts = grid%sp32 ; jte = grid%ep32 ;   ! note that tile is entire patch
147   kts = grid%sp33 ; kte = grid%ep33 ;   ! note that tile is entire patch
148         CASE ( DATA_ORDER_XZY )
149   ids = grid%sd31 ; ide = grid%ed31 ;
150   kds = grid%sd32 ; kde = grid%ed32 ;
151   jds = grid%sd33 ; jde = grid%ed33 ;
152
153   ims = grid%sm31 ; ime = grid%em31 ;
154   kms = grid%sm32 ; kme = grid%em32 ;
155   jms = grid%sm33 ; jme = grid%em33 ;
156
157   its = grid%sp31 ; ite = grid%ep31 ;   ! note that tile is entire patch
158   kts = grid%sp32 ; kte = grid%ep32 ;   ! note that tile is entire patch
159   jts = grid%sp33 ; jte = grid%ep33 ;   ! note that tile is entire patch
160
161   END SELECT
162
163
164!  stretch_grid = .true.
165!  FOR LES, set stretch to false
166   stretch_grid = .false.
167   delt = 3.
168!   z_scale = .50
169   z_scale = .40
170   pi = 2.*asin(1.0)
171   write(6,*) ' pi is ',pi
172   nxc = (ide-ids)/2
173   nyc = (jde-jds)/2
174
175   CALL model_to_grid_config_rec ( grid%id , model_config_rec , config_flags )
176
177! here we check to see if the boundary conditions are set properly
178
179   CALL boundary_condition_check( config_flags, bdyzone, error, grid%id )
180
181   moisture_init = .true.
182
183    grid%itimestep=0
184
185#ifdef DM_PARALLEL
186   CALL wrf_dm_bcast_bytes( icm , IWORDSIZE )
187   CALL wrf_dm_bcast_bytes( jcm , IWORDSIZE )
188#endif
189
190    CALL nl_set_mminlu(1, '    ')
191    CALL nl_set_iswater(1,0)
192    CALL nl_set_cen_lat(1,40.)
193    CALL nl_set_cen_lon(1,-105.)
194    CALL nl_set_truelat1(1,0.)
195    CALL nl_set_truelat2(1,0.)
196    CALL nl_set_moad_cen_lat (1,0.)
197    CALL nl_set_stand_lon (1,0.)
198    CALL nl_set_map_proj(1,0)
199
200
201!  here we initialize data we currently is not initialized
202!  in the input data
203
204    DO j = jts, jte
205      DO i = its, ite
206         grid%msftx(i,j)    = 1.
207         grid%msfty(i,j)    = 1.
208         grid%msfux(i,j)    = 1.
209         grid%msfuy(i,j)    = 1.
210         grid%msfvx(i,j)    = 1.
211         grid%msfvx_inv(i,j)= 1.
212         grid%msfvy(i,j)    = 1.
213         grid%sina(i,j)     = 0.
214         grid%cosa(i,j)     = 1.
215         grid%e(i,j)        = 0.
216!  for LES, include Coriolis force
217         grid%f(i,j)        = 0.  !!MARS MARS 1.e-4
218!!      grid%f(i,j)     = 2*EOMEG*SIN(grid%xlat(i,j)*degrad)
219      END DO
220   END DO
221
222    DO j = jts, jte
223    DO k = kts, kte
224      DO i = its, ite
225         grid%ww(i,k,j)     = 0.
226      END DO
227   END DO
228   END DO
229
230   grid%step_number = 0
231
232! set up the grid
233
234   IF (stretch_grid) THEN ! exponential stretch for eta (nearly constant dz)
235     DO k=1, kde
236      grid%znw(k) = (exp(-(k-1)/float(kde-1)/z_scale) - exp(-1./z_scale))/ &
237                                (1.-exp(-1./z_scale))
238     ENDDO
239   ELSE
240
241!!!MARS
242grid%znw(1)=1.000
243grid%znw(2)=0.9995 !5m
244grid%znw(3)=0.9980 !20m
245grid%znw(4)=0.9950 !55m
246DO k=5, kde
247   grid%znw(k) = grid%znw(4) * ( 1. - float(k-4)/float(kde-4) )
248ENDDO
249!!!!MARS
250!!
251!     DO k=1, kde
252!      grid%znw(k) = 1. - float(k-1)/float(kde-1)
253!     ENDDO
254
255   ENDIF
256
257   DO k=1, kde-1
258    grid%dnw(k) = grid%znw(k+1) - grid%znw(k)
259    grid%rdnw(k) = 1./grid%dnw(k)
260    grid%znu(k) = 0.5*(grid%znw(k+1)+grid%znw(k))
261   ENDDO
262   DO k=2, kde-1
263    grid%dn(k) = 0.5*(grid%dnw(k)+grid%dnw(k-1))
264    grid%rdn(k) = 1./grid%dn(k)
265    grid%fnp(k) = .5* grid%dnw(k  )/grid%dn(k)
266    grid%fnm(k) = .5* grid%dnw(k-1)/grid%dn(k)
267   ENDDO
268
269   cof1 = (2.*grid%dn(2)+grid%dn(3))/(grid%dn(2)+grid%dn(3))*grid%dnw(1)/grid%dn(2)
270   cof2 =     grid%dn(2)        /(grid%dn(2)+grid%dn(3))*grid%dnw(1)/grid%dn(3)
271   grid%cf1  = grid%fnp(2) + cof1
272   grid%cf2  = grid%fnm(2) - cof1 - cof2
273   grid%cf3  = cof2       
274
275   grid%cfn  = (.5*grid%dnw(kde-1)+grid%dn(kde-1))/grid%dn(kde-1)
276   grid%cfn1 = -.5*grid%dnw(kde-1)/grid%dn(kde-1)
277   grid%rdx = 1./config_flags%dx
278   grid%rdy = 1./config_flags%dy
279
280!  get the sounding from the ascii sounding file, first get dry sounding and
281!  calculate base state
282
283  dry_sounding = .true.
284  IF ( wrf_dm_on_monitor() ) THEN
285  write(6,*) ' getting dry sounding for base state '
286
287  CALL get_sounding( zk, p_in, pd_in, theta, rho, u, v, qv, dry_sounding, nl_max, nl_in )
288  ENDIF
289  CALL wrf_dm_bcast_real( zk , nl_max )
290  CALL wrf_dm_bcast_real( p_in , nl_max )
291  CALL wrf_dm_bcast_real( pd_in , nl_max )
292  CALL wrf_dm_bcast_real( theta , nl_max )
293  CALL wrf_dm_bcast_real( rho , nl_max )
294  CALL wrf_dm_bcast_real( u , nl_max )
295  CALL wrf_dm_bcast_real( v , nl_max )
296  CALL wrf_dm_bcast_real( qv , nl_max )
297  CALL wrf_dm_bcast_integer ( nl_in , 1 )
298
299  write(6,*) ' returned from reading sounding, nl_in is ',nl_in
300
301!!MARS
302!!MARS
303  open(unit=14,file='input_coord',form='formatted',status='old')
304  rewind(14)
305  read(14,*) lon_input
306  read(14,*) lat_input
307  close(14)
308  write(6,*) ' lon is ',lon_input
309  write(6,*) ' lat is ',lat_input
310!!MARS
311!!MARS
312
313!!MARS
314!!MARS
315  open(unit=18,file='input_more',form='formatted',status='old')
316  rewind(18)
317  read(18,*) alt_input, tsurf_input
318  close(18)
319  write(6,*) ' alt is ',alt_input
320  write(6,*) ' tsurf is ',tsurf_input
321!!MARS
322!!MARS
323
324!  find ptop for the desired ztop (ztop is input from the namelist),
325!  and find surface pressure
326
327  write(6,*) ' ztop above ground is ',config_flags%ztop
328  write(6,*) ' real ztop is ',config_flags%ztop + alt_input
329  grid%p_top = interp_0( p_in, zk, config_flags%ztop + alt_input, nl_in )
330
331  DO j=jts,jte
332  DO i=its,ite
333!!MARS
334    grid%ht(i,j) = alt_input
335    grid%tsk(i,j) = tsurf_input
336!!MARS
337    grid%xlat(i,j) = lat_input !+ float(j)*config_flags%dy/59000.
338    grid%xlong(i,j) = lon_input !+ float(i)*config_flags%dx/59000.
339    grid%mars_emiss(i,j)=0.95
340    grid%mars_cice(i,j)=0.
341    grid%slpx(i,j) = 0.
342    grid%slpy(i,j) = 0.
343   DO k=1,config_flags%num_soil_layers
344    grid%mars_tsoil(i,k,j) = 0.
345   ENDDO
346    grid%mars_gw(i,1,j) = 0.
347    grid%mars_gw(i,2,j) = 0.
348    grid%mars_gw(i,3,j) = 0.
349    grid%mars_gw(i,4,j) = 0.
350    grid%mars_gw(i,5,j) = 0.
351!!MARS
352  ENDDO
353  ENDDO
354
355  xs=ide/2 -3
356  xs=ids   -3
357  xe=xs + 6
358  ys=jde/2 -3
359  ye=ys + 6
360  mtn_ht = 500
361#ifdef MTN
362  DO j=max(ys,jds),min(ye,jde-1)
363  DO i=max(xs,ids),min(xe,ide-1)
364     grid%ht(i,j) = mtn_ht * 0.25 * &
365               ( 1. + COS ( 2*pi/(xe-xs) * ( i-xs ) + pi ) ) * &
366               ( 1. + COS ( 2*pi/(ye-ys) * ( j-ys ) + pi ) )
367  ENDDO
368  ENDDO
369#endif
370#ifdef EW_RIDGE
371  DO j=max(ys,jds),min(ye,jde-1)
372  DO i=ids,ide
373     grid%ht(i,j) = mtn_ht * 0.50 * &
374               ( 1. + COS ( 2*pi/(ye-ys) * ( j-ys ) + pi ) )
375  ENDDO
376  ENDDO
377#endif
378#ifdef NS_RIDGE
379  DO j=jds,jde
380  DO i=max(xs,ids),min(xe,ide-1)
381     grid%ht(i,j) = mtn_ht * 0.50 * &
382               ( 1. + COS ( 2*pi/(xe-xs) * ( i-xs ) + pi ) )
383  ENDDO
384  ENDDO
385#endif
386  DO j=jts,jte
387  DO i=its,ite
388    grid%phb(i,1,j) = g * grid%ht(i,j)
389    grid%ph0(i,1,j) = g * grid%ht(i,j)
390  ENDDO
391  ENDDO
392
393  DO J = jts, jte
394  DO I = its, ite
395
396    p_surf = interp_0( p_in, zk, grid%phb(i,1,j)/g, nl_in )
397    grid%mub(i,j) = p_surf-grid%p_top
398
399!  this is dry hydrostatic sounding (base state), so given grid%p (coordinate),
400!  interp theta (from interp) and compute 1/rho from eqn. of state
401
402    DO K = 1, kte-1
403      p_level = grid%znu(k)*(p_surf - grid%p_top) + grid%p_top
404      grid%pb(i,k,j) = p_level
405      grid%t_init(i,k,j) = interp_0( theta, p_in, p_level, nl_in ) - t0
406      grid%alb(i,k,j) = (r_d/p1000mb)*(grid%t_init(i,k,j)+t0)*(grid%pb(i,k,j)/p1000mb)**cvpm
407    ENDDO
408
409!  calc hydrostatic balance (alternatively we could interp the geopotential from the
410!  sounding, but this assures that the base state is in exact hydrostatic balance with
411!  respect to the model eqns.
412
413    DO k  = 2,kte
414      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)
415    ENDDO
416
417  ENDDO
418  ENDDO
419
420  IF ( wrf_dm_on_monitor() ) THEN
421    write(6,*) ' ptop is ',grid%p_top
422    write(6,*) ' base state grid%mub(1,1), p_surf is ',grid%mub(1,1),grid%mub(1,1)+grid%p_top
423  ENDIF
424
425!  calculate full state for each column - this includes moisture.
426
427  write(6,*) ' getting moist sounding for full state '
428  dry_sounding = .false.
429  CALL get_sounding( zk, p_in, pd_in, theta, rho, u, v, qv, dry_sounding, nl_max, nl_in )
430
431  DO J = jts, min(jde-1,jte)
432  DO I = its, min(ide-1,ite)
433
434!  At this point grid%p_top is already set. find the DRY mass in the column
435!  by interpolating the DRY pressure. 
436
437   pd_surf = interp_0( pd_in, zk, grid%phb(i,1,j)/g, nl_in )
438
439!  compute the perturbation mass and the full mass
440
441    grid%mu_1(i,j) = pd_surf-grid%p_top - grid%mub(i,j)
442    grid%mu_2(i,j) = grid%mu_1(i,j)
443    grid%mu0(i,j) = grid%mu_1(i,j) + grid%mub(i,j)
444
445! given the dry pressure and coordinate system, interp the potential
446! temperature and qv
447
448    do k=1,kde-1
449
450      p_level = grid%znu(k)*(pd_surf - grid%p_top) + grid%p_top
451
452      moist(i,k,j,P_QV) = interp_0( qv, pd_in, p_level, nl_in )
453      grid%t_1(i,k,j)          = interp_0( theta, pd_in, p_level, nl_in ) - t0
454      grid%t_2(i,k,j)          = grid%t_1(i,k,j)
455     
456
457    enddo
458
459!  integrate the hydrostatic equation (from the RHS of the bigstep
460!  vertical momentum equation) down from the top to get grid%p.
461!  first from the top of the model to the top pressure
462
463    k = kte-1  ! top level
464
465    qvf1 = 0.5*(moist(i,k,j,P_QV)+moist(i,k,j,P_QV))
466    qvf2 = 1./(1.+qvf1)
467    qvf1 = qvf1*qvf2
468
469!    grid%p(i,k,j) = - 0.5*grid%mu_1(i,j)/grid%rdnw(k)
470    grid%p(i,k,j) = - 0.5*(grid%mu_1(i,j)+qvf1*grid%mub(i,j))/grid%rdnw(k)/qvf2
471    qvf = 1. + rvovrd*moist(i,k,j,P_QV)
472    grid%alt(i,k,j) = (r_d/p1000mb)*(grid%t_1(i,k,j)+t0)*qvf* &
473                (((grid%p(i,k,j)+grid%pb(i,k,j))/p1000mb)**cvpm)
474    grid%al(i,k,j) = grid%alt(i,k,j) - grid%alb(i,k,j)
475
476!  down the column
477
478    do k=kte-2,1,-1
479      qvf1 = 0.5*(moist(i,k,j,P_QV)+moist(i,k+1,j,P_QV))
480      qvf2 = 1./(1.+qvf1)
481      qvf1 = qvf1*qvf2
482      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)
483      qvf = 1. + rvovrd*moist(i,k,j,P_QV)
484      grid%alt(i,k,j) = (r_d/p1000mb)*(grid%t_1(i,k,j)+t0)*qvf* &
485                  (((grid%p(i,k,j)+grid%pb(i,k,j))/p1000mb)**cvpm)
486      grid%al(i,k,j) = grid%alt(i,k,j) - grid%alb(i,k,j)
487    enddo
488
489!  this is the hydrostatic equation used in the model after the
490!  small timesteps.  In the model, grid%al (inverse density)
491!  is computed from the geopotential.
492
493
494    grid%ph_1(i,1,j) = 0.
495    DO k  = 2,kte
496      grid%ph_1(i,k,j) = grid%ph_1(i,k-1,j) - (1./grid%rdnw(k-1))*(       &
497                   (grid%mub(i,j)+grid%mu_1(i,j))*grid%al(i,k-1,j)+ &
498                    grid%mu_1(i,j)*grid%alb(i,k-1,j)  )
499                                                   
500      grid%ph_2(i,k,j) = grid%ph_1(i,k,j)
501      grid%ph0(i,k,j) = grid%ph_1(i,k,j) + grid%phb(i,k,j)
502    ENDDO
503
504    IF ( wrf_dm_on_monitor() ) THEN
505    if((i==2) .and. (j==2)) then
506     write(6,*) ' grid%ph_1 calc ',grid%ph_1(2,1,2),grid%ph_1(2,2,2),&
507                              grid%mu_1(2,2)+grid%mub(2,2),grid%mu_1(2,2), &
508                              grid%alb(2,1,2),grid%al(1,2,1),grid%rdnw(1)
509    endif
510    ENDIF
511
512  ENDDO
513  ENDDO
514
515!#if 0
516
517!  thermal perturbation to kick off convection
518
519  write(6,*) ' nxc, nyc for perturbation ',nxc,nyc
520  write(6,*) ' delt for perturbation ',delt
521
522! For LES, change the initial random perturbations
523! For 2D test, call randx outside I-loop
524! For 3D runs, call randx inside both I-J loops
525
526  DO J = jts, min(jde-1,jte)
527!   yrad = config_flags%dy*float(j-nyc)/10000.
528    yrad = 0.
529    DO I = its, min(ide-1,ite)
530!     xrad = config_flags%dx*float(i-nxc)/10000.
531      xrad = 0.
532      call random_number (randx)
533      randx = randx - 0.5
534!     DO K = 1, kte-1
535      DO K = 1, 4
536
537!  No bubbles for LES!
538!  put in preturbation theta (bubble) and recalc density.  note,
539!  the mass in the column is not changing, so when theta changes,
540!  we recompute density and geopotential
541
542!       zrad = 0.5*(grid%ph_1(i,k,j)+grid%ph_1(i,k+1,j)  &
543!                  +grid%phb(i,k,j)+grid%phb(i,k+1,j))/g
544!       zrad = (zrad-1500.)/1500.
545        zrad = 0.
546        RAD=SQRT(xrad*xrad+yrad*yrad+zrad*zrad)
547        IF(RAD <= 1.) THEN
548!          grid%t_1(i,k,j)=grid%t_1(i,k,j)+delt*COS(.5*PI*RAD)**2
549           grid%t_1(i,k,j)=grid%t_1(i,k,j)+ 0.1 *randx
550           grid%t_2(i,k,j)=grid%t_1(i,k,j)
551           qvf = 1. + rvovrd*moist(i,k,j,P_QV)
552           grid%alt(i,k,j) = (r_d/p1000mb)*(grid%t_1(i,k,j)+t0)*qvf* &
553                        (((grid%p(i,k,j)+grid%pb(i,k,j))/p1000mb)**cvpm)
554           grid%al(i,k,j) = grid%alt(i,k,j) - grid%alb(i,k,j)
555        ENDIF
556      ENDDO
557
558!  rebalance hydrostatically
559
560      DO k  = 2,kte
561        grid%ph_1(i,k,j) = grid%ph_1(i,k-1,j) - (1./grid%rdnw(k-1))*(       &
562                     (grid%mub(i,j)+grid%mu_1(i,j))*grid%al(i,k-1,j)+ &
563                      grid%mu_1(i,j)*grid%alb(i,k-1,j)  )
564                                                   
565        grid%ph_2(i,k,j) = grid%ph_1(i,k,j)
566        grid%ph0(i,k,j) = grid%ph_1(i,k,j) + grid%phb(i,k,j)
567      ENDDO
568
569    ENDDO
570  ENDDO
571
572!#endif
573
574   IF ( wrf_dm_on_monitor() ) THEN
575   write(6,*) ' grid%mu_1 from comp ', grid%mu_1(1,1)
576   write(6,*) ' full state sounding from comp, ph/g, grid%p, grid%al, grid%t_1, qv '
577   do k=1,kde-1
578     write(6,'(i3,1x,5(1x,1pe10.3))') k, (grid%ph_1(1,k,1)+grid%phb(1,k,1))/g, &
579                                      grid%p(1,k,1)+grid%pb(1,k,1), grid%alt(1,k,1), &
580                                      grid%t_1(1,k,1)+t0, moist(1,k,1,P_QV)
581   enddo
582
583   write(6,*) ' pert state sounding from comp, grid%ph_1, pp, alp, grid%t_1, qv '
584   do k=1,kde-1
585     write(6,'(i3,1x,5(1x,1pe10.3))') k, grid%ph_1(1,k,1), &
586                                      grid%p(1,k,1), grid%al(1,k,1), &
587                                      grid%t_1(1,k,1), moist(1,k,1,P_QV)
588   enddo
589   ENDIF
590
591! interp v
592
593  DO J = jts, jte
594  DO I = its, min(ide-1,ite)
595
596    IF (j == jds) THEN
597      z_at_v = grid%phb(i,1,j)/g
598    ELSE IF (j == jde) THEN
599      z_at_v = grid%phb(i,1,j-1)/g
600    ELSE
601      z_at_v = 0.5*(grid%phb(i,1,j)+grid%phb(i,1,j-1))/g
602    END IF
603    p_surf = interp_0( p_in, zk, z_at_v, nl_in )
604
605    DO K = 1, kte-1
606      p_level = grid%znu(k)*(p_surf - grid%p_top) + grid%p_top
607      grid%v_1(i,k,j) = interp_0( v, p_in, p_level, nl_in )
608      grid%v_2(i,k,j) = grid%v_1(i,k,j)
609    ENDDO
610
611  ENDDO
612  ENDDO
613
614! interp u
615
616  DO J = jts, min(jde-1,jte)
617  DO I = its, ite
618
619    IF (i == ids) THEN
620      z_at_u = grid%phb(i,1,j)/g
621    ELSE IF (i == ide) THEN
622      z_at_u = grid%phb(i-1,1,j)/g
623    ELSE
624      z_at_u = 0.5*(grid%phb(i,1,j)+grid%phb(i-1,1,j))/g
625    END IF
626
627    p_surf = interp_0( p_in, zk, z_at_u, nl_in )
628
629    DO K = 1, kte-1
630      p_level = grid%znu(k)*(p_surf - grid%p_top) + grid%p_top
631      grid%u_1(i,k,j) = interp_0( u, p_in, p_level, nl_in )
632      grid%u_2(i,k,j) = grid%u_1(i,k,j)
633    ENDDO
634
635  ENDDO
636  ENDDO
637
638!  set w
639
640  DO J = jts, min(jde-1,jte)
641  DO K = kts, kte
642  DO I = its, min(ide-1,ite)
643    grid%w_1(i,k,j) = 0.
644    grid%w_2(i,k,j) = 0.
645  ENDDO
646  ENDDO
647  ENDDO
648
649!!!MARS MARS
650IF (config_flags%init_MU .ne. 0.) THEN
651  grid%u_1 = grid%u_1*config_flags%init_MU
652  grid%u_2 = grid%u_2*config_flags%init_MU
653  print *, 'multiply zonal wind ', config_flags%init_MU
654ENDIF
655IF (config_flags%init_MV .ne. 0.) THEN
656  grid%v_1 = grid%v_1*config_flags%init_MV
657  grid%v_2 = grid%v_2*config_flags%init_MV
658  print *, 'multiply meridional wind ', config_flags%init_MV
659ENDIF
660IF (config_flags%init_U .ne. 0.) THEN
661  DO J = jts, min(jde-1,jte)
662  DO K = kts, kte-1
663  DO I = its, min(ide-1,ite)
664    grid%u_1(i,k,j) = config_flags%init_U
665    grid%u_2(i,k,j) = config_flags%init_U
666  ENDDO
667  ENDDO
668  ENDDO
669  print *, 'constant zonal wind ', config_flags%init_U
670  !!! ****** ou autre possibilité
671  !!! >   grid%u_1 = grid%u_1*0. + config_flags%init_U
672  !!! >   grid%u_2 = grid%u_2*0. + config_flags%init_U
673ENDIF
674IF (config_flags%init_V .ne. 0.) THEN
675  DO J = jts, min(jde-1,jte)
676  DO K = kts, kte-1
677  DO I = its, min(ide-1,ite)
678    grid%v_1(i,k,j) = config_flags%init_V
679    grid%v_2(i,k,j) = config_flags%init_V
680  ENDDO
681  ENDDO
682  ENDDO
683  print *, 'constant meridional wind ', config_flags%init_V
684ENDIF
685!!!MARS MARS
686
687
688!  set a few more things
689
690  DO J = jts, min(jde-1,jte)
691  DO K = kts, kte-1
692  DO I = its, min(ide-1,ite)
693    grid%h_diabatic(i,k,j) = 0.
694      !!!!! MARS NO WIND CASE
695      !grid%u_1(i,k,j) = 0.
696      !grid%u_2(i,k,j) = 0.
697      !grid%v_1(i,k,j) = 0.
698      !grid%v_2(i,k,j) = 0.
699      !!!!! MARS NO WIND CASE
700  ENDDO
701  ENDDO
702  ENDDO
703
704  IF ( wrf_dm_on_monitor() ) THEN
705  DO k=1,kte-1
706    grid%t_base(k) = grid%t_1(1,k,1)
707    grid%qv_base(k) = moist(1,k,1,P_QV)
708    grid%u_base(k) = grid%u_1(1,k,1)
709    grid%v_base(k) = grid%v_1(1,k,1)
710    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
711
712!!!!! MARS SIMPLE LES (PURE BUOYANCY)
713!!      grid%t_base(k)  = grid%t_init(its,k,jts)
714!      grid%t_base(k) = 0.
715!      grid%qv_base(k) = 0.
716!      grid%u_base(k)  = 0.
717!      grid%v_base(k)  = 0.
718!      grid%z_base(k) = 0.
719!!!!! MARS SIMPLE LES
720
721  ENDDO
722  ENDIF
723  CALL wrf_dm_bcast_real( grid%t_base , kte )
724  CALL wrf_dm_bcast_real( grid%qv_base , kte )
725  CALL wrf_dm_bcast_real( grid%u_base , kte )
726  CALL wrf_dm_bcast_real( grid%v_base , kte )
727  CALL wrf_dm_bcast_real( grid%z_base , kte )
728
729  DO J = jts, min(jde-1,jte)
730  DO I = its, min(ide-1,ite)
731     thtmp   = grid%t_2(i,1,j)+t0
732     ptmp    = grid%p(i,1,j)+grid%pb(i,1,j)
733     temp(1) = thtmp * (ptmp/p1000mb)**rcp
734     thtmp   = grid%t_2(i,2,j)+t0
735     ptmp    = grid%p(i,2,j)+grid%pb(i,2,j)
736     temp(2) = thtmp * (ptmp/p1000mb)**rcp
737     thtmp   = grid%t_2(i,3,j)+t0
738     ptmp    = grid%p(i,3,j)+grid%pb(i,3,j)
739     temp(3) = thtmp * (ptmp/p1000mb)**rcp
740
741!!    For LES-CBL, add 5 degrees to the surface temperature!
742!!
743!     grid%tsk(I,J)=grid%cf1*temp(1)+grid%cf2*temp(2)+grid%cf3*temp(3)
744!!     grid%tsk(I,J)=grid%cf1*temp(1)+grid%cf2*temp(2)+grid%cf3*temp(3)+5.
745     grid%tmn(I,J)=grid%tsk(I,J)-0.5
746
747  ENDDO
748  ENDDO
749
750 END SUBROUTINE init_domain_rk
751
752   SUBROUTINE init_module_initialize
753   END SUBROUTINE init_module_initialize
754
755!---------------------------------------------------------------------
756
757!  test driver for get_sounding
758!
759!      implicit none
760!      integer n
761!      parameter(n = 1000)
762!      real zk(n),p(n),theta(n),rho(n),u(n),v(n),qv(n),pd(n)
763!      logical dry
764!      integer nl,k
765!
766!      dry = .false.
767!      dry = .true.
768!      call get_sounding( zk, p, pd, theta, rho, u, v, qv, dry, n, nl )
769!      write(6,*) ' input levels ',nl
770!      write(6,*) ' sounding '
771!      write(6,*) '  k  height(m)  press (Pa) pd(Pa) theta (K) den(kg/m^3)  u(m/s)     v(m/s)    qv(g/g) '
772!      do k=1,nl
773!        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)
774!      enddo
775!      end
776!
777!---------------------------------------------------------------------------
778
779      subroutine get_sounding( zk, p, p_dry, theta, rho, &
780                               u, v, qv, dry, nl_max, nl_in )
781      implicit none
782
783      integer nl_max, nl_in
784      real zk(nl_max), p(nl_max), theta(nl_max), rho(nl_max), &
785           u(nl_max), v(nl_max), qv(nl_max), p_dry(nl_max)
786      logical dry
787
788      integer n
789      parameter(n=1000)
790      logical debug
791      parameter( debug = .true.)
792
793! input sounding data
794
795      real p_surf, th_surf, qv_surf
796      real pi_surf, pi(n)
797      real h_input(n), th_input(n), qv_input(n), u_input(n), v_input(n)
798
799! diagnostics
800
801      real rho_surf, p_input(n), rho_input(n)
802      real pm_input(n)  !  this are for full moist sounding
803
804! local data
805
806      real p1000mb,cv,cp,r,cvpm,g
807!      parameter (p1000mb = 1.e+05, r = 287, cp = 1003., cv = cp-r, cvpm = -cv/cp, g=9.81 )
808      parameter (p1000mb = 610., r = 192., cp = 844.6, cv = cp-r, cvpm = -cv/cp, g=3.72)
809      integer k, it, nl
810      real qvf, qvf1, dz
811
812!  first, read the sounding
813
814      call read_sounding( p_surf, th_surf, qv_surf, &
815                          h_input, th_input, qv_input, u_input, v_input,n, nl, debug )
816
817      if(dry) then
818       do k=1,nl
819         qv_input(k) = 0.
820       enddo
821      endif
822
823      if(debug) write(6,*) ' number of input levels = ',nl
824
825        nl_in = nl
826        if(nl_in .gt. nl_max ) then
827          write(6,*) ' too many levels for input arrays ',nl_in,nl_max
828          call wrf_error_fatal ( ' too many levels for input arrays ' )
829        end if
830
831!  compute diagnostics,
832!  first, convert qv(g/kg) to qv(g/g)
833
834      do k=1,nl
835        qv_input(k) = 0.001*qv_input(k)
836      enddo
837
838      p_surf = 100.*p_surf  ! convert to pascals
839      qvf = 1. + rvovrd*qv_input(1)
840      rho_surf = 1./((r/p1000mb)*th_surf*qvf*((p_surf/p1000mb)**cvpm))
841      pi_surf = (p_surf/p1000mb)**(r/cp)
842
843      if(debug) then
844        write(6,*) ' surface density is ',rho_surf
845        write(6,*) ' surface pi is      ',pi_surf
846      end if
847
848
849!  integrate moist sounding hydrostatically, starting from the
850!  specified surface pressure
851!  -> first, integrate from surface to lowest level
852
853          qvf = 1. + rvovrd*qv_input(1)
854          qvf1 = 1. + qv_input(1)
855          rho_input(1) = rho_surf
856          dz = h_input(1)
857          do it=1,10
858!            pm_input(1) = p_surf &
859!                    - 0.5*dz*(rho_surf+rho_input(1))*g*qvf1
860!!!MARS MARS MARS
861            pm_input(1) = p_surf
862            rho_input(1) = 1./((r/p1000mb)*th_input(1)*qvf*((pm_input(1)/p1000mb)**cvpm))
863          enddo
864
865! integrate up the column
866
867          do k=2,nl
868            rho_input(k) = rho_input(k-1)
869            dz = h_input(k)-h_input(k-1)
870            qvf1 = 0.5*(2.+(qv_input(k-1)+qv_input(k)))
871            qvf = 1. + rvovrd*qv_input(k)   ! qv is in g/kg here
872 
873            do it=1,10
874              pm_input(k) = pm_input(k-1) &
875                      - 0.5*dz*(rho_input(k)+rho_input(k-1))*g*qvf1
876              rho_input(k) = 1./((r/p1000mb)*th_input(k)*qvf*((pm_input(k)/p1000mb)**cvpm))
877            enddo
878          enddo
879
880!  we have the moist sounding
881
882!  next, compute the dry sounding using p at the highest level from the
883!  moist sounding and integrating down.
884
885        p_input(nl) = pm_input(nl)
886
887          do k=nl-1,1,-1
888            dz = h_input(k+1)-h_input(k)
889            p_input(k) = p_input(k+1) + 0.5*dz*(rho_input(k)+rho_input(k+1))*g
890          enddo
891
892
893        do k=1,nl
894
895          zk(k) = h_input(k)
896          p(k) = pm_input(k)
897          p_dry(k) = p_input(k)
898          theta(k) = th_input(k)
899          rho(k) = rho_input(k)
900          u(k) = u_input(k)
901          v(k) = v_input(k)
902          qv(k) = qv_input(k)
903
904        enddo
905
906     if(debug) then
907      write(6,*) ' sounding '
908      write(6,*) '  k  height(m)  press (Pa) pd(Pa) theta (K) den(kg/m^3)  u(m/s)     v(m/s)    qv(g/g) '
909      do k=1,nl
910        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)
911      enddo
912
913     end if
914
915      end subroutine get_sounding
916
917!-------------------------------------------------------
918
919      subroutine read_sounding( ps,ts,qvs,h,th,qv,u,v,n,nl,debug )
920      implicit none
921      integer n,nl
922      real ps,ts,qvs,h(n),th(n),qv(n),u(n),v(n)
923      logical end_of_file
924      logical debug
925
926      integer k
927
928      open(unit=10,file='input_sounding',form='formatted',status='old')
929      rewind(10)
930      read(10,*) ps, ts, qvs
931      if(debug) then
932        write(6,*) ' input sounding surface parameters '
933        write(6,*) ' surface pressure (mb) ',ps
934        write(6,*) ' surface pot. temp (K) ',ts
935        write(6,*) ' surface mixing ratio (g/kg) ',qvs
936      end if
937
938      end_of_file = .false.
939      k = 0
940
941      do while (.not. end_of_file)
942
943        read(10,*,end=100) h(k+1), th(k+1), qv(k+1), u(k+1), v(k+1)
944        k = k+1
945        if(debug) write(6,'(1x,i3,5(1x,e10.3))') k, h(k), th(k), qv(k), u(k), v(k)
946        go to 110
947 100    end_of_file = .true.
948 110    continue
949      enddo
950
951      nl = k
952
953      close(unit=10,status = 'keep')
954
955      end subroutine read_sounding
956
957END MODULE module_initialize_ideal
Note: See TracBrowser for help on using the repository browser.