MODFLOW 6  version 6.7.0.dev2
USGS Modular Hydrologic Model
MethodCellPollock.f90
Go to the documentation of this file.
2 
3  use kindmodule, only: dp, i4b
4  use constantsmodule, only: done, dzero
11  use particlemodule, only: particletype
12  implicit none
13 
14  private
15  public :: methodcellpollocktype
17 
19  contains
20  procedure, public :: apply => apply_mcp
21  procedure, public :: deallocate => destroy_mcp
22  procedure, public :: load => load_mcp
23  procedure, public :: load_subcell
24  procedure, public :: pass => pass_mcp
25  end type methodcellpollocktype
26 
27 contains
28 
29  !> @brief Create a tracking method
30  subroutine create_method_cell_pollock(method)
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
44  end subroutine create_method_cell_pollock
45 
46  !> @brief Destroy the tracking method
47  subroutine destroy_mcp(this)
48  class(methodcellpollocktype), intent(inout) :: this
49  deallocate (this%name)
50  end subroutine destroy_mcp
51 
52  !> @brief Load subcell tracking method
53  subroutine load_mcp(this, particle, next_level, submethod)
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  events=this%events, &
70  tracktimes=this%tracktimes)
71  submethod => method_subcell_plck
72  particle%itrdomain(next_level) = 1
73  end subroutine load_mcp
74 
75  !> @brief Having exited the lone subcell, pass the particle to the cell face
76  !! In this case the lone subcell is the cell.
77  subroutine pass_mcp(this, particle)
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(level_subfeature)
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(level_feature) = 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(level_feature) = entryface
115  end if
116  end subroutine pass_mcp
117 
118  !> @brief Apply Pollock's method to a rectangular cell
119  subroutine apply_mcp(this, particle, tmax)
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  call this%assess(particle, this%cell%defn, tmax)
134  if (.not. particle%advancing) return
135 
136  ! Transform model coordinates to local cell coordinates
137  ! (translated/rotated but not scaled relative to model)
138  xorigin = cell%xOrigin
139  yorigin = cell%yOrigin
140  zorigin = cell%zOrigin
141  sinrot = cell%sinrot
142  cosrot = cell%cosrot
143  call particle%transform(xorigin, yorigin, zorigin, &
144  sinrot, cosrot)
145 
146  ! Track the particle over the cell
147  call this%track(particle, 2, tmax)
148 
149  ! Transform cell coordinates back to model coordinates
150  call particle%transform(xorigin, yorigin, zorigin, &
151  sinrot, cosrot, invert=.true.)
152  call particle%reset_transform()
153  end select
154  end subroutine apply_mcp
155 
156  !> @brief Loads the lone rectangular subcell from the rectangular cell
157  subroutine load_subcell(this, particle, subcell) !
158  ! dummy
159  class(methodcellpollocktype), intent(inout) :: this
160  type(particletype), pointer, intent(inout) :: particle
161  type(subcellrecttype), intent(inout) :: subcell
162 
163  select type (cell => this%cell)
164  type is (cellrecttype)
165  ! Set subcell number to 1
166  subcell%isubcell = 1
167 
168  ! Subcell calculations will be done in local subcell coordinates
169  subcell%dx = cell%dx
170  subcell%dy = cell%dy
171  subcell%dz = cell%dz
172  subcell%sinrot = dzero
173  subcell%cosrot = done
174  subcell%xOrigin = dzero
175  subcell%yOrigin = dzero
176  subcell%zOrigin = dzero
177 
178  ! Set subcell edge velocities
179  subcell%vx1 = cell%vx1 ! cell velocities already account for retfactor and porosity
180  subcell%vx2 = cell%vx2
181  subcell%vy1 = cell%vy1
182  subcell%vy2 = cell%vy2
183  subcell%vz1 = cell%vz1
184  subcell%vz2 = cell%vz2
185  end select
186  end subroutine load_subcell
187 
188 end module methodcellpollockmodule
subroutine, public create_cell_rect(cell)
Create a new rectangular cell.
Definition: CellRect.f90:40
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
real(dp), parameter done
real constant 1
Definition: Constants.f90:76
This module defines variable data types.
Definition: kind.f90:8
subroutine pass_mcp(this, particle)
Having exited the lone subcell, pass the particle to the cell face In this case the lone subcell is t...
procedure subroutine, public create_method_cell_pollock(method)
Create a tracking method.
subroutine load_subcell(this, particle, subcell)
Loads the lone rectangular subcell from the rectangular cell.
subroutine destroy_mcp(this)
Destroy the tracking method.
subroutine load_mcp(this, particle, next_level, submethod)
Load subcell tracking method.
subroutine apply_mcp(this, particle, tmax)
Apply Pollock's method to a rectangular cell.
Particle tracking strategies.
Definition: Method.f90:2
@, public level_feature
Definition: Method.f90:37
@, public level_subfeature
Definition: Method.f90:38
Subcell-level tracking methods.
type(methodsubcellpollocktype), pointer, public method_subcell_plck
type(methodsubcellternarytype), pointer, public method_subcell_tern
subroutine, public create_subcell_rect(subcell)
Create a new rectangular subcell.
Definition: SubcellRect.f90:28
Base type for particle tracking methods.
Definition: Method.f90:54
Particle tracked by the PRT model.
Definition: Particle.f90:56
A subcell of a cell.
Definition: Subcell.f90:9