elastodynamicsx.pde

Module contents

The pde module contains tools for building a Partial Differential Equation from material and damping laws, as well as boundary conditions.

The PDE is of the form:
M*a + C*v + K(u) = b,
+ Boundary conditions.

The module also contains tools to formulate a time domain problem using implicit or explicit time schemes, although the preferred way to use these tools is through the elastodynamicsx.solvers.timestepper function.

— — —

Material laws:
from elastodynamicsx.pde import material
mat = material( (function_space, cell_tags, marker), type_, *args, **kwargs)
Body forces:
from elastodynamicsx.pde import BodyForce
bf  = BodyForce( (function_space, cell_tags, marker), value)
Boundary conditions:
from elastodynamicsx.pde import BoundaryCondition
bc  = BoundaryCondition( (function_space, facet_tags, marker), type_, value)
Assembling a PDE:
from elastodynamicsx.pde import PDE
pde = PDE(materials=[mat1, mat2, ...], bodyforces=[bf1, bf2, ...], bcs=[bc1, bc2, ...])
Getting the forms:
# Python functions
M_function = pde.M_fn
# Use as:
# u = ufl.TrialFunction(function_space)
# v = ufl.TestFunction(function_space)
# M_ufl = M_function(u, v)

# Compiled dolfinx forms
M_form = pde.M_form

# PETSc Matrices
M = pde.M()
class elastodynamicsx.pde.BodyForce(functionspace_tags_marker, value: Function | Constant | ndarray, **kwargs)[source]

Bases: object

Representation of the rhs term (the ‘b’ term) of a pde such as defined in the PDE class. An instance represents a single source.

Parameters:
  • functionspace_tags_marker – The function space, cell tags and marker(s) where the linear form is defined

  • value – The value of the body force as a ufl-compatible object

Keyword Arguments:

metadata – dict, compiler-specific parameters, passed to ufl.Measure(…, metadata=metadata)

Example

# V is a vector function space, dim=2

# Gaussian source distribution
X0_src = np.array([length/2, height/2, 0])  # center
R0_src = 0.1  # radius
F0_src = default_scalar_type([1, 0])  # amplitude of the source

def src_x(x):  # source(x): Gaussian
    nrm = 1 / (2 * np.pi * R0_src**2)  # normalize to int[src_x(x) dx]=1
    r = np.linalg.norm(x - X0_src[:, np.newaxis], axis=0)
    return nrm * np.exp(-1/2 * (r / R0_src)**2, dtype=default_scalar_type)

value = fem.Function(V)
value.interpolate(lambda x: F0_src[:, np.newaxis] * src_x(x)[np.newaxis,:])

bf = BodyForce(V, value)  # defined in the entire domain
bf = BodyForce((V, cell_tags, (1, 4)), value)  # definition restricted to tags 1 and 4
b_fn(v)[source]

The linear form function

property value: Function | Constant | ndarray
elastodynamicsx.pde.boundarycondition(functionspace_tags_marker, type_: str, *args, **kwargs) BoundaryConditionBase[source]

Builder function that instanciates the desired boundary condition (BC) type.

Possible BCs are:

‘Free’, ‘Clamp’ ‘Dirichlet’, ‘Neumann’, ‘Robin’ ‘Periodic’ ‘Dashpot’

Parameters:
  • functionspace_tags_marker – The function space, facet tags and marker(s) where the BC is applied

  • type – String identifier of the BC

  • *args – Passed to the required BC

Keyword Arguments:

**kwargs – Passed to the required BC

Example

General:

# Imposes u = u_D on boundary n°1
u_D = fem.Constant(V.mesh, [0.2,1.3])
bc = boundarycondition((V, facet_tags, 1), 'Dirichlet', u_D)

# Imposes sigma(u).n = T_N on boundary n°1
T_N = fem.Function(V)
T_N.interpolate(fancy_function)
bc = boundarycondition((V, facet_tags, 1), 'Neumann', T_N)

# Imposes sigma(u).n = z * (u - s) on boundary n°1
z = ufl.as_matrix([[1, 2], [3, 4]])
s = ufl.as_vector([0, 1])
bc = boundarycondition((V, facet_tags, 1), 'Robin', z, s)

Free or Clamp, on one or several tags:

# Imposes sigma(u).n=0 on boundary n°1
bc = boundarycondition((V, facet_tags, 1), 'Free')

# Imposes u=0 on boundaries n°1 and 2
bc = boundarycondition((V, facet_tags, (1, 2)), 'Clamp')

# Apply BC to all boundaries when facet_tags and marker are not specified
# or set to None
bc = boundarycondition(V, 'Clamp')

Periodic:

# Given:   x(2) = x(1) - P
# Imposes: u(2) = u(1)
# Where: x(i) are the coordinates on boundaries n°1, 2
#        P is a constant (translation) vector from slave to master
#        u(i) is the field on boundaries n°1,2
# Note that boundary n°2 is slave

Px, Py, Pz = length, 0, 0 #for x-periodic, and tag=left
Px, Py, Pz = 0, height, 0 #for y-periodic, and tag=bottom
Px, Py, Pz = 0, 0, depth  #for z-periodic, and tag=back
P  = [Px,Py,Pz]
bc = boundarycondition((V, facet_tags, 2), 'Periodic', P)

BCs involving the velocity:

# (for a scalar function_space)
# Imposes sigma(u).n = z * v on boundary n°1, with v=du/dt. Usually z=rho*c
bc = boundarycondition((V, facet_tags, 1), 'Dashpot', z)

# (for a vector function_space)
# Imposes sigma(u).n = z_N * v_N + z_T * v_T on boundary n°1,
# with v = du/dt, v_N = (v.n) n, v_T = v - v_N
# Usually z_N = rho * c_L and z_T = rho * c_S
bc = boundarycondition((V, facet_tags, 1), 'Dashpot', z_N, z_T)
Adapted from:

https://jsdokken.com/dolfinx-tutorial/chapter3/robin_neumann_dirichlet.html

elastodynamicsx.pde.material(functionspace_tags_marker, type_, *args, **kwargs) Material[source]

Builder method that instanciates the desired material type

Parameters:
  • functionspace_tags_marker

    Available possibilities are:

    • (function_space, cell_tags, marker): meaning application

      of the material in the cells whose tag correspond to marker

    • function_space: in this case cell_tags=None

      and marker=None, meaning application of the material in the entire domain

  • type

    Available options are:

    • linear (scalar): ‘scalar’

    • linear: ‘isotropic’, ‘cubic’, ‘hexagonal’, ‘trigonal’, ‘tetragonal’,

      ’orthotropic’, ‘monoclinic’, ‘triclinic’

    • nonlinear, hyperelastic: ‘murnaghan’, ‘saintvenant-kirchhoff’, ‘mooney-rivlin-incomp’

  • *args – Passed to the required material

Keyword Arguments:

**kwargs – Passed to the required material

Returns:

An instance of the desired material

Example

# ###
# Constant / variable material parameter
# ###
rho = 2.8
rho = fem.Constant(function_space.mesh, default_scalar_type(2.8))
rho = fem.Function(scalar_function_space); rho.interpolate(lambda x: 2.8 * np.ones(len(x)))

# ###
# Subdomain(s) or entire domain?
# ###

# restricted to subdomain number 1
aluminum = material((function_space, cell_tags, 1), 'isotropic',
                     rho=2.8, lambda_=58, mu=26)

# restricted to subdomains number 1, 2, 5
aluminum = material((function_space, cell_tags, (1,2,5)), 'isotropic',
                     rho=2.8, lambda_=58, mu=26)

# entire domain
aluminum = material( function_space, 'isotropic',
                     rho=2.8, lambda_=58, mu=26)

# ###
# Available laws
# ###

# Linear elasticity
mat = material( function_space, 'isotropic'  , rho, C12, C44)
mat = material( function_space, 'cubic'      , rho, C11, C12, C44)
mat = material( function_space, 'hexagonal'  , rho, C11, C13, C33, C44, C66)
mat = material( function_space, 'trigonal'   , rho, C11, C12, C13, C14, C25,
                                                    C33, C44)
mat = material( function_space, 'tetragonal' , rho, C11, C12, C13, C16, C33, C44, C66)
mat = material( function_space, 'orthotropic', rho, C11, C12, C13, C22, C23, C33, C44,
                                                    C55, C66)
mat = material( function_space, 'monoclinic' , rho, C11, C12, C13, C15, C22, C23, C25,
                                                    C33, C35, C44, C46, C55, C66)
mat = material( function_space, 'triclinic'  , rho, C11, C12, C13, C14, C15, C16, C22,
                                                    C23, C24, C25, C26, C33, C34, C35,
                                                    C36, C44, C45, C46, C55, C56, C66)

# Hyperelasticity
mat = material( function_space, 'saintvenant-kirchhoff', rho, C12, C44)
mat = material( function_space, 'murnaghan', rho, C12, C44, l, m, n)
mat = material( function_space, 'mooney-rivlin-incomp', rho, C1, C2)

# Arbitrary
mat = material( function_space,
                'custom',
                M_fn=lambda u,v: ufl.inner(u,v) * ufl.dx,
                K_fn=lambda u,v: ... * ufl.dx)
elastodynamicsx.pde.damping(type_: str, *args) Damping[source]

Builder method that instanciates the desired damping law type

Parameters:
  • type

    Available options are:

    • ’none’

    • ’rayleigh’

  • *args – Passed to the required damping law

Returns:

An instance of the desired damping law

Example

from elastodynamicsx.pde import damping, material

dmp = damping('none')
dmp = damping('rayleigh', eta_m, eta_k)

mat = material(V, 'isotropic', rho, lambda_, mu, damping=dmp)
class elastodynamicsx.pde.PDE(function_space: FunctionSpace, materials: Iterable[Material], **kwargs)[source]

Bases: object

Representation of a PDE of the kind:

M*a + C*v + K(u) = b + Boundary conditions

as an assembly of materials, forces and bcs defined over different subdomains.

Parameters:
  • functions_space

  • materials – a list of pde.Material instances

Keyword Arguments:
  • bodyforces – (default=[]) a list of pde.BodyForce instances

  • bcs – (default=[]) a list of fem.DirichletBCMetaClass and/or pde.BoundaryConditionBase instances

  • jit_options – (default=PDECONFIG.default_jit_options) options for the just-in-time compiler

  • finalize – (default=True) call self.finalize() on build

property is_linear: bool
property mpc
finalize() None[source]

Finalize dolfinx-mpc objects

M_fn(u, v)[source]

(bilinear) Mass form function

C_fn(u, v)[source]

(bilinear) Damping form function

K_fn(u, v)[source]

(bilinear) Stiffness form function

K0_fn(u, v)[source]

(bilinear) K0 stiffness form function (waveguide problems)

K1_fn(u, v)[source]

(bilinear) K1 stiffness form function (waveguide problems)

K2_fn(u, v)[source]

(bilinear) K2 stiffness form function (waveguide problems)

K_fn_CG(u, v)[source]

(bilinear) Stiffness form function (Continuous Galerkin)

K_fn_DG(u, v)[source]

(bilinear) Stiffness form function (Discontinuous Galerkin)

DG_numerical_flux(u, v)[source]

(bilinear) Numerical flux form function (Disontinuous Galerkin)

b_fn(v)[source]

Linear form function

property M_form: Form | None

Compiled mass bilinear form

property C_form: Form | None

Compiled damping bilinear form

property K_form: Form | None

Compiled stiffness bilinear form

property b_form: Form | None

Compiled linear form

M() Mat[source]

Mass matrix

C() Mat[source]

Damping matrix

K() Mat[source]

Stiffness matrix

b(omega=0) Vec[source]

Load vector

init_b() Vec[source]

Declares a zero vector compatible with the linear form

K0() Mat[source]

K0 stiffness matrix (waveguide problems)

K1() Mat[source]

K1 stiffness matrix (waveguide problems)

K2() Mat[source]

K2 stiffness matrix (waveguide problems)

class elastodynamicsx.pde.PDECONFIG[source]

Bases: object

Global configuration parameters for classes in the pde module

default_metadata: Dict | None = None

The default metadata used by all measures (dx, ds, dS, …) in the classes of the pde package: Material, BoundaryCondition, BodyForce

Example

Spectral Element Method with GLL elements of degree 6:

from elastodynamicsx.pde import PDECONFIG
from elastodynamicsx.utils import spectral_quadrature
specmd = spectral_quadrature("GLL", 6)
CONFIG.default_metadata = specmd
default_jit_options: Dict = {}

The default options for the just-in-time compiler used in the classes of the pde package: PDE, FEniCSxTimeScheme

See:

https://jsdokken.com/dolfinx-tutorial/chapter4/compiler_parameters.html

Example

from elastodynamicsx.pde import PDECONFIG
PDECONFIG.default_jit_options = {"cffi_extra_compile_args": ["-Ofast", "-march=native"],
                                 "cffi_libraries": ["m"]}

elastodynamicsx.pde.bodyforce

class elastodynamicsx.pde.bodyforce.BodyForce(functionspace_tags_marker, value: Function | Constant | ndarray, **kwargs)[source]

Bases: object

Representation of the rhs term (the ‘b’ term) of a pde such as defined in the PDE class. An instance represents a single source.

Parameters:
  • functionspace_tags_marker – The function space, cell tags and marker(s) where the linear form is defined

  • value – The value of the body force as a ufl-compatible object

Keyword Arguments:

metadata – dict, compiler-specific parameters, passed to ufl.Measure(…, metadata=metadata)

Example

# V is a vector function space, dim=2

# Gaussian source distribution
X0_src = np.array([length/2, height/2, 0])  # center
R0_src = 0.1  # radius
F0_src = default_scalar_type([1, 0])  # amplitude of the source

def src_x(x):  # source(x): Gaussian
    nrm = 1 / (2 * np.pi * R0_src**2)  # normalize to int[src_x(x) dx]=1
    r = np.linalg.norm(x - X0_src[:, np.newaxis], axis=0)
    return nrm * np.exp(-1/2 * (r / R0_src)**2, dtype=default_scalar_type)

value = fem.Function(V)
value.interpolate(lambda x: F0_src[:, np.newaxis] * src_x(x)[np.newaxis,:])

bf = BodyForce(V, value)  # defined in the entire domain
bf = BodyForce((V, cell_tags, (1, 4)), value)  # definition restricted to tags 1 and 4
b_fn(v)[source]

The linear form function

property value: Function | Constant | ndarray

elastodynamicsx.pde.boundaryconditions

Module for handling boundary conditions (BCs) of various types. The proper way to instantiate a BC is using the boundarycondition function.

class elastodynamicsx.pde.boundaryconditions.BoundaryConditionBase(functionspace_tags_marker, **kwargs)[source]

Bases: object

Base class for representating of a variety of boundary conditions (BCs). The proper way to instantiate a BC is using the boundarycondition function.

BCs split into three families:
  • BCStrongBase class: strong BCs -> map to fem.DirichletBC

  • BCMPCBase class: multi-point constraints BCs -> used in dolfinx_mpc

  • BCWeakBase class: weak BCs -> terms added to the weak forms

labels: List[str]

String identifier(s) of the class, recognized by the boundarycondition function

property bc

The boundary condition

property function_space
property facet_tags
property marker
class elastodynamicsx.pde.boundaryconditions.BCStrongBase(functionspace_tags_marker, **kwargs)[source]

Bases: BoundaryConditionBase

Base class for representating strong BCs. This class carries a fem.DirichletBC object.

class elastodynamicsx.pde.boundaryconditions.BCMPCBase(functionspace_tags_marker, **kwargs)[source]

Bases: BoundaryConditionBase

Base class for representating multi-point constraints BCs. Requires dolfinx_mpc to be installed.

class elastodynamicsx.pde.boundaryconditions.BCWeakBase(functionspace_tags_marker, **kwargs)[source]

Bases: BoundaryConditionBase

Base class for representating weak BCs.

property ds
C_fn(u, v)[source]
K_fn(u, v)[source]
b_fn(v)[source]
class elastodynamicsx.pde.boundaryconditions.BCDirichlet(functionspace_tags_marker, value: Function | Constant | ndarray, **kwargs)[source]

Bases: BCStrongBase

Representation of a Dirichlet BC. Set \(u = u_0\), with \(u_0\) the prescribed displacement.

labels: List[str] = ['dirichlet']

String identifier(s) of the class, recognized by the boundarycondition function

property value

Value of the prescribed displacement \(u_0\)

class elastodynamicsx.pde.boundaryconditions.BCClamp(functionspace_tags_marker, **kwargs)[source]

Bases: BCDirichlet

Same as a Dirichlet BC with value set to zero -> Set \(u = 0\).

labels: List[str] = ['clamp']

String identifier(s) of the class, recognized by the boundarycondition function

class elastodynamicsx.pde.boundaryconditions.BCNeumann(functionspace_tags_marker, value: Function | Constant | ndarray, **kwargs)[source]

Bases: BCWeakBase

Representation of a Neumann BC. Prescribes \(\sigma(u).n = t_n\), with \(t_n\) a boundary traction.

labels: List[str] = ['neumann']

String identifier(s) of the class, recognized by the boundarycondition function

b_fn(v)[source]
property value

Value of the \(t_n\) boundary traction

class elastodynamicsx.pde.boundaryconditions.BCFree(functionspace_tags_marker, **kwargs)[source]

Bases: BCNeumann

Same as a Neumann BC with value set to zero. This BC has actually no effect on the PDE.

labels: List[str] = ['free']

String identifier(s) of the class, recognized by the boundarycondition function

b_fn(v)[source]
class elastodynamicsx.pde.boundaryconditions.BCRobin(functionspace_tags_marker, value1: Function | Constant | ndarray, value2: Function | Constant | ndarray, **kwargs)[source]

Bases: BCWeakBase

Representation of a Robin BC. Prescribes \(\sigma(u).n = z * (u - s)\), with \(z\) an impedance matrix and \(s\) a given vector.

labels: List[str] = ['robin']

String identifier(s) of the class, recognized by the boundarycondition function

K_fn(u, v)[source]
b_fn(v)[source]
property value1
property value2
class elastodynamicsx.pde.boundaryconditions.BCDashpot(functionspace_tags_marker, *values: Function | Constant | ndarray, **kwargs)[source]

Bases: BCWeakBase

Representation of a Dashpot BC, i.e. an impedance condition involving the velocity. Prescribes \(\sigma(u).n = z_n * v_n + z_t * v_t\), with \(z_n\) and \(z_t\) the normal and transverse impedances, and with \(v_n\) and \(v_t\) the normal and transverse components of the velocity.

labels: List[str] = ['dashpot']

String identifier(s) of the class, recognized by the boundarycondition function

C_fn(u, v)[source]
property values
class elastodynamicsx.pde.boundaryconditions.BCCustom(functionspace_tags_marker, **kwargs)[source]

Bases: BCWeakBase

Representation of a boundary condition with user-defined damping (C), stiffness (K) and linear (b) forms

Keyword Arguments:
  • C_fn – (optional) callable returning the ufl (bilinear) damping form

  • K_fn – (optional) callable returning the ufl (bilinear) stiffness form

  • b_fn – (optional) callable returning the ufl linear form

labels: List[str] = ['custom']

String identifier(s) of the class, recognized by the boundarycondition function

class elastodynamicsx.pde.boundaryconditions.BCPeriodic(functionspace_tags_marker, *values: Function | Constant | ndarray, **kwargs)[source]

Bases: BCMPCBase

Representation of a Periodic BC (translation), such as \(u(x_2) = u(x_1)\), with \(x_2 = x_1 + P\) and \(P\) a translation vector from the slave to the master boundary

labels: List[str] = ['periodic']

String identifier(s) of the class, recognized by the boundarycondition function

elastodynamicsx.pde.boundaryconditions.get_dirichlet_BCs(bcs: Iterable[BoundaryConditionBase | DirichletBC] | Tuple[()]) List[DirichletBC][source]

Returns the BCs of Dirichlet type in bcs

elastodynamicsx.pde.boundaryconditions.get_mpc_BCs(bcs: Iterable[BoundaryConditionBase | DirichletBC] | Tuple[()]) List[BCMPCBase][source]

Returns the BCs to be built with dolfinx_mpc in bcs

elastodynamicsx.pde.boundaryconditions.get_weak_BCs(bcs: Iterable[BoundaryConditionBase | DirichletBC] | Tuple[()]) List[BCWeakBase][source]

Returns the weak BCs in bcs

elastodynamicsx.pde.boundaryconditions.boundarycondition(functionspace_tags_marker, type_: str, *args, **kwargs) BoundaryConditionBase[source]

Builder function that instanciates the desired boundary condition (BC) type.

Possible BCs are:

‘Free’, ‘Clamp’ ‘Dirichlet’, ‘Neumann’, ‘Robin’ ‘Periodic’ ‘Dashpot’

Parameters:
  • functionspace_tags_marker – The function space, facet tags and marker(s) where the BC is applied

  • type – String identifier of the BC

  • *args – Passed to the required BC

Keyword Arguments:

**kwargs – Passed to the required BC

Example

General:

# Imposes u = u_D on boundary n°1
u_D = fem.Constant(V.mesh, [0.2,1.3])
bc = boundarycondition((V, facet_tags, 1), 'Dirichlet', u_D)

# Imposes sigma(u).n = T_N on boundary n°1
T_N = fem.Function(V)
T_N.interpolate(fancy_function)
bc = boundarycondition((V, facet_tags, 1), 'Neumann', T_N)

# Imposes sigma(u).n = z * (u - s) on boundary n°1
z = ufl.as_matrix([[1, 2], [3, 4]])
s = ufl.as_vector([0, 1])
bc = boundarycondition((V, facet_tags, 1), 'Robin', z, s)

Free or Clamp, on one or several tags:

# Imposes sigma(u).n=0 on boundary n°1
bc = boundarycondition((V, facet_tags, 1), 'Free')

# Imposes u=0 on boundaries n°1 and 2
bc = boundarycondition((V, facet_tags, (1, 2)), 'Clamp')

# Apply BC to all boundaries when facet_tags and marker are not specified
# or set to None
bc = boundarycondition(V, 'Clamp')

Periodic:

# Given:   x(2) = x(1) - P
# Imposes: u(2) = u(1)
# Where: x(i) are the coordinates on boundaries n°1, 2
#        P is a constant (translation) vector from slave to master
#        u(i) is the field on boundaries n°1,2
# Note that boundary n°2 is slave

Px, Py, Pz = length, 0, 0 #for x-periodic, and tag=left
Px, Py, Pz = 0, height, 0 #for y-periodic, and tag=bottom
Px, Py, Pz = 0, 0, depth  #for z-periodic, and tag=back
P  = [Px,Py,Pz]
bc = boundarycondition((V, facet_tags, 2), 'Periodic', P)

BCs involving the velocity:

# (for a scalar function_space)
# Imposes sigma(u).n = z * v on boundary n°1, with v=du/dt. Usually z=rho*c
bc = boundarycondition((V, facet_tags, 1), 'Dashpot', z)

# (for a vector function_space)
# Imposes sigma(u).n = z_N * v_N + z_T * v_T on boundary n°1,
# with v = du/dt, v_N = (v.n) n, v_T = v - v_N
# Usually z_N = rho * c_L and z_T = rho * c_S
bc = boundarycondition((V, facet_tags, 1), 'Dashpot', z_N, z_T)
Adapted from:

https://jsdokken.com/dolfinx-tutorial/chapter3/robin_neumann_dirichlet.html

elastodynamicsx.pde.materials

elastodynamicsx.pde.pde

class elastodynamicsx.pde.pde.PDE(function_space: FunctionSpace, materials: Iterable[Material], **kwargs)[source]

Bases: object

Representation of a PDE of the kind:

M*a + C*v + K(u) = b + Boundary conditions

as an assembly of materials, forces and bcs defined over different subdomains.

Parameters:
  • functions_space

  • materials – a list of pde.Material instances

Keyword Arguments:
  • bodyforces – (default=[]) a list of pde.BodyForce instances

  • bcs – (default=[]) a list of fem.DirichletBCMetaClass and/or pde.BoundaryConditionBase instances

  • jit_options – (default=PDECONFIG.default_jit_options) options for the just-in-time compiler

  • finalize – (default=True) call self.finalize() on build

property is_linear: bool
property mpc
finalize() None[source]

Finalize dolfinx-mpc objects

M_fn(u, v)[source]

(bilinear) Mass form function

C_fn(u, v)[source]

(bilinear) Damping form function

K_fn(u, v)[source]

(bilinear) Stiffness form function

K0_fn(u, v)[source]

(bilinear) K0 stiffness form function (waveguide problems)

K1_fn(u, v)[source]

(bilinear) K1 stiffness form function (waveguide problems)

K2_fn(u, v)[source]

(bilinear) K2 stiffness form function (waveguide problems)

K_fn_CG(u, v)[source]

(bilinear) Stiffness form function (Continuous Galerkin)

K_fn_DG(u, v)[source]

(bilinear) Stiffness form function (Discontinuous Galerkin)

DG_numerical_flux(u, v)[source]

(bilinear) Numerical flux form function (Disontinuous Galerkin)

b_fn(v)[source]

Linear form function

property M_form: Form | None

Compiled mass bilinear form

property C_form: Form | None

Compiled damping bilinear form

property K_form: Form | None

Compiled stiffness bilinear form

property b_form: Form | None

Compiled linear form

M() Mat[source]

Mass matrix

C() Mat[source]

Damping matrix

K() Mat[source]

Stiffness matrix

b(omega=0) Vec[source]

Load vector

init_b() Vec[source]

Declares a zero vector compatible with the linear form

K0() Mat[source]

K0 stiffness matrix (waveguide problems)

K1() Mat[source]

K1 stiffness matrix (waveguide problems)

K2() Mat[source]

K2 stiffness matrix (waveguide problems)

elastodynamicsx.pde.timeschemes