WPS updated to 3.2.1
[wrffire.git] / WPS / metgrid / src / input_module.F
blobb4ee3e2493e8bc0a424bf8bd6f3e2f094d66d60f
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             call ext_int_get_var_ti_integer(handle, 'sr_x', &
179                                             trim(cname), temp(1), 1, temp(3), istatus)
180             call ext_int_get_var_ti_integer(handle, 'sr_y', &
181                                             trim(cname), temp(2), 1, temp(3), istatus)
182          end if
183 #endif
184 #ifdef IO_NETCDF
185          if (io_form_input == NETCDF) then
186             call ext_ncd_get_var_info(handle, cname, ndim, memorder, stagger, domain_start, domain_end, wrftype, istatus)
187             call ext_ncd_get_var_ti_integer(handle, 'sr_x', &
188                                             trim(cname), temp(1), 1, temp(3), istatus)
189             call ext_ncd_get_var_ti_integer(handle, 'sr_y', &
190                                             trim(cname), temp(2), 1, temp(3), istatus)
191          end if
192 #endif
193 #ifdef IO_GRIB1
194          if (io_form_input == GRIB1) then
195             call ext_gr1_get_var_info(handle, cname, ndim, memorder, stagger, domain_start, domain_end, wrftype, istatus)
196             call ext_gr1_get_var_ti_integer(handle, 'sr_x', &
197                                             trim(cname), temp(1), 1, temp(3), istatus)
198             call ext_gr1_get_var_ti_integer(handle, 'sr_y', &
199                                             trim(cname), temp(2), 1, temp(3), istatus)
200          end if
201 #endif
202          sr_x = temp(1)
203          sr_y = temp(2)
204      
205          call mprintf((istatus /= 0),ERROR,'In read_next_field(), problems with ext_pkg_get_var_info()')
207          start_patch_i = domain_start(1) 
208          start_patch_j = domain_start(2) 
209          end_patch_i = domain_end(1)
210          end_patch_j = domain_end(2)
211          if (ndim == 3) then
212             start_patch_k = domain_start(3) 
213             end_patch_k = domain_end(3) 
214          else
215             domain_start(3) = 1
216             domain_end(3) = 1
217             start_patch_k = 1
218             end_patch_k = 1
219          end if
220      
221          nullify(real_domain)
222      
223          allocate(real_domain(start_patch_i:end_patch_i, start_patch_j:end_patch_j, start_patch_k:end_patch_k))
224 #ifdef IO_BINARY
225          if (io_form_input == BINARY) then
226             call ext_int_read_field(handle, '0000-00-00_00:00:00', cname, real_domain, WRF_REAL, &
227                           1, 1, 0, memorder, stagger, &
228                           dimnames, domain_start, domain_end, domain_start, domain_end, &
229                           domain_start, domain_end, istatus)
230          end if
231 #endif
232 #ifdef IO_NETCDF
233          if (io_form_input == NETCDF) then
234             call ext_ncd_read_field(handle, '0000-00-00_00:00:00', cname, real_domain, WRF_REAL, &
235                           1, 1, 0, memorder, stagger, &
236                           dimnames, domain_start, domain_end, domain_start, domain_end, &
237                           domain_start, domain_end, istatus)
238          end if
239 #endif
240 #ifdef IO_GRIB1
241          if (io_form_input == GRIB1) then
242             call ext_gr1_read_field(handle, '0000-00-00_00:00:00', cname, real_domain, WRF_REAL, &
243                           1, 1, 0, memorder, stagger, &
244                           dimnames, domain_start, domain_end, domain_start, domain_end, &
245                           domain_start, domain_end, istatus)
246          end if
247 #endif
248      
249          call mprintf((istatus /= 0),ERROR,'In read_next_field(), got error code %i.', i1=istatus)
250      
251          if (io_form_input == BINARY) then
252             qd = q_remove(unit_desc)
253             cunits = qd%units
254             cdesc = qd%description
255             stagger = qd%stagger
256          else
257             cunits = ' '
258             cdesc = ' '
259             stagger = ' '
260         
261 #ifdef IO_NETCDF
262             if (io_form_input == NETCDF) then
263                call ext_ncd_get_var_ti_char(handle, 'units', cname, cunits, istatus)
264                call ext_ncd_get_var_ti_char(handle, 'description', cname, cdesc, istatus)
265                call ext_ncd_get_var_ti_char(handle, 'stagger', cname, stagger, istatus)
266             end if
267 #endif
268 #ifdef IO_GRIB1
269             if (io_form_input == GRIB1) then
270                call ext_gr1_get_var_ti_char(handle, 'units', cname, cunits, istatus)
271                call ext_gr1_get_var_ti_char(handle, 'description', cname, cdesc, istatus)
272                call ext_gr1_get_var_ti_char(handle, 'stagger', cname, stagger, istatus)
273             end if
274 #endif
275          end if
276   
277       end if ! (my_proc_id == IO_NODE .or. do_tiled_input)
279       if (nprocs > 1 .and. .not. do_tiled_input) then
280          call parallel_bcast_char(cname, len(cname))
281          call parallel_bcast_char(cunits, len(cunits))
282          call parallel_bcast_char(cdesc, len(cdesc))
283          call parallel_bcast_char(memorder, len(memorder))
284          call parallel_bcast_char(stagger, len(stagger))
285          call parallel_bcast_char(dimnames(1), 128)
286          call parallel_bcast_char(dimnames(2), 128)
287          call parallel_bcast_char(dimnames(3), 128)
288          call parallel_bcast_int(domain_start(3))
289          call parallel_bcast_int(domain_end(3))
290          call parallel_bcast_int(sr_x)
291          call parallel_bcast_int(sr_y)
292    
293          sp1 = my_minx
294          ep1 = my_maxx - 1
295          sp2 = my_miny
296          ep2 = my_maxy - 1
297          sp3 = domain_start(3)
298          ep3 = domain_end(3)
299    
300          if (internal_gridtype == 'C') then
301             if (my_x /= nproc_x - 1 .or. stagger == 'U' .or. sr_x > 1) then
302                ep1 = ep1 + 1
303             end if
304             if (my_y /= nproc_y - 1 .or. stagger == 'V' .or. sr_y > 1) then
305                ep2 = ep2 + 1
306             end if
307          else if (internal_gridtype == 'E') then
308             ep1 = ep1 + 1
309             ep2 = ep2 + 1
310          end if
311    
312          if (sr_x > 1) then
313             sp1 = (sp1-1)*sr_x+1
314             ep1 =  ep1   *sr_x
315          end if
316          if (sr_y > 1) then
317             sp2 = (sp2-1)*sr_y+1
318             ep2 =  ep2   *sr_y
319          end if
321          sm1 = sp1
322          em1 = ep1
323          sm2 = sp2
324          em2 = ep2
325          sm3 = sp3
326          em3 = ep3
327    
328          start_patch_i = sp1
329          end_patch_i   = ep1
330          start_patch_j = sp2
331          end_patch_j   = ep2
332          start_patch_k = sp3
333          end_patch_k   = ep3
334    
335          allocate(real_array(sm1:em1,sm2:em2,sm3:em3))
336          if (my_proc_id /= IO_NODE) then
337             allocate(real_domain(1,1,1))
338             domain_start(1) = 1
339             domain_start(2) = 1
340             domain_start(3) = 1
341             domain_end(1) = 1
342             domain_end(2) = 1
343             domain_end(3) = 1
344          end if
345          call scatter_whole_field_r(real_array, &
346                                    sm1, em1, sm2, em2, sm3, em3, &
347                                    sp1, ep1, sp2, ep2, sp3, ep3, &
348                                    real_domain, &
349                                    domain_start(1), domain_end(1), &
350                                    domain_start(2), domain_end(2), &
351                                    domain_start(3), domain_end(3))
352          deallocate(real_domain)
354       else
355   
356          real_array => real_domain
358       end if
360    end subroutine read_next_field
362    subroutine read_global_attrs(title, start_date, grid_type, dyn_opt,                                &
363                                 west_east_dim, south_north_dim, bottom_top_dim,                       &
364                                 we_patch_s, we_patch_e, we_patch_s_stag, we_patch_e_stag,             &
365                                 sn_patch_s, sn_patch_e, sn_patch_s_stag, sn_patch_e_stag,             &
366                                 map_proj, mminlu, num_land_cat, is_water, is_lake, is_ice, is_urban,  &
367                                 isoilwater, grid_id, parent_id, i_parent_start, j_parent_start,       &
368                                 i_parent_end, j_parent_end, dx, dy, cen_lat, moad_cen_lat, cen_lon,   &
369                                 stand_lon, truelat1, truelat2, pole_lat, pole_lon, parent_grid_ratio, &
370                                 corner_lats, corner_lons, sr_x, sr_y)
372       implicit none
373   
374       ! Arguments
375       integer, intent(out) :: dyn_opt, west_east_dim, south_north_dim, bottom_top_dim, map_proj,   &
376                  is_water, is_lake, we_patch_s, we_patch_e, we_patch_s_stag, we_patch_e_stag,      &
377                  sn_patch_s, sn_patch_e, sn_patch_s_stag, sn_patch_e_stag,                         &
378                  is_ice, is_urban, isoilwater, grid_id, parent_id, i_parent_start, j_parent_start, &
379                  i_parent_end, j_parent_end, parent_grid_ratio, sr_x, sr_y, num_land_cat
380       real, intent(out) :: dx, dy, cen_lat, moad_cen_lat, cen_lon, stand_lon, truelat1, truelat2,  &
381                  pole_lat, pole_lon
382       real, dimension(16), intent(out) :: corner_lats, corner_lons
383       character (len=128), intent(out) :: title, start_date, grid_type, mminlu
384   
385       ! Local variables
386       integer :: istatus, i
387       real :: wps_version
388       character (len=128) :: cunits, cdesc, cstagger
389       type (q_data) :: qd
390   
391       if (my_proc_id == IO_NODE .or. do_tiled_input) then
392   
393 #ifdef IO_BINARY
394          if (io_form_input == BINARY) then
395             istatus = 0
396             do while (istatus == 0) 
397                cunits = ' '
398                cdesc = ' '
399                cstagger = ' '
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)
407          
408                      qd%units = cunits
409                      qd%description = cdesc
410                      qd%stagger = cstagger
411                      call q_insert(unit_desc, qd)
412         
413                   end if
414        
415                end if
416             end do
417          end if
418 #endif
419      
420          call ext_get_dom_ti_char          ('TITLE', title)
421          if (index(title,'GEOGRID V3.2.1') /= 0) then
422             wps_version = 3.21
423          else if (index(title,'GEOGRID V3.2') /= 0) then
424             wps_version = 3.2
425          else if (index(title,'GEOGRID V3.1.1') /= 0) then
426             wps_version = 3.11
427          else if (index(title,'GEOGRID V3.1') /= 0) then
428             wps_version = 3.1
429          else if (index(title,'GEOGRID V3.0.1') /= 0) then
430             wps_version = 3.01
431          else
432             wps_version = 3.0
433          end if
434          call mprintf(.true.,DEBUG,'Reading static data from WPS version %f', f1=wps_version)
435          call ext_get_dom_ti_char          ('SIMULATION_START_DATE', start_date)
436          call ext_get_dom_ti_integer_scalar('WEST-EAST_GRID_DIMENSION', west_east_dim)
437          call ext_get_dom_ti_integer_scalar('SOUTH-NORTH_GRID_DIMENSION', south_north_dim)
438          call ext_get_dom_ti_integer_scalar('BOTTOM-TOP_GRID_DIMENSION', bottom_top_dim)
439          call ext_get_dom_ti_integer_scalar('WEST-EAST_PATCH_START_UNSTAG', we_patch_s)
440          call ext_get_dom_ti_integer_scalar('WEST-EAST_PATCH_END_UNSTAG', we_patch_e)
441          call ext_get_dom_ti_integer_scalar('WEST-EAST_PATCH_START_STAG', we_patch_s_stag)
442          call ext_get_dom_ti_integer_scalar('WEST-EAST_PATCH_END_STAG', we_patch_e_stag)
443          call ext_get_dom_ti_integer_scalar('SOUTH-NORTH_PATCH_START_UNSTAG', sn_patch_s)
444          call ext_get_dom_ti_integer_scalar('SOUTH-NORTH_PATCH_END_UNSTAG', sn_patch_e)
445          call ext_get_dom_ti_integer_scalar('SOUTH-NORTH_PATCH_START_STAG', sn_patch_s_stag)
446          call ext_get_dom_ti_integer_scalar('SOUTH-NORTH_PATCH_END_STAG', sn_patch_e_stag)
447          call ext_get_dom_ti_char          ('GRIDTYPE', grid_type)
448          call ext_get_dom_ti_real_scalar   ('DX', dx)
449          call ext_get_dom_ti_real_scalar   ('DY', dy)
450          call ext_get_dom_ti_integer_scalar('DYN_OPT', dyn_opt)
451          call ext_get_dom_ti_real_scalar   ('CEN_LAT', cen_lat)
452          call ext_get_dom_ti_real_scalar   ('CEN_LON', cen_lon)
453          call ext_get_dom_ti_real_scalar   ('TRUELAT1', truelat1)
454          call ext_get_dom_ti_real_scalar   ('TRUELAT2', truelat2)
455          call ext_get_dom_ti_real_scalar   ('MOAD_CEN_LAT', moad_cen_lat)
456          call ext_get_dom_ti_real_scalar   ('STAND_LON', stand_lon)
457          call ext_get_dom_ti_real_scalar   ('POLE_LAT', pole_lat)
458          call ext_get_dom_ti_real_scalar   ('POLE_LON', pole_lon)
459          call ext_get_dom_ti_real_vector   ('corner_lats', corner_lats, 16)
460          call ext_get_dom_ti_real_vector   ('corner_lons', corner_lons, 16)
461          call ext_get_dom_ti_integer_scalar('MAP_PROJ', map_proj)
462          call ext_get_dom_ti_char          ('MMINLU', mminlu)
463          if ( wps_version >= 3.01 ) then
464             call ext_get_dom_ti_integer_scalar('NUM_LAND_CAT', num_land_cat)
465          else
466             num_land_cat = 24
467          end if
468          call ext_get_dom_ti_integer_scalar('ISWATER', is_water)
469          if ( wps_version >= 3.01 ) then
470             call ext_get_dom_ti_integer_scalar('ISLAKE', is_lake)
471          else
472             is_lake = -1
473          end if
474          call ext_get_dom_ti_integer_scalar('ISICE', is_ice)
475          call ext_get_dom_ti_integer_scalar('ISURBAN', is_urban)
476          call ext_get_dom_ti_integer_scalar('ISOILWATER', isoilwater)
477          call ext_get_dom_ti_integer_scalar('grid_id', grid_id)
478          call ext_get_dom_ti_integer_scalar('parent_id', parent_id)
479          call ext_get_dom_ti_integer_scalar('i_parent_start', i_parent_start)
480          call ext_get_dom_ti_integer_scalar('j_parent_start', j_parent_start)
481          call ext_get_dom_ti_integer_scalar('i_parent_end', i_parent_end)
482          call ext_get_dom_ti_integer_scalar('j_parent_end', j_parent_end)
483          call ext_get_dom_ti_integer_scalar('parent_grid_ratio', parent_grid_ratio)
484          call ext_get_dom_ti_integer_scalar('sr_x', sr_x)
485          call ext_get_dom_ti_integer_scalar('sr_y', sr_y)
486    
487       end if
489   
490       if (nprocs > 1 .and. .not. do_tiled_input) then
491   
492          call parallel_bcast_char(title, len(title))
493          call parallel_bcast_char(start_date, len(start_date))
494          call parallel_bcast_char(grid_type, len(grid_type))
495          call parallel_bcast_int(west_east_dim)
496          call parallel_bcast_int(south_north_dim)
497          call parallel_bcast_int(bottom_top_dim)
498          call parallel_bcast_int(we_patch_s)
499          call parallel_bcast_int(we_patch_e)
500          call parallel_bcast_int(we_patch_s_stag)
501          call parallel_bcast_int(we_patch_e_stag)
502          call parallel_bcast_int(sn_patch_s)
503          call parallel_bcast_int(sn_patch_e)
504          call parallel_bcast_int(sn_patch_s_stag)
505          call parallel_bcast_int(sn_patch_e_stag)
506          call parallel_bcast_int(sr_x)
507          call parallel_bcast_int(sr_y)
509          ! Must figure out patch dimensions from info in parallel module
510 !         we_patch_s      = my_minx
511 !         we_patch_s_stag = my_minx
512 !         we_patch_e      = my_maxx - 1
513 !         sn_patch_s      = my_miny
514 !         sn_patch_s_stag = my_miny
515 !         sn_patch_e      = my_maxy - 1
517 !         if (trim(grid_type) == 'C') then
518 !            if (my_x /= nproc_x - 1) then
519 !               we_patch_e_stag = we_patch_e + 1
520 !            end if
521 !            if (my_y /= nproc_y - 1) then
522 !               sn_patch_e_stag = sn_patch_e + 1
523 !            end if
524 !         else if (trim(grid_type) == 'E') then
525 !            we_patch_e = we_patch_e + 1
526 !            sn_patch_e = sn_patch_e + 1
527 !            we_patch_e_stag = we_patch_e
528 !            sn_patch_e_stag = sn_patch_e
529 !         end if
531          call parallel_bcast_real(dx)
532          call parallel_bcast_real(dy)
533          call parallel_bcast_int(dyn_opt)
534          call parallel_bcast_real(cen_lat)
535          call parallel_bcast_real(cen_lon)
536          call parallel_bcast_real(truelat1)
537          call parallel_bcast_real(truelat2)
538          call parallel_bcast_real(pole_lat)
539          call parallel_bcast_real(pole_lon)
540          call parallel_bcast_real(moad_cen_lat)
541          call parallel_bcast_real(stand_lon)
542          do i=1,16
543             call parallel_bcast_real(corner_lats(i))
544             call parallel_bcast_real(corner_lons(i))
545          end do
546          call parallel_bcast_int(map_proj)
547          call parallel_bcast_char(mminlu, len(mminlu))
548          call parallel_bcast_int(is_water)
549          call parallel_bcast_int(is_lake)
550          call parallel_bcast_int(is_ice)
551          call parallel_bcast_int(is_urban)
552          call parallel_bcast_int(isoilwater)
553          call parallel_bcast_int(grid_id)
554          call parallel_bcast_int(parent_id)
555          call parallel_bcast_int(i_parent_start)
556          call parallel_bcast_int(i_parent_end)
557          call parallel_bcast_int(j_parent_start)
558          call parallel_bcast_int(j_parent_end)
559          call parallel_bcast_int(parent_grid_ratio)
560       end if
561   
562       internal_gridtype = grid_type
564    end subroutine read_global_attrs
567    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
568    ! Name: ext_get_dom_ti_integer
569    !
570    ! Purpose: Read a domain time-independent integer attribute from input.
571    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
572    subroutine ext_get_dom_ti_integer_scalar(var_name, var_value)
574       implicit none
576       ! Arguments
577       integer, intent(out) :: var_value
578       character (len=*), intent(in) :: var_name
580       ! Local variables
581       integer :: istatus, outcount
583 #ifdef IO_BINARY
584       if (io_form_input == BINARY) then
585          call ext_int_get_dom_ti_integer(handle, trim(var_name), &
586                                          var_value, &
587                                          1, outcount, istatus)
588       end if
589 #endif
590 #ifdef IO_NETCDF
591       if (io_form_input == NETCDF) then
592          call ext_ncd_get_dom_ti_integer(handle, trim(var_name), &
593                                          var_value, &
594                                          1, outcount, istatus)
595       end if
596 #endif
597 #ifdef IO_GRIB1
598       if (io_form_input == GRIB1) then
599          call ext_gr1_get_dom_ti_integer(handle, trim(var_name), &
600                                          var_value, &
601                                          1, outcount, istatus)
602       end if
603 #endif
605       call mprintf((istatus /= 0),ERROR,'Error while reading domain time-independent attribute.')
607    end subroutine ext_get_dom_ti_integer_scalar
610    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
611    ! Name: ext_get_dom_ti_integer
612    !
613    ! Purpose: Read a domain time-independent integer attribute from input.
614    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
615    subroutine ext_get_dom_ti_integer_vector(var_name, var_value, n)
617       implicit none
619       ! Arguments
620       integer, intent(in) :: n
621       integer, dimension(n), intent(out) :: var_value
622       character (len=*), intent(in) :: var_name
624       ! Local variables
625       integer :: istatus, outcount
627 #ifdef IO_BINARY
628       if (io_form_input == BINARY) then
629          call ext_int_get_dom_ti_integer(handle, trim(var_name), &
630                                          var_value, &
631                                          n, outcount, istatus)
632       end if
633 #endif
634 #ifdef IO_NETCDF
635       if (io_form_input == NETCDF) then
636          call ext_ncd_get_dom_ti_integer(handle, trim(var_name), &
637                                          var_value, &
638                                          n, outcount, istatus)
639       end if
640 #endif
641 #ifdef IO_GRIB1
642       if (io_form_input == GRIB1) then
643          call ext_gr1_get_dom_ti_integer(handle, trim(var_name), &
644                                          var_value, &
645                                          n, outcount, istatus)
646       end if
647 #endif
649       call mprintf((istatus /= 0),ERROR,'Error while reading domain time-independent attribute.')
651    end subroutine ext_get_dom_ti_integer_vector
654    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
655    ! Name: ext_get_dom_ti_real
656    !
657    ! Purpose: Read a domain time-independent real attribute from input.
658    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
659    subroutine ext_get_dom_ti_real_scalar(var_name, var_value)
661       implicit none
663       ! Arguments
664       real, intent(out) :: var_value
665       character (len=*), intent(in) :: var_name
667       ! Local variables
668       integer :: istatus, outcount
670 #ifdef IO_BINARY
671       if (io_form_input == BINARY) then
672          call ext_int_get_dom_ti_real(handle, trim(var_name), &
673                                          var_value, &
674                                          1, outcount, istatus)
675       end if
676 #endif
677 #ifdef IO_NETCDF
678       if (io_form_input == NETCDF) then
679          call ext_ncd_get_dom_ti_real(handle, trim(var_name), &
680                                          var_value, &
681                                          1, outcount, istatus)
682       end if
683 #endif
684 #ifdef IO_GRIB1
685       if (io_form_input == GRIB1) then
686          call ext_gr1_get_dom_ti_real(handle, trim(var_name), &
687                                          var_value, &
688                                          1, outcount, istatus)
689       end if
690 #endif
692       call mprintf((istatus /= 0),ERROR,'Error while reading domain time-independent attribute.')
694    end subroutine ext_get_dom_ti_real_scalar
697    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
698    ! Name: ext_get_dom_ti_real
699    !
700    ! Purpose: Read a domain time-independent real attribute from input.
701    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
702    subroutine ext_get_dom_ti_real_vector(var_name, var_value, n)
704       implicit none
706       ! Arguments
707       integer, intent(in) :: n
708       real, dimension(n), intent(out) :: var_value
709       character (len=*), intent(in) :: var_name
711       ! Local variables
712       integer :: istatus, outcount
714 #ifdef IO_BINARY
715       if (io_form_input == BINARY) then
716          call ext_int_get_dom_ti_real(handle, trim(var_name), &
717                                          var_value, &
718                                          n, outcount, istatus)
719       end if
720 #endif
721 #ifdef IO_NETCDF
722       if (io_form_input == NETCDF) then
723          call ext_ncd_get_dom_ti_real(handle, trim(var_name), &
724                                          var_value, &
725                                          n, outcount, istatus)
726       end if
727 #endif
728 #ifdef IO_GRIB1
729       if (io_form_input == GRIB1) then
730          call ext_gr1_get_dom_ti_real(handle, trim(var_name), &
731                                          var_value, &
732                                          n, outcount, istatus)
733       end if
734 #endif
736       call mprintf((istatus /= 0),ERROR,'Error while reading domain time-independent attribute.')
738    end subroutine ext_get_dom_ti_real_vector
741    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
742    ! Name: ext_get_dom_ti_char
743    !
744    ! Purpose: Read a domain time-independent character attribute from input.
745    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
746    subroutine ext_get_dom_ti_char(var_name, var_value)
748       implicit none
750       ! Arguments
751       character (len=*), intent(in) :: var_name
752       character (len=128), intent(out) :: var_value
754       ! Local variables
755       integer :: istatus
757 #ifdef IO_BINARY
758       if (io_form_input == BINARY) then
759          call ext_int_get_dom_ti_char(handle, trim(var_name), &
760                                          var_value, &
761                                          istatus)
762       end if
763 #endif
764 #ifdef IO_NETCDF
765       if (io_form_input == NETCDF) then
766          call ext_ncd_get_dom_ti_char(handle, trim(var_name), &
767                                          var_value, &
768                                          istatus)
769       end if
770 #endif
771 #ifdef IO_GRIB1
772       if (io_form_input == GRIB1) then
773          call ext_gr1_get_dom_ti_char(handle, trim(var_name), &
774                                          var_value, &
775                                          istatus)
776       end if
777 #endif
779       call mprintf((istatus /= 0),ERROR,'Error in reading domain time-independent attribute')
781    end subroutine ext_get_dom_ti_char
784    subroutine input_close()
786       implicit none
787   
788       ! Local variables
789       integer :: istatus
790   
791       istatus = 0
792       if (my_proc_id == IO_NODE .or. do_tiled_input) then
793 #ifdef IO_BINARY
794          if (io_form_input == BINARY) then
795             call ext_int_ioclose(handle, istatus)
796             call ext_int_ioexit(istatus)
797          end if
798 #endif
799 #ifdef IO_NETCDF
800          if (io_form_input == NETCDF) then
801             call ext_ncd_ioclose(handle, istatus)
802             call ext_ncd_ioexit(istatus)
803          end if
804 #endif
805 #ifdef IO_GRIB1
806          if (io_form_input == GRIB1) then
807             call ext_gr1_ioclose(handle, istatus)
808             call ext_gr1_ioexit(istatus)
809          end if
810 #endif
811       end if
813       call q_destroy(unit_desc)
815    end subroutine input_close
817 end module input_module