Eigenfunctions

This modules provides eigenfunctions for a certain set of parabolic problems. Therefore functions for the computation of the corresponding eigenvalues are included. The functions which compute the eigenvalues are deliberately separated from the predefined eigenfunctions in order to handle transformations and reduce effort by the controller implementation.

class AddMulFunction(function)[source]

Bases: object

(Temporary) Function class wich can multiplied with scalars and added with functions. Only needed to compute the matrix (of scalars) vector (of functions) product in FiniteTransformFunction. Will be no longer needed when pyinduct.core.Function is overloaded with __add__ and __mul__ operator.

Parameters:function (callable) –
class FiniteTransformFunction(function, M, l, scale_func=None, nested_lambda=False)[source]

Bases: pyinduct.core.Function

Provide a transformed pyinduct.core.Function \bar x(z) through the transformation \bar{\boldsymbol{\xi}} = T * \boldsymbol \xi, with the function vector \boldsymbol \xi\in\mathbb R^{2n} and with a given matrix T\in\mathbb R^{2n\times 2n}. The operator * denotes the matrix (of scalars) vector (of functions) product. The interim result \bar{\boldsymbol{\xi}} is a vector

\bar{\boldsymbol{\xi}} = (\bar\xi_{1,0},...,\bar\xi_{1,n-1},\bar\xi_{2,0},...,\bar\xi_{2,n-1})^T.

of functions

&\bar\xi_{1,j} = \bar x(jl_0 + z),\qquad j=0,...,n-1, \quad l_0=l/n, \quad z\in[0,l_0] \\
&\bar\xi_{2,j} = \bar x(l - jl_0 + z).

Finally, the provided function \bar x(z) is given through \bar\xi_{1,0},...,\bar\xi_{1,n-1}.

Note

For a more extensive documentation see section 4.2 in:

Parameters:
  • function (callable) –

    Function x(z) which will subdivided in 2n Functions

    &\bar\xi_{1,j} = x(jl_0 + z),\qquad j=0,...,n-1, \quad l_0=l/n, \quad z\in[0,l_0] \\
&\bar\xi_{2,j} = x(l - jl_0 + z).

    The vector of functions \boldsymbol\xi consist of these functions:

    \boldsymbol\xi = (\xi_{1,0},...,\xi_{1,n-1},\xi_{2,0},...,\xi_{2,n-1})^T) .

  • M (numpy.ndarray) – Matrix T\in\mathbb R^{2n\times 2n} of scalars.
  • l (numbers.Number) – Length of the domain (z\in[0,l]).
class SecondOrderDirichletEigenfunction(omega, param, spatial_domain, norm_fac=1.0)[source]

Bases: pyinduct.core.Function

Provide the eigenfunction \varphi(z) to an eigenvalue problem of the form

a_2\varphi''(z) + a_1&\varphi'(z) + a_0\varphi(z) = \lambda\varphi(z) \\
\varphi(0) &= 0 \\
\varphi(l) &= 0.

The eigenfrequency

\omega = \sqrt{-\frac{a_1^2}{4a_2^2}+\frac{a_0-\lambda}{a_2}}

must be provided.

Parameters:
  • om (numbers.Number) – eigenfrequency \omega
  • param (array_like) – \Big( a_2, a_1, a_0, None, None \Big)^T
  • spatial_domain (tuple) – Start point z_0 and end point z_1 of the spatial domain [z_0,z_1]\ni z.
  • norm_fac (numbers.Number) – Factor to scale the eigenfunctions.
class SecondOrderRobinEigenfunction(om, param, spatial_domain, phi_0=1)[source]

Bases: pyinduct.core.Function

Provide the eigenfunction \varphi(z) to an eigenvalue problem of the form

a_2\varphi''(z) + a_1&\varphi'(z) + a_0\varphi(z) = \lambda\varphi(z) \\
\varphi'(0) &= \alpha \varphi(0) \\
\varphi'(l) &= -\beta \varphi(l).

The eigenfrequency

\omega = \sqrt{-\frac{a_1^2}{4a_2^2}+\frac{a_0-\lambda}{a_2}}

must be provided (with compute_rad_robin_eigenfrequencies).

Parameters:
  • om (numbers.Number) – eigenfrequency \omega
  • param (array_like) – \Big( a_2, a_1, a_0, \alpha, \beta \Big)^T
  • spatial_domain (tuple) – Start point z_0 and end point z_1 of the spatial domain [z_0,z_1]\ni z.
  • phi_0 (numbers.Number) – Factor to scale the eigenfunctions (correspond \varphi(0)=\text{phi\_0}).
class TransformedSecondOrderEigenfunction(target_eigenvalue, init_state_vect, dgl_coefficients, domain)[source]

Bases: pyinduct.core.Function

Provide the eigenfunction \varphi(z) to an eigenvalue problem of the form

a_2(z)\varphi''(z) + a_1(z)\varphi'(z) + a_0(z)\varphi(z) = \lambda\varphi(z)

where \lambda is a predefined (potentially complex) eigenvalue and [z_0,z_1]\ni z is the domain.

Parameters:
  • target_eigenvalue (numbers.Number) – \lambda
  • init_state_vect (array_like) –

    \Big(\text{Re}\{\varphi(0)\}, \text{Re}\{\varphi'(0)\}, \text{Im}\{\varphi(0)\}, \text{Im}\{\varphi'(0)\}\Big)^T

  • dgl_coefficients (array_like) – \Big( a2(z), a1(z), a0(z) \Big)^T
  • domain (array_like) – \Big( z_0, ..... , z_1 \Big)
compute_rad_robin_eigenfrequencies(param, l, n_roots=10, show_plot=False)[source]

Return the first n_roots eigenfrequencies \omega (and eigenvalues \lambda)

\omega = \sqrt{-\frac{a_1^2}{4a_2^2}+\frac{a_0-\lambda}{a_2}}

to the eigenvalue problem

a_2\varphi''(z) + a_1&\varphi'(z) + a_0\varphi(z) = \lambda\varphi(z) \\
\varphi'(0) &= \alpha\varphi(0) \\
\varphi'(l) &= -\beta\varphi(l).

Parameters:
  • param (array_like) – \Big( a_2, a_1, a_0, \alpha, \beta \Big)^T
  • l (numbers.Number) – Right boundary value of the domain [0,l]\ni z.
  • n_roots (int) – Amount of eigenfrequencies to be compute.
  • show_plot (bool) – A plot window of the characteristic equation appears if it is True.
Returns:

\Big(\big[\omega_1,...,\omega_\text{n\_roots}\Big], \Big[\lambda_1,...,\lambda_\text{n\_roots}\big]\Big)

Return type:

tuple –> booth tuple elements are numpy.ndarrays of length nroots

get_adjoint_rad_evp_param(param)[source]

Return to the eigen value problem of the reaction-advection-diffusion equation with robin and/or dirichlet boundary conditions

a_2\varphi''(z) + a_1&\varphi'(z) + a_0\varphi(z) = \lambda\varphi(z) \\
\varphi(0) = 0 \quad &\text{or} \quad \varphi'(0) = \alpha\varphi(0) \\
\varphi`(l) = 0 \quad &\text{or} \quad \varphi'(l) = -\beta\varphi(l)

the parameters for the adjoint problem (with the same structure).

Parameters:param (array_like) – \Big( a_2, a_1, a_0, \alpha, \beta \Big)^T
Returns:Parameters \big(a_2, \tilde a_1=-a_1, a_0, \tilde \alpha, \tilde \beta \big) for the adjoint problem

a_2\psi''(z) + a_1&\psi'(z) + a_0\psi(z) = \lambda\psi(z) \\
\psi(0) = 0 \quad &\text{or} \quad \psi'(0) = \tilde\alpha\psi(0) \\
\psi`(l) = 0 \quad &\text{or} \quad \psi'(l) = -\tilde\beta\psi(l).

Return type:tuple
return_real_part(to_return)[source]

Check if the imaginary part of to_return vanishes and return the real part.

Parameters:to_return (numbers.Number or array_like) – Variable to check.
Raises:ValueError – If (all) imaginary part(s) not vanishes.
Returns:Real part of to_return.
Return type:numbers.Number or array_like
transform2intermediate(param, d_end=None)[source]

Transformation \tilde x(z,t)=x(z,t)e^{\int_0^z \frac{a_1(\bar z)}{2 a_2}\,d\bar z} which eliminate the advection term a_1 x(z,t) from the reaction-advection-diffusion equation

\dot x(z,t) = a_2 x''(z,t) + a_1(z) x'(z,t) + a_0(z) x(z,t)

with robin

x'(0,t) = \alpha x(0,t), \quad x'(l,t) = -\beta x(l,t)

or dirichlet

x(0,t) = 0, \quad x(l,t) = 0

or mixed boundary condition.

Parameters:param (array_like) – \Big( a_2, a_1, a_0, \alpha, \beta \Big)^T
Raises:TypeError – If a_1(z) is callable but no derivative handle is defined for it.
Returns:Parameters \big(a_2, \tilde a_1=0, \tilde a_0(z), \tilde \alpha, \tilde \beta \big) for the transformed system

\dot{\tilde{x}}(z,t) = a_2 \tilde x''(z,t) + \tilde a_0(z) \tilde x(z,t)

and the corresponding boundary conditions (\alpha and/or \beta set to None by dirichlet boundary condition).

Return type:tuple