MODFLOW 6  version 6.7.0.dev2
USGS Modular Hydrologic Model
binaryfilereadermodule Module Reference

Data Types

type  binaryfileheadertype
 
type  binaryfilereadertype
 
interface  read_header_if
 
interface  read_record_if
 

Functions/Subroutines

character(len=:) function, allocatable get_str (this)
 Get a string representation of the header. More...
 
subroutine build_index (this)
 Build an index of data sizes in the file. More...
 
subroutine peek_record (this)
 Peek to see if another record is available. More...
 
subroutine rewind (this)
 Rewind the file to the beginning. More...
 

Function/Subroutine Documentation

◆ build_index()

subroutine binaryfilereadermodule::build_index ( class(binaryfilereadertype), intent(inout)  this)

Definition at line 68 of file BinaryFileReader.f90.

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.

◆ get_str()

character(len=:) function, allocatable binaryfilereadermodule::get_str ( class(binaryfileheadertype), intent(in)  this)

Definition at line 53 of file BinaryFileReader.f90.

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)

◆ peek_record()

subroutine binaryfilereadermodule::peek_record ( class(binaryfilereadertype), intent(inout)  this)

Definition at line 101 of file BinaryFileReader.f90.

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
Here is the call graph for this function:

◆ rewind()

subroutine binaryfilereadermodule::rewind ( class(binaryfilereadertype), intent(inout)  this)

Definition at line 119 of file BinaryFileReader.f90.

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.