Skip to content

Commit

Permalink
fix(prt): avoid fpes, patch centroid calculation (#2119)
Browse files Browse the repository at this point in the history
Fix several related issues that can occur with the generalized tracking method on vertex grids.

1. The polygon centroid formula used when dividing cells into subcells can fail due to floating (im)precision when the vertex coordinates are all large and near to one another, i.e. big grids with really small cells. One option is to do a point-in-polygon check on the result of the centroid calculation and use a simple geometric mean if it fails, but the check is not cheap. Another option is to always use the geo mean instead of the polygon centroid — to guarantee that this is safe we will need to check that the cell vertices don't "wrap around" or contain duplicates I think. In the meantime this PR aims for an 80/20 solution and uses the geometric mean if the cell has 3 vertices, otherwise the full formula. This avoids the error in the common case that the very small cell is a triangle.

2. The tracking routine could encounter floating point exceptions due to a failure to terminate particles in subcells with no exit face.

3. Avoid floating point exceptions due to log of a negative number in the subroutine to calculate travel time to exit. This is done by taking the absolute value of the argument to `log`. Is this safe? I think so because the direction of travel is determined elsewhere.
  • Loading branch information
wpbonelli authored Dec 23, 2024
1 parent 298a4b2 commit 3f407b3
Show file tree
Hide file tree
Showing 3 changed files with 48 additions and 25 deletions.
37 changes: 29 additions & 8 deletions src/Solution/ParticleTracker/MethodCellTernary.f90
Original file line number Diff line number Diff line change
@@ -1,6 +1,6 @@
module MethodCellTernaryModule

use KindModule, only: DP, I4B
use KindModule, only: DP, I4B, LGP
use ErrorUtilModule, only: pstop
use MethodModule
use MethodSubcellPoolModule
Expand All @@ -9,7 +9,8 @@ module MethodCellTernaryModule
use SubcellTriModule, only: SubcellTriType, create_subcell_tri
use ParticleModule
use GeomUtilModule, only: area
use ConstantsModule, only: DZERO, DONE, DTWO
use ConstantsModule, only: DZERO, DONE, DTWO, DSIX
use GeomUtilModule, only: point_in_polygon
implicit none

private
Expand Down Expand Up @@ -86,6 +87,7 @@ subroutine load_mct(this, particle, next_level, submethod)
type is (SubcellTriType)
call this%load_subcell(particle, subcell)
end select

call method_subcell_tern%init( &
fmi=this%fmi, &
cell=this%cell, &
Expand Down Expand Up @@ -137,6 +139,7 @@ subroutine pass_mct(this, particle)
! Subcell top (cell top)
inface = this%nverts + 3
end select

if (inface .eq. -1) then
particle%iboundary(2) = 0
else if (inface .eq. 0) then
Expand Down Expand Up @@ -356,6 +359,8 @@ subroutine vertvelo(this)
real(DP), allocatable, dimension(:) :: vm0y
real(DP), allocatable, dimension(:) :: vm1x
real(DP), allocatable, dimension(:) :: vm1y
! real(DP), allocatable, dimension(:, :) :: poly
integer(I4B) :: nvert

select type (cell => this%cell)
type is (CellPolyType)
Expand Down Expand Up @@ -393,14 +398,30 @@ subroutine vertvelo(this)
unextx = -wk2 / le
unexty = wk1 / le

! Cell area
areacell = area(this%xvert, this%yvert)

! Cell centroid (in general, this is NOT the average of the vertex coordinates)
sixa = areacell * 6.d0
areacell = area(this%xvert, this%yvert)
sixa = areacell * DSIX
wk1 = -(this%xvert * this%yvertnext - this%xvertnext * this%yvert)
this%xctr = sum((this%xvert + this%xvertnext) * wk1) / sixa
this%yctr = sum((this%yvert + this%yvertnext) * wk1) / sixa
nvert = size(this%xvert)
if (nvert == 3) then
this%xctr = sum(this%xvert) / 3.0_DP
this%yctr = sum(this%yvert) / 3.0_DP
else
this%xctr = sum((this%xvert + this%xvertnext) * wk1) / sixa
this%yctr = sum((this%yvert + this%yvertnext) * wk1) / sixa
end if

! TODO: do we need to check if center is in polygon?
! only if using the centroid method which misbehaves
!
! allocate(poly(2, nvert))
! poly(1,:) = this%xvert
! poly(2,:) = this%yvert
! if (.not. point_in_polygon(this%xctr, this%yctr, poly)) then
! print *, "error -- centroid not in cell ", this%cell%defn%icell, this%xctr, this%yctr
! call pstop(1)
! end if
! deallocate(poly)

! Subcell areas
do i = 1, this%nverts
Expand Down
29 changes: 18 additions & 11 deletions src/Solution/ParticleTracker/MethodSubcellTernary.f90
Original file line number Diff line number Diff line change
Expand Up @@ -11,6 +11,7 @@ module MethodSubcellTernaryModule
use TernarySolveTrack, only: traverse_triangle, step_analytical, canonical
use PrtFmiModule, only: PrtFmiType
use BaseDisModule, only: DisBaseType
use MathUtilModule, only: is_close
implicit none

private
Expand All @@ -30,9 +31,9 @@ module MethodSubcellTernaryModule

!> @brief Create a new ternary subcell tracking method.
subroutine create_method_subcell_ternary(method)
! -- dummy
! dummy
type(MethodSubcellTernaryType), pointer :: method
! -- local
! local
type(SubcellTriType), pointer :: subcell

allocate (method)
Expand Down Expand Up @@ -218,14 +219,19 @@ subroutine track_subcell(this, subcell, particle, tmax)
! Exit through lateral subcell face
exitFace = itrifaceexit
dtexit = dtexitxy
else if (dtexitz < dtexitxy) then
else if (dtexitz < dtexitxy .and. dtexitz >= 0.0_DP) then
! Exit through top or bottom
if (itopbotexit == -1) then
exitFace = 4
else
exitFace = 5
end if
dtexit = dtexitz
else
particle%istatus = 9
particle%advancing = .false.
call this%save(particle, reason=3)
return
end if
texit = particle%ttrack + dtexit
t0 = particle%ttrack
Expand Down Expand Up @@ -257,17 +263,17 @@ subroutine track_subcell(this, subcell, particle, tmax)
! (local, unscaled) and other properties. The particle may at this
! point lie on a boundary of the subcell or may still be within it.
if (texit .gt. tmax) then
! -- The computed exit time is greater than the maximum time, so set
! -- final time for particle trajectory equal to maximum time.
! The computed exit time is greater than the maximum time, so set
! final time for particle trajectory equal to maximum time.
t = tmax
dt = t - t0
exitFace = 0
particle%istatus = 1
particle%advancing = .false.
reason = 2 ! timestep end
else
! -- The computed exit time is less than or equal to the maximum time,
! -- so set final time for particle trajectory equal to exit time.
! The computed exit time is less than or equal to the maximum time,
! so set final time for particle trajectory equal to exit time.
t = texit
dt = dtexit
reason = 1 ! (sub)cell transition
Expand All @@ -280,6 +286,7 @@ subroutine track_subcell(this, subcell, particle, tmax)
particle%z = z
particle%ttrack = t
particle%iboundary(3) = exitFace

call this%save(particle, reason=reason)
end subroutine track_subcell

Expand Down Expand Up @@ -417,7 +424,7 @@ subroutine calculate_dt(v1, v2, dx, xL, v, dvdx, &
end if

! Compute travel time to exit face. Return with status = 0
dt = log(vr) / dvdx
dt = log(abs(vr)) / dvdx
status = 0
end subroutine calculate_dt

Expand Down Expand Up @@ -481,13 +488,13 @@ subroutine calculate_xyz_position(dt, rxx, rxy, ryx, ryy, sxx, sxy, syy, &
else
! otherwise calculate z.
if (izstatus .eq. 2) then
! -- vz uniformly zero
! vz uniformly zero
z = zi
else if (izstatus .eq. 1) then
! -- vz uniform, nonzero
! vz uniform, nonzero
z = zi + vzi * dt
else
! -- vz nonuniform
! vz nonuniform
z = zbot + (vzi * dexp(az * dt) - vzbot) / az
end if
end if
Expand Down
7 changes: 1 addition & 6 deletions src/Solution/ParticleTracker/Particle.f90
Original file line number Diff line number Diff line change
Expand Up @@ -89,7 +89,6 @@ module ParticleModule
! state
integer(I4B), dimension(:, :), pointer, public, contiguous :: idomain !< array of indices for domains in the tracking domain hierarchy
integer(I4B), dimension(:, :), pointer, public, contiguous :: iboundary !< array of indices for tracking domain boundaries
integer(I4B), dimension(:), pointer, public, contiguous :: icp !< previous cell number (reduced)
integer(I4B), dimension(:), pointer, public, contiguous :: icu !< cell number (user)
integer(I4B), dimension(:), pointer, public, contiguous :: ilay !< layer
integer(I4B), dimension(:), pointer, public, contiguous :: izone !< current zone number
Expand Down Expand Up @@ -133,7 +132,6 @@ subroutine allocate_particle_store(this, np, mempath)
call mem_allocate(this%irpt, np, 'PLIRPT', mempath)
call mem_allocate(this%iprp, np, 'PLIPRP', mempath)
call mem_allocate(this%name, LENBOUNDNAME, np, 'PLNAME', mempath)
call mem_allocate(this%icp, np, 'PLICP', mempath)
call mem_allocate(this%icu, np, 'PLICU', mempath)
call mem_allocate(this%ilay, np, 'PLILAY', mempath)
call mem_allocate(this%izone, np, 'PLIZONE', mempath)
Expand Down Expand Up @@ -165,7 +163,6 @@ subroutine deallocate (this, mempath)
call mem_deallocate(this%iprp, 'PLIPRP', mempath)
call mem_deallocate(this%irpt, 'PLIRPT', mempath)
call mem_deallocate(this%name, 'PLNAME', mempath)
call mem_deallocate(this%icp, 'PLICP', mempath)
call mem_deallocate(this%icu, 'PLICU', mempath)
call mem_deallocate(this%ilay, 'PLILAY', mempath)
call mem_deallocate(this%izone, 'PLIZONE', mempath)
Expand Down Expand Up @@ -200,7 +197,6 @@ subroutine resize(this, np, mempath)
call mem_reallocate(this%iprp, np, 'PLIPRP', mempath)
call mem_reallocate(this%irpt, np, 'PLIRPT', mempath)
call mem_reallocate(this%name, LENBOUNDNAME, np, 'PLNAME', mempath)
call mem_reallocate(this%icp, np, 'PLICP', mempath)
call mem_reallocate(this%icu, np, 'PLICU', mempath)
call mem_reallocate(this%ilay, np, 'PLILAY', mempath)
call mem_reallocate(this%izone, np, 'PLIZONE', mempath)
Expand Down Expand Up @@ -244,7 +240,7 @@ subroutine load_particle(this, store, imdl, iprp, ip)
this%istopweaksink = store%istopweaksink(ip)
this%istopzone = store%istopzone(ip)
this%idrymeth = store%idrymeth(ip)
this%icp = store%icp(ip)
this%icp = 0
this%icu = store%icu(ip)
this%ilay = store%ilay(ip)
this%izone = store%izone(ip)
Expand Down Expand Up @@ -281,7 +277,6 @@ subroutine save_particle(this, particle, ip)
this%istopweaksink(ip) = particle%istopweaksink
this%istopzone(ip) = particle%istopzone
this%idrymeth(ip) = particle%idrymeth
this%icp(ip) = particle%icp
this%icu(ip) = particle%icu
this%ilay(ip) = particle%ilay
this%izone(ip) = particle%izone
Expand Down

0 comments on commit 3f407b3

Please sign in to comment.