R3BROOT
R3B analysis software
Loading...
Searching...
No Matches
R3BNeulandTacquilaMapped2Cal.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// ------------------------------------------------------------
15// ----- R3BNeulandTacquilaMapped2Cal -----
16// ----- Created 22-04-2014 by D.Kresan -----
17// ------------------------------------------------------------
18
20#include "FairLogger.h"
21#include "FairRootManager.h"
22#include "FairRunAna.h"
23#include "FairRuntimeDb.h"
24#include "R3BEventHeader.h"
25#include "R3BNeulandCalData.h"
26#include "R3BNeulandQCalPar.h"
28#include "R3BTCalEngine.h"
29#include "R3BTCalPar.h"
30#include "TClonesArray.h"
31#include "TH1F.h"
32#include "TMath.h"
33
34#define planes fNofPMTs / 100
35#define toID(x, y, z) (((x - 1) * 50 + (y - 1)) * 2 + (z - 1))
36
37Double_t wlk(Double_t x)
38{
39 Double_t y = 0;
40
41 Double_t par1 = 1500.; // +-0.2238
42 Double_t par2 = 0.00075; //+-2.355e-05
43 y = par1 * TMath::Power(x, par2) - (par1 * TMath::Power(400., par2)); // Michael's
44
45 // y=2.29083*log(x)-0.0870157*log(x)*log(x)-4.57824; // mine
46
47 return y;
48 // return 0.;
49}
50
52 : FairTask("NeulandTacquilaMapped2Cal", 1)
53 , fNEvents(0)
54 , fPulserMode(kFALSE)
55 , fWalkEnabled(kTRUE)
56 , fRawHit(NULL)
57 , fPmt(new TClonesArray("R3BNeulandCalData"))
58 , fNPmt(0)
59 , fTcalPar(NULL)
60 , fQCalPar(NULL)
61 , fTrigger(-1)
62 , fMap17Seen()
63 , fMapStopTime()
64 , fMapStopClock()
65 , fMapQdcOffset()
66 , fClockFreq(1. / TACQUILA_CLOCK_MHZ * 1000.)
67{
68}
69
71 : FairTask(name, iVerbose)
72 , fNEvents(0)
73 , fPulserMode(kFALSE)
74 , fWalkEnabled(kTRUE)
75 , fRawHit(NULL)
76 , fPmt(new TClonesArray("R3BNeulandCalData"))
77 , fNPmt(0)
78 , fTcalPar(NULL)
79 , fQCalPar(NULL)
80 , fTrigger(-1)
81 , fMap17Seen()
82 , fMapStopTime()
83 , fMapStopClock()
84 , fMapQdcOffset()
85 , fClockFreq(1. / TACQUILA_CLOCK_MHZ * 1000.)
86{
87}
88
90{
91 if (fPmt)
92 {
93 delete fPmt;
94 fPmt = NULL;
95 fNPmt = 0;
96 }
97}
98
100{
101 LOG(info) << "R3BNeulandTacquilaMapped2Cal::Init : read " << fTcalPar->GetNumModulePar() << " calibrated modules";
102 // fTcalPar->printParams();
103
104 FairRootManager* mgr = FairRootManager::Instance();
105 if (NULL == mgr)
106 {
107 LOG(fatal) << "FairRootManager not found";
108 }
109
110 header = dynamic_cast<R3BEventHeader*>(mgr->GetObject("EventHeader."));
111 if (NULL == header)
112 {
113 LOG(fatal) << "Branch R3BEventHeader not found";
114 }
115
116 fRawHit = dynamic_cast<TClonesArray*>(mgr->GetObject("NeulandTacquilaMappedData"));
117 if (NULL == fRawHit)
118 {
119 LOG(fatal) << "Branch NeulandTacquilaMappedData not found";
120 }
121
122 mgr->Register("NeulandCalData", "Neuland", fPmt, kTRUE);
123
124 fh_pulser_5_2 = new TH1F("h_pulser_5_2", "Single PMT resolution Bar 5 vs 2", 40000, -200., 200.);
125 fh_pulser_105_2 = new TH1F("h_pulser_105_2", "Single PMT resolution Bar 105 vs 2", 40000, -200., 200.);
126
127 SetParameter();
128
129 return kSUCCESS;
130}
131
133{
134 FairRunAna* ana = FairRunAna::Instance();
135 FairRuntimeDb* rtdb = ana->GetRuntimeDb();
136 fTcalPar = dynamic_cast<R3BTCalPar*>(rtdb->getContainer("LandTCalPar"));
137 fQCalPar = dynamic_cast<R3BNeulandQCalPar*>(rtdb->getContainer("NeulandQCalPar"));
138}
139
140void R3BNeulandTacquilaMapped2Cal::SetParameter()
141{
142
143 std::map<Int_t, Double_t> tempMapQdcOffset;
144 Int_t i = 0;
145 for (Int_t plane = 1; i <= planes; plane++)
146 for (Int_t bar = 1; bar <= 50; bar++)
147 for (Int_t side = 1; side <= 2; side++)
148 {
149 tempMapQdcOffset[i] = fQCalPar->GetParAt(plane, bar, side);
150 i++;
151 }
152
153 LOG(info) << "R3BNeulandTacquilaMapped2Cal::SetParameter : Number of Parameters: " << i;
154
155 fMapQdcOffset = tempMapQdcOffset;
156}
157
159{
161 SetParameter();
162 return kSUCCESS;
163}
164
166{
167 if (fTrigger >= 0)
168 {
169 if (header->GetTrigger() != fTrigger)
170 {
171 return;
172 }
173 }
174
175 Int_t nHits = fRawHit->GetEntriesFast();
176 if (fPulserMode)
177 {
178 if (nHits < fNofPMTs)
179 {
180 return;
181 }
182 }
183 else
184 {
185 if (nHits > (fNofPMTs / 2))
186 {
187 return;
188 }
189 }
190
191 if (nHits > 0)
192 {
193 MakeCal();
194 }
195
196 if (fPulserMode)
197 {
198 R3BNeulandCalData* pmt1;
199 Double_t time1 = nan("");
200 for (Int_t i = 0; i < fNPmt; i++)
201 {
202 pmt1 = dynamic_cast<R3BNeulandCalData*>(fPmt->At(i));
203 if (pmt1->GetBarId() == 2 && pmt1->GetSide() == 1)
204 {
205 time1 = pmt1->GetTime();
206 break;
207 }
208 }
209 for (Int_t i = 0; i < fNPmt; i++)
210 {
211 pmt1 = dynamic_cast<R3BNeulandCalData*>(fPmt->At(i));
212 if (pmt1->GetBarId() == 5 && pmt1->GetSide() == 1)
213 {
214 fh_pulser_5_2->Fill(pmt1->GetTime() - time1);
215 }
216 if (pmt1->GetBarId() == 105 && pmt1->GetSide() == 1)
217 {
218 fh_pulser_105_2->Fill(pmt1->GetTime() - time1);
219 }
220 }
221 }
222}
223
224void R3BNeulandTacquilaMapped2Cal::MakeCal()
225{
226 Int_t nHits = fRawHit->GetEntriesFast();
229 Int_t iPlane;
230 Int_t iPaddle;
231 Int_t iSide;
232 Int_t channel = 0;
233 Int_t tdc;
234 R3BTCalModulePar* par;
235 Double_t time;
236 Double_t time2;
237 Int_t qdc;
238
239 for (Int_t khit = 0; khit < nHits; khit++)
240 {
241 hit2 = dynamic_cast<R3BNeulandTacquilaMappedData*>(fRawHit->At(khit));
242 if (NULL == hit2)
243 {
244 continue;
245 }
246
247 iPlane = hit2->GetPlane();
248 iPaddle = hit2->GetPaddle();
249 iSide = hit2->GetSide();
250
251 if (!(par = fTcalPar->GetModuleParAt(iPlane, iPaddle, iSide)))
252 {
253 LOG(debug) << "R3BNeulandTacquilaMapped2Cal::Exec : Tcal par not found, channel: " << iPlane << " / "
254 << iPaddle << " / " << iSide;
255 continue;
256 }
257
258 tdc = hit2->GetTacData();
259 time = par->GetTimeTacquila(tdc);
260 if (time < 0. || time > fClockFreq)
261 {
262 LOG(error) << "R3BNeulandTacquilaMapped2Cal::Exec : error in time calibration: ch=" << channel
263 << ", tdc=" << tdc << ", time=" << time;
264 continue;
265 }
266
267 if (!(par = fTcalPar->GetModuleParAt(iPlane, iPaddle, iSide + 2)))
268 {
269 LOG(debug) << "R3BNeulandTacquilaMapped2Cal::Exec : Tcal par not found, channel: " << iPlane << " / "
270 << iPaddle << " / " << (iSide + 2);
271 continue;
272 }
273
274 tdc = hit2->GetStopT();
275 time2 = par->GetTimeTacquila(tdc);
276 if (time2 < 0. || time2 > fClockFreq)
277 {
278 LOG(error) << "R3BNeulandTacquilaMapped2Cal::Exec : error in time calibration: ch=" << channel
279 << ", tdc=" << tdc << ", time=" << time2;
280 continue;
281 }
282
283 qdc = hit2->GetQdcData() - fMapQdcOffset[toID(iPlane, iPaddle, iSide)];
284 qdc = std::max(qdc, 0);
285
286 time = time - time2 + hit2->GetClock() * fClockFreq;
287 if (fWalkEnabled)
288 {
289 time += wlk(qdc);
290 }
291 new ((*fPmt)[fNPmt]) R3BNeulandCalData((iPlane - 1) * 50 + iPaddle, iSide, time, NAN, qdc);
292 fNPmt += 1;
293 }
294}
295
297{
298 if (fVerbose && 0 == (fNEvents % 1000))
299 {
300 LOG(info) << "R3BNeulandTacquilaMapped2Cal::Exec : event=" << fNEvents << " nPMTs=" << fNPmt;
301 }
302
303 if (fPmt)
304 {
305 fPmt->Clear();
306 fNPmt = 0;
307 }
308 fMap17Seen.clear();
309 fMapStopTime.clear();
310 fMapStopClock.clear();
311
312 fNEvents += 1;
313}
314
316{
317 fh_pulser_5_2->Write();
318 fh_pulser_105_2->Write();
319}
320
ClassImp(R3B::Neuland::Cal2HitPar)
#define toID(x, y, z)
Double_t wlk(Double_t x)
Int_t GetSide() const
Double_t GetTime() const
Int_t GetBarId() const
Int_t GetParAt(Int_t plane, Int_t bar, Int_t side)
Method to get the pedestal offset of a PMT.
An analysis task to apply TCAL calibration for NeuLAND.
virtual void Exec(Option_t *option)
Method for event loop implementation.
virtual InitStatus ReInit()
Method for re-initialization of parameter containers in case the Run ID has changed.
virtual void FinishEvent()
A method for finish of processing of an event.
virtual void FinishTask()
Method for finish of the task execution.
virtual void SetParContainers()
Method for initialization of the parameter containers.
virtual InitStatus Init()
Method for task initialization.