Belle II Software  release-06-01-15
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 <TVector3.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 
32 using namespace Belle2;
33 using namespace std;
34 
35 
36 //-----------------------------------------------------------------
37 // Register the Module
38 //-----------------------------------------------------------------
39 REG_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.",
52  m_recoTracksStoreArrayName);
53 }
54 
56 {
57  StoreArray<RecoTrack> recoTracks(m_recoTracksStoreArrayName);
58  recoTracks.isOptional();
59 
60  m_rootFilePtr = new TFile(m_rootFileName.c_str(), "RECREATE");
61 
62  //Tree for SVD u overlapping clusters
63  m_t_U = new TTree("t_U", "Tree for SVD u-clusters");
64  m_t_U->Branch("svdClSNR", &m_svdClSNR, "svdClSNR/F");
65  m_t_U->Branch("svdClCharge", &m_svdClCharge, "svdClCharge/F");
66  m_t_U->Branch("svdStripCharge", &m_svdStripCharge);
67  m_t_U->Branch("svdStrip6Samples", &m_svdStrip6Samples);
68  m_t_U->Branch("svdClTime", &m_svdClTime, "svdClTime/F");
69  m_t_U->Branch("svdStripTime", &m_svdStripTime);
70  m_t_U->Branch("svdStripPosition", &m_svdStripPosition);
71  m_t_U->Branch("svdRes", &m_svdRes, "svdRes/F");
72  m_t_U->Branch("svdPitch", &m_svdPitch, "svdPitch/F");
73  m_t_U->Branch("svdWidth", &m_svdWidth, "svdWidth/F");
74  m_t_U->Branch("svdLength", &m_svdLength, "svdLength/F");
75  m_t_U->Branch("svdClIntStrPos", &m_svdClIntStrPos, "svdClIntStrPos/F");
76  m_t_U->Branch("svdClPos", &m_svdClPos, "svdClPos/F");
77  m_t_U->Branch("svdClPosErr", &m_svdClPosErr, "svdClPosErr/F");
78  m_t_U->Branch("svdTruePos", &m_svdTruePos, "svdTruePos/F");
79  m_t_U->Branch("svdClPhi", &m_svdClPhi, "svdClPhi/F");
80  m_t_U->Branch("svdClZ", &m_svdClZ, "svdClZ/F");
81  m_t_U->Branch("svdTrkd0", &m_svdTrkd0, "svdTrkd0/F");
82  m_t_U->Branch("svdTrkz0", &m_svdTrkz0, "svdTrkz0/F");
83  m_t_U->Branch("svdTrkpT", &m_svdTrkpT, "svdTrkpT/F");
84  m_t_U->Branch("svdTrkpCM", &m_svdTrkpCM, "svdTrkpCM/F");
85  m_t_U->Branch("svdTrkTraversedLength", &m_svdTrkTraversedLength, "svdTrkTraversedLength/F");
86  m_t_U->Branch("svdTrkPXDHits", &m_svdTrkPXDHits, "svdTrkPXDHits/i");
87  m_t_U->Branch("svdTrkSVDHits", &m_svdTrkSVDHits, "svdTrkSVDHits/i");
88  m_t_U->Branch("svdTrkCDCHits", &m_svdTrkCDCHits, "svdTrkCDCHits/i");
89  m_t_U->Branch("svdTrkPos", &m_svdTrkPos, "svdTrkPos/F");
90  m_t_U->Branch("svdTrkPosOS", &m_svdTrkPosOS, "svdTrkPosOS/F");
91  m_t_U->Branch("svdTrkPosErr", &m_svdTrkPosErr, "svdTrkPosErr/F");
92  m_t_U->Branch("svdTrkPosErrOS", &m_svdTrkPosErrOS, "svdTrkPosErrOS/F");
93  m_t_U->Branch("svdTrkQoP", &m_svdTrkQoP, "svdTrkQoP/F");
94  m_t_U->Branch("svdTrkPrime", &m_svdTrkPrime, "svdTrkPrime/F");
95  m_t_U->Branch("svdTrkPrimeOS", &m_svdTrkPrimeOS, "svdTrkPrimeOS/F");
96  m_t_U->Branch("svdTrkPosUnbiased", &m_svdTrkPosUnbiased, "svdTrkPosUnbiased/F");
97  m_t_U->Branch("svdTrkPosErrUnbiased", &m_svdTrkPosErrUnbiased, "svdTrkPosErrUnbiased/F");
98  m_t_U->Branch("svdTrkQoPUnbiased", &m_svdTrkQoPUnbiased, "svdTrkQoPUnbiased/F");
99  m_t_U->Branch("svdTrkPrimeUnbiased", &m_svdTrkPrimeUnbiased, "svdTrkPrimeUnbiased/F");
100  m_t_U->Branch("svdLayer", &m_svdLayer, "svdLayer/i");
101  m_t_U->Branch("svdLadder", &m_svdLadder, "svdLadder/i");
102  m_t_U->Branch("svdSensor", &m_svdSensor, "svdSensor/i");
103  m_t_U->Branch("svdSize", &m_svdSize, "svdSize/i");
104  m_t_U->Branch("svdTB", &m_svdTB, "svdTB/i");
105  //Tree for SVD v overlapping clusters
106  m_t_V = new TTree("t_V", "Tree for M_SVD v-clusters");
107  m_t_V->Branch("svdClSNR", &m_svdClSNR, "svdClSNR/F");
108  m_t_V->Branch("svdClCharge", &m_svdClCharge, "svdClCharge/F");
109  m_t_V->Branch("svdStripCharge", &m_svdStripCharge);
110  m_t_V->Branch("svdStrip6Samples", &m_svdStrip6Samples);
111  m_t_V->Branch("svdClTime", &m_svdClTime, "svdClTime/F");
112  m_t_V->Branch("svdStripTime", &m_svdStripTime);
113  m_t_V->Branch("svdStripPosition", &m_svdStripPosition);
114  m_t_V->Branch("svdRes", &m_svdRes, "svdRes/F");
115  m_t_V->Branch("svdPitch", &m_svdPitch, "svdPitch/F");
116  m_t_V->Branch("svdWidth", &m_svdWidth, "svdWidth/F");
117  m_t_V->Branch("svdLength", &m_svdLength, "svdLength/F");
118  m_t_V->Branch("svdClIntStrPos", &m_svdClIntStrPos, "svdClIntStrPos/F");
119  m_t_V->Branch("svdClPos", &m_svdClPos, "svdClPos/F");
120  m_t_V->Branch("svdClPosErr", &m_svdClPosErr, "svdClPosErr/F");
121  m_t_V->Branch("svdTruePos", &m_svdTruePos, "svdTruePos/F");
122  m_t_V->Branch("svdClPhi", &m_svdClPhi, "svdClPhi/F");
123  m_t_V->Branch("svdClZ", &m_svdClZ, "svdClZ/F");
124  m_t_V->Branch("svdTrkd0", &m_svdTrkd0, "svdTrkd0/F");
125  m_t_V->Branch("svdTrkz0", &m_svdTrkz0, "svdTrkz0/F");
126  m_t_V->Branch("svdTrkpT", &m_svdTrkpT, "svdTrkpT/F");
127  m_t_V->Branch("svdTrkpCM", &m_svdTrkpCM, "svdTrkpCM/F");
128  m_t_V->Branch("svdTrkTraversedLength", &m_svdTrkTraversedLength, "svdTrkTraversedLength/F");
129  m_t_V->Branch("svdTrkPXDHits", &m_svdTrkPXDHits, "svdTrkPXDHits/i");
130  m_t_V->Branch("svdTrkSVDHits", &m_svdTrkSVDHits, "svdTrkSVDHits/i");
131  m_t_V->Branch("svdTrkCDCHits", &m_svdTrkCDCHits, "svdTrkCDCHits/i");
132  m_t_V->Branch("svdTrkPos", &m_svdTrkPos, "svdTrkPos/F");
133  m_t_V->Branch("svdTrkPosOS", &m_svdTrkPosOS, "svdTrkPosOS/F");
134  m_t_V->Branch("svdTrkPosErr", &m_svdTrkPosErr, "svdTrkPosErr/F");
135  m_t_V->Branch("svdTrkPosErrOS", &m_svdTrkPosErrOS, "svdTrkPosErrOS/F");
136  m_t_V->Branch("svdTrkQoP", &m_svdTrkQoP, "svdTrkQoP/F");
137  m_t_V->Branch("svdTrkPrime", &m_svdTrkPrime, "svdTrkPrime/F");
138  m_t_V->Branch("svdTrkPrimeOS", &m_svdTrkPrimeOS, "svdTrkPrimeOS/F");
139  m_t_V->Branch("svdTrkPosUnbiased", &m_svdTrkPosUnbiased, "svdTrkPosUnbiased/F");
140  m_t_V->Branch("svdTrkPosErrUnbiased", &m_svdTrkPosErrUnbiased, "svdTrkPosErrUnbiased/F");
141  m_t_V->Branch("svdTrkQoPUnbiased", &m_svdTrkQoPUnbiased, "svdTrkQoPUnbiased/F");
142  m_t_V->Branch("svdTrkPrimeUnbiased", &m_svdTrkPrimeUnbiased, "svdTrkPrimeUnbiased/F");
143  m_t_V->Branch("svdLayer", &m_svdLayer, "svdLayer/i");
144  m_t_V->Branch("svdLadder", &m_svdLadder, "svdLadder/i");
145  m_t_V->Branch("svdSensor", &m_svdSensor, "svdSensor/i");
146  m_t_V->Branch("svdSize", &m_svdSize, "svdSize/i");
147  m_t_V->Branch("svdTB", &m_svdTB, "svdTB/i");
148 
149 }
150 
152 {
153 
154  //first check SVDEventInfo name
155  StoreObjPtr<SVDEventInfo> temp_eventinfo("SVDEventInfo");
156  std::string m_svdEventInfoName = "SVDEventInfo";
157  if (!temp_eventinfo.isValid())
158  m_svdEventInfoName = "SVDEventInfoSim";
159  StoreObjPtr<SVDEventInfo> eventinfo(m_svdEventInfoName);
160  if (!eventinfo) B2ERROR("No SVDEventInfo!");
161  m_svdTB = eventinfo->getModeByte().getTriggerBin();
162 
163  bool isMC = Environment::Instance().isMC();
164 
166  StoreArray<RecoTrack> recoTracks(m_recoTracksStoreArrayName);
167  for (const auto& trk : recoTracks) {
168  if (! trk.wasFitSuccessful()) {
169  continue;
170  }
171 
172 
173  RelationVector<Track> theTK = DataStore::getRelationsWithObj<Track>(&trk);
174 
175  if (theTK.size() == 0) {
176  continue;
177  }
178 
179 
180  const TrackFitResult* tfr = theTK[0]->getTrackFitResultWithClosestMass(Const::pion);
181 
182  if (tfr) {
183  m_svdTrkd0 = tfr->getD0();
184  m_svdTrkz0 = tfr->getZ0();
185  m_svdTrkpT = tfr->getMomentum().Perp();
186  TLorentzVector pStar = tfr->get4Momentum();
187  pStar.Boost(0, 0, 3. / 11);
188  m_svdTrkpCM = pStar.P();
189  }
190 
191 
192  const vector<SVDCluster* > svdClusters = trk.getSVDHitList();
193  B2DEBUG(40, "FITTED TRACK: NUMBER OF SVD HITS = " << svdClusters.size());
194 
195  m_svdTrkPXDHits = (trk.getPXDHitList()).size();
196  m_svdTrkSVDHits = (trk.getSVDHitList()).size();
197  m_svdTrkCDCHits = (trk.getCDCHitList()).size();
198 
199  for (unsigned int i = 0; i < svdClusters.size(); i++) {
200 
201  const SVDCluster* svd_1 = svdClusters[i];
202 
203  //get true hits, used only if isMC
204  RelationVector<SVDTrueHit> trueHit_1 = DataStore::getRelationsWithObj<SVDTrueHit>(svd_1);
205 
206  const RecoHitInformation* infoSVD_1 = trk.getRecoHitInformation(svd_1);
207  if (!infoSVD_1) {
208  continue;
209  }
210  const auto* hitTrackPoint_1 = trk.getCreatedTrackPoint(infoSVD_1);
211  const auto* fittedResult_1 = hitTrackPoint_1->getFitterInfo();
212  if (!fittedResult_1) {
213  continue;
214  }
215  const VxdID svd_id_1 = svd_1->getSensorID();
216  const unsigned short svd_Layer_1 = svd_id_1.getLayerNumber();
217  const unsigned short svd_Ladder_1 = svd_id_1.getLadderNumber();
218  const unsigned short svd_Sensor_1 = svd_id_1.getSensorNumber();
219 
220  try {
221  const TVectorD resUnBias_1 = fittedResult_1->getResidual(0, false).getState();
222  genfit::MeasuredStateOnPlane state_unbiased = fittedResult_1->getFittedState(false);
223  const TVectorD& svd_predIntersect_unbiased = state_unbiased.getState();
224  const TMatrixDSym& covMatrix_unbiased = state_unbiased.getCov();
225  genfit::MeasuredStateOnPlane state_1 = trk.getMeasuredStateOnPlaneFromRecoHit(infoSVD_1);
226  const TVectorD& svd_predIntersect_1 = state_1.getState();
227  const TMatrixDSym& covMatrix_1 = state_1.getCov();
228 
229  if (svd_1->isUCluster()) {
230 
231  const int strips_1 = svd_1->getSize();
232 
233  const double res_U_1 = resUnBias_1.GetMatrixArray()[0] * Unit::convertValueToUnit(1.0, "um");
234  const TVector3 svdLocal_1(svd_1->getPosition(), svd_predIntersect_1[4], 0.);
235  const VXD::SensorInfoBase& svdSensor_1 = geo.get(svd_id_1);
236  const TVector3& svdGlobal_1 = svdSensor_1.pointToGlobal(svdLocal_1);
237  double svdPhi_1 = atan2(svdGlobal_1(1), svdGlobal_1(0));
238  double svdZ_1 = svdGlobal_1(2);
239 
240  //Fill SVD tree for u-overlaps if required by the user
241  m_svdRes = res_U_1;
242  m_svdClTime = svd_1->getClsTime();
243  m_svdClSNR = svd_1->getSNR();
244  m_svdClCharge = svd_1->getCharge();
245  m_svdClPosErr = svd_1->getPositionSigma();
246  if (isMC && trueHit_1.size() > 0)
247  m_svdTruePos = trueHit_1[0]->getU();
248  else
249  m_svdTruePos = -99;
250  m_svdClPhi = svdPhi_1;
251  m_svdClZ = svdZ_1;
252  m_svdTrkPos = svd_predIntersect_1[3];
253  m_svdTrkPosOS = svd_predIntersect_1[4];
254  m_svdTrkPosErr = sqrt(covMatrix_1[3][3]);
255  m_svdTrkPosErrOS = sqrt(covMatrix_1[4][4]);
256  m_svdTrkQoP = svd_predIntersect_1[0];
257  m_svdTrkPrime = svd_predIntersect_1[1];
258  m_svdTrkPrimeOS = svd_predIntersect_1[2];
259  m_svdTrkTraversedLength = svdSensor_1.getThickness() * sqrt(1 + m_svdTrkPrimeOS * m_svdTrkPrimeOS + m_svdTrkPrime * m_svdTrkPrime);
260  m_svdTrkPosUnbiased = svd_predIntersect_unbiased[3];
261  m_svdClPos = m_svdRes / 1e4 + m_svdTrkPosUnbiased;
262  m_svdTrkPosErrUnbiased = sqrt(covMatrix_unbiased[3][3]);
263  m_svdTrkQoPUnbiased = svd_predIntersect_unbiased[0];
264  m_svdTrkPrimeUnbiased = svd_predIntersect_unbiased[1];
265  m_svdLayer = svd_Layer_1;
266  m_svdLadder = svd_Ladder_1;
267  m_svdSensor = svd_Sensor_1;
268  m_svdSize = strips_1;
269 
270  m_svdPitch = svdSensor_1.getUPitch(m_svdTrkPosOS);
271  m_svdWidth = svdSensor_1.getUSize(m_svdTrkPosOS);
272  m_svdLength = svdSensor_1.getVSize();
273 
274  m_svdClIntStrPos = fmod(m_svdClPos + 0.5 * m_svdWidth, m_svdPitch) / m_svdPitch;
275 
276  m_svdStripCharge.clear();
277  m_svdStripTime.clear();
278  m_svdStripPosition.clear();
279  m_svdStrip6Samples.clear();
280  //retrieve relations and set strip charges and times
281  RelationVector<SVDRecoDigit> theRecoDigits = DataStore::getRelationsWithObj<SVDRecoDigit>(svd_1);
282  if ((theRecoDigits.size() != m_svdSize) && (m_svdSize != 128)) //virtual cluster
283  B2ERROR(" Inconsistency with cluster size! # recoDigits = " << theRecoDigits.size() << " != " << m_svdSize << " cluster size");
284 
285  //skip clusters created beacuse of missing APV
286  if (m_svdSize < 128)
287  for (unsigned int d = 0; d < m_svdSize; d++) {
288 
289  SVDShaperDigit* ShaperDigit = theRecoDigits[d]->getRelated<SVDShaperDigit>();
290  array<float, 6> Samples = ShaperDigit->getSamples();
291 
292  m_svdStripCharge.push_back(theRecoDigits[d]->getCharge());
293  std::copy(std::begin(Samples), std::end(Samples), std::back_inserter(m_svdStrip6Samples));
294  m_svdStripTime.push_back(theRecoDigits[d]->getTime());
295  double misalignedStripPos = svdSensor_1.getUCellPosition(theRecoDigits[d]->getCellID());
296  //aligned strip pos = misaligned strip - ( misaligned cluster - aligned cluster)
297  m_svdStripPosition.push_back(misalignedStripPos - svd_1->getPosition() + m_svdClPos);
298  }
299 
300 
301 
302  m_t_U->Fill();
303 
304  } else {
305  const int strips_1 = svd_1->getSize();
306  const double res_V_1 = resUnBias_1.GetMatrixArray()[0] * Unit::convertValueToUnit(1.0, "um");
307  const TVector3 svdLocal_1(svd_predIntersect_1[3], svd_1->getPosition(), 0.);
308  const VXD::SensorInfoBase& svdSensor_1 = geo.get(svd_id_1);
309  const TVector3& svdGlobal_1 = svdSensor_1.pointToGlobal(svdLocal_1);
310  double svdPhi_1 = atan2(svdGlobal_1(1), svdGlobal_1(0));
311  double svdZ_1 = svdGlobal_1(2);
312 
313  m_svdRes = res_V_1;
314  m_svdClTime = svd_1->getClsTime();
315  m_svdClSNR = svd_1->getSNR();
316  m_svdClCharge = svd_1->getCharge();
317  m_svdClPosErr = svd_1->getPositionSigma();
318  if (isMC && trueHit_1.size() > 0)
319  m_svdTruePos = trueHit_1[0]->getV();
320  else
321  m_svdTruePos = -99;
322  m_svdClPhi = svdPhi_1;
323  m_svdClZ = svdZ_1;
324  m_svdTrkPos = svd_predIntersect_1[4];
325  m_svdTrkPosOS = svd_predIntersect_1[3];
326  m_svdTrkPosErr = sqrt(covMatrix_1[4][4]);
327  m_svdTrkPosErrOS = sqrt(covMatrix_1[3][3]);
328  m_svdTrkQoP = svd_predIntersect_1[0];
329  m_svdTrkPrime = svd_predIntersect_1[2];
330  m_svdTrkPrimeOS = svd_predIntersect_1[1];
331  m_svdTrkTraversedLength = svdSensor_1.getThickness() * sqrt(1 + m_svdTrkPrimeOS * m_svdTrkPrimeOS + m_svdTrkPrime * m_svdTrkPrime);
332  m_svdTrkPosUnbiased = svd_predIntersect_unbiased[4];
333  m_svdClPos = m_svdRes / 1e4 + m_svdTrkPosUnbiased;
334  m_svdTrkPosErrUnbiased = sqrt(covMatrix_unbiased[4][4]);
335  m_svdTrkQoPUnbiased = svd_predIntersect_unbiased[0];
336  m_svdTrkPrimeUnbiased = svd_predIntersect_unbiased[2];
337  m_svdLayer = svd_Layer_1;
338  m_svdLadder = svd_Ladder_1;
339  m_svdSensor = svd_Sensor_1;
340  m_svdSize = strips_1;
341 
342  m_svdPitch = svdSensor_1.getVPitch();
343  m_svdWidth = svdSensor_1.getUSize(m_svdTrkPos);
344  m_svdLength = svdSensor_1.getVSize();
345 
346  m_svdClIntStrPos = fmod(m_svdClPos + 0.5 * m_svdLength, m_svdPitch) / m_svdPitch;
347 
348  m_svdStripCharge.clear();
349  m_svdStripTime.clear();
350  m_svdStripPosition.clear();
351  m_svdStrip6Samples.clear();
352  //retrieve relations and set strip charges and times
353  RelationVector<SVDRecoDigit> theRecoDigits = DataStore::getRelationsWithObj<SVDRecoDigit>(svd_1);
354  if ((theRecoDigits.size() != m_svdSize) && (m_svdSize != 128)) //virtual cluster
355  B2ERROR(" Inconsistency with cluster size! # recoDigits = " << theRecoDigits.size() << " != " << m_svdSize << " cluster size");
356 
357  //skip clusters created beacuse of missing APV
358  if (m_svdSize < 128)
359  for (unsigned int d = 0; d < m_svdSize; d++) {
360  SVDShaperDigit* ShaperDigit = theRecoDigits[d]->getRelated<SVDShaperDigit>();
361  array<float, 6> Samples = ShaperDigit->getSamples();
362  m_svdStripCharge.push_back(theRecoDigits[d]->getCharge());
363  std::copy(std::begin(Samples), std::end(Samples), std::back_inserter(m_svdStrip6Samples));
364  m_svdStripTime.push_back(theRecoDigits[d]->getTime());
365  double misalignedStripPos = svdSensor_1.getVCellPosition(theRecoDigits[d]->getCellID());
366  //Aligned strip pos = misaligned strip - ( misaligned cluster - aligned cluster)
367  m_svdStripPosition.push_back(misalignedStripPos - svd_1->getPosition() + m_svdClPos);
368  }
369 
370  m_t_V->Fill();
371  }
372  } catch (...) {
373  B2INFO("oops...something went wrong in getting the unbiased state, skipping this cluster.");
374  continue;
375  }
376 
377  }
378  }
379 }
380 
381 
383 {
384 
385  if (m_rootFilePtr != nullptr) {
386  m_rootFilePtr->cd();
387  m_t_U->Write();
388  m_t_V->Write();
389  m_rootFilePtr->Close();
390  }
391 }
static const ChargedStable pion
charged pion particle
Definition: Const.h:542
bool isMC() const
Do we have generated, not real data?
Definition: Environment.cc:55
static Environment & Instance()
Static method to get a reference to the Environment instance.
Definition: Environment.cc:29
Base class for Modules.
Definition: Module.h:72
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:28
float getClsTime() const
Get average of waveform maximum times of cluster strip signals.
Definition: SVDCluster.h:133
float getSNR() const
Get cluster SNR.
Definition: SVDCluster.h:158
unsigned short getSize() const
Get cluster size.
Definition: SVDCluster.h:153
float getCharge() const
Get collected charge.
Definition: SVDCluster.h:143
VxdID getSensorID() const
Get the sensor ID.
Definition: SVDCluster.h:101
bool isUCluster() const
Get the direction of strips.
Definition: SVDCluster.h:109
float getPosition(double v=0) const
Get the coordinate of reconstructed hit.
Definition: SVDCluster.h:116
float getPositionSigma() const
Get the error of the reconstructed hit coordinate.
Definition: SVDCluster.h:128
The module is used to create a TTree to study SVD clusters, genfit unbiased residuals and many other ...
void initialize() override
Register input and output data.
void event() override
Compute the variables and fill the tree.
void terminate() override
Write the TTrees to the file.
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:95
bool isValid() const
Check whether the object was created.
Definition: StoreObjPtr.h:110
Values of the result of a track fit with a given particle hypothesis.
TLorentzVector get4Momentum() const
Getter for the 4Momentum at the closest approach of the track in the r/phi projection.
TVector3 getMomentum() const
Getter for vector of momentum at closest approach of track in r/phi projection.
double getD0() const
Getter for d0.
double getZ0() const
Getter for z0.
Class to faciliate easy access to sensor information of the VXD like coordinate transformations or pi...
Definition: GeoCache.h:39
static GeoCache & getInstance()
Return a reference to the singleton instance.
Definition: GeoCache.cc:213
static const SensorInfoBase & get(Belle2::VxdID id)
Return a reference to the SensorInfo of a given SensorID.
Definition: GeoCache.h:139
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.
double getVPitch(double v=0) const
Return the pitch of the sensor.
TVector3 pointToGlobal(const TVector3 &local, bool reco=false) const
Convert a point from local to global coordinates.
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
#StateOnPlane with additional covariance matrix.
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
Abstract base class for different kinds of events.