Belle II Software  release-05-02-19
PhaseSpaceAnalysisModule.cc
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2014 Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Thomas Madlener *
7  * *
8  **************************************************************************/
9 
10 #include <tracking/modules/spacePointCreator/PhaseSpaceAnalysisModule.h>
11 
12 #include <genfit/TrackCand.h>
13 #include <tracking/spacePointCreation/SpacePointTrackCand.h>
14 
15 #include <framework/datastore/StoreArray.h>
16 #include <framework/dataobjects/EventMetaData.h>
17 #include <framework/datastore/StoreObjPtr.h>
18 
19 #include <algorithm>
20 #include <numeric>
21 // #include <mdst/dataobjects/MCParticle.h> // already in header
22 
23 #include <TVector3.h> // COULDDO: B2Vector3.h
24 
25 using namespace std;
26 using namespace Belle2;
27 
28 REG_MODULE(PhaseSpaceAnalysis)
29 
31 {
32  setDescription("Module for analysing the phase space covered by TrackCands (resp. their related MCParticles)");
33 
34  std::vector<std::string> defaultNames = { std::string("") };
35  addParam("containerNames", m_PARAMcontainerNames, "Collection names of trackCands to be analized.", defaultNames);
36  addParam("trackCandTypes", m_PARAMtrackCandTypes, "Either 'GFTC' for genfit::TrackCand or 'SPTC' for SpacePointTrackCand." \
37  "If there are more 'containerNames' than 'trackCandTypes' the last type is assumed for all the names that cannot be matched"\
38  "one-to-one with a trackCandType. If there are as much names as types the i-th name is assumed to be of the i-th type");
39 
40  addParam("diffAnalysis", m_PARAMdiffAnalysis,
41  "Use the first container as reference and only collect data for those MCParticle that"\
42  " cannot be found in the other containers.", false);
43 
44  std::vector<std::string> defaultRootFName = { "PhaseSpaceAnalysis", "RECREATE" };
45  addParam("rootFileName", m_PARAMrootFileName, "Name of the output root file without '.root' file ending + write mode"\
46  " ('UPDATE' or 'RECREATE')", defaultRootFName);
47 
48  // initialize pointers to NULL (cppcheck)
49  m_rootFilePtr = NULL;
50  initializeCounters(0);
51 }
52 
53 // ===================================================== INITIALIZE ===============================================================
54 void PhaseSpaceAnalysisModule::initialize()
55 {
56  B2INFO("PhaseSpaceAnalysis ------------------------ initialize --------------------------");
57 
58  StoreArray<MCParticle> mcParticles;
59  mcParticles.isRequired();
60 
61  size_t nNames = m_PARAMcontainerNames.size();
62  size_t nTypes = m_PARAMtrackCandTypes.size();
63 
64  if (nTypes > nNames) { B2FATAL("Passed " << nTypes << " trackCandTypese but only " << nNames << " containerNames!"); }
65 
66  for (size_t iName = 0; iName < nNames; ++iName) {
67  size_t iType = iName < nTypes ? iName : nTypes - 1; // only acces values that are possible
68  string tcType = m_PARAMtrackCandTypes.at(iType);
69  if (tcType.compare(string("GFTC")) != 0 && tcType.compare(string("SPTC")) != 0) {
70  B2FATAL("Found id " << tcType << " in 'trackCandTypes' but only 'GFTC' and 'SPTC' are allowed!");
71  }
72 
73  string contName = m_PARAMcontainerNames.at(iName);
74  if (tcType.compare(string("GFTC")) == 0) {
75  StoreArray<genfit::TrackCand> trackCands(contName);
76  trackCands.isRequired();
77  m_tcStoreArrays.push_back(make_pair(trackCands, c_gftc));
78  } else {
79  StoreArray<SpacePointTrackCand> trackCands(contName);
80  trackCands.isRequired();
81  m_tcStoreArrays.push_back(make_pair(trackCands, c_sptc));
82  }
83  if (m_PARAMrootFileName.size() != 2 || (m_PARAMrootFileName[1] != "UPDATE" && m_PARAMrootFileName[1] != "RECREATE")) {
84  string output;
85  // cppcheck-suppress useStlAlgorithm
86  for (string id : m_PARAMrootFileName) { output += "'" + id + "' "; }
87  B2FATAL("PhaseSpaceAnalysis::initialize() : rootFileName is set wrong: entries are: " << output);
88  }
89 
90  m_treeNames.push_back(contName);
91  if (m_PARAMdiffAnalysis && iName != 0) m_treeNames.push_back(contName + "_diff");
92  }
93 
94  initializeRootFile(m_PARAMrootFileName[0] + string(".root"), m_PARAMrootFileName[1], m_treeNames);
95  initializeCounters(m_treeNames.size());
96 }
97 
98 
99 // ============================================================= EVENT ============================================================
100 void PhaseSpaceAnalysisModule::event()
101 {
102  // print out the number of the event in debug output
103  StoreObjPtr<EventMetaData> eventMetaDataPtr("EventMetaData", DataStore::c_Event);
104  const int eventCounter = eventMetaDataPtr->getEvent();
105  string arrayNames;
106  // cppcheck-suppress useStlAlgorithm
107  for (string name : m_PARAMcontainerNames) { arrayNames += " " + name; }
108  B2DEBUG(10, "PhaseSpaceAnalysis::event(). Processing event " << eventCounter << " for StoreArray names :" << arrayNames);
109 
110  StoreArray<MCParticle> mcParticles;
111  const int nMCParticles = mcParticles.getEntries();
112  B2DEBUG(25, "Found " << nMCParticles << " MCParticles for this event");
113 
114  vector<vector<int> > mcPartIds; // collect all mcParticle Ids of all trackCands
115 
116  for (pair<boost::any, e_trackCandType> storeArray : m_tcStoreArrays) {
117  // COULDDO: wrap this in try-catch, but since it is known what has to be casted, nothing should really happen
118  if (storeArray.second == c_gftc) {
119  mcPartIds.push_back(getMCParticleIDs(boost::any_cast<StoreArray<genfit::TrackCand> >(storeArray.first)));
120  } else if (storeArray.second == c_sptc) {
121  mcPartIds.push_back(getMCParticleIDs(boost::any_cast<StoreArray<SpacePointTrackCand> >(storeArray.first)));
122  }
123  B2DEBUG(25, mcPartIds.back().size() << " MCParticles of this event were set for this container");
124  }
125 
126  if (m_PARAMdiffAnalysis) { mcPartIds = getDiffIds(mcPartIds); }
127  // TODO: get the referee Status of a SPTC into the getValuesForRoot stuff
128 
129  // now loop over all MCParticle IDs
130  for (size_t iTree = 0; iTree < mcPartIds.size(); ++iTree) {
131  B2DEBUG(25, "Now collecting values for vector of MCParticle Ids " << iTree);
132  m_mcPartCtr[iTree]++; // increase counter
133 
134  m_rootVariables = RootVariables(); // clear root variables for each run
135  for (int id : mcPartIds[iTree]) {
136  if (id < 0) {
137  B2WARNING("Found a negative id in mcParticleId: " << id << \
138  ". It seems that it has not been set properly, I will skip this MC Particle");
139  m_skippedTCsCtr++;
140  continue;
141  }
142  MCParticle* mcParticle = mcParticles[id];
143  if (mcParticle == NULL) { // safety measure
144  m_noMcPartCtr++;
145  continue;
146  }
147  getValuesForRoot(mcParticle, m_rootVariables);
148  }
149  m_treePtrs[iTree]->Fill(); // write values to TTree
150  }
151 }
152 
153 // ================================================================= TERMINATE ====================================================
154 // TODO: update to new version and new counters
155 void PhaseSpaceAnalysisModule::terminate()
156 {
157  stringstream furtherInfo;
158  if (m_skippedTCsCtr || m_noMcPartCtr) {
159  furtherInfo << " There were " << m_skippedTCsCtr << " negative mcParticle IDs and " << m_noMcPartCtr <<
160  " NULL pointers to MCParticles";
161  }
162  unsigned int nMCParts = accumulate(m_mcPartCtr.begin(), m_mcPartCtr.end(), 0);
163  B2INFO("PhaseSpaceAnalysis::terminate(): Collected mcParticle info in " << m_PARAMcontainerNames.size() << " containers."\
164  " Collected information from " << nMCParts << " MCParticles and wrote them into " << m_treeNames.size() << " different trees."\
165  << furtherInfo.str());
166 
167  if (LogSystem::Instance().isLevelEnabled(LogConfig::c_Debug, 1, PACKAGENAME())) {
168  stringstream output;
169  output << "tree-wise summary (no of mcParticles):";
170  for (unsigned int ctr : m_mcPartCtr) output << " " << ctr;
171  B2DEBUG(1, output.str());
172  }
173 
174  // do ROOT stuff
175  if (m_rootFilePtr != NULL) {
176  m_rootFilePtr->cd(); //important! without this the famework root I/O could mix with the root I/O of this module
177  for (TTree* tree : m_treePtrs) { tree->Write(); }
178  m_rootFilePtr->Close();
179  }
180 }
181 
182 // ================================== initialize root ===============================
183 void PhaseSpaceAnalysisModule::initializeRootFile(std::string fileName, std::string writeOption, std::vector<std::string> treeNames)
184 {
185  B2DEBUG(10, "initializing root file. fileName: " << fileName << ", writeOption: " << writeOption);
186  m_rootFilePtr = new TFile(fileName.c_str(), writeOption.c_str());
187 
188  for (size_t i = 0; i < treeNames.size(); ++i) {
189  TTree* tree = new TTree(treeNames[i].c_str(), "data of mc particles");
190  tree->Branch("momentum_X", &m_rootVariables.MomX);
191  tree->Branch("momentum_Y", &m_rootVariables.MomY);
192  tree->Branch("momentum_Z", &m_rootVariables.MomZ);
193 
194  tree->Branch("eta", &m_rootVariables.Eta);
195  tree->Branch("p_T", &m_rootVariables.pT);
196 
197  tree->Branch("vertex_X", &m_rootVariables.VertX);
198  tree->Branch("vertex_Y", &m_rootVariables.VertY);
199  tree->Branch("vertex_Z", &m_rootVariables.VertZ);
200 
201  tree->Branch("energy", &m_rootVariables.Energy);
202  tree->Branch("mass", &m_rootVariables.Mass);
203 
204  tree->Branch("charge", &m_rootVariables.Charge);
205  tree->Branch("pdg_code", &m_rootVariables.PDG);
206 
207  m_treePtrs.push_back(tree);
208  }
209 }
210 
211 // =============================== collect values for root ====================================
212 void PhaseSpaceAnalysisModule::getValuesForRoot(Belle2::MCParticle* mcParticle, RootVariables& rootVariables)
213 {
214  B2DEBUG(100, "Collecting values for MCParticle " << mcParticle->getArrayIndex());
215  // collect all the momentum
216  const TVector3 momentum = mcParticle->getMomentum();
217  rootVariables.MomX.push_back(momentum.X());
218  rootVariables.MomY.push_back(momentum.Y());
219  rootVariables.MomZ.push_back(momentum.Z());
220 
221  rootVariables.pT.push_back(momentum.Pt());
222  rootVariables.Eta.push_back(momentum.Eta());
223 
224  B2DEBUG(499, "TVector3 momentum: (" << momentum.X() << "," << momentum.Y() << "," << momentum.Z() << \
225  "). This leads to p_T = " << momentum.Pt() << " and eta = " << momentum.Eta());
226 
227  const TVector3 vertex = mcParticle->getVertex();
228  rootVariables.VertX.push_back(vertex.Y());
229  rootVariables.VertY.push_back(vertex.Y());
230  rootVariables.VertZ.push_back(vertex.Z());
231 
232  B2DEBUG(499, "TVector3 vertex: (" << vertex.X() << "," << vertex.Y() << "," << vertex.Z() << ")");
233 
234  rootVariables.Charge.push_back(mcParticle->getCharge());
235  rootVariables.Energy.push_back(mcParticle->getEnergy());
236  rootVariables.PDG.push_back(mcParticle->getPDG());
237  rootVariables.Mass.push_back(mcParticle->getMass());
238 
239  B2DEBUG(499, "Charge " << mcParticle->getCharge() << ", PDG code" << mcParticle->getPDG() << ", Energy " << \
240  mcParticle->getEnergy() << ", Mass " << mcParticle->getMass());
241 }
242 
243 // ==================================================== GET MCPARTICLE IDS ========================================================
244 template<typename TrackCandType>
245 std::vector<int> PhaseSpaceAnalysisModule::getMCParticleIDs(Belle2::StoreArray<TrackCandType> trackCands)
246 {
247  const int nTrackCands = trackCands.getEntries();
248  B2DEBUG(25, "In StoreArray " << trackCands.getName() << ": " << nTrackCands << " track candidates for this event");
249  std::vector<int> mcPartIds;
250  for (int i = 0; i < nTrackCands; ++i) mcPartIds.push_back(trackCands[i]->getMcTrackId());
251  // COULDDO: more debug output
252  return mcPartIds;
253 }
254 
255 // ================================================== GET DIFF IDS ================================================================
256 // TODO: test
257 std::vector<std::vector<int> > PhaseSpaceAnalysisModule::getDiffIds(const std::vector<std::vector<int> >& allIDs)
258 {
259  if (allIDs.size() < 2) { // do nothing if there are less than two entries
260  B2DEBUG(50, "There are no more than 1 vectors passed to getDiffIds. No comparison possible. Returning passed vector");
261  return allIDs;
262  }
263 
264 // // COULDDO: remove (left here for checking if working properly, simply prints out all the entries in the vector of vectors)
265 // if(LogSystem::Instance().isLevelEnabled(LogConfig::c_Debug, 1999, PACKAGENAME())) { // very verbose dev output
266 // stringstream output;
267 // output << "content of allIDs: " << endl;
268 // for(vector<int> vec: allIDs) {
269 // output << "vector: ";
270 // for(int i : vec) output << i << " ";
271 // output << endl;
272 // }
273 // B2DEBUG(1999, output.str())
274 // }
275 
276  vector<vector<int> > diffIds;
277  vector<int> referenceIds = allIDs.at(0); // the first vector is the reference
278  diffIds.push_back(referenceIds);
279 
280  // loop over other vectors and look for missing IDs
281  for (size_t iVec = 1; iVec < allIDs.size(); ++iVec) {
282  vector<int> compareIds = allIDs[iVec];
283  diffIds.push_back(compareIds);
284  vector<int> notFoundIds;
285  for (int id : referenceIds) {
286  if (find(compareIds.begin(), compareIds.end(), id) == compareIds.end()) { // not found
287  notFoundIds.push_back(id);
288  }
289  }
290  diffIds.push_back(notFoundIds); // push_back empty vectors too! otherwise the filling of the TTrees will be wrong!
291  }
292 
293  return diffIds;
294 }
295 
296 // ================================================== GET REFEREE STATUSES ========================================================
297 // std::vector<unsigned short int>
298 // PhaseSpaceAnalysisModule::getRefereeStatuses(Belle2::StoreArray<Belle2::SpacePointTrackCand> trackCands)
299 // {
300 // std::vector<unsigned short int> statuses;
301 // for(int i = 0; i < trackCands.getEntries(); ++i) { statuses.push_back(trackCands[i]->getRefereeStatus()); }
302 // return statuses;
303 // }
Belle2::MCParticle::getEnergy
float getEnergy() const
Return particle energy in GeV.
Definition: MCParticle.h:158
Belle2::MCParticle::getCharge
float getCharge() const
Return the particle charge defined in TDatabasePDG.
Definition: MCParticle.cc:35
Belle2::PhaseSpaceAnalysisModule::RootVariables::VertZ
std::vector< double > VertZ
z position of vertex
Definition: PhaseSpaceAnalysisModule.h:66
REG_MODULE
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:652
Belle2::PhaseSpaceAnalysisModule::RootVariables
helper class to have all RootVariables assembled in one container
Definition: PhaseSpaceAnalysisModule.h:54
Belle2::PhaseSpaceAnalysisModule::RootVariables::MomY
std::vector< double > MomY
y momentum
Definition: PhaseSpaceAnalysisModule.h:56
Belle2::PhaseSpaceAnalysisModule::RootVariables::MomX
std::vector< double > MomX
x momentum
Definition: PhaseSpaceAnalysisModule.h:55
Belle2::MCParticle::getArrayIndex
int getArrayIndex() const
Get 0-based index of the particle in the corresponding MCParticle list.
Definition: MCParticle.h:255
Belle2::PhaseSpaceAnalysisModule::RootVariables::Eta
std::vector< double > Eta
pseudo rapidity.
Definition: PhaseSpaceAnalysisModule.h:60
Belle2::PhaseSpaceAnalysisModule::RootVariables::Mass
std::vector< double > Mass
mc particle mass
Definition: PhaseSpaceAnalysisModule.h:72
Belle2::Module
Base class for Modules.
Definition: Module.h:74
Belle2::PhaseSpaceAnalysisModule::RootVariables::pT
std::vector< double > pT
transverse momentum.
Definition: PhaseSpaceAnalysisModule.h:62
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::StoreObjPtr
Type-safe access to single objects in the data store.
Definition: ParticleList.h:33
Belle2::MCParticle::getVertex
TVector3 getVertex() const
Return production vertex position, shorthand for getProductionVertex().
Definition: MCParticle.h:194
Belle2::PhaseSpaceAnalysisModule::RootVariables::MomZ
std::vector< double > MomZ
z momentum
Definition: PhaseSpaceAnalysisModule.h:57
Belle2::PhaseSpaceAnalysisModule::RootVariables::Energy
std::vector< double > Energy
mc particle energy
Definition: PhaseSpaceAnalysisModule.h:71
Belle2::MCParticle::getPDG
int getPDG() const
Return PDG code of particle.
Definition: MCParticle.h:123
Belle2::MCParticle::getMass
float getMass() const
Return the particle mass in GeV.
Definition: MCParticle.h:146
Belle2::PhaseSpaceAnalysisModule::RootVariables::Charge
std::vector< double > Charge
mc particle charge
Definition: PhaseSpaceAnalysisModule.h:69
Belle2::PhaseSpaceAnalysisModule::RootVariables::VertY
std::vector< double > VertY
y position of vertex
Definition: PhaseSpaceAnalysisModule.h:65
Belle2::MCParticle::getMomentum
TVector3 getMomentum() const
Return momentum.
Definition: MCParticle.h:209
Belle2::MCParticle
A Class to store the Monte Carlo particle information.
Definition: MCParticle.h:43
Belle2::StoreArray< MCParticle >
Belle2::PhaseSpaceAnalysisModule::RootVariables::PDG
std::vector< int > PDG
mc particle pdg code
Definition: PhaseSpaceAnalysisModule.h:68
Belle2::StoreArray::getEntries
int getEntries() const
Get the number of objects in the array.
Definition: StoreArray.h:226
Belle2::PhaseSpaceAnalysisModule
Module for analysing the phase space of genfit::TrackCand(s) and SpacePointTrackCand(s) NOTE: this is...
Definition: PhaseSpaceAnalysisModule.h:44
Belle2::PhaseSpaceAnalysisModule::RootVariables::VertX
std::vector< double > VertX
x position of vertex
Definition: PhaseSpaceAnalysisModule.h:64