CbmRoot
Loading...
Searching...
No Matches
CbmCheckEvents.cxx
Go to the documentation of this file.
1/* Copyright (C) 2019-2021 GSI Helmholtzzentrum fuer Schwerionenforschung, Darmstadt
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Florian Uhlig [committer] */
4
5#include "CbmCheckEvents.h"
6
7//#include "CbmTofDigi.h"
8#include "CbmDefs.h"
9#include "CbmDigiManager.h"
10#include "CbmEvent.h"
11#include "CbmMuchBeamTimeDigi.h"
12#include "CbmStsDigi.h"
13#include "CbmTofDigi.h"
14
15#include "FairRootManager.h"
16#include "FairRunOnline.h"
17#include <Logger.h>
18
19#include "TClonesArray.h"
20#include "TH1.h"
21#include "TH2.h"
22#include "THttpServer.h"
23#include "TProfile.h"
24#include <TDirectory.h>
25#include <TFile.h>
26
27#include <iomanip>
28#include <typeinfo>
29
30using std::fixed;
31using std::setprecision;
32
33// ---- Default constructor -------------------------------------------
34CbmCheckEvents::CbmCheckEvents() : FairTask("CbmCheckEvents") {}
35
36// ---- Destructor ----------------------------------------------------
38
39// ---- Initialisation ----------------------------------------------
41{
42 // Load all necessary parameter containers from the runtime data base
43 /*
44 FairRunAna* ana = FairRunAna::Instance();
45 FairRuntimeDb* rtdb=ana->GetRuntimeDb();
46
47 <CbmCheckEventsDataMember> = (<ClassPointer>*)
48 (rtdb->getContainer("<ContainerName>"));
49 */
50}
51
52// ---- Init ----------------------------------------------------------
54{
55
56 // Get a handle from the IO manager
57 FairRootManager* ioman = FairRootManager::Instance();
58
59 // DigiManager
62 fDigiMan->Init();
63
64 // Get a pointer to the previous already existing data level
65
66 fBmonDigiVec = ioman->InitObjectAs<std::vector<CbmTofDigi> const*>("BmonDigi");
67 if (!fBmonDigiVec) {
68 fBmonDigiArr = dynamic_cast<TClonesArray*>(ioman->GetObject("BmonDigi"));
69 if (!fBmonDigiArr) { LOG(fatal) << "No TClonesArray with Bmon digis found."; }
70 }
71
72 if (!fDigiMan->IsPresent(ECbmModuleId::kSts)) { LOG(info) << "No TClonesArray with STS digis found."; }
73
74 if (!fDigiMan->IsPresent(ECbmModuleId::kMuch)) { LOG(info) << "No TClonesArray with MUCH digis found."; }
75
76 if (!fDigiMan->IsPresent(ECbmModuleId::kTof)) { LOG(info) << "No TClonesArray with TOF digis found."; }
77
78 fEvents = dynamic_cast<TClonesArray*>(ioman->GetObject("CbmEvent"));
79 if (nullptr == fEvents) {
80
81 if (nullptr != (ioman->GetObject("CbmEvent"))) {
82 LOG(error) << "Got pointer of type" << typeid(ioman->GetObject("CbmEvent")).name();
83 auto& tmp = *(ioman->GetObject("CbmEvent"));
84 LOG(error) << "Got Object of type" << typeid(tmp).name();
85 } // if( nullptr != (ioman->GetObject("CbmEvent") )
86 LOG(fatal) << "No TClonesArray with events found.";
87 } // if (nullptr == fEvents)
88
89
91
92 return kSUCCESS;
93}
94
95// ---- ReInit -------------------------------------------------------
96InitStatus CbmCheckEvents::ReInit() { return kSUCCESS; }
97
99{
100 fEventSize = new TH1F("fEventSize", "Event Size; # Digis; Counts", 1000, -0.5, 999.5);
101 fEventLength = new TH1F("fEventLength", "Event Length; time [ns]; Counts", 1000, -0.5, 999.5);
102 fEventsPerTS = new TH1F("fEventsPerTS", "Events per time slice; # Events; Counts", 1000, -0.5, 999.5);
103 fBmonInEvent = new TH1F("fBmonInEvent", "Number of Bmon digis in Event; # digis; Counts", 1000, -0.5, 999.5);
104 fStsInEvent = new TH1F("fStsInEvent", "Number of Sts digis in Event; # digis; Counts", 1000, -0.5, 999.5);
105 fMuchInEvent = new TH1F("fMuchInEvent", "Number of Much digis in Event; # digis; Counts", 1000, -0.5, 999.5);
106 fTofInEvent = new TH1F("fTofInEvent", "Number of Tof digis in Event; # digis; Counts", 1000, -0.5, 999.5);
108 new TH1F("fBmonDeltaT", "Time diff between first and last Bmon digi;dt [ns]; Counts", 1000, -0.5, 999.5);
109 fStsDeltaT = new TH1F("fStsDeltaT", "Time diff between first and last Sts digi;dt [ns]; Counts", 1000, -0.5, 999.5);
111 new TH1F("fMuchDeltaT", "Time diff between first and last Much digi;dt [ns]; Counts", 1000, -0.5, 999.5);
112 fTofDeltaT = new TH1F("fTofDeltaT", "Time diff between first and last Tof digi;dt [ns]; Counts", 1000, -0.5, 999.5);
113
114 fEventsvsTS = new TH2F("fEventsvsTS", "Nr. of events as fct. of TS", 10000, -0.5, 9999.5, 1000, -0.5, 999.5);
115 fLengthvsTS = new TProfile("fLengthvsTS", "Length of events as fct. of TS", 10000, -0.5, 9999.5, -0.5, 999.5);
116}
117// ---- Exec ----------------------------------------------------------
118void CbmCheckEvents::Exec(Option_t* /*option*/)
119{
120
121 LOG_IF(info, fNrTs % 1000 == 0) << "Analysing TS " << fNrTs;
122
123 Int_t nrEvents = fEvents->GetEntriesFast();
124
125 fEventsPerTS->Fill(nrEvents);
126 fEventsvsTS->Fill(fNrTs, nrEvents);
127
128 Int_t nrBmonDigis = -1;
129 if (fBmonDigiVec) nrBmonDigis = fBmonDigiVec->size();
130 else if (fBmonDigiArr)
131 nrBmonDigis = fBmonDigiArr->GetEntriesFast();
132 Int_t nrStsDigis = fDigiMan->GetNofDigis(ECbmModuleId::kSts);
133 Int_t nrMuchDigis = fDigiMan->GetNofDigis(ECbmModuleId::kMuch);
134 Int_t nrTofDigis = fDigiMan->GetNofDigis(ECbmModuleId::kTof);
135
136 LOG(debug) << "Events: " << nrEvents;
137 LOG(debug) << "BmonDigis: " << nrBmonDigis;
138 LOG(debug) << "StsDigis: " << nrStsDigis;
139 LOG(debug) << "MuchDigis: " << nrMuchDigis;
140 LOG(debug) << "TofDigis: " << nrTofDigis;
141
142 // Loop over all CbmEvents in the time slice
143 for (Int_t iEvent = 0; iEvent < nrEvents; iEvent++) {
144 CbmEvent* event = dynamic_cast<CbmEvent*>(fEvents->At(iEvent));
145 fEventSize->Fill(event->GetNofData());
146 fEventLength->Fill(event->GetEndTime() - event->GetStartTime());
147 fLengthvsTS->Fill(fNrTs, event->GetEndTime() - event->GetStartTime(), 1);
148 AnalyseEvent(event);
149 }
150
151 fNrTs++;
152}
153
155{
156 // Loop over the the digis and extract the maximum time
157 // difference between the digis
162}
163
164template<class Digi>
165void CbmCheckEvents::GetTimeDiff(CbmEvent* event, TH1* deltaT, TH1* size, ECbmDataType dataType)
166{
167 Double_t startTime {1.e18};
168 Double_t stopTime {0.};
169 Int_t nDigis = event->GetNofData(dataType);
170 size->Fill(nDigis);
171 for (Int_t iDigi = 0; iDigi < nDigis; iDigi++) {
172 UInt_t index = event->GetIndex(dataType, iDigi);
173 const Digi* digi = fDigiMan->Get<Digi>(index);
174 assert(digi);
175 if (digi->GetTime() < startTime) startTime = digi->GetTime();
176 if (digi->GetTime() > stopTime) stopTime = digi->GetTime();
177 }
178 deltaT->Fill(stopTime - startTime);
179}
180
181
182void CbmCheckEvents::GetTimeDiffBmon(CbmEvent* event, TH1* deltaT, TH1* size)
183{
184 Double_t startTime {1.e18};
185 Double_t stopTime {0.};
186 Int_t nDigis = event->GetNofData(ECbmDataType::kBmonDigi);
187 size->Fill(nDigis);
188 for (Int_t iDigi = 0; iDigi < nDigis; iDigi++) {
189 UInt_t index = event->GetIndex(ECbmDataType::kBmonDigi, iDigi);
190 //Double_t digiTime; (VF) not used
191 const CbmTofDigi* digi = nullptr;
192 if (fBmonDigiVec) digi = &(fBmonDigiVec->at(index));
193 else if (fBmonDigiArr)
194 digi = dynamic_cast<CbmTofDigi*>(fBmonDigiArr->At(index));
195 assert(digi);
196 if (digi->GetTime() < startTime) startTime = digi->GetTime();
197 if (digi->GetTime() > stopTime) stopTime = digi->GetTime();
198 }
199 deltaT->Fill(stopTime - startTime);
200}
201
202// ---- Finish --------------------------------------------------------
204{
205 TFile* oldFile = gFile;
206 TDirectory* oldDir = gDirectory;
207
208 TFile* outfile = TFile::Open("test2.root", "RECREATE");
209
210 fEventSize->Write();
211 fEventLength->Write();
212 fEventsPerTS->Write();
213
214 fBmonInEvent->Write();
215 fStsInEvent->Write();
216 fMuchInEvent->Write();
217 fTofInEvent->Write();
218
219 fBmonDeltaT->Write();
220 fStsDeltaT->Write();
221 fMuchDeltaT->Write();
222 fTofDeltaT->Write();
223
224 fEventsvsTS->Write();
225 fLengthvsTS->Write();
226
227 outfile->Close();
228 delete outfile;
229
230 gFile = oldFile;
231 gDirectory = oldDir;
232}
233
ClassImp(CbmConverterManager)
ECbmDataType
Definition CbmDefs.h:90
@ kTof
Time-of-flight Detector.
@ kSts
Silicon Tracking System.
@ kMuch
Muon detection system.
static constexpr size_t size()
Definition KfSimdPseudo.h:2
CbmDigiManager * fDigiMan
TProfile * fLengthvsTS
const std::vector< CbmTofDigi > * fBmonDigiVec
Interface to digi data.
TClonesArray * fEvents
virtual InitStatus Init()
TClonesArray * fBmonDigiArr
virtual void Exec(Option_t *)
virtual InitStatus ReInit()
void GetTimeDiff(CbmEvent *event, TH1 *, TH1 *, ECbmDataType type)
void GetTimeDiffBmon(CbmEvent *, TH1 *, TH1 *)
virtual void SetParContainers()
void AnalyseEvent(CbmEvent *event)
virtual void Finish()
static Int_t GetNofDigis(ECbmModuleId systemId)
static Bool_t IsPresent(ECbmModuleId systemId)
Presence of a digi branch.
void UseMuchBeamTimeDigi(Bool_t)
Use CbmMuchBeamTimeDigi instead of CbmMuchDigi for MUCH.
InitStatus Init()
Initialisation.
const Digi * Get(Int_t index) const
Get a digi object.
static CbmDigiManager * Instance()
Static instance.
Class characterising one event by a collection of links (indices) to data objects,...
Definition CbmEvent.h:34
Data class for expanded digital TOF information.
Definition CbmTofDigi.h:47
double GetTime() const
Inherited from CbmDigi.
Definition CbmTofDigi.h:131