CbmRoot
CbmCheckEvents.cxx
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 (LGPL) version 3, *
6  * copied verbatim in the file "LICENSE" *
7  ********************************************************************************/
8 #include "CbmCheckEvents.h"
9 
10 //#include "CbmTofDigi.h"
11 #include "CbmDefs.h"
12 #include "CbmDigiManager.h"
13 #include "CbmEvent.h"
14 #include "CbmMuchBeamTimeDigi.h"
15 #include "CbmStsDigi.h"
16 #include "CbmTofDigi.h"
17 
18 #include "FairLogger.h"
19 #include "FairRootManager.h"
20 #include "FairRunOnline.h"
21 
22 #include "TClonesArray.h"
23 #include "TH1.h"
24 #include "TH2.h"
25 #include "THttpServer.h"
26 #include "TProfile.h"
27 #include <TFile.h>
28 
29 #include <iomanip>
30 using std::fixed;
31 using std::setprecision;
32 
33 // ---- Default constructor -------------------------------------------
34 CbmCheckEvents::CbmCheckEvents() : FairTask("CbmCheckEvents") {}
35 
36 // ---- Destructor ----------------------------------------------------
38 
39 // ---- Initialisation ----------------------------------------------
41  // Load all necessary parameter containers from the runtime data base
42  /*
43  FairRunAna* ana = FairRunAna::Instance();
44  FairRuntimeDb* rtdb=ana->GetRuntimeDb();
45 
46  <CbmCheckEventsDataMember> = (<ClassPointer>*)
47  (rtdb->getContainer("<ContainerName>"));
48  */
49 }
50 
51 // ---- Init ----------------------------------------------------------
52 InitStatus CbmCheckEvents::Init() {
53 
54  // Get a handle from the IO manager
55  FairRootManager* ioman = FairRootManager::Instance();
56 
57  // DigiManager
60  fDigiMan->Init();
61 
62  // Get a pointer to the previous already existing data level
63 
64  fT0DigiVec = ioman->InitObjectAs<std::vector<CbmTofDigi> const*>("T0Digi");
65  if (!fT0DigiVec) {
66  fT0DigiArr = dynamic_cast<TClonesArray*>(ioman->GetObject("T0Digi"));
67  if (!fT0DigiArr) { LOG(fatal) << "No TClonesArray with T0 digis found."; }
68  }
69 
71  LOG(info) << "No TClonesArray with STS digis found.";
72  }
73 
75  LOG(info) << "No TClonesArray with MUCH digis found.";
76  }
77 
79  LOG(info) << "No TClonesArray with TOF digis found.";
80  }
81 
82  fEvents = static_cast<TClonesArray*>(ioman->GetObject("CbmEvent"));
83  if (!fEvents) { LOG(fatal) << "No TClonesArray with events found."; }
84 
85  CreateHistos();
86 
87  return kSUCCESS;
88 }
89 
90 // ---- ReInit -------------------------------------------------------
91 InitStatus CbmCheckEvents::ReInit() { return kSUCCESS; }
92 
94  fEventSize =
95  new TH1F("fEventSize", "Event Size; # Digis; Counts", 1000, -0.5, 999.5);
96  fEventLength = new TH1F(
97  "fEventLength", "Event Length; time [ns]; Counts", 1000, -0.5, 999.5);
98  fEventsPerTS = new TH1F("fEventsPerTS",
99  "Events per time slice; # Events; Counts",
100  1000,
101  -0.5,
102  999.5);
103  fT0InEvent = new TH1F("fT0InEvent",
104  "Number of T0 digis in Event; # digis; Counts",
105  1000,
106  -0.5,
107  999.5);
108  fStsInEvent = new TH1F("fStsInEvent",
109  "Number of Sts digis in Event; # digis; Counts",
110  1000,
111  -0.5,
112  999.5);
113  fMuchInEvent = new TH1F("fMuchInEvent",
114  "Number of Much digis in Event; # digis; Counts",
115  1000,
116  -0.5,
117  999.5);
118  fTofInEvent = new TH1F("fTofInEvent",
119  "Number of Tof digis in Event; # digis; Counts",
120  1000,
121  -0.5,
122  999.5);
123  fT0DeltaT =
124  new TH1F("fT0DeltaT",
125  "Time diff between first and last T0 digi;dt [ns]; Counts",
126  1000,
127  -0.5,
128  999.5);
129  fStsDeltaT =
130  new TH1F("fStsDeltaT",
131  "Time diff between first and last Sts digi;dt [ns]; Counts",
132  1000,
133  -0.5,
134  999.5);
135  fMuchDeltaT =
136  new TH1F("fMuchDeltaT",
137  "Time diff between first and last Much digi;dt [ns]; Counts",
138  1000,
139  -0.5,
140  999.5);
141  fTofDeltaT =
142  new TH1F("fTofDeltaT",
143  "Time diff between first and last Tof digi;dt [ns]; Counts",
144  1000,
145  -0.5,
146  999.5);
147 
148  fEventsvsTS = new TH2F("fEventsvsTS",
149  "Nr. of events as fct. of TS",
150  10000,
151  -0.5,
152  9999.5,
153  1000,
154  -0.5,
155  999.5);
156  fLengthvsTS = new TProfile("fLengthvsTS",
157  "Length of events as fct. of TS",
158  10000,
159  -0.5,
160  9999.5,
161  -0.5,
162  999.5);
163 }
164 // ---- Exec ----------------------------------------------------------
165 void CbmCheckEvents::Exec(Option_t* /*option*/) {
166 
167  LOG_IF(info, fNrTs % 1000 == 0) << "Analysing TS " << fNrTs;
168 
169  Int_t nrEvents = fEvents->GetEntriesFast();
170 
171  fEventsPerTS->Fill(nrEvents);
172  fEventsvsTS->Fill(fNrTs, nrEvents);
173 
174  Int_t nrT0Digis = -1;
175  if (fT0DigiVec)
176  nrT0Digis = fT0DigiVec->size();
177  else if (fT0DigiArr)
178  nrT0Digis = fT0DigiArr->GetEntriesFast();
179  Int_t nrStsDigis = fDigiMan->GetNofDigis(ECbmModuleId::kSts);
180  Int_t nrMuchDigis = fDigiMan->GetNofDigis(ECbmModuleId::kMuch);
181  Int_t nrTofDigis = fDigiMan->GetNofDigis(ECbmModuleId::kTof);
182 
183  LOG(debug) << "Events: " << nrEvents;
184  LOG(debug) << "T0Digis: " << nrT0Digis;
185  LOG(debug) << "StsDigis: " << nrStsDigis;
186  LOG(debug) << "MuchDigis: " << nrMuchDigis;
187  LOG(debug) << "TofDigis: " << nrTofDigis;
188 
189  // Loop over all CbmEvents in the time slice
190  for (Int_t iEvent = 0; iEvent < nrEvents; iEvent++) {
191  CbmEvent* event = static_cast<CbmEvent*>(fEvents->At(iEvent));
192  fEventSize->Fill(event->GetNofData());
193  fEventLength->Fill(event->GetEndTime() - event->GetStartTime());
194  fLengthvsTS->Fill(fNrTs, event->GetEndTime() - event->GetStartTime(), 1);
195  AnalyseEvent(event);
196  }
197 
198  fNrTs++;
199 }
200 
202  // Loop over the the digis and extract the maximum time
203  // difference between the digis
205  GetTimeDiff<CbmStsDigi>(
207  GetTimeDiff<CbmMuchBeamTimeDigi>(
209  GetTimeDiff<CbmTofDigi>(
211 }
212 
213 template<class Digi>
215  TH1* deltaT,
216  TH1* size,
217  ECbmDataType dataType) {
218  Double_t startTime {1.e18};
219  Double_t stopTime {0.};
220  Int_t nDigis = event->GetNofData(dataType);
221  size->Fill(nDigis);
222  for (Int_t iDigi = 0; iDigi < nDigis; iDigi++) {
223  UInt_t index = event->GetIndex(dataType, iDigi);
224  const Digi* digi = fDigiMan->Get<Digi>(index);
225  assert(digi);
226  if (digi->GetTime() < startTime) startTime = digi->GetTime();
227  if (digi->GetTime() > stopTime) stopTime = digi->GetTime();
228  }
229  deltaT->Fill(stopTime - startTime);
230 }
231 
232 
233 void CbmCheckEvents::GetTimeDiffT0(CbmEvent* event, TH1* deltaT, TH1* size) {
234  Double_t startTime {1.e18};
235  Double_t stopTime {0.};
236  Int_t nDigis = event->GetNofData(ECbmDataType::kT0Digi);
237  size->Fill(nDigis);
238  for (Int_t iDigi = 0; iDigi < nDigis; iDigi++) {
239  UInt_t index = event->GetIndex(ECbmDataType::kT0Digi, iDigi);
240  //Double_t digiTime; (VF) not used
241  const CbmTofDigi* digi = nullptr;
242  if (fT0DigiVec)
243  digi = &(fT0DigiVec->at(index));
244  else if (fT0DigiArr)
245  digi = dynamic_cast<CbmTofDigi*>(fT0DigiArr->At(index));
246  assert(digi);
247  if (digi->GetTime() < startTime) startTime = digi->GetTime();
248  if (digi->GetTime() > stopTime) stopTime = digi->GetTime();
249  }
250  deltaT->Fill(stopTime - startTime);
251 }
252 
253 // ---- Finish --------------------------------------------------------
255  TFile* old = gFile;
256  TFile* outfile = TFile::Open("test2.root", "RECREATE");
257 
258  fEventSize->Write();
259  fEventLength->Write();
260  fEventsPerTS->Write();
261 
262  fT0InEvent->Write();
263  fStsInEvent->Write();
264  fMuchInEvent->Write();
265  fTofInEvent->Write();
266 
267  fT0DeltaT->Write();
268  fStsDeltaT->Write();
269  fMuchDeltaT->Write();
270  fTofDeltaT->Write();
271 
272  fEventsvsTS->Write();
273  fLengthvsTS->Write();
274 
275  outfile->Close();
276  delete outfile;
277 
278  gFile = old;
279 }
280 
ECbmDataType::kStsDigi
@ kStsDigi
CbmCheckEvents::fTofInEvent
TH1 * fTofInEvent
Definition: CbmCheckEvents.h:87
CbmCheckEvents::fMuchDeltaT
TH1 * fMuchDeltaT
Definition: CbmCheckEvents.h:91
CbmDigiManager::Init
InitStatus Init()
Initialisation.
Definition: CbmDigiManager.cxx:71
CbmCheckEvents::CbmCheckEvents
CbmCheckEvents()
Definition: CbmCheckEvents.cxx:34
CbmCheckEvents::GetTimeDiffT0
void GetTimeDiffT0(CbmEvent *, TH1 *, TH1 *)
Definition: CbmCheckEvents.cxx:233
CbmCheckEvents::fEventsvsTS
TH2 * fEventsvsTS
Definition: CbmCheckEvents.h:94
CbmCheckEvents::ReInit
virtual InitStatus ReInit()
Definition: CbmCheckEvents.cxx:91
CbmCheckEvents::Exec
virtual void Exec(Option_t *)
Definition: CbmCheckEvents.cxx:165
ECbmDataType::kTofDigi
@ kTofDigi
CbmDigiManager::GetNofDigis
static Int_t GetNofDigis(ECbmModuleId systemId)
Definition: CbmDigiManager.cxx:62
CbmDigiManager::UseMuchBeamTimeDigi
void UseMuchBeamTimeDigi(Bool_t choice=kTRUE)
Use CbmMuchBeamTimeDigi instead of CbmMuchDigi for MUCH.
Definition: CbmDigiManager.h:130
CbmTofDigi.h
ECbmModuleId::kTof
@ kTof
Time-of-flight Detector.
CbmCheckEvents::fEventLength
TH1 * fEventLength
Definition: CbmCheckEvents.h:82
CbmCheckEvents::fDigiMan
CbmDigiManager * fDigiMan
Definition: CbmCheckEvents.h:68
CbmCheckEvents::GetTimeDiff
void GetTimeDiff(CbmEvent *event, TH1 *, TH1 *, ECbmDataType type)
Definition: CbmCheckEvents.cxx:214
CbmDigiManager::IsPresent
static Bool_t IsPresent(ECbmModuleId systemId)
Presence of a digi branch.
Definition: CbmDigiManager.cxx:112
CbmDigiManager::Instance
static CbmDigiManager * Instance()
Static instance.
Definition: CbmDigiManager.h:93
CbmCheckEvents.h
CbmEvent.h
CbmStsDigi.h
CbmCheckEvents::fMuchInEvent
TH1 * fMuchInEvent
Definition: CbmCheckEvents.h:86
CbmCheckEvents::fStsInEvent
TH1 * fStsInEvent
Definition: CbmCheckEvents.h:85
CbmCheckEvents::fT0DigiArr
TClonesArray * fT0DigiArr
Definition: CbmCheckEvents.h:71
ECbmDataType::kT0Digi
@ kT0Digi
CbmCheckEvents::fLengthvsTS
TProfile * fLengthvsTS
Definition: CbmCheckEvents.h:95
CbmDigiManager::Get
const Digi * Get(Int_t index) const
Get a digi object.
Definition: CbmDigiManager.h:52
ECbmDataType
ECbmDataType
Definition: CbmDefs.h:76
CbmCheckEvents::SetParContainers
virtual void SetParContainers()
Definition: CbmCheckEvents.cxx:40
CbmMuchBeamTimeDigi.h
CbmTofDigi::GetTime
Double_t GetTime() const
Inherited from CbmDigi.
Definition: CbmTofDigi.h:111
CbmCheckEvents::fTofDeltaT
TH1 * fTofDeltaT
Definition: CbmCheckEvents.h:92
CbmCheckEvents::fT0DigiVec
const std::vector< CbmTofDigi > * fT0DigiVec
Interface to digi data.
Definition: CbmCheckEvents.h:70
CbmCheckEvents::fEventsPerTS
TH1 * fEventsPerTS
Definition: CbmCheckEvents.h:83
CbmCheckEvents::fStsDeltaT
TH1 * fStsDeltaT
Definition: CbmCheckEvents.h:90
ClassImp
ClassImp(CbmConverterManager) InitStatus CbmConverterManager
Definition: CbmConverterManager.cxx:12
CbmCheckEvents::AnalyseEvent
void AnalyseEvent(CbmEvent *event)
Definition: CbmCheckEvents.cxx:201
CbmTofDigi
Data class for expanded digital TOF information.
Definition: CbmTofDigi.h:38
CbmCheckEvents
Definition: CbmCheckEvents.h:27
CbmCheckEvents::fNrTs
Int_t fNrTs
Definition: CbmCheckEvents.h:78
CbmDigiManager.h
CbmCheckEvents::CreateHistos
void CreateHistos()
Definition: CbmCheckEvents.cxx:93
CbmEvent
Class characterising one event by a collection of links (indices) to data objects,...
Definition: CbmEvent.h:30
CbmCheckEvents::Finish
virtual void Finish()
Definition: CbmCheckEvents.cxx:254
ECbmModuleId::kMuch
@ kMuch
Muon detection system.
CbmCheckEvents::fEventSize
TH1 * fEventSize
Definition: CbmCheckEvents.h:81
CbmCheckEvents::fT0InEvent
TH1 * fT0InEvent
Definition: CbmCheckEvents.h:84
ECbmModuleId::kSts
@ kSts
Silicon Tracking System.
CbmCheckEvents::~CbmCheckEvents
~CbmCheckEvents()
Definition: CbmCheckEvents.cxx:37
CbmCheckEvents::Init
virtual InitStatus Init()
Definition: CbmCheckEvents.cxx:52
ECbmDataType::kMuchDigi
@ kMuchDigi
CbmCheckEvents::fEvents
TClonesArray * fEvents
Definition: CbmCheckEvents.h:72
CbmCheckEvents::fT0DeltaT
TH1 * fT0DeltaT
Definition: CbmCheckEvents.h:89
CbmDefs.h