ATLAS Offline Software
Loading...
Searching...
No Matches
Muon::SegmentPlots Struct Reference

segment plots More...

#include <MuonInsideOutAnalysisPlots.h>

Collaboration diagram for Muon::SegmentPlots:

Public Member Functions

void book (TDirectory *dir, const TString &prefix)
void fill (int chIndex_, int sector_, int quality_)
void fill (const MuonValidationSegmentBlock &segments, int index, float p_, float betaTrack)

Public Attributes

TH1 * t0
TH1 * t0Trig
TH1 * beta
TH2 * beta2
TH1 * betaTrig
TH1 * t0Res
TH1 * t0ResTrig
TH1 * betaRes
TH1 * betaResTrig
TH1 * quality
TH2 * quality_chIndex
TH2 * quality_sector
ResPlots allx
ChamberResPlots chamberx
ResPlots ally
ChamberResPlots chambery
ResPlots allxz
ChamberResPlots chamberxz
ResPlots allyz
ChamberResPlots chamberyz
ResPlots allcy
ChamberResPlots chambercy

Detailed Description

segment plots

Definition at line 114 of file MuonInsideOutAnalysisPlots.h.

Member Function Documentation

◆ book()

void Muon::SegmentPlots::book ( TDirectory * dir,
const TString & prefix )

Definition at line 210 of file MuonInsideOutAnalysisPlots.cxx.

210 {
211 t0 = new TH1F(prefix+"t0","t0",100,-20.,50.);
212 t0Trig = new TH1F(prefix+"t0Trig","t0Trig",100,-20.,50.);
213 beta = new TH1F(prefix+"beta","beta",100,.0,2.);
214 beta2 = new TH2F(prefix+"betaTrack","betaTrack",100,.0,2.,100,.0,2.);
215 betaTrig = new TH1F(prefix+"betaTrig","betaTrig",100,.0,2.);
216 t0Res = new TH1F(prefix+"t0Res","t0Res",200,-15,15);
217 t0ResTrig = new TH1F(prefix+"t0TrigRes","t0TrigRes",200,-15,15);
218 betaRes = new TH1F(prefix+"betaRes","betaRes",200,-1.2,1.2);
219 betaResTrig = new TH1F(prefix+"betaTrigRes","betaTrigRes",200,-1.2,1.2);
220 quality = new TH1F(prefix+"quality","quality",6,-0.5,5.5);
221 quality_chIndex = new TH2F(prefix+"quality_chIndex","quality_chIndex",6,-0.5,5.5,toInt(ChIndex::ChIndexMax),-0.5,-0.5+toInt(ChIndex::ChIndexMax));
222 quality_sector = new TH2F(prefix+"quality_sector","quality_sector",6,-0.5,5.5,16,0.5,16.5);
223
224 allx.book(prefix+"x_");
225 chamberx.book(dir,prefix+"x_");
226
227 ally.book(prefix+"y_");
228 chambery.book(dir,prefix+"y_");
229
230 allxz.book(prefix+"xz_");
231 chamberxz.book(dir,prefix+"xz_");
232
233 allyz.book(prefix+"yz_");
234 chamberyz.book(dir,prefix+"yz_");
235
236 allcy.book(prefix+"cy_");
237 chambercy.book(dir,prefix+"cy_");
238 }
constexpr int toInt(const EnumType enumVal)
TH2F(name, title, nxbins, bins_par2, bins_par3, bins_par4, bins_par5=None, bins_par6=None, path='', **kwargs)
TH1F(name, title, nxbins, bins_par2, bins_par3=None, path='', **kwargs)

◆ fill() [1/2]

void Muon::SegmentPlots::fill ( const MuonValidationSegmentBlock & segments,
int index,
float p_,
float betaTrack )

Definition at line 241 of file MuonInsideOutAnalysisPlots.cxx.

241 {
242
243 allx.fill( (*segments.xresiduals.residual)[index],
244 (*segments.xresiduals.pull)[index], (*segments.xresiduals.expos_err)[index]);
245 chamberx.fill( (*segments.id.chIndex)[index], (*segments.xresiduals.residual)[index],
246 (*segments.xresiduals.pull)[index], (*segments.xresiduals.expos_err)[index], p_);
247
248 ally.fill( (*segments.yresiduals.residual)[index],
249 (*segments.yresiduals.pull)[index], (*segments.yresiduals.expos_err)[index]);
250 chambery.fill( (*segments.id.chIndex)[index], (*segments.yresiduals.residual)[index],
251 (*segments.yresiduals.pull)[index], (*segments.yresiduals.expos_err)[index], p_);
252
253 allxz.fill( (*segments.angleXZ.residual)[index],
254 (*segments.angleXZ.pull)[index], (*segments.angleXZ.expos_err)[index]);
255 chamberxz.fill( (*segments.id.chIndex)[index], (*segments.angleXZ.residual)[index],
256 (*segments.angleXZ.pull)[index], (*segments.angleXZ.expos_err)[index], p_);
257
258 allyz.fill( (*segments.angleYZ.residual)[index],
259 (*segments.angleYZ.pull)[index], (*segments.angleYZ.expos_err)[index]);
260 chamberyz.fill( (*segments.id.chIndex)[index], (*segments.angleYZ.residual)[index],
261 (*segments.angleYZ.pull)[index], (*segments.angleYZ.expos_err)[index], p_);
262
263 allcy.fill( (*segments.combinedYZ.residual)[index],
264 (*segments.combinedYZ.pull)[index], (*segments.combinedYZ.expos_err)[index]);
265 chambercy.fill( (*segments.id.chIndex)[index], (*segments.combinedYZ.residual)[index],
266 (*segments.combinedYZ.pull)[index], (*segments.combinedYZ.expos_err)[index], p_);
267
268 fill( (*segments.id.chIndex)[index], (*segments.id.sector)[index], (*segments.quality)[index] );
269
270 float dist = sqrt( (*segments.r)[index]*(*segments.r)[index] + (*segments.z)[index]*(*segments.z)[index] );
271 MuonBetaCalculationUtils betaUtils(0.);
272 float tof = betaUtils.calculateTof(1.,dist);
273 float t0_ = (*segments.t0)[index]+tof;
274 float t0Trig_ = (*segments.t0Trig)[index]+tof;
275 float t0Track = betaUtils.calculateTof(betaTrack,dist);
276 float betaSeg = betaUtils.calculateBeta(t0_,dist);
277 float betaSegTrig = betaUtils.calculateBeta(t0Trig_,dist);
278 t0->Fill((*segments.t0)[index]);
279 t0Trig->Fill((*segments.t0Trig)[index]);
280 beta->Fill(betaSeg);
281 beta2->Fill(betaSeg,betaTrack);
282 betaTrig->Fill(betaSegTrig);
283 t0Res->Fill(t0_-t0Track);
284 t0ResTrig->Fill(t0Trig_-t0Track);
285 betaRes->Fill(betaSeg-betaTrack);
286 betaResTrig->Fill(betaSegTrig-betaTrack);
287 }
void fill(int chIndex_, int sector_, int quality_)

◆ fill() [2/2]

void Muon::SegmentPlots::fill ( int chIndex_,
int sector_,
int quality_ )

Definition at line 289 of file MuonInsideOutAnalysisPlots.cxx.

289 {
290 quality->Fill(quality_);
291 quality_chIndex->Fill(quality_,chIndex_);
292 quality_sector->Fill(quality_,sector_);
293 }

Member Data Documentation

◆ allcy

ResPlots Muon::SegmentPlots::allcy

Definition at line 140 of file MuonInsideOutAnalysisPlots.h.

◆ allx

ResPlots Muon::SegmentPlots::allx

Definition at line 128 of file MuonInsideOutAnalysisPlots.h.

◆ allxz

ResPlots Muon::SegmentPlots::allxz

Definition at line 134 of file MuonInsideOutAnalysisPlots.h.

◆ ally

ResPlots Muon::SegmentPlots::ally

Definition at line 131 of file MuonInsideOutAnalysisPlots.h.

◆ allyz

ResPlots Muon::SegmentPlots::allyz

Definition at line 137 of file MuonInsideOutAnalysisPlots.h.

◆ beta

TH1* Muon::SegmentPlots::beta

Definition at line 117 of file MuonInsideOutAnalysisPlots.h.

◆ beta2

TH2* Muon::SegmentPlots::beta2

Definition at line 118 of file MuonInsideOutAnalysisPlots.h.

◆ betaRes

TH1* Muon::SegmentPlots::betaRes

Definition at line 122 of file MuonInsideOutAnalysisPlots.h.

◆ betaResTrig

TH1* Muon::SegmentPlots::betaResTrig

Definition at line 123 of file MuonInsideOutAnalysisPlots.h.

◆ betaTrig

TH1* Muon::SegmentPlots::betaTrig

Definition at line 119 of file MuonInsideOutAnalysisPlots.h.

◆ chambercy

ChamberResPlots Muon::SegmentPlots::chambercy

Definition at line 141 of file MuonInsideOutAnalysisPlots.h.

◆ chamberx

ChamberResPlots Muon::SegmentPlots::chamberx

Definition at line 129 of file MuonInsideOutAnalysisPlots.h.

◆ chamberxz

ChamberResPlots Muon::SegmentPlots::chamberxz

Definition at line 135 of file MuonInsideOutAnalysisPlots.h.

◆ chambery

ChamberResPlots Muon::SegmentPlots::chambery

Definition at line 132 of file MuonInsideOutAnalysisPlots.h.

◆ chamberyz

ChamberResPlots Muon::SegmentPlots::chamberyz

Definition at line 138 of file MuonInsideOutAnalysisPlots.h.

◆ quality

TH1* Muon::SegmentPlots::quality

Definition at line 124 of file MuonInsideOutAnalysisPlots.h.

◆ quality_chIndex

TH2* Muon::SegmentPlots::quality_chIndex

Definition at line 125 of file MuonInsideOutAnalysisPlots.h.

◆ quality_sector

TH2* Muon::SegmentPlots::quality_sector

Definition at line 126 of file MuonInsideOutAnalysisPlots.h.

◆ t0

TH1* Muon::SegmentPlots::t0

Definition at line 115 of file MuonInsideOutAnalysisPlots.h.

◆ t0Res

TH1* Muon::SegmentPlots::t0Res

Definition at line 120 of file MuonInsideOutAnalysisPlots.h.

◆ t0ResTrig

TH1* Muon::SegmentPlots::t0ResTrig

Definition at line 121 of file MuonInsideOutAnalysisPlots.h.

◆ t0Trig

TH1* Muon::SegmentPlots::t0Trig

Definition at line 116 of file MuonInsideOutAnalysisPlots.h.


The documentation for this struct was generated from the following files: