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

#include <TileCellSelector.h>

Inheritance diagram for TileCellSelector:

Public Member Functions

 TileCellSelector (const std::string &name, ISvcLocator *pSvcLocator)
virtual ~TileCellSelector ()
virtual StatusCode initialize () override
virtual StatusCode execute () override
virtual StatusCode finalize () override
virtual StatusCode sysInitialize () override
 Override sysInitialize.
virtual const DataObjIDColl & extraOutputDeps () const override
 Return the list of extra output dependencies.
ServiceHandle< StoreGateSvc > & evtStore ()
 The standard StoreGateSvc (event store) Returns (kind of) a pointer to the StoreGateSvc.
const ServiceHandle< StoreGateSvc > & detStore () const
 The standard StoreGateSvc/DetectorStore Returns (kind of) a pointer to the StoreGateSvc.
virtual StatusCode sysStart () override
 Handle START transition.
virtual std::vector< Gaudi::DataHandle * > inputHandles () const override
 Return this algorithm's input handles.
virtual std::vector< Gaudi::DataHandle * > outputHandles () const override
 Return this algorithm's output handles.
Gaudi::Details::PropertyBase & declareProperty (Gaudi::Property< T, V, H > &t)
void updateVHKA (Gaudi::Details::PropertyBase &)
MsgStream & msg () const
bool msgLvl (const MSG::Level lvl) const

Protected Member Functions

void renounceArray (SG::VarHandleKeyArray &handlesArray)
 remove all handles from I/O resolution
std::enable_if_t< std::is_void_v< std::result_of_t< decltype(&T::renounce)(T)> > &&!std::is_base_of_v< SG::VarHandleKeyArray, T > &&std::is_base_of_v< Gaudi::DataHandle, T >, void > renounce (T &h)
void extraDeps_update_handler (Gaudi::Details::PropertyBase &ExtraDeps)
 Add StoreName to extra input/output deps as needed.

Private Types

typedef ServiceHandle< StoreGateSvcStoreGateSvc_t

Private Member Functions

int Are3FF (std::vector< float > &OptFilterDigits, int OptFilterGain, int ch_type)
void printCell (const TileCell *cell)
Gaudi::Details::PropertyBase & declareGaudiProperty (Gaudi::Property< T, V, H > &hndl, const SG::VarHandleKeyType &)
 specialization for handling Gaudi::Property<SG::VarHandleKey>

Private Attributes

unsigned int m_counter
unsigned int m_accept
unsigned int m_minCell
unsigned int m_maxCell
unsigned int m_minChan
unsigned int m_maxChan
unsigned int m_jump
unsigned int m_const
unsigned int m_overLG
unsigned int m_overHG
unsigned int m_underLG
unsigned int m_underHG
unsigned int m_dqerr
unsigned int m_dmuerr
unsigned int m_warnerr
const TileIDm_tileID
const TileHWIDm_tileHWID
const TileCablingServicem_cabling
ToolHandle< ITileBadChanToolm_tileBadChanTool {this, "TileBadChanTool", "TileBadChanTool", "Tile bad channel tool"}
SG::ReadHandleKey< TileDQstatusm_dqStatusKey {this, "TileDQstatus", "TileDQstatus", "TileDQstatus key"}
ToolHandle< ITileDCSToolm_tileDCS {this, "TileDCSTool", "TileDCSTool", "Tile DCS tool"}
unsigned int m_runNum
unsigned int m_lumiBlock
unsigned int m_evtNum
unsigned int m_evtBCID
unsigned int m_tileFlag
unsigned int m_tileError
std::vector< bool > m_chanBad
std::vector< float > m_chanEne
std::vector< float > m_chanTime
std::vector< float > m_chanDsp
std::vector< float > m_chanTDsp
std::vector< float > m_chanQua
std::vector< bool > m_chanSel
std::vector< bool > m_chanToSkip
std::vector< bool > m_drawerToSkip
bool m_readCells
bool m_readRawChannels
bool m_readDigits
SG::ReadHandleKey< CaloCellContainerm_cellContainerKey
SG::ReadHandleKey< TileDigitsContainerm_digitsContainerKey
SG::ReadHandleKey< TileRawChannelContainerm_rawChannelContainerKey
SG::ReadHandleKey< xAOD::EventInfom_eventInfoKey
float m_minEneCell
float m_maxEneCell
float m_minEneChan [3] {}
float m_maxEneChan [3] {}
float m_minTimeCell
float m_maxTimeCell
float m_minTimeChan [3] {}
float m_maxTimeChan [3] {}
int m_ptnEneCell
int m_ptnEneChan [3] {}
int m_ptnTimeCell
int m_ptnTimeChan [3] {}
int m_selectGain
bool m_skipGain [2] {}
bool m_bitEneCell [ptnlength] {}
bool m_bitTimeCell [ptnlength] {}
bool m_bitEneChan [3][ptnlength] {}
bool m_bitTimeChan [3][ptnlength] {}
float m_secondMaxLevel
float m_jumpDeltaHG
float m_jumpDeltaLG
float m_pedDeltaHG
float m_pedDeltaLG
int m_constLength
int m_minBadDMU
int m_maxBadDMU
int m_minBadMB
bool m_skipEmpty
bool m_skipMasked
bool m_skipMBTS
bool m_checkDCS
bool m_checkJumps
bool m_checkDMUs
bool m_checkOverLG
bool m_checkOverHG
bool m_checkUnderLG
bool m_checkUnderHG
float m_overflowLG
float m_overflowHG
float m_underflowLG
float m_underflowHG
bool m_checkWarning
bool m_checkError
bool m_printOnly
std::vector< int > m_drawer
std::vector< int > m_drawerToCheck
std::vector< int > m_chanToCheck
int m_maxVerboseCnt
std::vector< int > m_nDrawerOff
std::string m_infoName
const TileInfom_tileInfo
float m_ADCmaxMinusEps = 0.0F
float m_ADCmaskValueMinusEps = 0.0F
DataObjIDColl m_extendedExtraObjects
StoreGateSvc_t m_evtStore
 Pointer to StoreGate (event store by default)
StoreGateSvc_t m_detStore
 Pointer to StoreGate (detector store by default)
std::vector< SG::VarHandleKeyArray * > m_vhka
bool m_varHandleArraysDeclared

Detailed Description

Definition at line 40 of file TileCellSelector.h.

Member Typedef Documentation

◆ StoreGateSvc_t

typedef ServiceHandle<StoreGateSvc> AthCommonDataStore< AthCommonMsg< Algorithm > >::StoreGateSvc_t
privateinherited

Definition at line 388 of file AthCommonDataStore.h.

Constructor & Destructor Documentation

◆ TileCellSelector()

TileCellSelector::TileCellSelector ( const std::string & name,
ISvcLocator * pSvcLocator )

Definition at line 36 of file TileCellSelector.cxx.

37 : AthAlgorithm(name, pSvcLocator)
38 , m_counter(0)
39 , m_accept(0)
40 , m_minCell(0)
41 , m_maxCell(0)
42 , m_minChan(0)
43 , m_maxChan(0)
44 , m_jump(0)
45 , m_const(0)
46 , m_overLG(0)
47 , m_overHG(0)
48 , m_underLG(0)
49 , m_underHG(0)
50 , m_dqerr(0)
51 , m_dmuerr(0)
52 , m_warnerr(0)
53 , m_tileID(0)
54 , m_tileHWID(0)
55 , m_cabling(0)
56 , m_runNum(0)
57 , m_lumiBlock(0)
58 , m_evtNum(0)
59 , m_evtBCID(0)
60 , m_tileFlag(0)
61 , m_tileError(0)
62 , m_readCells(true)
63 , m_readRawChannels(true)
64 , m_readDigits(true)
65 , m_tileInfo(0)
66{
67
68 declareProperty( "MinEnergyCell", m_minEneCell = -5000.); // cut on cell energy
69 declareProperty( "MaxEnergyCell", m_maxEneCell = 1000000.); // cut on cell energy
70 declareProperty( "PtnEnergyCell", m_ptnEneCell = 101); // cell energy pattern, accept events only below min (+1), between min-max (+10), above max (+100)
71 declareProperty( "MinEnergyChan", m_minEneChan[0] = -5000.); // cut on channel energy
72 declareProperty( "MaxEnergyChan", m_maxEneChan[0] = 500000.); // cut on channel energy
73 declareProperty( "PtnEnergyChan", m_ptnEneChan[0] = 101); // channel energy pattern
74 declareProperty( "MinEnergyGap", m_minEneChan[1] = -10000.); // cut on channel energy
75 declareProperty( "MaxEnergyGap", m_maxEneChan[1] = 500000.); // cut on channel energy
76 declareProperty( "PtnEnergyGap", m_ptnEneChan[1] = 101); // channel energy pattern
77 declareProperty( "MinEnergyMBTS", m_minEneChan[2] = -10000.); // cut on channel energy
78 declareProperty( "MaxEnergyMBTS", m_maxEneChan[2] = 500000.); // cut on channel energy
79 declareProperty( "PtnEnergyMBTS", m_ptnEneChan[2] = 101); // channel energy pattern
80
81 declareProperty( "MinTimeCell", m_minTimeCell = -100.); // cut on cell time
82 declareProperty( "MaxTimeCell", m_maxTimeCell = 100.); // cut on cell time
83 declareProperty( "PtnTimeCell", m_ptnTimeCell = 10); // cell time pattern, accept events only below min (+1), between min-max (+10), above max (+100)
84 declareProperty( "MinTimeChan", m_minTimeChan[0] = -100.); // cut on channel time
85 declareProperty( "MaxTimeChan", m_maxTimeChan[0] = 100.); // cut on channel time
86 declareProperty( "PtnTimeChan", m_ptnTimeChan[0] = 10); // channel time pattern
87 declareProperty( "MinTimeGap", m_minTimeChan[1] = -100.); // cut on channel time
88 declareProperty( "MaxTimeGap", m_maxTimeChan[1] = 100.); // cut on channel time
89 declareProperty( "PtnTimeGap", m_ptnTimeChan[1] = 10); // channel time pattern
90 declareProperty( "MinTimeMBTS", m_minTimeChan[2] = -100.); // cut on channel time
91 declareProperty( "MaxTimeMBTS", m_maxTimeChan[2] = 100.); // cut on channel time
92 declareProperty( "PtnTimeMBTS", m_ptnTimeChan[2] = 10); // channel time pattern
93
94 declareProperty( "SelectGain", m_selectGain = 2); // 0 - select LG only, 1 - HG only, 2 - both gains
97
98 // pattern - decimal number with up to 5 digits
99 // only values 1(=true) and 0(=false) for every digit are used
100 // digit 0 set to 1 - accept event if value < min
101 // digit 1 set to 1 - accept event if min < value < max
102 // digit 2 set to 1 - accept event if value > max
103 // digit 3 set to 1 - accept ene only if quality is good
104 // or accept time if time != 0
105 // digit 4 set to 1 - accept ene only if quality is bad
106 // or accept time if time == 0
107
108 declareProperty( "SecondMaxLevel",m_secondMaxLevel = 0.3); // sample below max should be above (max-min)*m_secondMax
109 declareProperty( "JumpDeltaHG", m_jumpDeltaHG = 50.0); // minimal jump in high gain
110 declareProperty( "JumpDeltaLG", m_jumpDeltaLG = 10.0); // minimal jump in low gain
111 declareProperty( "PedDetlaHG", m_pedDeltaHG = 4.1); // max variation of "const" value in high gain
112 declareProperty( "PedDetlaLG", m_pedDeltaLG = 4.1); // max variation of "const" value in low gain
113 declareProperty( "ConstLength", m_constLength = 6); // min number of consecutive samples of the same value
114 declareProperty( "MinBadDMU", m_minBadDMU = 4); // min number of bad DMUs to accept event
115 declareProperty( "MaxBadDMU", m_maxBadDMU = 15); // max number of bad DMUs to accept event
116 declareProperty( "MinBadMB", m_minBadMB = 4); // min number of bad motherboards in a drawer to accept event
117 declareProperty( "SkipEmpty", m_skipEmpty = true); // ignore empty channels in selection or not
118 declareProperty( "SkipMasked", m_skipMasked = true); // ignore masked channels in selection or not
119 declareProperty( "SkipMBTS", m_skipMBTS = true); // ignore MBTS channels in selection or not
120 declareProperty( "CheckDCS", m_checkDCS = true); // additional check for DCS status
121 declareProperty( "DrawerToDump", m_drawer); // for which drawer all channels should be printed
122 declareProperty( "DrawerToCheck",m_drawerToCheck); // for which drawer all checks should be performed
123 declareProperty( "ChannelToCheck",m_chanToCheck); // for which channels all checks should be performed
124
125 declareProperty( "CheckJumps", m_checkJumps = true); // global flag which allows to swithc on/off all checks in digits
126 declareProperty( "CheckDMUs", m_checkDMUs = true); // global flag which allows to swithc on/off DMU checks
127 declareProperty( "CheckOverLG" ,m_checkOverLG = true); // select events with overflow in low gain
128 declareProperty( "CheckOverHG", m_checkOverHG = false); // select events with overflow in high gain
129 declareProperty( "CheckUnderLG", m_checkUnderLG = false); // select events with underflow in low gain
130 declareProperty( "CheckUnderHG", m_checkUnderHG = false); // select events with underflow in high gain
131 declareProperty( "OverflowLG", m_overflowLG = -0.1); // threshold for overflow in low gain (smaller than ADCmax by this value)
132 declareProperty( "OverflowHG", m_overflowHG = -1.1); // threshold for overflow in high gain (smaller than ADCmax by this value)
133 declareProperty( "UnderflowLG", m_underflowLG = 0.1); // threshold for underflow in low gain
134 declareProperty( "UnderflowHG", m_underflowHG = 2.1); // threshold for underflow in high gain
135
136 declareProperty( "CheckWarning", m_checkWarning = false); // select events with warning status in TileCal status word
137 declareProperty( "CheckError", m_checkError = false); // select events with error status in TileCal status word
138 declareProperty( "PrintOnly", m_printOnly = false); // only print acccepted events, but do not accept anything
139
140 declareProperty( "MaxVerboseCnt",m_maxVerboseCnt=20); // max number of verbose output lines about drawer off
141
142 declareProperty("TileInfoName", m_infoName = "TileInfo");
143}
AthAlgorithm()
Default constructor:
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)
const TileCablingService * m_cabling
std::vector< int > m_chanToCheck
unsigned int m_minCell
unsigned int m_jump
unsigned int m_runNum
unsigned int m_lumiBlock
unsigned int m_evtBCID
unsigned int m_evtNum
unsigned int m_warnerr
unsigned int m_const
unsigned int m_underLG
unsigned int m_dmuerr
std::vector< int > m_drawerToCheck
unsigned int m_underHG
unsigned int m_counter
unsigned int m_maxCell
std::vector< int > m_drawer
unsigned int m_tileFlag
unsigned int m_accept
unsigned int m_overLG
const TileHWID * m_tileHWID
unsigned int m_maxChan
unsigned int m_minChan
const TileInfo * m_tileInfo
unsigned int m_overHG
const TileID * m_tileID
unsigned int m_tileError
unsigned int m_dqerr

◆ ~TileCellSelector()

TileCellSelector::~TileCellSelector ( )
virtual

Definition at line 146 of file TileCellSelector.cxx.

146 {
147}

Member Function Documentation

◆ Are3FF()

int TileCellSelector::Are3FF ( std::vector< float > & OptFilterDigits,
int OptFilterGain,
int ch_type )
private

Definition at line 2111 of file TileCellSelector.cxx.

2111 {
2112 bool allSaturated = true;
2113 int error = 0;
2114
2115 unsigned int nSamp = OptFilterDigits.size();
2116 if (nSamp) {
2117 float dmin = OptFilterDigits[0];
2118 float dmax = dmin;
2119
2120 for (unsigned int i = 1; i < nSamp; ++i) {
2121 float dig = OptFilterDigits[i];
2122 if (dig > dmax) dmax = dig;
2123 else if (dig < dmin) dmin = dig;
2124 }
2125 allSaturated = (dmin > m_ADCmaxMinusEps);
2126
2127 // FIXME:: set these parameters from JobOptions
2128 // FIXME:: move this method to base class
2129 const float epsilon = 4.1; // allow +/- 2 counts fluctuations around const value
2130 const float delta[4] = {29.9, 29.9, 49.9, 99.9}; // jump levels between constLG, constHG, non-constLG, non-constHG
2131 const float level0 = 29.9; // jump from this level to zero is bad
2132 const float level1 = 99.9; // jump from this level to m_tileInfo->ADCmax() is bad
2133 const float level2 = 199.9; // base line at this level is bad
2134 const float delt = std::min(std::min(std::min(delta[0], delta[1]), std::min(delta[2], delta[3])), level0);
2135
2136 if (!allSaturated && (dmax - dmin) > delt) {
2137 float abovemin = dmax;
2138 float belowmax = dmin;
2139 unsigned int nmin = 0;
2140 unsigned int nmax = 0;
2141 unsigned int pmin = nSamp;
2142 unsigned int pmax = nSamp;
2143 for (unsigned int i = 0; i < nSamp; ++i) {
2144 float smp = OptFilterDigits[i];
2145 if (smp - dmin < epsilon) {
2146 ++nmin;
2147 pmin = i;
2148 }
2149 if (dmax - smp < epsilon) {
2150 ++nmax;
2151 pmax = i;
2152 }
2153 if (smp < abovemin && smp > dmin) {
2154 abovemin = smp;
2155 }
2156 if (smp > belowmax && smp < dmax) {
2157 belowmax = smp;
2158 }
2159 }
2160
2161 if (abovemin != dmax || belowmax != dmin) { // more than two different values
2162 gain += 2; // shift index by 2, i.e. use thresholds for non-const levels
2163 }
2164
2165 if (dmin < 0.01 && dmax > m_ADCmaxMinusEps) { // jump from zero to saturation
2166 error = 1;
2167 } else if (dmin < 0.01 && abovemin > level0 && nmin > 1) { // at least two samples at zero, others - above pedestal
2168 error = 2;
2169 } else if (dmax > m_ADCmaxMinusEps && belowmax < level1 && nmax > 1) { // at least two saturated. others - close to pedestal
2170 error = 3;
2171 } else if (dmin>level2 && (gain==0 || ch_type<2) ) { // baseline above threshold is bad
2172 error = 9; // but should not apply that to MBTS
2173 } else if (nmax+nmin==nSamp && (dmax-dmin) > delta[gain]) {
2174 if (nmax>1 && nmin>1) { // at least 2 samples at two distinct levels
2175 error = 4;
2176 } else if (nmax==1) {
2177 if (pmax>0 && pmax<nSamp-1) { // jump up in one sample, but not at the edge
2178 error = 5;
2179 }
2180 } else if (nmin==1) { // jump down in one sample
2181 error = 6;
2182 }
2183 }
2184 if (!error && (dmax - dmin) > delta[gain]) {
2185 float secondMax = (dmax - dmin) * m_secondMaxLevel;
2186 if (pmax > 0 && pmax < nSamp - 1
2187 && std::max(OptFilterDigits[pmax - 1], OptFilterDigits[pmax + 1]) < dmin + secondMax) {
2188
2189 error = 7; // jump up in one sample in the middle. which is much higher than all others
2190 } else if (pmin > 0 && pmin < nSamp - 1
2191 && std::min(OptFilterDigits[pmin - 1], OptFilterDigits[pmin + 1]) > dmax - secondMax) {
2192
2193 error = 8; // jump down in one sample. which is much lower than all others
2194 }
2195 }
2196 }
2197 }
2198
2199 if (allSaturated)
2200 return 99;
2201 else
2202 return error;
2203}
const int nmax(200)

◆ declareGaudiProperty()

Gaudi::Details::PropertyBase & AthCommonDataStore< AthCommonMsg< Algorithm > >::declareGaudiProperty ( Gaudi::Property< T, V, H > & hndl,
const SG::VarHandleKeyType &  )
inlineprivateinherited

specialization for handling Gaudi::Property<SG::VarHandleKey>

Definition at line 156 of file AthCommonDataStore.h.

158 {
160 hndl.value(),
161 hndl.documentation());
162
163 }

◆ declareProperty()

Gaudi::Details::PropertyBase & AthCommonDataStore< AthCommonMsg< Algorithm > >::declareProperty ( Gaudi::Property< T, V, H > & t)
inlineinherited

Definition at line 145 of file AthCommonDataStore.h.

145 {
146 typedef typename SG::HandleClassifier<T>::type htype;
148 }
Gaudi::Details::PropertyBase & declareGaudiProperty(Gaudi::Property< T, V, H > &hndl, const SG::VarHandleKeyType &)
specialization for handling Gaudi::Property<SG::VarHandleKey>

◆ detStore()

const ServiceHandle< StoreGateSvc > & AthCommonDataStore< AthCommonMsg< Algorithm > >::detStore ( ) const
inlineinherited

The standard StoreGateSvc/DetectorStore Returns (kind of) a pointer to the StoreGateSvc.

Definition at line 95 of file AthCommonDataStore.h.

◆ evtStore()

ServiceHandle< StoreGateSvc > & AthCommonDataStore< AthCommonMsg< Algorithm > >::evtStore ( )
inlineinherited

The standard StoreGateSvc (event store) Returns (kind of) a pointer to the StoreGateSvc.

Definition at line 85 of file AthCommonDataStore.h.

◆ execute()

StatusCode TileCellSelector::execute ( )
overridevirtual

Definition at line 375 of file TileCellSelector.cxx.

375 {
376 //ATH_MSG_DEBUG ("execute()");
377
378 const EventContext& ctx = Gaudi::Hive::currentContext();
379
380 ++m_counter;
381
382 SG::ReadHandle<xAOD::EventInfo> eventInfo(m_eventInfoKey, ctx);
383 if ( eventInfo.isValid() ) {
384 m_runNum = eventInfo->runNumber();
385 m_lumiBlock = eventInfo->lumiBlock();
386 m_evtNum = eventInfo->eventNumber();
387 m_evtBCID = eventInfo->bcid();
388 m_tileFlag = eventInfo->eventFlags(EventInfo::Tile);
389 m_tileError = eventInfo->errorState(EventInfo::Tile);
390 } else {
391 m_runNum = 0;
392 m_evtNum = 0;
393 m_lumiBlock = 0;
394 m_evtBCID = 0;
395 m_tileFlag = 0;
396 m_tileError = 0;
397 }
398
399 std::ostringstream evState;
400 evState << "Run "<< std::setw(6) << m_runNum
401 <<" LB "<< std::setw(4) << m_lumiBlock
402 <<" Evt "<< std::setw(9) << m_evtNum
403 <<" ErrState " << m_tileError
404 <<" Flags 0x" << std::hex << m_tileFlag << std::dec;
405
406 std::ostringstream evtnum;
407 evtnum << "Run "<< std::setw(6) << m_runNum
408 <<" LB "<< std::setw(4) << m_lumiBlock
409 <<" Evt "<< std::setw(9) << m_evtNum
410 <<" BCID "<< std::setw(4) << m_evtBCID;
411
412 std::ostringstream nevtnum;
413 nevtnum << evtnum.str()
414 <<" nevt "<< std::setw(6) << m_counter;
415
416 bool emptyBad = true;
417 bool badFromCell = false;
418 m_chanBad.clear();
419 m_chanBad.resize(1+TileHWID::NOT_VALID_HASH,true);
420 m_chanEne.clear();
422 m_chanTime.clear();
424 m_chanQua.clear();
426 m_chanSel.clear();
427 m_chanSel.resize(1+TileHWID::NOT_VALID_HASH,false);
428
429 IdContext chan_context = m_tileHWID->channel_context();
430 IdentifierHash hash;
431 for (size_t i=0; i<m_drawer.size(); ++i) {
432 HWIdentifier drawer_id = m_tileHWID->drawer_id(m_drawer[i]);
433 HWIdentifier ch_id = m_tileHWID->channel_id(drawer_id,0);
434 m_tileHWID->get_hash(ch_id, hash, &chan_context);
435 auto itr = m_chanSel.begin() + hash;
436 std::fill(itr,itr+48,true);
437 }
438
439 bool statusOk = (m_checkWarning && m_tileError == EventInfo::Warning) ||
441 if (statusOk) {
442 ++m_warnerr;
443
444 using namespace boost::local_time;
445 using namespace boost::posix_time;
446 static const time_zone_ptr gva_tz(new posix_time_zone((std::string)"CET+01CEST01:00:00,M3.5.0/02:00:00,M10.5.0/03:00:00"));
447 local_date_time gva_time(from_time_t(eventInfo->timeStamp()),gva_tz);
448 evState << " " << gva_time << " ";
449
450
451 const char * part[5] = { "UNK", "LBA", "LBC", "EBA", "EBC" };
453 int dn = (m_tileFlag >> 16) & 0xF;
454 int n1 = (m_tileFlag >> 20) & 0x3F;
455 int n2 = (n1 + dn - 1) % 64;
456 int rr = ((m_tileFlag >> 26) & 0x3) + 1;
457 evState << " " << part[rr] <<std::setw(2)<<std::setfill('0')<<n1+1
458 << " - " << part[rr] <<std::setw(2)<<std::setfill('0')<<n2+1
459 << " " << dn << " consec bad ";
460 }
461 else if (m_tileError == EventInfo::Error) {
462 int dn = (m_tileFlag >> 16) & 0xF;
463 int p0 = m_tileFlag & 0xF; // sends any data above threshold
464 int p1 = (m_tileFlag >> 20) & 0xF; // 16 drawers masked
465 int p2 = (m_tileFlag >> 24) & 0xF; // 16 drawers off
466 int pp = 1;
467 for (int rr = 1; rr < 5; ++rr) {
468 if ((p2 & pp) || (p1 & pp)) {
469 evState << " " << part[rr];
470 if (p2 & pp) {
471 if (p0 & pp) evState << " off";
472 else evState << " OFF";
473 }
474 if (p1 & pp) {
475 if (p0 & pp) evState << " mask";
476 else evState << " MASK";
477 }
478 }
479 pp <<= 1;
480 }
481 evState << " " << dn << " consec bad ";
482 }
483
484 if (m_checkDCS) {
485 int n1 = -1;
486 int n2 = -1;
487 int dn = 0;
488 int rr = 0;
489 int m1 = -1;
490 int m2 = -1;
491 int dm = -1;
492 std::vector<int> allmod;
493 std::vector<int> consec;
494 for (int ros = 1; ros < 5; ++ros) {
495 int drmax = 65;
496 for (int dr = 0; dr < drmax; ++dr) {
497 int drawer = dr % 64;
498 if (m_tileDCS->getDCSStatus(ros, drawer) == TileDCSState::ALERT_DRAWER) {
499 if (m1 < 0) m1 = dr;
500 m2 = dr;
501 if (dr < 64) allmod.push_back((ros << 8) + dr);
502 } else if (m1 >= 0) {
503 dm = m2 - m1 + 1;
504 if (m1 == 0) drmax += dm;
505 if (dm > dn) {
506 n1 = m1;
507 n2 = m2;
508 dn = dm;
509 rr = ros;
510 consec.clear();
511 consec.push_back((ros << 8) + m1);
512 } else if (dm == dn) {
513 if (m1 < 64) consec.push_back((ros << 8) + m1);
514 }
515 m1 = m2 = -1;
516 }
517 }
518 }
519 evState << " DCS " << allmod.size() << " off ";
520 if (dn > 1) {
521 evState << dn;
522 if (consec.size() > 1) evState << "*" << (consec.size());
523 evState << " consec "
524 << part[rr] << std::setw(2) << std::setfill('0') << (n1 % 64) + 1 << " - "
525 << part[rr] << std::setw(2) << std::setfill('0') << (n2 % 64) + 1 << " ";
526 n1 += (rr << 8);
527 n2 += (rr << 8);
528 for (size_t n = 1; n < consec.size(); ++n) {
529 m1 = consec[n];
530 m2 = m1 + dn - 1;
531 evState << part[m1 >> 8] << std::setw(2) << std::setfill('0') << (m1 % 64) + 1 << " - "
532 << part[m2 >> 8] << std::setw(2) << std::setfill('0') << (m2 % 64) + 1 << " ";
533 for (size_t m = 0; m < allmod.size(); ++m) {
534 int mm = allmod[m];
535 if (mm >= m1 && mm <= m2) {
536 allmod[m] += n1 - m1;
537 }
538 }
539 }
540 } else {
541 n1 = n2 = dn = 0;
542 }
543 if (allmod.size() > (size_t) dn) {
544 for (size_t m = 0; m < allmod.size(); ++m) {
545 int mm = allmod[m];
546 if (!(mm >= n1 && mm <= n2)) {
547 evState << part[mm >> 8] << std::setw(2) << std::setfill('0') << (mm % 64) + 1 << " ";
548 }
549 }
550 }
551 }
552
553 ATH_MSG_DEBUG (evState.str() << " accepted");
554 }
555
556 int rawdata = -1;
557 const TileCell* cellminCh = 0;
558 const TileCell* cellmaxCh = 0;
559 const TileCell* tcellminCh = 0;
560 const TileCell* tcellmaxCh = 0;
561
562 if (m_readCells) {
563
564 // Get Calo cell container
565 SG::ReadHandle<CaloCellContainer> cellContainer(m_cellContainerKey);
566
567 if (!cellContainer.isValid()) {
568
569 ATH_MSG_WARNING("Unable to read CaloCellContainer from EventStore, disable reading of this container");
570 m_readCells = false;
571
572 } else {
573
574 float emin = 0.;
575 float emax = 0.;
576 float tmin = 0.;
577 float tmax = 0.;
578 float chmin = 0.;
579 float chmax = 0.;
580 float tcmin = 0.;
581 float tcmax = 0.;
582 const TileCell* cellmin = 0;
583 const TileCell* cellmax = 0;
584 const TileCell* tcellmin = 0;
585 const TileCell* tcellmax = 0;
586
587 // special case - check overflow here if digits container is not available
588 // should be careful here, because in TileCell overflow bit is set to 1
589 // both for overflow and underflow and underflow is HG are very often in gap cells
590 // also, overflow in HG might be masked if quality is too bad, so we'll not select all overflows...
591 // that's why only overflow in LG are checked
592 bool checkOver = (m_checkOverLG && m_digitsContainerKey.key().empty());
593
594 for (const CaloCell* cell : *cellContainer) {
595
596 Identifier id = cell->ID();
597 if ( m_tileID->is_tile(id) ) {
598 const TileCell* tile_cell = dynamic_cast<const TileCell*> (cell);
599 if (tile_cell==0) continue;
600 const CaloDetDescrElement * caloDDE = cell->caloDDE();
601 IdentifierHash hash1 = caloDDE->onl1();
602 IdentifierHash hash2 = caloDDE->onl2();
603 if ( m_chanToSkip[hash1] && m_chanToSkip[hash2] ) continue;
604 int ch_type = (hash2 == TileHWID::NOT_VALID_HASH) ? 1 : 0;
605 if (rawdata < 0) {
606 rawdata = (tile_cell->qbit1() & TileCell::MASK_CMPC) ? 0 : 1;
607 }
608
609 bool bad1 = tile_cell->badch1();
610 bool bad2 = tile_cell->badch2();
611 float ene1 = tile_cell->ene1();
612 float ene2 = tile_cell->ene2();
613 float time1 = tile_cell->time1();
614 float time2 = tile_cell->time2();
615 m_chanBad[hash1] = bad1;
616 m_chanBad[hash2] = bad2;
617 m_chanEne[hash1] = ene1;
618 m_chanEne[hash2] = ene2;
619 m_chanTime[hash1] = time1;
620 m_chanTime[hash2] = time2;
621 m_chanQua[hash1] = tile_cell->qual1();
622 m_chanQua[hash2] = tile_cell->qual2();
623
624 float ene = tile_cell->energy();
625 bool eneOk = false;
626 if (ene < m_minEneCell) {
627 eneOk = m_bitEneCell[0];
628 } else if (ene > m_maxEneCell) {
629 eneOk = m_bitEneCell[2];
630 } else {
631 eneOk = m_bitEneCell[1];
632 }
633
634 if (eneOk) {
635 if (bad1 && bad2) {
636 if (m_bitEneCell[3]) eneOk = false; // request good cells only, but cell is bad
637 } else {
638 if (m_bitEneCell[4]) eneOk = false; // request bad cells only, but cell is good
639 }
640 }
641
642 float time = tile_cell->time();
643 bool timeOk = false;
644 if (time < m_minTimeCell) {
645 timeOk = m_bitTimeCell[0];
646 } else if (time > m_maxTimeCell ) {
647 timeOk = m_bitTimeCell[2];
648 } else {
649 timeOk = m_bitTimeCell[1];
650 }
651
652 if (timeOk) {
653 if (time != 0.) {
654 if (m_bitTimeCell[4]) timeOk = false; // request time==0 only, but time!=0
655 } else {
656 if (m_bitTimeCell[3]) timeOk = false; // request time!=0 only, but time==0
657 }
658 }
659
660 if (timeOk && eneOk) {
661
662 ATH_MSG_VERBOSE( evtnum.str()
663 << " cell " << std::left << std::setw(14) << m_tileID->to_string(id,-2)
664 << " ene = " << ene << " time = " << time);
665
666 m_chanSel[hash1] = true;
667 m_chanSel[hash2] = true;
668
669 if (ene < emin) {
670 emin = ene;
671 cellmin = tile_cell;
672 } else if (ene > emax) {
673 emax = ene;
674 cellmax = tile_cell;
675 }
676
677 if (time<tmin) {
678 tmin = time;
679 tcellmin = tile_cell;
680 }
681 else if (time>tmax) {
682 tmax = time;
683 tcellmax = tile_cell;
684 }
685 }
686
687 if ( !(bad1 && bad2) ) {
688
689 bool ene1Ok = false;
690 bool time1Ok = false;
691
692 if ( !(bad1 || m_skipGain[tile_cell->gain1()]) ) {
693 if (time1 < m_minTimeChan[ch_type] ) {
694 time1Ok = m_bitTimeChan[ch_type][0];
695 } else if (time1 > m_maxTimeChan[ch_type] ) {
696 time1Ok = m_bitTimeChan[ch_type][2];
697 } else {
698 time1Ok = m_bitTimeChan[ch_type][1];
699 }
700
701 if (ene1 < m_minEneChan[ch_type] ) {
702 ene1Ok = m_bitEneChan[ch_type][0];
703 } else if (ene1 > m_maxEneChan[ch_type] ) {
704 ene1Ok = m_bitEneChan[ch_type][2];
705 } else {
706 ene1Ok = m_bitEneChan[ch_type][1];
707 }
708
709 if (ene1Ok) {
710 if (m_bitEneChan[ch_type][4]) ene1Ok = false; // request bad chan only, but chan is good
711 }
712
713 if (time1Ok) {
714 if (time1 != 0.) {
715 if (m_bitTimeChan[ch_type][4]) time1Ok = false; // request time==0 only, but time!=0
716 } else {
717 if (m_bitTimeChan[ch_type][3]) time1Ok = false; // request time!=0 only, but time==0
718 }
719 }
720 }
721
722 bool ene2Ok = false;
723 bool time2Ok = false;
724
725 if ( !(bad2 || m_skipGain[tile_cell->gain2()]) ) {
726 if (ene2 < m_minEneChan[ch_type] ) {
727 ene2Ok = m_bitEneChan[ch_type][0];
728 } else if (ene2 > m_maxEneChan[ch_type] ) {
729 ene2Ok = m_bitEneChan[ch_type][2];
730 } else {
731 ene2Ok = m_bitEneChan[ch_type][1];
732 }
733
734 if (time2 < m_minTimeChan[ch_type] ) {
735 time2Ok = m_bitTimeChan[ch_type][0];
736 } else if (time2 > m_maxTimeChan[ch_type] ) {
737 time2Ok = m_bitTimeChan[ch_type][2];
738 } else {
739 time2Ok = m_bitTimeChan[ch_type][1];
740 }
741
742 if (ene2Ok) {
743 if (m_bitEneChan[ch_type][4]) ene2Ok = false; // request bad chan only, but chan is good
744 }
745
746 if (time2Ok) {
747 if (time2 != 0.) {
748 if (m_bitTimeChan[ch_type][4]) time2Ok = false; // request time==0 only, but time!=0
749 } else {
750 if (m_bitTimeChan[ch_type][3]) time2Ok = false; // request time!=0 only, but time==0
751 }
752 }
753 }
754
755 bool over1=false;
756 bool over2=false;
757 if (checkOver) {
758 over1 = ( (!bad1) && (tile_cell->qbit1() & TileCell::MASK_OVER) && tile_cell->gain1()==TileID::LOWGAIN);
759 over2 = ( (!bad2) && (tile_cell->qbit2() & TileCell::MASK_OVER) && tile_cell->gain2()==TileID::LOWGAIN);
760 }
761
762 if ((ene1Ok && time1Ok) || over1) {
763
764 ATH_MSG_VERBOSE( evtnum.str()
765 << " cell " << std::left << std::setw(14) << m_tileID->to_string(id,-2)
766 << " ch_ene1 = " << ene1 << " ch_t1 = " << time1
767 << ((over1)?" overflow":""));
768
769 m_chanSel[hash1] = true;
770 m_chanSel[hash2] = true;
771
772 if (ene1 < chmin) {
773 chmin = ene1;
774 cellminCh = tile_cell;
775 } else if (ene1 > chmax) {
776 chmax = ene1;
777 cellmaxCh = tile_cell;
778 }
779
780 if (time1 < tcmin) {
781 tcmin = time1;
782 tcellminCh = tile_cell;
783 } else if (time1 > tcmax) {
784 tcmax = time1;
785 tcellmaxCh = tile_cell;
786 }
787 }
788
789 if ((ene2Ok && time2Ok) || over2) {
790
791 ATH_MSG_VERBOSE( evtnum.str()
792 << " cell " << std::left << std::setw(14) << m_tileID->to_string(id,-2)
793 << " ch_ene2 = " << ene2 << " ch_t2 = " << time2
794 << ((over2)?" overflow":""));
795
796 m_chanSel[hash1] = true;
797 m_chanSel[hash2] = true;
798
799 if (ene2 < chmin) {
800 chmin = ene2;
801 cellminCh = tile_cell;
802 } else if (ene2 > chmax) {
803 chmax = ene2;
804 cellmaxCh = tile_cell;
805 }
806
807 if (time2 < tcmin) {
808 tcmin = time2;
809 tcellminCh = tile_cell;
810 } else if (time2 > tcmax) {
811 tcmax = time2;
812 tcellmaxCh = tile_cell;
813 }
814 }
815
816 }
817 }
818 }
819
820 if (tcellmin && tcellmin != cellmin && tcellmin != cellmax) {
821 ATH_MSG_DEBUG( nevtnum.str()
822 << " cell " << std::left << std::setw(14) << m_tileID->to_string(tcellmin->ID(),-2)
823 << " ene = " << tcellmin->energy()
824 << " tmin = " << tcellmin->time());
825 }
826 if (tcellmax && tcellmax != cellmin && tcellmax != cellmax) {
827 ATH_MSG_DEBUG( nevtnum.str()
828 << " cell " << std::left << std::setw(14) << m_tileID->to_string(tcellmax->ID(),-2)
829 << " ene = " << tcellmax->energy()
830 << " tmax = " << tcellmax->energy());
831 }
832
833 if (tcellminCh && tcellminCh != cellminCh && tcellminCh != cellmaxCh) {
834 ATH_MSG_DEBUG( nevtnum.str()
835 << " cell " << std::left << std::setw(14) << m_tileID->to_string(tcellminCh->ID(),-2)
836 << " ch_ene = " << tcellminCh->ene1() << " " << tcellminCh->ene2()
837 << " ch_tmin = " << tcellminCh->time1() << " " << tcellminCh->time2());
838 }
839 if (tcellmaxCh && tcellmaxCh != cellminCh && tcellmaxCh != cellmaxCh) {
840 ATH_MSG_DEBUG( nevtnum.str()
841 << " cell " << std::left << std::setw(14) << m_tileID->to_string(tcellmaxCh->ID(),-2)
842 << " ch_ene = " << tcellmaxCh->ene1() << " " << tcellmaxCh->ene2()
843 << " ch_tmax = " << tcellmaxCh->time1() << " " << tcellmaxCh->time2());
844 }
845
846 if (cellmin) {
847 ++m_minCell;
848 statusOk = true;
849 const char * tit = (tcellmin == cellmin) ? " tmin = ": ((tcellmax == cellmin) ? " tmax = ": " t = ");
850 if (cellminCh!=cellmin) {
851 ATH_MSG_DEBUG( nevtnum.str()
852 << " cell " << std::left << std::setw(14) << m_tileID->to_string(cellmin->ID(),-2)
853 << " emin = " << emin
854 << tit << cellmin->time()
855 << " accepted");
856
857 } else {
858 ATH_MSG_DEBUG( nevtnum.str()
859 << " cell " << std::left << std::setw(14) << m_tileID->to_string(cellmin->ID(),-2)
860 << " emin = " << emin
861 << " ch_emin = " << chmin
862 << tit << cellmin->time()
863 << " accepted");
864 }
865 }
866
867 if (cellminCh) {
868 ++m_minChan;
869 statusOk = true;
870 const char * tit = (tcellminCh == cellminCh) ? " tmin = ": ((tcellmaxCh == cellminCh) ? " tmax = ": " t = ");
871 if (cellminCh!=cellmin) {
872 ATH_MSG_DEBUG( nevtnum.str()
873 << " cell " << std::left << std::setw(14) << m_tileID->to_string(cellminCh->ID(),-2)
874 << " ch_emin = " << chmin
875 << tit << cellminCh->time()
876 << " accepted");
877 }
878 }
879
880 if (cellmax) {
881 ++m_maxCell;
882 statusOk = true;
883 const char * tit = (tcellmin == cellmax) ? " tmin = ": ((tcellmax == cellmax) ? " tmax = ": " t = ");
884 if (cellmaxCh!=cellmax) {
885 ATH_MSG_DEBUG( nevtnum.str()
886 << " cell " << std::left << std::setw(14) << m_tileID->to_string(cellmax->ID(),-2)
887 << " emax = " << emax
888 << tit << cellmax->time()
889 << " accepted");
890
891 } else {
892 ATH_MSG_DEBUG( nevtnum.str()
893 << " cell " << std::left << std::setw(14) << m_tileID->to_string(cellmax->ID(),-2)
894 << " emax = " << emax
895 << " ch_emax = " << chmax
896 << tit << cellmax->time()
897 << " accepted");
898 }
899 }
900
901 if (cellmaxCh) {
902 ++m_maxChan;
903 statusOk = true;
904 const char * tit = (tcellminCh == cellmaxCh) ? " tmin = ": ((tcellmaxCh == cellmaxCh) ? " tmax = ": " t = ");
905 if (cellmaxCh!=cellmax) {
906 ATH_MSG_DEBUG( nevtnum.str()
907 << " cell " << std::left << std::setw(14) << m_tileID->to_string(cellmaxCh->ID(),-2)
908 << " ch_emax = " << chmax
909 << tit << cellmaxCh->time()
910 << " accepted");
911 }
912 }
913
914 emptyBad = false;
915 badFromCell = true;
916 }
917 }
918
919 const TileDQstatus* DQstatus(0);
920
921 if (m_readRawChannels) {
922
923 // Get Tile RawChannel container
924 SG::ReadHandle<TileRawChannelContainer> rawChannelContainer(m_rawChannelContainerKey);
925
926 if ( !rawChannelContainer.isValid() ) {
927 ATH_MSG_WARNING("Unable to read TileRawChannelContainer from EventStore, disable reading of this container");
928 m_readRawChannels = false;
929
930 } else {
931
932 float chmin = 0; // m_minEneChan[0];
933 float chmax = 0; // m_maxEneChan[0];
934 float tcmin = 0.;
935 float tcmax = 0.;
936 const TileRawChannel* minCh = 0;
937 const TileRawChannel* maxCh = 0;
938 const TileRawChannel* tminCh = 0;
939 const TileRawChannel* tmaxCh = 0;
940 TileRawChannelUnit::UNIT rChUnit = rawChannelContainer->get_unit();
941 bool allowAmpCheck = ( ( rChUnit == TileRawChannelUnit::MegaElectronVolts || // allow MeV only as units
943 bool fillChanEne = ( !m_readCells && allowAmpCheck ); // use amplitude from channel if cell container was not checked
944 if (!fillChanEne) {
945 m_chanDsp.clear();
947 m_chanTDsp.clear();
949 }
950
952 DQstatus = SG::makeHandle (m_dqStatusKey, ctx).get();
953 else
954 rawdata = 0;
955
956 IdContext chan_context = m_tileHWID->channel_context();
957 IdentifierHash hash;
958 int index=0, pmt;
959
960 int nbadMax = 0;
961 int nbadMBMax = 0;
962 const TileRawChannelCollection * collMax = 0;
963 const TileRawChannelCollection * collMBMax = 0;
964
965 bool someDQerrors = false;
966
967 for (const TileRawChannelCollection* rawChannelCollection : *rawChannelContainer) {
968
969 int frag = rawChannelCollection->identify();
970 bool eb = (frag > 0x2ff);
971 bool ebsp = (frag == 0x30e || frag == 0x411);
972
973 int ros = frag >> 8;
974 int drawer = frag & 0x3F;
975 unsigned int drawerIdx = TileCalibUtils::getDrawerIdx(ros,drawer);
976 if ( m_drawerToSkip[drawerIdx] ) continue;
977
978 int chMBTS = -1;
979 if (eb) {
980 for (int ch: {12,4,0}) {
981 m_cabling->h2s_cell_id_index(ros,drawer,ch,index,pmt);
982 if (index == -2) {
983 chMBTS = ch;
984 break;
985 }
986 }
987 }
988
989 HWIdentifier ch_id = m_tileHWID->channel_id(ros,drawer,0);
990 m_tileHWID->get_hash(ch_id, hash, &chan_context);
991 int hashNext = index = (int)hash + 48; // hash ID of the channel after current drawer
992
993 // all error words contains information in last 16 bits only
994 // but they are stored in collection as 32 bit numbers
995 uint32_t RODBCID = rawChannelCollection->getRODBCID();
996 uint32_t DSPBCID = rawChannelCollection->getFragDSPBCID();
997 uint32_t GlobalCRCErr = rawChannelCollection->getFragGlobalCRC() & 0x1;
998 uint32_t FE_DMUmask = rawChannelCollection->getFragFEChipMask();
999 uint32_t ROD_DMUmask = rawChannelCollection->getFragRODChipMask();
1000 uint32_t BCIDErr = rawChannelCollection->getFragBCID();
1001 uint32_t MemoryParityErr = rawChannelCollection->getFragMemoryPar();
1002 uint32_t HeaderFormatErr = rawChannelCollection->getFragHeaderBit();
1003 uint32_t HeaderParityErr = rawChannelCollection->getFragHeaderPar();
1004 uint32_t SampleFormatErr = rawChannelCollection->getFragSampleBit();
1005 uint32_t SampleParityErr = rawChannelCollection->getFragSamplePar();
1006 uint32_t SingleStrobeErr = rawChannelCollection->getFragSstrobe();
1007 uint32_t DoubleStrobeErr = rawChannelCollection->getFragDstrobe();
1008
1009 if (RODBCID!=0 && RODBCID != m_evtBCID ) {
1010 if (m_nDrawerOff[drawerIdx] < m_maxVerboseCnt) {
1011 ATH_MSG_VERBOSE( evtnum.str()
1012 << " drw " << drwname(rawChannelCollection->identify())
1013 << " ROD BCID " << RODBCID << " is wrong - skipping");
1014
1015 if (++m_nDrawerOff[drawerIdx] == m_maxVerboseCnt)
1016 ATH_MSG_VERBOSE( nevtnum.str()
1017 << " suppressing further messages about drawer 0x" << std::hex << rawChannelCollection->identify()
1018 << std::dec << " being bad");
1019 }
1020 someDQerrors = true;
1021 continue;
1022 }
1023
1024 if (DSPBCID >= 0x7FFF
1025 && GlobalCRCErr
1026 && FE_DMUmask == 0xFFFF
1027 && ROD_DMUmask == 0xFFFF
1028 && BCIDErr == 0xFFFF
1029 && MemoryParityErr == 0xFFFF
1030 && HeaderFormatErr == 0xFFFF
1031 && HeaderParityErr == 0xFFFF
1032 && SampleFormatErr == 0xFFFF
1033 && SampleParityErr == 0xFFFF
1034 && SingleStrobeErr == 0xFFFF
1035 && DoubleStrobeErr == 0xFFFF) {
1036
1037 if (m_nDrawerOff[drawerIdx] < m_maxVerboseCnt) {
1038 ATH_MSG_VERBOSE( evtnum.str()
1039 << " drw " << drwname(rawChannelCollection->identify())
1040 << " is OFF - skipping");
1041
1042 if (++m_nDrawerOff[drawerIdx] == m_maxVerboseCnt)
1043 ATH_MSG_VERBOSE( nevtnum.str()
1044 << " suppressing further messages about drawer 0x" << std::hex
1045 << rawChannelCollection->identify()
1046 << std::dec << " being bad");
1047 }
1048 continue;
1049 }
1050
1051 if (DSPBCID == 0
1052 && GlobalCRCErr == 0
1053 && FE_DMUmask == 0
1054 && ROD_DMUmask == 0
1055 && BCIDErr == 0
1056 && MemoryParityErr == 0
1057 && HeaderFormatErr == 0
1058 && HeaderParityErr == 0
1059 && SampleFormatErr == 0
1060 && SampleParityErr == 0
1061 && SingleStrobeErr == 0
1062 && DoubleStrobeErr == 0) {
1063
1064 if (m_nDrawerOff[drawerIdx] < m_maxVerboseCnt) {
1065 ATH_MSG_VERBOSE( evtnum.str()
1066 << " drw " << drwname(rawChannelCollection->identify())
1067 << " is MISSING - skipping");
1068
1069 if (++m_nDrawerOff[drawerIdx] == m_maxVerboseCnt)
1070 ATH_MSG_VERBOSE( nevtnum.str()
1071 << " suppressing further messages about drawer 0x" << std::hex
1072 << rawChannelCollection->identify() << std::dec << " being bad");
1073 }
1074 continue;
1075 }
1076
1077 if (GlobalCRCErr) {
1078 GlobalCRCErr = 0xFFFF; // global error - all wrong
1079 if (m_maxBadDMU<16) {
1080 if (m_nDrawerOff[drawerIdx] < m_maxVerboseCnt) {
1081 ATH_MSG_VERBOSE( evtnum.str()
1082 << " drw " << drwname(rawChannelCollection->identify())
1083 << " global CRC error - skipping");
1084
1085 if (++m_nDrawerOff[drawerIdx] == m_maxVerboseCnt)
1086 ATH_MSG_VERBOSE( nevtnum.str()
1087 << " suppressing further messages about drawer 0x" << std::hex
1088 << rawChannelCollection->identify() << std::dec << " being bad");
1089 }
1090 someDQerrors = true;
1091 continue;
1092 }
1093 }
1094
1095 if (HeaderFormatErr || HeaderParityErr || SampleFormatErr || SampleParityErr ) {
1096 FE_DMUmask = 0xFFFF; // can not trust FE mask, assume that all DMUs are good
1097 } else {
1098 if (eb) { // extended barrel
1099 if (ebsp) FE_DMUmask<<=1; // shift by one DMU in EBA15 EBC18
1100 FE_DMUmask = (FE_DMUmask & 0xFF) | ((FE_DMUmask & 0xF00)<<2); // shift upper half by two DMUs
1101 }
1102 }
1103
1104 FE_DMUmask = ~FE_DMUmask & 0xFFFF; // inversion for FE CRC
1105 ROD_DMUmask = ~ROD_DMUmask & 0xFFFF; // inversion for ROD CRC
1106
1107 if (BCIDErr & 0x2) { // DMU1 (second DMU) is bad - can not trust others
1108 BCIDErr = 0xFFFF; // assume that all DMUs are bad
1109 if (m_maxBadDMU < 16) {
1110 if (m_nDrawerOff[drawerIdx] < m_maxVerboseCnt) {
1111 ATH_MSG_VERBOSE( evtnum.str()
1112 << " drw " << drwname(rawChannelCollection->identify())
1113 << " BCID in DMU1 is bad - skipping");
1114
1115 if (++m_nDrawerOff[drawerIdx] == m_maxVerboseCnt)
1116 ATH_MSG_VERBOSE( nevtnum.str()
1117 << " suppressing further messages about drawer 0x"
1118 << std::hex << rawChannelCollection->identify() << std::dec << " being bad");
1119 }
1120 someDQerrors = true;
1121 continue;
1122 }
1123
1124 } else {
1125 // additional check if DQ frag BCID is the same as event BCID
1126 if ( DSPBCID!=0xDEAD && DSPBCID!=m_evtBCID ) { // DSP BCID doesn't match! all wrong
1127 BCIDErr = 0xFFFF;
1128 if (m_maxBadDMU < 16) {
1129 if (m_nDrawerOff[drawerIdx] < m_maxVerboseCnt) {
1130 ATH_MSG_VERBOSE( evtnum.str()
1131 << " drw " << drwname(rawChannelCollection->identify())
1132 << " DSP BCID is wrong - skipping");
1133
1134 if (++m_nDrawerOff[drawerIdx] == m_maxVerboseCnt)
1135 ATH_MSG_VERBOSE( nevtnum.str()
1136 << " suppressing further messages about drawer 0x"
1137 << std::hex << rawChannelCollection->identify() << std::dec << " being bad");
1138 }
1139 someDQerrors = true;
1140 continue;
1141 }
1142 }
1143 }
1144
1145 uint32_t error = GlobalCRCErr | FE_DMUmask | ROD_DMUmask | BCIDErr | MemoryParityErr |
1146 HeaderFormatErr | HeaderParityErr | SampleFormatErr | SampleParityErr;
1147
1148 if (error==0xFFFF && m_maxBadDMU<16) {
1149 if (m_nDrawerOff[drawerIdx] < m_maxVerboseCnt) {
1150 ATH_MSG_VERBOSE( evtnum.str()
1151 << " drw " << drwname(rawChannelCollection->identify())
1152 << " whole drawer is bad - skipping");
1153
1154 if (++m_nDrawerOff[drawerIdx] == m_maxVerboseCnt)
1155 ATH_MSG_VERBOSE( nevtnum.str()
1156 << " suppressing further messages about drawer 0x"
1157 << std::hex << rawChannelCollection->identify() << std::dec << " being bad");
1158 }
1159 someDQerrors = true;
1160 continue;
1161 }
1162
1163 // no global error detected - wait m_max_verbose_cnt good events and eventually enable error messages again
1164 if (m_nDrawerOff[drawerIdx]>=m_maxVerboseCnt) {
1165 if (++m_nDrawerOff[drawerIdx] == 2*m_maxVerboseCnt) {
1166 m_nDrawerOff[drawerIdx] = 0;
1167 ATH_MSG_VERBOSE( nevtnum.str()
1168 << " enabling messages about drawer 0x" << std::hex
1169 << rawChannelCollection->identify()
1170 << std::dec << " being bad after " << m_maxVerboseCnt << " good events");
1171 }
1172 }
1173
1174 uint32_t errMB = BCIDErr;
1175 if (eb) { // do not count non-existing DMUs in EB
1176 if (ebsp) {
1177 errMB &= 0x3cfe;
1178 } else {
1179 errMB &= 0x3cff;
1180 }
1181 }
1182 int nbadMB = 0;
1183 while (errMB) {
1184 if (errMB & 0xF) ++nbadMB;
1185 errMB >>= 4;
1186 }
1187 someDQerrors = (nbadMB >= m_minBadMB);
1188 if (nbadMB > nbadMBMax) {
1189 nbadMBMax = nbadMB;
1190 collMBMax = rawChannelCollection;
1191 }
1192
1193 int nerr = 0;
1194 for (uint32_t i = 0x8000; i != 0; i >>= 1) {
1195 if (error&i) {
1196 ++nerr;
1197 index-=3;
1198 } else {
1199 if (emptyBad && nbadMB < 4) {
1200 m_chanBad[--index] = false;
1201 m_chanBad[--index] = false;
1202 m_chanBad[--index] = false;
1203 } else {
1204 index -= 3;
1205 }
1206 }
1207 }
1208 //if (chMBTS>=0) m_chanBad[index+chMBTS] = true; // ignore completely MBTS channel
1209 int nbad = ((ebsp) ? nerr-5 : ((eb) ? nerr-4 : nerr));
1210
1211 if (nbad >= m_minBadDMU && nerr <= m_maxBadDMU) {
1212 someDQerrors = true;
1213 if (nbad > nbadMax) {
1214 nbadMax = nbad;
1215 collMax = rawChannelCollection;
1216 }
1217 }
1218 if (someDQerrors) { // will print later samples for all channels in a drawer
1219 for ( ; index<hashNext; ++index) { // but if drawer was completely bad, this loop will be skipped
1220 m_chanSel[index] = true;
1221 }
1222 }
1223
1224 if (allowAmpCheck || emptyBad) {
1225
1226 for (const TileRawChannel* rawChannel : *rawChannelCollection) {
1227
1228 HWIdentifier adcId = rawChannel->adc_HWID();
1229 HWIdentifier chId = m_tileHWID->channel_id(adcId);
1230 m_tileHWID->get_hash(chId, hash, &chan_context);
1231 if ( m_chanToSkip[hash] ) continue;
1232 int adc = m_tileHWID->adc(adcId);
1233 int channel = m_tileHWID->channel(adcId);
1234 int ch_type = 0;
1235 if (channel == chMBTS) {
1236 ch_type = 2;
1237 } else if ( (ebsp && (channel == 18 || channel == 19 || channel == 12 || channel == 13) )
1238 || (eb && (channel == 0 || channel == 1 || channel == 12 || channel == 13) ) ) {
1239 ch_type = 1;
1240 }
1241 if (emptyBad && !m_chanBad[hash] ) {
1242 m_chanBad[hash] = m_tileBadChanTool->getAdcStatus(drawerIdx,channel,adc).isBad() ||
1243 (DQstatus && !DQstatus->isAdcDQgood(ros,drawer,channel,adc)) ||
1244 (m_checkDCS && m_tileDCS->getDCSStatus(ros, drawer, channel) > TileDCSState::WARNING);
1245 }
1246
1247 if (allowAmpCheck) {
1248
1249 float amp = rawChannel->amplitude();
1250 float time = rawChannel->time();
1251 if (fillChanEne) {
1252 m_chanEne[hash] = amp;
1253 m_chanTime[hash] = time;
1254 m_chanQua[hash] = rawChannel->quality();
1255 } else {
1256 m_chanDsp[hash] = amp;
1257 m_chanTDsp[hash] = time;
1258 }
1259
1260 if ( (m_skipMasked && m_chanBad[hash]) ||
1261 (m_skipMBTS && channel == chMBTS) ||
1262 (m_skipEmpty && TileDQstatus::isChEmpty(ros, drawer, channel) > 0) ||
1263 m_skipGain[adc] )
1264 continue;
1265
1266 bool ampOk = false;
1267 if (amp < m_minEneChan[ch_type] ) {
1268 ampOk = m_bitEneChan[ch_type][0];
1269 } else if (amp > m_maxEneChan[ch_type] ) {
1270 ampOk = m_bitEneChan[ch_type][2];
1271 } else {
1272 ampOk = m_bitEneChan[ch_type][1];
1273 }
1274
1275 bool timeOk = false;
1276 if (time < m_minTimeChan[ch_type] ) {
1277 timeOk = m_bitTimeChan[ch_type][0];
1278 } else if (time > m_maxTimeChan[ch_type] ) {
1279 timeOk = m_bitTimeChan[ch_type][2];
1280 } else {
1281 timeOk = m_bitTimeChan[ch_type][1];
1282 }
1283
1284 if (ampOk && timeOk) {
1285
1286 ATH_MSG_VERBOSE(evtnum.str()
1287 << " chan " << std::left << std::setw(14) << m_tileHWID->to_string(adcId)
1288 << " ch_ene = " << amp << " ch_t = " << time);
1289
1290 m_chanSel[hash] = true;
1291
1292 if (amp < chmin) {
1293 chmin = amp;
1294 minCh = rawChannel;
1295 } else if (amp > chmax) {
1296 chmax = amp;
1297 maxCh = rawChannel;
1298 }
1299
1300 if (time<tcmin) {
1301 tcmin = time;
1302 tminCh = rawChannel;
1303 }
1304 else if (time>tcmax) {
1305 tcmax = time;
1306 tmaxCh = rawChannel;
1307 }
1308 }
1309 }
1310 }
1311 }
1312
1313 for (index = hashNext - 48; index < hashNext; ++index) {
1314 if ((m_chanSel[index] && rawdata) || someDQerrors) {
1315 ATH_MSG_VERBOSE(evtnum.str()
1316 << " drw " << drwname(rawChannelCollection->identify())
1317 << " nBadMB = " << nbadMB
1318 << " nBadDMU = " << nbad
1319 << " EvtBCID = " << m_evtBCID
1320 << " DSPBCID = " << rawChannelCollection->getFragDSPBCID()
1321 << " GlobCRC = " << rawChannelCollection->getFragGlobalCRC() << " " << GlobalCRCErr
1322 << " error = 0x" << std::hex << error
1323 << " FE_CRC = 0x" << rawChannelCollection->getFragFEChipMask() << " 0x" << FE_DMUmask
1324 << " ROD_CRC = 0x" << rawChannelCollection->getFragRODChipMask() << " 0x" << ROD_DMUmask
1325 << " BCIDErr = 0x" << rawChannelCollection->getFragBCID() << " 0x" << BCIDErr
1326 << " MemPar = 0x" << rawChannelCollection->getFragMemoryPar()
1327 << " HeadForm = 0x"<< rawChannelCollection->getFragHeaderBit()
1328 << " HeadPar = 0x" << rawChannelCollection->getFragHeaderPar()
1329 << " SampForm = 0x"<< rawChannelCollection->getFragSampleBit()
1330 << " SampPar = 0x" << rawChannelCollection->getFragSamplePar()
1331 << std::dec);
1332 break;
1333 }
1334 }
1335 }
1336
1337 if (nbadMBMax >= m_minBadMB) {
1338 ++m_dqerr;
1339 statusOk = true;
1340 ATH_MSG_DEBUG( nevtnum.str()
1341 << " drw " << drwname((collMBMax) ? collMBMax->identify() : 0)
1342 << " nBadMB = " << nbadMBMax
1343 << " accepted");
1344
1345 } else if (nbadMax >= m_minBadDMU && nbadMax <= m_maxBadDMU) {
1346 ++m_dqerr;
1347 statusOk = true;
1348 ATH_MSG_DEBUG( nevtnum.str()
1349 << " drw " << drwname((collMax)?collMax->identify():0)
1350 << " nBadDMU = " << nbadMax
1351 << " accepted");
1352 }
1353
1354 if (tminCh && tminCh != minCh && tminCh != maxCh) {
1355 ATH_MSG_DEBUG(nevtnum.str()
1356 << " chan " << std::left << std::setw(14) << m_tileHWID->to_string(tminCh->adc_HWID())
1357 << " ch_e = " << tminCh->amplitude()
1358 << " tmin =" << tminCh->time());
1359 }
1360
1361 if (tmaxCh && tmaxCh != minCh && tmaxCh != maxCh) {
1362 ATH_MSG_DEBUG(nevtnum.str()
1363 << " chan " << std::left << std::setw(14) << m_tileHWID->to_string(tmaxCh->adc_HWID())
1364 << " ch_e = " << tmaxCh->amplitude()
1365 << " tmax = " << tmaxCh->time());
1366 }
1367
1368 if (minCh) {
1369 if (!cellminCh) ++m_minChan;
1370 statusOk = true;
1371 const char * tit = (tminCh == minCh) ? " tmin = ": ((tmaxCh == minCh) ? " tmax = ": " t = ");
1372 ATH_MSG_DEBUG(nevtnum.str()
1373 << " chan " << std::left << std::setw(14) << m_tileHWID->to_string(minCh->adc_HWID())
1374 << " ch_emin = " << chmin
1375 << tit << minCh->time()
1376 << " accepted");
1377 }
1378 if (maxCh) {
1379 if (!cellmaxCh) ++m_maxChan;
1380 statusOk = true;
1381 const char * tit = (tminCh == maxCh) ? " tmin = ": ((tmaxCh == maxCh) ? " tmax = ": " t = ");
1382 ATH_MSG_DEBUG(nevtnum.str()
1383 << " chan " << std::left << std::setw(14) << m_tileHWID->to_string(maxCh->adc_HWID())
1384 << " ch_emax = " << chmax
1385 << tit << maxCh->time()
1386 << " accepted");
1387 }
1388 emptyBad = false;
1389 }
1390 }
1391
1392
1393 if (m_readDigits) {
1394
1395 // Pointer to a Tile digits container
1396 SG::ReadHandle<TileDigitsContainer> digitsContainer(m_digitsContainerKey);
1397
1398 if (!digitsContainer.isValid()) {
1399 ATH_MSG_WARNING("Unable to read TileDigitsContainer from EventStore, disable reading of this container");
1400 m_readDigits = false;
1401
1402 } else {
1403
1404 IdContext chan_context = m_tileHWID->channel_context();
1405 IdentifierHash hash;
1406 int index,pmt;
1407 int nConst = 0;
1408 int nJump = 0;
1409 int nDmuErr = 0;
1410 int nOverLG = 0;
1411 int nOverHG = 0;
1412 int nUnderLG = 0;
1413 int nUnderHG = 0;
1414
1415 for (const TileDigitsCollection * digitsCollection : *digitsContainer) {
1416
1417 int frag = digitsCollection->identify();
1418 bool eb = (frag > 0x2ff);
1419 bool ebsp = (frag == 0x30e || frag == 0x411);
1420
1421 int ros = frag >> 8;
1422 int drawer = frag & 0x3F;
1423 unsigned int drawerIdx = TileCalibUtils::getDrawerIdx(ros,drawer);
1424 if ( m_drawerToSkip[drawerIdx] ) continue;
1425
1426 int chMBTS = -1;
1427 if (eb) {
1428 for (int ch: {12,4,0}) {
1429 m_cabling->h2s_cell_id_index(ros,drawer,ch,index,pmt);
1430 if (index == -2) {
1431 chMBTS = ch;
1432 break;
1433 }
1434 }
1435 }
1436
1437 int nChBadDB = 0;
1438 int nChBadNC = 0;
1439 int nChBad = 0;
1440 int nChTot = 0;
1441 int nChDmu[16] = {0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0};
1442 int nChBadDmu[16] = {0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0};
1443
1444 for (const TileDigits* tile_digits : *digitsCollection) {
1445
1446 ++nChTot;
1447
1448 HWIdentifier adcId = tile_digits->adc_HWID();
1449 HWIdentifier chId = m_tileHWID->channel_id(adcId);
1450 m_tileHWID->get_hash(chId, hash, &chan_context);
1451 if ( m_chanToSkip[hash] ) continue;
1452 int channel = m_tileHWID->channel(adcId);
1453 int dmu = channel/3;
1454 ++nChDmu[dmu];
1455 int chEmpty = TileDQstatus::isChEmpty(ros, drawer, channel);
1456 bool isConnected = (chEmpty < 2);
1457 if (!isConnected) ++nChBadNC;
1458 int adc = m_tileHWID->adc(adcId);
1459 const char *cellname = "";
1460 int ch_type = 0;
1461 if (channel == chMBTS) {
1462 cellname = " MBTS";
1463 ch_type = 2;
1464 } else if ( (ebsp && (channel == 18 || channel == 19 || channel == 12 || channel == 13) )
1465 || (eb && (channel == 0 || channel == 1 || channel == 12 || channel == 13) ) ) {
1466 cellname = " GAP";
1467 ch_type = 1;
1468 } else if (chEmpty > 0) {
1469 cellname = " EMPTY";
1470 }
1471
1472 const char *badname = "";
1473 if (DQstatus && !DQstatus->isAdcDQgood(ros,drawer,channel,adc)) {
1474 badname = " BADDQ";
1475 if (isConnected) {
1476 ++nChBad;
1477 ++nChBadDmu[dmu];
1478 }
1479 } else if (m_checkDCS && m_tileDCS->getDCSStatus(ros, drawer, channel) > TileDCSState::WARNING) {
1480 badname = " BADDCS";
1481 } else if (m_tileBadChanTool->getAdcStatus(drawerIdx,channel,adc).isBad()) {
1482 badname = " BADDB";
1483 if (isConnected) {
1484 ++nChBadDB;
1485 }
1486 } else if (m_chanBad[hash]) {
1487 if (badFromCell) {
1488 if (ch_type != 2) badname = " BADQUAL";
1489 m_chanQua[hash] = 999;
1490 if (isConnected && ch_type!=2) {
1491 ++nChBad;
1492 ++nChBadDmu[dmu];
1493 }
1494 } else {
1495 badname = " BADUNKN"; // should never see this
1496 }
1497 } else if (badFromCell && m_chanEne[hash] == 0.0) {
1498 badname = " BADDIGI"; // temporary fix, (for May-2012 ESD), should never get this
1499 if (isConnected) {
1500 ++nChBad;
1501 ++nChBadDmu[dmu];
1502 }
1503 }
1504 const char *enename = " ene = ";
1505 const char *timename = " time = ";
1506 const char *qualname = " qual = ";
1507 if (badFromCell && badname[0] != 0) {
1508 enename = " BAD = ";
1509 if (m_chanDsp.size()) {
1510 qualname = " eDSP = ";
1512 }
1513 if (m_chanTDsp.size()) {
1514 timename = " tDSP = ";
1516 }
1517 }
1518
1519 char badnm[30];
1520 sprintf(badnm," BADDIGIEX%s",badname);
1521 float dmin,dmax;
1522
1523 std::vector<float> samples = tile_digits->samples();
1524 int nSamp = samples.size();
1525 if (nSamp > 6) {
1526
1527 bool useCh= !( (m_skipMBTS && channel == chMBTS) ||
1528 (m_skipEmpty && chEmpty > 0) );
1529 bool checkCh = !( (m_skipMasked && m_chanBad[hash]) ) && useCh && m_checkJumps;
1530
1531 int err = TileRawChannelBuilder::CorruptedData(ros,drawer,channel,adc,samples,dmin,dmax,m_ADCmaxMinusEps,m_ADCmaskValueMinusEps);
1532
1533 if (badname[0]==0) {
1534 if (err && err>-3) { // do not consider all zeros in empty samples as error
1535 if (isConnected || err != -2) {
1536 if (checkCh) m_chanSel[hash] = true;
1537 badname = badnm;
1538 ++nChBad;
1539 ++nChBadDmu[dmu];
1540 if (!isConnected) --nChBadNC;
1541 }
1542 if (err > 0) {
1543 if (err < 10) {
1544 badnm[9] = 48+err;
1545 } else if (err < 36) {
1546 badnm[9] = 55+err;
1547 }
1548 } else {
1549 badnm[9] = 48;
1550 }
1551 } else {
1552 // old error types, used for tests in August 2012
1553 // expect to see only warningE7 and warningE8 for gap/crack
1554 int warn = Are3FF(samples, adc, ch_type);
1555 if (warn) {
1556 if (checkCh) m_chanSel[hash] = true;
1557 sprintf(badnm," warningE%d%s",warn,badname);
1558 badname = badnm;
1559 }
1560 }
1561 }
1562
1563 if ((!err) // channel without bad patterns
1564 && (useCh) // normal connected channel
1565 && (badname[0] == 0 || badname[1] == 'w' // no digital error
1566 || (badname[4] == 'Q' && !m_skipMasked))) { // error from TileCell but it is ignored
1567
1568 if (adc) { // HG
1569
1570 if (dmax > m_overflowHG) {
1571 m_chanSel[hash] = true; // always print overflows
1572 if (m_checkOverHG){
1573 ++nOverHG;
1574 }
1575 }
1576 if (dmin < m_underflowHG) {
1577 m_chanSel[hash] = true; // always print underflows
1578 if (m_checkUnderHG){
1579 ++nUnderHG;
1580 }
1581 }
1582
1583 } else { // LG
1584
1585 if (dmax > m_overflowLG) {
1586 m_chanSel[hash] = true; // always print overflows
1587 if (m_checkOverLG){
1588 ++nOverLG;
1589 }
1590 }
1591 if (dmin < m_underflowLG) {
1592 m_chanSel[hash] = true; // always print underflows
1593 if (m_checkUnderLG){
1594 ++nUnderLG;
1595 }
1596 }
1597 }
1598 }
1599
1600 bool someSampErrors = false;
1601
1602 if (m_checkJumps && (checkCh || m_chanSel[hash])) {
1603
1604 float pedDelta = (adc ? m_pedDeltaHG : m_pedDeltaLG);
1605 float jumpDelta = (adc ? m_jumpDeltaHG : m_jumpDeltaLG);
1606
1607 float ped = samples[0];
1608 float dmin = ped;
1609 float dmax = ped;
1610 int nped = 1;
1611 int npedmax = 1;
1612 bool cnstPed = true;
1613 bool cnstPedmax = true;
1614 for (int i = 1; i < nSamp; ++i) {
1615 float smp = samples[i];
1616 float dped = smp - ped;
1617 if (fabs(dped) < pedDelta) {
1618 ++nped;
1619 if (dped != 0.0) {
1620 cnstPed = false;
1621 ped += dped/nped;
1622 }
1623 } else {
1624 if (nped>npedmax) {
1625 npedmax=nped;
1626 cnstPedmax = cnstPed;
1627 }
1628 cnstPed = true;
1629 ped = smp;
1630 nped = 1;
1631 }
1632 if (smp<dmin) {
1633 dmin = smp;
1634 } else if (smp>dmax) {
1635 dmax = smp;
1636 }
1637 }
1638 if (nped>npedmax) {
1639 npedmax=nped;
1640 cnstPedmax = cnstPed;
1641 }
1642
1643 if (dmax - dmin >= jumpDelta) {
1644 bool accEmin = (m_chanEne[hash]<m_minEneChan[ch_type]);
1645 bool accEmax = (m_chanEne[hash]>m_maxEneChan[ch_type]);
1646 bool accCnst = false;
1647 bool accJump = false;
1648 bool cnstMin = true;
1649 bool cnstMax = true;
1650 bool jumpNeg = false;
1651 bool jumpPos = false;
1652 bool jumpEnd = false;
1653 bool jumpZer = false;
1654 bool jumpOve = false;
1655 bool narrowUp = false;
1656 bool narrowDown = false;
1657 if (npedmax >= m_constLength && ((dmax-ped) >= jumpDelta || (ped-dmin) >= jumpDelta) ) {
1658 ++nConst;
1659 accCnst = true;
1660 }
1661 int nmin = 0;
1662 int nmax = 0;
1663 int pmin = -1;
1664 int pmax = -1;
1665 float abovemin = dmax;
1666 float belowmax = dmin;
1667 for (int i = 0; i < nSamp; ++i) {
1668 float smp = samples[i];
1669 if (smp - dmin < pedDelta) {
1670 ++nmin;
1671 pmin = i;
1672 if (smp != dmin) cnstMin = false;
1673 }
1674 if (dmax - smp < pedDelta) {
1675 ++nmax;
1676 pmax = i;
1677 if (smp != dmax) cnstMax = false;
1678 }
1679 if (smp < abovemin && smp > dmin) {
1680 abovemin = smp;
1681 }
1682 if (smp > belowmax && smp < dmax) {
1683 belowmax = smp;
1684 }
1685 }
1686 if (nmax + nmin == nSamp) {
1687 if (nmax > 1 && nmin > 1) {
1688 ++nJump;
1689 accJump = true;
1690 } else if (nmax == 1) {
1691 if (pmax < nSamp - 1) { // ignore jump in last sample
1692 ++nJump;
1693 accJump = true;
1694 jumpPos = true;
1695 cnstMax = false;
1696 }
1697 if (pmax == 0 || pmax == nSamp - 1) {
1698 jumpEnd = true;
1699 }
1700 } else if (nmin == 1) {
1701 ++nJump;
1702 accJump = true;
1703 jumpNeg = true;
1704 cnstMin = false;
1705 if (pmin == 0 || pmin == nSamp - 1) {
1706 jumpEnd = true;
1707 }
1708 }
1709 }
1710 if (dmin == 0.0) {
1711 if (!accJump) {
1712 ++nJump;
1713 accJump = true;
1714 cnstMin = false;
1715 cnstMax = false;
1716 }
1717 jumpZer = true;
1718 }
1719 if (dmax > m_ADCmaxMinusEps) {
1720 if (!accJump) {
1721 ++nJump;
1722 accJump = true;
1723 cnstMin = false;
1724 cnstMax = false;
1725 }
1726 jumpOve = true;
1727 }
1728 float secondMax = (dmax-dmin)*m_secondMaxLevel;
1729 if (pmax > 0 && pmax < nSamp-1 && std::max(samples[pmax-1], samples[pmax+1]) < dmin+secondMax) {
1730 if (!accJump) {
1731 ++nJump;
1732 accJump = true;
1733 cnstMax = false;
1734 if (nmin + nmax != nSamp) {
1735 cnstMin = false;
1736 }
1737 }
1738 narrowUp = true;
1739 }
1740 if (pmin > 0 && pmin < nSamp - 1 && std::min(samples[pmin - 1], samples[pmin + 1]) > dmax - secondMax) {
1741 if (!accJump) {
1742 ++nJump;
1743 accJump = true;
1744 cnstMin = false;
1745 if (nmin + nmax != nSamp) {
1746 cnstMax = false;
1747 }
1748 }
1749 narrowDown = true;
1750 }
1751
1752 if (accEmin || accEmax || accCnst || accJump) {
1753 someSampErrors = true;
1754 ATH_MSG_VERBOSE (evtnum.str()
1755 << " chan " << std::left << std::setw(14) << m_tileHWID->to_string(adcId)
1756 << enename << m_chanEne[hash] << " samp = " << samples[0]
1757 << " " << samples[1] << " " << samples[2] << " " << samples[3]
1758 << " " << samples[4] << " " << samples[5] << " " << samples[6]
1759 << timename << m_chanTime[hash]
1760 << qualname << m_chanQua[hash]
1761 << cellname << badname
1762 << ((accEmin) ? " neg_e" : "")
1763 << ((accEmax) ? " pos_e" : "")
1764 << ((accCnst) ? " const" : "")
1765 << ((accCnst&&cnstPedmax) ? "Const" : "")
1766 << ((accJump) ? " jump" : "")
1767 << ((accJump&&jumpZer) ? "Zero" : "")
1768 << ((accJump&&jumpOve) ? "Over" : "")
1769 << ((accJump&&jumpPos) ? "SingleUp" : ((narrowUp) ? "NarrowUp" : "") )
1770 << ((accJump&&jumpNeg) ? "SingleDown" : ((narrowDown) ? "NarrowDown" : "") )
1771 << ((accJump&&jumpEnd) ? "AtEdge" : "")
1772 << ((accJump&&cnstMin) ? "ConstMin" : "")
1773 << ((accJump&&cnstMax) ? "ConstMax" : "")
1774 << " " << dmax-dmin);
1775 }
1776 }
1777 }
1778
1779 if (someSampErrors) {
1780 m_chanSel[hash] = true;
1781 } else if (m_chanSel[hash]) {
1782 bool accEmin = (m_chanEne[hash] < m_minEneChan[ch_type]);
1783 bool accEmax = (m_chanEne[hash] > m_maxEneChan[ch_type]);
1784 bool jumpOve = (dmax>m_ADCmaxMinusEps);
1785 bool jumpZer = (dmin < 0.01);
1786 ATH_MSG_VERBOSE(evtnum.str()
1787 << " chan " << std::left << std::setw(14) << m_tileHWID->to_string(adcId)
1788 << enename << m_chanEne[hash] << " samp = " << samples[0]
1789 << " " << samples[1] << " " << samples[2] << " " << samples[3]
1790 << " " << samples[4] << " " << samples[5] << " " << samples[6]
1791 << timename << m_chanTime[hash]
1792 << qualname << m_chanQua[hash]
1793 << cellname << badname
1794 << ((accEmin) ? " neg_e" : "")
1795 << ((accEmax) ? " pos_e" : "")
1796 <<((jumpZer) ? " underflow" : "")
1797 <<((jumpOve) ? " overflow" : "") );
1798 }
1799
1800 }
1801 }
1802 if (m_checkDMUs && nChBad > 1 && nChBad + nChBadDB + nChBadNC < nChTot) {
1803 int nChBad1 = 0;
1804 int nChBad2 = 0;
1805 int nDmuBad1 = 0;
1806 int nDmuBad2 = 0;
1807 int nDmuBad = 0;
1808 int nDmuTot = 0;
1809 bool has23 = false;
1810 for (int dmu = 0; dmu < 16; ++dmu) {
1811 if (nChDmu[dmu] > 0) {
1812 ++nDmuTot;
1813 if (nChBadDmu[dmu] > 0) {
1814 ++nDmuBad;
1815 if (dmu < 8) {
1816 nChBad1 += nChBadDmu[dmu];
1817 ++nDmuBad1;
1818 } else {
1819 nChBad2 += nChBadDmu[dmu];
1820 ++nDmuBad2;
1821 }
1822 if (nChBadDmu[dmu] == 2 && nChDmu[dmu] == 3) has23 = true;
1823 }
1824 }
1825 }
1826 if (nDmuBad == 1 /* && (!has23) */ ) continue;
1827 if (nDmuBad == 2 && nChBad < 3) continue;
1828
1829 if (nDmuBad>2 || nChBad > 9 || nChTot > 19 || has23) {
1830
1831 ++nDmuErr;
1832
1833 for (const TileDigits* tile_digits : *digitsCollection) {
1834
1835 HWIdentifier adcId = tile_digits->adc_HWID();
1836 HWIdentifier chId = m_tileHWID->channel_id(adcId);
1837 m_tileHWID->get_hash(chId, hash, &chan_context);
1838 if ( m_chanToSkip[hash] ) continue;
1839 std::vector<float> samples = tile_digits->samples();
1840
1841 if (!m_chanSel[hash] && samples.size()>6) {
1842 int channel = m_tileHWID->channel(adcId);
1843 int adc = m_tileHWID->adc(adcId);
1844 int chEmpty = TileDQstatus::isChEmpty(ros, drawer, channel);
1845 const char *cellname = "";
1846 int ch_type = 0;
1847 if (channel == chMBTS) {
1848 cellname = " MBTS";
1849 ch_type = 2;
1850 } else if ( (ebsp && (channel == 18 || channel == 19 || channel == 12 || channel == 13) )
1851 || (eb && (channel == 0 || channel == 1 || channel == 12 || channel == 13) ) ) {
1852 cellname = " GAP";
1853 ch_type = 1;
1854 } else if (chEmpty > 0) {
1855 cellname = " EMPTY";
1856 }
1857 const char *badname = "";
1858 if (m_tileBadChanTool->getAdcStatus(drawerIdx, channel, adc).isBad()) {
1859 badname = " BADDB";
1860 } else if (DQstatus && !DQstatus->isAdcDQgood(ros, drawer, channel, adc)) {
1861 badname = " BADDQ";
1862 } else if (m_checkDCS && m_tileDCS->getDCSStatus(ros, drawer, channel) > TileDCSState::WARNING) {
1863 badname = " BADDCS";
1864 } else if (m_chanBad[hash]) {
1865 if (badFromCell) {
1866 if (ch_type != 2) badname = " BADQUAL";
1867 m_chanQua[hash] = 999;
1868 } else {
1869 badname = " BADUNKN"; // should never see this
1870 }
1871 } else if (badFromCell && m_chanEne[hash] == 0.0) {
1872 badname = " BADDIGI"; // temporary fix, (for May-2012 ESD), should never get this
1873 }
1874
1875 char badnm[30];
1876 sprintf(badnm, " BADDIGIEX%s", badname);
1877 float dmin, dmax;
1878
1879 int err = TileRawChannelBuilder::CorruptedData(ros, drawer,
1880 channel, adc, samples, dmin, dmax, m_ADCmaxMinusEps,m_ADCmaskValueMinusEps);
1881 if (err) {
1882 bool isConnected = (chEmpty < 2);
1883 if (isConnected || err != -2) {
1884 badname = badnm;
1885 }
1886 if (err > 0) {
1887 if (err < 10) {
1888 badnm[9] = 48 + err;
1889 } else if (err < 36) {
1890 badnm[9] = 55 + err;
1891 }
1892 } else {
1893 badnm[9] = 48;
1894 }
1895 }
1896
1897 const char *enename = " ene = ";
1898 const char *timename = " time = ";
1899 const char *qualname = " qual = ";
1900 if (badFromCell && badname[0] != 0) {
1901 enename = " BAD = ";
1902 if (m_chanDsp.size()) {
1903 qualname = " eDSP = ";
1905 }
1906 if (m_chanTDsp.size()) {
1907 timename = " tDSP = ";
1909 }
1910 }
1911
1912 bool accEmin = (m_chanEne[hash]<m_minEneChan[ch_type]);
1913 bool accEmax = (m_chanEne[hash]>m_maxEneChan[ch_type]);
1914 bool jumpOve = (dmax > m_ADCmaxMinusEps);
1915 bool jumpZer = (dmin < 0.01);
1916
1917 ATH_MSG_VERBOSE(evtnum.str()
1918 << " chan " << std::left << std::setw(14) << m_tileHWID->to_string(adcId)
1919 << enename << m_chanEne[hash] << " samp = " << samples[0]
1920 << " " << samples[1] << " " << samples[2] << " " << samples[3]
1921 << " " << samples[4] << " " << samples[5] << " " << samples[6]
1922 << timename << m_chanTime[hash]
1923 << qualname << m_chanQua[hash]
1924 << cellname << badname
1925 << ((accEmin) ? " neg_e" : "")
1926 << ((accEmax) ? " pos_e" : "")
1927 << ((jumpZer) ? " underflow" : "")
1928 << ((jumpOve) ? " overflow" : "") );
1929
1930 }
1931 }
1932 }
1933
1934 std::ostringstream badstr;
1935 badstr << " ch: " << nChBad1 << " + " << nChBad2 << " = " << nChBad << " / " << nChTot << " = " << 100*nChBad/nChTot
1936 << " % dmu: " << nDmuBad1 << " + " << nDmuBad2 << " = " << nDmuBad << " / " << nDmuTot << " ";
1937 for (int dmu=0; dmu<16; ++dmu) {
1938 if (nChDmu[dmu]>0) {
1939 badstr << " " << std::hex << dmu << "=" << nChBadDmu[dmu] << "/" << nChDmu[dmu];
1940 }
1941 }
1942 ATH_MSG_VERBOSE (evtnum.str()
1943 << " drw " << drwname(digitsCollection->identify()) << badstr.str());
1944 }
1945 }
1946
1947 if (nConst) {
1948 ++m_const;
1949 statusOk = true;
1950 ATH_MSG_DEBUG( nevtnum.str()
1951 << " n_const_sample_errors = " << nConst
1952 << " accepted");
1953 }
1954 if (nJump) {
1955 ++m_jump;
1956 statusOk = true;
1957 ATH_MSG_DEBUG( nevtnum.str()
1958 << " n_jump_sample_errors = " << nJump
1959 << " accepted");
1960 }
1961 if (nOverLG) {
1962 ++m_overLG;
1963 statusOk = true;
1964 ATH_MSG_DEBUG( nevtnum.str()
1965 << " n_overflow_LG = " << nOverLG
1966 << " accepted");
1967 }
1968 if (nOverHG) {
1969 ++m_overHG;
1970 statusOk = true;
1971 ATH_MSG_DEBUG(nevtnum.str()
1972 << " n_overflow_HG = " << nOverHG
1973 << " accepted");
1974 }
1975 if (nUnderLG) {
1976 ++m_underLG;
1977 statusOk = true;
1978 ATH_MSG_DEBUG( nevtnum.str()
1979 << " n_underflow_LG = " << nUnderLG
1980 << " accepted");
1981 }
1982 if (nUnderHG) {
1983 ++m_underHG;
1984 statusOk = true;
1985 ATH_MSG_DEBUG(nevtnum.str()
1986 << " n_underflow_HG = " << nUnderHG
1987 << " accepted");
1988 }
1989 if (nDmuErr) {
1990 ++m_dmuerr;
1991 statusOk = true;
1992 ATH_MSG_DEBUG( nevtnum.str()
1993 << " n_DMU_errors = " << nDmuErr
1994 << " accepted");
1995 }
1996 }
1997 }
1998
1999 if (m_printOnly)
2000 this->setFilterPassed (false);
2001 else
2002 this->setFilterPassed (statusOk);
2003
2004 if (statusOk) {
2005 ++m_accept;
2006 //ATH_MSG_VERBOSE (nevtnum.str() << " accepted");
2007 } else {
2008 //ATH_MSG_VERBOSE (nevtnum.str() << " rejected");
2009 }
2010
2011 return StatusCode::SUCCESS;
2012}
const boost::regex rr(r_r)
#define ATH_MSG_VERBOSE(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
static std::string drwname(int id)
static const Attributes_t empty
float time() const
get time (data member)
Definition CaloCell.h:368
double energy() const
get energy (data member)
Definition CaloCell.h:327
Identifier ID() const
get ID (from cached data member) non-virtual and inline for fast access
Definition CaloCell.h:295
IdentifierHash onl2() const
cell online identifier 2
IdentifierHash onl1() const
cell online identifier 1
static unsigned int getDrawerIdx(unsigned int ros, unsigned int drawer)
Returns a drawer hash.
std::vector< float > m_chanEne
bool m_bitTimeCell[ptnlength]
int Are3FF(std::vector< float > &OptFilterDigits, int OptFilterGain, int ch_type)
std::vector< float > m_chanDsp
bool m_bitEneChan[3][ptnlength]
ToolHandle< ITileBadChanTool > m_tileBadChanTool
SG::ReadHandleKey< TileDQstatus > m_dqStatusKey
std::vector< bool > m_chanSel
std::vector< bool > m_chanToSkip
ToolHandle< ITileDCSTool > m_tileDCS
SG::ReadHandleKey< TileRawChannelContainer > m_rawChannelContainerKey
std::vector< bool > m_chanBad
bool m_bitTimeChan[3][ptnlength]
SG::ReadHandleKey< TileDigitsContainer > m_digitsContainerKey
SG::ReadHandleKey< CaloCellContainer > m_cellContainerKey
SG::ReadHandleKey< xAOD::EventInfo > m_eventInfoKey
std::vector< int > m_nDrawerOff
std::vector< float > m_chanTime
std::vector< bool > m_drawerToSkip
std::vector< float > m_chanTDsp
std::vector< float > m_chanQua
bool m_bitEneCell[ptnlength]
uint8_t qual1(void) const
get quality of first PMT (data member)
Definition TileCell.h:197
float time1(void) const
get time of first PMT
Definition TileCell.h:192
int gain2(void) const
get gain of second PMT
Definition TileCell.cxx:175
bool badch1(void) const
check if first PMT is in bad channel list and masked
Definition TileCell.h:209
uint8_t qbit2(void) const
get quality bits of second PMT (data member)
Definition TileCell.h:206
int gain1(void) const
get gain of first PMT
Definition TileCell.cxx:168
uint8_t qual2(void) const
get quality of second PMT (data member)
Definition TileCell.h:200
float ene1(void) const
get energy of first PMT
Definition TileCell.h:187
bool badch2(void) const
check if second PMT is in bad channel list and masked
Definition TileCell.h:212
float time2(void) const
get time of second PMT
Definition TileCell.h:194
@ MASK_OVER
Definition TileCell.h:64
@ MASK_CMPC
Definition TileCell.h:66
uint8_t qbit1(void) const
get quality bits of first PMT (data member)
Definition TileCell.h:203
float ene2(void) const
get energy of second PMT
Definition TileCell.h:189
static int isChEmpty(int partition, int drawer, int ch)
True if channel is not fully implemented.
@ NOT_VALID_HASH
Definition TileHWID.h:314
static int CorruptedData(int ros, int drawer, int channel, int gain, const std::vector< float > &digits, float &dmin, float &dmax, float ADCmaxMinusEps, float ADCmaskValueMinusEps)
float time(int ind=0) const
float amplitude(int ind=0) const
HWIdentifier adc_HWID(void) const
Definition TileRawData.h:53
@ Tile
The Tile calorimeter.
@ Warning
The sub-detector issued a warning.
@ Error
The sub-detector issued an error.
time(flags, cells_name, *args, **kw)
str index
Definition DeMoScan.py:362
SG::ReadCondHandle< T > makeHandle(const SG::ReadCondHandleKey< T > &key, const EventContext &ctx=Gaudi::Hive::currentContext())
setEventNumber uint32_t

◆ extraDeps_update_handler()

void AthCommonDataStore< AthCommonMsg< Algorithm > >::extraDeps_update_handler ( Gaudi::Details::PropertyBase & ExtraDeps)
protectedinherited

Add StoreName to extra input/output deps as needed.

use the logic of the VarHandleKey to parse the DataObjID keys supplied via the ExtraInputs and ExtraOuputs Properties to add the StoreName if it's not explicitly given

◆ extraOutputDeps()

const DataObjIDColl & AthAlgorithm::extraOutputDeps ( ) const
overridevirtualinherited

Return the list of extra output dependencies.

This list is extended to include symlinks implied by inheritance relations.

Definition at line 50 of file AthAlgorithm.cxx.

51{
52 // If we didn't find any symlinks to add, just return the collection
53 // from the base class. Otherwise, return the extended collection.
54 if (!m_extendedExtraObjects.empty()) {
56 }
57 return Algorithm::extraOutputDeps();
58}
DataObjIDColl m_extendedExtraObjects

◆ finalize()

StatusCode TileCellSelector::finalize ( )
overridevirtual

Definition at line 2092 of file TileCellSelector.cxx.

2092 {
2093 //ATH_MSG_DEBUG ("finalize()");
2094
2095 ATH_MSG_INFO ("Processed " << m_counter << " events, accepted " << m_accept
2096 << " (" << m_minCell << "/" << m_minChan
2097 << "," << m_maxCell << "/" << m_maxChan
2098 << "," << m_jump << "/" << m_const
2099 << "," << m_overLG << "/" << m_overHG
2100 << "," << m_underLG << "/" << m_underHG
2101 << "," << m_dqerr << "/" << m_dmuerr
2102 << "," << m_warnerr << ") events.");
2103
2104 //ATH_MSG_DEBUG ("finalize() successful");
2105
2106 return StatusCode::SUCCESS;
2107}
#define ATH_MSG_INFO(x)

◆ initialize()

StatusCode TileCellSelector::initialize ( )
overridevirtual

Definition at line 152 of file TileCellSelector.cxx.

152 {
153 //ATH_MSG_DEBUG("in initialize()");
154
155 ATH_CHECK( m_eventInfoKey.initialize() );
156
157 ATH_CHECK(detStore()->retrieve(m_tileID, "TileID"));
158 ATH_CHECK(detStore()->retrieve(m_tileHWID, "TileHWID"));
159
161 m_ADCmaxMinusEps = m_tileInfo->ADCmax() - 0.01;
162 m_ADCmaskValueMinusEps = m_tileInfo->ADCmaskValue() - 0.01; // indicates channels which were masked in background dataset
163
165
166 ATH_CHECK(m_tileBadChanTool.retrieve());
167
168 if (m_checkDCS) {
169 ATH_CHECK(m_tileDCS.retrieve());
170 } else {
171 m_tileDCS.disable();
172 }
173
174 ATH_MSG_INFO( "Cell container "
175 << ((m_cellContainerKey.key().empty()) ? "NOT SET" : m_cellContainerKey.key()) );
176 ATH_MSG_INFO( "Digits container "
177 << ((m_digitsContainerKey.key().empty()) ? "NOT SET" : m_digitsContainerKey.key()));
178 ATH_MSG_INFO( "RawChannel container "
179 << ((m_rawChannelContainerKey.key().empty()) ? "NOT SET" : m_rawChannelContainerKey.key()));
180
181 ATH_MSG_INFO( "CheckJumps " << ((m_checkJumps) ? "true" : "false"));
182 ATH_MSG_INFO( "CheckDMUs " << ((m_checkDMUs) ? "true" : "false"));
183 ATH_MSG_INFO( "CheckOverLG " << ((m_checkOverLG) ? "true" : "false"));
184 ATH_MSG_INFO( "CheckOverHG " << ((m_checkOverHG) ? "true" : "false"));
185 ATH_MSG_INFO( "CheckUnderLG " << ((m_checkUnderLG) ? "true" : "false"));
186 ATH_MSG_INFO( "CheckUnderHG " << ((m_checkUnderHG) ? "true" : "false"));
187
188 if (m_overflowLG < 0) m_overflowLG += m_tileInfo->ADCmax();
189 if (m_overflowHG < 0) m_overflowHG += m_tileInfo->ADCmax();
190 ATH_MSG_INFO( "OverflowLG " << m_overflowLG);
191 ATH_MSG_INFO( "OverflowHG " << m_overflowHG);
192 ATH_MSG_INFO( "UnderflowLG " << m_underflowLG);
193 ATH_MSG_INFO( "UnderflowHG " << m_underflowHG);
194
195 ATH_MSG_INFO( "SkipEmpty " << ((m_skipEmpty) ? "true" : "false"));
196 ATH_MSG_INFO( "SkipMasked " << ((m_skipMasked) ? "true" : "false"));
197 ATH_MSG_INFO( "SkipMBTS " << ((m_skipMBTS) ? "true" : "false"));
198 ATH_MSG_INFO( "CheckDCS " << ((m_checkDCS) ? "true" : "false"));
199
200
201 m_readCells = !m_cellContainerKey.key().empty();
203 m_readDigits = !m_digitsContainerKey.key().empty();
204
205 if (m_readCells) {
206 ATH_MSG_INFO( "MinEnergyCell < " << m_minEneCell);
207 ATH_MSG_INFO( "MaxEnergyCell > " << m_maxEneCell);
208 ATH_MSG_INFO( "PtnEnergyCell = " << m_ptnEneCell);
209 ATH_MSG_INFO( "MinTimeCell < " << m_minTimeCell);
210 ATH_MSG_INFO( "MaxTimeCell > " << m_maxTimeCell);
211 ATH_MSG_INFO( "PtnTimeCell = " << m_ptnTimeCell);
212
213 ATH_CHECK( m_cellContainerKey.initialize() );
214 }
215
217 ATH_MSG_INFO( "MinEnergyChan < " << m_minEneChan[0]);
218 ATH_MSG_INFO( "MaxEnergyChan > " << m_maxEneChan[0]);
219 ATH_MSG_INFO( "PtnEnergyChan = " << m_ptnEneChan[0]);
220 ATH_MSG_INFO( "MinEnergyGap < " << m_minEneChan[1]);
221 ATH_MSG_INFO( "MaxEnergyGap > " << m_maxEneChan[1]);
222 ATH_MSG_INFO( "PtnEnergyGap = " << m_ptnEneChan[1]);
223 ATH_MSG_INFO( "MinTimeChan < " << m_minTimeChan[0]);
224 ATH_MSG_INFO( "MaxTimeChan > " << m_maxTimeChan[0]);
225 ATH_MSG_INFO( "PtnTimeChan = " << m_ptnTimeChan[0]);
226 ATH_MSG_INFO( "MinTimeGap < " << m_minTimeChan[1]);
227 ATH_MSG_INFO( "MaxTimeGap > " << m_maxTimeChan[1]);
228 ATH_MSG_INFO( "PtnTimeGap = " << m_ptnTimeChan[1]);
229 }
230
231 if (m_readRawChannels) {
232 ATH_MSG_INFO( "MinEnergyMBTS < " << m_minEneChan[2]);
233 ATH_MSG_INFO( "MaxEnergyMBTS > " << m_maxEneChan[2]);
234 ATH_MSG_INFO( "PtnEnergyMBTS = " << m_ptnEneChan[2]);
235 ATH_MSG_INFO( "MinTimeMBTS < " << m_minTimeChan[2]);
236 ATH_MSG_INFO( "MaxTimeMBTS > " << m_maxTimeChan[2]);
237 ATH_MSG_INFO( "PtnTimeMBTS = " << m_ptnTimeChan[2]);
238
239 ATH_CHECK( m_rawChannelContainerKey.initialize() );
240 }
241
242 switch (m_selectGain) {
243 case 0:
244 ATH_MSG_INFO( "Select Low gain channels only");
247 break;
248 case 1:
249 ATH_MSG_INFO( "Select High gain channels only");
252 break;
253 default:
254 ATH_MSG_INFO( "Select both gains");
255 break;
256 }
257
258 if (!m_digitsContainerKey.key().empty()) {
259 if (m_checkJumps) {
260 ATH_MSG_INFO( "JumpDeltaHG " << m_jumpDeltaHG);
261 ATH_MSG_INFO( "JumpDeltaLG " << m_jumpDeltaLG);
262 ATH_MSG_INFO( "PedDetlaHG " << m_pedDeltaHG);
263 ATH_MSG_INFO( "PedDetlaLG " << m_pedDeltaLG);
264 ATH_MSG_INFO( "ConstLength " << m_constLength);
265 }
266
267 ATH_CHECK( m_digitsContainerKey.initialize() );
268
269 }
270
271 if (!m_rawChannelContainerKey.key().empty()) {
272 if (m_checkDMUs) {
273 ATH_MSG_INFO( "MinBadDMU " << m_minBadDMU);
274 ATH_MSG_INFO( "MaxBadDMU " << m_maxBadDMU);
275 ATH_MSG_INFO( "MinBadMB " << m_minBadMB);
276 } else {
277 m_minBadDMU = 99;
278 m_maxBadDMU = -1;
279 m_minBadMB = 99;
280 }
281 ATH_MSG_INFO( "MaxVerboseCnt " << m_maxVerboseCnt);
282 }
283
284 ATH_MSG_INFO( "CheckWarning " << ((m_checkWarning)? "true" : "false"));
285 ATH_MSG_INFO( "CheckError " << ((m_checkError) ? "true" : "false"));
286 ATH_MSG_INFO( "PrintOnly " << ((m_printOnly) ? "true" : "false"));
287
288 if (m_drawer.size()>0) {
289 msg(MSG::INFO) << "Drawers which will be always printed:" << MSG::hex;
290 for (int frag: m_drawer) msg(MSG::INFO) << " 0x" << frag;
291 msg(MSG::INFO) << MSG::dec << endmsg;
292 }
293
294 if (m_drawerToCheck.size()>0) {
295 msg(MSG::INFO) << "Only those drawers will be checked:" << MSG::hex;
296 for (int frag: m_drawerToCheck) msg(MSG::INFO) << " 0x" << frag;
297 msg(MSG::INFO) << MSG::dec << endmsg;
298 }
299
300 if (m_chanToCheck.size()>0) {
301 msg(MSG::INFO) << "Only those channels will be checked:";
302 for (int ch: m_chanToCheck) msg(MSG::INFO) << " " << ch;
303 msg(MSG::INFO) << endmsg;
304 }
305
307
308 //ATH_MSG_DEBUG ("initialize() successful");
309
310 // convert patterns fo bolean arrays
311
312 int digit=1;
313 int scale=10;
314 for (int i=0; i<ptnlength; ++i) {
315
316 int bit=(m_ptnEneCell/digit)%scale;
317 m_bitEneCell[i] = (bit!=0);
318
319 bit=(m_ptnTimeCell/digit)%scale;
320 m_bitTimeCell[i] = (bit!=0);
321
322 for (int j=0; j<3; ++j) {
323
324 int bit=(m_ptnEneChan[j]/digit)%scale;
325 m_bitEneChan[j][i] = (bit!=0);
326
327 bit=(m_ptnTimeChan[j]/digit)%scale;
328 m_bitTimeChan[j][i] = (bit!=0);
329 }
330
331 digit *= scale;
332 }
333
334 m_chanToSkip.clear();
335 m_drawerToSkip.clear();
336 if (m_drawerToCheck.size()>0 || m_chanToCheck.size()>0 ) {
337 if (m_drawerToCheck.size()==0) {
338 m_drawerToCheck.resize(256);
339 auto itr = m_drawerToCheck.begin();
340 for (int frag : {0x100,0x200,0x300,0x400}) {
341 auto itr1 = itr+64;
342 std::iota(itr, itr1, frag);
343 itr = itr1;
344 }
345 } else if (m_chanToCheck.size()==0) {
346 m_chanToCheck.resize(48);
347 std::iota(m_chanToCheck.begin(), m_chanToCheck.end(), 0);
348 }
350 m_drawerToSkip.resize(1+TileCalibUtils::getDrawerIdx(4,63),true);
351 IdContext chan_context = m_tileHWID->channel_context();
352 IdentifierHash hash;
353 for (int frag : m_drawerToCheck) {
354 int ros = frag >> 8;
355 int drawer = frag & 0x3F;
356 unsigned int drawerIdx = TileCalibUtils::getDrawerIdx(ros,drawer);
357 m_drawerToSkip[drawerIdx] = false;
358 HWIdentifier ch_id = m_tileHWID->channel_id(ros,drawer,0);
359 m_tileHWID->get_hash(ch_id, hash, &chan_context);
360 for (int chan: m_chanToCheck)
361 m_chanToSkip[hash+chan] = false;
362 }
363 } else {
366 m_drawerToSkip.resize(1+TileCalibUtils::getDrawerIdx(4,63),false);
367 }
368
369 ATH_CHECK( m_dqStatusKey.initialize() );
370
371 return StatusCode::SUCCESS;
372}
#define endmsg
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ptnlength
const ServiceHandle< StoreGateSvc > & detStore() const
MsgStream & msg() const
static const TileCablingService * getInstance()
get pointer to service instance
retrieve(aClass, aKey=None)
Definition PyKernel.py:110

◆ inputHandles()

virtual std::vector< Gaudi::DataHandle * > AthCommonDataStore< AthCommonMsg< Algorithm > >::inputHandles ( ) const
overridevirtualinherited

Return this algorithm's input handles.

We override this to include handle instances from key arrays if they have not yet been declared. See comments on updateVHKA.

◆ msg()

MsgStream & AthCommonMsg< Algorithm >::msg ( ) const
inlineinherited

Definition at line 24 of file AthCommonMsg.h.

24 {
25 return this->msgStream();
26 }

◆ msgLvl()

bool AthCommonMsg< Algorithm >::msgLvl ( const MSG::Level lvl) const
inlineinherited

Definition at line 30 of file AthCommonMsg.h.

30 {
31 return this->msgLevel(lvl);
32 }

◆ outputHandles()

virtual std::vector< Gaudi::DataHandle * > AthCommonDataStore< AthCommonMsg< Algorithm > >::outputHandles ( ) const
overridevirtualinherited

Return this algorithm's output handles.

We override this to include handle instances from key arrays if they have not yet been declared. See comments on updateVHKA.

◆ printCell()

void TileCellSelector::printCell ( const TileCell * cell)
private

Definition at line 2016 of file TileCellSelector.cxx.

2016 {
2017 if (tile_cell==0) return;
2018
2019 // int drw = 0; // drawer number, range 0-63, the same for both channels
2020 // int ch1 = -1, ch2 = -1; // channel number, range 0-47 or -1 for unknown
2021 int ros1 = 0, ros2 = 0;
2022
2023 const CaloDetDescrElement * caloDDE = tile_cell->caloDDE();
2024 //Identifier id = tile_cell->ID();
2025
2026 IdentifierHash hash1 = caloDDE->onl1();
2027 if (hash1 != TileHWID::NOT_VALID_HASH) {
2028 HWIdentifier hw1 = m_tileHWID->channel_id(hash1);
2029 //ch1 = m_tileHWID->channel(hw1);
2030 //drw = m_tileHWID->drawer(hw1);
2031 ros1 = m_tileHWID->ros(hw1);
2032 }
2033
2034 IdentifierHash hash2 = caloDDE->onl2();
2035 if (hash2 != TileHWID::NOT_VALID_HASH) {
2036 HWIdentifier hw2 = m_tileHWID->channel_id(hash2);
2037 //ch2 = m_tileHWID->channel(hw2);
2038 //drw = m_tileHWID->drawer(hw2);
2039 ros2 = m_tileHWID->ros(hw2);
2040 if (ros1 == 0) ros1=ros2;
2041 } else {
2042 ros2 = ros1;
2043 }
2044
2045 // something is wrong
2046 if (ros1 == 0) return;
2047 /*
2048 int module = m_tileID->module(id);
2049
2050 int samp = m_tileID->sample(id);
2051
2052 bool single_PMT_scin = (samp == TileID::SAMP_E);
2053 bool single_PMT_C10 = (m_tileID->section(id) == TileID::GAPDET &&
2054 samp == TileID::SAMP_C &&
2055 (! m_cabling->C10_connected(m_tileID->module(id))) );
2056
2057 // distinguish cells with one or two PMTs
2058 bool single_PMT = single_PMT_C10 || single_PMT_scin;
2059
2060 // distinguish normal cells and fantoms (e.g. non-existing D4 in EBA15, EBC18
2061 // or non-existing E3/E4 - they might appear in CaloCellContainer)
2062 bool real_cell = single_PMT_C10 || m_cabling->TileGap_connected(id);
2063
2064 // note that in single PMT cell both badch1() and badch2() are changed together
2065 bool badch1 = (tile_cell->badch1());
2066 bool badch2 = (tile_cell->badch2());
2067
2068 // 0 = both PMTs are good; 1= 1 PMT is bad; 2= both PMTs are bad, or PMT is bad for single PMT cell
2069 int cell_isbad = (int)badch1 + (int)badch2;
2070
2071 int gn1 = tile_cell->gain1(); // gain of first PMT
2072 int gn2 = tile_cell->gain2(); // gain of second PMT
2073
2074 bool ch1Ok = (ch1>-1 && gn1 != CaloGain::INVALIDGAIN);
2075 bool ch2Ok = (ch2>-1 && gn2 != CaloGain::INVALIDGAIN);
2076
2077 // get the cell energy, time and position info
2078 double energy = tile_cell->energy();
2079 double time = tile_cell->time();
2080 double eta = tile_cell->eta();
2081 double phi = tile_cell->phi();
2082 double ene1 = tile_cell->ene1();
2083 double ene2 = tile_cell->ene2();
2084 double ediff = (single_PMT) ? 0.0 : tile_cell->eneDiff();
2085 double eratio = (energy!=0.0) ? ediff/energy : 0.0;
2086 double t1 = tile_cell->time1();
2087 double t2 = tile_cell->time2();
2088 double tdiff = (single_PMT) ? 0.0 : tile_cell->timeDiff();
2089 */
2090}

◆ renounce()

std::enable_if_t< std::is_void_v< std::result_of_t< decltype(&T::renounce)(T)> > &&!std::is_base_of_v< SG::VarHandleKeyArray, T > &&std::is_base_of_v< Gaudi::DataHandle, T >, void > AthCommonDataStore< AthCommonMsg< Algorithm > >::renounce ( T & h)
inlineprotectedinherited

Definition at line 380 of file AthCommonDataStore.h.

381 {
382 h.renounce();
384 }
std::enable_if_t< std::is_void_v< std::result_of_t< decltype(&T::renounce)(T)> > &&!std::is_base_of_v< SG::VarHandleKeyArray, T > &&std::is_base_of_v< Gaudi::DataHandle, T >, void > renounce(T &h)

◆ renounceArray()

void AthCommonDataStore< AthCommonMsg< Algorithm > >::renounceArray ( SG::VarHandleKeyArray & handlesArray)
inlineprotectedinherited

remove all handles from I/O resolution

Definition at line 364 of file AthCommonDataStore.h.

364 {
366 }

◆ sysInitialize()

StatusCode AthAlgorithm::sysInitialize ( )
overridevirtualinherited

Override sysInitialize.

Override sysInitialize from the base class.

Loop through all output handles, and if they're WriteCondHandles, automatically register them and this Algorithm with the CondSvc

Scan through all outputHandles, and if they're WriteCondHandles, register them with the CondSvc

Reimplemented from AthCommonDataStore< AthCommonMsg< Algorithm > >.

Reimplemented in AthAnalysisAlgorithm, AthFilterAlgorithm, AthHistogramAlgorithm, and PyAthena::Alg.

Definition at line 66 of file AthAlgorithm.cxx.

66 {
68
69 if (sc.isFailure()) {
70 return sc;
71 }
72 ServiceHandle<ICondSvc> cs("CondSvc",name());
73 for (auto h : outputHandles()) {
74 if (h->isCondition() && h->mode() == Gaudi::DataHandle::Writer) {
75 // do this inside the loop so we don't create the CondSvc until needed
76 if ( cs.retrieve().isFailure() ) {
77 ATH_MSG_WARNING("no CondSvc found: won't autoreg WriteCondHandles");
78 return StatusCode::SUCCESS;
79 }
80 if (cs->regHandle(this,*h).isFailure()) {
81 sc = StatusCode::FAILURE;
82 ATH_MSG_ERROR("unable to register WriteCondHandle " << h->fullKey()
83 << " with CondSvc");
84 }
85 }
86 }
87 return sc;
88}
#define ATH_MSG_ERROR(x)
static Double_t sc
virtual StatusCode sysInitialize() override
Override sysInitialize.
AthCommonDataStore(const std::string &name, T... args)
virtual std::vector< Gaudi::DataHandle * > outputHandles() const override
::StatusCode StatusCode
StatusCode definition for legacy code.

◆ sysStart()

virtual StatusCode AthCommonDataStore< AthCommonMsg< Algorithm > >::sysStart ( )
overridevirtualinherited

Handle START transition.

We override this in order to make sure that conditions handle keys can cache a pointer to the conditions container.

◆ updateVHKA()

void AthCommonDataStore< AthCommonMsg< Algorithm > >::updateVHKA ( Gaudi::Details::PropertyBase & )
inlineinherited

Definition at line 308 of file AthCommonDataStore.h.

308 {
309 // debug() << "updateVHKA for property " << p.name() << " " << p.toString()
310 // << " size: " << m_vhka.size() << endmsg;
311 for (auto &a : m_vhka) {
313 for (auto k : keys) {
314 k->setOwner(this);
315 }
316 }
317 }
std::vector< SG::VarHandleKeyArray * > m_vhka

Member Data Documentation

◆ m_accept

unsigned int TileCellSelector::m_accept
private

Definition at line 55 of file TileCellSelector.h.

◆ m_ADCmaskValueMinusEps

float TileCellSelector::m_ADCmaskValueMinusEps = 0.0F
private

Definition at line 180 of file TileCellSelector.h.

◆ m_ADCmaxMinusEps

float TileCellSelector::m_ADCmaxMinusEps = 0.0F
private

Definition at line 179 of file TileCellSelector.h.

◆ m_bitEneCell

bool TileCellSelector::m_bitEneCell[ptnlength] {}
private

Definition at line 139 of file TileCellSelector.h.

139{};

◆ m_bitEneChan

bool TileCellSelector::m_bitEneChan[3][ptnlength] {}
private

Definition at line 141 of file TileCellSelector.h.

141{};

◆ m_bitTimeCell

bool TileCellSelector::m_bitTimeCell[ptnlength] {}
private

Definition at line 140 of file TileCellSelector.h.

140{};

◆ m_bitTimeChan

bool TileCellSelector::m_bitTimeChan[3][ptnlength] {}
private

Definition at line 142 of file TileCellSelector.h.

142{};

◆ m_cabling

const TileCablingService* TileCellSelector::m_cabling
private

Definition at line 72 of file TileCellSelector.h.

◆ m_cellContainerKey

SG::ReadHandleKey<CaloCellContainer> TileCellSelector::m_cellContainerKey
private
Initial value:
{this,"CellContainerName",
"AllCalo", "Input Calo cell container key"}

Definition at line 111 of file TileCellSelector.h.

111 {this,"CellContainerName",
112 "AllCalo", "Input Calo cell container key"};

◆ m_chanBad

std::vector<bool> TileCellSelector::m_chanBad
private

Definition at line 97 of file TileCellSelector.h.

◆ m_chanDsp

std::vector<float> TileCellSelector::m_chanDsp
private

Definition at line 100 of file TileCellSelector.h.

◆ m_chanEne

std::vector<float> TileCellSelector::m_chanEne
private

Definition at line 98 of file TileCellSelector.h.

◆ m_chanQua

std::vector<float> TileCellSelector::m_chanQua
private

Definition at line 102 of file TileCellSelector.h.

◆ m_chanSel

std::vector<bool> TileCellSelector::m_chanSel
private

Definition at line 103 of file TileCellSelector.h.

◆ m_chanTDsp

std::vector<float> TileCellSelector::m_chanTDsp
private

Definition at line 101 of file TileCellSelector.h.

◆ m_chanTime

std::vector<float> TileCellSelector::m_chanTime
private

Definition at line 99 of file TileCellSelector.h.

◆ m_chanToCheck

std::vector<int> TileCellSelector::m_chanToCheck
private

Definition at line 172 of file TileCellSelector.h.

◆ m_chanToSkip

std::vector<bool> TileCellSelector::m_chanToSkip
private

Definition at line 104 of file TileCellSelector.h.

◆ m_checkDCS

bool TileCellSelector::m_checkDCS
private

Definition at line 155 of file TileCellSelector.h.

◆ m_checkDMUs

bool TileCellSelector::m_checkDMUs
private

Definition at line 157 of file TileCellSelector.h.

◆ m_checkError

bool TileCellSelector::m_checkError
private

Definition at line 167 of file TileCellSelector.h.

◆ m_checkJumps

bool TileCellSelector::m_checkJumps
private

Definition at line 156 of file TileCellSelector.h.

◆ m_checkOverHG

bool TileCellSelector::m_checkOverHG
private

Definition at line 159 of file TileCellSelector.h.

◆ m_checkOverLG

bool TileCellSelector::m_checkOverLG
private

Definition at line 158 of file TileCellSelector.h.

◆ m_checkUnderHG

bool TileCellSelector::m_checkUnderHG
private

Definition at line 161 of file TileCellSelector.h.

◆ m_checkUnderLG

bool TileCellSelector::m_checkUnderLG
private

Definition at line 160 of file TileCellSelector.h.

◆ m_checkWarning

bool TileCellSelector::m_checkWarning
private

Definition at line 166 of file TileCellSelector.h.

◆ m_const

unsigned int TileCellSelector::m_const
private

Definition at line 61 of file TileCellSelector.h.

◆ m_constLength

int TileCellSelector::m_constLength
private

Definition at line 148 of file TileCellSelector.h.

◆ m_counter

unsigned int TileCellSelector::m_counter
private

Definition at line 54 of file TileCellSelector.h.

◆ m_detStore

StoreGateSvc_t AthCommonDataStore< AthCommonMsg< Algorithm > >::m_detStore
privateinherited

Pointer to StoreGate (detector store by default)

Definition at line 393 of file AthCommonDataStore.h.

◆ m_digitsContainerKey

SG::ReadHandleKey<TileDigitsContainer> TileCellSelector::m_digitsContainerKey
private
Initial value:
{this,"DigitsContainerName",
"TileDigitsFlt", "Input Tile digits container key"}

Definition at line 114 of file TileCellSelector.h.

114 {this,"DigitsContainerName",
115 "TileDigitsFlt", "Input Tile digits container key"};

◆ m_dmuerr

unsigned int TileCellSelector::m_dmuerr
private

Definition at line 67 of file TileCellSelector.h.

◆ m_dqerr

unsigned int TileCellSelector::m_dqerr
private

Definition at line 66 of file TileCellSelector.h.

◆ m_dqStatusKey

SG::ReadHandleKey<TileDQstatus> TileCellSelector::m_dqStatusKey {this, "TileDQstatus", "TileDQstatus", "TileDQstatus key"}
private

Definition at line 74 of file TileCellSelector.h.

74{this, "TileDQstatus", "TileDQstatus", "TileDQstatus key"};

◆ m_drawer

std::vector<int> TileCellSelector::m_drawer
private

Definition at line 170 of file TileCellSelector.h.

◆ m_drawerToCheck

std::vector<int> TileCellSelector::m_drawerToCheck
private

Definition at line 171 of file TileCellSelector.h.

◆ m_drawerToSkip

std::vector<bool> TileCellSelector::m_drawerToSkip
private

Definition at line 105 of file TileCellSelector.h.

◆ m_eventInfoKey

SG::ReadHandleKey<xAOD::EventInfo> TileCellSelector::m_eventInfoKey
private
Initial value:
{this,"EventInfo",
"EventInfo", "Input event info key"}

Definition at line 121 of file TileCellSelector.h.

121 {this,"EventInfo",
122 "EventInfo", "Input event info key"};

◆ m_evtBCID

unsigned int TileCellSelector::m_evtBCID
private

Definition at line 80 of file TileCellSelector.h.

◆ m_evtNum

unsigned int TileCellSelector::m_evtNum
private

Definition at line 79 of file TileCellSelector.h.

◆ m_evtStore

StoreGateSvc_t AthCommonDataStore< AthCommonMsg< Algorithm > >::m_evtStore
privateinherited

Pointer to StoreGate (event store by default)

Definition at line 390 of file AthCommonDataStore.h.

◆ m_extendedExtraObjects

DataObjIDColl AthAlgorithm::m_extendedExtraObjects
privateinherited

Definition at line 79 of file AthAlgorithm.h.

◆ m_infoName

std::string TileCellSelector::m_infoName
private

Definition at line 177 of file TileCellSelector.h.

◆ m_jump

unsigned int TileCellSelector::m_jump
private

Definition at line 60 of file TileCellSelector.h.

◆ m_jumpDeltaHG

float TileCellSelector::m_jumpDeltaHG
private

Definition at line 144 of file TileCellSelector.h.

◆ m_jumpDeltaLG

float TileCellSelector::m_jumpDeltaLG
private

Definition at line 145 of file TileCellSelector.h.

◆ m_lumiBlock

unsigned int TileCellSelector::m_lumiBlock
private

Definition at line 78 of file TileCellSelector.h.

◆ m_maxBadDMU

int TileCellSelector::m_maxBadDMU
private

Definition at line 150 of file TileCellSelector.h.

◆ m_maxCell

unsigned int TileCellSelector::m_maxCell
private

Definition at line 57 of file TileCellSelector.h.

◆ m_maxChan

unsigned int TileCellSelector::m_maxChan
private

Definition at line 59 of file TileCellSelector.h.

◆ m_maxEneCell

float TileCellSelector::m_maxEneCell
private

Definition at line 125 of file TileCellSelector.h.

◆ m_maxEneChan

float TileCellSelector::m_maxEneChan[3] {}
private

Definition at line 127 of file TileCellSelector.h.

127{};

◆ m_maxTimeCell

float TileCellSelector::m_maxTimeCell
private

Definition at line 129 of file TileCellSelector.h.

◆ m_maxTimeChan

float TileCellSelector::m_maxTimeChan[3] {}
private

Definition at line 131 of file TileCellSelector.h.

131{};

◆ m_maxVerboseCnt

int TileCellSelector::m_maxVerboseCnt
private

Definition at line 174 of file TileCellSelector.h.

◆ m_minBadDMU

int TileCellSelector::m_minBadDMU
private

Definition at line 149 of file TileCellSelector.h.

◆ m_minBadMB

int TileCellSelector::m_minBadMB
private

Definition at line 151 of file TileCellSelector.h.

◆ m_minCell

unsigned int TileCellSelector::m_minCell
private

Definition at line 56 of file TileCellSelector.h.

◆ m_minChan

unsigned int TileCellSelector::m_minChan
private

Definition at line 58 of file TileCellSelector.h.

◆ m_minEneCell

float TileCellSelector::m_minEneCell
private

Definition at line 124 of file TileCellSelector.h.

◆ m_minEneChan

float TileCellSelector::m_minEneChan[3] {}
private

Definition at line 126 of file TileCellSelector.h.

126{};

◆ m_minTimeCell

float TileCellSelector::m_minTimeCell
private

Definition at line 128 of file TileCellSelector.h.

◆ m_minTimeChan

float TileCellSelector::m_minTimeChan[3] {}
private

Definition at line 130 of file TileCellSelector.h.

130{};

◆ m_nDrawerOff

std::vector<int> TileCellSelector::m_nDrawerOff
private

Definition at line 175 of file TileCellSelector.h.

◆ m_overflowHG

float TileCellSelector::m_overflowHG
private

Definition at line 163 of file TileCellSelector.h.

◆ m_overflowLG

float TileCellSelector::m_overflowLG
private

Definition at line 162 of file TileCellSelector.h.

◆ m_overHG

unsigned int TileCellSelector::m_overHG
private

Definition at line 63 of file TileCellSelector.h.

◆ m_overLG

unsigned int TileCellSelector::m_overLG
private

Definition at line 62 of file TileCellSelector.h.

◆ m_pedDeltaHG

float TileCellSelector::m_pedDeltaHG
private

Definition at line 146 of file TileCellSelector.h.

◆ m_pedDeltaLG

float TileCellSelector::m_pedDeltaLG
private

Definition at line 147 of file TileCellSelector.h.

◆ m_printOnly

bool TileCellSelector::m_printOnly
private

Definition at line 168 of file TileCellSelector.h.

◆ m_ptnEneCell

int TileCellSelector::m_ptnEneCell
private

Definition at line 132 of file TileCellSelector.h.

◆ m_ptnEneChan

int TileCellSelector::m_ptnEneChan[3] {}
private

Definition at line 133 of file TileCellSelector.h.

133{};

◆ m_ptnTimeCell

int TileCellSelector::m_ptnTimeCell
private

Definition at line 134 of file TileCellSelector.h.

◆ m_ptnTimeChan

int TileCellSelector::m_ptnTimeChan[3] {}
private

Definition at line 135 of file TileCellSelector.h.

135{};

◆ m_rawChannelContainerKey

SG::ReadHandleKey<TileRawChannelContainer> TileCellSelector::m_rawChannelContainerKey
private
Initial value:
{this,"RawChannelContainerName",
"TileRawChannelCnt",
"Input Tile raw channel container key"}

Definition at line 117 of file TileCellSelector.h.

117 {this,"RawChannelContainerName",
118 "TileRawChannelCnt",
119 "Input Tile raw channel container key"};

◆ m_readCells

bool TileCellSelector::m_readCells
private

Definition at line 107 of file TileCellSelector.h.

◆ m_readDigits

bool TileCellSelector::m_readDigits
private

Definition at line 109 of file TileCellSelector.h.

◆ m_readRawChannels

bool TileCellSelector::m_readRawChannels
private

Definition at line 108 of file TileCellSelector.h.

◆ m_runNum

unsigned int TileCellSelector::m_runNum
private

Definition at line 77 of file TileCellSelector.h.

◆ m_secondMaxLevel

float TileCellSelector::m_secondMaxLevel
private

Definition at line 143 of file TileCellSelector.h.

◆ m_selectGain

int TileCellSelector::m_selectGain
private

Definition at line 136 of file TileCellSelector.h.

◆ m_skipEmpty

bool TileCellSelector::m_skipEmpty
private

Definition at line 152 of file TileCellSelector.h.

◆ m_skipGain

bool TileCellSelector::m_skipGain[2] {}
private

Definition at line 137 of file TileCellSelector.h.

137{};

◆ m_skipMasked

bool TileCellSelector::m_skipMasked
private

Definition at line 153 of file TileCellSelector.h.

◆ m_skipMBTS

bool TileCellSelector::m_skipMBTS
private

Definition at line 154 of file TileCellSelector.h.

◆ m_tileBadChanTool

ToolHandle<ITileBadChanTool> TileCellSelector::m_tileBadChanTool {this, "TileBadChanTool", "TileBadChanTool", "Tile bad channel tool"}
private

Definition at line 73 of file TileCellSelector.h.

73{this, "TileBadChanTool", "TileBadChanTool", "Tile bad channel tool"};

◆ m_tileDCS

ToolHandle<ITileDCSTool> TileCellSelector::m_tileDCS {this, "TileDCSTool", "TileDCSTool", "Tile DCS tool"}
private

Definition at line 75 of file TileCellSelector.h.

75{this, "TileDCSTool", "TileDCSTool", "Tile DCS tool"};

◆ m_tileError

unsigned int TileCellSelector::m_tileError
private

Definition at line 95 of file TileCellSelector.h.

◆ m_tileFlag

unsigned int TileCellSelector::m_tileFlag
private

Definition at line 93 of file TileCellSelector.h.

◆ m_tileHWID

const TileHWID* TileCellSelector::m_tileHWID
private

Definition at line 71 of file TileCellSelector.h.

◆ m_tileID

const TileID* TileCellSelector::m_tileID
private

Definition at line 70 of file TileCellSelector.h.

◆ m_tileInfo

const TileInfo* TileCellSelector::m_tileInfo
private

Definition at line 178 of file TileCellSelector.h.

◆ m_underflowHG

float TileCellSelector::m_underflowHG
private

Definition at line 165 of file TileCellSelector.h.

◆ m_underflowLG

float TileCellSelector::m_underflowLG
private

Definition at line 164 of file TileCellSelector.h.

◆ m_underHG

unsigned int TileCellSelector::m_underHG
private

Definition at line 65 of file TileCellSelector.h.

◆ m_underLG

unsigned int TileCellSelector::m_underLG
private

Definition at line 64 of file TileCellSelector.h.

◆ m_varHandleArraysDeclared

bool AthCommonDataStore< AthCommonMsg< Algorithm > >::m_varHandleArraysDeclared
privateinherited

Definition at line 399 of file AthCommonDataStore.h.

◆ m_vhka

std::vector<SG::VarHandleKeyArray*> AthCommonDataStore< AthCommonMsg< Algorithm > >::m_vhka
privateinherited

Definition at line 398 of file AthCommonDataStore.h.

◆ m_warnerr

unsigned int TileCellSelector::m_warnerr
private

Definition at line 68 of file TileCellSelector.h.


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