Xmipp  v3.23.11-Nereus
movie_alignment_correlation_gpu.h
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 
26 #pragma once
27 
29 #include "data/fft_settings.h"
33 
37 template<typename T>
39 protected:
40  void defineParams();
41  void show();
42  void readParams();
43 
44 private:
45 
51  auto findGoodCorrelationSize(const Dimensions &ref);
52 
57  auto findGoodPatchSize();
58 
64  auto findGoodCropSize(const Dimensions &ref);
65 
71  auto getCorrelationHint(const Dimensions &d);
72 
78  std::vector<FramePatchMeta<T>> getPatchesLocation(const std::pair<T, T> &borders,
79  const Dimensions &patch);
80 
88  void getPatchData(const Rectangle<Point2D<T>> &patch,
89  const AlignmentResult<T> &globAlignment,
90  T *result);
91 
99  std::pair<T,T> getMovieBorders(const AlignmentResult<T> &globAlignment, int verbose);
100 
104  void LAOptimize();
105 
106  AlignmentResult<T> computeGlobalAlignment(const MetaData &movie,
107  const Image<T> &dark,
108  const Image<T> &igain) override;
109 
116  T* setFilter(float scale, const Dimensions &dims);
117 
119  struct AlignmentContext
120  {
121  int verbose;
122  float maxShift;
123  size_t N; // number of frames
124  std::pair<T, T> scale;
125  core::optional<size_t> refFrame;
126  Dimensions out = Dimensions(0);
127  size_t alignmentBytes() const {
128  return (N * (N-1) / 2) * 2 * sizeof(float); // 2D position for each correlation
129  }
130  };
131 
137  LocalAlignmentResult<T> localFromGlobal(const AlignmentResult<T> &globAlignment);
138 
139  void applyShiftsComputeAverage(
140  const MetaData& movieMD, const Image<T>& dark, const Image<T>& igain,
141  Image<T>& initialMic, size_t& Ninitial, Image<T>& averageMicrograph,
142  size_t& N, const AlignmentResult<T> &globAlignment) override;
143 
144 
148  auto getOutputStreamCount();
149 
156  void storeSizes(const Dimensions &orig, const Dimensions &opt, bool applyCrop);
157 
164  std::optional<Dimensions> getStoredSizes(const Dimensions &dim, bool applyCrop);
165 
169  void GAOptimize();
170 
179  void report(float scale,
180  typename CUDAFlexAlignScale<T>::Params &SP,
181  typename CUDAFlexAlignCorrelate<T>::Params &CP,
182  int streams, int threads);
183 
190  void getCroppedFrame(const Dimensions &settings, T *output, T *src);
191 
195  class Movie final
196  {
197  public:
204  auto set(const Dimensions &dim) {
205  mDim = dim;
206  mFrames.reserve(dim.n());
207  for (auto n = 0; n < dim.n(); ++n) {
208  // create multidim arrays, but without data
209  mFrames.emplace_back(1, 1, dim.y(), dim.x(), nullptr);
210  }
211  }
212 
216  auto &getDim() const {
217  return mDim;
218  }
219 
224  auto &getFrame(size_t index) const {
225  return mFrames[index];
226  }
227 
233  auto setFrameData(size_t index, T *ptr) {
234  mFrames[index].data = ptr;
235  }
236 
237  private:
238  std::vector<MultidimArray<T>> mFrames;
239  Dimensions mDim = Dimensions(0);
240  };
241 
250  T* loadFrame(const MetaData& movieMD, const Image<T>& dark, const Image<T>& igain, size_t index);
251 
258  auto computeAlignment(float *positions, const AlignmentContext &context);
259 
260  void releaseAll() override;
261 
262  LocalAlignmentResult<T> computeLocalAlignment(const MetaData &movie,
263  const Image<T> &dark, const Image<T> &igain,
264  const AlignmentResult<T> &globAlignment) override;
265 
273  std::string getKey(const std::string &keyword, const Dimensions &dim, bool crop) const {
274  std::stringstream ss;
275  ss << version << " " << mGpu.value().getUUID() << keyword << dim << " " << crop;
276  return ss.str();
277  }
278 
279  void applyShiftsComputeAverage(
280  const MetaData& movie, const Image<T>& dark, const Image<T>& igain,
281  Image<T>& initialMic, size_t& Ninitial, Image<T>& averageMicrograph,
282  size_t& N, const LocalAlignmentResult<T> &alignment) override;
283 
284  // holding the loaded movie
285  Movie movie;
286  // number of GPU streams for global alignment
287  int GAStreams = 1;
288  // parameters for global alignment scaling
289  typename CUDAFlexAlignScale<T>::Params GASP;
290  // parameters for global alignment correlation
291  typename CUDAFlexAlignCorrelate<T>::Params GACP;
292  // parameters for local alignment scaling
293  typename CUDAFlexAlignScale<T>::Params LASP;
294  // parameters for local alignment correlation
295  typename CUDAFlexAlignCorrelate<T>::Params LACP;
296 
298  int patchesAvg;
299 
301  std::string storage;
302 
303  core::optional<GPU> mGpu;
304 
308  std::string minMemoryStr = std::string("minMem");
309  std::string optSizeXStr = std::string("optSizeX");
310  std::string optSizeYStr = std::string("optSizeY");
311  std::string optBatchSizeStr = std::string("optBatchSize");
312 
313  static constexpr auto version = "2.0";
314 };
315 
316 //@}
int version() const
doublereal * d
CUDA_HD constexpr size_t x() const
Definition: dimensions.h:51
cudaStream_t * streams
void readParams()
Read argument from command line.
viol index
CUDA_HD constexpr size_t y() const
Definition: dimensions.h:60
int verbose
Verbosity level.
CUDA_HD constexpr size_t n() const
Definition: dimensions.h:78
int * n