Dear Navneet,

Great, I’m glad to hear that you were able to reconcile the two sets of 
equations. Please feel free to post any other questions that you may have to 
the forum. We’ll do my best to answer them!

Best,
Jean-Paul

> On 23 Apr 2020, at 12:36, navneet roshan <navneetrosh...@gmail.com> wrote:
> 
> Dear Jean-Paul
> 
> Thank you for the reply and suggestions I went through the modifications, I 
> could convert the step-44 code to respective one field and two field code.
> the results of two-field and the three-field(step-44) are matching closely.  
> I have one more question but it will take the discussion beyond the scope of 
> the question I asked earlier. I will ask that later. Thank you once again.
>  
> 
> Navneet R
> 
> 
> 
> On Wed, Mar 18, 2020 at 3:56 AM Jean-Paul Pelteret <jppelte...@gmail.com 
> <mailto:jppelte...@gmail.com>> wrote:
> Dear Navneet,
> 
> Yes, and no. The “p”s in the final result that you state are from the 
> independent field, but the “J” is actually the volumetric Jacobian, not the 
> independent field!
> 
> Remember that in the 3-field formulation this term no longer comes from the 
> linearisation of the constitutive law, but rather the linearisation of a 
> component of the residual associated with the displacement degrees of 
> freedom. In the end one can identify similarities in the structure of term 
> for the 1- and 3-field formulations, but I don't think that one should try to 
> move directly from the 1-field result to the 3-field result in the way that 
> you wish to. I encourage you to read the books and papers that are given in 
> the reference list for step-44. They might be able to highlight some key 
> parts of the (somewhat tedious) derivation that we didn’t cover in sufficient 
> detail in the tutorial documentation. 
> 
> To hyper-summarise what those references will show, if one pulls everything 
> back into the reference configuration then we end up with a displacement 
> residual contribution that looks something like this:
> R_{u} = \int \delta E : [S^{iso} + \tilde{p} J C^{-1}] dV
> where \delta E is the variation of the Green-Lagrange strain tensor, S^{iso} 
> is the isochoric part of the Piola-Kirchhoff stress tensor (derived from some 
> constitutive law), \tilde{p} is the independent pressure response field, J = 
> det(F) is the volumetric Jacobian, and C is the right Cauchy-Green 
> deformation tensor. This tangent contribution that you are seeking comes from 
> linearising the "\tilde{p} J C^{-1}” part with respect to the displacement, 
> i.e.
> \Delta_{u} R_{u} = \int \delta E : [H^{iso} + \tilde{p} 
> \frac{\partial}{\partial E}[J C^{-1}] ] : \Delta E dV + ….
> and pushing all of the tangents forward into the current configuration again.
> 
> Best,
> Jean-Paul
> 
>> On 06 Mar 2020, at 14:40, navneet roshan <navneetrosh...@gmail.com 
>> <mailto:navneetrosh...@gmail.com>> wrote:
>> 
>> Dear Jean-Paul,
>> 
>> Thank you so much for your timely and detailed response. If I repeat what I 
>> understood is
>> following, in the given expression    with  for the case of
>>  three variable formulation  and  are independent hence  , as .
>>  
>> Regards,
>> Navneet R
>> 
>> 
>> 
>> 
>> On Fri, Mar 6, 2020 at 2:47 AM Jean-Paul Pelteret <jppelte...@gmail.com 
>> <mailto:jppelte...@gmail.com>> wrote:
>> Dear Naveet,
>> 
>> I think that I see where things are not explained completely clearly. The 
>> sections “Hyperelastic materials”, “Neo-Hookean materials”, and “Elasticity 
>> tensors” all describe the theory for compressible materials. This is 
>> specifically mentioned in the “Neo-Hookean materials” part, and was probably 
>> done to give context for the modelling of quasi-incompressible models 
>> (although I don’t remember this with 100% clarity). So the “p” and “\hat{p}” 
>> that are mentioned there are the pressure-like terms that can be computed 
>> directly from the constitutive law for compressible materials. These, 
>> however, no longer apply for the 3-field formulation.
>> 
>> The strain energy function used in the section “Principle of stationary 
>> potential energy and the three-field formulation” (which derives the 
>> governing equations as implemented in the code) has the same form as that 
>> before, but different arguments. Specifically, the volumetric part now 
>> depends on the dilatation field “\tilde{J}" rather than the point-wise 
>> volumetric Jacobian “J". Similarly, the pressure response “p" is no longer 
>> computed from the energy density function, but is a field variable 
>> “\tilde{p}" itself. How these relate to the original volumetric/isochoric 
>> stresses and tangents is captured in the underpriced parts of the equations 
>> in “Principle of stationary potential energy and the three-field 
>> formulation”, namely the residual equation and its consistent linearisation.
>> 
>> What you had written is indeed correct for one-field elasticity. There is a 
>> code-gallery example (one-field elasticity) that use this exact expression 
>> for the volumetric material tangent:
>> https://github.com/dealii/code-gallery/blob/master/Quasi_static_Finite_strain_Compressible_Elasticity/cook_membrane.cc#L654-L666
>>  
>> <https://github.com/dealii/code-gallery/blob/master/Quasi_static_Finite_strain_Compressible_Elasticity/cook_membrane.cc#L654-L666>
>> and another 3-field elasticity example that uses the definition as given in 
>> step-44:
>> https://github.com/dealii/code-gallery/blob/master/Quasi_static_Finite_strain_Quasi_incompressible_ViscoElasticity/viscoelastic_strip_with_hole.cc#L561-L568
>>  
>> <https://github.com/dealii/code-gallery/blob/master/Quasi_static_Finite_strain_Quasi_incompressible_ViscoElasticity/viscoelastic_strip_with_hole.cc#L561-L568>
>> 
>> So, in summary, the definitions of the volumetric stress and tangent differ 
>> to one-field elasticity due to the introduction of the additional fields. 
>> This follows directly from the expression of the stationary point and 
>> linearisation, rather than from the constitutive law in the code. However, 
>> one can interpret various terms in the residual and linearisation as the 
>> equivalent of the volumetric stress and tangent (even though they only 
>> indirectly governed by the constitutive law, since they yield the pressure 
>> response to the added primary field). That’s what we tried to capture in the 
>> code.
>> 
>> I hope that this helps clarify any misunderstanding. Suggestions as to how 
>> this can be made more clear in the introductory text are always welcome, so 
>> please feel to open a pull request for them. 
>> 
>> Best,
>> Jean-Paul
>> 
>>> On 05 Mar 2020, at 12:51, navneet roshan <navneetrosh...@gmail.com 
>>> <mailto:navneetrosh...@gmail.com>> wrote:
>>> 
>>> Dear delii members
>>> 
>>> While modifying the material model in the step-44.cc <http://step-44.cc/>, 
>>> implementation of the function get_Jc_vol(). I found the implementation in 
>>> the code is different than the formula provided. The worse part is that the 
>>> solution stops converging after correcting the implementation. I will be 
>>> really grateful if some one can hint me as, I am stuck between convergence 
>>> and divergence of material models from quite some time.
>>> 
>>>    The given implementation is:
>>>                SymmetricTensor<4, dim> 
>>> <https://www.dealii.org/current/doxygen/deal.II/classSymmetricTensor.html> 
>>> get_Jc_vol() const
>>> {
>>>    return p_tilde * det_F *
>>>          (Physics::Elasticity::StandardTensors<dim>::IxI 
>>> <https://www.dealii.org/current/doxygen/deal.II/classPhysics_1_1Elasticity_1_1StandardTensors.html>
>>>  -
>>>          (2.0 * Physics::Elasticity::StandardTensors<dim>::S 
>>> <https://www.dealii.org/current/doxygen/deal.II/classPhysics_1_1Elasticity_1_1StandardTensors.html>))
>>> }
>>> 
>>> The implementation as per formula, should have been:
>>> 
>>>              SymmetricTensor<4, dim> get_Jc_vol() const
>>>              {
>>>                      return det_F  * (  (get_dPsi_vol_dJ() + 
>>> det_F*get_d2Psi_vol_dJ2() )
>>>                        *Physics::Elasticity::StandardTensors<dim>::IxI
>>>                       - (2.0 * get_dPsi_vol_dJ() 
>>> *Physics::Elasticity::StandardTensors<dim>::S) );
>>>     }
>>> 
>>> Thank you, 
>>> 
>>> Regards,
>>> Navneet R
>>> 
>>> 
>>> -- 
>>> The deal.II project is located at http://www.dealii.org/ 
>>> <http://www.dealii.org/>
>>> For mailing list/forum options, see 
>>> https://groups.google.com/d/forum/dealii?hl=en 
>>> <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>.
>>> To view this discussion on the web visit 
>>> https://groups.google.com/d/msgid/dealii/CAK9McD0Rm7GotVhecV4%3D2ExxGcZPUtPLVQNVi652skkXYG0mOQ%40mail.gmail.com
>>>  
>>> <https://groups.google.com/d/msgid/dealii/CAK9McD0Rm7GotVhecV4%3D2ExxGcZPUtPLVQNVi652skkXYG0mOQ%40mail.gmail.com?utm_medium=email&utm_source=footer>.
>> 
>> 
>> -- 
>> The deal.II project is located at http://www.dealii.org/ 
>> <http://www.dealii.org/>
>> For mailing list/forum options, see 
>> https://groups.google.com/d/forum/dealii?hl=en 
>> <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>.
>> To view this discussion on the web visit 
>> https://groups.google.com/d/msgid/dealii/111DFCE4-AA04-46BD-A091-481B15CDE17D%40gmail.com
>>  
>> <https://groups.google.com/d/msgid/dealii/111DFCE4-AA04-46BD-A091-481B15CDE17D%40gmail.com?utm_medium=email&utm_source=footer>.
>> 
>> -- 
>> The deal.II project is located at http://www.dealii.org/ 
>> <http://www.dealii.org/>
>> For mailing list/forum options, see 
>> https://groups.google.com/d/forum/dealii?hl=en 
>> <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>.
>> To view this discussion on the web visit 
>> https://groups.google.com/d/msgid/dealii/CAK9McD1dTA6%2BkeVhaCH9BjGdpjph--SH-ek%3D5GC2FWUip9cGvg%40mail.gmail.com
>>  
>> <https://groups.google.com/d/msgid/dealii/CAK9McD1dTA6%2BkeVhaCH9BjGdpjph--SH-ek%3D5GC2FWUip9cGvg%40mail.gmail.com?utm_medium=email&utm_source=footer>.
> 
> 
> -- 
> The deal.II project is located at http://www.dealii.org/ 
> <http://www.dealii.org/>
> For mailing list/forum options, see 
> https://groups.google.com/d/forum/dealii?hl=en 
> <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>.
> To view this discussion on the web visit 
> https://groups.google.com/d/msgid/dealii/71E8068B-C0EE-43CA-9A91-FED46C50F320%40gmail.com
>  
> <https://groups.google.com/d/msgid/dealii/71E8068B-C0EE-43CA-9A91-FED46C50F320%40gmail.com?utm_medium=email&utm_source=footer>.
> 
> -- 
> The deal.II project is located at http://www.dealii.org/ 
> <http://www.dealii.org/>
> For mailing list/forum options, see 
> https://groups.google.com/d/forum/dealii?hl=en 
> <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>.
> To view this discussion on the web visit 
> https://groups.google.com/d/msgid/dealii/CAK9McD1urJT4U%3DN-JdXp16OjYBARkVZXV8U9TLxHfGUahMxroQ%40mail.gmail.com
>  
> <https://groups.google.com/d/msgid/dealii/CAK9McD1urJT4U%3DN-JdXp16OjYBARkVZXV8U9TLxHfGUahMxroQ%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/6CD1A697-2E64-440C-9FD8-851B19E3DFB9%40gmail.com.

Reply via email to