21 integer(I4B) :: nverts
22 real(dp),
allocatable,
dimension(:) :: xvert
23 real(dp),
allocatable,
dimension(:) :: yvert
24 real(dp),
allocatable,
dimension(:) :: vne
25 real(dp),
allocatable,
dimension(:) :: vv0x
26 real(dp),
allocatable,
dimension(:) :: vv0y
27 real(dp),
allocatable,
dimension(:) :: vv1x
28 real(dp),
allocatable,
dimension(:) :: vv1y
38 integer(I4B),
allocatable,
dimension(:) :: iprev
39 real(dp),
allocatable,
dimension(:) :: xvertnext
40 real(dp),
allocatable,
dimension(:) :: yvertnext
41 integer(I4B),
public,
pointer :: zeromethod
65 method%type => method%cell%type
66 method%delegates = .true.
68 method%subcell => subcell
74 deallocate (this%type)
78 subroutine load_mct(this, particle, next_level, submethod)
82 integer(I4B),
intent(in) :: next_level
83 class(
methodtype),
pointer,
intent(inout) :: submethod
85 select type (subcell => this%subcell)
87 call this%load_subcell(particle, subcell)
92 subcell=this%subcell, &
93 trackctl=this%trackctl, &
94 tracktimes=this%tracktimes)
105 integer(I4B) :: exitFace
106 integer(I4B) :: inface
108 exitface = particle%iboundary(3)
109 isc = particle%idomain(3)
111 select case (exitface)
118 if (inface .eq. 0) inface = this%nverts
122 if (isc .gt. this%nverts) isc = 1
123 particle%idomain(3) = isc
124 particle%iboundary(3) = 3
129 if (isc .lt. 1) isc = this%nverts
130 particle%idomain(3) = isc
131 particle%iboundary(3) = 2
135 inface = this%nverts + 2
138 inface = this%nverts + 3
140 if (inface .eq. -1)
then
141 particle%iboundary(2) = 0
142 else if (inface .eq. 0)
then
143 particle%iboundary(2) = 0
145 particle%iboundary(2) = inface
155 real(DP),
intent(in) :: tmax
160 call this%update(particle, this%cell%defn)
161 if (.not. particle%advancing)
return
165 if (particle%z > this%cell%defn%top)
then
166 particle%z = this%cell%defn%top
167 call this%save(particle, reason=1)
171 select type (cell => this%cell)
174 this%nverts = cell%defn%npolyverts
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)
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))
199 do i = 1, this%nverts
200 this%xvert(i) = cell%defn%polyvert(1, i)
201 this%yvert(i) = cell%defn%polyvert(2, i)
204 this%ztop = cell%defn%top
205 this%zbot = cell%defn%bot
206 this%dz = this%ztop - this%zbot
208 do i = 1, this%nverts
211 this%iprev = cshift(this%iprev, -1)
212 this%xvertnext = cshift(this%xvert, 1)
213 this%yvertnext = cshift(this%yvert, 1)
220 call this%track(particle, 2, tmax)
254 select type (cell => this%cell)
258 isc = particle%idomain(3)
262 do iv0 = 1, this%nverts
265 x1 = this%xvertnext(iv0)
266 y1 = this%yvertnext(iv0)
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
281 if ((alphai .ge.
dzero) .and. &
282 (alphai + betai .le.
done))
then
288 print *,
"error -- initial triangle not found in cell ", ic, &
289 " for particle at ", particle%x, particle%y, particle%z
293 particle%idomain(3) = isc
296 subcell%isubcell = 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)
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
315 subcell%vzbot = this%vzbot
316 subcell%vztop = this%vztop
329 real(DP),
allocatable,
dimension(:) :: xvals
330 real(DP),
allocatable,
dimension(:) :: yvals
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
345 real(DP),
allocatable,
dimension(:) :: areasub
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
362 select type (cell => this%cell)
366 allocate (le(this%nverts))
367 allocate (unextx(this%nverts))
368 allocate (unexty(this%nverts))
369 allocate (areasub(this%nverts))
370 allocate (li(this%nverts))
371 allocate (unintx(this%nverts))
372 allocate (uninty(this%nverts))
373 allocate (xmid(this%nverts))
374 allocate (ymid(this%nverts))
375 allocate (lm(this%nverts))
376 allocate (umx(this%nverts))
377 allocate (umy(this%nverts))
378 allocate (kappax(this%nverts))
379 allocate (kappay(this%nverts))
380 allocate (vm0x(this%nverts))
381 allocate (vm0y(this%nverts))
382 allocate (vm1x(this%nverts))
383 allocate (vm1y(this%nverts))
384 allocate (unextxnext(this%nverts))
385 allocate (unextynext(this%nverts))
386 allocate (wk1(this%nverts))
387 allocate (wk2(this%nverts))
392 wk1 = this%xvertnext - this%xvert
393 wk2 = this%yvertnext - this%yvert
394 le = dsqrt(wk1 * wk1 + wk2 * wk2)
399 areacell =
area(this%xvert, this%yvert)
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
408 do i = 1, this%nverts
409 xvals(1) = this%xvert(i)
410 xvals(2) = this%xvertnext(i)
412 yvals(1) = this%yvert(i)
413 yvals(2) = this%yvertnext(i)
415 areasub(i) =
area(xvals, yvals)
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)
425 divcell = sum(le * this%vne) / areacell
428 wk1 = this%xvert - this%xctr
429 wk2 = this%yvert - this%yctr
430 li = dsqrt(wk1 * wk1 + wk2 * wk2)
434 unextxnext = cshift(unintx, 1)
435 unextynext = cshift(uninty, 1)
438 xmid = 5.d-1 * (this%xvert + this%xctr)
439 ymid = 5.d-1 * (this%yvert + this%yctr)
442 wk1 = cshift(xmid, 1) - xmid
443 wk2 = cshift(ymid, 1) - ymid
444 lm = dsqrt(wk1 * wk1 + wk2 * wk2)
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)
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)
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)
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
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
511 deallocate (unextxnext)
512 deallocate (unextynext)
522 le, li, lm, areasub, areacell, &
523 unintx, uninty, unextx, unexty, &
524 unintxnext, unintynext, &
525 kappax, kappay, vm0x, vm0y, vm1x, vm1y, hcsum)
531 real(DP),
dimension(:) :: le
532 real(DP),
dimension(:) :: li
533 real(DP),
dimension(:) :: lm
534 real(DP),
dimension(:) :: areasub
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
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
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))
599 do i = 2, this%nverts
601 vm0i(i) = (li(ip) * vm0i(ip) - le(ip) * this%vne(ip) &
602 + areasub(ip) * divcell) / li(i)
606 vm1i = cshift(vm0i, 1)
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
635 vm0e = 5.d-1 * (this%vne + unextx * this%vctrx + unexty * this%vctry)
641 vm0x = bi0x * vm0i + be0x * vm0e
642 vm0y = bi0y * vm0i + be0y * vm0e
643 vm1x = bi1x * vm1i + be1x * vm0e
644 vm1y = bi1y * vm1i + be1y * vm0e
648 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 done
real constant 1
subroutine pstop(status, message)
Stop the program, optionally specifying an error status code.
real(dp) function, public area(xv, yv, cw)
Calculate polygon area, with vertices given in CW or CCW order.
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.
subroutine load(this, particle, next_level, submethod)
Load subdomain tracking method (submethod)
subroutine pass(this, particle)
Pass a particle to the next subdomain, internal use only.
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.