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

#include <mpi_classify_CL2D.h>

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

Public Member Functions

 ProgClassifyCL2D (int argc, char **argv)
 MPI constructor. More...
 
 ~ProgClassifyCL2D ()
 Destructor. More...
 
 ProgClassifyCL2D (const ProgClassifyCL2D &)=delete
 
 ProgClassifyCL2D (const ProgClassifyCL2D &&)=delete
 
ProgClassifyCL2Doperator= (const ProgClassifyCL2D &)=delete
 
ProgClassifyCL2Doperator= (const ProgClassifyCL2D &&)=delete
 
void readParams ()
 Read. More...
 
void show () const
 Show. More...
 
void defineParams ()
 Usage. More...
 
void produceSideInfo ()
 Produce side info. More...
 
void run ()
 Run. 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 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
 Input selfile with the images to quantify. More...
 
FileName fnCodes0
 Input selfile with initial codes. More...
 
FileName fnOut
 Output rootname. More...
 
FileName fnODir
 Output directory. More...
 
int Niter
 Number of iterations. More...
 
int Ncodes0
 Initial number of code vectors. More...
 
int Ncodes
 Final number of code vectors. More...
 
int Nneighbours
 Number of neighbours. More...
 
double PminSize
 Minimum size of a node. More...
 
bool useCorrelation
 Use Correlation instead of Correntropy. More...
 
bool classicalMultiref
 Classical Multiref. More...
 
bool classifyAllImages
 Clasify all images. More...
 
bool classicalSplit
 Use ClassicalCriterion at split. More...
 
int NSplitTrials
 MaxTrials to split. More...
 
double maxShift
 Maximum shift. More...
 
bool normalizeImages
 Normalize input images. More...
 
bool mirrorImages
 Mirror. More...
 
bool useThresholdMask
 Use threshold mask. More...
 
double threshold
 Threshold to use. More...
 
bool alignImages
 Don't align images. More...
 
MetaDataDb SF
 
std::vector< size_t > objId
 
CL2D vq
 
MpiNodenode
 
double maxShift2
 
GaussianInterpolator gaussianInterpolator
 
size_t Ydim
 
size_t Xdim
 
MultidimArray< int > mask
 Mask for the background. More...
 
double sigma
 Noise in the images. 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 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

CL2D parameters.

Definition at line 224 of file mpi_classify_CL2D.h.

Constructor & Destructor Documentation

◆ ProgClassifyCL2D() [1/3]

ProgClassifyCL2D::ProgClassifyCL2D ( int  argc,
char **  argv 
)

MPI constructor.

Definition at line 1647 of file mpi_classify_CL2D.cpp.

1648 {
1649  node = new MpiNode(argc, argv);
1650  if (!node->isMaster())
1651  verbose = 0;
1652 }
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

◆ ~ProgClassifyCL2D()

ProgClassifyCL2D::~ProgClassifyCL2D ( )

Destructor.

Definition at line 1655 of file mpi_classify_CL2D.cpp.

1656 {
1657  delete node;
1658 }

◆ ProgClassifyCL2D() [2/3]

ProgClassifyCL2D::ProgClassifyCL2D ( const ProgClassifyCL2D )
delete

◆ ProgClassifyCL2D() [3/3]

ProgClassifyCL2D::ProgClassifyCL2D ( const ProgClassifyCL2D &&  )
delete

Member Function Documentation

◆ defineParams()

void ProgClassifyCL2D::defineParams ( )
virtual

Usage.

Reimplemented from XmippProgram.

Definition at line 1716 of file mpi_classify_CL2D.cpp.

1717 {
1718  addUsageLine("Divide a selfile into the desired number of classes. ");
1719  addUsageLine("+Vector quantization with correntropy and a probabilistic criterion is used for creating the subdivisions.");
1720  addUsageLine("+Correlation and the standard maximum correlation criterion can also be used and normally produce good results.");
1721  addUsageLine("+Correntropy and the probabilistic clustering criterion are recommended for images with very low SNR or cases in which the correlation have clear difficulties to converge.");
1722  addUsageLine("+");
1723  addUsageLine("+The algorithm is fully described in [[http://www.ncbi.nlm.nih.gov/pubmed/20362059][this article]].");
1724  addUsageLine("+");
1725  addUsageLine("+An interesting convergence criterion is the number of images changing classes between iterations. If a low percentage of the image change class, then the clustering is rather stable and clear.");
1726  addUsageLine("+If many images change class, it is likely that there is not enough SNR to determine so many classes. It is recommended to reduce the number of classes");
1727  addSeeAlsoLine("mpi_image_sort");
1728  addParamsLine(" -i <selfile> : Selfile with the input images");
1729  addParamsLine(" [--odir <dir=\".\">] : Output directory");
1730  addParamsLine(" [--oroot <root=class>] : Output rootname, by default, class");
1731  addParamsLine(" [--iter <N=20>] : Number of iterations");
1732  addParamsLine(" [--nref0 <N=2>] : Initial number of code vectors");
1733  addParamsLine("or --ref0 <selfile=\"\"> : Selfile with initial code vectors");
1734  addParamsLine(" [--nref <N=16>] : Final number of code vectors");
1735  addParamsLine(" [--neigh+ <N=4>] : Number of neighbour code vectors");
1736  addParamsLine(" : Set -1 for all");
1737  addParamsLine(" [--minsize+ <N=20>] : Percentage minimum node size");
1738  addParamsLine(" : Let's say that we have 1000 images and we are currently estimating 2 classes.");
1739  addParamsLine(" : On average each class should have about 500 images. ");
1740  addParamsLine(" : If one of the classes drops below 20% (by default) of 500,");
1741  addParamsLine(" : then that class is removed and the remaining class is splitted in two. ");
1742  addParamsLine(" : The same reasoning is applied when the number of classes is 4, but 20% ");
1743  addParamsLine(" : is calculated now on 250 (=1000/4).");
1744  addParamsLine(" [--distance <type=correntropy>] : Distance type");
1745  addParamsLine(" where <type>");
1746  addParamsLine(" correntropy correlation: See CL2D paper for the definition of correntropy");
1747  addParamsLine(" [--classicalMultiref] : Instead of enhanced clustering");
1748  addParamsLine(" [--classicalSplit] : Instead of enhanced clustering at the split iterations");
1749  addParamsLine(" [--maxSplitTrials <n=5>] : Maximum number of trials to split before giving up");
1750  addParamsLine(" [--maxShift <d=10>] : Maximum allowed shift");
1751  addParamsLine(" [--classifyAllImages] : By default, some images may not be classified. Use this option to classify them all.");
1752  addParamsLine(" [--dontNormalizeImages] : By default, input images are normalized to have 0 mean and standard deviation 1");
1753  addParamsLine(" [--dontMirrorImages] : By default, input images are studied unmirrored and mirrored");
1754  addParamsLine(" [--useThresholdMask <t>] : Use a mask to compare images. Remove pixels whose value is smaller or equal t");
1755  addParamsLine(" [--dontAlign] : Do not center the class representatives");
1756  addExampleLine("mpirun -np 3 `which xmipp_mpi_classify_CL2D` -i images.stk --nref 256 --oroot class --odir CL2Dresults --iter 10");
1757 }
void addSeeAlsoLine(const char *seeAlso)
void addExampleLine(const char *example, bool verbatim=true)
void addUsageLine(const char *line, bool verbatim=false)
void addParamsLine(const String &line)

◆ operator=() [1/2]

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

◆ operator=() [2/2]

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

◆ produceSideInfo()

void ProgClassifyCL2D::produceSideInfo ( )

Produce side info.

Definition at line 1759 of file mpi_classify_CL2D.cpp.

1760 {
1761  if (!useCorrelation)
1762  normalizeImages=false;
1763 
1765 
1766  gaussianInterpolator.initialize(6, 60000, false);
1767 
1768  // Get image dimensions
1769  SF.read(fnSel);
1770  SF.removeDisabled();
1771 
1772  size_t Zdim, Ndim;
1773  getImageSize(SF, Xdim, Ydim, Zdim, Ndim);
1774 
1775  // Prepare the Task distributor
1776  SF.findObjects(objId);
1777  // size_t Nimgs = objId.size();
1778 
1779  // Prepare mask for evaluating the noise outside
1780  mask.resize(prm->Ydim, prm->Xdim);
1781  mask.setXmippOrigin();
1783 
1784  // Read input code vectors if available
1785  std::vector<MultidimArray<double> > codes0;
1786  if (fnCodes0 != "")
1787  {
1788  Image<double> I;
1789  MetaDataVec SFCodes(fnCodes0);
1790 
1791  size_t Xdim0, Ydim0, Zdim0, Ndim0;
1792  getImageSize(SFCodes, Xdim0, Ydim0, Zdim0, Ndim0);
1793  if (Xdim0!=Xdim || Ydim0!=Ydim)
1794  REPORT_ERROR(ERR_MULTIDIM_SIZE,"Input reference and images are not of the same size");
1795 
1796  for (size_t objId : SFCodes.ids())
1797  {
1798  I.readApplyGeo(SFCodes, objId);
1799  I().setXmippOrigin();
1800  codes0.push_back(I());
1801  }
1802  Ncodes0 = codes0.size();
1803  }
1804  vq.initialize(SF, codes0);
1805 }
FileName fnSel
Input selfile with the images to quantify.
void resize(size_t Ndim, size_t Zdim, size_t Ydim, size_t Xdim, bool copy=true)
#define REPORT_ERROR(nerr, ErrormMsg)
Definition: xmipp_error.h:211
void findObjects(std::vector< size_t > &objectsOut, const MDQuery &query) const override
FileName fnCodes0
Input selfile with initial codes.
void BinaryCircularMask(MultidimArray< int > &mask, double radius, int mode, double x0, double y0, double z0)
Definition: mask.cpp:193
int Ncodes0
Initial number of code vectors.
void getImageSize(const MetaData &md, size_t &Xdim, size_t &Ydim, size_t &Zdim, size_t &Ndim, MDLabel image_label)
Incorrect MultidimArray size.
Definition: xmipp_error.h:174
int readApplyGeo(const FileName &name, const MDRow &row, const ApplyGeoParams &params=DefaultApplyGeoParams)
bool normalizeImages
Normalize input images.
std::vector< size_t > objId
double maxShift
Maximum shift.
bool useCorrelation
Use Correlation instead of Correntropy.
void initialize(MetaDataDb &_SF, std::vector< MultidimArray< double > > &_codes0)
Initialize.
virtual void removeDisabled()
void initialize(double xmax, int N, bool normalize=true)
void read(const FileName &inFile, const std::vector< MDLabel > *desiredLabels=NULL, bool decomposeStack=true) override
ProgClassifyCL2D * prm
MultidimArray< int > mask
Mask for the background.
GaussianInterpolator gaussianInterpolator
constexpr int INNER_MASK
Definition: mask.h:47

◆ readParams()

void ProgClassifyCL2D::readParams ( )
virtual

Read.

Reimplemented from XmippProgram.

Definition at line 1661 of file mpi_classify_CL2D.cpp.

1662 {
1663  fnSel = getParam("-i");
1664  fnOut = getParam("--oroot");
1665  fnODir = getParam("--odir");
1666  fnCodes0 = getParam("--ref0");
1667  Niter = getIntParam("--iter");
1668  Nneighbours = getIntParam("--neigh");
1669  Ncodes0 = getIntParam("--nref0");
1670  Ncodes = getIntParam("--nref");
1671  PminSize = getDoubleParam("--minsize");
1672  String aux;
1673  aux = getParam("--distance");
1674  useCorrelation = aux == "correlation";
1675  classicalMultiref = checkParam("--classicalMultiref");
1676  classicalSplit = checkParam("--classicalSplit");
1677  NSplitTrials = getIntParam("--maxSplitTrials");
1678  maxShift = getDoubleParam("--maxShift");
1679  classifyAllImages = checkParam("--classifyAllImages");
1680  normalizeImages = !checkParam("--dontNormalizeImages");
1681  mirrorImages = !checkParam("--dontMirrorImages");
1682  useThresholdMask = checkParam("--useThresholdMask");
1683  if (useThresholdMask)
1684  threshold=getDoubleParam("--useThresholdMask");
1685  alignImages = !checkParam("--dontAlign");
1686 
1687  prm = this; // FIXME HACK because of the global variable. Solve it properly
1688 }
FileName fnSel
Input selfile with the images to quantify.
int Nneighbours
Number of neighbours.
FileName fnODir
Output directory.
double getDoubleParam(const char *param, int arg=0)
int Niter
Number of iterations.
int NSplitTrials
MaxTrials to split.
FileName fnCodes0
Input selfile with initial codes.
int Ncodes0
Initial number of code vectors.
double PminSize
Minimum size of a node.
bool classicalMultiref
Classical Multiref.
bool classifyAllImages
Clasify all images.
int Ncodes
Final number of code vectors.
bool normalizeImages
Normalize input images.
const char * getParam(const char *param, int arg=0)
bool alignImages
Don&#39;t align images.
double maxShift
Maximum shift.
bool useCorrelation
Use Correlation instead of Correntropy.
FileName fnOut
Output rootname.
bool mirrorImages
Mirror.
std::string String
Definition: xmipp_strings.h:34
bool useThresholdMask
Use threshold mask.
bool classicalSplit
Use ClassicalCriterion at split.
bool checkParam(const char *param)
double threshold
Threshold to use.
ProgClassifyCL2D * prm
int getIntParam(const char *param, int arg=0)

◆ run()

void ProgClassifyCL2D::run ( )
virtual

Run.

Reimplemented from XmippProgram.

Definition at line 1807 of file mpi_classify_CL2D.cpp.

1808 {
1809  CREATE_LOG();
1810  show();
1811  produceSideInfo();
1812 
1813  // Run all iterations
1814  int level = 0;
1815  int Q = vq.P.size();
1816  bool originalClassicalMultiref=classicalMultiref;
1817  if (Q==1)
1818  classicalMultiref=true;
1819  vq.run(fnODir, fnOut, level);
1820 
1821  while (Q < Ncodes)
1822  {
1823  if (node->rank == 1)
1824  std::cout << "Spliting nodes ...\n";
1825 
1826  int Nclean = vq.cleanEmptyNodes();
1827  int Nsplits = XMIPP_MIN(Q,Ncodes-Q) + Nclean;
1828 
1829  for (int i = 0; i < Nsplits; i++)
1830  {
1831  vq.splitFirstNode();
1832  if (node->rank == 1)
1833  std::cout << "Currently there are " << vq.P.size() << " nodes"
1834  << std::endl;
1835  }
1836 
1837  Q = vq.P.size();
1838  level++;
1839  classicalMultiref=originalClassicalMultiref || Q==1;
1840  vq.run(fnODir, fnOut, level);
1841  }
1842  if (node->rank == 1)
1843  {
1844  std::sort(vq.P.begin(), vq.P.end(), SDescendingClusterSort());
1845  Q = vq.P.size();
1846  MetaDataDb SFq, SFclassified, SFaux, SFaux2;
1847  for (int q = 0; q < Q; q++)
1848  {
1849  SFq.read(formatString("class%06d_images@%s/level_%02d/%s_classes.xmd",q+1,fnODir.c_str(),level,fnOut.c_str()));
1850  SFq.fillConstant(MDL_REF, integerToString(q + 1));
1851  SFq.fillConstant(MDL_ENABLED, "1");
1852  SFclassified.unionAll(SFq);
1853  }
1854  SFaux = SF;
1855  SFaux.subtraction(SFclassified, MDL_IMAGE);
1856  SFaux.fillConstant(MDL_ENABLED, "-1");
1857  SFaux2.join1(SFclassified, SF, MDL_IMAGE, LEFT);
1858  SFclassified.clear();
1859  SFaux2.unionAll(SFaux);
1860  SFaux.clear();
1861  SFaux.sort(SFaux2, MDL_IMAGE);
1862  SFaux.write(fnODir+"/images.xmd");
1863  }
1864  CLOSE_LOG();
1865 }
void subtraction(const MetaDataDb &mdIn, const MDLabel label)
FileName fnODir
Output directory.
String integerToString(int I, int _width, char fill_with)
void show() const
Show.
bool classicalMultiref
Classical Multiref.
void run(const FileName &fnODir, const FileName &fnOut, int level)
int Ncodes
Final number of code vectors.
#define i
Is this image enabled? (int [-1 or 1])
void fillConstant(MDLabel label, const String &value) override
void clear() override
Definition: metadata_db.cpp:54
void write(const FileName &outFile, WriteModeMetaData mode=MD_OVERWRITE) const override
void sort(struct DCEL_T *dcel)
Definition: sorting.cpp:18
#define XMIPP_MIN(x, y)
Definition: xmipp_macros.h:181
size_t rank
Definition: xmipp_mpi.h:52
#define CREATE_LOG()
Class to which the image belongs (int)
void sort(MetaDataDb &MDin, const MDLabel sortLabel, bool asc=true, int limit=-1, int offset=0)
FileName fnOut
Output rootname.
String formatString(const char *format,...)
void read(const FileName &inFile, const std::vector< MDLabel > *desiredLabels=NULL, bool decomposeStack=true) override
void join1(const MetaDataDb &mdInLeft, const MetaDataDb &mdInRight, const MDLabel label, JoinType type=LEFT)
void unionAll(const MetaDataDb &mdIn)
void produceSideInfo()
Produce side info.
std::vector< CL2DClass * > P
List of nodes.
#define CLOSE_LOG()
void splitFirstNode()
int cleanEmptyNodes()
Name of an image (std::string)

◆ show()

void ProgClassifyCL2D::show ( ) const
virtual

Show.

Reimplemented from XmippProgram.

Definition at line 1690 of file mpi_classify_CL2D.cpp.

1690  {
1691  if (!verbose)
1692  return;
1693  std::cout << "Input images: " << fnSel << std::endl
1694  << "Output root: " << fnOut << std::endl
1695  << "Output dir: " << fnODir << std::endl
1696  << "Iterations: " << Niter << std::endl
1697  << "CodesSel0: " << fnCodes0 << std::endl
1698  << "Codes0: " << Ncodes0 << std::endl
1699  << "Codes: " << Ncodes << std::endl
1700  << "Neighbours: " << Nneighbours << std::endl
1701  << "Minimum node size: " << PminSize << std::endl
1702  << "Use Correlation: " << useCorrelation << std::endl
1703  << "Classical Multiref: " << classicalMultiref << std::endl
1704  << "Classical Split: " << classicalSplit << std::endl
1705  << "Max. Split trials: " << NSplitTrials << std::endl
1706  << "Maximum shift: " << maxShift << std::endl
1707  << "Classify all images: " << classifyAllImages << std::endl
1708  << "Normalize images: " << normalizeImages << std::endl
1709  << "Mirror images: " << mirrorImages << std::endl
1710  << "Align images: " << alignImages << std::endl
1711  ;
1712  if (useThresholdMask)
1713  std::cout << "Threshold mask: " << threshold << std::endl;
1714 }
FileName fnSel
Input selfile with the images to quantify.
int Nneighbours
Number of neighbours.
FileName fnODir
Output directory.
int Niter
Number of iterations.
int NSplitTrials
MaxTrials to split.
FileName fnCodes0
Input selfile with initial codes.
int Ncodes0
Initial number of code vectors.
double PminSize
Minimum size of a node.
bool classicalMultiref
Classical Multiref.
bool classifyAllImages
Clasify all images.
int Ncodes
Final number of code vectors.
bool normalizeImages
Normalize input images.
bool alignImages
Don&#39;t align images.
double maxShift
Maximum shift.
bool useCorrelation
Use Correlation instead of Correntropy.
int verbose
Verbosity level.
FileName fnOut
Output rootname.
bool mirrorImages
Mirror.
bool useThresholdMask
Use threshold mask.
bool classicalSplit
Use ClassicalCriterion at split.
double threshold
Threshold to use.

Member Data Documentation

◆ alignImages

bool ProgClassifyCL2D::alignImages

Don't align images.

Definition at line 284 of file mpi_classify_CL2D.h.

◆ classicalMultiref

bool ProgClassifyCL2D::classicalMultiref

Classical Multiref.

Definition at line 257 of file mpi_classify_CL2D.h.

◆ classicalSplit

bool ProgClassifyCL2D::classicalSplit

Use ClassicalCriterion at split.

Definition at line 263 of file mpi_classify_CL2D.h.

◆ classifyAllImages

bool ProgClassifyCL2D::classifyAllImages

Clasify all images.

Definition at line 260 of file mpi_classify_CL2D.h.

◆ fnCodes0

FileName ProgClassifyCL2D::fnCodes0

Input selfile with initial codes.

Definition at line 230 of file mpi_classify_CL2D.h.

◆ fnODir

FileName ProgClassifyCL2D::fnODir

Output directory.

Definition at line 236 of file mpi_classify_CL2D.h.

◆ fnOut

FileName ProgClassifyCL2D::fnOut

Output rootname.

Definition at line 233 of file mpi_classify_CL2D.h.

◆ fnSel

FileName ProgClassifyCL2D::fnSel

Input selfile with the images to quantify.

Definition at line 227 of file mpi_classify_CL2D.h.

◆ gaussianInterpolator

GaussianInterpolator ProgClassifyCL2D::gaussianInterpolator

Definition at line 329 of file mpi_classify_CL2D.h.

◆ mask

MultidimArray<int> ProgClassifyCL2D::mask

Mask for the background.

Definition at line 335 of file mpi_classify_CL2D.h.

◆ maxShift

double ProgClassifyCL2D::maxShift

Maximum shift.

Definition at line 269 of file mpi_classify_CL2D.h.

◆ maxShift2

double ProgClassifyCL2D::maxShift2

Definition at line 326 of file mpi_classify_CL2D.h.

◆ mirrorImages

bool ProgClassifyCL2D::mirrorImages

Mirror.

Definition at line 275 of file mpi_classify_CL2D.h.

◆ Ncodes

int ProgClassifyCL2D::Ncodes

Final number of code vectors.

Definition at line 245 of file mpi_classify_CL2D.h.

◆ Ncodes0

int ProgClassifyCL2D::Ncodes0

Initial number of code vectors.

Definition at line 242 of file mpi_classify_CL2D.h.

◆ Niter

int ProgClassifyCL2D::Niter

Number of iterations.

Definition at line 239 of file mpi_classify_CL2D.h.

◆ Nneighbours

int ProgClassifyCL2D::Nneighbours

Number of neighbours.

Definition at line 248 of file mpi_classify_CL2D.h.

◆ node

MpiNode* ProgClassifyCL2D::node

Definition at line 323 of file mpi_classify_CL2D.h.

◆ normalizeImages

bool ProgClassifyCL2D::normalizeImages

Normalize input images.

Definition at line 272 of file mpi_classify_CL2D.h.

◆ NSplitTrials

int ProgClassifyCL2D::NSplitTrials

MaxTrials to split.

Definition at line 266 of file mpi_classify_CL2D.h.

◆ objId

std::vector<size_t> ProgClassifyCL2D::objId

Definition at line 317 of file mpi_classify_CL2D.h.

◆ PminSize

double ProgClassifyCL2D::PminSize

Minimum size of a node.

Definition at line 251 of file mpi_classify_CL2D.h.

◆ SF

MetaDataDb ProgClassifyCL2D::SF

Definition at line 314 of file mpi_classify_CL2D.h.

◆ sigma

double ProgClassifyCL2D::sigma

Noise in the images.

Definition at line 338 of file mpi_classify_CL2D.h.

◆ threshold

double ProgClassifyCL2D::threshold

Threshold to use.

Definition at line 281 of file mpi_classify_CL2D.h.

◆ useCorrelation

bool ProgClassifyCL2D::useCorrelation

Use Correlation instead of Correntropy.

Definition at line 254 of file mpi_classify_CL2D.h.

◆ useThresholdMask

bool ProgClassifyCL2D::useThresholdMask

Use threshold mask.

Definition at line 278 of file mpi_classify_CL2D.h.

◆ vq

CL2D ProgClassifyCL2D::vq

Definition at line 320 of file mpi_classify_CL2D.h.

◆ Xdim

size_t ProgClassifyCL2D::Xdim

Definition at line 332 of file mpi_classify_CL2D.h.

◆ Ydim

size_t ProgClassifyCL2D::Ydim

Definition at line 332 of file mpi_classify_CL2D.h.


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