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

Public Member Functions

 ProgSortImages (int argc, char **argv)
 MPI constructor. More...
 
 ProgSortImages (const ProgSortImages &)=delete
 
 ProgSortImages (const ProgSortImages &&)=delete
 
 ~ProgSortImages ()
 MPI destructor. More...
 
ProgSortImagesoperator= (const ProgSortImages &)=delete
 
ProgSortImagesoperator= (const ProgSortImages &&)=delete
 
void readParams ()
 Read argument from command line. More...
 
void show ()
 Show. More...
 
void defineParams ()
 Usage. More...
 
void produceSideInfo (int rank)
 Produce side info. More...
 
void chooseNextImage (int rank, int Nproc)
 Choose next image. More...
 
void run ()
 Main routine. 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)
 
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

FileName fnSel
 
FileName fnRoot
 
MpiNodenode
 
FileName fnStack
 
bool center
 
MetaDataVec SFout
 
std::vector< MDRowVectoClassify
 
Matrix1D< int > stillToDo
 
Image< double > lastImage
 
MultidimArray< int > mask
 
- 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 35 of file mpi_image_sort.cpp.

Constructor & Destructor Documentation

◆ ProgSortImages() [1/3]

ProgSortImages::ProgSortImages ( int  argc,
char **  argv 
)
inline

MPI constructor.

Definition at line 68 of file mpi_image_sort.cpp.

69  {
70  node=new MpiNode(argc,argv);
71  if (!node->isMaster())
72  verbose=0;
73  }
int argc
Original command line arguments.
Definition: xmipp_program.h:86
const char ** argv
Definition: xmipp_program.h:87
int verbose
Verbosity level.
bool isMaster() const
Definition: xmipp_mpi.cpp:166

◆ ProgSortImages() [2/3]

ProgSortImages::ProgSortImages ( const ProgSortImages )
delete

◆ ProgSortImages() [3/3]

ProgSortImages::ProgSortImages ( const ProgSortImages &&  )
delete

◆ ~ProgSortImages()

ProgSortImages::~ProgSortImages ( )
inline

MPI destructor.

Definition at line 79 of file mpi_image_sort.cpp.

80  {
81  delete node;
82  }

Member Function Documentation

◆ chooseNextImage()

void ProgSortImages::chooseNextImage ( int  rank,
int  Nproc 
)
inline

Choose next image.

Definition at line 194 of file mpi_image_sort.cpp.

195  {
196  Image<double> bestImage, I;
198  double bestCorr=-1;
199  int bestIdx=-1;
200  int count=0;
201  FileName fnImageStack, fnImg;
202  AlignmentAux aux;
203  CorrelationAux aux2;
206  {
207  if (!VEC_ELEM(stillToDo,i))
208  continue;
209  if ((count+1)%Nproc!=rank)
210  {
211  ++count;
212  continue;
213  }
214  toClassify[i].getValue(MDL_IMAGE,fnImg);
215  I.read(fnImg);
216  I().setXmippOrigin();
217  double corr=alignImagesConsideringMirrors(lastImage(),I(),M,aux,aux2,aux3,&mask);
218  if (center)
219  centerImage(I(),aux2,aux3);
220  if (corr>bestCorr)
221  {
222  bestCorr=corr;
223  bestImage()=I();
224  bestIdx=i;
225  }
226  ++count;
227  }
228 
229  // Rank 0 receives from the other nodes their best image
230  double buffer[2];
231  if (rank==0)
232  {
233  MPI_Status status;
234  for (int n=1; n<Nproc; ++n)
235  {
236  MPI_Recv(buffer, 2, MPI_DOUBLE, n, 0, MPI_COMM_WORLD, &status);
237  if (buffer[1]>bestCorr)
238  {
239  bestIdx=(int)buffer[0];
240  bestCorr=buffer[1];
241  }
242  }
243  std::cout << "Images to go=" << stillToDo.sum()-1 << " current correlation= " << bestCorr << std::endl;
244  }
245  else
246  {
247  buffer[0]=bestIdx;
248  buffer[1]=bestCorr;
249  MPI_Send(buffer, 2, MPI_DOUBLE, 0, 0, MPI_COMM_WORLD);
250  }
251 
252  // Now rank 0 redistributes the best image
253  if (rank==0)
254  {
255  buffer[0]=bestIdx;
256  MPI_Bcast(buffer, 1, MPI_DOUBLE, 0, MPI_COMM_WORLD);
257  }
258  else
259  {
260  MPI_Bcast(buffer, 1, MPI_DOUBLE, 0, MPI_COMM_WORLD);
261  bestIdx=(int)buffer[0];
262  }
263 
264  // All compute the best image
265  toClassify[bestIdx].getValue(MDL_IMAGE,fnImg);
266  I.read(fnImg);
267  I().setXmippOrigin();
268  bestCorr=alignImagesConsideringMirrors(lastImage(),I(),M,aux,aux2,aux3,&mask);
269  if (center)
270  centerImage(I(),aux2,aux3);
271  lastImage()=bestImage()=I();
272 
273  if (rank==0)
274  {
275  int idxStack=SFout.size();
276  fnImageStack.compose(idxStack+1,fnStack);
277  bestImage.write(fnStack,idxStack,true,WRITE_APPEND);
278  MDRowVec &row = toClassify[bestIdx];
279  FileName fnImgOrig;
280  row.getValue(MDL_IMAGE,fnImgOrig);
281  row.setValue(MDL_IMAGE_ORIGINAL,fnImgOrig);
282  row.setValue(MDL_IMAGE,fnImageStack);
283  row.setValue(MDL_MAXCC,bestCorr);
284  SFout.addRow(row);
285  }
286  VEC_ELEM(stillToDo,bestIdx)=0;
287  }
#define VEC_ELEM(v, i)
Definition: matrix1d.h:245
HBITMAP buffer
Definition: svm-toy.cpp:37
void setValue(const MDObject &object) override
void write(const FileName &name="", size_t select_img=ALL_IMAGES, bool isStack=false, int mode=WRITE_OVERWRITE, CastWriteMode castMode=CW_CAST, int _swapWrite=0)
void compose(const String &str, const size_t no, const String &ext="")
size_t size() const override
#define i
size_t addRow(const MDRow &row) override
T & getValue(MDLabel label)
Matrix2D< double > centerImage(MultidimArray< double > &I, CorrelationAux &aux, RotationalCorrelationAux &aux2, int Niter, bool limitShift)
Definition: filters.cpp:3277
#define FOR_ALL_ELEMENTS_IN_MATRIX1D(v)
Definition: matrix1d.h:72
Image< double > lastImage
Name of an image from which MDL_IMAGE is coming from.
double sum(bool average=false) const
Definition: matrix1d.cpp:652
Maximum cross-correlation for the image (double)
Matrix1D< int > stillToDo
MultidimArray< int > mask
std::vector< MDRowVec > toClassify
int read(const FileName &name, DataMode datamode=DATA, size_t select_img=ALL_IMAGES, bool mapData=false, int mode=WRITE_READONLY)
double alignImagesConsideringMirrors(const MultidimArray< double > &Iref, const AlignmentTransforms &IrefTransforms, MultidimArray< double > &I, Matrix2D< double > &M, AlignmentAux &aux, CorrelationAux &aux2, RotationalCorrelationAux &aux3, bool wrap, const MultidimArray< int > *mask)
Definition: filters.cpp:2150
int * n
Name of an image (std::string)
MetaDataVec SFout

◆ defineParams()

void ProgSortImages::defineParams ( )
inlinevirtual

Usage.

Reimplemented from XmippProgram.

Definition at line 107 of file mpi_image_sort.cpp.

108  {
109  addUsageLine("Sort a set of images by local similarity");
110  addUsageLine("+The program takes the first image of the input selfile as reference.");
111  addUsageLine("+Then, it looks for the image among the remaining ones that is most similar to it once they have been aligned.");
112  addUsageLine("+Now, the most similar image is aligned to the reference image (the first one) and is added to the list of sorted images.");
113  addUsageLine("+The second image acts now as the reference, and the most similar image among the remaining images is added to the list of sorted images.");
114  addUsageLine("+This process is repeated until no image is left.");
115  addUsageLine("+");
116  addUsageLine("+Note that only the MPI version of this program exists");
117  addParamsLine(" -i <selfile> : selfile of images");
118  addParamsLine(" --oroot <rootname> : output rootname");
119  addParamsLine(" : rootname.stk contains the list of aligned images.");
120  addParamsLine(" : rootname.xmd contains the correspondence between aligned images ");
121  addParamsLine(" : and the original images as well as the correlation coefficient ");
122  addParamsLine(" : between the aligned image and its predecessor in the list of ");
123  addParamsLine(" : aligned images.");
124  addParamsLine(" [--dont_center] : Do not center images as they are sorted");
125  addExampleLine("mpirun -np 8 `which xmipp_mpi_image_sort` -i images.xmd --oroot sorted");
126  }
void addExampleLine(const char *example, bool verbatim=true)
void addUsageLine(const char *line, bool verbatim=false)
void addParamsLine(const String &line)

◆ operator=() [1/2]

ProgSortImages& ProgSortImages::operator= ( const ProgSortImages )
delete

◆ operator=() [2/2]

ProgSortImages& ProgSortImages::operator= ( const ProgSortImages &&  )
delete

◆ produceSideInfo()

void ProgSortImages::produceSideInfo ( int  rank)
inline

Produce side info.

Definition at line 129 of file mpi_image_sort.cpp.

130  {
131  fnStack = fnRoot + ".stk";
132 
133  if (rank == 0)
135 
136  // Read input selfile and reference
137  MetaDataVec SF;
138  SF.read(fnSel);
139  FileName fnImg;
140  toClassify.reserve(SF.size());
141  FileName fnImageStack;
142  CorrelationAux aux;
144  bool thereIsClassCount=SF.containsLabel(MDL_CLASS_COUNT);
145  bool firstSelected=false;
146 
147  for (auto& row : SF)
148  {
149  bool proceed=true;
150  if (thereIsClassCount)
151  {
152  size_t classCount;
153  row.getValue(MDL_CLASS_COUNT,classCount);
154  if (classCount==0)
155  proceed=false;
156  }
157  if (proceed)
158  {
159  row.detach(); // so that we can copy it outside of the metadata
160  toClassify.emplace_back(dynamic_cast<MDRowVec&>(row));
161  row.getValue(MDL_IMAGE,fnImg);
162  if (!firstSelected)
163  {
164  row.getValue(MDL_IMAGE,fnImg);
165  lastImage.read(fnImg);
166  if (center)
167  centerImage(lastImage(),aux,aux2);
168  if (rank==0)
170  fnImageStack.compose(1,fnStack);
171  if (rank==0)
172  {
173  row.setValue(MDL_IMAGE_ORIGINAL,fnImg);
174  row.setValue(MDL_MAXCC,1.0);
175  SFout.addRow(dynamic_cast<MDRowVec&>(row));
176  }
177  firstSelected=true;
178  }
179  }
180  }
181  if (toClassify.size()==0)
182  return;
185  VEC_ELEM(stillToDo,0)=0;
186 
187  // Prepare mask
188  mask.resize(lastImage());
191  }
#define VEC_ELEM(v, i)
Definition: matrix1d.h:245
void resize(size_t Ndim, size_t Zdim, size_t Ydim, size_t Xdim, bool copy=true)
void read(const FileName &inFile, const std::vector< MDLabel > *desiredLabels=nullptr, bool decomposeStack=true) override
void BinaryCircularMask(MultidimArray< int > &mask, double radius, int mode, double x0, double y0, double z0)
Definition: mask.cpp:193
void write(const FileName &name="", size_t select_img=ALL_IMAGES, bool isStack=false, int mode=WRITE_OVERWRITE, CastWriteMode castMode=CW_CAST, int _swapWrite=0)
void compose(const String &str, const size_t no, const String &ext="")
size_t size() const override
size_t addRow(const MDRow &row) override
Matrix2D< double > centerImage(MultidimArray< double > &I, CorrelationAux &aux, RotationalCorrelationAux &aux2, int Niter, bool limitShift)
Definition: filters.cpp:3277
Image< double > lastImage
Name of an image from which MDL_IMAGE is coming from.
#define XSIZE(v)
Maximum cross-correlation for the image (double)
Matrix1D< int > stillToDo
MultidimArray< int > mask
void deleteFile() const
std::vector< MDRowVec > toClassify
Number of images assigned to the same class as this image.
void resizeNoCopy(int Xdim)
Definition: matrix1d.h:458
int read(const FileName &name, DataMode datamode=DATA, size_t select_img=ALL_IMAGES, bool mapData=false, int mode=WRITE_READONLY)
bool containsLabel(const MDLabel label) const override
void initConstant(T val)
Definition: matrix1d.cpp:83
Name of an image (std::string)
MetaDataVec SFout
constexpr int INNER_MASK
Definition: mask.h:47

◆ readParams()

void ProgSortImages::readParams ( )
inlinevirtual

Read argument from command line.

Reimplemented from XmippProgram.

Definition at line 88 of file mpi_image_sort.cpp.

89  {
90  fnSel = getParam("-i");
91  fnRoot = getParam("--oroot");
92  center = !checkParam("--dont_center");
93  }
const char * getParam(const char *param, int arg=0)
bool checkParam(const char *param)

◆ run()

void ProgSortImages::run ( )
inlinevirtual

Main routine.

Reimplemented from XmippProgram.

Definition at line 290 of file mpi_image_sort.cpp.

291  {
292  show();
294  while (stillToDo.sum()>0)
296  if (node->rank==0)
297  SFout.write(fnRoot+".xmd");
298  }
size_t size
Definition: xmipp_mpi.h:52
void write(const FileName &outFile, WriteModeMetaData mode=MD_OVERWRITE) const
double sum(bool average=false) const
Definition: matrix1d.cpp:652
Matrix1D< int > stillToDo
void produceSideInfo(int rank)
Produce side info.
void show()
Show.
void chooseNextImage(int rank, int Nproc)
Choose next image.
size_t rank
Definition: xmipp_mpi.h:52
MetaDataVec SFout

◆ show()

void ProgSortImages::show ( )
inline

Show.

Definition at line 96 of file mpi_image_sort.cpp.

97  {
98  if (!verbose)
99  return;
100  std::cerr << "Input selfile: " << fnSel << std::endl
101  << "Output rootname: " << fnRoot << std::endl
102  << "Center: " << center << std::endl
103  ;
104  }
int verbose
Verbosity level.

Member Data Documentation

◆ center

bool ProgSortImages::center

Center

Definition at line 50 of file mpi_image_sort.cpp.

◆ fnRoot

FileName ProgSortImages::fnRoot

Filename output rootname

Definition at line 42 of file mpi_image_sort.cpp.

◆ fnSel

FileName ProgSortImages::fnSel

Filename selection file containing the images

Definition at line 39 of file mpi_image_sort.cpp.

◆ fnStack

FileName ProgSortImages::fnStack

Filename of the output stack

Definition at line 47 of file mpi_image_sort.cpp.

◆ lastImage

Image<double> ProgSortImages::lastImage

Definition at line 62 of file mpi_image_sort.cpp.

◆ mask

MultidimArray<int> ProgSortImages::mask

Definition at line 65 of file mpi_image_sort.cpp.

◆ node

MpiNode* ProgSortImages::node

Definition at line 44 of file mpi_image_sort.cpp.

◆ SFout

MetaDataVec ProgSortImages::SFout

Definition at line 53 of file mpi_image_sort.cpp.

◆ stillToDo

Matrix1D<int> ProgSortImages::stillToDo

Definition at line 59 of file mpi_image_sort.cpp.

◆ toClassify

std::vector< MDRowVec > ProgSortImages::toClassify

Definition at line 56 of file mpi_image_sort.cpp.


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