module_initialize_squall2d_y.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
121 #ifdef DM_PARALLEL
122 # include <em_data_calls.inc>
123 #endif
124
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
167 stretch_grid = .true.
168 delt = 3.
169 ! z_scale = .50
170 z_scale = .40
171 pi = 2.*asin(1.0)
172 write(6,*) ' pi is ',pi
173 nxc = (ide-ids)/2
174 nyc = (jde-jds)/2
175
176 CALL model_to_grid_config_rec ( grid%id , model_config_rec , config_flags )
177
178 ! here we check to see if the boundary conditions are set properly
179
180 CALL boundary_condition_check( config_flags, bdyzone, error, grid%id )
181
182 moisture_init = .true.
183
184 grid%itimestep=0
185
186 #ifdef DM_PARALLEL
187 CALL wrf_dm_bcast_bytes( icm , IWORDSIZE )
188 CALL wrf_dm_bcast_bytes( jcm , IWORDSIZE )
189 #endif
190
191 CALL nl_set_mminlu(1, ' ')
192 CALL nl_set_iswater(1,0)
193 CALL nl_set_cen_lat(1,40.)
194 CALL nl_set_cen_lon(1,-105.)
195 CALL nl_set_truelat1(1,0.)
196 CALL nl_set_truelat2(1,0.)
197 CALL nl_set_moad_cen_lat (1,0.)
198 CALL nl_set_stand_lon (1,0.)
199 CALL nl_set_map_proj(1,0)
200
201
202 ! here we initialize data we currently is not initialized
203 ! in the input data
204
205 DO j = jts, jte
206 DO i = its, ite
207 grid%msftx(i,j) = 1.
208 grid%msfty(i,j) = 1.
209 grid%msfux(i,j) = 1.
210 grid%msfuy(i,j) = 1.
211 grid%msfvx(i,j) = 1.
212 grid%msfvx_inv(i,j)= 1.
213 grid%msfvy(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) = 0.
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,*) ' getting dry sounding for base state '
272 dry_sounding = .true.
273 CALL get_sounding( zk, p_in, pd_in, theta, rho, u, v, qv, dry_sounding, nl_max, nl_in )
274
275 write(6,*) ' returned from reading sounding, nl_in is ',nl_in
276
277 ! find ptop for the desired ztop (ztop is input from the namelist),
278 ! and find surface pressure
279
280 grid%p_top = interp_0( p_in, zk, config_flags%ztop, nl_in )
281
282 DO j=jts,jte
283 DO i=its,ite ! flat surface
284 grid%em_phb(i,1,j) = 0.
285 grid%em_php(i,1,j) = 0.
286 grid%em_ph0(i,1,j) = 0.
287 grid%ht(i,j) = 0.
288 ENDDO
289 ENDDO
290
291 DO J = jts, jte
292 DO I = its, ite
293
294 p_surf = interp_0( p_in, zk, grid%em_phb(i,1,j)/g, nl_in )
295 grid%em_mub(i,j) = p_surf-grid%p_top
296
297 ! this is dry hydrostatic sounding (base state), so given grid%em_p (coordinate),
298 ! interp theta (from interp) and compute 1/rho from eqn. of state
299
300 DO K = 1, kte-1
301 p_level = grid%em_znu(k)*(p_surf - grid%p_top) + grid%p_top
302 grid%em_pb(i,k,j) = p_level
303 grid%em_t_init(i,k,j) = interp_0( theta, p_in, p_level, nl_in ) - t0
304 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
305 ENDDO
306
307 ! calc hydrostatic balance (alternatively we could interp the geopotential from the
308 ! sounding, but this assures that the base state is in exact hydrostatic balance with
309 ! respect to the model eqns.
310
311 DO k = 2,kte
312 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)
313 ENDDO
314
315 ENDDO
316 ENDDO
317
318 write(6,*) ' ptop is ',grid%p_top
319 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
320
321 ! calculate full state for each column - this includes moisture.
322
323 write(6,*) ' getting moist sounding for full state '
324 dry_sounding = .false.
325 CALL get_sounding( zk, p_in, pd_in, theta, rho, u, v, qv, dry_sounding, nl_max, nl_in )
326
327 DO J = jts, min(jde-1,jte)
328 DO I = its, min(ide-1,ite)
329
330 ! At this point grid%p_top is already set. find the DRY mass in the column
331 ! by interpolating the DRY pressure.
332
333 pd_surf = interp_0( pd_in, zk, grid%em_phb(i,1,j)/g, nl_in )
334
335 ! compute the perturbation mass and the full mass
336
337 grid%em_mu_1(i,j) = pd_surf-grid%p_top - grid%em_mub(i,j)
338 grid%em_mu_2(i,j) = grid%em_mu_1(i,j)
339 grid%em_mu0(i,j) = grid%em_mu_1(i,j) + grid%em_mub(i,j)
340
341 ! given the dry pressure and coordinate system, interp the potential
342 ! temperature and qv
343
344 do k=1,kde-1
345
346 p_level = grid%em_znu(k)*(pd_surf - grid%p_top) + grid%p_top
347
348 moist(i,k,j,P_QV) = interp_0( qv, pd_in, p_level, nl_in )
349 grid%em_t_1(i,k,j) = interp_0( theta, pd_in, p_level, nl_in ) - t0
350 grid%em_t_2(i,k,j) = grid%em_t_1(i,k,j)
351
352
353 enddo
354
355 ! integrate the hydrostatic equation (from the RHS of the bigstep
356 ! vertical momentum equation) down from the top to get grid%em_p.
357 ! first from the top of the model to the top pressure
358
359 k = kte-1 ! top level
360
361 qvf1 = 0.5*(moist(i,k,j,P_QV)+moist(i,k,j,P_QV))
362 qvf2 = 1./(1.+qvf1)
363 qvf1 = qvf1*qvf2
364
365 ! grid%em_p(i,k,j) = - 0.5*grid%em_mu_1(i,j)/grid%em_rdnw(k)
366 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
367 qvf = 1. + rvovrd*moist(i,k,j,P_QV)
368 grid%em_alt(i,k,j) = (r_d/p1000mb)*(grid%em_t_1(i,k,j)+t0)*qvf* &
369 (((grid%em_p(i,k,j)+grid%em_pb(i,k,j))/p1000mb)**cvpm)
370 grid%em_al(i,k,j) = grid%em_alt(i,k,j) - grid%em_alb(i,k,j)
371
372 ! down the column
373
374 do k=kte-2,1,-1
375 qvf1 = 0.5*(moist(i,k,j,P_QV)+moist(i,k+1,j,P_QV))
376 qvf2 = 1./(1.+qvf1)
377 qvf1 = qvf1*qvf2
378 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)
379 qvf = 1. + rvovrd*moist(i,k,j,P_QV)
380 grid%em_alt(i,k,j) = (r_d/p1000mb)*(grid%em_t_1(i,k,j)+t0)*qvf* &
381 (((grid%em_p(i,k,j)+grid%em_pb(i,k,j))/p1000mb)**cvpm)
382 grid%em_al(i,k,j) = grid%em_alt(i,k,j) - grid%em_alb(i,k,j)
383 enddo
384
385 ! this is the hydrostatic equation used in the model after the
386 ! small timesteps. In the model, grid%em_al (inverse density)
387 ! is computed from the geopotential.
388
389
390 grid%em_ph_1(i,1,j) = 0.
391 DO k = 2,kte
392 grid%em_ph_1(i,k,j) = grid%em_ph_1(i,k-1,j) - (1./grid%em_rdnw(k-1))*( &
393 (grid%em_mub(i,j)+grid%em_mu_1(i,j))*grid%em_al(i,k-1,j)+ &
394 grid%em_mu_1(i,j)*grid%em_alb(i,k-1,j) )
395
396 grid%em_ph_2(i,k,j) = grid%em_ph_1(i,k,j)
397 grid%em_ph0(i,k,j) = grid%em_ph_1(i,k,j) + grid%em_phb(i,k,j)
398 ENDDO
399
400 if((i==2) .and. (j==2)) then
401 write(6,*) ' grid%em_ph_1 calc ',grid%em_ph_1(2,1,2),grid%em_ph_1(2,2,2),&
402 grid%em_mu_1(2,2)+grid%em_mub(2,2),grid%em_mu_1(2,2), &
403 grid%em_alb(2,1,2),grid%em_al(1,2,1),grid%em_rdnw(1)
404 endif
405
406 ENDDO
407 ENDDO
408
409 !#if 0
410
411 ! thermal perturbation to kick off convection
412
413 write(6,*) ' nxc, nyc for perturbation ',nxc,nyc
414 write(6,*) ' delt for perturbation ',delt
415
416 DO J = jts, min(jde-1,jte)
417 yrad = config_flags%dy*float(j-nyc)/4000.
418 ! yrad = 0.
419 DO I = its, min(ide-1,ite)
420 ! xrad = config_flags%dx*float(i-nxc)/4000.
421 xrad = 0.
422 DO K = 1, kte-1
423
424 ! put in preturbation theta (bubble) and recalc density. note,
425 ! the mass in the column is not changing, so when theta changes,
426 ! we recompute density and geopotential
427
428 zrad = 0.5*(grid%em_ph_1(i,k,j)+grid%em_ph_1(i,k+1,j) &
429 +grid%em_phb(i,k,j)+grid%em_phb(i,k+1,j))/g
430 zrad = (zrad-1500.)/1500.
431 RAD=SQRT(xrad*xrad+yrad*yrad+zrad*zrad)
432 IF(RAD <= 1.) THEN
433 grid%em_t_1(i,k,j)=grid%em_t_1(i,k,j)+delt*COS(.5*PI*RAD)**2
434 grid%em_t_2(i,k,j)=grid%em_t_1(i,k,j)
435 qvf = 1. + rvovrd*moist(i,k,j,P_QV)
436 grid%em_alt(i,k,j) = (r_d/p1000mb)*(grid%em_t_1(i,k,j)+t0)*qvf* &
437 (((grid%em_p(i,k,j)+grid%em_pb(i,k,j))/p1000mb)**cvpm)
438 grid%em_al(i,k,j) = grid%em_alt(i,k,j) - grid%em_alb(i,k,j)
439 ENDIF
440 ENDDO
441
442 ! rebalance hydrostatically
443
444 DO k = 2,kte
445 grid%em_ph_1(i,k,j) = grid%em_ph_1(i,k-1,j) - (1./grid%em_rdnw(k-1))*( &
446 (grid%em_mub(i,j)+grid%em_mu_1(i,j))*grid%em_al(i,k-1,j)+ &
447 grid%em_mu_1(i,j)*grid%em_alb(i,k-1,j) )
448
449 grid%em_ph_2(i,k,j) = grid%em_ph_1(i,k,j)
450 grid%em_ph0(i,k,j) = grid%em_ph_1(i,k,j) + grid%em_phb(i,k,j)
451 ENDDO
452
453 ENDDO
454 ENDDO
455
456 !#endif
457
458 write(6,*) ' grid%em_mu_1 from comp ', grid%em_mu_1(1,1)
459 write(6,*) ' full state sounding from comp, ph, grid%em_p, grid%em_al, grid%em_t_1, qv '
460 do k=1,kde-1
461 write(6,'(i3,1x,5(1x,1pe10.3))') k, grid%em_ph_1(1,k,1)+grid%em_phb(1,k,1), &
462 grid%em_p(1,k,1)+grid%em_pb(1,k,1), grid%em_alt(1,k,1), &
463 grid%em_t_1(1,k,1)+t0, moist(1,k,1,P_QV)
464 enddo
465
466 write(6,*) ' pert state sounding from comp, grid%em_ph_1, pp, alp, grid%em_t_1, qv '
467 do k=1,kde-1
468 write(6,'(i3,1x,5(1x,1pe10.3))') k, grid%em_ph_1(1,k,1), &
469 grid%em_p(1,k,1), grid%em_al(1,k,1), &
470 grid%em_t_1(1,k,1), moist(1,k,1,P_QV)
471 enddo
472
473 ! interp v
474
475 DO J = jts, jte
476 DO I = its, min(ide-1,ite)
477
478 IF (j == jds) THEN
479 z_at_v = grid%em_phb(i,1,j)/g
480 ELSE IF (j == jde) THEN
481 z_at_v = grid%em_phb(i,1,j-1)/g
482 ELSE
483 z_at_v = 0.5*(grid%em_phb(i,1,j)+grid%em_phb(i,1,j-1))/g
484 END IF
485
486 p_surf = interp_0( p_in, zk, z_at_v, nl_in )
487
488 DO K = 1, kte
489 p_level = grid%em_znu(k)*(p_surf - grid%p_top) + grid%p_top
490 grid%em_v_1(i,k,j) = interp_0( v, p_in, p_level, nl_in )
491 grid%em_v_2(i,k,j) = grid%em_v_1(i,k,j)
492 ENDDO
493
494 ENDDO
495 ENDDO
496
497 ! interp u
498
499 DO J = jts, min(jde-1,jte)
500 DO I = its, ite
501
502 IF (i == ids) THEN
503 z_at_u = grid%em_phb(i,1,j)/g
504 ELSE IF (i == ide) THEN
505 z_at_u = grid%em_phb(i-1,1,j)/g
506 ELSE
507 z_at_u = 0.5*(grid%em_phb(i,1,j)+grid%em_phb(i-1,1,j))/g
508 END IF
509
510 p_surf = interp_0( p_in, zk, z_at_u, nl_in )
511
512 DO K = 1, kte
513 p_level = grid%em_znu(k)*(p_surf - grid%p_top) + grid%p_top
514 grid%em_u_1(i,k,j) = interp_0( u, p_in, p_level, nl_in )
515 grid%em_u_2(i,k,j) = grid%em_u_1(i,k,j)
516 ENDDO
517
518 ENDDO
519 ENDDO
520
521 ! set w
522
523 DO J = jts, min(jde-1,jte)
524 DO K = kts, kte
525 DO I = its, min(ide-1,ite)
526 grid%em_w_1(i,k,j) = 0.
527 grid%em_w_2(i,k,j) = 0.
528 ENDDO
529 ENDDO
530 ENDDO
531
532 ! set a few more things
533
534 DO J = jts, min(jde-1,jte)
535 DO K = kts, kte-1
536 DO I = its, min(ide-1,ite)
537 grid%h_diabatic(i,k,j) = 0.
538 ENDDO
539 ENDDO
540 ENDDO
541
542 DO k=1,kte-1
543 grid%em_t_base(k) = grid%em_t_1(1,k,1)
544 grid%qv_base(k) = moist(1,k,1,P_QV)
545 grid%u_base(k) = grid%em_u_1(1,k,1)
546 grid%v_base(k) = grid%em_v_1(1,k,1)
547 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
548 ENDDO
549
550 DO J = jts, min(jde-1,jte)
551 DO I = its, min(ide-1,ite)
552 thtmp = grid%em_t_2(i,1,j)+t0
553 ptmp = grid%em_p(i,1,j)+grid%em_pb(i,1,j)
554 temp(1) = thtmp * (ptmp/p1000mb)**rcp
555 thtmp = grid%em_t_2(i,2,j)+t0
556 ptmp = grid%em_p(i,2,j)+grid%em_pb(i,2,j)
557 temp(2) = thtmp * (ptmp/p1000mb)**rcp
558 thtmp = grid%em_t_2(i,3,j)+t0
559 ptmp = grid%em_p(i,3,j)+grid%em_pb(i,3,j)
560 temp(3) = thtmp * (ptmp/p1000mb)**rcp
561
562 grid%tsk(I,J)=grid%cf1*temp(1)+grid%cf2*temp(2)+grid%cf3*temp(3)
563 grid%tmn(I,J)=grid%tsk(I,J)-0.5
564 ENDDO
565 ENDDO
566
567 RETURN
568
569 END SUBROUTINE init_domain_rk
570
571 SUBROUTINE init_module_initialize
572 END SUBROUTINE init_module_initialize
573
574 !---------------------------------------------------------------------
575
576 ! test driver for get_sounding
577 !
578 ! implicit none
579 ! integer n
580 ! parameter(n = 1000)
581 ! real zk(n),p(n),theta(n),rho(n),u(n),v(n),qv(n),pd(n)
582 ! logical dry
583 ! integer nl,k
584 !
585 ! dry = .false.
586 ! dry = .true.
587 ! call get_sounding( zk, p, pd, theta, rho, u, v, qv, dry, n, nl )
588 ! write(6,*) ' input levels ',nl
589 ! write(6,*) ' sounding '
590 ! write(6,*) ' k height(m) press (Pa) pd(Pa) theta (K) den(kg/m^3) u(m/s) v(m/s) qv(g/g) '
591 ! do k=1,nl
592 ! 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)
593 ! enddo
594 ! end
595 !
596 !---------------------------------------------------------------------------
597
598 subroutine get_sounding( zk, p, p_dry, theta, rho, &
599 u, v, qv, dry, nl_max, nl_in )
600 implicit none
601
602 integer nl_max, nl_in
603 real zk(nl_max), p(nl_max), theta(nl_max), rho(nl_max), &
604 u(nl_max), v(nl_max), qv(nl_max), p_dry(nl_max)
605 logical dry
606
607 integer n
608 parameter(n=1000)
609 logical debug
610 parameter( debug = .true.)
611
612 ! input sounding data
613
614 real p_surf, th_surf, qv_surf
615 real pi_surf, pi(n)
616 real h_input(n), th_input(n), qv_input(n), u_input(n), v_input(n)
617
618 ! diagnostics
619
620 real rho_surf, p_input(n), rho_input(n)
621 real pm_input(n) ! this are for full moist sounding
622
623 ! local data
624
625 real p1000mb,cv,cp,r,cvpm,g
626 parameter (p1000mb = 1.e+05, r = 287, cp = 1003., cv = cp-r, cvpm = -cv/cp, g=9.81 )
627 integer k, it, nl
628 real qvf, qvf1, dz
629
630 ! first, read the sounding
631
632 call read_sounding( p_surf, th_surf, qv_surf, &
633 h_input, th_input, qv_input, u_input, v_input,n, nl, debug )
634
635 if(dry) then
636 do k=1,nl
637 qv_input(k) = 0.
638 enddo
639 endif
640
641 if(debug) write(6,*) ' number of input levels = ',nl
642
643 nl_in = nl
644 if(nl_in .gt. nl_max ) then
645 write(6,*) ' too many levels for input arrays ',nl_in,nl_max
646 call wrf_error_fatal ( ' too many levels for input arrays ' )
647 end if
648
649 ! compute diagnostics,
650 ! first, convert qv(g/kg) to qv(g/g)
651
652 do k=1,nl
653 qv_input(k) = 0.001*qv_input(k)
654 enddo
655
656 p_surf = 100.*p_surf ! convert to pascals
657 qvf = 1. + rvovrd*qv_input(1)
658 rho_surf = 1./((r/p1000mb)*th_surf*qvf*((p_surf/p1000mb)**cvpm))
659 pi_surf = (p_surf/p1000mb)**(r/cp)
660
661 if(debug) then
662 write(6,*) ' surface density is ',rho_surf
663 write(6,*) ' surface pi is ',pi_surf
664 end if
665
666
667 ! integrate moist sounding hydrostatically, starting from the
668 ! specified surface pressure
669 ! -> first, integrate from surface to lowest level
670
671 qvf = 1. + rvovrd*qv_input(1)
672 qvf1 = 1. + qv_input(1)
673 rho_input(1) = rho_surf
674 dz = h_input(1)
675 do it=1,10
676 pm_input(1) = p_surf &
677 - 0.5*dz*(rho_surf+rho_input(1))*g*qvf1
678 rho_input(1) = 1./((r/p1000mb)*th_input(1)*qvf*((pm_input(1)/p1000mb)**cvpm))
679 enddo
680
681 ! integrate up the column
682
683 do k=2,nl
684 rho_input(k) = rho_input(k-1)
685 dz = h_input(k)-h_input(k-1)
686 qvf1 = 0.5*(2.+(qv_input(k-1)+qv_input(k)))
687 qvf = 1. + rvovrd*qv_input(k) ! qv is in g/kg here
688
689 do it=1,10
690 pm_input(k) = pm_input(k-1) &
691 - 0.5*dz*(rho_input(k)+rho_input(k-1))*g*qvf1
692 rho_input(k) = 1./((r/p1000mb)*th_input(k)*qvf*((pm_input(k)/p1000mb)**cvpm))
693 enddo
694 enddo
695
696 ! we have the moist sounding
697
698 ! next, compute the dry sounding using p at the highest level from the
699 ! moist sounding and integrating down.
700
701 p_input(nl) = pm_input(nl)
702
703 do k=nl-1,1,-1
704 dz = h_input(k+1)-h_input(k)
705 p_input(k) = p_input(k+1) + 0.5*dz*(rho_input(k)+rho_input(k+1))*g
706 enddo
707
708
709 do k=1,nl
710
711 zk(k) = h_input(k)
712 p(k) = pm_input(k)
713 p_dry(k) = p_input(k)
714 theta(k) = th_input(k)
715 rho(k) = rho_input(k)
716 u(k) = u_input(k)
717 v(k) = v_input(k)
718 qv(k) = qv_input(k)
719
720 enddo
721
722 if(debug) then
723 write(6,*) ' sounding '
724 write(6,*) ' k height(m) press (Pa) pd(Pa) theta (K) den(kg/m^3) u(m/s) v(m/s) qv(g/g) '
725 do k=1,nl
726 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)
727 enddo
728
729 end if
730
731 end subroutine get_sounding
732
733 !-------------------------------------------------------
734
735 subroutine read_sounding( ps,ts,qvs,h,th,qv,u,v,n,nl,debug )
736 implicit none
737 integer n,nl
738 real ps,ts,qvs,h(n),th(n),qv(n),u(n),v(n)
739 logical end_of_file
740 logical debug
741
742 integer k
743
744 open(unit=10,file='input_sounding',form='formatted',status='old')
745 rewind(10)
746 read(10,*) ps, ts, qvs
747 if(debug) then
748 write(6,*) ' input sounding surface parameters '
749 write(6,*) ' surface pressure (mb) ',ps
750 write(6,*) ' surface pot. temp (K) ',ts
751 write(6,*) ' surface mixing ratio (g/kg) ',qvs
752 end if
753
754 end_of_file = .false.
755 k = 0
756
757 do while (.not. end_of_file)
758
759 read(10,*,end=100) h(k+1), th(k+1), qv(k+1), u(k+1), v(k+1)
760 k = k+1
761 if(debug) write(6,'(1x,i3,5(1x,e10.3))') k, h(k), th(k), qv(k), u(k), v(k)
762 go to 110
763 100 end_of_file = .true.
764 110 continue
765 enddo
766
767 nl = k
768
769 close(unit=10,status = 'keep')
770
771 end subroutine read_sounding
772
773 END MODULE module_initialize_ideal