Xmipp  v3.23.11-Nereus
reconstruct_ADMM.h
Go to the documentation of this file.
1 /***************************************************************************
2  *
3  * Authors: Masih Nilchian (masih.nilchian@epfl.ch)
4  * Carlos Oscar S. Sorzano (coss@cnb.csic.es)
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 #ifndef _PROG_RECONSTRUCT_ADMM_HH
27 #define _PROG_RECONSTRUCT_ADMM_HH
28 
29 #include <core/xmipp_program.h>
30 #include <core/matrix1d.h>
31 #include <core/xmipp_fftw.h>
32 #include <data/ctf.h>
33 #include <data/mask.h>
34 #include <core/symmetries.h>
35 
39 
41 {
42 public:
43  double supp;
45  double autocorrStep;
46  double alpha;
51 
52  void initializeKernel(double alpha, double a, double _projectionStep);
53 
54  double projectionValueAt(double u, double v);
55 
56  void convolveKernelWithItself(double _autocorrStep);
57 
58  void applyCTFToKernelAutocorrelation(CTFDescription &ctf, double Ts, MultidimArray<double> &autocorrelationWithCTF);
59 
60  void getKernelAutocorrelation(MultidimArray<double> &autocorrelation);
61 
62  /* Direction=x,y,z, mode=L (Lx), T (L^T x), 2 (L^TLx) */
63  void computeGradient(MultidimArray<double> &gradient, char direction, bool adjoint=false);
64 
67 };
68 
70 {
71 public:
72  FileName fnIn; // Set of input images
73  FileName fnRoot; // Rootname for output files
74  FileName fnFirst; // First reconstruction filename
75  FileName fnHtb; // Htb filename
76  FileName fnHtKH; // HtKH filename
78  double a; // Kaiser-Bessel maximum radius
79  double alpha; // Kaiser-Bessel shape
80  double mu; // Augmented Lagrange
81  double lambda; // Total variation regularization
82  double lambda1; // Tikhonov regularization
83 
84  double Ti; // Downsampling factor for the volume
85  double Tp; // Downsampling factor for the projections
86  bool useWeights; // Use weights if available in the input metadata
87  bool useCTF; // Use CTF if available in the input metadata
88  int Ncgiter; // Number of conjugate gradient iterations
89  int Nadmmiter; // Number of ADMM iterations
90  bool positivity; // Positivity constraint
91  bool applyMask;
92  double Ts; // Sampling rate
93  Mask mask; // Mask
96  size_t Nprocs;
97  size_t rank;
98 public:
100  void defineParams();
101  void readParams();
102  void show();
103  void run();
104 
105  void produceSideInfo();
106 
108  void project(double rot, double tilt, double psi, MultidimArray<double> &P, bool adjoint=false, double weight=1.);
109 
111  void project(const Matrix1D<double> &r1, const Matrix1D<double> &r2, MultidimArray<double> &P, bool adjoint=false, double weight=1.);
112 
115  void constructHtb();
116 
120  void computeHtKH(MultidimArray<double> &kernelV);
121 
123  void addRegularizationTerms();
124 
126  void applyKernel3D(MultidimArray<double> &x, MultidimArray<double> &AtAx);
127 
129  void applyLFilter(MultidimArray< std::complex<double> > &fourierL, bool adjoint=false);
130  void applyLtFilter(MultidimArray< std::complex<double> > &fourierL, MultidimArray<double> &u, MultidimArray<double> &d);
131 
133  void applyConjugateGradient();
134 
136  void doPOCSProjection();
137 
139  void updateUD();
140 
142  void produceVolume();
143 
144  // Share a volume among nodes
145  virtual void shareVolume(MultidimArray<double> &V) {}
146 
147  // Redefine how to synchronize
148  virtual void synchronize() {}
149 
150 public:
152  Image<double> CHtb; // First reconstructed volume
153  Image<double> Ck, Vk; // Reconstructed volume
154  MetaDataVec mdIn; // Set of images and angles
158 
159  MultidimArray<double> ux, uy, uz, dx, dy, dz, ud;
162 };
163 
165 #endif
MultidimArray< std::complex< double > > FourierProjectionAutocorr
MultidimArray< double > uz
MetaDataVec mdIn
Definition: mask.h:360
MultidimArray< double > projectionAutocorrWithCTF
double projectionValueAt(double u, double v)
void convolveKernelWithItself(double _autocorrStep)
MultidimArray< std::complex< double > > fourierLz
doublereal * x
void computeGradient(MultidimArray< double > &gradient, char direction, bool adjoint=false)
doublereal * d
void applyCTFToKernelAutocorrelation(CTFDescription &ctf, double Ts, MultidimArray< double > &autocorrelationWithCTF)
double dx
Image< double > CHtb
void getKernelAutocorrelation(MultidimArray< double > &autocorrelation)
MultidimArray< double > paddedx
void initializeKernel(double alpha, double a, double _projectionStep)
void direction(const MultidimArray< double > &orMap, MultidimArray< double > &qualityMap, double lambda, int size, MultidimArray< double > &dirMap, int x, int y)
MultidimArray< std::complex< double > > fourierKernelV
double projectionStep
std::string String
Definition: xmipp_strings.h:34
double psi(const double x)
doublereal * u
float r2
FourierTransformer transformerPaddedx
virtual void shareVolume(MultidimArray< double > &V)
Image< double > Vk
double autocorrStep
doublereal * a
double Ck
FourierTransformer transformer
Matrix1D< double > projectionProfile
void computeKernel3D(MultidimArray< double > &kernel)
virtual void synchronize()
float r1