[petsc-users] Memory leak at GPU when updating matrix of type mpiaijcusparse (CUDA)
jordic
jordic at cttc.upc.edu
Fri Feb 28 04:29:44 CST 2020
Dear all,
the following simple program:
//////////////////////////////////////////////////////////////////////////////////////
#include <petscmat.h>
PetscInt ierr=0;
int main(int argc,char **argv)
{
MPI_Comm comm;
PetscMPIInt rank,size;
PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
comm = PETSC_COMM_WORLD;
MPI_Comm_rank(comm,&rank);
MPI_Comm_size(comm,&size);
Mat A;
MatCreate(comm, &A);
MatSetSizes(A, 1, 1, PETSC_DETERMINE, PETSC_DETERMINE);
MatSetFromOptions(A);
PetscInt dnz=1, onz=0;
MatMPIAIJSetPreallocation(A, 0, &dnz, 0, &onz);
MatSetOption(A, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE);
MatSetOption(A, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE);
PetscInt igid=rank, jgid=rank;
PetscScalar value=rank+1.0;
// for(int i=0; i<10; ++i)
for(;;) //infinite loop
{
MatSetValue(A, igid, jgid, value, INSERT_VALUES);
MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);
MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);
}
MatDestroy(&A);
PetscFinalize();
return ierr;
}
//////////////////////////////////////////////////////////////////////////////////////
creates a simple diagonal matrix with one value per mpi-core. If the
type of the matrix is "mpiaij" (-mat_type mpiaij) there is no problem
but with "mpiaijcusparse" (-mat_type mpiaijcusparse) the memory usage at
the GPU grows with every iteration of the infinite loop. The only
solution that I found is to destroy and create the matrix every time
that it needs to be updated. Is there a better way to avoid this
problem?
I am using Petsc Release Version 3.12.2 with this configure options:
Configure options --package-prefix-hash=/home_nobck/user/petsc-hash-pkgs
--with-debugging=0 --with-fc=0 CC=gcc CXX=g++ --COPTFLAGS="-g -O3"
--CXXOPTFLAGS="-g -O3" --CUDAOPTFLAGS="-D_FORCE_INLINES -g -O3"
--with-mpi-include=/usr/lib/openmpi/include
--with-mpi-lib="-L/usr/lib/openmpi/lib -lmpi_cxx -lmpi" --with-cuda=1
--with-precision=double --with-cuda-include=/usr/include
--with-cuda-lib="-L/usr/lib/x86_64-linux-gnu -lcuda -lcudart -lcublas
-lcufft -lcusparse -lcusolver"
PETSC_ARCH=arch-ci-linux-opt-cxx-cuda-double
Thanks for your help,
Jorge
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20200228/e52df577/attachment.html>
More information about the petsc-users
mailing list