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#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
30using namespace std;
31using namespace Belle2;
32using namespace CDC;
33using namespace genfit;
34//-----------------------------------------------------------------
35// Register the Module
36//-----------------------------------------------------------------
37REG_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 transverse momentum smaller than this will not be recorded", 0.);
71}
72
74{
75 m_allHistos.clear();
76}
77//------------------------------------------------------------------
78// Function to define histograms
79//-----------------------------------------------------------------
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");
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 }
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
213{
214 REG_HISTOGRAM
215 m_Tracks.isRequired(m_trackArrayName);
218 m_CDCHits.isRequired(m_cdcHitArrayName);
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
235{
236}
237
239{
240 evtT0 = 0.;
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) {
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, "-------Fitted and Converged");
281
282 nfitted = nfitted + 1;
284 const Belle2::Track* b2track = track->getRelatedFrom<Belle2::Track>();
285 if (!b2track) {B2DEBUG(99, "No relation found"); continue;}
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;} // in case 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;
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();
312 phi0 = fitresult->getPhi0() * 180 / M_PI;
313 Pt = fitresult->getMomentum().Rho();
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
332 //try {
333 // DONT IGNORE THE EXCEPTION BEING THROWN HERE
335 // plotResults(track);
336 //} catch (...) {
337 //B2ERROR (" fatal 2! ");
338
339 //}
340 }
341 }
342
343 m_hNTracksPerEventFitted->Fill(nfitted);
344}
345
347{
348}
349
351{
352}
353
355{
360 bool hittrig = (sqrt((m_trigHitPos.X() - m_TriggerPos[0]) * (m_trigHitPos.X() - m_TriggerPos[0]) +
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();
379 const genfit::KalmanFitterInfo* kfi = tp->getKalmanFitterInfo();
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) {
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 B2Vector3D pocaOnWire = mop.getPlane()->getO();//Local wire position
399 const B2Vector3D pocaMom = mop.getMom();
400 alpha = cdcgeo.getAlpha(pocaOnWire, pocaMom) ;
401 theta = cdcgeo.getTheta(pocaMom);
402 //Convert to outgoing
403 x_mea = kfi->getMeasurementOnPlane(iMeas)->getState()(0);
404 x_b = kfi->getFittedState(true).getState()(3);// x fit biased
405 res_b = residual_b.getState()(0);
406 x_u = kfi->getFittedState(false).getState()(3);//x fit unbiased
407 res_u = residual_u.getState()(0);
408 if (x_u > 0) lr = 1;
409 else lr = 0;
410 if (fabs(alpha) > M_PI / 2) {
411 x_b *= -1;
412 res_b *= -1;
413 x_u *= -1;
414 res_u *= -1;
415 }
416 x_mea = copysign(x_mea, x_u);
417 lr = cdcgeo.getOutgoingLR(lr, alpha);
419 alpha = cdcgeo.getOutgoingAlpha(alpha);
420 B2DEBUG(199, "x_unbiased " << x_u << " |left_right " << lr);
421 if (m_calExpectedDriftTime) { t_fit = cdcgeo.getDriftTime(std::abs(x_u), lay, lr, alpha, theta);}
422 alpha *= 180 / M_PI;
423 theta *= 180 / M_PI;
424 m_hAlpha->Fill(alpha);
425 m_hTheta->Fill(theta);
426
427 //B2INFO("resi V " <<residual.getState()(1));
428 // weight_res = residual.getWeight();
429 absRes_b = std::abs(x_b + res_b) - std::abs(x_b);
430 absRes_u = std::abs(x_u + res_u) - std::abs(x_u);
431 weight = residual_u.getWeight();
432 res_b_err = std::sqrt(residual_b.getCov()(0, 0));
433 res_u_err = std::sqrt(residual_u.getCov()(0, 0));
434
435 t = tdcTrans->getDriftTime(tdc, wireid, mop.getTime(), pocaOnWire.Z(), adc);
436 z = pocaOnWire.Z();
437
438 // t = getCorrectedDriftTime(wireid, tdc, adc, z, z0);
440 CDCSimHit* simhit = rawCDC->getCDCHit()->getRelated<Belle2::CDCSimHit>();
441 if (simhit) {
442 x_sim = simhit->getDriftLength();
443 z_sim = simhit->getPosWire().Z();
444 dt_flight_sim = simhit->getFlightTime();
445 }
446 }
447 m_tree->Fill();
448 if (m_plotResidual) {
449 m_hDxDt[lay]->Fill(x_u, t);
450 m_hResidualU[lay]->Fill(res_b, weight);
451 m_hNormalizedResidualU[lay]->Fill(res_b / sqrt(residual_b.getCov()(0, 0)), weight);
452 }
453 if (m_fillExpertHistos) {
455 m_hNDFResidualU[lay]->Fill(ndf, res_b);
456 m_hNDFNormalizedResidualU[lay]->Fill(ndf, res_b / std::sqrt(residual_b.getCov()(0, 0)), weight);
457 }
458 } //NDF
459 // }//end of if isU
460 }//end of for
461 }//end of rawCDC
462 }//end of for tp
463}//end of func
464
466{
467 BOOST_FOREACH(const RecoHitInformation::UsedCDCHit * cdchit, track->getCDCHitList()) {
468 int iclay = getICLayer(cdchit->getISuperLayer(), cdchit->getILayer());
469 B2DEBUG(99, "In TrackCand: ICLayer: " << iclay << "IWire: " << cdchit->getIWire());
470 m_hHitDistInTrCand[iclay]->Fill(cdchit->getIWire());
471 m_h2DHitDistInTrCand->Fill(cdchit->getIWire(), iclay);
472 }
473}
474
476{
477 B2Vector3D trigpos(m_TriggerPos.at(0), m_TriggerPos.at(1), m_TriggerPos.at(2));
479 const genfit::AbsTrackRep* trackRepresentation = track->getCardinalRepresentation();
480 B2Vector3D pos(-200, 200, 200);
481 try {
482 genfit::MeasuredStateOnPlane mop = track->getMeasuredStateOnPlaneClosestTo(ROOT::Math::XYZVector(trigpos), trackRepresentation);
483 double l = mop.extrapolateToPlane(genfit::SharedPlanePtr(new genfit::DetPlane(trigpos, trigDir)));
484 if (fabs(l) < 1000) pos = mop.getPos();
485 } catch (const genfit::Exception& er) {
486 B2WARNING("extrapolate to Trigger counter failure" << er.what());
487 } catch (const std::runtime_error& er) {
488 B2WARNING("Runtime error encountered: " << er.what());
489 } catch (...) {
490 B2WARNING("Undefined exception encountered.");
491 }
492 return pos;
493}
495{
496 /* static CDCGeometryPar& cdcgeo = CDCGeometryPar::Instance();
497 for (int i = 0; i < 56; ++i) {
498 double rcell = (rinnerlayer[i] + routerlayer[i]) / 2;
499 double arcL = h.getArcLength2DAtCylindricalR(rcell);
500 const B2Vector3D hitpos = h.getPositionAtArcLength2D(arcL);
501 int cellID = cdcgeo.cellId(i, hitpos);
502 B2INFO("Hit at LayerID - CellID: " << i << "-" << cellID);
503 }
504 */
506 BOOST_FOREACH(const RecoHitInformation::UsedCDCHit * cdchit, track->getCDCHitList()) {
507 WireID Wid = WireID(cdchit->getID());
508 const genfit::TrackPoint* tp = track->getCreatedTrackPoint(track->getRecoHitInformation(cdchit));
509 //some hit didn't take account in fitting, so I use left/right info from track finding results.
510 int RLInfo = 0;
511 RecoHitInformation::RightLeftInformation rightLeftHitInformation = track->getRecoHitInformation(cdchit)->getRightLeftInformation();
512 if (rightLeftHitInformation == RecoHitInformation::RightLeftInformation::c_left) {
513 RLInfo = -1;
514 } else if (rightLeftHitInformation == RecoHitInformation::RightLeftInformation::c_right) {
515 RLInfo = 1;
516 } else continue;
517
518 if (!tp->hasRawMeasurements())
519 continue;
520 const genfit::KalmanFitterInfo* kfi = tp->getKalmanFitterInfo();
521 if (!kfi) continue;
522
523 // double max = std::max_element(kfi->getWeights(),kfi->getNumMeasurements());
524 double max = 0.;
525 unsigned short imea = 0;
526 for (unsigned int iMeas = 0; iMeas < kfi->getNumMeasurements(); ++iMeas) {
527 double ww = kfi->getWeights().at(iMeas);
528 if (ww > max) {max = ww; imea = iMeas;}
529 }
530 double xx = kfi->getMeasurementOnPlane(imea)->getState()(0);
531 m_hHitEff_soft[Wid.getICLayer()]->Fill(std::copysign(xx, RLInfo), max);
532 }
534}
535
537{
538 B2INFO("Start estimate residual for un-fitted layer");
539 B2INFO("position seed" << track->getPositionSeed().Y());
540 static CDCGeometryPar& cdcgeo = CDCGeometryPar::Instance();
542 const genfit::AbsTrackRep* trackRepresentation = track->getCardinalRepresentation();
543 x_b = 0; res_b = 0; res_b_err = 0;
544 x_u = 0; res_u = 0; res_u_err = 0;
545 z = 0; dt_flight_sim = 0; z_prop = 0; t = 0;
546 dt_prop = 0; dt_flight = 0; alpha = 0; theta = 0;
547 tdc = 0; adc = 0; lay = 0; IWire = 0;
548 // ndf =0; Pval=0; numhits=0; trigHitPos_x= 0; trigHitPos_z=0;
549 // trighit=1; lr=-1;
550 B2INFO("number of cdchit" << track->getCDCHitList().size());
551 B2INFO("number of point use int fit" << ndf + 4);
552
553 typedef std::pair<double, const RecoHitInformation*> SortingRecoHitPair;
554
555 for (const RecoHitInformation::UsedCDCHit* cdchit : track->getCDCHitList()) {
556 // RecoHitInformation* recoHitInfo = track->getRecoHitInformation(cdchit);
557 if (track->getRecoHitInformation(cdchit)->useInFit()) continue;
558 // yeah is true, but better to check for the above
559 //if ((recoHitInfo->getCreatedTrackPoint())) continue;
560 // This was wrong: the sorting parameter is not the hitID
561 int hitSortingParameter = track->getRecoHitInformation(cdchit)->getSortingParameter();
562
563 SortingRecoHitPair frontSideHit = std::make_pair(0, nullptr);;
564 SortingRecoHitPair backsideSideHit = std::make_pair(0, nullptr);;
565 SortingRecoHitPair hit4extraction; // = std::make_pair(0, nullptr); avoid cppcheck warning.
566
567 //find closest hit to hit which do not fit
568 // if (hitID < track->getNumberOfCDCHits() / 2) { //case for first part of track, searching forward, stop at first choice
569 for (const RecoHitInformation::UsedCDCHit* hit : track->getCDCHitList()) {
570 RecoHitInformation const* recoHitInfo_fw = track->getRecoHitInformation(hit);
571 if (recoHitInfo_fw->useInFit()) { //may be should check fit status of that hit, do it later.
572 frontSideHit = std::make_pair(recoHitInfo_fw->getSortingParameter(), recoHitInfo_fw);
573 break;
574 }
575 }
576 //}
577 // if (hitID > track->getNumberOfCDCHits() / 2) { //case for last part of track, searching backward, and stop at the first choice
578 auto hitListReverse = track->getCDCHitList();
579 std::reverse(hitListReverse.begin(), hitListReverse.end());
580 for (const RecoHitInformation::UsedCDCHit* hit : hitListReverse) {
581 RecoHitInformation const* recoHitInfo_bkw = track->getRecoHitInformation(hit);
582 if (recoHitInfo_bkw->useInFit()) {
583 // also get proper id here
584 backsideSideHit = std::make_pair(recoHitInfo_bkw->getSortingParameter(), recoHitInfo_bkw);
585 break;
586 }
587 }
588 B2DEBUG(99, "forward sorting parameter: " << frontSideHit.first << " |backward sorting parameter = " << backsideSideHit.first);
589 if (std::fabs(frontSideHit.first - hitSortingParameter) < std::fabs(backsideSideHit.first - hitSortingParameter)) {
590 hit4extraction = frontSideHit;
591 } else {
592 hit4extraction = backsideSideHit;
593 }
594
595 // no proper neighbouring hit found
596 if (hit4extraction.second == nullptr)
597 continue;
598
599 auto closestHitTrackPoint = track->getCreatedTrackPoint(hit4extraction.second);
600 // now we need to find the hit behind this sorting param !
601 // but easy: we have already the TrackPoint via the RecoHitInformation
602 genfit::MeasuredStateOnPlane meaOnPlane = closestHitTrackPoint->getFitterInfo(trackRepresentation)->getFittedState(
603 true /* biased version */);
604
605 //start to extrapolation
606 WireID wireid = WireID(cdchit->getID());
607 // double flightTime1 = meaOnPlane.getTime();
608 //Now reconstruct plane for hit
609 genfit::SharedPlanePtr plane = nullptr;
610 try {
611 plane = constructPlane(meaOnPlane, wireid);
612 } catch (const genfit::Exception& e) {
613 B2WARNING("Error happen, can not reconstruct plan for extrapolating" << e.what());
614 continue;
615 }
616 double segmentLength;
617 try {
618 segmentLength = meaOnPlane.extrapolateToPlane(plane);
619 } catch (const genfit::Exception& e) {
620 B2WARNING("Could not extrapolate the fit" << e.what());
621 continue;
622 }
623 IWire = wireid.getIWire();
624 lay = wireid.getICLayer();
625 const B2Vector3D pocaOnWire = meaOnPlane.getPlane()->getO();//Local wire position
626 const B2Vector3D pocaMom = meaOnPlane.getMom();
627 x_u = meaOnPlane.getState()(3);
628 alpha = cdcgeo.getAlpha(pocaOnWire, pocaMom) ;
629 theta = cdcgeo.getTheta(pocaMom);
630 z = pocaOnWire.Z();
631 z_prop = z - cdcgeo.wireBackwardPosition(wireid, CDCGeometryPar::c_Aligned).Z();
632 dt_prop = z_prop / 27.25;
633 //Convert to outgoing
634 if (x_u > 0) lr = 1;
635 else lr = 0;
636 if (fabs(alpha) > M_PI / 2) {
637 x_u *= -1;
638 }
639 lr = cdcgeo.getOutgoingLR(lr, alpha);
641 alpha = cdcgeo.getOutgoingAlpha(alpha);
642 dt_flight = meaOnPlane.getTime();
643 x_mea = tdcTrans->getDriftLength(cdchit->getTDCCount(), wireid, dt_flight, lr, pocaOnWire.Z(), alpha, theta, cdchit->getADCCount());
644 x_mea = std::copysign(x_mea, x_u);
645 res_u = x_mea - x_u;
646 absRes_u = fabs(x_mea) - fabs(x_u);
647 alpha *= 180 / M_PI;
648 theta *= 180 / M_PI;
649 m_hAlpha->Fill(alpha);
650 m_hTheta->Fill(theta);
652 CDCSimHit* simhit = cdchit->getRelated<Belle2::CDCSimHit>();
653 if (simhit) {
654 x_sim = simhit->getDriftLength();
655 z_sim = simhit->getPosWire().Z();
656 dt_flight_sim = simhit->getFlightTime();
657 }
658 }
659 B2DEBUG(199, "we calculate residua for lay - IWire: " << lay << " - " << IWire);
660 B2DEBUG(199, "distance between two hit" << segmentLength);
661 B2DEBUG(199, "Flight Time (extra | sim)" << dt_flight << " - " << dt_flight_sim);
662 B2DEBUG(199, "DriftLength (cal | sim)" << x_mea << " - " << x_sim);
663 m_tree->Fill();
664 }
665}
666
667const genfit::SharedPlanePtr CDCCRTestModule::constructPlane(const genfit::MeasuredStateOnPlane& state, WireID m_wireID)
668{
669 // We reconstruct plane from measuedStateOnPlane from one fitted hit.
670 // because I don't want to change state of this plane so I create other state to extrapolate.
671 // first: extrapolate to wire (ideal geometry, nosag) to find z pos than get virtual wire pos due to sag,
672 // extrapolate again to found z pos and U.
674 const StateOnPlane stateOnPlane = StateOnPlane(state.getState(), state.getPlane(), state.getRep());
675 genfit::StateOnPlane st(stateOnPlane);
676
677 const B2Vector3D& Wire1PosIdeal(cdcgeoTrans->getWireBackwardPosition(m_wireID));
678 const B2Vector3D& Wire2PosIdeal(cdcgeoTrans->getWireForwardPosition(m_wireID));
679
680 // unit vector of wire direction
681 B2Vector3D WireDirectionIdeal = Wire2PosIdeal - Wire1PosIdeal;
682 WireDirectionIdeal.SetMag(1.);//normalized
683
684 // extrapolate to find z
685 const genfit::AbsTrackRep* rep = state.getRep();
686 rep->extrapolateToLine(st, Wire1PosIdeal, WireDirectionIdeal);
687 const B2Vector3D& PocaIdeal = rep->getPos(st);
688
689 double zPOCA = (Wire1PosIdeal.Z()
690 + WireDirectionIdeal.Dot(PocaIdeal - Wire1PosIdeal) * WireDirectionIdeal.Z());
691
692 // Now re-extrapolate to new wire direction, wire sag was taking account.
693 const B2Vector3D& wire1(cdcgeoTrans->getWireBackwardPosition(m_wireID, zPOCA));
694 const B2Vector3D& wire2(cdcgeoTrans->getWireForwardPosition(m_wireID, zPOCA));
695
696 // unit vector of wire direction (include sag)
697 B2Vector3D wireDirection = wire2 - wire1;
698 wireDirection.SetMag(1.);
699
700 // extrapolate to find poca
701 rep->extrapolateToLine(st, wire1, wireDirection);
702 const B2Vector3D& poca = rep->getPos(st);
703 B2Vector3D dirInPoca = rep->getMom(st);
704 dirInPoca.SetMag(1.);
705 const B2Vector3D& pocaOnWire = wire1 + wireDirection.Dot(poca - wire1) * wireDirection;
706 if (fabs(wireDirection.Angle(dirInPoca)) < 0.01) {
707 B2WARNING("cannot construct det plane, track parallel with wire");
708 }
709 // construct orthogonal (unit) vector for plane
710 const B2Vector3D& U = wireDirection.Cross(dirInPoca);
711 genfit::SharedPlanePtr pl = genfit::SharedPlanePtr(new genfit::DetPlane(pocaOnWire, U, wireDirection));
712 return pl;
713}
714
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
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
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:184
B2Vector3D getPosWire() const
The method to get position on wire.
Definition: CDCSimHit.h:199
double getDriftLength() const
The method to get drift length.
Definition: CDCSimHit.h:181
StoreObjPtr< EventT0 > m_eventTimeStoreObject
Event timing.
std::string m_recoTrackArrayName
Belle2::RecoTrack StoreArray nam.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 efficience 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.
bool isRequired(const std::string &name="")
Ensure this array/object has been registered previously.
int getEntries() const
Get the number of objects in the array.
Definition: StoreArray.h:216
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:560
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28
Abstract base class for different kinds of events.
STL namespace.