Xmipp  v3.23.11-Nereus
angular_projection_matching.h
Go to the documentation of this file.
1 /***************************************************************************
2  *
3  * Authors: Sjors Scheres scheres@cnb.csic.es (2004)
4  * Authors: Roberto Marabini roberto@cnb.csic.es (2012)
5  *
6  * Unidad de Bioinformatica of Centro Nacional de Biotecnologia , CSIC
7  *
8  * This program is free software; you can redistribute it and/or modify
9  * it under the terms of the GNU General Public License as published by
10  * the Free Software Foundation; either version 2 of the License, or
11  * (at your option) any later version.
12  *
13  * This program is distributed in the hope that it will be useful,
14  * but WITHOUT ANY WARRANTY; without even the implied warranty of
15  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16  * GNU General Public License for more details.
17  *
18  * You should have received a copy of the GNU General Public License
19  * along with this program; if not, write to the Free Software
20  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA
21  * 02111-1307 USA
22  *
23  * All comments concerning this program package may be sent to the
24  * e-mail address 'xmipp@cnb.csic.es'
25  ***************************************************************************/
26 
27 #ifndef _angular_projection_matching_H
28 #define _angular_projection_matching_H
29 
30 #include "core/metadata_db.h"
31 #include "core/multidim_array.h"
32 #include "core/xmipp_program.h"
33 #include "core/xmipp_threads.h"
34 #include "data/polar.h"
35 #include "data/sampling.h"
36 
37 template<typename T>
38 class Image;
39 
40 constexpr int MY_OUPUT_SIZE = 10;
41 
43 
44 // Thread declaration
45 void * threadRotationallyAlignOneImage( void * data );
46 
47 // This structure is needed to pass parameters to threadRotationallyAlignOneImage
48 typedef struct{
49  size_t thread_id;
52  size_t this_image;
53  int * opt_refno;
54  double * opt_psi;
55  bool * opt_flip;
56  double * maxcorr;
59 
65 {
66 
67 public:
68 
76  size_t dim, paddim;
78  double pad;
80  double max_shift;
82  int Ri, Ro;
84  double avail_memory;
97  std::vector<int> pointer_allrefs2refsinmem; //order after removing redundant
98  std::vector<int> pointer_refsinmem2allrefs; //order in memory
100  std::vector <size_t> convert_refno_to_stack_position;
102  std::vector<size_t> ids;
121  std::vector<int> search5d_xoff, search5d_yoff;
127  int threads;
128  //Numbre of possible orientations
131  size_t nr_trans;
134 
136  bool do_scale;
138  double scale_step;
139  double scale_nsteps;
140 
143 public:
144 
146  virtual void readParams();
147 
149  virtual void defineParams();
150 
152  void show();
153 
155  void run();
156 
158  virtual void produceSideInfo();
159 
163  void rotationallyAlignOneImage(Matrix2D<double> &img, int imgno, int &opt_samplenr,
164  double &opt_psi, bool &opt_flip, double &maxcorr);
165 
170  const int &samplenr, const double &psi, const bool &opt_flip,
171  double &opt_xoff, double &opt_yoff, double &maxcorr);
172 
177  const int &samplenr, const double &psi, const bool &opt_flip,
178  const double &opt_xoff, const double &opt_yoff, const double &old_scale, double &opt_scale, double &maxcorr);
179 
183  void getCurrentReference(int refno, Polar_fftw_plans &local_plans);
184 
191  virtual void processAllImages();
192 
194  void processSomeImages(const std::vector<size_t> &imagesToProcess);
195 
198  void getCurrentImage(size_t imgid, Image<double> &img);
199 
203  virtual void writeOutputFiles();
204 
206  void destroyAndClean();
207 };
209 #endif
void scaleAlignOneImage(MultidimArray< double > &img, const int &samplenr, const double &psi, const bool &opt_flip, const double &opt_xoff, const double &opt_yoff, const double &old_scale, double &opt_scale, double &maxcorr)
void getCurrentReference(int refno, Polar_fftw_plans &local_plans)
virtual void readParams()
Read arguments from command line.
void processSomeImages(const std::vector< size_t > &imagesToProcess)
virtual void defineParams()
Define arguments accepted.
Polar< std::complex< double > > * fP_ref
void * threadRotationallyAlignOneImage(void *data)
void rotationallyAlignOneImage(Matrix2D< double > &img, int imgno, int &opt_samplenr, double &opt_psi, bool &opt_flip, double &maxcorr)
constexpr int MY_OUPUT_SIZE
void translationallyAlignOneImage(MultidimArray< double > &img, const int &samplenr, const double &psi, const bool &opt_flip, double &opt_xoff, double &opt_yoff, double &maxcorr)
Polar< std::complex< double > > * fPm_img
void getCurrentImage(size_t imgid, Image< double > &img)
Definition: polar.h:67
Polar< std::complex< double > > * fP_img
std::vector< size_t > convert_refno_to_stack_position
double psi(const double x)
WriteModeMetaData