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