Belle II Software  release-08-01-10
KLMExpertModule.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 <reconstruction/modules/KlId/KLMExpert/KLMExpertModule.h>
9 #include <mdst/dataobjects/KlId.h>
10 #include <framework/datastore/StoreArray.h>
11 #include <framework/logging/Logger.h>
12 
13 #include <mdst/dataobjects/ECLCluster.h>
14 #include <tracking/dataobjects/TrackClusterSeparation.h>
15 
16 #include <mva/interface/Interface.h>
17 #include <boost/algorithm/string/predicate.hpp>
18 
19 // here's where the functions are hidden
20 #include "reconstruction/modules/KlId/KLMExpert/KlId.h"
21 
22 using namespace std;
23 using namespace Belle2;
24 using namespace Belle2::KlongId;
25 
26 REG_MODULE(KLMExpert);
27 
28 KLMExpertModule::KLMExpertModule(): Module(), m_feature_variables(18, 0) //12
29 {
30  setDescription("Use to calculate KlId for each KLM cluster.");
31  addParam("classifierPath", m_identifier,
32  "path to the classifier you want to use. It is recommended to use the default classifiers and not to mess around with this.",
33  m_identifier);
35 }
36 
37 
38 
40 {
41 }
42 
43 
44 // --------------------------------------Module----------------------------------------------
46 {
47  // require existence of necessary datastore obj
48 
49  m_klmClusters.isRequired();
50 
51  m_klids.registerInDataStore();
52  m_klmClusters.registerRelationTo(m_klids);
53 
54 
55  if (not(boost::ends_with(m_identifier, ".root") or boost::ends_with(m_identifier, ".xml"))) {
56  m_weightfile_representation = std::unique_ptr<DBObjPtr<DatabaseRepresentationOfWeightfile>>(new
58  }
59 
61 
62 }
63 
64 
66 {
68  if (m_weightfile_representation->hasChanged()) {
69  std::stringstream ss((*m_weightfile_representation)->m_data);
70  auto weightfile = MVA::Weightfile::loadFromStream(ss);
71  init_mva(weightfile);
72  }
73  } else {
75  init_mva(weightfile);
76 
77  }
78 }
81 {
82 
83  auto supported_interfaces = MVA::AbstractInterface::getSupportedInterfaces();
84  MVA::GeneralOptions general_options;
85  weightfile.getOptions(general_options);
86 
87  m_expert = supported_interfaces[general_options.m_method]->getExpert();
88  m_expert->load(weightfile);
89 
90  std::vector<float> dummy;
91  dummy.resize(m_feature_variables.size(), 0);
92  m_dataset = std::unique_ptr<MVA::SingleDataset>(new MVA::SingleDataset(general_options, std::move(dummy), 0));
93 
94 }
95 
96 
98 {
99  // Use the neutralHadron hypothesis for the ECL
101 
102  //overwritten at the end of the cluster loop
103  KlId* klid = nullptr;
104 
105  // loop thru clusters in event and classify
106  for (KLMCluster& cluster : m_klmClusters) {
107 
108  const ROOT::Math::XYZVector& clusterPos = cluster.getClusterPosition();
109 
110  // get various KLMCluster vars
111  m_KLMglobalZ = clusterPos.Z();
112  m_KLMnCluster = m_klmClusters.getEntries();
113  m_KLMnLayer = cluster.getLayers();
114  m_KLMnInnermostLayer = cluster.getInnermostLayer();
115  m_KLMtime = cluster.getTime();
116  m_KLMenergy = cluster.getEnergy();
117  m_KLMhitDepth = cluster.getClusterPosition().R();
118 
119  // find nearest ecl cluster and calculate distance
120  pair<ECLCluster*, double> closestECLAndDist = findClosestECLCluster(clusterPos, eclHypothesis);
121  ECLCluster* closestECLCluster = get<0>(closestECLAndDist);
122  m_KLMECLDist = get<1>(closestECLAndDist);
123 
124  // get variables of the closest ECL cluster might be removed in future
125  if (!(closestECLCluster == nullptr)) {
126  m_KLMECLE = closestECLCluster -> getEnergy(eclHypothesis);
127  m_KLMECLE9oE25 = closestECLCluster -> getE9oE21();
128  m_KLMECLTerror = closestECLCluster -> getDeltaTime99();
129  m_KLMECLTiming = closestECLCluster -> getTime();
130  m_KLMECLEerror = closestECLCluster -> getUncertaintyEnergy();
131  m_KLMECLdeltaL = closestECLCluster -> getDeltaL();
132  m_KLMECLminTrackDist = closestECLCluster -> getMinTrkDistance();
133  m_KLMECLZMVA = closestECLCluster -> getZernikeMVA();
134  m_KLMECLZ40 = closestECLCluster -> getAbsZernike40();
135  m_KLMECLZ51 = closestECLCluster -> getAbsZernike51();
136  } else {
137  m_KLMECLdeltaL = -999;
138  m_KLMECLminTrackDist = -999;
139  m_KLMECLE = -999;
140  m_KLMECLE9oE25 = -999;
141  m_KLMECLTiming = -999;
142  m_KLMECLTerror = -999;
143  m_KLMECLEerror = -999;
144  m_KLMECLZMVA = -999;
145  m_KLMECLZ40 = -999;
146  m_KLMECLZ51 = -999;
147  }
148 
149  // calculate distance to next cluster
150  tuple<const KLMCluster*, double, double> closestKLMAndDist = findClosestKLMCluster(clusterPos);
151  m_KLMnextCluster = get<1>(closestKLMAndDist);
152  m_KLMavInterClusterDist = get<2>(closestKLMAndDist);
153 
154  m_KLMTrackSepDist = -999;
155  m_KLMTrackSepAngle = -999;
159 
160  auto trackSeperations = cluster.getRelationsTo<TrackClusterSeparation>();
161  float best_dist = 1e10;
162  for (auto trackSeperation : trackSeperations) {
163  float dist = trackSeperation.getDistance();
164  if (dist < best_dist) {
165  best_dist = dist;
166  m_KLMTrackSepDist = trackSeperation.getDistance();
167  m_KLMTrackSepAngle = trackSeperation.getTrackClusterAngle();
168  m_KLMInitialTrackSepAngle = trackSeperation.getTrackClusterInitialSeparationAngle();
169  m_KLMTrackRotationAngle = trackSeperation.getTrackRotationAngle();
170  m_KLMTrackClusterSepAngle = trackSeperation.getTrackClusterSeparationAngle();
171  }
172  }
173 
174 // reduced vars set
193 
194 
195  for (unsigned int i = 0; i < m_feature_variables.size(); ++i) {
196  if (!std::isfinite(m_feature_variables[i])) { m_feature_variables[i] = -999; }
197  m_dataset->m_input[i] = m_feature_variables[i];
198  }
199 
200  double IDMVAOut = m_expert->apply(*m_dataset)[0];
201  B2DEBUG(175, "KLM Expert classification: " << IDMVAOut);
202  klid = m_klids.appendNew();
203  cluster.addRelationTo(klid, IDMVAOut);
204 
205  }// for cluster in clusters
206 } // event
Class for accessing objects in the database.
Definition: DBObjPtr.h:21
ECL cluster data.
Definition: ECLCluster.h:27
EHypothesisBit
The hypothesis bits for this ECLCluster (Connected region (CR) is split using this hypothesis.
Definition: ECLCluster.h:31
@ c_neutralHadron
CR is reconstructed as a neutral hadron (N2)
KLM cluster data.
Definition: KLMCluster.h:28
float m_KLMECLE
energy measured in associated ECL cluster
float m_KLMECLdeltaL
distance between track entry pofloat and cluster center, might be removed
StoreArray< KlId > m_klids
storearray
float m_KLMInitialTrackSepAngle
angular distance from track to cluster at track starting point
float m_KLMECLZMVA
output of a BDT fitted on various Z-moments for the closest ECL cluster
float m_KLMTrackRotationAngle
angle between track at poca and trackbeginning
float m_KLMECLTiming
timing of associated ECL cluster
std::unique_ptr< MVA::SingleDataset > m_dataset
Pointer to the current dataset.
StoreArray< KLMCluster > m_klmClusters
storearray
virtual void initialize() override
init
float m_KLMECLTerror
uncertanty on time in associated ECL cluster
virtual void event() override
process event
float m_KLMtime
timing of KLM Cluster
float m_KLMECLDist
distance associated ECL <-> KLM cluster, extrapolated by genfit
float m_KLMnLayer
number of layers hit in KLM cluster
float m_KLMhitDepth
hit depth in KLM, distance to IP
float m_KLMECLminTrackDist
track distance between associated ECL cluster and track extrapolated into ECL
float m_KLMenergy
Energy deposit in KLM (0.2 GeV * nHitCells)
float m_KLMTrackSepAngle
angular distance from track separation object.
std::unique_ptr< MVA::Expert > m_expert
Pointer to the current MVA Expert.
float m_KLMECLZ40
zernike moment 4,0 of closest ECL
std::unique_ptr< DBObjPtr< DatabaseRepresentationOfWeightfile > > m_weightfile_representation
Database pointer to the Database representation of the weightfile.
virtual ~KLMExpertModule()
Destructor.
float m_KLMnInnermostLayer
number of innermost layers hit
float m_KLMECLEerror
uncertanty on E in associated ECL cluster
float m_KLMECLE9oE25
E in surrounding 9 crystals divided by surrounding 25 crydtalls.
std::vector< float > m_feature_variables
vars to be classified
float m_KLMECLZ51
zernike moment 5,1 of closest ECL
float m_KLMavInterClusterDist
average distance between all KLM clusters
virtual void beginRun() override
beginn run
float m_KLMTrackSepDist
distance from track separation object
float m_KLMglobalZ
global Z position in KLM
float m_KLMnCluster
varibales to write out.
float m_KLMnextCluster
distance to next KLM cluster
float m_KLMTrackClusterSepAngle
angle between trach momentum and cluster (measured from ip)
void init_mva(MVA::Weightfile &weightfile)
Initialize mva expert, dataset and features Called everytime the weightfile in the database changes i...
std::string m_identifier
mva identifier.
Klong identifcation (KlId) datastore object to store results from KlId calculations.
Definition: KlId.h:23
static std::map< std::string, AbstractInterface * > getSupportedInterfaces()
Returns interfaces supported by the MVA Interface.
Definition: Interface.h:53
static void initSupportedInterfaces()
Static function which initliazes all supported interfaces, has to be called once before getSupportedI...
Definition: Interface.cc:45
General options which are shared by all MVA trainings.
Definition: Options.h:62
Wraps the data of a single event into a Dataset.
Definition: Dataset.h:135
The Weightfile class serializes all information about a training into an xml tree.
Definition: Weightfile.h:38
static Weightfile loadFromStream(std::istream &stream)
Static function which deserializes a Weightfile from a stream.
Definition: Weightfile.cc:251
void getOptions(Options &options) const
Fills an Option object from the xml tree.
Definition: Weightfile.cc:67
static Weightfile loadFromFile(const std::string &filename)
Static function which loads a Weightfile from a file.
Definition: Weightfile.cc:206
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
Store one Track-KLMCluster separation as a ROOT object.
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
Helper functions for all klid modules to improve readability of the code.
Definition: KlId.h:27
std::tuple< const Belle2::KLMCluster *, double, double > findClosestKLMCluster(const ROOT::Math::XYZVector &klmClusterPosition)
find nearest KLMCluster, tis distance and the av intercluster distance
Definition: KlId.h:236
std::pair< Belle2::ECLCluster *, double > findClosestECLCluster(const ROOT::Math::XYZVector &klmClusterPosition, const Belle2::ECLCluster::EHypothesisBit eclhypothesis=Belle2::ECLCluster::EHypothesisBit::c_neutralHadron)
Find the closest ECLCluster with a neutral hadron hypothesis, and return it with its distance.
Definition: KlId.h:203
Abstract base class for different kinds of events.