dolfin team mailing list archive
-
dolfin team
-
Mailing list archive
-
Message #21782
Re: assemble of Matrix with Real spaces slow
On Thursday March 3 2011 11:20:03 Marie E. Rognes wrote:
> On 03/03/2011 08:03 PM, Johan Hake wrote:
> > 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.
>
> By replacing "CG", 1 by "R", 0 or?
OMG!! Yes, *flush*
That explains the memory usage :P
> > 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 can take a look on the FFC side, but not today.
Nice!
With the update correction from Marie the numbers now looks like:
With PETSc backend
Tensor without Mixed space | 0.11211 0.11211 1
With 1 global dofs | 1.9482 1.9482 1
With 2 global dofs | 2.8725 2.8725 1
With 4 global dofs | 5.1959 5.1959 1
With 8 global dofs | 10.524 10.524 1
With 16 global dofs | 25.574 25.574 1
With Epetra backend
Tensor without Mixed space | 0.87544 0.87544 1
With 1 global dofs | 1.7089 1.7089 1
With 2 global dofs | 2.6868 2.6868 1
With 4 global dofs | 4.28 4.28 1
With 8 global dofs | 8.123 8.123 1
With 16 global dofs | 17.394 17.394 1
Still a pretty big increase in time for just adding 16 scalar dofs to a system
of 274625 dofs in teh first place.
Johan
> --
> Marie
>
> > 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
> >
> >
> >
> > _______________________________________________
> > Mailing list: https://launchpad.net/~dolfin
> > Post to : dolfin@xxxxxxxxxxxxxxxxxxx
> > Unsubscribe : https://launchpad.net/~dolfin
> > More help : https://help.launchpad.net/ListHelp
>
> _______________________________________________
> Mailing list: https://launchpad.net/~dolfin
> Post to : dolfin@xxxxxxxxxxxxxxxxxxx
> Unsubscribe : https://launchpad.net/~dolfin
> More help : https://help.launchpad.net/ListHelp
# 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 = 64
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, "R", 0)]*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
References