Dear Simon,

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? The issue you are facing is that FEValues that you are using is
using a very abstract entry point that does precomputations that only pay
off if using the unit points many times. And even in the case of the same
unit points it is not really fast, it is a general-purpose baseline that I
would not recommend for high-performance purposes.

As a final note, I would mention that FEPointEvaluation falls back to
FEValues for complicated FiniteElement types, so it might be that you do
not get speedups in those cases. But we could work on it if you need it,
today we know much better what to do than a few years ago.

Best,
Martin

On Wed, 19 Oct 2022, 16:45 Simon Wiesheier, <simon.wieshe...@gmail.com>
wrote:

> " It's an environment variable. "
>
> I did
> $DEAL_II_NUM_THREADS
> and the variable is not set.
> But if it were set to one, why would this explain the gap between cpu and
> wall time?
>
> " My point is the constructor should not be called millions of times. You
> are not going to be able to get that function 100 times faster. It's best
> to find a way to call it less often. "
>
> 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.
>
> Do you a different suggestion to solve my problem, ie to evaluate the
> finite element field and its derivatives at 'p'?
>
> Best,
> Simon
>
>
> Am Mi., 19. Okt. 2022 um 16:17 Uhr schrieb Bruno Turcksin <
> bruno.turck...@gmail.com>:
>
>> Simon,
>>
>> Le mer. 19 oct. 2022 à 09:33, Simon Wiesheier <simon.wieshe...@gmail.com>
>> a écrit :
>>
>>> Thank you for your answer!
>>>
>>> " Did you set DEAL_II_NUM_THREADS=1?"
>>>
>>> How can I double-check that?
>>> ccmake .
>>> only shows my the variables CMAKE_BUILD_TYPE and deal.II_DIR .
>>> But I do  do knot if this is the right place to look for.
>>>
>> It's an environment variable. If you are using bash, you can do
>>
>> export DEAL_II_NUM_THREADS=1
>>
>>
>>>
>>> " That could explain why CPU and Wall time are different. Finally, if I
>>> understand correctly, you are calling the constructor of FEValues about 2.5
>>> million times. That means that the call to one FEValues constructor is
>>> 100/2.5e6 seconds about 40 microseconds. That doesn't seem too slow. "
>>>
>>> There was a typo in my post. It should be 160/2.5e6 seconds about 64
>>> microsecends.
>>>
>> My point is the constructor should not be called millions of times. You
>> are not going to be able to get that function 100 times faster. It's best
>> to find a way to call it less often.
>>
>> Best,
>>
>> Bruno
>>
>> --
>> 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/CAGVt9eMfVohOUToQOsBD_v%2BqU%3D0Em_XOMiwqFi2SM_0zLoy-sQ%40mail.gmail.com
>> <https://groups.google.com/d/msgid/dealii/CAGVt9eMfVohOUToQOsBD_v%2BqU%3D0Em_XOMiwqFi2SM_0zLoy-sQ%40mail.gmail.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/CAM50jEtyY576riC6yNqqMafXfGGvTXY8mhm%3Di7HMzr-U_LAxbQ%40mail.gmail.com
> <https://groups.google.com/d/msgid/dealii/CAM50jEtyY576riC6yNqqMafXfGGvTXY8mhm%3Di7HMzr-U_LAxbQ%40mail.gmail.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/CAG%3Dko9v8uURkxwW0Q8vLLpf8fXzQ_JE-6FwvP1hikNhPF6seyA%40mail.gmail.com.

Reply via email to