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 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  ! -- Update particle state, checking whether any reporting or
134  ! -- termination conditions apply
135  call this%update(particle, cell%defn)
136 
137  ! -- Return early if particle is done advancing
138  if (.not. particle%advancing) return
139 
140  ! -- If the particle is above the top of the cell (which is presumed to
141  ! -- represent a water table above the cell bottom), pass the particle
142  ! -- vertically and instantaneously to the cell top elevation and save
143  ! -- the particle state to output file(s).
144  if (particle%z > cell%defn%top) then
145  particle%z = cell%defn%top
146  call this%save(particle, reason=1) ! reason=1: cell transition
147  end if
148 
149  ! Transform particle location into local cell coordinates
150  ! (translated and rotated but not scaled relative to model).
151  ! Transform particle location back to model coordinates, then
152  ! reset transformation and eliminate accumulated roundoff error.
153  xorigin = cell%xOrigin
154  yorigin = cell%yOrigin
155  zorigin = cell%zOrigin
156  sinrot = cell%sinrot
157  cosrot = cell%cosrot
158  call particle%transform(xorigin, yorigin, zorigin, &
159  sinrot, cosrot)
160  call this%track(particle, 2, tmax)
161  call particle%transform(xorigin, yorigin, zorigin, &
162  sinrot, cosrot, invert=.true.)
163  call particle%transform(reset=.true.)
164  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%type => 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%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 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  trackfilectl=this%trackfilectl, &
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 168 of file MethodCellPollock.f90.

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