00001
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
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 "cuda/cuda_processor.h"
00047 #include "cuda/cuda_cmp.h"
00048 #endif
00049
00050 #ifdef SPARX_USING_CUDA
00051 #include <sparx/cuda/cuda_ccf.h>
00052 #endif
00053
00054 #define EMAN2_ALIGNER_DEBUG 0
00055
00056 using namespace EMAN;
00057
00058 const string TranslationalAligner::NAME = "translational";
00059 const string RotationalAligner::NAME = "rotational";
00060 const string RotationalAlignerIterative::NAME = "rotational_iterative";
00061 const string RotatePrecenterAligner::NAME = "rotate_precenter";
00062 const string RotateTranslateAligner::NAME = "rotate_translate";
00063 const string RotateTranslateAlignerIterative::NAME = "rotate_translate_iterative";
00064 const string RotateTranslateAlignerPawel::NAME = "rotate_translate_resample";
00065 const string RotateTranslateBestAligner::NAME = "rotate_translate_best";
00066 const string RotateFlipAligner::NAME = "rotate_flip";
00067 const string RotateFlipAlignerIterative::NAME = "rotate_flip.iterative";
00068 const string RotateTranslateFlipAligner::NAME = "rotate_translate_flip";
00069 const string RotateTranslateFlipAlignerIterative::NAME = "rotate_translate_flip_iterative";
00070 const string RotateTranslateFlipAlignerPawel::NAME = "rotate_translate_flip_resample";
00071 const string RTFExhaustiveAligner::NAME = "rtf_exhaustive";
00072 const string RTFSlowExhaustiveAligner::NAME = "rtf_slow_exhaustive";
00073 const string RefineAligner::NAME = "refine";
00074 const string Refine3DAlignerGrid::NAME = "refine_3d_grid";
00075 const string Refine3DAlignerQuaternion::NAME = "refine_3d";
00076 const string RT3DGridAligner::NAME = "rotate_translate_3d_grid";
00077 const string RT3DSphereAligner::NAME = "rotate_translate_3d";
00078 const string RT3DSymmetryAligner::NAME = "rotate_symmetry_3d";
00079 const string FRM2DAligner::NAME = "frm2d";
00080
00081
00082 template <> Factory < Aligner >::Factory()
00083 {
00084 force_add<TranslationalAligner>();
00085 force_add<RotationalAligner>();
00086 force_add<RotationalAlignerIterative>();
00087 force_add<RotatePrecenterAligner>();
00088 force_add<RotateTranslateAligner>();
00089 force_add<RotateTranslateAlignerIterative>();
00090 force_add<RotateTranslateAlignerPawel>();
00091 force_add<RotateFlipAligner>();
00092 force_add<RotateFlipAlignerIterative>();
00093 force_add<RotateTranslateFlipAligner>();
00094 force_add<RotateTranslateFlipAlignerIterative>();
00095 force_add<RotateTranslateFlipAlignerPawel>();
00096 force_add<RTFExhaustiveAligner>();
00097 force_add<RTFSlowExhaustiveAligner>();
00098 force_add<RefineAligner>();
00099 force_add<Refine3DAlignerGrid>();
00100 force_add<Refine3DAlignerQuaternion>();
00101 force_add<RT3DGridAligner>();
00102 force_add<RT3DSphereAligner>();
00103 force_add<RT3DSymmetryAligner>();
00104 force_add<FRM2DAligner>();
00105
00106 }
00107
00108 vector<Dict> Aligner::xform_align_nbest(EMData *, EMData *, const unsigned int, const string &, const Dict&) const
00109 {
00110 vector<Dict> solns;
00111 return solns;
00112 }
00113
00114
00115
00116
00117
00118 EMData *TranslationalAligner::align(EMData * this_img, EMData *to,
00119 const string&, const Dict&) const
00120 {
00121 if (!this_img) {
00122 return 0;
00123 }
00124
00125 if (to && !EMUtil::is_same_size(this_img, to))
00126 throw ImageDimensionException("Images must be the same size to perform translational alignment");
00127
00128 EMData *cf = 0;
00129 int nx = this_img->get_xsize();
00130 int ny = this_img->get_ysize();
00131 int nz = this_img->get_zsize();
00132
00133 int masked = params.set_default("masked",0);
00134 int useflcf = params.set_default("useflcf",0);
00135 bool use_cpu = true;
00136
00137 #ifdef EMAN2_USING_CUDA
00138 if(EMData::usecuda == 1) {
00139 if(!this_img->cudarwdata) this_img->copy_to_cuda();
00140 if(!to->cudarwdata) to->copy_to_cuda();
00141 if (masked) throw UnexpectedBehaviorException("Masked is not yet supported in CUDA");
00142 if (useflcf) throw UnexpectedBehaviorException("Useflcf is not yet supported in CUDA");
00143
00144 use_cpu = false;
00145 cf = this_img->calc_ccf(to);
00146 }
00147 #endif // EMAN2_USING_CUDA
00148
00149 if (use_cpu) {
00150 if (useflcf) cf = this_img->calc_flcf(to);
00151 else cf = this_img->calc_ccf(to);
00152 }
00153
00154
00155 if (masked) {
00156 EMData *msk=this_img->process("threshold.notzero");
00157 EMData *sqr=to->process("math.squared");
00158 EMData *cfn=msk->calc_ccf(sqr);
00159 cfn->process_inplace("math.sqrt");
00160 float *d1=cf->get_data();
00161 float *d2=cfn->get_data();
00162 for (size_t i=0; i<(size_t)nx*ny*nz; ++i) {
00163 if (d2[i]!=0) d1[i]/=d2[i];
00164 }
00165 cf->update();
00166 delete msk;
00167 delete sqr;
00168 delete cfn;
00169 }
00170
00171 int maxshiftx = params.set_default("maxshift",-1);
00172 int maxshifty = params["maxshift"];
00173 int maxshiftz = params["maxshift"];
00174 int nozero = params["nozero"];
00175
00176 if (maxshiftx <= 0) {
00177 maxshiftx = nx / 4;
00178 maxshifty = ny / 4;
00179 maxshiftz = nz / 4;
00180 }
00181
00182 if (maxshiftx > nx / 2 - 1) maxshiftx = nx / 2 - 1;
00183 if (maxshifty > ny / 2 - 1) maxshifty = ny / 2 - 1;
00184 if (maxshiftz > nz / 2 - 1) maxshiftz = nz / 2 - 1;
00185
00186 if (nx == 1) maxshiftx = 0;
00187 if (ny == 1) maxshifty = 0;
00188 if (nz == 1) maxshiftz = 0;
00189
00190
00191 if (nozero) {
00192 cf->zero_corner_circulant(1);
00193 }
00194
00195 IntPoint peak;
00196 #ifdef EMAN2_USING_CUDA
00197 if (!use_cpu) {
00198 if (nozero) throw UnexpectedBehaviorException("Nozero is not yet supported in CUDA");
00199 CudaPeakInfo* data = calc_max_location_wrap_cuda(cf->cudarwdata, cf->get_xsize(), cf->get_ysize(), cf->get_zsize(), maxshiftx, maxshifty, maxshiftz);
00200 peak = IntPoint(data->px,data->py,data->pz);
00201 free(data);
00202 }
00203 #endif // EMAN2_USING_CUDA
00204
00205 if (use_cpu) {
00206 peak = cf->calc_max_location_wrap(maxshiftx, maxshifty, maxshiftz);
00207 }
00208
00209 Vec3f cur_trans = Vec3f ( (float)-peak[0], (float)-peak[1], (float)-peak[2]);
00210
00211
00212 if (!to) {
00213 cur_trans /= 2.0f;
00214 int intonly = params.set_default("intonly",false);
00215 if (intonly) {
00216 cur_trans[0] = floor(cur_trans[0] + 0.5f);
00217 cur_trans[1] = floor(cur_trans[1] + 0.5f);
00218 cur_trans[2] = floor(cur_trans[2] + 0.5f);
00219 }
00220 }
00221
00222 if( cf ){
00223 delete cf;
00224 cf = 0;
00225 }
00226
00227 Dict params("trans",static_cast< vector<int> >(cur_trans));
00228 if (use_cpu){
00229 cf=this_img->process("math.translate.int",params);
00230 }
00231 Transform t;
00232 t.set_trans(cur_trans);
00233
00234 #ifdef EMAN2_USING_CUDA
00235 if (!use_cpu) {
00236
00237 cf = this_img->process("xform",Dict("transform",&t));
00238 }
00239 #endif // EMAN2_USING_CUDA
00240
00241 if ( nz != 1 ) {
00242
00243
00244 cf->set_attr("xform.align3d",&t);
00245 } else if ( ny != 1 ) {
00246
00247 cur_trans[2] = 0;
00248 t.set_trans(cur_trans);
00249 cf->set_attr("xform.align2d",&t);
00250 }
00251
00252 return cf;
00253 }
00254
00255 EMData * RotationalAligner::align_180_ambiguous(EMData * this_img, EMData * to, int rfp_mode) {
00256
00257
00258 EMData* this_img_rfp, * to_rfp;
00259 if (rfp_mode == 0) {
00260 this_img_rfp = this_img->make_rotational_footprint_e1();
00261 to_rfp = to->make_rotational_footprint_e1();
00262 } else if (rfp_mode == 1) {
00263 this_img_rfp = this_img->make_rotational_footprint();
00264 to_rfp = to->make_rotational_footprint();
00265 } else if (rfp_mode == 2) {
00266 this_img_rfp = this_img->make_rotational_footprint_cmc();
00267 to_rfp = to->make_rotational_footprint_cmc();
00268 } else {
00269 throw InvalidParameterException("rfp_mode must be 0,1 or 2");
00270 }
00271 int this_img_rfp_nx = this_img_rfp->get_xsize();
00272
00273
00274 EMData *cf = this_img_rfp->calc_ccfx(to_rfp, 0, this_img->get_ysize());
00275
00276
00277 delete this_img_rfp; this_img_rfp = 0;
00278 delete to_rfp; to_rfp = 0;
00279
00280
00281 float *data = cf->get_data();
00282
00283 float peak = 0;
00284 int peak_index = 0;
00285 Util::find_max(data, this_img_rfp_nx, &peak, &peak_index);
00286
00287 if( cf ) {
00288 delete cf;
00289 cf = 0;
00290 }
00291 float rot_angle = (float) (peak_index * 180.0f / this_img_rfp_nx);
00292
00293
00294 Transform tmp(Dict("type","2d","alpha",rot_angle));
00295 cf=this_img->process("xform",Dict("transform",(Transform*)&tmp));
00296
00297
00298
00299 cf->set_attr("xform.align2d",&tmp);
00300 return cf;
00301 }
00302
00303 EMData *RotationalAligner::align(EMData * this_img, EMData *to,
00304 const string& cmp_name, const Dict& cmp_params) const
00305 {
00306 if (!to) throw InvalidParameterException("Can not rotational align - the image to align to is NULL");
00307
00308 #ifdef EMAN2_USING_CUDA
00309 if(EMData::usecuda == 1) {
00310 if(!this_img->cudarwdata) this_img->copy_to_cuda();
00311 if(!to->cudarwdata) to->copy_to_cuda();
00312 }
00313 #endif
00314
00315
00316 int rfp_mode = params.set_default("rfp_mode",0);
00317 EMData* rot_aligned = RotationalAligner::align_180_ambiguous(this_img,to,rfp_mode);
00318 Transform * tmp = rot_aligned->get_attr("xform.align2d");
00319 Dict rot = tmp->get_rotation("2d");
00320 float rotate_angle_solution = rot["alpha"];
00321 delete tmp;
00322
00323 EMData *rot_align_180 = rot_aligned->process("math.rotate.180");
00324
00325
00326 float rot_cmp = rot_aligned->cmp(cmp_name, to, cmp_params);
00327 float rot_180_cmp = rot_align_180->cmp(cmp_name, to, cmp_params);
00328
00329
00330 float score = 0.0;
00331 EMData* result = NULL;
00332 if (rot_cmp < rot_180_cmp){
00333 result = rot_aligned;
00334 score = rot_cmp;
00335 delete rot_align_180; rot_align_180 = 0;
00336 } else {
00337 result = rot_align_180;
00338 score = rot_180_cmp;
00339 delete rot_aligned; rot_aligned = 0;
00340 rotate_angle_solution = rotate_angle_solution-180.0f;
00341 }
00342
00343
00344
00345 Transform tmp2(Dict("type","2d","alpha",rotate_angle_solution));
00346 result->set_attr("xform.align2d",&tmp2);
00347 return result;
00348 }
00349
00350
00351 EMData *RotatePrecenterAligner::align(EMData * this_img, EMData *to,
00352 const string&, const Dict&) const
00353 {
00354 if (!to) {
00355 return 0;
00356 }
00357
00358 int ny = this_img->get_ysize();
00359 int size = Util::calc_best_fft_size((int) (M_PI * ny * 1.5));
00360 EMData *e1 = this_img->unwrap(4, ny * 7 / 16, size, 0, 0, 1);
00361 EMData *e2 = to->unwrap(4, ny * 7 / 16, size, 0, 0, 1);
00362 EMData *cf = e1->calc_ccfx(e2, 0, ny);
00363
00364 float *data = cf->get_data();
00365
00366 float peak = 0;
00367 int peak_index = 0;
00368 Util::find_max(data, size, &peak, &peak_index);
00369 float a = (float) ((1.0f - 1.0f * peak_index / size) * 180. * 2);
00370
00371 Transform rot;
00372 rot.set_rotation(Dict("type","2d","alpha",(float)(a*180./M_PI)));
00373 EMData* rslt = this_img->process("xform",Dict("transform",&rot));
00374 rslt->set_attr("xform.align2d",&rot);
00375
00376
00377
00378
00379
00380
00381
00382
00383
00384
00385 if( e1 )
00386 {
00387 delete e1;
00388 e1 = 0;
00389 }
00390
00391 if( e2 )
00392 {
00393 delete e2;
00394 e2 = 0;
00395 }
00396
00397 if (cf) {
00398 delete cf;
00399 cf = 0;
00400 }
00401 return rslt;
00402 }
00403
00404 EMData *RotationalAlignerIterative::align(EMData * this_img, EMData *to,
00405 const string &, const Dict&) const
00406 {
00407 int r1 = params.set_default("r1",-1);
00408 int r2 = params.set_default("r2",-1);
00409
00410 EMData * to_polar = to->unwrap(r1,r2,-1,0,0,true);
00411 EMData * this_img_polar = this_img->unwrap(r1,r2,-1,0,0,true);
00412 int this_img_polar_nx = this_img_polar->get_xsize();
00413
00414 EMData *cf = this_img_polar->calc_ccfx(to_polar, 0, this_img->get_ysize());
00415
00416
00417 delete to_polar; to_polar = 0;
00418 delete this_img_polar; this_img_polar = 0;
00419
00420 float *data = cf->get_data();
00421 float peak = 0;
00422 int peak_index = 0;
00423 Util::find_max(data, this_img_polar_nx, &peak, &peak_index);
00424
00425 delete cf; cf = 0;
00426 float rot_angle = (float) (peak_index * 360.0f / this_img_polar_nx);
00427
00428
00429
00430 Transform tmp(Dict("type","2d","alpha",rot_angle));
00431 EMData * rotimg=this_img->process("xform",Dict("transform",(Transform*)&tmp));
00432 rotimg->set_attr("xform.align2d",&tmp);
00433
00434 return rotimg;
00435
00436 }
00437
00438 EMData *RotateTranslateAlignerIterative::align(EMData * this_img, EMData *to,
00439 const string & cmp_name, const Dict& cmp_params) const
00440 {
00441 int maxiter = params.set_default("maxiter", 3);
00442
00443 Dict trans_params;
00444 trans_params["intonly"] = 0;
00445 trans_params["maxshift"] = params.set_default("maxshift", -1);
00446 trans_params["useflcf"] = params.set_default("useflcf",0);
00447 trans_params["nozero"] = params.set_default("nozero",false);
00448
00449 Dict rot_params;
00450 rot_params["r1"] = params.set_default("r1", -1);
00451 rot_params["r2"] = params.set_default("r2", -1);
00452
00453 Transform t;
00454 EMData * moving_img = this_img;
00455 for(int it = 0; it < maxiter; it++){
00456
00457
00458 EMData * trans_align = moving_img->align("translational", to, trans_params, cmp_name, cmp_params);
00459 Transform * tt = trans_align->get_attr("xform.align2d");
00460 t = *tt*t;
00461 delete tt;
00462
00463
00464 EMData * rottrans_align = trans_align->align("rotational.iterative", to, rot_params, cmp_name, cmp_params);
00465 Transform * rt = rottrans_align->get_attr("xform.align2d");
00466 t = *rt*t;
00467 delete trans_align; trans_align = 0;
00468 delete rottrans_align; rottrans_align = 0;
00469 delete rt;
00470
00471
00472 if(it > 0){delete moving_img;}
00473
00474 moving_img = this_img->process("xform",Dict("transform",&t));
00475 }
00476
00477
00478 moving_img->set_attr("xform.align2d", &t);
00479
00480 return moving_img;
00481 }
00482
00483 EMData *RotateTranslateAlignerPawel::align(EMData * this_img, EMData *to,
00484 const string & cmp_name, const Dict& cmp_params) const
00485 {
00486 if (cmp_name != "dot" && cmp_name != "ccc") throw InvalidParameterException("Resample aligner only works for dot and ccc");
00487
00488 int maxtx = params.set_default("tx", 0);
00489 int maxty = params.set_default("ty", 0);
00490 int r1 = params.set_default("r1",-1);
00491 int r2 = params.set_default("r2",-1);
00492
00493 if(this_img->get_xsize()/2 - 1 - r2 - maxtx <= 0 || (r2 == -1 && maxtx > 0)) throw InvalidParameterException("nx/2 - 1 - r2 - tx must be greater than or = 0");
00494 if(this_img->get_ysize()/2 - 1 - r2 - maxty <= 0 || (r2 == -1 && maxty > 0)) throw InvalidParameterException("ny/2 - 1 - r2 - ty must be greater than or = 0");
00495
00496 float best_peak = -numeric_limits<float>::infinity();
00497 int best_peak_index = 0;
00498 int best_tx = 0;
00499 int best_ty = 0;
00500 int polarxsize = 0;
00501
00502 for(int x = -maxtx; x <= maxtx; x++){
00503 for(int y = -maxty; y <= maxty; y++){
00504
00505 EMData * to_polar = to->unwrap(r1,r2,-1,0,0,true);
00506 EMData * this_img_polar = this_img->unwrap(r1,r2,-1,x,y,true);
00507 EMData * cf = this_img_polar->calc_ccfx(to_polar, 0, this_img_polar->get_ysize());
00508
00509 polarxsize = this_img_polar->get_xsize();
00510
00511
00512 delete to_polar; to_polar = 0;
00513 delete this_img_polar; this_img_polar = 0;
00514
00515 float *data = cf->get_data();
00516 float peak = 0;
00517 int peak_index = 0;
00518 Util::find_max(data, polarxsize, &peak, &peak_index);
00519 delete cf; cf = 0;
00520
00521 if(peak > best_peak) {
00522 best_peak = peak;
00523 best_peak_index = peak_index;
00524 best_tx = x;
00525 best_ty = y;
00526 }
00527 }
00528 }
00529
00530 float rot_angle = (float) (best_peak_index * 360.0f / polarxsize);
00531
00532
00533 Transform tmp(Dict("type","2d","alpha",rot_angle,"tx",best_tx,"ty",best_ty));
00534 EMData* rotimg=this_img->process("xform",Dict("transform",(Transform*)&tmp));
00535 rotimg->set_attr("xform.align2d",&tmp);
00536
00537 return rotimg;
00538
00539 }
00540
00541 EMData *RotateTranslateAligner::align(EMData * this_img, EMData *to,
00542 const string & cmp_name, const Dict& cmp_params) const
00543 {
00544
00545 #ifdef EMAN2_USING_CUDA
00546 if(EMData::usecuda == 1) {
00547 if(!this_img->cudarwdata) this_img->copy_to_cuda();
00548 if(!to->cudarwdata) to->copy_to_cuda();
00549 }
00550 #endif
00551
00552
00553 int rfp_mode = params.set_default("rfp_mode",0);
00554 EMData *rot_align = RotationalAligner::align_180_ambiguous(this_img,to,rfp_mode);
00555 Transform * tmp = rot_align->get_attr("xform.align2d");
00556 Dict rot = tmp->get_rotation("2d");
00557 float rotate_angle_solution = rot["alpha"];
00558 delete tmp;
00559
00560 EMData *rot_align_180 = rot_align->copy();
00561 rot_align_180->process_inplace("math.rotate.180");
00562
00563 Dict trans_params;
00564 trans_params["intonly"] = 0;
00565 trans_params["maxshift"] = params.set_default("maxshift", -1);
00566 trans_params["useflcf"]=params.set_default("useflcf",0);
00567
00568
00569 trans_params["nozero"] = params.set_default("nozero",false);
00570 EMData* rot_trans = rot_align->align("translational", to, trans_params, cmp_name, cmp_params);
00571 if( rot_align ) {
00572 delete rot_align;
00573 rot_align = 0;
00574 }
00575
00576
00577 EMData* rot_180_trans = rot_align_180->align("translational", to, trans_params, cmp_name, cmp_params);
00578 if( rot_align_180 ) {
00579 delete rot_align_180;
00580 rot_align_180 = 0;
00581 }
00582
00583
00584 float cmp1 = rot_trans->cmp(cmp_name, to, cmp_params);
00585 float cmp2 = rot_180_trans->cmp(cmp_name, to, cmp_params);
00586
00587 EMData *result = 0;
00588 if (cmp1 < cmp2) {
00589 if( rot_180_trans ) {
00590 delete rot_180_trans;
00591 rot_180_trans = 0;
00592 }
00593 result = rot_trans;
00594 }
00595 else {
00596 if( rot_trans ) {
00597 delete rot_trans;
00598 rot_trans = 0;
00599 }
00600 result = rot_180_trans;
00601 rotate_angle_solution -= 180.f;
00602 }
00603
00604 Transform* t = result->get_attr("xform.align2d");
00605 t->set_rotation(Dict("type","2d","alpha",rotate_angle_solution));
00606 result->set_attr("xform.align2d",t);
00607 delete t;
00608
00609 return result;
00610 }
00611
00612
00613
00614
00615 EMData* RotateTranslateFlipAligner::align(EMData * this_img, EMData *to,
00616 const string & cmp_name, const Dict& cmp_params) const
00617 {
00618
00619 Dict rt_params("maxshift", params["maxshift"], "rfp_mode", params.set_default("rfp_mode",0),"useflcf",params.set_default("useflcf",0));
00620 EMData *rot_trans_align = this_img->align("rotate_translate",to,rt_params,cmp_name, cmp_params);
00621
00622
00623 EMData *flipped = params.set_default("flip", (EMData *) 0);
00624 bool delete_flag = false;
00625 if (flipped == 0) {
00626 flipped = to->process("xform.flip", Dict("axis", "x"));
00627 delete_flag = true;
00628 }
00629
00630 EMData * rot_trans_align_flip = this_img->align("rotate_translate", flipped, rt_params, cmp_name, cmp_params);
00631 Transform* t = rot_trans_align_flip->get_attr("xform.align2d");
00632 t->set_mirror(true);
00633 rot_trans_align_flip->set_attr("xform.align2d",t);
00634 delete t;
00635
00636
00637 float cmp1 = rot_trans_align->cmp(cmp_name, to, cmp_params);
00638 float cmp2 = rot_trans_align_flip->cmp(cmp_name, flipped, cmp_params);
00639
00640 if (delete_flag){
00641 if(flipped) {
00642 delete flipped;
00643 flipped = 0;
00644 }
00645 }
00646
00647 EMData *result = 0;
00648 if (cmp1 < cmp2 ) {
00649
00650 if( rot_trans_align_flip ) {
00651 delete rot_trans_align_flip;
00652 rot_trans_align_flip = 0;
00653 }
00654 result = rot_trans_align;
00655 }
00656 else {
00657 if( rot_trans_align ) {
00658 delete rot_trans_align;
00659 rot_trans_align = 0;
00660 }
00661 result = rot_trans_align_flip;
00662 result->process_inplace("xform.flip",Dict("axis","x"));
00663 }
00664
00665 return result;
00666 }
00667
00668 EMData* RotateTranslateFlipAlignerIterative::align(EMData * this_img, EMData *to,
00669 const string & cmp_name, const Dict& cmp_params) const
00670 {
00671
00672 Dict rt_params("maxshift", params["maxshift"],"r1",params.set_default("r1",-1),"r2",params.set_default("r2",-1));
00673 EMData *rot_trans_align = this_img->align("rotate_translate.iterative",to,rt_params,cmp_name, cmp_params);
00674
00675
00676 EMData *flipped = params.set_default("flip", (EMData *) 0);
00677 bool delete_flag = false;
00678 if (flipped == 0) {
00679 flipped = to->process("xform.flip", Dict("axis", "x"));
00680 delete_flag = true;
00681 }
00682
00683 EMData * rot_trans_align_flip = this_img->align("rotate_translate.iterative", flipped, rt_params, cmp_name, cmp_params);
00684 Transform* t = rot_trans_align_flip->get_attr("xform.align2d");
00685 t->set_mirror(true);
00686 rot_trans_align_flip->set_attr("xform.align2d",t);
00687 delete t;
00688
00689
00690 float cmp1 = rot_trans_align->cmp(cmp_name, to, cmp_params);
00691 float cmp2 = rot_trans_align_flip->cmp(cmp_name, flipped, cmp_params);
00692
00693 if (delete_flag){
00694 if(flipped) {
00695 delete flipped;
00696 flipped = 0;
00697 }
00698 }
00699
00700 EMData *result = 0;
00701 if (cmp1 < cmp2 ) {
00702
00703 if( rot_trans_align_flip ) {
00704 delete rot_trans_align_flip;
00705 rot_trans_align_flip = 0;
00706 }
00707 result = rot_trans_align;
00708 }
00709 else {
00710 if( rot_trans_align ) {
00711 delete rot_trans_align;
00712 rot_trans_align = 0;
00713 }
00714 result = rot_trans_align_flip;
00715 result->process_inplace("xform.flip",Dict("axis","x"));
00716 }
00717
00718 return result;
00719 }
00720
00721 EMData *RotateTranslateFlipAlignerPawel::align(EMData * this_img, EMData *to,
00722 const string & cmp_name, const Dict& cmp_params) const
00723 {
00724 if (cmp_name != "dot" && cmp_name != "ccc") throw InvalidParameterException("Resample aligner only works for dot and ccc");
00725
00726 int maxtx = params.set_default("tx", 0);
00727 int maxty = params.set_default("ty", 0);
00728 int r1 = params.set_default("r1",-1);
00729 int r2 = params.set_default("r2",-1);
00730
00731 if(this_img->get_xsize()/2 - 1 - r2 - maxtx <= 0 || (r2 == -1 && maxtx > 0)) throw InvalidParameterException("nx/2 - 1 - r2 - tx must be greater than or = 0");
00732 if(this_img->get_ysize()/2 - 1 - r2 - maxty <= 0 || (r2 == -1 && maxty > 0)) throw InvalidParameterException("ny/2 - 1 - r2 - ty must be greater than or = 0");
00733
00734 float best_peak = -numeric_limits<float>::infinity();
00735 int best_peak_index = 0;
00736 int best_tx = 0;
00737 int best_ty = 0;
00738 int polarxsize = 0;
00739 bool flip = false;
00740
00741 for(int x = -maxtx; x <= maxtx; x++){
00742 for(int y = -maxty; y <= maxty; y++){
00743
00744 EMData * to_polar = to->unwrap(r1,r2,-1,0,0,true);
00745 EMData * this_img_polar = this_img->unwrap(r1,r2,-1,x,y,true);
00746 EMData * cfflip = this_img_polar->calc_ccfx(to_polar, 0, this_img_polar->get_ysize(), false, true);
00747 EMData * cf = this_img_polar->calc_ccfx(to_polar, 0, this_img_polar->get_ysize());
00748
00749 polarxsize = this_img_polar->get_xsize();
00750
00751
00752 delete to_polar; to_polar = 0;
00753 delete this_img_polar; this_img_polar = 0;
00754
00755 float *data = cf->get_data();
00756 float peak = 0;
00757 int peak_index = 0;
00758 Util::find_max(data, polarxsize, &peak, &peak_index);
00759 delete cf; cf = 0;
00760
00761 if(peak > best_peak) {
00762 best_peak = peak;
00763 best_peak_index = peak_index;
00764 best_tx = x;
00765 best_ty = y;
00766 flip = false;
00767 }
00768
00769 data = cfflip->get_data();
00770 Util::find_max(data, polarxsize, &peak, &peak_index);
00771 delete cfflip; cfflip = 0;
00772
00773 if(peak > best_peak) {
00774 best_peak = peak;
00775 best_peak_index = peak_index;
00776 best_tx = x;
00777 best_ty = y;
00778 flip = true;
00779 }
00780 }
00781 }
00782
00783 float rot_angle = (float) (best_peak_index * 360.0f / polarxsize);
00784
00785
00786 Transform tmp(Dict("type","2d","alpha",rot_angle,"tx",best_tx,"ty",best_ty));
00787 EMData* rotimg=this_img->process("xform",Dict("transform",(Transform*)&tmp));
00788 rotimg->set_attr("xform.align2d",&tmp);
00789 if(flip == true) {
00790 rotimg->process_inplace("xform.flip",Dict("axis", "x"));
00791 }
00792
00793 return rotimg;
00794
00795 }
00796
00797 EMData *RotateFlipAligner::align(EMData * this_img, EMData *to,
00798 const string& cmp_name, const Dict& cmp_params) const
00799 {
00800 Dict rot_params("rfp_mode",params.set_default("rfp_mode",0));
00801 EMData *r1 = this_img->align("rotational", to, rot_params,cmp_name, cmp_params);
00802
00803
00804 EMData* flipped =to->process("xform.flip", Dict("axis", "x"));
00805 EMData *r2 = this_img->align("rotational", flipped,rot_params, cmp_name, cmp_params);
00806 Transform* t = r2->get_attr("xform.align2d");
00807 t->set_mirror(true);
00808 r2->set_attr("xform.align2d",t);
00809 delete t;
00810
00811 float cmp1 = r1->cmp(cmp_name, to, cmp_params);
00812 float cmp2 = r2->cmp(cmp_name, flipped, cmp_params);
00813
00814 delete flipped; flipped = 0;
00815
00816 EMData *result = 0;
00817
00818 if (cmp1 < cmp2) {
00819 if( r2 )
00820 {
00821 delete r2;
00822 r2 = 0;
00823 }
00824 result = r1;
00825 }
00826 else {
00827 if( r1 )
00828 {
00829 delete r1;
00830 r1 = 0;
00831 }
00832 result = r2;
00833 result->process_inplace("xform.flip",Dict("axis","x"));
00834 }
00835
00836 return result;
00837 }
00838
00839 EMData *RotateFlipAlignerIterative::align(EMData * this_img, EMData *to,
00840 const string& cmp_name, const Dict& cmp_params) const
00841 {
00842 Dict rot_params("r1",params.set_default("r1",-1),"r2",params.set_default("r2",-1));
00843 EMData *r1 = this_img->align("rotational.iterative", to, rot_params,cmp_name, cmp_params);
00844
00845 EMData* flipped =to->process("xform.flip", Dict("axis", "x"));
00846 EMData *r2 = this_img->align("rotational.iterative", flipped,rot_params, cmp_name, cmp_params);
00847 Transform* t = r2->get_attr("xform.align2d");
00848 t->set_mirror(true);
00849 r2->set_attr("xform.align2d",t);
00850 delete t;
00851
00852 float cmp1 = r1->cmp(cmp_name, to, cmp_params);
00853 float cmp2 = r2->cmp(cmp_name, flipped, cmp_params);
00854
00855 delete flipped; flipped = 0;
00856
00857 EMData *result = 0;
00858
00859 if (cmp1 < cmp2) {
00860 if( r2 )
00861 {
00862 delete r2;
00863 r2 = 0;
00864 }
00865 result = r1;
00866 }
00867 else {
00868 if( r1 )
00869 {
00870 delete r1;
00871 r1 = 0;
00872 }
00873 result = r2;
00874 result->process_inplace("xform.flip",Dict("axis","x"));
00875 }
00876
00877 return result;
00878 }
00879
00880
00881
00882
00883
00884
00885
00886
00887
00888
00889
00890
00891
00892
00893
00894 EMData *RTFExhaustiveAligner::align(EMData * this_img, EMData *to,
00895 const string & cmp_name, const Dict& cmp_params) const
00896 {
00897 EMData *flip = params.set_default("flip", (EMData *) 0);
00898 int maxshift = params.set_default("maxshift", this_img->get_xsize()/8);
00899 if (maxshift < 2) throw InvalidParameterException("maxshift must be greater than or equal to 2");
00900
00901 int ny = this_img->get_ysize();
00902 int xst = (int) floor(2 * M_PI * ny);
00903 xst = Util::calc_best_fft_size(xst);
00904
00905 Dict d("n",2);
00906 EMData *to_shrunk_unwrapped = to->process("math.medianshrink",d);
00907
00908 int to_copy_r2 = to_shrunk_unwrapped->get_ysize() / 2 - 2 - maxshift / 2;
00909 EMData *tmp = to_shrunk_unwrapped->unwrap(4, to_copy_r2, xst / 2, 0, 0, true);
00910 if( to_shrunk_unwrapped )
00911 {
00912 delete to_shrunk_unwrapped;
00913 to_shrunk_unwrapped = 0;
00914 }
00915 to_shrunk_unwrapped = tmp;
00916
00917 EMData *to_shrunk_unwrapped_copy = to_shrunk_unwrapped->copy();
00918 EMData* to_unwrapped = to->unwrap(4, to->get_ysize() / 2 - 2 - maxshift, xst, 0, 0, true);
00919 EMData *to_unwrapped_copy = to_unwrapped->copy();
00920
00921 bool delete_flipped = true;
00922 EMData *flipped = 0;
00923 if (flip) {
00924 delete_flipped = false;
00925 flipped = flip;
00926 }
00927 else {
00928 flipped = to->process("xform.flip", Dict("axis", "x"));
00929 }
00930 EMData *to_shrunk_flipped_unwrapped = flipped->process("math.medianshrink",d);
00931 tmp = to_shrunk_flipped_unwrapped->unwrap(4, to_copy_r2, xst / 2, 0, 0, true);
00932 if( to_shrunk_flipped_unwrapped )
00933 {
00934 delete to_shrunk_flipped_unwrapped;
00935 to_shrunk_flipped_unwrapped = 0;
00936 }
00937 to_shrunk_flipped_unwrapped = tmp;
00938 EMData *to_shrunk_flipped_unwrapped_copy = to_shrunk_flipped_unwrapped->copy();
00939 EMData* to_flip_unwrapped = flipped->unwrap(4, to->get_ysize() / 2 - 2 - maxshift, xst, 0, 0, true);
00940 EMData* to_flip_unwrapped_copy = to_flip_unwrapped->copy();
00941
00942 if (delete_flipped && flipped != 0) {
00943 delete flipped;
00944 flipped = 0;
00945 }
00946
00947 EMData *this_shrunk_2 = this_img->process("math.medianshrink",d);
00948
00949 float bestval = FLT_MAX;
00950 float bestang = 0;
00951 int bestflip = 0;
00952 float bestdx = 0;
00953 float bestdy = 0;
00954
00955 int half_maxshift = maxshift / 2;
00956
00957 int ur2 = this_shrunk_2->get_ysize() / 2 - 2 - half_maxshift;
00958 for (int dy = -half_maxshift; dy <= half_maxshift; dy += 1) {
00959 for (int dx = -half_maxshift; dx <= half_maxshift; dx += 1) {
00960 #ifdef _WIN32
00961 if (_hypot(dx, dy) <= half_maxshift) {
00962 #else
00963 if (hypot(dx, dy) <= half_maxshift) {
00964 #endif
00965 EMData *uw = this_shrunk_2->unwrap(4, ur2, xst / 2, dx, dy, true);
00966 EMData *uwc = uw->copy();
00967 EMData *a = uw->calc_ccfx(to_shrunk_unwrapped);
00968
00969 uwc->rotate_x(a->calc_max_index());
00970 float cm = uwc->cmp(cmp_name, to_shrunk_unwrapped_copy, cmp_params);
00971 if (cm < bestval) {
00972 bestval = cm;
00973 bestang = (float) (2.0 * M_PI * a->calc_max_index() / a->get_xsize());
00974 bestdx = (float)dx;
00975 bestdy = (float)dy;
00976 bestflip = 0;
00977 }
00978
00979
00980 if( a )
00981 {
00982 delete a;
00983 a = 0;
00984 }
00985 if( uw )
00986 {
00987 delete uw;
00988 uw = 0;
00989 }
00990 if( uwc )
00991 {
00992 delete uwc;
00993 uwc = 0;
00994 }
00995 uw = this_shrunk_2->unwrap(4, ur2, xst / 2, dx, dy, true);
00996 uwc = uw->copy();
00997 a = uw->calc_ccfx(to_shrunk_flipped_unwrapped);
00998
00999 uwc->rotate_x(a->calc_max_index());
01000 cm = uwc->cmp(cmp_name, to_shrunk_flipped_unwrapped_copy, cmp_params);
01001 if (cm < bestval) {
01002 bestval = cm;
01003 bestang = (float) (2.0 * M_PI * a->calc_max_index() / a->get_xsize());
01004 bestdx = (float)dx;
01005 bestdy = (float)dy;
01006 bestflip = 1;
01007 }
01008
01009 if( a )
01010 {
01011 delete a;
01012 a = 0;
01013 }
01014
01015 if( uw )
01016 {
01017 delete uw;
01018 uw = 0;
01019 }
01020 if( uwc )
01021 {
01022 delete uwc;
01023 uwc = 0;
01024 }
01025 }
01026 }
01027 }
01028 if( this_shrunk_2 )
01029 {
01030 delete this_shrunk_2;
01031 this_shrunk_2 = 0;
01032 }
01033 if( to_shrunk_unwrapped )
01034 {
01035 delete to_shrunk_unwrapped;
01036 to_shrunk_unwrapped = 0;
01037 }
01038 if( to_shrunk_unwrapped_copy )
01039 {
01040 delete to_shrunk_unwrapped_copy;
01041 to_shrunk_unwrapped_copy = 0;
01042 }
01043 if( to_shrunk_flipped_unwrapped )
01044 {
01045 delete to_shrunk_flipped_unwrapped;
01046 to_shrunk_flipped_unwrapped = 0;
01047 }
01048 if( to_shrunk_flipped_unwrapped_copy )
01049 {
01050 delete to_shrunk_flipped_unwrapped_copy;
01051 to_shrunk_flipped_unwrapped_copy = 0;
01052 }
01053 bestdx *= 2;
01054 bestdy *= 2;
01055 bestval = FLT_MAX;
01056
01057 float bestdx2 = bestdx;
01058 float bestdy2 = bestdy;
01059
01060
01061
01062 for (float dy = bestdy2 - 3; dy <= bestdy2 + 3; dy += 1.0 ) {
01063 for (float dx = bestdx2 - 3; dx <= bestdx2 + 3; dx += 1.0 ) {
01064
01065 #ifdef _WIN32
01066 if (_hypot(dx, dy) <= maxshift) {
01067 #else
01068 if (hypot(dx, dy) <= maxshift) {
01069 #endif
01070 EMData *uw = this_img->unwrap(4, this_img->get_ysize() / 2 - 2 - maxshift, xst, (int)dx, (int)dy, true);
01071 EMData *uwc = uw->copy();
01072 EMData *a = uw->calc_ccfx(to_unwrapped);
01073
01074 uwc->rotate_x(a->calc_max_index());
01075 float cm = uwc->cmp(cmp_name, to_unwrapped_copy, cmp_params);
01076
01077 if (cm < bestval) {
01078 bestval = cm;
01079 bestang = (float)(2.0 * M_PI * a->calc_max_index() / a->get_xsize());
01080 bestdx = dx;
01081 bestdy = dy;
01082 bestflip = 0;
01083 }
01084
01085 if( a )
01086 {
01087 delete a;
01088 a = 0;
01089 }
01090 if( uw )
01091 {
01092 delete uw;
01093 uw = 0;
01094 }
01095 if( uwc )
01096 {
01097 delete uwc;
01098 uwc = 0;
01099 }
01100 uw = this_img->unwrap(4, this_img->get_ysize() / 2 - 2 - maxshift, xst, (int)dx, (int)dy, true);
01101 uwc = uw->copy();
01102 a = uw->calc_ccfx(to_flip_unwrapped);
01103
01104 uwc->rotate_x(a->calc_max_index());
01105 cm = uwc->cmp(cmp_name, to_flip_unwrapped_copy, cmp_params);
01106
01107 if (cm < bestval) {
01108 bestval = cm;
01109 bestang = (float)(2.0 * M_PI * a->calc_max_index() / a->get_xsize());
01110 bestdx = dx;
01111 bestdy = dy;
01112 bestflip = 1;
01113 }
01114
01115 if( a )
01116 {
01117 delete a;
01118 a = 0;
01119 }
01120 if( uw )
01121 {
01122 delete uw;
01123 uw = 0;
01124 }
01125 if( uwc )
01126 {
01127 delete uwc;
01128 uwc = 0;
01129 }
01130 }
01131 }
01132 }
01133 if( to_unwrapped ) {delete to_unwrapped;to_unwrapped = 0;}
01134 if( to_shrunk_unwrapped ) { delete to_shrunk_unwrapped; to_shrunk_unwrapped = 0;}
01135 if (to_unwrapped_copy) { delete to_unwrapped_copy; to_unwrapped_copy = 0; }
01136 if (to_flip_unwrapped) { delete to_flip_unwrapped; to_flip_unwrapped = 0; }
01137 if (to_flip_unwrapped_copy) { delete to_flip_unwrapped_copy; to_flip_unwrapped_copy = 0;}
01138
01139 bestang *= (float)EMConsts::rad2deg;
01140 Transform t(Dict("type","2d","alpha",(float)bestang));
01141 t.set_pre_trans(Vec2f(-bestdx,-bestdy));
01142 if (bestflip) {
01143 t.set_mirror(true);
01144 }
01145
01146 EMData* ret = this_img->process("xform",Dict("transform",&t));
01147 ret->set_attr("xform.align2d",&t);
01148
01149 return ret;
01150 }
01151
01152
01153 EMData *RTFSlowExhaustiveAligner::align(EMData * this_img, EMData *to,
01154 const string & cmp_name, const Dict& cmp_params) const
01155 {
01156
01157 EMData *flip = params.set_default("flip", (EMData *) 0);
01158 int maxshift = params.set_default("maxshift", -1);
01159
01160 EMData *flipped = 0;
01161
01162 bool delete_flipped = true;
01163 if (flip) {
01164 delete_flipped = false;
01165 flipped = flip;
01166 }
01167 else {
01168 flipped = to->process("xform.flip", Dict("axis", "x"));
01169 }
01170
01171 int nx = this_img->get_xsize();
01172
01173 if (maxshift < 0) {
01174 maxshift = nx / 10;
01175 }
01176
01177 float angle_step = params.set_default("angstep", 0.0f);
01178 if ( angle_step == 0 ) angle_step = atan2(2.0f, (float)nx);
01179 else {
01180 angle_step *= (float)EMConsts::deg2rad;
01181 }
01182 float trans_step = params.set_default("transtep",1.0f);
01183
01184 if (trans_step <= 0) throw InvalidParameterException("transstep must be greater than 0");
01185 if (angle_step <= 0) throw InvalidParameterException("angstep must be greater than 0");
01186
01187
01188 Dict shrinkfactor("n",2);
01189 EMData *this_img_shrink = this_img->process("math.medianshrink",shrinkfactor);
01190 EMData *to_shrunk = to->process("math.medianshrink",shrinkfactor);
01191 EMData *flipped_shrunk = flipped->process("math.medianshrink",shrinkfactor);
01192
01193 int bestflip = 0;
01194 float bestdx = 0;
01195 float bestdy = 0;
01196
01197 float bestang = 0;
01198 float bestval = FLT_MAX;
01199
01200 int half_maxshift = maxshift / 2;
01201
01202
01203 for (int dy = -half_maxshift; dy <= half_maxshift; ++dy) {
01204 for (int dx = -half_maxshift; dx <= half_maxshift; ++dx) {
01205 if (hypot(dx, dy) <= maxshift) {
01206 for (float ang = -angle_step * 2.0f; ang <= (float)2 * M_PI; ang += angle_step * 4.0f) {
01207 EMData v(*this_img_shrink);
01208 Transform t(Dict("type","2d","alpha",static_cast<float>(ang*EMConsts::rad2deg)));
01209 t.set_trans((float)dx,(float)dy);
01210 v.transform(t);
01211
01212
01213 float lc = v.cmp(cmp_name, to_shrunk, cmp_params);
01214
01215 if (lc < bestval) {
01216 bestval = lc;
01217 bestang = ang;
01218 bestdx = (float)dx;
01219 bestdy = (float)dy;
01220 bestflip = 0;
01221 }
01222
01223 lc = v.cmp(cmp_name,flipped_shrunk , cmp_params);
01224 if (lc < bestval) {
01225 bestval = lc;
01226 bestang = ang;
01227 bestdx = (float)dx;
01228 bestdy = (float)dy;
01229 bestflip = 1;
01230 }
01231 }
01232 }
01233 }
01234 }
01235
01236 if( to_shrunk )
01237 {
01238 delete to_shrunk;
01239 to_shrunk = 0;
01240 }
01241 if( flipped_shrunk )
01242 {
01243 delete flipped_shrunk;
01244 flipped_shrunk = 0;
01245 }
01246 if( this_img_shrink )
01247 {
01248 delete this_img_shrink;
01249 this_img_shrink = 0;
01250 }
01251
01252 bestdx *= 2;
01253 bestdy *= 2;
01254 bestval = FLT_MAX;
01255
01256 float bestdx2 = bestdx;
01257 float bestdy2 = bestdy;
01258 float bestang2 = bestang;
01259
01260 for (float dy = bestdy2 - 3; dy <= bestdy2 + 3; dy += trans_step) {
01261 for (float dx = bestdx2 - 3; dx <= bestdx2 + 3; dx += trans_step) {
01262 if (hypot(dx, dy) <= maxshift) {
01263 for (float ang = bestang2 - angle_step * 6.0f; ang <= bestang2 + angle_step * 6.0f; ang += angle_step) {
01264 EMData v(*this_img);
01265 Transform t(Dict("type","2d","alpha",static_cast<float>(ang*EMConsts::rad2deg)));
01266 t.set_trans(dx,dy);
01267 v.transform(t);
01268
01269
01270 float lc = v.cmp(cmp_name, to, cmp_params);
01271
01272 if (lc < bestval) {
01273 bestval = lc;
01274 bestang = ang;
01275 bestdx = dx;
01276 bestdy = dy;
01277 bestflip = 0;
01278 }
01279
01280 lc = v.cmp(cmp_name, flipped, cmp_params);
01281
01282 if (lc < bestval) {
01283 bestval = lc;
01284 bestang = ang;
01285 bestdx = dx;
01286 bestdy = dy;
01287 bestflip = 1;
01288 }
01289 }
01290 }
01291 }
01292 }
01293
01294 if (delete_flipped) { delete flipped; flipped = 0; }
01295
01296 bestang *= (float)EMConsts::rad2deg;
01297 Transform t(Dict("type","2d","alpha",(float)bestang));
01298 t.set_trans(bestdx,bestdy);
01299
01300 if (bestflip) {
01301 t.set_mirror(true);
01302 }
01303
01304 EMData* rslt = this_img->process("xform",Dict("transform",&t));
01305 rslt->set_attr("xform.align2d",&t);
01306
01307 return rslt;
01308 }
01309
01310
01311
01312 static double refalifn(const gsl_vector * v, void *params)
01313 {
01314 Dict *dict = (Dict *) params;
01315
01316 double x = gsl_vector_get(v, 0);
01317 double y = gsl_vector_get(v, 1);
01318 double a = gsl_vector_get(v, 2);
01319
01320 EMData *this_img = (*dict)["this"];
01321 EMData *with = (*dict)["with"];
01322 bool mirror = (*dict)["mirror"];
01323
01324
01325
01326
01327
01328
01329 Transform t(Dict("type","2d","alpha",static_cast<float>(a)));
01330
01331
01332
01333 t.set_trans((float)x,(float)y);
01334 t.set_mirror(mirror);
01335 if (v->size>3) {
01336 float sca=(float)gsl_vector_get(v, 3);
01337 if (sca<.7 || sca>1.3) return 1.0e20;
01338 t.set_scale((float)gsl_vector_get(v, 3));
01339 }
01340 EMData *tmp = this_img->process("xform",Dict("transform",&t));
01341
01342
01343 Cmp* c = (Cmp*) ((void*)(*dict)["cmp"]);
01344 double result = c->cmp(tmp,with);
01345
01346
01347
01348
01349
01350
01351
01352
01353
01354
01355
01356
01357
01358
01360
01361
01362
01363
01364
01365
01366 if ( tmp != 0 ) delete tmp;
01367
01368 return result;
01369 }
01370
01371 static double refalifnfast(const gsl_vector * v, void *params)
01372 {
01373 Dict *dict = (Dict *) params;
01374 EMData *this_img = (*dict)["this"];
01375 EMData *img_to = (*dict)["with"];
01376 bool mirror = (*dict)["mirror"];
01377
01378 double x = gsl_vector_get(v, 0);
01379 double y = gsl_vector_get(v, 1);
01380 double a = gsl_vector_get(v, 2);
01381
01382 double r = this_img->dot_rotate_translate(img_to, (float)x, (float)y, (float)a, mirror);
01383 int nsec = this_img->get_xsize() * this_img->get_ysize();
01384 double result = 1.0 - r / nsec;
01385
01386
01387 return result;
01388 }
01389
01390
01391 EMData *RefineAligner::align(EMData * this_img, EMData *to,
01392 const string & cmp_name, const Dict& cmp_params) const
01393 {
01394
01395 if (!to) {
01396 return 0;
01397 }
01398
01399 EMData *result;
01400 int mode = params.set_default("mode", 0);
01401 float saz = 0.0;
01402 float sdx = 0.0;
01403 float sdy = 0.0;
01404 float sscale = 1.0;
01405 bool mirror = false;
01406 Transform* t;
01407 if (params.has_key("xform.align2d") ) {
01408 t = params["xform.align2d"];
01409 Dict params = t->get_params("2d");
01410 saz = params["alpha"];
01411 sdx = params["tx"];
01412 sdy = params["ty"];
01413 mirror = params["mirror"];
01414 sscale = params["scale"];
01415 } else {
01416 t = new Transform();
01417 }
01418
01419
01420 if ((float)(this_img->get_attr("sigma"))==0.0 || (float)(to->get_attr("sigma"))==0.0) {
01421 result = this_img->process("xform",Dict("transform",t));
01422 result->set_attr("xform.align2d",t);
01423 delete t;
01424 return result;
01425 }
01426
01427 float stepx = params.set_default("stepx",1.0f);
01428 float stepy = params.set_default("stepy",1.0f);
01429
01430 float stepaz = params.set_default("stepaz",5.0f);
01431 float stepscale = params.set_default("stepscale",0.0f);
01432
01433 int np = 3;
01434 if (stepscale!=0.0) np++;
01435 Dict gsl_params;
01436 gsl_params["this"] = this_img;
01437 gsl_params["with"] = to;
01438 gsl_params["snr"] = params["snr"];
01439 gsl_params["mirror"] = mirror;
01440
01441 const gsl_multimin_fminimizer_type *T = gsl_multimin_fminimizer_nmsimplex;
01442 gsl_vector *ss = gsl_vector_alloc(np);
01443
01444
01445 gsl_vector_set(ss, 0, stepx);
01446 gsl_vector_set(ss, 1, stepy);
01447 gsl_vector_set(ss, 2, stepaz);
01448 if (stepscale!=0.0) gsl_vector_set(ss,3,stepscale);
01449
01450 gsl_vector *x = gsl_vector_alloc(np);
01451 gsl_vector_set(x, 0, sdx);
01452 gsl_vector_set(x, 1, sdy);
01453 gsl_vector_set(x, 2, saz);
01454 if (stepscale!=0.0) gsl_vector_set(x,3,1.0);
01455
01456 Cmp *c = 0;
01457
01458 gsl_multimin_function minex_func;
01459 if (mode == 2) {
01460 minex_func.f = &refalifnfast;
01461 }
01462 else {
01463 c = Factory < Cmp >::get(cmp_name, cmp_params);
01464 gsl_params["cmp"] = (void *) c;
01465 minex_func.f = &refalifn;
01466 }
01467
01468 minex_func.n = np;
01469 minex_func.params = (void *) &gsl_params;
01470
01471 gsl_multimin_fminimizer *s = gsl_multimin_fminimizer_alloc(T, np);
01472 gsl_multimin_fminimizer_set(s, &minex_func, x, ss);
01473
01474 int rval = GSL_CONTINUE;
01475 int status = GSL_SUCCESS;
01476 int iter = 1;
01477
01478 float precision = params.set_default("precision",0.04f);
01479 int maxiter = params.set_default("maxiter",28);
01480
01481
01482
01483
01484 while (rval == GSL_CONTINUE && iter < maxiter) {
01485 iter++;
01486 status = gsl_multimin_fminimizer_iterate(s);
01487 if (status) {
01488 break;
01489 }
01490 rval = gsl_multimin_test_size(gsl_multimin_fminimizer_size(s), precision);
01491 }
01492
01493 int maxshift = params.set_default("maxshift",-1);
01494
01495 if (maxshift <= 0) {
01496 maxshift = this_img->get_xsize() / 4;
01497 }
01498 float fmaxshift = static_cast<float>(maxshift);
01499 if ( fmaxshift >= fabs((float)gsl_vector_get(s->x, 0)) && fmaxshift >= fabs((float)gsl_vector_get(s->x, 1)) && (stepscale==0 || (((float)gsl_vector_get(s->x, 3))<1.3 && ((float)gsl_vector_get(s->x, 3))<0.7)) )
01500 {
01501
01502 Transform tsoln(Dict("type","2d","alpha",(float)gsl_vector_get(s->x, 2)));
01503 tsoln.set_mirror(mirror);
01504 tsoln.set_trans((float)gsl_vector_get(s->x, 0),(float)gsl_vector_get(s->x, 1));
01505 if (stepscale!=0.0) tsoln.set_scale((float)gsl_vector_get(s->x, 3));
01506 result = this_img->process("xform",Dict("transform",&tsoln));
01507 result->set_attr("xform.align2d",&tsoln);
01508 } else {
01509
01510 result = this_img->process("xform",Dict("transform",t));
01511 result->set_attr("xform.align2d",t);
01512 }
01513
01514 delete t;
01515 t = 0;
01516
01517 gsl_vector_free(x);
01518 gsl_vector_free(ss);
01519 gsl_multimin_fminimizer_free(s);
01520
01521 if ( c != 0 ) delete c;
01522 return result;
01523 }
01524
01525 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)
01526 {
01527 Vec3f normal(n0,n1,n2);
01528 normal.normalize();
01529
01530 float omega = spincoeff*sqrt(n0*n0 + n1*n1 + n2*n2);
01531 Dict d;
01532 d["type"] = "spin";
01533 d["Omega"] = omega;
01534 d["n1"] = normal[0];
01535 d["n2"] = normal[1];
01536 d["n3"] = normal[2];
01537
01538
01539 Transform q(d);
01540 q.set_trans((float)x,(float)y,(float)z);
01541
01542 q = q*(*t);
01543
01544 return q;
01545 }
01546
01547 static double refalifn3dquat(const gsl_vector * v, void *params)
01548 {
01549 Dict *dict = (Dict *) params;
01550
01551 double n0 = gsl_vector_get(v, 0);
01552 double n1 = gsl_vector_get(v, 1);
01553 double n2 = gsl_vector_get(v, 2);
01554 double x = gsl_vector_get(v, 3);
01555 double y = gsl_vector_get(v, 4);
01556 double z = gsl_vector_get(v, 5);
01557
01558 EMData *this_img = (*dict)["this"];
01559 EMData *with = (*dict)["with"];
01560
01561
01562 Transform* t = (*dict)["transform"];
01563 float spincoeff = (*dict)["spincoeff"];
01564
01565 Transform soln = refalin3d_perturbquat(t,spincoeff,(float)n0,(float)n1,(float)n2,(float)x,(float)y,(float)z);
01566
01567 EMData *tmp = this_img->process("xform",Dict("transform",&soln));
01568 Cmp* c = (Cmp*) ((void*)(*dict)["cmp"]);
01569 double result = c->cmp(tmp,with);
01570 if ( tmp != 0 ) delete tmp;
01571 delete t; t = 0;
01572
01573 return result;
01574 }
01575
01576 EMData* Refine3DAlignerQuaternion::align(EMData * this_img, EMData *to,
01577 const string & cmp_name, const Dict& cmp_params) const
01578 {
01579
01580 if (!to || !this_img) throw NullPointerException("Input image is null");
01581
01582 if (to->get_ndim() != 3 || this_img->get_ndim() != 3) throw ImageDimensionException("The Refine3D aligner only works for 3D images");
01583
01584 #ifdef EMAN2_USING_CUDA
01585 if(EMData::usecuda == 1) {
01586 if(!this_img->cudarwdata) this_img->copy_to_cuda();
01587 if(!to->cudarwdata) to->copy_to_cuda();
01588 }
01589 #endif
01590
01591 float sdi = 0.0;
01592 float sdj = 0.0;
01593 float sdk = 0.0;
01594 float sdx = 0.0;
01595 float sdy = 0.0;
01596 float sdz = 0.0;
01597 bool mirror = false;
01598
01599 Transform* t;
01600 if (params.has_key("xform.align3d") ) {
01601
01602
01603
01604
01605 t = params["xform.align3d"];
01606 }else {
01607 t = new Transform();
01608 }
01609
01610 float spincoeff = params.set_default("spin_coeff",10.0f);
01611
01612 int np = 6;
01613 Dict gsl_params;
01614 gsl_params["this"] = this_img;
01615 gsl_params["with"] = to;
01616 gsl_params["snr"] = params["snr"];
01617 gsl_params["mirror"] = mirror;
01618 gsl_params["transform"] = t;
01619 gsl_params["spincoeff"] = spincoeff;
01620 Dict altered_cmp_params(cmp_params);
01621
01622 const gsl_multimin_fminimizer_type *T = gsl_multimin_fminimizer_nmsimplex;
01623 gsl_vector *ss = gsl_vector_alloc(np);
01624
01625 float stepi = params.set_default("stepi",1.0f);
01626 float stepj = params.set_default("stepj",1.0f);
01627 float stepk = params.set_default("stepk",1.0f);
01628 float stepx = params.set_default("stepx",1.0f);
01629 float stepy = params.set_default("stepy",1.0f);
01630 float stepz = params.set_default("stepz",1.0f);
01631
01632
01633 gsl_vector_set(ss, 0, stepi);
01634 gsl_vector_set(ss, 1, stepj);
01635 gsl_vector_set(ss, 2, stepk);
01636 gsl_vector_set(ss, 3, stepx);
01637 gsl_vector_set(ss, 4, stepy);
01638 gsl_vector_set(ss, 5, stepz);
01639
01640 gsl_vector *x = gsl_vector_alloc(np);
01641 gsl_vector_set(x, 0, sdi);
01642 gsl_vector_set(x, 1, sdj);
01643 gsl_vector_set(x, 2, sdk);
01644 gsl_vector_set(x, 3, sdx);
01645 gsl_vector_set(x, 4, sdy);
01646 gsl_vector_set(x, 5, sdz);
01647
01648 gsl_multimin_function minex_func;
01649 Cmp *c = Factory < Cmp >::get(cmp_name, altered_cmp_params);
01650
01651 gsl_params["cmp"] = (void *) c;
01652 minex_func.f = &refalifn3dquat;
01653
01654 minex_func.n = np;
01655 minex_func.params = (void *) &gsl_params;
01656
01657 gsl_multimin_fminimizer *s = gsl_multimin_fminimizer_alloc(T, np);
01658 gsl_multimin_fminimizer_set(s, &minex_func, x, ss);
01659
01660 int rval = GSL_CONTINUE;
01661 int status = GSL_SUCCESS;
01662 int iter = 1;
01663
01664 float precision = params.set_default("precision",0.01f);
01665 int maxiter = params.set_default("maxiter",100);
01666 while (rval == GSL_CONTINUE && iter < maxiter) {
01667 iter++;
01668 status = gsl_multimin_fminimizer_iterate(s);
01669 if (status) {
01670 break;
01671 }
01672 rval = gsl_multimin_test_size(gsl_multimin_fminimizer_size(s), precision);
01673 }
01674
01675 int maxshift = params.set_default("maxshift",-1);
01676
01677 if (maxshift <= 0) {
01678 maxshift = this_img->get_xsize() / 4;
01679 }
01680 float fmaxshift = static_cast<float>(maxshift);
01681
01682 EMData *result;
01683 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))
01684 {
01685 float n0 = (float)gsl_vector_get(s->x, 0);
01686 float n1 = (float)gsl_vector_get(s->x, 1);
01687 float n2 = (float)gsl_vector_get(s->x, 2);
01688 float x = (float)gsl_vector_get(s->x, 3);
01689 float y = (float)gsl_vector_get(s->x, 4);
01690 float z = (float)gsl_vector_get(s->x, 5);
01691
01692 Transform tsoln = refalin3d_perturbquat(t,spincoeff,n0,n1,n2,x,y,z);
01693
01694 result = this_img->process("xform",Dict("transform",&tsoln));
01695 result->set_attr("xform.align3d",&tsoln);
01696 result->set_attr("score", result->cmp(cmp_name,to,cmp_params));
01697
01698
01699 } else {
01700 result = this_img->process("xform",Dict("transform",t));
01701 result->set_attr("xform.align3d",t);
01702 }
01703
01704
01705 delete t;
01706 t = 0;
01707 gsl_vector_free(x);
01708 gsl_vector_free(ss);
01709 gsl_multimin_fminimizer_free(s);
01710
01711 if ( c != 0 ) delete c;
01712 return result;
01713 }
01714
01715 EMData*Refine3DAlignerGrid::align(EMData * this_img, EMData *to,
01716 const string & cmp_name, const Dict& cmp_params) const
01717 {
01718 if ( this_img->get_ndim() != 3 || to->get_ndim() != 3 ) {
01719 throw ImageDimensionException("This aligner only works for 3D images");
01720 }
01721
01722 Transform* t;
01723 if (params.has_key("xform.align3d") ) {
01724
01725
01726
01727
01728 t = params["xform.align3d"];
01729 }else {
01730 t = new Transform();
01731 }
01732
01733 int searchx = 0;
01734 int searchy = 0;
01735 int searchz = 0;
01736
01737 bool dotrans = params.set_default("dotrans",1);
01738 if (params.has_key("search")) {
01739 vector<string> check;
01740 check.push_back("searchx");
01741 check.push_back("searchy");
01742 check.push_back("searchz");
01743 for(vector<string>::const_iterator cit = check.begin(); cit != check.end(); ++cit) {
01744 if (params.has_key(*cit)) throw InvalidParameterException("The search parameter is mutually exclusive of the searchx, searchy, and searchz parameters");
01745 }
01746 int search = params["search"];
01747 searchx = search;
01748 searchy = search;
01749 searchz = search;
01750 } else {
01751 searchx = params.set_default("searchx",3);
01752 searchy = params.set_default("searchy",3);
01753 searchz = params.set_default("searchz",3);
01754 }
01755
01756 float delta = params.set_default("delta",1.0f);
01757 float range = params.set_default("range",10.0f);
01758 bool verbose = params.set_default("verbose",false);
01759
01760 bool tomography = (cmp_name == "ccc.tomo") ? 1 : 0;
01761 EMData * tofft = 0;
01762 if(dotrans || tomography){
01763 tofft = to->do_fft();
01764 }
01765
01766 #ifdef EMAN2_USING_CUDA
01767 if(EMData::usecuda == 1) {
01768 if(!this_img->isrodataongpu()) this_img->copy_to_cudaro();
01769 if(!to->cudarwdata) to->copy_to_cuda();
01770 if(to->cudarwdata){if(tofft) tofft->copy_to_cuda();}
01771 }
01772 #endif
01773
01774 Dict d;
01775 d["type"] = "eman";
01776 Dict best;
01777 best["score"] = 0.0f;
01778 bool use_cpu = true;
01779 Transform tran = Transform();
01780 for ( float alt = 0; alt < range; alt += delta) {
01781
01782 float saz = 360;
01783 if(alt != 0) saz = delta/sin(alt*M_PI/180.0f);
01784 for ( float az = 0; az < 360; az += saz ){
01785 if (verbose) {
01786 cout << "Trying angle alt " << alt << " az " << az << endl;
01787 }
01788
01789 for( float phi = -range-az; phi < range-az; phi += delta ) {
01790 d["alt"] = alt;
01791 d["phi"] = phi;
01792 d["az"] = az;
01793 Transform tr(d);
01794 tr = tr*(*t);
01795
01796 EMData* transformed = this_img->process("xform",Dict("transform",&tr));
01797
01798
01799 float score = 0.0f;
01800 if(dotrans || tomography){
01801 EMData* ccf = transformed->calc_ccf(tofft);
01802 #ifdef EMAN2_USING_CUDA
01803 if(to->cudarwdata){
01804 use_cpu = false;
01805 CudaPeakInfo* data = calc_max_location_wrap_cuda(ccf->cudarwdata, ccf->get_xsize(), ccf->get_ysize(), ccf->get_zsize(), searchx, searchy, searchz);
01806 tran.set_trans((float)-data->px, (float)-data->py, (float)-data->pz);
01807
01808
01809 tr = tran*tr;
01810 if (tomography) {
01811 float2 stats = get_stats_cuda(ccf->cudarwdata, ccf->get_xsize(), ccf->get_ysize(), ccf->get_zsize());
01812 score = -(data->peak - stats.x)/sqrt(stats.y);
01813 } else {
01814 score = -data->peak;
01815 }
01816 delete data;
01817 }
01818 #endif
01819 if(use_cpu){
01820 if(tomography) ccf->process_inplace("normalize");
01821
01822
01823
01824 IntPoint point = ccf->calc_max_location_wrap(searchx,searchy,searchz);
01825 tran.set_trans((float)-point[0], (float)-point[1], (float)-point[2]);
01826 score = -ccf->get_value_at_wrap(point[0], point[1], point[2]);
01827 tr = tran*tr;
01828
01829 }
01830 delete ccf; ccf =0;
01831 delete transformed; transformed = 0;
01832 }
01833
01834 if(!tomography){
01835 if(!transformed) transformed = this_img->process("xform",Dict("transform",&tr));
01836 score = transformed->cmp(cmp_name,to,cmp_params);
01837 delete transformed; transformed = 0;
01838 }
01839
01840 if(score < float(best["score"])) {
01841 best["score"] = score;
01842 best["xform.align3d"] = &tr;
01843 }
01844 }
01845 }
01846 }
01847
01848 if(tofft) {delete tofft; tofft = 0;}
01849
01850
01851 EMData* best_match = this_img->process("xform",Dict("transform", best["xform.align3d"]));
01852 best_match->set_attr("xform.align3d", best["xform.align3d"]);
01853 best_match->set_attr("score", float(best["score"]));
01854
01855
01856 Transform* zz = best_match->get_attr("xform.align3d");
01857 Vec3f zzz = zz->get_trans();
01858 cout << "x " << float(zzz[0]) << " y " << float(zzz[1]) << " z " << float(zzz[2]) << endl;
01859
01860 return best_match;
01861
01862 }
01863
01864 EMData* RT3DGridAligner::align(EMData * this_img, EMData *to, const string & cmp_name, const Dict& cmp_params) const
01865 {
01866
01867 vector<Dict> alis = xform_align_nbest(this_img,to,1,cmp_name,cmp_params);
01868
01869 Dict t;
01870 Transform* tr = (Transform*) alis[0]["xform.align3d"];
01871 t["transform"] = tr;
01872 EMData* soln = this_img->process("xform",t);
01873 soln->set_attr("xform.align3d",tr);
01874 delete tr; tr = 0;
01875
01876 return soln;
01877
01878 }
01879
01880 vector<Dict> RT3DGridAligner::xform_align_nbest(EMData * this_img, EMData * to, const unsigned int nsoln, const string & cmp_name, const Dict& cmp_params) const {
01881
01882 if ( this_img->get_ndim() != 3 || to->get_ndim() != 3 ) {
01883 throw ImageDimensionException("This aligner only works for 3D images");
01884 }
01885
01886 int searchx = 0;
01887 int searchy = 0;
01888 int searchz = 0;
01889
01890 bool dotrans = params.set_default("dotrans",1);
01891 if (params.has_key("search")) {
01892 vector<string> check;
01893 check.push_back("searchx");
01894 check.push_back("searchy");
01895 check.push_back("searchz");
01896 for(vector<string>::const_iterator cit = check.begin(); cit != check.end(); ++cit) {
01897 if (params.has_key(*cit)) throw InvalidParameterException("The search parameter is mutually exclusive of the searchx, searchy, and searchz parameters");
01898 }
01899 int search = params["search"];
01900 searchx = search;
01901 searchy = search;
01902 searchz = search;
01903 } else {
01904 searchx = params.set_default("searchx",3);
01905 searchy = params.set_default("searchy",3);
01906 searchz = params.set_default("searchz",3);
01907 }
01908
01909 float lalt = params.set_default("alt0",0.0f);
01910 float laz = params.set_default("az0",0.0f);
01911 float lphi = params.set_default("phi0",0.0f);
01912 float ualt = params.set_default("alt1",180.0f);
01913 float uphi = params.set_default("phi1",360.0f);
01914 float uaz = params.set_default("az1",360.0f);
01915 float dalt = params.set_default("dalt",10.f);
01916 float daz = params.set_default("daz",10.f);
01917 float dphi = params.set_default("dphi",10.f);
01918 bool verbose = params.set_default("verbose",false);
01919
01920
01921 Dict altered_cmp_params(cmp_params);
01922 if (cmp_name == "ccc.tomo") {
01923 altered_cmp_params.set_default("searchx", searchx);
01924 altered_cmp_params.set_default("searchy", searchy);
01925 altered_cmp_params.set_default("searchz", searchz);
01926 altered_cmp_params.set_default("norm", true);
01927 }
01928
01929 vector<Dict> solns;
01930 if (nsoln == 0) return solns;
01931 for (unsigned int i = 0; i < nsoln; ++i ) {
01932 Dict d;
01933 d["score"] = 1.e24;
01934 Transform t;
01935 d["xform.align3d"] = &t;
01936 solns.push_back(d);
01937 }
01938
01939 bool tomography = (cmp_name == "ccc.tomo") ? 1 : 0;
01940 EMData * tofft = 0;
01941 if(dotrans || tomography){
01942 tofft = to->do_fft();
01943 }
01944
01945 #ifdef EMAN2_USING_CUDA
01946 if(EMData::usecuda == 1) {
01947 if(!this_img->isrodataongpu()) this_img->copy_to_cudaro();
01948 if(!to->cudarwdata) to->copy_to_cuda();
01949 if(to->cudarwdata){if(tofft) tofft->copy_to_cuda();}
01950 }
01951 #endif
01952
01953 Dict d;
01954 d["type"] = "eman";
01955 Transform trans = Transform();
01956 bool use_cpu = true;
01957 for ( float alt = lalt; alt <= ualt; alt += dalt) {
01958
01959
01960 for ( float az = laz; az < uaz; az += daz ){
01961 if (verbose) {
01962 cout << "Trying angle alt " << alt << " az " << az << endl;
01963 }
01964 for( float phi = lphi; phi < uphi; phi += dphi ) {
01965 d["alt"] = alt;
01966 d["phi"] = phi;
01967 d["az"] = az;
01968 Transform t(d);
01969 EMData* transformed = this_img->process("xform",Dict("transform",&t));
01970
01971
01972 float best_score = 0.0f;
01973 if(dotrans || tomography){
01974 EMData* ccf = transformed->calc_ccf(tofft);
01975 #ifdef EMAN2_USING_CUDA
01976 if(to->cudarwdata){
01977 use_cpu = false;;
01978 CudaPeakInfo* data = calc_max_location_wrap_cuda(ccf->cudarwdata, ccf->get_xsize(), ccf->get_ysize(), ccf->get_zsize(), searchx, searchy, searchz);
01979 trans.set_trans((float)-data->px, (float)-data->py, (float)-data->pz);
01980 t = trans*t;
01981 if (tomography) {
01982 float2 stats = get_stats_cuda(ccf->cudarwdata, ccf->get_xsize(), ccf->get_ysize(), ccf->get_zsize());
01983 best_score = -(data->peak - stats.x)/sqrt(stats.y);
01984 } else {
01985 best_score = -data->peak;
01986 }
01987 delete data;
01988 }
01989 #endif
01990 if(use_cpu){
01991 if(tomography) ccf->process_inplace("normalize");
01992 IntPoint point = ccf->calc_max_location_wrap(searchx,searchy,searchz);
01993 trans.set_trans((float)-point[0], (float)-point[1], (float)-point[2]);
01994 t = trans*t;
01995 best_score = -ccf->get_value_at_wrap(point[0], point[1], point[2]);
01996 }
01997 delete ccf; ccf =0;
01998 delete transformed; transformed = 0;
01999 }
02000
02001 if(!tomography){
02002 if(!transformed) transformed = this_img->process("xform",Dict("transform",&t));
02003 best_score = transformed->cmp(cmp_name,to,cmp_params);
02004 delete transformed; transformed = 0;
02005 }
02006
02007 unsigned int j = 0;
02008 for ( vector<Dict>::iterator it = solns.begin(); it != solns.end(); ++it, ++j ) {
02009 if ( (float)(*it)["score"] > best_score ) {
02010 vector<Dict>::reverse_iterator rit = solns.rbegin();
02011 copy(rit+1,solns.rend()-j,rit);
02012 Dict& d = (*it);
02013 d["score"] = best_score;
02014 d["xform.align3d"] = &t;
02015 break;
02016 }
02017 }
02018 }
02019 }
02020 }
02021
02022 if(tofft) {delete tofft; tofft = 0;}
02023 return solns;
02024
02025 }
02026
02027 EMData* RT3DSphereAligner::align(EMData * this_img, EMData *to, const string & cmp_name, const Dict& cmp_params) const
02028 {
02029
02030 vector<Dict> alis = xform_align_nbest(this_img,to,1,cmp_name,cmp_params);
02031
02032 Dict t;
02033 Transform* tr = (Transform*) alis[0]["xform.align3d"];
02034 t["transform"] = tr;
02035 EMData* soln = this_img->process("xform",t);
02036 soln->set_attr("xform.align3d",tr);
02037 delete tr; tr = 0;
02038
02039 return soln;
02040
02041 }
02042
02043 vector<Dict> RT3DSphereAligner::xform_align_nbest(EMData * this_img, EMData * to, const unsigned int nsoln, const string & cmp_name, const Dict& cmp_params) const {
02044
02045 if ( this_img->get_ndim() != 3 || to->get_ndim() != 3 ) {
02046 throw ImageDimensionException("This aligner only works for 3D images");
02047 }
02048
02049 int searchx = 0;
02050 int searchy = 0;
02051 int searchz = 0;
02052
02053 bool dotrans = params.set_default("dotrans",1);
02054 if (params.has_key("search")) {
02055 vector<string> check;
02056 check.push_back("searchx");
02057 check.push_back("searchy");
02058 check.push_back("searchz");
02059 for(vector<string>::const_iterator cit = check.begin(); cit != check.end(); ++cit) {
02060 if (params.has_key(*cit)) throw InvalidParameterException("The search parameter is mutually exclusive of the searchx, searchy, and searchz parameters");
02061 }
02062 int search = params["search"];
02063 searchx = search;
02064 searchy = search;
02065 searchz = search;
02066 } else {
02067 searchx = params.set_default("searchx",3);
02068 searchy = params.set_default("searchy",3);
02069 searchz = params.set_default("searchz",3);
02070 }
02071
02072 float lphi = params.set_default("phi0",0.0f);
02073 float uphi = params.set_default("phi1",360.0f);
02074 float dphi = params.set_default("dphi",10.f);
02075 float threshold = params.set_default("threshold",0.f);
02076 if (threshold < 0.0f) throw InvalidParameterException("The threshold parameter must be greater than or equal to zero");
02077 bool verbose = params.set_default("verbose",false);
02078
02079
02080 Dict altered_cmp_params(cmp_params);
02081 if (cmp_name == "ccc.tomo") {
02082 altered_cmp_params.set_default("searchx", searchx);
02083 altered_cmp_params.set_default("searchy", searchy);
02084 altered_cmp_params.set_default("searchz", searchz);
02085 altered_cmp_params.set_default("norm", true);
02086 }
02087
02088 vector<Dict> solns;
02089 if (nsoln == 0) return solns;
02090 for (unsigned int i = 0; i < nsoln; ++i ) {
02091 Dict d;
02092 d["score"] = 1.e24;
02093 Transform t;
02094 d["xform.align3d"] = &t;
02095 solns.push_back(d);
02096 }
02097
02098 Dict d;
02099 d["inc_mirror"] = true;
02100 if ( params.has_key("delta") && params.has_key("n") ) {
02101 throw InvalidParameterException("The delta and n parameters are mutually exclusive in the RT3DSphereAligner aligner");
02102 } else if ( params.has_key("n") ) {
02103 d["n"] = params["n"];
02104 } else {
02105
02106 d["delta"] = params.set_default("delta",10.f);
02107 }
02108
02109 if ((string)params.set_default("orientgen","eman")=="eman") d["perturb"]=0;
02110 Symmetry3D* sym = Factory<Symmetry3D>::get((string)params.set_default("sym","c1"));
02111 vector<Transform> transforms = sym->gen_orientations((string)params.set_default("orientgen","eman"),d);
02112
02113 bool tomography = (cmp_name == "ccc.tomo") ? 1 : 0;
02114
02115
02116 EMData * this_imgfft = 0;
02117 if(dotrans || tomography){
02118 this_imgfft = this_img->do_fft();
02119 }
02120
02121
02122 #ifdef EMAN2_USING_CUDA
02123 if(EMData::usecuda == 1) {
02124 cout << "Using CUDA for 3D alignment" << endl;
02125 if(!to->isrodataongpu()) to->copy_to_cudaro();
02126 if(!this_img->cudarwdata) this_img->copy_to_cuda();
02127 if(this_imgfft) this_imgfft->copy_to_cuda();
02128 }
02129 #endif
02130
02131 Transform trans = Transform();
02132 bool use_cpu = true;
02133 for(vector<Transform>::const_iterator trans_it = transforms.begin(); trans_it != transforms.end(); trans_it++) {
02134 Dict params = trans_it->get_params("eman");
02135 Transform t(params);
02136 if (verbose) {
02137 float alt = params["alt"];
02138 float az = params["az"];
02139 cout << "Trying angle alt: " << alt << " az: " << az << endl;
02140 }
02141
02142 for( float phi = lphi; phi < uphi; phi += dphi ) {
02143
02144 params["phi"] = phi;
02145 t.set_rotation(params);
02146 EMData* transformed = to->process("xform",Dict("transform",&t));
02147
02148
02149 float best_score = 0.0f;
02150 if(dotrans || tomography){
02151 EMData* ccf = transformed->calc_ccf(this_imgfft);
02152 #ifdef EMAN2_USING_CUDA
02153 if(this_img->cudarwdata){
02154 use_cpu = false;;
02155 CudaPeakInfo* data = calc_max_location_wrap_cuda(ccf->cudarwdata, ccf->get_xsize(), ccf->get_ysize(), ccf->get_zsize(), searchx, searchy, searchz);
02156 trans.set_trans((float)-data->px, (float)-data->py, (float)-data->pz);
02157 t = trans*t;
02158 if (tomography) {
02159 float2 stats = get_stats_cuda(ccf->cudarwdata, ccf->get_xsize(), ccf->get_ysize(), ccf->get_zsize());
02160 best_score = -(data->peak - stats.x)/sqrt(stats.y);
02161 } else {
02162 best_score = -data->peak;
02163 }
02164 delete data;
02165 }
02166 #endif
02167 if(use_cpu){
02168 if(tomography) ccf->process_inplace("normalize");
02169 IntPoint point = ccf->calc_max_location_wrap(searchx,searchy,searchz);
02170 trans.set_trans((float)-point[0], (float)-point[1], (float)-point[2]);
02171 t = trans*t;
02172 best_score = -ccf->get_value_at_wrap(point[0], point[1], point[2]);
02173 }
02174 delete ccf; ccf =0;
02175 delete transformed; transformed = 0;
02176 }
02177
02178 if(!tomography){
02179 if(!transformed) transformed = to->process("xform",Dict("transform",&t));
02180 best_score = transformed->cmp(cmp_name,this_img,cmp_params);
02181 delete transformed; transformed = 0;
02182 }
02183
02184 unsigned int j = 0;
02185
02186 for ( vector<Dict>::iterator it = solns.begin(); it != solns.end(); ++it, ++j ) {
02187 if ( (float)(*it)["score"] > best_score ) {
02188 vector<Dict>::reverse_iterator rit = solns.rbegin();
02189 copy(rit+1,solns.rend()-j,rit);
02190 Dict& d = (*it);
02191 d["score"] = best_score;
02192 t.invert();
02193 d["xform.align3d"] = &t;
02194 break;
02195 }
02196 }
02197
02198 }
02199 }
02200 delete sym; sym = 0;
02201 if(this_imgfft) {delete this_imgfft; this_imgfft = 0;}
02202 return solns;
02203
02204 }
02205
02206
02207 EMData* RT3DSymmetryAligner::align(EMData * this_img, EMData *to, const string & cmp_name, const Dict& cmp_params) const
02208 {
02209
02210 vector<Dict> alis = xform_align_nbest(this_img,to,1,cmp_name,cmp_params);
02211
02212 Dict t;
02213 Transform* tr = (Transform*) alis[0]["xform.align3d"];
02214 t["transform"] = tr;
02215 EMData* soln = this_img->process("xform",t);
02216 soln->set_attr("xform.align3d",tr);
02217 delete tr; tr = 0;
02218
02219 return soln;
02220
02221 }
02222
02223 vector<Dict> RT3DSymmetryAligner::xform_align_nbest(EMData * this_img, EMData * to, const unsigned int nsoln, const string & cmp_name, const Dict& cmp_params) const
02224 {
02225
02226 bool verbose = params.set_default("verbose",false);
02227
02228 vector<Dict> solns;
02229 if (nsoln == 0) return solns;
02230 for (unsigned int i = 0; i < nsoln; ++i ) {
02231 Dict d;
02232 d["score"] = 1.e24;
02233 Transform t;
02234 d["xform.align3d"] = &t;
02235 solns.push_back(d);
02236 }
02237
02238
02239 vector<Transform> syms = Symmetry3D::get_symmetries((string)params.set_default("sym","icos"));
02240
02241 float score = 0.0f;
02242 for ( vector<Transform>::const_iterator symit = syms.begin(); symit != syms.end(); ++symit ) {
02243 Transform sympos = *symit;
02244
02245 EMData* transformed = this_img->process("xform",Dict("transform",&sympos));
02246 score = transformed->cmp(cmp_name,this_img,cmp_params);
02247 delete transformed; transformed = 0;
02248
02249 if (verbose) {
02250 Dict rots = sympos.get_rotation("eman");
02251 cout <<"Score is: " << score << " az " << float(rots["az"]) << " alt " << float(rots["alt"]) << " phi " << float(rots["phi"]) << endl;
02252 }
02253
02254 unsigned int j = 0;
02255 for ( vector<Dict>::iterator it = solns.begin(); it != solns.end(); ++it, ++j ) {
02256 if ( (float)(*it)["score"] > score ) {
02257 vector<Dict>::reverse_iterator rit = solns.rbegin();
02258 copy(rit+1,solns.rend()-j,rit);
02259 Dict& d = (*it);
02260 d["score"] = score;
02261 d["xform.align3d"] = &sympos;
02262 break;
02263 }
02264 }
02265 }
02266 return solns;
02267 }
02268
02269 namespace {
02270 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)
02271 {
02272 int size=rsize;
02273 float dx,dy;
02274 int bw=size/2;
02275 int MAXR=this_img->get_ysize()/2;
02276
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];
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();
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){
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++){
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){
02328 if(m!= 0){
02329 int ssize=tsize-2*m;
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
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){
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
02372
02373
02374
02375 Phi2=n_best*M_PI/bw;
02376 Phi=m_best*M_PI/bw;
02377 Vx=p*cos(Phi);
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;
02383 dy=-Ty;
02384
02385 EMData *this_tmp=this_img->copy();
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
02391 delete this_tmp;
02392 if(corre<=corre_fcs) {
02393 corre_fcs=corre;
02394 result[0+5*p] = float(p);
02395 result[1+5*p] = corre;
02396 result[2+5*p] = (Phi2-Phi)*180.0f;
02397 result[3+5*p] = Tx;
02398 result[4+5*p] = Ty;
02399 }
02400 maxcor[p]=corre_fcs;
02401 if(corre_fcs<maxcor_sofar) {
02402 maxcor_sofar=corre_fcs;
02403 rho_best=p;
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;
02409 }
02410 }
02411 }
02412
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;
02428 dx = -Tx;
02429 dy = -Ty;
02430 this_img->rotate(-ang_keep,0,0);
02431 this_img->translate(dx,dy,0.0);
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;
02442 }
02443 }
02444
02445
02446 EMData *FRM2DAligner::align(EMData * this_img, EMData * to,
02447 const string & cmp_name, const Dict& cmp_params) const
02448 {
02449 if (!this_img) {
02450 return 0;
02451 }
02452 if (to && !EMUtil::is_same_size(this_img, to))
02453 throw ImageDimensionException("Images must be the same size to perform translational alignment");
02454
02455 int nx=this_img->get_xsize();
02456 int ny=this_img->get_ysize();
02457 int size =(int)floor(M_PI*ny/4.0);
02458 size =Util::calc_best_fft_size(size);
02459 int MAXR=ny/2;
02460
02461 EMData *this_temp=this_img->copy();
02462 FloatPoint com_test,com_test1;
02463 com_test=this_temp->calc_center_of_mass();
02464 float com_this_x=com_test[0];
02465 float com_this_y=com_test[1];
02466 delete this_temp;
02467
02468
02469 EMData *that_temp=to->copy();
02470 com_test1=that_temp->calc_center_of_mass();
02471 float com_with_x=com_test1[0];
02472 float com_with_y=com_test1[1];
02473 delete that_temp;
02474
02475 EMData *avg_frm=to->copy();
02476 float dx,dy;
02477
02478
02479
02480 EMData *withpcs=avg_frm->unwrap_largerR(0,MAXR,size,float(MAXR));
02481
02482 EMData *withpcsfft=withpcs->oneDfftPolar(size, float(MAXR), float(MAXR));
02483
02484 float *sampl_fft=withpcsfft->get_data();
02485 delete avg_frm;
02486 delete withpcs;
02487
02488 int bw=size/2;
02489 unsigned long ind1=0, ind2=0, ind3=0, ind4=0, ind41=0;
02490 float pi2=2.0*M_PI, r2;
02491
02492 EMData *data_in=new EMData;
02493 data_in->set_complex(true);
02494 data_in->set_ri(1);
02495 data_in->set_size(2*size,1,1);
02496 float * comp_in=data_in->get_data();
02497
02498 int p_max=3;
02499 float *frm2dhhat=0;
02500
02501 if( (frm2dhhat=(float *)malloc((p_max+1)*(size+2)*bw*size*2* sizeof(float)))==NULL){
02502 cout <<"Error in allocating memory 13. \n";
02503 exit(1);
02504 }
02505
02506 float *sb=0, *cb=0;
02507 if((sb=new float[size])==NULL||(cb=new float[size])==NULL) {
02508 cout <<"can't allocate more memory, terminating. \n";
02509 exit(1);
02510 }
02511 for(int i=0;i<size;++i) {
02512 float beta=i*M_PI/bw;
02513 sb[i]=sin(beta);
02514 cb[i]=cos(beta);
02515 }
02516
02517 for(int p=0; p<=p_max; ++p){
02518 ind1=p*size*bw;
02519 float pp2=(float)(p*p);
02520 for(int n=0;n<bw;++n){
02521 ind2=ind1+n;
02522 for(int r=0;r<=MAXR;++r) {
02523 ind3=(ind2+r*bw)*size;
02524 float rr2=(float)(r*r);
02525 float rp2=(float)(r*p);
02526 for(int i=0;i<size;++i){
02527 r2=std::sqrt((float)(rr2+pp2-2.0*rp2*cb[i]));
02528 int r1=(int)floor(r2+0.5f);
02529 if(r1>MAXR){
02530 comp_in[2*i]=0.0f;
02531 comp_in[2*i+1]=0.0f;
02532 }
02533 else{
02534 float gn_r=sampl_fft[2*n+r1*(size+2)];
02535 float gn_i=sampl_fft[2*n+1+r1*(size+2)];
02536 float sb2, cb2, cn, sn;
02537 if(n!=0){
02538 if(r2 != 0.0){
02539 sb2=r*sb[i]/r2;
02540 cb2=(r*cb[i]-p)/r2;
02541 }
02542 else{
02543 sb2=0.0;
02544 cb2=1.0;
02545 }
02546 if(sb2>1.0) sb2= 1.0f;
02547 if(sb2<-1.0)sb2=-1.0f;
02548 if(cb2>1.0) cb2= 1.0f;
02549 if(cb2<-1.0)cb2=-1.0f;
02550 float beta2=atan2(sb2,cb2);
02551 if(beta2<0.0) beta2+=pi2;
02552 float nb2=n*beta2;
02553 cn=cos(nb2);
02554 sn=sin(nb2);
02555 }
02556 else{
02557 cn=1.0f; sn=0.0f;
02558 }
02559 comp_in[2*i]=cn*gn_r-sn*gn_i;
02560 comp_in[2*i+1]=-(cn*gn_i+sn*gn_r);
02561 }
02562 }
02563 EMData *data_out;
02564 data_out=data_in->do_fft();
02565 float * comp_out=data_out->get_data();
02566 for(int m=0;m<size;m++){
02567 ind4=(ind3+m)*2;
02568 ind41=ind4+1;
02569 frm2dhhat[ind4]=comp_out[2*m];
02570 frm2dhhat[ind41]=comp_out[2*m+1];
02571 }
02572 delete data_out;
02573 }
02574 }
02575 }
02576
02577 delete[] sb;
02578 delete[] cb;
02579 delete data_in;
02580 delete withpcsfft;
02581
02582 float dot_frm0=0.0f, dot_frm1=0.0f;
02583 EMData *da_nFlip=0, *da_yFlip=0, *dr_frm=0;
02584
02585 for (int iFlip=0;iFlip<=1;++iFlip){
02586 if (iFlip==0){dr_frm=this_img->copy(); da_nFlip=this_img->copy();}
02587 else {dr_frm=this_img->copy(); da_yFlip=this_img->copy();}
02588 if(iFlip==1) {com_this_x=nx-com_this_x; }
02589
02590 dx=-(com_this_x-nx/2);
02591 dy=-(com_this_y-ny/2);
02592 dr_frm->translate(dx,dy,0.0);
02593 EMData *selfpcs = dr_frm->unwrap_largerR(0,MAXR,size, (float)MAXR);
02594
02595 EMData *selfpcsfft = selfpcs->oneDfftPolar(size, (float)MAXR, (float)MAXR);
02596 delete selfpcs;
02597 delete dr_frm;
02598 if(iFlip==0)
02599 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);
02600 else
02601 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);
02602 delete selfpcsfft;
02603 }
02604
02605 delete[] frm2dhhat;
02606 if(dot_frm0 <=dot_frm1) {
02607 #ifdef DEBUG
02608 printf("best_corre=%f, no flip\n",dot_frm0);
02609 #endif
02610 delete da_yFlip;
02611 return da_nFlip;
02612 }
02613 else {
02614 #ifdef DEBUG
02615 printf("best_corre=%f, flipped\n",dot_frm1);
02616 #endif
02617 delete da_nFlip;
02618 return da_yFlip;
02619 }
02620 }
02621
02622 #ifdef SPARX_USING_CUDA
02623
02624 CUDA_Aligner::CUDA_Aligner(int id) {
02625 image_stack = NULL;
02626 image_stack_filtered = NULL;
02627 ccf = NULL;
02628 if (id != -1) cudasetup(id);
02629 }
02630
02631 void CUDA_Aligner::finish() {
02632 if (image_stack) free(image_stack);
02633 if (image_stack_filtered) free(image_stack_filtered);
02634 if (ccf) free(ccf);
02635 image_stack = NULL;
02636 image_stack_filtered = NULL;
02637 ccf = NULL;
02638 }
02639
02640 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) {
02641
02642 NIMA = nima;
02643 NX = nx;
02644 NY = ny;
02645 RING_LENGTH = ring_length;
02646 NRING = nring;
02647 STEP = step;
02648 KX = kx;
02649 KY = ky;
02650 OU = ou;
02651 CTF = ctf;
02652
02653 image_stack = (float *)malloc(NIMA*NX*NY*sizeof(float));
02654 if (CTF == 1) image_stack_filtered = (float *)malloc(NIMA*NX*NY*sizeof(float));
02655 ccf = (float *)malloc(2*(2*KX+1)*(2*KY+1)*NIMA*(RING_LENGTH+2)*sizeof(float));
02656 }
02657
02658 void CUDA_Aligner::insert_image(EMData *image, int num) {
02659
02660 int base_address = num*NX*NY;
02661
02662 for (int y=0; y<NY; y++)
02663 for (int x=0; x<NX; x++)
02664 image_stack[base_address+y*NX+x] = (*image)(x, y);
02665 }
02666
02667 void CUDA_Aligner::filter_stack(vector<float> ctf_params) {
02668
02669 float *params;
02670
02671 params = (float *)malloc(NIMA*6*sizeof(float));
02672
02673 for (int i=0; i<NIMA*6; i++) params[i] = ctf_params[i];
02674
02675 filter_image(image_stack, image_stack_filtered, NIMA, NX, NY, params);
02676
02677 free(params);
02678 }
02679
02680 void CUDA_Aligner::sum_oe(vector<float> ctf_params, vector<float> ali_params, EMData *ave1, EMData *ave2) {
02681
02682 float *ctf_p, *ali_p, *av1, *av2;
02683
02684 ctf_p = (float *)malloc(NIMA*6*sizeof(float));
02685 ali_p = (float *)malloc(NIMA*4*sizeof(float));
02686
02687 if (CTF == 1) {
02688 for (int i=0; i<NIMA*6; i++) ctf_p[i] = ctf_params[i];
02689 }
02690 for (int i=0; i<NIMA*4; i++) ali_p[i] = ali_params[i];
02691
02692 av1 = ave1->get_data();
02693 av2 = ave2->get_data();
02694
02695 rot_filt_sum(image_stack, NIMA, NX, NY, CTF, ctf_p, ali_p, av1, av2);
02696
02697 free(ctf_p);
02698 free(ali_p);
02699 }
02700
02701 vector<float> CUDA_Aligner::alignment_2d(EMData *ref_image_em, vector<float> sx_list, vector<float> sy_list, int silent) {
02702
02703 float *ref_image, max_ccf;
02704 int base_address, ccf_offset;
02705 float ts, tm;
02706 float ang, sx = 0, sy = 0, mirror, co, so, sxs, sys;
02707 float *sx2, *sy2;
02708 vector<float> align_result;
02709
02710 sx2 = (float *)malloc(NIMA*sizeof(float));
02711 sy2 = (float *)malloc(NIMA*sizeof(float));
02712
02713 ref_image = ref_image_em->get_data();
02714
02715 for (int i=0; i<NIMA; i++) {
02716 sx2[i] = sx_list[i];
02717 sy2[i] = sy_list[i];
02718 }
02719
02720 if (CTF == 1) {
02721 calculate_ccf(image_stack_filtered, ref_image, ccf, NIMA, NX, NY, RING_LENGTH, NRING, OU, STEP, KX, KY, sx2, sy2, silent);
02722 } else {
02723 calculate_ccf(image_stack, ref_image, ccf, NIMA, NX, NY, RING_LENGTH, NRING, OU, STEP, KX, KY, sx2, sy2, silent);
02724 }
02725
02726 ccf_offset = NIMA*(RING_LENGTH+2)*(2*KX+1)*(2*KY+1);
02727
02728 for (int im=0; im<NIMA; im++) {
02729 max_ccf = -1.0e22;
02730 for (int kx=-KX; kx<=KX; kx++) {
02731 for (int ky=-KY; ky<=KY; ky++) {
02732 base_address = (((ky+KY)*(2*KX+1)+(kx+KX))*NIMA+im)*(RING_LENGTH+2);
02733 for (int l=0; l<RING_LENGTH; l++) {
02734 ts = ccf[base_address+l];
02735 tm = ccf[base_address+l+ccf_offset];
02736 if (ts > max_ccf) {
02737 ang = float(l)/RING_LENGTH*360.0;
02738 sx = -kx*STEP;
02739 sy = -ky*STEP;
02740 mirror = 0;
02741 max_ccf = ts;
02742 }
02743 if (tm > max_ccf) {
02744 ang = float(l)/RING_LENGTH*360.0;
02745 sx = -kx*STEP;
02746 sy = -ky*STEP;
02747 mirror = 1;
02748 max_ccf = tm;
02749 }
02750 }
02751 }
02752 }
02753 co = cos(ang*M_PI/180.0);
02754 so = -sin(ang*M_PI/180.0);
02755 sxs = sx*co - sy*so;
02756 sys = sx*so + sy*co;
02757
02758 align_result.push_back(ang);
02759 align_result.push_back(sxs);
02760 align_result.push_back(sys);
02761 align_result.push_back(mirror);
02762 }
02763
02764 free(sx2);
02765 free(sy2);
02766
02767 return align_result;
02768 }
02769
02770
02771 vector<float> CUDA_Aligner::ali2d_single_iter(EMData *ref_image_em, vector<float> ali_params, float csx, float csy, int silent, float delta) {
02772
02773 float *ref_image, max_ccf;
02774 int base_address, ccf_offset;
02775 float ts, tm;
02776 float ang = 0.0, sx = 0.0, sy = 0.0, co, so, sxs, sys;
02777 int mirror;
02778 float *sx2, *sy2;
02779 vector<float> align_result;
02780
02781 sx2 = (float *)malloc(NIMA*sizeof(float));
02782 sy2 = (float *)malloc(NIMA*sizeof(float));
02783
02784 ref_image = ref_image_em->get_data();
02785
02786 for (int i=0; i<NIMA; i++) {
02787 ang = ali_params[i*4]/180.0*M_PI;
02788 sx = (ali_params[i*4+3] < 0.5)?(ali_params[i*4+1]-csx):(ali_params[i*4+1]+csx);
02789 sy = ali_params[i*4+2]-csy;
02790 co = cos(ang);
02791 so = sin(ang);
02792 sx2[i] = -(sx*co-sy*so);
02793 sy2[i] = -(sx*so+sy*co);
02794 }
02795
02796 if (CTF == 1) {
02797 calculate_ccf(image_stack_filtered, ref_image, ccf, NIMA, NX, NY, RING_LENGTH, NRING, OU, STEP, KX, KY, sx2, sy2, silent);
02798 } else {
02799 calculate_ccf(image_stack, ref_image, ccf, NIMA, NX, NY, RING_LENGTH, NRING, OU, STEP, KX, KY, sx2, sy2, silent);
02800 }
02801
02802 ccf_offset = NIMA*(RING_LENGTH+2)*(2*KX+1)*(2*KY+1);
02803
02804 float sx_sum = 0.0f;
02805 float sy_sum = 0.0f;
02806
02807 int dl;
02808 dl = static_cast<int>(delta/360.0*RING_LENGTH);
02809 if (dl<1) { dl = 1; }
02810
02811 for (int im=0; im<NIMA; im++) {
02812 max_ccf = -1.0e22;
02813 for (int kx=-KX; kx<=KX; kx++) {
02814 for (int ky=-KY; ky<=KY; ky++) {
02815 base_address = (((ky+KY)*(2*KX+1)+(kx+KX))*NIMA+im)*(RING_LENGTH+2);
02816 for (int l=0; l<RING_LENGTH; l+=dl) {
02817 ts = ccf[base_address+l];
02818 tm = ccf[base_address+l+ccf_offset];
02819 if (ts > max_ccf) {
02820 ang = float(l)/RING_LENGTH*360.0;
02821 sx = -kx*STEP;
02822 sy = -ky*STEP;
02823 mirror = 0;
02824 max_ccf = ts;
02825 }
02826 if (tm > max_ccf) {
02827 ang = float(l)/RING_LENGTH*360.0;
02828 sx = -kx*STEP;
02829 sy = -ky*STEP;
02830 mirror = 1;
02831 max_ccf = tm;
02832 }
02833 }
02834 }
02835 }
02836 co = cos(ang*M_PI/180.0);
02837 so = -sin(ang*M_PI/180.0);
02838
02839 sxs = (sx-sx2[im])*co-(sy-sy2[im])*so;
02840 sys = (sx-sx2[im])*so+(sy-sy2[im])*co;
02841
02842
02843
02844 align_result.push_back(ang);
02845 align_result.push_back(sxs);
02846 align_result.push_back(sys);
02847 align_result.push_back(mirror);
02848
02849 if (mirror == 0) { sx_sum += sxs; } else { sx_sum -= sxs; }
02850 sy_sum += sys;
02851 }
02852
02853 align_result.push_back(sx_sum);
02854 align_result.push_back(sy_sum);
02855
02856 free(sx2);
02857 free(sy2);
02858
02859 return align_result;
02860 }
02861
02862
02863 CUDA_multiref_aligner::CUDA_multiref_aligner(int id) {
02864 image_stack = NULL;
02865 ref_image_stack = NULL;
02866 ref_image_stack_filtered = NULL;
02867 ccf = NULL;
02868 ctf_params = NULL;
02869 ali_params = NULL;
02870 cudasetup(id);
02871 }
02872
02873
02874 void CUDA_multiref_aligner::finish() {
02875 if (image_stack) free(image_stack);
02876 if (ref_image_stack) free(ref_image_stack);
02877 if (ref_image_stack_filtered) free(ref_image_stack_filtered);
02878 if (ccf) free(ccf);
02879 if (ctf_params) free(ctf_params);
02880 if (ali_params) free(ali_params);
02881 image_stack = NULL;
02882 ref_image_stack = NULL;
02883 ref_image_stack_filtered = NULL;
02884 ccf = NULL;
02885 ctf_params = NULL;
02886 ali_params = NULL;
02887 }
02888
02889 void CUDA_multiref_aligner::setup(int nima, int nref, int nx, int ny, int ring_length, int nring, int ou, float step, int kx, int ky, bool ctf) {
02890
02891 NIMA = nima;
02892 NREF = nref;
02893 NX = nx;
02894 NY = ny;
02895 RING_LENGTH = ring_length;
02896 NRING = nring;
02897 STEP = step;
02898 KX = kx;
02899 KY = ky;
02900 OU = ou;
02901 CTF = ctf;
02902
02903
02904 MAX_IMAGE_BATCH = 10;
02905
02906 image_stack = (float *)malloc(NIMA*NX*NY*sizeof(float));
02907 ref_image_stack = (float *)malloc(NREF*NX*NY*sizeof(float));
02908 if (CTF == 1) ref_image_stack_filtered = (float *)malloc(NREF*NX*NY*sizeof(float));
02909 ccf = (float *)malloc(2*(2*KX+1)*(2*KY+1)*NREF*(RING_LENGTH+2)*MAX_IMAGE_BATCH*sizeof(float));
02910 }
02911
02912 void CUDA_multiref_aligner::setup_params(vector<float> all_ali_params, vector<float> all_ctf_params) {
02913
02914 ali_params = (float *)malloc(NIMA*4*sizeof(float));
02915 for (int i=0; i<NIMA*4; i++) ali_params[i] = all_ali_params[i];
02916 if (CTF == 1) {
02917 ctf_params = (float *)malloc(NIMA*6*sizeof(float));
02918 for (int i=0; i<NIMA*6; i++) ctf_params[i] = all_ctf_params[i];
02919 }
02920 }
02921
02922 void CUDA_multiref_aligner::insert_image(EMData *image, int num) {
02923
02924 int base_address = num*NX*NY;
02925
02926 for (int y=0; y<NY; y++)
02927 for (int x=0; x<NX; x++)
02928 image_stack[base_address+y*NX+x] = (*image)(x, y);
02929 }
02930
02931 void CUDA_multiref_aligner::insert_ref_image(EMData *image, int num) {
02932
02933 int base_address = num*NX*NY;
02934
02935 for (int y=0; y<NY; y++)
02936 for (int x=0; x<NX; x++)
02937 ref_image_stack[base_address+y*NX+x] = (*image)(x, y);
02938 }
02939
02940 vector<float> CUDA_multiref_aligner::multiref_ali2d(int silent) {
02941
02942 float *ctf_params_ref = (float *)malloc(NREF*6*sizeof(float));
02943 float *sx2 = (float *)malloc(NIMA*sizeof(float));
02944 float *sy2 = (float *)malloc(NIMA*sizeof(float));
02945 vector<float> align_results;
02946 int ccf_offset = NREF*(RING_LENGTH+2)*(2*KX+1)*(2*KY+1);
02947
02948 vector<int> batch_size;
02949 vector<int> batch_begin;
02950
02951 if (CTF == 1) {
02952 float previous_defocus = ctf_params[0];
02953 int current_size = 1;
02954 for (int i=1; i<NIMA; i++) {
02955 if (ctf_params[i*6] != previous_defocus || current_size >= MAX_IMAGE_BATCH) {
02956 batch_size.push_back(current_size);
02957 current_size = 1;
02958 previous_defocus = ctf_params[i*6];
02959 } else current_size++;
02960 }
02961 batch_size.push_back(current_size);
02962 } else {
02963 batch_size.resize(NIMA/MAX_IMAGE_BATCH, MAX_IMAGE_BATCH);
02964 if (NIMA%MAX_IMAGE_BATCH != 0) batch_size.push_back(NIMA%MAX_IMAGE_BATCH);
02965 }
02966 int num_batch = batch_size.size();
02967 batch_begin.resize(num_batch, 0);
02968 for (int i=1; i<num_batch; i++) batch_begin[i] = batch_size[i-1]+batch_begin[i-1];
02969 assert(batch_begin[num_batch-1]+batch_size[num_batch-1] == NIMA-1);
02970
02971 for (int i=0; i<NIMA; i++) {
02972 float ang = ali_params[i*4]/180.0*M_PI;
02973 float sx = ali_params[i*4+1];
02974 float sy = ali_params[i*4+2];
02975 float co = cos(ang);
02976 float so = sin(ang);
02977 sx2[i] = -(sx*co-sy*so);
02978 sy2[i] = -(sx*so+sy*co);
02979 }
02980
02981 for (int i=0; i<num_batch; i++) {
02982 if (CTF == 1) {
02983 for (int p=0; p<NREF; p++)
02984 for (int q=0; q<6; q++)
02985 ctf_params_ref[p*6+q] = ctf_params[batch_begin[i]*6+q];
02986 filter_image(ref_image_stack, ref_image_stack_filtered, NREF, NX, NY, ctf_params_ref);
02987 calculate_multiref_ccf(image_stack+batch_begin[i]*NX*NY, ref_image_stack_filtered, ccf, batch_size[i], NREF, NX, NY, RING_LENGTH, NRING, OU, STEP, KX, KY,
02988 sx2+batch_begin[i], sy2+batch_begin[i], silent);
02989 } else {
02990 calculate_multiref_ccf(image_stack+batch_begin[i]*NX*NY, ref_image_stack, ccf, batch_size[i], NREF, NX, NY, RING_LENGTH, NRING, OU, STEP, KX, KY,
02991 sx2+batch_begin[i], sy2+batch_begin[i], silent);
02992 }
02993
02994 for (int j=0; j<batch_size[i]; j++) {
02995 for (int im=0; im<NREF; im++) {
02996 float max_ccf = -1.0e22;
02997 float ang = 0.0, sx = 0.0, sy = 0.0;
02998 int mirror = 0;
02999 for (int kx=-KX; kx<=KX; kx++) {
03000 for (int ky=-KY; ky<=KY; ky++) {
03001 int base_address = (((ky+KY)*(2*KX+1)+(kx+KX))*NREF+im)*(RING_LENGTH+2)+ccf_offset*2*j;
03002 for (int l=0; l<RING_LENGTH; l++) {
03003 float ts = ccf[base_address+l];
03004 float tm = ccf[base_address+l+ccf_offset];
03005 if (ts > max_ccf) {
03006 ang = 360.0-float(l)/RING_LENGTH*360.0;
03007 sx = -kx*STEP;
03008 sy = -ky*STEP;
03009 mirror = 0;
03010 max_ccf = ts;
03011 }
03012 if (tm > max_ccf) {
03013 ang = float(l)/RING_LENGTH*360.0;
03014 sx = -kx*STEP;
03015 sy = -ky*STEP;
03016 mirror = 1;
03017 max_ccf = tm;
03018 }
03019 }
03020 }
03021 }
03022 float co = cos(ang*M_PI/180.0);
03023 float so = -sin(ang*M_PI/180.0);
03024
03025 int img_num = batch_begin[i]+j;
03026 float sxs = (sx-sx2[img_num])*co-(sy-sy2[img_num])*so;
03027 float sys = (sx-sx2[img_num])*so+(sy-sy2[img_num])*co;
03028
03029 align_results.push_back(max_ccf);
03030 align_results.push_back(ang);
03031 align_results.push_back(sxs);
03032 align_results.push_back(sys);
03033 align_results.push_back(mirror);
03034 }
03035 }
03036 }
03037
03038 free(ctf_params_ref);
03039 free(sx2);
03040 free(sy2);
03041
03042 return align_results;
03043 }
03044
03045 #endif
03046
03047
03048 void EMAN::dump_aligners()
03049 {
03050 dump_factory < Aligner > ();
03051 }
03052
03053 map<string, vector<string> > EMAN::dump_aligners_list()
03054 {
03055 return dump_factory_list < Aligner > ();
03056 }