ATLAS Offline Software
Loading...
Searching...
No Matches
DetailedIDNtupleTool.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3 */
4
6
8#include "TrkTrack/Track.h"
12
13//for Amg::error helper function:
15
18
22
27
29
30#include "TFile.h"
31#include "TTree.h"
32#include "TMath.h"
33#include <cmath>
34
35
36namespace InDet {
37//________________________________________________________________________
38 DetailedIDNtupleTool::DetailedIDNtupleTool(const std::string& type, const std::string& name, const IInterface* parent)
39 : AthAlgTool(type, name, parent)
40 , m_trackSumTool("Trk::TrackSummaryTool/TrackSummaryTool", this)
44 , m_storeTruth{true}
46 , m_file{}
47 , m_tree{}
48 , m_trackCollection{"Tracks"}
49 , m_tracksTruthName{"TrackTruthCollection"} {
50 declareInterface<Trk::IFillNtupleTool>(this);
51
52 declareProperty("FileName", m_filename = "IDAlign.root");
53 declareProperty("FilePath", m_filepath = "./");
54 declareProperty("TrackCollection", m_trackCollection);
55 declareProperty("TracksTruthName", m_tracksTruthName);
56 declareProperty("TrackSummaryTool", m_trackSumTool);
57 declareProperty("AlignModuleTool", m_alignModuleTool);
58 declareProperty("MatchProbability", m_matchProbability);
59 declareProperty("StoreTruth", m_storeTruth);
60 declareProperty("MatchProbability", m_matchProbability);
61 declareProperty("StoreNormalRefittedOnly", m_storeNormalRefittedOnly);
62 declareProperty("StoreConstrainedOnly", m_storeConstrainedOnly);
63 }
64
65//________________________________________________________________________
68
69//________________________________________________________________________
71 ATH_MSG_DEBUG("initialize() of DetailedIDNtupleTool");
72 ATH_MSG_DEBUG("creating file with name " << m_filepath << m_filename);
73
75 ATH_CHECK(m_truthToTrack.retrieve());
76
77 /* Retrieve extrapolator tool */
78 ATH_CHECK(m_extrapolator.retrieve());
79
80 // get AlignModuleTool
81 if (m_alignModuleTool.empty()) {
82 ATH_MSG_FATAL("m_alignModuleTool is empty in DetailedIDNtupleTool::initialize()");
83 return StatusCode::FAILURE;
84 } else {
85 ATH_CHECK(m_alignModuleTool.retrieve());
86 }
87 return StatusCode::SUCCESS;
88 }
89
90//________________________________________________________________________
92 if (m_file && m_file->IsOpen()) {
93 m_file->cd();
94 if (m_tree) m_tree->Write();
95 }
96 return StatusCode::SUCCESS;
97 }
98
99//________________________________________________________________________
101 ATH_MSG_DEBUG("finalize() of DetailedIDNtupleTool");
102 return StatusCode::SUCCESS;
103 }
104
105//________________________________________________________________________
106 void DetailedIDNtupleTool::dumpTrack(int /* itrk*/, const Trk::AlignTrack* alignTrack) {
107 ATH_MSG_DEBUG("In dumpTrack()");
109 const EventContext& ctx = Gaudi::Hive::currentContext();
110 // if hope to dump all the tracks, should set:
111 // m_storeNormalRefittedOnly = false, m_storeConstrainedOnly=false
114 return;
115 }
116
117 // get run and event numbers
118 ATH_MSG_DEBUG("Retrieving event info.");
119 const xAOD::EventInfo* eventInfo;
120 if (evtStore()->retrieve(eventInfo).isFailure()) {
121 ATH_MSG_ERROR("Could not retrieve event info.");
122 } else {
123 m_runNumber = eventInfo->runNumber();
124 m_evtNumber = eventInfo->eventNumber();
125 }
126
127 // initialize variables
128 constexpr double invalidParameterValue{-999.};
129 m_d0 = invalidParameterValue;
130 m_z0 = invalidParameterValue;
131 m_phi0 = invalidParameterValue;
132 m_theta = invalidParameterValue;
133 m_qoverp = invalidParameterValue;
134 m_pt = invalidParameterValue;
135 m_eta = invalidParameterValue;
136 constexpr double invalidChiSq{-1e12};
137 m_chi2 = invalidChiSq;
138 m_chi2prob = invalidChiSq;
139 constexpr int invalidDegreesOfFreedom{-999};
140
141 m_ndof = invalidDegreesOfFreedom;
142 m_xvtx = invalidParameterValue;
143 m_yvtx = invalidParameterValue;
144 m_zvtx = invalidParameterValue;
145
146 m_err_d0 = invalidParameterValue;
147 m_err_z0 = invalidParameterValue;
148 m_err_phi0 = invalidParameterValue;
149 m_err_theta = invalidParameterValue;
150 m_err_qoverp = invalidParameterValue;
151 m_toRef_d0 = invalidParameterValue;
152 m_toRef_z0 = invalidParameterValue;
153 m_toRef_phi0 = invalidParameterValue;
154 m_toRef_theta = invalidParameterValue;
155 m_toRef_qoverp = invalidParameterValue;
156
157
158 m_original_d0 = invalidParameterValue;
159 m_original_z0 = invalidParameterValue;
160 m_original_phi0 = invalidParameterValue;
161 m_original_theta = invalidParameterValue;
162 m_original_qoverp = invalidParameterValue;
163 m_original_pt =invalidParameterValue;
164 m_original_eta = -999.;
165 m_original_chi2 = invalidChiSq;
166 m_original_chi2prob = invalidChiSq;
167 m_original_ndof = invalidDegreesOfFreedom;
168 m_original_xvtx = invalidParameterValue;
169 m_original_yvtx = invalidParameterValue;
170 m_original_zvtx = invalidParameterValue;
171
172 m_original_err_d0 = invalidParameterValue;
173 m_original_err_z0 = invalidParameterValue;
174 m_original_err_phi0 = invalidParameterValue;
175 m_original_err_theta = invalidParameterValue;
176 m_original_err_qoverp = invalidParameterValue;;
177 m_original_toRef_d0 = invalidParameterValue;
178 m_original_toRef_z0 = invalidParameterValue;
179 m_original_toRef_phi0 = invalidParameterValue;
180 m_original_toRef_theta =invalidParameterValue;
181 m_original_toRef_qoverp = invalidParameterValue;;
182
183 m_truth_d0 = invalidParameterValue;
184 m_truth_z0 = invalidParameterValue;
185 m_truth_phi0 = invalidParameterValue;
186 m_truth_theta = invalidParameterValue;
187 m_truth_qoverp = invalidParameterValue;
188 m_truth_pt = invalidParameterValue;
189 m_truth_eta = invalidParameterValue;;
190 m_truth_prod_x = invalidParameterValue;
191 m_truth_prod_y = invalidParameterValue;
192 m_truth_prod_z = invalidParameterValue;;
193
194
195
196 const Trk::Perigee* aMeasPer = (alignTrack->perigeeParameters());
197
198 if (not aMeasPer) {
199 ATH_MSG_ERROR("Could not get Trk::Perigee of the alignTrack");
200 return;
201 } else {
202 m_d0 = aMeasPer->parameters()[Trk::d0];
203 m_z0 = aMeasPer->parameters()[Trk::z0];
204 m_phi0 = aMeasPer->parameters()[Trk::phi0];
205 m_theta = aMeasPer->parameters()[Trk::theta];
206 m_qoverp = aMeasPer->parameters()[Trk::qOverP];
207
208 const AmgSymMatrix(5) * locCov = aMeasPer->covariance();
209 m_err_d0 = Amg::error(*locCov, Trk::d0);
210 m_err_z0 = Amg::error(*locCov, Trk::z0);
211 m_err_phi0 = Amg::error(*locCov, Trk::phi0);
214
215 m_pt =
216 std::sqrt((aMeasPer->momentum().x()) * (aMeasPer->momentum().x()) + (aMeasPer->momentum().y()) *
217 (aMeasPer->momentum().y()));
218 m_eta = aMeasPer->eta();
219
220 m_xvtx = aMeasPer->position().x();
221 m_yvtx = aMeasPer->position().y();
222 m_zvtx = aMeasPer->position().z();
223
224
225
226 // get fit quality and chi2 probability of track
227 const Trk::FitQuality* fitQual = alignTrack->fitQuality();
228 if (not fitQual) ATH_MSG_ERROR("No fit quality assigned to the track");
229 else {
230 if (fitQual->chiSquared() > 0. && fitQual->numberDoF() > 0) {
231 m_chi2 = fitQual->chiSquared();
232 m_ndof = fitQual->numberDoF();
233 m_chi2prob = TMath::Prob(m_chi2, m_ndof);
234
235 ATH_MSG_DEBUG(" - chi2 : " << m_chi2);
236 ATH_MSG_DEBUG(" - ndof : " << m_ndof);
237 ATH_MSG_DEBUG(" - chi2/ndof : " << m_chi2 / (double) m_ndof);
238 ATH_MSG_DEBUG(" - chi2 propability : " << m_chi2prob);
239 }
240 }
241 }
242
243 const Trk::MeasurementBase* measbase = *(alignTrack->measurementsOnTrack()->begin());
244 Amg::Vector3D refPoint = measbase->associatedSurface().center();
245 const Trk::VertexOnTrack* vot = dynamic_cast<const Trk::VertexOnTrack*> (measbase);
246 if (!vot) {
247 ATH_MSG_ERROR(" Seems the pseudo-measuremnt in the alignTrack not exist!");
248 ATH_MSG_ERROR(" this pseudo-measurement has been rejected as outlier in the refitting!");
249 return;
250 } else {
251 ATH_MSG_DEBUG(" VoT get from the alignTrack:" << *vot);
252 }
253
254 ATH_MSG_DEBUG(" the pseudo-measurement position: " << refPoint);
255 const Trk::Track* originalTrack = alignTrack->originalTrack();
256
257 const Trk::Perigee* originalMeasPer = originalTrack->perigeeParameters();
258 if (!originalMeasPer) {
259 ATH_MSG_ERROR("No original track!");
260 } else {
261 m_original_d0 = originalMeasPer->parameters()[Trk::d0];
262 m_original_z0 = originalMeasPer->parameters()[Trk::z0];
263 m_original_phi0 = originalMeasPer->parameters()[Trk::phi0];
264 m_original_theta = originalMeasPer->parameters()[Trk::theta];
265 m_original_qoverp = originalMeasPer->parameters()[Trk::qOverP];
266
267 const AmgSymMatrix(5) * locError = originalMeasPer->covariance();
273
275 m_original_eta = -std::log(std::tan(0.5 * m_original_theta));
276
277 m_original_xvtx = originalMeasPer->position().x();
278 m_original_yvtx = originalMeasPer->position().y();
279 m_original_zvtx = originalMeasPer->position().z();
280
281 // get fit quality and chi2 probability of original Track
282 const Trk::FitQuality* originalFitQual = originalTrack->fitQuality();
283 if (not originalFitQual) ATH_MSG_ERROR("No fit quality assigned to the track");
284 else {
285 if (originalFitQual->chiSquared() > 0. && originalFitQual->numberDoF() > 0) {
286 m_original_chi2 = originalFitQual->chiSquared();
287 m_original_ndof = originalFitQual->numberDoF();
289
290 ATH_MSG_DEBUG(" - original chi2 : " << m_original_chi2);
291 ATH_MSG_DEBUG(" - original ndof : " << m_original_ndof);
292 ATH_MSG_DEBUG(" - original chi2/ndof : " << m_original_chi2 / (double) m_original_ndof);
293 ATH_MSG_DEBUG(" - original chi2 propability : " << m_original_chi2prob);
294 }
295 }
296
297
298 const Trk::PerigeeSurface persf(refPoint);
299 const Trk::Perigee* originalPerigeeAtRef =
300 dynamic_cast<const Trk::Perigee*>(m_extrapolator->extrapolateTrack(ctx, *originalTrack, persf).release());
301 if (!originalPerigeeAtRef) {
302 const Trk::Perigee* originalTrackPerigee = originalTrack->perigeeParameters();
303 if (originalTrackPerigee && ((originalTrackPerigee->associatedSurface())) == persf) {
304 ATH_MSG_DEBUG("Perigee of Track is already expressed to given vertex, a copy is returned.");
305 originalPerigeeAtRef = originalTrackPerigee->clone();
306 } else ATH_MSG_DEBUG("Extrapolation to Perigee failed, NULL pointer is returned.");
307 }
308
309 if (originalPerigeeAtRef) {
310 std::unique_ptr<const Trk::Perigee > originalMeasPerAtRef(originalPerigeeAtRef);
311 m_original_toRef_d0 = originalMeasPerAtRef->parameters()[Trk::d0];
312 m_original_toRef_z0 = originalMeasPerAtRef->parameters()[Trk::z0];
313 m_original_toRef_phi0 = originalMeasPerAtRef->parameters()[Trk::phi0];
314 m_original_toRef_theta = originalMeasPerAtRef->parameters()[Trk::theta];
315 m_original_toRef_qoverp = originalMeasPerAtRef->parameters()[Trk::qOverP];
316 }
317
318 const Trk::Perigee* PerigeeAtRef =
319 dynamic_cast<const Trk::Perigee*>(m_extrapolator->extrapolateTrack(ctx, *alignTrack, persf).release());
320 if (!PerigeeAtRef) {
321 const Trk::Perigee* alignTrackPerigee = alignTrack->perigeeParameters();
322 if (alignTrackPerigee && ((alignTrackPerigee->associatedSurface())) == persf) {
323 ATH_MSG_DEBUG("Perigee of AlignTrack is already expressed to given vertex, a copy is returned.");
324 PerigeeAtRef = alignTrackPerigee->clone();
325 } else ATH_MSG_DEBUG("Extrapolation to Perigee failed, NULL pointer is returned.");
326 }
327
328 //post-eigen, can simply use the TrackParameters * returned by m_extrapolator->extrapolate?
329 if (PerigeeAtRef) {
330 std::unique_ptr<const Trk::Perigee > MeasPerAtRef((PerigeeAtRef));
331 m_toRef_d0 = MeasPerAtRef->parameters()[Trk::d0];
332 m_toRef_z0 = MeasPerAtRef->parameters()[Trk::z0];
333 m_toRef_phi0 = MeasPerAtRef->parameters()[Trk::phi0];
334 m_toRef_theta = MeasPerAtRef->parameters()[Trk::theta];
335 m_toRef_qoverp = MeasPerAtRef->parameters()[Trk::qOverP];
336 }
337 }
338
340 if (m_storeTruth) retrieveTruthInfo(alignTrack);
341
342
343 m_file->cd();
344 m_tree->Fill();
345 ATH_MSG_DEBUG("tree filled");
346
347 return;
348 }
349
351 // although we select tracks using the TrackSelectionTool, we still need to get a complete TrackCollection
352 // from StoreGate for use in the track-truth map, otherwise the track-truth matching is screwed up
353
354 const TrackCollection* RecCollection = nullptr;
355
356 if (evtStore()->retrieve(RecCollection, m_trackCollection).isFailure()) {
357 ATH_MSG_WARNING("Track collection \"" << m_trackCollection << "\" not found.");
358 return false;
359 }
360 if (RecCollection) {
362 "Retrieved " << m_trackCollection << " with size " << RecCollection->size() <<
363 " reconstructed tracks from storegate");
364 } else {
365 ATH_MSG_WARNING("RecCollection is null pointer in DetailedIDNtupleTool::retrieveTruthInfo");
366 return false;
367 }
368
369 // get TrackTruthCollection
370 const TrackTruthCollection* TruthMap = nullptr;
371 if (StatusCode::SUCCESS != evtStore()->retrieve(TruthMap, m_tracksTruthName)) {
372 ATH_MSG_DEBUG("Cannot find " << m_tracksTruthName);
373 return false;
374 }
375
377 "Track Truth Collection with name " << m_tracksTruthName << " with size " << TruthMap->size() <<
378 " found in StoreGate");
379
380 bool flag = false;
381 // get fit quality and chi2 probability of track
382 const Trk::Perigee* startPerigee = alignTrack->perigeeParameters();
383 const Trk::Perigee* measPer = startPerigee;
384
385 if (measPer == nullptr) {
386 ATH_MSG_DEBUG("No measured perigee parameters assigned to the track");
387 return false;
388 }
389
391 tracklink.setElement(const_cast<Trk::Track*>(alignTrack->originalTrack()));
392 tracklink.setStorableObject(*RecCollection);
393 const ElementLink<TrackCollection> tracklink2 = tracklink;
394
395 TrackTruthCollection::const_iterator found = TruthMap->find(tracklink2);
396 if ((found != TruthMap->end()) && (found->second.probability() > m_matchProbability)) {
397 TrackTruth trtruth = found->second;
398 const HepMcParticleLink& HMPL = trtruth.particleLink();
399
400 if (HMPL.isValid()) {
401#ifdef HEPMC3
402 HepMC::ConstGenParticlePtr genparptr = HMPL.scptr();
403#else
404 const HepMC::GenParticle* genparptr = HMPL.cptr();
405#endif
406
407 if (genparptr) {
408 if (genparptr->production_vertex()) {
409 if (genparptr->pdg_id() == 0) {
410 ATH_MSG_INFO("PDG ID is zero in DetailedIDNtupleTool::retrieveTruthInfo");
411 } else {
412 const Trk::TrackParameters* generatedTrackPerigee = m_truthToTrack->makePerigeeParameters(genparptr);
413 if (!generatedTrackPerigee) ATH_MSG_WARNING("Unable to extrapolate genparticle to perigee!");
414 else {
415 flag = true;
416 m_truth_qoverpt = generatedTrackPerigee->parameters()[Trk::qOverP] / std::sin(
417 generatedTrackPerigee->parameters()[Trk::theta]);
418 m_truth_qoverp = generatedTrackPerigee->parameters()[Trk::qOverP];
419 m_truth_phi0 = generatedTrackPerigee->parameters()[Trk::phi0];
420 m_truth_d0 = generatedTrackPerigee->parameters()[Trk::d0];
421 m_truth_z0 = generatedTrackPerigee->parameters()[Trk::z0];
422 m_truth_theta = generatedTrackPerigee->parameters()[Trk::theta];
423 m_truth_eta = generatedTrackPerigee->eta();
424 m_truth_prod_x = genparptr->production_vertex()->position().x();
425 m_truth_prod_y = genparptr->production_vertex()->position().y();
426 m_truth_prod_z = genparptr->production_vertex()->position().z();
427
428 delete generatedTrackPerigee;
429 m_truth_pt = 1. / std::abs(m_truth_qoverpt);
430 m_truth_charge = 1;
431 if (m_truth_qoverpt < 0) m_truth_charge = -1;
432 if (m_truth_phi0 < 0) m_truth_phi0 += 2 * M_PI;
433 ATH_MSG_DEBUG("Found matched truth track with phi, PT = " << m_truth_phi0 << ", " << m_truth_pt);
434 }
435 }
436 }
437 }
438 }
439 }
440
441 return flag;
442 }
443
444//I think this is never used!
446 const Trk::TrackParameters* aMeasPer {};
447
448 for (const auto i:*(alignTrack->trackStateOnSurfaces())) {
450 aMeasPer = (i->trackParameters());
451 break;
452 }
453 }
454 if (not aMeasPer) {
455 ATH_MSG_ERROR("Could not get Trk::MeasuredPerigee of the alignTrack");
456 }
457 return aMeasPer;
458 }
459
462
463//________________________________________________________________________
466
469
470//________________________________________________________________________
473
474//________________________________________________________________________
476 //m_file = new TFile((m_filepath+m_filename).c_str(), "RECREATE");
477 m_file->cd();
478 m_tree = new TTree("IDAlign", "Inner Detector Alignment Ntuple");
479
480 m_tree->Branch("run", &m_runNumber, "run/I");
481 m_tree->Branch("evt", &m_evtNumber, "evt/I");
482
483 m_tree->Branch("xvtx", &m_xvtx, "xvtx/D");
484 m_tree->Branch("yvtx", &m_yvtx, "yvtx/D");
485 m_tree->Branch("zvtx", &m_zvtx, "zvtx/D");
486 m_tree->Branch("d0", &m_d0, "d0/D");
487 m_tree->Branch("z0", &m_z0, "z0/D");
488 m_tree->Branch("phi0", &m_phi0, "phi0/D");
489 m_tree->Branch("theta", &m_theta, "theta/D");
490 m_tree->Branch("qoverp", &m_qoverp, "qoverp/D");
491 m_tree->Branch("pt", &m_pt, "pt/D");
492 m_tree->Branch("eta", &m_eta, "eta/D");
493 m_tree->Branch("chi2", &m_chi2, "chi2/D");
494 m_tree->Branch("ndof", &m_ndof, "ndof/I");
495 m_tree->Branch("chi2prob", &m_chi2prob, "chi2prob/D");
496 m_tree->Branch("err_d0", &m_err_d0, "err_d0/D");
497 m_tree->Branch("err_z0", &m_err_z0, "err_z0/D");
498 m_tree->Branch("err_phi0", &m_err_phi0, "err_phi0/D");
499 m_tree->Branch("err_theta", &m_err_theta, "err_theta/D");
500 m_tree->Branch("err_qoverp", &m_err_qoverp, "err_qoverp/D");
501
502 if (m_storeTruth) {
503 m_tree->Branch("m_truth_prod_x", &m_truth_prod_x, "truth_prod_x/D");
504 m_tree->Branch("m_truth_prod_y", &m_truth_prod_y, "truth_prod_y/D");
505 m_tree->Branch("m_truth_prod_z", &m_truth_prod_z, "truth_prod_z/D");
506 m_tree->Branch("truth_d0", &m_truth_d0, "truth_d0/D");
507 m_tree->Branch("truth_z0", &m_truth_z0, "truth_z0/D");
508 m_tree->Branch("truth_phi0", &m_truth_phi0, "truth_phi0/D");
509 m_tree->Branch("truth_theta", &m_truth_theta, "truth_theta/D");
510 m_tree->Branch("truth_qoverp", &m_truth_qoverp, "truth_qoverp/D");
511 m_tree->Branch("truth_pt", &m_truth_pt, "truth_pt/D");
512 m_tree->Branch("truth_eta", &m_truth_eta, "truth_eta/D");
513 }
514
515 m_tree->Branch("original_xvtx", &m_original_xvtx, "original_xvtx/D");
516 m_tree->Branch("original_yvtx", &m_original_yvtx, "original_yvtx/D");
517 m_tree->Branch("original_zvtx", &m_original_zvtx, "original_zvtx/D");
518 m_tree->Branch("original_d0", &m_original_d0, "original_d0/D");
519 m_tree->Branch("original_z0", &m_original_z0, "original_z0/D");
520 m_tree->Branch("original_phi0", &m_original_phi0, "original_phi0/D");
521 m_tree->Branch("original_theta", &m_original_theta, "original_theta/D");
522 m_tree->Branch("original_qoverp", &m_original_qoverp, "original_qoverp/D");
523 m_tree->Branch("original_pt", &m_original_pt, "original_pt/D");
524 m_tree->Branch("original_eta", &m_original_eta, "original_eta/D");
525 m_tree->Branch("original_chi2", &m_original_chi2, "original_chi2/D");
526 m_tree->Branch("original_ndof", &m_original_ndof, "original_ndof/I");
527 m_tree->Branch("original_chi2prob", &m_original_chi2prob, "original_chi2prob/D");
528
529 m_tree->Branch("original_err_d0", &m_original_err_d0, "original_err_d0/D");
530 m_tree->Branch("original_err_z0", &m_original_err_z0, "original_err_z0/D");
531 m_tree->Branch("original_err_phi0", &m_original_err_phi0, "original_err_phi0/D");
532 m_tree->Branch("original_err_theta", &m_original_err_theta, "original_err_theta/D");
533 m_tree->Branch("original_err_qoverp", &m_original_err_qoverp, "original_err_qoverp/D");
534
535 m_tree->Branch("toRef_d0", &m_toRef_d0, "toRef_d0/D");
536 m_tree->Branch("toRef_z0", &m_toRef_z0, "toRef_z0/D");
537 m_tree->Branch("toRef_phi0", &m_toRef_phi0, "toRef_phi0/D");
538 m_tree->Branch("toRef_theta", &m_toRef_theta, "toRef_theta/D");
539 m_tree->Branch("toRef_qoverp", &m_toRef_qoverp, "toRef_qoverp/D");
540 m_tree->Branch("original_toRef_d0", &m_original_toRef_d0, "original_toRef_d0/D");
541 m_tree->Branch("original_toRef_z0", &m_original_toRef_z0, "original_toRef_z0/D");
542 m_tree->Branch("original_toRef_phi0", &m_original_toRef_phi0, "original_toRef_phi0/D");
543 m_tree->Branch("original_toRef_theta", &m_original_toRef_theta, "original_toRef_theta/D");
544 m_tree->Branch("original_toRef_qoverp", &m_original_toRef_qoverp, "original_toRef_qoverp/D");
545
546 return;
547 }
548} // end namespace
#define M_PI
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_FATAL(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
#define AmgSymMatrix(dim)
DataVector< Trk::Track > TrackCollection
This typedef represents a collection of Trk::Track objects.
AthAlgTool(const std::string &type, const std::string &name, const IInterface *parent)
Constructor with parameters:
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)
ServiceHandle< StoreGateSvc > & evtStore()
size_type size() const noexcept
Returns the number of elements in the collection.
bool m_storeTruth
retrieve the truth information
virtual void showStatistics() override
write statistics out to log file
std::string m_filename
name of ntuple file
bool m_storeConstrainedOnly
onlystore the AlignTrack which is normally refitted
DetailedIDNtupleTool(const std::string &type, const std::string &name, const IInterface *parent)
virtual StatusCode fillNtuple() override
writes trees and histograms to ntuple
PublicToolHandle< Trk::IExtrapolator > m_extrapolator
track extrapolator
virtual void storeHitmap() override
stores hitmap for writing to ntuple
ToolHandle< Trk::ITrackSummaryTool > m_trackSumTool
Pointer to track summary tool.
bool m_storeNormalRefittedOnly
only store the AlignTrack which is normally refitted
virtual StatusCode finalize() override
virtual void dumpTrack(int itrk, const Trk::AlignTrack *alignTrack) override
fills track information to ntuple
virtual void fillSummary() override
fills ntuple with event and track summary information
std::string m_filepath
path to ntuple file
PublicToolHandle< Trk::ITruthToTrack > m_truthToTrack
the truth to track Tool
virtual void fillHitmap() override
fills ntuple with hit information
ToolHandle< Trk::IAlignModuleTool > m_alignModuleTool
Pointer to AlignmModuleTool.
const Trk::TrackParameters * perigeeParameter(const Trk::AlignTrack *track) const
virtual StatusCode initialize() override
bool retrieveTruthInfo(const Trk::AlignTrack *trk)
double m_matchProbability
the probabililty cut in the truth matching
MC particle associated with a reco track + the quality of match.
Definition TrackTruth.h:14
const HepMcParticleLink & particleLink() const
Definition TrackTruth.h:26
const Track * originalTrack() const
retrieve pointer to original track
Definition AlignTrack.h:90
@ NormalRefitted
normally refitted, without adding any pseudo-measurement
Definition AlignTrack.h:48
@ BeamspotConstrained
refitted with beamspot constraint
Definition AlignTrack.h:49
AlignTrackType type() const
get and set the refit type
Definition AlignTrack.h:96
Class to represent and store fit qualities from track reconstruction in terms of and number of degre...
Definition FitQuality.h:97
int numberDoF() const
returns the number of degrees of freedom of the overall track or vertex fit as integer
Definition FitQuality.h:60
double chiSquared() const
returns the of the overall track fit
Definition FitQuality.h:56
This class is the pure abstract base class for all fittable tracking measurements.
virtual const Surface & associatedSurface() const =0
Interface method to get the associated Surface.
double eta() const
Access method for pseudorapidity - from momentum.
const Amg::Vector3D & momentum() const
Access method for the momentum.
const Amg::Vector3D & position() const
Access method for the position.
virtual ParametersT< DIM, T, S > * clone() const override final
Virtual clone.
virtual const S & associatedSurface() const override final
Access to the Surface method.
Class describing the Line to which the Perigee refers to.
const Amg::Vector3D & center() const
Returns the center position of the Surface.
@ Perigee
This represents a perigee, and so will contain a Perigee object only.
const Trk::TrackStates * trackStateOnSurfaces() const
return a pointer to a const DataVector of const TrackStateOnSurfaces.
const DataVector< const MeasurementBase > * measurementsOnTrack() const
return a pointer to a vector of MeasurementBase (NOT including any that come from outliers).
const Perigee * perigeeParameters() const
return Perigee.
const FitQuality * fitQuality() const
return a pointer to the fit quality const-overload
Class to handle Vertex On Tracks, it inherits from the common MeasurementBase.
uint32_t runNumber() const
The current event's run number.
uint64_t eventNumber() const
The current event's event number.
static std::string release
Definition computils.h:50
double error(const Amg::MatrixX &mat, int index)
return diagonal error of the matrix caller should ensure the matrix is symmetric and the index is in ...
Eigen::Matrix< double, 3, 1 > Vector3D
const GenParticle * ConstGenParticlePtr
Definition GenParticle.h:38
Primary Vertex Finder.
ParametersT< TrackParametersDim, Charged, PerigeeSurface > Perigee
@ phi0
Definition ParamDefs.h:65
@ theta
Definition ParamDefs.h:66
@ qOverP
perigee
Definition ParamDefs.h:67
@ d0
Definition ParamDefs.h:63
@ z0
Definition ParamDefs.h:64
ParametersBase< TrackParametersDim, Charged > TrackParameters
EventInfo_v1 EventInfo
Definition of the latest event info version.