Belle II Software  release-08-01-10
ParallelTrackFilterModule.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 <tracking/modules/trackFilter/ParallelTrackFilterModule.h>
10 #include <mdst/dataobjects/TrackFitResult.h>
11 #include <mdst/dataobjects/HitPatternVXD.h>
12 #include <mdst/dataobjects/HitPatternCDC.h>
13 #include <framework/datastore/StoreArray.h>
14 
15 using namespace Belle2;
16 
17 //-----------------------------------------------------------------
18 // Register the Module
19 //-----------------------------------------------------------------
20 REG_MODULE(ParallelTrackFilter);
21 
22 //-----------------------------------------------------------------
23 // Implementation
24 //-----------------------------------------------------------------
25 
27 {
28  // Set module properties
29  setDescription("Generates a new StoreArray from the input StoreArray which contains only tracks that meet the specified criteria.");
31 
32  // Parameter definitions
33  addParam("inputArrayName", m_inputArrayName, "StoreArray with the input tracks", m_inputArrayName);
34  addParam("outputINArrayName", m_outputINArrayName, "Output StoreArray with the tracks that pass the cuts", m_outputINArrayName);
35  addParam("outputOUTArrayName", m_outputOUTArrayName, "Output StoreArray with the tracks that do not pass the cuts",
37 
38  //selection parameter definition
39  addParam("min_d0", m_minD0, "minimum value of the d0", m_minD0);
40  addParam("max_d0", m_maxD0, "maximum value of the d0", m_maxD0);
41  addParam("min_z0", m_minZ0, "minimum value of the z0", m_minZ0);
42  addParam("max_z0", m_maxZ0, "maximum value of the z0", m_maxZ0);
43  addParam("min_pCM", m_minPCM, "minimum value of the center-of-mass-momentum", m_minPCM);
44  addParam("min_pT", m_minPT, "minimum value of the transverse momentum", m_minPT);
45  addParam("min_Pvalue", m_minPval, "minimum value of the P-Value of the track fit", m_minPval);
46  addParam("min_NumHitPXD", m_minNumHitsPXD, "minimum number of PXD hits associated to the trcak", m_minNumHitsPXD);
47  addParam("min_NumHitSVD", m_minNumHitsSVD, "minimum number of SVD hits associated to the trcak", m_minNumHitsSVD);
48  addParam("min_NumHitCDC", m_minNumHitsCDC, "minimum number of CDC hits associated to the trcak", m_minNumHitsCDC);
49 
50 }
51 
52 
54 {
55  B2DEBUG(22, "ParallelTrackFilterModule " + getName() + " parameters:"
56  << LogVar("inputArrayName", m_inputArrayName)
57  << LogVar("outputINArrayName", m_outputINArrayName)
58  << LogVar("outputOUTArrayName", m_outputOUTArrayName));
59 
60  // Can't initialize if the input array is not present
62  if (!inputArray.isOptional()) {
63  B2WARNING("Missing input tracks array, " + getName() + " is skipped"
64  << LogVar("inputArrayName", m_inputArrayName));
65  return;
66  }
67 
68  m_selectedTracks.registerSubset(inputArray, m_outputINArrayName);
69  m_selectedTracks.inheritAllRelations();
70 
71  m_notSelectedTracks.registerSubset(inputArray, m_outputOUTArrayName);
72  m_notSelectedTracks.inheritAllRelations();
73 }
74 
75 
77 {
78 
80  if (!inputArray.isOptional() || !m_selectedTracks.getSet()) {
81  B2DEBUG(22, "Missing Tracks array, " + getName() + " is skipped." << LogVar("inputArrayName", m_inputArrayName));
82  return;
83  }
84 
85  m_selectedTracks.select([this](const Track * track) {
86  return this->isSelected(track);
87  });
88 
89  m_notSelectedTracks.select([this](const Track * track) {
90  return !this->isSelected(track);
91  });
92 
93 }
94 
96 {
97  const TrackFitResult* tfr = track->getTrackFitResultWithClosestMass(Const::pion);
98  if (tfr == nullptr)
99  return false;
100 
101  if (tfr->getD0() < m_minD0 || tfr->getD0() > m_maxD0)
102  return false;
103 
104  if (tfr->getZ0() < m_minZ0 || tfr->getZ0() > m_maxZ0)
105  return false;
106 
107  if (tfr->getPValue() < m_minPval)
108  return false;
109 
110  if (tfr->getMomentum().Rho() < m_minPT)
111  return false;
112 
113  HitPatternVXD hitPatternVXD = tfr->getHitPatternVXD();
114  if (hitPatternVXD.getNSVDHits() < m_minNumHitsSVD || hitPatternVXD.getNPXDHits() < m_minNumHitsPXD)
115  return false;
116 
117  HitPatternCDC hitPatternCDC = tfr->getHitPatternCDC();
118  if (hitPatternCDC.getNHits() < m_minNumHitsCDC)
119  return false;
120 
121  return true;
122 }
static const ChargedStable pion
charged pion particle
Definition: Const.h:652
Hit pattern of CDC hits within a track.
Definition: HitPatternCDC.h:35
unsigned short getNHits() const
Get the total Number of CDC hits in the fit.
Hit pattern of the VXD within a track.
Definition: HitPatternVXD.h:37
unsigned short getNPXDHits() const
Get total number of hits in the PXD.
unsigned short getNSVDHits() const
Get total number of hits in the SVD.
Base class for Modules.
Definition: Module.h:72
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
const std::string & getName() const
Returns the name of the module.
Definition: Module.h:187
std::string m_outputOUTArrayName
StoreArray with the NOT selected output tracks.
SelectSubset< Track > m_selectedTracks
selected tracks
std::string m_outputINArrayName
StoreArray with the selected output tracks.
int m_minNumHitsCDC
miminum value of CDC hits
int m_minNumHitsPXD
miminum value of PXD hits
virtual void initialize() override
init the module
virtual void event() override
processes the event
bool isSelected(const Track *track)
determine if the track satisfies the selection criteria
double m_minPval
miminum P-value of the track fit
SelectSubset< Track > m_notSelectedTracks
not selected tracks
double m_minPCM
miminum value of the center of mass momentum
double m_minPT
miminum value of the transverse momentum
int m_minNumHitsSVD
miminum value of SVD hits
std::string m_inputArrayName
StoreArray with the input tracks.
ParallelTrackFilterModule()
Constructor: Sets the description, the properties and the parameters of the module.
bool isOptional(const std::string &name="")
Tell the DataStore about an optional input.
Accessor to arrays stored in the data store.
Definition: StoreArray.h:113
Values of the result of a track fit with a given particle hypothesis.
double getPValue() const
Getter for Chi2 Probability of the track fit.
double getD0() const
Getter for d0.
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.
HitPatternCDC getHitPatternCDC() const
Getter for the hit pattern in the CDC;.
HitPatternVXD getHitPatternVXD() const
Getter for the hit pattern in the VXD;.
Class that bundles various TrackFitResults.
Definition: Track.h:25
Class to store variables with their name which were sent to the logging service.
REG_MODULE(arichBtest)
Register the Module.
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
Abstract base class for different kinds of events.