Xmipp  v3.23.11-Nereus
Public Member Functions | Public Attributes | List of all members
ProgArtZernike3D Class Reference

#include <art_zernike3d.h>

Inheritance diagram for ProgArtZernike3D:
Inheritance graph
[legend]
Collaboration diagram for ProgArtZernike3D:
Collaboration graph
[legend]

Public Member Functions

 ProgArtZernike3D ()
 Empty constructor. More...
 
void readParams ()
 Read argument from command line. More...
 
void show ()
 Show. More...
 
void defineParams ()
 Define parameters. More...
 
void preProcess ()
 
void processImage (const FileName &fnImg, const FileName &fnImgOut, const MDRow &rowIn, MDRow &rowOut)
 
void numCoefficients (int l1, int l2)
 Length of coefficients vector. More...
 
void fillVectorTerms (int l1, int l2)
 Zernike and SPH coefficients allocation. More...
 
void deformVol (MultidimArray< float > &mP, MultidimArray< float > &mW, const MultidimArray< float > &mV, float rot, float tilt, float psi)
 Deform a volumen using Zernike-Spherical harmonic basis. More...
 
- Public Member Functions inherited from XmippMetadataProgram
MetaDatagetInputMd ()
 
MetaDataVecgetOutputMd ()
 
 XmippMetadataProgram ()
 Empty constructor. More...
 
virtual int tryRead (int argc, const char **argv, bool reportErrors=true)
 
virtual void init ()
 
virtual void setup (MetaData *md, const FileName &o="", const FileName &oroot="", bool applyGeo=false, MDLabel label=MDL_IMAGE)
 
virtual ~XmippMetadataProgram ()
 
void setMode (WriteModeMetaData _mode)
 
void setupRowOut (const FileName &fnImgIn, const MDRow &rowIn, const FileName &fnImgOut, MDRow &rowOut) const
 Prepare rowout. More...
 
virtual void wait ()
 Wait for the distributor to finish. More...
 
virtual void checkPoint ()
 For very long programs, it may be needed to write checkpoints. More...
 
- Public Member Functions inherited from XmippProgram
const char * getParam (const char *param, int arg=0)
 
const char * getParam (const char *param, const char *subparam, int arg=0)
 
int getIntParam (const char *param, int arg=0)
 
int getIntParam (const char *param, const char *subparam, int arg=0)
 
double getDoubleParam (const char *param, int arg=0)
 
double getDoubleParam (const char *param, const char *subparam, int arg=0)
 
float getFloatParam (const char *param, int arg=0)
 
float getFloatParam (const char *param, const char *subparam, int arg=0)
 
void getListParam (const char *param, StringVector &list)
 
int getCountParam (const char *param)
 
bool checkParam (const char *param)
 
bool existsParam (const char *param)
 
void addParamsLine (const String &line)
 
void addParamsLine (const char *line)
 
ParamDefgetParamDef (const char *param) const
 
virtual void quit (int exit_code=0) const
 
virtual int tryRun ()
 
void initProgress (size_t total, size_t stepBin=60)
 
void setProgress (size_t value=0)
 
void endProgress ()
 
void processDefaultComment (const char *param, const char *left)
 
void setDefaultComment (const char *param, const char *comment)
 
void setProgramName (const char *name)
 
void addUsageLine (const char *line, bool verbatim=false)
 
void clearUsage ()
 
void addExampleLine (const char *example, bool verbatim=true)
 
void addSeeAlsoLine (const char *seeAlso)
 
void addKeywords (const char *keywords)
 
const char * name () const
 
virtual void usage (int verb=0) const
 
virtual void usage (const String &param, int verb=2)
 
int version () const
 
virtual void read (int argc, const char **argv, bool reportErrors=true)
 
virtual void read (int argc, char **argv, bool reportErrors=true)
 
void read (const String &argumentsLine)
 
 XmippProgram ()
 
 XmippProgram (int argc, const char **argv)
 
virtual ~XmippProgram ()
 

Public Attributes

FileName fnVolR
 
FileName fnVolO
 
FileName fnOutDir
 Output directory. More...
 
int L1
 
int L2
 
Matrix1D< int > vL1
 
Matrix1D< int > vN
 
Matrix1D< int > vL2
 
Matrix1D< int > vM
 
double Ts
 
int RmaxDef
 
bool phaseFlipped
 
bool ignoreCTF
 
float lambda
 
int save_iter
 
bool useCTF
 
bool useZernike
 
int flagEnabled
 
bool resume
 
int niter
 
int sort_last_N
 
MultidimArray< int > mask2D
 
int Xdim
 
Image< float > V
 
Image< float > Vrefined
 
Image< float > Ifilteredp
 
Image< double > I
 
MultidimArray< int > Vmask
 
Image< float > P
 
Image< float > W
 
Image< float > Idiff
 
Matrix2D< float > A
 
float rot
 
float tilt
 
float psi
 
float shiftX
 
float shiftY
 
bool flip
 
bool hasCTF
 
float defocusU
 
float defocusV
 
float defocusAngle
 
CTFDescription ctf
 
FourierFilter FilterCTF
 
int vecSize
 
Matrix1D< float > clnm
 
bool showOptimization
 
MultidimArray< size_t > ordered_list
 
int current_save_iter
 
int num_images
 
int current_iter
 
int initX
 
int endX
 
int initY
 
int endY
 
int initZ
 
int endZ
 
- Public Attributes inherited from XmippMetadataProgram
FileName fn_in
 Filenames of input and output Metadata. More...
 
FileName fn_out
 
FileName baseName
 
FileName pathBaseName
 
FileName oextBaseName
 
bool apply_geo
 Apply geo. More...
 
size_t ndimOut
 Output dimensions. More...
 
size_t zdimOut
 
size_t ydimOut
 
size_t xdimOut
 
DataType datatypeOut
 
size_t mdInSize
 Number of input elements. More...
 
- Public Attributes inherited from XmippProgram
bool doRun
 
bool runWithoutArgs
 
int verbose
 Verbosity level. More...
 
int debug
 

Additional Inherited Members

- Protected Member Functions inherited from XmippMetadataProgram
virtual void initComments ()
 
virtual void postProcess ()
 
virtual bool getImageToProcess (size_t &objId, size_t &objIndex)
 
void show () const override
 
virtual void startProcessing ()
 
virtual void writeOutput ()
 
virtual void showProgress ()
 
virtual void defineLabelParam ()
 
- Protected Member Functions inherited from XmippProgram
void defineCommons ()
 
- Protected Attributes inherited from XmippMetadataProgram
WriteModeMetaData mode
 Metadata writing mode: OVERWRITE, APPEND. More...
 
FileName oext
 Output extension and root. More...
 
FileName oroot
 
MDLabel image_label
 MDLabel to be used to read/write images, usually will be MDL_IMAGE. More...
 
bool produces_an_output
 Indicate that a unique final output is produced. More...
 
bool produces_a_metadata
 Indicate that the unique final output file is a Metadata. More...
 
bool each_image_produces_an_output
 Indicate that an output is produced for each image in the input. More...
 
bool allow_apply_geo
 
bool decompose_stacks
 Input Metadata will treat a stack file as a set of images instead of a unique file. More...
 
bool delete_output_stack
 Delete previous output stack file prior to process images. More...
 
bool get_image_info
 Get the input image file dimensions to further operations. More...
 
bool save_metadata_stack
 Save the associated output metadata when output file is a stack. More...
 
bool track_origin
 Include the original input image filename in the output stack. More...
 
bool keep_input_columns
 Keep input metadata columns. More...
 
bool remove_disabled
 Remove disabled images from the input selfile. More...
 
bool allow_time_bar
 Show process time bar. More...
 
bool input_is_metadata
 Input is a metadata. More...
 
bool single_image
 Input is a single image. More...
 
bool input_is_stack
 Input is a stack. More...
 
bool output_is_stack
 Output is a stack. More...
 
bool create_empty_stackfile
 
bool delete_mdIn
 
size_t time_bar_step
 Some time bar related counters. More...
 
size_t time_bar_size
 
size_t time_bar_done
 
- Protected Attributes inherited from XmippProgram
int errorCode
 
ProgramDefprogDef
 Program definition and arguments parser. More...
 
std::map< String, CommentListdefaultComments
 
int argc
 Original command line arguments. More...
 
const char ** argv
 

Detailed Description

Predict Continuous Parameters.

Definition at line 39 of file art_zernike3d.h.

Constructor & Destructor Documentation

◆ ProgArtZernike3D()

ProgArtZernike3D::ProgArtZernike3D ( )

Empty constructor.

Definition at line 35 of file art_zernike3d.cpp.

36 {
37  resume = false;
38  produces_a_metadata = true;
40  showOptimization = false;
41 }
bool each_image_produces_an_output
Indicate that an output is produced for each image in the input.
bool produces_a_metadata
Indicate that the unique final output file is a Metadata.

Member Function Documentation

◆ defineParams()

void ProgArtZernike3D::defineParams ( )
virtual

Define parameters.

Reimplemented from XmippMetadataProgram.

Definition at line 88 of file art_zernike3d.cpp.

89 {
90  addUsageLine("Template-based canonical volume reconstruction through Zernike3D coefficients");
91  defaultComments["-i"].clear();
92  defaultComments["-i"].addComment("Metadata with initial alignment");
93  defaultComments["-o"].clear();
94  defaultComments["-o"].addComment("Refined volume");
96  addParamsLine(" [--ref <volume=\"\">] : Reference volume");
97  addParamsLine(" [--odir <outputDir=\".\">] : Output directory");
98  addParamsLine(" [--sampling <Ts=1>] : Sampling rate (A/pixel)");
99  addParamsLine(" [--RDef <r=-1>] : Maximum radius of the deformation (px). -1=Half of volume size");
100  addParamsLine(" [--l1 <l1=3>] : Degree Zernike Polynomials=1,2,3,...");
101  addParamsLine(" [--l2 <l2=2>] : Harmonical depth of the deformation=1,2,3,...");
102  addParamsLine(" [--useZernike] : Correct heterogeneity with Zernike3D coefficients");
103  addParamsLine(" [--useCTF] : Correct CTF during ART reconstruction");
104  addParamsLine(" [--phaseFlipped] : Input images have been phase flipped");
105  addParamsLine(" [--regularization <l=0.01>] : ART regularization weight");
106  addParamsLine(" [--niter <n=1>] : Number of ART iterations");
107  addParamsLine(" [--save_iter <s=0>] : Save intermidiate volume after #save_iter iterations");
108  addParamsLine(" [--sort_last <N=2>] : The algorithm sorts projections in the most orthogonally possible way. ");
109  addParamsLine(" : The most orthogonal way is defined as choosing the projection which maximizes the ");
110  addParamsLine(" : dot product with the N previous inserted projections. Use -1 to sort with all ");
111  addParamsLine(" : previous projections");
112  addParamsLine(" [--resume] : Resume processing");
113  addExampleLine("A typical use is:",false);
114  addExampleLine("xmipp_art_zernike3d -i anglesFromContinuousAssignment.xmd --ref reference.vol -o assigned_anglesAndDeformations.xmd --l1 3 --l2 2");
115 }
void addExampleLine(const char *example, bool verbatim=true)
void addUsageLine(const char *line, bool verbatim=false)
void addParamsLine(const String &line)
std::map< String, CommentList > defaultComments
Definition: xmipp_program.h:83

◆ deformVol()

void ProgArtZernike3D::deformVol ( MultidimArray< float > &  mP,
MultidimArray< float > &  mW,
const MultidimArray< float > &  mV,
float  rot,
float  tilt,
float  psi 
)

Deform a volumen using Zernike-Spherical harmonic basis.

◆ fillVectorTerms()

void ProgArtZernike3D::fillVectorTerms ( int  l1,
int  l2 
)

Zernike and SPH coefficients allocation.

Definition at line 263 of file art_zernike3d.cpp.

264 {
265  int idx = 0;
270  for (int h=0; h<=l2; h++) {
271  int totalSPH = 2*h+1;
272  int aux = totalSPH/2;
273  for (int l=h; l<=l1; l+=2) {
274  for (int m=0; m<totalSPH; m++) {
275  VEC_ELEM(vL1,idx) = l;
276  VEC_ELEM(vN,idx) = h;
277  VEC_ELEM(vL2,idx) = h;
278  VEC_ELEM(vM,idx) = m-aux;
279  idx++;
280  }
281  }
282  }
283 }
#define VEC_ELEM(v, i)
Definition: matrix1d.h:245
Matrix1D< int > vN
Definition: art_zernike3d.h:51
Matrix1D< int > vL2
Definition: art_zernike3d.h:51
Matrix1D< int > vM
Definition: art_zernike3d.h:51
void initZeros()
Definition: matrix1d.h:592
int m
Matrix1D< int > vL1
Definition: art_zernike3d.h:51

◆ numCoefficients()

void ProgArtZernike3D::numCoefficients ( int  l1,
int  l2 
)

Length of coefficients vector.

Definition at line 247 of file art_zernike3d.cpp.

248 {
249  for (int h=0; h<=l2; h++)
250  {
251  int numSPH = 2*h+1;
252  int count=l1-h+1;
253  int numEven=(count>>1)+(count&1 && !(h&1));
254  if (h%2 == 0) {
255  vecSize += numSPH*numEven;
256  }
257  else {
258  vecSize += numSPH*(l1-h+1-numEven);
259  }
260  }
261 }

◆ preProcess()

void ProgArtZernike3D::preProcess ( )
virtual

Produce side info. An exception is thrown if any of the files is not found

Reimplemented from XmippMetadataProgram.

Definition at line 117 of file art_zernike3d.cpp.

118 {
119 
120  // Check that metadata has all information neede
121  if (!getInputMd()->containsLabel(MDL_ANGLE_ROT) ||
122  !getInputMd()->containsLabel(MDL_ANGLE_TILT) ||
123  !getInputMd()->containsLabel(MDL_ANGLE_PSI))
124  {
125  REPORT_ERROR(ERR_MD_MISSINGLABEL,"Input metadata projection angles are missing. Exiting...");
126  }
127 
128  if (fnVolR != "")
129  {
130  V.read(fnVolR);
131  }
132  else
133  {
134  FileName fn_first_image;
135  Image<float> first_image;
136  getInputMd()->getRow(1)->getValue(MDL_IMAGE,fn_first_image);
137  first_image.read(fn_first_image);
138  int Xdim_first = static_cast<int>(XSIZE(first_image()));
139  V().initZeros(Xdim_first, Xdim_first, Xdim_first);
140 
141  }
142  V().setXmippOrigin();
143 
144  Xdim = static_cast<int>(XSIZE(V()));
145 
146  if (resume && fnVolO.exists()) {
148  } else {
149  Vrefined() = V();
150  }
151  Vrefined().setXmippOrigin();
152 
153  if (RmaxDef<0)
154  RmaxDef = Xdim/2;
155 
156  // Transformation matrix
157  A.initIdentity(3);
158 
159  // CTF Filter
162  FilterCTF.ctf.enable_CTFnoise = false;
163  FilterCTF.ctf.enable_CTF = true;
164 
165  // Area where Zernike3D basis is computed (and volume is updated)
166  Mask mask;
167  mask.type = BINARY_CIRCULAR_MASK;
168  mask.mode = INNER_MASK;
169  mask.R1 = RmaxDef;
170  mask.generate_mask(V());
171  Vmask = mask.get_binary_mask();
173 
174  // Area Zernike3D in 2D
175  mask.generate_mask(Xdim, Xdim);
176  mask2D = mask.get_binary_mask();
178 
179  vecSize = 0;
182 
183  initX = STARTINGX(Vrefined());
184  endX = FINISHINGX(Vrefined());
185  initY = STARTINGY(Vrefined());
186  endY = FINISHINGY(Vrefined());
187  initZ = STARTINGZ(Vrefined());
188  endZ = FINISHINGZ(Vrefined());
189 }
Rotation angle of an image (double,degrees)
CTFDescription ctf
#define FINISHINGX(v)
#define CTFINV
#define REPORT_ERROR(nerr, ErrormMsg)
Definition: xmipp_error.h:211
MultidimArray< int > mask2D
Definition: art_zernike3d.h:79
Definition: mask.h:360
Tilting angle of an image (double,degrees)
Special label to be used when gathering MDs in MpiMetadataPrograms.
#define FINISHINGZ(v)
#define STARTINGX(v)
bool enable_CTF
Enable CTF part.
Definition: ctf.h:275
#define STARTINGY(v)
FourierFilter FilterCTF
MultidimArray< int > Vmask
Definition: art_zernike3d.h:87
double R1
Definition: mask.h:413
#define XSIZE(v)
Missing expected label.
Definition: xmipp_error.h:158
int type
Definition: mask.h:402
virtual std::unique_ptr< MDRow > getRow(size_t id)=0
Image< float > V
Definition: art_zernike3d.h:83
Image< float > Vrefined
Definition: art_zernike3d.h:83
Matrix2D< float > A
Definition: art_zernike3d.h:95
bool exists() const
#define FINISHINGY(v)
void generate_mask(bool apply_geo=false)
Definition: mask.cpp:1577
#define BINARY_CIRCULAR_MASK
Definition: mask.h:365
int read(const FileName &name, DataMode datamode=DATA, size_t select_img=ALL_IMAGES, bool mapData=false, int mode=WRITE_READONLY)
void numCoefficients(int l1, int l2)
Length of coefficients vector.
bool enable_CTFnoise
Enable CTFnoise part.
Definition: ctf.h:273
const MultidimArray< int > & get_binary_mask() const
Definition: mask.h:707
#define STARTINGZ(v)
Name of an image (std::string)
void fillVectorTerms(int l1, int l2)
Zernike and SPH coefficients allocation.
void initIdentity()
Definition: matrix2d.h:673
int mode
Definition: mask.h:407
constexpr int INNER_MASK
Definition: mask.h:47

◆ processImage()

void ProgArtZernike3D::processImage ( const FileName fnImg,
const FileName fnImgOut,
const MDRow rowIn,
MDRow rowOut 
)
virtual

Predict angles and shift. At the input the pose parameters must have an initial guess of the parameters. At the output they have the estimated pose.

Implements XmippMetadataProgram.

Definition at line 196 of file art_zernike3d.cpp.

197 {
198  rowOut=rowIn; // FIXME have a look if this is needed. I don't think so, or see how to do this automatically in xmipp_metadata_program.cpp
199  flagEnabled=1;
200 
201  rowIn.getValue(MDL_ANGLE_ROT,rot);
203  rowIn.getValue(MDL_ANGLE_PSI,psi);
206  std::vector<double> vectortemp;
207  if (useZernike) {
208  rowIn.getValue(MDL_SPH_COEFFICIENTS,vectortemp);
209  clnm.initZeros(vectortemp.size()-8);
210  for(int i=0; i < vectortemp.size()-8; i++){
211  VEC_ELEM(clnm,i) = static_cast<float>(vectortemp[i]);
212  }
213  removeOverdeformation();
214  }
215  rowIn.getValueOrDefault(MDL_FLIP,flip, false);
216 
218  {
219  hasCTF=true;
220  FilterCTF.ctf.readFromMdRow(rowIn, false);
221  FilterCTF.ctf.Tm = Ts;
223  }
224  else
225  hasCTF=false;
226  MAT_ELEM(A,0,2)=shiftX;
227  MAT_ELEM(A,1,2)=shiftY;
228  MAT_ELEM(A,0,0)=1;
229  MAT_ELEM(A,0,1)=0;
230  MAT_ELEM(A,1,0)=0;
231  MAT_ELEM(A,1,1)=1;
232 
233  if (verbose>=2)
234  std::cout << "Processing " << fnImg << std::endl;
235 
236  I.read(fnImg);
237  I().setXmippOrigin();
238 
239  // Forward Model
240  artModel<Direction::Forward>();
241 
242  // ART update
243  artModel<Direction::Backward>();
244 
245 }
Rotation angle of an image (double,degrees)
#define VEC_ELEM(v, i)
Definition: matrix1d.h:245
Defocus U (Angstroms)
CTFDescription ctf
Tilting angle of an image (double,degrees)
Shift for the image in the X axis (double)
Special label to be used when gathering MDs in MpiMetadataPrograms.
Name for the CTF Model (std::string)
Image< double > I
Definition: art_zernike3d.h:85
#define i
#define MAT_ELEM(m, i, j)
Definition: matrix2d.h:116
FourierFilter FilterCTF
T & getValue(MDLabel label)
Flip the image? (bool)
Matrix1D< float > clnm
double Tm
Sampling rate (A/pixel)
Definition: ctf.h:240
void readFromMdRow(const MDRow &row, bool disable_if_not_K=true)
Definition: ctf.cpp:1172
int verbose
Verbosity level.
Matrix2D< float > A
Definition: art_zernike3d.h:95
void initZeros()
Definition: matrix1d.h:592
const T & getValueOrDefault(MDLabel label, const T &def) const
virtual bool containsLabel(MDLabel label) const =0
void produceSideInfo()
Produce Side information.
Definition: ctf.cpp:1392
Deformation coefficients.
int read(const FileName &name, DataMode datamode=DATA, size_t select_img=ALL_IMAGES, bool mapData=false, int mode=WRITE_READONLY)
Shift for the image in the Y axis (double)

◆ readParams()

void ProgArtZernike3D::readParams ( )
virtual

Read argument from command line.

Reimplemented from XmippMetadataProgram.

Definition at line 44 of file art_zernike3d.cpp.

45 {
47  fnVolR = getParam("--ref");
48  fnOutDir = getParam("--odir");
49  RmaxDef = getIntParam("--RDef");
50  phaseFlipped = checkParam("--phaseFlipped");
51  useCTF = checkParam("--useCTF");
52  Ts = getDoubleParam("--sampling");
53  L1 = getIntParam("--l1");
54  L2 = getIntParam("--l2");
55  useZernike = checkParam("--useZernike");
56  lambda = static_cast<float>(getDoubleParam("--regularization"));
57  resume = checkParam("--resume");
58  niter = getIntParam("--niter");
59  save_iter = getIntParam("--save_iter");
60  sort_last_N = getIntParam("--sort_last");
61  fnVolO = fnOutDir + "/Refined.vol";
62  keep_input_columns = true;
63 }
double getDoubleParam(const char *param, int arg=0)
FileName fnOutDir
Output directory.
Definition: art_zernike3d.h:47
const char * getParam(const char *param, int arg=0)
bool keep_input_columns
Keep input metadata columns.
bool checkParam(const char *param)
int getIntParam(const char *param, int arg=0)

◆ show()

void ProgArtZernike3D::show ( )

Show.

Definition at line 66 of file art_zernike3d.cpp.

67 {
68  if (!verbose)
69  return;
71  std::cout
72  << "Output directory: " << fnOutDir << std::endl
73  << "Reference volume: " << fnVolR << std::endl
74  << "Sampling: " << Ts << std::endl
75  << "Max. Radius Deform. " << RmaxDef << std::endl
76  << "Zernike Degree: " << L1 << std::endl
77  << "SH Degree: " << L2 << std::endl
78  << "Correct CTF: " << useCTF << std::endl
79  << "Correct heretogeneity: " << useZernike << std::endl
80  << "Phase flipped: " << phaseFlipped << std::endl
81  << "Regularization: " << lambda << std::endl
82  << "Number of iterations: " << niter << std::endl
83  << "Save every # iterations: " << save_iter << std::endl
84  ;
85 }
FileName fnOutDir
Output directory.
Definition: art_zernike3d.h:47
int verbose
Verbosity level.
void show() const override

Member Data Documentation

◆ A

Matrix2D<float> ProgArtZernike3D::A

Definition at line 95 of file art_zernike3d.h.

◆ clnm

Matrix1D<float> ProgArtZernike3D::clnm

Definition at line 113 of file art_zernike3d.h.

◆ ctf

CTFDescription ProgArtZernike3D::ctf

Definition at line 107 of file art_zernike3d.h.

◆ current_iter

int ProgArtZernike3D::current_iter

Definition at line 123 of file art_zernike3d.h.

◆ current_save_iter

int ProgArtZernike3D::current_save_iter

Definition at line 119 of file art_zernike3d.h.

◆ defocusAngle

float ProgArtZernike3D::defocusAngle

Definition at line 105 of file art_zernike3d.h.

◆ defocusU

float ProgArtZernike3D::defocusU

Definition at line 105 of file art_zernike3d.h.

◆ defocusV

float ProgArtZernike3D::defocusV

Definition at line 105 of file art_zernike3d.h.

◆ endX

int ProgArtZernike3D::endX

Definition at line 125 of file art_zernike3d.h.

◆ endY

int ProgArtZernike3D::endY

Definition at line 125 of file art_zernike3d.h.

◆ endZ

int ProgArtZernike3D::endZ

Definition at line 125 of file art_zernike3d.h.

◆ FilterCTF

FourierFilter ProgArtZernike3D::FilterCTF

Definition at line 109 of file art_zernike3d.h.

◆ flagEnabled

int ProgArtZernike3D::flagEnabled

Definition at line 69 of file art_zernike3d.h.

◆ flip

bool ProgArtZernike3D::flip

Definition at line 101 of file art_zernike3d.h.

◆ fnOutDir

FileName ProgArtZernike3D::fnOutDir

Output directory.

Definition at line 47 of file art_zernike3d.h.

◆ fnVolO

FileName ProgArtZernike3D::fnVolO

Filename of the refined volume

Definition at line 45 of file art_zernike3d.h.

◆ fnVolR

FileName ProgArtZernike3D::fnVolR

Filename of the reference volume

Definition at line 43 of file art_zernike3d.h.

◆ hasCTF

bool ProgArtZernike3D::hasCTF

Definition at line 103 of file art_zernike3d.h.

◆ I

Image<double> ProgArtZernike3D::I

Definition at line 85 of file art_zernike3d.h.

◆ Idiff

Image<float> ProgArtZernike3D::Idiff

Definition at line 93 of file art_zernike3d.h.

◆ Ifilteredp

Image<float> ProgArtZernike3D::Ifilteredp

Definition at line 83 of file art_zernike3d.h.

◆ ignoreCTF

bool ProgArtZernike3D::ignoreCTF

Definition at line 59 of file art_zernike3d.h.

◆ initX

int ProgArtZernike3D::initX

Definition at line 125 of file art_zernike3d.h.

◆ initY

int ProgArtZernike3D::initY

Definition at line 125 of file art_zernike3d.h.

◆ initZ

int ProgArtZernike3D::initZ

Definition at line 125 of file art_zernike3d.h.

◆ L1

int ProgArtZernike3D::L1

Degrees of Zernike polynomials and spherical harmonics

Definition at line 49 of file art_zernike3d.h.

◆ L2

int ProgArtZernike3D::L2

Definition at line 49 of file art_zernike3d.h.

◆ lambda

float ProgArtZernike3D::lambda

Definition at line 61 of file art_zernike3d.h.

◆ mask2D

MultidimArray<int> ProgArtZernike3D::mask2D

Definition at line 79 of file art_zernike3d.h.

◆ niter

int ProgArtZernike3D::niter

Definition at line 75 of file art_zernike3d.h.

◆ num_images

int ProgArtZernike3D::num_images

Definition at line 121 of file art_zernike3d.h.

◆ ordered_list

MultidimArray<size_t> ProgArtZernike3D::ordered_list

Definition at line 117 of file art_zernike3d.h.

◆ P

Image<float> ProgArtZernike3D::P

Definition at line 89 of file art_zernike3d.h.

◆ phaseFlipped

bool ProgArtZernike3D::phaseFlipped

Definition at line 57 of file art_zernike3d.h.

◆ psi

float ProgArtZernike3D::psi

Definition at line 97 of file art_zernike3d.h.

◆ resume

bool ProgArtZernike3D::resume

Resume computations

Definition at line 73 of file art_zernike3d.h.

◆ RmaxDef

int ProgArtZernike3D::RmaxDef

Maximum radius

Definition at line 55 of file art_zernike3d.h.

◆ rot

float ProgArtZernike3D::rot

Definition at line 97 of file art_zernike3d.h.

◆ save_iter

int ProgArtZernike3D::save_iter

Definition at line 63 of file art_zernike3d.h.

◆ shiftX

float ProgArtZernike3D::shiftX

Definition at line 99 of file art_zernike3d.h.

◆ shiftY

float ProgArtZernike3D::shiftY

Definition at line 99 of file art_zernike3d.h.

◆ showOptimization

bool ProgArtZernike3D::showOptimization

Definition at line 115 of file art_zernike3d.h.

◆ sort_last_N

int ProgArtZernike3D::sort_last_N

Definition at line 77 of file art_zernike3d.h.

◆ tilt

float ProgArtZernike3D::tilt

Definition at line 97 of file art_zernike3d.h.

◆ Ts

double ProgArtZernike3D::Ts

Sampling rate

Definition at line 53 of file art_zernike3d.h.

◆ useCTF

bool ProgArtZernike3D::useCTF

Definition at line 65 of file art_zernike3d.h.

◆ useZernike

bool ProgArtZernike3D::useZernike

Definition at line 67 of file art_zernike3d.h.

◆ V

Image<float> ProgArtZernike3D::V

Definition at line 83 of file art_zernike3d.h.

◆ vecSize

int ProgArtZernike3D::vecSize

Definition at line 111 of file art_zernike3d.h.

◆ vL1

Matrix1D<int> ProgArtZernike3D::vL1

Zernike and SPH coefficients vectors

Definition at line 51 of file art_zernike3d.h.

◆ vL2

Matrix1D<int> ProgArtZernike3D::vL2

Definition at line 51 of file art_zernike3d.h.

◆ vM

Matrix1D<int> ProgArtZernike3D::vM

Definition at line 51 of file art_zernike3d.h.

◆ Vmask

MultidimArray<int> ProgArtZernike3D::Vmask

Definition at line 87 of file art_zernike3d.h.

◆ vN

Matrix1D<int> ProgArtZernike3D::vN

Definition at line 51 of file art_zernike3d.h.

◆ Vrefined

Image<float> ProgArtZernike3D::Vrefined

Definition at line 83 of file art_zernike3d.h.

◆ W

Image<float> ProgArtZernike3D::W

Definition at line 91 of file art_zernike3d.h.

◆ Xdim

int ProgArtZernike3D::Xdim

Definition at line 81 of file art_zernike3d.h.


The documentation for this class was generated from the following files: