R3BROOT
R3B analysis software
Loading...
Searching...
No Matches
R3BNeulandCal2Hit.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
14#include "R3BNeulandCal2Hit.h"
15#include "FairLogger.h"
16#include "FairRootManager.h"
17#include "FairRunAna.h"
18#include "FairRuntimeDb.h"
19#include "R3BEventHeader.h"
20#include "R3BNeulandCommon.h"
21#include "R3BNeulandHitPar.h"
22#include "R3BTCalPar.h"
23#include "TCanvas.h"
24#include "TClonesArray.h"
25#include "TGraph2D.h"
26#include "TH1F.h"
27#include "TH2F.h"
28#include "TMath.h"
29#include "TVector.h"
30#include <R3BLogger.h>
31
32#include <array>
33
34namespace Neuland = R3B::Neuland;
35
36R3BNeulandCal2Hit::R3BNeulandCal2Hit(const char* name, const Int_t iVerbose)
37 : FairTask(name, iVerbose)
38 , fEventHeader(nullptr)
39 , fCalData("NeulandCalData")
40 , fHits("NeulandHits")
41 , fFirstPlaneHorizontal(true)
42 , fDistanceToTarget(Neuland::NaN)
43 , fGlobalTimeOffset(Neuland::NaN)
44 , fEnergyCutoff(Neuland::NaN)
45 , fWalkEnabled(kFALSE)
46{
47}
48
50{
51 auto ioman = FairRootManager::Instance();
52 if (ioman == nullptr)
53 {
54 throw std::runtime_error("R3BNeulandCal2Hit: No FairRootManager");
55 }
56
57 fEventHeader = dynamic_cast<R3BEventHeader*>(ioman->GetObject("EventHeader."));
58 if (fEventHeader == nullptr)
59 {
60 throw std::runtime_error("R3BNeulandCal2Hit: No R3BEventHeader");
61 }
62
63 fCalData.Init();
64 fHits.Init();
65
66 SetParameter();
67 return kSUCCESS;
68}
69
71{
72 fPar = dynamic_cast<R3BNeulandHitPar*>(FairRuntimeDb::instance()->getContainer("NeulandHitPar"));
73}
74
75void R3BNeulandCal2Hit::SetParameter()
76{
77 fParMap.clear();
78
79 if (std::isnan(fDistanceToTarget))
80 {
81 fDistanceToTarget = fPar->GetDistanceToTarget();
82 }
83 else
84 {
85 fPar->SetDistanceToTarget(fDistanceToTarget);
86 }
87
88 if (std::isnan(fGlobalTimeOffset))
89 {
90 fGlobalTimeOffset = fPar->GetGlobalTimeOffset();
91 }
92 else
93 {
94 fPar->SetGlobalTimeOffset(fGlobalTimeOffset);
95 }
96
97 if (std::isnan(fEnergyCutoff))
98 {
99 fEnergyCutoff = fPar->GetEnergyCutoff();
100 }
101 else
102 {
103 fPar->SetEnergyCutoff(fEnergyCutoff);
104 }
105
106 fNumberOfPlanes = fPar->GetNumberOfPlanes();
107 fDistancesToFirstPlane = fPar->GetDistancesToFirstPlane();
108 fAttenuationValues.resize(fNumberOfPlanes * Neuland::BarsPerPlane, 0.);
109 const auto nPars = fPar->GetNumModulePar();
110
111 for (auto i = 0; i < nPars; i++)
112 {
113 const auto id = fPar->GetModuleParAt(i)->GetModuleId() - 1;
114 fParMap[id] = *fPar->GetModuleParAt(i);
115 fAttenuationValues[id] = exp(Neuland::TotalBarLength / fParMap[id].GetLightAttenuationLength());
116 }
117
118 LOG(info) << "R3BNeulandCal2Hit::SetParameter : Number of Parameters: " << fPar->GetNumModulePar();
119}
120
122{
124 SetParameter();
125 return kSUCCESS;
126}
127
129{
130
131 if (fVerbose > 0 && ++fEventNumber % 10000 == 0)
132 std::cout << "\rR3BNeulandCal2Hit " << fEventNumber << " Events converted." << std::flush;
133
134 fHits.Reset();
135 fHitMap.clear();
136
137 auto calData = fCalData.Retrieve();
138
139 const auto start = fEventHeader->GetTStart();
140 const bool beam = !std::isnan(start);
141
142 for (auto calDataPtrPtr = calData.begin(); calDataPtrPtr != calData.end(); calDataPtrPtr++)
143 {
144 const auto barID = (*calDataPtrPtr)->GetBarId() - 1;
145
146 R3BLOG(debug2,
147 fmt::format("Input calData barID: {}, side: {}, time: {}, fTriggerTime: {}, fQdc: {}",
148 (*calDataPtrPtr)->GetBarId(),
149 (*calDataPtrPtr)->GetSide(),
150 (*calDataPtrPtr)->GetTime(),
151 (*calDataPtrPtr)->GetTriggerTime(),
152 (*calDataPtrPtr)->GetQdc()));
153
154 if (fParMap.find(barID) == fParMap.end())
155 continue; // We do not have parameters for this module
156
157 if (fHitMap.find(barID) == fHitMap.end()) // TODO what about multi hit?
158 {
159 // This is the first PMT of this bar, store the pointer and go to next cal entry
160 fHitMap[barID] = *calDataPtrPtr;
161 continue;
162 }
163
164 const auto& parameter = fParMap[barID];
165
166 std::array<const R3BNeulandCalData*, 2> cal;
167 if ((*calDataPtrPtr)->GetSide() == 1)
168 cal = { *calDataPtrPtr, fHitMap[barID] };
169 else
170 cal = { fHitMap[barID], *calDataPtrPtr };
171
172 const std::array<double, 2> qdc = { std::max(cal[0]->GetQdc() - parameter.GetPedestal(1), 1.),
173 std::max(cal[1]->GetQdc() - parameter.GetPedestal(2), 1.) };
174
175 const std::array<Double_t, 2> unsatEnergy = {
176 GetUnsaturatedEnergy(qdc[0], parameter.GetEnergyGain(1), parameter.GetPMTSaturation(1)),
177 GetUnsaturatedEnergy(qdc[1], parameter.GetEnergyGain(2), parameter.GetPMTSaturation(2))
178 };
179
180 const auto energy = TMath::Sqrt(fAttenuationValues[barID] * unsatEnergy[0] * unsatEnergy[1]);
181
182 if (energy < fEnergyCutoff)
183 continue;
184
185 std::array<Double_t, 2> tdc;
186
187 if (std::isnan(cal[0]->GetTriggerTime()) || std::isnan(cal[0]->GetTriggerTime()))
188 {
189 tdc = { cal[0]->GetTime() + parameter.GetTimeOffset(1), cal[1]->GetTime() + parameter.GetTimeOffset(2) };
190 }
191 else
192 {
193 tdc = { cal[0]->GetTime() - cal[0]->GetTriggerTime() + parameter.GetTimeOffset(1),
194 cal[1]->GetTime() - cal[1]->GetTriggerTime() + parameter.GetTimeOffset(2) };
195 // fmt::print("time0: {}, trigger0: {}, timeOffset: {}\n",
196 // cal[0]->GetTime(),
197 // cal[0]->GetTriggerTime(),
198 // parameter.GetTimeOffset(1));
199 // fmt::print("time1: {}, trigger1: {}, timeOffset: {}\n",
200 // cal[1]->GetTime(),
201 // cal[1]->GetTriggerTime(),
202 // parameter.GetTimeOffset(2));
203 }
204
205 // FIXME: this should be done in Mapped2Cal
206 // In Cal2Hit the difference between all bars should be checked
207 if (tdc[0] - tdc[1] < -0.5 * Neuland::MaxCalTime)
208 tdc[1] -= Neuland::MaxCalTime;
209 else if (tdc[0] - tdc[1] > 0.5 * Neuland::MaxCalTime)
210 tdc[0] -= Neuland::MaxCalTime;
211
212 auto time = (tdc[0] + tdc[1]) * 0.5 - fGlobalTimeOffset;
213
214 // NOTE: Walk correction has already been done in map2cal. Why here again?
215 if (fWalkEnabled)
216 time = time + WalkCorrection(energy);
217
218 if (beam)
219 {
220 // the shift is to get fmod to work as indented: 4 peaks -> 1 peak w/o stray data (e.g. at 5 * 2048)
221 // fmod gives the time between 0 and 5*2048
222 // remainder() gives the time between -0.5*5*2048 and +0.5*5*2048
223 fmt::print("time: {}, start: {}\n", time, start);
224 time = remainder(time - start, Neuland::MaxCalTime);
225 }
226 else
227 {
228 time = Neuland::NaN;
229 }
230
231 const auto plane = Neuland::ModuleID2PlaneID(barID);
232 const auto bar = (barID) % 50;
233
234 TVector3 pos;
235 TVector3 pixel;
236
237 // NOTE: First plane horizontal or not should be decided by IsPlaneHorizontal
238 if (Neuland::IsPlaneIDHorizontal(plane) == fFirstPlaneHorizontal)
239 {
240 pos[0] = parameter.GetEffectiveSpeed() * (tdc[1] - tdc[0]);
241 pos[1] = (bar + 0.5 - Neuland::BarsPerPlane * 0.5) * Neuland::BarSize_XY;
242
243 pixel[0] = std::min(std::max(0., pos[0] / 5. + 25), 49.);
244 pixel[1] = bar;
245 }
246 else
247 {
248 pos[0] = (bar + 0.5 - Neuland::BarsPerPlane * 0.5) * Neuland::BarSize_XY;
249 pos[1] = parameter.GetEffectiveSpeed() * (tdc[1] - tdc[0]);
250
251 pixel[0] = bar;
252 pixel[1] = std::min(std::max(0., pos[1] / 5. + 25), 49.);
253 }
254
255 pos[2] = fDistanceToTarget + fDistancesToFirstPlane[plane];
256 pixel[2] = plane;
257
258 auto hit = R3BNeulandHit{ barID + 1, tdc[0], tdc[1], time, unsatEnergy[0], unsatEnergy[1], energy, pos, pixel };
259 R3BLOG(debug1, fmt::format("Neuland hit: {}", hit));
260 fHits.Insert(hit);
261 }
262}
263
265{
266 if (fVerbose > 0)
267 std::cout << "\rR3BNeulandCal2Hit " << fEventNumber << " Events converted." << std::endl;
268}
269
270Double_t R3BNeulandCal2Hit::GetUnsaturatedEnergy(const Int_t qdc, const Double_t gain, const Double_t saturation) const
271{
272 return qdc / (gain - saturation * qdc);
273}
274
275Double_t R3BNeulandCal2Hit::WalkCorrection(Double_t x)
276{
277 Double_t y = 0;
278
279 // if (x<=12.) y=24.78;
280
281 if (x <= 12.)
282 y = 35. - 7.1897575 * log(x) + 1.2 * log(x) * log(x) + 0.23616;
283
284 if (x > 12.)
285 y = 25.5 - 0.055 * x;
286
287 return 24.78 - y;
288}
289
#define R3BLOG(severity, x)
Definition R3BLogger.h:35
ClassImp(R3B::Neuland::Cal2HitPar)
InitStatus Init() override
R3BNeulandCal2Hit(const char *name="R3BNeulandCal2Hit", const Int_t iVerbose=0)
void SetParContainers() override
InitStatus ReInit() override
void FinishTask() override
void Exec(Option_t *) override
Double_t GetDistanceToTarget(Int_t plane=0) const
Simulation of NeuLAND Bar/Paddle.
constexpr auto BarsPerPlane
constexpr auto ModuleID2PlaneID(int moduleID) -> int
constexpr auto BarSize_XY
constexpr auto TotalBarLength
constexpr auto NaN
constexpr auto IsPlaneIDHorizontal(int plane_id) -> bool
constexpr auto MaxCalTime