dolfin team mailing list archive
-
dolfin team
-
Mailing list archive
-
Message #22999
Re: evaluate_basis() vs evaluate_basis_all() in python dolfin
On Thursday May 5 2011 10:16:11 Neilen Marais wrote:
> Hi,
>
> On Thu, Apr 21, 2011 at 6:18 PM, Johan Hake <johan.hake@xxxxxxxxx> wrote:
> >> > We should probably ignore all the ufc::cell versions, and add unit
> >> > test for the one we do expose...
> >>
> >> Please excuse my naivete, but what is the essential difficulty in
> >> exposing the functions to the python side?
> >
> > We have not paid too much attention making methods taking generated code,
> > ufc objects, available from Python. This is too low level, and some of
> > the ufc objects gets reordered when the dolfin objects gets
> > instantitated. This is mostly the case when run in parrallel.
> >
> >> I ask since from a perspective of both testability of my own code and
> >> of learning what FEniCS does it is quite usefull to be able to call
> >> the lower level functions from python. For instance, I also want to
> >> use finite_element.map_to_reference_cell(), but it seems to give the
> >> same error that evaluate_basis_all() did before you fixed it :)
> >
> > Yes, this would be the same case here. We also need to consider the
> > interface. It would clutter the interface to duplicate too many
> > functions. But if Garth or Anders is not opposed to it we could expose
> > that one too.
>
> For me it would make sense to allow these things from python, since
> that is usually where I'd want to experiment with a new technique
> before trying to implement it fast in C++. I'm a little out of my
> depth here, but to me seems that "dolfin" is an easier to use wrapper
> around the lower level UFL stuff.
>
> I would potentially be happy to code to a slightly lower level, but
> let me sketch out what I'm looking for in more detail. I don't want to
> have to implement mesh interconnectivity handling and DOF mapping from
> scratch. At the same time I want the ability to formulate things at
> the level of basis function evaluations as needed. In part this is
> simply because the "standard" dolfin doesn't cater quite that well for
> stuff I need for electromagnetic simulation (yet?), and in part since
> I may want to try new and/or hare-brained things :)
>
> I'm quite happy to work with FEniCS developers to get the needed
> functionality into dolfin, and that is really first prize for me. In
> the meantime, however, my slave-driving overlords want a working EM
> FEM solver with academic outputs, so I can't neccesarily wait for that
> to happen first :)
>
> >> If it is not to difficult you could teach me how to fix such problems.
> >> Then I could go fishing rather than asking for fish periodically :)
> >
> > We need more fishermen! Just look at the last commit, where i added the
> > dolfin version. It should be pretty straight forward to implement a
> > similar method for the other methods which takes a ufc::cell.
>
> OK, as far as I can tell the trick is just to cast the dolfin cell to
> a UFC cell. I tried the patch below to get mapt_to_reference_cell()
> working, but receive:
>
> RuntimeError: map_to_reference_cell not yet implemented (introduced in UFC
> 2.0).
This you need to address to the FFC list, as FFC has not generated code for
that ufc method.
> Am I right in assuming that this is a limitation in the current dolfin
> C++ codebase rather than an issue with the code I added?
>
> On another note, if it simply down to casting a dolfin cell to a UFC
> cell, could one not define a python casting function that returns a
> "UFC compatible" python object, and then tell the user that if they
> want to access the lower level UFC functionality they need to do that
> first?
I think this should be OK for UFCCell. But we should typically not expose
ufc::dof_map, as this is might be reordered.
One way to handle this is to add a typemap at the Python side, so you can send
a dolfin::Cell to any method accepting a ufc::cell. This will fix the Python
side of this request.
But I wonder why we need a dolfin version of any of these methods as one in
C++ can just create a UFCCell and pass that one instead.
So a question to all C++ users. Would you mind if I removed all dolfin::Cell
versions of the dolfin::FiniteElement::foo and force you to call these methods
using a create UFCCell instead of dolfin::Cell. If not I think we need to add
dolfin::Cell version of all ufc::cell methods in dolfin::FiniteElement,
otherwise we will get new requests in the future.
Johan
>
> Cheers
> Neilen
>
>
> --- dolfin/fem/FiniteElement.h 2011-04-20 23:21:30 +0000
> +++ dolfin/fem/FiniteElement.h 2011-05-05 16:49:15 +0000
> @@ -192,6 +192,16 @@
> _ufc_element->evaluate_basis_all(values, coordinates, ufc_cell);
> }
>
> + /// Map from coordinate x in cell to coordinate xhat in reference cell
> + void map_to_reference_cell(double* xhat,
> + const double* x,
> + const Cell& cell) const
> + {
> + assert(_ufc_element);
> + UFCCell ufc_cell(cell);
> + _ufc_element->map_to_reference_cell(xhat, x, ufc_cell);
> + }
> +
>
> > Johan
> >
> >> Cheers
> >> Neilen
> >>
> >> _______________________________________________
> >> Mailing list: https://launchpad.net/~dolfin
> >> Post to : dolfin@xxxxxxxxxxxxxxxxxxx
> >> Unsubscribe : https://launchpad.net/~dolfin
> >> More help : https://help.launchpad.net/ListHelp
>
> _______________________________________________
> Mailing list: https://launchpad.net/~dolfin
> Post to : dolfin@xxxxxxxxxxxxxxxxxxx
> Unsubscribe : https://launchpad.net/~dolfin
> More help : https://help.launchpad.net/ListHelp
Follow ups
References