R3BROOT
R3B analysis software
Loading...
Searching...
No Matches
R3BNeulandCal2HitPar.cxx
Go to the documentation of this file.
1/******************************************************************************
2 * Copyright (C) 2019 GSI Helmholtzzentrum für Schwerionenforschung GmbH *
3 * Copyright (C) 2019-2025 Members of R3B Collaboration *
4 * *
5 * This software is distributed under the terms of the *
6 * GNU General Public Licence (GPL) version 3, *
7 * copied verbatim in the file "LICENSE". *
8 * *
9 * In applying this license GSI does not waive the privileges and immunities *
10 * granted to it by virtue of its status as an Intergovernmental Organization *
11 * or submit itself to any jurisdiction. *
12 ******************************************************************************/
13
15#include <FairRootManager.h>
16#include <FairTask.h>
17#include <Rtypes.h>
18#include <RtypesCore.h>
19#include <TDirectory.h>
20#include <TString.h>
21#include <cmath>
22#include <fairlogger/Logger.h>
23#include <iostream>
24#include <ostream>
25
26#include "R3BEventHeader.h"
27#include "R3BNeulandCalData.h"
28#include "R3BNeulandCommon.h"
32#include "R3BNeulandHitPar.h"
33
34#include "FairLogger.h"
35#include "FairRuntimeDb.h"
36
37#include "TCanvas.h"
38#include "TClonesArray.h"
39#include "TH2F.h"
40#include "TROOT.h"
41
42using namespace R3B::Neuland; // NOLINT
43
45 : R3BNeulandCal2HitPar("R3BNeulandCal2HitPar", 1)
46{
47}
48
49R3BNeulandCal2HitPar::R3BNeulandCal2HitPar(const char* name, const Int_t iVerbose)
50 : FairTask(name, iVerbose)
51 , fHitCalEngine(new Calibration::HitCalibrationEngine())
52 , fCosmicTracker(new Calibration::CosmicTracker())
53 , fNeulandHitPar(nullptr)
54 , fCalNeuland(nullptr)
55 , fMappedLos(nullptr)
56 , fCosmicTpat(0)
57 , fSavePlots(kFALSE)
58{
59}
60
62{
63 FairRootManager* mgr = FairRootManager::Instance();
64
65 if (!mgr)
66 {
67 LOG(fatal) << "R3BNeulandCal2HitPar::Init: FairRootManager not found";
68 return kFATAL;
69 }
70
71 fEventHeader = static_cast<R3BEventHeader*>(mgr->GetObject("EventHeader."));
72
73 if (!fEventHeader)
74 {
75 LOG(info) << "R3BNeulandCal2HitPar::Init Branch EventHeader not found, its ok";
76 // return kFATAL;
77 }
78
79 fCalNeuland = static_cast<TClonesArray*>(mgr->GetObject("NeulandCalData"));
80
81 if (!fCalNeuland)
82 {
83 LOG(fatal) << "R3BNeulandCal2HitPar::Init: Branch NeulandCalData not found";
84 return kFATAL;
85 }
86
87 fNeulandHitPar = static_cast<R3BNeulandHitPar*>(FairRuntimeDb::instance()->getContainer("NeulandHitPar"));
88 LOG(info) << "R3BNeulandCal2HitPar::Init: Number of Hit-Paramteres found: " << fNeulandHitPar->GetNumModulePar();
89
90 LOG(info) << "R3BNeulandCal2HitPar::Init: Initializing NeulandHitCalibrationEngine.";
92
94 fCosmicTracker->SetDistances(fNeulandHitPar->GetDistancesToFirstPlane());
95
96 return kSUCCESS;
97}
98
99void R3BNeulandCal2HitPar::Exec(Option_t* option)
100{
101 if (++fEventNumber % 100000 == 0 && fVerbose) // 10000
102 {
103 const auto msg = TString::Format(
104 "R3BNeulandCal2HitPar::Exec : Event: %10d, accepted Events: %10d", fEventNumber, fAcceptedEventNumber);
105 std::cout << msg << "\r" << std::flush;
106
107 std::cout << " cnts " << fIgorcnt0 << " ---- " << fIgorcnt1 << " ---- " << fIgorcnt2 << " " << std::endl;
108 }
109
110 // Skip if this is not a pure cosmic event
111 if (!IsCosmicEvent())
112 {
113 return;
114 }
115
116 ++fIgorcnt0;
117
118 LOG(debug) << "R3BNeulandCal2HitPar::Exec: Event " << fEventNumber - 1;
119 const auto nItems = fCalNeuland->GetEntriesFast();
120
121 if (nItems < 8)
122 {
123 LOG(debug) << " Event cannot be used: too few signals : " << nItems << "!";
124 return;
125 }
126
127 fHitCalEngine->Reset();
128 fCosmicTracker->Reset();
129
130 Int_t addedPoints = 0;
131
132 for (Int_t i = 0; i < nItems; i++)
133 {
134 const auto pmt = static_cast<R3BNeulandCalData*>(fCalNeuland->At(i));
135
136 const auto id = pmt->GetBarId() - 1;
137 const auto side = pmt->GetSide() - 1;
138
139 if (std::isnan(pmt->GetTriggerTime()))
140 {
141 fHitCalEngine->Set(id, side, pmt->GetTime(), pmt->GetQdc());
142 }
143 else
144 {
145 fHitCalEngine->Set(id, side, pmt->GetTime() - pmt->GetTriggerTime(), pmt->GetQdc());
146 }
147
148 if (fHitCalEngine->IsValid(id))
149 {
150 ++addedPoints;
151 fCosmicTracker->AddPoint(id, fHitCalEngine->GetPosition(id));
152 }
153 }
154
155 ++fIgorcnt1;
156
157 if (addedPoints < 6)
158 {
159 LOG(debug) << " Event cannot be used: too few Points : " << addedPoints << " (" << nItems << ")!";
160 return;
161 }
162
163 ++fIgorcnt2;
164
165 const auto& cosmicTrack = fCosmicTracker->GetTrack();
166
167 if (cosmicTrack.Interactions.size() == 0)
168 {
169 LOG(debug) << " Getting Cosmic Track : Failure!";
170 return;
171 }
172
173 LOG(debug) << " Getting Cosmic Track: Success!";
175
176 LOG(debug) << " Adding data to calibration";
177 fHitCalEngine->Add(cosmicTrack, fEventNumber);
178}
179
181{
182 LOG(info) << "R3BNeulandCal2HitPar::FinishTask: " << "Saved " << fAcceptedEventNumber << " Events.";
183
184 LOG(info) << "R3BNeulandCal2HitPar::FinishTask: " << "Starting Neuland Hit Calibration with "
185 << fAcceptedEventNumber << " Events.";
186
187 const auto batchMode = gROOT->IsBatch();
188 gROOT->SetBatch(kTRUE);
189 const auto defdir = gDirectory;
190
191 TDirectory* neulandDir = nullptr;
192 if (fSavePlots)
193 neulandDir = defdir->mkdir("NeuLAND");
194
195 auto parameters = fHitCalEngine->Calibrate(neulandDir);
196
197 defdir->cd();
198 gROOT->SetBatch(batchMode);
199
200 fNeulandHitPar->GetListOfModulePar()->Clear();
201
202 auto maxThreshold = 0.;
203 for (auto& parameter : parameters)
204 {
205 fNeulandHitPar->AddModulePar(new R3BNeulandHitModulePar(parameter));
206
207 for (auto side = 1; side <= 2; ++side)
208 {
209 const auto threshold =
210 parameter.GetPMTThreshold(side) * exp(TotalBarLength / parameter.GetLightAttenuationLength());
211 if (threshold > maxThreshold)
212 maxThreshold = threshold;
213 }
214 }
215
216 fNeulandHitPar->SetEnergyCutoff(maxThreshold);
217 LOG(info) << TString::Format("R3BNeulandCal2HitPar::FinishTask: Recommended Minimum Energy Cutoff : >%4.2f MeV\n",
218 maxThreshold);
219
220 fNeulandHitPar->setChanged();
221
222 LOG(info) << "R3BNeulandCal2HitPar::FinishTask: " << "Number of calibrated Bars: "
223 << fNeulandHitPar->GetNumModulePar();
224}
225
227
ClassImp(R3B::Neuland::Cal2HitPar)
std::unique_ptr< R3B::Neuland::Calibration::CosmicTracker > fCosmicTracker
R3BEventHeader * fEventHeader
R3BNeulandHitPar * fNeulandHitPar
std::unique_ptr< R3B::Neuland::Calibration::HitCalibrationEngine > fHitCalEngine
virtual void Exec(Option_t *option)
virtual InitStatus Init()
auto GetBarId() const -> int
Simulation of NeuLAND Bar/Paddle.
constexpr auto TotalBarLength