R3BROOT
R3B analysis software
Loading...
Searching...
No Matches
R3BNeulandTcal.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// ----- R3BNeuLandTcal -----
16// ----- Created 27-01-2015 by M.Heil -----
17// ------------------------------------------------------------
18
19#include "R3BNeulandTcal.h"
20
21#include "FairRootManager.h"
22#include "FairRunAna.h"
23#include "FairRuntimeDb.h"
24#include "FairTask.h"
25#include "R3BNeulandPmt.h"
27#include "R3BTCalEngine.h"
28#include "R3BTCalPar.h"
29
30#include "Rtypes.h"
31#include "RtypesCore.h"
32#include "TClonesArray.h"
33#include <cstddef>
34#include <fairlogger/Logger.h>
35
37 : FairTask("LandTcal", 1)
38 , fNEvents(0)
39 , fMappedHit(NULL)
40 , fPmt(new TClonesArray("R3BNeulandPmt"))
41 , fNPmt(0)
42 , fTcalPar(NULL)
43 , fTrigger(-1)
44 // , fMap17Seen()
45 // , fMapStopTime()
46 // , fMapStopClock()
47 , fClockFreq(1. / VFTX_CLOCK_MHZ * 1000.)
48{
49}
50
51R3BNeulandTcal::R3BNeulandTcal(const char* name, Int_t iVerbose)
52 : FairTask(name, iVerbose)
53 , fNEvents(0)
54 , fMappedHit(NULL)
55 , fPmt(new TClonesArray("R3BNeulandPmt"))
56 , fNPmt(0)
57 , fTcalPar(NULL)
58 , fTrigger(-1)
59 // , fMap17Seen()
60 // , fMapStopTime()
61 // , fMapStopClock()
62 , fClockFreq(1. / VFTX_CLOCK_MHZ * 1000.)
63{
64}
65
67{
68 if (fPmt)
69 {
70 delete fPmt;
71 fPmt = NULL;
72 fNPmt = 0;
73 }
74}
75
77{
78 LOG(info) << "R3BNeulandTcal::Init : read " << fTcalPar->GetNumModulePar() << " calibrated modules";
79 // fTcalPar->printParams();
80
81 FairRootManager* mgr = FairRootManager::Instance();
82 if (NULL == mgr)
83 {
84 LOG(fatal) << "FairRootManager not found";
85 }
86 /*
87 header = (R3BEventHeader*)mgr->GetObject("R3BEventHeader");
88 if (NULL == header)
89 {
90 FairLogger::GetLogger()->Fatal(MESSAGE_ORIGIN, "Branch R3BEventHeader not found");
91 }
92 */
93 fMappedHit = dynamic_cast<TClonesArray*>(mgr->GetObject("NeulandMappedData"));
94 if (NULL == fMappedHit)
95 {
96 LOG(fatal) << "Branch NeulandMappedData not found";
97 }
98
99 mgr->Register("NeulandPmt", "Land", fPmt, kTRUE);
100
101 return kSUCCESS;
102}
103
105{
106 FairRunAna* ana = FairRunAna::Instance();
107 FairRuntimeDb* rtdb = ana->GetRuntimeDb();
108 fTcalPar = dynamic_cast<R3BTCalPar*>(rtdb->getContainer("LandTCalPar"));
109}
110
112{
114 return kSUCCESS;
115}
116
117void R3BNeulandTcal::Exec(Option_t*)
118{
119 /*
120 if (fTrigger >= 0)
121 {
122 if (header->GetTrigger() != fTrigger)
123 {
124 return;
125 }
126 }
127*/
128 Int_t nHits = fMappedHit->GetEntriesFast();
129 /*
130 if (nHits > (fNofPMTs / 2))
131 {
132 return;
133 }
134 */
136 Int_t iPlane;
137 Int_t iBar;
138 const auto channel = 0;
139 Int_t tdc;
140 R3BTCalModulePar* par;
141 Double_t timeLE;
142 Double_t timeTE;
143
144 for (Int_t ihit = 0; ihit < nHits; ihit++)
145 {
146 hit = dynamic_cast<R3BPaddleTamexMappedData*>(fMappedHit->At(ihit));
147 if (NULL == hit)
148 {
149 continue;
150 }
151 iPlane = hit->GetPlaneId();
152 iBar = hit->GetBarId();
153 if (hit->Is17())
154 {
155 // 17-th channel
156 continue;
157 }
158
159 for (Int_t iSide = 0; iSide < 2; iSide++)
160 {
161
162 // !!! There is the Edge missing GetModulePar( ... ) !!!
163 // Convert TDC to [ns]
164 if (!(par = fTcalPar->GetModuleParAt(iPlane, iBar, iSide)))
165 {
166 LOG(debug) << "R3BNeulandTcal::Exec : Tcal par not found, barId: " << iBar << ", side: " << iSide;
167 continue;
168 }
169
170 tdc = hit->GetFineTime(iSide, 0); // PM, edge
171 timeLE = par->GetTimeVFTX(tdc);
172
173 if (!(par = fTcalPar->GetModuleParAt(iPlane, iBar, iSide + 2)))
174 {
175 LOG(debug) << "R3BNeulandTcal::Exec : Tcal par not found, barId: " << iBar << ", side: " << iSide;
176 continue;
177 }
178
179 tdc = hit->GetFineTime(iSide, 1);
180 timeTE = par->GetTimeVFTX(tdc);
181
182 /*
183 if (timeLE < -1000.)
184 {
185 continue;
186 }
187 */
188 if (timeLE < 0. || timeLE > fClockFreq || timeTE < 0. || timeTE > fClockFreq)
189 {
190 LOG(error) << "R3BNeulandTcal::Exec : error in time calibration: ch= " << channel << ", tdc= " << tdc
191 << ", time leading edge = " << timeLE << ", time trailing edge = " << timeTE;
192 continue;
193 }
194
195 timeLE = fClockFreq - timeLE + hit->GetCoarseTime(iSide, 0) * fClockFreq;
196 timeTE = fClockFreq - timeTE + hit->GetCoarseTime(iSide, 1) * fClockFreq;
197
198 new ((*fPmt)[fNPmt]) R3BNeulandPmt(iPlane, iBar, iSide, timeLE, timeTE - timeLE);
199 fNPmt += 1;
200 } // for Side
201 }
202}
203
205{
206 if (fVerbose && 0 == (fNEvents % 1))
207 {
208 LOG(info) << "R3BNeulandTcal::Exec : event=" << fNEvents << " nPMTs=" << fNPmt;
209 }
210
211 if (fPmt)
212 {
213 fPmt->Clear();
214 fNPmt = 0;
215 }
216
217 fNEvents += 1;
218}
219
221
ClassImp(R3B::Neuland::Cal2HitPar)
An analysis task to apply TCAL calibration for NeuLAND.
Int_t fNPmt
Number of produced time items per event.
virtual InitStatus Init()
Method for task initialization.
Int_t fTrigger
Trigger value.
TClonesArray * fMappedHit
Array with raw items - input data.
Int_t fNEvents
Event counter.
TClonesArray * fPmt
Array with time items - output data.
R3BTCalPar * fTcalPar
TCAL parameter container.
virtual void FinishEvent()
A method for finish of processing of an event.
R3BNeulandTcal()
Default constructor.
Double_t fClockFreq
Clock cycle in [ns].
virtual ~R3BNeulandTcal()
Destructor.
virtual void SetParContainers()
Method for initialization of the parameter containers.
virtual void FinishTask()
Method for finish of the task execution.
virtual InitStatus ReInit()
Method for re-initialization of parameter containers in case the Run ID has changed.
virtual void Exec(Option_t *option)
Method for event loop implementation.
const Int_t & GetFineTime(int t, int e) const
const Int_t & GetCoarseTime(int t, int e) const