Just two more comments then, to clarify some things:

What we really want to do from the form compiler side of things,
is to reuse the local element basis function tabulation.
That doesn't necessarily require using dolfin::Function if you
have a better idea from the dolfin side.

In the case where x is in a finite element space with the dofs
sent by dolfin to tabulate_tensor as w[0][...], computing the
value of a component of x is then
  sum_k w[0][k] FE0[ip][k]
and the Jacobian can then be computed as
  sum_k w[0][k] FE0_D010[ip][k] # derivative table values
because J = dx/dX.


To support user-defined maps, we can introduce a more generic
interface for user-defined Coefficents (i.e. the UFL equivalent of a
dolfin::GenericFunction). Both user-defined x, Jacobian and Coefficients
depend on defining a ufc interface for callbacks and generating calls
through
that interface so the work involved in UFL/FFC would be largely overlapping.
This can be implemented as a fake element space like the "Quadrature"
element, possibly by introducing better notation than "quadrature" element.
  x = Coefficient(CallbackSpace())
possibly with derivatives added to the callback interface:
  x = Coefficient(CallbackSpace(degree=1))
so we can use J=dx/dX and direct computation of J to
x.eval_reference_derivatives() or something.
The ufc callback interface can be designed to fit with GenericFunction.

So from the ufl/ffc side of things I consider the user-defined callbacks to
be
an orthogonal issue. From the dolfin side you're right these issues are
connected
through the choice of dolfin type for the Coordinate.

Have a good weekend :)

Martin


On 14 March 2014 16:40, Garth N. Wells <[email protected]> wrote:

>
>
> On Fri, 14 Mar, 2014 at 3:22 PM, Anders Logg <[email protected]> wrote:
>
>> On Fri, Mar 14, 2014 at 01:32:00PM +0000, Garth N. Wells wrote:
>>
>>>
>>>
>>>  On Fri, 14 Mar, 2014 at 12:21 PM, Martin Sandve Alnæs
>>>  <[email protected]> wrote:
>>>  >I'm already in the process of implementing from the
>>>  >UFL/UFLACS(FFC) side, and getting fundamentals fully in place
>>>  >there will be the leading priority. By using a regular Coefficient
>>>  >for x everything falls into place quite nicely, and I'm doing
>>>  >other symbolic geometry as well that simplifies a lot of things.
>>>
>>>  I don't think we've discussed how this should be done on the FFC
>>>  side. One extreme is that the map is provided as a Function, and FFC
>>>  handles the Jacobian computation, etc using the Function. Other
>>>  cases are that UFC provides an interface function that given a
>>>  parent coordinate returns the map and Jacobian at that point, or
>>>  ufc::cell stores 'higher-order node' coordinates that can be used to
>>>  construct a polynomial map in the generated code.
>>>
>>>  I'm not challenging (yet :)) what's being done - I'd just like to
>>>  know if there is a plan.
>>>
>>>
>>>  >Answering Garths comments:
>>>  >Doing boundary cells only will require generating an additional
>>>  >tabulate-tensor implementation as well as other details which will
>>>  >lead this discussion into information overload. We'll have to push
>>>  >that to a later stage, maybe it will be easier after we have
>>>  >better designed support for functions on
>>>  >submeshes/restrictions/whatever we call them.
>>>
>>>  I think it should be considered now since it will affect the choice
>>>  of data structure/implementation, e.g. it immediately discounts
>>>  using dolfin::Function without major changes to Function.
>>>
>>>
>>>  >Doing non-polynomial mappings would basically involve introducing
>>>  >callable functions to UFL, which can be useful in other contexts
>>>  >as well. Getting the framework in place for defining x as a
>>>  >Coefficient in UFL and a Function in DOLFIN will be the first step
>>>  >(which is already well underway). By later introducing the "space
>>>  >of functions that can be evaluated in quadrature points"
>>>  >(including derivatives) should then sort out the arbitrary
>>>  >mappings you want as a side effect.
>>>
>>>  >So I say getting in place coordinates represented by any vector
>>>  >valued finite element function is a pretty good first step.
>>>
>>>
>>>  In terms of implementation I agree, in terms of deciding how we
>>>  start to go about now is the time to think forward. Deliberating
>>>  implicit in the two cases that I presented is that they both
>>>  discount using dolfin::Function.
>>>
>>>  >
>>>  >Answering Anders comments:
>>>  >I think the right place to get the coordinate function from is the
>>>  >MeshGeometry.
>>>  >
>>>  >It should be accessible from the Mesh, because the UFL abstraction
>>>  >is   SpatialCoordinate(domain)
>>>  >  Jacobian(domain)
>>>  >  etc.
>>>  >and ufl.Domain carries the dolfin.Mesh as domain.data().
>>>  >
>>>  >So if we can add a member
>>>  >
>>>  >  shared_ptr<Function> MeshGeometry::coordinate_function()
>>>  >
>>>  >that returns nullptr for meshes without a coordinate function set,
>>>  >I think that's enough for now and I can try to get the basics
>>>  >working.
>>>
>>>
>>>  I don't get why we would want a dolfin::Function with all the extras
>>>  that come with it. For the isoparametric case, we just need need an
>>>  interface that supplies the coordinates of 'higher-order' points in
>>>  addition to the vertex coordinates. We could also associate a
>>>  dolfin::FiniteElement with mesh geometry to provide and complete
>>>  definition of the coordinate field.
>>>
>>
>> What functionality does Function have that we don't want for a
>> CoordinateFunction?
>>
>> - Perhaps the possibility to extract subfunctions?
>>
>> - What about dofmap? Would it be enough with a fixed scheme (for
>>   ordering the coordinate values) or will we at some point want to
>>   reorder the list of coordinates?
>>
>
>
> We don't want dof map re-ordering, map from re-ordered dofmap to ufc
> dofmap, detection of dofmap blocks, dofmap parallel re-numbering,
> complicated assignment operators, extrapolation, GenericVector-based axpy,
> GenericVector storage, rank, value_dim, ghost updates, . . . .
>
> There are important specialisations that can be performed for a coordinate
> map. We may want to order the coordinates (for data locality), but we don't
> need a dofmap for this.
>
>
>
>> - When I think of it now, the biggest difference is that we don't need
>>   to store the values in a GenericVector since I don't see a need for
>>   sending that data to a linear solver (ever).
>>
>
>
> The dofmap is a big difference too - we don't need something so
> heavy-weight for coordinates.
> Note that for some problems getting data out of a Function for assembly
> takes a noticeable part of the assembly time. The performance cost will be
> too great.
>
>
>
>> So in conclusion, perhaps we will need to implement a new class
>> CoordinateFunction
>>
>
>
> Maybe.
>
>
>  which implements the GenericFunction interface but
>> is otherwise different from Function.
>>
>
> Not necessarily. We could attach the data to ufc::cell/cell_geometry,
> which would avoid virtual function calls.
>
> In any case, can we focus a discussion CoordinateFunction versus other
> possibilities, like data attached to ufc::cell_geometry? We should discount
> dolfin::Function. It's too heavy-weight.
>
> In favour of CoordinateFunction, it could have a function to indicate the
> map type for a given cell, and via CoordinateFunction::eval functions could
> support user-defined maps in the future. For now it could just support the
> ufc::function interface.
>
> Garth
>
>
>> --
>> Anders
>>
>
>
_______________________________________________
fenics mailing list
[email protected]
http://fenicsproject.org/mailman/listinfo/fenics

Reply via email to