Xmipp  v3.23.11-Nereus
forward_zernike_volume.h
Go to the documentation of this file.
1 /***************************************************************************
2  *
3  * Authors: Carlos Oscar coss@cnb.csic.es
4  * David Herreros Calero dherreros@cnb.csic.es
5  *
6  * This program is free software; you can redistribute it and/or modify
7  * it under the terms of the GNU General Public License as published by
8  * the Free Software Foundation; either version 2 of the License, or
9  * (at your option) any later version.
10  *
11  * This program is distributed in the hope that it will be useful,
12  * but WITHOUT ANY WARRANTY; without even the implied warranty of
13  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14  * GNU General Public License for more details.
15  *
16  * You should have received a copy of the GNU General Public License
17  * along with this program; if not, write to the Free Software
18  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA
19  * 02111-1307 USA
20  *
21  * All comments concerning this program package may be sent to the
22  * e-mail address 'xmipp@cnb.uam.es'
23  ***************************************************************************/
24 
25 #ifndef _PROG_VOL_DEFORM_SPH
26 #define _PROG_VOL_DEFORM_SPH
27 
28 #include <vector>
29 #include "core/xmipp_program.h"
30 #include "core/xmipp_image.h"
31 #include <data/blobs.h>
32 
38 {
39 public:
42 
45 
48 
51 
54 
57 
60 
63 
65  int L1, L2;
66 
68  double Rmax;
69 
70 public:
72  int vecSize;
73 
76 
78  std::vector<double> absMaxR_vec;
79 
80  //Deformation in pixels, sumVI, sumVD
82 
83  // Regularization
84  double lambda;
85 
86  // Save output volume
88 
89  // Save the values of gx, gy and gz for local strain and rotation analysis
91 
92  // Blob
93  struct blobtype blob;
94  double blob_r;
95 
96  double sigma4;
97  // Gaussian projection table
99 
100  // Gaussian projection2 table
102 
103  // Loop step
105 
106  // Vector containing the degree of the spherical harmonics
108 
111 
112  //Copy of Optimizer steps
114 
115  std::vector<double>vec;
117 
118 public:
120  void defineParams();
121 
123  void readParams();
124 
126  void show();
127 
129  double distance(double *pclnm);
130 
132  void run();
133 
134  // /// Copy the coefficients from harmonical depth n-1 vector to harmonical depth n vector
135  // void copyvectors(Matrix1D<double> &oldvect,Matrix1D<double> &newvect);
136 
137  // /// Determine the positions to be minimize of a vector containing spherical harmonic coefficients
138  // void minimizepos(Matrix1D<double> &vectpos, Matrix1D<double> &prevpos);
139 
141  void minimizepos(int L1, int l2, Matrix1D<double> &steps);
142 
143  // ///Compute the number of spherical harmonics in l=0,1,...,depth
144  // void Numsph(Matrix1D<int> &sphD);
145 
147  void numCoefficients(int l1, int l2, int &vecSize);
148 
150  void fillVectorTerms(int l1, int l2);
151 
153  void computeStrain();
154 
156  void writeVector(std::string outPath, Matrix1D<double> v, bool append);
157 
158  void splattingAtPos(std::array<double, 3> r, double weight, MultidimArray<double> &mVO1, MultidimArray<double> &mVO2);
159 
160  template<bool SAVE_DEFORMATION>
161  void deformVolume();
162 
163  double splatVal(std::array<double, 3> r, double weight, const MultidimArray<double> &mV);
164 
165  std::string readNthLine(int N) const;
166 
167  std::vector<double> string2vector(std::string const &s) const;
168 
170 
171  void volume2Mask(MultidimArray<double> &vol, double thr);
172 
173  void rmsd(MultidimArray<double> vol1, MultidimArray<double> vol2, double &val);
174 };
175 
177 #endif
bool analyzeStrain
Save the deformation of each voxel for local strain and rotation analysis.
std::vector< double > absMaxR_vec
Maxima of reference volumes (in absolute value)
std::string readNthLine(int N) const
void volume2Blobs(MultidimArray< double > &vol, MultidimArray< double > &vol2, const MultidimArray< double > &mV, const MultidimArray< int > &mask)
void computeStrain()
Compute strain.
FileName fnMaskR
Filename of the reference volume mask.
FileName fnVolR
Reference volume.
int vecSize
Coefficient vector size.
Matrix1D< double > gaussianProjectionTable2
void minimizepos(int L1, int l2, Matrix1D< double > &steps)
Determine the positions to be minimize of a vector containing spherical harmonic coefficients.
Image< double > VI
Images.
double splatVal(std::array< double, 3 > r, double weight, const MultidimArray< double > &mV)
double Rmax
Maximum radius for the transformation.
MultidimArray< int > V_maskr
FileName fnVolI
Volume to deform.
void fillVectorTerms(int l1, int l2)
Zernike and SPH coefficients allocation.
void writeVector(std::string outPath, Matrix1D< double > v, bool append)
Save vector to file.
void defineParams()
Define params.
MultidimArray< int > V_maski
3D mask for reference volume
int L1
Degree of Zernike polynomials and spherical harmonics.
double distance(double *pclnm)
Distance.
std::vector< double > vec
double steps
std::vector< double > string2vector(std::string const &s) const
bool optimizeRadius
Radius optimization.
Matrix1D< double > gaussianProjectionTable
void rmsd(MultidimArray< double > vol1, MultidimArray< double > vol2, double &val)
void splattingAtPos(std::array< double, 3 > r, double weight, MultidimArray< double > &mVO1, MultidimArray< double > &mVO2)
FileName fnRoot
Root name for several output files.
MultidimArray< int > V_mask2
void readParams()
Read arguments from command line.
FileName fnVolOut
Output Volume (deformed input volume)
void numCoefficients(int l1, int l2, int &vecSize)
Length of coefficients vector.
void volume2Mask(MultidimArray< double > &vol, double thr)