Belle II Software development
SVDdEdxCollectorModule.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#include <svd/modules/svddEdxCalibrationCollector/SVDdEdxCollectorModule.h>
9
10#include <analysis/dataobjects/ParticleList.h>
11#include <analysis/VariableManager/Manager.h>
12#include <reconstruction/dataobjects/VXDDedxTrack.h>
13#include <mdst/dataobjects/Track.h>
14#include <mdst/dataobjects/TrackFitResult.h>
15#include <mdst/dataobjects/HitPatternVXD.h>
16
17#include <TTree.h>
18#include <string>
19
20using namespace std;
21using namespace Belle2;
22using namespace Belle2::Variable;
23
24//-----------------------------------------------------------------
25// Register the Module
26//-----------------------------------------------------------------
27REG_MODULE(SVDdEdxCollector);
28
29//-----------------------------------------------------------------
30// Implementation
31//-----------------------------------------------------------------
32
34{
35 // Set module properties
36
37 setDescription("Collector module used to create the ROOT ntuples used to produce dE/dx calibration payloads");
39
40 addParam("LambdaListName", m_LambdaListName, "Name of the Lambda particle list", std::string("Lambda0:cut"));
41 addParam("DstarListName", m_DstarListName, "Name of the Dstar particle list", std::string("D*+:cut"));
42 addParam("GammaListName", m_GammaListName, "Name of the Gamma particle list", std::string("gamma:cut"));
43 addParam("GenericListName", m_GenericListName, "Name of the genericTrack track particle list", std::string("pi+:cut"));
44}
45
47{
48 B2INFO("Initialisation of the trees");
49 std::string objectNameLambda = "Lambda";
50 std::string objectNameDstar = "Dstar";
51 std::string objectNameGamma = "Gamma";
52 std::string objectNameGeneric = "Generic";
53 // Data object creation --------------------------------------------------
54 TTree* LambdaTree = new TTree(objectNameLambda.c_str(), "");
55 TTree* DstarTree = new TTree(objectNameDstar.c_str(), "");
56 TTree* GammaTree = new TTree(objectNameGamma.c_str(), "");
57 TTree* GenericTree = new TTree(objectNameGeneric.c_str(), "");
58
59 // Event info for all trees
60 LambdaTree->Branch<int>("event", &m_evt);
61 LambdaTree->Branch<int>("exp", &m_exp);
62 LambdaTree->Branch<int>("run", &m_run);
63 LambdaTree->Branch<double>("time", &m_time);
64
65 DstarTree->Branch<int>("event", &m_evt);
66 DstarTree->Branch<int>("exp", &m_exp);
67 DstarTree->Branch<int>("run", &m_run);
68 DstarTree->Branch<double>("time", &m_time);
69
70 GammaTree->Branch<int>("event", &m_evt);
71 GammaTree->Branch<int>("exp", &m_exp);
72 GammaTree->Branch<int>("run", &m_run);
73 GammaTree->Branch<double>("time", &m_time);
74
75 // Specific decay info for all trees
76 LambdaTree->Branch<double>("InvM", &m_InvMLambda);
77 LambdaTree->Branch<double>("DIRA", &m_CosDirAngleLambda);
78 LambdaTree->Branch<double>("ProtonMomentum", &m_protonMomentum);
79 LambdaTree->Branch<double>("ProtonSVDdEdx", &m_protonSVDdEdx);
80 LambdaTree->Branch<double>("ProtonSVDdEdxErr", &m_protonSVDdEdxErr);
81 LambdaTree->Branch<double>("ProtonSVDdEdxTrackMomentum", &m_protondEdxTrackMomentum);
82 LambdaTree->Branch<int>("ProtonSVDdEdxTrackNHits", &m_protondEdxTrackNHits);
83 LambdaTree->Branch<double>("ProtonSVDdEdxTrackCosTheta", &m_protondEdxTrackCosTheta);
84 LambdaTree->Branch<double>("ProtonnSVDHits", &m_protonnSVDHits);
85
86 LambdaTree->Branch<double>("PionLambdaMomentum", &m_pionLambdap);
87 LambdaTree->Branch<double>("PionLambdaSVDdEdx", &m_pionLambdaSVDdEdx);
88 LambdaTree->Branch<double>("PionLambdaSVDdEdxErr", &m_pionLambdaSVDdEdxErr);
89 LambdaTree->Branch<double>("PionLambdaSVDdEdxTrackMomentum", &m_pionLambdadEdxTrackMomentum);
90 LambdaTree->Branch<int>("PionLambdaSVDdEdxTrackNHits", &m_pionLambdadEdxTrackNHits);
91 LambdaTree->Branch<double>("PionLambdaSVDdEdxTrackCosTheta", &m_pionLambdadEdxTrackCosTheta);
92 LambdaTree->Branch<double>("PionLambdanSVDHits", &m_pionLambdanSVDHits);
93
94 DstarTree->Branch<double>("InvM", &m_InvMDstar);
95 DstarTree->Branch<double>("D0InvM", &m_InvMD0);
96 DstarTree->Branch<double>("deltaM", &m_DeltaM);
97 DstarTree->Branch<double>("D0DIRA", &m_CosDirAngleD0);
98 DstarTree->Branch<double>("KaonMomentum", &m_kaonMomentum);
99 DstarTree->Branch<double>("KaonSVDdEdx", &m_kaonSVDdEdx);
100 DstarTree->Branch<double>("KaonSVDdEdxErr", &m_kaonSVDdEdxErr);
101 DstarTree->Branch<double>("KaonSVDdEdxTrackMomentum", &m_kaondEdxTrackMomentum);
102 DstarTree->Branch<int>("KaonSVDdEdxTrackNHits", &m_kaondEdxTrackNHits);
103 DstarTree->Branch<double>("KaonSVDdEdxTrackCosTheta", &m_kaondEdxTrackCosTheta);
104 DstarTree->Branch<double>("KaonnSVDHits", &m_kaonnSVDHits);
105
106 DstarTree->Branch<double>("PionDMomentum", &m_pionDMomentum);
107 DstarTree->Branch<double>("PionDSVDdEdx", &m_pionDSVDdEdx);
108 DstarTree->Branch<double>("PionDSVDdEdxErr", &m_pionDSVDdEdxErr);
109 DstarTree->Branch<double>("PionDSVDdEdxTrackMomentum", &m_pionDdEdxTrackMomentum);
110 DstarTree->Branch<int>("PionDSVDdEdxTrackNHits", &m_pionDdEdxTrackNHits);
111 DstarTree->Branch<double>("PionDSVDdEdxTrackCosTheta", &m_pionDdEdxTrackCosTheta);
112 DstarTree->Branch<double>("PionDnSVDHits", &m_pionDnSVDHits);
113
114 DstarTree->Branch<double>("SlowPionMomentum", &m_slowPionMomentum);
115 DstarTree->Branch<double>("SlowPionSVDdEdx", &m_slowPionSVDdEdx);
116 DstarTree->Branch<double>("SlowPionSVDdEdxErr", &m_slowPionSVDdEdxErr);
117 DstarTree->Branch<double>("SlowPionSVDdEdxTrackMomentum", &m_slowPiondEdxTrackMomentum);
118 DstarTree->Branch<int>("SlowPionSVDdEdxTrackNHits", &m_slowPiondEdxTrackNHits);
119 DstarTree->Branch<double>("SlowPionSVDdEdxTrackCosTheta", &m_slowPiondEdxTrackCosTheta);
120 DstarTree->Branch<double>("SlowPionnSVDHits", &m_slowPionnSVDHits);
121
122 GammaTree->Branch<double>("InvM", &m_InvMGamma);
123 GammaTree->Branch<double>("dr", &m_drGamma);
124 GammaTree->Branch<double>("DIRA", &m_CosDirAngleGamma);
125 GammaTree->Branch<double>("FirstElectronMomentum", &m_firstElectronMomentum);
126 GammaTree->Branch<double>("FirstElectronSVDdEdx", &m_firstElectronSVDdEdx);
127 GammaTree->Branch<double>("FirstElectronSVDdEdxErr", &m_firstElectronSVDdEdxErr);
128 GammaTree->Branch<double>("FirstElectronSVDdEdxTrackMomentum", &m_firstElectrondEdxTrackMomentum);
129 GammaTree->Branch<int>("FirstElectronSVDdEdxTrackNHits", &m_firstElectrondEdxTrackNHits);
130 GammaTree->Branch<double>("FirstElectronSVDdEdxTrackCosTheta", &m_firstElectrondEdxTrackCosTheta);
131 GammaTree->Branch<double>("FirstElectronnSVDHits", &m_firstElectronnSVDHits);
132
133 GammaTree->Branch<double>("SecondElectronMomentum", &m_secondElectronMomentum);
134 GammaTree->Branch<double>("SecondElectronSVDdEdx", &m_secondElectronSVDdEdx);
135 GammaTree->Branch<double>("SecondElectronSVDdEdxErr", &m_secondElectronSVDdEdxErr);
136 GammaTree->Branch<double>("SecondElectronSVDdEdxTrackMomentum", &m_secondElectrondEdxTrackMomentum);
137 GammaTree->Branch<int>("SecondElectronSVDdEdxTrackNHits", &m_secondElectrondEdxTrackNHits);
138 GammaTree->Branch<double>("SecondElectronSVDdEdxTrackCosTheta", &m_secondElectrondEdxTrackCosTheta);
139 GammaTree->Branch<double>("SecondElectronnSVDHits", &m_secondElectronnSVDHits);
140
141 GenericTree->Branch<double>("TrackMomentum", &m_genericTrackMomentum);
142 GenericTree->Branch<double>("TrackSVDdEdx", &m_genericTrackSVDdEdx);
143 GenericTree->Branch<double>("TrackSVDdEdxErr", &m_genericTrackSVDdEdxErr);
144 GenericTree->Branch<double>("TrackSVDdEdxTrackMomentum", &m_genericTrackdEdxTrackMomentum);
145 GenericTree->Branch<int>("TrackSVDdEdxTrackNHits", &m_genericTrackdEdxTrackNHits);
146 GenericTree->Branch<double>("TracknSVDHits", &m_genericTracknSVDHits);
147
148 // We register the objects so that our framework knows about them.
149 // Don't try and hold onto the pointers or fill these objects directly
150 // Use the getObjectPtr functions to access collector objects
151 registerObject<TTree>(objectNameLambda, LambdaTree);
152 registerObject<TTree>(objectNameDstar, DstarTree);
153 registerObject<TTree>(objectNameGamma, GammaTree);
154 registerObject<TTree>(objectNameGeneric, GenericTree);
155}
156
157VXDDedxTrack const* getSVDDedxFromParticle(Particle const* particle)
158{
159 const Track* track = particle->getTrack();
160 if (!track) {
161 return nullptr;
162 }
163
164 const VXDDedxTrack* dedxTrack = track->getRelatedTo<VXDDedxTrack>();
165 if (!dedxTrack) {
166 return nullptr;
167 }
168 return dedxTrack;
169}
170
172{
173
174 m_evt = m_emd->getEvent();
175 m_run = m_emd->getRun();
176 m_exp = m_emd->getExperiment();
177 m_time = m_emd->getTime() / 1e9 / 3600.; // from ns to hours
178
183
184 if (!LambdaParticles.isValid() && !DstarParticles.isValid() && !GammaParticles.isValid())
185 return;
186
187 if (LambdaParticles->getListSize() > 0) {
188 for (unsigned int iParticle = 0; iParticle < LambdaParticles->getListSize(); ++iParticle) {
189
190 std::vector<int> indicesLambda = LambdaParticles->getParticle(0)->getDaughterIndices();
191 if (indicesLambda.size() != 2)
192 return;
193 const Particle* partLambda = LambdaParticles->getParticle(0);
194 const Particle* partPFromLambda = LambdaParticles->getParticle(0)->getDaughter(0);
195 const Particle* partPiFromLambda = LambdaParticles->getParticle(0)->getDaughter(1);
196
197 const VXDDedxTrack* dedxTrackPFromLambda = getSVDDedxFromParticle(partPFromLambda);
198 const VXDDedxTrack* dedxTrackPiFromLambda = getSVDDedxFromParticle(partPiFromLambda);
199
200 m_InvMLambda = partLambda->getMass();
201 m_CosDirAngleLambda = std::get<double>(Variable::Manager::Instance().getVariable(
202 std::string("cosAngleBetweenMomentumAndVertexVector"))->function(partLambda));
203
204 m_protonMomentum = partPFromLambda->getMomentumMagnitude();
205 if (!dedxTrackPFromLambda) {
206 m_protonSVDdEdx = -999.0;
207 m_protonSVDdEdxErr = -999.0;
211 } else {
212 m_protonSVDdEdx = dedxTrackPFromLambda->getDedx(Const::EDetector::SVD);
213 m_protonSVDdEdxErr = dedxTrackPFromLambda->getDedxError(Const::EDetector::SVD);
214 m_protondEdxTrackMomentum = dedxTrackPFromLambda->getMomentum();
215 m_protondEdxTrackNHits = (int) dedxTrackPFromLambda->size();
216 m_protondEdxTrackCosTheta = dedxTrackPFromLambda->getCosTheta();
217 }
218 auto trackFitPFromLambda = partPFromLambda->getTrackFitResult();
219 if (!trackFitPFromLambda) {m_protonnSVDHits = Const::doubleNaN;}
220 else {m_protonnSVDHits = trackFitPFromLambda->getHitPatternVXD().getNSVDHits();}
221
222 m_pionLambdap = partPiFromLambda->getMomentumMagnitude();
223 if (!dedxTrackPiFromLambda) {
224 m_pionLambdaSVDdEdx = -999.0;
225 m_pionLambdaSVDdEdxErr = -999.0;
229 } else {
230 m_pionLambdaSVDdEdx = dedxTrackPiFromLambda->getDedx(Const::EDetector::SVD);
231 m_pionLambdaSVDdEdxErr = dedxTrackPiFromLambda->getDedxError(Const::EDetector::SVD);
232 m_pionLambdadEdxTrackMomentum = dedxTrackPiFromLambda->getMomentum();
233 m_pionLambdadEdxTrackNHits = (int) dedxTrackPiFromLambda->size();
234 m_pionLambdadEdxTrackCosTheta = dedxTrackPiFromLambda->getCosTheta();
235 }
236 auto trackFitPiFromLambda = partPiFromLambda->getTrackFitResult();
237 if (!trackFitPiFromLambda) {m_pionLambdanSVDHits = Const::doubleNaN;}
238 else {m_pionLambdanSVDHits = trackFitPiFromLambda->getHitPatternVXD().getNSVDHits();}
239
240 getObjectPtr<TTree>("Lambda")->Fill();
241 }
242 }
243
244 if (DstarParticles->getListSize() > 0) {
245 for (unsigned int iParticle = 0; iParticle < DstarParticles->getListSize(); ++iParticle) {
246
247 std::vector<int> indicesDstar = DstarParticles->getParticle(0)->getDaughterIndices();
248 if (indicesDstar.size() != 2)
249 return;
250
251 const Particle* partDstar = DstarParticles->getParticle(0);
252 const Particle* partD0 = DstarParticles->getParticle(0)->getDaughter(0);
253 const Particle* partPiS = DstarParticles->getParticle(0)->getDaughter(1);
254 const Particle* partKFromD = DstarParticles->getParticle(0)->getDaughter(0)->getDaughter(0);
255 const Particle* partPiFromD = DstarParticles->getParticle(0)->getDaughter(0)->getDaughter(1);
256
257 const VXDDedxTrack* dedxTrackPiS = getSVDDedxFromParticle(partPiS);
258 const VXDDedxTrack* dedxTrackKFromD = getSVDDedxFromParticle(partKFromD);
259 const VXDDedxTrack* dedxTrackPiFromD = getSVDDedxFromParticle(partPiFromD);
260
261
262 m_InvMDstar = partDstar->getMass();
263 m_InvMD0 = partD0->getMass();
265 m_CosDirAngleD0 = std::get<double>(Variable::Manager::Instance().getVariable(
266 std::string("cosAngleBetweenMomentumAndVertexVector"))->function(partD0));;
267
268 m_kaonMomentum = partKFromD->getMomentumMagnitude();
269 if (!dedxTrackKFromD) {
270 m_kaonSVDdEdx = -999.0;
271 m_kaonSVDdEdxErr = -999.0;
275 } else {
276 m_kaonSVDdEdx = dedxTrackKFromD->getDedx(Const::EDetector::SVD);
277 m_kaonSVDdEdxErr = dedxTrackKFromD->getDedxError(Const::EDetector::SVD);
278 m_kaondEdxTrackMomentum = dedxTrackKFromD->getMomentum();
279 m_kaondEdxTrackNHits = (int) dedxTrackKFromD->size();
280 m_kaondEdxTrackCosTheta = dedxTrackKFromD->getCosTheta();
281 }
282
283 auto trackFitKFromD = partKFromD->getTrackFitResult();
284 if (!trackFitKFromD) {m_kaonnSVDHits = Const::doubleNaN;}
285 else {m_kaonnSVDHits = trackFitKFromD->getHitPatternVXD().getNSVDHits();}
286
287 m_pionDMomentum = partPiFromD->getMomentumMagnitude();
288 if (!dedxTrackPiFromD) {
289 m_pionDSVDdEdx = -999.0;
290 m_pionDSVDdEdxErr = -999.0;
294 } else {
295 m_pionDSVDdEdx = dedxTrackPiFromD->getDedx(Const::EDetector::SVD);
296 m_pionDSVDdEdxErr = dedxTrackPiFromD->getDedxError(Const::EDetector::SVD);
297 m_pionDdEdxTrackMomentum = dedxTrackPiFromD->getMomentum();
298 m_pionDdEdxTrackNHits = (int) dedxTrackPiFromD->size();
299 m_pionDdEdxTrackCosTheta = dedxTrackPiFromD->getCosTheta();
300 }
301
302 auto trackFitPiFromD = partPiFromD->getTrackFitResult();
303 if (!trackFitPiFromD) {m_pionDnSVDHits = Const::doubleNaN;}
304 else {m_pionDnSVDHits = trackFitPiFromD->getHitPatternVXD().getNSVDHits();}
305
307 if (!dedxTrackPiS) {
308 m_slowPionSVDdEdx = -999.0;
309 m_slowPionSVDdEdxErr = -999.0;
313 } else {
314 m_slowPionSVDdEdx = dedxTrackPiS->getDedx(Const::EDetector::SVD);
315 m_slowPionSVDdEdxErr = dedxTrackPiS->getDedxError(Const::EDetector::SVD);
317 m_slowPiondEdxTrackNHits = (int) dedxTrackPiS->size();
319 }
320 auto trackFitPiS = partPiS->getTrackFitResult();
321 if (!trackFitPiS) {m_slowPionnSVDHits = Const::doubleNaN;}
322 else {m_slowPionnSVDHits = trackFitPiS->getHitPatternVXD().getNSVDHits();}
323
324 getObjectPtr<TTree>("Dstar")->Fill();
325 }
326
327 }
328
329 if (GammaParticles->getListSize() > 0) {
330 for (unsigned int iParticle = 0; iParticle < GammaParticles->getListSize(); ++iParticle) {
331 std::vector<int> indicesGamma = GammaParticles->getParticle(0)->getDaughterIndices();
332 if (indicesGamma.size() != 2)
333 return;
334
335 const Particle* partGamma = GammaParticles->getParticle(0);
336 const Particle* partE1FromGamma = GammaParticles->getParticle(0)->getDaughter(0);
337 const Particle* partE2FromGamma = GammaParticles->getParticle(0)->getDaughter(1);
338
339 const VXDDedxTrack* dedxTrackE1FromGamma = getSVDDedxFromParticle(partE1FromGamma);
340 const VXDDedxTrack* dedxTrackE2FromGamma = getSVDDedxFromParticle(partE2FromGamma);
341
342 m_InvMGamma = partGamma->getMass();
343 m_drGamma = std::get<double>(Variable::Manager::Instance().getVariable(std::string("dr"))->function(partGamma));
344 m_CosDirAngleGamma = std::get<double>(Variable::Manager::Instance().getVariable(
345 std::string("cosAngleBetweenMomentumAndVertexVector"))->function(partGamma));
346
347 auto trackFitE1FromGamma = partE1FromGamma->getTrackFitResult();
348 if (!trackFitE1FromGamma) {m_firstElectronnSVDHits = Const::doubleNaN;}
349 else {m_firstElectronnSVDHits = trackFitE1FromGamma->getHitPatternVXD().getNSVDHits();}
350
351 auto trackFitE2FromGamma = partE2FromGamma->getTrackFitResult();
352 if (!trackFitE2FromGamma) {m_secondElectronnSVDHits = Const::doubleNaN;}
353 else {m_secondElectronnSVDHits = trackFitE2FromGamma->getHitPatternVXD().getNSVDHits();}
354
355
357 if (!dedxTrackE1FromGamma) {
358 m_firstElectronSVDdEdx = -999.0;
363 } else {
364 m_firstElectronSVDdEdx = dedxTrackE1FromGamma->getDedx(Const::EDetector::SVD);
365 m_firstElectronSVDdEdxErr = dedxTrackE1FromGamma->getDedxError(Const::EDetector::SVD);
366 m_firstElectrondEdxTrackMomentum = dedxTrackE1FromGamma->getMomentum();
367 m_firstElectrondEdxTrackNHits = (int) dedxTrackE1FromGamma->size();
368 m_firstElectrondEdxTrackCosTheta = dedxTrackE1FromGamma->getCosTheta();
369 }
370
372 if (!dedxTrackE2FromGamma) {
378 } else {
379 m_secondElectronSVDdEdx = dedxTrackE2FromGamma->getDedx(Const::EDetector::SVD);
380 m_secondElectronSVDdEdxErr = dedxTrackE2FromGamma->getDedxError(Const::EDetector::SVD);
381 m_secondElectrondEdxTrackMomentum = dedxTrackE2FromGamma->getMomentum();
382 m_secondElectrondEdxTrackNHits = (int) dedxTrackE2FromGamma->size();
383 m_secondElectrondEdxTrackCosTheta = dedxTrackE2FromGamma->getCosTheta();
384 }
385
386 getObjectPtr<TTree>("Gamma")->Fill();
387 }
388 }
389
390 if (GenericParticles->getListSize() > 0) {
391 for (unsigned int iParticle = 0; iParticle < GenericParticles->getListSize(); ++iParticle) {
392
393 const Particle* partGeneric = GenericParticles->getParticle(0);
394
395 const VXDDedxTrack* dedxTrackGeneric = getSVDDedxFromParticle(partGeneric);
396
398
399 if (!dedxTrackGeneric) {
400 m_genericTrackSVDdEdx = -999.0;
405 } else {
406 m_genericTrackSVDdEdx = dedxTrackGeneric->getDedx(Const::EDetector::SVD);
407 m_genericTrackSVDdEdxErr = dedxTrackGeneric->getDedxError(Const::EDetector::SVD);
408 m_genericTrackdEdxTrackMomentum = dedxTrackGeneric->getMomentum();
409 m_genericTrackdEdxTrackNHits = (int) dedxTrackGeneric->size();
410 m_genericTrackdEdxTrackCosTheta = dedxTrackGeneric->getCosTheta();
411 }
412
413 auto trackFitGeneric = partGeneric->getTrackFitResult();
414 if (!trackFitGeneric) {m_genericTracknSVDHits = Const::doubleNaN;}
415 else {m_genericTracknSVDHits = trackFitGeneric->getHitPatternVXD().getNSVDHits();}
416
417 getObjectPtr<TTree>("Generic")->Fill();
418 }
419 }
420
421
422}
void registerObject(std::string name, T *obj)
Register object with a name, takes ownership, do not access the pointer beyond prepare()
StoreObjPtr< EventMetaData > m_emd
Current EventMetaData.
CalibrationCollectorModule()
Constructor. Sets the default prefix for calibration dataobjects.
T * getObjectPtr(std::string name)
Calls the CalibObjManager to get the requested stored collector data.
static const double doubleNaN
quiet_NaN
Definition Const.h:703
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
Class to store reconstructed particles.
Definition Particle.h:76
const TrackFitResult * getTrackFitResult() const
Returns the pointer to the TrackFitResult that was used to create this Particle (ParticleType == c_Tr...
Definition Particle.cc:925
double getMomentumMagnitude() const
Returns momentum magnitude.
Definition Particle.h:589
double getMass() const
Returns invariant mass (= nominal for FS particles)
Definition Particle.h:527
double m_pionDnSVDHits
nSVDHits for the pion from the D0
double m_secondElectronSVDdEdx
SVD dE/dx response for the second electron.
double m_CosDirAngleLambda
Cosine of the angle between momentum and vertex vectors of Lambda candidates.
double m_pionDMomentum
momentum for the pion from the D0
void prepare() override final
Initialize the module.
double m_kaonSVDdEdx
SVD dE/dx response for the kaon from the D0.
double m_genericTrackMomentum
momentum for the generic track
double m_firstElectronSVDdEdxErr
SVD dE/dx uncertainty for the first electron.
double m_firstElectronMomentum
momentum for the first electron
double m_pionLambdadEdxTrackMomentum
momentum for the pion from the Lambda, taken from the VXDDedxTrack
double m_protonMomentum
momentum for the proton from the Lambda
double m_pionLambdap
momentum for the pion from the Lambda
double m_pionLambdaSVDdEdxErr
SVD dE/dx uncertainty for the pion from the Lambda.
double m_pionLambdanSVDHits
nSVDHits for the pion from the Lambda
double m_kaondEdxTrackCosTheta
cosTheta for the kaon from the D0, taken from the VXDDedxTrack
std::string m_LambdaListName
Name of the Lambda particle list.
double m_drGamma
dr of converted photon candidates
double m_protonSVDdEdxErr
SVD dE/dx uncertainty for the proton from the Lambda.
double m_protonSVDdEdx
SVD dE/dx response for the proton from the Lambda.
double m_kaonSVDdEdxErr
SVD dE/dx uncertainty for the kaon from the D0.
double m_InvMDstar
Invariant mass of Dstar candidates.
double m_slowPionnSVDHits
nSVDHits for the pion from the Dstar
double m_firstElectrondEdxTrackCosTheta
cosTheta for the first electron, taken from the VXDDedxTrack
double m_genericTracknSVDHits
nSVDHits for the generic track
int m_genericTrackdEdxTrackNHits
number of hits for generic track, taken from the VXDDedxTrack
double m_InvMD0
Invariant mass of D0 candidates.
double m_protondEdxTrackMomentum
momentum for the proton from the Lambda, taken from the VXDDedxTrack
double m_slowPionMomentum
momentum for the pion from the Dstar
double m_secondElectronnSVDHits
nSVDHits for the second electron
double m_slowPionSVDdEdxErr
SVD dE/dx uncertainty for the pion from the Dstar.
double m_slowPiondEdxTrackMomentum
momentum for the pion from the Dstar, taken from the VXDDedxTrack
int m_pionLambdadEdxTrackNHits
number of hits for pion from the Lambda, taken from the VXDDedxTrack
double m_firstElectronnSVDHits
nSVDHits for the first electron
double m_slowPiondEdxTrackCosTheta
cosTheta for the pion from the Dstar, taken from the VXDDedxTrack
double m_slowPionSVDdEdx
SVD dE/dx response for the pion from the Dstar.
void collect() override final
Event processor.
double m_firstElectrondEdxTrackMomentum
momentum for the first electron, taken from the VXDDedxTrack
double m_secondElectronSVDdEdxErr
SVD dE/dx uncertainty for the second electron.
double m_CosDirAngleGamma
Cosine of the angle between momentum and vertex vectors of converted photon candidates.
double m_pionDdEdxTrackMomentum
momentum for the pion from the D0, taken from the VXDDedxTrack
std::string m_GenericListName
Name of the generic track particle list.
int m_slowPiondEdxTrackNHits
number of hits for pion from the Dstar, taken from the VXDDedxTrack
double m_pionLambdadEdxTrackCosTheta
cosTheta for the pion from the Lambda, taken from the VXDDedxTrack
double m_genericTrackdEdxTrackMomentum
momentum for the generic track, taken from the VXDDedxTrack
double m_DeltaM
deltaM = m(Dstar)-m(D0)
double m_InvMGamma
Invariant mass of converted photon candidates.
double m_pionDSVDdEdxErr
SVD dE/dx uncertainty for the pion from the D0.
double m_kaonnSVDHits
nSVDHits for the kaon from the D0
double m_firstElectronSVDdEdx
SVD dE/dx response for the first electron.
double m_CosDirAngleD0
Cosine of the angle between momentum and vertex vectors of D0 candidates.
double m_genericTrackSVDdEdxErr
SVD dE/dx uncertainty for the generic track.
double m_protonnSVDHits
nSVDHits for the proton from the Lambda
int m_pionDdEdxTrackNHits
number of hits for pion from the D0, taken from the VXDDedxTrack
double m_secondElectrondEdxTrackCosTheta
cosTheta for the second electron, taken from the VXDDedxTrack
double m_protondEdxTrackCosTheta
cosTheta for the proton from the Lambda, taken from the VXDDedxTrack
double m_genericTrackSVDdEdx
SVD dE/dx for the generic track.
double m_InvMLambda
Invariant mass of Lambda candidates.
int m_protondEdxTrackNHits
number of hits for proton from the Lambda, taken from the VXDDedxTrack
double m_pionLambdaSVDdEdx
SVD dE/dx response for the pion from the Lambda.
std::string m_DstarListName
Name of the Dstar particle list.
double m_pionDdEdxTrackCosTheta
cosTheta for the pion from the D0, taken from the VXDDedxTrack
double m_pionDSVDdEdx
SVD dE/dx response for the pion from the D0.
double m_secondElectrondEdxTrackMomentum
momentum for the second electron, taken from the VXDDedxTrack
double m_kaondEdxTrackMomentum
momentum for the kaon from the D0, taken from the VXDDedxTrack
int m_firstElectrondEdxTrackNHits
number of hits for first electron, taken from the VXDDedxTrack
int m_kaondEdxTrackNHits
number of hits for kaon from the D0, taken from the VXDDedxTrack
int m_secondElectrondEdxTrackNHits
number of hits for second electron, taken from the VXDDedxTrack
std::string m_GammaListName
Name of the Gamma particle list.
double m_secondElectronMomentum
momentum for the second electron
double m_kaonMomentum
momentum for the kaon from the D0
double m_genericTrackdEdxTrackCosTheta
cosTheta for the generic track, taken from the VXDDedxTrack
Type-safe access to single objects in the data store.
Definition StoreObjPtr.h:96
bool isValid() const
Check whether the object was created.
Class that bundles various TrackFitResults.
Definition Track.h:25
Debug output for VXDDedxPID module.
double getCosTheta() const
Return cos(theta) for this yrack.
int size() const
Return the number of hits for this track.
double getMomentum() const
Return the momentum valid at the IP.
static Manager & Instance()
get singleton instance.
Definition Manager.cc:26
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 getDedxError(Const::EDetector detector) const
Get the error on the dE/dx truncated mean for given detector.
double getDedx(Const::EDetector detector) const
Get dE/dx truncated mean for given detector.
VXDDedxTrack const * getSVDDedxFromParticle(Particle const *particle)
SVD dEdx value from particle.
Abstract base class for different kinds of events.
STL namespace.