Jpp 20.0.0-27-g39925593c-D
the software that should make you happy
Loading...
Searching...
No Matches
JRootToolkit.hh
Go to the documentation of this file.
1#ifndef __JROOT__JROOTTOOLKIT__
2#define __JROOT__JROOTTOOLKIT__
3
4#include <string>
5#include <istream>
6#include <ostream>
7#include <limits>
8#include <cctype>
9#include <vector>
10
11#pragma GCC diagnostic push
12#pragma GCC diagnostic ignored "-Wall"
13#include "TString.h"
14#include "TObjString.h"
15#include "TRegexp.h"
16#include "TPRegexp.h"
17#include "TFormula.h"
18#include "TFile.h"
19#include "TStreamerInfo.h"
20#include "TIterator.h"
21#include "TF1.h"
22#include "TH1.h"
23#include "TH2.h"
24#include "TH3.h"
25#include "TGraph.h"
26#include "TGraphErrors.h"
27#include "TGraphAsymmErrors.h"
28#include "TGraph2D.h"
29#include "TGraph2DErrors.h"
30#include "TMultiGraph.h"
31#include "TNtuple.h"
32#pragma GCC diagnostic pop
33
34#include "JLang/JType.hh"
35#include "JLang/JLangToolkit.hh"
36
37#include "JROOT/JRootFile.hh"
38
39/**
40 * \author mdejong, mlincett
41 */
42
43namespace JROOT {}
44namespace JPP { using namespace JROOT; }
45
46namespace JROOT {
47
48 using JLANG::JType;
49
50
51 /**
52 * Get ROOT name of given data type.
53 *
54 * \return name of object to be read
55 */
56 template<class T>
57 inline const char* getName()
58 {
59 return getName(JType<T>());
60 }
61
62
63 /**
64 * Get ROOT name of given data type.
65 *
66 * \param type data type
67 * \return name of object to be read
68 */
69 template<class T>
70 inline const char* getName(const JType<T>& type)
71 {
72 return T::Class_Name();
73 }
74
75
76 /**
77 * Reset TH3 object.
78 *
79 * \param h3 pointer to TH3 object
80 * \param reset reset contents if true
81 */
82 inline void resetObject(TH3* h3, const bool reset = false)
83 {
84 if (h3 != NULL) {
85
86 for (int ix = 1; ix <= h3->GetXaxis()->GetNbins(); ++ix) {
87 for (int iy = 1; iy <= h3->GetYaxis()->GetNbins(); ++iy) {
88 for (int iz = 1; iz <= h3->GetZaxis()->GetNbins(); ++iz) {
89
90 h3->SetBinError(ix, iy, iz, 0.0);
91
92 if (reset) {
93 h3->SetBinContent(ix, iy, iz, 0.0);
94 }
95 }
96 }
97 }
98 }
99 }
100
101
102 /**
103 * Reset TH2 object.
104 *
105 * \param h2 pointer to TH2 object
106 * \param reset reset contents if true
107 */
108 inline void resetObject(TH2* h2, const bool reset = false)
109 {
110 if (h2 != NULL) {
111
112 for (int ix = 1; ix <= h2->GetXaxis()->GetNbins(); ++ix) {
113 for (int iy = 1; iy <= h2->GetYaxis()->GetNbins(); ++iy) {
114
115 h2->SetBinError(ix, iy, 0.0);
116
117 if (reset) {
118 h2->SetBinContent(ix, iy, 0.0);
119 }
120 }
121 }
122 }
123 }
124
125
126 /**
127 * Reset TH1 object.
128 *
129 * \param h1 pointer to TH1 object
130 * \param reset reset contents if true
131 */
132 inline void resetObject(TH1* h1, const bool reset = false)
133 {
134#define RESET_OBJECT(T, P, O) if (dynamic_cast<T*>(P) != NULL) { resetObject(dynamic_cast<T*>(P), O); return; }
135
136 if (h1 != NULL) {
137
138 if (reset) {
139
140 h1->Reset();
141
142 } else {
143
144 RESET_OBJECT(TH3, h1, reset);
145 RESET_OBJECT(TH2, h1, reset);
146
147 for (int ix = 1; ix <= h1->GetXaxis()->GetNbins(); ++ix) {
148
149 h1->SetBinError(ix, 0.0);
150
151 if (reset) {
152 h1->SetBinContent(ix, 0.0);
153 }
154 }
155 }
156 }
157
158#undef RESET_OBJECT
159 }
160
161
162 /**
163 * Reset TGraph object.
164 *
165 * \param g1 pointer to TGraph object
166 * \param reset reset contents if true
167 */
168 inline void resetObject(TGraph* g1, const bool reset = false)
169 {
170 if (g1 != NULL) {
171
172 for (int i = 0; i != g1->GetN(); ++i) {
173
174 if (reset) {
175 g1->SetPoint(i, g1->GetX()[i], 0.0);
176 }
177 }
178 }
179 }
180
181
182 /**
183 * Reset TGraphErrors object.
184 *
185 * \param g1 pointer to TGraphErrors object
186 * \param reset reset contents if true
187 */
188 inline void resetObject(TGraphErrors* g1, const bool reset = false)
189 {
190 if (g1 != NULL) {
191
192 for (int i = 0; i != g1->GetN(); ++i) {
193
194 g1->SetPointError(i, 0.0, 0.0);
195
196 if (reset) {
197 g1->SetPoint(i, g1->GetX()[i], 0.0);
198 }
199 }
200 }
201 }
202
203
204 /**
205 * Reset TGraphAsymmErrors object.
206 *
207 * \param g1 pointer to TGraphErrors object
208 * \param reset reset contents if true
209 */
210 inline void resetObject(TGraphAsymmErrors* g1, const bool reset = false)
211 {
212 if (g1 != NULL) {
213
214 for (int i = 0; i != g1->GetN(); ++i) {
215
216 g1->SetPointError(i, 0.0, 0.0, 0.0, 0.0);
217
218 if (reset) {
219 g1->SetPoint(i, g1->GetX()[i], 0.0);
220 }
221 }
222 }
223 }
224
225
226 /**
227 * Reset TGraph2D object.
228 *
229 * \param g2 pointer to TGraph2D object
230 * \param reset reset contents if true
231 */
232 inline void resetObject(TGraph2D* g2, const bool reset = false)
233 {
234 if (g2 != NULL) {
235
236 for (int i = 0; i != g2->GetN(); ++i) {
237
238 if (reset) {
239 g2->SetPoint(i, g2->GetX()[i], g2->GetY()[i], 0.0);
240 }
241 }
242 }
243 }
244
245
246 /**
247 * Reset TGraph2DErrors object.
248 *
249 * \param g2 pointer to TGraph2DErrors object
250 * \param reset reset contents if true
251 */
252 inline void resetObject(TGraph2DErrors* g2, const bool reset = false)
253 {
254 if (g2 != NULL) {
255
256 for (int i = 0; i != g2->GetN(); ++i) {
257
258 g2->SetPointError(i, 0.0, 0.0, 0.0);
259
260 if (reset) {
261 g2->SetPoint(i, g2->GetX()[i], g2->GetY()[i], 0.0);
262 }
263 }
264 }
265 }
266
267
268 /**
269 * Reset TMultiGraph object.
270 *
271 * \param gs pointer to TMultiGraph object
272 * \param reset reset contents if true
273 */
274 inline void resetObject(TMultiGraph* gs, const bool reset = false)
275 {
276#define RESET_OBJECT(T, P, O) if (dynamic_cast<T*>(P) != NULL) { resetObject(dynamic_cast<T*>(P), O); continue; }
277
278 if (gs != NULL) {
279
280 for (TIter next(gs->GetListOfGraphs()); TGraph* graph = (TGraph*) next(); ) {
281 RESET_OBJECT(TGraphAsymmErrors, graph, reset);
282 RESET_OBJECT(TGraphErrors, graph, reset);
283 RESET_OBJECT(TGraph, graph, reset);
284 }
285 }
286
287#undef RESET_OBJECT
288 }
289
290
291 /**
292 * Add point to TGraph.
293 *
294 * \param g1 pointer to valid ROOT TGraph object
295 * \param x x value
296 * \param y y value
297 */
298 inline void AddPoint(TGraph* g1,
299 const Double_t x,
300 const Double_t y)
301 {
302 const Int_t n = g1->GetN();
303
304 g1->Set(n + 1);
305 g1->SetPoint(n, x, y);
306 }
307
308
309 /**
310 * Add point to TGraphErrors.
311 *
312 * \param g1 pointer to valid ROOT TGraph object
313 * \param x x value
314 * \param y y value
315 * \param ex x error
316 * \param ey y error
317 */
318 inline void AddPoint(TGraphErrors* g1,
319 const Double_t x,
320 const Double_t y,
321 const Double_t ex,
322 const Double_t ey)
323 {
324 const Int_t n = g1->GetN();
325
326 g1->Set(n + 1);
327 g1->SetPoint(n, x, y);
328 g1->SetPointError(n, ex, ey);
329 }
330
331
332 /**
333 * Add point to TGraphAsymmErrors.
334 *
335 * \param g1 pointer to valid ROOT TGraph object
336 * \param x x value
337 * \param y y value
338 * \param exl x error low
339 * \param exh x error high
340 * \param eyl y error low
341 * \param eyh y error high
342 */
343 inline void AddPoint(TGraphAsymmErrors* g1,
344 const Double_t x,
345 const Double_t y,
346 const Double_t exl,
347 const Double_t exh,
348 const Double_t eyl,
349 const Double_t eyh)
350 {
351 const Int_t n = g1->GetN();
352
353 g1->Set(n + 1);
354 g1->SetPoint(n, x, y);
355 g1->SetPointError(n, exl, exh, eyl, eyh);
356 }
357
358
359 /**
360 * Add point to TGraph2D.
361 *
362 * \param g1 pointer to valid ROOT TGraph object
363 * \param x x value
364 * \param y y value
365 * \param z z value
366 */
367 inline void AddPoint(TGraph2D* g1,
368 const Double_t x,
369 const Double_t y,
370 const Double_t z)
371 {
372 const Int_t n = g1->GetN();
373
374 g1->Set(n + 1);
375 g1->SetPoint(n, x, y, z);
376 }
377
378
379 /**
380 * Add point to TGraph2DErrors.
381 *
382 * \param g1 pointer to valid ROOT TGraph object
383 * \param x x value
384 * \param y y value
385 * \param z z value
386 * \param ex x error
387 * \param ey y error
388 * \param ez z error
389 */
390 inline void AddPoint(TGraph2DErrors* g1,
391 const Double_t x,
392 const Double_t y,
393 const Double_t z,
394 const Double_t ex,
395 const Double_t ey,
396 const Double_t ez)
397 {
398 const Int_t n = g1->GetN();
399
400 g1->Set(n + 1);
401 g1->SetPoint(n, x, y, z);
402 g1->SetPointError(n, ex, ey, ez);
403 }
404
405
406 /**
407 * Write object to ROOT file.
408 *
409 * \param file ROOT file
410 * \param object ROOT object
411 * \return ROOT file
412 */
413 inline TFile& operator<<(TFile& file, const TObject& object)
414 {
415 file.WriteTObject(&object);
416
417 return file;
418 }
419
420
421 /**
422 * Get ROOT streamer information of class with given name.
423 * Note that the class name should include the name space, if any.
424 *
425 * \param file pointer to ROOT file
426 * \param name class name
427 * \return pointer to TStreamerInfo (NULL in case of error)
428 */
429 inline const TStreamerInfo* getStreamerInfo(TFile* file, const char* const name)
430 {
431 if (file != NULL && file->IsOpen()) {
432
433 TList* ps = file->GetStreamerInfoList();
434
435 if (ps) {
436 return dynamic_cast<const TStreamerInfo*>(ps->FindObject(name));
437 }
438 }
439
440 return NULL;
441 }
442
443
444 /**
445 * Get ROOT streamer information of class with given name.
446 * Note that the class name should include the name space, if any.
447 *
448 * \param file_name file name
449 * \param name class name
450 * \return pointer to TStreamerInfo (NULL in case of error)
451 */
452 inline const TStreamerInfo* getStreamerInfo(const char* const file_name, const char* const name)
453 {
454 JRootInputFile file(file_name);
455
456 return getStreamerInfo(file.getFile(), name);
457 }
458
459
460 /**
461 * Get ROOT streamer version of class with given name.
462 * Note that the class name should include the name space, if any.
463 *
464 * \param file pointer to ROOT file
465 * \param name class name
466 * \return streamer version; (-1 in case of error)
467 */
468 inline int getStreamerVersion(TFile* file, const char* const name)
469 {
470 const TStreamerInfo* pStreamerInfo = getStreamerInfo(file, name);
471
472 if (pStreamerInfo != NULL)
473 return pStreamerInfo->GetClassVersion();
474 else
475 return -1;
476 }
477
478
479 /**
480 * Get ROOT streamer version of class with given name.
481 * Note that the class name should include the name space, if any.
482 *
483 * \param file_name file name
484 * \param name class name
485 * \return streamer version; (-1 in case of error)
486 */
487 inline int getStreamerVersion(const char* const file_name, const char* const name)
488 {
489 JRootInputFile file(file_name);
490
491 return getStreamerVersion(file.getFile(), name);
492 }
493
494
495 /**
496 * Match a regular expression with given string and return the specified matched parts.
497 *
498 * If no matches are found corresponding to the specified index, the original string is returned.
499 *
500 * \param regexp regular expression
501 * \param string input string
502 * \param index index of matched parts (starting at 1)
503 * \return matched part of string
504 */
505 inline TString parse(const TPRegexp& regexp, const TString& string, const int index = 1)
506 {
507 TPRegexp buffer = regexp;
508 TObjArray* array = buffer.MatchS(string);
509 TString result = string;
510
511 if (index - 1 < array->GetLast()) {
512 result = ((TObjString*) array->At(index))->GetName();
513 }
514
515 delete array;
516
517 return result;
518 }
519
520
521 /**
522 * Auxiliary data structure for a parameter index and its value.
523 */
525 /**
526 * Default constructor.
527 */
529 index(0),
530 value(0.0)
531 {}
532
533
534 /**
535 * Constructor.
536 *
537 * \param index parameter index
538 * \param value parameter value
539 */
541 const Double_t value) :
542 index(index),
543 value(value)
544 {}
545
546
547 /**
548 * Type conversion operator
549 *
550 *
551 * \return index
552 */
553 inline operator Int_t() const
554 {
555 return index;
556 }
557
558
559 Int_t index;
560 Double_t value;
561 };
562
563
564 /**
565 * Set fit parameter.
566 *
567 * \param f1 fit function
568 * \param parameter parameter index and value
569 */
570 inline bool setParameter(TF1& f1, const JFitParameter_t& parameter)
571 {
572 if (parameter.index >= 0 && parameter.index < f1.GetNpar()) {
573
574 f1.SetParameter(parameter.index, parameter.value);
575
576 return true;
577
578 } else {
579
580 return false;
581 }
582 }
583
584
585 /**
586 * Fix fit parameter.
587 *
588 * \param f1 fit function
589 * \param parameter parameter index and value
590 */
591 inline bool fixParameter(TF1& f1, const JFitParameter_t& parameter)
592 {
593 if (parameter.index >= 0 && parameter.index < f1.GetNpar()) {
594
595 f1.FixParameter(parameter.index, parameter.value);
596
597 return true;
598
599 } else {
600
601 return false;
602 }
603 }
604
605
606 /**
607 * Release fit parameter.
608 *
609 * \param f1 fit function
610 * \param index parameter index
611 */
612 inline bool releaseParameter(TF1& f1, const Int_t index)
613 {
614 if (index >= 0 && index < f1.GetNpar()) {
615
616 f1.ReleaseParameter(index);
617
618 return true;
619
620 } else {
621
622 return false;
623 }
624 }
625
626
627 /**
628 * Set fit parameter limits.
629 *
630 * \param f1 fit function
631 * \param index parameter index
632 * \param xmin lower limit
633 * \param xmax upper limit
634 */
635 inline bool setParLimits(TF1& f1, const Int_t index, Double_t xmin, Double_t xmax)
636 {
637 using namespace std;
638
639 if (index >= 0 && index < f1.GetNpar()) {
640
641 if (xmin == 0.0) { xmin = -numeric_limits<Double_t>::min(); }
642 if (xmax == 0.0) { xmax = +numeric_limits<Double_t>::min(); }
643
644 f1.SetParLimits(index, xmin, xmax);
645
646 return true;
647
648 } else {
649
650 return false;
651 }
652 }
653
654
655 /**
656 * Check if fit parameter is fixed.
657 *
658 * \param f1 fit function
659 * \param index parameter index
660 */
661 inline bool isParameterFixed(const TF1& f1, const Int_t index)
662 {
663 if (index >= 0 && index < f1.GetNpar()) {
664
665 Double_t xmin;
666 Double_t xmax;
667
668 f1.GetParLimits(index, xmin, xmax);
669
670 return (xmin != 0.0 && xmax != 0.0 && xmin >= xmax);
671
672 } else {
673
674 return false;
675 }
676 }
677 /**
678 * Helper method to convert a 1D histogram to a graph.
679 *
680 * The graph consists of:
681 * - bin centers as x values;
682 * - bin contents as y values.
683 *
684 * Note that underflow and overflow bins are ignored.
685 *
686 * \param h1 1D histogram
687 * \return pointer to newly create graph
688 */
689 inline TGraph* histogramToGraph(const TH1& h1)
690 {
691 const int N = h1.GetNbinsX();
692
693 std::vector<double> x(N), y(N);
694
695 for (int i = 0; i < N; i++) {
696 x[i] = h1.GetBinCenter (i + 1);
697 y[i] = h1.GetBinContent(i + 1);
698 }
699
700 return new TGraph(N, &x[0], &y[0]);
701 }
702
703
704 /**
705 * Helper method for ROOT histogram projections.
706 *
707 * \param h2 2D histogram
708 * \param xmin lower limit
709 * \param xmax upper limit
710 * \param projection projection (x|X|y|Y)
711 * \return pointer to newly created 1D histogram
712 */
713 inline TH1* projectHistogram(const TH2& h2, const Double_t xmin, const Double_t xmax, const char projection)
714 {
715 switch (projection) {
716
717 case 'x':
718 case 'X':
719 return h2.ProjectionX("_px", h2.GetYaxis()->FindBin(xmin), h2.GetYaxis()->FindBin(xmax));
720
721 case 'y':
722 case 'Y':
723 return h2.ProjectionY("_py", h2.GetXaxis()->FindBin(xmin), h2.GetXaxis()->FindBin(xmax));
724
725 default:
726 return NULL;
727 }
728 }
729}
730
731
732/**
733 * Read regular expression from input stream
734 *
735 * \param in input stream
736 * \param object regular expression
737 * \return output stream
738 */
739inline std::istream& operator>>(std::istream& in, TRegexp& object)
740{
741 std::string buffer;
742
743 if (in >> buffer) {
744 object = TRegexp(buffer.c_str());
745 }
746
747 return in;
748}
749
750
751/**
752 * Write regular expression to output stream
753 *
754 * \param out output stream
755 * \param object regular expression
756 * \return output stream
757 */
758inline std::ostream& operator<<(std::ostream& out, const TRegexp& object)
759{
760 return out;
761}
762
763
764/**
765 * Read formula from input stream.
766 *
767 * \param in input stream
768 * \param object formula
769 * \return input stream
770 */
771inline std::istream& operator>>(std::istream& in, TFormula& object)
772{
773 std::string buffer;
774
775 if (getline(in,buffer)) {
776 object = TFormula("", buffer.c_str());
777 }
778
779 return in;
780}
781
782
783/**
784 * Write formula to output stream.
785 *
786 * \param out output stream
787 * \param object formula
788 * \return output stream
789 */
790inline std::ostream& operator<<(std::ostream& out, const TFormula& object)
791{
792 return out << object.GetExpFormula();
793}
794
795#endif
796
Double_t g1(const Double_t x)
Function.
Definition JQuantiles.cc:25
std::istream & operator>>(std::istream &in, TRegexp &object)
Read regular expression from input stream.
#define RESET_OBJECT(T, P, O)
std::ostream & operator<<(std::ostream &out, const TRegexp &object)
Write regular expression to output stream.
TFile * getFile() const
Get file.
Definition JRootFile.hh:66
ROOT input file.
Definition JRootFile.hh:97
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
Auxiliary classes and methods for ROOT I/O.
const TStreamerInfo * getStreamerInfo(TFile *file, const char *const name)
Get ROOT streamer information of class with given name.
bool fixParameter(TF1 &f1, const JFitParameter_t &parameter)
Fix fit parameter.
TGraph * histogramToGraph(const TH1 &h1)
Helper method to convert a 1D histogram to a graph.
TFile & operator<<(TFile &file, const TObject &object)
Write object to ROOT file.
void AddPoint(TGraph *g1, const Double_t x, const Double_t y)
Add point to TGraph.
int getStreamerVersion(TFile *file, const char *const name)
Get ROOT streamer version of class with given name.
void resetObject(JManager< JKey_t, JValue_t > *object, const bool reset=false)
Reset JManager object.
Definition JManager.hh:394
TH1 * projectHistogram(const TH2 &h2, const Double_t xmin, const Double_t xmax, const char projection)
Helper method for ROOT histogram projections.
bool setParameter(TF1 &f1, const JFitParameter_t &parameter)
Set fit parameter.
bool isParameterFixed(const TF1 &f1, const Int_t index)
Check if fit parameter is fixed.
bool releaseParameter(TF1 &f1, const Int_t index)
Release fit parameter.
TString parse(const TPRegexp &regexp, const TString &string, const int index=1)
Match a regular expression with given string and return the specified matched parts.
bool setParLimits(TF1 &f1, const Int_t index, Double_t xmin, Double_t xmax)
Set fit parameter limits.
const char * getName()
Get ROOT name of given data type.
Auxiliary class for a type holder.
Definition JType.hh:19
Auxiliary data structure for a parameter index and its value.
JFitParameter_t(const Int_t index, const Double_t value)
Constructor.
JFitParameter_t()
Default constructor.