function

The function module defines the Evaluable class and derived objects, commonly referred to as nutils functions. They represent mappings from a nutils.topology onto Python space. The notabe class of Array objects map onto the space of Numpy arrays of predefined dimension and shape. Most functions used in nutils applicatons are of this latter type, including the geometry and function bases for analysis.

Nutils functions are essentially postponed python functions, stored in a tree structure of input/output dependencies. Many Array objects have directly recognizable numpy equivalents, such as Sin or Inverse. By not evaluating directly but merely stacking operations, complex operations can be defined prior to entering a quadrature loop, allowing for a higher level style programming. It also allows for automatic differentiation and code optimization.

It is important to realize that nutils functions do not map for a physical xy-domain but from a topology, where a point is characterized by the combination of an element and its local coordinate. This is a natural fit for typical finite element operations such as quadrature. Evaluation from physical coordinates is possible only via inverting of the geometry function, which is a fundamentally expensive and currently unsupported operation.

class nutils.function.Evaluable(args)

Bases: nutils.types.Singleton

Base class

dependencies

collection of all function arguments

ordereddeps

collection of all function arguments such that the arguments to dependencies[i] can be found in dependencies[:i]

dependencytree

lookup table of function arguments into ordereddeps, such that ordereddeps[i].__args[j] == ordereddeps[dependencytree[i][j]], and self.__args[j] == ordereddeps[dependencytree[-1][j]]

asciitree(self, richoutput=False)

string representation

graphviz(self, dotpath='dot', imgtype='png')

create function graph

stackstr(self, nlines=- 1)

print stack

prepare_eval(self, /, **kwargs)

Return a function tree suitable for evaluation.

exception nutils.function.EvaluationError(etype, evalue, evaluable, values)

Bases: Exception

evaluation error

__init__(self, etype, evalue, evaluable, values)

constructor

__str__(self)

string representation

__weakref__

list of weak references to the object (if defined)

class nutils.function.TransformChain(args, todims=None)

Bases: nutils.function.Evaluable

Chain of affine transformations.

Evaluates to a tuple of nutils.transform.TransformItem objects.

nutils.function.dot(a, b, axes=None)

Contract a and b along axes.

class nutils.function.Array(args, shape, dtype)

Bases: nutils.function.Evaluable

Base class for array valued functions.

shape

The shape of this array function.

Type

tuple of ints

ndim

The number of dimensions of this array array function. Equal to len(shape).

Type

int

dtype

The dtype of the array elements.

Type

int, float

dot(a, b, axes=None)

Contract a and b along axes.

class nutils.function.Normal(lgrad)

Bases: nutils.function.Array

normal

class nutils.function.Inverse(func)

Bases: nutils.function.Array

Matrix inverse of func over the last two axes. All other axes are treated element-wise.

class nutils.function.Interpolate(x, xp, fp, left=None, right=None)

Bases: nutils.function.Array

interpolate uniformly spaced data; stepwise for now

class nutils.function.Pointwise(*args)

Bases: nutils.function.Array

Abstract base class for pointwise array functions.

class nutils.function.Cos(*args)

Bases: nutils.function.Pointwise

Cosine, element-wise.

evalf = <ufunc 'cos'>
class nutils.function.Sin(*args)

Bases: nutils.function.Pointwise

Sine, element-wise.

evalf = <ufunc 'sin'>
class nutils.function.Tan(*args)

Bases: nutils.function.Pointwise

Tangent, element-wise.

evalf = <ufunc 'tan'>
class nutils.function.ArcSin(*args)

Bases: nutils.function.Pointwise

Inverse sine, element-wise.

evalf = <ufunc 'arcsin'>
class nutils.function.ArcCos(*args)

Bases: nutils.function.Pointwise

Inverse cosine, element-wise.

evalf = <ufunc 'arccos'>
class nutils.function.ArcTan(*args)

Bases: nutils.function.Pointwise

Inverse tangent, element-wise.

evalf = <ufunc 'arctan'>
class nutils.function.Sampled(points, expect)

Bases: nutils.function.Array

Basis-like identity operator.

Basis-like function that for every point in a predefined set evaluates to the unit vector corresponding to its index.

Parameters
  • points (1d Array) – Present point coordinates.

  • expect (2d Array) – Elementwise constant that evaluates to the predefined point coordinates; used for error checking and to inherit the shape.

class nutils.function.Zeros(shape, dtype)

Bases: nutils.function.Array

zero

class nutils.function.Guard(fun)

Bases: nutils.function.Array

bar all simplifications

class nutils.function.TrigNormal(angle)

Bases: nutils.function.Array

cos, sin

class nutils.function.TrigTangent(angle)

Bases: nutils.function.Array

-sin, cos

class nutils.function.Find(where)

Bases: nutils.function.Array

indices of boolean index vector

class nutils.function.DerivativeTargetBase(args, shape, dtype)

Bases: nutils.function.Array

base class for derivative targets

class nutils.function.Argument(name, shape, nderiv=0)

Bases: nutils.function.DerivativeTargetBase

Array argument, to be substituted before evaluation.

The Argument is an Array with a known shape, but whose values are to be defined later, before evaluation, e.g. using replace_arguments().

It is possible to take the derivative of an Array to an Argument:

>>> from nutils import function
>>> a = function.Argument('x', [])
>>> b = function.Argument('y', [])
>>> f = a**3 + b**2
>>> function.derivative(f, a).simplified == (3.*a**2).simplified
True

Furthermore, derivatives to the local cooardinates are remembered and applied to the replacement when using replace_arguments():

>>> from nutils import mesh
>>> domain, x = mesh.rectilinear([2,2])
>>> basis = domain.basis('spline', degree=2)
>>> c = function.Argument('c', basis.shape)
>>> replace_arguments(c.grad(x), dict(c=basis)) == basis.grad(x)
True
Parameters
  • name (str) – The Identifier of this argument.

  • shape (tuple of ints) – The shape of this argument.

  • nderiv (int, non-negative) – Number of times a derivative to the local coordinates is taken. Default: 0.

class nutils.function.LocalCoords(ndims)

Bases: nutils.function.DerivativeTargetBase

local coords derivative target

class nutils.function.DelayedJacobian(geom, *derivativestack)

Bases: nutils.function.Array

Placeholder for jacobian() until the dimension of the nutils.topology.Topology where this functions is being evaluated is known. The replacing is carried out by Evaluable.prepare_eval().

class nutils.function.Polyval(coeffs, points, ngrad=0)

Bases: nutils.function.Array

Computes the \(k\)-dimensional array

\[\begin{split}j_0,\dots,j_{k-1} \mapsto \sum_{\substack{i_0,\dots,i_{n-1}\in\mathbb{N}\\i_0+\cdots+i_{n-1}\le d}} p_0^{i_0} \cdots p_{n-1}^{i_{n-1}} c_{j_0,\dots,j_{k-1},i_0,\dots,i_{n-1}},\end{split}\]

where \(p\) are the \(n\)-dimensional local coordinates and \(c\) is the argument coeffs and \(d\) is the degree of the polynomial, where \(d\) is the length of the last \(n\) axes of coeffs.

Warning

All coefficients with a (combined) degree larger than \(d\) should be zero. Failing to do so won’t raise an Exception, but might give incorrect results.

class nutils.function.RevolutionAngle

Bases: nutils.function.Array

Pseudo coordinates of a nutils.topology.RevolutionTopology.

class nutils.function.Basis(ndofs, transforms, trans=<nutils.function.SelectChain object at 0x7f8005e46340>)

Bases: nutils.function.Array

Abstract base class for bases.

A basis is a sequence of elementwise polynomial functions.

Parameters

Notes

Subclasses must implement get_dofs() and get_coefficients() and if possible should redefine get_support().

get_support(self, dof)

Return the support of basis function dof.

If dof is an int, return the indices of elements that form the support of dof. If dof is an array, return the union of supports of the selected dofs as a unique array. The returned array is always unique, i.e. strict monotonic increasing.

Parameters

dof (int or array of int or bool) – Index or indices of basis function or a mask.

Returns

support – The elements (as indices) where function dof has support.

Return type

sorted and unique numpy.ndarray

abstract get_dofs(self, ielem)

Return an array of indices of basis functions with support on element ielem.

If ielem is an int, return the dofs on element ielem matching the coefficients array as returned by get_coefficients(). If ielem is an array, return the union of dofs on the selected elements as a unique array, i.e. a strict monotonic increasing array.

Parameters

ielem (int or array of int or bool) – Element number(s) or mask.

Returns

dofs – A 1D Array of indices.

Return type

numpy.ndarray

get_ndofs(self, ielem)

Return the number of basis functions with support on element ielem.

abstract get_coefficients(self, ielem)

Return an array of coefficients for all basis functions with support on element ielem.

Parameters

ielem (int) – Element number.

Returns

coefficients – Array of coefficients with shape (nlocaldofs,)+(degree,)*ndims, where the first axis corresponds to the dofs returned by get_dofs().

Return type

nutils.types.frozenarray

get_coeffshape(self, ielem)

Return the shape of the array of coefficients for basis functions with support on element ielem.

class nutils.function.PlainBasis(coefficients, dofs, ndofs, transforms, trans=<nutils.function.SelectChain object at 0x7f8005e46340>)

Bases: nutils.function.Basis

A general purpose implementation of a Basis.

Use this class only if there exists no specific implementation of Basis for the basis at hand.

Parameters
class nutils.function.DiscontBasis(coefficients, transforms, trans=<nutils.function.SelectChain object at 0x7f8005e46340>)

Bases: nutils.function.Basis

A discontinuous basis with monotonic increasing dofs.

Parameters
class nutils.function.MaskedBasis(parent, indices, trans=<nutils.function.SelectChain object at 0x7f8005e46340>)

Bases: nutils.function.Basis

An order preserving subset of another Basis.

Parameters
  • parent (Basis) – The basis to mask.

  • indices (array of ints) – The strict monotonic increasing indices of parent basis functions to keep.

  • trans (TransformChain) –

class nutils.function.StructuredBasis(coeffs, start_dofs, stop_dofs, dofs_shape, transforms, transforms_shape, trans=<nutils.function.SelectChain object at 0x7f8005e46340>)

Bases: nutils.function.Basis

A basis for class:nutils.transformseq.StructuredTransforms.

Parameters
  • coeffs (tuple of tuples of arrays) – Per dimension the coefficients of the basis functions per transform.

  • start_dofs (tuple of arrays of ints) – Per dimension the dof of the first entry in coeffs per transform.

  • stop_dofs (tuple of arrays of ints) – Per dimension one plus the dof of the last entry in coeffs per transform.

  • dofs_shape (tuple of ints) – The tensor shape of the dofs.

  • transforms (nutils.transformseq.Transforms) – The transforms on which this basis is defined.

  • transforms_shape (tuple of ints) – The tensor shape of the transforms.

  • trans (TransformChain) –

class nutils.function.PrunedBasis(parent, transmap, trans=<nutils.function.SelectChain object at 0x7f8005e46340>)

Bases: nutils.function.Basis

A subset of another Basis.

Parameters
  • parent (Basis) – The basis to prune.

  • transmap (one-dimensional array of ints) – The indices of transforms in parent that form this subset.

nutils.function.replace(func)

decorator for deep object replacement

Generates a deep replacement method for Immutable objects based on a callable that is applied (recursively) on individual constructor arguments.

Parameters

func – callable which maps (obj, …) onto replaced_obj

Returns

The method that searches the object to perform the replacements.

Return type

callable

nutils.function.partition(f, *levels)

Create a partition of unity for a scalar function f.

When n levels are specified, n+1 indicator functions are formed that evaluate to one if and only if the following condition holds:

indicator 0: f < levels[0]
indicator 1: levels[0] < f < levels[1]
...
indicator n-1: levels[n-2] < f < levels[n-1]
indicator n: f > levels[n-1]

At the interval boundaries the indicators evaluate to one half, in the remainder of the domain they evaluate to zero such that the whole forms a partition of unity. The partitions can be used to create a piecewise continuous function by means of multiplication and addition.

The following example creates a topology consiting of three elements, and a function f that is zero in the first element, parabolic in the second, and zero again in the third element.

>>> from nutils import mesh
>>> domain, x = mesh.rectilinear([3])
>>> left, center, right = partition(x[0], 1, 2)
>>> f = (1 - (2*x[0]-3)**2) * center
Parameters
  • f (Array) – Scalar-valued function

  • levels (scalar constants or Arrays) – The interval endpoints.

Returns

The indicator functions.

Return type

list of scalar Arrays

nutils.function.chain(funcs)
nutils.function.vectorize(args)

Combine scalar-valued bases into a vector-valued basis.

Parameters

args (iterable of 1-dimensional nutils.function.Array objects) –

Returns

Return type

Array

nutils.function.jacobian(geom, ndims)

Return \(\sqrt{|J^T J|}\) with \(J\) the gradient of geom to the local coordinate system with ndims dimensions (localgradient(geom, ndims)).

nutils.function.matmat(arg0, *args)

helper function, contracts last axis of arg0 with first axis of arg1, etc

nutils.function.derivative(func, var, seen=None)
nutils.function.localgradient(arg, ndims)

local derivative

nutils.function.outer(arg1, arg2=None, axis=0)

outer product

nutils.function.find(arg)
nutils.function.J(geometry, ndims=None)

Return \(\sqrt{|J^T J|}\) with \(J\) the gradient of geometry to the local coordinate system with ndims dimensions (localgradient(geom, ndims)).

nutils.function.replace_arguments(value, arguments)

Replace Argument objects in value.

Replace Argument objects in value according to the arguments map, taking into account derivatives to the local coordinates.

Parameters
Returns

The edited value.

Return type

Array

class nutils.function.Namespace(*, default_geometry_name='x', fallback_length=None, **kwargs)

Bases: object

Namespace for Array objects supporting assignments with tensor expressions.

The Namespace object is used to store Array objects.

>>> from nutils import function
>>> ns = function.Namespace()
>>> ns.A = function.zeros([3, 3])
>>> ns.x = function.zeros([3])
>>> ns.c = 2

In addition to the assignment of Array objects, it is also possible to specify an array using a tensor expression string — see nutils.expression.parse() for the syntax. All attributes defined in this namespace are available as variables in the expression. If the array defined by the expression has one or more dimensions the indices of the axes should be appended to the attribute name. Examples:

>>> ns.cAx_i = 'c A_ij x_j'
>>> ns.xAx = 'x_i A_ij x_j'

It is also possible to simply evaluate an expression without storing its value in the namespace by passing the expression to the method eval_ suffixed with appropriate indices:

>>> ns.eval_('2 c')
Array<>
>>> ns.eval_i('c A_ij x_j')
Array<3>
>>> ns.eval_ij('A_ij + A_ji')
Array<3,3>

For zero and one dimensional expressions the following shorthand can be used:

>>> '2 c' @ ns
Array<>
>>> 'A_ij x_j' @ ns
Array<3>

Sometimes the dimension of an expression cannot be determined, e.g. when evaluating the identity array:

>>> ns.eval_ij('δ_ij')
Traceback (most recent call last):
...
nutils.expression.ExpressionSyntaxError: Length of axis cannot be determined from the expression.
δ_ij
  ^

There are two ways to inform the namespace of the correct lengths. The first is to assign fixed lengths to certain indices via keyword argument length_<indices>:

>>> ns_fixed = function.Namespace(length_ij=2)
>>> ns_fixed.eval_ij('δ_ij')
Array<2,2>

Note that evaluating an expression with an incompatible length raises an exception:

>>> ns = function.Namespace(length_i=2)
>>> ns.a = numpy.array([1,2,3])
>>> 'a_i' @ ns
Traceback (most recent call last):
...
nutils.expression.ExpressionSyntaxError: Length of index i is fixed at 2 but the expression has length 3.
a_i
  ^

The second is to define a fallback length via the fallback_length argument:

>>> ns_fallback = function.Namespace(fallback_length=2)
>>> ns_fallback.eval_ij('δ_ij')
Array<2,2>

When evaluating an expression through this namespace the following functions are available: opposite, sin, cos, tan, sinh, cosh, tanh, arcsin, arccos, arctan2, arctanh, exp, abs, ln, log, log2, log10, sqrt and sign.

Parameters
  • default_geometry_name (str) – The name of the default geometry. This argument is passed to nutils.expression.parse(). Default: 'x'.

  • fallback_length (int, optional) – The fallback length of an axis if the length cannot be determined from the expression.

  • length_<indices> (int) – The fixed length of <indices>. All axes in the expression marked with one of the <indices> are asserted to have the specified length.

arg_shapes

A readonly map of argument names and shapes.

Type

dict

default_geometry_name

The name of the default geometry. See argument with the same name.

Type

str

__getstate__(self)

Pickle instructions

__setstate__(self, d)

Unpickle instructions

property default_geometry

The default geometry, shorthand for getattr(ns, ns.default_geometry_name).

Type

nutils.function.Array

__call__(*args, **subs)

Return a copy with arguments replaced by subs.

Return a copy of this namespace with Argument objects replaced according to subs.

Parameters

**subs (dict of str and nutils.function.Array objects) – Replacements of the Argument objects, identified by their names.

Returns

ns – The copy of this namespace with replaced Argument objects.

Return type

Namespace

copy_(self, *, default_geometry_name=None)

Return a copy of this namespace.

__getattr__(self, name)

Get attribute name.

__setattr__(self, name, value)

Set attribute name to value.

__delattr__(self, name)

Delete attribute name.

__rmatmul__(self, expr)

Evaluate zero or one dimensional expr or a list of expressions.