Isis 3 Programmer Reference
LinearAlgebra.cpp
Go to the documentation of this file.
1 
23 #include "LinearAlgebra.h"
24 
25 // std library
26 #include <cmath>
27 
28 // Qt Library
29 #include <QDebug>
30 #include <QtDebug>
31 #include <QtGlobal>
32 
33 // boost library
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>
38 
39 // other Isis library
40 #include "Angle.h"
41 #include "Constants.h"
42 #include "IException.h"
43 #include "IString.h"
44 
45 
46 namespace Isis {
53  }
54 
55 
60  }
61 
62 
74  bool LinearAlgebra::isIdentity(const Matrix &matrix) {
75  if (matrix.size1() != matrix.size2()) return false;
76 
77  for (unsigned int row = 0; row < matrix.size1(); row++) {
78  for (unsigned int col = 0; col < matrix.size2(); col++) {
79  if (row == col) {
80  if (!qFuzzyCompare(matrix(row, col), 1.0)) return false;
81  }
82  else {
83  // When using qFuzzy compare against a number that is zero
84  // you must offset the number per Qt documentation
85  if (!qFuzzyCompare(1.0 + matrix(row, col), 1.0)) return false;
86  }
87  }
88  }
89  return true;
90  }
91 
92 
102  bool LinearAlgebra::isOrthogonal(const Matrix &matrix) {
103  if (matrix.size1() != matrix.size2()) return false;
105  if (isIdentity(test)) {
106  // no need to test that transpose is both left and right inverse since the matrix is square
107  return true;
108  }
109  return false;
110  }
111 
112 
125  // derived from naif's isrot routine
127  // rotation matrices must be square
128  if (matrix.size1() != matrix.size2()) return false;
129 
130  /* An exerpt from NAIF's isrot() routine...
131  *
132  * A property of rotation matrices is that their columns form a
133  * right-handed, orthonormal basis in 3-dimensional space. The
134  * converse is true: all 3x3 matrices with this property are
135  * rotation matrices.
136  *
137  * An ordered set of three vectors V1, V2, V3 forms a right-handed,
138  * orthonormal basis if and only if
139  *
140  * 1) || V1 || = || V2 || = || V3 || = 1
141  *
142  * 2) V3 = V1 x V2. Since V1, V2, and V3 are unit vectors,
143  * we also have
144  *
145  * < V3, V1 x V2 > = 1.
146  *
147  * This quantity is the determinant of the matrix whose
148  * colums are V1, V2 and V3.
149  *
150  * When finite precision numbers are used, rotation matrices will
151  * usually fail to satisfy these criteria exactly. We must use
152  * criteria that indicate approximate conformance to the criteria
153  * listed above. We choose
154  *
155  * 1) | || Vi || - 1 | < NTOL, i = 1, 2, 3.
156  * -
157  *
158  * 2) Let
159  *
160  * Vi
161  * Ui = ------ , i = 1, 2, 3.
162  * ||Vi||
163  *
164  * Then we require
165  *
166  * | < U3, U1 x U2 > - 1 | < DTOL;
167  * -
168  *
169  * equivalently, letting U be the matrix whose columns
170  * are U1, U2, and U3, we insist on
171  *
172  * | det(U) - 1 | < DTOL.
173  * _
174  * The columns of M must resemble unit vectors. If the norms are
175  * outside of the allowed range, M is not a rotation matrix.
176  *
177  * Also, the columns of M are required to be pretty nearly orthogonal. The
178  * discrepancy is gauged by taking the determinant of the matrix UNIT,
179  * computed below, whose columns are the unitized columns of M.
180  */
181 
182  // create a matrix whose columns are unit vectors of the columns of the
183  // given matrix
184  Matrix unitMatrix(matrix.size1(), matrix.size2());
185 
186  for (unsigned int i = 0; i < unitMatrix.size1(); i++) {
187  // get the ith column vector from the given matrix
188  Vector columnVector = column(matrix, i);
189  // get the magnitude and if it is not near 1, this is not a rotation matrix
190  double columnMagnitude = magnitude(columnVector);
191  if ( columnMagnitude < 0.9 || columnMagnitude > 1.1 ) {
192  return false;
193  }
194 
195  // put the unitized row into the unitized matrix
196  setColumn(unitMatrix, columnVector / columnMagnitude, i);
197 
198  }
199 
200  try {
201  // get the determinant of the unitized matrix and if it is not near 1,
202  // this is not a rotation matrix
203  double det = determinant(unitMatrix);
204  if ( det >= 0.9 && det <= 1.1 ) {
205  return true;
206  }
207  }
208  catch (IException &e) {
209  // since we can only calculate the determinant for 2x2 or 3x3
210  QString msg = "Unable to determine whether the given matrix is a rotation matrix.";
212  }
213  return false;
214  }
215 
216 
224  bool LinearAlgebra::isZero(const Matrix &matrix) {
225  for (unsigned int i = 0; i < matrix.size1(); i++) {
226  for (unsigned int j = 0; j < matrix.size2(); j++) {
227  // When using qFuzzy compare against a number that is zero
228  // you must offset the number per Qt documentation
229  if (!qFuzzyCompare(matrix(i, j) + 1.0, 1.0)) return false;
230  }
231  }
232 
233  return true;
234  }
235 
236 
244  bool LinearAlgebra::isZero(const Vector &vector) {
245  for (unsigned int i = 0; i < vector.size(); i++) {
246  // When using qFuzzy compare against a number that is zero
247  // you must offset the number per Qt documentation
248  if (!qFuzzyCompare(vector(i) + 1.0, 1.0)) return false;
249  }
250 
251  return true;
252  }
253 
254 
264  return vector.empty();
265  }
266 
267 
275  bool LinearAlgebra::isUnit(const Vector &vector) {
276  if (qFuzzyCompare(dotProduct(vector, vector), 1.0)) {
277  return true;
278  }
279  return false;
280  }
281 
282 
296  try {
297  if (isOrthogonal(matrix)) {
298  return transpose(matrix);
299  }
300 
301  // the determinant method will verify the size of the matrix is 2x2 or 3x3
302  double det = determinant(matrix);
303 
304  if (qFuzzyCompare(det + 1.0, 1.0)) {
305  // the inverse exists <==> the determinant is not 0.0
306  QString msg = "The given matrix is not invertible. The determinant is 0.0.";
308  }
309 
310  // since the determinant is not zero, we can calculate the reciprocal
311  double scale = 1 / det;
312 
313  LinearAlgebra::Matrix inverse = identity(matrix.size1());
314  // find the inverse for 2x2
315  if (matrix.size1() == 2) {
316  inverse(0, 0) = scale * matrix(1, 1);
317  inverse(0, 1) = -scale * matrix(0, 1);
318 
319  inverse(1, 0) = -scale * matrix(1, 0);
320  inverse(1, 1) = scale * matrix(0, 0);
321  return inverse;
322  }
323  // else we have a 3x3
324  inverse(0, 0) = scale * ( matrix(1, 1) * matrix(2, 2) - matrix(2, 1) * matrix(1, 2) );
325  inverse(0, 1) = scale * ( matrix(0, 2) * matrix(2, 1) - matrix(2, 2) * matrix(0, 1) );
326  inverse(0, 2) = scale * ( matrix(0, 1) * matrix(1, 2) - matrix(1, 1) * matrix(0, 2) );
327 
328  inverse(1, 0) = scale * ( matrix(1, 2) * matrix(2, 0) - matrix(2, 2) * matrix(1, 0) );
329  inverse(1, 1) = scale * ( matrix(0, 0) * matrix(2, 2) - matrix(2, 0) * matrix(0, 2) );
330  inverse(1, 2) = scale * ( matrix(0, 2) * matrix(0, 1) - matrix(1, 2) * matrix(0, 0) );
331 
332  inverse(2, 0) = scale * ( matrix(1, 0) * matrix(2, 1) - matrix(2, 0) * matrix(1, 1) );
333  inverse(2, 1) = scale * ( matrix(0, 1) * matrix(2, 0) - matrix(2, 1) * matrix(0, 0) );
334  inverse(2, 2) = scale * ( matrix(0, 0) * matrix(1, 1) - matrix(1, 0) * matrix(0, 1) );
335 
336  return inverse;
337  }
338  catch (IException &e) {
339  QString msg = "Unable to invert the given matrix.";
341  }
342  }
343 
344 
353  // no error testing required so call the boost transpose function
354  return boost::numeric::ublas::trans(matrix);
355  }
356 
357 
368  if (size < 1) {
369  QString msg = "Can not create identity matrix of negative size ["
370  + toString((int) size) + "].";
372  }
373 
374  LinearAlgebra::Matrix m(size, size, 0.0);
375  for (int i = 0; i < size; i++) {
376  m(i, i) = 1.0;
377  }
378 
379  return m;
380  }
381 
382 
392  boost::numeric::ublas::zero_matrix<double> m(rows, columns);
393  return m;
394  }
395 
396 
405  boost::numeric::ublas::zero_vector<double> v(size);
406  return v;
407  }
408 
409 
423  double LinearAlgebra::determinant(const Matrix &matrix) {
424  if ( (matrix.size1() != matrix.size2()) || (matrix.size1() != 2 && matrix.size1() != 3) ) {
425  QString msg = "Unable to calculate the determinant for the given matrix. "
426  "This method only calculates the determinant for 2x2 or 3x3 matrices."
427  "The given matrix is [" + toString((int) matrix.size1()) + "x"
428  + toString((int) matrix.size2()) + "].";
430  }
431 
432  if (matrix.size1() == 2) {
433  return matrix(0, 0) * matrix(1, 1) - matrix(1, 0) * matrix(0, 1);
434  }
435  else {
436 
437  return matrix(0, 0) * (matrix(1, 1)*matrix(2, 2) - matrix(1, 2)*matrix(2, 1))
438  - matrix(0, 1) * (matrix(1, 0)*matrix(2, 2) - matrix(1, 2)*matrix(2, 0))
439  + matrix(0, 2) * (matrix(1, 0)*matrix(2, 1) - matrix(1, 1)*matrix(2, 0));
440  }
441 
442  }
443 
444 
461  // impossible to unitize the zero vector
463  QString msg = "Unable to normalize the zero vector.";
465  }
467  }
468 
469 
487  double LinearAlgebra::magnitude(const Vector &vector) {
488  double magnitude;
489  // avoid dividing by maxNorm = 0
490  // no stabilization is needed for the zero vector, magnitude is always 0.0
492  magnitude = 0.0;
493  }
494  else {
495  double maxNorm = LinearAlgebra::absoluteMaximum(vector);
496  magnitude = maxNorm * boost::numeric::ublas::norm_2(vector / maxNorm);
497  }
498  return magnitude;
499  }
500 
501 
514  double LinearAlgebra::absoluteMaximum(const Vector &vector) {
515  return boost::numeric::ublas::norm_inf(vector);
516  }
517 
518 
532  LinearAlgebra::Matrix LinearAlgebra::multiply(const Matrix &matrix1, const Matrix &matrix2) {
533  // Check to make sure we can multiply
534  if (matrix1.size2() != matrix2.size1()) {
535  QString msg = "Unable to multiply matrices with mismatched dimensions. "
536  "The left matrix has [" + toString((int) matrix1.size2())
537  + "] columns and the right matrix has [" + toString((int) matrix2.size1())
538  + "] rows.";
539 
541  }
542 
543  // Use boost product routine
544  return boost::numeric::ublas::prod(matrix1, matrix2);
545  }
546 
547 
564  // Check to make sure we can multiply
565  if (matrix.size2() != vector.size()) {
566  QString msg = "Unable to multiply matrix and vector with mismatched dimensions."
567  "The given vector has [" + toString((int) vector.size())
568  + "] components and the given matrix has ["
569  + toString((int) matrix.size2()) + "] columns.";
571  }
572 
573  // Use boost product routine
574  return boost::numeric::ublas::prod(matrix, vector);
575  }
576 
577 
586  LinearAlgebra::Vector LinearAlgebra::multiply(double scalar, const Vector &vector) {
587  // no error checks necessary (use the boost operator* method)
588  return scalar * vector;
589  }
590 
591 
600  LinearAlgebra::Matrix LinearAlgebra::multiply(double scalar, const Matrix &matrix) {
601  // no error checks necessary
602  return scalar * matrix;
603  }
604 
605 
617  LinearAlgebra::Vector LinearAlgebra::add(const Vector &vector1, const Vector &vector2) {
618  // Vectors best be the same size
619  if (vector1.size() != vector2.size()) {
620  QString msg = "Unable to add vectors with mismatched sizes."
621  "Vector1 has [" + toString((int) vector1.size())
622  + "] components and vector2 has ["
623  + toString((int) vector2.size()) + "] components.";
625  }
626 
627  return vector1 + vector2;
628  }
629 
630 
643  LinearAlgebra::Vector LinearAlgebra::subtract(const Vector &vector1, const Vector &vector2) {
644  // Vectors best be the same size
645  if (vector1.size() != vector2.size()) {
646  QString msg = "Unable to subtract vectors with mismatched sizes."
647  "Vector1 has [" + toString((int) vector1.size())
648  + "] components and vector2 has ["
649  + toString((int) vector2.size()) + "] components.";
651  }
652 
653  return vector1 - vector2;
654  }
655 
656 
672  if ((vector1.size() != 3) || (vector2.size() != 3)) {
673  QString msg = "Unable to calculate the cross product on vectors that are not size 3. "
674  "Vector1 has [" + toString((int) vector1.size())
675  + "] components and vector2 has ["
676  + toString((int) vector2.size()) + "] components.";
678  }
679 
680  Vector crossProd(3);
681  crossProd(0) = vector1(1)*vector2(2) - vector1(2)*vector2(1);
682  crossProd(1) = vector1(2)*vector2(0) - vector1(0)*vector2(2);
683  crossProd(2) = vector1(0)*vector2(1) - vector1(1)*vector2(0);
684 
685  return crossProd;
686  }
687 
688 
703  // is this derived from naif's ucrss routine???
705  const Vector &vector2) {
706 
707  double maxVector1 = LinearAlgebra::absoluteMaximum(vector1);
708  double maxVector2 = LinearAlgebra::absoluteMaximum(vector2);
709 
710  LinearAlgebra::Vector normalizedVector1 = vector1;
711  LinearAlgebra::Vector normalizedVector2 = vector2;
712 
713  if (maxVector1 != 0.0) {
714  normalizedVector1 /= maxVector1;
715  }
716 
717  if (maxVector2 != 0.0) {
718  normalizedVector2 /= maxVector2;
719  }
720 
721  LinearAlgebra::Vector vcross = crossProduct(normalizedVector1, normalizedVector2);
722 
723  return normalize(vcross);
724  }
725 
726 
741  if (vector1.size() != vector2.size()) {
742  QString msg = "Unable to compute the outer product for vectors with mismatched sizes."
743  "Vector1 has [" + toString((int) vector1.size())
744  + "] components and vector2 has ["
745  + toString((int) vector2.size()) + "] components.";
747  }
748  return boost::numeric::ublas::outer_prod(vector1, vector2);
749  }
750 
751 
764  double LinearAlgebra::dotProduct(const Vector &vector1, const Vector &vector2) {
765  // no error check needed - this is done by innerProduct
766  return LinearAlgebra::innerProduct(vector1, vector2);
767  }
768 
769 
784  double LinearAlgebra::innerProduct(const Vector &vector1, const Vector &vector2) {
785  if (vector1.size() != vector2.size()) {
786  QString msg = "Unable to compute the dot product for vectors with mismatched sizes."
787  "Vector1 has [" + toString((int) vector1.size())
788  + "] components and vector2 has ["
789  + toString((int) vector2.size()) + "] components.";
791  }
792  return boost::numeric::ublas::inner_prod(vector1, vector2);
793  }
794 
795 
815  // derived from naif's vproj routine
816  LinearAlgebra::Vector LinearAlgebra::project(const Vector &vector1, const Vector &vector2) {
817  if (vector1.size() != vector2.size()) {
818  QString msg = "Unable to project vector1 onto vector2 with mismatched sizes."
819  "Vector1 has [" + toString((int) vector1.size())
820  + "] components and vector2 has ["
821  + toString((int) vector2.size()) + "] components.";
823  }
824 
825  // if vector2 is the zero vector, then the projection of vector1 onto vector2
826  // is also the zero vector
827  if (LinearAlgebra::isZero(vector2)) return vector2;
828 
829  // find the dot product of the two vectors
830  double v1Dotv2 = LinearAlgebra::dotProduct(vector1, vector2);
831 
832  // find the square of the length of vector2 (this is the dot product of the vector2 with itself)
833  double v2Dotv2 = LinearAlgebra::dotProduct(vector2, vector2);
834 
835  return v1Dotv2/v2Dotv2 * vector2;
836  }
837 
838 
854  // derived from naif's vrotv routine
856  Angle angle) {
857  if ((vector.size() != 3) || (axis.size() != 3)) {
858  QString msg = "Unable to rotate vector about the given axis and angle. "
859  "Vectors must be of size 3 to perform rotation. "
860  "The given vector has [" + toString((int) vector.size())
861  + "] components and the given axis has ["
862  + toString((int) axis.size()) + "] components.";
864  }
865 
866  // if the given vector is the zero vector, then the rotation about the vector
867  // is also the zero vector
868  if (isZero(axis)) return vector;
869 
870  // Compute the unit vector that is codirectional with the given axis
871  Vector axisUnitVector = normalize(axis);
872 
873  // Compute the projection of the given vector onto the axis unit vector.
874  Vector projVectorOnAxis = project(vector, axisUnitVector);
875 
876  // Compute the component of the input orthogonal to the AXIS. Call it V1.
877  Vector v1 = vector - projVectorOnAxis;
878 
879  // Rotate V1 by 90 degrees about the AXIS and call the result V2.
880  Vector v2 = LinearAlgebra::crossProduct(axisUnitVector, v1);
881 
882  // Compute cos(angle)*v1 + sin(angle)*v2. This is v1 rotated about
883  // the axis in the plane normal to the axis, call the result rplane
884  double c = cos(angle.radians());
885  double s = sin(angle.radians());
886  Vector rplane = (c*v1) + (s*v2);
887 
888  // Add the rotated component in the normal plane to axis to the
889  // projection of v onto axis to obtain rotation.
890  return (rplane + projVectorOnAxis);
891  }
892 
893 
912  // Derived from NAIF's vperp routine
913  // possible that we may need to restrict size to dimesion 3, but doesn't seem
914  // necessary in the code.
916 
917  // if vector1 is the zero vector, just set p to zero and return.
918  if ( isZero(vector1) ) {
919  return vector1;
920  }
921 
922  // if vector2 is the zero vector, just set p to vector1 and return.
923  if ( isZero(vector2) ) {
924  return vector1;
925  }
926 
927  // normalize (using the max norm) the given vectors and project the first onto the second
928  double max1 = absoluteMaximum(vector1);
929  double max2 = absoluteMaximum(vector2);
930  Vector parallelVector = project(vector1 / max1, vector2 / max2);
931 
932  Vector perpendicularVector = vector1 - parallelVector * max1;
933 
934  return perpendicularVector;
935 
936  }
937 
938 
955  // Derived from NAIF's raxisa routine
957 
958  if ((rotationMatrix.size1() != 3) || (rotationMatrix.size2() != 3)) {
959  QString msg = "Unable to convert the given matrix to an axis of rotation "
960  "and a rotation angle. A 3x3 matrix is required. The given matrix is ["
961  + toString((int) rotationMatrix.size1()) + "x"
962  + toString((int) rotationMatrix.size2()) + "].";
964  }
965 
966  // if given matrix is not a rotation matrix, we can't proceed.
967  if ( !isRotationMatrix(rotationMatrix) ) {
968  QString msg = "Unable to convert the given matrix to an axis of rotation "
969  "and a rotation angle. The given matrix is not a rotation matrix.";
971  }
972 
973  Angle angle;
974  Vector axis(3);
975 
976  // Construct the quaternion corresponding to the input rotation matrix
977  Vector quaternion = LinearAlgebra::toQuaternion(rotationMatrix);
978 
979  LinearAlgebra::Vector subQuaternion(3);
980  subQuaternion[0] = quaternion[1];
981  subQuaternion[1] = quaternion[2];
982  subQuaternion[2] = quaternion[3];
983 
984 
985  // The quaternion we've just constructed is of the form:
986  //
987  // cos(ANGLE/2) + sin(ANGLE/2) * AXIS
988  //
989  // We take a few precautions to handle the case of an identity
990  // rotation.
991  if (isZero(subQuaternion)) {
992  axis(0) = 0.0;
993  axis(1) = 0.0;
994  axis(2) = 1.0;
995 
996  angle.setRadians(0.0);
997  }
998  else if (qFuzzyCompare(quaternion(0) + 1.0, 1.0)) {
999  axis = subQuaternion;
1000  angle.setRadians(Isis::PI);
1001  }
1002  else {
1003  axis = LinearAlgebra::normalize(subQuaternion);
1004  angle.setRadians(2.0 * atan2(magnitude(subQuaternion), quaternion(0)));
1005  }
1006 
1007  return qMakePair(axis, angle);
1008  }
1009 
1010 
1033  // Derived from NAIF's m2eul routine
1035  const QList<int> axes) {
1036 
1037  // check there are 3 axes in the set {1,2,3} with center axis not equal to first or last
1038  if (axes.size() != 3) {
1039  QString msg = "Unable to convert the given matrix to Euler angles. "
1040  "Exactly 3 axis codes are required. The given list has ["
1041  + toString((int) axes.size()) + "] axes.";
1043  }
1044 
1045  QSet<int> validAxes;
1046  validAxes << 1 << 2 << 3;
1047  if (!validAxes.contains(axes[0])
1048  || !validAxes.contains(axes[1])
1049  || !validAxes.contains(axes[2])) {
1050  QString msg = "Unable to convert the given matrix to Euler angles using the given axis codes "
1051  "[" + toString(axes[0]) + ", " + toString(axes[1]) + ", " + toString(axes[2])
1052  + "]. Axis codes must be 1, 2, or 3.";
1054  }
1055 
1056  if (axes[0] == axes[1] || axes[1] == axes[2]) {
1057  QString msg = "Unable to convert the given matrix to Euler angles using the given axis codes "
1058  "[" + toString(axes[0]) + ", " + toString(axes[1]) + ", " + toString(axes[2])
1059  + "]. The middle axis code must differ from its neighbors.";
1061  }
1062 
1063  // check the matrix is 3x3 rotation
1064  if ((rotationMatrix.size1() != 3) || (rotationMatrix.size2() != 3)) {
1065  QString msg = "Unable to convert the given matrix to Euler angles. A 3x3 matrix is required. "
1066  "The given matrix is ["
1067  + toString((int) rotationMatrix.size1()) + "x"
1068  + toString((int) rotationMatrix.size2()) + "].";
1070  }
1071 
1072  if ( !isRotationMatrix(rotationMatrix) ) {
1073  QString msg = "Unable to convert the given matrix to Euler angles. "
1074  "The given matrix is not a rotation matrix.";
1076  }
1077 
1078  // AXIS3, AXIS2, AXIS1 and R have passed their tests at this
1079  // point. We take the liberty of working with TEMPROT, a
1080  // version of R that has unitized columns.
1081  // DO I = 1, 3
1082  // CALL VHAT ( R(1,I), TEMPROT(1,I) )
1083  // END DO
1084 
1085  LinearAlgebra::Matrix tempRotation(3, 3);
1086  for (int i = 0; i < 3; i++) {
1087  Vector columnVector = column(rotationMatrix, i);
1088  // replace the ith column with a unitized column vector
1089  setColumn(tempRotation, normalize(columnVector), i);
1090  }
1091 
1092  QList<int> nextAxis;
1093  nextAxis << 2 << 3 << 1;
1094  double sign;
1095  // create a matrix of zeros and we will fill in the non zero components below
1096  LinearAlgebra::Matrix change = zeroMatrix(3, 3);
1097 
1098  Angle angle1, angle2, angle3;
1099  if ( axes[0] == axes[2] ) {
1100  if ( axes[1] == nextAxis[axes[0] - 1] ) {
1101  sign = 1.0;
1102  }
1103  else {
1104  sign = -1.0;
1105  }
1106  int c = 6 - axes[0] - axes[1];
1107 
1108  change( axes[0] - 1, 2 ) = 1.0;
1109  change( axes[1] - 1, 0 ) = 1.0;
1110  change( c - 1, 1 ) = sign * 1.0;
1111  LinearAlgebra::Matrix tempMatrix = multiply( tempRotation, change );
1112  tempRotation = multiply( transpose(change), tempMatrix );
1113 
1114  bool degen = ( qFuzzyCompare(tempRotation(0, 2) + 1.0, 1.0)
1115  && qFuzzyCompare(tempRotation(1, 2) + 1.0, 1.0) )
1116  || ( qFuzzyCompare(tempRotation(2, 0) + 1.0, 1.0)
1117  && qFuzzyCompare(tempRotation(2, 1) + 1.0, 1.0) )
1118  || qFuzzyCompare( qAbs(tempRotation(2, 2)), 1.0 );
1119 
1120  if ( degen ) {
1121 
1122  angle3.setRadians( 0.0 );
1123  angle2.setRadians( acos( tempRotation(2, 2) ) );
1124  angle1.setRadians( atan2( tempRotation(0, 1), tempRotation(0, 0) ) );
1125 
1126  }
1127  else {
1128  // the normal case.
1129  angle3.setRadians( atan2( tempRotation(0, 2), tempRotation(1, 2) ) );
1130  angle2.setRadians( acos( tempRotation(2, 2) ) );
1131  angle1.setRadians( atan2( tempRotation(2, 0), -tempRotation(2, 1) ) );
1132  }
1133  }
1134  else {
1135  // The axis order is c-b-a. We're going to find a matrix CHANGE
1136  // such that
1137  //
1138  // T
1139  // CHANGE R CHANGE
1140  if ( axes[1] == nextAxis[axes[0] - 1] ) {
1141  sign = 1.0;
1142  }
1143  else {
1144  sign = -1.0;
1145  }
1146  change( axes[0] - 1, 0 ) = 1.0;
1147  change( axes[1] - 1, 1 ) = 1.0;
1148  change( axes[2] - 1, 2 ) = sign * 1.0;
1149  LinearAlgebra::Matrix tempMatrix = multiply( tempRotation, change );
1150  tempRotation = multiply( transpose(change), tempMatrix );
1151  bool degen = ( qFuzzyCompare(tempRotation(0, 0) + 1.0, 1.0)
1152  && qFuzzyCompare(tempRotation(0, 1) + 1.0, 1.0) )
1153  || ( qFuzzyCompare(tempRotation(1, 2) + 1.0, 1.0)
1154  && qFuzzyCompare(tempRotation(2, 2) + 1.0, 1.0) )
1155  || qFuzzyCompare( qAbs(tempRotation(0, 2)), 1.0 );
1156 
1157  if ( degen ) {
1158 
1159  angle3.setRadians( 0.0 );
1160  angle2.setRadians( asin( -tempRotation(0, 2) ) );
1161  angle1.setRadians( sign * atan2( -tempRotation(1, 0), tempRotation(1, 1) ) );
1162 
1163  }
1164  else {
1165  // the normal case.
1166  angle3.setRadians( atan2( tempRotation(1, 2), tempRotation(2, 2) ) );
1167  angle2.setRadians( asin( -tempRotation(0, 2) ) );
1168  angle1.setRadians( sign * atan2( tempRotation(0, 1), tempRotation(0, 0) ) );
1169  }
1170  }
1171 
1173  eulerAngles << qMakePair(angle3, axes[0])
1174  << qMakePair(angle2, axes[1])
1175  << qMakePair(angle1, axes[2]);
1176  return eulerAngles;
1177  }
1178 
1179 
1193  // Derived from NAIF's m2q routine
1195 
1196  if ((rotationMatrix.size1() != 3) || (rotationMatrix.size2() != 3)) {
1197  QString msg = "Unable to convert the given matrix to a quaternion. "
1198  "A 3x3 matrix is required. The given matrix is ["
1199  + toString((int) rotationMatrix.size1()) + "x"
1200  + toString((int) rotationMatrix.size2()) + "].";
1202  }
1203 
1204  // if given matrix is not a rotation matrix, we can't proceed.
1205  if ( !isRotationMatrix(rotationMatrix) ) {
1206  QString msg = "Unable to convert the given matrix to an axis of rotation "
1207  "and a rotation angle. The given matrix is not a rotation matrix.";
1209  }
1210 
1211  // Let q = c + s1*i + s2*j + s3*k
1212  // where c is the constant and s1, s2, s3 are the coefficients of the imaginary parts
1213  // and |q| = 1, .
1214 
1215  // Denote the following for n,m in {1,2,3}
1216  // cc = c * c
1217  // csn = c * sn
1218  // snm = sn * sm
1219 
1220  // The rotation matrix corresponding to our quaternion is:
1221  //
1222  // | cc+s11-s22-s33 2*s12-2*cs3 2*s13+2*cs2 |
1223  // | 2*s12+2*cs3 cc-s11+s22-s33 2*s23-2*cs1 |
1224  // | 2*s13-2*cs2 2*s23+2*cs1 cc-s11-s22+s33 |
1225  //
1226  //
1227  // since |q| = cc + s11 + s22 + s33 = 1, we can use substitution on the diagonal entries to get
1228  //
1229  // | 1-2*s22-2*s33 2*s12-2*cs3 2*s13+2*cs2 |
1230  // | 2*s12+2*cs3 1-2*s11-2*s33 2*s23-2*cs1 |
1231  // | 2*s13-2*cs2 2*s23+2*cs1 1-2*s11-2*s22 |
1232  //
1233  //
1234  //
1235  // r(1,1) = 1.0 - 2*s22 - 2*s33
1236  // r(2,1) = 2*s12 + 2*cs3
1237  // r(3,1) = 2*s13 - 2*cs2
1238  //
1239  // r(1,2) = 2*s12 - 2*cs3
1240  // r(2,2) = 1.0 - 2*s11 - 2*s33
1241  // r(3,2) = 2*s23 + 2*cs1
1242  //
1243  // r(1,3) = 2*s13 + 2*cs2
1244  // r(2,3) = 2*s23 - 2*cs1
1245  // r(3,3) = 1.0 - 2*s11 - 2*s22
1246 
1247  // Using this matrix, we get the trace by summing the diagonal entries
1248  // trace = (1-2*s22-2*s33) + (1-2*s11-2*s33) + (1-2*s11-2*s22)
1249  // = 3 - 4*(s11 + s22 + s33)
1250 
1251  // Therefore
1252  // 1.0 + trace = 4 - 4*(s11 + s22 + s33)
1253  // = 4*(1 - s11 - s22 - s33) by factoring 4
1254  // = 4*(|q| - s11 - s22 - s33) since |q| = 1
1255  // = 4*cc since |q| = cc + s11 + s22 + s33
1256  //
1257  // Solving for c, we get that
1258  // 4*cc = 1 + trace
1259  // 2*c = +/- sqrt(1 + trace)
1260  // c = +/- 0.5 * sqrt( 1.0 + trace )
1261 
1262  // We also have the following where {n,m,p} = {1,2,3}
1263  // 1.0 + trace - 2.0*r(n,n) = 4.0 - 4.0(snn + smm + spp) - 2(1.0 - 2.0(smm + spp ))
1264  // = 4.0 - 4.0(snn + smm + spp) - 2.0 + 4.0(smm + spp )
1265  // = 2.0 - 4.0*snn
1266  //
1267  // Solving for snn, we get
1268  // 2.0 - 4.0*snn = 1.0 + trace - 2.0*r(n,n)
1269  // -2.0 + 4.0*snn = -1.0 - trace + 2.0*r(n,n)
1270  // 4.0*snn = 1.0 - trace + 2.0*r(n,n) ******** we will use this value *************
1271  // ******** in our code to minimize *************
1272  // ******** computations *************
1273  // 2.0*sn = +/- sqrt(1.0 - trace + 2.0*r(n,n))
1274  // sn = +/- 0.5 * sqrt(1.0 - trace + 2.0*r(n,n)) for any n in {1,2,3}
1275 
1276 
1277  // In addition to these observations, note that all of the product
1278  // pairs can easily be computed by simplifying the expressions
1279  // (r(n,m) - r(m,n)) / 4.0 and (r(n,m) + r(m,n)) / 4.0
1280  //
1281  // For example,
1282  // (r(3,2) - r(2,3)) / 4.0 = ( (2*s23 + 2*cs1) - (2*s23 - 2*cs1) ) / 4.0
1283  // = ( 4*cs1 ) / 4.0
1284  // = cs1
1285  //
1286  // So we have the following
1287  // cs1 = (r(3,2) - r(2,3))/4.0
1288  // cs2 = (r(1,3) - r(3,1))/4.0
1289  // cs3 = (r(2,1) - r(1,2))/4.0
1290  // s12 = (r(2,1) + r(1,2))/4.0
1291  // s13 = (r(3,1) + r(1,3))/4.0
1292  // s23 = (r(2,3) + r(3,2))/4.0
1293 
1294  //but taking sums or differences of numbers that are nearly equal
1295  //or nearly opposite results in a loss of precision. as a result
1296  //we should take some care in which terms to select when computing
1297  //c, s1, s2, s3. however, by simply starting with one of the
1298  //large quantities cc, s11, s22, or s33 we can make sure that we
1299  //use the best of the 6 quantities above when computing the
1300  //remaining components of the quaternion.
1301 
1302  double trace = rotationMatrix(0, 0) + rotationMatrix(1, 1) + rotationMatrix(2, 2);
1303  double mtrace = 1.0 - trace;
1304 
1305  // cc4 = 4 * c * c where c is the constant term (first) of the quaternion
1306  double cc4 = 1.0 + trace;
1307  // snn4 = 4 * sn * sn where sn is the coefficient for the n+1 term of the quaternion vector
1308  double s114 = mtrace + 2.0*rotationMatrix(0, 0);
1309  double s224 = mtrace + 2.0*rotationMatrix(1, 1);
1310  double s334 = mtrace + 2.0*rotationMatrix(2, 2);
1311 
1312  // Note that cc4 + s114 + s224 + s334 = 4|q| = 4
1313  // Thus, at least one of the 4 terms is greater than 1.
1314  double normalizingFactor;
1315  LinearAlgebra::Vector quaternion(4);
1316  if ( cc4 >= 1.0 ) { // true if the trace is non-negative
1317 
1318  quaternion[0] = sqrt( cc4 * 0.25 ); // note that q0 = c = sqrt(4*c*c*0.25)
1319  normalizingFactor = 1.0 / ( quaternion[0] * 4.0 );
1320 
1321  // multiply the difference of the entries symmetric about the diagonal by the factor 1/4c
1322  // to get s1, s2, and s3
1323  // (see calculations for cs1, cs2, cs3 in comments above)
1324  quaternion[1] = ( rotationMatrix(2, 1) - rotationMatrix(1, 2) ) * normalizingFactor;
1325  quaternion[2] = ( rotationMatrix(0, 2) - rotationMatrix(2, 0) ) * normalizingFactor;
1326  quaternion[3] = ( rotationMatrix(1, 0) - rotationMatrix(0, 1) ) * normalizingFactor;
1327 
1328  }
1329  else if ( s114 >= 1.0 ) {
1330 
1331  quaternion[1] = sqrt( s114 * 0.25 ); // note that q1 = s1 = sqrt(4*s1*s1*0.25)
1332  normalizingFactor = 1.0 /( quaternion[1] * 4.0 );
1333 
1334  quaternion[0] = ( rotationMatrix(2, 1) - rotationMatrix(1, 2) ) * normalizingFactor;
1335  quaternion[2] = ( rotationMatrix(0, 1) + rotationMatrix(1, 0) ) * normalizingFactor;
1336  quaternion[3] = ( rotationMatrix(0, 2) + rotationMatrix(2, 0) ) * normalizingFactor;
1337 
1338  }
1339  else if ( s224 >= 1.0 ) {
1340 
1341  quaternion[2] = sqrt( s224 * 0.25 );
1342  normalizingFactor = 1.0 /( quaternion[2] * 4.0 );
1343 
1344  quaternion[0] = ( rotationMatrix(0, 2) - rotationMatrix(2, 0) ) * normalizingFactor;
1345  quaternion[1] = ( rotationMatrix(0, 1) + rotationMatrix(1, 0) ) * normalizingFactor;
1346  quaternion[3] = ( rotationMatrix(1, 2) + rotationMatrix(2, 1) ) * normalizingFactor;
1347  }
1348  else { // s334 >= 1.0
1349 
1350  quaternion[3] = sqrt( s334 * 0.25 );
1351  normalizingFactor = 1.0 /( quaternion[3] * 4.0 );
1352 
1353  quaternion[0] = ( rotationMatrix(1, 0) - rotationMatrix(0, 1) ) * normalizingFactor;
1354  quaternion[1] = ( rotationMatrix(0, 2) + rotationMatrix(2, 0) ) * normalizingFactor;
1355  quaternion[2] = ( rotationMatrix(1, 2) + rotationMatrix(2, 1) ) * normalizingFactor;
1356 
1357  }
1358 
1359  // if the magnitude of this quaternion is not one, we polish it up a bit.
1360  if ( !isUnit(quaternion) ) {
1361  quaternion = normalize(quaternion);
1362  }
1363 
1364  // M2Q always returns a quaternion with scalar part greater than or equal to zero.
1365  if ( quaternion[0] < 0.0 ) {
1366  quaternion *= -1.0;
1367  }
1368 
1369  return quaternion;
1370  }
1371 
1372 
1387  // Derived from NAIF's axisar routine
1389  if (axis.size() != 3) {
1390  QString msg = "Unable to convert the given vector and angle to a rotation matrix. "
1391  "The given vector with size [" + toString((int) axis.size())
1392  + "] is not a 3D axis vector.";
1394  }
1395 
1396  // intialize the matrix with the 3x3 identity.
1397  LinearAlgebra::Matrix rotationMatrix = identity(3);
1398 
1399  // The matrix we want rotates every vector by angle about axis.
1400  // In particular, it does so to our basis vectors. The columns
1401  // of the output matrix are the images of the basis vectors
1402  // under this rotation.
1403  for (unsigned int i = 0; i < axis.size(); i++) {
1404  // rotate the ith column of the matrix (the ith standard basis vector)
1405  // about the given axis using the given angle
1406  Vector stdBasisVector = column(rotationMatrix, i);
1407 
1408  // replace the ith column with the rotated vector
1409  setColumn(rotationMatrix, rotate(stdBasisVector, axis, angle), i);
1410  }
1411 
1412  return rotationMatrix;
1413  }
1414 
1415 
1426  return LinearAlgebra::toMatrix(axisAngle.first, axisAngle.second);
1427  }
1428 
1429 
1444  // Derived from NAIF's eul2m routine
1446  const EulerAngle &angle2,
1447  const EulerAngle &angle1) {
1448  QSet<int> validAxes;
1449  validAxes << 1 << 2 << 3;
1450  if (!validAxes.contains(angle3.second)
1451  || !validAxes.contains(angle2.second)
1452  || !validAxes.contains(angle1.second)) {
1453  QString msg = "Unable to convert the given Euler angles to a matrix using the given axis "
1454  "codes [" + toString(angle3.second) + ", " + toString(angle2.second) + ", "
1455  + toString(angle1.second) + "]. Axis codes must be 1, 2, or 3.";
1457  }
1458 
1459  LinearAlgebra::Matrix m(3, 3);
1460  // implement eul2m
1461  // Calculate the 3x3 rotation matrix generated by a rotation
1462  // of a specified angle about a specified axis. This rotation
1463  // is thought of as rotating the coordinate system.
1464  //
1465  //ROTATE(angle1, axis1);
1466  double sinAngle = sin(angle1.first.radians());
1467  double cosAngle = cos(angle1.first.radians());
1468  // get the index offset based on the axis number corresponding to this angle
1469  double indexOffset = ((angle1.second % 3) + 3) % 3;
1470  QList<int> indexSet;
1471  indexSet << 2 << 0 << 1 << 2 << 0;
1472  int index1 = indexSet[indexOffset + 0];
1473  int index2 = indexSet[indexOffset + 1];
1474  int index3 = indexSet[indexOffset + 2];
1475  m(index1, index1) = 1.0; m(index1, index2) = 0.0; m(index1, index3) = 0.0;
1476  m(index2, index1) = 0.0; m(index2, index2) = cosAngle; m(index2, index3) = sinAngle;
1477  m(index3, index1) = 0.0; m(index3, index2) = -sinAngle; m(index3, index3) = cosAngle;
1478 
1479  //
1480  // C ROTMAT applies a rotation of ANGLE radians about axis IAXIS to a
1481  // C matrix. This rotation is thought of as rotating the coordinate
1482  // C system.
1483  // C
1484  //CALL ROTMAT ( R, ANGLE2, AXIS2, R1 )
1485  sinAngle = sin(angle2.first.radians());
1486  cosAngle = cos(angle2.first.radians());
1487  indexOffset = ((angle2.second % 3) + 3) % 3;
1488  index1 = indexSet[indexOffset + 0];
1489  index2 = indexSet[indexOffset + 1];
1490  index3 = indexSet[indexOffset + 2];
1491 
1492  LinearAlgebra::Matrix tempMatrix(3, 3);
1493  for (int i = 0; i < 3; i++) {
1494  tempMatrix(index1, i) = m(index1, i);
1495  tempMatrix(index2, i) = cosAngle*m(index2, i) + sinAngle*m(index3, i);
1496  tempMatrix(index3, i) = -sinAngle*m(index2, i) + cosAngle*m(index3, i);
1497  }
1498 
1499  //CALL ROTMAT ( tempMatrix, ANGLE3, AXIS3, m )
1500  sinAngle = sin(angle3.first.radians());
1501  cosAngle = cos(angle3.first.radians());
1502  indexOffset = ((angle3.second % 3) + 3) % 3;
1503  index1 = indexSet[indexOffset + 0];
1504  index2 = indexSet[indexOffset + 1];
1505  index3 = indexSet[indexOffset + 2];
1506  for (int i = 0; i < 3; i++) {
1507  m(index1, i) = tempMatrix(index1, i);
1508  m(index2, i) = cosAngle*tempMatrix(index2, i) + sinAngle*tempMatrix(index3, i);
1509  m(index3, i) = -sinAngle*tempMatrix(index2, i) + cosAngle*tempMatrix(index3, i);
1510  }
1511 
1512  return m;
1513  }
1514 
1515 
1529 
1530  if (eulerAngles.size() != 3) {
1531  QString msg = "Unable to convert the given Euler angles to a matrix. "
1532  "Exactly 3 Euler angles are required. The given list has ["
1533  + toString((int) eulerAngles.size()) + "] angles.";
1535  }
1536  return LinearAlgebra::toMatrix(eulerAngles[0], eulerAngles[1], eulerAngles[2]);
1537  }
1538 
1539 
1540 
1554  // Derived from NAIF's q2m routine
1556  if (quaternion.size() != 4) {
1557  QString msg = "Unable to convert the given vector to a rotation matrix. "
1558  "The given vector with [" + toString((int) quaternion.size())
1559  + "] components is not a quaternion.";
1561  }
1562 
1563  Vector q = quaternion;
1564  if ( !isUnit(quaternion) && !isZero(quaternion) ) {
1565  q = normalize(quaternion);
1566  }
1567 
1568  LinearAlgebra::Matrix matrix(3, 3);
1569  // For efficiency, we avoid duplicating calculations where possible.
1570  double q0q1 = q[0] * q[1];
1571  double q0q2 = q[0] * q[2];
1572  double q0q3 = q[0] * q[3];
1573 
1574  double q1q1 = q[1] * q[1];
1575  double q1q2 = q[1] * q[2];
1576  double q1q3 = q[1] * q[3];
1577 
1578  double q2q2 = q[2] * q[2];
1579  double q2q3 = q[2] * q[3];
1580 
1581  double q3q3 = q[3] * q[3];
1582 
1583  matrix(0, 0) = 1.0 - 2.0 * ( q2q2 + q3q3 );
1584  matrix(0, 1) = 2.0 * ( q1q2 - q0q3 );
1585  matrix(0, 2) = 2.0 * ( q1q3 + q0q2 );
1586 
1587  matrix(1, 0) = 2.0 * ( q1q2 + q0q3 );
1588  matrix(1, 1) = 1.0 - 2.0 * ( q1q1 + q3q3 );
1589  matrix(1, 2) = 2.0 * ( q2q3 - q0q1 );
1590 
1591  matrix(2, 0) = 2.0 * ( q1q3 - q0q2 );
1592  matrix(2, 1) = 2.0 * ( q2q3 + q0q1 );
1593  matrix(2, 2) = 1.0 - 2.0 * ( q1q1 + q2q2 );
1594 
1595  return matrix;
1596  }
1597 
1598 
1609  void LinearAlgebra::setRow(Matrix &matrix, const Vector &vector, int rowIndex) {
1610  if ( (rowIndex+1 > (int) matrix.size1()) ) {
1611  QString msg = "Unable to set the matrix row to the given vector. Row index "
1612  + toString(rowIndex) + " is out of bounds. The given matrix only has "
1613  + toString((int) matrix.size1()) + " rows." ;
1615  }
1616  boost::numeric::ublas::row(matrix, rowIndex) = vector;
1617  }
1618 
1619 
1630  void LinearAlgebra::setColumn(Matrix &matrix, const Vector &vector, int columnIndex) {
1631  if ( (columnIndex+1 > (int) matrix.size2()) ) {
1632  QString msg = "Unable to set the matrix column to the given vector. Column index "
1633  + toString(columnIndex) + " is out of bounds. The given matrix only has "
1634  + toString((int) matrix.size1()) + " columns." ;
1636  }
1637  boost::numeric::ublas::column(matrix, columnIndex) = vector;
1638  }
1639 
1640 
1653  LinearAlgebra::Vector LinearAlgebra::row(const Matrix &matrix, int rowIndex) {
1654  if ( (rowIndex+1 > (int) matrix.size1()) ) {
1655  QString msg = "Unable to get the matrix row to the given vector. Row index "
1656  + toString(rowIndex) + " is out of bounds. The given matrix only has "
1657  + toString((int) matrix.size1()) + " rows." ;
1659  }
1660  return boost::numeric::ublas::row(matrix, rowIndex);
1661  }
1662 
1663 
1676  LinearAlgebra::Vector LinearAlgebra::column(const Matrix &matrix, int columnIndex) {
1677  if ( (columnIndex+1 > (int) matrix.size2()) ) {
1678  QString msg = "Unable to get the matrix column to the given vector. Column index "
1679  + toString(columnIndex) + " is out of bounds. The given matrix only has "
1680  + toString((int) matrix.size1()) + " columns." ;
1682  }
1683  return boost::numeric::ublas::column(matrix, columnIndex);
1684  }
1685 
1686 
1697  LinearAlgebra::Vector LinearAlgebra::vector(double v0, double v1, double v2) {
1698  Vector v(3);
1699  v[0] = v0;
1700  v[1] = v1;
1701  v[2] = v2;
1702  return v;
1703  }
1704 
1705 
1717  LinearAlgebra::Vector LinearAlgebra::vector(double v0, double v1, double v2, double v3) {
1718  Vector v(4);
1719  v[0] = v0;
1720  v[1] = v1;
1721  v[2] = v2;
1722  v[3] = v3;
1723  return v;
1724  }
1725 
1726 
1735  void LinearAlgebra::setVec3(Vector *v, double v0, double v1, double v2) {
1736  (*v)[0] = v0;
1737  (*v)[1] = v1;
1738  (*v)[2] = v2;
1739  }
1740 
1741 
1751  void LinearAlgebra::setVec4(Vector *v, double v0, double v1, double v2, double v3) {
1752  (*v)[0] = v0;
1753  (*v)[1] = v1;
1754  (*v)[2] = v2;
1755  (*v)[3] = v3;
1756  }
1757 
1758 
1772  LinearAlgebra::Vector LinearAlgebra::subVector(const Vector &v, int startIndex, int size) {
1773  LinearAlgebra::Vector sub(size);
1774  for (int i = 0; i < size; i++) {
1775  sub[i] = v[i + startIndex];
1776  }
1777  return sub;
1778  }
1779 
1780 
1793  QDebug operator<<(QDebug dbg, const LinearAlgebra::Vector &vector) {
1794  QDebugStateSaver saver(dbg);
1795  dbg.noquote() << toString(vector);
1796  return dbg;
1797  }
1798 
1799 
1811  QDebug operator<<(QDebug dbg, const LinearAlgebra::Matrix &matrix) {
1812  QDebugStateSaver saver(dbg);
1813  for (unsigned int i = 0; i < matrix.size1(); i++) {
1814  dbg.noquote() << " ";
1815  for (unsigned int j = 0; j < matrix.size2(); j++) {
1816  dbg.noquote() << toString(matrix(i, j), 15) << " ";
1817  }
1818  dbg.noquote() << endl;
1819  }
1820  return dbg;
1821  }
1822 
1823 
1835  QString toString(const LinearAlgebra::Vector &vector, int precision) {
1836  QString result = "( ";
1837  for (unsigned int i = 0; i < vector.size(); i++) {
1838  result += toString(vector(i), precision);
1839  if (i != vector.size() - 1) result += ", ";
1840  }
1841  result += " )";
1842  return result;
1843  }
1844 }
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.
Definition: Process.h:32
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.
Definition: Angle.h:256
const double PI
The mathematical constant PI.
Definition: Constants.h:56
static bool isRotationMatrix(const Matrix &matrix)
Determines whether the given matrix is a rotation matrix.
double radians() const
Convert an angle to a double.
Definition: Angle.h:243
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)...
static Vector toQuaternion(const Matrix &rotationMatrix)
Converts a rotation&#39;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.
Definition: IString.cpp:226
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.
Definition: IException.h:162
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&#39;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&#39;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&#39;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.
Definition: IException.h:40
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.
Definition: IException.h:134
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.
Definition: Angle.h:62
Isis exception class.
Definition: IException.h:107
Namespace for ISIS/Bullet specific routines.
Definition: Apollo.h:31
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.
Definition: Hillshade.cpp:308
~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...