#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 |
| |
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 EMAN2_ALIGNER_DEBUG 0 |
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 2270 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().
02271 { 02272 int size=rsize; 02273 float dx,dy; 02274 int bw=size/2; 02275 int MAXR=this_img->get_ysize()/2; 02276 //int MAXR=size; 02277 unsigned long tsize=2*size; 02278 unsigned long ind1=0, ind2=0, ind3=0, ind4=0, ind41=0; 02279 unsigned long index0=0; 02280 int i=0, j=0, m=0, n=0, r=0; 02281 int loop_rho=0, rho_best=0; 02282 02283 float* gnr2 = new float[size*2]; 02284 float* maxcor = new float[size+1]; // MAXR need change 02285 02286 int p_max=p_max_input; 02287 float* result = new float[5*(p_max+1)]; 02288 float* cr=new float[size*(bw+1)]; 02289 float* ci=new float[size*(bw+1)]; 02290 EMData *data_in=new EMData; 02291 data_in->set_complex(true); 02292 data_in->set_fftodd(false); 02293 data_in->set_ri(true); 02294 data_in->set_size(size+2,size,1); 02295 float *in=data_in->get_data(); 02296 02297 float *self_sampl_fft = selfpcsfft->get_data(); // ming f(r) 02298 02299 float maxcor_sofar=0.0f; 02300 int p=0; 02301 02302 for(p=0; p<=p_max; ++p){ 02303 ind1=p*size*bw; 02304 for (i=0;i<size;++i) 02305 for (j=0;j<bw+1;++j){ 02306 cr[i*(bw+1)+j]=0.0; 02307 ci[i*(bw+1)+j]=0.0; 02308 } 02309 for(n=0;n<bw;++n){ // loop for n 02310 ind2=(ind1+n); 02311 index0=n*(bw+1); 02312 for(r=0;r<=MAXR;++r) { 02313 ind3=(ind2+r*bw)*size; 02314 for(m=0;m<size;m++){ // take back hat{h(n,r,p)}(m) 02315 ind4=(ind3+m)*2; 02316 ind41=ind4+1; 02317 gnr2[2*m]=frm2dhhat[ind4]; 02318 gnr2[2*m+1]=frm2dhhat[ind41]; 02319 } 02320 for(m=0;m<bw;++m){ 02321 float tempr=self_sampl_fft[2*m+r*(size+2)]*r; 02322 float tempi=self_sampl_fft[2*m+1+r*(size+2)]*r; 02323 float gnr2_r=gnr2[2*m]; 02324 float gnr2_i=gnr2[2*m+1]; 02325 cr[n*(bw+1)+m]+=gnr2_r*tempr+gnr2_i*tempi; 02326 ci[n*(bw+1)+m]+=gnr2_i*tempr-gnr2_r*tempi; 02327 if(n!=0){ // m,-n 02328 if(m!= 0){ 02329 int ssize=tsize-2*m; // ssize = 2*size-2m 02330 int ssize1=ssize+1; 02331 float gnr2_r=gnr2[ssize]; 02332 float gnr2_i=gnr2[ssize1]; 02333 cr[(size-n)*(bw+1)+m]+=gnr2_r*tempr-gnr2_i*tempi; 02334 ci[(size-n)*(bw+1)+m]-=gnr2_i*tempr+gnr2_r*tempi; 02335 } 02336 else{ 02337 cr[(size-n)*(bw+1)+m]+=*(gnr2)*tempr-*(gnr2+1)*tempi; 02338 ci[(size-n)*(bw+1)+m]-=*(gnr2+1)*tempr+*(gnr2)*tempi; 02339 } 02340 } 02341 } 02342 } 02343 } 02344 for (int cii=0; cii<size*(bw+1);++cii){ 02345 in[2*cii]=cr[cii]; 02346 in[2*cii+1]=ci[cii]; 02347 //printf("cii=%d,in[2i+1]=%f\n",cii, cr[cii]); 02348 } 02349 02350 EMData *data_out; 02351 data_out=data_in->do_ift(); 02352 float *c=data_out->get_data(); 02353 float tempr=0.0f, corre_fcs=999.0f; 02354 02355 int n_best=0, m_best=0; 02356 float temp=-100.0f; 02357 for(n=0;n<size;++n){// move Tri_2D to Tri = c(phi,phi';rho) 02358 for(m=0;m<size;++m){ 02359 temp=c[n*size+m]; 02360 if(temp>tempr) { 02361 tempr=temp; 02362 n_best=n; 02363 m_best=m; 02364 } 02365 } 02366 } 02367 delete data_out; 02368 02369 float corre,Phi2,Phi,Tx,Ty,Vx, Vy; 02370 02371 //for (n_best=0;n_best<bw;n_best++) 02372 // for (m_best=0;m_best<2*bw;m_best++){ 02373 //n_best=0; 02374 //m_best=70; 02375 Phi2=n_best*M_PI/bw; // ming this is reference image rotation angle 02376 Phi=m_best*M_PI/bw; // ming this is particle image rotation angle 02377 Vx=p*cos(Phi);//should use the angle of the centered one 02378 Vy=-p*sin(Phi); 02379 Tx=Vx+(floor(com_this_x+0.5f)-floor(com_with_x+0.5f)); 02380 Ty=Vy+(floor(com_this_y+0.5f)-floor(com_with_y+0.5f)); 02381 02382 dx=-Tx; // the Rota & Trans value (Tx,Ty, ang_keep) are for the projection image, 02383 dy=-Ty; // need to convert to raw image 02384 02385 EMData *this_tmp=this_img->copy();//ming change to to 02386 this_tmp->rotate(-(Phi2-Phi)*180.0f,0.0f,0.0f); 02387 this_tmp->translate(dx,dy,0.0); 02388 02389 corre=this_tmp->cmp(cmp_name,to,cmp_params); 02390 //printf("corre=%f\n",corre); 02391 delete this_tmp; 02392 if(corre<=corre_fcs) { //ming, cmp use smaller value stands for more similarity 02393 corre_fcs=corre; 02394 result[0+5*p] = float(p); // rho 02395 result[1+5*p] = corre; // correlation_fcs 02396 result[2+5*p] = (Phi2-Phi)*180.0f; // rotation angle 02397 result[3+5*p] = Tx; // Translation_x 02398 result[4+5*p] = Ty; // Translation_y 02399 } 02400 maxcor[p]=corre_fcs; // maximum correlation for current rho 02401 if(corre_fcs<maxcor_sofar) { 02402 maxcor_sofar=corre_fcs; // max correlation up to current rho 02403 rho_best=p; // the rho value with maxinum correlation value 02404 } 02405 if(p>=4){ 02406 if(maxcor[p] < maxcor[p-1] && maxcor[p-1] < maxcor[p-2]&& maxcor[p-2] < maxcor[p-3] && maxcor[p-3] < maxcor[p-4]){ 02407 loop_rho=1; 02408 break; //exit p loop 02409 } 02410 } 02411 } // end for p 02412 //}//test my method 02413 if(loop_rho == 1) 02414 p=p+1; 02415 int rb5=5*rho_best; 02416 float fsc = result[1+rb5]; 02417 float ang_keep = result[2+rb5]; 02418 float Tx = result[3+rb5]; 02419 float Ty = result[4+rb5]; 02420 delete[] gnr2; 02421 delete[] maxcor; 02422 delete[] result; 02423 delete cr; 02424 cr=0; 02425 delete ci; 02426 ci=0; 02427 delete data_in; // ming add 02428 dx = -Tx; // the Rota & Trans value (Tx,Ty, ang_keep) are for the projection image, 02429 dy = -Ty; // need to convert to raw image 02430 this_img->rotate(-ang_keep,0,0); // ming change this to this_img?? 02431 this_img->translate(dx,dy,0.0); // ming change this to this_img 02432 02433 02434 Transform tsoln(Dict("type","2d","alpha",ang_keep)); 02435 tsoln.set_trans(dx,dy); 02436 this_img->set_attr("xform.align2d",&tsoln); 02437 #ifdef DEBUG 02438 float fsc_best=this_img->cmp(cmp_name,to,cmp_params); 02439 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); 02440 #endif 02441 return fsc; // return the fsc coefficients 02442 } // FRM2D aligner sub_class
static double refalifn | ( | const gsl_vector * | v, | |
void * | params | |||
) | [static] |
Definition at line 1311 of file aligner.cpp.
References EMAN::Cmp::cmp(), EMAN::EMData::process(), t, x, and y.
Referenced by EMAN::RefineAligner::align().
01312 { 01313 Dict *dict = (Dict *) params; 01314 01315 double x = gsl_vector_get(v, 0); 01316 double y = gsl_vector_get(v, 1); 01317 double a = gsl_vector_get(v, 2); 01318 01319 EMData *this_img = (*dict)["this"]; 01320 EMData *with = (*dict)["with"]; 01321 bool mirror = (*dict)["mirror"]; 01322 01323 // float mean = (float)this_img->get_attr("mean"); 01324 // if ( Util::goodf(&mean) ) { 01325 // //cout << "tmps mean is nan even before rotation" << endl; 01326 // } 01327 01328 Transform t(Dict("type","2d","alpha",static_cast<float>(a))); 01329 // Transform3D t3d(Transform3D::EMAN, (float)a, 0.0f, 0.0f); 01330 // t3d.set_posttrans( (float) x, (float) y); 01331 // tmp->rotate_translate(t3d); 01332 t.set_trans((float)x,(float)y); 01333 t.set_mirror(mirror); 01334 if (v->size>3) { 01335 float sca=(float)gsl_vector_get(v, 3); 01336 if (sca<.7 || sca>1.3) return 1.0e20; 01337 t.set_scale((float)gsl_vector_get(v, 3)); 01338 } 01339 EMData *tmp = this_img->process("xform",Dict("transform",&t)); 01340 01341 // printf("GSL %f %f %f %d %f\n",x,y,a,mirror,(float)gsl_vector_get(v, 3)); 01342 Cmp* c = (Cmp*) ((void*)(*dict)["cmp"]); 01343 double result = c->cmp(tmp,with); 01344 01345 // DELETE AT SOME STAGE, USEFUL FOR PRERELEASE STUFF 01346 // float test_result = (float)result; 01347 // if ( Util::goodf(&test_result) ) { 01348 // cout << "result " << result << " " << x << " " << y << " " << a << endl; 01349 // cout << (float)this_img->get_attr("mean") << " " << (float)tmp->get_attr("mean") << " " << (float)with->get_attr("mean") << endl; 01350 // tmp->write_image("tmp.hdf"); 01351 // with->write_image("with.hdf"); 01352 // this_img->write_image("this_img.hdf"); 01353 // EMData* t = this_img->copy(); 01354 // cout << (float)t->get_attr("mean") << endl; 01355 // t->rotate_translate( t3d ); 01356 // cout << (float)t->get_attr("mean") << endl; 01357 // cout << "exit" << endl; 01359 // cout << (float)t->get_attr("mean") << endl; 01360 // cout << "now exit" << endl; 01361 // delete t; 01362 // } 01363 01364 01365 if ( tmp != 0 ) delete tmp; 01366 01367 return result; 01368 }
static double refalifn3dquat | ( | const gsl_vector * | v, | |
void * | params | |||
) | [static] |
Definition at line 1546 of file aligner.cpp.
References EMAN::Cmp::cmp(), EMAN::EMData::process(), refalin3d_perturbquat(), t, x, and y.
Referenced by EMAN::Refine3DAlignerQuaternion::align().
01547 { 01548 Dict *dict = (Dict *) params; 01549 01550 double n0 = gsl_vector_get(v, 0); 01551 double n1 = gsl_vector_get(v, 1); 01552 double n2 = gsl_vector_get(v, 2); 01553 double x = gsl_vector_get(v, 3); 01554 double y = gsl_vector_get(v, 4); 01555 double z = gsl_vector_get(v, 5); 01556 01557 EMData *this_img = (*dict)["this"]; 01558 EMData *with = (*dict)["with"]; 01559 // bool mirror = (*dict)["mirror"]; 01560 01561 Transform* t = (*dict)["transform"]; 01562 float spincoeff = (*dict)["spincoeff"]; 01563 01564 Transform soln = refalin3d_perturbquat(t,spincoeff,(float)n0,(float)n1,(float)n2,(float)x,(float)y,(float)z); 01565 01566 EMData *tmp = this_img->process("xform",Dict("transform",&soln)); 01567 Cmp* c = (Cmp*) ((void*)(*dict)["cmp"]); 01568 double result = c->cmp(tmp,with); 01569 if ( tmp != 0 ) delete tmp; 01570 delete t; t = 0; 01571 //cout << result << endl; 01572 return result; 01573 }
static double refalifnfast | ( | const gsl_vector * | v, | |
void * | params | |||
) | [static] |
Definition at line 1370 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().
01371 { 01372 Dict *dict = (Dict *) params; 01373 EMData *this_img = (*dict)["this"]; 01374 EMData *img_to = (*dict)["with"]; 01375 bool mirror = (*dict)["mirror"]; 01376 01377 double x = gsl_vector_get(v, 0); 01378 double y = gsl_vector_get(v, 1); 01379 double a = gsl_vector_get(v, 2); 01380 01381 double r = this_img->dot_rotate_translate(img_to, (float)x, (float)y, (float)a, mirror); 01382 int nsec = this_img->get_xsize() * this_img->get_ysize(); 01383 double result = 1.0 - r / nsec; 01384 01385 // cout << result << " x " << x << " y " << y << " az " << a << endl; 01386 return result; 01387 }
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 1524 of file aligner.cpp.
References EMAN::Vec3< Type >::normalize(), q, and sqrt().
Referenced by EMAN::Refine3DAlignerQuaternion::align(), and refalifn3dquat().
01525 { 01526 Vec3f normal(n0,n1,n2); 01527 normal.normalize(); 01528 01529 float omega = spincoeff*sqrt(n0*n0 + n1*n1 + n2*n2); // Here we compute the spin by the rotation axis vector length 01530 Dict d; 01531 d["type"] = "spin"; 01532 d["Omega"] = omega; 01533 d["n1"] = normal[0]; 01534 d["n2"] = normal[1]; 01535 d["n3"] = normal[2]; 01536 //cout << omega << " " << normal[0] << " " << normal[1] << " " << normal[2] << " " << n0 << " " << n1 << " " << n2 << endl; 01537 01538 Transform q(d); 01539 q.set_trans((float)x,(float)y,(float)z); 01540 01541 q = q*(*t); //compose transforms 01542 01543 return q; 01544 }