module_initialize_squall2d_x.F
References to this file elsewhere.
1 !IDEAL:MODEL_LAYER:INITIALIZATION
2 !
3
4 ! This MODULE holds the routines which are used to perform various initializations
5 ! for the individual domains.
6
7 ! This MODULE CONTAINS the following routines:
8
9 ! initialize_field_test - 1. Set different fields to different constant
10 ! values. This is only a test. If the correct
11 ! domain is not found (based upon the "id")
12 ! then a fatal error is issued.
13
14 !-----------------------------------------------------------------------
15
16 MODULE module_initialize_ideal
17
18 USE module_domain
19 USE module_io_domain
20 USE module_state_description
21 USE module_model_constants
22 USE module_bc
23 USE module_timing
24 USE module_configure
25 USE module_init_utilities
26 #ifdef DM_PARALLEL
27 USE module_dm
28 #endif
29
30
31 CONTAINS
32
33
34 !-------------------------------------------------------------------
35 ! this is a wrapper for the solver-specific init_domain routines.
36 ! Also dereferences the grid variables and passes them down as arguments.
37 ! This is crucial, since the lower level routines may do message passing
38 ! and this will get fouled up on machines that insist on passing down
39 ! copies of assumed-shape arrays (by passing down as arguments, the
40 ! data are treated as assumed-size -- ie. f77 -- arrays and the copying
41 ! business is avoided). Fie on the F90 designers. Fie and a pox.
42 ! NOTE: Modified to remove all but arrays of rank 4 or more from the
43 ! argument list. Arrays with rank>3 are still problematic due to the
44 ! above-noted fie- and pox-ities. TBH 20061129.
45
46 SUBROUTINE init_domain ( grid )
47
48 IMPLICIT NONE
49
50 ! Input data.
51 TYPE (domain), POINTER :: grid
52 ! Local data.
53 INTEGER :: dyn_opt
54 INTEGER :: idum1, idum2
55
56 CALL nl_get_dyn_opt( 1, dyn_opt )
57
58 CALL set_scalar_indices_from_config ( head_grid%id , idum1, idum2 )
59
60 IF ( dyn_opt .eq. 1 &
61 .or. dyn_opt .eq. 2 &
62 .or. dyn_opt .eq. 3 &
63 ) THEN
64 CALL init_domain_rk( grid &
65 !
66 #include <em_actual_new_args.inc>
67 !
68 )
69
70 ELSE
71 WRITE(0,*)' init_domain: unknown or unimplemented dyn_opt = ',dyn_opt
72 CALL wrf_error_fatal ( ' init_domain: unknown or unimplemented dyn_opt ' )
73 ENDIF
74
75 END SUBROUTINE init_domain
76
77 !-------------------------------------------------------------------
78
79 SUBROUTINE init_domain_rk ( grid &
80 !
81 # include <em_dummy_new_args.inc>
82 !
83 )
84 IMPLICIT NONE
85
86 ! Input data.
87 TYPE (domain), POINTER :: grid
88
89 # include <em_dummy_new_decl.inc>
90
91 TYPE (grid_config_rec_type) :: config_flags
92
93 ! Local data
94 INTEGER :: &
95 ids, ide, jds, jde, kds, kde, &
96 ims, ime, jms, jme, kms, kme, &
97 its, ite, jts, jte, kts, kte, &
98 i, j, k
99
100 ! Local data
101
102 INTEGER, PARAMETER :: nl_max = 1000
103 REAL, DIMENSION(nl_max) :: zk, p_in, theta, rho, u, v, qv, pd_in
104 INTEGER :: nl_in
105
106
107 INTEGER :: icm,jcm, ii, im1, jj, jm1, loop, error, fid, nxc, nyc
108 REAL :: u_mean,v_mean, f0, p_surf, p_level, qvf, z_at_v, z_at_u
109 REAL :: z_scale, xrad, yrad, zrad, rad, delt, cof1, cof2
110 ! REAL, EXTERNAL :: interp_0
111 REAL :: hm
112 REAL :: pi
113
114 ! stuff from original initialization that has been dropped from the Registry
115 REAL :: vnu, xnu, xnus, dinit0, cbh, p0_temp, t0_temp, zd, zt
116 REAL :: qvf1, qvf2, pd_surf
117 INTEGER :: it
118 real :: thtmp, ptmp, temp(3)
119
120 LOGICAL :: moisture_init
121 LOGICAL :: stretch_grid, dry_sounding
122
123
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/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 write(6,*) ' getting dry sounding for base state '
275 dry_sounding = .true.
276 CALL get_sounding( zk, p_in, pd_in, theta, rho, u, v, qv, dry_sounding, nl_max, nl_in )
277
278 write(6,*) ' returned from reading sounding, nl_in is ',nl_in
279
280 ! find ptop for the desired ztop (ztop is input from the namelist),
281 ! and find surface pressure
282
283 grid%p_top = interp_0( p_in, zk, config_flags%ztop, nl_in )
284
285 DO j=jts,jte
286 DO i=its,ite ! flat surface
287 grid%em_phb(i,1,j) = 0.
288 grid%em_php(i,1,j) = 0.
289 grid%em_ph0(i,1,j) = 0.
290 grid%ht(i,j) = 0.
291 ENDDO
292 ENDDO
293
294 DO J = jts, jte
295 DO I = its, ite
296
297 p_surf = interp_0( p_in, zk, grid%em_phb(i,1,j)/g, nl_in )
298 grid%em_mub(i,j) = p_surf-grid%p_top
299
300 ! this is dry hydrostatic sounding (base state), so given grid%em_p (coordinate),
301 ! interp theta (from interp) and compute 1/rho from eqn. of state
302
303 DO K = 1, kte-1
304 p_level = grid%em_znu(k)*(p_surf - grid%p_top) + grid%p_top
305 grid%em_pb(i,k,j) = p_level
306 grid%em_t_init(i,k,j) = interp_0( theta, p_in, p_level, nl_in ) - t0
307 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
308 ENDDO
309
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.
313
314 DO k = 2,kte
315 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)
316 ENDDO
317
318 ENDDO
319 ENDDO
320
321 write(6,*) ' ptop is ',grid%p_top
322 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
323
324 ! calculate full state for each column - this includes moisture.
325
326 write(6,*) ' getting moist sounding for full state '
327 dry_sounding = .false.
328 CALL get_sounding( zk, p_in, pd_in, theta, rho, u, v, qv, dry_sounding, nl_max, nl_in )
329
330 DO J = jts, min(jde-1,jte)
331 DO I = its, min(ide-1,ite)
332
333 ! At this point grid%p_top is already set. find the DRY mass in the column
334 ! by interpolating the DRY pressure.
335
336 pd_surf = interp_0( pd_in, zk, grid%em_phb(i,1,j)/g, nl_in )
337
338 ! compute the perturbation mass and the full mass
339
340 grid%em_mu_1(i,j) = pd_surf-grid%p_top - grid%em_mub(i,j)
341 grid%em_mu_2(i,j) = grid%em_mu_1(i,j)
342 grid%em_mu0(i,j) = grid%em_mu_1(i,j) + grid%em_mub(i,j)
343
344 ! given the dry pressure and coordinate system, interp the potential
345 ! temperature and qv
346
347 do k=1,kde-1
348
349 p_level = grid%em_znu(k)*(pd_surf - grid%p_top) + grid%p_top
350
351 moist(i,k,j,P_QV) = interp_0( qv, pd_in, p_level, nl_in )
352 grid%em_t_1(i,k,j) = interp_0( theta, pd_in, p_level, nl_in ) - t0
353 grid%em_t_2(i,k,j) = grid%em_t_1(i,k,j)
354
355
356 enddo
357
358 ! integrate the hydrostatic equation (from the RHS of the bigstep
359 ! vertical momentum equation) down from the top to get grid%em_p.
360 ! first from the top of the model to the top pressure
361
362 k = kte-1 ! top level
363
364 qvf1 = 0.5*(moist(i,k,j,P_QV)+moist(i,k,j,P_QV))
365 qvf2 = 1./(1.+qvf1)
366 qvf1 = qvf1*qvf2
367
368 ! grid%em_p(i,k,j) = - 0.5*grid%em_mu_1(i,j)/grid%em_rdnw(k)
369 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
370 qvf = 1. + rvovrd*moist(i,k,j,P_QV)
371 grid%em_alt(i,k,j) = (r_d/p1000mb)*(grid%em_t_1(i,k,j)+t0)*qvf* &
372 (((grid%em_p(i,k,j)+grid%em_pb(i,k,j))/p1000mb)**cvpm)
373 grid%em_al(i,k,j) = grid%em_alt(i,k,j) - grid%em_alb(i,k,j)
374
375 ! down the column
376
377 do k=kte-2,1,-1
378 qvf1 = 0.5*(moist(i,k,j,P_QV)+moist(i,k+1,j,P_QV))
379 qvf2 = 1./(1.+qvf1)
380 qvf1 = qvf1*qvf2
381 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)
382 qvf = 1. + rvovrd*moist(i,k,j,P_QV)
383 grid%em_alt(i,k,j) = (r_d/p1000mb)*(grid%em_t_1(i,k,j)+t0)*qvf* &
384 (((grid%em_p(i,k,j)+grid%em_pb(i,k,j))/p1000mb)**cvpm)
385 grid%em_al(i,k,j) = grid%em_alt(i,k,j) - grid%em_alb(i,k,j)
386 enddo
387
388 ! this is the hydrostatic equation used in the model after the
389 ! small timesteps. In the model, grid%em_al (inverse density)
390 ! is computed from the geopotential.
391
392
393 grid%em_ph_1(i,1,j) = 0.
394 DO k = 2,kte
395 grid%em_ph_1(i,k,j) = grid%em_ph_1(i,k-1,j) - (1./grid%em_rdnw(k-1))*( &
396 (grid%em_mub(i,j)+grid%em_mu_1(i,j))*grid%em_al(i,k-1,j)+ &
397 grid%em_mu_1(i,j)*grid%em_alb(i,k-1,j) )
398
399 grid%em_ph_2(i,k,j) = grid%em_ph_1(i,k,j)
400 grid%em_ph0(i,k,j) = grid%em_ph_1(i,k,j) + grid%em_phb(i,k,j)
401 ENDDO
402
403 if((i==2) .and. (j==2)) then
404 write(6,*) ' grid%em_ph_1 calc ',grid%em_ph_1(2,1,2),grid%em_ph_1(2,2,2),&
405 grid%em_mu_1(2,2)+grid%em_mub(2,2),grid%em_mu_1(2,2), &
406 grid%em_alb(2,1,2),grid%em_al(1,2,1),grid%em_rdnw(1)
407 endif
408
409 ENDDO
410 ENDDO
411
412 !#if 0
413
414 ! thermal perturbation to kick off convection
415
416 write(6,*) ' nxc, nyc for perturbation ',nxc,nyc
417 write(6,*) ' delt for perturbation ',delt
418
419 DO J = jts, min(jde-1,jte)
420 ! yrad = config_flags%dy*float(j-nyc)/4000.
421 yrad = 0.
422 DO I = its, min(ide-1,ite)
423 xrad = config_flags%dx*float(i-nxc)/4000.
424 ! xrad = 0.
425 DO K = 1, kte-1
426
427 ! put in preturbation theta (bubble) and recalc density. note,
428 ! the mass in the column is not changing, so when theta changes,
429 ! we recompute density and geopotential
430
431 zrad = 0.5*(grid%em_ph_1(i,k,j)+grid%em_ph_1(i,k+1,j) &
432 +grid%em_phb(i,k,j)+grid%em_phb(i,k+1,j))/g
433 zrad = (zrad-1500.)/1500.
434 RAD=SQRT(xrad*xrad+yrad*yrad+zrad*zrad)
435 IF(RAD <= 1.) THEN
436 grid%em_t_1(i,k,j)=grid%em_t_1(i,k,j)+delt*COS(.5*PI*RAD)**2
437 grid%em_t_2(i,k,j)=grid%em_t_1(i,k,j)
438 qvf = 1. + rvovrd*moist(i,k,j,P_QV)
439 grid%em_alt(i,k,j) = (r_d/p1000mb)*(grid%em_t_1(i,k,j)+t0)*qvf* &
440 (((grid%em_p(i,k,j)+grid%em_pb(i,k,j))/p1000mb)**cvpm)
441 grid%em_al(i,k,j) = grid%em_alt(i,k,j) - grid%em_alb(i,k,j)
442 ENDIF
443 ENDDO
444
445 ! rebalance hydrostatically
446
447 DO k = 2,kte
448 grid%em_ph_1(i,k,j) = grid%em_ph_1(i,k-1,j) - (1./grid%em_rdnw(k-1))*( &
449 (grid%em_mub(i,j)+grid%em_mu_1(i,j))*grid%em_al(i,k-1,j)+ &
450 grid%em_mu_1(i,j)*grid%em_alb(i,k-1,j) )
451
452 grid%em_ph_2(i,k,j) = grid%em_ph_1(i,k,j)
453 grid%em_ph0(i,k,j) = grid%em_ph_1(i,k,j) + grid%em_phb(i,k,j)
454 ENDDO
455
456 ENDDO
457 ENDDO
458
459 !#endif
460
461 write(6,*) ' grid%em_mu_1 from comp ', grid%em_mu_1(1,1)
462 write(6,*) ' full state sounding from comp, ph, grid%em_p, grid%em_al, grid%em_t_1, qv '
463 do k=1,kde-1
464 write(6,'(i3,1x,5(1x,1pe10.3))') k, grid%em_ph_1(1,k,1)+grid%em_phb(1,k,1), &
465 grid%em_p(1,k,1)+grid%em_pb(1,k,1), grid%em_alt(1,k,1), &
466 grid%em_t_1(1,k,1)+t0, moist(1,k,1,P_QV)
467 enddo
468
469 write(6,*) ' pert state sounding from comp, grid%em_ph_1, pp, alp, grid%em_t_1, qv '
470 do k=1,kde-1
471 write(6,'(i3,1x,5(1x,1pe10.3))') k, grid%em_ph_1(1,k,1), &
472 grid%em_p(1,k,1), grid%em_al(1,k,1), &
473 grid%em_t_1(1,k,1), moist(1,k,1,P_QV)
474 enddo
475
476 ! interp v
477
478 DO J = jts, jte
479 DO I = its, min(ide-1,ite)
480
481 IF (j == jds) THEN
482 z_at_v = grid%em_phb(i,1,j)/g
483 ELSE IF (j == jde) THEN
484 z_at_v = grid%em_phb(i,1,j-1)/g
485 ELSE
486 z_at_v = 0.5*(grid%em_phb(i,1,j)+grid%em_phb(i,1,j-1))/g
487 END IF
488
489 p_surf = interp_0( p_in, zk, z_at_v, nl_in )
490
491 DO K = 1, kte
492 p_level = grid%em_znu(k)*(p_surf - grid%p_top) + grid%p_top
493 grid%em_v_1(i,k,j) = interp_0( v, p_in, p_level, nl_in )
494 grid%em_v_2(i,k,j) = grid%em_v_1(i,k,j)
495 ENDDO
496
497 ENDDO
498 ENDDO
499
500 ! interp u
501
502 DO J = jts, min(jde-1,jte)
503 DO I = its, ite
504
505 IF (i == ids) THEN
506 z_at_u = grid%em_phb(i,1,j)/g
507 ELSE IF (i == ide) THEN
508 z_at_u = grid%em_phb(i-1,1,j)/g
509 ELSE
510 z_at_u = 0.5*(grid%em_phb(i,1,j)+grid%em_phb(i-1,1,j))/g
511 END IF
512
513 p_surf = interp_0( p_in, zk, z_at_u, nl_in )
514
515 DO K = 1, kte
516 p_level = grid%em_znu(k)*(p_surf - grid%p_top) + grid%p_top
517 grid%em_u_1(i,k,j) = interp_0( u, p_in, p_level, nl_in )
518 grid%em_u_2(i,k,j) = grid%em_u_1(i,k,j)
519 ENDDO
520
521 ENDDO
522 ENDDO
523
524 ! set w
525
526 DO J = jts, min(jde-1,jte)
527 DO K = kts, kte
528 DO I = its, min(ide-1,ite)
529 grid%em_w_1(i,k,j) = 0.
530 grid%em_w_2(i,k,j) = 0.
531 ENDDO
532 ENDDO
533 ENDDO
534
535 ! set a few more things
536
537 DO J = jts, min(jde-1,jte)
538 DO K = kts, kte-1
539 DO I = its, min(ide-1,ite)
540 grid%h_diabatic(i,k,j) = 0.
541 ENDDO
542 ENDDO
543 ENDDO
544
545 DO k=1,kte-1
546 grid%em_t_base(k) = grid%em_t_1(1,k,1)
547 grid%qv_base(k) = moist(1,k,1,P_QV)
548 grid%u_base(k) = grid%em_u_1(1,k,1)
549 grid%v_base(k) = grid%em_v_1(1,k,1)
550 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
551 ENDDO
552
553 DO J = jts, min(jde-1,jte)
554 DO I = its, min(ide-1,ite)
555 thtmp = grid%em_t_2(i,1,j)+t0
556 ptmp = grid%em_p(i,1,j)+grid%em_pb(i,1,j)
557 temp(1) = thtmp * (ptmp/p1000mb)**rcp
558 thtmp = grid%em_t_2(i,2,j)+t0
559 ptmp = grid%em_p(i,2,j)+grid%em_pb(i,2,j)
560 temp(2) = thtmp * (ptmp/p1000mb)**rcp
561 thtmp = grid%em_t_2(i,3,j)+t0
562 ptmp = grid%em_p(i,3,j)+grid%em_pb(i,3,j)
563 temp(3) = thtmp * (ptmp/p1000mb)**rcp
564
565 grid%tsk(I,J)=grid%cf1*temp(1)+grid%cf2*temp(2)+grid%cf3*temp(3)
566 grid%tmn(I,J)=grid%tsk(I,J)-0.5
567 ENDDO
568 ENDDO
569
570 RETURN
571
572 END SUBROUTINE init_domain_rk
573
574 SUBROUTINE init_module_initialize
575 END SUBROUTINE init_module_initialize
576
577 !---------------------------------------------------------------------
578
579 ! test driver for get_sounding
580 !
581 ! implicit none
582 ! integer n
583 ! parameter(n = 1000)
584 ! real zk(n),p(n),theta(n),rho(n),u(n),v(n),qv(n),pd(n)
585 ! logical dry
586 ! integer nl,k
587 !
588 ! dry = .false.
589 ! dry = .true.
590 ! call get_sounding( zk, p, pd, theta, rho, u, v, qv, dry, n, nl )
591 ! write(6,*) ' input levels ',nl
592 ! write(6,*) ' sounding '
593 ! write(6,*) ' k height(m) press (Pa) pd(Pa) theta (K) den(kg/m^3) u(m/s) v(m/s) qv(g/g) '
594 ! do k=1,nl
595 ! 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)
596 ! enddo
597 ! end
598 !
599 !---------------------------------------------------------------------------
600
601 subroutine get_sounding( zk, p, p_dry, theta, rho, &
602 u, v, qv, dry, nl_max, nl_in )
603 implicit none
604
605 integer nl_max, nl_in
606 real zk(nl_max), p(nl_max), theta(nl_max), rho(nl_max), &
607 u(nl_max), v(nl_max), qv(nl_max), p_dry(nl_max)
608 logical dry
609
610 integer n
611 parameter(n=1000)
612 logical debug
613 parameter( debug = .true.)
614
615 ! input sounding data
616
617 real p_surf, th_surf, qv_surf
618 real pi_surf, pi(n)
619 real h_input(n), th_input(n), qv_input(n), u_input(n), v_input(n)
620
621 ! diagnostics
622
623 real rho_surf, p_input(n), rho_input(n)
624 real pm_input(n) ! this are for full moist sounding
625
626 ! local data
627
628 real p1000mb,cv,cp,r,cvpm,g
629 parameter (p1000mb = 1.e+05, r = 287, cp = 1003., cv = cp-r, cvpm = -cv/cp, g=9.81 )
630 integer k, it, nl
631 real qvf, qvf1, dz
632
633 ! first, read the sounding
634
635 call read_sounding( p_surf, th_surf, qv_surf, &
636 h_input, th_input, qv_input, u_input, v_input,n, nl, debug )
637
638 if(dry) then
639 do k=1,nl
640 qv_input(k) = 0.
641 enddo
642 endif
643
644 if(debug) write(6,*) ' number of input levels = ',nl
645
646 nl_in = nl
647 if(nl_in .gt. nl_max ) then
648 write(6,*) ' too many levels for input arrays ',nl_in,nl_max
649 call wrf_error_fatal ( ' too many levels for input arrays ' )
650 end if
651
652 ! compute diagnostics,
653 ! first, convert qv(g/kg) to qv(g/g)
654
655 do k=1,nl
656 qv_input(k) = 0.001*qv_input(k)
657 enddo
658
659 p_surf = 100.*p_surf ! convert to pascals
660 qvf = 1. + rvovrd*qv_input(1)
661 rho_surf = 1./((r/p1000mb)*th_surf*qvf*((p_surf/p1000mb)**cvpm))
662 pi_surf = (p_surf/p1000mb)**(r/cp)
663
664 if(debug) then
665 write(6,*) ' surface density is ',rho_surf
666 write(6,*) ' surface pi is ',pi_surf
667 end if
668
669
670 ! integrate moist sounding hydrostatically, starting from the
671 ! specified surface pressure
672 ! -> first, integrate from surface to lowest level
673
674 qvf = 1. + rvovrd*qv_input(1)
675 qvf1 = 1. + qv_input(1)
676 rho_input(1) = rho_surf
677 dz = h_input(1)
678 do it=1,10
679 pm_input(1) = p_surf &
680 - 0.5*dz*(rho_surf+rho_input(1))*g*qvf1
681 rho_input(1) = 1./((r/p1000mb)*th_input(1)*qvf*((pm_input(1)/p1000mb)**cvpm))
682 enddo
683
684 ! integrate up the column
685
686 do k=2,nl
687 rho_input(k) = rho_input(k-1)
688 dz = h_input(k)-h_input(k-1)
689 qvf1 = 0.5*(2.+(qv_input(k-1)+qv_input(k)))
690 qvf = 1. + rvovrd*qv_input(k) ! qv is in g/kg here
691
692 do it=1,10
693 pm_input(k) = pm_input(k-1) &
694 - 0.5*dz*(rho_input(k)+rho_input(k-1))*g*qvf1
695 rho_input(k) = 1./((r/p1000mb)*th_input(k)*qvf*((pm_input(k)/p1000mb)**cvpm))
696 enddo
697 enddo
698
699 ! we have the moist sounding
700
701 ! next, compute the dry sounding using p at the highest level from the
702 ! moist sounding and integrating down.
703
704 p_input(nl) = pm_input(nl)
705
706 do k=nl-1,1,-1
707 dz = h_input(k+1)-h_input(k)
708 p_input(k) = p_input(k+1) + 0.5*dz*(rho_input(k)+rho_input(k+1))*g
709 enddo
710
711
712 do k=1,nl
713
714 zk(k) = h_input(k)
715 p(k) = pm_input(k)
716 p_dry(k) = p_input(k)
717 theta(k) = th_input(k)
718 rho(k) = rho_input(k)
719 u(k) = u_input(k)
720 v(k) = v_input(k)
721 qv(k) = qv_input(k)
722
723 enddo
724
725 if(debug) then
726 write(6,*) ' sounding '
727 write(6,*) ' k height(m) press (Pa) pd(Pa) theta (K) den(kg/m^3) u(m/s) v(m/s) qv(g/g) '
728 do k=1,nl
729 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)
730 enddo
731
732 end if
733
734 end subroutine get_sounding
735
736 !-------------------------------------------------------
737
738 subroutine read_sounding( ps,ts,qvs,h,th,qv,u,v,n,nl,debug )
739 implicit none
740 integer n,nl
741 real ps,ts,qvs,h(n),th(n),qv(n),u(n),v(n)
742 logical end_of_file
743 logical debug
744
745 integer k
746
747 open(unit=10,file='input_sounding',form='formatted',status='old')
748 rewind(10)
749 read(10,*) ps, ts, qvs
750 if(debug) then
751 write(6,*) ' input sounding surface parameters '
752 write(6,*) ' surface pressure (mb) ',ps
753 write(6,*) ' surface pot. temp (K) ',ts
754 write(6,*) ' surface mixing ratio (g/kg) ',qvs
755 end if
756
757 end_of_file = .false.
758 k = 0
759
760 do while (.not. end_of_file)
761
762 read(10,*,end=100) h(k+1), th(k+1), qv(k+1), u(k+1), v(k+1)
763 k = k+1
764 if(debug) write(6,'(1x,i3,5(1x,e10.3))') k, h(k), th(k), qv(k), u(k), v(k)
765 go to 110
766 100 end_of_file = .true.
767 110 continue
768 enddo
769
770 nl = k
771
772 close(unit=10,status = 'keep')
773
774 end subroutine read_sounding
775
776 END MODULE module_initialize_ideal