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

Data Types

type  methodcellpollocktype
 

Functions/Subroutines

procedure subroutine, public create_method_cell_pollock (method)
 Create a tracking method. More...
 
subroutine destroy_mcp (this)
 Destroy the tracking method. More...
 
subroutine load_mcp (this, particle, next_level, submethod)
 Load subcell tracking method. More...
 
subroutine pass_mcp (this, particle)
 Having exited the lone subcell, pass the particle to the cell face In this case the lone subcell is the cell. More...
 
subroutine apply_mcp (this, particle, tmax)
 Apply Pollock's method to a rectangular cell. More...
 
subroutine load_subcell (this, particle, subcell)
 Loads the lone rectangular subcell from the rectangular cell. More...
 

Function/Subroutine Documentation

◆ apply_mcp()

subroutine methodcellpollockmodule::apply_mcp ( class(methodcellpollocktype), intent(inout)  this,
type(particletype), intent(inout), pointer  particle,
real(dp), intent(in)  tmax 
)
private

Definition at line 118 of file MethodCellPollock.f90.

119  ! dummy
120  class(MethodCellPollockType), intent(inout) :: this
121  type(ParticleType), pointer, intent(inout) :: particle
122  real(DP), intent(in) :: tmax
123  ! local
124  real(DP) :: xOrigin
125  real(DP) :: yOrigin
126  real(DP) :: zOrigin
127  real(DP) :: sinrot
128  real(DP) :: cosrot
129 
130  select type (cell => this%cell)
131  type is (cellrecttype)
132  ! Update particle state, return early if done advancing
133  call this%update(particle, cell%defn)
134  if (.not. particle%advancing) return
135 
136  ! If the particle is above the top of the cell (presumed water table)
137  ! pass it vertically and instantaneously to the top
138  if (particle%z > cell%defn%top) then
139  particle%z = cell%defn%top
140  call this%save(particle, reason=1)
141  end if
142 
143  ! Transform particle location into local cell coordinates
144  ! (translated and rotated but not scaled relative to model).
145  xorigin = cell%xOrigin
146  yorigin = cell%yOrigin
147  zorigin = cell%zOrigin
148  sinrot = cell%sinrot
149  cosrot = cell%cosrot
150  call particle%transform(xorigin, yorigin, zorigin, &
151  sinrot, cosrot)
152 
153  ! Track the particle across the cell.
154  call this%track(particle, 2, tmax)
155 
156  ! Transform particle location back to model coordinates, then
157  ! reset transformation and eliminate accumulated roundoff error.
158  call particle%transform(xorigin, yorigin, zorigin, &
159  sinrot, cosrot, invert=.true.)
160  call particle%transform(reset=.true.)
161  end select

◆ create_method_cell_pollock()

procedure subroutine, public methodcellpollockmodule::create_method_cell_pollock ( type(methodcellpollocktype), pointer  method)

Definition at line 29 of file MethodCellPollock.f90.

30  ! dummy
31  type(MethodCellPollockType), pointer :: method
32  ! local
33  type(CellRectType), pointer :: cell
34  type(SubcellRectType), pointer :: subcell
35 
36  allocate (method)
37  call create_cell_rect(cell)
38  method%cell => cell
39  method%type => method%cell%type
40  method%delegates = .true.
41  call create_subcell_rect(subcell)
42  method%subcell => subcell
Here is the call graph for this function:
Here is the caller graph for this function:

◆ destroy_mcp()

subroutine methodcellpollockmodule::destroy_mcp ( class(methodcellpollocktype), intent(inout)  this)
private

Definition at line 46 of file MethodCellPollock.f90.

47  class(MethodCellPollockType), intent(inout) :: this
48  deallocate (this%type)

◆ load_mcp()

subroutine methodcellpollockmodule::load_mcp ( class(methodcellpollocktype), intent(inout)  this,
type(particletype), intent(inout), pointer  particle,
integer, intent(in)  next_level,
class(methodtype), intent(inout), pointer  submethod 
)
private

Definition at line 52 of file MethodCellPollock.f90.

53  ! modules
54  use subcellmodule, only: subcelltype
55  ! dummy
56  class(MethodCellPollockType), intent(inout) :: this
57  type(ParticleType), pointer, intent(inout) :: particle
58  integer, intent(in) :: next_level
59  class(MethodType), pointer, intent(inout) :: submethod
60 
61  select type (subcell => this%subcell)
62  type is (subcellrecttype)
63  call this%load_subcell(particle, subcell)
64  end select
65  call method_subcell_plck%init( &
66  cell=this%cell, &
67  subcell=this%subcell, &
68  trackctl=this%trackctl, &
69  tracktimes=this%tracktimes)
70  submethod => method_subcell_plck
71  particle%idomain(next_level) = 1
A subcell of a cell.
Definition: Subcell.f90:9

◆ load_subcell()

subroutine methodcellpollockmodule::load_subcell ( class(methodcellpollocktype), intent(inout)  this,
type(particletype), intent(inout), pointer  particle,
type(subcellrecttype), intent(inout)  subcell 
)
private

Definition at line 165 of file MethodCellPollock.f90.

166  ! dummy
167  class(MethodCellPollockType), intent(inout) :: this
168  type(ParticleType), pointer, intent(inout) :: particle
169  type(SubcellRectType), intent(inout) :: subcell
170 
171  select type (cell => this%cell)
172  type is (cellrecttype)
173  ! Set subcell number to 1
174  subcell%isubcell = 1
175 
176  ! Subcell calculations will be done in local subcell coordinates
177  subcell%dx = cell%dx
178  subcell%dy = cell%dy
179  subcell%dz = cell%dz
180  subcell%sinrot = dzero
181  subcell%cosrot = done
182  subcell%xOrigin = dzero
183  subcell%yOrigin = dzero
184  subcell%zOrigin = dzero
185 
186  ! Set subcell edge velocities
187  subcell%vx1 = cell%vx1 ! cell velocities already account for retfactor and porosity
188  subcell%vx2 = cell%vx2
189  subcell%vy1 = cell%vy1
190  subcell%vy2 = cell%vy2
191  subcell%vz1 = cell%vz1
192  subcell%vz2 = cell%vz2
193  end select

◆ pass_mcp()

subroutine methodcellpollockmodule::pass_mcp ( class(methodcellpollocktype), intent(inout)  this,
type(particletype), intent(inout), pointer  particle 
)

Definition at line 76 of file MethodCellPollock.f90.

77  ! dummy
78  class(MethodCellPollockType), intent(inout) :: this
79  type(ParticleType), pointer, intent(inout) :: particle
80  ! local
81  integer(I4B) :: exitface
82  integer(I4B) :: entryface
83 
84  exitface = particle%iboundary(3)
85  ! Map subcell exit face to cell face
86  select case (exitface) ! note: exitFace uses Dave's iface convention
87  case (0)
88  entryface = -1
89  case (1)
90  entryface = 1
91  case (2)
92  entryface = 3
93  case (3)
94  entryface = 4
95  case (4)
96  entryface = 2
97  case (5)
98  entryface = 6 ! note: inface=5 same as inface=1 due to wraparound
99  case (6)
100  entryface = 7
101  end select
102  if (entryface .eq. -1) then
103  particle%iboundary(2) = 0
104  else
105  if ((entryface .ge. 1) .and. (entryface .le. 4)) then
106  ! Account for local cell rotation
107  select type (cell => this%cell)
108  type is (cellrecttype)
109  entryface = entryface + cell%ipvOrigin - 1
110  end select
111  if (entryface .gt. 4) entryface = entryface - 4
112  end if
113  particle%iboundary(2) = entryface
114  end if