R3BROOT
R3B analysis software
Loading...
Searching...
No Matches
R3BNeulandMapped2CalPar.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 "FairLogger.h"
16#include "FairRootManager.h"
17#include "FairRtdbRun.h"
18#include "FairRunIdGenerator.h"
19#include "FairRuntimeDb.h"
20#include "R3BEventHeader.h"
22#include "R3BTCalEngine.h"
23#include "R3BTCalPar.h"
24#include "TClonesArray.h"
25#include "TF1.h"
26#include <iostream>
27#include <signal.h>
28#include <stdlib.h>
29
30using namespace std;
31
33 : FairTask("R3BNeulandMapped2TCalPar", 1)
34 , fMinStats(100000)
35 , fTrigger(-1)
36 // , fNofPMTs(0)
37 , fNEvents(0)
38 , fCal_Par(NULL)
39{
40}
41
42R3BNeulandMapped2CalPar::R3BNeulandMapped2CalPar(const char* name, Int_t iVerbose)
43 : FairTask(name, iVerbose)
44 , fMinStats(100000)
45 , fTrigger(-1)
46 //, fNofPMTs(0)
47 , fNEvents(0)
48 , fCal_Par(NULL)
49{
50}
51
53{
54 if (fEngine)
55 {
56 delete fEngine;
57 }
58}
59
61{
62 FairRootManager* rm = FairRootManager::Instance();
63 if (!rm)
64 {
65 return kFATAL;
66 }
67 header = dynamic_cast<R3BEventHeader*>(rm->GetObject("EventHeader."));
68 if (!header)
69 {
70 return kFATAL;
71 }
72 fHits = dynamic_cast<TClonesArray*>(rm->GetObject("NeulandMappedData"));
73 if (!fHits)
74 {
75 return kFATAL;
76 }
77 fHitsTrigger = dynamic_cast<TClonesArray*>(rm->GetObject("NeulandTrigMappedData"));
78 if (!fHitsTrigger)
79 {
80 LOG(info) << "Branch NeulandTrigMapped not found";
81 }
82
83 // container needs to be created in tcal/R3BTCalContFact.cxx AND R3BTCal needs
84 // to be set as dependency in CMakelists.txt (in this case in the land directory)
85 fCal_Par = dynamic_cast<R3BTCalPar*>(FairRuntimeDb::instance()->getContainer("LandTCalPar"));
86 fCal_Par->setChanged();
87
88 fEngine = new R3BTCalEngine(fCal_Par, fMinStats);
89
90 for (Int_t pln = 0; pln < fNofPlanes; pln++)
91 {
92 for (Int_t bar = 0; bar < fNofBarsPerPlane; bar++)
93 {
94 for (Int_t pmt = 0; pmt < 4; pmt++)
95 {
96 counts[pln][bar][pmt] = 0;
97 }
98 }
99 }
100
101 checkcounts = 0;
102
103 return kSUCCESS;
104}
105
106void R3BNeulandMapped2CalPar::Exec(Option_t* option)
107{
108
109 if (checkcounts == fNofPMTs)
110 {
111 // std::cout << "done " << std::endl;
112 // raise(SIGINT);
113 }
114
115 if (fTrigger >= 0)
116 {
117 if (header->GetTrigger() != fTrigger)
118 {
119 return;
120 }
121 }
122
123 Int_t nHits = fHits->GetEntries();
124 if (nHits > (fNofPMTs / 2))
125 {
126 // return;
127 }
128
129 // Loop over mapped hits
130 for (Int_t i = 0; i < nHits; i++)
131 {
132 auto hit = dynamic_cast<R3BPaddleTamexMappedData*>(fHits->At(i));
133 if (!hit)
134 {
135 continue;
136 }
137
138 // Check bar ID
139 Int_t iPlane = hit->GetPlaneId();
140 Int_t iBar = hit->GetBarId();
141 Int_t iSide = -1 == hit->fCoarseTime1LE ? 2 : 1;
142
143 // fill 1 1LE
144 // 2 1TE
145 // 3 2LE
146 // 4 2TE
147
148 Int_t iFine;
149 if (1 == iSide)
150 {
151 iFine = hit->fFineTime1LE;
152 }
153 else
154 {
155 iFine = hit->fFineTime2LE;
156 }
157
158 fEngine->Fill(iPlane, iBar, (iSide - 1) * 2 + 1, iFine);
159 counts[iPlane - 1][iBar - 1][(iSide - 1) * 2]++;
160 if (counts[iPlane - 1][iBar - 1][(iSide - 1) * 2] == fMinStats)
161 {
162 checkcounts++;
163 std::cout << iPlane << "a " << iBar << " " << iSide << std::endl;
164 std::cout << checkcounts << std::endl;
165 }
166
167 if (1 == iSide)
168 {
169 iFine = hit->fFineTime1TE;
170 }
171 else
172 {
173 iFine = hit->fFineTime2TE;
174 }
175 fEngine->Fill(iPlane, iBar, (iSide - 1) * 2 + 2, iFine);
176 counts[iPlane - 1][iBar - 1][(iSide - 1) * 2 + 1]++;
177 if (counts[iPlane - 1][iBar - 1][(iSide - 1) * 2 + 1] == fMinStats)
178 {
179 checkcounts++;
180 std::cout << iPlane << "b " << iBar << " " << iSide << std::endl;
181 std::cout << checkcounts << std::endl;
182 }
183 }
184
185 // Loop over mapped triggers
186 if (fHitsTrigger)
187 {
188 nHits = fHitsTrigger->GetEntriesFast();
189 for (Int_t i = 0; i < nHits; i++)
190 {
191 auto hit = dynamic_cast<R3BPaddleTamexMappedData*>(fHitsTrigger->At(i));
192 if (!hit)
193 {
194 continue;
195 }
196
197 // Check bar ID
198 auto iBar = hit->GetBarId();
199 auto iFine = hit->fFineTime1LE;
200
201 fEngine->Fill(100, iBar, 10, iFine);
202 }
203 }
204
205 // Increment events
206 fNEvents += 1;
207}
208
210
212{
213 fEngine->CalculateParamVFTX();
214 fCal_Par->printParams();
215}
216
ClassImp(R3B::Neuland::Cal2HitPar)
An analysis task for TCAL calibration of NeuLAND data.
virtual ~R3BNeulandMapped2CalPar()
Destructor.
virtual void FinishTask()
Method for finish of the task execution.
virtual void Exec(Option_t *option)
Method for event loop implementation.
R3BNeulandMapped2CalPar()
Default constructor.
virtual void FinishEvent()
A method for finish of processing of an event.
virtual InitStatus Init()
Method for task initialization.