23 integer(I4B) :: nverts
24 real(dp),
allocatable,
dimension(:) :: xvert
25 real(dp),
allocatable,
dimension(:) :: yvert
26 real(dp),
allocatable,
dimension(:) :: vne
27 real(dp),
allocatable,
dimension(:) :: vv0x
28 real(dp),
allocatable,
dimension(:) :: vv0y
29 real(dp),
allocatable,
dimension(:) :: vv1x
30 real(dp),
allocatable,
dimension(:) :: vv1y
40 integer(I4B),
allocatable,
dimension(:) :: iprev
41 real(dp),
allocatable,
dimension(:) :: xvertnext
42 real(dp),
allocatable,
dimension(:) :: yvertnext
43 integer(I4B),
public,
pointer :: zeromethod
49 procedure,
public :: pass =>
pass_mct
67 method%name => method%cell%type
68 method%delegates = .true.
70 method%subcell => subcell
76 deallocate (this%name)
80 subroutine load_mct(this, particle, next_level, submethod)
84 integer(I4B),
intent(in) :: next_level
85 class(
methodtype),
pointer,
intent(inout) :: submethod
87 select type (subcell => this%subcell)
89 call this%load_subcell(particle, subcell)
95 subcell=this%subcell, &
97 tracktimes=this%tracktimes)
108 integer(I4B) :: exitFace
109 integer(I4B) :: inface
114 select case (exitface)
121 if (inface .eq. 0) inface = this%nverts
125 if (isc .gt. this%nverts) isc = 1
132 if (isc .lt. 1) isc = this%nverts
138 inface = this%nverts + 2
141 inface = this%nverts + 3
144 if (inface .eq. -1)
then
146 else if (inface .eq. 0)
then
159 real(DP),
intent(in) :: tmax
162 real(DP) :: x, y, z, xO, yO
163 real(DP),
allocatable :: xs(:), ys(:)
166 select type (cell => this%cell)
168 call this%assess(particle, this%cell%defn, tmax)
169 if (.not. particle%advancing)
return
172 this%nverts = cell%defn%npolyverts
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)
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))
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))
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)
220 this%ztop = cell%defn%top
221 this%zbot = cell%defn%bot
222 this%dz = this%ztop - this%zbot
225 do i = 1, this%nverts
228 this%iprev = cshift(this%iprev, -1)
229 this%xvertnext = cshift(this%xvert, 1)
230 this%yvertnext = cshift(this%yvert, 1)
236 call particle%transform(xo, yo)
239 call this%track(particle, 2, tmax)
242 call particle%transform(xo, yo, invert=.true.)
243 call particle%reset_transform()
277 select type (cell => this%cell)
285 do iv0 = 1, this%nverts
288 x1 = this%xvertnext(iv0)
289 y1 = this%yvertnext(iv0)
296 di2 = xi * y2rel - yi * x2rel
297 d02 = x0 * y2rel - y0 * x2rel
298 d12 = x1rel * y2rel - y1rel * x2rel
299 di1 = xi * y1rel - yi * x1rel
300 d01 = x0 * y1rel - y0 * x1rel
301 alphai = (di2 - d02) / d12
302 betai = -(di1 - d01) / d12
304 if ((alphai .ge.
dzero) .and. &
305 (alphai + betai .le.
done))
then
311 print *,
"error -- initial triangle not found in cell ", ic, &
312 " for particle at ", particle%x, particle%y, particle%z
319 subcell%isubcell = isc
323 subcell%x0 = this%xvert(iv0)
324 subcell%y0 = this%yvert(iv0)
325 subcell%x1 = this%xvertnext(iv0)
326 subcell%y1 = this%yvertnext(iv0)
327 subcell%x2 = this%xctr
328 subcell%y2 = this%yctr
329 subcell%v0x = this%vv0x(iv0)
330 subcell%v0y = this%vv0y(iv0)
331 subcell%v1x = this%vv1x(iv0)
332 subcell%v1y = this%vv1y(iv0)
333 subcell%v2x = this%vctrx
334 subcell%v2y = this%vctry
335 subcell%ztop = this%ztop
336 subcell%zbot = this%zbot
338 subcell%vzbot = this%vzbot
339 subcell%vztop = this%vztop
352 real(DP),
allocatable,
dimension(:) :: xvals
353 real(DP),
allocatable,
dimension(:) :: yvals
360 real(DP),
allocatable,
dimension(:) :: wk1
361 real(DP),
allocatable,
dimension(:) :: wk2
362 real(DP),
allocatable,
dimension(:) :: unextxnext
363 real(DP),
allocatable,
dimension(:) :: unextynext
364 real(DP),
allocatable,
dimension(:) :: le
365 real(DP),
allocatable,
dimension(:) :: unextx
366 real(DP),
allocatable,
dimension(:) :: unexty
368 real(DP),
allocatable,
dimension(:) :: areasub
370 real(DP),
allocatable,
dimension(:) :: li
371 real(DP),
allocatable,
dimension(:) :: unintx
372 real(DP),
allocatable,
dimension(:) :: uninty
373 real(DP),
allocatable,
dimension(:) :: xmid
374 real(DP),
allocatable,
dimension(:) :: ymid
375 real(DP),
allocatable,
dimension(:) :: lm
376 real(DP),
allocatable,
dimension(:) :: umx
377 real(DP),
allocatable,
dimension(:) :: umy
378 real(DP),
allocatable,
dimension(:) :: kappax
379 real(DP),
allocatable,
dimension(:) :: kappay
380 real(DP),
allocatable,
dimension(:) :: vm0x
381 real(DP),
allocatable,
dimension(:) :: vm0y
382 real(DP),
allocatable,
dimension(:) :: vm1x
383 real(DP),
allocatable,
dimension(:) :: vm1y
385 integer(I4B) :: nvert
387 select type (cell => this%cell)
391 allocate (le(this%nverts))
392 allocate (unextx(this%nverts))
393 allocate (unexty(this%nverts))
394 allocate (areasub(this%nverts))
395 allocate (li(this%nverts))
396 allocate (unintx(this%nverts))
397 allocate (uninty(this%nverts))
398 allocate (xmid(this%nverts))
399 allocate (ymid(this%nverts))
400 allocate (lm(this%nverts))
401 allocate (umx(this%nverts))
402 allocate (umy(this%nverts))
403 allocate (kappax(this%nverts))
404 allocate (kappay(this%nverts))
405 allocate (vm0x(this%nverts))
406 allocate (vm0y(this%nverts))
407 allocate (vm1x(this%nverts))
408 allocate (vm1y(this%nverts))
409 allocate (unextxnext(this%nverts))
410 allocate (unextynext(this%nverts))
411 allocate (wk1(this%nverts))
412 allocate (wk2(this%nverts))
417 wk1 = this%xvertnext - this%xvert
418 wk2 = this%yvertnext - this%yvert
419 le = dsqrt(wk1 * wk1 + wk2 * wk2)
424 areacell = area(this%xvert, this%yvert)
425 sixa = areacell * dsix
426 wk1 = -(this%xvert * this%yvertnext - this%xvertnext * this%yvert)
427 nvert =
size(this%xvert)
428 this%xctr = sum((this%xvert + this%xvertnext) * wk1) / sixa
429 this%yctr = sum((this%yvert + this%yvertnext) * wk1) / sixa
435 do i = 1, this%nverts
436 xvals(1) = this%xvert(i)
437 xvals(2) = this%xvertnext(i)
439 yvals(1) = this%yvert(i)
440 yvals(2) = this%yvertnext(i)
442 areasub(i) = area(xvals, yvals)
446 term =
done / (cell%defn%porosity * cell%defn%retfactor * this%dz)
447 do i = 1, this%nverts
448 this%vne(i) = -cell%defn%faceflow(i) * term / le(i)
452 divcell = sum(le * this%vne) / areacell
455 wk1 = this%xvert - this%xctr
456 wk2 = this%yvert - this%yctr
457 li = dsqrt(wk1 * wk1 + wk2 * wk2)
461 unextxnext = cshift(unintx, 1)
462 unextynext = cshift(uninty, 1)
465 xmid = 5.d-1 * (this%xvert + this%xctr)
466 ymid = 5.d-1 * (this%yvert + this%yctr)
469 wk1 = cshift(xmid, 1) - xmid
470 wk2 = cshift(ymid, 1) - ymid
471 lm = dsqrt(wk1 * wk1 + wk2 * wk2)
488 call this%calc_thru_hcsum(vm0i0, divcell, le, li, lm, areasub, areacell, &
489 unintx, uninty, unextx, unexty, &
490 unextxnext, unextynext, &
491 kappax, kappay, vm0x, vm0y, vm1x, vm1y, hcsum0)
493 vm0ival = vm0i0 + perturb
494 call this%calc_thru_hcsum(vm0ival, divcell, le, li, lm, areasub, areacell, &
495 unintx, uninty, unextx, unexty, &
496 unextxnext, unextynext, &
497 kappax, kappay, vm0x, vm0y, vm1x, vm1y, hcsum)
499 jac = (hcsum - hcsum0) / perturb
500 vm0ival = vm0i0 - hcsum0 / jac
501 call this%calc_thru_hcsum(vm0ival, divcell, le, li, lm, areasub, areacell, &
502 unintx, uninty, unextx, unexty, &
503 unextxnext, unextynext, &
504 kappax, kappay, vm0x, vm0y, vm1x, vm1y, hcsum)
509 this%vv0x = 2.d0 * vm0x - this%vctrx
510 this%vv0y = 2.d0 * vm0y - this%vctry
511 this%vv1x = 2.d0 * vm1x - this%vctrx
512 this%vv1y = 2.d0 * vm1y - this%vctry
515 term =
done / (cell%defn%retfactor * cell%defn%porosity * areacell)
516 this%vzbot = cell%defn%faceflow(this%nverts + 2) * term
517 this%vztop = -cell%defn%faceflow(this%nverts + 3) * term
538 deallocate (unextxnext)
539 deallocate (unextynext)
549 le, li, lm, areasub, areacell, &
550 unintx, uninty, unextx, unexty, &
551 unintxnext, unintynext, &
552 kappax, kappay, vm0x, vm0y, vm1x, vm1y, hcsum)
558 real(DP),
dimension(:) :: le
559 real(DP),
dimension(:) :: li
560 real(DP),
dimension(:) :: lm
561 real(DP),
dimension(:) :: areasub
563 real(DP),
dimension(:) :: unintx
564 real(DP),
dimension(:) :: uninty
565 real(DP),
dimension(:) :: unextx
566 real(DP),
dimension(:) :: unexty
567 real(DP),
dimension(:) :: unintxnext
568 real(DP),
dimension(:) :: unintynext
569 real(DP),
dimension(:) :: kappax
570 real(DP),
dimension(:) :: kappay
571 real(DP),
dimension(:) :: vm0x
572 real(DP),
dimension(:) :: vm0y
573 real(DP),
dimension(:) :: vm1x
574 real(DP),
dimension(:) :: vm1y
576 real(DP),
allocatable,
dimension(:) :: vm0i
577 real(DP),
allocatable,
dimension(:) :: vm0e
578 real(DP),
allocatable,
dimension(:) :: vm1i
579 real(DP),
allocatable,
dimension(:) :: vm1e
580 real(DP),
allocatable,
dimension(:) :: uprod
581 real(DP),
allocatable,
dimension(:) :: det
582 real(DP),
allocatable,
dimension(:) :: wt
583 real(DP),
allocatable,
dimension(:) :: bi0x
584 real(DP),
allocatable,
dimension(:) :: be0x
585 real(DP),
allocatable,
dimension(:) :: bi0y
586 real(DP),
allocatable,
dimension(:) :: be0y
587 real(DP),
allocatable,
dimension(:) :: bi1x
588 real(DP),
allocatable,
dimension(:) :: be1x
589 real(DP),
allocatable,
dimension(:) :: bi1y
590 real(DP),
allocatable,
dimension(:) :: be1y
591 real(DP),
allocatable,
dimension(:) :: be01x
592 real(DP),
allocatable,
dimension(:) :: be01y
604 allocate (vm0i(this%nverts))
605 allocate (vm0e(this%nverts))
606 allocate (vm1i(this%nverts))
607 allocate (vm1e(this%nverts))
608 allocate (uprod(this%nverts))
609 allocate (det(this%nverts))
610 allocate (wt(this%nverts))
611 allocate (bi0x(this%nverts))
612 allocate (be0x(this%nverts))
613 allocate (bi0y(this%nverts))
614 allocate (be0y(this%nverts))
615 allocate (bi1x(this%nverts))
616 allocate (be1x(this%nverts))
617 allocate (bi1y(this%nverts))
618 allocate (be1y(this%nverts))
619 allocate (be01x(this%nverts))
620 allocate (be01y(this%nverts))
626 do i = 2, this%nverts
628 vm0i(i) = (li(ip) * vm0i(ip) - le(ip) * this%vne(ip) &
629 + areasub(ip) * divcell) / li(i)
633 vm1i = cshift(vm0i, 1)
636 uprod = unintx * unextx + uninty * unexty
637 det =
done - uprod * uprod
638 bi0x = (unintx - unextx * uprod) / det
639 be0x = (unextx - unintx * uprod) / det
640 bi0y = (uninty - unexty * uprod) / det
641 be0y = (unexty - uninty * uprod) / det
642 uprod = unintxnext * unextx + unintynext * unexty
643 det =
done - uprod * uprod
644 bi1x = (unintxnext - unextx * uprod) / det
645 be1x = (unextx - unintxnext * uprod) / det
646 bi1y = (unintynext - unexty * uprod) / det
647 be1y = (unexty - unintynext * uprod) / det
648 be01x = 5.d-1 * (be0x + be1x)
649 be01y = 5.d-1 * (be0y + be1y)
650 wt = areasub / areacell
651 emxx =
dtwo - sum(wt * be01x * unextx)
652 emxy = -sum(wt * be01x * unexty)
653 emyx = -sum(wt * be01y * unextx)
654 emyy =
dtwo - sum(wt * be01y * unexty)
655 rx = sum(wt * (bi0x * vm0i + bi1x * vm1i + be01x * this%vne))
656 ry = sum(wt * (bi0y * vm0i + bi1y * vm1i + be01y * this%vne))
657 emdet = emxx * emyy - emxy * emyx
658 this%vctrx = (emyy * rx - emxy * ry) / emdet
659 this%vctry = (emxx * ry - emyx * rx) / emdet
662 vm0e = 5.d-1 * (this%vne + unextx * this%vctrx + unexty * this%vctry)
668 vm0x = bi0x * vm0i + be0x * vm0e
669 vm0y = bi0y * vm0i + be0y * vm0e
670 vm1x = bi1x * vm1i + be1x * vm0e
671 vm1y = bi1y * vm1i + be1y * vm0e
675 hcsum = sum(lm * (kappax * (vm0x + vm1x) + kappay * (vm0y + vm1y)))
subroutine, public create_cell_poly(cell)
Create a new polygonal cell.
This module contains simulation constants.
real(dp), parameter dhalf
real constant 1/2
real(dp), parameter dzero
real constant zero
real(dp), parameter dtwo
real constant 2
real(dp), parameter dsix
real constant 6
real(dp), parameter done
real constant 1
subroutine pstop(status, message)
Stop the program, optionally specifying an error status code.
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.
real(dp) function, public area(xv, yv, cw)
Calculate polygon area, with vertices given in CW or CCW order.
logical function, public point_in_polygon(x, y, poly)
Check if a point is within a polygon.
This module defines variable data types.
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.
@, public level_subfeature
Subcell-level tracking methods.
type(methodsubcellternarytype), pointer, public method_subcell_tern
subroutine, public create_subcell_tri(subcell)
Create a new triangular subcell.
Base type for particle tracking methods.
Particle tracked by the PRT model.