Belle II Software development
SVDPerformanceTTreeModule.cc
1/**************************************************************************
2 * basf2 (Belle II Analysis Software Framework) *
3 * Author: The Belle II Collaboration *
4 * *
5 * See git log for contributors and copyright holders. *
6 * This file is licensed under LGPL-3.0, see LICENSE.md. *
7 **************************************************************************/
8
9#include <svd/modules/svdPerformance/SVDPerformanceTTreeModule.h>
10#include <framework/datastore/StoreArray.h>
11#include <framework/gearbox/Unit.h>
12#include <framework/core/Environment.h>
13#include <svd/dataobjects/SVDCluster.h>
14#include <svd/dataobjects/SVDRecoDigit.h>
15#include <svd/dataobjects/SVDShaperDigit.h>
16#include <svd/dataobjects/SVDTrueHit.h>
17#include <svd/dataobjects/SVDEventInfo.h>
18#include <vxd/dataobjects/VxdID.h>
19#include <vxd/geometry/GeoCache.h>
20#include <vxd/geometry/SensorInfoBase.h>
21#include <tracking/dataobjects/RecoTrack.h>
22#include <tracking/dataobjects/RecoHitInformation.h>
23#include <genfit/TrackPoint.h>
24#include <Math/Boost.h>
25#include <math.h>
26#include <iostream>
27#include <algorithm>
28#include <mdst/dataobjects/Track.h>
29#include <mdst/dataobjects/TrackFitResult.h>
30
31
32using namespace Belle2;
33using namespace std;
34
35
36//-----------------------------------------------------------------
37// Register the Module
38//-----------------------------------------------------------------
39REG_MODULE(SVDPerformanceTTree);
40
41//-----------------------------------------------------------------
42// Implementation
43//-----------------------------------------------------------------
44
46{
47 //Set module properties
48 setDescription("The module is used to create a TTree to study SVD clusters, genfit unbiased residuals and many other properties related to the track they belong to.");
49 //Parameter to take only specific RecoTracks as input
50 addParam("outputFileName", m_rootFileName, "Name of output root file.", std::string("SVDPerformanceTTree.root"));
51 addParam("recoTracksStoreArrayName", m_recoTracksStoreArrayName, "StoreArray name of the input and output RecoTracks.",
53}
54
56{
58 recoTracks.isOptional();
59 m_EventT0.isOptional();
60
61 m_rootFilePtr = new TFile(m_rootFileName.c_str(), "RECREATE");
62
63 //Tree for SVD u overlapping clusters
64 m_t_U = new TTree("t_U", "Tree for SVD u-clusters");
65 m_t_U->Branch("cdcEventT0", &m_cdcEventT0, "cdcEventT0/F");
66 m_t_U->Branch("cdcEventT0_6SRF", &m_cdcEventT0_6SRF, "cdcEventT0_6SRF/F");
67 m_t_U->Branch("cdcEventT0_3SRF", &m_cdcEventT0_3SRF, "cdcEventT0_3SRF/F");
68 m_t_U->Branch("cdcEventT0Err", &m_cdcEventT0Err, "cdcEventT0Err/F");
69 m_t_U->Branch("svdTB", &m_svdTB, "svdTB/i");
70 m_t_U->Branch("svdClSNR", &m_svdClSNR, "svdClSNR/F");
71 m_t_U->Branch("svdClCharge", &m_svdClCharge, "svdClCharge/F");
72 m_t_U->Branch("svdStripCharge", &m_svdStripCharge);
73 m_t_U->Branch("svdStrip6Samples", &m_svdStrip6Samples);
74 m_t_U->Branch("svdClTime", &m_svdClTime, "svdClTime/F");
75 m_t_U->Branch("svdClTimeErr", &m_svdClTimeErr, "svdClTimeErr/F");
76 m_t_U->Branch("svdClTime_6SRF", &m_svdClTime_6SRF, "svdClTime_6SRF/F");
77 m_t_U->Branch("svdClTime_3SRF", &m_svdClTime_3SRF, "svdClTime_3SRF/F");
78 m_t_U->Branch("svdFF", &m_svdFF, "svdFF/i");
79 m_t_U->Branch("svdStripTime", &m_svdStripTime);
80 m_t_U->Branch("svdStripPosition", &m_svdStripPosition);
81 m_t_U->Branch("svdRes", &m_svdRes, "svdRes/F");
82 m_t_U->Branch("svdPitch", &m_svdPitch, "svdPitch/F");
83 m_t_U->Branch("svdWidth", &m_svdWidth, "svdWidth/F");
84 m_t_U->Branch("svdLength", &m_svdLength, "svdLength/F");
85 m_t_U->Branch("svdClIntStrPos", &m_svdClIntStrPos, "svdClIntStrPos/F");
86 m_t_U->Branch("svdClPos", &m_svdClPos, "svdClPos/F");
87 m_t_U->Branch("svdClPosErr", &m_svdClPosErr, "svdClPosErr/F");
88 m_t_U->Branch("svdTruePos", &m_svdTruePos, "svdTruePos/F");
89 m_t_U->Branch("svdClPhi", &m_svdClPhi, "svdClPhi/F");
90 m_t_U->Branch("svdClZ", &m_svdClZ, "svdClZ/F");
91 m_t_U->Branch("svdTrkd0", &m_svdTrkd0, "svdTrkd0/F");
92 m_t_U->Branch("svdTrkz0", &m_svdTrkz0, "svdTrkz0/F");
93 m_t_U->Branch("svdTrkpT", &m_svdTrkpT, "svdTrkpT/F");
94 m_t_U->Branch("svdTrkpCM", &m_svdTrkpCM, "svdTrkpCM/F");
95 m_t_U->Branch("svdTrkTraversedLength", &m_svdTrkTraversedLength, "svdTrkTraversedLength/F");
96 m_t_U->Branch("svdTrkPXDHits", &m_svdTrkPXDHits, "svdTrkPXDHits/i");
97 m_t_U->Branch("svdTrkSVDHits", &m_svdTrkSVDHits, "svdTrkSVDHits/i");
98 m_t_U->Branch("svdTrkCDCHits", &m_svdTrkCDCHits, "svdTrkCDCHits/i");
99 m_t_U->Branch("svdTrkPos", &m_svdTrkPos, "svdTrkPos/F");
100 m_t_U->Branch("svdTrkPosOS", &m_svdTrkPosOS, "svdTrkPosOS/F");
101 m_t_U->Branch("svdTrkPosErr", &m_svdTrkPosErr, "svdTrkPosErr/F");
102 m_t_U->Branch("svdTrkPosErrOS", &m_svdTrkPosErrOS, "svdTrkPosErrOS/F");
103 m_t_U->Branch("svdTrkQoP", &m_svdTrkQoP, "svdTrkQoP/F");
104 m_t_U->Branch("svdTrkPrime", &m_svdTrkPrime, "svdTrkPrime/F");
105 m_t_U->Branch("svdTrkPrimeOS", &m_svdTrkPrimeOS, "svdTrkPrimeOS/F");
106 m_t_U->Branch("svdTrkPosUnbiased", &m_svdTrkPosUnbiased, "svdTrkPosUnbiased/F");
107 m_t_U->Branch("svdTrkPosErrUnbiased", &m_svdTrkPosErrUnbiased, "svdTrkPosErrUnbiased/F");
108 m_t_U->Branch("svdTrkQoPUnbiased", &m_svdTrkQoPUnbiased, "svdTrkQoPUnbiased/F");
109 m_t_U->Branch("svdTrkPrimeUnbiased", &m_svdTrkPrimeUnbiased, "svdTrkPrimeUnbiased/F");
110 m_t_U->Branch("svdLayer", &m_svdLayer, "svdLayer/i");
111 m_t_U->Branch("svdLadder", &m_svdLadder, "svdLadder/i");
112 m_t_U->Branch("svdSensor", &m_svdSensor, "svdSensor/i");
113 m_t_U->Branch("svdSize", &m_svdSize, "svdSize/i");
114 //Tree for SVD v overlapping clusters
115 m_t_V = new TTree("t_V", "Tree for SVD v-clusters");
116 m_t_V->Branch("cdcEventT0", &m_cdcEventT0, "cdcEventT0/F");
117 m_t_V->Branch("cdcEventT0_6SRF", &m_cdcEventT0_6SRF, "cdcEventT0_6SRF/F");
118 m_t_V->Branch("cdcEventT0_3SRF", &m_cdcEventT0_3SRF, "cdcEventT0_3SRF/F");
119 m_t_V->Branch("cdcEventT0Err", &m_cdcEventT0Err, "cdcEventT0Err/F");
120 m_t_V->Branch("svdTB", &m_svdTB, "svdTB/i");
121 m_t_V->Branch("svdClSNR", &m_svdClSNR, "svdClSNR/F");
122 m_t_V->Branch("svdClCharge", &m_svdClCharge, "svdClCharge/F");
123 m_t_V->Branch("svdStripCharge", &m_svdStripCharge);
124 m_t_V->Branch("svdStrip6Samples", &m_svdStrip6Samples);
125 m_t_V->Branch("svdClTime", &m_svdClTime, "svdClTime/F");
126 m_t_V->Branch("svdClTimeErr", &m_svdClTimeErr, "svdClTimeErr/F");
127 m_t_V->Branch("svdClTime_6SRF", &m_svdClTime_6SRF, "svdClTime_6SRF/F");
128 m_t_V->Branch("svdClTime_3SRF", &m_svdClTime_3SRF, "svdClTime_3SRF/F");
129 m_t_V->Branch("svdFF", &m_svdFF, "svdFF/i");
130 m_t_V->Branch("svdStripTime", &m_svdStripTime);
131 m_t_V->Branch("svdStripPosition", &m_svdStripPosition);
132 m_t_V->Branch("svdRes", &m_svdRes, "svdRes/F");
133 m_t_V->Branch("svdPitch", &m_svdPitch, "svdPitch/F");
134 m_t_V->Branch("svdWidth", &m_svdWidth, "svdWidth/F");
135 m_t_V->Branch("svdLength", &m_svdLength, "svdLength/F");
136 m_t_V->Branch("svdClIntStrPos", &m_svdClIntStrPos, "svdClIntStrPos/F");
137 m_t_V->Branch("svdClPos", &m_svdClPos, "svdClPos/F");
138 m_t_V->Branch("svdClPosErr", &m_svdClPosErr, "svdClPosErr/F");
139 m_t_V->Branch("svdTruePos", &m_svdTruePos, "svdTruePos/F");
140 m_t_V->Branch("svdClPhi", &m_svdClPhi, "svdClPhi/F");
141 m_t_V->Branch("svdClZ", &m_svdClZ, "svdClZ/F");
142 m_t_V->Branch("svdTrkd0", &m_svdTrkd0, "svdTrkd0/F");
143 m_t_V->Branch("svdTrkz0", &m_svdTrkz0, "svdTrkz0/F");
144 m_t_V->Branch("svdTrkpT", &m_svdTrkpT, "svdTrkpT/F");
145 m_t_V->Branch("svdTrkpCM", &m_svdTrkpCM, "svdTrkpCM/F");
146 m_t_V->Branch("svdTrkTraversedLength", &m_svdTrkTraversedLength, "svdTrkTraversedLength/F");
147 m_t_V->Branch("svdTrkPXDHits", &m_svdTrkPXDHits, "svdTrkPXDHits/i");
148 m_t_V->Branch("svdTrkSVDHits", &m_svdTrkSVDHits, "svdTrkSVDHits/i");
149 m_t_V->Branch("svdTrkCDCHits", &m_svdTrkCDCHits, "svdTrkCDCHits/i");
150 m_t_V->Branch("svdTrkPos", &m_svdTrkPos, "svdTrkPos/F");
151 m_t_V->Branch("svdTrkPosOS", &m_svdTrkPosOS, "svdTrkPosOS/F");
152 m_t_V->Branch("svdTrkPosErr", &m_svdTrkPosErr, "svdTrkPosErr/F");
153 m_t_V->Branch("svdTrkPosErrOS", &m_svdTrkPosErrOS, "svdTrkPosErrOS/F");
154 m_t_V->Branch("svdTrkQoP", &m_svdTrkQoP, "svdTrkQoP/F");
155 m_t_V->Branch("svdTrkPrime", &m_svdTrkPrime, "svdTrkPrime/F");
156 m_t_V->Branch("svdTrkPrimeOS", &m_svdTrkPrimeOS, "svdTrkPrimeOS/F");
157 m_t_V->Branch("svdTrkPosUnbiased", &m_svdTrkPosUnbiased, "svdTrkPosUnbiased/F");
158 m_t_V->Branch("svdTrkPosErrUnbiased", &m_svdTrkPosErrUnbiased, "svdTrkPosErrUnbiased/F");
159 m_t_V->Branch("svdTrkQoPUnbiased", &m_svdTrkQoPUnbiased, "svdTrkQoPUnbiased/F");
160 m_t_V->Branch("svdTrkPrimeUnbiased", &m_svdTrkPrimeUnbiased, "svdTrkPrimeUnbiased/F");
161 m_t_V->Branch("svdLayer", &m_svdLayer, "svdLayer/i");
162 m_t_V->Branch("svdLadder", &m_svdLadder, "svdLadder/i");
163 m_t_V->Branch("svdSensor", &m_svdSensor, "svdSensor/i");
164 m_t_V->Branch("svdSize", &m_svdSize, "svdSize/i");
165
166}
167
169{
171 m_apvClockPeriod = 1. / m_hwClock->getClockFrequency(Const::EDetector::SVD, "sampling");
172}
174{
175
176 m_cdcEventT0 = std::numeric_limits<float>::quiet_NaN();
177 m_cdcEventT0_6SRF = std::numeric_limits<float>::quiet_NaN();
178 m_cdcEventT0_3SRF = std::numeric_limits<float>::quiet_NaN();
179 m_cdcEventT0Err = std::numeric_limits<float>::quiet_NaN();
180
181 if (m_EventT0.isValid())
182 if (m_EventT0->hasEventT0()) {
183 if (m_EventT0->hasTemporaryEventT0(Const::EDetector::CDC)) {
184 const auto bestCDCEvtT0 = m_EventT0->getBestCDCTemporaryEventT0();
185 m_cdcEventT0 = bestCDCEvtT0->eventT0 ;
186 m_cdcEventT0Err = bestCDCEvtT0->eventT0Uncertainty;
187 }
188 }
189
190 //first check SVDEventInfo name
191 StoreObjPtr<SVDEventInfo> temp_eventinfo("SVDEventInfo");
192 std::string m_svdEventInfoName = "SVDEventInfo";
193 if (!temp_eventinfo.isValid())
194 m_svdEventInfoName = "SVDEventInfoSim";
195 StoreObjPtr<SVDEventInfo> eventinfo(m_svdEventInfoName);
196 if (!eventinfo) B2ERROR("No SVDEventInfo!");
197 m_svdTB = eventinfo->getModeByte().getTriggerBin();
198
199 bool isMC = Environment::Instance().isMC();
200
203 for (const auto& trk : recoTracks) {
204 if (! trk.wasFitSuccessful()) {
205 continue;
206 }
207
208
209 RelationVector<Track> theTK = DataStore::getRelationsWithObj<Track>(&trk);
210
211 if (theTK.size() == 0) {
212 continue;
213 }
214
215
216 const TrackFitResult* tfr = theTK[0]->getTrackFitResultWithClosestMass(Const::pion);
217
218 if (tfr) {
219 m_svdTrkd0 = tfr->getD0();
220 m_svdTrkz0 = tfr->getZ0();
221 m_svdTrkpT = tfr->getMomentum().Rho();
222 ROOT::Math::PxPyPzEVector pStar = tfr->get4Momentum();
223 ROOT::Math::BoostZ boost(3. / 11);
224 pStar = boost(pStar);
225 m_svdTrkpCM = pStar.P();
226 }
227
228
229 const vector<SVDCluster* > svdClusters = trk.getSVDHitList();
230 B2DEBUG(40, "FITTED TRACK: NUMBER OF SVD HITS = " << svdClusters.size());
231
232 m_svdTrkPXDHits = (trk.getPXDHitList()).size();
233 m_svdTrkSVDHits = (trk.getSVDHitList()).size();
234 m_svdTrkCDCHits = (trk.getCDCHitList()).size();
235
236 for (unsigned int i = 0; i < svdClusters.size(); i++) {
237
238 const SVDCluster* svd_1 = svdClusters[i];
239
240 //get true hits, used only if isMC
241 RelationVector<SVDTrueHit> trueHit_1 = DataStore::getRelationsWithObj<SVDTrueHit>(svd_1);
242
243 const RecoHitInformation* infoSVD_1 = trk.getRecoHitInformation(svd_1);
244 if (!infoSVD_1) {
245 continue;
246 }
247 const auto* hitTrackPoint_1 = trk.getCreatedTrackPoint(infoSVD_1);
248 const auto* fittedResult_1 = hitTrackPoint_1->getFitterInfo();
249 if (!fittedResult_1) {
250 continue;
251 }
252 const VxdID svd_id_1 = svd_1->getSensorID();
253 const unsigned short svd_Layer_1 = svd_id_1.getLayerNumber();
254 const unsigned short svd_Ladder_1 = svd_id_1.getLadderNumber();
255 const unsigned short svd_Sensor_1 = svd_id_1.getSensorNumber();
256
257 try {
258 const TVectorD resUnBias_1 = fittedResult_1->getResidual(0, false).getState();
259 genfit::MeasuredStateOnPlane state_unbiased = fittedResult_1->getFittedState(false);
260 const TVectorD& svd_predIntersect_unbiased = state_unbiased.getState();
261 const TMatrixDSym& covMatrix_unbiased = state_unbiased.getCov();
262 genfit::MeasuredStateOnPlane state_1 = trk.getMeasuredStateOnPlaneFromRecoHit(infoSVD_1);
263 const TVectorD& svd_predIntersect_1 = state_1.getState();
264 const TMatrixDSym& covMatrix_1 = state_1.getCov();
265
266 if (svd_1->isUCluster()) {
267
268 const int strips_1 = svd_1->getSize();
269
270 const double res_U_1 = resUnBias_1.GetMatrixArray()[0] * Unit::convertValueToUnit(1.0, "um");
271 const ROOT::Math::XYZVector svdLocal_1(svd_1->getPosition(), svd_predIntersect_1[4], 0.);
272 const VXD::SensorInfoBase& svdSensor_1 = geo.getSensorInfo(svd_id_1);
273 const ROOT::Math::XYZVector& svdGlobal_1 = svdSensor_1.pointToGlobal(svdLocal_1);
274 double svdPhi_1 = svdGlobal_1.Phi();
275 double svdZ_1 = svdGlobal_1.Z();
276
277 m_svdFF = svd_1->getFirstFrame();
278 //Fill SVD tree for u-overlaps if required by the user
279 m_svdRes = res_U_1;
280 m_svdClTime = svd_1->getClsTime();
282 m_svdClSNR = svd_1->getSNR();
283 m_svdClCharge = svd_1->getCharge();
285 if (isMC && trueHit_1.size() > 0)
286 m_svdTruePos = trueHit_1[0]->getU();
287 else
288 m_svdTruePos = -99;
289 m_svdClPhi = svdPhi_1;
290 m_svdClZ = svdZ_1;
291 m_svdTrkPos = svd_predIntersect_1[3];
292 m_svdTrkPosOS = svd_predIntersect_1[4];
293 m_svdTrkPosErr = sqrt(covMatrix_1[3][3]);
294 m_svdTrkPosErrOS = sqrt(covMatrix_1[4][4]);
295 m_svdTrkQoP = svd_predIntersect_1[0];
296 m_svdTrkPrime = svd_predIntersect_1[1];
297 m_svdTrkPrimeOS = svd_predIntersect_1[2];
299 m_svdTrkPosUnbiased = svd_predIntersect_unbiased[3];
301 m_svdTrkPosErrUnbiased = sqrt(covMatrix_unbiased[3][3]);
302 m_svdTrkQoPUnbiased = svd_predIntersect_unbiased[0];
303 m_svdTrkPrimeUnbiased = svd_predIntersect_unbiased[1];
304 m_svdLayer = svd_Layer_1;
305 m_svdLadder = svd_Ladder_1;
306 m_svdSensor = svd_Sensor_1;
307 m_svdSize = strips_1;
308
309 m_svdPitch = svdSensor_1.getUPitch(m_svdTrkPosOS);
310 m_svdWidth = svdSensor_1.getUSize(m_svdTrkPosOS);
311 m_svdLength = svdSensor_1.getVSize();
312
314
315 m_svdStripCharge.clear();
316 m_svdStripTime.clear();
317 m_svdStripPosition.clear();
318 m_svdStrip6Samples.clear();
319 //retrieve relations and set strip charges and times
320 RelationVector<SVDRecoDigit> theRecoDigits = DataStore::getRelationsWithObj<SVDRecoDigit>(svd_1);
321 if ((theRecoDigits.size() != m_svdSize) && (m_svdSize != 128)) //virtual cluster
322 B2ERROR(" Inconsistency with cluster size! # recoDigits = " << theRecoDigits.size() << " != " << m_svdSize << " cluster size");
323
324 //skip clusters created beacuse of missing APV
325 if (m_svdSize < 128)
326 for (unsigned int d = 0; d < m_svdSize; d++) {
327
328 SVDShaperDigit* ShaperDigit = theRecoDigits[d]->getRelated<SVDShaperDigit>();
329 array<float, 6> Samples = ShaperDigit->getSamples();
330
331 m_svdStripCharge.push_back(theRecoDigits[d]->getCharge());
332 std::copy(std::begin(Samples), std::end(Samples), std::back_inserter(m_svdStrip6Samples));
333 m_svdStripTime.push_back(theRecoDigits[d]->getTime());
334 double misalignedStripPos = svdSensor_1.getUCellPosition(theRecoDigits[d]->getCellID());
335 //aligned strip pos = misaligned strip - ( misaligned cluster - aligned cluster)
336 m_svdStripPosition.push_back(misalignedStripPos - svd_1->getPosition() + m_svdClPos);
337 }
338
339 m_cdcEventT0_3SRF = eventinfo->getTimeInSVDReference(m_cdcEventT0, m_svdFF);
341 m_svdClTime_3SRF = eventinfo->getTimeInSVDReference(m_svdClTime, m_svdFF);
343
344 m_t_U->Fill();
345
346 } else {
347 const int strips_1 = svd_1->getSize();
348 const double res_V_1 = resUnBias_1.GetMatrixArray()[0] * Unit::convertValueToUnit(1.0, "um");
349 const ROOT::Math::XYZVector svdLocal_1(svd_predIntersect_1[3], svd_1->getPosition(), 0.);
350 const VXD::SensorInfoBase& svdSensor_1 = geo.getSensorInfo(svd_id_1);
351 const ROOT::Math::XYZVector& svdGlobal_1 = svdSensor_1.pointToGlobal(svdLocal_1);
352 double svdPhi_1 = svdGlobal_1.Phi();
353 double svdZ_1 = svdGlobal_1.Z();
354
355 m_svdFF = svd_1->getFirstFrame();
356
357 m_svdRes = res_V_1;
358 m_svdClTime = svd_1->getClsTime();
360 m_svdClSNR = svd_1->getSNR();
361 m_svdClCharge = svd_1->getCharge();
363 if (isMC && trueHit_1.size() > 0)
364 m_svdTruePos = trueHit_1[0]->getV();
365 else
366 m_svdTruePos = -99;
367 m_svdClPhi = svdPhi_1;
368 m_svdClZ = svdZ_1;
369 m_svdTrkPos = svd_predIntersect_1[4];
370 m_svdTrkPosOS = svd_predIntersect_1[3];
371 m_svdTrkPosErr = sqrt(covMatrix_1[4][4]);
372 m_svdTrkPosErrOS = sqrt(covMatrix_1[3][3]);
373 m_svdTrkQoP = svd_predIntersect_1[0];
374 m_svdTrkPrime = svd_predIntersect_1[2];
375 m_svdTrkPrimeOS = svd_predIntersect_1[1];
377 m_svdTrkPosUnbiased = svd_predIntersect_unbiased[4];
379 m_svdTrkPosErrUnbiased = sqrt(covMatrix_unbiased[4][4]);
380 m_svdTrkQoPUnbiased = svd_predIntersect_unbiased[0];
381 m_svdTrkPrimeUnbiased = svd_predIntersect_unbiased[2];
382 m_svdLayer = svd_Layer_1;
383 m_svdLadder = svd_Ladder_1;
384 m_svdSensor = svd_Sensor_1;
385 m_svdSize = strips_1;
386
387 m_svdPitch = svdSensor_1.getVPitch();
388 m_svdWidth = svdSensor_1.getUSize(m_svdTrkPos);
389 m_svdLength = svdSensor_1.getVSize();
390
392
393 m_svdStripCharge.clear();
394 m_svdStripTime.clear();
395 m_svdStripPosition.clear();
396 m_svdStrip6Samples.clear();
397 //retrieve relations and set strip charges and times
398 RelationVector<SVDRecoDigit> theRecoDigits = DataStore::getRelationsWithObj<SVDRecoDigit>(svd_1);
399 if ((theRecoDigits.size() != m_svdSize) && (m_svdSize != 128)) //virtual cluster
400 B2ERROR(" Inconsistency with cluster size! # recoDigits = " << theRecoDigits.size() << " != " << m_svdSize << " cluster size");
401
402 //skip clusters created beacuse of missing APV
403 if (m_svdSize < 128)
404 for (unsigned int d = 0; d < m_svdSize; d++) {
405 SVDShaperDigit* ShaperDigit = theRecoDigits[d]->getRelated<SVDShaperDigit>();
406 array<float, 6> Samples = ShaperDigit->getSamples();
407 m_svdStripCharge.push_back(theRecoDigits[d]->getCharge());
408 std::copy(std::begin(Samples), std::end(Samples), std::back_inserter(m_svdStrip6Samples));
409 m_svdStripTime.push_back(theRecoDigits[d]->getTime());
410 double misalignedStripPos = svdSensor_1.getVCellPosition(theRecoDigits[d]->getCellID());
411 //Aligned strip pos = misaligned strip - ( misaligned cluster - aligned cluster)
412 m_svdStripPosition.push_back(misalignedStripPos - svd_1->getPosition() + m_svdClPos);
413 }
414
415 m_cdcEventT0_3SRF = eventinfo->getTimeInSVDReference(m_cdcEventT0, m_svdFF);
417 m_svdClTime_3SRF = eventinfo->getTimeInSVDReference(m_svdClTime, m_svdFF);
419
420
421 m_t_V->Fill();
422 }
423 } catch (...) {
424 B2INFO("oops...something went wrong in getting the unbiased state, skipping this cluster.");
425 continue;
426 }
427
428 }
429 }
430}
431
432
434{
435
436 if (m_rootFilePtr != nullptr) {
437 m_rootFilePtr->cd();
438 m_t_U->Write();
439 m_t_V->Write();
440 m_rootFilePtr->Close();
441 }
442}
static const ChargedStable pion
charged pion particle
Definition: Const.h:661
bool isMC() const
Do we have generated, not real data?
Definition: Environment.cc:54
static Environment & Instance()
Static method to get a reference to the Environment instance.
Definition: Environment.cc:28
Base class for Modules.
Definition: Module.h:72
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
This class stores additional information to every CDC/SVD/PXD hit stored in a RecoTrack.
Class for type safe access to objects that are referred to in relations.
size_t size() const
Get number of relations.
The SVD Cluster class This class stores all information about reconstructed SVD clusters.
Definition: SVDCluster.h:29
float getClsTime() const
Get average of waveform maximum times of cluster strip signals.
Definition: SVDCluster.h:134
float getSNR() const
Get cluster SNR.
Definition: SVDCluster.h:159
unsigned short getSize() const
Get cluster size.
Definition: SVDCluster.h:154
float getCharge() const
Get collected charge.
Definition: SVDCluster.h:144
VxdID getSensorID() const
Get the sensor ID.
Definition: SVDCluster.h:102
bool isUCluster() const
Get the direction of strips.
Definition: SVDCluster.h:110
float getPosition(double v=0) const
Get the coordinate of reconstructed hit.
Definition: SVDCluster.h:117
float getClsTimeSigma() const
Get the error of the reconstructed hit time.
Definition: SVDCluster.h:139
float getPositionSigma() const
Get the error of the reconstructed hit coordinate.
Definition: SVDCluster.h:129
int getFirstFrame() const
Get firstFrame of the MaxSum algorithm.
Definition: SVDCluster.h:169
float m_svdTrkPrimeOS
tan of incident angle projected on v/u,w (other side)
float m_svdClIntStrPos
cluster interstrip position
float m_svdTrkPosErrOS
track position error on the other side
TTree * m_t_V
tree containing info related to the V side clusters
void initialize() override
Register input and output data.
std::vector< float > m_svdStripPosition
absolute position of the strips of the cluster
float m_svdClPosErr
cluster position error
int m_svdTrkCDCHits
number of CDC hits on the track
float m_svdTrkPrime
tan of incident angle projected on u/v,w
void event() override
Compute the variables and fill the tree.
int m_svdTrkPXDHits
number of PXD hits on the track
void terminate() override
Write the TTrees to the file.
float m_cdcEventT0_3SRF
CDC event T0 in the 3-sample SVD ref frame.
std::vector< float > m_svdStripCharge
charge of the strips of the cluster
void beginRun() override
Compute the APV clock period.
float m_svdTrkTraversedLength
traversed length of the track in the sensor
float m_svdTrkPosUnbiased
unbiased track position
DBObjPtr< HardwareClockSettings > m_hwClock
Hardware Clocks.
std::vector< float > m_svdStrip6Samples
6 samples of the strips of the cluster
float m_svdRes
residual computed by genfit
StoreObjPtr< EventT0 > m_EventT0
event T0
float m_svdTrkPosErrUnbiased
unbiased track position error
TTree * m_t_U
tree containing info related to the U side clusters
float m_svdTrkPosOS
track position on the other side
float m_svdClTime_6SRF
cluster time in the 6-sample SVD ref frame
std::vector< float > m_svdStripTime
time of the strips of the cluster
float m_svdTrkPrimeUnbiased
unbiased tan of incident angle projected on u,w
float m_cdcEventT0_6SRF
CDC event T0 in the 6-sample SVD ref frame.
std::string m_recoTracksStoreArrayName
storeArray name of the input and output RecoTracks
TFile * m_rootFilePtr
pointer at root file used for storing histograms
float m_svdClTime_3SRF
cluster time in the 3-sample SVD ref frame
int m_svdTrkSVDHits
number of SVD hits on the track
The SVD ShaperDigit class.
APVFloatSamples getSamples() const
Get array of samples.
bool isOptional(const std::string &name="")
Tell the DataStore about an optional input.
Accessor to arrays stored in the data store.
Definition: StoreArray.h:113
Type-safe access to single objects in the data store.
Definition: StoreObjPtr.h:96
bool isValid() const
Check whether the object was created.
Definition: StoreObjPtr.h:111
Values of the result of a track fit with a given particle hypothesis.
ROOT::Math::PxPyPzEVector get4Momentum() const
Getter for the 4Momentum at the closest approach of the track in the r/phi projection.
double getD0() const
Getter for d0.
double getZ0() const
Getter for z0.
ROOT::Math::XYZVector getMomentum() const
Getter for vector of momentum at closest approach of track in r/phi projection.
Class to faciliate easy access to sensor information of the VXD like coordinate transformations or pi...
Definition: GeoCache.h:39
const SensorInfoBase & getSensorInfo(Belle2::VxdID id) const
Return a referecne to the SensorInfo of a given SensorID.
Definition: GeoCache.cc:67
static GeoCache & getInstance()
Return a reference to the singleton instance.
Definition: GeoCache.cc:214
Base class to provide Sensor Information for PXD and SVD.
double getVCellPosition(int vID) const
Return the position of a specific strip/pixel in v direction.
double getUPitch(double v=0) const
Return the pitch of the sensor.
double getVSize() const
Return the length of the sensor.
double getUCellPosition(int uID, int vID=-1) const
Return the position of a specific strip/pixel in u direction.
ROOT::Math::XYZVector pointToGlobal(const ROOT::Math::XYZVector &local, bool reco=false) const
Convert a point from local to global coordinates.
double getVPitch(double v=0) const
Return the pitch of the sensor.
double getThickness() const
Return the thickness of the sensor.
double getUSize(double v=0) const
Return the width of the sensor.
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:33
baseType getSensorNumber() const
Get the sensor id.
Definition: VxdID.h:100
baseType getLadderNumber() const
Get the ladder id.
Definition: VxdID.h:98
baseType getLayerNumber() const
Get the layer id.
Definition: VxdID.h:96
void addParam(const std::string &name, T &paramVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
Definition: Module.h:560
static double convertValueToUnit(double value, const std::string &unitString)
Converts a floating point value from the standard framework unit to the given unit.
Definition: UnitConst.cc:139
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28
Abstract base class for different kinds of events.
STL namespace.