module_initialize_hill2d_x.F

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