Belle II Software  release-06-02-00
SegmentNetworkAnalyzerModule.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 <vector>
10 
11 #include <tracking/modules/vxdtfRedesign/SegmentNetworkAnalyzerModule.h>
12 #include <tracking/spacePointCreation/SpacePoint.h>
13 #include <tracking/spacePointCreation/PurityCalculatorTools.h>
14 
15 
16 using namespace Belle2;
17 using namespace std;
18 
19 REG_MODULE(SegmentNetworkAnalyzer);
20 
22  Module(),
23  m_rFilePtr(nullptr),
24  m_treePtr(nullptr)
25 {
26  setDescription("Module for analyzing the SegmentNetwork.");
27 
28  addParam("networkInputName", m_PARAMnetworkName, "name of the StoreObjPtr to the DNN Container");
29  addParam("rootFileName", m_PARAMrootFileName, "name of the output root file name",
30  std::string("SegmentNetworkAnalyzer_output.root"));
31 
32 }
33 
35 {
36  B2INFO("SegmentNetworkAnalyzer::initialize() ------------------------------");
40 
41  m_rFilePtr = new TFile(m_PARAMrootFileName.c_str(), "recreate");
42  m_treePtr = new TTree("SegmentNetworkAnalyzer", "segment network analysis output");
43 
44  makeBranches();
45 }
46 
48 {
49  m_rootVariables = RootVariables(); // clear rootVariables
50 
51  auto& hitNetwork = m_network->accessHitNetwork();
52  auto& segmentNetwork = m_network->accessSegmentNetwork();
53  m_rootVariables.networkSize = segmentNetwork.size();
55 
56  for (const auto& outerHit : hitNetwork.getNodes()) {
57  for (const auto& centerHit : outerHit->getInnerNodes()) {
58  for (const auto& innerHit : centerHit->getInnerNodes()) {
59 
60  Segment<TrackNode>* innerSegment = new Segment<TrackNode>(centerHit->getEntry().m_sector->getFullSecID(),
61  innerHit->getEntry().m_sector->getFullSecID(),
62  &centerHit->getEntry(),
63  &innerHit->getEntry());
64  Segment<TrackNode>* outerSegment = new Segment<TrackNode>(outerHit->getEntry().m_sector->getFullSecID(),
65  centerHit->getEntry().m_sector->getFullSecID(),
66  &outerHit->getEntry(),
67  &centerHit->getEntry());
68 
69  bool passed = false;
70  // check if the outerSegment is in the network and then look if the inner is connected to it
71  if (auto outerNode = segmentNetwork.getNode(outerSegment->getID())) {
72  for (const auto& connectedNode : outerNode->getInnerNodes()) {
73  if (connectedNode->getEntry() == *innerSegment) {
74  passed = true;
75  break;
76  }
77  }
78  }
79  analyzeCombination(*innerSegment, *outerSegment, passed);
80 
81  } // end inner loop
82  } // end center loop
83  } // end outer loop
84 
85  m_treePtr->Fill();
86 }
87 
89 {
90  // write and close root file
91  if (m_rFilePtr && m_treePtr) {
92  m_rFilePtr->cd();
93  m_rFilePtr->Write();
94  m_rFilePtr->Close();
95  }
96 }
97 
99 {
100  m_treePtr->Branch("phi", &m_rootVariables.phi);
101  m_treePtr->Branch("theta", &m_rootVariables.theta);
102  m_treePtr->Branch("pT", &m_rootVariables.pT);
103  m_treePtr->Branch("signal", &m_rootVariables.signal);
104  m_treePtr->Branch("passed", &m_rootVariables.passed);
105  m_treePtr->Branch("pdg", &m_rootVariables.pdg);
106  m_treePtr->Branch("virtualIP", &m_rootVariables.virtualIP);
107  m_treePtr->Branch("networkSize", &m_rootVariables.networkSize);
108  m_treePtr->Branch("networkConns", &m_rootVariables.networkConnections);
109 }
110 
113  bool passed)
114 {
115  m_rootVariables.passed.push_back(passed);
116 
117  std::vector<const Belle2::SpacePoint*> combinationSPs; // order in which they are stored does not really matter
118  combinationSPs.push_back(outer.getOuterHit()->m_spacePoint);
119  combinationSPs.push_back(outer.getInnerHit()->m_spacePoint);
120  combinationSPs.push_back(inner.getInnerHit()->m_spacePoint);
121 
122  const vector<MCVXDPurityInfo> purityInfo = createPurityInfosVec(combinationSPs);
123  auto mcId = purityInfo[0].getPurity();
124  bool signal = mcId.first >= 0 && mcId.second == 1;
125  m_rootVariables.signal.push_back((int) signal);
126 
127  // define some default values for non-signal combinations
128  int pdg = 0;
129  double pT = -999;
130 
131  if (signal) {
132  const MCParticle* part = m_mcParticles[mcId.first];
133  pdg = part->getPDG();
134  pT = part->getMomentum().Pt();
135  }
136 
137  m_rootVariables.pdg.push_back(pdg);
138  m_rootVariables.pT.push_back(pT);
139 
140  // check if the innermost hit is virtualIP, if so get the outer hit to retrieve the position information from it
141  auto spacePoint = inner.getInnerHit()->m_spacePoint;
142  if (spacePoint->getType() == VXD::SensorInfoBase::VXD) {
143  m_rootVariables.virtualIP.push_back(1);
144  spacePoint = inner.getOuterHit()->m_spacePoint;
145  } else m_rootVariables.virtualIP.push_back(0);
146 
147  // angles are determined from inner most (accessible) hit in combination
148  auto position = spacePoint->getPosition();
149  m_rootVariables.phi.push_back(position.Phi());
150  m_rootVariables.theta.push_back(position.Theta());
151 
152  B2DEBUG(100, "Collected combination with: phi " << position.Phi() << ", theta " << position.Theta() <<
153  ", pdg " << pdg << ", pT " << pT << ", signal " << signal << ", passed " << passed << ", mcId " << mcId.first);
154 }
155 
156 template<typename EntryType, typename MetaInfoType >
158 {
159  size_t nLinks{};
160  for (const auto& outerNodes : network) {
161  nLinks += outerNodes->getInnerNodes().size();
162  }
163  return nLinks;
164 }
DirectedNodeNetwork< Belle2::TrackNode, Belle2::VoidMetaInfo > & accessHitNetwork()
Returns reference to the HitNetwork stored in this container, intended for read and write access.
DirectedNodeNetwork< Belle2::Segment< Belle2::TrackNode >, Belle2::CACell > & accessSegmentNetwork()
Returns reference to the SegmentNetwork stored in this container, intended for read and write access.
Network of directed nodes of the type EntryType.
A Class to store the Monte Carlo particle information.
Definition: MCParticle.h:32
Base class for Modules.
Definition: Module.h:72
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
std::string m_PARAMrootFileName
file name of the produced root file
void initialize() override
set up root file and check required Store Arrays
Belle2::StoreObjPtr< Belle2::DirectedNodeNetworkContainer > m_network
StoreObjPtr to the SegmentNetwork and TrackNode Network container.
RootVariables m_rootVariables
handle to collect all data for one event
void event() override
collect necessary data and put into TTree
std::string m_PARAMnetworkName
StoreArray name of the DirectedNodeNetworkContainer.
size_t getNConnections(Belle2::DirectedNodeNetwork< EntryType, MetaInfoType > &network) const
get the number of connections between the nodes of a network
void terminate() override
write and close root file
Belle2::StoreArray< Belle2::MCParticle > m_mcParticles
MCParticles StoreArray for obtaining MC Information.
void analyzeCombination(const Belle2::Segment< Belle2::TrackNode > &outer, const Belle2::Segment< Belle2::TrackNode > &inner, bool passed)
get necessary data from three hit combination and put them into the root variables
void makeBranches()
setup the Branches in the output TTree
The Segment class This class represents segments of track candidates needed for TrackFinderVXD-Module...
Definition: Segment.h:25
std::int64_t getID() const
************************* PUBLIC MEMBER FUNCTIONS *************************
Definition: Segment.h:84
const HitType * getInnerHit() const
returns inner hit of current Segment
Definition: Segment.h:96
const HitType * getOuterHit() const
returns outer hit of current Segment
Definition: Segment.h:100
bool isRequired(const std::string &name="")
Ensure this array/object has been registered previously.
@ VXD
Any type of VXD Sensor.
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 std::vector< Belle2::MCVXDPurityInfo > createPurityInfosVec(const std::vector< const Belle2::SpacePoint * > &spacePoints)
create a vector of MCVXDPurityInfos objects for a std::vector<Belle2::SpacePoints>.
Abstract base class for different kinds of events.
keep all the variables for rootoutput in one struct
std::vector< double > phi
phi of the innermost hit (that is not the virtual IP)
std::vector< int > pdg
pdg of the related MCParticle
std::vector< int > virtualIP
did the segment contain the virtual IP
std::vector< double > theta
theta of the innermost hit (that is not the virtual IP)
unsigned networkConnections
number of connections in network
std::vector< int > passed
did segment combination pass the three hit filters
std::vector< int > signal
was segment combination signal
std::vector< double > pT
pT of the related MCParticle
SpacePoint * m_spacePoint
Pointer to spacePoint.
Definition: TrackNode.h:42