ATLAS Offline Software
Loading...
Searching...
No Matches
TileHitsTestTool Class Reference

#include <TileHitsTestTool.h>

Inheritance diagram for TileHitsTestTool:
Collaboration diagram for TileHitsTestTool:

Public Member Functions

 TileHitsTestTool (const std::string &type, const std::string &name, const IInterface *parent)
StatusCode initialize ()
StatusCode processEvent ()

Protected Attributes

std::string m_path {"/truth/"}
ServiceHandle< ITHistSvc > m_histSvc {"THistSvc", "SimTestHisto"}

Private Attributes

const TileIDm_tileID
const TileTBIDm_tileTBID
const TileDetDescrManagerm_tileMgr
TH1 * m_tile_eta
TH1 * m_tile_phi
TH1 * m_tile_energy
TH1 * m_tile_log_energy
TH1 * m_tile_time
TH2 * m_tile_rz
TH2 * m_tile_etaphi
TH2 * m_tile_energyphi
TH2 * m_tile_energyeta
TH1 * m_mbts_side
TH1 * m_mbts_eta
TH1 * m_mbts_phi
TH2 * m_mbts_sidetaphi
TH1 * m_etot
bool m_testMBTS

structors and AlgTool implementation

std::string m_key
 The MC truth key.
HepMC::ConstGenParticlePtr getPrimary ()

Detailed Description

Definition at line 15 of file TileHitsTestTool.h.

Constructor & Destructor Documentation

◆ TileHitsTestTool()

TileHitsTestTool::TileHitsTestTool ( const std::string & type,
const std::string & name,
const IInterface * parent )

Definition at line 18 of file TileHitsTestTool.cxx.

19 : SimTestToolBase(type, name, parent),
20 m_tileID(0), m_tileTBID(0), m_tileMgr(0),
25 m_etot(0),
26 m_testMBTS(true)
27{
28 declareProperty("TestMBTS",m_testMBTS,"");
29 }
SimTestToolBase(const std::string &type, const std::string &name, const IInterface *parent)
const TileTBID * m_tileTBID
const TileID * m_tileID
const TileDetDescrManager * m_tileMgr

Member Function Documentation

◆ getPrimary()

HepMC::ConstGenParticlePtr SimTestToolBase::getPrimary ( )
protectedinherited

Definition at line 20 of file SimTestToolBase.cxx.

21{
22 const McEventCollection* mcCollection;
23 if (evtStore()->retrieve(mcCollection,m_key).isSuccess()) {
25 for (e=mcCollection->begin();e!=mcCollection->end(); ++e) {
26 for (auto p : (**e)) {
28 return p;
29 }
30 }
31 }
32 }
33 return 0;
34}
DataModel_detail::const_iterator< DataVector > const_iterator
Standard const_iterator.
Definition DataVector.h:838
const_iterator end() const noexcept
Return a const_iterator pointing past the end of the collection.
const_iterator begin() const noexcept
Return a const_iterator pointing at the beginning of the collection.
std::string m_key
The MC truth key.
bool is_simulation_particle(const T &p)
Method to establish if a particle (or barcode) was created during the simulation (TODO update to be s...
retrieve(aClass, aKey=None)
Definition PyKernel.py:110

◆ initialize()

StatusCode TileHitsTestTool::initialize ( )
virtual

Reimplemented from SimTestToolBase.

Definition at line 32 of file TileHitsTestTool.cxx.

32 {
36
37 std::string detName("Tile"); //could configure this in the future?
38 m_path+=detName+"/";
39
40 // TileCal histograms
41 _TH1D(m_tile_eta,"tile_eta",17,-1.65,1.65);
42 _TH1D(m_tile_phi,"tile_phi",16,-3.2,3.2);
43 _TH1D(m_tile_energy,"tile_energy",50,0,100);
44 _TH1D(m_tile_log_energy,"tile_log_energy",100,-14.,8.);
45 _SET_TITLE(m_tile_log_energy,"logarithm of energy deposits","log(energy [MeV])","dN/dlog(energy [MeV])");
46 _TH1D(m_tile_time,"tile_time",100,0,250);
47
48 _TH2D(m_tile_rz,"tile_rz",16,-6200.,6200.,25,1500.,4000.);
49 _TH2D(m_tile_etaphi,"tile_etaphi",65,-1.625,1.625,129,-3.175,3.175);
50 _TH2D(m_tile_energyphi,"tile_energyphi",100,0.,10.,129,-3.175,3.175);
51 _TH2D(m_tile_energyeta,"tile_energyeta",100,0.,10.,129,-3.175,3.175);
52
53 // MBTS histograms
54 _TH1D(m_mbts_side,"mbts_side",3,-1.5,1.5);
55 _TH1D(m_mbts_eta,"mbts_eta",2,-0.5,1.5);
56 _TH1D(m_mbts_phi,"mbts_phi",8,-0.5,7.5);
57 _TH2D(m_mbts_sidetaphi,"mbts_sidetaphi",5,-2.5,2.5,7,-0.5,7.5);
58
59 _TH1D(m_etot,(detName+"_etot").c_str(),100,0.,500.);
60
61 return StatusCode::SUCCESS;
62}
#define CHECK(...)
Evaluate an expression and check for errors.
#define _TH2D(var, name, nbinx, xmin, xmax, nbiny, ymin, ymax)
#define _TH1D(var, name, nbin, xmin, xmax)
#define _SET_TITLE(var, title, xaxis, yaxis)
std::string m_path

◆ processEvent()

StatusCode TileHitsTestTool::processEvent ( )

Definition at line 65 of file TileHitsTestTool.cxx.

65 {
66 // make sure we still have a TileMgr (SimFlags.ReleaseGeoModel=False)
68
69 const TileHitVector* hitVec = nullptr;
70
71 double etot=0.;
72 if (evtStore()->retrieve(hitVec,"TileHitVec") == StatusCode::SUCCESS) {
73 for (const TileHit& hit : *hitVec) {
74
75 Identifier pmt_id = hit.identify();
76 Identifier cell_id = m_tileID->cell_id(pmt_id);
77
78 const CaloDetDescrElement *ddElement = m_tileMgr->get_cell_element(cell_id);
79 if(ddElement) // skip E4' cells which are not yet described in TileDetDescrManager.
80 {
81 double eta = ddElement->eta();
82 double phi = ddElement->phi();
83 double radius = ddElement->r();
84 double z = ddElement->z();
85
86 int pmt = m_tileID->pmt(pmt_id);
87 if (pmt>0) phi += ddElement->dphi()/2.; // to see clearly two PMTs in a cell
88
89 // loop over subhits
90 double energy=0;
91 // ATH_MSG(INFO) << "TileHit size :" <<hit.size()<<endmsg;
92 for (int i=0; i<hit.size();++i) {
93 energy+=hit.energy(i);
94 m_tile_energy->Fill(hit.energy(i));
95 m_tile_log_energy->Fill( std::log(hit.energy(i)) );
96 m_tile_time->Fill(hit.time(i));
97 }
98 etot+=energy;
99
100 m_tile_eta->Fill(eta);
101 m_tile_phi->Fill(phi);
102 m_tile_rz->Fill(z,radius);
103 m_tile_etaphi->Fill(eta,phi);
104 m_tile_energyphi->Fill(energy,phi);
105 m_tile_energyeta->Fill(energy,eta);
106 }
107 }
108 }
109 //Check the Hit container with postfix _Fast which generated with FastCaloSim
110 const TileHitVector *hitVec_fast = nullptr;
111 if( evtStore()->contains<TileHitVector>("TileHitVec_Fast") &&
112 (evtStore()->retrieve(hitVec_fast,"TileHitVec_Fast")).isSuccess())
113 {
114 ATH_MSG_DEBUG ( "Retrieve FastCaloSim container TileHitVec_Fast." );
115 for (const TileHit& hit : *hitVec_fast)
116 {
117 Identifier pmt_id=hit.identify();
118 Identifier cell_id=m_tileID->cell_id(pmt_id);
119 const CaloDetDescrElement *ddElement=m_tileMgr->get_cell_element(cell_id);
120 if(ddElement) // skip E4' cells which are not yet described in TileDetDescrManager.
121 {
122 double eta=ddElement->eta();
123 double phi=ddElement->phi();
124 double radius=ddElement->r();
125 double z=ddElement->z();
126
127 int pmt=m_tileID->pmt(pmt_id);
128 if(pmt>0) phi+=ddElement->dphi()/2.;
129 double energy=0;
130 for(int i=0;i<hit.size();++i)
131 {
132 energy+=hit.energy(i);
133 m_tile_energy->Fill(hit.energy(i));
134 m_tile_log_energy->Fill( std::log(hit.energy(i))) ;
135 m_tile_time->Fill(hit.time(i));
136 }
137 etot+=energy;
138
139 m_tile_eta->Fill(eta);
140 m_tile_phi->Fill(phi);
141 m_tile_rz->Fill(z,radius);
142 m_tile_etaphi->Fill(eta,phi);
143 m_tile_energyphi->Fill(energy,phi);
144 m_tile_energyeta->Fill(energy,eta);
145 }
146 }
147 }
148
149 if(m_testMBTS) {
150 if (evtStore()->retrieve(hitVec,"MBTSHits") == StatusCode::SUCCESS) {
151 for (const TileHit& hit : *hitVec) {
152 Identifier mbts_id = hit.identify();
153 double side = m_tileTBID->side(mbts_id); // -1 or +1
154 double ieta = m_tileTBID->eta(mbts_id); // 0 for inner cell, 1 for outer cell
155 double iphi = m_tileTBID->phi(mbts_id); // 0-7, cell 0 at phi=0
156 m_mbts_side->Fill(side);
157 m_mbts_eta->Fill(ieta);
158 m_mbts_phi->Fill(iphi);
159 m_mbts_sidetaphi->Fill((ieta+1)*side,iphi);
160 }
161 }
162 }
163 m_etot->Fill(etot);
164
165
166
167 return StatusCode::SUCCESS;
168}
Scalar eta() const
pseudorapidity method
Scalar phi() const
phi method
#define ATH_MSG_DEBUG(x)
AtlasHitsVector< TileHit > TileHitVector
#define z
bool contains(const std::string &s, const std::string &regx)
does a string contain the substring
Definition hcg.cxx:114

Member Data Documentation

◆ m_etot

TH1* TileHitsTestTool::m_etot
private

Definition at line 38 of file TileHitsTestTool.h.

◆ m_histSvc

ServiceHandle<ITHistSvc> SimTestHisto::m_histSvc {"THistSvc", "SimTestHisto"}
protectedinherited

Definition at line 35 of file SimTestHisto.h.

35{"THistSvc", "SimTestHisto"};

◆ m_key

std::string SimTestToolBase::m_key
protectedinherited

The MC truth key.

Definition at line 34 of file SimTestToolBase.h.

◆ m_mbts_eta

TH1 * TileHitsTestTool::m_mbts_eta
private

Definition at line 35 of file TileHitsTestTool.h.

◆ m_mbts_phi

TH1 * TileHitsTestTool::m_mbts_phi
private

Definition at line 35 of file TileHitsTestTool.h.

◆ m_mbts_side

TH1* TileHitsTestTool::m_mbts_side
private

Definition at line 35 of file TileHitsTestTool.h.

◆ m_mbts_sidetaphi

TH2* TileHitsTestTool::m_mbts_sidetaphi
private

Definition at line 36 of file TileHitsTestTool.h.

◆ m_path

std::string SimTestHisto::m_path {"/truth/"}
protectedinherited

Definition at line 34 of file SimTestHisto.h.

34{"/truth/"};

◆ m_testMBTS

bool TileHitsTestTool::m_testMBTS
private

Definition at line 40 of file TileHitsTestTool.h.

◆ m_tile_energy

TH1 * TileHitsTestTool::m_tile_energy
private

Definition at line 32 of file TileHitsTestTool.h.

◆ m_tile_energyeta

TH2 * TileHitsTestTool::m_tile_energyeta
private

Definition at line 33 of file TileHitsTestTool.h.

◆ m_tile_energyphi

TH2 * TileHitsTestTool::m_tile_energyphi
private

Definition at line 33 of file TileHitsTestTool.h.

◆ m_tile_eta

TH1* TileHitsTestTool::m_tile_eta
private

Definition at line 32 of file TileHitsTestTool.h.

◆ m_tile_etaphi

TH2 * TileHitsTestTool::m_tile_etaphi
private

Definition at line 33 of file TileHitsTestTool.h.

◆ m_tile_log_energy

TH1 * TileHitsTestTool::m_tile_log_energy
private

Definition at line 32 of file TileHitsTestTool.h.

◆ m_tile_phi

TH1 * TileHitsTestTool::m_tile_phi
private

Definition at line 32 of file TileHitsTestTool.h.

◆ m_tile_rz

TH2* TileHitsTestTool::m_tile_rz
private

Definition at line 33 of file TileHitsTestTool.h.

◆ m_tile_time

TH1 * TileHitsTestTool::m_tile_time
private

Definition at line 32 of file TileHitsTestTool.h.

◆ m_tileID

const TileID* TileHitsTestTool::m_tileID
private

Definition at line 28 of file TileHitsTestTool.h.

◆ m_tileMgr

const TileDetDescrManager* TileHitsTestTool::m_tileMgr
private

Definition at line 30 of file TileHitsTestTool.h.

◆ m_tileTBID

const TileTBID* TileHitsTestTool::m_tileTBID
private

Definition at line 29 of file TileHitsTestTool.h.


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