← Back to team overview

dolfin team mailing list archive

assemble of Matrix with Real spaces slow

 

Hello!

I am using Mixed spaces with Reals quite alot. It turnes out that assemble 
forms with functions from MixedFunctionSpaces containing Real spaces are dead 
slow. The time spent also increase with the number of included Real spaces, 
even if none of them are included in form which is assembled. 

The attached test script illustrates this.

The test script also reviels that an unproportial time is spent in FFC 
generating the code. This time also increase with the number of Real spaces 
included. Turning of FErari helped a bit with this point.

I have not profiled any of this, but I just throw it out there. I do not 
recognize any difference between for example Epetra or PETSc backend as 
suggested in the fixed bug for building of sparsity pattern with global dofs.

My test has been done on a DOLFIN 0.9.9+. I haven't profiled it yet.

Output from summary:
  
  Tensor without Mixed space  |       0.11401     0.11401     1
  With 1 global dofs          |       0.40725     0.40725     1
  With 2 global dofs          |       0.94694     0.94694     1
  With 4 global dofs          |         2.763       2.763     1
  With 8 global dofs          |        9.6149      9.6149     1

Also the amount of memory used to build the sparsity patter seams to double 
for each step. The memory print for a 32x32x32 unit cube with 16 global dofs 
was 1.6 GB memory(!?).

Johan
# Test assemble time for Mixed space with Reals 
from dolfin import *

parameters["optimize"] = True
parameters["optimize_form"] = True
parameters["form_compiler"]["optimize"] = True
parameters["form_compiler"]["cpp_optimize"] = True
#parameters["linear_algebra_backend"] = "Epetra"

size = 32

mesh = UnitCube(size, size, size)

print "Size:", mesh.num_vertices()

V = FunctionSpace(mesh, "CG", 1)
v = TestFunction(V)
u = TrialFunction(V)

parameters["form_compiler"]["representation"] = "tensor"
A = assemble(u*v*dx)

timer = Timer("Tensor without Mixed space")
A = assemble(u*v*dx, tensor=A, reset_sparsity=False)
timer.stop()

for i in [1, 2, 4, 8]:#, 16]:
    print "Number of Global dofs:", i
    V = MixedFunctionSpace([FunctionSpace(mesh, "CG", 1)]+\
                           [FunctionSpace(mesh, "CG", 1)]*i)
    v = TestFunction(V)
    u = TrialFunction(V)

    A = assemble(u[0]*v[0]*dx)

    timer = Timer("With %d global dofs"%i)
    A = assemble(u[0]*v[0]*dx, tensor=A, reset_sparsity=False)
    timer.stop()

summary()

Follow ups