← Back to team overview

dolfin team mailing list archive

Re: tensor function space in Dolfin

 

On Monday 20 April 2009 22:59:19 Martin Sandve Alnæs wrote:
> See TensorElement in UFL, there should be a natural correspondence
>   FiniteElement <-> FunctionSpace
>   MixedElement <-> MixedFunctionSpace
>   VectorElement <-> VectorFunctionSpace
>   TensorElement <-> TensorFunctionSpace

I think we have a problem here.

The value_shape of the mixed element in your TensorFunctionSpace is probably 
not computed correct here Kent. This is done in the TensorElement construtor, 
which is not called in your case.

We could let TensorFunctionSpace inherit FunctionSpaceBase instead of 
MixedFunctionSpace. Then in the constructor of TensorFunctionSpace we can:

  element = TensorElement(...)

and pass the element to FunctionSpaceBase. Unfortunaltely we break the 
inheritance correlation with ufl. And some of the implemented PyDOLFIN 
functions such as split, check the provided functionspace for beeing a 
MixedFunctionSpace. this could be changed to check for the element of the 
provided FunctionSpace to be a MixedElement.

We still needs to construct the SubSpaces though, which could be done in the 
TensorFunctionSpace constructor similare to what is done in the 
MixedFunctionSpace constructor?

Not sure if this approach covers all bases? 


Johan


> Martin
>
> On Mon, Apr 20, 2009 at 10:30 PM,  <kent-and@xxxxxxxxx> wrote:
> > Are there any thoughts on how to implement TensorFunctionSpace in Dolfin
> > ?
> >
> > Copying code for VectorFunctionSpace it could be something like:
> >
> >
> >  class TensorFunctionSpace(MixedFunctionSpace):
> >    "VectorFunctionSpace represents a vector-valued finite element
> > function space."
> >
> >    def __init__(self, mesh, family, degree, dim=None):
> >        """Create tensor-valued finite element function space. The
> > function space
> >        may be created by
> >
> >            V = TensorFunctionSpace(mesh, family, domain, degree,
> > dim=None)
> >
> >            mesh    : a Mesh
> >            family  : a string specifying the element family
> >            degree  : the degree of the element
> >            dim     : an optional argument specifying the number of
> > components
> >
> >        If the dim argument is not provided, the dimension will be deduced
> > from
> >        the dimension of the mesh.
> >        """
> >
> >        # Create subspaces
> >        dim = dim or mesh.geometry().dim()
> >        spaces = dim*[VectorFunctionSpace(mesh, family, degree, dim)]
> >
> >        # Initialize base class
> >        MixedFunctionSpace.__init__(self, spaces)
> >
> >
> >
> > But then I get the following when in use:
> >
> >    condition or error(*message)
> >  File "/usr/local/lib/python2.5/site-packages/ufl/log.py", line 97, in
> > error raise UFLException(self._format_raw(*message))
> > ufl.log.UFLException: Transposed is only defined for rank 2 tensors.
> >
> > Any thoughts ?
> >
> >
> > Kent
> >
> >
> >
> > _______________________________________________
> > 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