sfepy.terms.terms_multilinear module

class sfepy.terms.terms_multilinear.ECauchyStressTerm(*args, **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].

Definition:

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

Call signature:

de_cauchy_stress

(material, parameter)

Arguments:
  • material : D_{ijkl}

  • parameter : \ul{w}

arg_shapes = {'material': 'S, S', 'parameter': 'D'}
arg_types = ('material', 'parameter')
get_function(mat, parameter, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
name = 'de_cauchy_stress'
class sfepy.terms.terms_multilinear.EConvectTerm(*args, **kwargs)[source]

Nonlinear convective term.

Definition:

\int_{\Omega} ((\ul{u} \cdot \nabla) \ul{u}) \cdot \ul{v}

Call signature:

de_convect

(virtual, state)

(parameter_1, parameter_2)

Arguments:
  • virtual/parameter_1: \ul{v}

  • state/parameter_2: \ul{u}

arg_shapes = {'parameter_1': 'D', 'parameter_2': 'D', 'state': 'D', 'virtual': ('D', 'state')}
arg_types = (('virtual', 'state'), ('parameter_1', 'parameter_2'))
get_function(virtual, state, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
modes = ('weak', 'eval')
name = 'de_convect'
class sfepy.terms.terms_multilinear.EDiffusionTerm(*args, **kwargs)[source]

General diffusion term.

Definition:

\int_{\Omega} K_{ij} \nabla_i q\, \nabla_j p

Call signature:

de_diffusion

(material, virtual, state)

(material, parameter_1, parameter_2)

Arguments:
  • material: K_{ij}

  • virtual/parameter_1: q

  • state/parameter_2: p

arg_shapes = {'material': 'D, D', 'parameter_1': 1, 'parameter_2': 1, 'state': 1, 'virtual': (1, 'state')}
arg_types = (('material', 'virtual', 'state'), ('material', 'parameter_1', 'parameter_2'))
get_function(mat, vvar, svar, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
modes = ('weak', 'eval')
name = 'de_diffusion'
class sfepy.terms.terms_multilinear.EDivGradTerm(*args, **kwargs)[source]

Vector field diffusion term.

Definition:

\int_{\Omega} \nabla \ul{v} : \nabla \ul{u} \mbox{ , }
\int_{\Omega} \nu\ \nabla \ul{v} : \nabla \ul{u}

Call signature:

de_div_grad

(opt_material, virtual, state)

(opt_material, parameter_1, parameter_2)

Arguments:
  • material: \nu (viscosity, optional)

  • virtual/parameter_1: \ul{v}

  • state/parameter_2: \ul{u}

arg_shapes = [{'opt_material': '1, 1', 'parameter_1': 'D', 'parameter_2': 'D', 'state': 'D', 'virtual': ('D', 'state')}, {'opt_material': None}]
arg_types = (('opt_material', 'virtual', 'state'), ('opt_material', 'parameter_1', 'parameter_2'))
get_function(mat, virtual, state, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
modes = ('weak', 'eval')
name = 'de_div_grad'
class sfepy.terms.terms_multilinear.EDivTerm(*args, **kwargs)[source]

Weighted divergence term.

Definition:

\int_{\Omega} \nabla \cdot \ul{v} \mbox { , }
\int_{\Omega} c \nabla \cdot \ul{v}

Call signature:

de_div

(opt_material, virtual)

(opt_material, parameter)

Arguments:
  • material: c (optional)

  • virtual/parameter: \ul{v}

arg_shapes = [{'opt_material': '1, 1', 'parameter': 'D', 'virtual': ('D', None)}, {'opt_material': None}]
arg_types = (('opt_material', 'virtual'), ('opt_material', 'parameter'))
get_function(mat, virtual, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
modes = ('weak', 'eval')
name = 'de_div'
class sfepy.terms.terms_multilinear.EDotTerm(*args, **kwargs)[source]

Volume and surface L^2(\Omega) weighted dot product for both scalar and vector fields. Can be evaluated. Can use derivatives.

Definition:

\int_{\cal{D}} q p \mbox{ , } \int_{\cal{D}} \ul{v} \cdot \ul{u}\\
\int_{\cal{D}} c q p \mbox{ , } \int_{\cal{D}} c \ul{v} \cdot \ul{u}\\
\int_{\cal{D}} \ul{v} \cdot (\ull{c}\, \ul{u})

Call signature:

de_dot

(opt_material, virtual, state)

(opt_material, parameter_1, parameter_2)

Arguments:
  • material: c or \ull{c} (optional)

  • virtual/parameter_1: q or \ul{v}

  • state/parameter_2: p or \ul{u}

arg_shapes = [{'opt_material': '1, 1', 'parameter_1': 1, 'parameter_2': 1, 'state': 1, 'virtual': (1, 'state')}, {'opt_material': None}, {'opt_material': '1, 1', 'parameter_1': 'D', 'parameter_2': 'D', 'state': 'D', 'virtual': ('D', 'state')}, {'opt_material': 'D, D'}, {'opt_material': None}, {'opt_material': '1, 1', 'parameter_1': 'N', 'parameter_2': 'N', 'state': 'N', 'virtual': ('N', 'state')}, {'opt_material': 'N, N'}, {'opt_material': None}]
arg_types = (('opt_material', 'virtual', 'state'), ('opt_material', 'parameter_1', 'parameter_2'))
get_function(mat, virtual, state, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
integration = ('cell', 'facet')
modes = ('weak', 'eval')
name = 'de_dot'
class sfepy.terms.terms_multilinear.EGradTerm(*args, **kwargs)[source]

Weighted gradient term.

Definition:

\int_{\Omega} \nabla \ul{v} \mbox { , }
\int_{\Omega} c \nabla \ul{v} \mbox { , }
\int_{\Omega} \ul{c} \cdot \nabla \ul{v} \mbox { , }
\int_{\Omega} \ull{c} \cdot \nabla \ul{v}

Call signature:

de_grad

(opt_material, parameter)

Arguments:
  • material: c (optional)

  • virtual/parameter: \ul{v}

arg_shapes = [{'opt_material': '1, 1', 'parameter': 'N'}, {'opt_material': 'N, 1'}, {'opt_material': 'N, N'}, {'opt_material': None}]
arg_types = ('opt_material', 'parameter')
get_function(mat, virtual, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
name = 'de_grad'
class sfepy.terms.terms_multilinear.EIntegrateOperatorTerm(*args, **kwargs)[source]

Volume and surface integral of a test function weighted by a scalar function c.

Definition:

\int_{\cal{D}} q \mbox{ or } \int_{\cal{D}} c q

Call signature:

de_integrate

(opt_material, virtual)

Arguments:
  • material : c (optional)

  • virtual : q

arg_shapes = [{'opt_material': '1, 1', 'virtual': (1, None)}, {'opt_material': None}]
arg_types = ('opt_material', 'virtual')
get_function(mat, virtual, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
integration = ('cell', 'facet')
name = 'de_integrate'
class sfepy.terms.terms_multilinear.ELaplaceTerm(*args, **kwargs)[source]

Laplace term with c coefficient. Can be evaluated. Can use derivatives.

Definition:

\int_{\Omega} \nabla q \cdot \nabla p \mbox{ , }
\int_{\Omega} c \nabla q \cdot \nabla p

Call signature:

de_laplace

(opt_material, virtual, state)

(opt_material, parameter_1, parameter_2)

Arguments:
  • material: c

  • virtual/parameter_1: q

  • state/parameter_2: p

arg_shapes = [{'opt_material': '1, 1', 'parameter_1': 1, 'parameter_2': 1, 'state': 1, 'virtual': (1, 'state')}, {'opt_material': None}]
arg_types = (('opt_material', 'virtual', 'state'), ('opt_material', 'parameter_1', 'parameter_2'))
get_function(mat, virtual, state, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
modes = ('weak', 'eval')
name = 'de_laplace'
class sfepy.terms.terms_multilinear.ELinearConvectTerm(*args, **kwargs)[source]

Linearized convective term.

Definition:

\int_{\Omega} ((\ul{w} \cdot \nabla) \ul{u}) \cdot \ul{v}

Call signature:

de_lin_convect

(virtual, parameter, state)

(parameter_1, parameter_2, parameter_3)

Arguments:
  • virtual/parameter_1: \ul{v}

  • parameter/parameter_2: \ul{w}

  • state/parameter_3: \ul{u}

arg_shapes = {'parameter': 'D', 'parameter_1': 'D', 'parameter_2': 'D', 'parameter_3': 'D', 'state': 'D', 'virtual': ('D', 'state')}
arg_types = (('virtual', 'parameter', 'state'), ('parameter_1', 'parameter_2', 'parameter_3'))
get_function(virtual, parameter, state, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
modes = ('weak', 'eval')
name = 'de_lin_convect'
class sfepy.terms.terms_multilinear.ELinearElasticTerm(*args, **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].

Definition:

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

Call signature:

de_lin_elastic

(material, virtual, state)

(material, parameter_1, parameter_2)

Arguments:
  • material: D_{ijkl}

  • virtual/parameter_1: \ul{v}

  • state/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_function(mat, virtual, state, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
modes = ('weak', 'eval')
name = 'de_lin_elastic'
class sfepy.terms.terms_multilinear.ELinearTractionTerm(*args, **kwargs)[source]

Linear traction term. The material parameter can have one of the following shapes:

  • 1 or (1, 1) - a given scalar pressure

  • (D, 1) - a traction vector

  • (S, 1) or (D, D) - a given stress in symmetric or non-symmetric tensor storage (in symmetric storage indicies are order as follows: 2D: [11, 22, 12], 3D: [11, 22, 33, 12, 13, 23])

Definition:

\int_{\Gamma} \ul{v} \cdot \ul{n}
    \mbox{ , }
\int_{\Gamma} c\, \ul{v} \cdot \ul{n}\\
\int_{\Gamma} \ul{v} \cdot (\ull{\sigma}\, \ul{n})
    \mbox{ , }
\int_{\Gamma} \ul{v} \cdot \ul{f}

Call signature:

de_surface_ltr

(opt_material, virtual)

(opt_material, parameter)

Arguments:
  • material: c, \ul{f}, \ul{\sigma} or \ull{\sigma}

  • virtual/parameter: \ul{v}

arg_shapes = [{'opt_material': 'S, 1', 'parameter': 'D', 'virtual': ('D', None)}, {'opt_material': None}, {'opt_material': '1, 1'}, {'opt_material': 'D, 1'}, {'opt_material': 'D, D'}]
arg_types = (('opt_material', 'virtual'), ('opt_material', 'parameter'))
get_function(traction, vvar, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
integration = 'facet'
modes = ('weak', 'eval')
name = 'de_surface_ltr'
class sfepy.terms.terms_multilinear.ENonPenetrationPenaltyTerm(*args, **kwargs)[source]

Non-penetration condition in the weak sense using a penalty.

Definition:

\int_{\Gamma} c (\ul{n} \cdot \ul{v}) (\ul{n} \cdot \ul{u})

Call signature:

de_non_penetration_p

(material, virtual, state)

Arguments:
  • material : c

  • virtual : \ul{v}

  • state : \ul{u}

arg_shapes = {'material': '1, 1', 'state': 'D', 'virtual': ('D', 'state')}
arg_types = ('material', 'virtual', 'state')
get_function(mat, virtual, state, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
integration = 'facet'
name = 'de_non_penetration_p'
class sfepy.terms.terms_multilinear.ENonSymElasticTerm(*args, **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{v} : \nabla \ul{u}

Call signature:

de_nonsym_elastic

(material, virtual, state)

(material, parameter_1, parameter_2)

Arguments:
  • material: \ull{D}

  • virtual/parameter_1: \ul{v}

  • state/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'))
get_function(mat, virtual, state, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
modes = ('weak', 'eval')
name = 'de_nonsym_elastic'
class sfepy.terms.terms_multilinear.EScalarDotMGradScalarTerm(*args, **kwargs)[source]

Volume dot product of a scalar gradient dotted with a material vector with a scalar.

Definition:

\int_{\Omega} q \ul{y} \cdot \nabla p \mbox{ , }
\int_{\Omega} p \ul{y} \cdot \nabla q

Call signature:

de_s_dot_mgrad_s

(material, virtual, state)

(material, state, virtual)

(material, parameter_1, parameter_2)

Arguments 1:
  • material : \ul{y}

  • virtual : q

  • state : p

Arguments 2:
  • material : \ul{y}

  • state : p

  • virtual : q

arg_shapes = [{'material': 'D, 1', 'parameter_1': 1, 'parameter_2': 1, 'state/grad_state': 1, 'state/grad_virtual': 1, 'virtual/grad_state': (1, None), 'virtual/grad_virtual': (1, None)}]
arg_types = (('material', 'virtual', 'state'), ('material', 'state', 'virtual'), ('material', 'parameter_1', 'parameter_2'))
get_function(mat, var1, var2, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
modes = ('grad_state', 'grad_virtual', 'eval')
name = 'de_s_dot_mgrad_s'
class sfepy.terms.terms_multilinear.EStokesTerm(*args, **kwargs)[source]

Stokes problem coupling term. Corresponds to weak forms of gradient and divergence terms.

Definition:

\int_{\Omega} p\, \nabla \cdot \ul{v} \mbox{ , }
\int_{\Omega} q\, \nabla \cdot \ul{u}\\
\int_{\Omega} c\, p\, \nabla \cdot \ul{v} \mbox{ , }
\int_{\Omega} c\, q\, \nabla \cdot \ul{u}

Call signature:

de_stokes

(opt_material, virtual, state)

(opt_material, state, virtual)

(opt_material, parameter_v, parameter_s)

Arguments 1:
  • material: c (optional)

  • virtual/parameter_v: \ul{v}

  • state/parameter_s: p

Arguments 2:
  • material : c (optional)

  • state : \ul{u}

  • virtual : q

arg_shapes = [{'opt_material': '1, 1', 'parameter_s': 1, 'parameter_v': 'D', 'state/div': 'D', 'state/grad': 1, 'virtual/div': (1, None), 'virtual/grad': ('D', None)}, {'opt_material': None}]
arg_types = (('opt_material', 'virtual', 'state'), ('opt_material', 'state', 'virtual'), ('opt_material', 'parameter_v', 'parameter_s'))
get_function(coef, vvar, svar, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
modes = ('grad', 'div', 'eval')
name = 'de_stokes'
class sfepy.terms.terms_multilinear.ETermBase(*args, **kwargs)[source]

Reserved letters:

c .. cells q .. quadrature points d-h .. DOFs axes r-z .. auxiliary axes

Layout specification letters:

c .. cells q .. quadrature points v .. variable component - matrix form (v, d) -> vector v*d g .. gradient component d .. local DOF (basis, node) 0 .. all material axes

build_expression(texpr, *eargs, mode=None, diff_var=None)[source]
can_backend = {'dask_single': <module 'dask.array' from '/home/eldaran/.local/lib/python3.10/site-packages/dask/array/__init__.py'>, 'dask_threads': <module 'dask.array' from '/home/eldaran/.local/lib/python3.10/site-packages/dask/array/__init__.py'>, 'jax': <module 'jax.numpy' from '/home/eldaran/.local/lib/python3.10/site-packages/jax/numpy/__init__.py'>, 'jax_vmap': <module 'jax.numpy' from '/home/eldaran/.local/lib/python3.10/site-packages/jax/numpy/__init__.py'>, 'numpy': <module 'numpy' from '/home/eldaran/.local/lib/python3.10/site-packages/numpy/__init__.py'>, 'numpy_loop': <module 'numpy' from '/home/eldaran/.local/lib/python3.10/site-packages/numpy/__init__.py'>, 'numpy_qloop': <module 'numpy' from '/home/eldaran/.local/lib/python3.10/site-packages/numpy/__init__.py'>, 'opt_einsum': <module 'opt_einsum' from '/home/eldaran/.local/lib/python3.10/site-packages/opt_einsum/__init__.py'>, 'opt_einsum_dask_single': <module 'dask.array' from '/home/eldaran/.local/lib/python3.10/site-packages/dask/array/__init__.py'>, 'opt_einsum_dask_threads': <module 'dask.array' from '/home/eldaran/.local/lib/python3.10/site-packages/dask/array/__init__.py'>, 'opt_einsum_loop': <module 'opt_einsum' from '/home/eldaran/.local/lib/python3.10/site-packages/opt_einsum/__init__.py'>, 'opt_einsum_qloop': <module 'opt_einsum' from '/home/eldaran/.local/lib/python3.10/site-packages/opt_einsum/__init__.py'>}
clear_cache()[source]
eval_complex(shape, fargs, mode='eval', term_mode=None, diff_var=None, **kwargs)[source]
eval_real(shape, fargs, mode='eval', term_mode=None, diff_var=None, **kwargs)[source]
static function_silent(out, eval_einsum, *args)[source]
static function_timer(out, eval_einsum, *args)[source]
get_eval_shape(*args, **kwargs)[source]
get_fargs(*args, **kwargs)[source]
get_normals(arg)[source]
get_operands(diff_var)[source]
get_paths(expressions, operands)[source]
layout_letters = 'cqgvd0'
make_function(texpr, *args, mode=None, diff_var=None)[source]
set_backend(backend='numpy', optimize=True, layout=None, **kwargs)[source]
set_verbosity(verbosity=None)[source]
verbosity = 0
class sfepy.terms.terms_multilinear.ExpressionArg(**kwargs)[source]
static from_term_arg(arg, term)[source]
get_bf(expr_cache)[source]
get_dofs(cache, expr_cache, oname)[source]
class sfepy.terms.terms_multilinear.ExpressionBuilder(n_add, cache)[source]
add_arg_dofs(iin, ein, name, n_components, iia=None)[source]
add_bf(iin, ein, name, cell_dependent=False)[source]
add_bfg(iin, ein, name)[source]
add_cell_scalar(name, cname)[source]
add_eye(iic, ein, name, iia=None)[source]
add_material_arg(arg, ii, ein)[source]
add_psg(iic, ein, name, iia=None)[source]
add_pvg(iic, ein, name, iia=None)[source]
add_qp_scalar(name, cname)[source]
add_state_arg(arg, ii, ein, modifier, diff_var)[source]
add_virtual_arg(arg, ii, ein, modifier)[source]
apply_layout(layout, operands, defaults=None, verbosity=0)[source]
build(texpr, *args, mode=None, diff_var=None)[source]
get_expressions(subscripts=None)[source]
static join_subscripts(subscripts, out_subscripts)[source]
letters = 'defgh'
make_eye(size)[source]
make_psg(dim)[source]
make_pvg(dim)[source]
print_shapes(subscripts, operands)[source]
transform(subscripts, operands, transformation='loop', **kwargs)[source]
class sfepy.terms.terms_multilinear.StokesTractionTerm(*args, **kwargs)[source]

Surface traction term for Stokes flow problem.

Definition:

\int_{\Gamma} \nu \ul{v}\cdot(\nabla \ul{u} \cdot \ul{n})

Call signature:

de_stokes_traction

(opt_material, virtual, state)

(opt_material, parameter_1, parameter_2)

Arguments:
  • material: \nu (viscosity, optional)

  • virtual/parameter_1: \ul{v}

  • state/parameter_2: \ul{u}

arg_shapes = [{'opt_material': '1, 1', 'parameter_1': 'D', 'parameter_2': 'D', 'state': 'D', 'virtual': ('D', 'state')}, {'opt_material': None}]
arg_types = (('opt_material', 'virtual', 'state'), ('opt_material', 'parameter_1', 'parameter_2'))
get_function(mat, virtual, state, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
integration = 'facet_extra'
modes = ('weak', 'eval')
name = 'de_stokes_traction'
class sfepy.terms.terms_multilinear.SurfaceFluxOperatorTerm(*args, **kwargs)[source]

Surface flux operator term.

Definition:

\int_{\Gamma} q \ul{n} \cdot \ull{K} \cdot \nabla p \mbox{ , }
\int_{\Gamma} p \ul{n} \cdot \ull{K} \cdot \nabla q

Call signature:

de_surface_flux

(material, virtual, state)

(material, state, virtual)

(material, parameter_1, parameter_2)

Arguments 1:
  • material : \ull{K}

  • virtual : q

  • state : p

Arguments 2:
  • material : \ull{K}

  • state : p

  • virtual : q

arg_shapes = [{'material': 'D, D', 'parameter_1': 1, 'parameter_2': 1, 'state/grad_state': 1, 'state/grad_virtual': 1, 'virtual/grad_state': (1, None), 'virtual/grad_virtual': (1, None)}]
arg_types = (('material', 'virtual', 'state'), ('material', 'state', 'virtual'), ('material', 'parameter_1', 'parameter_2'))
get_function(mat, var1, var2, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
integration = 'facet_extra'
modes = ('grad_state', 'grad_virtual', 'eval')
name = 'de_surface_flux'
class sfepy.terms.terms_multilinear.SurfacePiezoFluxOperatorTerm(*args, **kwargs)[source]

Surface piezoelectric flux operator term.

Corresponds to the electric flux due to mechanically induced electrical displacement.

Definition:

\int_{\Gamma} q g_{kij} e_{ij}(\ul{u}) n_k \mbox{ , }
\int_{\Gamma} p g_{kij} e_{ij}(\ul{v}) n_k

Call signature:

de_surface_piezo_flux

(material, virtual, state)

(material, state, virtual)

(material, parameter_1, parameter_2)

Arguments 1:
  • material : g_{kij}

  • virtual/parameter_1 : q

  • state/parameter_2 : \ul{u}

Arguments 2:
  • material : g_{kij}

  • state : p

  • virtual : \ul{v}

arg_shapes = [{'material': 'D, S', 'parameter_1': 1, 'parameter_2': 'D', 'state/grad_state': 'D', 'state/grad_virtual': 1, 'virtual/grad_state': (1, None), 'virtual/grad_virtual': ('D', None)}]
arg_types = (('material', 'virtual', 'state'), ('material', 'state', 'virtual'), ('material', 'parameter_1', 'parameter_2'))
get_function(mat, var1, var2, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
integration = 'facet_extra'
modes = ('grad_state', 'grad_virtual', 'eval')
name = 'de_surface_piezo_flux'
class sfepy.terms.terms_multilinear.SurfacePiezoFluxTerm(*args, **kwargs)[source]

Surface piezoelectric flux term.

Corresponds to the electric flux due to mechanically induced electrical displacement.

Definition:

\int_{\Gamma} g_{kij} e_{ij}(\ul{u}) n_k

Call signature:

ev_surface_piezo_flux

(material, parameter)

Arguments 1:
  • material : g_{kij}

  • parameter : \ul{u}

arg_shapes = {'material': 'D, S', 'parameter': 'D'}
arg_types = ('material', 'parameter')
get_function(mat, var1, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
integration = 'facet_extra'
modes = ('eval',)
name = 'ev_surface_piezo_flux'
sfepy.terms.terms_multilinear.append_all(seqs, item, ii=None)[source]
sfepy.terms.terms_multilinear.collect_modifiers(modifiers)[source]
sfepy.terms.terms_multilinear.find_free_indices(indices)[source]
sfepy.terms.terms_multilinear.get_einsum_ops(eargs, ebuilder, expr_cache)[source]
sfepy.terms.terms_multilinear.get_loop_indices(subs, loop_index)[source]
sfepy.terms.terms_multilinear.get_output_shape(out_subscripts, subscripts, operands)[source]
sfepy.terms.terms_multilinear.get_sizes(indices, operands)[source]
sfepy.terms.terms_multilinear.get_slice_ops(subs, ops, loop_index)[source]
sfepy.terms.terms_multilinear.parse_term_expression(texpr)[source]
sfepy.terms.terms_multilinear.sym2nonsym(sym_obj, axes=[3])[source]