Dear Praveen,

Regarding the performance of your code, there are several fundamental
problems:
- You run the constructor of Quadrature (allocates memory) and FEValues
(allocates lots of memory, evaluates a number of things you don't need)
that are both not made for use within the scenario of FEFieldFunction or
your manual loop.
- Since you're only using the shape values, you could skip all the
FEValues stuff and simply run the loop as
  solution += solution0(local_dof_indices[i]) * fe.shape_value(i, p_ref);
- At this point, the performance depends on the
transform_real_to_unit_cell (which might be slow due to Newton
iterations) and the polynomial evaluation of fe.shape_value and to some
extent the virtual function call of the latter. It depends on your
polynomial space whether this is expensive.
- You could speed up the polynomial evaluation significantly if the
polynomial basis had a way to evaluate all polynomials at once - one
could use the tensor product structure here. If you're interested, I
could give you some code to do that.

I am not aware of any generic performance optimization beyond that level
for generic mappings, but I'm also not an expert of semi-Lagrangian
methods. The works that I know of do usually work on Cartesian grids
where one can get tremendous speedups of the
transform_real_to_unit_cell. Then you spend most of the time in the
polynomial interpolation steps again - but probably a factor of 100
faster than with what you have now.

Best,
Martin


On 05.09.2017 15:18, Praveen C wrote:
> Hello Bruno
>
>> Are you sure that  FEFieldFunction.value() is the slowest part? In
>> 3D, point_inside could also be slow (you may need to do a few Newton
>> Iterations). It would be nice if you ran your code using a profiler
>> to see what is the bottleneck (building the quadrature in value(p1),
>> the evaluation itself, etc.).
>>
>
> If I comment out call to FEFieldFunction.value(), the code runs a LOT
> faster, so there is no doubt this is the slow part. I am in 2-d at
> present.
>
> I have set the correct cell by call to set_active_cell. Does the value
> function again check if the point is in this cell ?
>
> I implemented my own function to evaluate the solution as follows, but
> this is equally slow. I havent used a profiler but will timing parts
> of this below function be useful ?
>
> template<intdim>
> doubleProblem<dim>::evaluate_solution
> (typenameDoFHandler<dim>::active_cell_iterator cell,
>                                         constPoint<dim> &p) const
> {
>    // Transform p to unit cell
>    Point<dim> p_ref = mapping.transform_real_to_unit_cell(cell, p);
>
>    // Create quadrature rule
>    Quadrature<dim> quadrature(p_ref);
>
>    FEValues<dim> fe_values (mapping,
>                             fe,
>                             quadrature,
>                             update_values);
>    std::vector<types::global_dof_index> local_dof_indices
> (fe.dofs_per_cell);
>
>    fe_values.reinit(cell);
>    cell->get_dof_indices (local_dof_indices);
>    doublesolution = 0.0;
>    for(unsignedinti=0; i<fe.dofs_per_cell; ++i)
>       solution += solution0(local_dof_indices[i]) *
> fe_values.shape_value(i, 0);
>
>    returnsolution;
> }
>
> I am processing one point at a time which could be slow. I am now
> rewriting the code so that all points falling in a cell are evaluated
> at once. 
>
> Thanks
> praveen
> -- 
> The deal.II project is located at http://www.dealii.org/
> For mailing list/forum options, see
> https://groups.google.com/d/forum/dealii?hl=en
> ---
> You received this message because you are subscribed to the Google
> Groups "deal.II User Group" group.
> To unsubscribe from this group and stop receiving emails from it, send
> an email to dealii+unsubscr...@googlegroups.com
> <mailto:dealii+unsubscr...@googlegroups.com>.
> For more options, visit https://groups.google.com/d/optout.

-- 
The deal.II project is located at http://www.dealii.org/
For mailing list/forum options, see 
https://groups.google.com/d/forum/dealii?hl=en
--- 
You received this message because you are subscribed to the Google Groups 
"deal.II User Group" group.
To unsubscribe from this group and stop receiving emails from it, send an email 
to dealii+unsubscr...@googlegroups.com.
For more options, visit https://groups.google.com/d/optout.

Reply via email to