00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011 #ifndef EIGEN_HYPERPLANE_H
00012 #define EIGEN_HYPERPLANE_H
00013
00014 namespace Eigen {
00015
00033 template <typename _Scalar, int _AmbientDim, int _Options>
00034 class Hyperplane
00035 {
00036 public:
00037 EIGEN_MAKE_ALIGNED_OPERATOR_NEW_IF_VECTORIZABLE_FIXED_SIZE(_Scalar,_AmbientDim==Dynamic ? Dynamic : _AmbientDim+1)
00038 enum {
00039 AmbientDimAtCompileTime = _AmbientDim,
00040 Options = _Options
00041 };
00042 typedef _Scalar Scalar;
00043 typedef typename NumTraits<Scalar>::Real RealScalar;
00044 typedef DenseIndex Index;
00045 typedef Matrix<Scalar,AmbientDimAtCompileTime,1> VectorType;
00046 typedef Matrix<Scalar,Index(AmbientDimAtCompileTime)==Dynamic
00047 ? Dynamic
00048 : Index(AmbientDimAtCompileTime)+1,1,Options> Coefficients;
00049 typedef Block<Coefficients,AmbientDimAtCompileTime,1> NormalReturnType;
00050 typedef const Block<const Coefficients,AmbientDimAtCompileTime,1> ConstNormalReturnType;
00051
00053 inline explicit Hyperplane() {}
00054
00055 template<int OtherOptions>
00056 Hyperplane(const Hyperplane<Scalar,AmbientDimAtCompileTime,OtherOptions>& other)
00057 : m_coeffs(other.coeffs())
00058 {}
00059
00062 inline explicit Hyperplane(Index _dim) : m_coeffs(_dim+1) {}
00063
00067 inline Hyperplane(const VectorType& n, const VectorType& e)
00068 : m_coeffs(n.size()+1)
00069 {
00070 normal() = n;
00071 offset() = -n.dot(e);
00072 }
00073
00078 inline Hyperplane(const VectorType& n, Scalar d)
00079 : m_coeffs(n.size()+1)
00080 {
00081 normal() = n;
00082 offset() = d;
00083 }
00084
00088 static inline Hyperplane Through(const VectorType& p0, const VectorType& p1)
00089 {
00090 Hyperplane result(p0.size());
00091 result.normal() = (p1 - p0).unitOrthogonal();
00092 result.offset() = -p0.dot(result.normal());
00093 return result;
00094 }
00095
00099 static inline Hyperplane Through(const VectorType& p0, const VectorType& p1, const VectorType& p2)
00100 {
00101 EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(VectorType, 3)
00102 Hyperplane result(p0.size());
00103 result.normal() = (p2 - p0).cross(p1 - p0).normalized();
00104 result.offset() = -p0.dot(result.normal());
00105 return result;
00106 }
00107
00112
00113 explicit Hyperplane(const ParametrizedLine<Scalar, AmbientDimAtCompileTime>& parametrized)
00114 {
00115 normal() = parametrized.direction().unitOrthogonal();
00116 offset() = -parametrized.origin().dot(normal());
00117 }
00118
00119 ~Hyperplane() {}
00120
00122 inline Index dim() const { return AmbientDimAtCompileTime==Dynamic ? m_coeffs.size()-1 : Index(AmbientDimAtCompileTime); }
00123
00125 void normalize(void)
00126 {
00127 m_coeffs /= normal().norm();
00128 }
00129
00133 inline Scalar signedDistance(const VectorType& p) const { return normal().dot(p) + offset(); }
00134
00138 inline Scalar absDistance(const VectorType& p) const { return internal::abs(signedDistance(p)); }
00139
00142 inline VectorType projection(const VectorType& p) const { return p - signedDistance(p) * normal(); }
00143
00147 inline ConstNormalReturnType normal() const { return ConstNormalReturnType(m_coeffs,0,0,dim(),1); }
00148
00152 inline NormalReturnType normal() { return NormalReturnType(m_coeffs,0,0,dim(),1); }
00153
00157 inline const Scalar& offset() const { return m_coeffs.coeff(dim()); }
00158
00161 inline Scalar& offset() { return m_coeffs(dim()); }
00162
00166 inline const Coefficients& coeffs() const { return m_coeffs; }
00167
00171 inline Coefficients& coeffs() { return m_coeffs; }
00172
00179 VectorType intersection(const Hyperplane& other) const
00180 {
00181 EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(VectorType, 2)
00182 Scalar det = coeffs().coeff(0) * other.coeffs().coeff(1) - coeffs().coeff(1) * other.coeffs().coeff(0);
00183
00184
00185 if(internal::isMuchSmallerThan(det, Scalar(1)))
00186 {
00187 if(internal::abs(coeffs().coeff(1))>internal::abs(coeffs().coeff(0)))
00188 return VectorType(coeffs().coeff(1), -coeffs().coeff(2)/coeffs().coeff(1)-coeffs().coeff(0));
00189 else
00190 return VectorType(-coeffs().coeff(2)/coeffs().coeff(0)-coeffs().coeff(1), coeffs().coeff(0));
00191 }
00192 else
00193 {
00194 Scalar invdet = Scalar(1) / det;
00195 return VectorType(invdet*(coeffs().coeff(1)*other.coeffs().coeff(2)-other.coeffs().coeff(1)*coeffs().coeff(2)),
00196 invdet*(other.coeffs().coeff(0)*coeffs().coeff(2)-coeffs().coeff(0)*other.coeffs().coeff(2)));
00197 }
00198 }
00199
00206 template<typename XprType>
00207 inline Hyperplane& transform(const MatrixBase<XprType>& mat, TransformTraits traits = Affine)
00208 {
00209 if (traits==Affine)
00210 normal() = mat.inverse().transpose() * normal();
00211 else if (traits==Isometry)
00212 normal() = mat * normal();
00213 else
00214 {
00215 eigen_assert(0 && "invalid traits value in Hyperplane::transform()");
00216 }
00217 return *this;
00218 }
00219
00227 template<int TrOptions>
00228 inline Hyperplane& transform(const Transform<Scalar,AmbientDimAtCompileTime,Affine,TrOptions>& t,
00229 TransformTraits traits = Affine)
00230 {
00231 transform(t.linear(), traits);
00232 offset() -= normal().dot(t.translation());
00233 return *this;
00234 }
00235
00241 template<typename NewScalarType>
00242 inline typename internal::cast_return_type<Hyperplane,
00243 Hyperplane<NewScalarType,AmbientDimAtCompileTime,Options> >::type cast() const
00244 {
00245 return typename internal::cast_return_type<Hyperplane,
00246 Hyperplane<NewScalarType,AmbientDimAtCompileTime,Options> >::type(*this);
00247 }
00248
00250 template<typename OtherScalarType,int OtherOptions>
00251 inline explicit Hyperplane(const Hyperplane<OtherScalarType,AmbientDimAtCompileTime,OtherOptions>& other)
00252 { m_coeffs = other.coeffs().template cast<Scalar>(); }
00253
00258 template<int OtherOptions>
00259 bool isApprox(const Hyperplane<Scalar,AmbientDimAtCompileTime,OtherOptions>& other, typename NumTraits<Scalar>::Real prec = NumTraits<Scalar>::dummy_precision()) const
00260 { return m_coeffs.isApprox(other.m_coeffs, prec); }
00261
00262 protected:
00263
00264 Coefficients m_coeffs;
00265 };
00266
00267 }
00268
00269 #endif // EIGEN_HYPERPLANE_H