Tools & Helpers#

This page contains tools and helpers for constitutive material formulations.

View Force-Stretch Curves on Elementary Deformations

ViewMaterial(umat[, ux, ps, bx, statevars])

Create views on normal force per undeformed area vs stretch curves for the elementary homogeneous deformations uniaxial tension/compression, planar shear and biaxial tension of a given isotropic material formulation.

ViewMaterialIncompressible(umat[, ux, ps, ...])

Create views on normal force per undeformed area vs stretch curves for the elementary homogeneous incompressible deformations uniaxial tension/compression, planar shear and biaxial tension of a given isotropic material formulation.

Special Constitutive Materials

CompositeMaterial(material, other_material)

A composite material with two constitutive materials merged.

Volumetric(bulk[, parallel])

Neo-Hookean material formulation with deactivated shear modulus.

LineChange([parallel])

Line Change.

AreaChange([parallel])

Area Change.

VolumeChange([parallel])

Volume Change.

Other

constitution.lame_converter(E, nu)

Convert the pair of given material parameters Young's modulus \(E\) and Poisson ratio \(\nu\) to first and second Lamé - constants \(\lambda\) and \(\mu\).

constitution.lame_converter_orthotropic(E, nu, G)

Convert elastic orthotropic material parameters to Lamé constants.

Detailed API Reference

class felupe.ViewMaterial(umat, ux=array([0.7, 0.75, 0.8, 0.85, 0.9, 0.95, 1., 1.05, 1.1, 1.15, 1.2, 1.25, 1.3, 1.35, 1.4, 1.45, 1.5, 1.55, 1.6, 1.65, 1.7, 1.75, 1.8, 1.85, 1.9, 1.95, 2., 2.05, 2.1, 2.15, 2.2, 2.25, 2.3, 2.35, 2.4, 2.45, 2.5]), ps=array([1., 1.05, 1.1, 1.15, 1.2, 1.25, 1.3, 1.35, 1.4, 1.45, 1.5, 1.55, 1.6, 1.65, 1.7, 1.75, 1.8, 1.85, 1.9, 1.95, 2., 2.05, 2.1, 2.15, 2.2, 2.25, 2.3, 2.35, 2.4, 2.45, 2.5]), bx=array([1., 1.05, 1.1, 1.15, 1.2, 1.25, 1.3, 1.35, 1.4, 1.45, 1.5, 1.55, 1.6, 1.65, 1.7, 1.75]), statevars=None)[ソース]#

Create views on normal force per undeformed area vs stretch curves for the elementary homogeneous deformations uniaxial tension/compression, planar shear and biaxial tension of a given isotropic material formulation.

パラメータ:
  • umat (class) -- A class with methods for the gradient and hessian of the strain energy density function w.r.t. the deformation gradient. See Material for further details.

  • ux (ndarray, optional) -- Array with stretches for uniaxial tension/compression. Default is linsteps([0.7, 2.5], num=36).

  • ps (ndarray, optional) -- Array with stretches for planar shear. Default is linsteps([1.0, 2.5], num=30).

  • bx (ndarray, optional) -- Array with stretches for equi-biaxial tension. Default is linsteps([1.0, 1.75], num=15)`.

  • statevars (ndarray or None, optional) -- Array with state variables (default is None). If None, the state variables are assumed to be initially zero.

サンプル

>>> import felupe as fem
>>>
>>> umat = fem.OgdenRoxburgh(fem.NeoHooke(mu=1, bulk=2), r=3, m=1, beta=0)
>>> view = fem.ViewMaterial(
...     umat,
...     ux=fem.math.linsteps([1, 1.5, 1, 2, 1, 2.5, 1], num=15),
...     ps=None,
...     bx=None,
... )
>>> ax = view.plot(show_title=True, show_kwargs=True)
../../_images/tools-1.png
biaxial(stretches=None, **kwargs)[ソース]#

Normal force per undeformed area vs stretch curve for a equi-biaxial deformation.

パラメータ:

stretches (ndarray or None, optional) -- Array with stretches at which the forces are evaluated (default is None). If None, the stretches from initialization are used.

戻り値:

3-tuple with array of stretches and array of forces and the label.

戻り値の型:

tuple of ndarray

evaluate()#

Evaluate normal force per undeformed area vs stretch curves for the elementary homogeneous incompressible deformations uniaxial tension/compression, planar shear and biaxial tension. A load case is not included if its array of stretches (attribute ux, ps or bx) is None.

戻り値:

List with 3-tuple of stretch and force arrays and the label string for each load case.

戻り値の型:

list of 3-tuple

planar(stretches=None, **kwargs)[ソース]#

Normal force per undeformed area vs stretch curve for a planar shear deformation.

パラメータ:

stretches (ndarray or None, optional) -- Array with stretches at which the forces are evaluated (default is None). If None, the stretches from initialization are used.

戻り値:

3-tuple with array of stretches and array of forces and the label.

戻り値の型:

tuple of ndarray

plot(ax=None, show_title=True, show_kwargs=True, **kwargs)#

Plot normal force per undeformed area vs stretch curves for the elementary homogeneous incompressible deformations uniaxial tension/compression, planar shear and biaxial tension.

uniaxial(stretches=None, **kwargs)[ソース]#

Normal force per undeformed area vs. stretch curve for a uniaxial deformation.

パラメータ:

stretches (ndarray or None, optional) -- Array with stretches at which the forces are evaluated (default is None). If None, the stretches from initialization are used.

戻り値:

3-tuple with array of stretches and array of forces and the label.

戻り値の型:

tuple of ndarray

class felupe.ViewMaterialIncompressible(umat, ux=array([0.7, 0.75, 0.8, 0.85, 0.9, 0.95, 1., 1.05, 1.1, 1.15, 1.2, 1.25, 1.3, 1.35, 1.4, 1.45, 1.5, 1.55, 1.6, 1.65, 1.7, 1.75, 1.8, 1.85, 1.9, 1.95, 2., 2.05, 2.1, 2.15, 2.2, 2.25, 2.3, 2.35, 2.4, 2.45, 2.5]), ps=array([1., 1.05, 1.1, 1.15, 1.2, 1.25, 1.3, 1.35, 1.4, 1.45, 1.5, 1.55, 1.6, 1.65, 1.7, 1.75, 1.8, 1.85, 1.9, 1.95, 2., 2.05, 2.1, 2.15, 2.2, 2.25, 2.3, 2.35, 2.4, 2.45, 2.5]), bx=array([1., 1.05, 1.1, 1.15, 1.2, 1.25, 1.3, 1.35, 1.4, 1.45, 1.5, 1.55, 1.6, 1.65, 1.7, 1.75]), statevars=None)[ソース]#

Create views on normal force per undeformed area vs stretch curves for the elementary homogeneous incompressible deformations uniaxial tension/compression, planar shear and biaxial tension of a given isotropic material formulation.

パラメータ:
  • umat (class) -- A class with methods for the gradient and hessian of the strain energy density function w.r.t. the deformation gradient. See Material for further details.

  • ux (ndarray, optional) -- Array with stretches for incompressible uniaxial tension/compression. Default is linsteps([0.7, 2.5], num=36).

  • ps (ndarray, optional) -- Array with stretches for incompressible planar shear. Default is linsteps([1, 2.5], num=30).

  • bx (ndarray, optional) -- Array with stretches for incompressible equi-biaxial tension. Default is linsteps([1, 1.75], num=15).

  • statevars (ndarray or None, optional) -- Array with state variables (default is None). If None, the state variables are assumed to be initially zero.

サンプル

>>> import felupe as fem
>>>
>>> umat = fem.Hyperelastic(fem.extended_tube, Gc=0.2, Ge=0.2, beta=0.2, delta=0.1)
>>> preview = fem.ViewMaterialIncompressible(umat)
>>> ax = preview.plot(show_title=True, show_kwargs=True)
../../_images/tools-2.png
>>> umat = fem.OgdenRoxburgh(fem.NeoHooke(mu=1), r=3, m=1, beta=0)
>>> view = fem.ViewMaterialIncompressible(
...     umat,
...     ux=fem.math.linsteps([1, 1.5, 1, 2, 1, 2.5, 1], num=15),
...     ps=None,
...     bx=None,
... )
>>> ax = view.plot(show_title=True, show_kwargs=True)
../../_images/tools-3.png
biaxial(stretches=None)[ソース]#

Normal force per undeformed area vs stretch curve for a equi-biaxial incompressible deformation.

パラメータ:

stretches (ndarray or None, optional) -- Array with stretches at which the forces are evaluated (default is None). If None, the stretches from initialization are used.

戻り値:

3-tuple with array of stretches and array of forces and the label.

戻り値の型:

tuple of ndarray

evaluate()#

Evaluate normal force per undeformed area vs stretch curves for the elementary homogeneous incompressible deformations uniaxial tension/compression, planar shear and biaxial tension. A load case is not included if its array of stretches (attribute ux, ps or bx) is None.

戻り値:

List with 3-tuple of stretch and force arrays and the label string for each load case.

戻り値の型:

list of 3-tuple

planar(stretches=None)[ソース]#

Normal force per undeformed area vs stretch curve for a planar shear incompressible deformation.

パラメータ:

stretches (ndarray or None, optional) -- Array with stretches at which the forces are evaluated (default is None). If None, the stretches from initialization are used.

戻り値:

3-tuple with array of stretches and array of forces and the label.

戻り値の型:

tuple of ndarray

plot(ax=None, show_title=True, show_kwargs=True, **kwargs)#

Plot normal force per undeformed area vs stretch curves for the elementary homogeneous incompressible deformations uniaxial tension/compression, planar shear and biaxial tension.

uniaxial(stretches=None)[ソース]#

Normal force per undeformed area vs stretch curve for a uniaxial incompressible deformation.

パラメータ:

stretches (ndarray or None, optional) -- Array with stretches at which the forces are evaluated (default is None). If None, the stretches from initialization are used.

戻り値:

3-tuple with array of stretches and array of forces and the label.

戻り値の型:

tuple of ndarray

class felupe.CompositeMaterial(material, other_material)[ソース]#

A composite material with two constitutive materials merged. State variables are only considered for the first material.

パラメータ:

メモ

警告

Do not merge two constitutive materials with the same keys of material parameters. In this case, the values of these material parameters are taken from the first constitutive material.

サンプル

>>> import felupe as fem
>>>
>>> nh = fem.NeoHooke(mu=1.0)
>>> vol = fem.Volumetric(bulk=2.0)
>>> umat = nh & vol
>>> ax = umat.plot()
../../_images/tools-4.png
copy()#

Return a deep-copy of the constitutive material.

gradient(x, **kwargs)[ソース]#
hessian(x, **kwargs)[ソース]#
is_stable(x, hessian=None)#

Return a boolean mask for stability of isotropic material model formulations.

At a given deformation gradient, a normal force is applied on each principal stretch direction. If the resulting incremental stretches are positive, the material model formulation is considered to be stable at the given deformation gradient.

パラメータ:
  • x (list of ndarray) -- The list with input arguments. These contain the extracted fields of a FieldContainer.

  • hessian (ndarray or None, optional) -- Second partial derivative of the strain energy density function w.r.t. the deformation gradient. Default is None.

戻り値:

Boolean mask of stability.

戻り値の型:

ndarray

メモ

警告

This stability check will lead to a singular matrix for isotropic (hyperelastic) material model formulations without a volumetric part.

サンプル

First, let's check the stability of the Neo-Hookean material model formulation. The stability is evaluated on (valid) principal stretches of a biaxial deformation. All deformations are stable.

>>> import numpy as np
>>> import felupe as fem
>>>
>>> umat = fem.NeoHooke(mu=1.0, bulk=2.0)
>>> view = umat.view()
>>> λ = view.biaxial()[0]
>>>
>>> F = np.zeros((3, 3, 1, λ[0].size))
>>> for a in range(3):
...     F[a, a] = λ[a]
>>>
>>> umat.is_stable([F])
array([[ True,  True,  True,  True,  True,  True,  True,  True,  True,
         True,  True,  True,  True,  True,  True,  True]])

Now, let's check the stability of the Mooney-Rivlin material model formulation. The stability is evaluated on (valid) principal stretches of a biaxial deformation. Biaxial deformations are only stable up to a longitudinal stretch of 1.35.

>>> import numpy as np
>>> import felupe as fem
>>> import felupe.constitution.tensortrax as mat
>>>
>>> umat = fem.Hyperelastic(
...     mat.models.hyperelastic.mooney_rivlin,
...     C10=0.25,
...     C01=0.25,
... ) & fem.Volumetric(bulk=5000)
>>> view = umat.view()
>>> λ = view.biaxial()[0]
>>>
>>> F = np.zeros((3, 3, 1, λ[0].size))
>>> for a in range(3):
...     F[a, a] = λ[a]
>>>
>>> umat.is_stable([F])
array([[ True,  True,  True,  True,  True,  True,  True,  True, False,
        False, False, False, False, False, False, False]])
optimize(ux=None, ps=None, bx=None, incompressible=False, relative=False, **kwargs)#

Optimize the material parameters by a least-squares fit on experimental stretch-stress data.

パラメータ:
  • ux (array of shape (2, ...) or None, optional) -- Experimental uniaxial stretch and force-per-undeformed-area data (default is None).

  • ps (array of shape (2, ...) or None, optional) -- Experimental planar-shear stretch and force-per-undeformed-area data (default is None).

  • bx (array of shape (2, ...) or None, optional) -- Experimental biaxial stretch and force-per-undeformed-area data (default is None).

  • incompressible (bool, optional) -- A flag to enforce incompressible deformations (default is False).

  • relative (bool, optional) -- A flag to optimize relative instead of absolute residuals, i.e. (predicted - observed) / observed instead of predicted - observed (default is False).

  • **kwargs (dict, optional) -- Optional keyword arguments are passed to scipy.optimize.least_squares().

戻り値:

  • ConstitutiveMaterial -- A copy of the constitutive material with the optimized material parameters.

  • scipy.optimize.OptimizeResult -- Represents the optimization result.

メモ

警告

At least one load case, i.e. one of the arguments ux, ps or bx must not be None.

注釈

For JAX-based materials, double-precision is required to optimize material parameters.

import jax

jax.config.update("jax_enable_x64", True)

The vector of residuals is given in Eq. (3) in case of absolute residuals

(1)#\[ \begin{align}\begin{aligned}\begin{split}\boldsymbol{r} &= \begin{bmatrix} \boldsymbol{r}_\text{ux} \\ \boldsymbol{r}_\text{ps} \\ \boldsymbol{r}_\text{bx} \end{bmatrix}\end{split}\\r_\text{ux}(\lambda_i) &= P_\text{ux}(\lambda_i) - P_\text{ux, observed}(\lambda_i)\\r_\text{ps}(\lambda_i) &= P_\text{ps}(\lambda_i) - P_\text{ps, observed}(\lambda_i)\\r_\text{bx}(\lambda_i) &= P_\text{bx}(\lambda_i) - P_\text{bx, observed}(\lambda_i)\end{aligned}\end{align} \]

and in Eq. (4) in case of relative residuals.

(2)#\[ \begin{align}\begin{aligned}\begin{split}\boldsymbol{r} &= \begin{bmatrix} \boldsymbol{r}_\text{ux} \\ \boldsymbol{r}_\text{ps} \\ \boldsymbol{r}_\text{bx} \end{bmatrix}\end{split}\\r_\text{ux}(\lambda_i) &= \frac{ P_\text{ux}(\lambda_i) - P_\text{ux, observed}(\lambda_i)}{ P_\text{ux, observed}(\lambda_i) }\\r_\text{ps}(\lambda_i) &= \frac{ P_\text{ps}(\lambda_i) - P_\text{ps, observed}(\lambda_i)}{ P_\text{ps, observed}(\lambda_i) }\\r_\text{bx}(\lambda_i) &= \frac{ P_\text{bx}(\lambda_i) - P_\text{bx, observed}(\lambda_i)}{ P_\text{bx, observed}(\lambda_i) }\end{aligned}\end{align} \]

サンプル

The Anssari-Benam Bucchi material model formulation is best-fitted on Treloar's uniaxial and biaxial tension data [1]_.

>>> import numpy as np
>>> import felupe as fem
>>>
>>> λ, P = np.array(
...     [
...         [1.000, 0.00],
...         [1.240, 2.30],
...         [1.585, 4.16],
...         [2.180, 6.00],
...         [3.020, 8.80],
...         [4.030, 12.5],
...         [4.760, 16.2],
...         [5.750, 23.6],
...         [6.850, 38.5],
...         [7.250, 49.6],
...         [7.600, 64.4],
...     ]
... ).T * np.array([[1.0], [0.0980665]])
>>>
>>> umat = fem.Hyperelastic(fem.anssari_benam_bucchi)
>>> umat_new, res = umat.optimize(
...     ux=[λ, P], incompressible=True, relative=True
... )
>>>
>>> ux = np.linspace(λ.min(), λ.max(), num=50)
>>> ax = umat_new.plot(incompressible=True, ux=ux, bx=None, ps=None)
>>> ax.plot(λ, P, "C0x")
../../_images/tools-7.png

参考

scipy.optimize.least_squares

Solve a nonlinear least-squares problem with bounds on the variables.

参照

plot(incompressible=False, **kwargs)#

Return a plot with normal force per undeformed area vs. stretch curves for the elementary homogeneous deformations uniaxial tension/compression, planar shear and biaxial tension of a given isotropic material formulation.

パラメータ:
戻り値の型:

matplotlib.axes.Axes

参考

felupe.ViewMaterial

Create views on normal force per undeformed area vs. stretch curves for the elementary homogeneous deformations uniaxial tension/compression, planar shear and biaxial tension of a given isotropic material formulation.

felupe.ViewMaterialIncompressible

Create views on normal force per undeformed area vs. stretch curves for the elementary homogeneous incompressible deformations uniaxial tension/compression, planar shear and biaxial tension of a given isotropic material formulation.

screenshot(filename='umat.png', incompressible=False, **kwargs)#

Save a screenshot with normal force per undeformed area vs. stretch curves for the elementary homogeneous deformations uniaxial tension/compression, planar shear and biaxial tension of a given isotropic material formulation.

パラメータ:
  • filename (str, optional) -- The filename of the screenshot (default is "umat.png").

  • incompressible (bool, optional) -- A flag to enforce views on incompressible deformations (default is False).

  • **kwargs (dict, optional) -- Optional keyword-arguments for ViewMaterial or ViewMaterialIncompressible.

戻り値の型:

matplotlib.axes.Axes

参考

felupe.ViewMaterial

Create views on normal force per undeformed area vs. stretch curves for the elementary homogeneous deformations uniaxial tension/compression, planar shear and biaxial tension of a given isotropic material formulation.

felupe.ViewMaterialIncompressible

Create views on normal force per undeformed area vs. stretch curves for the elementary homogeneous incompressible deformations uniaxial tension/compression, planar shear and biaxial tension of a given isotropic material formulation.

view(incompressible=False, **kwargs)#

Create views on normal force per undeformed area vs. stretch curves for the elementary homogeneous deformations uniaxial tension/compression, planar shear and biaxial tension of a given isotropic material formulation.

パラメータ:
戻り値の型:

felupe.ViewMaterial or felupe.ViewMaterialIncompressible

参考

felupe.ViewMaterial

Create views on normal force per undeformed area vs. stretch curves for the elementary homogeneous deformations uniaxial tension/compression, planar shear and biaxial tension of a given isotropic material formulation.

felupe.ViewMaterialIncompressible

Create views on normal force per undeformed area vs. stretch curves for the elementary homogeneous incompressible deformations uniaxial tension/compression, planar shear and biaxial tension of a given isotropic material formulation.

class felupe.Volumetric(bulk, parallel=False)[ソース]#

Neo-Hookean material formulation with deactivated shear modulus.

copy()#

Return a deep-copy of the constitutive material.

function(x)#

Strain energy density function per unit undeformed volume of the Neo-Hookean material formulation.

パラメータ:

x (list of ndarray) -- List with the Deformation gradient F (3x3) as first item

gradient(x, out=None)#

Gradient of the strain energy density function per unit undeformed volume of the Neo-Hookean material formulation.

パラメータ:
  • x (list of ndarray) -- List with the Deformation gradient F (3x3) as first item

  • out (ndarray or None, optional) -- A location into which the result is stored (default is None).

hessian(x, out=None)#

Hessian of the strain energy density function per unit undeformed volume of the Neo-Hookean material formulation.

パラメータ:
  • x (list of ndarray) -- List with the Deformation gradient F (3x3) as first item

  • out (ndarray or None, optional) -- A location into which the result is stored (default is None).

is_stable(x, hessian=None)#

Return a boolean mask for stability of isotropic material model formulations.

At a given deformation gradient, a normal force is applied on each principal stretch direction. If the resulting incremental stretches are positive, the material model formulation is considered to be stable at the given deformation gradient.

パラメータ:
  • x (list of ndarray) -- The list with input arguments. These contain the extracted fields of a FieldContainer.

  • hessian (ndarray or None, optional) -- Second partial derivative of the strain energy density function w.r.t. the deformation gradient. Default is None.

戻り値:

Boolean mask of stability.

戻り値の型:

ndarray

メモ

警告

This stability check will lead to a singular matrix for isotropic (hyperelastic) material model formulations without a volumetric part.

サンプル

First, let's check the stability of the Neo-Hookean material model formulation. The stability is evaluated on (valid) principal stretches of a biaxial deformation. All deformations are stable.

>>> import numpy as np
>>> import felupe as fem
>>>
>>> umat = fem.NeoHooke(mu=1.0, bulk=2.0)
>>> view = umat.view()
>>> λ = view.biaxial()[0]
>>>
>>> F = np.zeros((3, 3, 1, λ[0].size))
>>> for a in range(3):
...     F[a, a] = λ[a]
>>>
>>> umat.is_stable([F])
array([[ True,  True,  True,  True,  True,  True,  True,  True,  True,
         True,  True,  True,  True,  True,  True,  True]])

Now, let's check the stability of the Mooney-Rivlin material model formulation. The stability is evaluated on (valid) principal stretches of a biaxial deformation. Biaxial deformations are only stable up to a longitudinal stretch of 1.35.

>>> import numpy as np
>>> import felupe as fem
>>> import felupe.constitution.tensortrax as mat
>>>
>>> umat = fem.Hyperelastic(
...     mat.models.hyperelastic.mooney_rivlin,
...     C10=0.25,
...     C01=0.25,
... ) & fem.Volumetric(bulk=5000)
>>> view = umat.view()
>>> λ = view.biaxial()[0]
>>>
>>> F = np.zeros((3, 3, 1, λ[0].size))
>>> for a in range(3):
...     F[a, a] = λ[a]
>>>
>>> umat.is_stable([F])
array([[ True,  True,  True,  True,  True,  True,  True,  True, False,
        False, False, False, False, False, False, False]])
optimize(ux=None, ps=None, bx=None, incompressible=False, relative=False, **kwargs)#

Optimize the material parameters by a least-squares fit on experimental stretch-stress data.

パラメータ:
  • ux (array of shape (2, ...) or None, optional) -- Experimental uniaxial stretch and force-per-undeformed-area data (default is None).

  • ps (array of shape (2, ...) or None, optional) -- Experimental planar-shear stretch and force-per-undeformed-area data (default is None).

  • bx (array of shape (2, ...) or None, optional) -- Experimental biaxial stretch and force-per-undeformed-area data (default is None).

  • incompressible (bool, optional) -- A flag to enforce incompressible deformations (default is False).

  • relative (bool, optional) -- A flag to optimize relative instead of absolute residuals, i.e. (predicted - observed) / observed instead of predicted - observed (default is False).

  • **kwargs (dict, optional) -- Optional keyword arguments are passed to scipy.optimize.least_squares().

戻り値:

  • ConstitutiveMaterial -- A copy of the constitutive material with the optimized material parameters.

  • scipy.optimize.OptimizeResult -- Represents the optimization result.

メモ

警告

At least one load case, i.e. one of the arguments ux, ps or bx must not be None.

注釈

For JAX-based materials, double-precision is required to optimize material parameters.

import jax

jax.config.update("jax_enable_x64", True)

The vector of residuals is given in Eq. (3) in case of absolute residuals

(3)#\[ \begin{align}\begin{aligned}\begin{split}\boldsymbol{r} &= \begin{bmatrix} \boldsymbol{r}_\text{ux} \\ \boldsymbol{r}_\text{ps} \\ \boldsymbol{r}_\text{bx} \end{bmatrix}\end{split}\\r_\text{ux}(\lambda_i) &= P_\text{ux}(\lambda_i) - P_\text{ux, observed}(\lambda_i)\\r_\text{ps}(\lambda_i) &= P_\text{ps}(\lambda_i) - P_\text{ps, observed}(\lambda_i)\\r_\text{bx}(\lambda_i) &= P_\text{bx}(\lambda_i) - P_\text{bx, observed}(\lambda_i)\end{aligned}\end{align} \]

and in Eq. (4) in case of relative residuals.

(4)#\[ \begin{align}\begin{aligned}\begin{split}\boldsymbol{r} &= \begin{bmatrix} \boldsymbol{r}_\text{ux} \\ \boldsymbol{r}_\text{ps} \\ \boldsymbol{r}_\text{bx} \end{bmatrix}\end{split}\\r_\text{ux}(\lambda_i) &= \frac{ P_\text{ux}(\lambda_i) - P_\text{ux, observed}(\lambda_i)}{ P_\text{ux, observed}(\lambda_i) }\\r_\text{ps}(\lambda_i) &= \frac{ P_\text{ps}(\lambda_i) - P_\text{ps, observed}(\lambda_i)}{ P_\text{ps, observed}(\lambda_i) }\\r_\text{bx}(\lambda_i) &= \frac{ P_\text{bx}(\lambda_i) - P_\text{bx, observed}(\lambda_i)}{ P_\text{bx, observed}(\lambda_i) }\end{aligned}\end{align} \]

サンプル

The Anssari-Benam Bucchi material model formulation is best-fitted on Treloar's uniaxial and biaxial tension data [1]_.

>>> import numpy as np
>>> import felupe as fem
>>>
>>> λ, P = np.array(
...     [
...         [1.000, 0.00],
...         [1.240, 2.30],
...         [1.585, 4.16],
...         [2.180, 6.00],
...         [3.020, 8.80],
...         [4.030, 12.5],
...         [4.760, 16.2],
...         [5.750, 23.6],
...         [6.850, 38.5],
...         [7.250, 49.6],
...         [7.600, 64.4],
...     ]
... ).T * np.array([[1.0], [0.0980665]])
>>>
>>> umat = fem.Hyperelastic(fem.anssari_benam_bucchi)
>>> umat_new, res = umat.optimize(
...     ux=[λ, P], incompressible=True, relative=True
... )
>>>
>>> ux = np.linspace(λ.min(), λ.max(), num=50)
>>> ax = umat_new.plot(incompressible=True, ux=ux, bx=None, ps=None)
>>> ax.plot(λ, P, "C0x")
../../_images/tools-10.png

参考

scipy.optimize.least_squares

Solve a nonlinear least-squares problem with bounds on the variables.

参照

plot(incompressible=False, **kwargs)#

Return a plot with normal force per undeformed area vs. stretch curves for the elementary homogeneous deformations uniaxial tension/compression, planar shear and biaxial tension of a given isotropic material formulation.

パラメータ:
戻り値の型:

matplotlib.axes.Axes

参考

felupe.ViewMaterial

Create views on normal force per undeformed area vs. stretch curves for the elementary homogeneous deformations uniaxial tension/compression, planar shear and biaxial tension of a given isotropic material formulation.

felupe.ViewMaterialIncompressible

Create views on normal force per undeformed area vs. stretch curves for the elementary homogeneous incompressible deformations uniaxial tension/compression, planar shear and biaxial tension of a given isotropic material formulation.

screenshot(filename='umat.png', incompressible=False, **kwargs)#

Save a screenshot with normal force per undeformed area vs. stretch curves for the elementary homogeneous deformations uniaxial tension/compression, planar shear and biaxial tension of a given isotropic material formulation.

パラメータ:
  • filename (str, optional) -- The filename of the screenshot (default is "umat.png").

  • incompressible (bool, optional) -- A flag to enforce views on incompressible deformations (default is False).

  • **kwargs (dict, optional) -- Optional keyword-arguments for ViewMaterial or ViewMaterialIncompressible.

戻り値の型:

matplotlib.axes.Axes

参考

felupe.ViewMaterial

Create views on normal force per undeformed area vs. stretch curves for the elementary homogeneous deformations uniaxial tension/compression, planar shear and biaxial tension of a given isotropic material formulation.

felupe.ViewMaterialIncompressible

Create views on normal force per undeformed area vs. stretch curves for the elementary homogeneous incompressible deformations uniaxial tension/compression, planar shear and biaxial tension of a given isotropic material formulation.

view(incompressible=False, **kwargs)#

Create views on normal force per undeformed area vs. stretch curves for the elementary homogeneous deformations uniaxial tension/compression, planar shear and biaxial tension of a given isotropic material formulation.

パラメータ:
戻り値の型:

felupe.ViewMaterial or felupe.ViewMaterialIncompressible

参考

felupe.ViewMaterial

Create views on normal force per undeformed area vs. stretch curves for the elementary homogeneous deformations uniaxial tension/compression, planar shear and biaxial tension of a given isotropic material formulation.

felupe.ViewMaterialIncompressible

Create views on normal force per undeformed area vs. stretch curves for the elementary homogeneous incompressible deformations uniaxial tension/compression, planar shear and biaxial tension of a given isotropic material formulation.

class felupe.LineChange(parallel=False)[ソース]#

Line Change.

\[d\boldsymbol{x} = \boldsymbol{F} d\boldsymbol{X}\]

Gradient:

\[\frac{\partial \boldsymbol{F}}{\partial \boldsymbol{F}} = \boldsymbol{I} \overset{ik}{\otimes} \boldsymbol{I}\]
function(extract)[ソース]#

Line change.

パラメータ:

extract (list of ndarray) -- List of extracted field values with Deformation gradient as first item.

戻り値:

F -- Deformation gradient

戻り値の型:

ndarray

gradient(extract, parallel=None)[ソース]#

Gradient of line change.

パラメータ:

extract (list of ndarray) -- List of extracted field values with Deformation gradient as first item.

戻り値:

Gradient of line change

戻り値の型:

ndarray

class felupe.AreaChange(parallel=False)[ソース]#

Area Change.

\[d\boldsymbol{a} = J \boldsymbol{F}^{-T} d\boldsymbol{A}\]

Gradient:

\[\frac{\partial J \boldsymbol{F}^{-T}}{\partial \boldsymbol{F}} = J \left( \boldsymbol{F}^{-T} \otimes \boldsymbol{F}^{-T} - \boldsymbol{F}^{-T} \overset{il}{\otimes} \boldsymbol{F}^{-T} \right)\]
function(extract, N=None, parallel=None)[ソース]#

Area change.

パラメータ:
  • extract (list of ndarray) -- List of extracted field values with Deformation gradient as first item.

  • N (ndarray or None, optional) -- Area normal vector (default is None)

戻り値:

Cofactor matrix of the deformation gradient

戻り値の型:

ndarray

gradient(extract, N=None, parallel=None)[ソース]#

Gradient of area change.

パラメータ:
  • extract (list of ndarray) -- List of extracted field values with Deformation gradient as first item.

  • N (ndarray or None, optional) -- Area normal vector (default is None)

戻り値:

Gradient of cofactor matrix of the deformation gradient

戻り値の型:

ndarray

class felupe.VolumeChange(parallel=False)[ソース]#

Volume Change.

\[d\boldsymbol{v} = \text{det}(\boldsymbol{F}) d\boldsymbol{V}\]

Gradient and hessian (equivalent to gradient of area change):

\[ \begin{align}\begin{aligned}\frac{\partial J}{\partial \boldsymbol{F}} &= J \boldsymbol{F}^{-T}\\\frac{\partial^2 J}{\partial \boldsymbol{F}\ \partial \boldsymbol{F}} &= J \left( \boldsymbol{F}^{-T} \otimes \boldsymbol{F}^{-T} - \boldsymbol{F}^{-T} \overset{il}{\otimes} \boldsymbol{F}^{-T} \right)\end{aligned}\end{align} \]
function(extract)[ソース]#

Gradient of volume change.

パラメータ:

extract (list of ndarray) -- List of extracted field values with Deformation gradient as first item.

戻り値:

J -- Determinant of the deformation gradient

戻り値の型:

ndarray

gradient(extract)[ソース]#

Gradient of volume change.

パラメータ:

F (ndarray) -- Deformation gradient

戻り値:

Gradient of the determinant of the deformation gradient

戻り値の型:

ndarray

hessian(extract, parallel=None)[ソース]#

Hessian of volume change.

パラメータ:

extract (list of ndarray) -- List of extracted field values with Deformation gradient as first item.

戻り値:

Hessian of the determinant of the deformation gradient

戻り値の型:

ndarray

felupe.constitution.lame_converter(E, nu)[ソース]#

Convert the pair of given material parameters Young's modulus \(E\) and Poisson ratio \(\nu\) to first and second Lamé - constants \(\lambda\) and \(\mu\).

パラメータ:
  • E (float) -- Young's modulus.

  • nu (float) -- Poisson ratio.

戻り値:

  • lmbda (float) -- First Lamé - constant.

  • mu (float) -- Second Lamé - constant (shear modulus).

メモ

\[ \begin{align}\begin{aligned}\lambda &= \frac{E \nu}{(1 + \nu) (1 - 2 \nu)}\\\mu &= \frac{E}{2 (1 + \nu)}\end{aligned}\end{align} \]
felupe.constitution.lame_converter_orthotropic(E, nu, G)[ソース]#

Convert elastic orthotropic material parameters to Lamé constants.

パラメータ:
  • E (list of float) -- List of the three elastic moduli \(E_1, E_2, E_3\).

  • nu (list of float) -- List of three poisson ratios \(\nu_{12}, \nu_{23}, \nu_{31}\).

  • G (list of float) -- List of three shear moduli \(G_{12}, G_{23}, G_{31}\).

戻り値:

  • lmbda (list of float) -- List of six (upper triangle) first Lamé parameters \(\lambda_{11}, \lambda_{12}, \lambda_{13}, \lambda_{22}, \lambda_{23}, \lambda_{33}\).

  • mu (list of float) -- List of the three second Lamé parameters \(\mu_1,\mu_2, \mu_3\).

メモ

The orthotropic material parameters are converted to orthotropic Lamé constants.

The compliance matrix as the inverse of the stiffness matrix with the parameters \(E_i\), \(\nu_{ij}\) and \(G_{ij}\) is given in Eq. (6).

(5)#\[\begin{split}\boldsymbol{C}^{-1} = \begin{bmatrix} \frac{1}{E_1} & -\frac{\nu_{21}}{E_2} & -\frac{\nu_{31}}{E_3} & 0 & 0 & 0 \\ -\frac{\nu_{12}}{E_1} & \frac{1}{E_2} & -\frac{\nu_{32}}{E_3} & 0 & 0 & 0 \\ -\frac{\nu_{13}}{E_1} & -\frac{\nu_{23}}{E_2} & \frac{1}{E_3} & 0 & 0 & 0 \\ 0 & 0 & 0 & \frac{1}{G_{12}} & 0 & 0 \\ 0 & 0 & 0 & 0 & \frac{1}{G_{23}} & 0 \\ 0 & 0 & 0 & 0 & 0 & \frac{1}{G_{31}} \end{bmatrix}\end{split}\]

The stiffness matrix with the Lamé constants is denoted in Eq. ortho-matrix.

(6)#\[\begin{split}\boldsymbol{C} = \begin{bmatrix} \lambda_{11} + 2 \mu_1 & \lambda_{12} & \lambda_{13} & 0 & 0 & 0 \\ \lambda_{11} & \lambda_{12} + 2 \mu_2 & \lambda_{13} & 0 & 0 & 0 \\ \lambda_{11} & \lambda_{12} & \lambda_{13} + 2 \mu_3 & 0 & 0 & 0 \\ 0 & 0 & 0 & \frac{\mu_1 + \mu_2}{2} & 0 & 0 \\ 0 & 0 & 0 & 0 & \frac{\mu_2 + \mu_3}{2} & 0 \\ 0 & 0 & 0 & 0 & 0 & \frac{\mu_3 + \mu_1}{2} \end{bmatrix}\end{split}\]

Eq. (6) is evaluated and inverted numerically to extract the Lamé constants.

参考

felupe.LinearElasticOrthotropic

Orthotropic linear-elastic material formulation.

felupe.constitution.tensortrax.models.hyperelastic.saint_venant_kirchhoff_orthotropic

Strain energy function of the orthotropic hyperelastic Saint-Venant Kirchhoff material formulation.