MODFLOW 6  version 6.7.0.dev0
USGS Modular Hydrologic Model
PackageMover.f90
Go to the documentation of this file.
2 
3  use kindmodule, only: dp, i4b
7 
8  implicit none
9  private
10  public :: packagemovertype
11  public :: set_packagemover_pointer
13 
15 
16  character(len=LENMEMPATH) :: memorypath !< the location in the memory manager where the variables are stored
17  integer(I4B), pointer :: nproviders
18  integer(I4B), pointer :: nreceivers
19  integer(I4B), dimension(:), pointer, contiguous :: iprmap => null() !< map between id1 and feature (needed for lake to map from outlet to lake number)
20  real(dp), dimension(:), pointer, contiguous :: qtformvr => null() !< total flow rate available for mover
21  real(dp), dimension(:), pointer, contiguous :: qformvr => null() !< currently available consumed water (changes during fc)
22  real(dp), dimension(:), pointer, contiguous :: qtomvr => null() !< actual amount of water sent to mover
23  real(dp), dimension(:), pointer, contiguous :: qfrommvr => null() !< actual amount of water received from mover
24  real(dp), dimension(:), pointer, contiguous :: qfrommvr0 => null() !< qfrommvr from previous iteration
25 
26  contains
27  procedure :: ar
28  procedure :: ad
29  procedure :: reset
30  procedure :: fc
31  procedure :: da
32  procedure :: allocate_scalars
33  procedure :: allocate_arrays
34  procedure :: get_qfrommvr
35  procedure :: get_qfrommvr0
36  procedure :: get_qtomvr
37  procedure :: accumulate_qformvr
38 
39  end type packagemovertype
40 
41 contains
42 
43  subroutine set_packagemover_pointer(packagemover, memPath)
44  type(packagemovertype), intent(inout) :: packagemover
45  character(len=*), intent(in) :: mempath
46  packagemover%memoryPath = mempath
47  call mem_setptr(packagemover%nproviders, 'NPROVIDERS', &
48  packagemover%memoryPath)
49  call mem_setptr(packagemover%nreceivers, 'NRECEIVERS', &
50  packagemover%memoryPath)
51  call mem_setptr(packagemover%iprmap, 'IPRMAP', packagemover%memoryPath)
52  call mem_setptr(packagemover%qtformvr, 'QTFORMVR', packagemover%memoryPath)
53  call mem_setptr(packagemover%qformvr, 'QFORMVR', packagemover%memoryPath)
54  call mem_setptr(packagemover%qtomvr, 'QTOMVR', packagemover%memoryPath)
55  call mem_setptr(packagemover%qfrommvr, 'QFROMMVR', packagemover%memoryPath)
56  call mem_setptr(packagemover%qfrommvr0, 'QFROMMVR0', packagemover%memoryPath)
57  end subroutine set_packagemover_pointer
58 
59  subroutine nulllify_packagemover_pointer(packagemover)
60  type(packagemovertype), intent(inout) :: packagemover
61  packagemover%memoryPath = ''
62  packagemover%nproviders => null()
63  packagemover%nreceivers => null()
64  packagemover%iprmap => null()
65  packagemover%qtformvr => null()
66  packagemover%qformvr => null()
67  packagemover%qtomvr => null()
68  packagemover%qfrommvr => null()
69  packagemover%qfrommvr0 => null()
70  end subroutine nulllify_packagemover_pointer
71 
72  subroutine ar(this, nproviders, nreceivers, memoryPath)
73  class(packagemovertype) :: this
74  integer, intent(in) :: nproviders
75  integer, intent(in) :: nreceivers
76  character(len=*), intent(in) :: memoryPath
77  !
78  this%memoryPath = memorypath
79  call this%allocate_scalars()
80  this%nproviders = nproviders
81  this%nreceivers = nreceivers
82  !
83  call this%allocate_arrays()
84  end subroutine ar
85 
86  subroutine ad(this)
87  class(packagemovertype) :: this
88  integer :: i
89  !
90  ! -- set qtomvr and qformvr to zero
91  do i = 1, this%nproviders
92  this%qtomvr(i) = dzero
93  this%qformvr(i) = dzero
94  end do
95  end subroutine ad
96 
97  subroutine reset(this)
98  class(packagemovertype) :: this
99  integer :: i
100  !
101  ! -- set frommvr and qtomvr to zero
102  do i = 1, this%nreceivers
103  this%qfrommvr0(i) = this%qfrommvr(i)
104  this%qfrommvr(i) = dzero
105  end do
106  do i = 1, this%nproviders
107  this%qtomvr(i) = dzero
108  this%qtformvr(i) = this%qformvr(i)
109  end do
110  end subroutine reset
111 
112  subroutine fc(this)
113  class(packagemovertype) :: this
114  integer :: i
115  !
116  ! -- set formvr to zero
117  do i = 1, this%nproviders
118  this%qformvr(i) = dzero
119  end do
120  end subroutine fc
121 
122  subroutine da(this)
123  class(packagemovertype) :: this
124  !
125  ! -- arrays
126  call mem_deallocate(this%iprmap)
127  call mem_deallocate(this%qtformvr)
128  call mem_deallocate(this%qformvr)
129  call mem_deallocate(this%qtomvr)
130  call mem_deallocate(this%qfrommvr)
131  call mem_deallocate(this%qfrommvr0)
132  !
133  ! -- scalars
134  call mem_deallocate(this%nproviders)
135  call mem_deallocate(this%nreceivers)
136  !
137  ! -- pointers
138  nullify (this%iprmap)
139  end subroutine da
140 
141  subroutine allocate_scalars(this)
142  class(packagemovertype) :: this
143  !
144  call mem_allocate(this%nproviders, 'NPROVIDERS', this%memoryPath)
145  call mem_allocate(this%nreceivers, 'NRECEIVERS', this%memoryPath)
146  !
147  this%nproviders = 0
148  this%nreceivers = 0
149  end subroutine allocate_scalars
150 
151  subroutine allocate_arrays(this)
152  class(packagemovertype) :: this
153  integer(I4B) :: i
154  !
155  call mem_allocate(this%iprmap, this%nproviders, 'IPRMAP', this%memoryPath)
156  call mem_allocate(this%qtformvr, this%nproviders, 'QTFORMVR', this%memoryPath)
157  call mem_allocate(this%qformvr, this%nproviders, 'QFORMVR', this%memoryPath)
158  call mem_allocate(this%qtomvr, this%nproviders, 'QTOMVR', this%memoryPath)
159  call mem_allocate(this%qfrommvr, this%nreceivers, 'QFROMMVR', this%memoryPath)
160  call mem_allocate(this%qfrommvr0, this%nreceivers, 'QFROMMVR0', &
161  this%memoryPath)
162  !
163  ! -- initialize
164  do i = 1, this%nproviders
165  this%iprmap(i) = i
166  this%qtformvr(i) = dzero
167  this%qformvr(i) = dzero
168  this%qtomvr(i) = dzero
169  end do
170  do i = 1, this%nreceivers
171  this%qfrommvr(i) = dzero
172  this%qfrommvr0(i) = dzero
173  end do
174  end subroutine allocate_arrays
175 
176  function get_qfrommvr(this, ireceiver) result(qfrommvr)
177  class(packagemovertype) :: this
178  real(dp) :: qfrommvr
179  integer, intent(in) :: ireceiver
180  qfrommvr = this%qfrommvr(ireceiver)
181  end function get_qfrommvr
182 
183  function get_qfrommvr0(this, ireceiver) result(qfrommvr0)
184  class(packagemovertype) :: this
185  real(dp) :: qfrommvr0
186  integer, intent(in) :: ireceiver
187  qfrommvr0 = this%qfrommvr0(ireceiver)
188  end function get_qfrommvr0
189 
190  function get_qtomvr(this, iprovider) result(qtomvr)
191  class(packagemovertype) :: this
192  real(dp) :: qtomvr
193  integer, intent(in) :: iprovider
194  qtomvr = this%qtomvr(iprovider)
195  end function get_qtomvr
196 
197  subroutine accumulate_qformvr(this, iprovider, qformvr)
198  class(packagemovertype) :: this
199  integer, intent(in) :: iprovider
200  real(DP), intent(in) :: qformvr
201  this%qformvr(iprovider) = this%qformvr(iprovider) + qformvr
202  end subroutine accumulate_qformvr
203 
204 end module packagemovermodule
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
integer(i4b), parameter lenmempath
maximum length of the memory path
Definition: Constants.f90:27
This module defines variable data types.
Definition: kind.f90:8
real(dp) function get_qtomvr(this, iprovider)
subroutine ar(this, nproviders, nreceivers, memoryPath)
subroutine ad(this)
subroutine reset(this)
subroutine da(this)
real(dp) function get_qfrommvr(this, ireceiver)
subroutine, public nulllify_packagemover_pointer(packagemover)
subroutine, public set_packagemover_pointer(packagemover, memPath)
real(dp) function get_qfrommvr0(this, ireceiver)
subroutine accumulate_qformvr(this, iprovider, qformvr)
subroutine allocate_arrays(this)
subroutine fc(this)
subroutine allocate_scalars(this)