aligner.cpp File Reference

#include "emfft.h"
#include "cmp.h"
#include "aligner.h"
#include "emdata.h"
#include "processor.h"
#include "util.h"
#include "symmetry.h"
#include <gsl/gsl_multimin.h>
#include "plugins/aligner_template.h"

Include dependency graph for aligner.cpp:

Go to the source code of this file.

Defines

#define EMAN2_ALIGNER_DEBUG   0
 
Id
aligner.cpp,v 1.247 2011/06/13 16:39:34 john Exp


Functions

static double refalifn (const gsl_vector *v, void *params)
static double refalifnfast (const gsl_vector *v, void *params)
static Transform refalin3d_perturbquat (const Transform *const t, const float &spincoeff, const float &n0, const float &n1, const float &n2, const float &x, const float &y, const float &z)
static double refalifn3dquat (const gsl_vector *v, void *params)
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)
void EMAN::dump_aligners ()
map< string, vector< string > > EMAN::dump_aligners_list ()


Define Documentation

#define EMAN2_ALIGNER_DEBUG   0

Id
aligner.cpp,v 1.247 2011/06/13 16:39:34 john Exp

Definition at line 54 of file aligner.cpp.


Function Documentation

float @0::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 
) [static]

Definition at line 2310 of file aligner.cpp.

References EMAN::EMData::cmp(), EMAN::EMData::copy(), EMAN::EMData::do_ift(), EMAN::EMData::get_data(), EMAN::EMData::get_ysize(), in, EMAN::EMData::rotate(), EMAN::EMData::set_attr(), EMAN::Transform::set_trans(), and EMAN::EMData::translate().

Referenced by EMAN::FRM2DAligner::align().

02311 {
02312         int size=rsize;
02313         float dx,dy;
02314         int bw=size/2;
02315         int MAXR=this_img->get_ysize()/2;
02316         //int MAXR=size;
02317         unsigned long tsize=2*size;
02318         unsigned long ind1=0, ind2=0, ind3=0, ind4=0, ind41=0;
02319         unsigned long index0=0;
02320         int i=0, j=0, m=0, n=0, r=0;
02321         int loop_rho=0, rho_best=0;
02322 
02323         float* gnr2   = new float[size*2];
02324         float* maxcor = new float[size+1];                  // MAXR need change
02325 
02326         int p_max=p_max_input;
02327         float* result = new float[5*(p_max+1)];
02328         float* cr=new float[size*(bw+1)];
02329         float* ci=new float[size*(bw+1)];
02330         EMData *data_in=new EMData;
02331         data_in->set_complex(true);
02332         data_in->set_fftodd(false);
02333         data_in->set_ri(true);
02334         data_in->set_size(size+2,size,1);
02335         float *in=data_in->get_data();
02336 
02337         float *self_sampl_fft = selfpcsfft->get_data(); // ming f(r)
02338 
02339         float maxcor_sofar=0.0f;
02340         int p=0;
02341 
02342         for(p=0; p<=p_max; ++p){
02343                 ind1=p*size*bw;
02344                 for (i=0;i<size;++i)
02345                         for (j=0;j<bw+1;++j){
02346                                 cr[i*(bw+1)+j]=0.0;
02347                                 ci[i*(bw+1)+j]=0.0;
02348                         }
02349         for(n=0;n<bw;++n){                                // loop for n
02350                 ind2=(ind1+n);
02351                 index0=n*(bw+1);
02352                         for(r=0;r<=MAXR;++r) {
02353                         ind3=(ind2+r*bw)*size;
02354                         for(m=0;m<size;m++){              // take back hat{h(n,r,p)}(m)
02355                                 ind4=(ind3+m)*2;
02356                                     ind41=ind4+1;
02357                                     gnr2[2*m]=frm2dhhat[ind4];
02358                                     gnr2[2*m+1]=frm2dhhat[ind41];
02359                                 }
02360                         for(m=0;m<bw;++m){
02361                                         float tempr=self_sampl_fft[2*m+r*(size+2)]*r;
02362                                 float tempi=self_sampl_fft[2*m+1+r*(size+2)]*r;
02363                                 float gnr2_r=gnr2[2*m];
02364                                 float gnr2_i=gnr2[2*m+1];
02365                                 cr[n*(bw+1)+m]+=gnr2_r*tempr+gnr2_i*tempi;
02366                                         ci[n*(bw+1)+m]+=gnr2_i*tempr-gnr2_r*tempi;
02367                                         if(n!=0){                                       // m,-n
02368                                         if(m!= 0){
02369                                                 int ssize=tsize-2*m;    // ssize = 2*size-2m
02370                                                 int ssize1=ssize+1;
02371                                                 float gnr2_r=gnr2[ssize];
02372                                                 float gnr2_i=gnr2[ssize1];
02373                                                         cr[(size-n)*(bw+1)+m]+=gnr2_r*tempr-gnr2_i*tempi;
02374                                                 ci[(size-n)*(bw+1)+m]-=gnr2_i*tempr+gnr2_r*tempi;
02375                                         }
02376                                                 else{
02377                                                         cr[(size-n)*(bw+1)+m]+=*(gnr2)*tempr-*(gnr2+1)*tempi;
02378                                                         ci[(size-n)*(bw+1)+m]-=*(gnr2+1)*tempr+*(gnr2)*tempi;
02379                                                 }
02380                                 }
02381                                 }
02382                         }
02383         }
02384         for (int cii=0; cii<size*(bw+1);++cii){
02385                         in[2*cii]=cr[cii];
02386                         in[2*cii+1]=ci[cii];
02387                         //printf("cii=%d,in[2i+1]=%f\n",cii, cr[cii]);
02388         }
02389 
02390         EMData *data_out;
02391                 data_out=data_in->do_ift();
02392                 float *c=data_out->get_data();
02393                 float tempr=0.0f, corre_fcs=999.0f;
02394 
02395             int n_best=0, m_best=0;
02396         float temp=-100.0f;
02397                 for(n=0;n<size;++n){// move Tri_2D to Tri = c(phi,phi';rho)
02398                         for(m=0;m<size;++m){
02399                                 temp=c[n*size+m];
02400                                 if(temp>tempr) {
02401                                         tempr=temp;
02402                                         n_best=n;
02403                                         m_best=m;
02404                                 }
02405                         }
02406                 }
02407                 delete data_out;
02408 
02409                 float corre,Phi2,Phi,Tx,Ty,Vx, Vy;
02410 
02411                 //for (n_best=0;n_best<bw;n_best++)
02412                   //  for (m_best=0;m_best<2*bw;m_best++){
02413                 //n_best=0;
02414                 //m_best=70;
02415                 Phi2=n_best*M_PI/bw;  // ming this is reference image rotation angle
02416                 Phi=m_best*M_PI/bw;   // ming this is particle image rotation angle
02417                 Vx=p*cos(Phi);//should use the angle of the centered one
02418                 Vy=-p*sin(Phi);
02419                 Tx=Vx+(floor(com_this_x+0.5f)-floor(com_with_x+0.5f));
02420                 Ty=Vy+(floor(com_this_y+0.5f)-floor(com_with_y+0.5f));
02421 
02422                 dx=-Tx; // the Rota & Trans value (Tx,Ty, ang_keep) are for the projection image,
02423                 dy=-Ty; // need to convert to raw image
02424 
02425                 EMData *this_tmp=this_img->copy();//ming change to to
02426                 this_tmp->rotate(-(Phi2-Phi)*180.0f,0.0f,0.0f);
02427                 this_tmp->translate(dx,dy,0.0);
02428 
02429                 corre=this_tmp->cmp(cmp_name,to,cmp_params);
02430                 //printf("corre=%f\n",corre);
02431                 delete this_tmp;
02432                 if(corre<=corre_fcs) { //ming, cmp use smaller value stands for more similarity
02433                         corre_fcs=corre;
02434                         result[0+5*p] = float(p);       // rho
02435                         result[1+5*p] = corre;          // correlation_fcs
02436                         result[2+5*p] = (Phi2-Phi)*180.0f;      // rotation angle
02437                         result[3+5*p] = Tx;             // Translation_x
02438                         result[4+5*p] = Ty;             // Translation_y
02439                 }
02440                 maxcor[p]=corre_fcs;                            //  maximum correlation for current rho
02441                 if(corre_fcs<maxcor_sofar) {
02442                         maxcor_sofar=corre_fcs;                 // max correlation up to current rho
02443                     rho_best=p;                         // the rho value with maxinum correlation value
02444                 }
02445                 if(p>=4){
02446                         if(maxcor[p] < maxcor[p-1] && maxcor[p-1] < maxcor[p-2]&& maxcor[p-2] < maxcor[p-3] && maxcor[p-3] < maxcor[p-4]){
02447                                 loop_rho=1;
02448                                 break; //exit p loop
02449                         }
02450                 }
02451         } // end for p
02452         //}//test my method
02453         if(loop_rho == 1)
02454                 p=p+1;
02455         int rb5=5*rho_best;
02456         float fsc      = result[1+rb5];
02457         float ang_keep = result[2+rb5];
02458         float Tx       = result[3+rb5];
02459         float Ty       = result[4+rb5];
02460         delete[] gnr2;
02461         delete[] maxcor;
02462         delete[] result;
02463         delete cr;
02464         cr=0;
02465         delete ci;
02466         ci=0;
02467         delete data_in; // ming add
02468         dx = -Tx;               // the Rota & Trans value (Tx,Ty, ang_keep) are for the projection image,
02469         dy = -Ty;               // need to convert to raw image
02470         this_img->rotate(-ang_keep,0,0); // ming change this to this_img??
02471         this_img->translate(dx,dy,0.0); // ming change this to this_img
02472 
02473 
02474         Transform  tsoln(Dict("type","2d","alpha",ang_keep));
02475         tsoln.set_trans(dx,dy);
02476         this_img->set_attr("xform.align2d",&tsoln);
02477 #ifdef DEBUG
02478         float fsc_best=this_img->cmp(cmp_name,to,cmp_params);
02479         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);
02480 #endif
02481         return fsc;     // return the fsc coefficients
02482 } // FRM2D aligner sub_class

static double refalifn ( const gsl_vector *  v,
void *  params 
) [static]

Definition at line 1323 of file aligner.cpp.

References EMAN::Cmp::cmp(), EMAN::EMData::process(), t, x, and y.

Referenced by EMAN::RefineAligner::align().

01324 {
01325         Dict *dict = (Dict *) params;
01326 
01327         double x = gsl_vector_get(v, 0);
01328         double y = gsl_vector_get(v, 1);
01329         double a = gsl_vector_get(v, 2);
01330 
01331         EMData *this_img = (*dict)["this"];
01332         EMData *with = (*dict)["with"];
01333         bool mirror = (*dict)["mirror"];
01334 
01335 //      float mean = (float)this_img->get_attr("mean");
01336 //      if ( Util::goodf(&mean) ) {
01337 //              //cout << "tmps mean is nan even before rotation" << endl;
01338 //      }
01339 
01340         Transform t(Dict("type","2d","alpha",static_cast<float>(a)));
01341 //      Transform3D t3d(Transform3D::EMAN, (float)a, 0.0f, 0.0f);
01342 //      t3d.set_posttrans( (float) x, (float) y);
01343 //      tmp->rotate_translate(t3d);
01344         t.set_trans((float)x,(float)y);
01345         t.set_mirror(mirror);
01346         if (v->size>3) {
01347                 float sca=(float)gsl_vector_get(v, 3);
01348                 if (sca<.7 || sca>1.3) return 1.0e20;
01349                 t.set_scale((float)gsl_vector_get(v, 3));
01350         }
01351         EMData *tmp = this_img->process("xform",Dict("transform",&t));
01352 
01353 //      printf("GSL %f %f %f %d %f\n",x,y,a,mirror,(float)gsl_vector_get(v, 3));
01354         Cmp* c = (Cmp*) ((void*)(*dict)["cmp"]);
01355         double result = c->cmp(tmp,with);
01356 
01357         // DELETE AT SOME STAGE, USEFUL FOR PRERELEASE STUFF
01358         //      float test_result = (float)result;
01359 //      if ( Util::goodf(&test_result) ) {
01360 //              cout << "result " << result << " " << x << " " << y << " " << a << endl;
01361 //              cout << (float)this_img->get_attr("mean") << " " << (float)tmp->get_attr("mean") << " " << (float)with->get_attr("mean") << endl;
01362 //              tmp->write_image("tmp.hdf");
01363 //              with->write_image("with.hdf");
01364 //              this_img->write_image("this_img.hdf");
01365 //              EMData* t = this_img->copy();
01366 //              cout << (float)t->get_attr("mean") << endl;
01367 //              t->rotate_translate( t3d );
01368 //              cout << (float)t->get_attr("mean") << endl;
01369 //              cout << "exit" << endl;
01371 //              cout << (float)t->get_attr("mean") << endl;
01372 //              cout << "now exit" << endl;
01373 //              delete t;
01374 //      }
01375 
01376 
01377         if ( tmp != 0 ) delete tmp;
01378 
01379         return result;
01380 }

static double refalifn3dquat ( const gsl_vector *  v,
void *  params 
) [static]

Definition at line 1558 of file aligner.cpp.

References EMAN::Cmp::cmp(), EMAN::EMData::process(), refalin3d_perturbquat(), t, x, and y.

Referenced by EMAN::Refine3DAlignerQuaternion::align().

01559 {
01560         Dict *dict = (Dict *) params;
01561 
01562         double n0 = gsl_vector_get(v, 0);
01563         double n1 = gsl_vector_get(v, 1);
01564         double n2 = gsl_vector_get(v, 2);
01565         double x = gsl_vector_get(v, 3);
01566         double y = gsl_vector_get(v, 4);
01567         double z = gsl_vector_get(v, 5);
01568 
01569         EMData *this_img = (*dict)["this"];
01570         EMData *with = (*dict)["with"];
01571 //      bool mirror = (*dict)["mirror"];
01572 
01573         Transform* t = (*dict)["transform"];
01574         float spincoeff = (*dict)["spincoeff"];
01575 
01576         Transform soln = refalin3d_perturbquat(t,spincoeff,(float)n0,(float)n1,(float)n2,(float)x,(float)y,(float)z);
01577 
01578         EMData *tmp = this_img->process("xform",Dict("transform",&soln));
01579         Cmp* c = (Cmp*) ((void*)(*dict)["cmp"]);
01580         double result = c->cmp(tmp,with);
01581         if ( tmp != 0 ) delete tmp;
01582         delete t; t = 0;
01583         //cout << result << endl;
01584         return result;
01585 }

static double refalifnfast ( const gsl_vector *  v,
void *  params 
) [static]

Definition at line 1382 of file aligner.cpp.

References EMAN::EMData::dot_rotate_translate(), EMAN::EMData::get_xsize(), EMAN::EMData::get_ysize(), x, and y.

Referenced by EMAN::RefineAligner::align().

01383 {
01384         Dict *dict = (Dict *) params;
01385         EMData *this_img = (*dict)["this"];
01386         EMData *img_to = (*dict)["with"];
01387         bool mirror = (*dict)["mirror"];
01388 
01389         double x = gsl_vector_get(v, 0);
01390         double y = gsl_vector_get(v, 1);
01391         double a = gsl_vector_get(v, 2);
01392 
01393         double r = this_img->dot_rotate_translate(img_to, (float)x, (float)y, (float)a, mirror);
01394         int nsec = this_img->get_xsize() * this_img->get_ysize();
01395         double result = 1.0 - r / nsec;
01396 
01397 //      cout << result << " x " << x << " y " << y << " az " << a <<  endl;
01398         return result;
01399 }

static Transform refalin3d_perturbquat ( const Transform *const   t,
const float &  spincoeff,
const float &  n0,
const float &  n1,
const float &  n2,
const float &  x,
const float &  y,
const float &  z 
) [static]

Definition at line 1536 of file aligner.cpp.

References EMAN::Vec3< Type >::normalize(), q, and sqrt().

Referenced by EMAN::Refine3DAlignerQuaternion::align(), and refalifn3dquat().

01537 {
01538         Vec3f normal(n0,n1,n2);
01539         normal.normalize();
01540         
01541         float omega = spincoeff*sqrt(n0*n0 + n1*n1 + n2*n2); // Here we compute the spin by the rotation axis vector length
01542         Dict d;
01543         d["type"] = "spin";
01544         d["Omega"] = omega;
01545         d["n1"] = normal[0];
01546         d["n2"] = normal[1];
01547         d["n3"] = normal[2];
01548         //cout << omega << " " << normal[0] << " " << normal[1] << " " << normal[2] << " " << n0 << " " << n1 << " " << n2 << endl;
01549         
01550         Transform q(d);
01551         q.set_trans((float)x,(float)y,(float)z);
01552         
01553         q = q*(*t); //compose transforms        
01554         
01555         return q;
01556 }


Generated on Tue Jul 12 13:45:55 2011 for EMAN2 by  doxygen 1.4.7