[petsc-users] How much memory does TS use (in terms of # of matrices)

Andrew Spott andrew.spott at gmail.com
Wed May 16 01:01:58 CDT 2012


I'm attempting to run some code with a rather large (>1GB) sparse matrix, and I keep getting kicked out of the que for what I believe is using too much memory.  The unfortunate thing is that as the steps go up, the memory usage seems to as well, which would normally indicate a memory leak.

But I'm not really doing anything that would lead to a memory leak, my jacobian is:

PetscErrorCode HamiltonianJ(TS timeStepContext, PetscReal t, Vec u, Mat *A, Mat *B, MatStructure *flag, void* context)
{
    Mat h;
    PetscErrorCode  err;
    Context         *cntx = (Context*)context;
    PetscScalar     ef = eField(t, cntx->params);
    //if (cntx->rank == 0) std::cout << "ef: " << ef;
    
    
    err = MatDuplicate(cntx->energy_eigenstates,MAT_COPY_VALUES,&h);
    err = MatAXPY(h, ef, cntx->dipole_matrix, DIFFERENT_NONZERO_PATTERN);

	err = MatAssemblyBegin(h, MAT_FINAL_ASSEMBLY);
	err = MatAssemblyEnd(h, MAT_FINAL_ASSEMBLY);
    
	*A = h;
    *flag = DIFFERENT_NONZERO_PATTERN;

    return err;
}

And my Monitor function is:

PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal t, Vec u, void *cntx)
{
    Context         *context = (Context*) cntx;
    PetscErrorCode  ierr = 0;
    PetscScalar     ef = eField(t, context->params);
    
    if (context->rank == 0 && context->progress) std::cout << step << ": " << t << " ef: " << ef << std::endl;
    
	return ierr;
}

I looked at it with valgrind, and got a summary of:

==24941== LEAK SUMMARY:
==24941==    definitely lost: 0 bytes in 0 blocks
==24941==    indirectly lost: 0 bytes in 0 blocks
==24941==      possibly lost: 1,291,626,491 bytes in 4,050 blocks
==24941==    still reachable: 1,795,550 bytes in 88 blocks
==24941==         suppressed: 0 bytes in 0 blocks

with a bunch of things like this:

==24937== 757,113,160 bytes in 10 blocks are possibly lost in loss record 1,834 of 1,834
==24937==    at 0x4A05306: memalign (vg_replace_malloc.c:532)
==24937==    by 0x43083C: PetscMallocAlign(unsigned long, int, char const*, char const*, char const*
, void**) (mal.c:30)
==24937==    by 0x434144: PetscTrMallocDefault(unsigned long, int, char const*, char const*, char co
nst*, void**) (mtr.c:196)
==24937==    by 0x8641BF: MatSeqAIJSetPreallocation_SeqAIJ (aij.c:3550)
==24937==    by 0x861DBA: MatSeqAIJSetPreallocation(_p_Mat*, int, int const*) (aij.c:3493)
==24937==    by 0x9201B5: MatMPIAIJSetPreallocation_MPIAIJ (mpiaij.c:3209)
==24937==    by 0x922600: MatMPIAIJSetPreallocation(_p_Mat*, int, int const*, int, int const*) (mpiaij.c:3897)
==24937==    by 0x9315BA: MatAXPY_MPIAIJ(_p_Mat*, std::complex<double>, _p_Mat*, MatStructure) (mpiaij.c:2254)
==24937==    by 0x5A9828: MatAXPY(_p_Mat*, std::complex<double>, _p_Mat*, MatStructure) (axpy.c:39)
==24937==    by 0x405EF7: HamiltonianJ(_p_TS*, double, _p_Vec*, _p_Mat**, _p_Mat**, MatStructure*, void*) (in /home/becker/ansp6066/code/EnergyBasisSchrodingerSolver/bin/propagate)
==24937==    by 0x736441: TSComputeRHSJacobian(_p_TS*, double, _p_Vec*, _p_Mat**, _p_Mat**, MatStructure*) (ts.c:186)
==24937==    by 0x736D68: TSComputeRHSFunctionLinear(_p_TS*, double, _p_Vec*, _p_Vec*, void*) (ts.c:2541)

I imagine these are not actual memory errors, but I thought I'm not sure what is leaking, so I thought I would bring it to your attention.

I can attach the full valgrind output, but at 2mb, it is rather large, so I figured I would wait for someone to ask for it.

Any idea what the memory leak could be?  Am I doing something stupid?

Thanks for the help, as always.

-Andrew


More information about the petsc-users mailing list