Belle II Software development
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
16using namespace Belle2;
17
18REG_MODULE(SegmentNetworkAnalyzer);
19
21 Module(),
22 m_rFilePtr(nullptr),
23 m_treePtr(nullptr)
24{
25 setDescription("Module for analyzing the SegmentNetwork.");
26
27 addParam("networkInputName", m_PARAMnetworkName, "name of the StoreObjPtr to the DNN Container");
28 addParam("rootFileName", m_PARAMrootFileName, "name of the output root file name",
29 std::string("SegmentNetworkAnalyzer_output.root"));
30
31}
32
34{
35 B2INFO("SegmentNetworkAnalyzer::initialize() ------------------------------");
39
40 m_rFilePtr = new TFile(m_PARAMrootFileName.c_str(), "recreate");
41 m_treePtr = new TTree("SegmentNetworkAnalyzer", "segment network analysis output");
42
44}
45
47{
48 m_rootVariables = RootVariables(); // clear rootVariables
49
50 auto& hitNetwork = m_network->accessHitNetwork();
51 auto& segmentNetwork = m_network->accessSegmentNetwork();
52 m_rootVariables.networkSize = segmentNetwork.size();
54
55 for (const auto& outerHit : hitNetwork.getNodes()) {
56 for (const auto& centerHit : outerHit->getInnerNodes()) {
57 Segment<TrackNode> outerSegment(outerHit->getEntry().m_sector->getFullSecID(),
58 centerHit->getEntry().m_sector->getFullSecID(),
59 &outerHit->getEntry(),
60 &centerHit->getEntry());
61
62 for (const auto& innerHit : centerHit->getInnerNodes()) {
63 Segment<TrackNode> innerSegment(centerHit->getEntry().m_sector->getFullSecID(),
64 innerHit->getEntry().m_sector->getFullSecID(),
65 &centerHit->getEntry(),
66 &innerHit->getEntry());
67
68 bool passed = false;
69 // check if the outerSegment is in the network and then look if the inner is connected to it
70 if (auto outerNode = segmentNetwork.getNode(outerSegment.getID())) {
71 for (const auto& connectedNode : outerNode->getInnerNodes()) {
72 if (connectedNode->getEntry() == innerSegment) {
73 passed = true;
74 break;
75 }
76 }
77 }
78 analyzeCombination(innerSegment, outerSegment, passed);
79
80 } // end inner loop
81 } // end center loop
82 } // end outer loop
83
84 m_treePtr->Fill();
85}
86
88{
89 // write and close root file
90 if (m_rFilePtr && m_treePtr) {
91 m_rFilePtr->cd();
92 m_rFilePtr->Write();
93 m_rFilePtr->Close();
94 }
95}
96
98{
99 m_treePtr->Branch("phi", &m_rootVariables.phi);
100 m_treePtr->Branch("theta", &m_rootVariables.theta);
101 m_treePtr->Branch("pT", &m_rootVariables.pT);
102 m_treePtr->Branch("signal", &m_rootVariables.signal);
103 m_treePtr->Branch("passed", &m_rootVariables.passed);
104 m_treePtr->Branch("pdg", &m_rootVariables.pdg);
105 m_treePtr->Branch("virtualIP", &m_rootVariables.virtualIP);
106 m_treePtr->Branch("networkSize", &m_rootVariables.networkSize);
107 m_treePtr->Branch("networkConns", &m_rootVariables.networkConnections);
108}
109
112 bool passed)
113{
114 m_rootVariables.passed.push_back(passed);
115
116 std::vector<const Belle2::SpacePoint*> combinationSPs; // order in which they are stored does not really matter
117 combinationSPs.push_back(outer.getOuterHit()->m_spacePoint);
118 combinationSPs.push_back(outer.getInnerHit()->m_spacePoint);
119 combinationSPs.push_back(inner.getInnerHit()->m_spacePoint);
120
121 const std::vector<MCVXDPurityInfo> purityInfo = createPurityInfosVec(combinationSPs);
122 auto mcId = purityInfo[0].getPurity();
123 bool signal = mcId.first >= 0 && mcId.second == 1;
124 m_rootVariables.signal.push_back((int) signal);
125
126 // define some default values for non-signal combinations
127 int pdg = 0;
128 double pT = -999;
129
130 if (signal) {
131 const MCParticle* part = m_mcParticles[mcId.first];
132 pdg = part->getPDG();
133 pT = part->getMomentum().Rho();
134 }
135
136 m_rootVariables.pdg.push_back(pdg);
137 m_rootVariables.pT.push_back(pT);
138
139 // check if the innermost hit is virtualIP, if so get the outer hit to retrieve the position information from it
140 auto spacePoint = inner.getInnerHit()->m_spacePoint;
141 if (spacePoint->getType() == VXD::SensorInfoBase::VXD) {
142 m_rootVariables.virtualIP.push_back(1);
143 spacePoint = inner.getOuterHit()->m_spacePoint;
144 } else m_rootVariables.virtualIP.push_back(0);
145
146 // angles are determined from inner most (accessible) hit in combination
147 auto position = spacePoint->getPosition();
148 m_rootVariables.phi.push_back(position.Phi());
149 m_rootVariables.theta.push_back(position.Theta());
150
151 B2DEBUG(22, "Collected combination with: phi " << position.Phi() << ", theta " << position.Theta() <<
152 ", pdg " << pdg << ", pT " << pT << ", signal " << signal << ", passed " << passed << ", mcId " << mcId.first);
153}
154
155template<typename EntryType, typename MetaInfoType >
157{
158 size_t nLinks{};
159 for (const auto& outerNodes : network) {
160 nLinks += outerNodes->getInnerNodes().size();
161 }
162 return nLinks;
163}
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
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
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 * getOuterHit() const
returns outer hit of current Segment
Definition: Segment.h:100
const HitType * getInnerHit() const
returns inner hit of current Segment
Definition: Segment.h:96
bool isRequired(const std::string &name="")
Ensure this array/object has been registered previously.
Accessor to arrays stored in the data store.
Definition: StoreArray.h:113
@ 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>.
B2Vector3D outerHit(0, 0, 0)
testing out of range behavior
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