Add new fields XLAT_C and XLONG_C
[WPS-merge.git] / metgrid / src / input_module.F
blobec9ec0f87857f76b8636c521d24820e9518fddf4
1 module input_module
3    use gridinfo_module
4    use misc_definitions_module
5    use module_debug
6 #ifdef IO_BINARY
7    use module_internal_header_util
8 #endif
9    use parallel_module
10    use queue_module
12    type (queue) :: unit_desc
14    ! WRF I/O API related variables
15    integer :: handle
17    integer :: num_calls
19    character (len=1) :: internal_gridtype
21    contains
24    subroutine input_init(nest_number, istatus)
26       implicit none
27   
28       ! Arguments
29       integer, intent(in) :: nest_number
30       integer, intent(out) :: istatus
31   
32 #include "wrf_io_flags.h"
33 #include "wrf_status_codes.h"
34   
35       ! Local variables
36       integer :: i
37       integer :: comm_1, comm_2
38       character (len=MAX_FILENAME_LEN) :: input_fname
39   
40       istatus = 0
41   
42       if (my_proc_id == IO_NODE .or. do_tiled_input) then
43   
44 #ifdef IO_BINARY
45          if (io_form_input == BINARY) call ext_int_ioinit('sysdep info', istatus)
46 #endif
47 #ifdef IO_NETCDF
48          if (io_form_input == NETCDF) call ext_ncd_ioinit('sysdep info', istatus)
49 #endif
50 #ifdef IO_GRIB1
51          if (io_form_input == GRIB1) call ext_gr1_ioinit('sysdep info', istatus)
52 #endif
53          call mprintf((istatus /= 0),ERROR,'Error in ext_pkg_ioinit')
54      
55          comm_1 = 1
56          comm_2 = 1
57          input_fname = ' '
58          if (gridtype == 'C') then
59 #ifdef IO_BINARY
60             if (io_form_input == BINARY) input_fname = trim(opt_output_from_geogrid_path)//'geo_em.d  .int'
61 #endif
62 #ifdef IO_NETCDF
63             if (io_form_input == NETCDF) input_fname = trim(opt_output_from_geogrid_path)//'geo_em.d  .nc'
64 #endif
65 #ifdef IO_GRIB1
66             if (io_form_input == GRIB1) input_fname = trim(opt_output_from_geogrid_path)//'geo_em.d  .grib'
67 #endif
68             i = len_trim(opt_output_from_geogrid_path)
69             write(input_fname(i+9:i+10),'(i2.2)') nest_number
70          else if (gridtype == 'E') then
71 #ifdef IO_BINARY
72             if (io_form_input == BINARY) input_fname = trim(opt_output_from_geogrid_path)//'geo_nmm.d  .int'
73 #endif
74 #ifdef IO_NETCDF
75             if (io_form_input == NETCDF) input_fname = trim(opt_output_from_geogrid_path)//'geo_nmm.d  .nc'
76 #endif
77 #ifdef IO_GRIB1
78             if (io_form_input == GRIB1) input_fname = trim(opt_output_from_geogrid_path)//'geo_nmm.d  .grib'
79 #endif
80             i = len_trim(opt_output_from_geogrid_path)
81             write(input_fname(i+10:i+11),'(i2.2)') nest_number
82          end if
84          if (nprocs > 1 .and. do_tiled_input) then
85             write(input_fname(len_trim(input_fname)+1:len_trim(input_fname)+5), '(a1,i4.4)') &
86                             '_', my_proc_id
87          end if
88      
89          istatus = 0
90 #ifdef IO_BINARY
91          if (io_form_input == BINARY) &
92             call ext_int_open_for_read(trim(input_fname), comm_1, comm_2, 'sysdep info', handle, istatus)
93 #endif
94 #ifdef IO_NETCDF
95          if (io_form_input == NETCDF) &
96             call ext_ncd_open_for_read(trim(input_fname), comm_1, comm_2, 'sysdep info', handle, istatus)
97 #endif
98 #ifdef IO_GRIB1
99          if (io_form_input == GRIB1) &
100             call ext_gr1_open_for_read(trim(input_fname), comm_1, comm_2, 'sysdep info', handle, istatus)
101 #endif
102          call mprintf((istatus /= 0),ERROR,'Couldn''t open file %s for input.',s1=input_fname)
103      
104          call q_init(unit_desc)
105   
106       end if ! (my_proc_id == IO_NODE .or. do_tiled_input)
107   
108       num_calls = 0
110    end subroutine input_init
113    subroutine read_next_field(start_patch_i, end_patch_i, &
114                               start_patch_j, end_patch_j, &
115                               start_patch_k, end_patch_k, &
116                               cname, cunits, cdesc, memorder, stagger, &
117                               dimnames, sr_x, sr_y, real_array, istatus)
119       implicit none
120   
121       ! Arguments
122       integer, intent(out) :: start_patch_i, end_patch_i, &
123                               start_patch_j, end_patch_j, &
124                               start_patch_k, end_patch_k, &
125                               sr_x, sr_y
126       real, pointer, dimension(:,:,:) :: real_array
127       character (len=*), intent(out) :: cname, memorder, stagger, cunits, cdesc
128       character (len=128), dimension(3), intent(inout) :: dimnames
129       integer, intent(inout) :: istatus
130   
131 #include "wrf_io_flags.h"
132 #include "wrf_status_codes.h"
133   
134       ! Local variables
135       integer :: ndim, wrftype
136       integer :: sm1, em1, sm2, em2, sm3, em3, sp1, ep1, sp2, ep2, sp3, ep3
137       integer, dimension(3) :: domain_start, domain_end, temp
138       real, pointer, dimension(:,:,:) :: real_domain
139       character (len=20) :: datestr
140       type (q_data) :: qd
141   
142       if (my_proc_id == IO_NODE .or. do_tiled_input) then
143   
144          if (num_calls == 0) then
145 #ifdef IO_BINARY
146             if (io_form_input == BINARY) call ext_int_get_next_time(handle, datestr, istatus)
147 #endif
148 #ifdef IO_NETCDF
149             if (io_form_input == NETCDF) call ext_ncd_get_next_time(handle, datestr, istatus)
150 #endif
151 #ifdef IO_GRIB1
152             if (io_form_input == GRIB1) call ext_gr1_get_next_time(handle, datestr, istatus)
153 #endif
154          end if
155      
156          num_calls = num_calls + 1
157    
158 #ifdef IO_BINARY
159          if (io_form_input == BINARY) call ext_int_get_next_var(handle, cname, istatus) 
160 #endif
161 #ifdef IO_NETCDF
162          if (io_form_input == NETCDF) call ext_ncd_get_next_var(handle, cname, istatus) 
163 #endif
164 #ifdef IO_GRIB1
165          if (io_form_input == GRIB1) call ext_gr1_get_next_var(handle, cname, istatus) 
166 #endif
167       end if
168   
169       if (nprocs > 1 .and. .not. do_tiled_input) call parallel_bcast_int(istatus)
170       if (istatus /= 0) return
171   
172       if (my_proc_id == IO_NODE .or. do_tiled_input) then
173   
174          istatus = 0
175 #ifdef IO_BINARY
176          if (io_form_input == BINARY) then
177             call ext_int_get_var_info(handle, cname, ndim, memorder, stagger, domain_start, domain_end, wrftype, istatus)
178          end if
179 #endif
180 #ifdef IO_NETCDF
181          if (io_form_input == NETCDF) then
182             call ext_ncd_get_var_info(handle, cname, ndim, memorder, stagger, domain_start, domain_end, wrftype, istatus)
183             call ext_ncd_get_var_ti_integer(handle, 'sr_x', &
184                                             trim(cname), temp(1), 1, temp(3), istatus)
185             call ext_ncd_get_var_ti_integer(handle, 'sr_y', &
186                                             trim(cname), temp(2), 1, temp(3), istatus)
187          end if
188 #endif
189 #ifdef IO_GRIB1
190          if (io_form_input == GRIB1) then
191             call ext_gr1_get_var_info(handle, cname, ndim, memorder, stagger, domain_start, domain_end, wrftype, istatus)
192             call ext_gr1_get_var_ti_integer(handle, 'sr_x', &
193                                             trim(cname), temp(1), 1, temp(3), istatus)
194             call ext_gr1_get_var_ti_integer(handle, 'sr_y', &
195                                             trim(cname), temp(2), 1, temp(3), istatus)
196          end if
197 #endif
198      
199          call mprintf((istatus /= 0),ERROR,'In read_next_field(), problems with ext_pkg_get_var_info()')
201          start_patch_i = domain_start(1) 
202          start_patch_j = domain_start(2) 
203          end_patch_i = domain_end(1)
204          end_patch_j = domain_end(2)
205          if (ndim == 3) then
206             start_patch_k = domain_start(3) 
207             end_patch_k = domain_end(3) 
208          else
209             domain_start(3) = 1
210             domain_end(3) = 1
211             start_patch_k = 1
212             end_patch_k = 1
213          end if
214      
215          nullify(real_domain)
216      
217          allocate(real_domain(start_patch_i:end_patch_i, start_patch_j:end_patch_j, start_patch_k:end_patch_k))
218 #ifdef IO_BINARY
219          if (io_form_input == BINARY) then
220             call ext_int_read_field(handle, '0000-00-00_00:00:00', cname, real_domain, WRF_REAL, &
221                           1, 1, 0, memorder, stagger, &
222                           dimnames, domain_start, domain_end, domain_start, domain_end, &
223                           domain_start, domain_end, istatus)
224          end if
225 #endif
226 #ifdef IO_NETCDF
227          if (io_form_input == NETCDF) then
228             call ext_ncd_read_field(handle, '0000-00-00_00:00:00', cname, real_domain, WRF_REAL, &
229                           1, 1, 0, memorder, stagger, &
230                           dimnames, domain_start, domain_end, domain_start, domain_end, &
231                           domain_start, domain_end, istatus)
232          end if
233 #endif
234 #ifdef IO_GRIB1
235          if (io_form_input == GRIB1) then
236             call ext_gr1_read_field(handle, '0000-00-00_00:00:00', cname, real_domain, WRF_REAL, &
237                           1, 1, 0, memorder, stagger, &
238                           dimnames, domain_start, domain_end, domain_start, domain_end, &
239                           domain_start, domain_end, istatus)
240          end if
241 #endif
242      
243          call mprintf((istatus /= 0),ERROR,'In read_next_field(), got error code %i.', i1=istatus)
245          if (io_form_input == BINARY) then
246             qd = q_remove(unit_desc)
247             cunits = qd%units
248             cdesc = qd%description
249             stagger = qd%stagger
250             sr_x = qd%sr_x
251             sr_y = qd%sr_y
252          else
253             cunits = ' '
254             cdesc = ' '
255             stagger = ' '
256             sr_x = temp(1)
257             sr_y = temp(2)
258         
259 #ifdef IO_NETCDF
260             if (io_form_input == NETCDF) then
261                call ext_ncd_get_var_ti_char(handle, 'units', cname, cunits, istatus)
262                call ext_ncd_get_var_ti_char(handle, 'description', cname, cdesc, istatus)
263                call ext_ncd_get_var_ti_char(handle, 'stagger', cname, stagger, istatus)
264             end if
265 #endif
266 #ifdef IO_GRIB1
267             if (io_form_input == GRIB1) then
268                call ext_gr1_get_var_ti_char(handle, 'units', cname, cunits, istatus)
269                call ext_gr1_get_var_ti_char(handle, 'description', cname, cdesc, istatus)
270                call ext_gr1_get_var_ti_char(handle, 'stagger', cname, stagger, istatus)
271             end if
272 #endif
273          end if
274   
275       end if ! (my_proc_id == IO_NODE .or. do_tiled_input)
277       if (nprocs > 1 .and. .not. do_tiled_input) then
278          call parallel_bcast_char(cname, len(cname))
279          call parallel_bcast_char(cunits, len(cunits))
280          call parallel_bcast_char(cdesc, len(cdesc))
281          call parallel_bcast_char(memorder, len(memorder))
282          call parallel_bcast_char(stagger, len(stagger))
283          call parallel_bcast_char(dimnames(1), 128)
284          call parallel_bcast_char(dimnames(2), 128)
285          call parallel_bcast_char(dimnames(3), 128)
286          call parallel_bcast_int(domain_start(3))
287          call parallel_bcast_int(domain_end(3))
288          call parallel_bcast_int(sr_x)
289          call parallel_bcast_int(sr_y)
290    
291          sp1 = my_minx
292          ep1 = my_maxx - 1
293          sp2 = my_miny
294          ep2 = my_maxy - 1
295          sp3 = domain_start(3)
296          ep3 = domain_end(3)
297    
298          if (internal_gridtype == 'C') then
299             if (my_x /= nproc_x - 1 .or. stagger == 'U' .or. stagger == 'CORNER' .or. sr_x > 1) then
300                ep1 = ep1 + 1
301             end if
302             if (my_y /= nproc_y - 1 .or. stagger == 'V' .or. stagger == 'CORNER' .or. sr_y > 1) then
303                ep2 = ep2 + 1
304             end if
305          else if (internal_gridtype == 'E') then
306             ep1 = ep1 + 1
307             ep2 = ep2 + 1
308          end if
309    
310          if (sr_x > 1) then
311             sp1 = (sp1-1)*sr_x+1
312             ep1 =  ep1   *sr_x
313          end if
314          if (sr_y > 1) then
315             sp2 = (sp2-1)*sr_y+1
316             ep2 =  ep2   *sr_y
317          end if
319          sm1 = sp1
320          em1 = ep1
321          sm2 = sp2
322          em2 = ep2
323          sm3 = sp3
324          em3 = ep3
325    
326          start_patch_i = sp1
327          end_patch_i   = ep1
328          start_patch_j = sp2
329          end_patch_j   = ep2
330          start_patch_k = sp3
331          end_patch_k   = ep3
332    
333          allocate(real_array(sm1:em1,sm2:em2,sm3:em3))
334          if (my_proc_id /= IO_NODE) then
335             allocate(real_domain(1,1,1))
336             domain_start(1) = 1
337             domain_start(2) = 1
338             domain_start(3) = 1
339             domain_end(1) = 1
340             domain_end(2) = 1
341             domain_end(3) = 1
342          end if
343          call scatter_whole_field_r(real_array, &
344                                    sm1, em1, sm2, em2, sm3, em3, &
345                                    sp1, ep1, sp2, ep2, sp3, ep3, &
346                                    real_domain, &
347                                    domain_start(1), domain_end(1), &
348                                    domain_start(2), domain_end(2), &
349                                    domain_start(3), domain_end(3))
350          deallocate(real_domain)
352       else
353   
354          real_array => real_domain
356       end if
358    end subroutine read_next_field
360    subroutine read_global_attrs(title, start_date, grid_type, dyn_opt,                                &
361                                 west_east_dim, south_north_dim, bottom_top_dim,                       &
362                                 we_patch_s, we_patch_e, we_patch_s_stag, we_patch_e_stag,             &
363                                 sn_patch_s, sn_patch_e, sn_patch_s_stag, sn_patch_e_stag,             &
364                                 map_proj, mminlu, num_land_cat, is_water, is_lake, is_ice, is_urban,  &
365                                 isoilwater, grid_id, parent_id, i_parent_start, j_parent_start,       &
366                                 i_parent_end, j_parent_end, dx, dy, cen_lat, moad_cen_lat, cen_lon,   &
367                                 stand_lon, truelat1, truelat2, pole_lat, pole_lon, parent_grid_ratio, &
368                                 corner_lats, corner_lons, sr_x, sr_y)
370       implicit none
371   
372       ! Arguments
373       integer, intent(out) :: dyn_opt, west_east_dim, south_north_dim, bottom_top_dim, map_proj,   &
374                  is_water, is_lake, we_patch_s, we_patch_e, we_patch_s_stag, we_patch_e_stag,      &
375                  sn_patch_s, sn_patch_e, sn_patch_s_stag, sn_patch_e_stag,                         &
376                  is_ice, is_urban, isoilwater, grid_id, parent_id, i_parent_start, j_parent_start, &
377                  i_parent_end, j_parent_end, parent_grid_ratio, sr_x, sr_y, num_land_cat
378       real, intent(out) :: dx, dy, cen_lat, moad_cen_lat, cen_lon, stand_lon, truelat1, truelat2,  &
379                  pole_lat, pole_lon
380       real, dimension(16), intent(out) :: corner_lats, corner_lons
381       character (len=128), intent(out) :: title, start_date, grid_type, mminlu
382   
383       ! Local variables
384       integer :: istatus, i
385       real :: wps_version
386       character (len=128) :: cunits, cdesc, cstagger
387       integer, dimension(3) :: sr
388       type (q_data) :: qd
389   
390       if (my_proc_id == IO_NODE .or. do_tiled_input) then
391   
392 #ifdef IO_BINARY
393          if (io_form_input == BINARY) then
394             istatus = 0
395             do while (istatus == 0) 
396                cunits = ' '
397                cdesc = ' '
398                cstagger = ' '
399                sr = 0
400                call ext_int_get_var_ti_char(handle, 'units', 'VAR', cunits, istatus)
401          
402                if (istatus == 0) then
403                   call ext_int_get_var_ti_char(handle, 'description', 'VAR', cdesc, istatus)
404           
405                   if (istatus == 0) then
406                      call ext_int_get_var_ti_char(handle, 'stagger', 'VAR', cstagger, istatus)
408                     if (istatus == 0) then
409                     call ext_int_get_var_ti_integer(handle, 'sr_x', 'VAR', sr(1), 1, sr(3), istatus)
411                         if (istatus == 0) then
412                         call ext_int_get_var_ti_integer(handle, 'sr_y', 'VAR', sr(2), 1, sr(3), istatus)
413          
414                              qd%units = cunits
415                              qd%description = cdesc
416                              qd%stagger = cstagger
417                              qd%sr_x = sr(1)
418                              qd%sr_y = sr(2)
419                              call q_insert(unit_desc, qd)
421                         end if
422                     end if
423                   end if
424                end if
425             end do
426          end if
427 #endif
428      
429          call ext_get_dom_ti_char          ('TITLE', title)
430          if (index(title,'GEOGRID V3.7.1') /= 0) then
431             wps_version = 3.71
432          else if (index(title,'GEOGRID V3.7') /= 0) then
433             wps_version = 3.7
434          else if (index(title,'GEOGRID V3.6.1') /= 0) then
435             wps_version = 3.61
436          else if (index(title,'GEOGRID V3.6') /= 0) then
437             wps_version = 3.6
438          else if (index(title,'GEOGRID V3.5.1') /= 0) then
439             wps_version = 3.51
440          else if (index(title,'GEOGRID V3.5') /= 0) then
441             wps_version = 3.5
442          else if (index(title,'GEOGRID V3.4.1') /= 0) then
443             wps_version = 3.41
444          else if (index(title,'GEOGRID V3.4') /= 0) then
445             wps_version = 3.4
446          else if (index(title,'GEOGRID V3.3.1') /= 0) then
447             wps_version = 3.31
448          else if (index(title,'GEOGRID V3.3') /= 0) then
449             wps_version = 3.3
450          else if (index(title,'GEOGRID V3.2.1') /= 0) then
451             wps_version = 3.21
452          else if (index(title,'GEOGRID V3.2') /= 0) then
453             wps_version = 3.2
454          else if (index(title,'GEOGRID V3.1.1') /= 0) then
455             wps_version = 3.11
456          else if (index(title,'GEOGRID V3.1') /= 0) then
457             wps_version = 3.1
458          else if (index(title,'GEOGRID V3.0.1') /= 0) then
459             wps_version = 3.01
460          else
461             wps_version = 3.0
462          end if
463          call mprintf(.true.,DEBUG,'Reading static data from WPS version %f', f1=wps_version)
464          call ext_get_dom_ti_char          ('SIMULATION_START_DATE', start_date)
465          call ext_get_dom_ti_integer_scalar('WEST-EAST_GRID_DIMENSION', west_east_dim)
466          call ext_get_dom_ti_integer_scalar('SOUTH-NORTH_GRID_DIMENSION', south_north_dim)
467          call ext_get_dom_ti_integer_scalar('BOTTOM-TOP_GRID_DIMENSION', bottom_top_dim)
468          call ext_get_dom_ti_integer_scalar('WEST-EAST_PATCH_START_UNSTAG', we_patch_s)
469          call ext_get_dom_ti_integer_scalar('WEST-EAST_PATCH_END_UNSTAG', we_patch_e)
470          call ext_get_dom_ti_integer_scalar('WEST-EAST_PATCH_START_STAG', we_patch_s_stag)
471          call ext_get_dom_ti_integer_scalar('WEST-EAST_PATCH_END_STAG', we_patch_e_stag)
472          call ext_get_dom_ti_integer_scalar('SOUTH-NORTH_PATCH_START_UNSTAG', sn_patch_s)
473          call ext_get_dom_ti_integer_scalar('SOUTH-NORTH_PATCH_END_UNSTAG', sn_patch_e)
474          call ext_get_dom_ti_integer_scalar('SOUTH-NORTH_PATCH_START_STAG', sn_patch_s_stag)
475          call ext_get_dom_ti_integer_scalar('SOUTH-NORTH_PATCH_END_STAG', sn_patch_e_stag)
476          call ext_get_dom_ti_char          ('GRIDTYPE', grid_type)
477          call ext_get_dom_ti_real_scalar   ('DX', dx)
478          call ext_get_dom_ti_real_scalar   ('DY', dy)
479          call ext_get_dom_ti_integer_scalar('DYN_OPT', dyn_opt)
480          call ext_get_dom_ti_real_scalar   ('CEN_LAT', cen_lat)
481          call ext_get_dom_ti_real_scalar   ('CEN_LON', cen_lon)
482          call ext_get_dom_ti_real_scalar   ('TRUELAT1', truelat1)
483          call ext_get_dom_ti_real_scalar   ('TRUELAT2', truelat2)
484          call ext_get_dom_ti_real_scalar   ('MOAD_CEN_LAT', moad_cen_lat)
485          call ext_get_dom_ti_real_scalar   ('STAND_LON', stand_lon)
486          call ext_get_dom_ti_real_scalar   ('POLE_LAT', pole_lat)
487          call ext_get_dom_ti_real_scalar   ('POLE_LON', pole_lon)
488          call ext_get_dom_ti_real_vector   ('corner_lats', corner_lats, 16)
489          call ext_get_dom_ti_real_vector   ('corner_lons', corner_lons, 16)
490          call ext_get_dom_ti_integer_scalar('MAP_PROJ', map_proj)
491          call ext_get_dom_ti_char          ('MMINLU', mminlu)
492          if ( wps_version >= 3.01 ) then
493             call ext_get_dom_ti_integer_scalar('NUM_LAND_CAT', num_land_cat)
494          else
495             num_land_cat = 24
496          end if
497          call ext_get_dom_ti_integer_scalar('ISWATER', is_water)
498          if ( wps_version >= 3.01 ) then
499             call ext_get_dom_ti_integer_scalar('ISLAKE', is_lake)
500          else
501             is_lake = -1
502          end if
503          call ext_get_dom_ti_integer_scalar('ISICE', is_ice)
504          call ext_get_dom_ti_integer_scalar('ISURBAN', is_urban)
505          call ext_get_dom_ti_integer_scalar('ISOILWATER', isoilwater)
506          call ext_get_dom_ti_integer_scalar('grid_id', grid_id)
507          call ext_get_dom_ti_integer_scalar('parent_id', parent_id)
508          call ext_get_dom_ti_integer_scalar('i_parent_start', i_parent_start)
509          call ext_get_dom_ti_integer_scalar('j_parent_start', j_parent_start)
510          call ext_get_dom_ti_integer_scalar('i_parent_end', i_parent_end)
511          call ext_get_dom_ti_integer_scalar('j_parent_end', j_parent_end)
512          call ext_get_dom_ti_integer_scalar('parent_grid_ratio', parent_grid_ratio)
513          call ext_get_dom_ti_integer_scalar('sr_x', sr_x)
514          call ext_get_dom_ti_integer_scalar('sr_y', sr_y)
515    
516       end if
518   
519       if (nprocs > 1 .and. .not. do_tiled_input) then
520   
521          call parallel_bcast_char(title, len(title))
522          call parallel_bcast_char(start_date, len(start_date))
523          call parallel_bcast_char(grid_type, len(grid_type))
524          call parallel_bcast_int(west_east_dim)
525          call parallel_bcast_int(south_north_dim)
526          call parallel_bcast_int(bottom_top_dim)
527          call parallel_bcast_int(we_patch_s)
528          call parallel_bcast_int(we_patch_e)
529          call parallel_bcast_int(we_patch_s_stag)
530          call parallel_bcast_int(we_patch_e_stag)
531          call parallel_bcast_int(sn_patch_s)
532          call parallel_bcast_int(sn_patch_e)
533          call parallel_bcast_int(sn_patch_s_stag)
534          call parallel_bcast_int(sn_patch_e_stag)
535          call parallel_bcast_int(sr_x)
536          call parallel_bcast_int(sr_y)
538          ! Must figure out patch dimensions from info in parallel module
539 !         we_patch_s      = my_minx
540 !         we_patch_s_stag = my_minx
541 !         we_patch_e      = my_maxx - 1
542 !         sn_patch_s      = my_miny
543 !         sn_patch_s_stag = my_miny
544 !         sn_patch_e      = my_maxy - 1
546 !         if (trim(grid_type) == 'C') then
547 !            if (my_x /= nproc_x - 1) then
548 !               we_patch_e_stag = we_patch_e + 1
549 !            end if
550 !            if (my_y /= nproc_y - 1) then
551 !               sn_patch_e_stag = sn_patch_e + 1
552 !            end if
553 !         else if (trim(grid_type) == 'E') then
554 !            we_patch_e = we_patch_e + 1
555 !            sn_patch_e = sn_patch_e + 1
556 !            we_patch_e_stag = we_patch_e
557 !            sn_patch_e_stag = sn_patch_e
558 !         end if
560          call parallel_bcast_real(dx)
561          call parallel_bcast_real(dy)
562          call parallel_bcast_int(dyn_opt)
563          call parallel_bcast_real(cen_lat)
564          call parallel_bcast_real(cen_lon)
565          call parallel_bcast_real(truelat1)
566          call parallel_bcast_real(truelat2)
567          call parallel_bcast_real(pole_lat)
568          call parallel_bcast_real(pole_lon)
569          call parallel_bcast_real(moad_cen_lat)
570          call parallel_bcast_real(stand_lon)
571          do i=1,16
572             call parallel_bcast_real(corner_lats(i))
573             call parallel_bcast_real(corner_lons(i))
574          end do
575          call parallel_bcast_int(map_proj)
576          call parallel_bcast_char(mminlu, len(mminlu))
577          call parallel_bcast_int(is_water)
578          call parallel_bcast_int(is_lake)
579          call parallel_bcast_int(is_ice)
580          call parallel_bcast_int(is_urban)
581          call parallel_bcast_int(isoilwater)
582          call parallel_bcast_int(grid_id)
583          call parallel_bcast_int(parent_id)
584          call parallel_bcast_int(i_parent_start)
585          call parallel_bcast_int(i_parent_end)
586          call parallel_bcast_int(j_parent_start)
587          call parallel_bcast_int(j_parent_end)
588          call parallel_bcast_int(parent_grid_ratio)
589       end if
590   
591       internal_gridtype = grid_type
593    end subroutine read_global_attrs
596    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
597    ! Name: ext_get_dom_ti_integer
598    !
599    ! Purpose: Read a domain time-independent integer attribute from input.
600    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
601    subroutine ext_get_dom_ti_integer_scalar(var_name, var_value, suppress_errors)
603       implicit none
605       ! Arguments
606       integer, intent(out) :: var_value
607       character (len=*), intent(in) :: var_name
608       logical, intent(in), optional :: suppress_errors
610       ! Local variables
611       integer :: istatus, outcount
613 #ifdef IO_BINARY
614       if (io_form_input == BINARY) then
615          call ext_int_get_dom_ti_integer(handle, trim(var_name), &
616                                          var_value, &
617                                          1, outcount, istatus)
618       end if
619 #endif
620 #ifdef IO_NETCDF
621       if (io_form_input == NETCDF) then
622          call ext_ncd_get_dom_ti_integer(handle, trim(var_name), &
623                                          var_value, &
624                                          1, outcount, istatus)
625       end if
626 #endif
627 #ifdef IO_GRIB1
628       if (io_form_input == GRIB1) then
629          call ext_gr1_get_dom_ti_integer(handle, trim(var_name), &
630                                          var_value, &
631                                          1, outcount, istatus)
632       end if
633 #endif
635       if (present(suppress_errors)) then
636          call mprintf((istatus /= 0 .and. .not.suppress_errors),ERROR,'Error while reading domain time-independent attribute.')
637       else
638          call mprintf((istatus /= 0),ERROR,'Error while reading domain time-independent attribute.')
639       end if
641    end subroutine ext_get_dom_ti_integer_scalar
644    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
645    ! Name: ext_get_dom_ti_integer
646    !
647    ! Purpose: Read a domain time-independent integer attribute from input.
648    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
649    subroutine ext_get_dom_ti_integer_vector(var_name, var_value, n)
651       implicit none
653       ! Arguments
654       integer, intent(in) :: n
655       integer, dimension(n), intent(out) :: var_value
656       character (len=*), intent(in) :: var_name
658       ! Local variables
659       integer :: istatus, outcount
661 #ifdef IO_BINARY
662       if (io_form_input == BINARY) then
663          call ext_int_get_dom_ti_integer(handle, trim(var_name), &
664                                          var_value, &
665                                          n, outcount, istatus)
666       end if
667 #endif
668 #ifdef IO_NETCDF
669       if (io_form_input == NETCDF) then
670          call ext_ncd_get_dom_ti_integer(handle, trim(var_name), &
671                                          var_value, &
672                                          n, outcount, istatus)
673       end if
674 #endif
675 #ifdef IO_GRIB1
676       if (io_form_input == GRIB1) then
677          call ext_gr1_get_dom_ti_integer(handle, trim(var_name), &
678                                          var_value, &
679                                          n, outcount, istatus)
680       end if
681 #endif
683       call mprintf((istatus /= 0),ERROR,'Error while reading domain time-independent attribute.')
685    end subroutine ext_get_dom_ti_integer_vector
688    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
689    ! Name: ext_get_dom_ti_real
690    !
691    ! Purpose: Read a domain time-independent real attribute from input.
692    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
693    subroutine ext_get_dom_ti_real_scalar(var_name, var_value)
695       implicit none
697       ! Arguments
698       real, intent(out) :: var_value
699       character (len=*), intent(in) :: var_name
701       ! Local variables
702       integer :: istatus, outcount
704 #ifdef IO_BINARY
705       if (io_form_input == BINARY) then
706          call ext_int_get_dom_ti_real(handle, trim(var_name), &
707                                          var_value, &
708                                          1, outcount, istatus)
709       end if
710 #endif
711 #ifdef IO_NETCDF
712       if (io_form_input == NETCDF) then
713          call ext_ncd_get_dom_ti_real(handle, trim(var_name), &
714                                          var_value, &
715                                          1, outcount, istatus)
716       end if
717 #endif
718 #ifdef IO_GRIB1
719       if (io_form_input == GRIB1) then
720          call ext_gr1_get_dom_ti_real(handle, trim(var_name), &
721                                          var_value, &
722                                          1, outcount, istatus)
723       end if
724 #endif
726       call mprintf((istatus /= 0),ERROR,'Error while reading domain time-independent attribute.')
728    end subroutine ext_get_dom_ti_real_scalar
731    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
732    ! Name: ext_get_dom_ti_real
733    !
734    ! Purpose: Read a domain time-independent real attribute from input.
735    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
736    subroutine ext_get_dom_ti_real_vector(var_name, var_value, n)
738       implicit none
740       ! Arguments
741       integer, intent(in) :: n
742       real, dimension(n), intent(out) :: var_value
743       character (len=*), intent(in) :: var_name
745       ! Local variables
746       integer :: istatus, outcount
748 #ifdef IO_BINARY
749       if (io_form_input == BINARY) then
750          call ext_int_get_dom_ti_real(handle, trim(var_name), &
751                                          var_value, &
752                                          n, outcount, istatus)
753       end if
754 #endif
755 #ifdef IO_NETCDF
756       if (io_form_input == NETCDF) then
757          call ext_ncd_get_dom_ti_real(handle, trim(var_name), &
758                                          var_value, &
759                                          n, outcount, istatus)
760       end if
761 #endif
762 #ifdef IO_GRIB1
763       if (io_form_input == GRIB1) then
764          call ext_gr1_get_dom_ti_real(handle, trim(var_name), &
765                                          var_value, &
766                                          n, outcount, istatus)
767       end if
768 #endif
770       call mprintf((istatus /= 0),ERROR,'Error while reading domain time-independent attribute.')
772    end subroutine ext_get_dom_ti_real_vector
775    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
776    ! Name: ext_get_dom_ti_char
777    !
778    ! Purpose: Read a domain time-independent character attribute from input.
779    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
780    subroutine ext_get_dom_ti_char(var_name, var_value)
782       implicit none
784       ! Arguments
785       character (len=*), intent(in) :: var_name
786       character (len=128), intent(out) :: var_value
788       ! Local variables
789       integer :: istatus
791 #ifdef IO_BINARY
792       if (io_form_input == BINARY) then
793          call ext_int_get_dom_ti_char(handle, trim(var_name), &
794                                          var_value, &
795                                          istatus)
796       end if
797 #endif
798 #ifdef IO_NETCDF
799       if (io_form_input == NETCDF) then
800          call ext_ncd_get_dom_ti_char(handle, trim(var_name), &
801                                          var_value, &
802                                          istatus)
803       end if
804 #endif
805 #ifdef IO_GRIB1
806       if (io_form_input == GRIB1) then
807          call ext_gr1_get_dom_ti_char(handle, trim(var_name), &
808                                          var_value, &
809                                          istatus)
810       end if
811 #endif
813       call mprintf((istatus /= 0),ERROR,'Error in reading domain time-independent attribute')
815    end subroutine ext_get_dom_ti_char
818    subroutine input_close()
820       implicit none
821   
822       ! Local variables
823       integer :: istatus
824   
825       istatus = 0
826       if (my_proc_id == IO_NODE .or. do_tiled_input) then
827 #ifdef IO_BINARY
828          if (io_form_input == BINARY) then
829             call ext_int_ioclose(handle, istatus)
830             call ext_int_ioexit(istatus)
831          end if
832 #endif
833 #ifdef IO_NETCDF
834          if (io_form_input == NETCDF) then
835             call ext_ncd_ioclose(handle, istatus)
836             call ext_ncd_ioexit(istatus)
837          end if
838 #endif
839 #ifdef IO_GRIB1
840          if (io_form_input == GRIB1) then
841             call ext_gr1_ioclose(handle, istatus)
842             call ext_gr1_ioexit(istatus)
843          end if
844 #endif
845       end if
847       call q_destroy(unit_desc)
849    end subroutine input_close
851 end module input_module