source: trunk/MESOSCALE_DEV/WORK/co2clouds_module_initialize_quarter_ss.F @ 357

Last change on this file since 357 was 207, checked in by aslmd, 13 years ago

MESOSCALE: A GENERAL CLEAN-UP FOLLOWING UPDATING THE USER MANUAL. EVERYTHING ESSENTIAL IS IN MESOSCALE (much lighter than before). EVERYTHING FOR DEVELOPPERS OR EXPERTS IS IN MESOSCALE_DEV.

  • Property svn:executable set to *
File size: 27.7 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  INTEGER :: xs , xe , ys , ye
121  REAL :: mtn_ht
122   LOGICAL, EXTERNAL :: wrf_dm_on_monitor
123
124!!MARS
125 REAL :: lon_input, lat_input, alt_input, tsurf_input
126!!MARS
127
128
129#ifdef DM_PARALLEL
130#    include <em_data_calls.inc>
131#endif
132
133
134   SELECT CASE ( model_data_order )
135         CASE ( DATA_ORDER_ZXY )
136   kds = grid%sd31 ; kde = grid%ed31 ;
137   ids = grid%sd32 ; ide = grid%ed32 ;
138   jds = grid%sd33 ; jde = grid%ed33 ;
139
140   kms = grid%sm31 ; kme = grid%em31 ;
141   ims = grid%sm32 ; ime = grid%em32 ;
142   jms = grid%sm33 ; jme = grid%em33 ;
143
144   kts = grid%sp31 ; kte = grid%ep31 ;   ! note that tile is entire patch
145   its = grid%sp32 ; ite = grid%ep32 ;   ! note that tile is entire patch
146   jts = grid%sp33 ; jte = grid%ep33 ;   ! note that tile is entire patch
147         CASE ( DATA_ORDER_XYZ )
148   ids = grid%sd31 ; ide = grid%ed31 ;
149   jds = grid%sd32 ; jde = grid%ed32 ;
150   kds = grid%sd33 ; kde = grid%ed33 ;
151
152   ims = grid%sm31 ; ime = grid%em31 ;
153   jms = grid%sm32 ; jme = grid%em32 ;
154   kms = grid%sm33 ; kme = grid%em33 ;
155
156   its = grid%sp31 ; ite = grid%ep31 ;   ! note that tile is entire patch
157   jts = grid%sp32 ; jte = grid%ep32 ;   ! note that tile is entire patch
158   kts = grid%sp33 ; kte = grid%ep33 ;   ! note that tile is entire patch
159         CASE ( DATA_ORDER_XZY )
160   ids = grid%sd31 ; ide = grid%ed31 ;
161   kds = grid%sd32 ; kde = grid%ed32 ;
162   jds = grid%sd33 ; jde = grid%ed33 ;
163
164   ims = grid%sm31 ; ime = grid%em31 ;
165   kms = grid%sm32 ; kme = grid%em32 ;
166   jms = grid%sm33 ; jme = grid%em33 ;
167
168   its = grid%sp31 ; ite = grid%ep31 ;   ! note that tile is entire patch
169   kts = grid%sp32 ; kte = grid%ep32 ;   ! note that tile is entire patch
170   jts = grid%sp33 ; jte = grid%ep33 ;   ! note that tile is entire patch
171
172   END SELECT
173
174
175   stretch_grid = .true.
176!   stretch_grid = .false.
177
178   delt = 3.
179   delt = 10.
180!   z_scale = .50
181   z_scale = .40
182
183!!!MARS
184!!!MARS
185!  open(unit=16,file='input_vert',form='formatted',status='old')
186!  rewind(16)
187!  read(16,*) delt, z_scale
188!  write(6,*) 'delt, z_scale are  ', delt, z_scale
189!  close(16)
190!!!MARS
191!!!MARS
192
193   pi = 2.*asin(1.0)
194   write(6,*) ' pi is ',pi
195   nxc = (ide-ids)/2
196   nyc = (jde-jds)/2
197
198   CALL model_to_grid_config_rec ( grid%id , model_config_rec , config_flags )
199
200! here we check to see if the boundary conditions are set properly
201
202   CALL boundary_condition_check( config_flags, bdyzone, error, grid%id )
203
204   moisture_init = .true.
205
206    grid%itimestep=0
207
208#ifdef DM_PARALLEL
209   CALL wrf_dm_bcast_bytes( icm , IWORDSIZE )
210   CALL wrf_dm_bcast_bytes( jcm , IWORDSIZE )
211#endif
212
213    CALL nl_set_mminlu(1, '    ')
214    CALL nl_set_iswater(1,0)
215    CALL nl_set_cen_lat(1,40.)
216    CALL nl_set_cen_lon(1,-105.)
217    CALL nl_set_truelat1(1,0.)
218    CALL nl_set_truelat2(1,0.)
219    CALL nl_set_moad_cen_lat (1,0.)
220    CALL nl_set_stand_lon (1,0.)
221    CALL nl_set_map_proj(1,0)
222
223
224!  here we initialize data we currently is not initialized
225!  in the input data
226
227    DO j = jts, jte
228      DO i = its, ite
229         grid%msft(i,j)     = 1.
230         grid%msfu(i,j)     = 1.
231         grid%msfv(i,j)     = 1.
232         grid%sina(i,j)     = 0.
233         grid%cosa(i,j)     = 1.
234         grid%e(i,j)        = 0.
235         grid%f(i,j)        = 0.         !! MARS: put coriolis here if needed
236
237      END DO
238   END DO
239
240    DO j = jts, jte
241    DO k = kts, kte
242      DO i = its, ite
243         grid%em_ww(i,k,j)     = 0.
244      END DO
245   END DO
246   END DO
247
248   grid%step_number = 0
249
250! set up the grid
251
252   IF (stretch_grid) THEN ! exponential stretch for eta (nearly constant dz)
253     DO k=1, kde
254      grid%em_znw(k) = (exp(-(k-1)/float(kde-1)/z_scale) - exp(-1./z_scale))/ &
255                                (1.-exp(-1./z_scale))
256     ENDDO
257   ELSE
258     DO k=1, kde
259      grid%em_znw(k) = 1. - float(k-1)/float(kde-1)
260     ENDDO
261   ENDIF
262
263!!MARS
264!!MARS
265  open(unit=12,file='levels',form='formatted',status='old')
266  rewind(12)
267  DO k=1, kde
268  read(12,*) grid%em_znw(k)
269  write(6,*) 'read level ', k,grid%em_znw(k)
270  ENDDO
271  close(12)
272!!MARS
273!!MARS
274 
275
276!   print *, 'coucou'
277!   stop
278
279   DO k=1, kde-1
280    grid%em_dnw(k) = grid%em_znw(k+1) - grid%em_znw(k)
281    grid%em_rdnw(k) = 1./grid%em_dnw(k)
282    grid%em_znu(k) = 0.5*(grid%em_znw(k+1)+grid%em_znw(k))
283   ENDDO
284   DO k=2, kde-1
285    grid%em_dn(k) = 0.5*(grid%em_dnw(k)+grid%em_dnw(k-1))
286    grid%em_rdn(k) = 1./grid%em_dn(k)
287    grid%em_fnp(k) = .5* grid%em_dnw(k  )/grid%em_dn(k)
288    grid%em_fnm(k) = .5* grid%em_dnw(k-1)/grid%em_dn(k)
289   ENDDO
290
291   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)
292   cof2 =     grid%em_dn(2)        /(grid%em_dn(2)+grid%em_dn(3))*grid%em_dnw(1)/grid%em_dn(3)
293   grid%cf1  = grid%em_fnp(2) + cof1
294   grid%cf2  = grid%em_fnm(2) - cof1 - cof2
295   grid%cf3  = cof2       
296
297   grid%cfn  = (.5*grid%em_dnw(kde-1)+grid%em_dn(kde-1))/grid%em_dn(kde-1)
298   grid%cfn1 = -.5*grid%em_dnw(kde-1)/grid%em_dn(kde-1)
299   grid%rdx = 1./config_flags%dx
300   grid%rdy = 1./config_flags%dy
301
302!  get the sounding from the ascii sounding file, first get dry sounding and
303!  calculate base state
304
305  dry_sounding = .true.
306  IF ( wrf_dm_on_monitor() ) THEN
307  write(6,*) ' getting dry sounding for base state '
308
309  CALL get_sounding( zk, p_in, pd_in, theta, rho, u, v, qv, dry_sounding, nl_max, nl_in )
310  ENDIF
311  CALL wrf_dm_bcast_real( zk , nl_max )
312  CALL wrf_dm_bcast_real( p_in , nl_max )
313  CALL wrf_dm_bcast_real( pd_in , nl_max )
314  CALL wrf_dm_bcast_real( theta , nl_max )
315  CALL wrf_dm_bcast_real( rho , nl_max )
316  CALL wrf_dm_bcast_real( u , nl_max )
317  CALL wrf_dm_bcast_real( v , nl_max )
318  CALL wrf_dm_bcast_real( qv , nl_max )
319  CALL wrf_dm_bcast_integer ( nl_in , 1 )
320
321  write(6,*) ' returned from reading sounding, nl_in is ',nl_in
322
323!  find ptop for the desired ztop (ztop is input from the namelist),
324!  and find surface pressure
325
326  grid%p_top = interp_0( p_in, zk, config_flags%ztop, nl_in )
327
328!!MARS
329!!MARS
330  open(unit=14,file='input_coord',form='formatted',status='old')
331  rewind(14)
332  read(14,*) lon_input
333  read(14,*) lat_input
334  close(14)
335  write(6,*) ' lon is ',lon_input
336  write(6,*) ' lat is ',lat_input
337!!MARS
338!!MARS
339
340!!MARS
341!!MARS
342  open(unit=18,file='input_more',form='formatted',status='old')
343  rewind(18)
344  read(18,*) alt_input, tsurf_input
345  close(18)
346  write(6,*) ' alt is ',alt_input
347  write(6,*) ' tsurf is ',tsurf_input
348!!MARS
349!!MARS
350
351  DO j=jts,jte
352  DO i=its,ite
353!!MARS
354    grid%ht(i,j) = alt_input
355    grid%tsk(i,j) = tsurf_input
356!!MARS
357    grid%xlat(i,j) = lat_input
358    grid%xlong(i,j) = lon_input
359    grid%mars_emiss(i,j)=0.95
360    grid%mars_cice(i,j)=0.
361    grid%slpx(i,j) = 0.
362    grid%slpy(i,j) = 0.
363!!MARS
364  ENDDO
365  ENDDO
366
367!   print *, 'coucou'
368!   stop
369
370
371  xs=ide/2 -3
372  xs=ids   -3
373  xe=xs + 6
374  ys=jde/2 -3
375  ye=ys + 6
376  mtn_ht = 500
377#ifdef MTN
378  DO j=max(ys,jds),min(ye,jde-1)
379  DO i=max(xs,ids),min(xe,ide-1)
380     grid%ht(i,j) = mtn_ht * 0.25 * &
381               ( 1. + COS ( 2*pi/(xe-xs) * ( i-xs ) + pi ) ) * &
382               ( 1. + COS ( 2*pi/(ye-ys) * ( j-ys ) + pi ) )
383  ENDDO
384  ENDDO
385#endif
386#ifdef EW_RIDGE
387  DO j=max(ys,jds),min(ye,jde-1)
388  DO i=ids,ide
389     grid%ht(i,j) = mtn_ht * 0.50 * &
390               ( 1. + COS ( 2*pi/(ye-ys) * ( j-ys ) + pi ) )
391  ENDDO
392  ENDDO
393#endif
394#ifdef NS_RIDGE
395  DO j=jds,jde
396  DO i=max(xs,ids),min(xe,ide-1)
397     grid%ht(i,j) = mtn_ht * 0.50 * &
398               ( 1. + COS ( 2*pi/(xe-xs) * ( i-xs ) + pi ) )
399  ENDDO
400  ENDDO
401#endif
402  DO j=jts,jte
403  DO i=its,ite
404    grid%em_phb(i,1,j) = g * grid%ht(i,j)
405    grid%em_ph0(i,1,j) = g * grid%ht(i,j)
406  ENDDO
407  ENDDO
408
409  DO J = jts, jte
410  DO I = its, ite
411
412    p_surf = interp_0( p_in, zk, grid%em_phb(i,1,j)/g, nl_in )
413    grid%em_mub(i,j) = p_surf-grid%p_top
414
415!  this is dry hydrostatic sounding (base state), so given grid%em_p (coordinate),
416!  interp theta (from interp) and compute 1/rho from eqn. of state
417
418    DO K = 1, kte-1
419      p_level = grid%em_znu(k)*(p_surf - grid%p_top) + grid%p_top
420      grid%em_pb(i,k,j) = p_level
421      grid%em_t_init(i,k,j) = interp_0( theta, p_in, p_level, nl_in ) - t0
422      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
423    ENDDO
424
425!  calc hydrostatic balance (alternatively we could interp the geopotential from the
426!  sounding, but this assures that the base state is in exact hydrostatic balance with
427!  respect to the model eqns.
428
429    DO k  = 2,kte
430      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)
431    ENDDO
432
433  ENDDO
434  ENDDO
435
436  IF ( wrf_dm_on_monitor() ) THEN
437    write(6,*) ' ptop is ',grid%p_top
438    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
439  ENDIF
440
441!  calculate full state for each column - this includes moisture.
442
443  write(6,*) ' getting moist sounding for full state '
444  dry_sounding = .false.
445  CALL get_sounding( zk, p_in, pd_in, theta, rho, u, v, qv, dry_sounding, nl_max, nl_in )
446
447  DO J = jts, min(jde-1,jte)
448  DO I = its, min(ide-1,ite)
449
450!  At this point grid%p_top is already set. find the DRY mass in the column
451!  by interpolating the DRY pressure. 
452
453   pd_surf = interp_0( pd_in, zk, grid%em_phb(i,1,j)/g, nl_in )
454
455!  compute the perturbation mass and the full mass
456
457    grid%em_mu_1(i,j) = pd_surf-grid%p_top - grid%em_mub(i,j)
458    grid%em_mu_2(i,j) = grid%em_mu_1(i,j)
459    grid%em_mu0(i,j) = grid%em_mu_1(i,j) + grid%em_mub(i,j)
460
461! given the dry pressure and coordinate system, interp the potential
462! temperature and qv
463
464    do k=1,kde-1
465
466      p_level = grid%em_znu(k)*(pd_surf - grid%p_top) + grid%p_top
467
468      moist(i,k,j,P_QV) = interp_0( qv, pd_in, p_level, nl_in )
469      grid%em_t_1(i,k,j)          = interp_0( theta, pd_in, p_level, nl_in ) - t0
470      grid%em_t_2(i,k,j)          = grid%em_t_1(i,k,j)
471     
472
473    enddo
474
475!  integrate the hydrostatic equation (from the RHS of the bigstep
476!  vertical momentum equation) down from the top to get grid%em_p.
477!  first from the top of the model to the top pressure
478
479    k = kte-1  ! top level
480
481    qvf1 = 0.5*(moist(i,k,j,P_QV)+moist(i,k,j,P_QV))
482    qvf2 = 1./(1.+qvf1)
483    qvf1 = qvf1*qvf2
484
485!    grid%em_p(i,k,j) = - 0.5*grid%em_mu_1(i,j)/grid%em_rdnw(k)
486    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
487    qvf = 1. + rvovrd*moist(i,k,j,P_QV)
488    grid%em_alt(i,k,j) = (r_d/p1000mb)*(grid%em_t_1(i,k,j)+t0)*qvf* &
489                (((grid%em_p(i,k,j)+grid%em_pb(i,k,j))/p1000mb)**cvpm)
490    grid%em_al(i,k,j) = grid%em_alt(i,k,j) - grid%em_alb(i,k,j)
491
492!  down the column
493
494    do k=kte-2,1,-1
495      qvf1 = 0.5*(moist(i,k,j,P_QV)+moist(i,k+1,j,P_QV))
496      qvf2 = 1./(1.+qvf1)
497      qvf1 = qvf1*qvf2
498      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)
499      qvf = 1. + rvovrd*moist(i,k,j,P_QV)
500      grid%em_alt(i,k,j) = (r_d/p1000mb)*(grid%em_t_1(i,k,j)+t0)*qvf* &
501                  (((grid%em_p(i,k,j)+grid%em_pb(i,k,j))/p1000mb)**cvpm)
502      grid%em_al(i,k,j) = grid%em_alt(i,k,j) - grid%em_alb(i,k,j)
503    enddo
504
505!  this is the hydrostatic equation used in the model after the
506!  small timesteps.  In the model, grid%em_al (inverse density)
507!  is computed from the geopotential.
508
509
510    grid%em_ph_1(i,1,j) = 0.
511    DO k  = 2,kte
512      grid%em_ph_1(i,k,j) = grid%em_ph_1(i,k-1,j) - (1./grid%em_rdnw(k-1))*(       &
513                   (grid%em_mub(i,j)+grid%em_mu_1(i,j))*grid%em_al(i,k-1,j)+ &
514                    grid%em_mu_1(i,j)*grid%em_alb(i,k-1,j)  )
515                                                   
516      grid%em_ph_2(i,k,j) = grid%em_ph_1(i,k,j)
517      grid%em_ph0(i,k,j) = grid%em_ph_1(i,k,j) + grid%em_phb(i,k,j)
518    ENDDO
519
520    IF ( wrf_dm_on_monitor() ) THEN
521    if((i==2) .and. (j==2)) then
522     write(6,*) ' grid%em_ph_1 calc ',grid%em_ph_1(2,1,2),grid%em_ph_1(2,2,2),&
523                              grid%em_mu_1(2,2)+grid%em_mub(2,2),grid%em_mu_1(2,2), &
524                              grid%em_alb(2,1,2),grid%em_al(1,2,1),grid%em_rdnw(1)
525    endif
526    ENDIF
527
528  ENDDO
529  ENDDO
530
531!#if 0
532
533!  thermal perturbation to kick off convection
534
535  write(6,*) ' nxc, nyc for perturbation ',nxc,nyc
536  write(6,*) ' delt for perturbation ',delt
537
538  DO J = jts, min(jde-1,jte)
539!    yrad = config_flags%dy*float(j-nyc)/10000.
540!! NO FUCKIN' BUBBLE !!
541yrad = 0.
542    DO I = its, min(ide-1,ite)
543!      xrad = config_flags%dx*float(i-nxc)/10000.
544!! NO FUCKIN' BUBBLE !!
545xrad = 0.
546      DO K = 1, kte-1
547
548!  put in preturbation theta (bubble) and recalc density.  note,
549!  the mass in the column is not changing, so when theta changes,
550!  we recompute density and geopotential
551
552!        zrad = 0.5*(grid%em_ph_1(i,k,j)+grid%em_ph_1(i,k+1,j)  &
553!                   +grid%em_phb(i,k,j)+grid%em_phb(i,k+1,j))/g
554!        zrad = (zrad-1500.)/1500.
555!! NO FUCKIN' BUBBLE !!
556zrad = 0.       
557        RAD=SQRT(xrad*xrad+yrad*yrad+zrad*zrad)
558        IF(RAD <= 1.) THEN
559!! NO FUCKIN' BUBBLE !!
560!           grid%em_t_1(i,k,j)=grid%em_t_1(i,k,j)+delt*COS(.5*PI*RAD)**2
561           grid%em_t_2(i,k,j)=grid%em_t_1(i,k,j)
562           qvf = 1. + rvovrd*moist(i,k,j,P_QV)
563           grid%em_alt(i,k,j) = (r_d/p1000mb)*(grid%em_t_1(i,k,j)+t0)*qvf* &
564                        (((grid%em_p(i,k,j)+grid%em_pb(i,k,j))/p1000mb)**cvpm)
565           grid%em_al(i,k,j) = grid%em_alt(i,k,j) - grid%em_alb(i,k,j)
566        ENDIF
567      ENDDO
568
569!  rebalance hydrostatically
570
571      DO k  = 2,kte
572        grid%em_ph_1(i,k,j) = grid%em_ph_1(i,k-1,j) - (1./grid%em_rdnw(k-1))*(       &
573                     (grid%em_mub(i,j)+grid%em_mu_1(i,j))*grid%em_al(i,k-1,j)+ &
574                      grid%em_mu_1(i,j)*grid%em_alb(i,k-1,j)  )
575                                                   
576        grid%em_ph_2(i,k,j) = grid%em_ph_1(i,k,j)
577        grid%em_ph0(i,k,j) = grid%em_ph_1(i,k,j) + grid%em_phb(i,k,j)
578      ENDDO
579
580    ENDDO
581  ENDDO
582
583!#endif
584
585   IF ( wrf_dm_on_monitor() ) THEN
586   write(6,*) ' grid%em_mu_1 from comp ', grid%em_mu_1(1,1)
587   write(6,*) ' full state sounding from comp, ph, grid%em_p, grid%em_al, grid%em_t_1, qv '
588   do k=1,kde-1
589     write(6,'(i3,1x,5(1x,1pe10.3))') k, grid%em_ph_1(1,k,1)+grid%em_phb(1,k,1), &
590                                      grid%em_p(1,k,1)+grid%em_pb(1,k,1), grid%em_alt(1,k,1), &
591                                      grid%em_t_1(1,k,1)+t0, moist(1,k,1,P_QV)
592   enddo
593
594   write(6,*) ' pert state sounding from comp, grid%em_ph_1, pp, alp, grid%em_t_1, qv '
595   do k=1,kde-1
596     write(6,'(i3,1x,5(1x,1pe10.3))') k, grid%em_ph_1(1,k,1), &
597                                      grid%em_p(1,k,1), grid%em_al(1,k,1), &
598                                      grid%em_t_1(1,k,1), moist(1,k,1,P_QV)
599   enddo
600   ENDIF
601
602! interp v
603
604  DO J = jts, jte
605  DO I = its, min(ide-1,ite)
606
607    IF (j == jds) THEN
608      z_at_v = grid%em_phb(i,1,j)/g
609    ELSE IF (j == jde) THEN
610      z_at_v = grid%em_phb(i,1,j-1)/g
611    ELSE
612      z_at_v = 0.5*(grid%em_phb(i,1,j)+grid%em_phb(i,1,j-1))/g
613    END IF
614    p_surf = interp_0( p_in, zk, z_at_v, nl_in )
615
616    DO K = 1, kte-1
617      p_level = grid%em_znu(k)*(p_surf - grid%p_top) + grid%p_top
618      grid%em_v_1(i,k,j) = interp_0( v, p_in, p_level, nl_in )
619      grid%em_v_2(i,k,j) = grid%em_v_1(i,k,j)
620    ENDDO
621
622  ENDDO
623  ENDDO
624
625! interp u
626
627  DO J = jts, min(jde-1,jte)
628  DO I = its, ite
629
630    IF (i == ids) THEN
631      z_at_u = grid%em_phb(i,1,j)/g
632    ELSE IF (i == ide) THEN
633      z_at_u = grid%em_phb(i-1,1,j)/g
634    ELSE
635      z_at_u = 0.5*(grid%em_phb(i,1,j)+grid%em_phb(i-1,1,j))/g
636    END IF
637
638    p_surf = interp_0( p_in, zk, z_at_u, nl_in )
639
640    DO K = 1, kte-1
641      p_level = grid%em_znu(k)*(p_surf - grid%p_top) + grid%p_top
642      grid%em_u_1(i,k,j) = interp_0( u, p_in, p_level, nl_in )
643      grid%em_u_2(i,k,j) = grid%em_u_1(i,k,j)
644    ENDDO
645
646  ENDDO
647  ENDDO
648
649!  set w
650
651  DO J = jts, min(jde-1,jte)
652  DO K = kts, kte
653  DO I = its, min(ide-1,ite)
654    grid%em_w_1(i,k,j) = 0.
655    grid%em_w_2(i,k,j) = 0.
656  ENDDO
657  ENDDO
658  ENDDO
659
660!  set a few more things
661
662  DO J = jts, min(jde-1,jte)
663  DO K = kts, kte-1
664  DO I = its, min(ide-1,ite)
665    grid%h_diabatic(i,k,j) = 0.
666  ENDDO
667  ENDDO
668  ENDDO
669
670  IF ( wrf_dm_on_monitor() ) THEN
671  DO k=1,kte-1
672    grid%em_t_base(k) = grid%em_t_1(1,k,1)
673    grid%qv_base(k) = moist(1,k,1,P_QV)
674    grid%u_base(k) = grid%em_u_1(1,k,1)
675    grid%v_base(k) = grid%em_v_1(1,k,1)
676    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
677  ENDDO
678  ENDIF
679  CALL wrf_dm_bcast_real( grid%em_t_base , kte )
680  CALL wrf_dm_bcast_real( grid%qv_base , kte )
681  CALL wrf_dm_bcast_real( grid%u_base , kte )
682  CALL wrf_dm_bcast_real( grid%v_base , kte )
683  CALL wrf_dm_bcast_real( grid%z_base , kte )
684
685  DO J = jts, min(jde-1,jte)
686  DO I = its, min(ide-1,ite)
687     thtmp   = grid%em_t_2(i,1,j)+t0
688     ptmp    = grid%em_p(i,1,j)+grid%em_pb(i,1,j)
689     temp(1) = thtmp * (ptmp/p1000mb)**rcp
690     thtmp   = grid%em_t_2(i,2,j)+t0
691     ptmp    = grid%em_p(i,2,j)+grid%em_pb(i,2,j)
692     temp(2) = thtmp * (ptmp/p1000mb)**rcp
693     thtmp   = grid%em_t_2(i,3,j)+t0
694     ptmp    = grid%em_p(i,3,j)+grid%em_pb(i,3,j)
695     temp(3) = thtmp * (ptmp/p1000mb)**rcp
696
697!!MARS
698!     grid%tsk(I,J)=grid%cf1*temp(1)+grid%cf2*temp(2)+grid%cf3*temp(3)
699     grid%tmn(I,J)=grid%tsk(I,J)-0.5
700!!!MARS
701!!TODO: passer la valeur a partir des donnees
702!grid%mars_tsoil(I,:,J)=grid%tsk(I,J)
703!!!MARS
704  ENDDO
705  ENDDO
706
707 END SUBROUTINE init_domain_rk
708
709   SUBROUTINE init_module_initialize
710   END SUBROUTINE init_module_initialize
711
712!---------------------------------------------------------------------
713
714!  test driver for get_sounding
715!
716!      implicit none
717!      integer n
718!      parameter(n = 1000)
719!      real zk(n),p(n),theta(n),rho(n),u(n),v(n),qv(n),pd(n)
720!      logical dry
721!      integer nl,k
722!
723!      dry = .false.
724!      dry = .true.
725!      call get_sounding( zk, p, pd, theta, rho, u, v, qv, dry, n, nl )
726!      write(6,*) ' input levels ',nl
727!      write(6,*) ' sounding '
728!      write(6,*) '  k  height(m)  press (Pa) pd(Pa) theta (K) den(kg/m^3)  u(m/s)     v(m/s)    qv(g/g) '
729!      do k=1,nl
730!        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)
731!      enddo
732!      end
733!
734!---------------------------------------------------------------------------
735
736      subroutine get_sounding( zk, p, p_dry, theta, rho, &
737                               u, v, qv, dry, nl_max, nl_in )
738      implicit none
739
740      integer nl_max, nl_in
741      real zk(nl_max), p(nl_max), theta(nl_max), rho(nl_max), &
742           u(nl_max), v(nl_max), qv(nl_max), p_dry(nl_max)
743      logical dry
744
745      integer n
746      parameter(n=1000)
747      logical debug
748      parameter( debug = .true.)
749
750! input sounding data
751
752      real p_surf, th_surf, qv_surf
753      real pi_surf, pi(n)
754      real h_input(n), th_input(n), qv_input(n), u_input(n), v_input(n)
755
756      !! special MARS
757      real r_input(n)
758      real cp_input(n)
759      real cv_input(n)
760      real cvpm_input(n)
761      real pfile_input(n)
762      real t_input(n)
763      real rhofile_input(n)
764      !! special MARS
765
766! diagnostics
767
768      real rho_surf, p_input(n), rho_input(n)
769      real pm_input(n)  !  this are for full moist sounding
770
771! local data
772
773      real p1000mb,cv,cp,r,cvpm,g
774!      parameter (p1000mb = 1.e+05, r = 287, cp = 1003., cv = cp-r, cvpm = -cv/cp, g=9.81 )
775      parameter (p1000mb = 610., r = 192., cp = 844.6, cv = cp-r, cvpm = -cv/cp, g=3.72)
776      integer k, it, nl
777      real qvf, qvf1, dz
778
779!  first, read the sounding
780
781      call read_sounding( p_surf, th_surf, qv_surf, &
782                          h_input, th_input, qv_input, u_input, v_input, r_input, cp_input, pfile_input, t_input, rhofile_input, n, nl, debug )
783
784
785      !! special MARS
786      do k=1,nl
787        cv_input(k) = cp_input(k) - r_input(k)
788        cvpm_input(k) = - cv_input(k) / cp_input(k)
789      enddo
790      !! special MARS
791
792      if(dry) then
793       do k=1,nl
794         qv_input(k) = 0.
795       enddo
796      endif
797
798      if(debug) write(6,*) ' number of input levels = ',nl
799
800        nl_in = nl
801        if(nl_in .gt. nl_max ) then
802          write(6,*) ' too many levels for input arrays ',nl_in,nl_max
803          call wrf_error_fatal ( ' too many levels for input arrays ' )
804        end if
805
806!  compute diagnostics,
807!  first, convert qv(g/kg) to qv(g/g)
808
809      do k=1,nl
810        qv_input(k) = 0.001*qv_input(k)
811      enddo
812
813      p_surf = 100.*p_surf  ! convert to pascals
814      qvf = 1. + rvovrd*qv_input(1)
815!      rho_surf = 1./((r/p1000mb)*th_surf*qvf*((p_surf/p1000mb)**cvpm))
816!      pi_surf = (p_surf/p1000mb)**(r/cp)
817      rho_surf =  1./((r_input(1)/p1000mb)*th_surf*qvf*((p_surf/p1000mb)**cvpm_input(1)))
818      pi_surf = (p_surf/p1000mb)**(r_input(1)/cp_input(1))
819
820
821      if(debug) then
822        write(6,*) ' surface density is ',rho_surf
823        write(6,*) ' surface pi is      ',pi_surf
824      end if
825
826
827!  integrate moist sounding hydrostatically, starting from the
828!  specified surface pressure
829!  -> first, integrate from surface to lowest level
830
831          qvf = 1. + rvovrd*qv_input(1)
832          qvf1 = 1. + qv_input(1)
833          rho_input(1) = rho_surf
834          dz = h_input(1)
835          do it=1,10
836!!MARS MARS
837            pm_input(1) = p_surf !&
838!                    - 0.5*dz*(rho_surf+rho_input(1))*g*qvf1
839!            rho_input(1) = 1./((r/p1000mb)*th_input(1)*qvf*((pm_input(1)/p1000mb)**cvpm))
840!! special MARS
841             rho_input(1) = 1./((r_input(1)/p1000mb)*th_input(1)*qvf*((pm_input(1)/p1000mb)**cvpm_input(1)))
842          enddo
843
844! integrate up the column
845
846          do k=2,nl
847            rho_input(k) = rho_input(k-1)
848!            dz = h_input(k)-h_input(k-1)
849! special MARS
850            dz = r_input(k) * t_input(k) * (- p_input(k) + p_input(k-1)) / p_input(k) / g
851!            dz = - cp_input(k) * (- t_input(k) + t_input(k-1)) / g
852            qvf1 = 0.5*(2.+(qv_input(k-1)+qv_input(k)))
853            qvf = 1. + rvovrd*qv_input(k)   ! qv is in g/kg here
854 
855            do it=1,10
856              pm_input(k) = pm_input(k-1) &
857                      - 0.5*dz*(rho_input(k)+rho_input(k-1))*g*qvf1
858!              rho_input(k) = 1./((r/p1000mb)*th_input(k)*qvf*((pm_input(k)/p1000mb)**cvpm))
859!! special MARS
860              rho_input(k) = 1./((r_input(k)/p1000mb)*th_input(k)*qvf*((pm_input(k)/p1000mb)**cvpm_input(k)))
861              !!
862      print *, p_input(k), pm_input(k), ((r_input(k)/p1000mb)*th_input(k)*qvf*((pm_input(k)/p1000mb)**cvpm_input(k))), k         
863            enddo
864          enddo
865
866!  we have the moist sounding
867
868!  next, compute the dry sounding using p at the highest level from the
869!  moist sounding and integrating down.
870
871        p_input(nl) = pm_input(nl)
872
873          do k=nl-1,1,-1
874            dz = h_input(k+1)-h_input(k)
875            p_input(k) = p_input(k+1) + 0.5*dz*(rho_input(k)+rho_input(k+1))*g
876          enddo
877
878
879        do k=1,nl
880
881          zk(k) = h_input(k)
882          !p(k) = pm_input(k)
883          !p_dry(k) = p_input(k)
884          p(k) = pfile_input(k)
885          p_dry(k) = pfile_input(k)
886          theta(k) = th_input(k)
887          !rho(k) = rho_input(k)
888          rho(k) = rhofile_input(k)
889          u(k) = u_input(k)
890          v(k) = v_input(k)
891          qv(k) = qv_input(k)
892
893        enddo
894
895     if(debug) then
896      write(6,*) ' sounding '
897      write(6,*) '  k  height(m)  press (Pa) pd(Pa) theta (K) den(kg/m^3)  u(m/s)     v(m/s)    qv(g/g) '
898      do k=1,nl
899        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)
900      enddo
901
902     end if
903
904      end subroutine get_sounding
905
906!-------------------------------------------------------
907
908      subroutine read_sounding( ps,ts,qvs,h,th,qv,u,v,r,cp,p,t,rho,n,nl,debug )
909      implicit none
910      integer n,nl
911      real ps,ts,qvs,h(n),th(n),qv(n),u(n),v(n),r(n),cp(n),p(n),t(n),rho(n)
912      logical end_of_file
913      logical debug
914
915      integer k
916
917      open(unit=10,file='input_sounding',form='formatted',status='old')
918      rewind(10)
919      read(10,*) ps, ts, qvs
920      if(debug) then
921        write(6,*) ' input sounding surface parameters '
922        write(6,*) ' surface pressure (mb) ',ps
923        write(6,*) ' surface pot. temp (K) ',ts
924        write(6,*) ' surface mixing ratio (g/kg) ',qvs
925      end if
926
927      end_of_file = .false.
928      k = 0
929
930      do while (.not. end_of_file)
931
932        read(10,*,end=100) h(k+1), th(k+1), qv(k+1), u(k+1), v(k+1)
933        k = k+1
934        if(debug) write(6,'(1x,i3,5(1x,e10.3))') k, h(k), th(k), qv(k), u(k), v(k)
935        go to 110
936 100    end_of_file = .true.
937 110    continue
938      enddo
939
940
941      !!! special MARS
942      open(unit=11,file='input_therm',form='formatted',status='old')
943      rewind(11)
944      end_of_file = .false.
945      k = 0
946      do while (.not. end_of_file)
947
948        read(11,*,end=101) r(k+1), cp(k+1), p(k+1), rho(k+1), t(k+1)
949        write(*,*) k, r(k+1), cp(k+1), p(k+1), rho(k+1), t(k+1)
950        k = k+1
951        go to 112
952 101    end_of_file = .true.
953 112    continue
954      enddo
955      !!! special MARS
956
957
958
959      nl = k
960
961      close(unit=10,status = 'keep')
962
963      end subroutine read_sounding
964
965END MODULE module_initialize
Note: See TracBrowser for help on using the repository browser.