MODFLOW 6  version 6.6.0.dev0
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
5  use methodmodule, only: methodtype
10  use particlemodule, only: particletype
11  implicit none
12 
13  private
14  public :: methodcellpollocktype
16 
18  contains
19  procedure, public :: apply => apply_mcp
20  procedure, public :: deallocate => destroy_mcp
21  procedure, public :: load => load_mcp
22  procedure, public :: load_subcell
23  procedure, public :: pass => pass_mcp
24  end type methodcellpollocktype
25 
26 contains
27 
28  !> @brief Create a tracking method
29  subroutine create_method_cell_pollock(method)
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%name => method%cell%type
40  method%delegates = .true.
41  call create_subcell_rect(subcell)
42  method%subcell => subcell
43  end subroutine create_method_cell_pollock
44 
45  !> @brief Destroy the tracking method
46  subroutine destroy_mcp(this)
47  class(methodcellpollocktype), intent(inout) :: this
48  deallocate (this%name)
49  end subroutine destroy_mcp
50 
51  !> @brief Load subcell tracking method
52  subroutine load_mcp(this, particle, next_level, submethod)
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
72  end subroutine load_mcp
73 
74  !> @brief Having exited the lone subcell, pass the particle to the cell face
75  !! In this case the lone subcell is the cell.
76  subroutine pass_mcp(this, particle)
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
115  end subroutine pass_mcp
116 
117  !> @brief Apply Pollock's method to a rectangular cell
118  subroutine apply_mcp(this, particle, tmax)
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  ! Check termination/reporting conditions
133  call this%check(particle, cell%defn)
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
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:31
Particle tracked by the PRT model.
Definition: Particle.f90:32
A subcell of a cell.
Definition: Subcell.f90:9