Belle II Software  release-06-02-00
CDCCRTestModule.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 "cdc/modules/cdcCRTest/CDCCRTestModule.h"
10 #include <TTree.h>
11 #include <framework/gearbox/Const.h>
12 #include <framework/core/HistoModule.h>
13 #include <framework/datastore/StoreObjPtr.h>
14 #include <framework/datastore/RelationArray.h>
15 
16 #include <genfit/TrackPoint.h>
17 #include <genfit/KalmanFitterInfo.h>
18 #include <genfit/MeasurementOnPlane.h>
19 #include <genfit/MeasuredStateOnPlane.h>
20 #include <genfit/StateOnPlane.h>
21 #include <boost/foreach.hpp>
22 
23 #include <cdc/translators/RealisticCDCGeometryTranslator.h>
24 #include <cdc/translators/RealisticTDCCountTranslator.h>
25 #include <cdc/dataobjects/CDCRecoHit.h>
26 #include <cdc/dataobjects/CDCSimHit.h>
27 #include "TDirectory.h"
28 #include <Math/ProbFuncMathCore.h>
29 
30 using namespace std;
31 using namespace Belle2;
32 using namespace CDC;
33 using namespace genfit;
34 //-----------------------------------------------------------------
35 // Register the Module
36 //-----------------------------------------------------------------
37 REG_MODULE(CDCCRTest)
38 
39 // Implementation
41 {
42  setDescription("CDC Cosmic ray test module");
43  setPropertyFlags(c_ParallelProcessingCertified); // specify this flag if you need parallel processing
44  addParam("CorrectToF", m_ToF, "if true, time of flight will take account in t", true);
45  addParam("CorrectToP", m_ToP, "if true, time of Propagation will take account in t", true);
46  addParam("RecoTracksColName", m_recoTrackArrayName, "Name of collection hold genfit::Track", std::string(""));
47  addParam("histogramDirectoryName", m_histogramDirectoryName,
48  "Track fit results histograms will be put into this directory", std::string("trackfit"));
49  addParam("NameOfTree", m_treeName, "name of tree in output file", string("tree"));
50  addParam("fillExpertHistograms", m_fillExpertHistos, "Fill additional histograms", true);
51  addParam("noBFit", m_noBFit, "If true -> #Params ==4, #params ==5 for calculate P-Val", false);
52  addParam("plotResidual", m_plotResidual, "plot biased residual, normalized res and xtplot for all layer", false);
53  addParam("calExpectedDriftTime", m_calExpectedDriftTime, "if true module will calculate expected drift time, it take a time",
54  true);
55  addParam("hitEfficiency", m_hitEfficiency, "calculate hit efficiency(Not work now) true:yes false:No", false);
56  addParam("TriggerPos", m_TriggerPos, "Trigger position use for cut and reconstruct Trigger image", std::vector<double> { -0.6, -13.25, 17.3});
57  addParam("NormTriggerPlaneDirection", m_TriggerPlaneDirection, "Normal trigger plane direction and reconstruct Trigger image",
58  std::vector<double> { 0, 1, 0});
59  addParam("TriggerSize", m_TriggerSize, "Trigger Size, (Width x length)", std::vector<double> {100, 50});
60  addParam("IwireLow", m_low, "Lower boundary of hit dist. Histogram", std::vector<int> {0, 0, 0, 0, 0, 0, 0, 0, 0});
61  addParam("IwireUpper", m_up, "Upper boundary of hit dist. Histogram", std::vector<int> {161, 161, 193, 225, 257, 289, 321, 355, 385});
62  addParam("StoreCDCSimHitInfo", m_StoreCDCSimHitInfo, "Store simulation info related to hit, driftLeng, flight time,z_onwire",
63  false);
64  addParam("EstimateResultForUnFittedLayer", m_EstimateResultForUnFittedLayer,
65  "Calculate residual for Layer that is set unUseInFit", true);
66  addParam("SmallerOutput", m_SmallerOutput, "If true, trigghit position, residual cov,absRes, will not be stored", false);
67  addParam("StoreTrackParams", m_StoreTrackParams, "Store Track Parameter or not, it will be multicount for each hit", true);
68  addParam("StoreHitDistribution", m_MakeHitDist, "Make hit distribution or not", true);
69  addParam("EventT0Extraction", m_EventT0Extraction, "use event t0 extract t0 or not", true);
70  addParam("MinimumPt", m_MinimumPt, "Tracks with tranverse momentum small than this will not recored", 0.);
71 }
72 
73 CDCCRTestModule::~CDCCRTestModule()
74 {
75  m_allHistos.clear();
76 }
77 //------------------------------------------------------------------
78 // Function to define histograms
79 //-----------------------------------------------------------------
80 void CDCCRTestModule::defineHisto()
81 {
82  m_tree = new TTree(m_treeName.c_str(), "tree");
83  m_tree->Branch("x_mea", &x_mea, "x_mea/D");
84  m_tree->Branch("x_u", &x_u, "x_u/D");
85  m_tree->Branch("x_b", &x_b, "x_b/D");
86  m_tree->Branch("z", &z, "z/D");
87  m_tree->Branch("alpha", &alpha, "alpha/D");
88  m_tree->Branch("theta", &theta, "theta/D");
89  m_tree->Branch("t", &t, "t/D");
90  m_tree->Branch("evtT0", &evtT0, "evtT0/D");
91  m_tree->Branch("adc", &adc, "adc/s");
92  m_tree->Branch("boardID", &boardID, "boardID/I");
93  m_tree->Branch("lay", &lay, "lay/I");
94  m_tree->Branch("weight", &weight, "weight/D");
95  m_tree->Branch("IWire", &IWire, "IWire/I");
96  m_tree->Branch("Pval", &Pval, "Pval/D");
97  m_tree->Branch("ndf", &ndf, "ndf/D");
98  // m_tree->Branch("trighit", &trighit, "trighit/I");
99  if (m_StoreTrackParams) {
100  m_tree->Branch("d0", &d0, "d0/D");
101  m_tree->Branch("z0", &z0, "z0/D");
102  m_tree->Branch("phi0", &phi0, "phi0/D");
103  m_tree->Branch("tanL", &tanL, "tanL/D");
104  m_tree->Branch("omega", &omega, "omega/D");
105  m_tree->Branch("Pt", &Pt, "Pt/D");
106  }
107  if (m_StoreCDCSimHitInfo) {
108  m_tree->Branch("z_sim", &z_sim, "z_sim/D");
109  m_tree->Branch("x_sim", &x_sim, "x_sim/D");
110  m_tree->Branch("dt_flight_sim", &dt_flight_sim, "dt_flight_sim/D");
111  }
112  if (m_calExpectedDriftTime) { // expected drift time, calculated form xfit
113  m_tree->Branch("t_fit", &t_fit, "t_fit/D");
114  }
115  if (!m_SmallerOutput) {
116  m_tree->Branch("tdc", &tdc, "tdc/I");
117  m_tree->Branch("z_prop", &z_prop, "z_prop/D");
118  m_tree->Branch("res_b", &res_b, "res_b/D");
119  m_tree->Branch("res_u", &res_u, "res_u/D");
120  m_tree->Branch("lr", &lr, "lr/I");
121  m_tree->Branch("trigHitPos_x", &trigHitPos_x, "trigHitPos_x/D");
122  m_tree->Branch("trigHitPos_z", &trigHitPos_z, "trigHitPos_z/D");
123  m_tree->Branch("numhits", &numhits, "numhits/I");
124  m_tree->Branch("res_b_err", &res_b_err, "res_b_err/D");
125  m_tree->Branch("res_u_err", &res_u_err, "res_u_err/D");
126  m_tree->Branch("absRes_u", &absRes_u, "absRes_u/D");
127  m_tree->Branch("absRes_b", &absRes_b, "absRes_b/D");
128  m_tree->Branch("dt_prop", &dt_prop, "dt_prop/D");
129  m_tree->Branch("dt_flight", &dt_flight, "dt_flight/D");
130  }
131 
132  // int N =m_Nchannel;//Number of Wire per Layer used;
133  TDirectory* oldDir = gDirectory;
134  TDirectory* histDir = oldDir->mkdir(m_histogramDirectoryName.c_str());
135  histDir->cd();
136  m_hNTracks = getHist("hNTracks", "number of tracks", 3, 0, 3);
137  m_hNTracks->GetXaxis()->SetBinLabel(1, "fitted, converged");
138  m_hNTracks->GetXaxis()->SetBinLabel(2, "fitted, not converged");
139  m_hNTracks->GetXaxis()->SetBinLabel(3, "TrackCand, but no Track");
140 
141  m_hNDF = getHist("hNDF", "NDF of fitted track;NDF;Tracks", 71, -1, 150);
142  m_hNHits = getHist("hNHits", "#hit of fitted track;#hit;Tracks", 61, -1, 150);
143  m_hNHits_trackcand = getHist("hNHits_trackcand", "#hit of track candidate;#hit;Tracks", 71, -1, 150);
144  m_hNTracksPerEvent = getHist("hNTracksPerEvent", "#tracks/Event;#Tracks;Event", 20, 0, 20);
145  m_hNTracksPerEventFitted = getHist("hNTracksPerEventFitted", "#tracks/Event After Fit;#Tracks;Event", 20, 0, 20);
146  m_hChi2 = getHist("hChi2", "#chi^{2} of tracks;#chi^{2};Tracks", 400, 0, 400);
147  m_hPhi0 = getHist("hPhi0", "#Phi_{0} of tracks;#phi_{0} (Degree);Tracks", 400, -190, 190);
148  m_hAlpha = getHist("hAlpha", "#alpha Dist.;#alpha (Degree);Hits", 360, -90, 90);
149  m_hTheta = getHist("hTheta", "#theta Dist.;#theta (Degree);Hits", 360, 0, 180);
150  m_hPval = getHist("hPval", "p-values of tracks;pVal;Tracks", 1000, 0, 1);
151  m_hEvtT0 = getHist("hEvtT0", "Event T0; EvtT0 (ns); #event", 200, -100, 100);
152 
153  m_hTriggerHitZX = getHist("TriggerHitZX", "Hit Position on trigger counter;z(cm);x(cm)", 300, -100, 100, 120, -15, 15);
154  if (m_MakeHitDist) {
155  m_h2DHitDistInCDCHit = getHist("2DHitDistInCDCHit", " CDCHit;WireID;LayerID",
156  m_up[8] - m_low[0], m_low[0], m_up[8], 56, 0, 56);
157  m_h2DHitDistInTrCand = getHist("2DHitDistInTrCand", "Track Cand ;WireID;LayerID",
158  m_up[8] - m_low[0], m_low[0], m_up[8], 56, 0, 56);
159  m_h2DHitDistInTrack = getHist("2DHitDistInTrack", "Fitted Track ;WireID;LayerID",
160  m_up[8] - m_low[0], m_low[0], m_up[8], 56, 0, 56);
161  }
162  if (m_fillExpertHistos) {
163  m_hNDFChi2 = getHist("hNDFChi2", "#chi^{2} of tracks;NDF;#chi^{2};Tracks", 8, 0, 8, 800, 0, 200);
164  m_hNDFPval = getHist("hNDFPval", "p-values of tracks;NDF;pVal;Tracks", 8, 0, 8, 100, 0, 1);
165  }
166  int sl;
167  for (int i = 0; i < 56; ++i) {
168  if (m_hitEfficiency) {
169  m_hHitEff_soft[i] = getHistProfile(Form("hHitEff_soft_L%d", i),
170  Form("hit efficiency(soft) of Layer %d ;Drift distance;Software Efficiency", i), 200, -1, 1);
171  }
172  if (m_MakeHitDist) {
173  if (i < 8) {sl = 0;} else { sl = floor((i - 8) / 6) + 1;}
174  m_hHitDistInCDCHit[i] = getHist(Form("hHitDistInCDCHit_layer%d", i), Form("Hit Dist. ICLayer_%d;WireID;#Hits", i),
175  m_up.at(sl) - m_low.at(sl), m_low.at(sl), m_up.at(sl));
176  m_hHitDistInCDCHit[i]->SetLineColor(kGreen);
177  m_hHitDistInTrCand[i] = getHist(Form("hHitDistInTrCand_layer%d", i), Form("Hit Dist. ICLayer_%d;WireID;#Hits", i),
178  m_up.at(sl) - m_low.at(sl), m_low.at(sl), m_up.at(sl));
179  m_hHitDistInTrCand[i]->SetLineColor(kRed);
180  m_hHitDistInTrack[i] = getHist(Form("hHitDistInTrack_layer%d", i), Form("Hit Dist. ICLayer_%d;WireID;#Hits", i),
181  m_up.at(sl) - m_low.at(sl), m_low.at(sl), m_up.at(sl));
182  }
183  const double normResRange = 20;
184  const double residualRange = 0.3;
185  std::string title, name;
186  if (m_plotResidual) {
187  name = (boost::format("hist_ResidualsU%1%") % i).str();
188  title = (boost::format("unnormalized, unbiased residuals in layer %1%;cm;Tracks") % i).str();
189  m_hResidualU[i] = getHist(name, title, 500, -residualRange, residualRange);
190 
191  name = (boost::format("hNormalizedResidualsU%1%") % i).str();
192  title = (boost::format("normalized, unbiased residuals in layer %1%;NDF;#sigma (cm);Tracks") % i).str();
193  m_hNormalizedResidualU[i] = getHist(name, title, 500, -normResRange, normResRange);
194 
195  name = (boost::format("DxDt%1%") % i).str();
196  title = (boost::format("Drift Length vs Drift time at Layer_%1%;Drift Length (cm);Drift time (ns)") % i).str();
197  m_hDxDt[i] = getHist(name, title, 200, -1, 1, 450, -50, 400);
198  }
199  if (m_fillExpertHistos) {
200  name = (boost::format("hNDFResidualsU%1%") % i).str();
201  title = (boost::format("unnormalized, unbiased residuals along U in layer %1%;NDF;cm;Tracks") % i).str();
202  m_hNDFResidualU[i] = getHist(name, title, 8, 0, 8, 1000, -residualRange, residualRange);
203 
204  name = (boost::format("hNDFNormalizedResidualsU%1%") % i).str();
205  title = (boost::format("normalized, unbiased residuals in layer %1%;NDF;#sigma (cm);Tracks") % i).str();
206  m_hNDFNormalizedResidualU[i] = getHist(name, title, 8, 0, 8, 1000, -normResRange, normResRange);
207  }
208  }
209  oldDir->cd();
210 }
211 
212 void CDCCRTestModule::initialize()
213 {
214  REG_HISTOGRAM
215  m_Tracks.isRequired(m_trackArrayName);
216  m_RecoTracks.isRequired(m_recoTrackArrayName);
217  m_TrackFitResults.isRequired(m_trackFitResultArrayName);
218  m_CDCHits.isRequired(m_cdcHitArrayName);
219  RelationArray relRecoTrackTrack(m_RecoTracks, m_Tracks, m_relRecoTrackTrackName);
220  //Store names to speed up creation later
221  m_relRecoTrackTrackName = relRecoTrackTrack.getName();
222 
223  for (size_t i = 0; i < m_allHistos.size(); ++i) {
224  m_allHistos[i]->Reset();
225  }
226  B2ASSERT("Trigger Position (TriggerPos) must be 3 components.", m_TriggerPos.size() == 3);
227  B2ASSERT("Normal vector of Trigger Plane (NormTriggerPlaneDirection) must be 3 components.", m_TriggerPlaneDirection.size() == 3);
228  B2ASSERT("Trigger size (TriggerSize) must be 2 component width and length(z direction)", m_TriggerSize.size() == 2);
229  B2ASSERT("List of Lower boundary (IWireLow) ( for histo must be 9 components, equivalent 9 supper layers", m_low.size() == 9);
230  B2ASSERT("List of Upper boundary (IWireUp) for histo must be 9 components, equivalent 9 supper layers", m_low.size() == 9);
231  B2INFO("Trigger Position (" << m_TriggerPos.at(0) << " ," << m_TriggerPos.at(1) << " ," << m_TriggerPos.at(2) << ")");
232 }
233 
234 void CDCCRTestModule::beginRun()
235 {
236 }
237 
238 void CDCCRTestModule::event()
239 {
240  evtT0 = 0.;
241  const RelationArray relTrackTrack(m_RecoTracks, m_Tracks, m_relRecoTrackTrackName);
242 
243  /* CDCHit distribution */
244  if (m_MakeHitDist) {
245  for (int i = 0; i < m_CDCHits.getEntries(); ++i) {
246  Belle2::CDCHit* hit = m_CDCHits[i];
247  m_hHitDistInCDCHit[getICLayer(hit->getISuperLayer(), hit->getILayer())]->Fill(hit->getIWire());
248  m_h2DHitDistInCDCHit->Fill(hit->getIWire(), getICLayer(hit->getISuperLayer(), hit->getILayer()));
249  }
250  }
251  // Loop over Recotracks
252  int nTr = m_RecoTracks.getEntries();
253  m_hNTracksPerEvent->Fill(nTr);
254 
255  int nfitted = 0;
256 
257  for (int i = 0; i < nTr; ++i) {
258  RecoTrack* track = m_RecoTracks[i];
259  if (track->getDirtyFlag()) {B2INFO("Dirty flag was set for track: " << track->getPositionSeed().Y()); continue;}
260  m_hNHits_trackcand->Fill(track->getNumberOfCDCHits());
261  if (m_MakeHitDist) {
262  getHitDistInTrackCand(track);
263  }
264  if (!track->hasTrackFitStatus()) {
265  m_hNTracks->Fill("Track not fitted", 1.0);
266  continue;
267  }
268  const genfit::FitStatus* fs = track->getTrackFitStatus();
269  if (!fs || !fs->isFitted()) {
270  m_hNTracks->Fill("Track not fitted", 1.0);
271  continue;
272  }
273  if (!fs->isFitConverged()) {//not fully convergence
274  m_hNTracks->Fill("fitted, not converged", 1.0);
275  B2DEBUG(99, "------Fitted but not converged");
276  continue;
277  }
278 
279  m_hNTracks->Fill("fitted, converged", 1.0);
280  B2DEBUG(99, "-------Fittted and Conveged");
281 
282  nfitted = nfitted + 1;
284  const Belle2::Track* b2track = track->getRelatedFrom<Belle2::Track>();
285  if (!b2track) {B2DEBUG(99, "No relation found"); continue;}
286  fitresult = b2track->getTrackFitResult(Const::muon);
287 
288  if (!fitresult) {
289  B2WARNING("track was fitted but Relation not found");
290  continue;
291  }
292 
293  if (m_noBFit) {ndf = fs->getNdf() + 1;} // incase no Magnetic field, NDF=4;
294  else {ndf = fs->getNdf();}
295  double Chi2 = fs->getChi2();
296  TrPval = std::max(0., ROOT::Math::chisquared_cdf_c(Chi2, ndf));
297  m_hPval->Fill(TrPval);
298  m_hNDF->Fill(ndf);
299  if (ndf < 15) continue;
300  if (m_EventT0Extraction) {
301  // event with is fail to extract t0 will be exclude from analysis
302  if (m_eventTimeStoreObject.isValid() && m_eventTimeStoreObject->hasEventT0()) {
303  evtT0 = m_eventTimeStoreObject->getEventT0();
304  m_hEvtT0->Fill(evtT0);
305  } else { continue;}
306  }
307 
308  d0 = fitresult->getD0();
309  z0 = fitresult->getZ0();
310  tanL = fitresult->getTanLambda();
311  omega = fitresult->getOmega();
312  phi0 = fitresult->getPhi0() * 180 / M_PI;
313  Pt = fitresult->getMomentum().Perp();
314  m_hPhi0->Fill(phi0);
315  m_hChi2->Fill(Chi2);
316  if (Pt < m_MinimumPt) continue;
317  if (m_hitEfficiency && track->getNumberOfCDCHits() > 30 && TrPval > 0.001) {
318  HitEfficiency(track);
319  }
320  if (m_fillExpertHistos) {
321  m_hNDFChi2->Fill(ndf, fs->getChi2());
322  m_hNDFPval->Fill(ndf, TrPval);
323  }
324  try {
325  plotResults(track);
326  } catch (const genfit::Exception& e) {
327  // at least log that there was something going very wrong
328  B2ERROR("Exception when calling the plotResults method" << e.what());
329  }
330 
331  if (m_EstimateResultForUnFittedLayer) {
332  //try {
333  // DONT IGNORE THE EXCEPTION BEING THROWN HERE
334  getResidualOfUnFittedLayer(track);
335  // plotResults(track);
336  //} catch (...) {
337  //B2ERROR (" fatal 2! ");
338 
339  //}
340  }
341  }
342 
343  m_hNTracksPerEventFitted->Fill(nfitted);
344 }
345 
346 void CDCCRTestModule::endRun()
347 {
348 }
349 
350 void CDCCRTestModule::terminate()
351 {
352 }
353 
354 void CDCCRTestModule::plotResults(Belle2::RecoTrack* track)
355 {
356  m_trigHitPos = getTriggerHitPosition(track);
357  trigHitPos_x = m_trigHitPos.X();
358  trigHitPos_z = m_trigHitPos.Z();
359  m_hTriggerHitZX->Fill(m_trigHitPos.Z(), m_trigHitPos.X());
360  bool hittrig = (sqrt((m_trigHitPos.X() - m_TriggerPos[0]) * (m_trigHitPos.X() - m_TriggerPos[0]) +
361  (m_trigHitPos.Y() - m_TriggerPos[1]) * (m_trigHitPos.Y() - m_TriggerPos[1])) < m_TriggerSize[0] / 2
362  && fabs(m_trigHitPos.Z() - m_TriggerPos[2]) < m_TriggerSize[1] / 2) ? true : false;
363  if (hittrig) {trighit = 1;}
364  else {trighit = 0;}
365  static CDCGeometryPar& cdcgeo = CDCGeometryPar::Instance();
367  m_hNHits->Fill(track->getNumberOfCDCHits());
368 
369  std::vector<genfit::TrackPoint*> tps = track->getHitPointsWithMeasurement();
370  numhits = tps.size();
371  BOOST_FOREACH(genfit::TrackPoint * tp, tps) {
372  if (!tp->hasRawMeasurements())
373  continue;
374 
375  const genfit::AbsMeasurement* raw = tp->getRawMeasurement(0);
376  const CDCRecoHit* rawCDC = dynamic_cast<const CDCRecoHit*>(raw);
377  if (rawCDC) {
378  WireID wireid = rawCDC->getWireID();
380  if (!kfi) {B2DEBUG(199, "No Fitter Info: Layer " << wireid.getICLayer()); continue;}
381 
382  for (unsigned int iMeas = 0; iMeas < kfi->getNumMeasurements(); ++iMeas) {
383  if ((kfi->getWeights().at(iMeas)) > 0.5) {
384  const genfit::MeasurementOnPlane& residual_b = kfi->getResidual(iMeas, true);
385  const genfit::MeasurementOnPlane& residual_u = kfi->getResidual(iMeas, false);
386  lay = wireid.getICLayer();
387  IWire = wireid.getIWire();
388  if (m_MakeHitDist) {
389  m_hHitDistInTrack[lay]->Fill(IWire);
390  m_h2DHitDistInTrack->Fill(IWire, lay);
391  }
392  boardID = cdcgeo.getBoardID(wireid);
393  Pval = TrPval;
394  tdc = rawCDC->getCDCHit()->getTDCCount();
395  adc = rawCDC->getCDCHit()->getADCCount();
396 
397  const genfit::MeasuredStateOnPlane& mop = kfi->getFittedState();
398  const TVector3 pocaOnWire = mop.getPlane()->getO();//Local wire position
399  const TVector3 pocaOnTrack = mop.getPlane()->getU();//residual direction
400  const TVector3 pocaMom = mop.getMom();
401  alpha = cdcgeo.getAlpha(pocaOnWire, pocaMom) ;
402  theta = cdcgeo.getTheta(pocaMom);
403  //Convert to outgoing
404  x_mea = kfi->getMeasurementOnPlane(iMeas)->getState()(0);
405  x_b = kfi->getFittedState(true).getState()(3);// x fit biased
406  res_b = residual_b.getState()(0);
407  x_u = kfi->getFittedState(false).getState()(3);//x fit unbiased
408  res_u = residual_u.getState()(0);
409  if (x_u > 0) lr = 1;
410  else lr = 0;
411  if (fabs(alpha) > M_PI / 2) {
412  x_b *= -1;
413  res_b *= -1;
414  x_u *= -1;
415  res_u *= -1;
416  }
417  x_mea = copysign(x_mea, x_u);
418  lr = cdcgeo.getOutgoingLR(lr, alpha);
419  theta = cdcgeo.getOutgoingTheta(alpha, theta);
420  alpha = cdcgeo.getOutgoingAlpha(alpha);
421  B2DEBUG(199, "x_unbiased " << x_u << " |left_right " << lr);
422  if (m_calExpectedDriftTime) { t_fit = cdcgeo.getDriftTime(std::abs(x_u), lay, lr, alpha , theta);}
423  alpha *= 180 / M_PI;
424  theta *= 180 / M_PI;
425  m_hAlpha->Fill(alpha);
426  m_hTheta->Fill(theta);
427 
428  //B2INFO("resi V " <<residual.getState()(1));
429  // weight_res = residual.getWeight();
430  absRes_b = std::abs(x_b + res_b) - std::abs(x_b);
431  absRes_u = std::abs(x_u + res_u) - std::abs(x_u);
432  weight = residual_u.getWeight();
433  res_b_err = std::sqrt(residual_b.getCov()(0, 0));
434  res_u_err = std::sqrt(residual_u.getCov()(0, 0));
435 
436  t = tdcTrans->getDriftTime(tdc, wireid, mop.getTime(), pocaOnWire.Z(), adc);
437  z = pocaOnWire.Z();
438 
439  // t = getCorrectedDriftTime(wireid, tdc, adc, z, z0);
440  if (m_StoreCDCSimHitInfo) {
441  CDCSimHit* simhit = rawCDC->getCDCHit()->getRelated<Belle2::CDCSimHit>();
442  if (simhit) {
443  x_sim = simhit->getDriftLength();
444  z_sim = simhit->getPosWire().Z();
445  dt_flight_sim = simhit->getFlightTime();
446  }
447  }
448  m_tree->Fill();
449  if (m_plotResidual) {
450  m_hDxDt[lay]->Fill(x_u, t);
451  m_hResidualU[lay]->Fill(res_b, weight);
452  m_hNormalizedResidualU[lay]->Fill(res_b / sqrt(residual_b.getCov()(0, 0)), weight);
453  }
454  if (m_fillExpertHistos) {
455  m_hNDFResidualU[lay]->Fill(ndf, res_b, weight);
456  m_hNDFResidualU[lay]->Fill(ndf, res_b);
457  m_hNDFNormalizedResidualU[lay]->Fill(ndf, res_b / std::sqrt(residual_b.getCov()(0, 0)), weight);
458  }
459  } //NDF
460  // }//end of if isU
461  }//end of for
462  }//end of rawCDC
463  }//end of for tp
464 }//end of func
465 
466 void CDCCRTestModule::getHitDistInTrackCand(const RecoTrack* track)
467 {
468  BOOST_FOREACH(const RecoHitInformation::UsedCDCHit * cdchit, track->getCDCHitList()) {
469  int iclay = getICLayer(cdchit->getISuperLayer(), cdchit->getILayer());
470  B2DEBUG(99, "In TrackCand: ICLayer: " << iclay << "IWire: " << cdchit->getIWire());
471  m_hHitDistInTrCand[iclay]->Fill(cdchit->getIWire());
472  m_h2DHitDistInTrCand->Fill(cdchit->getIWire(), iclay);
473  }
474 }
475 
476 TVector3 CDCCRTestModule::getTriggerHitPosition(RecoTrack* track)
477 {
478  TVector3 trigpos(m_TriggerPos.at(0), m_TriggerPos.at(1), m_TriggerPos.at(2));
479  TVector3 trigDir(m_TriggerPlaneDirection.at(0), m_TriggerPlaneDirection.at(1), m_TriggerPlaneDirection.at(2));
480  const genfit::AbsTrackRep* trackRepresentation = track->getCardinalRepresentation();
481  genfit::MeasuredStateOnPlane mop = track->getMeasuredStateOnPlaneClosestTo(trigpos, trackRepresentation);
482  TVector3 pos(-200, 200, 200);
483  try {
484  double l = mop.extrapolateToPlane(genfit::SharedPlanePtr(new genfit::DetPlane(trigpos, trigDir)));
485  if (fabs(l) < 1000) pos = mop.getPos();
486  } catch (const genfit::Exception& er) {
487  B2WARNING("extrapolate to Trigger counter failure" << er.what());
488  }
489  return pos;
490 }
491 void CDCCRTestModule::HitEfficiency(const Belle2::RecoTrack* track)
492 {
493  /* static CDCGeometryPar& cdcgeo = CDCGeometryPar::Instance();
494  for (int i = 0; i < 56; ++i) {
495  double rcell = (rinnerlayer[i] + routerlayer[i]) / 2;
496  double arcL = h.getArcLength2DAtCylindricalR(rcell);
497  const TVector3 hitpos = h.getPositionAtArcLength2D(arcL);
498  int cellID = cdcgeo.cellId(i, hitpos);
499  B2INFO("Hit at LayerID - CellID: " << i << "-" << cellID);
500  }
501  */
503  BOOST_FOREACH(const RecoHitInformation::UsedCDCHit * cdchit, track->getCDCHitList()) {
504  WireID Wid = WireID(cdchit->getID());
505  const genfit::TrackPoint* tp = track->getCreatedTrackPoint(track->getRecoHitInformation(cdchit));
506  //some hit didn't take account in fitting, so I use left/right info frm track finding results.
507  int RLInfo = 0;
508  RecoHitInformation::RightLeftInformation rightLeftHitInformation = track->getRecoHitInformation(cdchit)->getRightLeftInformation();
509  if (rightLeftHitInformation == RecoHitInformation::RightLeftInformation::c_left) {
510  RLInfo = -1;
511  } else if (rightLeftHitInformation == RecoHitInformation::RightLeftInformation::c_right) {
512  RLInfo = 1;
513  } else continue;
514 
515  if (!tp->hasRawMeasurements())
516  continue;
518  if (!kfi) continue;
519 
520  // double max = std::max_element(kfi->getWeights(),kfi->getNumMeasurements());
521  double max = 0.;
522  unsigned short imea = 0;
523  for (unsigned int iMeas = 0; iMeas < kfi->getNumMeasurements(); ++iMeas) {
524  double ww = kfi->getWeights().at(iMeas);
525  if (ww > max) {max = ww; imea = iMeas;}
526  }
527  double xx = kfi->getMeasurementOnPlane(imea)->getState()(0);
528  m_hHitEff_soft[Wid.getICLayer()]->Fill(std::copysign(xx, RLInfo), max);
529  }
531 }
532 
533 void CDCCRTestModule::getResidualOfUnFittedLayer(Belle2::RecoTrack* track)
534 {
535  B2INFO("Start estimate residual for un-fitted layer");
536  B2INFO("position seed" << track->getPositionSeed().Y());
537  static CDCGeometryPar& cdcgeo = CDCGeometryPar::Instance();
539  const genfit::AbsTrackRep* trackRepresentation = track->getCardinalRepresentation();
540  x_b = 0; res_b = 0; res_b_err = 0;
541  x_u = 0; res_u = 0; res_u_err = 0;
542  z = 0; dt_flight_sim = 0; z_prop = 0; t = 0;
543  dt_prop = 0; dt_flight = 0; alpha = 0; theta = 0;
544  tdc = 0; adc = 0; lay = 0; IWire = 0;
545  // ndf =0; Pval=0; numhits=0; trigHitPos_x= 0; trigHitPos_z=0;
546  // trighit=1; lr=-1;
547  B2INFO("number of cdchit" << track->getCDCHitList().size());
548  B2INFO("number of point use int fit" << ndf + 4);
549 
550  typedef std::pair<double, const RecoHitInformation*> SortingRecoHitPair;
551 
552  for (const RecoHitInformation::UsedCDCHit* cdchit : track->getCDCHitList()) {
553  // RecoHitInformation* recoHitInfo = track->getRecoHitInformation(cdchit);
554  if (track->getRecoHitInformation(cdchit)->useInFit()) continue;
555  // yeah is true, but better to check for the above
556  //if ((recoHitInfo->getCreatedTrackPoint())) continue;
557  // This was wrong: the sorting parameter is not the hitID
558  int hitSortingParameter = track->getRecoHitInformation(cdchit)->getSortingParameter();
559 
560  SortingRecoHitPair frontSideHit = std::make_pair(0, nullptr);;
561  SortingRecoHitPair backsideSideHit = std::make_pair(0, nullptr);;
562  SortingRecoHitPair hit4extraction; // = std::make_pair(0, nullptr); avoid cppcheck warning.
563 
564  //find closest hit to hit which do not fit
565  // if (hitID < track->getNumberOfCDCHits() / 2) { //case for first part of track, searching forward, stop at first choice
566  for (const RecoHitInformation::UsedCDCHit* hit : track->getCDCHitList()) {
567  RecoHitInformation const* recoHitInfo_fw = track->getRecoHitInformation(hit);
568  if (recoHitInfo_fw->useInFit()) { //may be should check fit status of that hit, do it later.
569  frontSideHit = std::make_pair(recoHitInfo_fw->getSortingParameter(), recoHitInfo_fw);
570  break;
571  }
572  }
573  //}
574  // if (hitID > track->getNumberOfCDCHits() / 2) { //case for last part of track, searching backward, and stop at the first choice
575  auto hitListReverse = track->getCDCHitList();
576  std::reverse(hitListReverse.begin() , hitListReverse.end());
577  for (const RecoHitInformation::UsedCDCHit* hit : hitListReverse) {
578  RecoHitInformation const* recoHitInfo_bkw = track->getRecoHitInformation(hit);
579  if (recoHitInfo_bkw->useInFit()) {
580  // also get proper id here
581  backsideSideHit = std::make_pair(recoHitInfo_bkw->getSortingParameter(), recoHitInfo_bkw);
582  break;
583  }
584  }
585  B2DEBUG(99, "forward sorting parameter: " << frontSideHit.first << " |backward sorting parameter = " << backsideSideHit.first);
586  if (std::fabs(frontSideHit.first - hitSortingParameter) < std::fabs(backsideSideHit.first - hitSortingParameter)) {
587  hit4extraction = frontSideHit;
588  } else {
589  hit4extraction = backsideSideHit;
590  }
591 
592  // no proper neighbouring hit found
593  if (hit4extraction.second == nullptr)
594  continue;
595 
596  auto closestHitTrackPoint = track->getCreatedTrackPoint(hit4extraction.second);
597  // now we need to find the hit behind this sorting param !
598  // but easy: we have already the TrackPoint via the RecoHitInformation
599  genfit::MeasuredStateOnPlane meaOnPlane = closestHitTrackPoint->getFitterInfo(trackRepresentation)->getFittedState(
600  true /* biased version */);
601 
602  //start to extrapolation
603  WireID wireid = WireID(cdchit->getID());
604  // double flightTime1 = meaOnPlane.getTime();
605  //Now reconstruct plane for hit
606  genfit::SharedPlanePtr plane = nullptr;
607  try {
608  plane = constructPlane(meaOnPlane, wireid);
609  } catch (const genfit::Exception& e) {
610  B2WARNING("Error happen, can not reconstruct plan for extrapolating" << e.what());
611  continue;
612  }
613  double segmentLength;
614  try {
615  segmentLength = meaOnPlane.extrapolateToPlane(plane);
616  } catch (const genfit::Exception& e) {
617  B2WARNING("Could not extrapolate the fit" << e.what());
618  continue;
619  }
620  IWire = wireid.getIWire();
621  lay = wireid.getICLayer();
622  const TVector3 pocaOnWire = meaOnPlane.getPlane()->getO();//Local wire position
623  const TVector3 pocaMom = meaOnPlane.getMom();
624  x_u = meaOnPlane.getState()(3);
625  alpha = cdcgeo.getAlpha(pocaOnWire, pocaMom) ;
626  theta = cdcgeo.getTheta(pocaMom);
627  z = pocaOnWire.Z();
628  z_prop = z - cdcgeo.wireBackwardPosition(wireid, CDCGeometryPar::c_Aligned).Z();
629  dt_prop = z_prop / 27.25;
630  //Convert to outgoing
631  if (x_u > 0) lr = 1;
632  else lr = 0;
633  if (fabs(alpha) > M_PI / 2) {
634  x_u *= -1;
635  }
636  lr = cdcgeo.getOutgoingLR(lr, alpha);
637  theta = cdcgeo.getOutgoingTheta(alpha, theta);
638  alpha = cdcgeo.getOutgoingAlpha(alpha);
639  dt_flight = meaOnPlane.getTime();
640  x_mea = tdcTrans->getDriftLength(cdchit->getTDCCount(), wireid, dt_flight, lr, pocaOnWire.Z(), alpha, theta, cdchit->getADCCount());
641  x_mea = std::copysign(x_mea, x_u);
642  res_u = x_mea - x_u;
643  absRes_u = fabs(x_mea) - fabs(x_u);
644  alpha *= 180 / M_PI;
645  theta *= 180 / M_PI;
646  m_hAlpha->Fill(alpha);
647  m_hTheta->Fill(theta);
648  if (m_StoreCDCSimHitInfo) {
649  CDCSimHit* simhit = cdchit->getRelated<Belle2::CDCSimHit>();
650  if (simhit) {
651  x_sim = simhit->getDriftLength();
652  z_sim = simhit->getPosWire().Z();
653  dt_flight_sim = simhit->getFlightTime();
654  }
655  }
656  B2DEBUG(199, "we calculate residua for lay - IWire: " << lay << " - " << IWire);
657  B2DEBUG(199, "distance between two hit" << segmentLength);
658  B2DEBUG(199, "Flight Time (extra | sim)" << dt_flight << " - " << dt_flight_sim);
659  B2DEBUG(199, "DriftLength (cal | sim)" << x_mea << " - " << x_sim);
660  m_tree->Fill();
661  }
662 }
663 
664 const genfit::SharedPlanePtr CDCCRTestModule::constructPlane(const genfit::MeasuredStateOnPlane& state, WireID m_wireID)
665 {
666  // We reconstruct plane from measuedStateOnPlane from one fitted hit.
667  // because I don't want to change state of this plane so I create other state to extrapolate.
668  // first: extrapolate to wire (ideal geometry, nosag) to find z pos than get virtual wire pos due to sag,
669  // extrapolate again to found z pos and U.
671  const StateOnPlane stateOnPlane = StateOnPlane(state.getState(), state.getPlane(), state.getRep());
672  genfit::StateOnPlane st(stateOnPlane);
673 
674  const TVector3& Wire1PosIdeal(cdcgeoTrans->getWireBackwardPosition(m_wireID));
675  const TVector3& Wire2PosIdeal(cdcgeoTrans->getWireForwardPosition(m_wireID));
676 
677  // unit vector of wire direction
678  TVector3 WireDirectionIdeal = Wire2PosIdeal - Wire1PosIdeal;
679  WireDirectionIdeal.SetMag(1.);//normalized
680 
681  // extraplate to find z
682  const genfit::AbsTrackRep* rep = state.getRep();
683  rep->extrapolateToLine(st, Wire1PosIdeal, WireDirectionIdeal);
684  const TVector3& PocaIdeal = rep->getPos(st);
685 
686  double zPOCA = (Wire1PosIdeal.Z()
687  + WireDirectionIdeal.Dot(PocaIdeal - Wire1PosIdeal) * WireDirectionIdeal.Z());
688 
689  // Now re-extrapolate to new wire direction, wire sag was taking account.
690  const TVector3& wire1(cdcgeoTrans->getWireBackwardPosition(m_wireID, zPOCA));
691  const TVector3& wire2(cdcgeoTrans->getWireForwardPosition(m_wireID, zPOCA));
692 
693  // unit vector of wire direction (include sag)
694  TVector3 wireDirection = wire2 - wire1;
695  wireDirection.SetMag(1.);
696 
697  // extraplate to find poca
698  rep->extrapolateToLine(st, wire1, wireDirection);
699  const TVector3& poca = rep->getPos(st);
700  TVector3 dirInPoca = rep->getMom(st);
701  dirInPoca.SetMag(1.);
702  const TVector3& pocaOnWire = wire1 + wireDirection.Dot(poca - wire1) * wireDirection;
703  if (fabs(wireDirection.Angle(dirInPoca)) < 0.01) {
704  B2WARNING("cannot construct det plane, track parallel with wire");
705  }
706  // construct orthogonal (unit) vector for plane
707  const TVector3& U = wireDirection.Cross(dirInPoca);
708  genfit::SharedPlanePtr pl = genfit::SharedPlanePtr(new genfit::DetPlane(pocaOnWire, U, wireDirection));
709  return pl;
710 }
711 
Class containing the result of the unpacker in raw data and the result of the digitizer in simulation...
Definition: CDCHit.h:40
unsigned short getIWire() const
Getter for iWire.
Definition: CDCHit.h:166
short getTDCCount() const
Getter for TDC count.
Definition: CDCHit.h:219
unsigned short getID() const
Getter for encoded wire number.
Definition: CDCHit.h:193
unsigned short getADCCount() const
Getter for integrated charge.
Definition: CDCHit.h:230
unsigned short getISuperLayer() const
Getter for iSuperLayer.
Definition: CDCHit.h:184
unsigned short getILayer() const
Getter for iLayer.
Definition: CDCHit.h:172
This class is used to transfer CDC information to the track fit.
Definition: CDCRecoHit.h:32
const CDCHit * getCDCHit() const
get the pointer to the CDCHit object that was used to create this CDCRecoHit object.
Definition: CDCRecoHit.h:112
WireID getWireID() const
Getter for WireID object.
Definition: CDCRecoHit.h:49
Example Detector.
Definition: CDCSimHit.h:22
TVector3 getPosWire() const
The method to get position on wire.
Definition: CDCSimHit.h:205
double getFlightTime() const
The method to get flight time.
Definition: CDCSimHit.h:190
double getDriftLength() const
The method to get drift length.
Definition: CDCSimHit.h:187
CDC Cosmic test calibration module.
The Class for CDC Geometry Parameters.
unsigned short getBoardID(const WireID &wID) const
Returns frontend board id. corresponding to the wire id.
const TVector3 wireBackwardPosition(int layerId, int cellId, EWirePosition set=c_Base) const
Returns the backward position of the input sense wire.
double getAlpha(const TVector3 &posOnWire, const TVector3 &momentum) const
Returns track incident angle in rphi plane (alpha in rad.).
double getOutgoingAlpha(const double alpha) const
Converts incoming- to outgoing-alpha.
unsigned short getOutgoingLR(const unsigned short lr, const double alpha) const
Converts incoming-lr to outgoing-lr.
double getOutgoingTheta(const double alpha, const double theta) const
Converts incoming- to outgoing-theta.
double getTheta(const TVector3 &momentum) const
Returns track incident angle (theta in rad.).
double getDriftTime(double dist, unsigned short layer, unsigned short lr, double alpha, double theta) const
Return the drift time to the sense wire.
This class uses the realistic detector geometry (the one after alignment procedure) for the translati...
const TVector3 getWireBackwardPosition(const WireID &wireID) override
Get wire position at backward end.
const TVector3 getWireForwardPosition(const WireID &wireID) override
Get wire position at forward end.
Translator mirroring the realistic Digitization.
double getDriftLength(unsigned short tdcCount, const WireID &wireID=WireID(), double timeOfFlightEstimator=0, bool leftRight=false, double z=0, double alpha=0, double theta=static_cast< double >(TMath::Pi()/2.), unsigned short adcCount=0) override
Get Drift length.
double getDriftTime(unsigned short tdcCount, const WireID &wireID, double timeOfFlightEstimator, double z, unsigned short adcCount) override
Get Drift time.
HistoModule.h is supposed to be used instead of Module.h for the modules with histogram definitions t...
Definition: HistoModule.h:29
This class stores additional information to every CDC/SVD/PXD hit stored in a RecoTrack.
RightLeftInformation
The RightLeft information of the hit which is only valid for CDC hits.
bool useInFit() const
Get the flag, whether this his should be used in a fit or not.
unsigned int getSortingParameter() const
Get the sorting parameter.
This is the Reconstruction Event-Data Model Track.
Definition: RecoTrack.h:76
Low-level class to create/modify relations between StoreArrays.
Definition: RelationArray.h:62
T * getRelated(const std::string &name="", const std::string &namedRelation="") const
Get the object to or from which this object has a relation.
const std::string & getName() const
Return name under which the object is saved in the DataStore.
Class that bundles various TrackFitResults.
Definition: Track.h:25
const TrackFitResult * getTrackFitResult(const Const::ChargedStable &chargedStable) const
Access to TrackFitResults.
Definition: Track.cc:17
Class to identify a wire inside the CDC.
Definition: WireID.h:34
unsigned short getICLayer() const
Getter for continuous layer numbering.
Definition: WireID.cc:24
Contains the measurement and covariance in raw detector coordinates.
Abstract base class for a track representation.
Definition: AbsTrackRep.h:66
virtual double extrapolateToLine(StateOnPlane &state, const TVector3 &linePoint, const TVector3 &lineDirection, bool stopAtBoundary=false, bool calcJacobianNoise=false) const =0
Extrapolates the state to the POCA to a line, and returns the extrapolation length and,...
virtual TVector3 getMom(const StateOnPlane &state) const =0
Get the cartesian momentum vector of a state.
virtual TVector3 getPos(const StateOnPlane &state) const =0
Get the cartesian position of a state.
Detector plane.
Definition: DetPlane.h:59
Exception class for error handling in GENFIT (provides storage for diagnostic information)
Definition: Exception.h:48
virtual const char * what() const noexcept
Standard error message handling for exceptions. use like "std::cerr << e.what();".
Definition: Exception.cc:52
Class where important numbers and properties of a fit can be stored.
Definition: FitStatus.h:80
bool isFitConverged(bool inAllPoints=true) const
Did the fit converge (in all Points or only partially)?
Definition: FitStatus.h:105
double getChi2() const
Get chi^2 of the fit.
Definition: FitStatus.h:120
bool isFitted() const
Has the track been fitted?
Definition: FitStatus.h:94
double getNdf() const
Get the degrees of freedom of the fit.
Definition: FitStatus.h:122
Collects information needed and produced by a AbsKalmanFitter implementations and is specific to one ...
MeasurementOnPlane getResidual(unsigned int iMeasurement=0, bool biased=false, bool onlyMeasurementErrors=true) const override
Get unbiased (default) or biased residual from ith measurement.
const MeasuredStateOnPlane & getFittedState(bool biased=true) const override
Get unbiased or biased (default) smoothed state.
std::vector< double > getWeights() const
Get weights of measurements.
#StateOnPlane with additional covariance matrix.
Measured coordinates on a plane.
A state with arbitrary dimension defined in a DetPlane.
Definition: StateOnPlane.h:47
Object containing AbsMeasurement and AbsFitterInfo objects.
Definition: TrackPoint.h:46
KalmanFitterInfo * getKalmanFitterInfo(const AbsTrackRep *rep=nullptr) const
Helper to avoid casting.
Definition: TrackPoint.cc:180
#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.
Defines for I/O streams used for error and debug printing.
std::shared_ptr< genfit::DetPlane > SharedPlanePtr
Shared Pointer to a DetPlane.