MODFLOW 6  version 6.6.0.dev0
USGS Modular Hydrologic Model
circulargeometrymodule Module Reference

Data Types

type  circulargeometrytype
 

Functions/Subroutines

real(dp) function area_sat (this)
 Return area as if geometry is fully saturated. More...
 
real(dp) function perimeter_sat (this)
 Return perimeter as if geometry is fully saturated. 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 circular object. More...
 
subroutine print_attributes (this, iout)
 Print the attributes for this object. More...
 

Function/Subroutine Documentation

◆ area_sat()

real(dp) function circulargeometrymodule::area_sat ( class(circulargeometrytype this)
private

Definition at line 28 of file CircularGeometry.f90.

29  ! -- modules
30  use constantsmodule, only: dtwo, dpi
31  ! -- return
32  real(DP) :: area_sat
33  ! -- dummy
34  class(CircularGeometryType) :: this
35  !
36  ! -- Calculate area
37  area_sat = dpi * this%radius**dtwo
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 circulargeometrymodule::area_wet ( class(circulargeometrytype this,
real(dp), intent(in)  depth 
)

Definition at line 56 of file CircularGeometry.f90.

57  ! -- modules
58  use constantsmodule, only: dtwo, dpi, dzero
59  ! -- return
60  real(DP) :: area_wet
61  ! -- dummy
62  class(CircularGeometryType) :: this
63  real(DP), intent(in) :: depth
64  !
65  ! -- Calculate area
66  if (depth <= dzero) then
67  area_wet = dzero
68  elseif (depth <= this%radius) then
69  area_wet = this%radius * this%radius * &
70  acos((this%radius - depth) / this%radius) - &
71  (this%radius - depth) * &
72  sqrt(this%radius * this%radius - (this%radius - depth)**dtwo)
73  elseif (depth <= dtwo * this%radius) then
74  area_wet = this%radius * this%radius * &
75  (dpi - acos((depth - this%radius) / this%radius)) - &
76  (this%radius - depth) * &
77  sqrt(this%radius * this%radius - (this%radius - depth)**dtwo)
78  else
79  area_wet = dpi * this%radius * this%radius
80  end if
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65

◆ perimeter_sat()

real(dp) function circulargeometrymodule::perimeter_sat ( class(circulargeometrytype this)

Definition at line 42 of file CircularGeometry.f90.

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

◆ perimeter_wet()

real(dp) function circulargeometrymodule::perimeter_wet ( class(circulargeometrytype this,
real(dp), intent(in)  depth 
)

Definition at line 85 of file CircularGeometry.f90.

86  ! -- modules
87  use constantsmodule, only: dtwo, dpi
88  ! -- return
89  real(DP) :: perimeter_wet
90  ! -- dummy
91  class(CircularGeometryType) :: this
92  real(DP), intent(in) :: depth
93  !
94  ! -- Calculate area
95  if (depth <= dzero) then
96  perimeter_wet = dzero
97  elseif (depth <= this%radius) then
98  perimeter_wet = dtwo * this%radius * acos((this%radius - depth) / &
99  this%radius)
100  elseif (depth <= dtwo * this%radius) then
101  perimeter_wet = dtwo * this%radius * (dpi - acos((depth - this%radius) / &
102  this%radius))
103  else
104  perimeter_wet = dtwo * dpi * this%radius
105  end if

◆ print_attributes()

subroutine circulargeometrymodule::print_attributes ( class(circulargeometrytype this,
integer(i4b), intent(in)  iout 
)

Definition at line 146 of file CircularGeometry.f90.

147  ! -- dummy
148  class(CircularGeometryType) :: this
149  ! -- local
150  integer(I4B), intent(in) :: iout
151  ! -- formats
152  character(len=*), parameter :: fmtnm = "(4x,a,a)"
153  character(len=*), parameter :: fmttd = "(4x,a,1(1PG15.6))"
154  !
155  ! -- call parent to print parent attributes
156  call this%BaseGeometryType%print_attributes(iout)
157  !
158  ! -- Print specifics of this geometry type
159  write (iout, fmttd) 'RADIUS = ', this%radius
160  write (iout, fmttd) 'SATURATED AREA = ', this%area_sat()
161  write (iout, fmttd) 'SATURATED WETTED PERIMETER = ', this%perimeter_sat()

◆ set_attribute()

subroutine circulargeometrymodule::set_attribute ( class(circulargeometrytype this,
character(len=*), intent(inout)  line 
)

Definition at line 110 of file CircularGeometry.f90.

111  ! -- module
112  use inputoutputmodule, only: urword
113  use constantsmodule, only: linelength
115  ! -- dummy
116  class(CircularGeometryType) :: this
117  character(len=LINELENGTH) :: errmsg
118  character(len=*), intent(inout) :: line
119  ! -- local
120  integer(I4B) :: lloc, istart, istop, ival
121  real(DP) :: rval
122  !
123  ! -- should change this and set id if uninitialized or store it
124  lloc = 1
125  call urword(line, lloc, istart, istop, 2, ival, rval, 0, 0)
126  this%id = ival
127 
128  ! -- Parse the attribute
129  call urword(line, lloc, istart, istop, 1, ival, rval, 0, 0)
130  select case (line(istart:istop))
131  case ('NAME')
132  call urword(line, lloc, istart, istop, 1, ival, rval, 0, 0)
133  this%name = line(istart:istop)
134  case ('RADIUS')
135  call urword(line, lloc, istart, istop, 3, ival, rval, 0, 0)
136  this%radius = rval
137  case default
138  write (errmsg, '(a,a)') &
139  'Unknown circular geometry attribute: ', line(istart:istop)
140  call store_error(errmsg, terminate=.true.)
141  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: