Hi Jean-Paul,

Thank you so much for the deliberate suggestions! It seems the code-gallery example is a good starting point for me to implement St-Venant Kirchhoff hyperelastic model with the help of “push-forward” transformations. I need to dig into this example to understand how the geometric stiffness tangent is embedded for free.

 

It’s good to know a tutorial that uses one-filed elasticity will be available. I believe that will help a lot for beginners. AD is an interesting tool for nonlinear FEM. In  an early discussion, Alex shared his similar codes using AD at Integrated material and spatial traction forces on boundary not equal (google.com). But since I’m a deal.II beginner and not familiar with AD, I have not studied his code further.

 

So I am going to study your code-gallery example first as it is well documented. I may ask for your further help and will let you know if I get the code roll.

 

Thank you again for your wonderful help!

 

Best,

Michael

 

 

From: Jean-Paul Pelteret
Sent: Thursday, June 24, 2021 11:43 PM
To: dealii@googlegroups.com
Subject: Re: [deal.II] Questions on Step-18

 

HI Michael.

 

To be honest, I have a hard time to understand the three-filed formulation in Step 44. 

 

That’s a fair comment, and a known deficit in our current tutorial list. I hope that by the end of the year we will have a tutorial that uses one-field elasticity, and would be a much better starting point for you and others (WIP: https://github.com/dealii/dealii/pull/10394).

 

However, I’m still interested in implementing the pure geometrical large finite deformation model. I believe that will make things simpler and help me understand more complex nonlinear models. Meanwhile, I try to compare the results with a reference which only consider the geometrical nonlinearity.

 

It sounds to me like you’re looking to implement finite strain elasticity with a St Venant-Kirchhoff constitutive law. That is, as far as I understand, the simplest constitutive law for finite strain elasticity.

In this model, the elastic tangent looks like that for linear elasticity, but the strain tensor of interest is no longer the linear/small strain tensor, but rather the (nonlinear) Green-Lagrange strain tensor. (Nevertheless, even though the material law is linear, you still have the geometric tangent contribution that appears after linearisation of the residual, which comes from the linearisation of the now nonlinearly deformation-dependent test function “de” in  \int de : \sigma dv .)

 

If this is indeed what you want, then it is possible to swap out the Neo-Hookean constitutive law that is implemented in the code-gallery example with the St-Venant Kirchhoff constitutive law. You’d get “for free” the geometric stiffness contribution to the tangent, and would just have to worry about the material tangent and definition of the stress. I think that the only complication is the default parameterisation for the constitutive law — it was parameterised in terms of the left Cauchy-Green tensor as the whole problem is set up with the Kirchhoff stress (that’s the Cauchy stress, but per unit reference volume). Nevertheless, with the help of these “push-forward” transformations you would be able to use whatever parameterisation you’d like for the constitutive law (e.g. in terms of the Green-Lagrange strain tensor) and transform the resulting stress and material tangent to the correct (i.e. the spatial) configuration as is required for that implementation of the weak forms. The important thing to remember is that, despite all of the transformations involved, you’re still solving the same BVP with the same constitutive law irrespective of which choice of stress-strain (energetic) conjugate pairs you use, and irrespective of whether its expressed in the total or updated Lagrangian form.

 

But anyway, this is just a suggestion. Irrespective of whether or not you choose to go that route, I wish you luck in implementing what you’re needing to.

 

Best,

Jean-Paul



On 25. Jun 2021, at 02:32, Michael Li <lianxi...@gmail.com> wrote:

 

Andrew, thanks for confirming that. The missing 1/2 does not affect the demonstration of functionalities of deal.II but it may change the results.

 

Jean-Paul, thanks for commenting on my second question. I want to study the pure geometrical nonlinear elasticity (large deformation with linear material). It should be the simplest nonlinear model in elasticity. Step 18 looks like a good start as an updated Lagrangian formulation; but it does not include the nonlinear part of the finite strain. I spent some on Step 44 because “ Quasi-Static Finite-Strain Compressible Elasticity” also mentions it is based on Step 44. To be honest, I have a hard time to understand the three-filed formulation in Step 44. So maybe I should get into “Quasi-Static … Elasticity” first which uses classical one-field formulation. However, I’m still interested in implementing the pure geometrical large finite deformation model. I believe that will make things simpler and help me understand more complex nonlinear models. Meanwhile, I try to compare the results with a reference which only consider the geometrical nonlinearity.

 

-Michael

 

 

From: Andrew McBride
Sent: Thursday, June 24, 2021 2:23 PM
To: deal.II User Group
Subject: Re: [deal.II] Questions on Step-18

 

Hi both,

 

Jean-Paul has addressed the second point nicely. On the first point, I think there is a 1/2 missing. The curl of the velocity gradient is the vorticity which is twice the angular velocity - hence I think you need a 1/2. Happy to be corrected on this.

 

Best,

Andrew




On 24 Jun 2021, at 21:03, Jean-Paul Pelteret <jppelte...@gmail.com> wrote:

 

Hi Michael,

 

I cannot comment on the first question, but might be able to assist a bit with the second. But may I first ask, what precisely are you trying to achieve with this extension? 

 

As interesting as it is, in the past I had found step-18 to deviate too significantly from the “classical” approach to elasticity to be a natural candidate extend towards finite strain elasticity, for example (this is kind-of implied by the caveat that you partially quoted). I’ve been looking at some of my textbooks (e.g. reviewing the topic of the updated Lagrangian formulation in Holzapfel’s “Nonlinear solid mechanics” and Wrigger’s “Nonlinear finite element methods”) to try to answer (2), but cannot trivially reconcile the two approaches. I think that there’s a bit too much going on to be able to correctly deduce by eye what the requisite changes are. I think that you’d need to reformulate the balance laws and consider their implications for the implemented weak forms — in particular, I think that you’d be missing a contribution that (for finite deformation) looks like the geometric stiffness, but there could be further differences that extend from the overall approach taken to the problem. 

 

If you’re interested in an examples that are more closely aligned with what you might see in the literature, then you can take a look at step-44 or the “Quasi-Static Finite-Strain Compressible Elasticity” code-gallery example, which is effectively step-44 reduced to the one-field total Lagrangian formulation that you’d find in many standard textbooks. It would be easy enough to rework this to use the updated Lagrangian approach, if that it what you desire.

 

I hope that this helps you a little.

 

Best,

Jean-Paul

 




On 22. Jun 2021, at 15:24, Michael Lee <lianxi...@gmail.com> wrote:

 

Hello,

I have two questions when studying Step-18. 

 

1) Should there be a factor 1/2 when calculating the rotation matrix angle (code in tutorial: angle = std::atan(curl) ) ? This is a mathematical problem. Can anyone give some material for the formula to clear out my doubt?

 

2) In the introduction, it says "we will consider a model in which we produce a small deformation, deform the physical coordinates of the body by this deformation, and then consider the next loading step again as a linear problem. This isn't consistent, since the assumption of linearity implies that deformations are infinitesimal and so moving around the vertices of our mesh by a finite amount before solving the next linear problem is an inconsistent approach." My question is how to make it consistent. Can I just add the nonlinear part of the strain tensor like the following (of course in both get_strain functions)? I just want to consider the geometrical nonlinearity.

 

  template <int dim>

  inline SymmetricTensor<2, dim>

  get_strain(const std::vector<Tensor<1, dim>> &grad)

  {

    Assert(grad.size() == dim, ExcInternalError());

 

    SymmetricTensor<2, dim> strain;

    for (unsigned int i = 0; i < dim; ++i)

      strain[i][i] = grad[i][i];

 

    for (unsigned int i = 0; i < dim; ++i)

      for (unsigned int j = i + 1; j < dim; ++j)

        strain[i][j] = (grad[i][j] + grad[j][i]       // linear part

                              grad[i][0] * grad[j][0] +    //nonlinear part

                              grad[i][1] * grad[j][1] +   

                              grad[i][2] * grad[j][2]) / 2;

 

    return strain;

  }

 

Thanks,

Michael

 

-- 
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/02f465df-3cab-4d84-bd48-92ca4a7981efn%40googlegroups.com.

 

 

-- 
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/C2745EC6-0BED-4AFC-A21C-C003D851C869%40gmail.com.

 

-- 
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/D7696104-6C14-43F2-8B03-CF7AA0E89391%40gmail.com.

 

 

-- 
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/5A41FB64-E568-44C2-82CA-2A7CAE923CD4%40hxcore.ol.

 

--
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/020CCBE9-E8BB-40AC-A205-41E8E2860553%40gmail.com.

 

--
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/56766414-40CB-47D3-A851-EB9583263D80%40hxcore.ol.

Reply via email to