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