Belle II Software development
KLMCalibrationChecker Class Reference

KLM calibration checker. More...

#include <KLMCalibrationChecker.h>

Public Member Functions

 KLMCalibrationChecker ()
 Constructor.
 
 ~KLMCalibrationChecker ()
 Destructor.
 
void setExperimentRun (int experiment, int run)
 Set experiment and run numbers.
 
void setTestingPayload (const std::string &testingPayloadName)
 Set testing payload name.
 
void setGlobalTag (const std::string &globalTagName)
 Set Global Tag name.
 
void setAlignmentResultsFile (const std::string &alignmentResultsFile)
 Set alignment results file.
 
void setStripEfficiencyResultsFile (const std::string &stripEfficiencyResultsFile)
 Set strip efficiency results file.
 
void setTimeCableDelayResultsFile (const std::string &timeCableDelayResultsFile)
 Set time cable delay results file.
 
void setTimeConstantsResultsFile (const std::string &timeConstantsResultsFile)
 Set time constants result file.
 
void checkAlignment ()
 Check alignment.
 
void checkStripEfficiency ()
 Check strip efficiency.
 
void createStripEfficiencyHistograms ()
 Create strip efficiency histograms.
 
void checkTimeCableDelay ()
 Check time cable delay.
 
void checkTimeConstants ()
 Check time constants.
 

Private Member Functions

void initializeDatabase ()
 Initialize the database.
 
void resetDatabase ()
 Reset the database.
 
template<class T >
void printPayloadInformation (DBObjPtr< T > &dbObject)
 Print payload information.
 

Private Attributes

int m_experiment
 Experiment number.
 
int m_run
 Run number.
 
std::string m_testingPayloadName = ""
 Testing payload location.
 
std::string m_GlobalTagName = ""
 Global Tag name.
 
std::string m_AlignmentResultsFile = "alignment.root"
 Output file for alignment results.
 
std::string m_StripEfficiencyResultsFile = "strip_efficiency.root"
 Output file for alignment results.
 
std::string m_TimeCableDelayResultsFile = "timeCableDelay.root"
 Output file for time cable delay results.
 
std::string m_TimeConstantsResultsFile = "timeConstants.root"
 Output file for time constants results.
 
const KLMElementNumbersm_ElementNumbers
 Element numbers.
 
StoreObjPtr< EventMetaDatam_EventMetaData
 Event metadata.
 

Detailed Description

KLM calibration checker.

Definition at line 32 of file KLMCalibrationChecker.h.

Constructor & Destructor Documentation

◆ KLMCalibrationChecker()

Constructor.

Definition at line 38 of file KLMCalibrationChecker.cc.

38 :
39 m_experiment(0),
40 m_run(0),
42{
43}
const KLMElementNumbers * m_ElementNumbers
Element numbers.
static const KLMElementNumbers & Instance()
Instantiation.

◆ ~KLMCalibrationChecker()

Destructor.

Definition at line 45 of file KLMCalibrationChecker.cc.

46{
47}

Member Function Documentation

◆ checkAlignment()

void checkAlignment ( )

Check alignment.

Definition at line 88 of file KLMCalibrationChecker.cc.

89{
90 /* Initialize the database. */
92 /* Now we can read the payload. */
93 DBObjPtr<BKLMAlignment> bklmAlignment;
94 DBObjPtr<EKLMAlignment> eklmAlignment;
95 DBObjPtr<EKLMSegmentAlignment> eklmSegmentAlignment;
96 DBObjPtr<BKLMAlignment> bklmAlignmentErrors("BKLMAlignment_ERRORS");
97 DBObjPtr<EKLMAlignment> eklmAlignmentErrors("EKLMAlignment_ERRORS");
98 DBObjPtr<EKLMSegmentAlignment> eklmSegmentAlignmentErrors("EKLMSegmentAlignment_ERRORS");
99 DBObjPtr<BKLMAlignment> bklmAlignmentCorrections("BKLMAlignment_CORRECTIONS");
100 DBObjPtr<EKLMAlignment> eklmAlignmentCorrections("EKLMAlignment_CORRECTIONS");
101 DBObjPtr<EKLMSegmentAlignment> eklmSegmentAlignmentCorrections("EKLMSegmentAlignment_CORRECTIONS");
102 if (!bklmAlignment.isValid() ||
103 !eklmAlignment.isValid() ||
104 !eklmSegmentAlignment.isValid() ||
105 !bklmAlignmentErrors.isValid() ||
106 !eklmAlignmentErrors.isValid() ||
107 !eklmSegmentAlignmentErrors.isValid() ||
108 !bklmAlignmentCorrections.isValid() ||
109 !eklmAlignmentCorrections.isValid() ||
110 !eklmSegmentAlignmentCorrections.isValid())
111 B2FATAL("Alignment data are not valid.");
112 if (m_GlobalTagName != "") {
113 printPayloadInformation(bklmAlignment);
114 printPayloadInformation(eklmAlignment);
115 printPayloadInformation(eklmSegmentAlignment);
116 printPayloadInformation(bklmAlignmentErrors);
117 printPayloadInformation(eklmAlignmentErrors);
118 printPayloadInformation(eklmSegmentAlignmentErrors);
119 printPayloadInformation(bklmAlignmentCorrections);
120 printPayloadInformation(eklmAlignmentCorrections);
121 printPayloadInformation(eklmSegmentAlignmentCorrections);
122 }
123 /* Create trees with alignment results. */
124 int section, sector, layer, plane, segment, param;
125 float value, correction, error;
126 TFile* alignmentResults = new TFile(m_AlignmentResultsFile.c_str(),
127 "recreate");
128 TTree* bklmModuleTree = new TTree("bklm_module",
129 "BKLM module alignment data.");
130 bklmModuleTree->Branch("experiment", &m_experiment, "experiment/I");
131 bklmModuleTree->Branch("run", &m_run, "run/I");
132 bklmModuleTree->Branch("section", &section, "section/I");
133 bklmModuleTree->Branch("sector", &sector, "sector/I");
134 bklmModuleTree->Branch("layer", &layer, "layer/I");
135 bklmModuleTree->Branch("param", &param, "param/I");
136 bklmModuleTree->Branch("value", &value, "value/F");
137 bklmModuleTree->Branch("correction", &correction, "correction/F");
138 bklmModuleTree->Branch("error", &error, "error/F");
139 TTree* eklmModuleTree = new TTree("eklm_module",
140 "EKLM module alignment data.");
141 eklmModuleTree->Branch("experiment", &m_experiment, "experiment/I");
142 eklmModuleTree->Branch("run", &m_run, "run/I");
143 eklmModuleTree->Branch("section", &section, "section/I");
144 eklmModuleTree->Branch("sector", &sector, "sector/I");
145 eklmModuleTree->Branch("layer", &layer, "layer/I");
146 eklmModuleTree->Branch("param", &param, "param/I");
147 eklmModuleTree->Branch("value", &value, "value/F");
148 eklmModuleTree->Branch("correction", &correction, "correction/F");
149 eklmModuleTree->Branch("error", &error, "error/F");
150 TTree* eklmSegmentTree = new TTree("eklm_segment",
151 "EKLM segment alignment data.");
152 eklmSegmentTree->Branch("experiment", &m_experiment, "experiment/I");
153 eklmSegmentTree->Branch("run", &m_run, "run/I");
154 eklmSegmentTree->Branch("section", &section, "section/I");
155 eklmSegmentTree->Branch("sector", &sector, "sector/I");
156 eklmSegmentTree->Branch("layer", &layer, "layer/I");
157 eklmSegmentTree->Branch("plane", &plane, "plane/I");
158 eklmSegmentTree->Branch("segment", &segment, "segment/I");
159 eklmSegmentTree->Branch("param", &param, "param/I");
160 eklmSegmentTree->Branch("value", &value, "value/F");
161 eklmSegmentTree->Branch("correction", &correction, "correction/F");
162 eklmSegmentTree->Branch("error", &error, "error/F");
163 const KLMAlignmentData* alignment, *alignmentError, *alignmentCorrection;
164 KLMAlignmentData zeroAlignment(0, 0, 0, 0, 0, 0);
166 for (KLMChannelIndex& klmModule : klmModules) {
167 KLMModuleNumber module = klmModule.getKLMModuleNumber();
168 if (klmModule.getSubdetector() == KLMElementNumbers::c_BKLM) {
169 alignment = bklmAlignment->getModuleAlignment(module);
170 alignmentError = bklmAlignmentErrors->getModuleAlignment(module);
171 alignmentCorrection =
172 bklmAlignmentCorrections->getModuleAlignment(module);
173 } else {
174 alignment = eklmAlignment->getModuleAlignment(module);
175 alignmentError = eklmAlignmentErrors->getModuleAlignment(module);
176 alignmentCorrection =
177 eklmAlignmentCorrections->getModuleAlignment(module);
178 }
179 if (alignment == nullptr)
180 B2FATAL("Incomplete KLM alignment data.");
181 if ((alignmentError == nullptr) && (alignmentCorrection == nullptr)) {
182 B2WARNING("Alignment is not determined for KLM module."
183 << LogVar("Module", module));
184 alignmentError = &zeroAlignment;
185 alignmentCorrection = &zeroAlignment;
186 } else if ((alignmentError == nullptr) ||
187 (alignmentCorrection == nullptr)) {
188 B2FATAL("Inconsistent undtermined parameters.");
189 }
190 section = klmModule.getSection();
191 sector = klmModule.getSector();
192 layer = klmModule.getLayer();
194 value = alignment->getDeltaU();
195 error = alignmentError->getDeltaU();
196 correction = alignmentCorrection->getDeltaU();
197 if (klmModule.getSubdetector() == KLMElementNumbers::c_BKLM)
198 bklmModuleTree->Fill();
199 else
200 eklmModuleTree->Fill();
201 /* cppcheck-suppress redundantAssignment */
203 /* cppcheck-suppress redundantAssignment */
204 value = alignment->getDeltaV();
205 /* cppcheck-suppress redundantAssignment */
206 error = alignmentError->getDeltaV();
207 /* cppcheck-suppress redundantAssignment */
208 correction = alignmentCorrection->getDeltaV();
209 if (klmModule.getSubdetector() == KLMElementNumbers::c_BKLM)
210 bklmModuleTree->Fill();
211 else
212 eklmModuleTree->Fill();
213 /* cppcheck-suppress redundantAssignment */
215 /* cppcheck-suppress redundantAssignment */
216 value = alignment->getDeltaW();
217 /* cppcheck-suppress redundantAssignment */
218 error = alignmentError->getDeltaW();
219 /* cppcheck-suppress redundantAssignment */
220 correction = alignmentCorrection->getDeltaW();
221 if (klmModule.getSubdetector() == KLMElementNumbers::c_BKLM)
222 bklmModuleTree->Fill();
223 else
224 eklmModuleTree->Fill();
225 /* cppcheck-suppress redundantAssignment */
227 /* cppcheck-suppress redundantAssignment */
228 value = alignment->getDeltaAlpha();
229 /* cppcheck-suppress redundantAssignment */
230 error = alignmentError->getDeltaAlpha();
231 /* cppcheck-suppress redundantAssignment */
232 correction = alignmentCorrection->getDeltaAlpha();
233 if (klmModule.getSubdetector() == KLMElementNumbers::c_BKLM)
234 bklmModuleTree->Fill();
235 else
236 eklmModuleTree->Fill();
237 /* cppcheck-suppress redundantAssignment */
239 /* cppcheck-suppress redundantAssignment */
240 value = alignment->getDeltaBeta();
241 /* cppcheck-suppress redundantAssignment */
242 error = alignmentError->getDeltaBeta();
243 /* cppcheck-suppress redundantAssignment */
244 correction = alignmentCorrection->getDeltaBeta();
245 if (klmModule.getSubdetector() == KLMElementNumbers::c_BKLM)
246 bklmModuleTree->Fill();
247 else
248 eklmModuleTree->Fill();
249 /* cppcheck-suppress redundantAssignment */
251 /* cppcheck-suppress redundantAssignment */
252 value = alignment->getDeltaGamma();
253 /* cppcheck-suppress redundantAssignment */
254 error = alignmentError->getDeltaGamma();
255 /* cppcheck-suppress redundantAssignment */
256 correction = alignmentCorrection->getDeltaGamma();
257 if (klmModule.getSubdetector() == KLMElementNumbers::c_BKLM)
258 bklmModuleTree->Fill();
259 else
260 eklmModuleTree->Fill();
261 }
264 eklmSegment = eklmSegments.beginEKLM();
265 eklmSegment.useEKLMSegments();
266 for (; eklmSegment != eklmSegments.endEKLM(); ++eklmSegment) {
267 int eklmSegmentNumber = eklmSegment.getEKLMSegmentNumber();
268 alignment = eklmSegmentAlignment->getSegmentAlignment(eklmSegmentNumber);
269 alignmentError =
270 eklmSegmentAlignmentErrors->getSegmentAlignment(eklmSegmentNumber);
271 alignmentCorrection =
272 eklmSegmentAlignmentCorrections->getSegmentAlignment(eklmSegmentNumber);
273 if (alignment == nullptr)
274 B2FATAL("Incomplete KLM alignment data.");
275 if ((alignmentError == nullptr) && (alignmentCorrection == nullptr)) {
276 /*
277 * The segment alignment is not determined now.
278 * TODO: If it will be turned on in the future, add a warning here.
279 */
280 alignmentError = &zeroAlignment;
281 alignmentCorrection = &zeroAlignment;
282 } else if ((alignmentError == nullptr) ||
283 (alignmentCorrection == nullptr)) {
284 B2FATAL("Inconsistent undtermined parameters.");
285 }
286 section = eklmSegment.getSection();
287 sector = eklmSegment.getSector();
288 layer = eklmSegment.getLayer();
289 plane = eklmSegment.getPlane();
290 segment = eklmSegment.getStrip();
292 value = alignment->getDeltaU();
293 error = alignmentError->getDeltaU();
294 correction = alignmentCorrection->getDeltaU();
295 eklmSegmentTree->Fill();
296 /* cppcheck-suppress redundantAssignment */
298 /* cppcheck-suppress redundantAssignment */
299 value = alignment->getDeltaV();
300 /* cppcheck-suppress redundantAssignment */
301 error = alignmentError->getDeltaV();
302 /* cppcheck-suppress redundantAssignment */
303 correction = alignmentCorrection->getDeltaV();
304 eklmSegmentTree->Fill();
305 /* cppcheck-suppress redundantAssignment */
307 /* cppcheck-suppress redundantAssignment */
308 value = alignment->getDeltaW();
309 /* cppcheck-suppress redundantAssignment */
310 error = alignmentError->getDeltaW();
311 /* cppcheck-suppress redundantAssignment */
312 correction = alignmentCorrection->getDeltaW();
313 eklmSegmentTree->Fill();
314 /* cppcheck-suppress redundantAssignment */
316 /* cppcheck-suppress redundantAssignment */
317 value = alignment->getDeltaAlpha();
318 /* cppcheck-suppress redundantAssignment */
319 error = alignmentError->getDeltaAlpha();
320 /* cppcheck-suppress redundantAssignment */
321 correction = alignmentCorrection->getDeltaAlpha();
322 eklmSegmentTree->Fill();
323 /* cppcheck-suppress redundantAssignment */
325 /* cppcheck-suppress redundantAssignment */
326 value = alignment->getDeltaBeta();
327 /* cppcheck-suppress redundantAssignment */
328 error = alignmentError->getDeltaBeta();
329 /* cppcheck-suppress redundantAssignment */
330 correction = alignmentCorrection->getDeltaBeta();
331 eklmSegmentTree->Fill();
332 /* cppcheck-suppress redundantAssignment */
334 /* cppcheck-suppress redundantAssignment */
335 value = alignment->getDeltaGamma();
336 /* cppcheck-suppress redundantAssignment */
337 error = alignmentError->getDeltaGamma();
338 /* cppcheck-suppress redundantAssignment */
339 correction = alignmentCorrection->getDeltaGamma();
340 eklmSegmentTree->Fill();
341 }
342 bklmModuleTree->Write();
343 eklmModuleTree->Write();
344 eklmSegmentTree->Write();
345 delete bklmModuleTree;
346 delete eklmModuleTree;
347 delete eklmSegmentTree;
348 delete alignmentResults;
349 /* Reset the database. Needed to avoid mess if we call this method multiple times with different GTs. */
351}
Class for accessing objects in the database.
Definition: DBObjPtr.h:21
KLM Alignment data.
float getDeltaU() const
Get shift in U.
float getDeltaV() const
Get shift in V.
@ c_DeltaAlpha
Rotation in alpha.
@ c_DeltaBeta
Rotation in beta.
@ c_DeltaU
Shift in U (EKLM: local X).
@ c_DeltaGamma
Rotation in gamma (EKLM: rotation in local plane).
@ c_DeltaV
Shift in V (EKLM: local Y).
float getDeltaW() const
Get shift in W.
float getDeltaGamma() const
Get rotation in alpha.
float getDeltaAlpha() const
Get rotation in alpha.
float getDeltaBeta() const
Get rotation in alpha.
void resetDatabase()
Reset the database.
void initializeDatabase()
Initialize the database.
std::string m_GlobalTagName
Global Tag name.
std::string m_AlignmentResultsFile
Output file for alignment results.
void printPayloadInformation(DBObjPtr< T > &dbObject)
Print payload information.
KLM channel index.
Class to store variables with their name which were sent to the logging service.
uint16_t KLMModuleNumber
Module number.

◆ checkStripEfficiency()

void checkStripEfficiency ( )

Check strip efficiency.

Definition at line 353 of file KLMCalibrationChecker.cc.

354{
355 /* Initialize the database. */
357 /* Now we can read the payload. */
358 DBObjPtr<KLMStripEfficiency> stripEfficiency;
359 if (!stripEfficiency.isValid())
360 B2FATAL("Strip efficiency data are not valid.");
361 if (m_GlobalTagName != "")
362 printPayloadInformation(stripEfficiency);
363 /* Create trees with strip efficiency measurement results. */
364 int subdetector, section, sector, layer, plane;
365 float efficiency, error;
366 TFile* stripEfficiencyResults =
367 new TFile(m_StripEfficiencyResultsFile.c_str(), "recreate");
368 TTree* efficiencyTree = new TTree("efficiency", "KLM strip efficiency data.");
369 efficiencyTree->Branch("experiment", &m_experiment, "experiment/I");
370 efficiencyTree->Branch("run", &m_run, "run/I");
371 efficiencyTree->Branch("subdetector", &subdetector, "subdetector/I");
372 efficiencyTree->Branch("section", &section, "section/I");
373 efficiencyTree->Branch("sector", &sector, "sector/I");
374 efficiencyTree->Branch("layer", &layer, "layer/I");
375 efficiencyTree->Branch("plane", &plane, "plane/I");
376 efficiencyTree->Branch("efficiency", &efficiency, "efficiency/F");
377 efficiencyTree->Branch("error", &error, "error/F");
379 for (KLMChannelIndex& klmPlane : klmPlanes) {
380 subdetector = klmPlane.getSubdetector();
381 section = klmPlane.getSection();
382 sector = klmPlane.getSector();
383 layer = klmPlane.getLayer();
384 plane = klmPlane.getPlane();
386 subdetector, section, sector, layer, plane, 1);
387 efficiency = stripEfficiency->getEfficiency(channel);
388 error = stripEfficiency->getEfficiencyError(channel);
389 efficiencyTree->Fill();
390 }
391 efficiencyTree->Write();
392 delete efficiencyTree;
393 delete stripEfficiencyResults;
394 /* Reset the database. Needed to avoid mess if we call this method multiple times with different GTs. */
396}
std::string m_StripEfficiencyResultsFile
Output file for alignment results.
KLMChannelNumber channelNumber(int subdetector, int section, int sector, int layer, int plane, int strip) const
Get channel number.
uint16_t KLMChannelNumber
Channel number.

◆ checkTimeCableDelay()

void checkTimeCableDelay ( )

Check time cable delay.

Definition at line 473 of file KLMCalibrationChecker.cc.

474{
475 /* Initialize the database. */
477 /* Now we can read the payload. */
478 DBObjPtr<KLMTimeCableDelay> timeCableDelay;
479 if (!timeCableDelay.isValid())
480 B2FATAL("Time Cable delay data are not valid.");
481 if (m_GlobalTagName != "")
482 printPayloadInformation(timeCableDelay);
483 /* Create trees with time cable delay measurement results. */
484 int subdetector, section, sector, layer, plane, strip, channelNumber;
485 double timeDelay;
486 TFile* timeCableDelayResults =
487 new TFile(m_TimeCableDelayResultsFile.c_str(), "recreate");
488 TTree* cableDelayTree = new TTree("cabledelay", "KLM timecabledelay data");
489 cableDelayTree->Branch("experiment", &m_experiment, "experiment/I");
490 cableDelayTree->Branch("run", &m_run, "run/I");
491 cableDelayTree->Branch("subdetector", &subdetector, "subdetector/I");
492 cableDelayTree->Branch("section", &section, "section/I");
493 cableDelayTree->Branch("sector", &sector, "sector/I");
494 cableDelayTree->Branch("layer", &layer, "layer/I");
495 cableDelayTree->Branch("plane", &plane, "plane/I");
496 cableDelayTree->Branch("strip", &strip, "strip/I");
497 cableDelayTree->Branch("channelNumber", &channelNumber, "channelNumber/I");
498 cableDelayTree->Branch("timeDelay", &timeDelay, "timeDelay/D");
500 for (KLMChannelIndex& klmStrip : klmStrips) {
501 subdetector = klmStrip.getSubdetector();
502 section = klmStrip.getSection();
503 sector = klmStrip.getSector();
504 layer = klmStrip.getLayer();
505 plane = klmStrip.getPlane();
506 strip = klmStrip.getStrip();
508 subdetector, section, sector, layer, plane, strip);
509 timeDelay = timeCableDelay->getTimeDelay(channel);
510 channelNumber = m_ElementNumbers->channelNumber(subdetector, section, sector, layer, plane, strip);
511 cableDelayTree->Fill();
512 }
513 cableDelayTree->Write();
514 delete cableDelayTree;
515 delete timeCableDelayResults;
516 /* Reset the database. Needed to avoid mess if we call this method multiple times with different GTs. */
518}
std::string m_TimeCableDelayResultsFile
Output file for time cable delay results.

◆ checkTimeConstants()

void checkTimeConstants ( )

Check time constants.

Definition at line 520 of file KLMCalibrationChecker.cc.

521{
522 /* Initialize the database. */
524 /* Now we can read the payload. */
525 DBObjPtr<KLMTimeConstants> timeConstants;
526 if (!timeConstants.isValid())
527 B2FATAL("Time Constants data are not valid.");
528 if (m_GlobalTagName != "")
529 printPayloadInformation(timeConstants);
530 /* Create trees with time constant measurement results. */
531 int subdetector, section, sector, layer, plane, strip, module, channelNumber;
532 float delayEKLM, delayBKLM, delayRPCPhi, delayRPCZ;
533 TFile* timeConstantsResults =
534 new TFile(m_TimeConstantsResultsFile.c_str(), "recreate");
535 TTree* constantsTree = new TTree("constants", "KLM timeConstants data");
536 constantsTree->Branch("experiment", &m_experiment, "experiment/I");
537 constantsTree->Branch("run", &m_run, "run/I");
538 constantsTree->Branch("subdetector", &subdetector, "subdetector/I");
539 constantsTree->Branch("section", &section, "section/I");
540 constantsTree->Branch("sector", &sector, "sector/I");
541 constantsTree->Branch("layer", &layer, "layer/I");
542 constantsTree->Branch("plane", &plane, "plane/I");
543 constantsTree->Branch("strip", &strip, "strip/I");
544 constantsTree->Branch("module", &module, "module/I");
545 constantsTree->Branch("channelNumber", &channelNumber, "channelNumber/I");
546 constantsTree->Branch("delayEKLM", &delayEKLM, "delayEKLM/F");
547 constantsTree->Branch("delayBKLM", &delayBKLM, "delayBKLM/F");
548 constantsTree->Branch("delayRPCPhi", &delayRPCPhi, "delayRPCPhi/F");
549 constantsTree->Branch("delayRPCZ", &delayRPCZ, "delayRPCZ/F");
551 for (KLMChannelIndex& klmStrip : klmStrips) {
552 subdetector = klmStrip.getSubdetector();
553 section = klmStrip.getSection();
554 sector = klmStrip.getSector();
555 layer = klmStrip.getLayer();
556 plane = klmStrip.getPlane();
557 strip = klmStrip.getStrip();
558 channelNumber = m_ElementNumbers->channelNumber(subdetector, section, sector, layer, plane, strip);
559 module = m_ElementNumbers->moduleNumber(subdetector, section, sector, layer);
560 delayEKLM = timeConstants->getDelay(KLMTimeConstants::c_EKLM);
561 delayBKLM = timeConstants->getDelay(KLMTimeConstants::c_BKLM);
562 delayRPCPhi = timeConstants->getDelay(KLMTimeConstants::c_RPCPhi);
563 delayRPCZ = timeConstants->getDelay(KLMTimeConstants::c_RPCZ);
564 constantsTree->Fill();
565 }
566 constantsTree->Write();
567 delete constantsTree;
568 delete timeConstantsResults;
569 /* Reset the database. Needed to avoid mess if we call this method multiple times with different GTs. */
571}
std::string m_TimeConstantsResultsFile
Output file for time constants results.
KLMModuleNumber moduleNumber(int subdetector, int section, int sector, int layer) const
Get module number.
@ c_BKLM
BKLM scintillator.
@ c_EKLM
EKLM scintillator.

◆ createStripEfficiencyHistograms()

void createStripEfficiencyHistograms ( )

Create strip efficiency histograms.

Definition at line 398 of file KLMCalibrationChecker.cc.

399{
400 /* Initialize the database. */
402 /* Now we can read the payload. */
403 DBObjPtr<KLMStripEfficiency> stripEfficiency;
404 if (not stripEfficiency.isValid())
405 B2FATAL("Strip efficiency data are not valid.");
406 if (m_GlobalTagName != "")
407 printPayloadInformation(stripEfficiency);
408 /* Finally, loop over KLM sectors to check the efficiency. */
410 TCanvas* canvas = new TCanvas();
411 for (KLMChannelIndex& klmSector : klmSectors) {
412 int subdetector = klmSector.getSubdetector();
413 int section = klmSector.getSection();
414 int sector = klmSector.getSector();
415 /* Setup the histogram. */
416 TH1F* hist = new TH1F("plane_histogram", "", 30, 0.5, 30.5);
417 hist->GetYaxis()->SetTitle("Efficiency");
418 hist->SetMinimum(0.4);
419 hist->SetMaximum(1.);
420 hist->SetMarkerStyle(20);
421 hist->SetMarkerSize(0.5);
422 TString title;
423 if (subdetector == KLMElementNumbers::c_BKLM) {
425 title.Form("BKLM backward sector %d", sector);
426 else
427 title.Form("BKLM forward sector %d", sector);
428 hist->SetTitle(title.Data());
429 hist->GetXaxis()->SetTitle("(Layer - 1) * 2 + plane + 1");
430 for (int layer = 1; layer <= BKLMElementNumbers::getMaximalLayerNumber(); layer++) {
431 for (int plane = 0; plane <= BKLMElementNumbers::getMaximalPlaneNumber(); plane++) {
432 int bin = (layer - 1) * 2 + plane + 1;
433 float efficiency = stripEfficiency->getBarrelEfficiency(section, sector, layer, plane, 2);
434 float efficiencyError = stripEfficiency->getBarrelEfficiencyError(section, sector, layer, plane, 2);
435 hist->SetBinContent(bin, efficiency);
436 hist->SetBinError(bin, efficiencyError);
437 }
438 }
439 } else {
441 hist->SetBins(24, 0.5, 24.5);
442 title.Form("EKLM backward sector %d", sector);
443 } else {
444 hist->SetBins(28, 0.5, 28.5);
445 title.Form("EKLM forward sector %d", sector);
446 }
447 hist->SetTitle(title.Data());
448 hist->GetXaxis()->SetTitle("(Layer - 1) * 2 + plane");
449 const EKLMElementNumbers* elementNumbersEKLM = &(EKLMElementNumbers::Instance());
450 for (int layer = 1; layer <= elementNumbersEKLM->getMaximalDetectorLayerNumber(section); layer++) {
451 for (int plane = 1; plane <= EKLMElementNumbers::getMaximalPlaneNumber(); plane++) {
452 int bin = (layer - 1) * 2 + plane;
453 float efficiency = stripEfficiency->getEndcapEfficiency(section, sector, layer, plane, 2);
454 float efficiencyError = stripEfficiency->getEndcapEfficiencyError(section, sector, layer, plane, 2);
455 hist->SetBinContent(bin, efficiency);
456 hist->SetBinError(bin, efficiencyError);
457 }
458 }
459 }
460 hist->Draw("e");
461 TString name;
462 name.Form("efficiency_subdetector_%d_section_%d_sector_%d.pdf", subdetector, section, sector);
463 canvas->Print(name.Data());
464 canvas->Update();
465 delete hist;
466 }
467 delete canvas;
468 /* Reset the database. Needed to avoid mess if we call this method multiple times with different GTs. */
470}
static constexpr int getMaximalLayerNumber()
Get maximal layer number (1-based).
static constexpr int getMaximalPlaneNumber()
Get maximal plane number (0-based).
EKLM element numbers.
int getMaximalDetectorLayerNumber(int section) const
Get maximal detector layer number.
static const EKLMElementNumbers & Instance()
Instantiation.
static constexpr int getMaximalPlaneNumber()
Get maximal plane number.

◆ initializeDatabase()

void initializeDatabase ( )
private

Initialize the database.

Definition at line 59 of file KLMCalibrationChecker.cc.

60{
61 /* Mimic a module initialization. */
63 m_EventMetaData.registerInDataStore();
65 if (!m_EventMetaData.isValid())
66 m_EventMetaData.construct(1, m_run, m_experiment);
67 /* Database instance and configuration. */
68 DBStore& dbStore = DBStore::Instance();
69 dbStore.update();
70 dbStore.updateEvent();
71 auto& dbConfiguration = Conditions::Configuration::getInstance();
72 if ((m_testingPayloadName != "") and (m_GlobalTagName == ""))
73 dbConfiguration.prependTestingPayloadLocation(m_testingPayloadName);
74 else if ((m_testingPayloadName == "") and (m_GlobalTagName != ""))
75 dbConfiguration.prependGlobalTag(m_GlobalTagName);
76 else
77 B2FATAL("Setting both testing payload and Global Tag or setting no one of them.");
78}
static Configuration & getInstance()
Get a reference to the instance which will be used when the Database is initialized.
Singleton class to cache database objects.
Definition: DBStore.h:31
static DataStore & Instance()
Instance of singleton Store.
Definition: DataStore.cc:54
void setInitializeActive(bool active)
Setter for m_initializeActive.
Definition: DataStore.cc:94
std::string m_testingPayloadName
Testing payload location.
StoreObjPtr< EventMetaData > m_EventMetaData
Event metadata.
static DBStore & Instance()
Instance of a singleton DBStore.
Definition: DBStore.cc:28
void updateEvent()
Updates all intra-run dependent objects.
Definition: DBStore.cc:142
void update()
Updates all objects that are outside their interval of validity.
Definition: DBStore.cc:79

◆ printPayloadInformation()

void printPayloadInformation ( DBObjPtr< T > &  dbObject)
inlineprivate

Print payload information.

Definition at line 144 of file KLMCalibrationChecker.h.

145 {
146 B2INFO("Analyzing the following payload:"
147 << LogVar("Global Tag", m_GlobalTagName.c_str())
148 << LogVar("Name", dbObject.getName())
149 << LogVar("Revision", dbObject.getRevision())
150 << LogVar("IoV", dbObject.getIoV()));
151 }

◆ resetDatabase()

void resetDatabase ( )
private

Reset the database.

Definition at line 80 of file KLMCalibrationChecker.cc.

81{
82 /* Reset both DataStore and Database. */
85 DBStore::Instance().reset(false);
86}
void reset(EDurability durability)
Frees memory occupied by data store items and removes all objects from the map.
Definition: DataStore.cc:86
void reset(bool keepEntries=false)
Invalidate all payloads.
Definition: DBStore.cc:177
static Database & Instance()
Instance of a singleton Database.
Definition: Database.cc:42
static void reset(bool keepConfig=false)
Reset the database instance.
Definition: Database.cc:50

◆ setAlignmentResultsFile()

void setAlignmentResultsFile ( const std::string &  alignmentResultsFile)
inline

Set alignment results file.

Definition at line 70 of file KLMCalibrationChecker.h.

71 {
72 m_AlignmentResultsFile = alignmentResultsFile;
73 }

◆ setExperimentRun()

void setExperimentRun ( int  experiment,
int  run 
)

Set experiment and run numbers.

Definition at line 49 of file KLMCalibrationChecker.cc.

50{
51 m_experiment = experiment;
52 m_run = run;
53 if (m_EventMetaData.isValid()) {
54 m_EventMetaData->setExperiment(experiment);
55 m_EventMetaData->setRun(run);
56 }
57}

◆ setGlobalTag()

void setGlobalTag ( const std::string &  globalTagName)
inline

Set Global Tag name.

Definition at line 62 of file KLMCalibrationChecker.h.

63 {
64 m_GlobalTagName = globalTagName;
65 }

◆ setStripEfficiencyResultsFile()

void setStripEfficiencyResultsFile ( const std::string &  stripEfficiencyResultsFile)
inline

Set strip efficiency results file.

Definition at line 78 of file KLMCalibrationChecker.h.

80 {
81 m_StripEfficiencyResultsFile = stripEfficiencyResultsFile;
82 }

◆ setTestingPayload()

void setTestingPayload ( const std::string &  testingPayloadName)
inline

Set testing payload name.

Definition at line 54 of file KLMCalibrationChecker.h.

55 {
56 m_testingPayloadName = testingPayloadName;
57 }

◆ setTimeCableDelayResultsFile()

void setTimeCableDelayResultsFile ( const std::string &  timeCableDelayResultsFile)
inline

Set time cable delay results file.

Definition at line 87 of file KLMCalibrationChecker.h.

89 {
90 m_TimeCableDelayResultsFile = timeCableDelayResultsFile;
91 }

◆ setTimeConstantsResultsFile()

void setTimeConstantsResultsFile ( const std::string &  timeConstantsResultsFile)
inline

Set time constants result file.

Definition at line 96 of file KLMCalibrationChecker.h.

98 {
99 m_TimeConstantsResultsFile = timeConstantsResultsFile;
100 }

Member Data Documentation

◆ m_AlignmentResultsFile

std::string m_AlignmentResultsFile = "alignment.root"
private

Output file for alignment results.

Definition at line 166 of file KLMCalibrationChecker.h.

◆ m_ElementNumbers

const KLMElementNumbers* m_ElementNumbers
private

Element numbers.

Definition at line 178 of file KLMCalibrationChecker.h.

◆ m_EventMetaData

StoreObjPtr<EventMetaData> m_EventMetaData
private

Event metadata.

Definition at line 181 of file KLMCalibrationChecker.h.

◆ m_experiment

int m_experiment
private

Experiment number.

Definition at line 154 of file KLMCalibrationChecker.h.

◆ m_GlobalTagName

std::string m_GlobalTagName = ""
private

Global Tag name.

Definition at line 163 of file KLMCalibrationChecker.h.

◆ m_run

int m_run
private

Run number.

Definition at line 157 of file KLMCalibrationChecker.h.

◆ m_StripEfficiencyResultsFile

std::string m_StripEfficiencyResultsFile = "strip_efficiency.root"
private

Output file for alignment results.

Definition at line 169 of file KLMCalibrationChecker.h.

◆ m_testingPayloadName

std::string m_testingPayloadName = ""
private

Testing payload location.

Definition at line 160 of file KLMCalibrationChecker.h.

◆ m_TimeCableDelayResultsFile

std::string m_TimeCableDelayResultsFile = "timeCableDelay.root"
private

Output file for time cable delay results.

Definition at line 172 of file KLMCalibrationChecker.h.

◆ m_TimeConstantsResultsFile

std::string m_TimeConstantsResultsFile = "timeConstants.root"
private

Output file for time constants results.

Definition at line 175 of file KLMCalibrationChecker.h.


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