MODFLOW 6  version 6.8.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 128 of file MethodCellPollock.f90.

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

◆ create_method_cell_pollock()

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

Definition at line 31 of file MethodCellPollock.f90.

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

52  class(MethodCellPollockType), intent(inout) :: this
53  deallocate (this%name)
54 
55  ! Deallocate owned subcell method
56  call this%method_subcell_plck%deallocate()
57  deallocate (this%method_subcell_plck)

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

62  ! modules
63  use subcellmodule, only: subcelltype
64  ! dummy
65  class(MethodCellPollockType), intent(inout) :: this
66  type(ParticleType), pointer, intent(inout) :: particle
67  integer, intent(in) :: next_level
68  class(MethodType), pointer, intent(inout) :: submethod
69 
70  select type (subcell => this%subcell)
71  type is (subcellrecttype)
72  call this%load_subcell(particle, subcell)
73  end select
74  call this%method_subcell_plck%init( &
75  fmi=this%fmi, &
76  cell=this%cell, &
77  subcell=this%subcell, &
78  events=this%events, &
79  tracktimes=this%tracktimes)
80  submethod => this%method_subcell_plck
81  particle%itrdomain(next_level) = 1
A subcell of a cell.
Definition: Subcell.f90:10

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

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

◆ pass_mcp()

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

Definition at line 86 of file MethodCellPollock.f90.

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