MODFLOW 6  version 6.6.0.dev0
USGS Modular Hydrologic Model
MethodCellTernary.f90
Go to the documentation of this file.
2 
3  use kindmodule, only: dp, i4b
4  use errorutilmodule, only: pstop
5  use methodmodule
10  use particlemodule
12  use geomutilmodule, only: area
13  use constantsmodule, only: dzero, done, dtwo
14  implicit none
15 
16  private
17  public :: methodcellternarytype
19 
21  private
22  integer(I4B) :: nverts !< number of vertices
23  real(dp), allocatable, dimension(:) :: xvert
24  real(dp), allocatable, dimension(:) :: yvert !< cell vertex coordinates
25  real(dp), allocatable, dimension(:) :: vne !< cell edge normal velocities
26  real(dp), allocatable, dimension(:) :: vv0x
27  real(dp), allocatable, dimension(:) :: vv0y
28  real(dp), allocatable, dimension(:) :: vv1x
29  real(dp), allocatable, dimension(:) :: vv1y !< cell vertex velocities
30  real(dp) :: xctr
31  real(dp) :: yctr !< cell center coordinates
32  real(dp) :: vctrx
33  real(dp) :: vctry !< cell center velocities
34  real(dp) :: ztop
35  real(dp) :: zbot !< cell top and bottom elevations
36  real(dp) :: dz !< cell thickness
37  real(dp) :: vztop
38  real(dp) :: vzbot !< cell top and bottom velocities
39  integer(I4B), allocatable, dimension(:) :: iprev !< array of shifted indices
40  real(dp), allocatable, dimension(:) :: xvertnext
41  real(dp), allocatable, dimension(:) :: yvertnext !< arrays of shifted cell vertex coordinates
42  integer(I4B), public, pointer :: zeromethod
43  contains
44  procedure, public :: apply => apply_mct
45  procedure, public :: deallocate => destroy_mct
46  procedure, public :: load => load_mct
47  procedure, public :: load_subcell
48  procedure, public :: pass => pass_mct
49  procedure :: vertvelo
50  procedure :: calc_thru_hcsum
51  end type methodcellternarytype
52 
53 contains
54 
55  !> @brief Create a tracking method
56  subroutine create_method_cell_ternary(method)
57  ! dummy
58  type(methodcellternarytype), pointer :: method
59  ! local
60  type(cellpolytype), pointer :: cell
61  type(subcelltritype), pointer :: subcell
62 
63  allocate (method)
64  call create_cell_poly(cell)
65  method%cell => cell
66  method%type => method%cell%type
67  method%delegates = .true.
68  call create_subcell_tri(subcell)
69  method%subcell => subcell
70  end subroutine create_method_cell_ternary
71 
72  !> @brief Destroy the tracking method
73  subroutine destroy_mct(this)
74  class(methodcellternarytype), intent(inout) :: this
75  deallocate (this%type)
76  end subroutine destroy_mct
77 
78  !> @brief Load subcell into tracking method
79  subroutine load_mct(this, particle, next_level, submethod)
80  ! dummy
81  class(methodcellternarytype), intent(inout) :: this
82  type(particletype), pointer, intent(inout) :: particle
83  integer(I4B), intent(in) :: next_level
84  class(methodtype), pointer, intent(inout) :: submethod
85 
86  select type (subcell => this%subcell)
87  type is (subcelltritype)
88  call this%load_subcell(particle, subcell)
89  end select
90  call method_subcell_tern%init( &
91  cell=this%cell, &
92  subcell=this%subcell, &
93  trackfilectl=this%trackfilectl, &
94  tracktimes=this%tracktimes)
95  submethod => method_subcell_tern
96  end subroutine load_mct
97 
98  !> @brief Pass particle to next subcell if there is one, or to the cell face
99  subroutine pass_mct(this, particle)
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
147  end subroutine pass_mct
148 
149  !> @brief Apply the ternary method to a polygonal cell
150  subroutine apply_mct(this, particle, tmax)
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, checking whether any reporting or
160  ! termination conditions apply
161  call this%update(particle, this%cell%defn)
162 
163  ! Return early if particle is done advancing
164  if (.not. particle%advancing) return
165 
166  ! If the particle is above the top of the cell (presumed water table)
167  ! pass it vertically and instantaneously to the cell top and save the
168  ! particle state to file
169  if (particle%z > this%cell%defn%top) then
170  particle%z = this%cell%defn%top
171  call this%save(particle, reason=1) ! reason=1: cell transition
172  end if
173 
174  select type (cell => this%cell)
175  type is (cellpolytype)
176  ! Number of vertices
177  this%nverts = cell%defn%npolyverts
178  ! (Re)allocate type-bound arrays
179  if (allocated(this%xvert)) then
180  deallocate (this%xvert)
181  deallocate (this%yvert)
182  deallocate (this%vne)
183  deallocate (this%vv0x)
184  deallocate (this%vv0y)
185  deallocate (this%vv1x)
186  deallocate (this%vv1y)
187  deallocate (this%iprev)
188  deallocate (this%xvertnext)
189  deallocate (this%yvertnext)
190  end if
191  allocate (this%xvert(this%nverts))
192  allocate (this%yvert(this%nverts))
193  allocate (this%vne(this%nverts))
194  allocate (this%vv0x(this%nverts))
195  allocate (this%vv0y(this%nverts))
196  allocate (this%vv1x(this%nverts))
197  allocate (this%vv1y(this%nverts))
198  allocate (this%iprev(this%nverts))
199  allocate (this%xvertnext(this%nverts))
200  allocate (this%yvertnext(this%nverts))
201  ! Cell vertices
202  do i = 1, this%nverts
203  this%xvert(i) = cell%defn%polyvert(1, i)
204  this%yvert(i) = cell%defn%polyvert(2, i)
205  end do
206  ! Top, bottom, and thickness
207  this%ztop = cell%defn%top
208  this%zbot = cell%defn%bot
209  this%dz = this%ztop - this%zbot
210  ! Shifted arrays
211  do i = 1, this%nverts
212  this%iprev(i) = i
213  end do
214  this%iprev = cshift(this%iprev, -1)
215  this%xvertnext = cshift(this%xvert, 1)
216  this%yvertnext = cshift(this%yvert, 1)
217  end select
218 
219  ! Calculate vertex velocities
220  call this%vertvelo()
221 
222  ! Track across subcells
223  call this%track(particle, 2, tmax)
224 
225  end subroutine apply_mct
226 
227  !> @brief Loads a triangular subcell from the polygonal cell
228  subroutine load_subcell(this, particle, subcell)
229  ! dummy
230  class(methodcellternarytype), intent(inout) :: this
231  type(particletype), pointer, intent(inout) :: particle
232  class(subcelltritype), intent(inout) :: subcell
233  ! local
234  integer(I4B) :: ic
235  integer(I4B) :: isc
236  integer(I4B) :: iv0
237  real(DP) :: x0
238  real(DP) :: y0
239  real(DP) :: x1
240  real(DP) :: y1
241  real(DP) :: x2
242  real(DP) :: y2
243  real(DP) :: x1rel
244  real(DP) :: y1rel
245  real(DP) :: x2rel
246  real(DP) :: y2rel
247  real(DP) :: xi
248  real(DP) :: yi
249  real(DP) :: di2
250  real(DP) :: d02
251  real(DP) :: d12
252  real(DP) :: di1
253  real(DP) :: d01
254  real(DP) :: alphai
255  real(DP) :: betai
256 
257  select type (cell => this%cell)
258  type is (cellpolytype)
259  ic = cell%defn%icell
260  subcell%icell = ic
261  isc = particle%idomain(3)
262  if (isc .le. 0) then
263  xi = particle%x
264  yi = particle%y
265  do iv0 = 1, this%nverts
266  x0 = this%xvert(iv0)
267  y0 = this%yvert(iv0)
268  x1 = this%xvertnext(iv0)
269  y1 = this%yvertnext(iv0)
270  x2 = this%xctr
271  y2 = this%yctr
272  x1rel = x1 - x0
273  y1rel = y1 - y0
274  x2rel = x2 - x0
275  y2rel = y2 - y0
276  di2 = xi * y2rel - yi * x2rel
277  d02 = x0 * y2rel - y0 * x2rel
278  d12 = x1rel * y2rel - y1rel * x2rel
279  di1 = xi * y1rel - yi * x1rel
280  d01 = x0 * y1rel - y0 * x1rel
281  alphai = (di2 - d02) / d12
282  betai = -(di1 - d01) / d12
283  ! assumes particle is in the cell, so no check needed for beta
284  if ((alphai .ge. dzero) .and. &
285  (alphai + betai .le. done)) then
286  isc = iv0
287  exit
288  end if
289  end do
290  if (isc .le. 0) then
291  print *, "error -- initial triangle not found in cell ", ic, &
292  " for particle at ", particle%x, particle%y, particle%z
293 
294  call pstop(1)
295  else
296  particle%idomain(3) = isc
297  end if
298  end if
299  subcell%isubcell = isc
300 
301  ! Set coordinates and velocities at vertices of triangular subcell
302  iv0 = isc
303  subcell%x0 = this%xvert(iv0)
304  subcell%y0 = this%yvert(iv0)
305  subcell%x1 = this%xvertnext(iv0)
306  subcell%y1 = this%yvertnext(iv0)
307  subcell%x2 = this%xctr
308  subcell%y2 = this%yctr
309  subcell%v0x = this%vv0x(iv0)
310  subcell%v0y = this%vv0y(iv0)
311  subcell%v1x = this%vv1x(iv0) ! the indices here actually refer to subcells, not vertices
312  subcell%v1y = this%vv1y(iv0)
313  subcell%v2x = this%vctrx
314  subcell%v2y = this%vctry
315  subcell%ztop = this%ztop
316  subcell%zbot = this%zbot
317  subcell%dz = this%dz
318  subcell%vzbot = this%vzbot
319  subcell%vztop = this%vztop
320  end select
321  end subroutine load_subcell
322 
323  !> @brief Calculate vertex velocities
324  subroutine vertvelo(this)
325  use constantsmodule, only: dzero, done, dhalf
326  ! dummy
327  class(methodcellternarytype), intent(inout) :: this
328  ! local
329  real(DP) :: term
330  integer(I4B) :: i
331  real(DP) :: perturb
332  real(DP), allocatable, dimension(:) :: xvals
333  real(DP), allocatable, dimension(:) :: yvals
334  real(DP) :: sixa
335  real(DP) :: vm0i0
336  real(DP) :: vm0ival
337  real(DP) :: hcsum0
338  real(DP) :: hcsum
339  real(DP) :: jac
340  real(DP), allocatable, dimension(:) :: wk1
341  real(DP), allocatable, dimension(:) :: wk2
342  real(DP), allocatable, dimension(:) :: unextxnext
343  real(DP), allocatable, dimension(:) :: unextynext
344  real(DP), allocatable, dimension(:) :: le
345  real(DP), allocatable, dimension(:) :: unextx
346  real(DP), allocatable, dimension(:) :: unexty
347  real(DP) :: areacell
348  real(DP), allocatable, dimension(:) :: areasub
349  real(DP) :: divcell
350  real(DP), allocatable, dimension(:) :: li
351  real(DP), allocatable, dimension(:) :: unintx
352  real(DP), allocatable, dimension(:) :: uninty
353  real(DP), allocatable, dimension(:) :: xmid
354  real(DP), allocatable, dimension(:) :: ymid
355  real(DP), allocatable, dimension(:) :: lm
356  real(DP), allocatable, dimension(:) :: umx
357  real(DP), allocatable, dimension(:) :: umy
358  real(DP), allocatable, dimension(:) :: kappax
359  real(DP), allocatable, dimension(:) :: kappay
360  real(DP), allocatable, dimension(:) :: vm0x
361  real(DP), allocatable, dimension(:) :: vm0y
362  real(DP), allocatable, dimension(:) :: vm1x
363  real(DP), allocatable, dimension(:) :: vm1y
364 
365  select type (cell => this%cell)
366  type is (cellpolytype)
367 
368  ! Allocate local arrays
369  allocate (le(this%nverts)) ! lengths of exterior (cell) edges
370  allocate (unextx(this%nverts)) ! x components of unit normals to exterior edges
371  allocate (unexty(this%nverts)) ! y components of unit normals to exterior edges
372  allocate (areasub(this%nverts)) ! subcell areas
373  allocate (li(this%nverts)) ! lengths of interior edges ("spokes")
374  allocate (unintx(this%nverts)) ! x components of unit normals to interior edges
375  allocate (uninty(this%nverts)) ! y components of unit normals to interior edges
376  allocate (xmid(this%nverts)) ! x coordinates of midpoints
377  allocate (ymid(this%nverts)) ! y coordinates of midpoints
378  allocate (lm(this%nverts)) ! lengths of midpoint connectors
379  allocate (umx(this%nverts)) ! x components of midpoint-connector (cw) unit vectors
380  allocate (umy(this%nverts)) ! y components of midpoint-connector (cw) unit vectors
381  allocate (kappax(this%nverts)) ! x components of kappa vectors
382  allocate (kappay(this%nverts)) ! y components of kappa vectors
383  allocate (vm0x(this%nverts)) ! x component of vm0
384  allocate (vm0y(this%nverts)) ! y component of vm0
385  allocate (vm1x(this%nverts)) ! x component of vm1
386  allocate (vm1y(this%nverts)) ! y component of vm1
387  allocate (unextxnext(this%nverts)) ! vector of "next" interior unit-normal x coordinates defined for convenience
388  allocate (unextynext(this%nverts)) ! vector of "next" interior unit-normal y coordinates defined for convenience
389  allocate (wk1(this%nverts))
390  allocate (wk2(this%nverts))
391  allocate (xvals(3))
392  allocate (yvals(3))
393 
394  ! Exterior edge unit normals (outward) and lengths
395  wk1 = this%xvertnext - this%xvert
396  wk2 = this%yvertnext - this%yvert
397  le = dsqrt(wk1 * wk1 + wk2 * wk2)
398  unextx = -wk2 / le
399  unexty = wk1 / le
400 
401  ! Cell area
402  areacell = area(this%xvert, this%yvert)
403 
404  ! Cell centroid (in general, this is NOT the average of the vertex coordinates)
405  sixa = areacell * 6.d0
406  wk1 = -(this%xvert * this%yvertnext - this%xvertnext * this%yvert)
407  this%xctr = sum((this%xvert + this%xvertnext) * wk1) / sixa
408  this%yctr = sum((this%yvert + this%yvertnext) * wk1) / sixa
409 
410  ! Subcell areas
411  do i = 1, this%nverts
412  xvals(1) = this%xvert(i)
413  xvals(2) = this%xvertnext(i)
414  xvals(3) = this%xctr
415  yvals(1) = this%yvert(i)
416  yvals(2) = this%yvertnext(i)
417  yvals(3) = this%yctr
418  areasub(i) = area(xvals, yvals)
419  end do
420 
421  ! Cell-edge normal velocities (outward)
422  term = done / (cell%defn%porosity * cell%defn%retfactor * this%dz)
423  do i = 1, this%nverts
424  this%vne(i) = -cell%defn%faceflow(i) * term / le(i)
425  end do
426 
427  ! Cell divergence (2D)
428  divcell = sum(le * this%vne) / areacell
429 
430  ! Interior edge (cw) unit normals and lengths
431  wk1 = this%xvert - this%xctr
432  wk2 = this%yvert - this%yctr
433  li = dsqrt(wk1 * wk1 + wk2 * wk2)
434  unintx = wk2 / li
435  uninty = -wk1 / li
436  ! Shifted arrays for convenience
437  unextxnext = cshift(unintx, 1)
438  unextynext = cshift(uninty, 1)
439 
440  ! Midpoints of interior edges
441  xmid = 5.d-1 * (this%xvert + this%xctr)
442  ymid = 5.d-1 * (this%yvert + this%yctr)
443 
444  ! Unit midpoint-connector (cw) vectors and lengths
445  wk1 = cshift(xmid, 1) - xmid
446  wk2 = cshift(ymid, 1) - ymid
447  lm = dsqrt(wk1 * wk1 + wk2 * wk2)
448  umx = wk1 / lm
449  umy = wk2 / lm
450 
451  ! Kappa vectors (K tensor times unit midpoint-connector vectors) do not
452  ! account for anisotropy, which is consistent with the way internal face
453  ! flow calculations are done in MP7. The isotropic value of K does not
454  ! matter in this case because it cancels out of the calculations, so
455  ! K = 1 is assumed for simplicity.
456  kappax = umx
457  kappay = umy
458 
459  ! Use linearity to find vm0i[0] such that curl of the head gradient
460  ! is zero
461  perturb = 1.d-2
462  ! Calculations at base value
463  vm0i0 = 0.d0
464  call this%calc_thru_hcsum(vm0i0, divcell, le, li, lm, areasub, areacell, &
465  unintx, uninty, unextx, unexty, &
466  unextxnext, unextynext, &
467  kappax, kappay, vm0x, vm0y, vm1x, vm1y, hcsum0)
468  ! Calculations at perturbed value
469  vm0ival = vm0i0 + perturb
470  call this%calc_thru_hcsum(vm0ival, divcell, le, li, lm, areasub, areacell, &
471  unintx, uninty, unextx, unexty, &
472  unextxnext, unextynext, &
473  kappax, kappay, vm0x, vm0y, vm1x, vm1y, hcsum)
474  ! Calculations at root value
475  jac = (hcsum - hcsum0) / perturb
476  vm0ival = vm0i0 - hcsum0 / jac
477  call this%calc_thru_hcsum(vm0ival, divcell, le, li, lm, areasub, areacell, &
478  unintx, uninty, unextx, unexty, &
479  unextxnext, unextynext, &
480  kappax, kappay, vm0x, vm0y, vm1x, vm1y, hcsum)
481 
482  ! Project linearly to get corner (vertex) velocities. Note that velocity
483  ! vv1 is at the next vertex cw from vv0, so vv0(i) and vv1(i) are the
484  ! two vertex velocities used by triangular subcell i.
485  this%vv0x = 2.d0 * vm0x - this%vctrx
486  this%vv0y = 2.d0 * vm0y - this%vctry
487  this%vv1x = 2.d0 * vm1x - this%vctrx
488  this%vv1y = 2.d0 * vm1y - this%vctry
489 
490  ! Set top and bottom velocities
491  term = done / (cell%defn%retfactor * cell%defn%porosity * areacell)
492  this%vzbot = cell%defn%faceflow(this%nverts + 2) * term
493  this%vztop = -cell%defn%faceflow(this%nverts + 3) * term
494 
495  ! Deallocate local arrays
496  deallocate (le)
497  deallocate (unextx)
498  deallocate (unexty)
499  deallocate (areasub)
500  deallocate (li)
501  deallocate (unintx)
502  deallocate (uninty)
503  deallocate (xmid)
504  deallocate (ymid)
505  deallocate (lm)
506  deallocate (umx)
507  deallocate (umy)
508  deallocate (kappax)
509  deallocate (kappay)
510  deallocate (vm0x)
511  deallocate (vm0y)
512  deallocate (vm1x)
513  deallocate (vm1y)
514  deallocate (unextxnext)
515  deallocate (unextynext)
516  deallocate (wk1)
517  deallocate (wk2)
518  deallocate (xvals)
519  deallocate (yvals)
520 
521  end select
522  end subroutine vertvelo
523 
524  subroutine calc_thru_hcsum(this, vm0ival, divcell, &
525  le, li, lm, areasub, areacell, &
526  unintx, uninty, unextx, unexty, &
527  unintxnext, unintynext, &
528  kappax, kappay, vm0x, vm0y, vm1x, vm1y, hcsum)
529  ! dummy
530  class(methodcellternarytype), intent(inout) :: this
531  real(DP) :: vm0ival
532  real(DP) :: divcell
533  real(DP) :: hcsum
534  real(DP), dimension(:) :: le
535  real(DP), dimension(:) :: li
536  real(DP), dimension(:) :: lm
537  real(DP), dimension(:) :: areasub
538  real(DP) :: areacell
539  real(DP), dimension(:) :: unintx
540  real(DP), dimension(:) :: uninty
541  real(DP), dimension(:) :: unextx
542  real(DP), dimension(:) :: unexty
543  real(DP), dimension(:) :: unintxnext
544  real(DP), dimension(:) :: unintynext
545  real(DP), dimension(:) :: kappax
546  real(DP), dimension(:) :: kappay
547  real(DP), dimension(:) :: vm0x
548  real(DP), dimension(:) :: vm0y
549  real(DP), dimension(:) :: vm1x
550  real(DP), dimension(:) :: vm1y
551  ! local
552  real(DP), allocatable, dimension(:) :: vm0i
553  real(DP), allocatable, dimension(:) :: vm0e
554  real(DP), allocatable, dimension(:) :: vm1i
555  real(DP), allocatable, dimension(:) :: vm1e
556  real(DP), allocatable, dimension(:) :: uprod
557  real(DP), allocatable, dimension(:) :: det
558  real(DP), allocatable, dimension(:) :: wt
559  real(DP), allocatable, dimension(:) :: bi0x
560  real(DP), allocatable, dimension(:) :: be0x
561  real(DP), allocatable, dimension(:) :: bi0y
562  real(DP), allocatable, dimension(:) :: be0y
563  real(DP), allocatable, dimension(:) :: bi1x
564  real(DP), allocatable, dimension(:) :: be1x
565  real(DP), allocatable, dimension(:) :: bi1y
566  real(DP), allocatable, dimension(:) :: be1y
567  real(DP), allocatable, dimension(:) :: be01x
568  real(DP), allocatable, dimension(:) :: be01y
569  real(DP) :: emxx
570  real(DP) :: emxy
571  real(DP) :: emyx
572  real(DP) :: emyy
573  real(DP) :: rx
574  real(DP) :: ry
575  real(DP) :: emdet
576  integer(I4B) :: i
577  integer(I4B) :: ip
578 
579  ! Allocate local arrays
580  allocate (vm0i(this%nverts))
581  allocate (vm0e(this%nverts))
582  allocate (vm1i(this%nverts))
583  allocate (vm1e(this%nverts))
584  allocate (uprod(this%nverts))
585  allocate (det(this%nverts))
586  allocate (wt(this%nverts))
587  allocate (bi0x(this%nverts))
588  allocate (be0x(this%nverts))
589  allocate (bi0y(this%nverts))
590  allocate (be0y(this%nverts))
591  allocate (bi1x(this%nverts))
592  allocate (be1x(this%nverts))
593  allocate (bi1y(this%nverts))
594  allocate (be1y(this%nverts))
595  allocate (be01x(this%nverts))
596  allocate (be01y(this%nverts))
597 
598  ! Set vm0i(1)
599  vm0i(1) = vm0ival
600 
601  ! Get remaining vm0i's sequentially using divergence conditions
602  do i = 2, this%nverts
603  ip = this%iprev(i)
604  vm0i(i) = (li(ip) * vm0i(ip) - le(ip) * this%vne(ip) &
605  + areasub(ip) * divcell) / li(i)
606  end do
607 
608  ! Get vm1i's from vm0i's using continuity conditions
609  vm1i = cshift(vm0i, 1)
610 
611  ! Get centroid velocity by setting up and solving 2x2 linear system
612  uprod = unintx * unextx + uninty * unexty
613  det = done - uprod * uprod
614  bi0x = (unintx - unextx * uprod) / det
615  be0x = (unextx - unintx * uprod) / det
616  bi0y = (uninty - unexty * uprod) / det
617  be0y = (unexty - uninty * uprod) / det
618  uprod = unintxnext * unextx + unintynext * unexty
619  det = done - uprod * uprod
620  bi1x = (unintxnext - unextx * uprod) / det
621  be1x = (unextx - unintxnext * uprod) / det
622  bi1y = (unintynext - unexty * uprod) / det
623  be1y = (unexty - unintynext * uprod) / det
624  be01x = 5.d-1 * (be0x + be1x)
625  be01y = 5.d-1 * (be0y + be1y)
626  wt = areasub / areacell
627  emxx = dtwo - sum(wt * be01x * unextx)
628  emxy = -sum(wt * be01x * unexty)
629  emyx = -sum(wt * be01y * unextx)
630  emyy = dtwo - sum(wt * be01y * unexty)
631  rx = sum(wt * (bi0x * vm0i + bi1x * vm1i + be01x * this%vne))
632  ry = sum(wt * (bi0y * vm0i + bi1y * vm1i + be01y * this%vne))
633  emdet = emxx * emyy - emxy * emyx
634  this%vctrx = (emyy * rx - emxy * ry) / emdet
635  this%vctry = (emxx * ry - emyx * rx) / emdet
636 
637  ! Get vm0e's using "known" conditions
638  vm0e = 5.d-1 * (this%vne + unextx * this%vctrx + unexty * this%vctry)
639 
640  ! Get vm1e's from uniformity along exterior edges
641  vm1e = vm0e
642 
643  ! Transform vm0 and vm1 to (x, y) coordinates
644  vm0x = bi0x * vm0i + be0x * vm0e
645  vm0y = bi0y * vm0i + be0y * vm0e
646  vm1x = bi1x * vm1i + be1x * vm0e
647  vm1y = bi1y * vm1i + be1y * vm0e
648 
649  ! Calculate head-cycle summation (which is proportional to
650  ! the curl of the head gradient)
651  hcsum = sum(lm * (kappax * (vm0x + vm1x) + kappay * (vm0y + vm1y)))
652 
653  ! Deallocate local arrays
654  deallocate (vm0i)
655  deallocate (vm0e)
656  deallocate (vm1i)
657  deallocate (vm1e)
658  deallocate (uprod)
659  deallocate (det)
660  deallocate (wt)
661  deallocate (bi0x)
662  deallocate (be0x)
663  deallocate (bi0y)
664  deallocate (be0y)
665  deallocate (bi1x)
666  deallocate (be1x)
667  deallocate (bi1y)
668  deallocate (be1y)
669  deallocate (be01x)
670  deallocate (be01y)
671 
672  end subroutine calc_thru_hcsum
673 
674 end module methodcellternarymodule
675 
subroutine, public create_cell_poly(cell)
Create a new polygonal cell.
Definition: CellPoly.f90:20
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 dtwo
real constant 2
Definition: Constants.f90:79
real(dp), parameter done
real constant 1
Definition: Constants.f90:76
subroutine pstop(status, message)
Stop the program, optionally specifying an error status code.
Definition: ErrorUtil.f90:24
real(dp) function, public area(xv, yv, cw)
Calculate polygon area, with vertices given in CW or CCW order.
Definition: GeomUtil.f90:370
This module defines variable data types.
Definition: kind.f90:8
subroutine apply_mct(this, particle, tmax)
Apply the ternary method to a polygonal cell.
subroutine load_mct(this, particle, next_level, submethod)
Load subcell into tracking method.
subroutine destroy_mct(this)
Destroy the tracking method.
subroutine, public create_method_cell_ternary(method)
Create a tracking method.
subroutine load_subcell(this, particle, subcell)
Loads a triangular subcell from the polygonal cell.
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)
subroutine vertvelo(this)
Calculate vertex velocities.
subroutine pass_mct(this, particle)
Pass particle to next subcell if there is one, or to the cell face.
Particle tracking strategies.
Definition: Method.f90:2
subroutine load(this, particle, next_level, submethod)
Load subdomain tracking method (submethod)
Definition: Method.f90:141
subroutine pass(this, particle)
Pass a particle to the next subdomain, internal use only.
Definition: Method.f90:150
Subcell-level tracking methods.
type(methodsubcellternarytype), pointer, public method_subcell_tern
subroutine, public create_subcell_tri(subcell)
Create a new triangular subcell.
Definition: SubcellTri.f90:26
Base type for particle tracking methods.
Definition: Method.f90:29
Particle tracked by the PRT model.
Definition: Particle.f90:32
Manages particle track (i.e. pathline) files.
Definition: TrackData.f90:38