source: trunk/mesoscale/LMD_LES_MARS/modif_mars/module_initialize_les.F.hill @ 126

Last change on this file since 126 was 115, checked in by aslmd, 14 years ago

LMD_LES_MARS: LES avec topo. Registry OK avec LES simple.

File size: 28.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, xa
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             icm = ide/2
332             jcm = jde/2
333             !!3D hill
334             !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
335             !!MARS : mountain
336             !!MARS : mountain ex. hm = 2000. xa = 6.0
337               open(unit=22,file='ze_hill',form='formatted',status='old')
338               rewind(22)
339               read(22,*) hm, xa
340               write(6,*) 'height, width ', hm, xa
341               close(22)
342             !!MARS
343             !!MARS
344             !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
345
346  DO j=jts,jte
347  DO i=its,ite
348!!MARS
349    grid%ht(i,j) = alt_input
350
351             !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
352             grid%ht(i,j) = alt_input + hm/(1.+(float(i-icm)/xa)**2+(float(j-jcm)/xa)**2)
353             !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
354
355    grid%tsk(i,j) = tsurf_input
356!!MARS
357    grid%xlat(i,j) = lat_input !+ float(j)*config_flags%dy/59000.
358    grid%xlong(i,j) = lon_input !+ float(i)*config_flags%dx/59000.
359    grid%mars_emiss(i,j)=0.95
360    grid%mars_cice(i,j)=0.
361    grid%mars_wice(i,j)=0.
362    grid%slpx(i,j) = 0.
363    grid%slpy(i,j) = 0.
364   DO k=1,config_flags%num_soil_layers
365    grid%mars_tsoil(i,k,j) = 0.
366   ENDDO
367    grid%mars_gw(i,1,j) = 0.
368    grid%mars_gw(i,2,j) = 0.
369    grid%mars_gw(i,3,j) = 0.
370    grid%mars_gw(i,4,j) = 0.
371    grid%mars_gw(i,5,j) = 0.
372!!MARS
373  ENDDO
374  ENDDO
375
376  xs=ide/2 -3
377  xs=ids   -3
378  xe=xs + 6
379  ys=jde/2 -3
380  ye=ys + 6
381  mtn_ht = 500
382#ifdef MTN
383  DO j=max(ys,jds),min(ye,jde-1)
384  DO i=max(xs,ids),min(xe,ide-1)
385     grid%ht(i,j) = mtn_ht * 0.25 * &
386               ( 1. + COS ( 2*pi/(xe-xs) * ( i-xs ) + pi ) ) * &
387               ( 1. + COS ( 2*pi/(ye-ys) * ( j-ys ) + pi ) )
388  ENDDO
389  ENDDO
390#endif
391#ifdef EW_RIDGE
392  DO j=max(ys,jds),min(ye,jde-1)
393  DO i=ids,ide
394     grid%ht(i,j) = mtn_ht * 0.50 * &
395               ( 1. + COS ( 2*pi/(ye-ys) * ( j-ys ) + pi ) )
396  ENDDO
397  ENDDO
398#endif
399#ifdef NS_RIDGE
400  DO j=jds,jde
401  DO i=max(xs,ids),min(xe,ide-1)
402     grid%ht(i,j) = mtn_ht * 0.50 * &
403               ( 1. + COS ( 2*pi/(xe-xs) * ( i-xs ) + pi ) )
404  ENDDO
405  ENDDO
406#endif
407  DO j=jts,jte
408  DO i=its,ite
409    grid%phb(i,1,j) = g * grid%ht(i,j)
410    grid%ph0(i,1,j) = g * grid%ht(i,j)
411  ENDDO
412  ENDDO
413
414  DO J = jts, jte
415  DO I = its, ite
416
417    p_surf = interp_0( p_in, zk, grid%phb(i,1,j)/g, nl_in )
418    grid%mub(i,j) = p_surf-grid%p_top
419
420!  this is dry hydrostatic sounding (base state), so given grid%p (coordinate),
421!  interp theta (from interp) and compute 1/rho from eqn. of state
422
423    DO K = 1, kte-1
424      p_level = grid%znu(k)*(p_surf - grid%p_top) + grid%p_top
425      grid%pb(i,k,j) = p_level
426      grid%t_init(i,k,j) = interp_0( theta, p_in, p_level, nl_in ) - t0
427      grid%alb(i,k,j) = (r_d/p1000mb)*(grid%t_init(i,k,j)+t0)*(grid%pb(i,k,j)/p1000mb)**cvpm
428    ENDDO
429
430!  calc hydrostatic balance (alternatively we could interp the geopotential from the
431!  sounding, but this assures that the base state is in exact hydrostatic balance with
432!  respect to the model eqns.
433
434    DO k  = 2,kte
435      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)
436    ENDDO
437
438  ENDDO
439  ENDDO
440
441  IF ( wrf_dm_on_monitor() ) THEN
442    write(6,*) ' ptop is ',grid%p_top
443    write(6,*) ' base state grid%mub(1,1), p_surf is ',grid%mub(1,1),grid%mub(1,1)+grid%p_top
444  ENDIF
445
446!  calculate full state for each column - this includes moisture.
447
448  write(6,*) ' getting moist sounding for full state '
449  dry_sounding = .false.
450  CALL get_sounding( zk, p_in, pd_in, theta, rho, u, v, qv, dry_sounding, nl_max, nl_in )
451
452  DO J = jts, min(jde-1,jte)
453  DO I = its, min(ide-1,ite)
454
455!  At this point grid%p_top is already set. find the DRY mass in the column
456!  by interpolating the DRY pressure. 
457
458   pd_surf = interp_0( pd_in, zk, grid%phb(i,1,j)/g, nl_in )
459
460!  compute the perturbation mass and the full mass
461
462    grid%mu_1(i,j) = pd_surf-grid%p_top - grid%mub(i,j)
463    grid%mu_2(i,j) = grid%mu_1(i,j)
464    grid%mu0(i,j) = grid%mu_1(i,j) + grid%mub(i,j)
465
466! given the dry pressure and coordinate system, interp the potential
467! temperature and qv
468
469    do k=1,kde-1
470
471      p_level = grid%znu(k)*(pd_surf - grid%p_top) + grid%p_top
472
473      moist(i,k,j,P_QV) = interp_0( qv, pd_in, p_level, nl_in )
474      grid%t_1(i,k,j)          = interp_0( theta, pd_in, p_level, nl_in ) - t0
475      grid%t_2(i,k,j)          = grid%t_1(i,k,j)
476     
477
478    enddo
479
480!  integrate the hydrostatic equation (from the RHS of the bigstep
481!  vertical momentum equation) down from the top to get grid%p.
482!  first from the top of the model to the top pressure
483
484    k = kte-1  ! top level
485
486    qvf1 = 0.5*(moist(i,k,j,P_QV)+moist(i,k,j,P_QV))
487    qvf2 = 1./(1.+qvf1)
488    qvf1 = qvf1*qvf2
489
490!    grid%p(i,k,j) = - 0.5*grid%mu_1(i,j)/grid%rdnw(k)
491    grid%p(i,k,j) = - 0.5*(grid%mu_1(i,j)+qvf1*grid%mub(i,j))/grid%rdnw(k)/qvf2
492    qvf = 1. + rvovrd*moist(i,k,j,P_QV)
493    grid%alt(i,k,j) = (r_d/p1000mb)*(grid%t_1(i,k,j)+t0)*qvf* &
494                (((grid%p(i,k,j)+grid%pb(i,k,j))/p1000mb)**cvpm)
495    grid%al(i,k,j) = grid%alt(i,k,j) - grid%alb(i,k,j)
496
497!  down the column
498
499    do k=kte-2,1,-1
500      qvf1 = 0.5*(moist(i,k,j,P_QV)+moist(i,k+1,j,P_QV))
501      qvf2 = 1./(1.+qvf1)
502      qvf1 = qvf1*qvf2
503      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)
504      qvf = 1. + rvovrd*moist(i,k,j,P_QV)
505      grid%alt(i,k,j) = (r_d/p1000mb)*(grid%t_1(i,k,j)+t0)*qvf* &
506                  (((grid%p(i,k,j)+grid%pb(i,k,j))/p1000mb)**cvpm)
507      grid%al(i,k,j) = grid%alt(i,k,j) - grid%alb(i,k,j)
508    enddo
509
510!  this is the hydrostatic equation used in the model after the
511!  small timesteps.  In the model, grid%al (inverse density)
512!  is computed from the geopotential.
513
514
515    grid%ph_1(i,1,j) = 0.
516    DO k  = 2,kte
517      grid%ph_1(i,k,j) = grid%ph_1(i,k-1,j) - (1./grid%rdnw(k-1))*(       &
518                   (grid%mub(i,j)+grid%mu_1(i,j))*grid%al(i,k-1,j)+ &
519                    grid%mu_1(i,j)*grid%alb(i,k-1,j)  )
520                                                   
521      grid%ph_2(i,k,j) = grid%ph_1(i,k,j)
522      grid%ph0(i,k,j) = grid%ph_1(i,k,j) + grid%phb(i,k,j)
523    ENDDO
524
525    IF ( wrf_dm_on_monitor() ) THEN
526    if((i==2) .and. (j==2)) then
527     write(6,*) ' grid%ph_1 calc ',grid%ph_1(2,1,2),grid%ph_1(2,2,2),&
528                              grid%mu_1(2,2)+grid%mub(2,2),grid%mu_1(2,2), &
529                              grid%alb(2,1,2),grid%al(1,2,1),grid%rdnw(1)
530    endif
531    ENDIF
532
533  ENDDO
534  ENDDO
535
536!#if 0
537
538!  thermal perturbation to kick off convection
539
540  write(6,*) ' nxc, nyc for perturbation ',nxc,nyc
541  write(6,*) ' delt for perturbation ',delt
542
543! For LES, change the initial random perturbations
544! For 2D test, call randx outside I-loop
545! For 3D runs, call randx inside both I-J loops
546
547  DO J = jts, min(jde-1,jte)
548!   yrad = config_flags%dy*float(j-nyc)/10000.
549    yrad = 0.
550    DO I = its, min(ide-1,ite)
551!     xrad = config_flags%dx*float(i-nxc)/10000.
552      xrad = 0.
553      call random_number (randx)
554      randx = randx - 0.5
555!     DO K = 1, kte-1
556      DO K = 1, 4
557
558!  No bubbles for LES!
559!  put in preturbation theta (bubble) and recalc density.  note,
560!  the mass in the column is not changing, so when theta changes,
561!  we recompute density and geopotential
562
563!       zrad = 0.5*(grid%ph_1(i,k,j)+grid%ph_1(i,k+1,j)  &
564!                  +grid%phb(i,k,j)+grid%phb(i,k+1,j))/g
565!       zrad = (zrad-1500.)/1500.
566        zrad = 0.
567        RAD=SQRT(xrad*xrad+yrad*yrad+zrad*zrad)
568        IF(RAD <= 1.) THEN
569!          grid%t_1(i,k,j)=grid%t_1(i,k,j)+delt*COS(.5*PI*RAD)**2
570           grid%t_1(i,k,j)=grid%t_1(i,k,j)+ 0.1 *randx
571           grid%t_2(i,k,j)=grid%t_1(i,k,j)
572           qvf = 1. + rvovrd*moist(i,k,j,P_QV)
573           grid%alt(i,k,j) = (r_d/p1000mb)*(grid%t_1(i,k,j)+t0)*qvf* &
574                        (((grid%p(i,k,j)+grid%pb(i,k,j))/p1000mb)**cvpm)
575           grid%al(i,k,j) = grid%alt(i,k,j) - grid%alb(i,k,j)
576        ENDIF
577      ENDDO
578
579!  rebalance hydrostatically
580
581      DO k  = 2,kte
582        grid%ph_1(i,k,j) = grid%ph_1(i,k-1,j) - (1./grid%rdnw(k-1))*(       &
583                     (grid%mub(i,j)+grid%mu_1(i,j))*grid%al(i,k-1,j)+ &
584                      grid%mu_1(i,j)*grid%alb(i,k-1,j)  )
585                                                   
586        grid%ph_2(i,k,j) = grid%ph_1(i,k,j)
587        grid%ph0(i,k,j) = grid%ph_1(i,k,j) + grid%phb(i,k,j)
588      ENDDO
589
590    ENDDO
591  ENDDO
592
593!#endif
594
595   IF ( wrf_dm_on_monitor() ) THEN
596   write(6,*) ' grid%mu_1 from comp ', grid%mu_1(1,1)
597   write(6,*) ' full state sounding from comp, ph/g, grid%p, grid%al, grid%t_1, qv '
598   do k=1,kde-1
599     write(6,'(i3,1x,5(1x,1pe10.3))') k, (grid%ph_1(1,k,1)+grid%phb(1,k,1))/g, &
600                                      grid%p(1,k,1)+grid%pb(1,k,1), grid%alt(1,k,1), &
601                                      grid%t_1(1,k,1)+t0, moist(1,k,1,P_QV)
602   enddo
603
604   write(6,*) ' pert state sounding from comp, grid%ph_1, pp, alp, grid%t_1, qv '
605   do k=1,kde-1
606     write(6,'(i3,1x,5(1x,1pe10.3))') k, grid%ph_1(1,k,1), &
607                                      grid%p(1,k,1), grid%al(1,k,1), &
608                                      grid%t_1(1,k,1), moist(1,k,1,P_QV)
609   enddo
610   ENDIF
611
612! interp v
613
614  DO J = jts, jte
615  DO I = its, min(ide-1,ite)
616
617    IF (j == jds) THEN
618      z_at_v = grid%phb(i,1,j)/g
619    ELSE IF (j == jde) THEN
620      z_at_v = grid%phb(i,1,j-1)/g
621    ELSE
622      z_at_v = 0.5*(grid%phb(i,1,j)+grid%phb(i,1,j-1))/g
623    END IF
624    p_surf = interp_0( p_in, zk, z_at_v, nl_in )
625
626    DO K = 1, kte-1
627      p_level = grid%znu(k)*(p_surf - grid%p_top) + grid%p_top
628      grid%v_1(i,k,j) = interp_0( v, p_in, p_level, nl_in )
629      grid%v_2(i,k,j) = grid%v_1(i,k,j)
630    ENDDO
631
632  ENDDO
633  ENDDO
634
635! interp u
636
637  DO J = jts, min(jde-1,jte)
638  DO I = its, ite
639
640    IF (i == ids) THEN
641      z_at_u = grid%phb(i,1,j)/g
642    ELSE IF (i == ide) THEN
643      z_at_u = grid%phb(i-1,1,j)/g
644    ELSE
645      z_at_u = 0.5*(grid%phb(i,1,j)+grid%phb(i-1,1,j))/g
646    END IF
647
648    p_surf = interp_0( p_in, zk, z_at_u, nl_in )
649
650    DO K = 1, kte-1
651      p_level = grid%znu(k)*(p_surf - grid%p_top) + grid%p_top
652      grid%u_1(i,k,j) = interp_0( u, p_in, p_level, nl_in )
653      grid%u_2(i,k,j) = grid%u_1(i,k,j)
654    ENDDO
655
656  ENDDO
657  ENDDO
658
659!  set w
660
661  DO J = jts, min(jde-1,jte)
662  DO K = kts, kte
663  DO I = its, min(ide-1,ite)
664    grid%w_1(i,k,j) = 0.
665    grid%w_2(i,k,j) = 0.
666  ENDDO
667  ENDDO
668  ENDDO
669
670!!!MARS MARS
671IF (config_flags%init_MU .ne. 0.) THEN
672  grid%u_1 = grid%u_1*config_flags%init_MU
673  grid%u_2 = grid%u_2*config_flags%init_MU
674  print *, 'multiply zonal wind ', config_flags%init_MU
675ENDIF
676IF (config_flags%init_MV .ne. 0.) THEN
677  grid%v_1 = grid%v_1*config_flags%init_MV
678  grid%v_2 = grid%v_2*config_flags%init_MV
679  print *, 'multiply meridional wind ', config_flags%init_MV
680ENDIF
681IF (config_flags%init_U .ne. 0.) THEN
682  DO J = jts, min(jde-1,jte)
683  DO K = kts, kte-1
684  DO I = its, min(ide-1,ite)
685    grid%u_1(i,k,j) = config_flags%init_U
686    grid%u_2(i,k,j) = config_flags%init_U
687  ENDDO
688  ENDDO
689  ENDDO
690  print *, 'constant zonal wind ', config_flags%init_U
691  !!! ****** ou autre possibilité
692  !!! >   grid%u_1 = grid%u_1*0. + config_flags%init_U
693  !!! >   grid%u_2 = grid%u_2*0. + config_flags%init_U
694ENDIF
695IF (config_flags%init_V .ne. 0.) THEN
696  DO J = jts, min(jde-1,jte)
697  DO K = kts, kte-1
698  DO I = its, min(ide-1,ite)
699    grid%v_1(i,k,j) = config_flags%init_V
700    grid%v_2(i,k,j) = config_flags%init_V
701  ENDDO
702  ENDDO
703  ENDDO
704  print *, 'constant meridional wind ', config_flags%init_V
705ENDIF
706!!!MARS MARS
707
708
709!  set a few more things
710
711  DO J = jts, min(jde-1,jte)
712  DO K = kts, kte-1
713  DO I = its, min(ide-1,ite)
714    grid%h_diabatic(i,k,j) = 0.
715      !!!!! MARS NO WIND CASE
716      !grid%u_1(i,k,j) = 0.
717      !grid%u_2(i,k,j) = 0.
718      !grid%v_1(i,k,j) = 0.
719      !grid%v_2(i,k,j) = 0.
720      !!!!! MARS NO WIND CASE
721  ENDDO
722  ENDDO
723  ENDDO
724
725  IF ( wrf_dm_on_monitor() ) THEN
726  DO k=1,kte-1
727    grid%t_base(k) = grid%t_1(1,k,1)
728    grid%qv_base(k) = moist(1,k,1,P_QV)
729    grid%u_base(k) = grid%u_1(1,k,1)
730    grid%v_base(k) = grid%v_1(1,k,1)
731    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
732
733!!!!! MARS SIMPLE LES (PURE BUOYANCY)
734!!      grid%t_base(k)  = grid%t_init(its,k,jts)
735!      grid%t_base(k) = 0.
736!      grid%qv_base(k) = 0.
737!      grid%u_base(k)  = 0.
738!      grid%v_base(k)  = 0.
739!      grid%z_base(k) = 0.
740!!!!! MARS SIMPLE LES
741
742  ENDDO
743  ENDIF
744  CALL wrf_dm_bcast_real( grid%t_base , kte )
745  CALL wrf_dm_bcast_real( grid%qv_base , kte )
746  CALL wrf_dm_bcast_real( grid%u_base , kte )
747  CALL wrf_dm_bcast_real( grid%v_base , kte )
748  CALL wrf_dm_bcast_real( grid%z_base , kte )
749
750  DO J = jts, min(jde-1,jte)
751  DO I = its, min(ide-1,ite)
752     thtmp   = grid%t_2(i,1,j)+t0
753     ptmp    = grid%p(i,1,j)+grid%pb(i,1,j)
754     temp(1) = thtmp * (ptmp/p1000mb)**rcp
755     thtmp   = grid%t_2(i,2,j)+t0
756     ptmp    = grid%p(i,2,j)+grid%pb(i,2,j)
757     temp(2) = thtmp * (ptmp/p1000mb)**rcp
758     thtmp   = grid%t_2(i,3,j)+t0
759     ptmp    = grid%p(i,3,j)+grid%pb(i,3,j)
760     temp(3) = thtmp * (ptmp/p1000mb)**rcp
761
762!!    For LES-CBL, add 5 degrees to the surface temperature!
763!!
764!     grid%tsk(I,J)=grid%cf1*temp(1)+grid%cf2*temp(2)+grid%cf3*temp(3)
765!!     grid%tsk(I,J)=grid%cf1*temp(1)+grid%cf2*temp(2)+grid%cf3*temp(3)+5.
766     grid%tmn(I,J)=grid%tsk(I,J)-0.5
767
768  ENDDO
769  ENDDO
770
771 END SUBROUTINE init_domain_rk
772
773   SUBROUTINE init_module_initialize
774   END SUBROUTINE init_module_initialize
775
776!---------------------------------------------------------------------
777
778!  test driver for get_sounding
779!
780!      implicit none
781!      integer n
782!      parameter(n = 1000)
783!      real zk(n),p(n),theta(n),rho(n),u(n),v(n),qv(n),pd(n)
784!      logical dry
785!      integer nl,k
786!
787!      dry = .false.
788!      dry = .true.
789!      call get_sounding( zk, p, pd, theta, rho, u, v, qv, dry, n, nl )
790!      write(6,*) ' input levels ',nl
791!      write(6,*) ' sounding '
792!      write(6,*) '  k  height(m)  press (Pa) pd(Pa) theta (K) den(kg/m^3)  u(m/s)     v(m/s)    qv(g/g) '
793!      do k=1,nl
794!        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)
795!      enddo
796!      end
797!
798!---------------------------------------------------------------------------
799
800      subroutine get_sounding( zk, p, p_dry, theta, rho, &
801                               u, v, qv, dry, nl_max, nl_in )
802      implicit none
803
804      integer nl_max, nl_in
805      real zk(nl_max), p(nl_max), theta(nl_max), rho(nl_max), &
806           u(nl_max), v(nl_max), qv(nl_max), p_dry(nl_max)
807      logical dry
808
809      integer n
810      parameter(n=1000)
811      logical debug
812      parameter( debug = .true.)
813
814! input sounding data
815
816      real p_surf, th_surf, qv_surf
817      real pi_surf, pi(n)
818      real h_input(n), th_input(n), qv_input(n), u_input(n), v_input(n)
819
820! diagnostics
821
822      real rho_surf, p_input(n), rho_input(n)
823      real pm_input(n)  !  this are for full moist sounding
824
825! local data
826
827      real p1000mb,cv,cp,r,cvpm,g
828!      parameter (p1000mb = 1.e+05, r = 287, cp = 1003., cv = cp-r, cvpm = -cv/cp, g=9.81 )
829!      parameter (p1000mb = 610., r = 192., cp = 844.6, cv = cp-r, cvpm = -cv/cp, g=3.72)
830      parameter (p1000mb = 610., r = 191., cp = 744.5, cv = cp-r, cvpm = -cv/cp, g=3.72)
831      integer k, it, nl
832      real qvf, qvf1, dz
833
834!  first, read the sounding
835
836      call read_sounding( p_surf, th_surf, qv_surf, &
837                          h_input, th_input, qv_input, u_input, v_input,n, nl, debug )
838
839      if(dry) then
840       do k=1,nl
841         qv_input(k) = 0.
842       enddo
843      endif
844
845      if(debug) write(6,*) ' number of input levels = ',nl
846
847        nl_in = nl
848        if(nl_in .gt. nl_max ) then
849          write(6,*) ' too many levels for input arrays ',nl_in,nl_max
850          call wrf_error_fatal ( ' too many levels for input arrays ' )
851        end if
852
853!  compute diagnostics,
854!  first, convert qv(g/kg) to qv(g/g)
855
856      do k=1,nl
857        qv_input(k) = 0.001*qv_input(k)
858      enddo
859
860      p_surf = 100.*p_surf  ! convert to pascals
861      qvf = 1. + rvovrd*qv_input(1)
862      rho_surf = 1./((r/p1000mb)*th_surf*qvf*((p_surf/p1000mb)**cvpm))
863      pi_surf = (p_surf/p1000mb)**(r/cp)
864
865      if(debug) then
866        write(6,*) ' surface density is ',rho_surf
867        write(6,*) ' surface pi is      ',pi_surf
868      end if
869
870
871!  integrate moist sounding hydrostatically, starting from the
872!  specified surface pressure
873!  -> first, integrate from surface to lowest level
874
875          qvf = 1. + rvovrd*qv_input(1)
876          qvf1 = 1. + qv_input(1)
877          rho_input(1) = rho_surf
878          dz = h_input(1)
879          do it=1,10
880!            pm_input(1) = p_surf &
881!                    - 0.5*dz*(rho_surf+rho_input(1))*g*qvf1
882!!!MARS MARS MARS
883            pm_input(1) = p_surf
884            rho_input(1) = 1./((r/p1000mb)*th_input(1)*qvf*((pm_input(1)/p1000mb)**cvpm))
885          enddo
886
887! integrate up the column
888
889          do k=2,nl
890            rho_input(k) = rho_input(k-1)
891            dz = h_input(k)-h_input(k-1)
892            qvf1 = 0.5*(2.+(qv_input(k-1)+qv_input(k)))
893            qvf = 1. + rvovrd*qv_input(k)   ! qv is in g/kg here
894 
895            do it=1,10
896              pm_input(k) = pm_input(k-1) &
897                      - 0.5*dz*(rho_input(k)+rho_input(k-1))*g*qvf1
898              rho_input(k) = 1./((r/p1000mb)*th_input(k)*qvf*((pm_input(k)/p1000mb)**cvpm))
899            enddo
900          enddo
901
902!  we have the moist sounding
903
904!  next, compute the dry sounding using p at the highest level from the
905!  moist sounding and integrating down.
906
907        p_input(nl) = pm_input(nl)
908
909          do k=nl-1,1,-1
910            dz = h_input(k+1)-h_input(k)
911            p_input(k) = p_input(k+1) + 0.5*dz*(rho_input(k)+rho_input(k+1))*g
912          enddo
913
914
915        do k=1,nl
916
917          zk(k) = h_input(k)
918          p(k) = pm_input(k)
919          p_dry(k) = p_input(k)
920          theta(k) = th_input(k)
921          rho(k) = rho_input(k)
922          u(k) = u_input(k)
923          v(k) = v_input(k)
924          qv(k) = qv_input(k)
925
926        enddo
927
928     if(debug) then
929      write(6,*) ' sounding '
930      write(6,*) '  k  height(m)  press (Pa) pd(Pa) theta (K) den(kg/m^3)  u(m/s)     v(m/s)    qv(g/g) '
931      do k=1,nl
932        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)
933      enddo
934
935     end if
936
937      end subroutine get_sounding
938
939!-------------------------------------------------------
940
941      subroutine read_sounding( ps,ts,qvs,h,th,qv,u,v,n,nl,debug )
942      implicit none
943      integer n,nl
944      real ps,ts,qvs,h(n),th(n),qv(n),u(n),v(n)
945      logical end_of_file
946      logical debug
947
948      integer k
949
950      open(unit=10,file='input_sounding',form='formatted',status='old')
951      rewind(10)
952      read(10,*) ps, ts, qvs
953      if(debug) then
954        write(6,*) ' input sounding surface parameters '
955        write(6,*) ' surface pressure (mb) ',ps
956        write(6,*) ' surface pot. temp (K) ',ts
957        write(6,*) ' surface mixing ratio (g/kg) ',qvs
958      end if
959
960      end_of_file = .false.
961      k = 0
962
963      do while (.not. end_of_file)
964
965        read(10,*,end=100) h(k+1), th(k+1), qv(k+1), u(k+1), v(k+1)
966        k = k+1
967        if(debug) write(6,'(1x,i3,5(1x,e10.3))') k, h(k), th(k), qv(k), u(k), v(k)
968        go to 110
969 100    end_of_file = .true.
970 110    continue
971      enddo
972
973      nl = k
974
975      close(unit=10,status = 'keep')
976
977      end subroutine read_sounding
978
979END MODULE module_initialize_ideal
Note: See TracBrowser for help on using the repository browser.