sfepy.terms.terms_elastic module

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

Evaluate Cauchy strain tensor.

It is given in the usual vector form exploiting symmetry: in 3D it has 6 components with the indices ordered as [11, 22, 33, 12, 13, 23], in 2D it has 3 components with the indices ordered as [11, 22,
12]. The last three (non-diagonal) components are doubled so that it is energetically conjugate to the Cauchy stress tensor with the same storage.

Supports ‘eval’, ‘el_avg’ and ‘qp’ evaluation modes.

Definition:

\int_{\cal{D}} \ull{e}(\ul{w})

Call signature:

ev_cauchy_strain

(parameter)

Arguments:
  • parameter : \ul{w}

arg_shapes = {'parameter': 'D'}
arg_types = ('parameter',)
static function(out, strain, vg, fmode)[source]
get_eval_shape(parameter, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
get_fargs(parameter, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
integration = ('cell', 'facet_extra')
name = 'ev_cauchy_strain'
class sfepy.terms.terms_elastic.CauchyStressETHTerm(name, arg_str, integral, region, **kwargs)[source]

Evaluate fading memory Cauchy stress tensor.

It is given in the usual vector form exploiting symmetry: in 3D it has 6 components with the indices ordered as [11, 22, 33, 12, 13, 23], in 2D it has 3 components with the indices ordered as [11, 22, 12].

Assumes an exponential approximation of the convolution kernel resulting in much higher efficiency.

Supports ‘eval’, ‘el_avg’ and ‘qp’ evaluation modes.

Definition:

\int_{\Omega} \int_0^t \Hcal_{ijkl}(t-\tau)\,e_{kl}(\ul{w}(\tau))
\difd{\tau}

Call signature:

ev_cauchy_stress_eth

(ts, material_0, material_1, parameter)

Arguments:
  • ts : TimeStepper instance

  • material_0 : \Hcal_{ijkl}(0)

  • material_1 : \exp(-\lambda \Delta t) (decay at t_1)

  • parameter : \ul{w}

arg_shapes = {'material_0': 'S, S', 'material_1': '1, 1', 'parameter': 'D'}
arg_types = ('ts', 'material_0', 'material_1', 'parameter')
get_eval_shape(ts, mat0, mat1, parameter, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
get_fargs(ts, mat0, mat1, state, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
name = 'ev_cauchy_stress_eth'
class sfepy.terms.terms_elastic.CauchyStressTHTerm(name, arg_str, integral, region, **kwargs)[source]

Evaluate fading memory Cauchy stress tensor.

It is given in the usual vector form exploiting symmetry: in 3D it has 6 components with the indices ordered as [11, 22, 33, 12, 13, 23], in 2D it has 3 components with the indices ordered as [11, 22, 12].

Supports ‘eval’, ‘el_avg’ and ‘qp’ evaluation modes.

Definition:

\int_{\Omega} \int_0^t \Hcal_{ijkl}(t-\tau)\,e_{kl}(\ul{w}(\tau))
\difd{\tau}

Call signature:

ev_cauchy_stress_th

(ts, material, parameter)

Arguments:
  • ts : TimeStepper instance

  • material : \Hcal_{ijkl}(\tau)

  • parameter : \ul{w}

arg_shapes = {'material': '.: N, S, S', 'parameter': 'D'}
arg_types = ('ts', 'material', 'parameter')
get_eval_shape(ts, mats, parameter, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
get_fargs(ts, mats, state, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
name = 'ev_cauchy_stress_th'
class sfepy.terms.terms_elastic.CauchyStressTerm(name, arg_str, integral, region, **kwargs)[source]

Evaluate Cauchy stress tensor.

It is given in the usual vector form exploiting symmetry: in 3D it has 6 components with the indices ordered as [11, 22, 33, 12, 13, 23], in 2D it has 3 components with the indices ordered as [11, 22, 12].

Supports ‘eval’, ‘el_avg’ and ‘qp’ evaluation modes.

Definition:

\int_{\cal{D}} D_{ijkl} e_{kl}(\ul{w})

Call signature:

ev_cauchy_stress

(material, parameter)

Arguments:
  • material : D_{ijkl}

  • parameter : \ul{w}

arg_shapes = {'material': 'S, S', 'parameter': 'D'}
arg_types = ('material', 'parameter')
static function(out, coef, strain, mat, vg, fmode)[source]
get_eval_shape(mat, parameter, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
get_fargs(mat, parameter, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
integration = ('cell', 'facet_extra')
name = 'ev_cauchy_stress'
class sfepy.terms.terms_elastic.ElasticWaveCauchyTerm(name, arg_str, integral, region, **kwargs)[source]

Elastic dispersion term involving the wave strain g_{ij}, g_{ij}(\ul{u}) = \frac{1}{2}(u_i \kappa_j + \kappa_i u_j), with the wave vector \ul{\kappa} and the elastic strain e_{ij}. D_{ijkl} is given in the usual matrix form exploiting symmetry: in 3D it is 6\times6 with the indices ordered as [11, 22, 33,
12, 13, 23], in 2D it is 3\times3 with the indices ordered as [11, 22, 12].

Definition:

\int_{\Omega} D_{ijkl}\ g_{ij}(\ul{v}) e_{kl}(\ul{u})\\
\int_{\Omega} D_{ijkl}\ g_{ij}(\ul{u}) e_{kl}(\ul{v})

Call signature:

dw_elastic_wave_cauchy

(material_1, material_2, virtual, state)

(material_1, material_2, state, virtual)

Arguments 1:
  • material_1 : D_{ijkl}

  • material_2 : \ul{\kappa}

  • virtual : \ul{v}

  • state : \ul{u}

Arguments 2:
  • material_1 : D_{ijkl}

  • material_2 : \ul{\kappa}

  • state : \ul{u}

  • virtual : \ul{v}

arg_shapes = {'material_1': 'S, S', 'material_2': '.: D', 'state': 'D', 'virtual': ('D', 'state')}
arg_types = (('material_1', 'material_2', 'virtual', 'state'), ('material_1', 'material_2', 'state', 'virtual'))
static function(out, out_qp, geo, fmode)[source]
geometries = ['2_3', '2_4', '3_4', '3_8']
get_fargs(mat, kappa, gvar, evar, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
modes = ('ge', 'eg')
name = 'dw_elastic_wave_cauchy'
class sfepy.terms.terms_elastic.ElasticWaveTerm(name, arg_str, integral, region, **kwargs)[source]

Elastic dispersion term involving the wave strain g_{ij}, g_{ij}(\ul{u}) = \frac{1}{2}(u_i \kappa_j + \kappa_i u_j), with the wave vector \ul{\kappa}. D_{ijkl} is given in the usual matrix form exploiting symmetry: in 3D it is 6\times6 with the indices ordered as [11, 22, 33, 12, 13, 23], in 2D it is 3\times3 with the indices ordered as [11, 22, 12].

Definition:

\int_{\Omega} D_{ijkl}\ g_{ij}(\ul{v}) g_{kl}(\ul{u})

Call signature:

dw_elastic_wave

(material_1, material_2, virtual, state)

Arguments:
  • material_1 : D_{ijkl}

  • material_2 : \ul{\kappa}

  • virtual : \ul{v}

  • state : \ul{u}

arg_shapes = {'material_1': 'S, S', 'material_2': '.: D', 'state': 'D', 'virtual': ('D', 'state')}
arg_types = ('material_1', 'material_2', 'virtual', 'state')
static function(out, out_qp, geo, fmode)[source]
geometries = ['2_3', '2_4', '3_4', '3_8']
get_fargs(mat, kappa, virtual, state, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
name = 'dw_elastic_wave'
class sfepy.terms.terms_elastic.LinearDRotSpringTerm(name, arg_str, integral, region, **kwargs)[source]
Call signature:

dw_lin_dspring_rot

(opt_material, material, virtual, state)

arg_shapes = [{'material': 'S, 1', 'opt_material': 'D, 1', 'state': 'S', 'virtual': ('S', 'state')}, {'material': 'S, S'}, {'opt_material': None}]
arg_types = ('opt_material', 'material', 'virtual', 'state')
name = 'dw_lin_dspring_rot'
class sfepy.terms.terms_elastic.LinearDSpringTerm(name, arg_str, integral, region, **kwargs)[source]

Linear spring element with the stiffness transformed into the element direction.

Definition:

f^{(i)}_k = -f^{(j)}_k = K_{kl} (u^{(j)}_l - u^{(i)}_l)\\
\quad \forall \mbox{ elements } T_K^{i,j}\\
\mbox{ in a region connecting nodes } i, j

Call signature:

dw_lin_dspring

(opt_material, material, virtual, state)

Arguments:
  • opt_material : \ul{d}

  • material : \ul{k}

  • virtual: \ul{v}

  • state: \ul{u}

Stiffness matrix \ul{K} = \ul{T(\ul{d})}^T \ul{K(\ul{k})} \ul{T(\ul{d})} is defined by 6 components \ul{k} = [k_{u1}, k_{u2}, k_{u3}, k_{r1}, k_{r2}, k_{r3}] in 3D and by 3 components \ul{k} = [k_{u1}, k_{u2}, k_{r1}], where k_{ui} is the stiffness for the displacement DOF and r_{ui} is for the rotational DOF. Note that the components of \ul{k} are in the local coordinates system specified by a given direction \ul{d} or by the vector \ul{d} = \ul{x}^{(j)} - \ul{x}^{(i)} for non-coincidental end nodes. The stiffness parameter \ul{K} can also be defined as a 6x6 matrix in 3D or a 3x3 matrix in 2D.

arg_shapes = [{'material': 'D, 1', 'opt_material': 'D, 1', 'state': 'D', 'virtual': ('D', 'state')}, {'material': 'D, D'}, {'opt_material': None}]
arg_types = ('opt_material', 'material', 'virtual', 'state')
static function(out, mat, vec, mtx_t, diff_var)[source]
geometries = ['1_2', '2_1_2', '3_1_2']
get_fargs(dvec, mat, virtual, state, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
integration_order = 0
name = 'dw_lin_dspring'
class sfepy.terms.terms_elastic.LinearElasticETHTerm(name, arg_str, integral, region, **kwargs)[source]

This term has the same definition as dw_lin_elastic_th, but assumes an exponential approximation of the convolution kernel resulting in much higher efficiency. Can use derivatives.

Definition:

\int_{\Omega} \left [\int_0^t
\Hcal_{ijkl}(t-\tau)\,e_{kl}(\ul{u}(\tau)) \difd{\tau}
\right]\,e_{ij}(\ul{v})

Call signature:

dw_lin_elastic_eth

(ts, material_0, material_1, virtual, state)

Arguments:
  • ts : TimeStepper instance

  • material_0 : \Hcal_{ijkl}(0)

  • material_1 : \exp(-\lambda \Delta t) (decay at t_1)

  • virtual : \ul{v}

  • state : \ul{u}

arg_shapes = {'material_0': 'S, S', 'material_1': '1, 1', 'state': 'D', 'virtual': ('D', 'state')}
arg_types = ('ts', 'material_0', 'material_1', 'virtual', 'state')
static function(out, coef, strain, mtx_d, cmap, is_diff)
get_fargs(ts, mat0, mat1, virtual, state, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
name = 'dw_lin_elastic_eth'
class sfepy.terms.terms_elastic.LinearElasticIsotropicTerm(name, arg_str, integral, region, **kwargs)[source]

Isotropic linear elasticity term.

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_iso

(material_1, material_2, virtual, state)

(material_1, material_2, parameter_1, parameter_2)

Arguments:
  • material_1: \lambda

  • material_2: \mu

  • virtual/parameter_1: \ul{v}

  • state/parameter_2: \ul{u}

arg_shapes = {'material_1': '1, 1', 'material_2': '1, 1', 'parameter_1': 'D', 'parameter_2': 'D', 'state': 'D', 'virtual': ('D', 'state')}
arg_types = (('material_1', 'material_2', 'virtual', 'state'), ('material_1', 'material_2', 'parameter_1', 'parameter_2'))
geometries = ['2_3', '2_4', '3_4', '3_8']
get_eval_shape(mat1, mat2, virtual, state, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
get_fargs(lam, mu, virtual, state, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
name = 'dw_lin_elastic_iso'
class sfepy.terms.terms_elastic.LinearElasticTHTerm(name, arg_str, integral, region, **kwargs)[source]

Fading memory linear elastic (viscous) term. Can use derivatives.

Definition:

\int_{\Omega} \left [\int_0^t
\Hcal_{ijkl}(t-\tau)\,e_{kl}(\ul{u}(\tau)) \difd{\tau}
\right]\,e_{ij}(\ul{v})

Call signature:

dw_lin_elastic_th

(ts, material, virtual, state)

Arguments:
  • ts : TimeStepper instance

  • material : \Hcal_{ijkl}(\tau)

  • virtual : \ul{v}

  • state : \ul{u}

arg_shapes = {'material': '.: N, S, S', 'state': 'D', 'virtual': ('D', 'state')}
arg_types = ('ts', 'material', 'virtual', 'state')
static function(out, coef, strain, mtx_d, cmap, is_diff)
get_fargs(ts, mats, virtual, state, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
name = 'dw_lin_elastic_th'
class sfepy.terms.terms_elastic.LinearElasticTerm(name, arg_str, integral, region, **kwargs)[source]

General linear elasticity term, with D_{ijkl} given in the usual matrix form exploiting symmetry: in 3D it is 6\times6 with the indices ordered as [11, 22, 33, 12, 13, 23], in 2D it is 3\times3 with the indices ordered as [11, 22, 12]. Can be evaluated. Can use derivatives.

Definition:

\int_{\Omega} D_{ijkl}\ e_{ij}(\ul{v}) e_{kl}(\ul{u})

Call signature:

dw_lin_elastic

(material, virtual, state)

(material, parameter_1, parameter_2)

Arguments 1:
  • material : D_{ijkl}

  • virtual : \ul{v}

  • state : \ul{u}

Arguments 2:
  • material : D_{ijkl}

  • parameter_1 : \ul{w}

  • parameter_2 : \ul{u}

arg_shapes = {'material': 'S, S', 'parameter_1': 'D', 'parameter_2': 'D', 'state': 'D', 'virtual': ('D', 'state')}
arg_types = (('material', 'virtual', 'state'), ('material', 'parameter_1', 'parameter_2'))
get_eval_shape(mat, virtual, state, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
get_fargs(mat, virtual, state, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
modes = ('weak', 'eval')
name = 'dw_lin_elastic'
set_arg_types()[source]
class sfepy.terms.terms_elastic.LinearPrestressTerm(name, arg_str, integral, region, **kwargs)[source]

Linear prestress term, with the prestress \sigma_{ij} given either in the usual vector form exploiting symmetry: in 3D it has 6 components with the indices ordered as [11, 22, 33, 12, 13, 23], in 2D it has 3 components with the indices ordered as [11, 22, 12], or in the matrix (possibly non-symmetric) form. Can be evaluated.

Definition:

\int_{\Omega} \sigma_{ij} e_{ij}(\ul{v})

Call signature:

dw_lin_prestress

(material, virtual)

(material, parameter)

Arguments 1:
  • material : \sigma_{ij}

  • virtual : \ul{v}

Arguments 2:
  • material : \sigma_{ij}

  • parameter : \ul{u}

arg_shapes = [{'material': 'S, 1', 'parameter': 'D', 'virtual': ('D', None)}, {'material': 'D, D'}]
arg_types = (('material', 'virtual'), ('material', 'parameter'))
d_lin_prestress(out, strain, mat, vg, fmode)[source]
get_eval_shape(mat, virtual, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
get_fargs(mat, virtual, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
modes = ('weak', 'eval')
name = 'dw_lin_prestress'
set_arg_types()[source]
class sfepy.terms.terms_elastic.LinearSpringTerm(name, arg_str, integral, region, **kwargs)[source]

Linear spring element.

Definition:

\ul{f}^{(i)} = - \ul{f}^{(j)} = k (\ul{u}^{(j)} - \ul{u}^{(i)})\\
\quad \forall \mbox{ elements } T_K^{i,j}\\
\mbox{ in a region connecting nodes } i, j

Call signature:

dw_lin_spring

(material, virtual, state)

Arguments 1:
  • material : k

  • virtual : \ul{v}

  • state : \ul{u}

arg_shapes = {'material': '1, 1', 'state': 'D', 'virtual': ('D', 'state')}
arg_types = ('material', 'virtual', 'state')
static function(out, stiffness, vec, diff_var)[source]
geometries = ['1_2', '2_1_2', '3_1_2']
get_fargs(mat, virtual, state, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
integration_order = 0
name = 'dw_lin_spring'
class sfepy.terms.terms_elastic.LinearStrainFiberTerm(name, arg_str, integral, region, **kwargs)[source]

Linear (pre)strain fiber term with the unit direction vector \ul{d}.

Definition:

\int_{\Omega} D_{ijkl} e_{ij}(\ul{v}) \left(d_k d_l\right)

Call signature:

dw_lin_strain_fib

(material_1, material_2, virtual)

Arguments:
  • material_1 : D_{ijkl}

  • material_2 : \ul{d}

  • virtual : \ul{v}

arg_shapes = {'material_1': 'S, S', 'material_2': 'D, 1', 'virtual': ('D', None)}
arg_types = ('material_1', 'material_2', 'virtual')
static function(out, mtx_d, mat, cmap)
get_fargs(mat1, mat2, virtual, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
name = 'dw_lin_strain_fib'
class sfepy.terms.terms_elastic.LinearTrussInternalForceTerm(name, arg_str, integral, region, **kwargs)[source]

Evaluate internal force in the element direction. To be used with ‘el_avg’ or ‘qp’ evaluation modes which give the same results. The material parameter EA is equal to Young modulus times element coss-section. The internal force is given by F^{(i)} = -F^{(j)} = EA / l (U^{(j)} - U^{(i)}), where l is the element length and U, F are the nodal displacements and the nodal forces in the element direction.

Definition:

F = EA / l (U^{(j)} - U^{(i)})\\
\quad \forall \mbox{ elements } T_K^{i,j}\\
\mbox{ in a region connecting nodes } i, j

Call signature:

ev_lin_truss_force

(material, parameter)

Arguments:
  • material : EA

  • parameter : \ul{w}

arg_shapes = {'material': '1, 1', 'parameter': 'D'}
arg_types = ('material', 'parameter')
static function(out, mat, vec, mtx_t)[source]
geometries = ['1_2', '2_1_2', '3_1_2']
get_eval_shape(mat, parameter, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
get_fargs(mat, parameter, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
integration_order = 0
name = 'ev_lin_truss_force'
class sfepy.terms.terms_elastic.LinearTrussTerm(name, arg_str, integral, region, **kwargs)[source]

Evaluate internal force in the element direction. To be used with ‘el_avg’ or ‘qp’ evaluation modes which give the same results. The material parameter EA is equal to Young modulus times element coss-section. The internal force is given by F^{(i)} = -F^{(j)} = EA / l (U^{(j)} - U^{(i)}), where l is the element length and U, F are the nodal displacements and the nodal forces in the element direction.

Definition:

F^{(i)} = -F^{(j)} = EA / l (U^{(j)} - U^{(i)})\\
\quad \forall \mbox{ elements } T_K^{i,j}\\
\mbox{ in a region connecting nodes } i, j

Call signature:

dw_lin_truss

(material, virtual, state)

Arguments:
  • material : EA

  • parameter : \ul{w}

arg_shapes = {'material': '1, 1', 'state': 'D', 'virtual': ('D', 'state')}
arg_types = ('material', 'virtual', 'state')
static function(out, mat, vec, mtx_t, length, diff_var)[source]
geometries = ['1_2', '2_1_2', '3_1_2']
get_fargs(mat, virtual, state, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
static get_mtx_t_and_length(coors, dx=None)[source]
integration_order = 0
name = 'dw_lin_truss'
class sfepy.terms.terms_elastic.NonsymElasticTerm(name, arg_str, integral, region, **kwargs)[source]

Elasticity term with non-symmetric gradient. The indices of matrix D_{ijkl} are ordered as [11, 12, 13, 21, 22, 23, 31, 32, 33] in 3D and as [11, 12, 21, 22] in 2D.

Definition:

\int_{\Omega} \ull{D} \nabla\ul{u} : \nabla\ul{v}

Call signature:

dw_nonsym_elastic

(material, virtual, state)

(material, parameter_1, parameter_2)

Arguments 1:
  • material : \ull{D}

  • virtual : \ul{v}

  • state : \ul{u}

Arguments 2:
  • material : \ull{D}

  • parameter_1 : \ul{w}

  • parameter_2 : \ul{u}

arg_shapes = {'material': 'D2, D2', 'parameter_1': 'D', 'parameter_2': 'D', 'state': 'D', 'virtual': ('D', 'state')}
arg_types = (('material', 'virtual', 'state'), ('material', 'parameter_1', 'parameter_2'))
geometries = ['2_3', '2_4', '3_4', '3_8']
get_eval_shape(mat, virtual, state, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
get_fargs(mat, virtual, state, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
modes = ('weak', 'eval')
name = 'dw_nonsym_elastic'
set_arg_types()[source]
class sfepy.terms.terms_elastic.SDLinearElasticTerm(name, arg_str, integral, region, **kwargs)[source]

Sensitivity analysis of the linear elastic term.

Definition:

\int_{\Omega} \hat{D}_{ijkl}\ e_{ij}(\ul{v}) e_{kl}(\ul{u})

\hat{D}_{ijkl} = D_{ijkl}(\nabla \cdot \ul{\Vcal})
- D_{ijkq}{\partial \Vcal_l \over \partial x_q}
- D_{iqkl}{\partial \Vcal_j \over \partial x_q}

Call signature:

ev_sd_lin_elastic

(material, parameter_w, parameter_u, parameter_mv)

Arguments:
  • material : D_{ijkl}

  • parameter_w : \ul{w}

  • parameter_u : \ul{u}

  • parameter_mv : \ul{\Vcal}

arg_shapes = {'material': 'S, S', 'parameter_mv': 'D', 'parameter_u': 'D', 'parameter_w': 'D'}
arg_types = ('material', 'parameter_w', 'parameter_u', 'parameter_mv')
static function(out, coef, grad_v, grad_u, grad_w, mtx_d, cmap)
geometries = ['2_3', '2_4', '3_4', '3_8']
get_eval_shape(mat, par_w, par_u, par_mv, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
get_fargs(mat, par_w, par_u, par_mv, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
name = 'ev_sd_lin_elastic'