1 !IDEAL:MODEL_LAYER:INITIALIZATION
3 ! This MODULE holds the routines which are used to perform various initializations
4 ! for the individual domains.
6 !-----------------------------------------------------------------------
8 MODULE module_initialize_ideal
12 USE module_state_description
13 USE module_model_constants
17 USE module_init_utilities
26 !-------------------------------------------------------------------
27 ! this is a wrapper for the solver-specific init_domain routines.
28 ! Also dereferences the grid variables and passes them down as arguments.
29 ! This is crucial, since the lower level routines may do message passing
30 ! and this will get fouled up on machines that insist on passing down
31 ! copies of assumed-shape arrays (by passing down as arguments, the
32 ! data are treated as assumed-size -- ie. f77 -- arrays and the copying
33 ! business is avoided). Fie on the F90 designers. Fie and a pox.
35 SUBROUTINE init_domain ( grid )
40 TYPE (domain), POINTER :: grid
42 INTEGER :: idum1, idum2
44 CALL set_scalar_indices_from_config ( head_grid%id , idum1, idum2 )
46 CALL init_domain_rk( grid &
48 #include <actual_new_args.inc>
52 END SUBROUTINE init_domain
54 !-------------------------------------------------------------------
56 SUBROUTINE init_domain_rk ( grid &
58 # include <dummy_new_args.inc>
64 TYPE (domain), POINTER :: grid
66 # include <dummy_decl.inc>
68 TYPE (grid_config_rec_type) :: config_flags
72 ids, ide, jds, jde, kds, kde, &
73 ims, ime, jms, jme, kms, kme, &
74 its, ite, jts, jte, kts, kte, &
79 INTEGER, PARAMETER :: nl_max = 1000
80 REAL, DIMENSION(nl_max) :: zk, p_in, theta, rho, u, v, qv, pd_in
83 INTEGER :: icm,jcm, ii, im1, jj, jm1, loop, error, fid, nxc, nyc
84 REAL :: u_mean,v_mean, f0, p_surf, p_level, qvf, z_at_v, z_at_u
85 REAL :: z_scale, xrad, yrad, zrad, rad, delt, cof1, cof2
86 ! REAL, EXTERNAL :: interp_0
90 ! stuff from original initialization that has been dropped from the Registry
91 REAL :: vnu, xnu, xnus, dinit0, cbh, p0_temp, t0_temp, zd, zt
92 REAL :: qvf1, qvf2, pd_surf
95 LOGICAL :: moisture_init
96 LOGICAL :: stretch_grid, dry_sounding, debug
98 ! kludge space for initial jet
100 INTEGER, parameter :: nz_jet=64, ny_jet=80
101 REAL, DIMENSION(nz_jet, ny_jet) :: u_jet, rho_jet, th_jet, z_jet
103 ! perturbation parameters
105 REAL, PARAMETER :: htbub=8000., radbub=2000000., radz=8000., tpbub=1.0
107 INTEGER :: icen, jcen
108 real :: thtmp, ptmp, temp(3)
110 SELECT CASE ( model_data_order )
111 CASE ( DATA_ORDER_ZXY )
112 kds = grid%sd31 ; kde = grid%ed31 ;
113 ids = grid%sd32 ; ide = grid%ed32 ;
114 jds = grid%sd33 ; jde = grid%ed33 ;
116 kms = grid%sm31 ; kme = grid%em31 ;
117 ims = grid%sm32 ; ime = grid%em32 ;
118 jms = grid%sm33 ; jme = grid%em33 ;
120 kts = grid%sp31 ; kte = grid%ep31 ; ! note that tile is entire patch
121 its = grid%sp32 ; ite = grid%ep32 ; ! note that tile is entire patch
122 jts = grid%sp33 ; jte = grid%ep33 ; ! note that tile is entire patch
123 CASE ( DATA_ORDER_XYZ )
124 ids = grid%sd31 ; ide = grid%ed31 ;
125 jds = grid%sd32 ; jde = grid%ed32 ;
126 kds = grid%sd33 ; kde = grid%ed33 ;
128 ims = grid%sm31 ; ime = grid%em31 ;
129 jms = grid%sm32 ; jme = grid%em32 ;
130 kms = grid%sm33 ; kme = grid%em33 ;
132 its = grid%sp31 ; ite = grid%ep31 ; ! note that tile is entire patch
133 jts = grid%sp32 ; jte = grid%ep32 ; ! note that tile is entire patch
134 kts = grid%sp33 ; kte = grid%ep33 ; ! note that tile is entire patch
135 CASE ( DATA_ORDER_XZY )
136 ids = grid%sd31 ; ide = grid%ed31 ;
137 kds = grid%sd32 ; kde = grid%ed32 ;
138 jds = grid%sd33 ; jde = grid%ed33 ;
140 ims = grid%sm31 ; ime = grid%em31 ;
141 kms = grid%sm32 ; kme = grid%em32 ;
142 jms = grid%sm33 ; jme = grid%em33 ;
144 its = grid%sp31 ; ite = grid%ep31 ; ! note that tile is entire patch
145 kts = grid%sp32 ; kte = grid%ep32 ; ! note that tile is entire patch
146 jts = grid%sp33 ; jte = grid%ep33 ; ! note that tile is entire patch
154 stretch_grid = .true.
158 write(6,*) ' pi is ',pi
162 CALL model_to_grid_config_rec ( grid%id , model_config_rec , config_flags )
164 ! here we check to see if the boundary conditions are set properly
166 CALL boundary_condition_check( config_flags, bdyzone, error, grid%id )
168 moisture_init = .true.
173 CALL wrf_dm_bcast_bytes( icm , IWORDSIZE )
174 CALL wrf_dm_bcast_bytes( jcm , IWORDSIZE )
177 CALL nl_set_mminlu(1,' ')
178 CALL nl_set_iswater(1,0)
179 CALL nl_set_cen_lat(1,40.)
180 CALL nl_set_cen_lon(1,-105.)
181 CALL nl_set_truelat1(1,0.)
182 CALL nl_set_truelat2(1,0.)
183 CALL nl_set_moad_cen_lat (1,0.)
184 CALL nl_set_stand_lon (1,0.)
185 CALL nl_set_pole_lon (1,0.)
186 CALL nl_set_pole_lat (1,90.)
187 CALL nl_set_map_proj(1,0)
190 ! here we initialize data we currently is not initialized
202 grid%msfvx_inv(i,j)= 1.
224 IF (stretch_grid) THEN ! exponential stretch for eta (nearly constant dz)
226 grid%znw(k) = (exp(-(k-1)/float(kde-1)/z_scale) - exp(-1./z_scale))/ &
227 (1.-exp(-1./z_scale))
231 grid%znw(k) = 1. - float(k-1)/float(kde-1)
236 grid%dnw(k) = grid%znw(k+1) - grid%znw(k)
237 grid%rdnw(k) = 1./grid%dnw(k)
238 grid%znu(k) = 0.5*(grid%znw(k+1)+grid%znw(k))
241 grid%dn(k) = 0.5*(grid%dnw(k)+grid%dnw(k-1))
242 grid%rdn(k) = 1./grid%dn(k)
243 grid%fnp(k) = .5* grid%dnw(k )/grid%dn(k)
244 grid%fnm(k) = .5* grid%dnw(k-1)/grid%dn(k)
247 cof1 = (2.*grid%dn(2)+grid%dn(3))/(grid%dn(2)+grid%dn(3))*grid%dnw(1)/grid%dn(2)
248 cof2 = grid%dn(2) /(grid%dn(2)+grid%dn(3))*grid%dnw(1)/grid%dn(3)
249 grid%cf1 = grid%fnp(2) + cof1
250 grid%cf2 = grid%fnm(2) - cof1 - cof2
253 grid%cfn = (.5*grid%dnw(kde-1)+grid%dn(kde-1))/grid%dn(kde-1)
254 grid%cfn1 = -.5*grid%dnw(kde-1)/grid%dn(kde-1)
255 grid%rdx = 1./config_flags%dx
256 grid%rdy = 1./config_flags%dy
258 ! get the sounding from the ascii sounding file, first get dry sounding and
259 ! calculate base state
261 write(6,*) ' reading input jet sounding '
262 call read_input_jet( u_jet, rho_jet, th_jet, z_jet, nz_jet, ny_jet )
264 write(6,*) ' getting dry sounding for base state '
265 write(6,*) ' using middle column in jet sounding, j = ',ny_jet/2
266 dry_sounding = .true.
268 dry_sounding = .true.
269 debug = .true. ! this will produce print of the sounding
270 CALL get_sounding( zk, p_in, pd_in, theta, rho, u, v, qv, dry_sounding, &
271 nl_max, nl_in, u_jet, rho_jet, th_jet, z_jet, &
272 nz_jet, ny_jet, ny_jet/2, debug )
274 write(6,*) ' returned from reading sounding, nl_in is ',nl_in
276 ! find ptop for the desired ztop (ztop is input from the namelist),
277 ! and find surface pressure
279 ! For the jet, using the middle column for the base state means that
280 ! we will be extrapolating above the highest height data to the south
283 grid%p_top = interp_0( p_in, zk, config_flags%ztop, nl_in )
286 DO i=its,ite ! flat surface
297 p_surf = interp_0( p_in, zk, grid%phb(i,1,j)/g, nl_in )
298 grid%mub(i,j) = p_surf-grid%p_top
300 ! this is dry hydrostatic sounding (base state), so given grid%p (coordinate),
301 ! interp theta (from interp) and compute 1/rho from eqn. of state
304 p_level = grid%znu(k)*(p_surf - grid%p_top) + grid%p_top
305 grid%pb(i,k,j) = p_level
306 grid%t_init(i,k,j) = interp_0( theta, p_in, p_level, nl_in ) - t0
307 grid%alb(i,k,j) = (r_d/p1000mb)*(grid%t_init(i,k,j)+t0)*(grid%pb(i,k,j)/p1000mb)**cvpm
310 ! calc hydrostatic balance (alternatively we could interp the geopotential from the
311 ! sounding, but this assures that the base state is in exact hydrostatic balance with
312 ! respect to the model eqns.
315 grid%phb(i,k,j) = grid%phb(i,k-1,j) - grid%dnw(k-1)*grid%mub(i,j)*grid%alb(i,k-1,j)
321 write(6,*) ' ptop is ',grid%p_top
322 write(6,*) ' base state grid%mub(1,1), p_surf is ',grid%mub(1,1),grid%mub(1,1)+grid%p_top
324 ! calculate full state for each column - this includes moisture.
326 write(6,*) ' getting grid%moist sounding for full state '
328 dry_sounding = .true.
329 IF (config_flags%mp_physics /= 0) dry_sounding = .false.
331 DO J = jts, min(jde-1,jte)
333 ! get sounding for this point
335 debug = .false. ! this will turn off print of the sounding
336 CALL get_sounding( zk, p_in, pd_in, theta, rho, u, v, qv, dry_sounding, &
337 nl_max, nl_in, u_jet, rho_jet, th_jet, z_jet, &
338 nz_jet, ny_jet, j, debug )
340 DO I = its, min(ide-1,ite)
342 ! we could just do the first point in "i" and copy from there, but we'll
343 ! be lazy and do all the points as if they are all, independent
345 ! At this point grid%p_top is already set. find the DRY mass in the column
346 ! by interpolating the DRY pressure.
348 pd_surf = interp_0( pd_in, zk, grid%phb(i,1,j)/g, nl_in )
350 ! compute the perturbation mass and the full mass
352 grid%mu_1(i,j) = pd_surf-grid%p_top - grid%mub(i,j)
353 grid%mu_2(i,j) = grid%mu_1(i,j)
354 grid%mu0(i,j) = grid%mu_1(i,j) + grid%mub(i,j)
356 ! given the dry pressure and coordinate system, interp the potential
361 p_level = grid%znu(k)*(pd_surf - grid%p_top) + grid%p_top
363 grid%moist(i,k,j,P_QV) = interp_0( qv, pd_in, p_level, nl_in )
364 grid%t_1(i,k,j) = interp_0( theta, pd_in, p_level, nl_in ) - t0
365 grid%t_2(i,k,j) = grid%t_1(i,k,j)
370 ! integrate the hydrostatic equation (from the RHS of the bigstep
371 ! vertical momentum equation) down from the top to get grid%p.
372 ! first from the top of the model to the top pressure
374 k = kte-1 ! top level
376 qvf1 = 0.5*(grid%moist(i,k,j,P_QV)+grid%moist(i,k,j,P_QV))
380 ! grid%p(i,k,j) = - 0.5*grid%mu_1(i,j)/grid%rdnw(k)
381 grid%p(i,k,j) = - 0.5*(grid%mu_1(i,j)+qvf1*grid%mub(i,j))/grid%rdnw(k)/qvf2
382 qvf = 1. + rvovrd*grid%moist(i,k,j,P_QV)
383 grid%alt(i,k,j) = (r_d/p1000mb)*(grid%t_1(i,k,j)+t0)*qvf* &
384 (((grid%p(i,k,j)+grid%pb(i,k,j))/p1000mb)**cvpm)
385 grid%al(i,k,j) = grid%alt(i,k,j) - grid%alb(i,k,j)
390 qvf1 = 0.5*(grid%moist(i,k,j,P_QV)+grid%moist(i,k+1,j,P_QV))
393 grid%p(i,k,j) = grid%p(i,k+1,j) - (grid%mu_1(i,j) + qvf1*grid%mub(i,j))/qvf2/grid%rdn(k+1)
394 qvf = 1. + rvovrd*grid%moist(i,k,j,P_QV)
395 grid%alt(i,k,j) = (r_d/p1000mb)*(grid%t_1(i,k,j)+t0)*qvf* &
396 (((grid%p(i,k,j)+grid%pb(i,k,j))/p1000mb)**cvpm)
397 grid%al(i,k,j) = grid%alt(i,k,j) - grid%alb(i,k,j)
400 ! this is the hydrostatic equation used in the model after the
401 ! small timesteps. In the model, grid%al (inverse density)
402 ! is computed from the geopotential.
405 grid%ph_1(i,1,j) = 0.
407 grid%ph_1(i,k,j) = grid%ph_1(i,k-1,j) - (1./grid%rdnw(k-1))*( &
408 (grid%mub(i,j)+grid%mu_1(i,j))*grid%al(i,k-1,j)+ &
409 grid%mu_1(i,j)*grid%alb(i,k-1,j) )
411 grid%ph_2(i,k,j) = grid%ph_1(i,k,j)
412 grid%ph0(i,k,j) = grid%ph_1(i,k,j) + grid%phb(i,k,j)
418 p_level = grid%znu(k)*(p_surf - grid%p_top) + grid%p_top
419 grid%u_1(i,k,j) = interp_0( u, p_in, p_level, nl_in )
420 grid%u_2(i,k,j) = grid%u_1(i,k,j)
426 ! thermal perturbation to kick off convection
428 write(6,*) ' nxc, nyc for perturbation ',nxc,nyc
429 write(6,*) ' delt for perturbation ',tpbub
431 DO J = jts, min(jde-1,jte)
432 yrad = config_flags%dy*float(j-jde/2-1)/radbub
433 DO I = its, min(ide-1,ite)
434 xrad = float(i-1)/float(ide-ids)
438 ! put in preturbation theta (bubble) and recalc density. note,
439 ! the mass in the column is not changing, so when theta changes,
440 ! we recompute density and geopotential
442 zrad = 0.5*(grid%ph_1(i,k,j)+grid%ph_1(i,k+1,j) &
443 +grid%phb(i,k,j)+grid%phb(i,k+1,j))/g
444 zrad = (zrad-htbub)/radz
445 RAD=SQRT(yrad*yrad+zrad*zrad)
447 tp = tpbub*cos(rad*piov2)*cos(rad*piov2)*cos(xrad*2*pi+pi)
448 grid%t_1(i,k,j)=grid%t_1(i,k,j)+tp
449 grid%t_2(i,k,j)=grid%t_1(i,k,j)
450 qvf = 1. + rvovrd*grid%moist(i,k,j,P_QV)
451 grid%alt(i,k,j) = (r_d/p1000mb)*(grid%t_1(i,k,j)+t0)*qvf* &
452 (((grid%p(i,k,j)+grid%pb(i,k,j))/p1000mb)**cvpm)
453 grid%al(i,k,j) = grid%alt(i,k,j) - grid%alb(i,k,j)
457 ! rebalance hydrostatically
460 grid%ph_1(i,k,j) = grid%ph_1(i,k-1,j) - (1./grid%rdnw(k-1))*( &
461 (grid%mub(i,j)+grid%mu_1(i,j))*grid%al(i,k-1,j)+ &
462 grid%mu_1(i,j)*grid%alb(i,k-1,j) )
464 grid%ph_2(i,k,j) = grid%ph_1(i,k,j)
465 grid%ph0(i,k,j) = grid%ph_1(i,k,j) + grid%phb(i,k,j)
473 write(6,*) ' grid%mu_1 from comp ', grid%mu_1(1,1)
474 write(6,*) ' pert state sounding from comp, grid%ph_1, pp, grid%al, grid%t_1, qv '
476 write(6,'(i3,1x,5(1x,1pe10.3))') k, grid%ph_1(1,k,1),grid%p(1,k,1), grid%al(1,k,1), &
477 grid%t_1(1,k,1), grid%moist(1,k,1,P_QV)
480 write(6,*) ' grid%mu_1 from comp ', grid%mu_1(1,1)
481 write(6,*) ' full state sounding from comp, ph, grid%p, grid%al, grid%t_1, qv '
482 write(6,*) ' at j = 1 '
484 write(6,'(i3,1x,5(1x,1pe10.3))') k, grid%ph_1(1,k,1)+grid%phb(1,k,1), &
485 grid%p(1,k,1)+grid%pb(1,k,1), grid%al(1,k,1)+grid%alb(1,k,1), &
486 grid%t_1(1,k,1)+t0, grid%moist(1,k,1,P_QV)
490 write(6,*) ' full state sounding from comp, ph, grid%p, grid%al, grid%t_1, qv '
491 write(6,*) ' at j = jde/2 '
493 write(6,'(i3,1x,5(1x,1pe10.3))') k, grid%ph_1(1,k,jde/2)+grid%phb(1,k,jde/2), &
494 grid%p(1,k,jde/2)+grid%pb(1,k,jde/2), grid%al(1,k,jde/2)+grid%alb(1,k,jde/2), &
495 grid%t_1(1,k,jde/2)+t0, grid%moist(1,k,jde/2,P_QV)
498 write(6,*) ' full state sounding from comp, ph, grid%p, grid%al, grid%t_1, qv '
499 write(6,*) ' at j = jde-1 '
501 write(6,'(i3,1x,5(1x,1pe10.3))') k, grid%ph_1(1,k,jde-1)+grid%phb(1,k,jde-1), &
502 grid%p(1,k,jde-1)+grid%pb(1,k,jde-1), grid%al(1,k,jde-1)+grid%alb(1,k,jde-1), &
503 grid%t_1(1,k,jde-1)+t0, grid%moist(1,k,jde-1,P_QV)
509 DO I = its, min(ide-1,ite)
513 grid%v_2(i,k,j) = grid%v_1(i,k,j)
519 ! fill out last i row for u
521 DO J = jts, min(jde-1,jte)
525 grid%u_1(i,k,j) = grid%u_1(its,k,j)
526 grid%u_2(i,k,j) = grid%u_2(its,k,j)
534 DO J = jts, min(jde-1,jte)
536 DO I = its, min(ide-1,ite)
543 ! set a few more things
545 DO J = jts, min(jde-1,jte)
547 DO I = its, min(ide-1,ite)
548 grid%h_diabatic(i,k,j) = 0.
554 grid%t_base(k) = grid%t_1(1,k,1)
555 grid%qv_base(k) = grid%moist(1,k,1,P_QV)
556 grid%u_base(k) = grid%u_1(1,k,1)
557 grid%v_base(k) = grid%v_1(1,k,1)
560 DO J = jts, min(jde-1,jte)
561 DO I = its, min(ide-1,ite)
562 thtmp = grid%t_2(i,1,j)+t0
563 ptmp = grid%p(i,1,j)+grid%pb(i,1,j)
564 temp(1) = thtmp * (ptmp/p1000mb)**rcp
565 thtmp = grid%t_2(i,2,j)+t0
566 ptmp = grid%p(i,2,j)+grid%pb(i,2,j)
567 temp(2) = thtmp * (ptmp/p1000mb)**rcp
568 thtmp = grid%t_2(i,3,j)+t0
569 ptmp = grid%p(i,3,j)+grid%pb(i,3,j)
570 temp(3) = thtmp * (ptmp/p1000mb)**rcp
572 grid%tsk(I,J)=grid%cf1*temp(1)+grid%cf2*temp(2)+grid%cf3*temp(3)
573 if (i .eq. 1) print*,'sfctem',j,temp(1),temp(2),temp(3),grid%tsk(I,J)
574 grid%tmn(I,J)=grid%tsk(I,J)-0.5
580 END SUBROUTINE init_domain_rk
582 !---------------------------------------------------------------------
584 SUBROUTINE init_module_initialize
585 END SUBROUTINE init_module_initialize
587 !---------------------------------------------------------------------
589 ! TEST DRIVER FOR "read_input_jet" and "get_sounding"
591 integer, parameter :: nz_jet=64, ny_jet=80
592 real, dimension(nz_jet,ny_jet) :: u_jet, rho_jet, &
595 real, dimension(nz_jet,ny_jet) :: zk,p,p_dry,theta,rho,u,v,qv
596 logical :: dry, debug
599 call read_input_jet( u_jet, rho_jet, th_jet, z_jet, nz_jet, ny_jet )
602 call parray( u_jet, nz_jet, ny_jet)
603 call parray( rho_jet, nz_jet, ny_jet)
604 call parray( th_jet, nz_jet, ny_jet)
613 call get_sounding( zk(:,j),p(:,j),p_dry(:,j),theta(:,j), &
614 rho(:,j),u(:,j), v(:,j), qv(:,j), &
615 dry, nz_jet, nl, u_jet, rho_jet, th_jet, &
616 z_jet, nz_jet, ny_jet, j, debug )
621 write(6,*) ' lowest level p, th, and rho, highest level p '
624 write(6,*) j, p(1,j),theta(1,j),rho(1,j), p(nz_jet,j)
625 ! write(6,*) j, p(1,j),theta(1,j)-th_jet(1,j),rho(1,j)-rho_jet(1,j)
628 call parray( p, nz_jet, ny_jet)
629 call parray( p_dry, nz_jet, ny_jet)
630 call parray( theta, nz_jet, ny_jet)
636 !---------------------------------
638 subroutine parray(a,m,n)
648 write(6,'('' dimensions m,n '',2i6)')m,n
649 call set(.05,.95,.05,.95,0.,1.,0.,1.,1)
651 call setusv('LW',2000)
652 ! CALL CONREC(a,m,m,n,cmax,cmin,cinc,-1,-638,-922)
653 CALL CONREC(b,n,n,m,0.,0.,0.,-1,-638,-922)
658 ! END TEST DRIVER FOR "read_input_jet" and "get_sounding"
661 !------------------------------------------------------------------
663 subroutine get_sounding( zk, p, p_dry, theta, rho, &
664 u, v, qv, dry, nl_max, nl_in, &
665 u_jet, rho_jet, th_jet, z_jet, &
666 nz_jet, ny_jet, j_point, debug )
669 integer nl_max, nl_in
670 real zk(nl_max), p(nl_max), theta(nl_max), rho(nl_max), &
671 u(nl_max), v(nl_max), qv(nl_max), p_dry(nl_max)
674 integer nz_jet, ny_jet, j_point
675 real, dimension(nz_jet, ny_jet) :: u_jet, rho_jet, th_jet, z_jet
681 ! input sounding data
683 real p_surf, th_surf, qv_surf
685 real h_input(n), th_input(n), qv_input(n), u_input(n), v_input(n)
689 real rho_surf, p_input(n), rho_input(n)
690 real pm_input(n) ! this are for full moist sounding
699 ! first, read the sounding
701 ! call read_sounding( p_surf, th_surf, qv_surf, &
702 ! h_input, th_input, qv_input, u_input, v_input,n, nl, debug )
704 call calc_jet_sounding( p_surf, th_surf, qv_surf, &
705 h_input, th_input, qv_input, u_input, v_input, &
706 n, nl, debug, u_jet, rho_jet, th_jet, z_jet, j_point, &
707 nz_jet, ny_jet, dry )
717 if(debug) write(6,*) ' number of input levels = ',nl
720 if(nl_in .gt. nl_max ) then
721 write(6,*) ' too many levels for input arrays ',nl_in,nl_max
722 call wrf_error_fatal ( ' too many levels for input arrays ' )
725 ! compute diagnostics,
726 ! first, convert qv(g/kg) to qv(g/g)
729 ! qv_input(k) = 0.001*qv_input(k)
731 ! p_surf = 100.*p_surf ! convert to pascals
733 qvf = 1. + rvovrd*qv_input(1)
734 rho_surf = 1./((r/p1000mb)*th_surf*qvf*((p_surf/p1000mb)**cvpm))
735 pi_surf = (p_surf/p1000mb)**(r/cp)
738 write(6,*) ' surface density is ',rho_surf
739 write(6,*) ' surface pi is ',pi_surf
743 ! integrate moist sounding hydrostatically, starting from the
744 ! specified surface pressure
745 ! -> first, integrate from surface to lowest level
747 qvf = 1. + rvovrd*qv_input(1)
748 qvf1 = 1. + qv_input(1)
749 rho_input(1) = rho_surf
752 pm_input(1) = p_surf &
753 - 0.5*dz*(rho_surf+rho_input(1))*g*qvf1
754 rho_input(1) = 1./((r/p1000mb)*th_input(1)*qvf*((pm_input(1)/p1000mb)**cvpm))
757 ! integrate up the column
760 rho_input(k) = rho_input(k-1)
761 dz = h_input(k)-h_input(k-1)
762 qvf1 = 0.5*(2.+(qv_input(k-1)+qv_input(k)))
763 qvf = 1. + rvovrd*qv_input(k) ! qv is in g/kg here
766 pm_input(k) = pm_input(k-1) &
767 - 0.5*dz*(rho_input(k)+rho_input(k-1))*g*qvf1
768 rho_input(k) = 1./((r/p1000mb)*th_input(k)*qvf*((pm_input(k)/p1000mb)**cvpm))
772 ! we have the moist sounding
774 ! next, compute the dry sounding using p at the highest level from the
775 ! moist sounding and integrating down.
777 p_input(nl) = pm_input(nl)
780 dz = h_input(k+1)-h_input(k)
781 p_input(k) = p_input(k+1) + 0.5*dz*(rho_input(k)+rho_input(k+1))*g
789 p_dry(k) = p_input(k)
790 theta(k) = th_input(k)
791 rho(k) = rho_input(k)
799 write(6,*) ' sounding '
800 write(6,*) ' k height(m) press (Pa) pd(Pa) theta (K) den(kg/m^3) u(m/s) v(m/s) qv(g/g) '
802 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)
807 end subroutine get_sounding
809 !------------------------------------------------------------------
811 subroutine calc_jet_sounding( p_surf, th_surf, qv_surf, &
812 h, th, qv, u, v, n, nl, debug, &
813 u_jet, rho_jet, th_jet, z_jet, &
814 jp, nz_jet, ny_jet, dry )
816 integer :: n, nl, jp, nz_jet, ny_jet
818 real, dimension(nz_jet, ny_jet) :: u_jet, rho_jet, th_jet, z_jet
819 real, dimension(n) :: h,th,qv,u,v
820 real :: p_surf, th_surf, qv_surf
821 logical :: debug, dry
823 real, dimension(1:nz_jet) :: rho, rel_hum, p
828 real :: tmppi, es, qvs, temperature
830 ! get sounding from column jp
836 rho(k) = rho_jet(k,jp)
843 if(h(k) .gt. 8000.) then
846 rel_hum(k)=(1.-0.90*(h(k)/8000.)**1.25)
848 rel_hum(k) = min(0.7,rel_hum(k))
856 ! next, compute pressure
859 p(k) = p1000mb*(R_d*rho(k)*th(k)/p1000mb)**cpovcv
862 ! here we adjust for fixed moisture profile
866 ! here we assume the input theta is th_v, so we reset theta accordingly
869 tmppi=(p(k)/p1000mb)**rcp
870 temperature = tmppi*th(k)
871 if (temperature .gt. svpt0) then
872 es = 1000.*svp1*exp(svp2*(temperature-svpt0)/(temperature-svp3))
873 qvs = ep_2*es/(p(k)-es)
875 es = 1000.*svp1*exp( 21.8745584*(temperature-273.16)/(temperature-7.66) )
876 qvs = ep_2*es/(p(k)-es)
878 qv(k) = rel_hum(k)*qvs
879 th(k) = th(k)/(1.+.61*qv(k))
884 ! finally, set the surface data. We'll just do a simple extrapolation
886 p_surf = 1.5*p(1) - 0.5*p(2)
887 th_surf = 1.5*th(1) - 0.5*th(2)
888 qv_surf = 1.5*qv(1) - 0.5*qv(2)
890 end subroutine calc_jet_sounding
892 !---------------------------------------------------------------------
894 SUBROUTINE read_input_jet( u, r, t, zk, nz, ny )
897 integer, intent(in) :: nz,ny
898 real, dimension(nz,ny), intent(out) :: u,r,t,zk
899 integer :: ny_in, nz_in, j,k
900 real, dimension(ny,nz) :: field_in
902 ! this code assumes it is called on processor 0 only
904 OPEN(unit=10, file='input_jet', form='unformatted', status='old' )
907 if((ny_in /= ny ) .or. (nz_in /= nz)) then
908 write(0,*) ' error in input jet dimensions '
909 write(0,*) ' ny, ny_input, nz, nz_input ', ny, ny_in, nz,nz_in
910 write(0,*) ' error exit '
911 call wrf_error_fatal ( ' error in input jet dimensions ' )
916 u(k,j) = field_in(j,k)
922 t(k,j) = field_in(j,k)
929 r(k,j) = field_in(j,k)
935 zk(k,j) = 125. + 250.*float(k-1)
939 end subroutine read_input_jet
941 END MODULE module_initialize_ideal