MODFLOW 6  version 6.7.0.dev0
USGS Modular Hydrologic Model
MethodCellTernary.f90
Go to the documentation of this file.
2 
3  use kindmodule, only: dp, i4b, lgp
4  use errorutilmodule, only: pstop
5  use methodmodule
10  use particlemodule
11  use geomutilmodule, only: area, transform
12  use constantsmodule, only: dzero, done, dtwo, dsix
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%name => 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%name)
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 
91  call method_subcell_tern%init( &
92  fmi=this%fmi, &
93  cell=this%cell, &
94  subcell=this%subcell, &
95  trackctl=this%trackctl, &
96  tracktimes=this%tracktimes)
97  submethod => method_subcell_tern
98  end subroutine load_mct
99 
100  !> @brief Pass particle to next subcell if there is one, or to the cell face
101  subroutine pass_mct(this, particle)
102  ! dummy
103  class(methodcellternarytype), intent(inout) :: this
104  type(particletype), pointer, intent(inout) :: particle
105  ! local
106  integer(I4B) :: isc
107  integer(I4B) :: exitFace
108  integer(I4B) :: inface
109 
110  exitface = particle%iboundary(3)
111  isc = particle%idomain(3)
112 
113  select case (exitface)
114  case (0)
115  ! Subcell interior (cell interior)
116  inface = -1
117  case (1)
118  ! Subcell face 1 (cell face)
119  inface = isc
120  if (inface .eq. 0) inface = this%nverts
121  case (2)
122  ! Subcell face --> next subcell in "cycle" (cell interior)
123  isc = isc + 1
124  if (isc .gt. this%nverts) isc = 1
125  particle%idomain(3) = isc
126  particle%iboundary(3) = 3
127  inface = 0
128  case (3)
129  ! Subcell face --> preceding subcell in "cycle" (cell interior)
130  isc = isc - 1
131  if (isc .lt. 1) isc = this%nverts
132  particle%idomain(3) = isc
133  particle%iboundary(3) = 2
134  inface = 0
135  case (4)
136  ! Subcell bottom (cell bottom)
137  inface = this%nverts + 2
138  case (5)
139  ! Subcell top (cell top)
140  inface = this%nverts + 3
141  end select
142 
143  if (inface .eq. -1) then
144  particle%iboundary(2) = 0
145  else if (inface .eq. 0) then
146  particle%iboundary(2) = 0
147  else
148  particle%iboundary(2) = inface
149  end if
150  end subroutine pass_mct
151 
152  !> @brief Apply the ternary method to a polygonal cell
153  subroutine apply_mct(this, particle, tmax)
154  use constantsmodule, only: dzero, done, dhalf
155  ! dummy
156  class(methodcellternarytype), intent(inout) :: this
157  type(particletype), pointer, intent(inout) :: particle
158  real(DP), intent(in) :: tmax
159  ! local
160  integer(I4B) :: i
161  real(DP) :: x, y, z, xO, yO
162  real(DP), allocatable :: xs(:), ys(:)
163 
164  ! (Re)allocate type-bound arrays
165  select type (cell => this%cell)
166  type is (cellpolytype)
167  ! Check termination/reporting conditions
168  call this%check(particle, this%cell%defn)
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 
245  end select
246 
247  end subroutine apply_mct
248 
249  !> @brief Loads a triangular subcell from the polygonal cell
250  subroutine load_subcell(this, particle, subcell)
251  ! dummy
252  class(methodcellternarytype), intent(inout) :: this
253  type(particletype), pointer, intent(inout) :: particle
254  class(subcelltritype), intent(inout) :: subcell
255  ! local
256  integer(I4B) :: ic
257  integer(I4B) :: isc
258  integer(I4B) :: iv0
259  real(DP) :: x0
260  real(DP) :: y0
261  real(DP) :: x1
262  real(DP) :: y1
263  real(DP) :: x2
264  real(DP) :: y2
265  real(DP) :: x1rel
266  real(DP) :: y1rel
267  real(DP) :: x2rel
268  real(DP) :: y2rel
269  real(DP) :: xi
270  real(DP) :: yi
271  real(DP) :: di2
272  real(DP) :: d02
273  real(DP) :: d12
274  real(DP) :: di1
275  real(DP) :: d01
276  real(DP) :: alphai
277  real(DP) :: betai
278 
279  select type (cell => this%cell)
280  type is (cellpolytype)
281  ic = cell%defn%icell
282  subcell%icell = ic
283  isc = particle%idomain(3)
284  if (isc .le. 0) then
285  xi = particle%x
286  yi = particle%y
287  do iv0 = 1, this%nverts
288  x0 = this%xvert(iv0)
289  y0 = this%yvert(iv0)
290  x1 = this%xvertnext(iv0)
291  y1 = this%yvertnext(iv0)
292  x2 = this%xctr
293  y2 = this%yctr
294  x1rel = x1 - x0
295  y1rel = y1 - y0
296  x2rel = x2 - x0
297  y2rel = y2 - y0
298  di2 = xi * y2rel - yi * x2rel
299  d02 = x0 * y2rel - y0 * x2rel
300  d12 = x1rel * y2rel - y1rel * x2rel
301  di1 = xi * y1rel - yi * x1rel
302  d01 = x0 * y1rel - y0 * x1rel
303  alphai = (di2 - d02) / d12
304  betai = -(di1 - d01) / d12
305  ! assumes particle is in the cell, so no check needed for beta
306  if ((alphai .ge. dzero) .and. &
307  (alphai + betai .le. done)) then
308  isc = iv0
309  exit
310  end if
311  end do
312  if (isc .le. 0) then
313  print *, "error -- initial triangle not found in cell ", ic, &
314  " for particle at ", particle%x, particle%y, particle%z
315 
316  call pstop(1)
317  else
318  particle%idomain(3) = isc
319  end if
320  end if
321  subcell%isubcell = isc
322 
323  ! Set coordinates and velocities at vertices of triangular subcell
324  iv0 = isc
325  subcell%x0 = this%xvert(iv0)
326  subcell%y0 = this%yvert(iv0)
327  subcell%x1 = this%xvertnext(iv0)
328  subcell%y1 = this%yvertnext(iv0)
329  subcell%x2 = this%xctr
330  subcell%y2 = this%yctr
331  subcell%v0x = this%vv0x(iv0)
332  subcell%v0y = this%vv0y(iv0)
333  subcell%v1x = this%vv1x(iv0) ! the indices here actually refer to subcells, not vertices
334  subcell%v1y = this%vv1y(iv0)
335  subcell%v2x = this%vctrx
336  subcell%v2y = this%vctry
337  subcell%ztop = this%ztop
338  subcell%zbot = this%zbot
339  subcell%dz = this%dz
340  subcell%vzbot = this%vzbot
341  subcell%vztop = this%vztop
342  end select
343  end subroutine load_subcell
344 
345  !> @brief Calculate vertex velocities
346  subroutine vertvelo(this)
347  use constantsmodule, only: dzero, done, dhalf
348  ! dummy
349  class(methodcellternarytype), intent(inout) :: this
350  ! local
351  real(DP) :: term
352  integer(I4B) :: i
353  real(DP) :: perturb
354  real(DP), allocatable, dimension(:) :: xvals
355  real(DP), allocatable, dimension(:) :: yvals
356  real(DP) :: sixa
357  real(DP) :: vm0i0
358  real(DP) :: vm0ival
359  real(DP) :: hcsum0
360  real(DP) :: hcsum
361  real(DP) :: jac
362  real(DP), allocatable, dimension(:) :: wk1
363  real(DP), allocatable, dimension(:) :: wk2
364  real(DP), allocatable, dimension(:) :: unextxnext
365  real(DP), allocatable, dimension(:) :: unextynext
366  real(DP), allocatable, dimension(:) :: le
367  real(DP), allocatable, dimension(:) :: unextx
368  real(DP), allocatable, dimension(:) :: unexty
369  real(DP) :: areacell
370  real(DP), allocatable, dimension(:) :: areasub
371  real(DP) :: divcell
372  real(DP), allocatable, dimension(:) :: li
373  real(DP), allocatable, dimension(:) :: unintx
374  real(DP), allocatable, dimension(:) :: uninty
375  real(DP), allocatable, dimension(:) :: xmid
376  real(DP), allocatable, dimension(:) :: ymid
377  real(DP), allocatable, dimension(:) :: lm
378  real(DP), allocatable, dimension(:) :: umx
379  real(DP), allocatable, dimension(:) :: umy
380  real(DP), allocatable, dimension(:) :: kappax
381  real(DP), allocatable, dimension(:) :: kappay
382  real(DP), allocatable, dimension(:) :: vm0x
383  real(DP), allocatable, dimension(:) :: vm0y
384  real(DP), allocatable, dimension(:) :: vm1x
385  real(DP), allocatable, dimension(:) :: vm1y
386  ! real(DP), allocatable, dimension(:, :) :: poly
387  integer(I4B) :: nvert
388 
389  select type (cell => this%cell)
390  type is (cellpolytype)
391 
392  ! Allocate local arrays
393  allocate (le(this%nverts)) ! lengths of exterior (cell) edges
394  allocate (unextx(this%nverts)) ! x components of unit normals to exterior edges
395  allocate (unexty(this%nverts)) ! y components of unit normals to exterior edges
396  allocate (areasub(this%nverts)) ! subcell areas
397  allocate (li(this%nverts)) ! lengths of interior edges ("spokes")
398  allocate (unintx(this%nverts)) ! x components of unit normals to interior edges
399  allocate (uninty(this%nverts)) ! y components of unit normals to interior edges
400  allocate (xmid(this%nverts)) ! x coordinates of midpoints
401  allocate (ymid(this%nverts)) ! y coordinates of midpoints
402  allocate (lm(this%nverts)) ! lengths of midpoint connectors
403  allocate (umx(this%nverts)) ! x components of midpoint-connector (cw) unit vectors
404  allocate (umy(this%nverts)) ! y components of midpoint-connector (cw) unit vectors
405  allocate (kappax(this%nverts)) ! x components of kappa vectors
406  allocate (kappay(this%nverts)) ! y components of kappa vectors
407  allocate (vm0x(this%nverts)) ! x component of vm0
408  allocate (vm0y(this%nverts)) ! y component of vm0
409  allocate (vm1x(this%nverts)) ! x component of vm1
410  allocate (vm1y(this%nverts)) ! y component of vm1
411  allocate (unextxnext(this%nverts)) ! vector of "next" interior unit-normal x coordinates defined for convenience
412  allocate (unextynext(this%nverts)) ! vector of "next" interior unit-normal y coordinates defined for convenience
413  allocate (wk1(this%nverts))
414  allocate (wk2(this%nverts))
415  allocate (xvals(3))
416  allocate (yvals(3))
417 
418  ! Exterior edge unit normals (outward) and lengths
419  wk1 = this%xvertnext - this%xvert
420  wk2 = this%yvertnext - this%yvert
421  le = dsqrt(wk1 * wk1 + wk2 * wk2)
422  unextx = -wk2 / le
423  unexty = wk1 / le
424 
425  ! Cell centroid (in general, this is NOT the average of the vertex coordinates)
426  areacell = area(this%xvert, this%yvert)
427  sixa = areacell * dsix
428  wk1 = -(this%xvert * this%yvertnext - this%xvertnext * this%yvert)
429  nvert = size(this%xvert)
430  this%xctr = sum((this%xvert + this%xvertnext) * wk1) / sixa
431  this%yctr = sum((this%yvert + this%yvertnext) * wk1) / sixa
432 
433  ! TODO: can we use some of the terms of the centroid calculation
434  ! to do a cheap point in polygon check?
435  !
436  ! allocate(poly(2, nvert))
437  ! poly(1,:) = this%xvert
438  ! poly(2,:) = this%yvert
439  ! if (.not. point_in_polygon(this%xctr, this%yctr, poly)) then
440  ! print *, "error -- centroid not in cell ", this%cell%defn%icell, this%xctr, this%yctr
441  ! call pstop(1)
442  ! end if
443  ! deallocate(poly)
444 
445  ! Subcell areas
446  do i = 1, this%nverts
447  xvals(1) = this%xvert(i)
448  xvals(2) = this%xvertnext(i)
449  xvals(3) = this%xctr
450  yvals(1) = this%yvert(i)
451  yvals(2) = this%yvertnext(i)
452  yvals(3) = this%yctr
453  areasub(i) = area(xvals, yvals)
454  end do
455 
456  ! Cell-edge normal velocities (outward)
457  term = done / (cell%defn%porosity * cell%defn%retfactor * this%dz)
458  do i = 1, this%nverts
459  this%vne(i) = -cell%defn%faceflow(i) * term / le(i)
460  end do
461 
462  ! Cell divergence (2D)
463  divcell = sum(le * this%vne) / areacell
464 
465  ! Interior edge (cw) unit normals and lengths
466  wk1 = this%xvert - this%xctr
467  wk2 = this%yvert - this%yctr
468  li = dsqrt(wk1 * wk1 + wk2 * wk2)
469  unintx = wk2 / li
470  uninty = -wk1 / li
471  ! Shifted arrays for convenience
472  unextxnext = cshift(unintx, 1)
473  unextynext = cshift(uninty, 1)
474 
475  ! Midpoints of interior edges
476  xmid = 5.d-1 * (this%xvert + this%xctr)
477  ymid = 5.d-1 * (this%yvert + this%yctr)
478 
479  ! Unit midpoint-connector (cw) vectors and lengths
480  wk1 = cshift(xmid, 1) - xmid
481  wk2 = cshift(ymid, 1) - ymid
482  lm = dsqrt(wk1 * wk1 + wk2 * wk2)
483  umx = wk1 / lm
484  umy = wk2 / lm
485 
486  ! Kappa vectors (K tensor times unit midpoint-connector vectors) do not
487  ! account for anisotropy, which is consistent with the way internal face
488  ! flow calculations are done in MP7. The isotropic value of K does not
489  ! matter in this case because it cancels out of the calculations, so
490  ! K = 1 is assumed for simplicity.
491  kappax = umx
492  kappay = umy
493 
494  ! Use linearity to find vm0i[0] such that curl of the head gradient
495  ! is zero
496  perturb = 1.d-2
497  ! Calculations at base value
498  vm0i0 = 0.d0
499  call this%calc_thru_hcsum(vm0i0, divcell, le, li, lm, areasub, areacell, &
500  unintx, uninty, unextx, unexty, &
501  unextxnext, unextynext, &
502  kappax, kappay, vm0x, vm0y, vm1x, vm1y, hcsum0)
503  ! Calculations at perturbed value
504  vm0ival = vm0i0 + perturb
505  call this%calc_thru_hcsum(vm0ival, divcell, le, li, lm, areasub, areacell, &
506  unintx, uninty, unextx, unexty, &
507  unextxnext, unextynext, &
508  kappax, kappay, vm0x, vm0y, vm1x, vm1y, hcsum)
509  ! Calculations at root value
510  jac = (hcsum - hcsum0) / perturb
511  vm0ival = vm0i0 - hcsum0 / jac
512  call this%calc_thru_hcsum(vm0ival, divcell, le, li, lm, areasub, areacell, &
513  unintx, uninty, unextx, unexty, &
514  unextxnext, unextynext, &
515  kappax, kappay, vm0x, vm0y, vm1x, vm1y, hcsum)
516 
517  ! Project linearly to get corner (vertex) velocities. Note that velocity
518  ! vv1 is at the next vertex cw from vv0, so vv0(i) and vv1(i) are the
519  ! two vertex velocities used by triangular subcell i.
520  this%vv0x = 2.d0 * vm0x - this%vctrx
521  this%vv0y = 2.d0 * vm0y - this%vctry
522  this%vv1x = 2.d0 * vm1x - this%vctrx
523  this%vv1y = 2.d0 * vm1y - this%vctry
524 
525  ! Set top and bottom velocities
526  term = done / (cell%defn%retfactor * cell%defn%porosity * areacell)
527  this%vzbot = cell%defn%faceflow(this%nverts + 2) * term
528  this%vztop = -cell%defn%faceflow(this%nverts + 3) * term
529 
530  ! Deallocate local arrays
531  deallocate (le)
532  deallocate (unextx)
533  deallocate (unexty)
534  deallocate (areasub)
535  deallocate (li)
536  deallocate (unintx)
537  deallocate (uninty)
538  deallocate (xmid)
539  deallocate (ymid)
540  deallocate (lm)
541  deallocate (umx)
542  deallocate (umy)
543  deallocate (kappax)
544  deallocate (kappay)
545  deallocate (vm0x)
546  deallocate (vm0y)
547  deallocate (vm1x)
548  deallocate (vm1y)
549  deallocate (unextxnext)
550  deallocate (unextynext)
551  deallocate (wk1)
552  deallocate (wk2)
553  deallocate (xvals)
554  deallocate (yvals)
555 
556  end select
557  end subroutine vertvelo
558 
559  subroutine calc_thru_hcsum(this, vm0ival, divcell, &
560  le, li, lm, areasub, areacell, &
561  unintx, uninty, unextx, unexty, &
562  unintxnext, unintynext, &
563  kappax, kappay, vm0x, vm0y, vm1x, vm1y, hcsum)
564  ! dummy
565  class(methodcellternarytype), intent(inout) :: this
566  real(DP) :: vm0ival
567  real(DP) :: divcell
568  real(DP) :: hcsum
569  real(DP), dimension(:) :: le
570  real(DP), dimension(:) :: li
571  real(DP), dimension(:) :: lm
572  real(DP), dimension(:) :: areasub
573  real(DP) :: areacell
574  real(DP), dimension(:) :: unintx
575  real(DP), dimension(:) :: uninty
576  real(DP), dimension(:) :: unextx
577  real(DP), dimension(:) :: unexty
578  real(DP), dimension(:) :: unintxnext
579  real(DP), dimension(:) :: unintynext
580  real(DP), dimension(:) :: kappax
581  real(DP), dimension(:) :: kappay
582  real(DP), dimension(:) :: vm0x
583  real(DP), dimension(:) :: vm0y
584  real(DP), dimension(:) :: vm1x
585  real(DP), dimension(:) :: vm1y
586  ! local
587  real(DP), allocatable, dimension(:) :: vm0i
588  real(DP), allocatable, dimension(:) :: vm0e
589  real(DP), allocatable, dimension(:) :: vm1i
590  real(DP), allocatable, dimension(:) :: vm1e
591  real(DP), allocatable, dimension(:) :: uprod
592  real(DP), allocatable, dimension(:) :: det
593  real(DP), allocatable, dimension(:) :: wt
594  real(DP), allocatable, dimension(:) :: bi0x
595  real(DP), allocatable, dimension(:) :: be0x
596  real(DP), allocatable, dimension(:) :: bi0y
597  real(DP), allocatable, dimension(:) :: be0y
598  real(DP), allocatable, dimension(:) :: bi1x
599  real(DP), allocatable, dimension(:) :: be1x
600  real(DP), allocatable, dimension(:) :: bi1y
601  real(DP), allocatable, dimension(:) :: be1y
602  real(DP), allocatable, dimension(:) :: be01x
603  real(DP), allocatable, dimension(:) :: be01y
604  real(DP) :: emxx
605  real(DP) :: emxy
606  real(DP) :: emyx
607  real(DP) :: emyy
608  real(DP) :: rx
609  real(DP) :: ry
610  real(DP) :: emdet
611  integer(I4B) :: i
612  integer(I4B) :: ip
613 
614  ! Allocate local arrays
615  allocate (vm0i(this%nverts))
616  allocate (vm0e(this%nverts))
617  allocate (vm1i(this%nverts))
618  allocate (vm1e(this%nverts))
619  allocate (uprod(this%nverts))
620  allocate (det(this%nverts))
621  allocate (wt(this%nverts))
622  allocate (bi0x(this%nverts))
623  allocate (be0x(this%nverts))
624  allocate (bi0y(this%nverts))
625  allocate (be0y(this%nverts))
626  allocate (bi1x(this%nverts))
627  allocate (be1x(this%nverts))
628  allocate (bi1y(this%nverts))
629  allocate (be1y(this%nverts))
630  allocate (be01x(this%nverts))
631  allocate (be01y(this%nverts))
632 
633  ! Set vm0i(1)
634  vm0i(1) = vm0ival
635 
636  ! Get remaining vm0i's sequentially using divergence conditions
637  do i = 2, this%nverts
638  ip = this%iprev(i)
639  vm0i(i) = (li(ip) * vm0i(ip) - le(ip) * this%vne(ip) &
640  + areasub(ip) * divcell) / li(i)
641  end do
642 
643  ! Get vm1i's from vm0i's using continuity conditions
644  vm1i = cshift(vm0i, 1)
645 
646  ! Get centroid velocity by setting up and solving 2x2 linear system
647  uprod = unintx * unextx + uninty * unexty
648  det = done - uprod * uprod
649  bi0x = (unintx - unextx * uprod) / det
650  be0x = (unextx - unintx * uprod) / det
651  bi0y = (uninty - unexty * uprod) / det
652  be0y = (unexty - uninty * uprod) / det
653  uprod = unintxnext * unextx + unintynext * unexty
654  det = done - uprod * uprod
655  bi1x = (unintxnext - unextx * uprod) / det
656  be1x = (unextx - unintxnext * uprod) / det
657  bi1y = (unintynext - unexty * uprod) / det
658  be1y = (unexty - unintynext * uprod) / det
659  be01x = 5.d-1 * (be0x + be1x)
660  be01y = 5.d-1 * (be0y + be1y)
661  wt = areasub / areacell
662  emxx = dtwo - sum(wt * be01x * unextx)
663  emxy = -sum(wt * be01x * unexty)
664  emyx = -sum(wt * be01y * unextx)
665  emyy = dtwo - sum(wt * be01y * unexty)
666  rx = sum(wt * (bi0x * vm0i + bi1x * vm1i + be01x * this%vne))
667  ry = sum(wt * (bi0y * vm0i + bi1y * vm1i + be01y * this%vne))
668  emdet = emxx * emyy - emxy * emyx
669  this%vctrx = (emyy * rx - emxy * ry) / emdet
670  this%vctry = (emxx * ry - emyx * rx) / emdet
671 
672  ! Get vm0e's using "known" conditions
673  vm0e = 5.d-1 * (this%vne + unextx * this%vctrx + unexty * this%vctry)
674 
675  ! Get vm1e's from uniformity along exterior edges
676  vm1e = vm0e
677 
678  ! Transform vm0 and vm1 to (x, y) coordinates
679  vm0x = bi0x * vm0i + be0x * vm0e
680  vm0y = bi0y * vm0i + be0y * vm0e
681  vm1x = bi1x * vm1i + be1x * vm0e
682  vm1y = bi1y * vm1i + be1y * vm0e
683 
684  ! Calculate head-cycle summation (which is proportional to
685  ! the curl of the head gradient)
686  hcsum = sum(lm * (kappax * (vm0x + vm1x) + kappay * (vm0y + vm1y)))
687 
688  ! Deallocate local arrays
689  deallocate (vm0i)
690  deallocate (vm0e)
691  deallocate (vm1i)
692  deallocate (vm1e)
693  deallocate (uprod)
694  deallocate (det)
695  deallocate (wt)
696  deallocate (bi0x)
697  deallocate (be0x)
698  deallocate (bi0y)
699  deallocate (be0y)
700  deallocate (bi1x)
701  deallocate (be1x)
702  deallocate (bi1y)
703  deallocate (be1y)
704  deallocate (be01x)
705  deallocate (be01y)
706 
707  end subroutine calc_thru_hcsum
708 
709 end module methodcellternarymodule
710 
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 dsix
real constant 6
Definition: Constants.f90:82
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
subroutine, public transform(xin, yin, zin, xout, yout, zout, xorigin, yorigin, zorigin, sinrot, cosrot, invert)
Apply a 3D translation and optional 2D rotation to coordinates.
Definition: GeomUtil.f90:183
real(dp) function, public area(xv, yv, cw)
Calculate polygon area, with vertices given in CW or CCW order.
Definition: GeomUtil.f90:370
logical function, public point_in_polygon(x, y, poly)
Check if a point is within a polygon.
Definition: GeomUtil.f90:27
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 the subdomain tracking method (submethod).
Definition: Method.f90:144
subroutine pass(this, particle)
Pass the particle to the next subdomain.
Definition: Method.f90:153
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:31
Particle tracked by the PRT model.
Definition: Particle.f90:32