#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 2178 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().
02179 { 02180 int size=rsize; 02181 float dx,dy; 02182 int bw=size/2; 02183 int MAXR=this_img->get_ysize()/2; 02184 //int MAXR=size; 02185 unsigned long tsize=2*size; 02186 unsigned long ind1=0, ind2=0, ind3=0, ind4=0, ind41=0; 02187 unsigned long index0=0; 02188 int i=0, j=0, m=0, n=0, r=0; 02189 int loop_rho=0, rho_best=0; 02190 02191 float* gnr2 = new float[size*2]; 02192 float* maxcor = new float[size+1]; // MAXR need change 02193 02194 int p_max=p_max_input; 02195 float* result = new float[5*(p_max+1)]; 02196 float* cr=new float[size*(bw+1)]; 02197 float* ci=new float[size*(bw+1)]; 02198 EMData *data_in=new EMData; 02199 data_in->set_complex(true); 02200 data_in->set_fftodd(false); 02201 data_in->set_ri(true); 02202 data_in->set_size(size+2,size,1); 02203 float *in=data_in->get_data(); 02204 02205 float *self_sampl_fft = selfpcsfft->get_data(); // ming f(r) 02206 02207 float maxcor_sofar=0.0f; 02208 int p=0; 02209 02210 for(p=0; p<=p_max; ++p){ 02211 ind1=p*size*bw; 02212 for (i=0;i<size;++i) 02213 for (j=0;j<bw+1;++j){ 02214 cr[i*(bw+1)+j]=0.0; 02215 ci[i*(bw+1)+j]=0.0; 02216 } 02217 for(n=0;n<bw;++n){ // loop for n 02218 ind2=(ind1+n); 02219 index0=n*(bw+1); 02220 for(r=0;r<=MAXR;++r) { 02221 ind3=(ind2+r*bw)*size; 02222 for(m=0;m<size;m++){ // take back hat{h(n,r,p)}(m) 02223 ind4=(ind3+m)*2; 02224 ind41=ind4+1; 02225 gnr2[2*m]=frm2dhhat[ind4]; 02226 gnr2[2*m+1]=frm2dhhat[ind41]; 02227 } 02228 for(m=0;m<bw;++m){ 02229 float tempr=self_sampl_fft[2*m+r*(size+2)]*r; 02230 float tempi=self_sampl_fft[2*m+1+r*(size+2)]*r; 02231 float gnr2_r=gnr2[2*m]; 02232 float gnr2_i=gnr2[2*m+1]; 02233 cr[n*(bw+1)+m]+=gnr2_r*tempr+gnr2_i*tempi; 02234 ci[n*(bw+1)+m]+=gnr2_i*tempr-gnr2_r*tempi; 02235 if(n!=0){ // m,-n 02236 if(m!= 0){ 02237 int ssize=tsize-2*m; // ssize = 2*size-2m 02238 int ssize1=ssize+1; 02239 float gnr2_r=gnr2[ssize]; 02240 float gnr2_i=gnr2[ssize1]; 02241 cr[(size-n)*(bw+1)+m]+=gnr2_r*tempr-gnr2_i*tempi; 02242 ci[(size-n)*(bw+1)+m]-=gnr2_i*tempr+gnr2_r*tempi; 02243 } 02244 else{ 02245 cr[(size-n)*(bw+1)+m]+=*(gnr2)*tempr-*(gnr2+1)*tempi; 02246 ci[(size-n)*(bw+1)+m]-=*(gnr2+1)*tempr+*(gnr2)*tempi; 02247 } 02248 } 02249 } 02250 } 02251 } 02252 for (int cii=0; cii<size*(bw+1);++cii){ 02253 in[2*cii]=cr[cii]; 02254 in[2*cii+1]=ci[cii]; 02255 //printf("cii=%d,in[2i+1]=%f\n",cii, cr[cii]); 02256 } 02257 02258 EMData *data_out; 02259 data_out=data_in->do_ift(); 02260 float *c=data_out->get_data(); 02261 float tempr=0.0f, corre_fcs=999.0f; 02262 02263 int n_best=0, m_best=0; 02264 float temp=-100.0f; 02265 for(n=0;n<size;++n){// move Tri_2D to Tri = c(phi,phi';rho) 02266 for(m=0;m<size;++m){ 02267 temp=c[n*size+m]; 02268 if(temp>tempr) { 02269 tempr=temp; 02270 n_best=n; 02271 m_best=m; 02272 } 02273 } 02274 } 02275 delete data_out; 02276 02277 float corre,Phi2,Phi,Tx,Ty,Vx, Vy; 02278 02279 //for (n_best=0;n_best<bw;n_best++) 02280 // for (m_best=0;m_best<2*bw;m_best++){ 02281 //n_best=0; 02282 //m_best=70; 02283 Phi2=n_best*M_PI/bw; // ming this is reference image rotation angle 02284 Phi=m_best*M_PI/bw; // ming this is particle image rotation angle 02285 Vx=p*cos(Phi);//should use the angle of the centered one 02286 Vy=-p*sin(Phi); 02287 Tx=Vx+(floor(com_this_x+0.5f)-floor(com_with_x+0.5f)); 02288 Ty=Vy+(floor(com_this_y+0.5f)-floor(com_with_y+0.5f)); 02289 02290 dx=-Tx; // the Rota & Trans value (Tx,Ty, ang_keep) are for the projection image, 02291 dy=-Ty; // need to convert to raw image 02292 02293 EMData *this_tmp=this_img->copy();//ming change to to 02294 this_tmp->rotate(-(Phi2-Phi)*180.0f,0.0f,0.0f); 02295 this_tmp->translate(dx,dy,0.0); 02296 02297 corre=this_tmp->cmp(cmp_name,to,cmp_params); 02298 //printf("corre=%f\n",corre); 02299 delete this_tmp; 02300 if(corre<=corre_fcs) { //ming, cmp use smaller value stands for more similarity 02301 corre_fcs=corre; 02302 result[0+5*p] = float(p); // rho 02303 result[1+5*p] = corre; // correlation_fcs 02304 result[2+5*p] = (Phi2-Phi)*180.0f; // rotation angle 02305 result[3+5*p] = Tx; // Translation_x 02306 result[4+5*p] = Ty; // Translation_y 02307 } 02308 maxcor[p]=corre_fcs; // maximum correlation for current rho 02309 if(corre_fcs<maxcor_sofar) { 02310 maxcor_sofar=corre_fcs; // max correlation up to current rho 02311 rho_best=p; // the rho value with maxinum correlation value 02312 } 02313 if(p>=4){ 02314 if(maxcor[p] < maxcor[p-1] && maxcor[p-1] < maxcor[p-2]&& maxcor[p-2] < maxcor[p-3] && maxcor[p-3] < maxcor[p-4]){ 02315 loop_rho=1; 02316 break; //exit p loop 02317 } 02318 } 02319 } // end for p 02320 //}//test my method 02321 if(loop_rho == 1) 02322 p=p+1; 02323 int rb5=5*rho_best; 02324 float fsc = result[1+rb5]; 02325 float ang_keep = result[2+rb5]; 02326 float Tx = result[3+rb5]; 02327 float Ty = result[4+rb5]; 02328 delete[] gnr2; 02329 delete[] maxcor; 02330 delete[] result; 02331 delete cr; 02332 cr=0; 02333 delete ci; 02334 ci=0; 02335 delete data_in; // ming add 02336 dx = -Tx; // the Rota & Trans value (Tx,Ty, ang_keep) are for the projection image, 02337 dy = -Ty; // need to convert to raw image 02338 this_img->rotate(-ang_keep,0,0); // ming change this to this_img?? 02339 this_img->translate(dx,dy,0.0); // ming change this to this_img 02340 02341 02342 Transform tsoln(Dict("type","2d","alpha",ang_keep)); 02343 tsoln.set_trans(dx,dy); 02344 this_img->set_attr("xform.align2d",&tsoln); 02345 #ifdef DEBUG 02346 float fsc_best=this_img->cmp(cmp_name,to,cmp_params); 02347 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); 02348 #endif 02349 return fsc; // return the fsc coefficients 02350 } // FRM2D aligner sub_class
static double refalifn | ( | const gsl_vector * | v, | |
void * | params | |||
) | [static] |
Definition at line 1310 of file aligner.cpp.
References EMAN::Cmp::cmp(), EMAN::EMData::process(), t, x, and y.
Referenced by EMAN::RefineAligner::align().
01311 { 01312 Dict *dict = (Dict *) params; 01313 01314 double x = gsl_vector_get(v, 0); 01315 double y = gsl_vector_get(v, 1); 01316 double a = gsl_vector_get(v, 2); 01317 01318 EMData *this_img = (*dict)["this"]; 01319 EMData *with = (*dict)["with"]; 01320 bool mirror = (*dict)["mirror"]; 01321 01322 // float mean = (float)this_img->get_attr("mean"); 01323 // if ( Util::goodf(&mean) ) { 01324 // //cout << "tmps mean is nan even before rotation" << endl; 01325 // } 01326 01327 Transform t(Dict("type","2d","alpha",static_cast<float>(a))); 01328 // Transform3D t3d(Transform3D::EMAN, (float)a, 0.0f, 0.0f); 01329 // t3d.set_posttrans( (float) x, (float) y); 01330 // tmp->rotate_translate(t3d); 01331 t.set_trans((float)x,(float)y); 01332 t.set_mirror(mirror); 01333 if (v->size>3) { 01334 float sca=(float)gsl_vector_get(v, 3); 01335 if (sca<.7 || sca>1.3) return 1.0e20; 01336 t.set_scale((float)gsl_vector_get(v, 3)); 01337 } 01338 EMData *tmp = this_img->process("xform",Dict("transform",&t)); 01339 01340 // printf("GSL %f %f %f %d %f\n",x,y,a,mirror,(float)gsl_vector_get(v, 3)); 01341 Cmp* c = (Cmp*) ((void*)(*dict)["cmp"]); 01342 double result = c->cmp(tmp,with); 01343 01344 // DELETE AT SOME STAGE, USEFUL FOR PRERELEASE STUFF 01345 // float test_result = (float)result; 01346 // if ( Util::goodf(&test_result) ) { 01347 // cout << "result " << result << " " << x << " " << y << " " << a << endl; 01348 // cout << (float)this_img->get_attr("mean") << " " << (float)tmp->get_attr("mean") << " " << (float)with->get_attr("mean") << endl; 01349 // tmp->write_image("tmp.hdf"); 01350 // with->write_image("with.hdf"); 01351 // this_img->write_image("this_img.hdf"); 01352 // EMData* t = this_img->copy(); 01353 // cout << (float)t->get_attr("mean") << endl; 01354 // t->rotate_translate( t3d ); 01355 // cout << (float)t->get_attr("mean") << endl; 01356 // cout << "exit" << endl; 01358 // cout << (float)t->get_attr("mean") << endl; 01359 // cout << "now exit" << endl; 01360 // delete t; 01361 // } 01362 01363 01364 if ( tmp != 0 ) delete tmp; 01365 01366 return result; 01367 }
static double refalifn3dquat | ( | const gsl_vector * | v, | |
void * | params | |||
) | [static] |
Definition at line 1545 of file aligner.cpp.
References EMAN::Cmp::cmp(), EMAN::EMData::process(), refalin3d_perturbquat(), t, x, and y.
Referenced by EMAN::Refine3DAlignerQuaternion::align().
01546 { 01547 Dict *dict = (Dict *) params; 01548 01549 double n0 = gsl_vector_get(v, 0); 01550 double n1 = gsl_vector_get(v, 1); 01551 double n2 = gsl_vector_get(v, 2); 01552 double x = gsl_vector_get(v, 3); 01553 double y = gsl_vector_get(v, 4); 01554 double z = gsl_vector_get(v, 5); 01555 01556 EMData *this_img = (*dict)["this"]; 01557 EMData *with = (*dict)["with"]; 01558 // bool mirror = (*dict)["mirror"]; 01559 01560 Transform* t = (*dict)["transform"]; 01561 float spincoeff = (*dict)["spincoeff"]; 01562 01563 Transform soln = refalin3d_perturbquat(t,spincoeff,(float)n0,(float)n1,(float)n2,(float)x,(float)y,(float)z); 01564 01565 EMData *tmp = this_img->process("xform",Dict("transform",&soln)); 01566 Cmp* c = (Cmp*) ((void*)(*dict)["cmp"]); 01567 double result = c->cmp(tmp,with); 01568 if ( tmp != 0 ) delete tmp; 01569 delete t; t = 0; 01570 //cout << result << endl; 01571 return result; 01572 }
static double refalifnfast | ( | const gsl_vector * | v, | |
void * | params | |||
) | [static] |
Definition at line 1369 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().
01370 { 01371 Dict *dict = (Dict *) params; 01372 EMData *this_img = (*dict)["this"]; 01373 EMData *img_to = (*dict)["with"]; 01374 bool mirror = (*dict)["mirror"]; 01375 01376 double x = gsl_vector_get(v, 0); 01377 double y = gsl_vector_get(v, 1); 01378 double a = gsl_vector_get(v, 2); 01379 01380 double r = this_img->dot_rotate_translate(img_to, (float)x, (float)y, (float)a, mirror); 01381 int nsec = this_img->get_xsize() * this_img->get_ysize(); 01382 double result = 1.0 - r / nsec; 01383 01384 // cout << result << " x " << x << " y " << y << " az " << a << endl; 01385 return result; 01386 }
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 1523 of file aligner.cpp.
References EMAN::Vec3< Type >::normalize(), q, and sqrt().
Referenced by EMAN::Refine3DAlignerQuaternion::align(), and refalifn3dquat().
01524 { 01525 Vec3f normal(n0,n1,n2); 01526 normal.normalize(); 01527 01528 float omega = spincoeff*sqrt(n0*n0 + n1*n1 + n2*n2); // Here we compute the spin by the rotation axis vector length 01529 Dict d; 01530 d["type"] = "spin"; 01531 d["Omega"] = omega; 01532 d["n1"] = normal[0]; 01533 d["n2"] = normal[1]; 01534 d["n3"] = normal[2]; 01535 //cout << omega << " " << normal[0] << " " << normal[1] << " " << normal[2] << " " << n0 << " " << n1 << " " << n2 << endl; 01536 01537 Transform q(d); 01538 q.set_trans((float)x,(float)y,(float)z); 01539 01540 q = q*(*t); //compose transforms 01541 01542 return q; 01543 }