Hi,

I am trying to use symbolic differentiation to obtain the Kirchhoff stress 
tensor using the material model given in step 44. I have attached the 
minimal example herewith, which reproduces the following error:

In line 63, the differentiation of the energy function (psi_sd) with 
respect to left cauchy green tensor (b_sd) results in  the Kirchhoff stress 
tensor tau_sd = zero tensor. I believe this is due to the fact that tau_sd 
is stored as a function of F_sd (symbolic deformation gradient) and is 
therefore unable to recognise the dependence with b_sd (symbolic left 
cauchy green tensor).

How can i resolve this issue? Thanks in advance.

Regards
Vinayak

-- 
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/61d019b1-626b-4c3c-9c17-b5c323d1e77en%40googlegroups.com.
//-----------------------------------------------------------
//
//    Copyright (C) 2020 by the deal.II authors
//
//    This file is subject to LGPL and may not be distributed
//    without copyright and license information. Please refer
//    to the file deal.II/doc/license.html for the  text  and
//    further information on this license.
//
//-----------------------------------------------------------

#include <deal.II/base/symmetric_tensor.h>
#include <deal.II/base/tensor.h>

#include <deal.II/differentiation/sd.h>

#include <deal.II/physics/elasticity/kinematics.h>
#include <deal.II/physics/elasticity/standard_tensors.h>

using namespace dealii;

int
main()
{
  const Differentiation::SD::Expression nu_e_sd("nu_e");
  const Differentiation::SD::Expression lambda_sd("lambda");
  const Differentiation::SD::Expression mu_e_sd("mu_e");
  const Differentiation::SD::Expression kappa_e_sd("kappa_e");
  const Differentiation::SD::Expression c_1_e_sd("c_1_e");

  const int dim      = 3;
  const int spacedim = 3;

  const Tensor<2, dim, Differentiation::SD::Expression> F_sd(
    Differentiation::SD::make_tensor_of_symbols<2, dim>("F"));
  std::cout << F_sd << std::endl;

  const Differentiation::SD::Expression det_F_sd = determinant(F_sd);

  const Tensor<2, 3, dealii::Differentiation::SD::Expression> F_bar_sd =
    std::pow(det_F_sd, -1 / 3) * F_sd;

  SymmetricTensor<2, spacedim, Differentiation::SD::Expression> b_bar_sd =
    symmetrize(F_bar_sd * transpose(F_bar_sd));
  SymmetricTensor<2, spacedim, Differentiation::SD::Expression> b_sd =
    std::pow(det_F_sd, 1 / 3) * b_bar_sd;
  std::cout << "b_sd:  " << b_sd << std::endl;


  const Differentiation::SD::Expression psi_vol_sd =
    (kappa_e_sd / 4.0) *
    (det_F_sd * det_F_sd - 1.0 - 2.0 * Differentiation::SD::log(det_F_sd));

  const Differentiation::SD::Expression psi_iso_sd =
    c_1_e_sd * (trace(b_bar_sd) - dim);

  Differentiation::SD::Expression psi_sd;

  psi_sd = psi_vol_sd + psi_iso_sd;

  SymmetricTensor<2, dim, Differentiation::SD::Expression> tau_sd;

  tau_sd = 2 * b_sd * Differentiation::SD::differentiate(psi_sd, b_sd);

  Tensor<2, dim, Differentiation::SD::Expression> first_PK_sd;
  first_PK_sd = 2 * b_sd * Differentiation::SD::differentiate(psi_sd, F_sd);

  std::cout << "psi_sd:  " << psi_sd << std::endl;
  std::cout << "first_PK_sd:  " << first_PK_sd << std::endl;

  std::cout << "tau_sd:  " << tau_sd << std::endl;
}

Reply via email to