This file uses 4 large data sets from the H1 collaboration at DESY Hamburg. One can access these data sets (277 MBytes) from the standard Root web site at: https://root.cern.ch/files/h1/ The Physics plots below generated by this example cannot be produced when using smaller data sets.
A chain of 4 files (originally converted from PAW ntuples) is used to illustrate the various ways to loop on Root data sets. Each data set contains a Root Tree named "h42" The class definition in h1analysis.h has been generated automatically by the Root utility TTree::MakeSelector using one of the files with the following statement:
The entry list is saved to a file "elist.root" by the Terminate function. To see the list of selected events, you can do elist->Print("all")
. The selection function has selected 7525 events out of the 283813 events in the chain of files. (2.65 per cent)
You can repeat the steps B, C and D using the script compiler by replacing "h1analysis.C" by "h1analysis.C+" or "h1analysis.C++" in a new session (see F).
The same analysis can be run on PROOF. For a quick try start a PROOF-Lite session
You can then repeat step B above. Step C can also be executed in PROOF. However, step D cannot be executed in PROOF as in the local session (i.e. just passing option 'useList'): to use the entry list you have to
call Process as in step B. Of course this works also for local processing.
{
if (
x <= 0.13957)
return 0;
+ par[2] / 2.5066/par[4]*
TMath::Exp(-xp3/2/par[4]/par[4]));
return res;
}
{
if (
x <= 0.13957)
return 0;
return res;
}
{
Info(
"Begin",
"starting h1analysis with process option: %s",
option.Data());
delete gDirectory->GetList()->FindObject(
"elist");
if (
option.Contains(
"fillList")) {
}
Info(
"Begin",
"creating an entry-list");
}
if (
option.Contains(
"useList")) {
} else {
}
}
}
{
"starting h1analysis with process option: %s (tree: %p)",
option.Data(),
tree);
hdmd =
new TH1F(
"hdmd",
"dm_d",40,0.13,0.17);
h2 =
new TH2F(
"h2",
"ptD0 vs dm_d",30,0.135,0.165,30,-3,6);
if (
option.Contains(
"fillList")) {
else
}
}
}
{
}
}
{
}
{
if (
hdmd ==
nullptr ||
h2 ==
nullptr) {
return;
}
c1->SetBottomMargin(0.15);
if (
gROOT->GetListOfFunctions()->FindObject(
"f5"))
delete gROOT->GetFunction(
"f5");
c2->SetBottomMargin(0.15);
if (
gROOT->GetListOfFunctions()->FindObject(
"f2"))
delete gROOT->GetFunction(
"f2");
Printf(
"Entry list 'elist' created:");
TFile efile(
"elist.root",
"recreate");
} else {
Error(
"Terminate",
"entry list requested but not found in output");
}
}
}
void Printf(const char *fmt,...)
Formats a string in a circular formatting buffer and prints the string.
R__EXTERN TStyle * gStyle
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title.
virtual void SetMarkerStyle(Style_t mstyle=1)
Set the marker style.
virtual Int_t GetEntry(Long64_t entry=0, Int_t getall=0)
Read all leaves of entry and return total number of bytes read.
virtual Bool_t Enter(Long64_t entry, TTree *tree=nullptr)
Add entry #entry to the list.
virtual void SetDirectory(TDirectory *dir)
Add reference to directory dir. dir can be 0.
void Print(const Option_t *option="") const override
Print this list.
virtual void SetParameters(const Double_t *params)
1-D histogram with a double per channel (see TH1 documentation)}
1-D histogram with a float per channel (see TH1 documentation)}
virtual TFitResultPtr Fit(const char *formula, Option_t *option="", Option_t *goption="", Double_t xmin=0, Double_t xmax=0)
Fit histogram with function fname.
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
void Draw(Option_t *option="") override
Draw this histogram with options.
TList * GetListOfFunctions() const
2-D histogram with a float per channel (see TH1 documentation)}
Int_t Fill(Double_t) override
Invalid Fill method.
virtual void FitSlicesX(TF1 *f1=nullptr, Int_t firstybin=0, Int_t lastybin=-1, Int_t cut=0, Option_t *option="QNR", TObjArray *arr=nullptr)
Project slices along X in case of a 2-D histogram, then fit each slice with function f1 and make a hi...
TObject * FindObject(const char *name) const override
Find object using its name.
Use the TLine constructor to create a simple line.
TObject * FindObject(const char *name) const override
Find an object in this list using its name.
void Add(TObject *obj) override
The TNamed class is the base class for all named ROOT classes.
TObject * Clone(const char *newname="") const override
Make a clone of an object using the Streamer facility.
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
virtual Int_t Write(const char *name=nullptr, Int_t option=0, Int_t bufsize=0)
Write this object to the current directory.
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
virtual void Draw(Option_t *option="")
Default Draw method for all objects.
virtual void Info(const char *method, const char *msgfmt,...) const
Issue info message.
The histogram statistics painter class.
void SetOptStat(Int_t stat=1)
Set the stat option.
TList * fInput
List of objects available during processing.
TSelectorList * fOutput
! List of objects created during processing
Long64_t fStatus
Selector status.
const char * GetOption() const override
void SetOptStat(Int_t stat=1)
The type of information printed in the histogram statistics box can be selected via the parameter mod...
void SetOptFit(Int_t fit=1)
The type of information about fit parameters printed in the histogram statistics box can be selected ...
A TTree represents a columnar dataset.
virtual void SetEntryList(TEntryList *list, Option_t *opt="")
Set an EntryList.
void SlaveTerminate() override
void Init(TTree *tree) override
Bool_t Process(Long64_t entry) override
The Process() function is called for each entry in the tree (or possibly keyed object in the case of ...
void Begin(TTree *tree) override
void SlaveBegin(TTree *tree) override
void Terminate() override
Double_t fdm5(Double_t *xx, Double_t *par)
Double_t fdm2(Double_t *xx, Double_t *par)
Double_t Exp(Double_t x)
Returns the base-e exponential function of x, which is e raised to the power x.
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Returns x raised to the power y.
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.