Belle II Software  release-08-01-10
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/Boost.h>
26 #include <math.h>
27 #include <iostream>
28 #include <algorithm>
29 #include <mdst/dataobjects/Track.h>
30 #include <mdst/dataobjects/TrackFitResult.h>
31 
32 
33 using namespace Belle2;
34 using namespace std;
35 
36 
37 //-----------------------------------------------------------------
38 // Register the Module
39 //-----------------------------------------------------------------
40 REG_MODULE(SVDPerformanceTTree);
41 
42 //-----------------------------------------------------------------
43 // Implementation
44 //-----------------------------------------------------------------
45 
47 {
48  //Set module properties
49  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.");
50  //Parameter to take only specific RecoTracks as input
51  addParam("outputFileName", m_rootFileName, "Name of output root file.", std::string("SVDPerformanceTTree.root"));
52  addParam("recoTracksStoreArrayName", m_recoTracksStoreArrayName, "StoreArray name of the input and output RecoTracks.",
54 }
55 
57 {
59  recoTracks.isOptional();
60  m_EventT0.isOptional();
61 
62  m_rootFilePtr = new TFile(m_rootFileName.c_str(), "RECREATE");
63 
64  //Tree for SVD u overlapping clusters
65  m_t_U = new TTree("t_U", "Tree for SVD u-clusters");
66  m_t_U->Branch("cdcEventT0", &m_cdcEventT0, "cdcEventT0/F");
67  m_t_U->Branch("cdcEventT0_6SRF", &m_cdcEventT0_6SRF, "cdcEventT0_6SRF/F");
68  m_t_U->Branch("cdcEventT0_3SRF", &m_cdcEventT0_3SRF, "cdcEventT0_3SRF/F");
69  m_t_U->Branch("cdcEventT0Err", &m_cdcEventT0Err, "cdcEventT0Err/F");
70  m_t_U->Branch("svdTB", &m_svdTB, "svdTB/i");
71  m_t_U->Branch("svdClSNR", &m_svdClSNR, "svdClSNR/F");
72  m_t_U->Branch("svdClCharge", &m_svdClCharge, "svdClCharge/F");
73  m_t_U->Branch("svdStripCharge", &m_svdStripCharge);
74  m_t_U->Branch("svdStrip6Samples", &m_svdStrip6Samples);
75  m_t_U->Branch("svdClTime", &m_svdClTime, "svdClTime/F");
76  m_t_U->Branch("svdClTimeErr", &m_svdClTimeErr, "svdClTimeErr/F");
77  m_t_U->Branch("svdClTime_6SRF", &m_svdClTime_6SRF, "svdClTime_6SRF/F");
78  m_t_U->Branch("svdClTime_3SRF", &m_svdClTime_3SRF, "svdClTime_3SRF/F");
79  m_t_U->Branch("svdFF", &m_svdFF, "svdFF/i");
80  m_t_U->Branch("svdStripTime", &m_svdStripTime);
81  m_t_U->Branch("svdStripPosition", &m_svdStripPosition);
82  m_t_U->Branch("svdRes", &m_svdRes, "svdRes/F");
83  m_t_U->Branch("svdPitch", &m_svdPitch, "svdPitch/F");
84  m_t_U->Branch("svdWidth", &m_svdWidth, "svdWidth/F");
85  m_t_U->Branch("svdLength", &m_svdLength, "svdLength/F");
86  m_t_U->Branch("svdClIntStrPos", &m_svdClIntStrPos, "svdClIntStrPos/F");
87  m_t_U->Branch("svdClPos", &m_svdClPos, "svdClPos/F");
88  m_t_U->Branch("svdClPosErr", &m_svdClPosErr, "svdClPosErr/F");
89  m_t_U->Branch("svdTruePos", &m_svdTruePos, "svdTruePos/F");
90  m_t_U->Branch("svdClPhi", &m_svdClPhi, "svdClPhi/F");
91  m_t_U->Branch("svdClZ", &m_svdClZ, "svdClZ/F");
92  m_t_U->Branch("svdTrkd0", &m_svdTrkd0, "svdTrkd0/F");
93  m_t_U->Branch("svdTrkz0", &m_svdTrkz0, "svdTrkz0/F");
94  m_t_U->Branch("svdTrkpT", &m_svdTrkpT, "svdTrkpT/F");
95  m_t_U->Branch("svdTrkpCM", &m_svdTrkpCM, "svdTrkpCM/F");
96  m_t_U->Branch("svdTrkTraversedLength", &m_svdTrkTraversedLength, "svdTrkTraversedLength/F");
97  m_t_U->Branch("svdTrkPXDHits", &m_svdTrkPXDHits, "svdTrkPXDHits/i");
98  m_t_U->Branch("svdTrkSVDHits", &m_svdTrkSVDHits, "svdTrkSVDHits/i");
99  m_t_U->Branch("svdTrkCDCHits", &m_svdTrkCDCHits, "svdTrkCDCHits/i");
100  m_t_U->Branch("svdTrkPos", &m_svdTrkPos, "svdTrkPos/F");
101  m_t_U->Branch("svdTrkPosOS", &m_svdTrkPosOS, "svdTrkPosOS/F");
102  m_t_U->Branch("svdTrkPosErr", &m_svdTrkPosErr, "svdTrkPosErr/F");
103  m_t_U->Branch("svdTrkPosErrOS", &m_svdTrkPosErrOS, "svdTrkPosErrOS/F");
104  m_t_U->Branch("svdTrkQoP", &m_svdTrkQoP, "svdTrkQoP/F");
105  m_t_U->Branch("svdTrkPrime", &m_svdTrkPrime, "svdTrkPrime/F");
106  m_t_U->Branch("svdTrkPrimeOS", &m_svdTrkPrimeOS, "svdTrkPrimeOS/F");
107  m_t_U->Branch("svdTrkPosUnbiased", &m_svdTrkPosUnbiased, "svdTrkPosUnbiased/F");
108  m_t_U->Branch("svdTrkPosErrUnbiased", &m_svdTrkPosErrUnbiased, "svdTrkPosErrUnbiased/F");
109  m_t_U->Branch("svdTrkQoPUnbiased", &m_svdTrkQoPUnbiased, "svdTrkQoPUnbiased/F");
110  m_t_U->Branch("svdTrkPrimeUnbiased", &m_svdTrkPrimeUnbiased, "svdTrkPrimeUnbiased/F");
111  m_t_U->Branch("svdLayer", &m_svdLayer, "svdLayer/i");
112  m_t_U->Branch("svdLadder", &m_svdLadder, "svdLadder/i");
113  m_t_U->Branch("svdSensor", &m_svdSensor, "svdSensor/i");
114  m_t_U->Branch("svdSize", &m_svdSize, "svdSize/i");
115  //Tree for SVD v overlapping clusters
116  m_t_V = new TTree("t_V", "Tree for SVD v-clusters");
117  m_t_V->Branch("cdcEventT0", &m_cdcEventT0, "cdcEventT0/F");
118  m_t_V->Branch("cdcEventT0_6SRF", &m_cdcEventT0_6SRF, "cdcEventT0_6SRF/F");
119  m_t_V->Branch("cdcEventT0_3SRF", &m_cdcEventT0_3SRF, "cdcEventT0_3SRF/F");
120  m_t_V->Branch("cdcEventT0Err", &m_cdcEventT0Err, "cdcEventT0Err/F");
121  m_t_V->Branch("svdTB", &m_svdTB, "svdTB/i");
122  m_t_V->Branch("svdClSNR", &m_svdClSNR, "svdClSNR/F");
123  m_t_V->Branch("svdClCharge", &m_svdClCharge, "svdClCharge/F");
124  m_t_V->Branch("svdStripCharge", &m_svdStripCharge);
125  m_t_V->Branch("svdStrip6Samples", &m_svdStrip6Samples);
126  m_t_V->Branch("svdClTime", &m_svdClTime, "svdClTime/F");
127  m_t_V->Branch("svdClTimeErr", &m_svdClTimeErr, "svdClTimeErr/F");
128  m_t_V->Branch("svdClTime_6SRF", &m_svdClTime_6SRF, "svdClTime_6SRF/F");
129  m_t_V->Branch("svdClTime_3SRF", &m_svdClTime_3SRF, "svdClTime_3SRF/F");
130  m_t_V->Branch("svdFF", &m_svdFF, "svdFF/i");
131  m_t_V->Branch("svdStripTime", &m_svdStripTime);
132  m_t_V->Branch("svdStripPosition", &m_svdStripPosition);
133  m_t_V->Branch("svdRes", &m_svdRes, "svdRes/F");
134  m_t_V->Branch("svdPitch", &m_svdPitch, "svdPitch/F");
135  m_t_V->Branch("svdWidth", &m_svdWidth, "svdWidth/F");
136  m_t_V->Branch("svdLength", &m_svdLength, "svdLength/F");
137  m_t_V->Branch("svdClIntStrPos", &m_svdClIntStrPos, "svdClIntStrPos/F");
138  m_t_V->Branch("svdClPos", &m_svdClPos, "svdClPos/F");
139  m_t_V->Branch("svdClPosErr", &m_svdClPosErr, "svdClPosErr/F");
140  m_t_V->Branch("svdTruePos", &m_svdTruePos, "svdTruePos/F");
141  m_t_V->Branch("svdClPhi", &m_svdClPhi, "svdClPhi/F");
142  m_t_V->Branch("svdClZ", &m_svdClZ, "svdClZ/F");
143  m_t_V->Branch("svdTrkd0", &m_svdTrkd0, "svdTrkd0/F");
144  m_t_V->Branch("svdTrkz0", &m_svdTrkz0, "svdTrkz0/F");
145  m_t_V->Branch("svdTrkpT", &m_svdTrkpT, "svdTrkpT/F");
146  m_t_V->Branch("svdTrkpCM", &m_svdTrkpCM, "svdTrkpCM/F");
147  m_t_V->Branch("svdTrkTraversedLength", &m_svdTrkTraversedLength, "svdTrkTraversedLength/F");
148  m_t_V->Branch("svdTrkPXDHits", &m_svdTrkPXDHits, "svdTrkPXDHits/i");
149  m_t_V->Branch("svdTrkSVDHits", &m_svdTrkSVDHits, "svdTrkSVDHits/i");
150  m_t_V->Branch("svdTrkCDCHits", &m_svdTrkCDCHits, "svdTrkCDCHits/i");
151  m_t_V->Branch("svdTrkPos", &m_svdTrkPos, "svdTrkPos/F");
152  m_t_V->Branch("svdTrkPosOS", &m_svdTrkPosOS, "svdTrkPosOS/F");
153  m_t_V->Branch("svdTrkPosErr", &m_svdTrkPosErr, "svdTrkPosErr/F");
154  m_t_V->Branch("svdTrkPosErrOS", &m_svdTrkPosErrOS, "svdTrkPosErrOS/F");
155  m_t_V->Branch("svdTrkQoP", &m_svdTrkQoP, "svdTrkQoP/F");
156  m_t_V->Branch("svdTrkPrime", &m_svdTrkPrime, "svdTrkPrime/F");
157  m_t_V->Branch("svdTrkPrimeOS", &m_svdTrkPrimeOS, "svdTrkPrimeOS/F");
158  m_t_V->Branch("svdTrkPosUnbiased", &m_svdTrkPosUnbiased, "svdTrkPosUnbiased/F");
159  m_t_V->Branch("svdTrkPosErrUnbiased", &m_svdTrkPosErrUnbiased, "svdTrkPosErrUnbiased/F");
160  m_t_V->Branch("svdTrkQoPUnbiased", &m_svdTrkQoPUnbiased, "svdTrkQoPUnbiased/F");
161  m_t_V->Branch("svdTrkPrimeUnbiased", &m_svdTrkPrimeUnbiased, "svdTrkPrimeUnbiased/F");
162  m_t_V->Branch("svdLayer", &m_svdLayer, "svdLayer/i");
163  m_t_V->Branch("svdLadder", &m_svdLadder, "svdLadder/i");
164  m_t_V->Branch("svdSensor", &m_svdSensor, "svdSensor/i");
165  m_t_V->Branch("svdSize", &m_svdSize, "svdSize/i");
166 
167 }
168 
170 {
172  m_apvClockPeriod = 1. / m_hwClock->getClockFrequency(Const::EDetector::SVD, "sampling");
173 }
175 {
176 
177  m_cdcEventT0 = std::numeric_limits<float>::quiet_NaN();
178  m_cdcEventT0_6SRF = std::numeric_limits<float>::quiet_NaN();
179  m_cdcEventT0_3SRF = std::numeric_limits<float>::quiet_NaN();
180  m_cdcEventT0Err = std::numeric_limits<float>::quiet_NaN();
181 
182  if (m_EventT0.isValid())
183  if (m_EventT0->hasEventT0()) {
184  if (m_EventT0->hasTemporaryEventT0(Const::EDetector::CDC)) {
185  const auto bestCDCEvtT0 = m_EventT0->getBestCDCTemporaryEventT0();
186  m_cdcEventT0 = bestCDCEvtT0->eventT0 ;
187  m_cdcEventT0Err = bestCDCEvtT0->eventT0Uncertainty;
188  }
189  }
190 
191  //first check SVDEventInfo name
192  StoreObjPtr<SVDEventInfo> temp_eventinfo("SVDEventInfo");
193  std::string m_svdEventInfoName = "SVDEventInfo";
194  if (!temp_eventinfo.isValid())
195  m_svdEventInfoName = "SVDEventInfoSim";
196  StoreObjPtr<SVDEventInfo> eventinfo(m_svdEventInfoName);
197  if (!eventinfo) B2ERROR("No SVDEventInfo!");
198  m_svdTB = eventinfo->getModeByte().getTriggerBin();
199 
200  bool isMC = Environment::Instance().isMC();
201 
204  for (const auto& trk : recoTracks) {
205  if (! trk.wasFitSuccessful()) {
206  continue;
207  }
208 
209 
210  RelationVector<Track> theTK = DataStore::getRelationsWithObj<Track>(&trk);
211 
212  if (theTK.size() == 0) {
213  continue;
214  }
215 
216 
217  const TrackFitResult* tfr = theTK[0]->getTrackFitResultWithClosestMass(Const::pion);
218 
219  if (tfr) {
220  m_svdTrkd0 = tfr->getD0();
221  m_svdTrkz0 = tfr->getZ0();
222  m_svdTrkpT = tfr->getMomentum().Rho();
223  ROOT::Math::PxPyPzEVector pStar = tfr->get4Momentum();
224  ROOT::Math::BoostZ boost(3. / 11);
225  pStar = boost(pStar);
226  m_svdTrkpCM = pStar.P();
227  }
228 
229 
230  const vector<SVDCluster* > svdClusters = trk.getSVDHitList();
231  B2DEBUG(40, "FITTED TRACK: NUMBER OF SVD HITS = " << svdClusters.size());
232 
233  m_svdTrkPXDHits = (trk.getPXDHitList()).size();
234  m_svdTrkSVDHits = (trk.getSVDHitList()).size();
235  m_svdTrkCDCHits = (trk.getCDCHitList()).size();
236 
237  for (unsigned int i = 0; i < svdClusters.size(); i++) {
238 
239  const SVDCluster* svd_1 = svdClusters[i];
240 
241  //get true hits, used only if isMC
242  RelationVector<SVDTrueHit> trueHit_1 = DataStore::getRelationsWithObj<SVDTrueHit>(svd_1);
243 
244  const RecoHitInformation* infoSVD_1 = trk.getRecoHitInformation(svd_1);
245  if (!infoSVD_1) {
246  continue;
247  }
248  const auto* hitTrackPoint_1 = trk.getCreatedTrackPoint(infoSVD_1);
249  const auto* fittedResult_1 = hitTrackPoint_1->getFitterInfo();
250  if (!fittedResult_1) {
251  continue;
252  }
253  const VxdID svd_id_1 = svd_1->getSensorID();
254  const unsigned short svd_Layer_1 = svd_id_1.getLayerNumber();
255  const unsigned short svd_Ladder_1 = svd_id_1.getLadderNumber();
256  const unsigned short svd_Sensor_1 = svd_id_1.getSensorNumber();
257 
258  try {
259  const TVectorD resUnBias_1 = fittedResult_1->getResidual(0, false).getState();
260  genfit::MeasuredStateOnPlane state_unbiased = fittedResult_1->getFittedState(false);
261  const TVectorD& svd_predIntersect_unbiased = state_unbiased.getState();
262  const TMatrixDSym& covMatrix_unbiased = state_unbiased.getCov();
263  genfit::MeasuredStateOnPlane state_1 = trk.getMeasuredStateOnPlaneFromRecoHit(infoSVD_1);
264  const TVectorD& svd_predIntersect_1 = state_1.getState();
265  const TMatrixDSym& covMatrix_1 = state_1.getCov();
266 
267  if (svd_1->isUCluster()) {
268 
269  const int strips_1 = svd_1->getSize();
270 
271  const double res_U_1 = resUnBias_1.GetMatrixArray()[0] * Unit::convertValueToUnit(1.0, "um");
272  const ROOT::Math::XYZVector svdLocal_1(svd_1->getPosition(), svd_predIntersect_1[4], 0.);
273  const VXD::SensorInfoBase& svdSensor_1 = geo.get(svd_id_1);
274  const ROOT::Math::XYZVector& svdGlobal_1 = svdSensor_1.pointToGlobal(svdLocal_1);
275  double svdPhi_1 = atan2(svdGlobal_1.Y(), svdGlobal_1.X()); // maybe use svdGlobal_1.Phi()
276  double svdZ_1 = svdGlobal_1.Z();
277 
278  m_svdFF = svd_1->getFirstFrame();
279  //Fill SVD tree for u-overlaps if required by the user
280  m_svdRes = res_U_1;
281  m_svdClTime = svd_1->getClsTime();
282  m_svdClTimeErr = svd_1->getClsTimeSigma();
283  m_svdClSNR = svd_1->getSNR();
284  m_svdClCharge = svd_1->getCharge();
285  m_svdClPosErr = svd_1->getPositionSigma();
286  if (isMC && trueHit_1.size() > 0)
287  m_svdTruePos = trueHit_1[0]->getU();
288  else
289  m_svdTruePos = -99;
290  m_svdClPhi = svdPhi_1;
291  m_svdClZ = svdZ_1;
292  m_svdTrkPos = svd_predIntersect_1[3];
293  m_svdTrkPosOS = svd_predIntersect_1[4];
294  m_svdTrkPosErr = sqrt(covMatrix_1[3][3]);
295  m_svdTrkPosErrOS = sqrt(covMatrix_1[4][4]);
296  m_svdTrkQoP = svd_predIntersect_1[0];
297  m_svdTrkPrime = svd_predIntersect_1[1];
298  m_svdTrkPrimeOS = svd_predIntersect_1[2];
300  m_svdTrkPosUnbiased = svd_predIntersect_unbiased[3];
302  m_svdTrkPosErrUnbiased = sqrt(covMatrix_unbiased[3][3]);
303  m_svdTrkQoPUnbiased = svd_predIntersect_unbiased[0];
304  m_svdTrkPrimeUnbiased = svd_predIntersect_unbiased[1];
305  m_svdLayer = svd_Layer_1;
306  m_svdLadder = svd_Ladder_1;
307  m_svdSensor = svd_Sensor_1;
308  m_svdSize = strips_1;
309 
310  m_svdPitch = svdSensor_1.getUPitch(m_svdTrkPosOS);
311  m_svdWidth = svdSensor_1.getUSize(m_svdTrkPosOS);
312  m_svdLength = svdSensor_1.getVSize();
313 
315 
316  m_svdStripCharge.clear();
317  m_svdStripTime.clear();
318  m_svdStripPosition.clear();
319  m_svdStrip6Samples.clear();
320  //retrieve relations and set strip charges and times
321  RelationVector<SVDRecoDigit> theRecoDigits = DataStore::getRelationsWithObj<SVDRecoDigit>(svd_1);
322  if ((theRecoDigits.size() != m_svdSize) && (m_svdSize != 128)) //virtual cluster
323  B2ERROR(" Inconsistency with cluster size! # recoDigits = " << theRecoDigits.size() << " != " << m_svdSize << " cluster size");
324 
325  //skip clusters created beacuse of missing APV
326  if (m_svdSize < 128)
327  for (unsigned int d = 0; d < m_svdSize; d++) {
328 
329  SVDShaperDigit* ShaperDigit = theRecoDigits[d]->getRelated<SVDShaperDigit>();
330  array<float, 6> Samples = ShaperDigit->getSamples();
331 
332  m_svdStripCharge.push_back(theRecoDigits[d]->getCharge());
333  std::copy(std::begin(Samples), std::end(Samples), std::back_inserter(m_svdStrip6Samples));
334  m_svdStripTime.push_back(theRecoDigits[d]->getTime());
335  double misalignedStripPos = svdSensor_1.getUCellPosition(theRecoDigits[d]->getCellID());
336  //aligned strip pos = misaligned strip - ( misaligned cluster - aligned cluster)
337  m_svdStripPosition.push_back(misalignedStripPos - svd_1->getPosition() + m_svdClPos);
338  }
339 
340  m_cdcEventT0_3SRF = eventinfo->getTimeInSVDReference(m_cdcEventT0, m_svdFF);
342  m_svdClTime_3SRF = eventinfo->getTimeInSVDReference(m_svdClTime, m_svdFF);
344 
345  m_t_U->Fill();
346 
347  } else {
348  const int strips_1 = svd_1->getSize();
349  const double res_V_1 = resUnBias_1.GetMatrixArray()[0] * Unit::convertValueToUnit(1.0, "um");
350  const ROOT::Math::XYZVector svdLocal_1(svd_predIntersect_1[3], svd_1->getPosition(), 0.);
351  const VXD::SensorInfoBase& svdSensor_1 = geo.get(svd_id_1);
352  const ROOT::Math::XYZVector& svdGlobal_1 = svdSensor_1.pointToGlobal(svdLocal_1);
353  double svdPhi_1 = atan2(svdGlobal_1.Y(), svdGlobal_1.X()); // maybe use svdGlobal_1.Phi()
354  double svdZ_1 = svdGlobal_1.Z();
355 
356  m_svdFF = svd_1->getFirstFrame();
357 
358  m_svdRes = res_V_1;
359  m_svdClTime = svd_1->getClsTime();
360  m_svdClTimeErr = svd_1->getClsTimeSigma();
361  m_svdClSNR = svd_1->getSNR();
362  m_svdClCharge = svd_1->getCharge();
363  m_svdClPosErr = svd_1->getPositionSigma();
364  if (isMC && trueHit_1.size() > 0)
365  m_svdTruePos = trueHit_1[0]->getV();
366  else
367  m_svdTruePos = -99;
368  m_svdClPhi = svdPhi_1;
369  m_svdClZ = svdZ_1;
370  m_svdTrkPos = svd_predIntersect_1[4];
371  m_svdTrkPosOS = svd_predIntersect_1[3];
372  m_svdTrkPosErr = sqrt(covMatrix_1[4][4]);
373  m_svdTrkPosErrOS = sqrt(covMatrix_1[3][3]);
374  m_svdTrkQoP = svd_predIntersect_1[0];
375  m_svdTrkPrime = svd_predIntersect_1[2];
376  m_svdTrkPrimeOS = svd_predIntersect_1[1];
378  m_svdTrkPosUnbiased = svd_predIntersect_unbiased[4];
380  m_svdTrkPosErrUnbiased = sqrt(covMatrix_unbiased[4][4]);
381  m_svdTrkQoPUnbiased = svd_predIntersect_unbiased[0];
382  m_svdTrkPrimeUnbiased = svd_predIntersect_unbiased[2];
383  m_svdLayer = svd_Layer_1;
384  m_svdLadder = svd_Ladder_1;
385  m_svdSensor = svd_Sensor_1;
386  m_svdSize = strips_1;
387 
388  m_svdPitch = svdSensor_1.getVPitch();
389  m_svdWidth = svdSensor_1.getUSize(m_svdTrkPos);
390  m_svdLength = svdSensor_1.getVSize();
391 
393 
394  m_svdStripCharge.clear();
395  m_svdStripTime.clear();
396  m_svdStripPosition.clear();
397  m_svdStrip6Samples.clear();
398  //retrieve relations and set strip charges and times
399  RelationVector<SVDRecoDigit> theRecoDigits = DataStore::getRelationsWithObj<SVDRecoDigit>(svd_1);
400  if ((theRecoDigits.size() != m_svdSize) && (m_svdSize != 128)) //virtual cluster
401  B2ERROR(" Inconsistency with cluster size! # recoDigits = " << theRecoDigits.size() << " != " << m_svdSize << " cluster size");
402 
403  //skip clusters created beacuse of missing APV
404  if (m_svdSize < 128)
405  for (unsigned int d = 0; d < m_svdSize; d++) {
406  SVDShaperDigit* ShaperDigit = theRecoDigits[d]->getRelated<SVDShaperDigit>();
407  array<float, 6> Samples = ShaperDigit->getSamples();
408  m_svdStripCharge.push_back(theRecoDigits[d]->getCharge());
409  std::copy(std::begin(Samples), std::end(Samples), std::back_inserter(m_svdStrip6Samples));
410  m_svdStripTime.push_back(theRecoDigits[d]->getTime());
411  double misalignedStripPos = svdSensor_1.getVCellPosition(theRecoDigits[d]->getCellID());
412  //Aligned strip pos = misaligned strip - ( misaligned cluster - aligned cluster)
413  m_svdStripPosition.push_back(misalignedStripPos - svd_1->getPosition() + m_svdClPos);
414  }
415 
416  m_cdcEventT0_3SRF = eventinfo->getTimeInSVDReference(m_cdcEventT0, m_svdFF);
418  m_svdClTime_3SRF = eventinfo->getTimeInSVDReference(m_svdClTime, m_svdFF);
420 
421 
422  m_t_V->Fill();
423  }
424  } catch (...) {
425  B2INFO("oops...something went wrong in getting the unbiased state, skipping this cluster.");
426  continue;
427  }
428 
429  }
430  }
431 }
432 
433 
435 {
436 
437  if (m_rootFilePtr != nullptr) {
438  m_rootFilePtr->cd();
439  m_t_U->Write();
440  m_t_V->Write();
441  m_rootFilePtr->Close();
442  }
443 }
static const ChargedStable pion
charged pion particle
Definition: Const.h:652
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.
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
static GeoCache & getInstance()
Return a reference to the singleton instance.
Definition: GeoCache.cc:214
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.
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
#StateOnPlane with additional covariance matrix.
REG_MODULE(arichBtest)
Register the Module.
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
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28
Abstract base class for different kinds of events.