FairRoot/PandaRoot
PndTimeStructureAnaTask.h
Go to the documentation of this file.
1 /********************************************************************************
2  * Copyright (C) 2014 GSI Helmholtzzentrum fuer Schwerionenforschung GmbH *
3  * *
4  * This software is distributed under the terms of the *
5  * GNU Lesser General Public Licence version 3 (LGPL) version 3, *
6  * copied verbatim in the file "LICENSE" *
7  ********************************************************************************/
8 
13 #ifndef PndTimeStructureAnaTask_H
14 #define PndTimeStructureAnaTask_H
15 
16 #include "FairTask.h" // for FairTask, InitStatus
17 #include "FairTSBufferFunctional.h"
18 
19 #include "TH1.h"
20 #include "TGraph.h"
21 
22 #include "Rtypes.h" // for Bool_t, Int_t, kTRUE, etc
23 #include "TString.h" // for TString
24 
25 class FairTimeStamp;
26 class TClonesArray;
27 
28 struct DataObject
29 {
31 
32  DataObject(TString branchName): fBranchName(branchName), fBranch(0), fOldEventNr(-1), fOldTimeStamp(0), fEventMixture(kFALSE)
33  {
34  TString histoName("hBetweenEvents");
35  histoName += fBranchName;
36  TString eventHistoName("hInsideEvent");
37  eventHistoName += fBranchName;
38  TString eventDiffHistoName("hEventLength");
39  eventDiffHistoName += fBranchName;
40  fTimeHisto = new TH1D(histoName, histoName, 10000, 0, 10000);
41  fEventHisto = new TH1D(eventHistoName, eventHistoName, 100, 0, 100);
42  fEventDiffHisto = new TH1D(eventDiffHistoName, eventDiffHistoName, 1000, 0, 1000);
43  }
44  std::vector<Double_t> CalcIntegral(TH1* h1){
45  std::vector<Double_t> result;
46  Double_t sum = 0;
47  for (int i = 0; i < h1->GetNbinsX(); i++){
48  sum += h1->GetBinContent(i);
49  result.push_back(sum);
50  }
51  return result;
52  }
53 
54  std::vector<Double_t> GetBinCenters(TH1* h1){
55  std::vector<Double_t> result;
56  for (int i = 0; i < h1->GetNbinsX(); i++){
57  result.push_back(h1->GetBinCenter(i));
58  }
59  return result;
60  }
61  void FillGraphs(){
62 
63  TString eventGapName = fBranchName + "EventGap";
64  fEventGap = new TH1D(eventGapName.Data(), eventGapName.Data(), fTimeHisto->GetNbinsX(), fTimeHisto->GetXaxis()->GetXmin(), fTimeHisto->GetXaxis()->GetXmax());
66 
67  TString eventGapPercName = fBranchName + "EventGapPerc";
68  fEventGapPerc = new TH1D(eventGapPercName.Data(),eventGapPercName.Data(), fTimeHisto->GetNbinsX(), fTimeHisto->GetXaxis()->GetXmin(), fTimeHisto->GetXaxis()->GetXmax());
70 
71  TString overlapName = fBranchName + "Overlap";
72  fOverlap = new TH1D(overlapName.Data(), overlapName.Data(),fEventHisto->GetNbinsX(), fEventHisto->GetXaxis()->GetXmin(), fEventHisto->GetXaxis()->GetXmax());
74 
75  TString overlapPercName = fBranchName + "OverlapPerc";
76  fOverlapPerc = new TH1D(overlapPercName.Data(), overlapPercName.Data(), fEventHisto->GetNbinsX(), fEventHisto->GetXaxis()->GetXmin(), fEventHisto->GetXaxis()->GetXmax());
78  }
79 
80  void FillGraph(TH1D* theGraph, TH1* h1){
81  std::vector<Double_t> values = CalcIntegral(h1);
82  std::vector<Double_t> binCenters = GetBinCenters(h1);
83  for (size_t i = 0; i < values.size(); i++)
84  theGraph->Fill(binCenters[i], values[i]);
85  }
86 
87  void FillGraphPercent(TH1D* theGraph, TH1* h1){
88  std::vector<Double_t> values = CalcIntegral(h1);
89  std::vector<Double_t> binCenters = GetBinCenters(h1);
90  for (size_t i = 0; i < values.size(); i++)
91  theGraph->Fill(binCenters[i], values[i]/values.back() * 100.0);
92  }
93 
95  TClonesArray* fBranch;
96  TH1D* fTimeHisto;
97  TH1D* fEventHisto;
99  TH1D* fEventGap;
101  TH1D* fOverlap;
103  Int_t fOldEventNr;
106  std::map<Int_t, std::pair<Double_t, Double_t> > fEventStartStopMap;
107 };
108 
109 class PndTimeStructureAnaTask : public FairTask
110 {
111  public:
112 
115  FairTask("TimeStructureAnaTask"), fEntryNr(0)
116  {
117  SetVerbose(0);
118  }
119 
122  FairTask(name), fEntryNr(0)
123  {
124  SetVerbose(0);
125  }
126 
129  }
130 
132  virtual InitStatus Init();
133  virtual InitStatus ReInit();
134 
136  virtual void Exec(Option_t* opt);
137  virtual void FinishEvent();
138  virtual void FinishTask();
139 
140  virtual void SetParContainers() {};
141 
143  fData.push_back(DataObject(name));
144  name += "Prim";
145  fDataPrim.push_back(DataObject(name));
146  }
147 
148 
149  protected:
150 
151  std::vector<DataObject> fData; //data for all type of particles
152  std::vector<DataObject> fDataPrim; //data for primary particles
155 
156  Int_t fEntryNr;
157 
159 
160 };
161 
162 #endif
std::vector< DataObject > fData
Int_t i
Definition: run_full.C:25
DataObject(TString branchName)
std::vector< DataObject > fDataPrim
std::vector< Double_t > GetBinCenters(TH1 *h1)
void FillGraph(TH1D *theGraph, TH1 *h1)
ClassDef(PndTimeStructureAnaTask, 1)
std::map< Int_t, std::pair< Double_t, Double_t > > fEventStartStopMap
Double_t
PndTimeStructureAnaTask(const char *name)
TString name
virtual void Exec(Option_t *opt)
std::vector< Double_t > CalcIntegral(TH1 *h1)
drchit SetVerbose(iVerbose)
TClonesArray * fBranch
void FillGraphPercent(TH1D *theGraph, TH1 *h1)