Xmipp  v3.23.11-Nereus
Functions
xmipp_fft.cpp File Reference
#include "bilib/dft.h"
#include "xmipp_fft.h"
#include "args.h"
#include "histogram.h"
Include dependency graph for xmipp_fft.cpp:

Go to the source code of this file.

Functions

void Whole2Half (const MultidimArray< std::complex< double > > &in, MultidimArray< std::complex< double > > &out)
 
void Half2Whole (const MultidimArray< std::complex< double > > &in, MultidimArray< std::complex< double > > &out, size_t oridim)
 
void Complex2RealImag (const MultidimArray< std::complex< double > > &in, MultidimArray< double > &real, MultidimArray< double > &imag)
 
void RealImag2Complex (const MultidimArray< double > &real, const MultidimArray< double > &imag, MultidimArray< std::complex< double > > &out)
 
void FourierTransform (const MultidimArray< double > &in, MultidimArray< std::complex< double > > &out)
 
void InverseFourierTransform (const MultidimArray< std::complex< double > > &in, MultidimArray< double > &out)
 
void FourierTransformHalf (const MultidimArray< double > &in, MultidimArray< std::complex< double > > &out)
 
void InverseFourierTransformHalf (const MultidimArray< std::complex< double > > &in, MultidimArray< double > &out, int oridim)
 
void FourierTransform (const MultidimArray< std::complex< double > > &in, MultidimArray< std::complex< double > > &out)
 
void InverseFourierTransform (const MultidimArray< std::complex< double > > &in, MultidimArray< std::complex< double > > &out)
 
void FourierTransformHalf (const MultidimArray< std::complex< double > > &in, MultidimArray< std::complex< double > > &out)
 
void InverseFourierTransformHalf (const MultidimArray< std::complex< double > > &in, MultidimArray< std::complex< double > > &out, int orixdim)
 
void centerFFT2 (MultidimArray< double > &v)
 
void ShiftFFT (MultidimArray< std::complex< double > > &v, double xshift)
 
void ShiftFFT (MultidimArray< std::complex< double > > &v, double xshift, double yshift)
 
void ShiftFFT (MultidimArray< std::complex< double > > &v, double xshift, double yshift, double zshift)
 
void CenterOriginFFT (MultidimArray< std::complex< double > > &v, bool forward)
 
template<typename T >
void xmipp2PSD (const MultidimArray< T > &input, MultidimArray< T > &output, bool takeLog)
 
template void xmipp2PSD< float > (const MultidimArray< float > &, MultidimArray< float > &, bool)
 
template void xmipp2PSD< double > (const MultidimArray< double > &, MultidimArray< double > &, bool)
 
void xmipp2CTF (const MultidimArray< double > &input, MultidimArray< double > &output)
 

Function Documentation

◆ FourierTransform()

void FourierTransform ( const MultidimArray< std::complex< double > > &  in,
MultidimArray< std::complex< double > > &  out 
)

Complex Direct Fourier Transform 1D ----------------------------------------—

Definition at line 210 of file xmipp_fft.cpp.

212 {
213  // Only implemented for 1D transforms
214  in.checkDimension(1);
215 
216  int N = XSIZE(in);
217  MultidimArray<double> re(N), tmpre(N), tmpim(N), im(N), cas(N);
218  out.resizeNoCopy(N);
219 
220  GetCaS(MULTIDIM_ARRAY(cas), N);
222  MULTIDIM_ARRAY(im), N);
224  MULTIDIM_ARRAY(tmpre), MULTIDIM_ARRAY(tmpim),
225  MULTIDIM_ARRAY(cas), N);
227  MULTIDIM_ARRAY(out), N);
228 }
void resizeNoCopy(const MultidimArray< T1 > &v)
#define MULTIDIM_ARRAY(v)
int GetCaS(double CaS[], long SignalLength)
void RealImag2Complex(const MultidimArray< double > &real, const MultidimArray< double > &imag, MultidimArray< std::complex< double > > &out)
Definition: xmipp_fft.cpp:114
#define XSIZE(v)
int DftRealImaginaryToRealImaginary(double Re2Re[], double Im2Im[], double *TmpRe, double *TmpIm, double CaS[], long SignalLength)
void Complex2RealImag(const MultidimArray< std::complex< double > > &in, MultidimArray< double > &real, MultidimArray< double > &imag)
Definition: xmipp_fft.cpp:103

◆ FourierTransformHalf()

void FourierTransformHalf ( const MultidimArray< std::complex< double > > &  in,
MultidimArray< std::complex< double > > &  out 
)

Complex Fourier Transform 1D, output half of (centro-symmetric) transform -—

Definition at line 252 of file xmipp_fft.cpp.

254 {
255  // Only implemented for 1D transforms
256  in.checkDimension(1);
257 
259  FourierTransform(in, aux);
260  Whole2Half(aux, out);
261 }
void Whole2Half(const MultidimArray< std::complex< double > > &in, MultidimArray< std::complex< double > > &out)
Definition: xmipp_fft.cpp:34
void FourierTransform(const MultidimArray< double > &in, MultidimArray< std::complex< double > > &out)
Definition: xmipp_fft.cpp:124

◆ InverseFourierTransform()

void InverseFourierTransform ( const MultidimArray< std::complex< double > > &  in,
MultidimArray< std::complex< double > > &  out 
)

Complex Inverse Fourier Transform 1D. --------------------------------------—

Definition at line 231 of file xmipp_fft.cpp.

233 {
234  // Only implemented for 1D transforms
235  in.checkDimension(1);
236 
237  int N = XSIZE(in);
238  MultidimArray<double> tmpre(N), tmpim(N), re(N), im(N), cas(N);
239  out.resizeNoCopy(N);
240 
241  GetCaS(MULTIDIM_ARRAY(cas), N);
243  MULTIDIM_ARRAY(im), N);
245  MULTIDIM_ARRAY(tmpre), MULTIDIM_ARRAY(tmpim),
246  MULTIDIM_ARRAY(cas), N);
248  MULTIDIM_ARRAY(out), N);
249 }
void resizeNoCopy(const MultidimArray< T1 > &v)
int InvDftRealImaginaryToRealImaginary(double Re2Re[], double Im2Im[], double *TmpRe, double *TmpIm, double CaS[], long SignalLength)
#define MULTIDIM_ARRAY(v)
int GetCaS(double CaS[], long SignalLength)
void RealImag2Complex(const MultidimArray< double > &real, const MultidimArray< double > &imag, MultidimArray< std::complex< double > > &out)
Definition: xmipp_fft.cpp:114
#define XSIZE(v)
void Complex2RealImag(const MultidimArray< std::complex< double > > &in, MultidimArray< double > &real, MultidimArray< double > &imag)
Definition: xmipp_fft.cpp:103

◆ InverseFourierTransformHalf()

void InverseFourierTransformHalf ( const MultidimArray< std::complex< double > > &  in,
MultidimArray< std::complex< double > > &  out,
int  orixdim 
)

Complex Inverse Fourier Transform 1D, input half of (centro-symmetric) transform -—

Definition at line 264 of file xmipp_fft.cpp.

266 {
267  // Only implemented for 1D transforms
268  in.checkDimension(1);
269 
271  Half2Whole(in, aux, orixdim);
272  InverseFourierTransform(aux, out);
273  out.setXmippOrigin();
274 }
void InverseFourierTransform(const MultidimArray< std::complex< double > > &in, MultidimArray< double > &out)
Definition: xmipp_fft.cpp:155
void Half2Whole(const MultidimArray< std::complex< double > > &in, MultidimArray< std::complex< double > > &out, size_t oridim)
Definition: xmipp_fft.cpp:66

◆ xmipp2PSD< double >()

template void xmipp2PSD< double > ( const MultidimArray< double > &  ,
MultidimArray< double > &  ,
bool   
)

◆ xmipp2PSD< float >()

template void xmipp2PSD< float > ( const MultidimArray< float > &  ,
MultidimArray< float > &  ,
bool   
)