ATLAS Offline Software
Loading...
Searching...
No Matches
SCT_CalibBsErrorTool.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
12#include "SCT_CalibUtilities.h"
13#include "SCT_CalibNumbers.h"
14
15#include "Identifier/Identifier.h"
17
18// RAW data access
21
22#include "GaudiKernel/ITHistSvc.h"
23
24#include "TFile.h"
25#include "TH1I.h"
26
27#include <set>
28
29using namespace SCT_CalibAlgs;
30
31static const std::string pathRoot{"/BSErrors/"};
32static const std::string detectorPaths[] {"SCTEC/", "SCTB/","SCTEA/"};
33static const std::string detectorNames[] {"negativeEndcap", "barrel", "positiveEndcap"};
34static const int n_BSErrorType{15};
35static const int firstBSErrorType{0};
36static const int lastBSErrorType{14};
37
38SCT_CalibBsErrorTool::SCT_CalibBsErrorTool(const std::string& type, const std::string& name, const IInterface* parent) :
39 base_class(type, name, parent)
40{
41}
42
43StatusCode
45 ATH_CHECK( (m_thistSvc = service("THistSvc")).isValid() );
46 ATH_CHECK(detStore()->retrieve(m_pSCTHelper, "SCT_ID"));
48
49 m_maxHash = m_pSCTHelper->wafer_hash_max();
50 m_waferItrBegin = m_pSCTHelper->wafer_begin();
51 m_waferItrEnd = m_pSCTHelper->wafer_end();
52
53 return StatusCode::SUCCESS;
54}
55
56StatusCode
58 ATH_MSG_INFO("Finalize of SCT_CalibBsErrorTool");
59
60 return StatusCode::SUCCESS;
61}
62
63bool
65 bool result{true};
66 m_phistoVector.clear();
67 std::string histoName{pathRoot+"GENERAL/"};
68 //histogram for numbers of events
69 m_numberOfEventsHisto = new TH1I{"events", "Events", 1, 0.5, 1.5};
70 if (m_thistSvc->regHist(histoName.c_str(), m_numberOfEventsHisto ).isFailure()) {
71 ATH_MSG_ERROR("Error in booking BSErrors histogram");
72 }
73 //--- BSErrors for each wafer
76 for (; waferItr not_eq waferItrE; ++waferItr) {
77 Identifier waferId{*waferItr};
78 const int bec{m_pSCTHelper->barrel_ec(waferId)};
79 const std::string formattedPosition{formatPosition(waferId, m_pSCTHelper)};
80 std::string histotitle{std::string{"SCT "} + detectorNames[bec2Index(bec)] + std::string{" BSErrors : plane "} + formattedPosition};
81 const std::string name{pathRoot+detectorPaths[bec2Index(m_pSCTHelper->barrel_ec(waferId))] + formattedPosition};
82 TH1F* hitmapHisto_tmp{new TH1F{TString{formattedPosition}, TString{histotitle}, n_BSErrorType, firstBSErrorType-0.5, lastBSErrorType+0.5}};
83 if (m_thistSvc->regHist(name.c_str(), hitmapHisto_tmp).isFailure()) {
84 ATH_MSG_ERROR("Error in booking BSErrors histogram");
85 }
86 m_phistoVector.push_back(hitmapHisto_tmp);
87 }
88 return result;
89}
90
91bool
92SCT_CalibBsErrorTool::read(const std::string& fileName) {
93 bool result{true};
94 //pointers to the histos are deleted by m_thistSvc methods
95 m_phistoVector.clear();
96 TFile* fileHitmap{TFile::Open(fileName.c_str())};
97 ATH_MSG_INFO("opening Hitmap file : " << fileName.c_str());
98
99 if (fileHitmap==nullptr) {
100 ATH_MSG_ERROR("can not open Hitmap file : " << fileName.c_str());
101 return result;
102 }
103 //histogram for numbers of events
104 m_numberOfEventsHisto = static_cast<TH1I*>(fileHitmap->Get("GENERAL/events"));
105 if (m_numberOfEventsHisto==nullptr) {
106 ATH_MSG_ERROR("Error in reading EventNumber histogram");
107 }
108 //histograms for each wafer
110 for (; waferItr not_eq m_waferItrEnd; ++waferItr) {
111 const Identifier& waferId{*waferItr};
112 const std::string formattedPosition{formatPosition(waferId, m_pSCTHelper)};
113 std::string name{detectorPaths[bec2Index(m_pSCTHelper->barrel_ec(waferId))] + formattedPosition};
114 TH1F* hitmapHisto_tmp{static_cast<TH1F*>(fileHitmap->Get(name.c_str()))};
115 if (hitmapHisto_tmp==nullptr) {
116 ATH_MSG_ERROR("Error in reading BSErrors histogram");
117 } else {
118 m_phistoVector.push_back(hitmapHisto_tmp);
119 }
120 }
121 return result;
122}
123
124/*
125bool
126SCT_CalibBsErrorTool::read(const std::string& fileName) {
127 ATH_MSG_ERROR("Reading BsError histograms from " << fileName.c_str() << " is not supported!");
128 return false;
129}
130*/
131
132bool
133SCT_CalibBsErrorTool::fill(const bool fromData) {
134 if (fromData) {
135 return fillFromData();
136 }
137 bool result{true};
138 //--- Number of event
139 m_numberOfEventsHisto->Fill( 1 );
140 //--- Fill BSErrors
141 const int maxEntry{static_cast<int>(m_scterr_type->size())};
142 for (int i{0}; i != maxEntry; ++i ) {
143 int bec{(*m_scterr_bec)[i]};
144 int layer{(*m_scterr_layer)[i]};
145 int phi{(*m_scterr_phi)[i]};
146 int eta{(*m_scterr_eta)[i]};
147 int side{(*m_scterr_side)[i]};
148 int type{(*m_scterr_type)[i]};
149 Identifier waferId{m_pSCTHelper->wafer_id( bec, layer, phi, eta, side )};
150 fillBsErrorsForWafer(waferId, type);
151 }
152 return result;
153}
154
155bool
157 bool result{true};
158 //--- Number of event
159 m_numberOfEventsHisto->Fill( 1 );
160 //--- Loop over BSErrors
161 const EventContext& ctx{Gaudi::Hive::currentContext()};
162 for (int type = 0; type < SCT_ByteStreamErrors::NUM_ERROR_TYPES; ++type) {
163 const std::set<IdentifierHash> errorSet{m_bytestreamErrorsTool->getErrorSet(type, ctx)};
164 for(const auto& idHash : errorSet) {
165 Identifier waferId{m_pSCTHelper->wafer_id(idHash)};
166 fillBsErrorsForWafer(waferId, type);
167 }
168 }
169 return result;
170}
171
172void
174 int iWaferHash{static_cast<int>(m_pSCTHelper->wafer_hash(waferId))};
175 const std::string osWafer{formatPosition(waferId, m_pSCTHelper,".")};
176 //--- Protection for wrong waferID
177 if ( iWaferHash < 0 || iWaferHash >= m_maxHash ) {
178 ATH_MSG_WARNING("WaferHash " << iWaferHash << " is out of range : [ bec.layer.eta.phi.side, BSErrorType ] = [ " << osWafer << ", " << type << " ]");
179 } else {
180 ATH_MSG_DEBUG("BSError : [ bec.layer.eta.phi.side, Type ] = [ " << osWafer<< ", " << type << " ]");
181 m_phistoVector[ iWaferHash ]->Fill( type );
182 }
183}
Scalar eta() const
pseudorapidity method
Scalar phi() const
phi method
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
bool isValid(const T &p)
Av: we implement here an ATLAS-sepcific convention: all particles which are 99xxxxx are fine.
Definition AtlasPID.h:878
static const std::string detectorPaths[]
static const int n_BSErrorType
static const std::string pathRoot
static const int lastBSErrorType
static const int firstBSErrorType
static const std::string detectorNames[]
Header file for the SCT_CalibBsErrorTool class.
header file for the SCTCalibUtilities
virtual bool read(const std::string &fileName)
SCT_ID::const_id_iterator m_waferItrBegin
SCT_CalibBsErrorTool(const std::string &, const std::string &, const IInterface *)
void fillBsErrorsForWafer(const Identifier &waferId, const int type)
virtual StatusCode finalize()
SCT_ID::const_id_iterator m_waferItrEnd
ToolHandle< ISCT_ByteStreamErrorsTool > m_bytestreamErrorsTool
virtual StatusCode initialize()
virtual bool fill(const bool fromData=false)
std::vector< Identifier >::const_iterator const_id_iterator
Definition SCT_ID.h:73
std::string formatPosition(const Identifier &waferId, const SCT_ID *helper, const std::string &delimiter, const bool includeSide)
unsigned int bec2Index(const int bec)