Xmipp  v3.23.11-Nereus
mpi_ml_align2d.cpp
Go to the documentation of this file.
1 /***************************************************************************
2  * Authors: J.M. de la Rosa Trevin (jmdelarosa@cnb.csic.es)
3  *
4  *
5  * Unidad de Bioinformatica of Centro Nacional de Biotecnologia , CSIC
6  *
7  * This program is free software; you can redistribute it and/or modify
8  * it under the terms of the GNU General Public License as published by
9  * the Free Software Foundation; either version 2 of the License, or
10  * (at your option) any later version.
11  *
12  * This program is distributed in the hope that it will be useful,
13  * but WITHOUT ANY WARRANTY; without even the implied warranty of
14  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15  * GNU General Public License for more details.
16  *
17  * You should have received a copy of the GNU General Public License
18  * along with this program; if not, write to the Free Software
19  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA
20  * 02111-1307 USA
21  *
22  * All comments concerning this program package may be sent to the
23  * e-mail address 'xmipp@cnb.csic.es'
24  ***************************************************************************/
25 
26 #include "mpi_ml_align2d.h"
27 
28 /* Some constast to message passing tags */
29 constexpr int TAG_SEED = 1;
30 constexpr int TAG_DOCFILESIZE = 2;
31 constexpr int TAG_DOCFILE = 3;
32 
33 template<typename T>
35 {
36  node = nullptr;
37  program = prm;
38 }
39 
40 template<typename T>
42 {
43  node = mpinode;
44  created_node = false;
45  program = prm;
46 }
47 
48 template<typename T>
50 {
51 
52  // Write intermediate files
53  if (!node->isMaster())
54  {
55  // All slaves send docfile data to the master
56  int s_size = MULTIDIM_SIZE(docfiledata);
57  MPI_Send(&s_size, 1, MPI_INT, 0, TAG_DOCFILESIZE,
58  MPI_COMM_WORLD);
59  MPI_Send(MULTIDIM_ARRAY(docfiledata), s_size, MPI_DOUBLE,
60  0, TAG_DOCFILE, MPI_COMM_WORLD);
61  }
62  else
63  {
64  // Master fills docfile
65  // Master's own contribution
66  auto * ml2d = (ML2DBaseProgram*)program;
67  ml2d->addPartialDocfileData(docfiledata, ml2d->myFirstImg, ml2d->myLastImg);
68  int s_size;
69  size_t first_img, last_img;
70  MPI_Status status;
71 
72  for (size_t docCounter = 1; docCounter < node->size; ++docCounter)
73  {
74  // receive in order
75  MPI_Recv(&s_size, 1, MPI_INT, docCounter, TAG_DOCFILESIZE,
76  MPI_COMM_WORLD, &status);
77  MPI_Recv(MULTIDIM_ARRAY(docfiledata), s_size,
78  MPI_DOUBLE, docCounter, TAG_DOCFILE,
79  MPI_COMM_WORLD, &status);
80  divide_equally(ml2d->nr_images_global, node->size, docCounter, first_img, last_img);
81  ml2d->addPartialDocfileData(docfiledata, first_img, last_img);
82  }
83  }
84 }
85 
87 {}
88 
90 {}
91 
92 
93 
95 {
97  myLastImg);
98 }
99 
101 {
102  node->barrierWait();
104  //Also sync after finishing produceSideInfo2
105  node->barrierWait();
106 }
107 
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
154 
156 {
157  // Write output files
160 }
161 
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 }
174 
175 //Just for debuging
177 {
178  if (node->isMaster())
179  ProgML2D::printModel(msg, model);
180 }
181 
182 void MpiProgML2D::usage(int verb) const
183 {
184  if (node->isMaster())
185  ProgML2D::usage();
186 }
187 
188 #define SET_RANK_AND_SIZE() rank = node->rank; size = node->size;
189 
191 {}
192 
194 {}
195 
197 {
198  MultidimArray<double> Maux, Vaux;
199  double aux;
200  Maux.resize(dim, dim);
201  Maux.setXmippOrigin();
202 
204  MPI_Allreduce(&LL, &aux, 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
205  LL = aux;
206  MPI_Allreduce(&sumcorr, &aux, 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
207  sumcorr = aux;
208  MPI_Allreduce(&wsum_sigma_offset, &aux, 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
209  wsum_sigma_offset = aux;
210  if (do_kstest)
211  {
212  Vaux.resize(sumhist);
213  MPI_Allreduce(MULTIDIM_ARRAY(sumhist), MULTIDIM_ARRAY(Vaux),
214  MULTIDIM_SIZE(sumhist), MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
217  for (size_t ires = 0; ires < hdim; ires++)
218  {
219  Vaux.resize(sumhist);
220  MPI_Allreduce(MULTIDIM_ARRAY(resolhist[ires]), MULTIDIM_ARRAY(Vaux),
221  MULTIDIM_SIZE(resolhist[ires]), MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
224  }
225  }
226  for (int refno = 0;refno < model.n_ref; refno++)
227  {
228  MPI_Allreduce(MULTIDIM_ARRAY(wsum_Mref[refno]), MULTIDIM_ARRAY(Maux),
229  MULTIDIM_SIZE(wsum_Mref[refno]), MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
230  wsum_Mref[refno] = Maux;
231  if (do_ctf_correction)
232  {
233  MPI_Allreduce(MULTIDIM_ARRAY(wsum_ctfMref[refno]), MULTIDIM_ARRAY(Maux),
234  MULTIDIM_SIZE(wsum_ctfMref[refno]), MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
235  wsum_ctfMref[refno] = Maux;
236  }
237  MPI_Allreduce(&sumw[refno], &aux, 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
238  sumw[refno] = aux;
239  MPI_Allreduce(&sumw2[refno], &aux, 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
240  sumw2[refno] = aux;
241  MPI_Allreduce(&sumwsc2[refno], &aux, 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
242  sumwsc2[refno] = aux;
243  MPI_Allreduce(&sumwsc[refno], &aux, 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
244  sumwsc[refno] = aux;
245  MPI_Allreduce(&sumw_mirror[refno], &aux, 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
246  sumw_mirror[refno] = aux;
247  }
248  for (size_t ifocus = 0;ifocus < nr_focus;ifocus++)
249  {
250  for (size_t ii = 0; ii < Mwsum_sigma2[ifocus].size(); ii++)
251  {
252  MPI_Allreduce(&Mwsum_sigma2[ifocus][ii], &aux, 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
253  Mwsum_sigma2[ifocus][ii] = aux;
254  }
255  if (do_student)
256  {
257  MPI_Allreduce(&sumw_defocus[ifocus], &aux, 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
258  sumw_defocus[ifocus] = aux;
259  }
260  }
261 }//end of expectation
262 
264 {
265  node->barrierWait();
267  //Also sync after finishing produceSideInfo2
268  node->barrierWait();
269 }
270 
272 {
275 }
276 
277 
279 {
280  // Write output files
281  LOG("MpiProgMLF2D::endIteration : before sendDocfile");
283  LOG("MpiProgMLF2D::endIteration : before writeOutputFiles");
285  LOG("MpiProgMLF2D::endIteration : before updateWienerFilters");
287  LOG("MpiProgMLF2D::endIteration : before barrierWait");
288  node->barrierWait();
289 }
290 
292 {
293  //All nodes should arrive to writeOutput files at same time
294  LOG("MpiProgMLF2D::writeOutputFiles : waiting before writing");
295  node->barrierWait();
296  //Only master write files
297  if (node->isMaster())
298  {
299  ProgMLF2D::writeOutputFiles(model, outputType);
300  LOG("MpiProgMLF2D::writeOutputFiles : master writing ");
301  }
302  else if (outputType == OUT_REFS)
303  {
305  LOG(formatString("MpiProgMLF2D::writeOutputFiles : slave setting outRefsMd = %s", outRefsMd.c_str()).c_str());
306  }
307  //All nodes wait until files are written
308  LOG("MpiProgMLF2D::writeOutputFiles : waiting after writing");
309  node->barrierWait();
310 }
311 
312 
virtual void usage(int verb=0) const
bool do_student
IN DEVELOPMENT.
Definition: mlf_align2d.h:135
size_t size
Definition: xmipp_mpi.h:52
virtual void printModel(const String &msg, const ModelML2D &model)
Definition: ml_align2d.cpp:281
double wsum_sigma_offset
Definition: mlf_align2d.h:153
Definition: ml2d.h:76
Definition: ml2d.h:76
void resize(size_t Ndim, size_t Zdim, size_t Ydim, size_t Xdim, bool copy=true)
constexpr int TAG_DOCFILESIZE
std::vector< double > sumw_mirror
Definition: ml_align2d.h:75
#define MULTIDIM_SIZE(v)
void sendDocfile(const MultidimArray< double > &data)
double wsum_sigma_offset
Definition: ml_align2d.h:74
ModelML2D model
Definition: ml2d.h:224
double LL
Definition: mlf_align2d.h:153
size_t divide_equally(size_t N, size_t size, size_t rank, size_t &first, size_t &last)
constexpr int TAG_SEED
#define MULTIDIM_ARRAY(v)
size_t myLastImg
Definition: ml2d.h:168
int n_ref
Definition: ml2d.h:83
std::vector< double > sumw_defocus
Definition: mlf_align2d.h:156
OutputType
Definition: ml2d.h:76
constexpr int TAG_DOCFILE
void expectation()
Integrate over all experimental images.
virtual void produceSideInfo2()
Try to merge produceSideInfo1 and 2.
Definition: ml_align2d.cpp:357
void expectation()
After normal ML2D expectation, data must be collected from nodes.
size_t hdim
Definition: ml2d.h:151
FileName fn_root
Definition: ml2d.h:132
virtual void expectation()
Integrate over all experimental images.
void setNumberOfLocalImages()
double sumfracweight
Definition: ml_align2d.h:73
void barrierWait()
Definition: xmipp_mpi.cpp:171
#define FN_CLASSES_MD(base)
Definition: ml2d.h:55
size_t nr_focus
Definition: mlf_align2d.h:114
std::vector< double > sumwsc2
Definition: mlf_align2d.h:156
double LL
Definition: ml_align2d.h:73
size_t nr_images_local
Definition: ml2d.h:166
void endIteration()
Redefine endIteration for some synchronization.
MultidimArray< double > spectral_signal
Definition: ml2d.h:253
size_t myFirstImg
Definition: ml2d.h:168
MpiNode * node
void produceSideInfo2()
bool do_kstest
Statistical analysis of the noise distributions.
Definition: mlf_align2d.h:141
std::vector< double > sumw
Definition: ml_align2d.h:75
double wsum_sigma_noise
Definition: ml_align2d.h:74
std::vector< MultidimArray< double > > wsum_Mref
Definition: mlf_align2d.h:154
FileName getIterExtraPath(const FileName &fn_root, int iter, bool makePath)
Definition: ml2d.cpp:305
#define FOR_ALL_DIRECT_ELEMENTS_IN_MULTIDIMARRAY(v)
MultidimArray< double > docfiledata
Definition: ml2d.h:234
virtual void produceSideInfo2()
#define DIRECT_MULTIDIM_ELEM(v, n)
size_t dim
Definition: ml2d.h:151
int factor_nref
Definition: ml2d.h:215
bool do_ctf_correction
Definition: mlf_align2d.h:100
std::vector< Histogram1D > resolhist
Definition: mlf_align2d.h:146
virtual void writeOutputFiles(const ModelML2D &model, OutputType outputType=OUT_FINAL)
Write model parameters.
std::vector< double > sumwsc2
Definition: ml_align2d.h:75
std::vector< double > sumwsc
Definition: mlf_align2d.h:156
size_t rank
Definition: xmipp_mpi.h:52
#define LOG(msg)
std::vector< double > sumw
Definition: mlf_align2d.h:156
void produceSideInfo2()
std::vector< double > sumw_mirror
Definition: mlf_align2d.h:156
std::vector< double > sumw2
Definition: mlf_align2d.h:156
void writeOutputFiles(const ModelML2D &model, OutputType outputType=OUT_FINAL)
Write model parameters.
MpiML2DBase(T *prm)
std::vector< double > sumwsc
Definition: ml_align2d.h:75
void produceSideInfo()
std::vector< MultidimArray< double > > wsum_Mref
Definition: ml_align2d.h:76
std::string String
Definition: xmipp_strings.h:34
std::vector< std::vector< double > > Mwsum_sigma2
Definition: mlf_align2d.h:155
String formatString(const char *format,...)
Definition: ml2d.h:79
void expectation()
After normal ML2D expectation, data must be collected from nodes.
virtual void usage(int verb=0) const
void endIteration()
Redefine endIteration for some synchronization.
String outRefsMd
Definition: ml2d.h:250
bool isMaster() const
Definition: xmipp_mpi.cpp:166
size_t nr_images_global
Definition: ml2d.h:164
ProgClassifyCL2D * prm
void writeOutputFiles(const ModelML2D &model, OutputType outputType=OUT_FINAL)
Write model parameters.
Histogram1D sumhist
Definition: mlf_align2d.h:145
void updateWienerFilters(const MultidimArray< double > &spectral_signal, std::vector< double > &sumw_defocus, int iter)
#define SET_RANK_AND_SIZE()
int * n
void printModel(const String &msg, const ModelML2D &model)
virtual void writeOutputFiles(const ModelML2D &model, OutputType outputType)
Write out reference images, selfile and logfile.
double sumcorr
Definition: mlf_align2d.h:153
virtual void produceSideInfo()
Setup lots of stuff.
std::vector< double > sumw2
Definition: ml_align2d.h:75
std::vector< MultidimArray< double > > wsum_ctfMref
Definition: mlf_align2d.h:154