Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
THnBase.cxx
Go to the documentation of this file.
1// @(#)root/hist:$Id$
2// Author: Axel Naumann (2011-12-20)
3
4/*************************************************************************
5 * Copyright (C) 1995-2012, Rene Brun and Fons Rademakers. *
6 * All rights reserved. *
7 * *
8 * For the licensing terms see $ROOTSYS/LICENSE. *
9 * For the list of contributors see $ROOTSYS/README/CREDITS. *
10 *************************************************************************/
11
12#include "THnBase.h"
13
14#include "TAxis.h"
15#include "TBrowser.h"
16#include "TError.h"
17#include "TClass.h"
18#include "TF1.h"
19#include "TH1D.h"
20#include "TH2D.h"
21#include "TH3D.h"
22#include "THn.h"
23#include "THnSparse.h"
24#include "TMath.h"
25#include "TRandom.h"
26#include "TVirtualPad.h"
27
28#include "HFitInterface.h"
29#include "Fit/DataRange.h"
30#include "Fit/SparseData.h"
33
34
35/** \class THnBase
36 \ingroup Hist
37Multidimensional histogram base.
38Defines common functionality and interfaces for THn, THnSparse.
39*/
40
42
43////////////////////////////////////////////////////////////////////////////////
44/// Construct a THnBase with "dim" dimensions,
45/// "nbins" holds the number of bins for each dimension;
46/// "xmin" and "xmax" the minimal and maximal value for each dimension.
47/// The arrays "xmin" and "xmax" can be NULL; in that case SetBinEdges()
48/// must be called for each dimension.
49
50THnBase::THnBase(const char* name, const char* title, Int_t dim,
51 const Int_t* nbins, const Double_t* xmin, const Double_t* xmax):
52TNamed(name, title), fNdimensions(dim), fAxes(dim), fBrowsables(dim),
53fEntries(0), fTsumw(0), fTsumw2(-1.), fTsumwx(dim), fTsumwx2(dim),
54fIntegral(0), fIntegralStatus(kNoInt)
55{
56 for (Int_t i = 0; i < fNdimensions; ++i) {
57 TAxis* axis = new TAxis(nbins[i], xmin ? xmin[i] : 0., xmax ? xmax[i] : 1.);
58 axis->SetName(TString::Format("axis%d", i));
59 fAxes.AddAtAndExpand(axis, i);
60 }
61 SetTitle(title);
63}
64
65THnBase::THnBase(const char *name, const char *title, Int_t dim, const Int_t *nbins,
66 const std::vector<std::vector<double>> &xbins)
67 : TNamed(name, title), fNdimensions(dim), fAxes(dim), fBrowsables(dim), fEntries(0), fTsumw(0), fTsumw2(-1.),
68 fTsumwx(dim), fTsumwx2(dim), fIntegral(0), fIntegralStatus(kNoInt)
69{
70 if (Int_t(xbins.size()) != fNdimensions) {
71 Error("THnBase", "Mismatched number of dimensions %d with number of bin edge vectors %zu", fNdimensions,
72 xbins.size());
73 }
74 for (Int_t i = 0; i < fNdimensions; ++i) {
75 if (Int_t(xbins[i].size()) != (nbins[i] + 1)) {
76 Error("THnBase", "Mismatched number of bins %d with number of bin edges %zu", nbins[i], xbins[i].size());
77 }
78 TAxis *axis = new TAxis(nbins[i], xbins[i].data());
79 axis->SetName(TString::Format("axis%d", i));
80 fAxes.AddAtAndExpand(axis, i);
81 }
82 SetTitle(title);
84}
85
87 : TNamed(other), fNdimensions(other.fNdimensions), fAxes(fNdimensions), fBrowsables(fNdimensions),
88 fEntries(other.fEntries), fTsumw(other.fTsumw), fTsumw2(other.fTsumw2), fTsumwx(other.fTsumwx),
89 fTsumwx2(other.fTsumwx2), fIntegral(other.fIntegral), fIntegralStatus(other.fIntegralStatus)
90{
91
92 for (Int_t i = 0; i < fNdimensions; ++i) {
93 TAxis *axis = new TAxis(*static_cast<TAxis *>(other.fAxes[i]));
94 fAxes.AddAtAndExpand(axis, i);
95 }
97}
98
100{
101
102 if (this == &other)
103 return *this;
104
105 TNamed::operator=(other);
109 fEntries = other.fEntries;
110 fTsumw = other.fTsumw;
111 fTsumw2 = other.fTsumw2;
112 fTsumwx = other.fTsumwx;
113 fTsumwx2 = other.fTsumwx2;
114 fIntegral = other.fIntegral;
116
117 for (Int_t i = 0; i < fNdimensions; ++i) {
118 TAxis *axis = new TAxis(*static_cast<TAxis *>(other.fAxes[i]));
119 fAxes.AddAtAndExpand(axis, i);
120 }
121 fAxes.SetOwner();
122
123 return *this;
124}
125
127 : TNamed(std::move(other)), fNdimensions(other.fNdimensions), fAxes(other.fAxes), fBrowsables(fNdimensions),
128 fEntries(other.fEntries), fTsumw(other.fTsumw), fTsumw2(other.fTsumw2), fTsumwx(std::move(other.fTsumwx)),
129 fTsumwx2(std::move(other.fTsumwx2)), fIntegral(std::move(other.fIntegral)), fIntegralStatus(other.fIntegralStatus)
130{
131
132 other.fAxes.SetOwner(false);
133 other.fAxes.Clear();
134 fAxes.SetOwner();
135}
136
138{
139
140 if (this == &other)
141 return *this;
142
143 TNamed::operator=(std::move(other));
144 fNdimensions = other.fNdimensions;
145 fAxes = other.fAxes;
147 fEntries = other.fEntries;
148 fTsumw = other.fTsumw;
149 fTsumw2 = other.fTsumw2;
150 fTsumwx = std::move(other.fTsumwx);
151 fTsumwx2 = std::move(other.fTsumwx2);
152 fIntegral = std::move(other.fIntegral);
153 fIntegralStatus = other.fIntegralStatus;
154
155 other.fAxes.SetOwner(false);
156 other.fAxes.Clear();
157 fAxes.SetOwner();
158
159 return *this;
160}
161
162////////////////////////////////////////////////////////////////////////////////
163/// Destruct a THnBase
164
166 if (fIntegralStatus != kNoInt)
167 fIntegral.clear();
168}
169
170
171////////////////////////////////////////////////////////////////////////////////
172/// Create a new THnBase object that is of the same type as *this,
173/// but with dimensions and bins given by axes.
174/// If keepTargetAxis is true, the axes will keep their original xmin / xmax,
175/// else they will be restricted to the range selected (first / last).
176
177THnBase* THnBase::CloneEmpty(const char* name, const char* title,
178 const TObjArray* axes, Bool_t keepTargetAxis) const
179{
180 THnBase* ret = (THnBase*)IsA()->New();
181 Int_t chunkSize = 1024 * 16;
183 chunkSize = ((const THnSparse*)this)->GetChunkSize();
184 }
185 ret->Init(name, title, axes, keepTargetAxis, chunkSize);
186 return ret;
187}
188
189
190////////////////////////////////////////////////////////////////////////////////
191/// Initialize axes and name.
192
193void THnBase::Init(const char* name, const char* title,
194 const TObjArray* axes, Bool_t keepTargetAxis,
195 Int_t chunkSize /*= 1024 * 16*/)
196{
197 SetNameTitle(name, title);
198
199 TIter iAxis(axes);
200 const TAxis* axis = 0;
201 Int_t pos = 0;
202 Int_t *nbins = new Int_t[axes->GetEntriesFast()];
203 while ((axis = (TAxis*)iAxis())) {
204 TAxis* reqaxis = new TAxis(*axis);
205 if (!keepTargetAxis && axis->TestBit(TAxis::kAxisRange)) {
206 Int_t binFirst = axis->GetFirst();
207 // The lowest egde of the underflow is meaningless.
208 if (binFirst == 0)
209 binFirst = 1;
210 Int_t binLast = axis->GetLast();
211 // The overflow edge is implicit.
212 if (binLast > axis->GetNbins())
213 binLast = axis->GetNbins();
214 Int_t nBins = binLast - binFirst + 1;
215 if (axis->GetXbins()->GetSize()) {
216 // non-uniform bins:
217 reqaxis->Set(nBins, axis->GetXbins()->GetArray() + binFirst - 1);
218 } else {
219 // uniform bins:
220 reqaxis->Set(nBins, axis->GetBinLowEdge(binFirst), axis->GetBinUpEdge(binLast));
221 }
222 reqaxis->ResetBit(TAxis::kAxisRange);
223 }
224
225 nbins[pos] = reqaxis->GetNbins();
226 fAxes.AddAtAndExpand(reqaxis, pos++);
227 }
228 fAxes.SetOwner();
229
231 InitStorage(nbins, chunkSize);
234 delete [] nbins;
235}
236
237
238////////////////////////////////////////////////////////////////////////////////
239/// Create an empty histogram with name and title with a given
240/// set of axes. Create a TH1D/TH2D/TH3D, depending on the number
241/// of elements in axes.
242
243TH1* THnBase::CreateHist(const char* name, const char* title,
244 const TObjArray* axes,
245 Bool_t keepTargetAxis ) const {
246 const int ndim = axes->GetSize();
247
248 TH1* hist = 0;
249 // create hist with dummy axes, we fix them later.
250 if (ndim == 1)
251 hist = new TH1D(name, title, 1, 0., 1.);
252 else if (ndim == 2)
253 hist = new TH2D(name, title, 1, 0., 1., 1, 0., 1.);
254 else if (ndim == 3)
255 hist = new TH3D(name, title, 1, 0., 1., 1, 0., 1., 1, 0., 1.);
256 else {
257 Error("CreateHist", "Cannot create histogram %s with %d dimensions!", name, ndim);
258 return 0;
259 }
260
261 TAxis* hax[3] = {hist->GetXaxis(), hist->GetYaxis(), hist->GetZaxis()};
262 for (Int_t d = 0; d < ndim; ++d) {
263 TAxis* reqaxis = (TAxis*)(*axes)[d];
264 hax[d]->SetTitle(reqaxis->GetTitle());
265 if (!keepTargetAxis && reqaxis->TestBit(TAxis::kAxisRange)) {
266 // axis cannot extend to underflow/overflows (fix ROOT-8781)
267 Int_t binFirst = std::max(reqaxis->GetFirst(),1);
268 Int_t binLast = std::min(reqaxis->GetLast(), reqaxis->GetNbins() );
269 Int_t nBins = binLast - binFirst + 1;
270 if (reqaxis->GetXbins()->GetSize()) {
271 // non-uniform bins:
272 hax[d]->Set(nBins, reqaxis->GetXbins()->GetArray() + binFirst - 1);
273 } else {
274 // uniform bins:
275 hax[d]->Set(nBins, reqaxis->GetBinLowEdge(binFirst), reqaxis->GetBinUpEdge(binLast));
276 }
277 } else {
278 if (reqaxis->GetXbins()->GetSize()) {
279 // non-uniform bins:
280 hax[d]->Set(reqaxis->GetNbins(), reqaxis->GetXbins()->GetArray());
281 } else {
282 // uniform bins:
283 hax[d]->Set(reqaxis->GetNbins(), reqaxis->GetXmin(), reqaxis->GetXmax());
284 }
285 }
286 }
287
288 hist->Rebuild();
289
290 return hist;
291}
292
293////////////////////////////////////////////////////////////////////////////////
294/// Create a THn / THnSparse object from a histogram deriving from TH1.
295
296THnBase* THnBase::CreateHnAny(const char* name, const char* title,
297 const TH1* h, Bool_t sparse, Int_t chunkSize)
298{
299 // Get the dimension of the TH1
300 int ndim = h->GetDimension();
301
302 // Axis properties
303 int nbins[3] = {0,0,0};
304 double minRange[3] = {0.,0.,0.};
305 double maxRange[3] = {0.,0.,0.};
306 const TAxis* axis[3] = { h->GetXaxis(), h->GetYaxis(), h->GetZaxis() };
307 for (int i = 0; i < ndim; ++i) {
308 nbins[i] = axis[i]->GetNbins();
309 minRange[i] = axis[i]->GetXmin();
310 maxRange[i] = axis[i]->GetXmax();
311 }
312
313 // Create the corresponding THnSparse, depending on the storage
314 // type of the TH1. The class name will be "TH??\0" where the first
315 // ? is 1,2 or 3 and the second ? indicates the storage as C, S,
316 // I, F or D.
317 THnBase* s = 0;
318 const char* cname( h->ClassName() );
319 if (cname[0] == 'T' && cname[1] == 'H'
320 && cname[2] >= '1' && cname[2] <= '3' && cname[4] == 0) {
321
322#define R__THNBCASE(TAG) \
323if (sparse) { \
324s = new _NAME2_(THnSparse,TAG)(name, title, ndim, nbins, \
325minRange, maxRange, chunkSize); \
326} else { \
327s = new _NAME2_(THn,TAG)(name, title, ndim, nbins, \
328minRange, maxRange); \
329} \
330break;
331
332 switch (cname[3]) {
333 case 'F': R__THNBCASE(F);
334 case 'D': R__THNBCASE(D);
335 case 'I': R__THNBCASE(I);
336 case 'S': R__THNBCASE(S);
337 case 'C': R__THNBCASE(C);
338 }
339#undef R__THNBCASE
340 }
341 if (!s) {
342 ::Warning("THnSparse::CreateHnAny", "Unknown Type of Histogram");
343 return 0;
344 }
345
346 for (int i = 0; i < ndim; ++i) {
347 s->GetAxis(i)->SetTitle(axis[i]->GetTitle());
348 }
349
350 // Get the array to know the number of entries of the TH1
351 const TArray *array = dynamic_cast<const TArray*>(h);
352 if (!array) {
353 ::Warning("THnSparse::CreateHnAny", "Unknown Type of Histogram");
354 return 0;
355 }
356
357 s->Add(h);
358 return s;
359}
360
361
362////////////////////////////////////////////////////////////////////////////////
363/// Create a THnSparse (if "sparse") or THn from "hn", possibly
364/// converting THn <-> THnSparse.
365
366THnBase* THnBase::CreateHnAny(const char* name, const char* title,
367 const THnBase* hn, Bool_t sparse,
368 Int_t chunkSize /*= 1024 * 16*/)
369{
370 TClass* type = 0;
371 if (hn->InheritsFrom(THnSparse::Class())) {
372 if (sparse) type = hn->IsA();
373 else {
374 char bintype;
375 if (hn->InheritsFrom(THnSparseD::Class())) bintype = 'D';
376 else if (hn->InheritsFrom(THnSparseF::Class())) bintype = 'F';
377 else if (hn->InheritsFrom(THnSparseL::Class())) bintype = 'L';
378 else if (hn->InheritsFrom(THnSparseI::Class())) bintype = 'I';
379 else if (hn->InheritsFrom(THnSparseS::Class())) bintype = 'S';
380 else if (hn->InheritsFrom(THnSparseC::Class())) bintype = 'C';
381 else {
382 hn->Error("CreateHnAny", "Type %s not implemented; please inform the ROOT team!",
383 hn->IsA()->GetName());
384 return 0;
385 }
386 type = TClass::GetClass(TString::Format("THn%c", bintype));
387 }
388 } else if (hn->InheritsFrom(THn::Class())) {
389 if (!sparse) type = hn->IsA();
390 else {
391 char bintype = 0;
392 if (hn->InheritsFrom(THnD::Class())) bintype = 'D';
393 else if (hn->InheritsFrom(THnF::Class())) bintype = 'F';
394 else if (hn->InheritsFrom(THnC::Class())) bintype = 'C';
395 else if (hn->InheritsFrom(THnS::Class())) bintype = 'S';
396 else if (hn->InheritsFrom(THnI::Class())) bintype = 'I';
397 else if (hn->InheritsFrom(THnL::Class())) bintype = 'L';
398 else if (hn->InheritsFrom(THnL64::Class())) {
399 hn->Error("CreateHnAny", "Type THnSparse with Long64_t bins is not available!");
400 return 0;
401 }
402 if (bintype) {
403 type = TClass::GetClass(TString::Format("THnSparse%c", bintype));
404 }
405 }
406 } else {
407 hn->Error("CreateHnAny", "Unhandled type %s, not deriving from THn nor THnSparse!",
408 hn->IsA()->GetName());
409 return 0;
410 }
411 if (!type) {
412 hn->Error("CreateHnAny", "Unhandled type %s, please inform the ROOT team!",
413 hn->IsA()->GetName());
414 return 0;
415 }
416
417 THnBase* ret = (THnBase*)type->New();
418 ret->Init(name, title, hn->GetListOfAxes(),
419 kFALSE /*keepTargetAxes*/, chunkSize);
420
421 ret->Add(hn);
422 return ret;
423}
424
425
426////////////////////////////////////////////////////////////////////////////////
427/// Fill the THnBase with the bins of hist that have content
428/// or error != 0.
429
430void THnBase::Add(const TH1* hist, Double_t c /*=1.*/)
431{
432 Long64_t nbins = hist->GetNcells();
433 int x[3] = {0,0,0};
434 for (int i = 0; i < nbins; ++i) {
435 double value = hist->GetBinContent(i);
436 double error = hist->GetBinError(i);
437 if (!value && !error) continue;
438 hist->GetBinXYZ(i, x[0], x[1], x[2]);
440 SetBinError(x, error * c);
441 }
442}
443
444////////////////////////////////////////////////////////////////////////////////
445/// Fit a THnSparse with function f
446///
447/// since the data is sparse by default a likelihood fit is performed
448/// merging all the regions with empty bins for better performance efficiency
449///
450/// Since the THnSparse is not drawn no graphics options are passed
451/// Here is the list of possible options
452///
453/// = "I" Use integral of function in bin instead of value at bin center
454/// = "X" Use chi2 method (default is log-likelihood method)
455/// = "U" Use a User specified fitting algorithm (via SetFCN)
456/// = "Q" Quiet mode (minimum printing)
457/// = "V" Verbose mode (default is between Q and V)
458/// = "E" Perform better Errors estimation using Minos technique
459/// = "B" Use this option when you want to fix one or more parameters
460/// and the fitting function is like "gaus", "expo", "poln", "landau".
461/// = "M" More. Improve fit results
462/// = "R" Use the Range specified in the function range
463
465{
466
467 Foption_t fitOption;
468
469 if (!TH1::FitOptionsMake(option,fitOption)) return 0;
470
471 // The function used to fit cannot be stored in a THnSparse. It
472 // cannot be drawn either. Perhaps in the future.
473 fitOption.Nostore = true;
474 // Use likelihood fit if not specified
475 if (!fitOption.Chi2) fitOption.Like = true;
476 // create range and minimizer options with default values
478 for ( int i = 0; i < GetNdimensions(); ++i ) {
479 TAxis *axis = GetAxis(i);
480 range.AddRange(i, axis->GetXmin(), axis->GetXmax());
481 }
483
484 return ROOT::Fit::FitObject(this, f , fitOption , minOption, goption, range);
485}
486
487////////////////////////////////////////////////////////////////////////////////
488/// Generate an n-dimensional random tuple based on the histogrammed
489/// distribution. If subBinRandom, the returned tuple will be additionally
490/// randomly distributed within the randomized bin, using a flat
491/// distribution.
492
493void THnBase::GetRandom(Double_t *rand, Bool_t subBinRandom /* = kTRUE */)
494{
495 // check whether the integral array is valid
498
499 // generate a random bin
500 Double_t p = gRandom->Rndm();
501 Long64_t idx = TMath::BinarySearch(GetNbins() + 1, fIntegral.data(), p);
502 const Int_t nStaticBins = 40;
503 Int_t bin[nStaticBins];
504 Int_t* pBin = bin;
505 if (GetNdimensions() > nStaticBins) {
506 pBin = new Int_t[GetNdimensions()];
507 }
508 GetBinContent(idx, pBin);
509
510 // convert bin coordinates to real values
511 for (Int_t i = 0; i < fNdimensions; i++) {
512 rand[i] = GetAxis(i)->GetBinCenter(pBin[i]);
513
514 // randomize the vector within a bin
515 if (subBinRandom)
516 rand[i] += (gRandom->Rndm() - 0.5) * GetAxis(i)->GetBinWidth(pBin[i]);
517 }
518 if (pBin != bin) {
519 delete [] pBin;
520 }
521
522 return;
523}
524
525////////////////////////////////////////////////////////////////////////////////
526/// Check whether bin coord is in range, as defined by TAxis::SetRange().
527
529{
530 Int_t min = 0;
531 Int_t max = 0;
532 for (Int_t i = 0; i < fNdimensions; ++i) {
533 TAxis *axis = GetAxis(i);
534 if (!axis->TestBit(TAxis::kAxisRange)) continue;
535 min = axis->GetFirst();
536 max = axis->GetLast();
537 if (coord[i] < min || coord[i] > max)
538 return kFALSE;
539 }
540 return kTRUE;
541}
542
543////////////////////////////////////////////////////////////////////////////////
544/// Project all bins into a ndim-dimensional THn / THnSparse (whatever
545/// *this is) or if (ndim < 4 and !wantNDim) a TH1/2/3 histogram,
546/// keeping only axes in dim (specifying ndim dimensions).
547/// If "option" contains:
548/// - "E" errors will be calculated.
549/// - "A" ranges of the target axes will be ignored.
550/// - "O" original axis range of the target axes will be
551/// kept, but only bins inside the selected range
552/// will be filled.
553
555 Bool_t wantNDim,
556 Option_t* option /*= ""*/) const
557{
559 name +="_proj";
560
561 for (Int_t d = 0; d < ndim; ++d) {
562 name += "_";
563 name += dim[d];
564 }
565
566 TString title(GetTitle());
567 Ssiz_t posInsert = title.First(';');
568 if (posInsert == kNPOS) {
569 title += " projection ";
570 for (Int_t d = 0; d < ndim; ++d)
571 title += GetAxis(dim[d])->GetTitle();
572 } else {
573 for (Int_t d = ndim - 1; d >= 0; --d) {
574 title.Insert(posInsert, GetAxis(d)->GetTitle());
575 if (d)
576 title.Insert(posInsert, ", ");
577 }
578 title.Insert(posInsert, " projection ");
579 }
580
581 TObjArray newaxes(ndim);
582 for (Int_t d = 0; d < ndim; ++d) {
583 newaxes.AddAt(GetAxis(dim[d]),d);
584 }
585
586 THnBase* hn = 0;
587 TH1* hist = 0;
588 TObject* ret = 0;
589
590 Bool_t* hadRange = 0;
591 Bool_t ignoreTargetRange = (option && (strchr(option, 'A') || strchr(option, 'a')));
592 Bool_t keepTargetAxis = ignoreTargetRange || (option && (strchr(option, 'O') || strchr(option, 'o')));
593 if (ignoreTargetRange) {
594 hadRange = new Bool_t[ndim];
595 for (Int_t d = 0; d < ndim; ++d){
596 TAxis *axis = GetAxis(dim[d]);
597 hadRange[d] = axis->TestBit(TAxis::kAxisRange);
599 }
600 }
601
602 if (wantNDim)
603 ret = hn = CloneEmpty(name, title, &newaxes, keepTargetAxis);
604 else
605 ret = hist = CreateHist(name, title, &newaxes, keepTargetAxis);
606
607 if (keepTargetAxis) {
608 // make the whole axes visible, i.e. unset the range
609 if (wantNDim) {
610 for (Int_t d = 0; d < ndim; ++d) {
611 hn->GetAxis(d)->SetRange(0, 0);
612 }
613 } else {
614 hist->GetXaxis()->SetRange(0, 0);
615 hist->GetYaxis()->SetRange(0, 0);
616 hist->GetZaxis()->SetRange(0, 0);
617 }
618 }
619
620 Bool_t haveErrors = GetCalculateErrors();
621 Bool_t wantErrors = haveErrors || (option && (strchr(option, 'E') || strchr(option, 'e')));
622
623 Int_t* bins = new Int_t[ndim];
624 Long64_t myLinBin = 0;
625
626 THnIter iter(this, kTRUE /*use axis range*/);
627
628 while ((myLinBin = iter.Next()) >= 0) {
629 Double_t v = GetBinContent(myLinBin);
630
631 for (Int_t d = 0; d < ndim; ++d) {
632 bins[d] = iter.GetCoord(dim[d]);
633 if (!keepTargetAxis && GetAxis(dim[d])->TestBit(TAxis::kAxisRange)) {
634 Int_t binOffset = GetAxis(dim[d])->GetFirst();
635 // Don't subtract even more if underflow is alreday included:
636 if (binOffset > 0) --binOffset;
637 bins[d] -= binOffset;
638 }
639 }
640
641 Long64_t targetLinBin = -1;
642 if (!wantNDim) {
643 if (ndim == 1) targetLinBin = bins[0];
644 else if (ndim == 2) targetLinBin = hist->GetBin(bins[0], bins[1]);
645 else if (ndim == 3) targetLinBin = hist->GetBin(bins[0], bins[1], bins[2]);
646 } else {
647 targetLinBin = hn->GetBin(bins, kTRUE /*allocate*/);
648 }
649
650 if (wantErrors) {
651 Double_t err2 = 0.;
652 if (haveErrors) {
653 err2 = GetBinError2(myLinBin);
654 } else {
655 err2 = v;
656 }
657 if (wantNDim) {
658 hn->AddBinError2(targetLinBin, err2);
659 } else {
660 Double_t preverr = hist->GetBinError(targetLinBin);
661 hist->SetBinError(targetLinBin, TMath::Sqrt(preverr * preverr + err2));
662 }
663 }
664
665 // only _after_ error calculation, or sqrt(v) is taken into account!
666 if (wantNDim)
667 hn->AddBinContent(targetLinBin, v);
668 else
669 hist->AddBinContent(targetLinBin, v);
670 }
671
672 delete [] bins;
673
674 if (wantNDim) {
675 hn->SetEntries(fEntries);
676 } else {
677 if (!iter.HaveSkippedBin()) {
678 hist->SetEntries(fEntries);
679 } else {
680 // re-compute the entries
681 // in case of error calculation (i.e. when Sumw2() is set)
682 // use the effective entries for the entries
683 // since this is the only way to estimate them
684 hist->ResetStats();
685 Double_t entries = hist->GetEffectiveEntries();
686 if (!wantErrors) {
687 // to avoid numerical rounding
688 entries = TMath::Floor(entries + 0.5);
689 }
690 hist->SetEntries(entries);
691 }
692 }
693
694 if (hadRange) {
695 // reset kAxisRange bit:
696 for (Int_t d = 0; d < ndim; ++d)
697 GetAxis(dim[d])->SetBit(TAxis::kAxisRange, hadRange[d]);
698
699 delete [] hadRange;
700 }
701
702 return ret;
703}
704
705////////////////////////////////////////////////////////////////////////////////
706/// Scale contents and errors of this histogram by c:
707/// this = this * c
708/// It does not modify the histogram's number of entries.
709
711{
712
713 Double_t nEntries = GetEntries();
714 // Scale the contents & errors
715 Bool_t haveErrors = GetCalculateErrors();
716 Long64_t i = 0;
717 THnIter iter(this);
718 while ((i = iter.Next()) >= 0) {
719 // Get the content of the bin from the current histogram
721 SetBinContent(i, c * v);
722 if (haveErrors) {
723 Double_t err2 = GetBinError2(i);
724 SetBinError2(i, c * c * err2);
725 }
726 }
727 SetEntries(nEntries);
728}
729
730////////////////////////////////////////////////////////////////////////////////
731/// Add() implementation for both rebinned histograms and those with identical
732/// binning. See THnBase::Add().
733
735{
736 if (fNdimensions != h->GetNdimensions()) {
737 Warning("RebinnedAdd", "Different number of dimensions, cannot carry out operation on the histograms");
738 return;
739 }
740
741 // Trigger error calculation if h has it
742 if (!GetCalculateErrors() && h->GetCalculateErrors())
743 Sumw2();
744 Bool_t haveErrors = GetCalculateErrors();
745
746 Double_t* x = 0;
747 if (rebinned) {
748 x = new Double_t[fNdimensions];
749 }
750 Int_t* coord = new Int_t[fNdimensions];
751
752 // Expand the exmap if needed, to reduce collisions
753 Long64_t numTargetBins = GetNbins() + h->GetNbins();
754 Reserve(numTargetBins);
755
756 Long64_t i = 0;
757 THnIter iter(h);
758 // Add to this whatever is found inside the other histogram
759 while ((i = iter.Next(coord)) >= 0) {
760 // Get the content of the bin from the second histogram
761 Double_t v = h->GetBinContent(i);
762
763 Long64_t mybinidx = -1;
764 if (rebinned) {
765 // Get the bin center given a coord
766 for (Int_t j = 0; j < fNdimensions; ++j)
767 x[j] = h->GetAxis(j)->GetBinCenter(coord[j]);
768
769 mybinidx = GetBin(x, kTRUE /* allocate*/);
770 } else {
771 mybinidx = GetBin(coord, kTRUE /*allocate*/);
772 }
773
774 if (haveErrors) {
775 Double_t err2 = h->GetBinError2(i) * c * c;
776 AddBinError2(mybinidx, err2);
777 }
778 // only _after_ error calculation, or sqrt(v) is taken into account!
779 AddBinContent(mybinidx, c * v);
780 }
781
782 delete [] coord;
783 delete [] x;
784
785 // add also the statistics
786 fTsumw += c * h->fTsumw;
787 if (haveErrors) {
788 fTsumw2 += c * c * h->fTsumw2;
789 if (h->fTsumwx.fN == fNdimensions && h->fTsumwx2.fN == fNdimensions) {
790 for (Int_t d = 0; d < fNdimensions; ++d) {
791 fTsumwx[d] += c * h->fTsumwx[d];
792 fTsumwx2[d] += c * h->fTsumwx2[d];
793 }
794 }
795 }
796
797 Double_t nEntries = GetEntries() + c * h->GetEntries();
798 SetEntries(nEntries);
799}
800
801////////////////////////////////////////////////////////////////////////////////
802/// Add contents of h scaled by c to this histogram:
803/// this = this + c * h
804/// Note that if h has Sumw2 set, Sumw2 is automatically called for this
805/// if not already set.
806
808{
809 // Check consistency of the input
810 if (!CheckConsistency(h, "Add")) return;
811
813}
814
815////////////////////////////////////////////////////////////////////////////////
816/// Add contents of h scaled by c to this histogram:
817/// this = this + c * h
818/// Note that if h has Sumw2 set, Sumw2 is automatically called for this
819/// if not already set.
820/// In contrast to Add(), RebinnedAdd() does not require consistent binning of
821/// this and h; instead, each bin's center is used to determine the target bin.
822
824{
825 AddInternal(h, c, kTRUE);
826}
827
828
829////////////////////////////////////////////////////////////////////////////////
830/// Merge this with a list of THnBase's. All THnBase's provided
831/// in the list must have the same bin layout!
832
834{
835 if (!list) return 0;
836 if (list->IsEmpty()) return (Long64_t)GetEntries();
837
838 Long64_t sumNbins = GetNbins();
839 TIter iter(list);
840 const TObject* addMeObj = 0;
841 while ((addMeObj = iter())) {
842 const THnBase* addMe = dynamic_cast<const THnBase*>(addMeObj);
843 if (addMe) {
844 sumNbins += addMe->GetNbins();
845 }
846 }
847 Reserve(sumNbins);
848
849 iter.Reset();
850 while ((addMeObj = iter())) {
851 const THnBase* addMe = dynamic_cast<const THnBase*>(addMeObj);
852 if (!addMe)
853 Error("Merge", "Object named %s is not THnBase! Skipping it.",
854 addMeObj->GetName());
855 else
856 Add(addMe);
857 }
858 return (Long64_t)GetEntries();
859}
860
861
862////////////////////////////////////////////////////////////////////////////////
863/// Multiply this histogram by histogram h
864/// this = this * h
865/// Note that if h has Sumw2 set, Sumw2 is automatically called for this
866/// if not already set.
867
869{
870 // Check consistency of the input
871 if(!CheckConsistency(h, "Multiply"))return;
872
873 // Trigger error calculation if h has it
874 Bool_t wantErrors = kFALSE;
875 if (GetCalculateErrors() || h->GetCalculateErrors())
876 wantErrors = kTRUE;
877
878 if (wantErrors) Sumw2();
879
880 Double_t nEntries = GetEntries();
881 // Now multiply the contents: in this case we have the intersection of the sets of bins
882 Int_t* coord = new Int_t[fNdimensions];
883 Long64_t i = 0;
884 THnIter iter(this);
885 // Add to this whatever is found inside the other histogram
886 while ((i = iter.Next(coord)) >= 0) {
887 // Get the content of the bin from the current histogram
889 // Now look at the bin with the same coordinates in h
890 Long64_t idxh = h->GetBin(coord);
891 Double_t v2 = 0.;
892 if (idxh >= 0) v2 = h->GetBinContent(idxh);
893 SetBinContent(i, v1 * v2);
894 if (wantErrors) {
895 Double_t err1 = GetBinError(i) * v2;
896 Double_t err2 = 0.;
897 if (idxh >= 0) err2 = h->GetBinError(idxh) * v1;
898 SetBinError(i, TMath::Sqrt((err2 * err2 + err1 * err1)));
899 }
900 }
901 SetEntries(nEntries);
902
903 delete [] coord;
904}
905
906////////////////////////////////////////////////////////////////////////////////
907/// Performs the operation: this = this*c*f1
908/// if errors are defined, errors are also recalculated.
909///
910/// Only bins inside the function range are recomputed.
911/// IMPORTANT NOTE: If you intend to use the errors of this histogram later
912/// you should call Sumw2 before making this operation.
913/// This is particularly important if you fit the histogram after
914/// calling Multiply()
915
917{
918 Int_t* coord = new Int_t[fNdimensions];
920
921 Bool_t wantErrors = GetCalculateErrors();
922 if (wantErrors) Sumw2();
923
924 Long64_t i = 0;
925 THnIter iter(this);
926 // Add to this whatever is found inside the other histogram
927 while ((i = iter.Next(coord)) >= 0) {
929
930 // Get the bin coordinates given an index array
931 for (Int_t j = 0; j < fNdimensions; ++j)
932 x[j] = GetAxis(j)->GetBinCenter(coord[j]);
933
934 if (!f->IsInside(x))
935 continue;
937
938 // Evaluate function at points
939 Double_t fvalue = f->EvalPar(x, NULL);
940
941 SetBinContent(i, c * fvalue * value);
942 if (wantErrors) {
943 Double_t error = GetBinError(i);
944 SetBinError(i, c * fvalue * error);
945 }
946 }
947
948 delete [] x;
949 delete [] coord;
950}
951
952////////////////////////////////////////////////////////////////////////////////
953/// Divide this histogram by h
954/// this = this/(h)
955/// Note that if h has Sumw2 set, Sumw2 is automatically called for
956/// this if not already set.
957/// The resulting errors are calculated assuming uncorrelated content.
958
960{
961 // Check consistency of the input
962 if (!CheckConsistency(h, "Divide"))return;
963
964 // Trigger error calculation if h has it
965 Bool_t wantErrors=GetCalculateErrors();
966 if (!GetCalculateErrors() && h->GetCalculateErrors())
967 wantErrors=kTRUE;
968
969 // Remember original histogram statistics
970 Double_t nEntries = fEntries;
971
972 if (wantErrors) Sumw2();
973 Bool_t didWarn = kFALSE;
974
975 // Now divide the contents: also in this case we have the intersection of the sets of bins
976 Int_t* coord = new Int_t[fNdimensions];
977 Long64_t i = 0;
978 THnIter iter(this);
979 // Add to this whatever is found inside the other histogram
980 while ((i = iter.Next(coord)) >= 0) {
981 // Get the content of the bin from the first histogram
983 // Now look at the bin with the same coordinates in h
984 Long64_t hbin = h->GetBin(coord);
985 Double_t v2 = h->GetBinContent(hbin);
986 if (!v2) {
987 v1 = 0.;
988 v2 = 1.;
989 if (!didWarn) {
990 Warning("Divide(h)", "Histogram h has empty bins - division by zero! Setting bin to 0.");
991 didWarn = kTRUE;
992 }
993 }
994 SetBinContent(i, v1 / v2);
995 if (wantErrors) {
996 Double_t err1 = GetBinError(i) * v2;
997 Double_t err2 = h->GetBinError(hbin) * v1;
998 Double_t b22 = v2 * v2;
999 Double_t err = (err1 * err1 + err2 * err2) / (b22 * b22);
1000 SetBinError2(i, err);
1001 }
1002 }
1003 delete [] coord;
1004 SetEntries(nEntries);
1005}
1006
1007////////////////////////////////////////////////////////////////////////////////
1008/// Replace contents of this histogram by multiplication of h1 by h2
1009/// this = (c1*h1)/(c2*h2)
1010/// Note that if h1 or h2 have Sumw2 set, Sumw2 is automatically called for
1011/// this if not already set.
1012/// The resulting errors are calculated assuming uncorrelated content.
1013/// However, if option ="B" is specified, Binomial errors are computed.
1014/// In this case c1 and c2 do not make real sense and they are ignored.
1015
1017{
1018
1019 TString opt = option;
1020 opt.ToLower();
1021 Bool_t binomial = kFALSE;
1022 if (opt.Contains("b")) binomial = kTRUE;
1023
1024 // Check consistency of the input
1025 if (!CheckConsistency(h1, "Divide") || !CheckConsistency(h2, "Divide"))return;
1026 if (!c2) {
1027 Error("Divide","Coefficient of dividing histogram cannot be zero");
1028 return;
1029 }
1030
1031 Reset();
1032
1033 // Trigger error calculation if h1 or h2 have it
1034 if (!GetCalculateErrors() && (h1->GetCalculateErrors()|| h2->GetCalculateErrors() != 0))
1035 Sumw2();
1036
1037 // Count filled bins
1038 Long64_t nFilledBins=0;
1039
1040 // Now divide the contents: we have the intersection of the sets of bins
1041
1042 Int_t* coord = new Int_t[fNdimensions];
1043 memset(coord, 0, sizeof(Int_t) * fNdimensions);
1044 Bool_t didWarn = kFALSE;
1045
1046 Long64_t i = 0;
1047 THnIter iter(h1);
1048 // Add to this whatever is found inside the other histogram
1049 while ((i = iter.Next(coord)) >= 0) {
1050 // Get the content of the bin from the first histogram
1052 // Now look at the bin with the same coordinates in h2
1053 Long64_t h2bin = h2->GetBin(coord);
1054 Double_t v2 = h2->GetBinContent(h2bin);
1055 if (!v2) {
1056 v1 = 0.;
1057 v2 = 1.;
1058 if (!didWarn) {
1059 Warning("Divide(h1, h2)", "Histogram h2 has empty bins - division by zero! Setting bin to 0.");
1060 didWarn = kTRUE;
1061 }
1062 }
1063 nFilledBins++;
1064 Long64_t myBin = GetBin(coord);
1065 SetBinContent(myBin, c1 * v1 / c2 / v2);
1066 if(GetCalculateErrors()){
1067 Double_t err1 = h1->GetBinError(i);
1068 Double_t err2 = h2->GetBinError(h2bin);
1069 Double_t errSq = 0.;
1070 if (binomial) {
1071 if (v1 != v2) {
1072 Double_t w = v1 / v2;
1073 err2 *= w;
1074 errSq = TMath::Abs( ( (1. - 2.*w) * err1 * err1 + err2 * err2 ) / (v2 * v2) );
1075 }
1076 } else {
1077 c1 *= c1;
1078 c2 *= c2;
1079 Double_t b22 = v2 * v2 * c2;
1080 err1 *= v2;
1081 err2 *= v1;
1082 errSq = c1 * c2 * (err1 * err1 + err2 * err2) / (b22 * b22);
1083 }
1084 SetBinError2(myBin, errSq);
1085 }
1086 }
1087
1088 delete [] coord;
1089 SetFilledBins(nFilledBins);
1090
1091 // Set as entries in the result histogram the entries in the numerator
1093}
1094
1095////////////////////////////////////////////////////////////////////////////////
1096/// Consistency check on (some of) the parameters of two histograms (for operations).
1097
1098Bool_t THnBase::CheckConsistency(const THnBase *h, const char *tag) const
1099{
1100 if (fNdimensions != h->GetNdimensions()) {
1101 Warning(tag, "Different number of dimensions, cannot carry out operation on the histograms");
1102 return kFALSE;
1103 }
1104 for (Int_t dim = 0; dim < fNdimensions; dim++){
1105 if (GetAxis(dim)->GetNbins() != h->GetAxis(dim)->GetNbins()) {
1106 Warning(tag, "Different number of bins on axis %i, cannot carry out operation on the histograms", dim);
1107 return kFALSE;
1108 }
1109 }
1110 return kTRUE;
1111}
1112
1113////////////////////////////////////////////////////////////////////////////////
1114/// Set the axis # of bins and bin limits on dimension idim
1115
1116void THnBase::SetBinEdges(Int_t idim, const Double_t* bins)
1117{
1118 TAxis* axis = (TAxis*) fAxes[idim];
1119 axis->Set(axis->GetNbins(), bins);
1120}
1121
1122////////////////////////////////////////////////////////////////////////////////
1123/// Change (i.e. set) the title.
1124///
1125/// If title is in the form "stringt;string0;string1;string2 ..."
1126/// the histogram title is set to stringt, the title of axis0 to string0,
1127/// of axis1 to string1, of axis2 to string2, etc, just like it is done
1128/// for TH1/TH2/TH3.
1129/// To insert the character ";" in one of the titles, one should use "#;"
1130/// or "#semicolon".
1131
1132void THnBase::SetTitle(const char *title)
1133{
1134 fTitle = title;
1135 fTitle.ReplaceAll("#;",2,"#semicolon",10);
1136
1137 Int_t endHistTitle = fTitle.First(';');
1138 if (endHistTitle >= 0) {
1139 // title contains a ';' so parse the axis titles
1140 Int_t posTitle = endHistTitle + 1;
1141 Int_t lenTitle = fTitle.Length();
1142 Int_t dim = 0;
1143 while (posTitle > 0 && posTitle < lenTitle && dim < fNdimensions){
1144 Int_t endTitle = fTitle.Index(";", posTitle);
1145 TString axisTitle = fTitle(posTitle, endTitle - posTitle);
1146 axisTitle.ReplaceAll("#semicolon", 10, ";", 1);
1147 GetAxis(dim)->SetTitle(axisTitle);
1148 dim++;
1149 if (endTitle > 0)
1150 posTitle = endTitle + 1;
1151 else
1152 posTitle = -1;
1153 }
1154 // Remove axis titles from histogram title
1155 fTitle.Remove(endHistTitle, lenTitle - endHistTitle);
1156 }
1157
1158 fTitle.ReplaceAll("#semicolon", 10, ";", 1);
1159
1160}
1161
1162////////////////////////////////////////////////////////////////////////////////
1163/// Combine the content of "group" neighboring bins into
1164/// a new bin and return the resulting THnBase.
1165/// For group=2 and a 3 dimensional histogram, all "blocks"
1166/// of 2*2*2 bins will be put into a bin.
1167
1169{
1170 Int_t* ngroup = new Int_t[GetNdimensions()];
1171 for (Int_t d = 0; d < GetNdimensions(); ++d)
1172 ngroup[d] = group;
1173 THnBase* ret = RebinBase(ngroup);
1174 delete [] ngroup;
1175 return ret;
1176}
1177
1178////////////////////////////////////////////////////////////////////////////////
1179/// Combine the content of "group" neighboring bins for each dimension
1180/// into a new bin and return the resulting THnBase.
1181/// For group={2,1,1} and a 3 dimensional histogram, pairs of x-bins
1182/// will be grouped.
1183
1185{
1186 Int_t ndim = GetNdimensions();
1187 TString name(GetName());
1188 for (Int_t d = 0; d < ndim; ++d)
1189 name += Form("_%d", group[d]);
1190
1191
1192 TString title(GetTitle());
1193 Ssiz_t posInsert = title.First(';');
1194 if (posInsert == kNPOS) {
1195 title += " rebin ";
1196 for (Int_t d = 0; d < ndim; ++d)
1197 title += Form("{%d}", group[d]);
1198 } else {
1199 for (Int_t d = ndim - 1; d >= 0; --d)
1200 title.Insert(posInsert, Form("{%d}", group[d]));
1201 title.Insert(posInsert, " rebin ");
1202 }
1203
1204 TObjArray newaxes(ndim);
1205 newaxes.SetOwner();
1206 for (Int_t d = 0; d < ndim; ++d) {
1207 newaxes.AddAt(new TAxis(*GetAxis(d) ),d);
1208 if (group[d] > 1) {
1209 TAxis* newaxis = (TAxis*) newaxes.At(d);
1210 Int_t newbins = (newaxis->GetNbins() + group[d] - 1) / group[d];
1211 if (newaxis->GetXbins() && newaxis->GetXbins()->GetSize()) {
1212 // variable bins
1213 Double_t *edges = new Double_t[newbins + 1];
1214 for (Int_t i = 0; i < newbins + 1; ++i)
1215 if (group[d] * i <= newaxis->GetNbins())
1216 edges[i] = newaxis->GetXbins()->At(group[d] * i);
1217 else edges[i] = newaxis->GetXmax();
1218 newaxis->Set(newbins, edges);
1219 delete [] edges;
1220 } else {
1221 newaxis->Set(newbins, newaxis->GetXmin(), newaxis->GetXmax());
1222 }
1223 }
1224 }
1225
1226 THnBase* h = CloneEmpty(name.Data(), title.Data(), &newaxes, kTRUE);
1227 Bool_t haveErrors = GetCalculateErrors();
1228 Bool_t wantErrors = haveErrors;
1229
1230 Int_t* bins = new Int_t[ndim];
1231 Int_t* coord = new Int_t[fNdimensions];
1232
1233 Long64_t i = 0;
1234 THnIter iter(this);
1235 while ((i = iter.Next(coord)) >= 0) {
1237 for (Int_t d = 0; d < ndim; ++d) {
1238 bins[d] = TMath::CeilNint( (double) coord[d]/group[d] );
1239 }
1240 Long64_t idxh = h->GetBin(bins, kTRUE /*allocate*/);
1241
1242 if (wantErrors) {
1243 // wantErrors == haveErrors, thus:
1244 h->AddBinError2(idxh, GetBinError2(i));
1245 }
1246
1247 // only _after_ error calculation, or sqrt(v) is taken into account!
1248 h->AddBinContent(idxh, v);
1249 }
1250
1251 delete [] bins;
1252 delete [] coord;
1253
1254 h->SetEntries(fEntries);
1255
1256 return h;
1257
1258}
1259
1260////////////////////////////////////////////////////////////////////////////////
1261/// Clear the histogram
1262
1263void THnBase::ResetBase(Option_t * /*option = ""*/)
1264{
1265 fEntries = 0.;
1266 fTsumw = 0.;
1267 fTsumw2 = -1.;
1268 if (fIntegralStatus != kNoInt) {
1269 fIntegral.clear();
1271 }
1272}
1273
1274////////////////////////////////////////////////////////////////////////////////
1275/// Calculate the integral of the histogram
1276
1278{
1279 // delete old integral
1280 if (fIntegralStatus != kNoInt) {
1281 fIntegral.clear();
1283 }
1284
1285 // check number of bins
1286 if (GetNbins() == 0) {
1287 Error("ComputeIntegral", "The histogram must have at least one bin.");
1288 return 0.;
1289 }
1290
1291 // allocate integral array
1292 fIntegral.resize(GetNbins() + 1);
1293 fIntegral[0] = 0.;
1294
1295 // fill integral array with contents of regular bins (non over/underflow)
1296 Int_t* coord = new Int_t[fNdimensions];
1297 Long64_t i = 0;
1298 THnIter iter(this);
1299 while ((i = iter.Next(coord)) >= 0) {
1301
1302 // check whether the bin is regular
1303 bool regularBin = true;
1304 for (Int_t dim = 0; dim < fNdimensions; dim++) {
1305 if (coord[dim] < 1 || coord[dim] > GetAxis(dim)->GetNbins()) {
1306 regularBin = false;
1307 break;
1308 }
1309 }
1310
1311 // if outlayer, count it with zero weight
1312 if (!regularBin) v = 0.;
1313
1314 fIntegral[i + 1] = fIntegral[i] + v;
1315 }
1316 delete [] coord;
1317
1318 // check sum of weights
1319 if (fIntegral[GetNbins()] == 0.) {
1320 Error("ComputeIntegral", "No hits in regular bins (non over/underflow).");
1321 fIntegral.clear();
1322 return 0.;
1323 }
1324
1325 // normalize the integral array
1326 for (Long64_t j = 0; j <= GetNbins(); ++j)
1328
1329 // set status to valid
1331 return fIntegral[GetNbins()];
1332}
1333
1334////////////////////////////////////////////////////////////////////////////////
1335/// Print bin with linex index "idx".
1336/// For valid options see PrintBin(Long64_t idx, Int_t* bin, Option_t* options).
1337
1338void THnBase::PrintBin(Long64_t idx, Option_t* options) const
1339{
1340 Int_t* coord = new Int_t[fNdimensions];
1341 PrintBin(idx, coord, options);
1342 delete [] coord;
1343}
1344
1345////////////////////////////////////////////////////////////////////////////////
1346/// Print one bin. If "idx" is != -1 use that to determine the bin,
1347/// otherwise (if "idx" == -1) use the coordinate in "bin".
1348/// If "options" contains:
1349/// - '0': only print bins with an error or content != 0
1350/// Return whether the bin was printed (depends on options)
1351
1353{
1354 Double_t v = -42;
1355 if (idx == -1) {
1356 idx = GetBin(bin);
1357 v = GetBinContent(idx);
1358 } else {
1359 v = GetBinContent(idx, bin);
1360 }
1361
1362 Double_t err = 0.;
1363 if (GetCalculateErrors()) {
1364 if (idx != -1) {
1365 err = GetBinError(idx);
1366 }
1367 }
1368
1369 if (v == 0. && err == 0. && options && strchr(options, '0')) {
1370 // suppress zeros, and we have one.
1371 return kFALSE;
1372 }
1373
1374 TString coord;
1375 for (Int_t dim = 0; dim < fNdimensions; ++dim) {
1376 coord += bin[dim];
1377 coord += ',';
1378 }
1379 coord.Remove(coord.Length() - 1);
1380
1381 if (GetCalculateErrors()) {
1382 Printf("Bin at (%s) = %g (+/- %g)", coord.Data(), v, err);
1383 } else {
1384 Printf("Bin at (%s) = %g", coord.Data(), v);
1385 }
1386
1387 return kTRUE;
1388}
1389
1390////////////////////////////////////////////////////////////////////////////////
1391/// Print "howmany" entries starting at "from". If "howmany" is -1, print all.
1392/// If "options" contains:
1393/// - 'x': print in the order of axis bins, i.e. (0,0,...,0), (0,0,...,1),...
1394/// - '0': only print bins with content != 0
1395
1396void THnBase::PrintEntries(Long64_t from /*=0*/, Long64_t howmany /*=-1*/,
1397 Option_t* options /*=0*/) const
1398{
1399 if (from < 0) from = 0;
1400 if (howmany == -1) howmany = GetNbins();
1401
1402 Int_t* bin = new Int_t[fNdimensions];
1403
1404 if (options && (strchr(options, 'x') || strchr(options, 'X'))) {
1405 Int_t* nbins = new Int_t[fNdimensions];
1406 for (Int_t dim = fNdimensions - 1; dim >= 0; --dim) {
1407 nbins[dim] = GetAxis(dim)->GetNbins();
1408 bin[dim] = from % nbins[dim];
1409 from /= nbins[dim];
1410 }
1411
1412 for (Long64_t i = 0; i < howmany; ++i) {
1413 if (!PrintBin(-1, bin, options))
1414 ++howmany;
1415 // Advance to next bin:
1416 ++bin[fNdimensions - 1];
1417 for (Int_t dim = fNdimensions - 1; dim >= 0; --dim) {
1418 if (bin[dim] >= nbins[dim]) {
1419 bin[dim] = 0;
1420 if (dim > 0) {
1421 ++bin[dim - 1];
1422 } else {
1423 howmany = -1; // aka "global break"
1424 }
1425 }
1426 }
1427 }
1428 delete [] nbins;
1429 } else {
1430 for (Long64_t i = from; i < from + howmany; ++i) {
1431 if (!PrintBin(i, bin, options))
1432 ++howmany;
1433 }
1434 }
1435 delete [] bin;
1436}
1437
1438////////////////////////////////////////////////////////////////////////////////
1439/// Print a THnBase. If "option" contains:
1440/// - 'a': print axis details
1441/// - 'm': print memory usage
1442/// - 's': print statistics
1443/// - 'c': print its content, too (this can generate a LOT of output!)
1444/// Other options are forwarded to PrintEntries().
1445
1446void THnBase::Print(Option_t* options) const
1447{
1448 Bool_t optAxis = options && (strchr(options, 'A') || (strchr(options, 'a')));
1449 Bool_t optMem = options && (strchr(options, 'M') || (strchr(options, 'm')));
1450 Bool_t optStat = options && (strchr(options, 'S') || (strchr(options, 's')));
1451 Bool_t optContent = options && (strchr(options, 'C') || (strchr(options, 'c')));
1452
1453 Printf("%s (*0x%zx): \"%s\" \"%s\"", IsA()->GetName(), (size_t)this, GetName(), GetTitle());
1454 Printf(" %d dimensions, %g entries in %lld filled bins", GetNdimensions(), GetEntries(), GetNbins());
1455
1456 if (optAxis) {
1457 for (Int_t dim = 0; dim < fNdimensions; ++dim) {
1458 TAxis* axis = GetAxis(dim);
1459 Printf(" axis %d \"%s\": %d bins (%g..%g), %s bin sizes", dim,
1460 axis->GetTitle(), axis->GetNbins(), axis->GetXmin(), axis->GetXmax(),
1461 (axis->GetXbins() ? "variable" : "fixed"));
1462 }
1463 }
1464
1465 if (optStat) {
1466 Printf(" %s error calculation", (GetCalculateErrors() ? "with" : "without"));
1467 if (GetCalculateErrors()) {
1468 Printf(" Sum(w)=%g, Sum(w^2)=%g", GetSumw(), GetSumw2());
1469 for (Int_t dim = 0; dim < fNdimensions; ++dim) {
1470 Printf(" axis %d: Sum(w*x)=%g, Sum(w*x^2)=%g", dim, GetSumwx(dim), GetSumwx2(dim));
1471 }
1472 }
1473 }
1474
1475 if (optMem && InheritsFrom(THnSparse::Class())) {
1476 const THnSparse* hsparse = dynamic_cast<const THnSparse*>(this);
1477 Printf(" coordinates stored in %d chunks of %d entries\n %g of bins filled using %g of memory compared to an array",
1478 hsparse->GetNChunks(), hsparse->GetChunkSize(),
1479 hsparse->GetSparseFractionBins(), hsparse->GetSparseFractionMem());
1480 }
1481
1482 if (optContent) {
1483 Printf(" BIN CONTENT:");
1484 PrintEntries(0, -1, options);
1485 }
1486}
1487
1488
1489////////////////////////////////////////////////////////////////////////////////
1490/// Browse a THnSparse: create an entry (ROOT::THnSparseBrowsable) for each
1491/// dimension.
1492
1494{
1495 if (fBrowsables.IsEmpty()) {
1496 for (Int_t dim = 0; dim < fNdimensions; ++dim) {
1498 }
1500 }
1501
1502 for (Int_t dim = 0; dim < fNdimensions; ++dim) {
1503 b->Add(fBrowsables[dim]);
1504 }
1505}
1506
1507
1508
1509
1510/** \class ROOT::Internal::THnBaseBinIter
1511 Iterator over THnBase bins (internal implementation).
1512*/
1513
1514/// Destruct a bin iterator.
1515
1517 // Not much to do, but pin vtable
1518}
1519
1520
1521
1522/** \class THnIter
1523 Iterator over THnBase bins
1524*/
1525
1527
1529 // Destruct a bin iterator.
1530 delete fIter;
1531}
1532
1533
1534/** \class ROOT::Internal::THnBaseBrowsable
1535 TBrowser helper for THnBase.
1536*/
1537
1538
1540
1541////////////////////////////////////////////////////////////////////////////////
1542/// Construct a THnBaseBrowsable.
1543
1545fHist(hist), fAxis(axis), fProj(0)
1546{
1547 TString axisName = hist->GetAxis(axis)->GetName();
1548 if (axisName.IsNull()) {
1549 axisName = TString::Format("axis%d", axis);
1550 }
1551
1552 SetNameTitle(axisName,
1553 TString::Format("Projection on %s of %s", axisName.Data(),
1554 hist->IsA()->GetName()).Data());
1555}
1556
1557////////////////////////////////////////////////////////////////////////////////
1558/// Destruct a THnBaseBrowsable.
1559
1561{
1562 delete fProj;
1563}
1564
1565////////////////////////////////////////////////////////////////////////////////
1566/// Browse an axis of a THnBase, i.e. draw its projection.
1567
1569{
1570 if (!fProj) {
1571 fProj = fHist->Projection(fAxis);
1572 }
1573 fProj->Draw(b ? b->GetDrawOption() : "");
1574 gPad->Update();
1575}
1576
#define d(i)
Definition RSha256.hxx:102
#define b(i)
Definition RSha256.hxx:100
#define f(i)
Definition RSha256.hxx:104
#define c(i)
Definition RSha256.hxx:101
#define h(i)
Definition RSha256.hxx:106
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
int Int_t
Definition RtypesCore.h:45
constexpr Bool_t kFALSE
Definition RtypesCore.h:101
constexpr Ssiz_t kNPOS
Definition RtypesCore.h:124
long long Long64_t
Definition RtypesCore.h:80
constexpr Bool_t kTRUE
Definition RtypesCore.h:100
const char Option_t
Definition RtypesCore.h:66
#define ClassImp(name)
Definition Rtypes.h:377
winID h TVirtualViewer3D TVirtualGLPainter p
Option_t Option_t option
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char cname
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void value
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t type
char name[80]
Definition TGX11.cxx:110
float xmin
float xmax
#define R__THNBCASE(TAG)
R__EXTERN TRandom * gRandom
Definition TRandom.h:62
char * Form(const char *fmt,...)
Formats a string in a circular formatting buffer.
Definition TString.cxx:2467
void Printf(const char *fmt,...)
Formats a string in a circular formatting buffer and prints the string.
Definition TString.cxx:2481
#define gPad
class describing the range in the coordinates it supports multiple range in a coordinate.
Definition DataRange.h:35
void AddRange(unsigned int icoord, double xmin, double xmax)
add a range [xmin,xmax] for the new coordinate icoord Adding a range does not delete existing one,...
Definition DataRange.cxx:94
virtual ~THnBaseBinIter()
Destruct a bin iterator.
Definition THnBase.cxx:1516
TBrowser helper for THnBase.
Definition THnBase.h:297
~THnBaseBrowsable() override
Destruct a THnBaseBrowsable.
Definition THnBase.cxx:1560
THnBaseBrowsable(THnBase *hist, Int_t axis)
Construct a THnBaseBrowsable.
Definition THnBase.cxx:1544
void Browse(TBrowser *b) override
Browse an axis of a THnBase, i.e. draw its projection.
Definition THnBase.cxx:1568
Double_t At(Int_t i) const
Definition TArrayD.h:79
void Set(Int_t n) override
Set size of this array to n doubles.
Definition TArrayD.cxx:106
const Double_t * GetArray() const
Definition TArrayD.h:43
Abstract array base class.
Definition TArray.h:31
Int_t GetSize() const
Definition TArray.h:47
Class to manage histogram axis.
Definition TAxis.h:30
const char * GetTitle() const override
Returns title of object.
Definition TAxis.h:130
virtual Double_t GetBinCenter(Int_t bin) const
Return center of bin.
Definition TAxis.cxx:478
const TArrayD * GetXbins() const
Definition TAxis.h:131
Double_t GetXmax() const
Definition TAxis.h:135
@ kAxisRange
Definition TAxis.h:61
virtual Double_t GetBinLowEdge(Int_t bin) const
Return low edge of bin.
Definition TAxis.cxx:518
virtual void Set(Int_t nbins, Double_t xmin, Double_t xmax)
Initialize axis with fix bins.
Definition TAxis.cxx:759
Int_t GetLast() const
Return last bin on the axis i.e.
Definition TAxis.cxx:469
Double_t GetXmin() const
Definition TAxis.h:134
Int_t GetNbins() const
Definition TAxis.h:121
virtual void SetRange(Int_t first=0, Int_t last=0)
Set the viewing range for the axis using bin numbers.
Definition TAxis.cxx:952
virtual Double_t GetBinWidth(Int_t bin) const
Return bin width.
Definition TAxis.cxx:540
virtual Double_t GetBinUpEdge(Int_t bin) const
Return up edge of bin.
Definition TAxis.cxx:528
Int_t GetFirst() const
Return first bin on the axis i.e.
Definition TAxis.cxx:458
Using a TBrowser one can browse all ROOT objects.
Definition TBrowser.h:37
TClass instances represent classes, structs and namespaces in the ROOT type system.
Definition TClass.h:81
void * New(ENewType defConstructor=kClassNew, Bool_t quiet=kFALSE) const
Return a pointer to a newly allocated object of this class.
Definition TClass.cxx:4978
static TClass * GetClass(const char *name, Bool_t load=kTRUE, Bool_t silent=kFALSE)
Static method returning pointer to TClass of the specified class name.
Definition TClass.cxx:2968
Collection abstract base class.
Definition TCollection.h:65
virtual void SetOwner(Bool_t enable=kTRUE)
Set whether this collection is the owner (enable==true) of its content.
virtual Bool_t IsEmpty() const
virtual Int_t GetSize() const
Return the capacity of the collection, i.e.
1-Dim function class
Definition TF1.h:213
static void RejectPoint(Bool_t reject=kTRUE)
Static function to set the global flag to reject points the fgRejectPoint global flag is tested by al...
Definition TF1.cxx:3646
Provides an indirection to the TFitResult class and with a semantics identical to a TFitResult pointe...
1-D histogram with a double per channel (see TH1 documentation)}
Definition TH1.h:620
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:58
virtual Double_t GetEffectiveEntries() const
Number of effective entries of the histogram.
Definition TH1.cxx:4444
virtual void Rebuild(Option_t *option="")
Using the current bin info, recompute the arrays for contents and errors.
Definition TH1.cxx:7075
TAxis * GetZaxis()
Definition TH1.h:324
virtual void AddBinContent(Int_t bin)
Increment bin content by 1.
Definition TH1.cxx:1242
virtual Double_t GetBinError(Int_t bin) const
Return value of error associated to bin number bin.
Definition TH1.cxx:8929
static Int_t FitOptionsMake(Option_t *option, Foption_t &Foption)
Decode string choptin and fill fitOption structure.
Definition TH1.cxx:4616
TAxis * GetXaxis()
Definition TH1.h:322
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 abo...
Definition TH1.cxx:4938
virtual Int_t GetNcells() const
Definition TH1.h:298
virtual Int_t GetBin(Int_t binx, Int_t biny=0, Int_t binz=0) const
Return Global bin number corresponding to binx,y,z.
Definition TH1.cxx:4925
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-empt...
Definition TH1.cxx:9072
TAxis * GetYaxis()
Definition TH1.h:323
virtual Double_t GetEntries() const
Return the current number of entries.
Definition TH1.cxx:4419
virtual void ResetStats()
Reset the statistics including the number of entries and replace with values calculated from bin cont...
Definition TH1.cxx:7812
virtual Double_t GetBinContent(Int_t bin) const
Return content of bin number bin.
Definition TH1.cxx:5025
virtual void SetEntries(Double_t n)
Definition TH1.h:387
2-D histogram with a double per channel (see TH1 documentation)}
Definition TH2.h:300
3-D histogram with a double per channel (see TH1 documentation)}
Definition TH3.h:307
Multidimensional histogram base.
Definition THnBase.h:43
virtual void Sumw2()=0
void SetEntries(Double_t entries)
Definition THnBase.h:192
virtual void SetFilledBins(Long64_t)
Definition THnBase.h:105
void Browse(TBrowser *b) override
Browse a THnSparse: create an entry (ROOT::THnSparseBrowsable) for each dimension.
Definition THnBase.cxx:1493
virtual void InitStorage(Int_t *nbins, Int_t chunkSize)=0
Double_t GetBinError(const Int_t *idx) const
Definition THnBase.h:187
void SetBinError(const Int_t *idx, Double_t e)
Definition THnBase.h:189
enum THnBase::@84 fIntegralStatus
! status of integral
Double_t fEntries
Number of entries, spread over chunks.
Definition THnBase.h:48
Double_t GetSumw2() const
Definition THnBase.h:212
Bool_t IsInRange(Int_t *coord) const
Check whether bin coord is in range, as defined by TAxis::SetRange().
Definition THnBase.cxx:528
TFitResultPtr Fit(TF1 *f1, Option_t *option="", Option_t *goption="")
Fit a THnSparse with function f.
Definition THnBase.cxx:464
void AddBinContent(const Int_t *x, Double_t v=1.)
Definition THnBase.h:191
void Scale(Double_t c)
Scale contents and errors of this histogram by c: this = this * c It does not modify the histogram's ...
Definition THnBase.cxx:710
virtual void SetBinError2(Long64_t bin, Double_t e2)=0
TObjArray * GetListOfAxes()
Definition THnBase.h:128
TH1 * CreateHist(const char *name, const char *title, const TObjArray *axes, Bool_t keepTargetAxis) const
Create an empty histogram with name and title with a given set of axes.
Definition THnBase.cxx:243
Bool_t PrintBin(Long64_t idx, Int_t *coord, Option_t *options) const
Print one bin.
Definition THnBase.cxx:1352
Double_t GetSumwx(Int_t dim) const
Definition THnBase.h:213
Bool_t CheckConsistency(const THnBase *h, const char *tag) const
Consistency check on (some of) the parameters of two histograms (for operations).
Definition THnBase.cxx:1098
Int_t GetNdimensions() const
Definition THnBase.h:140
~THnBase() override
Destruct a THnBase.
Definition THnBase.cxx:165
static THnBase * CreateHnAny(const char *name, const char *title, const TH1 *h1, Bool_t sparse, Int_t chunkSize=1024 *16)
Create a THn / THnSparse object from a histogram deriving from TH1.
Definition THnBase.cxx:296
TArrayD fTsumwx2
Total sum of weight*X*X for each dimension.
Definition THnBase.h:52
void ResetBase(Option_t *option="")
Clear the histogram.
Definition THnBase.cxx:1263
virtual Long64_t GetNbins() const =0
THnBase * RebinBase(Int_t group) const
Combine the content of "group" neighboring bins into a new bin and return the resulting THnBase.
Definition THnBase.cxx:1168
TClass * IsA() const override
Definition THnBase.h:289
TObjArray fAxes
Axes of the histogram.
Definition THnBase.h:46
void SetBinEdges(Int_t idim, const Double_t *bins)
Set the axis # of bins and bin limits on dimension idim.
Definition THnBase.cxx:1116
Bool_t GetCalculateErrors() const
Definition THnBase.h:141
void AddInternal(const THnBase *h, Double_t c, Bool_t rebinned)
Add() implementation for both rebinned histograms and those with identical binning.
Definition THnBase.cxx:734
void SetTitle(const char *title) override
Change (i.e.
Definition THnBase.cxx:1132
void PrintEntries(Long64_t from=0, Long64_t howmany=-1, Option_t *options=nullptr) const
Print "howmany" entries starting at "from".
Definition THnBase.cxx:1396
Double_t GetBinContent(const Int_t *idx) const
Definition THnBase.h:195
void Add(const THnBase *h, Double_t c=1.)
Add contents of h scaled by c to this histogram: this = this + c * h Note that if h has Sumw2 set,...
Definition THnBase.cxx:807
Double_t fTsumw2
Total sum of weights squared; -1 if no errors are calculated.
Definition THnBase.h:50
virtual void Reserve(Long64_t)
Definition THnBase.h:104
THnBase()
Definition THnBase.h:61
Double_t GetEntries() const
Definition THnBase.h:138
virtual Double_t GetBinError2(Long64_t linidx) const =0
virtual void Reset(Option_t *option="")=0
void SetBinContent(const Int_t *idx, Double_t v)
Definition THnBase.h:205
void Multiply(const THnBase *h)
Multiply this histogram by histogram h this = this * h Note that if h has Sumw2 set,...
Definition THnBase.cxx:868
void Divide(const THnBase *h)
Divide this histogram by h this = this/(h) Note that if h has Sumw2 set, Sumw2 is automatically calle...
Definition THnBase.cxx:959
void Print(Option_t *option="") const override
Print a THnBase.
Definition THnBase.cxx:1446
Double_t GetSumwx2(Int_t dim) const
Definition THnBase.h:214
TObjArray fBrowsables
! Browser-helpers for each axis
Definition THnBase.h:47
Double_t GetSumw() const
Definition THnBase.h:211
TObject * ProjectionAny(Int_t ndim, const Int_t *dim, Bool_t wantNDim, Option_t *option="") const
Project all bins into a ndim-dimensional THn / THnSparse (whatever *this is) or if (ndim < 4 and !...
Definition THnBase.cxx:554
void GetRandom(Double_t *rand, Bool_t subBinRandom=kTRUE)
Generate an n-dimensional random tuple based on the histogrammed distribution.
Definition THnBase.cxx:493
TAxis * GetAxis(Int_t dim) const
Definition THnBase.h:130
THnBase * CloneEmpty(const char *name, const char *title, const TObjArray *axes, Bool_t keepTargetAxis) const
Create a new THnBase object that is of the same type as *this, but with dimensions and bins given by ...
Definition THnBase.cxx:177
Long64_t Merge(TCollection *list)
Merge this with a list of THnBase's.
Definition THnBase.cxx:833
Double_t fTsumw
Total sum of weights.
Definition THnBase.h:49
virtual Long64_t GetBin(const Int_t *idx) const =0
Double_t ComputeIntegral()
Calculate the integral of the histogram.
Definition THnBase.cxx:1277
@ kNoInt
Definition THnBase.h:55
@ kValidInt
Definition THnBase.h:56
virtual void AddBinError2(Long64_t bin, Double_t e2)=0
void RebinnedAdd(const THnBase *h, Double_t c=1.)
Add contents of h scaled by c to this histogram: this = this + c * h Note that if h has Sumw2 set,...
Definition THnBase.cxx:823
Int_t fNdimensions
Number of dimensions.
Definition THnBase.h:45
void Init(const char *name, const char *title, const TObjArray *axes, Bool_t keepTargetAxis, Int_t chunkSize=1024 *16)
Initialize axes and name.
Definition THnBase.cxx:193
THnBase & operator=(const THnBase &other)
Definition THnBase.cxx:99
std::vector< Double_t > fIntegral
! vector with bin weight sums
Definition THnBase.h:53
TArrayD fTsumwx
Total sum of weight*X for each dimension.
Definition THnBase.h:51
Iterator over THnBase bins.
Definition THnBase.h:330
Bool_t HaveSkippedBin() const
Definition THnBase.h:345
Int_t GetCoord(Int_t dim) const
Definition THnBase.h:344
~THnIter() override
Definition THnBase.cxx:1528
Long64_t Next(Int_t *coord=nullptr)
Return the next bin's index.
Definition THnBase.h:340
static TClass * Class()
Efficient multidimensional histogram.
Definition THnSparse.h:36
Double_t GetSparseFractionBins() const
Return the amount of filled bins over all bins.
Double_t GetSparseFractionMem() const
Return the amount of used memory over memory that would be used by a non-sparse n-dimensional histogr...
Int_t GetNChunks() const
Definition THnSparse.h:88
Int_t GetChunkSize() const
Definition THnSparse.h:87
static TClass * Class()
static TClass * Class()
static TClass * Class()
void Reset()
The TNamed class is the base class for all named ROOT classes.
Definition TNamed.h:29
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
Definition TNamed.cxx:164
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
const char * GetTitle() const override
Returns title of object.
Definition TNamed.h:48
TString fTitle
Definition TNamed.h:33
void Clear(Option_t *option="") override
Set name and title to empty strings ("").
Definition TNamed.cxx:64
virtual void SetName(const char *name)
Set the name of the TNamed.
Definition TNamed.cxx:140
TNamed & operator=(const TNamed &rhs)
TNamed assignment operator.
Definition TNamed.cxx:51
virtual void SetNameTitle(const char *name, const char *title)
Set all the TNamed parameters (name and title).
Definition TNamed.cxx:154
An array of TObjects.
Definition TObjArray.h:31
Int_t GetEntriesFast() const
Definition TObjArray.h:58
void AddAt(TObject *obj, Int_t idx) override
Add object at position ids.
virtual void AddAtAndExpand(TObject *obj, Int_t idx)
Add object at position idx.
TObject * At(Int_t idx) const override
Definition TObjArray.h:164
Bool_t IsEmpty() const override
Definition TObjArray.h:65
Mother of all ROOT objects.
Definition TObject.h:41
virtual const char * GetName() const
Returns name of object.
Definition TObject.cxx:439
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
Definition TObject.h:201
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition TObject.cxx:956
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition TObject.cxx:774
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition TObject.cxx:525
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:970
void ResetBit(UInt_t f)
Definition TObject.h:200
Double_t Rndm() override
Machine independent random number generator.
Definition TRandom.cxx:552
Basic string class.
Definition TString.h:139
Ssiz_t Length() const
Definition TString.h:421
void ToLower()
Change string to lower-case.
Definition TString.cxx:1170
TString & Insert(Ssiz_t pos, const char *s)
Definition TString.h:661
Ssiz_t First(char c) const
Find first occurrence of a character c.
Definition TString.cxx:531
const char * Data() const
Definition TString.h:380
TString & ReplaceAll(const TString &s1, const TString &s2)
Definition TString.h:704
Bool_t IsNull() const
Definition TString.h:418
TString & Remove(Ssiz_t pos)
Definition TString.h:685
static TString Format(const char *fmt,...)
Static method which formats a string using a printf style format descriptor and return a TString.
Definition TString.cxx:2356
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition TString.h:636
Ssiz_t Index(const char *pat, Ssiz_t i=0, ECaseCompare cmp=kExact) const
Definition TString.h:651
return c1
Definition legend1.C:41
Double_t x[n]
Definition legend1.C:17
TH1F * h1
Definition legend1.C:5
return c2
Definition legend2.C:14
#define F(x, y, z)
#define I(x, y, z)
TFitResultPtr FitObject(TH1 *h1, TF1 *f1, Foption_t &option, const ROOT::Math::MinimizerOptions &moption, const char *goption, ROOT::Fit::DataRange &range)
fitting function for a TH1 (called from TH1::Fit)
Definition HFitImpl.cxx:965
Double_t Floor(Double_t x)
Rounds x downward, returning the largest integral value that is not greater than x.
Definition TMath.h:678
Double_t Sqrt(Double_t x)
Returns the square root of x.
Definition TMath.h:660
Int_t CeilNint(Double_t x)
Returns the nearest integer of TMath::Ceil(x).
Definition TMath.h:672
Long64_t BinarySearch(Long64_t n, const T *array, T value)
Binary search in an array of n values to locate value.
Definition TMathBase.h:347
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
Definition TMathBase.h:123
int Like
Definition Foption.h:34
int Nostore
Definition Foption.h:41
int Chi2
Definition Foption.h:32