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
00036 #include "reconstructor.h"
00037 #include "plugins/reconstructor_template.h"
00038 #include "ctf.h"
00039 #include "emassert.h"
00040 #include "symmetry.h"
00041 #include <cstring>
00042 #include <fstream>
00043 #include <iomanip>
00044 #include <boost/bind.hpp>
00045 #include <boost/format.hpp>
00046
00047 #include <gsl/gsl_statistics_double.h>
00048 #include <gsl/gsl_fit.h>
00049
00050 using namespace EMAN;
00051 using std::complex;
00052
00053
00054 #include <iostream>
00055 using std::cerr;
00056 using std::endl;
00057 using std::cout;
00058
00059 #include <algorithm>
00060
00061
00062 #include <iomanip>
00063 using std::setprecision;
00064
00065 template < typename T > void checked_delete( T*& x )
00066 {
00067 typedef char type_must_be_complete[ sizeof(T)? 1: -1 ];
00068 (void) sizeof(type_must_be_complete);
00069 delete x;
00070 x = NULL;
00071 }
00072
00073 const string FourierReconstructor::NAME = "fourier";
00074 const string FourierReconstructorSimple2D::NAME = "fouriersimple2D";
00075 const string WienerFourierReconstructor::NAME = "wiener_fourier";
00076 const string BackProjectionReconstructor::NAME = "back_projection";
00077 const string nn4Reconstructor::NAME = "nn4";
00078 const string nnSSNR_Reconstructor::NAME = "nnSSNR";
00079 const string nn4_ctfReconstructor::NAME = "nn4_ctf";
00080 const string nnSSNR_ctfReconstructor::NAME = "nnSSNR_ctf";
00081
00082 template <> Factory < Reconstructor >::Factory()
00083 {
00084 force_add<FourierReconstructor>();
00085 force_add<FourierReconstructorSimple2D>();
00086
00087 force_add<WienerFourierReconstructor>();
00088 force_add<BackProjectionReconstructor>();
00089 force_add<nn4Reconstructor>();
00090 force_add<nnSSNR_Reconstructor>();
00091 force_add<nn4_ctfReconstructor>();
00092 force_add<nnSSNR_ctfReconstructor>();
00093
00094 }
00095
00096 void FourierReconstructorSimple2D::setup()
00097 {
00098 nx = params.set_default("nx",0);
00099
00100 if ( nx < 0 ) throw InvalidValueException(nx, "nx must be positive");
00101
00102 bool is_fftodd = (nx % 2 == 1);
00103
00104 ny = nx;
00105 nx += 2-is_fftodd;
00106
00107 image = new EMData();
00108 image->set_size(nx, ny);
00109 image->set_complex(true);
00110 image->set_fftodd(is_fftodd);
00111 image->set_ri(true);
00112
00113 tmp_data = new EMData();
00114 tmp_data->set_size(nx/2, nx);
00115 }
00116
00117 int FourierReconstructorSimple2D::insert_slice(const EMData* const slice, const Transform & euler, const float weight)
00118 {
00119
00120
00121 if (!slice) throw NullPointerException("EMData pointer (input image) is NULL");
00122
00123 if ( slice->get_ndim() != 1 ) throw ImageDimensionException("Image dimension must be 1");
00124
00125
00126 if (slice->is_complex()) throw ImageFormatException("The image is complex, expecting real");
00127
00128 EMData* working_slice = slice->process("xform.phaseorigin.tocorner");
00129
00130
00131 working_slice->do_fft_inplace();
00132
00133 float* rdata = image->get_data();
00134 float* norm = tmp_data->get_data();
00135 float* dat = working_slice->get_data();
00136
00137 float g[4];
00138 int offset[4];
00139 float dt[2];
00140 offset[0] = 0; offset[1] = 2; offset[2] = nx; offset[3] = nx+2;
00141
00142 float alt = -((float)(euler.get_rotation("2d"))["alpha"])*M_PI/180.0f;
00143 for (int x = 0; x < working_slice->get_xsize() / 2; x++) {
00144
00145 float rx = (float) x;
00146
00147 float xx = rx*cos(alt);
00148 float yy = rx*sin(alt);
00149 float cc = 1.0;
00150
00151 if (xx < 0) {
00152 xx = -xx;
00153 yy = -yy;
00154 cc = -1.0;
00155 }
00156
00157 yy += ny / 2;
00158
00159
00160 dt[0] = dat[2*x];
00161 dt[1] = cc * dat[2*x+1];
00162
00163
00164 int x0 = (int) floor(xx);
00165 int y0 = (int) floor(yy);
00166
00167 int i = 2*x0 + y0*nx;
00168
00169 float dx = xx - x0;
00170 float dy = yy - y0;
00171
00172 g[0] = Util::agauss(1, dx, dy, 0, EMConsts::I2G);
00173 g[1] = Util::agauss(1, 1 - dx, dy, 0, EMConsts::I2G);
00174 g[2] = Util::agauss(1, dx, 1 - dy, 0, EMConsts::I2G);
00175 g[3] = Util::agauss(1, 1 - dx, 1 - dy, 0, EMConsts::I2G);
00176
00177
00178 if ( x0 == nx-2 ) {
00179 int k = i + offset[0];
00180 rdata[k] += g[0] * dt[0];
00181 rdata[k + 1] += g[0] * dt[1];
00182 norm[k/2] += g[0];
00183
00184 k = i + offset[2];
00185 rdata[k] += g[2] * dt[0];
00186 rdata[k + 1] += g[2] * dt[1];
00187 norm[k/2] += g[2];
00188 continue;
00189
00190 }
00191
00192 if ( x0 > nx-2 ) {
00193 int dif = x0 - (nx-2);
00194 x0 -= dif;
00195 }
00196
00197 if ( y0 == ny -1 ) {
00198 int k = i + offset[0];
00199 rdata[k] += g[0] * dt[0];
00200 rdata[k + 1] += g[0] * dt[1];
00201 norm[k/2] += g[0];
00202
00203 k = i + offset[1];
00204 rdata[k] += g[1] * dt[0];
00205 rdata[k + 1] += g[1] * dt[1];
00206 norm[k/2] += g[1];
00207 continue;
00208 }
00209
00210 if ( y0 > ny-1) {
00211 int dif = y0 - (ny-1);
00212 y0 -= dif;
00213 }
00214
00215 if (x0 >= nx - 2 || y0 >= ny - 1) continue;
00216
00217
00218
00219
00220 for (int j = 0; j < 4; j++)
00221 {
00222 int k = i + offset[j];
00223 rdata[k] += g[j] * dt[0];
00224 rdata[k + 1] += g[j] * dt[1];
00225 norm[k/2] += g[j];
00226
00227 }
00228 }
00229
00230 return 0;
00231
00232 }
00233
00234 EMData *FourierReconstructorSimple2D::finish(bool doift)
00235 {
00236 normalize_threed();
00237
00238 image->process_inplace("xform.fourierorigin.tocorner");
00239 image->do_ift_inplace();
00240 image->depad();
00241 image->process_inplace("xform.phaseorigin.tocenter");
00242
00243 return image;
00244 }
00245
00246 void ReconstructorVolumeData::normalize_threed(const bool sqrt_damp,const bool wiener)
00247
00248 {
00249 float* norm = tmp_data->get_data();
00250 float* rdata = image->get_data();
00251
00252
00253
00254
00255 if ( 0 == norm ) throw NullPointerException("The normalization volume was null!");
00256 if ( 0 == rdata ) throw NullPointerException("The complex reconstruction volume was null!");
00257
00258
00259
00260
00261
00262 float wfilt=0.0;
00263 if (wiener) wfilt=1.0;
00264
00265 for (int i = 0; i < subnx * subny * subnz; i += 2) {
00266 float d = norm[i/2];
00267 if (sqrt_damp) d*=sqrt(d);
00268 if (d == 0) {
00269 rdata[i] = 0;
00270 rdata[i + 1] = 0;
00271 }
00272 else {
00273
00274
00275 rdata[i] /= d+wfilt;
00276 rdata[i + 1] /= d+wfilt;
00277 }
00278 }
00279
00280
00281 if (subx0==0 && subnx>1 && subny==ny && subnz==nz) {
00282 for (int z=0; z<=nz/2; z++) {
00283 for (int y=1; y<=ny; y++) {
00284 if (y==0 && z==0) continue;
00285
00286 int i=(y%ny)*subnx+(z%nz)*subnx*subny;
00287 int i2=(ny-y)*subnx+((nz-z)%nz)*subnx*subny;
00288 float ar=(rdata[i]+rdata[i2])/2.0f;
00289 float ai=(rdata[i+1]-rdata[i2+1])/2.0f;
00290 rdata[i]=ar;
00291 rdata[i2]=ar;
00292 rdata[i+1]=ai;
00293 rdata[i2+1]=-ai;
00294 }
00295 }
00296 }
00297
00298 if (subx0+subnx==nx && subnx>1 && subny==ny && subnz==nz) {
00299 for (int z=0; z<=nz/2; z++) {
00300 for (int y=1; y<=ny; y++) {
00301 if (y==0 && z==0) continue;
00302
00303 int i=(y%ny)*subnx+(z%nz)*subnx*subny+subnx-2;
00304 int i2=(ny-y)*subnx+((nz-z)%nz)*subnx*subny+subnx-2;
00305 float ar=(rdata[i]+rdata[i2])/2.0f;
00306 float ai=(rdata[i+1]-rdata[i2+1])/2.0f;
00307 rdata[i]=ar;
00308 rdata[i2]=ar;
00309 rdata[i+1]=ai;
00310 rdata[i2+1]=-ai;
00311 }
00312 }
00313 }
00314 }
00315
00316 void FourierReconstructor::load_default_settings()
00317 {
00318 inserter=0;
00319 image=0;
00320 tmp_data=0;
00321 }
00322
00323 void FourierReconstructor::free_memory()
00324 {
00325 if (image) { delete image; image=0; }
00326 if (tmp_data) { delete tmp_data; tmp_data=0; }
00327 if ( inserter != 0 )
00328 {
00329 delete inserter;
00330 inserter = 0;
00331 }
00332 }
00333
00334 #include <sstream>
00335
00336 void FourierReconstructor::load_inserter()
00337 {
00338
00339
00340 Dict parms;
00341 parms["data"] = image;
00342 parms["norm"] = tmp_data->get_data();
00343
00344
00345
00346
00347
00348
00349
00350
00351 if ( inserter != 0 )
00352 {
00353 delete inserter;
00354 }
00355
00356 inserter = Factory<FourierPixelInserter3D>::get((string)params["mode"], parms);
00357 inserter->init();
00358 }
00359
00360 void FourierReconstructor::setup()
00361 {
00362
00363 params.set_default("mode","gauss_2");
00364
00365 vector<int> size=params["size"];
00366
00367 nx = size[0];
00368 ny = size[1];
00369 nz = size[2];
00370 nx2=nx/2-1;
00371 ny2=ny/2;
00372 nz2=nz/2;
00373
00374
00375
00376 bool is_fftodd = (nx % 2 == 1);
00377
00378
00379
00380 nx += 2-is_fftodd;
00381
00382 if (params.has_key("subvolume")) {
00383 vector<int> sub=params["subvolume"];
00384 subx0=sub[0];
00385 suby0=sub[1];
00386 subz0=sub[2];
00387 subnx=sub[3];
00388 subny=sub[4];
00389 subnz=sub[5];
00390
00391 if (subx0<0 || suby0<0 || subz0<0 || subx0+subnx>nx || suby0+subny>ny || subz0+subnz>nz)
00392 throw ImageDimensionException("The subvolume cannot extend outside the reconstructed volume");
00393
00394 }
00395 else {
00396 subx0=suby0=subz0=0;
00397 subnx=nx;
00398 subny=ny;
00399 subnz=nz;
00400 }
00401
00402
00403
00404 if (image) delete image;
00405 image = new EMData();
00406 image->set_size(subnx, subny, subnz);
00407 image->set_complex(true);
00408 image->set_fftodd(is_fftodd);
00409 image->set_ri(true);
00410 image->to_zero();
00411
00412 if (params.has_key("subvolume")) {
00413 image->set_attr("subvolume_x0",subx0);
00414 image->set_attr("subvolume_y0",suby0);
00415 image->set_attr("subvolume_z0",subz0);
00416 image->set_attr("subvolume_full_nx",nx);
00417 image->set_attr("subvolume_full_ny",ny);
00418 image->set_attr("subvolume_full_nz",nz);
00419 }
00420
00421 if (tmp_data) delete tmp_data;
00422 tmp_data = new EMData();
00423 tmp_data->set_size(subnx/2, subny, subnz);
00424 tmp_data->to_zero();
00425 tmp_data->update();
00426
00427 load_inserter();
00428
00429 if ( (bool) params["quiet"] == false )
00430 {
00431 cout << "3D Fourier dimensions are " << nx << " " << ny << " " << nz << endl;
00432 cout << "3D Fourier subvolume is " << subnx << " " << subny << " " << subnz << endl;
00433 cout << "You will require approximately " << setprecision(3) << (subnx*subny*subnz*sizeof(float)*1.5)/1000000000.0 << "GB of memory to reconstruct this volume" << endl;
00434 }
00435 }
00436
00437 void FourierReconstructor::setup_seed(EMData* seed,float seed_weight) {
00438
00439 params.set_default("mode","gauss_2");
00440
00441 vector<int> size=params["size"];
00442
00443 nx = size[0];
00444 ny = size[1];
00445 nz = size[2];
00446 nx2=nx/2-1;
00447 ny2=ny/2;
00448 nz2=nz/2;
00449
00450
00451
00452 bool is_fftodd = (nx % 2 == 1);
00453
00454
00455
00456 nx += 2-is_fftodd;
00457
00458 if (params.has_key("subvolume")) {
00459 vector<int> sub=params["subvolume"];
00460 subx0=sub[0];
00461 suby0=sub[1];
00462 subz0=sub[2];
00463 subnx=sub[3];
00464 subny=sub[4];
00465 subnz=sub[5];
00466
00467 if (subx0<0 || suby0<0 || subz0<0 || subx0+subnx>nx || suby0+subny>ny || subz0+subnz>nz)
00468 throw ImageDimensionException("The subvolume cannot extend outside the reconstructed volume");
00469
00470 }
00471 else {
00472 subx0=suby0=subz0=0;
00473 subnx=nx;
00474 subny=ny;
00475 subnz=nz;
00476 }
00477
00478 if (seed->get_xsize()!=subnx || seed->get_ysize()!=subny || seed->get_zsize()!=subnz || !seed->is_complex())
00479 throw ImageDimensionException("The dimensions of the seed volume do not match the reconstruction size");
00480
00481
00482 image = seed;
00483 if (params.has_key("subvolume")) {
00484 image->set_attr("subvolume_x0",subx0);
00485 image->set_attr("subvolume_y0",suby0);
00486 image->set_attr("subvolume_z0",subz0);
00487 image->set_attr("subvolume_full_nx",nx);
00488 image->set_attr("subvolume_full_ny",ny);
00489 image->set_attr("subvolume_full_nz",nz);
00490 }
00491
00492 if (tmp_data) delete tmp_data;
00493 tmp_data = new EMData();
00494 tmp_data->set_size(subnx/2, subny, subnz);
00495 tmp_data->to_value(seed_weight);
00496
00497 load_inserter();
00498
00499 if ( (bool) params["quiet"] == false )
00500 {
00501 cout << "Seeded direct Fourier inversion";
00502 cout << "3D Fourier dimensions are " << nx << " " << ny << " " << nz << endl;
00503 cout << "3D Fourier subvolume is " << subnx << " " << subny << " " << subnz << endl;
00504 cout << "You will require approximately " << setprecision(3) << (subnx*subny*subnz*sizeof(float)*1.5)/1000000000.0 << "GB of memory to reconstruct this volume" << endl;
00505 }
00506 }
00507
00508
00509 EMData* FourierReconstructor::preprocess_slice( const EMData* const slice, const Transform& t )
00510 {
00511
00512 EMData* return_slice = 0;
00513 Transform tmp(t);
00514 tmp.set_rotation(Dict("type","eman"));
00515
00516 if (tmp.is_identity()) return_slice=slice->copy();
00517 else return_slice = slice->process("xform",Dict("transform",&tmp));
00518
00519 return_slice->process_inplace("xform.phaseorigin.tocorner");
00520
00521
00522
00523 return_slice->do_fft_inplace();
00524
00525
00526
00527 return_slice->mult((float)sqrt(1.0f/(return_slice->get_ysize())*return_slice->get_xsize()));
00528
00529
00530
00531
00532 return_slice->set_attr("reconstruct_preproc",(int)1);
00533 return return_slice;
00534 }
00535
00536
00537 int FourierReconstructor::insert_slice(const EMData* const input_slice, const Transform & arg, const float weight)
00538 {
00539
00540 if (!input_slice) throw NullPointerException("EMData pointer (input image) is NULL");
00541
00542 Transform * rotation;
00543
00544
00545
00546 rotation = new Transform(arg);
00547
00548
00549 EMData *slice;
00550 if (input_slice->get_attr_default("reconstruct_preproc",(int) 0)) slice=input_slice->copy();
00551 else slice = preprocess_slice( input_slice, *rotation);
00552
00553
00554
00555
00556 rotation->set_scale(1.0);
00557 rotation->set_mirror(false);
00558 rotation->set_trans(0,0,0);
00559
00560
00561 do_insert_slice_work(slice, *rotation, weight);
00562
00563 delete rotation; rotation=0;
00564 delete slice;
00565
00566
00567 return 0;
00568 }
00569
00570 void FourierReconstructor::do_insert_slice_work(const EMData* const input_slice, const Transform & arg,const float weight)
00571 {
00572
00573
00574
00575
00576
00577
00578
00579
00580
00581
00582 vector<Transform> syms = Symmetry3D::get_symmetries((string)params["sym"]);
00583
00584 float inx=(float)(input_slice->get_xsize());
00585 float iny=(float)(input_slice->get_ysize());
00586
00587 for ( vector<Transform>::const_iterator it = syms.begin(); it != syms.end(); ++it ) {
00588 Transform t3d = arg*(*it);
00589 for (int y = -iny/2; y < iny/2; y++) {
00590 for (int x = 0; x <= inx/2; x++) {
00591
00592 float rx = (float) x/(inx-2.0f);
00593 float ry = (float) y/iny;
00594
00595 Vec3f coord(rx,ry,0);
00596 coord = coord*t3d;
00597 float xx = coord[0];
00598 float yy = coord[1];
00599 float zz = coord[2];
00600
00601
00602 xx=xx*(nx-2);
00603 yy=yy*ny;
00604 zz=zz*nz;
00605
00606
00607
00608
00609
00610
00611
00612
00613
00614 inserter->insert_pixel(xx,yy,zz,input_slice->get_complex_at(x,y),weight);
00615 }
00616 }
00617 }
00618 }
00619
00620 int FourierReconstructor::determine_slice_agreement(EMData* input_slice, const Transform & arg, const float weight,bool sub)
00621 {
00622
00623 if (!input_slice) throw NullPointerException("EMData pointer (input image) is NULL");
00624
00625 Transform * rotation;
00626 rotation = new Transform(arg);
00627
00628 EMData *slice;
00629 if (input_slice->get_attr_default("reconstruct_preproc",(int) 0)) slice=input_slice->copy();
00630 else slice = preprocess_slice( input_slice, *rotation);
00631
00632
00633
00634
00635 rotation->set_scale(1.0);
00636 rotation->set_mirror(false);
00637 rotation->set_trans(0,0,0);
00638
00639
00640 if (sub) do_insert_slice_work(slice, *rotation, -weight);
00641
00642
00643 do_compare_slice_work(slice, *rotation,weight);
00644
00645 input_slice->set_attr("reconstruct_norm",slice->get_attr("reconstruct_norm"));
00646 input_slice->set_attr("reconstruct_absqual",slice->get_attr("reconstruct_absqual"));
00647
00648 input_slice->set_attr("reconstruct_weight",slice->get_attr("reconstruct_weight"));
00649
00650
00651 if (sub) do_insert_slice_work(slice, *rotation, weight);
00652
00653
00654 delete rotation;
00655 delete slice;
00656
00657
00658 return 0;
00659
00660 }
00661
00662 void FourierReconstructor::do_compare_slice_work(EMData* input_slice, const Transform & arg,float weight)
00663 {
00664
00665 float dt[3];
00666 float dt2[2];
00667 float *dat = input_slice->get_data();
00668 vector<Transform> syms = Symmetry3D::get_symmetries((string)params["sym"]);
00669
00670 float inx=(float)(input_slice->get_xsize());
00671 float iny=(float)(input_slice->get_ysize());
00672
00673 double dot=0;
00674 double vweight=0;
00675 double power=0;
00676 double power2=0;
00677 for ( vector<Transform>::const_iterator it = syms.begin(); it != syms.end(); ++it ) {
00678 Transform t3d = arg*(*it);
00679 for (int y = -iny/2; y < iny/2; y++) {
00680 for (int x = 0; x <= inx/2; x++) {
00681 if (x==0 && y==0) continue;
00682
00683 float rx = (float) x/(inx-2);
00684 float ry = (float) y/iny;
00685
00686
00687
00688
00689 Vec3f coord(rx,ry,0);
00690 coord = coord*t3d;
00691 float xx = coord[0];
00692 float yy = coord[1];
00693 float zz = coord[2];
00694
00695
00696 if (fabs(xx)>0.5 || fabs(yy)>=0.5 || fabs(zz)>=0.5) continue;
00697
00698
00699 xx=xx*(nx-2);
00700 yy=yy*ny;
00701 zz=zz*nz;
00702
00703
00704 int idx = (int)(x * 2 + inx*(y<0?iny+y:y));
00705 dt2[0] = dat[idx];
00706 dt2[1] = dat[idx+1];
00707
00708
00709 if (!pixel_at(xx,yy,zz,dt) || dt[2]==0) continue;
00710
00711
00712 dot+=(dt[0]*dt2[0]+dt[1]*dt2[1])*dt[2];
00713 vweight+=dt[2];
00714 power+=(dt[0]*dt[0]+dt[1]*dt[1])*dt[2];
00715 power2+=(dt2[0]*dt2[0]+dt2[1]*dt2[1])*dt[2];
00716 }
00717 }
00718 }
00719
00720 dot/=sqrt(power*power2);
00721
00722 input_slice->set_attr("reconstruct_norm",(float)(power2<=0?1.0:sqrt(power/power2)));
00723 input_slice->set_attr("reconstruct_absqual",(float)dot);
00724 float rw=weight<=0?1.0f:1.0f/weight;
00725 input_slice->set_attr("reconstruct_qual",(float)(dot*rw/((rw-1.0)*dot+1.0)));
00726 input_slice->set_attr("reconstruct_weight",(float)vweight/(float)(subnx*subny*subnz));
00727
00728
00729 }
00730
00731 bool FourierReconstructor::pixel_at(const float& xx, const float& yy, const float& zz, float *dt)
00732 {
00733 int x0 = (int) floor(xx);
00734 int y0 = (int) floor(yy);
00735 int z0 = (int) floor(zz);
00736
00737 float *rdata=image->get_data();
00738 float *norm=tmp_data->get_data();
00739 float normsum=0,normsum2=0;
00740
00741 dt[0]=dt[1]=dt[2]=0.0;
00742
00743 if (nx==subnx) {
00744 if (x0<-nx2-1 || y0<-ny2-1 || z0<-nz2-1 || x0>nx2 || y0>ny2 || z0>nz2 ) return false;
00745
00746
00747 int x1=x0+1;
00748 int y1=y0+1;
00749 int z1=z0+1;
00750 if (x0<-nx2) x0=-nx2;
00751 if (x1>nx2) x1=nx2;
00752 if (y0<-ny2) y0=-ny2;
00753 if (y1>ny2) y1=ny2;
00754 if (z0<-nz2) z0=-nz2;
00755 if (z1>nz2) z1=nz2;
00756
00757 size_t idx=0;
00758 float r, gg;
00759 for (int k = z0 ; k <= z1; k++) {
00760 for (int j = y0 ; j <= y1; j++) {
00761 for (int i = x0; i <= x1; i ++) {
00762 r = Util::hypot3sq((float) i - xx, j - yy, k - zz);
00763 idx=image->get_complex_index_fast(i,j,k);
00764 gg = Util::fast_exp(-r / EMConsts::I2G);
00765
00766 dt[0]+=gg*rdata[idx];
00767 dt[1]+=(i<0?-1.0f:1.0f)*gg*rdata[idx+1];
00768 dt[2]+=norm[idx/2]*gg;
00769 normsum2+=gg;
00770 normsum+=gg*norm[idx/2];
00771 }
00772 }
00773 }
00774 if (normsum==0) return false;
00775 dt[0]/=normsum;
00776 dt[1]/=normsum;
00777 dt[2]/=normsum2;
00778
00779 return true;
00780 }
00781 else {
00782 size_t idx;
00783 float r, gg;
00784 for (int k = z0 ; k <= z0 + 1; k++) {
00785 for (int j = y0 ; j <= y0 + 1; j++) {
00786 for (int i = x0; i <= x0 + 1; i ++) {
00787 r = Util::hypot3sq((float) i - xx, j - yy, k - zz);
00788 idx=image->get_complex_index(i,j,k,subx0,suby0,subz0,nx,ny,nz);
00789 gg = Util::fast_exp(-r / EMConsts::I2G)*norm[idx/2];
00790
00791 dt[0]+=gg*rdata[idx];
00792 dt[1]+=(i<0?-1.0f:1.0f)*gg*rdata[idx+1];
00793 dt[2]+=norm[idx/2];
00794 normsum+=gg;
00795 }
00796 }
00797 }
00798
00799 if (normsum==0) return false;
00800 return true;
00801 }
00802 }
00803
00804
00805 EMData *FourierReconstructor::finish(bool doift)
00806 {
00807
00808
00809
00810 bool sqrtnorm=params.set_default("sqrtnorm",false);
00811 normalize_threed(sqrtnorm);
00812
00813
00814
00815
00816
00817
00818
00819
00820
00821
00822
00823
00824
00825
00826
00827 if (doift) {
00828 image->do_ift_inplace();
00829 image->depad();
00830 image->process_inplace("xform.phaseorigin.tocenter");
00831 }
00832
00833
00834
00835
00836
00837
00838
00839
00840
00841
00842
00843
00844
00845
00846
00847
00848 image->update();
00849
00850 if (params.has_key("savenorm") && strlen((const char *)params["savenorm"])>0) {
00851 if (tmp_data->get_ysize()%2==0 && tmp_data->get_zsize()%2==0) tmp_data->process_inplace("xform.fourierorigin.tocenter");
00852 tmp_data->write_image((const char *)params["savenorm"]);
00853 }
00854
00855 delete tmp_data;
00856 tmp_data=0;
00857 EMData *ret=image;
00858 image=0;
00859
00860 return ret;
00861 }
00862
00863 int WienerFourierReconstructor::insert_slice(const EMData* const input_slice, const Transform & arg, const float weight)
00864 {
00865
00866 if (!input_slice) throw NullPointerException("EMData pointer (input image) is NULL");
00867
00868 Transform * rotation;
00869
00870
00871
00872 rotation = new Transform(arg);
00873
00874
00875 if (!input_slice->has_attr("ctf_snr_total"))
00876 throw NotExistingObjectException("ctf_snr_total","No SNR information present in class-average. Must use the ctf.auto or ctfw.auto averager.");
00877
00878 EMData *slice;
00879 if (input_slice->get_attr_default("reconstruct_preproc",(int) 0)) slice=input_slice->copy();
00880 else slice = preprocess_slice( input_slice, *rotation);
00881
00882
00883
00884
00885 rotation->set_scale(1.0);
00886 rotation->set_mirror(false);
00887 rotation->set_trans(0,0,0);
00888
00889
00890 do_insert_slice_work(slice, *rotation, weight);
00891
00892 delete rotation; rotation=0;
00893 delete slice;
00894
00895
00896 return 0;
00897 }
00898
00899 void WienerFourierReconstructor::do_insert_slice_work(const EMData* const input_slice, const Transform & arg,const float inweight)
00900 {
00901
00902 vector<Transform> syms = Symmetry3D::get_symmetries((string)params["sym"]);
00903
00904 float inx=(float)(input_slice->get_xsize());
00905 float iny=(float)(input_slice->get_ysize());
00906
00907 int undo_wiener=(int)input_slice->get_attr_default("ctf_wiener_filtered",0);
00908
00909
00910 vector<float> snr=input_slice->get_attr("ctf_snr_total");
00911 float sub=1.0;
00912 if (inweight<0) sub=-1.0;
00913 float weight;
00914
00915 for ( vector<Transform>::const_iterator it = syms.begin(); it != syms.end(); ++it ) {
00916 Transform t3d = arg*(*it);
00917 for (int y = -iny/2; y < iny/2; y++) {
00918 for (int x = 0; x <= inx/2; x++) {
00919
00920 float rx = (float) x/(inx-2.0f);
00921 float ry = (float) y/iny;
00922
00923
00924 float rn = hypot(rx,ry);
00925 if (rn>=.5) continue;
00926 rn*=snr.size()*2.0;
00927 int rni=(int)floor(rn);
00928 if (rni>=snr.size()-1) weight=snr[snr.size()-1]*sub;
00929 else {
00930 rn-=rni;
00931 weight=Util::linear_interpolate(snr[rni],snr[rni+1],rn);
00932 }
00933
00934
00935 Vec3f coord(rx,ry,0);
00936 coord = coord*t3d;
00937 float xx = coord[0];
00938 float yy = coord[1];
00939 float zz = coord[2];
00940
00941
00942 xx=xx*(nx-2);
00943 yy=yy*ny;
00944 zz=zz*nz;
00945
00946
00947 if (undo_wiener) inserter->insert_pixel(xx,yy,zz,(input_slice->get_complex_at(x,y))*((weight+1.0f)/weight),weight*sub);
00948 else inserter->insert_pixel(xx,yy,zz,input_slice->get_complex_at(x,y),weight*sub);
00949 }
00950 }
00951 }
00952 }
00953
00954 int WienerFourierReconstructor::determine_slice_agreement(EMData* input_slice, const Transform & arg, const float weight,bool sub)
00955 {
00956
00957 if (!input_slice) throw NullPointerException("EMData pointer (input image) is NULL");
00958
00959 Transform * rotation;
00960 rotation = new Transform(arg);
00961
00962 EMData *slice;
00963 if (input_slice->get_attr_default("reconstruct_preproc",(int) 0)) slice=input_slice->copy();
00964 else slice = preprocess_slice( input_slice, *rotation);
00965
00966
00967
00968
00969 rotation->set_scale(1.0);
00970 rotation->set_mirror(false);
00971 rotation->set_trans(0,0,0);
00972
00973 tmp_data->write_image("dbug.hdf",0);
00974
00975
00976 if (sub) do_insert_slice_work(slice, *rotation, -weight);
00977
00978
00979 do_compare_slice_work(slice, *rotation,weight);
00980
00981 input_slice->set_attr("reconstruct_norm",slice->get_attr("reconstruct_norm"));
00982 input_slice->set_attr("reconstruct_absqual",slice->get_attr("reconstruct_absqual"));
00983
00984 input_slice->set_attr("reconstruct_weight",slice->get_attr("reconstruct_weight"));
00985
00986
00987 if (sub) do_insert_slice_work(slice, *rotation, weight);
00988
00989
00990 delete rotation;
00991 delete slice;
00992
00993
00994 return 0;
00995
00996 }
00997
00998 void WienerFourierReconstructor::do_compare_slice_work(EMData* input_slice, const Transform & arg,float weight)
00999 {
01000
01001 float dt[3];
01002 float dt2[2];
01003 float *dat = input_slice->get_data();
01004 vector<Transform> syms = Symmetry3D::get_symmetries((string)params["sym"]);
01005
01006 float inx=(float)(input_slice->get_xsize());
01007 float iny=(float)(input_slice->get_ysize());
01008
01009 double dot=0;
01010 double vweight=0;
01011 double power=0;
01012 double power2=0;
01013 for ( vector<Transform>::const_iterator it = syms.begin(); it != syms.end(); ++it ) {
01014 Transform t3d = arg*(*it);
01015 for (int y = -iny/2; y < iny/2; y++) {
01016 for (int x = 0; x <= inx/2; x++) {
01017 if (x==0 && y==0) continue;
01018
01019 float rx = (float) x/(inx-2);
01020 float ry = (float) y/iny;
01021
01022
01023
01024
01025 Vec3f coord(rx,ry,0);
01026 coord = coord*t3d;
01027 float xx = coord[0];
01028 float yy = coord[1];
01029 float zz = coord[2];
01030
01031
01032 if (fabs(xx)>0.5 || fabs(yy)>=0.5 || fabs(zz)>=0.5) continue;
01033
01034
01035 xx=xx*(nx-2);
01036 yy=yy*ny;
01037 zz=zz*nz;
01038
01039
01040 int idx = (int)(x * 2 + inx*(y<0?iny+y:y));
01041 dt2[0] = dat[idx];
01042 dt2[1] = dat[idx+1];
01043
01044
01045 if (!pixel_at(xx,yy,zz,dt) || dt[2]<=0) continue;
01046
01047
01048 dot+=(dt[0]*dt2[0]+dt[1]*dt2[1])*dt[2];
01049 vweight+=dt[2];
01050 power+=(dt[0]*dt[0]+dt[1]*dt[1])*dt[2];
01051 power2+=(dt2[0]*dt2[0]+dt2[1]*dt2[1])*dt[2];
01052 }
01053 }
01054 }
01055
01056 dot/=sqrt(power*power2);
01057
01058 input_slice->set_attr("reconstruct_norm",(float)(power2<=0?1.0:sqrt(power/power2)));
01059 input_slice->set_attr("reconstruct_absqual",(float)dot);
01060 float rw=weight<=0?1.0f:1.0f/weight;
01061 input_slice->set_attr("reconstruct_qual",(float)(dot*rw/((rw-1.0)*dot+1.0)));
01062 input_slice->set_attr("reconstruct_weight",(float)vweight/(float)(subnx*subny*subnz));
01063
01064
01065 }
01066
01067 bool WienerFourierReconstructor::pixel_at(const float& xx, const float& yy, const float& zz, float *dt)
01068 {
01069 int x0 = (int) floor(xx);
01070 int y0 = (int) floor(yy);
01071 int z0 = (int) floor(zz);
01072
01073 float *rdata=image->get_data();
01074 float *norm=tmp_data->get_data();
01075 float normsum=0,normsum2=0;
01076
01077 dt[0]=dt[1]=dt[2]=0.0;
01078
01079 if (nx==subnx) {
01080 if (x0<-nx2-1 || y0<-ny2-1 || z0<-nz2-1 || x0>nx2 || y0>ny2 || z0>nz2 ) return false;
01081
01082
01083 int x1=x0+1;
01084 int y1=y0+1;
01085 int z1=z0+1;
01086 if (x0<-nx2) x0=-nx2;
01087 if (x1>nx2) x1=nx2;
01088 if (y0<-ny2) y0=-ny2;
01089 if (y1>ny2) y1=ny2;
01090 if (z0<-nz2) z0=-nz2;
01091 if (z1>nz2) z1=nz2;
01092
01093 size_t idx=0;
01094 float r, gg;
01095 for (int k = z0 ; k <= z1; k++) {
01096 for (int j = y0 ; j <= y1; j++) {
01097 for (int i = x0; i <= x1; i ++) {
01098 r = Util::hypot3sq((float) i - xx, j - yy, k - zz);
01099 idx=image->get_complex_index_fast(i,j,k);
01100 gg = Util::fast_exp(-r / EMConsts::I2G);
01101
01102 dt[0]+=gg*rdata[idx];
01103 dt[1]+=(i<0?-1.0f:1.0f)*gg*rdata[idx+1];
01104 dt[2]+=norm[idx/2]*gg;
01105 normsum2+=gg;
01106 normsum+=gg*norm[idx/2];
01107 }
01108 }
01109 }
01110 if (normsum==0) return false;
01111 dt[0]/=normsum;
01112 dt[1]/=normsum;
01113 dt[2]/=normsum2;
01114
01115 return true;
01116 }
01117 else {
01118 size_t idx;
01119 float r, gg;
01120 for (int k = z0 ; k <= z0 + 1; k++) {
01121 for (int j = y0 ; j <= y0 + 1; j++) {
01122 for (int i = x0; i <= x0 + 1; i ++) {
01123 r = Util::hypot3sq((float) i - xx, j - yy, k - zz);
01124 idx=image->get_complex_index(i,j,k,subx0,suby0,subz0,nx,ny,nz);
01125 gg = Util::fast_exp(-r / EMConsts::I2G)*norm[idx/2];
01126
01127 dt[0]+=gg*rdata[idx];
01128 dt[1]+=(i<0?-1.0f:1.0f)*gg*rdata[idx+1];
01129 dt[2]+=norm[idx/2];
01130 normsum+=gg;
01131 }
01132 }
01133 }
01134
01135 if (normsum==0) return false;
01136 return true;
01137 }
01138 }
01139
01140
01141 EMData *WienerFourierReconstructor::finish(bool doift)
01142 {
01143
01144 bool sqrtnorm=params.set_default("sqrtnorm",false);
01145 normalize_threed(sqrtnorm,true);
01146
01147 if (doift) {
01148 image->do_ift_inplace();
01149 image->depad();
01150 image->process_inplace("xform.phaseorigin.tocenter");
01151 }
01152
01153 image->update();
01154
01155 if (params.has_key("savenorm") && strlen((const char *)params["savenorm"])>0) {
01156 if (tmp_data->get_ysize()%2==0 && tmp_data->get_zsize()%2==0) tmp_data->process_inplace("xform.fourierorigin.tocenter");
01157 tmp_data->write_image((const char *)params["savenorm"]);
01158 }
01159
01160 delete tmp_data;
01161 tmp_data=0;
01162 EMData *ret=image;
01163 image=0;
01164
01165 return ret;
01166 }
01167
01168
01169
01170
01171
01172
01173
01174
01175
01176
01177
01178
01179
01180
01181
01182
01183
01184
01185
01186
01187
01188
01189
01190
01191
01192
01193
01194
01195
01196
01197
01198
01199
01200
01201
01202
01203
01204
01205
01206
01207
01208
01209
01210
01211
01212
01213
01214
01215
01216
01217
01218
01219
01220
01221
01222
01223
01224
01225
01226
01227
01228
01229
01230
01231
01232
01233
01234
01235
01236
01237
01238
01239
01240
01241
01242
01243
01244
01245
01246
01247
01248
01249
01250
01251
01252
01253
01254
01255
01256
01257
01258
01259
01260
01261
01262
01263
01264
01265
01266
01267
01268
01269
01270
01271
01272
01273
01274
01275
01276
01277
01278
01279
01280
01281
01282
01283
01284
01285
01286
01287
01288
01289
01290
01291
01292
01293
01294
01295
01296
01297
01298
01299
01300
01301
01302
01303
01304
01305
01306
01307
01308
01309
01310
01311
01312
01313
01314
01315
01316
01317
01318
01319
01320
01321
01322
01323
01324
01325
01326
01327
01328
01329
01330
01331
01332
01333
01334
01335
01336
01337
01338
01339
01340
01341
01342
01343
01344
01345
01346
01347
01348
01349
01350
01351
01352
01353
01354
01355
01356
01357
01358
01359
01360
01361
01362
01363
01364
01365
01366
01367
01368
01369
01370
01371
01372
01373
01374
01375
01376
01377
01378
01379
01380
01381
01382
01383
01384
01385
01386
01387
01388
01389
01390
01391
01392
01393
01394
01395
01396
01397
01398
01399
01400
01401
01402
01403
01404
01405
01406
01407
01408
01409
01410
01411
01412
01413
01414
01415
01416
01417
01418
01419
01420
01421
01422
01423
01424
01425
01426
01427
01428
01429
01430
01431
01432
01433
01434
01435
01436
01437
01438
01439
01440
01441
01442
01443
01444
01445
01446
01447
01448
01449
01450
01451
01452
01453
01454
01455
01456
01457
01458
01459
01460
01461
01462
01463
01464
01465
01466
01467
01468
01469
01470
01471
01472
01473
01474
01475
01476
01477
01478
01479
01480
01481
01482
01483
01484
01485
01486
01487
01488
01489
01490
01491
01492
01493
01494
01495
01496
01497
01498
01499
01500
01501
01502
01503
01504
01505
01506
01507
01508
01509
01510
01511
01512
01513
01514
01515
01516
01517
01518
01519
01520
01521
01522
01523
01524
01525
01526
01527
01528
01529
01530
01531
01532
01533
01534
01535
01536
01537
01538
01539
01540
01541
01542
01543
01544
01545
01546
01547
01548
01549
01550
01551
01552
01553
01554
01555
01556
01557
01558
01559
01560
01561
01562
01563
01564
01565
01566
01567
01568
01569
01570
01571
01572
01573
01574
01575
01576
01577
01578
01579
01580
01581
01582
01583
01584
01585
01586
01587
01588
01589
01590
01591
01592
01593
01594
01595
01596
01597
01598
01599
01600
01601
01602
01603
01604
01605
01606
01607
01608
01609
01610
01611
01612
01613
01614
01615
01616
01617
01618
01619
01620
01621
01622
01623
01624
01625
01626
01627
01628
01629
01630
01631
01632
01633
01634
01635
01636
01637
01638
01639
01640
01641
01642
01643
01644
01645
01646
01647
01648
01649
01650
01651
01652
01653
01654
01655
01656
01657
01658
01659
01660
01661
01662
01663
01664
01665
01666
01667
01668
01669
01670
01671
01672
01673
01674
01675
01676
01677
01678
01679
01680
01681
01682
01683
01684
01685
01686
01687
01688
01689
01690
01691
01692
01693
01694
01695
01696
01697
01698
01699
01700
01701
01702
01703
01704
01705
01706
01707
01708
01709
01710
01711
01712
01713
01714
01715
01716
01717
01718
01719
01720
01721
01722
01723
01724
01725
01726
01727
01728
01729 void BackProjectionReconstructor::setup()
01730 {
01731 int size = params["size"];
01732 image = new EMData();
01733 nx = size;
01734 ny = size;
01735 if ( (int) params["zsample"] != 0 ) nz = params["zsample"];
01736 else nz = size;
01737 image->set_size(nx, ny, nz);
01738 }
01739
01740 EMData* BackProjectionReconstructor::preprocess_slice(const EMData* const slice, const Transform& t)
01741 {
01742
01743 EMData* return_slice = slice->process("normalize.edgemean");
01744 return_slice->process_inplace("filter.linearfourier");
01745
01746 Transform tmp(t);
01747 tmp.set_rotation(Dict("type","eman"));
01748 Vec2f trans = tmp.get_trans_2d();
01749 float scale = tmp.get_scale();
01750 bool mirror = tmp.get_mirror();
01751 if (trans[0] != 0 || trans[1] != 0 || scale != 1.0 ) {
01752 return_slice->transform(tmp);
01753 } else if ( mirror == true ) {
01754 return_slice = slice->process("xform.flip",Dict("axis","x"));
01755 }
01756
01757 return return_slice;
01758 }
01759
01760 int BackProjectionReconstructor::insert_slice(const EMData* const input, const Transform &t, const float sliceweight)
01761 {
01762 if (!input) {
01763 LOGERR("try to insert NULL slice");
01764 return 1;
01765 }
01766
01767 if (input->get_xsize() != input->get_ysize() || input->get_xsize() != nx) {
01768 LOGERR("tried to insert image that was not correction dimensions");
01769 return 1;
01770 }
01771
01772 Transform * transform;
01773 if ( input->has_attr("xform.projection") ) {
01774 transform = (Transform*) (input->get_attr("xform.projection"));
01775 } else {
01776 transform = new Transform(t);
01777 }
01778 EMData* slice = preprocess_slice(input, t);
01779
01780 float weight = params["weight"];
01781 slice->mult(weight);
01782
01783 EMData *tmp = new EMData();
01784 tmp->set_size(nx, ny, nz);
01785
01786 float *slice_data = slice->get_data();
01787 float *tmp_data = tmp->get_data();
01788
01789 size_t nxy = nx * ny;
01790 size_t nxy_size = nxy * sizeof(float);;
01791 for (int i = 0; i < nz; ++i) {
01792 memcpy(&tmp_data[nxy * i], slice_data, nxy_size);
01793 }
01794
01795 transform->set_scale(1.0);
01796 transform->set_mirror(false);
01797 transform->set_trans(0,0,0);
01798 transform->invert();
01799
01800 tmp->transform(*transform);
01801 image->add(*tmp);
01802
01803 if(transform) {delete transform; transform=0;}
01804 delete tmp;
01805 delete slice;
01806
01807 return 0;
01808 }
01809
01810 EMData *BackProjectionReconstructor::finish(bool doift)
01811 {
01812
01813 Symmetry3D* sym = Factory<Symmetry3D>::get((string)params["sym"]);
01814 vector<Transform> syms = sym->get_syms();
01815
01816 for ( vector<Transform>::const_iterator it = syms.begin(); it != syms.end(); ++it ) {
01817
01818 EMData tmpcopy(*image);
01819 tmpcopy.transform(*it);
01820 image->add(tmpcopy);
01821 }
01822
01823 image->mult(1.0f/(float)sym->get_nsym());
01824 delete sym;
01825 return image;
01826 }
01827
01828 EMData* EMAN::padfft_slice( const EMData* const slice, const Transform& t, int npad )
01829 {
01830 int nx = slice->get_xsize();
01831 int ny = slice->get_ysize();
01832 int ndim = (ny==1) ? 1 : 2;
01833
01834 if( ndim==2 && nx!=ny )
01835 {
01836
01837 throw std::runtime_error("Tried to padfft a 2D slice which is not square.");
01838 }
01839
01840
01841 EMData* temp = slice->average_circ_sub();
01842
01843 Assert( temp != NULL );
01844 EMData* zeropadded = temp->norm_pad( false, npad );
01845 Assert( zeropadded != NULL );
01846 checked_delete( temp );
01847
01848 zeropadded->do_fft_inplace();
01849 EMData* padfftslice = zeropadded;
01850
01851
01852 Vec2f trans = t.get_trans_2d();
01853 float sx = -trans[0];
01854 float sy = -trans[1];
01855 if(sx != 0.0f || sy != 0.0)
01856 padfftslice->process_inplace("filter.shift", Dict("x_shift", sx, "y_shift", sy, "z_shift", 0.0f));
01857
01858 int remove = slice->get_attr_default("remove", 0);
01859 padfftslice->set_attr( "remove", remove );
01860
01861
01862
01863 padfftslice->center_origin_fft();
01864 return padfftslice;
01865 }
01866
01867 nn4Reconstructor::nn4Reconstructor()
01868 {
01869 m_volume = NULL;
01870 m_wptr = NULL;
01871 m_result = NULL;
01872 }
01873
01874 nn4Reconstructor::nn4Reconstructor( const string& symmetry, int size, int npad )
01875 {
01876 m_volume = NULL;
01877 m_wptr = NULL;
01878 m_result = NULL;
01879 setup( symmetry, size, npad );
01880 load_default_settings();
01881 print_params();
01882 }
01883
01884 nn4Reconstructor::~nn4Reconstructor()
01885 {
01886 if( m_delete_volume )
01887 checked_delete(m_volume);
01888
01889 if( m_delete_weight )
01890 checked_delete( m_wptr );
01891
01892 checked_delete( m_result );
01893 }
01894
01895 enum weighting_method { NONE, ESTIMATE, VORONOI };
01896
01897 float max2d( int kc, const vector<float>& pow_a )
01898 {
01899 float max = 0.0;
01900 for( int i=-kc; i <= kc; ++i ) {
01901 for( int j=-kc; j <= kc; ++j ) {
01902 if( i==0 && j==0 ) continue;
01903 {
01904 int c = 2*kc+1 - std::abs(i) - std::abs(j);
01905 max = max + pow_a[c];
01906 }
01907 }
01908 }
01909 return max;
01910 }
01911
01912 float max3d( int kc, const vector<float>& pow_a )
01913 {
01914 float max = 0.0;
01915 for( int i=-kc; i <= kc; ++i ) {
01916 for( int j=-kc; j <= kc; ++j ) {
01917 for( int k=-kc; k <= kc; ++k ) {
01918 if( i==0 && j==0 && k==0 ) continue;
01919
01920 {
01921 int c = 3*kc+1 - std::abs(i) - std::abs(j) - std::abs(k);
01922 max = max + pow_a[c];
01923
01924
01925 }
01926 }
01927 }
01928 }
01929 return max;
01930 }
01931
01932
01933 void nn4Reconstructor::setup()
01934 {
01935 int size = params["size"];
01936 int npad = params["npad"];
01937
01938
01939 string symmetry;
01940 if( params.has_key("symmetry") ) {
01941 symmetry = params["symmetry"].to_str();
01942 } else {
01943 symmetry = "c1";
01944 }
01945
01946 if( params.has_key("ndim") ) {
01947 m_ndim = params["ndim"];
01948 } else {
01949 m_ndim = 3;
01950 }
01951
01952 if( params.has_key( "snr" ) ) {
01953 m_osnr = 1.0f/float( params["snr"] );
01954 } else {
01955 m_osnr = 0.0;
01956 }
01957
01958 setup( symmetry, size, npad );
01959 }
01960
01961 void nn4Reconstructor::setup( const string& symmetry, int size, int npad )
01962 {
01963 m_weighting = ESTIMATE;
01964 m_wghta = 0.2f;
01965
01966 m_symmetry = symmetry;
01967 m_npad = npad;
01968 m_nsym = Transform::get_nsym(m_symmetry);
01969
01970 m_vnx = size;
01971 m_vny = size;
01972 m_vnz = (m_ndim==3) ? size : 1;
01973
01974 m_vnxp = size*npad;
01975 m_vnyp = size*npad;
01976 m_vnzp = (m_ndim==3) ? size*npad : 1;
01977
01978 m_vnxc = m_vnxp/2;
01979 m_vnyc = m_vnyp/2;
01980 m_vnzc = (m_ndim==3) ? m_vnzp/2 : 1;
01981
01982 buildFFTVolume();
01983 buildNormVolume();
01984
01985 }
01986
01987
01988 void nn4Reconstructor::buildFFTVolume() {
01989 int offset = 2 - m_vnxp%2;
01990
01991 if( params.has_key("fftvol") ) {
01992 m_volume = params["fftvol"];
01993 m_delete_volume = false;
01994 } else {
01995 m_volume = new EMData();
01996 m_delete_volume = true;
01997 }
01998
01999 if( m_volume->get_xsize() != m_vnxp+offset &&
02000 m_volume->get_ysize() != m_vnyp &&
02001 m_volume->get_zsize() != m_vnzp ) {
02002 m_volume->set_size(m_vnxp+offset,m_vnyp,m_vnzp);
02003 m_volume->to_zero();
02004 }
02005
02006
02007
02008
02009 if ( m_vnxp % 2 == 0 ) m_volume->set_fftodd(0);
02010 else m_volume->set_fftodd(1);
02011
02012
02013 m_volume->set_nxc(m_vnxp/2);
02014 m_volume->set_complex(true);
02015 m_volume->set_ri(true);
02016 m_volume->set_fftpad(true);
02017 m_volume->set_attr("npad", m_npad);
02018 m_volume->set_array_offsets(0,1,1);
02019 }
02020
02021 void nn4Reconstructor::buildNormVolume() {
02022
02023 if( params.has_key("weight") ) {
02024 m_wptr = params["weight"];
02025 m_delete_weight = false;
02026 } else {
02027 m_wptr = new EMData();
02028 m_delete_weight = true;
02029 }
02030
02031 if( m_wptr->get_xsize() != m_vnxc+1 &&
02032 m_wptr->get_ysize() != m_vnyp &&
02033 m_wptr->get_zsize() != m_vnzp ) {
02034 m_wptr->set_size(m_vnxc+1,m_vnyp,m_vnzp);
02035 m_wptr->to_zero();
02036 }
02037
02038 m_wptr->set_array_offsets(0,1,1);
02039 }
02040
02041 void printImage( const EMData* line )
02042 {
02043 Assert( line->get_zsize()==1 );
02044
02045
02046 int nx = line->get_xsize();
02047 int ny = line->get_ysize();
02048 for( int j=0; j < ny; ++j ) {
02049 for( int i=0; i < nx; ++i ) printf( "%10.3f ", line->get_value_at(i,j) );
02050 printf( "\n" );
02051 }
02052 }
02053
02054
02055
02056 int nn4Reconstructor::insert_slice(const EMData* const slice, const Transform& t, const float weight) {
02057
02058 if (!slice) {
02059 LOGERR("try to insert NULL slice");
02060 return 1;
02061 }
02062
02063 int padffted= slice->get_attr_default( "padffted", 0 );
02064 if( m_ndim==3 ) {
02065 if ( padffted==0 && (slice->get_xsize()!=slice->get_ysize() || slice->get_xsize()!=m_vnx) ) {
02066
02067 LOGERR("Tried to insert a slice that is the wrong size.");
02068 return 1;
02069 }
02070 } else {
02071 Assert( m_ndim==2 );
02072 if( slice->get_ysize() !=1 ) {
02073 LOGERR( "for 2D reconstruction, a line is excepted" );
02074 return 1;
02075 }
02076 }
02077
02078 EMData* padfft = NULL;
02079
02080 if( padffted != 0 ) padfft = new EMData(*slice);
02081 else padfft = padfft_slice( slice, t, m_npad );
02082
02083 int mult= slice->get_attr_default( "mult", 1 );
02084 Assert( mult > 0 );
02085
02086 if( m_ndim==3 ) {
02087 insert_padfft_slice( padfft, t, mult );
02088 } else {
02089 float alpha = padfft->get_attr( "alpha" );
02090 alpha = alpha/180.0f*M_PI;
02091 for(int i=0; i < m_vnxc+1; ++i ) {
02092 float xnew = i*cos(alpha);
02093 float ynew = -i*sin(alpha);
02094 float btqr = padfft->get_value_at( 2*i, 0, 0 );
02095 float btqi = padfft->get_value_at( 2*i+1, 0, 0 );
02096 if( xnew < 0.0 ) {
02097 xnew *= -1;
02098 ynew *= -1;
02099 btqi *= -1;
02100 }
02101
02102 int ixn = int(xnew+0.5+m_vnxp) - m_vnxp;
02103 int iyn = int(ynew+0.5+m_vnyp) - m_vnyp;
02104
02105 if(iyn < 0 ) iyn += m_vnyp;
02106
02107 (*m_volume)( 2*ixn, iyn+1, 1 ) += btqr *float(mult);
02108 (*m_volume)( 2*ixn+1, iyn+1, 1 ) += btqi * float(mult);
02109 (*m_wptr)(ixn,iyn+1, 1) += float(mult);
02110 }
02111
02112 }
02113 checked_delete( padfft );
02114 return 0;
02115 }
02116
02117 int nn4Reconstructor::insert_padfft_slice( EMData* padfft, const Transform& t, int mult )
02118 {
02119 Assert( padfft != NULL );
02120
02121 for (int isym=0; isym < m_nsym; isym++) {
02122 Transform tsym = t.get_sym(m_symmetry, isym);
02123 m_volume->nn( m_wptr, padfft, tsym, mult);
02124 }
02125 return 0;
02126 }
02127
02128 #define tw(i,j,k) tw[ i-1 + (j-1+(k-1)*iy)*ix ]
02129
02130 void circumference( EMData* win )
02131 {
02132 float *tw = win->get_data();
02133
02134 int ix = win->get_xsize();
02135 int iy = win->get_ysize();
02136 int iz = win->get_zsize();
02137 int L2 = (ix/2)*(ix/2);
02138 int L2P = (ix/2-1)*(ix/2-1);
02139
02140 int IP = ix/2+1;
02141 int JP = iy/2+1;
02142 int KP = iz/2+1;
02143
02144 float TNR = 0.0f;
02145 size_t m = 0;
02146 for (int k = 1; k <= iz; ++k) {
02147 for (int j = 1; j <= iy; ++j) {
02148 for (int i = 1; i <= ix; ++i) {
02149 size_t LR = (k-KP)*(k-KP)+(j-JP)*(j-JP)+(i-IP)*(i-IP);
02150 if (LR<=(size_t)L2) {
02151 if(LR >= (size_t)L2P && LR <= (size_t)L2) {
02152 TNR += tw(i,j,k);
02153 ++m;
02154 }
02155 }
02156 }
02157 }
02158 }
02159
02160 TNR /=float(m);
02161 for (int k = 1; k <= iz; ++k) {
02162 for (int j = 1; j <= iy; ++j) {
02163 for (int i = 1; i <= ix; ++i) {
02164 size_t LR = (k-KP)*(k-KP)+(j-JP)*(j-JP)+(i-IP)*(i-IP);
02165 if (LR<=(size_t)L2) tw(i,j,k) -= TNR; else tw(i,j,k) = 0.0f;
02166 }
02167 }
02168 }
02169
02170 }
02171
02172
02173 EMData* nn4Reconstructor::finish(bool doift)
02174 {
02175 if( m_ndim==3 ) {
02176 m_volume->symplane0(m_wptr);
02177 } else {
02178 for( int i=1; i <= m_vnyp; ++i ) {
02179
02180 if( (*m_wptr)(0, i, 1)==0.0 ) {
02181 int j = m_vnyp + 1 - i;
02182 (*m_wptr)(0, i, 1) = (*m_wptr)(0, j, 1);
02183 (*m_volume)(0, i, 1) = (*m_volume)(0, j, 1);
02184 (*m_volume)(1, i, 1) = (*m_volume)(1, j, 1);
02185 }
02186 }
02187 }
02188
02189
02190 int box = 7;
02191 int kc = (box-1)/2;
02192 vector< float > pow_a( m_ndim*kc+1, 1.0 );
02193 for( unsigned int i=1; i < pow_a.size(); ++i ) pow_a[i] = pow_a[i-1] * exp(m_wghta);
02194 pow_a.back()=0.0;
02195
02196 float alpha = 0.0;
02197 if( m_ndim==3) {
02198 int vol = box*box*box;
02199 float max = max3d( kc, pow_a );
02200 alpha = ( 1.0f - 1.0f/(float)vol ) / max;
02201 } else {
02202 int ara = box*box;
02203 float max = max2d( kc, pow_a );
02204 alpha = ( 1.0f - 1.0f/(float)ara ) / max;
02205 }
02206
02207 int ix,iy,iz;
02208 for (iz = 1; iz <= m_vnzp; iz++) {
02209 for (iy = 1; iy <= m_vnyp; iy++) {
02210 for (ix = 0; ix <= m_vnxc; ix++) {
02211 if ( (*m_wptr)(ix,iy,iz) > 0) {
02212 float tmp;
02213 tmp = (-2*((ix+iy+iz)%2)+1)/((*m_wptr)(ix,iy,iz)+m_osnr);
02214
02215 if( m_weighting == ESTIMATE ) {
02216 int cx = ix;
02217 int cy = (iy<=m_vnyc) ? iy - 1 : iy - 1 - m_vnyp;
02218 int cz = (iz<=m_vnzc) ? iz - 1 : iz - 1 - m_vnzp;
02219 float sum = 0.0;
02220 for( int ii = -kc; ii <= kc; ++ii ) {
02221 int nbrcx = cx + ii;
02222 if( nbrcx >= m_vnxc ) continue;
02223 for( int jj= -kc; jj <= kc; ++jj ) {
02224 int nbrcy = cy + jj;
02225 if( nbrcy <= -m_vnyc || nbrcy >= m_vnyc ) continue;
02226
02227 int kcz = (m_ndim==3) ? kc : 0;
02228 for( int kk = -kcz; kk <= kcz; ++kk ) {
02229 int nbrcz = cz + kk;
02230 if( nbrcz <= -m_vnyc || nbrcz >= m_vnyc ) continue;
02231 if( nbrcx < 0 ) {
02232 nbrcx = -nbrcx;
02233 nbrcy = -nbrcy;
02234 nbrcz = -nbrcz;
02235 }
02236 int nbrix = nbrcx;
02237 int nbriy = nbrcy >= 0 ? nbrcy + 1 : nbrcy + 1 + m_vnyp;
02238 int nbriz = nbrcz >= 0 ? nbrcz + 1 : nbrcz + 1 + m_vnzp;
02239 if( (*m_wptr)( nbrix, nbriy, nbriz ) == 0 ) {
02240 int c = m_ndim*kc+1 - std::abs(ii) - std::abs(jj) - std::abs(kk);
02241 sum = sum + pow_a[c];
02242 }
02243 }
02244 }
02245 }
02246 float wght = 1.0f / ( 1.0f - alpha * sum );
02247 tmp = tmp * wght;
02248 }
02249 (*m_volume)(2*ix,iy,iz) *= tmp;
02250 (*m_volume)(2*ix+1,iy,iz) *= tmp;
02251 }
02252 }
02253 }
02254 }
02255
02256
02257
02258
02259 m_volume->do_ift_inplace();
02260
02261
02262 m_volume->depad();
02263 circumference( m_volume );
02264 m_volume->set_array_offsets( 0, 0, 0 );
02265
02266 m_result = m_volume->copy();
02267 return m_result;
02268 }
02269 #undef tw
02270
02271
02272
02273
02274
02275 nnSSNR_Reconstructor::nnSSNR_Reconstructor()
02276 {
02277 m_volume = NULL;
02278 m_wptr = NULL;
02279 m_wptr2 = NULL;
02280 m_result = NULL;
02281 }
02282
02283 nnSSNR_Reconstructor::nnSSNR_Reconstructor( const string& symmetry, int size, int npad)
02284 {
02285 m_volume = NULL;
02286 m_wptr = NULL;
02287 m_wptr2 = NULL;
02288 m_result = NULL;
02289
02290 setup( symmetry, size, npad );
02291 }
02292
02293 nnSSNR_Reconstructor::~nnSSNR_Reconstructor()
02294 {
02295 if( m_delete_volume ) checked_delete(m_volume);
02296
02297 if( m_delete_weight ) checked_delete( m_wptr );
02298
02299 if( m_delete_weight2 ) checked_delete( m_wptr2 );
02300
02301 checked_delete( m_result );
02302 }
02303
02304 void nnSSNR_Reconstructor::setup()
02305 {
02306 int size = params["size"];
02307 int npad = params["npad"];
02308
02309 string symmetry;
02310 if( params.has_key("symmetry") ) symmetry = params["symmetry"].to_str();
02311 else symmetry = "c1";
02312
02313 setup( symmetry, size, npad );
02314 }
02315
02316 void nnSSNR_Reconstructor::setup( const string& symmetry, int size, int npad )
02317 {
02318
02319 m_weighting = ESTIMATE;
02320 m_wghta = 0.2f;
02321 m_wghtb = 0.004f;
02322
02323 m_symmetry = symmetry;
02324 m_npad = npad;
02325 m_nsym = Transform::get_nsym(m_symmetry);
02326
02327 m_vnx = size;
02328 m_vny = size;
02329 m_vnz = size;
02330
02331 m_vnxp = size*npad;
02332 m_vnyp = size*npad;
02333 m_vnzp = size*npad;
02334
02335 m_vnxc = m_vnxp/2;
02336 m_vnyc = m_vnyp/2;
02337 m_vnzc = m_vnzp/2;
02338
02339 buildFFTVolume();
02340 buildNormVolume();
02341 buildNorm2Volume();
02342
02343 }
02344
02345 void nnSSNR_Reconstructor::buildFFTVolume() {
02346
02347 if( params.has_key("fftvol") ) {
02348 m_volume = params["fftvol"];
02349 m_delete_volume = false;
02350 } else {
02351 m_volume = new EMData();
02352 m_delete_volume = true;
02353 }
02354 m_volume->set_size(m_vnxp+ 2 - m_vnxp%2,m_vnyp,m_vnzp);
02355 m_volume->to_zero();
02356 if ( m_vnxp % 2 == 0 ) m_volume->set_fftodd(0);
02357 else m_volume->set_fftodd(1);
02358
02359 m_volume->set_nxc(m_vnxc);
02360 m_volume->set_complex(true);
02361 m_volume->set_ri(true);
02362 m_volume->set_fftpad(true);
02363 m_volume->set_attr("npad", m_npad);
02364 m_volume->set_array_offsets(0,1,1);
02365 }
02366
02367 void nnSSNR_Reconstructor::buildNormVolume() {
02368 if( params.has_key("weight") ) {
02369 m_wptr = params["weight"];
02370 m_delete_weight = false;
02371 } else {
02372 m_wptr = new EMData();
02373 m_delete_weight = true;
02374 }
02375
02376 m_wptr->set_size(m_vnxc+1,m_vnyp,m_vnzp);
02377 m_wptr->to_zero();
02378
02379 m_wptr->set_array_offsets(0,1,1);
02380 }
02381
02382 void nnSSNR_Reconstructor::buildNorm2Volume() {
02383
02384 if( params.has_key("weight2") ) {
02385 m_wptr2 = params["weight2"];
02386 m_delete_weight2 = false;
02387 } else {
02388 m_wptr2 = new EMData();
02389 m_delete_weight2 = true;
02390 }
02391 m_wptr2->set_size(m_vnxc+1,m_vnyp,m_vnzp);
02392 m_wptr2->to_zero();
02393 m_wptr2->set_array_offsets(0,1,1);
02394 }
02395
02396
02397 int nnSSNR_Reconstructor::insert_slice(const EMData* const slice, const Transform& t, const float weight) {
02398
02399 if (!slice) {
02400 LOGERR("try to insert NULL slice");
02401 return 1;
02402 }
02403
02404 int padffted=slice->get_attr_default( "padffted", 0 );
02405
02406 if ( padffted==0 && (slice->get_xsize()!=slice->get_ysize() || slice->get_xsize()!=m_vnx) ) {
02407
02408 LOGERR("Tried to insert a slice that has wrong size.");
02409 return 1;
02410 }
02411
02412 EMData* padfft = NULL;
02413
02414 if( padffted != 0 ) padfft = new EMData(*slice);
02415 else padfft = padfft_slice( slice, t, m_npad );
02416
02417 int mult = slice->get_attr_default("mult", 1);
02418
02419 Assert( mult > 0 );
02420 insert_padfft_slice( padfft, t, mult );
02421
02422 checked_delete( padfft );
02423 return 0;
02424 }
02425
02426 int nnSSNR_Reconstructor::insert_padfft_slice( EMData* padfft, const Transform& t, int mult )
02427 {
02428 Assert( padfft != NULL );
02429
02430 for (int isym=0; isym < m_nsym; isym++) {
02431 Transform tsym = t.get_sym(m_symmetry, isym);
02432 m_volume->nn_SSNR( m_wptr, m_wptr2, padfft, tsym, mult);
02433 }
02434 return 0;
02435 }
02436
02437
02438 #define tw(i,j,k) tw[ i-1 + (j-1+(k-1)*iy)*ix ]
02439 EMData* nnSSNR_Reconstructor::finish(bool doift)
02440 {
02441
02442
02443
02444
02445
02446
02447 int kz, ky;
02448
02449 int box = 7;
02450 int kc = (box-1)/2;
02451 float alpha = 0.0;
02452 float argx, argy, argz;
02453 vector< float > pow_a( 3*kc+1, 1.0 );
02454 float w = params["w"];
02455 EMData* SSNR = params["SSNR"];
02456
02457
02458
02459
02460 EMData* vol_ssnr = new EMData();
02461 vol_ssnr->set_size(m_vnxp+ 2 - m_vnxp%2, m_vnyp ,m_vnzp);
02462 vol_ssnr->to_zero();
02463 if ( m_vnxp % 2 == 0 ) vol_ssnr->set_fftodd(0);
02464 else vol_ssnr->set_fftodd(1);
02465 vol_ssnr->set_nxc(m_vnxc);
02466 vol_ssnr->set_complex(true);
02467 vol_ssnr->set_ri(true);
02468 vol_ssnr->set_fftpad(false);
02469
02470
02471 float dx2 = 1.0f/float(m_vnxc)/float(m_vnxc);
02472 float dy2 = 1.0f/float(m_vnyc)/float(m_vnyc);
02473 #ifdef _WIN32
02474 float dz2 = 1.0f/_cpp_max(float(m_vnzc),1.0f)/_cpp_max(float(m_vnzc),1.0f);
02475 int inc = Util::round(float(_cpp_max(_cpp_max(m_vnxc,m_vnyc),m_vnzc))/w);
02476 #else
02477 float dz2 = 1.0f/std::max(float(m_vnzc),1.0f)/std::max(float(m_vnzc),1.0f);
02478 int inc = Util::round(float(std::max(std::max(m_vnxc,m_vnyc),m_vnzc))/w);
02479 #endif //_WIN32
02480 SSNR->set_size(inc+1,4,1);
02481
02482 float *nom = new float[inc+1];
02483 float *denom = new float[inc+1];
02484 int *nn = new int[inc+1];
02485 int *ka = new int[inc+1];
02486 float wght = 1.0f;
02487 for (int i = 0; i <= inc; i++) {
02488 nom[i] = 0.0f;
02489 denom[i] = 0.0f;
02490 nn[i] = 0;
02491 ka[i] = 0;
02492 }
02493
02494 m_volume->symplane1(m_wptr, m_wptr2);
02495
02496 if ( m_weighting == ESTIMATE ) {
02497 int vol = box*box*box;
02498 for( unsigned int i=1; i < pow_a.size(); ++i ) pow_a[i] = pow_a[i-1] * exp(m_wghta);
02499 pow_a[3*kc] = 0.0;
02500 float max = max3d( kc, pow_a );
02501 alpha = ( 1.0f - 1.0f/(float)vol ) / max;
02502 }
02503
02504 for (int iz = 1; iz <= m_vnzp; iz++) {
02505 if ( iz-1 > m_vnzc ) kz = iz-1-m_vnzp; else kz = iz-1;
02506 argz = float(kz*kz)*dz2;
02507 for (int iy = 1; iy <= m_vnyp; iy++) {
02508 if ( iy-1 > m_vnyc ) ky = iy-1-m_vnyp; else ky = iy-1;
02509 argy = argz + float(ky*ky)*dy2;
02510 for (int ix = 0; ix <= m_vnxc; ix++) {
02511 float Kn = (*m_wptr)(ix,iy,iz);
02512 argx = std::sqrt(argy + float(ix*ix)*dx2);
02513 int r = Util::round(float(inc)*argx);
02514 if ( r >= 0 && Kn > 4.5f ) {
02515 if ( m_weighting == ESTIMATE ) {
02516 int cx = ix;
02517 int cy = (iy<=m_vnyc) ? iy - 1 : iy - 1 - m_vnyp;
02518 int cz = (iz<=m_vnzc) ? iz - 1 : iz - 1 - m_vnzp;
02519
02520 float sum = 0.0;
02521 for( int ii = -kc; ii <= kc; ++ii ) {
02522 int nbrcx = cx + ii;
02523 if( nbrcx >= m_vnxc ) continue;
02524 for ( int jj= -kc; jj <= kc; ++jj ) {
02525 int nbrcy = cy + jj;
02526 if( nbrcy <= -m_vnyc || nbrcy >= m_vnyc ) continue;
02527 for( int kk = -kc; kk <= kc; ++kk ) {
02528 int nbrcz = cz + jj;
02529 if ( nbrcz <= -m_vnyc || nbrcz >= m_vnyc ) continue;
02530 if( nbrcx < 0 ) {
02531 nbrcx = -nbrcx;
02532 nbrcy = -nbrcy;
02533 nbrcz = -nbrcz;
02534 }
02535 int nbrix = nbrcx;
02536 int nbriy = nbrcy >= 0 ? nbrcy + 1 : nbrcy + 1 + m_vnyp;
02537 int nbriz = nbrcz >= 0 ? nbrcz + 1 : nbrcz + 1 + m_vnzp;
02538 if( (*m_wptr)( nbrix, nbriy, nbriz ) == 0 ) {
02539 int c = 3*kc+1 - std::abs(ii) - std::abs(jj) - std::abs(kk);
02540 sum = sum + pow_a[c];
02541 }
02542 }
02543 }
02544 }
02545 wght = 1.0f / ( 1.0f - alpha * sum );
02546 }
02547 float nominator = std::norm(m_volume->cmplx(ix,iy,iz))/Kn;
02548 float denominator = ((*m_wptr2)(ix,iy,iz)-nominator)/(Kn-1.0f);
02549
02550 if( (ix>0 || (kz>=0 && (ky>=0 || kz!=0)))) {
02551 if ( r <= inc ) {
02552 nom[r] += nominator*wght;
02553 denom[r] += denominator/Kn*wght;
02554 nn[r] += 2;
02555 ka[r] += int(Kn);
02556 }
02557
02558
02559
02560
02561
02562
02563
02564
02565
02566
02567
02568
02569
02570
02571
02572
02573 }
02574 (*vol_ssnr)(2*ix, iy-1, iz-1) = denominator*wght;
02575
02576
02577 }
02578 }
02579 }
02580 }
02581
02582 for (int i = 0; i <= inc; i++) {
02583 (*SSNR)(i,0,0) = nom[i];
02584 (*SSNR)(i,1,0) = denom[i];
02585 (*SSNR)(i,2,0) = static_cast<float>(nn[i]);
02586 (*SSNR)(i,3,0) = static_cast<float>(ka[i]);
02587 }
02588 vol_ssnr->update();
02589 return vol_ssnr;
02590 }
02591 #undef tw
02592
02593
02594
02595
02596
02597
02598
02599 nn4_ctfReconstructor::nn4_ctfReconstructor()
02600 {
02601 m_volume = NULL;
02602 m_wptr = NULL;
02603 m_result = NULL;
02604 }
02605
02606 nn4_ctfReconstructor::nn4_ctfReconstructor( const string& symmetry, int size, int npad, float snr, int sign )
02607 {
02608 setup( symmetry, size, npad, snr, sign );
02609 }
02610
02611 nn4_ctfReconstructor::~nn4_ctfReconstructor()
02612 {
02613 if( m_delete_volume ) checked_delete(m_volume);
02614
02615 if( m_delete_weight ) checked_delete( m_wptr );
02616
02617 checked_delete( m_result );
02618 }
02619
02620 void nn4_ctfReconstructor::setup()
02621 {
02622 if( ! params.has_key("size") ) throw std::logic_error("Error: image size is not given");
02623
02624 int size = params["size"];
02625 int npad = params.has_key("npad") ? int(params["npad"]) : 4;
02626
02627 int sign = 1;
02628 string symmetry = params.has_key("symmetry")? params["symmetry"].to_str() : "c1";
02629
02630 float snr = params["snr"];
02631
02632 m_varsnr = params.has_key("varsnr") ? int(params["varsnr"]) : 0;
02633 setup( symmetry, size, npad, snr, sign );
02634
02635 }
02636
02637 void nn4_ctfReconstructor::setup( const string& symmetry, int size, int npad, float snr, int sign )
02638 {
02639 m_weighting = ESTIMATE;
02640 if( params.has_key("weighting") ) {
02641 int tmp = int( params["weighting"] );
02642 if( tmp==0 ) m_weighting = NONE;
02643 }
02644
02645
02646
02647 m_wghta = 0.2f;
02648 m_wghtb = 0.004f;
02649
02650 m_symmetry = symmetry;
02651 m_npad = npad;
02652 m_sign = sign;
02653 m_nsym = Transform::get_nsym(m_symmetry);
02654
02655 m_snr = snr;
02656
02657 m_vnx = size;
02658 m_vny = size;
02659 m_vnz = size;
02660
02661 m_vnxp = size*npad;
02662 m_vnyp = size*npad;
02663 m_vnzp = size*npad;
02664
02665 m_vnxc = m_vnxp/2;
02666 m_vnyc = m_vnyp/2;
02667 m_vnzc = m_vnzp/2;
02668
02669 buildFFTVolume();
02670 buildNormVolume();
02671 }
02672
02673 void nn4_ctfReconstructor::buildFFTVolume() {
02674 int offset = 2 - m_vnxp%2;
02675 if( params.has_key("fftvol") ) {
02676 m_volume = params["fftvol"];
02677 m_delete_volume = false;
02678 } else {
02679 m_volume = new EMData();
02680 m_delete_volume = true;
02681 }
02682
02683 if( m_volume->get_xsize() != m_vnxp+offset &&
02684 m_volume->get_ysize() != m_vnyp &&
02685 m_volume->get_zsize() != m_vnzp )
02686 {
02687 m_volume->set_size(m_vnxp+offset,m_vnyp,m_vnzp);
02688 m_volume->to_zero();
02689 }
02690
02691 m_volume->set_nxc(m_vnxp/2);
02692 m_volume->set_complex(true);
02693 m_volume->set_ri(true);
02694 m_volume->set_fftpad(true);
02695 m_volume->set_attr("npad", m_npad);
02696 m_volume->set_array_offsets(0,1,1);
02697 }
02698
02699 void nn4_ctfReconstructor::buildNormVolume()
02700 {
02701 if( params.has_key("weight") )
02702 {
02703 m_wptr = params["weight"];
02704 m_delete_weight = false;
02705 } else {
02706 m_wptr = new EMData();
02707 m_delete_weight = true;
02708 }
02709
02710 if( m_wptr->get_xsize() != m_vnxc+1 &&
02711 m_wptr->get_ysize() != m_vnyp &&
02712 m_wptr->get_zsize() != m_vnzp )
02713 {
02714 m_wptr->set_size(m_vnxc+1,m_vnyp,m_vnzp);
02715 m_wptr->to_zero();
02716 }
02717
02718 m_wptr->set_array_offsets(0,1,1);
02719
02720 }
02721
02722 int nn4_ctfReconstructor::insert_slice(const EMData* const slice, const Transform& t, const float weight)
02723 {
02724
02725 if (!slice) {
02726 LOGERR("try to insert NULL slice");
02727 return 1;
02728 }
02729
02730 int buffed = slice->get_attr_default( "buffed", 0 );
02731 if( buffed > 0 )
02732 {
02733 int mult = slice->get_attr_default( "mult", 1 );
02734 insert_buffed_slice( slice, mult );
02735 return 0;
02736 }
02737
02738 int padffted= slice->get_attr_default("padffted", 0);
02739 if( padffted==0 && (slice->get_xsize()!=slice->get_ysize() || slice->get_xsize()!=m_vnx) )
02740 {
02741
02742 LOGERR("Tried to insert a slice that is the wrong size.");
02743 return 1;
02744 }
02745
02746 EMData* padfft = NULL;
02747
02748 if( padffted != 0 ) padfft = new EMData(*slice);
02749 else padfft = padfft_slice( slice, t, m_npad );
02750
02751 int mult= slice->get_attr_default("mult", 1);
02752
02753 Assert( mult > 0 );
02754 insert_padfft_slice( padfft, t, mult );
02755
02756 checked_delete( padfft );
02757
02758 return 0;
02759 }
02760
02761 int nn4_ctfReconstructor::insert_buffed_slice( const EMData* buffed, int mult )
02762 {
02763 const float* bufdata = buffed->get_data();
02764 float* cdata = m_volume->get_data();
02765 float* wdata = m_wptr->get_data();
02766
02767 int npoint = buffed->get_xsize()/4;
02768 for( int i=0; i < npoint; ++i ) {
02769
02770 int pos2 = int( bufdata[4*i] );
02771 int pos1 = pos2 * 2;
02772 cdata[pos1 ] += bufdata[4*i+1]*mult;
02773 cdata[pos1+1] += bufdata[4*i+2]*mult;
02774 wdata[pos2 ] += bufdata[4*i+3]*mult;
02775
02776
02777
02778
02779
02780 }
02781 return 0;
02782 }
02783
02784 int nn4_ctfReconstructor::insert_padfft_slice( EMData* padfft, const Transform& t, int mult )
02785 {
02786 Assert( padfft != NULL );
02787 float tmp = padfft->get_attr("ctf_applied");
02788 int ctf_applied = (int) tmp;
02789
02790 for( int isym=0; isym < m_nsym; isym++) {
02791 Transform tsym = t.get_sym( m_symmetry, isym );
02792
02793 if(ctf_applied) m_volume->nn_ctf_applied(m_wptr, padfft, tsym, mult);
02794 else m_volume->nn_ctf(m_wptr, padfft, tsym, mult);
02795 }
02796
02797 return 0;
02798
02799 }
02800
02801 #define tw(i,j,k) tw[ i-1 + (j-1+(k-1)*iy)*ix ]
02802 EMData* nn4_ctfReconstructor::finish(bool doift)
02803 {
02804 m_volume->set_array_offsets(0, 1, 1);
02805 m_wptr->set_array_offsets(0, 1, 1);
02806 m_volume->symplane0_ctf(m_wptr);
02807
02808 int box = 7;
02809 int vol = box*box*box;
02810 int kc = (box-1)/2;
02811 vector< float > pow_a( 3*kc+1, 1.0 );
02812 for( unsigned int i=1; i < pow_a.size(); ++i ) pow_a[i] = pow_a[i-1] * exp(m_wghta);
02813 pow_a[3*kc]=0.0;
02814
02815
02816 float max = max3d( kc, pow_a );
02817 float alpha = ( 1.0f - 1.0f/(float)vol ) / max;
02818 float osnr = 1.0f/m_snr;
02819
02820
02821 int ix,iy,iz;
02822 for (iz = 1; iz <= m_vnzp; iz++) {
02823 for (iy = 1; iy <= m_vnyp; iy++) {
02824 for (ix = 0; ix <= m_vnxc; ix++) {
02825 if ( (*m_wptr)(ix,iy,iz) > 0.0f) {
02826 int iyp = (iy<=m_vnyc) ? iy - 1 : iy-m_vnyp-1;
02827 int izp = (iz<=m_vnzc) ? iz - 1 : iz-m_vnzp-1;
02828 float tmp=0.0;
02829 if( m_varsnr )
02830 {
02831 float freq = sqrt( (float)(ix*ix+iyp*iyp+izp*izp) );
02832 tmp = (-2*((ix+iy+iz)%2)+1)/((*m_wptr)(ix,iy,iz)+freq*osnr)*m_sign;
02833 }
02834 else
02835 {
02836 tmp = (-2*((ix+iy+iz)%2)+1)/((*m_wptr)(ix,iy,iz)+osnr)*m_sign;
02837 }
02838
02839 if( m_weighting == ESTIMATE ) {
02840 int cx = ix;
02841 int cy = (iy<=m_vnyc) ? iy - 1 : iy - 1 - m_vnyp;
02842 int cz = (iz<=m_vnzc) ? iz - 1 : iz - 1 - m_vnzp;
02843 float sum = 0.0;
02844 for( int ii = -kc; ii <= kc; ++ii ) {
02845 int nbrcx = cx + ii;
02846 if( nbrcx >= m_vnxc ) continue;
02847 for( int jj= -kc; jj <= kc; ++jj ) {
02848 int nbrcy = cy + jj;
02849 if( nbrcy <= -m_vnyc || nbrcy >= m_vnyc ) continue;
02850 for( int kk = -kc; kk <= kc; ++kk ) {
02851 int nbrcz = cz + jj;
02852 if( nbrcz <= -m_vnyc || nbrcz >= m_vnyc ) continue;
02853 if( nbrcx < 0 ) {
02854 nbrcx = -nbrcx;
02855 nbrcy = -nbrcy;
02856 nbrcz = -nbrcz;
02857 }
02858
02859 int nbrix = nbrcx;
02860 int nbriy = nbrcy >= 0 ? nbrcy + 1 : nbrcy + 1 + m_vnyp;
02861 int nbriz = nbrcz >= 0 ? nbrcz + 1 : nbrcz + 1 + m_vnzp;
02862 if( (*m_wptr)( nbrix, nbriy, nbriz ) == 0.0 ) {
02863 int c = 3*kc+1 - std::abs(ii) - std::abs(jj) - std::abs(kk);
02864 sum = sum + pow_a[c];
02865
02866
02867 }
02868 }
02869 }
02870 }
02871 float wght = 1.0f / ( 1.0f - alpha * sum );
02872
02873
02874
02875
02876
02877
02878
02879
02880
02881 tmp = tmp * wght;
02882 }
02883 (*m_volume)(2*ix,iy,iz) *= tmp;
02884 (*m_volume)(2*ix+1,iy,iz) *= tmp;
02885 }
02886 }
02887 }
02888 }
02889
02890
02891 m_volume->do_ift_inplace();
02892 m_volume->depad();
02893 circumference( m_volume );
02894 m_volume->set_array_offsets( 0, 0, 0 );
02895 m_result = m_volume->copy();
02896 return m_result;
02897 }
02898 #undef tw
02899
02900
02901
02902
02903
02904
02905
02906
02907 nnSSNR_ctfReconstructor::nnSSNR_ctfReconstructor()
02908 {
02909 m_volume = NULL;
02910 m_wptr = NULL;
02911 m_wptr2 = NULL;
02912 m_wptr3 = NULL;
02913 m_result = NULL;
02914 }
02915
02916 nnSSNR_ctfReconstructor::nnSSNR_ctfReconstructor( const string& symmetry, int size, int npad, float snr, int sign)
02917 {
02918 m_volume = NULL;
02919 m_wptr = NULL;
02920 m_wptr2 = NULL;
02921 m_wptr3 = NULL;
02922 m_result = NULL;
02923
02924 setup( symmetry, size, npad, snr, sign );
02925 }
02926
02927 nnSSNR_ctfReconstructor::~nnSSNR_ctfReconstructor()
02928 {
02929
02930 if( m_delete_volume )
02931 checked_delete(m_volume);
02932 if( m_delete_weight )
02933 checked_delete( m_wptr );
02934 if ( m_delete_weight2 )
02935 checked_delete( m_wptr2 );
02936 if ( m_delete_weight3 )
02937 checked_delete( m_wptr3 );
02938 checked_delete( m_result );
02939 }
02940
02941 void nnSSNR_ctfReconstructor::setup()
02942 {
02943 int size = params["size"];
02944 int npad = params["npad"];
02945 int sign = params["sign"];
02946 float snr = params["snr"];
02947 string symmetry;
02948 if( params.has_key("symmetry") ) {
02949 symmetry = params["symmetry"].to_str();
02950 } else {
02951 symmetry = "c1";
02952 }
02953 setup( symmetry, size, npad, snr, sign );
02954 }
02955 void nnSSNR_ctfReconstructor::setup( const string& symmetry, int size, int npad, float snr, int sign )
02956 {
02957
02958 m_weighting = ESTIMATE;
02959 m_wghta = 0.2f;
02960 m_wghtb = 0.004f;
02961 wiener = 1;
02962
02963 m_symmetry = symmetry;
02964 m_npad = npad;
02965 m_nsym = Transform::get_nsym(m_symmetry);
02966
02967 m_sign = sign;
02968 m_snr = snr;
02969
02970 m_vnx = size;
02971 m_vny = size;
02972 m_vnz = size;
02973
02974 m_vnxp = size*npad;
02975 m_vnyp = size*npad;
02976 m_vnzp = size*npad;
02977
02978 m_vnxc = m_vnxp/2;
02979 m_vnyc = m_vnyp/2;
02980 m_vnzc = m_vnzp/2;
02981
02982 buildFFTVolume();
02983 buildNormVolume();
02984 buildNorm2Volume();
02985 buildNorm3Volume();
02986 }
02987
02988 void nnSSNR_ctfReconstructor::buildFFTVolume() {
02989
02990 int offset = 2 - m_vnxp%2;
02991 if( params.has_key("fftvol") ) {
02992 m_volume = params["fftvol"];
02993 m_delete_volume = false;
02994 } else {
02995 m_volume = new EMData();
02996 m_delete_volume = true;
02997 }
02998
02999 m_volume->set_size(m_vnxp+offset,m_vnyp,m_vnzp);
03000 m_volume->to_zero();
03001
03002 if ( m_vnxp % 2 == 0 ) {
03003 m_volume->set_fftodd(0);
03004 } else {
03005 m_volume->set_fftodd(1);
03006 }
03007
03008 m_volume->set_nxc(m_vnxc);
03009 m_volume->set_complex(true);
03010 m_volume->set_ri(true);
03011 m_volume->set_fftpad(true);
03012 m_volume->set_attr("npad", m_npad);
03013 m_volume->set_array_offsets(0,1,1);
03014 }
03015
03016
03017 void nnSSNR_ctfReconstructor::buildNormVolume()
03018 {
03019 if( params.has_key("weight") ) {
03020 m_wptr = params["weight"];
03021 m_delete_weight = false;
03022 } else {
03023 m_wptr = new EMData();
03024 m_delete_weight = true;
03025 }
03026 m_wptr->set_size(m_vnxc+1,m_vnyp,m_vnzp);
03027 m_wptr->to_zero();
03028 m_wptr->set_array_offsets(0,1,1);
03029 }
03030
03031 void nnSSNR_ctfReconstructor::buildNorm2Volume() {
03032
03033 if( params.has_key("weight2") ) {
03034 m_wptr2 = params["weight2"];
03035 m_delete_weight2 = false;
03036 } else {
03037 m_wptr2 = new EMData();
03038 m_delete_weight2 = true;
03039 }
03040 m_wptr2->set_size(m_vnxc+1,m_vnyp,m_vnzp);
03041 m_wptr2->to_zero();
03042 m_wptr2->set_array_offsets(0,1,1);
03043 }
03044
03045 void nnSSNR_ctfReconstructor::buildNorm3Volume() {
03046
03047 if( params.has_key("weight3") ) {
03048 m_wptr3 = params["weight3"];
03049 m_delete_weight3 = false;
03050 } else {
03051 m_wptr3 = new EMData();
03052 m_delete_weight3 = true;
03053 }
03054 m_wptr3->set_size(m_vnxc+1,m_vnyp,m_vnzp);
03055 m_wptr3->to_zero();
03056 m_wptr3->set_array_offsets(0,1,1);
03057 }
03058
03059 int nnSSNR_ctfReconstructor::insert_slice(const EMData *const slice, const Transform& t, const float weight) {
03060
03061 if (!slice) {
03062 LOGERR("try to insert NULL slice");
03063 return 1;
03064 }
03065 int padffted= slice->get_attr_default("padffted", 0);
03066 if ( padffted==0 && (slice->get_xsize()!=slice->get_ysize() || slice->get_xsize()!=m_vnx) ) {
03067
03068 LOGERR("Tried to insert a slice that is the wrong size.");
03069 return 1;
03070 }
03071 EMData* padfft = NULL;
03072
03073 if( padffted != 0 ) {
03074 padfft = new EMData(*slice);
03075 } else {
03076 padfft = padfft_slice( slice, t, m_npad );
03077 }
03078
03079 int mult= slice->get_attr_default("mult", 1);
03080
03081 Assert( mult > 0 );
03082 insert_padfft_slice( padfft, t, mult );
03083
03084 checked_delete( padfft );
03085 return 0;
03086 }
03087 int nnSSNR_ctfReconstructor::insert_padfft_slice( EMData* padfft, const Transform& t, int mult )
03088 {
03089 Assert( padfft != NULL );
03090
03091
03092 for (int isym=0; isym < m_nsym; isym++) {
03093 Transform tsym = t.get_sym(m_symmetry, isym);
03094 m_volume->nn_SSNR_ctf(m_wptr, m_wptr2, m_wptr3, padfft, tsym, mult);
03095 }
03096
03097 return 0;
03098 }
03099
03100 #define tw(i,j,k) tw[ i-1 + (j-1+(k-1)*iy)*ix ]
03101 EMData* nnSSNR_ctfReconstructor::finish(bool doift)
03102 {
03103
03104
03105
03106
03107
03108
03109
03110
03111
03112
03113
03114
03115
03116
03117 int kz, ky;
03118 int box = 7;
03119 int kc = (box-1)/2;
03120 float alpha = 0.0;
03121 float argx, argy, argz;
03122 vector< float > pow_a( 3*kc+1, 1.0 );
03123 float w = params["w"];
03124 float dx2 = 1.0f/float(m_vnxc)/float(m_vnxc);
03125 float dy2 = 1.0f/float(m_vnyc)/float(m_vnyc);
03126 #ifdef _WIN32
03127 float dz2 = 1.0f/_cpp_max(float(m_vnzc),1.0f)/_cpp_max(float(m_vnzc),1.0f);
03128 int inc = Util::round(float(_cpp_max(_cpp_max(m_vnxc,m_vnyc),m_vnzc))/w);
03129 #else
03130 float dz2 = 1.0f/std::max(float(m_vnzc),1.0f)/std::max(float(m_vnzc),1.0f);
03131 int inc = Util::round(float(std::max(std::max(m_vnxc,m_vnyc),m_vnzc))/w);
03132 #endif //_WIN32
03133
03134 EMData* SSNR = params["SSNR"];
03135 SSNR->set_size(inc+1,4,1);
03136
03137
03138
03139
03140 EMData* vol_ssnr = new EMData();
03141 vol_ssnr->set_size(m_vnxp+ 2 - m_vnxp%2, m_vnyp ,m_vnzp);
03142 vol_ssnr->to_zero();
03143 if ( m_vnxp % 2 == 0 ) vol_ssnr->set_fftodd(0);
03144 else vol_ssnr->set_fftodd(1);
03145 vol_ssnr->set_nxc(m_vnxc);
03146 vol_ssnr->set_complex(true);
03147 vol_ssnr->set_ri(true);
03148 vol_ssnr->set_fftpad(false);
03149
03150 float *nom = new float[inc+1];
03151 float *denom = new float[inc+1];
03152 int *ka = new int[inc+1];
03153 int *nn = new int[inc+1];
03154 float wght = 1.f;
03155 for (int i = 0; i <= inc; i++) {
03156 nom[i] = 0.0f;
03157 denom[i] = 0.0f;
03158 nn[i] = 0;
03159 ka[i] = 0;
03160 }
03161 m_volume->symplane2(m_wptr, m_wptr2, m_wptr3);
03162 if ( m_weighting == ESTIMATE ) {
03163 int vol = box*box*box;
03164 for( unsigned int i=1; i < pow_a.size(); ++i ) pow_a[i] = pow_a[i-1] * exp(m_wghta);
03165 pow_a[3*kc] = 0.0;
03166 float max = max3d( kc, pow_a );
03167 alpha = ( 1.0f - 1.0f/(float)vol ) / max;
03168 }
03169 for (int iz = 1; iz <= m_vnzp; iz++) {
03170 if ( iz-1 > m_vnzc ) kz = iz-1-m_vnzp; else kz = iz-1;
03171 argz = float(kz*kz)*dz2;
03172 for (int iy = 1; iy <= m_vnyp; iy++) {
03173 if ( iy-1 > m_vnyc ) ky = iy-1-m_vnyp; else ky = iy-1;
03174 argy = argz + float(ky*ky)*dy2;
03175 for (int ix = 0; ix <= m_vnxc; ix++) {
03176 float Kn = (*m_wptr3)(ix,iy,iz);
03177 argx = std::sqrt(argy + float(ix*ix)*dx2);
03178 int r = Util::round(float(inc)*argx);
03179 if ( r >= 0 && Kn > 4.5f ) {
03180 if ( m_weighting == ESTIMATE ) {
03181 int cx = ix;
03182 int cy = (iy<=m_vnyc) ? iy - 1 : iy - 1 - m_vnyp;
03183 int cz = (iz<=m_vnzc) ? iz - 1 : iz - 1 - m_vnzp;
03184 float sum = 0.0;
03185 for( int ii = -kc; ii <= kc; ++ii ) {
03186 int nbrcx = cx + ii;
03187 if( nbrcx >= m_vnxc ) continue;
03188 for ( int jj= -kc; jj <= kc; ++jj ) {
03189 int nbrcy = cy + jj;
03190 if( nbrcy <= -m_vnyc || nbrcy >= m_vnyc ) continue;
03191 for( int kk = -kc; kk <= kc; ++kk ) {
03192 int nbrcz = cz + jj;
03193 if ( nbrcz <= -m_vnyc || nbrcz >= m_vnyc ) continue;
03194 if( nbrcx < 0 ) {
03195 nbrcx = -nbrcx;
03196 nbrcy = -nbrcy;
03197 nbrcz = -nbrcz;
03198 }
03199 int nbrix = nbrcx;
03200 int nbriy = nbrcy >= 0 ? nbrcy + 1 : nbrcy + 1 + m_vnyp;
03201 int nbriz = nbrcz >= 0 ? nbrcz + 1 : nbrcz + 1 + m_vnzp;
03202 if( (*m_wptr)( nbrix, nbriy, nbriz ) == 0 ) {
03203 int c = 3*kc+1 - std::abs(ii) - std::abs(jj) - std::abs(kk);
03204 sum = sum + pow_a[c];
03205 }
03206 }
03207 }
03208 }
03209
03210 wght = 1.0f / ( 1.0f - alpha * sum );
03211 }
03212 float nominator = std::norm(m_volume->cmplx(ix,iy,iz))/(*m_wptr)(ix,iy,iz);
03213 float denominator = ((*m_wptr2)(ix,iy,iz)-std::norm(m_volume->cmplx(ix,iy,iz))/(*m_wptr)(ix,iy,iz))/(Kn-1.0f);
03214
03215 if( (ix>0 || (kz>=0 && (ky>=0 || kz!=0)))) {
03216 if ( r <= inc ) {
03217 nom[r] += nominator*wght;
03218 denom[r] += denominator/(*m_wptr)(ix,iy,iz)*wght;
03219 nn[r] += 2;
03220 ka[r] += int(Kn);
03221 }
03222
03223
03224
03225
03226
03227
03228
03229
03230
03231
03232
03233
03234
03235
03236
03237 }
03238 (*vol_ssnr)(2*ix, iy-1, iz-1) = denominator*wght;
03239 }
03240 }
03241 }
03242 }
03243 for (int i = 0; i <= inc; i++) {
03244 (*SSNR)(i,0,0) = nom[i];
03245 (*SSNR)(i,1,0) = denom[i];
03246 (*SSNR)(i,2,0) = static_cast<float>(nn[i]);
03247 (*SSNR)(i,3,0) = static_cast<float>(ka[i]);
03248 }
03249 vol_ssnr->update();
03250 return vol_ssnr;
03251
03252 }
03253 #undef tw
03254
03255
03256
03257
03258 void EMAN::dump_reconstructors()
03259 {
03260 dump_factory < Reconstructor > ();
03261 }
03262
03263 map<string, vector<string> > EMAN::dump_reconstructors_list()
03264 {
03265 return dump_factory_list < Reconstructor > ();
03266 }
03267
03268
03269
03270
03271
03272 using std::ofstream;
03273 using std::ifstream;
03274
03275
03276 newfile_store::newfile_store( const string& filename, int npad, bool ctf )
03277 : m_bin_file( filename + ".bin" ),
03278 m_txt_file( filename + ".txt" )
03279 {
03280 m_npad = npad;
03281 m_ctf = ctf;
03282 }
03283
03284 newfile_store::~newfile_store( )
03285 {
03286 }
03287
03288 void newfile_store::add_image( EMData* emdata, const Transform& tf )
03289 {
03290 if( m_bin_of == NULL )
03291 {
03292 m_bin_of = shared_ptr<ofstream>( new ofstream(m_bin_file.c_str(), std::ios::out|std::ios::binary) );
03293 m_txt_of = shared_ptr<ofstream>( new ofstream(m_txt_file.c_str()) );
03294 }
03295
03296
03297 EMData* padfft = padfft_slice( emdata, tf, m_npad );
03298
03299 int nx = padfft->get_xsize();
03300 int ny = padfft->get_ysize();
03301 int n2 = ny / 2;
03302 int n = ny;
03303
03304 float voltage=0.0f, pixel=0.0f, Cs=0.0f, ampcont=0.0f, bfactor=0.0f, defocus=0.0f;
03305
03306 if( m_ctf )
03307 {
03308 Ctf* ctf = emdata->get_attr( "ctf" );
03309 Dict params = ctf->to_dict();
03310 voltage = params["voltage"];
03311 pixel = params["apix"];
03312 Cs = params["cs"];
03313 ampcont = params["ampcont"];
03314 bfactor = params["bfactor"];
03315 defocus = params["defocus"];
03316 if(ctf) {delete ctf; ctf=0;}
03317 }
03318
03319 vector<point_t> points;
03320 for( int j=-ny/2+1; j <= ny/2; j++ )
03321 {
03322 int jp = (j>=0) ? j+1 : ny+j+1;
03323 for( int i=0; i <= n2; ++i )
03324 {
03325 int r2 = i*i + j*j;
03326 if( (r2<ny*ny/4) && !( (i==0) && (j<0) ) )
03327 {
03328 float ctf;
03329 if( m_ctf )
03330 {
03331 float ak = std::sqrt( r2/float(ny*ny) )/pixel;
03332 ctf = Util::tf( defocus, ak, voltage, Cs, ampcont, bfactor, 1);
03333 }
03334 else
03335 {
03336 ctf = 1.0;
03337 }
03338
03339 float xnew = i*tf[0][0] + j*tf[1][0];
03340 float ynew = i*tf[0][1] + j*tf[1][1];
03341 float znew = i*tf[0][2] + j*tf[1][2];
03342 std::complex<float> btq;
03343 if (xnew < 0.)
03344 {
03345 xnew = -xnew;
03346 ynew = -ynew;
03347 znew = -znew;
03348 btq = conj(padfft->cmplx(i,jp-1));
03349 }
03350 else
03351 {
03352 btq = padfft->cmplx(i,jp-1);
03353 }
03354
03355 int ixn = int(xnew + 0.5 + n) - n;
03356 int iyn = int(ynew + 0.5 + n) - n;
03357 int izn = int(znew + 0.5 + n) - n;
03358 if ((ixn <= n2) && (iyn >= -n2) && (iyn <= n2) && (izn >= -n2) && (izn <= n2))
03359 {
03360 int ixf, iyf, izf;
03361 if (ixn >= 0)
03362 {
03363 int iza, iya;
03364 if (izn >= 0)
03365 iza = izn + 1;
03366 else
03367 iza = n + izn + 1;
03368
03369 if (iyn >= 0)
03370 iya = iyn + 1;
03371 else
03372 iya = n + iyn + 1;
03373
03374 ixf = ixn;
03375 iyf = iya;
03376 izf = iza;
03377 }
03378 else
03379 {
03380 int izt, iyt;
03381 if (izn > 0)
03382 izt = n - izn + 1;
03383 else
03384 izt = -izn + 1;
03385
03386 if (iyn > 0)
03387 iyt = n - iyn + 1;
03388 else
03389 iyt = -iyn + 1;
03390
03391 ixf = -ixn;
03392 iyf = iyt;
03393 izf = izt;
03394 }
03395
03396
03397 int pos2 = ixf + (iyf-1)*nx/2 + (izf-1)*ny*nx/2;
03398 float ctfv1 = btq.real() * ctf;
03399 float ctfv2 = btq.imag() * ctf;
03400 float ctf2 = ctf*ctf;
03401
03402 point_t p;
03403 p.pos2 = pos2;
03404 p.real = ctfv1;
03405 p.imag = ctfv2;
03406 p.ctf2 = ctf2;
03407
03408 points.push_back( p );
03409 }
03410 }
03411 }
03412 }
03413
03414
03415 int npoint = points.size();
03416 std::istream::off_type offset = (m_offsets.size()==0) ? 0 : m_offsets.back();
03417 offset += npoint*sizeof(point_t);
03418 m_offsets.push_back( offset );
03419
03420 *m_txt_of << m_offsets.back() << std::endl;
03421 m_bin_of->write( (char*)(&points[0]), sizeof(point_t)*npoint );
03422 checked_delete( padfft );
03423 }
03424
03425 void newfile_store::get_image( int id, EMData* buf )
03426 {
03427 if( m_offsets.size()==0 )
03428 {
03429 ifstream is( m_txt_file.c_str() );
03430 std::istream::off_type off;
03431 while( is >> off )
03432 {
03433 m_offsets.push_back( off );
03434 }
03435
03436 m_bin_if = shared_ptr<std::ifstream>( new ifstream(m_bin_file.c_str(), std::ios::in|std::ios::binary) );
03437 }
03438
03439 Assert( m_bin_if != NULL );
03440
03441 std::istream::off_type offset = (id==0) ? 0 : m_offsets[id-1];
03442 Assert( offset >= 0 );
03443 m_bin_if->seekg( offset, std::ios::beg );
03444
03445
03446 if( m_bin_if->bad() || m_bin_if->fail() || m_bin_if->eof() )
03447 {
03448 std::cout << "bad or fail or eof while fetching id, offset: " << id << " " << offset << std::endl;
03449 throw std::logic_error( "bad happen" );
03450 }
03451
03452 int bufsize = (m_offsets[id] - offset)/sizeof(float);
03453 if( buf->get_xsize() != bufsize )
03454 {
03455 buf->set_size( bufsize, 1, 1 );
03456 }
03457
03458 char* data = (char*)(buf->get_data());
03459 m_bin_if->read( data, sizeof(float)*bufsize );
03460 buf->update();
03461 }
03462
03463 void newfile_store::read( int nprj )
03464 {
03465 if( m_offsets.size()==0 )
03466 {
03467 ifstream is( m_txt_file.c_str() );
03468 std::istream::off_type off;
03469 while( is >> off )
03470 {
03471 m_offsets.push_back( off );
03472 }
03473 }
03474
03475 if( m_bin_if==NULL )
03476 {
03477 m_bin_if = shared_ptr< ifstream>( new ifstream(m_bin_file.c_str(), std::ios::in|std::ios::binary) );
03478 }
03479
03480
03481 int npoint = m_offsets[0]/sizeof(point_t);
03482 std::ios::off_type prjsize = m_offsets[0];
03483
03484 try
03485 {
03486 m_points.resize(nprj * npoint);
03487 }
03488 catch( std::exception& e )
03489 {
03490 std::cout << "Error: " << e.what() << std::endl;
03491 }
03492
03493 int ip = 0;
03494 for( int i=0; i < nprj; ++i )
03495 {
03496 m_bin_if->read( (char*)(&m_points[ip]), prjsize );
03497 if( m_bin_if->bad() || m_bin_if->fail() || m_bin_if->eof() )
03498 {
03499 std::cout << "Error: file hander bad or fail or eof" << std::endl;
03500 return;
03501 }
03502 ip += npoint;
03503 }
03504 }
03505
03506 void newfile_store::add_tovol( EMData* fftvol, EMData* wgtvol, const vector<int>& mults, int pbegin, int pend )
03507 {
03508 float* vdata = fftvol->get_data();
03509 float* wdata = wgtvol->get_data();
03510
03511 int npoint = m_offsets[0]/sizeof(point_t);
03512
03513 Assert( int(m_points.size())== (pend - pbegin)*npoint );
03514
03515 for( int iprj=pbegin; iprj < pend; ++iprj )
03516 {
03517 int m = mults[iprj];
03518 if( m==0 ) continue;
03519
03520 int ipt = (iprj-pbegin)*npoint;
03521 for( int i=0; i < npoint; ++i )
03522 {
03523 int pos2 = m_points[ipt].pos2;
03524 int pos1 = pos2*2;
03525
03526 wdata[pos2] += m_points[ipt].ctf2*m;
03527 vdata[pos1] += m_points[ipt].real*m;
03528 vdata[pos1+1]+= m_points[ipt].imag*m;
03529 ++ipt;
03530 }
03531 }
03532 }
03533
03534 void newfile_store::restart()
03535 {
03536 m_bin_if = shared_ptr< ifstream>( new ifstream(m_bin_file.c_str(), std::ios::in|std::ios::binary) );
03537 }
03538
03539 file_store::file_store(const string& filename, int npad, int write, bool ctf)
03540 : m_bin_file(filename + ".bin"),
03541 m_txt_file(filename + ".txt")
03542 {
03543 m_ctf = ctf;
03544 m_prev = -1;
03545 m_npad = npad;
03546 m_write = write;
03547 }
03548
03549 file_store::~file_store()
03550 {
03551 }
03552
03553 void file_store::add_image( EMData* emdata, const Transform& tf )
03554 {
03555
03556 EMData* padfft = padfft_slice( emdata, tf, m_npad );
03557
03558 float* data = padfft->get_data();
03559
03560 if( m_write && m_bin_ohandle == NULL )
03561 {
03562 m_bin_ohandle = shared_ptr< ofstream >( new ofstream(m_bin_file.c_str(), std::ios::out | std::ios::binary) );
03563 m_txt_ohandle = shared_ptr< ofstream >( new ofstream(m_txt_file.c_str() ) );
03564 if( m_ctf )
03565 *m_txt_ohandle << "Cs pixel voltage ctf_applied amp_contrast defocus ";
03566
03567 *m_txt_ohandle << "phi theta psi" << std::endl;
03568 }
03569
03570 m_x_out = padfft->get_xsize();
03571 m_y_out = padfft->get_ysize();
03572 m_z_out = padfft->get_zsize();
03573 m_totsize = m_x_out*m_y_out*m_z_out;
03574
03575 if( m_ctf )
03576 {
03577 Ctf* ctf = padfft->get_attr( "ctf" );
03578 Dict ctf_params = ctf->to_dict();
03579
03580 m_ctf_applied = padfft->get_attr( "ctf_applied" );
03581
03582 m_Cs = ctf_params["cs"];
03583 m_pixel = ctf_params["apix"];
03584 m_voltage = ctf_params["voltage"];
03585 m_amp_contrast = ctf_params["ampcont"];
03586 m_defocuses.push_back( ctf_params["defocus"] );
03587 if(ctf) {delete ctf; ctf=0;}
03588 }
03589
03590 Dict params = tf.get_rotation( "spider" );
03591 float phi = params.get( "phi" );
03592 float tht = params.get( "theta" );
03593 float psi = params.get( "psi" );
03594
03595
03596 m_phis.push_back( phi );
03597 m_thetas.push_back( tht );
03598 m_psis.push_back( psi );
03599
03600 if( m_write )
03601 {
03602 m_bin_ohandle->write( (char*)data, sizeof(float)*m_totsize );
03603
03604 if( m_ctf )
03605 {
03606 *m_txt_ohandle << m_Cs << " ";
03607 *m_txt_ohandle << m_pixel << " ";
03608 *m_txt_ohandle << m_voltage << " ";
03609 *m_txt_ohandle << m_ctf_applied << " ";
03610 *m_txt_ohandle << m_amp_contrast << " ";
03611 *m_txt_ohandle << m_defocuses.back() << " ";
03612 }
03613 *m_txt_ohandle << m_phis.back() << " ";
03614 *m_txt_ohandle << m_thetas.back() << " ";
03615 *m_txt_ohandle << m_psis.back() << " ";
03616 *m_txt_ohandle << m_x_out << " ";
03617 *m_txt_ohandle << m_y_out << " ";
03618 *m_txt_ohandle << m_z_out << " ";
03619 *m_txt_ohandle << m_totsize << std::endl;
03620 }
03621
03622 checked_delete(padfft);
03623
03624 }
03625
03626 void file_store::get_image( int id, EMData* padfft )
03627 {
03628
03629 if( m_phis.size() == 0 ) {
03630 ifstream m_txt_ifs( m_txt_file.c_str() );
03631
03632 if( !m_txt_ifs )
03633 {
03634 std::cerr << "Error: file " << m_txt_file << " does not exist" << std::endl;
03635 }
03636
03637 string line;
03638 std::getline( m_txt_ifs, line );
03639
03640 float first, defocus, phi, theta, psi;
03641
03642
03643
03644 while( m_txt_ifs >> first ) {
03645
03646 if( m_ctf )
03647 {
03648 m_Cs = first;
03649 m_txt_ifs >> m_pixel >> m_voltage;
03650 m_txt_ifs >> m_ctf_applied >> m_amp_contrast;
03651 m_txt_ifs >> defocus >> phi >> theta >> psi;
03652 m_defocuses.push_back( defocus );
03653 }
03654 else
03655 {
03656 phi = first;
03657 m_txt_ifs >> theta >> psi;
03658 }
03659
03660 m_txt_ifs >> m_x_out >> m_y_out >> m_z_out >> m_totsize;
03661 m_phis.push_back( phi );
03662 m_thetas.push_back( theta );
03663 m_psis.push_back( psi );
03664 }
03665 }
03666
03667 Assert( m_ihandle != NULL );
03668
03669 std::istream::off_type offset = id*sizeof(float)*m_totsize;
03670 Assert( offset >= 0 );
03671
03672 if( offset > 0 )
03673 {
03674 m_ihandle->seekg(offset, std::ios::beg);
03675 }
03676
03677 if( m_ihandle->bad() )
03678 {
03679 std::cout << "bad while fetching id, offset: " << id << " " << offset << std::endl;
03680 throw std::logic_error( "bad happen" );
03681 }
03682
03683 if( m_ihandle->fail() )
03684 {
03685 std::cout << "fail while fetching id, offset, curoff: " << id << " " << offset << std::endl;
03686 throw std::logic_error( "fail happen" );
03687 }
03688
03689 if( m_ihandle->eof() )
03690 {
03691 std::cout << "eof while fetching id, offset: " << id << " " << offset << std::endl;
03692 throw std::logic_error( "eof happen" );
03693 }
03694
03695 if( padfft->get_xsize() != m_x_out ||
03696 padfft->get_ysize() != m_y_out ||
03697 padfft->get_zsize() != m_z_out )
03698 {
03699 padfft->set_size(m_x_out, m_y_out, m_z_out);
03700 }
03701
03702 char* data = (char*)(padfft->get_data());
03703 m_ihandle->read( data, sizeof(float)*m_totsize );
03704 padfft->update();
03705
03706 if( m_ctf )
03707 {
03708 padfft->set_attr( "Cs", m_Cs );
03709 padfft->set_attr( "Pixel_size", m_pixel );
03710 padfft->set_attr( "voltage", m_voltage );
03711 padfft->set_attr( "ctf_applied", m_ctf_applied );
03712 padfft->set_attr( "amp_contrast", m_amp_contrast );
03713 padfft->set_attr( "defocus", m_defocuses[id] );
03714 }
03715
03716 padfft->set_attr( "padffted", 1 );
03717 padfft->set_attr( "phi", m_phis[id] );
03718 padfft->set_attr( "theta", m_thetas[id] );
03719 padfft->set_attr( "psi", m_psis[id] );
03720
03721 }
03722
03723 void file_store::restart( )
03724 {
03725 if( m_ihandle == NULL )
03726 {
03727 m_ihandle = shared_ptr< ifstream >( new ifstream(m_bin_file.c_str(), std::ios::in | std::ios::binary) );
03728 }
03729
03730 if( m_ihandle->bad() || m_ihandle->fail() || m_ihandle->eof() )
03731 {
03732 m_ihandle->open( m_bin_file.c_str(), std::ios::binary );
03733 }
03734
03735 m_ihandle->seekg( 0, std::ios::beg );
03736 }
03737
03738