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

#include <histogram.h>

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

Public Member Functions

 Histogram2D ()
 
 Histogram2D (const Histogram2D &H)
 
void clear ()
 
Histogram2Doperator= (const Histogram2D &H)
 
void assign (const Histogram2D &H)
 
void init (double imin_val, double imax_val, int in_steps, double jmin_val, double jmax_val, int jn_steps)
 
void insert_value (double v, double u)
 
void write (const FileName &fn)
 
void val2index (double v, double u, int &i, int &j) const
 
void index2val (double i, double j, double &v, double &u) const
 
double Ihist_min () const
 
double Ihist_max () const
 
double Istep () const
 
int IstepNo () const
 
double Jhist_min () const
 
double Jhist_max () const
 
double Jstep () const
 
int JstepNo () const
 
int sampleNo () 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 imin
 
double imax
 
double istep_size
 
double jmin
 
double jmax
 
double jstep_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 Histogram2D &hist)
 

Detailed Description

Histograms with 2 parameters

The histogram with 2 parameters can be regarded as an approximation to the joint probability density function of two variables (just dividing the histogram by its total mass). Ie, the 2D histogram is the count of times that a certain combination of values in two variables has occurred. For example, this 2D histograms can be used to plot the projection distribution over the topological sphere, in such situation only the first two Euler angles are interesting and we could plot how many projections are there with first angle equal to 45 degrees and second equal to 0 degrees, and so on covering the whole range for each angle.

The 2D histogram is defined, as in the 1D case, by the respective ranges for the two variables and the number of intervals on each range. The distribution and limits of intervals are just the 2D extension of the graph shown in histograms 1D.

Definition at line 699 of file histogram.h.

Constructor & Destructor Documentation

◆ Histogram2D() [1/2]

Histogram2D::Histogram2D ( )
inline

Empty constructor

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

Definition at line 720 of file histogram.h.

721  {
722  clear();
723  }
void clear()
Definition: histogram.cpp:476

◆ Histogram2D() [2/2]

Histogram2D::Histogram2D ( const Histogram2D H)
inline

Copy constructor

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

Histogram2D hist2(hist1);

Definition at line 733 of file histogram.h.

734  {
735  *this = H;
736  }

Member Function Documentation

◆ assign()

void Histogram2D::assign ( const Histogram2D H)

Another function for assigment.

Definition at line 506 of file histogram.cpp.

507 {
508  *this = H;
509 }

◆ clear()

void Histogram2D::clear ( )
virtual

Empties an histogram

Remember to initialise it before using it again.

hist.clear();

Implements MultidimArrayBase.

Definition at line 476 of file histogram.cpp.

477 {
478  imin = 0;
479  imax = 0;
480  istep_size = 0;
481  jmin = 0;
482  jmax = 0;
483  jstep_size = 0;
484  no_samples = 0;
486 }
double istep_size
Definition: histogram.h:705
double jmax
Definition: histogram.h:707
double jstep_size
Definition: histogram.h:708
double imax
Definition: histogram.h:704
int no_samples
Definition: histogram.h:709
double jmin
Definition: histogram.h:706
double imin
Definition: histogram.h:703

◆ Ihist_max()

double Histogram2D::Ihist_max ( ) const
inline

Maximum i value where the histogram is defined

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

Definition at line 866 of file histogram.h.

867  {
868  return imax;
869  }
double imax
Definition: histogram.h:704

◆ Ihist_min()

double Histogram2D::Ihist_min ( ) const
inline

Minimum i value where the histogram is defined

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

Definition at line 855 of file histogram.h.

856  {
857  return imin;
858  }
double imin
Definition: histogram.h:703

◆ index2val()

void Histogram2D::index2val ( double  i,
double  j,
double &  v,
double &  u 
) const
inline

Index –> Value

Given the code of one interval, this function returns the value of its starting point (its left-top border point, ie, its lowest corner). If the intervals are defined as [v0,vF) and [u0,uF), this function returns the point (v0,u0)

hist.index2val(5, 1, v, u);

Definition at line 843 of file histogram.h.

844  {
845  v = imin + i * istep_size;
846  u = jmin + j * jstep_size;
847  }
double istep_size
Definition: histogram.h:705
#define i
double jstep_size
Definition: histogram.h:708
#define j
doublereal * u
double jmin
Definition: histogram.h:706
double imin
Definition: histogram.h:703

◆ init()

void Histogram2D::init ( double  imin_val,
double  imax_val,
int  in_steps,
double  jmin_val,
double  jmax_val,
int  jn_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 for each variable within which the values will be counted, and the number of steps (discrete bars) in these ranges.If the value is outside the defined ranges it will not be taken into account although we have asked for its insertion in the histogram.

hist.init(0, 90, 100, 0, 360, 200);
// 100 steps in the range V=0...90 and 200 steps for U=0...360

Definition at line 512 of file histogram.cpp.

514 {
515  // V axis
516  imin = imin_val;
517  imax = imax_val;
518  istep_size = (double) (imax_val - imin_val) / (double) in_steps;
519 
520  // U axis
521  jmin = jmin_val;
522  jmax = jmax_val;
523  jstep_size = (double) (jmax_val - jmin_val) / (double) jn_steps;
524 
525  initZeros(in_steps, jn_steps);
526  no_samples = 0;
527 }
double istep_size
Definition: histogram.h:705
double jmax
Definition: histogram.h:707
double jstep_size
Definition: histogram.h:708
double imax
Definition: histogram.h:704
int no_samples
Definition: histogram.h:709
double jmin
Definition: histogram.h:706
double imin
Definition: histogram.h:703

◆ insert_value()

void Histogram2D::insert_value ( double  v,
double  u 
)

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. Notice that two values are needed, these are the two values of the two variables in our example of the projection at with first Euler angle=45 and second=0, the insertion of this projection in the 2D histogram would be like in the following example.

hist.insert_value(45, 0);

Definition at line 530 of file histogram.cpp.

531 {
532  int i, j;
533  val2index(v, u, i, j);
534  if (i == -1 || j == -1)
535  return; // it is outside our scope
536  int Xdim=(int)XSIZE(*this);
537  int Ydim=(int)YSIZE(*this);
538  i = CLIP(i, 0, Ydim);
539  j = CLIP(j, 0, Xdim);
540  A2D_ELEM(*this, i, j)++;
541  no_samples++;
542 }
#define YSIZE(v)
#define A2D_ELEM(v, i, j)
#define i
#define CLIP(x, x0, xF)
Definition: xmipp_macros.h:260
void val2index(double v, double u, int &i, int &j) const
Definition: histogram.h:813
#define XSIZE(v)
#define j
int no_samples
Definition: histogram.h:709
doublereal * u

◆ Istep()

double Histogram2D::Istep ( ) const
inline

Step size in i for the histogram

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

Definition at line 877 of file histogram.h.

878  {
879  return istep_size;
880  }
double istep_size
Definition: histogram.h:705

◆ IstepNo()

int Histogram2D::IstepNo ( ) const
inline

Number of steps in i in the histogram

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

Definition at line 888 of file histogram.h.

889  {
890  return YSIZE(*this);
891  }
#define YSIZE(v)

◆ Jhist_max()

double Histogram2D::Jhist_max ( ) const
inline

Maximum j value where the histogram is defined

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

Definition at line 910 of file histogram.h.

911  {
912  return jmax;
913  }
double jmax
Definition: histogram.h:707

◆ Jhist_min()

double Histogram2D::Jhist_min ( ) const
inline

Minimum j value where the histogram is defined

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

Definition at line 899 of file histogram.h.

900  {
901  return jmin;
902  }
double jmin
Definition: histogram.h:706

◆ Jstep()

double Histogram2D::Jstep ( ) const
inline

Step size in j for the histogram

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

Definition at line 921 of file histogram.h.

922  {
923  return jstep_size;
924  }
double jstep_size
Definition: histogram.h:708

◆ JstepNo()

int Histogram2D::JstepNo ( ) const
inline

Number of steps in j in the histogram

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

Definition at line 932 of file histogram.h.

933  {
934  return XSIZE(*this);
935  }
#define XSIZE(v)

◆ operator=()

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

Assignment.

Definition at line 489 of file histogram.cpp.

490 {
491  if (this != &H)
492  {
494  imin = H.imin;
495  imax = H.imax;
497  jmin = H.jmin;
498  jmax = H.jmax;
501  }
502  return *this;
503 }
double istep_size
Definition: histogram.h:705
double jmax
Definition: histogram.h:707
double jstep_size
Definition: histogram.h:708
double imax
Definition: histogram.h:704
MultidimArray< T > & operator=(const MultidimArray< T > &op1)
int no_samples
Definition: histogram.h:709
double jmin
Definition: histogram.h:706
double imin
Definition: histogram.h:703

◆ sampleNo()

int Histogram2D::sampleNo ( ) const
inline

Number of samples introduced in the histogram

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

Definition at line 943 of file histogram.h.

944  {
945  return no_samples;
946  }
int no_samples
Definition: histogram.h:709

◆ val2index()

void Histogram2D::val2index ( double  v,
double  u,
int &  i,
int &  j 
) const
inline

Value –> Index

Given two values for the two variables 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 for that variable (i or j). The following example tells you that the interval corresponding to (45,0) is that with code (i,j), i.e, you could access to hist()(i,j) to know how many projections are there in the same interval as this projection.

hist.val2index(45, 0, i, j);

Definition at line 813 of file histogram.h.

814  {
815  if (v == imax)
816  i = IstepNo() - 1;
817  else
818  i = (int) FLOOR((v - imin) / istep_size);
819 
820  if (u == jmax)
821  j = JstepNo() - 1;
822 
823  j = (int) FLOOR((u - jmin) / jstep_size);
824 
825  if (i < 0 || i >= IstepNo())
826  i = -1;
827 
828  if (j < 0 || j >= JstepNo())
829  j = -1;
830  }
double istep_size
Definition: histogram.h:705
int IstepNo() const
Definition: histogram.h:888
int JstepNo() const
Definition: histogram.h:932
double jmax
Definition: histogram.h:707
#define i
double jstep_size
Definition: histogram.h:708
#define FLOOR(x)
Definition: xmipp_macros.h:240
double imax
Definition: histogram.h:704
#define j
doublereal * u
double jmin
Definition: histogram.h:706
double imin
Definition: histogram.h:703

◆ write()

void Histogram2D::write ( const FileName fn)

Write an histogram to disk

Definition at line 561 of file histogram.cpp.

562 {
563  std::ofstream fh;
564  fh.open(fn.c_str(), std::ios::out);
565  if (!fh)
566  REPORT_ERROR(ERR_IO_NOTOPEN, (std::string)"histogram2D::write: File " + fn + " cannot be openned for output");
567  fh << *this;
568  fh.close();
569 }
#define REPORT_ERROR(nerr, ErrormMsg)
Definition: xmipp_error.h:211
File cannot be open.
Definition: xmipp_error.h:137

Friends And Related Function Documentation

◆ operator<<

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

Show an histogram

The first column and second column are the (X,Y) coordinates of each histogram measure. The third one is the histogram measure.

Definition at line 545 of file histogram.cpp.

546 {
548  aux.resize(hist.IstepNo() * hist.JstepNo(), 3);
549  int n = 0;
551  {
552  hist.index2val(i, j, A2D_ELEM(aux, n, 0), A2D_ELEM(aux, n, 1));
553  A2D_ELEM(aux, n, 2) = A2D_ELEM(hist, i, j);
554  n++;
555  }
556  o << aux;
557  return o;
558 }
#define A2D_ELEM(v, i, j)
void resize(size_t Ndim, size_t Zdim, size_t Ydim, size_t Xdim, bool copy=true)
int IstepNo() const
Definition: histogram.h:888
int JstepNo() const
Definition: histogram.h:932
void index2val(double i, double j, double &v, double &u) const
Definition: histogram.h:843
#define i
#define FOR_ALL_ELEMENTS_IN_ARRAY2D(m)
#define j
int * n

Member Data Documentation

◆ imax

double Histogram2D::imax

Definition at line 704 of file histogram.h.

◆ imin

double Histogram2D::imin

Definition at line 703 of file histogram.h.

◆ istep_size

double Histogram2D::istep_size

Definition at line 705 of file histogram.h.

◆ jmax

double Histogram2D::jmax

Definition at line 707 of file histogram.h.

◆ jmin

double Histogram2D::jmin

Definition at line 706 of file histogram.h.

◆ jstep_size

double Histogram2D::jstep_size

Definition at line 708 of file histogram.h.

◆ no_samples

int Histogram2D::no_samples

Definition at line 709 of file histogram.h.


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