sfepy.terms.terms_sensitivity module

class sfepy.terms.terms_sensitivity.ESDDiffusionTerm(*args, **kwargs)[source]

Diffusion sensitivity analysis term.

Definition:

\int_{\Omega} \hat{K}_{ij} \nabla_i q\, \nabla_j p

\hat{K}_{ij} = K_{ij}\left(
    \delta_{ik}\delta_{jl} \nabla \cdot \ul{\Vcal}
  - \delta_{ik}{\partial \Vcal_j \over \partial x_l}
  - \delta_{jl}{\partial \Vcal_i \over \partial x_k}\right)

Call signature:

de_sd_diffusion

(material, virtual, state, parameter_mv)

(material, parameter_1, parameter_2, parameter_mv)

Arguments:
  • material: K_{ij}

  • virtual/parameter_1: q

  • state/parameter_2: p

  • parameter_mv: \ul{\Vcal}

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

Sensitivity (shape derivative) of diffusion term de_div_grad.

Definition:

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

\hat{I}_{ijkl} =
    \delta_{ik}\delta_{jl} \nabla \cdot \ul{\Vcal}
  - \delta_{ik}\delta_{js} {\partial \Vcal_l \over \partial x_s}
  - \delta_{is}\delta_{jl} {\partial \Vcal_k \over \partial x_s}

Call signature:

de_sd_div_grad

(opt_material, virtual, state, parameter_mv)

(opt_material, parameter_1, parameter_2, parameter_mv)

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

  • virtual/parameter_1: \ul{v}

  • state/parameter_2: \ul{u}

  • parameter_mv: \ul{\Vcal}

arg_shapes = [{'opt_material': '1, 1', 'virtual': ('D', 'state'), 'state': 'D', 'parameter_1': 'D', 'parameter_2': 'D', 'parameter_mv': 'D'}, {'opt_material': None}]
arg_types = (('opt_material', 'virtual', 'state', 'parameter_mv'), ('opt_material', 'parameter_1', 'parameter_2', 'parameter_mv'))
get_function(mat, vvar, svar, par_mv, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
modes = ('weak', 'eval')
name = 'de_sd_div_grad'
class sfepy.terms.terms_sensitivity.ESDDotTerm(*args, **kwargs)[source]

Sensitivity (shape derivative) of dot product of scalars or vectors.

Definition:

\int_\Omega q p (\nabla \cdot \ul{\Vcal}) \mbox{ , }
\int_\Omega (\ul{v} \cdot \ul{u}) (\nabla \cdot \ul{\Vcal})\\
\int_\Omega c q p (\nabla \cdot \ul{\Vcal}) \mbox{ , }
\int_\Omega c (\ul{v} \cdot \ul{u}) (\nabla \cdot \ul{\Vcal})\\
\int_\Omega \ul{v} \cdot (\ull{M}\, \ul{u}) (\nabla \cdot \ul{\Vcal})

Call signature:

de_sd_dot

(opt_material, virtual, state, parameter_mv)

(opt_material, parameter_1, parameter_2, parameter_mv)

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

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

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

  • parameter_mv : \ul{\Vcal}

arg_shapes = [{'opt_material': '1, 1', 'virtual': (1, 'state'), 'state': 1, 'parameter_1': 1, 'parameter_2': 1, 'parameter_mv': 'D'}, {'opt_material': None}, {'opt_material': '1, 1', 'virtual': ('D', 'state'), 'state': 'D', 'parameter_1': 'D', 'parameter_2': 'D', 'parameter_mv': 'D'}, {'opt_material': 'D, D'}, {'opt_material': None}]
arg_types = (('opt_material', 'virtual', 'state', 'parameter_mv'), ('opt_material', 'parameter_1', 'parameter_2', 'parameter_mv'))
get_function(mat, vvar, svar, par_mv, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
modes = ('weak', 'eval')
name = 'de_sd_dot'
class sfepy.terms.terms_sensitivity.ESDLinearElasticTerm(*args, **kwargs)[source]

Sensitivity analysis of the linear elastic term. D_{ijkl} can be given in symmetric or non-symmetric form.

Definition:

\int_{\Omega} \hat{D}_{ijkl}
{\partial v_i \over \partial x_j}
{\partial u_k \over \partial x_l}

\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:

de_sd_lin_elastic

(material, virtual, state, parameter_mv)

(material, parameter_1, parameter_2, parameter_mv)

Arguments 1:
  • material : \ull{D}

  • virtual/parameter_v : \ul{v}

  • state/parameter_s : \ul{u}

  • parameter_mv : \ul{\Vcal}

arg_shapes = [{'material': 'S, S', 'virtual': ('D', 'state'), 'state': 'D', 'parameter_1': 'D', 'parameter_2': 'D', 'parameter_mv': 'D'}, {'material': 'D2, D2'}]
arg_types = (('material', 'virtual', 'state', 'parameter_mv'), ('material', 'parameter_1', 'parameter_2', 'parameter_mv'))
geometries = ['2_3', '2_4', '3_4', '3_8']
get_function(mat, vvar, svar, par_mv, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
modes = ('weak', 'eval')
name = 'de_sd_lin_elastic'
class sfepy.terms.terms_sensitivity.ESDLinearTractionTerm(*args, **kwargs)[source]

Sensitivity of the linear traction term.

Definition:

\int_{\Gamma} \ul{v} \cdot \left[\left(\ull{\hat{\sigma}}\,
\nabla \cdot \ul{\cal{V}} - \ull{{\hat\sigma}}\, \nabla \ul{\cal{V}}
\right)\ul{n}\right]

\ull{\hat\sigma} = \ull{I} \mbox{ , }
\ull{\hat\sigma} = c\,\ull{I} \mbox{ or }
\ull{\hat\sigma} = \ull{\sigma}

Call signature:

de_sd_surface_ltr

(opt_material, virtual, parameter_mv)

(opt_material, parameter, parameter_mv)

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

  • virtual/parameter: \ul{v}

  • parameter_mv: \ul{\Vcal}

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

Sensitivity (shape derivative) of the piezoelectric coupling term.

Definition:

\int_{\Omega} \hat{g}_{kij}\ e_{ij}(\ul{v}) \nabla_k p \mbox{ , }
\int_{\Omega} \hat{g}_{kij}\ e_{ij}(\ul{u}) \nabla_k q

\hat{g}_{kij} = g_{kij}(\nabla \cdot \ul{\Vcal})
- g_{kil}{\partial \Vcal_j \over \partial x_l}
- g_{lij}{\partial \Vcal_k \over \partial x_l}

Call signature:

de_sd_piezo_coupling

(material, virtual, state, parameter_mv)

(material, state, virtual, parameter_mv)

(material, parameter_v, parameter_s, parameter_mv)

Arguments 1:
  • material : g_{kij}

  • virtual/parameter_v : \ul{v}

  • state/parameter_s : p

  • parameter_mv : \ul{\Vcal}

Arguments 2:
  • material : g_{kij}

  • state : \ul{u}

  • virtual : q

  • parameter_mv : \ul{\Vcal}

arg_shapes = {'material': 'D, S', 'parameter_mv': 'D', 'parameter_s': 1, 'parameter_v': 'D', 'state/div': 'D', 'state/grad': 1, 'virtual/div': (1, None), 'virtual/grad': ('D', None)}
arg_types = (('material', 'virtual', 'state', 'parameter_mv'), ('material', 'state', 'virtual', 'parameter_mv'), ('material', 'parameter_v', 'parameter_s', 'parameter_mv'))
geometries = ['2_3', '2_4', '3_4', '3_8']
get_function(mat, vvar, svar, par_mv, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
modes = ('grad', 'div', 'eval')
name = 'de_sd_piezo_coupling'
class sfepy.terms.terms_sensitivity.ESDStokesTerm(*args, **kwargs)[source]

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

Definition:

\int_{\Omega} p\, \hat{I}_{ij} {\partial v_i \over \partial x_j} \mbox{ , }
\int_{\Omega} q\, \hat{I}_{ij} {\partial u_i \over \partial x_j}

\hat{I}_{ij} = \delta_{ij} \nabla \cdot \Vcal
  - {\partial \Vcal_j \over \partial x_i}

Call signature:

de_sd_stokes

(opt_material, virtual, state, parameter_mv)

(opt_material, state, virtual, parameter_mv)

(opt_material, parameter_v, parameter_s, parameter_mv)

Arguments 1:
  • virtual/parameter_v: \ul{v}

  • state/parameter_s: p

  • parameter_mv: \ul{\Vcal}

Arguments 2:
  • state : \ul{u}

  • virtual : q

  • parameter_mv: \ul{\Vcal}

arg_shapes = [{'opt_material': '1, 1', 'virtual/grad': ('D', None), 'state/grad': 1, 'virtual/div': (1, None), 'state/div': 'D', 'parameter_v': 'D', 'parameter_s': 1, 'parameter_mv': 'D'}, {'opt_material': None}]
arg_types = (('opt_material', 'virtual', 'state', 'parameter_mv'), ('opt_material', 'state', 'virtual', 'parameter_mv'), ('opt_material', 'parameter_v', 'parameter_s', 'parameter_mv'))
get_function(coef, vvar, svar, par_mv, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
modes = ('grad', 'div', 'eval')
name = 'de_sd_stokes'
texpr = 'ij,i.j,0'
class sfepy.terms.terms_sensitivity.ESDVectorDotGradScalarTerm(*args, **kwargs)[source]

Sensitivity of volume dot product of a vector and a gradient of scalar.

Definition:

\int_{\Omega} \hat{I}_{ij} {\partial p \over \partial x_j}\, v_i
\mbox{ , }
\int_{\Omega} \hat{I}_{ij} {\partial q \over \partial x_j}\, u_i

\hat{I}_{ij} = \delta_{ij} \nabla \cdot \Vcal
  - {\partial \Vcal_j \over \partial x_i}

Call signature:

de_sd_v_dot_grad_s

(opt_material, virtual, state, parameter_mv)

(opt_material, state, virtual, parameter_mv)

(opt_material, parameter_v, parameter_s, parameter_mv)

Arguments 1:
  • virtual/parameter_v: \ul{v}

  • state/parameter_s: p

  • parameter_mv: \ul{\Vcal}

Arguments 2:
  • state : \ul{u}

  • virtual : q

  • parameter_mv: \ul{\Vcal}

name = 'de_sd_v_dot_grad_s'
texpr = 'ij,i,0.j'
sfepy.terms.terms_sensitivity.get_nonsym_grad_op(sgrad)[source]