MODFLOW 6  version 6.8.0.dev0
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 163 of file MethodCellTernary.f90.

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

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

◆ create_method_cell_ternary()

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

Definition at line 59 of file MethodCellTernary.f90.

60  ! dummy
61  type(MethodCellTernaryType), pointer :: method
62  ! local
63  type(CellPolyType), pointer :: cell
64  type(SubcellTriType), pointer :: subcell
65 
66  allocate (method)
67  call create_cell_poly(cell)
68  method%cell => cell
69  method%name => method%cell%type
70  method%delegates = .true.
71  call create_subcell_tri(subcell)
72  method%subcell => subcell
73 
74  ! Create subcell method this method can delegate to
75  call create_method_subcell_ternary(method%method_subcell_tern)
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 79 of file MethodCellTernary.f90.

80  class(MethodCellTernaryType), intent(inout) :: this
81  deallocate (this%name)
82 
83  ! Deallocate owned subcell method
84  call this%method_subcell_tern%deallocate()
85  deallocate (this%method_subcell_tern)

◆ 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 89 of file MethodCellTernary.f90.

90  ! dummy
91  class(MethodCellTernaryType), intent(inout) :: this
92  type(ParticleType), pointer, intent(inout) :: particle
93  integer(I4B), intent(in) :: next_level
94  class(MethodType), pointer, intent(inout) :: submethod
95 
96  select type (subcell => this%subcell)
97  type is (subcelltritype)
98  call this%load_subcell(particle, subcell)
99  end select
100 
101  call this%method_subcell_tern%init( &
102  fmi=this%fmi, &
103  cell=this%cell, &
104  subcell=this%subcell, &
105  events=this%events, &
106  tracktimes=this%tracktimes)
107  submethod => this%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 255 of file MethodCellTernary.f90.

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

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

◆ vertvelo()

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

Definition at line 334 of file MethodCellTernary.f90.

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