dolfin team mailing list archive
-
dolfin team
-
Mailing list archive
-
Message #13143
Re: [FFC-dev] lagrange multiplier on boundary
On ma., 2009-04-20 at 17:17 -0400, Shawn Walker wrote:
>
> On Mon, 20 Apr 2009, kent-and@xxxxxxxxx wrote:
>
> >>
> >> On Mon, 20 Apr 2009, Anders Logg wrote:
> >>
> >>> On Thu, Apr 16, 2009 at 09:19:51PM -0400, Shawn Walker wrote:
> >>>>
> >>>> On Fri, 17 Apr 2009, Garth N. Wells wrote:
> >>>>
> >>>>> kent-and@xxxxxxxxx wrote:
> >>>>>> I would also like this capability! It is something that often shows
> >>>>>> up
> >>>>>> in inverse/optimal control problems.
> >>>>>>
> >>>>>> Written in FFC/UFL your first equation reads:
> >>>>>>
> >>>>>> dot(u,v)*dx - p*div(v)*dx + lmbda*dot(v,n)*ds
> >>>>>>
> >>>>>> where u, p, lmbda are trial functions.
> >>>>>>
> >>>>>> You could form one system or create a block matrix. Anyhow
> >>>>>> the term
> >>>>>> lmbda*dot(v,n)*ds
> >>>>>> would lead to a matrix with a very big kernel since you are not able
> >>>>>> to
> >>>>>> restrict the dofs of lmbda only to the boundary.
> >>>>>>
> >>>>>> What you can currently do is to restrict the functionspace for lmbda
> >>>>>> to
> >>>>>> all the cells
> >>>>>> associated with the boundary.
> >>>>>>
> >>>>>> Using restricted functionspaces (in a simpler fashion) can be found
> >>>>>> in
> >>>>>> demo/function/restriction.
> >>>>>>
> >>>>>> The restriction does only work on cells for now.
> >>>>>>
> >>>>>> We could discuss Uzawa and/or block matrices for this problem but I
> >>>>>> think
> >>>>>> the simplest start is to create one system to begin with.
> >>>>>>
> >>>>>> Whether it makes sense that lmbda lives on the whole cell associated
> >>>>>> with
> >>>>>> the boundary, I don't know.
> >>>>>>
> >>>>>
> >>>>> It should live only on the boundary. In practice this only becomes an
> >>>>> issue for higher-order elements with internal dofs.
> >>>>>
> >>>>> Garth
> >>>>
> >>>> Yes, I agree.
> >>>>
> >>>> So how ridiculous is it to enable FFC/DOLFIN to have finite element
> >>>> functions that are only defined on the boundary of the domain? I'm
> >>>> guessing there would be some special DoFmappings to go from the global
> >>>> domain numbering to a boundary numbering only. This would be really
> >>>> nice
> >>>> to have. There are lots of cases in practice that have these kinds of
> >>>> boundary functions.
> >>>>
> >>>> - Shawn
> >>>
> >>> It's not impossible but it requires some thought. I think Garth has
> >>> asked about this for a long time as well (to have function spaces that
> >>> only live on facets). I don't really know how to best handle it.
> >>>
> >>> --
> >>> Anders
> >>
> >> ok. I just implemented what I needed in MATLAB and that formulation
> >> works. But it would certainly be great to have it in FENICS.
> >>
> >> - Shawn
> >
> > A possible way to do it with not to much work (?) in FEniCS would be to
> > to create the matrix on the space with to many degrees of freedom and
> > a large kernel. Then create a projection matrix based on the degrees
> > of freedom you want. You may then project the matrix onto the space you
> > want. This is similar to what is currently done in the function restriction
> > (allthought it is between the lines here). Both PETSc and Trilinos support
> > matrix matrix multiplication (I think).
> >
> > This is maybe not the most elegant solution, but it is not all bad.
> >
> > Kent
>
> Is it really that bad to put this directly into FENICS? It seems that all
> you need is a separate DoFmap for the variables that live on the boundary
> only. Is it not possible to setup a DoF numbering on the boundary only?
> Maybe, in addition, have a way of mapping the DoF's on the boundary to the
> corresponding DoF's in the global numbering?
>
> Wasn't something like this done for the parallel stuff? I seem to
> remember something about having mapping between DoF's on different
> sub-domains, but I wasn't really paying attention.
>
> If someone could summarize the difficulties here, that would be great.
> I am not that familiar with this aspect of FENICS.
>
> - Shawn
In Dolfin, the dof_map is very much tied to the mesh.
As I understand it, it is bascially two ways to create
dof_maps on parts of the mesh.
1) Create a Mesh on the subset and provide a mapping between the
two meshes. This might be the cleanest solution (?).
2) Having a dof_map and an array of indicators of which of the dofs
that should be used. In this case an auxillary dofmap must be made.
This dof_map is a collapsed version of the old dof_map. ie.
If we have a dofmap
(0,1,2,3) and and indicator (1,-1,1,-1) then the new
dof_map should be (0,1) and the mapping between them should
be { 0:0, 1:-1, 2:1, 3:-1 }
Hence, two auxillary arrays of ints are needed.
In addition to this, some editing of the assemble function might be
needed, but I don't think this should be hard.
Kent
Follow ups
References