Interesting problem. Memory profiling has always been challenging, but I've not
seen signs of leaks when I've looked. Your PySparse w/ GC curve seems to
reinforce that. There's a lot more going on with Trilinos, but I don't know why
it would be leaky, even with garbage collection. Out of curiosity, what are you
using to track memory consumption?
A couple of things that might help the non-gc case, but probably won't make
much difference with gc:
* I would declare your equation only once and then reset the value of psi,
tcoeff_e, energy_e, etc. as needed
* Overwriting psi with the residual of the sweep and then renormalizing seems
fraught.
Something like:
def normalize(psi):
a = 1.0/sqrt( inner_product(psi,psi) ) # normalizing factor#
return psi*a
psi.equation = (0 == DiffusionTerm(-tcoeff_e) + (V - energy_e) * psi )
psi_solver.iterations = 5
while abs(energy_old-energy_e) > cutoff:
res = psi.equation.sweep(var = psi, solver = psi_solver )
psi.updateOld()
psi.value = normalize(res)
energy_old = energy_e
energy_e = H_expect_e(psi)
psi.setValue( (alpha_mix*psi + (1.0-alpha_mix)*psi.old).value )
On Mar 21, 2016, at 2:14 PM, Michael Waters <[email protected]> wrote:
> Hello,
>
> I have a large memory leak while solving Schroedinger's equation through
> self-consistent iteration. I limit the number of solver iterations with my
> sweeps so that I can update the energy eigenvalues. I also mix my old and new
> solutions. Using both of these gives me very stable convergence. I've tried
> making a simple version of my code for demonstration:
>
>
> psi = CellVariable(name = 'Electron Wavefunction', mesh=mesh, value = 0.0,
> hasOld = True)
> V = CellVariable(name = 'External Potential', mesh=mesh, value = 0.0)
> tcoeff_e = CellVariable(name = 'Electron Kinetic Coeff', mesh=mesh, value =
> tcoeff)
>
> def inner_product(psi_a,psi_b):
> ans = ((psi_a * psi_b).cellVolumeAverage * simulation_volume).value
> return ans
>
> def normalize(psi):
> a = 1.0/sqrt( inner_product(psi,psi) ) # normalizing factor#
> psi.setValue( psi*a)
>
>
> def sandwich_H_e(psi_a, psi_b):
> Hpsi = (-tcoeff_e.arithmeticFaceValue* psi_b.faceGrad).divergence + V *
> psi_b
> return ((psi_a*Hpsi).cellVolumeAverage * simulation_volume).value
>
> def H_expect_e(psi):
> return sandwich_H_e(psi, psi)
>
>
> while abs(energy_old-energy_e) > cutoff:
> psi.equation = (0 == DiffusionTerm(-tcoeff_e) + (V - energy_e) *
> psi )
> psi_solver.iterations = 5
> psi = psi.equation.sweep(var = psi, solver = psi_solver )
>
> normalize(psi)
> energy_old = energy_e
> energy_e = H_expect_e(psi)
>
> psi.setValue( (alpha_mix*psi + (1.0-alpha_mix)*psi.old).value )
>
>
> My actual code is very long because there are many CellVariables that sum
> together to form V. As for the solvers, I am using Pysparse's LinearPCGSolver
> w/ SsorPreconditioner or Trilinos' LinearPCGSolver w/ JacobiPreconditioner. I
> tried adding gc.collect() statement at the end of each iteration in the
> while loop, which improves the behavior for both solvers. Here is a graph of
> the memory usage versus time for the solvers using the same system:
> <abbcifch..png>
> The long pauses in the cases without garbage collection correspond to writing
> to disk. With Trilinos, I find the leaks to behave the same in parallel as
> serial, so I am just showing serial results. Looking at the results of using
> Pysparse, it seems likely that my code is cause of the leak there.
>
> I am new to this kind of debugging, so any help would be greatly appreciated.
>
> Thanks,
> -Michael Waters
>
>
> _______________________________________________
> fipy mailing list
> [email protected]
> http://www.ctcms.nist.gov/fipy
> [ NIST internal ONLY: https://email.nist.gov/mailman/listinfo/fipy ]
_______________________________________________
fipy mailing list
[email protected]
http://www.ctcms.nist.gov/fipy
[ NIST internal ONLY: https://email.nist.gov/mailman/listinfo/fipy ]