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