← Back to team overview

dolfin team mailing list archive

Re: Matrix concatenation

 

On Tue, Sep 30, 2008 at 3:22 PM, Matthew Knepley <knepley@xxxxxxxxx> wrote:

> On Tue, Sep 30, 2008 at 4:46 AM, Evan Lezar <evanlezar@xxxxxxxxx> wrote:
> > On Tue, Sep 30, 2008 at 10:53 AM, Jed Brown <jed@xxxxxxxx> wrote:
> >>
> >> On Tue 2008-09-30 10:48, Evan Lezar wrote:
> >> > I was wondering what the best way is to concatenate PETSc matrices?  I
> >> > am
> >> > currently using the set() method, but this takes quite a bit of time.
> >>
> >> What are you trying to do?  You should probably preallocate the full
> >> matrix from the beginning.
> >
> > I am assembling a number of submatrices.  And then trying to concatenate
> > them to form one large eigensystem.  I do pre-allocate the matrix that I
> > want to place the entries in.
>
> I agree with Jed. The best way to do this is to assemble the full system
> matrix,
> and then if submatrices are required, you pull them out using
> MatGetSubmatrix().
>
> The only situation I can see that is different is in some algorithms
> like FETI in
> which submatrices are unassembled. Then you would have to use a MATIS
> structure, or something else depending on your algorithm.
>
>  Matt
>
> > It may be that I would be able to use mixed elements to do the same
> thing.
> > I will see if I can get a code snippet together.
> >
> > Evan
> >
>

I am trying to model the transverse and axial components of the electric
field in a waveguide - using Nededec vector elements and Lagrange elements
respectively.

Consider the code snippet below that describes the matrices:

====================================

    vector_element = FiniteElement("Nedelec", "triangle", 1)
    nodal_element = FiniteElement("Lagrange", "triangle", 1)

    N_v = TestFunction(vector_element);
    N_u = TrialFunction(vector_element);

    L_v = TestFunction(nodal_element);
    L_u = TrialFunction(nodal_element);

    """ define the forms """
    a_tt = dot(curl_t(N_v), curl_t(N_u))*dx - dot(N_v, N_u)*dx

    b_tt = dot(N_v, N_u)*dx
    b_tz = dot(N_v, grad(L_u))*dx
    b_zt = dot(grad(L_v), N_u)*dx

    b_zz = dot(grad(L_v), grad(L_u))*dx - L_v*L_u*dx

    """ define the mesh """
    width = 22.86e-3
    height = 10.16e-3
    mesh = Rectangle(0, width, 0, height, 1, 2, 0)
    mesh.refine()
    mesh.refine()
    c = 2.99792458e8

    """ Assemble the sub-matrices """

    B_tt = PETScMatrix()
    assemble(b_tt, mesh, tensor=B_tt)
    B_tz = assemble(b_tz, mesh)
    B_zt = assemble(b_zt, mesh)
    B_zz = assemble(b_zz, mesh)

    A_tt = PETScMatrix()
    assemble((a_tt, mesh, tensor=A_tt)

=============================
I then need to construct the following matrices to solve an eigensystem:

L = [ B_tt,  B_tz ]
       [ B_zt,  B_zz ]

and

R = [ B_tt + k*A_tt, B_tz  ]
       [ B_zt,                B_zz ]



So, is someone knows how I can define a mixed element that achieves the same
it would be much appreciated, or a means to concatenate the matrices
effectively.

Evan



> >>
> >> Jed
> >>
> >> _______________________________________________
> >> 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
> >
> >
>
>
>
> --
> What most experimenters take for granted before they begin their
> experiments is infinitely more interesting than any results to which
> their experiments lead.
> -- Norbert Wiener
>

Follow ups

References