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)

string representation

graphviz(self)

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

tuple of ints – The shape of this array function.

ndim

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

dtype

int, float – The dtype of the array elements.

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.BlockAdd(funcs)

Bases: nutils.function.Array

block addition (used for DG)

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.OldSampled(data, trans=<nutils.function.SelectChain object at 0x7f89ead160c8>)

Bases: nutils.function.Array

sampled

class nutils.function.Sampled(sample, array, trans=<nutils.function.SelectChain object at 0x7f89ead160c8>)

Bases: nutils.function.Array

Convert sampled data to evaluable array.

Using the result of nutils.sample.Sample.eval(), create an evaluable array that upon evaluation recovers the original function in the set of points matching the original sampling.

Parameters:
  • sample (nutils.sample.Sample) – The set of points that the data was sampled on.
  • array – The sampled data.
  • trans (TransformChain (optional)) – The transformation chain that is used to locate the sample points.
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.

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.bringforward(arg, axis)

bring axis forward

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.polyfunc(coeffs, dofs, ndofs, transforms, *, issorted=True)

Create an inflated Polyval with coefficients coeffs and corresponding dofs dofs. The arguments coeffs, dofs and transforms are assumed to have matching order. In addition, if issorted is true, the transforms argument is assumed to be sorted.

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')

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>

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'.
arg_shapes

view of dict – A readonly map of argument names and shapes.

default_geometry_name

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

default_geometry

nutils.function.Array – The default geometry, shorthand for getattr(ns, ns.default_geometry_name).

__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.