34 #include <boost/assign/std/vector.hpp>
35 #include <boost/numeric/ublas/io.hpp>
36 #include <boost/numeric/ublas/matrix.hpp>
37 #include <boost/numeric/ublas/matrix_proxy.hpp>
74 if (matrix.size1() != matrix.size2())
return false;
76 for (
unsigned int row = 0;
row < matrix.size1();
row++) {
77 for (
unsigned int col = 0; col < matrix.size2(); col++) {
79 if (!qFuzzyCompare(matrix(
row, col), 1.0))
return false;
84 if (!qFuzzyCompare(1.0 + matrix(
row, col), 1.0))
return false;
102 if (matrix.size1() != matrix.size2())
return false;
127 if (matrix.size1() != matrix.size2())
return false;
183 Matrix unitMatrix(matrix.size1(), matrix.size2());
185 for (
unsigned int i = 0; i < unitMatrix.size1(); i++) {
189 double columnMagnitude =
magnitude(columnVector);
190 if ( columnMagnitude < 0.9 || columnMagnitude > 1.1 ) {
195 setColumn(unitMatrix, columnVector / columnMagnitude, i);
203 if ( det >= 0.9 && det <= 1.1 ) {
209 QString msg =
"Unable to determine whether the given matrix is a rotation matrix.";
224 for (
unsigned int i = 0; i < matrix.size1(); i++) {
225 for (
unsigned int j = 0; j < matrix.size2(); j++) {
228 if (!qFuzzyCompare(matrix(i, j) + 1.0, 1.0))
return false;
244 for (
unsigned int i = 0; i < vector.size(); i++) {
247 if (!qFuzzyCompare(
vector(i) + 1.0, 1.0))
return false;
263 return vector.empty();
275 if (qFuzzyCompare(
dotProduct(vector, vector), 1.0)) {
303 if (qFuzzyCompare(det + 1.0, 1.0)) {
305 QString msg =
"The given matrix is not invertible. The determinant is 0.0.";
310 double scale = 1 / det;
314 if (matrix.size1() == 2) {
315 inverse(0, 0) = scale * matrix(1, 1);
316 inverse(0, 1) = -scale * matrix(0, 1);
318 inverse(1, 0) = -scale * matrix(1, 0);
319 inverse(1, 1) = scale * matrix(0, 0);
323 inverse(0, 0) = scale * ( matrix(1, 1) * matrix(2, 2) - matrix(2, 1) * matrix(1, 2) );
324 inverse(0, 1) = scale * ( matrix(0, 2) * matrix(2, 1) - matrix(2, 2) * matrix(0, 1) );
325 inverse(0, 2) = scale * ( matrix(0, 1) * matrix(1, 2) - matrix(1, 1) * matrix(0, 2) );
327 inverse(1, 0) = scale * ( matrix(1, 2) * matrix(2, 0) - matrix(2, 2) * matrix(1, 0) );
328 inverse(1, 1) = scale * ( matrix(0, 0) * matrix(2, 2) - matrix(2, 0) * matrix(0, 2) );
329 inverse(1, 2) = scale * ( matrix(0, 2) * matrix(0, 1) - matrix(1, 2) * matrix(0, 0) );
331 inverse(2, 0) = scale * ( matrix(1, 0) * matrix(2, 1) - matrix(2, 0) * matrix(1, 1) );
332 inverse(2, 1) = scale * ( matrix(0, 1) * matrix(2, 0) - matrix(2, 1) * matrix(0, 0) );
333 inverse(2, 2) = scale * ( matrix(0, 0) * matrix(1, 1) - matrix(1, 0) * matrix(0, 1) );
338 QString msg =
"Unable to invert the given matrix.";
353 return boost::numeric::ublas::trans(matrix);
368 QString msg =
"Can not create identity matrix of negative size ["
374 for (
int i = 0; i < size; i++) {
391 boost::numeric::ublas::zero_matrix<double> m(rows, columns);
404 boost::numeric::ublas::zero_vector<double> v(size);
423 if ( (matrix.size1() != matrix.size2()) || (matrix.size1() != 2 && matrix.size1() != 3) ) {
424 QString msg =
"Unable to calculate the determinant for the given matrix. "
425 "This method only calculates the determinant for 2x2 or 3x3 matrices."
426 "The given matrix is [" +
toString((
int) matrix.size1()) +
"x"
427 +
toString((
int) matrix.size2()) +
"].";
431 if (matrix.size1() == 2) {
432 return matrix(0, 0) * matrix(1, 1) - matrix(1, 0) * matrix(0, 1);
436 return matrix(0, 0) * (matrix(1, 1)*matrix(2, 2) - matrix(1, 2)*matrix(2, 1))
437 - matrix(0, 1) * (matrix(1, 0)*matrix(2, 2) - matrix(1, 2)*matrix(2, 0))
438 + matrix(0, 2) * (matrix(1, 0)*matrix(2, 1) - matrix(1, 1)*matrix(2, 0));
462 QString msg =
"Unable to normalize the zero vector.";
495 magnitude = maxNorm * boost::numeric::ublas::norm_2(vector / maxNorm);
514 return boost::numeric::ublas::norm_inf(vector);
533 if (matrix1.size2() != matrix2.size1()) {
534 QString msg =
"Unable to multiply matrices with mismatched dimensions. "
535 "The left matrix has [" +
toString((
int) matrix1.size2())
536 +
"] columns and the right matrix has [" +
toString((
int) matrix2.size1())
543 return boost::numeric::ublas::prod(matrix1, matrix2);
564 if (matrix.size2() != vector.size()) {
565 QString msg =
"Unable to multiply matrix and vector with mismatched dimensions."
566 "The given vector has [" +
toString((
int) vector.size())
567 +
"] components and the given matrix has ["
568 +
toString((
int) matrix.size2()) +
"] columns.";
573 return boost::numeric::ublas::prod(matrix, vector);
601 return scalar * matrix;
618 if (vector1.size() != vector2.size()) {
619 QString msg =
"Unable to add vectors with mismatched sizes."
620 "Vector1 has [" +
toString((
int) vector1.size())
621 +
"] components and vector2 has ["
622 +
toString((
int) vector2.size()) +
"] components.";
626 return vector1 + vector2;
644 if (vector1.size() != vector2.size()) {
645 QString msg =
"Unable to subtract vectors with mismatched sizes."
646 "Vector1 has [" +
toString((
int) vector1.size())
647 +
"] components and vector2 has ["
648 +
toString((
int) vector2.size()) +
"] components.";
652 return vector1 - vector2;
671 if ((vector1.size() != 3) || (vector2.size() != 3)) {
672 QString msg =
"Unable to calculate the cross product on vectors that are not size 3. "
673 "Vector1 has [" +
toString((
int) vector1.size())
674 +
"] components and vector2 has ["
675 +
toString((
int) vector2.size()) +
"] components.";
680 crossProd(0) = vector1(1)*vector2(2) - vector1(2)*vector2(1);
681 crossProd(1) = vector1(2)*vector2(0) - vector1(0)*vector2(2);
682 crossProd(2) = vector1(0)*vector2(1) - vector1(1)*vector2(0);
712 if (maxVector1 != 0.0) {
713 normalizedVector1 /= maxVector1;
716 if (maxVector2 != 0.0) {
717 normalizedVector2 /= maxVector2;
740 if (vector1.size() != vector2.size()) {
741 QString msg =
"Unable to compute the outer product for vectors with mismatched sizes."
742 "Vector1 has [" +
toString((
int) vector1.size())
743 +
"] components and vector2 has ["
744 +
toString((
int) vector2.size()) +
"] components.";
747 return boost::numeric::ublas::outer_prod(vector1, vector2);
784 if (vector1.size() != vector2.size()) {
785 QString msg =
"Unable to compute the dot product for vectors with mismatched sizes."
786 "Vector1 has [" +
toString((
int) vector1.size())
787 +
"] components and vector2 has ["
788 +
toString((
int) vector2.size()) +
"] components.";
791 return boost::numeric::ublas::inner_prod(vector1, vector2);
816 if (vector1.size() != vector2.size()) {
817 QString msg =
"Unable to project vector1 onto vector2 with mismatched sizes."
818 "Vector1 has [" +
toString((
int) vector1.size())
819 +
"] components and vector2 has ["
820 +
toString((
int) vector2.size()) +
"] components.";
834 return v1Dotv2/v2Dotv2 * vector2;
856 if ((vector.size() != 3) || (axis.size() != 3)) {
857 QString msg =
"Unable to rotate vector about the given axis and angle. "
858 "Vectors must be of size 3 to perform rotation. "
859 "The given vector has [" +
toString((
int) vector.size())
860 +
"] components and the given axis has ["
861 +
toString((
int) axis.size()) +
"] components.";
876 Vector v1 = vector - projVectorOnAxis;
883 double c = cos(angle.
radians());
884 double s = sin(angle.
radians());
885 Vector rplane = (c*v1) + (s*v2);
889 return (rplane + projVectorOnAxis);
929 Vector parallelVector =
project(vector1 / max1, vector2 / max2);
931 Vector perpendicularVector = vector1 - parallelVector * max1;
933 return perpendicularVector;
957 if ((rotationMatrix.size1() != 3) || (rotationMatrix.size2() != 3)) {
958 QString msg =
"Unable to convert the given matrix to an axis of rotation "
959 "and a rotation angle. A 3x3 matrix is required. The given matrix is ["
960 +
toString((
int) rotationMatrix.size1()) +
"x"
961 +
toString((
int) rotationMatrix.size2()) +
"].";
967 QString msg =
"Unable to convert the given matrix to an axis of rotation "
968 "and a rotation angle. The given matrix is not a rotation matrix.";
979 subQuaternion[0] = quaternion[1];
980 subQuaternion[1] = quaternion[2];
981 subQuaternion[2] = quaternion[3];
990 if (
isZero(subQuaternion)) {
997 else if (qFuzzyCompare(quaternion(0) + 1.0, 1.0)) {
998 axis = subQuaternion;
1006 return qMakePair(axis, angle);
1037 if (axes.size() != 3) {
1038 QString msg =
"Unable to convert the given matrix to Euler angles. "
1039 "Exactly 3 axis codes are required. The given list has ["
1040 +
toString((
int) axes.size()) +
"] axes.";
1045 validAxes << 1 << 2 << 3;
1046 if (!validAxes.contains(axes[0])
1047 || !validAxes.contains(axes[1])
1048 || !validAxes.contains(axes[2])) {
1049 QString msg =
"Unable to convert the given matrix to Euler angles using the given axis codes "
1051 +
"]. Axis codes must be 1, 2, or 3.";
1055 if (axes[0] == axes[1] || axes[1] == axes[2]) {
1056 QString msg =
"Unable to convert the given matrix to Euler angles using the given axis codes "
1058 +
"]. The middle axis code must differ from its neighbors.";
1063 if ((rotationMatrix.size1() != 3) || (rotationMatrix.size2() != 3)) {
1064 QString msg =
"Unable to convert the given matrix to Euler angles. A 3x3 matrix is required. "
1065 "The given matrix is ["
1066 +
toString((
int) rotationMatrix.size1()) +
"x"
1067 +
toString((
int) rotationMatrix.size2()) +
"].";
1072 QString msg =
"Unable to convert the given matrix to Euler angles. "
1073 "The given matrix is not a rotation matrix.";
1085 for (
int i = 0; i < 3; i++) {
1092 nextAxis << 2 << 3 << 1;
1097 Angle angle1, angle2, angle3;
1098 if ( axes[0] == axes[2] ) {
1099 if ( axes[1] == nextAxis[axes[0] - 1] ) {
1105 int c = 6 - axes[0] - axes[1];
1107 change( axes[0] - 1, 2 ) = 1.0;
1108 change( axes[1] - 1, 0 ) = 1.0;
1109 change( c - 1, 1 ) = sign * 1.0;
1113 bool degen = ( qFuzzyCompare(tempRotation(0, 2) + 1.0, 1.0)
1114 && qFuzzyCompare(tempRotation(1, 2) + 1.0, 1.0) )
1115 || ( qFuzzyCompare(tempRotation(2, 0) + 1.0, 1.0)
1116 && qFuzzyCompare(tempRotation(2, 1) + 1.0, 1.0) )
1117 || qFuzzyCompare( qAbs(tempRotation(2, 2)), 1.0 );
1122 angle2.
setRadians( acos( tempRotation(2, 2) ) );
1123 angle1.
setRadians( atan2( tempRotation(0, 1), tempRotation(0, 0) ) );
1128 angle3.
setRadians( atan2( tempRotation(0, 2), tempRotation(1, 2) ) );
1129 angle2.
setRadians( acos( tempRotation(2, 2) ) );
1130 angle1.
setRadians( atan2( tempRotation(2, 0), -tempRotation(2, 1) ) );
1139 if ( axes[1] == nextAxis[axes[0] - 1] ) {
1145 change( axes[0] - 1, 0 ) = 1.0;
1146 change( axes[1] - 1, 1 ) = 1.0;
1147 change( axes[2] - 1, 2 ) = sign * 1.0;
1150 bool degen = ( qFuzzyCompare(tempRotation(0, 0) + 1.0, 1.0)
1151 && qFuzzyCompare(tempRotation(0, 1) + 1.0, 1.0) )
1152 || ( qFuzzyCompare(tempRotation(1, 2) + 1.0, 1.0)
1153 && qFuzzyCompare(tempRotation(2, 2) + 1.0, 1.0) )
1154 || qFuzzyCompare( qAbs(tempRotation(0, 2)), 1.0 );
1159 angle2.
setRadians( asin( -tempRotation(0, 2) ) );
1160 angle1.
setRadians( sign * atan2( -tempRotation(1, 0), tempRotation(1, 1) ) );
1165 angle3.
setRadians( atan2( tempRotation(1, 2), tempRotation(2, 2) ) );
1166 angle2.
setRadians( asin( -tempRotation(0, 2) ) );
1167 angle1.
setRadians( sign * atan2( tempRotation(0, 1), tempRotation(0, 0) ) );
1172 eulerAngles << qMakePair(angle3, axes[0])
1173 << qMakePair(angle2, axes[1])
1174 << qMakePair(angle1, axes[2]);
1195 if ((rotationMatrix.size1() != 3) || (rotationMatrix.size2() != 3)) {
1196 QString msg =
"Unable to convert the given matrix to a quaternion. "
1197 "A 3x3 matrix is required. The given matrix is ["
1198 +
toString((
int) rotationMatrix.size1()) +
"x"
1199 +
toString((
int) rotationMatrix.size2()) +
"].";
1205 QString msg =
"Unable to convert the given matrix to an axis of rotation "
1206 "and a rotation angle. The given matrix is not a rotation matrix.";
1301 double trace = rotationMatrix(0, 0) + rotationMatrix(1, 1) + rotationMatrix(2, 2);
1302 double mtrace = 1.0 - trace;
1305 double cc4 = 1.0 + trace;
1307 double s114 = mtrace + 2.0*rotationMatrix(0, 0);
1308 double s224 = mtrace + 2.0*rotationMatrix(1, 1);
1309 double s334 = mtrace + 2.0*rotationMatrix(2, 2);
1313 double normalizingFactor;
1317 quaternion[0] = sqrt( cc4 * 0.25 );
1318 normalizingFactor = 1.0 / ( quaternion[0] * 4.0 );
1323 quaternion[1] = ( rotationMatrix(2, 1) - rotationMatrix(1, 2) ) * normalizingFactor;
1324 quaternion[2] = ( rotationMatrix(0, 2) - rotationMatrix(2, 0) ) * normalizingFactor;
1325 quaternion[3] = ( rotationMatrix(1, 0) - rotationMatrix(0, 1) ) * normalizingFactor;
1328 else if ( s114 >= 1.0 ) {
1330 quaternion[1] = sqrt( s114 * 0.25 );
1331 normalizingFactor = 1.0 /( quaternion[1] * 4.0 );
1333 quaternion[0] = ( rotationMatrix(2, 1) - rotationMatrix(1, 2) ) * normalizingFactor;
1334 quaternion[2] = ( rotationMatrix(0, 1) + rotationMatrix(1, 0) ) * normalizingFactor;
1335 quaternion[3] = ( rotationMatrix(0, 2) + rotationMatrix(2, 0) ) * normalizingFactor;
1338 else if ( s224 >= 1.0 ) {
1340 quaternion[2] = sqrt( s224 * 0.25 );
1341 normalizingFactor = 1.0 /( quaternion[2] * 4.0 );
1343 quaternion[0] = ( rotationMatrix(0, 2) - rotationMatrix(2, 0) ) * normalizingFactor;
1344 quaternion[1] = ( rotationMatrix(0, 1) + rotationMatrix(1, 0) ) * normalizingFactor;
1345 quaternion[3] = ( rotationMatrix(1, 2) + rotationMatrix(2, 1) ) * normalizingFactor;
1349 quaternion[3] = sqrt( s334 * 0.25 );
1350 normalizingFactor = 1.0 /( quaternion[3] * 4.0 );
1352 quaternion[0] = ( rotationMatrix(1, 0) - rotationMatrix(0, 1) ) * normalizingFactor;
1353 quaternion[1] = ( rotationMatrix(0, 2) + rotationMatrix(2, 0) ) * normalizingFactor;
1354 quaternion[2] = ( rotationMatrix(1, 2) + rotationMatrix(2, 1) ) * normalizingFactor;
1359 if ( !
isUnit(quaternion) ) {
1364 if ( quaternion[0] < 0.0 ) {
1388 if (axis.size() != 3) {
1389 QString msg =
"Unable to convert the given vector and angle to a rotation matrix. "
1390 "The given vector with size [" +
toString((
int) axis.size())
1391 +
"] is not a 3D axis vector.";
1402 for (
unsigned int i = 0; i < axis.size(); i++) {
1411 return rotationMatrix;
1448 validAxes << 1 << 2 << 3;
1449 if (!validAxes.contains(angle3.second)
1450 || !validAxes.contains(angle2.second)
1451 || !validAxes.contains(angle1.second)) {
1452 QString msg =
"Unable to convert the given Euler angles to a matrix using the given axis "
1453 "codes [" +
toString(angle3.second) +
", " +
toString(angle2.second) +
", "
1454 +
toString(angle1.second) +
"]. Axis codes must be 1, 2, or 3.";
1465 double sinAngle = sin(angle1.first.radians());
1466 double cosAngle = cos(angle1.first.radians());
1468 double indexOffset = ((angle1.second % 3) + 3) % 3;
1470 indexSet << 2 << 0 << 1 << 2 << 0;
1471 int index1 = indexSet[indexOffset + 0];
1472 int index2 = indexSet[indexOffset + 1];
1473 int index3 = indexSet[indexOffset + 2];
1474 m(index1, index1) = 1.0; m(index1, index2) = 0.0; m(index1, index3) = 0.0;
1475 m(index2, index1) = 0.0; m(index2, index2) = cosAngle; m(index2, index3) = sinAngle;
1476 m(index3, index1) = 0.0; m(index3, index2) = -sinAngle; m(index3, index3) = cosAngle;
1484 sinAngle = sin(angle2.first.radians());
1485 cosAngle = cos(angle2.first.radians());
1486 indexOffset = ((angle2.second % 3) + 3) % 3;
1487 index1 = indexSet[indexOffset + 0];
1488 index2 = indexSet[indexOffset + 1];
1489 index3 = indexSet[indexOffset + 2];
1492 for (
int i = 0; i < 3; i++) {
1493 tempMatrix(index1, i) = m(index1, i);
1494 tempMatrix(index2, i) = cosAngle*m(index2, i) + sinAngle*m(index3, i);
1495 tempMatrix(index3, i) = -sinAngle*m(index2, i) + cosAngle*m(index3, i);
1499 sinAngle = sin(angle3.first.radians());
1500 cosAngle = cos(angle3.first.radians());
1501 indexOffset = ((angle3.second % 3) + 3) % 3;
1502 index1 = indexSet[indexOffset + 0];
1503 index2 = indexSet[indexOffset + 1];
1504 index3 = indexSet[indexOffset + 2];
1505 for (
int i = 0; i < 3; i++) {
1506 m(index1, i) = tempMatrix(index1, i);
1507 m(index2, i) = cosAngle*tempMatrix(index2, i) + sinAngle*tempMatrix(index3, i);
1508 m(index3, i) = -sinAngle*tempMatrix(index2, i) + cosAngle*tempMatrix(index3, i);
1529 if (eulerAngles.size() != 3) {
1530 QString msg =
"Unable to convert the given Euler angles to a matrix. "
1531 "Exactly 3 Euler angles are required. The given list has ["
1532 +
toString((
int) eulerAngles.size()) +
"] angles.";
1555 if (quaternion.size() != 4) {
1556 QString msg =
"Unable to convert the given vector to a rotation matrix. "
1557 "The given vector with [" +
toString((
int) quaternion.size())
1558 +
"] components is not a quaternion.";
1569 double q0q1 = q[0] * q[1];
1570 double q0q2 = q[0] * q[2];
1571 double q0q3 = q[0] * q[3];
1573 double q1q1 = q[1] * q[1];
1574 double q1q2 = q[1] * q[2];
1575 double q1q3 = q[1] * q[3];
1577 double q2q2 = q[2] * q[2];
1578 double q2q3 = q[2] * q[3];
1580 double q3q3 = q[3] * q[3];
1582 matrix(0, 0) = 1.0 - 2.0 * ( q2q2 + q3q3 );
1583 matrix(0, 1) = 2.0 * ( q1q2 - q0q3 );
1584 matrix(0, 2) = 2.0 * ( q1q3 + q0q2 );
1586 matrix(1, 0) = 2.0 * ( q1q2 + q0q3 );
1587 matrix(1, 1) = 1.0 - 2.0 * ( q1q1 + q3q3 );
1588 matrix(1, 2) = 2.0 * ( q2q3 - q0q1 );
1590 matrix(2, 0) = 2.0 * ( q1q3 - q0q2 );
1591 matrix(2, 1) = 2.0 * ( q2q3 + q0q1 );
1592 matrix(2, 2) = 1.0 - 2.0 * ( q1q1 + q2q2 );
1609 if ( (rowIndex+1 > (
int) matrix.size1()) ) {
1610 QString msg =
"Unable to set the matrix row to the given vector. Row index "
1611 +
toString(rowIndex) +
" is out of bounds. The given matrix only has "
1612 +
toString((
int) matrix.size1()) +
" rows." ;
1615 boost::numeric::ublas::row(matrix, rowIndex) =
vector;
1630 if ( (columnIndex+1 > (
int) matrix.size2()) ) {
1631 QString msg =
"Unable to set the matrix column to the given vector. Column index "
1632 +
toString(columnIndex) +
" is out of bounds. The given matrix only has "
1633 +
toString((
int) matrix.size1()) +
" columns." ;
1636 boost::numeric::ublas::column(matrix, columnIndex) =
vector;
1653 if ( (rowIndex+1 > (
int) matrix.size1()) ) {
1654 QString msg =
"Unable to get the matrix row to the given vector. Row index "
1655 +
toString(rowIndex) +
" is out of bounds. The given matrix only has "
1656 +
toString((
int) matrix.size1()) +
" rows." ;
1659 return boost::numeric::ublas::row(matrix, rowIndex);
1676 if ( (columnIndex+1 > (
int) matrix.size2()) ) {
1677 QString msg =
"Unable to get the matrix column to the given vector. Column index "
1678 +
toString(columnIndex) +
" is out of bounds. The given matrix only has "
1679 +
toString((
int) matrix.size1()) +
" columns." ;
1682 return boost::numeric::ublas::column(matrix, columnIndex);
1773 for (
int i = 0; i < size; i++) {
1774 sub[i] = v[i + startIndex];
1793 QDebugStateSaver saver(dbg);
1811 QDebugStateSaver saver(dbg);
1812 for (
unsigned int i = 0; i < matrix.size1(); i++) {
1813 dbg.noquote() <<
" ";
1814 for (
unsigned int j = 0; j < matrix.size2(); j++) {
1815 dbg.noquote() <<
toString(matrix(i, j), 15) <<
" ";
1817 dbg.noquote() << endl;
1835 QString result =
"( ";
1836 for (
unsigned int i = 0; i < vector.size(); i++) {
1837 result +=
toString(vector(i), precision);
1838 if (i != vector.size() - 1) result +=
", ";
static Vector project(const Vector &vector1, const Vector &vector2)
Compute the vector projection of vector1 onto vector2.
static Vector vector(double v0, double v1, double v2)
Constructs a 3 dimensional vector with the given component values.
static double innerProduct(const Vector &vector1, const Vector &vector2)
Computes the inner product of the given vectors.
static Matrix transpose(const Matrix &matrix)
Returns the transpose of the given matrix.
static Vector subtract(const Vector &vector1, const Vector &vector2)
Subtracts the right vector from the left vector.
void setRadians(double radians)
Set the angle in units of Radians.
double radians() const
Convert an angle to a double.
static bool isRotationMatrix(const Matrix &matrix)
Determines whether the given matrix is a rotation matrix.
static bool isIdentity(const Matrix &matrix)
Determines whether the given matrix is the identity.
static Vector column(const Matrix &matrix, int columnIndex)
Returns a vector whose components match those of the given matrix column.
static Matrix multiply(const Matrix &matrix1, const Matrix &matrix2)
Returns the product of two matrices.
static double magnitude(const Vector &vector)
Computes the magnitude (i.e., the length) of the given vector using the Euclidean norm (L2 norm)...
const double PI(3.14159265358979323846)
The mathematical constant PI.
static Vector toQuaternion(const Matrix &rotationMatrix)
Converts a rotation's representation from a matrix to a unit quaternion.
static Vector rotate(const Vector &vector, const Vector &axis, Angle angle)
Rotates a vector about an axis vector given a specified angle.
static Vector subVector(const Vector &v, int start, int size)
Constructs a vector of given size using the given vector and starting index.
LinearAlgebra()
Default constructor for a LinearAlgebra object.
static void setColumn(Matrix &matrix, const Vector &vector, int columnIndex)
Sets the column of the given matrix to the values of the given vector.
boost::numeric::ublas::matrix< double > Matrix
Definition for an Isis::LinearAlgebra::Matrix of doubles.
QString toString(bool boolToConvert)
Global function to convert a boolean to a string.
static Vector normalize(const Vector &vector)
Returns a unit vector that is codirectional with the given vector by dividing each component of the v...
This error is for when a programmer made an API call that was illegal.
static Matrix zeroMatrix(int rows, int columns)
Returns a matrix with given dimensions that is filled with zeroes.
static double dotProduct(const Vector &vector1, const Vector &vector2)
Computes the dot product of the given vectors.
static Matrix outerProduct(const Vector &vector1, const Vector &vector2)
Computes the outer product of the given vectors.
static Matrix toMatrix(const AxisAngle &axisAngle)
Converts a rotation's representation from an axis of rotation and its corresponding rotation angle to...
static Vector add(const Vector &vector1, const Vector &vector2)
Adds the two given vectors.
boost::numeric::ublas::vector< double > Vector
Definition for an Isis::LinearAlgebra::Vector of doubles.
static QList< EulerAngle > toEulerAngles(const Matrix &rotationMatrix, const QList< int > axes)
Converts a rotation's representation from a matrix to a set of Euler angles with corresponding axes...
static bool isEmpty(const Vector &vector)
Determines whether the given vector is empty (i.e.
static Vector row(const Matrix &matrix, int rowIndex)
Returns a vector whose components match those of the given matrix row.
static Matrix inverse(const Matrix &matrix)
Returns the inverse of a 2x2 or 3x3 matrix.
static AxisAngle toAxisAngle(const Matrix &rotationMatrix)
Converts a rotation's representation from a matrix to a axis of rotation and its corresponding rotati...
static bool isZero(const Matrix &matrix)
Determines whether the given matrix is filled with zereos.
static double determinant(const Matrix &matrix)
Returns the determinant of the given 3x3 matrix.
#define _FILEINFO_
Macro for the filename and line number.
static Matrix identity(int size)
Returns the identity matrix of size NxN.
static Vector zeroVector(int size)
Returns a vector of given length that is filled with zeroes.
A type of error that cannot be classified as any of the other error types.
static void setRow(Matrix &matrix, const Vector &vector, int rowIndex)
Sets the row of the given matrix to the values of the given vector.
static void setVec4(Vector *v, double v0, double v1, double v2, double v3)
Fills the first four elements of the given vector with the given values.
static Vector normalizedCrossProduct(const Vector &vector1, const Vector &vector2)
Divides each vector by its corresponding absolute maximum, computes the cross product of the new vect...
static bool isUnit(const Vector &vector)
Determines whether the given vector is a unit vector.
Defines an angle and provides unit conversions.
static bool isOrthogonal(const Matrix &matrix)
Determines whether the given matrix is orthogonal by verifying that the matrix and its tranpose are i...
static Vector crossProduct(const Vector &vector1, const Vector &vector2)
Returns the cross product of two vectors.
QDebug operator<<(QDebug debug, const Hillshade &hillshade)
Print this class out to a QDebug object.
~LinearAlgebra()
Destructor for a LinearAlgebra object.
static void setVec3(Vector *v, double v0, double v1, double v2)
Fills the first three elements of the given vector with the given values.
static double absoluteMaximum(const Vector &vector)
Returns the maximum norm (L-infinity norm) for the given vector.
static Vector perpendicular(const Vector &vector1, const Vector &vector2)
Finds the unique vector P such that A = V + P, V is parallel to B and P is perpendicular to B...