MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
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
38  !
39  ! -- Return
40  return
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dpi
real constant
Definition: Constants.f90:127
real(dp), parameter dtwo
real constant 2
Definition: Constants.f90:78

◆ area_wet()

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

Definition at line 62 of file CircularGeometry.f90.

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

◆ perimeter_sat()

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

Definition at line 45 of file CircularGeometry.f90.

46  ! -- modules
47  use constantsmodule, only: dtwo, dpi
48  ! -- return
49  real(DP) :: perimeter_sat
50  ! -- dummy
51  class(CircularGeometryType) :: this
52  !
53  ! -- Calculate area
54  perimeter_sat = dtwo * dpi * this%radius
55  !
56  ! -- Return
57  return

◆ perimeter_wet()

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

Definition at line 94 of file CircularGeometry.f90.

95  ! -- modules
96  use constantsmodule, only: dtwo, dpi
97  ! -- return
98  real(DP) :: perimeter_wet
99  ! -- dummy
100  class(CircularGeometryType) :: this
101  real(DP), intent(in) :: depth
102  !
103  ! -- Calculate area
104  if (depth <= dzero) then
105  perimeter_wet = dzero
106  elseif (depth <= this%radius) then
107  perimeter_wet = dtwo * this%radius * acos((this%radius - depth) / &
108  this%radius)
109  elseif (depth <= dtwo * this%radius) then
110  perimeter_wet = dtwo * this%radius * (dpi - acos((depth - this%radius) / &
111  this%radius))
112  else
113  perimeter_wet = dtwo * dpi * this%radius
114  end if
115  !
116  ! -- Return
117  return

◆ print_attributes()

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

Definition at line 161 of file CircularGeometry.f90.

162  ! -- dummy
163  class(CircularGeometryType) :: this
164  ! -- local
165  integer(I4B), intent(in) :: iout
166  ! -- formats
167  character(len=*), parameter :: fmtnm = "(4x,a,a)"
168  character(len=*), parameter :: fmttd = "(4x,a,1(1PG15.6))"
169  !
170  ! -- call parent to print parent attributes
171  call this%BaseGeometryType%print_attributes(iout)
172  !
173  ! -- Print specifics of this geometry type
174  write (iout, fmttd) 'RADIUS = ', this%radius
175  write (iout, fmttd) 'SATURATED AREA = ', this%area_sat()
176  write (iout, fmttd) 'SATURATED WETTED PERIMETER = ', this%perimeter_sat()
177  !
178  ! -- Return
179  return

◆ set_attribute()

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

Definition at line 122 of file CircularGeometry.f90.

123  ! -- module
124  use inputoutputmodule, only: urword
125  use constantsmodule, only: linelength
127  ! -- dummy
128  class(CircularGeometryType) :: this
129  character(len=LINELENGTH) :: errmsg
130  character(len=*), intent(inout) :: line
131  ! -- local
132  integer(I4B) :: lloc, istart, istop, ival
133  real(DP) :: rval
134  !
135  ! -- should change this and set id if uninitialized or store it
136  lloc = 1
137  call urword(line, lloc, istart, istop, 2, ival, rval, 0, 0)
138  this%id = ival
139 
140  ! -- Parse the attribute
141  call urword(line, lloc, istart, istop, 1, ival, rval, 0, 0)
142  select case (line(istart:istop))
143  case ('NAME')
144  call urword(line, lloc, istart, istop, 1, ival, rval, 0, 0)
145  this%name = line(istart:istop)
146  case ('RADIUS')
147  call urword(line, lloc, istart, istop, 3, ival, rval, 0, 0)
148  this%radius = rval
149  case default
150  write (errmsg, '(a,a)') &
151  'Unknown circular geometry attribute: ', line(istart:istop)
152  call store_error(errmsg, terminate=.true.)
153  end select
154  !
155  ! -- Return
156  return
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:44
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: