Memory Leak in SNES Solver

Issue #1030 new
Moritz Huck created an issue

Hi, while using multiphenics I found a memory leak. Together with the developer of multiphenics (Francesco Ballarin) we could track the problem to dolfin. The following code produces growing memory usage: (The used dolfin version is 2018.2.dev.0)

from dolfin import *
class DirichletBoundary(SubDomain):
    def inside(self, x, on_boundary):
        return abs(x[0] - 1.0) < DOLFIN_EPS and on_boundary

mesh = UnitSquareMesh(32, 32)

V = FunctionSpace(mesh, "CG", 1)

g = Constant(1.0)
bc = DirichletBC(V, g, DirichletBoundary())

# Define variational problem
u = Function(V)
v = TestFunction(V)
f = Expression("x[0]*sin(x[1])", degree=2)
F = inner((1 + u**2)*grad(u), grad(v))*dx - f*v*dx


#Leaks
J=derivative(F,u)
P=NonlinearVariationalProblem(F,u,bc,J)
solver=NonlinearVariationalSolver(P)
solver.parameters.update({"nonlinear_solver": "snes","snes_solver": {"linear_solver": "mumps"}})
for i in range(10000):       
        solver.solve()
## Doenst leak
#for i in range(10000):
#       solve(F == 0, u, bc,solver_parameters={"nonlinear_solver": "snes","snes_solver": {"linear_solver": "mumps"})

dolfinnonlinear.png

If the code piece #Doenst leak is used the memory consumption stays constant.

Francesco found a valgrind report which might be related:

dolfin::PETScSNESSolver::init(dolfin::NonlinearProblem&, dolfin::GenericVector&) (PETScSNESSolver.cpp:184)

which might be a VecDuplicate without a corresponding VecDestroy.

Best Regards, Moritz

Comments (6)

  1. Francesco Ballarin

    Thanks Moritz for filing the issue we found.

    If our analysis of the valgrind output is correct, we can prepare a pull request that destroys the f_tmp vector in ~PETScSNESSolver(). Nonetheless, we are a bit puzzled by the fact that solving with solve(F == 0, ...) doesn't show the issue. What could be the reason?

    Thanks,

    Francesco & Moritz

  2. Lawrence Mitchell

    The _snes_ctx object has a destructor which calls VecDestroy(&f_tmp) (c.f. PetscSNESSolver.h).

    So when you do solve(F == 0, ...), then that duplicated vector is cleaned up.

    But, when you reuse the solver object, this->init is called every time before the solve method is called, and every time it duplicates the vector anew.

    So I suppose you could guard the VecDuplicate, and the subsequence error checking, by if(!_snes_ctx.f_tmp) { ... }

    But perhaps it should be the case that PETScSNESSolver::init is only ever called once, rather than once per solve.

  3. Francesco Ballarin

    Ah, thanks! I didn't notice that the destructor was implemented in the header file.

    Let's hear from (you and) others which option they prefer among the two you suggested, also considering that (if I understood correctly) PetscSNESSolver will not be part of dolfinx anymore.

  4. Log in to comment