Parallel3dEulerGeometry

class odl.tomo.geometry.parallel.Parallel3dEulerGeometry(apart, dpart, det_pos_init=(0, 1, 0), **kwargs)[source]

Bases: ParallelBeamGeometry

Parallel beam geometry in 3d.

The motion parameters are two or three Euler angles, and the detector is flat and two-dimensional.

In the standard configuration, the detector reference point starts at (0, 1, 0), and the initial detector axes are [(1, 0, 0), (0, 0, 1)].

For details, check the online docs.

Attributes:
angles

All angles of this geometry as an array.

check_bounds

If True, methods computing vectors check input arguments.

det_axes_init

Initial axes of the detector.

det_grid

Sampling grid of det_params.

det_params

Continuous detector parameter range, an IntervalProd.

det_partition

Partition of the detector parameter set into subsets.

det_pos_init

Initial position of the detector reference point.

detector

Detector representation of this geometry.

grid

Joined sampling grid for motion and detector.

implementation_cache

Dictionary acting as a cache for this geometry.

motion_grid

Sampling grid of motion_params.

motion_params

Continuous motion parameter range, an IntervalProd.

motion_partition

Partition of the motion parameter set into subsets.

ndim

Number of dimensions of the geometry.

params

Joined parameter set for motion and detector.

partition

Joined parameter set partition for motion and detector.

translation

Shift of the origin of this geometry.

Methods

det_axes(angles)

Return the detector axes tuple at angles.

det_point_position(mparam, dparam)

Return the detector point at (mparam, dparam).

det_refpoint(angle)

Return the position(s) of the detector ref.

det_to_src(angle, dparam)

Direction from a detector location to the source.

frommatrix(apart, dpart, init_matrix, **kwargs)

Create an instance of Parallel3dEulerGeometry using a matrix.

rotation_matrix(angles)

Return the rotation matrix to the system state at angles.

__init__(apart, dpart, det_pos_init=(0, 1, 0), **kwargs)[source]

Initialize a new instance.

Parameters:
apart2- or 3-dim. RectPartition

Partition of the angle parameter set.

dpart2-dim. RectPartition

Partition of the detector parameter set.

det_pos_initarray-like, shape (3,), optional

Initial position of the detector reference point.

Other Parameters:
det_axes_init2-tuple of array-like's (shape (3,)), optional

Initial axes defining the detector orientation. The default depends on det_pos_init, see Notes.

translationarray-like, shape (3,), optional

Global translation of the geometry. This is added last in any method that computes an absolute vector, e.g., det_refpoint, and also shifts the center of rotation. Default: (0, 0, 0)

check_boundsbool, optional

If True, methods computing vectors check input arguments. Checks are vectorized and add only a small overhead. Default: True

Notes

In the default configuration, the initial detector reference point is (0, 1, 0), and the initial detector axes are [(1, 0, 0), (0, 0, 1)]. If a different det_pos_init is chosen, the new default axes are given as a rotation of the original ones by a matrix that transforms (0, 1, 0) to the new (normalized) det_pos_init. This matrix is calculated with the rotation_matrix_from_to function. Expressed in code, we have

init_rot = rotation_matrix_from_to((0, 1, 0), det_pos_init)
det_axes_init[0] = init_rot.dot((1, 0, 0))
det_axes_init[1] = init_rot.dot((0, 0, 1))

Examples

Initialization with default parameters and 2 Euler angles:

>>> apart = odl.uniform_partition([0, 0], [np.pi, 2 * np.pi],
...                               (10, 20))
>>> dpart = odl.uniform_partition([-1, -1], [1, 1], (20, 20))
>>> geom = Parallel3dEulerGeometry(apart, dpart)
>>> geom.det_refpoint([0, 0])
array([ 0.,  1.,  0.])
>>> geom.det_point_position([0, 0], [-1, 1])
array([-1.,  1.,  1.])

Checking the default orientation:

>>> e_x, e_y, e_z = np.eye(3)  # standard unit vectors
>>> np.allclose(geom.det_pos_init, e_y)
True
>>> np.allclose(geom.det_axes_init, (e_x, e_z))
True

Specifying an initial detector position by default rotates the standard configuration to this position:

>>> geom = Parallel3dEulerGeometry(apart, dpart,
...                                det_pos_init=(1, 0, 0))
>>> np.allclose(geom.det_pos_init, e_x)
True
>>> np.allclose(geom.det_axes_init, (-e_y, e_z))
True
>>> geom = Parallel3dEulerGeometry(apart, dpart,
...                                det_pos_init=(0, 0, 1))
>>> np.allclose(geom.det_pos_init, e_z)
True
>>> np.allclose(geom.det_axes_init, (e_x, -e_y))
True

The initial detector axes can also be set explicitly:

>>> geom = Parallel3dEulerGeometry(
...     apart, dpart, det_pos_init=(-1, 0, 0),
...     det_axes_init=((0, 1, 0), (0, 0, 1)))
>>> np.allclose(geom.det_pos_init, -e_x)
True
>>> np.allclose(geom.det_axes_init, (e_y, e_z))
True