Skip to content

Commit

Permalink
- fix mem leak from KSPBuildResidual in convergence check (#1435)
Browse files Browse the repository at this point in the history
  • Loading branch information
mjr-deltares authored Nov 9, 2023
1 parent 1b220d4 commit 4921c1b
Show file tree
Hide file tree
Showing 3 changed files with 30 additions and 4 deletions.
4 changes: 4 additions & 0 deletions src/Solution/NumericalSolution.f90
Original file line number Diff line number Diff line change
Expand Up @@ -1222,6 +1222,10 @@ subroutine sln_da(this)
call this%cnvg_summary%destroy()
deallocate (this%cnvg_summary)
!
! -- linear solver
call this%linear_solver%destroy()
deallocate (this%linear_solver)
!
! -- linear solver settings
call this%linear_settings%destroy()
deallocate (this%linear_settings)
Expand Down
25 changes: 25 additions & 0 deletions src/Solution/PETSc/PetscConvergence.F90
Original file line number Diff line number Diff line change
Expand Up @@ -18,6 +18,8 @@ module PetscConvergenceModule
real(DP) :: dvclose
integer(I4B) :: max_its
type(ConvergenceSummaryType), pointer :: cnvg_summary => null()
contains
procedure :: destroy
end type PetscContextType

! passing our context into PETSc requires an explicit interface
Expand Down Expand Up @@ -73,6 +75,8 @@ subroutine petsc_check_convergence(ksp, n, rnorm, flag, context, ierr)

summary => context%cnvg_summary

! NB: KSPBuildResidual needs to have its vector destroyed
! to avoid a memory leak, KSPBuildSolution doesn't...
call KSPBuildSolution(ksp, PETSC_NULL_VEC, x, ierr)
CHKERRQ(ierr)
call KSPBuildResidual(ksp, PETSC_NULL_VEC, PETSC_NULL_VEC, res, ierr)
Expand All @@ -88,6 +92,8 @@ subroutine petsc_check_convergence(ksp, n, rnorm, flag, context, ierr)
CHKERRQ(ierr)
call VecCopy(res, context%res_old, ierr)
CHKERRQ(ierr)
call VecDestroy(res, ierr)
CHKERRQ(ierr)
flag = KSP_CONVERGED_ITERATING
end if
return
Expand Down Expand Up @@ -122,6 +128,9 @@ subroutine petsc_check_convergence(ksp, n, rnorm, flag, context, ierr)
call VecCopy(res, context%res_old, ierr)
CHKERRQ(ierr)

call VecDestroy(res, ierr)
CHKERRQ(ierr)

! get dv and dr per local model
call VecGetArrayF90(context%delta_x, local_dx, ierr)
CHKERRQ(ierr)
Expand Down Expand Up @@ -171,4 +180,20 @@ subroutine petsc_check_convergence(ksp, n, rnorm, flag, context, ierr)

end subroutine petsc_check_convergence

subroutine destroy(this)
class(PetscContextType) :: this
! local
integer(I4B) :: ierr

call VecDestroy(this%x_old, ierr)
CHKERRQ(ierr)
call VecDestroy(this%res_old, ierr)
CHKERRQ(ierr)
call VecDestroy(this%delta_x, ierr)
CHKERRQ(ierr)
call VecDestroy(this%delta_res, ierr)
CHKERRQ(ierr)

end subroutine destroy

end module PetscConvergenceModule
5 changes: 1 addition & 4 deletions src/Solution/PETSc/PetscSolver.F90
Original file line number Diff line number Diff line change
Expand Up @@ -419,10 +419,7 @@ subroutine petsc_destroy(this)
deallocate (this%vec_residual)

! delete context
call VecDestroy(this%petsc_ctx%delta_x, ierr)
CHKERRQ(ierr)
call VecDestroy(this%petsc_ctx%x_old, ierr)
CHKERRQ(ierr)
call this%petsc_ctx%destroy()
deallocate (this%petsc_ctx)

call VecDestroy(this%pc_context%diag, ierr)
Expand Down

0 comments on commit 4921c1b

Please sign in to comment.