transforms

The BigFeta solver uses AlignerTransform objects to construct linear least squares elements from tilespecs. The final transform in each tilespec transform list is converted into an AlignerTransform object, which has a renderapi.transform.transform object as its base class. The AlignerTransform object has a few other methods to enable constructing the matrix, setting regularizations, extracting the starting transform values and setting the solved transform values.

class bigfeta.transform.transform.AlignerTransform(name=None, transform=None, fullsize=False, order=2)

Bases: object

general transform object that the solver expects

__init__(name=None, transform=None, fullsize=False, order=2)
Parameters:
  • name (str) – specifies the intended transform for the type of solve
  • transform (renderapi.transform.Transform) – The new AlignerTransform will inherit from this transform, if possible.
  • fullsize (bool) – only applies to affine transform. Remains for legacy reason as an explicit demonstration of the equivalence of fullsize and halfsize transforms.
  • order (int) – used in Polynomial2DTransform
class bigfeta.transform.affine_model.AlignerAffineModel(transform=None, fullsize=False)

Bases: renderapi.transform.leaf.affine_models.AffineModel

Object for implementing full or half-size affine transforms

__init__(transform=None, fullsize=False)
Parameters:
  • transform (renderapi.transform.Transform) – The new AlignerTransform will inherit from this transform, if possible.
  • fullsize (bool) – only applies to affine transform. Remains for legacy reason as an explicit demonstration of the equivalence of fullsize and halfsize transforms.
block_from_pts(pts, w, col_ind, col_max)

partial sparse block for a transform/match

Parameters:
  • pts (numpy.ndarray) – N x 2, the x, y values of the match (either p or q)
  • w (numpy.ndarray) – size N, the weights associated with the pts
  • col_ind (int) – the starting column index for this tile
  • col_max (int) – total number of columns in the matrix
Returns:

  • block (scipy.sparse.csr_matrix) – the partial block for this transform
  • w (numpy.ndarray) – the weights associated with the rows of this block
  • rhs (numpy.ndarray) – N/2 x 2 (halfsize) or N x 1 (fullsize) right hand side for this transform. generally all zeros. could implement fixed tiles in rhs later.

from_solve_vec(vec)

reads values from solution and sets transform parameters

Parameters:vec (numpy.ndarray) – input to this function is sliced so that vec[0] is the first harvested value for this transform
Returns:n – number of rows read from vec. Used to increment vec slice for next transform
Return type:int
regularization(regdict)

regularization vector from this transform

Parameters:regdict (dict) – bigfeta.schemas.regularization. controls regularization values
Returns:reg – array of regularization values of length DOF_per_tile
Return type:numpy.ndarray
to_solve_vec()

sets solve vector values from transform parameters

Returns:vec – N/2 x 2 for halfsize, N x 2 for fullsize
Return type:numpy.ndarray
class bigfeta.transform.polynomial_model.AlignerPolynomial2DTransform(transform=None, order=2)

Bases: renderapi.transform.leaf.polynomial_models.Polynomial2DTransform

Object for implementing half-size polynomial transforms

__init__(transform=None, order=2)
Parameters:
  • transform (renderapi.transform.Transform) – The new AlignerTransform will inherit from this transform, if possible.
  • order (int) – order of the intended polynomial.
block_from_pts(pts, w, col_ind, col_max)

partial sparse block for a transform/match

Parameters:
  • pts (numpy.ndarray) – N x 2, the x, y values of the match (either p or q)
  • w (numpy.ndarray) – the weights associated with the pts
  • col_ind (int) – the starting column index for this tile
  • col_max (int) – number of columns in the matrix
Returns:

  • block (scipy.sparse.csr_matrix) – the partial block for this transform
  • w (numpy.ndarray) – the weights associated with the rows of this block
  • rhs (numpy.ndarray) – N/2 x 2 (halfsize) or N x 1 (fullsize) right hand side for this transform. generally all zeros. could implement fixed tiles in rhs later.

from_solve_vec(vec)

reads values from solution and sets transform parameters

Parameters:vec (numpy.ndarray) – input to this function is sliced so that vec[0] is the first harvested value for this transform
Returns:n – number of rows read from vec. Used to increment vec slice for next transform
Return type:int
regularization(regdict)

regularization vector

Parameters:regdict (dict) – bigfeta.schemas.regularization. controls regularization values
Returns:reg – array of regularization values of length DOF_per_tile
Return type:numpy.ndarray
to_solve_vec()

sets solve vector values from transform parameters

Returns:vec – N x 2 transform parameters in solve form
Return type:numpy.ndarray
class bigfeta.transform.rotation_model.AlignerRotationModel(transform=None)

Bases: renderapi.transform.leaf.affine_models.AffineModel

Object for implementing rotation transform

__init__(transform=None)
Parameters:transform (renderapi.transform.Transform) – The new AlignerTransform will inherit from this transform, if possible.
block_from_pts(pts, w, col_ind, col_max)
partial sparse block for a transform/match.
Note: for rotation, a pre-processing step is called at the tilepair level.
Parameters:
  • pts (numpy.ndarray) – N x 1, preprocessed from preprocess()
  • w (numpy.ndarray) – the weights associated with the pts
  • col_ind (int) – the starting column index for this tile
  • col_max (int) – number of columns in the matrix
Returns:

from_solve_vec(vec)

reads values from solution and sets transform parameters

Parameters:vec (numpy.ndarray) – input to this function is sliced so that vec[0] is the first harvested value for this transform
Returns:n – number of rows read from vec. Used to increment vec slice for next transform
Return type:int
static preprocess(ppts, qpts, w, npts_max=None, choose_random=True, cutoff_distance=15)
tilepair-level preprocessing step for rotation transform.
derives the relative center-of-mass angles between all p’s and q’s to avoid angular discontinuity. Will filter out points very close to center-of-mass. Tilepairs with relative rotations near 180deg will not avoid the discontinuity.
Parameters:
  • ppts (numpy.ndarray) – N x 2. The p tile correspondence coordinates
  • qpts (numpy.ndarray) – N x 2. The q tile correspondence coordinates
  • w (numpy.ndarray) – size N. The weights.
  • npts_max (int) – maximum points to include after processing
  • choose_random (bool) – whether to randomly reduce the input points to npts_max
  • cutoff_distance (float) – distance from center of mass of each point cloud in which to exclude points
Returns:

  • pa (numpy.ndarray) – M x 1 preprocessed angular distances. -0.5 x delta angle M <= N depending on filter
  • qa (numpy.ndarray) – M x 1 preprocessed angular distances. 0.5 x delta angle M <= N depending on filter
  • w (numpy.ndarray) – size M. filtered weights.

regularization(regdict)

regularization vector

Parameters:regdict (dict) – bigfeta.schemas.regularization. controls regularization values
Returns:reg – array of regularization values of length DOF_per_tile
Return type:numpy.ndarray
to_solve_vec()

sets solve vector values from transform parameters

Returns:vec – N x 1 transform parameters in solve form
Return type:numpy.ndarray
class bigfeta.transform.similarity_model.AlignerSimilarityModel(transform=None)

Bases: renderapi.transform.leaf.affine_models.AffineModel

Object for implementing similarity transform.

__init__(transform=None)
Parameters:transform (renderapi.transform.Transform) – The new AlignerTransform will inherit from this transform, if possible.
block_from_pts(pts, w, col_ind, col_max)
partial sparse block for a transform/match.
similarity constrains the center-of-mass coordinates to transform according to the same affine transform as the coordinates, save translation.
Parameters:
  • pts (numpy.ndarray) – N x 2, the x, y values of the match (either p or q)
  • w (numpy.ndarray) – the weights associated with the pts
  • col_ind (int) – the starting column index for this tile
  • col_max (int) – number of columns in the matrix
Returns:

  • block (scipy.sparse.csr_matrix) – the partial block for this transform
  • w (numpy.ndarray) – the weights associated with the rows of this block
  • rhs (numpy.ndarray) – N x 1 (fullsize) right hand side for this transform. generally all zeros. could implement fixed tiles in rhs later.

from_solve_vec(vec)

reads values from solution and sets transform parameters

Parameters:vec (numpy.ndarray) – input to this function is sliced so that vec[0] is the first harvested value for this transform
Returns:n – number of rows read from vec. Used to increment vec slice for next transform
Return type:int
regularization(regdict)

regularization vector

Parameters:regdict (dict) – bigfeta.schemas.regularization. controls regularization values
Returns:reg – array of regularization values of length DOF_per_tile
Return type:numpy.ndarray
to_solve_vec()

sets solve vector values from transform parameters

Returns:vec – N x 1 transform parameters in solve form
Return type:numpy.ndarray
class bigfeta.transform.thinplatespline_model.AlignerThinPlateSplineTransform(transform=None)

Bases: renderapi.transform.leaf.thin_plate_spline.ThinPlateSplineTransform

Object for implementing thin plate spline transform

__init__(transform=None)
Parameters:transform (renderapi.transform.Transform) – The new AlignerTransform will inherit from this transform, if possible.
block_from_pts(pts, w, col_ind, col_max)

partial sparse block for a tilepair/match

Parameters:
  • pts (numpy.ndarray) – N x 2, the x, y values of the match (either p or q)
  • w (numpy.ndarray) – the weights associated with the pts
  • col_ind (int) – the starting column index for this tile
  • col_max (int) – number of columns in the matrix
Returns:

from_solve_vec(vec)

reads values from solution and sets transform parameters

Parameters:vec (numpy.ndarray) – input to this function is sliced so that vec[0] is the first harvested value for this transform
Returns:n – number of rows read from vec. Used to increment vec slice for next transform
Return type:int
regularization(regdict)

regularization vector

Parameters:regdict (dict) – bigfeta.schemas.regularization. controls regularization values
Returns:reg – array of regularization values of length DOF_per_tile
Return type:numpy.ndarray
scale

tuple of scale for x, y. For setting regularization, it is useful to watch scale (logged output for the solver) to look for unwanted distortions and shrinking. Other transforms have scale implemented inside of renderapi.

to_solve_vec()

sets solve vector values from transform parameters

Returns:vec – N x 2 transform parameters in solve form
Return type:numpy.ndarray
class bigfeta.transform.translation_model.AlignerTranslationModel(transform=None)

Bases: renderapi.transform.leaf.affine_models.AffineModel

Object for implementing translation transform

__init__(transform=None)
Parameters:transform (renderapi.transform.Transform) – The new AlignerTransform will inherit from this transform, if possible.
block_from_pts(pts, w, col_ind, col_max)

partial sparse block for a tilepair/match

Parameters:
  • pts (numpy.ndarray) – N x 2, the x, y values of the match (either p or q)
  • w (numpy.ndarray) – the weights associated with the pts
  • col_ind (int) – the starting column index for this tile
  • col_max (int) – number of columns in the matrix
Returns:

from_solve_vec(vec)

reads values from solution and sets transform parameters

Parameters:vec (numpy.ndarray) – input to this function is sliced so that vec[0] is the first harvested value for this transform
Returns:n – number of rows read from vec. Used to increment vec slice for next transform
Return type:int
regularization(regdict)

regularization vector

Parameters:regdict (dict) – bigfeta.schemas.regularization. controls regularization values
Returns:reg – array of regularization values of length DOF_per_tile
Return type:numpy.ndarray
to_solve_vec()

sets solve vector values from transform parameters

Returns:vec – 1 x 2 transform parameters in solve form
Return type:numpy.ndarray
exception bigfeta.transform.utils.AlignerTransformException

Bases: Exception

Exception class for AlignerTransforms

bigfeta.transform.utils.aff_matrices(thetas, offs=None)

affine matrices from thetas

bigfeta.transform.utils.aff_matrix(theta, offs=None)

affine matrix or augmented affine matrix given a rotation angle.

Parameters:
  • theta (float) – rotation angle in radians
  • offs (numpy.ndarray) – the translations to include
Returns:

M – 2 x 2 (for offs=None) affine matrix or 3 x 3 augmented matrix

Return type:

numpy.ndarray