Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TProfile2D Class Reference

Profile2D histograms are used to display the mean value of Z and its error for each cell in X,Y.

Profile2D histograms are in many cases an elegant replacement of three-dimensional histograms : the inter-relation of three measured quantities X, Y and Z can always be visualized by a three-dimensional histogram or scatter-plot; its representation on the line-printer is not particularly satisfactory, except for sparse data. If Z is an unknown (but single-valued) approximate function of X,Y this function is displayed by a profile2D histogram with much better precision than by a scatter-plot.

The following formulae show the cumulated contents (capital letters) and the values displayed by the printing or plotting routines (small letters) of the elements for cell i, j.

\[ \begin{align} H(i,j) &= \sum w \cdot Z \\ E(i,j) &= \sum w \cdot Z^2 \\ W(i,j) &= \sum w \\ h(i,j) &= \frac{H(i,j)}{W(i,j)} \\ s(i,j) &= \sqrt{E(i,j)/W(i,j)- h(i,j)^2} \\ e(i,j) &= \frac{s(i,j)}{\sqrt{W(i,j)}} \end{align} \]

The bin content is always the mean of the Z values, but errors change depending on options:

\[ \begin{align} \text{GetBinContent}(i,j) &= h(i,j) \\ \text{GetBinError}(i,j) &= \begin{cases} e(i,j) &\text{if option="" (default). Error of the mean of all z values.} \\ s(i,j) &\text{if option="s". Standard deviation of z values.} \\ \begin{cases} e(j) &\text{if } h(i,j) \ne 0 \\ 1/\sqrt{12 N} &\text{if } h(i,j)=0 \end{cases} &\text{if option="i". This is useful for storing integers such as ADC counts.} \\ 1/\sqrt{W(i,j)} &\text{if option="g". Error of a weighted mean when combining measurements with variances of } w. \\ \end{cases} \end{align} \]

In the special case where s(I,J) is zero (eg, case of 1 entry only in one cell) the bin error e(I,J) is computed from the average of the s(I,J) for all cells if the static function TProfile2D::Approximate has been called. This simple/crude approximation was suggested in order to keep the cell during a fit operation. But note that this approximation is not the default behaviour.

Creating and drawing a 2D profile

{
auto c1 = new TCanvas("c1","Profile histogram example",200,10,700,500);
auto hprof2d = new TProfile2D("hprof2d","Profile of pz versus px and py",40,-4,4,40,-4,4,0,20);
Float_t px, py, pz;
for ( Int_t i=0; i<25000; i++) {
gRandom->Rannor(px,py);
pz = px*px + py*py;
hprof2d->Fill(px,py,pz,1);
}
hprof2d->Draw();
}
float Float_t
Definition RtypesCore.h:57
R__EXTERN TRandom * gRandom
Definition TRandom.h:62
The Canvas class.
Definition TCanvas.h:23
TProfile2D()
Default constructor for Profile2D histograms.
virtual void Rannor(Float_t &a, Float_t &b)
Return 2 numbers distributed following a gaussian with mean=0 and sigma=1.
Definition TRandom.cxx:500
return c1
Definition legend1.C:41

Definition at line 27 of file TProfile2D.h.

Public Member Functions

 TProfile2D ()
 Default constructor for Profile2D histograms.
 
 TProfile2D (const char *name, const char *title, Int_t nbinsx, const Double_t *xbins, Int_t nbinsy, const Double_t *ybins, Option_t *option="")
 Create a 2-D Profile with variable bins in X and variable bins in Y.
 
 TProfile2D (const char *name, const char *title, Int_t nbinsx, const Double_t *xbins, Int_t nbinsy, Double_t ylow, Double_t yup, Option_t *option="")
 Create a 2-D Profile with variable bins in X and fix bins in Y.
 
 TProfile2D (const char *name, const char *title, Int_t nbinsx, Double_t xlow, Double_t xup, Int_t nbinsy, const Double_t *ybins, Option_t *option="")
 Create a 2-D Profile with fix bins in X and variable bins in Y.
 
 TProfile2D (const char *name, const char *title, Int_t nbinsx, Double_t xlow, Double_t xup, Int_t nbinsy, Double_t ylow, Double_t yup, Double_t zlow, Double_t zup, Option_t *option="")
 Constructor for Profile2D histograms with range in z.
 
 TProfile2D (const char *name, const char *title, Int_t nbinsx, Double_t xlow, Double_t xup, Int_t nbinsy, Double_t ylow, Double_t yup, Option_t *option="")
 Normal Constructor for Profile histograms.
 
 TProfile2D (const TProfile2D &profile)
 Copy constructor.
 
 ~TProfile2D () override
 Default destructor for Profile2D histograms.
 
Bool_t Add (const TH1 *h1, const TH1 *h2, Double_t c1=1, Double_t c2=1) override
 Replace contents of this profile2D by the addition of h1 and h2.
 
Bool_t Add (const TH1 *h1, Double_t c1=1) override
 Performs the operation: this = this + c1*h1 .
 
Bool_t Add (TF1 *h1, Double_t c1=1, Option_t *option="") override
 Performs the operation: this = this + c1*f1 .
 
Int_t BufferEmpty (Int_t action=0) override
 Fill histogram with all entries in the buffer.
 
void BuildOptions (Double_t zmin, Double_t zmax, Option_t *option)
 Set Profile2D histogram structure and options.
 
void Copy (TObject &hnew) const override
 Copy a Profile2D histogram to a new profile2D histogram.
 
Bool_t Divide (const TH1 *h1) override
 Divide this profile2D by h1.
 
Bool_t Divide (const TH1 *h1, const TH1 *h2, Double_t c1=1, Double_t c2=1, Option_t *option="") override
 Replace contents of this profile2D by the division of h1 by h2.
 
Bool_t Divide (TF1 *h1, Double_t c1=1) override
 Performs the operation: this = this/(c1*f1) .
 
void ExtendAxis (Double_t x, TAxis *axis) override
 Profile histogram is resized along axis such that x is in the axis range.
 
virtual Int_t Fill (const char *namex, const char *namey, Double_t z, Double_t w=1.)
 Fill a Profile2D histogram (no weights).
 
virtual Int_t Fill (const char *namex, Double_t y, Double_t z, Double_t w=1.)
 Fill a Profile2D histogram (no weights).
 
virtual Int_t Fill (Double_t x, const char *namey, Double_t z, Double_t w=1.)
 Fill a Profile2D histogram (no weights).
 
Int_t Fill (Double_t x, Double_t y, Double_t z) override
 Fill a Profile2D histogram (no weights).
 
virtual Int_t Fill (Double_t x, Double_t y, Double_t z, Double_t w)
 Fill a Profile2D histogram with weights.
 
Double_t GetBinContent (Int_t bin) const override
 Return bin content of a Profile2D histogram.
 
Double_t GetBinContent (Int_t binx, Int_t biny) const override
 
Double_t GetBinContent (Int_t binx, Int_t biny, Int_t) const override
 
virtual Double_t GetBinEffectiveEntries (Int_t bin)
 Return bin effective entries for a weighted filled Profile histogram.
 
virtual Double_t GetBinEntries (Int_t bin) const
 Return bin entries of a Profile2D histogram.
 
Double_t GetBinError (Int_t bin) const override
 Return bin error of a Profile2D histogram.
 
Double_t GetBinError (Int_t binx, Int_t biny) const override
 
Double_t GetBinError (Int_t binx, Int_t biny, Int_t) const override
 
virtual TArrayDGetBinSumw2 ()
 
virtual const TArrayDGetBinSumw2 () const
 
Option_tGetErrorOption () const
 Return option to compute profile2D errors.
 
Double_t GetNumberOfBins ()
 
void GetStats (Double_t *stats) const override
 Fill the array stats from the contents of this profile.
 
virtual Double_t GetZmax () const
 
virtual Double_t GetZmin () const
 
TClassIsA () const override
 
void LabelsDeflate (Option_t *axis="X") override
 Reduce the number of bins for this axis to the number of bins having a label.
 
void LabelsInflate (Option_t *axis="X") override
 Double the number of bins for axis.
 
void LabelsOption (Option_t *option="h", Option_t *axis="X") override
 Set option(s) to draw axis with labels.
 
Long64_t Merge (TCollection *list) override
 Merge all histograms in the collection in this histogram.
 
Bool_t Multiply (const TH1 *h1) override
 Multiply this profile2D by h1.
 
Bool_t Multiply (const TH1 *h1, const TH1 *h2, Double_t c1=1, Double_t c2=1, Option_t *option="") override
 Replace contents of this profile2D by multiplication of h1 by h2.
 
Bool_t Multiply (TF1 *h1, Double_t c1=1) override
 Performs the operation: this = this*c1*f1.
 
TProfile2Doperator= (const TProfile2D &profile)
 
TProfileProfileX (const char *name="_pfx", Int_t firstybin=0, Int_t lastybin=-1, Option_t *option="") const
 Project a 2-D histogram into a profile histogram along X.
 
TProfileProfileY (const char *name="_pfy", Int_t firstxbin=0, Int_t lastxbin=-1, Option_t *option="") const
 Project a 2-D histogram into a profile histogram along X.
 
TH2DProjectionXY (const char *name="_pxy", Option_t *option="e") const
 Project this profile2D into a 2-D histogram along X,Y.
 
void PutStats (Double_t *stats) override
 Replace current statistics with the values in array stats.
 
TProfile2DRebin2D (Int_t nxgroup=2, Int_t nygroup=2, const char *newname="") override
 Rebin this histogram grouping nxgroup/nygroup bins along the xaxis/yaxis together.
 
TProfile2DRebinX (Int_t ngroup=2, const char *newname="") override
 Rebin only the X axis.
 
TProfile2DRebinY (Int_t ngroup=2, const char *newname="") override
 Rebin only the Y axis.
 
void Reset (Option_t *option="") override
 Reset contents of a Profile2D histogram.
 
void SavePrimitive (std::ostream &out, Option_t *option="") override
 Save primitive as a C++ statement(s) on output stream out.
 
void Scale (Double_t c1=1, Option_t *option="") override
 Multiply this profile2D by a constant c1.
 
virtual void SetBinEntries (Int_t bin, Double_t w)
 Set the number of entries in bin.
 
void SetBins (Int_t nbinsx, Double_t xmin, Double_t xmax, Int_t nbinsy, Double_t ymin, Double_t ymax) override
 Redefine x and y axis parameters.
 
void SetBins (Int_t nx, const Double_t *xBins, Int_t ny, const Double_t *yBins) override
 Redefine x and y axis parameters for variable bin sizes.
 
void SetBinsLength (Int_t n=-1) override
 Set total number of bins including under/overflow.
 
void SetBuffer (Int_t buffersize, Option_t *option="") override
 Set the buffer size in units of 8 bytes (double).
 
virtual void SetErrorOption (Option_t *option="")
 Set option to compute profile2D errors.
 
void Streamer (TBuffer &) override
 Stream an object of class TProfile2D.
 
void StreamerNVirtual (TBuffer &ClassDef_StreamerNVirtual_b)
 
void Sumw2 (Bool_t flag=kTRUE) override
 Create/Delete structure to store sum of squares of weights per bin.
 
- Public Member Functions inherited from TH2D
 TH2D ()
 Constructor.
 
 TH2D (const char *name, const char *title, Int_t nbinsx, const Double_t *xbins, Int_t nbinsy, const Double_t *ybins)
 Constructor (see TH2::TH2 for explanation of parameters)
 
 TH2D (const char *name, const char *title, Int_t nbinsx, const Double_t *xbins, Int_t nbinsy, Double_t ylow, Double_t yup)
 Constructor (see TH2::TH2 for explanation of parameters)
 
 TH2D (const char *name, const char *title, Int_t nbinsx, const Float_t *xbins, Int_t nbinsy, const Float_t *ybins)
 Constructor (see TH2::TH2 for explanation of parameters)
 
 TH2D (const char *name, const char *title, Int_t nbinsx, Double_t xlow, Double_t xup, Int_t nbinsy, const Double_t *ybins)
 Constructor (see TH2::TH2 for explanation of parameters)
 
 TH2D (const char *name, const char *title, Int_t nbinsx, Double_t xlow, Double_t xup, Int_t nbinsy, Double_t ylow, Double_t yup)
 Constructor (see TH2::TH2 for explanation of parameters)
 
 TH2D (const TH2D &h2d)
 Copy constructor.
 
 TH2D (const TMatrixDBase &m)
 Constructor Construct a 2-D histogram from a TMatrixDBase.
 
 ~TH2D () override
 Destructor.
 
void AddBinContent (Int_t bin) override
 Increment bin content by 1.
 
void AddBinContent (Int_t bin, Double_t w) override
 Increment bin content by a weight w.
 
TH2Doperator= (const TH2D &h1)
 Operator =.
 
void StreamerNVirtual (TBuffer &ClassDef_StreamerNVirtual_b)
 
- Public Member Functions inherited from TH2
 ~TH2 () override
 Destructor.
 
void FillN (Int_t ntimes, const Double_t *x, const Double_t *y, const Double_t *w, Int_t stride=1) override
 Fill a 2-D histogram with an array of values and weights.
 
void FillN (Int_t, const Double_t *, const Double_t *, Int_t) override
 Fill this histogram with an array x and weights w.
 
void FillRandom (const char *fname, Int_t ntimes=5000, TRandom *rng=nullptr) override
 Fill histogram following distribution in function fname.
 
void FillRandom (TH1 *h, Int_t ntimes=5000, TRandom *rng=nullptr) override
 Fill histogram following distribution in histogram h.
 
virtual void FitSlicesX (TF1 *f1=nullptr, Int_t firstybin=0, Int_t lastybin=-1, Int_t cut=0, Option_t *option="QNR", TObjArray *arr=nullptr)
 Project slices along X in case of a 2-D histogram, then fit each slice with function f1 and make a histogram for each fit parameter Only bins along Y between firstybin and lastybin are considered.
 
virtual void FitSlicesY (TF1 *f1=nullptr, Int_t firstxbin=0, Int_t lastxbin=-1, Int_t cut=0, Option_t *option="QNR", TObjArray *arr=nullptr)
 Project slices along Y in case of a 2-D histogram, then fit each slice with function f1 and make a histogram for each fit parameter Only bins along X between firstxbin and lastxbin are considered.
 
Int_t GetBin (Int_t binx, Int_t biny, Int_t binz=0) const override
 Return Global bin number corresponding to binx,y,z.
 
virtual Double_t GetBinErrorLow (Int_t bin) const
 Return lower error associated to bin number bin.
 
virtual Double_t GetBinErrorLow (Int_t binx, Int_t biny)
 
virtual Double_t GetBinErrorUp (Int_t bin) const
 Return upper error associated to bin number bin.
 
virtual Double_t GetBinErrorUp (Int_t binx, Int_t biny)
 
virtual Double_t GetBinWithContent2 (Double_t c, Int_t &binx, Int_t &biny, Int_t firstxbin=1, Int_t lastxbin=-1, Int_t firstybin=1, Int_t lastybin=-1, Double_t maxdiff=0) const
 compute first cell (binx,biny) in the range [firstxbin,lastxbin][firstybin,lastybin] for which diff = abs(cell_content-c) <= maxdiff In case several cells in the specified range with diff=0 are found the first cell found is returned in binx,biny.
 
virtual Double_t GetCorrelationFactor (Int_t axis1=1, Int_t axis2=2) const
 Return correlation factor between axis1 and axis2.
 
virtual Double_t GetCovariance (Int_t axis1=1, Int_t axis2=2) const
 Return covariance between axis1 and axis2.
 
virtual void GetRandom2 (Double_t &x, Double_t &y, TRandom *rng=nullptr)
 Return 2 random numbers along axis x and y distributed according to the cell-contents of this 2-D histogram return a NaN if the histogram has a bin with negative content.
 
virtual Double_t Integral (Int_t binx1, Int_t binx2, Int_t biny1, Int_t biny2, Option_t *option="") const
 Return integral of bin contents in range [firstxbin,lastxbin],[firstybin,lastybin] for a 2-D histogram By default the integral is computed as the sum of bin contents in the range.
 
virtual Double_t Integral (Int_t, Int_t, Int_t, Int_t, Int_t, Int_t, Option_t *="") const
 
Double_t Integral (Option_t *option="") const override
 Return integral of bin contents.
 
virtual Double_t IntegralAndError (Int_t binx1, Int_t binx2, Int_t biny1, Int_t biny2, Double_t &err, Option_t *option="") const
 Return integral of bin contents in range [firstxbin,lastxbin],[firstybin,lastybin] for a 2-D histogram.
 
Double_t Interpolate (Double_t x, Double_t y) const override
 Given a point P(x,y), Interpolate approximates the value via bilinear interpolation based on the four nearest bin centers see Wikipedia, Bilinear Interpolation Andy Mastbaum 10/8/2008 vaguely based on R.Raja 6-Sep-2008.
 
Double_t Interpolate (Double_t x, Double_t y, Double_t z) const override
 illegal for a TH2
 
Double_t KolmogorovTest (const TH1 *h2, Option_t *option="") const override
 Statistical test of compatibility in shape between THIS histogram and h2, using Kolmogorov test.
 
TProfileProfileX (const char *name="_pfx", Int_t firstybin=1, Int_t lastybin=-1, Option_t *option="") const
 Project a 2-D histogram into a profile histogram along X.
 
TProfileProfileY (const char *name="_pfy", Int_t firstxbin=1, Int_t lastxbin=-1, Option_t *option="") const
 Project a 2-D histogram into a profile histogram along Y.
 
TH1DProjectionX (const char *name="_px", Int_t firstybin=0, Int_t lastybin=-1, Option_t *option="") const
 Project a 2-D histogram into a 1-D histogram along X.
 
TH1DProjectionY (const char *name="_py", Int_t firstxbin=0, Int_t lastxbin=-1, Option_t *option="") const
 Project a 2-D histogram into a 1-D histogram along Y.
 
TH1DQuantilesX (Double_t prob=0.5, const char *name="_qx") const
 Compute the X distribution of quantiles in the other variable Y name is the name of the returned histogram prob is the probability content for the quantile (0.5 is the default for the median) An approximate error for the quantile is computed assuming that the distribution in the other variable is normal.
 
TH1DQuantilesY (Double_t prob=0.5, const char *name="_qy") const
 Compute the Y distribution of quantiles in the other variable X name is the name of the returned histogram prob is the probability content for the quantile (0.5 is the default for the median) An approximate error for the quantile is computed assuming that the distribution in the other variable is normal.
 
TH2Rebin (Int_t ngroup=2, const char *newname="", const Double_t *xbins=nullptr) override
 Override TH1::Rebin as TH2::RebinX Rebinning in variable binning as for TH1 is not allowed If a non-null pointer is given an error is flagged see RebinX and Rebin2D.
 
void SetBinContent (Int_t bin, Double_t content) override
 Set bin content.
 
void SetBinContent (Int_t binx, Int_t biny, Double_t content) override
 
void SetBinContent (Int_t binx, Int_t biny, Int_t, Double_t content) override
 
virtual void SetShowProjectionX (Int_t nbins=1)
 When the mouse is moved in a pad containing a 2-d view of this histogram a second canvas shows the projection along X corresponding to the mouse position along Y.
 
virtual void SetShowProjectionXY (Int_t nbinsY=1, Int_t nbinsX=1)
 When the mouse is moved in a pad containing a 2-d view of this histogram two canvases show the projection along X and Y corresponding to the mouse position along Y and X, respectively.
 
virtual void SetShowProjectionY (Int_t nbins=1)
 When the mouse is moved in a pad containing a 2-d view of this histogram a second canvas shows the projection along Y corresponding to the mouse position along X.
 
TH1ShowBackground (Int_t niter=20, Option_t *option="same") override
 This function calculates the background spectrum in this histogram.
 
Int_t ShowPeaks (Double_t sigma=2, Option_t *option="", Double_t threshold=0.05) override
 Interface to TSpectrum2::Search the function finds peaks in this histogram where the width is > sigma and the peak maximum greater than threshold*maximum bin content of this.
 
void Smooth (Int_t ntimes=1, Option_t *option="") override
 Smooth bin contents of this 2-d histogram using kernel algorithms similar to the ones used in the raster graphics community.
 
void StreamerNVirtual (TBuffer &ClassDef_StreamerNVirtual_b)
 
- Public Member Functions inherited from TH1
 ~TH1 () override
 Histogram default destructor.
 
virtual Double_t AndersonDarlingTest (const TH1 *h2, Double_t &advalue) const
 Same function as above but returning also the test statistic value.
 
virtual Double_t AndersonDarlingTest (const TH1 *h2, Option_t *option="") const
 Statistical test of compatibility in shape between this histogram and h2, using the Anderson-Darling 2 sample test.
 
void Browse (TBrowser *b) override
 Browse the Histogram object.
 
virtual Bool_t CanExtendAllAxes () const
 Returns true if all axes are extendable.
 
virtual Double_t Chi2Test (const TH1 *h2, Option_t *option="UU", Double_t *res=nullptr) const
 \( \chi^{2} \) test for comparing weighted and unweighted histograms
 
virtual Double_t Chi2TestX (const TH1 *h2, Double_t &chi2, Int_t &ndf, Int_t &igood, Option_t *option="UU", Double_t *res=nullptr) const
 The computation routine of the Chisquare test.
 
virtual Double_t Chisquare (TF1 *f1, Option_t *option="") const
 Compute and return the chisquare of this histogram with respect to a function The chisquare is computed by weighting each histogram point by the bin error By default the full range of the histogram is used.
 
virtual void ClearUnderflowAndOverflow ()
 Remove all the content from the underflow and overflow bins, without changing the number of entries After calling this method, every undeflow and overflow bins will have content 0.0 The Sumw2 is also cleared, since there is no more content in the bins.
 
TObjectClone (const char *newname="") const override
 Make a complete copy of the underlying object.
 
virtual Double_t ComputeIntegral (Bool_t onlyPositive=false)
 Compute integral (cumulative sum of bins) The result stored in fIntegral is used by the GetRandom functions.
 
virtual void DirectoryAutoAdd (TDirectory *)
 Perform the automatic addition of the histogram to the given directory.
 
Int_t DistancetoPrimitive (Int_t px, Int_t py) override
 Compute distance from point px,py to a line.
 
void Draw (Option_t *option="") override
 Draw this histogram with options.
 
virtual TH1DrawCopy (Option_t *option="", const char *name_postfix="_copy") const
 Copy this histogram and Draw in the current pad.
 
virtual TH1DrawNormalized (Option_t *option="", Double_t norm=1) const
 Draw a normalized copy of this histogram.
 
virtual void DrawPanel ()
 Display a panel with all histogram drawing options.
 
virtual void Eval (TF1 *f1, Option_t *option="")
 Evaluate function f1 at the center of bins of this histogram.
 
void ExecuteEvent (Int_t event, Int_t px, Int_t py) override
 Execute action corresponding to one event.
 
virtual TH1FFT (TH1 *h_output, Option_t *option)
 This function allows to do discrete Fourier transforms of TH1 and TH2.
 
virtual Int_t FindBin (Double_t x, Double_t y=0, Double_t z=0)
 Return Global bin number corresponding to x,y,z.
 
virtual Int_t FindFirstBinAbove (Double_t threshold=0, Int_t axis=1, Int_t firstBin=1, Int_t lastBin=-1) const
 Find first bin with content > threshold for axis (1=x, 2=y, 3=z) if no bins with content > threshold is found the function returns -1.
 
virtual Int_t FindFixBin (Double_t x, Double_t y=0, Double_t z=0) const
 Return Global bin number corresponding to x,y,z.
 
virtual Int_t FindLastBinAbove (Double_t threshold=0, Int_t axis=1, Int_t firstBin=1, Int_t lastBin=-1) const
 Find last bin with content > threshold for axis (1=x, 2=y, 3=z) if no bins with content > threshold is found the function returns -1.
 
TObjectFindObject (const char *name) const override
 Search object named name in the list of functions.
 
TObjectFindObject (const TObject *obj) const override
 Search object obj in the list of functions.
 
virtual TFitResultPtr Fit (const char *formula, Option_t *option="", Option_t *goption="", Double_t xmin=0, Double_t xmax=0)
 Fit histogram with function fname.
 
virtual TFitResultPtr Fit (TF1 *f1, Option_t *option="", Option_t *goption="", Double_t xmin=0, Double_t xmax=0)
 Fit histogram with the function pointer f1.
 
virtual void FitPanel ()
 Display a panel with all histogram fit options.
 
TH1GetAsymmetry (TH1 *h2, Double_t c2=1, Double_t dc2=0)
 Return a histogram containing the asymmetry of this histogram with h2, where the asymmetry is defined as:
 
virtual Color_t GetAxisColor (Option_t *axis="X") const
 Return the number of divisions for "axis".
 
virtual Float_t GetBarOffset () const
 
virtual Float_t GetBarWidth () const
 
virtual Double_t GetBinCenter (Int_t bin) const
 Return bin center for 1D histogram.
 
virtual EBinErrorOpt GetBinErrorOption () const
 
virtual Double_t GetBinLowEdge (Int_t bin) const
 Return bin lower edge for 1D histogram.
 
virtual Double_t GetBinWidth (Int_t bin) const
 Return bin width for 1D histogram.
 
virtual Double_t GetBinWithContent (Double_t c, Int_t &binx, Int_t firstx=0, Int_t lastx=0, Double_t maxdiff=0) const
 Compute first binx in the range [firstx,lastx] for which diff = abs(bin_content-c) <= maxdiff.
 
virtual void GetBinXYZ (Int_t binglobal, Int_t &binx, Int_t &biny, Int_t &binz) const
 Return binx, biny, binz corresponding to the global bin number globalbin see TH1::GetBin function above.
 
const Double_tGetBuffer () const
 
Int_t GetBufferLength () const
 
Int_t GetBufferSize () const
 
virtual Double_t GetCellContent (Int_t binx, Int_t biny) const
 
virtual Double_t GetCellError (Int_t binx, Int_t biny) const
 
virtual void GetCenter (Double_t *center) const
 Fill array with center of bins for 1D histogram Better to use h1.GetXaxis()->GetCenter(center)
 
virtual Int_t GetContour (Double_t *levels=nullptr)
 Return contour values into array levels if pointer levels is non zero.
 
virtual Double_t GetContourLevel (Int_t level) const
 Return value of contour number level.
 
virtual Double_t GetContourLevelPad (Int_t level) const
 Return the value of contour number "level" in Pad coordinates.
 
TH1GetCumulative (Bool_t forward=kTRUE, const char *suffix="_cumulative") const
 Return a pointer to a histogram containing the cumulative content.
 
virtual Int_t GetDimension () const
 
TDirectoryGetDirectory () const
 
virtual Double_t GetEffectiveEntries () const
 Number of effective entries of the histogram.
 
virtual Double_t GetEntries () const
 Return the current number of entries.
 
virtual TF1GetFunction (const char *name) const
 Return pointer to function with name.
 
virtual Double_tGetIntegral ()
 Return a pointer to the array of bins integral.
 
virtual Double_t GetKurtosis (Int_t axis=1) const
 
virtual Color_t GetLabelColor (Option_t *axis="X") const
 Return the "axis" label color.
 
virtual Style_t GetLabelFont (Option_t *axis="X") const
 Return the "axis" label font.
 
virtual Float_t GetLabelOffset (Option_t *axis="X") const
 Return the "axis" label offset.
 
virtual Float_t GetLabelSize (Option_t *axis="X") const
 Return the "axis" label size.
 
TListGetListOfFunctions () const
 
virtual void GetLowEdge (Double_t *edge) const
 Fill array with low edge of bins for 1D histogram Better to use h1.GetXaxis()->GetLowEdge(edge)
 
virtual Double_t GetMaximum (Double_t maxval=FLT_MAX) const
 Return maximum value smaller than maxval of bins in the range, unless the value has been overridden by TH1::SetMaximum, in which case it returns that value.
 
virtual Int_t GetMaximumBin () const
 Return location of bin with maximum value in the range.
 
virtual Int_t GetMaximumBin (Int_t &locmax, Int_t &locmay, Int_t &locmaz) const
 Return location of bin with maximum value in the range.
 
virtual Double_t GetMaximumStored () const
 
virtual Double_t GetMean (Int_t axis=1) const
 For axis = 1,2 or 3 returns the mean value of the histogram along X,Y or Z axis.
 
virtual Double_t GetMeanError (Int_t axis=1) const
 Return standard error of mean of this histogram along the X axis.
 
virtual Double_t GetMinimum (Double_t minval=-FLT_MAX) const
 Return minimum value larger than minval of bins in the range, unless the value has been overridden by TH1::SetMinimum, in which case it returns that value.
 
virtual void GetMinimumAndMaximum (Double_t &min, Double_t &max) const
 Retrieve the minimum and maximum values in the histogram.
 
virtual Int_t GetMinimumBin () const
 Return location of bin with minimum value in the range.
 
virtual Int_t GetMinimumBin (Int_t &locmix, Int_t &locmiy, Int_t &locmiz) const
 Return location of bin with minimum value in the range.
 
virtual Double_t GetMinimumStored () const
 
virtual Int_t GetNbinsX () const
 
virtual Int_t GetNbinsY () const
 
virtual Int_t GetNbinsZ () const
 
virtual Int_t GetNcells () const
 
virtual Int_t GetNdivisions (Option_t *axis="X") const
 Return the number of divisions for "axis".
 
virtual Double_t GetNormFactor () const
 
char * GetObjectInfo (Int_t px, Int_t py) const override
 Redefines TObject::GetObjectInfo.
 
Option_tGetOption () const override
 
TVirtualHistPainterGetPainter (Option_t *option="")
 Return pointer to painter.
 
virtual Int_t GetQuantiles (Int_t nprobSum, Double_t *q, const Double_t *probSum=nullptr)
 Compute Quantiles for this histogram Quantile x_q of a probability distribution Function F is defined as.
 
virtual Double_t GetRandom (TRandom *rng=nullptr) const
 Return a random number distributed according the histogram bin contents.
 
Double_t GetRMS (Int_t axis=1) const
 This function returns the Standard Deviation (Sigma) of the distribution not the Root Mean Square (RMS).
 
Double_t GetRMSError (Int_t axis=1) const
 
virtual Double_t GetSkewness (Int_t axis=1) const
 
EStatOverflows GetStatOverflows () const
 Get the behaviour adopted by the object about the statoverflows. See EStatOverflows for more information.
 
virtual Double_t GetStdDev (Int_t axis=1) const
 Returns the Standard Deviation (Sigma).
 
virtual Double_t GetStdDevError (Int_t axis=1) const
 Return error of standard deviation estimation for Normal distribution.
 
virtual Double_t GetSumOfWeights () const
 Return the sum of weights excluding under/overflows.
 
virtual TArrayDGetSumw2 ()
 
virtual const TArrayDGetSumw2 () const
 
virtual Int_t GetSumw2N () const
 
virtual Float_t GetTickLength (Option_t *axis="X") const
 Return the "axis" tick length.
 
virtual Style_t GetTitleFont (Option_t *axis="X") const
 Return the "axis" title font.
 
virtual Float_t GetTitleOffset (Option_t *axis="X") const
 Return the "axis" title offset.
 
virtual Float_t GetTitleSize (Option_t *axis="X") const
 Return the "axis" title size.
 
TAxisGetXaxis ()
 
const TAxisGetXaxis () const
 
TAxisGetYaxis ()
 
const TAxisGetYaxis () const
 
TAxisGetZaxis ()
 
const TAxisGetZaxis () const
 
Bool_t IsBinOverflow (Int_t bin, Int_t axis=0) const
 Return true if the bin is overflow.
 
Bool_t IsBinUnderflow (Int_t bin, Int_t axis=0) const
 Return true if the bin is underflow.
 
virtual Bool_t IsHighlight () const
 
Long64_t Merge (TCollection *list, Option_t *option)
 Add all histograms in the collection to this histogram.
 
void Paint (Option_t *option="") override
 Control routine to paint any kind of histograms.
 
void Print (Option_t *option="") const override
 Print some global quantities for this histogram.
 
virtual void RebinAxis (Double_t x, TAxis *axis)
 
virtual void Rebuild (Option_t *option="")
 Using the current bin info, recompute the arrays for contents and errors.
 
void RecursiveRemove (TObject *obj) override
 Recursively remove object from the list of functions.
 
virtual void ResetStats ()
 Reset the statistics including the number of entries and replace with values calculated from bin content.
 
virtual void SetAxisColor (Color_t color=1, Option_t *axis="X")
 Set color to draw the axis line and tick marks.
 
virtual void SetAxisRange (Double_t xmin, Double_t xmax, Option_t *axis="X")
 Set the "axis" range.
 
virtual void SetBarOffset (Float_t offset=0.25)
 Set the bar offset as fraction of the bin width for drawing mode "B".
 
virtual void SetBarWidth (Float_t width=0.5)
 Set the width of bars as fraction of the bin width for drawing mode "B".
 
virtual void SetBinError (Int_t bin, Double_t error)
 Set the bin Error Note that this resets the bin eror option to be of Normal Type and for the non-empty bin the bin error is set by default to the square root of their content.
 
virtual void SetBinError (Int_t binx, Int_t biny, Double_t error)
 See convention for numbering bins in TH1::GetBin.
 
virtual void SetBinError (Int_t binx, Int_t biny, Int_t binz, Double_t error)
 See convention for numbering bins in TH1::GetBin.
 
virtual void SetBinErrorOption (EBinErrorOpt type)
 
virtual UInt_t SetCanExtend (UInt_t extendBitMask)
 Make the histogram axes extendable / not extendable according to the bit mask returns the previous bit mask specifying which axes are extendable.
 
virtual void SetCellContent (Int_t binx, Int_t biny, Double_t content)
 
virtual void SetCellError (Int_t binx, Int_t biny, Double_t content)
 
virtual void SetContent (const Double_t *content)
 Replace bin contents by the contents of array content.
 
virtual void SetContour (Int_t nlevels, const Double_t *levels=nullptr)
 Set the number and values of contour levels.
 
virtual void SetContourLevel (Int_t level, Double_t value)
 Set value for one contour level.
 
virtual void SetDirectory (TDirectory *dir)
 By default, when a histogram is created, it is added to the list of histogram objects in the current directory in memory.
 
virtual void SetEntries (Double_t n)
 
virtual void SetError (const Double_t *error)
 Replace bin errors by values in array error.
 
virtual void SetHighlight (Bool_t set=kTRUE)
 Set highlight (enable/disable) mode for the histogram by default highlight mode is disable.
 
virtual void SetLabelColor (Color_t color=1, Option_t *axis="X")
 Set axis labels color.
 
virtual void SetLabelFont (Style_t font=62, Option_t *axis="X")
 Set font number used to draw axis labels.
 
virtual void SetLabelOffset (Float_t offset=0.005, Option_t *axis="X")
 Set offset between axis and axis' labels.
 
virtual void SetLabelSize (Float_t size=0.02, Option_t *axis="X")
 Set size of axis' labels.
 
virtual void SetMaximum (Double_t maximum=-1111)
 
virtual void SetMinimum (Double_t minimum=-1111)
 
void SetName (const char *name) override
 Change the name of this histogram.
 
void SetNameTitle (const char *name, const char *title) override
 Change the name and title of this histogram.
 
virtual void SetNdivisions (Int_t n=510, Option_t *axis="X")
 Set the number of divisions to draw an axis.
 
virtual void SetNormFactor (Double_t factor=1)
 
virtual void SetOption (Option_t *option=" ")
 
void SetStatOverflows (EStatOverflows statOverflows)
 See GetStatOverflows for more information.
 
virtual void SetStats (Bool_t stats=kTRUE)
 Set statistics option on/off.
 
virtual void SetTickLength (Float_t length=0.02, Option_t *axis="X")
 Set the axis' tick marks length.
 
void SetTitle (const char *title) override
 Change/set the title.
 
virtual void SetTitleFont (Style_t font=62, Option_t *axis="X")
 Set the axis' title font.
 
virtual void SetTitleOffset (Float_t offset=1, Option_t *axis="X")
 Specify a parameter offset to control the distance between the axis and the axis' title.
 
virtual void SetTitleSize (Float_t size=0.02, Option_t *axis="X")
 Set the axis' title size.
 
virtual void SetXTitle (const char *title)
 
virtual void SetYTitle (const char *title)
 
virtual void SetZTitle (const char *title)
 
void StreamerNVirtual (TBuffer &ClassDef_StreamerNVirtual_b)
 
void UseCurrentStyle () override
 Copy current attributes from/to current style.
 
- Public Member Functions inherited from TNamed
 TNamed ()
 
 TNamed (const char *name, const char *title)
 
 TNamed (const TNamed &named)
 TNamed copy ctor.
 
 TNamed (const TString &name, const TString &title)
 
virtual ~TNamed ()
 TNamed destructor.
 
void Clear (Option_t *option="") override
 Set name and title to empty strings ("").
 
TObjectClone (const char *newname="") const override
 Make a clone of an object using the Streamer facility.
 
Int_t Compare (const TObject *obj) const override
 Compare two TNamed objects.
 
virtual void FillBuffer (char *&buffer)
 Encode TNamed into output buffer.
 
const char * GetName () const override
 Returns name of object.
 
const char * GetTitle () const override
 Returns title of object.
 
ULong_t Hash () const override
 Return hash value for this object.
 
Bool_t IsSortable () const override
 
void ls (Option_t *option="") const override
 List TNamed name and title.
 
TNamedoperator= (const TNamed &rhs)
 TNamed assignment operator.
 
void Print (Option_t *option="") const override
 Print TNamed name and title.
 
virtual Int_t Sizeof () const
 Return size of the TNamed part of the TObject.
 
void StreamerNVirtual (TBuffer &ClassDef_StreamerNVirtual_b)
 
- Public Member Functions inherited from TObject
 TObject ()
 TObject constructor.
 
 TObject (const TObject &object)
 TObject copy ctor.
 
virtual ~TObject ()
 TObject destructor.
 
void AbstractMethod (const char *method) const
 Use this method to implement an "abstract" method that you don't want to leave purely abstract.
 
virtual void AppendPad (Option_t *option="")
 Append graphics object to current pad.
 
ULong_t CheckedHash ()
 Check and record whether this class has a consistent Hash/RecursiveRemove setup (*) and then return the regular Hash value for this object.
 
virtual const char * ClassName () const
 Returns name of class to which the object belongs.
 
virtual void Delete (Option_t *option="")
 Delete this object.
 
virtual void DrawClass () const
 Draw class inheritance tree of the class to which this object belongs.
 
virtual TObjectDrawClone (Option_t *option="") const
 Draw a clone of this object in the current selected pad with: gROOT->SetSelectedPad(c1).
 
virtual void Dump () const
 Dump contents of object on stdout.
 
virtual void Error (const char *method, const char *msgfmt,...) const
 Issue error message.
 
virtual void Execute (const char *method, const char *params, Int_t *error=nullptr)
 Execute method on this object with the given parameter string, e.g.
 
virtual void Execute (TMethod *method, TObjArray *params, Int_t *error=nullptr)
 Execute method on this object with parameters stored in the TObjArray.
 
virtual void Fatal (const char *method, const char *msgfmt,...) const
 Issue fatal error message.
 
virtual Option_tGetDrawOption () const
 Get option used by the graphics system to draw this object.
 
virtual const char * GetIconName () const
 Returns mime type name of object.
 
virtual UInt_t GetUniqueID () const
 Return the unique object id.
 
virtual Bool_t HandleTimer (TTimer *timer)
 Execute action in response of a timer timing out.
 
Bool_t HasInconsistentHash () const
 Return true is the type of this object is known to have an inconsistent setup for Hash and RecursiveRemove (i.e.
 
virtual void Info (const char *method, const char *msgfmt,...) const
 Issue info message.
 
virtual Bool_t InheritsFrom (const char *classname) const
 Returns kTRUE if object inherits from class "classname".
 
virtual Bool_t InheritsFrom (const TClass *cl) const
 Returns kTRUE if object inherits from TClass cl.
 
virtual void Inspect () const
 Dump contents of this object in a graphics canvas.
 
void InvertBit (UInt_t f)
 
Bool_t IsDestructed () const
 IsDestructed.
 
virtual Bool_t IsEqual (const TObject *obj) const
 Default equal comparison (objects are equal if they have the same address in memory).
 
virtual Bool_t IsFolder () const
 Returns kTRUE in case object contains browsable objects (like containers or lists of other objects).
 
R__ALWAYS_INLINE Bool_t IsOnHeap () const
 
R__ALWAYS_INLINE Bool_t IsZombie () const
 
void MayNotUse (const char *method) const
 Use this method to signal that a method (defined in a base class) may not be called in a derived class (in principle against good design since a child class should not provide less functionality than its parent, however, sometimes it is necessary).
 
virtual Bool_t Notify ()
 This method must be overridden to handle object notification (the base implementation is no-op).
 
void Obsolete (const char *method, const char *asOfVers, const char *removedFromVers) const
 Use this method to declare a method obsolete.
 
void operator delete (void *ptr)
 Operator delete.
 
void operator delete (void *ptr, void *vp)
 Only called by placement new when throwing an exception.
 
void operator delete[] (void *ptr)
 Operator delete [].
 
void operator delete[] (void *ptr, void *vp)
 Only called by placement new[] when throwing an exception.
 
void * operator new (size_t sz)
 
void * operator new (size_t sz, void *vp)
 
void * operator new[] (size_t sz)
 
void * operator new[] (size_t sz, void *vp)
 
TObjectoperator= (const TObject &rhs)
 TObject assignment operator.
 
virtual void Pop ()
 Pop on object drawn in a pad to the top of the display list.
 
virtual Int_t Read (const char *name)
 Read contents of object with specified name from the current directory.
 
void ResetBit (UInt_t f)
 
virtual void SaveAs (const char *filename="", Option_t *option="") const
 Save this object in the file specified by filename.
 
void SetBit (UInt_t f)
 
void SetBit (UInt_t f, Bool_t set)
 Set or unset the user status bits as specified in f.
 
virtual void SetDrawOption (Option_t *option="")
 Set drawing option for object.
 
virtual void SetUniqueID (UInt_t uid)
 Set the unique object id.
 
void StreamerNVirtual (TBuffer &ClassDef_StreamerNVirtual_b)
 
virtual void SysError (const char *method, const char *msgfmt,...) const
 Issue system error message.
 
R__ALWAYS_INLINE Bool_t TestBit (UInt_t f) const
 
Int_t TestBits (UInt_t f) const
 
virtual void Warning (const char *method, const char *msgfmt,...) const
 Issue warning message.
 
virtual Int_t Write (const char *name=nullptr, Int_t option=0, Int_t bufsize=0)
 Write this object to the current directory.
 
virtual Int_t Write (const char *name=nullptr, Int_t option=0, Int_t bufsize=0) const
 Write this object to the current directory.
 
- Public Member Functions inherited from TAttLine
 TAttLine ()
 AttLine default constructor.
 
 TAttLine (Color_t lcolor, Style_t lstyle, Width_t lwidth)
 AttLine normal constructor.
 
virtual ~TAttLine ()
 AttLine destructor.
 
void Copy (TAttLine &attline) const
 Copy this line attributes to a new TAttLine.
 
Int_t DistancetoLine (Int_t px, Int_t py, Double_t xp1, Double_t yp1, Double_t xp2, Double_t yp2)
 Compute distance from point px,py to a line.
 
virtual Color_t GetLineColor () const
 Return the line color.
 
virtual Style_t GetLineStyle () const
 Return the line style.
 
virtual Width_t GetLineWidth () const
 Return the line width.
 
virtual void Modify ()
 Change current line attributes if necessary.
 
virtual void ResetAttLine (Option_t *option="")
 Reset this line attributes to default values.
 
virtual void SaveLineAttributes (std::ostream &out, const char *name, Int_t coldef=1, Int_t stydef=1, Int_t widdef=1)
 Save line attributes as C++ statement(s) on output stream out.
 
virtual void SetLineAttributes ()
 Invoke the DialogCanvas Line attributes.
 
virtual void SetLineColor (Color_t lcolor)
 Set the line color.
 
virtual void SetLineColorAlpha (Color_t lcolor, Float_t lalpha)
 Set a transparent line color.
 
virtual void SetLineStyle (Style_t lstyle)
 Set the line style.
 
virtual void SetLineWidth (Width_t lwidth)
 Set the line width.
 
void StreamerNVirtual (TBuffer &ClassDef_StreamerNVirtual_b)
 
- Public Member Functions inherited from TAttFill
 TAttFill ()
 AttFill default constructor.
 
 TAttFill (Color_t fcolor, Style_t fstyle)
 AttFill normal constructor.
 
virtual ~TAttFill ()
 AttFill destructor.
 
void Copy (TAttFill &attfill) const
 Copy this fill attributes to a new TAttFill.
 
virtual Color_t GetFillColor () const
 Return the fill area color.
 
virtual Style_t GetFillStyle () const
 Return the fill area style.
 
virtual Bool_t IsTransparent () const
 
virtual void Modify ()
 Change current fill area attributes if necessary.
 
virtual void ResetAttFill (Option_t *option="")
 Reset this fill attributes to default values.
 
virtual void SaveFillAttributes (std::ostream &out, const char *name, Int_t coldef=1, Int_t stydef=1001)
 Save fill attributes as C++ statement(s) on output stream out.
 
virtual void SetFillAttributes ()
 Invoke the DialogCanvas Fill attributes.
 
virtual void SetFillColor (Color_t fcolor)
 Set the fill area color.
 
virtual void SetFillColorAlpha (Color_t fcolor, Float_t falpha)
 Set a transparent fill color.
 
virtual void SetFillStyle (Style_t fstyle)
 Set the fill area style.
 
void StreamerNVirtual (TBuffer &ClassDef_StreamerNVirtual_b)
 
- Public Member Functions inherited from TAttMarker
 TAttMarker ()
 TAttMarker default constructor.
 
 TAttMarker (Color_t color, Style_t style, Size_t msize)
 TAttMarker normal constructor.
 
virtual ~TAttMarker ()
 TAttMarker destructor.
 
void Copy (TAttMarker &attmarker) const
 Copy this marker attributes to a new TAttMarker.
 
virtual Color_t GetMarkerColor () const
 Return the marker color.
 
virtual Size_t GetMarkerSize () const
 Return the marker size.
 
virtual Style_t GetMarkerStyle () const
 Return the marker style.
 
virtual void Modify ()
 Change current marker attributes if necessary.
 
virtual void ResetAttMarker (Option_t *toption="")
 Reset this marker attributes to the default values.
 
virtual void SaveMarkerAttributes (std::ostream &out, const char *name, Int_t coldef=1, Int_t stydef=1, Int_t sizdef=1)
 Save line attributes as C++ statement(s) on output stream out.
 
virtual void SetMarkerAttributes ()
 Invoke the DialogCanvas Marker attributes.
 
virtual void SetMarkerColor (Color_t mcolor=1)
 Set the marker color.
 
virtual void SetMarkerColorAlpha (Color_t mcolor, Float_t malpha)
 Set a transparent marker color.
 
virtual void SetMarkerSize (Size_t msize=1)
 Set the marker size.
 
virtual void SetMarkerStyle (Style_t mstyle=1)
 Set the marker style.
 
void StreamerNVirtual (TBuffer &ClassDef_StreamerNVirtual_b)
 
- Public Member Functions inherited from TArrayD
 TArrayD ()
 Default TArrayD ctor.
 
 TArrayD (const TArrayD &array)
 Copy constructor.
 
 TArrayD (Int_t n)
 Create TArrayD object and set array size to n doubles.
 
 TArrayD (Int_t n, const Double_t *array)
 Create TArrayD object and initialize it with values of array.
 
virtual ~TArrayD ()
 Delete TArrayD object.
 
void AddAt (Double_t c, Int_t i)
 Set the double c value at position i in the array.
 
void Adopt (Int_t n, Double_t *array)
 Adopt array arr into TArrayD, i.e.
 
Double_t At (Int_t i) const
 
void Copy (TArrayD &array) const
 
Double_tGetArray ()
 
const Double_tGetArray () const
 
Double_t GetAt (Int_t i) const override
 
Stat_t GetSum () const
 
TArrayDoperator= (const TArrayD &rhs)
 TArrayD assignment operator.
 
Double_toperator[] (Int_t i)
 
Double_t operator[] (Int_t i) const
 
void Reset ()
 
void Reset (Double_t val)
 
void Set (Int_t n) override
 Set size of this array to n doubles.
 
void Set (Int_t n, const Double_t *array)
 Set size of this array to n doubles and set the contents This function should not be called if the array was declared via Adopt.
 
void SetAt (Double_t v, Int_t i) override
 
void StreamerNVirtual (TBuffer &ClassDef_StreamerNVirtual_b)
 
- Public Member Functions inherited from TArray
 TArray ()
 
 TArray (const TArray &a)
 
 TArray (Int_t n)
 
virtual ~TArray ()
 
Int_t GetSize () const
 
TArrayoperator= (const TArray &rhs)
 
void StreamerNVirtual (TBuffer &ClassDef_StreamerNVirtual_b)
 

Static Public Member Functions

static void Approximate (Bool_t approx=kTRUE)
 Static function, set the fgApproximate flag.
 
static TClassClass ()
 
static const char * Class_Name ()
 
static constexpr Version_t Class_Version ()
 
static const char * DeclFileName ()
 
- Static Public Member Functions inherited from TH2D
static TClassClass ()
 
static const char * Class_Name ()
 
static constexpr Version_t Class_Version ()
 
static const char * DeclFileName ()
 
- Static Public Member Functions inherited from TH2
static TClassClass ()
 
static const char * Class_Name ()
 
static constexpr Version_t Class_Version ()
 
static const char * DeclFileName ()
 
- Static Public Member Functions inherited from TH1
static void AddDirectory (Bool_t add=kTRUE)
 Sets the flag controlling the automatic add of histograms in memory.
 
static Bool_t AddDirectoryStatus ()
 Static function: cannot be inlined on Windows/NT.
 
static TClassClass ()
 
static const char * Class_Name ()
 
static constexpr Version_t Class_Version ()
 
static const char * DeclFileName ()
 
static Int_t FitOptionsMake (Option_t *option, Foption_t &Foption)
 Decode string choptin and fill fitOption structure.
 
static Int_t GetDefaultBufferSize ()
 Static function return the default buffer size for automatic histograms the parameter fgBufferSize may be changed via SetDefaultBufferSize.
 
static Bool_t GetDefaultSumw2 ()
 Return kTRUE if TH1::Sumw2 must be called when creating new histograms.
 
static void SetDefaultBufferSize (Int_t buffersize=1000)
 Static function to set the default buffer size for automatic histograms.
 
static void SetDefaultSumw2 (Bool_t sumw2=kTRUE)
 When this static function is called with sumw2=kTRUE, all new histograms will automatically activate the storage of the sum of squares of errors, ie TH1::Sumw2 is automatically called.
 
static void SmoothArray (Int_t NN, Double_t *XX, Int_t ntimes=1)
 Smooth array xx, translation of Hbook routine hsmoof.F.
 
static void StatOverflows (Bool_t flag=kTRUE)
 if flag=kTRUE, underflows and overflows are used by the Fill functions in the computation of statistics (mean value, StdDev).
 
static TH1TransformHisto (TVirtualFFT *fft, TH1 *h_output, Option_t *option)
 For a given transform (first parameter), fills the histogram (second parameter) with the transform output data, specified in the third parameter If the 2nd parameter h_output is empty, a new histogram (TH1D or TH2D) is created and the user is responsible for deleting it.
 
- Static Public Member Functions inherited from TNamed
static TClassClass ()
 
static const char * Class_Name ()
 
static constexpr Version_t Class_Version ()
 
static const char * DeclFileName ()
 
- Static Public Member Functions inherited from TObject
static TClassClass ()
 
static const char * Class_Name ()
 
static constexpr Version_t Class_Version ()
 
static const char * DeclFileName ()
 
static Longptr_t GetDtorOnly ()
 Return destructor only flag.
 
static Bool_t GetObjectStat ()
 Get status of object stat flag.
 
static void SetDtorOnly (void *obj)
 Set destructor only flag.
 
static void SetObjectStat (Bool_t stat)
 Turn on/off tracking of objects in the TObjectTable.
 
- Static Public Member Functions inherited from TAttLine
static TClassClass ()
 
static const char * Class_Name ()
 
static constexpr Version_t Class_Version ()
 
static const char * DeclFileName ()
 
- Static Public Member Functions inherited from TAttFill
static TClassClass ()
 
static const char * Class_Name ()
 
static constexpr Version_t Class_Version ()
 
static const char * DeclFileName ()
 
- Static Public Member Functions inherited from TAttMarker
static TClassClass ()
 
static const char * Class_Name ()
 
static constexpr Version_t Class_Version ()
 
static const char * DeclFileName ()
 
static Width_t GetMarkerLineWidth (Style_t style)
 Internal helper function that returns the line width of the given marker style (0 = filled marker)
 
static Style_t GetMarkerStyleBase (Style_t style)
 Internal helper function that returns the corresponding marker style with line width 1 for the given style.
 
- Static Public Member Functions inherited from TArrayD
static TClassClass ()
 
static const char * Class_Name ()
 
static constexpr Version_t Class_Version ()
 
static const char * DeclFileName ()
 
- Static Public Member Functions inherited from TArray
static TClassClass ()
 
static const char * Class_Name ()
 
static constexpr Version_t Class_Version ()
 
static const char * DeclFileName ()
 
static TArrayReadArray (TBuffer &b, const TClass *clReq)
 Read TArray object from buffer.
 
static void WriteArray (TBuffer &b, const TArray *a)
 Write TArray object to buffer.
 

Protected Member Functions

virtual Int_t BufferFill (Double_t x, Double_t y, Double_t z, Double_t w)
 Accumulate arguments in buffer.
 
Int_t BufferFill (Double_t, Double_t) override
 accumulate arguments in buffer.
 
Int_t BufferFill (Double_t, Double_t, Double_t) override
 accumulate arguments in buffer.
 
TProfileDoProfile (bool onX, const char *name, Int_t firstbin, Int_t lastbin, Option_t *option) const override
 Implementation of ProfileX or ProfileY for a TProfile2D.
 
Int_t Fill (const char *, Double_t) override
 Increment bin with namex with a weight w.
 
virtual Int_t Fill (const char *namex, const char *namey, Double_t w)
 Increment cell defined by namex,namey by a weight w.
 
virtual Int_t Fill (const char *namex, Double_t y, Double_t w)
 Increment cell defined by namex,y by a weight w.
 
Int_t Fill (const Double_t *v)
 
virtual Int_t Fill (Double_t x, const char *namey, Double_t w)
 Increment cell defined by x,namey by a weight w.
 
Int_t Fill (Double_t x, Double_t y) override
 Increment cell defined by x,y by 1.
 
virtual Int_t Fill (Double_t x, Double_t y, Double_t w)
 Increment cell defined by x,y by a weight w.
 
Int_t Fill (Double_t) override
 Invalid Fill method.
 
Int_t Fill (Double_t, Double_t) override
 Increment cell defined by x,y by 1.
 
Double_t GetBinErrorSqUnchecked (Int_t bin) const override
 
Double_t RetrieveBinContent (Int_t bin) const override
 Raw retrieval of bin content on internal data structure see convention for numbering bins in TH1::GetBin.
 
void SetBins (const Int_t *nbins, const Double_t *range)
 
- Protected Member Functions inherited from TH2D
void UpdateBinContent (Int_t bin, Double_t content) override
 Raw update of bin content on internal data structure see convention for numbering bins in TH1::GetBin.
 
- Protected Member Functions inherited from TH2
 TH2 ()
 2-D histogram default constructor.
 
 TH2 (const char *name, const char *title, Int_t nbinsx, const Double_t *xbins, Int_t nbinsy, const Double_t *ybins)
 Constructor for Double_t variable bin size 2-D histograms.
 
 TH2 (const char *name, const char *title, Int_t nbinsx, const Double_t *xbins, Int_t nbinsy, Double_t ylow, Double_t yup)
 Constructor for variable bin size (along X axis) 2-D histograms using an input array of type double.
 
 TH2 (const char *name, const char *title, Int_t nbinsx, const Float_t *xbins, Int_t nbinsy, const Float_t *ybins)
 Constructor for variable bin size (along X and Y axis) 2-D histograms using input arrays of type float.
 
 TH2 (const char *name, const char *title, Int_t nbinsx, Double_t xlow, Double_t xup, Int_t nbinsy, const Double_t *ybins)
 Constructor for Double_t variable bin size (along Y axis) 2-D histograms.
 
 TH2 (const char *name, const char *title, Int_t nbinsx, Double_t xlow, Double_t xup, Int_t nbinsy, Double_t ylow, Double_t yup)
 Constructor for fix bin size 2-D histograms.
 
virtual void DoFitSlices (bool onX, TF1 *f1, Int_t firstbin, Int_t lastbin, Int_t cut, Option_t *option, TObjArray *arr)
 
virtual TH1DDoProjection (bool onX, const char *name, Int_t firstbin, Int_t lastbin, Option_t *option) const
 Internal (protected) method for performing projection on the X or Y axis called by ProjectionX or ProjectionY.
 
virtual TH1DDoQuantiles (bool onX, const char *name, Double_t prob) const
 Implementation of quantiles for x or y.
 
Double_t Interpolate (Double_t x) const override
 illegal for a TH2
 
- Protected Member Functions inherited from TH1
 TH1 ()
 Histogram default constructor.
 
 TH1 (const char *name, const char *title, Int_t nbinsx, const Double_t *xbins)
 Constructor for variable bin size histograms using an input array of type double.
 
 TH1 (const char *name, const char *title, Int_t nbinsx, const Float_t *xbins)
 Constructor for variable bin size histograms using an input array of type float.
 
 TH1 (const char *name, const char *title, Int_t nbinsx, Double_t xlow, Double_t xup)
 Constructor for fix bin size histograms.
 
virtual Int_t AutoP2FindLimits (Double_t min, Double_t max)
 Buffer-based estimate of the histogram range using the power of 2 algorithm.
 
Int_t AxisChoice (Option_t *axis) const
 Choose an axis according to "axis".
 
virtual void DoFillN (Int_t ntimes, const Double_t *x, const Double_t *w, Int_t stride=1)
 Internal method to fill histogram content from a vector called directly by TH1::BufferEmpty.
 
virtual Double_t DoIntegral (Int_t ix1, Int_t ix2, Int_t iy1, Int_t iy2, Int_t iz1, Int_t iz2, Double_t &err, Option_t *opt, Bool_t doerr=kFALSE) const
 Internal function compute integral and optionally the error between the limits specified by the bin number values working for all histograms (1D, 2D and 3D)
 
virtual Bool_t FindNewAxisLimits (const TAxis *axis, const Double_t point, Double_t &newMin, Double_t &newMax)
 finds new limits for the axis so that point is within the range and the limits are compatible with the previous ones (see TH1::Merge).
 
UInt_t GetAxisLabelStatus () const
 Internal function used in TH1::Fill to see which axis is full alphanumeric, i.e.
 
Bool_t GetStatOverflowsBehaviour () const
 
Bool_t IsEmpty () const
 Check if a histogram is empty (this is a protected method used mainly by TH1Merger )
 
virtual void SavePrimitiveHelp (std::ostream &out, const char *hname, Option_t *option="")
 Helper function for the SavePrimitive functions from TH1 or classes derived from TH1, eg TProfile, TProfile2D.
 
- Protected Member Functions inherited from TObject
virtual void DoError (int level, const char *location, const char *fmt, va_list va) const
 Interface to ErrorHandler (protected).
 
void MakeZombie ()
 
- Protected Member Functions inherited from TArray
Bool_t BoundsOk (const char *where, Int_t at) const
 
Bool_t OutOfBoundsError (const char *where, Int_t i) const
 Generate an out-of-bounds error. Always returns false.
 

Protected Attributes

TArrayD fBinEntries
 Number of entries per bin.
 
TArrayD fBinSumw2
 Array of sum of squares of weights per bin.
 
EErrorType fErrorMode
 Option to compute errors.
 
Bool_t fScaling
 ! True when TProfile2D::Scale is called
 
Double_t fTsumwz
 Total Sum of weight*Z.
 
Double_t fTsumwz2
 Total Sum of weight*Z*Z.
 
Double_t fZmax
 Upper limit in Z (if set)
 
Double_t fZmin
 Lower limit in Z (if set)
 
- Protected Attributes inherited from TH2
Double_t fScalefactor
 Scale factor.
 
Double_t fTsumwxy
 Total Sum of weight*X*Y.
 
Double_t fTsumwy
 Total Sum of weight*Y.
 
Double_t fTsumwy2
 Total Sum of weight*Y*Y.
 
- Protected Attributes inherited from TH1
Short_t fBarOffset
 (1000*offset) for bar charts or legos
 
Short_t fBarWidth
 (1000*width) for bar charts or legos
 
EBinErrorOpt fBinStatErrOpt
 Option for bin statistical errors.
 
Double_tfBuffer
 [fBufferSize] entry buffer
 
Int_t fBufferSize
 fBuffer size
 
TArrayD fContour
 Array to display contour levels.
 
Int_t fDimension
 ! Histogram dimension (1, 2 or 3 dim)
 
TDirectoryfDirectory
 ! Pointer to directory holding this histogram
 
Double_t fEntries
 Number of entries.
 
TListfFunctions
 ->Pointer to list of functions (fits and user)
 
Double_tfIntegral
 ! Integral of bins used by GetRandom
 
Double_t fMaximum
 Maximum value for plotting.
 
Double_t fMinimum
 Minimum value for plotting.
 
Int_t fNcells
 Number of bins(1D), cells (2D) +U/Overflows.
 
Double_t fNormFactor
 Normalization factor.
 
TString fOption
 Histogram options.
 
TVirtualHistPainterfPainter
 ! Pointer to histogram painter
 
EStatOverflows fStatOverflows
 Per object flag to use under/overflows in statistics.
 
TArrayD fSumw2
 Array of sum of squares of weights.
 
Double_t fTsumw
 Total Sum of weights.
 
Double_t fTsumw2
 Total Sum of squares of weights.
 
Double_t fTsumwx
 Total Sum of weight*X.
 
Double_t fTsumwx2
 Total Sum of weight*X*X.
 
TAxis fXaxis
 X axis descriptor.
 
TAxis fYaxis
 Y axis descriptor.
 
TAxis fZaxis
 Z axis descriptor.
 
- Protected Attributes inherited from TNamed
TString fName
 
TString fTitle
 
- Protected Attributes inherited from TAttLine
Color_t fLineColor
 Line color.
 
Style_t fLineStyle
 Line style.
 
Width_t fLineWidth
 Line width.
 
- Protected Attributes inherited from TAttFill
Color_t fFillColor
 Fill area color.
 
Style_t fFillStyle
 Fill area style.
 
- Protected Attributes inherited from TAttMarker
Color_t fMarkerColor
 Marker color.
 
Size_t fMarkerSize
 Marker size.
 
Style_t fMarkerStyle
 Marker style.
 

Static Protected Attributes

static Bool_t fgApproximate = kFALSE
 Bin error approximation option.
 
- Static Protected Attributes inherited from TH1
static Bool_t fgAddDirectory = kTRUE
 ! Flag to add histograms to the directory
 
static Int_t fgBufferSize = 1000
 ! Default buffer size for automatic histograms
 
static Bool_t fgDefaultSumw2 = kFALSE
 ! Flag to call TH1::Sumw2 automatically at histogram creation time
 
static Bool_t fgStatOverflows = kFALSE
 ! Flag to use under/overflows in statistics
 

Private Member Functions

Double_tGetB ()
 
Double_tGetB2 ()
 
Double_tGetW ()
 
Double_tGetW2 ()
 
void SetBins (Int_t, const Double_t *) override
 Redefine x axis parameters with variable bin sizes.
 
void SetBins (Int_t, const Double_t *, Int_t, const Double_t *, Int_t, const Double_t *) override
 Redefine x, y and z axis parameters with variable bin sizes.
 
void SetBins (Int_t, Double_t, Double_t) override
 Redefine x axis parameters.
 
void SetBins (Int_t, Double_t, Double_t, Int_t, Double_t, Double_t, Int_t, Double_t, Double_t) override
 Redefine x, y and z axis parameters.
 

Friends

class TH1Merger
 
class TProfileHelper
 

Additional Inherited Members

- Public Types inherited from TH1
enum  {
  kNoAxis = 0 , kXaxis = (1ULL << ( 0 )) , kYaxis = (1ULL << ( 1 )) , kZaxis = (1ULL << ( 2 )) ,
  kAllAxes = kXaxis | kYaxis | kZaxis
}
 Enumeration specifying which axes can be extended. More...
 
enum  { kNstat = 13 }
 Size of statistics data (size of array used in GetStats()/ PutStats ) More...
 
enum  EBinErrorOpt { kNormal = 0 , kPoisson = 1 , kPoisson2 = 2 }
 Enumeration specifying type of statistics for bin errors. More...
 
enum  EStatOverflows { kIgnore = 0 , kConsider = 1 , kNeutral = 2 }
 Enumeration specifying the way to treat statoverflow. More...
 
enum  EStatusBits {
  kNoStats = (1ULL << ( 9 )) , kUserContour = (1ULL << ( 10 )) , kLogX = (1ULL << ( 15 )) , kIsZoomed = (1ULL << ( 16 )) ,
  kNoTitle = (1ULL << ( 17 )) , kIsAverage = (1ULL << ( 18 )) , kIsNotW = (1ULL << ( 19 )) , kAutoBinPTwo = (1ULL << ( 20 )) ,
  kIsHighlight = (1ULL << ( 21 ))
}
 TH1 status bits. More...
 
- Public Types inherited from TObject
enum  {
  kIsOnHeap = 0x01000000 , kNotDeleted = 0x02000000 , kZombie = 0x04000000 , kInconsistent = 0x08000000 ,
  kBitMask = 0x00ffffff
}
 
enum  { kSingleKey = (1ULL << ( 0 )) , kOverwrite = (1ULL << ( 1 )) , kWriteDelete = (1ULL << ( 2 )) }
 
enum  EDeprecatedStatusBits { kObjInCanvas = (1ULL << ( 3 )) }
 
enum  EStatusBits {
  kCanDelete = (1ULL << ( 0 )) , kMustCleanup = (1ULL << ( 3 )) , kIsReferenced = (1ULL << ( 4 )) , kHasUUID = (1ULL << ( 5 )) ,
  kCannotPick = (1ULL << ( 6 )) , kNoContextMenu = (1ULL << ( 8 )) , kInvalidObject = (1ULL << ( 13 ))
}
 
- Public Attributes inherited from TArrayD
Double_tfArray
 
- Public Attributes inherited from TArray
Int_t fN
 
- Protected Types inherited from TObject
enum  { kOnlyPrepStep = (1ULL << ( 3 )) }
 
- Static Protected Member Functions inherited from TH1
static Int_t AutoP2GetBins (Int_t n)
 Auxiliary function to get the next power of 2 integer value larger then n.
 
static Double_t AutoP2GetPower2 (Double_t x, Bool_t next=kTRUE)
 Auxiliary function to get the power of 2 next (larger) or previous (smaller) a given x.
 
static bool CheckAxisLimits (const TAxis *a1, const TAxis *a2)
 Check that the axis limits of the histograms are the same.
 
static bool CheckBinLabels (const TAxis *a1, const TAxis *a2)
 Check that axis have same labels.
 
static bool CheckBinLimits (const TAxis *a1, const TAxis *a2)
 Check bin limits.
 
static bool CheckConsistency (const TH1 *h1, const TH1 *h2)
 Check histogram compatibility.
 
static bool CheckConsistentSubAxes (const TAxis *a1, Int_t firstBin1, Int_t lastBin1, const TAxis *a2, Int_t firstBin2=0, Int_t lastBin2=0)
 Check that two sub axis are the same.
 
static bool CheckEqualAxes (const TAxis *a1, const TAxis *a2)
 Check that the axis are the same.
 
static Bool_t RecomputeAxisLimits (TAxis &destAxis, const TAxis &anAxis)
 Finds new limits for the axis for the Merge function.
 
static Bool_t SameLimitsAndNBins (const TAxis &axis1, const TAxis &axis2)
 Same limits and bins.
 

#include <TProfile2D.h>

Inheritance diagram for TProfile2D:
[legend]

Constructor & Destructor Documentation

◆ TProfile2D() [1/7]

TProfile2D::TProfile2D ( )

Default constructor for Profile2D histograms.

Definition at line 88 of file TProfile2D.cxx.

◆ TProfile2D() [2/7]

TProfile2D::TProfile2D ( const char *  name,
const char *  title,
Int_t  nx,
Double_t  xlow,
Double_t  xup,
Int_t  ny,
Double_t  ylow,
Double_t  yup,
Double_t  zlow,
Double_t  zup,
Option_t option = "" 
)

Constructor for Profile2D histograms with range in z.

The first eight parameters are similar to TH2D::TH2D. Only the values of Z between ZMIN and ZMAX will be considered at filling time. zmin and zmax will also be the maximum and minimum values on the z scale when drawing the profile2D.

See TProfile2D::BuildOptions for more explanations on errors

Definition at line 169 of file TProfile2D.cxx.

◆ TProfile2D() [3/7]

TProfile2D::TProfile2D ( const char *  name,
const char *  title,
Int_t  nx,
Double_t  xlow,
Double_t  xup,
Int_t  ny,
Double_t  ylow,
Double_t  yup,
Option_t option = "" 
)

Normal Constructor for Profile histograms.

The first eight parameters are similar to TH2D::TH2D. All values of z are accepted at filling time. To fill a profile2D histogram, one must use TProfile2D::Fill function.

Note that when filling the profile histogram the function Fill checks if the variable z is between fZmin and fZmax. If a minimum or maximum value is set for the Z scale before filling, then all values below zmin or above zmax will be discarded. Setting the minimum or maximum value for the Z scale before filling has the same effect as calling the special TProfile2D constructor below where zmin and zmax are specified.

H(I,J) is printed as the cell contents. The errors computed are s(I,J) if CHOPT='S' (spread option), or e(I,J) if CHOPT=' ' (error on mean).

See TProfile2D::BuildOptions for explanation of parameters

see other constructors below with all possible combinations of fix and variable bin size like in TH2D.

Definition at line 125 of file TProfile2D.cxx.

◆ TProfile2D() [4/7]

TProfile2D::TProfile2D ( const char *  name,
const char *  title,
Int_t  nbinsx,
const Double_t xbins,
Int_t  nbinsy,
Double_t  ylow,
Double_t  yup,
Option_t option = "" 
)

Create a 2-D Profile with variable bins in X and fix bins in Y.

Definition at line 135 of file TProfile2D.cxx.

◆ TProfile2D() [5/7]

TProfile2D::TProfile2D ( const char *  name,
const char *  title,
Int_t  nbinsx,
Double_t  xlow,
Double_t  xup,
Int_t  nbinsy,
const Double_t ybins,
Option_t option = "" 
)

Create a 2-D Profile with fix bins in X and variable bins in Y.

Definition at line 144 of file TProfile2D.cxx.

◆ TProfile2D() [6/7]

TProfile2D::TProfile2D ( const char *  name,
const char *  title,
Int_t  nbinsx,
const Double_t xbins,
Int_t  nbinsy,
const Double_t ybins,
Option_t option = "" 
)

Create a 2-D Profile with variable bins in X and variable bins in Y.

Definition at line 153 of file TProfile2D.cxx.

◆ TProfile2D() [7/7]

TProfile2D::TProfile2D ( const TProfile2D profile)

Copy constructor.

Definition at line 206 of file TProfile2D.cxx.

◆ ~TProfile2D()

TProfile2D::~TProfile2D ( )
override

Default destructor for Profile2D histograms.

Definition at line 98 of file TProfile2D.cxx.

Member Function Documentation

◆ Add() [1/3]

Bool_t TProfile2D::Add ( const TH1 h1,
const TH1 h2,
Double_t  c1 = 1,
Double_t  c2 = 1 
)
overridevirtual

Replace contents of this profile2D by the addition of h1 and h2.

this = c1*h1 + c2*h2

Reimplemented from TH1.

Definition at line 249 of file TProfile2D.cxx.

◆ Add() [2/3]

Bool_t TProfile2D::Add ( const TH1 h1,
Double_t  c1 = 1 
)
overridevirtual

Performs the operation: this = this + c1*h1 .

Reimplemented from TH1.

Definition at line 230 of file TProfile2D.cxx.

◆ Add() [3/3]

Bool_t TProfile2D::Add ( TF1 h1,
Double_t  c1 = 1,
Option_t option = "" 
)
overridevirtual

Performs the operation: this = this + c1*f1 .

Reimplemented from TH1.

Definition at line 221 of file TProfile2D.cxx.

◆ Approximate()

void TProfile2D::Approximate ( Bool_t  approx = kTRUE)
static

Static function, set the fgApproximate flag.

When the flag is true, the function GetBinError will approximate the bin error with the average profile error on all bins in the following situation only

  • the number of bins in the profile2D is less than 10404 (eg 100x100)
  • the bin number of entries is small ( <5)
  • the estimated bin error is extremely small compared to the bin content (see TProfile2D::GetBinError)

Definition at line 277 of file TProfile2D.cxx.

◆ BufferEmpty()

Int_t TProfile2D::BufferEmpty ( Int_t  action = 0)
overridevirtual

Fill histogram with all entries in the buffer.

  • action = -1 histogram is reset and refilled from the buffer (called by THistPainter::Paint)
  • action = 0 histogram is filled from the buffer
  • action = 1 histogram is filled and buffer is deleted The buffer is automatically deleted when the number of entries in the buffer is greater than the number of entries in the histogram

Reimplemented from TH2.

Definition at line 291 of file TProfile2D.cxx.

◆ BufferFill() [1/3]

Int_t TProfile2D::BufferFill ( Double_t  x,
Double_t  y,
Double_t  z,
Double_t  w 
)
protectedvirtual

Accumulate arguments in buffer.

When buffer is full, empty the buffer.

  • fBuffer[0] = number of entries in buffer
  • fBuffer[1] = w of first entry
  • fBuffer[2] = x of first entry
  • fBuffer[3] = y of first entry
  • fBuffer[4] = z of first entry

Definition at line 358 of file TProfile2D.cxx.

◆ BufferFill() [2/3]

Int_t TProfile2D::BufferFill ( Double_t  x,
Double_t  w 
)
inlineoverrideprotectedvirtual

accumulate arguments in buffer.

When buffer is full, empty the buffer

  • fBuffer[0] = number of entries in buffer
  • fBuffer[1] = w of first entry
  • fBuffer[2] = x of first entry

Reimplemented from TH2.

Definition at line 44 of file TProfile2D.h.

◆ BufferFill() [3/3]

Int_t TProfile2D::BufferFill ( Double_t  x,
Double_t  y,
Double_t  w 
)
inlineoverrideprotectedvirtual

accumulate arguments in buffer.

When buffer is full, empty the buffer

fBuffer[0] = number of entries in buffer
fBuffer[1] = w of first entry
fBuffer[2] = x of first entry
fBuffer[3] = y of first entry
Double_t * fBuffer
[fBufferSize] entry buffer
Definition TH1.h:107
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17
Definition first.py:1

Reimplemented from TH2.

Definition at line 45 of file TProfile2D.h.

◆ BuildOptions()

void TProfile2D::BuildOptions ( Double_t  zmin,
Double_t  zmax,
Option_t option 
)

Set Profile2D histogram structure and options.

  • zmin: minimum value allowed for z
  • zmax: maximum value allowed for z if (zmin = zmax = 0) there are no limits on the allowed z values (zmin = -inf, zmax = +inf)
  • option: this is the option for the computation of the t error of the profile ( TProfile2D::GetBinError ) possible values for the options are documented in TProfile2D::SetErrorOption

    See TProfile::BuildOptions for a detailed description

Definition at line 189 of file TProfile2D.cxx.

◆ Class()

static TClass * TProfile2D::Class ( )
static
Returns
TClass describing this class

◆ Class_Name()

static const char * TProfile2D::Class_Name ( )
static
Returns
Name of this class

◆ Class_Version()

static constexpr Version_t TProfile2D::Class_Version ( )
inlinestaticconstexpr
Returns
Version of this class

Definition at line 148 of file TProfile2D.h.

◆ Copy()

void TProfile2D::Copy ( TObject hnew) const
overridevirtual

Copy a Profile2D histogram to a new profile2D histogram.

Reimplemented from TH2D.

Definition at line 386 of file TProfile2D.cxx.

◆ DeclFileName()

static const char * TProfile2D::DeclFileName ( )
inlinestatic
Returns
Name of the file containing the class declaration

Definition at line 148 of file TProfile2D.h.

◆ Divide() [1/3]

Bool_t TProfile2D::Divide ( const TH1 h1)
overridevirtual

Divide this profile2D by h1.

this = this/h1

This function return kFALSE if the divide operation failed

Reimplemented from TH1.

Definition at line 428 of file TProfile2D.cxx.

◆ Divide() [2/3]

Bool_t TProfile2D::Divide ( const TH1 h1,
const TH1 h2,
Double_t  c1 = 1,
Double_t  c2 = 1,
Option_t option = "" 
)
overridevirtual

Replace contents of this profile2D by the division of h1 by h2.

this = c1*h1/(c2*h2)

This function return kFALSE if the divide operation failed

Reimplemented from TH1.

Definition at line 511 of file TProfile2D.cxx.

◆ Divide() [3/3]

Bool_t TProfile2D::Divide ( TF1 h1,
Double_t  c1 = 1 
)
overridevirtual

Performs the operation: this = this/(c1*f1) .

This function is not implemented

Reimplemented from TH1.

Definition at line 415 of file TProfile2D.cxx.

◆ DoProfile()

TProfile * TProfile2D::DoProfile ( bool  onX,
const char *  name,
Int_t  firstbin,
Int_t  lastbin,
Option_t option 
) const
overrideprotectedvirtual

Implementation of ProfileX or ProfileY for a TProfile2D.

Do correctly the combination of the bin averages when doing the projection

Reimplemented from TH2.

Definition at line 1397 of file TProfile2D.cxx.

◆ ExtendAxis()

void TProfile2D::ExtendAxis ( Double_t  x,
TAxis axis 
)
overridevirtual

Profile histogram is resized along axis such that x is in the axis range.

The new axis limits are recomputed by doubling iteratively the current axis range until the specified value x is within the limits. The algorithm makes a copy of the histogram, then loops on all bins of the old histogram to fill the extended histogram. Takes into account errors (Sumw2) if any. The axis must be extendable before invoking this function.

Ex: h->GetXaxis()->SetCanExtend(kTRUE)

Reimplemented from TH1.

Definition at line 1512 of file TProfile2D.cxx.

◆ Fill() [1/14]

Int_t TH2::Fill ( const char *  namex,
Double_t  w 
)
inlineoverrideprotectedvirtual

Increment bin with namex with a weight w.

if x is less than the low-edge of the first bin, the Underflow bin is incremented if x is equal to or greater than the upper edge of last bin, the Overflow bin is incremented

If the weight is not equal to 1, the storage of the sum of squares of weights is automatically triggered and the sum of the squares of weights is incremented by \( w^2 \) in the bin corresponding to x.

The function returns the corresponding bin number which has its content incremented by w.

Reimplemented from TH2.

Definition at line 58 of file TH2.h.

◆ Fill() [2/14]

Int_t TH2::Fill ( const char *  namex,
const char *  namey,
Double_t  w 
)
protectedvirtual

Increment cell defined by namex,namey by a weight w.

  • if x or/and y is less than the low-edge of the corresponding axis first bin, the Underflow cell is incremented.
  • if x or/and y is equal to or greater than the upper edge of corresponding axis last bin, the Overflow cell is incremented.
    • If the weight is not equal to 1, the storage of the sum of squares of weights is automatically triggered and the sum of the squares of weights is incremented by w^2 in the bin corresponding to namex,namey

The function returns the corresponding global bin number which has its content incremented by w

Reimplemented from TH2.

Definition at line 79 of file TH2.cxx.

◆ Fill() [3/14]

Int_t TProfile2D::Fill ( const char *  namex,
const char *  namey,
Double_t  z,
Double_t  w = 1. 
)
virtual

Fill a Profile2D histogram (no weights).

Definition at line 736 of file TProfile2D.cxx.

◆ Fill() [4/14]

Int_t TH2::Fill ( const char *  namex,
Double_t  y,
Double_t  w 
)
protectedvirtual

Increment cell defined by namex,y by a weight w.

  • if x or/and y is less than the low-edge of the corresponding axis first bin, the Underflow cell is incremented.
  • if x or/and y is equal to or greater than the upper edge of corresponding axis last bin, the Overflow cell is incremented.
  • If the weight is not equal to 1, the storage of the sum of squares of weights is automatically triggered and the sum of the squares of weights is incremented by w^2 in the bin corresponding to namex,y

The function returns the corresponding global bin number which has its content incremented by w

Reimplemented from TH2.

Definition at line 78 of file TH2.cxx.

◆ Fill() [5/14]

Int_t TProfile2D::Fill ( const char *  namex,
Double_t  y,
Double_t  z,
Double_t  w = 1. 
)
virtual

Fill a Profile2D histogram (no weights).

Definition at line 781 of file TProfile2D.cxx.

◆ Fill() [6/14]

Int_t TProfile2D::Fill ( const Double_t v)
inlineprotected

Definition at line 51 of file TProfile2D.h.

◆ Fill() [7/14]

Int_t TH2::Fill ( Double_t  x,
const char *  namey,
Double_t  w 
)
protectedvirtual

Increment cell defined by x,namey by a weight w.

  • if x or/and y is less than the low-edge of the corresponding axis first bin, the Underflow cell is incremented.
  • if x or/and y is equal to or greater than the upper edge of corresponding axis last bin, the Overflow cell is incremented.
  • If the weight is not equal to 1, the storage of the sum of squares of weights is automatically triggered and the sum of the squares of weights is incremented by w^2 in the bin corresponding to x,y.

The function returns the corresponding global bin number which has its content incremented by w

Reimplemented from TH2.

Definition at line 77 of file TH2.cxx.

◆ Fill() [8/14]

Int_t TProfile2D::Fill ( Double_t  x,
const char *  namey,
Double_t  z,
Double_t  w = 1. 
)
virtual

Fill a Profile2D histogram (no weights).

Definition at line 691 of file TProfile2D.cxx.

◆ Fill() [9/14]

Int_t TH2::Fill ( Double_t  x,
Double_t  y 
)
overrideprotectedvirtual

Increment cell defined by x,y by 1.

  • if x or/and y is less than the low-edge of the corresponding axis first bin, the Underflow cell is incremented.
  • if x or/and y is equal to or greater than the upper edge of corresponding axis last bin, the Overflow cell is incremented.
  • If the storage of the sum of squares of weights has been triggered, via the function Sumw2, then the sum of the squares of weights is incremented by 1 in the cell corresponding to x,y.

The function returns the corresponding global bin number which has its content incremented by 1

Reimplemented from TH2.

Definition at line 75 of file TH2.cxx.

◆ Fill() [10/14]

Int_t TH2::Fill ( Double_t  x,
Double_t  y,
Double_t  w 
)
protectedvirtual

Increment cell defined by x,y by a weight w.

  • if x or/and y is less than the low-edge of the corresponding axis first bin, the Underflow cell is incremented.
  • if x or/and y is equal to or greater than the upper edge of corresponding axis last bin, the Overflow cell is incremented.
    • If the weight is not equal to 1, the storage of the sum of squares of weights is automatically triggered and the sum of the squares of weights is incremented by w^2 in the bin corresponding to x,y

The function returns the corresponding global bin number which has its content incremented by w

Reimplemented from TH2.

Definition at line 76 of file TH2.cxx.

◆ Fill() [11/14]

Int_t TProfile2D::Fill ( Double_t  x,
Double_t  y,
Double_t  z 
)
overridevirtual

Fill a Profile2D histogram (no weights).

Reimplemented from TH2.

Definition at line 608 of file TProfile2D.cxx.

◆ Fill() [12/14]

Int_t TProfile2D::Fill ( Double_t  x,
Double_t  y,
Double_t  z,
Double_t  w 
)
virtual

Fill a Profile2D histogram with weights.

Definition at line 648 of file TProfile2D.cxx.

◆ Fill() [13/14]

Int_t TH2::Fill ( Double_t  )
overrideprotectedvirtual

Invalid Fill method.

Reimplemented from TH2.

Definition at line 57 of file TH2.cxx.

◆ Fill() [14/14]

Int_t TProfile2D::Fill ( Double_t  x,
Double_t  y 
)
inlineoverrideprotectedvirtual

Increment cell defined by x,y by 1.

  • if x or/and y is less than the low-edge of the corresponding axis first bin, the Underflow cell is incremented.
  • if x or/and y is equal to or greater than the upper edge of corresponding axis last bin, the Overflow cell is incremented.
  • If the storage of the sum of squares of weights has been triggered, via the function Sumw2, then the sum of the squares of weights is incremented by 1 in the cell corresponding to x,y.

The function returns the corresponding global bin number which has its content incremented by 1

Reimplemented from TH2.

Definition at line 56 of file TProfile2D.h.

◆ GetB()

Double_t * TProfile2D::GetB ( )
inlineprivate

Definition at line 63 of file TProfile2D.h.

◆ GetB2()

Double_t * TProfile2D::GetB2 ( )
inlineprivate

Definition at line 64 of file TProfile2D.h.

◆ GetBinContent() [1/3]

Double_t TProfile2D::GetBinContent ( Int_t  bin) const
overridevirtual

Return bin content of a Profile2D histogram.

Reimplemented from TH2.

Definition at line 829 of file TProfile2D.cxx.

◆ GetBinContent() [2/3]

Double_t TProfile2D::GetBinContent ( Int_t  binx,
Int_t  biny 
) const
inlineoverridevirtual

Reimplemented from TH2.

Definition at line 109 of file TProfile2D.h.

◆ GetBinContent() [3/3]

Double_t TProfile2D::GetBinContent ( Int_t  binx,
Int_t  biny,
Int_t   
) const
inlineoverridevirtual

Reimplemented from TH2.

Definition at line 110 of file TProfile2D.h.

◆ GetBinEffectiveEntries()

Double_t TProfile2D::GetBinEffectiveEntries ( Int_t  bin)
virtual

Return bin effective entries for a weighted filled Profile histogram.

In case of an unweighted profile, it is equivalent to the number of entries per bin The effective entries is defined as the square of the sum of the weights divided by the sum of the weights square. TProfile::Sumw2() must be called before filling the profile with weights. Only by calling this method the sum of the square of the weights per bin is stored.

Definition at line 858 of file TProfile2D.cxx.

◆ GetBinEntries()

Double_t TProfile2D::GetBinEntries ( Int_t  bin) const
virtual

Return bin entries of a Profile2D histogram.

Definition at line 842 of file TProfile2D.cxx.

◆ GetBinError() [1/3]

Double_t TProfile2D::GetBinError ( Int_t  bin) const
overridevirtual

Return bin error of a Profile2D histogram.

Computing errors: A moving field

The computation of errors for a TProfile2D has evolved with the versions of ROOT. The difficulty is in computing errors for bins with low statistics.

  • prior to version 3.10, we had no special treatment of low statistic bins. As a result, these bins had huge errors. The reason is that the expression eprim2 is very close to 0 (rounding problems) or 0.
  • The algorithm is modified/protected for the case when a TProfile2D is projected (ProjectionX). The previous algorithm generated a N^2 problem when projecting a TProfile2D with a large number of bins (eg 100000).
  • in version 3.10/02, a new static function TProfile::Approximate is introduced to enable or disable (default) the approximation. (see also comments in TProfile::GetBinError)

Reimplemented from TH1.

Definition at line 881 of file TProfile2D.cxx.

◆ GetBinError() [2/3]

Double_t TProfile2D::GetBinError ( Int_t  binx,
Int_t  biny 
) const
inlineoverridevirtual

Reimplemented from TH1.

Definition at line 112 of file TProfile2D.h.

◆ GetBinError() [3/3]

Double_t TProfile2D::GetBinError ( Int_t  binx,
Int_t  biny,
Int_t   
) const
inlineoverridevirtual

Reimplemented from TH1.

Definition at line 113 of file TProfile2D.h.

◆ GetBinErrorSqUnchecked()

Double_t TProfile2D::GetBinErrorSqUnchecked ( Int_t  bin) const
inlineoverrideprotectedvirtual

Reimplemented from TH1.

Definition at line 60 of file TProfile2D.h.

◆ GetBinSumw2() [1/2]

virtual TArrayD * TProfile2D::GetBinSumw2 ( )
inlinevirtual

Definition at line 116 of file TProfile2D.h.

◆ GetBinSumw2() [2/2]

virtual const TArrayD * TProfile2D::GetBinSumw2 ( ) const
inlinevirtual

Definition at line 117 of file TProfile2D.h.

◆ GetErrorOption()

Option_t * TProfile2D::GetErrorOption ( ) const

Return option to compute profile2D errors.

Definition at line 889 of file TProfile2D.cxx.

◆ GetNumberOfBins()

Double_t TProfile2D::GetNumberOfBins ( )
inline

Definition at line 146 of file TProfile2D.h.

◆ GetStats()

void TProfile2D::GetStats ( Double_t stats) const
overridevirtual

Fill the array stats from the contents of this profile.

The array stats must be correctly dimensioned in the calling program.

  • stats[0] = sumw
  • stats[1] = sumw2
  • stats[2] = sumwx
  • stats[3] = sumwx2
  • stats[4] = sumwy
  • stats[5] = sumwy2
  • stats[6] = sumwxy
  • stats[7] = sumwz
  • stats[8] = sumwz2

If no axis-subrange is specified (via TAxis::SetRange), the array stats is simply a copy of the statistics quantities computed at filling time. If a sub-range is specified, the function recomputes these quantities from the bin contents in the current axis range.

Reimplemented from TH2.

Definition at line 916 of file TProfile2D.cxx.

◆ GetW()

Double_t * TProfile2D::GetW ( )
inlineprivate

Definition at line 65 of file TProfile2D.h.

◆ GetW2()

Double_t * TProfile2D::GetW2 ( )
inlineprivate

Definition at line 66 of file TProfile2D.h.

◆ GetZmax()

virtual Double_t TProfile2D::GetZmax ( ) const
inlinevirtual

Definition at line 121 of file TProfile2D.h.

◆ GetZmin()

virtual Double_t TProfile2D::GetZmin ( ) const
inlinevirtual

Definition at line 120 of file TProfile2D.h.

◆ IsA()

TClass * TProfile2D::IsA ( ) const
inlineoverridevirtual
Returns
TClass describing current object

Reimplemented from TH2D.

Definition at line 148 of file TProfile2D.h.

◆ LabelsDeflate()

void TProfile2D::LabelsDeflate ( Option_t axis = "X")
overridevirtual

Reduce the number of bins for this axis to the number of bins having a label.

Reimplemented from TH1.

Definition at line 980 of file TProfile2D.cxx.

◆ LabelsInflate()

void TProfile2D::LabelsInflate ( Option_t ax = "X")
overridevirtual

Double the number of bins for axis.

Refill histogram This function is called by TAxis::FindBin(const char *label)

Reimplemented from TH1.

Definition at line 990 of file TProfile2D.cxx.

◆ LabelsOption()

void TProfile2D::LabelsOption ( Option_t option = "h",
Option_t ax = "X" 
)
overridevirtual

Set option(s) to draw axis with labels.

option might have the following values:

  • "a" sort by alphabetic order
  • ">" sort by decreasing values
  • "<" sort by increasing values
  • "h" draw labels horizontal
  • "v" draw labels vertical
  • "u" draw labels up (end of label right adjusted)
  • "d" draw labels down (start of label left adjusted)

Reimplemented from TH1.

Definition at line 1008 of file TProfile2D.cxx.

◆ Merge()

Long64_t TProfile2D::Merge ( TCollection li)
overridevirtual

Merge all histograms in the collection in this histogram.

This function computes the min/max for the axes, compute a new number of bins, if necessary, add bin contents, errors and statistics. If overflows are present and limits are different the function will fail. The function returns the total number of entries in the result histogram if the merge is successful, -1 otherwise.

IMPORTANT remark. The 2 axis x and y may have different number of bins and different limits, BUT the largest bin width must be a multiple of the smallest bin width and the upper limit must also be a multiple of the bin width.

Reimplemented from TH1.

Definition at line 1227 of file TProfile2D.cxx.

◆ Multiply() [1/3]

Bool_t TProfile2D::Multiply ( const TH1 h1)
overridevirtual

Multiply this profile2D by h1.

this = this*h1

Reimplemented from TH1.

Definition at line 1246 of file TProfile2D.cxx.

◆ Multiply() [2/3]

Bool_t TProfile2D::Multiply ( const TH1 h1,
const TH1 h2,
Double_t  c1 = 1,
Double_t  c2 = 1,
Option_t option = "" 
)
overridevirtual

Replace contents of this profile2D by multiplication of h1 by h2.

this = (c1*h1)*(c2*h2)

Reimplemented from TH1.

Definition at line 1257 of file TProfile2D.cxx.

◆ Multiply() [3/3]

Bool_t TProfile2D::Multiply ( TF1 h1,
Double_t  c1 = 1 
)
overridevirtual

Performs the operation: this = this*c1*f1.

Reimplemented from TH1.

Definition at line 1235 of file TProfile2D.cxx.

◆ operator=()

TProfile2D & TProfile2D::operator= ( const TProfile2D profile)

Definition at line 211 of file TProfile2D.cxx.

◆ ProfileX()

TProfile * TProfile2D::ProfileX ( const char *  name = "_pfx",
Int_t  firstybin = 0,
Int_t  lastybin = -1,
Option_t option = "" 
) const

Project a 2-D histogram into a profile histogram along X.

The projection is made from the channels along the Y axis ranging from firstybin to lastybin included. The result is a 1D profile which contains the combination of all the considered bins along Y By default, bins 1 to ny are included When all bins are included, the number of entries in the projection is set to the number of entries of the 2-D histogram, otherwise the number of entries is incremented by 1 for all non empty cells.

The option can also be used to specify the projected profile error type. Values which can be used are 's', 'i', or 'g'. See TProfile::BuildOptions for details

Definition at line 1368 of file TProfile2D.cxx.

◆ ProfileY()

TProfile * TProfile2D::ProfileY ( const char *  name = "_pfy",
Int_t  firstxbin = 0,
Int_t  lastxbin = -1,
Option_t option = "" 
) const

Project a 2-D histogram into a profile histogram along X.

The projection is made from the channels along the X axis ranging from firstybin to lastybin included. The result is a 1D profile which contains the combination of all the considered bins along X By default, bins 1 to ny are included When all bins are included, the number of entries in the projection is set to the number of entries of the 2-D histogram, otherwise the number of entries is incremented by 1 for all non empty cells.

The option can also be used to specify the projected profile error type. Values which can be used are 's', 'i', or 'g'. See TProfile::BuildOptions for details

Definition at line 1387 of file TProfile2D.cxx.

◆ ProjectionXY()

TH2D * TProfile2D::ProjectionXY ( const char *  name = "_pxy",
Option_t option = "e" 
) const

Project this profile2D into a 2-D histogram along X,Y.

The projection is always of the type TH2D.

  • if option "E" is specified the errors of the projected histogram are computed and set to be equal to the errors of the profile. Option "E" is defined as the default one in the header file.
  • if option "" is specified the histogram errors are simply the sqrt of its content
  • if option "B" is specified, the content of bin of the returned histogram will be equal to the GetBinEntries(bin) of the profile,
  • if option "C=E" the bin contents of the projection are set to the bin errors of the profile
  • if option "W" is specified the bin content of the projected histogram is set to the product of the bin content of the profile and the entries. With this option the returned histogram will be equivalent to the one obtained by filling directly a TH2D using the 3-rd value as a weight. This option makes sense only for profile filled with all weights =1. When the profile is weighted (filled with weights different than 1) the bin error of the projected histogram (obtained using this option "W") cannot be correctly computed from the information stored in the profile. In that case the obtained histogram contains as bin error square the weighted sum of the square of the profiled observable (TProfile2D::fSumw2[bin] )

Definition at line 1287 of file TProfile2D.cxx.

◆ PutStats()

void TProfile2D::PutStats ( Double_t stats)
overridevirtual

Replace current statistics with the values in array stats.

Reimplemented from TH2.

Definition at line 1472 of file TProfile2D.cxx.

◆ Rebin2D()

TProfile2D * TProfile2D::Rebin2D ( Int_t  nxgroup = 2,
Int_t  nygroup = 2,
const char *  newname = "" 
)
overridevirtual

Rebin this histogram grouping nxgroup/nygroup bins along the xaxis/yaxis together.

if newname is not blank a new profile hnew is created. else the current histogram is modified (default) The parameter nxgroup/nygroup indicate how many bins along the xaxis/yaxis of this have to be merged into one bin of hnew If the original profile has errors stored (via Sumw2), the resulting profile has new errors correctly calculated.

examples: if hpxpy is an existing TProfile2D profile with 40 x 40 bins

hpxpy->Rebin2D(); // merges two bins along the xaxis and yaxis in one
// Carefull: previous contents of hpxpy are lost
hpxpy->Rebin2D(3,5); // merges 3 bins along the xaxis and 5 bins along the yaxis in one
// Carefull: previous contents of hpxpy are lost
hpxpy->RebinX(5); //merges five bins along the xaxis in one in hpxpy
TProfile2D *hnew = hpxpy->RebinY(5,"hnew"); // creates a new profile hnew
// merging 5 bins of hpxpy along the yaxis in one bin
Profile2D histograms are used to display the mean value of Z and its error for each cell in X,...
Definition TProfile2D.h:27
TProfile2D * RebinY(Int_t ngroup=2, const char *newname="") override
Rebin only the Y axis.

NOTE : If nxgroup/nygroup is not an exact divider of the number of bins, along the xaxis/yaxis the top limit(s) of the rebinned profile is changed to the upper edge of the xbin=newxbins*nxgroup resp. ybin=newybins*nygroup and the remaining bins are added to the overflow bin. Statistics will be recomputed from the new bin contents.

Reimplemented from TH2.

Definition at line 1550 of file TProfile2D.cxx.

◆ RebinX()

TProfile2D * TProfile2D::RebinX ( Int_t  ngroup = 2,
const char *  newname = "" 
)
overridevirtual

Rebin only the X axis.

see Rebin2D

Reimplemented from TH2.

Definition at line 1835 of file TProfile2D.cxx.

◆ RebinY()

TProfile2D * TProfile2D::RebinY ( Int_t  ngroup = 2,
const char *  newname = "" 
)
overridevirtual

Rebin only the Y axis.

see Rebin2D

Reimplemented from TH2.

Definition at line 1843 of file TProfile2D.cxx.

◆ Reset()

void TProfile2D::Reset ( Option_t option = "")
overridevirtual

Reset contents of a Profile2D histogram.

Reimplemented from TH2D.

Definition at line 1488 of file TProfile2D.cxx.

◆ RetrieveBinContent()

Double_t TProfile2D::RetrieveBinContent ( Int_t  bin) const
inlineoverrideprotectedvirtual

Raw retrieval of bin content on internal data structure see convention for numbering bins in TH1::GetBin.

Reimplemented from TH2D.

Definition at line 58 of file TProfile2D.h.

◆ SavePrimitive()

void TProfile2D::SavePrimitive ( std::ostream &  out,
Option_t option = "" 
)
overridevirtual

Save primitive as a C++ statement(s) on output stream out.

Note the following restrictions in the code generated:

  • variable bin size not implemented
  • SetErrorOption not implemented

Reimplemented from TH1.

Definition at line 1854 of file TProfile2D.cxx.

◆ Scale()

void TProfile2D::Scale ( Double_t  c1 = 1,
Option_t option = "" 
)
overridevirtual

Multiply this profile2D by a constant c1.

`this = c1*this

This function uses the services of TProfile2D::Add

Reimplemented from TH1.

Definition at line 1908 of file TProfile2D.cxx.

◆ SetBinEntries()

void TProfile2D::SetBinEntries ( Int_t  bin,
Double_t  w 
)
virtual

Set the number of entries in bin.

Definition at line 1916 of file TProfile2D.cxx.

◆ SetBins() [1/7]

void TProfile2D::SetBins ( const Int_t nbins,
const Double_t range 
)
inlineprotected

Definition at line 49 of file TProfile2D.h.

◆ SetBins() [2/7]

void TProfile2D::SetBins ( Int_t  nbinsx,
Double_t  xmin,
Double_t  xmax,
Int_t  nbinsy,
Double_t  ymin,
Double_t  ymax 
)
overridevirtual

Redefine x and y axis parameters.

Reimplemented from TH1.

Definition at line 1924 of file TProfile2D.cxx.

◆ SetBins() [3/7]

void TProfile2D::SetBins ( Int_t  nx,
const Double_t xBins,
Int_t  ny,
const Double_t yBins 
)
overridevirtual

Redefine x and y axis parameters for variable bin sizes.

Reimplemented from TH1.

Definition at line 1934 of file TProfile2D.cxx.

◆ SetBins() [4/7]

void TProfile2D::SetBins ( Int_t  nx,
const Double_t xBins 
)
inlineoverrideprivatevirtual

Redefine x axis parameters with variable bin sizes.

The X axis parameters are modified. The bins content array is resized if errors (Sumw2) the errors array is resized The previous bin contents are lost To change only the axis limits, see TAxis::SetRange xBins is supposed to be of length nx+1

Reimplemented from TH1.

Definition at line 69 of file TProfile2D.h.

◆ SetBins() [5/7]

void TProfile2D::SetBins ( Int_t  nx,
const Double_t xBins,
Int_t  ny,
const Double_t yBins,
Int_t  nz,
const Double_t zBins 
)
inlineoverrideprivatevirtual

Redefine x, y and z axis parameters with variable bin sizes.

The X, Y and Z axis parameters are modified. The bins content array is resized if errors (Sumw2) the errors array is resized The previous bin contents are lost To change only the axis limits, see TAxis::SetRange xBins is supposed to be of length nx+1, yBins is supposed to be of length ny+1, zBins is supposed to be of length nz+1

Reimplemented from TH1.

Definition at line 73 of file TProfile2D.h.

◆ SetBins() [6/7]

void TProfile2D::SetBins ( Int_t  nx,
Double_t  xmin,
Double_t  xmax 
)
inlineoverrideprivatevirtual

Redefine x axis parameters.

The X axis parameters are modified. The bins content array is resized if errors (Sumw2) the errors array is resized The previous bin contents are lost To change only the axis limits, see TAxis::SetRange

Reimplemented from TH1.

Definition at line 67 of file TProfile2D.h.

◆ SetBins() [7/7]

void TProfile2D::SetBins ( Int_t  nx,
Double_t  xmin,
Double_t  xmax,
Int_t  ny,
Double_t  ymin,
Double_t  ymax,
Int_t  nz,
Double_t  zmin,
Double_t  zmax 
)
inlineoverrideprivatevirtual

Redefine x, y and z axis parameters.

The X, Y and Z axis parameters are modified. The bins content array is resized if errors (Sumw2) the errors array is resized The previous bin contents are lost To change only the axis limits, see TAxis::SetRange

Reimplemented from TH1.

Definition at line 71 of file TProfile2D.h.

◆ SetBinsLength()

void TProfile2D::SetBinsLength ( Int_t  n = -1)
overridevirtual

Set total number of bins including under/overflow.

Reallocate bin contents array

Reimplemented from TH2D.

Definition at line 1945 of file TProfile2D.cxx.

◆ SetBuffer()

void TProfile2D::SetBuffer ( Int_t  buffersize,
Option_t option = "" 
)
overridevirtual

Set the buffer size in units of 8 bytes (double).

Reimplemented from TH1.

Definition at line 1954 of file TProfile2D.cxx.

◆ SetErrorOption()

void TProfile2D::SetErrorOption ( Option_t option = "")
virtual

Set option to compute profile2D errors.

The computation of the bin errors is based on the parameter option:

  • ' ' (Default) The bin errors are the standard error on the mean of the bin profiled values (Z), i.e. the standard error of the bin contents. Note that if TProfile::Approximate() is called, an approximation is used when the spread in Z is 0 and the number of bin entries is > 0
  • 's' The bin errors are the standard deviations of the Z bin values Note that if TProfile::Approximate() is called, an approximation is used when the spread in Z is 0 and the number of bin entries is > 0
  • 'i' Errors are as in default case (standard errors of the bin contents) The only difference is for the case when the spread in Z is zero. In this case for N > 0 the error is 1./SQRT(12.*N)
  • 'g' Errors are 1./SQRT(W) for W not equal to 0 and 0 for W = 0. W is the sum in the bin of the weights of the profile. This option is for combining measurements z +/- dz, and the profile is filled with values y and weights z = 1/dz**2

See TProfile::BuildOptions for a detailed explanation of all options

Definition at line 1992 of file TProfile2D.cxx.

◆ Streamer()

void TProfile2D::Streamer ( TBuffer R__b)
overridevirtual

Stream an object of class TProfile2D.

Reimplemented from TH2D.

Definition at line 2000 of file TProfile2D.cxx.

◆ StreamerNVirtual()

void TProfile2D::StreamerNVirtual ( TBuffer ClassDef_StreamerNVirtual_b)
inline

Definition at line 148 of file TProfile2D.h.

◆ Sumw2()

void TProfile2D::Sumw2 ( Bool_t  flag = kTRUE)
overridevirtual

Create/Delete structure to store sum of squares of weights per bin.

This is needed to compute the correct statistical quantities of a profile filled with weights

This function is automatically called when the histogram is created if the static function TH1::SetDefaultSumw2 has been called before. If flag is false the structure is deleted

Reimplemented from TH1.

Definition at line 2041 of file TProfile2D.cxx.

Friends And Related Symbol Documentation

◆ TH1Merger

friend class TH1Merger
friend

Definition at line 31 of file TProfile2D.h.

◆ TProfileHelper

friend class TProfileHelper
friend

Definition at line 30 of file TProfile2D.h.

Member Data Documentation

◆ fBinEntries

TArrayD TProfile2D::fBinEntries
protected

Number of entries per bin.

Definition at line 34 of file TProfile2D.h.

◆ fBinSumw2

TArrayD TProfile2D::fBinSumw2
protected

Array of sum of squares of weights per bin.

Definition at line 41 of file TProfile2D.h.

◆ fErrorMode

EErrorType TProfile2D::fErrorMode
protected

Option to compute errors.

Definition at line 35 of file TProfile2D.h.

◆ fgApproximate

Bool_t TProfile2D::fgApproximate = kFALSE
staticprotected

Bin error approximation option.

Definition at line 42 of file TProfile2D.h.

◆ fScaling

Bool_t TProfile2D::fScaling
protected

! True when TProfile2D::Scale is called

Definition at line 38 of file TProfile2D.h.

◆ fTsumwz

Double_t TProfile2D::fTsumwz
protected

Total Sum of weight*Z.

Definition at line 39 of file TProfile2D.h.

◆ fTsumwz2

Double_t TProfile2D::fTsumwz2
protected

Total Sum of weight*Z*Z.

Definition at line 40 of file TProfile2D.h.

◆ fZmax

Double_t TProfile2D::fZmax
protected

Upper limit in Z (if set)

Definition at line 37 of file TProfile2D.h.

◆ fZmin

Double_t TProfile2D::fZmin
protected

Lower limit in Z (if set)

Definition at line 36 of file TProfile2D.h.

Libraries for TProfile2D:

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