Belle II Software  release-06-02-00
FastBDTClassifierTrainingModule.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/vxdtfRedesign/FastBDTClassifierTrainingModule.h>
10 #include <tracking/trackFindingVXD/filterTools/FBDTClassifier.h>
11 #include <tracking/trackFindingVXD/segmentNetwork/DirectedNodeNetwork.h>
12 #include <tracking/trackFindingVXD/segmentNetwork/TrackNode.h>
13 
14 #include <tracking/spacePointCreation/PurityCalculatorTools.h>
15 #include <tracking/spacePointCreation/MCVXDPurityInfo.h>
16 
17 #include <framework/logging/Logger.h>
18 
19 #include <fstream>
20 #include <sstream>
21 
22 using namespace Belle2;
23 using namespace std;
24 
25 REG_MODULE(FastBDTClassifierTraining);
26 
28 {
29  setDescription("TODO");
30 
31  addParam("outputFileName",
33  "Output file name to which the trained FBDTClassifier will be stored",
34  std::string("FBDTClassifier.dat"));
35 
36  addParam("networkInputName", m_PARAMnetworkInputName,
37  "Name of the StoreObjPtr where the network container used in this module is stored", std::string(""));
38 
39  addParam("train", m_PARAMdoTrain, "Set if the module should train a classifier after collecting or not", true);
40  addParam("nTrees", m_PARAMnTrees, "Number of Trees used in the FastBDT", 100);
41  addParam("treeDepth", m_PARAMtreeDepth, "Tree depth of the trees used in the FastBDT", 3);
42  addParam("shrinkage", m_PARAMshrinkage, "Shrinkage parameter used in the FastBDT", 0.15);
43  addParam("randRatio", m_PARAMrandRatio, "ratio of randomly chosen samples for training of one tree", 0.5);
44  addParam("storeSamples", m_PARAMstoreSamples, "store the collected samples into a file", false);
45  addParam("samplesFileName", m_PARAMsamplesFileName, "the file name into which/from whicht the collected samples are stored/read",
46  std::string("FBDTClassifier_samples.dat"));
47  addParam("useSamples", m_PARAMuseSamples,
48  "use samples for training that have been collected previously and bypass the collection of samples", false);
49 }
50 
52 {
53 
55  B2ERROR("storeSamples and useSamples are both set to true. However, only one option can be set at a time");
56  }
57 
58  if (m_PARAMnTrees < 1) {
59  B2WARNING("nTrees was set to " << m_PARAMnTrees << ". Has to be at least 1. Setting to 1.");
60  m_PARAMnTrees = 1;
61  }
62 
63  if (m_PARAMtreeDepth < 0) {
64  B2WARNING("Trees have to be at least a stump, but treeDepth was set to " << m_PARAMtreeDepth << ". Setting to 3 (default).");
65  m_PARAMtreeDepth = 3;
66  }
67 
68  if (m_PARAMshrinkage < 0 || m_PARAMshrinkage > 1) { // TODO: check this
69  B2WARNING("shrinkage has to be in [0,1] but was " << m_PARAMrandRatio << ". Setting to 0.15 (default).");
70  m_PARAMshrinkage = .15;
71  }
72 
73  if (m_PARAMrandRatio < 0 || m_PARAMrandRatio > 1) {
74  B2WARNING("randRatio has to be in [0,1] but was " << m_PARAMrandRatio << ". Setting to 0.5 (default).");
75  m_PARAMrandRatio = 0.5;
76  }
77 
78  if (m_PARAMuseSamples) {
79  ifstream sampFile(m_PARAMsamplesFileName);
80  if (!sampFile.is_open() || !sampFile.good()) {
81  B2ERROR("Was not able to open the samples file: " << m_PARAMsamplesFileName);
82  }
83 
85  } else { // only if no samples are provided the collection from the DNN is necessary
87  }
88 }
89 
91 {
92  if (m_PARAMuseSamples) return; // don't collect anything during event if samples are provided
93 
95 
96  // B2DEBUG(1, "size of hitNetwork " << hitNetwork.getNodes().size());
97 
98  size_t samplePriorEvent = m_samples.size();
99 
100  // XXXHit is of type DirectedNode<TrackNode, VoidMetaInfo>
101  for (const auto& outerHit : hitNetwork.getNodes()) { // loop over all outer nodes
102  for (const auto& centerHit : outerHit->getInnerNodes()) { // loop over all center nodes attached to outer node
103  for (const auto& innerHit : centerHit->getInnerNodes()) { // loop over all inner nodes attached to center node
104  m_samples.push_back(makeTrainSample(outerHit->getEntry().m_spacePoint,
105  centerHit->getEntry().m_spacePoint,
106  innerHit->getEntry().m_spacePoint));
107  } // inner node loop
108  } // center node loop
109  } // outer node loop
110 
111  B2DEBUG(10, "collected " << m_samples.size() - samplePriorEvent << " training samples in this event");
112 
113 }
114 
116 {
117  if (m_PARAMstoreSamples) {
118  B2DEBUG(1, "Storing the collected samples to file: " << m_PARAMsamplesFileName);
119  ofstream sampStream(m_PARAMsamplesFileName);
120  sampStream.precision(16); // increase precision for sample writeout
121  writeSamplesToStream(sampStream, m_samples);
122  sampStream.close();
123  }
124  if (m_PARAMdoTrain) {
125  FBDTClassifier<9> classifier{};
126  B2DEBUG(1, "Training a FBDTClassifier with " << m_samples.size() << " input samples. Training Parameters: \n" <<
127  "nTrees: " << m_PARAMnTrees << "\n" <<
128  "treeDetph: " << m_PARAMtreeDepth << "\n" <<
129  "shrinkage: " << m_PARAMshrinkage << "\n" <<
130  "randRatio: " << m_PARAMrandRatio << "\n");
132 
133  ofstream ofs(m_PARAMfbdtOutFileName);
134  classifier.writeToStream(ofs);
135  ofs.close();
136  }
137 }
138 
141  const Belle2::SpacePoint* inner)
142 {
143  vector<MCVXDPurityInfo> purityInfos = createPurityInfosVec({outer, center, inner});
144  auto mcId = purityInfos[0].getPurity(); // there is at least one entry in this vector!
145  bool signal = (mcId.first >= 0 && mcId.second == 1); // only assign true for actual MCParticle and purity 1
146 
147  std::array<double, 9> coords {{
148  inner->X(), inner->Y(), inner->Z(),
149  center->X(), center->Y(), center->Z(),
150  outer->X(), outer->Y(), outer->Z()
151  }};
152 
153  TrainSample sample(coords, signal);
154 
155  if (LogSystem::Instance().isLevelEnabled(LogConfig::c_Debug, 499, PACKAGENAME())) {
156  std::stringstream coordOutput;
157  for (double d : sample.hits) coordOutput << d << " ";
158 
159  B2DEBUG(499, "Created TrainingsSample with coordinates: ( " << coordOutput.str() << " ) " << sample.signal);
160  }
161 
162  return sample;
163 }
164 
165 // void FastBDTClassifierTrainingModule::readSamplesFromStream(std::istream& is)
166 // {
167 // std::string line;
168 // while(!is.eof()) {
169 // getline(is, line);
170 // if(line.empty()) break;
171 // stringstream ss(line);
172 // std::array<double, 9> coords;
173 // for(double& c : coords) ss >> c;
174 // bool sig; ss >> sig;
175 
176 // m_samples.push_back(FBDTTrainSample<9>(coords, sig));
177 // }
178 
179 // B2INFO("Read in " << m_samples.size() << " training samples.");
180 // }
181 
182 // void FastBDTClassifierTrainingModule::writeSamplesToStream(std::ostream& os) const
183 // {
184 // for (const auto& event : m_samples) {
185 // for (const auto& val : event.hits) {
186 // os << val << " ";
187 // }
188 // os << event.signal << std::endl;
189 // }
190 // B2INFO("Wrote out " << m_samples.size() << " training samples.");
191 // }
DirectedNodeNetwork< Belle2::TrackNode, Belle2::VoidMetaInfo > & accessHitNetwork()
Returns reference to the HitNetwork stored in this container, intended for read and write access.
Network of directed nodes of the type EntryType.
std::vector< Node * > & getNodes()
Returns all nodes of the network.
double m_PARAMrandRatio
ratio of samples to be used for training one tree in the FastBDT.
std::string m_PARAMsamplesFileName
filename to be used to store / read collect samples
void initialize() override
initialize the module
Belle2::StoreObjPtr< Belle2::DirectedNodeNetworkContainer > m_network
StoreObjPtr to access the DNNs that are used in this module.
void event() override
collect all possible combinations and store them
std::vector< TrainSample > m_samples
vector in which all samples are collected on the fly in event.
const TrainSample makeTrainSample(const Belle2::SpacePoint *outerHit, const Belle2::SpacePoint *centerHit, const Belle2::SpacePoint *innerHit)
create a trainings sample from the three hit combination
void terminate() override
take the collected data and train a FBDTClassifier and store it in the given output file
bool m_PARAMdoTrain
actually train a classifier or only do collection
std::string m_PARAMfbdtOutFileName
output file name into which the FBDTClassifier is stored.
double m_PARAMshrinkage
shrinkage parameter of FastBDT.
bool m_PARAMstoreSamples
store the collected samples into a file
std::string m_PARAMnetworkInputName
name of the StoreObjPtr in which the network container is stored which contains the network that is u...
bool m_PARAMuseSamples
use pre-collected samples for training and bypass the collection step
int m_PARAMnTrees
number of trees in the FastBDT.
@ c_Debug
Debug: for code development.
Definition: LogConfig.h:26
static LogSystem & Instance()
Static method to get a reference to the LogSystem instance.
Definition: LogSystem.cc:31
Base class for Modules.
Definition: Module.h:72
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
SpacePoint typically is build from 1 PXDCluster or 1-2 SVDClusters.
Definition: SpacePoint.h:42
double Z() const
return the z-value of the global position of the SpacePoint
Definition: SpacePoint.h:129
double X() const
return the x-value of the global position of the SpacePoint
Definition: SpacePoint.h:123
double Y() const
return the y-value of the global position of the SpacePoint
Definition: SpacePoint.h:126
bool isRequired(const std::string &name="")
Ensure this array/object has been registered previously.
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
static void readSamplesFromStream(std::istream &is, std::vector< FBDTTrainSample< Ndims > > &samples)
read samples from stream and append them to samples
static std::vector< Belle2::MCVXDPurityInfo > createPurityInfosVec(const std::vector< const Belle2::SpacePoint * > &spacePoints)
create a vector of MCVXDPurityInfos objects for a std::vector<Belle2::SpacePoints>.
static void writeSamplesToStream(std::ostream &os, const std::vector< FBDTTrainSample< Ndims > > &samples)
write all samples to stream
Abstract base class for different kinds of events.
bundle together the classifier input and the target value into one struct for easier passing around.