Eigenfunctions

This modules provides eigenfunctions for a certain set of second order spatial operators. 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 within the controller implementation.

class SecondOrderOperator(a2=0, a1=0, a0=0, alpha1=0, alpha0=0, beta1=0, beta0=0, domain=(-inf, inf))[source]

Bases: object

Interface class to collect all important parameters that describe a second order ordinary differential equation.

Parameters:
  • a2 (Number or callable) – coefficient a_2.
  • a1 (Number or callable) – coefficient a_1.
  • a0 (Number or callable) – coefficient a_0.
  • alpha1 (Number) – coefficient \alpha_1.
  • alpha0 (Number) – coefficient \alpha_0.
  • beta1 (Number) – coefficient \beta_1.
  • beta0 (Number) – coefficient \beta_0.
static from_dict(param_dict, domain=None)[source]
static from_list(param_list, domain=None)[source]
get_adjoint_problem()[source]

Return the parameters of the operator A^* describing the the problem

(\textup{A}^*\psi)(z) = \bar a_2 \partial_z^2 \psi(z)
+ \bar a_1 \partial_z \psi(z)
+ \bar a_0 \psi(z) \:,

where the \bar a_i are constant and whose boundary conditions are given by

\bar\alpha_1 \partial_z \psi(z_1) + \bar\alpha_0 \psi(z_1)
  &= 0 \\
\bar\beta_1 \partial_z \psi(z_2) + \bar\beta_0 \psi(z_2)
  &= 0 .

The following mapping is used:

\bar a_2 = a_2,  \quad\bar a_1 = -a_1, \quad\bar a_0 = a_0, \\
\bar\alpha_1 = -1, \quad
     \bar\alpha_0 = \frac{a_1}{a_2} - \frac{\alpha_0}{\alpha_1}
     , \\
\bar\beta_1 = -1, \quad
    \bar\beta_0 = \frac{a_1}{a_2} - \frac{\beta_0}{\beta_1} \:.

Returns:Parameter set describing A^* .
Return type:SecondOrderOperator
class SecondOrderEigenVector(char_pair, coefficients, domain, derivative_order)[source]

Bases: pyinduct.core.Function

This class provides eigenvectors of the form

\varphi(z) =  e^{\eta z} \big(\kappa_1 \cos(\nu z)
+ \sin(\nu z) \big) \:,

of a linear second order spatial operator \textup{A} denoted by

(\textup{A}\varphi)(z) = a_2 \partial_z^2 \varphi(z)
+ a_1 \partial_z \varphi(z)
+ a_0 \varphi(z)

where the a_i are constant and whose boundary conditions are given by

\alpha_1 \partial_z x(z_1) + \alpha_0 x(z_1) &= 0 \\
\beta_1 \partial_z x(z_2) + \beta_0 x(z_2) &= 0 .

To calculate the corresponding eigenvectors, the problem

(\textup{A}\varphi)(z) = \lambda \varphi(z)

is solved for the eigenvalues \lambda , making use of the characteristic roots p given by

p = \underbrace{-\frac{a_1}{a_2}}_{=:\eta}
+ j\underbrace{\sqrt{ \frac{a_0 - \lambda}{a_2}
        - \left(\frac{a_1}{2a_2}\right)^2 }}_{=:\nu}

Note

To easily instantiate a set of eigenvectors for a certain system, use the cure_hint() of this class or even better the helper-function cure_interval() .

Warn:
Since an eigenvalue corresponds to a pair of conjugate complex characteristic roots, latter are only calculated for the positive half-plane since the can be mirrored. To obtain the orthonormal properties of the generated eigenvectors, the eigenvalue corresponding to the characteristic root 0+0j is ignored, since it leads to the zero function.
Parameters:
  • char_pair (tuple of complex) – Characteristic root, corresponding to the eigenvalue \lambda for which the eigenvector is to be determined. (Can be obtained by convert_to_characteristic_root())
  • coefficients (tuple) – Constants of the exponential ansatz solution.
Returns:

The eigenvector.

Return type:

SecondOrderEigenVector

static calculate_eigenvalues(domain, params, count, extended_output=False, **kwargs)[source]

Determine the eigenvalues of the problem given by parameters defined on domain .

Parameters:
  • domain (Domain) – Domain of the spatial problem.
  • params (bunch-like) – Parameters of the system, see __init__() for details on their definition. Long story short, it must contain a_2, a_1, a_0, \alpha_0, \alpha_1,
\beta_0 \text{ and } \beta_1 .
  • count (int) – Amount of eigenvalues to generate.
  • extended_output (bool) – If true, not only eigenvalues but also the corresponding characteristic roots and coefficients of the eigenvectors are returned. Defaults to False.
Keyword Arguments:
 

debug (bool) – If provided, this parameter will cause several debug windows to open.

Returns:

\lambda , ordered in increasing order or tuple of (\lambda, p, \boldsymbol{\kappa} ) if extended_output is True.

Return type:

array or tuple of arrays

static convert_to_characteristic_root(params, eigenvalue)[source]

Converts a given eigenvalue \lambda into a characteristic root p by using the provided parameters. The relation is given by

p = -\frac{a_1}{a_2}
+ j\sqrt{ \frac{a_0 - \lambda}{a_2}
        - \left(\frac{a_1}{2a_2}\right)^2 }

Parameters:
  • params (bunch) – system parameters, see cure_hint() .
  • eigenvalue (real) – eigenvalue \lambda
Returns:

characteristic root p

Return type:

complex number

static convert_to_eigenvalue(params, char_roots)[source]

Converts a pair of characteristic roots p_{1,2} into an eigenvalue \lambda by using the provided parameters. The relation is given by

\lambda = a_2 p^2 + a_1 p + a_0

Parameters:
  • params (SecondOrderOperator) – System parameters.
  • char_roots (tuple or array of tuples) – Characteristic roots
static cure_hint(domain, params, count, derivative_order, **kwargs)[source]

Helper to cure an interval with eigenvectors.

Parameters:
  • domain (Domain) – Domain of the spatial problem.
  • params (SecondOrderOperator) – Parameters of the system, see __init__() for details on their definition. Long story short, it must contain a_2, a_1, a_0, \alpha_0, \alpha_1, \beta_0 \text{ and }
\beta_1 .
  • count (int) – Amount of eigenvectors to generate.
  • derivative_order (int) – Amount of derivative handles to provide.
  • kwargs – will be passed to calculate_eigenvalues()
Keyword Arguments:
 

debug (bool) – If provided, this parameter will cause several debug windows to open.

Returns:

An array holding the eigenvalues paired with a basis spanned by the eigenvectors.

Return type:

tuple of (array, Base)

class SecondOrderEigenfunction[source]

Bases: object

Wrapper for all eigenvalue problems of the form

a_2\varphi''(z) + a_1&\varphi'(z) + a_0\varphi(z) =
\lambda\varphi(z), \qquad a_2, a_1, a_0, \lambda \in \mathbb R

with eigenfunctions \varphi and eigenvalues \lambda. The roots of the characteristic equation (belonging to the ode) are denoted by

p = \eta \pm j\omega, \qquad \eta \in \mathbb R,
\quad \omega \in \mathbb C

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

In the following the variable \omega is called an eigenfrequency.

static eigfreq_eigval_hint(param, l, n_roots)[source]
Parameters:
  • param (array_like) – Parameters (a_2, a_1, a_0, None, None).
  • l – End of the domain z\in[0, 1].
  • n_roots (int) – Number of eigenfrequencies/eigenvalues to be compute.
Returns:

Booth tuple elements are numpy.ndarrays of the same length, one for eigenfrequencies and one for eigenvalues.

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

Return type:

tuple

static eigval_tf_eigfreq(param, eig_val=None, eig_freq=None)[source]

Provide corresponding of eigenvalues/eigenfrequencies for given eigenfreqeuncies/eigenvalues, depending on which type is given.

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

respectively

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

Parameters:
  • param (array_like) – Parameters (a_2, a_1, a_0, None, None).
  • eig_val (array_like) – Eigenvalues \lambda.
  • eig_freq (array_like) – Eigenfrequencies \omega.
Returns:

Eigenfrequencies \omega or eigenvalues \lambda.

Return type:

numpy.array

static get_adjoint_problem(param)[source]

Return the parameters of the adjoint eigenvalue problem for the given parameter set. Hereby, dirichlet or robin boundary condition at z=0

\varphi(0) = 0 \quad &\text{or} \quad
\varphi'(0) = \alpha\varphi(0)

and dirichlet or robin boundary condition at z=l

\varphi`(l) = 0 \quad &\text{or} \quad
\varphi'(l) = -\beta\varphi(l)

can be imposed.

Parameters:param (array_like) –

To define a homogeneous dirichlet boundary condition set alpha or beta to None at the corresponding side. Possibilities:

  • \Big( a_2, a_1, a_0, \alpha, \beta \Big)^T,
  • \Big( a_2, a_1, a_0, None, \beta \Big)^T,
  • \Big( a_2, a_1, a_0, \alpha, None \Big)^T or
  • \Big( a_2, a_1, a_0, None, None \Big)^T.
Returns:Parameters \big(a_2, \tilde a_1, a_0, \tilde \alpha, \tilde \beta \big) for the adjoint problem

a_2\psi''(z) + \tilde 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)

with

\tilde a_1 = -a_1, \quad \tilde\alpha = \frac{a_1}{a_2}\alpha, \quad \tilde\beta = -\frac{a_1}{a_2}\beta.

Return type:tuple
classmethod solve_evp_hint(param, l, n=None, eig_val=None, eig_freq=None, max_order=2, scale=None)[source]

Provide the first n eigenvalues and eigenfunctions (wraped inside a pyinduct base). For the exact formulation of the considered eigenvalue problem, have a look at the docstring from the eigenfunction class from which you will call this method.

You must call this classmethod with one and only one of the kwargs:

or (and) pass the kwarg scale (then n is set to len(scale)). If you have the kwargs eig_val and eig_freq already calculated then these are preferable, in the sense of performance.

Parameters:
  • param – Parameters (a_2, a_1, a_0, ...) see evp_class.__doc__.
  • l – End of the eigenfunction domain (start is 0).
  • n – Number of eigenvalues/eigenfunctions to compute.
  • eig_freq (array_like) – Pass your own choice of eigenfrequencies here.
  • eig_val (array_like) – Pass your own choice of eigenvalues here.
  • max_order – Maximum derivative order which must provided by the eigenfunctions.
  • scale (array_like) – Here you can pass a list of values to scale the eigenfunctions.
Returns:

Tuple with one list for the eigenvalues and one base which fractions are the eigenfunctions.

class SecondOrderDirichletEigenfunction(om, param, l, scale=1, max_der_order=2)[source]

Bases: pyinduct.eigenfunctions.LambdifiedSympyExpression, pyinduct.eigenfunctions.SecondOrderEigenfunction

This class provides an eigenfunction \varphi(z) to eigenvalue problems 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 (for example with the eigfreq_eigval_hint() of this class).

Parameters:
  • om (numbers.Number) – eigenfrequency \omega
  • param (array_like) – \Big( a_2, a_1, a_0, None, None \Big)^T
  • l (numbers.Number) – End of the domain z\in [0,l].
  • scale (numbers.Number) – Factor to scale the eigenfunctions.
  • max_der_order (int) – Number of derivative handles that are needed.
static eigfreq_eigval_hint(param, l, n_roots)[source]

Return the first n_roots eigenfrequencies \omega and eigenvalues \lambda.

\omega_i = \sqrt{-\frac{a_1^2}{4a_2^2}+
\frac{a_0-\lambda_i}{a_2}} \quad i = 1,...,\text{n\_roots}

to the considered eigenvalue problem.

Parameters:
  • param (array_like) – \Big( a_2, a_1, a_0, None, None \Big)^T
  • l (numbers.Number) – Right boundary value of the domain [0,l]\ni z.
  • n_roots (int) – Amount of eigenfrequencies to be compute.
Returns:

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

Return type:

tuple –> two numpy.ndarrays of length n_roots

class SecondOrderRobinEigenfunction(om, param, l, scale=1, max_der_order=2)[source]

Bases: pyinduct.core.Function, pyinduct.eigenfunctions.SecondOrderEigenfunction

This class provides an eigenfunction \varphi(z) to the eigenvalue problem given by

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 (for example with the eigfreq_eigval_hint() of this class).

Parameters:
  • om (numbers.Number) – eigenfrequency \omega
  • param (array_like) – \Big( a_2, a_1, a_0, \alpha, \beta \Big)^T
  • l (numbers.Number) – End of the domain z\in [0,l].
  • scale (numbers.Number) – Factor to scale the eigenfunctions (corresponds to \varphi(0)=\text{phi\_0}).
  • max_der_order (int) – Number of derivative handles that are needed.
static eigfreq_eigval_hint(param, l, n_roots, show_plot=False)[source]

Return the first n_roots eigenfrequencies \omega and eigenvalues \lambda.

\omega_i = \sqrt{
- \frac{a_1^2}{4a_2^2}
+ \frac{a_0 - \lambda_i}{a_2}}
\quad i = 1, \dotsc, \text{n\_roots}

to the considered eigenvalue problem.

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 compute.
  • show_plot (bool) – Show a plot window of the characteristic equation.
Returns:

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

Return type:

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

class TransformedSecondOrderEigenfunction(target_eigenvalue, init_state_vector, dgl_coefficients, domain)[source]

Bases: pyinduct.core.Function

This class provides an eigenfunction \varphi(z) to the eigenvalue problem given by

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

where \lambda \in \mathbb{C} denotes an eigenvalue and z \in [z_0, \dotsc, z_n] the domain.

Parameters:
  • target_eigenvalue (numbers.Number) – \lambda
  • init_state_vector (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) – Function handles \Big( a2(z), a1(z), a0(z) \Big)^T .
  • domain (Domain) – Spatial domain of the problem.
class AddMulFunction(function)[source]

Bases: object

(Temporary) Function class which 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 Function is overloaded with __add__ and __mul__ operator.

Parameters:function (callable) –
class LambdifiedSympyExpression(sympy_funcs, spat_symbol, spatial_domain)[source]

Bases: pyinduct.core.Function

This class provides a Function \varphi(z) based on a lambdified sympy expression. The sympy expressions for the function and it’s spatial derivatives must be provided as the list sympy_funcs. The expressions must be provided with increasing derivative order, starting with order 0.

Parameters:
  • sympy_funcs (array_like) – Sympy expressions for the function and the derivatives: \varphi(z), \varphi'(z), ....
  • spat_symbol – Sympy symbol for the spatial variable z.
  • spatial_domain (tuple) – Domain on which \varphi(z) is defined (e.g.: spatial_domain=(0, 1)).
class FiniteTransformFunction(function, M, l, scale_func=None, nested_lambda=False)[source]

Bases: pyinduct.core.Function

This class provides a transformed 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) that will act as start for the generation of 2n Functions \xi_{i,j}

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

    The vector of functions \boldsymbol\xi will then be constituted as follows:

    \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]).