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

Last change on this file since 3026 was 2089, checked in by mlefevre, 6 years ago

MESOSCALE. Add reading of water tracer for ideal mode

File size: 37.8 KB
Line 
1!IDEAL:MODEL_LAYER:INITIALIZATION
2!
3
4!  This MODULE holds the routines which are used to perform various initializations
5!  for the individual domains. 
6
7!  This MODULE CONTAINS the following routines:
8
9!  initialize_field_test - 1. Set different fields to different constant
10!                             values.  This is only a test.  If the correct
11!                             domain is not found (based upon the "id")
12!                             then a fatal error is issued.               
13
14!-----------------------------------------------------------------------
15
16MODULE module_initialize
17
18   USE module_domain
19   USE module_io_domain
20   USE module_state_description
21   USE module_model_constants
22   USE module_bc
23   USE module_timing
24   USE module_configure
25   USE module_init_utilities
26#ifdef DM_PARALLEL
27   USE module_dm
28#endif
29
30
31CONTAINS
32
33
34!-------------------------------------------------------------------
35! this is a wrapper for the solver-specific init_domain routines.
36! Also dereferences the grid variables and passes them down as arguments.
37! This is crucial, since the lower level routines may do message passing
38! and this will get fouled up on machines that insist on passing down
39! copies of assumed-shape arrays (by passing down as arguments, the
40! data are treated as assumed-size -- ie. f77 -- arrays and the copying
41! business is avoided).  Fie on the F90 designers.  Fie and a pox.
42
43   SUBROUTINE init_domain ( grid )
44
45   IMPLICIT NONE
46
47   !  Input data.
48   TYPE (domain), POINTER :: grid
49   !  Local data.
50   INTEGER                :: dyn_opt
51   INTEGER :: idum1, idum2
52
53   CALL nl_get_dyn_opt( 1,dyn_opt )
54   
55   CALL set_scalar_indices_from_config ( head_grid%id , idum1, idum2 )
56
57   IF (      dyn_opt .eq. 1 &
58        .or. dyn_opt .eq. 2 &
59        .or. dyn_opt .eq. 3 &
60                                       ) THEN
61     CALL init_domain_rk( grid &
62!
63#include <em_actual_new_args.inc>
64!
65                        )
66
67   ELSE
68     WRITE(0,*)' init_domain: unknown or unimplemented dyn_opt = ',dyn_opt
69     CALL wrf_error_fatal ( ' init_domain: unknown or unimplemented dyn_opt ' )
70   ENDIF
71
72   END SUBROUTINE init_domain
73
74!-------------------------------------------------------------------
75
76   SUBROUTINE init_domain_rk ( grid &
77!
78# include <em_dummy_new_args.inc>
79!
80)
81   IMPLICIT NONE
82
83   !  Input data.
84   TYPE (domain), POINTER :: grid
85
86# include <em_dummy_new_decl.inc>
87
88   TYPE (grid_config_rec_type)              :: config_flags
89
90   !  Local data
91   INTEGER                             ::                       &
92                                  ids, ide, jds, jde, kds, kde, &
93                                  ims, ime, jms, jme, kms, kme, &
94                                  its, ite, jts, jte, kts, kte, &
95                                  i, j, k
96
97   ! Local data
98!****Mars
99    REAL    :: x_param,y_param,rho_param,dilat
100    REAL    :: mulu, mulv, addu, addv
101!****Mars
102   INTEGER, PARAMETER :: nl_max = 1000
103   REAL, DIMENSION(nl_max) :: zk, p_in, theta, tk, rho, u, v, qv, pd_in
104   INTEGER :: nl_in
105
106
107   INTEGER :: icm,jcm, ii, im1, jj, jm1, loop, error, fid, nxc, nyc
108   REAL    :: u_mean,v_mean, f0, p_surf, p_level, qvf, z_at_v, z_at_u
109   REAL    :: xrad, yrad, zrad, rad, delt, cof1, cof2
110!   REAL, EXTERNAL :: interp_0
111   REAL    :: hm, xa
112   REAL    :: pi
113
114!  stuff from original initialization that has been dropped from the Registry
115   REAL    :: vnu, xnu, xnus, dinit0, cbh, p0_temp, t0_temp, zd, zt
116   REAL    :: qvf1, qvf2, pd_surf
117   INTEGER :: it
118   real :: thtmp, ptmp, temp(3)
119
120   LOGICAL :: moisture_init
121   LOGICAL :: stretch_grid, dry_sounding
122
123  INTEGER :: xs , xe , ys , ye
124  REAL :: mtn_ht
125   LOGICAL, EXTERNAL :: wrf_dm_on_monitor
126
127!!MARS
128 REAL :: lon_input, lat_input, alt_input, tsurf_input
129 INTEGER :: ierr
130!!MARS
131 REAL, DIMENSION(nl_max) :: profdustq,profdustn
132 REAL, DIMENSION(nl_max) :: prescribed_sw,prescribed_lw
133
134      REAL :: pfu, pfd, phm
135      INTEGER :: hypsometric_opt = 1 ! classic
136      !INTEGER :: hypsometric_opt = 2 ! Wee et al. 2012 correction
137
138#ifdef DM_PARALLEL
139#    include <em_data_calls.inc>
140#endif
141
142   call init_module_model_constants
143
144   SELECT CASE ( model_data_order )
145         CASE ( DATA_ORDER_ZXY )
146   kds = grid%sd31 ; kde = grid%ed31 ;
147   ids = grid%sd32 ; ide = grid%ed32 ;
148   jds = grid%sd33 ; jde = grid%ed33 ;
149
150   kms = grid%sm31 ; kme = grid%em31 ;
151   ims = grid%sm32 ; ime = grid%em32 ;
152   jms = grid%sm33 ; jme = grid%em33 ;
153
154   kts = grid%sp31 ; kte = grid%ep31 ;   ! note that tile is entire patch
155   its = grid%sp32 ; ite = grid%ep32 ;   ! note that tile is entire patch
156   jts = grid%sp33 ; jte = grid%ep33 ;   ! note that tile is entire patch
157         CASE ( DATA_ORDER_XYZ )
158   ids = grid%sd31 ; ide = grid%ed31 ;
159   jds = grid%sd32 ; jde = grid%ed32 ;
160   kds = grid%sd33 ; kde = grid%ed33 ;
161
162   ims = grid%sm31 ; ime = grid%em31 ;
163   jms = grid%sm32 ; jme = grid%em32 ;
164   kms = grid%sm33 ; kme = grid%em33 ;
165
166   its = grid%sp31 ; ite = grid%ep31 ;   ! note that tile is entire patch
167   jts = grid%sp32 ; jte = grid%ep32 ;   ! note that tile is entire patch
168   kts = grid%sp33 ; kte = grid%ep33 ;   ! note that tile is entire patch
169         CASE ( DATA_ORDER_XZY )
170   ids = grid%sd31 ; ide = grid%ed31 ;
171   kds = grid%sd32 ; kde = grid%ed32 ;
172   jds = grid%sd33 ; jde = grid%ed33 ;
173
174   ims = grid%sm31 ; ime = grid%em31 ;
175   kms = grid%sm32 ; kme = grid%em32 ;
176   jms = grid%sm33 ; jme = grid%em33 ;
177
178   its = grid%sp31 ; ite = grid%ep31 ;   ! note that tile is entire patch
179   kts = grid%sp32 ; kte = grid%ep32 ;   ! note that tile is entire patch
180   jts = grid%sp33 ; jte = grid%ep33 ;   ! note that tile is entire patch
181
182   END SELECT
183
184!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
185!!MARS : NOFILE      === no mountain
186!!MARS : FILE xa=0.  === a linear slope with elevation hm
187!!MARS : FILE        === mountain height hm, width xa
188  open(unit=22,file='ze_hill',form='formatted',status='old',iostat=ierr)
189  IF (ierr .eq. 0) THEN
190    rewind(22)
191    read(22,*) hm, xa
192    write(6,*) 'hm, xa ', hm, xa
193    close(22)
194  ENDIF
195!!MARS
196!!MARS
197!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
198
199
200   icm = ide/2
201!****Mars
202   jcm = jde/2
203!
204!   xa1  = 5000./500.
205!   xal1 = 4000./500.
206!   pii  = 2.*asin(1.0)
207!   hm1  = 250.
208!!   hm1  = 1000.
209!****Mars
210
211   delt = 3.
212!   delt = 10.
213
214!****Mars
215   stretch_grid = .true.
216!   stretch_grid = .false.
217!****Mars
218!   z_scale = .50
219   z_scale = .40
220   pi = 2.*asin(1.0)
221   write(6,*) ' pi is ',pi
222   nxc = (ide-ids)/2
223   nyc = (jde-jds)/2
224
225!!!MARS
226!!!MARS
227!  open(unit=16,file='input_vert',form='formatted',status='old')
228!  rewind(16)
229!  read(16,*) delt, z_scale
230!  write(6,*) 'delt, z_scale are  ', delt, z_scale
231!  close(16)
232!!!MARS
233!!!MARS
234
235   CALL model_to_grid_config_rec ( grid%id , model_config_rec , config_flags )
236
237! here we check to see if the boundary conditions are set properly
238
239   CALL boundary_condition_check( config_flags, bdyzone, error, grid%id )
240
241   moisture_init = .true.
242
243    grid%itimestep=0
244
245#ifdef DM_PARALLEL
246   CALL wrf_dm_bcast_bytes( icm , IWORDSIZE )
247   CALL wrf_dm_bcast_bytes( jcm , IWORDSIZE )
248#endif
249
250    CALL nl_set_mminlu(1,'    ')
251    CALL nl_set_iswater(1,0)
252    CALL nl_set_cen_lat(1,40.)
253    CALL nl_set_cen_lon(1,-105.)
254    CALL nl_set_truelat1(1,0.)
255    CALL nl_set_truelat2(1,0.)
256    CALL nl_set_moad_cen_lat (1,0.)
257    CALL nl_set_stand_lon (1,0.)
258    CALL nl_set_map_proj(1,0)
259
260
261!  here we initialize data we currently is not initialized
262!  in the input data
263
264    DO j = jts, jte
265      DO i = its, ite
266         grid%msft(i,j)     = 1.
267         grid%msfu(i,j)     = 1.
268         grid%msfv(i,j)     = 1.
269         grid%sina(i,j)     = 0.
270         grid%cosa(i,j)     = 1.
271         grid%e(i,j)        = 0.
272         grid%f(i,j)        = 0.
273
274      END DO
275   END DO
276
277    DO j = jts, jte
278    DO k = kts, kte
279      DO i = its, ite
280         grid%em_ww(i,k,j)     = 0.
281      END DO
282   END DO
283   END DO
284
285   grid%step_number = 0
286
287!! set up the grid
288!
289!   IF (stretch_grid) THEN ! exponential stretch for eta (nearly constant dz)
290!     DO k=1, kde
291!      grid%em_znw(k) = (exp(-(k-1)/float(kde-1)/z_scale) - exp(-1./z_scale))/ &
292!                                (1.-exp(-1./z_scale))
293!     ENDDO
294!   ELSE
295!     DO k=1, kde
296!      grid%em_znw(k) = 1. - float(k-1)/float(kde-1)
297!     ENDDO
298!   ENDIF
299
300!!MARS
301!!MARS
302  open(unit=12,file='levels',form='formatted',status='old')
303  rewind(12)
304  DO k=1, kde
305  read(12,*) grid%em_znw(k)
306  write(6,*) 'read level ', k,grid%em_znw(k)
307  ENDDO
308  close(12)
309!!MARS
310!!MARS
311   
312
313   DO k=1, kde-1
314    grid%em_dnw(k) = grid%em_znw(k+1) - grid%em_znw(k)
315    grid%em_rdnw(k) = 1./grid%em_dnw(k)
316    grid%em_znu(k) = 0.5*(grid%em_znw(k+1)+grid%em_znw(k))
317   ENDDO
318   DO k=2, kde-1
319    grid%em_dn(k) = 0.5*(grid%em_dnw(k)+grid%em_dnw(k-1))
320    grid%em_rdn(k) = 1./grid%em_dn(k)
321    grid%em_fnp(k) = .5* grid%em_dnw(k  )/grid%em_dn(k)
322    grid%em_fnm(k) = .5* grid%em_dnw(k-1)/grid%em_dn(k)
323   ENDDO
324
325   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)
326   cof2 =     grid%em_dn(2)        /(grid%em_dn(2)+grid%em_dn(3))*grid%em_dnw(1)/grid%em_dn(3)
327   grid%cf1  = grid%em_fnp(2) + cof1
328   grid%cf2  = grid%em_fnm(2) - cof1 - cof2
329   grid%cf3  = cof2       
330
331   grid%cfn  = (.5*grid%em_dnw(kde-1)+grid%em_dn(kde-1))/grid%em_dn(kde-1)
332   grid%cfn1 = -.5*grid%em_dnw(kde-1)/grid%em_dn(kde-1)
333   grid%rdx = 1./config_flags%dx
334   grid%rdy = 1./config_flags%dy
335
336!  get the sounding from the ascii sounding file, first get dry sounding and
337!  calculate base state
338
339     !!!!
340     !!!! user-modified wind speed
341     !!!!
342     mulu = 1.  !! default
343     mulv = 1.  !! default
344     addu = 0.  !! default
345     addv = 0.  !! default
346     IF (config_flags%init_MU .ne. 0.) mulu = config_flags%init_MU
347     IF (config_flags%init_MV .ne. 0.) mulv = config_flags%init_MV
348     IF (config_flags%init_U  .ne. 0.) addu = config_flags%init_U
349     IF (config_flags%init_V  .ne. 0.) addv = config_flags%init_V
350     write(6,*) ' coeff for winds: ', mulu, mulv, addu, addv
351
352  dry_sounding = .true.
353  IF ( wrf_dm_on_monitor() ) THEN
354  write(6,*) ' getting dry sounding for base state '
355
356  CALL get_sounding( zk, p_in, pd_in, theta, tk, rho, u, v, qv, dry_sounding, nl_max, nl_in, &
357                        mulu, mulv, addu, addv )
358  ENDIF
359  CALL wrf_dm_bcast_real( zk , nl_max )
360  CALL wrf_dm_bcast_real( p_in , nl_max )
361  CALL wrf_dm_bcast_real( pd_in , nl_max )
362  CALL wrf_dm_bcast_real( theta , nl_max )
363  CALL wrf_dm_bcast_real( tk , nl_max )
364  CALL wrf_dm_bcast_real( rho , nl_max )
365  CALL wrf_dm_bcast_real( u , nl_max )
366  CALL wrf_dm_bcast_real( v , nl_max )
367  CALL wrf_dm_bcast_real( qv , nl_max )
368  CALL wrf_dm_bcast_integer ( nl_in , 1 )
369
370  write(6,*) ' returned from reading sounding, nl_in is ',nl_in
371
372!  find ptop for the desired ztop (ztop is input from the namelist),
373!  and find surface pressure
374
375  grid%p_top = interp_0( p_in, zk, config_flags%ztop, nl_in )
376
377!!MARS
378!!MARS
379  open(unit=14,file='input_coord',form='formatted',status='old')
380  rewind(14)
381  read(14,*) lon_input
382  read(14,*) lat_input
383  close(14)
384  write(6,*) ' lon is ',lon_input
385  write(6,*) ' lat is ',lat_input
386!!MARS
387!!MARS
388
389!!MARS
390!!MARS
391  open(unit=18,file='input_more',form='formatted',status='old')
392  rewind(18)
393  read(18,*) alt_input, tsurf_input
394  close(18)
395  write(6,*) ' alt is ',alt_input
396  write(6,*) ' tsurf is ',tsurf_input
397!!MARS
398!!MARS
399
400  DO j=jts,jte
401  DO i=its,ite
402!!MARS
403IF (ierr .eq. 0) THEN
404  write(6,*) ' IDEALIZED TOPOGRAPHY '
405  IF (xa .ne. 0.) THEN
406      !!!2D hill
407      !grid%ht(i,j) = alt_input + hm/(1.+(float(i-icm)/xa)**2)
408      !    grid%ht(i,j) = hm1*exp(-(( float(i-icm)/xa1)**2))   &
409      !               *( (cos(pii*float(i-icm)/xal1))**2 )
410    IF (hm .gt. 0.) THEN
411      write(6,*) '3D hill. height, width: ',hm,xa
412      write(6,*) 'input sounding is out of the mountain'
413      grid%ht(i,j) = alt_input + hm/(1.+(float(i-icm)/xa)**2+(float(j-jcm)/xa)**2)
414    ELSE IF (hm .lt. 0.) THEN
415      write(6,*) '3D crater. height, width: ',hm,xa
416      write(6,*) 'input sounding is at the bottom of crater'
417      grid%ht(i,j) = (alt_input - hm) + hm/(1.+(float(i-icm)/xa)**2+(float(j-jcm)/xa)**2)
418      !! AS: cannot use same formula as hill because it would force ideal.exe to extrapolate
419      !!       which is not possible given how interp_0 is written
420    ELSE
421      write(6,*) 'Nothing. Height is 0. Flat topography'
422      grid%ht(i,j) = alt_input
423    ENDIF
424  ELSE
425    write(6,*) 'linear slope '
426    write(6,*) 'height ',hm
427    IF (hm .gt. 0.) THEN 
428      grid%ht(i,j) = alt_input + hm * float(i) 
429    ELSE IF (hm .lt. 0.) THEN
430      !! see above, crater case
431      grid%ht(i,j) = (alt_input - hm) + hm * float(i)
432    ELSE
433      write(6,*) 'Nothing. Height is 0. Flat topography'
434      grid%ht(i,j) = alt_input
435    ENDIF
436  ENDIF
437      !!!3D crater   
438      !!        grid%ht(i,j) = hm - hm/(1.+(float(i-icm)/xa)**2+(float(j-jcm)/xa)**2)
439      !!3D crater w/ rims
440      !  x_param = float(i-icm)
441      !  y_param = float(j-jcm)
442      !  dilat = xa/2
443      !  rho_param = sqrt(x_param**2 + y_param**2)
444      !  ! revolution surface ; seed is a fourth order polynom
445      !  grid%ht(i,j) = (rho_param+6*dilat)*(rho_param+10*dilat)
446      !  grid%ht(i,j) = (rho_param-6*dilat)*(rho_param-10*dilat)*grid%ht(i,j)
447      !  ! flat terrain elsewhere - smooth gradient (no abrupt fall)
448      !  grid%ht(i,j) = grid%ht(i,j)*(tanh(rho_param+7*dilat)/2 - tanh(rho_param-7*dilat)/2)
449      !  grid%ht(i,j) = hm - (hm*.4/1500)*grid%ht(i,j)/(dilat**4)
450      !           !NONONONONON 
451      !           !grid%ht(i,j) = grid%ht(i,j) + alt_input 
452      !   !if (rho_param .GE. dilat*10) ht(i,j) = hm
453ELSE
454    write(6,*) ' FLAT SURFACE '
455    grid%ht(i,j) = alt_input
456    !grid%ht(i,j) = 0.
457ENDIF
458    grid%tsk(i,j) = tsurf_input
459    grid%m_tsurf(i,j) = tsurf_input
460!!MARS
461    grid%xlat(i,j) = lat_input
462    grid%xlong(i,j) = lon_input!+float(i)*config_flags%dx/59000.
463    grid%m_emiss(i,j)=0.95
464    grid%m_co2ice(i,j)=0.
465    grid%m_h2oice(i,j)=0.
466!! >> Used for restarts only:
467    grid%m_q2(i,:,j)=0.
468    grid%m_fluxrad(i,j)=0.
469    grid%m_wstar(i,j)=0.
470!! <<
471    write(6,*) 'NOTE TO SELF. slpx and slpy set to 0 which means no slope insolation.'
472    grid%slpx(i,j) = 0.
473    grid%slpy(i,j) = 0.
474   DO k=1,config_flags%num_soil_layers
475    grid%m_tsoil(i,k,j) = 0.
476   ENDDO
477    !!! COMMENT THE LINES BELOW IF YOU DON'T WANT CORIOLIS TERMS
478    !!! cf. doc WRF2008 page 11 for e and f expressions
479    grid%e(i,j) = 2. * EOMEG * COS(pi*lat_input/180.)
480    grid%f(i,j) = 2. * EOMEG * SIN(pi*lat_input/180.)
481    write(6,*) 'CALCULATE CORIOLIS TERM', grid%f(i,j),grid%e(i,j)
482!!MARS
483  ENDDO
484  ENDDO
485
486  xs=ide/2 -3
487  xs=ids   -3
488  xe=xs + 6
489  ys=jde/2 -3
490  ye=ys + 6
491  mtn_ht = 500
492#ifdef MTN
493  DO j=max(ys,jds),min(ye,jde-1)
494  DO i=max(xs,ids),min(xe,ide-1)
495     grid%ht(i,j) = mtn_ht * 0.25 * &
496               ( 1. + COS ( 2*pi/(xe-xs) * ( i-xs ) + pi ) ) * &
497               ( 1. + COS ( 2*pi/(ye-ys) * ( j-ys ) + pi ) )
498  ENDDO
499  ENDDO
500#endif
501#ifdef EW_RIDGE
502  DO j=max(ys,jds),min(ye,jde-1)
503  DO i=ids,ide
504     grid%ht(i,j) = mtn_ht * 0.50 * &
505               ( 1. + COS ( 2*pi/(ye-ys) * ( j-ys ) + pi ) )
506  ENDDO
507  ENDDO
508#endif
509#ifdef NS_RIDGE
510  DO j=jds,jde
511  DO i=max(xs,ids),min(xe,ide-1)
512     grid%ht(i,j) = mtn_ht * 0.50 * &
513               ( 1. + COS ( 2*pi/(xe-xs) * ( i-xs ) + pi ) )
514  ENDDO
515  ENDDO
516#endif
517  DO j=jts,jte
518  DO i=its,ite
519    grid%em_phb(i,1,j) = g * grid%ht(i,j)
520    grid%em_ph0(i,1,j) = g * grid%ht(i,j)
521  ENDDO
522  ENDDO
523
524!!!dans hill_2d OK
525!    grid%em_phb(i,1,j) = g*grid%ht(i,j)
526!    grid%em_php(i,1,j) = 0.
527!    grid%em_ph0(i,1,j) = grid%em_phb(i,1,j)
528
529  DO J = jts, jte
530  DO I = its, ite
531
532    p_surf = interp_0( p_in, zk, grid%em_phb(i,1,j)/g, nl_in )
533    grid%em_mub(i,j) = p_surf-grid%p_top
534
535!  this is dry hydrostatic sounding (base state), so given grid%em_p (coordinate),
536!  interp theta (from interp) and compute 1/rho from eqn. of state
537
538    DO K = 1, kte-1
539      p_level = grid%em_znu(k)*(p_surf - grid%p_top) + grid%p_top
540      grid%em_pb(i,k,j) = p_level
541
542! OLD METHOD
543!      grid%em_t_init(i,k,j) = interp_0( theta, p_in, p_level, nl_in ) - t0
544! NEW METHOD: Wee et al. 2012
545! interpolate temperature. then convert to potential temperature.
546      grid%em_t_init(i,k,j) = interp_0( tk, p_in, p_level, nl_in )
547      !! l un ou l autre pareil
548      grid%em_t_init(i,k,j) = - t0 + (grid%em_t_init(i,k,j) * ((p1000mb/p_level)**rcp))
549      !grid%em_t_init(i,k,j) = - t0 + (grid%em_t_init(i,k,j) * ((610./p_level)**(1.0/3.9)))
550      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
551    ENDDO
552
553!  calc hydrostatic balance (alternatively we could interp the geopotential from the
554!  sounding, but this assures that the base state is in exact hydrostatic balance with
555!  respect to the model eqns.
556
557   IF (hypsometric_opt == 1) THEN
558    DO k  = 2,kte
559      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)
560    ENDDO
561   ELSE IF (hypsometric_opt == 2) THEN
562    DO k = 2,kte
563      pfu = grid%em_mub(i,j)*grid%em_znw(k)   + grid%p_top
564      pfd = grid%em_mub(i,j)*grid%em_znw(k-1)   + grid%p_top
565      phm = grid%em_mub(i,j)*grid%em_znu(k-1)   + grid%p_top
566      grid%em_phb(i,k,j) = grid%em_phb(i,k-1,j) + grid%em_alb(i,k-1,j)*phm*LOG(pfd/pfu)
567    END DO
568   END IF
569
570
571  ENDDO
572  ENDDO
573
574  IF ( wrf_dm_on_monitor() ) THEN
575    write(6,*) ' ptop is ',grid%p_top
576    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
577  ENDIF
578
579!  calculate full state for each column - this includes moisture.
580
581!!!!!MARS MARS
582!  write(6,*) ' getting moist sounding for full state '
583!  dry_sounding = .false.
584  dry_sounding = .true.
585  CALL get_sounding( zk, p_in, pd_in, theta, tk, rho, u, v, qv, dry_sounding, nl_max, nl_in, &
586                              mulu, mulv, addu, addv )
587
588  DO J = jts, min(jde-1,jte)
589  DO I = its, min(ide-1,ite)
590
591!  At this point grid%p_top is already set. find the DRY mass in the column
592!  by interpolating the DRY pressure. 
593
594   pd_surf = interp_0( pd_in, zk, grid%em_phb(i,1,j)/g, nl_in )
595
596!  compute the perturbation mass and the full mass
597
598    grid%em_mu_1(i,j) = pd_surf-grid%p_top - grid%em_mub(i,j)
599    grid%em_mu_2(i,j) = grid%em_mu_1(i,j)
600    grid%em_mu0(i,j) = grid%em_mu_1(i,j) + grid%em_mub(i,j)
601
602! given the dry pressure and coordinate system, interp the potential
603! temperature and qv
604
605    do k=1,kde-1
606
607      p_level = grid%em_znu(k)*(pd_surf - grid%p_top) + grid%p_top
608
609      moist(i,k,j,P_QV) = interp_0( qv, pd_in, p_level, nl_in )
610      grid%em_t_1(i,k,j)          = interp_0( theta, pd_in, p_level, nl_in ) - t0
611      grid%em_t_2(i,k,j)          = grid%em_t_1(i,k,j)
612     
613
614    enddo
615
616!  integrate the hydrostatic equation (from the RHS of the bigstep
617!  vertical momentum equation) down from the top to get grid%em_p.
618!  first from the top of the model to the top pressure
619
620    k = kte-1  ! top level
621
622    qvf1 = 0.5*(moist(i,k,j,P_QV)+moist(i,k,j,P_QV))
623    qvf2 = 1./(1.+qvf1)
624    qvf1 = qvf1*qvf2
625
626!    grid%em_p(i,k,j) = - 0.5*grid%em_mu_1(i,j)/grid%em_rdnw(k)
627    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
628    qvf = 1. + rvovrd*moist(i,k,j,P_QV)
629    grid%em_alt(i,k,j) = (r_d/p1000mb)*(grid%em_t_1(i,k,j)+t0)*qvf* &
630                (((grid%em_p(i,k,j)+grid%em_pb(i,k,j))/p1000mb)**cvpm)
631    grid%em_al(i,k,j) = grid%em_alt(i,k,j) - grid%em_alb(i,k,j)
632
633!  down the column
634
635    do k=kte-2,1,-1
636      qvf1 = 0.5*(moist(i,k,j,P_QV)+moist(i,k+1,j,P_QV))
637      qvf2 = 1./(1.+qvf1)
638      qvf1 = qvf1*qvf2
639      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)
640      qvf = 1. + rvovrd*moist(i,k,j,P_QV)
641      grid%em_alt(i,k,j) = (r_d/p1000mb)*(grid%em_t_1(i,k,j)+t0)*qvf* &
642                  (((grid%em_p(i,k,j)+grid%em_pb(i,k,j))/p1000mb)**cvpm)
643      grid%em_al(i,k,j) = grid%em_alt(i,k,j) - grid%em_alb(i,k,j)
644    enddo
645
646!  this is the hydrostatic equation used in the model after the
647!  small timesteps.  In the model, grid%em_al (inverse density)
648!  is computed from the geopotential.
649
650
651    grid%em_ph_1(i,1,j) = 0.
652   IF (hypsometric_opt == 1) THEN
653    DO k  = 2,kte
654      grid%em_ph_1(i,k,j) = grid%em_ph_1(i,k-1,j) - (1./grid%em_rdnw(k-1))*(       &
655                   (grid%em_mub(i,j)+grid%em_mu_1(i,j))*grid%em_al(i,k-1,j)+ &
656                    grid%em_mu_1(i,j)*grid%em_alb(i,k-1,j)  )
657                                                   
658      grid%em_ph_2(i,k,j) = grid%em_ph_1(i,k,j)
659      grid%em_ph0(i,k,j) = grid%em_ph_1(i,k,j) + grid%em_phb(i,k,j)
660    ENDDO
661   ELSE IF (hypsometric_opt == 2) THEN
662
663             ! Alternative hydrostatic eq.: dZ = -al*p*dLOG(p), where p is dry pressure.
664             ! Note that al*p approximates Rd*T and dLOG(p) does z.
665             ! Here T varies mostly linear with z, the first-order integration produces better result.
666
667               grid%em_ph_2(i,1,j) = grid%em_phb(i,1,j)
668               DO k = 2,kte
669                  pfu = grid%em_mu0(i,j)*grid%em_znw(k)   + grid%p_top
670                  pfd = grid%em_mu0(i,j)*grid%em_znw(k-1) + grid%p_top
671                  phm = grid%em_mu0(i,j)*grid%em_znu(k-1) + grid%p_top
672                  grid%em_ph_2(i,k,j) = grid%em_ph_2(i,k-1,j) + grid%em_alt(i,k-1,j)*phm*LOG(pfd/pfu)
673               END DO
674
675               DO k = 1,kte
676                  grid%em_ph_2(i,k,j) = grid%em_ph_2(i,k,j) - grid%em_phb(i,k,j)
677                  grid%em_ph_1(i,k,j) = grid%em_ph_2(i,k,j)
678               END DO
679
680   END IF
681
682
683
684    IF ( wrf_dm_on_monitor() ) THEN
685    if((i==2) .and. (j==2)) then
686     write(6,*) ' grid%em_ph_1 calc ',grid%em_ph_1(2,1,2),grid%em_ph_1(2,2,2),&
687                              grid%em_mu_1(2,2)+grid%em_mub(2,2),grid%em_mu_1(2,2), &
688                              grid%em_alb(2,1,2),grid%em_al(1,2,1),grid%em_rdnw(1)
689    endif
690    ENDIF
691
692  ENDDO
693  ENDDO
694
695!#if 0
696
697!  thermal perturbation to kick off convection
698
699  write(6,*) ' nxc, nyc for perturbation ',nxc,nyc
700  write(6,*) ' delt for perturbation ',delt
701
702  DO J = jts, min(jde-1,jte)
703    yrad = config_flags%dy*float(j-nyc)/10000.
704yrad = 0.
705    DO I = its, min(ide-1,ite)
706      xrad = config_flags%dx*float(i-nxc)/10000.
707xrad = 0.
708      DO K = 1, kte-1
709
710!  put in preturbation theta (bubble) and recalc density.  note,
711!  the mass in the column is not changing, so when theta changes,
712!  we recompute density and geopotential
713
714        zrad = 0.5*(grid%em_ph_1(i,k,j)+grid%em_ph_1(i,k+1,j)  &
715                   +grid%em_phb(i,k,j)+grid%em_phb(i,k+1,j))/g
716        zrad = (zrad-1500.)/1500.
717        RAD=SQRT(xrad*xrad+yrad*yrad+zrad*zrad)
718        IF(RAD <= 1.) THEN
719           grid%em_t_1(i,k,j)=grid%em_t_1(i,k,j)+delt*COS(.5*PI*RAD)**2
720           grid%em_t_2(i,k,j)=grid%em_t_1(i,k,j)
721           qvf = 1. + rvovrd*moist(i,k,j,P_QV)
722           grid%em_alt(i,k,j) = (r_d/p1000mb)*(grid%em_t_1(i,k,j)+t0)*qvf* &
723                        (((grid%em_p(i,k,j)+grid%em_pb(i,k,j))/p1000mb)**cvpm)
724           grid%em_al(i,k,j) = grid%em_alt(i,k,j) - grid%em_alb(i,k,j)
725        ENDIF
726      ENDDO
727
728!  rebalance hydrostatically
729
730   IF (hypsometric_opt == 1) THEN
731
732      DO k  = 2,kte
733        grid%em_ph_1(i,k,j) = grid%em_ph_1(i,k-1,j) - (1./grid%em_rdnw(k-1))*(       &
734                     (grid%em_mub(i,j)+grid%em_mu_1(i,j))*grid%em_al(i,k-1,j)+ &
735                      grid%em_mu_1(i,j)*grid%em_alb(i,k-1,j)  )
736                                                   
737        grid%em_ph_2(i,k,j) = grid%em_ph_1(i,k,j)
738        grid%em_ph0(i,k,j) = grid%em_ph_1(i,k,j) + grid%em_phb(i,k,j)
739      ENDDO
740
741   ELSE IF (hypsometric_opt == 2) THEN
742
743             ! Alternative hydrostatic eq.: dZ = -al*p*dLOG(p), where p is dry pressure.
744             ! Note that al*p approximates Rd*T and dLOG(p) does z.
745             ! Here T varies mostly linear with z, the first-order integration produces better result.
746
747               grid%em_ph_2(i,1,j) = grid%em_phb(i,1,j)
748               DO k = 2,kte
749                  pfu = grid%em_mu0(i,j)*grid%em_znw(k)   + grid%p_top
750                  pfd = grid%em_mu0(i,j)*grid%em_znw(k-1) + grid%p_top
751                  phm = grid%em_mu0(i,j)*grid%em_znu(k-1) + grid%p_top
752                  grid%em_ph_2(i,k,j) = grid%em_ph_2(i,k-1,j) + grid%em_alt(i,k-1,j)*phm*LOG(pfd/pfu)
753               END DO
754
755               DO k = 1,kte
756                  grid%em_ph_2(i,k,j) = grid%em_ph_2(i,k,j) - grid%em_phb(i,k,j)
757                  grid%em_ph_1(i,k,j) = grid%em_ph_2(i,k,j)
758               END DO
759
760   END IF
761
762
763    ENDDO
764  ENDDO
765
766!#endif
767
768   IF ( wrf_dm_on_monitor() ) THEN
769   write(6,*) ' grid%em_mu_1 from comp ', grid%em_mu_1(1,1)
770   write(6,*) ' full state sounding from comp, ph, grid%em_p, grid%em_al, grid%em_t_1, qv '
771   do k=1,kde-1
772     write(6,'(i3,1x,5(1x,1pe10.3))') k, grid%em_ph_1(1,k,1)+grid%em_phb(1,k,1), &
773                                      grid%em_p(1,k,1)+grid%em_pb(1,k,1), grid%em_alt(1,k,1), &
774                                      grid%em_t_1(1,k,1)+t0, moist(1,k,1,P_QV)
775   enddo
776
777   write(6,*) ' pert state sounding from comp, grid%em_ph_1, pp, alp, grid%em_t_1, qv '
778   do k=1,kde-1
779     write(6,'(i3,1x,5(1x,1pe10.3))') k, grid%em_ph_1(1,k,1), &
780                                      grid%em_p(1,k,1), grid%em_al(1,k,1), &
781                                      grid%em_t_1(1,k,1), moist(1,k,1,P_QV)
782   enddo
783   ENDIF
784
785! interp v
786
787  DO J = jts, jte
788  DO I = its, min(ide-1,ite)
789
790    IF (j == jds) THEN
791      z_at_v = grid%em_phb(i,1,j)/g
792    ELSE IF (j == jde) THEN
793      z_at_v = grid%em_phb(i,1,j-1)/g
794    ELSE
795      z_at_v = 0.5*(grid%em_phb(i,1,j)+grid%em_phb(i,1,j-1))/g
796    END IF
797
798    p_surf = interp_0( p_in, zk, z_at_v, nl_in )
799
800    DO K = 1, kte-1
801      p_level = grid%em_znu(k)*(p_surf - grid%p_top) + grid%p_top
802      grid%em_v_1(i,k,j) = interp_0( v, p_in, p_level, nl_in )
803      grid%em_v_2(i,k,j) = grid%em_v_1(i,k,j)
804    ENDDO
805
806  ENDDO
807  ENDDO
808
809! interp u
810
811  DO J = jts, min(jde-1,jte)
812  DO I = its, ite
813
814    IF (i == ids) THEN
815      z_at_u = grid%em_phb(i,1,j)/g
816    ELSE IF (i == ide) THEN
817      z_at_u = grid%em_phb(i-1,1,j)/g
818    ELSE
819      z_at_u = 0.5*(grid%em_phb(i,1,j)+grid%em_phb(i-1,1,j))/g
820    END IF
821
822    p_surf = interp_0( p_in, zk, z_at_u, nl_in )
823
824    DO K = 1, kte-1
825      p_level = grid%em_znu(k)*(p_surf - grid%p_top) + grid%p_top
826      grid%em_u_1(i,k,j) = interp_0( u, p_in, p_level, nl_in )
827      grid%em_u_2(i,k,j) = grid%em_u_1(i,k,j)
828    ENDDO
829
830  ENDDO
831  ENDDO
832
833!  set w
834
835  DO J = jts, min(jde-1,jte)
836  DO K = kts, kte
837  DO I = its, min(ide-1,ite)
838    grid%em_w_1(i,k,j) = 0.
839    grid%em_w_2(i,k,j) = 0.
840  ENDDO
841  ENDDO
842  ENDDO
843
844!  set a few more things
845
846  DO J = jts, min(jde-1,jte)
847  DO K = kts, kte-1
848  DO I = its, min(ide-1,ite)
849    grid%h_diabatic(i,k,j) = 0.
850  ENDDO
851  ENDDO
852  ENDDO
853
854  IF ( wrf_dm_on_monitor() ) THEN
855  DO k=1,kte-1
856    grid%em_t_base(k) = grid%em_t_1(1,k,1)
857    grid%qv_base(k) = moist(1,k,1,P_QV)
858    grid%u_base(k) = grid%em_u_1(1,k,1)
859    grid%v_base(k) = grid%em_v_1(1,k,1)
860    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
861  ENDDO
862  ENDIF
863  CALL wrf_dm_bcast_real( grid%em_t_base , kte )
864  CALL wrf_dm_bcast_real( grid%qv_base , kte )
865  CALL wrf_dm_bcast_real( grid%u_base , kte )
866  CALL wrf_dm_bcast_real( grid%v_base , kte )
867  CALL wrf_dm_bcast_real( grid%z_base , kte )
868
869  DO J = jts, min(jde-1,jte)
870  DO I = its, min(ide-1,ite)
871     thtmp   = grid%em_t_2(i,1,j)+t0
872     ptmp    = grid%em_p(i,1,j)+grid%em_pb(i,1,j)
873     temp(1) = thtmp * (ptmp/p1000mb)**rcp
874     thtmp   = grid%em_t_2(i,2,j)+t0
875     ptmp    = grid%em_p(i,2,j)+grid%em_pb(i,2,j)
876     temp(2) = thtmp * (ptmp/p1000mb)**rcp
877     thtmp   = grid%em_t_2(i,3,j)+t0
878     ptmp    = grid%em_p(i,3,j)+grid%em_pb(i,3,j)
879     temp(3) = thtmp * (ptmp/p1000mb)**rcp
880
881!!MARS
882!     grid%tsk(I,J)=grid%cf1*temp(1)+grid%cf2*temp(2)+grid%cf3*temp(3)
883     grid%tmn(I,J)=grid%tsk(I,J)-0.5
884!!!MARS
885!!TODO: passer la valeur a partir des donnees
886!grid%mars_tsoil(I,:,J)=grid%tsk(I,J)
887!!!MARS
888  ENDDO
889  ENDDO
890    IF (planet.eq."prescribed") Then
891      call read_hr(profdustq,profdustn,nl_in)
892      open(unit=17,file="prescribed_sw.txt",action="write")
893      open(unit=18,file="prescribed_lw.txt",action="write")
894      DO k=1,kte!-1
895        p_level = grid%em_znu(k)*(pd_surf - grid%p_top) + grid%p_top
896        prescribed_sw(k) = interp_0( profdustq, pd_in, p_level, nl_in )
897        prescribed_lw(k) = interp_0( profdustn, pd_in, p_level, nl_in )
898        write (17,*) prescribed_sw(k)
899        write (18,*) prescribed_lw(k)
900      ENDDO
901      close(unit=17)
902      close(unit=18)
903    ENDIF
904
905    if (      ( config_flags%mars == 1  ) &
906         .OR. ( config_flags%mars == 11 ) &
907         .OR. ( config_flags%mars == 12 ) ) then
908      print *, '**** INTERPOLATE HV **** RANK 2 in SCALAR'
909      DO k=1,kte-1
910         p_level = grid%em_znu(k)*(pd_surf - grid%p_top) + grid%p_top
911         scalar(its:ite,k,jts:jte,2) = interp_0( qv, pd_in, p_level, nl_in )
912         scalar(its:ite,k,jts:jte,3) = 0.
913           !! water ice is set to 0 (was put into water vapor when building prof
914           !from MCD)
915      ENDDO
916      print *, "WATER VAPOR",scalar(its,:,jts,2)
917    endif
918
919 END SUBROUTINE init_domain_rk
920
921   SUBROUTINE init_module_initialize
922   END SUBROUTINE init_module_initialize
923
924!---------------------------------------------------------------------
925
926!  test driver for get_sounding
927!
928!      implicit none
929!      integer n
930!      parameter(n = 1000)
931!      real zk(n),p(n),theta(n),rho(n),u(n),v(n),qv(n),pd(n)
932!      logical dry
933!      integer nl,k
934!
935!      dry = .false.
936!      dry = .true.
937!      call get_sounding( zk, p, pd, theta, rho, u, v, qv, dry, n, nl )
938!      write(6,*) ' input levels ',nl
939!      write(6,*) ' sounding '
940!      write(6,*) '  k  height(m)  press (Pa) pd(Pa) theta (K) den(kg/m^3)  u(m/s)     v(m/s)    qv(g/g) '
941!      do k=1,nl
942!        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)
943!      enddo
944!      end
945!
946!---------------------------------------------------------------------------
947
948      subroutine get_sounding( zk, p, p_dry, theta, tk, rho, &
949                               u, v, qv, dry, nl_max, nl_in, &
950                               mulu, mulv, addu, addv )
951      implicit none
952
953      integer nl_max, nl_in
954      real zk(nl_max), p(nl_max), theta(nl_max), tk(nl_max), rho(nl_max), &
955           u(nl_max), v(nl_max), qv(nl_max), p_dry(nl_max)
956      logical dry
957
958      integer n
959      parameter(n=1000)
960      logical debug
961
962!      parameter( debug = .false.)
963!****Mars
964      parameter( debug = .true.)
965      real mulu, mulv, addu, addv
966
967
968! input sounding data
969
970      real p_surf, th_surf, qv_surf
971      real pi_surf, pi(n)
972      real h_input(n), th_input(n), qv_input(n), u_input(n), v_input(n)
973
974      !! special MARS
975      real r_input(n)
976      real cp_input(n)
977      real cv_input(n)
978      real cvpm_input(n)
979      real pfile_input(n)
980      real t_input(n)
981      real rhofile_input(n)
982      !! special MARS
983
984! diagnostics
985
986      real rho_surf, p_input(n), rho_input(n)
987      real pm_input(n)  !  this are for full moist sounding
988
989! local data
990
991      !real p1000mb,cv,cp,r,cvpm,g
992!****Mars
993!      parameter (p1000mb = 1.e+05, r = 287, cp = 1003., cv = cp-r, cvpm = -cv/cp, g=9.81 )
994!      parameter (p1000mb = 610., r = 192., cp = 844.6, cv = cp-r, cvpm = -cv/cp, g=3.72)
995!      parameter (p1000mb = 610., r = 191., cp = 744.5, cv = cp-r, cvpm = -cv/cp, g=3.72)
996!****Mars
997      integer k, it, nl
998      real qvf, qvf1, dz
999
1000!  first, read the sounding
1001
1002      call read_sounding( p_surf, th_surf, qv_surf, &
1003                          h_input, th_input, qv_input, &
1004                          u_input, v_input, r_input, cp_input, &
1005                          pfile_input, t_input, rhofile_input, n, nl, debug )
1006
1007
1008      !! special MARS
1009      do k=1,nl
1010        cv_input(k) = cp_input(k) - r_input(k)
1011        cvpm_input(k) = - cv_input(k) / cp_input(k)
1012      enddo
1013      !! special MARS
1014
1015      if(dry) then
1016       do k=1,nl
1017         qv_input(k) = 0.
1018       enddo
1019      endif
1020
1021      if(debug) write(6,*) ' number of input levels = ',nl
1022
1023        nl_in = nl
1024        if(nl_in .gt. nl_max ) then
1025          write(6,*) ' too many levels for input arrays ',nl_in,nl_max
1026          call wrf_error_fatal ( ' too many levels for input arrays ' )
1027        end if
1028
1029!  compute diagnostics,
1030!  first, convert qv(g/kg) to qv(g/g)
1031
1032      do k=1,nl
1033        qv_input(k) = 0.001*qv_input(k)
1034      enddo
1035
1036      p_surf = 100.*p_surf  ! convert to pascals
1037      qvf = 1. + rvovrd*qv_input(1)
1038      rho_surf = 1./((r_d/p1000mb)*th_surf*qvf*((p_surf/p1000mb)**cvpm))
1039      pi_surf = (p_surf/p1000mb)**(rcp)
1040          !!!!!! rcp variable
1041          !rho_surf =  1./((r_input(1)/p1000mb)*th_surf*qvf*((p_surf/p1000mb)**cvpm_input(1)))
1042          !pi_surf = (p_surf/p1000mb)**(r_input(1)/cp_input(1))
1043
1044
1045      if(debug) then
1046        write(6,*) ' surface density is ',rho_surf
1047        write(6,*) ' surface pi is      ',pi_surf
1048      end if
1049
1050
1051!  integrate moist sounding hydrostatically, starting from the
1052!  specified surface pressure
1053!  -> first, integrate from surface to lowest level
1054
1055          qvf = 1. + rvovrd*qv_input(1)
1056          qvf1 = 1. + qv_input(1)
1057          rho_input(1) = rho_surf
1058          dz = h_input(1)
1059          do it=1,10
1060!!MARS MARS
1061            pm_input(1) = p_surf !&
1062!                    - dz*(0.25*rho_surf+0.75*rho_input(1))*g*qvf1  !!! BEURK
1063!                    - 0.5*dz*(rho_surf+rho_input(1))*g*qvf1  !! parce que couche 1 tres proche
1064            rho_input(1) = 1./((r_d/p1000mb)*th_input(1)*qvf*((pm_input(1)/p1000mb)**cvpm))
1065               !!!!!!! rcp variable
1066               !rho_input(1) = 1./((r_input(1)/p1000mb)*th_input(1)*qvf*((pm_input(1)/p1000mb)**cvpm_input(1)))
1067          enddo
1068
1069! integrate up the column
1070
1071          do k=2,nl
1072            rho_input(k) = rho_input(k-1)
1073            dz = h_input(k)-h_input(k-1)
1074               !!!!!!! rcp variable
1075               !dz = r_input(k) * t_input(k) * (- p_input(k) + p_input(k-1)) / p_input(k) / g
1076               !!dz = - cp_input(k) * (- t_input(k) + t_input(k-1)) / g
1077            qvf1 = 0.5*(2.+(qv_input(k-1)+qv_input(k)))
1078            qvf = 1. + rvovrd*qv_input(k)   ! qv is in g/kg here
1079
1080print *, 'input', pfile_input(k), rhofile_input(k) 
1081
1082            do it=1,10 !!ou moins??? non. !! trop de rho(k-1) donne une pression trop faible puis crash
1083                                          !! mais coeff ci-dessous vont varier la pression calculÃe
1084              pm_input(k) = pm_input(k-1) &
1085                      - dz*(0.75*rho_input(k)+0.25*rho_input(k-1))*g*qvf1
1086                      !- 0.5*dz*(rho_input(k)+rho_input(k-1))*g*qvf1
1087              rho_input(k) = 1./((r_d/p1000mb)*th_input(k)*qvf*((pm_input(k)/p1000mb)**cvpm))
1088!!
1089!! marche pas
1090!!
1091!pm_input(k) = pm_input(k-1) &
1092!        - 0.5*dz*(1./rho_input(k)+1./rho_input(k-1))*g*qvf1
1093!rho_input(k) = (r/p1000mb)*th_input(k)*qvf*((pm_input(k)/p1000mb)**cvpm)
1094              !!!!!!! rcp variable
1095              !rho_input(k) = 1./((r_input(k)/p1000mb)*th_input(k)*qvf*((pm_input(k)/p1000mb)**cvpm_input(k)))
1096              !print *, p_input(k), pm_input(k),((r_input(k)/p1000mb)*th_input(k)*qvf*((pm_input(k)/p1000mb)**cvpm_input(k))),k
1097print *, it, pm_input(k), rho_input(k), dz
1098            enddo
1099          enddo
1100
1101
1102!  we have the moist sounding
1103
1104!  next, compute the dry sounding using p at the highest level from the
1105!  moist sounding and integrating down.
1106
1107        p_input(nl) = pm_input(nl)
1108
1109          do k=nl-1,1,-1
1110            dz = h_input(k+1)-h_input(k)
1111            p_input(k) = p_input(k+1) + 0.5*dz*(rho_input(k)+rho_input(k+1))*g
1112          enddo
1113
1114
1115        do k=1,nl
1116
1117          zk(k) = h_input(k)
1118          p(k) = pm_input(k)
1119          p_dry(k) = p_input(k)
1120          theta(k) = th_input(k)
1121          tk(k) = t_input(k)
1122          rho(k) = rho_input(k)
1123          u(k) = mulu*u_input(k) + addu
1124          v(k) = mulv*v_input(k) + addv
1125          qv(k) = qv_input(k)
1126
1127         !!!! direct input from file
1128         write(6,*) '*** DIRECT INPUT FROM FILE ***'
1129         p(k) = pfile_input(k)
1130         p_dry(k) = pfile_input(k)
1131         rho(k) = rhofile_input(k)
1132
1133        enddo
1134
1135     if(debug) then
1136      write(6,*) ' sounding '
1137      write(6,*) '  k  height(m)  press (Pa) pd(Pa) theta (K) den(kg/m^3)  u(m/s)     v(m/s)    qv(g/g) '
1138      do k=1,nl
1139        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)
1140      enddo
1141
1142     end if
1143
1144      end subroutine get_sounding
1145
1146!-------------------------------------------------------
1147
1148      subroutine read_sounding( ps,ts,qvs,h,th,qv,u,v,r,cp,p,t,rho,n,nl,debug )
1149      implicit none
1150      integer n,nl
1151      real ps,ts,qvs,h(n),th(n),qv(n),u(n),v(n),r(n),cp(n),p(n),t(n),rho(n)
1152      logical end_of_file
1153      logical debug
1154
1155      integer k
1156
1157      open(unit=10,file='input_sounding',form='formatted',status='old')
1158      rewind(10)
1159      read(10,*) ps, ts, qvs
1160      if(debug) then
1161        write(6,*) ' input sounding surface parameters '
1162        write(6,*) ' surface pressure (mb) ',ps
1163        write(6,*) ' surface pot. temp (K) ',ts
1164        write(6,*) ' surface mixing ratio (g/kg) ',qvs
1165      end if
1166
1167      end_of_file = .false.
1168      k = 0
1169
1170      do while (.not. end_of_file)
1171
1172        read(10,*,end=100) h(k+1), th(k+1), qv(k+1), u(k+1), v(k+1)
1173        k = k+1
1174        if(debug) write(6,'(1x,i3,5(1x,e10.3))') k, h(k), th(k), qv(k), u(k), v(k)
1175        go to 110
1176 100    end_of_file = .true.
1177 110    continue
1178      enddo
1179
1180
1181      !!! special MARS
1182      open(unit=11,file='input_therm',form='formatted',status='old')
1183      rewind(11)
1184      end_of_file = .false.
1185      k = 0
1186      do while (.not. end_of_file)
1187
1188        read(11,*,end=101) r(k+1), cp(k+1), p(k+1), rho(k+1), t(k+1)
1189        write(*,*) k, r(k+1), cp(k+1), p(k+1), rho(k+1), t(k+1)
1190        k = k+1
1191        go to 112
1192 101    end_of_file = .true.
1193 112    continue
1194      enddo
1195      !!! special MARS
1196
1197
1198
1199      nl = k
1200
1201      close(unit=10,status = 'keep')
1202
1203      end subroutine read_sounding
1204
1205      subroutine read_hr(hr_sw,hr_lw,n)
1206      implicit none
1207      integer n
1208      real hr_sw(n),hr_lw(n)
1209      logical end_of_file
1210
1211      integer k
1212
1213! first element is the surface
1214
1215      open(unit=11,file='input_hr',form='formatted',status='old')
1216      rewind(11)
1217      end_of_file = .false.
1218      k = 0
1219      do while (.not. end_of_file)
1220
1221        read(11,*,end=102) hr_sw(k+1),hr_lw(k+1)
1222        write(*,*) k,hr_sw(k+1),hr_lw(k+1)
1223        k = k+1
1224        go to 113
1225 102    end_of_file = .true.
1226 113    continue
1227      enddo
1228
1229      close(unit=11,status = 'keep')
1230
1231      end subroutine read_hr
1232
1233END MODULE module_initialize
Note: See TracBrowser for help on using the repository browser.