Belle II Software development
KLMChannelStatus.cc
1/**************************************************************************
2 * basf2 (Belle II Analysis Software Framework) *
3 * Author: The Belle II Collaboration *
4 * *
5 * See git log for contributors and copyright holders. *
6 * This file is licensed under LGPL-3.0, see LICENSE.md. *
7 **************************************************************************/
8
9/* Own header. */
10#include <klm/dbobjects/KLMChannelStatus.h>
11
12/* KLM headers. */
13#include <klm/dataobjects/KLMChannelIndex.h>
14#include <klm/dataobjects/KLMElementNumbers.h>
15
16/* Basf2 headers. */
17#include <framework/logging/Logger.h>
18
19using namespace Belle2;
20
22{
23}
24
26{
27}
28
31 std::map<KLMChannelNumber, enum ChannelStatus>::const_iterator it;
32 it = m_ChannelStatus.find(channel);
33 if (it == m_ChannelStatus.end())
34 return c_Unknown;
35 return it->second;
36}
37
39 enum ChannelStatus status)
40{
41 std::map<KLMChannelNumber, enum ChannelStatus>::iterator it;
42 it = m_ChannelStatus.find(channel);
43 if (it == m_ChannelStatus.end()) {
44 m_ChannelStatus.insert(
45 std::pair<KLMChannelNumber, enum ChannelStatus>(channel, status));
46 } else {
47 it->second = status;
48 }
49}
50
52{
53 KLMChannelIndex klmChannels;
54 for (KLMChannelIndex& klmChannel : klmChannels)
55 setChannelStatus(klmChannel.getKLMChannelNumber(), status);
56}
57
59{
60 int active;
61 int subdetector, section, sector, layer;
62 const KLMElementNumbers* elementNumbers =
64 elementNumbers->moduleNumberToElementNumbers(
65 module, &subdetector, &section, &sector, &layer);
66 KLMChannelIndex klmModule(subdetector, section, sector, layer, 1, 1);
67 /* Plane number is 0-based for BKLM. */
68 if (subdetector == KLMElementNumbers::c_BKLM) {
69 KLMChannelIndex klmModuleTemp(subdetector, section, sector, layer, 0, 1);
70 klmModule = klmModuleTemp;
71 }
73 KLMChannelIndex klmNextModule(klmModule);
74 ++klmNextModule;
75 KLMChannelIndex klmChannel(klmModule);
77 active = 0;
78 for (; klmChannel != klmNextModule; ++klmChannel) {
79 KLMChannelNumber channel = klmChannel.getKLMChannelNumber();
80 ChannelStatus status = getChannelStatus(channel);
81 if (status == c_Unknown)
82 B2FATAL("Incomplete KLM channel data.");
83 if (status != c_Dead)
84 active++;
85 }
86 return active;
87}
88
90{
91 if (this->m_ChannelStatus.size() != status.m_ChannelStatus.size())
92 return false;
93 std::map<KLMChannelNumber, enum ChannelStatus>::iterator it, it2;
94 it = this->m_ChannelStatus.begin();
95 it2 = status.m_ChannelStatus.begin();
96 while (it != this->m_ChannelStatus.end()) {
97 if (it->first != it2->first)
98 return false;
99 if (it->second != it2->second)
100 return false;
101 ++it;
102 ++it2;
103 }
104 return true;
105}
106
108{
109 unsigned int channels = 0;
110 if (this->m_ChannelStatus.size() != status.m_ChannelStatus.size())
111 return 0;
112 std::map<KLMChannelNumber, enum ChannelStatus>::iterator it, it2;
113 it = this->m_ChannelStatus.begin();
114 it2 = status.m_ChannelStatus.begin();
115 while (it != this->m_ChannelStatus.end()) {
116 if (it->first != it2->first)
117 return 0;
118 if ((it->second == c_Normal) && (it2->second != c_Normal))
119 ++channels;
120 ++it;
121 ++it2;
122 }
123 return channels;
124}
KLM channel index.
KLMChannelNumber getKLMChannelNumber() const
Get KLM channel number.
void setIndexLevel(enum IndexLevel indexLevel)
Set index level.
KLM channel status.
void setChannelStatus(KLMChannelNumber channel, enum ChannelStatus status)
Set channel status.
enum ChannelStatus getChannelStatus(KLMChannelNumber channel) const
Get channel status.
void setStatusAllChannels(enum ChannelStatus status)
Set status for all channels.
int getActiveStripsInModule(KLMChannelNumber module) const
Get number of active strips in the specified KLM module.
bool operator==(KLMChannelStatus &status)
Operator ==.
std::map< KLMChannelNumber, enum ChannelStatus > m_ChannelStatus
Channel data.
ChannelStatus
Channel status.
@ c_Normal
Normally operating channel.
@ c_Dead
Dead channel (no signal).
@ c_Unknown
Unknown status (no data).
unsigned int newNormalChannels(KLMChannelStatus &status)
Number of new channels with status c_Normal that have a different status in another channel-status da...
KLM element numbers.
static const KLMElementNumbers & Instance()
Instantiation.
void moduleNumberToElementNumbers(KLMModuleNumber module, int *subdetector, int *section, int *sector, int *layer) const
Get element numbers by module number.
uint16_t KLMChannelNumber
Channel number.
Abstract base class for different kinds of events.