MODFLOW 6  version 6.6.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 150 of file MethodCellTernary.f90.

151  use constantsmodule, only: dzero, done, dhalf
152  ! dummy
153  class(MethodCellTernaryType), intent(inout) :: this
154  type(ParticleType), pointer, intent(inout) :: particle
155  real(DP), intent(in) :: tmax
156  ! local
157  integer(I4B) :: i
158 
159  ! Update particle state, return early if done advancing
160  call this%update(particle, this%cell%defn)
161  if (.not. particle%advancing) return
162 
163  ! If the particle is above the top of the cell (presumed water table)
164  ! pass it vertically and instantaneously to the top
165  if (particle%z > this%cell%defn%top) then
166  particle%z = this%cell%defn%top
167  call this%save(particle, reason=1)
168  end if
169 
170  ! (Re)allocate type-bound arrays
171  select type (cell => this%cell)
172  type is (cellpolytype)
173  ! Number of vertices
174  this%nverts = cell%defn%npolyverts
175  ! (Re)allocate type-bound arrays
176  if (allocated(this%xvert)) then
177  deallocate (this%xvert)
178  deallocate (this%yvert)
179  deallocate (this%vne)
180  deallocate (this%vv0x)
181  deallocate (this%vv0y)
182  deallocate (this%vv1x)
183  deallocate (this%vv1y)
184  deallocate (this%iprev)
185  deallocate (this%xvertnext)
186  deallocate (this%yvertnext)
187  end if
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  ! Cell vertices
199  do i = 1, this%nverts
200  this%xvert(i) = cell%defn%polyvert(1, i)
201  this%yvert(i) = cell%defn%polyvert(2, i)
202  end do
203  ! Top, bottom, and thickness
204  this%ztop = cell%defn%top
205  this%zbot = cell%defn%bot
206  this%dz = this%ztop - this%zbot
207  ! Shifted arrays
208  do i = 1, this%nverts
209  this%iprev(i) = i
210  end do
211  this%iprev = cshift(this%iprev, -1)
212  this%xvertnext = cshift(this%xvert, 1)
213  this%yvertnext = cshift(this%yvert, 1)
214  end select
215 
216  ! Calculate vertex velocities.
217  call this%vertvelo()
218 
219  ! Track the particle across the cell.
220  call this%track(particle, 2, tmax)
221 
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 521 of file MethodCellTernary.f90.

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

◆ create_method_cell_ternary()

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

Definition at line 55 of file MethodCellTernary.f90.

56  ! dummy
57  type(MethodCellTernaryType), pointer :: method
58  ! local
59  type(CellPolyType), pointer :: cell
60  type(SubcellTriType), pointer :: subcell
61 
62  allocate (method)
63  call create_cell_poly(cell)
64  method%cell => cell
65  method%type => method%cell%type
66  method%delegates = .true.
67  call create_subcell_tri(subcell)
68  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 72 of file MethodCellTernary.f90.

73  class(MethodCellTernaryType), intent(inout) :: this
74  deallocate (this%type)

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

79  ! dummy
80  class(MethodCellTernaryType), intent(inout) :: this
81  type(ParticleType), pointer, intent(inout) :: particle
82  integer(I4B), intent(in) :: next_level
83  class(MethodType), pointer, intent(inout) :: submethod
84 
85  select type (subcell => this%subcell)
86  type is (subcelltritype)
87  call this%load_subcell(particle, subcell)
88  end select
89  call method_subcell_tern%init( &
90  fmi=this%fmi, &
91  cell=this%cell, &
92  subcell=this%subcell, &
93  trackctl=this%trackctl, &
94  tracktimes=this%tracktimes)
95  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 225 of file MethodCellTernary.f90.

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

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

◆ vertvelo()

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

Definition at line 321 of file MethodCellTernary.f90.

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