ATLAS Offline Software
Loading...
Searching...
No Matches
SCT_CalibHitmapTool.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3*/
4
10
11#include "SCT_CalibHitmapTool.h"
12#include "SCT_CalibUtilities.h"
13#include "SCT_CalibNumbers.h"
14
15// RAW data access
19
20#include "Identifier/Identifier.h"
22
23#include "GaudiKernel/ITHistSvc.h"
24
25#include "TH1I.h"
26#include "TH1F.h"
27#include "TFile.h"
28#include "TFileCollection.h"
29#include "THashList.h"
30
31using namespace SCT_CalibAlgs;
32static const std::string pathRoot{"/HitMaps/"};
33static const std::string detectorNames[] {"negativeEndcap", "barrel", "positiveEndcap"};
34static const std::string detectorPaths[] {"SCTEC/", "SCTB/","SCTEA/"};
35
36SCT_CalibHitmapTool::SCT_CalibHitmapTool(const std::string& type, const std::string& name, const IInterface* parent):
37 base_class(type, name, parent)
38{
39}
40
41StatusCode
43 m_thistSvc = service("THistSvc");
44 ATH_CHECK( m_thistSvc.isValid() );
45
46 ATH_CHECK(detStore()->retrieve(m_pSCTHelper, "SCT_ID"));
47 //
48 m_waferItrBegin = m_pSCTHelper->wafer_begin();
49 m_waferItrEnd = m_pSCTHelper->wafer_end();
50 //
51
52 // Read Handle Key
53 ATH_CHECK(m_rdoContainerKey.initialize());
54
55 return StatusCode::SUCCESS;
56}
57
58StatusCode
60 ATH_MSG_VERBOSE("SCT_CalibHitmapTool::finalize()");
61
62 return StatusCode::SUCCESS;
63}
64
65bool
67 bool result{true};
68 //pointers to the histos are deleted by m_thistSvc methods
69 m_phistoVector.clear();
70 std::string histoName{pathRoot + "GENERAL/"};
71 //histogram for numbers of events
72 m_numberOfEventsHisto = new TH1I{"events", "Events", 1, 0.5, 1.5};
73 if (m_thistSvc->regHist(histoName.c_str(), m_numberOfEventsHisto).isFailure()) {
74 ATH_MSG_ERROR("Error in booking EventNumber histogram");
75 }
76 //histograms for each wafer
78 std::string hitmapPaths[3];
79 for (int i{0}; i<3; ++i) {
80 hitmapPaths[i]=pathRoot+detectorPaths[i];
81 }
82 for (; waferItr not_eq m_waferItrEnd; ++waferItr) {
83 const Identifier& waferId{*waferItr};
84 const int bec{m_pSCTHelper->barrel_ec(waferId)};
85 const std::string formattedPosition{formatPosition(waferId, m_pSCTHelper)};
86 std::string histotitle{std::string{"SCT "} + detectorNames[bec2Index(bec)] + std::string{" Hitmap: plane "} + formattedPosition};
87 std::string name{hitmapPaths[bec2Index(m_pSCTHelper->barrel_ec(waferId))] + formattedPosition};
88 TH1F* hitmapHisto_tmp{new TH1F{TString{formattedPosition}, TString{histotitle}, nbins, firstStrip-0.5, lastStrip+0.5}};
89
90 if (m_thistSvc->regHist( name.c_str(), hitmapHisto_tmp ).isFailure()) {
91 ATH_MSG_ERROR("Error in booking Hitmap histogram");
92 } else {
93 m_phistoVector.push_back(hitmapHisto_tmp);
94 }
95
96 }
97 return result;
98}
99
100bool
101SCT_CalibHitmapTool::read(const std::string& fileName) {
102 bool result{true};
103 //pointers to the histos are deleted by m_thistSvc methods
104 m_phistoVector.clear();
105 TFile* fileHitmap{TFile::Open(fileName.c_str())};
106 ATH_MSG_INFO("opening Hitmap file : " << fileName.c_str());
107
108 if (fileHitmap==nullptr) {
109 ATH_MSG_ERROR("can not open Hitmap file : " << fileName.c_str());
110 return result;
111 }
112 //histogram for numbers of events
113 m_numberOfEventsHisto = static_cast<TH1I*>(fileHitmap->Get("GENERAL/events"));
114 if (m_numberOfEventsHisto==nullptr) {
115 ATH_MSG_ERROR("Error in reading EventNumber histogram");
116 }
117 //histograms for each wafer
119 for (; waferItr not_eq m_waferItrEnd; ++waferItr) {
120 const Identifier& waferId{*waferItr};
121 const std::string formattedPosition{formatPosition(waferId, m_pSCTHelper)};
122 std::string name{detectorPaths[bec2Index(m_pSCTHelper->barrel_ec(waferId))] + formattedPosition};
123 TH1F* hitmapHisto_tmp{static_cast<TH1F*>(fileHitmap->Get(name.c_str()))};
124 if (hitmapHisto_tmp==nullptr) {
125 ATH_MSG_ERROR("Error in reading Hitmap histogram");
126 } else {
127 m_phistoVector.push_back(hitmapHisto_tmp);
128 }
129 }
130 return result;
131}
132
133bool
134SCT_CalibHitmapTool::fill(const bool fromData) {
135 //cout<<"fromData "<<fromData<<endl;
136 if (fromData) {
137 return fillFromData();
138 }
139 bool result{true};
140 //--- Number of events
141 m_numberOfEventsHisto->Fill(1);
142
143 //--- Fill hitmap
144 const int MaxEntry{static_cast<int>(m_sct_waferHash.size())};
145 for (int i{0}; i != MaxEntry; ++i) {
146 const int theFirstStrip{m_sct_firstStrip[i]};
147 const int endStrip{m_sct_rdoGroupSize[i] + theFirstStrip};
148 const int index{m_sct_waferHash[i]};
149 TH1F* pThisHisto{m_phistoVector[index]};
150 for (int strip{theFirstStrip}; strip!=endStrip; ++strip) {
151 pThisHisto->Fill(strip);
152 }
153 }
154 return result;
155}
156
157bool
159 bool result{true};
160 m_numberOfEventsHisto->Fill(1);
162 if (not prdoContainer.isValid()) ATH_MSG_ERROR("Failed to retrieve the SCT RDO container");
163 SCT_RDO_Container::const_iterator itr{prdoContainer->begin()};
164 const SCT_RDO_Container::const_iterator end{prdoContainer->end()};
165 for (; itr !=end; ++itr) {
166 const InDetRawDataCollection<SCT_RDORawData>* SCT_Collection{*itr};
167 if (not SCT_Collection) continue;
168 const Identifier waferId{SCT_Collection->identify()};
169 const IdentifierHash waferHash{m_pSCTHelper->wafer_hash(waferId)};
170 TH1F* pThisHisto{m_phistoVector[static_cast<int>(waferHash)]};
171 DataVector<SCT_RDORawData>::const_iterator rdoItr{SCT_Collection->begin()};
172 const DataVector<SCT_RDORawData>::const_iterator rdoEnd{SCT_Collection->end()};
173 for (; rdoItr != rdoEnd; ++rdoItr) {
174 int strip{m_pSCTHelper->strip((*rdoItr)->identify())};
175 const int endStrip{(*rdoItr)->getGroupSize() + strip};
176 for (; strip != endStrip; ++strip) {
177 pThisHisto->Fill(strip);
178 }
179 }
180 }
181 return result;
182}
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_VERBOSE(x)
static const std::string detectorPaths[]
static const std::string pathRoot
static const std::string detectorNames[]
Header file for the SCT_CalibHitmapTool class.
header file for the SCTCalibUtilities
Handle class for reading from StoreGate.
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.
This is a "hash" representation of an Identifier.
virtual Identifier identify() const override final
SCT_ID::const_id_iterator m_waferItrEnd
SG::ReadHandleKey< SCT_RDO_Container > m_rdoContainerKey
virtual bool read(const std::string &fileName)
virtual bool fill(const bool fromData=false)
virtual StatusCode initialize()
const SCT_ID * m_pSCTHelper
SCT_CalibHitmapTool(const std::string &, const std::string &, const IInterface *)
virtual StatusCode finalize()
SCT_ID::const_id_iterator m_waferItrBegin
std::vector< Identifier >::const_iterator const_id_iterator
Definition SCT_ID.h:73
virtual bool isValid() override final
Can the handle be successfully dereferenced?
std::string formatPosition(const Identifier &waferId, const SCT_ID *helper, const std::string &delimiter, const bool includeSide)
unsigned int bec2Index(const int bec)
Definition index.py:1