Xmipp  v3.23.11-Nereus
resolution_monotomo.h
Go to the documentation of this file.
1 /***************************************************************************
2  *
3  * Authors: Jose Luis Vilas, jlvilas@cnb.csic.es
4  * Oier Lauzirika olauzirika@cnb.csic.es
5  * Carlos Oscar S. Sorzano coss@cnb.csic.es (2016)
6  *
7  * Unidad de Bioinformatica of Centro Nacional de Biotecnologia , CSIC
8  *
9  * This program is free software; you can redistribute it and/or modify
10  * it under the terms of the GNU General Public License as published by
11  * the Free Software Foundation; either version 2 of the License, or
12  * (at your option) any later version.
13  *
14  * This program is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17  * GNU General Public License for more details.
18  *
19  * You should have received a copy of the GNU General Public License
20  * along with this program; if not, write to the Free Software
21  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA
22  * 02111-1307 USA
23  *
24  * All comments concerning this program package may be sent to the
25  * e-mail address 'xmipp@cnb.csic.es'
26  ***************************************************************************/
27 
28 #ifndef _PROG_MONOTOMO
29 #define _PROG_MONOTOMO
30 
31 #include <iostream>
32 #include <core/xmipp_program.h>
33 #include <core/xmipp_image.h>
34 #include <core/xmipp_fft.h>
35 #include <core/xmipp_fftw.h>
36 #include <math.h>
37 #include <limits>
38 #include <complex>
39 #include <data/fourier_filter.h>
40 #include <data/filters.h>
41 #include <string>
42 #include <data/aft.h>
43 #include <memory>
44 
45 
46 class ProgMonoTomo : public XmippProgram
47 {
48 public:
51 
53  double sampling, minRes, maxRes;
54 
58 
61 
64 
65  std::vector<std::complex<float>> fourierSignal, fourierNoise;
66 
67 public:
68 
69  void defineParams();
70  void readParams();
71  void produceSideInfo();
72 
73  /* Mogonogenid amplitud of a volume, given an input volume,
74  * the monogenic amplitud is calculated and low pass filtered at frequency w1*/
75  void amplitudeMonogenicSignal3D(const std::vector<std::complex<float>> &myfftV, float freq, float freqH, float freqL, MultidimArray<float> &amplitude, int count, FileName fnDebug);
76 
77  void firstMonoResEstimation(MultidimArray< std::complex<double> > &myfftV,
78  double freq, double freqH, double freqL, MultidimArray<double> &amplitude,
79  int count, FileName fnDebug, double &mean_Signal,
80  double &mean_noise, double &thresholdFirstEstimation);
81 
83 
84  //Computes the noise distribution inside a box with size boxsize, of a given map, and determines the percentile 95
85  // which is stored in thresholdMatrix.
86  void localNoise(MultidimArray<float> &noiseMap, Matrix2D<double> &noiseMatrix, int boxsize, Matrix2D<double> &thresholdMatrix);
87 
89  std::vector<float> &list);
90 
91  void resolution2eval(int &count_res, double step,
92  double &resolution, double &last_resolution,
93  double &freq, double &freqL,
94  int &last_fourier_idx,
95  bool &continueIter, bool &breakIter);
96 
98 
100  std::vector<float> &list, float &cut_value, float &resolutionThreshold);
101 
102  void getFilteringResolution(size_t idx, float freq, float lastResolution, float freqL, float &resolution);
103 
104  void gaussFilter(const MultidimArray<float> &vol, const float, MultidimArray<float> &VRiesz);
105 
106  void run();
107 
108 public:
110  std::vector<std::complex<float>> fftVRiesz, fftVRiesz_aux;
113  std::unique_ptr<AFT<float>> forward_transformer, backward_transformer;
114 };
116 #endif
void lowestResolutionbyPercentile(MultidimArray< double > &resolutionVol, std::vector< double > &list, double &cut_value, double &resolutionThreshold)
void resolution2eval(int &count_res, double step, double &resolution, double &last_resolution, double &freq, double &freqL, int &last_fourier_idx, bool &continueIter, bool &breakIter)
std::unique_ptr< AFT< float > > forward_transformer
void localNoise(MultidimArray< double > &noiseMap, Matrix2D< double > &noiseMatrix, int boxsize, Matrix2D< double > &thresholdMatrix)
MultidimArray< double > VRiesz
void getFilteringResolution(size_t idx, float freq, float lastResolution, float freqL, float &resolution)
Matrix2D< double > resolutionMatrix
void smoothBorders(MultidimArray< float > &vol, MultidimArray< int > &pMask)
MultidimArray< float > VRiesz
Image< int > mask
std::vector< std::complex< float > > fourierSignal
std::unique_ptr< AFT< float > > backward_transformer
void amplitudeMonogenicSignal3D(MultidimArray< std::complex< double > > &myfftV, double freq, double freqH, double freqL, MultidimArray< double > &amplitude, int count, FileName fnDebug)
std::vector< std::complex< float > > fftVRiesz_aux
MultidimArray< std::complex< double > > fftVRiesz
std::vector< std::complex< float > > fourierNoise
Matrix2D< double > maskMatrix
void gaussFilter(const MultidimArray< float > &vol, const float, MultidimArray< float > &VRiesz)
void firstMonoResEstimation(MultidimArray< std::complex< double > > &myfftV, double freq, double freqH, double freqL, MultidimArray< double > &amplitude, int count, FileName fnDebug, double &mean_Signal, double &mean_noise, double &thresholdFirstEstimation)
void median3x3x3(MultidimArray< double > vol, MultidimArray< double > &filtered)
void postProcessingLocalResolutions(MultidimArray< double > &resolutionVol, std::vector< double > &list)