Xmipp  v3.23.11-Nereus
polar_rotation_estimator.cpp
Go to the documentation of this file.
1 /***************************************************************************
2  *
3  * Authors: David Strelak (davidstrelak@gmail.com)
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 USA
21  *
22  * All comments concerning this program package may be sent to the
23  * e-mail address 'xmipp@cnb.csic.es'
24  ***************************************************************************/
25 
27 
28 namespace Alignment {
29 
30 template<typename T>
31 void PolarRotationEstimator<T>::init2D() {
32  release();
33  auto s = this->getSettings();
34  if (1 != s.hw.size()) {
35  REPORT_ERROR(ERR_ARG_INCORRECT, "Only one CPU thread expected");
36  }
37  try {
38  m_cpu = dynamic_cast<CPU*>(s.hw.at(0));
39  } catch (std::bad_cast&) {
40  REPORT_ERROR(ERR_ARG_INCORRECT, "Instance of CPU expected");
41  }
42 
43  if (std::is_same<T, float>()) {
44  m_dataAux.resize(s.refDims.y(), s.refDims.x());
45  }
46 }
47 
48 template<typename T>
49 void PolarRotationEstimator<T>::load2DReferenceOneToN(const T *ref) {
50  auto isReady = this->isInitialized();
51  if ( ! isReady) {
52  REPORT_ERROR(ERR_LOGIC_ERROR, "Not ready to load a reference signal");
53  }
54 
55  MultidimArray<double> tmp = convert(const_cast<T*>(ref));
56  tmp.setXmippOrigin();
57  auto s = this->getSettings();
58  polarFourierTransform<false>(tmp, m_refPolarFourierI, false,
59  s.firstRing, s.lastRing, m_refPlans, 1);
60  m_rotCorrAux.resize(2 * m_refPolarFourierI.getSampleNoOuterRing() - 1);
61  m_aux.local_transformer.setReal(m_rotCorrAux);
62 }
63 
64 template<>
65 MultidimArray<double> PolarRotationEstimator<float>::convert(float *data) {
66  const size_t s = this->getSettings().otherDims.sizeSingle();
67  for (size_t i = 0; i < s; ++i) {
68  m_dataAux.data[i] = data[i];
69  }
70  return m_dataAux;
71 }
72 
73 template<>
74 MultidimArray<double> PolarRotationEstimator<double>::convert(double *data) {
75  const auto s = this->getSettings().otherDims.copyForN(1);
76  return MultidimArray<double>(
77  s.n(), s.z(),
78  s.y(), s.x(),
79  data);
80 }
81 
82 template<typename T>
83 void PolarRotationEstimator<T>::computeRotation2DOneToN(T *others) {
84  bool isReady = this->isInitialized() && this->isRefLoaded();
85 
86  if ( ! isReady) {
87  REPORT_ERROR(ERR_LOGIC_ERROR, "Not ready to execute. Call init() and load reference");
88  }
89  auto s = this->getSettings();
90  for (size_t n = 0; n < s.otherDims.n(); ++n) {
91  size_t offset = n * s.otherDims.sizeSingle();
92  MultidimArray<double> tmp = convert(others + offset);
93  tmp.setXmippOrigin();
94  polarFourierTransform<false>(tmp, m_polarFourierI, true,
95  s.firstRing, s.lastRing, m_plans, 1);
96  this->getRotations2D().emplace_back(
97  best_rotation(m_refPolarFourierI, m_polarFourierI, m_aux));
98  }
99 }
100 
101 template<typename T>
102 void PolarRotationEstimator<T>::release() {
103  delete m_plans;
104  delete m_refPlans;
105  m_rotCorrAux.clear();
106  m_dataAux.clear();
107 
108  setDefault();
109 }
110 
111 template<typename T>
112 void PolarRotationEstimator<T>::setDefault() {
113  m_plans = nullptr;
114  m_refPlans = nullptr;
115  m_polarFourierI = Polar<std::complex<double>>();
116  m_refPolarFourierI = Polar<std::complex<double>>();
117 }
118 
119 template<typename T>
120 void PolarRotationEstimator<T>::check() {
121  const auto s = this->getSettings();
122  if (s.batch != 1) {
123  std::cerr << "Batch processing is not supported. Signals will be processed one by one.\n";
124  }
125  if (s.allowDataOverwrite) {
126  std::cerr << "allowDataOverwrite flag is ignored, as it's not yet supported.\n";
127  }
128  if (s.refDims.x() != s.refDims.y()) {
129  REPORT_ERROR(ERR_ARG_INCORRECT, "This estimator can work only with square signal");
130  }
131  if (s.refDims.isPadded()) {
132  REPORT_ERROR(ERR_ARG_INCORRECT, "Padded signal is not supported");
133  }
134  if (s.refDims.x() < 6) {
135  // we need some edge around the biggest ring, to avoid invalid memory access
136  REPORT_ERROR(ERR_ARG_INCORRECT, "The input signal is too small.");
137  }
138 }
139 
140 // explicit instantiation
141 template class PolarRotationEstimator<float>;
142 template class PolarRotationEstimator<double>;
143 
144 } /* namespace Alignment */
#define REPORT_ERROR(nerr, ErrormMsg)
Definition: xmipp_error.h:211
template void polarFourierTransform< false >(const MultidimArray< double > &in, Polar< std::complex< double > > &out, bool conjugated, int first_ring, int last_ring, Polar_fftw_plans *&plans, int BsplineOrder)
#define i
double best_rotation(const Polar< std::complex< double > > &I1, const Polar< std::complex< double > > &I2, RotationalCorrelationAux &aux)
Definition: polar.cpp:212
Incorrect argument received.
Definition: xmipp_error.h:113
Definition: polar.h:67
Definition: cpu.h:37
int * n
Some logical error in the pipeline.
Definition: xmipp_error.h:147