sfepy.terms.terms_jax module

Proof-of-concept JAX-based terms supporting automatic differentiation.

class sfepy.terms.terms_jax.LinearElasticLADTerm(name, arg_str, integral, region, **kwargs)[source]

Homogeneous isotropic linear elasticity term differentiable w.r.t. material parameters \lambda, \mu (Lamé’s parameters).

Definition:

\int_{\Omega} D_{ijkl}\ e_{ij}(\ul{v}) e_{kl}(\ul{u})\\ \mbox{ with } \\
D_{ijkl} = \mu (\delta_{ik} \delta_{jl}+\delta_{il} \delta_{jk}) +
\lambda \ \delta_{ij} \delta_{kl}

Call signature:

dw_lin_elastic_l_ad

(material_1, material_2, virtual, state)

Arguments:
  • material_1: \lambda (Lamé’s first parameter)

  • material_2: \mu (Lamé’s second parameter, shear modulus)

  • virtual/parameter_1: \ul{v}

  • state/parameter_2: \ul{u}

arg_shapes = {'material_1': '1, 1', 'material_2': '1, 1', 'state': 'D', 'virtual': ('D', 'state')}
arg_types = (('material_1', 'material_2', 'virtual', 'state'),)
diff_info = {'material_1': 1, 'material_2': 1}
static function(out, fun, fargs)[source]
get_fargs(material1, material2, virtual, state, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
modes = ('weak',)
name = 'dw_lin_elastic_l_ad'
class sfepy.terms.terms_jax.LinearElasticYPADTerm(name, arg_str, integral, region, **kwargs)[source]

Homogeneous isotropic linear elasticity term differentiable w.r.t. material parameters E (Young’s modulus), \nu (Poisson’s ratio).

Definition:

\int_{\Omega} D_{ijkl}\ e_{ij}(\ul{v}) e_{kl}(\ul{u})\\ \mbox{ with } \\
D_{ijkl} = \mu (\delta_{ik} \delta_{jl}+\delta_{il} \delta_{jk}) +
\lambda \ \delta_{ij} \delta_{kl}, \\ \mbox{ where } \\
\lambda = E \nu / ((1 + \nu)(1 - 2\nu)), \\ \mu = E / 2(1 + \nu)

Call signature:

dw_lin_elastic_yp_ad

(material_1, material_2, virtual, state)

Arguments:
  • material_1: E

  • material_2: \nu

  • virtual/parameter_1: \ul{v}

  • state/parameter_2: \ul{u}

arg_shapes = {'material_1': '1, 1', 'material_2': '1, 1', 'state': 'D', 'virtual': ('D', 'state')}
arg_types = (('material_1', 'material_2', 'virtual', 'state'),)
diff_info = {'material_1': 1, 'material_2': 1}
static function(out, fun, fargs)[source]
get_fargs(material1, material2, virtual, state, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
modes = ('weak',)
name = 'dw_lin_elastic_yp_ad'
class sfepy.terms.terms_jax.MassADTerm(name, arg_str, integral, region, **kwargs)[source]

Homogeneous mass term differentiable w.r.t. the material parameter.

Definition:

\int_{\cal{D}} \rho \ul{v} \cdot \ul{u}

Call signature:

dw_mass_ad

(material, virtual, state)

Arguments:
  • material_1: \rho

  • virtual: \ul{v}

  • state: \ul{u}

arg_shapes = {'material': '1, 1', 'state': 'D', 'virtual': ('D', 'state')}
arg_types = (('material', 'virtual', 'state'),)
diff_info = {'material': 1}
static function(out, fun, fargs)[source]
get_fargs(material_density, virtual, state, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
modes = ('weak',)
name = 'dw_mass_ad'
class sfepy.terms.terms_jax.NeoHookeanTLADTerm(name, arg_str, integral, region, **kwargs)[source]

Homogeneous Hyperelastic neo-Hookean term differentiable w.r.t. the material parameter. Effective stress S_{ij} = \mu
J^{-\frac{2}{3}}(\delta_{ij} - \frac{1}{3}C_{kk}C_{ij}^{-1}).

Definition:

\int_{\Omega} S_{ij}(\ul{u}) \delta E_{ij}(\ul{u};\ul{v})

Call signature:

dw_tl_he_neohook_ad

(material, virtual, state)

Arguments:
  • material : \mu

  • virtual : \ul{v}

  • state : \ul{u}

arg_shapes = {'material': '1, 1', 'state': 'D', 'virtual': ('D', 'state')}
arg_types = ('material', 'virtual', 'state')
diff_info = {'material': 1}
static function(out, fun, fargs)[source]
geometries = ['2_3', '2_4', '3_4', '3_8']
get_fargs(material, virtual, state, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
modes = ('weak',)
name = 'dw_tl_he_neohook_ad'
class sfepy.terms.terms_jax.OgdenTLADTerm(name, arg_str, integral, region, **kwargs)[source]

Homogeneous hyperelastic Ogden model term differentiable w.r.t. the material parameters, with the strain energy density

W = \frac{\mu}{\alpha} \, \left(
    \bar\lambda_1^{\alpha} + \bar\lambda_2^{\alpha}
     + \bar\lambda_3^{\alpha} - 3 \right) \; ,

where \lambda_k, k=1, 2, 3 are the principal stretches, whose squares are the principal values of the right Cauchy-Green deformation tensor \mathbf{C}. For more details see OgdenTLTerm.

WARNING: The current implementation fails to compute the tangent matrix when \mathbf{C} has multiple eigenvalues (e.g. zero deformation). In that case nans are returned, as a result of dividing by zero. See [1], Section 11.2.3, page 385.

[1] Borst, R. de, Crisfield, M.A., Remmers, J.J.C., Verhoosel, C.V., 2012. Nonlinear Finite Element Analysis of Solids and Structures, 2nd edition. ed. Wiley, Hoboken, NJ.

Definition:

\int_{\Omega} S_{ij}(\ul{u}) \delta E_{ij}(\ul{u};\ul{v})

Call signature:

dw_tl_he_ogden_ad

(material_mu, material_alpha, virtual, state)

Arguments:
  • material_1 : \mu

  • material_2 : \alpha

  • virtual : \ul{v}

  • state : \ul{u}

arg_shapes = {'material_alpha': '1, 1', 'material_mu': '1, 1', 'state': 'D', 'virtual': ('D', 'state')}
arg_types = ('material_mu', 'material_alpha', 'virtual', 'state')
diff_info = {'material_alpha': 1, 'material_mu': 1}
static function(out, fun, fargs)[source]
geometries = ['3_4', '3_8']
get_fargs(material_mu, material_alpha, virtual, state, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
modes = ('weak',)
name = 'dw_tl_he_ogden_ad'
sfepy.terms.terms_jax.ceval_elasticity_l(lam, mu, vbfg, ubfg, det, cu)[source]
sfepy.terms.terms_jax.ceval_elasticity_yp(young, poisson, plane, vbfg, ubfg, det, cu)[source]
sfepy.terms.terms_jax.ceval_mass(density, vbf, ubf, det, cu)[source]
sfepy.terms.terms_jax.ceval_neohook(mu, vbfg, ubfg, det, cu)[source]
sfepy.terms.terms_jax.ceval_neohook0(mu, vbfg, ubfg, det, cu)[source]
sfepy.terms.terms_jax.ceval_ogden(mu, alpha, vbfg, ubfg, det, cu)[source]
sfepy.terms.terms_jax.eval_alpha_ogden(mu, alpha, vbfg, ubfg, det, cu)

Vectorized version of ceval_ogden. Takes similar arguments as ceval_ogden but with additional array axes over which ceval_ogden is mapped.

Original documentation:

Jacobian of ceval_ogden with respect to positional argument(s) 1. Takes the same arguments as ceval_ogden but returns the jacobian of the output with respect to the arguments at positions 1.

sfepy.terms.terms_jax.eval_density_mass(density, vbf, ubf, det, cu)

Vectorized version of ceval_mass. Takes similar arguments as ceval_mass but with additional array axes over which ceval_mass is mapped.

Original documentation:

Jacobian of ceval_mass with respect to positional argument(s) 0. Takes the same arguments as ceval_mass but returns the jacobian of the output with respect to the arguments at positions 0.

sfepy.terms.terms_jax.eval_elasticity_l(lam, mu, vbfg, ubfg, det, cu)

Vectorized version of ceval_elasticity_l. Takes similar arguments as ceval_elasticity_l but with additional array axes over which ceval_elasticity_l is mapped.

sfepy.terms.terms_jax.eval_elasticity_yp(young, poisson, plane, vbfg, ubfg, det, cu)

Vectorized version of ceval_elasticity_yp. Takes similar arguments as ceval_elasticity_yp but with additional array axes over which ceval_elasticity_yp is mapped.

sfepy.terms.terms_jax.eval_jac_elasticity_l(lam, mu, vbfg, ubfg, det, cu)

Vectorized version of ceval_elasticity_l. Takes similar arguments as ceval_elasticity_l but with additional array axes over which ceval_elasticity_l is mapped.

Original documentation:

Jacobian of ceval_elasticity_l with respect to positional argument(s) -1. Takes the same arguments as ceval_elasticity_l but returns the jacobian of the output with respect to the arguments at positions -1.

sfepy.terms.terms_jax.eval_jac_elasticity_yp(young, poisson, plane, vbfg, ubfg, det, cu)

Vectorized version of ceval_elasticity_yp. Takes similar arguments as ceval_elasticity_yp but with additional array axes over which ceval_elasticity_yp is mapped.

Original documentation:

Jacobian of ceval_elasticity_yp with respect to positional argument(s) -1. Takes the same arguments as ceval_elasticity_yp but returns the jacobian of the output with respect to the arguments at positions -1.

sfepy.terms.terms_jax.eval_jac_mass(density, vbf, ubf, det, cu)

Vectorized version of ceval_mass. Takes similar arguments as ceval_mass but with additional array axes over which ceval_mass is mapped.

Original documentation:

Jacobian of ceval_mass with respect to positional argument(s) -1. Takes the same arguments as ceval_mass but returns the jacobian of the output with respect to the arguments at positions -1.

sfepy.terms.terms_jax.eval_jac_neohook(mu, vbfg, ubfg, det, cu)

Vectorized version of ceval_neohook. Takes similar arguments as ceval_neohook but with additional array axes over which ceval_neohook is mapped.

Original documentation:

Jacobian of ceval_neohook with respect to positional argument(s) -1. Takes the same arguments as ceval_neohook but returns the jacobian of the output with respect to the arguments at positions -1.

sfepy.terms.terms_jax.eval_jac_ogden(mu, alpha, vbfg, ubfg, det, cu)

Vectorized version of ceval_ogden. Takes similar arguments as ceval_ogden but with additional array axes over which ceval_ogden is mapped.

Original documentation:

Jacobian of ceval_ogden with respect to positional argument(s) -1. Takes the same arguments as ceval_ogden but returns the jacobian of the output with respect to the arguments at positions -1.

sfepy.terms.terms_jax.eval_lam_elasticity_l(lam, mu, vbfg, ubfg, det, cu)

Vectorized version of ceval_elasticity_l. Takes similar arguments as ceval_elasticity_l but with additional array axes over which ceval_elasticity_l is mapped.

Original documentation:

Jacobian of ceval_elasticity_l with respect to positional argument(s) 0. Takes the same arguments as ceval_elasticity_l but returns the jacobian of the output with respect to the arguments at positions 0.

sfepy.terms.terms_jax.eval_mass(density, vbf, ubf, det, cu)

Vectorized version of ceval_mass. Takes similar arguments as ceval_mass but with additional array axes over which ceval_mass is mapped.

sfepy.terms.terms_jax.eval_mu_elasticity_l(lam, mu, vbfg, ubfg, det, cu)

Vectorized version of ceval_elasticity_l. Takes similar arguments as ceval_elasticity_l but with additional array axes over which ceval_elasticity_l is mapped.

Original documentation:

Jacobian of ceval_elasticity_l with respect to positional argument(s) 1. Takes the same arguments as ceval_elasticity_l but returns the jacobian of the output with respect to the arguments at positions 1.

sfepy.terms.terms_jax.eval_mu_neohook(mu, vbfg, ubfg, det, cu)

Vectorized version of ceval_neohook. Takes similar arguments as ceval_neohook but with additional array axes over which ceval_neohook is mapped.

Original documentation:

Jacobian of ceval_neohook with respect to positional argument(s) 0. Takes the same arguments as ceval_neohook but returns the jacobian of the output with respect to the arguments at positions 0.

sfepy.terms.terms_jax.eval_mu_ogden(mu, alpha, vbfg, ubfg, det, cu)

Vectorized version of ceval_ogden. Takes similar arguments as ceval_ogden but with additional array axes over which ceval_ogden is mapped.

Original documentation:

Jacobian of ceval_ogden with respect to positional argument(s) 0. Takes the same arguments as ceval_ogden but returns the jacobian of the output with respect to the arguments at positions 0.

sfepy.terms.terms_jax.eval_neohook(mu, vbfg, ubfg, det, cu)

Vectorized version of ceval_neohook. Takes similar arguments as ceval_neohook but with additional array axes over which ceval_neohook is mapped.

sfepy.terms.terms_jax.eval_ogden(mu, alpha, vbfg, ubfg, det, cu)

Vectorized version of ceval_ogden. Takes similar arguments as ceval_ogden but with additional array axes over which ceval_ogden is mapped.

sfepy.terms.terms_jax.eval_poisson_elasticity_yp(young, poisson, plane, vbfg, ubfg, det, cu)

Vectorized version of ceval_elasticity_yp. Takes similar arguments as ceval_elasticity_yp but with additional array axes over which ceval_elasticity_yp is mapped.

Original documentation:

Jacobian of ceval_elasticity_yp with respect to positional argument(s) 1. Takes the same arguments as ceval_elasticity_yp but returns the jacobian of the output with respect to the arguments at positions 1.

sfepy.terms.terms_jax.eval_young_elasticity_yp(young, poisson, plane, vbfg, ubfg, det, cu)

Vectorized version of ceval_elasticity_yp. Takes similar arguments as ceval_elasticity_yp but with additional array axes over which ceval_elasticity_yp is mapped.

Original documentation:

Jacobian of ceval_elasticity_yp with respect to positional argument(s) 0. Takes the same arguments as ceval_elasticity_yp but returns the jacobian of the output with respect to the arguments at positions 0.

sfepy.terms.terms_jax.get_neohook_strain_energy(mu, C)[source]
sfepy.terms.terms_jax.get_neohook_strain_energy_f(mu, F)[source]
sfepy.terms.terms_jax.get_neohook_stress_1pk(mu, gu)[source]
sfepy.terms.terms_jax.get_neohook_stress_2pk(mu, gu)[source]
sfepy.terms.terms_jax.get_ogden_strain_energy_f(mu, alpha, F)[source]
sfepy.terms.terms_jax.get_ogden_stress_1pk(mu, alpha, gu)[source]
sfepy.terms.terms_jax.get_strain(gu)[source]
sfepy.terms.terms_jax.get_stress(lam, mu, gu)[source]