aligner.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 #include "emfft.h"
00036 #include "cmp.h"
00037 #include "aligner.h"
00038 #include "emdata.h"
00039 #include "processor.h"
00040 #include "util.h"
00041 #include "symmetry.h"
00042 #include <gsl/gsl_multimin.h>
00043 #include "plugins/aligner_template.h"
00044 
00045 #ifdef EMAN2_USING_CUDA
00046         #include <sparx/cuda/cuda_ccf.h>
00047 #endif
00048 
00049 #define EMAN2_ALIGNER_DEBUG 0
00050 
00051 using namespace EMAN;
00052 
00053 const string TranslationalAligner::NAME = "translational";
00054 const string RotationalAligner::NAME = "rotational";
00055 const string RotatePrecenterAligner::NAME = "rotate_precenter";
00056 const string RotateTranslateAligner::NAME = "rotate_translate";
00057 const string RotateTranslateBestAligner::NAME = "rotate_translate_best";
00058 const string RotateFlipAligner::NAME = "rotate_flip";
00059 const string RotateTranslateFlipAligner::NAME = "rotate_translate_flip";
00060 const string RTFExhaustiveAligner::NAME = "rtf_exhaustive";
00061 const string RTFSlowExhaustiveAligner::NAME = "rtf_slow_exhaustive";
00062 const string RefineAligner::NAME = "refine";
00063 const string Refine3DAligner::NAME = "refine.3d";
00064 const string RT3DGridAligner::NAME = "rt.3d.grid";
00065 const string RT3DSphereAligner::NAME = "rt.3d.sphere";
00066 #ifdef FFTW2
00067 const string FRM2DAligner::NAME = "FRM2D";
00068 #endif  //FFTW2
00069 
00070 template <> Factory < Aligner >::Factory()
00071 {
00072         force_add<TranslationalAligner>();
00073         force_add<RotationalAligner>();
00074         force_add<RotatePrecenterAligner>();
00075         force_add<RotateTranslateAligner>();
00076         force_add<RotateFlipAligner>();
00077         force_add<RotateTranslateFlipAligner>();
00078         force_add<RTFExhaustiveAligner>();
00079         force_add<RTFSlowExhaustiveAligner>();
00080         force_add<RefineAligner>();
00081         force_add<Refine3DAligner>();
00082         force_add<RT3DGridAligner>();
00083         force_add<RT3DSphereAligner>();
00084 #ifdef  FFTW2
00085         force_add<FRM2DAligner>();
00086 #endif  //FFTW2
00087 //      force_add<XYZAligner>();
00088 }
00089 
00090 // Note, the translational aligner assumes that the correlation image
00091 // generated by the calc_ccf function is centered on the bottom left corner
00092 // That is, if you did at calc_cff using identical images, the
00093 // peak would be at 0,0
00094 EMData *TranslationalAligner::align(EMData * this_img, EMData *to,
00095                                         const string&, const Dict&) const
00096 {
00097         if (!this_img) {
00098                 return 0;
00099         }
00100 
00101         if (to && !EMUtil::is_same_size(this_img, to))
00102                 throw ImageDimensionException("Images must be the same size to perform translational alignment");
00103 
00104         EMData *cf = 0;
00105         int nx = this_img->get_xsize();
00106         int ny = this_img->get_ysize();
00107         int nz = this_img->get_zsize();
00108 
00109         int masked = params.set_default("masked",0);
00110         int useflcf = params.set_default("useflcf",0);
00111         bool use_cpu = true;
00112 #ifdef EMAN2_USING_CUDA
00113         if (this_img->gpu_operation_preferred() ) {
00114 //              cout << "Translate on GPU" << endl;
00115                 use_cpu = false;
00116                 cf = this_img->calc_ccf_cuda(to,false,false);
00117         }
00118 #endif // EMAN2_USING_CUDA
00119         if (use_cpu) {
00120                 if (useflcf) cf = this_img->calc_flcf(to);
00121                 else cf = this_img->calc_ccf(to);
00122         }
00123 
00124         // This is too expensive
00125         if (masked) {
00126                 EMData *msk=this_img->process("threshold.notzero");
00127                 EMData *sqr=to->process("math.squared");
00128                 EMData *cfn=msk->calc_ccf(sqr);
00129                 cfn->process_inplace("math.sqrt");
00130                 float *d1=cf->get_data();
00131                 float *d2=cfn->get_data();
00132                 for (int i=0; i<nx*ny*nz; i++) {
00133                         if (d2[i]!=0) d1[i]/=d2[i];
00134                 }
00135                 cf->update();
00136                 delete msk;
00137                 delete sqr;
00138                 delete cfn;
00139         }
00140 
00141 //
00142 
00143         int maxshiftx = params.set_default("maxshift",-1);
00144         int maxshifty = params["maxshift"];
00145         int maxshiftz = params["maxshift"];
00146         int nozero = params["nozero"];
00147 
00148         if (maxshiftx <= 0) {
00149                 maxshiftx = nx / 4;
00150                 maxshifty = ny / 4;
00151                 maxshiftz = nz / 4;
00152         }
00153 
00154         if (maxshiftx > nx / 2 - 1) maxshiftx = nx / 2 - 1;
00155         if (maxshifty > ny / 2 - 1)     maxshifty = ny / 2 - 1;
00156         if (maxshiftz > nz / 2 - 1) maxshiftz = nz / 2 - 1;
00157 
00158         if (nx == 1) maxshiftx = 0; // This is justhere for completeness really... plus it saves errors
00159         if (ny == 1) maxshifty = 0;
00160         if (nz == 1) maxshiftz = 0;
00161 
00162         // If nozero the portion of the image in the center (and its 8-connected neighborhood) is zeroed
00163         if (nozero) {
00164                 cf->zero_corner_circulant(1);
00165         }
00166 
00167         IntPoint peak;
00168 #ifdef EMAN2_USING_CUDA
00169         if (!use_cpu) {
00170                 EMDataForCuda tmp = cf->get_data_struct_for_cuda();
00171                 int* p = calc_max_location_wrap_cuda(&tmp,maxshiftx, maxshifty, maxshiftz);
00172                 peak = IntPoint(p[0],p[1],p[2]);
00173                 free(p);
00174         }
00175 #endif // EMAN2_USING_CUDA
00176         if (use_cpu) {
00177                 peak = cf->calc_max_location_wrap(maxshiftx, maxshifty, maxshiftz);
00178         }
00179         Vec3f cur_trans = Vec3f ( (float)-peak[0], (float)-peak[1], (float)-peak[2]);
00180 
00181         if (!to) {
00182                 cur_trans /= 2.0f; // If aligning theimage to itself then only go half way -
00183                 int intonly = params.set_default("intonly",false);
00184                 if (intonly) {
00185                         cur_trans[0] = floor(cur_trans[0] + 0.5f);
00186                         cur_trans[1] = floor(cur_trans[1] + 0.5f);
00187                         cur_trans[2] = floor(cur_trans[2] + 0.5f);
00188                 }
00189         }
00190 
00191         if( cf ){
00192                 delete cf;
00193                 cf = 0;
00194         }
00195         Dict params("trans",static_cast< vector<int> >(cur_trans));
00196         cf=this_img->process("math.translate.int",params);
00197         Transform t;
00198         t.set_trans(cur_trans);
00199         if ( nz != 1 ) {
00200 //              Transform* t = get_set_align_attr("xform.align3d",cf,this_img);
00201 //              t->set_trans(cur_trans);
00202                 cf->set_attr("xform.align3d",&t);
00203         } else if ( ny != 1 ) {
00204                 //Transform* t = get_set_align_attr("xform.align2d",cf,this_img);
00205                 cur_trans[2] = 0; // just make sure of it
00206                 t.set_trans(cur_trans);
00207                 cf->set_attr("xform.align2d",&t);
00208         }
00209 
00210         return cf;
00211 }
00212 
00213 EMData * RotationalAligner::align_180_ambiguous(EMData * this_img, EMData * to, int rfp_mode) {
00214 
00215         // Make translationally invariant rotational footprints
00216         EMData* this_img_rfp, * to_rfp;
00217         if (rfp_mode == 0) {
00218                 this_img_rfp = this_img->make_rotational_footprint_e1();
00219                 to_rfp = to->make_rotational_footprint_e1();
00220         } else if (rfp_mode == 1) {
00221                 this_img_rfp = this_img->make_rotational_footprint();
00222                 to_rfp = to->make_rotational_footprint();
00223         } else if (rfp_mode == 2) {
00224                 this_img_rfp = this_img->make_rotational_footprint_cmc();
00225                 to_rfp = to->make_rotational_footprint_cmc();
00226         } else {
00227                 throw InvalidParameterException("rfp_mode must be 0,1 or 2");
00228         }
00229         int this_img_rfp_nx = this_img_rfp->get_xsize();
00230 
00231         // Do row-wise correlation, returning a sum.
00232         EMData *cf = this_img_rfp->calc_ccfx(to_rfp, 0, this_img->get_ysize());
00233 
00234         // Delete them, they're no longer needed
00235         delete this_img_rfp; this_img_rfp = 0;
00236         delete to_rfp; to_rfp = 0;
00237 
00238         // Now solve the rotational alignment by finding the max in the column sum
00239         float *data = cf->get_data();
00240         float peak = 0;
00241         int peak_index = 0;
00242         Util::find_max(data, this_img_rfp_nx, &peak, &peak_index);
00243 
00244         if( cf ) {
00245                 delete cf;
00246                 cf = 0;
00247         }
00248         float rot_angle = (float) (peak_index * 180.0f / this_img_rfp_nx);
00249 
00250         // Return the result
00251         Transform tmp(Dict("type","2d","alpha",rot_angle));
00252         cf=this_img->process("xform",Dict("transform",(Transform*)&tmp));
00253 //      Transform* t = get_set_align_attr("xform.align2d",cf,this_img);
00254 //      Dict d("type","2d","alpha",rot_angle);
00255 //      t->set_rotation(d);
00256         cf->set_attr("xform.align2d",&tmp);
00257         return cf;
00258 }
00259 
00260 EMData *RotationalAligner::align(EMData * this_img, EMData *to,
00261                         const string& cmp_name, const Dict& cmp_params) const
00262 {
00263         if (!to) throw InvalidParameterException("Can not rotational align - the image to align to is NULL");
00264 
00265         // Perform 180 ambiguous alignment
00266         int rfp_mode = params.set_default("rfp_mode",0);
00267         EMData* rot_aligned = RotationalAligner::align_180_ambiguous(this_img,to,rfp_mode);
00268         Transform * tmp = rot_aligned->get_attr("xform.align2d");
00269         Dict rot = tmp->get_rotation("2d");
00270         float rotate_angle_solution = rot["alpha"];
00271         delete tmp;
00272 
00273         EMData *rot_align_180 = rot_aligned->process("math.rotate.180");
00274 
00275         // Generate the comparison metrics for both rotational candidates
00276         float rot_cmp = rot_aligned->cmp(cmp_name, to, cmp_params);
00277         float rot_180_cmp = rot_align_180->cmp(cmp_name, to, cmp_params);
00278         // Decide on the result
00279         float score = 0.0;
00280         EMData* result = NULL;
00281         if (rot_cmp < rot_180_cmp){
00282                 result = rot_aligned;
00283                 score = rot_cmp;
00284                 delete rot_align_180; rot_align_180 = 0;
00285         } else {
00286                 result = rot_align_180;
00287                 score = rot_180_cmp;
00288                 delete rot_aligned; rot_aligned = 0;
00289                 rotate_angle_solution = rotate_angle_solution-180.0f;
00290         }
00291 
00292 //      Transform* t = get_align_attr("xform.align2d",result);
00293 //      t->set_rotation(Dict("type","2d","alpha",rotate_angle_solution));
00294         Transform tmp2(Dict("type","2d","alpha",rotate_angle_solution));
00295         result->set_attr("xform.align2d",&tmp2);
00296         return result;
00297 }
00298 
00299 
00300 EMData *RotatePrecenterAligner::align(EMData * this_img, EMData *to,
00301                         const string&, const Dict&) const
00302 {
00303         if (!to) {
00304                 return 0;
00305         }
00306 
00307         int ny = this_img->get_ysize();
00308         int size = Util::calc_best_fft_size((int) (M_PI * ny * 1.5));
00309         EMData *e1 = this_img->unwrap(4, ny * 7 / 16, size, 0, 0, 1);
00310         EMData *e2 = to->unwrap(4, ny * 7 / 16, size, 0, 0, 1);
00311         EMData *cf = e1->calc_ccfx(e2, 0, ny);
00312 
00313         float *data = cf->get_data();
00314 
00315         float peak = 0;
00316         int peak_index = 0;
00317         Util::find_max(data, size, &peak, &peak_index);
00318         float a = (float) ((1.0f - 1.0f * peak_index / size) * 180. * 2);
00319 
00320         Transform rot;
00321         rot.set_rotation(Dict("type","2d","alpha",(float)(a*180./M_PI)));
00322         EMData* rslt = this_img->process("xform",Dict("transform",&rot));
00323         rslt->set_attr("xform.align2d",&rot);
00324 //
00325 //      Transform* t = get_set_align_attr("xform.align2d",rslt,this_img);
00326 //      t->set_rotation(Dict("type","2d","alpha",-a));
00327 //
00328 //      EMData* result this_img->transform(Dict("type","2d","alpha",(float)(a*180./M_PI)));
00329 //
00330 //      cf->set_attr("xform.align2d",t);
00331 //      delete t;
00332 //      cf->update();
00333 
00334         if( e1 )
00335         {
00336                 delete e1;
00337                 e1 = 0;
00338         }
00339 
00340         if( e2 )
00341         {
00342                 delete e2;
00343                 e2 = 0;
00344         }
00345 
00346         if (cf) {
00347                 delete cf;
00348                 cf = 0;
00349         }
00350         return rslt;
00351 }
00352 
00353 EMData *RotateTranslateAligner::align(EMData * this_img, EMData *to,
00354                         const string & cmp_name, const Dict& cmp_params) const
00355 {
00356         // Get the 180 degree ambiguously rotationally aligned and its 180 degree rotation counterpart
00357         int rfp_mode = params.set_default("rfp_mode",0);
00358         EMData *rot_align  =  RotationalAligner::align_180_ambiguous(this_img,to,rfp_mode);
00359         Transform * tmp = rot_align->get_attr("xform.align2d");
00360         Dict rot = tmp->get_rotation("2d");
00361         float rotate_angle_solution = rot["alpha"];
00362         delete tmp;
00363 
00364         EMData *rot_align_180 = rot_align->copy();
00365         rot_align_180->process_inplace("math.rotate.180");
00366 
00367         Dict trans_params;
00368         trans_params["intonly"]  = 0;
00369         trans_params["maxshift"] = params.set_default("maxshift", -1);
00370         trans_params["useflcf"]=params.set_default("useflcf",0);
00371 
00372         // Do the first case translational alignment
00373         trans_params["nozero"]   = params.set_default("nozero",false);
00374         EMData* rot_trans = rot_align->align("translational", to, trans_params, cmp_name, cmp_params);
00375         if( rot_align ) { // Clean up
00376                 delete rot_align;
00377                 rot_align = 0;
00378         }
00379 
00380         // Do the second case translational alignment
00381         EMData*  rot_180_trans = rot_align_180->align("translational", to, trans_params, cmp_name, cmp_params);
00382         if( rot_align_180 )     { // Clean up
00383                 delete rot_align_180;
00384                 rot_align_180 = 0;
00385         }
00386 
00387         // Finally decide on the result
00388         float cmp1 = rot_trans->cmp(cmp_name, to, cmp_params);
00389         float cmp2 = rot_180_trans->cmp(cmp_name, to, cmp_params);
00390 
00391         EMData *result = 0;
00392         if (cmp1 < cmp2) { // Assumes smaller is better - thus all comparitors should support "smaller is better"
00393                 if( rot_180_trans )     {
00394                         delete rot_180_trans;
00395                         rot_180_trans = 0;
00396                 }
00397                 result = rot_trans;
00398         }
00399         else {
00400                 if( rot_trans ) {
00401                         delete rot_trans;
00402                         rot_trans = 0;
00403                 }
00404                 result = rot_180_trans;
00405                 rotate_angle_solution -= 180.f;
00406         }
00407 
00408         Transform* t = result->get_attr("xform.align2d");
00409         t->set_rotation(Dict("type","2d","alpha",rotate_angle_solution));
00410         result->set_attr("xform.align2d",t);
00411         delete t;
00412 
00413         return result;
00414 }
00415 
00416 
00417 
00418 
00419 EMData* RotateTranslateFlipAligner::align(EMData * this_img, EMData *to,
00420                                                                                   const string & cmp_name, const Dict& cmp_params) const
00421 {
00422         // Get the non flipped rotational, tranlsationally aligned image
00423         Dict rt_params("maxshift", params["maxshift"], "rfp_mode", params.set_default("rfp_mode",0),"useflcf",params.set_default("useflcf",0));
00424         EMData *rot_trans_align = this_img->align("rotate_translate",to,rt_params,cmp_name, cmp_params);
00425 
00426         // Do the same alignment, but using the flipped version of the image
00427         EMData *flipped = params.set_default("flip", (EMData *) 0);
00428         bool delete_flag = false;
00429         if (flipped == 0) {
00430                 flipped = to->process("xform.flip", Dict("axis", "x"));
00431                 delete_flag = true;
00432         }
00433 
00434         EMData * rot_trans_align_flip = this_img->align("rotate_translate", flipped, rt_params, cmp_name, cmp_params);
00435         Transform* t = rot_trans_align_flip->get_attr("xform.align2d");
00436         t->set_mirror(true);
00437         rot_trans_align_flip->set_attr("xform.align2d",t);
00438         delete t;
00439 
00440         // Now finally decide on what is the best answer
00441         float cmp1 = rot_trans_align->cmp(cmp_name, to, cmp_params);
00442         float cmp2 = rot_trans_align_flip->cmp(cmp_name, flipped, cmp_params);
00443 
00444         if (delete_flag){
00445                 if(flipped) {
00446                         delete flipped;
00447                         flipped = 0;
00448                 }
00449         }
00450 
00451         EMData *result = 0;
00452         if (cmp1 < cmp2 )  {
00453 
00454                 if( rot_trans_align_flip ) {
00455                         delete rot_trans_align_flip;
00456                         rot_trans_align_flip = 0;
00457                 }
00458                 result = rot_trans_align;
00459         }
00460         else {
00461                 if( rot_trans_align ) {
00462                         delete rot_trans_align;
00463                         rot_trans_align = 0;
00464                 }
00465                 result = rot_trans_align_flip;
00466                 result->process_inplace("xform.flip",Dict("axis","x"));
00467         }
00468 
00469         return result;
00470 }
00471 
00472 
00473 
00474 
00475 EMData *RotateFlipAligner::align(EMData * this_img, EMData *to,
00476                         const string& cmp_name, const Dict& cmp_params) const
00477 {
00478         Dict rot_params("rfp_mode",params.set_default("rfp_mode",0));
00479         EMData *r1 = this_img->align("rotational", to, rot_params,cmp_name, cmp_params);
00480 
00481 
00482         EMData* flipped =to->process("xform.flip", Dict("axis", "x"));
00483         EMData *r2 = this_img->align("rotational", flipped,rot_params, cmp_name, cmp_params);
00484         Transform* t = r2->get_attr("xform.align2d");
00485         t->set_mirror(true);
00486         r2->set_attr("xform.align2d",t);
00487         delete t;
00488 
00489         float cmp1 = r1->cmp(cmp_name, to, cmp_params);
00490         float cmp2 = r2->cmp(cmp_name, flipped, cmp_params);
00491 
00492         delete flipped; flipped = 0;
00493 
00494         EMData *result = 0;
00495 
00496         if (cmp1 < cmp2) {
00497                 if( r2 )
00498                 {
00499                         delete r2;
00500                         r2 = 0;
00501                 }
00502                 result = r1;
00503         }
00504         else {
00505                 if( r1 )
00506                 {
00507                         delete r1;
00508                         r1 = 0;
00509                 }
00510                 result = r2;
00511                 result->process_inplace("xform.flip",Dict("axis","x"));
00512         }
00513 
00514         return result;
00515 }
00516 
00517 
00518 // David Woolford says FIXME
00519 // You will note the excessive amount of EMData copying that's going in this function
00520 // This is because functions that are operating on the EMData objects are changing them
00521 // and if you do not use copies the whole algorithm breaks. I did not have time to go
00522 // through and rectify this situation.
00523 // David Woolford says - this problem is related to the fact that many functions that
00524 // take EMData pointers as arguments do not take them as constant pointers to constant
00525 // objects, instead they are treated as raw (completely changeable) pointers. This means
00526 // it's hard to track down which functions are changing the EMData objects, because they
00527 // all do (in name). If this behavior is unavoidable then ignore this comment, however if possible it would
00528 // be good to make things const as much as possible. For example in alignment, technically
00529 // the argument EMData objects (raw pointers) should not be altered... should they?
00530 //
00531 // But const can be very annoying sometimes...
00532 EMData *RTFExhaustiveAligner::align(EMData * this_img, EMData *to,
00533                         const string & cmp_name, const Dict& cmp_params) const
00534 {
00535         EMData *flip = params.set_default("flip", (EMData *) 0);
00536         int maxshift = params.set_default("maxshift", this_img->get_xsize()/8);
00537         if (maxshift < 2) throw InvalidParameterException("maxshift must be greater than or equal to 2");
00538 
00539         int ny = this_img->get_ysize();
00540         int xst = (int) floor(2 * M_PI * ny);
00541         xst = Util::calc_best_fft_size(xst);
00542 
00543         Dict d("n",2);
00544         EMData *to_shrunk_unwrapped = to->process("math.medianshrink",d);
00545 
00546         int to_copy_r2 = to_shrunk_unwrapped->get_ysize() / 2 - 2 - maxshift / 2;
00547         EMData *tmp = to_shrunk_unwrapped->unwrap(4, to_copy_r2, xst / 2, 0, 0, true);
00548         if( to_shrunk_unwrapped )
00549         {
00550                 delete to_shrunk_unwrapped;
00551                 to_shrunk_unwrapped = 0;
00552         }
00553         to_shrunk_unwrapped = tmp;
00554 
00555         EMData *to_shrunk_unwrapped_copy = to_shrunk_unwrapped->copy();
00556         EMData* to_unwrapped = to->unwrap(4, to->get_ysize() / 2 - 2 - maxshift, xst, 0, 0, true);
00557         EMData *to_unwrapped_copy = to_unwrapped->copy();
00558 
00559         bool delete_flipped = true;
00560         EMData *flipped = 0;
00561         if (flip) {
00562                 delete_flipped = false;
00563                 flipped = flip;
00564         }
00565         else {
00566                 flipped = to->process("xform.flip", Dict("axis", "x"));
00567         }
00568         EMData *to_shrunk_flipped_unwrapped = flipped->process("math.medianshrink",d);
00569         tmp = to_shrunk_flipped_unwrapped->unwrap(4, to_copy_r2, xst / 2, 0, 0, true);
00570         if( to_shrunk_flipped_unwrapped )
00571         {
00572                 delete to_shrunk_flipped_unwrapped;
00573                 to_shrunk_flipped_unwrapped = 0;
00574         }
00575         to_shrunk_flipped_unwrapped = tmp;
00576         EMData *to_shrunk_flipped_unwrapped_copy = to_shrunk_flipped_unwrapped->copy();
00577         EMData* to_flip_unwrapped = flipped->unwrap(4, to->get_ysize() / 2 - 2 - maxshift, xst, 0, 0, true);
00578         EMData* to_flip_unwrapped_copy = to_flip_unwrapped->copy();
00579 
00580         if (delete_flipped && flipped != 0) {
00581                 delete flipped;
00582                 flipped = 0;
00583         }
00584 
00585         EMData *this_shrunk_2 = this_img->process("math.medianshrink",d);
00586 
00587         float bestval = FLT_MAX;
00588         float bestang = 0;
00589         int bestflip = 0;
00590         float bestdx = 0;
00591         float bestdy = 0;
00592 
00593         int half_maxshift = maxshift / 2;
00594 
00595         int ur2 = this_shrunk_2->get_ysize() / 2 - 2 - half_maxshift;
00596         for (int dy = -half_maxshift; dy <= half_maxshift; dy += 1) {
00597                 for (int dx = -half_maxshift; dx <= half_maxshift; dx += 1) {
00598 #ifdef  _WIN32
00599                         if (_hypot(dx, dy) <= half_maxshift) {
00600 #else
00601                         if (hypot(dx, dy) <= half_maxshift) {
00602 #endif
00603                                 EMData *uw = this_shrunk_2->unwrap(4, ur2, xst / 2, dx, dy, true);
00604                                 EMData *uwc = uw->copy();
00605                                 EMData *a = uw->calc_ccfx(to_shrunk_unwrapped);
00606 
00607                                 uwc->rotate_x(a->calc_max_index());
00608                                 float cm = uwc->cmp(cmp_name, to_shrunk_unwrapped_copy, cmp_params);
00609                                 if (cm < bestval) {
00610                                         bestval = cm;
00611                                         bestang = (float) (2.0 * M_PI * a->calc_max_index() / a->get_xsize());
00612                                         bestdx = (float)dx;
00613                                         bestdy = (float)dy;
00614                                         bestflip = 0;
00615                                 }
00616 
00617 
00618                                 if( a )
00619                                 {
00620                                         delete a;
00621                                         a = 0;
00622                                 }
00623                                 if( uw )
00624                                 {
00625                                         delete uw;
00626                                         uw = 0;
00627                                 }
00628                                 if( uwc )
00629                                 {
00630                                         delete uwc;
00631                                         uwc = 0;
00632                                 }
00633                                 uw = this_shrunk_2->unwrap(4, ur2, xst / 2, dx, dy, true);
00634                                 uwc = uw->copy();
00635                                 a = uw->calc_ccfx(to_shrunk_flipped_unwrapped);
00636 
00637                                 uwc->rotate_x(a->calc_max_index());
00638                                 cm = uwc->cmp(cmp_name, to_shrunk_flipped_unwrapped_copy, cmp_params);
00639                                 if (cm < bestval) {
00640                                         bestval = cm;
00641                                         bestang = (float) (2.0 * M_PI * a->calc_max_index() / a->get_xsize());
00642                                         bestdx = (float)dx;
00643                                         bestdy = (float)dy;
00644                                         bestflip = 1;
00645                                 }
00646 
00647                                 if( a )
00648                                 {
00649                                         delete a;
00650                                         a = 0;
00651                                 }
00652 
00653                                 if( uw )
00654                                 {
00655                                         delete uw;
00656                                         uw = 0;
00657                                 }
00658                                 if( uwc )
00659                                 {
00660                                         delete uwc;
00661                                         uwc = 0;
00662                                 }
00663                         }
00664                 }
00665         }
00666         if( this_shrunk_2 )
00667         {
00668                 delete this_shrunk_2;
00669                 this_shrunk_2 = 0;
00670         }
00671         if( to_shrunk_unwrapped )
00672         {
00673                 delete to_shrunk_unwrapped;
00674                 to_shrunk_unwrapped = 0;
00675         }
00676         if( to_shrunk_unwrapped_copy )
00677         {
00678                 delete to_shrunk_unwrapped_copy;
00679                 to_shrunk_unwrapped_copy = 0;
00680         }
00681         if( to_shrunk_flipped_unwrapped )
00682         {
00683                 delete to_shrunk_flipped_unwrapped;
00684                 to_shrunk_flipped_unwrapped = 0;
00685         }
00686         if( to_shrunk_flipped_unwrapped_copy )
00687         {
00688                 delete to_shrunk_flipped_unwrapped_copy;
00689                 to_shrunk_flipped_unwrapped_copy = 0;
00690         }
00691         bestdx *= 2;
00692         bestdy *= 2;
00693         bestval = FLT_MAX;
00694 
00695         float bestdx2 = bestdx;
00696         float bestdy2 = bestdy;
00697         // Note I tried steps less than 1.0 (sub pixel precision) and it actually appeared detrimental
00698         // So my advice is to stick with dx += 1.0 etc unless you really are looking to fine tune this
00699         // algorithm
00700         for (float dy = bestdy2 - 3; dy <= bestdy2 + 3; dy += 1.0 ) {
00701                 for (float dx = bestdx2 - 3; dx <= bestdx2 + 3; dx += 1.0 ) {
00702 
00703 #ifdef  _WIN32
00704                         if (_hypot(dx, dy) <= maxshift) {
00705 #else
00706                         if (hypot(dx, dy) <= maxshift) {
00707 #endif
00708                                 EMData *uw = this_img->unwrap(4, this_img->get_ysize() / 2 - 2 - maxshift, xst, (int)dx, (int)dy, true);
00709                                 EMData *uwc = uw->copy();
00710                                 EMData *a = uw->calc_ccfx(to_unwrapped);
00711 
00712                                 uwc->rotate_x(a->calc_max_index());
00713                                 float cm = uwc->cmp(cmp_name, to_unwrapped_copy, cmp_params);
00714 
00715                                 if (cm < bestval) {
00716                                         bestval = cm;
00717                                         bestang = (float)(2.0 * M_PI * a->calc_max_index() / a->get_xsize());
00718                                         bestdx = dx;
00719                                         bestdy = dy;
00720                                         bestflip = 0;
00721                                 }
00722 
00723                                 if( a )
00724                                 {
00725                                         delete a;
00726                                         a = 0;
00727                                 }
00728                                 if( uw )
00729                                 {
00730                                         delete uw;
00731                                         uw = 0;
00732                                 }
00733                                 if( uwc )
00734                                 {
00735                                         delete uwc;
00736                                         uwc = 0;
00737                                 }
00738                                 uw = this_img->unwrap(4, this_img->get_ysize() / 2 - 2 - maxshift, xst, (int)dx, (int)dy, true);
00739                                 uwc = uw->copy();
00740                                 a = uw->calc_ccfx(to_flip_unwrapped);
00741 
00742                                 uwc->rotate_x(a->calc_max_index());
00743                                 cm = uwc->cmp(cmp_name, to_flip_unwrapped_copy, cmp_params);
00744 
00745                                 if (cm < bestval) {
00746                                         bestval = cm;
00747                                         bestang = (float)(2.0 * M_PI * a->calc_max_index() / a->get_xsize());
00748                                         bestdx = dx;
00749                                         bestdy = dy;
00750                                         bestflip = 1;
00751                                 }
00752 
00753                                 if( a )
00754                                 {
00755                                         delete a;
00756                                         a = 0;
00757                                 }
00758                                 if( uw )
00759                                 {
00760                                         delete uw;
00761                                         uw = 0;
00762                                 }
00763                                 if( uwc )
00764                                 {
00765                                         delete uwc;
00766                                         uwc = 0;
00767                                 }
00768                         }
00769                 }
00770         }
00771         if( to_unwrapped ) {delete to_unwrapped;to_unwrapped = 0;}
00772         if( to_shrunk_unwrapped ) {     delete to_shrunk_unwrapped;     to_shrunk_unwrapped = 0;}
00773         if (to_unwrapped_copy) { delete to_unwrapped_copy; to_unwrapped_copy = 0; }
00774         if (to_flip_unwrapped) { delete to_flip_unwrapped; to_flip_unwrapped = 0; }
00775         if (to_flip_unwrapped_copy) { delete to_flip_unwrapped_copy; to_flip_unwrapped_copy = 0;}
00776 
00777         bestang *= (float)EMConsts::rad2deg;
00778         Transform t(Dict("type","2d","alpha",(float)bestang));
00779         t.set_pre_trans(Vec2f(-bestdx,-bestdy));
00780         if (bestflip) {
00781                 t.set_mirror(true);
00782         }
00783 
00784         EMData* ret = this_img->process("xform",Dict("transform",&t));
00785         ret->set_attr("xform.align2d",&t);
00786 
00787         return ret;
00788 }
00789 
00790 
00791 EMData *RTFSlowExhaustiveAligner::align(EMData * this_img, EMData *to,
00792                         const string & cmp_name, const Dict& cmp_params) const
00793 {
00794 
00795         EMData *flip = params.set_default("flip", (EMData *) 0);
00796         int maxshift = params.set_default("maxshift", -1);
00797 
00798         EMData *flipped = 0;
00799 
00800         bool delete_flipped = true;
00801         if (flip) {
00802                 delete_flipped = false;
00803                 flipped = flip;
00804         }
00805         else {
00806                 flipped = to->process("xform.flip", Dict("axis", "x"));
00807         }
00808 
00809         int nx = this_img->get_xsize();
00810 
00811         if (maxshift < 0) {
00812                 maxshift = nx / 10;
00813         }
00814 
00815         float angle_step =  params.set_default("angstep", 0.0f);
00816         if ( angle_step == 0 ) angle_step = atan2(2.0f, (float)nx);
00817         else {
00818                 angle_step *= (float)EMConsts::deg2rad; //convert to radians
00819         }
00820         float trans_step =  params.set_default("transtep",1.0f);
00821 
00822         if (trans_step <= 0) throw InvalidParameterException("transstep must be greater than 0");
00823         if (angle_step <= 0) throw InvalidParameterException("angstep must be greater than 0");
00824 
00825 
00826         Dict shrinkfactor("n",2);
00827         EMData *this_img_shrink = this_img->process("math.medianshrink",shrinkfactor);
00828         EMData *to_shrunk = to->process("math.medianshrink",shrinkfactor);
00829         EMData *flipped_shrunk = flipped->process("math.medianshrink",shrinkfactor);
00830 
00831         int bestflip = 0;
00832         float bestdx = 0;
00833         float bestdy = 0;
00834 
00835         float bestang = 0;
00836         float bestval = FLT_MAX;
00837 
00838         int half_maxshift = maxshift / 2;
00839 
00840 
00841         for (int dy = -half_maxshift; dy <= half_maxshift; ++dy) {
00842                 for (int dx = -half_maxshift; dx <= half_maxshift; ++dx) {
00843                         if (hypot(dx, dy) <= maxshift) {
00844                                 for (float ang = -angle_step * 2.0f; ang <= (float)2 * M_PI; ang += angle_step * 4.0f) {
00845                                         EMData v(*this_img_shrink);
00846                                         Transform t(Dict("type","2d","alpha",static_cast<float>(ang*EMConsts::rad2deg)));
00847                                         t.set_trans((float)dx,(float)dy);
00848                                         v.transform(t);
00849 //                                      v.rotate_translate(ang*EMConsts::rad2deg, 0.0f, 0.0f, (float)dx, (float)dy, 0.0f);
00850 
00851                                         float lc = v.cmp(cmp_name, to_shrunk, cmp_params);
00852 
00853                                         if (lc < bestval) {
00854                                                 bestval = lc;
00855                                                 bestang = ang;
00856                                                 bestdx = (float)dx;
00857                                                 bestdy = (float)dy;
00858                                                 bestflip = 0;
00859                                         }
00860 
00861                                         lc = v.cmp(cmp_name,flipped_shrunk , cmp_params);
00862                                         if (lc < bestval) {
00863                                                 bestval = lc;
00864                                                 bestang = ang;
00865                                                 bestdx = (float)dx;
00866                                                 bestdy = (float)dy;
00867                                                 bestflip = 1;
00868                                         }
00869                                 }
00870                         }
00871                 }
00872         }
00873 
00874         if( to_shrunk )
00875         {
00876                 delete to_shrunk;
00877                 to_shrunk = 0;
00878         }
00879         if( flipped_shrunk )
00880         {
00881                 delete flipped_shrunk;
00882                 flipped_shrunk = 0;
00883         }
00884         if( this_img_shrink )
00885         {
00886                 delete this_img_shrink;
00887                 this_img_shrink = 0;
00888         }
00889 
00890         bestdx *= 2;
00891         bestdy *= 2;
00892         bestval = FLT_MAX;
00893 
00894         float bestdx2 = bestdx;
00895         float bestdy2 = bestdy;
00896         float bestang2 = bestang;
00897 
00898         for (float dy = bestdy2 - 3; dy <= bestdy2 + 3; dy += trans_step) {
00899                 for (float dx = bestdx2 - 3; dx <= bestdx2 + 3; dx += trans_step) {
00900                         if (hypot(dx, dy) <= maxshift) {
00901                                 for (float ang = bestang2 - angle_step * 6.0f; ang <= bestang2 + angle_step * 6.0f; ang += angle_step) {
00902                                         EMData v(*this_img);
00903                                         Transform t(Dict("type","2d","alpha",static_cast<float>(ang*EMConsts::rad2deg)));
00904                                         t.set_trans(dx,dy);
00905                                         v.transform(t);
00906 //                                      v.rotate_translate(ang*EMConsts::rad2deg, 0.0f, 0.0f, (float)dx, (float)dy, 0.0f);
00907 
00908                                         float lc = v.cmp(cmp_name, to, cmp_params);
00909 
00910                                         if (lc < bestval) {
00911                                                 bestval = lc;
00912                                                 bestang = ang;
00913                                                 bestdx = dx;
00914                                                 bestdy = dy;
00915                                                 bestflip = 0;
00916                                         }
00917 
00918                                         lc = v.cmp(cmp_name, flipped, cmp_params);
00919 
00920                                         if (lc < bestval) {
00921                                                 bestval = lc;
00922                                                 bestang = ang;
00923                                                 bestdx = dx;
00924                                                 bestdy = dy;
00925                                                 bestflip = 1;
00926                                         }
00927                                 }
00928                         }
00929                 }
00930         }
00931 
00932         if (delete_flipped) { delete flipped; flipped = 0; }
00933 
00934         bestang *= (float)EMConsts::rad2deg;
00935         Transform t(Dict("type","2d","alpha",(float)bestang));
00936         t.set_trans(bestdx,bestdy);
00937 
00938         if (bestflip) {
00939                 t.set_mirror(true);
00940         }
00941 
00942         EMData* rslt = this_img->process("xform",Dict("transform",&t));
00943         rslt->set_attr("xform.align2d",&t);
00944 
00945         return rslt;
00946 }
00947 
00948 
00949 
00950 static double refalifn(const gsl_vector * v, void *params)
00951 {
00952         Dict *dict = (Dict *) params;
00953 
00954         double x = gsl_vector_get(v, 0);
00955         double y = gsl_vector_get(v, 1);
00956         double a = gsl_vector_get(v, 2);
00957 
00958         EMData *this_img = (*dict)["this"];
00959         EMData *with = (*dict)["with"];
00960         bool mirror = (*dict)["mirror"];
00961 
00962 //      float mean = (float)this_img->get_attr("mean");
00963 //      if ( Util::goodf(&mean) ) {
00964 //              //cout << "tmps mean is nan even before rotation" << endl;
00965 //      }
00966 
00967         Transform t(Dict("type","2d","alpha",static_cast<float>(a)));
00968 //      Transform3D t3d(Transform3D::EMAN, (float)a, 0.0f, 0.0f);
00969 //      t3d.set_posttrans( (float) x, (float) y);
00970 //      tmp->rotate_translate(t3d);
00971         t.set_trans((float)x,(float)y);
00972         t.set_mirror(mirror);
00973         EMData *tmp = this_img->process("xform",Dict("transform",&t));
00974 
00975         Cmp* c = (Cmp*) ((void*)(*dict)["cmp"]);
00976         double result = c->cmp(tmp,with);
00977 
00978         // DELETE AT SOME STAGE, USEFUL FOR PRERELEASE STUFF
00979         //      float test_result = (float)result;
00980 //      if ( Util::goodf(&test_result) ) {
00981 //              cout << "result " << result << " " << x << " " << y << " " << a << endl;
00982 //              cout << (float)this_img->get_attr("mean") << " " << (float)tmp->get_attr("mean") << " " << (float)with->get_attr("mean") << endl;
00983 //              tmp->write_image("tmp.hdf");
00984 //              with->write_image("with.hdf");
00985 //              this_img->write_image("this_img.hdf");
00986 //              EMData* t = this_img->copy();
00987 //              cout << (float)t->get_attr("mean") << endl;
00988 //              t->rotate_translate( t3d );
00989 //              cout << (float)t->get_attr("mean") << endl;
00990 //              cout << "exit" << endl;
00992 //              cout << (float)t->get_attr("mean") << endl;
00993 //              cout << "now exit" << endl;
00994 //              delete t;
00995 //      }
00996 
00997 
00998         if ( tmp != 0 ) delete tmp;
00999 
01000         return result;
01001 }
01002 
01003 static double refalifnfast(const gsl_vector * v, void *params)
01004 {
01005         Dict *dict = (Dict *) params;
01006         EMData *this_img = (*dict)["this"];
01007         EMData *img_to = (*dict)["with"];
01008         bool mirror = (*dict)["mirror"];
01009 
01010         double x = gsl_vector_get(v, 0);
01011         double y = gsl_vector_get(v, 1);
01012         double a = gsl_vector_get(v, 2);
01013 
01014         double r = this_img->dot_rotate_translate(img_to, (float)x, (float)y, (float)a, mirror);
01015         int nsec = this_img->get_xsize() * this_img->get_ysize();
01016         double result = 1.0 - r / nsec;
01017 
01018 //      cout << result << " x " << x << " y " << y << " az " << a <<  endl;
01019         return result;
01020 }
01021 
01022 
01023 EMData *RefineAligner::align(EMData * this_img, EMData *to,
01024         const string & cmp_name, const Dict& cmp_params) const
01025 {
01026 
01027         if (!to) {
01028                 return 0;
01029         }
01030 
01031         EMData *result;
01032         int mode = params.set_default("mode", 0);
01033         float saz = 0.0;
01034         float sdx = 0.0;
01035         float sdy = 0.0;
01036         bool mirror = false;
01037         Transform* t;
01038         if (params.has_key("xform.align2d") ) {
01039                 t = params["xform.align2d"];
01040                 Dict params = t->get_params("2d");
01041                 saz = params["alpha"];
01042                 sdx = params["tx"];
01043                 sdy = params["ty"];
01044                 mirror = params["mirror"];
01045 
01046         } else {
01047                 t = new Transform(); // is the identity
01048         }
01049 
01050         // We do this to prevent the GSL routine from crashing on an invalid alignment
01051         if ((float)(this_img->get_attr("sigma"))==0.0 || (float)(to->get_attr("sigma"))==0.0) {
01052                 result = this_img->process("xform",Dict("transform",t));
01053                 result->set_attr("xform.align2d",t);
01054                 delete t;
01055                 return result;
01056         }
01057 
01058         int np = 3;
01059         Dict gsl_params;
01060         gsl_params["this"] = this_img;
01061         gsl_params["with"] = to;
01062         gsl_params["snr"]  = params["snr"];
01063         gsl_params["mirror"] = mirror;
01064 
01065 
01066 
01067         const gsl_multimin_fminimizer_type *T = gsl_multimin_fminimizer_nmsimplex;
01068         gsl_vector *ss = gsl_vector_alloc(np);
01069 
01070         float stepx = params.set_default("stepx",1.0f);
01071         float stepy = params.set_default("stepy",1.0f);
01072         // Default step is 5 degree - note in EMAN1 it was 0.1 radians
01073         float stepaz = params.set_default("stepaz",5.0f);
01074 
01075         gsl_vector_set(ss, 0, stepx);
01076         gsl_vector_set(ss, 1, stepy);
01077         gsl_vector_set(ss, 2, stepaz);
01078 
01079         gsl_vector *x = gsl_vector_alloc(np);
01080         gsl_vector_set(x, 0, sdx);
01081         gsl_vector_set(x, 1, sdy);
01082         gsl_vector_set(x, 2, saz);
01083 
01084         Cmp *c = 0;
01085 
01086         gsl_multimin_function minex_func;
01087         if (mode == 2) {
01088                 minex_func.f = &refalifnfast;
01089         }
01090         else {
01091                 c = Factory < Cmp >::get(cmp_name, cmp_params);
01092                 gsl_params["cmp"] = (void *) c;
01093                 minex_func.f = &refalifn;
01094         }
01095 
01096         minex_func.n = np;
01097         minex_func.params = (void *) &gsl_params;
01098 
01099         gsl_multimin_fminimizer *s = gsl_multimin_fminimizer_alloc(T, np);
01100         gsl_multimin_fminimizer_set(s, &minex_func, x, ss);
01101 
01102         int rval = GSL_CONTINUE;
01103         int status = GSL_SUCCESS;
01104         int iter = 1;
01105 
01106         float precision = params.set_default("precision",0.04f);
01107         int maxiter = params.set_default("maxiter",28);
01108 
01109 //      printf("Refine sx=%1.2f sy=%1.2f sa=%1.2f prec=%1.4f maxit=%d\n",stepx,stepy,stepaz,precision,maxiter);
01110 //      printf("%1.2f %1.2f %1.1f  ->",(float)gsl_vector_get(s->x, 0),(float)gsl_vector_get(s->x, 1),(float)gsl_vector_get(s->x, 2));
01111 
01112         while (rval == GSL_CONTINUE && iter < maxiter) {
01113                 iter++;
01114                 status = gsl_multimin_fminimizer_iterate(s);
01115                 if (status) {
01116                         break;
01117                 }
01118                 rval = gsl_multimin_test_size(gsl_multimin_fminimizer_size(s), precision);
01119         }
01120 
01121         int maxshift = params.set_default("maxshift",-1);
01122 
01123         if (maxshift <= 0) {
01124                 maxshift = this_img->get_xsize() / 4;
01125         }
01126         float fmaxshift = static_cast<float>(maxshift);
01127         if ( fmaxshift >= fabs((float)gsl_vector_get(s->x, 0)) && fmaxshift >= fabs((float)gsl_vector_get(s->x, 1))  )
01128         {
01129 //              printf(" Refine good %1.2f %1.2f %1.1f\n",(float)gsl_vector_get(s->x, 0),(float)gsl_vector_get(s->x, 1),(float)gsl_vector_get(s->x, 2));
01130                 Transform  tsoln(Dict("type","2d","alpha",(float)gsl_vector_get(s->x, 2)));
01131                 tsoln.set_mirror(mirror);
01132                 tsoln.set_trans((float)gsl_vector_get(s->x, 0),(float)gsl_vector_get(s->x, 1));
01133                 result = this_img->process("xform",Dict("transform",&tsoln));
01134                 result->set_attr("xform.align2d",&tsoln);
01135         } else { // The refine aligner failed - this shift went beyond the max shift
01136 //              printf(" Refine Failed %1.2f %1.2f %1.1f\n",(float)gsl_vector_get(s->x, 0),(float)gsl_vector_get(s->x, 1),(float)gsl_vector_get(s->x, 2));
01137                 result = this_img->process("xform",Dict("transform",t));
01138                 result->set_attr("xform.align2d",t);
01139         }
01140 
01141         delete t;
01142         t = 0;
01143 
01144         gsl_vector_free(x);
01145         gsl_vector_free(ss);
01146         gsl_multimin_fminimizer_free(s);
01147 
01148         if ( c != 0 ) delete c;
01149         return result;
01150 }
01151 
01152 static Transform refalin3d_perturb(const Transform*const t, const float& delta, const float& arc, const float& phi, const float& x, const float& y, const float& z)
01153 {
01154         Dict orig_params = t->get_params("eman");
01155         float orig_phi = orig_params["phi"];
01156         float orig_x = orig_params["tx"];
01157         float orig_y = orig_params["ty"];
01158         float orig_z = orig_params["tz"];
01159         orig_params["phi"] = 0;
01160         orig_params["tx"] = 0;
01161         orig_params["ty"] = 0;
01162         orig_params["tz"] = 0;
01163         Transform t_no_phi(orig_params);
01164 
01165         Vec3f zz(0,0,1);
01166 
01167         Vec3f vv  = t_no_phi.transpose()*zz;
01168         Vec3f normal = Vec3f(-vv[2],0,-vv[0]);
01169 
01170         normal.normalize();
01171 
01172         Dict d;
01173         d["type"] = "spin";
01174         d["Omega"] = arc;
01175         d["n1"] = vv[0];
01176         d["n2"] = vv[1];
01177         d["n3"] = vv[2];
01178 
01179         Transform q(d);
01180 
01181         Vec3f  rot_vec = q*normal;
01182         rot_vec.normalize();
01183 
01184         Dict e;
01185         e["type"] = "spin";
01186         e["Omega"] = delta;
01187         e["n1"] = rot_vec[0];
01188         e["n2"] = rot_vec[1];
01189         e["n3"] = rot_vec[2];
01190 
01191         Transform perturb(e);
01192 
01193         Dict g;
01194         g["type"] = "eman";
01195         g["alt"] = 0;
01196         g["az"] = 0;
01197         g["phi"] = 0*phi+orig_phi;
01198 
01199         Transform phi_rot(g);
01200         Transform soln = t_no_phi*perturb*phi_rot;
01201         soln.set_trans(x+orig_x,y+orig_y,z+orig_z);
01202 
01203         Dict params = soln.get_params("eman");
01204         return soln;
01205 }
01206 
01207 static double refalifn3d(const gsl_vector * v, void *params)
01208 {
01209         Dict *dict = (Dict *) params;
01210         double x = gsl_vector_get(v, 0);
01211         double y = gsl_vector_get(v, 1);
01212         double z = gsl_vector_get(v, 2);
01213         double arc = gsl_vector_get(v, 3);
01214         double delta = gsl_vector_get(v, 4);
01215         double phi = gsl_vector_get(v, 5);
01216         EMData *this_img = (*dict)["this"];
01217         EMData *with = (*dict)["with"];
01218 //      bool mirror = (*dict)["mirror"];
01219 
01220         Transform* t = (*dict)["transform"];
01221 
01222         Transform soln = refalin3d_perturb(t,(float)delta,(float)arc,(float)phi,(float)x,(float)y,(float)z);
01223 
01224         EMData *tmp = this_img->process("xform",Dict("transform",&soln));
01225         Cmp* c = (Cmp*) ((void*)(*dict)["cmp"]);
01226         double result = c->cmp(tmp,with);
01227         if ( tmp != 0 ) delete tmp;
01228         delete t; t = 0;
01229 //      cout << result << " " << az << " " << alt << " " << phi << " " << x << " " << y << " " << z << endl;
01230         return result;
01231 }
01232 
01233 /*
01234 static double refalifn3d(const gsl_vector * v, void *params)
01235 {
01236         Dict *dict = (Dict *) params;
01237 
01238         double x = gsl_vector_get(v, 0);
01239         double y = gsl_vector_get(v, 1);
01240         double z = gsl_vector_get(v, 2);
01241         double az = gsl_vector_get(v, 3);
01242         double alt = gsl_vector_get(v, 4);
01243         double phi = gsl_vector_get(v, 5);
01244         EMData *this_img = (*dict)["this"];
01245         EMData *with = (*dict)["with"];
01246         bool mirror = (*dict)["mirror"];
01247 
01248         Dict d("type","eman","az",static_cast<float>(az));
01249         d["alt"] = static_cast<float>(alt);
01250         d["phi"] = static_cast<float>(phi);
01251         Transform t(d);
01252         t.set_trans((float)x,(float)y,(float)z);
01253         t.set_mirror(mirror);
01254         EMData *tmp = this_img->process("xform",Dict("transform",&t));
01255 
01256         Cmp* c = (Cmp*) ((void*)(*dict)["cmp"]);
01257         double result = c->cmp(tmp,with);
01258         if ( tmp != 0 ) delete tmp;
01259 
01260         return result;
01261 }*/
01262 
01263 EMData* Refine3DAligner::align(EMData * this_img, EMData *to,
01264         const string & cmp_name, const Dict& cmp_params) const
01265 {
01266 
01267         if (!to || !this_img) throw NullPointerException("Input image is null"); // not sure if this is necessary, it was there before I started
01268 
01269         if (to->get_ndim() != 3 || this_img->get_ndim() != 3) throw ImageDimensionException("The Refine3D aligner only works for 3D images");
01270 
01271         float saz = 0.0;
01272         float sphi = 0.0;
01273         float salt = 0.0;
01274         float sdx = 0.0;
01275         float sdy = 0.0;
01276         float sdz = 0.0;
01277         bool mirror = false;
01278         Transform* t;
01279         if (params.has_key("xform.align3d") ) {
01280                 // Unlike the 2d refine aligner, this class doesn't require the starting transform's
01281                 // parameters to form the starting guess. Instead the Transform itself
01282                 // is perturbed carefully (using quaternion rotation) to overcome problems that arise
01283                 // when you use orthogonally-based Euler angles
01284                 t = params["xform.align3d"];
01285         }else {
01286                 t = new Transform(); // is the identity
01287         }
01288 
01289         int np = 6; // the number of dimensions
01290         Dict gsl_params;
01291         gsl_params["this"] = this_img;
01292         gsl_params["with"] = to;
01293         gsl_params["snr"]  = params["snr"];
01294         gsl_params["mirror"] = mirror;
01295         gsl_params["transform"] = t;
01296 
01297         const gsl_multimin_fminimizer_type *T = gsl_multimin_fminimizer_nmsimplex;
01298         gsl_vector *ss = gsl_vector_alloc(np);
01299 
01300         float stepx = params.set_default("stepx",1.0f);
01301         float stepy = params.set_default("stepy",1.0f);
01302         float stepz = params.set_default("stepz",1.0f);
01303         // Default step is 5 degree - note in EMAN1 it was 0.1 radians
01304         float half_circle_step = 180.0f; // This shouldn't be altered
01305         float stepphi = params.set_default("stephi",5.0f);
01306         float stepdelta = params.set_default("stepdelta",5.0f);
01307 
01308         gsl_vector_set(ss, 0, stepx);
01309         gsl_vector_set(ss, 1, stepy);
01310         gsl_vector_set(ss, 2, stepz);
01311         gsl_vector_set(ss, 3, half_circle_step);
01312         gsl_vector_set(ss, 4, stepdelta);
01313         gsl_vector_set(ss, 5, stepphi);
01314 
01315         gsl_vector *x = gsl_vector_alloc(np);
01316         gsl_vector_set(x, 0, sdx);
01317         gsl_vector_set(x, 1, sdy);
01318         gsl_vector_set(x, 2, sdz);
01319         gsl_vector_set(x, 3, saz);
01320         gsl_vector_set(x, 4, salt);
01321         gsl_vector_set(x, 5, sphi);
01322 
01323         gsl_multimin_function minex_func;
01324         Cmp *c = Factory < Cmp >::get(cmp_name, cmp_params);
01325         gsl_params["cmp"] = (void *) c;
01326         minex_func.f = &refalifn3d;
01327 
01328         minex_func.n = np;
01329         minex_func.params = (void *) &gsl_params;
01330 
01331         gsl_multimin_fminimizer *s = gsl_multimin_fminimizer_alloc(T, np);
01332         gsl_multimin_fminimizer_set(s, &minex_func, x, ss);
01333 
01334         int rval = GSL_CONTINUE;
01335         int status = GSL_SUCCESS;
01336         int iter = 1;
01337 
01338         float precision = params.set_default("precision",0.04f);
01339         int maxiter = params.set_default("maxiter",60);
01340         while (rval == GSL_CONTINUE && iter < maxiter) {
01341                 iter++;
01342                 status = gsl_multimin_fminimizer_iterate(s);
01343                 if (status) {
01344                         break;
01345                 }
01346                 rval = gsl_multimin_test_size(gsl_multimin_fminimizer_size(s), precision);
01347         }
01348 
01349         int maxshift = params.set_default("maxshift",-1);
01350 
01351         if (maxshift <= 0) {
01352                 maxshift = this_img->get_xsize() / 4;
01353         }
01354         float fmaxshift = static_cast<float>(maxshift);
01355         EMData *result;
01356         if ( fmaxshift >= (float)gsl_vector_get(s->x, 0) && fmaxshift >= (float)gsl_vector_get(s->x, 1)  && fmaxshift >= (float)gsl_vector_get(s->x, 2))
01357         {
01358                 float x = (float)gsl_vector_get(s->x, 0);
01359                 float y = (float)gsl_vector_get(s->x, 1);
01360                 float z = (float)gsl_vector_get(s->x, 2);
01361 
01362                 float arc = (float)gsl_vector_get(s->x, 3);
01363                 float delta = (float)gsl_vector_get(s->x, 4);
01364                 float phi = (float)gsl_vector_get(s->x, 5);
01365 
01366                 Transform tsoln = refalin3d_perturb(t,delta,arc,phi,x,y,z);
01367 
01368                 result = this_img->process("xform",Dict("transform",&tsoln));
01369                 result->set_attr("xform.align3d",&tsoln);
01370 
01371         } else { // The refine aligner failed - this shift went beyond the max shift
01372                 result = this_img->process("xform",Dict("transform",t));
01373                 result->set_attr("xform.align3d",t);
01374         }
01375 
01376         delete t;
01377         t = 0;
01378 
01379         gsl_vector_free(x);
01380         gsl_vector_free(ss);
01381         gsl_multimin_fminimizer_free(s);
01382 
01383         if ( c != 0 ) delete c;
01384         return result;
01385 }
01386 
01387 EMData* RT3DGridAligner::align(EMData * this_img, EMData *to, const string & cmp_name, const Dict& cmp_params) const
01388 {
01389 
01390         vector<Dict> alis = xform_align_nbest(this_img,to,1,cmp_name,cmp_params);
01391 
01392         Dict t;
01393         Transform* tr = (Transform*) alis[0]["xform.align3d"];
01394         t["transform"] = tr;
01395         EMData* soln = this_img->process("xform",t);
01396         soln->set_attr("xform.align3d",tr);
01397         delete tr; tr = 0;
01398 
01399         return soln;
01400 
01401 }
01402 
01403 vector<Dict> RT3DGridAligner::xform_align_nbest(EMData * this_img, EMData * to, const unsigned int nsoln, const string & cmp_name, const Dict& cmp_params) const {
01404 
01405         if ( this_img->get_ndim() != 3 || to->get_ndim() != 3 ) {
01406                 throw ImageDimensionException("This aligner only works for 3D images");
01407         }
01408 
01409         int searchx = 0;
01410         int searchy = 0;
01411         int searchz = 0;
01412 
01413         if (params.has_key("search")) {
01414                 vector<string> check;
01415                 check.push_back("searchx");
01416                 check.push_back("searchy");
01417                 check.push_back("searchz");
01418                 for(vector<string>::const_iterator cit = check.begin(); cit != check.end(); ++cit) {
01419                         if (params.has_key(*cit)) throw InvalidParameterException("The search parameter is mutually exclusive of the searchx, searchy, and searchz parameters");
01420                 }
01421                 int search  = params["search"];
01422                 searchx = search;
01423                 searchy = search;
01424                 searchz = search;
01425         } else {
01426                 searchx = params.set_default("searchx",3);
01427                 searchy = params.set_default("searchy",3);
01428                 searchz = params.set_default("searchz",3);
01429         }
01430 
01431         float ralt = params.set_default("ralt",180.f);
01432         float rphi = params.set_default("rphi",180.f);
01433         float raz = params.set_default("raz",180.f);
01434         float dalt = params.set_default("dalt",10.f);
01435         float daz = params.set_default("daz",10.f);
01436         float dphi = params.set_default("dphi",10.f);
01437         float threshold = params.set_default("threshold",0.f);
01438         if (threshold < 0.0f) throw InvalidParameterException("The threshold parameter must be greater than or equal to zero");
01439         bool verbose = params.set_default("verbose",false);
01440 
01441         vector<Dict> solns;
01442         if (nsoln == 0) return solns; // What was the user thinking?
01443         for (unsigned int i = 0; i < nsoln; ++i ) {
01444                 Dict d;
01445                 d["score"] = 1.e24;
01446                 Transform t; // identity by default
01447                 d["xform.align3d"] = &t; // deep copy is going on here
01448                 solns.push_back(d);
01449         }
01450         Dict d;
01451         d["type"] = "eman"; // d is used in the loop below
01452         for ( float alt = 0.0f; alt <= ralt; alt += dalt) {
01453                 // An optimization for the range of az is made at the top of the sphere
01454                 // If you think about it, this is just a coarse way of making this approach slightly more efficient
01455                 if (verbose) {
01456                         cout << "Trying angle " << alt << endl;
01457                 }
01458 
01459                 float begin_az = -raz;
01460                 float end_az = raz;
01461                 if (alt == 0.0f) {
01462                         begin_az = 0.0f;
01463                         end_az = 0.0f;
01464                 }
01465 
01466                 for ( float az = begin_az; az <= end_az; az += daz ){
01467                         for( float phi = -rphi-az; phi <= rphi-az; phi += dphi ) {
01468                                 d["alt"] = alt;
01469                                 d["phi"] = phi;
01470                                 d["az"] = az;
01471                                 Transform t(d);
01472                                 EMData* transformed = this_img->process("xform",Dict("transform",&t));
01473                                 EMData* ccf = transformed->calc_ccf(to);
01474 
01475                                 IntPoint point = ccf->calc_max_location_wrap(searchx,searchy,searchz);
01476                                 Dict altered_cmp_params(cmp_params);
01477                                 if (cmp_name == "dot.tomo") {
01478                                         altered_cmp_params["ccf"] = ccf;
01479                                         altered_cmp_params["tx"] = point[0];
01480                                         altered_cmp_params["ty"] = point[1];
01481                                         altered_cmp_params["tz"] = point[2];
01482                                 }
01483 
01484                                 float best_score = transformed->cmp(cmp_name,to,altered_cmp_params);
01485                                 delete transformed; transformed =0;
01486                                 delete ccf; ccf = 0;
01487 
01488                                 unsigned int j = 0;
01489                                 for ( vector<Dict>::iterator it = solns.begin(); it != solns.end(); ++it, ++j ) {
01490                                         if ( (float)(*it)["score"] > best_score ) {  // Note greater than - EMAN2 preferes minimums as a matter of policy
01491                                                 vector<Dict>::reverse_iterator rit = solns.rbegin();
01492                                                 copy(rit+1,solns.rend()-j,rit);
01493                                                 Dict& d = (*it);
01494                                                 d["score"] = best_score;
01495                                                 d["xform.align3d"] = &t;
01496                                                 break;
01497                                         }
01498                                 }
01499                         }
01500                 }
01501         }
01502 
01503         return solns;
01504 
01505 }
01506 
01507 EMData* RT3DSphereAligner::align(EMData * this_img, EMData *to, const string & cmp_name, const Dict& cmp_params) const
01508 {
01509 
01510         vector<Dict> alis = xform_align_nbest(this_img,to,1,cmp_name,cmp_params);
01511 
01512         Dict t;
01513         Transform* tr = (Transform*) alis[0]["xform.align3d"];
01514         t["transform"] = tr;
01515         EMData* soln = this_img->process("xform",t);
01516         soln->set_attr("xform.align3d",tr);
01517         delete tr; tr = 0;
01518 
01519         return soln;
01520 
01521 }
01522 
01523 vector<Dict> RT3DSphereAligner::xform_align_nbest(EMData * this_img, EMData * to, const unsigned int nsoln, const string & cmp_name, const Dict& cmp_params) const {
01524 
01525         if ( this_img->get_ndim() != 3 || to->get_ndim() != 3 ) {
01526                 throw ImageDimensionException("This aligner only works for 3D images");
01527         }
01528 
01529         int searchx = 0;
01530         int searchy = 0;
01531         int searchz = 0;
01532 
01533         if (params.has_key("search")) {
01534                 vector<string> check;
01535                 check.push_back("searchx");
01536                 check.push_back("searchy");
01537                 check.push_back("searchz");
01538                 for(vector<string>::const_iterator cit = check.begin(); cit != check.end(); ++cit) {
01539                         if (params.has_key(*cit)) throw InvalidParameterException("The search parameter is mutually exclusive of the searchx, searchy, and searchz parameters");
01540                 }
01541                 int search  = params["search"];
01542                 searchx = search;
01543                 searchy = search;
01544                 searchz = search;
01545         } else {
01546                 searchx = params.set_default("searchx",3);
01547                 searchy = params.set_default("searchy",3);
01548                 searchz = params.set_default("searchz",3);
01549         }
01550 
01551         float rphi = params.set_default("rphi",180.f);
01552         float dphi = params.set_default("dphi",10.f);
01553         float threshold = params.set_default("threshold",0.f);
01554         if (threshold < 0.0f) throw InvalidParameterException("The threshold parameter must be greater than or equal to zero");
01555         bool verbose = params.set_default("verbose",false);
01556 
01557         vector<Dict> solns;
01558         if (nsoln == 0) return solns; // What was the user thinking?
01559         for (unsigned int i = 0; i < nsoln; ++i ) {
01560                 Dict d;
01561                 d["score"] = 1.e24;
01562                 Transform t; // identity by default
01563                 d["xform.align3d"] = &t; // deep copy is going on here
01564                 solns.push_back(d);
01565         }
01566 
01567         Dict d;
01568         d["inc_mirror"] = true; // This should probably always be true for 3D case. If it ever changes we might have to make inc_mirror a parameter
01569         if ( params.has_key("delta") && params.has_key("n") ) {
01570                 throw InvalidParameterException("The delta and n parameters are mutually exclusive in the RT3DSphereAligner aligner");
01571         } else if ( params.has_key("n") ) {
01572                 d["n"] = params["n"];
01573         } else {
01574                 // If they didn't specify either then grab the default delta - if they did supply delta we're still safe doing this
01575                 d["delta"] = params.set_default("delta",10.f);
01576         }
01577 
01578         Symmetry3D* sym = Factory<Symmetry3D>::get((string)params.set_default("sym","c1"));
01579         vector<Transform> transforms = sym->gen_orientations((string)params.set_default("orientgen","eman"),d);
01580 
01581         float verbose_alt = -1.0f;;
01582         for(vector<Transform>::const_iterator trans_it = transforms.begin(); trans_it != transforms.end(); trans_it++) {
01583                 Dict params = trans_it->get_params("eman");
01584                 float az = params["az"];
01585                 if (verbose) {
01586                         float alt = params["alt"];
01587                         if ( alt != verbose_alt ) {
01588                                 verbose_alt = alt;
01589                                 cout << "Trying angle " << alt << endl;
01590                         }
01591                 }
01592                 for( float phi = -rphi-az; phi <= rphi-az; phi += dphi ) {
01593                         params["phi"] = phi;
01594                         Transform t(params);
01595                         EMData* transformed = this_img->process("xform",Dict("transform",&t));
01596                         EMData* ccf = transformed->calc_ccf(to);
01597                         IntPoint point = ccf->calc_max_location_wrap(searchx,searchy,searchz);
01598                         Dict altered_cmp_params(cmp_params);
01599                         if (cmp_name == "dot.tomo") {
01600                                 altered_cmp_params["ccf"] = ccf;
01601                                 altered_cmp_params["tx"] = point[0];
01602                                 altered_cmp_params["ty"] = point[1];
01603                                 altered_cmp_params["tz"] = point[2];
01604                         }
01605 
01606                         float best_score = transformed->cmp(cmp_name,to,altered_cmp_params);
01607                         delete transformed; transformed =0;
01608                         delete ccf; ccf =0;
01609 
01610                         unsigned int j = 0;
01611                         for ( vector<Dict>::iterator it = solns.begin(); it != solns.end(); ++it, ++j ) {
01612                                 if ( (float)(*it)["score"] > best_score ) { // Note greater than - EMAN2 preferes minimums as a matter of policy
01613                                         vector<Dict>::reverse_iterator rit = solns.rbegin();
01614                                         copy(rit+1,solns.rend()-j,rit);
01615                                         Dict& d = (*it);
01616                                         d["score"] = best_score;
01617                                         d["xform.align3d"] = &t; // deep copy is going on here
01618                                         break;
01619                                 }
01620                         }
01621 
01622                 }
01623         }
01624         delete sym; sym = 0;
01625         return solns;
01626 
01627 }
01628 
01629 #ifdef  FFTW2
01630 namespace {
01631 float frm_2d_Align(EMData *this_img, EMData *to, float *frm2dhhat, EMData* selfpcsfft, int p_max_input,int rsize, float &com_this_x, float &com_this_y, float &com_with_x, float &com_with_y,const string & cmp_name, const Dict& cmp_params)
01632 {
01633         int size=rsize;
01634         float dx,dy;
01635         int bw=size/2;
01636 
01637         unsigned long tsize=2*size;
01638         unsigned long ind1=0, ind2=0, ind3=0, ind4=0, ind41=0;
01639         unsigned long index0=0;
01640         int i=0, j=0, m=0, n=0, r=0;
01641         int tm=0, tm1=0, ipm=0, tipm=0, tipm1=0, loop_rho=0, rho_best=0;
01642 
01643         float* gnr2   = new float[size*2];
01644         float* maxcor = new float[size+1];                  // MAXR need change
01645 
01646         int p_max=p_max_input;
01647         float* result = new float[5*(p_max+1)];
01648         fftw_real * c = new fftw_real[size * size];   // ming, correlation of real space
01649         fftw_complex * C = new fftw_complex[size * (bw+1)]; // ming, correlation in fourier space which is complex number.
01650         rfftwnd_plan pinv;                           // Inverse 2D FFW plan; see FFTW menu p.7-11
01651         pinv = rfftw2d_create_plan(size, size, FFTW_COMPLEX_TO_REAL, FFTW_ESTIMATE);  //out-of-place translation
01652         float *self_sampl_fft = selfpcsfft->get_data(); // ming f(r)
01653         float maxcor_sofar=0.0f;
01654         int p=0;
01655 
01656         for(p=0; p<=p_max; p++){
01657                 ind1=p*size*bw;
01658                 int r_va1=1;
01659                 int r_va2=size;
01660                 for (i=0;i<size;i++)
01661                         for (j=0;j<bw+1;j++){
01662                                 C[i*(bw+1) + j].im = 0.0f;
01663                         }
01664 
01665         for(n=0;n<bw;n++){                                // loop for n
01666                 ind2 = (ind1+n);
01667                 index0=n*(bw+1);
01668                         for(r=r_va1;r<=r_va2;r++) {
01669                         ind3 = (ind2 + r*bw)*size;
01670                         for(m=0;m<size;m++){              // take back hat{h(n,r,p)}(m)
01671                                 ind4 = (ind3+m)*2;
01672                                     ind41= ind4+1;
01673                                     tm=m*2;
01674                                     tm1=tm+1;
01675                                     gnr2[tm]  = frm2dhhat[ind4];
01676                                     gnr2[tm1] = frm2dhhat[ind41];
01677                                 }
01678                         for(m=0;m<bw;m++){
01679                             tm=m*2; tm1=tm+1;
01680                                         ipm=index0+m; tipm=ipm*2; tipm1=tipm+1;               // index0=n*(bw+1)
01681                                 float tempr = self_sampl_fft[tm + r*(size+2)] * r;
01682                                 float tempi = self_sampl_fft[tm1+ r*(size+2)] * r;
01683                                 float gnr2_r = gnr2[tm];
01684                                 float gnr2_i = gnr2[tm1];
01685                                         C[n*(bw+1) + m].re += gnr2_r * tempr + gnr2_i * tempi;     //  m,+n
01686                                         C[n*(bw+1) + m].im += gnr2_i * tempr - gnr2_r * tempi;     // (a+ib)(c-id)=(ac+bd)+i(bc-ad)
01687                                         if(n != 0){                                     // m,-n
01688                                         if(m != 0){
01689                                                 int ssize = tsize-tm;   // ssize = 2*size-2m
01690                                                 int ssize1= ssize+1;
01691                                                 float gnr2_r = gnr2[ssize];
01692                                                 float gnr2_i = gnr2[ssize1];
01693                                                         C[(size-n)*(bw+1) + m].re += gnr2_r * tempr - gnr2_i * tempi;   // (a-ib)(c-id)=(ac-bd)-i(ad+bc)
01694                                                 C[(size-n)*(bw+1) + m].im -= gnr2_i * tempr + gnr2_r * tempi;// hhat*f(r)
01695                                         }
01696                                                 else{
01697                                                         C[(size-n)*(bw+1) + m].re += *(gnr2)  * tempr - *(gnr2+1)* tempi;
01698                                                         C[(size-n)*(bw+1) + m].im -= *(gnr2+1)* tempr + *(gnr2)  * tempi;
01699                                                 }
01700                                 }
01701                                 }
01702                         }
01703         }
01704 
01705                 rfftwnd_one_complex_to_real(pinv, C, c); // ming, this c has valuce at m n and rho
01706                 //c=C->do_ift();
01707                 float tempr=0.0f, corre_fcs=999.0f;
01708 
01709             int n_best=0, m_best=0;
01710                 for(n=0;n<size;n++){                    // move Tri_2D to Tri = c(phi,phi';rho)
01711                         float temp=0.0f;
01712                         for(m=0;m<size;m++){
01713                                 temp=c[n*size + m];
01714                                 if(temp>tempr) {
01715                                         tempr=temp;
01716                                         n_best=n;
01717                                         m_best=m;
01718                                 }
01719                         }
01720                 }
01721 
01722                 float corre,Phi2,Phi,Tx,Ty,Vx, Vy;;
01723                 //EMData *this_tmp;
01724         //      bw=size;
01725                 //for (n_best=0;n_best<bw;n_best++)
01726                         //for (m_best=0;m_best<2*bw;m_best++){
01727                 //n_best=0;
01728                 //m_best=70;
01729                 Phi2= n_best*M_PI/bw;  // ming this is reference image rotation angle
01730                 Phi = m_best*M_PI/bw;   // ming this is particle image rotation angle
01731                 Vx =  p*cos(Phi);//should use the angle of the centered one
01732                 Vy = -p*sin(Phi);
01733 
01734                 Tx=Vx+(floor(com_this_x+0.5) - floor(com_with_x+0.5));
01735                 Ty=Vy+(floor(com_this_y+0.5) - floor(com_with_y+0.5));
01736 
01737                 dx = -Tx;       // the Rota & Trans value (Tx,Ty, ang_keep) are for the projection image,
01738                 dy = -Ty;       // need to convert to raw image
01739 
01740                 EMData *this_tmp=this_img->copy();//ming change to to
01741                 this_tmp->rotate(-(Phi2-Phi)*180.0,0.0,0.0);
01742                 this_tmp->translate(dx,dy,0.0);
01743 
01744                 corre=this_tmp->cmp(cmp_name,to,cmp_params);
01745                 printf("corre=%f n_best=%d m_best=%d p=%d\n",corre,n_best,m_best,p);
01746                         //}
01747                 delete this_tmp;
01748                 //float corre_fcs=0.0f;
01749                 if(corre<=corre_fcs) { //ming, cmp use smaller value stands for more similarity
01750                         corre_fcs=corre;
01751                         result[0+5*p] = float(p);       // rho
01752                         result[1+5*p] = corre;          // correlation_fcs
01753                         result[2+5*p] = (Phi2-Phi)*180.0;       // rotation angle
01754                         result[3+5*p] = Tx;             // Translation_x
01755                         result[4+5*p] = Ty;             // Translation_y
01756                 }
01757                 maxcor[p]=corre_fcs;                            //  maximum correlation for current rho
01758                 if(corre_fcs<maxcor_sofar) {
01759                         maxcor_sofar=corre_fcs;                 // max correlation up to current rho
01760                     rho_best=p;                         // the rho value with maxinum correlation value
01761                 }
01762                 if(p>=4){
01763                         if(maxcor[p] < maxcor[p-1] && maxcor[p-1] < maxcor[p-2]&& maxcor[p-2] < maxcor[p-3] && maxcor[p-3] < maxcor[p-4]){
01764                                 loop_rho=1;
01765                                 break; //exit p loop
01766                         }
01767                 }
01768         }
01769         //}
01770         if(loop_rho == 1)
01771                 p=p+1;
01772         int rb5=5*rho_best;
01773         float fsc      = result[1+rb5];
01774         float ang_keep = result[2+rb5];
01775         //float ang_keep=150.;
01776         float Tx       = result[3+rb5];
01777         float Ty       = result[4+rb5];
01778         delete[] gnr2;
01779         delete[] maxcor;
01780         delete[] result;
01781         delete c;
01782         c = 0;
01783         delete C;
01784         C = 0;
01785         rfftwnd_destroy_plan(pinv);
01786         dx = -Tx;               // the Rota & Trans value (Tx,Ty, ang_keep) are for the projection image,
01787         dy = -Ty;               // need to convert to raw image
01788         this_img->rotate(-ang_keep,0,0); // ming change this to this_img??
01789         this_img->translate(dx,dy,0.0); // ming change this to this_img
01790 
01791         float fsc_best=this_img->cmp(cmp_name,to,cmp_params);
01792         printf("rho_best=%d fsc=%f fsc_best=%f dx=%f dy=%f ang_keep=%f com_withx=%f com_selfx=%f com_selfy=%f\n",rho_best,fsc,fsc_best,dx,dy,ang_keep,com_with_x,com_this_x,com_this_y);
01793         return fsc;     // return the fsc coefficients
01794 } // FRM2D aligner sub_class
01795         }
01796 
01797 
01798 EMData *FRM2DAligner::align(EMData * this_img, EMData * to,
01799                         const string & cmp_name, const Dict& cmp_params) const
01800 {
01801         if (!this_img) {
01802                 return 0;
01803         }
01804         if (to && !EMUtil::is_same_size(this_img, to))
01805                 throw ImageDimensionException("Images must be the same size to perform translational alignment");
01806 
01807         int nx=this_img->get_xsize();
01808         int ny=this_img->get_ysize();
01809         int size = (int) floor(M_PI*ny/4.0);
01810         size =Util::calc_best_fft_size(size);//ming   bestfftsize(size);
01811         int MAXR=ny/2;
01812 
01813         EMData *this_temp=this_img->copy(); // ming change avg to to
01814         FloatPoint com_test;
01815         com_test=this_temp->calc_center_of_mass();//ming add
01816         float com_this_x=com_test[0];
01817         float com_this_y=com_test[1];
01818         delete this_temp;
01819         this_temp=to->copy();
01820         com_test=this_temp->calc_center_of_mass();
01821         float com_with_x=com_test[0];
01822         float com_with_y=com_test[1];
01823         delete this_temp;
01824 
01825 
01826         float dx=-(com_with_x - nx/2); //ming
01827         float dy=-(com_with_y - ny/2); //ming
01828         //avg_frm->translate(-dx,-dy,0.0);
01829         EMData *avg_frm=to->copy();
01830         EMData *withpcs=avg_frm->unwrap_largerR(0,MAXR,size,float(MAXR)); // ming, something wrong inside this subroutine
01831         EMData *withpcsfft=withpcs->oneDfftPolar(size, MAXR, MAXR);
01832         //EMData *withpcsfft=withpcs->do_fft();
01833         //withpcsfft->ap2ri();
01834         float *sampl_fft = withpcsfft->get_data(); // ming change getDataRO() to get_data
01835         delete avg_frm;
01836         delete withpcs;
01837 
01838         int bw=size/2.;
01839         unsigned long ind1=0, ind2=0, ind3=0, ind4=0, ind41=0;
01840         float pi2=2.0*M_PI, r2;
01841         fftw_complex * gnr2_in = new fftw_complex[size];
01842         fftw_complex * gnr2 = new fftw_complex[size];
01843         fftw_plan planc_fftw;                   // plan for 1D FFTW (Forward), complex->complex  to get gnr2
01844         planc_fftw = fftw_create_plan(size, FFTW_FORWARD, FFTW_ESTIMATE );    // out-of-plalce transform, see FFTW menu p.3-4
01845         int p_max=params.set_default("p_max",5);
01846         float *frm2dhhat=0;
01847         if( (frm2dhhat = (float *) malloc( (p_max+1)*(size+2)*bw*size*2* sizeof(float)) ) == NULL ) {
01848                 cout <<"Error in allocating memory 13. \n";
01849                 exit(1);
01850         }
01851         float *sb=0, *cb=0;             // sin(beta) and cos(beta) for get h_hat, 300>size
01852         if( (sb = new float[size]) == NULL || (cb = new float[size]) == NULL) {
01853                 cout <<"can't allocate more memory, terminating. \n";
01854                 exit(1);
01855         }
01856         for(int i=0;i<size;i++) {        // beta sampling, to calculate beta' and r'
01857                 float beta=i*M_PI/bw;
01858                 sb[i]=sin(beta);
01859                 cb[i]=cos(beta);
01860         }
01861         for(int p=0; p<=p_max; p++){
01862                 ind1=p*size*bw;
01863         float pp2 = p*p;
01864                 int r_va1 = 1;
01865                 int r_va2 = size;
01866 
01867                 for(int n=0;n<bw;n++){         /* loop for n */
01868                 int tn=n*2;
01869                 int tn1=tn+1;
01870                 ind2 = ind1+n;
01871                 for(int r=r_va1;r<=r_va2;r++) {
01872                                 ind3 = (ind2+r*bw)*size;
01873                         float rr2 = r*r;
01874                                 float rp2 = r*p;
01875                         for(int i=0;i<size;i++){                            // beta sampling, to get beta' and r'
01876                                 if(p==0) r2 = (float) r;
01877                                 else     r2 = std::sqrt((float)(rr2+pp2-2.0*rp2*cb[i]));   // r2->r'
01878                                 int r1=floor(r2+0.5);                        // for computing gn(r')
01879                                 if(r1>MAXR){
01880                                         gnr2_in[i].re  = 0.0f;
01881                                         gnr2_in[i].im  = 0.0f;
01882                                 }
01883                                 else{
01884                                         float gn_r = sampl_fft[tn +r1*(size+2)];           // real part of gn(r')
01885                                         float gn_i = sampl_fft[tn1+r1*(size+2)];           // imaginary part of gn(r')
01886                                                 float sb2, cb2, cn, sn;
01887                                                 if(n != 0){
01888                                                         if(r2 != 0.0){
01889                                                                 if(p==0) {
01890                                                                         sb2=sb[i];
01891                                                                         cb2=cb[i];
01892                                                                 }
01893                                                                 else{
01894                                                                         sb2=r*sb[i]/r2;
01895                                                                         cb2=(r*cb[i]-p)/r2;
01896                                                                 }
01897                                                         }
01898                                                 else{
01899                                                                 sb2=0.0;
01900                                                                 cb2=1.0;
01901                                                         }
01902                                                 if(sb2>1.0) sb2= 1.0f;
01903                                                 if(sb2<-1.0)sb2=-1.0f;
01904                                                 if(cb2>1.0) cb2= 1.0f;
01905                                                 if(cb2<-1.0)cb2=-1.0f;
01906                                                 float beta2=atan2(sb2,cb2);
01907                                                 if(beta2<0.0) beta2 += pi2;
01908                                                 float nb2=n*beta2;
01909                                                 cn=cos(nb2);
01910                                                         sn=sin(nb2);
01911                                                 }
01912                                         else{
01913                                                         cn=1.0f; sn=0.0f;
01914                                                 }
01915                                                 gnr2_in[i].re =   cn*gn_r - sn*gn_i;         // exp(-i*n*beta')*gn(r')
01916                                                 gnr2_in[i].im = -(cn*gn_i + sn*gn_r);
01917                                 }
01918                         }
01919 
01920                         fftw_one(planc_fftw, gnr2_in, gnr2); //complex to complex
01921                                 for(int m=0;m<size;m++){                                     // store hat{h(n,r,p)}(m)
01922                                         ind4 = (ind3+m)*2;
01923                                         ind41 = ind4+1;
01924                                         frm2dhhat[ind4]  = gnr2[m].re;
01925                                         frm2dhhat[ind41] = gnr2[m].im;
01926                                 }
01927                         }
01928                 }
01929         }
01930         delete[] sb;
01931         delete[] cb;
01932         delete [] gnr2_in;
01933         gnr2_in = 0;
01934         delete [] gnr2;
01935         gnr2 = 0;
01936         fftw_destroy_plan(planc_fftw);
01937         delete withpcsfft;
01938 
01939         float dot_frm0=0.0f, dot_frm1=0.0f;
01940         EMData *da_nFlip=0, *da_yFlip=0;
01941         for (int iFlip=0;iFlip<=1;iFlip++){
01942                 EMData *dr_frm=0;
01943                 if (iFlip==0) {
01944                         dr_frm=this_img->copy();
01945                         da_nFlip=this_img->copy();
01946                 }   //ming, this actually refer to this_img  //delete dr_frm somewhere
01947                 else {
01948                         dr_frm=this_img->copy();
01949                         da_yFlip=this_img->copy();
01950                 }
01951                 if(iFlip==1) com_this_x=nx-com_this_x;   //ming   // image mirror about Y axis, so y keeps the same
01952 
01953                 dx=-(com_this_x - nx/2); //ming
01954                 dy=-(com_this_y - ny/2); //ming
01955                 dr_frm->translate(dx,dy,0.0); // this
01956                 EMData *selfpcs = dr_frm->unwrap_largerR(0,MAXR,size, (float)MAXR);
01957                 EMData *selfpcsfft = selfpcs->oneDfftPolar(size, MAXR, MAXR);
01958                 //EMData *selfpcsfft = selfpcs->do_fft();       // 1DFFT of polar sampling
01959                 //selfpcsfft->ap2ri();
01960                 delete dr_frm;
01961                 delete selfpcs;
01962                 if(iFlip==0)
01963                         dot_frm0=frm_2d_Align(da_nFlip,to, frm2dhhat, selfpcsfft, p_max, size, com_this_x, com_this_y, com_with_x, com_with_y,cmp_name,cmp_params);
01964                 else
01965                         dot_frm1=frm_2d_Align(da_yFlip,to, frm2dhhat, selfpcsfft, p_max, size, com_this_x, com_this_y, com_with_x, com_with_y,cmp_name,cmp_params);
01966                 delete selfpcsfft;
01967         }
01968         if(dot_frm0 <= dot_frm1) {
01969                 delete da_yFlip;
01970                 printf("best_corre=%f, no flip\n",dot_frm0);
01971                 return da_nFlip;
01972         }
01973         else {
01974                 delete da_nFlip;
01975                 printf("best_corre=%f, flipped\n",dot_frm1);
01976                 return da_yFlip;
01977         }
01978 }
01979 #endif  //FFTW2
01980 
01981 CUDA_Aligner::CUDA_Aligner() {
01982         image_stack = NULL;
01983         image_stack_filtered = NULL;
01984         ccf = NULL;
01985 }
01986 
01987 #ifdef EMAN2_USING_CUDA
01988 void CUDA_Aligner::finish() {
01989         if (image_stack) delete image_stack;
01990         if (image_stack_filtered) delete image_stack_filtered;
01991         if (ccf) delete ccf;
01992 }
01993 
01994 void CUDA_Aligner::setup(int nima, int nx, int ny, int ring_length, int nring, int ou, float step, int kx, int ky, bool ctf) {
01995 
01996         NIMA = nima;
01997         NX = nx;
01998         NY = ny;
01999         RING_LENGTH = ring_length;
02000         NRING = nring;
02001         STEP = step;
02002         KX = kx;
02003         KY = ky;
02004         OU = ou;
02005         CTF = ctf;
02006         
02007         image_stack = (float *)malloc(NIMA*NX*NY*sizeof(float));
02008         if (CTF == 1) image_stack_filtered = (float *)malloc(NIMA*NX*NY*sizeof(float));
02009         ccf = (float *)malloc(2*(2*KX+1)*(2*KY+1)*NIMA*(RING_LENGTH+2)*sizeof(float));
02010 }
02011 
02012 void CUDA_Aligner::insert_image(EMData *image, int num) {
02013 
02014         int base_address = num*NX*NY;
02015 
02016         for (int y=0; y<NY; y++)
02017                 for (int x=0; x<NX; x++)
02018                         image_stack[base_address+y*NX+x] = (*image)(x, y);
02019 }
02020 
02021 void CUDA_Aligner::filter_stack(vector<float> ctf_params, int id) {
02022         
02023         float *params;
02024         
02025         params = (float *)malloc(NIMA*6*sizeof(float)); 
02026         
02027         for (int i=0; i<NIMA*6; i++) params[i] = ctf_params[i];
02028 
02029         filter_image(image_stack, image_stack_filtered, NIMA, NX, NY, params, id);
02030 
02031         free(params);
02032 }
02033 
02034 void CUDA_Aligner::sum_oe(vector<float> ctf_params, vector<float> ali_params, EMData *ave1, EMData *ave2, int id) {
02035         
02036         float *ctf_p, *ali_p, *av1, *av2;
02037         
02038         ctf_p = (float *)malloc(NIMA*6*sizeof(float));
02039         ali_p = (float *)malloc(NIMA*4*sizeof(float));
02040         
02041         if (CTF == 1) {
02042                 for (int i=0; i<NIMA*6; i++)  ctf_p[i] = ctf_params[i];
02043         }
02044         for (int i=0; i<NIMA*4; i++)   ali_p[i] = ali_params[i];
02045         
02046         av1 = ave1->get_data();
02047         av2 = ave2->get_data();
02048         
02049         rot_filt_sum(image_stack, NIMA, NX, NY, CTF, ctf_p, ali_p, av1, av2, id);
02050         
02051         free(ctf_p);
02052         free(ali_p);
02053 }
02054 
02055 vector<float> CUDA_Aligner::alignment_2d(EMData *ref_image_em, vector<float> sx_list, vector<float> sy_list, int id, int silent) {
02056 
02057         float *ref_image, max_ccf;
02058         int base_address, ccf_offset;
02059         float ts, tm;
02060         float ang, sx = 0, sy = 0, mirror, co, so, sxs, sys;
02061         float *sx2, *sy2;
02062         vector<float> align_result;
02063 
02064         sx2 = (float *)malloc(NIMA*sizeof(float));
02065         sy2 = (float *)malloc(NIMA*sizeof(float));
02066 
02067         ref_image = ref_image_em->get_data();
02068         
02069         for (int i=0; i<NIMA; i++) {
02070                 sx2[i] = sx_list[i];
02071                 sy2[i] = sy_list[i];
02072         }
02073         
02074         if (CTF == 1) {
02075                 calculate_ccf(image_stack_filtered, ref_image, ccf, NIMA, NX, NY, RING_LENGTH, NRING, OU, STEP, KX, KY, sx2, sy2, id, silent);
02076         } else {
02077                 calculate_ccf(image_stack, ref_image, ccf, NIMA, NX, NY, RING_LENGTH, NRING, OU, STEP, KX, KY, sx2, sy2, id, silent);
02078         }
02079 
02080         ccf_offset = NIMA*(RING_LENGTH+2)*(2*KX+1)*(2*KY+1);
02081 
02082         for (int im=0; im<NIMA; im++) {
02083                 max_ccf = -1.0e22;
02084                 for (int kx=-KX; kx<=KX; kx++) {
02085                         for (int ky=-KY; ky<=KY; ky++) {
02086                                 base_address = (((ky+KY)*(2*KX+1)+(kx+KX))*NIMA+im)*(RING_LENGTH+2);
02087                                 for (int l=0; l<RING_LENGTH; l++) {
02088                                         ts = ccf[base_address+l];
02089                                         tm = ccf[base_address+l+ccf_offset];
02090                                         if (ts > max_ccf) {
02091                                                 ang = float(l)/RING_LENGTH*360.0;
02092                                                 sx = -kx*STEP;
02093                                                 sy = -ky*STEP;
02094                                                 mirror = 0;
02095                                                 max_ccf = ts;
02096                                         }
02097                                         if (tm > max_ccf) {
02098                                                 ang = float(l)/RING_LENGTH*360.0; 
02099                                                 sx = -kx*STEP;
02100                                                 sy = -ky*STEP;
02101                                                 mirror = 1;
02102                                                 max_ccf = tm;
02103                                         }
02104                                 }
02105                         }
02106                 }
02107                 co =  cos(ang*M_PI/180.0);
02108                 so = -sin(ang*M_PI/180.0);
02109                 sxs = sx*co - sy*so;
02110                 sys = sx*so + sy*co;
02111 
02112                 align_result.push_back(ang);
02113                 align_result.push_back(sxs);
02114                 align_result.push_back(sys);
02115                 align_result.push_back(mirror);
02116         }
02117         
02118         free(sx2);
02119         free(sy2);
02120         
02121         return align_result;
02122 }
02123 
02124 
02125 vector<float> CUDA_Aligner::ali2d_single_iter(EMData *ref_image_em, vector<float> ali_params, float csx, float csy, int id, int silent, float delta) {
02126 
02127         float *ref_image, max_ccf;
02128         int base_address, ccf_offset;
02129         float ts, tm;
02130         float ang, sx, sy, mirror, co, so, sxs, sys;
02131         float *sx2, *sy2;
02132         vector<float> align_result;
02133 
02134         sx2 = (float *)malloc(NIMA*sizeof(float));
02135         sy2 = (float *)malloc(NIMA*sizeof(float));
02136 
02137         ref_image = ref_image_em->get_data();
02138         
02139         for (int i=0; i<NIMA; i++) {
02140                 ang = ali_params[i*4]/180.0*M_PI;
02141                 sx = (ali_params[i*4+3] < 0.5)?(ali_params[i*4+1]-csx):(ali_params[i*4+1]+csx);
02142                 sy = ali_params[i*4+2]-csy;
02143                 co = cos(ang);
02144                 so = sin(ang);
02145                 sx2[i] = -(sx*co-sy*so);
02146                 sy2[i] = -(sx*so+sy*co);
02147         }
02148         
02149         if (CTF == 1) {
02150                 calculate_ccf(image_stack_filtered, ref_image, ccf, NIMA, NX, NY, RING_LENGTH, NRING, OU, STEP, KX, KY, sx2, sy2, id, silent);
02151         } else {
02152                 calculate_ccf(image_stack, ref_image, ccf, NIMA, NX, NY, RING_LENGTH, NRING, OU, STEP, KX, KY, sx2, sy2, id, silent);
02153         }
02154 
02155         ccf_offset = NIMA*(RING_LENGTH+2)*(2*KX+1)*(2*KY+1);
02156 
02157         float sx_sum = 0.0f;
02158         float sy_sum = 0.0f;
02159 
02160         int dl;
02161         dl = static_cast<int>(delta/360.0*RING_LENGTH);
02162         if (dl<1) { dl = 1; }   
02163         
02164         for (int im=0; im<NIMA; im++) {
02165                 max_ccf = -1.0e22;
02166                 for (int kx=-KX; kx<=KX; kx++) {
02167                         for (int ky=-KY; ky<=KY; ky++) {
02168                                 base_address = (((ky+KY)*(2*KX+1)+(kx+KX))*NIMA+im)*(RING_LENGTH+2);
02169                                 for (int l=0; l<RING_LENGTH; l+=dl) {
02170                                         ts = ccf[base_address+l];
02171                                         tm = ccf[base_address+l+ccf_offset];
02172                                         if (ts > max_ccf) {
02173                                                 ang = float(l)/RING_LENGTH*360.0;
02174                                                 sx = -kx*STEP;
02175                                                 sy = -ky*STEP;
02176                                                 mirror = 0;
02177                                                 max_ccf = ts;
02178                                         }
02179                                         if (tm > max_ccf) {
02180                                                 ang = float(l)/RING_LENGTH*360.0; 
02181                                                 sx = -kx*STEP;
02182                                                 sy = -ky*STEP;
02183                                                 mirror = 1;
02184                                                 max_ccf = tm;
02185                                         }
02186                                 }
02187                         }
02188                 }
02189                 co =  cos(ang*M_PI/180.0);
02190                 so = -sin(ang*M_PI/180.0);
02191                 
02192                 sxs = (sx-sx2[im])*co-(sy-sy2[im])*so;
02193                 sys = (sx-sx2[im])*so+(sy-sy2[im])*co;
02194 
02195                 align_result.push_back(ang);
02196                 align_result.push_back(sxs);
02197                 align_result.push_back(sys);
02198                 align_result.push_back(mirror);
02199                 
02200                 if (mirror == 0)  { sx_sum += sxs; }  else { sx_sum -= sxs; }
02201                 sy_sum += sys;
02202         }
02203         
02204         align_result.push_back(sx_sum);
02205         align_result.push_back(sy_sum);
02206         
02207         free(sx2);
02208         free(sy2);
02209         
02210         return align_result;
02211 }
02212 
02213 
02214 #endif
02215 
02216 void EMAN::dump_aligners()
02217 {
02218         dump_factory < Aligner > ();
02219 }
02220 
02221 map<string, vector<string> > EMAN::dump_aligners_list()
02222 {
02223         return dump_factory_list < Aligner > ();
02224 }

Generated on Mon Jul 19 12:40:09 2010 for EMAN2 by  doxygen 1.4.7