Xmipp  v3.23.11-Nereus
fourier_projection.h
Go to the documentation of this file.
1 /***************************************************************************
2  *
3  * Authors: Carlos Oscar S. Sorzano (coss@cnb.csic.es)
4  *
5  * Unidad de Bioinformatica of Centro Nacional de Biotecnologia , CSIC
6  *
7  * This program is free software; you can redistribute it and/or modify
8  * it under the terms of the GNU General Public License as published by
9  * the Free Software Foundation; either version 2 of the License, or
10  * (at your option) any later version.
11  *
12  * This program is distributed in the hope that it will be useful,
13  * but WITHOUT ANY WARRANTY; without even the implied warranty of
14  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15  * GNU General Public License for more details.
16  *
17  * You should have received a copy of the GNU General Public License
18  * along with this program; if not, write to the Free Software
19  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA
20  * 02111-1307 USAcd ..
21  * All comments concerning this program package may be sent to the
22  * e-mail address 'xmipp@cnb.csic.es'
23  ***************************************************************************/
24 
25 #ifndef _CORE_FOURIER_FOURIER_PROJECTION_H
26 #define _CORE_FOURIER_FOURIER_PROJECTION_H
27 
28 #include "core/matrix2d.h"
29 #include "core/xmipp_fftw.h"
30 #include "core/xmipp_image.h"
31 
35 
38 
39 
56 class Projection: public Image<double>
57 {
58 public:
63 
72 
80 
87  void reset(int Ydim, int Xdim);
88 
95  void setAngles(double _rot, double _tilt, double _psi);
96 
102  void read(const FileName& fn, const bool only_apply_shifts = false,
103  DataMode datamode = DATA, MDRow * row = nullptr);
104 
107  void assign(const Projection& P);
108 };
109 
112 {
113 public:
117  double maxFrequency;
119  double BSplineDeg;
120 
121 public:
122  // Auxiliary FFT transformer
124 
125  // Volume to project
127 
128  // Real and imaginary B-spline coefficients for Fourier of the volume
131 
132  // Projection in Fourier space
134 
135  // Projection in real space
137 
138  // Phase shift image
141 
142  // Original volume size
144 
145  // Volume padded size
147 
148  // Euler matrix
150 public:
151  /* Empty constructor */
152  FourierProjector(double paddFactor, double maxFreq, int degree);
153 
154  /*
155  * The constructor of the class
156  */
157  FourierProjector(MultidimArray<double> &V, double paddFactor, double maxFreq, int BSplinedegree);
158 
162  void project(double rot, double tilt, double psi, const MultidimArray<double> *ctf=nullptr);
163 
165  void updateVolume(MultidimArray<double> &V);
166 public:
168  void produceSideInfo();
169 
171  void produceSideInfoProjection();
172 };
173 
174 /*
175  * This function gets an object form the FourierProjection class and makes the desired projection in Fourier space
176  */
177 void projectVolume(FourierProjector &projector, Projection &P, int Ydim, int Xdim,
178  double rot, double tilt, double psi, const MultidimArray<double> *ctf=nullptr);
179 
181 
182 #endif
Matrix2D< double > eulert
void reset(int Ydim, int Xdim)
double psi(const size_t n=0) const
Matrix2D< double > euler
MultidimArray< double > phaseShiftImgB
MultidimArray< double > phaseShiftImgA
double rot(const size_t n=0) const
DataMode
double tilt(const size_t n=0) const
MultidimArray< double > VfourierRealCoefs
double maxFrequency
Maximum Frequency for pixels.
FourierTransformer transformer2D
MultidimArray< double > VfourierImagCoefs
double BSplineDeg
The order of B-Spline for interpolation.
void projectVolume(FourierProjector &projector, Projection &P, int Ydim, int Xdim, double rot, double tilt, double psi, const MultidimArray< double > *ctf=nullptr)
Matrix2D< double > E
Matrix1D< double > direction
MultidimArray< std::complex< double > > projectionFourier
Image< double > projection
void read(const FileName &fn, const bool only_apply_shifts=false, DataMode datamode=DATA, MDRow *row=nullptr)
MultidimArray< double > * volume
void setAngles(double _rot, double _tilt, double _psi)
void assign(const Projection &P)
double paddingFactor
Padding factor.