Hi Simon,

Sorry that I didn't get back to you on this before now. I think that I get 
the gist of what you're trying to do, but I must admit that I'm a bit 
skeptical as to whether or not its going to work (the AD framework was not 
really conceptualised to be used in this manner, so it will be interesting 
to hear whether you end up finding success in this endeavor). Basically, 
what you need to do is somehow link the degrees of freedom in the found 
cell with those that you created C_ad with.  I guess suggests that you need 
to express how these values might change as your invariant point changes. 
If you *somehow* manage to form an association between these two sets of 
DoFs then you can use one of the various FEValuesViews 
get_function_gradients_from_local_dof_values() 
<https://dealii.org/current/doxygen/deal.II/classFEValuesViews_1_1Scalar.html#abdc8d5d7b846cc7e0d1e15e8b0cf66b5>
 
functions to compute gradients that have the required sensitivities encoded 
in them. This is actually what we already do in ScratchData 
<https://github.com/dealii/dealii/blob/master/include/deal.II/meshworker/scratch_data.h#L1712-L1713>
:: 
<https://github.com/dealii/dealii/blob/master/include/deal.II/meshworker/scratch_data.h#L1712-L1713>
get_gradients() 
<https://github.com/dealii/dealii/blob/master/include/deal.II/meshworker/scratch_data.h#L1712-L1713>,
 
which makes this class AD-capable. 

I hope that this helps, if you haven't succeeded with this already.

Best,
Jean-Paul

On Tuesday, July 19, 2022 at 11:07:50 AM UTC+2 Simon wrote:

> Hi Jean-Paul,
>
> thanks for pointing that out.
>
> However, I believe there has been a misunderstanding:
> I am working here with a FE discretized energy on a seperate 
> triangulation, say, 'triangulation_energy', with the coordinates being the 
> first and third invariant of 'C'; The corresponding nodal energy values are 
> stored in the vector 'solution_energy' and are computed once based on a 
> analytical energy function. 
> So for each quadrature point in my geometry (first triangulation) I am 
> computing 'C=F^T*F', thence the first and third invariant of 'C' and store 
> them in a Point<2> object denoted as 'invariant_point'. 
> Then, I employ 
> 'GridTools::find_active_cell_around_point(triangulation_energy, 
> invariant_point)' to finally create a second FEValues object with the only 
> purpose to compute the gradient of the energy ('solution_energy') at  
> 'invariant_point'.
> 'get_function_gradients(solution_energy, grad_energy_at_ref_point ) gives 
> me the gradient of 'solution_energy' with respect to the first and third 
> invariant of 'C' because the real coordinates on the associated 
> triangulation are the first and third invariant of 'C'.
> Now, I have all the quantities to compute the stress tensor as
> S_ad = 2*(grad_energy_at_ref_point[0][0]*invariant1_wrt_C_ad 
>                                                                 + 
> grad_energy_at_ref_point[0][1]*invariant3_wrt_C_ad);
>
> My goal is to compute the 4th order tangent \frac{\partial S_ad}{\partial 
> C_ad}.
>
> All that said, is it not correct to start like this?:
> SymmetricTensor<2,dim> C = symmetrize(F^t *F);
> ad_helper.register_independent_variable(C, C_dofs);
> const SymmetricTensor<2, dim, ADNumberType> C_ad = 
> ad_helper.get_sensitive_variables(C_dofs);
>
> Best,
> Simon
>
>
>
> Am Mo., 18. Juli 2022 um 22:46 Uhr schrieb Jean-Paul Pelteret <
> jppel...@gmail.com>:
>
>> Hi Simon,
>>
>> Unfortunately I don’t have the time at this moment to give you a full 
>> explanation as to why the logic of your code is wrong, but in essence you 
>> have the sequence of operations inverted. You need to compute your energy 
>> based on the “sensitive” DoF values that would come from the reinit’d 
>> FEValues operation. You cannot compute something like C_ad ahead of time, 
>> like it appears that you have.
>>
>> Take a look at these lines in the (unfinished) tutorial that will 
>> illustrate the steps required to linearise an energy functional. This 
>> should hopefully help steer you in the right direction.
>>
>> https://github.com/dealii/dealii/pull/10394/files#diff-fc8b83d1370bfd7eb558ea76175bfc0e8d6305023d54b17ec9cccb0fba9b1e02R1758-R1831
>>
>> Best,
>> Jean-Paul
>>
>> On 18. Jul 2022, at 19:46, Simon <simon.w...@gmail.com> wrote:
>>
>> Dear all,
>>
>> I want to retrieve a tangent at quadrature point level using AD - 
>> attached is a snippet of my code:
>>
>> The computation of the dependent variable 'S_ad' in line 36 involves the 
>> call 
>> get_function_gradients (solution_energy, grad_energy_at_ref_point), 
>> see line 35. 
>> 'solution_energy' is a scalar function which depends on the independent 
>> variable 'C_ad' as can be seen in the derivations from line 8 on.
>> My issue is that I see no way to pass the NumberTypes type 
>> (=Differentiation::AD::NumberTypes::sacado_dfad_dfad) to the corrosponding 
>> FEValues object, and, as a consequence, the AD framework is not aware that 
>> 'solution_energy' is a function of 'C_ad'. 
>> In particular, 'grad_energy_at_ref_point' will be treated as a constant 
>> with respect to 'C_ad' in line 36f. 
>>
>> How can I incorporate this dependency in the AD framework?
>>
>> Thank you,
>> Simon
>>
>> -- 
>> 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/990f7f7f-37a2-4f9d-b6c2-bbf7dda21e8dn%40googlegroups.com
>>  
>> <https://groups.google.com/d/msgid/dealii/990f7f7f-37a2-4f9d-b6c2-bbf7dda21e8dn%40googlegroups.com?utm_medium=email&utm_source=footer>
>> .
>> <Code_AD.PNG>
>>
>>
>> -- 
>> 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/D175C472-E708-4129-A15E-3FA67666AE9A%40gmail.com
>>  
>> <https://groups.google.com/d/msgid/dealii/D175C472-E708-4129-A15E-3FA67666AE9A%40gmail.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/db66463a-0283-4808-bb38-8004511aad91n%40googlegroups.com.

Reply via email to