EMAN2
transform.h
Go to the documentation of this file.
1/*
2 * Author: Steven Ludtke (sludtke@bcm.edu)
3 * Copyright (c) 2000-2006 Baylor College of Medicine
4 *
5 * This software is issued under a joint BSD/GNU license. You may use the
6 * source code in this file under either license. However, note that the
7 * complete EMAN2 and SPARX software packages have some GPL dependencies,
8 * so you are responsible for compliance with the licenses of these packages
9 * if you opt to use BSD licensing. The warranty disclaimer below holds
10 * in either instance.
11 *
12 * This complete copyright notice must be included in any revised version of the
13 * source code. Additional authorship citations may be added, but existing
14 * author citations must be preserved.
15 *
16 * This program is free software; you can redistribute it and/or modify
17 * it under the terms of the GNU General Public License as published by
18 * the Free Software Foundation; either version 2 of the License, or
19 * (at your option) any later version.
20 *
21 * This program is distributed in the hope that it will be useful,
22 * but WITHOUT ANY WARRANTY; without even the implied warranty of
23 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
24 * GNU General Public License for more details.
25 *
26 * You should have received a copy of the GNU General Public License
27 * along with this program; if not, write to the Free Software
28 * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
29 *
30 * */
31
32
33#ifndef eman__transform_h__
34#define eman__transform_h__ 1
35
36#include "vec3.h"
37#include "emobject.h"
38#include <cstdio>
39// #include <vector>
40// using std::vector;
41//
42// #include <map>
43// using std::map;
44//
45// #include <string>
46// using std::string;
47
48
49
50namespace EMAN
51{
75 class Transform {
76// friend Transform EMAN::operator*(const Transform & M2, const Transform & M1);
77 public:
78 static const float ERR_LIMIT;
79
83 Transform();
84
88 Transform( const Transform& rhs );
89
93 Transform& operator=(const Transform& that);
94
99 bool operator==(const Transform& rhs) const;
100
105 bool operator!=(const Transform& rhs) const;
106
110 Transform(const Dict& d);
111
115 Transform(const float array[12]);
116
120 Transform(const vector<float> array);
121
122
124
129 void set_rotation( const Dict &rotation );
130
136 void set_rotation(const Vec3f & v);
137
141 void rotate_origin(const Transform& by);
142
149 void rotate_origin_newBasis(const Transform& tcs, const float& omega, const float& n1, const float& n2, const float& n3);
150
154 void rotate(const Transform& by);
155
160 Dict get_rotation(const string& euler_type = "eman") const;
161
162
167
172
177
183 void set_params(const Dict& d);
184
191 void set_params_inverse(const Dict& d);
192
197 Dict get_params(const string& euler_type) const;
198
203 Dict get_params_inverse(const string& euler_type) const;
204
205 //=============== set and get post trans =============
211 void set_trans(const float& x, const float& y, const float& z=0);
212
216 inline void set_trans(const Vec3f& v) { set_trans(v[0],v[1],v[2]); }
217
221 inline void set_trans(const Vec2f& v) { set_trans(v[0],v[1]); }
222
226 Vec3f get_trans() const;
227
233 void translate(const float& tx, const float& ty, const float& tz=0);
234
238 inline void translate(const Vec3f& v) { translate(v[0],v[1],v[2]); }
239
243 inline void translate(const Vec2f& v) { translate(v[0],v[1]); }
244
253 void translate_newBasis(const Transform& tcs, const float& tx, const float& ty, const float& tz=0);
254
258 inline void translate(const Transform& tcs, const Vec3f& v) { translate_newBasis(tcs, v[0],v[1],v[2]); }
259
263 Vec2f get_trans_2d() const;
264
265 //================= get pre trans is supported =========
266
271 Vec3f get_pre_trans() const;
272
277 Vec2f get_pre_trans_2d() const;
278
283 template<typename type>
284 void set_pre_trans(const type& v);
285
286 //=============== set and get scale =============
290 void set_scale(const float& scale);
291
295 float get_scale() const;
296
300 void scale(const float& scale);
301
302 //=============== set and get post x mirror =============
306 bool get_mirror() const;
307
311 void set_mirror(const bool x_mirror);
312
313 //=============== other stuff ============================
319 void get_scale_and_mirror(float& scale, bool& x_mirror) const;
320
323 void to_identity();
324
327 bool is_identity() const;
328
331 bool is_rot_identity() const;
332
339 void orthogonalize();
340
344 float get_determinant() const;
345
348 void printme() const {
349 printf("%8.6f %8.6f %8.6f %8.6f\n",matrix[0][0],matrix[0][1],matrix[0][2],matrix[0][3]);
350 printf("%8.6f %8.6f %8.6f %8.6f\n",matrix[1][0],matrix[1][1],matrix[1][2],matrix[1][3]);
351 printf("%8.6f %8.6f %8.6f %8.6f\n",matrix[2][0],matrix[2][1],matrix[2][2],matrix[2][3]);
352 printf("%8.6f %8.6f %8.6f %8.6f\n",0.0,0.0,0.0,1.0);
353
354 }
355
356 //=============== get set matrix ============================
360 void set_matrix(const vector<float>& v);
361
362 // returns a Python/JSON string representation of the matrix form
363 string get_matrix_string(int precision);
364
368 vector<float> get_matrix() const;
369
373 vector<float> get_matrix_4x4() const;
374
378 void invert();
379
383 Transform inverse() const;
384
387 void transpose_inplace();
388
392 Transform transpose() const;
393
396 inline float at(int r,int c) const { return matrix[r][c]; }
397
400 inline void set(int r, int c, float value) { matrix[r][c] = value; }
401
405 inline float * operator[] (int i) { return matrix[i]; }
406
410 inline const float * operator[] (int i) const { return matrix[i]; }
411
417 inline Vec2f transform(const float& x, const float& y) const {
418// assert_valid_2d();
419 Vec2f ret;
420 ret[0] = matrix[0][0]*x + matrix[0][1]*y + matrix[0][3];
421 ret[1] = matrix[1][0]*x + matrix[1][1]*y + matrix[1][3];
422 return ret;
423 }
424
429 template<typename Type>
430 inline Vec2f transform(const Vec2<Type>& v) const {
431 return transform(v[0],v[1]);
432 }
433
440 inline Vec3f transform(const float& x, const float& y, const float& z) const {
441// assert_consistent_type(THREED);
442 Vec3f ret;
443 ret[0] = matrix[0][0] * x + matrix[0][1] * y + matrix[0][2] * z + matrix[0][3];
444 ret[1] = matrix[1][0] * x + matrix[1][1] * y + matrix[1][2] * z + matrix[1][3];
445 ret[2] = matrix[2][0] * x + matrix[2][1] * y + matrix[2][2] * z + matrix[2][3];
446 return ret;
447 }
448
453 template<typename Type>
454 inline Vec3f transform(const Vec3<Type>& v) const {
455// assert_consistent_type(THREED); // Transform does the assertion
456 return transform(v[0],v[1],v[2]);
457 }
458
459
465 inline Vec3f get_matrix3_row(int i) const {
466 return Vec3f(matrix[i][0], matrix[i][1], matrix[i][2]);
467 }
468
471 static int get_nsym(const string & sym);
472
477 Transform get_sym(const string & sym, int n) const;
478 // The next two are used solely in sparx, mainly to deal with tet.
479 Transform get_sym_sparx(const string & sym, int n) const;
480 vector<Transform > get_sym_proj(const string & sym) const;
481
484 void copy_matrix_into_array(float* const) const;
485
490 Transform negate() const;
491
495 static Transform tet_3_to_2();
496
500 static Transform icos_5_to_2();
501
502 private:
503 float matrix[3][4];
504
505 void assert_valid_2d() const;
506
508 static vector<string> permissable_2d_not_rot;
509 static vector<string> permissable_3d_not_rot;
510 static map<string,vector<string> > permissable_rot_keys;
511
515
523 void detect_problem_keys(const Dict& d);
524
525 };
527 Transform operator*(const Transform & M2, const Transform & M1);
528
530 template<typename Type>
531 Vec3f operator*( const Transform& M, const Vec3<Type> & v)
532 {
533 return M.transform(v);
534 }
535
537 template<typename Type>
538 Vec2f operator*( const Transform& M, const Vec2<Type> & v)
539 {
540 return M.transform(v);
541 }
542
546 template<typename Type>
547 Vec3f operator*(const Vec3<Type> & v, const Transform & M)
548 {
549 float x = v[0] * M[0][0] + v[1] * M[1][0] + v[2] * M[2][0] ;
550 float y = v[0] * M[0][1] + v[1] * M[1][1] + v[2] * M[2][1];
551 float z = v[0] * M[0][2] + v[1] * M[1][2] + v[2] * M[2][2];
552 return Vec3f(x, y, z);
553 }
554
555 template<typename type>
556 void Transform::set_pre_trans(const type& v) {
557
558 Transform tmp;
559 Dict rot = get_rotation("eman");
560 tmp.set_rotation(rot);
561
562 float scale = get_scale();
563 if (scale != 1.0 ) tmp.set_scale(scale);
564
565 Transform trans;
566 trans.set_trans(v);
567
568 trans = tmp*trans;
569
570 Transform tmp2;
571 tmp2.set_rotation(rot);
572 tmp2.invert(); // invert
573 if (scale != 1.0 ) tmp2.set_scale(1.0f/scale);
574
575
576 trans = trans*tmp2;
577
578 set_trans(trans.get_trans());
579 }
580
624// class Transform3D
625// {
626// public:
627// static const float ERR_LIMIT;
628// enum EulerType
629// {
630// UNKNOWN,
631// EMAN,
632// IMAGIC,
633// SPIN,
634// QUATERNION,
635// SGIROT,
636// SPIDER,
637// MRC,
638// XYZ,
639// MATRIX
640// };
641//
642// /** Default constructor
643// * Internal matrix is the identity
644// */
645// Transform3D();
646//
647// /** Copy constructor
648// * @param rhs the object to be copied
649// */
650// Transform3D( const Transform3D& rhs );
651//
652// /** Construct a Transform3D object describing a rotation, assuming the EMAN Euler type
653// * @param az EMAN - az
654// * @param alt EMAN - alt
655// * @param phi EMAN - phi
656// */
657// Transform3D(const float& az,const float& alt,const float& phi); // EMAN by default
658//
659// /** Construct a Transform3D object describing a rotation (assuming the EMAN Euler type) and
660// * a post translation
661// * @param az EMAN - az
662// * @param alt EMAN - alt
663// * @param phi EMAN - phi
664// * @param posttrans the post translation vector
665// */
666// Transform3D(const float& az, const float& alt, const float& phi, const Vec3f& posttrans);
667//
668// /** Construct a Transform3D object describing a pre trans, a rotation
669// * assuming the EMAN Euler type) and a post translation
670// * @param pretrans the pre translation vector
671// * @param az EMAN - az
672// * @param alt EMAN - alt
673// * @param phi EMAN - phi
674// * @param posttrans the post translation vector
675// */
676// Transform3D(const Vec3f & pretrans, const float& az,const float& alt,const float& phi, const Vec3f& posttrans);
677//
678// /** Construct a Transform3D object describing a rotation, using a specific Euler type.
679// * works for EMAN, SPIDER, MRC, IMAGIC and XYZ (Euler types required 3 parameters)
680// * @param euler_type the Euler type either EMAN, SPIDER, MRC, IMAGIC and XYZ
681// * @param a1 EMAN - az, SPIDER - phi, MRC - phi, IMAGIC - alpha, XYZ - xtilt
682// * @param a2 EMAN - alt, SPIDER - theta, MRC - theta, IMAGIC - beta, XYZ - ytilt
683// * @param a3 EMAN - phi, SPIDER - psi, MRC - omega, IMAGIC - gamma, XYZ - ztilt
684// */
685// Transform3D(EulerType euler_type, const float& a1, const float& a2, const float& a3) ;
686//
687// /** Construct a Transform3D object describing a rotation, using a specific Euler type.
688// * Works for Euler types that require 4 parameters
689// * @param euler_type the Euler type either QUATERNION, SPIN or SGIROT
690// * @param a1 QUATERNION - e0, SPN and SGIROT - omega
691// * @param a2 QUATERNION - e1, SPN and SGIROT - n1
692// * @param a3 QUATERNION - e2, SPN and SGIROT - n2
693// * @param a4 QUATERNION - e3, SPN and SGIROT - n3
694// */
695// Transform3D(EulerType euler_type, const float& a1, const float& a2, const float& a3, const float& a4) ; // o
696//
697// /** Construct a Transform3D object consisting only of a rotation, using a specific Euler type.
698// * Works for all Euler types
699// * @param euler_type any Euler type
700// * @param rotation a dictionary containing all key-entry pair required of the associated Euler type
701// */
702// Transform3D(EulerType euler_type, const Dict& rotation);
703//
704// /** Construct a Transform3D object consisting only of a rotation by initializing the internal
705// * rotation matrix component wise
706// * @param mij the value to be placed in the internal transformation matrix at coordinate (i-1,j-1)
707// */
708// Transform3D(const float& m11, const float& m12, const float& m13,
709// const float& m21, const float& m22, const float& m23,
710// const float& m31, const float& m32, const float& m33);
711//
712// /** Destructor
713// */
714// ~Transform3D();
715//
716// /** FIXME insert comments
717// */
718// void apply_scale(const float& scale);
719//
720// /** FIXME insert comments
721// */
722// void set_scale(const float& scale);
723//
724// /** Reorthogonalize the matrix
725// */
726// void orthogonalize();
727//
728// /** create the transpose in place
729// */
730// void transpose();
731//
732// void set_rotation(const float& az, const float& alt,const float& phi);
733//
734// /** Sets the rotation as defined by the EulerType
735// * works for EMAN, SPIDER, MRC, IMAGIC and XYZ
736// * @param euler_type the Euler type either EMAN, SPIDER, MRC, IMAGIC and XYZ
737// * @param a1 EMAN - az, SPIDER - phi, MRC - phi, IMAGIC - alpha, XYZ - xtilt
738// * @param a2 EMAN - alt, SPIDER - theta, MRC - theta, IMAGIC - beta, XYZ - ytilt
739// * @param a3 EMAN - phi, SPIDER - psi, MRC - omega, IMAGIC - gamma, XYZ - ztilt
740// */
741// void set_rotation(EulerType euler_type,const float& a1,const float& a2,const float& a3);
742//
743// /** Set quaternion-based rotations
744// * Works for QUATERNION, SPIN and SGIROT
745// * @param euler_type the Euler type either QUATERNION, SPIN or SGIROT
746// * @param a1 QUATERNION - e0, SPN and SGIROT - omega
747// * @param a2 QUATERNION - e1, SPN and SGIROT - n1
748// * @param a3 QUATERNION - e2, SPN and SGIROT - n2
749// * @param a4 QUATERNION - e3, SPN and SGIROT - n3
750// */
751// void set_rotation(EulerType euler_type,const float& a1,const float& a2,const float& a3, const float& a4);
752//
753// /** set the internal rotation matrix component wise
754// * @param mij the value to be placed in the internal transformation matrix at coordinate (i-1,j-1)
755// */
756// void set_rotation(const float& m11, const float& m12, const float& m13,
757// const float& m21, const float& m22, const float& m23,
758// const float& m31, const float& m32, const float& m33);
759//
763//
764// /** Set a rotation using a specific Euler type and the dictionary interface
765// * Works for all Euler types
766// * @param euler_type any Euler type
767// * @param rotation a dictionary containing all key-entry pair required of the associated Euler type
768// */
769// void set_rotation(EulerType euler_type, const Dict &rotation );
770//
771// /** Get a rotation in any Euler format
772// * @param euler_type the requested Euler type
773// * @return a dictionary containing the key-entry pairs describing the rotations in terms of the requested Euler type
774// */
775// Dict get_rotation(EulerType euler_type=EMAN) const;
776//
777// /** returns a rotation that maps a pair of unit vectors, a,b to a second pair A,B
778// * @param eahat, ebhat, eAhat, eBhat are all unit vectors
779// * @return a transform3D rotation
780// */
781// void set_rotation(const Vec3f & eahat, const Vec3f & ebhat,
782// const Vec3f & eAhat, const Vec3f & eBhat);
783// /** returns the magnitude of the rotation
784// */
785// float get_mag() const;
786//
787// /** returns the spin-axis (or finger) of the rotation
788// */
789// Vec3f get_finger() const;
790//
791// /** Gets one of two pre translation vectors
792// * @param flag if 0 returns the pre translation vector, if 1 all translation is treated as pre
793// * @return the translation vector
794// */
795// Vec3f get_pretrans( int flag=0) const; // flag=1 => all trans is pre
796//
797// /** Gets one of two post translation vectors.
798// * when the flag is 1 then the contents of the Transform3D matrix right column are returned
799// * @param flag if 0 returns the post translation vector, if 1 all translation is treated as post
800// * @return the translation vector
801// */
802// Vec3f get_posttrans(int flag=0) const; // flag=1 => all trans is post
803//
804// /** Get the total translation as a post translation. Calls get_postrans(1)
805// * @return the translation vector
806// */
807// Vec3f get_total_posttrans() const;
808//
809// /** Get the total translation as a pre translation. Calls get_pretrans(1)
810// * @return the translation vector
811// */
812// Vec3f get_total_pretrans() const;
813//
814// /** This doesn't do anything, it returns an empty vector. Why is it being used?
815// */
816// Vec3f get_center() const; // This doesn't do anything
817//
818// /** Get a matrix column as a Vec3f
819// * @param i the column number (starting at 0)
820// * @return the ith column
821// */
822// Vec3f get_matrix3_col(int i) const;
823//
824// /** Get a matrix row as a Vec3f
825// * @param i the row number (starting at 0)
826// * @return the ith row
827// */
828// Vec3f get_matrix3_row(int i) const;
829//
830// /** Perform a full transform a Vec3f using the internal transformation matrix
831// * @param v3f the vector to be transformed
832// * @return the transformed vector
833// */
834// Vec3f transform(const Vec3f & v3f) const;
835//
836// /** Rotate a Vec3f using the internal rotation matrix
837// * @param v3f the vector to be rotated
838// * @return the rotated vector
839// */
840// Vec3f rotate(const Vec3f & v3f) const;
841//
842// /** FIXME insert comments
843// */
844// Transform3D inverseUsingAngs() const;
845//
846// /** FIXME insert comments
847// */
848// Transform3D inverse() const;
849//
850//
851// /** Print the Transform3D matrix
852// */
853// void printme() const {
854// for (int i=0; i<3; i++) {
855// printf("%6.15f\t%6.15f\t%6.15f\t%6.1f\n",
856// matrix[i][0],matrix[i][1],matrix[i][2],matrix[i][3]);
857// }
858// printf("%6.3f\t%6.3f\t%6.3f\t%6.3f\n",0.0,0.0,0.0,1.0);
859// printf("\n");
860// }
861//
862// /** Get the value stored in the internal transformation matrix at at coordinate (r,c)
863// */
864// inline float at(int r,int c) const { return matrix[r][c]; }
865//
866// /** Set the value stored in the internal transformation matrix at at coordinate (r,c) to value
867// */
868// void set(int r, int c, float value) { matrix[r][c] = value; }
869//
870// /** Operator[] convenience
871// * so Transform3D[2][2] etc terminology can be used
872// */
873// inline float * operator[] (int i) { return matrix[i]; }
874//
875// /** Operator[] convenience
876// * so Transform3D[2][2] etc terminology can be used
877// */
878// inline const float * operator[] (int i) const { return matrix[i]; }
879//
880// static int get_nsym(const string & sym);
881// Transform3D get_sym(const string & sym, int n) const;
882//
883// /** Set functions FIXME insert more comments from here down
884// */
885// void set_center(const Vec3f & center);
886// void set_pretrans(const Vec3f & pretrans);
887// void set_pretrans(const float& dx, const float& dy, const float& dz);
888// void set_pretrans(const float& dx, const float& dy);
889// void set_pretrans(const Vec2f& pretrans);
890// void set_posttrans(const Vec3f & posttrans);
891// void set_posttrans(const float& dx, const float& dy, const float& dz);
892// void set_posttrans(const float& dx, const float& dy);
893// void set_posttrans(const Vec2f& posttrans);
894//
895// void set_post_x_mirror(const bool b) { post_x_mirror = b; }
896// bool get_post_x_mirror() const { return post_x_mirror; }
897//
898// float get_scale() const;
899//
900// void to_identity();
901// bool is_identity();
902//
903// /** Convert a list of euler angles to a vector of Transform3D objects.
904// *
905// * @param[in] eulertype The type of Euler angles that is being passed in.
906// * @param[in] angles A flat vector of angles.
907// *
908// * @return Vector of pointers to Transform3D objects.
909// */
910// static vector<Transform3D*>
911// angles2tfvec(EulerType eulertype, const vector<float> angles);
912//
913// /** This added by d.woolford, will eventually be removed by author
914// */
915// void dsaw_zero_hack(){
916// for (int j=0; j<4; ++j) {
917// for (int i=0; i<4; i++) {
918// if ( fabs(matrix[j][i]) < 0.000001 )
919// matrix[j][i] = 0.0;
920// }
921// }
922//
923// }
924//
925// protected:
926// enum SymType
927// { CSYM,
928// DSYM,
929// TET_SYM,
930// ICOS_SYM,
931// OCT_SYM,
932// ISYM,
933// UNKNOWN_SYM
934// };
935//
936// void init();
937//
938// static SymType get_sym_type(const string & symname);
939//
940// float matrix[4][4];
941//
942// bool post_x_mirror;
943//
944// Transform3D::EulerType s;
945// }; // ends Class
946//
947// Transform3D operator*(const Transform3D & M1, const Transform3D & M2);
948// Vec3f operator*(const Vec3f & v , const Transform3D & M);
949// Vec3f operator*(const Transform3D & M, const Vec3f & v );
950
951// template<typename Type>
952// Vec3f operator*(const Vec3<Type> & v, const Transform3D & M) // YYY
953// {
955// float x = v[0] * M[0][0] + v[1] * M[1][0] + v[2] * M[2][0] ;
956// float y = v[0] * M[0][1] + v[1] * M[1][1] + v[2] * M[2][1];
957// float z = v[0] * M[0][2] + v[1] * M[1][2] + v[2] * M[2][2];
958// return Vec3f(x, y, z);
959// }
960//
961// template<typename Type>
962// Vec3f operator*( const Transform3D & M, const Vec3<Type> & v) // YYY
963// {
965// float x = M[0][0] * v[0] + M[0][1] * v[1] + M[0][2] * v[2] + M[0][3];
966// float y = M[1][0] * v[0] + M[1][1] * v[1] + M[1][2] * v[2] + M[1][3];
967// float z = M[2][0] * v[0] + M[2][1] * v[1] + M[2][2] * v[2] + M[2][3];
968// return Vec3f(x, y, z);
969// }
970//
971//
972// template<typename Type>
973// Vec2f operator*( const Transform3D & M, const Vec2<Type> & v) // YYY
974// {
976// float x = M[0][0] * v[0] + M[0][1] * v[1] + M[0][3] ;
977// float y = M[1][0] * v[0] + M[1][1] * v[1] + M[1][3];
978// return Vec2f(x, y);
979// }
980
981} // ends NameSpace EMAN
982
983
984
985#endif
Dict is a dictionary to store <string, EMObject> pair.
Definition: emobject.h:385
A Transform object is a somewhat specialized object designed specifically for EMAN2/Sparx storage of ...
Definition: transform.h:75
Dict get_params(const string &euler_type) const
Get the parameters of the entire transform, using a specific euler convention.
Definition: transform.cpp:479
Transform get_hflip_transform() const
How do I get the transform that will yield the horizontally flipped projection?
Definition: transform.cpp:792
bool operator!=(const Transform &rhs) const
Unequality comparision operator.
Definition: transform.cpp:128
float at(int r, int c) const
Get the value stored in the internal transformation matrix at at coordinate (r,c)
Definition: transform.h:396
bool get_mirror() const
Query whether x_mirroring is occurring.
Definition: transform.cpp:1250
Vec2f transform(const float &x, const float &y) const
Transform 2D coordinates using the internal transformation matrix.
Definition: transform.h:417
void copy_matrix_into_array(float *const) const
Definition: transform.cpp:158
void rotate_origin_newBasis(const Transform &tcs, const float &omega, const float &n1, const float &n2, const float &n3)
Increment the rotation by multipling the rotation bit of the argument transfrom by the rotation part ...
Definition: transform.cpp:738
Transform get_sym(const string &sym, int n) const
Apply the symmetry deduced from the function arguments to this Transform and return the result.
Definition: transform.cpp:1389
bool is_rot_identity() const
Returns whethers or this matrix rotation is the identity.
Definition: transform.cpp:238
void translate(const Transform &tcs, const Vec3f &v)
Increment the current translation using vec3f& v and a non standard basis.
Definition: transform.h:258
void assert_valid_2d() const
Definition: transform.cpp:1366
float matrix[3][4]
Definition: transform.h:503
Vec3f get_matrix3_row(int i) const
Get a matrix row as a Vec3f required for back compatibility with Tranform3D - see PawelProjector.
Definition: transform.h:465
static const float ERR_LIMIT
Definition: transform.h:78
void translate(const Vec2f &v)
Increment the current translation using vec2f& v.
Definition: transform.h:243
Vec2f get_pre_trans_2d() const
2D version of getting the translation vector as though this object was MSRT_ not MTSR,...
Definition: transform.cpp:1111
Vec2f transform(const Vec2< Type > &v) const
Transform a 2D vector using the internal transformation matrix.
Definition: transform.h:430
string get_matrix_string(int precision)
Definition: transform.cpp:169
Transform get_rotation_transform() const
Get the rotation part of the tranformation matrix as a Transform object.
Definition: transform.cpp:688
bool operator==(const Transform &rhs) const
Equality comparision operator.
Definition: transform.cpp:119
Vec3f transform(const Vec3< Type > &v) const
Transform a 3D vector using the internal transformation matrix.
Definition: transform.h:454
void set_mirror(const bool x_mirror)
Set whether or not x_mirroring is occurring.
Definition: transform.cpp:1238
void scale(const float &scale)
Increment the scale.
Definition: transform.cpp:1159
void set_rotation(const Dict &rotation)
Set a rotation using a specific Euler type and the dictionary interface Works for all Euler types.
Definition: transform.cpp:519
Dict get_params_inverse(const string &euler_type) const
Get the parameters of the inverse of the transform as though it were in RSMT order not MTSR.
Definition: transform.cpp:499
Vec3f get_pre_trans() const
Get the translation vector as though this object was MSRT_ not MTSR, where T_ is what you want Note M...
Definition: transform.cpp:1100
static Transform icos_5_to_2()
Get the transform that moves any icosahedron generated by eman2 so that it matches the 2-2-2 (MRC,...
Definition: transform.cpp:63
Transform get_vflip_transform() const
How do I get the transform that will yield the vertically flipped projection?
Definition: transform.cpp:811
static vector< string > permissable_3d_not_rot
Definition: transform.h:509
vector< Transform > get_sym_proj(const string &sym) const
Definition: transform.cpp:1407
Dict get_rotation(const string &euler_type="eman") const
Get a rotation in any Euler format.
Definition: transform.cpp:829
Transform get_sym_sparx(const string &sym, int n) const
Definition: transform.cpp:1398
void translate(const float &tx, const float &ty, const float &tz=0)
Increment the current translation by tx, ty, tz.
Definition: transform.cpp:1063
static map< string, vector< string > > permissable_rot_keys
Definition: transform.h:510
Vec3f get_trans() const
Get the post trans as a vec3f.
Definition: transform.cpp:1046
void to_identity()
Force the internal matrix to become the identity.
Definition: transform.cpp:207
static vector< string > permissable_2d_not_rot
This map is used to validate keys in the argument maps - e.g. if the type is 2d and the angle is not ...
Definition: transform.h:508
void printme() const
Print the contents of the internal matrix verbatim to standard out.
Definition: transform.h:348
void set(int r, int c, float value)
Set the value stored in the internal transformation matrix at at coordinate (r,c) to value.
Definition: transform.h:400
void detect_problem_keys(const Dict &d)
Test to ensure the parametes in the given dictionary are valid Throws if an error is detected Generic...
Definition: transform.cpp:369
vector< float > get_matrix_4x4() const
Get the 4x4 transformation matrix using a vector.
Definition: transform.cpp:192
void translate(const Vec3f &v)
Increment the current translation using vec3f& v.
Definition: transform.h:238
Transform negate() const
Negates the Transform - a useful way, for example, for getting an orientation on the opposite side of...
Definition: transform.cpp:781
void set_params(const Dict &d)
Set the parameters of the entire transform.
Definition: transform.cpp:254
Transform()
Default constructor Internal matrix is the identity.
Definition: transform.cpp:100
void rotate(const Transform &by)
Increment the rotation by multipling the rotation bit of the argument transfrom by the current transf...
Definition: transform.cpp:763
void set_trans(const Vec3f &v)
Set the post translation component using a Vec3f.
Definition: transform.h:216
float * operator[](int i)
Operator[] convenience so Transform3D[2][2] etc terminology can be used.
Definition: transform.h:405
Vec3f transform(const float &x, const float &y, const float &z) const
Transform 3D coordinates using the internal transformation matrix.
Definition: transform.h:440
float get_scale() const
Get the scale that was applied.
Definition: transform.cpp:1145
Transform inverse() const
Get the inverse of this transformation matrix.
Definition: transform.cpp:1327
void get_scale_and_mirror(float &scale, bool &x_mirror) const
Get scale and x_mirror with 1 function call.
Definition: transform.cpp:1260
void rotate_origin(const Transform &by)
Increment the rotation by multipling the rotation bit of the argument transfrom by the rotation part ...
Definition: transform.cpp:720
void set_scale(const float &scale)
Set the scale.
Definition: transform.cpp:1123
void set_pre_trans(const type &v)
Set the translational component of the matrix as though it was MSRT_ not MTSR, where T_ is the pre tr...
Definition: transform.h:556
Transform & operator=(const Transform &that)
Assignment operator.
Definition: transform.cpp:110
void set_matrix(const vector< float > &v)
Set the transformation matrix using a vector.
Definition: transform.cpp:147
float get_determinant() const
Get the determinant of the matrix.
Definition: transform.cpp:1279
bool is_identity() const
Returns whethers or this matrix is the identity.
Definition: transform.cpp:222
void orthogonalize()
Reorthogonalize the rotation part of the matrix in place.
Definition: transform.cpp:1179
void set_trans(const float &x, const float &y, const float &z=0)
Set the post translation component.
Definition: transform.cpp:1036
void init_permissable_keys()
Called internally to initialize permissable_2d_not_rot, permissable_3d_not_rot, and permissable_rot_k...
Definition: transform.cpp:285
void set_trans(const Vec2f &v)
Set the post translation component using a Vec2f.
Definition: transform.h:221
void translate_newBasis(const Transform &tcs, const float &tx, const float &ty, const float &tz=0)
Increment the current translation by tx, ty, tz using a non standard basis Actualy what it does is re...
Definition: transform.cpp:1072
void invert()
Get the inverse of this transformation matrix.
Definition: transform.cpp:1293
vector< float > get_matrix() const
Get the transformation matrix using a vector.
Definition: transform.cpp:181
void transpose_inplace()
Get the transpose of this transformation matrix.
Definition: transform.cpp:1333
Vec2f get_trans_2d() const
Get the degenerant 2D post trans as a vec2f.
Definition: transform.cpp:1088
void set_params_inverse(const Dict &d)
Set the parameters of the entire transform as though they there in the inverse format.
Definition: transform.cpp:425
static Transform tet_3_to_2()
Get the transform that moves any tetrahedron generated by eman2 so that it matches the 2-2-2 (MRC,...
Definition: transform.cpp:85
static int get_nsym(const string &sym)
get the number of symmetries associated with the given symmetry name
Definition: transform.cpp:1570
Transform transpose() const
Get the transpose of this transformation matrix.
Definition: transform.cpp:1346
The Vec2 is precisely the same as Vec3 except it works exclusively in 2D Note there are convenient ty...
Definition: vec3.h:710
E2Exception class.
Definition: aligner.h:40
Vec3< float > Vec3f
Definition: vec3.h:693
EMData * operator*(const EMData &em, float n)
Definition: emdata.cpp:3201
#define y(i, j)
Definition: projector.cpp:1516
#define x(i)
Definition: projector.cpp:1517