" What type of Mapping are you using? If you take a look at
https://github.com/dealii/dealii/blob/ad13824e599601ee170cb2fd1c7c3099d3d5b0f7/source/matrix_free/fe_point_evaluation.cc#L40-L95
you can see when the fast path of FEPointEvaluation is taken. Indeed, the
slow path is (FEValues). One question: are you running in release or debug
mode? "

I use FE_Q<1>(1) with a MappingQ<1>(1) and
FE_Q<2>(1) with a MappingQ<2>(1).

I am running in release mode.

Best,
Simon

Am Do., 20. Okt. 2022 um 16:53 Uhr schrieb Peter Munch <
peterrmue...@gmail.com>:

> > FEPointEvaluation creates an FEValues object along with a quadrature
> object under the hood.
> Closer inspection revealed that all constructors, destructors,...
> associated with FEPointEvaluation
> need roughly 5000 instructions more (per call!).
> That said, FEValues is indeed the faster approach, at least for FE_Q
> elements.
>
> What type of Mapping are you using? If you take a look at
> https://github.com/dealii/dealii/blob/ad13824e599601ee170cb2fd1c7c3099d3d5b0f7/source/matrix_free/fe_point_evaluation.cc#L40-L95
> you can see when the fast path of FEPointEvaluation is taken. Indeed, the
> slow path is (FEValues). One question: are you running in release or debug
> mode?
>
> Hope this brings us closer to the issue,
> Peter
>
> On Thursday, 20 October 2022 at 16:47:17 UTC+2 Simon wrote:
>
>> Update:
>>
>> I profiled my program with valgrind --tool=callgrind and could figure out
>> that
>> FEPointEvaluation creates an FEValues object along with a quadrature
>> object under the hood.
>> Closer inspection revealed that all constructors, destructors,...
>> associated with FEPointEvaluation
>> need roughly 5000 instructions more (per call!).
>> That said, FEValues is indeed the faster approach, at least for FE_Q
>> elements.
>>
>> export DEAL_II_NUM_THREADS=1
>> eliminated the gap between cpu and wall time.
>> Using FEValues directly, I get cpu time 19.8 seconds
>> and in the case of FEPointEvaluation cpu time = 21.9 seconds;
>> Wall times are in the same ballpark.
>> Out of curiosity, why produces multi-threading such high wall times (200
>> seconds) in my case?.
>>
>> These times are far too big given that the solution of the linear system
>> takes only about 13 seconds.
>> But based on what all of you have said, there is probably no other to way
>> to implement my problem.
>>
>> Best
>> Simon
>>
>> Am Do., 20. Okt. 2022 um 11:55 Uhr schrieb Simon Wiesheier <
>> simon.w...@gmail.com>:
>>
>>> Dear Martin and Wolfgang,
>>>
>>> " You seem to be looking for FEPointEvaluation. That class is shown in
>>> step-19 and provides, for simple FiniteElement types, a much faster way to
>>> evaluate solutions at arbitrary points within a cell. Do you want to give
>>> it a try? "
>>>
>>> I implemented the FEPointEvaluation approach like this:
>>>
>>> FEPointEvaluation<1,1> fe_eval(mapping,
>>>                                         FE_Q<1>(1),
>>>                                         update_gradients |
>>> update_values);
>>> fe_eval.reinit(cell,
>>> make_array_view(std::vector<Point<1>>{ref_point_energy_vol}));
>>> Vector<double> p_dofs(2);
>>> cell->get_dof_values(solution_global, p_dofs);
>>> fe_eval.evaluate(make_array_view(p_dofs),
>>>                                     EvaluationFlags::values |
>>> EvaluationFlags::gradients);
>>> double val = fe_eval.get_value(0);
>>> Tensor<1,1> grad = fe_eval.get_gradient(0);
>>>
>>> I am using FE_Q elements of degree one and a MappingQ object also of
>>> degree one.
>>>
>>> Frankly, I do not really understand the measured computation times.
>>> My program has several loadsteps with nested Newton iterations:
>>> Loadstep 1:
>>> Assembly 1: cpu time 12.8 sec  wall time 268.7 sec
>>> Assembly 2: cpu time 17.7 sec  wall time 275.2 sec
>>> Assembly 3: cpu time 22.3 sec  wall time 272.6 sec
>>> Assembly 4: cpu time 23.8 sec  wall time 271.3sec
>>> Loadstep 2:
>>> Assembly 1: cpu time 14.3 sec  wall time 260.0 sec
>>> Assembly 2: cpu time 16.9 sec  wall time 262.1 sec
>>> Assembly 3: cpu time 18.5 sec  wall time 270.6 sec
>>> Assembly 4: cpu time 17.1 sec  wall time 262.2 sec
>>> ...
>>>
>>> Using FEValues instead of FEPointEvaluation, the results are:
>>> Loadstep 1:
>>> Assembly 1: cpu time 23.9 sec  wall time 171.0 sec
>>> Assembly 2: cpu time 32.5 sec  wall time 168.9 sec
>>> Assembly 3: cpu time 33.2 sec  wall time 168.0 sec
>>> Assembly 4: cpu time 32.7 sec  wall time 166.9 sec
>>> Loadstep 2:
>>> Assembly 1: cpu time 24.9 sec  wall time 168.0 sec
>>> Assembly 2: cpu time 34.7 sec  wall time 167.3 sec
>>> Assembly 3: cpu time 33.9 sec  wall time 167.8 sec
>>> Assembly 4: cpu time 34.3 sec  wall time 167.7 sec
>>> ...
>>>
>>> Clearly, the fluctuations using FEValues are smaller than in case of
>>> FEPointEvaluation.
>>> Anyway, using FEPointEvaluation the cpu time is smaller but the wall
>>> time substantially bigger.
>>> If I am not mistaken, the values cpu time 34.3 sec and wall time 167.7
>>> sec mean that
>>> the cpu needs 34.3 sec to execute my assembly routine and has to wait in
>>> the
>>> remaining 167.7-34.3 seconds.
>>> This huge gap between cpu and wall time has to be related to what I do
>>> with FEValues or FEPointEvaluation
>>> as cpu and wall time are nearly balanced if I use either neither of them.
>>> What might be the problem?
>>>
>>> Best
>>> Simon
>>>
>>>
>>>
>>>
>>>
>>> Am Mi., 19. Okt. 2022 um 22:34 Uhr schrieb Wolfgang Bangerth <
>>> bang...@colostate.edu>:
>>>
>>>> On 10/19/22 08:45, Simon Wiesheier wrote:
>>>> >
>>>> > What I want to do boils down to the following:
>>>> > Given the reference co-ordinates of a point 'p', along with the cell
>>>> on
>>>> > which 'p' lives,
>>>> > give me the value and gradient of a finite element function evaluated
>>>> at
>>>> > 'p'.
>>>> >
>>>> > My idea was to create a quadrature object with 'p' being the only
>>>> > quadrature point and pass this
>>>> > quadrature object to the FEValues object and finally do the
>>>> > .reinit(cell) call (then, of course, get_function_values()...)
>>>> > 'p' is different for all (2.5 million) quadrature points, which is
>>>> why I
>>>> > create the FEValues object so many times.
>>>>
>>>> It's worth pointing out that is exactly what
>>>> VectorTools::point_values()
>>>> does.
>>>>
>>>> (As others have already mentioned, if you want to do that many many
>>>> times over, this is too expensive and you should be using
>>>> FEPointEvaluation instead.)
>>>>
>>>> Best
>>>>   W.
>>>>
>>>> --
>>>> ------------------------------------------------------------------------
>>>> Wolfgang Bangerth          email:                 bang...@colostate.edu
>>>>                             www:
>>>> http://www.math.colostate.edu/~bangerth/
>>>>
>>>> --
>>>> 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+un...@googlegroups.com.
>>>> To view this discussion on the web visit
>>>> https://groups.google.com/d/msgid/dealii/cd1c8fa0-443d-b7bf-b433-f5ab033a247c%40colostate.edu
>>>> .
>>>>
>>> --
> 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.
> To view this discussion on the web visit
> https://groups.google.com/d/msgid/dealii/cb9c0d1a-77f5-4920-86dc-bc38ee6ab1dfn%40googlegroups.com
> <https://groups.google.com/d/msgid/dealii/cb9c0d1a-77f5-4920-86dc-bc38ee6ab1dfn%40googlegroups.com?utm_medium=email&utm_source=footer>
> .
>

-- 
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.
To view this discussion on the web visit 
https://groups.google.com/d/msgid/dealii/CAM50jEub3%3DYtqx09moMdCP-eJuhgPA-hjtqtPfmWtBFq%2Bd%3DyKw%40mail.gmail.com.

Reply via email to