evaluable
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.evaluable.Evaluable
Bases:
DataClass
Base class
- property arguments
a frozenset of all arguments of this evaluable
- asciitree(self, richoutput=False)
string representation
- eval(self, /, **evalargs)
Evaluate function on a specified element, point set.
- argument_degree(self, argument)
return the highest power of argument of self is polynomial, or raise NotPolynomal otherwise.
- nutils.evaluable.sum(arg, axis=None)
Sum array elements over a given axis.
- nutils.evaluable.dot(a, b, axes)
Contract
a
andb
alongaxes
.
- nutils.evaluable.align(arg, where, shape)
Align array to target shape.
The align operation can be considered the opposite of transpose: instead of specifying for each axis of the return value the original position in the argument, align specifies for each axis of the argument the new position in the return value. In addition, the return value may be of higher dimension, with new axes being inserted according to the
shape
argument.
- nutils.evaluable.unalign(*args, naxes=None)
Remove (joint) inserted axes.
Given one or more array arguments, return the shortest common axis vector along with function arguments such that the original arrays can be recovered by
align()
. Axes beyond the firstnaxes
are not considered for removal, keep their position (as seen from the right), and are not part of the common axis vector. Those axes should be added to the axis vector before callingalign()
.If
naxes
isNone
(the default), all arguments must have the same number of axes andnaxes
is set to this number.
- class nutils.evaluable.AsEvaluableArray
Bases:
object
Protocol for conversion into an
Array
.- property as_evaluable_array: Array
Lower this object to a
nutils.evaluable.Array
.
- __weakref__
list of weak references to the object
- class nutils.evaluable.Array
Bases:
Evaluable
Base class for array valued functions.
- sum(arg, axis=None)
Sum array elements over a given axis.
- dot(a, b, axes)
Contract
a
andb
alongaxes
.
- property as_evaluable_array
return self
- class nutils.evaluable.Orthonormal(basis, vector)
Bases:
Array
make a vector orthonormal to a subspace
- class nutils.evaluable.Inverse(func)
Bases:
Array
Matrix inverse of
func
over the last two axes. All other axes are treated element-wise.
- class nutils.evaluable.Pointwise
Bases:
Array
Abstract base class for pointwise array functions.
- classmethod outer(*args)
Alternative constructor that outer-aligns the arguments.
The output shape of this pointwise function is the sum of all shapes of its arguments. When called with multiple arguments, the first argument will be appended with singleton axes to match the output shape, the second argument will be prepended with as many singleton axes as the dimension of the original first argument and appended to match the output shape, and so forth and so on.
- class nutils.evaluable.Holomorphic(arg)
Bases:
Pointwise
Abstract base class for holomorphic array functions.
- class nutils.evaluable.Cos(arg)
Bases:
Holomorphic
Cosine, element-wise.
- class nutils.evaluable.Sin(arg)
Bases:
Holomorphic
Sine, element-wise.
- class nutils.evaluable.Tan(arg)
Bases:
Holomorphic
Tangent, element-wise.
- class nutils.evaluable.ArcSin(arg)
Bases:
Holomorphic
Inverse sine, element-wise.
- class nutils.evaluable.ArcCos(arg)
Bases:
Holomorphic
Inverse cosine, element-wise.
- class nutils.evaluable.ArcTan(arg)
Bases:
Holomorphic
Inverse tangent, element-wise.
- class nutils.evaluable.CosH(arg)
Bases:
Holomorphic
Hyperbolic cosine, element-wise.
- class nutils.evaluable.SinH(arg)
Bases:
Holomorphic
Hyperbolic sine, element-wise.
- class nutils.evaluable.TanH(arg)
Bases:
Holomorphic
Hyperbolic tangent, element-wise.
- class nutils.evaluable.ArcTanH(arg)
Bases:
Holomorphic
Inverse hyperbolic tangent, element-wise.
- class nutils.evaluable.Sampled(points, target, interpolation)
Bases:
Array
Basis-like identity operator.
Basis-like function that for every point evaluates to a partition of unity of a predefined set based on the selected interpolation scheme.
- Parameters:
- class nutils.evaluable.WithDerivative(func, var, derivative)
Bases:
Array
Wrap the given function and define the derivative to a target.
The wrapper is typically used together with a virtual derivative target like
IdentifierDerivativeTarget
. The wrapper is removed in the simplified form.- Parameters:
func (
Array
) – The function to wrap.var (
DerivativeTargetBase
) – The derivative target.derivative (
Array
) – The derivative with shapefunc.shape + var.shape
.
See also
IdentifierDerivativeTarget
a virtual derivative target
- class nutils.evaluable.Argument(name, shape, dtype=<class 'float'>)
Bases:
DerivativeTargetBase
Array argument, to be substituted before evaluation.
The
Argument
is anArray
with a known shape, but whose values are to be defined later, before evaluation, e.g. usingreplace_arguments()
.It is possible to take the derivative of an
Array
to anArgument
:>>> from nutils import evaluable >>> a = evaluable.Argument('x', ()) >>> b = evaluable.Argument('y', ()) >>> f = a**3. + b**2. >>> evaluable.derivative(f, b).simplified == 2.*b True
- Parameters:
- class nutils.evaluable.IdentifierDerivativeTarget(identifier, shape)
Bases:
DerivativeTargetBase
Virtual derivative target distinguished by an identifier.
- Parameters:
See also
WithDerivative
Array
wrapper with additional derivative
- class nutils.evaluable.Polyval(coeffs, points)
Bases:
Array
Evaluate a polynomial
The polynomials are of the form
\[Σ_{k ∈ ℤ^n | Σ_i k_i ≤ p} c_k ∏_i x_i^(k_i)\]where \(c\) is a vector of coefficients, \(x\) a vector of \(n\) variables and \(p\) a nonnegative integer degree. The coefficients are assumed to be in reverse [lexicographic order]: the coefficient for powers \(j ∈ ℤ^n\) comes before the coefficient for powers \(k ∈ ℤ^n / {j}\) iff \(j_i > k_i\), where \(i = max_l(j_l ≠ k_l)\), the index of the last non-matching power.
- Parameters:
- class nutils.evaluable.PolyDegree(ncoeffs, nvars)
Bases:
Array
Returns the degree of a polynomial given the number of coefficients and number of variables
- Parameters:
Notes
See
Polyval
for a definition of the polynomial.
- class nutils.evaluable.PolyNCoeffs(nvars, degree)
Bases:
Array
Returns the number of coefficients for a polynomial of given degree and number of variables
- Parameters:
Notes
See
Polyval
for a definition of the polynomial.
- class nutils.evaluable.PolyMul(coeffs_left, coeffs_right, vars)
Bases:
Array
Compute the coefficients for the product of two polynomials
Return the coefficients such that calling
Polyval
on this result is equal to the product ofPolyval
called on the individual arrays of coefficients (with the appropriate selection of the variables as described by parametervars
).- Parameters:
coeffs_left (
Array
) – The coefficients for the left operand. The last axis is treated as the coefficients axis.coeffs_right (
Array
) – The coefficients for the right operand. The last axis is treated as the coefficients axis.vars (
tuple
ofnutils_poly.MulVar
) – For each variable of this product,var
defines if the variable exists in the left polynomial, the right or both.
Notes
See
Polyval
for a definition of the polynomial.
- class nutils.evaluable.PolyGrad(coeffs, nvars)
Bases:
Array
Compute the coefficients for the gradient of a polynomial
The last two axes of this array are the axis of variables and the axis of coefficients.
- Parameters:
Notes
See
Polyval
for a definition of the polynomial.
- class nutils.evaluable.Legendre(x, degree)
Bases:
Array
Series of Legendre polynomial up to and including the given degree.
- Parameters:
- class nutils.evaluable.Choose(index, choices)
Bases:
Array
Function equivalent of
numpy.choose()
.
- class nutils.evaluable.TransformCoords(target, source, index, coords)
Bases:
Array
Transform coordinates from one coordinate system to another (spatial part)
- Parameters:
target (
nutils.transformseq.Transforms
, optional) – The target coordinate system. If None the target is root coordinate system.source (
nutils.transformseq.Transforms
) – The source coordinate system.index (scalar, integer
Array
) – The index part of the source coordinates.coords (
Array
) – The spatial part of the source coordinates.
- class nutils.evaluable.TransformIndex(target, source, index)
Bases:
Array
Transform coordinates from one coordinate system to another (index part)
- Parameters:
target (
nutils.transformseq.Transforms
) – The target coordinate system.source (
nutils.transformseq.Transforms
) – The source coordinate system.index (scalar, integer
Array
) – The index part of the source coordinates.
- class nutils.evaluable.TransformLinear(target, source, index)
Bases:
Array
Linear part of a coordinate transformation
- Parameters:
target (
nutils.transformseq.Transforms
, optional) – The target coordinate system. If None the target is the root coordinate system.source (
nutils.transformseq.Transforms
) – The source coordinate system.index (scalar, integer
Array
) – The index part of the source coordinates.
- class nutils.evaluable.TransformBasis(source, index)
Bases:
Array
Vector basis for the root and a source coordinate system
The columns of this matrix form a vector basis for the space of root coordinates. The first n vectors also span the space of source coordinates mapped to the root, where n is the dimension of the source coordinate system. The remainder is not a span of the complement space in general.
No additional properties are guaranteed beyond the above. In particular, if the source coordinate system has the same dimension as the root, the basis is not necessarily the same as
TransformLinear(None, source, index)
.- Parameters:
source (
nutils.transformseq.Transforms
) – The source coordinate system.index (scalar, integer
Array
) – The index part of the source coordinates.
- class nutils.evaluable.Loop(loop_id, length)
Bases:
Array
Base class for evaluable loops.
Subclasses must implement
method
evalf_loop_init(init_arg)
andmethod
evalf_loop_body(output, body_arg)
.
- class nutils.evaluable.SearchSorted(arg, array, sorter, side)
Bases:
Array
Find index of evaluable array into sorted numpy array.
- nutils.evaluable.zero_all_arguments(value)
Replace all function arguments by zeros.
- class nutils.evaluable.Monomial(values, args, indices, powers)
Bases:
Array
Helper object for factor.
Performs a sparse tensor multiplication, without summation, and returns the result as a dense vector; inflation and reshape is the responsibility of factor. The factors of the multiplication are
values
and all theargs
, which are scattered into values via theindices
. Thepowers
argument contains multiplicities, for the following reason.With reference to the example of the factor doc string, derivative will generate an evaluable of the form array’(arg) darg = darray_darg(0) darg + .5 arg d2array_darg2 darg + .5 darg d2array_darg2 arg. By Schwarz’s theorem d2array_darg2 is symmetric, and the latter two terms are equal. However, since the row and column indices of d2array_darg2 differ, we cannot detect this equality but rather need to embed the information explicitly.
In this situation, the
powers
argument contains the value 2 to indicate that its position is symmetric with the next, and the first integral can therefore be doubled. With that, the derivative takes the desired form of array’(arg) == darray_darg(0) + d2array_darg2 arg.
- nutils.evaluable.factor(array)
Convert array to a sparse polynomial.
This function forms the equivalent polynomial of an evaluable array, of which the coefficients are already evaluated, with the aim to speed up repeated evaluations in e.g. time or Newton loops.
For example, if
array
is a quadratic function ofarg
, then thefactor
function returns the equivalent Taylor series:array(arg) -> array(0) + darray_darg(0) arg + .5 arg d2array_darg2 arg
The coefficients
array(0)
,darray_darg(0)
andd2array_darg2
are evaluated as sparse tensors. As such, the new representation retains all sparsity of the original.
- nutils.evaluable.derivative(func, var, seen=None)
- nutils.evaluable.prependaxes(func, shape)
Prepend axes with specified shape to func.
- nutils.evaluable.appendaxes(func, shape)
Append axes with specified shape to func.
- nutils.evaluable.replace_arguments(value, arguments)
Replace
Argument
objects invalue
.Replace
Argument
objects invalue
according to thearguments
map, taking into account derivatives to the local coordinates.
- nutils.evaluable.einsum(fmt, *args, **dims)
Multiply and/or contract arrays via format string.
The format string consists of a comma separated list of axis labels, followed by
->
and the axis labels of the return value. For example, the following swaps the axes of a matrix:>>> a45 = ones(tuple(map(constant, [4,5]))) # 4x5 matrix >>> einsum('ij->ji', a45) nutils.evaluable.Transpose<f:(5),(4)>
Axis labels that do not occur in the return value are summed. For example, the following performs a matrix-vector product:
>>> a5 = ones(tuple(map(constant, [5]))) # vector with length 5 >>> einsum('ij,j->i', a45, a5) nutils.evaluable.Sum<f:4>
The following contracts a third order tensor, a matrix, and a vector, and transposes the result:
>>> a234 = ones(tuple(map(constant, [2,3,4]))) # 2x3x4 tensor >>> einsum('ijk,kl,l->ji', a234, a45, a5) nutils.evaluable.Sum<f:3,2>
In case the dimension of the input and output arrays may vary, a variable length axes group can be denoted by a capital. Its length is automatically established based on the dimension of the input arrays. The following example performs a tensor product of an array and a vector:
>>> einsum('A,i->Ai', a234, a5) nutils.evaluable.Multiply<f:2,3,4,5>
The format string may contain multiple variable length axes groups, but their lengths must be resolvable from left to right. In case this is not possible, lengths may be specified as keyword arguments.
>>> einsum('AjB,i->AijB', a234, a5, B=1) nutils.evaluable.Multiply<f:2,5,3,4>
- nutils.evaluable.eval_sparse(funcs, **arguments)
Evaluate one or several Array objects as sparse data.
- nutils.evaluable.eval_once(func, *, stats=None, arguments={}, _simplify=True, _optimize=True)
Evaluate one or several Array objects by compiling it for single use.
- Parameters:
func (
Evaluable
or (possibly nested) tuples ofEvaluable
s) – The function or functions to compile.stats (
'log'
orNone
) – If'log'
the compiled function will log the durations of individualEvaluable
s referenced byfunc
.arguments (
dict
(default: None)) – Optional arguments for function evaluation.
- Returns:
results
- Return type:
tuple
of (values, indices, shape) triplets
- nutils.evaluable.compile(func, /, *, stats=None, cache_const_intermediates=True, _simplify=True, _optimize=True)
Returns a callable that evaluates
func
.The return value of the callable matches the structure of
func
. Iffunc
is a singleEvaluable
, the returned callable returns a single array. Iffunc
is a tuple, the returned callable returns a tuple as well.If
func
depends onArgument
s, their values must be passed to the callable using keywords arguments. Arguments that are unknown tofunc
are silently ignored.- Parameters:
func (
Evaluable
or (possibly nested) tuples ofEvaluable
s) – The function or functions to compile.stats (
'log'
orNone
) – If'log'
the compiled function will log the durations of individualEvaluable
s referenced byfunc
.cache_const_intermediates (
bool
) – If true, the returned callable caches parts offunc
that can be reused for a second call.
- Returns:
The callable that evaluates
func
.- Return type:
Examples
Compiling a single function with argument
arg
.>>> arg = Argument('arg', (), int) >>> f = arg + constant(1) >>> compiled = compile(f) >>> compiled({'arg': 1}) 2
The return value of the compiled function matches the structure of the functions passed to
compile()
:>>> g = arg * constant(3) >>> h = arg * arg >>> compiled = compile((f, (g, h))) >>> compiled({'arg': 2}) (3, (6, 4))