Xmipp  v3.23.11-Nereus
Public Member Functions | Public Attributes | List of all members
ProgMpiAngularProjectLibrary Class Reference
Inheritance diagram for ProgMpiAngularProjectLibrary:
Inheritance graph
[legend]
Collaboration diagram for ProgMpiAngularProjectLibrary:
Collaboration graph
[legend]

Public Member Functions

 ProgMpiAngularProjectLibrary ()
 
void readParams ()
 
void defineParams ()
 
void show ()
 
void preRun ()
 
void process ()
 
void createGroupSamplingFiles (void)
 
void error_exit (const char *msg)
 
void run ()
 
- Public Member Functions inherited from ProgAngularProjectLibrary
 ProgAngularProjectLibrary ()
 
 ProgAngularProjectLibrary (const ProgAngularProjectLibrary &)=delete
 
 ProgAngularProjectLibrary (const ProgAngularProjectLibrary &&)=delete
 
 ~ProgAngularProjectLibrary ()
 
ProgAngularProjectLibraryoperator= (const ProgAngularProjectLibrary &)=delete
 
ProgAngularProjectLibraryoperator= (const ProgAngularProjectLibrary &&)=delete
 
void readParams ()
 
void defineParams ()
 
void show ()
 
void run ()
 
void createGroupSamplingFiles (void)
 
void get_sym_vectors (std::vector< Matrix1D< double > > &sym_points)
 
void project_angle_vector (int my_init, int my_end, bool verbose=true)
 
- 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 (int verb=0) 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 Attributes

int nProcs
 
int mpi_job_size
 
int numberOfJobs
 
int rank
 
MPI_Status status
 
- Public Attributes inherited from ProgAngularProjectLibrary
Sampling mysampling
 
double sampling
 
double perturb_projection_vector
 
double psi_sampling
 
double max_tilt_angle
 
double min_tilt_angle
 
FileName output_file
 
FileName output_file_root
 
FileName fn_groups
 
projectionType projType
 Type of projection algorithm: More...
 
double paddFactor
 The padding factor for Fourier projection. More...
 
double maxFrequency
 The maximum frequency for Fourier projection. More...
 
int BSplineDeg
 The type of interpolation (NEAR. More...
 
FileName FnexperimentalImages
 
bool angular_distance_bool
 
bool remove_points_far_away_from_experimental_data_bool
 
double angular_distance
 
bool compute_closer_sampling_point_bool
 
bool compute_neighbors_bool
 
int symmetry
 
FileName fn_sym
 symmetry file for sampling More...
 
FileName fn_sym_neigh
 symmetry file for heighbors computation More...
 
int sym_order
 
FileName input_volume
 
Image< double > inputVol
 
int Xdim
 
int Ydim
 
bool only_winner
 
RealShearsInfoVshears
 
FourierProjectorVfourier
 
- Public Attributes inherited from XmippProgram
bool doRun
 
bool runWithoutArgs
 
int verbose
 Verbosity level. More...
 
int debug
 

Additional Inherited Members

- 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
 

Detailed Description

Definition at line 46 of file mpi_angular_project_library.cpp.

Constructor & Destructor Documentation

◆ ProgMpiAngularProjectLibrary()

ProgMpiAngularProjectLibrary::ProgMpiAngularProjectLibrary ( )
inline

Definition at line 68 of file mpi_angular_project_library.cpp.

69  {
70  //parent class constructor will be called by deault without parameters
71  MPI_Comm_size(MPI_COMM_WORLD, &nProcs);
72  MPI_Comm_rank(MPI_COMM_WORLD, &rank);
73  //Blocks until all process have reached this routine.
74  //very likelly this is
75  MPI_Barrier(MPI_COMM_WORLD);
76  }

Member Function Documentation

◆ createGroupSamplingFiles()

void ProgMpiAngularProjectLibrary::createGroupSamplingFiles ( void  )
inline

◆ defineParams()

void ProgMpiAngularProjectLibrary::defineParams ( )
inlinevirtual

Function in which the param of each Program are defined.

Reimplemented from XmippProgram.

Definition at line 89 of file mpi_angular_project_library.cpp.

90  {
92  addParamsLine(" [--mpi_job_size <size=10>]: Number of images sent to a cpu in a single job");
93  addParamsLine(" : 10 may be a good value");
94  addParamsLine(" : if -1 the computer will put the maximum");
95  addParamsLine(" : posible value that may not be the best option");
96  }
void addParamsLine(const String &line)

◆ error_exit()

void ProgMpiAngularProjectLibrary::error_exit ( const char *  msg)
inline

Definition at line 488 of file mpi_angular_project_library.cpp.

489  {
490  fprintf(stderr, "%s", msg);
491  MPI_Abort(MPI_COMM_WORLD, EXIT_FAILURE);
492  }
fprintf(glob_prnt.io, "\)

◆ preRun()

void ProgMpiAngularProjectLibrary::preRun ( )
inline

Definition at line 107 of file mpi_angular_project_library.cpp.

108  {
109  //do not use rank 0 here, since rank 1 will manipulate this file later
110  if (rank == 1)
111  {
112  unlink(output_file.c_str());
113  }
114 
115 
116  //#define DEBUGTIME
117 #ifdef DEBUGTIME
118 #include <ctime>
119 
120  time_t start,end;
121  double time_dif;
122  time (&start);
123  time (&end);
124  time_dif = difftime (end,start);
125  start=end;
126  std::cerr<<" starting prerun rank= "<<rank<<std::endl;
127 #endif
128 
129  int my_seed;
130  if (rank == 0)
131  {
132  show();
133  //random numbers must be the same in all nodes
135  {
136  srand(getpid());
137  my_seed = rand();
138  }
139  }
140 #ifdef DEBUGTIME
141  time (&end);
142  time_dif = difftime (end,start);
143  start=end;
144  std::cerr<<" set rand seed rank= "<<rank<<std::endl;
145 #endif
146  //Bcast must be seem by all processors
148  {
149  MPI_Bcast (&my_seed, 1, MPI_INT, 0, MPI_COMM_WORLD);
151 #ifdef DEBUGTIME
152 
153  time (&end);
154  time_dif = difftime (end,start);
155  start=end;
156  std::cerr<<" after perturb rank= "<<rank<<std::endl;
157 #endif
158 
159  }
160  //all ranks
162  //symmetry for sampling may be different from neighbourhs
164  REPORT_ERROR(ERR_NUMERICAL, (std::string)"angular_project_library::run Invalid symmetry" + fn_sym);//set sampling must go before set noise
165  if(angular_distance_bool!=0)
167 
168 #ifdef DEBUGTIME
169 
170  time (&end);
171  time_dif = difftime (end,start);
172  start=end;
173  std::cerr<<" setsampling rank= "<<rank<<std::endl;
174 #endif
175  //true -> half_sphere
178  //store symmetry matrices, this is faster than computing them each time
180 
181 #ifdef DEBUGTIME
182 
183  time (&end);
184  time_dif = difftime (end,start);
185  start=end;
186  std::cerr<<" compute sampling points rank= "<<rank<<std::endl;
187 #endif
188 
189  // We first sample The whole sphere
190  // Then we remove point redundant due to sampling symmetry
191  // use old symmetry, this is geometric does not use L_R
193 
194  //=========================
195  //======================
196  //recompute symmetry with neigh symmetry
198  REPORT_ERROR(ERR_NUMERICAL, (std::string)"angular_project_library::run Invalid neig symmetry" + fn_sym_neigh);
201  //precompute product between symmetry matrices and experimental data
202  if (FnexperimentalImages.size() > 0)
204 #ifdef DEBUGTIME
205 
206  time (&end);
207  time_dif = difftime (end,start);
208  start=end;
209  std::cerr<<" remove redundant rank= "<<rank<<std::endl;
210 #endif
211 
212  if (FnexperimentalImages.size() > 0 &&
214  {
216 #ifdef DEBUGTIME
217 
218  time (&end);
219  time_dif = difftime (end,start);
220  start=end;
221  std::cerr<<" remove points far away rank= "<<rank<<std::endl;
222 #endif
223 
224  }
225 
226  /* save files */
227  if (rank == 0)
228  {
230  {
231  //find sampling point closer to experimental point (only 0) and bool
232  //and save docfile with this information
234  }
235  //mysampling.createSymFile(symmetry, sym_order);
237  }
238 
239  if (rank != 0)
240  {
242  inputVol().setXmippOrigin();
243  Xdim = XSIZE(inputVol());
244  Ydim = YSIZE(inputVol());
245  }
246  if (rank == 0)
247  {
249  {
252  }
253  }
254  //release some memory
256 
257  if (mpi_job_size != -1)
258  {
260  }
261  else
262  {
263  numberOfJobs=nProcs-1;//one node is the master
265  }
266 
267  verbose=false;
268  #define DEBUG
269 #ifdef DEBUG
270 
271  if (rank == 1)
272  {
273  std::cerr << "numberOfJobs: " << numberOfJobs << std::endl
274  << "number of projections to be created: " << mysampling.no_redundant_sampling_points_angles.size()
275  <<std::endl;
276  }
277 #endif
278  #undef DEBUG
279 
280  //create blanck outfile so main header does not need to be re-writen
281  //rank 1 has already xdim and ydim.
282  if (rank == 1)
283  {
284  int numberProjections = std::ceil(360.0 / psi_sampling);
285  numberProjections *= mysampling.no_redundant_sampling_points_angles.size();
286  std::cerr << "creating Blank file: "
287  << output_file << " for "
288  << numberProjections << " projections of size: " << Xdim << " " << Ydim
289  << std::endl;
290  createEmptyFile(output_file,Xdim,Ydim,1,numberProjections,true,WRITE_REPLACE);
291  }
292 
293  MPI_Barrier(MPI_COMM_WORLD);
294  }
#define YSIZE(v)
void setSampling(double sampling)
Definition: sampling.cpp:121
void removeRedundantPoints(const int symmetry, int sym_order)
Definition: sampling.cpp:691
void setNoise(double deviation, int my_seed=-1)
Definition: sampling.cpp:134
#define REPORT_ERROR(nerr, ErrormMsg)
Definition: xmipp_error.h:211
void createAsymUnitFile(const FileName &docfilename)
Definition: sampling.cpp:1440
void setNeighborhoodRadius(double neighborhood)
Definition: sampling.cpp:140
void findClosestSamplingPoint(const MetaData &DFi, const FileName &output_file_root)
Definition: sampling.cpp:1991
int readSymmetryFile(FileName fn_sym, double accuracy=SYM_ACCURACY)
Definition: symmetries.cpp:33
bool isSymmetryGroup(FileName fn_sym, int &pgGroup, int &pgOrder)
Definition: symmetries.cpp:601
FileName fn_sym
symmetry file for sampling
void computeNeighbors(bool only_winner=false)
Definition: sampling.cpp:1715
void saveSamplingFile(const FileName &fn_base, bool write_vectors=true, bool write_sampling_sphere=false)
Definition: sampling.cpp:1495
FileName fn_sym_neigh
symmetry file for heighbors computation
std::vector< Matrix1D< double > > exp_data_projection_direction_by_L_R
Definition: sampling.h:108
void computeSamplingPoints(bool only_half_sphere=true, double max_tilt=180, double min_tilt=0)
Definition: sampling.cpp:155
#define XSIZE(v)
Error related to numerical calculation.
Definition: xmipp_error.h:179
void createEmptyFile(const FileName &filename, int xdim, int ydim, int Zdim, size_t select_img, bool isStack, int mode, int _swapWrite, const MDRowVec *md)
int verbose
Verbosity level.
void fillLRRepository(void)
Definition: sampling.cpp:2216
std::vector< Matrix1D< double > > no_redundant_sampling_points_angles
Definition: sampling.h:121
int read(const FileName &name, DataMode datamode=DATA, size_t select_img=ALL_IMAGES, bool mapData=false, int mode=WRITE_READONLY)
void fillExpDataProjectionDirectionByLR(const MetaData &DFi)
Definition: sampling.cpp:2256
void removePointsFarAwayFromExperimentalData()
Definition: sampling.cpp:1928
SymList SL
Definition: sampling.h:138

◆ process()

void ProgMpiAngularProjectLibrary::process ( )
inline

Definition at line 296 of file mpi_angular_project_library.cpp.

297  {
298  if (rank == 0)
299  {
300  int stopTagsSent =0;
301  int c = XMIPP_MAX(1, numberOfJobs / 60);
303  for (int i=0;i<numberOfJobs;)
304  {
305  //collect data if available
306  //be aware that mpi_Probe will block the program until a message is received
307  MPI_Probe(MPI_ANY_SOURCE, MPI_ANY_TAG, MPI_COMM_WORLD, &status);
308  // worker is free
309  if (status.MPI_TAG == TAG_FREEWORKER)
310  {
311  MPI_Recv(nullptr, 0, MPI_INT, MPI_ANY_SOURCE, TAG_FREEWORKER,
312  MPI_COMM_WORLD, &status);
313  //#define DEBUG
314 #ifdef DEBUG
315 
316  std::cerr << "Mr_f received TAG_FREEWORKER from worker " << status.MPI_SOURCE << std::endl;
317 #endif
318  //send work
319  MPI_Send(&i,
320  1,
321  MPI_INT,
322  status.MPI_SOURCE,
324  MPI_COMM_WORLD);
325  i++; //increase job number
326  //#define DEBUG
327 #ifdef DEBUG
328 
329  std::cerr << "Ms_f sent TAG_WORKFORWORKER to worker " << status.MPI_SOURCE << std::endl;
330  std::cerr << "Sent jobNo " << i << std::endl;
331 #endif
332  //#undef DEBUG
333  }
334  else
335  {
336  std::cerr << "M_f Received unknown TAG" << std::endl;
337  exit(0);
338  }
339 
340  if (i % c == 0)
341  progress_bar(i);
342  }
343  progress_bar(numberOfJobs);
344 
345  //send TAG_STOP
346  while (stopTagsSent < (nProcs-1))
347  {
348  MPI_Recv(nullptr, 0, MPI_INT, MPI_ANY_SOURCE, TAG_FREEWORKER,
349  MPI_COMM_WORLD, &status);
350 #ifdef DEBUG
351 
352  std::cerr << "Mr received TAG_FREEWORKER from worker " << status.MPI_SOURCE << std::endl;
353  std::cerr << "Ms sent TAG_STOP to worker" << status.MPI_SOURCE << std::endl;
354 #endif
355  #undef DEBUG
356 
357  MPI_Send(nullptr, 0, MPI_INT, status.MPI_SOURCE, TAG_STOP, MPI_COMM_WORLD);
358  stopTagsSent++;
359  }
360  //only rank 0 create sel file
361  if(rank==0)
362  {
363 
364  MetaDataVec mySFin;
365  mySFin.read(output_file_root + "_angles.doc");
366 #define ANGLESDOC
367 #ifdef ANGLESDOC
368 // std::cerr << "DEBUG_ROB, output_file_root + angle.doc:"
369 // << output_file_root + "_angles.doc" << std::endl;
370 #endif
371  MetaDataVec mySF;
372  FileName fn_temp;
373 
374  size_t myCounter = 0;
375  size_t id;
376  int ref;
377 
378  for (size_t objId : mySFin.ids())
379  {
380  // FIXME: read whole row at once
381  double x,y,z, rot, tilt, psi;
382  mySFin.getValue(MDL_ANGLE_ROT,rot,objId);
383  mySFin.getValue(MDL_ANGLE_TILT,tilt,objId);
384  mySFin.getValue(MDL_ANGLE_PSI,psi,objId);
385  mySFin.getValue(MDL_X,x,objId);
386  mySFin.getValue(MDL_Y,y,objId);
387  mySFin.getValue(MDL_Z,z,objId);
388  mySFin.getValue(MDL_REF,ref,objId);
389 
390  for (int i = 0; i < std::ceil(360.0 / psi_sampling); ++i) {
391  //FIXME, do I have order?
392  fn_temp.compose(++myCounter,output_file);
393  id = mySF.addObject();
394  mySF.setValue(MDL_IMAGE,fn_temp, id);
395  mySF.setValue(MDL_ENABLED,1, id);
396 
397  mySF.setValue(MDL_ANGLE_ROT,rot, id);
398  mySF.setValue(MDL_ANGLE_TILT,tilt, id);
399  mySF.setValue(MDL_ANGLE_PSI,psi + i * psi_sampling, id);
400  mySF.setValue(MDL_X,x, id);
401  mySF.setValue(MDL_Y,y, id);
402  mySF.setValue(MDL_Z,z, id);
403  mySF.setValue(MDL_SCALE,1.0,id);
404  mySF.setValue(MDL_REF,ref,id);
405  }
406  }
407  fn_temp=output_file_root+".doc";
408  mySF.setComment("x,y,z refer to the coordinates of the unitary vector at direction given by the euler angles");
409  mySF.write(fn_temp);
410 #ifndef ANGLESDOC
411  unlink((output_file_root+"_angles.doc").c_str());
412 #endif
413 #undef ANGLESDOC
414  }
415  }
416  else
417  {
418  // Select only relevant part of selfile for this rank
419  // job number
420  // job size
421  // aux variable
422  while (1)
423  {
424  int jobNumber;
425  //I am free
426  MPI_Send(nullptr, 0, MPI_INT, 0, TAG_FREEWORKER, MPI_COMM_WORLD);
427  //#define DEBUG
428 #ifdef DEBUG
429 
430  std::cerr << "W" << rank << " " << "sent TAG_FREEWORKER to master " << std::endl;
431 #endif
432 #undef DEBUG
433  //get yor next task
434  MPI_Probe(0, MPI_ANY_TAG, MPI_COMM_WORLD, &status);
435 #ifdef DEBUG
436 
437  std::cerr << "W" << rank << " " << "probe MPI_ANY_TAG " << std::endl;
438 #endif
439 
440  if (status.MPI_TAG == TAG_STOP)//no more jobs exit
441  {
442  //If I do not read this tag
443  //master will no further process
444  //a posibility is a non-blocking send
445  MPI_Recv(nullptr, 0, MPI_INT, 0, TAG_STOP,
446  MPI_COMM_WORLD, &status);
447 #ifdef DEBUG
448 
449  std::cerr << "Wr" << rank << " " << "TAG_STOP" << std::endl;
450 #endif
451 
452  break;
453  }
454  if (status.MPI_TAG == TAG_WORKFORWORKER)
455  //there is still some work to be done
456  {
457  //get the jobs number
458  MPI_Recv(&jobNumber, 1, MPI_INT, 0, TAG_WORKFORWORKER, MPI_COMM_WORLD, &status);
459 #ifdef DEBUG
460 
461  std::cerr << "Wr" << rank << " " << "TAG_WORKFORWORKER" << std::endl;
462  std::cerr << "jobNumber " << jobNumber << std::endl;
463 #endif
464 #undef DEBUG
465  // Process all images
467  std::min((size_t)((jobNumber+1)* mpi_job_size -1),
469  verbose);
470  //get yor next task
471  }
472  else
473  {
474  std::cerr << "3) Received unknown TAG I quit" << std::endl;
475  exit(0);
476  }
477  }
478  }
479  }
void init_progress_bar(long total)
Rotation angle of an image (double,degrees)
void min(Image< double > &op1, const Image< double > &op2)
#define XMIPP_MAX(x, y)
Definition: xmipp_macros.h:193
void read(const FileName &inFile, const std::vector< MDLabel > *desiredLabels=nullptr, bool decomposeStack=true) override
doublereal * c
Z component (double)
Tilting angle of an image (double,degrees)
static double * y
constexpr int TAG_WORKFORWORKER
void compose(const String &str, const size_t no, const String &ext="")
Special label to be used when gathering MDs in MpiMetadataPrograms.
void write(const FileName &outFile, WriteModeMetaData mode=MD_OVERWRITE) const
virtual IdIteratorProxy< false > ids()
doublereal * x
#define i
Is this image enabled? (int [-1 or 1])
virtual void setComment(const String &newComment="No comment")
constexpr int TAG_STOP
bool setValue(const MDObject &mdValueIn, size_t id)
size_t addObject() override
void project_angle_vector(int my_init, int my_end, bool verbose=true)
void progress_bar(long rlen)
double z
scaling factor for an image or volume (double)
int verbose
Verbosity level.
bool getValue(MDObject &mdValueOut, size_t id) const override
Class to which the image belongs (int)
std::vector< Matrix1D< double > > no_redundant_sampling_points_angles
Definition: sampling.h:121
Y component (double)
double psi(const double x)
X component (double)
Name of an image (std::string)
constexpr int TAG_FREEWORKER

◆ readParams()

void ProgMpiAngularProjectLibrary::readParams ( )
inlinevirtual

Function in which each program will read parameters that it need. If some error occurs the usage will be printed out.

Reimplemented from XmippProgram.

Definition at line 80 of file mpi_angular_project_library.cpp.

81  {
82  if (nProcs < 2)
83  error_exit("This program cannot be executed in a single working node");
85  mpi_job_size = getIntParam("--mpi_job_size");
86  }
int getIntParam(const char *param, int arg=0)

◆ run()

void ProgMpiAngularProjectLibrary::run ( )
inlinevirtual

This function will be start running the program. it also should be implemented by derived classes.

Reimplemented from XmippProgram.

Definition at line 494 of file mpi_angular_project_library.cpp.

495  {
496  try
497  {
498  preRun();
499  process();
501  }
502  catch (XmippError &XE)
503  {
504  std::cerr << "Error!" <<std::endl;
505  std::cerr << XE.what();
506  MPI_Finalize();
507  exit(1);
508  }
509 
510  if (rank == 0)
511  std::cout << "Done!" <<std::endl;
512  }

◆ show()

void ProgMpiAngularProjectLibrary::show ( )
inline

Definition at line 100 of file mpi_angular_project_library.cpp.

101  {
103  std::cerr << " Size of mpi jobs " << mpi_job_size <<std::endl;
104  }

Member Data Documentation

◆ mpi_job_size

int ProgMpiAngularProjectLibrary::mpi_job_size

Dvide the job in this number block with this number of images

Definition at line 56 of file mpi_angular_project_library.cpp.

◆ nProcs

int ProgMpiAngularProjectLibrary::nProcs

Number of Procesors

Definition at line 53 of file mpi_angular_project_library.cpp.

◆ numberOfJobs

int ProgMpiAngularProjectLibrary::numberOfJobs

Number of jobs

Definition at line 59 of file mpi_angular_project_library.cpp.

◆ rank

int ProgMpiAngularProjectLibrary::rank

computing node number. Master=0

Definition at line 62 of file mpi_angular_project_library.cpp.

◆ status

MPI_Status ProgMpiAngularProjectLibrary::status

status after am MPI call

Definition at line 65 of file mpi_angular_project_library.cpp.


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