MODFLOW 6  version 6.7.0.dev2
USGS Modular Hydrologic Model
BinaryFileReader.f90
Go to the documentation of this file.
2 
3  use kindmodule, only: i4b, i8b, dp, lgp
4  use errorutilmodule, only: pstop
6 
8 
10  integer(I4B) :: pos, size
11  integer(I4B) :: kper, kstp
12  real(dp) :: pertim, totim
13  contains
14  procedure :: get_str
15  end type binaryfileheadertype
16 
17  type, abstract :: binaryfilereadertype
18  integer(I4B) :: inunit
19  integer(I4B) :: nrecords
20  logical(LGP) :: indexed
21  logical(LGP) :: endoffile
22  integer(I4B), allocatable :: record_sizes(:)
23  class(binaryfileheadertype), allocatable :: header
24  class(binaryfileheadertype), allocatable :: headernext
25  contains
26  procedure(read_header_if), deferred :: read_header
27  procedure(read_record_if), deferred :: read_record
28  procedure :: build_index
29  procedure :: peek_record
30  procedure :: rewind
31  end type binaryfilereadertype
32 
33  abstract interface
34  subroutine read_header_if(this, success, iout)
36  import i4b, lgp
37  class(binaryfilereadertype), intent(inout) :: this
38  logical(LGP), intent(out) :: success
39  integer(I4B), intent(in), optional :: iout
40  end subroutine read_header_if
41 
42  subroutine read_record_if(this, success, iout)
44  import i4b, lgp
45  class(binaryfilereadertype), intent(inout) :: this
46  logical(LGP), intent(out) :: success
47  integer(I4B), intent(in), optional :: iout
48  end subroutine read_record_if
49  end interface
50 contains
51 
52  !> @brief Get a string representation of the header.
53  function get_str(this) result(str)
54  class(binaryfileheadertype), intent(in) :: this
55  character(len=:), allocatable :: str
56 
57  write (str, '(*(G0))') &
58  'Binary file header (pos: ', this%pos, &
59  ', kper: ', this%kper, &
60  ', kstp: ', this%kstp, &
61  ', pertim: ', this%pertim, &
62  ', totim: ', this%totim, &
63  ')'
64  str = trim(str)
65  end function get_str
66 
67  !> @brief Build an index of data sizes in the file
68  subroutine build_index(this)
69  class(binaryfilereadertype), intent(inout) :: this
70  ! local
71  logical(LGP) :: success
72  integer(I4B) :: pos, record_size
73 
74  this%indexed = .false.
75  this%nrecords = 0
76  call this%rewind()
77  if (allocated(this%record_sizes)) deallocate (this%record_sizes)
78 
79  ! first pass: count
80  do
81  call this%read_record(success)
82  if (.not. success) exit
83  this%nrecords = this%nrecords + 1
84  end do
85  call this%rewind()
86  allocate (this%record_sizes(this%nrecords))
87  ! second pass: save
88  do
89  call this%read_record(success)
90  inquire (this%inunit, pos=pos)
91  record_size = pos - this%header%pos
92  if (.not. success) exit
93  this%record_sizes(this%nrecords) = this%header%size + record_size
94  end do
95  call this%rewind()
96 
97  this%indexed = .true.
98  end subroutine build_index
99 
100  !> @brief Peek to see if another record is available.
101  subroutine peek_record(this)
102  class(binaryfilereadertype), intent(inout) :: this
103  ! local
104  integer(I4B) :: iostat
105 
106  if (.not. this%endoffile) then
107  inquire (this%inunit, pos=this%headernext%pos)
108  read (this%inunit, iostat=iostat) this%headernext%kstp, this%headernext%kper
109  if (iostat == 0) then
110  call fseek_stream(this%inunit, -2 * i4b, 1, iostat)
111  else if (iostat < 0) then
112  this%endoffile = .true.
113  if (allocated(this%headernext)) deallocate (this%headernext)
114  end if
115  end if
116  end subroutine peek_record
117 
118  !> @brief Rewind the file to the beginning.
119  subroutine rewind (this)
120  class(binaryfilereadertype), intent(inout) :: this
121 
122  rewind(this%inunit)
123  if (allocated(this%header)) deallocate (this%header)
124  if (allocated(this%headernext)) deallocate (this%headernext)
125  allocate (binaryfileheadertype :: this%header)
126  allocate (binaryfileheadertype :: this%headernext)
127  this%header%pos = 1
128  this%headernext%pos = 1
129  this%endoffile = .false.
130  end subroutine rewind
131 
132 end module binaryfilereadermodule
character(len=:) function, allocatable get_str(this)
Get a string representation of the header.
subroutine rewind(this)
Rewind the file to the beginning.
subroutine build_index(this)
Build an index of data sizes in the file.
subroutine peek_record(this)
Peek to see if another record is available.
subroutine pstop(status, message)
Stop the program, optionally specifying an error status code.
Definition: ErrorUtil.f90:24
subroutine, public fseek_stream(iu, offset, whence, status)
Move the file pointer.
This module defines variable data types.
Definition: kind.f90:8