Belle II Software development
KLMElectronicsMapImporter Class Reference

KLM database importer. More...

#include <KLMElectronicsMapImporter.h>

Public Member Functions

 KLMElectronicsMapImporter ()
 Constructor.
 
 ~KLMElectronicsMapImporter ()
 Destructor.
 
void setIOV (int experimentLow, int runLow, int experimentHigh, int runHigh)
 Set interval of validity.
 
void clearElectronicsMap ()
 Clear electronics map (to be able to import its multiple versions).
 
void loadBKLMElectronicsMap (int version)
 Load BKLM electronics map.
 
void loadEKLMElectronicsMap (int version, bool mc)
 Load EKLM electronics map.
 
void setLane (int subdetector, int section, int sector, int layer, int lane)
 Set non-default lane for all channels in a module.
 
void setLane (int subdetector, int section, int sector, int layer, int plane, int lane)
 Set non-default lane for all channels in a plane.
 
void importElectronicsMap ()
 Import electronics map.
 

Private Member Functions

int getEKLMStripFirmwareBySoftware (int stripSoftware) const
 Get EKLM firmware strip number by software strip number.
 
void setChannelsEKLMSegment (int section, int sector, int layer, int plane, int segment, int firmwareSegment)
 Set channels for EKLM segment.
 
void addEKLMLane (int section, int sector, int layer, int copper, int slot, int lane)
 Add EKLM electronics map lane.
 

Private Attributes

const KLMElementNumbersm_ElementNumbers
 Element numbers.
 
std::map< KLMChannelNumber, KLMElectronicsChannelm_ChannelMap
 Data for creation of the electronics map.
 
int m_ExperimentLow = 0
 Low experiment.
 
int m_RunLow = 0
 Low run.
 
int m_ExperimentHigh = -1
 High experiment.
 
int m_RunHigh = -1
 High run.
 

Detailed Description

KLM database importer.

Definition at line 27 of file KLMElectronicsMapImporter.h.

Constructor & Destructor Documentation

◆ KLMElectronicsMapImporter()

Constructor.

Definition at line 24 of file KLMElectronicsMapImporter.cc.

24 :
26{
27}
const KLMElementNumbers * m_ElementNumbers
Element numbers.
static const KLMElementNumbers & Instance()
Instantiation.

◆ ~KLMElectronicsMapImporter()

Destructor.

Definition at line 29 of file KLMElectronicsMapImporter.cc.

30{
31}

Member Function Documentation

◆ addEKLMLane()

void addEKLMLane ( int  section,
int  sector,
int  layer,
int  copper,
int  slot,
int  lane 
)
private

Add EKLM electronics map lane.

Parameters
[in]sectionSection.
[in]sectorSector.
[in]layerLayer.
[in]copperCopper.
[in]slotSlot.
[in]laneLane.

Definition at line 187 of file KLMElectronicsMapImporter.cc.

189{
190 for (int plane = 1; plane <= EKLMElementNumbers::getMaximalPlaneNumber(); ++plane) {
191 for (int strip = 1; strip <= EKLMElementNumbers::getMaximalStripNumber(); ++strip) {
192 int axis = plane - 1;
193 int channel = getEKLMStripFirmwareBySoftware(strip);
195 section, sector, layer, plane, strip);
196 m_ChannelMap.insert(
197 std::pair<KLMChannelNumber, KLMElectronicsChannel>(
198 detectorChannel,
199 KLMElectronicsChannel(EKLM_ID + copper, slot, lane, axis, channel)));
200 }
201 }
202}
static constexpr int getMaximalStripNumber()
Get maximal strip number.
static constexpr int getMaximalPlaneNumber()
Get maximal plane number.
BKLM electronics channel.
std::map< KLMChannelNumber, KLMElectronicsChannel > m_ChannelMap
Data for creation of the electronics map.
int getEKLMStripFirmwareBySoftware(int stripSoftware) const
Get EKLM firmware strip number by software strip number.
KLMChannelNumber channelNumberEKLM(int section, int sector, int layer, int plane, int strip) const
Get channel number for EKLM.
uint16_t KLMChannelNumber
Channel number.

◆ clearElectronicsMap()

void clearElectronicsMap ( )

Clear electronics map (to be able to import its multiple versions).

Definition at line 42 of file KLMElectronicsMapImporter.cc.

43{
44 m_ChannelMap.clear();
45}

◆ getEKLMStripFirmwareBySoftware()

int getEKLMStripFirmwareBySoftware ( int  stripSoftware) const
private

Get EKLM firmware strip number by software strip number.

Parameters
[in]stripSoftwareSoftware strip number.

Definition at line 175 of file KLMElectronicsMapImporter.cc.

176{
177 int segment, strip;
178 /* Segment is 0-based here. */
179 segment = (stripSoftware - 1) / EKLMElementNumbers::getNStripsSegment();
180 /* Order of segment readout boards in the firmware is opposite. */
181 segment = 4 - segment;
182 strip = segment * EKLMElementNumbers::getNStripsSegment() +
183 (stripSoftware - 1) % EKLMElementNumbers::getNStripsSegment() + 1;
184 return strip;
185}
static constexpr int getNStripsSegment()
Get number of strips in a segment.

◆ importElectronicsMap()

void importElectronicsMap ( )

Import electronics map.

Definition at line 428 of file KLMElectronicsMapImporter.cc.

429{
431 electronicsMap.construct();
432 std::map<KLMChannelNumber, KLMElectronicsChannel>::iterator it;
433 for (it = m_ChannelMap.begin(); it != m_ChannelMap.end(); ++it) {
434 electronicsMap->addChannel(
435 it->first,
436 it->second.getCopper(),
437 it->second.getSlot(),
438 it->second.getLane(),
439 it->second.getAxis(),
440 it->second.getChannel());
441 }
444 electronicsMap.import(iov);
445}
bool import(const IntervalOfValidity &iov)
Import the object to database.
Definition: DBImportBase.cc:36
Class for importing a single object to the database.
void construct(Args &&... params)
Construct an object of type T in this DBImportObjPtr using the provided constructor arguments.
A class that describes the interval of experiments/runs for which an object in the database is valid.

◆ loadBKLMElectronicsMap()

void loadBKLMElectronicsMap ( int  version)

Load BKLM electronics map.

Parameters
[in]versionmap version: 1) Before experiment 10. 2) Experiment 10 and later (mapping in chimney sector has changed)

Definition at line 47 of file KLMElectronicsMapImporter.cc.

48{
49 const int minimalVersion = 1;
50 const int maximalVersion = 2;
51 if (version < minimalVersion || version > maximalVersion) {
52 B2FATAL("Incorrect version (" << version << ") of BKLM electronics map. "
53 "It must be from " << minimalVersion << " to " << maximalVersion);
54 }
55 int copperId = 0;
56 int slotId = 0;
57 int laneId;
58 int axisId = 0;
60 for (KLMChannelIndex bklmPlane = bklmPlanes.beginBKLM();
61 bklmPlane != bklmPlanes.endBKLM(); ++bklmPlane) {
62 int section = bklmPlane.getSection();
63 int sector = bklmPlane.getSector();
64 int layer = bklmPlane.getLayer();
65 int plane = bklmPlane.getPlane();
66
68 if (sector == 3 || sector == 4 || sector == 5 || sector == 6)
69 copperId = 1 + BKLM_ID;
70 if (sector == 1 || sector == 2 || sector == 7 || sector == 8)
71 copperId = 2 + BKLM_ID;
72 }
74 if (sector == 3 || sector == 4 || sector == 5 || sector == 6)
75 copperId = 3 + BKLM_ID;
76 if (sector == 1 || sector == 2 || sector == 7 || sector == 8)
77 copperId = 4 + BKLM_ID;
78 }
79
80 if (sector == 3 || sector == 4 || sector == 5 || sector == 6)
81 slotId = sector - 2;
82 if (sector == 1 || sector == 2)
83 slotId = sector + 2;
84 if (sector == 7 || sector == 8)
85 slotId = sector - 6;
86
88 laneId = layer + 5;
89 axisId = plane;
90 } else {
91 laneId = layer;
93 axisId = 1;
95 axisId = 0;
96 }
97
98 int MaxiChannel = BKLMElementNumbers::getNStrips(
99 section, sector, layer, plane);
100
101 bool dontFlip = false;
103 (sector == 7 || sector == 8 || sector == 1 || sector == 2))
104 dontFlip = true;
106 (sector == 4 || sector == 5 || sector == 6 || sector == 7))
107 dontFlip = true;
108
109 for (int iStrip = 1; iStrip <= MaxiChannel; iStrip++) {
110 int channelId = iStrip;
111
112 if (!(dontFlip && layer >= BKLMElementNumbers::c_FirstRPCLayer && plane == BKLMElementNumbers::c_PhiPlane))
113 channelId = MaxiChannel - iStrip + 1;
114
115 if (plane == BKLMElementNumbers::c_PhiPlane) {
116 if (layer == 1)
117 channelId += 4;
118 if (layer == 2)
119 channelId += 2;
120 }
121
122 if (plane == BKLMElementNumbers::c_ZPlane) {
124 int channelCheck = channelId;
127 if (layer == 1) {
128 if (version == 1) {
129 if (channelCheck > 0 && channelCheck < 9)
130 channelId = 9 - channelId;
131 if (channelCheck > 8 && channelCheck < 24)
132 channelId = 54 - channelId;
133 if (channelCheck > 23 && channelCheck < 39)
134 channelId = 54 - channelId;
135 } else {
136 if (channelCheck > 0 && channelCheck < 9)
137 channelId = 9 - channelId; // 8 : 1
138 if (channelCheck > 8 && channelCheck < 24)
139 channelId = 39 - channelId; // 30 : 16
140 if (channelCheck > 23 && channelCheck < 39)
141 channelId = 69 - channelId; // 45 : 31
142 }
143 }
144 if (layer == 2) {
145 if (channelCheck > 0 && channelCheck < 10)
146 channelId = 10 - channelId;
147 if (channelCheck > 9 && channelCheck < 24)
148 channelId = 40 - channelId;
149 if (channelCheck > 23 && channelCheck < 39)
150 channelId = 69 - channelId;
151 }
152 } else { // All the sectors except the chimney one
153 if (channelCheck > 0 && channelCheck < 10)
154 channelId = 10 - channelId;
155 if (channelCheck > 9 && channelCheck < 25)
156 channelId = 40 - channelId;
157 if (channelCheck > 24 && channelCheck < 40)
158 channelId = 70 - channelId;
159 if (channelCheck > 39 && channelCheck < 55)
160 channelId = 100 - channelId;
161 }
162 }
163 }
164
166 section, sector, layer, plane, iStrip);
167 m_ChannelMap.insert(
168 std::pair<KLMChannelNumber, KLMElectronicsChannel>(
169 detectorChannel,
170 KLMElectronicsChannel(copperId, slotId, laneId, axisId, channelId)));
171 }
172 }
173}
@ c_ChimneySector
Chimney sector: BB3 in 1-based notation; BB2 in 0-based notation.
@ c_FirstRPCLayer
First RPC layer.
static int getNStrips(int section, int sector, int layer, int plane)
Get number of strips.
KLM channel index.
KLMChannelIndex beginBKLM()
First channel for BKLM.
KLMChannelIndex & endBKLM()
Last channel for BKLM.
KLMChannelNumber channelNumberBKLM(int section, int sector, int layer, int plane, int strip) const
Get channel number for BKLM.

◆ loadEKLMElectronicsMap()

void loadEKLMElectronicsMap ( int  version,
bool  mc 
)

Load EKLM electronics map.

Parameters
[in]versionMap version: 1) Phase 2, experiment 3. Connection of cables was wrong for backward sectors 2 and 3. 2) Phase 3, starting from experiment 4.
[in]mcMC or data (MC does not have occasional cable switches that exist in the data).

Definition at line 222 of file KLMElectronicsMapImporter.cc.

223{
224 const int minimalVersion = 1;
225 const int maximalVersion = 2;
226 if (version < minimalVersion || version > maximalVersion) {
227 B2FATAL("Incorrect version (" << version << ") of EKLM electronics map. "
228 "It must be from " << minimalVersion << " to " << maximalVersion);
229 }
230 /* Backward section. */
231 /* Sector 1. */
241 /* Wrong connection was fixed between phase 2 and phase 3. */
242 if (mc || (version >= 2)) {
245 } else {
248 }
250 /* Sector 2. */
251 if (version == 1) {
264 } else {
277 }
278 /* Sector 3. */
279 if (version == 1) {
292 } else {
305 }
306 /* Sector 4. */
319 /* Forward section. */
320 /* Sector 1. */
335 /* Sector 2. */
350 /* Sector 3. */
365 /* Sector 4. */
380 /* Switch of internal cables. */
381 if (!mc) {
382 setChannelsEKLMSegment(1, 1, 5, 1, 1, 4);
383 setChannelsEKLMSegment(1, 1, 5, 1, 2, 5);
384 setChannelsEKLMSegment(1, 1, 5, 2, 1, 4);
385 setChannelsEKLMSegment(1, 1, 5, 2, 2, 5);
386 setChannelsEKLMSegment(2, 1, 10, 2, 3, 2);
387 setChannelsEKLMSegment(2, 1, 10, 2, 4, 3);
388 }
389}
void setChannelsEKLMSegment(int section, int sector, int layer, int plane, int segment, int firmwareSegment)
Set channels for EKLM segment.
void addEKLMLane(int section, int sector, int layer, int copper, int slot, int lane)
Add EKLM electronics map lane.

◆ setChannelsEKLMSegment()

void setChannelsEKLMSegment ( int  section,
int  sector,
int  layer,
int  plane,
int  segment,
int  firmwareSegment 
)
private

Set channels for EKLM segment.

Parameters
[in]sectionSection.
[in]sectorSector.
[in]layerLayer.
[in]planePlane.
[in]segmentSegment.
[in]firmwareSegmentSegment number in firmware.

Definition at line 204 of file KLMElectronicsMapImporter.cc.

207{
208 std::map<KLMChannelNumber, KLMElectronicsChannel>:: iterator it;
209 for (int strip = 1; strip <= EKLMElementNumbers::getNStripsSegment(); ++strip) {
210 int stripPlane = strip + EKLMElementNumbers::getNStripsSegment() * (segment - 1);
212 section, sector, layer, plane, stripPlane);
213 it = m_ChannelMap.find(detectorChannel);
214 if (it == m_ChannelMap.end())
215 B2FATAL("The KLM electronics map is not loaded or incomplete.");
216 int channel = (firmwareSegment - 1) * EKLMElementNumbers::getNStripsSegment() +
217 strip;
218 it->second.setChannel(channel);
219 }
220}

◆ setIOV()

void setIOV ( int  experimentLow,
int  runLow,
int  experimentHigh,
int  runHigh 
)

Set interval of validity.

Definition at line 33 of file KLMElectronicsMapImporter.cc.

35{
36 m_ExperimentLow = experimentLow;
37 m_RunLow = runLow;
38 m_ExperimentHigh = experimentHigh;
39 m_RunHigh = runHigh;
40}

◆ setLane() [1/2]

void setLane ( int  subdetector,
int  section,
int  sector,
int  layer,
int  lane 
)

Set non-default lane for all channels in a module.

Parameters
[in]subdetectorSubdetector.
[in]sectionSection.
[in]sectorSector.
[in]layerLayer.
[in]laneLane.

Definition at line 391 of file KLMElectronicsMapImporter.cc.

393{
394 std::map<KLMChannelNumber, KLMElectronicsChannel>::iterator it;
395 int minimalPlane = m_ElementNumbers->getMinimalPlaneNumber(subdetector);
396 KLMChannelIndex klmChannel(subdetector, section, sector, layer, minimalPlane, 1);
397 KLMChannelIndex klmModule(klmChannel);
398 klmModule.setIndexLevel(KLMChannelIndex::c_IndexLevelLayer);
399 KLMChannelIndex klmNextModule(klmModule);
400 ++klmNextModule;
401 for (; klmChannel != klmNextModule; ++klmChannel) {
402 KLMChannelNumber channel = klmChannel.getKLMChannelNumber();
403 it = m_ChannelMap.find(channel);
404 if (it == m_ChannelMap.end())
405 B2FATAL("The KLM electronics map is not loaded or incomplete.");
406 it->second.setLane(lane);
407 }
408}
int getMinimalPlaneNumber(int subdetector) const
Get minimal plane number.

◆ setLane() [2/2]

void setLane ( int  subdetector,
int  section,
int  sector,
int  layer,
int  plane,
int  lane 
)

Set non-default lane for all channels in a plane.

Parameters
[in]subdetectorSubdetector.
[in]sectionSection.
[in]sectorSector.
[in]layerLayer.
[in]planePlane.
[in]laneLane.

Definition at line 410 of file KLMElectronicsMapImporter.cc.

412{
413 std::map<KLMChannelNumber, KLMElectronicsChannel>::iterator it;
414 KLMChannelIndex klmChannel(subdetector, section, sector, layer, plane, 1);
415 KLMChannelIndex klmPlane(klmChannel);
416 klmPlane.setIndexLevel(KLMChannelIndex::c_IndexLevelPlane);
417 KLMChannelIndex klmNextPlane(klmPlane);
418 ++klmNextPlane;
419 for (; klmChannel != klmNextPlane; ++klmChannel) {
420 KLMChannelNumber channel = klmChannel.getKLMChannelNumber();
421 it = m_ChannelMap.find(channel);
422 if (it == m_ChannelMap.end())
423 B2FATAL("The KLM electronics map is not loaded or incomplete.");
424 it->second.setLane(lane);
425 }
426}

Member Data Documentation

◆ m_ChannelMap

std::map<KLMChannelNumber, KLMElectronicsChannel> m_ChannelMap
private

Data for creation of the electronics map.

Definition at line 141 of file KLMElectronicsMapImporter.h.

◆ m_ElementNumbers

const KLMElementNumbers* m_ElementNumbers
private

Element numbers.

Definition at line 138 of file KLMElectronicsMapImporter.h.

◆ m_ExperimentHigh

int m_ExperimentHigh = -1
private

High experiment.

Definition at line 150 of file KLMElectronicsMapImporter.h.

◆ m_ExperimentLow

int m_ExperimentLow = 0
private

Low experiment.

Definition at line 144 of file KLMElectronicsMapImporter.h.

◆ m_RunHigh

int m_RunHigh = -1
private

High run.

Definition at line 153 of file KLMElectronicsMapImporter.h.

◆ m_RunLow

int m_RunLow = 0
private

Low run.

Definition at line 147 of file KLMElectronicsMapImporter.h.


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