MODFLOW 6  version 6.7.0.dev1
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 119 of file MethodCellPollock.f90.

120  ! dummy
121  class(MethodCellPollockType), intent(inout) :: this
122  type(ParticleType), pointer, intent(inout) :: particle
123  real(DP), intent(in) :: tmax
124  ! local
125  real(DP) :: xOrigin
126  real(DP) :: yOrigin
127  real(DP) :: zOrigin
128  real(DP) :: sinrot
129  real(DP) :: cosrot
130 
131  select type (cell => this%cell)
132  type is (cellrecttype)
133  ! Check termination/reporting conditions
134  call this%check(particle, cell%defn, tmax)
135  if (.not. particle%advancing) return
136 
137  ! Transform model coordinates to local cell coordinates
138  ! (translated/rotated but not scaled relative to model)
139  xorigin = cell%xOrigin
140  yorigin = cell%yOrigin
141  zorigin = cell%zOrigin
142  sinrot = cell%sinrot
143  cosrot = cell%cosrot
144  call particle%transform(xorigin, yorigin, zorigin, &
145  sinrot, cosrot)
146 
147  ! Track the particle over the cell
148  call this%track(particle, 2, tmax)
149 
150  ! Transform cell coordinates back to model coordinates
151  call particle%transform(xorigin, yorigin, zorigin, &
152  sinrot, cosrot, invert=.true.)
153  call particle%reset_transform()
154  end select

◆ create_method_cell_pollock()

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

Definition at line 30 of file MethodCellPollock.f90.

31  ! dummy
32  type(MethodCellPollockType), pointer :: method
33  ! local
34  type(CellRectType), pointer :: cell
35  type(SubcellRectType), pointer :: subcell
36 
37  allocate (method)
38  call create_cell_rect(cell)
39  method%cell => cell
40  method%name => method%cell%type
41  method%delegates = .true.
42  call create_subcell_rect(subcell)
43  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 47 of file MethodCellPollock.f90.

48  class(MethodCellPollockType), intent(inout) :: this
49  deallocate (this%name)

◆ 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 53 of file MethodCellPollock.f90.

54  ! modules
55  use subcellmodule, only: subcelltype
56  ! dummy
57  class(MethodCellPollockType), intent(inout) :: this
58  type(ParticleType), pointer, intent(inout) :: particle
59  integer, intent(in) :: next_level
60  class(MethodType), pointer, intent(inout) :: submethod
61 
62  select type (subcell => this%subcell)
63  type is (subcellrecttype)
64  call this%load_subcell(particle, subcell)
65  end select
66  call method_subcell_plck%init( &
67  cell=this%cell, &
68  subcell=this%subcell, &
69  trackctl=this%trackctl, &
70  tracktimes=this%tracktimes)
71  submethod => method_subcell_plck
72  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 158 of file MethodCellPollock.f90.

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

◆ pass_mcp()

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

Definition at line 77 of file MethodCellPollock.f90.

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