On 3 September 2015 at 16:50, Garth N. Wells <gn...@cam.ac.uk> wrote:

>
> On 3 September 2015 at 14:42, Chris Richardson <ch...@bpi.cam.ac.uk>
> wrote:
>
>> On 03/09/2015 14:10, Martin Sandve Alnæs wrote:
>>
>>> Part of the complexity in this chain steps from the numerous dolfin
>>> signatures. As discussed before, these also need cleanup to clarify
>>> collective/non-collective operations. While at it, we could also
>>> vectorize on x by allowing a number of coordinates wherever there's an
>>> x argument in the code below.
>>>
>>> For a moment I'd like to ignore the ufc/ffc parts and figure out how
>>> many eval signatures we really need in dolfin? In particular to
>>> resolve the collective/non-collective part. Then we can design the
>>> callback mechanism to match the ideal dolfin design.
>>>
>>> Here's one take (dodging the type of 'cell'):
>>>
>>
>> Is it as simple as just having a "collective" and "non_collective" eval?
>>
>> Typically, processes will iterate over some set of "local" points, asking
>> for evaluations ( - OK, they may be "local" on one mesh, but not another).
>> When the point is not "on process", the process needs to get it, somehow.
>> But all the other processes are unaware - they are thinking about
>> different points.
>>
>> It is quite unusual for all processes to want to evaluate the same point
>> at the same time, for a simple collective operation.
>>
>> My thought was that we would have to store up a list of points, and then
>> do a collective operation to resolve all the cross-process evaluations.
>>
>> Or maybe I have missed something...
>>
>>
> My impression is that Martin is just suggesting the we split evals into
> non-collective and collective 'groups', and have as few versions of eval in
> each group. With fewer evals, we have few code paths to think about.
>
> The interface for collective evals is orthogonal to this, and something we
> need to look closely at.
>
> Garth
>


Exactly. The collective variants of eval will always end up calling
non-collective variants of eval but never the other way around.
Also the collective eval variants won't be virtual functions.

The evaluation of a GenericFunction on cells that are known to be
on the local process differs between Function and Expression, but
can be handled with a single signature "virtual void eval(v,x,cell) const =
0;"
in GenericFunction, implemented by Function and by the users subclass
of Expression. The only thing we need to clean up there in the ufc::function
is the ufc::cell vs dolfin::cell.

For a collective eval over multiple points, I currently see only two
variants:
a) each process pass a different set of point(s) and wants the
corresponding results
b) each process actually pass the same point(s) and wants all results
In both cases there is first communication of 'who owns which points',
then each process makes calls to the non-collective evals for its own
points,
then the results are communicated back to the right place.
I believe a) is what Chris was talking about as the most common operation.

Martin



>
>
>
>> Chris
>>
>>
>>> class GenericFunction
>>> {
>>> public:
>>>   // Known cell on process, non-collective call:
>>>
>>>   // Note: If we allow multiple x within cell,
>>>   //   the implementation can be vectorized
>>>
>>>   virtual void noncollective_eval(v, x, cell) const = 0;
>>>
>>>   // Unknown cell, collective call:
>>>
>>>   // Note: there could be a non-collective version with
>>>
>>>   //   unknown cell, which would fail if the point(s) are not local
>>>
>>>   void collective_eval(v, x) const // not virtual!
>>>   {
>>>
>>>       // Note: If we allow multiple x, we can reduce
>>>       //   the number of mpi calls significantly in here
>>>
>>>       cell = find cell x is on
>>>
>>>       if on this process:
>>>           noncollective_eval(v, x, cell)
>>>
>>>           send v
>>>
>>>       else:
>>>
>>>           receive v
>>>
>>>       // vectorized would be something like his:
>>>
>>>       find all (x,cell) on this process
>>>
>>>       non_collective_eval(v,x,cell) for all local points
>>>
>>>       send and receive all v
>>>   }
>>> };
>>>
>>> class Function: public GenericFunction
>>> {
>>> public:
>>>     virtual void noncollective_eval(v, x, cell) const
>>>     {
>>>         dofs = restrict to cell
>>>         fe = evaluate basis in x on cell
>>>         v = sum dofs*fe
>>>     }
>>> };
>>>
>>> class Expression: public GenericFunction
>>> {
>>> public:
>>>     virtual void noncollective_eval(v, x, cell) const = 0;
>>> };
>>>
>>> class MyExpression: public Expression
>>> {
>>> public:
>>>     virtual void noncollective_eval(v, x, cell)
>>>     {
>>>         v = user defined code (may involve e.g. cell number lookup)
>>>     }
>>> }
>>>
>>> I don't really see why we need more than one collective and one
>>> non-collective eval.
>>>
>>> Martin
>>>
>>> On 3 September 2015 at 11:54, Chris Richardson <ch...@bpi.cam.ac.uk>
>>> wrote:
>>>
>>> This is a mailing list continuation of a bitbucket discussion:
>>>>
>>>>
>>> https://bitbucket.org/fenics-project/dolfin/pull-requests/240/interpolation-improvements
>>>
>>>> [1]
>>>>
>>>>
>>>> I have been trying to clean up Function::interpolate() a bit.
>>>> The following issues have come up:
>>>>
>>>> 1. The chain of function calls is long and complex:
>>>>
>>>> e.g. to interpolate in the general case from another Function v:
>>>>
>>>> - this is just an illustrative sketch, with details removed:
>>>>
>>>> Function::interpolate(GenericFunction& v)
>>>> FunctionSpace::interpolate(Vector& coeffs, const GenericFunction&
>>>> v)
>>>> Iterate over cells {
>>>> v.restrict(cell_coeffs, element, cell, coords, ufc_cell)
>>>> Function::restrict(w, element, cell, coords, ufc_cell)
>>>> ...
>>>> restrict_as_ufc_function(w, element, dolfin_cell, coords,
>>>> ufc_cell)
>>>> element.evaluate_dofs()
>>>> ufc_element->evaluate_dofs()
>>>> ufc_function::evaluate(vals, coords, ufc_cell)
>>>> GenericFunction::evaluate(vals, coords, ufc_cell)
>>>> GenericFunction::eval(vals, x, ufc_cell)
>>>> GenericFunction::eval(vals, x)
>>>> Function::eval(vals, x)
>>>> {
>>>> id =
>>>> mesh.bb_tree.compute_first_entity_collision(point)
>>>> Cell (mesh, id)
>>>> Function::eval(val, x, cell, ufc_cell)
>>>> restrict(coeffs, element, cell, coords, ufc_cell)
>>>> element.evaluate_basis()
>>>> values += coeffs * basis
>>>> }
>>>>
>>>> 2. Many calls use ufc_cell and dolfin_cell with the same data. This
>>>> looks bad, and is bad. Can we find a way to transition away from
>>>> needing both?
>>>>
>>>> 3. In parallel, this model totally breaks down, as it is only at
>>>> the very bottom that we realise the point is 'off process' and we
>>>> need a collective operation to get it.
>>>>
>>>> _______________________________________________
>>>> fenics mailing list
>>>> fenics@fenicsproject.org
>>>> http://fenicsproject.org/mailman/listinfo/fenics [2]
>>>>
>>>
>>>
>>>
>>> Links:
>>> ------
>>> [1]
>>>
>>> https://bitbucket.org/fenics-project/dolfin/pull-requests/240/interpolation-improvements
>>> [2] http://fenicsproject.org/mailman/listinfo/fenics
>>>
>>
>> --
>> Chris Richardson
>> BP Institute
>> Madingley Road
>> Cambridge CB3 0EZ
>>
>> _______________________________________________
>> fenics mailing list
>> fenics@fenicsproject.org
>> http://fenicsproject.org/mailman/listinfo/fenics
>>
>
>
_______________________________________________
fenics mailing list
fenics@fenicsproject.org
http://fenicsproject.org/mailman/listinfo/fenics

Reply via email to