Hi,

I am currently comparing a fenics based solver to some other open source 
Navier-Stokes solvers and I just realized that I was using much more memory, 
perhaps as much as a factor 3. I always thought it was the matrices, vectors 
and Krylov solvers that were consuming memory, but on inspection I see that 
memory is totally dominated by the boundary conditions (mesh connectivity?). 
This was quite surprising to me, so I guess my question is - is this a well 
known issue?

I have made a small test for memory use of a Poisson solver to illustrate. 
Running it on one cpu returns the table showing the increment in RSS memory use 
and total memory use for each new function down the script. 

In [1]: run memorytest.py
                     RSS memory Increment       Total
Starting weight of dolfin        75988 KB      75988 KB
Mesh                             16760 KB      92748 KB
FunctionSpace                    41600 KB     134348 KB
Matrix assembly                   7792 KB     142140 KB
Function                           304 KB     142444 KB
Solve                              868 KB     143312 KB
Boundary condition              321428 KB     464740 KB
New solve                            0 KB     464740 KB

Creating and applying the Dirichlet boundary condition calls mesh.init(D) and 
mesh.init(D-1, D). I think this is were the extra memory comes in. The memory 
cost is nearly 5 times more than for all the other stuff! Any comments? 

Mikael



from dolfin import *
from commands import getoutput
from os import getpid
set_log_active(False)

def getRSSMemoryUsage():
    mypid = getpid()
    return eval(getoutput("ps -o rss %s" % mypid).split()[1])

class MemoryUsage:
    def __init__(self, s):
        self.memory = 0
        info_blue(' '*21+'RSS memory Increment       Total')
        self(s)
        
    def __call__(self, s):
        self.prev = self.memory
        self.memory = getRSSMemoryUsage()
        #self.memory = memory_usage(False)[0]
        info_blue('{0:26s}  {1:10d} KB {2:10d} KB'.format(s, 
                   self.memory-self.prev, self.memory))

memoryusage = MemoryUsage('Starting weight of dolfin')

mesh = UnitCubeMesh(40, 40, 40)

memoryusage("Mesh")

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

memoryusage('FunctionSpace')

u, v = TrialFunction(V), TestFunction(V)
A = assemble(inner(grad(u), grad(v))*dx)

memoryusage("Matrix assembly")

b = assemble(Constant(6)*v*dx)
u = Function(V)

memoryusage("Function")

solve(A, u.vector(), b, 'gmres', 'hypre_amg')

memoryusage('Solve')

bc = DirichletBC(V, Constant(0), DomainBoundary())
bc.apply(A, b)

memoryusage('Boundary condition')

solve(A, u.vector(), b, 'gmres', 'hypre_amg')

memoryusage('New solve')

_______________________________________________
fenics mailing list
[email protected]
http://fenicsproject.org/mailman/listinfo/fenics

Reply via email to