sfepy.terms.terms_hyperelastic_base module

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

Deformation gradient \ull{F} in quadrature points for term_mode=’def_grad’ (default) or the jacobian J if term_mode=’jacobian’.

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

Definition:

\ull{F} = \pdiff{\ul{x}}{\ul{X}}|_{qp}
= \ull{I} + \pdiff{\ul{u}}{\ul{X}}|_{qp} \;, \\
\ul{x} = \ul{X} + \ul{u} \;, J = \det{(\ull{F})}

Call signature:

ev_def_grad

(parameter)

Arguments:
  • parameter : \ul{u}

arg_shapes = {'parameter': 'D'}
arg_types = ('parameter',)
static function(out, vec, vg, econn, term_mode, 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_def_grad'
class sfepy.terms.terms_hyperelastic_base.HyperElasticBase(*args, **kwargs)[source]

Base class for all hyperelastic terms in TL/UL formulation.

HyperElasticBase.__call__() computes element contributions given either stress (-> residual) or tangent modulus (-> tangent sitffnes matrix), i.e. constitutive relation type (CRT) related data. The CRT data are computed in subclasses implementing particular CRT (e.g. neo-Hookean material), in self.compute_crt_data().

Modes:

  • 0: total formulation

  • 1: updated formulation

Notes

This is not a proper Term!

arg_shapes = {'material': '1, 1', 'state': 'D', 'virtual': ('D', 'state')}
arg_types = ('material', 'virtual', 'state')
compute_stress(mat, family_data, **kwargs)[source]
compute_tan_mod(mat, family_data, **kwargs)[source]
static function(out, fun, *args)[source]
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]
static integrate(out, val_qp, vg, fmode)[source]
integration = ('cell', 'facet_extra')
class sfepy.terms.terms_hyperelastic_base.HyperElasticFamilyData(**kwargs)[source]

Base class for hyperelastic family data.

The common (family) data are cached in the evaluate cache of state variable.

data_shapes = {'det_f': ('n_el', 'n_qp', 1, 1), 'green_strain': ('n_el', 'n_qp', 'sym', 1), 'in2_b': ('n_el', 'n_qp', 1, 1), 'in2_c': ('n_el', 'n_qp', 1, 1), 'inv_f': ('n_el', 'n_qp', 'dim', 'dim'), 'mtx_f': ('n_el', 'n_qp', 'dim', 'dim'), 'sym_b': ('n_el', 'n_qp', 'sym', 1), 'sym_c': ('n_el', 'n_qp', 'sym', 1), 'sym_inv_c': ('n_el', 'n_qp', 'sym', 1), 'tr_b': ('n_el', 'n_qp', 1, 1), 'tr_c': ('n_el', 'n_qp', 1, 1)}
init_data_struct(state_shape, name='family_data')[source]