PointwiseInner¶
-
class
odl.operator.tensor_ops.
PointwiseInner
(*args, **kwargs)[source]¶ Bases:
odl.operator.tensor_ops.PointwiseInnerBase
Take the point-wise inner product with a given vector field.
This operator takes the (weighted) inner product
<F(x), G(x)> = sum_j ( w_j * F_j(x) * conj(G_j(x)) )
for a given vector field
G
, whereF
is the vector field acting as a variable to this operator.This implies that the
Operator.domain
is a power space of a discretized function space. For example, ifX
is aDiscretizedSpace
space, thenProductSpace(X, d)
is a valid domain for any positive integerd
.- Attributes
adjoint
Adjoint of this operator.
base_space
Base space
X
of this operator’s domain and range.domain
Set of objects on which this operator can be evaluated.
inverse
Return the operator inverse.
is_functional
True
if this operator’s range is aField
.is_linear
True
if this operator is linear.is_weighted
True
if weighting is not 1 or all ones.range
Set in which the result of an evaluation of this operator lies.
vecfield
Fixed vector field
G
of this inner product.weights
Weighting array of this operator.
Methods
_call
(self, vf, out)Implement
self(vf, out)
.derivative
(self, point)Return the operator derivative at
point
.norm
(self[, estimate])Return the operator norm of this operator.
-
__init__
(self, vfspace, vecfield, weighting=None)[source]¶ Initialize a new instance.
- Parameters
- vfspace
ProductSpace
Space of vector fields on which the operator acts. It has to be a product space of identical spaces, i.e. a power space.
- vecfield
vfspace
element-like
Vector field with which to calculate the point-wise inner product of an input vector field
- weighting
array-like
or float, optional Weighting array or constant for the norm. If an array is given, its length must be equal to
len(domain)
, and all entries must be positive. By default, the weights are is taken fromdomain.weighting
. Note that this excludes unusual weightings with custom inner product, norm or dist.
- vfspace
Examples
We make a tiny vector field space in 2D and create the point-wise inner product operator with a fixed vector field. The operator maps a vector field to a scalar function:
>>> spc = odl.uniform_discr([-1, -1], [1, 1], (1, 2)) >>> vfspace = odl.ProductSpace(spc, 2) >>> fixed_vf = np.array([[[0, 1]], ... [[1, -1]]]) >>> pw_inner = PointwiseInner(vfspace, fixed_vf) >>> pw_inner.range == spc True
Now we can calculate the inner product in each point:
>>> x = vfspace.element([[[1, -4]], ... [[0, 3]]]) >>> print(pw_inner(x)) [[ 0., -7.]]