ATLAS Offline Software
Loading...
Searching...
No Matches
MMHitsTestTool.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2023 CERN for the benefit of the ATLAS collaboration
3*/
4
5#include "MMHitsTestTool.h"
6
7#include "Identifier/Identifier.h"
8
10
13
15
18
20#include "CLHEP/Vector/LorentzVector.h"
21
22#include "GaudiKernel/NTuple.h"
23#include "GaudiKernel/SmartDataPtr.h"
24#include "GaudiKernel/IDataProviderSvc.h"
25#include "GaudiKernel/ITHistSvc.h"
26#include "GaudiKernel/INTupleSvc.h"
27
28#include "TH2D.h"
29#include "TTree.h"
30#include "TROOT.h"
31#include "TFile.h"
32#include "TF1.h"
33#include "TH1F.h"
34
35using namespace MuonGM;
36using namespace std;
37
38
41
42 if (m_DoMMTest) {
43 const MMSimHitCollection* p_collection = nullptr;
44 CHECK(evtStore()->retrieve(p_collection,"MM_Hits"));
45 for (const MMSimHit& hit : *p_collection) {
46 Amg::Vector3D u = (hit).globalPosition();
48 //Useful link on how to retrieve variables: http://acode-browser.usatlas.bnl.gov/lxr/source/atlas/MuonSpectrometer/MuonValidation/MuonPRDTest/src/MMSimHitVariables.cxx
49 //Get station names and make plots for each wedge
51 int simId = (hit).MMId();
52 std::string sim_stationName = hitHelper->GetStationName(simId);
53 //Declare station name strings
54 static const std::string s_m1s1("M1S1");
55 static const std::string s_m2s1("M2S1");
56 static const std::string s_m1l1("M1L1");
57 static const std::string s_m2l1("M2L1");
58 static const std::string s_m1s2("M1S2");
59 static const std::string s_m2s2("M2S2");
60 static const std::string s_m1l2("M1L2");
61 static const std::string s_m2l2("M2L2");
62
63 //----------------------------------Wedge 1 Histos begin-------------------------------------------------------------------------
64 //M1S1 (Note: M1->Module 1, S1->Small sector, wedge 1)
65 if (sim_stationName==s_m1s1 && u.z()>0){
66 m_MMTransverseEta1SmallWedge1->Fill(u.x(),u.y());
67 }
68
69 //M2S1 (Note: M2->Module 2, S1->Small sector, wedge 1)
70 if (sim_stationName==s_m2s1 && u.z()>0){
71 m_MMTransverseEta2SmallWedge1->Fill(u.x(),u.y());
72 }
73
74 //M1L1 (Note: M1->Module 1, S1->Large sector, wedge 1)
75 if (sim_stationName==s_m1l1 && u.z()>0){
76 m_MMTransverseEta1LargeWedge1->Fill(u.x(),u.y());
77 }
78
79 //M2L1 (Note: M2->Module 2, S1->Large sector, wedge 1.)
80 if (sim_stationName==s_m2l1 && u.z()>0){
81 m_MMTransverseEta2LargeWedge1->Fill(u.x(),u.y());
82 }
83
84 //Plots of transverse and rZ view for wedge 1 in either +ve Z or -ve Z and done separately for large and small sectors
85 double rwedge1 = sqrt(u.x()*u.x()+u.y()*u.y()); //Evaluate r
86
87 //Small sectors at +ve Z
88 if ((sim_stationName==s_m1s1 || sim_stationName==s_m2s1) && u.z() > 0){
90 m_MM_SmallWedge1_rZview_positiveZ->Fill(u.z(), rwedge1);
91 m_MM_rPlot_S1_posZ->Fill(rwedge1);
92 }
93
94 //Small sectors at -ve Z
95 if ((sim_stationName==s_m1s1 || sim_stationName==s_m2s1) && u.z() < 0){
97 m_MM_SmallWedge1_rZview_negativeZ->Fill(u.z(), rwedge1);
98 }
99
100 //Large sectors at +ve Z
101 if ((sim_stationName==s_m1l1 || sim_stationName==s_m2l1) && u.z() > 0){
103 m_MM_LargeWedge1_rZview_positiveZ->Fill(u.z(), rwedge1);
104 m_MM_rPlot_L1_posZ->Fill(rwedge1);
105 }
106
107 //Large sectors at -ve Z
108 if ((sim_stationName==s_m1l1 || sim_stationName==s_m2l1) && u.z() < 0){
110 m_MM_LargeWedge1_rZview_negativeZ->Fill(u.z(), rwedge1);
111 }
112 //----------------------------------Wedge 1 Histos end-------------------------------------------------------------------------
113
114 //----------------------------------Wedge 2 Histos begin---------------------------------------------------------------------
115 //M1S2 (Note: M1->Module 1, S2->Small sector, wedge 2)
116 if (sim_stationName==s_m1s2 && u.z()>0){
117 m_MMTransverseEta1SmallWedge2->Fill(u.x(),u.y());
118 }
119
120 //M2S2 (Note: M1->Module 2, S2->Small sector, wedge 2)
121 if (sim_stationName==s_m2s2 && u.z()>0){
122 m_MMTransverseEta2SmallWedge2->Fill(u.x(),u.y());
123 }
124
125 //M1L2 (Note: M1->Module 1, L2->Large sector, wedge 2)
126 if (sim_stationName==s_m1l2 && u.z()>0){
127 m_MMTransverseEta1LargeWedge2->Fill(u.x(),u.y());
128 }
129
130 //M2L2 (Note: M2->Module 2, L2->Large sector, wedge 2)
131 if (sim_stationName==s_m2l2 && u.z()>0){
132 m_MMTransverseEta2LargeWedge2->Fill(u.x(),u.y());
133 }
134
135 //Plots of transverse and rZ view for wedge 1 in either +ve Z or -ve Z and done separately for large and small sectors
136 double rwedge2 = sqrt(u.x()*u.x()+u.y()*u.y()); //Evaluate r
137
138 //Small sectors at +ve Z
139 if ((sim_stationName==s_m1s2 || sim_stationName==s_m2s2) && u.z() > 0){
141 m_MM_SmallWedge2_rZview_positiveZ->Fill(u.z(), rwedge2);
142 m_MM_rPlot_S2_posZ->Fill(rwedge2);
143 }
144
145 //Small sectors at -ve Z
146 if ((sim_stationName==s_m1s2 || sim_stationName==s_m2s2) && u.z() < 0){
148 m_MM_SmallWedge2_rZview_negativeZ->Fill(u.z(), rwedge2);
149 }
150
151 //Large sectors at +ve Z
152 if ((sim_stationName==s_m1l2 || sim_stationName==s_m2l2) && u.z() > 0){
154 m_MM_LargeWedge2_rZview_positiveZ->Fill(u.z(), rwedge2);
155 m_MM_rPlot_L2_posZ->Fill(rwedge2);
156 }
157
158 //Large sectors at -ve Z
159 if ((sim_stationName==s_m1l2 || sim_stationName==s_m2l2) && u.z() < 0){
161 m_MM_LargeWedge2_rZview_negativeZ->Fill(u.z(), rwedge2);
162 }
163 //-----------------------------------------------Wedge 2 Histos end-----------------------------------------------------------
164
165
166 // GeoMMHit ghit(hit);
167 // if (!ghit) continue;
168 // Amg::Vector3D u = ghit.getGlobalPosition();
169 // CHECK(executeFillHistos(u));
170 }
171 }
172
173 return StatusCode::SUCCESS;
174}
175
176
179
180 //-------------------------------Wedge 1 begin-------------------------------------------------------------------------------
181 _TH1D(m_MM_rPlot_S1_posZ,"MM_rPlot_S1_posZ",10000,0.,14000.);
182 _TH1D(m_MM_rPlot_L1_posZ,"MM_rPlot_L1_posZ",5000,0.,5000.);
183
184 _TH2D(m_MM_SmallWedge1_rZview_positiveZ,"MM_rZView_S1_posZ",1000,7000.,8000.,5000,0.,5000.);
185 _TH2D(m_MM_SmallWedge1_rZview_negativeZ,"MM_rZView_S1_negZ",1000,-8000.,-7000.,5000,0.,5000.);
186 _TH2D(m_MM_LargeWedge1_rZview_positiveZ,"MM_rZView_L1_posZ",1000,7000.,8000.,5000,0.,5000.);
187 _TH2D(m_MM_LargeWedge1_rZview_negativeZ,"MM_rZView_L1_negZ",1000,-8000.,-7000.,5000,0.,5000.);
188
189 _TH2D(m_MMTransverseEta1SmallWedge1,"MM_TransverseView_M1S1_posZ",1200,-6000.,6000.,1200,-6000.,6000.);
190 _TH2D(m_MMTransverseEta2SmallWedge1,"MM_TransverseView_M2S1_posZ",1200,-6000.,6000.,1200,-6000.,6000.);
191 _TH2D(m_MMTransverseEta1LargeWedge1,"MM_TransverseView_M1L1_posZ",1200,-6000.,6000.,1200,-6000.,6000.);
192 _TH2D(m_MMTransverseEta2LargeWedge1,"MM_TransverseView_M2L1_posZ",1200,-6000.,6000.,1200,-6000.,6000.);
193
194 _TH2D(m_MM_SmallWedge1_TransverseView_positiveZ,"MM_TransverseView_S1_posZ",1200,-6000.,6000.,1200,-6000.,6000.);
195 _TH2D(m_MM_SmallWedge1_TransverseView_negativeZ,"MM_TransverseView_S1_negZ",1200,-6000.,6000.,1200,-6000.,6000.);
196 _TH2D(m_MM_LargeWedge1_TransverseView_positiveZ,"MM_TransverseView_L1_posZ",1200,-6000.,6000.,1200,-6000.,6000.);
197 _TH2D(m_MM_LargeWedge1_TransverseView_negativeZ,"MM_TransverseView_L1_negZ",1200,-6000.,6000.,1200,-6000.,6000.);
198 //--------------------------------Wedge 1 end--------------------------------------------------------------------------------
199
200 //--------------------------------Wedge 2 begin------------------------------------------------------------------------------
201 _TH1D(m_MM_rPlot_S2_posZ,"MM_rPlot_S2_posZ",10000,0.,14000.);
202 _TH1D(m_MM_rPlot_L2_posZ,"MM_rPlot_L2_posZ",10000,0.,14000.);
203
204 _TH2D(m_MM_SmallWedge2_rZview_positiveZ,"MM_rZView_S2_posZ",1000,7000.,8000.,5000,0.,5000.);
205 _TH2D(m_MM_SmallWedge2_rZview_negativeZ,"MM_rZView_S2_negZ",1000,-8000.,-7000.,5000,0.,5000.);
206 _TH2D(m_MM_LargeWedge2_rZview_positiveZ,"MM_rZView_L2_posZ",1000,7000.,8000.,5000,0.,5000.);
207 _TH2D(m_MM_LargeWedge2_rZview_negativeZ,"MM_rZView_L2_negZ",1000,-8000.,-7000.,5000,0.,5000.);
208
209 _TH2D(m_MMTransverseEta1SmallWedge2,"MM_TransverseView_M1S2_posZ",1200,-6000.,6000.,1200,-6000.,6000.);
210 _TH2D(m_MMTransverseEta2SmallWedge2,"MM_TransverseView_M2S2_posZ",1200,-6000.,6000.,1200,-6000.,6000.);
211 _TH2D(m_MMTransverseEta1LargeWedge2,"MM_TransverseView_M1L2_posZ",1200,-6000.,6000.,1200,-6000.,6000.);
212 _TH2D(m_MMTransverseEta2LargeWedge2,"MM_TransverseView_M2L2_posZ",1200,-6000.,6000.,1200,-6000.,6000.);
213
214 _TH2D(m_MM_SmallWedge2_TransverseView_positiveZ,"MM_TransverseView_S2_posZ",1200,-6000.,6000.,1200,-6000.,6000.);
215 _TH2D(m_MM_SmallWedge2_TransverseView_negativeZ,"MM_TransverseView_S2_negZ",1200,-6000.,6000.,1200,-6000.,6000.);
216 _TH2D(m_MM_LargeWedge2_TransverseView_positiveZ,"MM_TransverseView_L2_posZ",1200,-6000.,6000.,1200,-6000.,6000.);
217 _TH2D(m_MM_LargeWedge2_TransverseView_negativeZ,"MM_TransverseView_L2_negZ",1200,-6000.,6000.,1200,-6000.,6000.);
218 //--------------------------------Wedge 2 end---------------------------------------------------------------------------------------
219
220
221 return StatusCode::SUCCESS;
222}
#define CHECK(...)
Evaluate an expression and check for errors.
AtlasHitsVector< MMSimHit > MMSimHitCollection
#define _TH2D(var, name, nbinx, xmin, xmax, nbiny, ymin, ymax)
#define _TH1D(var, name, nbin, xmin, xmax)
TH2 * m_MMTransverseEta1SmallWedge1
TH2 * m_MMTransverseEta2LargeWedge1
TH2 * m_MM_SmallWedge1_rZview_negativeZ
TH2 * m_MM_SmallWedge1_TransverseView_negativeZ
TH2 * m_MMTransverseEta1LargeWedge1
TH2 * m_MM_LargeWedge1_TransverseView_negativeZ
TH2 * m_MM_SmallWedge2_rZview_positiveZ
TH2 * m_MM_LargeWedge1_TransverseView_positiveZ
TH2 * m_MMTransverseEta2SmallWedge2
TH2 * m_MM_LargeWedge1_rZview_positiveZ
TH2 * m_MM_SmallWedge2_TransverseView_negativeZ
TH2 * m_MM_SmallWedge1_rZview_positiveZ
TH2 * m_MM_LargeWedge2_TransverseView_positiveZ
TH2 * m_MMTransverseEta1LargeWedge2
TH2 * m_MM_LargeWedge2_rZview_negativeZ
TH2 * m_MM_SmallWedge2_TransverseView_positiveZ
virtual StatusCode processEvent() override final
TH2 * m_MMTransverseEta2LargeWedge2
TH2 * m_MM_SmallWedge2_rZview_negativeZ
TH2 * m_MM_SmallWedge1_TransverseView_positiveZ
TH2 * m_MM_LargeWedge2_rZview_positiveZ
TH2 * m_MMTransverseEta2SmallWedge1
TH2 * m_MM_LargeWedge2_TransverseView_negativeZ
TH2 * m_MMTransverseEta1SmallWedge2
virtual StatusCode initialize() override final
TH2 * m_MM_LargeWedge1_rZview_negativeZ
static const MicromegasHitIdHelper * GetHelper()
std::string GetStationName(const int &hid) const
StatusCode executeFillHistos(const Amg::Vector3D &)
virtual StatusCode initialize() override
Eigen::Matrix< double, 3, 1 > Vector3D
Ensure that the Athena extensions are properly loaded.
Definition GeoMuonHits.h:27
STL namespace.