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 <[email protected]> wrote:
>
> Dear delii members
>
> While modifying the material model in the 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 [email protected]
> <mailto:[email protected]>.
> 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/
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 [email protected].
To view this discussion on the web visit
https://groups.google.com/d/msgid/dealii/111DFCE4-AA04-46BD-A091-481B15CDE17D%40gmail.com.