MODFLOW 6  version 6.7.0.dev2
USGS Modular Hydrologic Model
methodcellternarymodule Module Reference

Data Types

type  methodcellternarytype
 

Functions/Subroutines

subroutine, public create_method_cell_ternary (method)
 Create a tracking method. More...
 
subroutine destroy_mct (this)
 Destroy the tracking method. More...
 
subroutine load_mct (this, particle, next_level, submethod)
 Load subcell into tracking method. More...
 
subroutine pass_mct (this, particle)
 Pass particle to next subcell if there is one, or to the cell face. More...
 
subroutine apply_mct (this, particle, tmax)
 Apply the ternary method to a polygonal cell. More...
 
subroutine load_subcell (this, particle, subcell)
 Loads a triangular subcell from the polygonal cell. More...
 
subroutine vertvelo (this)
 Calculate vertex velocities. More...
 
subroutine calc_thru_hcsum (this, vm0ival, divcell, le, li, lm, areasub, areacell, unintx, uninty, unextx, unexty, unintxnext, unintynext, kappax, kappay, vm0x, vm0y, vm1x, vm1y, hcsum)
 

Function/Subroutine Documentation

◆ apply_mct()

subroutine methodcellternarymodule::apply_mct ( class(methodcellternarytype), intent(inout)  this,
type(particletype), intent(inout), pointer  particle,
real(dp), intent(in)  tmax 
)
private

Definition at line 154 of file MethodCellTernary.f90.

155  use constantsmodule, only: dzero, done, dhalf
156  ! dummy
157  class(MethodCellTernaryType), intent(inout) :: this
158  type(ParticleType), pointer, intent(inout) :: particle
159  real(DP), intent(in) :: tmax
160  ! local
161  integer(I4B) :: i
162  real(DP) :: x, y, z, xO, yO
163  real(DP), allocatable :: xs(:), ys(:)
164 
165  ! (Re)allocate type-bound arrays
166  select type (cell => this%cell)
167  type is (cellpolytype)
168  call this%assess(particle, this%cell%defn, tmax)
169  if (.not. particle%advancing) return
170 
171  ! Number of vertices
172  this%nverts = cell%defn%npolyverts
173 
174  ! (Re)allocate type-bound arrays
175  if (allocated(this%xvert)) then
176  deallocate (this%xvert)
177  deallocate (this%yvert)
178  deallocate (this%vne)
179  deallocate (this%vv0x)
180  deallocate (this%vv0y)
181  deallocate (this%vv1x)
182  deallocate (this%vv1y)
183  deallocate (this%iprev)
184  deallocate (this%xvertnext)
185  deallocate (this%yvertnext)
186  end if
187 
188  allocate (this%xvert(this%nverts))
189  allocate (this%yvert(this%nverts))
190  allocate (this%vne(this%nverts))
191  allocate (this%vv0x(this%nverts))
192  allocate (this%vv0y(this%nverts))
193  allocate (this%vv1x(this%nverts))
194  allocate (this%vv1y(this%nverts))
195  allocate (this%iprev(this%nverts))
196  allocate (this%xvertnext(this%nverts))
197  allocate (this%yvertnext(this%nverts))
198 
199  ! New origin is the corner of the cell's
200  ! bounding box closest to the old origin
201  allocate (xs(this%nverts))
202  allocate (ys(this%nverts))
203  xs = cell%defn%polyvert(1, :)
204  ys = cell%defn%polyvert(2, :)
205  xo = xs(minloc(abs(xs), dim=1))
206  yo = ys(minloc(abs(ys), dim=1))
207  deallocate (xs)
208  deallocate (ys)
209 
210  ! Cell vertices
211  do i = 1, this%nverts
212  x = cell%defn%polyvert(1, i)
213  y = cell%defn%polyvert(2, i)
214  call transform(x, y, dzero, x, y, z, xo, yo)
215  this%xvert(i) = x
216  this%yvert(i) = y
217  end do
218 
219  ! Top, bottom, and thickness
220  this%ztop = cell%defn%top
221  this%zbot = cell%defn%bot
222  this%dz = this%ztop - this%zbot
223 
224  ! Shifted arrays
225  do i = 1, this%nverts
226  this%iprev(i) = i
227  end do
228  this%iprev = cshift(this%iprev, -1)
229  this%xvertnext = cshift(this%xvert, 1)
230  this%yvertnext = cshift(this%yvert, 1)
231 
232  ! Calculate vertex velocities.
233  call this%vertvelo()
234 
235  ! Transform particle coordinates
236  call particle%transform(xo, yo)
237 
238  ! Track the particle across the cell.
239  call this%track(particle, 2, tmax)
240 
241  ! Transform particle coordinates back
242  call particle%transform(xo, yo, invert=.true.)
243  call particle%reset_transform()
244  end select
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dhalf
real constant 1/2
Definition: Constants.f90:68
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
real(dp), parameter done
real constant 1
Definition: Constants.f90:76

◆ calc_thru_hcsum()

subroutine methodcellternarymodule::calc_thru_hcsum ( class(methodcellternarytype), intent(inout)  this,
real(dp)  vm0ival,
real(dp)  divcell,
real(dp), dimension(:)  le,
real(dp), dimension(:)  li,
real(dp), dimension(:)  lm,
real(dp), dimension(:)  areasub,
real(dp)  areacell,
real(dp), dimension(:)  unintx,
real(dp), dimension(:)  uninty,
real(dp), dimension(:)  unextx,
real(dp), dimension(:)  unexty,
real(dp), dimension(:)  unintxnext,
real(dp), dimension(:)  unintynext,
real(dp), dimension(:)  kappax,
real(dp), dimension(:)  kappay,
real(dp), dimension(:)  vm0x,
real(dp), dimension(:)  vm0y,
real(dp), dimension(:)  vm1x,
real(dp), dimension(:)  vm1y,
real(dp)  hcsum 
)

Definition at line 548 of file MethodCellTernary.f90.

553  ! dummy
554  class(MethodCellTernaryType), intent(inout) :: this
555  real(DP) :: vm0ival
556  real(DP) :: divcell
557  real(DP) :: hcsum
558  real(DP), dimension(:) :: le
559  real(DP), dimension(:) :: li
560  real(DP), dimension(:) :: lm
561  real(DP), dimension(:) :: areasub
562  real(DP) :: areacell
563  real(DP), dimension(:) :: unintx
564  real(DP), dimension(:) :: uninty
565  real(DP), dimension(:) :: unextx
566  real(DP), dimension(:) :: unexty
567  real(DP), dimension(:) :: unintxnext
568  real(DP), dimension(:) :: unintynext
569  real(DP), dimension(:) :: kappax
570  real(DP), dimension(:) :: kappay
571  real(DP), dimension(:) :: vm0x
572  real(DP), dimension(:) :: vm0y
573  real(DP), dimension(:) :: vm1x
574  real(DP), dimension(:) :: vm1y
575  ! local
576  real(DP), allocatable, dimension(:) :: vm0i
577  real(DP), allocatable, dimension(:) :: vm0e
578  real(DP), allocatable, dimension(:) :: vm1i
579  real(DP), allocatable, dimension(:) :: vm1e
580  real(DP), allocatable, dimension(:) :: uprod
581  real(DP), allocatable, dimension(:) :: det
582  real(DP), allocatable, dimension(:) :: wt
583  real(DP), allocatable, dimension(:) :: bi0x
584  real(DP), allocatable, dimension(:) :: be0x
585  real(DP), allocatable, dimension(:) :: bi0y
586  real(DP), allocatable, dimension(:) :: be0y
587  real(DP), allocatable, dimension(:) :: bi1x
588  real(DP), allocatable, dimension(:) :: be1x
589  real(DP), allocatable, dimension(:) :: bi1y
590  real(DP), allocatable, dimension(:) :: be1y
591  real(DP), allocatable, dimension(:) :: be01x
592  real(DP), allocatable, dimension(:) :: be01y
593  real(DP) :: emxx
594  real(DP) :: emxy
595  real(DP) :: emyx
596  real(DP) :: emyy
597  real(DP) :: rx
598  real(DP) :: ry
599  real(DP) :: emdet
600  integer(I4B) :: i
601  integer(I4B) :: ip
602 
603  ! Allocate local arrays
604  allocate (vm0i(this%nverts))
605  allocate (vm0e(this%nverts))
606  allocate (vm1i(this%nverts))
607  allocate (vm1e(this%nverts))
608  allocate (uprod(this%nverts))
609  allocate (det(this%nverts))
610  allocate (wt(this%nverts))
611  allocate (bi0x(this%nverts))
612  allocate (be0x(this%nverts))
613  allocate (bi0y(this%nverts))
614  allocate (be0y(this%nverts))
615  allocate (bi1x(this%nverts))
616  allocate (be1x(this%nverts))
617  allocate (bi1y(this%nverts))
618  allocate (be1y(this%nverts))
619  allocate (be01x(this%nverts))
620  allocate (be01y(this%nverts))
621 
622  ! Set vm0i(1)
623  vm0i(1) = vm0ival
624 
625  ! Get remaining vm0i's sequentially using divergence conditions
626  do i = 2, this%nverts
627  ip = this%iprev(i)
628  vm0i(i) = (li(ip) * vm0i(ip) - le(ip) * this%vne(ip) &
629  + areasub(ip) * divcell) / li(i)
630  end do
631 
632  ! Get vm1i's from vm0i's using continuity conditions
633  vm1i = cshift(vm0i, 1)
634 
635  ! Get centroid velocity by setting up and solving 2x2 linear system
636  uprod = unintx * unextx + uninty * unexty
637  det = done - uprod * uprod
638  bi0x = (unintx - unextx * uprod) / det
639  be0x = (unextx - unintx * uprod) / det
640  bi0y = (uninty - unexty * uprod) / det
641  be0y = (unexty - uninty * uprod) / det
642  uprod = unintxnext * unextx + unintynext * unexty
643  det = done - uprod * uprod
644  bi1x = (unintxnext - unextx * uprod) / det
645  be1x = (unextx - unintxnext * uprod) / det
646  bi1y = (unintynext - unexty * uprod) / det
647  be1y = (unexty - unintynext * uprod) / det
648  be01x = 5.d-1 * (be0x + be1x)
649  be01y = 5.d-1 * (be0y + be1y)
650  wt = areasub / areacell
651  emxx = dtwo - sum(wt * be01x * unextx)
652  emxy = -sum(wt * be01x * unexty)
653  emyx = -sum(wt * be01y * unextx)
654  emyy = dtwo - sum(wt * be01y * unexty)
655  rx = sum(wt * (bi0x * vm0i + bi1x * vm1i + be01x * this%vne))
656  ry = sum(wt * (bi0y * vm0i + bi1y * vm1i + be01y * this%vne))
657  emdet = emxx * emyy - emxy * emyx
658  this%vctrx = (emyy * rx - emxy * ry) / emdet
659  this%vctry = (emxx * ry - emyx * rx) / emdet
660 
661  ! Get vm0e's using "known" conditions
662  vm0e = 5.d-1 * (this%vne + unextx * this%vctrx + unexty * this%vctry)
663 
664  ! Get vm1e's from uniformity along exterior edges
665  vm1e = vm0e
666 
667  ! Transform vm0 and vm1 to (x, y) coordinates
668  vm0x = bi0x * vm0i + be0x * vm0e
669  vm0y = bi0y * vm0i + be0y * vm0e
670  vm1x = bi1x * vm1i + be1x * vm0e
671  vm1y = bi1y * vm1i + be1y * vm0e
672 
673  ! Calculate head-cycle summation (which is proportional to
674  ! the curl of the head gradient)
675  hcsum = sum(lm * (kappax * (vm0x + vm1x) + kappay * (vm0y + vm1y)))
676 
677  ! Deallocate local arrays
678  deallocate (vm0i)
679  deallocate (vm0e)
680  deallocate (vm1i)
681  deallocate (vm1e)
682  deallocate (uprod)
683  deallocate (det)
684  deallocate (wt)
685  deallocate (bi0x)
686  deallocate (be0x)
687  deallocate (bi0y)
688  deallocate (be0y)
689  deallocate (bi1x)
690  deallocate (be1x)
691  deallocate (bi1y)
692  deallocate (be1y)
693  deallocate (be01x)
694  deallocate (be01y)
695 

◆ create_method_cell_ternary()

subroutine, public methodcellternarymodule::create_method_cell_ternary ( type(methodcellternarytype), pointer  method)

Definition at line 57 of file MethodCellTernary.f90.

58  ! dummy
59  type(MethodCellTernaryType), pointer :: method
60  ! local
61  type(CellPolyType), pointer :: cell
62  type(SubcellTriType), pointer :: subcell
63 
64  allocate (method)
65  call create_cell_poly(cell)
66  method%cell => cell
67  method%name => method%cell%type
68  method%delegates = .true.
69  call create_subcell_tri(subcell)
70  method%subcell => subcell
Here is the call graph for this function:
Here is the caller graph for this function:

◆ destroy_mct()

subroutine methodcellternarymodule::destroy_mct ( class(methodcellternarytype), intent(inout)  this)
private

Definition at line 74 of file MethodCellTernary.f90.

75  class(MethodCellTernaryType), intent(inout) :: this
76  deallocate (this%name)

◆ load_mct()

subroutine methodcellternarymodule::load_mct ( class(methodcellternarytype), intent(inout)  this,
type(particletype), intent(inout), pointer  particle,
integer(i4b), intent(in)  next_level,
class(methodtype), intent(inout), pointer  submethod 
)
private

Definition at line 80 of file MethodCellTernary.f90.

81  ! dummy
82  class(MethodCellTernaryType), intent(inout) :: this
83  type(ParticleType), pointer, intent(inout) :: particle
84  integer(I4B), intent(in) :: next_level
85  class(MethodType), pointer, intent(inout) :: submethod
86 
87  select type (subcell => this%subcell)
88  type is (subcelltritype)
89  call this%load_subcell(particle, subcell)
90  end select
91 
92  call method_subcell_tern%init( &
93  fmi=this%fmi, &
94  cell=this%cell, &
95  subcell=this%subcell, &
96  events=this%events, &
97  tracktimes=this%tracktimes)
98  submethod => method_subcell_tern

◆ load_subcell()

subroutine methodcellternarymodule::load_subcell ( class(methodcellternarytype), intent(inout)  this,
type(particletype), intent(inout), pointer  particle,
class(subcelltritype), intent(inout)  subcell 
)

Definition at line 248 of file MethodCellTernary.f90.

249  ! dummy
250  class(MethodCellTernaryType), intent(inout) :: this
251  type(ParticleType), pointer, intent(inout) :: particle
252  class(SubcellTriType), intent(inout) :: subcell
253  ! local
254  integer(I4B) :: ic
255  integer(I4B) :: isc
256  integer(I4B) :: iv0
257  real(DP) :: x0
258  real(DP) :: y0
259  real(DP) :: x1
260  real(DP) :: y1
261  real(DP) :: x2
262  real(DP) :: y2
263  real(DP) :: x1rel
264  real(DP) :: y1rel
265  real(DP) :: x2rel
266  real(DP) :: y2rel
267  real(DP) :: xi
268  real(DP) :: yi
269  real(DP) :: di2
270  real(DP) :: d02
271  real(DP) :: d12
272  real(DP) :: di1
273  real(DP) :: d01
274  real(DP) :: alphai
275  real(DP) :: betai
276 
277  select type (cell => this%cell)
278  type is (cellpolytype)
279  ic = cell%defn%icell
280  subcell%icell = ic
281  isc = particle%itrdomain(level_subfeature)
282  if (isc .le. 0) then
283  xi = particle%x
284  yi = particle%y
285  do iv0 = 1, this%nverts
286  x0 = this%xvert(iv0)
287  y0 = this%yvert(iv0)
288  x1 = this%xvertnext(iv0)
289  y1 = this%yvertnext(iv0)
290  x2 = this%xctr
291  y2 = this%yctr
292  x1rel = x1 - x0
293  y1rel = y1 - y0
294  x2rel = x2 - x0
295  y2rel = y2 - y0
296  di2 = xi * y2rel - yi * x2rel
297  d02 = x0 * y2rel - y0 * x2rel
298  d12 = x1rel * y2rel - y1rel * x2rel
299  di1 = xi * y1rel - yi * x1rel
300  d01 = x0 * y1rel - y0 * x1rel
301  alphai = (di2 - d02) / d12
302  betai = -(di1 - d01) / d12
303  ! assumes particle is in the cell, so no check needed for beta
304  if ((alphai .ge. dzero) .and. &
305  (alphai + betai .le. done)) then
306  isc = iv0
307  exit
308  end if
309  end do
310  if (isc .le. 0) then
311  print *, "error -- initial triangle not found in cell ", ic, &
312  " for particle at ", particle%x, particle%y, particle%z
313 
314  call pstop(1)
315  else
316  particle%itrdomain(level_subfeature) = isc
317  end if
318  end if
319  subcell%isubcell = isc
320 
321  ! Set coordinates and velocities at vertices of triangular subcell
322  iv0 = isc
323  subcell%x0 = this%xvert(iv0)
324  subcell%y0 = this%yvert(iv0)
325  subcell%x1 = this%xvertnext(iv0)
326  subcell%y1 = this%yvertnext(iv0)
327  subcell%x2 = this%xctr
328  subcell%y2 = this%yctr
329  subcell%v0x = this%vv0x(iv0)
330  subcell%v0y = this%vv0y(iv0)
331  subcell%v1x = this%vv1x(iv0) ! the indices here actually refer to subcells, not vertices
332  subcell%v1y = this%vv1y(iv0)
333  subcell%v2x = this%vctrx
334  subcell%v2y = this%vctry
335  subcell%ztop = this%ztop
336  subcell%zbot = this%zbot
337  subcell%dz = this%dz
338  subcell%vzbot = this%vzbot
339  subcell%vztop = this%vztop
340  end select
Here is the call graph for this function:

◆ pass_mct()

subroutine methodcellternarymodule::pass_mct ( class(methodcellternarytype), intent(inout)  this,
type(particletype), intent(inout), pointer  particle 
)
private

Definition at line 102 of file MethodCellTernary.f90.

103  ! dummy
104  class(MethodCellTernaryType), intent(inout) :: this
105  type(ParticleType), pointer, intent(inout) :: particle
106  ! local
107  integer(I4B) :: isc
108  integer(I4B) :: exitFace
109  integer(I4B) :: inface
110 
111  exitface = particle%iboundary(level_subfeature)
112  isc = particle%itrdomain(level_subfeature)
113 
114  select case (exitface)
115  case (0)
116  ! Subcell interior (cell interior)
117  inface = -1
118  case (1)
119  ! Subcell face 1 (cell face)
120  inface = isc
121  if (inface .eq. 0) inface = this%nverts
122  case (2)
123  ! Subcell face --> next subcell in "cycle" (cell interior)
124  isc = isc + 1
125  if (isc .gt. this%nverts) isc = 1
126  particle%itrdomain(level_subfeature) = isc
127  particle%iboundary(level_subfeature) = 3
128  inface = 0
129  case (3)
130  ! Subcell face --> preceding subcell in "cycle" (cell interior)
131  isc = isc - 1
132  if (isc .lt. 1) isc = this%nverts
133  particle%itrdomain(level_subfeature) = isc
134  particle%iboundary(level_subfeature) = 2
135  inface = 0
136  case (4)
137  ! Subcell bottom (cell bottom)
138  inface = this%nverts + 2
139  case (5)
140  ! Subcell top (cell top)
141  inface = this%nverts + 3
142  end select
143 
144  if (inface .eq. -1) then
145  particle%iboundary(level_feature) = 0
146  else if (inface .eq. 0) then
147  particle%iboundary(level_feature) = 0
148  else
149  particle%iboundary(level_feature) = inface
150  end if

◆ vertvelo()

subroutine methodcellternarymodule::vertvelo ( class(methodcellternarytype), intent(inout)  this)
private

Definition at line 344 of file MethodCellTernary.f90.

345  use constantsmodule, only: dzero, done, dhalf
346  ! dummy
347  class(MethodCellTernaryType), intent(inout) :: this
348  ! local
349  real(DP) :: term
350  integer(I4B) :: i
351  real(DP) :: perturb
352  real(DP), allocatable, dimension(:) :: xvals
353  real(DP), allocatable, dimension(:) :: yvals
354  real(DP) :: sixa
355  real(DP) :: vm0i0
356  real(DP) :: vm0ival
357  real(DP) :: hcsum0
358  real(DP) :: hcsum
359  real(DP) :: jac
360  real(DP), allocatable, dimension(:) :: wk1
361  real(DP), allocatable, dimension(:) :: wk2
362  real(DP), allocatable, dimension(:) :: unextxnext
363  real(DP), allocatable, dimension(:) :: unextynext
364  real(DP), allocatable, dimension(:) :: le
365  real(DP), allocatable, dimension(:) :: unextx
366  real(DP), allocatable, dimension(:) :: unexty
367  real(DP) :: areacell
368  real(DP), allocatable, dimension(:) :: areasub
369  real(DP) :: divcell
370  real(DP), allocatable, dimension(:) :: li
371  real(DP), allocatable, dimension(:) :: unintx
372  real(DP), allocatable, dimension(:) :: uninty
373  real(DP), allocatable, dimension(:) :: xmid
374  real(DP), allocatable, dimension(:) :: ymid
375  real(DP), allocatable, dimension(:) :: lm
376  real(DP), allocatable, dimension(:) :: umx
377  real(DP), allocatable, dimension(:) :: umy
378  real(DP), allocatable, dimension(:) :: kappax
379  real(DP), allocatable, dimension(:) :: kappay
380  real(DP), allocatable, dimension(:) :: vm0x
381  real(DP), allocatable, dimension(:) :: vm0y
382  real(DP), allocatable, dimension(:) :: vm1x
383  real(DP), allocatable, dimension(:) :: vm1y
384  ! real(DP), allocatable, dimension(:, :) :: poly
385  integer(I4B) :: nvert
386 
387  select type (cell => this%cell)
388  type is (cellpolytype)
389 
390  ! Allocate local arrays
391  allocate (le(this%nverts)) ! lengths of exterior (cell) edges
392  allocate (unextx(this%nverts)) ! x components of unit normals to exterior edges
393  allocate (unexty(this%nverts)) ! y components of unit normals to exterior edges
394  allocate (areasub(this%nverts)) ! subcell areas
395  allocate (li(this%nverts)) ! lengths of interior edges ("spokes")
396  allocate (unintx(this%nverts)) ! x components of unit normals to interior edges
397  allocate (uninty(this%nverts)) ! y components of unit normals to interior edges
398  allocate (xmid(this%nverts)) ! x coordinates of midpoints
399  allocate (ymid(this%nverts)) ! y coordinates of midpoints
400  allocate (lm(this%nverts)) ! lengths of midpoint connectors
401  allocate (umx(this%nverts)) ! x components of midpoint-connector (cw) unit vectors
402  allocate (umy(this%nverts)) ! y components of midpoint-connector (cw) unit vectors
403  allocate (kappax(this%nverts)) ! x components of kappa vectors
404  allocate (kappay(this%nverts)) ! y components of kappa vectors
405  allocate (vm0x(this%nverts)) ! x component of vm0
406  allocate (vm0y(this%nverts)) ! y component of vm0
407  allocate (vm1x(this%nverts)) ! x component of vm1
408  allocate (vm1y(this%nverts)) ! y component of vm1
409  allocate (unextxnext(this%nverts)) ! vector of "next" interior unit-normal x coordinates defined for convenience
410  allocate (unextynext(this%nverts)) ! vector of "next" interior unit-normal y coordinates defined for convenience
411  allocate (wk1(this%nverts))
412  allocate (wk2(this%nverts))
413  allocate (xvals(3))
414  allocate (yvals(3))
415 
416  ! Exterior edge unit normals (outward) and lengths
417  wk1 = this%xvertnext - this%xvert
418  wk2 = this%yvertnext - this%yvert
419  le = dsqrt(wk1 * wk1 + wk2 * wk2)
420  unextx = -wk2 / le
421  unexty = wk1 / le
422 
423  ! Cell centroid (in general, this is NOT the average of the vertex coordinates)
424  areacell = area(this%xvert, this%yvert)
425  sixa = areacell * dsix
426  wk1 = -(this%xvert * this%yvertnext - this%xvertnext * this%yvert)
427  nvert = size(this%xvert)
428  this%xctr = sum((this%xvert + this%xvertnext) * wk1) / sixa
429  this%yctr = sum((this%yvert + this%yvertnext) * wk1) / sixa
430 
431  ! TODO: can we use some of the terms of the centroid calculation
432  ! to do a cheap point in polygon check?
433 
434  ! Subcell areas
435  do i = 1, this%nverts
436  xvals(1) = this%xvert(i)
437  xvals(2) = this%xvertnext(i)
438  xvals(3) = this%xctr
439  yvals(1) = this%yvert(i)
440  yvals(2) = this%yvertnext(i)
441  yvals(3) = this%yctr
442  areasub(i) = area(xvals, yvals)
443  end do
444 
445  ! Cell-edge normal velocities (outward)
446  term = done / (cell%defn%porosity * cell%defn%retfactor * this%dz)
447  do i = 1, this%nverts
448  this%vne(i) = -cell%defn%faceflow(i) * term / le(i)
449  end do
450 
451  ! Cell divergence (2D)
452  divcell = sum(le * this%vne) / areacell
453 
454  ! Interior edge (cw) unit normals and lengths
455  wk1 = this%xvert - this%xctr
456  wk2 = this%yvert - this%yctr
457  li = dsqrt(wk1 * wk1 + wk2 * wk2)
458  unintx = wk2 / li
459  uninty = -wk1 / li
460  ! Shifted arrays for convenience
461  unextxnext = cshift(unintx, 1)
462  unextynext = cshift(uninty, 1)
463 
464  ! Midpoints of interior edges
465  xmid = 5.d-1 * (this%xvert + this%xctr)
466  ymid = 5.d-1 * (this%yvert + this%yctr)
467 
468  ! Unit midpoint-connector (cw) vectors and lengths
469  wk1 = cshift(xmid, 1) - xmid
470  wk2 = cshift(ymid, 1) - ymid
471  lm = dsqrt(wk1 * wk1 + wk2 * wk2)
472  umx = wk1 / lm
473  umy = wk2 / lm
474 
475  ! Kappa vectors (K tensor times unit midpoint-connector vectors) do not
476  ! account for anisotropy, which is consistent with the way internal face
477  ! flow calculations are done in MP7. The isotropic value of K does not
478  ! matter in this case because it cancels out of the calculations, so
479  ! K = 1 is assumed for simplicity.
480  kappax = umx
481  kappay = umy
482 
483  ! Use linearity to find vm0i[0] such that curl of the head gradient
484  ! is zero
485  perturb = 1.d-2
486  ! Calculations at base value
487  vm0i0 = 0.d0
488  call this%calc_thru_hcsum(vm0i0, divcell, le, li, lm, areasub, areacell, &
489  unintx, uninty, unextx, unexty, &
490  unextxnext, unextynext, &
491  kappax, kappay, vm0x, vm0y, vm1x, vm1y, hcsum0)
492  ! Calculations at perturbed value
493  vm0ival = vm0i0 + perturb
494  call this%calc_thru_hcsum(vm0ival, divcell, le, li, lm, areasub, areacell, &
495  unintx, uninty, unextx, unexty, &
496  unextxnext, unextynext, &
497  kappax, kappay, vm0x, vm0y, vm1x, vm1y, hcsum)
498  ! Calculations at root value
499  jac = (hcsum - hcsum0) / perturb
500  vm0ival = vm0i0 - hcsum0 / jac
501  call this%calc_thru_hcsum(vm0ival, divcell, le, li, lm, areasub, areacell, &
502  unintx, uninty, unextx, unexty, &
503  unextxnext, unextynext, &
504  kappax, kappay, vm0x, vm0y, vm1x, vm1y, hcsum)
505 
506  ! Project linearly to get corner (vertex) velocities. Note that velocity
507  ! vv1 is at the next vertex cw from vv0, so vv0(i) and vv1(i) are the
508  ! two vertex velocities used by triangular subcell i.
509  this%vv0x = 2.d0 * vm0x - this%vctrx
510  this%vv0y = 2.d0 * vm0y - this%vctry
511  this%vv1x = 2.d0 * vm1x - this%vctrx
512  this%vv1y = 2.d0 * vm1y - this%vctry
513 
514  ! Set top and bottom velocities
515  term = done / (cell%defn%retfactor * cell%defn%porosity * areacell)
516  this%vzbot = cell%defn%faceflow(this%nverts + 2) * term
517  this%vztop = -cell%defn%faceflow(this%nverts + 3) * term
518 
519  ! Deallocate local arrays
520  deallocate (le)
521  deallocate (unextx)
522  deallocate (unexty)
523  deallocate (areasub)
524  deallocate (li)
525  deallocate (unintx)
526  deallocate (uninty)
527  deallocate (xmid)
528  deallocate (ymid)
529  deallocate (lm)
530  deallocate (umx)
531  deallocate (umy)
532  deallocate (kappax)
533  deallocate (kappay)
534  deallocate (vm0x)
535  deallocate (vm0y)
536  deallocate (vm1x)
537  deallocate (vm1y)
538  deallocate (unextxnext)
539  deallocate (unextynext)
540  deallocate (wk1)
541  deallocate (wk2)
542  deallocate (xvals)
543  deallocate (yvals)
544 
545  end select