1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176
|
/*=========================================================================
*
* Copyright NumFOCUS
*
* Licensed under the Apache License, Version 2.0 (the "License");
* you may not use this file except in compliance with the License.
* You may obtain a copy of the License at
*
* https://www.apache.org/licenses/LICENSE-2.0.txt
*
* Unless required by applicable law or agreed to in writing, software
* distributed under the License is distributed on an "AS IS" BASIS,
* WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
* See the License for the specific language governing permissions and
* limitations under the License.
*
*=========================================================================*/
#ifndef itkEuler3DTransform_h
#define itkEuler3DTransform_h
#include <iostream>
#include "itkRigid3DTransform.h"
namespace itk
{
/** \class Euler3DTransform
*
* \brief Euler3DTransform of a vector space (e.g. space coordinates)
*
* This transform applies a rotation and translation to the space given 3 euler
* angles and a 3D translation. Rotation is about a user specified center.
*
* The parameters for this transform can be set either using individual Set
* methods or in serialized form using SetParameters() and SetFixedParameters().
*
* The serialization of the optimizable parameters is an array of 6 elements.
* The first 3 represents three euler angle of rotation respectively about
* the X, Y and Z axis. The last 3 parameters defines the translation in each
* dimension.
*
* The serialization of the fixed parameters is an array of 4 elements.
* The first 3 elements define the center of rotation. The final
* element indicates the status of ComputeZYX (i.e. 1.0 if true, 0.0 if false).
*
* \ingroup ITKTransform
*/
template <typename TParametersValueType = double>
class ITK_TEMPLATE_EXPORT Euler3DTransform : public Rigid3DTransform<TParametersValueType>
{
public:
ITK_DISALLOW_COPY_AND_MOVE(Euler3DTransform);
/** Standard class type aliases. */
using Self = Euler3DTransform;
using Superclass = Rigid3DTransform<TParametersValueType>;
using Pointer = SmartPointer<Self>;
using ConstPointer = SmartPointer<const Self>;
/** New macro for creation of through a Smart Pointer. */
itkNewMacro(Self);
/** \see LightObject::GetNameOfClass() */
itkOverrideGetNameOfClassMacro(Euler3DTransform);
/** Dimension of the space. */
static constexpr unsigned int SpaceDimension = 3;
static constexpr unsigned int InputSpaceDimension = 3;
static constexpr unsigned int OutputSpaceDimension = 3;
static constexpr unsigned int ParametersDimension = 6;
using typename Superclass::ParametersType;
using typename Superclass::ParametersValueType;
using typename Superclass::FixedParametersType;
using typename Superclass::FixedParametersValueType;
using typename Superclass::JacobianType;
using typename Superclass::JacobianPositionType;
using typename Superclass::InverseJacobianPositionType;
using typename Superclass::ScalarType;
using typename Superclass::InputVectorType;
using typename Superclass::OutputVectorType;
using typename Superclass::InputCovariantVectorType;
using typename Superclass::OutputCovariantVectorType;
using typename Superclass::InputVnlVectorType;
using typename Superclass::OutputVnlVectorType;
using typename Superclass::InputPointType;
using typename Superclass::OutputPointType;
using typename Superclass::MatrixType;
using typename Superclass::InverseMatrixType;
using typename Superclass::CenterType;
using typename Superclass::TranslationType;
using typename Superclass::OffsetType;
using AngleType = typename Superclass::ScalarType;
/** Set/Get the transformation from a container of parameters
* This is typically used by optimizers. There are 6 parameters. The first
* three represent the angles to rotate around the coordinate axis, and the
* last three represents the offset. */
void
SetParameters(const ParametersType & parameters) override;
const ParametersType &
GetParameters() const override;
const FixedParametersType &
GetFixedParameters() const override;
void
SetFixedParameters(const FixedParametersType & parameters) override;
/** Set the rotational part of the transform. */
void
SetRotation(ScalarType angleX, ScalarType angleY, ScalarType angleZ);
itkGetConstMacro(AngleX, ScalarType);
itkGetConstMacro(AngleY, ScalarType);
itkGetConstMacro(AngleZ, ScalarType);
/** This method computes the Jacobian matrix of the transformation.
* given point or vector, returning the transformed point or
* vector. The rank of the Jacobian will also indicate if the
* transform is invertible at this point. */
void
ComputeJacobianWithRespectToParameters(const InputPointType & p, JacobianType & jacobian) const override;
using Superclass::ComputeJacobianWithRespectToPosition;
/** The Euler angle representation of a rotation is not unique and
* depends on the order of rotations. In general there are 12
* options. This class supports two of them, ZXY and ZYX. The
* default is ZXY. These functions set and get the value which
* indicates whether the rotation is ZYX or ZXY.
*/
virtual void
SetComputeZYX(const bool flag);
itkGetConstMacro(ComputeZYX, bool);
/** Set the state to the identity.
*
* Sets the angles to a 0 value.
*/
void
SetIdentity() override;
protected:
Euler3DTransform(const MatrixType & matrix, const OutputPointType & offset);
Euler3DTransform(unsigned int parametersDimension);
Euler3DTransform();
~Euler3DTransform() override = default;
void
PrintSelf(std::ostream & os, Indent indent) const override;
/** Set values of angles directly without recomputing other parameters. */
void
SetVarRotation(ScalarType angleX, ScalarType angleY, ScalarType angleZ);
/** Compute the components of the rotation matrix in the superclass. */
void
ComputeMatrix() override;
/** Compute angles from the rotation matrix. */
void
ComputeMatrixParameters() override;
private:
ScalarType m_AngleX{};
ScalarType m_AngleY{};
ScalarType m_AngleZ{};
bool m_ComputeZYX{};
}; // class Euler3DTransform
} // namespace itk
#ifndef ITK_MANUAL_INSTANTIATION
# include "itkEuler3DTransform.hxx"
#endif
#endif /* itkEuler3DTransform_h */
|