MODFLOW 6  version 6.7.0.dev0
USGS Modular Hydrologic Model
tabletermmodule Module Reference

Data Types

type  tabletermtype
 

Functions/Subroutines

subroutine initialize (this, tag, width, alignment)
 
integer(i4b) function get_width (this)
 
integer(i4b) function get_alignment (this)
 
integer(i4b) function get_header_lines (this)
 
subroutine allocate_scalars (this)
 
subroutine da (this)
 
subroutine set_header (this, nlines)
 
subroutine get_header (this, iline, cval)
 

Function/Subroutine Documentation

◆ allocate_scalars()

subroutine tabletermmodule::allocate_scalars ( class(tabletermtype this)

Definition at line 177 of file TableTerm.f90.

178  ! -- modules
179  ! -- dummy
180  class(TableTermType) :: this
181  !
182  ! -- allocate scalars
183  allocate (this%tag)
184  allocate (this%alignment)
185  allocate (this%width)
186  allocate (this%nheader_lines)
187  !
188  ! -- initialize scalars
189  this%nheader_lines = 0

◆ da()

subroutine tabletermmodule::da ( class(tabletermtype this)

Definition at line 194 of file TableTerm.f90.

195  ! -- modules
196  ! -- dummy
197  class(TableTermType) :: this
198  ! -- local
199  !integer(I4B) :: n
200  !
201  ! -- deallocate scalars
202  deallocate (this%tag)
203  deallocate (this%alignment)
204  deallocate (this%width)
205  deallocate (this%nheader_lines)
206  deallocate (this%header_lines)
207  !
208  ! -- return

◆ get_alignment()

integer(i4b) function tabletermmodule::get_alignment ( class(tabletermtype this)

Definition at line 153 of file TableTerm.f90.

154  ! -- return variable
155  integer(I4B) :: get_alignment
156  ! -- modules
157  ! -- dummy
158  class(TableTermType) :: this
159  ! -- local
160  get_alignment = this%alignment

◆ get_header()

subroutine tabletermmodule::get_header ( class(tabletermtype this,
integer(i4b), intent(in)  iline,
character(len=*), intent(inout)  cval 
)

Definition at line 256 of file TableTerm.f90.

257  ! -- modules
258  ! -- dummy
259  class(TableTermType) :: this
260  integer(I4B), intent(in) :: iline
261  character(len=*), intent(inout) :: cval
262  ! -- return variable
263  ! -- local
264  !
265  ! -- set return value
266  cval = this%header_lines(iline) (1:this%width)
267  !
268  ! -- return

◆ get_header_lines()

integer(i4b) function tabletermmodule::get_header_lines ( class(tabletermtype this)

Definition at line 165 of file TableTerm.f90.

166  ! -- return variable
167  integer(I4B) :: get_header_lines
168  ! -- modules
169  ! -- dummy
170  class(TableTermType) :: this
171  ! -- local
172  get_header_lines = this%nheader_lines

◆ get_width()

integer(i4b) function tabletermmodule::get_width ( class(tabletermtype this)

Definition at line 141 of file TableTerm.f90.

142  ! -- return variable
143  integer(I4B) :: get_width
144  ! -- modules
145  ! -- dummy
146  class(TableTermType) :: this
147  ! -- local
148  get_width = this%width

◆ initialize()

subroutine tabletermmodule::initialize ( class(tabletermtype this,
character(len=*), intent(in)  tag,
integer(i4b), intent(in)  width,
integer(i4b), intent(in), optional  alignment 
)

Definition at line 43 of file TableTerm.f90.

44  ! -- modules
45  ! -- dummy
46  class(TableTermType) :: this
47  character(len=*), intent(in) :: tag
48  integer(I4B), intent(in) :: width
49  integer(I4B), intent(in), optional :: alignment
50  ! -- local
51  character(len=LINELENGTH) :: string
52  character(len=LINELENGTH) :: tstring
53  character(len=LINELENGTH), allocatable, dimension(:) :: words
54  integer(I4B) :: nwords
55  integer(I4B) :: ilen
56  integer(I4B) :: i
57  integer(I4B) :: j
58 
59  !
60  ! -- allocate scalars
61  call this%allocate_scalars()
62 
63  ! -- process dummy variables
64  this%tag = tag
65 
66  if (present(alignment)) then
67  this%alignment = alignment
68  else
69  this%alignment = tabcenter
70  end if
71 
72  this%width = width
73  !
74  ! -- parse tag into words
75  call parseline(tag, nwords, words, 0)
76  !
77  ! -- abbreviate any words that exceed the specified width
78  ! and trim trailing characters
79  do i = 1, nwords
80  ilen = len(trim(words(i)))
81  if (ilen > width) then
82  words(i) (width:width) = '.'
83  do j = width + 1, ilen
84  words(i) (j:j) = ' '
85  end do
86  end if
87  end do
88  !
89  ! -- combine words that fit into width
90  i = 0
91  do
92  i = i + 1
93  if (i > nwords) then
94  exit
95  end if
96  string = trim(adjustl(words(i)))
97  tstring = string
98  do j = i + 1, nwords
99  if (len(trim(adjustl(string))) > 0) then
100  tstring = trim(adjustl(tstring))//' '//trim(adjustl(words(j)))
101  else
102  tstring = trim(adjustl(words(j)))
103  end if
104  ilen = len(trim(adjustl(tstring)))
105  if (ilen == 0) then
106  continue
107  else if (ilen <= width) then
108  words(j) = ' '
109  string = tstring
110  else
111  exit
112  end if
113  end do
114  words(i) = string
115  end do
116  !
117  ! -- calculate the number of header lines
118  do i = 1, nwords
119  ilen = len(trim(adjustl(words(i))))
120  if (ilen > 0) then
121  this%nheader_lines = this%nheader_lines + 1
122  end if
123  end do
124  !
125  ! allocate initial_lines and fill with words
126  allocate (this%initial_lines(this%nheader_lines))
127  do i = 1, this%nheader_lines
128  this%initial_lines(i) = words(i) (1:width)
129  end do
130  !
131  ! -- deallocate words
132  deallocate (words)
133  !
134  ! -- return
135  return
136 
Here is the call graph for this function:

◆ set_header()

subroutine tabletermmodule::set_header ( class(tabletermtype this,
integer(i4b), intent(in)  nlines 
)

Definition at line 213 of file TableTerm.f90.

214  ! -- modules
215  ! -- dummy
216  class(TableTermType) :: this
217  integer(I4B), intent(in) :: nlines
218  ! -- local
219  character(len=this%width) :: string
220  integer(I4B) :: idiff
221  integer(I4B) :: i0
222  integer(I4B) :: i
223  integer(I4B) :: j
224  !
225  ! -- initialize variables
226  string = ' '
227  !
228  ! allocate header_lines
229  allocate (this%header_lines(nlines))
230  !
231  ! -- initialize header lines
232  do i = 1, nlines
233  this%header_lines(i) = string
234  end do
235  !
236  ! -- fill header_lines with initial_lines from
237  ! bottom to top
238  idiff = nlines - this%nheader_lines
239  i0 = 1 - idiff
240  do i = this%nheader_lines, 1, -1
241  j = i + idiff
242  this%header_lines(j) = this%initial_lines(i)
243  end do
244  !
245  ! -- deallocate temporary header lines
246  deallocate (this%initial_lines)
247  !
248  ! -- reinitialize nheader_lines
249  this%nheader_lines = nlines
250  !
251  ! -- return