ATLAS Offline Software
Loading...
Searching...
No Matches
TileHitToNtuple.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//*****************************************************************************
6// Filename : TileHitToNtuple.cxx
7// Author : Gia
8// Created : March, 2003
9//
10// DESCRIPTION:
11// Implement the algorithm
12//
13// HISTORY:
14//
15// BUGS:
16//
17//*****************************************************************************
18
19//Gaudi includes
20#include "GaudiKernel/INTupleSvc.h"
21#include "GaudiKernel/IDataProviderSvc.h"
22#include "GaudiKernel/SmartDataPtr.h"
23
24//Atlas include
27
28// Calo include
30
31//TileCalo include
34
35TileHitToNtuple::TileHitToNtuple(const std::string& name, ISvcLocator* pSvcLocator)
36 : AthAlgorithm(name, pSvcLocator)
37 , m_ntuplePtr(0)
38 , m_ntupleID("h30")
39 , m_ntupleLoc("/FILE1/TileRec")
40 , m_hitContainer("TileHitCnt")
41 , m_tileID(0)
42{
43 declareProperty("TileHitContainer", m_hitContainer);
44 declareProperty("NTupleLoc", m_ntupleLoc);
45 declareProperty("NTupleID", m_ntupleID);
46 declareProperty("CommitNtuple", m_commitNtuple = true);
47 declareProperty("MaxLength", m_maxLength = 4999);
48
49}
50
54
55//****************************************************************************
56//* Initialization
57//****************************************************************************
58
60{
61
62 ATH_MSG_INFO( "Initialization started" );
63
64 // retrieve TileID helper from det store
65 CHECK( detStore()->retrieve(m_tileID) );
66
67 m_ntupleLoc="/NTUPLES" + m_ntupleLoc;
68
69 SmartDataPtr<NTuple::Directory> DirPtr(ntupleSvc(), m_ntupleLoc);
70 if(!DirPtr) DirPtr=ntupleSvc()->createDirectory(m_ntupleLoc);
71 if(!DirPtr) {
72 ATH_MSG_ERROR( "Invalid Ntuple Directory: " );
73 return StatusCode::FAILURE;
74 }
75 m_ntuplePtr=ntupleSvc()->book(DirPtr.ptr(), m_ntupleID,
76 CLID_ColumnWiseTuple, "TileHit-Ntuple");
77 if(!m_ntuplePtr) {
78
79 std::string ntupleCompleteID=m_ntupleLoc+"/"+m_ntupleID;
80
81 NTuplePtr nt(ntupleSvc(),ntupleCompleteID);
82 if (!nt) {
83 ATH_MSG_ERROR( "Failed to book or to retrieve ntuple "
84 << ntupleCompleteID );
85 return StatusCode::FAILURE;
86 } else {
87 ATH_MSG_INFO( "Reaccessing ntuple " << ntupleCompleteID );
88 m_ntuplePtr = nt;
89 }
90 }
91
92 CHECK( m_ntuplePtr->addItem("TileHit/nhit",m_nchan,0,m_maxLength) );
93 CHECK( m_ntuplePtr->addItem("TileHit/totalE",m_tolE) );
94
95 CHECK( m_ntuplePtr->addItem("TileHit/energy",m_nchan,m_energy) );
96 CHECK( m_ntuplePtr->addItem("TileHit/time",m_nchan,m_time) );
97
98 CHECK( m_ntuplePtr->addItem("TileHit/detector",m_nchan,m_detector,0,3) );
99 CHECK( m_ntuplePtr->addItem("TileHit/side",m_nchan,m_side,-1,1) );
100 CHECK( m_ntuplePtr->addItem("TileHit/sample",m_nchan,m_sample,0,3) );
101 CHECK( m_ntuplePtr->addItem("TileHit/pmt", m_nchan, m_pmt,0,1) );
102 CHECK( m_ntuplePtr->addItem("TileHit/eta",m_nchan,m_eta,0,15) );
103 CHECK( m_ntuplePtr->addItem("TileHit/phi",m_nchan,m_phi,0,63) );
104
105 ATH_MSG_INFO( "Initialization completed" );
106 return StatusCode::SUCCESS;
107}
108
109//****************************************************************************
110//* Execution
111//****************************************************************************
112
114{
115
116 // step1: read TileHits from TDS
117 const TileHitContainer* HitCnt = nullptr;
118 CHECK( evtStore()->retrieve(HitCnt, m_hitContainer) );
119
120 // step2: put items in ntuple
121 SelectAllObject<TileHitContainer> selHits(HitCnt);
122 SelectAllObject<TileHitContainer>::const_iterator it=selHits.begin();
123 SelectAllObject<TileHitContainer>::const_iterator end=selHits.end();
124
125 m_nchan=0;
126 m_tolE=0.0;
127 int n_hit = 0;
128 for(; it != end; ++it) {
129
130 const TileHit * cinp = (*it);
131 Identifier id=cinp->identify();
132
133 int size = cinp->size();
134 for(int i=0;i<size;++i)
135 {
136 m_tolE+=cinp->energy(i);
137 m_energy[m_nchan]=cinp->energy(i);
138 m_time[m_nchan]=cinp->time(i);
139
140 m_detector[m_nchan]=m_tileID->section(id);
141 m_side[m_nchan]=m_tileID->side(id);
142 m_sample[m_nchan]=m_tileID->sample(id);
143 m_pmt[m_nchan]=m_tileID->pmt(id);
144 m_eta[m_nchan]=m_tileID->tower(id);
145 m_phi[m_nchan]=m_tileID->module(id);
146
147 m_nchan++;
148 if (m_nchan >= m_maxLength) break;
149 }
150
151 if (msgLvl(MSG::VERBOSE)) {
152 msg(MSG::VERBOSE) << " iHit=" << n_hit << " id="
153 << m_tileID->to_string(id, -1) << " ener=";
154
155 for (int i = 0; i < size; ++i)
156 msg(MSG::VERBOSE) << cinp->energy(i) << " ";
157
158 msg(MSG::VERBOSE) << "time=";
159 for (int i = 0; i < size; ++i)
160 msg(MSG::VERBOSE) << cinp->time() << " ";
161
162 msg(MSG::VERBOSE) << endmsg;
163 }
164
165 ++n_hit;
166 if (m_nchan >= m_maxLength) {
167 ATH_MSG_DEBUG( "Number of hits exceeds maximum (" << m_maxLength
168 << "), ignore all the rest" );
169 break;
170 }
171 }
172
173 // step3: commit ntuple
174 if ( m_commitNtuple ) {
175 ATH_MSG_DEBUG( "Committing Ntuple" );
176 CHECK( ntupleSvc()->writeRecord(m_ntuplePtr) );
177 }
178
179 // Execution completed.
180 ATH_MSG_DEBUG( "execute() completed successfully" );
181 return StatusCode::SUCCESS;
182}
183
184//****************************************************************************
185//* Finalize
186//****************************************************************************
187
189{
190 ATH_MSG_INFO( "finalize() completed successfully" );
191 return StatusCode::SUCCESS;
192}
193
#define endmsg
#define ATH_MSG_ERROR(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_DEBUG(x)
Helpers for checking error return status codes and reporting errors.
#define CHECK(...)
Evaluate an expression and check for errors.
SelectAllObjectMT< DCC, OBJECT > SelectAllObject
INTupleSvc * ntupleSvc()
AthAlgorithm(const std::string &name, ISvcLocator *pSvcLocator)
Constructor with parameters:
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)
const ServiceHandle< StoreGateSvc > & detStore() const
bool msgLvl(const MSG::Level lvl) const
MsgStream & msg() const
const_iterator end()
const_iterator begin()
NTuple::Array< int > m_detector
NTuple::Array< float > m_energy
NTuple::Array< int > m_sample
NTuple::Item< double > m_tolE
NTuple::Array< int > m_side
NTuple::Array< int > m_eta
virtual ~TileHitToNtuple()
StatusCode finalize()
std::string m_hitContainer
TileHitToNtuple(const std::string &name, ISvcLocator *pSvcLocator)
const TileID * m_tileID
NTuple::Array< int > m_phi
std::string m_ntupleID
StatusCode initialize()
NTuple::Tuple * m_ntuplePtr
NTuple::Array< float > m_time
NTuple::Array< int > m_pmt
NTuple::Item< int > m_nchan
std::string m_ntupleLoc
float time(int ind=0) const
Return time of ind-th sub-hit.
float energy(int ind=0) const
Return energy of ind-th sub-hit.
Identifier identify(void) const
Return logical ID of the pmt.
int size(void) const
Return length of energy/time vectors.