sfepy.terms.terms_dg module

Discontinous Galekrin method specific terms

Note

In einsum calls the following convention is used:

i represents iterating over all cells of a region;

n represents iterating over selected cells of a region, for example over cells on boundary;

b represents iterating over basis functions of state variable;

d represents iterating over basis functions of test variable;

k, l , m represent iterating over geometric dimensions, for example coordinates of velocity or facet normal vector or rows and columns of diffusion tensor;

q represents iterating over quadrature points;

f represents iterating over facets of cell;

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

Lax-Friedrichs flux term for advection of scalar quantity p with the advection velocity \ul{a} given as a material parameter (a known function of space and time).

Definition:

\int_{\partial{T_K}} \ul{n} \cdot \ul{f}^{*} (p_{in}, p_{out})q

where

\ul{f}^{*}(p_{in}, p_{out}) =  \ul{a}  \frac{p_{in} + p_{out}}{2}  +
(1 - \alpha) \ul{n} C \frac{ p_{in} - p_{out}}{2},

\alpha \in [0, 1]; \alpha = 0 for upwind scheme, \alpha = 1 for central scheme, and

C = \max_{p \in [?, ?]}\left\lvert n_x a_1 +
n_y a_2 \right\rvert =
\max_{p \in [?, ?]} \left\lvert  \ul{n} \cdot \ul{a} \right\rvert

the p_{in} resp. p_{out} is solution on the boundary of the element provided by element itself resp. its neighbor and \ul{a} is advection velocity.

Call signature:

dw_dg_advect_laxfrie_flux

(opt_material, material_advelo, virtual, state)

Arguments 1:
  • material : \ul{a}

  • virtual : q

  • state : p

Arguments 3:
  • material : \ul{a}

  • virtual : q

  • state : p

  • opt_material : \alpha

alpha = 0
arg_shapes = [{'material_advelo': 'D, 1', 'opt_material': '.: 1', 'state': 1, 'virtual': (1, 'state')}, {'opt_material': None}]
arg_types = ('opt_material', 'material_advelo', 'virtual', 'state')
function(out, state, diff_var, field, region, advelo)[source]
get_fargs(alpha, advelo, test, state, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
integration = 'cell'
modes = ('weak',)
name = 'dw_dg_advect_laxfrie_flux'
symbolic = {'expression': 'div(a*p)*w', 'map': {'a': 'material', 'p': 'state', 'v': 'virtual'}}
class sfepy.terms.terms_dg.DGTerm(name, arg_str, integral, region, **kwargs)[source]

Abstract base class for DG terms, provides alternative call_function and eval_real methods to accommodate returning iels and vals.

call_function(out, fargs)[source]
eval_real(shape, fargs, mode='eval', term_mode=None, diff_var=None, **kwargs)[source]
poly_space_base = 'legendre'
class sfepy.terms.terms_dg.DiffusionDGFluxTerm(name, arg_str, integral, region, **kwargs)[source]

Basic DG diffusion flux term for scalar quantity.

Definition:

\int_{\partial{T_K}} D \langle \nabla p \rangle [q] \mbox{ , }
\int_{\partial{T_K}} D \langle \nabla q \rangle [p]

where

\langle \nabla \phi \rangle = \frac{\nabla\phi_{in} + \nabla\phi_{out}}{2}

[\phi] = \phi_{in} - \phi_{out}

Math:

The p_{in} resp. p_{out} is solution on the boundary of the element provided by element itself resp. its neighbour.

Call signature:

dw_dg_diffusion_flux

(material, state, virtual)

(material, virtual, state)

Arguments 1:
  • material : D

  • state : p

  • virtual : q

Arguments 2:
  • material : D

  • virtual : q

  • state : p

arg_shapes = [{'material': '1, 1', 'state/avg_state': 1, 'state/avg_virtual': 1, 'virtual/avg_state': (1, None), 'virtual/avg_virtual': (1, None)}]
arg_types = (('material', 'state', 'virtual'), ('material', 'virtual', 'state'))
function(out, state, diff_var, field, region, D)[source]
get_fargs(diff_tensor, test, state, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
integration = 'cell'
modes = ('avg_state', 'avg_virtual')
name = 'dw_dg_diffusion_flux'
class sfepy.terms.terms_dg.DiffusionInteriorPenaltyTerm(name, arg_str, integral, region, **kwargs)[source]

Penalty term used to counteract discontinuity arising when modeling diffusion using Discontinuous Galerkin schemes.

Definition:

\int_{\partial{T_K}} \bar{D} C_w \frac{Ord^2}{d(\partial{T_K})}[p][q]

where

[\phi] = \phi_{in} - \phi_{out}

Math:

the p_{in} resp. p_{out} is solution on the boundary of the element provided by element itself resp. its neighbour.

Call signature:

dw_dg_interior_penalty

(material, material_Cw, virtual, state)

Arguments:
  • material : D

  • material : C_w

  • state : p

  • virtual : q

arg_shapes = [{'material': '1, 1', 'material_Cw': '.: 1', 'state': 1, 'virtual': (1, 'state')}]
arg_types = ('material', 'material_Cw', 'virtual', 'state')
function(out, state, diff_var, field, region, Cw, diff_tensor)[source]
get_fargs(diff_tensor, Cw, test, state, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
modes = ('weak',)
name = 'dw_dg_interior_penalty'
class sfepy.terms.terms_dg.NonlinearHyperbolicDGFluxTerm(name, arg_str, integral, region, **kwargs)[source]

Lax-Friedrichs flux term for nonlinear hyperpolic term of scalar quantity p with the vector function \ul{f} given as a material parameter.

Definition:

\int_{\partial{T_K}} \ul{n} \cdot f^{*} (p_{in}, p_{out})q

where

\ul{f}^{*}(p_{in}, p_{out}) =  \frac{\ul{f}(p_{in})
    + \ul{f}(p_{out})}{2}  +
(1 - \alpha) \ul{n} C \frac{ p_{in} - p_{out}}{2},

\alpha \in [0, 1]; \alpha = 0 for upwind scheme, \alpha = 1 for central scheme, and

C =
\max_{p \in [?, ?]}\left\lvert
n_x \frac{d f_1}{d p} + n_y \frac{d f_2}{d p}
+ \cdots
\right\rvert =

\max_{p \in [?, ?]} \left\lvert
\vec{n}\cdot\frac{d\ul{f}}{dp}(p)
 \right\rvert

the p_{in} resp. p_{out} is solution on the boundary of the element provided by element itself resp. its neighbor.

Call signature:

dw_dg_nonlinear_laxfrie_flux

(opt_material, fun, fun_d, virtual, state)

Arguments 1:
  • material : \ul{f}

  • material : \frac{d\ul{f}}{d p}

  • virtual : q

  • state : p

Arguments 3:
  • material : \ul{f}

  • material : \frac{d\ul{f}}{d p}

  • virtual : q

  • state : p

  • opt_material : \alpha

alf = 0
arg_shapes = [{'material_fun': '.: 1', 'material_fun_d': '.: 1', 'opt_material': '.: 1', 'state': 1, 'virtual': (1, 'state')}, {'opt_material': None}]
arg_types = ('opt_material', 'fun', 'fun_d', 'virtual', 'state')
function(out, state, field, region, f, df)[source]
get_fargs(alpha, fun, dfun, test, state, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
integration = 'cell'
modes = ('weak',)
name = 'dw_dg_nonlinear_laxfrie_flux'
symbolic = {'expression': 'div(f(p))*w', 'map': {'f': 'function', 'p': 'state', 'v': 'virtual'}}
class sfepy.terms.terms_dg.NonlinearScalarDotGradTerm(name, arg_str, integral, region, **kwargs)[source]

Product of virtual and divergence of vector function of state or volume dot product of vector function of state and gradient of scalar virtual.

Definition:

\int_{\Omega} q \cdot \nabla \cdot \ul{f}(p) = \int_{\Omega} q \cdot
   \text{div} \ul{f}(p)  \mbox{ , }
\int_{\Omega} \ul{f}(p) \cdot \nabla q

Call signature:

dw_ns_dot_grad_s

(fun, fun_d, virtual, state)

(fun, fun_d, state, virtual)

Arguments 1:
  • function : \ul{f}

  • virtual : q

  • state : p

Arguments 2:
  • function : \ul{f}

  • state : p

  • virtual : q

TODO maybe this term would fit better to terms_dot?

arg_shapes = [{'material_fun': '.: 1', 'material_fun_d': '.: 1', 'state/grad_state': 1, 'state/grad_virtual': 1, 'virtual/grad_state': (1, None), 'virtual/grad_virtual': (1, None)}]
arg_types = (('fun', 'fun_d', 'virtual', 'state'), ('fun', 'fun_d', 'state', 'virtual'))
static function(out, out_qp, geo, fmode)[source]
get_fargs(fun, dfun, var1, var2, mode=None, term_mode=None, diff_var=None, **kwargs)[source]
modes = ('grad_state', 'grad_virtual')
name = 'dw_ns_dot_grad_s'