1#ifndef __JROOT__JROOTTOOLKIT__
2#define __JROOT__JROOTTOOLKIT__
11#pragma GCC diagnostic push
12#pragma GCC diagnostic ignored "-Wall"
14#include "TObjString.h"
19#include "TStreamerInfo.h"
26#include "TGraphErrors.h"
27#include "TGraphAsymmErrors.h"
29#include "TGraph2DErrors.h"
30#include "TMultiGraph.h"
32#pragma GCC diagnostic pop
46namespace JPP {
using namespace JROOT; }
74 return T::Class_Name();
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) {
92 h3->SetBinError(ix, iy, iz, 0.0);
95 h3->SetBinContent(ix, iy, iz, 0.0);
114 for (
int ix = 1; ix <= h2->GetXaxis()->GetNbins(); ++ix) {
115 for (
int iy = 1; iy <= h2->GetYaxis()->GetNbins(); ++iy) {
117 h2->SetBinError(ix, iy, 0.0);
120 h2->SetBinContent(ix, iy, 0.0);
136#define RESET_OBJECT(T, P, O) if (dynamic_cast<T*>(P) != NULL) { resetObject(dynamic_cast<T*>(P), O); return; }
149 for (
int ix = 1; ix <= h1->GetXaxis()->GetNbins(); ++ix) {
151 h1->SetBinError(ix, 0.0);
154 h1->SetBinContent(ix, 0.0);
174 for (
int i = 0; i !=
g1->GetN(); ++i) {
177 g1->SetPoint(i,
g1->GetX()[i], 0.0);
194 for (
int i = 0; i !=
g1->GetN(); ++i) {
196 g1->SetPointError(i, 0.0, 0.0);
199 g1->SetPoint(i,
g1->GetX()[i], 0.0);
216 for (
int i = 0; i !=
g1->GetN(); ++i) {
218 g1->SetPointError(i, 0.0, 0.0, 0.0, 0.0);
221 g1->SetPoint(i,
g1->GetX()[i], 0.0);
238 for (
int i = 0; i != g2->GetN(); ++i) {
241 g2->SetPoint(i, g2->GetX()[i], g2->GetY()[i], 0.0);
254 inline void resetObject(TGraph2DErrors* g2,
const bool reset =
false)
258 for (
int i = 0; i != g2->GetN(); ++i) {
260 g2->SetPointError(i, 0.0, 0.0, 0.0);
263 g2->SetPoint(i, g2->GetX()[i], g2->GetY()[i], 0.0);
276 inline void resetObject(TMultiGraph* gs,
const bool reset =
false)
278#define RESET_OBJECT(T, P, O) if (dynamic_cast<T*>(P) != NULL) { resetObject(dynamic_cast<T*>(P), O); continue; }
282 for (TIter next(gs->GetListOfGraphs()); TGraph* graph = (TGraph*) next(); ) {
304 const Int_t n =
g1->GetN();
307 g1->SetPoint(n, x, y);
326 const Int_t n =
g1->GetN();
329 g1->SetPoint(n, x, y);
330 g1->SetPointError(n, ex, ey);
353 const Int_t n =
g1->GetN();
356 g1->SetPoint(n, x, y);
357 g1->SetPointError(n, exl, exh, eyl, eyh);
374 const Int_t n =
g1->GetN();
377 g1->SetPoint(n, x, y, z);
400 const Int_t n =
g1->GetN();
403 g1->SetPoint(n, x, y, z);
404 g1->SetPointError(n, ex, ey, ez);
417 file.WriteTObject(&
object);
433 if (file != NULL && file->IsOpen()) {
434 return dynamic_cast<const TStreamerInfo*
>(file->GetStreamerInfoList()->FindObject(name));
449 inline const TStreamerInfo*
getStreamerInfo(
const char*
const file_name,
const char*
const name)
469 if (pStreamerInfo != NULL)
470 return pStreamerInfo->GetClassVersion();
502 inline TString
parse(
const TPRegexp& regexp,
const TString&
string,
const int index = 1)
504 TPRegexp buffer = regexp;
505 TObjArray* array = buffer.MatchS(
string);
506 TString result = string;
508 if (index - 1 < array->GetLast()) {
509 result = ((TObjString*) array->At(index))->GetName();
538 const Double_t
value) :
550 inline operator Int_t()
const
569 if (parameter.
index >= 0 && parameter.
index < f1.GetNpar()) {
571 f1.SetParameter(parameter.
index, parameter.
value);
590 if (parameter.
index >= 0 && parameter.
index < f1.GetNpar()) {
592 f1.FixParameter(parameter.
index, parameter.
value);
611 if (index >= 0 && index < f1.GetNpar()) {
613 f1.ReleaseParameter(index);
632 inline bool setParLimits(TF1& f1,
const Int_t index, Double_t xmin, Double_t xmax)
636 if (index >= 0 && index < f1.GetNpar()) {
638 if (xmin == 0.0) { xmin = -numeric_limits<Double_t>::min(); }
639 if (xmax == 0.0) { xmax = +numeric_limits<Double_t>::min(); }
641 f1.SetParLimits(index, xmin, xmax);
660 if (index >= 0 && index < f1.GetNpar()) {
665 f1.GetParLimits(index, xmin, xmax);
667 return (xmin != 0.0 && xmax != 0.0 && xmin >= xmax);
688 const int N = h1.GetNbinsX();
692 for (
int i = 0; i < N; i++) {
693 x[i] = h1.GetBinCenter (i + 1);
694 y[i] = h1.GetBinContent(i + 1);
697 return new TGraph(N, &x[0], &y[0]);
710 inline TH1*
projectHistogram(
const TH2& h2,
const Double_t xmin,
const Double_t xmax,
const char projection)
712 switch (projection) {
716 return h2.ProjectionX(
"_px", h2.GetYaxis()->FindBin(xmin), h2.GetYaxis()->FindBin(xmax));
720 return h2.ProjectionY(
"_py", h2.GetXaxis()->FindBin(xmin), h2.GetXaxis()->FindBin(xmax));
736inline std::istream&
operator>>(std::istream& in, TRegexp&
object)
741 object = TRegexp(buffer.c_str());
755inline std::ostream&
operator<<(std::ostream& out,
const TRegexp&
object)
768inline std::istream&
operator>>(std::istream& in, TFormula&
object)
772 if (getline(in,buffer)) {
773 object = TFormula(
"", buffer.c_str());
787inline std::ostream&
operator<<(std::ostream& out,
const TFormula&
object)
789 return out <<
object.GetExpFormula();
I/O formatting auxiliaries.
Double_t g1(const Double_t x)
Function.
TFile * getFile() const
Get file.
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 ¶meter)
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.
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 ¶meter)
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 ®exp, 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.
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.