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 :: idum1, idum2 |
---|
54 | |
---|
55 | CALL set_scalar_indices_from_config ( head_grid%id , idum1, idum2 ) |
---|
56 | |
---|
57 | CALL init_domain_rk( grid & |
---|
58 | ! |
---|
59 | #include <actual_new_args.inc> |
---|
60 | ! |
---|
61 | ) |
---|
62 | END SUBROUTINE init_domain |
---|
63 | |
---|
64 | !------------------------------------------------------------------- |
---|
65 | |
---|
66 | SUBROUTINE init_domain_rk ( grid & |
---|
67 | ! |
---|
68 | # include <dummy_new_args.inc> |
---|
69 | ! |
---|
70 | ) |
---|
71 | IMPLICIT NONE |
---|
72 | |
---|
73 | ! Input data. |
---|
74 | TYPE (domain), POINTER :: grid |
---|
75 | |
---|
76 | # include <dummy_new_decl.inc> |
---|
77 | |
---|
78 | TYPE (grid_config_rec_type) :: config_flags |
---|
79 | |
---|
80 | ! Local data |
---|
81 | INTEGER :: & |
---|
82 | ids, ide, jds, jde, kds, kde, & |
---|
83 | ims, ime, jms, jme, kms, kme, & |
---|
84 | its, ite, jts, jte, kts, kte, & |
---|
85 | i, j, k |
---|
86 | |
---|
87 | ! Local data |
---|
88 | |
---|
89 | INTEGER, PARAMETER :: nl_max = 1000 |
---|
90 | REAL, DIMENSION(nl_max) :: zk, p_in, theta, rho, u, v, qv, pd_in |
---|
91 | INTEGER :: nl_in |
---|
92 | |
---|
93 | |
---|
94 | INTEGER :: icm,jcm, ii, im1, jj, jm1, loop, error, fid, nxc, nyc |
---|
95 | REAL :: u_mean,v_mean, f0, p_surf, p_level, qvf, z_at_v, z_at_u |
---|
96 | REAL :: z_scale, xrad, yrad, zrad, rad, delt, cof1, cof2 |
---|
97 | ! REAL, EXTERNAL :: interp_0 |
---|
98 | REAL :: hm |
---|
99 | REAL :: pi |
---|
100 | |
---|
101 | ! stuff from original initialization that has been dropped from the Registry |
---|
102 | REAL :: vnu, xnu, xnus, dinit0, cbh, p0_temp, t0_temp, zd, zt |
---|
103 | REAL :: qvf1, qvf2, pd_surf |
---|
104 | INTEGER :: it |
---|
105 | real :: thtmp, ptmp, temp(3) |
---|
106 | |
---|
107 | LOGICAL :: moisture_init |
---|
108 | LOGICAL :: stretch_grid, dry_sounding |
---|
109 | |
---|
110 | SELECT CASE ( model_data_order ) |
---|
111 | CASE ( DATA_ORDER_ZXY ) |
---|
112 | kds = grid%sd31 ; kde = grid%ed31 ; |
---|
113 | ids = grid%sd32 ; ide = grid%ed32 ; |
---|
114 | jds = grid%sd33 ; jde = grid%ed33 ; |
---|
115 | |
---|
116 | kms = grid%sm31 ; kme = grid%em31 ; |
---|
117 | ims = grid%sm32 ; ime = grid%em32 ; |
---|
118 | jms = grid%sm33 ; jme = grid%em33 ; |
---|
119 | |
---|
120 | kts = grid%sp31 ; kte = grid%ep31 ; ! note that tile is entire patch |
---|
121 | its = grid%sp32 ; ite = grid%ep32 ; ! note that tile is entire patch |
---|
122 | jts = grid%sp33 ; jte = grid%ep33 ; ! note that tile is entire patch |
---|
123 | CASE ( DATA_ORDER_XYZ ) |
---|
124 | ids = grid%sd31 ; ide = grid%ed31 ; |
---|
125 | jds = grid%sd32 ; jde = grid%ed32 ; |
---|
126 | kds = grid%sd33 ; kde = grid%ed33 ; |
---|
127 | |
---|
128 | ims = grid%sm31 ; ime = grid%em31 ; |
---|
129 | jms = grid%sm32 ; jme = grid%em32 ; |
---|
130 | kms = grid%sm33 ; kme = grid%em33 ; |
---|
131 | |
---|
132 | its = grid%sp31 ; ite = grid%ep31 ; ! note that tile is entire patch |
---|
133 | jts = grid%sp32 ; jte = grid%ep32 ; ! note that tile is entire patch |
---|
134 | kts = grid%sp33 ; kte = grid%ep33 ; ! note that tile is entire patch |
---|
135 | CASE ( DATA_ORDER_XZY ) |
---|
136 | ids = grid%sd31 ; ide = grid%ed31 ; |
---|
137 | kds = grid%sd32 ; kde = grid%ed32 ; |
---|
138 | jds = grid%sd33 ; jde = grid%ed33 ; |
---|
139 | |
---|
140 | ims = grid%sm31 ; ime = grid%em31 ; |
---|
141 | kms = grid%sm32 ; kme = grid%em32 ; |
---|
142 | jms = grid%sm33 ; jme = grid%em33 ; |
---|
143 | |
---|
144 | its = grid%sp31 ; ite = grid%ep31 ; ! note that tile is entire patch |
---|
145 | kts = grid%sp32 ; kte = grid%ep32 ; ! note that tile is entire patch |
---|
146 | jts = grid%sp33 ; jte = grid%ep33 ; ! note that tile is entire patch |
---|
147 | |
---|
148 | END SELECT |
---|
149 | |
---|
150 | |
---|
151 | stretch_grid = .true. |
---|
152 | delt = 3. |
---|
153 | ! z_scale = .50 |
---|
154 | z_scale = .40 |
---|
155 | pi = 2.*asin(1.0) |
---|
156 | write(6,*) ' pi is ',pi |
---|
157 | nxc = (ide-ids)/2 |
---|
158 | nyc = jde/2 |
---|
159 | |
---|
160 | CALL model_to_grid_config_rec ( grid%id , model_config_rec , config_flags ) |
---|
161 | |
---|
162 | ! here we check to see if the boundary conditions are set properly |
---|
163 | |
---|
164 | CALL boundary_condition_check( config_flags, bdyzone, error, grid%id ) |
---|
165 | |
---|
166 | moisture_init = .true. |
---|
167 | |
---|
168 | grid%itimestep=0 |
---|
169 | |
---|
170 | #ifdef DM_PARALLEL |
---|
171 | CALL wrf_dm_bcast_bytes( icm , IWORDSIZE ) |
---|
172 | CALL wrf_dm_bcast_bytes( jcm , IWORDSIZE ) |
---|
173 | #endif |
---|
174 | |
---|
175 | CALL nl_set_mminlu(1, ' ') |
---|
176 | CALL nl_set_iswater(1,0) |
---|
177 | CALL nl_set_cen_lat(1,40.) |
---|
178 | CALL nl_set_cen_lon(1,-105.) |
---|
179 | CALL nl_set_truelat1(1,0.) |
---|
180 | CALL nl_set_truelat2(1,0.) |
---|
181 | CALL nl_set_moad_cen_lat (1,0.) |
---|
182 | CALL nl_set_stand_lon (1,0.) |
---|
183 | CALL nl_set_pole_lon (1,0.) |
---|
184 | CALL nl_set_pole_lat (1,90.) |
---|
185 | CALL nl_set_map_proj(1,0) |
---|
186 | |
---|
187 | |
---|
188 | ! here we initialize data we currently is not initialized |
---|
189 | ! in the input data |
---|
190 | |
---|
191 | DO j = jts, jte |
---|
192 | DO i = its, ite |
---|
193 | grid%msftx(i,j) = 1. |
---|
194 | grid%msfty(i,j) = 1. |
---|
195 | grid%msfux(i,j) = 1. |
---|
196 | grid%msfuy(i,j) = 1. |
---|
197 | grid%msfvx(i,j) = 1. |
---|
198 | grid%msfvx_inv(i,j)= 1. |
---|
199 | grid%msfvy(i,j) = 1. |
---|
200 | grid%sina(i,j) = 0. |
---|
201 | grid%cosa(i,j) = 1. |
---|
202 | grid%e(i,j) = 0. |
---|
203 | grid%f(i,j) = 0. |
---|
204 | |
---|
205 | END DO |
---|
206 | END DO |
---|
207 | |
---|
208 | DO j = jts, jte |
---|
209 | DO k = kts, kte |
---|
210 | DO i = its, ite |
---|
211 | grid%ww(i,k,j) = 0. |
---|
212 | END DO |
---|
213 | END DO |
---|
214 | END DO |
---|
215 | |
---|
216 | grid%step_number = 0 |
---|
217 | |
---|
218 | ! set up the grid |
---|
219 | |
---|
220 | IF (stretch_grid) THEN ! exponential stretch for eta (nearly constant dz) |
---|
221 | DO k=1, kde |
---|
222 | grid%znw(k) = (exp(-(k-1)/float(kde-1)/z_scale) - exp(-1./z_scale))/ & |
---|
223 | (1.-exp(-1./z_scale)) |
---|
224 | ENDDO |
---|
225 | ELSE |
---|
226 | DO k=1, kde |
---|
227 | grid%znw(k) = 1. - float(k-1)/float(kde-1) |
---|
228 | ENDDO |
---|
229 | ENDIF |
---|
230 | |
---|
231 | DO k=1, kde-1 |
---|
232 | grid%dnw(k) = grid%znw(k+1) - grid%znw(k) |
---|
233 | grid%rdnw(k) = 1./grid%dnw(k) |
---|
234 | grid%znu(k) = 0.5*(grid%znw(k+1)+grid%znw(k)) |
---|
235 | ENDDO |
---|
236 | DO k=2, kde-1 |
---|
237 | grid%dn(k) = 0.5*(grid%dnw(k)+grid%dnw(k-1)) |
---|
238 | grid%rdn(k) = 1./grid%dn(k) |
---|
239 | grid%fnp(k) = .5* grid%dnw(k )/grid%dn(k) |
---|
240 | grid%fnm(k) = .5* grid%dnw(k-1)/grid%dn(k) |
---|
241 | ENDDO |
---|
242 | |
---|
243 | cof1 = (2.*grid%dn(2)+grid%dn(3))/(grid%dn(2)+grid%dn(3))*grid%dnw(1)/grid%dn(2) |
---|
244 | cof2 = grid%dn(2) /(grid%dn(2)+grid%dn(3))*grid%dnw(1)/grid%dn(3) |
---|
245 | grid%cf1 = grid%fnp(2) + cof1 |
---|
246 | grid%cf2 = grid%fnm(2) - cof1 - cof2 |
---|
247 | grid%cf3 = cof2 |
---|
248 | |
---|
249 | grid%cfn = (.5*grid%dnw(kde-1)+grid%dn(kde-1))/grid%dn(kde-1) |
---|
250 | grid%cfn1 = -.5*grid%dnw(kde-1)/grid%dn(kde-1) |
---|
251 | grid%rdx = 1./config_flags%dx |
---|
252 | grid%rdy = 1./config_flags%dy |
---|
253 | |
---|
254 | ! get the sounding from the ascii sounding file, first get dry sounding and |
---|
255 | ! calculate base state |
---|
256 | |
---|
257 | write(6,*) ' getting dry sounding for base state ' |
---|
258 | dry_sounding = .true. |
---|
259 | CALL get_sounding( zk, p_in, pd_in, theta, rho, u, v, qv, dry_sounding, nl_max, nl_in ) |
---|
260 | |
---|
261 | write(6,*) ' returned from reading sounding, nl_in is ',nl_in |
---|
262 | |
---|
263 | ! find ptop for the desired ztop (ztop is input from the namelist), |
---|
264 | ! and find surface pressure |
---|
265 | |
---|
266 | grid%p_top = interp_0( p_in, zk, config_flags%ztop, nl_in ) |
---|
267 | |
---|
268 | DO j=jts,jte |
---|
269 | DO i=its,ite ! flat surface |
---|
270 | grid%phb(i,1,j) = 0. |
---|
271 | grid%php(i,1,j) = 0. |
---|
272 | grid%ph0(i,1,j) = 0. |
---|
273 | grid%ht(i,j) = 0. |
---|
274 | ENDDO |
---|
275 | ENDDO |
---|
276 | |
---|
277 | DO J = jts, jte |
---|
278 | DO I = its, ite |
---|
279 | |
---|
280 | p_surf = interp_0( p_in, zk, grid%phb(i,1,j)/g, nl_in ) |
---|
281 | grid%mub(i,j) = p_surf-grid%p_top |
---|
282 | |
---|
283 | ! this is dry hydrostatic sounding (base state), so given grid%p (coordinate), |
---|
284 | ! interp theta (from interp) and compute 1/rho from eqn. of state |
---|
285 | |
---|
286 | DO K = 1, kte-1 |
---|
287 | p_level = grid%znu(k)*(p_surf - grid%p_top) + grid%p_top |
---|
288 | grid%pb(i,k,j) = p_level |
---|
289 | grid%t_init(i,k,j) = interp_0( theta, p_in, p_level, nl_in ) - t0 |
---|
290 | grid%alb(i,k,j) = (r_d/p1000mb)*(grid%t_init(i,k,j)+t0)*(grid%pb(i,k,j)/p1000mb)**cvpm |
---|
291 | ENDDO |
---|
292 | |
---|
293 | ! calc hydrostatic balance (alternatively we could interp the geopotential from the |
---|
294 | ! sounding, but this assures that the base state is in exact hydrostatic balance with |
---|
295 | ! respect to the model eqns. |
---|
296 | |
---|
297 | DO k = 2,kte |
---|
298 | grid%phb(i,k,j) = grid%phb(i,k-1,j) - grid%dnw(k-1)*grid%mub(i,j)*grid%alb(i,k-1,j) |
---|
299 | ENDDO |
---|
300 | |
---|
301 | ENDDO |
---|
302 | ENDDO |
---|
303 | |
---|
304 | write(6,*) ' ptop is ',grid%p_top |
---|
305 | write(6,*) ' base state grid%mub(1,1), p_surf is ',grid%mub(1,1),grid%mub(1,1)+grid%p_top |
---|
306 | |
---|
307 | ! calculate full state for each column - this includes moisture. |
---|
308 | |
---|
309 | write(6,*) ' getting moist sounding for full state ' |
---|
310 | dry_sounding = .false. |
---|
311 | CALL get_sounding( zk, p_in, pd_in, theta, rho, u, v, qv, dry_sounding, nl_max, nl_in ) |
---|
312 | |
---|
313 | DO J = jts, min(jde-1,jte) |
---|
314 | DO I = its, min(ide-1,ite) |
---|
315 | |
---|
316 | ! At this point grid%p_top is already set. find the DRY mass in the column |
---|
317 | ! by interpolating the DRY pressure. |
---|
318 | |
---|
319 | pd_surf = interp_0( pd_in, zk, grid%phb(i,1,j)/g, nl_in ) |
---|
320 | |
---|
321 | ! compute the perturbation mass and the full mass |
---|
322 | |
---|
323 | grid%mu_1(i,j) = pd_surf-grid%p_top - grid%mub(i,j) |
---|
324 | grid%mu_2(i,j) = grid%mu_1(i,j) |
---|
325 | grid%mu0(i,j) = grid%mu_1(i,j) + grid%mub(i,j) |
---|
326 | |
---|
327 | ! given the dry pressure and coordinate system, interp the potential |
---|
328 | ! temperature and qv |
---|
329 | |
---|
330 | do k=1,kde-1 |
---|
331 | |
---|
332 | p_level = grid%znu(k)*(pd_surf - grid%p_top) + grid%p_top |
---|
333 | |
---|
334 | moist(i,k,j,P_QV) = interp_0( qv, pd_in, p_level, nl_in ) |
---|
335 | grid%t_1(i,k,j) = interp_0( theta, pd_in, p_level, nl_in ) - t0 |
---|
336 | grid%t_2(i,k,j) = grid%t_1(i,k,j) |
---|
337 | |
---|
338 | |
---|
339 | enddo |
---|
340 | |
---|
341 | ! integrate the hydrostatic equation (from the RHS of the bigstep |
---|
342 | ! vertical momentum equation) down from the top to get grid%p. |
---|
343 | ! first from the top of the model to the top pressure |
---|
344 | |
---|
345 | k = kte-1 ! top level |
---|
346 | |
---|
347 | qvf1 = 0.5*(moist(i,k,j,P_QV)+moist(i,k,j,P_QV)) |
---|
348 | qvf2 = 1./(1.+qvf1) |
---|
349 | qvf1 = qvf1*qvf2 |
---|
350 | |
---|
351 | ! grid%p(i,k,j) = - 0.5*grid%mu_1(i,j)/grid%rdnw(k) |
---|
352 | grid%p(i,k,j) = - 0.5*(grid%mu_1(i,j)+qvf1*grid%mub(i,j))/grid%rdnw(k)/qvf2 |
---|
353 | qvf = 1. + rvovrd*moist(i,k,j,P_QV) |
---|
354 | grid%alt(i,k,j) = (r_d/p1000mb)*(grid%t_1(i,k,j)+t0)*qvf* & |
---|
355 | (((grid%p(i,k,j)+grid%pb(i,k,j))/p1000mb)**cvpm) |
---|
356 | grid%al(i,k,j) = grid%alt(i,k,j) - grid%alb(i,k,j) |
---|
357 | |
---|
358 | ! down the column |
---|
359 | |
---|
360 | do k=kte-2,1,-1 |
---|
361 | qvf1 = 0.5*(moist(i,k,j,P_QV)+moist(i,k+1,j,P_QV)) |
---|
362 | qvf2 = 1./(1.+qvf1) |
---|
363 | qvf1 = qvf1*qvf2 |
---|
364 | grid%p(i,k,j) = grid%p(i,k+1,j) - (grid%mu_1(i,j) + qvf1*grid%mub(i,j))/qvf2/grid%rdn(k+1) |
---|
365 | qvf = 1. + rvovrd*moist(i,k,j,P_QV) |
---|
366 | grid%alt(i,k,j) = (r_d/p1000mb)*(grid%t_1(i,k,j)+t0)*qvf* & |
---|
367 | (((grid%p(i,k,j)+grid%pb(i,k,j))/p1000mb)**cvpm) |
---|
368 | grid%al(i,k,j) = grid%alt(i,k,j) - grid%alb(i,k,j) |
---|
369 | enddo |
---|
370 | |
---|
371 | ! this is the hydrostatic equation used in the model after the |
---|
372 | ! small timesteps. In the model, grid%al (inverse density) |
---|
373 | ! is computed from the geopotential. |
---|
374 | |
---|
375 | |
---|
376 | grid%ph_1(i,1,j) = 0. |
---|
377 | DO k = 2,kte |
---|
378 | grid%ph_1(i,k,j) = grid%ph_1(i,k-1,j) - (1./grid%rdnw(k-1))*( & |
---|
379 | (grid%mub(i,j)+grid%mu_1(i,j))*grid%al(i,k-1,j)+ & |
---|
380 | grid%mu_1(i,j)*grid%alb(i,k-1,j) ) |
---|
381 | |
---|
382 | grid%ph_2(i,k,j) = grid%ph_1(i,k,j) |
---|
383 | grid%ph0(i,k,j) = grid%ph_1(i,k,j) + grid%phb(i,k,j) |
---|
384 | ENDDO |
---|
385 | |
---|
386 | if((i==2) .and. (j==2)) then |
---|
387 | write(6,*) ' grid%ph_1 calc ',grid%ph_1(2,1,2),grid%ph_1(2,2,2),& |
---|
388 | grid%mu_1(2,2)+grid%mub(2,2),grid%mu_1(2,2), & |
---|
389 | grid%alb(2,1,2),grid%al(1,2,1),grid%rdnw(1) |
---|
390 | endif |
---|
391 | |
---|
392 | ENDDO |
---|
393 | ENDDO |
---|
394 | |
---|
395 | !#if 0 |
---|
396 | |
---|
397 | ! thermal perturbation to kick off convection |
---|
398 | |
---|
399 | write(6,*) ' nxc, nyc for perturbation ',nxc,nyc |
---|
400 | write(6,*) ' delt for perturbation ',delt |
---|
401 | |
---|
402 | DO J = jts, min(jde-1,jte) |
---|
403 | ! yrad = config_flags%dy*float(j-nyc)/4000. |
---|
404 | yrad = 0. |
---|
405 | DO I = its, min(ide-1,ite) |
---|
406 | xrad = config_flags%dx*float(i-nxc)/4000. |
---|
407 | ! xrad = 0. |
---|
408 | DO K = 1, kte-1 |
---|
409 | |
---|
410 | ! put in preturbation theta (bubble) and recalc density. note, |
---|
411 | ! the mass in the column is not changing, so when theta changes, |
---|
412 | ! we recompute density and geopotential |
---|
413 | |
---|
414 | zrad = 0.5*(grid%ph_1(i,k,j)+grid%ph_1(i,k+1,j) & |
---|
415 | +grid%phb(i,k,j)+grid%phb(i,k+1,j))/g |
---|
416 | zrad = (zrad-1500.)/1500. |
---|
417 | RAD=SQRT(xrad*xrad+yrad*yrad+zrad*zrad) |
---|
418 | IF(RAD <= 1.) THEN |
---|
419 | grid%t_1(i,k,j)=grid%t_1(i,k,j)+delt*COS(.5*PI*RAD)**2 |
---|
420 | grid%t_2(i,k,j)=grid%t_1(i,k,j) |
---|
421 | qvf = 1. + rvovrd*moist(i,k,j,P_QV) |
---|
422 | grid%alt(i,k,j) = (r_d/p1000mb)*(grid%t_1(i,k,j)+t0)*qvf* & |
---|
423 | (((grid%p(i,k,j)+grid%pb(i,k,j))/p1000mb)**cvpm) |
---|
424 | grid%al(i,k,j) = grid%alt(i,k,j) - grid%alb(i,k,j) |
---|
425 | ENDIF |
---|
426 | ENDDO |
---|
427 | |
---|
428 | ! rebalance hydrostatically |
---|
429 | |
---|
430 | DO k = 2,kte |
---|
431 | grid%ph_1(i,k,j) = grid%ph_1(i,k-1,j) - (1./grid%rdnw(k-1))*( & |
---|
432 | (grid%mub(i,j)+grid%mu_1(i,j))*grid%al(i,k-1,j)+ & |
---|
433 | grid%mu_1(i,j)*grid%alb(i,k-1,j) ) |
---|
434 | |
---|
435 | grid%ph_2(i,k,j) = grid%ph_1(i,k,j) |
---|
436 | grid%ph0(i,k,j) = grid%ph_1(i,k,j) + grid%phb(i,k,j) |
---|
437 | ENDDO |
---|
438 | |
---|
439 | ENDDO |
---|
440 | ENDDO |
---|
441 | |
---|
442 | !#endif |
---|
443 | |
---|
444 | write(6,*) ' grid%mu_1 from comp ', grid%mu_1(1,1) |
---|
445 | write(6,*) ' full state sounding from comp, ph, grid%p, grid%al, grid%t_1, qv ' |
---|
446 | do k=1,kde-1 |
---|
447 | write(6,'(i3,1x,5(1x,1pe10.3))') k, grid%ph_1(1,k,1)+grid%phb(1,k,1), & |
---|
448 | grid%p(1,k,1)+grid%pb(1,k,1), grid%alt(1,k,1), & |
---|
449 | grid%t_1(1,k,1)+t0, moist(1,k,1,P_QV) |
---|
450 | enddo |
---|
451 | |
---|
452 | write(6,*) ' pert state sounding from comp, grid%ph_1, pp, alp, grid%t_1, qv ' |
---|
453 | do k=1,kde-1 |
---|
454 | write(6,'(i3,1x,5(1x,1pe10.3))') k, grid%ph_1(1,k,1), & |
---|
455 | grid%p(1,k,1), grid%al(1,k,1), & |
---|
456 | grid%t_1(1,k,1), moist(1,k,1,P_QV) |
---|
457 | enddo |
---|
458 | |
---|
459 | ! interp v |
---|
460 | |
---|
461 | DO J = jts, jte |
---|
462 | DO I = its, min(ide-1,ite) |
---|
463 | |
---|
464 | IF (j == jds) THEN |
---|
465 | z_at_v = grid%phb(i,1,j)/g |
---|
466 | ELSE IF (j == jde) THEN |
---|
467 | z_at_v = grid%phb(i,1,j-1)/g |
---|
468 | ELSE |
---|
469 | z_at_v = 0.5*(grid%phb(i,1,j)+grid%phb(i,1,j-1))/g |
---|
470 | END IF |
---|
471 | |
---|
472 | p_surf = interp_0( p_in, zk, z_at_v, nl_in ) |
---|
473 | |
---|
474 | DO K = 1, kte |
---|
475 | p_level = grid%znu(k)*(p_surf - grid%p_top) + grid%p_top |
---|
476 | grid%v_1(i,k,j) = interp_0( v, p_in, p_level, nl_in ) |
---|
477 | grid%v_2(i,k,j) = grid%v_1(i,k,j) |
---|
478 | ENDDO |
---|
479 | |
---|
480 | ENDDO |
---|
481 | ENDDO |
---|
482 | |
---|
483 | ! interp u |
---|
484 | |
---|
485 | DO J = jts, min(jde-1,jte) |
---|
486 | DO I = its, ite |
---|
487 | |
---|
488 | IF (i == ids) THEN |
---|
489 | z_at_u = grid%phb(i,1,j)/g |
---|
490 | ELSE IF (i == ide) THEN |
---|
491 | z_at_u = grid%phb(i-1,1,j)/g |
---|
492 | ELSE |
---|
493 | z_at_u = 0.5*(grid%phb(i,1,j)+grid%phb(i-1,1,j))/g |
---|
494 | END IF |
---|
495 | |
---|
496 | p_surf = interp_0( p_in, zk, z_at_u, nl_in ) |
---|
497 | |
---|
498 | DO K = 1, kte |
---|
499 | p_level = grid%znu(k)*(p_surf - grid%p_top) + grid%p_top |
---|
500 | grid%u_1(i,k,j) = interp_0( u, p_in, p_level, nl_in ) |
---|
501 | grid%u_2(i,k,j) = grid%u_1(i,k,j) |
---|
502 | ENDDO |
---|
503 | |
---|
504 | ENDDO |
---|
505 | ENDDO |
---|
506 | |
---|
507 | ! set w |
---|
508 | |
---|
509 | DO J = jts, min(jde-1,jte) |
---|
510 | DO K = kts, kte |
---|
511 | DO I = its, min(ide-1,ite) |
---|
512 | grid%w_1(i,k,j) = 0. |
---|
513 | grid%w_2(i,k,j) = 0. |
---|
514 | ENDDO |
---|
515 | ENDDO |
---|
516 | ENDDO |
---|
517 | |
---|
518 | ! set a few more things |
---|
519 | |
---|
520 | DO J = jts, min(jde-1,jte) |
---|
521 | DO K = kts, kte-1 |
---|
522 | DO I = its, min(ide-1,ite) |
---|
523 | grid%h_diabatic(i,k,j) = 0. |
---|
524 | ENDDO |
---|
525 | ENDDO |
---|
526 | ENDDO |
---|
527 | |
---|
528 | DO k=1,kte-1 |
---|
529 | grid%t_base(k) = grid%t_1(1,k,1) |
---|
530 | grid%qv_base(k) = moist(1,k,1,P_QV) |
---|
531 | grid%u_base(k) = grid%u_1(1,k,1) |
---|
532 | grid%v_base(k) = grid%v_1(1,k,1) |
---|
533 | grid%z_base(k) = 0.5*(grid%phb(1,k,1)+grid%phb(1,k+1,1)+grid%ph_1(1,k,1)+grid%ph_1(1,k+1,1))/g |
---|
534 | ENDDO |
---|
535 | |
---|
536 | DO J = jts, min(jde-1,jte) |
---|
537 | DO I = its, min(ide-1,ite) |
---|
538 | thtmp = grid%t_2(i,1,j)+t0 |
---|
539 | ptmp = grid%p(i,1,j)+grid%pb(i,1,j) |
---|
540 | temp(1) = thtmp * (ptmp/p1000mb)**rcp |
---|
541 | thtmp = grid%t_2(i,2,j)+t0 |
---|
542 | ptmp = grid%p(i,2,j)+grid%pb(i,2,j) |
---|
543 | temp(2) = thtmp * (ptmp/p1000mb)**rcp |
---|
544 | thtmp = grid%t_2(i,3,j)+t0 |
---|
545 | ptmp = grid%p(i,3,j)+grid%pb(i,3,j) |
---|
546 | temp(3) = thtmp * (ptmp/p1000mb)**rcp |
---|
547 | |
---|
548 | grid%tsk(I,J)=grid%cf1*temp(1)+grid%cf2*temp(2)+grid%cf3*temp(3) |
---|
549 | grid%tmn(I,J)=grid%tsk(I,J)-0.5 |
---|
550 | ENDDO |
---|
551 | ENDDO |
---|
552 | |
---|
553 | RETURN |
---|
554 | |
---|
555 | END SUBROUTINE init_domain_rk |
---|
556 | |
---|
557 | SUBROUTINE init_module_initialize |
---|
558 | END SUBROUTINE init_module_initialize |
---|
559 | |
---|
560 | !--------------------------------------------------------------------- |
---|
561 | |
---|
562 | ! test driver for get_sounding |
---|
563 | ! |
---|
564 | ! implicit none |
---|
565 | ! integer n |
---|
566 | ! parameter(n = 1000) |
---|
567 | ! real zk(n),p(n),theta(n),rho(n),u(n),v(n),qv(n),pd(n) |
---|
568 | ! logical dry |
---|
569 | ! integer nl,k |
---|
570 | ! |
---|
571 | ! dry = .false. |
---|
572 | ! dry = .true. |
---|
573 | ! call get_sounding( zk, p, pd, theta, rho, u, v, qv, dry, n, nl ) |
---|
574 | ! write(6,*) ' input levels ',nl |
---|
575 | ! write(6,*) ' sounding ' |
---|
576 | ! write(6,*) ' k height(m) press (Pa) pd(Pa) theta (K) den(kg/m^3) u(m/s) v(m/s) qv(g/g) ' |
---|
577 | ! do k=1,nl |
---|
578 | ! 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) |
---|
579 | ! enddo |
---|
580 | ! end |
---|
581 | ! |
---|
582 | !--------------------------------------------------------------------------- |
---|
583 | |
---|
584 | subroutine get_sounding( zk, p, p_dry, theta, rho, & |
---|
585 | u, v, qv, dry, nl_max, nl_in ) |
---|
586 | implicit none |
---|
587 | |
---|
588 | integer nl_max, nl_in |
---|
589 | real zk(nl_max), p(nl_max), theta(nl_max), rho(nl_max), & |
---|
590 | u(nl_max), v(nl_max), qv(nl_max), p_dry(nl_max) |
---|
591 | logical dry |
---|
592 | |
---|
593 | integer n |
---|
594 | parameter(n=1000) |
---|
595 | logical debug |
---|
596 | parameter( debug = .true.) |
---|
597 | |
---|
598 | ! input sounding data |
---|
599 | |
---|
600 | real p_surf, th_surf, qv_surf |
---|
601 | real pi_surf, pi(n) |
---|
602 | real h_input(n), th_input(n), qv_input(n), u_input(n), v_input(n) |
---|
603 | |
---|
604 | ! diagnostics |
---|
605 | |
---|
606 | real rho_surf, p_input(n), rho_input(n) |
---|
607 | real pm_input(n) ! this are for full moist sounding |
---|
608 | |
---|
609 | ! local data |
---|
610 | |
---|
611 | real r |
---|
612 | parameter (r = r_d) |
---|
613 | integer k, it, nl |
---|
614 | real qvf, qvf1, dz |
---|
615 | |
---|
616 | ! first, read the sounding |
---|
617 | |
---|
618 | call read_sounding( p_surf, th_surf, qv_surf, & |
---|
619 | h_input, th_input, qv_input, u_input, v_input,n, nl, debug ) |
---|
620 | |
---|
621 | if(dry) then |
---|
622 | do k=1,nl |
---|
623 | qv_input(k) = 0. |
---|
624 | enddo |
---|
625 | endif |
---|
626 | |
---|
627 | if(debug) write(6,*) ' number of input levels = ',nl |
---|
628 | |
---|
629 | nl_in = nl |
---|
630 | if(nl_in .gt. nl_max ) then |
---|
631 | write(6,*) ' too many levels for input arrays ',nl_in,nl_max |
---|
632 | call wrf_error_fatal ( ' too many levels for input arrays ' ) |
---|
633 | end if |
---|
634 | |
---|
635 | ! compute diagnostics, |
---|
636 | ! first, convert qv(g/kg) to qv(g/g) |
---|
637 | |
---|
638 | do k=1,nl |
---|
639 | qv_input(k) = 0.001*qv_input(k) |
---|
640 | enddo |
---|
641 | |
---|
642 | p_surf = 100.*p_surf ! convert to pascals |
---|
643 | qvf = 1. + rvovrd*qv_input(1) |
---|
644 | rho_surf = 1./((r/p1000mb)*th_surf*qvf*((p_surf/p1000mb)**cvpm)) |
---|
645 | pi_surf = (p_surf/p1000mb)**(r/cp) |
---|
646 | |
---|
647 | if(debug) then |
---|
648 | write(6,*) ' surface density is ',rho_surf |
---|
649 | write(6,*) ' surface pi is ',pi_surf |
---|
650 | end if |
---|
651 | |
---|
652 | |
---|
653 | ! integrate moist sounding hydrostatically, starting from the |
---|
654 | ! specified surface pressure |
---|
655 | ! -> first, integrate from surface to lowest level |
---|
656 | |
---|
657 | qvf = 1. + rvovrd*qv_input(1) |
---|
658 | qvf1 = 1. + qv_input(1) |
---|
659 | rho_input(1) = rho_surf |
---|
660 | dz = h_input(1) |
---|
661 | do it=1,10 |
---|
662 | pm_input(1) = p_surf & |
---|
663 | - 0.5*dz*(rho_surf+rho_input(1))*g*qvf1 |
---|
664 | rho_input(1) = 1./((r/p1000mb)*th_input(1)*qvf*((pm_input(1)/p1000mb)**cvpm)) |
---|
665 | enddo |
---|
666 | |
---|
667 | ! integrate up the column |
---|
668 | |
---|
669 | do k=2,nl |
---|
670 | rho_input(k) = rho_input(k-1) |
---|
671 | dz = h_input(k)-h_input(k-1) |
---|
672 | qvf1 = 0.5*(2.+(qv_input(k-1)+qv_input(k))) |
---|
673 | qvf = 1. + rvovrd*qv_input(k) ! qv is in g/kg here |
---|
674 | |
---|
675 | do it=1,10 |
---|
676 | pm_input(k) = pm_input(k-1) & |
---|
677 | - 0.5*dz*(rho_input(k)+rho_input(k-1))*g*qvf1 |
---|
678 | rho_input(k) = 1./((r/p1000mb)*th_input(k)*qvf*((pm_input(k)/p1000mb)**cvpm)) |
---|
679 | enddo |
---|
680 | enddo |
---|
681 | |
---|
682 | ! we have the moist sounding |
---|
683 | |
---|
684 | ! next, compute the dry sounding using p at the highest level from the |
---|
685 | ! moist sounding and integrating down. |
---|
686 | |
---|
687 | p_input(nl) = pm_input(nl) |
---|
688 | |
---|
689 | do k=nl-1,1,-1 |
---|
690 | dz = h_input(k+1)-h_input(k) |
---|
691 | p_input(k) = p_input(k+1) + 0.5*dz*(rho_input(k)+rho_input(k+1))*g |
---|
692 | enddo |
---|
693 | |
---|
694 | |
---|
695 | do k=1,nl |
---|
696 | |
---|
697 | zk(k) = h_input(k) |
---|
698 | p(k) = pm_input(k) |
---|
699 | p_dry(k) = p_input(k) |
---|
700 | theta(k) = th_input(k) |
---|
701 | rho(k) = rho_input(k) |
---|
702 | u(k) = u_input(k) |
---|
703 | v(k) = v_input(k) |
---|
704 | qv(k) = qv_input(k) |
---|
705 | |
---|
706 | enddo |
---|
707 | |
---|
708 | if(debug) then |
---|
709 | write(6,*) ' sounding ' |
---|
710 | write(6,*) ' k height(m) press (Pa) pd(Pa) theta (K) den(kg/m^3) u(m/s) v(m/s) qv(g/g) ' |
---|
711 | do k=1,nl |
---|
712 | 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) |
---|
713 | enddo |
---|
714 | |
---|
715 | end if |
---|
716 | |
---|
717 | end subroutine get_sounding |
---|
718 | |
---|
719 | !------------------------------------------------------- |
---|
720 | |
---|
721 | subroutine read_sounding( ps,ts,qvs,h,th,qv,u,v,n,nl,debug ) |
---|
722 | implicit none |
---|
723 | integer n,nl |
---|
724 | real ps,ts,qvs,h(n),th(n),qv(n),u(n),v(n) |
---|
725 | logical end_of_file |
---|
726 | logical debug |
---|
727 | |
---|
728 | integer k |
---|
729 | |
---|
730 | open(unit=10,file='input_sounding',form='formatted',status='old') |
---|
731 | rewind(10) |
---|
732 | read(10,*) ps, ts, qvs |
---|
733 | if(debug) then |
---|
734 | write(6,*) ' input sounding surface parameters ' |
---|
735 | write(6,*) ' surface pressure (mb) ',ps |
---|
736 | write(6,*) ' surface pot. temp (K) ',ts |
---|
737 | write(6,*) ' surface mixing ratio (g/kg) ',qvs |
---|
738 | end if |
---|
739 | |
---|
740 | end_of_file = .false. |
---|
741 | k = 0 |
---|
742 | |
---|
743 | do while (.not. end_of_file) |
---|
744 | |
---|
745 | read(10,*,end=100) h(k+1), th(k+1), qv(k+1), u(k+1), v(k+1) |
---|
746 | k = k+1 |
---|
747 | if(debug) write(6,'(1x,i3,5(1x,e10.3))') k, h(k), th(k), qv(k), u(k), v(k) |
---|
748 | go to 110 |
---|
749 | 100 end_of_file = .true. |
---|
750 | 110 continue |
---|
751 | enddo |
---|
752 | |
---|
753 | nl = k |
---|
754 | |
---|
755 | close(unit=10,status = 'keep') |
---|
756 | |
---|
757 | end subroutine read_sounding |
---|
758 | |
---|
759 | END MODULE module_initialize_ideal |
---|