36 module time_interp_external2_mod
45 use platform_mod,
only : r8_kind, r4_kind
55 use time_interp_mod,
only :
time_interp, time_interp_init
58 use platform_mod,
only: r8_kind, fms_path_len, fms_file_len
59 use horiz_interp_mod,
only :
horiz_interp, horiz_interp_type
60 use fms2_io_mod,
only : valid_t, fmsnetcdfdomainfile_t,
open_file, get_unlimited_dimension_name, &
61 variable_att_exists, fmsnetcdffile_t, &
62 variable_exists, get_valid, get_variable_num_dimensions,
read_data, &
63 is_valid,
close_file, get_dimension_size, get_variable_dimension_names, &
64 get_variable_size, get_time_calendar, get_variable_missing, get_variable_units
70 #include<file_version.h>
72 integer,
parameter,
public :: NO_REGION=0, inside_region=1, outside_region=2
73 integer,
parameter,
private :: modulo_year= 0001
74 integer,
parameter,
private :: LINEAR_TIME_INTERP = 1
75 integer,
parameter,
public :: SUCCESS = 0, err_field_not_found = 1
76 real(r8_kind),
parameter,
private :: DEFAULT_MISSING_VALUE = -1e20_r8_kind
77 integer,
private :: max_fields = 100, max_files= 40
78 integer,
private :: num_fields = 0, num_files=0
80 integer,
private :: num_io_buffers = 2
81 logical,
private :: module_initialized = .false.
82 logical,
private :: debug_this_module = .false.
87 public get_external_fileobj
97 type(fmsnetcdffile_t),
pointer :: fileobj=>null()
98 character(len=128) :: name, units
99 integer :: siz(4), ndim
100 character(len=32) :: axisname(4)
103 type(
time_type),
dimension(:),
pointer :: start_time =>null(), end_time =>null()
104 type(
time_type),
dimension(:),
pointer :: period =>null()
105 logical :: modulo_time
106 real(r8_kind),
dimension(:,:,:,:),
pointer :: domain_data =>null()
107 logical,
dimension(:,:,:,:),
pointer :: mask =>null()
108 integer,
dimension(:),
pointer :: ibuf =>null()
109 real(r8_kind),
dimension(:,:,:,:),
pointer :: src_data =>null()
110 type(valid_t) :: valid
112 logical :: domain_present
113 real(r8_kind) :: slope, intercept
114 integer :: isc,iec,jsc,jec
115 type(
time_type) :: modulo_time_beg, modulo_time_end
116 logical :: have_modulo_times, correct_leap_year_inconsistency
117 integer :: region_type
118 integer :: is_region, ie_region, js_region, je_region
119 integer :: is_src, ie_src, js_src, je_src
121 integer :: numwindows
122 logical,
dimension(:,:),
pointer :: need_compute=>null()
123 real(r8_kind) :: missing
129 character(len=FMS_FILE_LEN) :: filename =
''
130 type(FmsNetcdfFile_t),
pointer :: fileobj => null()
146 module procedure time_interp_external_0d_r4
147 module procedure time_interp_external_2d_r4
148 module procedure time_interp_external_3d_r4
149 module procedure time_interp_external_0d_r8
150 module procedure time_interp_external_2d_r8
151 module procedure time_interp_external_3d_r8
155 module procedure time_interp_external_bridge_0d_r4
156 module procedure time_interp_external_bridge_2d_r4
157 module procedure time_interp_external_bridge_3d_r4
158 module procedure time_interp_external_bridge_0d_r8
159 module procedure time_interp_external_bridge_2d_r8
160 module procedure time_interp_external_bridge_3d_r8
168 type(ext_fieldtype),
save,
private,
pointer :: loaded_fields(:) => null()
169 type(filetype),
save,
private,
pointer :: opened_files(:) => null()
171 real(r8_kind),
private,
parameter :: time_interp_missing=-1e99_r8_kind
183 integer :: io_status, logunit, ierr
185 namelist /time_interp_external_nml/ num_io_buffers, debug_this_module, &
186 max_fields, max_files
190 if(module_initialized)
return
194 call write_version_number(
"TIME_INTERP_EXTERNAL_MOD", version)
196 read (input_nml_file, time_interp_external_nml, iostat=io_status)
199 write(logunit,time_interp_external_nml)
203 module_initialized = .true.
205 call time_interp_init()
271 verbose,axis_names, axis_sizes,override,correct_leap_year_inconsistency,&
272 permit_calendar_conversion,use_comp_domain,ierr, nwindows, ignore_axis_atts, ongrid )
274 character(len=*),
intent(in) :: file,fieldname
275 logical,
intent(in),
optional :: verbose
276 character(len=*),
intent(in),
optional :: desired_units
277 type(
domain2d),
intent(in),
optional :: domain
278 integer,
intent(inout),
optional :: axis_sizes(4)
279 character(len=*),
intent(out),
optional :: axis_names(4)
280 logical,
intent(in),
optional :: override, correct_leap_year_inconsistency,&
281 permit_calendar_conversion,use_comp_domain
282 integer,
intent(out),
optional :: ierr
283 integer,
intent(in),
optional :: nwindows
284 logical,
optional :: ignore_axis_atts
285 logical,
optional :: ongrid
287 logical :: ongrid_local
291 real(r8_kind) :: slope, intercept
292 integer :: ndim,ntime,i,j
293 integer :: iscomp,iecomp,jscomp,jecomp,isglobal,ieglobal,jsglobal,jeglobal
294 integer :: isdata,iedata,jsdata,jedata, dxsize, dysize,dxsize_max,dysize_max
295 logical :: verb, transpose_xy,use_comp_domain1
296 real(r8_kind),
dimension(:),
allocatable :: tstamp, tstart, tend, tavg
297 character(len=1) :: cart
298 character(len=1),
dimension(4) :: cart_dir
299 character(len=128) :: units, fld_units
300 character(len=128) :: msg, calendar_type, timebeg, timeend
301 character(len=128) :: timename, timeunits
302 character(len=128),
allocatable :: axisname(:)
303 integer,
allocatable :: axislen(:)
304 integer :: siz(4), siz_in(4), gxsize, gysize,gxsize_max, gysize_max
306 integer :: yr, mon, day, hr, minu, sec
307 integer :: len, nfile, nfields_orig, nbuf, nx,ny
308 integer :: numwindows
309 logical :: ignore_axatts
310 logical :: have_modulo_time
311 type(fmsnetcdffile_t),
pointer :: fileobj=>null()
312 integer,
dimension(:),
allocatable :: pes
314 if (.not. module_initialized)
call mpp_error(fatal,
'Must call time_interp_external_init first')
315 if(
present(ierr)) ierr = success
316 ignore_axatts=.false.
317 cart_dir(1)=
'X';cart_dir(2)=
'Y';cart_dir(3)=
'Z';cart_dir(4)=
'T'
318 if(
present(ignore_axis_atts)) ignore_axatts = ignore_axis_atts
319 use_comp_domain1 = .false.
320 if(
PRESENT(use_comp_domain)) use_comp_domain1 = use_comp_domain
322 if (
PRESENT(verbose)) verb=verbose
323 if (debug_this_module) verb = .true.
325 if(
present(nwindows)) numwindows = nwindows
328 if (
PRESENT(desired_units))
then
329 units = desired_units
330 call mpp_error(fatal,
'==> Unit conversion via time_interp_external &
331 &has been temporarily deprecated. Previous versions of&
332 &this module used udunits_mod to perform unit conversion.&
333 & Udunits_mod is in the process of being replaced since &
334 &there were portability issues associated with this code.&
335 & Please remove the desired_units argument from calls to &
340 if(trim(opened_files(i)%filename) == trim(file))
then
346 num_files = num_files + 1
347 if(num_files > max_files)
then
353 call mpp_error(fatal,
"time_interp_external: num_files is greater than max_files, "// &
354 "increase time_interp_external_nml max_files")
356 opened_files(num_files)%filename = trim(file)
357 allocate(opened_files(num_files)%fileobj)
358 fileobj => opened_files(num_files)%fileobj
360 if(.not.
open_file(opened_files(num_files)%fileobj, trim(file),
'read')) &
361 call mpp_error(fatal,
'time_interp_external_mod: Error in opening file '//trim(file))
363 fileobj => opened_files(nfile)%fileobj
366 call get_unlimited_dimension_name(fileobj, timename)
367 call get_dimension_size(fileobj, timename, ntime)
370 write(msg,
'(a15,a,a58)')
'external field ',trim(fieldname),&
371 ' does not have an associated record dimension (REQUIRED) '
376 call get_time_calendar(fileobj, timename, calendar_type)
379 call get_variable_units(fileobj, timename, timeunits)
385 call mpp_get_current_pelist(pes)
386 allocate(tstamp(ntime),tstart(ntime),tend(ntime),tavg(ntime))
389 if (mpp_root_pe() .eq.
mpp_pe())
call read_data(fileobj, timename, tstamp)
390 call mpp_broadcast(tstamp,
size(tstamp), mpp_root_pe(), pelist=pes)
393 transpose_xy = .false.
394 isdata=1; iedata=1; jsdata=1; jedata=1
398 if (
PRESENT(domain))
then
400 nx = iecomp-iscomp+1; ny = jecomp-jscomp+1
401 call mpp_get_data_domain(domain,isdata,iedata,jsdata,jedata,dxsize,dxsize_max,dysize,dysize_max)
402 call mpp_get_global_domain(domain,isglobal,ieglobal,jsglobal,jeglobal,gxsize,gxsize_max,gysize,gysize_max)
403 ongrid_local = .false.
404 if (
present(ongrid)) ongrid_local = ongrid
407 if (ongrid_local)
then
413 elseif(use_comp_domain1)
then
414 call mpp_error(fatal,
"init_external_field:"//&
415 " use_comp_domain=true but domain is not present")
419 nfields_orig = num_fields
421 if (.not. variable_exists(fileobj, fieldname) )
then
422 if (
present(ierr))
then
423 ierr = err_field_not_found
426 call mpp_error(fatal,
'external field "'//trim(fieldname)//
'" not found in file "'//trim(file)//
'"')
433 if(variable_att_exists(fileobj, fieldname,
'time_avg_info'))
then
434 if(variable_exists(fileobj,
'average_T1'))
call read_data(fileobj,
'average_T1', tstart)
435 if(variable_exists(fileobj,
'average_T2'))
call read_data(fileobj,
'average_T2', tend)
436 if(variable_exists(fileobj,
'average_DT'))
call read_data(fileobj,
'average_DT', tavg)
439 num_fields = num_fields + 1
440 if(num_fields > max_fields)
then
446 call mpp_error(fatal,
"time_interp_external: num_fields is greater than max_fields, "// &
447 "increase time_interp_external_nml max_fields")
451 call get_variable_units(fileobj, fieldname, fld_units)
454 loaded_fields(num_fields)%fileobj => fileobj
455 loaded_fields(num_fields)%name = trim(fieldname)
456 loaded_fields(num_fields)%units = trim(fld_units)
457 loaded_fields(num_fields)%isc = 1
458 loaded_fields(num_fields)%iec = 1
459 loaded_fields(num_fields)%jsc = 1
460 loaded_fields(num_fields)%jec = 1
461 loaded_fields(num_fields)%region_type = no_region
462 loaded_fields(num_fields)%is_region = 0
463 loaded_fields(num_fields)%ie_region = -1
464 loaded_fields(num_fields)%js_region = 0
465 loaded_fields(num_fields)%je_region = -1
466 if (
PRESENT(domain))
then
467 loaded_fields(num_fields)%domain_present = .true.
468 loaded_fields(num_fields)%domain = domain
469 loaded_fields(num_fields)%isc=iscomp;loaded_fields(num_fields)%iec = iecomp
470 loaded_fields(num_fields)%jsc=jscomp;loaded_fields(num_fields)%jec = jecomp
472 loaded_fields(num_fields)%domain_present = .false.
475 loaded_fields(num_fields)%valid = get_valid(fileobj, fieldname)
476 ndim = get_variable_num_dimensions(fileobj, fieldname)
478 'invalid array rank <=4d fields supported')
480 loaded_fields(num_fields)%siz = 1
481 loaded_fields(num_fields)%ndim = ndim
482 loaded_fields(num_fields)%tdim = 4
484 loaded_fields(num_fields)%missing = get_variable_missing(fileobj, fieldname)
486 allocate(axisname(ndim), axislen(ndim))
488 call get_variable_dimension_names(fileobj, fieldname, axisname)
489 call get_variable_size(fileobj, fieldname, axislen)
490 do j=1,loaded_fields(num_fields)%ndim
493 if (cart ==
'N' .and. .not. ignore_axatts)
then
494 write(msg,
'(a,"/",a)') trim(file),trim(fieldname)
495 call mpp_error(fatal,
'file/field '//trim(msg)// &
496 ' couldnt recognize axis atts in time_interp_external')
497 else if (cart ==
'N' .and. ignore_axatts)
then
502 if (j.eq.2) transpose_xy = .true.
503 if (.not.
PRESENT(domain) .and. .not.
PRESENT(override))
then
508 loaded_fields(num_fields)%isc=iscomp;loaded_fields(num_fields)%iec=iecomp
509 elseif (
PRESENT(override))
then
511 if (
PRESENT(axis_sizes)) axis_sizes(1) = len
513 loaded_fields(num_fields)%axisname(1) = axisname(j)
514 if(use_comp_domain1)
then
515 loaded_fields(num_fields)%siz(1) = nx
517 loaded_fields(num_fields)%siz(1) = dxsize
519 if (len /= gxsize)
then
520 write(msg,
'(a,"/",a)') trim(file),trim(fieldname)
521 call mpp_error(fatal,
'time_interp_ext, file/field '//trim(msg)//
' x dim doesnt match model')
524 loaded_fields(num_fields)%axisname(2) = axisname(j)
525 if (.not.
PRESENT(domain) .and. .not.
PRESENT(override))
then
530 loaded_fields(num_fields)%jsc=jscomp;loaded_fields(num_fields)%jec=jecomp
531 elseif (
PRESENT(override))
then
533 if (
PRESENT(axis_sizes)) axis_sizes(2) = len
535 if(use_comp_domain1)
then
536 loaded_fields(num_fields)%siz(2) = ny
538 loaded_fields(num_fields)%siz(2) = dysize
540 if (len /= gysize)
then
541 write(msg,
'(a,"/",a)') trim(file),trim(fieldname)
542 call mpp_error(fatal,
'time_interp_ext, file/field '//trim(msg)//
' y dim doesnt match model')
545 loaded_fields(num_fields)%axisname(3) = axisname(j)
546 loaded_fields(num_fields)%siz(3) = len
548 loaded_fields(num_fields)%axisname(4) = axisname(j)
549 loaded_fields(num_fields)%siz(4) = ntime
550 loaded_fields(num_fields)%tdim = j
553 siz = loaded_fields(num_fields)%siz
554 if(
PRESENT(axis_names)) axis_names = loaded_fields(num_fields)%axisname
555 if (
PRESENT(axis_sizes) .and. .not.
PRESENT(override))
then
556 axis_sizes = loaded_fields(num_fields)%siz
559 if (verb)
write(outunit,
'(a,4i6)')
'field x,y,z,t local size= ',siz
560 if (verb)
write(outunit,*)
'field contains data in units = ',trim(loaded_fields(num_fields)%units)
561 if (transpose_xy)
call mpp_error(fatal,
'axis ordering not supported')
562 if (num_io_buffers .le. 1)
call mpp_error(fatal,
'time_interp_ext:num_io_buffers should be at least 2')
563 nbuf = min(num_io_buffers,siz(4))
565 loaded_fields(num_fields)%numwindows = numwindows
566 allocate(loaded_fields(num_fields)%need_compute(nbuf, numwindows))
567 loaded_fields(num_fields)%need_compute = .true.
569 allocate(loaded_fields(num_fields)%domain_data(isdata:iedata,jsdata:jedata,siz(3),nbuf),&
570 loaded_fields(num_fields)%mask(isdata:iedata,jsdata:jedata,siz(3),nbuf) )
571 loaded_fields(num_fields)%mask = .false.
572 loaded_fields(num_fields)%domain_data = 0.0_r8_kind
573 slope=1.0_r8_kind;intercept=0.0_r8_kind
579 loaded_fields(num_fields)%slope = slope
580 loaded_fields(num_fields)%intercept = intercept
581 allocate(loaded_fields(num_fields)%ibuf(nbuf))
582 loaded_fields(num_fields)%ibuf = -1
583 loaded_fields(num_fields)%nbuf = 0
584 if(
PRESENT(override))
then
585 loaded_fields(num_fields)%is_src = 1
586 loaded_fields(num_fields)%ie_src = gxsize
587 loaded_fields(num_fields)%js_src = 1
588 loaded_fields(num_fields)%je_src = gysize
589 allocate(loaded_fields(num_fields)%src_data(gxsize,gysize,siz(3),nbuf))
591 loaded_fields(num_fields)%is_src = isdata
592 loaded_fields(num_fields)%ie_src = iedata
593 loaded_fields(num_fields)%js_src = jsdata
594 loaded_fields(num_fields)%je_src = jedata
595 allocate(loaded_fields(num_fields)%src_data(isdata:iedata,jsdata:jedata,siz(3),nbuf))
598 allocate(loaded_fields(num_fields)%time(ntime))
599 allocate(loaded_fields(num_fields)%period(ntime))
600 allocate(loaded_fields(num_fields)%start_time(ntime))
601 allocate(loaded_fields(num_fields)%end_time(ntime))
604 loaded_fields(num_fields)%time(j) =
get_cal_time(tstamp(j),trim(timeunits),trim(calendar_type), &
605 & permit_calendar_conversion)
606 loaded_fields(num_fields)%start_time(j) =
get_cal_time(tstart(j),trim(timeunits),trim(calendar_type), &
607 & permit_calendar_conversion)
608 loaded_fields(num_fields)%end_time(j) =
get_cal_time( tend(j),trim(timeunits),trim(calendar_type), &
609 & permit_calendar_conversion)
612 if (loaded_fields(num_fields)%modulo_time)
then
613 call set_time_modulo(loaded_fields(num_fields)%Time)
614 call set_time_modulo(loaded_fields(num_fields)%start_time)
615 call set_time_modulo(loaded_fields(num_fields)%end_time)
618 if(
present(correct_leap_year_inconsistency))
then
619 loaded_fields(num_fields)%correct_leap_year_inconsistency = correct_leap_year_inconsistency
621 loaded_fields(num_fields)%correct_leap_year_inconsistency = .false.
624 if(have_modulo_time)
then
626 loaded_fields(num_fields)%modulo_time_beg =
set_time(timebeg)
627 loaded_fields(num_fields)%modulo_time_end =
set_time(timeend)
629 loaded_fields(num_fields)%modulo_time_beg =
set_date(timebeg)
630 loaded_fields(num_fields)%modulo_time_end =
set_date(timeend)
632 loaded_fields(num_fields)%have_modulo_times = .true.
634 loaded_fields(num_fields)%have_modulo_times = .false.
637 call mpp_error(note,
'time_interp_external_mod: file '//trim(file)//
' has only one time level')
640 loaded_fields(num_fields)%period(j) = loaded_fields(num_fields)%end_time(j) &
641 - loaded_fields(num_fields)%start_time(j)
642 if (loaded_fields(num_fields)%period(j) >
set_time(0,0))
then
643 call get_time(loaded_fields(num_fields)%period(j), sec, day)
644 sec = sec/2+mod(day,2)*43200
646 loaded_fields(num_fields)%time(j) = loaded_fields(num_fields)%start_time(j)+&
649 if (j > 1 .and. j < ntime)
then
650 tdiff = loaded_fields(num_fields)%time(j+1) - loaded_fields(num_fields)%time(j-1)
652 sec = sec/2+mod(day,2)*43200
654 loaded_fields(num_fields)%period(j) =
set_time(sec,day)
655 sec = sec/2+mod(day,2)*43200
657 loaded_fields(num_fields)%start_time(j) = loaded_fields(num_fields)%time(j) -
set_time(sec,day)
658 loaded_fields(num_fields)%end_time(j) = loaded_fields(num_fields)%time(j) +
set_time(sec,day)
659 elseif ( j == 1)
then
660 tdiff = loaded_fields(num_fields)%time(2) - loaded_fields(num_fields)%time(1)
662 loaded_fields(num_fields)%period(j) =
set_time(sec,day)
663 sec = sec/2+mod(day,2)*43200
665 loaded_fields(num_fields)%start_time(j) = loaded_fields(num_fields)%time(j) -
set_time(sec,day)
666 loaded_fields(num_fields)%end_time(j) = loaded_fields(num_fields)%time(j) +
set_time(sec,day)
668 tdiff = loaded_fields(num_fields)%time(ntime) - loaded_fields(num_fields)%time(ntime-1)
670 loaded_fields(num_fields)%period(j) =
set_time(sec,day)
671 sec = sec/2+mod(day,2)*43200
673 loaded_fields(num_fields)%start_time(j) = loaded_fields(num_fields)%time(j) -
set_time(sec,day)
674 loaded_fields(num_fields)%end_time(j) = loaded_fields(num_fields)%time(j) +
set_time(sec,day)
681 if (loaded_fields(num_fields)%time(j) >= loaded_fields(num_fields)%time(j+1))
then
682 write(msg,
'(A,i20)')
"times not monotonically increasing. Filename: " &
683 //trim(file)//
" field: "//trim(fieldname)//
" timeslice: ", j
688 loaded_fields(num_fields)%modulo_time =
get_axis_modulo(fileobj, timename)
691 if (loaded_fields(num_fields)%modulo_time)
write(outunit,*)
'data are being treated as modulo in time'
693 write(outunit,*)
'time index, ', j
694 call get_date(loaded_fields(num_fields)%start_time(j),yr,mon,day,hr,minu,sec)
695 write(outunit,
'(a,i4,a,i2,a,i2,1x,i2,a,i2,a,i2)') &
696 'start time: yyyy/mm/dd hh:mm:ss= ',yr,
'/',mon,
'/',day,hr,
':',minu,
':',sec
697 call get_date(loaded_fields(num_fields)%time(j),yr,mon,day,hr,minu,sec)
698 write(outunit,
'(a,i4,a,i2,a,i2,1x,i2,a,i2,a,i2)') &
699 'mid time: yyyy/mm/dd hh:mm:ss= ',yr,
'/',mon,
'/',day,hr,
':',minu,
':',sec
700 call get_date(loaded_fields(num_fields)%end_time(j),yr,mon,day,hr,minu,sec)
701 write(outunit,
'(a,i4,a,i2,a,i2,1x,i2,a,i2,a,i2)') &
702 'end time: yyyy/mm/dd hh:mm:ss= ',yr,
'/',mon,
'/',day,hr,
':',minu,
':',sec
706 deallocate(axisname, axislen)
707 deallocate(tstamp, tstart, tend, tavg)
714 function get_external_fileobj(filename, fileobj)
715 character(len=*),
intent(in) :: filename
716 type(fmsnetcdffile_t),
intent(out) :: fileobj
717 logical :: get_external_fileobj
720 get_external_fileobj = .false.
722 if(trim(opened_files(i)%filename) == trim(filename))
then
723 fileobj = opened_files(i)%fileobj
724 get_external_fileobj = .true.
729 end function get_external_fileobj
732 subroutine set_time_modulo(Time)
734 type(
time_type),
intent(inout),
dimension(:) :: time
737 integer :: yr, mon, dy, hr, minu, sec
739 ntime =
size(time(:))
742 call get_date(time(n), yr, mon, dy, hr, minu, sec)
744 time(n) =
set_date(yr, mon, dy, hr, minu, sec)
748 end subroutine set_time_modulo
753 subroutine load_record(field, rec, interp, is_in, ie_in, js_in, je_in, window_id_in)
755 integer ,
intent(in) :: rec
756 type(horiz_interp_type),
intent(in),
optional :: interp
757 integer,
intent(in),
optional :: is_in, ie_in, js_in, je_in
758 integer,
intent(in),
optional :: window_id_in
762 integer :: isw,iew,jsw,jew
763 integer :: is_region, ie_region, js_region, je_region, i, j
764 integer :: start(4), nread(4)
765 logical :: need_compute
767 real(r8_kind) :: mask_in(size(field%src_data,1),size(field%src_data,2),size(field%src_data,3))
768 real(r8_kind),
allocatable :: mask_out(:,:,:)
769 real(r4_kind),
allocatable :: hi_tmp_data(:,:,:,:)
770 real(r4_kind),
allocatable :: hi_tmp_msk_out(:,:,:)
771 real(r4_kind),
allocatable :: hi_tmp_src_data(:,:,:,:)
774 if(
PRESENT(window_id_in) ) window_id = window_id_in
775 need_compute = .true.
782 need_compute = .false.
785 field%nbuf = field%nbuf + 1
786 if(field%nbuf >
size(field%domain_data,4).or.field%nbuf <= 0) field%nbuf = 1
789 field%need_compute(ib,:) = .true.
791 if (debug_this_module)
write(outunit,*)
'reading record without domain for field ',trim(field%name)
793 start(1) = field%is_src; nread(1) = field%ie_src - field%is_src + 1
794 start(2) = field%js_src; nread(2) = field%je_src - field%js_src + 1
795 start(3) = 1; nread(3) =
size(field%src_data,3)
796 start(field%tdim) = rec; nread(field%tdim) = 1
797 call read_data(field%fileobj,field%name,field%src_data(:,:,:,ib:ib),corner=start,edge_lengths=nread)
800 isw=field%isc;iew=field%iec
801 jsw=field%jsc;jew=field%jec
803 if( field%numwindows > 1)
then
804 if( .NOT.
PRESENT(is_in) .OR. .NOT.
PRESENT(ie_in) .OR. .NOT.
PRESENT(js_in) .OR. .NOT.
PRESENT(je_in) )
then
806 &
'time_interp_external(load_record): is_in, ie_in, js_in, je_in must be present when numwindows>1')
808 isw = isw + is_in - 1
809 iew = isw + ie_in - is_in
810 jsw = jsw + js_in - 1
811 jew = jsw + je_in - js_in
816 need_compute = field%need_compute(ib, window_id)
817 if(need_compute)
then
818 if(
PRESENT(interp))
then
819 is_region = field%is_region; ie_region = field%ie_region
820 js_region = field%js_region; je_region = field%je_region
821 mask_in = 0.0_r8_kind
822 where (is_valid(field%src_data(:,:,:,ib), field%valid)) mask_in = 1.0_r8_kind
823 if ( field%region_type .NE. no_region )
then
824 if( any(mask_in == 0.0_r8_kind) )
then
825 call mpp_error(fatal,
"time_interp_external: mask_in should be all 1 when region_type is not NO_REGION")
827 if( field%region_type == outside_region)
then
828 do j = js_region, je_region
829 do i = is_region, ie_region
830 mask_in(i,j,:) = 0.0_r8_kind
834 do j = 1,
size(mask_in,2)
835 do i = 1,
size(mask_in,1)
836 if( j<js_region .OR. j>je_region .OR. i<is_region .OR. i>ie_region ) mask_in(i,j,:) = 0.0_r8_kind
844 if (interp%horizInterpReals4_type%is_allocated)
then
846 allocate(hi_tmp_msk_out(isw:iew,jsw:jew,
SIZE(field%src_data,3)))
847 allocate(hi_tmp_data(lbound(field%domain_data,1):ubound(field%domain_data,1), &
848 lbound(field%domain_data,2):ubound(field%domain_data,2), &
849 lbound(field%domain_data,3):ubound(field%domain_data,3), &
850 lbound(field%domain_data,4):ubound(field%domain_data,4)))
851 allocate(hi_tmp_src_data(lbound(field%src_data,1):ubound(field%src_data,1), &
852 lbound(field%src_data,2):ubound(field%src_data,2), &
853 lbound(field%src_data,3):ubound(field%src_data,3), &
854 lbound(field%src_data,4):ubound(field%src_data,4)))
856 hi_tmp_data = real(field%domain_data, r4_kind)
857 hi_tmp_src_data = real(field%src_data, r4_kind)
859 call horiz_interp(interp, hi_tmp_src_data(:,:,:,ib), hi_tmp_data(isw:iew,jsw:jew,:,ib), &
860 mask_in=real(mask_in,r4_kind), mask_out=hi_tmp_msk_out)
862 field%domain_data = real(hi_tmp_data, r8_kind)
863 field%mask(isw:iew,jsw:jew,:,ib) = hi_tmp_msk_out(isw:iew,jsw:jew,:) > 0.0_r4_kind
865 if(
allocated(hi_tmp_data))
deallocate(hi_tmp_data)
866 if(
allocated(hi_tmp_msk_out))
deallocate(hi_tmp_msk_out)
867 if(
allocated(hi_tmp_src_data))
deallocate(hi_tmp_src_data)
869 allocate(mask_out(isw:iew,jsw:jew,
size(field%src_data,3)))
870 call horiz_interp(interp, field%src_data(:,:,:,ib),field%domain_data(isw:iew,jsw:jew,:,ib), &
873 field%mask(isw:iew,jsw:jew,:,ib) = mask_out(isw:iew,jsw:jew,:) > 0.0_r8_kind
878 if ( field%region_type .NE. no_region )
then
879 call mpp_error(fatal,
"time_interp_external: region_type should be NO_REGION when interp is not present")
881 field%domain_data(isw:iew,jsw:jew,:,ib) = field%src_data(isw:iew,jsw:jew,:,ib)
882 field%mask(isw:iew,jsw:jew,:,ib) = is_valid(field%domain_data(isw:iew,jsw:jew,:,ib),field%valid)
885 where(field%mask(isw:iew,jsw:jew,:,ib)) field%domain_data(isw:iew,jsw:jew,:,ib) = &
886 field%domain_data(isw:iew,jsw:jew,:,ib)*field%slope + field%intercept
887 field%need_compute(ib, window_id) = .false.
896 integer ,
intent(in) :: rec
899 integer :: start(4), nread(4)
907 field%nbuf = field%nbuf + 1
908 if(field%nbuf >
size(field%domain_data,4).or.field%nbuf <= 0) field%nbuf = 1
912 if (debug_this_module)
write(outunit,*)
'reading record without domain for field ',trim(field%name)
914 start(3) = 1; nread(3) =
size(field%src_data,3)
915 start(field%tdim) = rec; nread(field%tdim) = 1
916 call read_data(field%fileobj,field%name,field%src_data(:,:,:,ib),corner=start,edge_lengths=nread)
917 if ( field%region_type .NE. no_region )
then
918 call mpp_error(fatal,
"time_interp_external: region_type should be NO_REGION when field is scalar")
920 field%domain_data(1,1,:,ib) = field%src_data(1,1,:,ib)
921 field%mask(1,1,:,ib) = is_valid(field%domain_data(1,1,:,ib),field%valid)
923 where(field%mask(1,1,:,ib)) field%domain_data(1,1,:,ib) = &
924 field%domain_data(1,1,:,ib)*field%slope + field%intercept
931 integer,
intent(in) :: index
932 integer,
intent(in) :: is, ie, js, je
935 if( is == loaded_fields(index)%is_src .AND. ie == loaded_fields(index)%ie_src .AND. &
936 js == loaded_fields(index)%js_src .AND. ie == loaded_fields(index)%je_src )
return
938 if( .NOT.
ASSOCIATED(loaded_fields(index)%src_data) )
call mpp_error(fatal, &
939 "time_interp_external: field(index)%src_data is not associated")
940 nk =
size(loaded_fields(index)%src_data,3)
941 nbuf =
size(loaded_fields(index)%src_data,4)
942 deallocate(loaded_fields(index)%src_data)
943 allocate(loaded_fields(index)%src_data(is:ie,js:je,nk,nbuf))
944 loaded_fields(index)%is_src = is
945 loaded_fields(index)%ie_src = ie
946 loaded_fields(index)%js_src = js
947 loaded_fields(index)%je_src = je
952 subroutine set_override_region(index, region_type, is_region, ie_region, js_region, je_region)
953 integer,
intent(in) :: index, region_type
954 integer,
intent(in) :: is_region, ie_region, js_region, je_region
956 loaded_fields(index)%region_type = region_type
957 loaded_fields(index)%is_region = is_region
958 loaded_fields(index)%ie_region = ie_region
959 loaded_fields(index)%js_region = js_region
960 loaded_fields(index)%je_region = je_region
964 end subroutine set_override_region
968 integer,
intent(in) :: n
973 if (
associated(opened_files))
then
974 if (n <=
size(opened_files))
return
980 if(
Associated(ptr(i)%fileobj))
then
982 deallocate(ptr(i)%fileobj)
986 if (
associated(opened_files))
then
987 ptr(1:
size(opened_files)) = opened_files(:)
988 deallocate(opened_files)
996 integer,
intent(in) :: n
1001 if (
associated(loaded_fields))
then
1002 if (n <=
size(loaded_fields))
return
1007 ptr(i)%fileobj => null()
1012 ptr(i)%domain = null_domain2d
1013 if (
ASSOCIATED(ptr(i)%time))
DEALLOCATE(ptr(i)%time, stat=ier)
1014 if (
ASSOCIATED(ptr(i)%start_time))
DEALLOCATE(ptr(i)%start_time, stat=ier)
1015 if (
ASSOCIATED(ptr(i)%end_time))
DEALLOCATE(ptr(i)%end_time, stat=ier)
1016 if (
ASSOCIATED(ptr(i)%period))
DEALLOCATE(ptr(i)%period, stat=ier)
1017 ptr(i)%modulo_time=.false.
1018 if (
ASSOCIATED(ptr(i)%domain_data))
DEALLOCATE(ptr(i)%domain_data, stat=ier)
1019 if (
ASSOCIATED(ptr(i)%ibuf))
DEALLOCATE(ptr(i)%ibuf, stat=ier)
1020 if (
ASSOCIATED(ptr(i)%src_data))
DEALLOCATE(ptr(i)%src_data, stat=ier)
1022 ptr(i)%domain_present=.false.
1023 ptr(i)%slope=1.0_r8_kind
1024 ptr(i)%intercept=0.0_r8_kind
1025 ptr(i)%isc=-1;ptr(i)%iec=-1
1026 ptr(i)%jsc=-1;ptr(i)%jec=-1
1028 if (
associated(loaded_fields))
then
1029 ptr(1:
size(loaded_fields)) = loaded_fields(:)
1030 deallocate(loaded_fields)
1040 integer,
dimension(:) :: buf
1050 if (buf(i) == indx)
then
1066 if (index .lt. 1 .or. index .gt. num_fields) &
1067 call mpp_error(fatal,
'invalid index in call to get_external_field_size')
1083 if (index .lt. 1 .or. index .gt. num_fields) &
1084 call mpp_error(fatal,
'invalid index in call to get_external_field_size')
1092 integer ,
intent(in) :: index
1097 if (index < 1.or.index > num_fields) &
1098 call mpp_error(fatal,
'invalid index in call to get_time_axis')
1100 n = min(
size(time),
size(loaded_fields(index)%time))
1102 time(1:n) = loaded_fields(index)%time(1:n)
1115 deallocate(loaded_fields(i)%time,loaded_fields(i)%start_time,loaded_fields(i)%end_time,&
1116 loaded_fields(i)%period,loaded_fields(i)%domain_data,loaded_fields(i)%mask,loaded_fields(i)%ibuf)
1117 if (
ASSOCIATED(loaded_fields(i)%src_data))
deallocate(loaded_fields(i)%src_data)
1118 loaded_fields(i)%domain = null_domain2d
1119 loaded_fields(i)%nbuf = 0
1120 loaded_fields(i)%slope = 0.0_r8_kind
1121 loaded_fields(i)%intercept = 0.0_r8_kind
1127 deallocate(opened_files(i)%fileobj)
1130 deallocate(loaded_fields)
1131 deallocate(opened_files)
1135 module_initialized = .false.
1139 #include "time_interp_external2_r4.fh"
1140 #include "time_interp_external2_r8.fh"
1142 #include "time_interp_external2_bridge_r4.fh"
1143 #include "time_interp_external2_bridge_r8.fh"
1145 end module time_interp_external2_mod
subroutine, public get_axis_cart(fileobj, axisname, cart)
Returns X,Y,Z or T cartesian attribute.
logical function, public get_axis_modulo(fileobj, axisname)
Checks if 'modulo' variable exists for a given axis.
logical function, public get_axis_modulo_times(fileobj, axisname, tbeg, tend)
Close a netcdf or domain file opened with open_file or open_virtual_file.
Opens a given netcdf or domain file.
Read data from a defined field in a file.
integer function, public check_nml_error(IOSTAT, NML_NAME)
Checks the iostat argument that is returned after reading a namelist and determines if the error code...
subroutine, public write_version_number(version, tag, unit)
Prints to the log file (or a specified unit) the version id string and tag name.
Added for mixed precision support. Updates force time_manager math to be done with kind=8 reals _wrap...
Subroutine for performing the horizontal interpolation between two grids.
These routines retrieve the axis specifications associated with the compute domains....
These routines retrieve the axis specifications associated with the data domains. The domain is a der...
These routines retrieve the axis specifications associated with the global domains....
The domain2D type contains all the necessary information to define the global, compute and data domai...
integer function stdout()
This function returns the current standard fortran unit numbers for output.
integer function stdlog()
This function returns the current standard fortran unit numbers for log messages. Log messages,...
integer function mpp_npes()
Returns processor count for current pelist.
integer function mpp_pe()
Returns processor ID.
Perform parallel broadcasts.
subroutine, public time_interp_external_exit()
exit time_interp_external_mod. Close all open files and release storage
integer function, dimension(4), public get_external_field_size(index)
Returns size of field after call to init_external_field. Ordering is X/Y/Z/T. This call only makes se...
subroutine, public reset_src_data_region(index, is, ie, js, je)
Reallocates src_data for field from module level loaded_fields array.
subroutine load_record(field, rec, interp, is_in, ie_in, js_in, je_in, window_id_in)
load specified record from file Always loads in r8, casts down for horiz_interp if interp argument is...
subroutine, public time_interp_external_init()
Initialize the time_interp_external_mod module.
real(r8_kind) function, public get_external_field_missing(index)
return missing value
subroutine realloc_files(n)
reallocates array of fields, increasing its size
integer function, public init_external_field(file, fieldname, domain, desired_units, verbose, axis_names, axis_sizes, override, correct_leap_year_inconsistency, permit_calendar_conversion, use_comp_domain, ierr, nwindows, ignore_axis_atts, ongrid)
Initialize an external field. Buffer "num_io_buffers" (default=2) in memory to reduce memory allocati...
integer function, private find_buf_index(indx, buf)
simple linear search for given value in given list TODO should use better search if this list is bigg...
subroutine, public get_time_axis(index, time)
subroutine realloc_fields(n)
reallocates array of fields,increasing its size
subroutine load_record_0d(field, rec)
Given a initialized ext_fieldtype and record number, loads the given index into fieldsrc_data.
Provide data from external file interpolated to current model time. Data may be local to current proc...
Represents external fields.
Holds filename and file object.
Returns a weight and dates or indices for interpolating between two dates. The interface fraction_of_...
subroutine, public get_time(Time, seconds, days, ticks, err_msg)
Returns days and seconds ( < 86400 ) corresponding to a time. err_msg should be checked for any error...
subroutine, public get_date(time, year, month, day, hour, minute, second, tick, err_msg)
Gets the date for different calendar types. Given a time_interval, returns the corresponding date und...
integer function, public days_in_year(Time)
Returns the number of days in the calendar year corresponding to the date represented by time for the...
type(time_type) function, public increment_time(Time, seconds, days, ticks, err_msg, allow_neg_inc)
Increments a time by seconds and days.
integer function, public days_in_month(Time, err_msg)
Given a time, computes the corresponding date given the selected date time mapping algorithm.
integer function, public get_calendar_type()
Returns default calendar type for mapping from time to date.
Given an input date in year, month, days, etc., creates a time_type that represents this time interva...
Given some number of seconds and days, returns the corresponding time_type.
Type to represent amounts of time. Implemented as seconds and days to allow for larger intervals.