transform.cpp

Go to the documentation of this file.
00001 
00005 /*
00006  * Author: Steven Ludtke, 04/10/2003 (sludtke@bcm.edu)
00007  * Copyright (c) 2000-2006 Baylor College of Medicine
00008  *
00009  * This software is issued under a joint BSD/GNU license. You may use the
00010  * source code in this file under either license. However, note that the
00011  * complete EMAN2 and SPARX software packages have some GPL dependencies,
00012  * so you are responsible for compliance with the licenses of these packages
00013  * if you opt to use BSD licensing. The warranty disclaimer below holds
00014  * in either instance.
00015  *
00016  * This complete copyright notice must be included in any revised version of the
00017  * source code. Additional authorship citations may be added, but existing
00018  * author citations must be preserved.
00019  *
00020  * This program is free software; you can redistribute it and/or modify
00021  * it under the terms of the GNU General Public License as published by
00022  * the Free Software Foundation; either version 2 of the License, or
00023  * (at your option) any later version.
00024  *
00025  * This program is distributed in the hope that it will be useful,
00026  * but WITHOUT ANY WARRANTY; without even the implied warranty of
00027  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
00028  * GNU General Public License for more details.
00029  *
00030  * You should have received a copy of the GNU General Public License
00031  * along with this program; if not, write to the Free Software
00032  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
00033  *
00034  * */
00035 
00036 #include "transform.h"
00037 #include "util.h"
00038 #include "emobject.h"
00039 #include <cctype> // for std::tolower
00040 #include <cstring>  // for memcpy
00041 #include "symmetry.h"
00042 using namespace EMAN;
00043 
00044 #ifdef WIN32
00045 #ifndef M_PI
00046 #define M_PI 3.14159265358979323846
00047 #endif
00048 #endif
00049 
00050 #include <algorithm> // for std::transform
00051 
00052 #include <gsl_matrix.h>
00053 #include <gsl_blas.h>
00054 #include <gsl_linalg.h>
00055 
00056 #include <ostream>
00057 using std::ostream_iterator;
00058 
00059 //const float Transform3D::ERR_LIMIT = 0.000001f;
00060 
00061 const float Transform::ERR_LIMIT = 0.000001f;
00062 
00063 vector<string>  Transform::permissable_2d_not_rot;
00064 vector<string>  Transform::permissable_3d_not_rot;
00065 map<string,vector<string> >  Transform::permissable_rot_keys;
00066 
00067 Transform Transform::icos_5_to_2() {
00068         Transform t;
00069         Dict d;
00070         d["type"] = "eman";
00071         d["phi"] = 0;
00072         d["az"] = 90.0f;
00073         d["alt"] = 31.717474; // 5 fold to the nearest 2-fold
00085         t.set_rotation(d);
00086         return t;
00087 }
00088 
00089 Transform Transform::tet_3_to_2() {
00090         Transform t;
00091         Dict d;
00092         d["type"] = "eman";
00093         d["phi"] = 45.0f;
00094         d["az"] = 0.0f;
00095         d["alt"] = 54.73561f; // 3 fold to a 2 fold
00099         t.set_rotation(d);
00100         return t;
00101 }
00102 
00103 
00104 Transform::Transform()
00105 {
00106         to_identity();
00107 }
00108 
00109 Transform::Transform( const Transform& that )
00110 {
00111         *this = that;
00112 }
00113 
00114 Transform& Transform::operator=(const Transform& that ) {
00115         if (this != &that ) {
00116                 memcpy(matrix,that.matrix,12*sizeof(float));
00117 //              transform_type = that.transform_type;
00118         }
00119         return *this;
00120 }
00121 
00122 bool Transform::operator==(const Transform& rhs) const{
00123         if (memcmp(this->matrix, rhs.matrix, 3*4*sizeof(float)) == 0) {
00124                 return true;
00125         }
00126         else {
00127                 return false;
00128         }
00129 }
00130 
00131 bool Transform::operator!=(const Transform& rhs) const{
00132         return !(operator==(rhs));
00133 }
00134 
00135 Transform::Transform(const Dict& d)  {
00136         to_identity();
00137         set_params(d);
00138 }
00139 
00140 
00141 Transform::Transform(const float array[12]) {
00142         memcpy(matrix,array,12*sizeof(float));
00143 }
00144 
00145 Transform::Transform(const vector<float> array)
00146 {
00147         set_matrix(array);
00148 }
00149 
00150 void Transform::set_matrix(const vector<float>& v)
00151 {
00152         if (v.size() != 12 ) throw InvalidParameterException("The construction array must be of size 12");
00153 
00154         for(int i=0; i<3; ++i) {
00155                 for(int j=0; j<4; ++j) {
00156                         matrix[i][j] = v[i*4+j];
00157                 }
00158         }
00159 }
00160 
00161 void Transform::copy_matrix_into_array(float* const array) const {
00162 
00163         int idx = 0;
00164         for(int i=0; i<3; ++i) {
00165                 for(int j=0; j<4; ++j) {
00166                         array[idx] = matrix[i][j];
00167                         idx ++;
00168                 }
00169         }
00170 }
00171 
00172 vector<float> Transform::get_matrix() const
00173 {
00174         vector<float> ret(12);
00175         for(int i=0; i<3; ++i) {
00176                 for(int j=0; j<4; ++j) {
00177                         ret[i*4+j] = matrix[i][j];
00178                 }
00179         }
00180         return ret;
00181 }
00182 
00183 vector<float> Transform::get_matrix_4x4() const
00184 {
00185         vector<float> ret(16);
00186         for(int i=0; i<3; ++i) {
00187                 for(int j=0; j<4; ++j) {
00188                         ret[i*4+j] = matrix[i][j];
00189                 }
00190         }
00191         ret[12] = 0.0; 
00192         ret[13] = 0.0;
00193         ret[14] = 0.0;
00194         ret[15] = 1.0;
00195         
00196         return ret;
00197 }
00198 void Transform::to_identity()
00199 {
00200 //      transform_type = UNKNOWN;
00201         for(int i=0; i<3; ++i) {
00202                 for(int j=0; j<4; ++j) {
00203                         if(i==j) {
00204                                 matrix[i][j] = 1;
00205                         }
00206                         else {
00207                                 matrix[i][j] = 0;
00208                         }
00209                 }
00210         }
00211 }
00212 
00213 bool Transform::is_identity() const {
00214         for(int i=0; i<3; ++i) {
00215                 for(int j=0; j<4; ++j) {
00216                         float c = matrix[i][j];
00217                         Util::apply_precision(c,ERR_LIMIT);
00218                         if(i==j) {
00219                                 if (c != 1.0) return false;
00220                         }
00221                         else {
00222                                 if (c != 0.0) return false;
00223                         }
00224                 }
00225         }
00226         return true;
00227 }
00228 
00229 bool Transform::is_rot_identity() const {
00230         for(int i=0; i<3; ++i) {
00231                 for(int j=0; j<3; ++j) {
00232                         float c = matrix[i][j];
00233                         Util::apply_precision(c,ERR_LIMIT);
00234                         if(i==j) {
00235                                 if (c != 1.0) return false;
00236                         }
00237                         else {
00238                                 if (c != 0.0) return false;
00239                         }
00240                 }
00241         }
00242         return true;
00243 }
00244 
00245 void Transform::set_params(const Dict& d) {
00246         detect_problem_keys(d);
00247 
00248         if (d.has_key_ci("type") ) set_rotation(d);
00249 
00250         if (d.has_key_ci("scale")) {
00251                 float scale = static_cast<float>(d.get_ci("scale"));
00252                 set_scale(scale);
00253         }
00254 
00255         float dx=0,dy=0,dz=0;
00256 
00257         if (d.has_key_ci("tx")) dx = static_cast<float>(d.get_ci("tx"));
00258         if (d.has_key_ci("ty")) dy = static_cast<float>(d.get_ci("ty"));
00259         if (d.has_key_ci("tz")) dz = static_cast<float>(d.get_ci("tz"));
00260 
00261         if ( dx != 0.0 || dy != 0.0 || dz != 0.0 ) {
00262                 set_trans(dx,dy,dz);
00263         }
00264 
00265         if (d.has_key_ci("mirror")) {
00266                 EMObject e = d.get_ci("mirror");
00267                 if ( (e.get_type() != EMObject::BOOL ) && (e.get_type() != EMObject::INT ) && (e.get_type() != EMObject::UNSIGNEDINT ) )
00268                         throw InvalidParameterException("Error, mirror must be a bool or an int");
00269 
00270                 bool mirror = static_cast<bool>(e);
00271                 set_mirror(mirror);
00272         }
00273 }
00274 
00275 
00276 void Transform::init_permissable_keys()
00277 {
00278 
00279         permissable_2d_not_rot.push_back("tx");
00280         permissable_2d_not_rot.push_back("ty");
00281         permissable_2d_not_rot.push_back("scale");
00282         permissable_2d_not_rot.push_back("mirror");
00283         permissable_2d_not_rot.push_back("type");
00284 
00285         permissable_3d_not_rot.push_back("tx");
00286         permissable_3d_not_rot.push_back("ty");
00287         permissable_3d_not_rot.push_back("tz");
00288         permissable_3d_not_rot.push_back("scale");
00289         permissable_3d_not_rot.push_back("mirror");
00290         permissable_3d_not_rot.push_back("type");
00291 
00292         vector<string> tmp;
00293         tmp.push_back("alpha");
00294         permissable_rot_keys["2d"] = tmp;
00295 
00296         tmp.clear();
00297         tmp.push_back("alt");
00298         tmp.push_back("az");
00299         tmp.push_back("phi");
00300         permissable_rot_keys["eman"] = tmp;
00301 
00302         tmp.clear();
00303         tmp.push_back("psi");
00304         tmp.push_back("theta");
00305         tmp.push_back("phi");
00306         permissable_rot_keys["spider"] = tmp;
00307 
00308         tmp.clear();
00309         tmp.push_back("alpha");
00310         tmp.push_back("beta");
00311         tmp.push_back("gamma");
00312         permissable_rot_keys["imagic"] = tmp;
00313 
00314         tmp.clear();
00315         tmp.push_back("ztilt");
00316         tmp.push_back("xtilt");
00317         tmp.push_back("ytilt");
00318         permissable_rot_keys["xyz"] = tmp;
00319 
00320         tmp.clear();
00321         tmp.push_back("phi");
00322         tmp.push_back("theta");
00323         tmp.push_back("omega");
00324         permissable_rot_keys["mrc"] = tmp;
00325 
00326         tmp.clear();
00327         tmp.push_back("e0");
00328         tmp.push_back("e1");
00329         tmp.push_back("e2");
00330         tmp.push_back("e3");
00331         permissable_rot_keys["quaternion"] = tmp;
00332 
00333         tmp.clear();
00334         tmp.push_back("n1");
00335         tmp.push_back("n2");
00336         tmp.push_back("n3");
00337         tmp.push_back("Omega");
00338         permissable_rot_keys["spin"] = tmp;
00339 
00340         tmp.clear();
00341         tmp.push_back("n1");
00342         tmp.push_back("n2");
00343         tmp.push_back("n3");
00344         tmp.push_back("q");
00345         permissable_rot_keys["sgirot"] = tmp;
00346 
00347         tmp.clear();
00348         tmp.push_back("m11");
00349         tmp.push_back("m12");
00350         tmp.push_back("m13");
00351         tmp.push_back("m21");
00352         tmp.push_back("m22");
00353         tmp.push_back("m23");
00354         tmp.push_back("m31");
00355         tmp.push_back("m32");
00356         tmp.push_back("m33");
00357         permissable_rot_keys["matrix"] = tmp;
00358 }
00359 
00360 
00361 void Transform::detect_problem_keys(const Dict& d) {
00362         if (permissable_rot_keys.size() == 0 ) {
00363                 init_permissable_keys();
00364         }
00365 
00366         vector<string> verification;
00367         vector<string> problem_keys;
00368         bool is_2d = false;
00369         if (d.has_key_ci("type") ) {
00370                 string type = Util::str_to_lower((string)d["type"]);
00371                 bool problem = false;
00372                 if (permissable_rot_keys.find(type) == permissable_rot_keys.end() ) {
00373                         problem_keys.push_back(type);
00374                         problem = true;
00375                 }
00376                 if ( !problem ) {
00377                         vector<string> perm = permissable_rot_keys[type];
00378                         std::copy(perm.begin(),perm.end(),back_inserter(verification));
00379 
00380                         if ( type == "2d" ) {
00381                                 is_2d = true;
00382                                 std::copy(permissable_2d_not_rot.begin(),permissable_2d_not_rot.end(),back_inserter(verification));
00383                         }
00384                 }
00385         }
00386         if ( !is_2d ) {
00387                 std::copy(permissable_3d_not_rot.begin(),permissable_3d_not_rot.end(),back_inserter(verification));
00388         }
00389 
00390         for (Dict::const_iterator it = d.begin(); it != d.end();  ++it) {
00391                 if ( std::find(verification.begin(),verification.end(), it->first) == verification.end() ) {
00392                         problem_keys.push_back(it->first);
00393                 }
00394         }
00395 
00396         if (problem_keys.size() != 0 ) {
00397                 string error;
00398                 if (problem_keys.size() == 1) {
00399                         error = "Transform Error: The \"" +problem_keys[0]+ "\" key is unsupported";
00400                 } else {
00401                         error = "Transform Error: The ";
00402                         for(vector<string>::const_iterator cit = problem_keys.begin(); cit != problem_keys.end(); ++cit ) {
00403                                 if ( cit != problem_keys.begin() ) {
00404                                         if (cit == (problem_keys.end() -1) ) error += " and ";
00405                                         else error += ", ";
00406                                 }
00407                                 error += "\"";
00408                                 error += *cit;
00409                                 error += "\"";
00410                         }
00411                         error += " keys are unsupported";
00412                 }
00413                 throw InvalidParameterException(error);
00414         }
00415 }
00416 
00417 void Transform::set_params_inverse(const Dict& d) {
00418         detect_problem_keys(d);
00419 
00420         if (d.has_key_ci("type") ) set_rotation(d);
00421 
00422         float dx=0,dy=0,dz=0;
00423         if (d.has_key_ci("tx")) dx = static_cast<float>(d.get_ci("tx"));
00424         if (d.has_key_ci("ty")) dy = static_cast<float>(d.get_ci("ty"));
00425         if (d.has_key_ci("tz")) dz = static_cast<float>(d.get_ci("tz"));
00426 
00427         if ( (dx != 0.0 || dy != 0.0 || dz != 0.0) && d.has_key_ci("type") ) {
00428                 Transform pre_trans;
00429                 pre_trans.set_trans(dx,dy,dz);
00430 
00431                 Transform tmp;
00432                 tmp.set_rotation(d);
00433 
00434                 if (d.has_key_ci("scale")) {
00435                         float scale = static_cast<float>(d.get_ci("scale"));
00436                         tmp.set_scale(scale);
00437                 }
00438 
00439                 Transform solution_trans = tmp*pre_trans;
00440 
00441                 if (d.has_key_ci("scale")) {
00442                         Transform tmp;
00443                         float scale = static_cast<float>(d.get_ci("scale"));
00444                         tmp.set_scale(scale);
00445                         solution_trans = solution_trans*tmp;
00446                 }
00447 
00448                 tmp = Transform();
00449                 tmp.set_rotation(d);
00450                 solution_trans = solution_trans*tmp;
00451                 set_trans(solution_trans.get_trans());
00452         }
00453 
00454         if (d.has_key_ci("scale")) {
00455                 float scale = static_cast<float>(d.get_ci("scale"));
00456                 set_scale(scale);
00457         }
00458 
00459         if (d.has_key_ci("mirror")) {
00460                 EMObject e = d.get_ci("mirror");
00461                 if ( (e.get_type() != EMObject::BOOL ) && (e.get_type() != EMObject::INT ) && (e.get_type() != EMObject::UNSIGNEDINT ) )
00462                         throw InvalidParameterException("Error, mirror must be a bool or an int");
00463 
00464                 bool mirror = static_cast<bool>(e);
00465                 set_mirror(mirror);
00466         }
00467         invert();
00468 }
00469 
00470 
00471 Dict Transform::get_params(const string& euler_type) const {
00472         Dict params = get_rotation(euler_type);
00473 
00474         Vec3f v = get_trans();
00475         params["tx"] = v[0]; params["ty"] = v[1];
00476 
00477         string type = Util::str_to_lower(euler_type);
00478         if ( type != "2d") params["tz"] = v[2];
00479 
00480         float scale = get_scale();
00481         params["scale"] = scale;
00482 
00483         bool mirror = get_mirror();
00484         params["mirror"] = mirror;
00485 
00486         return params;
00487 }
00488 
00489 
00490 
00491 Dict Transform::get_params_inverse(const string& euler_type) const {
00492         Transform inv(inverse());
00493 
00494         Dict params = inv.get_rotation(euler_type);
00495         Vec3f v = inv.get_pre_trans();
00496         params["tx"] = v[0]; params["ty"] = v[1];
00497 
00498         string type = Util::str_to_lower(euler_type);
00499         if ( type != "2d") params["tz"] = v[2];
00500 
00501         float scale = inv.get_scale();
00502         params["scale"] = scale;
00503 
00504         bool mirror = inv.get_mirror();
00505         params["mirror"] = mirror;
00506 
00507         return params;
00508 }
00509 
00510 
00511 void Transform::set_rotation(const Dict& rotation)
00512 {
00513         detect_problem_keys(rotation);
00514         string euler_type;
00515 
00516         if (!rotation.has_key_ci("type") ){
00517                         throw InvalidParameterException("argument dictionary does not contain the type key");
00518         }
00519 
00520         euler_type = static_cast<string>(rotation.get_ci("type"));// Warning, will throw
00521 
00522 
00523         double e0=0; double e1=0; double e2=0; double e3=0;
00524         double Omega=0;
00525         double az  = 0;
00526         double alt = 0;
00527         double phi = 0;
00528         double cxtilt = 0;
00529         double sxtilt = 0;
00530         double cytilt = 0;
00531         double sytilt = 0;
00532         double cztilt = 0;
00533         double sztilt = 0;
00534         bool is_quaternion = 0;
00535         bool is_matrix = 0;
00536         bool is_xyz = 0;
00537 
00538         bool x_mirror;
00539         float scale;
00540         // Get these before anything changes so we can apply them again after the rotation is set
00541         get_scale_and_mirror(scale,x_mirror);
00542         if (scale == 0) throw UnexpectedBehaviorException("The determinant of the Transform is 0. This is unexpected.");
00543 
00544         string type = Util::str_to_lower(euler_type);
00545         if (type == "2d") {
00546                 assert_valid_2d();
00547                 az  = 0;
00548                 alt = 0;
00549                 phi = (double)rotation["alpha"] ;
00550         } else if ( type == "eman" ) {
00551 //              validate_and_set_type(THREED);
00552                 az  = (double)rotation["az"] ;
00553                 alt = (double)rotation["alt"]  ;
00554                 phi = (double)rotation["phi"] ;
00555         } else if ( type == "imagic" ) {
00556 //              validate_and_set_type(THREED);
00557                 az  = (double)rotation["alpha"] ;
00558                 alt = (double)rotation["beta"]  ;
00559                 phi = (double)rotation["gamma"] ;
00560         } else if ( type == "spider" ) {
00561 //              validate_and_set_type(THREED);
00562                 az =  (double)rotation["phi"]    + 90.0;
00563                 alt = (double)rotation["theta"] ;
00564                 phi = (double)rotation["psi"]    - 90.0;
00565         } else if ( type == "xyz" ) {
00566 //              validate_and_set_type(THREED);
00567                 is_xyz = 1;
00568                 cxtilt = cos(EMConsts::deg2rad*(double)rotation["xtilt"]);
00569                 sxtilt = sin(EMConsts::deg2rad*(double)rotation["xtilt"]);
00570                 cytilt = cos(EMConsts::deg2rad*(double)rotation["ytilt"]);
00571                 sytilt = sin(EMConsts::deg2rad*(double)rotation["ytilt"]);
00572                 cztilt = cos(EMConsts::deg2rad*(double)rotation["ztilt"]);
00573                 sztilt = sin(EMConsts::deg2rad*(double)rotation["ztilt"]);
00574         } else if ( type == "mrc" ) {
00575 //              validate_and_set_type(THREED);
00576                 az  = (double)rotation["phi"]   + 90.0f ;
00577                 alt = (double)rotation["theta"] ;
00578                 phi = (double)rotation["omega"] - 90.0f ;
00579         } else if ( type == "quaternion" ) {
00580 //              validate_and_set_type(THREED);
00581                 is_quaternion = 1;
00582                 e0 = (double)rotation["e0"];
00583                 e1 = (double)rotation["e1"];
00584                 e2 = (double)rotation["e2"];
00585                 e3 = (double)rotation["e3"];
00586         } else if ( type == "spin" ) {
00587 //              validate_and_set_type(THREED);
00588                 is_quaternion = 1;
00589                 Omega = (double)rotation["Omega"];
00590                 e0 = cos(Omega*EMConsts::deg2rad/2.0);
00591                 e1 = sin(Omega*EMConsts::deg2rad/2.0) * (double)rotation["n1"];
00592                 e2 = sin(Omega*EMConsts::deg2rad/2.0) * (double)rotation["n2"];
00593                 e3 = sin(Omega*EMConsts::deg2rad/2.0) * (double)rotation["n3"];
00594         } else if ( type == "sgirot" ) {
00595 //              validate_and_set_type(THREED);
00596                 is_quaternion = 1;
00597                 Omega = (double)rotation["q"] ;
00598                 e0 = cos(Omega*EMConsts::deg2rad/2.0);
00599                 e1 = sin(Omega*EMConsts::deg2rad/2.0) * (double)rotation["n1"];
00600                 e2 = sin(Omega*EMConsts::deg2rad/2.0) * (double)rotation["n2"];
00601                 e3 = sin(Omega*EMConsts::deg2rad/2.0) * (double)rotation["n3"];
00602         } else if ( type == "matrix" ) {
00603                 is_matrix = 1;
00604                 matrix[0][0] = (float)rotation["m11"];
00605                 matrix[0][1] = (float)rotation["m12"];
00606                 matrix[0][2] = (float)rotation["m13"];
00607                 matrix[1][0] = (float)rotation["m21"];
00608                 matrix[1][1] = (float)rotation["m22"];
00609                 matrix[1][2] = (float)rotation["m23"];
00610                 matrix[2][0] = (float)rotation["m31"];
00611                 matrix[2][1] = (float)rotation["m32"];
00612                 matrix[2][2] = (float)rotation["m33"];
00613         } else {
00614 //              transform_type = UNKNOWN;
00615                 throw InvalidStringException(euler_type, "unknown Euler Type");
00616         }
00617 
00618         double azp  =  az*EMConsts::deg2rad;
00619         double altp = alt*EMConsts::deg2rad;
00620         double phip = phi*EMConsts::deg2rad;
00621 
00622         if (!is_quaternion && !is_matrix && !is_xyz) {
00623                 matrix[0][0] =  (float)(cos(phip)*cos(azp) - cos(altp)*sin(azp)*sin(phip));
00624                 matrix[0][1] =  (float)(cos(phip)*sin(azp) + cos(altp)*cos(azp)*sin(phip));
00625                 matrix[0][2] =  (float)(sin(altp)*sin(phip));
00626                 matrix[1][0] =  (float)(-sin(phip)*cos(azp) - cos(altp)*sin(azp)*cos(phip));
00627                 matrix[1][1] =  (float)(-sin(phip)*sin(azp) + cos(altp)*cos(azp)*cos(phip));
00628                 matrix[1][2] =  (float)(sin(altp)*cos(phip));
00629                 matrix[2][0] =  (float)(sin(altp)*sin(azp));
00630                 matrix[2][1] =  (float)(-sin(altp)*cos(azp));
00631                 matrix[2][2] =  (float)cos(altp);
00632         }
00633         if (is_quaternion){
00634                 matrix[0][0] = (float)(e0 * e0 + e1 * e1 - e2 * e2 - e3 * e3);
00635                 matrix[0][1] = (float)(2.0f * (e1 * e2 + e0 * e3));
00636                 matrix[0][2] = (float)(2.0f * (e1 * e3 - e0 * e2));
00637                 matrix[1][0] = (float)(2.0f * (e2 * e1 - e0 * e3));
00638                 matrix[1][1] = (float)(e0 * e0 - e1 * e1 + e2 * e2 - e3 * e3);
00639                 matrix[1][2] = (float)(2.0f * (e2 * e3 + e0 * e1));
00640                 matrix[2][0] = (float)(2.0f * (e3 * e1 + e0 * e2));
00641                 matrix[2][1] = (float)(2.0f * (e3 * e2 - e0 * e1));
00642                 matrix[2][2] = (float)(e0 * e0 - e1 * e1 - e2 * e2 + e3 * e3);
00643                 // keep in mind matrix[0][2] is M13 gives an e0 e2 piece, etc
00644         }
00645         if (is_xyz){
00646                 matrix[0][0] =  (float)(cytilt*cztilt);
00647                 matrix[0][1] =  (float)(cxtilt*sztilt+sxtilt*sytilt*cztilt);
00648                 matrix[0][2] =  (float)(sxtilt*sztilt-cxtilt*sytilt*cztilt);
00649                 matrix[1][0] =  (float)(-cytilt*sztilt);
00650                 matrix[1][1] =  (float)(cxtilt*cztilt-sxtilt*sytilt*sztilt);
00651                 matrix[1][2] =  (float)(sxtilt*cztilt+cxtilt*sytilt*sztilt);
00652                 matrix[2][0] =  (float)(sytilt);
00653                 matrix[2][1] =  (float)(-sxtilt*cytilt);
00654                 matrix[2][2] =  (float)(cxtilt*cytilt);
00655         }
00656         
00657         // Apply scale if it existed previously
00658         if (scale != 1.0f) {
00659                 for(int i=0; i<3; ++i) {
00660                         for(int j=0; j<3; ++j) {
00661                                 matrix[i][j] *= scale;
00662                         }
00663                 }
00664         }
00665 
00666         // Apply post x mirroring if it was applied previously
00667         if ( x_mirror ) {
00668                 for(int j=0; j<3; ++j) {
00669                         matrix[0][j] *= -1.0f;
00670                 }
00671         }
00672 }
00673 
00674 Transform Transform::get_rotation_transform() const
00675 {
00676         Transform ret(*this);
00677         ret.set_scale(1.0);
00678         ret.set_mirror(false);
00679         ret.set_trans(0,0,0);
00680         //ret.orthogonalize(); // ?
00681         return ret;
00682 }
00683 
00684 void Transform::set_rotation(const Vec3f & v)
00685 {
00686         if ( v[0] == 0 && v[1] == 0 && v[2] == 0 )
00687                 throw UnexpectedBehaviorException("Can't set rotation for the null vector");
00688 
00689         Vec3f v1(v);
00690         v1.normalize();
00691 
00692         double theta = acos(v1[2]); // in radians
00693         double psi = atan2(v1[1],-v1[0]);
00694 
00695         Dict d;
00696         d["theta"] = (double)EMConsts::rad2deg*theta;
00697         d["psi"] = (double)EMConsts::rad2deg*psi;
00698         d["phi"] = (double)0.0;
00699         d["type"] = "spider";
00700 
00701         set_rotation(d);
00702 
00703 
00704 }
00705 
00706 void Transform::rotate_origin(const Transform& by)
00707 {
00708         vector<float> multmatrix = by.get_matrix();
00709         // First Multiply and put the result in a temp matrix
00710         Transform result;
00711         for (int i=0; i<3; i++) {
00712                 for (int j=0; j<3; j++) {
00713                         result[i][j] = multmatrix[i*4]*matrix[0][j] +  multmatrix[i*4+1]*matrix[1][j] + multmatrix[i*4+2]*matrix[2][j];
00714                 }
00715         }
00716         //Then put the result from the tmep matrix in the original one
00717         for (int i=0; i<3; i++) {
00718                 for (int j=0; j<3; j++) {
00719                         matrix[i][j] = result[i][j];
00720                 }
00721         }
00722 }
00723 
00724 void Transform::rotate_origin_newBasis(const Transform& tcs, const float& Omega, const float& n1, const float& n2, const float& n3)
00725 {
00726         //Get the rotational inverse
00727         Transform tcsinv = Transform(tcs);
00728         tcsinv.set_trans(0.0, 0.0, 0.0);
00729         tcsinv.set_scale(1.0);
00730         tcsinv.invert();
00731         
00732         //Get the current rotation
00733         Transform temp = Transform();
00734         temp.set_trans(n1, n2, n3);
00735         Transform cc = tcsinv*temp;
00736         Vec3f cctrans = cc.get_trans();
00737         
00738         //set the right rotation
00739         Dict spinrot = Dict();
00740         spinrot["type"] = "spin";
00741         spinrot["Omega"] = Omega;
00742         spinrot["n1"] = cctrans[0];
00743         spinrot["n2"] = cctrans[1];
00744         spinrot["n3"] = cctrans[2];
00745         Transform rightrot = Transform(spinrot);
00746         rotate_origin(rightrot);
00747 }
00748 
00749 void Transform::rotate(const Transform& by)
00750 {
00751         vector<float> multmatrix = by.get_matrix();
00752         // First Multiply and put the result in a temp matrix
00753         Transform result;
00754         for (int i=0; i<3; i++) {
00755                 for (int j=0; j<4; j++) {
00756                         result[i][j] = multmatrix[i*4]*matrix[0][j] +  multmatrix[i*4+1]*matrix[1][j] + multmatrix[i*4+2]*matrix[2][j];
00757                 }
00758         }
00759         //Then put the result from the tmep matrix in the original one
00760         for (int i=0; i<3; i++) {
00761                 for (int j=0; j<4; j++) {
00762                         matrix[i][j] = result[i][j];
00763                 }
00764         }
00765 }
00766 
00767 Transform Transform::negate() const
00768 {
00769         Transform t(*this);
00770         for(unsigned int i = 0; i < 3; ++i) {
00771                 for(unsigned int j = 0; j < 4; ++j) {
00772                         t.set(i,j,t[i][j]*-1);
00773                 }
00774         }
00775         return t;
00776 }
00777 
00778 Transform Transform::get_hflip_transform() const {
00779 
00780         Dict rot = get_rotation("eman");
00781         rot["alt"] = 180.0f + static_cast<float>(rot["alt"]);
00782         rot["phi"] = 180.0f - static_cast<float>(rot["phi"]);
00783 
00784         Transform ret(*this); // Is the identity
00785         ret.set_rotation(rot);
00786 
00787         Vec3f trans = get_trans();
00788         trans[0] = -trans[0];
00789         ret.set_trans(trans);
00790 
00791 //      ret.set_mirror(self.get_mirror());
00792 
00793         return ret;
00794 }
00795 
00796 Transform Transform::get_vflip_transform() const {
00797 
00798         Dict rot = get_rotation("eman");
00799         rot["alt"] = 180.0f + static_cast<float>(rot["alt"]);
00800         rot["phi"] = - static_cast<float>(rot["phi"]);
00801 
00802         Transform ret(*this);
00803         ret.set_rotation(rot);
00804 
00805         Vec3f trans = get_trans();
00806         trans[1] = -trans[1];
00807         ret.set_trans(trans);
00808 
00809         return ret;
00810 }
00811 
00812 Dict Transform::get_rotation(const string& euler_type) const
00813 {
00814         Dict result;
00815 
00816         //float max = 1 - ERR_LIMIT;
00817         float scale;
00818         bool x_mirror;
00819         get_scale_and_mirror(scale,x_mirror);
00820         if (scale == 0) throw UnexpectedBehaviorException("The determinant of the Transform is 0. This is unexpected.");
00821 
00822         double cosalt = matrix[2][2]/scale;
00823         double x_mirror_scale = (x_mirror ? -1.0f : 1.0f);
00824         double inv_scale = 1.0f/scale;
00825 
00826         double az  = 0;
00827         double alt = 0;
00828         double phi = 0;
00829         double phiS = 0;  // like az  (but in SPIDER ZYZ)
00830         double psiS = 0;  // like phi (but in SPIDER ZYZ)
00831 
00832         // get alt, az, phi in EMAN convention
00833 
00834         if (cosalt >= 1) {  // that is, alt close to 0
00835                         alt = 0;
00836                         az = 0;
00837                         phi = (double)EMConsts::rad2deg * atan2(x_mirror_scale*matrix[0][1], x_mirror_scale*matrix[0][0]);
00838         } else if (cosalt <= -1) {  // that is, alt close to 180
00839                         alt = 180;
00840                         az = 0;
00841                         phi = (double)EMConsts::rad2deg * atan2(-x_mirror_scale*matrix[0][1], x_mirror_scale*matrix[0][0]);
00842         } else {   // for non exceptional cases:  0 < alt < 180
00843 
00844                 az  = (double)EMConsts::rad2deg * atan2(scale*matrix[2][0], -scale*matrix[2][1]);
00845 
00846                 if (matrix[2][2]==0.0)
00847                         alt = 90.0;
00848                 else
00849                         alt = (double)EMConsts::rad2deg * atan(sqrt((double)matrix[2][0]*matrix[2][0]+(double)matrix[2][1]*matrix[2][1])/fabs(matrix[2][2]));
00850 
00851                 if (matrix[2][2] * scale < 0)
00852                         alt = 180.0f-alt;
00853                 
00854                 phi = (double)EMConsts::rad2deg * atan2(x_mirror_scale*(double)matrix[0][2], (double)matrix[1][2]);
00855 
00856         } // ends separate cases: alt close to 0, 180, or neither
00857 
00858         phi = phi-360.0*floor(phi/360.0);
00859         az  = az -360.0*floor(az/360.0);
00860 
00861 //  get phiS, psiS (SPIDER)
00862         if (cosalt >= 1) {  // that is, alt close to 0
00863                 phiS = 0;
00864                 psiS = phi;
00865         } else if (cosalt <= -1) {  // that is, alt close to 180
00866                 phiS = 0;
00867                 psiS = phi + 180.0;
00868         } else {
00869                 phiS = az  - 90.0;
00870                 psiS = phi + 90.0;
00871         }
00872 
00873         phiS = phiS-360.0*floor(phiS/360.0);
00874         psiS = psiS-360.0*floor(psiS/360.0);
00875 
00876 //   do some quaternionic stuff here
00877 
00878         double nphi = (az-phi)/2.0;
00879     // The next is also e0
00880         double cosOover2 = cos((az+phi)*EMConsts::deg2rad/2.0) * cos(alt*EMConsts::deg2rad/2.0);
00881         double sinOover2 = sqrt(1 -cosOover2*cosOover2);
00882         double cosnTheta = sin((az+phi)*EMConsts::deg2rad/2.0) * cos(alt*EMConsts::deg2rad/2.0) / sqrt(1-cosOover2*cosOover2);
00883         double sinnTheta = sqrt(1-cosnTheta*cosnTheta);
00884         double n1 = sinnTheta*cos(nphi*EMConsts::deg2rad);
00885         double n2 = sinnTheta*sin(nphi*EMConsts::deg2rad);
00886         double n3 = cosnTheta;
00887         double xtilt = 0;
00888         double ytilt = 0;
00889         double ztilt = 0;
00890 
00891 
00892         if (cosOover2<0) {
00893                 cosOover2*=-1; n1 *=-1; n2*=-1; n3*=-1;
00894         }
00895 
00896         string type = Util::str_to_lower(euler_type);
00897 
00898         result["type"] = type;
00899         if (type == "2d") {
00900                 assert_valid_2d();
00901                 result["alpha"]  = phi;
00902         } else if (type == "eman") {
00903 //              assert_consistent_type(THREED);
00904                 result["az"]  = az;
00905                 result["alt"] = alt;
00906                 result["phi"] = phi;
00907         } else if (type == "imagic") {
00908 //              assert_consistent_type(THREED);
00909                 result["alpha"] = az;
00910                 result["beta"]  = alt;
00911                 result["gamma"] = phi;
00912         } else if (type == "spider") {
00913 //              assert_consistent_type(THREED);
00914                 result["phi"]   = phiS;  // The first Euler like az
00915                 result["theta"] = alt;
00916                 result["psi"]   = psiS;
00917         } else if (type == "mrc") {
00918 //              assert_consistent_type(THREED);
00919                 result["phi"]   = phiS;
00920                 result["theta"] = alt;
00921                 result["omega"] = psiS;
00922         } else if (type == "xyz") {               // need to double-check these 3 equations ********
00923 //              assert_consistent_type(THREED);
00924                 xtilt = atan2(-sin(EMConsts::deg2rad*phiS)*sin(EMConsts::deg2rad*alt),cos(EMConsts::deg2rad*alt));
00925                 ytilt = asin(  cos(EMConsts::deg2rad*phiS)*sin(EMConsts::deg2rad*alt));
00926                 ztilt = psiS*EMConsts::deg2rad - atan2(sin(xtilt), cos(xtilt) *sin(ytilt));
00927 
00928                 xtilt *= EMConsts::rad2deg; ytilt *= EMConsts::rad2deg; ztilt *= EMConsts::rad2deg;
00929                 xtilt = xtilt-360*.0*floor((xtilt+180.0)/360.0);
00930                 ytilt = ytilt-360*.0*floor((ytilt+180.0)/360.0);  //already in range [-90,90] but anyway...
00931                 ztilt = ztilt-360*.0*floor((ztilt+180.0)/360.0);
00932 
00933                 result["xtilt"]  = xtilt;
00934                 result["ytilt"]  = ytilt;
00935                 result["ztilt"]  = ztilt;
00936         } else if (type == "quaternion") {
00937 //              assert_consistent_type(THREED);
00938                 result["e0"] = cosOover2 ;
00939                 result["e1"] = sinOover2 * n1 ;
00940                 result["e2"] = sinOover2 * n2;
00941                 result["e3"] = sinOover2 * n3;
00942         } else if (type == "spin") {
00943 //              assert_consistent_type(THREED);
00944                 result["Omega"] = 2.0*EMConsts::rad2deg * acos(cosOover2);
00945                 result["n1"] = n1;
00946                 result["n2"] = n2;
00947                 result["n3"] = n3;
00948         } else if (type == "sgirot") {
00949 //              assert_consistent_type(THREED);
00950                 result["q"] = 2.0*EMConsts::rad2deg * acos(cosOover2);
00951                 result["n1"] = n1;
00952                 result["n2"] = n2;
00953                 result["n3"] = n3;
00954         } else if (type == "matrix") {
00955 //              assert_consistent_type(THREED);
00956                 result["m11"] = x_mirror_scale*matrix[0][0]*inv_scale;
00957                 result["m12"] = x_mirror_scale*matrix[0][1]*inv_scale;
00958                 result["m13"] = x_mirror_scale*matrix[0][2]*inv_scale;
00959                 result["m21"] = matrix[1][0]*inv_scale;
00960                 result["m22"] = matrix[1][1]*inv_scale;
00961                 result["m23"] = matrix[1][2]*inv_scale;
00962                 result["m31"] = matrix[2][0]*inv_scale;
00963                 result["m32"] = matrix[2][1]*inv_scale;
00964                 result["m33"] = matrix[2][2]*inv_scale;
00965         } else {
00966                 throw InvalidStringException(euler_type, "unknown Euler Type");
00967         }
00968 
00969         return result;
00970 }
00971 
00972 void Transform::set_trans(const float& x, const float& y, const float& z)
00973 {
00974         bool x_mirror = get_mirror();
00975 
00976         if (x_mirror) matrix[0][3] = -x;
00977         else matrix[0][3] = x;
00978         matrix[1][3] = y;
00979         matrix[2][3] = z;
00980 }
00981 
00982 Vec3f Transform::get_trans() const
00983 {
00984         // No type asserted
00985         bool x_mirror = get_mirror();
00986         Vec3f v;
00987         if (x_mirror) v[0] = -matrix[0][3];
00988         else v[0] = matrix[0][3];
00989         v[1] = matrix[1][3];
00990         v[2] = matrix[2][3];
00991 
00992         Util::apply_precision(v[0],ERR_LIMIT);
00993         Util::apply_precision(v[1],ERR_LIMIT);
00994         Util::apply_precision(v[2],ERR_LIMIT);
00995 
00996         return v;
00997 }
00998 
00999 void Transform::translate(const float& tx, const float& ty, const float& tz)
01000 {
01001         bool x_mirror = get_mirror();
01002         if (x_mirror) matrix[0][3] = -matrix[0][3] + tx;
01003         else matrix[0][3] = matrix[0][3] + tx;
01004         matrix[1][3] = matrix[1][3] + ty;
01005         matrix[2][3] = matrix[2][3] + tz;
01006 }
01007 
01008 void Transform::translate_newBasis(const Transform& tcs, const float& tx, const float& ty, const float& tz)
01009 {
01010         //Get the rotational inverse
01011         Transform tcsinv = Transform(tcs);
01012         tcsinv.set_trans(0.0, 0.0, 0.0);
01013         tcsinv.invert();
01014         
01015         //Now move the coordinate system
01016         Transform temp = Transform();
01017         temp.set_trans(tx, ty, tz);
01018         Transform nb_trans = tcsinv*temp;
01019         
01020         translate(nb_trans.get_trans());
01021         
01022 }
01023 
01024 Vec2f Transform::get_trans_2d() const
01025 {
01026         bool x_mirror = get_mirror();
01027         Vec2f v;
01028         if (x_mirror) v[0] = -matrix[0][3];
01029         else v[0] = matrix[0][3];
01030         v[1] = matrix[1][3];
01031         return v;
01032 }
01033 
01034 
01035 
01036 Vec3f Transform::get_pre_trans() const
01037 {
01038         Transform T(*this);
01039         T.set_trans(0,0,0);
01040         T.invert();
01041 
01042         Transform soln  = T*(*this);
01043 //      soln.printme();
01044         return soln.get_trans();
01045 }
01046 
01047 Vec2f Transform::get_pre_trans_2d() const
01048 {
01049         Transform T(*this);
01050         T.set_trans(0,0,0);
01051         T.invert();
01052 
01053         Transform soln  = T*(*this);
01054 //      soln.printme();
01055         return soln.get_trans_2d();
01056 }
01057 
01058 
01059 void Transform::set_scale(const float& new_scale) {
01060         if (new_scale <= 0) {
01061                 throw InvalidValueException(new_scale,"The scale factor in a Transform object must be positive and non zero");
01062         }
01063         // Transform = MTSR (Mirroring, Translation, Scaling, Rotate)
01064         // So changing the scale boils down to this....
01065 
01066         float old_scale = get_scale();
01067 
01068         float n_scale = new_scale;
01069         Util::apply_precision(n_scale,ERR_LIMIT);
01070 
01071         float corrected_scale = n_scale/old_scale;
01072         if ( corrected_scale != 1.0 ) {
01073                 for(int i = 0; i < 3;  ++i ) {
01074                         for(int j = 0; j < 3; ++j ) {
01075                                 matrix[i][j] *= corrected_scale;
01076                         }
01077                 }
01078         }
01079 }
01080 
01081 float Transform::get_scale() const {
01082         float determinant = get_determinant();
01083         if (determinant < 0 ) determinant *= -1;
01084 
01085         float scale = std::pow(determinant,1.0f/3.0f);
01086         int int_scale = static_cast<int>(scale);
01087         float scale_residual = scale-static_cast<float>(int_scale);
01088         if  ( scale_residual < ERR_LIMIT ) { scale = static_cast<float>(int_scale); };
01089 
01090         Util::apply_precision(scale, ERR_LIMIT);
01091 
01092         return scale;
01093 }
01094 
01095 void Transform::scale(const float& scale)
01096 {
01097         float determinant = get_determinant();
01098         if (determinant < 0) determinant *= -1.0f;
01099         float newscale = std::pow(determinant,1.0f/3.0f) + scale;
01100         //IF scale is < 0.01 The det can sometime be computed to be 0 (The precission is a bit Rubbish here)
01101         if(newscale > 0.01) 
01102         {
01103                 //cout << newscale << endl;
01104                 set_scale(newscale); // If scale ~ 0 things blowup, so we need a little fudge factor
01105         }
01106 }
01107 
01108 void print_matrix(gsl_matrix* M, unsigned int r, unsigned int c, const string& message ) {
01109         cout << "Message is " << message << endl;
01110         for ( unsigned int i = 0; i < r; ++i )
01111         {
01112                 for ( unsigned int j = 0; j < c; ++j )
01113                 {
01114                         cout << gsl_matrix_get(M,i,j) << " ";
01115                 }
01116                 cout << endl;
01117         }
01118 }
01119 
01120 void Transform::orthogonalize()
01121 {
01122         float scale;
01123         bool x_mirror;
01124         get_scale_and_mirror(scale,x_mirror);
01125         if (scale == 0) throw UnexpectedBehaviorException("The determinant of the Transform is 0. This is unexpected.");
01126         double inv_scale = 1.0/static_cast<double>(scale);
01127         double mirror_scale = (x_mirror == true ? -1.0:1.0);
01128 
01129         gsl_matrix * R = gsl_matrix_calloc(3,3);
01130         for ( unsigned int i = 0; i < 3; ++i )
01131         {
01132                 for ( unsigned int j = 0; j < 3; ++j )
01133                 {
01134                         if (i == 0 && mirror_scale != 1.0 ) {
01135                                 gsl_matrix_set( R, i, j, static_cast<double>(matrix[i][j])*mirror_scale*inv_scale );
01136                         }
01137                         else {
01138                                 gsl_matrix_set( R, i, j, static_cast<double>(matrix[i][j])*inv_scale );
01139                         }
01140                 }
01141         }
01142 
01143         gsl_matrix * V = gsl_matrix_calloc(3,3);
01144         gsl_vector * S = gsl_vector_calloc(3);
01145         gsl_vector * work = gsl_vector_calloc(3);
01146         gsl_linalg_SV_decomp (R, V, S, work); // Now R is U of the SVD R = USV^T
01147 
01148         gsl_matrix * Soln = gsl_matrix_calloc(3,3);
01149         gsl_blas_dgemm (CblasNoTrans, CblasTrans, 1.0, R, V, 0.0, Soln);
01150 
01151         for ( unsigned int i = 0; i < 3; ++i )
01152         {
01153                 for ( unsigned int j = 0; j < 3; ++j )
01154                 {
01155                         matrix[i][j] = static_cast<float>( gsl_matrix_get(Soln,i,j) );
01156                 }
01157         }
01158 
01159         // Apply scale if it existed previously
01160         if (scale != 1.0f) {
01161                 for(int i=0; i<3; ++i) {
01162                         for(int j=0; j<3; ++j) {
01163                                 matrix[i][j] *= scale;
01164                         }
01165                 }
01166         }
01167 
01168         // Apply post x mirroring if it was applied previouslys
01169         if ( x_mirror ) {
01170                 for(int j=0; j<3; ++j) {
01171                         matrix[0][j] *= -1.0f;
01172                 }
01173         }
01174 
01175         gsl_matrix_free(V); gsl_matrix_free(R); gsl_matrix_free(Soln);
01176         gsl_vector_free(S); gsl_vector_free(work);
01177 }
01178 
01179 void Transform::set_mirror(const bool x_mirror ) {
01180 
01181         bool old_x_mirror = get_mirror();
01182         if (old_x_mirror == x_mirror) return; // The user is setting the same value
01183         else {
01184                 // Toggle the mirroring operation
01185                 for (int j = 0; j < 4; ++j ) {
01186                         matrix[0][j] *= -1;
01187                 }
01188         }
01189 }
01190 
01191 bool Transform::get_mirror() const {
01192         float determinant = get_determinant();
01193 
01194         bool x_mirror = false;
01195         if ( determinant < 0 ) x_mirror = true;
01196 
01197         return x_mirror;
01198 
01199 }
01200 
01201 void Transform::get_scale_and_mirror(float& scale, bool& x_mirror) const {
01202 
01203         float determinant = get_determinant();
01204         x_mirror = false;
01205         if ( determinant < 0 ) {
01206                 x_mirror = true;
01207                 determinant *= -1;
01208         }
01209         if (determinant != 1 ) {
01210                 scale = std::pow(determinant,1.0f/3.0f);
01211                 int int_scale = static_cast<int>(scale);
01212                 float scale_residual = scale-static_cast<float>(int_scale);
01213                 if  ( scale_residual < ERR_LIMIT ) { scale = static_cast<float>(int_scale); };
01214         }
01215         else scale = 1;
01216 
01217         Util::apply_precision(scale,ERR_LIMIT);
01218 }
01219 
01220 float Transform::get_determinant() const
01221 {
01222         float det;
01223         double det2;
01224         det2  = matrix[0][0]*((double)matrix[1][1]*matrix[2][2]-(double)matrix[2][1]*matrix[1][2]);
01225         det2 -= matrix[0][1]*((double)matrix[1][0]*matrix[2][2]-(double)matrix[2][0]*matrix[1][2]);
01226         det2 += matrix[0][2]*((double)matrix[1][0]*matrix[2][1]-(double)matrix[2][0]*matrix[1][1]);
01227 
01228         det = (float)det2;
01229         Util::apply_precision(det,ERR_LIMIT);
01230 
01231         return det;
01232 }
01233 
01234 void Transform::invert() {
01235 
01236         double m00 = matrix[0][0]; double m01=matrix[0][1]; double m02=matrix[0][2];
01237         double m10 = matrix[1][0]; double m11=matrix[1][1]; double m12=matrix[1][2];
01238         double m20 = matrix[2][0]; double m21=matrix[2][1]; double m22=matrix[2][2];
01239         double v0  = matrix[0][3]; double v1 =matrix[1][3]; double v2 =matrix[2][3];
01240 
01241         double cof00 = m11*m22-m12*m21;
01242         double cof11 = m22*m00-m20*m02;
01243         double cof22 = m00*m11-m01*m10;
01244         double cof01 = m10*m22-m20*m12;
01245         double cof02 = m10*m21-m20*m11;
01246         double cof12 = m00*m21-m01*m20;
01247         double cof10 = m01*m22-m02*m21;
01248         double cof20 = m01*m12-m02*m11;
01249         double cof21 = m00*m12-m10*m02;
01250 
01251         double det = m00* cof00 + m02* cof02 -m01*cof01;
01252 
01253         matrix[0][0] =   (float)(cof00/det);
01254         matrix[0][1] = - (float)(cof10/det);
01255         matrix[0][2] =   (float)(cof20/det);
01256         matrix[1][0] = - (float)(cof01/det);
01257         matrix[1][1] =   (float)(cof11/det);
01258         matrix[1][2] = - (float)(cof21/det);
01259         matrix[2][0] =   (float)(cof02/det);
01260         matrix[2][1] = - (float)(cof12/det);
01261         matrix[2][2] =   (float)(cof22/det);
01262 
01263         matrix[0][3] =  (float)((- cof00*v0 + cof10*v1 - cof20*v2)/det);
01264         matrix[1][3] =  (float)((  cof01*v0 - cof11*v1 + cof21*v2)/det);
01265         matrix[2][3] =  (float)((- cof02*v0 + cof12*v1 - cof22*v2)/det);
01266 }
01267 
01268 Transform Transform::inverse() const {
01269         Transform t(*this);
01270         t.invert();
01271         return t;
01272 }
01273 
01274 void Transform::transpose_inplace() {
01275         float tempij;
01276         for (int i = 0; i < 3; i++) {
01277                 for (int j = 0; j < i; j++) {
01278                         if (i != j) {
01279                                 tempij= matrix[i][j];
01280                                 matrix[i][j] = matrix[j][i];
01281                                 matrix[j][i] = tempij;
01282                         }
01283                 }
01284         }
01285 }
01286 
01287 Transform Transform::transpose() const {
01288         Transform t(*this);
01289         t.transpose_inplace();
01290         return t;
01291 }
01292 
01293 
01294 Transform EMAN::operator*(const Transform & M2, const Transform & M1)     // YYY
01295 {
01296         Transform result;
01297         for (int i=0; i<3; i++) {
01298                 for (int j=0; j<4; j++) {
01299                         result[i][j] = M2[i][0] * M1[0][j] +  M2[i][1] * M1[1][j] + M2[i][2] * M1[2][j];
01300                 }
01301                 result[i][3] += M2[i][3];
01302         }
01303 
01304         return result;
01305 }
01306 
01307 void Transform::assert_valid_2d() const {
01308         int rotation_error = 0;
01309         int translation_error = 0;
01310         if (fabs(matrix[2][0]) > ERR_LIMIT) rotation_error++;
01311         if (fabs(matrix[2][1]) > ERR_LIMIT) rotation_error++;
01312         if (fabs(matrix[2][3]) > ERR_LIMIT) translation_error++;
01313         if (fabs(matrix[0][2]) > ERR_LIMIT) rotation_error++;
01314         if (fabs(matrix[1][2]) > ERR_LIMIT) rotation_error++;
01315 //      if (fabs(matrix[2][2]-1.0) >ERR_LIMIT) rotation_error++; 
01316         if ( translation_error && rotation_error ) {
01317                 throw UnexpectedBehaviorException("Error, the internal matrix contains 3D rotations and 3D translations. This object can not be considered 2D");
01318         } else if ( translation_error ) {
01319                 throw UnexpectedBehaviorException("Error, the internal matrix contains a non zero z component for a 3D translation. This object can not be considered 2D");
01320         }
01321         else if ( rotation_error ) {
01322                 throw UnexpectedBehaviorException("Error, the internal matrix contains 3D rotations and this object can not be considered 2D");
01323         }
01324 
01325 }
01326 
01327 
01328 
01329 Transform Transform::get_sym(const string & sym_name, int n) const
01330 {
01331         Symmetry3D* sym = Factory<Symmetry3D>::get(sym_name);
01332         Transform ret;
01333         ret = (*this) * sym->get_sym(n);
01334         delete sym;
01335         return ret;
01336 }
01337 
01338 vector<Transform > Transform::get_sym_proj(const string & sym_name) const
01339 {
01340         vector<Transform> ret;
01341         Transform t;
01342         Symmetry3D* sym = Factory<Symmetry3D>::get(sym_name);
01343         int nsym = sym->get_nsym();
01344         int n = nsym;
01345         
01346         if ((sym_name[0] == 'c' || sym_name[0] == 'd' ) &&  fabs(matrix[2][2]) < 1.e-6){
01347                 
01348                 Dict d1,d2;                             
01349                 d2["theta"] = (double)90.0;
01350                 d2["psi"] = (double)0.0;
01351                 d2["phi"] = (double)0.0;
01352                 d2["type"] = "spider";
01353                 d1 = this->get_rotation("spider");
01354                 
01355                 if (sym_name[0] == 'c') {
01356                         if( nsym%2 == 0)        n = nsym/2;
01357                         
01358                         for (int k=0;k<n;k++) {                         
01359                                 d2["phi"] = (double)d1["phi"] + k*double(360.0)/ nsym;
01360                                 d2["psi"] = d1["psi"];
01361                                 t.set_rotation(d2);
01362                                 ret.push_back( t );
01363                         }
01364                                 
01365                 }
01366                 else {
01367                         nsym = nsym/2;
01368                         
01369                         if (nsym%2 == 0) {
01370                                 n = nsym;
01371                                 float cos_phi = cos( EMConsts::deg2rad*360.0/2/nsym );
01372                         
01373                                 for (int k=0;k<n;k++){
01374                                         
01375                                         if(k%2==0)      {
01376                                         
01377                                                 d2["phi"] = (double)d1["phi"] + k/2*double(360.0)/ nsym;
01378                                                 d2["psi"] = d1["psi"];
01379                                                 t.set_rotation(d2);
01380                                                 ret.push_back( t );     
01381                                         }
01382                                         else    {
01383                                                         
01384                                                 if( ( fabs(1.0-matrix[2][0])>1.0e-6 )&& fabs( matrix[2][0]-cos_phi)>1.0e-6  ){
01385                                                         //cout<<"jumped into"<<endl;
01386                                                         d2["phi"] = k/2*double(360.0)/ nsym +180 - (double)d1["phi"];
01387                                                         d2["psi"] = (double)d1["psi"] + 180;
01388                                                         t.set_rotation(d2);
01389                                                         ret.push_back( t );
01390                                                 }
01391                                         }
01392                                 
01393                                 }
01394                         }
01395                         
01396                         
01397                         
01398                         else    {
01399                                 n = nsym*2;
01400                                 float cos_phi = cos( EMConsts::deg2rad*360.0/4/nsym );
01401                                 for (int k=0;k<n;k++){
01402                                         
01403                                         if(k%4==0)      {
01404                                         
01405                                                 d2["phi"] = (double)d1["phi"] + k/4*360.0/ nsym;
01406                                                 d2["psi"] = (double)d1["psi"];
01407                                                 t.set_rotation(d2);
01408                                                 ret.push_back( t );     
01409                                         }
01410                                         else if( k%4 ==1)       {
01411                                                 if( ( fabs(1.0-matrix[2][0])>1.0e-6 )&& fabs( matrix[2][0]-cos_phi)>1.0e-6  ){
01412                                                 
01413                                                         d2["phi"] = k/4*360.0/nsym + 360.0/2/nsym+180 - (double)d1["phi"];
01414                                                         d2["psi"] = (double)d1["psi"] + 180;
01415                                                         t.set_rotation(d2);
01416                                                         ret.push_back( t );
01417                                                 }
01418                                 
01419                                         }
01420                                         
01421                                         else if( k%4 ==2)       {
01422                                         
01423                                                 d2["phi"] =  k/4*360.0/ nsym+360.0/2/nsym+180 + (double)d1["phi"];
01424                                                 d2["psi"] = (double)d1["psi"];
01425                                                 t.set_rotation(d2);
01426                                                 ret.push_back( t );
01427                                 
01428                                         }
01429                                         
01430                                         else if( k%4 ==3)       {
01431                                                 if( ( fabs(1.0-matrix[2][0])>1.0e-6 )&& fabs( matrix[2][0]-cos_phi)>1.0e-6  ) {
01432                                                         d2["phi"] = k/4*360.0/nsym+ 2.0*360.0/2/nsym - (double)d1["phi"];
01433                                                         d2["psi"] = (double)d1["psi"] + 180;
01434                                                         t.set_rotation(d2);
01435                                                         ret.push_back( t );
01436                                                 }
01437                                         }
01438                                 
01439                                 }
01440                         }
01441                         
01442                 }
01443                 
01444         }
01445         else {
01446                 for (int k=0;k<nsym;k++) {
01447                         t =  sym->get_sym(k);
01448                         ret.push_back( (*this) * t );
01449                 }
01450                 
01451         }
01452         delete sym;
01453         return ret;
01454 }
01455 
01456 
01457 int Transform::get_nsym(const string & sym_name)
01458 {
01459         Symmetry3D* sym = Factory<Symmetry3D>::get(sym_name);
01460         int nsym = sym->get_nsym();
01461         delete sym;
01462         return nsym;
01463 }
01464 
01465 //
01466 //Transform3D::Transform3D()  //    C1
01467 //{
01468 //      init();
01469 //}
01470 //
01471 //Transform3D::Transform3D( const Transform3D& rhs )
01472 //{
01473 //    for( int i=0; i < 4; ++i )
01474 //    {
01475 //        for( int j=0; j < 4; ++j )
01476 //      {
01477 //          matrix[i][j] = rhs.matrix[i][j];
01478 //      }
01479 //    }
01480 //}
01481 //
01483 //Transform3D::Transform3D(const float& az, const float& alt, const float& phi)
01484 //{
01485 //      init();
01486 //      set_rotation(az,alt,phi);
01487 //}
01488 //
01489 //
01491 //Transform3D::Transform3D(const float& az, const float& alt, const float& phi, const Vec3f& posttrans )
01492 //{
01493 //      init(); // This is called in set_rotation
01494 //      set_rotation(az,alt,phi);
01495 //      set_posttrans(posttrans);
01496 //}
01497 //
01498 //Transform3D::Transform3D(const float& m11, const float& m12, const float& m13,
01499 //                                               const float& m21, const float& m22, const float& m23,
01500 //                                               const float& m31, const float& m32, const float& m33)
01501 //{
01502 //      init();
01503 //      set_rotation(m11,m12,m13,m21,m22,m23,m31,m32,m33);
01504 //}
01505 //
01507 //Transform3D::Transform3D(EulerType euler_type, const float& a1, const float& a2, const float& a3)
01508 //{
01509 //      init();
01510 //      set_rotation(euler_type,a1,a2,a3);
01511 //}
01512 //
01513 //Transform3D::Transform3D(EulerType euler_type, const float& a1, const float& a2, const float& a3, const float& a4)
01514 //{
01515 //      init();
01516 //      set_rotation(euler_type,a1,a2,a3,a4);
01517 //}
01518 //
01519 //
01521 //Transform3D::Transform3D(EulerType euler_type, const Dict& rotation)  //YYY
01522 //{
01523 //      init();
01524 //      set_rotation(euler_type,rotation);
01525 //}
01526 //
01527 //
01529 //
01530 //Transform3D::Transform3D(  const Vec3f& pretrans,  const float& az, const float& alt, const float& phi, const Vec3f& posttrans )  //YYY  by default EMAN
01531 //{
01532 //      init();
01533 //      set_pretrans(pretrans);
01534 //      set_rotation(az,alt,phi);
01535 //      set_posttrans(posttrans);
01536 //}
01537 //
01538 //
01539 //
01540 //
01541 //Transform3D::~Transform3D()
01542 //{
01543 //}
01544 //
01545 //
01546 //
01547 //void Transform3D::to_identity()
01548 //{
01552 //
01553 //      for(int i=0; i<4; ++i) {
01554 //              for(int j=0; j<4; ++j) {
01555 //                      if(i==j) {
01556 //                              matrix[i][j] = 1;
01557 //                      }
01558 //                      else {
01559 //                              matrix[i][j] = 0;
01560 //                      }
01561 //              }
01562 //      }
01563 //      post_x_mirror = false;
01564 //      set_center(Vec3f(0,0,0));
01565 //}
01566 //
01567 //
01568 //
01569 //bool Transform3D::is_identity()  // YYY
01570 //{
01571 //      for (int i=0; i<4; i++) {
01572 //              for (int j=0; j<4; j++) {
01573 //                      if (i==j && matrix[i][j]!=1.0) return 0;
01574 //                      if (i!=j && matrix[i][j]!=0.0) return 0;
01575 //              }
01576 //      }
01577 //      return 1;
01578 //}
01579 //
01580 //
01581 //void Transform3D::set_center(const Vec3f & center) //YYN
01582 //{
01583 //      set_pretrans( Vec3f(0,0,0)-center);
01584 //      for (int i = 0; i < 3; i++) {
01585 //              matrix[i][3]=center[i];
01586 //      }
01587 //}
01588 //
01591 //void Transform3D::init()  // M1
01592 //{
01593 //      to_identity();
01594 //}
01595 //
01597 //
01598 //void Transform3D::set_pretrans(const float& dx, const float& dy, const float& dz) // YYY
01599 //{    set_pretrans( Vec3f(dx,dy,dz)); }
01600 //
01601 //
01602 //void Transform3D::set_pretrans(const float& dx, const float& dy) // YYY
01603 //{    set_pretrans( Vec3f(dx,dy,0)); }
01604 //
01605 //void Transform3D::set_pretrans(const Vec2f& pretrans) // YYY
01606 //{    set_pretrans( Vec3f(pretrans[0],pretrans[1],0)); }
01607 //
01608 //void Transform3D::set_pretrans(const Vec3f & preT)  // flag=1 means keep the old value of total trans
01609 //{
01610 //              int flag=0;
01611 //
01614 //    if (flag==0){
01615 //              matrix[0][3] = matrix[3][0] + matrix[0][0]*preT[0] + matrix[0][1]*preT[1] + matrix[0][2]*preT[2]  ;
01616 //              matrix[1][3] = matrix[3][1] + matrix[1][0]*preT[0] + matrix[1][1]*preT[1] + matrix[1][2]*preT[2]  ;
01617 //              matrix[2][3] = matrix[3][2] + matrix[2][0]*preT[0] + matrix[2][1]*preT[1] + matrix[2][2]*preT[2]  ;
01618 //      }
01620 //    if (flag==1){
01621 //              matrix[3][0] = matrix[0][3] - (matrix[0][0]*preT[0] + matrix[0][1]*preT[1] + matrix[0][2]*preT[2])  ;
01622 //              matrix[3][1] = matrix[1][3] - (matrix[1][0]*preT[0] + matrix[1][1]*preT[1] + matrix[1][2]*preT[2])  ;
01623 //              matrix[3][2] = matrix[2][3] - (matrix[2][0]*preT[0] + matrix[2][1]*preT[1] + matrix[2][2]*preT[2])  ;
01624 //      }
01625 //}
01626 //
01627 //
01628 //void Transform3D::set_posttrans(const float& dx, const float& dy, const float& dz) // YYY
01629 //{    set_posttrans( Vec3f(dx,dy,dz)); }
01630 //
01631 //
01632 //void Transform3D::set_posttrans(const float& dx, const float& dy) // YYY
01633 //{    set_posttrans( Vec3f(dx,dy,0)); }
01634 //
01635 //void Transform3D::set_posttrans(const Vec2f& posttrans) // YYY
01636 //{    set_pretrans( Vec3f(posttrans[0],posttrans[1],0)); }
01637 //
01638 //void Transform3D::set_posttrans(const Vec3f & posttrans) // flag=1 means keep the old value of total trans
01639 //{
01640 //      int flag=0;
01641 //    Vec3f preT   = get_pretrans(0) ;
01642 //      for (int i = 0; i < 3; i++) {
01643 //              matrix[3][i] = posttrans[i];
01644 //      }
01647 //      if (flag==0) {
01648 //              matrix[0][3] = matrix[3][0] + matrix[0][0]*preT[0] + matrix[0][1]*preT[1] + matrix[0][2]*preT[2]  ;
01649 //              matrix[1][3] = matrix[3][1] + matrix[1][0]*preT[0] + matrix[1][1]*preT[1] + matrix[1][2]*preT[2]  ;
01650 //              matrix[2][3] = matrix[3][2] + matrix[2][0]*preT[0] + matrix[2][1]*preT[1] + matrix[2][2]*preT[2]  ;
01651 //      }
01653 //      if (flag==1) { // Don't do anything
01654 //      }
01655 //}
01656 //
01657 //
01658 //
01659 //
01660 //void Transform3D::apply_scale(const float& scale)    // YYY
01661 //{
01662 //      for (int i = 0; i < 3; i++) {
01663 //              for (int j = 0; j < 4; j++) {
01664 //                      matrix[i][j] *= scale;
01665 //              }
01666 //      }
01667 //      for (int j = 0; j < 3; j++) {
01668 //              matrix[3][j] *= scale;
01669 //      }
01670 //}
01671 //
01672 //void Transform3D::orthogonalize()  // YYY
01673 //{
01674 //      //EulerType EMAN;
01675 //      float scale = get_scale() ;
01676 //      float inverseScale= 1/scale ;
01677 //      apply_scale(inverseScale);
01680 //}
01681 //
01682 //
01683 //void Transform3D::transpose()  // YYY
01684 //{
01685 //      float tempij;
01686 //      for (int i = 0; i < 3; i++) {
01687 //              for (int j = 0; j < i; j++) {
01688 //                      tempij= matrix[i][j];
01689 //                      matrix[i][j] = matrix[j][i];
01690 //                      matrix[j][i] = tempij;
01691 //              }
01692 //      }
01693 //}
01694 //
01695 //void Transform3D::set_scale(const float& scale)    // YYY
01696 //{
01697 //      float OldScale= get_scale();
01698 //      float Scale2Apply = scale/OldScale;
01699 //      apply_scale(Scale2Apply);
01700 //}
01701 //
01702 //float Transform3D::get_mag() const //
01703 //{
01704 //      EulerType eulertype= SPIN ;
01705 //      Dict AA= get_rotation(eulertype);
01706 //      return AA["Omega"];
01707 //}
01708 //
01709 //Vec3f Transform3D::get_finger() const //
01710 //{
01711 //      EulerType eulertype= SPIN ;
01712 //      Dict AA= get_rotation(eulertype);
01713 //      return Vec3f(AA["n1"],AA["n2"],AA["n3"]);
01714 //}
01715 //
01716 //Vec3f Transform3D::get_posttrans(int flag) const    //
01717 //{
01718 //      if (flag==0){
01719 //              return Vec3f(matrix[3][0], matrix[3][1], matrix[3][2]);
01720 //      }
01721 //      // otherwise as if all the translation was post
01722 //      return Vec3f(matrix[0][3], matrix[1][3], matrix[2][3]);
01723 //}
01724 //
01725 //Vec3f Transform3D::get_total_posttrans() const {
01726 //      return get_posttrans(1);
01727 //}
01728 //
01729 //Vec3f Transform3D::get_total_pretrans() const {
01730 //      return get_pretrans(1);
01731 //}
01732 //
01733 //
01734 //Vec3f Transform3D::get_pretrans(int flag) const    // Fix Me
01735 //{
01737 //
01738 //      Vec3f pretrans;
01739 //      Vec3f posttrans(matrix[3][0], matrix[3][1], matrix[3][2]);
01740 //      Vec3f tottrans(matrix[0][3], matrix[1][3], matrix[2][3]);
01741 //      Vec3f totminuspost;
01742 //
01743 //      totminuspost = tottrans;
01744 //      if (flag==0) {
01745 //              totminuspost = tottrans-posttrans;
01746 //      }
01747 //
01748 //      Transform3D Rinv = inverse();
01749 //      for (int i=0; i<3; i++) {
01750 //                float ptnow=0;
01751 //              for (int j=0; j<3; j++) {
01752 //                      ptnow +=   Rinv.matrix[i][j]* totminuspost[j] ;
01753 //              }
01754 //              pretrans.set_value_at(i,ptnow) ;  //
01755 //      }
01756 //      return pretrans;
01757 //}
01758 //
01759 //
01760 // Vec3f Transform3D::get_center() const  // YYY
01761 // {
01762 //      return Vec3f();
01763 // }
01764 //
01765 //
01766 //
01767 //Vec3f Transform3D::get_matrix3_col(int i) const     // YYY
01768 //{
01769 //      return Vec3f(matrix[0][i], matrix[1][i], matrix[2][i]);
01770 //}
01771 //
01772 //
01773 //Vec3f Transform3D::get_matrix3_row(int i) const     // YYY
01774 //{
01775 //      return Vec3f(matrix[i][0], matrix[i][1], matrix[i][2]);
01776 //}
01777 //
01778 //Vec3f Transform3D::transform(const Vec3f & v3f) const     // YYY
01779 //{
01781 //      float x = matrix[0][0] * v3f[0] + matrix[0][1] * v3f[1] + matrix[0][2] * v3f[2] + matrix[0][3] ;
01782 //      float y = matrix[1][0] * v3f[0] + matrix[1][1] * v3f[1] + matrix[1][2] * v3f[2] + matrix[1][3] ;
01783 //      float z = matrix[2][0] * v3f[0] + matrix[2][1] * v3f[1] + matrix[2][2] * v3f[2] + matrix[2][3] ;
01784 //      return Vec3f(x, y, z);
01785 //}
01786 //
01787 //
01788 //Vec3f Transform3D::rotate(const Vec3f & v3f) const     // YYY
01789 //{
01791 //      float x = matrix[0][0] * v3f[0] + matrix[0][1] * v3f[1] + matrix[0][2] * v3f[2]  ;
01792 //      float y = matrix[1][0] * v3f[0] + matrix[1][1] * v3f[1] + matrix[1][2] * v3f[2]  ;
01793 //      float z = matrix[2][0] * v3f[0] + matrix[2][1] * v3f[1] + matrix[2][2] * v3f[2]  ;
01794 //      return Vec3f(x, y, z);
01795 //}
01796 //
01797 //
01798 //Transform3D EMAN::operator*(const Transform3D & M2, const Transform3D & M1)     // YYY
01799 //{
01802 //      Transform3D resultant;
01803 //      for (int i=0; i<3; i++) {
01804 //              for (int j=0; j<4; j++) {
01805 //                      resultant[i][j] = M2[i][0] * M1[0][j] +  M2[i][1] * M1[1][j] + M2[i][2] * M1[2][j];
01806 //              }
01807 //              resultant[i][3] += M2[i][3];  // add on the new translation (not included above)
01808 //      }
01809 //
01810 //      for (int j=0; j<3; j++) {
01811 //              resultant[3][j] = M2[3][j];
01812 //      }
01813 //
01814 //      return resultant; // This will have the post_trans of M2
01815 //}
01816 
01817 /*             Here starts the pure rotation stuff */
01818 
01888 //
01889 //void Transform3D::set_rotation(const float& az, const float& alt, const float& phi )
01890 //{
01891 //      EulerType euler_type=EMAN;
01892 //      Dict rot;
01893 //      rot["az"]  = az;
01894 //      rot["alt"] = alt;
01895 //      rot["phi"] = phi;
01896 //      set_rotation(euler_type, rot);
01897 //}
01898 //
01900 //void Transform3D::set_rotation(EulerType euler_type, const float& a1, const float& a2, const float& a3)
01901 //{
01902 //      init();
01903 //      Dict rot;
01904 //      switch(euler_type) {
01905 //              case EMAN:
01906 //                      rot["az"]  = a1;
01907 //                      rot["alt"] = a2;
01908 //                      rot["phi"] = a3;
01909 //                      break;
01910 //              case SPIDER:
01911 //                      rot["phi"]   = a1;
01912 //                      rot["theta"] = a2;
01913 //                      rot["psi"]   = a3;
01914 //                      break;
01915 //              case IMAGIC:
01916 //                      rot["alpha"]   = a1;
01917 //                      rot["beta"] = a2;
01918 //                      rot["gamma"]   = a3;
01919 //                      break;
01920 //              case MRC:
01921 //                      rot["phi"]   = a1;
01922 //                      rot["theta"] = a2;
01923 //                      rot["omega"]   = a3;
01924 //                      break;
01925 //              case XYZ:
01926 //                      rot["xtilt"]   = a1;
01927 //                      rot["ytilt"] = a2;
01928 //                      rot["ztilt"]   = a3;
01929 //                      break;
01930 //              default:
01931 //              throw InvalidValueException(euler_type, "cannot instantiate this Euler Type");
01932 //      }  // ends switch euler_type
01933 //      set_rotation(euler_type, rot);
01934 //}
01935 //
01937 //void Transform3D::set_rotation(EulerType euler_type, const float& a1, const float& a2, const float& a3, const float& a4)
01938 //{
01939 //      init();
01940 //      Dict rot;
01941 //      switch(euler_type) {
01942 //              case QUATERNION:
01943 //                      rot["e0"]  = a1;
01944 //                      rot["e1"] = a2;
01945 //                      rot["e2"] = a3;
01946 //                      rot["e3"] = a4;
01947 //                      break;
01948 //              case SGIROT:
01949 //                      rot["q"]  = a1;
01950 //                      rot["n1"] = a2;
01951 //                      rot["n2"] = a3;
01952 //                      rot["n3"] = a4;
01953 //              case SPIN:
01954 //                      rot["Omega"]  = a1;
01955 //                      rot["n1"] = a2;
01956 //                      rot["n2"] = a3;
01957 //                      rot["n3"] = a4;
01958 //                      break;
01959 //              default:
01960 //                      throw InvalidValueException(euler_type, "cannot instantiate this Euler Type");
01961 //      }  // ends switch euler_type
01962 //      set_rotation(euler_type, rot);
01963 //}
01964 //
01965 //void Transform3D::set_rotation(EulerType euler_type, const Dict& rotation)
01966 //{
01967 //      double e0  = 0; double e1=0; double e2=0; double e3=0;
01968 //      double Omega=0;
01969 //      double az  = 0;
01970 //      double alt = 0;
01971 //      double phi = 0;
01972 //      double cxtilt = 0;
01973 //      double sxtilt = 0;
01974 //      double cytilt = 0;
01975 //      double sytilt = 0;
01976 //      bool is_quaternion = 0;
01977 //      bool is_matrix = 0;
01978 //
01979 //      switch(euler_type) {
01980 //      case EMAN:
01981 //              az  = (double)rotation["az"] ;
01982 //              alt = (double)rotation["alt"]  ;
01983 //              phi = (double)rotation["phi"] ;
01984 //              break;
01985 //      case IMAGIC:
01986 //              az  = (double)rotation["alpha"] ;
01987 //              alt = (double)rotation["beta"]  ;
01988 //              phi = (double)rotation["gamma"] ;
01989 //              break;
01990 //
01991 //      case SPIDER:
01992 //              az =  (double)rotation["phi"]    + 90.0;
01993 //              alt = (double)rotation["theta"] ;
01994 //              phi = (double)rotation["psi"]    - 90.0;
01995 //              break;
01996 //
01997 //      case XYZ:
01998 //              cxtilt = cos( (M_PI/180.0f)*(double)rotation["xtilt"]);
01999 //              sxtilt = sin( (M_PI/180.0f)*(double)rotation["xtilt"]);
02000 //              cytilt = cos( (M_PI/180.0f)*(double)rotation["ytilt"]);
02001 //              sytilt = sin( (M_PI/180.0f)*(double)rotation["ytilt"]);
02002 //              az =  (180.0f/M_PI)*atan2(-cytilt*sxtilt,sytilt)   + 90.0f ;
02003 //              alt = (180.0f/M_PI)*acos(cytilt*cxtilt)  ;
02004 //              phi = (double)rotation["ztilt"] +(180.0f/M_PI)*atan2(sxtilt,cxtilt*sytilt)   - 90.0f ;
02005 //              break;
02006 //
02007 //      case MRC:
02008 //              az  = (double)rotation["phi"]   + 90.0f ;
02009 //              alt = (double)rotation["theta"] ;
02010 //              phi = (double)rotation["omega"] - 90.0f ;
02011 //              break;
02012 //
02013 //      case QUATERNION:
02014 //              is_quaternion = 1;
02015 //              e0 = (double)rotation["e0"];
02016 //              e1 = (double)rotation["e1"];
02017 //              e2 = (double)rotation["e2"];
02018 //              e3 = (double)rotation["e3"];
02019 //              break;
02020 //
02021 //      case SPIN:
02022 //              is_quaternion = 1;
02023 //              Omega = (double)rotation["Omega"];
02024 //              e0 = cos(Omega*M_PI/360.0f);
02025 //              e1 = sin(Omega*M_PI/360.0f)* (double)rotation["n1"];
02026 //              e2 = sin(Omega*M_PI/360.0f)* (double)rotation["n2"];
02027 //              e3 = sin(Omega*M_PI/360.0f)* (double)rotation["n3"];
02028 //              break;
02029 //
02030 //      case SGIROT:
02031 //              is_quaternion = 1;
02032 //              Omega = (double)rotation["q"]  ;
02033 //              e0 = cos(Omega*M_PI/360.0f);
02034 //              e1 = sin(Omega*M_PI/360.0f)* (double)rotation["n1"];
02035 //              e2 = sin(Omega*M_PI/360.0f)* (double)rotation["n2"];
02036 //              e3 = sin(Omega*M_PI/360.0f)* (double)rotation["n3"];
02037 //              break;
02038 //
02039 //      case MATRIX:
02040 //              is_matrix = 1;
02041 //              matrix[0][0] = (float)rotation["m11"]  ;
02042 //              matrix[0][1] = (float)rotation["m12"]  ;
02043 //              matrix[0][2] = (float)rotation["m13"]  ;
02044 //              matrix[1][0] = (float)rotation["m21"]  ;
02045 //              matrix[1][1] = (float)rotation["m22"]  ;
02046 //              matrix[1][2] = (float)rotation["m23"]  ;
02047 //              matrix[2][0] = (float)rotation["m31"]  ;
02048 //              matrix[2][1] = (float)rotation["m32"]  ;
02049 //              matrix[2][2] = (float)rotation["m33"]  ;
02050 //              break;
02051 //
02052 //      default:
02053 //              throw InvalidValueException(euler_type, "unknown Euler Type");
02054 //      }  // ends switch euler_type
02055 //
02056 //
02057 //      Vec3f postT  = get_posttrans( ) ;
02058 //      Vec3f preT   = get_pretrans( ) ;
02059 //
02060 //
02061 //      double azp  = fmod(az,360.0)*M_PI/180.0;
02062 //      double altp  = alt*M_PI/180.0;
02063 //      double phip = fmod(phi,360.0)*M_PI/180.0;
02064 //
02065 //      if (!is_quaternion && !is_matrix) {
02066 //              matrix[0][0] =  (float)(cos(phip)*cos(azp) - cos(altp)*sin(azp)*sin(phip));
02067 //              matrix[0][1] =  (float)(cos(phip)*sin(azp) + cos(altp)*cos(azp)*sin(phip));
02068 //              matrix[0][2] =  (float)(sin(altp)*sin(phip));
02069 //              matrix[1][0] =  (float)(-sin(phip)*cos(azp) - cos(altp)*sin(azp)*cos(phip));
02070 //              matrix[1][1] =  (float)(-sin(phip)*sin(azp) + cos(altp)*cos(azp)*cos(phip));
02071 //              matrix[1][2] =  (float)(sin(altp)*cos(phip));
02072 //              matrix[2][0] =  (float)(sin(altp)*sin(azp));
02073 //              matrix[2][1] =  (float)(-sin(altp)*cos(azp));
02074 //              matrix[2][2] =  (float)cos(altp);
02075 //      }
02076 //      if (is_quaternion){
02077 //              matrix[0][0] = (float)(e0 * e0 + e1 * e1 - e2 * e2 - e3 * e3);
02078 //              matrix[0][1] = (float)(2.0f * (e1 * e2 + e0 * e3));
02079 //              matrix[0][2] = (float)(2.0f * (e1 * e3 - e0 * e2));
02080 //              matrix[1][0] = (float)(2.0f * (e2 * e1 - e0 * e3));
02081 //              matrix[1][1] = (float)(e0 * e0 - e1 * e1 + e2 * e2 - e3 * e3);
02082 //              matrix[1][2] = (float)(2.0f * (e2 * e3 + e0 * e1));
02083 //              matrix[2][0] = (float)(2.0f * (e3 * e1 + e0 * e2));
02084 //              matrix[2][1] = (float)(2.0f * (e3 * e2 - e0 * e1));
02085 //              matrix[2][2] = (float)(e0 * e0 - e1 * e1 - e2 * e2 + e3 * e3);
02086 //              // keep in mind matrix[0][2] is M13 gives an e0 e2 piece, etc
02087 //      }
02088 //      // Now do post and pretrans: vfinal = vpost + R vpre;
02089 //
02090 //      matrix[0][3] = postT[0] + matrix[0][0]*preT[0] + matrix[0][1]*preT[1] + matrix[0][2]*preT[2]  ;
02091 //      matrix[1][3] = postT[1] + matrix[1][0]*preT[0] + matrix[1][1]*preT[1] + matrix[1][2]*preT[2]  ;
02092 //      matrix[2][3] = postT[2] + matrix[2][0]*preT[0] + matrix[2][1]*preT[1] + matrix[2][2]*preT[2]  ;
02093 //}
02094 //
02095 //
02096 //void Transform3D::set_rotation(const float& m11, const float& m12, const float& m13,
02097 //                                                         const float& m21, const float& m22, const float& m23,
02098 //                const float& m31, const float& m32, const float& m33)
02099 //{
02100 //      EulerType euler_type = MATRIX;
02101 //      Dict rot;
02102 //      rot["m11"]  = m11;
02103 //      rot["m12"]  = m12;
02104 //      rot["m13"]  = m13;
02105 //      rot["m21"]  = m21;
02106 //      rot["m22"]  = m22;
02107 //      rot["m23"]  = m23;
02108 //      rot["m31"]  = m31;
02109 //      rot["m32"]  = m32;
02110 //      rot["m33"]  = m33;
02111 //      set_rotation(euler_type, rot);  // Or should it be &rot ?
02112 //}
02113 //
02114 //void Transform3D::set_rotation(const Vec3f & eahat, const Vec3f & ebhat,
02115 //                                    const Vec3f & eAhat, const Vec3f & eBhat)
02116 //{// this rotation rotates unit vectors a,b into A,B;
02118 //      Vec3f eahatcp(eahat);
02119 //      Vec3f ebhatcp(ebhat);
02120 //      Vec3f eAhatcp(eAhat);
02121 //      Vec3f eBhatcp(eBhat);
02122 //
02123 //      eahatcp.normalize();
02124 //      ebhatcp.normalize();
02125 //      eAhatcp.normalize();
02126 //      eBhatcp.normalize();
02127 //
02128 //      Vec3f aMinusA(eahatcp);
02129 //      aMinusA  -= eAhatcp;
02130 //      Vec3f bMinusB(ebhatcp);
02131 //      bMinusB  -= eBhatcp;
02132 //
02133 //
02134 //      Vec3f  nhat;
02135 //      float aAlength = aMinusA.length();
02136 //      float bBlength = bMinusB.length();
02137 //      if (aAlength==0){
02138 //              nhat=eahatcp;
02139 //      }else if (bBlength==0){
02140 //              nhat=ebhatcp;
02141 //      }else{
02142 //              nhat= aMinusA.cross(bMinusB);
02143 //              nhat.normalize();
02144 //      }
02145 //
02147 //
02148 //      Vec3f neahat  = eahatcp.cross(nhat);
02149 //      Vec3f nebhat  = ebhatcp.cross(nhat);
02150 //      Vec3f neAhat  = eAhatcp.cross(nhat);
02151 //      Vec3f neBhat  = eBhatcp.cross(nhat);
02152 //
02153 //      double cosOmegaA = (neahat.dot(neAhat))  / (neahat.dot(neahat));
02155 //      double sinOmegaA = (neahat.dot(eAhatcp)) / (neahat.dot(neahat));
02157 //
02158 //      double OmegaA = atan2(sinOmegaA,cosOmegaA);
02160 //
02161 //      EulerType euler_type=SPIN;
02162 //      Dict rotation;
02163 //      rotation["n1"]= nhat[0];
02164 //      rotation["n2"]= nhat[1];
02165 //      rotation["n3"]= nhat[2];
02166 //      rotation["Omega"] = (float)(OmegaA*180.0/M_PI);
02167 //      set_rotation(euler_type,  rotation);
02168 //}
02169 //
02170 //
02171 //float Transform3D::get_scale() const     // YYY
02172 //{
02173 //      // Assumes uniform scaling, calculation uses Z only.
02174 //      float scale =0;
02175 //      for (int i=0; i<3; i++) {
02176 //              for (int j=0; j<3; j++) {
02177 //                      scale = scale + matrix[i][j]*matrix[i][j];
02178 //              }
02179 //      }
02180 //
02181 //      return sqrt(scale/3);
02182 //}
02183 //
02184 //
02185 //
02186 //Dict Transform3D::get_rotation(EulerType euler_type) const
02187 //{
02188 //      Dict result;
02189 //
02190 //      double max = 1 - ERR_LIMIT;
02191 //      double sca=get_scale();
02192 //      double cosalt=matrix[2][2]/sca;
02193 //
02194 //
02195 //      double az=0;
02196 //      double alt = 0;
02197 //      double phi=0;
02198 //      double phiS = 0; // like az     (but in SPIDER ZXZ)
02199 //      double psiS =0;  // like phi  (but in SPIDER ZYZ)
02200 //
02201 //
02203 //
02204 //      if (cosalt > max) {  // that is, alt close to 0
02205 //              alt = 0;
02206 //              az=0;
02207 //              phi = (double)EMConsts::rad2deg * atan2(matrix[0][1], matrix[0][0]);
02208 //      }
02209 //      else if (cosalt < -max) { // alt close to pi
02210 //              alt = 180;
02211 //              az=0;
02212 //              phi=360.0f-(double)EMConsts::rad2deg * atan2(matrix[0][1], matrix[0][0]);
02213 //      }
02214 //      else {
02215 //              alt = (double)EMConsts::rad2deg * acos(cosalt);
02216 //              az  = 360.0f+(double)EMConsts::rad2deg * atan2(matrix[2][0], -matrix[2][1]);
02217 //              phi = 360.0f+(double)EMConsts::rad2deg * atan2(matrix[0][2], matrix[1][2]);
02218 //      }
02219 //      az =fmod(az+180.0,360.0)-180.0;
02220 //      phi=fmod(phi+180.0,360.0)-180.0;
02221 //
02223 //      if (fabs(cosalt) > max) {  // that is, alt close to 0
02224 //              phiS=0;
02225 //              psiS = az+phi;
02226 //      }
02227 //      else {
02228 //              phiS = az   - 90.0;
02229 //              psiS = phi  + 90.0;
02230 //      }
02231 //      phiS = fmod((phiS   + 360.0 ), 360.0) ;
02232 //      psiS = fmod((psiS   + 360.0 ), 360.0) ;
02233 //
02235 //
02236 //      double nphi = (az-phi)/2.0;
02237 //    // The next is also e0
02238 //      double cosOover2 = (cos((az+phi)*M_PI/360) * cos(alt*M_PI/360)) ;
02239 //      double sinOover2 = sqrt(1 -cosOover2*cosOover2);
02240 //      double cosnTheta = sin((az+phi)*M_PI/360) * cos(alt*M_PI/360) / sqrt(1-cosOover2*cosOover2) ;
02241 //      double sinnTheta = sqrt(1-cosnTheta*cosnTheta);
02242 //      double n1 = sinnTheta*cos(nphi*M_PI/180);
02243 //      double n2 = sinnTheta*sin(nphi*M_PI/180);
02244 //      double n3 = cosnTheta;
02245 //      double xtilt = 0;
02246 //      double ytilt = 0;
02247 //      double ztilt = 0;
02248 //
02249 //
02250 //      if (cosOover2<0) {
02251 //              cosOover2*=-1; n1 *=-1; n2*=-1; n3*=-1;
02252 //      }
02253 //
02254 //
02255 //      switch (euler_type) {
02256 //      case EMAN:
02257 //              result["az"]  = az;
02258 //              result["alt"] = alt;
02259 //              result["phi"] = phi;
02260 //              break;
02261 //
02262 //      case IMAGIC:
02263 //              result["alpha"] = az;
02264 //              result["beta"] = alt;
02265 //              result["gamma"] = phi;
02266 //              break;
02267 //
02268 //      case SPIDER:
02269 //              result["phi"]   = phiS;  // The first Euler like az
02270 //              result["theta"] = alt;
02271 //              result["psi"]   = psiS;
02272 //              break;
02273 //
02274 //      case MRC:
02275 //              result["phi"]   = phiS;
02276 //              result["theta"] = alt;
02277 //              result["omega"] = psiS;
02278 //              break;
02279 //
02280 //      case XYZ:
02281 //              xtilt = atan2(-sin((M_PI/180.0f)*phiS)*sin((M_PI/180.0f)*alt),cos((M_PI/180.0f)*alt));
02282 //              ytilt = asin(  cos((M_PI/180.0f)*phiS)*sin((M_PI/180.0f)*alt));
02283 //              ztilt = psiS*M_PI/180.0f - atan2(sin(xtilt), cos(xtilt) *sin(ytilt));
02284 //
02285 //              xtilt=fmod(xtilt*180/M_PI+540.0,360.0) -180.0;
02286 //              ztilt=fmod(ztilt*180/M_PI+540.0,360.0) -180.0;
02287 //
02288 //              result["xtilt"]  = xtilt;
02289 //              result["ytilt"]  = ytilt*180/M_PI;
02290 //              result["ztilt"]  = ztilt;
02291 //              break;
02292 //
02293 //      case QUATERNION:
02294 //              result["e0"] = cosOover2 ;
02295 //              result["e1"] = sinOover2 * n1 ;
02296 //              result["e2"] = sinOover2 * n2;
02297 //              result["e3"] = sinOover2 * n3;
02298 //              break;
02299 //
02300 //      case SPIN:
02301 //              result["Omega"] =360.0f* acos(cosOover2)/ M_PI ;
02302 //              result["n1"] = n1;
02303 //              result["n2"] = n2;
02304 //              result["n3"] = n3;
02305 //              break;
02306 //
02307 //      case SGIROT:
02308 //              result["q"] = 360.0f*acos(cosOover2)/M_PI ;
02309 //              result["n1"] = n1;
02310 //              result["n2"] = n2;
02311 //              result["n3"] = n3;
02312 //              break;
02313 //
02314 //      case MATRIX:
02315 //              result["m11"] = matrix[0][0] ;
02316 //              result["m12"] = matrix[0][1] ;
02317 //              result["m13"] = matrix[0][2] ;
02318 //              result["m21"] = matrix[1][0] ;
02319 //              result["m22"] = matrix[1][1] ;
02320 //              result["m23"] = matrix[1][2] ;
02321 //              result["m31"] = matrix[2][0] ;
02322 //              result["m32"] = matrix[2][1] ;
02323 //              result["m33"] = matrix[2][2] ;
02324 //              break;
02325 //
02326 //      default:
02327 //              throw InvalidValueException(euler_type, "unknown Euler Type");
02328 //      }
02329 //
02330 //      return result;
02331 //}
02332 //
02333 //Transform3D Transform3D::inverseUsingAngs() const    //   YYN need to test it for sure
02334 //{
02335 //      // First Find the scale
02336 //      EulerType eE=EMAN;
02337 //
02338 //
02339 //      float scale   = get_scale();
02340 //      Vec3f preT   = get_pretrans( ) ;
02341 //      Vec3f postT   = get_posttrans( ) ;
02342 //      Dict angs     = get_rotation(eE);
02343 //      Dict invAngs  ;
02344 //
02345 //      invAngs["phi"]   = 180.0f - (float) angs["az"] ;
02346 //      invAngs["az"]    = 180.0f - (float) angs["phi"] ;
02347 //      invAngs["alt"]   = angs["alt"] ;
02348 //
02355 //
02356 //      float inverseScale= 1/scale ;
02357 //
02358 //      Transform3D invM;
02359 //
02360 //      invM.set_rotation(EMAN, invAngs);
02361 //      invM.apply_scale(inverseScale);
02362 //      invM.set_pretrans(-postT );
02363 //      invM.set_posttrans(-preT );
02364 //
02365 //
02366 //      return invM;
02367 //
02368 //}
02369 //
02370 //Transform3D Transform3D::inverse() const    //   YYN need to test it for sure
02371 //{
02372 //      // This assumes the matrix is 4 by 4 and the last row reads [0 0 0 1]
02373 //
02374 //      double m00 = matrix[0][0]; double m01=matrix[0][1]; double m02=matrix[0][2];
02375 //      double m10 = matrix[1][0]; double m11=matrix[1][1]; double m12=matrix[1][2];
02376 //      double m20 = matrix[2][0]; double m21=matrix[2][1]; double m22=matrix[2][2];
02377 //      double v0  = matrix[0][3]; double v1 =matrix[1][3]; double v2 =matrix[2][3];
02378 //
02379 //    double cof00 = m11*m22-m12*m21;
02380 //    double cof11 = m22*m00-m20*m02;
02381 //    double cof22 = m00*m11-m01*m10;
02382 //    double cof01 = m10*m22-m20*m12;
02383 //    double cof02 = m10*m21-m20*m11;
02384 //    double cof12 = m00*m21-m01*m20;
02385 //    double cof10 = m01*m22-m02*m21;
02386 //    double cof20 = m01*m12-m02*m11;
02387 //    double cof21 = m00*m12-m10*m02;
02388 //
02389 //    double Det = m00* cof00 + m02* cof02 -m01*cof01;
02390 //
02391 //    Transform3D invM;
02392 //
02393 //    invM.matrix[0][0] =  (float)(cof00/Det);
02394 //    invM.matrix[0][1] =  (float)(-cof10/Det);
02395 //    invM.matrix[0][2] =  (float)(cof20/Det);
02396 //    invM.matrix[1][0] =  (float)(-cof01/Det);
02397 //    invM.matrix[1][1] =  (float)(cof11/Det);
02398 //    invM.matrix[1][2] =  (float)(-cof21/Det);
02399 //    invM.matrix[2][0] =  (float)(cof02/Det);
02400 //    invM.matrix[2][1] =  (float)(-cof12/Det);
02401 //    invM.matrix[2][2] =  (float)(cof22/Det);
02402 //
02403 //    invM.matrix[0][3] =  (float)((-cof00*v0 + cof10*v1 - cof20*v2)/Det);
02404 //    invM.matrix[1][3] =  (float)(( cof01*v0 - cof11*v1 + cof21*v2)/Det);
02405 //    invM.matrix[2][3] =  (float)((-cof02*v0 + cof12*v1 - cof22*v2)/Det);
02406 //
02407 //      Vec3f postT   = get_posttrans( ) ;
02408 //      Vec3f invMpre   = - postT;
02409 //      Vec3f invMpost   ;
02410 //      for ( int i = 0; i < 3; i++) {
02411 //              invMpost[i] = invM.matrix[i][3];
02412 //              for ( int j = 0; j < 3; j++) {
02413 //                      invMpost[i] += - invM.matrix[i][j]*invMpre[j];
02414 //              }
02415 //              invM.matrix[3][i] = invMpost[i];
02416 //      }
02417 //
02418 //      return invM;
02419 //}
02420 //
02421 //
02422 //
02424 //
02425 //Transform3D Transform3D::get_sym(const string & symname, int n) const
02426 //{
02427 //      int nsym = get_nsym(symname);
02428 //
02431 //
02432 //      // see www.math.utah.edu/~alfeld/math/polyhedra/polyhedra.html for pictures
02433 //      // By default we will put largest symmetry along z-axis.
02434 //
02435 //      // Each Platonic Solid has 2E symmetry elements.
02436 //
02437 //
02438 //      // An icosahedron has   m=5, n=3, F=20 E=30=nF/2, V=12=nF/m,since vertices shared by 5 triangles;
02439 //      // It is composed of 20 triangles. E=3*20/2;
02440 //
02441 //
02442 //      // An dodecahedron has m=3, n=5   F=12 E=30  V=20
02443 //      // It is composed of 12 pentagons. E=5*12/2;   V= 5*12/3, since vertices shared by 3 pentagons;
02444 //
02445 //
02446 //
02447 //    // The ICOS symmetry group has the face along z-axis
02448 //
02449 //      float lvl0=0;                             //  there is one pentagon on top; five-fold along z
02450 //      float lvl1= 63.4349f; // that is atan(2)  // there are 5 pentagons with centers at this height (angle)
02451 //      float lvl2=116.5651f; //that is 180-lvl1  // there are 5 pentagons with centers at this height (angle)
02452 //      float lvl3=180.f;                           // there is one pentagon on the bottom
02453 //             // Notice that 63.439 is the angle between two faces of the dual object
02454 //
02455 //      static double ICOS[180] = { // This is with a pentagon normal to z
02456 //                0,lvl0,0,    0,lvl0,288,   0,lvl0,216,   0,lvl0,144,  0,lvl0,72,
02457 //                0,lvl1,36,   0,lvl1,324,   0,lvl1,252,   0,lvl1,180,  0,lvl1,108,
02458 //               72,lvl1,36,  72,lvl1,324,  72,lvl1,252,  72,lvl1,180,  72,lvl1,108,
02459 //              144,lvl1,36, 144,lvl1,324, 144,lvl1,252, 144,lvl1,180, 144,lvl1,108,
02460 //              216,lvl1,36, 216,lvl1,324, 216,lvl1,252, 216,lvl1,180, 216,lvl1,108,
02461 //              288,lvl1,36, 288,lvl1,324, 288,lvl1,252, 288,lvl1,180, 288,lvl1,108,
02462 //               36,lvl2,0,   36,lvl2,288,  36,lvl2,216,  36,lvl2,144,  36,lvl2,72,
02463 //              108,lvl2,0,  108,lvl2,288, 108,lvl2,216, 108,lvl2,144, 108,lvl2,72,
02464 //              180,lvl2,0,  180,lvl2,288, 180,lvl2,216, 180,lvl2,144, 180,lvl2,72,
02465 //              252,lvl2,0,  252,lvl2,288, 252,lvl2,216, 252,lvl2,144, 252,lvl2,72,
02466 //              324,lvl2,0,  324,lvl2,288, 324,lvl2,216, 324,lvl2,144, 324,lvl2,72,
02467 //                0,lvl3,0,    0,lvl3,288,   0,lvl3,216,   0,lvl3,144,   0,lvl3,72
02468 //      };
02469 //
02470 //
02471 //      // A cube has   m=3, n=4, F=6 E=12=nF/2, V=8=nF/m,since vertices shared by 3 squares;
02472 //      // It is composed of 6 squares.
02473 //
02474 //
02475 //      // An octahedron has   m=4, n=3, F=8 E=12=nF/2, V=6=nF/m,since vertices shared by 4 triangles;
02476 //      // It is composed of 8 triangles.
02477 //
02478 //    // We have placed the OCT symmetry group with a face along the z-axis
02479 //        lvl0=0;
02480 //      lvl1=90;
02481 //      lvl2=180;
02482 //
02483 //      static float OCT[72] = {// This is with a face of a cube along z
02484 //                    0,lvl0,0,   0,lvl0,90,    0,lvl0,180,    0,lvl0,270,
02485 //                    0,lvl1,0,   0,lvl1,90,    0,lvl1,180,    0,lvl1,270,
02486 //                   90,lvl1,0,  90,lvl1,90,   90,lvl1,180,   90,lvl1,270,
02487 //                  180,lvl1,0, 180,lvl1,90,  180,lvl1,180,  180,lvl1,270,
02488 //                  270,lvl1,0, 270,lvl1,90,  270,lvl1,180,  270,lvl1,270,
02489 //                    0,lvl2,0,   0,lvl2,90,    0,lvl2,180,    0,lvl2,270
02490 //      };
02491 //      // B^4=A^3=1;  BABA=1; implies   AA=BAB, ABA=B^3 , AB^2A = BBBABBB and
02492 //      //   20 words with at most a single A
02493 //    //   1 B BB BBB A  BA AB BBA BAB ABB BBBA BBAB BABB ABBB BBBAB BBABB BABBB
02494 //    //                        BBBABB BBABBB BBBABBB
02495 //     // also     ABBBA is distinct yields 4 more words
02496 //     //    ABBBA   BABBBA BBABBBA BBBABBBA
02497 //     // for a total of 24 words
02498 //     // Note A BBB A BBB A  reduces to BBABB
02499 //     //  and  B A BBB A is the same as A BBB A BBB etc.
02500 //
02501 //    // The TET symmetry group has a face along the z-axis
02502 //    // It has n=m=3; F=4, E=6=nF/2, V=4=nF/m
02503 //        lvl0=0;         // There is a face along z
02504 //      lvl1=109.4712f;  //  that is acos(-1/3)  // There  are 3 faces at this angle
02505 //
02506 //      static float TET[36] = {// This is with the face along z
02507 //            0,lvl0,0,   0,lvl0,120,    0,lvl0,240,
02508 //            0,lvl1,60,   0,lvl1,180,    0,lvl1,300,
02509 //          120,lvl1,60, 120,lvl1,180,  120,lvl1,300,
02510 //          240,lvl1,60, 240,lvl1,180,  240,lvl1,300
02511 //      };
02512 //      // B^3=A^3=1;  BABA=1; implies   A^2=BAB, ABA=B^2 , AB^2A = B^2AB^2 and
02513 //      //   12 words with at most a single A
02514 //    //   1 B BB  A  BA AB BBA BAB ABB BBAB BABB BBABB
02515 //    // at most one A is necessary
02516 //
02517 //      Transform3D ret;
02518 //      SymType type = get_sym_type(symname);
02519 //
02520 //      switch (type) {
02521 //      case CSYM:
02522 //              ret.set_rotation( n * 360.0f / nsym, 0, 0);
02523 //              break;
02524 //      case DSYM:
02525 //              if (n >= nsym / 2) {
02526 //                      ret.set_rotation((n - nsym/2) * 360.0f / (nsym / 2),180.0f, 0);
02527 //              }
02528 //              else {
02529 //                      ret.set_rotation( n * 360.0f / (nsym / 2),0, 0);
02530 //              }
02531 //              break;
02532 //      case ICOS_SYM:
02533 //              ret.set_rotation((float)ICOS[n * 3 ],
02534 //                               (float)ICOS[n * 3 + 1],
02535 //                               (float)ICOS[n * 3 + 2] );
02536 //              break;
02537 //      case OCT_SYM:
02538 //              ret.set_rotation((float)OCT[n * 3],
02539 //                               (float)OCT[n * 3 + 1],
02540 //                               (float)OCT[n * 3 + 2] );
02541 //              break;
02542 //      case TET_SYM:
02543 //              ret.set_rotation((float)TET[n * 3 ],
02544 //                               (float)TET[n * 3 + 1] ,
02545 //                               (float)TET[n * 3 + 2] );
02546 //              break;
02547 //      case ISYM:
02548 //              ret.set_rotation(0, 0, 0);
02549 //              break;
02550 //      default:
02551 //              throw InvalidValueException(type, symname);
02552 //      }
02553 //
02554 //      ret = (*this) * ret;
02555 //
02556 //      return ret;
02557 //}
02558 //
02559 //int Transform3D::get_nsym(const string & name)
02560 //{
02561 //      string symname = name;
02562 //
02563 //      for (size_t i = 0; i < name.size(); i++) {
02564 //              if (isalpha(name[i])) {
02565 //                      symname[i] = (char)tolower(name[i]);
02566 //              }
02567 //      }
02568 //
02569 //      SymType type = get_sym_type(symname);
02570 //      int nsym = 0;
02571 //
02572 //      switch (type) {
02573 //      case CSYM:
02574 //              nsym = atoi(symname.c_str() + 1);
02575 //              break;
02576 //      case DSYM:
02577 //              nsym = atoi(symname.c_str() + 1) * 2;
02578 //              break;
02579 //      case ICOS_SYM:
02580 //              nsym = 60;
02581 //              break;
02582 //      case OCT_SYM:
02583 //              nsym = 24;
02584 //              break;
02585 //      case TET_SYM:
02586 //              nsym = 12;
02587 //              break;
02588 //      case ISYM:
02589 //              nsym = 1;
02590 //              break;
02591 //      case UNKNOWN_SYM:
02592 //      default:
02593 //              throw InvalidValueException(type, name);
02594 //      }
02595 //      return nsym;
02596 //}
02597 //
02598 //
02599 //
02600 //Transform3D::SymType Transform3D::get_sym_type(const string & name)
02601 //{
02602 //      SymType t = UNKNOWN_SYM;
02603 //
02604 //      if (name[0] == 'c') {
02605 //              t = CSYM;
02606 //      }
02607 //      else if (name[0] == 'd') {
02608 //              t = DSYM;
02609 //      }
02610 //      else if (name == "icos") {
02611 //              t = ICOS_SYM;
02612 //      }
02613 //      else if (name == "oct") {
02614 //              t = OCT_SYM;
02615 //      }
02616 //      else if (name == "tet") {
02617 //              t = TET_SYM;
02618 //      }
02619 //      else if (name == "i" || name == "") {
02620 //              t = ISYM;
02621 //      }
02622 //      return t;
02623 //}
02624 //
02625 //vector<Transform3D*>
02626 //Transform3D::angles2tfvec(EulerType eulertype, const vector<float> ang) {
02627 //      int nangles = ang.size() / 3;
02628 //      vector<Transform3D*> tfvec;
02629 //      for (int i = 0; i < nangles; i++) {
02630 //              tfvec.push_back(new Transform3D(eulertype,ang[3*i],ang[3*i+1],ang[3*i+2]));
02631 //      }
02632 //      return tfvec;
02633 //}
02634 //
02635 
02636 
02637 
02638 /* vim: set ts=4 noet: */
02639 
02640 
02641 /*    Rotation stuff */
02642 

Generated on Thu Nov 17 12:43:03 2011 for EMAN2 by  doxygen 1.4.7