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