EMAN2
emdata_transform.cpp
Go to the documentation of this file.
00001 
00005 /*
00006  * Author: Steven Ludtke, 04/10/2003 (sludtke@bcm.edu)
00007  * Copyright (c) 2000-2006 Baylor College of Medicine
00008  *
00009  * This software is issued under a joint BSD/GNU license. You may use the
00010  * source code in this file under either license. However, note that the
00011  * complete EMAN2 and SPARX software packages have some GPL dependencies,
00012  * so you are responsible for compliance with the licenses of these packages
00013  * if you opt to use BSD licensing. The warranty disclaimer below holds
00014  * in either instance.
00015  *
00016  * This complete copyright notice must be included in any revised version of the
00017  * source code. Additional authorship citations may be added, but existing
00018  * author citations must be preserved.
00019  *
00020  * This program is free software; you can redistribute it and/or modify
00021  * it under the terms of the GNU General Public License as published by
00022  * the Free Software Foundation; either version 2 of the License, or
00023  * (at your option) any later version.
00024  *
00025  * This program is distributed in the hope that it will be useful,
00026  * but WITHOUT ANY WARRANTY; without even the implied warranty of
00027  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
00028  * GNU General Public License for more details.
00029  *
00030  * You should have received a copy of the GNU General Public License
00031  * along with this program; if not, write to the Free Software
00032  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
00033  *
00034  * */
00035 
00036 #include "emdata.h"
00037 #include "emfft.h"
00038 
00039 #include <cstring>
00040 #include <cstdio>
00041 
00042 #include  "gsl_sf_result.h"
00043 #include  "gsl_sf_bessel.h"
00044 #include <iostream>
00045 #include <algorithm>
00046 #include <vector>
00047 #include <utility>
00048 #include <cmath>
00049 #include "util.h"
00050 
00051 //#ifdef EMAN2_USING_CUDA
00052 //#include "cuda/cuda_processor.h"
00053 //#endif
00054 
00055 using namespace EMAN;
00056 using namespace std;
00057 typedef vector< pair<float,int> > vp;
00058 
00059 EMData *EMData::do_fft() const
00060 {
00061         ENTERFUNC;
00062 #ifdef FFT_CACHING
00063         if (fftcache!=0) {
00064                 return fftcache->copy();
00065         }
00066 #endif //FFT_CACHING
00067 
00068         if (is_complex() ) { // ming add 08/17/2010
00069 #ifdef NATIVE_FFT
00070                 LOGERR(" NATIVE_FFT does not support complex to complex.");  // PAP
00071                 throw ImageFormatException("real image expected. Input image is complex image.");
00072 #else
00073                 EMData *temp_in=copy();
00074                 EMData *dat= copy_head();
00075                 int offset;
00076                 if(is_fftpadded()) {
00077                         offset = is_fftodd() ? 1 : 2;
00078                 }
00079                 else offset=0;
00080                 //printf("offset=%d\n",offset);
00081                 EMfft::complex_to_complex_nd(temp_in->get_data(),dat->get_data(),nx-offset,ny,nz);
00082 
00083                 if(dat->get_ysize()==1 && dat->get_zsize()==1) dat->set_complex_x(true);
00084 
00085                 dat->update();
00086                 delete temp_in;
00087                 EXITFUNC;
00088                 return dat;
00089 #endif // NATIVE_FFT
00090         } else {
00091                 int nxreal = nx;
00092                 int offset = 2 - nx%2;
00093                 int nx2 = nx + offset;
00094                 EMData* dat = copy_head();
00095                 dat->set_size(nx2, ny, nz);
00096                 //dat->to_zero();  // do not need it, real_to_complex will do it right anyway
00097                 if (offset == 1) dat->set_fftodd(true);
00098                 else             dat->set_fftodd(false);
00099 
00100                 float *d = dat->get_data();
00101                 //std::cout<<" do_fft "<<rdata[5]<<"  "<<d[5]<<std::endl;
00102                 EMfft::real_to_complex_nd(get_data(), d, nxreal, ny, nz);
00103 
00104                 dat->update();
00105                 dat->set_fftpad(true);
00106                 dat->set_complex(true);
00107                 if(dat->get_ysize()==1 && dat->get_zsize()==1) dat->set_complex_x(true);
00108                 dat->set_ri(true);
00109 
00110                 EXITFUNC;
00111 #ifdef FFT_CACHING
00112 //              printf("%p %d\n",this,nxyz);
00113                 if (nxyz<80000000) {
00114                         fftcache=dat->copy();
00115                 }
00116 #endif //FFT_CACHING
00117                 return dat;
00118         }
00119 }
00120 
00121 void EMData::do_fft_inplace()
00122 {
00123         ENTERFUNC;
00124 
00125         if ( is_complex() ) {
00126                 LOGERR("real image expected. Input image is complex image.");
00127                 throw ImageFormatException("real image expected. Input image is complex image.");
00128         }
00129 
00130         size_t offset;
00131         int nxreal;
00132         get_data(); // Required call if GPU caching is being used. Otherwise harmless
00133         if (!is_fftpadded()) {
00134                 // need to extend the matrix along x
00135                 // meaning nx is the un-fftpadded size
00136                 nxreal = nx;
00137                 offset = 2 - nx%2;
00138                 if (1 == offset) set_fftodd(true);
00139                 else             set_fftodd(false);
00140                 int nxnew = nx + offset;
00141                 set_size(nxnew, ny, nz);
00142                 for (int iz = nz-1; iz >= 0; iz--) {
00143                         for (int iy = ny-1; iy >= 0; iy--) {
00144                                 for (int ix = nxreal-1; ix >= 0; ix--) {
00145                                         size_t oldxpos = ix + (iy + iz*ny)*(size_t)nxreal;
00146                                         size_t newxpos = ix + (iy + iz*ny)*(size_t)nxnew;
00147                                         (*this)(newxpos) = (*this)(oldxpos);
00148                                 }
00149                         }
00150                 }
00151                 set_fftpad(true);
00152         } else {
00153                 offset = is_fftodd() ? 1 : 2;
00154                 nxreal = nx - offset;
00155         }
00156         EMfft::real_to_complex_nd(rdata, rdata, nxreal, ny, nz);
00157 
00158         set_complex(true);
00159         if(ny==1 && nz==1)  set_complex_x(true);
00160         set_ri(true);
00161 
00162         update();
00163 
00164         EXITFUNC;
00165         //return this;
00166 }
00167 
00168 #ifdef EMAN2_USING_CUDA
00169 
00170 #include "cuda/cuda_emfft.h"
00171 
00172 EMData *EMData::do_fft_cuda()
00173 {
00174         ENTERFUNC;
00175 
00176         if ( is_complex() ) {
00177                 LOGERR("real image expected. Input image is complex image.");
00178                 throw ImageFormatException("real image expected. Input image is complex image.");
00179         }
00180 
00181         int offset = 2 - nx%2;
00182         EMData* dat = new EMData(0,0,nx+offset,ny,nz,attr_dict);
00183         if(!dat->rw_alloc()) throw UnexpectedBehaviorException("Bad alloc");
00184         //cout << "Doing CUDA FFT " << cudarwdata << endl;
00185         if(cudarwdata == 0){copy_to_cuda();}
00186         cuda_dd_fft_real_to_complex_nd(cudarwdata, dat->cudarwdata, nx, ny,nz, 1);
00187 
00188         if (offset == 1) dat->set_fftodd(true);
00189         else             dat->set_fftodd(false);
00190 
00191         dat->set_fftpad(true);
00192         dat->set_complex(true);
00193         if(dat->get_ysize()==1 && dat->get_zsize()==1) dat->set_complex_x(true);
00194         dat->set_ri(true);
00195         dat->update();
00196 
00197         EXITFUNC;
00198         return dat;
00199 }
00200 
00201 void EMData::do_fft_inplace_cuda()
00202 {
00203         ENTERFUNC;
00204 
00205         if ( is_complex() ) {
00206                 LOGERR("real image expected. Input image is complex image.");
00207                 throw ImageFormatException("real image expected. Input image is complex image.");
00208         }
00209 
00210         int offset = 2 - nx%2;
00211         float* tempcudadata = 0;
00212         cudaError_t error = cudaMalloc((void**)&tempcudadata,(nx + offset)*ny*nz*sizeof(float));
00213         if( error != cudaSuccess) throw ImageFormatException("Couldn't allocate memory.");
00214         
00215         //cout << "Doing CUDA FFT inplace" << cudarwdata << endl;
00216         if(cudarwdata == 0){copy_to_cuda();}
00217         cuda_dd_fft_real_to_complex_nd(cudarwdata, tempcudadata, nx, ny,nz, 1);
00218         // this section is a bit slight of hand it actually does the FFT out of place but this avoids and EMData object creation and detruction...
00219         cudaError_t ferror = cudaFree(cudarwdata);
00220         if ( ferror != cudaSuccess) throw UnexpectedBehaviorException( "CudaFree failed:" + string(cudaGetErrorString(error)));
00221         cudarwdata = tempcudadata;
00222         num_bytes = (nx + offset)*ny*nz*sizeof(float);
00223 
00224         if (offset == 1) set_fftodd(true);
00225         else             set_fftodd(false);
00226 
00227         nx = nx + offset; // don't want to call set_size b/c that will delete my cudadata, remember what I am doing is a bit slignt of hand....
00228         set_fftpad(true);
00229         set_complex(true);
00230         if(get_ysize()==1 && get_zsize()==1) set_complex_x(true);
00231         set_ri(true);
00232         update();
00233 
00234         EXITFUNC;
00235 //      return this;
00236 }
00237 
00238 EMData *EMData::do_ift_cuda()
00239 {
00240         ENTERFUNC;
00241 
00242         if (!is_complex()) {
00243                 LOGERR("complex image expected. Input image is real image.");
00244                 throw ImageFormatException("complex image expected. Input image is real image.");
00245         }
00246 
00247         if (!is_ri()) {
00248                 throw ImageFormatException("complex ri expected. Got amplitude/phase.");
00249         }
00250 
00251         int offset = is_fftodd() ? 1 : 2;
00252         EMData* dat = new EMData(0,0,nx-offset,ny,nz,attr_dict);
00253         if(!dat->rw_alloc()) throw UnexpectedBehaviorException("Bad alloc");
00254         
00255         if(cudarwdata == 0){copy_to_cuda();}
00256 
00257         
00258         int ndim = get_ndim();
00259         if ( ndim == 1 ) {
00260                 cuda_dd_fft_complex_to_real_nd(cudarwdata,dat->cudarwdata, nx-offset,1,1,1);
00261         } else if (ndim == 2) {
00262                 cuda_dd_fft_complex_to_real_nd(cudarwdata,dat->cudarwdata, ny,nx-offset,1,1);
00263         } else if (ndim == 3) {
00264                 cuda_dd_fft_complex_to_real_nd(cudarwdata,dat->cudarwdata, nz,ny,nx-offset,1);
00265         } else throw ImageDimensionException("No cuda FFT support of images with dimensions exceeding 3");
00266         
00267         // SCALE the inverse FFT
00268         float scale = 1.0f/static_cast<float>((dat->get_size()));
00269         dat->mult(scale); 
00270 
00271         dat->set_fftpad(false);
00272         dat->set_fftodd(false);
00273         dat->set_complex(false);
00274         if(dat->get_ysize()==1 && dat->get_zsize()==1)  dat->set_complex_x(false);
00275         dat->set_ri(false);
00276 //      dat->gpu_update();
00277         dat->update(); 
00278         
00279         EXITFUNC;
00280         return dat;
00281 }
00282 
00283 /*
00284    FFT in place does not depad, hence this routine is of limited use b/c mem operations on the device are quite SLOW, JFF
00285    use
00286 */
00287 
00288 void EMData::do_ift_inplace_cuda()
00289 {
00290         ENTERFUNC;
00291 
00292         if (!is_complex()) {
00293                 LOGERR("complex image expected. Input image is real image.");
00294                 throw ImageFormatException("complex image expected. Input image is real image.");
00295         }
00296 
00297         if (!is_ri()) {
00298                 LOGWARN("run IFT on AP data, only RI should be used. ");
00299         }
00300 
00301         int offset = is_fftodd() ? 1 : 2;
00302         
00303         if(cudarwdata == 0){copy_to_cuda();}
00304         
00305         int ndim = get_ndim();
00306         if ( ndim == 1 ) {
00307                 cuda_dd_fft_complex_to_real_nd(cudarwdata,cudarwdata, nx-offset,1,1,1);
00308         } else if (ndim == 2) {
00309                 cuda_dd_fft_complex_to_real_nd(cudarwdata,cudarwdata, ny,nx-offset,1,1);
00310         } else if (ndim == 3) {
00311                 cuda_dd_fft_complex_to_real_nd(cudarwdata,cudarwdata, nz,ny,nx-offset,1);
00312         } else throw ImageDimensionException("No cuda FFT support of images with dimensions exceeding 3");
00313 #if defined     FFTW2 || defined FFTW3 //native fft and ACML already done normalization
00314         // SCALE the inverse FFT
00315         int nxo = nx - offset;
00316         float scale = 1.0f / (nxo * ny * nz);
00317         mult(scale); //if we are just doing a CCF, this is a waste!
00318 #endif //FFTW2 || FFTW3
00319 
00320         set_fftpad(true);
00321         set_complex(false);
00322 
00323         if(ny==1 && nz==1) set_complex_x(false);
00324         set_ri(false);
00325         update();
00326         
00327         EXITFUNC;
00328 //      return this;
00329 }
00330 
00331 #endif //EMAN2_USING_CUDA
00332 
00333 EMData *EMData::do_ift()
00334 {
00335         ENTERFUNC;
00336 
00337         if (!is_complex()) {
00338                 LOGERR("complex image expected. Input image is real image.");
00339                 throw ImageFormatException("complex image expected. Input image is real image.");
00340         }
00341 
00342         if (!is_ri()) {
00343                 LOGWARN("run IFT on AP data, only RI should be used. Converting.");
00344         }
00345 
00346         get_data(); // Required call if GPU caching is being used. Otherwise harmless
00347         EMData* dat = copy_head();
00348         dat->set_size(nx, ny, nz);
00349         ap2ri();
00350 
00351         float *d = dat->get_data();
00352         int ndim = get_ndim();
00353 
00354         /* Do inplace IFT on a image copy, because the complex to real transform of
00355          * nd will destroy its input array even for out-of-place transforms.
00356          */
00357         memcpy((char *) d, (char *) rdata, (size_t)nx * ny * nz * sizeof(float));
00358 
00359         int offset = is_fftodd() ? 1 : 2;
00360         //cout << "Sending offset " << offset << " " << nx-offset << endl;
00361         if (ndim == 1) {
00362                 EMfft::complex_to_real_nd(d, d, nx - offset, ny, nz);
00363         } else {
00364                 EMfft::complex_to_real_nd(d, d, nx - offset, ny, nz);
00365 
00366                 size_t row_size = (nx - offset) * sizeof(float);
00367                 for (size_t i = 1; i < (size_t)ny * nz; i++) {
00368                         memmove((char *) &d[i * (nx - offset)], (char *) &d[i * nx], row_size);
00369                 }
00370         }
00371 
00372         dat->set_size(nx - offset, ny, nz);     //remove the padding
00373 #if defined     FFTW2 || defined FFTW3 //native fft and ACML already done normalization
00374         // SCALE the inverse FFT
00375         float scale = 1.0f / ((nx - offset) * ny * nz);
00376         dat->mult(scale);
00377 #endif  //FFTW2 || FFTW3
00378         dat->set_fftodd(false);
00379         dat->set_fftpad(false);
00380         dat->set_complex(false);
00381         if(dat->get_ysize()==1 && dat->get_zsize()==1)  dat->set_complex_x(false);
00382         dat->set_ri(false);
00383         dat->update();
00384 
00385 
00386         EXITFUNC;
00387         return dat;
00388 }
00389 
00390 /*
00391    FFT in place does not depad, return real x-extended image (needs to be depadded before use as PAP does in CCF routines)
00392    use
00393 */
00394 void EMData::do_ift_inplace()
00395 {
00396         ENTERFUNC;
00397 
00398         if (!is_complex()) {
00399                 LOGERR("complex image expected. Input image is real image.");
00400                 throw ImageFormatException("complex image expected. Input image is real image.");
00401         }
00402 
00403         if (!is_ri()) {
00404                 LOGWARN("run IFT on AP data, only RI should be used. ");
00405         }
00406         ap2ri();
00407 
00408         int offset = is_fftodd() ? 1 : 2;
00409         float* data = get_data();
00410         EMfft::complex_to_real_nd(data, data, nx - offset, ny, nz);
00411 
00412 #if defined     FFTW2 || defined FFTW3  //native fft and ACML already done normalization
00413         // SCALE the inverse FFT
00414         int nxo = nx - offset;
00415         float scale = 1.0f / ((size_t)nxo * ny * nz);
00416         mult(scale);
00417 #endif //FFTW2 || FFTW3
00418 
00419         set_fftpad(true);
00420         set_complex(false);
00421         if(ny==1 && nz==1) set_complex_x(false);
00422         set_ri(false);
00423         update();
00424 
00425         EXITFUNC;
00426 //      return this;
00427 }
00428 #undef rdata
00429 
00430 
00431 std::string EMData::render_ap24(int x0, int y0, int ixsize, int iysize,
00432                                                  int bpl, float scale, int mingray, int maxgray,
00433                                                  float render_min, float render_max,float gamma,int flags)
00434 {
00435         ENTERFUNC;
00436 
00437         int asrgb;
00438         int hist=(flags&2)/2;
00439         int invy=(flags&4)?1:0;
00440 
00441         if (!is_complex()) throw ImageDimensionException("complex only");
00442 
00443         if (get_ndim() != 2) {
00444                 throw ImageDimensionException("2D only");
00445         }
00446 
00447         if (is_complex()) ri2ap();
00448 
00449         if (render_max <= render_min) {
00450                 render_max = render_min + 0.01f;
00451         }
00452 
00453         if (gamma<=0) gamma=1.0;
00454 
00455         // Calculating a full floating point gamma for
00456         // each pixel in the image slows rendering unacceptably
00457         // however, applying a gamma-mapping to an 8 bit colorspace
00458         // has unaccepable accuracy. So, we oversample the 8 bit colorspace
00459         // as a 12 bit colorspace and apply the gamma mapping to that
00460         // This should produce good accuracy for gamma values
00461         // larger than 0.5 (and a high upper limit)
00462         static int smg0=0,smg1=0;       // while this destroys threadsafety in the rendering process
00463         static float sgam=0;            // it is necessary for speed when rendering large numbers of small images
00464         static unsigned char gammamap[4096];
00465         if (gamma!=1.0 && (smg0!=mingray || smg1!=maxgray || sgam!=gamma)) {
00466                 for (int i=0; i<4096; i++) {
00467                         if (mingray<maxgray) gammamap[i]=(unsigned char)(mingray+(maxgray-mingray+0.999)*pow(((float)i/4096.0f),gamma));
00468                         else gammamap[4095-i]=(unsigned char)(mingray+(maxgray-mingray+0.999)*pow(((float)i/4096.0f),gamma));
00469                 }
00470         }
00471         smg0=mingray;   // so we don't recompute the map unless something changes
00472         smg1=maxgray;
00473         sgam=gamma;
00474 
00475         if (flags&8) asrgb=4;
00476         else if (flags&1) asrgb=3;
00477         else throw ImageDimensionException("must set flag 1 or 8");
00478 
00479         std::string ret=std::string();
00480 //      ret.resize(iysize*bpl);
00481         ret.assign(iysize*bpl+hist*1024,char(mingray));
00482         unsigned char *data=(unsigned char *)ret.data();
00483         unsigned int *histd=(unsigned int *)(data+iysize*bpl);
00484         if (hist) {
00485                 for (int i=0; i<256; i++) histd[i]=0;
00486         }
00487 
00488         float rm = render_min;
00489         float inv_scale = 1.0f / scale;
00490         int ysize = iysize;
00491         int xsize = ixsize;
00492 
00493         int ymin = 0;
00494         if (iysize * inv_scale > ny) {
00495                 ymin = (int) (iysize - ny / inv_scale);
00496         }
00497 
00498         float gs = (maxgray - mingray) / (render_max - render_min);
00499         float gs2 = 4095.999f / (render_max - render_min);
00500 //      float gs2 = 1.0 / (render_max - render_min);
00501         if (render_max < render_min) {
00502                 gs = 0;
00503                 rm = FLT_MAX;
00504         }
00505 
00506         int dsx = -1;
00507         int dsy = 0;
00508         int remx = 0;
00509         int remy = 0;
00510         const int scale_n = 100000;
00511 
00512         int addi = 0;
00513         int addr = 0;
00514         if (inv_scale == floor(inv_scale)) {
00515                 dsx = (int) inv_scale;
00516                 dsy = (int) (inv_scale * nx);
00517         }
00518         else {
00519                 addi = (int) floor(inv_scale);
00520                 addr = (int) (scale_n * (inv_scale - floor(inv_scale)));
00521         }
00522 
00523         int xmin = 0;
00524         if (x0 < 0) {
00525                 xmin = (int) (-x0 / inv_scale);
00526                 xsize -= (int) floor(x0 / inv_scale);
00527                 x0 = 0;
00528         }
00529 
00530         if ((xsize - xmin) * inv_scale > (nx - x0)) {
00531                 xsize = (int) ((nx - x0) / inv_scale + xmin);
00532         }
00533         int ymax = ysize - 1;
00534         if (y0 < 0) {
00535                 ymax = (int) (ysize + y0 / inv_scale - 1);
00536                 ymin += (int) floor(y0 / inv_scale);
00537                 y0 = 0;
00538         }
00539 
00540         if (xmin < 0) xmin = 0;
00541         if (ymin < 0) ymin = 0;
00542         if (xsize > ixsize) xsize = ixsize;
00543         if (ymax > iysize) ymax = iysize;
00544 
00545         int lmax = nx * ny - 1;
00546 
00547         int mid=nx*ny/2;
00548         float* image_data = get_data();
00549         if (dsx != -1) {
00550                 int l = y0 * nx;
00551                 for (int j = ymax; j >= ymin; j--) {
00552                         int ll = x0;
00553                         for (int i = xmin; i < xsize; i++) {
00554                                 if (l + ll > lmax || ll >= nx - 2) break;
00555 
00556                                 int k = 0;
00557                                 unsigned char p;
00558                                 int ph;
00559                                 if (ll >= nx / 2) {
00560                                         if (l >= (ny - inv_scale) * nx) k = 2 * (ll - nx / 2) + 2;
00561                                         else k = 2 * (ll - nx / 2) + l + 2 + nx;
00562                                         if (k>=mid) k-=mid;             // These 2 lines handle the Fourier origin being in the corner, not the middle
00563                                         else k+=mid;
00564                                         ph = (int)(image_data[k+1]*768/(2.0*M_PI))+384; // complex phase as integer 0-767
00565                                 }
00566                                 else {
00567                                         k = nx * ny - (l + 2 * ll) - 2;
00568                                         ph = (int)(-image_data[k+1]*768/(2.0*M_PI))+384;        // complex phase as integer 0-767
00569                                         if (k>=mid) k-=mid;             // These 2 lines handle the Fourier origin being in the corner, not the middle
00570                                         else k+=mid;
00571                                 }
00572                                 float t = image_data[k];
00573                                 if (t <= rm)  p = mingray;
00574                                 else if (t >= render_max) p = maxgray;
00575                                 else if (gamma!=1.0) {
00576                                         k=(int)(gs2 * (t-render_min));          // map float value to 0-4096 range
00577                                         p = gammamap[k];                                        // apply gamma using precomputed gamma map
00578                                 }
00579                                 else {
00580                                         p = (unsigned char) (gs * (t - render_min));
00581                                         p += mingray;
00582                                 }
00583                                 if (ph<256) {
00584                                         data[i * asrgb + j * bpl] = p*(255-ph)/256;
00585                                         data[i * asrgb + j * bpl+1] = p*ph/256;
00586                                         data[i * asrgb + j * bpl+2] = 0;
00587                                 }
00588                                 else if (ph<512) {
00589                                         data[i * asrgb + j * bpl+1] = p*(511-ph)/256;
00590                                         data[i * asrgb + j * bpl+2] = p*(ph-256)/256;
00591                                         data[i * asrgb + j * bpl] = 0;
00592                                 }
00593                                 else {
00594                                         data[i * asrgb + j * bpl+2] = p*(767-ph)/256;
00595                                         data[i * asrgb + j * bpl] = p*(ph-512)/256;
00596                                         data[i * asrgb + j * bpl+1] = 0;
00597                                 }
00598                                 if (hist) histd[p]++;
00599                                 ll += dsx;
00600                         }
00601                         l += dsy;
00602                 }
00603         }
00604         else {
00605                 remy = 10;
00606                 int l = y0 * nx;
00607                 for (int j = ymax; j >= ymin; j--) {
00608                         int br = l;
00609                         remx = 10;
00610                         int ll = x0;
00611                         for (int i = xmin; i < xsize - 1; i++) {
00612                                 if (l + ll > lmax || ll >= nx - 2) {
00613                                         break;
00614                                 }
00615                                 int k = 0;
00616                                 unsigned char p;
00617                                 int ph;
00618                                 if (ll >= nx / 2) {
00619                                         if (l >= (ny * nx - nx)) k = 2 * (ll - nx / 2) + 2;
00620                                         else k = 2 * (ll - nx / 2) + l + 2 + nx;
00621                                         if (k>=mid) k-=mid;             // These 2 lines handle the Fourier origin being in the corner, not the middle
00622                                         else k+=mid;
00623                                         ph = (int)(image_data[k+1]*768/(2.0*M_PI))+384; // complex phase as integer 0-767
00624                                 }
00625                                 else {
00626                                         k = nx * ny - (l + 2 * ll) - 2;
00627                                         if (k>=mid) k-=mid;             // These 2 lines handle the Fourier origin being in the corner, not the middle
00628                                         else k+=mid;
00629                                         ph = (int)(-image_data[k+1]*768/(2.0*M_PI))+384;        // complex phase as integer 0-767
00630                                 }
00631 
00632                                 float t = image_data[k];
00633                                 if (t <= rm)
00634                                         p = mingray;
00635                                 else if (t >= render_max) {
00636                                         p = maxgray;
00637                                 }
00638                                 else if (gamma!=1.0) {
00639                                         k=(int)(gs2 * (t-render_min));          // map float value to 0-4096 range
00640                                         p = gammamap[k];                                        // apply gamma using precomputed gamma map
00641                                 }
00642                                 else {
00643                                         p = (unsigned char) (gs * (t - render_min));
00644                                         p += mingray;
00645                                 }
00646                                 if (ph<256) {
00647                                         data[i * asrgb + j * bpl] = p*(255-ph)/256;
00648                                         data[i * asrgb + j * bpl+1] = p*ph/256;
00649                                         data[i * asrgb + j * bpl+2] = 0;
00650                                 }
00651                                 else if (ph<512) {
00652                                         data[i * asrgb + j * bpl+1] = p*(511-ph)/256;
00653                                         data[i * asrgb + j * bpl+2] = p*(ph-256)/256;
00654                                         data[i * asrgb + j * bpl] = 0;
00655                                 }
00656                                 else {
00657                                         data[i * asrgb + j * bpl+2] = p*(767-ph)/256;
00658                                         data[i * asrgb + j * bpl] = p*(ph-512)/256;
00659                                         data[i * asrgb + j * bpl+1] = 0;
00660                                 }
00661                                 if (hist) histd[p]++;
00662                                 ll += addi;
00663                                 remx += addr;
00664                                 if (remx > scale_n) {
00665                                         remx -= scale_n;
00666                                         ll++;
00667                                 }
00668                         }
00669                         l = br + addi * nx;
00670                         remy += addr;
00671                         if (remy > scale_n) {
00672                                 remy -= scale_n;
00673                                 l += nx;
00674                         }
00675                 }
00676         }
00677 
00678         // this replicates r -> g,b
00679         if (asrgb==4) {
00680                 for (int j=ymin*bpl; j<=ymax*bpl; j+=bpl) {
00681                         for (int i=xmin; i<xsize*4; i+=4) {
00682                                 data[i+j+3]=255;
00683                         }
00684                 }
00685         }
00686 
00687         EXITFUNC;
00688 
00689         // ok, ok, not the most efficient place to do this, but it works
00690         if (invy) {
00691                 int x,y;
00692                 char swp;
00693                 for (y=0; y<iysize/2; y++) {
00694                         for (x=0; x<ixsize; x++) {
00695                                 swp=ret[y*bpl+x];
00696                                 ret[y*bpl+x]=ret[(iysize-y-1)*bpl+x];
00697                                 ret[(iysize-y-1)*bpl+x]=swp;
00698                         }
00699                 }
00700         }
00701 
00702     //  return PyString_FromStringAndSize((const char*) data,iysize*bpl);
00703         return ret;
00704 }
00705 
00706 
00707 void EMData::render_amp24( int x0, int y0, int ixsize, int iysize,
00708                                                   int bpl, float scale, int mingray, int maxgray,
00709                                                   float render_min, float render_max, void *ref,
00710                                                   void cmap(void *, int coord, unsigned char *tri))
00711 {
00712         ENTERFUNC;
00713 
00714         if (get_ndim() != 2) {
00715                 throw ImageDimensionException("2D only");
00716         }
00717 
00718         if (is_complex()) {
00719                 ri2ap();
00720         }
00721 
00722         if (render_max <= render_min) {
00723                 render_max = render_min + 0.01f;
00724         }
00725 
00726         std::string ret=std::string();
00727         ret.resize(iysize*bpl);
00728         unsigned char *data=(unsigned char *)ret.data();
00729 
00730         float rm = render_min;
00731         float inv_scale = 1.0f / scale;
00732         int ysize = iysize;
00733         int xsize = ixsize;
00734         const int scale_n = 100000;
00735 
00736         int ymin = 0;
00737         if ( iysize * inv_scale > ny) {
00738                 ymin = (int) (iysize - ny / inv_scale);
00739         }
00740         float gs = (maxgray - mingray) / (render_max - render_min);
00741         if (render_max < render_min) {
00742                 gs = 0;
00743                 rm = FLT_MAX;
00744         }
00745         int dsx = -1;
00746         int dsy = 0;
00747         if (inv_scale == floor(inv_scale)) {
00748                 dsx = (int) inv_scale;
00749                 dsy = (int) (inv_scale * nx);
00750         }
00751         int addi = 0;
00752         int addr = 0;
00753 
00754         if (dsx == -1) {
00755                 addi = (int) floor(inv_scale);
00756                 addr = (int) (scale_n * (inv_scale - floor(inv_scale)));
00757         }
00758 
00759         int remx = 0;
00760         int remy = 0;
00761         int xmin = 0;
00762         if (x0 < 0) {
00763                 xmin = (int) (-x0 / inv_scale);
00764                 xsize -= (int) floor(x0 / inv_scale);
00765                 x0 = 0;
00766         }
00767 
00768         if ((xsize - xmin) * inv_scale > (nx - x0)) {
00769                 xsize = (int) ((nx - x0) / inv_scale + xmin);
00770         }
00771         int ymax = ysize - 1;
00772         if (y0 < 0) {
00773                 ymax = (int) (ysize + y0 / inv_scale - 1);
00774                 ymin += (int) floor(y0 / inv_scale);
00775                 y0 = 0;
00776         }
00777 
00778 
00779         if (xmin < 0) {
00780                 xmin = 0;
00781         }
00782 
00783         if (ymin < 0) {
00784                 ymin = 0;
00785         }
00786         if (xsize > ixsize) {
00787                 xsize = ixsize;
00788         }
00789         if (ymax > iysize) {
00790                 ymax = iysize;
00791         }
00792 
00793         int lmax = nx * ny - 1;
00794         unsigned char tri[3];
00795         float* image_data = get_data();
00796         if (is_complex()) {
00797                 if (dsx != -1) {
00798                         int l = y0 * nx;
00799                         for (int j = ymax; j >= ymin; j--) {
00800                                 int ll = x0;
00801                                 for (int i = xmin; i < xsize; i++, ll += dsx) {
00802                                         if (l + ll > lmax || ll >= nx - 2) {
00803                                                 break;
00804                                         }
00805                                         int kk = 0;
00806                                         if (ll >= nx / 2) {
00807                                                 if (l >= (ny - inv_scale) * nx) {
00808                                                         kk = 2 * (ll - nx / 2) + 2;
00809                                                 }
00810                                                 else {
00811                                                         kk = 2 * (ll - nx / 2) + l + 2 + nx;
00812                                                 }
00813                                         }
00814                                         else {
00815                                                 kk = nx * ny - (l + 2 * ll) - 2;
00816                                         }
00817                                         int k = 0;
00818                                         float t = image_data[kk];
00819                                         if (t <= rm) {
00820                                                 k = mingray;
00821                                         }
00822                                         else if (t >= render_max) {
00823                                                 k = maxgray;
00824                                         }
00825                                         else {
00826                                                 k = (int) (gs * (t - render_min));
00827                                                 k += mingray;
00828                                         }
00829                                         tri[0] = static_cast < unsigned char >(k);
00830                                         cmap(ref, kk, tri);
00831                                         data[i * 3 + j * bpl] = tri[0];
00832                                         data[i * 3 + 1 + j * bpl] = tri[1];
00833                                         data[i * 3 + 2 + j * bpl] = tri[2];
00834                                 }
00835                                 l += dsy;
00836                         }
00837                 }
00838                 else {
00839                         remy = 10;
00840                         for (int j = ymax, l = y0 * nx; j >= ymin; j--) {
00841                                 int br = l;
00842                                 remx = 10;
00843                                 for (int i = xmin, ll = x0; i < xsize - 1; i++) {
00844                                         if (l + ll > lmax || ll >= nx - 2) {
00845                                                 break;
00846                                         }
00847                                         int kk = 0;
00848                                         if (ll >= nx / 2) {
00849                                                 if (l >= (ny * nx - nx)) {
00850                                                         kk = 2 * (ll - nx / 2) + 2;
00851                                                 }
00852                                                 else {
00853                                                         kk = 2 * (ll - nx / 2) + l + 2 + nx;
00854                                                 }
00855                                         }
00856                                         else {
00857                                                 kk = nx * ny - (l + 2 * ll) - 2;
00858                                         }
00859                                         int k = 0;
00860                                         float t = image_data[kk];
00861                                         if (t <= rm) {
00862                                                 k = mingray;
00863                                         }
00864                                         else if (t >= render_max) {
00865                                                 k = maxgray;
00866                                         }
00867                                         else {
00868                                                 k = (int) (gs * (t - render_min));
00869                                                 k += mingray;
00870                                         }
00871                                         tri[0] = static_cast < unsigned char >(k);
00872                                         cmap(ref, kk, tri);
00873                                         data[i * 3 + j * bpl] = tri[0];
00874                                         data[i * 3 + 1 + j * bpl] = tri[1];
00875                                         data[i * 3 + 2 + j * bpl] = tri[2];
00876                                         ll += addi;
00877                                         remx += addr;
00878                                         if (remx > scale_n) {
00879                                                 remx -= scale_n;
00880                                                 ll++;
00881                                         }
00882                                 }
00883                                 l = br + addi * nx;
00884                                 remy += addr;
00885                                 if (remy > scale_n) {
00886                                         remy -= scale_n;
00887                                         l += nx;
00888                                 }
00889                         }
00890                 }
00891         }
00892         else {
00893                 if (dsx != -1) {
00894                         for (int j = ymax, l = x0 + y0 * nx; j >= ymin; j--) {
00895                                 int br = l;
00896                                 for (int i = xmin; i < xsize; i++, l += dsx) {
00897                                         if (l > lmax) {
00898                                                 break;
00899                                         }
00900                                         float t = image_data[l];
00901                                         int k = 0;
00902                                         if (t <= rm) {
00903                                                 k = mingray;
00904                                         }
00905                                         else if (t >= render_max) {
00906                                                 k = maxgray;
00907                                         }
00908                                         else {
00909                                                 k = (int) (gs * (t - render_min));
00910                                                 k += mingray;
00911                                         }
00912                                         tri[0] = static_cast < unsigned char >(k);
00913                                         cmap(ref, l, tri);
00914                                         data[i * 3 + j * bpl] = tri[0];
00915                                         data[i * 3 + 1 + j * bpl] = tri[1];
00916                                         data[i * 3 + 2 + j * bpl] = tri[2];
00917                                 }
00918                                 l = br + dsy;
00919                         }
00920                 }
00921                 else {
00922                         remy = 10;
00923                         for (int j = ymax, l = x0 + y0 * nx; j >= ymin; j--) {
00924                                 int br = l;
00925                                 remx = 10;
00926                                 for (int i = xmin; i < xsize; i++) {
00927                                         if (l > lmax) {
00928                                                 break;
00929                                         }
00930                                         float t = image_data[l];
00931                                         int k = 0;
00932                                         if (t <= rm) {
00933                                                 k = mingray;
00934                                         }
00935                                         else if (t >= render_max) {
00936                                                 k = maxgray;
00937                                         }
00938                                         else {
00939                                                 k = (int) (gs * (t - render_min));
00940                                                 k += mingray;
00941                                         }
00942                                         tri[0] = static_cast < unsigned char >(k);
00943                                         cmap(ref, l, tri);
00944                                         data[i * 3 + j * bpl] = tri[0];
00945                                         data[i * 3 + 1 + j * bpl] = tri[1];
00946                                         data[i * 3 + 2 + j * bpl] = tri[2];
00947                                         l += addi;
00948                                         remx += addr;
00949                                         if (remx > scale_n) {
00950                                                 remx -= scale_n;
00951                                                 l++;
00952                                         }
00953                                 }
00954                                 l = br + addi * nx;
00955                                 remy += addr;
00956                                 if (remy > scale_n) {
00957                                         remy -= scale_n;
00958                                         l += nx;
00959                                 }
00960                         }
00961                 }
00962         }
00963 
00964         EXITFUNC;
00965 }
00966 
00967 void EMData::ap2ri()
00968 {
00969         ENTERFUNC;
00970 
00971         if (!is_complex() || is_ri()) {
00972                 return;
00973         }
00974         
00975         Util::ap2ri(get_data(), (size_t)nx * ny * nz);
00976         set_ri(true);
00977         update();
00978         EXITFUNC;
00979 }
00980 
00981 void EMData::ri2inten()
00982 {
00983         ENTERFUNC;
00984 
00985         if (!is_complex()) return;
00986         if (!is_ri()) ap2ri();
00987 
00988         float * data = get_data();
00989         size_t size = (size_t)nx * ny * nz;
00990         for (size_t i = 0; i < size; i += 2) {
00991                 data[i]=data[i]*data[i]+data[i+1]*data[i+1];
00992                 data[i+1]=0;
00993         }
00994 
00995         set_attr("is_intensity", int(1));
00996         update();
00997         EXITFUNC;
00998 }
00999 
01000 
01001 void EMData::ri2ap()
01002 {
01003         ENTERFUNC;
01004 
01005         if (!is_complex() || !is_ri()) {
01006                 return;
01007         }
01008 
01009         float * data = get_data();
01010 
01011         size_t size = (size_t)nx * ny * nz;
01012         for (size_t i = 0; i < size; i += 2) {
01013 #ifdef  _WIN32
01014                 float f = (float)_hypot(data[i], data[i + 1]);
01015 #else
01016                 float f = (float)hypot(data[i], data[i + 1]);
01017 #endif
01018                 if (data[i] == 0 && data[i + 1] == 0) {
01019                         data[i + 1] = 0;
01020                 }
01021                 else {
01022                         data[i + 1] = atan2(data[i + 1], data[i]);
01023                 }
01024                 data[i] = f;
01025         }
01026 
01027         set_ri(false);
01028         update();
01029         EXITFUNC;
01030 }
01031 
01032 
01033 float calc_bessel(const int n, const float& x) {
01034         gsl_sf_result result;
01035 //      int success =
01036         gsl_sf_bessel_Jn_e(n,(double)x, &result);
01037         return (float)result.val;
01038 }
01039 
01040 EMData*   EMData::bispecRotTransInvN(int N, int NK)
01041 {
01042 
01043         int EndP = this -> get_xsize(); // length(fTrueVec);
01044         int Mid  = (int) ((1+EndP)/2);
01045         int End = 2*Mid-1;
01046 
01047         int CountxyMax = End*End;
01048 
01049         int   *SortfkInds       = new    int[CountxyMax];
01050         int   *kVecX            = new    int[CountxyMax];
01051         int   *kVecY            = new    int[CountxyMax];
01052         float *fkVecR           = new  float[CountxyMax];
01053         float *fkVecI           = new  float[CountxyMax];
01054         float *absD1fkVec       = new  float[CountxyMax];
01055         float *absD1fkVecSorted = new  float[CountxyMax];
01056 
01057         float *jxjyatan2         = new  float[End*End]; //  jxjyatan2[jy*End + jx]  = atan2(jy+1-Mid , jx +1 -Mid);
01058 
01059         EMData * ThisCopy = new EMData(End,End);
01060 
01061         for (int jx=0; jx <End ; jx++) {
01062                 for (int jy=0; jy <End ; jy++) {
01063                         float ValNow = this -> get_value_at(jx,jy);
01064                         ThisCopy -> set_value_at(jx,jy,ValNow);
01065 //              cout<< " jxM= " << jx+1<<" jyM= " << jy+1<< "ValNow" << ValNow << endl; //    Works
01066         }}
01067 
01068 
01069         EMData* fk = ThisCopy -> do_fft();
01070         fk          ->process_inplace("xform.fourierorigin.tocenter");
01071 
01072 //      EMData* fk
01073         EMData* fkRCopy = new EMData(End,End);
01074         EMData* fkICopy = new EMData(End,End);
01075         EMData* fkCopy  = new EMData(End,End);
01076 
01077 
01078         for  (int jCount= 0; jCount<End*End; jCount++) {
01079 //              jCount = jy*End + jx;
01080                 int jx             = jCount%End ;
01081                 int jy             = (jCount-jx)/End ;
01082                 jxjyatan2[jCount]  = atan2((float)(jy+1-Mid) , (float)(jx +1-Mid));
01083         }
01084 
01085 
01086         for (int kEx= 0; kEx<2*Mid; kEx=kEx+2) { // kEx twice the value of the Fourier
01087                                                 // x variable: EMAN index for real, imag
01088                 int kx    = kEx/2;              // kx  is  the value of the Fourier variable
01089                 int kIx   = kx+Mid-1; // This is the value of the index for a matlab image (-1)
01090                 int kCx   =  -kx ;
01091                 int kCIx  = kCx+ Mid-1 ;
01092                 for (int kEy= 0 ; kEy<End; kEy++) { // This is the value of the EMAN index
01093                         int kIy              =  kEy       ; //  This is the value of the index for a matlab image (-1)
01094                         int ky               =  kEy+1-Mid; // (kEy+ Mid-1)%End - Mid+1 ;  // This is the actual value of the Fourier variable
01095                         float realVal        =  fk -> get_value_at(kEx  ,kEy) ;
01096                         float imagVal        =  fk -> get_value_at(kEx+1,kEy) ;
01097                         float absVal         =  ::sqrt(realVal*realVal+imagVal*imagVal);
01098                         float fkAng          =  atan2(imagVal,realVal);
01099 
01100                         float NewRealVal   ;
01101                         float NewImagVal   ;
01102                         float AngMatlab    ;
01103 
01104                         if (kIx==Mid-1) {
01105 //                              AngMatlab = -fkAng - 2.*M_PI*(kIy+ 1-Mid)*(Mid)/End;
01106 //                      cout<< "i= " << i << " kIx= " << kIx << " kIy=" << kIy << " fkVecR[i] =" << fkVecR[i]<< " fkVecI[i]="  << fkVecI[i] <<"  angle[i]= "  << AngMatlab << endl;
01107                         }
01108 
01109                         if (kIx>Mid-1){
01110 //                      cout<< "i= " << i << " kIx= " << kIx << " kIy=" << kIy << " fkVecR[i] =" << fkVecR[i]<< " fkVecI[i]="  << fkVecI[i] <<"  angle[i]= "  << AngMatlab << endl;
01111                         }
01112 
01113                         AngMatlab = fkAng - 2.0f*M_PI*(kx +ky)*(Mid)/End;
01114                         NewRealVal  =   absVal*cos(AngMatlab);
01115                         NewImagVal  =   absVal*sin(AngMatlab);
01116 
01117 
01118                         fkVecR[kIy+kIx *End] =  NewRealVal ;
01119                         fkVecR[kIy+kCIx*End] =  NewRealVal ;
01120                         fkVecI[kIy+kIx *End] =  NewImagVal ;
01121                         fkVecI[kIy+kCIx*End] = -NewImagVal ;
01122                         absD1fkVec[kIy + kIx  *End] = absVal;
01123                         absD1fkVec[kIy + kCIx *End] = absVal;
01124                         kVecX[kIy+kIx  *End] =  kx      ;
01125                         kVecX[kIy+kCIx *End] =  kCx    ;
01126                         kVecY[kIy+kIx  *End] =  ky     ;
01127                         kVecY[kIy+kCIx *End] =  ky     ;
01128 //                      printf("kx=%d,ky=%d,tempVal =%f+ i %4.2f \n",kx,ky,realVal,imagVal );
01129 //                      cout << "kx = " << kx << "; ky = "<< ky << "; val is" << realVal<<"+ i "<<imagVal<< endl;
01130 
01131 //                      cout << "kIMx = "<< kIx+1 << "; kIMy = "<< kIy+1 <<"; fkAng*9/ 2pi is " << fkAng*9/2/M_PI<<  endl;
01132 //                      cout << "kIMx = "<< kIx+1 << "; kIMy = "<< kIy+1 <<"; absval is " << absVal<<  "; realval is " << NewRealVal<< "; imagval is " << NewImagVal<< endl;
01133                         fkCopy  -> set_value_at(kIx ,kIy, absVal);
01134                         fkCopy  -> set_value_at(kCIx,kIy, absVal);
01135                         fkRCopy -> set_value_at(kIx, kIy, NewRealVal);
01136                         fkRCopy -> set_value_at(kCIx,kIy, NewRealVal);
01137                         fkICopy -> set_value_at(kIx, kIy, NewImagVal);
01138                         fkICopy -> set_value_at(kCIx,kIy,-NewImagVal);
01139 
01140                 }
01141         }
01142 #ifdef _WIN32
01143         _unlink("fkCopy.???");
01144         _unlink("fk?Copy.???");
01145 #else
01146         system("rm -f fkCopy.???");
01147         system("rm -f fk?Copy.???");
01148 #endif  //_WIN32
01149         fkCopy  -> write_image("fkCopy.img");
01150         fkRCopy -> write_image("fkRCopy.img");
01151         fkICopy -> write_image("fkICopy.img");
01152 
01153         cout << "Starting the sort "<< endl;
01154 
01155         vector< pair<float, int> > absInds;
01156         for(int i  = 0; i < CountxyMax; ++i ) {
01157                 pair<float,int> p;
01158                 p = make_pair(absD1fkVec[i],i); // p = make_pair(rand(),i);
01159                 absInds.push_back( p);
01160         }
01161 
01162         std::sort(absInds.begin(),absInds.end());
01163 
01164         for(int i  = 0; i < CountxyMax; ++i ) {
01165                 pair<float,int> p   ;
01166                 p = absInds[i]         ;
01167                 absD1fkVecSorted[CountxyMax-1-i] =  p.first ;
01168                 SortfkInds[CountxyMax-1-i]       =  p.second ;
01169         }
01170 
01171         cout << "Ending the sort "<< endl;
01172 
01173 // float AngsMat[] ={2.8448, -0.3677,-0.2801,-1.0494,-1.7836,-2.5179, 2.9959, 3.0835,-0.1290,-0.8876,2.1829, 2.2705,1.5011,0.7669,0.0327,-0.7366,-0.6489,2.4215,-1.6029,1.4676,1.5552,0.7859,0.0517,-0.6825,-1.4518,-1.3642,1.7063,-1.7845,1.2859,1.3736,0.6043,-0.1299,-0.8642,-1.6335,-1.5459,1.5247,-1.6546,1.4159,1.5036,0.7342,0,-0.7342,-1.5036,-1.4159,1.6546,-1.5247,1.5459,1.6335,0.8642,0.1299,-0.6043,-1.3736,-1.286,1.7846,-1.7063,1.3642,1.4519,0.6825,-0.0517,-0.7859,-1.5553,-1.4676,1.6029,-2.4216,0.649,0.7366,-0.0327,-0.767,-1.5012,-2.2705,-2.1829,0.8877,0.1291,-3.0836,-2.9959,2.5179,1.7837,1.0495,0.2801,0.3677,-2.8449};
01174 
01175 
01176         for(int i  = 0; i < CountxyMax; ++i ) {  // creates a new fkVec
01177                 int Si  = SortfkInds[i];
01178                 int kIx = (int)  Si/End;  kIx = (int)  i/End; // i = kIx*End+kIy
01179 //              int kIy = Si  - kIx*End;  kIy = i  - kIx*End;
01180 //              int iC = (End-1-kIx)*End + (End-1-kIy);
01181 //              if (i<30) { cout<< "i= " << i << " kIx= " << kIx << " kIy=" << kIy << " valAft=" << absD1fkVecSorted[i]<< " valBef="  <<     absD1fkVec[Si] << "  SortfkInds = " << Si <<endl; }// This worked
01182 //              cout<< "i= " << i << " kIx= " << kIx << " kIy=" << kIy << " fkVecR[i] =" << fkVecR[i]<< " fkVecI[i]="  << fkVecI[i] <<"  angle[i]= "  << fkAng << endl;
01183         }
01184         cout<< "Ratio of Last Amplitude to First Amplitude= " << absD1fkVecSorted[NK] /absD1fkVecSorted[0]  << endl;
01185 
01186 //      pause;
01187 
01188 //      for(int i  = 0; i < NK; ++i ) { // Prints out the new fkVec ,  CountxyMax
01189 //              int Si= SortfkInds[i];
01190 //              int kIx = (int)  Si/End; // i = kIx*End+kIy
01191 //              int kIy = Si  - kIx*End;
01192 //              cout << " kIxM= " << kIx+1 << " kIyM=" << kIy+1 << " fkVecAbs=" << ::sqrt(fkVecR[Si]*fkVecR[Si] +  fkVecI[Si]* fkVecI[Si]) << " fkVecAbs=" << absD1fkVecSorted[i] << " kx= " << kVecX[Si] <<  " ky=" << kVecY[Si] <<  endl;
01193 //      }
01194 
01195 //       angEMAN+angMat+angDiff    =0  mod 2 pi
01196 
01197 //      angDiff=  2*pi*((-4):4)*(Mid)/End; angEMAN+angMat+angDiff= integer*2 *pi
01198 //              [  absD1fkVecSorted, SortfkInds] =sort( absD1fkVec,'descend') ;
01199 //      Util::sort_mat(&absD1fkVec[0],&absD1fkVec[Countxy],&SortfkInds[0],&SortfkInds[Countxy]);
01200 
01201 
01202 //      Let radial sampling be 0:0.5:(Mid-1)
01203 
01204  //     int NK=  min(12,CountxyMax) ;
01205 
01206 
01207 
01208         cout << "NK = " << NK << endl;
01209         float frR= 3.0/4.0;
01210         int LradRange= (int) (floor(Mid/frR)) ;
01211 
01212         float *radRange = new float[LradRange]; //= 0:.75:(Mid-1);
01213         radRange[0]=0;
01214         for (int irad=1; irad < LradRange; irad++){
01215                         radRange[irad] = radRange[irad-1] + frR; }
01216 
01217 
01218 
01219          // should equal to (2*Mid-1)
01220         cout << "Starting the calculation of invariants for N= " << N << endl;
01221 
01222 /*      int NMax=5;            */
01223 
01224         EMData* RotTransInv = new EMData();
01225         RotTransInv -> set_size(LradRange,LradRange);
01226 
01227 
01228 //      float  *RotTransInv       = new float[LradRange*LradRange ] ;
01229 //      float  *RotTransInvN      = new float[LradRange*LradRange*(NMax+1) ] ;
01230 
01231 //      for (int N=0 ; N<NMax; N++) {
01232 
01233         for (int jr1=0; jr1 < LradRange ; jr1++ ) { // LradRange
01234                 float r1= radRange[jr1];
01235 //              cout << "Pre jr2 "<< endl;
01236                 for (int jr2=0;  jr2<LradRange;  jr2++ ) { //LradRange
01237                         float r2= radRange[jr2];
01238                         float RotTransInvTemp=0;
01239                         for (int jCountkxy =0; jCountkxy<NK; jCountkxy++){
01240                                 int Countkxy =SortfkInds[jCountkxy] ;   // 1: CountxyMax
01241                                 int kx = kVecX[Countkxy] ;
01242                                 int ky = kVecY[Countkxy] ;
01243                                 float k2 = (float)(kx*kx+ky*ky);
01244                                 if (k2==0) { continue;}
01245                                 float phiK =0;
01246                                 if (k2>0) phiK= jxjyatan2[ (ky+Mid-1)*End + kx+Mid-1];  phiK= atan2((float)ky,(float)kx);
01247 
01248                                 float fkR     = fkVecR[Countkxy] ;
01249                                 float fkI     = fkVecI[Countkxy]  ;
01250 /*                              printf("jCountkxy=%d, Countkxy=%d,absD1fkVec(Countkxy)=%f,\t\t kx=%d, ky=%d \n", jCountkxy, Countkxy, absD1fkVec[Countkxy], kx, ky);*/
01251 
01252                                 for (int jCountqxy =0; jCountqxy<NK; jCountqxy++){
01253                                         int Countqxy =SortfkInds[jCountqxy] ;   // Countqxy is the index for absD1fkVec
01254                                         int qx   = kVecX[Countqxy] ;
01255                                         int qy   = kVecY[Countqxy] ;
01256                                         int q2   = qx*qx+qy*qy;
01257                                         if (q2==0) {continue;}
01258                                         float phiQ =0;
01259                                         if (q2>0) phiQ = jxjyatan2[ (qy+Mid-1)*End + qx+Mid-1];   phiQ=atan2((float)qy,(float)qx);
01260                                         float fqR     = fkVecR[Countqxy]  ;
01261                                         float fqI     = fkVecI[Countqxy]  ;
01262                                         int kCx  = (-kx-qx);
01263                                         int kCy  = (-ky-qy);
01264                                         int kCIx = ((kCx+Mid+2*End)%End);// labels of the image in C
01265                                         int kCIy = ((kCy+Mid+2*End)%End);
01266                                         kCx  = kCIx-Mid; // correct
01267                                         kCy  = kCIy-Mid; // correct
01268                                         int CountCxy = kCIx*End+kCIy;
01269                                         float fCR     = fkVecR[CountCxy];
01270                                         float fCI     = fkVecI[CountCxy];
01271                                         if (jr1+jr2==-1) {
01272                                         printf("jCountqxy=%d , Countqxy=%d, absD1fkVec(Countqxy)=%f,qx=%d, qy=%d \n", jCountqxy, Countqxy, absD1fkVec[Countqxy],qx, qy);
01273                                         printf(" CountCxy=%d,absD1fkVec[CountCxy]=%f,  kCx=%d,     kCy=%d \n",CountCxy, absD1fkVec[CountCxy], kCx, kCy );
01274                                         }
01275                                         for (int p=0; p<NK; p++){
01276 //                                              printf("p=%d, SortfkInds[p]=%d, CountCxy =%d \n", p,SortfkInds[p], CountCxy);
01277                                                 if (SortfkInds[p]==CountCxy){
01278                                                         float Arg1 = 2.0f*M_PI*r1*::sqrt((float) q2)/End;
01279                                                         float Arg2 = 2.0f*M_PI*r2*::sqrt((float) k2)/End;
01280 //                                                      printf("Arg1=%4.2f, Arg2=%4.2f,  \n",Arg1, Arg2 );
01281 //                                                      if (Arg1+ Arg2<15) {
01282                                                                 float bispectemp  = (fkR*(fqR*fCR -fqI*fCI) -fkI*(fqI*fCR  +fqR*fCI))
01283                                                                 * cos(N*(phiK-phiQ+M_PI));
01284                                                                 bispectemp  -= (fkR*(fqR*fCI + fqI*fCR) +fkI*(fqR*fCR - fqI*fCI))
01285                                                                 * sin(N*(phiK-phiQ+M_PI));
01286                                                                 float bess1 = calc_bessel(N, Arg1 );
01287                                                                 float bess2 = calc_bessel(N, Arg2 );
01288 //                      printf("fkr=%4.2f, fqr=%4.2f, bess1=%4.2f,bess2=%4.2f \n",fkR, fqR, bess1, bess2);
01289 /*                      printf("p =%d, SortfkInds[p]=%d, CountCxy=%d, Arg1 =%4.2f, bess1=%4.2f,  \n",
01290                                 p, SortfkInds[p],CountCxy, Arg1, bess1);*/
01291                                                                 RotTransInvTemp   = RotTransInvTemp  + bispectemp  * bess1*bess2 ;
01292 //                                                      }
01293                                                 }
01294                                         }
01295                                 } // jCountqxy
01296                         } // jCountkxy
01297                         RotTransInv -> set_value_at(jr1,jr2, RotTransInvTemp)   ;
01298 /*              RotTransInvN[jr1 + LradRange*jr2+LradRange*LradRange*N] = RotTransInvTemp  ;*/
01299                 } //jr2
01300         } //jr1
01301 // }//N
01302 
01303         return  RotTransInv ;
01304 
01305 
01306 }
01307 
01308 
01309 
01310 /*
01311 // find example
01312 #include <iostream>
01313 #include <algorithm>
01314 #include <vector>
01315 using namespace std;
01316 
01317 int main () {
01318   int myints[] = { 10, 20, 30 ,40 };
01319   int * p;
01320 
01321   // pointer to array element:
01322   p = find(myints,myints+4,30);
01323   ++p;
01324   cout << "The element following 30 is " << *p << endl;
01325 
01326   vector<int> myvector (myints,myints+4);
01327   vector<int>::iterator it;
01328 
01329   // iterator to vector element:
01330   it = find (myvector.begin(), myvector.end(), 30);
01331   ++it;
01332   cout << "The element following 30 is " << *it << endl;
01333 
01334   return 0;
01335 }*/
01336 
01337 EMData*   EMData::bispecRotTransInvDirect(int type)
01338 {
01339 
01340         int EndP = this -> get_xsize(); // length(fTrueVec);
01341         int Mid  = (int) ((1+EndP)/2);
01342         int End = 2*Mid-1;
01343 
01344         int CountxyMax = End*End;
01345 
01346 //      int   *SortfkInds       = new    int[CountxyMax];
01347         int   *kVecX            = new    int[CountxyMax];
01348         int   *kVecY            = new    int[CountxyMax];
01349         float *fkVecR           = new  float[CountxyMax];
01350         float *fkVecI           = new  float[CountxyMax];
01351         float *absD1fkVec       = new  float[CountxyMax];
01352 //      float *absD1fkVecSorted = new  float[CountxyMax];
01353 
01354 
01355         float *jxjyatan2         = new  float[End*End];
01356 
01357 
01358         EMData * ThisCopy = new EMData(End,End);
01359 
01360         for (int jx=0; jx <End ; jx++) {  // create jxjyatan2
01361                 for (int jy=0; jy <End ; jy++) {
01362                         float ValNow = this -> get_value_at(jx,jy);
01363                         ThisCopy -> set_value_at(jx,jy,ValNow);
01364                         jxjyatan2[jy*End + jx]  = atan2((float)(jy+1-Mid) , (float)(jx +1 -Mid));
01365 //              cout<< " jxM= " << jx+1<<" jyM= " << jy+1<< "ValNow" << ValNow << endl; //    Works
01366         }}
01367 
01368 
01369         EMData* fk = ThisCopy -> do_fft();
01370         fk          ->process_inplace("xform.fourierorigin.tocenter");
01371 
01372 //      Create kVecX , kVecy etc
01373 
01374         for (int kEx= 0; kEx<2*Mid; kEx=kEx+2) { // kEx twice the value of the Fourier
01375                                                 // x variable: EMAN index for real, imag
01376                 int kx    = kEx/2;              // kx  is  the value of the Fourier variable
01377                 int kIx   = kx+Mid-1; // This is the value of the index for a matlab image (-1)
01378                 int kCx   = -kx ;
01379                 int kCIx  = kCx+ Mid-1 ;
01380                 for (int kEy= 0 ; kEy<End; kEy++) { // This is the value of the EMAN index
01381                         int kIy              =  kEy       ; //  This is the value of the index for a matlab image (-1)
01382                         int ky               =  kEy+1-Mid; // (kEy+ Mid-1)%End - Mid+1 ;  // This is the actual value of the Fourier variable
01383                         float realVal        =  fk -> get_value_at(kEx  ,kEy) ;
01384                         float imagVal        =  fk -> get_value_at(kEx+1,kEy) ;
01385                         float absVal         =  ::sqrt(realVal*realVal+imagVal*imagVal);
01386                         float fkAng          =  atan2(imagVal,realVal);
01387 
01388                         float NewRealVal   ;
01389                         float NewImagVal   ;
01390                         float AngMatlab    ;
01391 
01392                         if (kIx==Mid-1) {
01393 //                              AngMatlab = -fkAng - 2.*M_PI*(kIy+ 1-Mid)*(Mid)/End;
01394                         }
01395 
01396                         if (kIx>Mid-1){
01397 //                      cout<< "i= " << i << " kIx= " << kIx << " kIy=" << kIy << " fkVecR[i] =" << fkVecR[i]<< " fkVecI[i]="  << fkVecI[i] <<"  angle[i]= "  << AngMatlab << endl;
01398                         }
01399 
01400                         AngMatlab = fkAng - 2.0f*M_PI*(kx +ky)*(Mid)/End;
01401                         NewRealVal  =   absVal*cos(AngMatlab);
01402                         NewImagVal  =   absVal*sin(AngMatlab);
01403 
01404 
01405                         fkVecR[ kIy +kIx *End] =  NewRealVal ;
01406                         fkVecR[(End-1-kIy)+kCIx*End] =  NewRealVal ;
01407                         fkVecI[ kIy +kIx *End] =  NewImagVal ;
01408                         fkVecI[(End-1-kIy)+kCIx*End] = -NewImagVal ;
01409                         absD1fkVec[(End-1-kIy) + kIx  *End] = absVal;
01410                         absD1fkVec[(End-1-kIy) + kCIx *End] = absVal;
01411                         kVecX[kIy+kIx  *End] =  kx      ;
01412                         kVecX[kIy+kCIx *End] =  kCx    ;
01413                         kVecY[kIy+kIx  *End] =  ky     ;
01414                         kVecY[kIy+kCIx *End] =  ky     ;
01415 
01416  //                     cout << " kIxM= " << kIx+1 << " kIy=" << kIy+1 << " fkVecR[i] =" << NewRealVal << " fkVecI[i]="  << NewImagVal <<"  angle[i]= "  << AngMatlab << " Total Index" << kIy+kIx *End << endl;
01417 
01418 //                      printf("kx=%d,ky=%d,tempVal =%f+ i %4.2f \n",kx,ky,realVal,imagVal );
01419 //                      cout << "kx = " << kx << "; ky = "<< ky << "; val is" << realVal<<"+ i "<<imagVal<< endl;
01420 
01421 //                      cout << "kIMx = "<< kIx+1 << "; kIMy = "<< kIy+1 <<"; fkAng*9/ 2pi is " << fkAng*9/2/M_PI<<  endl;
01422 //                      cout << "kIMx = "<< kIx+1 << "; kIMy = "<< kIy+1 <<"; absval is " << absVal<<  "; realval is " << NewRealVal<< "; imagval is " << NewImagVal<< endl;
01423                 }
01424         }
01425 
01426 
01427 //      for (int TotalInd = 0 ;  TotalInd < CountxyMax ; TotalInd++){
01428 //              int kx     = kVecX[TotalInd]; // This is the value of the index for a matlab image (-1)
01429 //              int kIx    = kx+Mid-1; // This is the value of the index for a matlab image (-1)
01430 //              int ky     = kVecY[TotalInd];
01431 //              int kIy    = ky+Mid-1; // This is the value of the index for a matlab image (-1)
01432                 //float fkR  = fkVecR[kIy+kIx *End]  ;
01433                 //float fkI  = fkVecI[kIy+kIx *End]  ;
01434 //      }
01435 
01436         float frR= 3.0/4.0;
01437         frR= 1;
01438         int LradRange= (int) (1+floor(Mid/frR -.1)) ;
01439 
01440         float *radRange = new float[LradRange]; //= 0:.75:(Mid-1);
01441         for (int irad=0; irad < LradRange; irad++){
01442                         radRange[irad] =  frR*irad;
01443 //                      cout << " irad = " << irad << " radRange[irad]= " <<  radRange[irad] <<  " LradRange= " << LradRange << endl;
01444         }
01445 
01446         cout << "Starting the calculation of invariants" << endl;
01447 
01448 
01449         if (type==0) {
01450                 int LthetaRange  = 59;
01451                 float ftR        = (2.0f*M_PI/LthetaRange );
01452                 float *thetaRange = new float[LthetaRange]; //= 0:.75:(Mid-1);
01453 
01454                 for (int ith=0; ith < LthetaRange; ith++){
01455                                 thetaRange[ith] =  ftR*ith; }
01456 
01457                 int TotalVol = LradRange*LradRange*LthetaRange;
01458 
01459                 float *RotTransInv   = new  float[TotalVol];
01460                 float *WeightInv     = new  float[TotalVol];
01461 
01462                 for (int jW=0; jW<TotalVol; jW++) {
01463                         RotTransInv[jW] = 0;
01464                         WeightInv[jW]   = 0;
01465                 }
01466 
01467                 for (int jW=0; jW<TotalVol; jW++) {
01468                         RotTransInv[jW] = 0;
01469                         WeightInv[jW]   = 0;
01470                 }
01471         //      float  *RotTransInv       = new float[LradRange*LradRange ] ;
01472         //      float  *RotTransInvN      = new float[LradRange*LradRange*(NMax+1) ] ;
01473 
01474                 for (int Countkxy =0; Countkxy<CountxyMax; Countkxy++){  // Main Section for type 0
01475                         int kx = kVecX[Countkxy] ;
01476                         int ky = kVecY[Countkxy] ;
01477                         float k2 = ::sqrt((float)(kx*kx+ky*ky));
01478                         float phiK =0;
01479                         if (k2>0)    phiK = jxjyatan2[ (ky+Mid-1)*End + kx+Mid-1]; //  phiK=atan2(ky,kx);
01480                         float fkR     = fkVecR[(ky+Mid-1) + (kx+Mid-1) *End] ;
01481                         float fkI     = fkVecI[(ky+Mid-1) + (kx+Mid-1) *End]  ;
01482         //              printf("Countkxy=%d,\t kx=%d, ky=%d, fkR=%3.2f,fkI=%3.2f \n", Countkxy, kx, ky, fkR, fkI);
01483 
01484                         if ((k2==0)|| (k2>Mid) ) { continue;}
01485 
01486                         for (int Countqxy =0; Countqxy<CountxyMax; Countqxy++){   // This is the innermost loop
01487                                 int qx   = kVecX[Countqxy] ;
01488                                 int qy   = kVecY[Countqxy] ;
01489                                 float q2   = ::sqrt((float)(qx*qx+qy*qy));
01490                                 if ((q2==0)|| (q2>Mid) ) {continue;}
01491                                 float phiQ =0;
01492                                 if (q2>0)   phiQ = jxjyatan2[ (qy+Mid-1)*End + qx+Mid-1]; // phiQ=atan2(qy,qx);
01493                                 float fqR     = fkVecR[(qy+Mid-1) + (qx+Mid-1) *End] ;
01494                                 float fqI     = fkVecI[(qy+Mid-1) + (qx+Mid-1) *End]  ;
01495                                 int kCx  = (-kx-qx);
01496                                 int kCy  = (-ky-qy);
01497                                 int kCIx = ((kCx+Mid+2*End)%End);// labels of the image in C
01498                                 int kCIy = ((kCy+Mid+2*End)%End);
01499                                 kCx  = ((kCIx+End-1)%End)+1-Mid; // correct
01500                                 kCy  = ((kCIy+End-1)%End)+1-Mid ; // correct
01501 
01502 //                              float C2   = ::sqrt((float)(kCx*kCx+ kCy*kCy));
01503                                 int CountCxy  = (kCx+Mid-1)*End+(kCy+Mid-1);
01504                                 float fCR     = fkVecR[CountCxy];
01505                                 float fCI     = fkVecI[CountCxy];
01506         /*                      if (Countkxy==1) {
01507                                         printf(" Countqxy=%d, absD1fkVec(Countqxy)=%f,qx=%d, qy=%d \n", Countqxy, absD1fkVec[Countqxy],qx, qy);
01508                                         printf(" CountCxy=%d, absD1fkVec[CountCxy]=%f,kCx=%d,kCy=%d \n",CountCxy, absD1fkVec[CountCxy], kCx, kCy );
01509                                 }*/
01510 //                              float   phiC = jxjyatan2[ (kCy+Mid-1)*End + kCx+Mid-1];
01511                                 float   phiQK = (4*M_PI+phiQ-phiK);
01512                                 while (phiQK> (2*M_PI)) phiQK -= (2*M_PI);
01513 
01514 
01515 
01516                                 float bispectemp  = (fkR*(fqR*fCR -fqI*fCI) -fkI*(fqI*fCR  +fqR*fCI));
01517 
01518                                 if  ((q2<k2) )  continue;
01519 //                              if  ((q2<k2) || (C2<k2) || (C2<q2))  continue;
01520 
01521         //                              printf(" CountCxy=%d, absD1fkVec[CountCxy]=%f,kCx=%d,kCy=%d \n",CountCxy, absD1fkVec[CountCxy], kCx, kCy );
01522 
01523         //                      up to here, matched perfectly with Matlab
01524 
01525                                 int k2IndLo  = 0; while ((k2>=radRange[k2IndLo+1]) && (k2IndLo+1 < LradRange ) ) k2IndLo +=1;
01526                                 int k2IndHi = k2IndLo;
01527                                 float k2Lo= radRange[k2IndLo];
01528                                 if (k2IndLo+1< LradRange) {
01529                                         k2IndHi   = k2IndLo+1;
01530                                 }
01531 //                              float k2Hi= radRange[k2IndHi];
01532 
01533                                 float kCof =k2-k2Lo;
01534 
01535                                 int q2IndLo  = 0; while ((q2>=radRange[q2IndLo+1]) && (q2IndLo+1 < LradRange ) ) q2IndLo +=1;
01536                                 int q2IndHi=q2IndLo;
01537                                 float q2Lo= radRange[q2IndLo];
01538                                 if (q2IndLo+1 < LradRange)  {
01539                                         q2IndHi   = q2IndLo+1 ;
01540                                 }
01541                                 float qCof = q2-q2Lo;
01542 
01543                                 if ((qCof<0) || (qCof >1) ) {
01544                                         cout<< "Weird! qCof="<< qCof <<  " q2="<< q2 << " q2IndLo="<< q2IndLo << endl ;
01545                                         int x    ;
01546                                         cin >> x ;
01547                                 }
01548 
01549                                 int thetaIndLo = 0; while ((phiQK>=thetaRange[thetaIndLo+1])&& (thetaIndLo+1<LthetaRange)) thetaIndLo +=1;
01550                                 int thetaIndHi = thetaIndLo;
01551 
01552                                 float thetaLo  = thetaRange[thetaIndLo];
01553                                 float thetaHi = thetaLo;
01554                                 float thetaCof = 0;
01555 
01556                                 if (thetaIndLo+1< LthetaRange) {
01557                                         thetaIndHi = thetaIndLo +1;
01558                                 }else{
01559                                         thetaIndHi=0;
01560                                 }
01561 
01562                                 thetaHi    = thetaRange[thetaIndHi];
01563 
01564                                 if (thetaHi==thetaLo) {
01565                                         thetaCof =0 ;
01566                                 } else {
01567                                         thetaCof   = (phiQK-thetaLo)/(thetaHi-thetaLo);
01568                                 }
01569 
01570                                 if ((thetaCof>2*M_PI)  ) {
01571                                         cout<< "Weird! thetaCof="<< thetaCof <<endl ;
01572                                         thetaCof=0;
01573                                 }
01574 
01575 
01576         //                      if ((thetaIndLo>=58) || (k2IndLo >= LradRange-1) || (q2IndLo >= LradRange-1) ) {
01577 
01578 
01579                                 for (int jk =1; jk<=2; jk++){
01580                                 for (int jq =1; jq<=2; jq++){
01581                                 for (int jtheta =1; jtheta<=2; jtheta++){
01582 
01583                                         float Weight = (kCof+(1-2*kCof)*(jk==1))*(qCof+(1-2*qCof)*(jq==1))
01584                                                         * (thetaCof+(1-2*thetaCof)*(jtheta==1));
01585 
01586 
01587                                         int k2Ind      =  k2IndLo*(jk==1)      +   k2IndHi*(jk==2);
01588                                         int q2Ind      =  q2IndLo*(jq==1)      +   q2IndHi*(jq==2);
01589                                         int thetaInd   =  thetaIndLo*(jtheta==1)  + thetaIndHi*(jtheta ==2);
01590                                         int TotalInd   = thetaInd*LradRange*LradRange+q2Ind*LradRange+k2Ind;
01591         /*                              if (TotalInd+1 >=  LthetaRange*LradRange*LradRange) {
01592                                                 cout << "Weird!!! TotalInd="<< TotalInd << " IndMax" << LthetaRange*LradRange*LradRange << " LradRange=" << LradRange << endl;
01593                                                 cout << "k2Ind= "<< k2Ind  << " q2Ind="<< q2Ind  << " thetaInd="<< thetaInd  << " q2IndLo="<< q2IndLo  << " q2IndHi="<< q2IndHi  <<  endl;
01594                                                 cout << "k2=" << k2 << "q2=" << q2 << " phiQK=" << phiQK*180.0/M_PI<< endl;
01595                                         }*/
01596 
01597                                         RotTransInv[TotalInd] += Weight*bispectemp;
01598                                         WeightInv[TotalInd]   +=  Weight;
01599         //                              cout << "k2Ind= "<< k2Ind  << " q2Ind="<< q2Ind  << "Weight=" << Weight << endl;
01600                                 }}}
01601                         } // Countqxy
01602                 } // Countkxy
01603 
01604                 cout << "Finished Main Section " << endl;
01605 
01606         /*              RotTransInvN[jr1 + LradRange*jr2+LradRange*LradRange*N] = RotTransInvTemp  ;*/
01607 
01608                 cout << " LradRange " <<LradRange <<" LthetaRange " << LthetaRange << endl;
01609                 EMData *RotTransInvF  = new  EMData(LradRange,LradRange,LthetaRange);
01610                 EMData *WeightImage   = new  EMData(LradRange,LradRange,LthetaRange);
01611 
01612         //      cout << "FFFFFFF" << endl;
01613         //
01614         //      RotTransInvF -> set_size(LradRange,LradRange,LthetaRange);
01615         //
01616         //      cout << "GGGG" << endl;
01617 
01618                 for (int jtheta =0; jtheta < LthetaRange; jtheta++){    // write out to RotTransInvF
01619                 for (int jq =0; jq<LradRange; jq++){ // LradRange
01620                 for (int jk =0; jk<LradRange ; jk++){// LradRange
01621         //              cout << "Hi There" << endl;
01622                         int TotalInd   = jtheta*LradRange*LradRange+jq*LradRange+jk;
01623                         float Weight = WeightInv[TotalInd];
01624                         WeightImage    -> set_value_at(jk,jq,jtheta,Weight);
01625                         RotTransInvF   -> set_value_at(jk,jq,jtheta,0);
01626                         if (Weight <= 0) continue;
01627                         RotTransInvF -> set_value_at(jk,jq,jtheta,RotTransInv[TotalInd] / Weight);//  include /Weight
01628                         int newjtheta = (LthetaRange- jtheta)%LthetaRange;
01629                         RotTransInvF -> set_value_at(jq,jk,newjtheta,RotTransInv[TotalInd]/Weight );//  include /Weight
01630                                 }
01631                         }
01632                 }
01633 
01634                 cout << " Almost Done " << endl;
01635 #ifdef  _WIN32
01636                 _unlink("WeightImage.???");
01637 #else
01638                 system("rm -f WeightImage.???");
01639 #endif  //_WIN32
01640                 WeightImage  -> write_image("WeightImage.img");
01641 
01642                 return  RotTransInvF ;
01643         } // Finish type 0
01644 
01645         if (type==1) {
01646                 int TotalVol = LradRange*LradRange;
01647 
01648                 float *RotTransInv   = new  float[TotalVol];
01649                 float *WeightInv     = new  float[TotalVol];
01650 
01651                 for (int jW=0; jW<TotalVol; jW++) {
01652                         RotTransInv[jW] = 0;
01653                         WeightInv[jW]   = 0;
01654                 }
01655 
01656 
01657                 for (int Countkxy =0; Countkxy<CountxyMax; Countkxy++){
01658                         int kx = kVecX[Countkxy] ;
01659                         int ky = kVecY[Countkxy] ;
01660                         float k2 = ::sqrt((float)(kx*kx+ky*ky));
01661                         float fkR     = fkVecR[(ky+Mid-1) + (kx+Mid-1) *End] ;
01662                         float fkI     = fkVecI[(ky+Mid-1) + (kx+Mid-1) *End]  ;
01663         //              printf("Countkxy=%d,\t kx=%d, ky=%d, fkR=%3.2f,fkI=%3.2f \n", Countkxy, kx, ky, fkR, fkI);
01664 
01665                         if ((k2==0)|| (k2>Mid) ) { continue;}
01666 
01667                         for (int Countqxy =0; Countqxy<CountxyMax; Countqxy++){   // This is the innermost loop
01668 
01669 //                      up to here, matched perfectly with Matlab
01670                                 int qx   = kVecX[Countqxy] ;
01671                                 int qy   = kVecY[Countqxy] ;
01672                                 float q2   = ::sqrt((float)(qx*qx+qy*qy));
01673                                 if ((q2==0)|| (q2>Mid) ) {continue;}
01674                                 if  ((q2<k2) )   continue;
01675 
01676                                 float fqR     = fkVecR[(qy+Mid-1) + (qx+Mid-1) *End] ;
01677                                 float fqI     = fkVecI[(qy+Mid-1) + (qx+Mid-1) *End]  ;
01678 
01679                                 int kCx  = (-kx-qx);
01680                                 int kCy  = (-ky-qy);
01681                                 int kCIx = ((kCx+Mid+2*End)%End);// labels of the image in C
01682                                 int kCIy = ((kCy+Mid+2*End)%End);
01683                                 kCx  = ((kCIx+End-1)%End)+1-Mid; // correct
01684                                 kCy  = ((kCIy+End-1)%End)+1-Mid ; // correct
01685 
01686 //                              float C2   = ::sqrt((float)(kCx*kCx+ kCy*kCy));
01687                                 int CountCxy  = (kCx+Mid-1)*End+(kCy+Mid-1);
01688                                 float fCR     = fkVecR[CountCxy];
01689                                 float fCI     = fkVecI[CountCxy];
01690 
01691 
01692                                 float bispectemp  = (fkR*(fqR*fCR -fqI*fCI) -fkI*(fqI*fCR  +fqR*fCI));
01693 
01694 
01695                                 int k2IndLo  = 0; while ((k2>=radRange[k2IndLo+1]) && (k2IndLo+1 < LradRange ) ) k2IndLo +=1;
01696                                 int k2IndHi = k2IndLo;
01697                                 float k2Lo= radRange[k2IndLo];
01698                                 if (k2IndLo+1< LradRange) {
01699                                         k2IndHi   = k2IndLo+1;
01700                                 }
01701 //                              float k2Hi= radRange[k2IndHi];
01702 
01703                                 float kCof =k2-k2Lo;
01704 
01705                                 int q2IndLo  = 0; while ((q2>=radRange[q2IndLo+1]) && (q2IndLo+1 < LradRange ) ) q2IndLo +=1;
01706                                 int q2IndHi=q2IndLo;
01707                                 float q2Lo= radRange[q2IndLo];
01708                                 if (q2IndLo+1 < LradRange)  {
01709                                         q2IndHi   = q2IndLo+1 ;
01710                                 }
01711                                 float qCof = q2-q2Lo;
01712 
01713 
01714                                 for (int jk =1; jk<=2; jk++){
01715                                 for (int jq =1; jq<=2; jq++){
01716 
01717                                         float Weight = (kCof+(1-2*kCof)*(jk==1))*(qCof+(1-2*qCof)*(jq==1));
01718 
01719                                         int k2Ind      =  k2IndLo*(jk==1)      +   k2IndHi*(jk==2);
01720                                         int q2Ind      =  q2IndLo*(jq==1)      +   q2IndHi*(jq==2);
01721                                         int TotalInd   = q2Ind*LradRange+k2Ind;
01722                                         RotTransInv[TotalInd] += Weight*bispectemp;
01723                                         WeightInv[TotalInd]   +=  Weight;
01724         //                              cout << "k2Ind= "<< k2Ind  << " q2Ind="<< q2Ind  << "Weight=" << Weight << endl;
01725                                 }}
01726                         } // Countqxy
01727                 } // Countkxy
01728 
01729 //              cout << "Finished Main Section " << endl;
01730 //              cout << " LradRange " <<LradRange <<  endl;
01731 
01732 
01733                 EMData *RotTransInvF  = new  EMData(LradRange,LradRange);
01734                 EMData *WeightImage   = new  EMData(LradRange,LradRange);
01735 
01736                 for (int jk =0; jk<LradRange ; jk++){// LradRange
01737                 for (int jq =jk; jq<LradRange; jq++){ // LradRange
01738                         int TotalInd      = jq*LradRange+jk;
01739                         int TotalIndBar   = jq*LradRange+jk;
01740                         float Weight = WeightInv[TotalInd] + WeightInv[TotalIndBar];
01741                         if (Weight <=0) continue;
01742                         WeightImage    -> set_value_at(jk,jq,Weight);
01743                         WeightImage    -> set_value_at(jq,jk,Weight);
01744 #ifdef _WIN32
01745                         float ValNow  = pow( (RotTransInv[TotalInd] + RotTransInv[TotalIndBar]) / Weight, 1.0f/3.0f )  ;
01746 #else
01747                         float ValNow  = cbrt( (RotTransInv[TotalInd] + RotTransInv[TotalIndBar]) / Weight )  ;
01748 #endif  //_WIN32
01749                         RotTransInvF -> set_value_at(jk,jq,ValNow);//  include /Weight
01750                         RotTransInvF -> set_value_at(jq,jk,ValNow );//  include /Weight
01751                 }}
01752 
01753 #ifdef  _WIN32
01754                 _unlink("WeightImage.???");
01755 #else
01756                 system("rm -f WeightImage.???");
01757 #endif  //_WIN32
01758                 WeightImage  -> write_image("WeightImage.img");
01759 
01760                 return  RotTransInvF ;
01761         }
01762         return 0;
01763 }
01764 
01765 
01766 void EMData::insert_clip(const EMData * const block, const IntPoint &origin) {
01767         int nx1 = block->get_xsize();
01768         int ny1 = block->get_ysize();
01769         int nz1 = block->get_zsize();
01770 
01771         Region area(origin[0], origin[1], origin[2],nx1, ny1, nz1);
01772 
01773         //make sure the block fits in EMData 
01774         int x0 = (int) area.origin[0];
01775         x0 = x0 < 0 ? 0 : x0;
01776 
01777         int y0 = (int) area.origin[1];
01778         y0 = y0 < 0 ? 0 : y0;
01779 
01780         int z0 = (int) area.origin[2];
01781         z0 = z0 < 0 ? 0 : z0;
01782 
01783         int x1 = (int) (area.origin[0] + area.size[0]);
01784         x1 = x1 > nx ? nx : x1;
01785 
01786         int y1 = (int) (area.origin[1] + area.size[1]);
01787         y1 = y1 > ny ? ny : y1;
01788 
01789         int z1 = (int) (area.origin[2] + area.size[2]);
01790         z1 = z1 > nz ? nz : z1;
01791         if (z1 <= 0) {
01792                 z1 = 1;
01793         }
01794 
01795         int xd0 = (int) (area.origin[0] < 0 ? -area.origin[0] : 0);
01796         int yd0 = (int) (area.origin[1] < 0 ? -area.origin[1] : 0);
01797         int zd0 = (int) (area.origin[2] < 0 ? -area.origin[2] : 0);
01798 
01799         if (x1 < x0 || y1 < y0 || z1 < z0) return; // out of bounds, this is fine, nothing happens
01800 
01801         size_t clipped_row_size = (x1-x0) * sizeof(float);
01802         size_t src_secsize =  (size_t)(nx1 * ny1);
01803         size_t dst_secsize = (size_t)(nx * ny);
01804 
01805 /*
01806 #ifdef EMAN2_USING_CUDA
01807         if(block->cudarwdata){
01808                 // this is VERY slow.....
01809                 float *cudasrc = block->cudarwdata + zd0 * src_secsize + yd0 * nx1 + xd0;
01810                 if(!cudarwdata) rw_alloc();
01811                 float *cudadst = cudarwdata + z0 * dst_secsize + y0 * nx + x0;
01812                 for (int i = z0; i < z1; i++) {
01813                         for (int j = y0; j < y1; j++) {
01814                                 //printf("%x %x %d\n", cudadst, cudasrc, clipped_row_size);
01815                                 cudaMemcpy(cudadst,cudasrc,clipped_row_size,cudaMemcpyDeviceToDevice);
01816                                 cudasrc += nx1;
01817                                 cudadst += nx;
01818                         }
01819                         cudasrc += src_gap;
01820                         cudadst += dst_gap;
01821                 }
01822                 return;
01823         }
01824 #endif
01825 */
01826         float *src = block->get_data() + zd0 * src_secsize + yd0 * nx1 + xd0;
01827         float *dst = get_data() + z0 * dst_secsize + y0 * nx + x0;
01828         
01829         size_t src_gap = src_secsize - (y1-y0) * nx1;
01830         size_t dst_gap = dst_secsize - (y1-y0) * nx;
01831         
01832         for (int i = z0; i < z1; i++) {
01833                 for (int j = y0; j < y1; j++) {
01834                         EMUtil::em_memcpy(dst, src, clipped_row_size);
01835                         src += nx1;
01836                         dst += nx;
01837                 }
01838                 src += src_gap;
01839                 dst += dst_gap;
01840         }
01841         
01842 #ifdef EMAN2_USING_CUDA 
01843         if(block->cudarwdata){
01844                 copy_to_cuda(); // copy back to device as padding is faster on the host
01845         }
01846 #endif
01847 
01848         update();
01849         EXITFUNC;
01850 }
01851 
01852 
01853 void EMData::insert_scaled_sum(EMData *block, const FloatPoint &center,
01854                                                    float scale, float)
01855 {
01856         ENTERFUNC;
01857         float * data = get_data();
01858         if (get_ndim()==3) {
01859                 // Start by determining the region to operate on
01860                 int xs=(int)floor(block->get_xsize()*scale/2.0);
01861                 int ys=(int)floor(block->get_ysize()*scale/2.0);
01862                 int zs=(int)floor(block->get_zsize()*scale/2.0);
01863                 int x0=(int)center[0]-xs;
01864                 int x1=(int)center[0]+xs;
01865                 int y0=(int)center[1]-ys;
01866                 int y1=(int)center[1]+ys;
01867                 int z0=(int)center[2]-zs;
01868                 int z1=(int)center[2]+zs;
01869 
01870                 if (x1<0||y1<0||z1<0||x0>get_xsize()||y0>get_ysize()||z0>get_zsize()) return;   // object is completely outside the target volume
01871 
01872                 // make sure we stay inside the volume
01873                 if (x0<0) x0=0;
01874                 if (y0<0) y0=0;
01875                 if (z0<0) z0=0;
01876                 if (x1>=get_xsize()) x1=get_xsize()-1;
01877                 if (y1>=get_ysize()) y1=get_ysize()-1;
01878                 if (z1>=get_zsize()) z1=get_zsize()-1;
01879 
01880                 float bx=block->get_xsize()/2.0f;
01881                 float by=block->get_ysize()/2.0f;
01882                 float bz=block->get_zsize()/2.0f;
01883 
01884                 size_t idx;
01885                 for (int x=x0; x<=x1; x++) {
01886                         for (int y=y0; y<=y1; y++) {
01887                                 for (int z=z0; z<=z1; z++) {
01888                                         idx = x + y * nx + (size_t)z * nx * ny;
01889                                         data[idx] +=
01890                                                 block->sget_value_at_interp((x-center[0])/scale+bx,(y-center[1])/scale+by,(z-center[2])/scale+bz);
01891                                 }
01892                         }
01893                 }
01894                 update();
01895         }
01896         else if (get_ndim()==2) {
01897                 // Start by determining the region to operate on
01898                 int xs=(int)floor(block->get_xsize()*scale/2.0);
01899                 int ys=(int)floor(block->get_ysize()*scale/2.0);
01900                 int x0=(int)center[0]-xs;
01901                 int x1=(int)center[0]+xs;
01902                 int y0=(int)center[1]-ys;
01903                 int y1=(int)center[1]+ys;
01904 
01905                 if (x1<0||y1<0||x0>get_xsize()||y0>get_ysize()) return; // object is completely outside the target volume
01906 
01907                 // make sure we stay inside the volume
01908                 if (x0<0) x0=0;
01909                 if (y0<0) y0=0;
01910                 if (x1>=get_xsize()) x1=get_xsize()-1;
01911                 if (y1>=get_ysize()) y1=get_ysize()-1;
01912 
01913                 float bx=block->get_xsize()/2.0f;
01914                 float by=block->get_ysize()/2.0f;
01915 
01916                 for (int x=x0; x<=x1; x++) {
01917                         for (int y=y0; y<=y1; y++) {
01918                                 data[x + y * nx] +=
01919                                         block->sget_value_at_interp((x-center[0])/scale+bx,(y-center[1])/scale+by);
01920                         }
01921                 }
01922                 update();
01923         }
01924         else {
01925                 LOGERR("insert_scaled_sum supports only 2D and 3D data");
01926                 throw ImageDimensionException("2D/3D only");
01927         }
01928 
01929         EXITFUNC;
01930 }
01931 //                      else if ( m == 0 )
01932 //                      {
01933 //                              if ( n_f == -ny/2 )
01934 //                              {
01935 //                                      t2++;
01936 // //                                   continue;
01937 //                                      for (int y = 0; y < return_slice->get_ysize(); ++y) {
01938 //                                              for (int x = 0; x < return_slice->get_xsize(); ++x) {
01939 //                                                      double cur_val = return_slice->get_value_at(x,y);
01940 //                                                      return_slice->set_value_at(x,y,cur_val+dat[idx]*std::pow(-1.0f,y));
01941 //                                              }
01942 //                                      }
01943 //                                      if (phase > 0.01 ) cout << "foo 2 " << phase << " " << amp << " " << dat[idx] << endl;
01944 //                              }
01945 //                              else
01946 //                              {
01947 //                                      if ( n_f < 1 ) continue;
01948 //                                      t3++;
01949 //                                      for (int y = 0; y < return_slice->get_ysize(); ++y) {
01950 //                                              for (int x = 0; x < return_slice->get_xsize(); ++x) {
01951 //                                                      double cur_val = return_slice->get_value_at(x,y);
01952 //                                                      return_slice->set_value_at(x,y,cur_val+2*amp*cos(ndash*y+phase));
01953 //                                              }
01954 //                                      }
01955 //                              }
01956 //                      }
01957 //                      else if ( n_f == -ny/2 )
01958 //                      {
01959 //                              if ( m == ((nx-2)/2) )
01960 //                              {
01961 //                                      t4++;
01962 //                                      for (int y = 0; y < return_slice->get_ysize(); ++y) {
01963 //                                              for (int x = 0; x < return_slice->get_xsize(); ++x) {
01964 //                                                      double cur_val = return_slice->get_value_at(x,y);
01965 //                                                      return_slice->set_value_at(x,y,cur_val+dat[idx]*std::pow(-1.0f,x+y));
01966 //                                              }
01967 //                                      }
01968 //                                      if (phase > 0.01 ) cout << "foo 4 " << phase << " " << amp << " " << dat[idx] << endl;
01969 //                              }
01970 //                              else
01971 //                              {
01972 //                                      t5++;
01973 //                                      for (int y = 0; y < return_slice->get_ysize(); ++y) {
01974 //                                              for (int x = 0; x < return_slice->get_xsize(); ++x) {
01975 //                                                      double cur_val = return_slice->get_value_at(x,y);
01976 //                                                      return_slice->set_value_at(x,y,cur_val+2*amp*cos(mdash*x+phase));
01977 //                                              }
01978 //                                      }
01979 //                              }
01980 //                      }
01981 //                      else if ( n_f == 0 )
01982 //                      {
01983 //                              if ( m == ((nx-2)/2) )
01984 //                              {
01985 //                                      t6++;
01986 //                                      for (int y = 0; y < return_slice->get_ysize(); ++y) {
01987 //                                              for (int x = 0; x < return_slice->get_xsize(); ++x) {
01988 //                                                      double cur_val = return_slice->get_value_at(x,y);
01989 //                                                      return_slice->set_value_at(x,y,cur_val+dat[idx]*std::pow(-1.0f,x));
01990 //                                              }
01991 //                                      }
01992 //                                      if (phase > 0.01 ) cout << "foo 3 " << phase << " " << amp << " " << dat[idx] << endl;
01993 //                              }
01994 //                              else
01995 //                              {
01996 //                                      t7++;
01997 //                                      for (int y = 0; y < return_slice->get_ysize(); ++y) {
01998 //                                              for (int x = 0; x < return_slice->get_xsize(); ++x) {
01999 //                                                      double cur_val = return_slice->get_value_at(x,y);
02000 //                                                      return_slice->set_value_at(x,y,cur_val+2*amp*cos(mdash*x+M_PI*y + phase));
02001 //                                              }
02002 //                                      }
02003 //                              }
02004 //                      }
02005 //                      else if ( m == ((nx-2)/2) )
02006 //                      {
02007 //                              if ( n_f < 1 ) continue;
02008 //                              t8++;
02009 //                              for (int y = 0; y < return_slice->get_ysize(); ++y) {
02010 //                                      for (int x = 0; x < return_slice->get_xsize(); ++x) {
02011 //                                              double cur_val = return_slice->get_value_at(x,y);
02012 //                                              return_slice->set_value_at(x,y,cur_val+2*amp*cos(ndash*y+M_PI*x+phase));
02013 //                                      }
02014 //                              }
02015 //                      }
02016 // }