1 | !WRF:DRIVER_LAYER:MAIN |
---|
2 | ! |
---|
3 | |
---|
4 | PROGRAM ndown_em |
---|
5 | |
---|
6 | USE module_machine |
---|
7 | USE module_domain, ONLY : domain, head_grid, alloc_and_configure_domain, & |
---|
8 | domain_clock_set, domain_clock_get, get_ijk_from_grid |
---|
9 | USE module_domain_type, ONLY : program_name |
---|
10 | USE module_initialize_real, ONLY : wrfu_initialize, rebalance_driver |
---|
11 | USE module_integrate |
---|
12 | USE module_driver_constants |
---|
13 | USE module_configure, ONLY : grid_config_rec_type, model_config_rec |
---|
14 | USE module_io_domain |
---|
15 | USE module_utility |
---|
16 | USE module_check_a_mundo |
---|
17 | |
---|
18 | USE module_timing |
---|
19 | USE module_wrf_error |
---|
20 | #ifdef DM_PARALLEL |
---|
21 | USE module_dm |
---|
22 | #endif |
---|
23 | |
---|
24 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
25 | !new for bc |
---|
26 | USE module_bc |
---|
27 | USE module_big_step_utilities_em |
---|
28 | USE module_get_file_names |
---|
29 | #ifdef WRF_CHEM |
---|
30 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
31 | ! for chemistry |
---|
32 | USE module_input_chem_data |
---|
33 | ! USE module_input_chem_bioemiss |
---|
34 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
35 | #endif |
---|
36 | |
---|
37 | IMPLICIT NONE |
---|
38 | ! interface |
---|
39 | INTERFACE |
---|
40 | ! mediation-supplied |
---|
41 | SUBROUTINE med_read_wrf_chem_bioemiss ( grid , config_flags) |
---|
42 | USE module_domain |
---|
43 | TYPE (domain) grid |
---|
44 | TYPE (grid_config_rec_type) config_flags |
---|
45 | END SUBROUTINE med_read_wrf_chem_bioemiss |
---|
46 | |
---|
47 | SUBROUTINE init_domain_constants_em_ptr ( parent , nest ) |
---|
48 | USE module_domain |
---|
49 | USE module_configure |
---|
50 | TYPE(domain), POINTER :: parent , nest |
---|
51 | END SUBROUTINE init_domain_constants_em_ptr |
---|
52 | |
---|
53 | SUBROUTINE vertical_interp (nested_grid,znw_c,znu_c,cf1_c,cf2_c,cf3_c,cfn_c,cfn1_c,k_dim_c) |
---|
54 | USE module_domain |
---|
55 | USE module_configure |
---|
56 | TYPE(domain), POINTER :: nested_grid |
---|
57 | INTEGER , INTENT (IN) :: k_dim_c |
---|
58 | REAL , INTENT (IN) :: cf1_c,cf2_c,cf3_c,cfn_c,cfn1_c |
---|
59 | REAL , DIMENSION(k_dim_c) , INTENT (IN) :: znw_c,znu_c |
---|
60 | END SUBROUTINE vertical_interp |
---|
61 | |
---|
62 | |
---|
63 | END INTERFACE |
---|
64 | INTEGER :: ids , ide , jds , jde , kds , kde |
---|
65 | INTEGER :: ims , ime , jms , jme , kms , kme |
---|
66 | INTEGER :: ips , ipe , jps , jpe , kps , kpe |
---|
67 | INTEGER :: its , ite , jts , jte , kts , kte |
---|
68 | INTEGER :: nids, nide, njds, njde, nkds, nkde, & |
---|
69 | nims, nime, njms, njme, nkms, nkme, & |
---|
70 | nips, nipe, njps, njpe, nkps, nkpe |
---|
71 | INTEGER :: spec_bdy_width |
---|
72 | INTEGER :: i , j , k , nvchem |
---|
73 | INTEGER :: time_loop_max , time_loop |
---|
74 | INTEGER :: total_time_sec , file_counter |
---|
75 | INTEGER :: julyr , julday , iswater , map_proj |
---|
76 | INTEGER :: icnt |
---|
77 | |
---|
78 | REAL :: dt , new_bdy_frq |
---|
79 | REAL :: gmt , cen_lat , cen_lon , dx , dy , truelat1 , truelat2 , moad_cen_lat , stand_lon |
---|
80 | |
---|
81 | REAL , DIMENSION(:,:,:) , ALLOCATABLE :: ubdy3dtemp1 , vbdy3dtemp1 , tbdy3dtemp1 , pbdy3dtemp1 , qbdy3dtemp1 |
---|
82 | REAL , DIMENSION(:,:,:) , ALLOCATABLE :: mbdy2dtemp1 |
---|
83 | REAL , DIMENSION(:,:,:) , ALLOCATABLE :: ubdy3dtemp2 , vbdy3dtemp2 , tbdy3dtemp2 , pbdy3dtemp2 , qbdy3dtemp2 |
---|
84 | REAL , DIMENSION(:,:,:) , ALLOCATABLE :: mbdy2dtemp2 |
---|
85 | REAL , DIMENSION(:,:,:) , ALLOCATABLE :: cbdy3dtemp1 , cbdy3dtemp2 |
---|
86 | REAL , DIMENSION(:,:,:,:) , ALLOCATABLE :: cbdy3dtemp0 |
---|
87 | |
---|
88 | CHARACTER(LEN=19) :: start_date_char , current_date_char , end_date_char |
---|
89 | CHARACTER(LEN=19) :: stopTimeStr |
---|
90 | |
---|
91 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
92 | |
---|
93 | INTEGER :: num_veg_cat , num_soil_top_cat , num_soil_bot_cat |
---|
94 | |
---|
95 | REAL :: time |
---|
96 | INTEGER :: rc |
---|
97 | |
---|
98 | INTEGER :: loop , levels_to_process |
---|
99 | INTEGER , PARAMETER :: max_sanity_file_loop = 100 |
---|
100 | |
---|
101 | TYPE (domain) , POINTER :: keep_grid, grid_ptr, null_domain, parent_grid , nested_grid |
---|
102 | TYPE (domain) :: dummy |
---|
103 | TYPE (grid_config_rec_type) :: config_flags |
---|
104 | INTEGER :: number_at_same_level |
---|
105 | INTEGER :: time_step_begin_restart |
---|
106 | |
---|
107 | INTEGER :: max_dom , domain_id , fid , fido, fidb , oid , idum1 , idum2 , ierr |
---|
108 | INTEGER :: status_next_var |
---|
109 | INTEGER :: debug_level |
---|
110 | LOGICAL :: input_from_file , need_new_file |
---|
111 | CHARACTER (LEN=19) :: date_string |
---|
112 | |
---|
113 | #ifdef DM_PARALLEL |
---|
114 | INTEGER :: nbytes |
---|
115 | INTEGER, PARAMETER :: configbuflen = 4* CONFIG_BUF_LEN |
---|
116 | INTEGER :: configbuf( configbuflen ) |
---|
117 | LOGICAL , EXTERNAL :: wrf_dm_on_monitor |
---|
118 | #endif |
---|
119 | |
---|
120 | INTEGER :: idsi |
---|
121 | CHARACTER (LEN=80) :: inpname , outname , bdyname |
---|
122 | CHARACTER (LEN=80) :: si_inpname |
---|
123 | character *19 :: temp19 |
---|
124 | character *24 :: temp24 , temp24b |
---|
125 | character(len=24) :: start_date_hold |
---|
126 | |
---|
127 | CHARACTER (LEN=256) :: message |
---|
128 | integer :: ii |
---|
129 | |
---|
130 | #include "version_decl" |
---|
131 | |
---|
132 | !!!!!!!!!!!!!!!!!!!!! mousta |
---|
133 | integer :: n_ref_m,k_dim_c,k_dim_n |
---|
134 | real :: cf1_c,cf2_c,cf3_c,cfn_c,cfn1_c |
---|
135 | REAL , DIMENSION(:) , ALLOCATABLE :: znw_c,znu_c |
---|
136 | !!!!!!!!!!!!!!!!!!!!!!!!!!11 |
---|
137 | |
---|
138 | ! Interface block for routine that passes pointers and needs to know that they |
---|
139 | ! are receiving pointers. |
---|
140 | |
---|
141 | INTERFACE |
---|
142 | |
---|
143 | SUBROUTINE med_interp_domain ( parent_grid , nested_grid ) |
---|
144 | USE module_domain |
---|
145 | USE module_configure |
---|
146 | TYPE(domain), POINTER :: parent_grid , nested_grid |
---|
147 | END SUBROUTINE med_interp_domain |
---|
148 | |
---|
149 | SUBROUTINE Setup_Timekeeping( parent_grid ) |
---|
150 | USE module_domain |
---|
151 | TYPE(domain), POINTER :: parent_grid |
---|
152 | END SUBROUTINE Setup_Timekeeping |
---|
153 | |
---|
154 | SUBROUTINE vert_cor(parent_grid,znw_c,k_dim_c) |
---|
155 | USE module_domain |
---|
156 | TYPE(domain), POINTER :: parent_grid |
---|
157 | integer , intent(in) :: k_dim_c |
---|
158 | real , dimension(k_dim_c), INTENT(IN) :: znw_c |
---|
159 | END SUBROUTINE vert_cor |
---|
160 | END INTERFACE |
---|
161 | |
---|
162 | ! Define the name of this program (program_name defined in module_domain) |
---|
163 | |
---|
164 | program_name = "NDOWN_EM " // TRIM(release_version) // " PREPROCESSOR" |
---|
165 | |
---|
166 | #ifdef DM_PARALLEL |
---|
167 | CALL disable_quilting |
---|
168 | #endif |
---|
169 | |
---|
170 | ! Initialize the modules used by the WRF system. Many of the CALLs made from the |
---|
171 | ! init_modules routine are NO-OPs. Typical initializations are: the size of a |
---|
172 | ! REAL, setting the file handles to a pre-use value, defining moisture and |
---|
173 | ! chemistry indices, etc. |
---|
174 | |
---|
175 | CALL init_modules(1) ! Phase 1 returns after MPI_INIT() (if it is called) |
---|
176 | #ifdef NO_LEAP_CALENDAR |
---|
177 | CALL WRFU_Initialize( defaultCalendar=WRFU_CAL_NOLEAP, rc=rc ) |
---|
178 | #else |
---|
179 | CALL WRFU_Initialize( defaultCalendar=WRFU_CAL_GREGORIAN, rc=rc ) |
---|
180 | #endif |
---|
181 | CALL init_modules(2) ! Phase 2 resumes after MPI_INIT() (if it is called) |
---|
182 | |
---|
183 | ! Get the NAMELIST data. This is handled in the initial_config routine. All of the |
---|
184 | ! NAMELIST input variables are assigned to the model_config_rec structure. Below, |
---|
185 | ! note for parallel processing, only the monitor processor handles the raw Fortran |
---|
186 | ! I/O, and then broadcasts the info to each of the other nodes. |
---|
187 | |
---|
188 | #ifdef DM_PARALLEL |
---|
189 | IF ( wrf_dm_on_monitor() ) THEN |
---|
190 | CALL initial_config |
---|
191 | ENDIF |
---|
192 | CALL get_config_as_buffer( configbuf, configbuflen, nbytes ) |
---|
193 | CALL wrf_dm_bcast_bytes( configbuf, nbytes ) |
---|
194 | CALL set_config_as_buffer( configbuf, configbuflen ) |
---|
195 | CALL wrf_dm_initialize |
---|
196 | #else |
---|
197 | CALL initial_config |
---|
198 | #endif |
---|
199 | |
---|
200 | CALL check_nml_consistency |
---|
201 | CALL set_physics_rconfigs |
---|
202 | |
---|
203 | ! If we are running ndown, and that is WHERE we are now, make sure that we account |
---|
204 | ! for folks forgetting to say that the aux_input2 stream is in place. |
---|
205 | |
---|
206 | IF ( model_config_rec%io_form_auxinput2 .EQ. 0 ) THEN |
---|
207 | CALL wrf_error_fatal( 'ndown: Please set io_form_auxinput2 in the time_control namelist record (probably to 2).') |
---|
208 | END IF |
---|
209 | |
---|
210 | !!!!!!!!!!!!!!! mousta |
---|
211 | n_ref_m = model_config_rec % vert_refine_fact |
---|
212 | k_dim_c = model_config_rec % e_vert(1) |
---|
213 | k_dim_n = k_dim_c |
---|
214 | if (n_ref_m .ge. 2) k_dim_n = (k_dim_c - 1) * n_ref_m + 1 |
---|
215 | model_config_rec % e_vert(1) = k_dim_n |
---|
216 | model_config_rec % e_vert(2) = k_dim_n |
---|
217 | ALLOCATE(znw_c(k_dim_c)) |
---|
218 | ALLOCATE(znu_c(k_dim_c)) |
---|
219 | WRITE ( message , FMT = '(A,3I5)' ) 'KDIM_C', k_dim_c , model_config_rec % e_vert(1) , model_config_rec % e_vert(2) |
---|
220 | CALL wrf_debug ( 99,message ) |
---|
221 | !!!!!!!!!!!!!!! mousta |
---|
222 | |
---|
223 | ! And here is an instance of using the information in the NAMELIST. |
---|
224 | |
---|
225 | CALL nl_get_debug_level ( 1, debug_level ) |
---|
226 | CALL set_wrf_debug_level ( debug_level ) |
---|
227 | |
---|
228 | ! Allocated and configure the mother domain. Since we are in the nesting down |
---|
229 | ! mode, we know a) we got a nest, and b) we only got 1 nest. |
---|
230 | |
---|
231 | NULLIFY( null_domain ) |
---|
232 | |
---|
233 | CALL wrf_message ( program_name ) |
---|
234 | CALL wrf_debug ( 100 , 'ndown_em: calling alloc_and_configure_domain coarse ' ) |
---|
235 | CALL alloc_and_configure_domain ( domain_id = 1 , & |
---|
236 | grid = head_grid , & |
---|
237 | parent = null_domain , & |
---|
238 | kid = -1 ) |
---|
239 | |
---|
240 | parent_grid => head_grid |
---|
241 | |
---|
242 | ! Set up time initializations. |
---|
243 | |
---|
244 | CALL Setup_Timekeeping ( parent_grid ) |
---|
245 | |
---|
246 | CALL domain_clock_set( head_grid, & |
---|
247 | time_step_seconds=model_config_rec%interval_seconds ) |
---|
248 | CALL wrf_debug ( 100 , 'ndown_em: calling model_to_grid_config_rec ' ) |
---|
249 | CALL model_to_grid_config_rec ( parent_grid%id , model_config_rec , config_flags ) |
---|
250 | CALL wrf_debug ( 100 , 'ndown_em: calling set_scalar_indices_from_config ' ) |
---|
251 | CALL set_scalar_indices_from_config ( parent_grid%id , idum1, idum2 ) |
---|
252 | |
---|
253 | ! Initialize the I/O for WRF. |
---|
254 | |
---|
255 | CALL wrf_debug ( 100 , 'ndown_em: calling init_wrfio' ) |
---|
256 | CALL init_wrfio |
---|
257 | |
---|
258 | ! Some of the configuration values may have been modified from the initial READ |
---|
259 | ! of the NAMELIST, so we re-broadcast the configuration records. |
---|
260 | |
---|
261 | #ifdef DM_PARALLEL |
---|
262 | CALL get_config_as_buffer( configbuf, configbuflen, nbytes ) |
---|
263 | CALL wrf_dm_bcast_bytes( configbuf, nbytes ) |
---|
264 | CALL set_config_as_buffer( configbuf, configbuflen ) |
---|
265 | #endif |
---|
266 | |
---|
267 | ! We need to current and starting dates for the output files. The times need to be incremented |
---|
268 | ! so that the lateral BC files are not overwritten. |
---|
269 | |
---|
270 | #ifdef PLANET |
---|
271 | WRITE ( start_date_char , FMT = '(I4.4,"-",I5.5,"_",I2.2,":",I2.2,":",I2.2)' ) & |
---|
272 | model_config_rec%start_year (parent_grid%id) , & |
---|
273 | model_config_rec%start_day (parent_grid%id) , & |
---|
274 | model_config_rec%start_hour (parent_grid%id) , & |
---|
275 | model_config_rec%start_minute(parent_grid%id) , & |
---|
276 | model_config_rec%start_second(parent_grid%id) |
---|
277 | |
---|
278 | WRITE ( end_date_char , FMT = '(I4.4,"-",I5.5,"_",I2.2,":",I2.2,":",I2.2)' ) & |
---|
279 | model_config_rec% end_year (parent_grid%id) , & |
---|
280 | model_config_rec% end_day (parent_grid%id) , & |
---|
281 | model_config_rec% end_hour (parent_grid%id) , & |
---|
282 | model_config_rec% end_minute(parent_grid%id) , & |
---|
283 | model_config_rec% end_second(parent_grid%id) |
---|
284 | #else |
---|
285 | WRITE ( start_date_char , FMT = '(I4.4,"-",I2.2,"-",I2.2,"_",I2.2,":",I2.2,":",I2.2)' ) & |
---|
286 | model_config_rec%start_year (parent_grid%id) , & |
---|
287 | model_config_rec%start_month (parent_grid%id) , & |
---|
288 | model_config_rec%start_day (parent_grid%id) , & |
---|
289 | model_config_rec%start_hour (parent_grid%id) , & |
---|
290 | model_config_rec%start_minute(parent_grid%id) , & |
---|
291 | model_config_rec%start_second(parent_grid%id) |
---|
292 | |
---|
293 | WRITE ( end_date_char , FMT = '(I4.4,"-",I2.2,"-",I2.2,"_",I2.2,":",I2.2,":",I2.2)' ) & |
---|
294 | model_config_rec% end_year (parent_grid%id) , & |
---|
295 | model_config_rec% end_month (parent_grid%id) , & |
---|
296 | model_config_rec% end_day (parent_grid%id) , & |
---|
297 | model_config_rec% end_hour (parent_grid%id) , & |
---|
298 | model_config_rec% end_minute(parent_grid%id) , & |
---|
299 | model_config_rec% end_second(parent_grid%id) |
---|
300 | #endif |
---|
301 | |
---|
302 | ! Override stop time with value computed above. |
---|
303 | CALL domain_clock_set( parent_grid, stop_timestr=end_date_char ) |
---|
304 | |
---|
305 | CALL geth_idts ( end_date_char , start_date_char , total_time_sec ) |
---|
306 | |
---|
307 | new_bdy_frq = model_config_rec%interval_seconds |
---|
308 | time_loop_max = total_time_sec / model_config_rec%interval_seconds + 1 |
---|
309 | |
---|
310 | start_date = start_date_char // '.0000' |
---|
311 | current_date = start_date_char // '.0000' |
---|
312 | start_date_hold = start_date_char // '.0000' |
---|
313 | current_date_char = start_date_char |
---|
314 | |
---|
315 | ! Get a list of available file names to try. This fills up the eligible_file_name |
---|
316 | ! array with number_of_eligible_files entries. This routine issues a nonstandard |
---|
317 | ! call (system). |
---|
318 | |
---|
319 | file_counter = 1 |
---|
320 | need_new_file = .FALSE. |
---|
321 | CALL unix_ls ( 'wrfout' , parent_grid%id ) |
---|
322 | |
---|
323 | ! Open the input data (wrfout_d01_xxxxxx) for reading. |
---|
324 | |
---|
325 | CALL wrf_debug ( 100 , 'ndown_em main: calling open_r_dataset for ' // TRIM(eligible_file_name(file_counter)) ) |
---|
326 | CALL open_r_dataset ( fid, TRIM(eligible_file_name(file_counter)) , head_grid , config_flags , "DATASET=AUXINPUT1", ierr ) |
---|
327 | IF ( ierr .NE. 0 ) THEN |
---|
328 | WRITE( wrf_err_message , FMT='(A,A,A,I8)' ) 'program ndown: error opening ',TRIM(eligible_file_name(file_counter)), & |
---|
329 | ' for reading ierr=',ierr |
---|
330 | CALL WRF_ERROR_FATAL ( wrf_err_message ) |
---|
331 | ENDIF |
---|
332 | |
---|
333 | ! We know how many time periods to process, so we begin. |
---|
334 | |
---|
335 | big_time_loop_thingy : DO time_loop = 1 , time_loop_max |
---|
336 | |
---|
337 | ! Which date are we currently soliciting? |
---|
338 | |
---|
339 | CALL geth_newdate ( date_string , start_date_char , ( time_loop - 1 ) * NINT ( new_bdy_frq) ) |
---|
340 | WRITE ( message , FMT = '(A,I5," ",A,A)' ) '-------->>> Processing data: loop=',time_loop,' date/time = ',date_string |
---|
341 | CALL wrf_debug ( 99,message ) |
---|
342 | current_date_char = date_string |
---|
343 | current_date = date_string // '.0000' |
---|
344 | start_date = date_string // '.0000' |
---|
345 | WRITE ( message , FMT = '(A,I5," ",A,A)' ) 'loopmax = ', time_loop_max, ' ending date = ',end_date_char |
---|
346 | CALL wrf_debug ( 99,message ) |
---|
347 | CALL domain_clock_set( parent_grid, & |
---|
348 | current_timestr=current_date(1:19) ) |
---|
349 | |
---|
350 | ! Which times are in this file, and more importantly, are any of them the |
---|
351 | ! ones that we want? We need to loop over times in each files, loop |
---|
352 | ! over files. |
---|
353 | |
---|
354 | get_the_right_time : DO |
---|
355 | |
---|
356 | CALL wrf_get_next_time ( fid , date_string , status_next_var ) |
---|
357 | WRITE ( message , FMT = '(A,A,A,A,A,I5)' ) 'file date/time = ',date_string,' desired date = ',& |
---|
358 | current_date_char,' status = ', status_next_var |
---|
359 | CALL wrf_debug ( 99,message ) |
---|
360 | |
---|
361 | IF ( status_next_var .NE. 0 ) THEN |
---|
362 | CALL wrf_debug ( 100 , 'ndown_em main: calling close_dataset for ' // TRIM(eligible_file_name(file_counter)) ) |
---|
363 | CALL close_dataset ( fid , config_flags , "DATASET=INPUT" ) |
---|
364 | file_counter = file_counter + 1 |
---|
365 | IF ( file_counter .GT. number_of_eligible_files ) THEN |
---|
366 | WRITE( wrf_err_message , FMT='(A,A,A,I8)' ) 'program ndown: opening too many files' |
---|
367 | CALL WRF_ERROR_FATAL ( wrf_err_message ) |
---|
368 | END IF |
---|
369 | CALL wrf_debug ( 100 , 'ndown_em main: calling open_r_dataset for ' // TRIM(eligible_file_name(file_counter)) ) |
---|
370 | CALL open_r_dataset ( fid, TRIM(eligible_file_name(file_counter)) , head_grid , config_flags , "DATASET=INPUT", ierr ) |
---|
371 | IF ( ierr .NE. 0 ) THEN |
---|
372 | WRITE( wrf_err_message , FMT='(A,A,A,I8)' ) 'program ndown: error opening ',TRIM(eligible_file_name(file_counter)), & |
---|
373 | ' for reading ierr=',ierr |
---|
374 | CALL WRF_ERROR_FATAL ( wrf_err_message ) |
---|
375 | ENDIF |
---|
376 | CYCLE get_the_right_time |
---|
377 | ELSE IF ( TRIM(date_string) .LT. TRIM(current_date_char) ) THEN |
---|
378 | CYCLE get_the_right_time |
---|
379 | ELSE IF ( TRIM(date_string) .EQ. TRIM(current_date_char) ) THEN |
---|
380 | EXIT get_the_right_time |
---|
381 | ELSE IF ( TRIM(date_string) .GT. TRIM(current_date_char) ) THEN |
---|
382 | WRITE( wrf_err_message , FMT='(A,A,A,A,A)' ) 'Found ',TRIM(date_string),' before I found ',TRIM(current_date_char),'.' |
---|
383 | CALL WRF_ERROR_FATAL ( wrf_err_message ) |
---|
384 | END IF |
---|
385 | END DO get_the_right_time |
---|
386 | |
---|
387 | CALL wrf_debug ( 100 , 'wrf: calling input_history' ) |
---|
388 | CALL wrf_get_previous_time ( fid , date_string , status_next_var ) |
---|
389 | WRITE ( message , * ) 'CFB' ,cf1_c,cf2_c,cf3_c,cfn_c,cfn1_c,znw_c(1),znu_c(1) |
---|
390 | CALL wrf_debug ( 99,message ) |
---|
391 | CALL input_history ( fid , head_grid , config_flags, ierr) |
---|
392 | !!!!!!!!!!!!!1 mousta |
---|
393 | cf1_c = head_grid%cf1 |
---|
394 | cf2_c = head_grid%cf2 |
---|
395 | cf3_c = head_grid%cf3 |
---|
396 | |
---|
397 | cfn_c = head_grid%cfn |
---|
398 | cfn1_c = head_grid%cfn1 |
---|
399 | |
---|
400 | do k = 1,k_dim_c |
---|
401 | znw_c(k) = head_grid%znw(k) |
---|
402 | znu_c(k) = head_grid%znu(k) |
---|
403 | enddo |
---|
404 | call vert_cor(head_grid,znw_c,k_dim_c) |
---|
405 | WRITE ( message , * ) 'CFA' ,cf1_c,cf2_c,cf3_c,cfn_c,cfn1_c,znw_c(1),znu_c(1) |
---|
406 | CALL wrf_debug ( 99,message ) |
---|
407 | WRITE ( message , * ) 'CFV' ,head_grid%cf1,head_grid%cf2,head_grid%cf3,head_grid%cfn,head_grid%cfn1,& |
---|
408 | head_grid%znw(1),head_grid%znu(1) , head_grid%e_vert , parent_grid%cf1 , parent_grid%znw(1) , parent_grid%znu(1) |
---|
409 | CALL wrf_debug ( 99,message ) |
---|
410 | !!!!!!!!!!!!!1 mousta |
---|
411 | CALL wrf_debug ( 100 , 'wrf: back from input_history' ) |
---|
412 | |
---|
413 | ! Get the coarse grid info for later transfer to the fine grid domain. |
---|
414 | |
---|
415 | CALL wrf_get_dom_ti_integer ( fid , 'MAP_PROJ' , map_proj , 1 , icnt , ierr ) |
---|
416 | CALL wrf_get_dom_ti_real ( fid , 'DX' , dx , 1 , icnt , ierr ) |
---|
417 | CALL wrf_get_dom_ti_real ( fid , 'DY' , dy , 1 , icnt , ierr ) |
---|
418 | CALL wrf_get_dom_ti_real ( fid , 'CEN_LAT' , cen_lat , 1 , icnt , ierr ) |
---|
419 | CALL wrf_get_dom_ti_real ( fid , 'CEN_LON' , cen_lon , 1 , icnt , ierr ) |
---|
420 | CALL wrf_get_dom_ti_real ( fid , 'TRUELAT1' , truelat1 , 1 , icnt , ierr ) |
---|
421 | CALL wrf_get_dom_ti_real ( fid , 'TRUELAT2' , truelat2 , 1 , icnt , ierr ) |
---|
422 | CALL wrf_get_dom_ti_real ( fid , 'MOAD_CEN_LAT' , moad_cen_lat , 1 , icnt , ierr ) |
---|
423 | CALL wrf_get_dom_ti_real ( fid , 'STAND_LON' , stand_lon , 1 , icnt , ierr ) |
---|
424 | ! CALL wrf_get_dom_ti_real ( fid , 'GMT' , gmt , 1 , icnt , ierr ) |
---|
425 | ! CALL wrf_get_dom_ti_integer ( fid , 'JULYR' , julyr , 1 , icnt , ierr ) |
---|
426 | ! CALL wrf_get_dom_ti_integer ( fid , 'JULDAY' , julday , 1 , icnt , ierr ) |
---|
427 | CALL wrf_get_dom_ti_integer ( fid , 'ISWATER' , iswater , 1 , icnt , ierr ) |
---|
428 | |
---|
429 | ! First time in, do this: allocate sapce for the fine grid, get the config flags, open the |
---|
430 | ! wrfinput and wrfbdy files. This COULD be done outside the time loop, I think, so check it |
---|
431 | ! out and move it up if you can. |
---|
432 | |
---|
433 | IF ( time_loop .EQ. 1 ) THEN |
---|
434 | |
---|
435 | CALL wrf_message ( program_name ) |
---|
436 | CALL wrf_debug ( 100 , 'wrf: calling alloc_and_configure_domain fine ' ) |
---|
437 | CALL alloc_and_configure_domain ( domain_id = 2 , & |
---|
438 | grid = nested_grid , & |
---|
439 | parent = parent_grid , & |
---|
440 | kid = 1 ) |
---|
441 | |
---|
442 | CALL wrf_debug ( 100 , 'wrf: calling model_to_grid_config_rec ' ) |
---|
443 | CALL model_to_grid_config_rec ( nested_grid%id , model_config_rec , config_flags ) |
---|
444 | CALL wrf_debug ( 100 , 'wrf: calling set_scalar_indices_from_config ' ) |
---|
445 | CALL set_scalar_indices_from_config ( nested_grid%id , idum1, idum2 ) |
---|
446 | |
---|
447 | ! Set up time initializations for the fine grid. |
---|
448 | |
---|
449 | CALL Setup_Timekeeping ( nested_grid ) |
---|
450 | ! Strictly speaking, nest stop time should come from model_config_rec... |
---|
451 | CALL domain_clock_get( parent_grid, stop_timestr=stopTimeStr ) |
---|
452 | CALL domain_clock_set( nested_grid, & |
---|
453 | current_timestr=current_date(1:19), & |
---|
454 | stop_timestr=stopTimeStr , & |
---|
455 | time_step_seconds= & |
---|
456 | model_config_rec%interval_seconds ) |
---|
457 | |
---|
458 | ! Generate an output file from this program, which will be an input file to WRF. |
---|
459 | |
---|
460 | CALL nl_set_bdyfrq ( nested_grid%id , new_bdy_frq ) |
---|
461 | config_flags%bdyfrq = new_bdy_frq |
---|
462 | |
---|
463 | #ifdef WRF_CHEM |
---|
464 | nested_grid%chem_opt = parent_grid%chem_opt |
---|
465 | nested_grid%chem_in_opt = parent_grid%chem_in_opt |
---|
466 | #endif |
---|
467 | |
---|
468 | ! Initialize constants and 1d arrays in fine grid from the parent. |
---|
469 | |
---|
470 | CALL init_domain_constants_em_ptr ( parent_grid , nested_grid ) |
---|
471 | |
---|
472 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
473 | |
---|
474 | CALL wrf_debug ( 100 , 'ndown_em main: calling open_w_dataset for wrfinput' ) |
---|
475 | CALL construct_filename1( outname , 'wrfinput' , nested_grid%id , 2 ) |
---|
476 | CALL open_w_dataset ( fido, TRIM(outname) , nested_grid , config_flags , output_input , "DATASET=INPUT", ierr ) |
---|
477 | IF ( ierr .NE. 0 ) THEN |
---|
478 | WRITE( wrf_err_message , FMT='(A,A,A,I8)' ) 'program ndown: error opening ',TRIM(outname),' for reading ierr=',ierr |
---|
479 | CALL WRF_ERROR_FATAL ( wrf_err_message ) |
---|
480 | ENDIF |
---|
481 | |
---|
482 | ! Various sizes that we need to be concerned about. |
---|
483 | |
---|
484 | ids = nested_grid%sd31 |
---|
485 | ide = nested_grid%ed31 |
---|
486 | kds = nested_grid%sd32 |
---|
487 | kde = nested_grid%ed32 |
---|
488 | jds = nested_grid%sd33 |
---|
489 | jde = nested_grid%ed33 |
---|
490 | |
---|
491 | ims = nested_grid%sm31 |
---|
492 | ime = nested_grid%em31 |
---|
493 | kms = nested_grid%sm32 |
---|
494 | kme = nested_grid%em32 |
---|
495 | jms = nested_grid%sm33 |
---|
496 | jme = nested_grid%em33 |
---|
497 | |
---|
498 | ips = nested_grid%sp31 |
---|
499 | ipe = nested_grid%ep31 |
---|
500 | kps = nested_grid%sp32 |
---|
501 | kpe = nested_grid%ep32 |
---|
502 | jps = nested_grid%sp33 |
---|
503 | jpe = nested_grid%ep33 |
---|
504 | |
---|
505 | |
---|
506 | print *, ids , ide , jds , jde , kds , kde |
---|
507 | print *, ims , ime , jms , jme , kms , kme |
---|
508 | print *, ips , ipe , jps , jpe , kps , kpe |
---|
509 | |
---|
510 | spec_bdy_width = model_config_rec%spec_bdy_width |
---|
511 | print *,'spec_bdy_width=',spec_bdy_width |
---|
512 | |
---|
513 | ! This is the space needed to save the current 3d data for use in computing |
---|
514 | ! the lateral boundary tendencies. |
---|
515 | |
---|
516 | ALLOCATE ( ubdy3dtemp1(ims:ime,kms:kme,jms:jme) ) |
---|
517 | ALLOCATE ( vbdy3dtemp1(ims:ime,kms:kme,jms:jme) ) |
---|
518 | ALLOCATE ( tbdy3dtemp1(ims:ime,kms:kme,jms:jme) ) |
---|
519 | ALLOCATE ( pbdy3dtemp1(ims:ime,kms:kme,jms:jme) ) |
---|
520 | ALLOCATE ( qbdy3dtemp1(ims:ime,kms:kme,jms:jme) ) |
---|
521 | ALLOCATE ( mbdy2dtemp1(ims:ime,1:1, jms:jme) ) |
---|
522 | ALLOCATE ( ubdy3dtemp2(ims:ime,kms:kme,jms:jme) ) |
---|
523 | ALLOCATE ( vbdy3dtemp2(ims:ime,kms:kme,jms:jme) ) |
---|
524 | ALLOCATE ( tbdy3dtemp2(ims:ime,kms:kme,jms:jme) ) |
---|
525 | ALLOCATE ( pbdy3dtemp2(ims:ime,kms:kme,jms:jme) ) |
---|
526 | ALLOCATE ( qbdy3dtemp2(ims:ime,kms:kme,jms:jme) ) |
---|
527 | ALLOCATE ( mbdy2dtemp2(ims:ime,1:1, jms:jme) ) |
---|
528 | ALLOCATE ( cbdy3dtemp0(ims:ime,kms:kme,jms:jme,1:num_chem) ) |
---|
529 | ALLOCATE ( cbdy3dtemp1(ims:ime,kms:kme,jms:jme) ) |
---|
530 | ALLOCATE ( cbdy3dtemp2(ims:ime,kms:kme,jms:jme) ) |
---|
531 | |
---|
532 | END IF |
---|
533 | |
---|
534 | CALL domain_clock_set( nested_grid, & |
---|
535 | current_timestr=current_date(1:19), & |
---|
536 | time_step_seconds= & |
---|
537 | model_config_rec%interval_seconds ) |
---|
538 | |
---|
539 | ! Do the horizontal interpolation. |
---|
540 | |
---|
541 | nested_grid%imask_nostag = 1 |
---|
542 | nested_grid%imask_xstag = 1 |
---|
543 | nested_grid%imask_ystag = 1 |
---|
544 | nested_grid%imask_xystag = 1 |
---|
545 | |
---|
546 | CALL med_interp_domain ( head_grid , nested_grid ) |
---|
547 | WRITE ( message , * ) 'MOUSTA_L', nested_grid%mu_2(ips,jps) , nested_grid%u_2(ips,kde-2,jps) |
---|
548 | CALL wrf_debug ( 99,message ) |
---|
549 | CALL vertical_interp (nested_grid,znw_c,znu_c,cf1_c,cf2_c,cf3_c,cfn_c,cfn1_c,k_dim_c) |
---|
550 | WRITE ( message , * ) 'MOUSTA_V', nested_grid%mu_2(ips,jps) , nested_grid%u_2(ips,kde-2,jps) |
---|
551 | CALL wrf_debug ( 99,message ) |
---|
552 | nested_grid%ht_int = nested_grid%ht |
---|
553 | |
---|
554 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
555 | |
---|
556 | IF ( time_loop .EQ. 1 ) THEN |
---|
557 | |
---|
558 | ! Dimension info for fine grid. |
---|
559 | |
---|
560 | CALL get_ijk_from_grid ( nested_grid , & |
---|
561 | nids, nide, njds, njde, nkds, nkde, & |
---|
562 | nims, nime, njms, njme, nkms, nkme, & |
---|
563 | nips, nipe, njps, njpe, nkps, nkpe ) |
---|
564 | |
---|
565 | ! Store horizontally interpolated terrain in temp location |
---|
566 | |
---|
567 | CALL copy_3d_field ( nested_grid%ht_fine , nested_grid%ht , & |
---|
568 | nids , nide , njds , njde , 1 , 1 , & |
---|
569 | nims , nime , njms , njme , 1 , 1 , & |
---|
570 | nips , nipe , njps , njpe , 1 , 1 ) |
---|
571 | |
---|
572 | ! Open the fine grid SI static file. |
---|
573 | |
---|
574 | CALL construct_filename1( si_inpname , 'wrfndi' , nested_grid%id , 2 ) |
---|
575 | CALL wrf_debug ( 100 , 'med_sidata_input: calling open_r_dataset for ' // TRIM(si_inpname) ) |
---|
576 | CALL open_r_dataset ( idsi, TRIM(si_inpname) , nested_grid , config_flags , "DATASET=INPUT", ierr ) |
---|
577 | IF ( ierr .NE. 0 ) THEN |
---|
578 | CALL wrf_error_fatal( 'real: error opening FG input for reading: ' // TRIM (si_inpname) ) |
---|
579 | END IF |
---|
580 | |
---|
581 | ! Input data. |
---|
582 | |
---|
583 | CALL wrf_debug ( 100 , 'ndown_em: calling input_auxinput2' ) |
---|
584 | CALL input_auxinput2 ( idsi , nested_grid , config_flags , ierr ) |
---|
585 | nested_grid%ht_input = nested_grid%ht |
---|
586 | |
---|
587 | ! Close this fine grid static input file. |
---|
588 | |
---|
589 | CALL wrf_debug ( 100 , 'ndown_em: closing fine grid static input' ) |
---|
590 | CALL close_dataset ( idsi , config_flags , "DATASET=INPUT" ) |
---|
591 | |
---|
592 | ! Blend parent and nest field of terrain. |
---|
593 | |
---|
594 | CALL blend_terrain ( nested_grid%ht_fine , nested_grid%ht , & |
---|
595 | nids , nide , njds , njde , 1 , 1 , & |
---|
596 | nims , nime , njms , njme , 1 , 1 , & |
---|
597 | nips , nipe , njps , njpe , 1 , 1 ) |
---|
598 | |
---|
599 | nested_grid%ht_input = nested_grid%ht |
---|
600 | nested_grid%ht_int = nested_grid%ht_fine |
---|
601 | |
---|
602 | ! We need a fine grid landuse in the interpolation. So we need to generate |
---|
603 | ! that field now. |
---|
604 | |
---|
605 | IF ( ( nested_grid%ivgtyp(ips,jps) .GT. 0 ) .AND. & |
---|
606 | ( nested_grid%isltyp(ips,jps) .GT. 0 ) ) THEN |
---|
607 | DO j = jps, MIN(jde-1,jpe) |
---|
608 | DO i = ips, MIN(ide-1,ipe) |
---|
609 | nested_grid% vegcat(i,j) = nested_grid%ivgtyp(i,j) |
---|
610 | nested_grid%soilcat(i,j) = nested_grid%isltyp(i,j) |
---|
611 | END DO |
---|
612 | END DO |
---|
613 | |
---|
614 | ELSE IF ( ( nested_grid% vegcat(ips,jps) .GT. 0.5 ) .AND. & |
---|
615 | ( nested_grid%soilcat(ips,jps) .GT. 0.5 ) ) THEN |
---|
616 | DO j = jps, MIN(jde-1,jpe) |
---|
617 | DO i = ips, MIN(ide-1,ipe) |
---|
618 | nested_grid%ivgtyp(i,j) = NINT(nested_grid% vegcat(i,j)) |
---|
619 | nested_grid%isltyp(i,j) = NINT(nested_grid%soilcat(i,j)) |
---|
620 | END DO |
---|
621 | END DO |
---|
622 | |
---|
623 | ELSE |
---|
624 | num_veg_cat = SIZE ( nested_grid%landusef , DIM=2 ) |
---|
625 | num_soil_top_cat = SIZE ( nested_grid%soilctop , DIM=2 ) |
---|
626 | num_soil_bot_cat = SIZE ( nested_grid%soilcbot , DIM=2 ) |
---|
627 | |
---|
628 | CALL land_percentages ( nested_grid%xland , & |
---|
629 | nested_grid%landusef , nested_grid%soilctop , nested_grid%soilcbot , & |
---|
630 | nested_grid%isltyp , nested_grid%ivgtyp , & |
---|
631 | num_veg_cat , num_soil_top_cat , num_soil_bot_cat , & |
---|
632 | ids , ide , jds , jde , kds , kde , & |
---|
633 | ims , ime , jms , jme , kms , kme , & |
---|
634 | ips , ipe , jps , jpe , kps , kpe , & |
---|
635 | model_config_rec%iswater(nested_grid%id) ) |
---|
636 | |
---|
637 | END IF |
---|
638 | |
---|
639 | DO j = jps, MIN(jde-1,jpe) |
---|
640 | DO i = ips, MIN(ide-1,ipe) |
---|
641 | nested_grid%lu_index(i,j) = nested_grid%ivgtyp(i,j) |
---|
642 | END DO |
---|
643 | END DO |
---|
644 | |
---|
645 | #ifndef PLANET |
---|
646 | CALL check_consistency ( nested_grid%ivgtyp , nested_grid%isltyp , nested_grid%landmask , & |
---|
647 | ids , ide , jds , jde , kds , kde , & |
---|
648 | ims , ime , jms , jme , kms , kme , & |
---|
649 | ips , ipe , jps , jpe , kps , kpe , & |
---|
650 | model_config_rec%iswater(nested_grid%id) ) |
---|
651 | |
---|
652 | CALL check_consistency2( nested_grid%ivgtyp , nested_grid%isltyp , nested_grid%landmask , & |
---|
653 | nested_grid%tmn , nested_grid%tsk , nested_grid%sst , nested_grid%xland , & |
---|
654 | nested_grid%tslb , nested_grid%smois , nested_grid%sh2o , & |
---|
655 | config_flags%num_soil_layers , nested_grid%id , & |
---|
656 | ids , ide , jds , jde , kds , kde , & |
---|
657 | ims , ime , jms , jme , kms , kme , & |
---|
658 | ips , ipe , jps , jpe , kps , kpe , & |
---|
659 | model_config_rec%iswater(nested_grid%id) ) |
---|
660 | #endif |
---|
661 | |
---|
662 | END IF |
---|
663 | |
---|
664 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
665 | |
---|
666 | ! We have 2 terrain elevations. One is from input and the other is from the |
---|
667 | ! the horizontal interpolation. |
---|
668 | |
---|
669 | nested_grid%ht_fine = nested_grid%ht_input |
---|
670 | nested_grid%ht = nested_grid%ht_int |
---|
671 | |
---|
672 | ! We have both the interpolated fields and the higher-resolution static fields. From these |
---|
673 | ! the rebalancing is now done. Note also that the field nested_grid%ht is now from the |
---|
674 | ! fine grid input file (after this call is completed). |
---|
675 | |
---|
676 | CALL rebalance_driver ( nested_grid ) |
---|
677 | |
---|
678 | ! Different things happen during the different time loops: |
---|
679 | ! first loop - write wrfinput file, close data set, copy files to holder arrays |
---|
680 | ! middle loops - diff 3d/2d arrays, compute and output bc |
---|
681 | ! last loop - diff 3d/2d arrays, compute and output bc, write wrfbdy file, close wrfbdy file |
---|
682 | |
---|
683 | IF ( time_loop .EQ. 1 ) THEN |
---|
684 | |
---|
685 | ! Set the time info. |
---|
686 | |
---|
687 | print *,'current_date = ',current_date |
---|
688 | CALL domain_clock_set( nested_grid, & |
---|
689 | current_timestr=current_date(1:19) ) |
---|
690 | #ifdef WRF_CHEM |
---|
691 | ! |
---|
692 | ! Put in chemistry data |
---|
693 | ! |
---|
694 | IF( nested_grid%chem_opt .NE. 0 ) then |
---|
695 | ! IF( nested_grid%chem_in_opt .EQ. 0 ) then |
---|
696 | ! Read the chemistry data from a previous wrf forecast (wrfout file) |
---|
697 | ! Generate chemistry data from a idealized vertical profile |
---|
698 | ! message = 'STARTING WITH BACKGROUND CHEMISTRY ' |
---|
699 | CALL wrf_message ( message ) |
---|
700 | |
---|
701 | ! CALL input_chem_profile ( nested_grid ) |
---|
702 | |
---|
703 | if(nested_grid%biomass_burn_opt == BIOMASSB) THEN |
---|
704 | message = 'READING BIOMASS BURNING EMISSIONS DATA ' |
---|
705 | CALL wrf_message ( message ) |
---|
706 | CALL med_read_wrf_chem_emissopt3 ( nested_grid , config_flags) |
---|
707 | end if |
---|
708 | |
---|
709 | if(nested_grid%dust_opt == 1 .or. nested_grid%dmsemis_opt == 1 & |
---|
710 | .or. nested_grid%chem_opt == 300 .or. nested_grid%chem_opt == 301) then |
---|
711 | message = 'READING GOCART BG AND/OR DUST and DMS REF FIELDS' |
---|
712 | CALL wrf_message ( message ) |
---|
713 | CALL med_read_wrf_chem_gocart_bg ( nested_grid , config_flags) |
---|
714 | end if |
---|
715 | |
---|
716 | if( nested_grid%bio_emiss_opt .eq. 2 )then |
---|
717 | message = 'READING BEIS3.11 EMISSIONS DATA' |
---|
718 | CALL wrf_message ( message ) |
---|
719 | CALL med_read_wrf_chem_bioemiss ( nested_grid , config_flags) |
---|
720 | else if( nested_grid%bio_emiss_opt == 3 ) THEN |
---|
721 | message = 'READING MEGAN 2 EMISSIONS DATA' |
---|
722 | CALL wrf_message ( message ) |
---|
723 | CALL med_read_wrf_chem_bioemiss ( nested_grid , config_flags) |
---|
724 | endif |
---|
725 | ! ELSE |
---|
726 | ! message = 'RUNNING WITHOUT CHEMISTRY INITIALIZATION' |
---|
727 | ! CALL wrf_message ( message ) |
---|
728 | ! ENDIF |
---|
729 | ENDIF |
---|
730 | #endif |
---|
731 | |
---|
732 | ! Output the first time period of the data. |
---|
733 | |
---|
734 | CALL output_input ( fido , nested_grid , config_flags , ierr ) |
---|
735 | |
---|
736 | CALL wrf_put_dom_ti_integer ( fido , 'MAP_PROJ' , map_proj , 1 , ierr ) |
---|
737 | ! CALL wrf_put_dom_ti_real ( fido , 'DX' , dx , 1 , ierr ) |
---|
738 | ! CALL wrf_put_dom_ti_real ( fido , 'DY' , dy , 1 , ierr ) |
---|
739 | CALL wrf_put_dom_ti_real ( fido , 'CEN_LAT' , cen_lat , 1 , ierr ) |
---|
740 | CALL wrf_put_dom_ti_real ( fido , 'CEN_LON' , cen_lon , 1 , ierr ) |
---|
741 | CALL wrf_put_dom_ti_real ( fido , 'TRUELAT1' , truelat1 , 1 , ierr ) |
---|
742 | CALL wrf_put_dom_ti_real ( fido , 'TRUELAT2' , truelat2 , 1 , ierr ) |
---|
743 | CALL wrf_put_dom_ti_real ( fido , 'MOAD_CEN_LAT' , moad_cen_lat , 1 , ierr ) |
---|
744 | CALL wrf_put_dom_ti_real ( fido , 'STAND_LON' , stand_lon , 1 , ierr ) |
---|
745 | CALL wrf_put_dom_ti_integer ( fido , 'ISWATER' , iswater , 1 , ierr ) |
---|
746 | |
---|
747 | ! These change if the initial time for the nest is not the same as the |
---|
748 | ! first time period in the WRF output file. |
---|
749 | ! Now that we know the starting date, we need to set the GMT, JULYR, and JULDAY |
---|
750 | ! values for the global attributes. This call is based on the setting of the |
---|
751 | ! current_date string. |
---|
752 | |
---|
753 | CALL geth_julgmt ( julyr , julday , gmt) |
---|
754 | CALL nl_set_julyr ( nested_grid%id , julyr ) |
---|
755 | CALL nl_set_julday ( nested_grid%id , julday ) |
---|
756 | CALL nl_set_gmt ( nested_grid%id , gmt ) |
---|
757 | CALL wrf_put_dom_ti_real ( fido , 'GMT' , gmt , 1 , ierr ) |
---|
758 | CALL wrf_put_dom_ti_integer ( fido , 'JULYR' , julyr , 1 , ierr ) |
---|
759 | CALL wrf_put_dom_ti_integer ( fido , 'JULDAY' , julday , 1 , ierr ) |
---|
760 | print *,'current_date =',current_date |
---|
761 | print *,'julyr=',julyr |
---|
762 | print *,'julday=',julday |
---|
763 | print *,'gmt=',gmt |
---|
764 | |
---|
765 | ! Close the input (wrfout_d01_000000, for example) file. That's right, the |
---|
766 | ! input is an output file. Who'd've thunk. |
---|
767 | |
---|
768 | CALL close_dataset ( fido , config_flags , "DATASET=INPUT" ) |
---|
769 | |
---|
770 | ! We need to save the 3d/2d data to compute a difference during the next loop. Couple the |
---|
771 | ! 3d fields with total mu (mub + mu_2) and the stagger-specific map scale factor. |
---|
772 | |
---|
773 | ! u, theta, h, scalars coupled with my, v coupled with mx |
---|
774 | CALL couple ( nested_grid%mu_2 , nested_grid%mub , ubdy3dtemp1 , nested_grid%u_2 , & |
---|
775 | 'u' , nested_grid%msfuy , & |
---|
776 | ids, ide, jds, jde, kds, kde, ims, ime, jms, jme, kms, kme, ips, ipe, jps, jpe, kps, kpe ) |
---|
777 | CALL couple ( nested_grid%mu_2 , nested_grid%mub , vbdy3dtemp1 , nested_grid%v_2 , & |
---|
778 | 'v' , nested_grid%msfvx , & |
---|
779 | ids, ide, jds, jde, kds, kde, ims, ime, jms, jme, kms, kme, ips, ipe, jps, jpe, kps, kpe ) |
---|
780 | CALL couple ( nested_grid%mu_2 , nested_grid%mub , tbdy3dtemp1 , nested_grid%t_2 , & |
---|
781 | 't' , nested_grid%msfty , & |
---|
782 | ids, ide, jds, jde, kds, kde, ims, ime, jms, jme, kms, kme, ips, ipe, jps, jpe, kps, kpe ) |
---|
783 | CALL couple ( nested_grid%mu_2 , nested_grid%mub , pbdy3dtemp1 , nested_grid%ph_2 , & |
---|
784 | 'h' , nested_grid%msfty , & |
---|
785 | ids, ide, jds, jde, kds, kde, ims, ime, jms, jme, kms, kme, ips, ipe, jps, jpe, kps, kpe ) |
---|
786 | CALL couple ( nested_grid%mu_2 , nested_grid%mub , qbdy3dtemp1 , nested_grid%moist(ims:ime,kms:kme,jms:jme,P_QV) , & |
---|
787 | 't' , nested_grid%msfty , & |
---|
788 | ids, ide, jds, jde, kds, kde, ims, ime, jms, jme, kms, kme, ips, ipe, jps, jpe, kps, kpe ) |
---|
789 | |
---|
790 | DO j = jps , jpe |
---|
791 | DO i = ips , ipe |
---|
792 | mbdy2dtemp1(i,1,j) = nested_grid%mu_2(i,j) |
---|
793 | END DO |
---|
794 | END DO |
---|
795 | |
---|
796 | ! There are 2 components to the lateral boundaries. First, there is the starting |
---|
797 | ! point of this time period - just the outer few rows and columns. |
---|
798 | |
---|
799 | CALL stuff_bdy ( ubdy3dtemp1 , nested_grid%u_bxs, nested_grid%u_bxe, & |
---|
800 | nested_grid%u_bys, nested_grid%u_bye, & |
---|
801 | 'U' , spec_bdy_width , & |
---|
802 | ids , ide , jds , jde , kds , kde , & |
---|
803 | ims , ime , jms , jme , kms , kme , & |
---|
804 | ips , ipe , jps , jpe , kps , kpe ) |
---|
805 | CALL stuff_bdy ( vbdy3dtemp1 , nested_grid%v_bxs, nested_grid%v_bxe, & |
---|
806 | nested_grid%v_bys, nested_grid%v_bye, & |
---|
807 | 'V' , spec_bdy_width , & |
---|
808 | ids , ide , jds , jde , kds , kde , & |
---|
809 | ims , ime , jms , jme , kms , kme , & |
---|
810 | ips , ipe , jps , jpe , kps , kpe ) |
---|
811 | CALL stuff_bdy ( tbdy3dtemp1 , nested_grid%t_bxs, nested_grid%t_bxe, & |
---|
812 | nested_grid%t_bys, nested_grid%t_bye, & |
---|
813 | 'T' , spec_bdy_width , & |
---|
814 | ids , ide , jds , jde , kds , kde , & |
---|
815 | ims , ime , jms , jme , kms , kme , & |
---|
816 | ips , ipe , jps , jpe , kps , kpe ) |
---|
817 | CALL stuff_bdy ( pbdy3dtemp1 , nested_grid%ph_bxs, nested_grid%ph_bxe, & |
---|
818 | nested_grid%ph_bys, nested_grid%ph_bye, & |
---|
819 | 'W' , spec_bdy_width , & |
---|
820 | ids , ide , jds , jde , kds , kde , & |
---|
821 | ims , ime , jms , jme , kms , kme , & |
---|
822 | ips , ipe , jps , jpe , kps , kpe ) |
---|
823 | CALL stuff_bdy ( qbdy3dtemp1 , nested_grid%moist_bxs(jms:jme,kms:kme,1:spec_bdy_width,P_QV), & |
---|
824 | nested_grid%moist_bxe(jms:jme,kms:kme,1:spec_bdy_width,P_QV), & |
---|
825 | nested_grid%moist_bys(ims:ime,kms:kme,1:spec_bdy_width,P_QV), & |
---|
826 | nested_grid%moist_bye(ims:ime,kms:kme,1:spec_bdy_width,P_QV), & |
---|
827 | 'T' , spec_bdy_width , & |
---|
828 | ids , ide , jds , jde , kds , kde , & |
---|
829 | ims , ime , jms , jme , kms , kme , & |
---|
830 | ips , ipe , jps , jpe , kps , kpe ) |
---|
831 | CALL stuff_bdy ( mbdy2dtemp1 , nested_grid%mu_bxs, nested_grid%mu_bxe, & |
---|
832 | nested_grid%mu_bys, nested_grid%mu_bye, & |
---|
833 | 'M' , spec_bdy_width , & |
---|
834 | ids , ide , jds , jde , 1 , 1 , & |
---|
835 | ims , ime , jms , jme , 1 , 1 , & |
---|
836 | ips , ipe , jps , jpe , 1 , 1 ) |
---|
837 | #ifdef WRF_CHEM |
---|
838 | do nvchem=1,num_chem |
---|
839 | ! if(nvchem.eq.p_o3)then |
---|
840 | ! write(0,*)'fill ch_b',cbdy3dtemp1(5,1,5),nvchem |
---|
841 | ! endif |
---|
842 | cbdy3dtemp1(ips:ipe,kps:kpe,jps:jpe)=nested_grid%chem(ips:ipe,kps:kpe,jps:jpe,nvchem) |
---|
843 | ! if(nvchem.eq.p_o3)then |
---|
844 | ! write(0,*)'fill ch_b',cbdy3dtemp1(5,1,5) |
---|
845 | ! endif |
---|
846 | CALL stuff_bdy ( cbdy3dtemp1 , nested_grid%chem_bxs(jms:jme,kms:kme,1:spec_bdy_width,nvchem), & |
---|
847 | nested_grid%chem_bxe(jms:jme,kms:kme,1:spec_bdy_width,nvchem), & |
---|
848 | nested_grid%chem_bys(ims:ime,kms:kme,1:spec_bdy_width,nvchem), & |
---|
849 | nested_grid%chem_bye(ims:ime,kms:kme,1:spec_bdy_width,nvchem), & |
---|
850 | 'T' , spec_bdy_width , & |
---|
851 | ids , ide , jds , jde , kds , kde , & |
---|
852 | ims , ime , jms , jme , kms , kme , & |
---|
853 | ips , ipe , jps , jpe , kps , kpe ) |
---|
854 | cbdy3dtemp0(ips:ipe,kps:kpe,jps:jpe,nvchem)=cbdy3dtemp1(ips:ipe,kps:kpe,jps:jpe) |
---|
855 | ! if(nvchem.eq.p_o3)then |
---|
856 | ! write(0,*)'filled ch_b',time_loop,cbdy3dtemp1(ips,kps,jps),cbdy3dtemp0(ips,kps,jps,nvchem) |
---|
857 | ! endif |
---|
858 | enddo |
---|
859 | #endif |
---|
860 | ELSE IF ( ( time_loop .GT. 1 ) .AND. ( time_loop .LT. time_loop_max ) ) THEN |
---|
861 | |
---|
862 | ! u, theta, h, scalars coupled with my, v coupled with mx |
---|
863 | CALL couple ( nested_grid%mu_2 , nested_grid%mub , ubdy3dtemp2 , nested_grid%u_2 , & |
---|
864 | 'u' , nested_grid%msfuy , & |
---|
865 | ids, ide, jds, jde, kds, kde, ims, ime, jms, jme, kms, kme, ips, ipe, jps, jpe, kps, kpe ) |
---|
866 | CALL couple ( nested_grid%mu_2 , nested_grid%mub , vbdy3dtemp2 , nested_grid%v_2 , & |
---|
867 | 'v' , nested_grid%msfvx , & |
---|
868 | ids, ide, jds, jde, kds, kde, ims, ime, jms, jme, kms, kme, ips, ipe, jps, jpe, kps, kpe ) |
---|
869 | CALL couple ( nested_grid%mu_2 , nested_grid%mub , tbdy3dtemp2 , nested_grid%t_2 , & |
---|
870 | 't' , nested_grid%msfty , & |
---|
871 | ids, ide, jds, jde, kds, kde, ims, ime, jms, jme, kms, kme, ips, ipe, jps, jpe, kps, kpe ) |
---|
872 | CALL couple ( nested_grid%mu_2 , nested_grid%mub , pbdy3dtemp2 , nested_grid%ph_2 , & |
---|
873 | 'h' , nested_grid%msfty , & |
---|
874 | ids, ide, jds, jde, kds, kde, ims, ime, jms, jme, kms, kme, ips, ipe, jps, jpe, kps, kpe ) |
---|
875 | CALL couple ( nested_grid%mu_2 , nested_grid%mub , qbdy3dtemp2 , nested_grid%moist(ims:ime,kms:kme,jms:jme,P_QV) , & |
---|
876 | 't' , nested_grid%msfty , & |
---|
877 | ids, ide, jds, jde, kds, kde, ims, ime, jms, jme, kms, kme, ips, ipe, jps, jpe, kps, kpe ) |
---|
878 | |
---|
879 | DO j = jps , jpe |
---|
880 | DO i = ips , ipe |
---|
881 | mbdy2dtemp2(i,1,j) = nested_grid%mu_2(i,j) |
---|
882 | END DO |
---|
883 | END DO |
---|
884 | |
---|
885 | ! During all of the loops after the first loop, we first compute the boundary |
---|
886 | ! tendencies with the current data values and the previously save information |
---|
887 | ! stored in the *bdy3dtemp1 arrays. |
---|
888 | |
---|
889 | CALL stuff_bdytend ( ubdy3dtemp2 , ubdy3dtemp1 , new_bdy_frq , & |
---|
890 | nested_grid%u_btxs, nested_grid%u_btxe , & |
---|
891 | nested_grid%u_btys, nested_grid%u_btye , & |
---|
892 | 'U' , & |
---|
893 | spec_bdy_width , & |
---|
894 | ids , ide , jds , jde , kds , kde , & |
---|
895 | ims , ime , jms , jme , kms , kme , & |
---|
896 | ips , ipe , jps , jpe , kps , kpe ) |
---|
897 | CALL stuff_bdytend ( vbdy3dtemp2 , vbdy3dtemp1 , new_bdy_frq , & |
---|
898 | nested_grid%v_btxs, nested_grid%v_btxe , & |
---|
899 | nested_grid%v_btys, nested_grid%v_btye , & |
---|
900 | 'V' , & |
---|
901 | spec_bdy_width , & |
---|
902 | ids , ide , jds , jde , kds , kde , & |
---|
903 | ims , ime , jms , jme , kms , kme , & |
---|
904 | ips , ipe , jps , jpe , kps , kpe ) |
---|
905 | CALL stuff_bdytend ( tbdy3dtemp2 , tbdy3dtemp1 , new_bdy_frq , & |
---|
906 | nested_grid%t_btxs, nested_grid%t_btxe , & |
---|
907 | nested_grid%t_btys, nested_grid%t_btye , & |
---|
908 | 'T' , & |
---|
909 | spec_bdy_width , & |
---|
910 | ids , ide , jds , jde , kds , kde , & |
---|
911 | ims , ime , jms , jme , kms , kme , & |
---|
912 | ips , ipe , jps , jpe , kps , kpe ) |
---|
913 | CALL stuff_bdytend ( pbdy3dtemp2 , pbdy3dtemp1 , new_bdy_frq , & |
---|
914 | nested_grid%ph_btxs, nested_grid%ph_btxe , & |
---|
915 | nested_grid%ph_btys, nested_grid%ph_btye , & |
---|
916 | 'W' , & |
---|
917 | spec_bdy_width , & |
---|
918 | ids , ide , jds , jde , kds , kde , & |
---|
919 | ims , ime , jms , jme , kms , kme , & |
---|
920 | ips , ipe , jps , jpe , kps , kpe ) |
---|
921 | CALL stuff_bdytend ( qbdy3dtemp2 , qbdy3dtemp1 , new_bdy_frq , & |
---|
922 | nested_grid%moist_btxs(jms:jme,kms:kme,1:spec_bdy_width,P_QV), & |
---|
923 | nested_grid%moist_btxe(jms:jme,kms:kme,1:spec_bdy_width,P_QV), & |
---|
924 | nested_grid%moist_btys(ims:ime,kms:kme,1:spec_bdy_width,P_QV), & |
---|
925 | nested_grid%moist_btye(ims:ime,kms:kme,1:spec_bdy_width,P_QV), & |
---|
926 | 'T' , & |
---|
927 | spec_bdy_width , & |
---|
928 | ids , ide , jds , jde , kds , kde , & |
---|
929 | ims , ime , jms , jme , kms , kme , & |
---|
930 | ips , ipe , jps , jpe , kps , kpe ) |
---|
931 | CALL stuff_bdytend ( mbdy2dtemp2 , mbdy2dtemp1 , new_bdy_frq , & |
---|
932 | nested_grid%mu_btxs, nested_grid%mu_btxe , & |
---|
933 | nested_grid%mu_btys, nested_grid%mu_btye , & |
---|
934 | 'M' , & |
---|
935 | spec_bdy_width , & |
---|
936 | ids , ide , jds , jde , 1 , 1 , & |
---|
937 | ims , ime , jms , jme , 1 , 1 , & |
---|
938 | ips , ipe , jps , jpe , 1 , 1 ) |
---|
939 | #ifdef WRF_CHEM |
---|
940 | do nvchem=1,num_chem |
---|
941 | cbdy3dtemp1(ips:ipe,kps:kpe,jps:jpe)=cbdy3dtemp0(ips:ipe,kps:kpe,jps:jpe,nvchem) |
---|
942 | cbdy3dtemp2(ips:ipe,kps:kpe,jps:jpe)=nested_grid%chem(ips:ipe,kps:kpe,jps:jpe,nvchem) |
---|
943 | ! if(nvchem.eq.p_o3)then |
---|
944 | ! write(0,*)'fill 1ch_b2',time_loop,cbdy3dtemp1(ips,kps,jps),cbdy3dtemp0(ips,kps,jps,nvchem),cbdy3dtemp2(ips,kps,jps) |
---|
945 | ! endif |
---|
946 | CALL stuff_bdytend ( cbdy3dtemp2 , cbdy3dtemp1 , new_bdy_frq , & |
---|
947 | nested_grid%chem_btxs(jms:jme,kms:kme,1:spec_bdy_width,nvchem), & |
---|
948 | nested_grid%chem_btxe(jms:jme,kms:kme,1:spec_bdy_width,nvchem), & |
---|
949 | nested_grid%chem_btys(ims:ime,kms:kme,1:spec_bdy_width,nvchem), & |
---|
950 | nested_grid%chem_btye(ims:ime,kms:kme,1:spec_bdy_width,nvchem), & |
---|
951 | 'T' , & |
---|
952 | spec_bdy_width , & |
---|
953 | ids , ide , jds , jde , kds , kde , & |
---|
954 | ims , ime , jms , jme , kms , kme , & |
---|
955 | ips , ipe , jps , jpe , kps , kpe ) |
---|
956 | cbdy3dtemp0(ips:ipe,kps:kpe,jps:jpe,nvchem)=cbdy3dtemp2(ips:ipe,kps:kpe,jps:jpe) |
---|
957 | ! if(nvchem.eq.p_o3)then |
---|
958 | ! write(0,*)'fill 2ch_b2',cbdy3dtemp1(ips,kps,jps),cbdy3dtemp0(ips,kps,jps,nvchem),cbdy3dtemp2(ips,kps,jps) |
---|
959 | ! endif |
---|
960 | enddo |
---|
961 | #endif |
---|
962 | IF ( time_loop .EQ. 2 ) THEN |
---|
963 | |
---|
964 | ! Generate an output file from this program, which will be an input file to WRF. |
---|
965 | |
---|
966 | CALL wrf_debug ( 100 , 'ndown_em main: calling open_w_dataset for wrfbdy' ) |
---|
967 | CALL construct_filename1( bdyname , 'wrfbdy' , nested_grid%id , 2 ) |
---|
968 | CALL open_w_dataset ( fidb, TRIM(bdyname) , nested_grid , config_flags , output_boundary , & |
---|
969 | "DATASET=BOUNDARY", ierr ) |
---|
970 | IF ( ierr .NE. 0 ) THEN |
---|
971 | WRITE( wrf_err_message , FMT='(A,A,A,I8)' ) 'program ndown: error opening ',TRIM(bdyname),' for reading ierr=',ierr |
---|
972 | CALL WRF_ERROR_FATAL ( wrf_err_message ) |
---|
973 | ENDIF |
---|
974 | |
---|
975 | END IF |
---|
976 | |
---|
977 | ! Both pieces of the boundary data are now available to be written. |
---|
978 | |
---|
979 | CALL domain_clock_set( nested_grid, & |
---|
980 | current_timestr=current_date(1:19) ) |
---|
981 | temp24= current_date |
---|
982 | temp24b=start_date_hold |
---|
983 | start_date = start_date_hold |
---|
984 | CALL geth_newdate ( temp19 , temp24b(1:19) , (time_loop-2) * model_config_rec%interval_seconds ) |
---|
985 | current_date = temp19 // '.0000' |
---|
986 | CALL geth_julgmt ( julyr , julday , gmt) |
---|
987 | CALL nl_set_julyr ( nested_grid%id , julyr ) |
---|
988 | CALL nl_set_julday ( nested_grid%id , julday ) |
---|
989 | CALL nl_set_gmt ( nested_grid%id , gmt ) |
---|
990 | CALL wrf_put_dom_ti_real ( fidb , 'GMT' , gmt , 1 , ierr ) |
---|
991 | CALL wrf_put_dom_ti_integer ( fidb , 'JULYR' , julyr , 1 , ierr ) |
---|
992 | CALL wrf_put_dom_ti_integer ( fidb , 'JULDAY' , julday , 1 , ierr ) |
---|
993 | CALL domain_clock_set( nested_grid, & |
---|
994 | current_timestr=current_date(1:19) ) |
---|
995 | print *,'bdy time = ',time_loop-1,' bdy date = ',current_date,' ',start_date |
---|
996 | CALL output_boundary ( fidb , nested_grid , config_flags , ierr ) |
---|
997 | current_date = temp24 |
---|
998 | start_date = temp24b |
---|
999 | CALL domain_clock_set( nested_grid, & |
---|
1000 | current_timestr=current_date(1:19) ) |
---|
1001 | |
---|
1002 | IF ( time_loop .EQ. 2 ) THEN |
---|
1003 | CALL wrf_put_dom_ti_real ( fidb , 'BDYFRQ' , new_bdy_frq , 1 , ierr ) |
---|
1004 | END IF |
---|
1005 | |
---|
1006 | ! We need to save the 3d data to compute a difference during the next loop. Couple the |
---|
1007 | ! 3d fields with total mu (mub + mu_2) and the stagger-specific map scale factor. |
---|
1008 | ! We load up the boundary data again for use in the next loop. |
---|
1009 | |
---|
1010 | DO j = jps , jpe |
---|
1011 | DO k = kps , kpe |
---|
1012 | DO i = ips , ipe |
---|
1013 | ubdy3dtemp1(i,k,j) = ubdy3dtemp2(i,k,j) |
---|
1014 | vbdy3dtemp1(i,k,j) = vbdy3dtemp2(i,k,j) |
---|
1015 | tbdy3dtemp1(i,k,j) = tbdy3dtemp2(i,k,j) |
---|
1016 | pbdy3dtemp1(i,k,j) = pbdy3dtemp2(i,k,j) |
---|
1017 | qbdy3dtemp1(i,k,j) = qbdy3dtemp2(i,k,j) |
---|
1018 | END DO |
---|
1019 | END DO |
---|
1020 | END DO |
---|
1021 | |
---|
1022 | DO j = jps , jpe |
---|
1023 | DO i = ips , ipe |
---|
1024 | mbdy2dtemp1(i,1,j) = mbdy2dtemp2(i,1,j) |
---|
1025 | END DO |
---|
1026 | END DO |
---|
1027 | |
---|
1028 | ! There are 2 components to the lateral boundaries. First, there is the starting |
---|
1029 | ! point of this time period - just the outer few rows and columns. |
---|
1030 | |
---|
1031 | CALL stuff_bdy ( ubdy3dtemp1 , & |
---|
1032 | nested_grid%u_bxs, nested_grid%u_bxe , & |
---|
1033 | nested_grid%u_bys, nested_grid%u_bye , & |
---|
1034 | 'U' , spec_bdy_width , & |
---|
1035 | ids , ide , jds , jde , kds , kde , & |
---|
1036 | ims , ime , jms , jme , kms , kme , & |
---|
1037 | ips , ipe , jps , jpe , kps , kpe ) |
---|
1038 | CALL stuff_bdy ( vbdy3dtemp1 , & |
---|
1039 | nested_grid%v_bxs, nested_grid%v_bxe , & |
---|
1040 | nested_grid%v_bys, nested_grid%v_bye , & |
---|
1041 | 'V' , spec_bdy_width , & |
---|
1042 | ids , ide , jds , jde , kds , kde , & |
---|
1043 | ims , ime , jms , jme , kms , kme , & |
---|
1044 | ips , ipe , jps , jpe , kps , kpe ) |
---|
1045 | CALL stuff_bdy ( tbdy3dtemp1 , & |
---|
1046 | nested_grid%t_bxs, nested_grid%t_bxe , & |
---|
1047 | nested_grid%t_bys, nested_grid%t_bye , & |
---|
1048 | 'T' , spec_bdy_width , & |
---|
1049 | ids , ide , jds , jde , kds , kde , & |
---|
1050 | ims , ime , jms , jme , kms , kme , & |
---|
1051 | ips , ipe , jps , jpe , kps , kpe ) |
---|
1052 | CALL stuff_bdy ( pbdy3dtemp1 , & |
---|
1053 | nested_grid%ph_bxs, nested_grid%ph_bxe , & |
---|
1054 | nested_grid%ph_bys, nested_grid%ph_bye , & |
---|
1055 | 'W' , spec_bdy_width , & |
---|
1056 | ids , ide , jds , jde , kds , kde , & |
---|
1057 | ims , ime , jms , jme , kms , kme , & |
---|
1058 | ips , ipe , jps , jpe , kps , kpe ) |
---|
1059 | CALL stuff_bdy ( qbdy3dtemp1 , & |
---|
1060 | nested_grid%moist_bxs(jms:jme,kms:kme,1:spec_bdy_width,P_QV), & |
---|
1061 | nested_grid%moist_bxe(jms:jme,kms:kme,1:spec_bdy_width,P_QV), & |
---|
1062 | nested_grid%moist_bys(ims:ime,kms:kme,1:spec_bdy_width,P_QV), & |
---|
1063 | nested_grid%moist_bye(ims:ime,kms:kme,1:spec_bdy_width,P_QV), & |
---|
1064 | 'T' , spec_bdy_width , & |
---|
1065 | ids , ide , jds , jde , kds , kde , & |
---|
1066 | ims , ime , jms , jme , kms , kme , & |
---|
1067 | ips , ipe , jps , jpe , kps , kpe ) |
---|
1068 | #ifdef WRF_CHEM |
---|
1069 | do nvchem=1,num_chem |
---|
1070 | cbdy3dtemp1(ips:ipe,kps:kpe,jps:jpe)=cbdy3dtemp0(ips:ipe,kps:kpe,jps:jpe,nvchem) |
---|
1071 | ! if(nvchem.eq.p_o3)then |
---|
1072 | ! write(0,*)'fill 2ch_b3',cbdy3dtemp1(ips,kps,jps),cbdy3dtemp0(ips,kps,jps,nvchem) |
---|
1073 | ! endif |
---|
1074 | CALL stuff_bdy ( cbdy3dtemp1 , & |
---|
1075 | nested_grid%chem_bxs(jms:jme,kms:kme,1:spec_bdy_width,nvchem), & |
---|
1076 | nested_grid%chem_bxe(jms:jme,kms:kme,1:spec_bdy_width,nvchem), & |
---|
1077 | nested_grid%chem_bys(ims:ime,kms:kme,1:spec_bdy_width,nvchem), & |
---|
1078 | nested_grid%chem_bye(ims:ime,kms:kme,1:spec_bdy_width,nvchem), & |
---|
1079 | 'T' , spec_bdy_width , & |
---|
1080 | ids , ide , jds , jde , kds , kde , & |
---|
1081 | ims , ime , jms , jme , kms , kme , & |
---|
1082 | ips , ipe , jps , jpe , kps , kpe ) |
---|
1083 | ! cbdy3dtemp0(ips:ipe,kps:kpe,jps:jpe,nvchem)=cbdy3dtemp1(ips:ipe,kps:kpe,jps:jpe) |
---|
1084 | ! if(nvchem.eq.p_o3)then |
---|
1085 | ! write(0,*)'fill 2ch_b3',cbdy3dtemp1(ips,kps,jps),cbdy3dtemp0(ips,kps,jps,nvchem) |
---|
1086 | ! endif |
---|
1087 | enddo |
---|
1088 | #endif |
---|
1089 | CALL stuff_bdy ( mbdy2dtemp1 , & |
---|
1090 | nested_grid%mu_bxs, nested_grid%mu_bxe , & |
---|
1091 | nested_grid%mu_bys, nested_grid%mu_bye , & |
---|
1092 | 'M' , spec_bdy_width , & |
---|
1093 | ids , ide , jds , jde , 1 , 1 , & |
---|
1094 | ims , ime , jms , jme , 1 , 1 , & |
---|
1095 | ips , ipe , jps , jpe , 1 , 1 ) |
---|
1096 | |
---|
1097 | ELSE IF ( time_loop .EQ. time_loop_max ) THEN |
---|
1098 | |
---|
1099 | ! u, theta, h, scalars coupled with my, v coupled with mx |
---|
1100 | CALL couple ( nested_grid%mu_2 , nested_grid%mub , ubdy3dtemp2 , nested_grid%u_2 , & |
---|
1101 | 'u' , nested_grid%msfuy , & |
---|
1102 | ids, ide, jds, jde, kds, kde, ims, ime, jms, jme, kms, kme, ips, ipe, jps, jpe, kps, kpe ) |
---|
1103 | CALL couple ( nested_grid%mu_2 , nested_grid%mub , vbdy3dtemp2 , nested_grid%v_2 , & |
---|
1104 | 'v' , nested_grid%msfvx , & |
---|
1105 | ids, ide, jds, jde, kds, kde, ims, ime, jms, jme, kms, kme, ips, ipe, jps, jpe, kps, kpe ) |
---|
1106 | CALL couple ( nested_grid%mu_2 , nested_grid%mub , tbdy3dtemp2 , nested_grid%t_2 , & |
---|
1107 | 't' , nested_grid%msfty , & |
---|
1108 | ids, ide, jds, jde, kds, kde, ims, ime, jms, jme, kms, kme, ips, ipe, jps, jpe, kps, kpe ) |
---|
1109 | CALL couple ( nested_grid%mu_2 , nested_grid%mub , pbdy3dtemp2 , nested_grid%ph_2 , & |
---|
1110 | 'h' , nested_grid%msfty , & |
---|
1111 | ids, ide, jds, jde, kds, kde, ims, ime, jms, jme, kms, kme, ips, ipe, jps, jpe, kps, kpe ) |
---|
1112 | CALL couple ( nested_grid%mu_2 , nested_grid%mub , qbdy3dtemp2 , nested_grid%moist(ims:ime,kms:kme,jms:jme,P_QV) , & |
---|
1113 | 't' , nested_grid%msfty , & |
---|
1114 | ids, ide, jds, jde, kds, kde, ims, ime, jms, jme, kms, kme, ips, ipe, jps, jpe, kps, kpe ) |
---|
1115 | mbdy2dtemp2(:,1,:) = nested_grid%mu_2(:,:) |
---|
1116 | |
---|
1117 | ! During all of the loops after the first loop, we first compute the boundary |
---|
1118 | ! tendencies with the current data values and the previously save information |
---|
1119 | ! stored in the *bdy3dtemp1 arrays. |
---|
1120 | #ifdef WRF_CHEM |
---|
1121 | do nvchem=1,num_chem |
---|
1122 | cbdy3dtemp1(ips:ipe,kps:kpe,jps:jpe)=cbdy3dtemp0(ips:ipe,kps:kpe,jps:jpe,nvchem) |
---|
1123 | cbdy3dtemp2(ips:ipe,kps:kpe,jps:jpe)=nested_grid%chem(ips:ipe,kps:kpe,jps:jpe,nvchem) |
---|
1124 | ! if(nvchem.eq.p_o3)then |
---|
1125 | ! write(0,*)'fill 1ch_b4',cbdy3dtemp1(ips,kps,jps),cbdy3dtemp0(ips,kps,jps,nvchem),cbdy3dtemp2(ips,kps,jps) |
---|
1126 | ! endif |
---|
1127 | CALL stuff_bdytend ( cbdy3dtemp2 , cbdy3dtemp1 , new_bdy_frq , & |
---|
1128 | nested_grid%chem_btxs(jms:jme,kms:kme,1:spec_bdy_width,nvchem), & |
---|
1129 | nested_grid%chem_btxe(jms:jme,kms:kme,1:spec_bdy_width,nvchem), & |
---|
1130 | nested_grid%chem_btys(ims:ime,kms:kme,1:spec_bdy_width,nvchem), & |
---|
1131 | nested_grid%chem_btye(ims:ime,kms:kme,1:spec_bdy_width,nvchem), & |
---|
1132 | 'T' , & |
---|
1133 | spec_bdy_width , & |
---|
1134 | ids , ide , jds , jde , kds , kde , & |
---|
1135 | ims , ime , jms , jme , kms , kme , & |
---|
1136 | ips , ipe , jps , jpe , kps , kpe ) |
---|
1137 | cbdy3dtemp0(ips:ipe,kps:kpe,jps:jpe,nvchem)=cbdy3dtemp2(ips:ipe,kps:kpe,jps:jpe) |
---|
1138 | ! if(nvchem.eq.p_o3)then |
---|
1139 | ! write(0,*)'fill 2ch_b4',cbdy3dtemp1(ips,kps,jps),cbdy3dtemp0(ips,kps,jps,nvchem),cbdy3dtemp2(ips,kps,jps) |
---|
1140 | ! endif |
---|
1141 | enddo |
---|
1142 | #endif |
---|
1143 | |
---|
1144 | CALL stuff_bdytend ( ubdy3dtemp2 , ubdy3dtemp1 , new_bdy_frq , & |
---|
1145 | nested_grid%u_btxs , nested_grid%u_btxe , & |
---|
1146 | nested_grid%u_btys , nested_grid%u_btye , & |
---|
1147 | 'U' , & |
---|
1148 | spec_bdy_width , & |
---|
1149 | ids , ide , jds , jde , kds , kde , & |
---|
1150 | ims , ime , jms , jme , kms , kme , & |
---|
1151 | ips , ipe , jps , jpe , kps , kpe ) |
---|
1152 | CALL stuff_bdytend ( vbdy3dtemp2 , vbdy3dtemp1 , new_bdy_frq , & |
---|
1153 | nested_grid%v_btxs , nested_grid%v_btxe , & |
---|
1154 | nested_grid%v_btys , nested_grid%v_btye , & |
---|
1155 | 'V' , & |
---|
1156 | spec_bdy_width , & |
---|
1157 | ids , ide , jds , jde , kds , kde , & |
---|
1158 | ims , ime , jms , jme , kms , kme , & |
---|
1159 | ips , ipe , jps , jpe , kps , kpe ) |
---|
1160 | CALL stuff_bdytend ( tbdy3dtemp2 , tbdy3dtemp1 , new_bdy_frq , & |
---|
1161 | nested_grid%t_btxs , nested_grid%t_btxe , & |
---|
1162 | nested_grid%t_btys , nested_grid%t_btye , & |
---|
1163 | 'T' , & |
---|
1164 | spec_bdy_width , & |
---|
1165 | ids , ide , jds , jde , kds , kde , & |
---|
1166 | ims , ime , jms , jme , kms , kme , & |
---|
1167 | ips , ipe , jps , jpe , kps , kpe ) |
---|
1168 | CALL stuff_bdytend ( pbdy3dtemp2 , pbdy3dtemp1 , new_bdy_frq , & |
---|
1169 | nested_grid%ph_btxs , nested_grid%ph_btxe , & |
---|
1170 | nested_grid%ph_btys , nested_grid%ph_btye , & |
---|
1171 | 'W' , & |
---|
1172 | spec_bdy_width , & |
---|
1173 | ids , ide , jds , jde , kds , kde , & |
---|
1174 | ims , ime , jms , jme , kms , kme , & |
---|
1175 | ips , ipe , jps , jpe , kps , kpe ) |
---|
1176 | CALL stuff_bdytend ( qbdy3dtemp2 , qbdy3dtemp1 , new_bdy_frq , & |
---|
1177 | nested_grid%moist_btxs(jms:jme,kms:kme,1:spec_bdy_width,P_QV) , & |
---|
1178 | nested_grid%moist_btxe(jms:jme,kms:kme,1:spec_bdy_width,P_QV) , & |
---|
1179 | nested_grid%moist_btys(ims:ime,kms:kme,1:spec_bdy_width,P_QV) , & |
---|
1180 | nested_grid%moist_btye(ims:ime,kms:kme,1:spec_bdy_width,P_QV) , & |
---|
1181 | 'T' , & |
---|
1182 | spec_bdy_width , & |
---|
1183 | ids , ide , jds , jde , kds , kde , & |
---|
1184 | ims , ime , jms , jme , kms , kme , & |
---|
1185 | ips , ipe , jps , jpe , kps , kpe ) |
---|
1186 | CALL stuff_bdytend ( mbdy2dtemp2 , mbdy2dtemp1 , new_bdy_frq , & |
---|
1187 | nested_grid%mu_btxs , nested_grid%mu_btxe , & |
---|
1188 | nested_grid%mu_btys , nested_grid%mu_btye , & |
---|
1189 | 'M' , & |
---|
1190 | spec_bdy_width , & |
---|
1191 | ids , ide , jds , jde , 1 , 1 , & |
---|
1192 | ims , ime , jms , jme , 1 , 1 , & |
---|
1193 | ips , ipe , jps , jpe , 1 , 1 ) |
---|
1194 | |
---|
1195 | IF ( time_loop .EQ. 2 ) THEN |
---|
1196 | |
---|
1197 | ! Generate an output file from this program, which will be an input file to WRF. |
---|
1198 | |
---|
1199 | CALL wrf_debug ( 100 , 'ndown_em main: calling open_w_dataset for wrfbdy' ) |
---|
1200 | CALL construct_filename1( bdyname , 'wrfbdy' , nested_grid%id , 2 ) |
---|
1201 | CALL open_w_dataset ( fidb, TRIM(bdyname) , nested_grid , config_flags , output_boundary , & |
---|
1202 | "DATASET=BOUNDARY", ierr ) |
---|
1203 | IF ( ierr .NE. 0 ) THEN |
---|
1204 | WRITE( wrf_err_message , FMT='(A,A,A,I8)' ) 'program ndown: error opening ',TRIM(bdyname),' for reading ierr=',ierr |
---|
1205 | CALL WRF_ERROR_FATAL ( wrf_err_message ) |
---|
1206 | ENDIF |
---|
1207 | |
---|
1208 | END IF |
---|
1209 | |
---|
1210 | ! Both pieces of the boundary data are now available to be written. |
---|
1211 | |
---|
1212 | CALL domain_clock_set( nested_grid, & |
---|
1213 | current_timestr=current_date(1:19) ) |
---|
1214 | temp24= current_date |
---|
1215 | temp24b=start_date_hold |
---|
1216 | start_date = start_date_hold |
---|
1217 | CALL geth_newdate ( temp19 , temp24b(1:19) , (time_loop-2) * model_config_rec%interval_seconds ) |
---|
1218 | current_date = temp19 // '.0000' |
---|
1219 | CALL geth_julgmt ( julyr , julday , gmt) |
---|
1220 | CALL nl_set_julyr ( nested_grid%id , julyr ) |
---|
1221 | CALL nl_set_julday ( nested_grid%id , julday ) |
---|
1222 | CALL nl_set_gmt ( nested_grid%id , gmt ) |
---|
1223 | CALL wrf_put_dom_ti_real ( fidb , 'GMT' , gmt , 1 , ierr ) |
---|
1224 | CALL wrf_put_dom_ti_integer ( fidb , 'JULYR' , julyr , 1 , ierr ) |
---|
1225 | CALL wrf_put_dom_ti_integer ( fidb , 'JULDAY' , julday , 1 , ierr ) |
---|
1226 | CALL domain_clock_set( nested_grid, & |
---|
1227 | current_timestr=current_date(1:19) ) |
---|
1228 | CALL output_boundary ( fidb , nested_grid , config_flags , ierr ) |
---|
1229 | current_date = temp24 |
---|
1230 | start_date = temp24b |
---|
1231 | CALL domain_clock_set( nested_grid, & |
---|
1232 | current_timestr=current_date(1:19) ) |
---|
1233 | |
---|
1234 | IF ( time_loop .EQ. 2 ) THEN |
---|
1235 | CALL wrf_put_dom_ti_real ( fidb , 'BDYFRQ' , new_bdy_frq , 1 , ierr ) |
---|
1236 | END IF |
---|
1237 | |
---|
1238 | ! Since this is the last time through here, we need to close the boundary file. |
---|
1239 | |
---|
1240 | CALL model_to_grid_config_rec ( nested_grid%id , model_config_rec , config_flags ) |
---|
1241 | CALL close_dataset ( fidb , config_flags , "DATASET=BOUNDARY" ) |
---|
1242 | |
---|
1243 | |
---|
1244 | END IF |
---|
1245 | |
---|
1246 | ! Process which time now? |
---|
1247 | |
---|
1248 | END DO big_time_loop_thingy |
---|
1249 | DEALLOCATE(znw_c) |
---|
1250 | DEALLOCATE(znu_c) |
---|
1251 | CALL model_to_grid_config_rec ( parent_grid%id , model_config_rec , config_flags ) |
---|
1252 | CALL med_shutdown_io ( parent_grid , config_flags ) |
---|
1253 | |
---|
1254 | CALL wrf_debug ( 0 , 'ndown_em: SUCCESS COMPLETE NDOWN_EM INIT' ) |
---|
1255 | |
---|
1256 | CALL wrf_shutdown |
---|
1257 | |
---|
1258 | CALL WRFU_Finalize( rc=rc ) |
---|
1259 | |
---|
1260 | END PROGRAM ndown_em |
---|
1261 | |
---|
1262 | SUBROUTINE land_percentages ( xland , & |
---|
1263 | landuse_frac , soil_top_cat , soil_bot_cat , & |
---|
1264 | isltyp , ivgtyp , & |
---|
1265 | num_veg_cat , num_soil_top_cat , num_soil_bot_cat , & |
---|
1266 | ids , ide , jds , jde , kds , kde , & |
---|
1267 | ims , ime , jms , jme , kms , kme , & |
---|
1268 | its , ite , jts , jte , kts , kte , & |
---|
1269 | iswater ) |
---|
1270 | USE module_soil_pre |
---|
1271 | |
---|
1272 | IMPLICIT NONE |
---|
1273 | |
---|
1274 | INTEGER , INTENT(IN) :: ids , ide , jds , jde , kds , kde , & |
---|
1275 | ims , ime , jms , jme , kms , kme , & |
---|
1276 | its , ite , jts , jte , kts , kte , & |
---|
1277 | iswater |
---|
1278 | |
---|
1279 | INTEGER , INTENT(IN) :: num_veg_cat , num_soil_top_cat , num_soil_bot_cat |
---|
1280 | REAL , DIMENSION(ims:ime,1:num_veg_cat,jms:jme) , INTENT(INOUT):: landuse_frac |
---|
1281 | REAL , DIMENSION(ims:ime,1:num_soil_top_cat,jms:jme) , INTENT(IN):: soil_top_cat |
---|
1282 | REAL , DIMENSION(ims:ime,1:num_soil_bot_cat,jms:jme) , INTENT(IN):: soil_bot_cat |
---|
1283 | INTEGER , DIMENSION(ims:ime,jms:jme), INTENT(OUT) :: isltyp , ivgtyp |
---|
1284 | REAL , DIMENSION(ims:ime,jms:jme) , INTENT(OUT) :: xland |
---|
1285 | |
---|
1286 | CALL process_percent_cat_new ( xland , & |
---|
1287 | landuse_frac , soil_top_cat , soil_bot_cat , & |
---|
1288 | isltyp , ivgtyp , & |
---|
1289 | num_veg_cat , num_soil_top_cat , num_soil_bot_cat , & |
---|
1290 | ids , ide , jds , jde , kds , kde , & |
---|
1291 | ims , ime , jms , jme , kms , kme , & |
---|
1292 | its , ite , jts , jte , kts , kte , & |
---|
1293 | iswater ) |
---|
1294 | |
---|
1295 | END SUBROUTINE land_percentages |
---|
1296 | |
---|
1297 | SUBROUTINE check_consistency ( ivgtyp , isltyp , landmask , & |
---|
1298 | ids , ide , jds , jde , kds , kde , & |
---|
1299 | ims , ime , jms , jme , kms , kme , & |
---|
1300 | its , ite , jts , jte , kts , kte , & |
---|
1301 | iswater ) |
---|
1302 | |
---|
1303 | IMPLICIT NONE |
---|
1304 | |
---|
1305 | INTEGER , INTENT(IN) :: ids , ide , jds , jde , kds , kde , & |
---|
1306 | ims , ime , jms , jme , kms , kme , & |
---|
1307 | its , ite , jts , jte , kts , kte , & |
---|
1308 | iswater |
---|
1309 | INTEGER , DIMENSION(ims:ime,jms:jme), INTENT(INOUT) :: isltyp , ivgtyp |
---|
1310 | REAL , DIMENSION(ims:ime,jms:jme), INTENT(INOUT) :: landmask |
---|
1311 | |
---|
1312 | LOGICAL :: oops |
---|
1313 | INTEGER :: oops_count , i , j |
---|
1314 | |
---|
1315 | oops = .FALSE. |
---|
1316 | oops_count = 0 |
---|
1317 | |
---|
1318 | DO j = jts, MIN(jde-1,jte) |
---|
1319 | DO i = its, MIN(ide-1,ite) |
---|
1320 | IF ( ( ( landmask(i,j) .LT. 0.5 ) .AND. ( ivgtyp(i,j) .NE. iswater ) ) .OR. & |
---|
1321 | ( ( landmask(i,j) .GT. 0.5 ) .AND. ( ivgtyp(i,j) .EQ. iswater ) ) ) THEN |
---|
1322 | print *,'mismatch in landmask and veg type' |
---|
1323 | print *,'i,j=',i,j, ' landmask =',NINT(landmask(i,j)),' ivgtyp=',ivgtyp(i,j) |
---|
1324 | oops = .TRUE. |
---|
1325 | oops_count = oops_count + 1 |
---|
1326 | landmask(i,j) = 0 |
---|
1327 | ivgtyp(i,j)=16 |
---|
1328 | isltyp(i,j)=14 |
---|
1329 | END IF |
---|
1330 | END DO |
---|
1331 | END DO |
---|
1332 | |
---|
1333 | IF ( oops ) THEN |
---|
1334 | CALL wrf_debug( 0, 'mismatch in check_consistency, turned to water points, be careful' ) |
---|
1335 | END IF |
---|
1336 | |
---|
1337 | END SUBROUTINE check_consistency |
---|
1338 | |
---|
1339 | SUBROUTINE check_consistency2( ivgtyp , isltyp , landmask , & |
---|
1340 | tmn , tsk , sst , xland , & |
---|
1341 | tslb , smois , sh2o , & |
---|
1342 | num_soil_layers , id , & |
---|
1343 | ids , ide , jds , jde , kds , kde , & |
---|
1344 | ims , ime , jms , jme , kms , kme , & |
---|
1345 | its , ite , jts , jte , kts , kte , & |
---|
1346 | iswater ) |
---|
1347 | |
---|
1348 | USE module_configure |
---|
1349 | USE module_optional_input |
---|
1350 | |
---|
1351 | INTEGER , INTENT(IN) :: ids , ide , jds , jde , kds , kde , & |
---|
1352 | ims , ime , jms , jme , kms , kme , & |
---|
1353 | its , ite , jts , jte , kts , kte |
---|
1354 | INTEGER , INTENT(IN) :: num_soil_layers , id |
---|
1355 | |
---|
1356 | INTEGER , DIMENSION(ims:ime,jms:jme) :: ivgtyp , isltyp |
---|
1357 | REAL , DIMENSION(ims:ime,jms:jme) :: landmask , tmn , tsk , sst , xland |
---|
1358 | REAL , DIMENSION(ims:ime,num_soil_layers,jms:jme) :: tslb , smois , sh2o |
---|
1359 | |
---|
1360 | INTEGER :: oops1 , oops2 |
---|
1361 | INTEGER :: i , j , k |
---|
1362 | |
---|
1363 | fix_tsk_tmn : SELECT CASE ( model_config_rec%sf_surface_physics(id) ) |
---|
1364 | |
---|
1365 | CASE ( SLABSCHEME , LSMSCHEME , RUCLSMSCHEME ) |
---|
1366 | DO j = jts, MIN(jde-1,jte) |
---|
1367 | DO i = its, MIN(ide-1,ite) |
---|
1368 | IF ( ( landmask(i,j) .LT. 0.5 ) .AND. ( flag_sst .EQ. 1 ) ) THEN |
---|
1369 | tmn(i,j) = sst(i,j) |
---|
1370 | tsk(i,j) = sst(i,j) |
---|
1371 | ELSE IF ( landmask(i,j) .LT. 0.5 ) THEN |
---|
1372 | tmn(i,j) = tsk(i,j) |
---|
1373 | END IF |
---|
1374 | END DO |
---|
1375 | END DO |
---|
1376 | END SELECT fix_tsk_tmn |
---|
1377 | |
---|
1378 | ! Is the TSK reasonable? |
---|
1379 | |
---|
1380 | DO j = jts, MIN(jde-1,jte) |
---|
1381 | DO i = its, MIN(ide-1,ite) |
---|
1382 | IF ( tsk(i,j) .LT. 170 .or. tsk(i,j) .GT. 400. ) THEN |
---|
1383 | print *,'error in the TSK' |
---|
1384 | print *,'i,j=',i,j |
---|
1385 | print *,'landmask=',landmask(i,j) |
---|
1386 | print *,'tsk, sst, tmn=',tsk(i,j),sst(i,j),tmn(i,j) |
---|
1387 | if(tmn(i,j).gt.170. .and. tmn(i,j).lt.400.)then |
---|
1388 | tsk(i,j)=tmn(i,j) |
---|
1389 | else if(sst(i,j).gt.170. .and. sst(i,j).lt.400.)then |
---|
1390 | tsk(i,j)=sst(i,j) |
---|
1391 | else |
---|
1392 | CALL wrf_error_fatal ( 'TSK unreasonable' ) |
---|
1393 | end if |
---|
1394 | END IF |
---|
1395 | END DO |
---|
1396 | END DO |
---|
1397 | |
---|
1398 | ! Is the TMN reasonable? |
---|
1399 | |
---|
1400 | DO j = jts, MIN(jde-1,jte) |
---|
1401 | DO i = its, MIN(ide-1,ite) |
---|
1402 | IF ( ( ( tmn(i,j) .LT. 170. ) .OR. ( tmn(i,j) .GT. 400. ) ) .AND. ( landmask(i,j) .GT. 0.5 ) ) THEN |
---|
1403 | print *,'error in the TMN' |
---|
1404 | print *,'i,j=',i,j |
---|
1405 | print *,'landmask=',landmask(i,j) |
---|
1406 | print *,'tsk, sst, tmn=',tsk(i,j),sst(i,j),tmn(i,j) |
---|
1407 | if(tsk(i,j).gt.170. .and. tsk(i,j).lt.400.)then |
---|
1408 | tmn(i,j)=tsk(i,j) |
---|
1409 | else if(sst(i,j).gt.170. .and. sst(i,j).lt.400.)then |
---|
1410 | tmn(i,j)=sst(i,j) |
---|
1411 | else |
---|
1412 | CALL wrf_error_fatal ( 'TMN unreasonable' ) |
---|
1413 | endif |
---|
1414 | END IF |
---|
1415 | END DO |
---|
1416 | END DO |
---|
1417 | |
---|
1418 | ! Is the TSLB reasonable? |
---|
1419 | |
---|
1420 | DO j = jts, MIN(jde-1,jte) |
---|
1421 | DO i = its, MIN(ide-1,ite) |
---|
1422 | IF ( ( ( tslb(i,1,j) .LT. 170. ) .OR. ( tslb(i,1,j) .GT. 400. ) ) .AND. ( landmask(i,j) .GT. 0.5 ) ) THEN |
---|
1423 | print *,'error in the TSLB' |
---|
1424 | print *,'i,j=',i,j |
---|
1425 | print *,'landmask=',landmask(i,j) |
---|
1426 | print *,'tsk, sst, tmn=',tsk(i,j),sst(i,j),tmn(i,j) |
---|
1427 | print *,'tslb = ',tslb(i,:,j) |
---|
1428 | print *,'old smois = ',smois(i,:,j) |
---|
1429 | DO l = 1 , num_soil_layers |
---|
1430 | sh2o(i,l,j) = 0.0 |
---|
1431 | END DO |
---|
1432 | DO l = 1 , num_soil_layers |
---|
1433 | smois(i,l,j) = 0.3 |
---|
1434 | END DO |
---|
1435 | if(tsk(i,j).gt.170. .and. tsk(i,j).lt.400.)then |
---|
1436 | DO l = 1 , num_soil_layers |
---|
1437 | tslb(i,l,j)=tsk(i,j) |
---|
1438 | END DO |
---|
1439 | else if(sst(i,j).gt.170. .and. sst(i,j).lt.400.)then |
---|
1440 | DO l = 1 , num_soil_layers |
---|
1441 | tslb(i,l,j)=sst(i,j) |
---|
1442 | END DO |
---|
1443 | else if(tmn(i,j).gt.170. .and. tmn(i,j).lt.400.)then |
---|
1444 | DO l = 1 , num_soil_layers |
---|
1445 | tslb(i,l,j)=tmn(i,j) |
---|
1446 | END DO |
---|
1447 | else |
---|
1448 | CALL wrf_error_fatal ( 'TSLB unreasonable' ) |
---|
1449 | endif |
---|
1450 | END IF |
---|
1451 | END DO |
---|
1452 | END DO |
---|
1453 | |
---|
1454 | ! Let us make sure (again) that the landmask and the veg/soil categories match. |
---|
1455 | |
---|
1456 | oops1=0 |
---|
1457 | oops2=0 |
---|
1458 | DO j = jts, MIN(jde-1,jte) |
---|
1459 | DO i = its, MIN(ide-1,ite) |
---|
1460 | IF ( ( ( landmask(i,j) .LT. 0.5 ) .AND. ( ivgtyp(i,j) .NE. iswater .OR. isltyp(i,j) .NE. 14 ) ) .OR. & |
---|
1461 | ( ( landmask(i,j) .GT. 0.5 ) .AND. ( ivgtyp(i,j) .EQ. iswater .OR. isltyp(i,j) .EQ. 14 ) ) ) THEN |
---|
1462 | IF ( tslb(i,1,j) .GT. 1. ) THEN |
---|
1463 | oops1=oops1+1 |
---|
1464 | ivgtyp(i,j) = 5 |
---|
1465 | isltyp(i,j) = 8 |
---|
1466 | landmask(i,j) = 1 |
---|
1467 | xland(i,j) = 1 |
---|
1468 | ELSE IF ( sst(i,j) .GT. 1. ) THEN |
---|
1469 | oops2=oops2+1 |
---|
1470 | ivgtyp(i,j) = iswater |
---|
1471 | isltyp(i,j) = 14 |
---|
1472 | landmask(i,j) = 0 |
---|
1473 | xland(i,j) = 2 |
---|
1474 | ELSE |
---|
1475 | print *,'the landmask and soil/veg cats do not match' |
---|
1476 | print *,'i,j=',i,j |
---|
1477 | print *,'landmask=',landmask(i,j) |
---|
1478 | print *,'ivgtyp=',ivgtyp(i,j) |
---|
1479 | print *,'isltyp=',isltyp(i,j) |
---|
1480 | print *,'iswater=', iswater |
---|
1481 | print *,'tslb=',tslb(i,:,j) |
---|
1482 | print *,'sst=',sst(i,j) |
---|
1483 | CALL wrf_error_fatal ( 'mismatch_landmask_ivgtyp' ) |
---|
1484 | END IF |
---|
1485 | END IF |
---|
1486 | END DO |
---|
1487 | END DO |
---|
1488 | if (oops1.gt.0) then |
---|
1489 | print *,'points artificially set to land : ',oops1 |
---|
1490 | endif |
---|
1491 | if(oops2.gt.0) then |
---|
1492 | print *,'points artificially set to water: ',oops2 |
---|
1493 | endif |
---|
1494 | |
---|
1495 | END SUBROUTINE check_consistency2 |
---|
1496 | !!!!!!!!!!!!!!!!!!!!!!!!!!!11 |
---|
1497 | SUBROUTINE vert_cor(parent,znw_c,k_dim_c) |
---|
1498 | USE module_domain |
---|
1499 | IMPLICIT NONE |
---|
1500 | TYPE(domain), POINTER :: parent |
---|
1501 | integer , intent(in) :: k_dim_c |
---|
1502 | real , dimension(k_dim_c), INTENT(IN) :: znw_c |
---|
1503 | |
---|
1504 | integer :: kde_c , kde_n ,n_refine,ii,kkk,k |
---|
1505 | real :: dznw_m,cof1,cof2 |
---|
1506 | |
---|
1507 | kde_c = k_dim_c |
---|
1508 | kde_n = parent%e_vert |
---|
1509 | n_refine = parent % vert_refine_fact |
---|
1510 | |
---|
1511 | |
---|
1512 | kkk = 0 |
---|
1513 | do k = 1 , kde_c-1 |
---|
1514 | dznw_m = znw_c(k+1) - znw_c(k) |
---|
1515 | do ii = 1,n_refine |
---|
1516 | kkk = kkk + 1 |
---|
1517 | parent%znw(kkk) = znw_c(k) + float(ii-1)/float(n_refine)*dznw_m |
---|
1518 | enddo |
---|
1519 | enddo |
---|
1520 | parent%znw(kde_n) = znw_c(kde_c) |
---|
1521 | parent%znw(1) = znw_c(1) |
---|
1522 | |
---|
1523 | DO k=1, kde_n-1 |
---|
1524 | parent%dnw(k) = parent%znw(k+1) - parent%znw(k) |
---|
1525 | parent%rdnw(k) = 1./parent%dnw(k) |
---|
1526 | parent%znu(k) = 0.5*(parent%znw(k+1)+parent%znw(k)) |
---|
1527 | END DO |
---|
1528 | |
---|
1529 | DO k=2, kde_n-1 |
---|
1530 | parent%dn(k) = 0.5*(parent%dnw(k)+parent%dnw(k-1)) |
---|
1531 | parent%rdn(k) = 1./parent%dn(k) |
---|
1532 | parent%fnp(k) = .5* parent%dnw(k )/parent%dn(k) |
---|
1533 | parent%fnm(k) = .5* parent%dnw(k-1)/parent%dn(k) |
---|
1534 | END DO |
---|
1535 | |
---|
1536 | cof1 = (2.*parent%dn(2)+parent%dn(3))/(parent%dn(2)+parent%dn(3))*parent%dnw(1)/parent%dn(2) |
---|
1537 | cof2 = parent%dn(2) /(parent%dn(2)+parent%dn(3))*parent%dnw(1)/parent%dn(3) |
---|
1538 | |
---|
1539 | parent%cf1 = parent%fnp(2) + cof1 |
---|
1540 | parent%cf2 = parent%fnm(2) - cof1 - cof2 |
---|
1541 | parent%cf3 = cof2 |
---|
1542 | |
---|
1543 | parent%cfn = (.5*parent%dnw(kde_n-1)+parent%dn(kde_n-1))/parent%dn(kde_n-1) |
---|
1544 | parent%cfn1 = -.5*parent%dnw(kde_n-1)/parent%dn(kde_n-1) |
---|
1545 | |
---|
1546 | |
---|
1547 | |
---|
1548 | END SUBROUTINE vert_cor |
---|
1549 | |
---|
1550 | |
---|
1551 | SUBROUTINE init_domain_constants_em_ptr ( parent , nest ) |
---|
1552 | USE module_domain |
---|
1553 | USE module_configure |
---|
1554 | IMPLICIT NONE |
---|
1555 | TYPE(domain), POINTER :: parent , nest |
---|
1556 | INTERFACE |
---|
1557 | SUBROUTINE init_domain_constants_em ( parent , nest ) |
---|
1558 | USE module_domain |
---|
1559 | USE module_configure |
---|
1560 | TYPE(domain) :: parent , nest |
---|
1561 | END SUBROUTINE init_domain_constants_em |
---|
1562 | END INTERFACE |
---|
1563 | CALL init_domain_constants_em ( parent , nest ) |
---|
1564 | END SUBROUTINE init_domain_constants_em_ptr |
---|
1565 | |
---|
1566 | |
---|
1567 | SUBROUTINE vertical_interp (parent_grid,znw_c,znu_c,cf1_c,cf2_c,cf3_c,cfn_c,cfn1_c,k_dim_c) |
---|
1568 | USE module_domain |
---|
1569 | USE module_configure |
---|
1570 | IMPLICIT NONE |
---|
1571 | TYPE(domain), POINTER :: parent_grid |
---|
1572 | INTEGER , INTENT (IN) :: k_dim_c |
---|
1573 | REAL , INTENT (IN) :: cf1_c,cf2_c,cf3_c,cfn_c,cfn1_c |
---|
1574 | REAL , DIMENSION(k_dim_c) , INTENT (IN) :: znw_c,znu_c |
---|
1575 | |
---|
1576 | integer :: kde_c , kde_n ,n_refine,ii,kkk |
---|
1577 | integer :: i , j, k , itrace |
---|
1578 | real :: p_top_m , p_surf_m , mu_m , hsca_m , pre_c ,pre_n |
---|
1579 | |
---|
1580 | real, allocatable, dimension(:) :: alt_w_c , alt_u_c ,pro_w_c , pro_u_c |
---|
1581 | real, allocatable, dimension(:) :: alt_w_n , alt_u_n ,pro_w_n , pro_u_n |
---|
1582 | |
---|
1583 | INTEGER :: nids, nide, njds, njde, nkds, nkde, & |
---|
1584 | nims, nime, njms, njme, nkms, nkme, & |
---|
1585 | nips, nipe, njps, njpe, nkps, nkpe |
---|
1586 | |
---|
1587 | |
---|
1588 | hsca_m = 6.7 |
---|
1589 | n_refine = model_config_rec % vert_refine_fact |
---|
1590 | kde_c = k_dim_c |
---|
1591 | kde_n = parent_grid%e_vert |
---|
1592 | |
---|
1593 | CALL get_ijk_from_grid ( parent_grid , & |
---|
1594 | nids, nide, njds, njde, nkds, nkde, & |
---|
1595 | nims, nime, njms, njme, nkms, nkme, & |
---|
1596 | nips, nipe, njps, njpe, nkps, nkpe ) |
---|
1597 | |
---|
1598 | print * , 'MOUSTA_VER ',parent_grid%e_vert , kde_c , kde_n |
---|
1599 | print *, nids , nide , njds , njde , nkds , nkde |
---|
1600 | print *, nims , nime , njms , njme , nkms , nkme |
---|
1601 | print *, nips , nipe , njps , njpe , nkps , nkpe |
---|
1602 | |
---|
1603 | |
---|
1604 | |
---|
1605 | allocate (alt_w_c(kde_c)) |
---|
1606 | allocate (alt_u_c(kde_c+1)) |
---|
1607 | allocate (pro_w_c(kde_c)) |
---|
1608 | allocate (pro_u_c(kde_c+1)) |
---|
1609 | |
---|
1610 | allocate (alt_w_n(kde_n)) |
---|
1611 | allocate (alt_u_n(kde_n+1)) |
---|
1612 | allocate (pro_w_n(kde_n)) |
---|
1613 | allocate (pro_u_n(kde_n+1)) |
---|
1614 | |
---|
1615 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!11111 |
---|
1616 | p_top_m = parent_grid%p_top |
---|
1617 | p_surf_m = 1.e5 |
---|
1618 | mu_m = p_surf_m - p_top_m |
---|
1619 | print * , 'p_top_m', p_top_m |
---|
1620 | ! parent |
---|
1621 | do k = 1,kde_c |
---|
1622 | pre_c = mu_m * znw_c(k) + p_top_m |
---|
1623 | alt_w_c(k) = -hsca_m * alog(pre_c/p_surf_m) |
---|
1624 | enddo |
---|
1625 | do k = 1,kde_c-1 |
---|
1626 | pre_c = mu_m * znu_c(k) + p_top_m |
---|
1627 | alt_u_c(k+1) = -hsca_m * alog(pre_c/p_surf_m) |
---|
1628 | enddo |
---|
1629 | alt_u_c(1) = alt_w_c(1) |
---|
1630 | alt_u_c(kde_c+1) = alt_w_c(kde_c) |
---|
1631 | ! nest |
---|
1632 | do k = 1,kde_n |
---|
1633 | pre_n = mu_m * parent_grid%znw(k) + p_top_m |
---|
1634 | alt_w_n(k) = -hsca_m * alog(pre_n/p_surf_m) |
---|
1635 | enddo |
---|
1636 | do k = 1,kde_n-1 |
---|
1637 | pre_n = mu_m * parent_grid%znu(k) + p_top_m |
---|
1638 | alt_u_n(k+1) = -hsca_m * alog(pre_n/p_surf_m) |
---|
1639 | enddo |
---|
1640 | alt_u_n(1) = alt_w_n(1) |
---|
1641 | alt_u_n(kde_n+1) = alt_w_n(kde_n) |
---|
1642 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
1643 | |
---|
1644 | IF ( SIZE( parent_grid%u_2, 1 ) * SIZE( parent_grid%u_2, 3 ) .GT. 1 ) THEN |
---|
1645 | do j = njms , njme |
---|
1646 | do i = nims , nime |
---|
1647 | |
---|
1648 | do k = 1,kde_c-1 |
---|
1649 | pro_u_c(k+1) = parent_grid%u_2(i,k,j) |
---|
1650 | enddo |
---|
1651 | pro_u_c(1 ) = cf1_c*parent_grid%u_2(i,1,j) & |
---|
1652 | + cf2_c*parent_grid%u_2(i,2,j) & |
---|
1653 | + cf3_c*parent_grid%u_2(i,3,j) |
---|
1654 | |
---|
1655 | pro_u_c(kde_c+1) = cfn_c *parent_grid%u_2(i,kde_c-1,j) & |
---|
1656 | + cfn1_c*parent_grid%u_2(i,kde_c-2,j) |
---|
1657 | |
---|
1658 | call inter(pro_u_c,alt_u_c,kde_c+1,pro_u_n,alt_u_n,kde_n+1) |
---|
1659 | |
---|
1660 | do k = 1,kde_n-1 |
---|
1661 | parent_grid%u_2(i,k,j) = pro_u_n(k+1) |
---|
1662 | enddo |
---|
1663 | |
---|
1664 | enddo |
---|
1665 | enddo |
---|
1666 | ENDIF |
---|
1667 | |
---|
1668 | IF ( SIZE( parent_grid%v_2, 1 ) * SIZE( parent_grid%v_2, 3 ) .GT. 1 ) THEN |
---|
1669 | do j = njms , njme |
---|
1670 | do i = nims , nime |
---|
1671 | |
---|
1672 | do k = 1,kde_c-1 |
---|
1673 | pro_u_c(k+1) = parent_grid%v_2(i,k,j) |
---|
1674 | enddo |
---|
1675 | pro_u_c(1 ) = cf1_c*parent_grid%v_2(i,1,j) & |
---|
1676 | + cf2_c*parent_grid%v_2(i,2,j) & |
---|
1677 | + cf3_c*parent_grid%v_2(i,3,j) |
---|
1678 | |
---|
1679 | pro_u_c(kde_c+1) = cfn_c *parent_grid%v_2(i,kde_c-1,j) & |
---|
1680 | + cfn1_c*parent_grid%v_2(i,kde_c-2,j) |
---|
1681 | |
---|
1682 | call inter(pro_u_c,alt_u_c,kde_c+1,pro_u_n,alt_u_n,kde_n+1) |
---|
1683 | |
---|
1684 | do k = 1,kde_n-1 |
---|
1685 | parent_grid%v_2(i,k,j) = pro_u_n(k+1) |
---|
1686 | enddo |
---|
1687 | |
---|
1688 | enddo |
---|
1689 | enddo |
---|
1690 | ENDIF |
---|
1691 | |
---|
1692 | IF ( SIZE( parent_grid%w_2, 1 ) * SIZE( parent_grid%w_2, 3 ) .GT. 1 ) THEN |
---|
1693 | do j = njms , njme |
---|
1694 | do i = nims , nime |
---|
1695 | |
---|
1696 | do k = 1,kde_c |
---|
1697 | pro_w_c(k) = parent_grid%w_2(i,k,j) |
---|
1698 | enddo |
---|
1699 | call inter(pro_w_c,alt_w_c,kde_c,pro_w_n,alt_w_n,kde_n) |
---|
1700 | |
---|
1701 | do k = 1,kde_n |
---|
1702 | parent_grid%w_2(i,k,j) = pro_w_n(k) |
---|
1703 | enddo |
---|
1704 | enddo |
---|
1705 | enddo |
---|
1706 | ENDIF |
---|
1707 | |
---|
1708 | IF ( SIZE( parent_grid%t_2, 1 ) * SIZE( parent_grid%t_2, 3 ) .GT. 1 ) THEN |
---|
1709 | do j = njms , njme |
---|
1710 | do i = nims , nime |
---|
1711 | |
---|
1712 | do k = 1,kde_c-1 |
---|
1713 | pro_u_c(k+1) = parent_grid%t_2(i,k,j) |
---|
1714 | enddo |
---|
1715 | pro_u_c(1 ) = cf1_c*parent_grid%t_2(i,1,j) & |
---|
1716 | + cf2_c*parent_grid%t_2(i,2,j) & |
---|
1717 | + cf3_c*parent_grid%t_2(i,3,j) |
---|
1718 | |
---|
1719 | pro_u_c(kde_c+1) = cfn_c *parent_grid%t_2(i,kde_c-1,j) & |
---|
1720 | + cfn1_c*parent_grid%t_2(i,kde_c-2,j) |
---|
1721 | |
---|
1722 | call inter(pro_u_c,alt_u_c,kde_c+1,pro_u_n,alt_u_n,kde_n+1) |
---|
1723 | |
---|
1724 | do k = 1,kde_n-1 |
---|
1725 | parent_grid%t_2(i,k,j) = pro_u_n(k+1) |
---|
1726 | enddo |
---|
1727 | |
---|
1728 | enddo |
---|
1729 | enddo |
---|
1730 | ENDIF |
---|
1731 | |
---|
1732 | DO itrace = PARAM_FIRST_SCALAR, num_moist |
---|
1733 | IF ( SIZE( parent_grid%moist, 1 ) * SIZE( parent_grid%moist, 3 ) .GT. 1 ) THEN |
---|
1734 | do j = njms , njme |
---|
1735 | do i = nims , nime |
---|
1736 | |
---|
1737 | do k = 1,kde_c-1 |
---|
1738 | pro_u_c(k+1) = parent_grid%moist(i,k,j,itrace) |
---|
1739 | enddo |
---|
1740 | pro_u_c(1 ) = cf1_c*parent_grid%moist(i,1,j,itrace) & |
---|
1741 | + cf2_c*parent_grid%moist(i,2,j,itrace) & |
---|
1742 | + cf3_c*parent_grid%moist(i,3,j,itrace) |
---|
1743 | |
---|
1744 | pro_u_c(kde_c+1) = cfn_c *parent_grid%moist(i,kde_c-1,j,itrace) & |
---|
1745 | + cfn1_c*parent_grid%moist(i,kde_c-2,j,itrace) |
---|
1746 | |
---|
1747 | call inter(pro_u_c,alt_u_c,kde_c+1,pro_u_n,alt_u_n,kde_n+1) |
---|
1748 | |
---|
1749 | do k = 1,kde_n-1 |
---|
1750 | parent_grid%moist(i,k,j,itrace) = pro_u_n(k+1) |
---|
1751 | enddo |
---|
1752 | |
---|
1753 | enddo |
---|
1754 | enddo |
---|
1755 | ENDIF |
---|
1756 | ENDDO |
---|
1757 | |
---|
1758 | DO itrace = PARAM_FIRST_SCALAR, num_dfi_moist |
---|
1759 | IF ( SIZE( parent_grid%dfi_moist, 1 ) * SIZE( parent_grid%dfi_moist, 3 ) .GT. 1 ) THEN |
---|
1760 | do j = njms , njme |
---|
1761 | do i = nims , nime |
---|
1762 | |
---|
1763 | do k = 1,kde_c-1 |
---|
1764 | pro_u_c(k+1) = parent_grid%dfi_moist(i,k,j,itrace) |
---|
1765 | enddo |
---|
1766 | pro_u_c(1 ) = cf1_c*parent_grid%dfi_moist(i,1,j,itrace) & |
---|
1767 | + cf2_c*parent_grid%dfi_moist(i,2,j,itrace) & |
---|
1768 | + cf3_c*parent_grid%dfi_moist(i,3,j,itrace) |
---|
1769 | |
---|
1770 | pro_u_c(kde_c+1) = cfn_c *parent_grid%dfi_moist(i,kde_c-1,j,itrace) & |
---|
1771 | + cfn1_c*parent_grid%dfi_moist(i,kde_c-2,j,itrace) |
---|
1772 | |
---|
1773 | call inter(pro_u_c,alt_u_c,kde_c+1,pro_u_n,alt_u_n,kde_n+1) |
---|
1774 | |
---|
1775 | do k = 1,kde_n-1 |
---|
1776 | parent_grid%dfi_moist(i,k,j,itrace) = pro_u_n(k+1) |
---|
1777 | enddo |
---|
1778 | |
---|
1779 | enddo |
---|
1780 | enddo |
---|
1781 | ENDIF |
---|
1782 | ENDDO |
---|
1783 | |
---|
1784 | DO itrace = PARAM_FIRST_SCALAR, num_scalar |
---|
1785 | IF ( SIZE( parent_grid%scalar, 1 ) * SIZE( parent_grid%scalar, 3 ) .GT. 1 ) THEN |
---|
1786 | do j = njms , njme |
---|
1787 | do i = nims , nime |
---|
1788 | |
---|
1789 | do k = 1,kde_c-1 |
---|
1790 | pro_u_c(k+1) = parent_grid%scalar(i,k,j,itrace) |
---|
1791 | enddo |
---|
1792 | pro_u_c(1 ) = cf1_c*parent_grid%scalar(i,1,j,itrace) & |
---|
1793 | + cf2_c*parent_grid%scalar(i,2,j,itrace) & |
---|
1794 | + cf3_c*parent_grid%scalar(i,3,j,itrace) |
---|
1795 | |
---|
1796 | pro_u_c(kde_c+1) = cfn_c *parent_grid%scalar(i,kde_c-1,j,itrace) & |
---|
1797 | + cfn1_c*parent_grid%scalar(i,kde_c-2,j,itrace) |
---|
1798 | |
---|
1799 | call inter(pro_u_c,alt_u_c,kde_c+1,pro_u_n,alt_u_n,kde_n+1) |
---|
1800 | |
---|
1801 | do k = 1,kde_n-1 |
---|
1802 | parent_grid%scalar(i,k,j,itrace) = pro_u_n(k+1) |
---|
1803 | enddo |
---|
1804 | |
---|
1805 | enddo |
---|
1806 | enddo |
---|
1807 | ENDIF |
---|
1808 | ENDDO |
---|
1809 | |
---|
1810 | DO itrace = PARAM_FIRST_SCALAR, num_dfi_scalar |
---|
1811 | IF ( SIZE( parent_grid%dfi_scalar, 1 ) * SIZE( parent_grid%dfi_scalar, 3 ) .GT. 1 ) THEN |
---|
1812 | do j = njms , njme |
---|
1813 | do i = nims , nime |
---|
1814 | |
---|
1815 | do k = 1,kde_c-1 |
---|
1816 | pro_u_c(k+1) = parent_grid%dfi_scalar(i,k,j,itrace) |
---|
1817 | enddo |
---|
1818 | pro_u_c(1 ) = cf1_c*parent_grid%dfi_scalar(i,1,j,itrace) & |
---|
1819 | + cf2_c*parent_grid%dfi_scalar(i,2,j,itrace) & |
---|
1820 | + cf3_c*parent_grid%dfi_scalar(i,3,j,itrace) |
---|
1821 | |
---|
1822 | pro_u_c(kde_c+1) = cfn_c *parent_grid%dfi_scalar(i,kde_c-1,j,itrace) & |
---|
1823 | + cfn1_c*parent_grid%dfi_scalar(i,kde_c-2,j,itrace) |
---|
1824 | |
---|
1825 | call inter(pro_u_c,alt_u_c,kde_c+1,pro_u_n,alt_u_n,kde_n+1) |
---|
1826 | |
---|
1827 | do k = 1,kde_n-1 |
---|
1828 | parent_grid%dfi_scalar(i,k,j,itrace) = pro_u_n(k+1) |
---|
1829 | enddo |
---|
1830 | |
---|
1831 | enddo |
---|
1832 | enddo |
---|
1833 | ENDIF |
---|
1834 | ENDDO |
---|
1835 | |
---|
1836 | |
---|
1837 | |
---|
1838 | IF ( SIZE( parent_grid%f_ice_phy, 1 ) * SIZE( parent_grid%f_ice_phy, 3 ) .GT. 1 ) THEN |
---|
1839 | do j = njms , njme |
---|
1840 | do i = nims , nime |
---|
1841 | |
---|
1842 | do k = 1,kde_c-1 |
---|
1843 | pro_u_c(k+1) = parent_grid%f_ice_phy(i,k,j) |
---|
1844 | enddo |
---|
1845 | pro_u_c(1 ) = cf1_c*parent_grid%f_ice_phy(i,1,j) & |
---|
1846 | + cf2_c*parent_grid%f_ice_phy(i,2,j) & |
---|
1847 | + cf3_c*parent_grid%f_ice_phy(i,3,j) |
---|
1848 | |
---|
1849 | pro_u_c(kde_c+1) = cfn_c *parent_grid%f_ice_phy(i,kde_c-1,j) & |
---|
1850 | + cfn1_c*parent_grid%f_ice_phy(i,kde_c-2,j) |
---|
1851 | |
---|
1852 | call inter(pro_u_c,alt_u_c,kde_c+1,pro_u_n,alt_u_n,kde_n+1) |
---|
1853 | |
---|
1854 | do k = 1,kde_n-1 |
---|
1855 | parent_grid%f_ice_phy(i,k,j) = pro_u_n(k+1) |
---|
1856 | enddo |
---|
1857 | |
---|
1858 | enddo |
---|
1859 | enddo |
---|
1860 | ENDIF |
---|
1861 | |
---|
1862 | IF ( SIZE( parent_grid%f_rain_phy, 1 ) * SIZE( parent_grid%f_rain_phy, 3 ) .GT. 1 ) THEN |
---|
1863 | do j = njms , njme |
---|
1864 | do i = nims , nime |
---|
1865 | |
---|
1866 | do k = 1,kde_c-1 |
---|
1867 | pro_u_c(k+1) = parent_grid%f_rain_phy(i,k,j) |
---|
1868 | enddo |
---|
1869 | pro_u_c(1 ) = cf1_c*parent_grid%f_rain_phy(i,1,j) & |
---|
1870 | + cf2_c*parent_grid%f_rain_phy(i,2,j) & |
---|
1871 | + cf3_c*parent_grid%f_rain_phy(i,3,j) |
---|
1872 | |
---|
1873 | pro_u_c(kde_c+1) = cfn_c *parent_grid%f_rain_phy(i,kde_c-1,j) & |
---|
1874 | + cfn1_c*parent_grid%f_rain_phy(i,kde_c-2,j) |
---|
1875 | |
---|
1876 | call inter(pro_u_c,alt_u_c,kde_c+1,pro_u_n,alt_u_n,kde_n+1) |
---|
1877 | |
---|
1878 | do k = 1,kde_n-1 |
---|
1879 | parent_grid%f_rain_phy(i,k,j) = pro_u_n(k+1) |
---|
1880 | enddo |
---|
1881 | |
---|
1882 | enddo |
---|
1883 | enddo |
---|
1884 | ENDIF |
---|
1885 | |
---|
1886 | |
---|
1887 | IF ( SIZE( parent_grid%f_rimef_phy, 1 ) * SIZE( parent_grid%f_rimef_phy, 3 ) .GT. 1 ) THEN |
---|
1888 | do j = njms , njme |
---|
1889 | do i = nims , nime |
---|
1890 | |
---|
1891 | do k = 1,kde_c-1 |
---|
1892 | pro_u_c(k+1) = parent_grid%f_rimef_phy(i,k,j) |
---|
1893 | enddo |
---|
1894 | pro_u_c(1 ) = cf1_c*parent_grid%f_rimef_phy(i,1,j) & |
---|
1895 | + cf2_c*parent_grid%f_rimef_phy(i,2,j) & |
---|
1896 | + cf3_c*parent_grid%f_rimef_phy(i,3,j) |
---|
1897 | |
---|
1898 | pro_u_c(kde_c+1) = cfn_c *parent_grid%f_rimef_phy(i,kde_c-1,j) & |
---|
1899 | + cfn1_c*parent_grid%f_rimef_phy(i,kde_c-2,j) |
---|
1900 | |
---|
1901 | call inter(pro_u_c,alt_u_c,kde_c+1,pro_u_n,alt_u_n,kde_n+1) |
---|
1902 | |
---|
1903 | do k = 1,kde_n-1 |
---|
1904 | parent_grid%f_rimef_phy(i,k,j) = pro_u_n(k+1) |
---|
1905 | enddo |
---|
1906 | |
---|
1907 | enddo |
---|
1908 | enddo |
---|
1909 | ENDIF |
---|
1910 | |
---|
1911 | IF ( SIZE( parent_grid%h_diabatic, 1 ) * SIZE( parent_grid%h_diabatic, 3 ) .GT. 1 ) THEN |
---|
1912 | do j = njms , njme |
---|
1913 | do i = nims , nime |
---|
1914 | |
---|
1915 | do k = 1,kde_c-1 |
---|
1916 | pro_u_c(k+1) = parent_grid%h_diabatic(i,k,j) |
---|
1917 | enddo |
---|
1918 | pro_u_c(1 ) = cf1_c*parent_grid%h_diabatic(i,1,j) & |
---|
1919 | + cf2_c*parent_grid%h_diabatic(i,2,j) & |
---|
1920 | + cf3_c*parent_grid%h_diabatic(i,3,j) |
---|
1921 | |
---|
1922 | pro_u_c(kde_c+1) = cfn_c *parent_grid%h_diabatic(i,kde_c-1,j) & |
---|
1923 | + cfn1_c*parent_grid%h_diabatic(i,kde_c-2,j) |
---|
1924 | |
---|
1925 | call inter(pro_u_c,alt_u_c,kde_c+1,pro_u_n,alt_u_n,kde_n+1) |
---|
1926 | |
---|
1927 | do k = 1,kde_n-1 |
---|
1928 | parent_grid%h_diabatic(i,k,j) = pro_u_n(k+1) |
---|
1929 | enddo |
---|
1930 | |
---|
1931 | enddo |
---|
1932 | enddo |
---|
1933 | ENDIF |
---|
1934 | |
---|
1935 | IF ( SIZE( parent_grid%rthraten, 1 ) * SIZE( parent_grid%rthraten, 3 ) .GT. 1 ) THEN |
---|
1936 | do j = njms , njme |
---|
1937 | do i = nims , nime |
---|
1938 | |
---|
1939 | do k = 1,kde_c-1 |
---|
1940 | pro_u_c(k+1) = parent_grid%rthraten(i,k,j) |
---|
1941 | enddo |
---|
1942 | pro_u_c(1 ) = cf1_c*parent_grid%rthraten(i,1,j) & |
---|
1943 | + cf2_c*parent_grid%rthraten(i,2,j) & |
---|
1944 | + cf3_c*parent_grid%rthraten(i,3,j) |
---|
1945 | |
---|
1946 | pro_u_c(kde_c+1) = cfn_c *parent_grid%rthraten(i,kde_c-1,j) & |
---|
1947 | + cfn1_c*parent_grid%rthraten(i,kde_c-2,j) |
---|
1948 | |
---|
1949 | call inter(pro_u_c,alt_u_c,kde_c+1,pro_u_n,alt_u_n,kde_n+1) |
---|
1950 | |
---|
1951 | do k = 1,kde_n-1 |
---|
1952 | parent_grid%rthraten(i,k,j) = pro_u_n(k+1) |
---|
1953 | enddo |
---|
1954 | |
---|
1955 | enddo |
---|
1956 | enddo |
---|
1957 | ENDIF |
---|
1958 | |
---|
1959 | |
---|
1960 | |
---|
1961 | |
---|
1962 | |
---|
1963 | |
---|
1964 | |
---|
1965 | |
---|
1966 | |
---|
1967 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
1968 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
1969 | deallocate (alt_w_c) |
---|
1970 | deallocate (alt_u_c) |
---|
1971 | deallocate (pro_w_c) |
---|
1972 | deallocate (pro_u_c) |
---|
1973 | |
---|
1974 | deallocate (alt_w_n) |
---|
1975 | deallocate (alt_u_n) |
---|
1976 | deallocate (pro_w_n) |
---|
1977 | deallocate (pro_u_n) |
---|
1978 | |
---|
1979 | |
---|
1980 | END SUBROUTINE vertical_interp |
---|
1981 | |
---|
1982 | !!!!!!!!!!!!!!!!!!!!!!!!11 |
---|
1983 | SUBROUTINE inter(pro_c,alt_c,kde_c,pro_n,alt_n,kde_n) |
---|
1984 | |
---|
1985 | IMPLICIT NONE |
---|
1986 | INTEGER , INTENT(IN) :: kde_c,kde_n |
---|
1987 | REAL , DIMENSION(kde_c) , INTENT(IN ) :: pro_c,alt_c |
---|
1988 | REAL , DIMENSION(kde_n) , INTENT(IN ) :: alt_n |
---|
1989 | REAL , DIMENSION(kde_n) , INTENT(OUT) :: pro_n |
---|
1990 | |
---|
1991 | real ,dimension(kde_c) :: a,b,c,d |
---|
1992 | real :: p |
---|
1993 | integer :: i,j |
---|
1994 | |
---|
1995 | |
---|
1996 | call coeff_mon(alt_c,pro_c,a,b,c,d,kde_c) |
---|
1997 | |
---|
1998 | do i = 1,kde_n-1 |
---|
1999 | |
---|
2000 | do j=1,kde_c-1 |
---|
2001 | |
---|
2002 | if ( (alt_n(i) .ge. alt_c(j)).and.(alt_n(i) .lt. alt_c(j+1)) ) then |
---|
2003 | p = alt_n(i)-alt_c(j) |
---|
2004 | pro_n(i) = p*( p*(a(j)*p+b(j))+c(j)) + d(j) |
---|
2005 | goto 20 |
---|
2006 | endif |
---|
2007 | enddo |
---|
2008 | 20 continue |
---|
2009 | enddo |
---|
2010 | |
---|
2011 | pro_n(kde_n) = pro_c(kde_c) |
---|
2012 | |
---|
2013 | |
---|
2014 | END SUBROUTINE inter |
---|
2015 | |
---|
2016 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!11 |
---|
2017 | |
---|
2018 | subroutine coeff_mon(x,y,a,b,c,d,n) |
---|
2019 | |
---|
2020 | implicit none |
---|
2021 | |
---|
2022 | integer :: n |
---|
2023 | real ,dimension(n) :: x,y,a,b,c,d |
---|
2024 | real ,dimension(n) :: h,s,p,yp |
---|
2025 | |
---|
2026 | integer :: i |
---|
2027 | |
---|
2028 | |
---|
2029 | do i=1,n-1 |
---|
2030 | h(i) = (x(i+1)-x(i)) |
---|
2031 | s(i) = (y(i+1)-y(i)) / h(i) |
---|
2032 | enddo |
---|
2033 | |
---|
2034 | do i=2,n-1 |
---|
2035 | p(i) = (s(i-1)*h(i)+s(i)*h(i-1)) / (h(i-1)+h(i)) |
---|
2036 | enddo |
---|
2037 | |
---|
2038 | p(1) = s(1) |
---|
2039 | p(n) = s(n-1) |
---|
2040 | |
---|
2041 | do i=1,n |
---|
2042 | yp(i) = p(i) |
---|
2043 | enddo |
---|
2044 | !!!!!!!!!!!!!!!!!!!!! |
---|
2045 | |
---|
2046 | do i=2,n-1 |
---|
2047 | yp(i) = (sign(1.,s(i-1))+sign(1.,s(i)))* min( abs(s(i-1)),abs(s(i)),0.5*abs(p(i))) |
---|
2048 | enddo |
---|
2049 | |
---|
2050 | do i = 1,n-1 |
---|
2051 | a(i) = (yp(i)+yp(i+1)-2.*s(i))/(h(i)*h(i)) |
---|
2052 | b(i) = (3.*s(i)-2.*yp(i)-yp(i+1))/h(i) |
---|
2053 | c(i) = yp(i) |
---|
2054 | d(i) = y(i) |
---|
2055 | enddo |
---|
2056 | |
---|
2057 | end subroutine coeff_mon |
---|
2058 | |
---|
2059 | |
---|