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

#include <histogram.h>

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

Public Member Functions

 Histogram1D ()
 
 Histogram1D (const Histogram1D &H)
 
void clear ()
 
Histogram1Doperator= (const Histogram1D &H)
 
void assign (const Histogram1D &H)
 
void init (double min_val, double max_val, int n_steps)
 
void insert_value (double val)
 
double percentil (double percent_mass)
 
double mass_below (double value)
 
double mass_above (double value)
 
void write (const FileName &fn, MDLabel=MDL_X, MDLabel=MDL_COUNT)
 
void val2index (double v, int &i) const
 
void index2val (double i, double &v) const
 
double hist_min () const
 
double hist_max () const
 
double step () const
 
int stepNo () const
 
double sampleNo () const
 
double entropy () const
 
- Public Member Functions inherited from MultidimArray< double >
 MultidimArray ()
 
 MultidimArray (size_t Ndim, size_t Zdim, size_t Ydim, size_t Xdim, double *data)
 
 MultidimArray (size_t Ndim, int Zdim, int Ydim, int Xdim)
 
 MultidimArray (int Zdim, int Ydim, int Xdim)
 
 MultidimArray (int Ydim, int Xdim)
 
 MultidimArray (int Xdim)
 
 MultidimArray (const MultidimArray< double > &V)
 
 MultidimArray (MultidimArray< double > &&V) noexcept
 
 MultidimArray (const Matrix1D< double > &V)
 
 MultidimArray (const std::vector< double > &vector)
 
virtual ~MultidimArray ()
 
void clear ()
 
void coreInit () noexcept
 
void swap (MultidimArray< double > &other) noexcept
 
void coreAllocate (size_t _ndim, int _zdim, int _ydim, int _xdim)
 
void coreAllocate ()
 
void coreAllocateReuse ()
 
FILE * mmapFile (double *&_data, size_t nzyxDim) const
 
void coreDeallocate () noexcept
 
void alias (const MultidimArray< double > &m)
 
void aliasRow (const MultidimArray< double > &m, size_t select_row)
 
void aliasSlice (const MultidimArray< double > &m, size_t select_slice)
 
void aliasImageInStack (const MultidimArray< double > &m, size_t select_image)
 
void resize (size_t Ndim, size_t Zdim, size_t Ydim, size_t Xdim, bool copy=true)
 
void resize (const MultidimArray< T1 > &v, bool copy=true)
 
void resizeNoCopy (const MultidimArray< T1 > &v)
 
void checkDimensionWithDebug (int dim, const char *file, int line) const
 
void selfWindow (int n0, int z0, int y0, int x0, int nF, int zF, int yF, int xF, double init_value=0)
 
void selfWindow (int z0, int y0, int x0, int zF, int yF, int xF, double init_value=0)
 
void selfWindow (int y0, int x0, int yF, int xF, double init_value=0)
 
void selfWindow (int x0, int xF, double init_value=0)
 
void window (MultidimArray< T1 > &result, int n0, int z0, int y0, int x0, int nF, int zF, int yF, int xF, T1 init_value=0) const
 
void window (MultidimArray< T1 > &result, int z0, int y0, int x0, int zF, int yF, int xF, T1 init_value=0) const
 
void window (MultidimArray< T1 > &result, int y0, int x0, int yF, int xF, T1 init_value=0) const
 
void window (MultidimArray< T1 > &result, int x0, int xF, T1 init_value=0) const
 
void patch (MultidimArray< double > patchArray, int x, int y)
 
double & operator() (const Matrix1D< double > &v) const
 
double & operator() (const Matrix1D< int > &v) const
 
double & operator() (size_t n, int k, int i, int j) const
 
double & operator() (int k, int i, int j) const
 
double & operator() (int i, int j) const
 
double & operator() (int i) const
 
double & operator[] (size_t i) const
 
void * getArrayPointer () const
 
void getImage (size_t n, MultidimArray< double > &M, size_t n2=0) const
 
void getSlice (int k, MultidimArray< T1 > &M, char axis='Z', bool reverse=false, size_t n=0) const
 
void getSliceAsMatrix (size_t k, Matrix2D< double > &m) const
 
void getAliasAsRowVector (Matrix1D< double > &m) const
 
void setSlice (int k, const MultidimArray< T1 > &v, size_t n=0)
 
void reslice (MultidimArray< T1 > &out, AxisView face, bool flip=false, size_t n=0) const
 
void reslice (AxisView face, bool flip=false, size_t n=0)
 
void getCol (size_t j, MultidimArray< double > &v) const
 
void setCol (size_t j, const MultidimArray< double > &v)
 
void getRow (size_t i, MultidimArray< double > &v) const
 
void setRow (int i, const MultidimArray< double > &v)
 
void getReal (MultidimArray< double > &realImg) const
 
void getImag (MultidimArray< double > &imagImg) const
 
void toPhysical (int k_log, int i_log, int j_log, int &k_phys, int &i_phys, int &j_phys) const
 
void toPhysical (int i_log, int j_log, int &i_phys, int &j_phys) const
 
void toPhysical (int i_log, int &i_phys) const
 
void toLogical (int k_phys, int i_phys, int j_phys, int &k_log, int &i_log, int &j_log) const
 
void toLogical (int i_phys, int j_phys, int &i_log, int &j_log) const
 
void toLogical (int i_phys, int &i_log) const
 
double interpolatedElement3D (double x, double y, double z, double outside_value=(double) 0) const
 
double interpolatedElement2D (double x, double y, double outside_value=(double) 0) const
 
double interpolatedElement2DOutsideZero (double x, double y) const
 
double interpolatedElement1D (double x, double outside_value=(double) 0) const
 
double interpolatedElementBSpline3D (double x, double y, double z, int SplineDegree=3) const
 
double interpolatedElementBSpline2D (double x, double y, int SplineDegree=3) const
 
double interpolatedElementBSpline2D_Degree3 (double x, double y) const
 
double interpolatedElementBSpline1D (double x, int SplineDegree=3) const
 
void printStats (std::ostream &out=std::cout) const
 
double computeMax () const
 
void maxIndex (int &jmax) const
 
void maxIndex (size_t &lmax, int &kmax, int &imax, int &jmax) const
 
double computeMin () const
 
void minIndex (int &lmin, int &kmin, int &imin, int &jmin) const
 
void minIndex (int &kmin, int &imin, int &jmin) const
 
void minIndex (int &imin, int &jmin) const
 
void minIndex (int &jmin) const
 
void computeDoubleMinMax (double &minval, double &maxval) const
 
void computeDoubleMinMaxRange (double &minval, double &maxval, size_t offset, size_t size) const
 
double computeAvg () const
 
double computeStddev () const
 
void computeStats (double &avg, double &stddev, double &minval, double &maxval) const
 
void computeStats (double &avg, double &stddev, double &min_val, double &max_val, Matrix1D< int > &corner1, Matrix1D< int > &corner2, size_t n=0)
 
void computeAvgStdev (U &avg, U &stddev) const
 
void computeAvgStdev_within_binary_mask (const MultidimArray< int > &mask, double &avg, double &stddev) const
 
void computeMedian_within_binary_mask (const MultidimArray< int > &mask, double &median) const
 
double computeMedian () const
 
void rangeAdjust (double minF, double maxF)
 
void rangeAdjust (double minF, double maxF, MultidimArray< int > &mask)
 
void rangeAdjust (const MultidimArray< double > &example, const MultidimArray< int > *mask=NULL)
 
void statisticsAdjust (U avgF, U stddevF)
 
void initConstant (double val)
 
void initZeros (const MultidimArray< T1 > &op)
 
void initZeros ()
 
void initZeros (size_t Ndim, size_t Zdim, size_t Ydim, size_t Xdim)
 
void initZeros (int Xdim)
 
void initZeros (int Ydim, int Xdim)
 
void initZeros (int Zdim, int Ydim, int Xdim)
 
void initLinear (double minF, double maxF, int n=1, const String &mode="incr")
 
void initRandom (double op1, double op2, RandomMode mode=RND_UNIFORM)
 
void addNoise (double op1, double op2, const String &mode="uniform", double df=3.) const
 
MultidimArray< double > operator+ (const MultidimArray< double > &op1) const
 
MultidimArray< double > operator- (const MultidimArray< double > &op1) const
 
MultidimArray< double > operator* (const MultidimArray< double > &op1) const
 
MultidimArray< double > operator/ (const MultidimArray< double > &op1) const
 
void operator+= (const MultidimArray< double > &op1)
 
void operator-= (const MultidimArray< double > &op1)
 
void operator*= (const MultidimArray< double > &op1)
 
void operator/= (const MultidimArray< double > &op1)
 
double dotProduct (const MultidimArray< double > &op1)
 
MultidimArray< double > operator+ (double op1) const
 
MultidimArray< double > operator- (double op1) const
 
MultidimArray< double > operator* (double op1) const
 
MultidimArray< double > operator/ (double op1) const
 
void operator+= (const double &op1)
 
void operator-= (const double &op1)
 
void operator*= (const double &op1)
 
void operator/= (const double &op1)
 
void equal (double op1, MultidimArray< char > &result) const
 
MultidimArray< double > operator- () const
 
bool equal (const MultidimArray< double > &op, double accuracy=XMIPP_EQUAL_ACCURACY) const
 
MultidimArray< double > & operator= (const MultidimArray< double > &op1)
 
MultidimArray< double > & operator= (MultidimArray< double > &&other) noexcept
 
MultidimArray< double > & operator= (const Matrix2D< double > &op1)
 
void copy (Matrix2D< double > &op1) const
 
double *** adaptForNumericalRecipes3D (size_t n=0) const
 
void killAdaptationForNumericalRecipes3D (double ***m) const
 
double ** adaptForNumericalRecipes2D (size_t n=0) const
 
double * adaptForNumericalRecipes22D () const
 
void loadFromNumericalRecipes2D (double **m, int Ydim, int Xdim)
 
void killAdaptationForNumericalRecipes2D (double **m) const
 
void killAdaptationForNumericalRecipes22D (double **m) const
 
double * adaptForNumericalRecipes1D () const
 
void killAdaptationForNumericalRecipes1D (double *m) const
 
void centerOfMass (Matrix1D< double > &center, void *mask=NULL, size_t n=0)
 
void sort (MultidimArray< double > &result) const
 
void indexSort (MultidimArray< int > &indx) const
 
void cumlativeDensityFunction (MultidimArray< double > &cdf)
 
void threshold (const String &type, double a, double b=0, MultidimArray< int > *mask=NULL)
 
size_t countThreshold (const String &type, double a, double b, MultidimArray< int > *mask=NULL)
 
void substitute (double oldv, double newv, double accuracy=XMIPP_EQUAL_ACCURACY, MultidimArray< int > *mask=NULL)
 
void randomSubstitute (double oldv, double avgv, double sigv, double accuracy=XMIPP_EQUAL_ACCURACY, MultidimArray< int > *mask=NULL)
 
void binarize (double val=0, double accuracy=XMIPP_EQUAL_ACCURACY, MultidimArray< int > *mask=NULL)
 
void binarizeRange (double valMin=0, double valMax=255, MultidimArray< int > *mask=NULL)
 
void selfROUND ()
 
void selfCEIL ()
 
void selfFLOOR ()
 
void selfABS ()
 
void selfNormalizeInterval (double minPerc=0.25, double maxPerc=0.75, int Npix=1000)
 
void selfSQRT ()
 
double sum () const
 
double sum2 () const
 
void selfLog10 ()
 
void selfLog ()
 
void selfReverseX ()
 
void selfReverseY ()
 
void selfReverseZ ()
 
void profile (int x0, int y0, int xF, int yF, int N, MultidimArray< double > &profile) const
 
void showWithGnuPlot (const String &xlabel, const String &title)
 
void edit ()
 
void write (const FileName &fn) const
 
- Public Member Functions inherited from MultidimArrayBase
virtual ~MultidimArrayBase ()
 
int getDim () const
 
void setMmap (bool mmap)
 
void maxIndex (ArrayCoord &pos) const
 
void maxIndex (int &kmax, int &imax, int &jmax) const
 
void maxIndex (int &imax, int &jmax) const
 
void maxIndex (int &jmax) const
 
void printShape (std::ostream &out=std::cout) const
 
void setNdim (int Ndim)
 
void setZdim (int Zdim)
 
void setYdim (int Ydim)
 
void setXdim (int Xdim)
 
void setDimensions (int Xdim, int Ydim, int Zdim, size_t Ndim)
 
void setDimensions (ArrayDim &newDim)
 
void getDimensions (size_t &Xdim, size_t &Ydim, size_t &Zdim, size_t &Ndim) const
 
void getDimensions (ArrayDim &idim) const
 
ArrayDim getDimensions () const
 
void getDimensions (int *size) const
 
size_t getSize () const
 
void resize (size_t Zdim, size_t Ydim, size_t Xdim)
 
void resize (size_t Ydim, size_t Xdim)
 
void resize (size_t Xdim)
 
void resize (ArrayDim &adim, bool copy=true)
 
void resizeNoCopy (size_t Ndim, size_t Zdim, size_t Ydim, size_t Xdim)
 
void resizeNoCopy (size_t Zdim, size_t Ydim, size_t Xdim)
 
void resizeNoCopy (size_t Ydim, size_t Xdim)
 
void resizeNoCopy (size_t Xdim)
 
size_t rowNumber () const
 
size_t colNumber () const
 
void copyShape (const MultidimArrayBase &m)
 
bool sameShape (const MultidimArrayBase &op) const
 
void setXmippOrigin ()
 
void resetOrigin ()
 
void moveOriginTo (int k, int i, int j)
 
void moveOriginTo (int i, int j)
 
int startingZ () const
 
int finishingZ () const
 
int startingY () const
 
int finishingY () const
 
int startingX () const
 
int finishingX () const
 
bool isCorner (const Matrix1D< double > &v) const
 
bool outside (int k, int i, int j) const
 
bool outside (int i, int j) const
 
bool outside (int i) const
 
bool outside (const Matrix1D< double > &r) const
 

Public Attributes

double hmin
 
double hmax
 
double step_size
 
double istep_size
 
int no_samples
 
- Public Attributes inherited from MultidimArray< double >
double * data
 
- Public Attributes inherited from MultidimArrayBase
bool destroyData
 
size_t ndim
 
size_t zdim
 
size_t ydim
 
size_t xdim
 
size_t yxdim
 
size_t zyxdim
 
size_t nzyxdim
 
int zinit
 
int yinit
 
int xinit
 
bool mmapOn
 
FILE * mFd
 
size_t nzyxdimAlloc
 

Friends

std::ostream & operator<< (std::ostream &o, const Histogram1D &hist)
 

Detailed Description

Histograms with 1 parameter

This class of histograms are the usual ones where we want to count the number of elements within a certain range of a variable, then we make the histogram of that variable. The range is divided into small subranges within which the values will be grouped. Any value outside the global range will not be counted in the histogram.

To see exactly which is the division between subranges let's have a look on the following example where an histogram between 0 and 2 is computed with 5 steps.

[......) [.......]
[ ) [ ]
[ )[......) [ ]
[ )[ ) [......)[ ]
[ )[ )[......)[ )[ ]
[ )[ )[ )[ )[ ]
[ 0 )[ 1 )[ 2 )[ 3 )[ 4 ]
[ )[ )[ )[ )[ ]
|---------|---------|---------|---------|
0.0 0.5 1.0 1.5 2.0

The border points are 0.0, 0.4, 0.8, 1.2, 1.6 and 2.0. The brackets and parenthesis try to represent where the border point belongs to, and the numbers within the bars are the index of each bar whithin the histogram. The height of each bar is the number of times that a value within that subrange has been inserted. Be careful that this is not a probability density function (pdf), to be so it should be divided by the total number of values inserted.

The simplest way of computing a histograms is the following:

// Variable definition
// Matrix initialisation
A.init_random(0, 100);
// Histogram calculation with 200 bins
compute_hist(A, hist, 200);
// Effective range computation
double eff0 = hist.percentil(2.5);
double effF = hist.percentil(97.5);

The following example shows how to work with the histograms. In it we will compute which is the central range within which the 95% of the values of a matrix are comprised. This example could be even simplified by using the function compute_hist but it has been kept like this to show the idea behind the histograms

// Variable definition
double min_val, max_val;
double eff0, effF;
// Matrix initialisation
A.init_random(0, 100);
// Histogram initialisation
min_val = A.min();
max_val = A.max();
hist.init(min_val, max_val, 200);
// Histogram computation
for (int i=STARTINGY(A); i<=FINISHINGY(A); i++)
for (int j=STARTINGX(A); j<=FINISHINGX(A); j++)
// Effective range computation
eff0 = hist.percentil(2.5);
effF = hist.percentil(97.5);
std::cout << "The effective range goes from " << eff0
<< " to " << effF << std::endl;

Definition at line 121 of file histogram.h.

Constructor & Destructor Documentation

◆ Histogram1D() [1/2]

Histogram1D::Histogram1D ( )
inline

Empty constructor

Creates an empty histogram. Before using it you must initialise it with init.

Definition at line 140 of file histogram.h.

141  {
142  clear();
143  }
void clear()
Definition: histogram.cpp:40

◆ Histogram1D() [2/2]

Histogram1D::Histogram1D ( const Histogram1D H)
inline

Copy constructor

Makes an exact copy of the given histogram into another histogram.

Histogram1D hist2(hist1);

Definition at line 153 of file histogram.h.

154  {
155  clear();
156  *this = H;
157  }
void clear()
Definition: histogram.cpp:40

Member Function Documentation

◆ assign()

void Histogram1D::assign ( const Histogram1D H)

Another function for assignament.

Definition at line 66 of file histogram.cpp.

67 {
68  *this = H;
69 }

◆ clear()

void Histogram1D::clear ( )
virtual

Empties an histogram

Remember to initialise it before using it again.

hist.clear();

Implements MultidimArrayBase.

Definition at line 40 of file histogram.cpp.

41 {
42  hmin = 0;
43  hmax = 0;
44  step_size = 0;
45  istep_size = 0;
46  no_samples = 0;
48 }
double hmin
Definition: histogram.h:125
double hmax
Definition: histogram.h:126
int no_samples
Definition: histogram.h:129
double step_size
Definition: histogram.h:127
double istep_size
Definition: histogram.h:128

◆ entropy()

double Histogram1D::entropy ( ) const

Measure the entropy of this histogram.

Before computing the entropy, the histogram is corrected with a Laplace correction. The entropy is computed as sum(-p*log(p))

Definition at line 238 of file histogram.cpp.

239 {
241  p.initZeros(XSIZE(*this));
242  double pSum = 0;
244  {
245  A1D_ELEM(p,i) = A1D_ELEM(*this,i) + 1;
246  pSum += A1D_ELEM(p,i);
247  }
248  double entropy = 0;
249  double ipSum = 1.0 / pSum;
251  {
252  double pi = A1D_ELEM(p,i) * ipSum;
253  entropy -= pi * log(pi);
254  }
255  return entropy;
256 }
#define A1D_ELEM(v, i)
#define i
void log(Image< double > &op)
#define XSIZE(v)
#define FOR_ALL_ELEMENTS_IN_ARRAY1D(v)
double entropy() const
Definition: histogram.cpp:238
#define pi
void initZeros(const MultidimArray< T1 > &op)

◆ hist_max()

double Histogram1D::hist_max ( ) const
inline

Maximum value where the histogram is defined

std::cout << "Maximum value for histogram " << hist.max() << std::endl;

Definition at line 317 of file histogram.h.

318  {
319  return hmax;
320  }
double hmax
Definition: histogram.h:126

◆ hist_min()

double Histogram1D::hist_min ( ) const
inline

Minimum value where the histogram is defined

std::cout << "Minimum value for histogram " << hist.min() << std::endl;

Definition at line 306 of file histogram.h.

307  {
308  return hmin;
309  }
double hmin
Definition: histogram.h:125

◆ index2val()

void Histogram1D::index2val ( double  i,
double &  v 
) const
inline

Index –> Value

Given the code of one interval, this function returns the value of its starting point (its left border point). If the intervals are defined as [a,b), this function returns a.

hist.index2val(0, beginning_of_interval_0);

Definition at line 295 of file histogram.h.

296  {
297  v = hmin + i * step_size;
298  }
double hmin
Definition: histogram.h:125
#define i
double step_size
Definition: histogram.h:127

◆ init()

void Histogram1D::init ( double  min_val,
double  max_val,
int  n_steps 
)

Initialisation of the histogram

This is the operation which allows the histogram to be used. This should be performed before inserting any value in it. The information given to this initialisation is the range within which the values will be counted, and the number of steps (discrete bars) in this range. If the value is outside this range it will not be taken into account although we have asked for its insertion in the histogram.

hist.init(-3, 3, 100);
// 100 steps in the range -3...3

Definition at line 71 of file histogram.cpp.

72 {
73  hmin = min_val;
74  hmax = max_val;
75  step_size = (double) (max_val - min_val) / (double) n_steps; // CO: n_steps-1->n_steps
76  istep_size = 1.0 / step_size;
78  no_samples = 0;
79 }
double hmin
Definition: histogram.h:125
double hmax
Definition: histogram.h:126
int no_samples
Definition: histogram.h:129
double step_size
Definition: histogram.h:127
double istep_size
Definition: histogram.h:128

◆ insert_value()

void Histogram1D::insert_value ( double  val)

Insert a value within histogram

The right interval is chosen according to the initialisation of the histogram and the count of elements in that interval is incremented by 1. If the value lies outside the global range of the histogram nothing is done.

hist.insert_value(3);

Definition at line 83 of file histogram.cpp.

84 {
85  int i;
86  int Xdim=(int)XSIZE(*this);
87 
88  // The following code is equivalent to val2index(val, i);
89  if (val == hmax)
90  {
91  i = Xdim - 1;
92  ++DIRECT_A1D_ELEM(*this, i);
93  ++no_samples;
94  }
95  else
96  {
97  double aux = (val - hmin) * istep_size;
98  i = (int) aux;
99 
100  if (i < 0 || i >= Xdim)
101  return; // the value is outside our scope
102 
103  ++DIRECT_A1D_ELEM(*this, i);
104  ++no_samples;
105  }
106 #ifdef DEBUG
107 
108  std::cout << " hmin " << hmin << " hmax " << hmax << " value " << val
109  << " index " << i << " (step_size= " << step_size << ")" << std::endl;
110 #endif
111 }
double hmin
Definition: histogram.h:125
double hmax
Definition: histogram.h:126
#define i
#define DIRECT_A1D_ELEM(v, i)
#define XSIZE(v)
int no_samples
Definition: histogram.h:129
double step_size
Definition: histogram.h:127
double istep_size
Definition: histogram.h:128

◆ mass_above()

double Histogram1D::mass_above ( double  value)
inline

Mass above

Returns the number of points which are above a certain value

Definition at line 245 of file histogram.h.

246  {
247  return no_samples - mass_below(value);
248  }
double mass_below(double value)
Definition: histogram.cpp:215
int no_samples
Definition: histogram.h:129

◆ mass_below()

double Histogram1D::mass_below ( double  value)

Mass below

Returns the number of points which are below a certain value

Definition at line 215 of file histogram.cpp.

216 {
217  // Trivial cases
218  if (value <= hmin)
219  return 0;
220  if (value >= hmax)
221  return no_samples;
222 
223  // Any other case, find index of corresponding piece
224  int i = 0;
225  double acc = 0;
226  double current_value;
227  index2val(i, current_value);
228  while (current_value <= value)
229  {
230  acc += A1D_ELEM(*this, i);
231  i++;
232  index2val(i, current_value);
233  }
234  return acc;
235 }
#define A1D_ELEM(v, i)
double hmin
Definition: histogram.h:125
double hmax
Definition: histogram.h:126
#define i
void index2val(double i, double &v) const
Definition: histogram.h:295
int current_value
Definition: svm-toy.cpp:54
int no_samples
Definition: histogram.h:129

◆ operator=()

Histogram1D & Histogram1D::operator= ( const Histogram1D H)

Assignment

Definition at line 51 of file histogram.cpp.

52 {
53  if (this != &H)
54  {
56  hmin = H.hmin;
57  hmax = H.hmax;
58  step_size = H.step_size;
61  }
62  return *this;
63 }
double hmin
Definition: histogram.h:125
double hmax
Definition: histogram.h:126
MultidimArray< T > & operator=(const MultidimArray< T > &op1)
int no_samples
Definition: histogram.h:129
double step_size
Definition: histogram.h:127
double istep_size
Definition: histogram.h:128

◆ percentil()

double Histogram1D::percentil ( double  percent_mass)

Returns the percentil value

This function returns the value within the range for which a given percent mass of the total number of elements are below it. For instance, if we have 120 values distributed between 0 and 45, and we ask for the percentil of 60%, the returned value is that within 0 and 45 for which 120 * 0.6 = 72 elements are smaller than it.

percentil60=hist.percentil(60);

Definition at line 160 of file histogram.cpp.

161 {
162  int i = 0;
163  double acc = 0;
164  double required_mass;
165  double percentil_i;
166  double ret_val;
167 
168  // Check it is a correct mass
169  if (percent_mass > 100)
170  REPORT_ERROR(ERR_VALUE_INCORRECT, "Asked for a percentil greater than 100");
171 
172  // Trivial cases
173  if (percent_mass == 0)
174  return (hmin);
175  if (percent_mass == 100)
176  return (hmax);
177 
178  // Any other case, find index of corresponding piece
179  required_mass = (double) no_samples * percent_mass / 100.0;
180  int N_diff_from_0 = 0;
181  while (acc < required_mass)
182  {
183  acc += A1D_ELEM(*this, i);
184  if (A1D_ELEM(*this, i) > 0)
185  N_diff_from_0++;
186  i++;
187  }
188 
189  // If the sum is just the one we want OK
190  if (acc == required_mass)
191  percentil_i = i;
192  // If there is only one sample different from 0
193  // then there is no way of setting the threshold in the middle
194  // Let's put it at the beginning of the bin
195  else if (N_diff_from_0 == 1)
196  percentil_i = i - 1;
197  // If not, then go back a step and compute which fraction of the
198  // bar is needed to finish the required mass
199  else
200  {
201  /* CO: We cannot assure that there is at least what is supposed to be
202  above this threshold. Let's move to the safe side
203  i--;
204  acc -= A1D_ELEM(*this,i);
205  percentil_i=i+(required_mass-acc)/(double) A1D_ELEM(*this,i); */
206  percentil_i = i - 1;
207  }
208 
209  // Now translate from index to range
210  index2val(percentil_i, ret_val);
211  return ret_val;
212 }
#define REPORT_ERROR(nerr, ErrormMsg)
Definition: xmipp_error.h:211
#define A1D_ELEM(v, i)
double hmin
Definition: histogram.h:125
double hmax
Definition: histogram.h:126
#define i
void index2val(double i, double &v) const
Definition: histogram.h:295
int no_samples
Definition: histogram.h:129
Incorrect value received.
Definition: xmipp_error.h:195

◆ sampleNo()

double Histogram1D::sampleNo ( ) const
inline

Number of samples introduced in the histogram

std::cout << "No. Samples in the histogram " << hist.sampleNo() << std::endl;

Definition at line 350 of file histogram.h.

351  {
352  return no_samples;
353  }
int no_samples
Definition: histogram.h:129

◆ step()

double Histogram1D::step ( ) const
inline

Step size for the histogram

std::cout << "Step size of the histogram " << hist.step() << std::endl;

Definition at line 328 of file histogram.h.

329  {
330  return step_size;
331  }
double step_size
Definition: histogram.h:127

◆ stepNo()

int Histogram1D::stepNo ( ) const
inline

Number of steps in the histogram

std::cout << "No. Steps in the histogram " << hist.stepNo() << std::endl;

Definition at line 339 of file histogram.h.

340  {
341  return XSIZE(*this);
342  }
#define XSIZE(v)

◆ val2index()

void Histogram1D::val2index ( double  v,
int &  i 
) const
inline

Value –> Index

Given a value it returns the code of the interval where it should be counted. If it is outside the global range of the histogram it returns -1

hist.val2index(12.3, interval_for_it);

Definition at line 271 of file histogram.h.

272  {
273  if (v == hmax)
274  i = XSIZE(*this) - 1;
275  else
276  {
277  double aux=(v - hmin) * istep_size;
278  i = (int) FLOOR(aux);
279  }
280 
281  if (i < 0 || i >= (int)XSIZE(*this))
282  i = -1;
283  }
double hmin
Definition: histogram.h:125
double hmax
Definition: histogram.h:126
#define i
#define FLOOR(x)
Definition: xmipp_macros.h:240
#define XSIZE(v)
double istep_size
Definition: histogram.h:128

◆ write()

void Histogram1D::write ( const FileName fn,
MDLabel  mdlValue = MDL_X,
MDLabel  mdlCount = MDL_COUNT 
)

Write an histogram to disk.

Definition at line 129 of file histogram.cpp.

132 {
133  MetaDataVec auxMD;
134  MDRowVec row;
135  double auxD;
136  size_t auxT;
138  {
139  index2val(i, auxD);
140  this->index2val(i, auxD);
141  row.setValue(mdlValue, auxD);
142  auxT=(size_t)A1D_ELEM(*this, i);
143  row.setValue(mdlCount, auxT);
144  auxMD.addRow(row);
145  }
146  auxMD.write(fn);
147  // std::ofstream fh;
148  // fh.open(fn.c_str(), std::ios::out);
149  // if (!fh)
150  // REPORT_ERROR(ERR_IO_NOTOPEN, (std::string)"Histogram1D::write: File " + fn +
151  // " cannot be openned for output");
152  // fh << "; Value Count\n";
153  // fh << *this;
154  // fh.close();
155 }
void setValue(const MDObject &object) override
#define A1D_ELEM(v, i)
void write(const FileName &outFile, WriteModeMetaData mode=MD_OVERWRITE) const
#define i
size_t addRow(const MDRow &row) override
void index2val(double i, double &v) const
Definition: histogram.h:295
#define FOR_ALL_ELEMENTS_IN_ARRAY1D(v)

Friends And Related Function Documentation

◆ operator<<

std::ostream& operator<< ( std::ostream &  o,
const Histogram1D hist 
)
friend

Show an histogram

The first column is the value associated to each histogram measure. The second one is the histogram measure.

Definition at line 115 of file histogram.cpp.

116 {
118  aux.resize(hist.stepNo(), 2);
120  {
121  hist.index2val(i, A2D_ELEM(aux, i, 0));
122  A2D_ELEM(aux, i, 1) = A1D_ELEM(hist, i);
123  }
124  o << aux;
125  return o;
126 }
#define A2D_ELEM(v, i, j)
void resize(size_t Ndim, size_t Zdim, size_t Ydim, size_t Xdim, bool copy=true)
#define A1D_ELEM(v, i)
#define i
int stepNo() const
Definition: histogram.h:339
void index2val(double i, double &v) const
Definition: histogram.h:295
#define FOR_ALL_ELEMENTS_IN_ARRAY1D(v)

Member Data Documentation

◆ hmax

double Histogram1D::hmax

Definition at line 126 of file histogram.h.

◆ hmin

double Histogram1D::hmin

Definition at line 125 of file histogram.h.

◆ istep_size

double Histogram1D::istep_size

Definition at line 128 of file histogram.h.

◆ no_samples

int Histogram1D::no_samples

Definition at line 129 of file histogram.h.

◆ step_size

double Histogram1D::step_size

Definition at line 127 of file histogram.h.


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