Belle II Software  release-06-02-00
StudyMaterialEffectsModule.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 #include <tracking/modules/trackFinderVXDTests/StudyMaterialEffectsModule.h>
9 
10 #include <framework/pcore/ProcHandler.h>
11 #include <framework/logging/Logger.h>
12 
13 #include <mdst/dataobjects/MCParticle.h>
14 #include <svd/dataobjects/SVDCluster.h>
15 
16 #include <vxd/geometry/SensorInfoBase.h>
17 
18 #include <TVector3.h>
19 #include "TMath.h"
20 
21 
22 using namespace std;
23 using namespace Belle2;
24 
25 
26 REG_MODULE(StudyMaterialEffects)
27 
28 StudyMaterialEffectsModule::StudyMaterialEffectsModule() : Module(), m_tree("materialEffectsStudyOutput" , DataStore::c_Persistent)
29 {
30  setDescription("StudyMaterialEffects- should be used with single track pGuns and without magnetic field.");
31 
32  setPropertyFlags(c_ParallelProcessingCertified | c_TerminateInAllProcesses);
33 }
34 
35 
36 void StudyMaterialEffectsModule::initialize()
37 {
38  m_spacePoints.isRequired("spacePoints");
39 
40  m_file = new TFile("materialStudy.root", "RECREATE");
41 
42  m_file->cd();
43 
44  m_tree.registerInDataStore();
45  m_tree.construct("materialEffectsStudyTree", "Raw data of two-hit-combinations for a sectorMap");
46 
47  m_tree->get().Branch("ResidualPhiL3L3", &(m_PhiL3L3));
48  m_tree->get().Branch("ResidualPhiL3L4", &(m_PhiL3L4));
49  m_tree->get().Branch("ResidualPhiL3L6", &(m_PhiL3L6));
50  m_tree->get().Branch("ResidualThetaL3L3", &(m_ThetaL3L3));
51  m_tree->get().Branch("ResidualThetaL3L4", &(m_ThetaL3L4));
52  m_tree->get().Branch("ResidualThetaL3L6", &(m_ThetaL3L6));
53  m_tree->get().Branch("ResidualScatterAngleL3L3", &(m_ScatterAngleL3L3));
54  m_tree->get().Branch("ResidualScatterAngleL3L4", &(m_ScatterAngleL3L4));
55  m_tree->get().Branch("ResidualcatterAngleL3L6", &(m_ScatterAngleL3L6));
56  m_tree->get().Branch("ResidualScatterAngleV2GradL3L3", &(m_ScatterAngleGradL3L3));
57  m_tree->get().Branch("ResidualScatterAngleV2GradL3L4", &(m_ScatterAngleGradL3L4));
58  m_tree->get().Branch("ResidualScatterAngleV2GradL3L6", &(m_ScatterAngleGradL3L6));
59  m_tree->get().Branch("ResidualScatterAngleV3GradL3L3", &(m_ScatterAngleV3GradL3L3));
60  m_tree->get().Branch("ResidualScatterAngleV3GradL3L4", &(m_ScatterAngleV3GradL3L4));
61  m_tree->get().Branch("ResidualScatterAngleV3GradL3L6", &(m_ScatterAngleV3GradL3L6));
62  m_tree->get().Branch("ResidualcatterAngleL3L6", &(m_ScatterAngleL3L6));
63  m_tree->get().Branch("ResidualXYL3L4", &(m_distXY));
64  m_tree->get().Branch("ResidualMomentumL3bL3e", &(m_deltaPL3L3));
65  m_tree->get().Branch("ResidualMomentumL3bL4b", &(m_deltaPL3L4));
66  m_tree->get().Branch("ResidualMomentumL3bL6e", &(m_deltaPL3L6));
67 
68  B2WARNING("StudyMaterialEffectsModule::initialize: nBranches: " << m_tree->get().GetNbranches());
69 }
70 
71 
72 
73 void StudyMaterialEffectsModule::event()
74 {
75  // position and momentum in global coordinates:
76  B2Vector3D l3HitPosBegin;
77  B2Vector3D l3HitPosEnd;
78  B2Vector3D l4HitPosBegin;
79  B2Vector3D l6HitPosEnd;
80  B2Vector3D l3MomentumBegin;
81  B2Vector3D l3MomentumEnd;
82  B2Vector3D l4MomentumBegin;
83  B2Vector3D l6MomentumEnd;
84 
85  bool wasFoundL3(false);
86  bool wasFoundL4(false);
87  bool wasFoundL6(false);
88 
89 
90  for (const auto& aSP : m_spacePoints) {
91  VxdID vxdID = aSP.getVxdID();
92  const SVDTrueHit* trueHit = getTrueHit(aSP);
93  if (trueHit == nullptr) continue;
94 
95  B2Vector3D entryHitPos = getGlobalPosition(trueHit, vxdID, true);
96  B2Vector3D exitHitPos = getGlobalPosition(trueHit, vxdID, false);
97  B2Vector3D entryMomentum = getGlobalMomentumVector(trueHit, vxdID, true);
98  B2Vector3D exitMomentum = getGlobalMomentumVector(trueHit, vxdID, false);
99  if (vxdID.getLayerNumber() == 3) {
100  l3HitPosBegin = entryHitPos;
101  l3HitPosEnd = exitHitPos;
102  l3MomentumBegin = entryMomentum;
103  l3MomentumEnd = exitMomentum;
104  wasFoundL3 = true;
105  }
106  if (vxdID.getLayerNumber() == 4) {
107  l4HitPosBegin = entryHitPos;
108  l4MomentumBegin = entryMomentum;
109  wasFoundL4 = true;
110  }
111  if (vxdID.getLayerNumber() == 6) {
112  l6HitPosEnd = exitHitPos;
113  l6MomentumEnd = exitMomentum;
114  wasFoundL6 = true;
115  }
116  }
117 
118  if (!wasFoundL3 or !wasFoundL4 or !wasFoundL6) return;
119  m_COUNTERsuccessfullEvents++;
120 
121  m_PhiL3L3 = l3HitPosBegin.Phi() - l3HitPosEnd.Phi();
122  m_PhiL3L4 = l3HitPosBegin.Phi() - l4HitPosBegin.Phi();
123  m_PhiL3L6 = l3HitPosBegin.Phi() - l6HitPosEnd.Phi();
124 
125  m_ThetaL3L3 = l3HitPosBegin.Theta() - l3HitPosEnd.Theta();
126  m_ThetaL3L4 = l3HitPosBegin.Theta() - l4HitPosBegin.Theta();
127  m_ThetaL3L6 = l3HitPosBegin.Theta() - l6HitPosEnd.Theta();
128 
129  m_ScatterAngleL3L3 = sqrt(pow(m_ThetaL3L3, 2) + pow(m_PhiL3L3, 2));
130  m_ScatterAngleL3L4 = sqrt(pow(m_ThetaL3L4, 2) + pow(m_PhiL3L4, 2));
131  m_ScatterAngleL3L6 = sqrt(pow(m_ThetaL3L6, 2) + pow(m_PhiL3L6, 2));
132 
133  m_ScatterAngleGradL3L3 = (l3HitPosEnd - l3HitPosBegin).Angle(l3HitPosBegin) * 180. / TMath::Pi();
134  m_ScatterAngleGradL3L4 = (l4HitPosBegin - l3HitPosBegin).Angle(l3HitPosBegin) * 180. / TMath::Pi();
135  m_ScatterAngleGradL3L6 = (l6HitPosEnd - l3HitPosBegin).Angle(l3HitPosBegin) * 180. / TMath::Pi();
136 
137 
138  m_ScatterAngleV3GradL3L3 = l3MomentumEnd.Angle(l3MomentumBegin) * 180. / TMath::Pi();
139  m_ScatterAngleV3GradL3L4 = l4MomentumBegin.Angle(l3MomentumBegin) * 180. / TMath::Pi();
140  m_ScatterAngleV3GradL3L6 = l6MomentumEnd.Angle(l3MomentumBegin) * 180. / TMath::Pi();
141 
142 
143  m_distXY = (l3HitPosBegin - l4HitPosBegin).Perp();
144 
145  m_deltaPL3L3 = l3MomentumBegin.Mag() - l3MomentumEnd.Mag();
146  m_deltaPL3L4 = l3MomentumBegin.Mag() - l4MomentumBegin.Mag();
147  m_deltaPL3L6 = l3MomentumBegin.Mag() - l6MomentumEnd.Mag();
148 
149  m_tree->get().Fill();
150 }
151 
152 
153 
154 void StudyMaterialEffectsModule::terminate()
155 {
156  B2WARNING("StudyMaterialEffects: there were " << m_COUNTERsuccessfullEvents << " events with hits resulting in root output");
157  m_file->cd();
158 
159  if (!ProcHandler::parallelProcessingUsed() or ProcHandler::isOutputProcess()) {
160  //use TFile you created in initialize()
161  m_tree->write(m_file);
162  }
163 }
164 
165 
166 
167 B2Vector3D StudyMaterialEffectsModule::getGlobalPosition(const SVDTrueHit* trueHit, VxdID vxdID, bool useEntry)
168 {
169  const Belle2::VXD::SensorInfoBase* aSensorInfo = &VXD::GeoCache::getInstance().getSensorInfo(vxdID);
170 
171  if (useEntry) {
172  B2Vector3D pos = aSensorInfo->pointToGlobal(TVector3(trueHit->getEntryU(), trueHit->getEntryV(), 0), true);
173  return pos;
174  }
175  B2Vector3D pos = aSensorInfo->pointToGlobal(TVector3(trueHit->getExitU(), trueHit->getExitV(), 0), true);
176  return pos;
177 }
178 
179 
180 
181 const SVDTrueHit* StudyMaterialEffectsModule::getTrueHit(const SpacePoint& aSP)
182 {
183 // const auto* aTrueHit = aSP.getRelationsTo<SVDCluster>("ALL")[0]->getRelationsTo<SVDTrueHit>("ALL")[0];;
184  const SVDTrueHit* aTrueHit = nullptr;
185  for (const SVDCluster& aCluster : aSP.getRelationsTo<SVDCluster>("ALL")) {
186  for (const MCParticle& aMcParticle : aCluster.getRelationsTo<MCParticle>("ALL")) {
187  if (aMcParticle.hasStatus(MCParticle::c_PrimaryParticle) == false) continue;
188 
189  aTrueHit = aCluster.getRelationsTo<SVDTrueHit>("ALL")[0];
190 // for (const auto* aTrueHit : aCluster.getRelationsTo<SVDTrueHit>("ALL")) {
191 //
192 // }
193  }
194  }
195  return aTrueHit;
196 }
197 
198 
199 B2Vector3D StudyMaterialEffectsModule::getGlobalMomentumVector(const SVDTrueHit* trueHit, VxdID vxdID, bool useEntry)
200 {
201  const Belle2::VXD::SensorInfoBase* aSensorInfo = &VXD::GeoCache::getInstance().getSensorInfo(vxdID);
202 
203  if (useEntry) {
204 // B2Vector3D mom = aSensorInfo->pointToGlobal(trueHit->getEntryMomentum());
205  return aSensorInfo->vectorToGlobal(trueHit->getEntryMomentum(), true);
206  }
207 // B2Vector3D pos = aSensorInfo->pointToGlobal(trueHit->getExitMomentum());
208  return aSensorInfo->vectorToGlobal(trueHit->getExitMomentum(), true);
209 }
DataType Phi() const
The azimuth angle.
Definition: B2Vector3.h:136
DataType Theta() const
The polar angle.
Definition: B2Vector3.h:138
DataType Mag() const
The magnitude (rho in spherical coordinate system).
Definition: B2Vector3.h:144
DataType Angle(const B2Vector3< DataType > &q) const
The angle w.r.t.
Definition: B2Vector3.h:287
In the store you can park objects that have to be accessed by various modules.
Definition: DataStore.h:51
A Class to store the Monte Carlo particle information.
Definition: MCParticle.h:32
Base class for Modules.
Definition: Module.h:72
RelationVector< TO > getRelationsTo(const std::string &name="", const std::string &namedRelation="") const
Get the relations that point from this object to another store array.
The SVD Cluster class This class stores all information about reconstructed SVD clusters.
Definition: SVDCluster.h:28
Class SVDTrueHit - Records of tracks that either enter or leave the sensitive volume.
Definition: SVDTrueHit.h:33
SpacePoint typically is build from 1 PXDCluster or 1-2 SVDClusters.
Definition: SpacePoint.h:42
float getEntryU() const
Return local u coordinate of hit when entering silicon.
Definition: VXDTrueHit.h:81
TVector3 getEntryMomentum() const
Return momentum at the start point of the track.
Definition: VXDTrueHit.h:99
TVector3 getExitMomentum() const
Return momentum at the endpoint of the track.
Definition: VXDTrueHit.h:101
float getExitU() const
Return local u coordinate of hit at the endpoint of the track.
Definition: VXDTrueHit.h:87
float getExitV() const
Return local v coordinate of hit at the endpoint of the track.
Definition: VXDTrueHit.h:89
float getEntryV() const
Return local v coordinate of the start point of the track.
Definition: VXDTrueHit.h:83
Base class to provide Sensor Information for PXD and SVD.
TVector3 pointToGlobal(const TVector3 &local, bool reco=false) const
Convert a point from local to global coordinates.
TVector3 vectorToGlobal(const TVector3 &local, bool reco=false) const
Convert a vector from local to global coordinates.
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:33
baseType getLayerNumber() const
Get the layer id.
Definition: VxdID.h:96
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
Abstract base class for different kinds of events.