Placeholder

In pyinduct.placeholder you find placeholders for symbolic Term definitions.

class FieldVariable(function_label, order=0, 0, weight_label=None, location=None, exponent=1, raised_spatially=False)[source]

Bases: pyinduct.placeholder.Placeholder

Class that represents terms of the systems field variable x(z, t).

Parameters
  • function_label (str) – Label of shapefunctions to use for approximation, see register_base() for more information about how to register an approximation basis.

  • tuple of int (order) – Tuple of temporal_order and spatial_order derivation order.

  • weight_label (str) – Label of weights for which coefficients are to be calculated (defaults to function_label).

  • location – Where the expression is to be evaluated.

  • exponent – Exponent of the term.

Examples

Assuming some shapefunctions have been registered under the label "phi" the following expressions hold:

  • \frac{\partial^{3}}{\partial t \partial z^2}x(z, t)

>>> x_dt_dzz = FieldVariable("phi", order=(1, 2))
  • \frac{\partial^2}{\partial t^2}x(3, t)

>>> x_dtt_at_3 = FieldVariable("phi", order=(2, 0), location=3)
class TestFunction(function_label, order=0, location=None, approx_label=None)[source]

Bases: pyinduct.placeholder.SpatialPlaceholder

Class that works as a placeholder for test functions in an equation.

Parameters
  • function_label (str) – Label of the function test base.

  • order (int) – Spatial derivative order.

  • location (Number) – Point of evaluation / argument of the function.

  • approx_label (str) – Label of the approximation test base.

class Base(fractions, matching_base_lbls=None, intermediate_base_lbls=None)[source]

Bases: pyinduct.core.ApproximationBasis

Base class for approximation bases.

In general, a Base is formed by a certain amount of BaseFractions and therefore forms finite-dimensional subspace of the distributed problem’s domain. Most of the time, the user does not need to interact with this class.

Parameters
  • fractions (iterable of BaseFraction) – List, array or dict of BaseFraction’s

  • matching_base_lbls (list of str) – List of labels from exactly matching bases, for which no transformation is necessary. Useful for transformations from bases that ‘live’ in different function spaces but evolve with the same time dynamic/coefficients (e.g. modal bases).

  • intermediate_base_lbls (list of str) – If it is certain that this base instance will be asked (as destination base) to return a transformation to a source base, whose implementation is cumbersome, its label can be provided here. This will trigger the generation of the transformation using build-in features. The algorithm, implemented in get_weights_transformation is then called again with the intermediate base as destination base and the ‘old’ source base. With this technique arbitrary long transformation chains are possible, if the provided intermediate bases again define intermediate bases.

derive(self, order)[source]

Basic implementation of derive function. Empty implementation, overwrite to use this functionality.

Parameters

order (numbers.Number) – derivative order

Returns

derived object

Return type

Base

function_space_hint(self)[source]

Hint that returns properties that characterize the functional space of the fractions. It can be used to determine if function spaces match.

Note

Overwrite to implement custom functionality.

get_attribute(self, attr)[source]

Retrieve an attribute from the fractions of the base.

Parameters

attr (str) – Attribute to query the fractions for.

Returns

Array of len(fractions) holding the attributes. With None entries if the attribute is missing.

Return type

np.ndarray

raise_to(self, power)[source]

Factory method to obtain instances of this base, raised by the given power.

Parameters

power – power to raise the basis onto.

scalar_product_hint(self)[source]

Hint that returns steps for scalar product calculation with elements of this base.

Note

Overwrite to implement custom functionality.

scale(self, factor)[source]

Factory method to obtain instances of this base, scaled by the given factor.

Parameters

factor – factor or function to scale this base with.

transformation_hint(self, info)[source]

Method that provides a information about how to transform weights from one BaseFraction into another.

In Detail this function has to return a callable, which will take the weights of the source- and return the weights of the target system. It may have keyword arguments for other data which is required to perform the transformation. Information about these extra keyword arguments should be provided in form of a dictionary whose keys are keyword arguments of the returned transformation handle.

Note

This implementation covers the most basic case, where the two BaseFraction’s are of same type. For any other case it will raise an exception. Overwrite this Method in your implementation to support conversion between bases that differ from yours.

Parameters

infoTransformationInfo

Raises

NotImplementedError

Returns

Transformation handle

class ConstantFunction(constant, **kwargs)[source]

Bases: pyinduct.core.Function

A Function that returns a constant value.

This function can be differentiated without limits.

Parameters

constant (number) – value to return

Keyword Arguments

**kwargs – All other kwargs get passed to Function.

derive(self, order=1)[source]

Spatially derive this Function.

This is done by neglecting order derivative handles and to select handle \text{order} - 1 as the new evaluation_handle.

Parameters

order (int) – the amount of derivations to perform

Raises
  • TypeError – If order is not of type int.

  • ValueError – If the requested derivative order is higher than the provided one.

Returns

Function the derived function.

class EquationTerm(scale, arg)[source]

Bases: object

Base class for all accepted terms in a weak formulation.

Parameters
  • scale

  • arg

class Function(eval_handle, domain=- np.inf, np.inf, nonzero=- np.inf, np.inf, derivative_handles=None)[source]

Bases: pyinduct.core.BaseFraction

Most common instance of a BaseFraction. This class handles all tasks concerning derivation and evaluation of functions. It is used broad across the toolbox and therefore incorporates some very specific attributes. For example, to ensure the accurateness of numerical handling functions may only evaluated in areas where they provide nonzero return values. Also their domain has to be taken into account. Therefore the attributes domain and nonzero are provided.

To save implementation time, ready to go version like LagrangeFirstOrder are provided in the pyinduct.simulation module.

For the implementation of new shape functions subclass this implementation or directly provide a callable eval_handle and callable derivative_handles if spatial derivatives are required for the application.

Parameters
  • eval_handle (callable) – Callable object that can be evaluated.

  • domain ((list of) tuples) – Domain on which the eval_handle is defined.

  • nonzero (tuple) – Region in which the eval_handle will return

  • output. Must be a subset of domain (nonzero) –

  • derivative_handles (list) – List of callable(s) that contain

  • of eval_handle (derivatives) –

add_neutral_element(self)[source]

Return the neutral element of addition for this object.

In other words: self + ret_val == self.

derivative_handles(self)
derive(self, order=1)[source]

Spatially derive this Function.

This is done by neglecting order derivative handles and to select handle \text{order} - 1 as the new evaluation_handle.

Parameters

order (int) – the amount of derivations to perform

Raises
  • TypeError – If order is not of type int.

  • ValueError – If the requested derivative order is higher than the provided one.

Returns

Function the derived function.

static from_data(x, y, **kwargs)[source]

Create a Function based on discrete data by interpolating.

The interpolation is done by using interp1d from scipy, the kwargs will be passed.

Parameters
  • x (array-like) – Places where the function has been evaluated .

  • y (array-like) – Function values at x.

  • **kwargs – all kwargs get passed to Function .

Returns

An interpolating function.

Return type

Function

function_handle(self)
function_space_hint(self)[source]

Return the hint that this function is an element of the an scalar product space which is uniquely defined by the scalar product scalar_product_hint().

Note

If you are working on different function spaces, you have to overwrite this hint in order to provide more properties which characterize your specific function space. For example the domain of the functions.

get_member(self, idx)[source]

Implementation of the abstract parent method.

Since the Function has only one member (itself) the parameter idx is ignored and self is returned.

Parameters

idx – ignored.

Returns

self

mul_neutral_element(self)[source]

Return the neutral element of multiplication for this object.

In other words: self * ret_val == self.

raise_to(self, power)[source]

Raises the function to the given power.

Warning

Derivatives are lost after this action is performed.

Parameters

power (numbers.Number) – power to raise the function to

Returns

raised function

scalar_product_hint(self)[source]

Return the hint that the _dot_product_l2() has to calculated to gain the scalar product.

scale(self, factor)[source]

Factory method to scale a Function.

Parameters

factornumbers.Number or a callable.

class Input(function_handle, index=0, order=0, exponent=1)[source]

Bases: pyinduct.placeholder.Placeholder

Class that works as a placeholder for an input of the system.

Parameters
  • function_handle (callable) – Handle that will be called by the simulation unit.

  • index (int) – If the system’s input is vectorial, specify the element to be used.

  • order (int) – temporal derivative order of this term (See Placeholder).

  • exponent (numbers.Number) – See FieldVariable.

Note

if order is nonzero, the callable is expected to return the temporal derivatives of the input signal by returning an array of len(order)+1.

class IntegralTerm(integrand, limits, scale=1.0)[source]

Bases: pyinduct.placeholder.EquationTerm

Class that represents an integral term in a weak equation.

Parameters
  • integrand

  • limits (tuple) –

  • scale

class ObserverGain(observer_feedback)[source]

Bases: pyinduct.placeholder.Placeholder

Class that works as a placeholder for the observer error gain.

Parameters

observer_feedback (ObserverFeedback) – Handle that will be called by the simulation unit.

class Placeholder(data, order=0, 0, location=None)[source]

Bases: object

Base class that works as a placeholder for terms that are later parsed into a canonical form.

Parameters
  • data (arbitrary) – data to store in the placeholder.

  • order (tuple) – (temporal_order, spatial_order) derivative orders that are to be applied before evaluation.

  • location (numbers.Number) – Location to evaluate at before further computation.

Todo

convert order and location into attributes with setter and getter methods. This will close the gap of unchecked values for order and location that can be sneaked in by the copy constructors by circumventing code doubling.

derivative(self, temp_order=0, spat_order=0)[source]

Mimics a copy constructor and adds the given derivative orders.

Note

The desired derivative order order is added to the original order.

Parameters
  • temp_order – Temporal derivative order to be added.

  • spat_order – Spatial derivative order to be added.

Returns

New Placeholder instance with the desired derivative order.

class Product(a, b=None)[source]

Bases: object

Represents a product.

Parameters
  • a

  • b

get_arg_by_class(self, cls)[source]

Extract element from product that is an instance of cls.

Parameters

cls

Returns

Return type

list

class ScalarFunction(function_label, order=0, location=None)[source]

Bases: pyinduct.placeholder.SpatialPlaceholder

Class that works as a placeholder for spatial functions in an equation. An example could be spatial dependent coefficients.

Parameters
  • function_label (str) – label under which the function is registered

  • order (int) – spatial derivative order to use

  • location – location to evaluate at

Warns
  • There seems to be a problem when this function is used in combination

  • with the :py:class:`.Product` class. Make sure to provide this class as

  • first argument to any product you define.

Todo

see warning.

static from_scalar(scalar, label, **kwargs)[source]

create a ScalarFunction from scalar values.

Parameters
  • scalar (array like) – Input that is used to generate the placeholder. If a number is given, a constant function will be created, if it is callable it will be wrapped in a Function and registered.

  • label (string) – Label to register the created base.

  • **kwargs – All kwargs that are not mentioned below will be passed to Function.

Keyword Arguments
  • order (int) – See constructor.

  • location (int) – See constructor.

  • overwrite (bool) – See register_base()

Returns

Placeholder object that can be used in a weak formulation.

Return type

ScalarFunction

class ScalarProductTerm(arg1, arg2, scale=1.0)[source]

Bases: pyinduct.placeholder.EquationTerm

Class that represents a scalar product in a weak equation.

Parameters
  • arg1 – Fieldvariable (Shapefunctions) to be projected.

  • arg2 – Testfunctions to project on.

  • scale (Number) – Scaling of expression.

class ScalarTerm(argument, scale=1.0)[source]

Bases: pyinduct.placeholder.EquationTerm

Class that represents a scalar term in a weak equation.

Parameters
  • argument

  • scale

class Scalars(values, target_term=None, target_form=None, test_func_lbl=None)[source]

Bases: pyinduct.placeholder.Placeholder

Placeholder for scalar values that scale the equation system, gained by the projection of the pde onto the test basis.

Note

The arguments target_term and target_form are used inside the parser. For frontend use, just specify the values.

Parameters
  • values – Iterable object containing the scalars for every k-th equation.

  • target_term – Coefficient matrix to add_to().

  • target_form – Desired weight set.

class SpatialDerivedFieldVariable(function_label, order, weight_label=None, location=None)[source]

Bases: pyinduct.placeholder.FieldVariable

Class that represents terms of the systems field variable x(z, t).

Parameters
  • function_label (str) – Label of shapefunctions to use for approximation, see register_base() for more information about how to register an approximation basis.

  • tuple of int (order) – Tuple of temporal_order and spatial_order derivation order.

  • weight_label (str) – Label of weights for which coefficients are to be calculated (defaults to function_label).

  • location – Where the expression is to be evaluated.

  • exponent – Exponent of the term.

Examples

Assuming some shapefunctions have been registered under the label "phi" the following expressions hold:

  • \frac{\partial^{3}}{\partial t \partial z^2}x(z, t)

>>> x_dt_dzz = FieldVariable("phi", order=(1, 2))
  • \frac{\partial^2}{\partial t^2}x(3, t)

>>> x_dtt_at_3 = FieldVariable("phi", order=(2, 0), location=3)
class SpatialPlaceholder(data, order=0, location=None)[source]

Bases: pyinduct.placeholder.Placeholder

Base class for all spatially-only dependent placeholders. The deeper meaning of this abstraction layer is to offer an easier to use interface.

derive(self, order=1)[source]

Take the (spatial) derivative of this object. :param order: Derivative order.

Returns

The derived expression.

Return type

Placeholder

class TemporalDerivedFieldVariable(function_label, order, weight_label=None, location=None)[source]

Bases: pyinduct.placeholder.FieldVariable

Class that represents terms of the systems field variable x(z, t).

Parameters
  • function_label (str) – Label of shapefunctions to use for approximation, see register_base() for more information about how to register an approximation basis.

  • tuple of int (order) – Tuple of temporal_order and spatial_order derivation order.

  • weight_label (str) – Label of weights for which coefficients are to be calculated (defaults to function_label).

  • location – Where the expression is to be evaluated.

  • exponent – Exponent of the term.

Examples

Assuming some shapefunctions have been registered under the label "phi" the following expressions hold:

  • \frac{\partial^{3}}{\partial t \partial z^2}x(z, t)

>>> x_dt_dzz = FieldVariable("phi", order=(1, 2))
  • \frac{\partial^2}{\partial t^2}x(3, t)

>>> x_dtt_at_3 = FieldVariable("phi", order=(2, 0), location=3)
evaluate_placeholder_function(placeholder, input_values)[source]

Evaluate a given placeholder object, that contains functions.

Parameters
Returns

numpy.ndarray of results.

get_base(label)[source]

Retrieve registered set of initial functions by their label.

Parameters

label (str) – String, label of functions to retrieve.

Returns

initial_functions

get_common_form(placeholders)[source]

Extracts the common target form from a list of scalars while making sure that the given targets are equivalent.

Parameters

placeholders – Placeholders with possibly differing target forms.

Returns

Common target form.

Return type

str

get_common_target(scalars)[source]

Extracts the common target from list of scalars while making sure that targets are equivalent.

Parameters

scalars (Scalars) –

Returns

Common target.

Return type

dict

is_registered(label)[source]

Checks whether a specific label has already been registered.

Args: label (str): Label to check for.

Returns

True if registered, False if not.

Return type

bool

register_base(label, base, overwrite=False)[source]

Register a basis to make it accessible all over the pyinduct framework.

Parameters
  • base (ApproximationBase) – base to register

  • label (str) – String that will be used as label.

  • overwrite – Force overwrite if a basis is already registered under this label.

sanitize_input(input_object, allowed_type)[source]

Sanitizes input data by testing if input_object is an array of type allowed_type.

Parameters
  • input_object – Object which is to be checked.

  • allowed_type – desired type

Returns

input_object