FanBeamGeometry¶
- class odl.tomo.geometry.conebeam.FanBeamGeometry(apart, dpart, src_radius, det_radius, det_curvature_radius=None, src_to_det_init=(0, 1), src_shift_func=None, det_shift_func=None, **kwargs)[source]¶
Bases:
DivergentBeamGeometry
Fan beam (2d cone beam) geometry.
The source moves on a circle with radius
src_radius
, and the detector reference point is opposite to the source, i.e. at maximum distance, on a circle with radiusdet_radius
. One of the two radii can be chosen as 0, which corresponds to a stationary source or detector, respectively.The motion parameter is the 1d rotation angle parameterizing source and detector positions simultaneously.
In the standard configuration, the detector is perpendicular to the ray direction, its reference point is initially at
(0, 1)
, and the initial detector axis is(1, 0)
.For details, check the online docs.
- Attributes:
angles
Discrete angles given in this geometry.
check_bounds
If
True
, methods computing vectors check input arguments.det_axis_init
Detector axis at angle 0.
det_curvature_radius
Detector curve radius of this geometry.
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_radius
Detector circle radius of this geometry.
det_shift_func
Detector shifts in the geometry.
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.
src_radius
Source circle radius of this geometry.
src_shift_func
Source shifts in the geometry.
src_to_det_init
Initial source-to-detector unit vector.
translation
Shift of the origin of this geometry.
Methods
det_axis
(angle)Return the detector axis at
angle
.det_point_position
(mparam, dparam)Return the detector point at
(mparam, dparam)
.det_refpoint
(angle)Return the detector reference point position at
angle
.det_to_src
(angle, dparam[, normalized])Vector or direction from a detector location to the source.
frommatrix
(apart, dpart, src_radius, ...[, ...])Create an instance of
FanBeamGeometry
using a matrix.rotation_matrix
(angle)Return the rotation matrix for
angle
.src_position
(angle)Return the source position at
angle
.- __init__(apart, dpart, src_radius, det_radius, det_curvature_radius=None, src_to_det_init=(0, 1), src_shift_func=None, det_shift_func=None, **kwargs)[source]¶
Initialize a new instance.
- Parameters:
- apart1-dim.
RectPartition
Partition of the angle interval.
- dpart1-dim.
RectPartition
Partition of the detector parameter interval.
- src_radiusnonnegative float
Radius of the source circle.
- det_radiusnonnegative float
Radius of the detector circle. Must be nonzero if
src_radius
is zero.- det_curvature_radiusnonnegative float, optional
Radius of the detector curvature. If
None
, a flat detector is used, otherwise must be positive.- src_to_det_init
array-like
(shape(2,)
), optional Initial state of the vector pointing from source to detector reference point. The zero vector is not allowed.
- src_shift_funccallable, optional
Function with signature
src_shift_func(angle) -> shift
returning a source shift for a given angle. Each shift is interpreted as a vector[shift_d, shift_t]
, where "d" and "t" denote shifts in the detector-to-source and tangent directions, respectively.- det_shift_funccallable, optional
Function with signature
det_shift_func(angle) -> shift
returning a detector shift for a given angle. Each shift is interpreted as a vector[shift_d, shift_t]
, where "d" and "t" denote shifts in the source-to-detector and tangent directions, respectively.
- apart1-dim.
- Other Parameters:
- det_axis_init
array-like
(shape(2,)
), optional Initial axis defining the detector orientation. The default depends on
src_to_det_init
, see Notes.- translation
array-like
, shape(2,)
, 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.- check_boundsbool, optional
If
True
, methods computing vectors check input arguments. Checks are vectorized and add only a small overhead. Default:True
- det_axis_init
Notes
In the default configuration, the initial source-to-detector vector is
(0, 1)
, and the initial detector axis is(1, 0)
. If a differentsrc_to_det_init
is chosen, the new default axis is given as a rotation of the original one by a matrix that transforms(0, 1)
to the new (normalized)src_to_det_init
. This matrix is calculated with therotation_matrix_from_to
function. Expressed in code, we haveinit_rot = rotation_matrix_from_to((0, 1), src_to_det_init) det_axis_init = init_rot.dot((1, 0))
Examples
Initialization with default parameters and some radii:
>>> apart = odl.uniform_partition(0, 2 * np.pi, 10) >>> dpart = odl.uniform_partition(-1, 1, 20) >>> geom = FanBeamGeometry(apart, dpart, src_radius=1, det_radius=5) >>> geom.src_position(0) array([ 0., -1.]) >>> geom.det_refpoint(0) array([ 0., 5.]) >>> geom.det_point_position(0, 1) # (0, 5) + 1 * (1, 0) array([ 1., 5.])
Checking the default orientation:
>>> geom.src_to_det_init array([ 0., 1.]) >>> geom.det_axis_init array([ 1., 0.])
Specifying curvature of the detector:
>>> apart = odl.uniform_partition(0, 2 * np.pi, 10) >>> dpart = odl.uniform_partition(-np.pi / 2, np.pi / 2, 10) >>> geom = FanBeamGeometry(apart, dpart, src_radius=1, det_radius=5, ... det_curvature_radius=10) >>> geom.src_position(0) array([ 0., -1.]) >>> geom.det_refpoint(0) array([ 0., 5.]) >>> # (0, 5) + 10 * (sin(pi/6), cos(pi/6) - 1) >>> np.round(geom.det_point_position(0, np.pi / 6), 2) array([ 5. , 3.66])
Specifying an initial detector position by default rotates the standard configuration to this position:
>>> e_x, e_y = np.eye(2) # standard unit vectors >>> geom = FanBeamGeometry(apart, dpart, src_radius=1, det_radius=5, ... src_to_det_init=(1, 0)) >>> np.allclose(geom.src_to_det_init, e_x) True >>> np.allclose(geom.det_axis_init, -e_y) True >>> geom = FanBeamGeometry(apart, dpart, src_radius=1, det_radius=5, ... src_to_det_init=(0, -1)) >>> np.allclose(geom.src_to_det_init, -e_y) True >>> np.allclose(geom.det_axis_init, -e_x) True
The initial detector axis can also be set explicitly:
>>> geom = FanBeamGeometry( ... apart, dpart, src_radius=1, det_radius=5, ... src_to_det_init=(1, 0), det_axis_init=(0, 1)) >>> np.allclose(geom.src_to_det_init, e_x) True >>> np.allclose(geom.det_axis_init, e_y) True
Specifying a flying focal spot and detector offset:
>>> apart = odl.uniform_partition(0, 2 * np.pi, 4) >>> geom = FanBeamGeometry( ... apart, dpart, ... src_radius=1, det_radius=5, ... src_shift_func=lambda angle: odl.tomo.flying_focal_spot( ... angle, apart=apart, shifts=[(0.1, 0), (0, 0.1)]), ... det_shift_func=lambda angle: [0.0, 0.05]) >>> geom.src_shift_func(geom.angles) array([[ 0.1, 0. ], [ 0. , 0.1], [ 0.1, 0. ], [ 0. , 0.1]]) >>> geom.det_shift_func(geom.angles) [0.0, 0.05]