mirror of
https://gitlab.com/libeigen/eigen.git
synced 2025-07-09 06:31:47 +08:00
backporting bugfix in Quaternion::setFromTwoVectors()
This commit is contained in:
parent
9bff5e4f67
commit
c6eb9ef60e
@ -346,7 +346,6 @@ inline Quaternion<Scalar>& Quaternion<Scalar>::setFromTwoVectors(const MatrixBas
|
|||||||
{
|
{
|
||||||
Vector3 v0 = a.normalized();
|
Vector3 v0 = a.normalized();
|
||||||
Vector3 v1 = b.normalized();
|
Vector3 v1 = b.normalized();
|
||||||
Vector3 axis = v0.cross(v1);
|
|
||||||
Scalar c = v0.dot(v1);
|
Scalar c = v0.dot(v1);
|
||||||
|
|
||||||
// if dot == 1, vectors are the same
|
// if dot == 1, vectors are the same
|
||||||
@ -354,7 +353,17 @@ inline Quaternion<Scalar>& Quaternion<Scalar>::setFromTwoVectors(const MatrixBas
|
|||||||
{
|
{
|
||||||
// set to identity
|
// set to identity
|
||||||
this->w() = 1; this->vec().setZero();
|
this->w() = 1; this->vec().setZero();
|
||||||
|
return *this;
|
||||||
}
|
}
|
||||||
|
// if dot == -1, vectors are opposites
|
||||||
|
if (ei_isApprox(c,Scalar(-1)))
|
||||||
|
{
|
||||||
|
this->vec() = v0.unitOrthogonal();
|
||||||
|
this->w() = 0;
|
||||||
|
return *this;
|
||||||
|
}
|
||||||
|
|
||||||
|
Vector3 axis = v0.cross(v1);
|
||||||
Scalar s = ei_sqrt((Scalar(1)+c)*Scalar(2));
|
Scalar s = ei_sqrt((Scalar(1)+c)*Scalar(2));
|
||||||
Scalar invs = Scalar(1)/s;
|
Scalar invs = Scalar(1)/s;
|
||||||
this->vec() = axis * invs;
|
this->vec() = axis * invs;
|
||||||
|
Loading…
x
Reference in New Issue
Block a user