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

A Transform object is a somewhat specialized object designed specifically for EMAN2/Sparx storage of alignment parameters and euler orientations. More...

#include <transform.h>

Public Member Functions

 Transform ()
 Default constructor Internal matrix is the identity. More...
 
 Transform (const Transform &rhs)
 Copy constructor. More...
 
Transformoperator= (const Transform &that)
 Assignment operator. More...
 
bool operator== (const Transform &rhs) const
 Equality comparision operator. More...
 
bool operator!= (const Transform &rhs) const
 Unequality comparision operator. More...
 
 Transform (const Dict &d)
 Construction using a dictionary. More...
 
 Transform (const float array[12])
 Construction using an array of floats. More...
 
 Transform (const vector< float > array)
 Construction using a vector of size 12. More...
 
 ~Transform ()
 
void set_rotation (const Dict &rotation)
 Set a rotation using a specific Euler type and the dictionary interface Works for all Euler types. More...
 
void set_rotation (const Vec3f &v)
 Determine the rotation that would transform a vector pointing in the Z direction so that it points in the direction of the argument vector Automatically normalizes the vector. More...
 
void rotate_origin (const Transform &by)
 Increment the rotation by multipling the rotation bit of the argument transfrom by the rotation part of the current transfrom. More...
 
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 of the current transfrom This version rotates in the standard coordinate system, even after it have been modified by tcs. More...
 
void rotate (const Transform &by)
 Increment the rotation by multipling the rotation bit of the argument transfrom by the current transfrom. More...
 
Dict get_rotation (const string &euler_type="eman") const
 Get a rotation in any Euler format. More...
 
Transform get_rotation_transform () const
 Get the rotation part of the tranformation matrix as a Transform object. More...
 
Transform get_hflip_transform () const
 How do I get the transform that will yield the horizontally flipped projection? More...
 
Transform get_vflip_transform () const
 How do I get the transform that will yield the vertically flipped projection? More...
 
void set_params (const Dict &d)
 Set the parameters of the entire transform. More...
 
void set_params_inverse (const Dict &d)
 Set the parameters of the entire transform as though they there in the inverse format. More...
 
Dict get_params (const string &euler_type) const
 Get the parameters of the entire transform, using a specific euler convention. More...
 
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. More...
 
void set_trans (const float &x, const float &y, const float &z=0)
 Set the post translation component. More...
 
void set_trans (const Vec3f &v)
 Set the post translation component using a Vec3f. More...
 
void set_trans (const Vec2f &v)
 Set the post translation component using a Vec2f. More...
 
Vec3f get_trans () const
 Get the post trans as a vec3f. More...
 
void translate (const float &tx, const float &ty, const float &tz=0)
 Increment the current translation by tx, ty, tz. More...
 
void translate (const Vec3f &v)
 Increment the current translation using vec3f& v. More...
 
void translate (const Vec2f &v)
 Increment the current translation using vec2f& v. More...
 
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 remove the effect of tcs when a composite transfrom tcs*t (where t is the current transform) This function is used in the scenegraph. More...
 
void translate (const Transform &tcs, const Vec3f &v)
 Increment the current translation using vec3f& v and a non standard basis. More...
 
Vec2f get_trans_2d () const
 Get the degenerant 2D post trans as a vec2f. More...
 
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 means post x mirror, T means translation, S means scale, and R means rotaiton. More...
 
Vec2f get_pre_trans_2d () const
 2D version of getting the translation vector as though this object was MSRT_ not MTSR, where T_ is what you want Note M means post x mirror, T means translation, S means scale, and R means rotation More...
 
template<typename type >
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 translation. More...
 
void set_scale (const float &scale)
 Set the scale. More...
 
float get_scale () const
 Get the scale that was applied. More...
 
void scale (const float &scale)
 Increment the scale. More...
 
bool get_mirror () const
 Query whether x_mirroring is occurring. More...
 
void set_mirror (const bool x_mirror)
 Set whether or not x_mirroring is occurring. More...
 
void get_scale_and_mirror (float &scale, bool &x_mirror) const
 Get scale and x_mirror with 1 function call. More...
 
void to_identity ()
 Force the internal matrix to become the identity. More...
 
bool is_identity () const
 Returns whethers or this matrix is the identity. More...
 
bool is_rot_identity () const
 Returns whethers or this matrix rotation is the identity. More...
 
void orthogonalize ()
 Reorthogonalize the rotation part of the matrix in place. More...
 
float get_determinant () const
 Get the determinant of the matrix. More...
 
void printme () const
 Print the contents of the internal matrix verbatim to standard out. More...
 
void set_matrix (const vector< float > &v)
 Set the transformation matrix using a vector. More...
 
string get_matrix_string (int precision)
 
vector< float > get_matrix () const
 Get the transformation matrix using a vector. More...
 
vector< float > get_matrix_4x4 () const
 Get the 4x4 transformation matrix using a vector. More...
 
void invert ()
 Get the inverse of this transformation matrix. More...
 
Transform inverse () const
 Get the inverse of this transformation matrix. More...
 
void transpose_inplace ()
 Get the transpose of this transformation matrix. More...
 
Transform transpose () const
 Get the transpose of this transformation matrix. More...
 
float at (int r, int c) const
 Get the value stored in the internal transformation matrix at at coordinate (r,c) More...
 
void set (int r, int c, float value)
 Set the value stored in the internal transformation matrix at at coordinate (r,c) to value. More...
 
float * operator[] (int i)
 Operator[] convenience so Transform3D[2][2] etc terminology can be used. More...
 
const float * operator[] (int i) const
 Operator[] convenience so Transform3D[2][2] etc terminology can be used. More...
 
Vec2f transform (const float &x, const float &y) const
 Transform 2D coordinates using the internal transformation matrix. More...
 
template<typename Type >
Vec2f transform (const Vec2< Type > &v) const
 Transform a 2D vector using the internal transformation matrix. More...
 
Vec3f transform (const float &x, const float &y, const float &z) const
 Transform 3D coordinates using the internal transformation matrix. More...
 
template<typename Type >
Vec3f transform (const Vec3< Type > &v) const
 Transform a 3D vector using the internal transformation matrix. More...
 
Vec3f get_matrix3_row (int i) const
 Get a matrix row as a Vec3f required for back compatibility with Tranform3D - see PawelProjector. More...
 
Transform get_sym (const string &sym, int n) const
 Apply the symmetry deduced from the function arguments to this Transform and return the result. More...
 
Transform get_sym_sparx (const string &sym, int n) const
 
vector< Transformget_sym_proj (const string &sym) const
 
void copy_matrix_into_array (float *const) const
 
Transform negate () const
 Negates the Transform - a useful way, for example, for getting an orientation on the opposite side of the sphere. More...
 

Static Public Member Functions

static int get_nsym (const string &sym)
 get the number of symmetries associated with the given symmetry name More...
 
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, FREALIGN) convention. More...
 
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, FREALIGN) convention. More...
 

Static Public Attributes

static const float ERR_LIMIT = 0.000001f
 

Private Member Functions

void assert_valid_2d () const
 
void init_permissable_keys ()
 Called internally to initialize permissable_2d_not_rot, permissable_3d_not_rot, and permissable_rot_keys static members. More...
 
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 - works in every circumstance (set_params, set_rotation, set_params_inv) Uses static members permissable_2d_not_rot, permissable_3d_not_rot, and permissable_rot_keys as basis of decision. More...
 

Private Attributes

float matrix [3][4]
 

Static Private Attributes

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 "alpha" then we should throw. More...
 
static vector< string > permissable_3d_not_rot
 
static map< string, vector< string > > permissable_rot_keys
 

Detailed Description

A Transform object is a somewhat specialized object designed specifically for EMAN2/Sparx storage of alignment parameters and euler orientations.

It's designed to store four transformations in a specific order, namely Transform = MTSR Where M is a mirroring operation (about the x-axis) or the identity T is a Translation matrix S is a uniform scaling matrix R is a rotation matrix This means you can call set_scale, set_trans, set_rotation in any order but still have the operations arranged internally in the order of MTSR. This is somewhat restrictive, for example in the context of how OpenGL handles transformations, but in practice is nicely suited to the situations that arise in EMAN2 - namely, alignment and projection orientation characterization.

Note that you can fool the Transform object into storing any matrix by using the constructors that take array arguments. This can useful, for example, for shearing your image.

See http://blake.bcm.tmc.edu/emanwiki/Eman2TransformInPython for using it from Python and detailed discussion See test_transform.py for examples of the way it is unit tested See http://blake.bcm.tmc.edu/emanwiki/EMAN2/Tutorials/RotateTranslate for examples showing how to transform EMDatas with it.

Author
David Woolford (based on Philip Baldwin's original Transform3D code)
Date
September 2008

Definition at line 75 of file transform.h.

Constructor & Destructor Documentation

◆ Transform() [1/5]

Transform::Transform ( )

Default constructor Internal matrix is the identity.

Definition at line 100 of file transform.cpp.

101{
102 to_identity();
103}
void to_identity()
Force the internal matrix to become the identity.
Definition: transform.cpp:207

References to_identity().

Referenced by rotate_origin_newBasis(), set_params_inverse(), and translate_newBasis().

◆ Transform() [2/5]

Transform::Transform ( const Transform rhs)

Copy constructor.

Parameters
rhsthe object to be copied

Definition at line 105 of file transform.cpp.

106{
107 *this = that;
108}

◆ Transform() [3/5]

Transform::Transform ( const Dict d)

Construction using a dictionary.

Parameters
dthe dictionary containing the parameters

Definition at line 132 of file transform.cpp.

132 {
133 to_identity();
134 set_params(d);
135}
void set_params(const Dict &d)
Set the parameters of the entire transform.
Definition: transform.cpp:254

References set_params(), and to_identity().

◆ Transform() [4/5]

Transform::Transform ( const float  array[12])

Construction using an array of floats.

Parameters
arraythe array of values that will become the internal matrix. row order (3 rows of 4)

Definition at line 138 of file transform.cpp.

138 {
139 memcpy(matrix,array,12*sizeof(float));
140}
float matrix[3][4]
Definition: transform.h:503

References matrix.

◆ Transform() [5/5]

Transform::Transform ( const vector< float >  array)

Construction using a vector of size 12.

Parameters
arraythe array of values that will become the internal matrix. row order (3 rows of 4)

Definition at line 142 of file transform.cpp.

143{
144 set_matrix(array);
145}
void set_matrix(const vector< float > &v)
Set the transformation matrix using a vector.
Definition: transform.cpp:147

References set_matrix().

◆ ~Transform()

EMAN::Transform::~Transform ( )
inline

Definition at line 123 of file transform.h.

123{ }

Member Function Documentation

◆ assert_valid_2d()

void Transform::assert_valid_2d ( ) const
private

Definition at line 1366 of file transform.cpp.

1366 {
1367 int rotation_error = 0;
1368 int translation_error = 0;
1369 if (fabs(matrix[2][0]) > ERR_LIMIT) rotation_error++;
1370 if (fabs(matrix[2][1]) > ERR_LIMIT) rotation_error++;
1371 if (fabs(matrix[2][3]) > ERR_LIMIT) translation_error++;
1372 if (fabs(matrix[0][2]) > ERR_LIMIT) rotation_error++;
1373 if (fabs(matrix[1][2]) > ERR_LIMIT) rotation_error++;
1374// if (fabs(matrix[2][2]-1.0) >ERR_LIMIT) rotation_error++;
1375 if (matrix[2][2] <=0) rotation_error++; // previous line commented out due to addition of scaling. This line will insure we don't have alt=180.0
1376 if ( translation_error && rotation_error ) {
1377 throw UnexpectedBehaviorException("Error, the internal matrix contains 3D rotations and 3D translations. This object can not be considered 2D");
1378 } else if ( translation_error ) {
1379 throw UnexpectedBehaviorException("Error, the internal matrix contains a non zero z component for a 3D translation. This object can not be considered 2D");
1380 }
1381 else if ( rotation_error ) {
1382 throw UnexpectedBehaviorException("Error, the internal matrix contains 3D rotations and this object can not be considered 2D");
1383 }
1384
1385}
static const float ERR_LIMIT
Definition: transform.h:78
#define UnexpectedBehaviorException(desc)
Definition: exception.h:400

References ERR_LIMIT, matrix, and UnexpectedBehaviorException.

Referenced by get_rotation(), and set_rotation().

◆ at()

float EMAN::Transform::at ( int  r,
int  c 
) const
inline

Get the value stored in the internal transformation matrix at at coordinate (r,c)

Definition at line 396 of file transform.h.

396{ return matrix[r][c]; }

References matrix.

Referenced by EMAN::PawelProjector::backproject3d().

◆ copy_matrix_into_array()

void Transform::copy_matrix_into_array ( float * const  array) const

Definition at line 158 of file transform.cpp.

158 {
159
160 int idx = 0;
161 for(int i=0; i<3; ++i) {
162 for(int j=0; j<4; ++j) {
163 array[idx] = matrix[i][j];
164 idx ++;
165 }
166 }
167}

References matrix.

Referenced by EMAN::FourierReconstructor::do_compare_slice_work(), EMAN::FourierReconstructor::do_insert_slice_work(), EMAN::TransformProcessor::process(), EMAN::TransformProcessor::process_inplace(), and EMAN::StandardProjector::project3d().

◆ detect_problem_keys()

void Transform::detect_problem_keys ( const Dict d)
private

Test to ensure the parametes in the given dictionary are valid Throws if an error is detected Generic - works in every circumstance (set_params, set_rotation, set_params_inv) Uses static members permissable_2d_not_rot, permissable_3d_not_rot, and permissable_rot_keys as basis of decision.

Parameters
dthe dictionary that was the function argument of the set_params, set_rotation or the set_params_inv function
Exceptions
InvalidParameterExceptionif the dictionary is invalid in anyway

Definition at line 369 of file transform.cpp.

369 {
370 if (permissable_rot_keys.size() == 0 ) {
372 }
373
374 vector<string> verification;
375 vector<string> problem_keys;
376 bool is_2d = false;
377 if (d.has_key_ci("type") ) {
378 string type = Util::str_to_lower((string)d["type"]);
379 bool problem = false;
380 if (permissable_rot_keys.find(type) == permissable_rot_keys.end() ) {
381 problem_keys.push_back(type);
382 problem = true;
383 }
384 if ( !problem ) {
385 vector<string> perm = permissable_rot_keys[type];
386 std::copy(perm.begin(),perm.end(),back_inserter(verification));
387
388 if ( type == "2d" ) {
389 is_2d = true;
390 std::copy(permissable_2d_not_rot.begin(),permissable_2d_not_rot.end(),back_inserter(verification));
391 }
392 }
393 }
394 if ( !is_2d ) {
395 std::copy(permissable_3d_not_rot.begin(),permissable_3d_not_rot.end(),back_inserter(verification));
396 }
397
398 for (Dict::const_iterator it = d.begin(); it != d.end(); ++it) {
399 if ( std::find(verification.begin(),verification.end(), it->first) == verification.end() ) {
400 problem_keys.push_back(it->first);
401 }
402 }
403
404 if (problem_keys.size() != 0 ) {
405 string error;
406 if (problem_keys.size() == 1) {
407 error = "Transform Error: The \"" +problem_keys[0]+ "\" key is unsupported";
408 } else {
409 error = "Transform Error: The ";
410 for(vector<string>::const_iterator cit = problem_keys.begin(); cit != problem_keys.end(); ++cit ) {
411 if ( cit != problem_keys.begin() ) {
412 if (cit == (problem_keys.end() -1) ) error += " and ";
413 else error += ", ";
414 }
415 error += "\"";
416 error += *cit;
417 error += "\"";
418 }
419 error += " keys are unsupported";
420 }
421 throw InvalidParameterException(error);
422 }
423}
Const iterator support for the Dict object This is just a wrapper, everything is inherited from the m...
Definition: emobject.h:674
bool has_key_ci(const string &key) const
Ask the Dictionary if it as a particular key in a case insensitive way.
Definition: emobject.cpp:1140
iterator end()
Definition: emobject.cpp:1061
iterator begin()
Definition: emobject.cpp:1045
static vector< string > permissable_3d_not_rot
Definition: transform.h:509
static map< string, vector< string > > permissable_rot_keys
Definition: transform.h:510
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 init_permissable_keys()
Called internally to initialize permissable_2d_not_rot, permissable_3d_not_rot, and permissable_rot_k...
Definition: transform.cpp:285
static string str_to_lower(const string &s)
Return a lower case version of the argument string.
Definition: util.cpp:283
EMData * copy() const
This file is a part of "emdata.h", to use functions in this file, you should "#include "emdata....
#define InvalidParameterException(desc)
Definition: exception.h:361

References EMAN::Dict::begin(), copy(), EMAN::Dict::end(), EMAN::Dict::has_key_ci(), init_permissable_keys(), InvalidParameterException, permissable_2d_not_rot, permissable_3d_not_rot, permissable_rot_keys, and EMAN::Util::str_to_lower().

Referenced by set_params(), set_params_inverse(), and set_rotation().

◆ get_determinant()

float Transform::get_determinant ( ) const

Get the determinant of the matrix.

Returns
the determinant

Definition at line 1279 of file transform.cpp.

1280{
1281 float det;
1282 double det2;
1283 det2 = matrix[0][0]*((double)matrix[1][1]*matrix[2][2]-(double)matrix[2][1]*matrix[1][2]);
1284 det2 -= matrix[0][1]*((double)matrix[1][0]*matrix[2][2]-(double)matrix[2][0]*matrix[1][2]);
1285 det2 += matrix[0][2]*((double)matrix[1][0]*matrix[2][1]-(double)matrix[2][0]*matrix[1][1]);
1286
1287 det = (float)det2;
1289
1290 return det;
1291}
static void apply_precision(float &value, const float &precision)
Definition: util.h:1309

References EMAN::Util::apply_precision(), ERR_LIMIT, and matrix.

Referenced by get_mirror(), get_scale(), get_scale_and_mirror(), and scale().

◆ get_hflip_transform()

Transform Transform::get_hflip_transform ( ) const

How do I get the transform that will yield the horizontally flipped projection?

Returns
the transform that will yield the horizontally flipped projection

Definition at line 792 of file transform.cpp.

792 {
793
794 Dict rot = get_rotation("eman");
795 rot["alt"] = 180.0f + static_cast<float>(rot["alt"]);
796 rot["phi"] = 180.0f - static_cast<float>(rot["phi"]);
797 // This is the same as new_alt= 180-alt, new_phi=-phi, new_az=180+az
798
799 Transform ret(*this); // Is the identity
800 ret.set_rotation(rot);
801
802 Vec3f trans = get_trans();
803 trans[0] = -trans[0];
804 ret.set_trans(trans);
805
806// ret.set_mirror(self.get_mirror());
807
808 return ret;
809}
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_rotation(const string &euler_type="eman") const
Get a rotation in any Euler format.
Definition: transform.cpp:829
Vec3f get_trans() const
Get the post trans as a vec3f.
Definition: transform.cpp:1046

References get_rotation(), get_trans(), set_rotation(), and set_trans().

◆ get_matrix()

vector< float > Transform::get_matrix ( ) const

Get the transformation matrix using a vector.

Returns
a vector - 3 rows of 4 - that stores the values of the transformation matrix

Definition at line 181 of file transform.cpp.

182{
183 vector<float> ret(12);
184 for(int i=0; i<3; ++i) {
185 for(int j=0; j<4; ++j) {
186 ret[i*4+j] = matrix[i][j];
187 }
188 }
189 return ret;
190}

References matrix.

Referenced by EMAN::EMData::extract_box(), rotate(), and rotate_origin().

◆ get_matrix3_row()

Vec3f EMAN::Transform::get_matrix3_row ( int  i) const
inline

Get a matrix row as a Vec3f required for back compatibility with Tranform3D - see PawelProjector.

Parameters
ithe row number (starting at 0)
Returns
the ith row

Definition at line 465 of file transform.h.

465 {
466 return Vec3f(matrix[i][0], matrix[i][1], matrix[i][2]);
467 }
Vec3< float > Vec3f
Definition: vec3.h:693

References matrix.

Referenced by EMAN::PawelProjector::project3d().

◆ get_matrix_4x4()

vector< float > Transform::get_matrix_4x4 ( ) const

Get the 4x4 transformation matrix using a vector.

Returns
a vector - 4 rows of 4 - that stores the values of the transformation matrix

Definition at line 192 of file transform.cpp.

193{
194 vector<float> ret(16);
195 for(int i=0; i<3; ++i) {
196 for(int j=0; j<4; ++j) {
197 ret[i*4+j] = matrix[i][j];
198 }
199 }
200 ret[12] = 0.0;
201 ret[13] = 0.0;
202 ret[14] = 0.0;
203 ret[15] = 1.0;
204
205 return ret;
206}

References matrix.

◆ get_matrix_string()

string Transform::get_matrix_string ( int  precision)

Definition at line 169 of file transform.cpp.

169 {
170 if (precision<2||precision>9) precision=9;
171 char fmt[256];
172 char buf[256]; // sufficient for any potential result
173 sprintf(fmt,"[%%1.%0dg,%%1.%0dg,%%1.%0dg,%%1.%0dg,%%1.%0dg,%%1.%0dg,%%1.%0dg,%%1.%0dg,%%1.%0dg,%%1.%0dg,%%1.%0dg,%%1.%0dg]",
174 precision,precision,precision,precision,precision,precision,precision,precision,precision,precision,precision,precision);
175 sprintf(buf,fmt,matrix[0][0],matrix[0][1],matrix[0][2],matrix[0][3],matrix[1][0],matrix[1][1],matrix[1][2],matrix[1][3],
176 matrix[2][0],matrix[2][1],matrix[2][2],matrix[2][3]);
177 return string(buf);
178}

References matrix.

◆ get_mirror()

bool Transform::get_mirror ( ) const

Query whether x_mirroring is occurring.

Returns
whether x_mirroring is occurring

Definition at line 1250 of file transform.cpp.

1250 {
1251 float determinant = get_determinant();
1252
1253 bool x_mirror = false;
1254 if ( determinant < 0 ) x_mirror = true;
1255
1256 return x_mirror;
1257
1258}
float get_determinant() const
Get the determinant of the matrix.
Definition: transform.cpp:1279

References get_determinant().

Referenced by get_params(), get_params_inverse(), get_trans(), get_trans_2d(), EMAN::GaussFFTProjector::project3d(), set_mirror(), set_trans(), EMAN::TransformProcessor::transform(), translate(), and EMAN::RT2DTreeAligner::xform_align_nbest().

◆ get_nsym()

int Transform::get_nsym ( const string &  sym)
static

get the number of symmetries associated with the given symmetry name

Definition at line 1570 of file transform.cpp.

1571{
1572 Symmetry3D* sym = Factory<Symmetry3D>::get(sym_name);
1573 int nsym = sym->get_nsym();
1574 delete sym;
1575 return nsym;
1576}
static T * get(const string &instance_name)
Definition: emobject.h:781
Symmetry3D - A base class for 3D Symmetry objects.
Definition: symmetry.h:57
virtual int get_nsym() const =0
The total number of unique symmetry operations that will be return by this object when a calling prog...

References EMAN::Factory< T >::get(), and EMAN::Symmetry3D::get_nsym().

Referenced by EMAN::PointArray::set_from(), EMAN::nn4Reconstructor::setup(), EMAN::nn4_rectReconstructor::setup(), EMAN::nnSSNR_Reconstructor::setup(), EMAN::nn4_ctfReconstructor::setup(), EMAN::nnSSNR_ctfReconstructor::setup(), EMAN::nn4_ctfwReconstructor::setup(), EMAN::nn4_ctfwsReconstructor::setup(), and EMAN::nn4_ctf_rectReconstructor::setup().

◆ get_params()

Dict Transform::get_params ( const string &  euler_type) const

Get the parameters of the entire transform, using a specific euler convention.

Parameters
euler_typethe euler type of the retrieved rotation
Returns
a dictionary containing the parameters

Definition at line 479 of file transform.cpp.

479 {
480 Dict params = get_rotation(euler_type);
481
482 Vec3f v = get_trans();
483 params["tx"] = v[0]; params["ty"] = v[1];
484
485 string type = Util::str_to_lower(euler_type);
486 if ( type != "2d") params["tz"] = v[2];
487
488 float scale = get_scale();
489 params["scale"] = scale;
490
491 bool mirror = get_mirror();
492 params["mirror"] = mirror;
493
494 return params;
495}
bool get_mirror() const
Query whether x_mirroring is occurring.
Definition: transform.cpp:1250
void scale(const float &scale)
Increment the scale.
Definition: transform.cpp:1159
float get_scale() const
Get the scale that was applied.
Definition: transform.cpp:1145

References get_mirror(), get_rotation(), get_scale(), get_trans(), scale(), and EMAN::Util::str_to_lower().

Referenced by EMAN::TomoTiltEdgeMaskProcessor::process_inplace(), EMAN::RT3DLocalTreeAligner::testort(), EMAN::RT2Dto3DTreeAligner::testort(), EMAN::RT3DTreeAligner::testort(), EMAN::SpiderIO::write_single_header(), EMAN::RT2Dto3DTreeAligner::xform_align_nbest(), EMAN::RT3DTreeAligner::xform_align_nbest(), and EMAN::RT3DLocalTreeAligner::xform_align_nbest().

◆ get_params_inverse()

Dict Transform::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.

Parameters
euler_typethe euler type of the retrieved rotation
Returns
a dictionary containing the parameters

Definition at line 499 of file transform.cpp.

499 {
500 Transform inv(inverse());
501
502 Dict params = inv.get_rotation(euler_type);
503 Vec3f v = inv.get_pre_trans();
504 params["tx"] = v[0]; params["ty"] = v[1];
505
506 string type = Util::str_to_lower(euler_type);
507 if ( type != "2d") params["tz"] = v[2];
508
509 float scale = inv.get_scale();
510 params["scale"] = scale;
511
512 bool mirror = inv.get_mirror();
513 params["mirror"] = mirror;
514
515 return params;
516}
Transform inverse() const
Get the inverse of this transformation matrix.
Definition: transform.cpp:1327

References get_mirror(), get_pre_trans(), get_rotation(), get_scale(), inverse(), scale(), and EMAN::Util::str_to_lower().

◆ get_pre_trans()

Vec3f Transform::get_pre_trans ( ) const

Get the translation vector as though this object was MSRT_ not MTSR, where T_ is what you want Note M means post x mirror, T means translation, S means scale, and R means rotaiton.

Returns
the pre translation vector

Definition at line 1100 of file transform.cpp.

1101{
1102 Transform T(*this);
1103 T.set_trans(0,0,0);
1104 T.invert();
1105
1106 Transform soln = T*(*this);
1107// soln.printme();
1108 return soln.get_trans();
1109}

References get_trans(), invert(), and set_trans().

Referenced by get_params_inverse().

◆ get_pre_trans_2d()

Vec2f Transform::get_pre_trans_2d ( ) const

2D version of getting the translation vector as though this object was MSRT_ not MTSR, where T_ is what you want Note M means post x mirror, T means translation, S means scale, and R means rotation

Returns
the pre translation vector

Definition at line 1111 of file transform.cpp.

1112{
1113 Transform T(*this);
1114 T.set_trans(0,0,0);
1115 T.invert();
1116
1117 Transform soln = T*(*this);
1118// soln.printme();
1119 return soln.get_trans_2d();
1120}
Vec2f get_trans_2d() const
Get the degenerant 2D post trans as a vec2f.
Definition: transform.cpp:1088

References get_trans_2d(), invert(), and set_trans().

◆ get_rotation()

Dict Transform::get_rotation ( const string &  euler_type = "eman") const

Get a rotation in any Euler format.

Parameters
euler_typethe requested Euler type
Returns
a dictionary containing the key-entry pairs describing the rotations in terms of the requested Euler type

Definition at line 829 of file transform.cpp.

830{
831 Dict result;
832 //printf(" Hello from Transform \n");
833
834 //float max = 1 - ERR_LIMIT;
835 float scale;
836 bool x_mirror;
837 get_scale_and_mirror(scale,x_mirror);
838 if (scale == 0) throw UnexpectedBehaviorException("The determinant of the Transform is 0. This is unexpected.");
839
840 double cosalt = matrix[2][2]/scale;
841 double sinalt = sqrt(matrix[2][0]*matrix[2][0]+matrix[2][1]*matrix[2][1])/scale; // PRB May 2021
842 double x_mirror_scale = (x_mirror ? -1.0f : 1.0f);
843 double inv_scale = 1.0f/scale;
844
845 double az = 0;
846 double alt = 0;
847 double phi = 0;
848 double phiS = 0; // like az (but in SPIDER ZYZ)
849 double psiS = 0; // like phi (but in SPIDER ZYZ)
850
851 // get alt, az, phi in EMAN convention
852
853 if (cosalt >= 1) { // that is, alt close to 0
854 alt = 0;
855 az = 0;
856 phi = (double)EMConsts::rad2deg * atan2(x_mirror_scale*matrix[0][1], x_mirror_scale*matrix[0][0]);
857 } else if (cosalt <= -1) { // that is, alt close to 180
858 alt = 180;
859 az = 0;
860 phi = (double)EMConsts::rad2deg * atan2(-x_mirror_scale*matrix[0][1], x_mirror_scale*matrix[0][0]);
861 } else { // for non exceptional cases: 0 < alt < 180
862
863 az = (double)EMConsts::rad2deg * atan2(scale*matrix[2][0], -scale*matrix[2][1]);
864
865 if (matrix[2][2]==0.0)
866 alt = 90.0;
867 else
868 alt = (double)EMConsts::rad2deg * atan(sqrt((double)matrix[2][0]*matrix[2][0]+(double)matrix[2][1]*matrix[2][1])/fabs(matrix[2][2]));
869
870 if (matrix[2][2] * scale < 0)
871 alt = 180.0f-alt;
872
873 phi = (double)EMConsts::rad2deg * atan2(x_mirror_scale*(double)matrix[0][2], (double)matrix[1][2]);
874
875 } // ends separate cases: alt close to 0, 180, or neither
876
877 //phi = phi-360.0*floor(phi/360.0);
878 //az = az -360.0*floor(az/360.0);
879
880// get phiS, psiS (SPIDER)
881 if (cosalt >= 1) { // that is, alt close to 0
882 phiS = 0;
883 psiS = phi;
884 } else if (cosalt <= -1) { // that is, alt close to 180
885 phiS = 0;
886 psiS = phi + 180.0;
887 } else {
888 phiS = az - 90.0;
889 psiS = phi + 90.0;
890 }
891
892 phiS = phiS-360.0*floor(phiS/360.0);
893 psiS = psiS-360.0*floor(psiS/360.0);
894
895// do some quaternionic stuff here
896 double xtilt = 0;
897 double ytilt = 0;
898 double ztilt = 0;
899
900
901 string type = Util::str_to_lower(euler_type);
902
903 result["type"] = type;
904 if (type == "2d") {
906 result["alpha"] = phi;
907 } else if (type == "eman") {
908// assert_consistent_type(THREED);
909 result["az"] = az;
910 result["alt"] = alt;
911 result["phi"] = phi;
912 } else if (type == "imagic") {
913// assert_consistent_type(THREED);
914 result["alpha"] = az;
915 result["beta"] = alt;
916 result["gamma"] = phi;
917 } else if (type == "spider") {
918// assert_consistent_type(THREED);
919 result["phi"] = phiS; // The first Euler like az
920 result["theta"] = alt;
921 result["psi"] = psiS;
922 } else if (type == "mrc") {
923// assert_consistent_type(THREED);
924 result["phi"] = phiS;
925 result["theta"] = alt;
926 result["omega"] = psiS;
927 } else if (type == "xyz") { // need to double-check these 3 equations ********
928// assert_consistent_type(THREED);
929 xtilt = atan2(-sin(EMConsts::deg2rad*phiS)*sin(EMConsts::deg2rad*alt),cos(EMConsts::deg2rad*alt));
930 ytilt = asin( cos(EMConsts::deg2rad*phiS)*sin(EMConsts::deg2rad*alt));
931 ztilt = psiS*EMConsts::deg2rad - atan2(sin(xtilt), cos(xtilt) *sin(ytilt));
932
933 xtilt *= EMConsts::rad2deg; ytilt *= EMConsts::rad2deg; ztilt *= EMConsts::rad2deg;
934 xtilt = xtilt-360*.0*floor((xtilt+180.0)/360.0);
935 ytilt = ytilt-360*.0*floor((ytilt+180.0)/360.0); //already in range [-90,90] but anyway...
936 ztilt = ztilt-360*.0*floor((ztilt+180.0)/360.0);
937
938 result["xtilt"] = xtilt;
939 result["ytilt"] = ytilt;
940 result["ztilt"] = ztilt;
941 } else if ((type == "quaternion") || (type == "spin") || (type == "sgirot")) {
942
943 // The cosOover2 is also e0
944// double nphi = (az-phi)/2.0;
945// double cosOover2 = cos((az+phi)*EMConsts::deg2rad/2.0) * cos(alt*EMConsts::deg2rad/2.0);
946// printf("%f %f %f",matrix[0][0],matrix[1][1],matrix[2][2]);
947 double traceR = matrix[0][0]+matrix[1][1]+matrix[2][2]; // This should be 1 + 2 cos omega
948 double cosomega = (traceR-1.0)/2.0;
949
950 double A0 =(matrix[1][2]-matrix[2][1])/2.0;
951 double A1 =(matrix[2][0]-matrix[0][2])/2.0;
952 double A2 =(matrix[0][1]-matrix[1][0])/2.0;
953 double sinomega = sqrt(A0*A0+A1*A1+A2*A2);
954
955 if (cosomega> 1.0) cosomega= 1.0;
956 if (cosomega<-1.0) cosomega=-1.0;
957
958 // matrix(x,y)-matrix(y,x) = 2 n_z sin(omega) etc
959 // trace matrix = 1 + 2 cos(omega)
960 double cosOover2= sqrt((1.0 +cosomega)/2.0);
961 //double sinOover2= sinomega/cosOover2/2.0;
962 double sinOover2= sqrt((1.0 -cosomega)/2.0);
963 //double sinomega = 2* sinOover2*cosOover2; // Not a great formula. See above
964
965 double n1 = A0/sinomega;
966 double n2 = A1/sinomega;
967 double n3 = A2/sinomega;
968
969
970
971 if (sinOover2==1) {// This will also mean sinomega=0, omega =pi,
972 n1 = sqrt((matrix[0][0]+1)/2.0) ;
973 n2 = sqrt((matrix[1][1]+1)/2.0) ;
974 n3 = sqrt((matrix[2][2]+1)/2.0) ;
975 if (n1 <= n2){ // n1 is the smallest
976 if (n2<=n3) {
977 n3 = sqrt((matrix[2][2]+1)/2.0) ; n1=matrix[0][2]/n3/2.0; n2=matrix[1][2]/n3/2.0; }
978 else { n2 = sqrt((matrix[1][1]+1)/2.0) ; n1=matrix[0][1]/n2/2.0; n3=matrix[2][1]/n2/2.0; }}
979 else { // n2 is the smallest
980 if (n1<=n3) {
981 n3 = sqrt((matrix[2][2]+1)/2.0) ; n1=matrix[0][2]/n3/2.0; n2=matrix[1][2]/n3/2.0;}
982 else { n1 = sqrt((matrix[0][0]+1)/2.0) ; n2=matrix[1][0]/n1/2.0; n3=matrix[2][0]/n1/2.0;}}
983 }
984
985
986 if (sinOover2==0) {// This will also mean omega =0,
987 n1 = 0;
988 n2 = 0;
989 n3 = 1;
990 }
991
992
993// printf("traceR=%lf,OneMinusCosomega=%lf,sinOover2=%lf,cosOover2=%lf,sinomega=%lf,cosomega=%lf,n3=%lf \n",traceR,1-cosomega,sinOover2,cosOover2,sinomega,cosomega,n3);
994
995 if (type == "quaternion"){
996 result["e0"] = cosOover2 ;
997 result["e1"] = sinOover2 * n1 ;
998 result["e2"] = sinOover2 * n2;
999 result["e3"] = sinOover2 * n3;
1000 }
1001
1002 if (type == "spin"){
1003 double Omega = EMConsts::rad2deg * asin(sinomega);
1004 if (cosomega< 0) { Omega = 180-Omega;}
1005 result["omega"] = Omega; //Changed by PRB
1006 result["n1"] = n1;
1007 result["n2"] = n2;
1008 result["n3"] = n3;
1009 }
1010
1011 if (type == "sgirot"){
1012 result["q"] = EMConsts::rad2deg * acos(cosomega);
1013 result["n1"] = n1;
1014 result["n2"] = n2;
1015 result["n3"] = n3;
1016 }
1017
1018 } else if (type == "matrix") {
1019// assert_consistent_type(THREED);
1020 result["m11"] = x_mirror_scale*matrix[0][0]*inv_scale;
1021 result["m12"] = x_mirror_scale*matrix[0][1]*inv_scale;
1022 result["m13"] = x_mirror_scale*matrix[0][2]*inv_scale;
1023 result["m21"] = matrix[1][0]*inv_scale;
1024 result["m22"] = matrix[1][1]*inv_scale;
1025 result["m23"] = matrix[1][2]*inv_scale;
1026 result["m31"] = matrix[2][0]*inv_scale;
1027 result["m32"] = matrix[2][1]*inv_scale;
1028 result["m33"] = matrix[2][2]*inv_scale;
1029 } else {
1030 throw InvalidStringException(euler_type, "unknown Euler Type");
1031 }
1032
1033 return result;
1034}
static const double rad2deg
Definition: emobject.h:77
static const double deg2rad
Definition: emobject.h:78
void assert_valid_2d() const
Definition: transform.cpp:1366
void get_scale_and_mirror(float &scale, bool &x_mirror) const
Get scale and x_mirror with 1 function call.
Definition: transform.cpp:1260
EMData * sqrt() const
return square root of current image
#define InvalidStringException(str, desc)
Definition: exception.h:306

References assert_valid_2d(), EMAN::EMConsts::deg2rad, get_scale_and_mirror(), InvalidStringException, matrix, EMAN::EMConsts::rad2deg, scale(), sqrt(), EMAN::Util::str_to_lower(), and UnexpectedBehaviorException.

Referenced by EMAN::file_store::add_image(), EMAN::PawelProjector::backproject3d(), EMAN::ChaoProjector::backproject3d(), get_hflip_transform(), get_params(), get_params_inverse(), get_vflip_transform(), EMAN::FourierReconstructorSimple2D::insert_slice(), EMAN::FourierGriddingProjector::project3d(), EMAN::ChaoProjector::project3d(), set_pre_trans(), EMAN::TransformProcessor::transform(), EMAN::RT2Dto3DTreeAligner::xform_align_nbest(), EMAN::RT3DTreeAligner::xform_align_nbest(), EMAN::RT3DLocalTreeAligner::xform_align_nbest(), and EMAN::RT3DSymmetryAligner::xform_align_nbest().

◆ get_rotation_transform()

Transform Transform::get_rotation_transform ( ) const

Get the rotation part of the tranformation matrix as a Transform object.

Returns
a Transform describing the rotation only

Definition at line 688 of file transform.cpp.

689{
690 Transform ret(*this);
691 ret.set_scale(1.0);
692 ret.set_mirror(false);
693 ret.set_trans(0,0,0);
694 //ret.orthogonalize(); // ?
695 return ret;
696}

References set_mirror(), set_scale(), and set_trans().

Referenced by EMAN::GaussFFTProjector::project3d().

◆ get_scale()

float Transform::get_scale ( ) const

Get the scale that was applied.

Returns
the scale factor

Definition at line 1145 of file transform.cpp.

1145 {
1146 float determinant = get_determinant();
1147 if (determinant < 0 ) determinant *= -1;
1148
1149 float scale = std::pow(determinant,1.0f/3.0f);
1150 int int_scale = static_cast<int>(scale);
1151 float scale_residual = scale-static_cast<float>(int_scale);
1152 if ( scale_residual < ERR_LIMIT ) { scale = static_cast<float>(int_scale); };
1153
1155
1156 return scale;
1157}

References EMAN::Util::apply_precision(), ERR_LIMIT, get_determinant(), and scale().

Referenced by get_params(), get_params_inverse(), EMAN::TransformProcessor::process(), EMAN::TransformProcessor::process_inplace(), EMAN::GaussFFTProjector::project3d(), set_pre_trans(), and set_scale().

◆ get_scale_and_mirror()

void Transform::get_scale_and_mirror ( float &  scale,
bool &  x_mirror 
) const

Get scale and x_mirror with 1 function call.

More efficient than calling get_scale and get_x_mirror separately

Parameters
scalea reference to the value that will be assigned the scale value
x_mirrora reference to the value that will be assigned the x_mirror value

Definition at line 1260 of file transform.cpp.

1260 {
1261
1262 float determinant = get_determinant();
1263 x_mirror = false;
1264 if ( determinant < 0 ) {
1265 x_mirror = true;
1266 determinant *= -1;
1267 }
1268 if (determinant != 1 ) {
1269 scale = std::pow(determinant,1.0f/3.0f);
1270 int int_scale = static_cast<int>(scale);
1271 float scale_residual = scale-static_cast<float>(int_scale);
1272 if ( scale_residual < ERR_LIMIT ) { scale = static_cast<float>(int_scale); };
1273 }
1274 else scale = 1;
1275
1277}

References EMAN::Util::apply_precision(), ERR_LIMIT, get_determinant(), and scale().

Referenced by get_rotation(), orthogonalize(), and set_rotation().

◆ get_sym()

Transform 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 at line 1389 of file transform.cpp.

1390{
1391 Symmetry3D* sym = Factory<Symmetry3D>::get(sym_name);
1392 Transform ret;
1393 ret = (*this) * sym->get_sym(n);
1394 delete sym;
1395 return ret;
1396}
virtual Transform get_sym(const int n) const =0
Every Symmetry3D object must provide access to the full set of its symmetry operators via this functi...

References EMAN::Factory< T >::get(), and EMAN::Symmetry3D::get_sym().

Referenced by EMAN::PointArray::set_from().

◆ get_sym_proj()

vector< Transform > Transform::get_sym_proj ( const string &  sym) const

Definition at line 1407 of file transform.cpp.

1408{
1409 vector<Transform> ret;
1410 int nsym;
1411 string lstr = Util::str_to_lower(sym_name);
1412 if( lstr == "tet" ) {
1413 nsym = 12;
1414 static float TET[108] = {
1415 1.0f, 0.0f, 0.0f, 0.0f, 1.0f, 0.0f, 0.0f, 0.0f, 1.0f,
1416 -0.5f, 0.86602540378f, 0.0f, -0.86602540378f, -0.5f, 0.0f, 0.0f, 0.0f, 1.0f,
1417 -0.5f, -0.86602540378f, 0.0f, 0.86602540378f, -0.5f, -0.0f, 0.0f, 0.0f, 1.0f,
1418 -0.16666666667f, 0.86602540378f, -0.47140452079f, 0.28867513459f, 0.5f, 0.81649658093f, 0.94280904158f, 0.0f, -0.33333333333f,
1419 0.33333333333f, 0.0f, 0.94280904158f, 0.0f, -1.0f, 0.0f, 0.94280904158f, 0.0f, -0.33333333333f,
1420 -0.16666666667f, -0.86602540378f, -0.47140452079f, -0.28867513459f, 0.5f, -0.81649658093f, 0.94280904158f, 0.0f, -0.33333333333f,
1421 -0.66666666667f, -0.57735026919f, -0.47140452079f, -0.57735026919f, 0.0f, 0.81649658093f, -0.47140452079f, 0.81649658093f, -0.33333333333f,
1422 -0.16666666667f, 0.28867513459f, 0.94280904158f, 0.86602540378f, 0.5f, 0.0f, -0.47140452079f, 0.81649658093f, -0.33333333333f,
1423 0.83333333333f, 0.28867513459f, -0.47140452079f, -0.28867513459f, -0.5f, -0.81649658093f, -0.47140452079f, 0.81649658093f, -0.33333333333f,
1424 0.83333333333f, -0.28867513459f, -0.47140452079f, 0.28867513459f, -0.5f, 0.81649658093f, -0.47140452079f, -0.81649658093f, -0.33333333333f,
1425 -0.16666666667f, -0.28867513459f, 0.94280904158f, -0.86602540378f, 0.5f, 0.0f, -0.47140452079f, -0.81649658093f, -0.33333333333f,
1426 -0.66666666667f, 0.57735026919f, -0.47140452079f, 0.57735026919f, -0.0f, -0.81649658093f, -0.47140452079f, -0.81649658093f, -0.33333333333f
1427 };
1428
1429 for (int k=0;k<nsym;k++) {
1430 Transform t;
1431 for (int i=0; i<3; i++) {
1432 for (int j=0; j<3; j++) {
1433 t.matrix[i][j] = TET[9*k + i*3 +j];
1434 }
1435 }
1436 //vector<float> z = t.get_matrix();
1437 //for (int i=0; i<12; i++) cout<<z[i]<<endl;
1438 ret.push_back( (*this) * t );
1439 }
1440 } else if( lstr == "oct" ) {
1441 nsym = 24;
1442 static float TET[216] = {
1443 1.0f, 0.0f, 0.0f, 0.0f, 1.0f, 0.0f, 0.0f, 0.0f, 1.0f,
1444 0.0f, 1.0f, 0.0f, -1.0f, 0.0f, 0.0f, 0.0f, 0.0f, 1.0f,
1445 -1.0f, 0.0f, 0.0f, 0.0f, -1.0f, 0.0f, 0.0f, 0.0f, 1.0f,
1446 0.0f, -1.0f, 0.0f, 1.0f, 0.0f, 0.0f, 0.0f, 0.0f, 1.0f,
1447 0.0f, 0.0f, -1.0f, 0.0f, 1.0f, 0.0f, 1.0f, 0.0f, 0.0f,
1448 0.0f, 0.0f, -1.0f, -1.0f, 0.0f, 0.0f, 0.0f, 1.0f, 0.0f,
1449 0.0f, 0.0f, -1.0f, 0.0f, -1.0f, 0.0f, -1.0f, 0.0f, 0.0f,
1450 0.0f, 0.0f, -1.0f, 1.0f, 0.0f, 0.0f, 0.0f, -1.0f, 0.0f,
1451 0.0f, 1.0f, 0.0f, 0.0f, 0.0f, 1.0f, 1.0f, 0.0f, 0.0f,
1452 -1.0f, 0.0f, 0.0f, 0.0f, 0.0f, 1.0f, 0.0f, 1.0f, 0.0f,
1453 0.0f, -1.0f, 0.0f, 0.0f, 0.0f, 1.0f, -1.0f, 0.0f, 0.0f,
1454 1.0f, 0.0f, 0.0f, 0.0f, 0.0f, 1.0f, 0.0f, -1.0f, 0.0f,
1455 0.0f, 0.0f, 1.0f, 0.0f, -1.0f, 0.0f, 1.0f, 0.0f, 0.0f,
1456 0.0f, 0.0f, 1.0f, 1.0f, 0.0f, 0.0f, 0.0f, 1.0f, 0.0f,
1457 0.0f, 0.0f, 1.0f, 0.0f, 1.0f, 0.0f, -1.0f, 0.0f, 0.0f,
1458 0.0f, 0.0f, 1.0f, -1.0f, 0.0f, 0.0f, 0.0f, -1.0f, 0.0f,
1459 0.0f, -1.0f, 0.0f, 0.0f, 0.0f, -1.0f, 1.0f, 0.0f, 0.0f,
1460 1.0f, 0.0f, 0.0f, 0.0f, 0.0f, -1.0f, 0.0f, 1.0f, 0.0f,
1461 0.0f, 1.0f, 0.0f, 0.0f, 0.0f, -1.0f, -1.0f, 0.0f, 0.0f,
1462 -1.0f, 0.0f, 0.0f, 0.0f, 0.0f, -1.0f, 0.0f, -1.0f, 0.0f,
1463 -1.0f, 0.0f, 0.0f, 0.0f, 1.0f, 0.0f, 0.0f, 0.0f, -1.0f,
1464 0.0f, 1.0f, 0.0f, 1.0f, 0.0f, 0.0f, 0.0f, 0.0f, -1.0f,
1465 1.0f, 0.0f, 0.0f, 0.0f, -1.0f, 0.0f, 0.0f, 0.0f, -1.0f,
1466 0.0f, -1.0f, 0.0f, -1.0f, 0.0f, 0.0f, 0.0f, 0.0f, -1.0f
1467 };
1468
1469 for (int k=0;k<nsym;k++) {
1470 Transform t;
1471 for (int i=0; i<3; i++) {
1472 for (int j=0; j<3; j++) {
1473 t.matrix[i][j] = TET[9*k + i*3 +j];
1474 }
1475 }
1476 //vector<float> z = t.get_matrix();
1477 //for (int i=0; i<12; i++) cout<<z[i]<<endl;
1478 ret.push_back( (*this) * t );
1479 }
1480 } else if( lstr == "icos" ) {
1481 nsym = 60;
1482 static float TET[540] = {
1483 1.0f, 0.0f, 0.0f, 0.0f, 1.0f, 0.0f, 0.0f, 0.0f, 1.0f,
1484 0.30901699437f, 0.9510565163f, 0.0f, -0.9510565163f, 0.30901699437f, 0.0f, 0.0f, 0.0f, 1.0f,
1485 -0.80901699437f, 0.58778525229f, 0.0f, -0.58778525229f, -0.80901699437f, 0.0f, 0.0f, 0.0f, 1.0f,
1486 -0.80901699437f, -0.58778525229f, 0.0f, 0.58778525229f, -0.80901699437f, 0.0f, 0.0f, 0.0f, 1.0f,
1487 0.30901699437f, -0.9510565163f, 0.0f, 0.9510565163f, 0.30901699437f, 0.0f, 0.0f, 0.0f, 1.0f,
1488 0.36180339887f, 0.58778525229f, -0.72360679775f, -0.26286555606f, 0.80901699437f, 0.52573111212f, 0.894427191f, 0.0f, 0.4472135955f,
1489 -0.13819660113f, 0.9510565163f, 0.27639320225f, -0.42532540418f, -0.30901699437f, 0.85065080835f, 0.894427191f, 0.0f, 0.4472135955f,
1490 -0.4472135955f, 0.0f, 0.894427191f, 0.0f, -1.0f, 0.0f, 0.894427191f, 0.0f, 0.4472135955f,
1491 -0.13819660113f, -0.9510565163f, 0.27639320225f, 0.42532540418f, -0.30901699437f, -0.85065080835f, 0.894427191f, 0.0f, 0.4472135955f,
1492 0.36180339887f, -0.58778525229f, -0.72360679775f, 0.26286555606f, 0.80901699437f, -0.52573111212f, 0.894427191f, 0.0f, 0.4472135955f,
1493 -0.4472135955f, 0.52573111212f, -0.72360679775f, -0.85065080835f, 0.0f, 0.52573111212f, 0.27639320225f, 0.85065080835f, 0.4472135955f,
1494 -0.9472135955f, 0.16245984812f, 0.27639320225f, 0.16245984812f, -0.5f, 0.85065080835f, 0.27639320225f, 0.85065080835f, 0.4472135955f,
1495 -0.13819660113f, -0.42532540418f, 0.894427191f, 0.9510565163f, -0.30901699437f, 0.0f, 0.27639320225f, 0.85065080835f, 0.4472135955f,
1496 0.86180339887f, -0.42532540418f, 0.27639320225f, 0.42532540418f, 0.30901699437f, -0.85065080835f, 0.27639320225f, 0.85065080835f, 0.4472135955f,
1497 0.67082039325f, 0.16245984812f, -0.72360679775f, -0.68819096024f, 0.5f, -0.52573111212f, 0.27639320225f, 0.85065080835f, 0.4472135955f,
1498 -0.63819660113f, -0.26286555606f, -0.72360679775f, -0.26286555606f, -0.80901699437f, 0.52573111212f, -0.72360679775f, 0.52573111212f, 0.4472135955f,
1499 -0.4472135955f, -0.85065080835f, 0.27639320225f, 0.52573111212f, 0.0f, 0.85065080835f, -0.72360679775f, 0.52573111212f, 0.4472135955f,
1500 0.36180339887f, -0.26286555606f, 0.894427191f, 0.58778525229f, 0.80901699437f, 0.0f, -0.72360679775f, 0.52573111212f, 0.4472135955f,
1501 0.67082039325f, 0.68819096024f, 0.27639320225f, -0.16245984812f, 0.5f, -0.85065080835f, -0.72360679775f, 0.52573111212f, 0.4472135955f,
1502 0.0527864045f, 0.68819096024f, -0.72360679775f, -0.68819096024f, -0.5f, -0.52573111212f, -0.72360679775f, 0.52573111212f, 0.4472135955f,
1503 0.0527864045f, -0.68819096024f, -0.72360679775f, 0.68819096024f, -0.5f, 0.52573111212f, -0.72360679775f, -0.52573111212f, 0.4472135955f,
1504 0.67082039325f, -0.68819096024f, 0.27639320225f, 0.16245984812f, 0.5f, 0.85065080835f, -0.72360679775f, -0.52573111212f, 0.4472135955f,
1505 0.36180339887f, 0.26286555606f, 0.894427191f, -0.58778525229f, 0.80901699437f, 0.0f, -0.72360679775f, -0.52573111212f, 0.4472135955f,
1506 -0.4472135955f, 0.85065080835f, 0.27639320225f, -0.52573111212f, 0.0f, -0.85065080835f, -0.72360679775f, -0.52573111212f, 0.4472135955f,
1507 -0.63819660113f, 0.26286555606f, -0.72360679775f, 0.26286555606f, -0.80901699437f, -0.52573111212f, -0.72360679775f, -0.52573111212f, 0.4472135955f,
1508 0.67082039325f, -0.16245984812f, -0.72360679775f, 0.68819096024f, 0.5f, 0.52573111212f, 0.27639320225f, -0.85065080835f, 0.4472135955f,
1509 0.86180339887f, 0.42532540418f, 0.27639320225f, -0.42532540418f, 0.30901699437f, 0.85065080835f, 0.27639320225f, -0.85065080835f, 0.4472135955f,
1510 -0.13819660113f, 0.42532540418f, 0.894427191f, -0.9510565163f, -0.30901699437f, 0.0f, 0.27639320225f, -0.85065080835f, 0.4472135955f,
1511 -0.9472135955f, -0.16245984812f, 0.27639320225f, -0.16245984812f, -0.5f, -0.85065080835f, 0.27639320225f, -0.85065080835f, 0.4472135955f,
1512 -0.4472135955f, -0.52573111212f, -0.72360679775f, 0.85065080835f, 0.0f, -0.52573111212f, 0.27639320225f, -0.85065080835f, 0.4472135955f,
1513 -0.36180339887f, -0.26286555606f, -0.894427191f, -0.58778525229f, 0.80901699437f, 0.0f, 0.72360679775f, 0.52573111212f, -0.4472135955f,
1514 -0.67082039325f, 0.68819096024f, -0.27639320225f, 0.16245984812f, 0.5f, 0.85065080835f, 0.72360679775f, 0.52573111212f, -0.4472135955f,
1515 -0.0527864045f, 0.68819096024f, 0.72360679775f, 0.68819096024f, -0.5f, 0.52573111212f, 0.72360679775f, 0.52573111212f, -0.4472135955f,
1516 0.63819660113f, -0.26286555606f, 0.72360679775f, 0.26286555606f, -0.80901699437f, -0.52573111212f, 0.72360679775f, 0.52573111212f, -0.4472135955f,
1517 0.4472135955f, -0.85065080835f, -0.27639320225f, -0.52573111212f, 0.0f, -0.85065080835f, 0.72360679775f, 0.52573111212f, -0.4472135955f,
1518 0.13819660113f, -0.42532540418f, -0.894427191f, -0.9510565163f, -0.30901699437f, 0.0f, -0.27639320225f, 0.85065080835f, -0.4472135955f,
1519 -0.86180339887f, -0.42532540418f, -0.27639320225f, -0.42532540418f, 0.30901699437f, 0.85065080835f, -0.27639320225f, 0.85065080835f, -0.4472135955f,
1520 -0.67082039325f, 0.16245984812f, 0.72360679775f, 0.68819096024f, 0.5f, 0.52573111212f, -0.27639320225f, 0.85065080835f, -0.4472135955f,
1521 0.4472135955f, 0.52573111212f, 0.72360679775f, 0.85065080835f, 0.0f, -0.52573111212f, -0.27639320225f, 0.85065080835f, -0.4472135955f,
1522 0.9472135955f, 0.16245984812f, -0.27639320225f, -0.16245984812f, -0.5f, -0.85065080835f, -0.27639320225f, 0.85065080835f, -0.4472135955f,
1523 0.4472135955f, 0.0f, -0.894427191f, 0.0f, -1.0f, 0.0f, -0.894427191f, 0.0f, -0.4472135955f,
1524 0.13819660113f, -0.9510565163f, -0.27639320225f, -0.42532540418f, -0.30901699437f, 0.85065080835f, -0.894427191f, 0.0f, -0.4472135955f,
1525 -0.36180339887f, -0.58778525229f, 0.72360679775f, -0.26286555606f, 0.80901699437f, 0.52573111212f, -0.894427191f, 0.0f, -0.4472135955f,
1526 -0.36180339887f, 0.58778525229f, 0.72360679775f, 0.26286555606f, 0.80901699437f, -0.52573111212f, -0.894427191f, 0.0f, -0.4472135955f,
1527 0.13819660113f, 0.9510565163f, -0.27639320225f, 0.42532540418f, -0.30901699437f, -0.85065080835f, -0.894427191f, 0.0f, -0.4472135955f,
1528 0.13819660113f, 0.42532540418f, -0.894427191f, 0.9510565163f, -0.30901699437f, 0.0f, -0.27639320225f, -0.85065080835f, -0.4472135955f,
1529 0.9472135955f, -0.16245984812f, -0.27639320225f, 0.16245984812f, -0.5f, 0.85065080835f, -0.27639320225f, -0.85065080835f, -0.4472135955f,
1530 0.4472135955f, -0.52573111212f, 0.72360679775f, -0.85065080835f, 0.0f, 0.52573111212f, -0.27639320225f, -0.85065080835f, -0.4472135955f,
1531 -0.67082039325f, -0.16245984812f, 0.72360679775f, -0.68819096024f, 0.5f, -0.52573111212f, -0.27639320225f, -0.85065080835f, -0.4472135955f,
1532 -0.86180339887f, 0.42532540418f, -0.27639320225f, 0.42532540418f, 0.30901699437f, -0.85065080835f, -0.27639320225f, -0.85065080835f, -0.4472135955f,
1533 -0.36180339887f, 0.26286555606f, -0.894427191f, 0.58778525229f, 0.80901699437f, 0.0f, 0.72360679775f, -0.52573111212f, -0.4472135955f,
1534 0.4472135955f, 0.85065080835f, -0.27639320225f, 0.52573111212f, 0.0f, 0.85065080835f, 0.72360679775f, -0.52573111212f, -0.4472135955f,
1535 0.63819660113f, 0.26286555606f, 0.72360679775f, -0.26286555606f, -0.80901699437f, 0.52573111212f, 0.72360679775f, -0.52573111212f, -0.4472135955f,
1536 -0.0527864045f, -0.68819096024f, 0.72360679775f, -0.68819096024f, -0.5f, -0.52573111212f, 0.72360679775f, -0.52573111212f, -0.4472135955f,
1537 -0.67082039325f, -0.68819096024f, -0.27639320225f, -0.16245984812f, 0.5f, -0.85065080835f, 0.72360679775f, -0.52573111212f, -0.4472135955f,
1538 -1.0f, 0.0f, 0.0f, 0.0f, 1.0f, 0.0f, 0.0f, 0.0f, -1.0f,
1539 -0.30901699437f, 0.9510565163f, 0.0f, 0.9510565163f, 0.30901699437f, 0.0f, 0.0f, 0.0f, -1.0f,
1540 0.80901699437f, 0.58778525229f, 0.0f, 0.58778525229f, -0.80901699437f, 0.0f, 0.0f, 0.0f, -1.0f,
1541 0.80901699437f, -0.58778525229f, 0.0f, -0.58778525229f, -0.80901699437f, 0.0f, 0.0f, 0.0f, -1.0f,
1542 -0.30901699437f, -0.9510565163f, 0.0f, -0.9510565163f, 0.30901699437f, 0.0f, 0.0f, 0.0f, -1.0f
1543 };
1544
1545 for (int k=0;k<nsym;k++) {
1546 Transform t;
1547 for (int i=0; i<3; i++) {
1548 for (int j=0; j<3; j++) {
1549 t.matrix[i][j] = TET[9*k + i*3 +j];
1550 }
1551 }
1552 //vector<float> z = t.get_matrix();
1553 //for (int i=0; i<12; i++) cout<<z[i]<<endl;
1554 ret.push_back( (*this) * t );
1555 }
1556 } else {
1558 nsym = sym->get_nsym();
1559 for (int k=0;k<nsym;k++) {
1560 Transform t;
1561 t = sym->get_sym(k);
1562 ret.push_back( (*this) * t );
1563 }
1564 delete sym;
1565 }
1566 return ret;
1567}

References EMAN::Factory< T >::get(), EMAN::Symmetry3D::get_nsym(), EMAN::Symmetry3D::get_sym(), matrix, and EMAN::Util::str_to_lower().

Referenced by get_sym_sparx(), EMAN::nn4Reconstructor::insert_padfft_slice(), EMAN::nn4_rectReconstructor::insert_padfft_slice(), EMAN::nnSSNR_Reconstructor::insert_padfft_slice(), EMAN::nnSSNR_ctfReconstructor::insert_padfft_slice(), EMAN::nn4_ctf_rectReconstructor::insert_padfft_slice(), EMAN::nn4_ctfReconstructor::insert_padfft_slice(), EMAN::nn4_ctfwReconstructor::insert_padfft_slice_weighted(), and EMAN::nn4_ctfwsReconstructor::insert_padfft_slice_weighted().

◆ get_sym_sparx()

Transform Transform::get_sym_sparx ( const string &  sym,
int  n 
) const

Definition at line 1398 of file transform.cpp.

1399{
1400 Transform ret;
1401 vector<Transform> tut;
1402 tut = ret.get_sym_proj(sym_name);
1403 ret = (*this) * tut[n];
1404 return ret;
1405}
vector< Transform > get_sym_proj(const string &sym) const
Definition: transform.cpp:1407

References get_sym_proj().

◆ get_trans()

Vec3f Transform::get_trans ( ) const

◆ get_trans_2d()

Vec2f Transform::get_trans_2d ( ) const

Get the degenerant 2D post trans as a vec2f.

Returns
the 2D translation vector

Definition at line 1088 of file transform.cpp.

1089{
1090 bool x_mirror = get_mirror();
1091 Vec2f v;
1092 if (x_mirror) v[0] = -matrix[0][3];
1093 else v[0] = matrix[0][3];
1094 v[1] = matrix[1][3];
1095 return v;
1096}
The Vec2 is precisely the same as Vec3 except it works exclusively in 2D Note there are convenient ty...
Definition: vec3.h:710

References get_mirror(), and matrix.

Referenced by get_pre_trans_2d(), and EMAN::padfft_slice().

◆ get_vflip_transform()

Transform Transform::get_vflip_transform ( ) const

How do I get the transform that will yield the vertically flipped projection?

Returns
the transform that will yield the vertically flipped projection

Definition at line 811 of file transform.cpp.

811 {
812
813 Dict rot = get_rotation("eman");
814 rot["alt"] = 180.0f + static_cast<float>(rot["alt"]);
815 rot["phi"] = - static_cast<float>(rot["phi"]);
816
817 // This is the same as new_alt= 180-alt, new_phi=180-phi, new_az=180+az
818
819 Transform ret(*this);
820 ret.set_rotation(rot);
821
822 Vec3f trans = get_trans();
823 trans[1] = -trans[1];
824 ret.set_trans(trans);
825
826 return ret;
827}

References get_rotation(), get_trans(), set_rotation(), and set_trans().

◆ icos_5_to_2()

Transform Transform::icos_5_to_2 ( )
static

Get the transform that moves any icosahedron generated by eman2 so that it matches the 2-2-2 (MRC, FREALIGN) convention.

Returns
a Transforms, alt=58.282523, az=270

Doctor Steve says Phil's answer put the 2-fold in the wrong place based on the standard Virus convention (empirically). It was also 2 to 5 not 5 to 2. It is possible to rotate to a 2-fold directly from a 5-fold, though there are 2 possible orientations for the 2-2-2 convention, and this finds only one of them :^(

Doctor Phil says: alt = (acos(cos(pi/5)/sqrt(3)/sin(pi/5)) + acos(2*cos(pi/5)/ sqrt(3) ) )*180/pi This is the angle between a 5 and a 3 plus the angle between a 3 and a 2

Definition at line 63 of file transform.cpp.

63 {
64 Transform t;
65 Dict d;
66 d["type"] = "eman";
67 d["phi"] = 0;
68 d["az"] = 90.0f;
69 d["alt"] = 31.717474; // 5 fold to the nearest 2-fold
81 t.set_rotation(d);
82 return t;
83}
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

References set_rotation().

◆ init_permissable_keys()

void Transform::init_permissable_keys ( )
private

Called internally to initialize permissable_2d_not_rot, permissable_3d_not_rot, and permissable_rot_keys static members.

Definition at line 285 of file transform.cpp.

286{
287
288 permissable_2d_not_rot.push_back("tx");
289 permissable_2d_not_rot.push_back("ty");
290 permissable_2d_not_rot.push_back("scale");
291 permissable_2d_not_rot.push_back("mirror");
292 permissable_2d_not_rot.push_back("type");
293
294 permissable_3d_not_rot.push_back("tx");
295 permissable_3d_not_rot.push_back("ty");
296 permissable_3d_not_rot.push_back("tz");
297 permissable_3d_not_rot.push_back("scale");
298 permissable_3d_not_rot.push_back("mirror");
299 permissable_3d_not_rot.push_back("type");
300
301 vector<string> tmp;
302 tmp.push_back("alpha");
303 permissable_rot_keys["2d"] = tmp;
304
305 tmp.clear();
306 tmp.push_back("alt");
307 tmp.push_back("az");
308 tmp.push_back("phi");
309 permissable_rot_keys["eman"] = tmp;
310
311 tmp.clear();
312 tmp.push_back("psi");
313 tmp.push_back("theta");
314 tmp.push_back("phi");
315 permissable_rot_keys["spider"] = tmp;
316
317 tmp.clear();
318 tmp.push_back("alpha");
319 tmp.push_back("beta");
320 tmp.push_back("gamma");
321 permissable_rot_keys["imagic"] = tmp;
322
323 tmp.clear();
324 tmp.push_back("ztilt");
325 tmp.push_back("xtilt");
326 tmp.push_back("ytilt");
327 permissable_rot_keys["xyz"] = tmp;
328
329 tmp.clear();
330 tmp.push_back("phi");
331 tmp.push_back("theta");
332 tmp.push_back("omega");
333 permissable_rot_keys["mrc"] = tmp;
334
335 tmp.clear();
336 tmp.push_back("e0");
337 tmp.push_back("e1");
338 tmp.push_back("e2");
339 tmp.push_back("e3");
340 permissable_rot_keys["quaternion"] = tmp;
341
342 tmp.clear();
343 tmp.push_back("n1");
344 tmp.push_back("n2");
345 tmp.push_back("n3");
346 tmp.push_back("omega");
347 permissable_rot_keys["spin"] = tmp;
348
349 tmp.clear();
350 tmp.push_back("n1");
351 tmp.push_back("n2");
352 tmp.push_back("n3");
353 tmp.push_back("q");
354 permissable_rot_keys["sgirot"] = tmp;
355
356 tmp.clear();
357 tmp.push_back("m11");
358 tmp.push_back("m12");
359 tmp.push_back("m13");
360 tmp.push_back("m21");
361 tmp.push_back("m22");
362 tmp.push_back("m23");
363 tmp.push_back("m31");
364 tmp.push_back("m32");
365 tmp.push_back("m33");
366 permissable_rot_keys["matrix"] = tmp;
367}

References permissable_2d_not_rot, permissable_3d_not_rot, and permissable_rot_keys.

Referenced by detect_problem_keys().

◆ inverse()

Transform Transform::inverse ( ) const

◆ invert()

void Transform::invert ( )

Get the inverse of this transformation matrix.

Returns
the inverse of this transformation matrix

Definition at line 1293 of file transform.cpp.

1293 {
1294
1295 double m00 = matrix[0][0]; double m01=matrix[0][1]; double m02=matrix[0][2];
1296 double m10 = matrix[1][0]; double m11=matrix[1][1]; double m12=matrix[1][2];
1297 double m20 = matrix[2][0]; double m21=matrix[2][1]; double m22=matrix[2][2];
1298 double v0 = matrix[0][3]; double v1 =matrix[1][3]; double v2 =matrix[2][3];
1299
1300 double cof00 = m11*m22-m12*m21;
1301 double cof11 = m22*m00-m20*m02;
1302 double cof22 = m00*m11-m01*m10;
1303 double cof01 = m10*m22-m20*m12;
1304 double cof02 = m10*m21-m20*m11;
1305 double cof12 = m00*m21-m01*m20;
1306 double cof10 = m01*m22-m02*m21;
1307 double cof20 = m01*m12-m02*m11;
1308 double cof21 = m00*m12-m10*m02;
1309
1310 double det = m00* cof00 + m02* cof02 -m01*cof01;
1311
1312 matrix[0][0] = (float)(cof00/det);
1313 matrix[0][1] = - (float)(cof10/det);
1314 matrix[0][2] = (float)(cof20/det);
1315 matrix[1][0] = - (float)(cof01/det);
1316 matrix[1][1] = (float)(cof11/det);
1317 matrix[1][2] = - (float)(cof21/det);
1318 matrix[2][0] = (float)(cof02/det);
1319 matrix[2][1] = - (float)(cof12/det);
1320 matrix[2][2] = (float)(cof22/det);
1321
1322 matrix[0][3] = (float)((- cof00*v0 + cof10*v1 - cof20*v2)/det);
1323 matrix[1][3] = (float)(( cof01*v0 - cof11*v1 + cof21*v2)/det);
1324 matrix[2][3] = (float)((- cof02*v0 + cof12*v1 - cof22*v2)/det);
1325}
#define v0(i)
Definition: analyzer.cpp:698

References matrix, and v0.

Referenced by get_pre_trans(), get_pre_trans_2d(), EMAN::Symmetry3D::in_which_asym_unit(), inverse(), EMAN::GaussFFTProjector::project3d(), EMAN::Symmetry3D::reduce(), rotate_origin_newBasis(), set_params_inverse(), set_pre_trans(), translate_newBasis(), EMAN::RT3DSphereAligner::xform_align_nbest(), EMAN::RT3DTreeAligner::xform_align_nbest(), and EMAN::RT3DLocalTreeAligner::xform_align_nbest().

◆ is_identity()

bool Transform::is_identity ( ) const

Returns whethers or this matrix is the identity.

Definition at line 222 of file transform.cpp.

222 {
223 for(int i=0; i<3; ++i) {
224 for(int j=0; j<4; ++j) {
225 float c = matrix[i][j];
227 if(i==j) {
228 if (c != 1.0) return false;
229 }
230 else {
231 if (c != 0.0) return false;
232 }
233 }
234 }
235 return true;
236}

References EMAN::Util::apply_precision(), ERR_LIMIT, and matrix.

Referenced by EMAN::TestTomoImage::insert_rectangle(), EMAN::FourierReconstructor::preprocess_slice(), and EMAN::FourierIterReconstructor::preprocess_slice().

◆ is_rot_identity()

bool Transform::is_rot_identity ( ) const

Returns whethers or this matrix rotation is the identity.

Definition at line 238 of file transform.cpp.

238 {
239 for(int i=0; i<3; ++i) {
240 for(int j=0; j<3; ++j) {
241 float c = matrix[i][j];
243 if(i==j) {
244 if (c != 1.0) return false;
245 }
246 else {
247 if (c != 0.0) return false;
248 }
249 }
250 }
251 return true;
252}

References EMAN::Util::apply_precision(), ERR_LIMIT, and matrix.

◆ negate()

Transform Transform::negate ( ) const

Negates the Transform - a useful way, for example, for getting an orientation on the opposite side of the sphere.

Returns
a transform that has been negated

Definition at line 781 of file transform.cpp.

782{
783 Transform t(*this);
784 for(unsigned int i = 0; i < 3; ++i) {
785 for(unsigned int j = 0; j < 4; ++j) {
786 t.set(i,j,t[i][j]*-1);
787 }
788 }
789 return t;
790}

References set().

◆ operator!=()

bool Transform::operator!= ( const Transform rhs) const

Unequality comparision operator.

Parameters
rhsthe Transform object compared to
Returns
a boolean indicate if the pass in transform is not equal this object

Definition at line 128 of file transform.cpp.

128 {
129 return !(operator==(rhs));
130}
bool operator==(const Transform &rhs) const
Equality comparision operator.
Definition: transform.cpp:119

References operator==().

◆ operator=()

Transform & Transform::operator= ( const Transform that)

Assignment operator.

Parameters
thatthat which this will become

Definition at line 110 of file transform.cpp.

110 {
111 if (this != &that ) {
112// for (int i=0; i<3; i++)
113// for (int j=0; j<4; j++) matrix[i][j]=that.matrix[i][j];
114 memcpy(matrix,that.matrix,12*sizeof(float));
115 }
116 return *this;
117}

References matrix.

◆ operator==()

bool Transform::operator== ( const Transform rhs) const

Equality comparision operator.

Parameters
rhsthe Transform object compared to
Returns
a boolean indicate if the pass in transform is equal this object

Definition at line 119 of file transform.cpp.

119 {
120 if (memcmp(this->matrix, rhs.matrix, 3*4*sizeof(float)) == 0) {
121 return true;
122 }
123 else {
124 return false;
125 }
126}

References matrix.

Referenced by operator!=().

◆ operator[]() [1/2]

float * EMAN::Transform::operator[] ( int  i)
inline

Operator[] convenience so Transform3D[2][2] etc terminology can be used.

Definition at line 405 of file transform.h.

405{ return matrix[i]; }

References matrix.

◆ operator[]() [2/2]

const float * EMAN::Transform::operator[] ( int  i) const
inline

Operator[] convenience so Transform3D[2][2] etc terminology can be used.

Definition at line 410 of file transform.h.

410{ return matrix[i]; }

References matrix.

◆ orthogonalize()

void Transform::orthogonalize ( )

Reorthogonalize the rotation part of the matrix in place.

Does this by performing the SVD decomposition of the rotation matrix R such that R = USV^T - since the eigenvalues of a rotation matrix are all 1 we enforce that S should be the identity and produce a corrected matrix R' = UV^T

Definition at line 1179 of file transform.cpp.

1180{
1181 float scale;
1182 bool x_mirror;
1183 get_scale_and_mirror(scale,x_mirror);
1184 if (scale == 0) throw UnexpectedBehaviorException("The determinant of the Transform is 0. This is unexpected.");
1185 double inv_scale = 1.0/static_cast<double>(scale);
1186 double mirror_scale = (x_mirror == true ? -1.0:1.0);
1187
1188 gsl_matrix * R = gsl_matrix_calloc(3,3);
1189 for ( unsigned int i = 0; i < 3; ++i )
1190 {
1191 for ( unsigned int j = 0; j < 3; ++j )
1192 {
1193 if (i == 0 && mirror_scale != 1.0 ) {
1194 gsl_matrix_set( R, i, j, static_cast<double>(matrix[i][j])*mirror_scale*inv_scale );
1195 }
1196 else {
1197 gsl_matrix_set( R, i, j, static_cast<double>(matrix[i][j])*inv_scale );
1198 }
1199 }
1200 }
1201
1202 gsl_matrix * V = gsl_matrix_calloc(3,3);
1203 gsl_vector * S = gsl_vector_calloc(3);
1204 gsl_vector * work = gsl_vector_calloc(3);
1205 gsl_linalg_SV_decomp (R, V, S, work); // Now R is U of the SVD R = USV^T
1206
1207 gsl_matrix * Soln = gsl_matrix_calloc(3,3);
1208 gsl_blas_dgemm (CblasNoTrans, CblasTrans, 1.0, R, V, 0.0, Soln);
1209
1210 for ( unsigned int i = 0; i < 3; ++i )
1211 {
1212 for ( unsigned int j = 0; j < 3; ++j )
1213 {
1214 matrix[i][j] = static_cast<float>( gsl_matrix_get(Soln,i,j) );
1215 }
1216 }
1217
1218 // Apply scale if it existed previously
1219 if (scale != 1.0f) {
1220 for(int i=0; i<3; ++i) {
1221 for(int j=0; j<3; ++j) {
1222 matrix[i][j] *= scale;
1223 }
1224 }
1225 }
1226
1227 // Apply post x mirroring if it was applied previouslys
1228 if ( x_mirror ) {
1229 for(int j=0; j<3; ++j) {
1230 matrix[0][j] *= -1.0f;
1231 }
1232 }
1233
1234 gsl_matrix_free(V); gsl_matrix_free(R); gsl_matrix_free(Soln);
1235 gsl_vector_free(S); gsl_vector_free(work);
1236}
#define V(i, j)
Definition: analyzer.cpp:697

References get_scale_and_mirror(), matrix, scale(), UnexpectedBehaviorException, and V.

◆ printme()

void EMAN::Transform::printme ( ) const
inline

Print the contents of the internal matrix verbatim to standard out.

Definition at line 348 of file transform.h.

348 {
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 }

References matrix.

◆ rotate()

void Transform::rotate ( const Transform by)

Increment the rotation by multipling the rotation bit of the argument transfrom by the current transfrom.

Parameters
rotationmultiplican, a tranform, R'', by which to multiply the current one, R' == R''R'

Definition at line 763 of file transform.cpp.

764{
765 vector<float> multmatrix = by.get_matrix();
766 // First Multiply and put the result in a temp matrix
767 Transform result;
768 for (int i=0; i<3; i++) {
769 for (int j=0; j<4; j++) {
770 result[i][j] = multmatrix[i*4]*matrix[0][j] + multmatrix[i*4+1]*matrix[1][j] + multmatrix[i*4+2]*matrix[2][j];
771 }
772 }
773 //Then put the result from the tmep matrix in the original one
774 for (int i=0; i<3; i++) {
775 for (int j=0; j<4; j++) {
776 matrix[i][j] = result[i][j];
777 }
778 }
779}
vector< float > get_matrix() const
Get the transformation matrix using a vector.
Definition: transform.cpp:181

References get_matrix(), and matrix.

◆ rotate_origin()

void Transform::rotate_origin ( const Transform by)

Increment the rotation by multipling the rotation bit of the argument transfrom by the rotation part of the current transfrom.

Parameters
bymultiplican, a tranform, R'', by which to multiply the current one, R' == R''R'

Definition at line 720 of file transform.cpp.

721{
722 vector<float> multmatrix = by.get_matrix();
723 // First Multiply and put the result in a temp matrix
724 Transform result;
725 for (int i=0; i<3; i++) {
726 for (int j=0; j<3; j++) {
727 result[i][j] = multmatrix[i*4]*matrix[0][j] + multmatrix[i*4+1]*matrix[1][j] + multmatrix[i*4+2]*matrix[2][j];
728 }
729 }
730 //Then put the result from the tmep matrix in the original one
731 for (int i=0; i<3; i++) {
732 for (int j=0; j<3; j++) {
733 matrix[i][j] = result[i][j];
734 }
735 }
736}

References get_matrix(), and matrix.

Referenced by rotate_origin_newBasis().

◆ rotate_origin_newBasis()

void Transform::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 of the current transfrom This version rotates in the standard coordinate system, even after it have been modified by tcs.

The effect is to undo the distortion casued by dcs. Useful in the scenegraph

Parameters
tcs,thestansfrom that moves us to a non standard coordinate system
bymultiplican, a tranform, R'', by which to multiply the current one, R' == R''R'

Definition at line 738 of file transform.cpp.

739{
740 //Get the rotational inverse
741 Transform tcsinv = Transform(tcs);
742 tcsinv.set_trans(0.0, 0.0, 0.0);
743 tcsinv.set_scale(1.0);
744 tcsinv.invert();
745
746 //Get the current rotation
747 Transform temp = Transform();
748 temp.set_trans(n1, n2, n3);
749 Transform cc = tcsinv*temp;
750 Vec3f cctrans = cc.get_trans();
751
752 //set the right rotation
753 Dict spinrot = Dict();
754 spinrot["type"] = "spin";
755 spinrot["omega"] = omega;
756 spinrot["n1"] = cctrans[0];
757 spinrot["n2"] = cctrans[1];
758 spinrot["n3"] = cctrans[2];
759 Transform rightrot = Transform(spinrot);
760 rotate_origin(rightrot);
761}
Transform()
Default constructor Internal matrix is the identity.
Definition: transform.cpp:100
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_trans(const float &x, const float &y, const float &z=0)
Set the post translation component.
Definition: transform.cpp:1036
void invert()
Get the inverse of this transformation matrix.
Definition: transform.cpp:1293

References get_trans(), invert(), rotate_origin(), set_scale(), set_trans(), and Transform().

◆ scale()

void Transform::scale ( const float &  scale)

Increment the scale.

Parameters
scalethe amount increment by

Definition at line 1159 of file transform.cpp.

1160{
1161 float determinant = get_determinant();
1162 if (determinant < 0) determinant *= -1.0f;
1163 float newscale = std::pow(determinant,1.0f/3.0f) + scale;
1164 if(newscale > 0.0001) set_scale(newscale); // If scale ~ 0 things blowup, so we need a little fudge factor
1165}

References get_determinant(), scale(), and set_scale().

Referenced by get_params(), get_params_inverse(), get_rotation(), get_scale(), get_scale_and_mirror(), orthogonalize(), scale(), set_params(), set_params_inverse(), set_pre_trans(), and set_rotation().

◆ set()

void EMAN::Transform::set ( int  r,
int  c,
float  value 
)
inline

Set the value stored in the internal transformation matrix at at coordinate (r,c) to value.

Definition at line 400 of file transform.h.

400{ matrix[r][c] = value; }

References matrix.

Referenced by EMAN::PointArray::align_2d(), and negate().

◆ set_matrix()

void Transform::set_matrix ( const vector< float > &  v)

Set the transformation matrix using a vector.

Must be of length 12.

Parameters
vthe transformation matrix stored as a vector - 3 rows of 4.

Definition at line 147 of file transform.cpp.

148{
149 if (v.size() != 12 ) throw InvalidParameterException("The construction array must be of size 12");
150
151 for(int i=0; i<3; ++i) {
152 for(int j=0; j<4; ++j) {
153 matrix[i][j] = v[i*4+j];
154 }
155 }
156}

References InvalidParameterException, and matrix.

Referenced by EMAN::EMObject::operator Transform *(), and Transform().

◆ set_mirror()

void Transform::set_mirror ( const bool  x_mirror)

Set whether or not x_mirroring is occurring.

Parameters
x_mirrorwhether x_mirroring should be applied

Definition at line 1238 of file transform.cpp.

1238 {
1239
1240 bool old_x_mirror = get_mirror();
1241 if (old_x_mirror == x_mirror) return; // The user is setting the same value
1242 else {
1243 // Toggle the mirroring operation
1244 for (int j = 0; j < 4; ++j ) {
1245 matrix[0][j] *= -1;
1246 }
1247 }
1248}

References get_mirror(), and matrix.

Referenced by EMAN::FourierReconstructor::determine_slice_agreement(), EMAN::WienerFourierReconstructor::determine_slice_agreement(), get_rotation_transform(), EMAN::FourierReconstructor::insert_slice(), EMAN::WienerFourierReconstructor::insert_slice(), EMAN::FourierIterReconstructor::insert_slice(), EMAN::FourierReconstructor::projection(), refalifn(), set_params(), and set_params_inverse().

◆ set_params()

void Transform::set_params ( const Dict d)

Set the parameters of the entire transform.

keys acted upon are "type" - if this exists then the correct euler angles need to be included - also "tx","ty","tz", "scale", and "mirror"

Parameters
dthe dictionary containing the parameters

Definition at line 254 of file transform.cpp.

254 {
256
257 if (d.has_key_ci("type") ) set_rotation(d);
258
259 if (d.has_key_ci("scale")) {
260 float scale = static_cast<float>(d.get_ci("scale"));
262 }
263
264 float dx=0,dy=0,dz=0;
265
266 if (d.has_key_ci("tx")) dx = static_cast<float>(d.get_ci("tx"));
267 if (d.has_key_ci("ty")) dy = static_cast<float>(d.get_ci("ty"));
268 if (d.has_key_ci("tz")) dz = static_cast<float>(d.get_ci("tz"));
269
270 if ( dx != 0.0 || dy != 0.0 || dz != 0.0 ) {
271 set_trans(dx,dy,dz);
272 }
273
274 if (d.has_key_ci("mirror")) {
275 EMObject e = d.get_ci("mirror");
276 if ( (e.get_type() != EMObject::BOOL ) && (e.get_type() != EMObject::INT ) && (e.get_type() != EMObject::UNSIGNEDINT ) )
277 throw InvalidParameterException("Error, mirror must be a bool or an int");
278
279 bool mirror = static_cast<bool>(e);
280 set_mirror(mirror);
281 }
282}
EMObject get_ci(const string &key) const
Get the EMObject corresponding to the particular key using case insensitivity.
Definition: emobject.cpp:1128
EMObject is a wrapper class for types including int, float, double, etc as defined in ObjectType.
Definition: emobject.h:123
ObjectType get_type() const
Get the ObjectType This is an enumerated type first declared in the class EMObjectTypes.
Definition: emobject.h:224
void set_mirror(const bool x_mirror)
Set whether or not x_mirroring is occurring.
Definition: transform.cpp:1238
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

References EMAN::EMObject::BOOL, detect_problem_keys(), EMAN::Dict::get_ci(), EMAN::EMObject::get_type(), EMAN::Dict::has_key_ci(), EMAN::EMObject::INT, InvalidParameterException, scale(), set_mirror(), set_rotation(), set_scale(), set_trans(), and EMAN::EMObject::UNSIGNEDINT.

Referenced by EMAN::RT3DLocalTreeAligner::testort(), EMAN::RT2Dto3DTreeAligner::testort(), EMAN::RT3DTreeAligner::testort(), Transform(), EMAN::RT2DTreeAligner::xform_align_nbest(), EMAN::RT2Dto3DTreeAligner::xform_align_nbest(), EMAN::RT3DTreeAligner::xform_align_nbest(), and EMAN::RT3DLocalTreeAligner::xform_align_nbest().

◆ set_params_inverse()

void Transform::set_params_inverse ( const Dict d)

Set the parameters of the entire transform as though they there in the inverse format.

in other words, calling set_params_inverse(get_params_inverse()) should essentially leave the object unchanged.

Parameters
dthe dictionary containing the inverse parameters

Definition at line 425 of file transform.cpp.

425 {
427
428 if (d.has_key_ci("type") ) set_rotation(d);
429
430 float dx=0,dy=0,dz=0;
431 if (d.has_key_ci("tx")) dx = static_cast<float>(d.get_ci("tx"));
432 if (d.has_key_ci("ty")) dy = static_cast<float>(d.get_ci("ty"));
433 if (d.has_key_ci("tz")) dz = static_cast<float>(d.get_ci("tz"));
434
435 if ( (dx != 0.0 || dy != 0.0 || dz != 0.0) && d.has_key_ci("type") ) {
436 Transform pre_trans;
437 pre_trans.set_trans(dx,dy,dz);
438
439 Transform tmp;
440 tmp.set_rotation(d);
441
442 if (d.has_key_ci("scale")) {
443 float scale = static_cast<float>(d.get_ci("scale"));
444 tmp.set_scale(scale);
445 }
446
447 Transform solution_trans = tmp*pre_trans;
448
449 if (d.has_key_ci("scale")) {
450 Transform tmp;
451 float scale = static_cast<float>(d.get_ci("scale"));
452 tmp.set_scale(scale);
453 solution_trans = solution_trans*tmp;
454 }
455
456 tmp = Transform();
457 tmp.set_rotation(d);
458 solution_trans = solution_trans*tmp;
459 set_trans(solution_trans.get_trans());
460 }
461
462 if (d.has_key_ci("scale")) {
463 float scale = static_cast<float>(d.get_ci("scale"));
465 }
466
467 if (d.has_key_ci("mirror")) {
468 EMObject e = d.get_ci("mirror");
469 if ( (e.get_type() != EMObject::BOOL ) && (e.get_type() != EMObject::INT ) && (e.get_type() != EMObject::UNSIGNEDINT ) )
470 throw InvalidParameterException("Error, mirror must be a bool or an int");
471
472 bool mirror = static_cast<bool>(e);
473 set_mirror(mirror);
474 }
475 invert();
476}

References EMAN::EMObject::BOOL, detect_problem_keys(), EMAN::Dict::get_ci(), get_trans(), EMAN::EMObject::get_type(), EMAN::Dict::has_key_ci(), EMAN::EMObject::INT, InvalidParameterException, invert(), scale(), set_mirror(), set_rotation(), set_scale(), set_trans(), Transform(), and EMAN::EMObject::UNSIGNEDINT.

◆ set_pre_trans()

template<typename type >
void EMAN::Transform::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 translation.

Internally the correct form of MTSR is computed.

Parameters
vthe vector (Vec3f or Vec2f) that is the pre trans

Definition at line 556 of file transform.h.

556 {
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 }

References get_rotation(), get_scale(), get_trans(), invert(), scale(), set_rotation(), set_scale(), and set_trans().

Referenced by EMAN::PointArray::align_2d(), and EMAN::EMData::rotate_translate().

◆ set_rotation() [1/2]

void Transform::set_rotation ( const Dict rotation)

Set a rotation using a specific Euler type and the dictionary interface Works for all Euler types.

Parameters
rotationa dictionary containing all key-entry pair required of the associated Euler type

Definition at line 519 of file transform.cpp.

520{
521 detect_problem_keys(rotation);
522 string euler_type;
523
524 if (!rotation.has_key_ci("type") ){
525 throw InvalidParameterException("argument dictionary does not contain the type key");
526 }
527
528 euler_type = static_cast<string>(rotation.get_ci("type"));// Warning, will throw
529
530
531 double e0=0; double e1=0; double e2=0; double e3=0;
532 double omega=0;
533 double az = 0;
534 double alt = 0;
535 double phi = 0;
536 double cxtilt = 0;
537 double sxtilt = 0;
538 double cytilt = 0;
539 double sytilt = 0;
540 double cztilt = 0;
541 double sztilt = 0;
542 bool is_quaternion = 0;
543 bool is_matrix = 0;
544 bool is_xyz = 0;
545
546 bool x_mirror;
547 float scale;
548 // Get these before anything changes so we can apply them again after the rotation is set
549 get_scale_and_mirror(scale,x_mirror);
550 if (scale == 0) throw UnexpectedBehaviorException("The determinant of the Transform is 0. This is unexpected.");
551
552 string type = Util::str_to_lower(euler_type);
553 if (type == "2d") {
555 az = 0;
556 alt = 0;
557 phi = (double)rotation["alpha"] ;
558 } else if ( type == "eman" ) {
559// validate_and_set_type(THREED);
560 az = (double)rotation["az"] ;
561 alt = (double)rotation["alt"] ;
562 phi = (double)rotation["phi"] ;
563 } else if ( type == "imagic" ) {
564// validate_and_set_type(THREED);
565 az = (double)rotation["alpha"] ;
566 alt = (double)rotation["beta"] ;
567 phi = (double)rotation["gamma"] ;
568 } else if ( type == "spider" ) {
569// validate_and_set_type(THREED);
570 az = (double)rotation["phi"] + 90.0;
571 alt = (double)rotation["theta"] ;
572 phi = (double)rotation["psi"] - 90.0;
573 } else if ( type == "xyz" ) {
574// validate_and_set_type(THREED);
575 is_xyz = 1;
576 cxtilt = cos(EMConsts::deg2rad*(double)rotation["xtilt"]);
577 sxtilt = sin(EMConsts::deg2rad*(double)rotation["xtilt"]);
578 cytilt = cos(EMConsts::deg2rad*(double)rotation["ytilt"]);
579 sytilt = sin(EMConsts::deg2rad*(double)rotation["ytilt"]);
580 cztilt = cos(EMConsts::deg2rad*(double)rotation["ztilt"]);
581 sztilt = sin(EMConsts::deg2rad*(double)rotation["ztilt"]);
582 } else if ( type == "mrc" ) {
583// validate_and_set_type(THREED);
584 az = (double)rotation["phi"] + 90.0f ;
585 alt = (double)rotation["theta"] ;
586 phi = (double)rotation["omega"] - 90.0f ;
587 } else if ( type == "quaternion" ) {
588// validate_and_set_type(THREED);
589 is_quaternion = 1;
590 e0 = (double)rotation["e0"];
591 e1 = (double)rotation["e1"];
592 e2 = (double)rotation["e2"];
593 e3 = (double)rotation["e3"];
594 } else if ( type == "spin" ) {
595// validate_and_set_type(THREED);
596 is_quaternion = 1;
597 omega = (double)rotation["omega"];
598 double norm=Util::hypot3((double)rotation["n1"],(double)rotation["n2"],(double)rotation["n3"]);
599 if (norm==0.0) {
600 e0=1.0;
601 e1=e2=e3=0.0;
602 } else {
603 e0 = cos(omega*EMConsts::deg2rad/2.0);
604 e1 = sin(omega*EMConsts::deg2rad/2.0) * (double)rotation["n1"]/norm;
605 e2 = sin(omega*EMConsts::deg2rad/2.0) * (double)rotation["n2"]/norm;
606 e3 = sin(omega*EMConsts::deg2rad/2.0) * (double)rotation["n3"]/norm;
607 }
608 } else if ( type == "sgirot" ) {
609// validate_and_set_type(THREED);
610 is_quaternion = 1;
611 omega = (double)rotation["q"] ;
612 e0 = cos(omega*EMConsts::deg2rad/2.0);
613 e1 = sin(omega*EMConsts::deg2rad/2.0) * (double)rotation["n1"];
614 e2 = sin(omega*EMConsts::deg2rad/2.0) * (double)rotation["n2"];
615 e3 = sin(omega*EMConsts::deg2rad/2.0) * (double)rotation["n3"];
616 } else if ( type == "matrix" ) {
617 is_matrix = 1;
618 matrix[0][0] = (float)rotation["m11"];
619 matrix[0][1] = (float)rotation["m12"];
620 matrix[0][2] = (float)rotation["m13"];
621 matrix[1][0] = (float)rotation["m21"];
622 matrix[1][1] = (float)rotation["m22"];
623 matrix[1][2] = (float)rotation["m23"];
624 matrix[2][0] = (float)rotation["m31"];
625 matrix[2][1] = (float)rotation["m32"];
626 matrix[2][2] = (float)rotation["m33"];
627 } else {
628// transform_type = UNKNOWN;
629 throw InvalidStringException(euler_type, "unknown Euler Type");
630 }
631
632 double azp = az*EMConsts::deg2rad;
633 double altp = alt*EMConsts::deg2rad;
634 double phip = phi*EMConsts::deg2rad;
635
636 if (!is_quaternion && !is_matrix && !is_xyz) {
637 matrix[0][0] = (float)(cos(phip)*cos(azp) - cos(altp)*sin(azp)*sin(phip));
638 matrix[0][1] = (float)(cos(phip)*sin(azp) + cos(altp)*cos(azp)*sin(phip));
639 matrix[0][2] = (float)(sin(altp)*sin(phip));
640 matrix[1][0] = (float)(-sin(phip)*cos(azp) - cos(altp)*sin(azp)*cos(phip));
641 matrix[1][1] = (float)(-sin(phip)*sin(azp) + cos(altp)*cos(azp)*cos(phip));
642 matrix[1][2] = (float)(sin(altp)*cos(phip));
643 matrix[2][0] = (float)(sin(altp)*sin(azp));
644 matrix[2][1] = (float)(-sin(altp)*cos(azp));
645 matrix[2][2] = (float)cos(altp);
646 }
647 if (is_quaternion){
648 matrix[0][0] = (float)(e0 * e0 + e1 * e1 - e2 * e2 - e3 * e3);
649 matrix[0][1] = (float)(2.0f * (e1 * e2 + e0 * e3));
650 matrix[0][2] = (float)(2.0f * (e1 * e3 - e0 * e2));
651 matrix[1][0] = (float)(2.0f * (e2 * e1 - e0 * e3));
652 matrix[1][1] = (float)(e0 * e0 - e1 * e1 + e2 * e2 - e3 * e3);
653 matrix[1][2] = (float)(2.0f * (e2 * e3 + e0 * e1));
654 matrix[2][0] = (float)(2.0f * (e3 * e1 + e0 * e2));
655 matrix[2][1] = (float)(2.0f * (e3 * e2 - e0 * e1));
656 matrix[2][2] = (float)(e0 * e0 - e1 * e1 - e2 * e2 + e3 * e3);
657 // keep in mind matrix[0][2] is M13 gives an e0 e2 piece, etc
658 }
659 if (is_xyz){
660 matrix[0][0] = (float)(cytilt*cztilt);
661 matrix[0][1] = (float)(cxtilt*sztilt+sxtilt*sytilt*cztilt);
662 matrix[0][2] = (float)(sxtilt*sztilt-cxtilt*sytilt*cztilt);
663 matrix[1][0] = (float)(-cytilt*sztilt);
664 matrix[1][1] = (float)(cxtilt*cztilt-sxtilt*sytilt*sztilt);
665 matrix[1][2] = (float)(sxtilt*cztilt+cxtilt*sytilt*sztilt);
666 matrix[2][0] = (float)(sytilt);
667 matrix[2][1] = (float)(-sxtilt*cytilt);
668 matrix[2][2] = (float)(cxtilt*cytilt);
669 }
670
671 // Apply scale if it existed previously
672 if (scale != 1.0f) {
673 for(int i=0; i<3; ++i) {
674 for(int j=0; j<3; ++j) {
675 matrix[i][j] *= scale;
676 }
677 }
678 }
679
680 // Apply post x mirroring if it was applied previously
681 if ( x_mirror ) {
682 for(int j=0; j<3; ++j) {
683 matrix[0][j] *= -1.0f;
684 }
685 }
686}
static float hypot3(int x, int y, int z)
Euclidean distance function in 3D: f(x,y,z) = sqrt(x*x + y*y + z*z);.
Definition: util.h:827

References assert_valid_2d(), EMAN::EMConsts::deg2rad, detect_problem_keys(), EMAN::Dict::get_ci(), get_scale_and_mirror(), EMAN::Dict::has_key_ci(), EMAN::Util::hypot3(), InvalidParameterException, InvalidStringException, matrix, scale(), EMAN::Util::str_to_lower(), and UnexpectedBehaviorException.

Referenced by EMAN::PointArray::calc_transform(), get_hflip_transform(), get_vflip_transform(), icos_5_to_2(), EMAN::FourierReconstructor::preprocess_slice(), EMAN::FourierIterReconstructor::preprocess_slice(), EMAN::TransformProcessor::process(), EMAN::TransformProcessor::process_inplace(), EMAN::TestTomoImage::process_inplace(), EMAN::EMData::rotate_translate(), set_params(), set_params_inverse(), set_pre_trans(), set_rotation(), and tet_3_to_2().

◆ set_rotation() [2/2]

void Transform::set_rotation ( const Vec3f v)

Determine the rotation that would transform a vector pointing in the Z direction so that it points in the direction of the argument vector Automatically normalizes the vector.

Parameters
vthe direction you want to solve for

Definition at line 698 of file transform.cpp.

699{
700 if ( v[0] == 0 && v[1] == 0 && v[2] == 0 )
701 throw UnexpectedBehaviorException("Can't set rotation for the null vector");
702
703 Vec3f v1(v);
704 v1.normalize();
705
706 double theta = acos(v1[2]); // in radians
707 double psi = atan2(v1[1],-v1[0]);
708
709 Dict d;
710 d["theta"] = (double)EMConsts::rad2deg*theta;
711 d["psi"] = (double)EMConsts::rad2deg*psi;
712 d["phi"] = (double)0.0;
713 d["type"] = "spider";
714
715 set_rotation(d);
716
717
718}

References EMAN::Vec3< Type >::normalize(), EMAN::EMConsts::rad2deg, set_rotation(), and UnexpectedBehaviorException.

◆ set_scale()

void Transform::set_scale ( const float &  scale)

Set the scale.

Parameters
scalethe amount to scale by

Definition at line 1123 of file transform.cpp.

1123 {
1124 if (new_scale <= 0) {
1125 throw InvalidValueException(new_scale,"The scale factor in a Transform object must be positive and non zero");
1126 }
1127 // Transform = MTSR (Mirroring, Translation, Scaling, Rotate)
1128 // So changing the scale boils down to this....
1129
1130 float old_scale = get_scale();
1131
1132 float n_scale = new_scale;
1134
1135 float corrected_scale = n_scale/old_scale;
1136 if ( corrected_scale != 1.0 ) {
1137 for(int i = 0; i < 3; ++i ) {
1138 for(int j = 0; j < 3; ++j ) {
1139 matrix[i][j] *= corrected_scale;
1140 }
1141 }
1142 }
1143}
#define InvalidValueException(val, desc)
Definition: exception.h:285

References EMAN::Util::apply_precision(), ERR_LIMIT, get_scale(), InvalidValueException, and matrix.

Referenced by EMAN::ScaleAlignerABS::align_using_base(), EMAN::FourierReconstructor::determine_slice_agreement(), EMAN::WienerFourierReconstructor::determine_slice_agreement(), get_rotation_transform(), EMAN::FourierReconstructor::insert_slice(), EMAN::WienerFourierReconstructor::insert_slice(), EMAN::FourierIterReconstructor::insert_slice(), EMAN::ScaleTransformProcessor::process(), EMAN::ScaleTransformProcessor::process_inplace(), EMAN::FourierReconstructor::projection(), refalifn(), rotate_origin_newBasis(), scale(), EMAN::EMData::scale(), set_params(), set_params_inverse(), and set_pre_trans().

◆ set_trans() [1/3]

void Transform::set_trans ( const float &  x,
const float &  y,
const float &  z = 0 
)

Set the post translation component.

Parameters
xthe x translation
ythe y translation
zthe z translation

Definition at line 1036 of file transform.cpp.

1037{
1038 bool x_mirror = get_mirror();
1039
1040 if (x_mirror) matrix[0][3] = -x;
1041 else matrix[0][3] = x;
1042 matrix[1][3] = y;
1043 matrix[2][3] = z;
1044}
#define y(i, j)
Definition: projector.cpp:1516
#define x(i)
Definition: projector.cpp:1517

References get_mirror(), matrix, x, and y.

Referenced by EMAN::PointArray::calc_transform(), EMAN::FourierReconstructor::determine_slice_agreement(), EMAN::WienerFourierReconstructor::determine_slice_agreement(), EMAN::TestUtil::emobject_transformarray_to_py(), EMAN::TestUtil::get_debug_transform(), get_hflip_transform(), get_pre_trans(), get_pre_trans_2d(), get_rotation_transform(), EMAN::HSym::get_sym(), get_vflip_transform(), EMAN::FourierReconstructor::insert_slice(), EMAN::WienerFourierReconstructor::insert_slice(), EMAN::FourierIterReconstructor::insert_slice(), EMAN::TransformProcessor::process(), EMAN::TransformProcessor::process_inplace(), EMAN::PhaseToMassCenterProcessor::process_inplace(), EMAN::ToCenterProcessor::process_inplace(), EMAN::ToMassCenterProcessor::process_inplace(), EMAN::TestTomoImage::process_inplace(), EMAN::FourierReconstructor::projection(), refalifn(), refalin3d_perturbquat(), rotate_origin_newBasis(), EMAN::EMData::rotate_translate(), set_params(), set_params_inverse(), set_pre_trans(), EMAN::RT3DLocalTreeAligner::testort(), EMAN::RT3DTreeAligner::testort(), EMAN::EMData::translate(), translate_newBasis(), EMAN::RT3DGridAligner::xform_align_nbest(), EMAN::RT3DSphereAligner::xform_align_nbest(), EMAN::RT2DTreeAligner::xform_align_nbest(), EMAN::RT2Dto3DTreeAligner::xform_align_nbest(), EMAN::RT3DTreeAligner::xform_align_nbest(), and EMAN::RT3DLocalTreeAligner::xform_align_nbest().

◆ set_trans() [2/3]

void EMAN::Transform::set_trans ( const Vec2f v)
inline

Set the post translation component using a Vec2f.

Parameters
vthe 2D translation vector

Definition at line 221 of file transform.h.

221{ set_trans(v[0],v[1]); }

References set_trans().

Referenced by set_trans().

◆ set_trans() [3/3]

void EMAN::Transform::set_trans ( const Vec3f v)
inline

Set the post translation component using a Vec3f.

Parameters
vthe 3D translation vector

Definition at line 216 of file transform.h.

216{ set_trans(v[0],v[1],v[2]); }

References set_trans().

Referenced by set_trans().

◆ tet_3_to_2()

Transform Transform::tet_3_to_2 ( )
static

Get the transform that moves any tetrahedron generated by eman2 so that it matches the 2-2-2 (MRC, FREALIGN) convention.

Returns
a Transforms, alt=54.73561, phi=45

Doctor Phil says: AltAngle= acos(-1/3.0)*90/pi;

Definition at line 85 of file transform.cpp.

85 {
86 Transform t;
87 Dict d;
88 d["type"] = "eman";
89 d["phi"] = 45.0f;
90 d["az"] = 0.0f;
91 d["alt"] = 54.73561f; // 3 fold to a 2 fold
95 t.set_rotation(d);
96 return t;
97}

References set_rotation().

◆ to_identity()

void Transform::to_identity ( )

Force the internal matrix to become the identity.

Definition at line 207 of file transform.cpp.

208{
209// transform_type = UNKNOWN;
210 for(int i=0; i<3; ++i) {
211 for(int j=0; j<4; ++j) {
212 if(i==j) {
213 matrix[i][j] = 1;
214 }
215 else {
216 matrix[i][j] = 0;
217 }
218 }
219 }
220}

References matrix.

Referenced by Transform().

◆ transform() [1/4]

Vec2f EMAN::Transform::transform ( const float &  x,
const float &  y 
) const
inline

Transform 2D coordinates using the internal transformation matrix.

Parameters
xthe x coordinate of the transformed point
ythe y coordinate of the transformed point
Returns
the transformed vector

Definition at line 417 of file transform.h.

417 {
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 }
Vec2< float > Vec2f
Definition: vec3.h:1071

References matrix, x, and y.

Referenced by EMAN::EMData::get_rotated_clip(), EMAN::operator*(), EMAN::SegmentSubunitProcessor::process_inplace(), and transform().

◆ transform() [2/4]

Vec3f EMAN::Transform::transform ( const float &  x,
const float &  y,
const float &  z 
) const
inline

Transform 3D coordinates using the internal transformation matrix.

Parameters
xthe x coordinate of the transformed point
ythe y coordinate of the transformed point
zthe z coordinate of the transformed point
Returns
the transformed vector

Definition at line 440 of file transform.h.

440 {
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 }

References matrix, x, and y.

◆ transform() [3/4]

template<typename Type >
Vec2f EMAN::Transform::transform ( const Vec2< Type > &  v) const
inline

Transform a 2D vector using the internal transformation matrix.

Parameters
va two dimensional vector to be transformed
Returns
the transformed vector

Definition at line 430 of file transform.h.

430 {
431 return transform(v[0],v[1]);
432 }
Vec2f transform(const float &x, const float &y) const
Transform 2D coordinates using the internal transformation matrix.
Definition: transform.h:417

References transform().

◆ transform() [4/4]

template<typename Type >
Vec3f EMAN::Transform::transform ( const Vec3< Type > &  v) const
inline

Transform a 3D vector using the internal transformation matrix.

Parameters
va three dimensional vector to be transformed
Returns
the transformed vector

Definition at line 454 of file transform.h.

454 {
455// assert_consistent_type(THREED); // Transform does the assertion
456 return transform(v[0],v[1],v[2]);
457 }

References transform().

◆ translate() [1/4]

void Transform::translate ( const float &  tx,
const float &  ty,
const float &  tz = 0 
)

Increment the current translation by tx, ty, tz.

Parameters
txthe x incrementation
tythe y incrementation
tzthe z incrementation

Definition at line 1063 of file transform.cpp.

1064{
1065 bool x_mirror = get_mirror();
1066 if (x_mirror) matrix[0][3] = -matrix[0][3] + tx;
1067 else matrix[0][3] = matrix[0][3] + tx;
1068 matrix[1][3] = matrix[1][3] + ty;
1069 matrix[2][3] = matrix[2][3] + tz;
1070}

References get_mirror(), and matrix.

Referenced by translate_newBasis().

◆ translate() [2/4]

void EMAN::Transform::translate ( const Transform tcs,
const Vec3f v 
)
inline

Increment the current translation using vec3f& v and a non standard basis.

Parameters
vthe 3D translation vector

Definition at line 258 of file transform.h.

258{ translate_newBasis(tcs, v[0],v[1],v[2]); }
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

References translate_newBasis().

◆ translate() [3/4]

void EMAN::Transform::translate ( const Vec2f v)
inline

Increment the current translation using vec2f& v.

Parameters
vthe @D translation vector

Definition at line 243 of file transform.h.

243{ translate(v[0],v[1]); }
void translate(const float &tx, const float &ty, const float &tz=0)
Increment the current translation by tx, ty, tz.
Definition: transform.cpp:1063

References translate().

Referenced by translate().

◆ translate() [4/4]

void EMAN::Transform::translate ( const Vec3f v)
inline

Increment the current translation using vec3f& v.

Parameters
vthe 3D translation vector

Definition at line 238 of file transform.h.

238{ translate(v[0],v[1],v[2]); }

References translate().

Referenced by translate().

◆ translate_newBasis()

void Transform::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 remove the effect of tcs when a composite transfrom tcs*t (where t is the current transform) This function is used in the scenegraph.

Parameters
tcsthe transform specifing the new basis vectors
txthe x incrementation
tythe y incrementation
tzthe z incrementation

Definition at line 1072 of file transform.cpp.

1073{
1074 //Get the rotational inverse
1075 Transform tcsinv = Transform(tcs);
1076 tcsinv.set_trans(0.0, 0.0, 0.0);
1077 tcsinv.invert();
1078
1079 //Now move the coordinate system
1080 Transform temp = Transform();
1081 temp.set_trans(tx, ty, tz);
1082 Transform nb_trans = tcsinv*temp;
1083
1084 translate(nb_trans.get_trans());
1085
1086}

References get_trans(), invert(), set_trans(), Transform(), and translate().

Referenced by translate().

◆ transpose()

Transform Transform::transpose ( ) const

Get the transpose of this transformation matrix.

Returns
the transpose of this transformation matrix

Definition at line 1346 of file transform.cpp.

1346 {
1347 Transform t(*this);
1348 t.transpose_inplace();
1349 return t;
1350}

References transpose_inplace().

Referenced by EMAN::PDBReader::right_transform().

◆ transpose_inplace()

void Transform::transpose_inplace ( )

Get the transpose of this transformation matrix.

Definition at line 1333 of file transform.cpp.

1333 {
1334 float tempij;
1335 for (int i = 0; i < 3; i++) {
1336 for (int j = 0; j < i; j++) {
1337 if (i != j) {
1338 tempij= matrix[i][j];
1339 matrix[i][j] = matrix[j][i];
1340 matrix[j][i] = tempij;
1341 }
1342 }
1343 }
1344}

References matrix.

Referenced by transpose().

Member Data Documentation

◆ ERR_LIMIT

const float Transform::ERR_LIMIT = 0.000001f
static

◆ matrix

float EMAN::Transform::matrix[3][4]
private

◆ permissable_2d_not_rot

vector< string > Transform::permissable_2d_not_rot
staticprivate

This map is used to validate keys in the argument maps - e.g. if the type is 2d and the angle is not "alpha" then we should throw.

Definition at line 508 of file transform.h.

Referenced by detect_problem_keys(), and init_permissable_keys().

◆ permissable_3d_not_rot

vector< string > Transform::permissable_3d_not_rot
staticprivate

Definition at line 509 of file transform.h.

Referenced by detect_problem_keys(), and init_permissable_keys().

◆ permissable_rot_keys

map< string, vector< string > > Transform::permissable_rot_keys
staticprivate

Definition at line 510 of file transform.h.

Referenced by detect_problem_keys(), and init_permissable_keys().


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