On 6 September 2013 12:45, Jan Blechta <[email protected]> wrote:
> It seems that for P1*R function space DofMapBuilder::reorder_local
> scales quadratically with number of DOFs. Is it inevitable?
>

Yes. The 'Real' dof is coupled to all dofs, destroying graph sparsity.

The solutions are (a) extend DofMapBuilder to re-order just sub-dofmap
or (b) do not use Real for anything other than small problems (for
nullspace, they can be handled by an iterative solve). It doesn't play
well with linear solvers, it's bad for assembly and its bad in
parallel.

Garth


> =================================================================
> from dolfin import *
> import time
> import numpy as np
> import matplotlib.pyplot as plt
>
> #set_log_level(DEBUG)
> #parameters['reorder_dofs_serial'] = False
>
> for n in range(25, 33):
>
>     mesh = UnitCubeMesh(n, n, n)
>     #mesh = UnitSquareMesh(10*n, 10*n)
>     P1 = FunctionSpace(mesh, 'CG', 1)
>     R  = FunctionSpace(mesh, 'R', 0)
>
>     # Measure time for creation of P1*R space
>     t = time.time()
>     V = P1*R
>     t = time.time() - t
>
>     # Report
>     print 'DOFs:', V.dim(), 'TIME:', t, 'TIME/DOF:', t/V.dim()
>     plt.loglog(V.dim(), t, 'o')
>
> # Plot some quadratic function for comparison
> xl = plt.xlim()
> yl = plt.ylim()
> x  = np.linspace(*xl)
> y  = map(lambda x:yl[0]/(xl[0]*xl[0])*x*x, x)
> plt.loglog(x, y)
>
> # Show plot
> plt.show()
> =================================================================
>
> Jan
> _______________________________________________
> fenics mailing list
> [email protected]
> http://fenicsproject.org/mailman/listinfo/fenics
_______________________________________________
fenics mailing list
[email protected]
http://fenicsproject.org/mailman/listinfo/fenics

Reply via email to