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

#include <mpi_ml_align2d.h>

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

Public Member Functions

 MpiProgML2D ()
 
 MpiProgML2D (MpiNode *node)
 
void setNumberOfLocalImages ()
 
void produceSideInfo2 ()
 
void writeOutputFiles (const ModelML2D &model, OutputType outputType=OUT_FINAL)
 Write model parameters. More...
 
void expectation ()
 After normal ML2D expectation, data must be collected from nodes. More...
 
void endIteration ()
 Redefine endIteration for some synchronization. More...
 
void printModel (const String &msg, const ModelML2D &model)
 
virtual void usage (int verb=0) const
 
- Public Member Functions inherited from ProgML2D
void readParams ()
 Read arguments from command line. More...
 
void defineParams ()
 Params definition. More...
 
 ProgML2D ()
 
virtual void show ()
 Show info at starting program. More...
 
virtual void produceSideInfo ()
 Try to merge produceSideInfo1 and 2. More...
 
void calculatePdfInplane ()
 Calculate probability density distribution for in-plane transformations. More...
 
void rotateReference ()
 Fill vector of matrices with all rotations of reference. More...
 
void reverseRotateReference ()
 Apply reverse rotations to all matrices in vector and fill new matrix with their sum. More...
 
void preselectLimitedDirections (double &phi, double &theta)
 
void preselectFastSignificant ()
 
void expectationSingleImage (Matrix1D< double > &opt_offsets)
 ML-integration over all (or -fast) translations. More...
 
void createThreads ()
 Create working threads. More...
 
void destroyThreads ()
 Exit threads and free memory. More...
 
int getThreadRefnoJob (int &refno)
 Assign refno jobs to threads. More...
 
void awakeThreads (ThreadTask task, int start_refno, int load=1)
 Awake threads for different tasks. More...
 
void doThreadRotateReferenceRefno ()
 Thread code to parallelize refno loop in rotateReference. More...
 
void doThreadReverseRotateReferenceRefno ()
 Thread code to parallelize refno loop in reverseRotateReference. More...
 
void doThreadPreselectFastSignificantRefno ()
 Thread code to parallelize refno loop in preselectFastSignificant. More...
 
void doThreadExpectationSingleImageRefno ()
 Thread code to parallelize refno loop in expectationSingleImage. More...
 
void doThreadESIUpdateRefno ()
 Thread code to parallelize update loop in ESI. More...
 
virtual void iteration ()
 Perform an iteration. More...
 
void maximizeModel (ModelML2D &model)
 Update all model parameters. More...
 
virtual void maximization ()
 Update all model parameters, adapted for IEM blocks use. More...
 
void correctScaleAverage ()
 Correct references scale. More...
 
virtual void addPartialDocfileData (const MultidimArray< double > &data, size_t first, size_t last)
 Add docfiledata to docfile. More...
 
virtual void readModel (ModelML2D &model, int block)
 Read model from file. More...
 
FileName getBaseName (String suffix="", int number=-1)
 Get base name based on fn_root and some number. More...
 
- Public Member Functions inherited from ML2DBaseProgram
void initSamplingStuff ()
 
 ML2DBaseProgram ()
 
virtual void randomizeImagesOrder ()
 Randomize initial images order, only once. More...
 
virtual bool checkConvergence ()
 
virtual void run ()
 Main function of the program. More...
 
virtual void defineBasicParams (XmippProgram *prog)
 
virtual void defineAdditionalParams (XmippProgram *prog, const char *sectionLine)
 
virtual void defineHiddenParams (XmippProgram *prog)
 
- 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)
 
virtual void initComments ()
 
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 (const String &param, int verb=2)
 
int version () const
 
virtual void show () 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 Member Functions inherited from MpiML2DBase< MpiProgML2D >
void readMpi (int argc, char **argv)
 
 MpiML2DBase (MpiProgML2D *prm)
 
 MpiML2DBase (MpiProgML2D *prm, MpiNode *node)
 
 MpiML2DBase (const MpiML2DBase &)=delete
 
 MpiML2DBase (const MpiML2DBase &&)=delete
 
 ~MpiML2DBase ()
 
MpiML2DBaseoperator= (const MpiML2DBase &)=delete
 
MpiML2DBaseoperator= (const MpiML2DBase &&)=delete
 
void sendDocfile (const MultidimArray< double > &data)
 

Additional Inherited Members

- Public Attributes inherited from ProgML2D
bool no_iem
 
MultidimArray< int > mask
 
MultidimArray< int > omask
 
int threadTask
 
barrier_t barrier
 
barrier_t barrier2
 
barrier_t barrier3
 
pthread_t * th_ids
 
structThreadTasksthreads_d
 
double LL
 
double sumfracweight
 
double wsum_sigma_noise
 
double wsum_sigma_offset
 
double sumw_allrefs
 
std::vector< double > sumw
 
std::vector< double > sumw2
 
std::vector< double > sumwsc
 
std::vector< double > sumwsc2
 
std::vector< double > sumw_mirror
 
std::vector< MultidimArray< double > > wsum_Mref
 
MultidimArray< int > Msignificant
 
MultidimArray< double > Mimg
 
std::vector< double > allref_offsets
 
std::vector< double > pdf_directions
 
std::vector< MultidimArray< double > > mref
 
std::vector< MultidimArray< std::complex< double > > > fref
 
std::vector< MultidimArray< std::complex< double > > > wsumimgs
 
int opt_refno
 
int iopt_psi
 
int iopt_flip
 
double trymindiff
 
double opt_scale
 
double bgmean
 
double opt_psi
 
double fracweight
 
double maxweight2
 
double dLL
 
std::vector< MultidimArray< std::complex< double > > > Fimg_flip
 
std::vector< MultidimArray< std::complex< double > > > mysumimgs
 
std::vector< double > refw
 
std::vector< double > refw2
 
std::vector< double > refwsc2
 
std::vector< double > refw_mirror
 
std::vector< double > sumw_refpsi
 
double wsum_corr
 
double sum_refw
 
double sum_refw2
 
double maxweight
 
double wsum_sc
 
double wsum_sc2
 
double wsum_offset
 
double old_bgmean
 
double mindiff
 
size_t sigdim
 
int ioptx
 
int iopty
 
std::vector< int > ioptx_ref
 
std::vector< int > iopty_ref
 
std::vector< int > ioptflip_ref
 
double pfs_mindiff
 
MultidimArray< double > pfs_maxweight
 
MultidimArray< double > pfs_weight
 
int pfs_count
 
std::vector< double > A2
 
double Xi2
 
size_t refno_index
 
size_t refno_load_param
 
int refno_load
 
int refno_count
 
int mygroup
 
- Public Attributes inherited from ML2DBaseProgram
FileName fn_img
 
FileName fn_ref
 
FileName fn_root
 
FileName fn_frac
 
FileName fn_sig
 
FileName fn_doc
 
FileName fn_oext
 
FileName fn_scratch
 
FileName fn_control
 
String cline
 
bool do_mirror
 
bool fix_fractions
 
bool fix_sigma_offset
 
bool fix_sigma_noise
 
int istart
 
int iter
 
int Niter
 
size_t dim
 
size_t dim2
 
size_t hdim
 
double ddim2
 
size_t nr_psi
 
size_t nr_flip
 
double psi_step
 
double psi_max
 
size_t nr_nomirror_flips
 
size_t nr_images_global
 
size_t nr_images_local
 
size_t myFirstImg
 
size_t myLastImg
 
double eps
 
MetaDataDb MDimg
 
MetaDataDb MDref
 
MetaDataDb MDlog
 
std::vector< Matrix2D< double > > F
 
std::vector< Image< double > > Iold
 
MultidimArray< double > P_phi
 
MultidimArray< double > Mr2
 
bool fast_mode
 
double C_fast
 
std::vector< Matrix1D< double > > Vtrans
 
bool zero_offsets
 
bool limit_rot
 
double search_rot
 
bool save_mem1
 
bool save_mem2
 
bool save_mem3
 
std::vector< double > imgs_oldphi
 
std::vector< double > imgs_oldtheta
 
bool do_ML3D
 
bool do_generate_refs
 
std::vector< std::vector< double > > imgs_offsets
 
double trymindiff_factor
 
double df
 
double df2
 
double dfsigma2
 
bool do_norm
 Re-normalize internally. More...
 
std::vector< double > imgs_scale
 
std::vector< double > imgs_bgmean
 
std::vector< double > imgs_trymindiff
 
std::vector< int > imgs_optrefno
 
double average_scale
 
int factor_nref
 
int refs_per_class
 
std::vector< double > conv
 
size_t current_image
 
ModelML2D model
 
ModelML2Dcurrent_model
 
size_t blocks
 
size_t current_block
 
std::vector< size_t > img_id
 
MultidimArray< double > docfiledata
 
bool do_restart
 
bool referenceExclusive
 
bool allowFastOption
 
bool allowThreads
 
bool allowIEM
 
bool allowRestart
 
String defaultRoot
 
int defaultNiter
 
int defaultStartingIter
 
int threads
 
String outRefsMd
 
MultidimArray< double > spectral_signal
 
int seed
 
- Public Attributes inherited from XmippProgram
bool doRun
 
bool runWithoutArgs
 
int verbose
 Verbosity level. More...
 
int debug
 
- Protected Member Functions inherited from XmippProgram
void defineCommons ()
 
- 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
 
- Protected Attributes inherited from MpiML2DBase< MpiProgML2D >
MpiNodenode
 
bool created_node
 
MpiProgML2Dprogram
 

Detailed Description

Class to parallelize the ML 2D alignment program

Definition at line 94 of file mpi_ml_align2d.h.

Constructor & Destructor Documentation

◆ MpiProgML2D() [1/2]

MpiProgML2D::MpiProgML2D ( )

Default constructor

Definition at line 86 of file mpi_ml_align2d.cpp.

86  :MpiML2DBase(this)
87 {}
MpiML2DBase(MpiProgML2D *prm)

◆ MpiProgML2D() [2/2]

MpiProgML2D::MpiProgML2D ( MpiNode node)

Constructor passing the MpiNode

Definition at line 89 of file mpi_ml_align2d.cpp.

89  :MpiML2DBase(this, node)
90 {}
MpiML2DBase(MpiProgML2D *prm)

Member Function Documentation

◆ endIteration()

void MpiProgML2D::endIteration ( )
virtual

Redefine endIteration for some synchronization.

Reimplemented from ML2DBaseProgram.

Definition at line 155 of file mpi_ml_align2d.cpp.

156 {
157  // Write output files
160 }
Definition: ml2d.h:76
void sendDocfile(const MultidimArray< double > &data)
ModelML2D model
Definition: ml2d.h:224
MultidimArray< double > docfiledata
Definition: ml2d.h:234
void writeOutputFiles(const ModelML2D &model, OutputType outputType=OUT_FINAL)
Write model parameters.

◆ expectation()

void MpiProgML2D::expectation ( )
virtual

After normal ML2D expectation, data must be collected from nodes.

Reimplemented from ProgML2D.

Definition at line 108 of file mpi_ml_align2d.cpp.

109 {
111  double aux;
112  Maux.resize(dim, dim);
113  Maux.setXmippOrigin();
114 
116  //After expectation, collect data from all nodes
117  // Here MPI_allreduce of all wsums,LL and sumfracweight !!!
118  MPI_Allreduce(&LL, &aux, 1, MPI_DOUBLE, MPI_SUM,
119  MPI_COMM_WORLD);
120  LL = aux;
121  MPI_Allreduce(&sumfracweight, &aux, 1, MPI_DOUBLE, MPI_SUM,
122  MPI_COMM_WORLD);
123  sumfracweight = aux;
124  MPI_Allreduce(&wsum_sigma_noise, &aux, 1, MPI_DOUBLE,
125  MPI_SUM, MPI_COMM_WORLD);
126  wsum_sigma_noise = aux;
127  MPI_Allreduce(&wsum_sigma_offset, &aux, 1, MPI_DOUBLE,
128  MPI_SUM, MPI_COMM_WORLD);
129  wsum_sigma_offset = aux;
130  for (int refno = 0; refno < model.n_ref * factor_nref; refno++)
131  {
132  MPI_Allreduce(MULTIDIM_ARRAY(wsum_Mref[refno]),
133  MULTIDIM_ARRAY(Maux),
134  MULTIDIM_SIZE(wsum_Mref[refno]), MPI_DOUBLE,
135  MPI_SUM, MPI_COMM_WORLD);
136  wsum_Mref[refno] = Maux;
137  MPI_Allreduce(&sumw[refno], &aux, 1, MPI_DOUBLE,
138  MPI_SUM, MPI_COMM_WORLD);
139  sumw[refno] = aux;
140  MPI_Allreduce(&sumwsc2[refno], &aux, 1, MPI_DOUBLE,
141  MPI_SUM, MPI_COMM_WORLD);
142  sumwsc2[refno] = aux;
143  MPI_Allreduce(&sumw_mirror[refno], &aux, 1, MPI_DOUBLE,
144  MPI_SUM, MPI_COMM_WORLD);
145  sumw_mirror[refno] = aux;
146  MPI_Allreduce(&sumw2[refno], &aux, 1, MPI_DOUBLE,
147  MPI_SUM, MPI_COMM_WORLD);
148  sumw2[refno] = aux;
149  MPI_Allreduce(&sumwsc[refno], &aux, 1, MPI_DOUBLE,
150  MPI_SUM, MPI_COMM_WORLD);
151  sumwsc[refno] = aux;
152  }
153 }//end of expectation
void resize(size_t Ndim, size_t Zdim, size_t Ydim, size_t Xdim, bool copy=true)
std::vector< double > sumw_mirror
Definition: ml_align2d.h:75
#define MULTIDIM_SIZE(v)
double wsum_sigma_offset
Definition: ml_align2d.h:74
ModelML2D model
Definition: ml2d.h:224
#define MULTIDIM_ARRAY(v)
int n_ref
Definition: ml2d.h:83
virtual void expectation()
Integrate over all experimental images.
double sumfracweight
Definition: ml_align2d.h:73
double LL
Definition: ml_align2d.h:73
std::vector< double > sumw
Definition: ml_align2d.h:75
double wsum_sigma_noise
Definition: ml_align2d.h:74
size_t dim
Definition: ml2d.h:151
int factor_nref
Definition: ml2d.h:215
std::vector< double > sumwsc2
Definition: ml_align2d.h:75
std::vector< double > sumwsc
Definition: ml_align2d.h:75
std::vector< MultidimArray< double > > wsum_Mref
Definition: ml_align2d.h:76
std::vector< double > sumw2
Definition: ml_align2d.h:75

◆ printModel()

void MpiProgML2D::printModel ( const String msg,
const ModelML2D model 
)
virtual

Reimplemented from ProgML2D.

Definition at line 176 of file mpi_ml_align2d.cpp.

177 {
178  if (node->isMaster())
179  ProgML2D::printModel(msg, model);
180 }
virtual void printModel(const String &msg, const ModelML2D &model)
Definition: ml_align2d.cpp:281
bool isMaster() const
Definition: xmipp_mpi.cpp:166

◆ produceSideInfo2()

void MpiProgML2D::produceSideInfo2 ( )
virtual

All mpi nodes should syncronize at this point that's why the need of override the implementation.

Reimplemented from ProgML2D.

Definition at line 100 of file mpi_ml_align2d.cpp.

101 {
102  node->barrierWait();
104  //Also sync after finishing produceSideInfo2
105  node->barrierWait();
106 }
virtual void produceSideInfo2()
Try to merge produceSideInfo1 and 2.
Definition: ml_align2d.cpp:357
void barrierWait()
Definition: xmipp_mpi.cpp:171

◆ setNumberOfLocalImages()

void MpiProgML2D::setNumberOfLocalImages ( )
virtual

Only take a part of images for process

Reimplemented from ML2DBaseProgram.

Definition at line 94 of file mpi_ml_align2d.cpp.

95 {
97  myLastImg);
98 }
size_t size
Definition: xmipp_mpi.h:52
size_t divide_equally(size_t N, size_t size, size_t rank, size_t &first, size_t &last)
size_t myLastImg
Definition: ml2d.h:168
size_t nr_images_local
Definition: ml2d.h:166
size_t myFirstImg
Definition: ml2d.h:168
size_t rank
Definition: xmipp_mpi.h:52
size_t nr_images_global
Definition: ml2d.h:164

◆ usage()

void MpiProgML2D::usage ( int  verb = 0) const
virtual

Print the usage of the program, reduced version

Reimplemented from XmippProgram.

Definition at line 182 of file mpi_ml_align2d.cpp.

183 {
184  if (node->isMaster())
185  ProgML2D::usage();
186 }
virtual void usage(int verb=0) const
bool isMaster() const
Definition: xmipp_mpi.cpp:166

◆ writeOutputFiles()

void MpiProgML2D::writeOutputFiles ( const ModelML2D model,
OutputType  outputType = OUT_FINAL 
)
virtual

Write model parameters.

Reimplemented from ProgML2D.

Definition at line 162 of file mpi_ml_align2d.cpp.

163 {
164  //All nodes should arrive to writeOutput files at same time
165  node->barrierWait();
166  //Only master write files
167  if (node->isMaster())
168  ProgML2D::writeOutputFiles(model, outputType);
169  else if (outputType == OUT_REFS)
171  //All nodes wait until files are written
172  node->barrierWait();
173 }
Definition: ml2d.h:76
FileName fn_root
Definition: ml2d.h:132
void barrierWait()
Definition: xmipp_mpi.cpp:171
#define FN_CLASSES_MD(base)
Definition: ml2d.h:55
FileName getIterExtraPath(const FileName &fn_root, int iter, bool makePath)
Definition: ml2d.cpp:305
virtual void writeOutputFiles(const ModelML2D &model, OutputType outputType=OUT_FINAL)
Write model parameters.
String outRefsMd
Definition: ml2d.h:250
bool isMaster() const
Definition: xmipp_mpi.cpp:166

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