module_initialize_b_wave.F
References to this file elsewhere.
1 !IDEAL:MODEL_LAYER:INITIALIZATION
2
3 ! This MODULE holds the routines which are used to perform various initializations
4 ! for the individual domains.
5
6 !-----------------------------------------------------------------------
7
8 MODULE module_initialize_ideal
9
10 USE module_domain
11 USE module_io_domain
12 USE module_state_description
13 USE module_model_constants
14 USE module_bc
15 USE module_timing
16 USE module_configure
17 USE module_init_utilities
18 #ifdef DM_PARALLEL
19 USE module_dm
20 #endif
21
22
23 CONTAINS
24
25
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.
34
35 SUBROUTINE init_domain ( grid )
36
37 IMPLICIT NONE
38
39 ! Input data.
40 TYPE (domain), POINTER :: grid
41 ! Local data.
42 INTEGER :: dyn_opt
43 INTEGER :: idum1, idum2
44
45 CALL nl_get_dyn_opt( 1, dyn_opt )
46
47 CALL set_scalar_indices_from_config ( head_grid%id , idum1, idum2 )
48
49 IF ( dyn_opt .eq. 1 &
50 .or. dyn_opt .eq. 2 &
51 .or. dyn_opt .eq. 3 &
52 ) THEN
53 CALL init_domain_rk( grid &
54 !
55 #include <em_actual_new_args.inc>
56 !
57 )
58
59 ELSE
60 WRITE(0,*)' init_domain: unknown or unimplemented dyn_opt = ',dyn_opt
61 call wrf_error_fatal ( ' init_domain: unknown or unimplemented dyn_opt ' )
62 ENDIF
63
64 END SUBROUTINE init_domain
65
66 !-------------------------------------------------------------------
67
68 SUBROUTINE init_domain_rk ( grid &
69 !
70 # include <em_dummy_new_args.inc>
71 !
72 )
73 IMPLICIT NONE
74
75 ! Input data.
76 TYPE (domain), POINTER :: grid
77
78 # include <em_dummy_decl.inc>
79
80 TYPE (grid_config_rec_type) :: config_flags
81
82 ! Local data
83 INTEGER :: &
84 ids, ide, jds, jde, kds, kde, &
85 ims, ime, jms, jme, kms, kme, &
86 its, ite, jts, jte, kts, kte, &
87 i, j, k
88
89 ! Local data
90
91 INTEGER, PARAMETER :: nl_max = 1000
92 REAL, DIMENSION(nl_max) :: zk, p_in, theta, rho, u, v, qv, pd_in
93 INTEGER :: nl_in
94
95 INTEGER :: icm,jcm, ii, im1, jj, jm1, loop, error, fid, nxc, nyc
96 REAL :: u_mean,v_mean, f0, p_surf, p_level, qvf, z_at_v, z_at_u
97 REAL :: z_scale, xrad, yrad, zrad, rad, delt, cof1, cof2
98 ! REAL, EXTERNAL :: interp_0
99 REAL :: hm
100 REAL :: pi
101
102 ! stuff from original initialization that has been dropped from the Registry
103 REAL :: vnu, xnu, xnus, dinit0, cbh, p0_temp, t0_temp, zd, zt
104 REAL :: qvf1, qvf2, pd_surf
105 INTEGER :: it
106
107 LOGICAL :: moisture_init
108 LOGICAL :: stretch_grid, dry_sounding, debug
109
110 ! kludge space for initial jet
111
112 INTEGER, parameter :: nz_jet=64, ny_jet=80
113 REAL, DIMENSION(nz_jet, ny_jet) :: u_jet, rho_jet, th_jet, z_jet
114
115 ! perturbation parameters
116
117 REAL, PARAMETER :: htbub=8000., radbub=2000000., radz=8000., tpbub=1.0
118 REAL :: piov2, tp
119 INTEGER :: icen, jcen
120 real :: thtmp, ptmp, temp(3)
121
122 #ifdef DM_PARALLEL
123 # include <em_data_calls.inc>
124 #endif
125
126 SELECT CASE ( model_data_order )
127 CASE ( DATA_ORDER_ZXY )
128 kds = grid%sd31 ; kde = grid%ed31 ;
129 ids = grid%sd32 ; ide = grid%ed32 ;
130 jds = grid%sd33 ; jde = grid%ed33 ;
131
132 kms = grid%sm31 ; kme = grid%em31 ;
133 ims = grid%sm32 ; ime = grid%em32 ;
134 jms = grid%sm33 ; jme = grid%em33 ;
135
136 kts = grid%sp31 ; kte = grid%ep31 ; ! note that tile is entire patch
137 its = grid%sp32 ; ite = grid%ep32 ; ! note that tile is entire patch
138 jts = grid%sp33 ; jte = grid%ep33 ; ! note that tile is entire patch
139 CASE ( DATA_ORDER_XYZ )
140 ids = grid%sd31 ; ide = grid%ed31 ;
141 jds = grid%sd32 ; jde = grid%ed32 ;
142 kds = grid%sd33 ; kde = grid%ed33 ;
143
144 ims = grid%sm31 ; ime = grid%em31 ;
145 jms = grid%sm32 ; jme = grid%em32 ;
146 kms = grid%sm33 ; kme = grid%em33 ;
147
148 its = grid%sp31 ; ite = grid%ep31 ; ! note that tile is entire patch
149 jts = grid%sp32 ; jte = grid%ep32 ; ! note that tile is entire patch
150 kts = grid%sp33 ; kte = grid%ep33 ; ! note that tile is entire patch
151 CASE ( DATA_ORDER_XZY )
152 ids = grid%sd31 ; ide = grid%ed31 ;
153 kds = grid%sd32 ; kde = grid%ed32 ;
154 jds = grid%sd33 ; jde = grid%ed33 ;
155
156 ims = grid%sm31 ; ime = grid%em31 ;
157 kms = grid%sm32 ; kme = grid%em32 ;
158 jms = grid%sm33 ; jme = grid%em33 ;
159
160 its = grid%sp31 ; ite = grid%ep31 ; ! note that tile is entire patch
161 kts = grid%sp32 ; kte = grid%ep32 ; ! note that tile is entire patch
162 jts = grid%sp33 ; jte = grid%ep33 ; ! note that tile is entire patch
163
164 END SELECT
165
166 piov2 = 2.*atan(1.0)
167 icen = ide/4
168 jcen = jde/2
169
170 stretch_grid = .true.
171 delt = 0.
172 z_scale = .50
173 pi = 2.*asin(1.0)
174 write(6,*) ' pi is ',pi
175 nxc = (ide-ids)/4
176 nyc = (jde-jds)/2
177
178 CALL model_to_grid_config_rec ( grid%id , model_config_rec , config_flags )
179
180 ! here we check to see if the boundary conditions are set properly
181
182 CALL boundary_condition_check( config_flags, bdyzone, error, grid%id )
183
184 moisture_init = .true.
185
186 grid%itimestep=0
187
188 #ifdef DM_PARALLEL
189 CALL wrf_dm_bcast_bytes( icm , IWORDSIZE )
190 CALL wrf_dm_bcast_bytes( jcm , IWORDSIZE )
191 #endif
192
193 CALL nl_set_mminlu(1,' ')
194 CALL nl_set_iswater(1,0)
195 CALL nl_set_cen_lat(1,40.)
196 CALL nl_set_cen_lon(1,-105.)
197 CALL nl_set_truelat1(1,0.)
198 CALL nl_set_truelat2(1,0.)
199 CALL nl_set_moad_cen_lat (1,0.)
200 CALL nl_set_stand_lon (1,0.)
201 CALL nl_set_map_proj(1,0)
202
203
204 ! here we initialize data we currently is not initialized
205 ! in the input data
206
207 DO j = jts, jte
208 DO i = its, ite
209
210 grid%ht(i,j) = 0.
211 grid%msft(i,j) = 1.
212 grid%msfu(i,j) = 1.
213 grid%msfv(i,j) = 1.
214 grid%sina(i,j) = 0.
215 grid%cosa(i,j) = 1.
216 grid%e(i,j) = 0.
217 grid%f(i,j) = 1.e-04
218
219 END DO
220 END DO
221
222 DO j = jts, jte
223 DO k = kts, kte
224 DO i = its, ite
225 grid%em_ww(i,k,j) = 0.
226 END DO
227 END DO
228 END DO
229
230 grid%step_number = 0
231
232 ! set up the grid
233
234 IF (stretch_grid) THEN ! exponential stretch for eta (nearly constant dz)
235 DO k=1, kde
236 grid%em_znw(k) = (exp(-(k-1)/float(kde-1)/z_scale) - exp(-1./z_scale))/ &
237 (1.-exp(-1./z_scale))
238 ENDDO
239 ELSE
240 DO k=1, kde
241 grid%em_znw(k) = 1. - float(k-1)/float(kde-1)
242 ENDDO
243 ENDIF
244
245 DO k=1, kde-1
246 grid%em_dnw(k) = grid%em_znw(k+1) - grid%em_znw(k)
247 grid%em_rdnw(k) = 1./grid%em_dnw(k)
248 grid%em_znu(k) = 0.5*(grid%em_znw(k+1)+grid%em_znw(k))
249 ENDDO
250 DO k=2, kde-1
251 grid%em_dn(k) = 0.5*(grid%em_dnw(k)+grid%em_dnw(k-1))
252 grid%em_rdn(k) = 1./grid%em_dn(k)
253 grid%em_fnp(k) = .5* grid%em_dnw(k )/grid%em_dn(k)
254 grid%em_fnm(k) = .5* grid%em_dnw(k-1)/grid%em_dn(k)
255 ENDDO
256
257 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)
258 cof2 = grid%em_dn(2) /(grid%em_dn(2)+grid%em_dn(3))*grid%em_dnw(1)/grid%em_dn(3)
259 grid%cf1 = grid%em_fnp(2) + cof1
260 grid%cf2 = grid%em_fnm(2) - cof1 - cof2
261 grid%cf3 = cof2
262
263 grid%cfn = (.5*grid%em_dnw(kde-1)+grid%em_dn(kde-1))/grid%em_dn(kde-1)
264 grid%cfn1 = -.5*grid%em_dnw(kde-1)/grid%em_dn(kde-1)
265 grid%rdx = 1./config_flags%dx
266 grid%rdy = 1./config_flags%dy
267
268 ! get the sounding from the ascii sounding file, first get dry sounding and
269 ! calculate base state
270
271 write(6,*) ' reading input jet sounding '
272 call read_input_jet( u_jet, rho_jet, th_jet, z_jet, nz_jet, ny_jet )
273
274 write(6,*) ' getting dry sounding for base state '
275 write(6,*) ' using middle column in jet sounding, j = ',ny_jet/2
276 dry_sounding = .true.
277
278 dry_sounding = .true.
279 debug = .true. ! this will produce print of the sounding
280 CALL get_sounding( zk, p_in, pd_in, theta, rho, u, v, qv, dry_sounding, &
281 nl_max, nl_in, u_jet, rho_jet, th_jet, z_jet, &
282 nz_jet, ny_jet, ny_jet/2, debug )
283
284 write(6,*) ' returned from reading sounding, nl_in is ',nl_in
285
286 ! find ptop for the desired ztop (ztop is input from the namelist),
287 ! and find surface pressure
288
289 ! For the jet, using the middle column for the base state means that
290 ! we will be extrapolating above the highest height data to the south
291 ! of the centerline.
292
293 grid%p_top = interp_0( p_in, zk, config_flags%ztop, nl_in )
294
295 DO j=jts,jte
296 DO i=its,ite ! flat surface
297 grid%em_phb(i,1,j) = 0.
298 grid%em_php(i,1,j) = 0.
299 grid%em_ph0(i,1,j) = 0.
300 grid%ht(i,j) = 0.
301 ENDDO
302 ENDDO
303
304 DO J = jts, jte
305 DO I = its, ite
306
307 p_surf = interp_0( p_in, zk, grid%em_phb(i,1,j)/g, nl_in )
308 grid%em_mub(i,j) = p_surf-grid%p_top
309
310 ! this is dry hydrostatic sounding (base state), so given grid%em_p (coordinate),
311 ! interp theta (from interp) and compute 1/rho from eqn. of state
312
313 DO K = 1, kte-1
314 p_level = grid%em_znu(k)*(p_surf - grid%p_top) + grid%p_top
315 grid%em_pb(i,k,j) = p_level
316 grid%em_t_init(i,k,j) = interp_0( theta, p_in, p_level, nl_in ) - t0
317 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
318 ENDDO
319
320 ! calc hydrostatic balance (alternatively we could interp the geopotential from the
321 ! sounding, but this assures that the base state is in exact hydrostatic balance with
322 ! respect to the model eqns.
323
324 DO k = 2,kte
325 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)
326 ENDDO
327
328 ENDDO
329 ENDDO
330
331 write(6,*) ' ptop is ',grid%p_top
332 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
333
334 ! calculate full state for each column - this includes moisture.
335
336 write(6,*) ' getting grid%moist sounding for full state '
337
338 dry_sounding = .true.
339 IF (config_flags%mp_physics /= 0) dry_sounding = .false.
340
341 DO J = jts, min(jde-1,jte)
342
343 ! get sounding for this point
344
345 debug = .false. ! this will turn off print of the sounding
346 CALL get_sounding( zk, p_in, pd_in, theta, rho, u, v, qv, dry_sounding, &
347 nl_max, nl_in, u_jet, rho_jet, th_jet, z_jet, &
348 nz_jet, ny_jet, j, debug )
349
350 DO I = its, min(ide-1,ite)
351
352 ! we could just do the first point in "i" and copy from there, but we'll
353 ! be lazy and do all the points as if they are all, independent
354
355 ! At this point grid%p_top is already set. find the DRY mass in the column
356 ! by interpolating the DRY pressure.
357
358 pd_surf = interp_0( pd_in, zk, grid%em_phb(i,1,j)/g, nl_in )
359
360 ! compute the perturbation mass and the full mass
361
362 grid%em_mu_1(i,j) = pd_surf-grid%p_top - grid%em_mub(i,j)
363 grid%em_mu_2(i,j) = grid%em_mu_1(i,j)
364 grid%em_mu0(i,j) = grid%em_mu_1(i,j) + grid%em_mub(i,j)
365
366 ! given the dry pressure and coordinate system, interp the potential
367 ! temperature and qv
368
369 do k=1,kde-1
370
371 p_level = grid%em_znu(k)*(pd_surf - grid%p_top) + grid%p_top
372
373 grid%moist(i,k,j,P_QV) = interp_0( qv, pd_in, p_level, nl_in )
374 grid%em_t_1(i,k,j) = interp_0( theta, pd_in, p_level, nl_in ) - t0
375 grid%em_t_2(i,k,j) = grid%em_t_1(i,k,j)
376
377
378 enddo
379
380 ! integrate the hydrostatic equation (from the RHS of the bigstep
381 ! vertical momentum equation) down from the top to get grid%em_p.
382 ! first from the top of the model to the top pressure
383
384 k = kte-1 ! top level
385
386 qvf1 = 0.5*(grid%moist(i,k,j,P_QV)+grid%moist(i,k,j,P_QV))
387 qvf2 = 1./(1.+qvf1)
388 qvf1 = qvf1*qvf2
389
390 ! grid%em_p(i,k,j) = - 0.5*grid%em_mu_1(i,j)/grid%em_rdnw(k)
391 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
392 qvf = 1. + rvovrd*grid%moist(i,k,j,P_QV)
393 grid%em_alt(i,k,j) = (r_d/p1000mb)*(grid%em_t_1(i,k,j)+t0)*qvf* &
394 (((grid%em_p(i,k,j)+grid%em_pb(i,k,j))/p1000mb)**cvpm)
395 grid%em_al(i,k,j) = grid%em_alt(i,k,j) - grid%em_alb(i,k,j)
396
397 ! down the column
398
399 do k=kte-2,1,-1
400 qvf1 = 0.5*(grid%moist(i,k,j,P_QV)+grid%moist(i,k+1,j,P_QV))
401 qvf2 = 1./(1.+qvf1)
402 qvf1 = qvf1*qvf2
403 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)
404 qvf = 1. + rvovrd*grid%moist(i,k,j,P_QV)
405 grid%em_alt(i,k,j) = (r_d/p1000mb)*(grid%em_t_1(i,k,j)+t0)*qvf* &
406 (((grid%em_p(i,k,j)+grid%em_pb(i,k,j))/p1000mb)**cvpm)
407 grid%em_al(i,k,j) = grid%em_alt(i,k,j) - grid%em_alb(i,k,j)
408 enddo
409
410 ! this is the hydrostatic equation used in the model after the
411 ! small timesteps. In the model, grid%em_al (inverse density)
412 ! is computed from the geopotential.
413
414
415 grid%em_ph_1(i,1,j) = 0.
416 DO k = 2,kte
417 grid%em_ph_1(i,k,j) = grid%em_ph_1(i,k-1,j) - (1./grid%em_rdnw(k-1))*( &
418 (grid%em_mub(i,j)+grid%em_mu_1(i,j))*grid%em_al(i,k-1,j)+ &
419 grid%em_mu_1(i,j)*grid%em_alb(i,k-1,j) )
420
421 grid%em_ph_2(i,k,j) = grid%em_ph_1(i,k,j)
422 grid%em_ph0(i,k,j) = grid%em_ph_1(i,k,j) + grid%em_phb(i,k,j)
423 ENDDO
424
425 ! interp u
426
427 DO K = 1, kte
428 p_level = grid%em_znu(k)*(p_surf - grid%p_top) + grid%p_top
429 grid%em_u_1(i,k,j) = interp_0( u, p_in, p_level, nl_in )
430 grid%em_u_2(i,k,j) = grid%em_u_1(i,k,j)
431 ENDDO
432
433 ENDDO
434 ENDDO
435
436 ! thermal perturbation to kick off convection
437
438 write(6,*) ' nxc, nyc for perturbation ',nxc,nyc
439 write(6,*) ' delt for perturbation ',tpbub
440
441 DO J = jts, min(jde-1,jte)
442 yrad = config_flags%dy*float(j-jde/2-1)/radbub
443 DO I = its, min(ide-1,ite)
444 xrad = float(i-1)/float(ide-ids)
445
446 DO K = 1, kte-1
447
448 ! put in preturbation theta (bubble) and recalc density. note,
449 ! the mass in the column is not changing, so when theta changes,
450 ! we recompute density and geopotential
451
452 zrad = 0.5*(grid%em_ph_1(i,k,j)+grid%em_ph_1(i,k+1,j) &
453 +grid%em_phb(i,k,j)+grid%em_phb(i,k+1,j))/g
454 zrad = (zrad-htbub)/radz
455 RAD=SQRT(yrad*yrad+zrad*zrad)
456 IF(RAD <= 1.) THEN
457 tp = tpbub*cos(rad*piov2)*cos(rad*piov2)*cos(xrad*2*pi+pi)
458 grid%em_t_1(i,k,j)=grid%em_t_1(i,k,j)+tp
459 grid%em_t_2(i,k,j)=grid%em_t_1(i,k,j)
460 qvf = 1. + rvovrd*grid%moist(i,k,j,P_QV)
461 grid%em_alt(i,k,j) = (r_d/p1000mb)*(grid%em_t_1(i,k,j)+t0)*qvf* &
462 (((grid%em_p(i,k,j)+grid%em_pb(i,k,j))/p1000mb)**cvpm)
463 grid%em_al(i,k,j) = grid%em_alt(i,k,j) - grid%em_alb(i,k,j)
464 ENDIF
465 ENDDO
466
467 ! rebalance hydrostatically
468
469 DO k = 2,kte
470 grid%em_ph_1(i,k,j) = grid%em_ph_1(i,k-1,j) - (1./grid%em_rdnw(k-1))*( &
471 (grid%em_mub(i,j)+grid%em_mu_1(i,j))*grid%em_al(i,k-1,j)+ &
472 grid%em_mu_1(i,j)*grid%em_alb(i,k-1,j) )
473
474 grid%em_ph_2(i,k,j) = grid%em_ph_1(i,k,j)
475 grid%em_ph0(i,k,j) = grid%em_ph_1(i,k,j) + grid%em_phb(i,k,j)
476 ENDDO
477
478 ENDDO
479 ENDDO
480
481 !#endif
482
483 write(6,*) ' grid%em_mu_1 from comp ', grid%em_mu_1(1,1)
484 write(6,*) ' pert state sounding from comp, grid%em_ph_1, pp, grid%em_al, grid%em_t_1, qv '
485 do k=1,kde-1
486 write(6,'(i3,1x,5(1x,1pe10.3))') k, grid%em_ph_1(1,k,1),grid%em_p(1,k,1), grid%em_al(1,k,1), &
487 grid%em_t_1(1,k,1), grid%moist(1,k,1,P_QV)
488 enddo
489
490 write(6,*) ' grid%em_mu_1 from comp ', grid%em_mu_1(1,1)
491 write(6,*) ' full state sounding from comp, ph, grid%em_p, grid%em_al, grid%em_t_1, qv '
492 write(6,*) ' at j = 1 '
493 do k=1,kde-1
494 write(6,'(i3,1x,5(1x,1pe10.3))') k, grid%em_ph_1(1,k,1)+grid%em_phb(1,k,1), &
495 grid%em_p(1,k,1)+grid%em_pb(1,k,1), grid%em_al(1,k,1)+grid%em_alb(1,k,1), &
496 grid%em_t_1(1,k,1)+t0, grid%moist(1,k,1,P_QV)
497 enddo
498
499
500 write(6,*) ' full state sounding from comp, ph, grid%em_p, grid%em_al, grid%em_t_1, qv '
501 write(6,*) ' at j = jde/2 '
502 do k=1,kde-1
503 write(6,'(i3,1x,5(1x,1pe10.3))') k, grid%em_ph_1(1,k,jde/2)+grid%em_phb(1,k,jde/2), &
504 grid%em_p(1,k,jde/2)+grid%em_pb(1,k,jde/2), grid%em_al(1,k,jde/2)+grid%em_alb(1,k,jde/2), &
505 grid%em_t_1(1,k,jde/2)+t0, grid%moist(1,k,jde/2,P_QV)
506 enddo
507
508 write(6,*) ' full state sounding from comp, ph, grid%em_p, grid%em_al, grid%em_t_1, qv '
509 write(6,*) ' at j = jde-1 '
510 do k=1,kde-1
511 write(6,'(i3,1x,5(1x,1pe10.3))') k, grid%em_ph_1(1,k,jde-1)+grid%em_phb(1,k,jde-1), &
512 grid%em_p(1,k,jde-1)+grid%em_pb(1,k,jde-1), grid%em_al(1,k,jde-1)+grid%em_alb(1,k,jde-1), &
513 grid%em_t_1(1,k,jde-1)+t0, grid%moist(1,k,jde-1,P_QV)
514 enddo
515
516 ! set v
517
518 DO J = jts, jte
519 DO I = its, min(ide-1,ite)
520
521 DO K = 1, kte
522 grid%em_v_1(i,k,j) = 0.
523 grid%em_v_2(i,k,j) = grid%em_v_1(i,k,j)
524 ENDDO
525
526 ENDDO
527 ENDDO
528
529 ! fill out last i row for u
530
531 DO J = jts, min(jde-1,jte)
532 DO I = ite, ite
533
534 DO K = 1, kte
535 grid%em_u_1(i,k,j) = grid%em_u_1(its,k,j)
536 grid%em_u_2(i,k,j) = grid%em_u_2(its,k,j)
537 ENDDO
538
539 ENDDO
540 ENDDO
541
542 ! set w
543
544 DO J = jts, min(jde-1,jte)
545 DO K = kts, kte
546 DO I = its, min(ide-1,ite)
547 grid%em_w_1(i,k,j) = 0.
548 grid%em_w_2(i,k,j) = 0.
549 ENDDO
550 ENDDO
551 ENDDO
552
553 ! set a few more things
554
555 DO J = jts, min(jde-1,jte)
556 DO K = kts, kte-1
557 DO I = its, min(ide-1,ite)
558 grid%h_diabatic(i,k,j) = 0.
559 ENDDO
560 ENDDO
561 ENDDO
562
563 DO k=1,kte-1
564 grid%em_t_base(k) = grid%em_t_1(1,k,1)
565 grid%qv_base(k) = grid%moist(1,k,1,P_QV)
566 grid%u_base(k) = grid%em_u_1(1,k,1)
567 grid%v_base(k) = grid%em_v_1(1,k,1)
568 ENDDO
569
570 DO J = jts, min(jde-1,jte)
571 DO I = its, min(ide-1,ite)
572 thtmp = grid%em_t_2(i,1,j)+t0
573 ptmp = grid%em_p(i,1,j)+grid%em_pb(i,1,j)
574 temp(1) = thtmp * (ptmp/p1000mb)**rcp
575 thtmp = grid%em_t_2(i,2,j)+t0
576 ptmp = grid%em_p(i,2,j)+grid%em_pb(i,2,j)
577 temp(2) = thtmp * (ptmp/p1000mb)**rcp
578 thtmp = grid%em_t_2(i,3,j)+t0
579 ptmp = grid%em_p(i,3,j)+grid%em_pb(i,3,j)
580 temp(3) = thtmp * (ptmp/p1000mb)**rcp
581
582 grid%tsk(I,J)=grid%cf1*temp(1)+grid%cf2*temp(2)+grid%cf3*temp(3)
583 if (i .eq. 1) print*,'sfctem',j,temp(1),temp(2),temp(3),grid%tsk(I,J)
584 grid%tmn(I,J)=grid%tsk(I,J)-0.5
585 ENDDO
586 ENDDO
587
588 RETURN
589
590 END SUBROUTINE init_domain_rk
591
592 !---------------------------------------------------------------------
593
594 SUBROUTINE init_module_initialize
595 END SUBROUTINE init_module_initialize
596
597 !---------------------------------------------------------------------
598 #if 0
599 ! TEST DRIVER FOR "read_input_jet" and "get_sounding"
600 implicit none
601 integer, parameter :: nz_jet=64, ny_jet=80
602 real, dimension(nz_jet,ny_jet) :: u_jet, rho_jet, &
603 th_jet, z_jet
604
605 real, dimension(nz_jet,ny_jet) :: zk,p,p_dry,theta,rho,u,v,qv
606 logical :: dry, debug
607 integer :: j, nl
608
609 call read_input_jet( u_jet, rho_jet, th_jet, z_jet, nz_jet, ny_jet )
610
611 call opngks
612 call parray( u_jet, nz_jet, ny_jet)
613 call parray( rho_jet, nz_jet, ny_jet)
614 call parray( th_jet, nz_jet, ny_jet)
615 ! call clsgks
616
617 ! set up initial jet
618
619 debug = .true.
620 dry = .true.
621 do j=1,ny_jet
622
623 call get_sounding( zk(:,j),p(:,j),p_dry(:,j),theta(:,j), &
624 rho(:,j),u(:,j), v(:,j), qv(:,j), &
625 dry, nz_jet, nl, u_jet, rho_jet, th_jet, &
626 z_jet, nz_jet, ny_jet, j, debug )
627 debug = .false.
628
629 enddo
630
631 write(6,*) ' lowest level p, th, and rho, highest level p '
632
633 do j=1,ny_jet
634 write(6,*) j, p(1,j),theta(1,j),rho(1,j), p(nz_jet,j)
635 ! write(6,*) j, p(1,j),theta(1,j)-th_jet(1,j),rho(1,j)-rho_jet(1,j)
636 enddo
637
638 call parray( p, nz_jet, ny_jet)
639 call parray( p_dry, nz_jet, ny_jet)
640 call parray( theta, nz_jet, ny_jet)
641
642 call clsgks
643
644 end
645
646 !---------------------------------
647
648 subroutine parray(a,m,n)
649 dimension a(m,n)
650 dimension b(n,m)
651
652 do i=1,m
653 do j=1,n
654 b(j,i) = a(i,j)
655 enddo
656 enddo
657
658 write(6,'('' dimensions m,n '',2i6)')m,n
659 call set(.05,.95,.05,.95,0.,1.,0.,1.,1)
660 call perim(4,5,4,5)
661 call setusv('LW',2000)
662 ! CALL CONREC(a,m,m,n,cmax,cmin,cinc,-1,-638,-922)
663 CALL CONREC(b,n,n,m,0.,0.,0.,-1,-638,-922)
664 call frame
665 return
666 end
667
668 ! END TEST DRIVER FOR "read_input_jet" and "get_sounding"
669 #endif
670
671 !------------------------------------------------------------------
672
673 subroutine get_sounding( zk, p, p_dry, theta, rho, &
674 u, v, qv, dry, nl_max, nl_in, &
675 u_jet, rho_jet, th_jet, z_jet, &
676 nz_jet, ny_jet, j_point, debug )
677 implicit none
678
679 integer nl_max, nl_in
680 real zk(nl_max), p(nl_max), theta(nl_max), rho(nl_max), &
681 u(nl_max), v(nl_max), qv(nl_max), p_dry(nl_max)
682 logical dry
683
684 integer nz_jet, ny_jet, j_point
685 real, dimension(nz_jet, ny_jet) :: u_jet, rho_jet, th_jet, z_jet
686
687 integer n
688 parameter(n=1000)
689 logical debug
690
691 ! input sounding data
692
693 real p_surf, th_surf, qv_surf
694 real pi_surf, pi(n)
695 real h_input(n), th_input(n), qv_input(n), u_input(n), v_input(n)
696
697 ! diagnostics
698
699 real rho_surf, p_input(n), rho_input(n)
700 real pm_input(n) ! this are for full moist sounding
701
702 ! local data
703
704 real p1000mb,cv,cp,r,cvpm,g
705 parameter (p1000mb = 1.e+05, r = 287, cp = 1003., cv = cp-r, cvpm = -cv/cp, g=9.81 )
706 integer k, it, nl
707 real qvf, qvf1, dz
708
709 ! first, read the sounding
710
711 ! call read_sounding( p_surf, th_surf, qv_surf, &
712 ! h_input, th_input, qv_input, u_input, v_input,n, nl, debug )
713
714 call calc_jet_sounding( p_surf, th_surf, qv_surf, &
715 h_input, th_input, qv_input, u_input, v_input, &
716 n, nl, debug, u_jet, rho_jet, th_jet, z_jet, j_point, &
717 nz_jet, ny_jet, dry )
718
719 nl = nz_jet
720
721 if(dry) then
722 do k=1,nl
723 qv_input(k) = 0.
724 enddo
725 endif
726
727 if(debug) write(6,*) ' number of input levels = ',nl
728
729 nl_in = nl
730 if(nl_in .gt. nl_max ) then
731 write(6,*) ' too many levels for input arrays ',nl_in,nl_max
732 call wrf_error_fatal ( ' too many levels for input arrays ' )
733 end if
734
735 ! compute diagnostics,
736 ! first, convert qv(g/kg) to qv(g/g)
737 !
738 ! do k=1,nl
739 ! qv_input(k) = 0.001*qv_input(k)
740 ! enddo
741 ! p_surf = 100.*p_surf ! convert to pascals
742
743 qvf = 1. + rvovrd*qv_input(1)
744 rho_surf = 1./((r/p1000mb)*th_surf*qvf*((p_surf/p1000mb)**cvpm))
745 pi_surf = (p_surf/p1000mb)**(r/cp)
746
747 if(debug) then
748 write(6,*) ' surface density is ',rho_surf
749 write(6,*) ' surface pi is ',pi_surf
750 end if
751
752
753 ! integrate moist sounding hydrostatically, starting from the
754 ! specified surface pressure
755 ! -> first, integrate from surface to lowest level
756
757 qvf = 1. + rvovrd*qv_input(1)
758 qvf1 = 1. + qv_input(1)
759 rho_input(1) = rho_surf
760 dz = h_input(1)
761 do it=1,10
762 pm_input(1) = p_surf &
763 - 0.5*dz*(rho_surf+rho_input(1))*g*qvf1
764 rho_input(1) = 1./((r/p1000mb)*th_input(1)*qvf*((pm_input(1)/p1000mb)**cvpm))
765 enddo
766
767 ! integrate up the column
768
769 do k=2,nl
770 rho_input(k) = rho_input(k-1)
771 dz = h_input(k)-h_input(k-1)
772 qvf1 = 0.5*(2.+(qv_input(k-1)+qv_input(k)))
773 qvf = 1. + rvovrd*qv_input(k) ! qv is in g/kg here
774
775 do it=1,10
776 pm_input(k) = pm_input(k-1) &
777 - 0.5*dz*(rho_input(k)+rho_input(k-1))*g*qvf1
778 rho_input(k) = 1./((r/p1000mb)*th_input(k)*qvf*((pm_input(k)/p1000mb)**cvpm))
779 enddo
780 enddo
781
782 ! we have the moist sounding
783
784 ! next, compute the dry sounding using p at the highest level from the
785 ! moist sounding and integrating down.
786
787 p_input(nl) = pm_input(nl)
788
789 do k=nl-1,1,-1
790 dz = h_input(k+1)-h_input(k)
791 p_input(k) = p_input(k+1) + 0.5*dz*(rho_input(k)+rho_input(k+1))*g
792 enddo
793
794
795 do k=1,nl
796
797 zk(k) = h_input(k)
798 p(k) = pm_input(k)
799 p_dry(k) = p_input(k)
800 theta(k) = th_input(k)
801 rho(k) = rho_input(k)
802 u(k) = u_input(k)
803 v(k) = v_input(k)
804 qv(k) = qv_input(k)
805
806 enddo
807
808 if(debug) then
809 write(6,*) ' sounding '
810 write(6,*) ' k height(m) press (Pa) pd(Pa) theta (K) den(kg/m^3) u(m/s) v(m/s) qv(g/g) '
811 do k=1,nl
812 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)
813 enddo
814
815 end if
816
817 end subroutine get_sounding
818
819 !------------------------------------------------------------------
820
821 subroutine calc_jet_sounding( p_surf, th_surf, qv_surf, &
822 h, th, qv, u, v, n, nl, debug, &
823 u_jet, rho_jet, th_jet, z_jet, &
824 jp, nz_jet, ny_jet, dry )
825 implicit none
826 integer :: n, nl, jp, nz_jet, ny_jet
827
828 real, dimension(nz_jet, ny_jet) :: u_jet, rho_jet, th_jet, z_jet
829 real, dimension(n) :: h,th,qv,u,v
830 real :: p_surf, th_surf, qv_surf
831 logical :: debug, dry
832
833 real, dimension(1:nz_jet) :: rho, rel_hum, p
834 integer :: k
835
836 ! some local stuff
837
838 real :: tmppi, es, qvs, temperature
839 real, parameter :: p1000mb=1.e+05, rcp=287./1004.5, svpt0=273.15, &
840 svp3 = 29.65, ep_2=287./461.6, r_d = 287., &
841 cpovcv = 1004./(1004.-287.), &
842 svp1 = 0.6112, svp2 = 17.67
843
844 ! get sounding from column jp
845
846 do k=1,nz_jet
847 h(k) = z_jet(k,jp)
848 th(k) = th_jet(k,jp)
849 qv(k) = 0.
850 rho(k) = rho_jet(k,jp)
851 u(k) = u_jet(k,jp)
852 v(k) = 0.
853 enddo
854
855 if (.not.dry) then
856 DO k=1,nz_jet
857 if(h(k) .gt. 8000.) then
858 rel_hum(k)=0.1
859 else
860 rel_hum(k)=(1.-0.90*(h(k)/8000.)**1.25)
861 end if
862 rel_hum(k) = min(0.7,rel_hum(k))
863 ENDDO
864 else
865 do k=1,nz_jet
866 rel_hum(k) = 0.
867 enddo
868 endif
869
870 ! next, compute pressure
871
872 do k=1,nz_jet
873 p(k) = p1000mb*(R_d*rho(k)*th(k)/p1000mb)**cpovcv
874 enddo
875
876 ! here we adjust for fixed moisture profile
877
878 IF (.not.dry) THEN
879
880 ! here we assume the input theta is th_v, so we reset theta accordingly
881
882 DO k=1,nz_jet
883 tmppi=(p(k)/p1000mb)**rcp
884 temperature = tmppi*th(k)
885 if (temperature .gt. svpt0) then
886 es = 1000.*svp1*exp(svp2*(temperature-svpt0)/(temperature-svp3))
887 qvs = ep_2*es/(p(k)-es)
888 else
889 es = 1000.*svp1*exp( 21.8745584*(temperature-273.16)/(temperature-7.66) )
890 qvs = ep_2*es/(p(k)-es)
891 endif
892 qv(k) = rel_hum(k)*qvs
893 th(k) = th(k)/(1.+.61*qv(k))
894 ENDDO
895
896 ENDIF
897
898 ! finally, set the surface data. We'll just do a simple extrapolation
899
900 p_surf = 1.5*p(1) - 0.5*p(2)
901 th_surf = 1.5*th(1) - 0.5*th(2)
902 qv_surf = 1.5*qv(1) - 0.5*qv(2)
903
904 end subroutine calc_jet_sounding
905
906 !---------------------------------------------------------------------
907
908 SUBROUTINE read_input_jet( u, r, t, zk, nz, ny )
909 implicit none
910
911 integer, intent(in) :: nz,ny
912 real, dimension(nz,ny), intent(out) :: u,r,t,zk
913 integer :: ny_in, nz_in, j,k
914 real, dimension(ny,nz) :: field_in
915
916 ! this code assumes it is called on processor 0 only
917
918 OPEN(unit=10, file='input_jet', form='unformatted', status='old' )
919 REWIND(10)
920 read(10) ny_in,nz_in
921 if((ny_in /= ny ) .or. (nz_in /= nz)) then
922 write(0,*) ' error in input jet dimensions '
923 write(0,*) ' ny, ny_input, nz, nz_input ', ny, ny_in, nz,nz_in
924 write(0,*) ' error exit '
925 call wrf_error_fatal ( ' error in input jet dimensions ' )
926 end if
927 read(10) field_in
928 do j=1,ny
929 do k=1,nz
930 u(k,j) = field_in(j,k)
931 enddo
932 enddo
933 read(10) field_in
934 do j=1,ny
935 do k=1,nz
936 t(k,j) = field_in(j,k)
937 enddo
938 enddo
939
940 read(10) field_in
941 do j=1,ny
942 do k=1,nz
943 r(k,j) = field_in(j,k)
944 enddo
945 enddo
946
947 do j=1,ny
948 do k=1,nz
949 zk(k,j) = 125. + 250.*float(k-1)
950 enddo
951 enddo
952
953 end subroutine read_input_jet
954
955 END MODULE module_initialize_ideal