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