dolfin team mailing list archive
-
dolfin team
-
Mailing list archive
-
Message #04422
Re: [HG dolfin] Implement assembly over interior facets. The result seems to
On Fri, Mar 09, 2007 at 07:10:31PM +0100, Garth N. Wells wrote:
> Anders Logg wrote:
> > On Fri, Mar 09, 2007 at 03:07:43PM +0100, Garth N. Wells wrote:
> >> Why do we need a class AssemblyMatrix? GenericMatrix/GenericVector
> >> should provide the necessary virtual functions for assembly.
> >>
> >> Garth
> >
> > Yes, definitely. The assembler assumes it gets a GenericTensor that
> > has a function
> >
> > void add(real* block, uint* size, uint** entries)
> >
> > (entries is probably a bad choice of name here, let's change it)
> > This makes it possible to implement just one assembly function and
> > keep the Assembler class small.
> >
>
> Then GenericVector and GenericMatrix should then derive from GenericTensor?
Yes. So far it is just two functions that are needed, init() and add().
> > Anyway, the class AssemblyMatrix is just a temporary implementation of
> > something that implements the GenericTensor interface. My suggestion
> > is that we overload the add() function above in both GenericMatrix and
> > GenericVector. The overloaded functions will map to the correct add()
> > function in the respective classes. For GenericMatrix, we'll have
> >
> > void add(real* block, uint* size, uint** entries)
> > {
> > add(block, size[0], entries[0], size[1], entries[1]);
> > }
> >
> > So, even if AssemblyMatrix could be of some use (it seems to be a fast
> > option for assembly, the STL vector of maps beat PETSc in at least one
> > benchmark I made)
>
> Don't draw conclusions too quickly. STL vector of maps is fast to
> assemble, but is close to useless for linear algebra. The current
> implementation for the assembly uBLAS matrices use something like STL
> vector of maps (fast to assemble but slow for linear algebra) and
> converts it to compressed row storage (fast for linear algebra). The
> conversion can take some time. A STL vector of integer maps should be
> created and then used to initialise a particular sparse storage format
> (before first use only). Then we wouldn't have to convert matrix formats.
Do you mean this should be done in DofMap?
/Anders
> Garth
>
> it is just a temporary implementation until we decide
> > to implement the GenericTensor interface in GenericMatrix and
> > GenericVector.
> >
> > /Anders
> >
> >> DOLFIN wrote:
> >>> One or more new changesets pushed to the primary DOLFIN repository.
> >>> A short summary of the last three changesets is included below.
> >>>
> >>> changeset: 2698:eb995ba168f6249aa061bf16b6b094d7d697090b
> >>> tag: tip
> >>> user: "Anders Logg <logg@xxxxxxxxx>"
> >>> date: Fri Mar 09 13:40:42 2007 +0100
> >>> files: TODO src/kernel/fem/Assembler.cpp src/kernel/fem/UFC.cpp src/kernel/fem/dolfin/UFC.h src/sandbox/assembly/Poisson.form src/sandbox/assembly/Poisson.h src/sandbox/assembly/PoissonOld.form src/sandbox/assembly/PoissonOld.h src/sandbox/assembly/main.cpp
> >>> description:
> >>> Implement assembly over interior facets. The result seems to
> >>> be wrong (compared to the old assembly) but the new assembly
> >>> is considerably faster... :-)
> >>>
> >>>
> >>> changeset: 2697:d9af2a3e540db91dc68936ff188c588b0cc259e9
> >>> parent: 2696:e7e54e87c948b580979f1041facb491e5dacd102
> >>> parent: 2695:1099832e1b156f53450783a81243609beab71fff
> >>> user: "Anders Logg <logg@xxxxxxxxx>"
> >>> date: Fri Mar 09 09:17:50 2007 +0100
> >>> files: src/kernel/la/PETScEigenvalueSolver.cpp src/kernel/la/dolfin/PETScEigenvalueSolver.h
> >>> description:
> >>> Merge assembly over exterior facets
> >>>
> >>>
> >>> changeset: 2696:e7e54e87c948b580979f1041facb491e5dacd102
> >>> parent: 2690:24a6d6a566cdbdccb2c253cc62522fd265dad1af
> >>> user: "Anders Logg <logg@xxxxxxxxx>"
> >>> date: Thu Mar 08 23:30:09 2007 +0100
> >>> files: src/kernel/fem/Assembler.cpp src/sandbox/assembly/Poisson.form src/sandbox/assembly/Poisson.h src/sandbox/assembly/main.cpp
> >>> description:
> >>> Implement assembly over exterior facets (interior facets remain).
> >>>
> >>>
> >>> ----------------------------------------------------------------------
> >>> For more details, visit http://www.fenics.org/hg/dolfin
> >>> _______________________________________________
> >>> DOLFIN-dev mailing list
> >>> DOLFIN-dev@xxxxxxxxxx
> >>> http://www.fenics.org/mailman/listinfo/dolfin-dev
> >>>
> >> _______________________________________________
> >> DOLFIN-dev mailing list
> >> DOLFIN-dev@xxxxxxxxxx
> >> http://www.fenics.org/mailman/listinfo/dolfin-dev
> > _______________________________________________
> > DOLFIN-dev mailing list
> > DOLFIN-dev@xxxxxxxxxx
> > http://www.fenics.org/mailman/listinfo/dolfin-dev
> >
>
>
> _______________________________________________
> DOLFIN-dev mailing list
> DOLFIN-dev@xxxxxxxxxx
> http://www.fenics.org/mailman/listinfo/dolfin-dev
Follow ups
References