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 setEventT0HitResolutionResultsFile (const std::string &eventT0HitResolutionResultsFile)
 Set EventT0 hit resolution results 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.
 
void checkEventT0HitResolution ()
 Check EventT0 hit resolution.
 

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.
 
std::string m_EventT0HitResolutionResultsFile = "eventT0HitResolution.root"
 Output file for EventT0 hit resolution 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 39 of file KLMCalibrationChecker.cc.

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

◆ ~KLMCalibrationChecker()

Destructor.

Definition at line 46 of file KLMCalibrationChecker.cc.

47{
48}

Member Function Documentation

◆ checkAlignment()

void checkAlignment ( )

Check alignment.

Definition at line 89 of file KLMCalibrationChecker.cc.

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

◆ checkEventT0HitResolution()

void checkEventT0HitResolution ( )

Check EventT0 hit resolution.

Definition at line 574 of file KLMCalibrationChecker.cc.

575{
576 /* Initialize the database. */
578 /* Now we can read the payload. */
579 DBObjPtr<KLMEventT0HitResolution> eventT0HitResolution;
580 if (!eventT0HitResolution.isValid())
581 B2FATAL("EventT0 Hit Resolution data are not valid.");
582 if (m_GlobalTagName != "")
583 printPayloadInformation(eventT0HitResolution);
584 /* Create tree with EventT0 hit resolution (one row per detector category). */
585 int category;
586 float sigma, sigmaErr;
587 TFile* eventT0HitResolutionResults =
588 new TFile(m_EventT0HitResolutionResultsFile.c_str(), "recreate");
589 TTree* resolutionTree = new TTree("eventT0HitResolution", "KLM EventT0 hit resolution data");
590 resolutionTree->Branch("experiment", &m_experiment, "experiment/I");
591 resolutionTree->Branch("run", &m_run, "run/I");
592 resolutionTree->Branch("category", &category, "category/I");
593 resolutionTree->Branch("sigma", &sigma, "sigma/F");
594 resolutionTree->Branch("sigmaErr", &sigmaErr, "sigmaErr/F");
595 /* Category order must match KLMEventT0HitResolution::Category enum. */
596 const int categories[5] = {
602 };
603 for (int cat : categories) {
604 category = cat;
605 sigma = eventT0HitResolution->getSigma(cat);
606 sigmaErr = eventT0HitResolution->getSigmaErr(cat);
607 resolutionTree->Fill();
608 }
609 resolutionTree->Write();
610 delete resolutionTree;
611 delete eventT0HitResolutionResults;
612 /* Reset the database. Needed to avoid mess if we call this method multiple times with different GTs. */
614}
bool isValid() const
Check whether a valid object was obtained from the database.
std::string m_EventT0HitResolutionResultsFile
Output file for EventT0 hit resolution results.
@ c_RPC
BKLM RPC (combined phi/z).

◆ checkStripEfficiency()

void checkStripEfficiency ( )

Check strip efficiency.

Definition at line 354 of file KLMCalibrationChecker.cc.

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

◆ checkTimeCableDelay()

void checkTimeCableDelay ( )

Check time cable delay.

Definition at line 474 of file KLMCalibrationChecker.cc.

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

◆ checkTimeConstants()

void checkTimeConstants ( )

Check time constants.

Definition at line 521 of file KLMCalibrationChecker.cc.

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

◆ createStripEfficiencyHistograms()

void createStripEfficiencyHistograms ( )

Create strip efficiency histograms.

Definition at line 399 of file KLMCalibrationChecker.cc.

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

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

◆ printPayloadInformation()

template<class T>
void printPayloadInformation ( DBObjPtr< T > & dbObject)
inlineprivate

Print payload information.

Definition at line 158 of file KLMCalibrationChecker.h.

159 {
160 B2INFO("Analyzing the following payload:"
161 << LogVar("Global Tag", m_GlobalTagName.c_str())
162 << LogVar("Name", dbObject.getName())
163 << LogVar("Revision", dbObject.getRevision())
164 << LogVar("IoV", dbObject.getIoV()));
165 }

◆ resetDatabase()

void resetDatabase ( )
private

Reset the database.

Definition at line 81 of file KLMCalibrationChecker.cc.

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

◆ 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 }

◆ setEventT0HitResolutionResultsFile()

void setEventT0HitResolutionResultsFile ( const std::string & eventT0HitResolutionResultsFile)
inline

Set EventT0 hit resolution results file.

Definition at line 105 of file KLMCalibrationChecker.h.

107 {
108 m_EventT0HitResolutionResultsFile = eventT0HitResolutionResultsFile;
109 }

◆ setExperimentRun()

void setExperimentRun ( int experiment,
int run )

Set experiment and run numbers.

Definition at line 50 of file KLMCalibrationChecker.cc.

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

◆ 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 180 of file KLMCalibrationChecker.h.

◆ m_ElementNumbers

const KLMElementNumbers* m_ElementNumbers
private

Element numbers.

Definition at line 195 of file KLMCalibrationChecker.h.

◆ m_EventMetaData

StoreObjPtr<EventMetaData> m_EventMetaData
private

Event metadata.

Definition at line 198 of file KLMCalibrationChecker.h.

◆ m_EventT0HitResolutionResultsFile

std::string m_EventT0HitResolutionResultsFile = "eventT0HitResolution.root"
private

Output file for EventT0 hit resolution results.

Definition at line 192 of file KLMCalibrationChecker.h.

◆ m_experiment

int m_experiment
private

Experiment number.

Definition at line 168 of file KLMCalibrationChecker.h.

◆ m_GlobalTagName

std::string m_GlobalTagName = ""
private

Global Tag name.

Definition at line 177 of file KLMCalibrationChecker.h.

◆ m_run

int m_run
private

Run number.

Definition at line 171 of file KLMCalibrationChecker.h.

◆ m_StripEfficiencyResultsFile

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

Output file for alignment results.

Definition at line 183 of file KLMCalibrationChecker.h.

◆ m_testingPayloadName

std::string m_testingPayloadName = ""
private

Testing payload location.

Definition at line 174 of file KLMCalibrationChecker.h.

◆ m_TimeCableDelayResultsFile

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

Output file for time cable delay results.

Definition at line 186 of file KLMCalibrationChecker.h.

◆ m_TimeConstantsResultsFile

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

Output file for time constants results.

Definition at line 189 of file KLMCalibrationChecker.h.


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