ATLAS Offline Software
Loading...
Searching...
No Matches
ALFAHitAnalysis.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
5#include "ALFAHitAnalysis.h"
6
7// Section of includes for Pixel and SCT tests
10
11
12
14 ATH_MSG_DEBUG( "Initializing ALFAHitAnalysis" );
15
16 ATH_CHECK(m_readKey.initialize());
18 std::stringstream s;
19 for (unsigned int j=0; j<8; j++) {
20 s.str("");
21 s << "edep_in_det_no." << j+1;
22 float Emax = 5;
23 if (j==3) Emax = 150;
24 m_h_E_full_sum_h[j] = new TH1D(s.str().c_str(), s.str().c_str(), 100, 0, Emax);
25 m_h_E_full_sum_h[j]->StatOverflows();
26 ATH_CHECK( histSvc()->regHist( m_path + m_h_E_full_sum_h[j]->GetName(), m_h_E_full_sum_h[j] ) );
27
28 s.str("");
29 s << "edep_per_layer_det_no." << j+1;
30 m_h_E_layer_sum_h[j] = new TH1D(s.str().c_str(), s.str().c_str(), 100, 0, 15);
31 m_h_E_layer_sum_h[j]->StatOverflows();
32 ATH_CHECK( histSvc()->regHist( m_path + m_h_E_layer_sum_h[j]->GetName(), m_h_E_layer_sum_h[j] ) );
33
34 s.str("");
35 s << "hit_layer_det_no." << j+1;
36 m_h_hit_layer[j] = new TH1D(s.str().c_str(), s.str().c_str(), 50, 0, 50);
37 m_h_hit_layer[j]->StatOverflows();
38 ATH_CHECK( histSvc()->regHist( m_path + m_h_hit_layer[j]->GetName(), m_h_hit_layer[j] ) );
39
40 s.str("");
41 s << "hit_fiber_det_no." << j+1;
42 m_h_hit_fiber[j] = new TH1D(s.str().c_str(), s.str().c_str(), 100, 0, 60);
43 m_h_hit_fiber[j]->StatOverflows();
44 ATH_CHECK( histSvc()->regHist( m_path + m_h_hit_fiber[j]->GetName(), m_h_hit_fiber[j] ) );
45 }
46
47 m_tree = new TTree("ALFA", "ALFA");
48 std::string fullNtupleName = "/" + m_ntupleFileName + "/" ;
49 ATH_CHECK(histSvc()->regTree(fullNtupleName,m_tree));
50
51 m_tree->Branch("plate", &m_plate);
52 m_tree->Branch("station", &m_station);
53 m_tree->Branch("energy", &m_energy);
54 m_tree->Branch("fiber", &m_fiber);
55
56 return StatusCode::SUCCESS;
57}
58
59
61 ATH_MSG_DEBUG( "In ALFAHitAnalysis::execute()" );
62
63 m_station->clear();
64 m_plate->clear();
65 m_fiber->clear();
66 m_sign->clear();
67 m_energy->clear();
68
69 //cleaning
70 int fiber, plate, sign, station;
71 double E_fiber_sum[8][10][64][2], E_full_sum[8], E_layer_sum[8][20];
72 for (int l= 0; l<8; l++) {
73 E_full_sum[l] = 0.;
74 for (int i = 0; i < 10; i++) {
75 E_layer_sum[l][i] = 0.;
76 E_layer_sum[l][i+10] = 0.;
77 for (int j = 0; j < 64; j++) {
78 for (int k = 0; k < 2; k++) {
79 E_fiber_sum[l][i][j][k] = 0.;
80 }
81 }
82 }
83 }
84
86 const EventContext& ctx{Gaudi::Hive::currentContext()};
87 const ALFA_HitCollection* col_alfa{nullptr};
88 ATH_CHECK(SG::get(col_alfa, m_readKey, ctx));
89 for (iter = (*col_alfa).begin(); iter != (*col_alfa).end(); ++iter) {
90 station = (*iter).GetStationNumber();
91 plate = (*iter).GetPlateNumber();
92 fiber = (*iter).GetFiberNumber();
93 sign = (*iter).GetSignFiber();
94 E_fiber_sum[station-1][plate-1][fiber-1][(1-sign)/2] += ((*iter).GetEnergyDeposit());
95
96 m_station->push_back(station);
97 m_plate->push_back(plate);
98 m_fiber->push_back(fiber);
99 m_sign->push_back(sign);
100 m_energy->push_back((*iter).GetEnergyDeposit());
101 } // End iteration
102
103 for (int l=0; l<8; l++){
104 for (int i=0; i<10; i++){
105 for (int j=0; j<64; j++){
106 for (int k=0; k<2; k++){
107 E_full_sum[l] += E_fiber_sum[l][i][j][k];
108 E_layer_sum[l][2*i+k] += E_fiber_sum[l][i][j][k];
109 if (E_fiber_sum[l][i][j][k] > 0.) {
110 m_h_hit_layer[l]->Fill(2*i+k+1);
111 m_h_hit_fiber[l]->Fill(j+1);
112 }
113 }
114 }
115 }
116 }
117 for (int l=0; l<8; l++) {
118 m_h_E_full_sum_h[l]->Fill(E_full_sum[l]);
119 for (int i = 0; i< 20; i++) {
120 m_h_E_layer_sum_h[l]->Fill(E_layer_sum[l][i]);
121 }
122 }
123
124 if (m_tree) m_tree->Fill();
125
126 return StatusCode::SUCCESS;
127}
AtlasHitsVector< ALFA_Hit >::const_iterator ALFA_HitConstIter
AtlasHitsVector< ALFA_Hit > ALFA_HitCollection
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_DEBUG(x)
Handle class for reading from StoreGate.
int sign(int a)
std::array< TH1 *, 8 > m_h_hit_layer
std::array< TH1 *, 8 > m_h_E_full_sum_h
Some variables.
virtual StatusCode execute()
std::array< TH1 *, 8 > m_h_hit_fiber
std::vector< int > * m_fiber
std::vector< int > * m_station
Gaudi::Property< std::string > m_ntupleFileName
std::vector< double > * m_energy
std::vector< int > * m_sign
std::vector< int > * m_plate
SG::ReadHandleKey< ALFA_HitCollection > m_readKey
Gaudi::Property< std::string > m_path
std::array< TH1 *, 8 > m_h_E_layer_sum_h
virtual StatusCode initialize()
const ServiceHandle< ITHistSvc > & histSvc() const
The standard THistSvc (for writing histograms and TTrees and more to a root file) Returns (kind of) a...
const T * get(const ReadCondHandleKey< T > &key, const EventContext &ctx)
Convenience function to retrieve an object given a ReadCondHandleKey.