7#include "LinearAlgebra.h"
18#include <boost/assign/std/vector.hpp>
19#include <boost/numeric/ublas/io.hpp>
20#include <boost/numeric/ublas/matrix.hpp>
21#include <boost/numeric/ublas/matrix_proxy.hpp>
29#include "IException.h"
62 if (matrix.size1() != matrix.size2())
return false;
64 for (
unsigned int row = 0;
row < matrix.size1();
row++) {
65 for (
unsigned int col = 0; col < matrix.size2(); col++) {
67 if (!qFuzzyCompare(matrix(
row, col), 1.0))
return false;
72 if (!qFuzzyCompare(1.0 + matrix(
row, col), 1.0))
return false;
90 if (matrix.size1() != matrix.size2())
return false;
115 if (matrix.size1() != matrix.size2())
return false;
171 Matrix unitMatrix(matrix.size1(), matrix.size2());
173 for (
unsigned int i = 0; i < unitMatrix.size1(); i++) {
177 double columnMagnitude =
magnitude(columnVector);
178 if ( columnMagnitude < 0.9 || columnMagnitude > 1.1 ) {
183 setColumn(unitMatrix, columnVector / columnMagnitude, i);
191 if ( det >= 0.9 && det <= 1.1 ) {
197 QString msg =
"Unable to determine whether the given matrix is a rotation matrix.";
212 for (
unsigned int i = 0; i < matrix.size1(); i++) {
213 for (
unsigned int j = 0; j < matrix.size2(); j++) {
216 if (!qFuzzyCompare(matrix(i, j) + 1.0, 1.0))
return false;
232 for (
unsigned int i = 0; i <
vector.size(); i++) {
235 if (!qFuzzyCompare(
vector(i) + 1.0, 1.0))
return false;
291 if (qFuzzyCompare(det + 1.0, 1.0)) {
293 QString msg =
"The given matrix is not invertible. The determinant is 0.0.";
298 double scale = 1 / det;
302 if (matrix.size1() == 2) {
303 inverse(0, 0) = scale * matrix(1, 1);
304 inverse(0, 1) = -scale * matrix(0, 1);
306 inverse(1, 0) = -scale * matrix(1, 0);
307 inverse(1, 1) = scale * matrix(0, 0);
311 inverse(0, 0) = scale * ( matrix(1, 1) * matrix(2, 2) - matrix(2, 1) * matrix(1, 2) );
312 inverse(0, 1) = scale * ( matrix(0, 2) * matrix(2, 1) - matrix(2, 2) * matrix(0, 1) );
313 inverse(0, 2) = scale * ( matrix(0, 1) * matrix(1, 2) - matrix(1, 1) * matrix(0, 2) );
315 inverse(1, 0) = scale * ( matrix(1, 2) * matrix(2, 0) - matrix(2, 2) * matrix(1, 0) );
316 inverse(1, 1) = scale * ( matrix(0, 0) * matrix(2, 2) - matrix(2, 0) * matrix(0, 2) );
317 inverse(1, 2) = scale * ( matrix(0, 2) * matrix(0, 1) - matrix(1, 2) * matrix(0, 0) );
319 inverse(2, 0) = scale * ( matrix(1, 0) * matrix(2, 1) - matrix(2, 0) * matrix(1, 1) );
320 inverse(2, 1) = scale * ( matrix(0, 1) * matrix(2, 0) - matrix(2, 1) * matrix(0, 0) );
321 inverse(2, 2) = scale * ( matrix(0, 0) * matrix(1, 1) - matrix(1, 0) * matrix(0, 1) );
326 QString msg =
"Unable to invert the given matrix.";
341 arma::mat arMat(matrix.size1(), matrix.size2());
342 for (
size_t i = 0; i < matrix.size1(); i++) {
343 for (
size_t j = 0; j < matrix.size2(); j++) {
344 arMat(i, j) = matrix(i, j);
348 arma::mat invArMat = arma::pinv(arMat);
352 for (
size_t i = 0; i < invArMat.n_rows; i++) {
353 for (
size_t j = 0; j < invArMat.n_cols; j++) {
354 inverse(i, j) = invArMat(i, j);
371 return boost::numeric::ublas::trans(matrix);
386 QString msg =
"Can not create identity matrix of negative size ["
392 for (
int i = 0; i < size; i++) {
409 boost::numeric::ublas::zero_matrix<double> m(rows, columns);
422 boost::numeric::ublas::zero_vector<double> v(size);
441 if ( (matrix.size1() != matrix.size2()) || (matrix.size1() != 2 && matrix.size1() != 3) ) {
442 QString msg =
"Unable to calculate the determinant for the given matrix. "
443 "This method only calculates the determinant for 2x2 or 3x3 matrices."
444 "The given matrix is [" +
toString((
int) matrix.size1()) +
"x"
445 +
toString((
int) matrix.size2()) +
"].";
449 if (matrix.size1() == 2) {
450 return matrix(0, 0) * matrix(1, 1) - matrix(1, 0) * matrix(0, 1);
454 return matrix(0, 0) * (matrix(1, 1)*matrix(2, 2) - matrix(1, 2)*matrix(2, 1))
455 - matrix(0, 1) * (matrix(1, 0)*matrix(2, 2) - matrix(1, 2)*matrix(2, 0))
456 + matrix(0, 2) * (matrix(1, 0)*matrix(2, 1) - matrix(1, 1)*matrix(2, 0));
480 QString msg =
"Unable to normalize the zero vector.";
532 return boost::numeric::ublas::norm_inf(
vector);
551 if (matrix1.size2() != matrix2.size1()) {
552 QString msg =
"Unable to multiply matrices with mismatched dimensions. "
553 "The left matrix has [" +
toString((
int) matrix1.size2())
554 +
"] columns and the right matrix has [" +
toString((
int) matrix2.size1())
561 return boost::numeric::ublas::prod(matrix1, matrix2);
582 if (matrix.size2() !=
vector.size()) {
583 QString msg =
"Unable to multiply matrix and vector with mismatched dimensions."
585 +
"] components and the given matrix has ["
586 +
toString((
int) matrix.size2()) +
"] columns.";
591 return boost::numeric::ublas::prod(matrix,
vector);
619 return scalar * matrix;
636 if (vector1.size() != vector2.size()) {
637 QString msg =
"Unable to add vectors with mismatched sizes."
638 "Vector1 has [" +
toString((
int) vector1.size())
639 +
"] components and vector2 has ["
640 +
toString((
int) vector2.size()) +
"] components.";
644 return vector1 + vector2;
662 if (vector1.size() != vector2.size()) {
663 QString msg =
"Unable to subtract vectors with mismatched sizes."
664 "Vector1 has [" +
toString((
int) vector1.size())
665 +
"] components and vector2 has ["
666 +
toString((
int) vector2.size()) +
"] components.";
670 return vector1 - vector2;
689 if ((vector1.size() != 3) || (vector2.size() != 3)) {
690 QString msg =
"Unable to calculate the cross product on vectors that are not size 3. "
691 "Vector1 has [" +
toString((
int) vector1.size())
692 +
"] components and vector2 has ["
693 +
toString((
int) vector2.size()) +
"] components.";
698 crossProd(0) = vector1(1)*vector2(2) - vector1(2)*vector2(1);
699 crossProd(1) = vector1(2)*vector2(0) - vector1(0)*vector2(2);
700 crossProd(2) = vector1(0)*vector2(1) - vector1(1)*vector2(0);
730 if (maxVector1 != 0.0) {
731 normalizedVector1 /= maxVector1;
734 if (maxVector2 != 0.0) {
735 normalizedVector2 /= maxVector2;
758 if (vector1.size() != vector2.size()) {
759 QString msg =
"Unable to compute the outer product for vectors with mismatched sizes."
760 "Vector1 has [" +
toString((
int) vector1.size())
761 +
"] components and vector2 has ["
762 +
toString((
int) vector2.size()) +
"] components.";
765 return boost::numeric::ublas::outer_prod(vector1, vector2);
802 if (vector1.size() != vector2.size()) {
803 QString msg =
"Unable to compute the dot product for vectors with mismatched sizes."
804 "Vector1 has [" +
toString((
int) vector1.size())
805 +
"] components and vector2 has ["
806 +
toString((
int) vector2.size()) +
"] components.";
809 return boost::numeric::ublas::inner_prod(vector1, vector2);
834 if (vector1.size() != vector2.size()) {
835 QString msg =
"Unable to project vector1 onto vector2 with mismatched sizes."
836 "Vector1 has [" +
toString((
int) vector1.size())
837 +
"] components and vector2 has ["
838 +
toString((
int) vector2.size()) +
"] components.";
852 return v1Dotv2/v2Dotv2 * vector2;
874 if ((
vector.size() != 3) || (axis.size() != 3)) {
875 QString msg =
"Unable to rotate vector about the given axis and angle. "
876 "Vectors must be of size 3 to perform rotation. "
878 +
"] components and the given axis has ["
879 +
toString((
int) axis.size()) +
"] components.";
901 double c = cos(angle.
radians());
902 double s = sin(angle.
radians());
903 Vector rplane = (c*v1) + (s*v2);
907 return (rplane + projVectorOnAxis);
947 Vector parallelVector =
project(vector1 / max1, vector2 / max2);
949 Vector perpendicularVector = vector1 - parallelVector * max1;
951 return perpendicularVector;
975 if ((rotationMatrix.size1() != 3) || (rotationMatrix.size2() != 3)) {
976 QString msg =
"Unable to convert the given matrix to an axis of rotation "
977 "and a rotation angle. A 3x3 matrix is required. The given matrix is ["
978 +
toString((
int) rotationMatrix.size1()) +
"x"
979 +
toString((
int) rotationMatrix.size2()) +
"].";
985 QString msg =
"Unable to convert the given matrix to an axis of rotation "
986 "and a rotation angle. The given matrix is not a rotation matrix.";
997 subQuaternion[0] = quaternion[1];
998 subQuaternion[1] = quaternion[2];
999 subQuaternion[2] = quaternion[3];
1008 if (
isZero(subQuaternion)) {
1015 else if (qFuzzyCompare(quaternion(0) + 1.0, 1.0)) {
1016 axis = subQuaternion;
1024 return qMakePair(axis, angle);
1055 if (axes.size() != 3) {
1056 QString msg =
"Unable to convert the given matrix to Euler angles. "
1057 "Exactly 3 axis codes are required. The given list has ["
1058 +
toString((
int) axes.size()) +
"] axes.";
1063 validAxes << 1 << 2 << 3;
1064 if (!validAxes.contains(axes[0])
1065 || !validAxes.contains(axes[1])
1066 || !validAxes.contains(axes[2])) {
1067 QString msg =
"Unable to convert the given matrix to Euler angles using the given axis codes "
1069 +
"]. Axis codes must be 1, 2, or 3.";
1073 if (axes[0] == axes[1] || axes[1] == axes[2]) {
1074 QString msg =
"Unable to convert the given matrix to Euler angles using the given axis codes "
1076 +
"]. The middle axis code must differ from its neighbors.";
1081 if ((rotationMatrix.size1() != 3) || (rotationMatrix.size2() != 3)) {
1082 QString msg =
"Unable to convert the given matrix to Euler angles. A 3x3 matrix is required. "
1083 "The given matrix is ["
1084 +
toString((
int) rotationMatrix.size1()) +
"x"
1085 +
toString((
int) rotationMatrix.size2()) +
"].";
1090 QString msg =
"Unable to convert the given matrix to Euler angles. "
1091 "The given matrix is not a rotation matrix.";
1103 for (
int i = 0; i < 3; i++) {
1110 nextAxis << 2 << 3 << 1;
1115 Angle angle1, angle2, angle3;
1116 if ( axes[0] == axes[2] ) {
1117 if ( axes[1] == nextAxis[axes[0] - 1] ) {
1123 int c = 6 - axes[0] - axes[1];
1125 change( axes[0] - 1, 2 ) = 1.0;
1126 change( axes[1] - 1, 0 ) = 1.0;
1127 change( c - 1, 1 ) = sign * 1.0;
1131 bool degen = ( qFuzzyCompare(tempRotation(0, 2) + 1.0, 1.0)
1132 && qFuzzyCompare(tempRotation(1, 2) + 1.0, 1.0) )
1133 || ( qFuzzyCompare(tempRotation(2, 0) + 1.0, 1.0)
1134 && qFuzzyCompare(tempRotation(2, 1) + 1.0, 1.0) )
1135 || qFuzzyCompare( qAbs(tempRotation(2, 2)), 1.0 );
1140 angle2.
setRadians( acos( tempRotation(2, 2) ) );
1141 angle1.
setRadians( atan2( tempRotation(0, 1), tempRotation(0, 0) ) );
1146 angle3.
setRadians( atan2( tempRotation(0, 2), tempRotation(1, 2) ) );
1147 angle2.
setRadians( acos( tempRotation(2, 2) ) );
1148 angle1.
setRadians( atan2( tempRotation(2, 0), -tempRotation(2, 1) ) );
1157 if ( axes[1] == nextAxis[axes[0] - 1] ) {
1163 change( axes[0] - 1, 0 ) = 1.0;
1164 change( axes[1] - 1, 1 ) = 1.0;
1165 change( axes[2] - 1, 2 ) = sign * 1.0;
1168 bool degen = ( qFuzzyCompare(tempRotation(0, 0) + 1.0, 1.0)
1169 && qFuzzyCompare(tempRotation(0, 1) + 1.0, 1.0) )
1170 || ( qFuzzyCompare(tempRotation(1, 2) + 1.0, 1.0)
1171 && qFuzzyCompare(tempRotation(2, 2) + 1.0, 1.0) )
1172 || qFuzzyCompare( qAbs(tempRotation(0, 2)), 1.0 );
1177 angle2.
setRadians( asin( -tempRotation(0, 2) ) );
1178 angle1.
setRadians( sign * atan2( -tempRotation(1, 0), tempRotation(1, 1) ) );
1183 angle3.
setRadians( atan2( tempRotation(1, 2), tempRotation(2, 2) ) );
1184 angle2.
setRadians( asin( -tempRotation(0, 2) ) );
1185 angle1.
setRadians( sign * atan2( tempRotation(0, 1), tempRotation(0, 0) ) );
1190 eulerAngles << qMakePair(angle3, axes[0])
1191 << qMakePair(angle2, axes[1])
1192 << qMakePair(angle1, axes[2]);
1213 if ((rotationMatrix.size1() != 3) || (rotationMatrix.size2() != 3)) {
1214 QString msg =
"Unable to convert the given matrix to a quaternion. "
1215 "A 3x3 matrix is required. The given matrix is ["
1216 +
toString((
int) rotationMatrix.size1()) +
"x"
1217 +
toString((
int) rotationMatrix.size2()) +
"].";
1223 QString msg =
"Unable to convert the given matrix to an axis of rotation "
1224 "and a rotation angle. The given matrix is not a rotation matrix.";
1319 double trace = rotationMatrix(0, 0) + rotationMatrix(1, 1) + rotationMatrix(2, 2);
1320 double mtrace = 1.0 - trace;
1323 double cc4 = 1.0 + trace;
1325 double s114 = mtrace + 2.0*rotationMatrix(0, 0);
1326 double s224 = mtrace + 2.0*rotationMatrix(1, 1);
1327 double s334 = mtrace + 2.0*rotationMatrix(2, 2);
1331 double normalizingFactor;
1335 quaternion[0] = sqrt( cc4 * 0.25 );
1336 normalizingFactor = 1.0 / ( quaternion[0] * 4.0 );
1341 quaternion[1] = ( rotationMatrix(2, 1) - rotationMatrix(1, 2) ) * normalizingFactor;
1342 quaternion[2] = ( rotationMatrix(0, 2) - rotationMatrix(2, 0) ) * normalizingFactor;
1343 quaternion[3] = ( rotationMatrix(1, 0) - rotationMatrix(0, 1) ) * normalizingFactor;
1346 else if ( s114 >= 1.0 ) {
1348 quaternion[1] = sqrt( s114 * 0.25 );
1349 normalizingFactor = 1.0 /( quaternion[1] * 4.0 );
1351 quaternion[0] = ( rotationMatrix(2, 1) - rotationMatrix(1, 2) ) * normalizingFactor;
1352 quaternion[2] = ( rotationMatrix(0, 1) + rotationMatrix(1, 0) ) * normalizingFactor;
1353 quaternion[3] = ( rotationMatrix(0, 2) + rotationMatrix(2, 0) ) * normalizingFactor;
1356 else if ( s224 >= 1.0 ) {
1358 quaternion[2] = sqrt( s224 * 0.25 );
1359 normalizingFactor = 1.0 /( quaternion[2] * 4.0 );
1361 quaternion[0] = ( rotationMatrix(0, 2) - rotationMatrix(2, 0) ) * normalizingFactor;
1362 quaternion[1] = ( rotationMatrix(0, 1) + rotationMatrix(1, 0) ) * normalizingFactor;
1363 quaternion[3] = ( rotationMatrix(1, 2) + rotationMatrix(2, 1) ) * normalizingFactor;
1367 quaternion[3] = sqrt( s334 * 0.25 );
1368 normalizingFactor = 1.0 /( quaternion[3] * 4.0 );
1370 quaternion[0] = ( rotationMatrix(1, 0) - rotationMatrix(0, 1) ) * normalizingFactor;
1371 quaternion[1] = ( rotationMatrix(0, 2) + rotationMatrix(2, 0) ) * normalizingFactor;
1372 quaternion[2] = ( rotationMatrix(1, 2) + rotationMatrix(2, 1) ) * normalizingFactor;
1377 if ( !
isUnit(quaternion) ) {
1382 if ( quaternion[0] < 0.0 ) {
1406 if (axis.size() != 3) {
1407 QString msg =
"Unable to convert the given vector and angle to a rotation matrix. "
1408 "The given vector with size [" +
toString((
int) axis.size())
1409 +
"] is not a 3D axis vector.";
1420 for (
unsigned int i = 0; i < axis.size(); i++) {
1429 return rotationMatrix;
1466 validAxes << 1 << 2 << 3;
1467 if (!validAxes.contains(angle3.second)
1468 || !validAxes.contains(angle2.second)
1469 || !validAxes.contains(angle1.second)) {
1470 QString msg =
"Unable to convert the given Euler angles to a matrix using the given axis "
1471 "codes [" +
toString(angle3.second) +
", " +
toString(angle2.second) +
", "
1472 +
toString(angle1.second) +
"]. Axis codes must be 1, 2, or 3.";
1483 double sinAngle = sin(angle1.first.radians());
1484 double cosAngle = cos(angle1.first.radians());
1486 double indexOffset = ((angle1.second % 3) + 3) % 3;
1488 indexSet << 2 << 0 << 1 << 2 << 0;
1489 int index1 = indexSet[indexOffset + 0];
1490 int index2 = indexSet[indexOffset + 1];
1491 int index3 = indexSet[indexOffset + 2];
1492 m(index1, index1) = 1.0; m(index1, index2) = 0.0; m(index1, index3) = 0.0;
1493 m(index2, index1) = 0.0; m(index2, index2) = cosAngle; m(index2, index3) = sinAngle;
1494 m(index3, index1) = 0.0; m(index3, index2) = -sinAngle; m(index3, index3) = cosAngle;
1502 sinAngle = sin(angle2.first.radians());
1503 cosAngle = cos(angle2.first.radians());
1504 indexOffset = ((angle2.second % 3) + 3) % 3;
1505 index1 = indexSet[indexOffset + 0];
1506 index2 = indexSet[indexOffset + 1];
1507 index3 = indexSet[indexOffset + 2];
1510 for (
int i = 0; i < 3; i++) {
1511 tempMatrix(index1, i) = m(index1, i);
1512 tempMatrix(index2, i) = cosAngle*m(index2, i) + sinAngle*m(index3, i);
1513 tempMatrix(index3, i) = -sinAngle*m(index2, i) + cosAngle*m(index3, i);
1517 sinAngle = sin(angle3.first.radians());
1518 cosAngle = cos(angle3.first.radians());
1519 indexOffset = ((angle3.second % 3) + 3) % 3;
1520 index1 = indexSet[indexOffset + 0];
1521 index2 = indexSet[indexOffset + 1];
1522 index3 = indexSet[indexOffset + 2];
1523 for (
int i = 0; i < 3; i++) {
1524 m(index1, i) = tempMatrix(index1, i);
1525 m(index2, i) = cosAngle*tempMatrix(index2, i) + sinAngle*tempMatrix(index3, i);
1526 m(index3, i) = -sinAngle*tempMatrix(index2, i) + cosAngle*tempMatrix(index3, i);
1547 if (eulerAngles.size() != 3) {
1548 QString msg =
"Unable to convert the given Euler angles to a matrix. "
1549 "Exactly 3 Euler angles are required. The given list has ["
1550 +
toString((
int) eulerAngles.size()) +
"] angles.";
1573 if (quaternion.size() != 4) {
1574 QString msg =
"Unable to convert the given vector to a rotation matrix. "
1575 "The given vector with [" +
toString((
int) quaternion.size())
1576 +
"] components is not a quaternion.";
1587 double q0q1 = q[0] * q[1];
1588 double q0q2 = q[0] * q[2];
1589 double q0q3 = q[0] * q[3];
1591 double q1q1 = q[1] * q[1];
1592 double q1q2 = q[1] * q[2];
1593 double q1q3 = q[1] * q[3];
1595 double q2q2 = q[2] * q[2];
1596 double q2q3 = q[2] * q[3];
1598 double q3q3 = q[3] * q[3];
1600 matrix(0, 0) = 1.0 - 2.0 * ( q2q2 + q3q3 );
1601 matrix(0, 1) = 2.0 * ( q1q2 - q0q3 );
1602 matrix(0, 2) = 2.0 * ( q1q3 + q0q2 );
1604 matrix(1, 0) = 2.0 * ( q1q2 + q0q3 );
1605 matrix(1, 1) = 1.0 - 2.0 * ( q1q1 + q3q3 );
1606 matrix(1, 2) = 2.0 * ( q2q3 - q0q1 );
1608 matrix(2, 0) = 2.0 * ( q1q3 - q0q2 );
1609 matrix(2, 1) = 2.0 * ( q2q3 + q0q1 );
1610 matrix(2, 2) = 1.0 - 2.0 * ( q1q1 + q2q2 );
1627 if ( (rowIndex+1 > (
int) matrix.size1()) ) {
1628 QString msg =
"Unable to set the matrix row to the given vector. Row index "
1629 +
toString(rowIndex) +
" is out of bounds. The given matrix only has "
1630 +
toString((
int) matrix.size1()) +
" rows." ;
1633 boost::numeric::ublas::row(matrix, rowIndex) =
vector;
1648 if ( (columnIndex+1 > (
int) matrix.size2()) ) {
1649 QString msg =
"Unable to set the matrix column to the given vector. Column index "
1650 +
toString(columnIndex) +
" is out of bounds. The given matrix only has "
1651 +
toString((
int) matrix.size1()) +
" columns." ;
1654 boost::numeric::ublas::column(matrix, columnIndex) =
vector;
1671 if ( (rowIndex+1 > (
int) matrix.size1()) ) {
1672 QString msg =
"Unable to get the matrix row to the given vector. Row index "
1673 +
toString(rowIndex) +
" is out of bounds. The given matrix only has "
1674 +
toString((
int) matrix.size1()) +
" rows." ;
1677 return boost::numeric::ublas::row(matrix, rowIndex);
1694 if ( (columnIndex+1 > (
int) matrix.size2()) ) {
1695 QString msg =
"Unable to get the matrix column to the given vector. Column index "
1696 +
toString(columnIndex) +
" is out of bounds. The given matrix only has "
1697 +
toString((
int) matrix.size1()) +
" columns." ;
1700 return boost::numeric::ublas::column(matrix, columnIndex);
1791 for (
int i = 0; i < size; i++) {
1792 sub[i] = v[i + startIndex];
1811 QDebugStateSaver saver(dbg);
1829 QDebugStateSaver saver(dbg);
1830 for (
unsigned int i = 0; i < matrix.size1(); i++) {
1831 dbg.noquote() <<
" ";
1832 for (
unsigned int j = 0; j < matrix.size2(); j++) {
1833 dbg.noquote() <<
toString(matrix(i, j), 15) <<
" ";
1835 dbg.noquote() << Qt::endl;
1853 QString result =
"( ";
1854 for (
unsigned int i = 0; i < vector.size(); i++) {
1855 result +=
toString(vector(i), precision);
1856 if (i != vector.size() - 1) result +=
", ";
Defines an angle and provides unit conversions.
double radians() const
Convert an angle to a double.
void setRadians(double radians)
Set the angle in units of Radians.
@ Unknown
A type of error that cannot be classified as any of the other error types.
@ Programmer
This error is for when a programmer made an API call that was illegal.
static AxisAngle toAxisAngle(const Matrix &rotationMatrix)
Converts a rotation's representation from a matrix to a axis of rotation and its corresponding rotati...
static Matrix identity(int size)
Returns the identity matrix of size NxN.
~LinearAlgebra()
Destructor for a LinearAlgebra object.
static Vector subtract(const Vector &vector1, const Vector &vector2)
Subtracts the right vector from the left vector.
static Vector crossProduct(const Vector &vector1, const Vector &vector2)
Returns the cross product of two vectors.
static Vector zeroVector(int size)
Returns a vector of given length that is filled with zeroes.
static Vector vector(double v0, double v1, double v2)
Constructs a 3 dimensional vector with the given component values.
LinearAlgebra()
Default constructor for a LinearAlgebra object.
static Vector add(const Vector &vector1, const Vector &vector2)
Adds the two given vectors.
static Vector row(const Matrix &matrix, int rowIndex)
Returns a vector whose components match those of the given matrix row.
static Vector subVector(const Vector &v, int start, int size)
Constructs a vector of given size using the given vector and starting index.
static Matrix toMatrix(const AxisAngle &axisAngle)
Converts a rotation's representation from an axis of rotation and its corresponding rotation angle to...
static bool isUnit(const Vector &vector)
Determines whether the given vector is a unit vector.
static Vector column(const Matrix &matrix, int columnIndex)
Returns a vector whose components match those of the given matrix column.
static bool isEmpty(const Vector &vector)
Determines whether the given vector is empty (i.e.
static Vector project(const Vector &vector1, const Vector &vector2)
Compute the vector projection of vector1 onto vector2.
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 setColumn(Matrix &matrix, const Vector &vector, int columnIndex)
Sets the column 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 Matrix multiply(const Matrix &matrix1, const Matrix &matrix2)
Returns the product of two matrices.
static double innerProduct(const Vector &vector1, const Vector &vector2)
Computes the inner product of the given vectors.
static Matrix zeroMatrix(int rows, int columns)
Returns a matrix with given dimensions that is filled with zeroes.
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 Vector rotate(const Vector &vector, const Vector &axis, Angle angle)
Rotates a vector about an axis vector given a specified angle.
QPair< Vector, Angle > AxisAngle
Definition for an Axis-Angle pair.
static double magnitude(const Vector &vector)
Computes the magnitude (i.e., the length) of the given vector using the Euclidean norm (L2 norm).
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,...
static Matrix pseudoinverse(const Matrix &matrix)
Returns the pseudoinverse of a matrix.
static double determinant(const Matrix &matrix)
Returns the determinant of the given 3x3 matrix.
static Matrix inverse(const Matrix &matrix)
Returns the inverse of a 2x2 or 3x3 matrix.
boost::numeric::ublas::vector< double > Vector
Definition for an Isis::LinearAlgebra::Vector of doubles.
static bool isIdentity(const Matrix &matrix)
Determines whether the given matrix is the identity.
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 bool isZero(const Matrix &matrix)
Determines whether the given matrix is filled with zereos.
QPair< Angle, int > EulerAngle
Definition for an EulerAngle pair.
static Matrix outerProduct(const Vector &vector1, const Vector &vector2)
Computes the outer product of the given vectors.
boost::numeric::ublas::matrix< double > Matrix
Definition for an Isis::LinearAlgebra::Matrix of doubles.
static bool isOrthogonal(const Matrix &matrix)
Determines whether the given matrix is orthogonal by verifying that the matrix and its tranpose are i...
static double absoluteMaximum(const Vector &vector)
Returns the maximum norm (L-infinity norm) for the given vector.
static Vector toQuaternion(const Matrix &rotationMatrix)
Converts a rotation's representation from a matrix to a unit quaternion.
static Matrix transpose(const Matrix &matrix)
Returns the transpose of the given matrix.
static double dotProduct(const Vector &vector1, const Vector &vector2)
Computes the dot product of the given vectors.
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 Vector normalize(const Vector &vector)
Returns a unit vector that is codirectional with the given vector by dividing each component of the v...
static bool isRotationMatrix(const Matrix &matrix)
Determines whether the given matrix is a rotation matrix.
This is free and unencumbered software released into the public domain.
This is free and unencumbered software released into the public domain.
This is free and unencumbered software released into the public domain.
QString toString(bool boolToConvert)
Global function to convert a boolean to a string.
QDebug operator<<(QDebug debug, const Hillshade &hillshade)
Print this class out to a QDebug object.
const double PI
The mathematical constant PI.