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