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 :: idum1, idum2 |
---|
51 | |
---|
52 | CALL set_scalar_indices_from_config ( head_grid%id , idum1, idum2 ) |
---|
53 | |
---|
54 | CALL init_domain_rk( grid & |
---|
55 | ! |
---|
56 | #include <actual_new_args.inc> |
---|
57 | ! |
---|
58 | ) |
---|
59 | |
---|
60 | END SUBROUTINE init_domain |
---|
61 | |
---|
62 | !------------------------------------------------------------------- |
---|
63 | |
---|
64 | SUBROUTINE init_domain_rk ( grid & |
---|
65 | ! |
---|
66 | # include <dummy_new_args.inc> |
---|
67 | ! |
---|
68 | ) |
---|
69 | IMPLICIT NONE |
---|
70 | |
---|
71 | ! Input data. |
---|
72 | TYPE (domain), POINTER :: grid |
---|
73 | |
---|
74 | # include <dummy_new_decl.inc> |
---|
75 | |
---|
76 | TYPE (grid_config_rec_type) :: config_flags |
---|
77 | |
---|
78 | ! Local data |
---|
79 | INTEGER :: & |
---|
80 | ids, ide, jds, jde, kds, kde, & |
---|
81 | ims, ime, jms, jme, kms, kme, & |
---|
82 | its, ite, jts, jte, kts, kte, & |
---|
83 | i, j, k |
---|
84 | |
---|
85 | ! Local data |
---|
86 | |
---|
87 | INTEGER, PARAMETER :: nl_max = 1000 |
---|
88 | REAL, DIMENSION(nl_max) :: zk, p_in, theta, rho, u, v, qv, pd_in |
---|
89 | REAL*8, DIMENSION(nl_max) :: pd_in8 |
---|
90 | INTEGER :: nl_in |
---|
91 | |
---|
92 | INTEGER :: icm,jcm, ii, im1, jj, jm1, loop, error, fid, nxc, nyc |
---|
93 | REAL :: u_mean,v_mean, f0, p_surf, p_level, qvf, z_at_v, z_at_u |
---|
94 | REAL*8 :: p_level8 |
---|
95 | REAL :: xrad, yrad, zrad, rad, delt, cof1, cof2 |
---|
96 | ! REAL, EXTERNAL :: interp_0 |
---|
97 | REAL :: hm |
---|
98 | REAL :: pi |
---|
99 | |
---|
100 | ! stuff from original initialization that has been dropped from the Registry |
---|
101 | REAL :: vnu, xnu, xnus, dinit0, cbh, p0_temp, t0_temp, zd, zt |
---|
102 | REAL :: qvf1, qvf2, pd_surf |
---|
103 | INTEGER :: it |
---|
104 | real :: thtmp, ptmp, temp(3) |
---|
105 | |
---|
106 | LOGICAL :: moisture_init |
---|
107 | LOGICAL :: stretch_grid, dry_sounding |
---|
108 | |
---|
109 | INTEGER :: xs , xe , ys , ye |
---|
110 | REAL :: mtn_ht |
---|
111 | LOGICAL, EXTERNAL :: wrf_dm_on_monitor |
---|
112 | ! For LES, add randx |
---|
113 | real :: randx |
---|
114 | |
---|
115 | !!MARS |
---|
116 | REAL :: lon_input, lat_input, alt_input, tsurf_input |
---|
117 | ! for mode 3 |
---|
118 | REAL, DIMENSION(nl_max) :: profdustq,profdustn |
---|
119 | REAL, DIMENSION(nl_max) :: prescribed_sw,prescribed_lw,prescribed_dyn |
---|
120 | REAL, DIMENSION(nl_max) :: hrsw,hrlw,hrdyn |
---|
121 | REAL, DIMENSION(nl_max) :: lsf_dt,lsf_dq,lsfdt,lsfdq |
---|
122 | REAL, DIMENSION(nl_max) :: venus_hrdyn |
---|
123 | REAL, DIMENSION(nl_max) :: altitude |
---|
124 | REAL*8, DIMENSION(nl_max) :: trac |
---|
125 | !!MARS |
---|
126 | |
---|
127 | REAL :: pfu, pfd, phm |
---|
128 | INTEGER :: hypsometric_opt = 1 ! classic |
---|
129 | !INTEGER :: hypsometric_opt = 2 ! Wee et al. 2012 correction |
---|
130 | |
---|
131 | LOGICAL :: logp = .true. ! use logp to interpolate (and not p) |
---|
132 | |
---|
133 | #ifdef DM_PARALLEL |
---|
134 | # include <data_calls.inc> |
---|
135 | #endif |
---|
136 | |
---|
137 | |
---|
138 | SELECT CASE ( model_data_order ) |
---|
139 | CASE ( DATA_ORDER_ZXY ) |
---|
140 | kds = grid%sd31 ; kde = grid%ed31 ; |
---|
141 | ids = grid%sd32 ; ide = grid%ed32 ; |
---|
142 | jds = grid%sd33 ; jde = grid%ed33 ; |
---|
143 | |
---|
144 | kms = grid%sm31 ; kme = grid%em31 ; |
---|
145 | ims = grid%sm32 ; ime = grid%em32 ; |
---|
146 | jms = grid%sm33 ; jme = grid%em33 ; |
---|
147 | |
---|
148 | kts = grid%sp31 ; kte = grid%ep31 ; ! note that tile is entire patch |
---|
149 | its = grid%sp32 ; ite = grid%ep32 ; ! note that tile is entire patch |
---|
150 | jts = grid%sp33 ; jte = grid%ep33 ; ! note that tile is entire patch |
---|
151 | CASE ( DATA_ORDER_XYZ ) |
---|
152 | ids = grid%sd31 ; ide = grid%ed31 ; |
---|
153 | jds = grid%sd32 ; jde = grid%ed32 ; |
---|
154 | kds = grid%sd33 ; kde = grid%ed33 ; |
---|
155 | |
---|
156 | ims = grid%sm31 ; ime = grid%em31 ; |
---|
157 | jms = grid%sm32 ; jme = grid%em32 ; |
---|
158 | kms = grid%sm33 ; kme = grid%em33 ; |
---|
159 | |
---|
160 | its = grid%sp31 ; ite = grid%ep31 ; ! note that tile is entire patch |
---|
161 | jts = grid%sp32 ; jte = grid%ep32 ; ! note that tile is entire patch |
---|
162 | kts = grid%sp33 ; kte = grid%ep33 ; ! note that tile is entire patch |
---|
163 | CASE ( DATA_ORDER_XZY ) |
---|
164 | ids = grid%sd31 ; ide = grid%ed31 ; |
---|
165 | kds = grid%sd32 ; kde = grid%ed32 ; |
---|
166 | jds = grid%sd33 ; jde = grid%ed33 ; |
---|
167 | |
---|
168 | ims = grid%sm31 ; ime = grid%em31 ; |
---|
169 | kms = grid%sm32 ; kme = grid%em32 ; |
---|
170 | jms = grid%sm33 ; jme = grid%em33 ; |
---|
171 | |
---|
172 | its = grid%sp31 ; ite = grid%ep31 ; ! note that tile is entire patch |
---|
173 | kts = grid%sp32 ; kte = grid%ep32 ; ! note that tile is entire patch |
---|
174 | jts = grid%sp33 ; jte = grid%ep33 ; ! note that tile is entire patch |
---|
175 | |
---|
176 | END SELECT |
---|
177 | |
---|
178 | IF (planet == "mars" .or. planet == "titan") THEN |
---|
179 | stretch_grid = .false. |
---|
180 | !! FOR LES, set stretch to false |
---|
181 | ELSE |
---|
182 | stretch_grid = .true. !! VENUS |
---|
183 | ENDIF |
---|
184 | delt = 3. |
---|
185 | ! z_scale = .50 |
---|
186 | ! z_scale = .10 |
---|
187 | ! z_scale = .25 |
---|
188 | ! z_scale = .15 |
---|
189 | pi = 2.*asin(1.0) |
---|
190 | write(6,*) ' pi is ',pi |
---|
191 | nxc = (ide-ids)/2 |
---|
192 | nyc = (jde-jds)/2 |
---|
193 | |
---|
194 | CALL model_to_grid_config_rec ( grid%id , model_config_rec , config_flags ) |
---|
195 | |
---|
196 | ! here we check to see if the boundary conditions are set properly |
---|
197 | |
---|
198 | CALL boundary_condition_check( config_flags, bdyzone, error, grid%id ) |
---|
199 | |
---|
200 | moisture_init = .true. |
---|
201 | |
---|
202 | grid%itimestep=0 |
---|
203 | |
---|
204 | #ifdef DM_PARALLEL |
---|
205 | CALL wrf_dm_bcast_bytes( icm , IWORDSIZE ) |
---|
206 | CALL wrf_dm_bcast_bytes( jcm , IWORDSIZE ) |
---|
207 | #endif |
---|
208 | |
---|
209 | CALL nl_set_mminlu(1, ' ') |
---|
210 | CALL nl_set_iswater(1,0) |
---|
211 | CALL nl_set_cen_lat(1,40.) |
---|
212 | CALL nl_set_cen_lon(1,-105.) |
---|
213 | CALL nl_set_truelat1(1,0.) |
---|
214 | CALL nl_set_truelat2(1,0.) |
---|
215 | CALL nl_set_moad_cen_lat (1,0.) |
---|
216 | CALL nl_set_stand_lon (1,0.) |
---|
217 | CALL nl_set_map_proj(1,0) |
---|
218 | |
---|
219 | |
---|
220 | ! here we initialize data we currently is not initialized |
---|
221 | ! in the input data |
---|
222 | |
---|
223 | DO j = jts, jte |
---|
224 | DO i = its, ite |
---|
225 | grid%msftx(i,j) = 1. |
---|
226 | grid%msfty(i,j) = 1. |
---|
227 | grid%msfux(i,j) = 1. |
---|
228 | grid%msfuy(i,j) = 1. |
---|
229 | grid%msfvx(i,j) = 1. |
---|
230 | grid%msfvx_inv(i,j)= 1. |
---|
231 | grid%msfvy(i,j) = 1. |
---|
232 | grid%sina(i,j) = 0. |
---|
233 | grid%cosa(i,j) = 1. |
---|
234 | grid%e(i,j) = 0. |
---|
235 | ! for LES, include Coriolis force |
---|
236 | grid%f(i,j) = 0. !!MARS MARS 1.e-4 |
---|
237 | !! grid%f(i,j) = 2*EOMEG*SIN(grid%xlat(i,j)*degrad) |
---|
238 | END DO |
---|
239 | END DO |
---|
240 | |
---|
241 | DO j = jts, jte |
---|
242 | DO k = kts, kte |
---|
243 | DO i = its, ite |
---|
244 | grid%ww(i,k,j) = 0. |
---|
245 | END DO |
---|
246 | END DO |
---|
247 | END DO |
---|
248 | |
---|
249 | grid%step_number = 0 |
---|
250 | |
---|
251 | ! set up the grid |
---|
252 | |
---|
253 | IF (stretch_grid) THEN ! exponential stretch for eta (nearly constant dz) |
---|
254 | DO k=1, kde |
---|
255 | grid%znw(k) = (exp(-(k-1)/float(kde-1)/z_scale) - exp(-1./z_scale))/ & |
---|
256 | (1.-exp(-1./z_scale)) |
---|
257 | ENDDO |
---|
258 | ELSE |
---|
259 | !!MARS |
---|
260 | !!MARS |
---|
261 | IF (planet .ne. "mars" .and. planet .ne. "titan") THEN |
---|
262 | open(unit=12,file='levels',form='formatted',status='old') |
---|
263 | rewind(12) |
---|
264 | DO k=1, kde |
---|
265 | read(12,*) grid%znw(k) |
---|
266 | write(6,*) 'read level ', k,grid%znw(k) |
---|
267 | ENDDO |
---|
268 | close(12) |
---|
269 | ENDIF |
---|
270 | !!MARS |
---|
271 | !!MARS |
---|
272 | ! !DO k=1, kde |
---|
273 | ! ! grid%znw(k) = 1. - float(k-1)/float(kde-1) |
---|
274 | ! !ENDDO |
---|
275 | ENDIF |
---|
276 | |
---|
277 | !! SPECIFIC FOR LES PBL MARS |
---|
278 | IF (planet == "mars" .or. planet == "titan") THEN |
---|
279 | !!!MARS |
---|
280 | grid%znw(1)=1.000 |
---|
281 | grid%znw(2)=0.9995 !5m |
---|
282 | grid%znw(3)=0.9980 !20m |
---|
283 | grid%znw(4)=0.9950 !55m |
---|
284 | DO k=5, kde |
---|
285 | grid%znw(k) = grid%znw(4) * ( 1. - float(k-4)/float(kde-4) ) |
---|
286 | ENDDO |
---|
287 | ENDIF |
---|
288 | |
---|
289 | open(unit=12,file='levels',form='formatted',status='old') |
---|
290 | rewind(12) |
---|
291 | DO k=1, kde |
---|
292 | write(12,*) grid%znw(k) |
---|
293 | write(6,*) 'for update_inputs_physiq_mod (e.g. Titan, generic)',k,grid%znw(k) |
---|
294 | ENDDO |
---|
295 | close(12) |
---|
296 | |
---|
297 | |
---|
298 | |
---|
299 | DO k=1, kde-1 |
---|
300 | grid%dnw(k) = grid%znw(k+1) - grid%znw(k) |
---|
301 | grid%rdnw(k) = 1./grid%dnw(k) |
---|
302 | grid%znu(k) = 0.5*(grid%znw(k+1)+grid%znw(k)) |
---|
303 | ENDDO |
---|
304 | DO k=2, kde-1 |
---|
305 | grid%dn(k) = 0.5*(grid%dnw(k)+grid%dnw(k-1)) |
---|
306 | grid%rdn(k) = 1./grid%dn(k) |
---|
307 | grid%fnp(k) = .5* grid%dnw(k )/grid%dn(k) |
---|
308 | grid%fnm(k) = .5* grid%dnw(k-1)/grid%dn(k) |
---|
309 | ENDDO |
---|
310 | |
---|
311 | cof1 = (2.*grid%dn(2)+grid%dn(3))/(grid%dn(2)+grid%dn(3))*grid%dnw(1)/grid%dn(2) |
---|
312 | cof2 = grid%dn(2) /(grid%dn(2)+grid%dn(3))*grid%dnw(1)/grid%dn(3) |
---|
313 | grid%cf1 = grid%fnp(2) + cof1 |
---|
314 | grid%cf2 = grid%fnm(2) - cof1 - cof2 |
---|
315 | grid%cf3 = cof2 |
---|
316 | |
---|
317 | grid%cfn = (.5*grid%dnw(kde-1)+grid%dn(kde-1))/grid%dn(kde-1) |
---|
318 | grid%cfn1 = -.5*grid%dnw(kde-1)/grid%dn(kde-1) |
---|
319 | grid%rdx = 1./config_flags%dx |
---|
320 | grid%rdy = 1./config_flags%dy |
---|
321 | |
---|
322 | ! get the sounding from the ascii sounding file, first get dry sounding and |
---|
323 | ! calculate base state |
---|
324 | |
---|
325 | dry_sounding = .true. |
---|
326 | IF ( wrf_dm_on_monitor() ) THEN |
---|
327 | write(6,*) ' getting dry sounding for base state ' |
---|
328 | |
---|
329 | CALL get_sounding( zk, p_in, pd_in, theta, rho, u, v, qv, dry_sounding, nl_max, nl_in ) |
---|
330 | ENDIF |
---|
331 | CALL wrf_dm_bcast_real( zk , nl_max ) |
---|
332 | CALL wrf_dm_bcast_real( p_in , nl_max ) |
---|
333 | CALL wrf_dm_bcast_real( pd_in , nl_max ) |
---|
334 | CALL wrf_dm_bcast_real( theta , nl_max ) |
---|
335 | CALL wrf_dm_bcast_real( rho , nl_max ) |
---|
336 | CALL wrf_dm_bcast_real( u , nl_max ) |
---|
337 | CALL wrf_dm_bcast_real( v , nl_max ) |
---|
338 | CALL wrf_dm_bcast_real( qv , nl_max ) |
---|
339 | CALL wrf_dm_bcast_integer ( nl_in , 1 ) |
---|
340 | |
---|
341 | write(6,*) ' returned from reading sounding, nl_in is ',nl_in |
---|
342 | |
---|
343 | !!MARS |
---|
344 | !!MARS |
---|
345 | open(unit=14,file='input_coord',form='formatted',status='old') |
---|
346 | rewind(14) |
---|
347 | read(14,*) lon_input |
---|
348 | read(14,*) lat_input |
---|
349 | close(14) |
---|
350 | write(6,*) ' lon is ',lon_input |
---|
351 | write(6,*) ' lat is ',lat_input |
---|
352 | !!MARS |
---|
353 | !!MARS |
---|
354 | |
---|
355 | !!MARS |
---|
356 | !!MARS |
---|
357 | open(unit=18,file='input_more',form='formatted',status='old') |
---|
358 | rewind(18) |
---|
359 | read(18,*) alt_input, tsurf_input |
---|
360 | close(18) |
---|
361 | write(6,*) ' alt is ',alt_input |
---|
362 | write(6,*) ' tsurf is ',tsurf_input |
---|
363 | !!MARS |
---|
364 | !!MARS |
---|
365 | |
---|
366 | ! find ptop for the desired ztop (ztop is input from the namelist), |
---|
367 | ! and find surface pressure |
---|
368 | |
---|
369 | write(6,*) ' ztop above ground is ',config_flags%ztop |
---|
370 | write(6,*) ' real ztop is ',config_flags%ztop + alt_input |
---|
371 | grid%p_top = interp_0( p_in, zk, config_flags%ztop + alt_input, nl_in ) |
---|
372 | |
---|
373 | DO j=jts,jte |
---|
374 | DO i=its,ite |
---|
375 | !!MARS |
---|
376 | grid%ht(i,j) = alt_input |
---|
377 | grid%m_tsurf(i,j) = tsurf_input |
---|
378 | !!MARS |
---|
379 | grid%xlat(i,j) = lat_input !+ float(j)*config_flags%dy/59000. |
---|
380 | grid%xlong(i,j) = lon_input !+ float(i)*config_flags%dx/59000. |
---|
381 | grid%m_emiss(i,j)=0.95 |
---|
382 | grid%m_co2ice(i,j)=0. |
---|
383 | grid%m_h2oice(i,j)=0. |
---|
384 | !! >> Used for restarts only: |
---|
385 | grid%m_q2(i,:,j)=0. |
---|
386 | grid%m_fluxrad(i,j)=0. |
---|
387 | grid%m_wstar(i,j)=0. |
---|
388 | !! << |
---|
389 | grid%slpx(i,j) = 0. |
---|
390 | grid%slpy(i,j) = 0. |
---|
391 | DO k=1,config_flags%num_soil_layers |
---|
392 | grid%m_tsoil(i,k,j) = 0. |
---|
393 | ENDDO |
---|
394 | grid%m_gw(i,1,j) = 0. |
---|
395 | grid%m_gw(i,2,j) = 0. |
---|
396 | grid%m_gw(i,3,j) = 0. |
---|
397 | grid%m_gw(i,4,j) = 0. |
---|
398 | grid%m_gw(i,5,j) = 0. |
---|
399 | !!MARS |
---|
400 | ENDDO |
---|
401 | ENDDO |
---|
402 | |
---|
403 | xs=ide/2 -3 |
---|
404 | xs=ids -3 |
---|
405 | xe=xs + 6 |
---|
406 | ys=jde/2 -3 |
---|
407 | ye=ys + 6 |
---|
408 | mtn_ht = 500 |
---|
409 | #ifdef MTN |
---|
410 | DO j=max(ys,jds),min(ye,jde-1) |
---|
411 | DO i=max(xs,ids),min(xe,ide-1) |
---|
412 | grid%ht(i,j) = alt_input + mtn_ht * 0.25 * & |
---|
413 | ( 1. + COS ( 2*pi/(xe-xs) * ( i-xs ) + pi ) ) * & |
---|
414 | ( 1. + COS ( 2*pi/(ye-ys) * ( j-ys ) + pi ) ) |
---|
415 | ENDDO |
---|
416 | ENDDO |
---|
417 | #endif |
---|
418 | #ifdef EW_RIDGE |
---|
419 | DO j=max(ys,jds),min(ye,jde-1) |
---|
420 | DO i=ids,ide |
---|
421 | grid%ht(i,j) = mtn_ht * 0.50 * & |
---|
422 | ( 1. + COS ( 2*pi/(ye-ys) * ( j-ys ) + pi ) ) |
---|
423 | ENDDO |
---|
424 | ENDDO |
---|
425 | #endif |
---|
426 | #ifdef NS_RIDGE |
---|
427 | DO j=jds,jde |
---|
428 | DO i=max(xs,ids),min(xe,ide-1) |
---|
429 | grid%ht(i,j) = mtn_ht * 0.50 * & |
---|
430 | ( 1. + COS ( 2*pi/(xe-xs) * ( i-xs ) + pi ) ) |
---|
431 | ENDDO |
---|
432 | ENDDO |
---|
433 | #endif |
---|
434 | DO j=jts,jte |
---|
435 | DO i=its,ite |
---|
436 | grid%phb(i,1,j) = g * grid%ht(i,j) |
---|
437 | grid%ph0(i,1,j) = g * grid%ht(i,j) |
---|
438 | ENDDO |
---|
439 | ENDDO |
---|
440 | |
---|
441 | IF (.not.logp) THEN |
---|
442 | write(6,*) 'interpolate in p' |
---|
443 | ELSE |
---|
444 | write(6,*) 'interpolate in logp' |
---|
445 | ENDIF |
---|
446 | |
---|
447 | DO J = jts, jte |
---|
448 | DO I = its, ite |
---|
449 | |
---|
450 | p_surf = interp_0( p_in, zk, grid%phb(i,1,j)/g, nl_in ) |
---|
451 | grid%mub(i,j) = p_surf-grid%p_top |
---|
452 | |
---|
453 | ! this is dry hydrostatic sounding (base state), so given grid%p (coordinate), |
---|
454 | ! interp theta (from interp) and compute 1/rho from eqn. of state |
---|
455 | |
---|
456 | DO K = 1, kte-1 |
---|
457 | p_level = grid%znu(k)*(p_surf - grid%p_top) + grid%p_top |
---|
458 | grid%pb(i,k,j) = p_level |
---|
459 | IF (.not.logp) THEN |
---|
460 | grid%t_init(i,k,j) = interp_0( theta, p_in, p_level, nl_in ) - t0 |
---|
461 | ELSE |
---|
462 | grid%t_init(i,k,j) = interp_0_log( theta, p_in, p_level, nl_in ) - t0 |
---|
463 | ENDIF |
---|
464 | grid%alb(i,k,j) = (r_d/p1000mb)*(grid%t_init(i,k,j)+t0)*(grid%pb(i,k,j)/p1000mb)**cvpm |
---|
465 | ENDDO |
---|
466 | |
---|
467 | ! calc hydrostatic balance (alternatively we could interp the geopotential from the |
---|
468 | ! sounding, but this assures that the base state is in exact hydrostatic balance with |
---|
469 | ! respect to the model eqns. |
---|
470 | |
---|
471 | IF (hypsometric_opt == 1) THEN |
---|
472 | DO k = 2,kte |
---|
473 | 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) |
---|
474 | ENDDO |
---|
475 | ELSE IF (hypsometric_opt == 2) THEN |
---|
476 | DO k = 2,kte |
---|
477 | pfu = grid%mub(i,j)*grid%znw(k) + grid%p_top |
---|
478 | pfd = grid%mub(i,j)*grid%znw(k-1) + grid%p_top |
---|
479 | phm = grid%mub(i,j)*grid%znu(k-1) + grid%p_top |
---|
480 | grid%phb(i,k,j) = grid%phb(i,k-1,j) + grid%alb(i,k-1,j)*phm*LOG(pfd/pfu) |
---|
481 | END DO |
---|
482 | END IF |
---|
483 | |
---|
484 | |
---|
485 | ENDDO |
---|
486 | ENDDO |
---|
487 | IF ( wrf_dm_on_monitor() ) THEN |
---|
488 | write(6,*) ' ptop is ',grid%p_top |
---|
489 | write(6,*) ' base state grid%mub(1,1), p_surf is ',grid%mub(1,1),grid%mub(1,1)+grid%p_top |
---|
490 | ENDIF |
---|
491 | |
---|
492 | ! calculate full state for each column - this includes moisture. |
---|
493 | |
---|
494 | write(6,*) ' getting moist sounding for full state ' |
---|
495 | dry_sounding = .false. |
---|
496 | CALL get_sounding( zk, p_in, pd_in, theta, rho, u, v, qv, dry_sounding, nl_max, nl_in ) |
---|
497 | |
---|
498 | DO J = jts, min(jde-1,jte) |
---|
499 | DO I = its, min(ide-1,ite) |
---|
500 | |
---|
501 | ! At this point grid%p_top is already set. find the DRY mass in the column |
---|
502 | ! by interpolating the DRY pressure. |
---|
503 | |
---|
504 | pd_surf = interp_0( pd_in, zk, grid%phb(i,1,j)/g, nl_in ) |
---|
505 | |
---|
506 | ! compute the perturbation mass and the full mass |
---|
507 | |
---|
508 | grid%mu_1(i,j) = pd_surf-grid%p_top - grid%mub(i,j) |
---|
509 | grid%mu_2(i,j) = grid%mu_1(i,j) |
---|
510 | grid%mu0(i,j) = grid%mu_1(i,j) + grid%mub(i,j) |
---|
511 | |
---|
512 | ! given the dry pressure and coordinate system, interp the potential |
---|
513 | ! temperature and qv |
---|
514 | |
---|
515 | do k=1,kde-1 |
---|
516 | |
---|
517 | p_level = grid%znu(k)*(pd_surf - grid%p_top) + grid%p_top |
---|
518 | IF (.not.logp) THEN |
---|
519 | moist(i,k,j,P_QV) = interp_0( qv, pd_in, p_level, nl_in ) |
---|
520 | grid%t_1(i,k,j) = interp_0( theta, pd_in, p_level, nl_in ) - t0 |
---|
521 | ELSE |
---|
522 | moist(i,k,j,P_QV) = interp_0_log( qv, pd_in, p_level, nl_in ) |
---|
523 | grid%t_1(i,k,j) = interp_0_log( theta, pd_in, p_level, nl_in ) - t0 |
---|
524 | ENDIF |
---|
525 | grid%t_2(i,k,j) = grid%t_1(i,k,j) |
---|
526 | |
---|
527 | |
---|
528 | enddo |
---|
529 | |
---|
530 | ! integrate the hydrostatic equation (from the RHS of the bigstep |
---|
531 | ! vertical momentum equation) down from the top to get grid%p. |
---|
532 | ! first from the top of the model to the top pressure |
---|
533 | |
---|
534 | k = kte-1 ! top level |
---|
535 | |
---|
536 | qvf1 = 0.5*(moist(i,k,j,P_QV)+moist(i,k,j,P_QV)) |
---|
537 | qvf2 = 1./(1.+qvf1) |
---|
538 | qvf1 = qvf1*qvf2 |
---|
539 | |
---|
540 | ! grid%p(i,k,j) = - 0.5*grid%mu_1(i,j)/grid%rdnw(k) |
---|
541 | grid%p(i,k,j) = - 0.5*(grid%mu_1(i,j)+qvf1*grid%mub(i,j))/grid%rdnw(k)/qvf2 |
---|
542 | qvf = 1. + rvovrd*moist(i,k,j,P_QV) |
---|
543 | grid%alt(i,k,j) = (r_d/p1000mb)*(grid%t_1(i,k,j)+t0)*qvf* & |
---|
544 | (((grid%p(i,k,j)+grid%pb(i,k,j))/p1000mb)**cvpm) |
---|
545 | grid%al(i,k,j) = grid%alt(i,k,j) - grid%alb(i,k,j) |
---|
546 | |
---|
547 | ! down the column |
---|
548 | |
---|
549 | do k=kte-2,1,-1 |
---|
550 | qvf1 = 0.5*(moist(i,k,j,P_QV)+moist(i,k+1,j,P_QV)) |
---|
551 | qvf2 = 1./(1.+qvf1) |
---|
552 | qvf1 = qvf1*qvf2 |
---|
553 | 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) |
---|
554 | qvf = 1. + rvovrd*moist(i,k,j,P_QV) |
---|
555 | grid%alt(i,k,j) = (r_d/p1000mb)*(grid%t_1(i,k,j)+t0)*qvf* & |
---|
556 | (((grid%p(i,k,j)+grid%pb(i,k,j))/p1000mb)**cvpm) |
---|
557 | grid%al(i,k,j) = grid%alt(i,k,j) - grid%alb(i,k,j) |
---|
558 | enddo |
---|
559 | |
---|
560 | ! this is the hydrostatic equation used in the model after the |
---|
561 | ! small timesteps. In the model, grid%al (inverse density) |
---|
562 | ! is computed from the geopotential. |
---|
563 | |
---|
564 | |
---|
565 | grid%ph_1(i,1,j) = 0. |
---|
566 | IF (hypsometric_opt == 1) THEN |
---|
567 | DO k = 2,kte |
---|
568 | grid%ph_1(i,k,j) = grid%ph_1(i,k-1,j) - (1./grid%rdnw(k-1))*( & |
---|
569 | (grid%mub(i,j)+grid%mu_1(i,j))*grid%al(i,k-1,j)+ & |
---|
570 | grid%mu_1(i,j)*grid%alb(i,k-1,j) ) |
---|
571 | |
---|
572 | grid%ph_2(i,k,j) = grid%ph_1(i,k,j) |
---|
573 | grid%ph0(i,k,j) = grid%ph_1(i,k,j) + grid%phb(i,k,j) |
---|
574 | ENDDO |
---|
575 | ELSE IF (hypsometric_opt == 2) THEN |
---|
576 | |
---|
577 | ! Alternative hydrostatic eq.: dZ = -al*p*dLOG(p), where p is dry pressure. |
---|
578 | ! Note that al*p approximates Rd*T and dLOG(p) does z. |
---|
579 | ! Here T varies mostly linear with z, the first-order integration produces better result. |
---|
580 | |
---|
581 | grid%ph_2(i,1,j) = grid%phb(i,1,j) |
---|
582 | DO k = 2,kte |
---|
583 | pfu = grid%mu0(i,j)*grid%znw(k) + grid%p_top |
---|
584 | pfd = grid%mu0(i,j)*grid%znw(k-1) + grid%p_top |
---|
585 | phm = grid%mu0(i,j)*grid%znu(k-1) + grid%p_top |
---|
586 | grid%ph_2(i,k,j) = grid%ph_2(i,k-1,j) + grid%alt(i,k-1,j)*phm*LOG(pfd/pfu) |
---|
587 | END DO |
---|
588 | |
---|
589 | DO k = 1,kte |
---|
590 | grid%ph_2(i,k,j) = grid%ph_2(i,k,j) - grid%phb(i,k,j) |
---|
591 | grid%ph_1(i,k,j) = grid%ph_2(i,k,j) |
---|
592 | END DO |
---|
593 | |
---|
594 | END IF |
---|
595 | |
---|
596 | |
---|
597 | |
---|
598 | IF ( wrf_dm_on_monitor() ) THEN |
---|
599 | if((i==2) .and. (j==2)) then |
---|
600 | write(6,*) ' grid%ph_1 calc ',grid%ph_1(2,1,2),grid%ph_1(2,2,2),& |
---|
601 | grid%mu_1(2,2)+grid%mub(2,2),grid%mu_1(2,2), & |
---|
602 | grid%alb(2,1,2),grid%al(1,2,1),grid%rdnw(1) |
---|
603 | endif |
---|
604 | ENDIF |
---|
605 | |
---|
606 | ENDDO |
---|
607 | ENDDO |
---|
608 | |
---|
609 | !#if 0 |
---|
610 | |
---|
611 | ! thermal perturbation to kick off convection |
---|
612 | |
---|
613 | write(6,*) ' nxc, nyc for perturbation ',nxc,nyc |
---|
614 | write(6,*) ' delt for perturbation ',delt |
---|
615 | |
---|
616 | ! For LES, change the initial random perturbations |
---|
617 | ! For 2D test, call randx outside I-loop |
---|
618 | ! For 3D runs, call randx inside both I-J loops |
---|
619 | |
---|
620 | DO J = jts, min(jde-1,jte) |
---|
621 | ! yrad = config_flags%dy*float(j-nyc)/10000. |
---|
622 | yrad = 0. |
---|
623 | DO I = its, min(ide-1,ite) |
---|
624 | ! xrad = config_flags%dx*float(i-nxc)/10000. |
---|
625 | xrad = 0. |
---|
626 | call random_number (randx) |
---|
627 | randx = randx - 0.5 |
---|
628 | ! DO K = 1, kte-1 |
---|
629 | DO K = 1, 4 |
---|
630 | |
---|
631 | ! No bubbles for LES! |
---|
632 | ! put in preturbation theta (bubble) and recalc density. note, |
---|
633 | ! the mass in the column is not changing, so when theta changes, |
---|
634 | ! we recompute density and geopotential |
---|
635 | |
---|
636 | ! zrad = 0.5*(grid%ph_1(i,k,j)+grid%ph_1(i,k+1,j) & |
---|
637 | ! +grid%phb(i,k,j)+grid%phb(i,k+1,j))/g |
---|
638 | ! zrad = (zrad-1500.)/1500. |
---|
639 | zrad = 0. |
---|
640 | RAD=SQRT(xrad*xrad+yrad*yrad+zrad*zrad) |
---|
641 | IF(RAD <= 1.) THEN |
---|
642 | ! grid%t_1(i,k,j)=grid%t_1(i,k,j)+delt*COS(.5*PI*RAD)**2 |
---|
643 | grid%t_1(i,k,j)=grid%t_1(i,k,j)+ 0.1 *randx |
---|
644 | grid%t_2(i,k,j)=grid%t_1(i,k,j) |
---|
645 | qvf = 1. + rvovrd*moist(i,k,j,P_QV) |
---|
646 | grid%alt(i,k,j) = (r_d/p1000mb)*(grid%t_1(i,k,j)+t0)*qvf* & |
---|
647 | (((grid%p(i,k,j)+grid%pb(i,k,j))/p1000mb)**cvpm) |
---|
648 | grid%al(i,k,j) = grid%alt(i,k,j) - grid%alb(i,k,j) |
---|
649 | ENDIF |
---|
650 | ENDDO |
---|
651 | |
---|
652 | ! rebalance hydrostatically |
---|
653 | |
---|
654 | IF (hypsometric_opt == 1) THEN |
---|
655 | |
---|
656 | DO k = 2,kte |
---|
657 | grid%ph_1(i,k,j) = grid%ph_1(i,k-1,j) - (1./grid%rdnw(k-1))*( & |
---|
658 | (grid%mub(i,j)+grid%mu_1(i,j))*grid%al(i,k-1,j)+ & |
---|
659 | grid%mu_1(i,j)*grid%alb(i,k-1,j) ) |
---|
660 | |
---|
661 | grid%ph_2(i,k,j) = grid%ph_1(i,k,j) |
---|
662 | grid%ph0(i,k,j) = grid%ph_1(i,k,j) + grid%phb(i,k,j) |
---|
663 | ENDDO |
---|
664 | |
---|
665 | ELSE IF (hypsometric_opt == 2) THEN |
---|
666 | |
---|
667 | ! Alternative hydrostatic eq.: dZ = -al*p*dLOG(p), where p is dry pressure. |
---|
668 | ! Note that al*p approximates Rd*T and dLOG(p) does z. |
---|
669 | ! Here T varies mostly linear with z, the first-order integration produces better result. |
---|
670 | |
---|
671 | grid%ph_2(i,1,j) = grid%phb(i,1,j) |
---|
672 | DO k = 2,kte |
---|
673 | pfu = grid%mu0(i,j)*grid%znw(k) + grid%p_top |
---|
674 | pfd = grid%mu0(i,j)*grid%znw(k-1) + grid%p_top |
---|
675 | phm = grid%mu0(i,j)*grid%znu(k-1) + grid%p_top |
---|
676 | grid%ph_2(i,k,j) = grid%ph_2(i,k-1,j) + grid%alt(i,k-1,j)*phm*LOG(pfd/pfu) |
---|
677 | END DO |
---|
678 | |
---|
679 | DO k = 1,kte |
---|
680 | grid%ph_2(i,k,j) = grid%ph_2(i,k,j) - grid%phb(i,k,j) |
---|
681 | grid%ph_1(i,k,j) = grid%ph_2(i,k,j) |
---|
682 | END DO |
---|
683 | |
---|
684 | END IF |
---|
685 | |
---|
686 | |
---|
687 | ENDDO |
---|
688 | ENDDO |
---|
689 | |
---|
690 | !#endif |
---|
691 | |
---|
692 | IF ( wrf_dm_on_monitor() ) THEN |
---|
693 | write(6,*) ' grid%mu_1 from comp ', grid%mu_1(1,1) |
---|
694 | write(6,*) ' full state sounding from comp, ph/g, grid%p, grid%al, grid%t_1, qv ' |
---|
695 | do k=1,kde-1 |
---|
696 | write(6,'(i3,1x,5(1x,1pe10.3))') k, (grid%ph_1(1,k,1)+grid%phb(1,k,1))/g, & |
---|
697 | grid%p(1,k,1)+grid%pb(1,k,1), grid%alt(1,k,1), & |
---|
698 | grid%t_1(1,k,1)+t0, moist(1,k,1,P_QV) |
---|
699 | enddo |
---|
700 | |
---|
701 | write(6,*) ' pert state sounding from comp, grid%ph_1, pp, alp, grid%t_1, qv ' |
---|
702 | do k=1,kde-1 |
---|
703 | write(6,'(i3,1x,5(1x,1pe10.3))') k, grid%ph_1(1,k,1), & |
---|
704 | grid%p(1,k,1), grid%al(1,k,1), & |
---|
705 | grid%t_1(1,k,1), moist(1,k,1,P_QV) |
---|
706 | enddo |
---|
707 | ENDIF |
---|
708 | |
---|
709 | ! interp v |
---|
710 | |
---|
711 | DO J = jts, jte |
---|
712 | DO I = its, min(ide-1,ite) |
---|
713 | |
---|
714 | IF (j == jds) THEN |
---|
715 | z_at_v = grid%phb(i,1,j)/g |
---|
716 | ELSE IF (j == jde) THEN |
---|
717 | z_at_v = grid%phb(i,1,j-1)/g |
---|
718 | ELSE |
---|
719 | z_at_v = 0.5*(grid%phb(i,1,j)+grid%phb(i,1,j-1))/g |
---|
720 | END IF |
---|
721 | |
---|
722 | p_surf = interp_0( p_in, zk, z_at_v, nl_in ) |
---|
723 | |
---|
724 | DO K = 1, kte-1 |
---|
725 | p_level = grid%znu(k)*(p_surf - grid%p_top) + grid%p_top |
---|
726 | IF (.not.logp) THEN |
---|
727 | grid%v_1(i,k,j) = interp_0( v, p_in, p_level, nl_in ) |
---|
728 | ELSE |
---|
729 | grid%v_1(i,k,j) = interp_0_log( v, p_in, p_level, nl_in ) |
---|
730 | ENDIF |
---|
731 | grid%v_2(i,k,j) = grid%v_1(i,k,j) |
---|
732 | ENDDO |
---|
733 | |
---|
734 | ENDDO |
---|
735 | ENDDO |
---|
736 | |
---|
737 | ! interp u |
---|
738 | |
---|
739 | DO J = jts, min(jde-1,jte) |
---|
740 | DO I = its, ite |
---|
741 | |
---|
742 | IF (i == ids) THEN |
---|
743 | z_at_u = grid%phb(i,1,j)/g |
---|
744 | ELSE IF (i == ide) THEN |
---|
745 | z_at_u = grid%phb(i-1,1,j)/g |
---|
746 | ELSE |
---|
747 | z_at_u = 0.5*(grid%phb(i,1,j)+grid%phb(i-1,1,j))/g |
---|
748 | END IF |
---|
749 | |
---|
750 | p_surf = interp_0( p_in, zk, z_at_u, nl_in ) |
---|
751 | |
---|
752 | DO K = 1, kte-1 |
---|
753 | p_level = grid%znu(k)*(p_surf - grid%p_top) + grid%p_top |
---|
754 | IF (.not.logp) THEN |
---|
755 | grid%u_1(i,k,j) = interp_0( u, p_in, p_level, nl_in ) |
---|
756 | ELSE |
---|
757 | grid%u_1(i,k,j) = interp_0_log( u, p_in, p_level, nl_in ) |
---|
758 | ENDIF |
---|
759 | grid%u_2(i,k,j) = grid%u_1(i,k,j) |
---|
760 | ENDDO |
---|
761 | |
---|
762 | ENDDO |
---|
763 | ENDDO |
---|
764 | |
---|
765 | ! set w |
---|
766 | |
---|
767 | DO J = jts, min(jde-1,jte) |
---|
768 | DO K = kts, kte |
---|
769 | DO I = its, min(ide-1,ite) |
---|
770 | grid%w_1(i,k,j) = 0. |
---|
771 | grid%w_2(i,k,j) = 0. |
---|
772 | ENDDO |
---|
773 | ENDDO |
---|
774 | ENDDO |
---|
775 | |
---|
776 | !!!MARS MARS |
---|
777 | IF (config_flags%init_MU .ne. 0.) THEN |
---|
778 | grid%u_1 = grid%u_1*config_flags%init_MU |
---|
779 | grid%u_2 = grid%u_2*config_flags%init_MU |
---|
780 | print *, 'multiply zonal wind ', config_flags%init_MU |
---|
781 | ENDIF |
---|
782 | IF (config_flags%init_MV .ne. 0.) THEN |
---|
783 | grid%v_1 = grid%v_1*config_flags%init_MV |
---|
784 | grid%v_2 = grid%v_2*config_flags%init_MV |
---|
785 | print *, 'multiply meridional wind ', config_flags%init_MV |
---|
786 | ENDIF |
---|
787 | IF (config_flags%init_U .ne. 0.) THEN |
---|
788 | DO J = jts, min(jde-1,jte) |
---|
789 | DO K = kts, kte-1 |
---|
790 | DO I = its, min(ide-1,ite) |
---|
791 | grid%u_1(i,k,j) = config_flags%init_U |
---|
792 | grid%u_2(i,k,j) = config_flags%init_U |
---|
793 | ENDDO |
---|
794 | ENDDO |
---|
795 | ENDDO |
---|
796 | print *, 'constant zonal wind ', config_flags%init_U |
---|
797 | !!! ****** ou autre possibilité |
---|
798 | !!! > grid%u_1 = grid%u_1*0. + config_flags%init_U |
---|
799 | !!! > grid%u_2 = grid%u_2*0. + config_flags%init_U |
---|
800 | ENDIF |
---|
801 | IF (config_flags%init_V .ne. 0.) THEN |
---|
802 | DO J = jts, min(jde-1,jte) |
---|
803 | DO K = kts, kte-1 |
---|
804 | DO I = its, min(ide-1,ite) |
---|
805 | grid%v_1(i,k,j) = config_flags%init_V |
---|
806 | grid%v_2(i,k,j) = config_flags%init_V |
---|
807 | ENDDO |
---|
808 | ENDDO |
---|
809 | ENDDO |
---|
810 | print *, 'constant meridional wind ', config_flags%init_V |
---|
811 | ENDIF |
---|
812 | |
---|
813 | !!!MARS MARS |
---|
814 | |
---|
815 | |
---|
816 | ! set a few more things |
---|
817 | |
---|
818 | DO J = jts, min(jde-1,jte) |
---|
819 | DO K = kts, kte-1 |
---|
820 | DO I = its, min(ide-1,ite) |
---|
821 | grid%h_diabatic(i,k,j) = 0. |
---|
822 | !!!!! MARS NO WIND CASE |
---|
823 | !grid%u_1(i,k,j) = 0. |
---|
824 | !grid%u_2(i,k,j) = 0. |
---|
825 | !grid%v_1(i,k,j) = 0. |
---|
826 | !grid%v_2(i,k,j) = 0. |
---|
827 | !!!!! MARS NO WIND CASE |
---|
828 | ENDDO |
---|
829 | ENDDO |
---|
830 | ENDDO |
---|
831 | |
---|
832 | IF ( wrf_dm_on_monitor() ) THEN |
---|
833 | DO k=1,kte-1 |
---|
834 | grid%t_base(k) = grid%t_1(1,k,1) |
---|
835 | grid%qv_base(k) = moist(1,k,1,P_QV) |
---|
836 | grid%u_base(k) = grid%u_1(1,k,1) |
---|
837 | grid%v_base(k) = grid%v_1(1,k,1) |
---|
838 | 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 |
---|
839 | |
---|
840 | !!!!! MARS SIMPLE LES (PURE BUOYANCY) |
---|
841 | !! grid%t_base(k) = grid%t_init(its,k,jts) |
---|
842 | ! grid%t_base(k) = 0. |
---|
843 | ! grid%qv_base(k) = 0. |
---|
844 | ! grid%u_base(k) = 0. |
---|
845 | ! grid%v_base(k) = 0. |
---|
846 | ! grid%z_base(k) = 0. |
---|
847 | !!!!! MARS SIMPLE LES |
---|
848 | |
---|
849 | ENDDO |
---|
850 | ENDIF |
---|
851 | CALL wrf_dm_bcast_real( grid%t_base , kte ) |
---|
852 | CALL wrf_dm_bcast_real( grid%qv_base , kte ) |
---|
853 | CALL wrf_dm_bcast_real( grid%u_base , kte ) |
---|
854 | CALL wrf_dm_bcast_real( grid%v_base , kte ) |
---|
855 | CALL wrf_dm_bcast_real( grid%z_base , kte ) |
---|
856 | |
---|
857 | DO J = jts, min(jde-1,jte) |
---|
858 | DO I = its, min(ide-1,ite) |
---|
859 | thtmp = grid%t_2(i,1,j)+t0 |
---|
860 | ptmp = grid%p(i,1,j)+grid%pb(i,1,j) |
---|
861 | temp(1) = thtmp * (ptmp/p1000mb)**rcp |
---|
862 | thtmp = grid%t_2(i,2,j)+t0 |
---|
863 | ptmp = grid%p(i,2,j)+grid%pb(i,2,j) |
---|
864 | temp(2) = thtmp * (ptmp/p1000mb)**rcp |
---|
865 | thtmp = grid%t_2(i,3,j)+t0 |
---|
866 | ptmp = grid%p(i,3,j)+grid%pb(i,3,j) |
---|
867 | temp(3) = thtmp * (ptmp/p1000mb)**rcp |
---|
868 | |
---|
869 | !! For LES-CBL, add 5 degrees to the surface temperature! |
---|
870 | !! |
---|
871 | ! grid%tsk(I,J)=grid%cf1*temp(1)+grid%cf2*temp(2)+grid%cf3*temp(3) |
---|
872 | !! grid%tsk(I,J)=grid%cf1*temp(1)+grid%cf2*temp(2)+grid%cf3*temp(3)+5. |
---|
873 | grid%tmn(I,J)=grid%tsk(I,J)-0.5 |
---|
874 | |
---|
875 | ENDDO |
---|
876 | ENDDO |
---|
877 | |
---|
878 | !!!!! MARS |
---|
879 | |
---|
880 | ! interpolate water vapor |
---|
881 | if ( ( config_flags%mars == 1 ) & |
---|
882 | .OR. ( config_flags%mars == 11 ) & |
---|
883 | .OR. ( config_flags%mars == 12 ) ) then |
---|
884 | print *, '**** INTERPOLATE HV **** RANK 2 in SCALAR' |
---|
885 | DO k=1,kte-1 |
---|
886 | p_level = grid%znu(k)*(pd_surf - grid%p_top) + grid%p_top |
---|
887 | IF (.not.logp) THEN |
---|
888 | scalar(its:ite,k,jts:jte,2) = interp_0( qv, pd_in, p_level, nl_in ) |
---|
889 | ELSE |
---|
890 | scalar(its:ite,k,jts:jte,2) = interp_0_log( qv, pd_in, p_level, nl_in ) |
---|
891 | ENDIF |
---|
892 | scalar(its:ite,k,jts:jte,3) = 0. |
---|
893 | !! water ice is set to 0 (was put into water vapor when building prof from MCD) |
---|
894 | ENDDO |
---|
895 | print *, "WATER VAPOR",scalar(its,:,jts,2) |
---|
896 | endif |
---|
897 | |
---|
898 | ! interpolate qdust |
---|
899 | if ( ( config_flags%mars == 11 ) & |
---|
900 | .OR. ( config_flags%mars == 12 ) ) then |
---|
901 | call read_dust(profdustq,profdustn,nl_in) |
---|
902 | print *, '**** INTERPOLATE DUSTQ **** RANK 4 in SCALAR' |
---|
903 | print *, '**** INTERPOLATE DUSTN **** RANK 5 in SCALAR' |
---|
904 | DO k=1,kte-1 |
---|
905 | p_level = grid%znu(k)*(pd_surf - grid%p_top) + grid%p_top |
---|
906 | IF (.not.logp) THEN |
---|
907 | scalar(its:ite,k,jts:jte,4) = interp_0( profdustq, pd_in, p_level,nl_in ) |
---|
908 | scalar(its:ite,k,jts:jte,5) = interp_0( profdustn, pd_in, p_level,nl_in ) |
---|
909 | ELSE |
---|
910 | scalar(its:ite,k,jts:jte,4) = interp_0_log( profdustq, pd_in, p_level, nl_in ) |
---|
911 | scalar(its:ite,k,jts:jte,5) = interp_0_log( profdustn, pd_in, p_level, nl_in ) |
---|
912 | ENDIF |
---|
913 | ENDDO |
---|
914 | print *, "DUST Q", scalar(its,:,jts,4) |
---|
915 | print *, "DUST N", scalar(its,:,jts,5) |
---|
916 | endif |
---|
917 | |
---|
918 | if ( config_flags%mars == 12 ) then |
---|
919 | scalar(its:ite,1:kte-1,jts:jte,6) = 0. |
---|
920 | scalar(its:ite,1:kte-1,jts:jte,7) = 0. |
---|
921 | endif |
---|
922 | |
---|
923 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
924 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
925 | IF (planet .ne. "mars") THEN |
---|
926 | call read_dust(profdustq,profdustn,nl_in) |
---|
927 | DO k=1,kte!-1 |
---|
928 | p_level = grid%znu(k)*(pd_surf - grid%p_top) + grid%p_top |
---|
929 | DO j = jts, jte |
---|
930 | DO i = its, ite |
---|
931 | !!! we use Q2 as a vehicle for heating rates! sick! |
---|
932 | grid%m_q2(i,k,j) = interp_0_log( profdustq, pd_in, p_level, nl_in ) & |
---|
933 | + interp_0_log( profdustn, pd_in, p_level, nl_in ) |
---|
934 | ENDDO |
---|
935 | ENDDO |
---|
936 | !print*,'grid%m_q2' |
---|
937 | !print*,k,grid%m_q2(1,k,1) |
---|
938 | ENDDO |
---|
939 | ENDIF |
---|
940 | |
---|
941 | IF (planet.eq."prescribed") Then |
---|
942 | call read_hr(hrsw,hrlw,hrdyn,nl_in) |
---|
943 | open(unit=17,file="prescribed_sw.txt",action="write") |
---|
944 | open(unit=18,file="prescribed_lw.txt",action="write") |
---|
945 | open(unit=19,file="prescribed_dyn.txt",action="write") |
---|
946 | DO k=1,kte!-1 |
---|
947 | p_level = grid%znu(k)*(pd_surf - grid%p_top) + grid%p_top |
---|
948 | prescribed_sw(k) = interp_0_log( hrsw, pd_in, p_level, nl_in ) |
---|
949 | prescribed_lw(k) = interp_0_log( hrlw, pd_in, p_level, nl_in ) |
---|
950 | prescribed_dyn(k) = interp_0_log( hrdyn, pd_in, p_level, nl_in ) |
---|
951 | write (17,*) prescribed_sw(k) |
---|
952 | write (18,*) prescribed_lw(k) |
---|
953 | write (19,*) prescribed_dyn(k) |
---|
954 | ENDDO |
---|
955 | close(unit=19) |
---|
956 | close(unit=18) |
---|
957 | close(unit=17) |
---|
958 | ENDIF |
---|
959 | |
---|
960 | IF (planet.eq."venus") Then |
---|
961 | call read_hr(hrsw,hrlw,hrdyn,nl_in) |
---|
962 | open(unit=20,file="venus_hrdyn.txt",action="write") |
---|
963 | DO k=1,kte!-1 |
---|
964 | p_level = grid%znu(k)*(pd_surf - grid%p_top) + grid%p_top |
---|
965 | venus_hrdyn(k) = interp_0_log( hrdyn, pd_in, p_level, nl_in ) |
---|
966 | write (20,*) venus_hrdyn(k) |
---|
967 | ENDDO |
---|
968 | close(unit=20) |
---|
969 | ENDIF |
---|
970 | |
---|
971 | IF (planet.eq."generic") THEN |
---|
972 | call read_lsf(lsfdt,lsfdq,nl_in) |
---|
973 | open(unit=17,file="lsf.txt",action="write") |
---|
974 | DO k=1,kte!-1 |
---|
975 | p_level = grid%znu(k)*(pd_surf - grid%p_top) + grid%p_top |
---|
976 | lsf_dt = interp_0_log( lsfdt, pd_in, p_level, nl_in ) |
---|
977 | lsf_dq = interp_0_log( lsfdq, pd_in, p_level, nl_in ) |
---|
978 | write (17,*) lsf_dt(k),lsf_dq(k) |
---|
979 | ENDDO |
---|
980 | ENDIF |
---|
981 | |
---|
982 | IF ((planet.eq."venus") .AND. ( config_flags%mars == 34 )) Then |
---|
983 | pd_in8(:)=pd_in(:) |
---|
984 | do i = 1,34 |
---|
985 | call read_tracer(trac,num_scalar,i,nl_in) |
---|
986 | DO k=1,kte-1 |
---|
987 | p_level = grid%znu(k)*(pd_surf - grid%p_top) + grid%p_top |
---|
988 | p_level8=p_level |
---|
989 | scalar(its:ite,k,jts:jte,i+1) = interp_0_log2( trac, pd_in8, p_level8, nl_in ) |
---|
990 | ENDDO |
---|
991 | ENDDO |
---|
992 | !close(unit=22) |
---|
993 | ENDIF |
---|
994 | |
---|
995 | open(unit=21,file="altitude.txt",action="write") |
---|
996 | DO k=1,kte!-1 |
---|
997 | p_level = grid%znu(k)*(pd_surf - grid%p_top) + grid%p_top |
---|
998 | altitude(k) = interp_0_log( zk, pd_in, p_level, nl_in ) |
---|
999 | write (21,*) altitude(k) |
---|
1000 | ENDDO |
---|
1001 | close(unit=21) |
---|
1002 | |
---|
1003 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
1004 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
1005 | |
---|
1006 | |
---|
1007 | !!!!! MARS |
---|
1008 | |
---|
1009 | |
---|
1010 | END SUBROUTINE init_domain_rk |
---|
1011 | |
---|
1012 | SUBROUTINE init_module_initialize |
---|
1013 | END SUBROUTINE init_module_initialize |
---|
1014 | |
---|
1015 | !--------------------------------------------------------------------- |
---|
1016 | |
---|
1017 | ! test driver for get_sounding |
---|
1018 | ! |
---|
1019 | ! implicit none |
---|
1020 | ! integer n |
---|
1021 | ! parameter(n = 1000) |
---|
1022 | ! real zk(n),p(n),theta(n),rho(n),u(n),v(n),qv(n),pd(n) |
---|
1023 | ! logical dry |
---|
1024 | ! integer nl,k |
---|
1025 | ! |
---|
1026 | ! dry = .false. |
---|
1027 | ! dry = .true. |
---|
1028 | ! call get_sounding( zk, p, pd, theta, rho, u, v, qv, dry, n, nl ) |
---|
1029 | ! write(6,*) ' input levels ',nl |
---|
1030 | ! write(6,*) ' sounding ' |
---|
1031 | ! write(6,*) ' k height(m) press (Pa) pd(Pa) theta (K) den(kg/m^3) u(m/s) v(m/s) qv(g/g) ' |
---|
1032 | ! do k=1,nl |
---|
1033 | ! 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) |
---|
1034 | ! enddo |
---|
1035 | ! end |
---|
1036 | ! |
---|
1037 | !--------------------------------------------------------------------------- |
---|
1038 | |
---|
1039 | subroutine get_sounding( zk, p, p_dry, theta, rho, & |
---|
1040 | u, v, qv, dry, nl_max, nl_in ) |
---|
1041 | implicit none |
---|
1042 | |
---|
1043 | integer nl_max, nl_in |
---|
1044 | real zk(nl_max), p(nl_max), theta(nl_max), rho(nl_max), & |
---|
1045 | u(nl_max), v(nl_max), qv(nl_max), p_dry(nl_max) |
---|
1046 | logical dry |
---|
1047 | |
---|
1048 | integer n |
---|
1049 | parameter(n=1000) |
---|
1050 | logical debug |
---|
1051 | parameter( debug = .true.) |
---|
1052 | |
---|
1053 | ! input sounding data |
---|
1054 | |
---|
1055 | real p_surf, th_surf, qv_surf |
---|
1056 | real pi_surf, pi(n) |
---|
1057 | real h_input(n), th_input(n), qv_input(n), u_input(n), v_input(n) |
---|
1058 | |
---|
1059 | ! input therm data (element 0 is the ground so it's n+1 but n is 1000 anyway so...) |
---|
1060 | |
---|
1061 | real r_therm(n),cp_therm(n),p_therm(n),rho_therm(n),t_therm(n) |
---|
1062 | |
---|
1063 | ! diagnostics |
---|
1064 | |
---|
1065 | real rho_surf, p_input(n), rho_input(n) |
---|
1066 | real pm_input(n) ! this are for full moist sounding |
---|
1067 | |
---|
1068 | ! local data |
---|
1069 | |
---|
1070 | ! real p1000mb,cv,cp,r,cvpm,g |
---|
1071 | ! parameter (p1000mb = 1.e+05, r = 287, cp = 1003., cv = cp-r, cvpm = -cv/cp, g=9.81 ) |
---|
1072 | ! parameter (p1000mb = 610., r = 192., cp = 844.6, cv = cp-r, cvpm = -cv/cp, g=3.72) |
---|
1073 | ! parameter (p1000mb = 610., r = 191., cp = 744.5, cv = cp-r, cvpm = -cv/cp, g=3.72) |
---|
1074 | !parameter (p1000mb = 92.e+05, r = 192., cp = 900.0, cv = cp-r, cvpm = -cv/cp, g=8.87) |
---|
1075 | integer k, it, nl |
---|
1076 | real qvf, qvf1, dz |
---|
1077 | |
---|
1078 | LOGICAL :: direct_from_file |
---|
1079 | |
---|
1080 | IF (planet == "mars") THEN |
---|
1081 | direct_from_file = .false. |
---|
1082 | ELSE |
---|
1083 | direct_from_file = .true. |
---|
1084 | ENDIF |
---|
1085 | |
---|
1086 | ! first, read the sounding |
---|
1087 | |
---|
1088 | call read_sounding( p_surf, th_surf, qv_surf, & |
---|
1089 | h_input, th_input, qv_input, & |
---|
1090 | u_input, v_input,n, nl, debug ) |
---|
1091 | |
---|
1092 | ! and the therm : |
---|
1093 | |
---|
1094 | call read_therm(r_therm,cp_therm,p_therm,rho_therm,t_therm,n) |
---|
1095 | |
---|
1096 | if(debug) write(6,*) ' number of input levels = ',nl |
---|
1097 | nl_in = nl |
---|
1098 | if(nl_in .gt. nl_max ) then |
---|
1099 | write(6,*) ' too many levels for input arrays ',nl_in,nl_max |
---|
1100 | call wrf_error_fatal ( ' too many levels for input arrays ' ) |
---|
1101 | end if |
---|
1102 | |
---|
1103 | IF (.NOT. direct_from_file) THEN |
---|
1104 | |
---|
1105 | ! To use r/cp as defined above, one has to recompute teta from T (default MCD computes |
---|
1106 | ! teta for a variable r/cp) |
---|
1107 | |
---|
1108 | do k=1,nl |
---|
1109 | th_input(k) = t_therm(k)*(p1000mb/p_therm(k))**(rcp) |
---|
1110 | enddo |
---|
1111 | th_surf = t_therm(1)*(p1000mb/p_therm(1))**(rcp) |
---|
1112 | ! ----- |
---|
1113 | |
---|
1114 | !if(dry) then |
---|
1115 | ! do k=1,nl |
---|
1116 | ! qv_input(k) = 0. |
---|
1117 | ! enddo |
---|
1118 | !endif |
---|
1119 | |
---|
1120 | ! compute diagnostics, |
---|
1121 | ! first, convert qv(g/kg) to qv(g/g) |
---|
1122 | |
---|
1123 | do k=1,nl |
---|
1124 | !!!!!!!!!!!!!! MARS |
---|
1125 | !! from mol/mol to kg/kg |
---|
1126 | qv_input(k) = qv_input(k)*18./mwdry |
---|
1127 | qv_input(k) = 0.001*qv_input(k) |
---|
1128 | enddo |
---|
1129 | |
---|
1130 | p_surf = 100.*p_surf ! convert to pascals |
---|
1131 | qvf = 1. + rvovrd*qv_input(1) |
---|
1132 | !!MARS |
---|
1133 | qvf = 1. |
---|
1134 | rho_surf = 1./((r_d/p1000mb)*th_surf*qvf*((p_surf/p1000mb)**cvpm)) |
---|
1135 | pi_surf = (p_surf/p1000mb)**(rcp) |
---|
1136 | |
---|
1137 | if(debug) then |
---|
1138 | write(6,*) ' surface density is ',rho_surf |
---|
1139 | write(6,*) ' surface pi is ',pi_surf |
---|
1140 | end if |
---|
1141 | |
---|
1142 | |
---|
1143 | ! integrate moist sounding hydrostatically, starting from the |
---|
1144 | ! specified surface pressure |
---|
1145 | ! -> first, integrate from surface to lowest level |
---|
1146 | |
---|
1147 | qvf = 1. + rvovrd*qv_input(1) |
---|
1148 | qvf1 = 1. + qv_input(1) |
---|
1149 | !!MARS |
---|
1150 | qvf = 1. |
---|
1151 | qvf1 = 1. |
---|
1152 | rho_input(1) = rho_surf |
---|
1153 | dz = h_input(1) |
---|
1154 | do it=1,10 |
---|
1155 | ! pm_input(1) = p_surf & |
---|
1156 | ! - 0.5*dz*(rho_surf+rho_input(1))*g*qvf1 |
---|
1157 | !!!MARS MARS MARS |
---|
1158 | pm_input(1) = p_surf |
---|
1159 | rho_input(1) = 1./((r_d/p1000mb)*th_input(1)*qvf*((pm_input(1)/p1000mb)**cvpm)) |
---|
1160 | enddo |
---|
1161 | |
---|
1162 | ! integrate up the column |
---|
1163 | |
---|
1164 | do k=2,nl |
---|
1165 | rho_input(k) = rho_input(k-1) |
---|
1166 | dz = h_input(k)-h_input(k-1) |
---|
1167 | qvf1 = 0.5*(2.+(qv_input(k-1)+qv_input(k))) |
---|
1168 | qvf = 1. + rvovrd*qv_input(k) ! qv is in g/kg here |
---|
1169 | !!MARS |
---|
1170 | qvf = 1. |
---|
1171 | qvf1 = 1. |
---|
1172 | !!MARS |
---|
1173 | |
---|
1174 | do it=1,10 |
---|
1175 | pm_input(k) = pm_input(k-1) & |
---|
1176 | - 0.5*dz*(rho_input(k)+rho_input(k-1))*g*qvf1 |
---|
1177 | rho_input(k) = 1./((r_d/p1000mb)*th_input(k)*qvf*((pm_input(k)/p1000mb)**cvpm)) |
---|
1178 | enddo |
---|
1179 | enddo |
---|
1180 | |
---|
1181 | ! we have the moist sounding |
---|
1182 | |
---|
1183 | ! next, compute the dry sounding using p at the highest level from the |
---|
1184 | ! moist sounding and integrating down. |
---|
1185 | |
---|
1186 | p_input(nl) = pm_input(nl) |
---|
1187 | |
---|
1188 | do k=nl-1,1,-1 |
---|
1189 | dz = h_input(k+1)-h_input(k) |
---|
1190 | p_input(k) = p_input(k+1) + 0.5*dz*(rho_input(k)+rho_input(k+1))*g |
---|
1191 | enddo |
---|
1192 | |
---|
1193 | ELSE !IF (.NOT. direct_from_file) THEN |
---|
1194 | |
---|
1195 | do k=1,nl |
---|
1196 | !!!! direct input from file |
---|
1197 | write(6,*) '*** DIRECT INPUT FROM FILE ***' |
---|
1198 | pm_input(k) = p_therm(k) |
---|
1199 | p_input(k) = p_therm(k) |
---|
1200 | rho_input(k) = rho_therm(k) |
---|
1201 | enddo |
---|
1202 | |
---|
1203 | ENDIF |
---|
1204 | |
---|
1205 | |
---|
1206 | do k=1,nl |
---|
1207 | |
---|
1208 | zk(k) = h_input(k) |
---|
1209 | p(k) = pm_input(k) |
---|
1210 | p_dry(k) = p_input(k) |
---|
1211 | theta(k) = th_input(k) |
---|
1212 | rho(k) = rho_input(k) |
---|
1213 | u(k) = u_input(k) |
---|
1214 | v(k) = v_input(k) |
---|
1215 | qv(k) = qv_input(k) |
---|
1216 | |
---|
1217 | enddo |
---|
1218 | |
---|
1219 | if(debug) then |
---|
1220 | write(6,*) ' sounding ' |
---|
1221 | write(6,*) ' k height(m) press (Pa) pd(Pa) theta (K) den(kg/m^3) u(m/s) v(m/s) qv(g/g) ' |
---|
1222 | do k=1,nl |
---|
1223 | 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) |
---|
1224 | enddo |
---|
1225 | |
---|
1226 | end if |
---|
1227 | |
---|
1228 | end subroutine get_sounding |
---|
1229 | |
---|
1230 | !------------------------------------------------------- |
---|
1231 | |
---|
1232 | subroutine read_sounding( ps,ts,qvs,h,th,qv,u,v,n,nl,debug ) |
---|
1233 | implicit none |
---|
1234 | integer n,nl |
---|
1235 | real ps,ts,qvs,h(n),th(n),qv(n),u(n),v(n) |
---|
1236 | logical end_of_file |
---|
1237 | logical debug |
---|
1238 | |
---|
1239 | integer k |
---|
1240 | |
---|
1241 | open(unit=10,file='input_sounding',form='formatted',status='old') |
---|
1242 | rewind(10) |
---|
1243 | read(10,*) ps, ts, qvs |
---|
1244 | if(debug) then |
---|
1245 | write(6,*) ' input sounding surface parameters ' |
---|
1246 | write(6,*) ' surface pressure (mb) ',ps |
---|
1247 | write(6,*) ' surface pot. temp (K) ',ts |
---|
1248 | write(6,*) ' surface mixing ratio (g/kg) ',qvs |
---|
1249 | end if |
---|
1250 | |
---|
1251 | end_of_file = .false. |
---|
1252 | k = 0 |
---|
1253 | |
---|
1254 | do while (.not. end_of_file) |
---|
1255 | |
---|
1256 | read(10,*,end=100) h(k+1), th(k+1), qv(k+1), u(k+1), v(k+1) |
---|
1257 | k = k+1 |
---|
1258 | if(debug) write(6,'(1x,i3,5(1x,e10.3))') k, h(k), th(k), qv(k), u(k), v(k) |
---|
1259 | go to 110 |
---|
1260 | 100 end_of_file = .true. |
---|
1261 | 110 continue |
---|
1262 | enddo |
---|
1263 | |
---|
1264 | nl = k |
---|
1265 | |
---|
1266 | close(unit=10,status = 'keep') |
---|
1267 | |
---|
1268 | end subroutine read_sounding |
---|
1269 | |
---|
1270 | subroutine read_therm(r,cp,p,rho,t,n) |
---|
1271 | implicit none |
---|
1272 | integer n |
---|
1273 | real r(n),cp(n),p(n),rho(n),t(n) |
---|
1274 | logical end_of_file |
---|
1275 | |
---|
1276 | integer k |
---|
1277 | |
---|
1278 | ! first element is the surface |
---|
1279 | |
---|
1280 | open(unit=11,file='input_therm',form='formatted',status='old') |
---|
1281 | rewind(11) |
---|
1282 | end_of_file = .false. |
---|
1283 | k = 0 |
---|
1284 | do while (.not. end_of_file) |
---|
1285 | |
---|
1286 | read(11,*,end=101) r(k+1), cp(k+1), p(k+1), rho(k+1), t(k+1) |
---|
1287 | write(*,*) k, r(k+1), cp(k+1), p(k+1), rho(k+1), t(k+1) |
---|
1288 | k = k+1 |
---|
1289 | go to 112 |
---|
1290 | 101 end_of_file = .true. |
---|
1291 | 112 continue |
---|
1292 | enddo |
---|
1293 | |
---|
1294 | close(unit=11,status = 'keep') |
---|
1295 | |
---|
1296 | end subroutine read_therm |
---|
1297 | |
---|
1298 | !!cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
1299 | subroutine read_dust(pdustq,pdustn,n) |
---|
1300 | implicit none |
---|
1301 | integer n |
---|
1302 | real pdustq(n+1),pdustn(n+1) |
---|
1303 | logical end_of_file |
---|
1304 | |
---|
1305 | integer k |
---|
1306 | |
---|
1307 | ! first element is the surface |
---|
1308 | |
---|
1309 | open(unit=11,file='input_dust',form='formatted',status='old') |
---|
1310 | rewind(11) |
---|
1311 | end_of_file = .false. |
---|
1312 | k = 0 |
---|
1313 | do while (.not. end_of_file) |
---|
1314 | |
---|
1315 | read(11,*,end=102) pdustq(k+1),pdustn(k+1) |
---|
1316 | write(*,*) k,pdustq(k+1),pdustn(k+1) |
---|
1317 | k = k+1 |
---|
1318 | go to 113 |
---|
1319 | 102 end_of_file = .true. |
---|
1320 | 113 continue |
---|
1321 | enddo |
---|
1322 | |
---|
1323 | close(unit=11,status = 'keep') |
---|
1324 | |
---|
1325 | end subroutine read_dust |
---|
1326 | !!cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
1327 | subroutine read_hr(hr_sw,hr_lw,hr_dyn,n) |
---|
1328 | implicit none |
---|
1329 | integer n |
---|
1330 | real hr_sw(n+1),hr_lw(n+1),hr_dyn(n+1) |
---|
1331 | logical end_of_file |
---|
1332 | |
---|
1333 | integer k |
---|
1334 | |
---|
1335 | ! first element is the surface |
---|
1336 | |
---|
1337 | open(unit=11,file='input_hr',form='formatted',status='old') |
---|
1338 | rewind(11) |
---|
1339 | end_of_file = .false. |
---|
1340 | k = 0 |
---|
1341 | do while (.not. end_of_file) |
---|
1342 | |
---|
1343 | read(11,*,end=103) hr_sw(k+1),hr_lw(k+1),hr_dyn(k+1) |
---|
1344 | write(*,*) k,hr_sw(k+1),hr_lw(k+1),hr_dyn(k+1) |
---|
1345 | k = k+1 |
---|
1346 | go to 114 |
---|
1347 | 103 end_of_file = .true. |
---|
1348 | 114 continue |
---|
1349 | enddo |
---|
1350 | |
---|
1351 | close(unit=11,status = 'keep') |
---|
1352 | |
---|
1353 | end subroutine read_hr |
---|
1354 | |
---|
1355 | subroutine read_lsf(dt,dq,n) |
---|
1356 | implicit none |
---|
1357 | integer n |
---|
1358 | real dt(n+1),dq(n+1) |
---|
1359 | logical end_of_file |
---|
1360 | |
---|
1361 | integer k |
---|
1362 | |
---|
1363 | ! first element is the surface |
---|
1364 | |
---|
1365 | open(unit=12,file='input_lsf',form='formatted',status='old') |
---|
1366 | rewind(12) |
---|
1367 | end_of_file = .false. |
---|
1368 | k = 0 |
---|
1369 | do while (.not. end_of_file) |
---|
1370 | |
---|
1371 | read(12,*,end=103) dt(k+1),dq(k+1) |
---|
1372 | write(*,*) k,dt(k+1),dq(k+1) |
---|
1373 | k = k+1 |
---|
1374 | go to 114 |
---|
1375 | 103 end_of_file = .true. |
---|
1376 | 114 continue |
---|
1377 | enddo |
---|
1378 | |
---|
1379 | close(unit=12,status = 'keep') |
---|
1380 | |
---|
1381 | end subroutine read_lsf |
---|
1382 | |
---|
1383 | subroutine read_tracer(trace,nq,qn,n) |
---|
1384 | implicit none |
---|
1385 | integer n,qn,nq ! qn : number of the tracer |
---|
1386 | real*8 tra(nq-1,n+1) |
---|
1387 | real*8 trace(n+1) !output |
---|
1388 | logical end_of_file |
---|
1389 | |
---|
1390 | integer k,j |
---|
1391 | |
---|
1392 | ! first element is the surface |
---|
1393 | open(unit=14,file='input_tracer',form='formatted',status='old') |
---|
1394 | rewind(14) |
---|
1395 | end_of_file = .false. |
---|
1396 | DO k=1,n |
---|
1397 | read(14,*) tra(:,k) |
---|
1398 | write(*,*) k,tra(qn,k) |
---|
1399 | ENDDO |
---|
1400 | |
---|
1401 | close(14) |
---|
1402 | trace(:)=tra(qn,:) |
---|
1403 | end subroutine read_tracer |
---|
1404 | |
---|
1405 | |
---|
1406 | END MODULE module_initialize_ideal |
---|