EMAN2
Public Member Functions | Static Public Member Functions | Private Attributes | List of all members
EMAN::Matrix4 Class Reference

#include <vecmath.h>

Public Member Functions

 Matrix4 ()
 
 Matrix4 (const Vector4 &row0, const Vector4 &row1, const Vector4 &row2, const Vector4 &row3)
 
 Matrix4 (const Vector3 &row0, const Vector3 &row1, const Vector3 &row2)
 
 Matrix4 (const Matrix4 &m)
 
Matrix4operator= (const Matrix4 &m)
 
int index (int row, int col) const
 
const double & operator() (int row, int col) const
 
double & operator() (int row, int col)
 
Vector4 row (int r) const
 
Vector4 column (int c) const
 
Matrix4operator*= (const Matrix4 &m)
 
Matrix4operator+= (const Matrix4 &m)
 
Matrix4operator-= (const Matrix4 &m)
 
Matrix4 transpose () const
 
Matrix4 operator+ (const Matrix4 &m) const
 
Matrix4 operator- (const Matrix4 &m) const
 
Vector3 operator* (const Vector3 &v) const
 
Point3 operator* (const Point3 &p) const
 
Vector4 operator* (const Vector4 &v) const
 
Matrix4 operator* (const Matrix4 &b) const
 
Matrix4 inverse () const
 
bool operator== (const Matrix4 &m) const
 
bool approxEqual (const Matrix4 &m, double eps=1e-12) const
 
void print () const
 

Static Public Member Functions

static Matrix4 identity ()
 
static Matrix4 translation (const Point3 &p)
 
static Matrix4 translation (const Vector3 &v)
 
static Matrix4 rotation (const Vector3 &u, const Vector3 &v, const Vector3 &w)
 
static Matrix4 rotation (const Vector3 &axis, double angle)
 
static Matrix4 xrotation (double angle)
 
static Matrix4 yrotation (double angle)
 
static Matrix4 zrotation (double angle)
 
static Matrix4 scaling (const Vector3 &s)
 
static Matrix4 scaling (double x, double y, double z)
 
static Matrix4 scaling (double scale)
 

Private Attributes

double mat [16]
 

Detailed Description

Definition at line 691 of file vecmath.h.

Constructor & Destructor Documentation

◆ Matrix4() [1/4]

EMAN::Matrix4::Matrix4 ( )
inline

Definition at line 693 of file vecmath.h.

693 {
694 for ( int i = 0; i < 4; i++ )
695 for ( int j = 0; j < 4; j++ )
696 mat[ index(i,j) ] = (i == j) ? 1.0 : 0.0;
697 }
int index(int row, int col) const
Definition: vecmath.h:729
double mat[16]
Definition: vecmath.h:938

References index(), and mat.

Referenced by identity(), rotation(), scaling(), translation(), xrotation(), yrotation(), and zrotation().

◆ Matrix4() [2/4]

EMAN::Matrix4::Matrix4 ( const Vector4 row0,
const Vector4 row1,
const Vector4 row2,
const Vector4 row3 
)
inline

Definition at line 699 of file vecmath.h.

699 {
700 for ( int i = 0; i < 4; i++ ) {
701 mat[ index( 0, i ) ] = row0[i];
702 mat[ index( 1, i ) ] = row1[i];
703 mat[ index( 2, i ) ] = row2[i];
704 mat[ index( 3, i ) ] = row3[i];
705 }
706 }

References index(), and mat.

◆ Matrix4() [3/4]

EMAN::Matrix4::Matrix4 ( const Vector3 row0,
const Vector3 row1,
const Vector3 row2 
)
inline

Definition at line 708 of file vecmath.h.

708 {
709 for ( int i = 0; i < 3; i++ ) {
710 mat[ index( 0, i ) ] = row0[i];
711 mat[ index( 1, i ) ] = row1[i];
712 mat[ index( 2, i ) ] = row2[i];
713 mat[ index(i,3) ] = 0.0;
714 mat[ index(3,i) ] = 0.0;
715 }
716 mat[ index(3,3) ] = 1.0;
717 }

References index(), and mat.

◆ Matrix4() [4/4]

EMAN::Matrix4::Matrix4 ( const Matrix4 m)
inline

Definition at line 719 of file vecmath.h.

719 {
720 (*this) = m;
721 }

Member Function Documentation

◆ approxEqual()

bool EMAN::Matrix4::approxEqual ( const Matrix4 m,
double  eps = 1e-12 
) const
inline

Definition at line 825 of file vecmath.h.

825 {
826 for ( int i = 0; i < 16; i++ )
827 if ( isZero( mat[i] - m.mat[i], eps ) )
828 return false;
829 return true;
830 }
bool isZero(double in_d, double in_dEps=1e-16)
Definition: vecmath.h:44

References EMAN::isZero(), and mat.

◆ column()

Vector4 EMAN::Matrix4::column ( int  c) const
inline

Definition at line 738 of file vecmath.h.

738 {
739 return Vector4( mat[index(0,c)], mat[index(1,c)], mat[index(2,c)], mat[index(3,c)] );
740 }

References index(), and mat.

◆ identity()

static Matrix4 EMAN::Matrix4::identity ( )
inlinestatic

Definition at line 839 of file vecmath.h.

839 {
840 return Matrix4(Vector4(1, 0, 0, 0),
841 Vector4(0, 1, 0, 0),
842 Vector4(0, 0, 1, 0),
843 Vector4(0, 0, 0, 1));
844 }

References Matrix4().

Referenced by inverse().

◆ index()

int EMAN::Matrix4::index ( int  row,
int  col 
) const
inline

Definition at line 729 of file vecmath.h.

729{ Assert( row >= 0 && row < 4 ); Assert( col >= 0 && col < 4 ); return col * 4 + row; }
Vector4 row(int r) const
Definition: vecmath.h:734
#define Assert(s)
Define Assert() function that is effective only when -DDEBUG is used.
Definition: emassert.h:42

References Assert, and row().

Referenced by column(), Matrix4(), operator()(), and row().

◆ inverse()

Matrix4 EMAN::Matrix4::inverse ( ) const
inline

Definition at line 956 of file vecmath.h.

956 {
957 // Compute inverse using Gauss-Jordan elimination; caller is responsible
958 // for ensuring that the matrix isn't singular.
959 Matrix4 a(*this);
961 int i, j;
962 int p;
963
964 for (j = 0; j < 4; j++) {
965 p = j;
966 for (i = j + 1; i < 4; i++) {
967 if (fabs(a(i,j)) > fabs(a(p,j)))
968 p = i;
969 }
970 // Swap rows p and j
971 if ( p != j ) {
972 for ( i = 0; i < 4; i++ ) {
973 const double ta = a(p,i);
974 a(p,i) = a(j,i);
975 a(j,i) = ta;
976
977 const double tb = b(p,i);
978 b(p,i) = b(j,i);
979 b(j,i) = tb;
980 }
981 }
982
983 const double s = a(j,j); // if s == 0, the matrix is singular
984 Assert( isZero( s ) == false );
985 for ( i = 0; i < 4; i++ ) {
986 a(j,i) *= (1.0 / s);
987 b(j,i) *= (1.0 / s);
988 }
989 // Eliminate off-diagonal elements
990 for (i = 0; i < 4; i++) {
991 if (i != j) {
992 for ( int k = 0; k < 4; k++ ) {
993 b(i,k) -= a(i,j) * b(j,k);
994 a(i,k) -= a(i,j) * a(j,k);
995 }
996 }
997 }
998 }
999 return b;
1000 }
static Matrix4 identity()
Definition: vecmath.h:839

References Assert, identity(), and EMAN::isZero().

◆ operator()() [1/2]

double & EMAN::Matrix4::operator() ( int  row,
int  col 
)
inline

Definition at line 732 of file vecmath.h.

732{ return mat[ index(row,col) ]; }

References index(), mat, and row().

◆ operator()() [2/2]

const double & EMAN::Matrix4::operator() ( int  row,
int  col 
) const
inline

Definition at line 731 of file vecmath.h.

731{ return mat[ index(row,col) ]; }

References index(), mat, and row().

◆ operator*() [1/4]

Matrix4 EMAN::Matrix4::operator* ( const Matrix4 b) const
inline

Definition at line 804 of file vecmath.h.

804 {
805 Matrix4 matRet;
806 for ( int i = 0; i < 4; i++ ) {
807 for ( int j = 0; j < 4; j++ ) {
808 matRet(i,j) = 0.0;
809 for ( int k = 0; k < 4; k++ )
810 matRet(i,j) += (*this)(i,k) * b(k,j);
811 }
812 }
813 return matRet;
814 }

◆ operator*() [2/4]

Point3 EMAN::Matrix4::operator* ( const Point3 p) const
inline

Definition at line 788 of file vecmath.h.

788 {
789 const Point3 pt((*this)(0,0) * p[0] + (*this)(0,1) * p[1] + (*this)(0,2) * p[2] + (*this)(0,3),
790 (*this)(1,0) * p[0] + (*this)(1,1) * p[1] + (*this)(1,2) * p[2] + (*this)(1,3),
791 (*this)(2,0) * p[0] + (*this)(2,1) * p[1] + (*this)(2,2) * p[2] + (*this)(2,3));
792 const double w = (*this)(3,0) * p[0] + (*this)(3,1) * p[1] + (*this)(3,2) * p[2] + (*this)(3,3);
793 Assert( isZero( w ) == false );
794 return Point3( pt[0] / w, pt[1] / w, pt[2] / w );
795 }

References Assert, and EMAN::isZero().

◆ operator*() [3/4]

Vector3 EMAN::Matrix4::operator* ( const Vector3 v) const
inline

Definition at line 782 of file vecmath.h.

782 {
783 return Vector3((*this)(0,0) * v[0] + (*this)(0,1) * v[1] + (*this)(0,2) * v[2],
784 (*this)(1,0) * v[0] + (*this)(1,1) * v[1] + (*this)(1,2) * v[2],
785 (*this)(2,0) * v[0] + (*this)(2,1) * v[1] + (*this)(2,2) * v[2]);
786 }

◆ operator*() [4/4]

Vector4 EMAN::Matrix4::operator* ( const Vector4 v) const
inline

Definition at line 797 of file vecmath.h.

797 {
798 return Vector4((*this)(0,0) * v[0] + (*this)(0,1) * v[1] + (*this)(0,2) * v[2] + (*this)(0,3) * v[3],
799 (*this)(1,0) * v[0] + (*this)(1,1) * v[1] + (*this)(1,2) * v[2] + (*this)(1,3) * v[3],
800 (*this)(2,0) * v[0] + (*this)(2,1) * v[1] + (*this)(2,2) * v[2] + (*this)(2,3) * v[3],
801 (*this)(3,0) * v[0] + (*this)(3,1) * v[1] + (*this)(3,2) * v[2] + (*this)(3,3) * v[3]);
802 }

◆ operator*=()

Matrix4 & EMAN::Matrix4::operator*= ( const Matrix4 m)
inline

Definition at line 742 of file vecmath.h.

742 {
743 const Matrix4 matRet = (*this) * m;
744 (*this) = matRet;
745 return *this;
746 }

◆ operator+()

Matrix4 EMAN::Matrix4::operator+ ( const Matrix4 m) const
inline

Definition at line 768 of file vecmath.h.

768 {
769 Matrix4 matRet;
770 for ( int i = 0; i < 16; i++ )
771 matRet.mat[i] = mat[i] + m.mat[i];
772 return matRet;
773 }

References mat.

◆ operator+=()

Matrix4 & EMAN::Matrix4::operator+= ( const Matrix4 m)
inline

Definition at line 748 of file vecmath.h.

748 {
749 const Matrix4 matRet = (*this) + m;
750 (*this) = matRet;
751 return *this;
752 }

◆ operator-()

Matrix4 EMAN::Matrix4::operator- ( const Matrix4 m) const
inline

Definition at line 775 of file vecmath.h.

775 {
776 Matrix4 matRet;
777 for ( int i = 0; i < 16; i++ )
778 matRet.mat[i] = mat[i] - m.mat[i];
779 return matRet;
780 }

References mat.

◆ operator-=()

Matrix4 & EMAN::Matrix4::operator-= ( const Matrix4 m)
inline

Definition at line 754 of file vecmath.h.

754 {
755 const Matrix4 matRet = (*this) - m;
756 (*this) = matRet;
757 return *this;
758 }

◆ operator=()

Matrix4 & EMAN::Matrix4::operator= ( const Matrix4 m)
inline

Definition at line 723 of file vecmath.h.

723 {
724 memcpy( &mat[0], &m.mat[0], sizeof(double) * 16 );
725 return *this;
726 }

References mat.

◆ operator==()

bool EMAN::Matrix4::operator== ( const Matrix4 m) const
inline

Definition at line 818 of file vecmath.h.

818 {
819 for ( int i = 0; i < 16; i++ )
820 if ( mat[i] != m.mat[i] )
821 return false;
822 return true;
823 }

References mat.

◆ print()

void EMAN::Matrix4::print ( ) const
inline

Definition at line 832 of file vecmath.h.

832 {
833 std::cout << "( " << (*this)(0,0) << ", " << (*this)(0,1) << ", " << (*this)(0,2) << ", " << (*this)(0,3) << "\n";
834 std::cout << " " << (*this)(1,0) << ", " << (*this)(1,1) << ", " << (*this)(1,2) << ", " << (*this)(1,3) << "\n";
835 std::cout << " " << (*this)(2,0) << ", " << (*this)(2,1) << ", " << (*this)(2,2) << ", " << (*this)(2,3) << "\n";
836 std::cout << " " << (*this)(3,0) << ", " << (*this)(3,1) << ", " << (*this)(3,2) << ", " << (*this)(3,3) << ")\n";
837 }

◆ rotation() [1/2]

static Matrix4 EMAN::Matrix4::rotation ( const Vector3 axis,
double  angle 
)
inlinestatic

Definition at line 867 of file vecmath.h.

867 {
868 Vector3 a = axis;
869 a.normalize();
870 const double c = cos(angle);
871 const double s = sin(angle);
872 const double t = 1 - c;
873
874 return Matrix4(Vector4(t * a[0] * a[0] + c,
875 t * a[0] * a[1] - s * a[2],
876 t * a[0] * a[2] + s * a[1],
877 0),
878 Vector4(t * a[0] * a[1] + s * a[2],
879 t * a[1] * a[1] + c,
880 t * a[1] * a[2] - s * a[0],
881 0),
882 Vector4(t * a[0] * a[2] - s * a[1],
883 t * a[1] * a[2] + s * a[0],
884 t * a[2] * a[2] + c,
885 0),
886 Vector4(0, 0, 0, 1));
887 }

References Matrix4(), and EMAN::Vector3::normalize().

◆ rotation() [2/2]

static Matrix4 EMAN::Matrix4::rotation ( const Vector3 u,
const Vector3 v,
const Vector3 w 
)
inlinestatic

Definition at line 860 of file vecmath.h.

860 {
861 return Matrix4(Vector4(u[0], u[1], u[2], 0),
862 Vector4(v[0], v[1], v[2], 0),
863 Vector4(w[0], w[1], w[2], 0),
864 Vector4(0 , 0 , 0 , 1));
865 }

References Matrix4().

◆ row()

Vector4 EMAN::Matrix4::row ( int  r) const
inline

Definition at line 734 of file vecmath.h.

734 {
735 return Vector4( mat[index(r,0)], mat[index(r,1)], mat[index(r,2)], mat[index(r,3)] );
736 }

References index(), and mat.

Referenced by index(), operator()(), and EMAN::operator<<().

◆ scaling() [1/3]

static Matrix4 EMAN::Matrix4::scaling ( const Vector3 s)
inlinestatic

Definition at line 919 of file vecmath.h.

919 {
920 return Matrix4(Vector4(s[0], 0 , 0 , 0),
921 Vector4(0 , s[1], 0 , 0),
922 Vector4(0 , 0 , s[2], 0),
923 Vector4(0 , 0 , 0 , 1));
924 }

References Matrix4().

Referenced by scaling().

◆ scaling() [2/3]

static Matrix4 EMAN::Matrix4::scaling ( double  scale)
inlinestatic

Definition at line 933 of file vecmath.h.

933 {
934 return scaling(Vector3(scale, scale, scale));
935 }
static Matrix4 scaling(const Vector3 &s)
Definition: vecmath.h:919

References scaling().

◆ scaling() [3/3]

static Matrix4 EMAN::Matrix4::scaling ( double  x,
double  y,
double  z 
)
inlinestatic

Definition at line 926 of file vecmath.h.

926 {
927 return Matrix4(Vector4(x, 0, 0, 0),
928 Vector4(0, y, 0, 0),
929 Vector4(0, 0, z, 0),
930 Vector4(0, 0, 0, 1));
931 }
#define y(i, j)
Definition: projector.cpp:1516
#define x(i)
Definition: projector.cpp:1517

References Matrix4(), x, and y.

◆ translation() [1/2]

static Matrix4 EMAN::Matrix4::translation ( const Point3 p)
inlinestatic

Definition at line 846 of file vecmath.h.

846 {
847 return Matrix4(Vector4(1, 0, 0, p[0]),
848 Vector4(0, 1, 0, p[1]),
849 Vector4(0, 0, 1, p[2]),
850 Vector4(0, 0, 0, 1));
851 }

References Matrix4().

◆ translation() [2/2]

static Matrix4 EMAN::Matrix4::translation ( const Vector3 v)
inlinestatic

Definition at line 853 of file vecmath.h.

853 {
854 return Matrix4(Vector4(1, 0, 0, v[0]),
855 Vector4(0, 1, 0, v[1]),
856 Vector4(0, 0, 1, v[2]),
857 Vector4(0, 0, 0, 1));
858 }

References Matrix4().

◆ transpose()

Matrix4 EMAN::Matrix4::transpose ( ) const
inline

Definition at line 760 of file vecmath.h.

760 {
761 Matrix4 matRet;
762 for ( int i = 0; i < 4; i++ )
763 for ( int j = 0; j < 4; j++ )
764 matRet(i,j) = (*this)(j,i);
765 return matRet;
766 }

◆ xrotation()

static Matrix4 EMAN::Matrix4::xrotation ( double  angle)
inlinestatic

Definition at line 889 of file vecmath.h.

889 {
890 const double c = cos(angle);
891 const double s = sin(angle);
892
893 return Matrix4(Vector4(1, 0, 0, 0),
894 Vector4(0, c, -s, 0),
895 Vector4(0, s, c, 0),
896 Vector4(0, 0, 0, 1));
897 }

References Matrix4().

◆ yrotation()

static Matrix4 EMAN::Matrix4::yrotation ( double  angle)
inlinestatic

Definition at line 899 of file vecmath.h.

899 {
900 double c = cos(angle);
901 double s = sin(angle);
902
903 return Matrix4(Vector4( c, 0, s, 0),
904 Vector4( 0, 1, 0, 0),
905 Vector4(-s, 0, c, 0),
906 Vector4( 0, 0, 0, 1));
907 }

References Matrix4().

◆ zrotation()

static Matrix4 EMAN::Matrix4::zrotation ( double  angle)
inlinestatic

Definition at line 909 of file vecmath.h.

909 {
910 const double c = cos(angle);
911 const double s = sin(angle);
912
913 return Matrix4(Vector4(c, -s, 0, 0),
914 Vector4(s, c, 0, 0),
915 Vector4(0, 0, 1, 0),
916 Vector4(0, 0, 0, 1));
917 }

References Matrix4().

Member Data Documentation

◆ mat

double EMAN::Matrix4::mat[16]
private

The documentation for this class was generated from the following file: