51 THnBase::THnBase(
const char* name,
const char* title, Int_t dim,
52 const Int_t* nbins,
const Double_t* xmin,
const Double_t* xmax):
53 TNamed(name, title), fNdimensions(dim), fAxes(dim), fBrowsables(dim),
54 fEntries(0), fTsumw(0), fTsumw2(-1.), fTsumwx(dim), fTsumwx2(dim),
55 fIntegral(0), fIntegralStatus(kNoInt)
57 for (Int_t i = 0; i < fNdimensions; ++i) {
58 TAxis* axis =
new TAxis(nbins[i], xmin ? xmin[i] : 0., xmax ? xmax[i] : 1.);
59 axis->SetName(TString::Format(
"axis%d", i));
60 fAxes.AddAtAndExpand(axis, i);
70 if (fIntegralStatus != kNoInt)
delete [] fIntegral;
80 THnBase* THnBase::CloneEmpty(
const char* name,
const char* title,
81 const TObjArray* axes, Bool_t keepTargetAxis)
const
83 THnBase* ret = (THnBase*)IsA()->New();
84 Int_t chunkSize = 1024 * 16;
85 if (InheritsFrom(THnSparse::Class())) {
86 chunkSize = ((
const THnSparse*)
this)->GetChunkSize();
88 ret->Init(name, title, axes, keepTargetAxis, chunkSize);
96 void THnBase::Init(
const char* name,
const char* title,
97 const TObjArray* axes, Bool_t keepTargetAxis,
100 SetNameTitle(name, title);
103 const TAxis* axis = 0;
105 Int_t *nbins =
new Int_t[axes->GetEntriesFast()];
106 while ((axis = (TAxis*)iAxis())) {
107 TAxis* reqaxis =
new TAxis(*axis);
108 if (!keepTargetAxis && axis->TestBit(TAxis::kAxisRange)) {
109 Int_t binFirst = axis->GetFirst();
113 Int_t binLast = axis->GetLast();
115 if (binLast > axis->GetNbins())
116 binLast = axis->GetNbins();
117 Int_t nBins = binLast - binFirst + 1;
118 if (axis->GetXbins()->GetSize()) {
120 reqaxis->Set(nBins, axis->GetXbins()->GetArray() + binFirst - 1);
123 reqaxis->Set(nBins, axis->GetBinLowEdge(binFirst), axis->GetBinUpEdge(binLast));
125 reqaxis->ResetBit(TAxis::kAxisRange);
128 nbins[pos] = reqaxis->GetNbins();
129 fAxes.AddAtAndExpand(
new TAxis(*reqaxis), pos++);
133 fNdimensions = axes->GetEntriesFast();
134 InitStorage(nbins, chunkSize);
144 TH1* THnBase::CreateHist(
const char* name,
const char* title,
145 const TObjArray* axes,
146 Bool_t keepTargetAxis )
const {
147 const int ndim = axes->GetSize();
152 hist =
new TH1D(name, title, 1, 0., 1.);
154 hist =
new TH2D(name, title, 1, 0., 1., 1, 0., 1.);
156 hist =
new TH3D(name, title, 1, 0., 1., 1, 0., 1., 1, 0., 1.);
158 Error(
"CreateHist",
"Cannot create histogram %s with %d dimensions!", name, ndim);
162 TAxis* hax[3] = {hist->GetXaxis(), hist->GetYaxis(), hist->GetZaxis()};
163 for (Int_t d = 0; d < ndim; ++d) {
164 TAxis* reqaxis = (TAxis*)(*axes)[d];
165 hax[d]->SetTitle(reqaxis->GetTitle());
166 if (!keepTargetAxis && reqaxis->TestBit(TAxis::kAxisRange)) {
168 Int_t binFirst = std::max(reqaxis->GetFirst(),1);
169 Int_t binLast = std::min(reqaxis->GetLast(), reqaxis->GetNbins() );
170 Int_t nBins = binLast - binFirst + 1;
171 if (reqaxis->GetXbins()->GetSize()) {
173 hax[d]->Set(nBins, reqaxis->GetXbins()->GetArray() + binFirst - 1);
176 hax[d]->Set(nBins, reqaxis->GetBinLowEdge(binFirst), reqaxis->GetBinUpEdge(binLast));
179 if (reqaxis->GetXbins()->GetSize()) {
181 hax[d]->Set(reqaxis->GetNbins(), reqaxis->GetXbins()->GetArray());
184 hax[d]->Set(reqaxis->GetNbins(), reqaxis->GetXmin(), reqaxis->GetXmax());
197 THnBase* THnBase::CreateHnAny(
const char* name,
const char* title,
198 const TH1* h, Bool_t sparse, Int_t chunkSize)
201 int ndim = h->GetDimension();
204 int nbins[3] = {0,0,0};
205 double minRange[3] = {0.,0.,0.};
206 double maxRange[3] = {0.,0.,0.};
207 const TAxis* axis[3] = { h->GetXaxis(), h->GetYaxis(), h->GetZaxis() };
208 for (
int i = 0; i < ndim; ++i) {
209 nbins[i] = axis[i]->GetNbins();
210 minRange[i] = axis[i]->GetXmin();
211 maxRange[i] = axis[i]->GetXmax();
219 const char* cname( h->ClassName() );
220 if (cname[0] ==
'T' && cname[1] ==
'H'
221 && cname[2] >=
'1' && cname[2] <=
'3' && cname[4] == 0) {
223 #define R__THNBCASE(TAG) \
225 s = new _NAME2_(THnSparse,TAG)(name, title, ndim, nbins, \
226 minRange, maxRange, chunkSize); \
228 s = new _NAME2_(THn,TAG)(name, title, ndim, nbins, \
229 minRange, maxRange); \
234 case 'F': R__THNBCASE(F);
235 case 'D': R__THNBCASE(D);
236 case 'I': R__THNBCASE(I);
237 case 'S': R__THNBCASE(S);
238 case 'C': R__THNBCASE(C);
243 ::Warning(
"THnSparse::CreateHnAny",
"Unknown Type of Histogram");
247 for (
int i = 0; i < ndim; ++i) {
248 s->GetAxis(i)->SetTitle(axis[i]->GetTitle());
252 const TArray *array =
dynamic_cast<const TArray*
>(h);
254 ::Warning(
"THnSparse::CreateHnAny",
"Unknown Type of Histogram");
267 THnBase* THnBase::CreateHnAny(
const char* name,
const char* title,
268 const THnBase* hn, Bool_t sparse,
272 if (hn->InheritsFrom(THnSparse::Class())) {
273 if (sparse) type = hn->IsA();
276 if (hn->InheritsFrom(THnSparseD::Class())) bintype =
'D';
277 else if (hn->InheritsFrom(THnSparseF::Class())) bintype =
'F';
278 else if (hn->InheritsFrom(THnSparseL::Class())) bintype =
'L';
279 else if (hn->InheritsFrom(THnSparseI::Class())) bintype =
'I';
280 else if (hn->InheritsFrom(THnSparseS::Class())) bintype =
'S';
281 else if (hn->InheritsFrom(THnSparseC::Class())) bintype =
'C';
283 hn->Error(
"CreateHnAny",
"Type %s not implemented; please inform the ROOT team!",
284 hn->IsA()->GetName());
287 type = TClass::GetClass(TString::Format(
"THn%c", bintype));
289 }
else if (hn->InheritsFrom(THn::Class())) {
290 if (!sparse) type = hn->IsA();
293 if (hn->InheritsFrom(THnD::Class())) bintype =
'D';
294 else if (hn->InheritsFrom(THnF::Class())) bintype =
'F';
295 else if (hn->InheritsFrom(THnC::Class())) bintype =
'C';
296 else if (hn->InheritsFrom(THnS::Class())) bintype =
'S';
297 else if (hn->InheritsFrom(THnI::Class())) bintype =
'I';
298 else if (hn->InheritsFrom(THnL::Class())) bintype =
'L';
299 else if (hn->InheritsFrom(THnL64::Class())) {
300 hn->Error(
"CreateHnAny",
"Type THnSparse with Long64_t bins is not available!");
304 type = TClass::GetClass(TString::Format(
"THnSparse%c", bintype));
308 hn->Error(
"CreateHnAny",
"Unhandled type %s, not deriving from THn nor THnSparse!",
309 hn->IsA()->GetName());
313 hn->Error(
"CreateHnAny",
"Unhandled type %s, please inform the ROOT team!",
314 hn->IsA()->GetName());
318 THnBase* ret = (THnBase*)type->New();
319 ret->Init(name, title, hn->GetListOfAxes(),
331 void THnBase::Add(
const TH1* hist, Double_t c )
333 Long64_t nbins = hist->GetNcells();
335 for (
int i = 0; i < nbins; ++i) {
336 double value = hist->GetBinContent(i);
337 double error = hist->GetBinError(i);
338 if (!value && !error)
continue;
339 hist->GetBinXYZ(i, x[0], x[1], x[2]);
340 SetBinContent(x, value * c);
341 SetBinError(x, error * c);
365 TFitResultPtr THnBase::Fit(TF1 *f ,Option_t *option ,Option_t *goption)
370 if (!TH1::FitOptionsMake(option,fitOption))
return 0;
374 fitOption.Nostore =
true;
376 if (!fitOption.Chi2) fitOption.Like =
true;
378 ROOT::Fit::DataRange range(GetNdimensions());
379 for (
int i = 0; i < GetNdimensions(); ++i ) {
380 TAxis *axis = GetAxis(i);
381 range.AddRange(i, axis->GetXmin(), axis->GetXmax());
383 ROOT::Math::MinimizerOptions minOption;
385 return ROOT::Fit::FitObject(
this, f , fitOption , minOption, goption, range);
394 void THnBase::GetRandom(Double_t *rand, Bool_t subBinRandom )
397 if (fIntegralStatus != kValidInt)
401 Double_t p = gRandom->Rndm();
402 Long64_t idx = TMath::BinarySearch(GetNbins() + 1, fIntegral, p);
403 const Int_t nStaticBins = 40;
404 Int_t bin[nStaticBins];
406 if (GetNdimensions() > nStaticBins) {
407 pBin =
new Int_t[GetNdimensions()];
409 GetBinContent(idx, pBin);
412 for (Int_t i = 0; i < fNdimensions; i++) {
413 rand[i] = GetAxis(i)->GetBinCenter(pBin[i]);
417 rand[i] += (gRandom->Rndm() - 0.5) * GetAxis(i)->GetBinWidth(pBin[i]);
429 Bool_t THnBase::IsInRange(Int_t *coord)
const
433 for (Int_t i = 0; i < fNdimensions; ++i) {
434 TAxis *axis = GetAxis(i);
435 if (!axis->TestBit(TAxis::kAxisRange))
continue;
436 min = axis->GetFirst();
437 max = axis->GetLast();
438 if (coord[i] < min || coord[i] > max)
455 TObject* THnBase::ProjectionAny(Int_t ndim,
const Int_t* dim,
457 Option_t* option )
const
459 TString name(GetName());
462 for (Int_t d = 0; d < ndim; ++d) {
467 TString title(GetTitle());
468 Ssiz_t posInsert = title.First(
';');
469 if (posInsert == kNPOS) {
470 title +=
" projection ";
471 for (Int_t d = 0; d < ndim; ++d)
472 title += GetAxis(dim[d])->GetTitle();
474 for (Int_t d = ndim - 1; d >= 0; --d) {
475 title.Insert(posInsert, GetAxis(d)->GetTitle());
477 title.Insert(posInsert,
", ");
479 title.Insert(posInsert,
" projection ");
482 TObjArray newaxes(ndim);
483 for (Int_t d = 0; d < ndim; ++d) {
484 newaxes.AddAt(GetAxis(dim[d]),d);
491 Bool_t* hadRange = 0;
492 Bool_t ignoreTargetRange = (option && (strchr(option,
'A') || strchr(option,
'a')));
493 Bool_t keepTargetAxis = ignoreTargetRange || (option && (strchr(option,
'O') || strchr(option,
'o')));
494 if (ignoreTargetRange) {
495 hadRange =
new Bool_t[ndim];
496 for (Int_t d = 0; d < ndim; ++d){
497 TAxis *axis = GetAxis(dim[d]);
498 hadRange[d] = axis->TestBit(TAxis::kAxisRange);
499 axis->SetBit(TAxis::kAxisRange, kFALSE);
504 ret = hn = CloneEmpty(name, title, &newaxes, keepTargetAxis);
506 ret = hist = CreateHist(name, title, &newaxes, keepTargetAxis);
508 if (keepTargetAxis) {
511 for (Int_t d = 0; d < ndim; ++d) {
512 hn->GetAxis(d)->SetRange(0, 0);
515 hist->GetXaxis()->SetRange(0, 0);
516 hist->GetYaxis()->SetRange(0, 0);
517 hist->GetZaxis()->SetRange(0, 0);
521 Bool_t haveErrors = GetCalculateErrors();
522 Bool_t wantErrors = haveErrors || (option && (strchr(option,
'E') || strchr(option,
'e')));
524 Int_t* bins =
new Int_t[ndim];
525 Long64_t myLinBin = 0;
527 THnIter iter(
this, kTRUE );
529 while ((myLinBin = iter.Next()) >= 0) {
530 Double_t v = GetBinContent(myLinBin);
532 for (Int_t d = 0; d < ndim; ++d) {
533 bins[d] = iter.GetCoord(dim[d]);
534 if (!keepTargetAxis && GetAxis(dim[d])->TestBit(TAxis::kAxisRange)) {
535 Int_t binOffset = GetAxis(dim[d])->GetFirst();
537 if (binOffset > 0) --binOffset;
538 bins[d] -= binOffset;
542 Long64_t targetLinBin = -1;
544 if (ndim == 1) targetLinBin = bins[0];
545 else if (ndim == 2) targetLinBin = hist->GetBin(bins[0], bins[1]);
546 else if (ndim == 3) targetLinBin = hist->GetBin(bins[0], bins[1], bins[2]);
548 targetLinBin = hn->GetBin(bins, kTRUE );
554 err2 = GetBinError2(myLinBin);
559 hn->AddBinError2(targetLinBin, err2);
561 Double_t preverr = hist->GetBinError(targetLinBin);
562 hist->SetBinError(targetLinBin, TMath::Sqrt(preverr * preverr + err2));
568 hn->AddBinContent(targetLinBin, v);
570 hist->AddBinContent(targetLinBin, v);
576 hn->SetEntries(fEntries);
578 if (!iter.HaveSkippedBin()) {
579 hist->SetEntries(fEntries);
586 Double_t entries = hist->GetEffectiveEntries();
589 entries = TMath::Floor(entries + 0.5);
591 hist->SetEntries(entries);
597 for (Int_t d = 0; d < ndim; ++d)
598 GetAxis(dim[d])->SetBit(TAxis::kAxisRange, hadRange[d]);
611 void THnBase::Scale(Double_t c)
614 Double_t nEntries = GetEntries();
616 Bool_t haveErrors = GetCalculateErrors();
619 while ((i = iter.Next()) >= 0) {
621 Double_t v = GetBinContent(i);
622 SetBinContent(i, c * v);
624 Double_t err2 = GetBinError2(i);
625 SetBinError2(i, c * c * err2);
628 SetEntries(nEntries);
635 void THnBase::AddInternal(
const THnBase* h, Double_t c, Bool_t rebinned)
637 if (fNdimensions != h->GetNdimensions()) {
638 Warning(
"RebinnedAdd",
"Different number of dimensions, cannot carry out operation on the histograms");
643 if (!GetCalculateErrors() && h->GetCalculateErrors())
645 Bool_t haveErrors = GetCalculateErrors();
649 x =
new Double_t[fNdimensions];
651 Int_t* coord =
new Int_t[fNdimensions];
654 Long64_t numTargetBins = GetNbins() + h->GetNbins();
655 Reserve(numTargetBins);
660 while ((i = iter.Next(coord)) >= 0) {
662 Double_t v = h->GetBinContent(i);
664 Long64_t mybinidx = -1;
667 for (Int_t j = 0; j < fNdimensions; ++j)
668 x[j] = h->GetAxis(j)->GetBinCenter(coord[j]);
670 mybinidx = GetBin(x, kTRUE );
672 mybinidx = GetBin(coord, kTRUE );
676 Double_t err2 = h->GetBinError2(i) * c * c;
677 AddBinError2(mybinidx, err2);
680 AddBinContent(mybinidx, c * v);
686 Double_t nEntries = GetEntries() + c * h->GetEntries();
687 SetEntries(nEntries);
696 void THnBase::Add(
const THnBase* h, Double_t c)
699 if (!CheckConsistency(h,
"Add"))
return;
701 AddInternal(h, c, kFALSE);
712 void THnBase::RebinnedAdd(
const THnBase* h, Double_t c)
714 AddInternal(h, c, kTRUE);
722 Long64_t THnBase::Merge(TCollection* list)
725 if (list->IsEmpty())
return (Long64_t)GetEntries();
727 Long64_t sumNbins = GetNbins();
729 const TObject* addMeObj = 0;
730 while ((addMeObj = iter())) {
731 const THnBase* addMe =
dynamic_cast<const THnBase*
>(addMeObj);
733 sumNbins += addMe->GetNbins();
739 while ((addMeObj = iter())) {
740 const THnBase* addMe =
dynamic_cast<const THnBase*
>(addMeObj);
742 Error(
"Merge",
"Object named %s is not THnBase! Skipping it.",
743 addMeObj->GetName());
747 return (Long64_t)GetEntries();
757 void THnBase::Multiply(
const THnBase* h)
760 if(!CheckConsistency(h,
"Multiply"))
return;
763 Bool_t wantErrors = kFALSE;
764 if (GetCalculateErrors() || h->GetCalculateErrors())
767 if (wantErrors) Sumw2();
769 Double_t nEntries = GetEntries();
771 Int_t* coord =
new Int_t[fNdimensions];
775 while ((i = iter.Next(coord)) >= 0) {
777 Double_t v1 = GetBinContent(i);
779 Long64_t idxh = h->GetBin(coord);
781 if (idxh >= 0) v2 = h->GetBinContent(idxh);
782 SetBinContent(i, v1 * v2);
784 Double_t err1 = GetBinError(i) * v2;
786 if (idxh >= 0) err2 = h->GetBinError(idxh) * v1;
787 SetBinError(i, TMath::Sqrt((err2 * err2 + err1 * err1)));
790 SetEntries(nEntries);
805 void THnBase::Multiply(TF1* f, Double_t c)
807 Int_t* coord =
new Int_t[fNdimensions];
808 Double_t* x =
new Double_t[fNdimensions];
810 Bool_t wantErrors = GetCalculateErrors();
811 if (wantErrors) Sumw2();
816 while ((i = iter.Next(coord)) >= 0) {
817 Double_t value = GetBinContent(i);
820 for (Int_t j = 0; j < fNdimensions; ++j)
821 x[j] = GetAxis(j)->GetBinCenter(coord[j]);
825 TF1::RejectPoint(kFALSE);
828 Double_t fvalue = f->EvalPar(x, NULL);
830 SetBinContent(i, c * fvalue * value);
832 Double_t error = GetBinError(i);
833 SetBinError(i, c * fvalue * error);
848 void THnBase::Divide(
const THnBase *h)
851 if (!CheckConsistency(h,
"Divide"))
return;
854 Bool_t wantErrors=GetCalculateErrors();
855 if (!GetCalculateErrors() && h->GetCalculateErrors())
859 Double_t nEntries = fEntries;
861 if (wantErrors) Sumw2();
862 Bool_t didWarn = kFALSE;
865 Int_t* coord =
new Int_t[fNdimensions];
869 while ((i = iter.Next(coord)) >= 0) {
871 Double_t v1 = GetBinContent(i);
873 Long64_t hbin = h->GetBin(coord);
874 Double_t v2 = h->GetBinContent(hbin);
879 Warning(
"Divide(h)",
"Histogram h has empty bins - division by zero! Setting bin to 0.");
883 SetBinContent(i, v1 / v2);
885 Double_t err1 = GetBinError(i) * v2;
886 Double_t err2 = h->GetBinError(hbin) * v1;
887 Double_t b22 = v2 * v2;
888 Double_t err = (err1 * err1 + err2 * err2) / (b22 * b22);
889 SetBinError2(i, err);
893 SetEntries(nEntries);
905 void THnBase::Divide(
const THnBase *h1,
const THnBase *h2, Double_t c1, Double_t c2, Option_t *option)
908 TString opt = option;
910 Bool_t binomial = kFALSE;
911 if (opt.Contains(
"b")) binomial = kTRUE;
914 if (!CheckConsistency(h1,
"Divide") || !CheckConsistency(h2,
"Divide"))
return;
916 Error(
"Divide",
"Coefficient of dividing histogram cannot be zero");
923 if (!GetCalculateErrors() && (h1->GetCalculateErrors()|| h2->GetCalculateErrors() != 0))
927 Long64_t nFilledBins=0;
931 Int_t* coord =
new Int_t[fNdimensions];
932 memset(coord, 0,
sizeof(Int_t) * fNdimensions);
933 Bool_t didWarn = kFALSE;
938 while ((i = iter.Next(coord)) >= 0) {
940 Double_t v1 = h1->GetBinContent(i);
942 Long64_t h2bin = h2->GetBin(coord);
943 Double_t v2 = h2->GetBinContent(h2bin);
948 Warning(
"Divide(h1, h2)",
"Histogram h2 has empty bins - division by zero! Setting bin to 0.");
953 Long64_t myBin = GetBin(coord);
954 SetBinContent(myBin, c1 * v1 / c2 / v2);
955 if(GetCalculateErrors()){
956 Double_t err1 = h1->GetBinError(i);
957 Double_t err2 = h2->GetBinError(h2bin);
961 Double_t w = v1 / v2;
963 errSq = TMath::Abs( ( (1. - 2.*w) * err1 * err1 + err2 * err2 ) / (v2 * v2) );
968 Double_t b22 = v2 * v2 * c2;
971 errSq = c1 * c2 * (err1 * err1 + err2 * err2) / (b22 * b22);
973 SetBinError2(myBin, errSq);
978 SetFilledBins(nFilledBins);
981 SetEntries(h1->GetEntries());
987 Bool_t THnBase::CheckConsistency(
const THnBase *h,
const char *tag)
const
989 if (fNdimensions != h->GetNdimensions()) {
990 Warning(tag,
"Different number of dimensions, cannot carry out operation on the histograms");
993 for (Int_t dim = 0; dim < fNdimensions; dim++){
994 if (GetAxis(dim)->GetNbins() != h->GetAxis(dim)->GetNbins()) {
995 Warning(tag,
"Different number of bins on axis %i, cannot carry out operation on the histograms", dim);
1005 void THnBase::SetBinEdges(Int_t idim,
const Double_t* bins)
1007 TAxis* axis = (TAxis*) fAxes[idim];
1008 axis->Set(axis->GetNbins(), bins);
1021 void THnBase::SetTitle(
const char *title)
1024 fTitle.ReplaceAll(
"#;",2,
"#semicolon",10);
1026 Int_t endHistTitle = fTitle.First(
';');
1027 if (endHistTitle >= 0) {
1029 Int_t posTitle = endHistTitle + 1;
1030 Int_t lenTitle = fTitle.Length();
1032 while (posTitle > 0 && posTitle < lenTitle && dim < fNdimensions){
1033 Int_t endTitle = fTitle.Index(
";", posTitle);
1034 TString axisTitle = fTitle(posTitle, endTitle - posTitle);
1035 axisTitle.ReplaceAll(
"#semicolon", 10,
";", 1);
1036 GetAxis(dim)->SetTitle(axisTitle);
1039 posTitle = endTitle + 1;
1044 fTitle.Remove(endHistTitle, lenTitle - endHistTitle);
1047 fTitle.ReplaceAll(
"#semicolon", 10,
";", 1);
1057 THnBase* THnBase::RebinBase(Int_t group)
const
1059 Int_t* ngroup =
new Int_t[GetNdimensions()];
1060 for (Int_t d = 0; d < GetNdimensions(); ++d)
1062 THnBase* ret = RebinBase(ngroup);
1073 THnBase* THnBase::RebinBase(
const Int_t* group)
const
1075 Int_t ndim = GetNdimensions();
1076 TString name(GetName());
1077 for (Int_t d = 0; d < ndim; ++d)
1078 name += Form(
"_%d", group[d]);
1081 TString title(GetTitle());
1082 Ssiz_t posInsert = title.First(
';');
1083 if (posInsert == kNPOS) {
1085 for (Int_t d = 0; d < ndim; ++d)
1086 title += Form(
"{%d}", group[d]);
1088 for (Int_t d = ndim - 1; d >= 0; --d)
1089 title.Insert(posInsert, Form(
"{%d}", group[d]));
1090 title.Insert(posInsert,
" rebin ");
1093 TObjArray newaxes(ndim);
1095 for (Int_t d = 0; d < ndim; ++d) {
1096 newaxes.AddAt(
new TAxis(*GetAxis(d) ),d);
1098 TAxis* newaxis = (TAxis*) newaxes.At(d);
1099 Int_t newbins = (newaxis->GetNbins() + group[d] - 1) / group[d];
1100 if (newaxis->GetXbins() && newaxis->GetXbins()->GetSize()) {
1102 Double_t *edges =
new Double_t[newbins + 1];
1103 for (Int_t i = 0; i < newbins + 1; ++i)
1104 if (group[d] * i <= newaxis->GetNbins())
1105 edges[i] = newaxis->GetXbins()->At(group[d] * i);
1106 else edges[i] = newaxis->GetXmax();
1107 newaxis->Set(newbins, edges);
1110 newaxis->Set(newbins, newaxis->GetXmin(), newaxis->GetXmax());
1115 THnBase* h = CloneEmpty(name.Data(), title.Data(), &newaxes, kTRUE);
1116 Bool_t haveErrors = GetCalculateErrors();
1117 Bool_t wantErrors = haveErrors;
1119 Int_t* bins =
new Int_t[ndim];
1120 Int_t* coord =
new Int_t[fNdimensions];
1124 while ((i = iter.Next(coord)) >= 0) {
1125 Double_t v = GetBinContent(i);
1126 for (Int_t d = 0; d < ndim; ++d) {
1127 bins[d] = TMath::CeilNint( (
double) coord[d]/group[d] );
1129 Long64_t idxh = h->GetBin(bins, kTRUE );
1133 h->AddBinError2(idxh, GetBinError2(i));
1137 h->AddBinContent(idxh, v);
1143 h->SetEntries(fEntries);
1152 void THnBase::ResetBase(Option_t * )
1157 if (fIntegralStatus != kNoInt) {
1158 delete [] fIntegral;
1159 fIntegralStatus = kNoInt;
1166 Double_t THnBase::ComputeIntegral()
1169 if (fIntegralStatus != kNoInt) {
1170 delete [] fIntegral;
1171 fIntegralStatus = kNoInt;
1175 if (GetNbins() == 0) {
1176 Error(
"ComputeIntegral",
"The histogram must have at least one bin.");
1181 fIntegral =
new Double_t [GetNbins() + 1];
1185 Int_t* coord =
new Int_t[fNdimensions];
1188 while ((i = iter.Next(coord)) >= 0) {
1189 Double_t v = GetBinContent(i);
1192 bool regularBin =
true;
1193 for (Int_t dim = 0; dim < fNdimensions; dim++) {
1194 if (coord[dim] < 1 || coord[dim] > GetAxis(dim)->GetNbins()) {
1201 if (!regularBin) v = 0.;
1203 fIntegral[i + 1] = fIntegral[i] + v;
1208 if (fIntegral[GetNbins()] == 0.) {
1209 Error(
"ComputeIntegral",
"No hits in regular bins (non over/underflow).");
1210 delete [] fIntegral;
1215 for (Long64_t j = 0; j <= GetNbins(); ++j)
1216 fIntegral[j] = fIntegral[j] / fIntegral[GetNbins()];
1219 fIntegralStatus = kValidInt;
1220 return fIntegral[GetNbins()];
1227 void THnBase::PrintBin(Long64_t idx, Option_t* options)
const
1229 Int_t* coord =
new Int_t[fNdimensions];
1230 PrintBin(idx, coord, options);
1241 Bool_t THnBase::PrintBin(Long64_t idx, Int_t* bin, Option_t* options)
const
1246 v = GetBinContent(idx);
1248 v = GetBinContent(idx, bin);
1252 if (GetCalculateErrors()) {
1254 err = GetBinError(idx);
1258 if (v == 0. && err == 0. && options && strchr(options,
'0')) {
1264 for (Int_t dim = 0; dim < fNdimensions; ++dim) {
1268 coord.Remove(coord.Length() - 1);
1270 if (GetCalculateErrors()) {
1271 Printf(
"Bin at (%s) = %g (+/- %g)", coord.Data(), v, err);
1273 Printf(
"Bin at (%s) = %g", coord.Data(), v);
1285 void THnBase::PrintEntries(Long64_t from , Long64_t howmany ,
1286 Option_t* options )
const
1288 if (from < 0) from = 0;
1289 if (howmany == -1) howmany = GetNbins();
1291 Int_t* bin =
new Int_t[fNdimensions];
1293 if (options && (strchr(options,
'x') || strchr(options,
'X'))) {
1294 Int_t* nbins =
new Int_t[fNdimensions];
1295 for (Int_t dim = fNdimensions - 1; dim >= 0; --dim) {
1296 nbins[dim] = GetAxis(dim)->GetNbins();
1297 bin[dim] = from % nbins[dim];
1301 for (Long64_t i = 0; i < howmany; ++i) {
1302 if (!PrintBin(-1, bin, options))
1305 ++bin[fNdimensions - 1];
1306 for (Int_t dim = fNdimensions - 1; dim >= 0; --dim) {
1307 if (bin[dim] >= nbins[dim]) {
1319 for (Long64_t i = from; i < from + howmany; ++i) {
1320 if (!PrintBin(i, bin, options))
1335 void THnBase::Print(Option_t* options)
const
1337 Bool_t optAxis = options && (strchr(options,
'A') || (strchr(options,
'a')));
1338 Bool_t optMem = options && (strchr(options,
'M') || (strchr(options,
'm')));
1339 Bool_t optStat = options && (strchr(options,
'S') || (strchr(options,
's')));
1340 Bool_t optContent = options && (strchr(options,
'C') || (strchr(options,
'c')));
1342 Printf(
"%s (*0x%lx): \"%s\" \"%s\"", IsA()->GetName(), (
unsigned long)
this, GetName(), GetTitle());
1343 Printf(
" %d dimensions, %g entries in %lld filled bins", GetNdimensions(), GetEntries(), GetNbins());
1346 for (Int_t dim = 0; dim < fNdimensions; ++dim) {
1347 TAxis* axis = GetAxis(dim);
1348 Printf(
" axis %d \"%s\": %d bins (%g..%g), %s bin sizes", dim,
1349 axis->GetTitle(), axis->GetNbins(), axis->GetXmin(), axis->GetXmax(),
1350 (axis->GetXbins() ?
"variable" :
"fixed"));
1355 Printf(
" %s error calculation", (GetCalculateErrors() ?
"with" :
"without"));
1356 if (GetCalculateErrors()) {
1357 Printf(
" Sum(w)=%g, Sum(w^2)=%g", GetSumw(), GetSumw2());
1358 for (Int_t dim = 0; dim < fNdimensions; ++dim) {
1359 Printf(
" axis %d: Sum(w*x)=%g, Sum(w*x^2)=%g", dim, GetSumwx(dim), GetSumwx2(dim));
1364 if (optMem && InheritsFrom(THnSparse::Class())) {
1365 const THnSparse* hsparse =
dynamic_cast<const THnSparse*
>(
this);
1366 Printf(
" coordinates stored in %d chunks of %d entries\n %g of bins filled using %g of memory compared to an array",
1367 hsparse->GetNChunks(), hsparse->GetChunkSize(),
1368 hsparse->GetSparseFractionBins(), hsparse->GetSparseFractionMem());
1372 Printf(
" BIN CONTENT:");
1373 PrintEntries(0, -1, options);
1382 void THnBase::Browse(TBrowser *b)
1384 if (fBrowsables.IsEmpty()) {
1385 for (Int_t dim = 0; dim < fNdimensions; ++dim) {
1386 fBrowsables.AddAtAndExpand(
new ROOT::Internal::THnBaseBrowsable(
this, dim), dim);
1388 fBrowsables.SetOwner();
1391 for (Int_t dim = 0; dim < fNdimensions; ++dim) {
1392 b->Add(fBrowsables[dim]);
1405 ROOT::Internal::THnBaseBinIter::~THnBaseBinIter() {
1417 THnIter::~THnIter() {
1428 ClassImp(ROOT::Internal::THnBaseBrowsable);
1433 ROOT::Internal::THnBaseBrowsable::THnBaseBrowsable(THnBase* hist, Int_t axis):
1434 fHist(hist), fAxis(axis), fProj(0)
1436 TString axisName = hist->GetAxis(axis)->GetName();
1437 if (axisName.IsNull()) {
1438 axisName = TString::Format(
"axis%d", axis);
1441 SetNameTitle(axisName,
1442 TString::Format(
"Projection on %s of %s", axisName.Data(),
1443 hist->IsA()->GetName()).Data());
1449 ROOT::Internal::THnBaseBrowsable::~THnBaseBrowsable()
1457 void ROOT::Internal::THnBaseBrowsable::Browse(TBrowser* b)
1460 fProj = fHist->Projection(fAxis);
1462 fProj->Draw(b ? b->GetDrawOption() :
"");