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

Data Types

type  rectangulargeometrytype
 

Functions/Subroutines

real(dp) function area_sat (this)
 Return saturated area. More...
 
real(dp) function perimeter_sat (this)
 Return saturated perimeter. More...
 
real(dp) function area_wet (this, depth)
 Return wetted area. More...
 
real(dp) function perimeter_wet (this, depth)
 Return wetted perimeter. More...
 
subroutine set_attribute (this, line)
 Set a parameter for this rectangular object. More...
 
subroutine print_attributes (this, iout)
 Print the attributes for this object. More...
 

Function/Subroutine Documentation

◆ area_sat()

real(dp) function rectangulargeometrymodule::area_sat ( class(rectangulargeometrytype this)
private

Definition at line 27 of file RectangularGeometry.f90.

28  ! -- modules
29  use constantsmodule, only: dtwo, dpi
30  ! -- return
31  real(DP) :: area_sat
32  ! -- dummy
33  class(RectangularGeometryType) :: this
34  !
35  ! -- Calculate area
36  area_sat = this%height * this%width
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dpi
real constant
Definition: Constants.f90:128
real(dp), parameter dtwo
real constant 2
Definition: Constants.f90:79

◆ area_wet()

real(dp) function rectangulargeometrymodule::area_wet ( class(rectangulargeometrytype this,
real(dp), intent(in)  depth 
)

Definition at line 55 of file RectangularGeometry.f90.

56  ! -- modules
57  use constantsmodule, only: dtwo, dpi, dzero
58  ! -- return
59  real(DP) :: area_wet
60  ! -- dummy
61  class(RectangularGeometryType) :: this
62  real(DP), intent(in) :: depth
63  !
64  ! -- Calculate area
65  if (depth <= dzero) then
66  area_wet = dzero
67  elseif (depth <= this%height) then
68  area_wet = depth * this%width
69  else
70  area_wet = this%width * this%height
71  end if
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65

◆ perimeter_sat()

real(dp) function rectangulargeometrymodule::perimeter_sat ( class(rectangulargeometrytype this)

Definition at line 41 of file RectangularGeometry.f90.

42  ! -- modules
43  use constantsmodule, only: dtwo, dpi
44  ! -- return
45  real(DP) :: perimeter_sat
46  ! -- dummy
47  class(RectangularGeometryType) :: this
48  !
49  ! -- Calculate area
50  perimeter_sat = dtwo * (this%height + this%width)

◆ perimeter_wet()

real(dp) function rectangulargeometrymodule::perimeter_wet ( class(rectangulargeometrytype this,
real(dp), intent(in)  depth 
)

Definition at line 76 of file RectangularGeometry.f90.

77  ! -- modules
78  use constantsmodule, only: dtwo, dpi
79  ! -- return
80  real(DP) :: perimeter_wet
81  ! -- dummy
82  class(RectangularGeometryType) :: this
83  real(DP), intent(in) :: depth
84  !
85  ! -- Calculate area
86  if (depth <= dzero) then
87  perimeter_wet = dzero
88  elseif (depth <= this%height) then
89  perimeter_wet = dtwo * (depth + this%width)
90  else
91  perimeter_wet = dtwo * (this%height + this%width)
92  end if

◆ print_attributes()

subroutine rectangulargeometrymodule::print_attributes ( class(rectangulargeometrytype this,
integer(i4b), intent(in)  iout 
)

Definition at line 136 of file RectangularGeometry.f90.

137  ! -- dummy
138  class(RectangularGeometryType) :: this
139  ! -- local
140  integer(I4B), intent(in) :: iout
141  ! -- formats
142  character(len=*), parameter :: fmtnm = "(4x,a,a)"
143  character(len=*), parameter :: fmttd = "(4x,a,1(1PG15.6))"
144  !
145  ! -- call parent to print parent attributes
146  call this%BaseGeometryType%print_attributes(iout)
147  !
148  ! -- Print specifics of this geometry type
149  write (iout, fmttd) 'HEIGHT = ', this%height
150  write (iout, fmttd) 'WIDTH = ', this%width
151  write (iout, fmttd) 'SATURATED AREA = ', this%area_sat()
152  write (iout, fmttd) 'SATURATED WETTED PERIMETER = ', this%perimeter_sat()

◆ set_attribute()

subroutine rectangulargeometrymodule::set_attribute ( class(rectangulargeometrytype this,
character(len=*), intent(inout)  line 
)

Definition at line 97 of file RectangularGeometry.f90.

98  ! -- module
99  use inputoutputmodule, only: urword
100  use constantsmodule, only: linelength
102  ! -- dummy
103  class(RectangularGeometryType) :: this
104  character(len=LINELENGTH) :: errmsg
105  character(len=*), intent(inout) :: line
106  ! -- local
107  integer(I4B) :: lloc, istart, istop, ival
108  real(DP) :: rval
109  !
110  ! -- should change this and set id if uninitialized or store it
111  lloc = 1
112  call urword(line, lloc, istart, istop, 2, ival, rval, 0, 0)
113  this%id = ival
114  !
115  ! -- Parse the attribute
116  call urword(line, lloc, istart, istop, 1, ival, rval, 0, 0)
117  select case (line(istart:istop))
118  case ('NAME')
119  call urword(line, lloc, istart, istop, 1, ival, rval, 0, 0)
120  this%name = line(istart:istop)
121  case ('HEIGHT')
122  call urword(line, lloc, istart, istop, 3, ival, rval, 0, 0)
123  this%height = rval
124  case ('WIDTH')
125  call urword(line, lloc, istart, istop, 3, ival, rval, 0, 0)
126  this%width = rval
127  case default
128  write (errmsg, '(a,a)') &
129  'Unknown rectangular geometry attribute: ', line(istart:istop)
130  call store_error(errmsg, terminate=.true.)
131  end select
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:45
subroutine, public urword(line, icol, istart, istop, ncode, n, r, iout, in)
Extract a word from a string.
This module contains simulation methods.
Definition: Sim.f90:10
subroutine, public store_error(msg, terminate)
Store an error message.
Definition: Sim.f90:92
integer(i4b) function, public count_errors()
Return number of errors.
Definition: Sim.f90:59
Here is the call graph for this function: