1 module interp_option_module
5 use misc_definitions_module
9 integer, parameter :: BUFSIZE=128
11 integer :: num_entries
12 integer, pointer, dimension(:) :: output_stagger
13 real, pointer, dimension(:) :: masked, fill_missing, missing_value, &
14 interp_mask_val, interp_land_mask_val, interp_water_mask_val
15 logical, pointer, dimension(:) :: output_this_field, is_u_field, is_v_field, is_derived_field, is_mandatory
16 character (len=128), pointer, dimension(:) :: fieldname, interp_method, v_interp_method, &
17 interp_mask, interp_land_mask, interp_water_mask, &
18 flag_in_output, output_name, from_input, z_dim_name, level_template
19 character (len=1), pointer, dimension(:) :: interp_mask_relational, interp_land_mask_relational, interp_water_mask_relational
20 type (list), pointer, dimension(:) :: fill_lev_list
21 type (list) :: flag_in_output_list
25 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
26 ! Name: read_interp_table
29 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
30 subroutine read_interp_table()
33 integer :: i, p1, p2, idx, eos, ispace, funit, istatus, nparams, s1, s2
34 logical :: is_used, have_specification
35 character (len=128) :: lev_string, fill_string, flag_string, flag_val
36 character (len=BUFSIZE) :: buffer
39 inquire(unit=funit, opened=is_used)
40 if (.not. is_used) exit
46 open(funit, file=trim(opt_metgrid_tbl_path)//'METGRID.TBL', form='formatted', status='old', err=1001)
48 do while (istatus == 0)
49 read(funit, '(a)', iostat=istatus) buffer
50 if (istatus == 0) then
53 ! Is this line a comment?
54 if (buffer(1:1) == '#') then
56 ! Are we beginning a new field specification?
57 else if (index(buffer,'=====') /= 0) then
58 if (nparams > 0) num_entries = num_entries + 1
62 eos = index(buffer,'#')
63 if (eos /= 0) buffer(eos:BUFSIZE) = ' '
65 ! Does this line contain at least one parameter specification?
66 if (index(buffer,'=') /= 0) then
76 ! Allocate one extra array element to act as the default
77 ! BUG: Maybe this will not be necessary if we move to a module with query routines for
78 ! parsing the METGRID.TBL
79 num_entries = num_entries + 1
81 allocate(fieldname(num_entries))
82 allocate(interp_method(num_entries))
83 allocate(v_interp_method(num_entries))
84 allocate(masked(num_entries))
85 allocate(fill_missing(num_entries))
86 allocate(missing_value(num_entries))
87 allocate(fill_lev_list(num_entries))
88 allocate(interp_mask(num_entries))
89 allocate(interp_land_mask(num_entries))
90 allocate(interp_water_mask(num_entries))
91 allocate(interp_mask_val(num_entries))
92 allocate(interp_land_mask_val(num_entries))
93 allocate(interp_water_mask_val(num_entries))
94 allocate(interp_mask_relational(num_entries))
95 allocate(interp_land_mask_relational(num_entries))
96 allocate(interp_water_mask_relational(num_entries))
97 allocate(level_template(num_entries))
98 allocate(flag_in_output(num_entries))
99 allocate(output_name(num_entries))
100 allocate(from_input(num_entries))
101 allocate(z_dim_name(num_entries))
102 allocate(output_stagger(num_entries))
103 allocate(output_this_field(num_entries))
104 allocate(is_u_field(num_entries))
105 allocate(is_v_field(num_entries))
106 allocate(is_derived_field(num_entries))
107 allocate(is_mandatory(num_entries))
114 flag_in_output(i) = ' '
117 z_dim_name(i) = 'num_metgrid_levels'
118 interp_method(i) = 'nearest_neighbor'
119 v_interp_method(i) = 'linear_log_p'
120 masked(i) = NOT_MASKED
121 fill_missing(i) = NAN
122 missing_value(i) = NAN
123 call list_init(fill_lev_list(i))
125 interp_land_mask(i) = ' '
126 interp_water_mask(i) = ' '
127 interp_mask_val(i) = NAN
128 interp_land_mask_val(i) = NAN
129 interp_water_mask_val(i) = NAN
130 interp_mask_relational(i) = ' '
131 interp_land_mask_relational(i) = ' '
132 interp_water_mask_relational(i) = ' '
133 level_template(i) = ' '
134 if (gridtype == 'C') then
135 output_stagger(i) = M
136 else if (gridtype == 'E') then
137 output_stagger(i) = HH
139 output_this_field(i) = .true.
140 is_u_field(i) = .false.
141 is_v_field(i) = .false.
142 is_derived_field(i) = .false.
143 is_mandatory(i) = .false.
145 call list_init(flag_in_output_list)
151 do while (istatus == 0)
153 read(funit, '(a)', iostat=istatus) buffer
154 if (istatus == 0) then
157 ! Is this line a comment?
158 if (buffer(1:1) == '#') then
161 ! Are we beginning a new field specification?
162 else if (index(buffer,'=====') /= 0) then !{
163 if (nparams > 0) i = i + 1
167 ! Check whether the current line is a comment
168 if (buffer(1:1) /= '#') then
169 have_specification = .true.
171 have_specification = .false.
174 ! If only part of the line is a comment, just turn the comment into spaces
175 eos = index(buffer,'#')
176 if (eos /= 0) buffer(eos:BUFSIZE) = ' '
178 do while (have_specification) !{
180 ! If this line has no semicolon, it may contain a single specification,
181 ! so we set have_specification = .false. to prevent the line from being
182 ! processed again and "pretend" that the last character was a semicolon
183 eos = index(buffer,';')
185 have_specification = .false.
189 idx = index(buffer(1:eos-1),'=')
191 if (idx /= 0) then !{
192 nparams = nparams + 1
194 if (index('name',trim(buffer(1:idx-1))) /= 0 .and. &
195 len_trim('name') == len_trim(buffer(1:idx-1))) then
197 do while ((ispace < eos) .and. (buffer(ispace:ispace) /= ' '))
201 fieldname(i)(1:ispace-idx) = buffer(idx+1:ispace-1)
203 else if (index('from_input',trim(buffer(1:idx-1))) /= 0 .and. &
204 len_trim('from_input') == len_trim(buffer(1:idx-1))) then
206 do while ((ispace < eos) .and. (buffer(ispace:ispace) /= ' '))
210 from_input(i)(1:ispace-idx) = buffer(idx+1:ispace-1)
212 else if (index('z_dim_name',trim(buffer(1:idx-1))) /= 0 .and. &
213 len_trim('z_dim_name') == len_trim(buffer(1:idx-1))) then
215 do while ((ispace < eos) .and. (buffer(ispace:ispace) /= ' '))
219 z_dim_name(i)(1:ispace-idx) = buffer(idx+1:ispace-1)
221 else if (index('output_stagger',trim(buffer(1:idx-1))) /= 0 .and. &
222 len_trim('output_stagger') == len_trim(buffer(1:idx-1))) then
223 if (index('M',trim(buffer(idx+1:eos-1))) /= 0) then
224 output_stagger(i) = M
225 else if (index('U',trim(buffer(idx+1:eos-1))) /= 0) then
226 output_stagger(i) = U
227 else if (index('V',trim(buffer(idx+1:eos-1))) /= 0) then
228 output_stagger(i) = V
229 else if (index('HH',trim(buffer(idx+1:eos-1))) /= 0) then
230 output_stagger(i) = HH
231 else if (index('VV',trim(buffer(idx+1:eos-1))) /= 0) then
232 output_stagger(i) = VV
235 else if (index('output',trim(buffer(1:idx-1))) /= 0 .and. &
236 len_trim('output') == len_trim(buffer(1:idx-1))) then
237 if (index('yes',trim(buffer(idx+1:eos-1))) /= 0) then
238 output_this_field(i) = .true.
239 else if (index('no',trim(buffer(idx+1:eos-1))) /= 0) then
240 output_this_field(i) = .false.
243 else if (index('is_u_field',trim(buffer(1:idx-1))) /= 0 .and. &
244 len_trim('is_u_field') == len_trim(buffer(1:idx-1))) then
245 if (index('yes',trim(buffer(idx+1:eos-1))) /= 0) then
246 is_u_field(i) = .true.
247 else if (index('no',trim(buffer(idx+1:eos-1))) /= 0) then
248 is_u_field(i) = .false.
251 else if (index('is_v_field',trim(buffer(1:idx-1))) /= 0 .and. &
252 len_trim('is_v_field') == len_trim(buffer(1:idx-1))) then
253 if (index('yes',trim(buffer(idx+1:eos-1))) /= 0) then
254 is_v_field(i) = .true.
255 else if (index('no',trim(buffer(idx+1:eos-1))) /= 0) then
256 is_v_field(i) = .false.
259 else if (index('derived',trim(buffer(1:idx-1))) /= 0 .and. &
260 len_trim('derived') == len_trim(buffer(1:idx-1))) then
261 if (index('yes',trim(buffer(idx+1:eos-1))) /= 0) then
262 is_derived_field(i) = .true.
263 else if (index('no',trim(buffer(idx+1:eos-1))) /= 0) then
264 is_derived_field(i) = .false.
267 else if (index('mandatory',trim(buffer(1:idx-1))) /= 0 .and. &
268 len_trim('mandatory') == len_trim(buffer(1:idx-1))) then
269 if (index('yes',trim(buffer(idx+1:eos-1))) /= 0) then
270 is_mandatory(i) = .true.
271 else if (index('no',trim(buffer(idx+1:eos-1))) /= 0) then
272 is_mandatory(i) = .false.
275 else if (index('interp_option',trim(buffer(1:idx-1))) /= 0 .and. &
276 len_trim('interp_option') == len_trim(buffer(1:idx-1))) then
278 do while ((ispace < eos) .and. (buffer(ispace:ispace) /= ' '))
281 interp_method(i) = ' '
282 interp_method(i)(1:ispace-idx) = buffer(idx+1:ispace-1)
284 else if (index('vertical_interp_option',trim(buffer(1:idx-1))) /= 0 .and. &
285 len_trim('vertical_interp_option') == len_trim(buffer(1:idx-1))) then
287 do while ((ispace < eos) .and. (buffer(ispace:ispace) /= ' '))
290 v_interp_method(i) = ' '
291 v_interp_method(i)(1:ispace-idx) = buffer(idx+1:ispace-1)
293 else if (index('level_template',trim(buffer(1:idx-1))) /= 0 .and. &
294 len_trim('level_template') == len_trim(buffer(1:idx-1))) then
296 do while ((ispace < eos) .and. (buffer(ispace:ispace) /= ' '))
299 level_template(i)(1:ispace-idx) = buffer(idx+1:ispace-1)
301 else if (index('interp_mask',trim(buffer(1:idx-1))) /= 0 .and. &
302 len_trim('interp_mask') == len_trim(buffer(1:idx-1))) then
304 do while ((ispace < eos) .and. (buffer(ispace:ispace) /= ' '))
307 p1 = index(buffer(idx+1:ispace-1),'(')
308 p2 = index(buffer(idx+1:ispace-1),')')
309 s1 = index(buffer(idx+1:ispace-1),'<')
310 s2 = index(buffer(idx+1:ispace-1),'>')
311 if (p1 == 0 .or. p2 == 0) then
312 call mprintf(.true.,WARN, &
313 'Problem in specifying interp_mask flag. Setting masked flag to 0.')
315 interp_mask(i)(1:ispace-idx) = buffer(idx+1:ispace-1)
316 interp_mask_val(i) = 0
318 ! Parenthesis found; additionally, there may be a relational symbol
319 if ((s1 /= 0) .OR. (s2 /= 0)) then
321 interp_mask_relational(i) = buffer(idx+s1:idx+s1)
322 else if (s2 > 0) then
323 interp_mask_relational(i) = buffer(idx+s2:idx+s2)
326 interp_mask(i)(1:p1) = buffer(idx+1:idx+p1-1)
327 read(buffer(idx+p1+2:idx+p2-1),*,err=1000) interp_mask_val(i)
329 ! No relational symbol
331 interp_mask(i)(1:p1) = buffer(idx+1:idx+p1-1)
332 read(buffer(idx+p1+1:idx+p2-1),*,err=1000) interp_mask_val(i)
336 else if (index('interp_land_mask',trim(buffer(1:idx-1))) /= 0 .and. &
337 len_trim('interp_land_mask') == len_trim(buffer(1:idx-1))) then
339 do while ((ispace < eos) .and. (buffer(ispace:ispace) /= ' '))
342 p1 = index(buffer(idx+1:ispace-1),'(')
343 p2 = index(buffer(idx+1:ispace-1),')')
344 s1 = index(buffer(idx+1:ispace-1),'<')
345 s2 = index(buffer(idx+1:ispace-1),'>')
346 if (p1 == 0 .or. p2 == 0) then
347 call mprintf(.true.,WARN, &
348 'Problem in specifying interp_land_mask flag. Setting masked flag to 0.')
349 interp_land_mask(i) = ' '
350 interp_land_mask(i)(1:ispace-idx) = buffer(idx+1:ispace-1)
351 interp_land_mask_val(i) = 0
353 ! Parenthesis found; additionally, there may be a relational symbol
354 if ((s1 /= 0) .OR. (s2 /= 0)) then
356 interp_land_mask_relational(i) = buffer(idx+s1:idx+s1)
357 else if (s2 > 0) then
358 interp_land_mask_relational(i) = buffer(idx+s2:idx+s2)
360 interp_land_mask(i) = ' '
361 interp_land_mask(i)(1:p1) = buffer(idx+1:idx+p1-1)
362 read(buffer(idx+p1+2:idx+p2-1),*,err=1000) interp_land_mask_val(i)
364 ! No relational symbol
365 interp_land_mask(i) = ' '
366 interp_land_mask(i)(1:p1) = buffer(idx+1:idx+p1-1)
367 read(buffer(idx+p1+1:idx+p2-1),*,err=1000) interp_land_mask_val(i)
371 else if (index('interp_water_mask',trim(buffer(1:idx-1))) /= 0 .and. &
372 len_trim('interp_water_mask') == len_trim(buffer(1:idx-1))) then
374 do while ((ispace < eos) .and. (buffer(ispace:ispace) /= ' '))
377 p1 = index(buffer(idx+1:ispace-1),'(')
378 p2 = index(buffer(idx+1:ispace-1),')')
379 s1 = index(buffer(idx+1:ispace-1),'<')
380 s2 = index(buffer(idx+1:ispace-1),'>')
381 if (p1 == 0 .or. p2 == 0) then
382 call mprintf(.true.,WARN, &
383 'Problem in specifying interp_water_mask flag. Setting masked flag to 0.')
384 interp_water_mask(i) = ' '
385 interp_water_mask(i)(1:ispace-idx) = buffer(idx+1:ispace-1)
386 interp_water_mask_val(i) = 0
388 ! Parenthesis found; additionally, there may be a relational symbol
389 if ((s1 /= 0) .OR. (s2 /= 0)) then
391 interp_water_mask_relational(i) = buffer(idx+s1:idx+s1)
392 else if (s2 > 0) then
393 interp_water_mask_relational(i) = buffer(idx+s2:idx+s2)
395 interp_water_mask(i) = ' '
396 interp_water_mask(i)(1:p1) = buffer(idx+1:idx+p1-1)
397 read(buffer(idx+p1+2:idx+p2-1),*,err=1000) interp_water_mask_val(i)
399 ! No relational symbol
400 interp_water_mask(i) = ' '
401 interp_water_mask(i)(1:p1) = buffer(idx+1:idx+p1-1)
402 read(buffer(idx+p1+1:idx+p2-1),*,err=1000) interp_water_mask_val(i)
406 else if (index('masked',trim(buffer(1:idx-1))) /= 0 .and. &
407 len_trim('masked') == len_trim(buffer(1:idx-1))) then
408 if (index('water',trim(buffer(idx+1:eos-1))) /= 0) then
409 masked(i) = MASKED_WATER
410 else if (index('land',trim(buffer(idx+1:eos-1))) /= 0) then
411 masked(i) = MASKED_LAND
412 else if (index('both',trim(buffer(idx+1:eos-1))) /= 0) then
413 masked(i) = MASKED_BOTH
416 else if (index('flag_in_output',trim(buffer(1:idx-1))) /= 0 .and. &
417 len_trim('flag_in_output') == len_trim(buffer(1:idx-1))) then
419 flag_string(1:eos-idx-1) = buffer(idx+1:eos-1)
420 if (list_search(flag_in_output_list, ckey=flag_string, cvalue=flag_val)) then
421 call mprintf(.true.,WARN, 'In METGRID.TBL, %s is given as a flag more than once.', &
423 flag_in_output(i)(1:eos-idx-1) = buffer(idx+1:eos-1)
425 flag_in_output(i)(1:eos-idx-1) = buffer(idx+1:eos-1)
426 write(flag_val,'(i1)') 1
427 call list_insert(flag_in_output_list, ckey=flag_string, cvalue=flag_val)
430 else if (index('output_name',trim(buffer(1:idx-1))) /= 0 .and. &
431 len_trim('output_name') == len_trim(buffer(1:idx-1))) then
433 do while ((ispace < eos) .and. (buffer(ispace:ispace) /= ' '))
437 output_name(i)(1:ispace-idx) = buffer(idx+1:ispace-1)
439 else if (index('fill_missing',trim(buffer(1:idx-1))) /= 0 .and. &
440 len_trim('fill_missing') == len_trim(buffer(1:idx-1))) then
441 read(buffer(idx+1:eos-1),*) fill_missing(i)
443 else if (index('missing_value',trim(buffer(1:idx-1))) /= 0 .and. &
444 len_trim('missing_value') == len_trim(buffer(1:idx-1))) then
445 read(buffer(idx+1:eos-1),*) missing_value(i)
447 else if (index('fill_lev',trim(buffer(1:idx-1))) /= 0 .and. &
448 len_trim('fill_lev') == len_trim(buffer(1:idx-1))) then
450 do while ((ispace < eos) .and. (buffer(ispace:ispace) /= ' '))
454 fill_string(1:ispace-idx-1) = buffer(idx+1:ispace-1)
455 ispace = index(fill_string,':')
456 if (ispace /= 0) then
457 write(lev_string,'(a)') fill_string(1:ispace-1)
459 write(lev_string,'(a)') 'all'
461 write(fill_string,'(a)') trim(fill_string(ispace+1:128))
462 fill_string(128-ispace:128) = ' '
463 call list_insert(fill_lev_list(i), ckey=lev_string, cvalue=fill_string)
466 call mprintf(.true.,WARN, 'In METGRID.TBL, unrecognized option %s in entry %i.', s1=buffer(1:idx-1), i1=idx)
469 end if !} index(buffer(1:eos-1),'=') /= 0
471 ! BUG: If buffer has non-whitespace characters but no =, then maybe a wrong specification?
473 buffer = buffer(eos+1:BUFSIZE)
474 end do ! while eos /= 0 }
476 end if !} index(buffer, '=====') /= 0
481 call check_table_specs()
487 1000 call mprintf(.true.,ERROR,'The mask value of the interp_mask specification must '// &
488 'be a real value, enclosed in parentheses immediately after the field name.')
490 1001 call mprintf(.true.,ERROR,'Could not open file METGRID.TBL')
491 1002 call mprintf(.true.,ERROR,'Symbol expected < >. Check METGRID.TBL for missing symbol or erroreous entry')
493 end subroutine read_interp_table
496 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
497 ! Name: check_table_specs
499 ! Pupose: Perform basic consistency and sanity checks on the METGRID.TBL
500 ! entries supplied by the user.
501 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
502 subroutine check_table_specs()
511 ! For C grid, U field must be on U staggering, and V field must be on
512 ! V staggering; for E grid, U and V must be on VV staggering.
513 if (gridtype == 'C') then
514 if (is_u_field(i) .and. output_stagger(i) /= U) then
515 call mprintf(.true.,ERROR,'In entry %i of METGRID.TBL, the wind U-component field '// &
516 'must be interpolated to the U staggered grid points.',i1=i)
517 else if (is_v_field(i) .and. output_stagger(i) /= V) then
518 call mprintf(.true.,ERROR,'In entry %i of METGRID.TBL, the wind V-component field '// &
519 'must be interpolated to the V staggered grid points.',i1=i)
522 if (output_stagger(i) == VV) then
523 call mprintf(.true.,ERROR,'In entry %i of METGRID.TBL, VV is not a valid output staggering for ARW.',i1=i)
524 else if (output_stagger(i) == HH) then
525 call mprintf(.true.,ERROR,'In entry %i of METGRID.TBL, HH is not a valid output staggering for ARW.',i1=i)
528 if (masked(i) /= NOT_MASKED .and. output_stagger(i) /= M) then
529 call mprintf(.true.,ERROR,'In entry %i of METGRID.TBL, staggered output field '// &
530 'cannot use the ''masked'' option.',i1=i)
533 else if (gridtype == 'E') then
534 if (is_u_field(i) .and. output_stagger(i) /= VV) then
535 call mprintf(.true.,ERROR,'In entry %i of METGRID.TBL, the wind U-component field '// &
536 'must be interpolated to the V staggered grid points.',i1=i)
537 else if (is_v_field(i) .and. output_stagger(i) /= VV) then
538 call mprintf(.true.,ERROR,'In entry %i of METGRID.TBL, the wind V-component field '// &
539 'must be interpolated to the V staggered grid points.',i1=i)
542 if (output_stagger(i) == M) then
543 call mprintf(.true.,ERROR,'In entry %i of METGRID.TBL, M is not a valid output staggering for NMM.',i1=i)
544 else if (output_stagger(i) == U) then
545 call mprintf(.true.,ERROR,'In entry %i of METGRID.TBL, U is not a valid output staggering for NMM.',i1=i)
546 else if (output_stagger(i) == V) then
547 call mprintf(.true.,ERROR,'In entry %i of METGRID.TBL, V is not a valid output staggering for NMM.',i1=i)
550 if (masked(i) /= NOT_MASKED .and. output_stagger(i) /= HH) then
551 call mprintf(.true.,ERROR,'In entry %i of METGRID.TBL, staggered output field '// &
552 'cannot use the ''masked'' option.',i1=i)
558 end subroutine check_table_specs
561 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
562 ! Name: get_z_dim_name
565 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
566 subroutine get_z_dim_name(fldname, zdim_name)
571 character (len=*), intent(in) :: fldname
572 character (len=32), intent(out) :: zdim_name
577 zdim_name = z_dim_name(num_entries)(1:32)
579 if (trim(fldname) == trim(fieldname(i))) then
580 zdim_name = z_dim_name(i)(1:32)
585 end subroutine get_z_dim_name
588 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
589 ! Name: get_gcell_threshold
592 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
593 subroutine get_gcell_threshold(interp_opt, threshold, istatus)
598 integer, intent(out) :: istatus
599 real, intent(out) :: threshold
600 character (len=128), intent(in) :: interp_opt
608 i = index(interp_opt,'average_gcell')
611 ! Check for a threshold
612 p1 = index(interp_opt(i:128),'(')
613 p2 = index(interp_opt(i:128),')')
614 if (p1 /= 0 .and. p2 /= 0) then
615 read(interp_opt(p1+1:p2-1),*,err=1000) threshold
617 call mprintf(.true.,WARN, 'Problem in specifying threshold for average_gcell interp option. Setting threshold to 1.0')
625 1000 call mprintf(.true.,ERROR, &
626 'Threshold option to average_gcell interpolator must be a real number, '// &
627 'enclosed in parentheses immediately after keyword "average_gcell"')
629 end subroutine get_gcell_threshold
632 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
633 ! Name: get_constant_fill_lev
636 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
637 subroutine get_constant_fill_lev(fill_opt, fill_const, istatus)
642 integer, intent(out) :: istatus
643 real, intent(out) :: fill_const
644 character (len=128), intent(in) :: fill_opt
652 i = index(fill_opt,'const')
655 ! Check for a threshold
656 p1 = index(fill_opt(i:128),'(')
657 p2 = index(fill_opt(i:128),')')
658 if (p1 /= 0 .and. p2 /= 0) then
659 read(fill_opt(p1+1:p2-1),*,err=1000) fill_const
661 call mprintf(.true.,WARN, 'Problem in specifying fill_lev constant. Setting fill_const to %f', f1=NAN)
669 1000 call mprintf(.true.,ERROR, &
670 'Constant option to fill_lev must be a real number, enclosed in parentheses '// &
671 'immediately after keyword "const"')
673 end subroutine get_constant_fill_lev
676 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
677 ! Name: get_fill_src_level
680 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
681 subroutine get_fill_src_level(fill_opt, fill_src, fill_src_level)
686 integer, intent(out) :: fill_src_level
687 character (len=128), intent(in) :: fill_opt
688 character (len=128), intent(out) :: fill_src
693 ! Check for a level in parentheses
694 p1 = index(fill_opt,'(')
695 p2 = index(fill_opt,')')
696 if (p1 /= 0 .and. p2 /= 0) then
697 read(fill_opt(p1+1:p2-1),*,err=1000) fill_src_level
699 write(fill_src,'(a)') fill_opt(1:p1-1)
707 1000 call mprintf(.true.,ERROR, &
708 'For fill_lev specification, level in source field must be an integer, '// &
709 'enclosed in parentheses immediately after the fieldname')
711 end subroutine get_fill_src_level
714 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
715 ! Name: interp_option_destroy
718 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
719 subroutine interp_option_destroy()
726 deallocate(fieldname)
727 deallocate(from_input)
728 deallocate(z_dim_name)
729 deallocate(interp_method)
730 deallocate(v_interp_method)
732 deallocate(fill_missing)
733 deallocate(missing_value)
735 call list_destroy(fill_lev_list(i))
737 deallocate(fill_lev_list)
738 deallocate(interp_mask)
739 deallocate(interp_land_mask)
740 deallocate(interp_water_mask)
741 deallocate(interp_mask_val)
742 deallocate(interp_land_mask_val)
743 deallocate(interp_water_mask_val)
744 deallocate(interp_mask_relational)
745 deallocate(interp_land_mask_relational)
746 deallocate(interp_water_mask_relational)
747 deallocate(level_template)
748 deallocate(flag_in_output)
749 deallocate(output_name)
750 deallocate(output_stagger)
751 deallocate(output_this_field)
752 deallocate(is_u_field)
753 deallocate(is_v_field)
754 deallocate(is_derived_field)
755 deallocate(is_mandatory)
756 call list_destroy(flag_in_output_list)
758 end subroutine interp_option_destroy
760 end module interp_option_module