[Trilinos-Users] Epetra_FECrsGraph question

Kent-Andre Mardal kent-and at simula.no
Wed Aug 13 06:20:49 MDT 2008


On ti., 2008-08-12 at 08:52 -0600, Williams, Alan B wrote:
> Hello,
> 
> When NumIndicesPerRow is not set at construction time, the crsgraph does memory allocations each time InsertGlobalIndices is called. That is probably the cause of your slow performance.
> 
> We have finite-element applications which use CrsGraph and CrsMatrix, and the typical execution sequence is to first perform an initialization step where the matrix-graph structure is computed and the CrsGraph and CrsMatrix objects are constructed, and then repeatedly fill the linear system and solve it. (In a typical application there are many fill/solve steps per initialization step.) This produces good performance.
> If the structure of the matrix changes, due to h-adaptivity or load-balancing etc., then a new initialization step must be performed before new fill/solve steps are performed.
> 
> Alan

Hello, 

But what I don't understand is that there is such a great difference in
the case where the finite element method contains vector type elements
in 3D. While in the scalar case, Epetra is comparable to the other
backends. 

Here are the numbers (the difference is mostly due to the EpetraGraph): 

Assemble total  |  Poisson2DP1  Poisson2DP2  Poisson2DP3  THStokes2D
StabStokes2D  Elasticity3D  NSEMomentum3D
---------------------------------------------------------------------------------------------------------------
uBLAS           |         0.14         1.53         2.79        5.13           1.2          3.31           3.35
PETSc           |         0.24         1.45         2.52        4.81           1.1          3.14           3.28
Epetra          |         0.21         1.44         2.42        4.33          1.28          7.05           7.26
MTL4            |         0.11         1.09         1.71        1.81          0.46          1.19           1.38
STL             |          0.1         1.19         2.27        3.85          0.73          1.86           2.09

Kent-Andre Mardal
Simula Research Laboratory

> 
> 
> > -----Original Message-----
> > From: trilinos-users-bounces at software.sandia.gov
> > [mailto:trilinos-users-bounces at software.sandia.gov] On Behalf
> > Of Kent-Andre Mardal
> > Sent: Monday, August 11, 2008 6:45 AM
> > To: trilinos-users at software.sandia.gov
> > Subject: [Trilinos-Users] Epetra_FECrsGraph question
> >
> > Hi,
> >
> >
> > Lately we have been doing a set of benchmarks for comparing the
> > efficiency of the (finite element) assembly routine in Dolfin with
> > respect to its different backends, MTL, PETSC, uBlas, and Epetra.
> > The benchmark results can be found here:
> > http://www.fenics.org/wiki/Benchmark
> >
> > It turns out that for equations using vector elements, Epetra is
> > somewhat slower than the other backends,  and most of this is
> > due to the
> > creation of the Epetra_FECrsGraph.
> >
> > We have implemented it roughly as follows:
> >
> > Epetra_Map row_map(dims[0], 0, Comm);
> > Epetra_FECrsGraph* epetra_graph;
> > epetra_graph = new Epetra_FECrsGraph(Copy, row_map, 0);
> > int numRows, numCols;
> > int *rows;
> > int *cols;
> > for each element:
> >     epetra_graph->InsertGlobalIndices(numRows, rows, numCols, cols);
> >
> > I assume that the problem here is that NumIndicesPerRow is not set
> > appropriately.
> > However, this code is 'library code' and we do not know (with the
> > current design in Dolfin) the number of indices per row when creating
> > the Epetra_FECrsGraph.
> >
> > Do you know any ways to speed up the code (without knowing
> > NumIndicesPerRow) ?
> >
> > Kent-Andre Mardal
> > Simula Research Laboratory
> >
> >
> >
> > _______________________________________________
> > Trilinos-Users mailing list
> > Trilinos-Users at software.sandia.gov
> > http://software.sandia.gov/mailman/listinfo/trilinos-users
> >
> 




More information about the Trilinos-Users mailing list