ATLAS Offline Software
Loading...
Searching...
No Matches
LumiCalculator.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
6
7// GoodRunsLists
10
11//
12#include "LumiCalc/CoolQuery.h"
14
15// ROOT
16#include "TH1F.h"
17#include "TTree.h"
18#include "TFile.h"
19#include "TString.h"
20#include "TObjString.h"
21#include "TROOT.h"
22
23// stl includes
24#include <memory>
25#include <iomanip>
26#include <iostream>
27#include <set>
28#include <regex.h>
29#include <stdexcept>
30
32 : m_LumiTree (nullptr)
33 , m_recordTTree (false)
34 , m_State (true)
35
36 , m_trigger ("COOLONL_TRIGGER/")
37 , m_lumioff ("COOLOFL_TRIGGER/")
38 , m_lumionl ("COOLONL_TRIGGER/")
39 , m_lumitag ("OflLumi-8TeV-002")// for offline: OflLumi_CosmicFake, OflLumi_TopMix
40 , m_lumimethod ("")// offline channels: ATLAS_PREFERRED, OflLumi_Fake0, OflLumi_Fake:, TopMixLumi
41 , m_laroff ("COOLOFL_LAR/")
42 , m_bsonl ("COOLONL_INDET/")
43 , m_bstag ("IndetBeamposOnl-HLT-UPD1-001-00")
44
45 , m_parlvl1menufolder ("/TRIGGER/LVL1/Menu")
46 , m_parhltmenufolder ("/TRIGGER/HLT/Menu")// ChainCounter is here for COOLONL_TRIGGER/COMP20
47 , m_parhltprescalesfolder ("/TRIGGER/HLT/Prescales")// ChainCounter is here for COOLONL_TRIGGER/COMP20
48 , m_parlumilvl1folder ("/TRIGGER/LUMI/LVL1COUNTERS")
49 , m_parlumihltfolder ("/TRIGGER/LUMI/HLTCOUNTERS")
50 , m_parlvl1prescalesfolder ("/TRIGGER/LVL1/Prescales")
51 , m_parlvl1lblbfolder ("/TRIGGER/LUMI/LBLB")// for time information
52 , m_parlareventvetofolder ("/LAR/BadChannelsOfl/EventVeto") // For LAr event veto
53 , m_paronlbeamspotfolder ("/Indet/Onl/Beampos") // For invalid online beamspot
54
55 , m_logger( "LumiCalculator" )
56 , m_lbcollname ("LumiBlocks")
57 , m_uselivetrigger (false)
58 , m_verbose (false)
59 , m_lbstarttime(0.)
60 , m_lbendtime(0.)
61
62 , m_totalDelL(0.)
63 , m_totalL(0.)
64 , m_totalLRun(0.)
65 , m_totaltime(0.)
66 , m_instLumi(0.)
67 , m_AvEvtsPerBX(0.)
68 , m_delLumi(0.)
69 , m_intLumi(0.)
70 , m_deltaT(0.)
71 , m_TotaldeltaT(0.)
72 , m_livefrac(0.)
73 , m_livetime(0.)
74 , m_l1acc(0)
76 , m_l1accof(false)
77 , m_l2acc(0)
78 , m_l3acc(0)
79 , m_totall1acc(0)
82 , m_totall2acc(0)
83 , m_totall3acc(0)
84 , m_l1prescale(0.)
85 , m_l2prescale(0.)
86 , m_l3prescale(0.)
90 , m_afterprescaleof(false)
92 , m_beforeprescaleof(false)
93 , m_runnbr(0)
94 , m_lbstart(0)
95 , m_lbstop(0)
97 , m_lbstop_prev(0)
98 , m_runnbr_prev(0)
100 , m_totalbadblock(0)
101 , m_clumiblocknbr(0)
103 , m_triglevel(0)
106 , m_totalPrescale(0.)
108 , m_lumiWOPrescale(0.)
110 , m_lumiLAr(0.)
111 , m_t_lumiLAr(0.)
112 , m_t_totalDelL(0.)
113 , m_t_totalL(0.)
114 , m_t_totalLRun(0.)
115 , m_t_totaltime(0.)
116 , m_t_deltaT(0.)
117 , m_t_l1acc(0)
118 , m_t_l2acc(0)
119 , m_t_l3acc(0)
123 , m_lartime(0.)
124 , m_larfrac(0.)
125 , m_bsvalid(0.)
126 , m_effxsec(1.)
127 , m_l1rate(0.)
128 , m_l2rate(0.)
129 , m_l3rate(0.)
136 , m_mintrigrate(5./120.)
137 , m_collsgrl(0)
138 , m_ntrigplb(0)
139 , m_trigrateplb(0)
140 , m_lumiplb(0)
142 , m_intlumi(0)
146 , m_intlumiruns(0)
149 , m_avgintperbx(0)
150 , m_makePlots(false)
151 , m_makecollList(false)
152 , m_Lumiid(0)
153 , m_L3id(0)
154 , m_L2id(0)
155 , m_L1id(0)
156 , m_LiveL1id(0)
157 , m_L1Valid(false)
158 , m_L2Valid(false)
159 , m_L3Valid(false)
160 , m_LiveValid(false)
161 , m_onlinelumi(false)
162 , m_uselar(false)
163 , m_usebs(false)
164 , m_minrun(0)
165 , m_maxrun(0)
166{
167
168 // by default we use the "offline data" database name
173
174 m_lumichannel = 0;
175
176}
177
178
180
181 // delete collisions grl
182 if (m_collsgrl!=0) { delete m_collsgrl; m_collsgrl=0; }
183
184 std::vector<TH1F*>::iterator itr;
185 for (itr=m_ntrigplbVec.begin(); itr!=m_ntrigplbVec.end(); ++itr) { delete *itr; }
186 for (itr=m_trigrateplbVec.begin(); itr!=m_trigrateplbVec.end(); ++itr) { delete *itr; }
187 for (itr=m_lumiplbVec.begin(); itr!=m_lumiplbVec.end(); ++itr) { delete *itr; }
188 for (itr=m_lumitrigrateplbVec.begin(); itr!=m_lumitrigrateplbVec.end(); ++itr) { delete *itr; }
189 for (itr=m_intlumiVec.begin(); itr!=m_intlumiVec.end(); ++itr) { delete *itr; }
190 for (itr=m_intlumitrigrateVec.begin(); itr!=m_intlumitrigrateVec.end(); ++itr) { delete *itr; }
191 for (itr=m_lumitrigrateplb_recordedVec.begin(); itr!=m_lumitrigrateplb_recordedVec.end(); ++itr) { delete *itr; }
192 for (itr=m_intlumitrigrate_recordedVec.begin(); itr!=m_intlumitrigrate_recordedVec.end(); ++itr) { delete *itr; }
193
194 m_ntrigplbVec.clear();
195 m_trigrateplbVec.clear();
196 m_lumiplbVec.clear();
197 m_lumitrigrateplbVec.clear();
198 m_intlumiVec.clear();
199 m_intlumitrigrateVec.clear();
202
203 if (m_intlumiruns!=0) { delete m_intlumiruns; m_intlumiruns=0; }
206 if (m_avgintperbx != 0) { delete m_avgintperbx; m_avgintperbx = 0; }
207
208}
209
211 // Register branches
213 if(m_LumiTree != 0){
214
215 m_LumiTree->Branch("LBStartTime/D", &m_lbstarttime);
216 m_LumiTree->Branch("LBEndTime/D", &m_lbendtime);
217 m_LumiTree->Branch("Trigger", &m_triggerchain);
218 m_LumiTree->Branch("LBCollection", &m_lbcollname);
219 m_LumiTree->Branch("RunNbr", &m_runnbr);
220 m_LumiTree->Branch("IOVRStart", &m_lbstart);
221 m_LumiTree->Branch("IOVREnd", &m_lbstop);
222 m_LumiTree->Branch("LBStart", &m_clumiblocknbr);
223 m_LumiTree->Branch("Inst_m_Lumi", &m_instLumi);
224 m_LumiTree->Branch("LiveTime", &m_livetime);
225 m_LumiTree->Branch("L1Presc", &m_l1prescale);
226 m_LumiTree->Branch("L2Presc", &m_l2prescale);
227 m_LumiTree->Branch("L3Presc", &m_l3prescale);
228 m_LumiTree->Branch("L1Count", &m_l1acc);
229 m_LumiTree->Branch("L1CountOverFlow", &m_l1accof);
230 m_LumiTree->Branch("L2Count", &m_l2acc);
231 m_LumiTree->Branch("L3Count", &m_l3acc);
232 m_LumiTree->Branch("L1AfterPrescale", &m_afterprescale);
233 m_LumiTree->Branch("L1AfterPrescaleOverFlow", &m_afterprescaleof);
234 m_LumiTree->Branch("L1BeforePrescale", &m_beforeprescale);
235 m_LumiTree->Branch("L1BeforePrescaleOverFlow", &m_beforeprescaleof);
236 m_LumiTree->Branch("Livefrac", &m_livefrac);
237 m_LumiTree->Branch("LArfrac", &m_larfrac);
238 m_LumiTree->Branch("DeltaT", &m_deltaT);
239 m_LumiTree->Branch("L1Rate", &m_l1rate);
240 m_LumiTree->Branch("IntLumi",&m_intLumi);
241 m_LumiTree->Branch("L1Ratediveffxsec",&m_l1ratediveffxsec);
242 m_LumiTree->Branch("TotalLumi",&m_totalL);
243 m_LumiTree->Branch("Total_L1Ratediveffxsec",&m_total_l1ratediveffxsec);
244 m_LumiTree->Branch("TotalLumiRun",&m_totalLRun);
245 m_LumiTree->Branch("Total_L1RatediveffxsecRun",&m_total_l1ratediveffxsecRun);
246 m_LumiTree->Branch("L1RatediveffxsecRecorded",&m_l1ratediveffxsec_recorded);
247 m_LumiTree->Branch("Total_L1RatediveffxsecRecorded",&m_total_l1ratediveffxsec_recorded);
248 m_LumiTree->Branch("Total_L1RatediveffxsecRunRecorded",&m_total_l1ratediveffxsecRun_recorded);
249 m_LumiTree->Branch("AvergeInteractionPerXing",&m_AvEvtsPerBX);
250 m_LumiTree->Branch("BSValid", &m_bsvalid);
251 }
252
253}
254
255
256void LumiCalculator::SetCollName(const std::string& lbcollname){
257 m_lbcollname = lbcollname;
258}
259
260
264
266
267 // Print warning
268 if (mc)
269 m_logger << Root::kWARNING << "Monte Carlo mode no longer supported!" << Root::GEndl;
270
271}
272
274
275 m_onlinelumi = online;
276
277}
278
279void LumiCalculator::UseLumiTag(const std::string& tag){
280 m_lumitag = tag;
281}
282
283void LumiCalculator::UseLumiMethod(const std::string& method){
284 m_lumimethod = method;
285}
286
288 m_lumichannel = chan;
289 m_lumimethod = "";
290}
291
292
293void LumiCalculator::UseLiveTrigger(bool live, std::string& livetrigger){
294 m_uselivetrigger = live;
296}
297
298void LumiCalculator::UseLArNoiseDB(bool lar, const std::string& lartag) {
299 m_uselar = lar;
300 m_lartag = lartag;
301}
302
303void LumiCalculator::UseBeamspot(bool bs, const std::string& bstag) {
304 m_usebs = bs;
305 m_bstag = bstag;
306}
307
309 if(m_LumiTree != 0)return m_LumiTree;
310 return 0;
311}
312
313//______________________________________________________________________________
314void LumiCalculator::IntegrateLumi ATLAS_NOT_THREAD_SAFE (const xAOD::LumiBlockRangeContainer * iovc, const std::string& triggerchain){
315
316
317 CoolQuery * cq_lumi = NULL;
318 CoolQuery * cq_trigger = NULL;
319 CoolQuery* cq_lar = NULL;
320 CoolQuery* cq_bs = NULL;
321
322 // Peek at first run in iovc to check if we want Run1 or Run2 data
324 // const IOVRange * iovr = (*it);
325 IOVRange iovr2 (IOVTime((*it)->startRunNumber(),(*it)->startLumiBlockNumber()),
326 IOVTime((*it)->stopRunNumber(),(*it)->stopLumiBlockNumber()));
327
328 bool isrun2 = false;
329 std::string onlfolder;
330 std::string oflfolder;
331 if (iovr2.start().run() > 222222) {
332 isrun2 = true;
333 m_data_db="CONDBR2";
334 onlfolder = "/TRIGGER/LUMI/OnlPrefLumi";
335 oflfolder = "/TRIGGER/OFLLUMI/OflPrefLumi";
336 } else {
337 isrun2 = false;
338 m_data_db="COMP200";
339 onlfolder = "/TRIGGER/LUMI/LBLESTONL";
340 oflfolder = "/TRIGGER/OFLLUMI/LBLESTOFL";
341 }
342
343 // define all connection strings
344
345 // check for online DB
346 if(m_onlinelumi){
347 m_lumi_database = m_lumionl + m_data_db;
348 m_parlumiestfolder = std::move(onlfolder);
349 } else {
350 m_lumi_database = m_lumioff + m_data_db;
351 m_parlumiestfolder = std::move(oflfolder);
352 }
353
354 m_trig_database = m_trigger + m_data_db;
355 m_lar_database = m_laroff + m_data_db;
356 m_bs_database = m_bsonl + m_data_db;
357
358
359 m_logger << Root::kINFO << "Luminosity database: " << m_lumi_database << Root::GEndl;
360 m_logger << Root::kINFO << "Trigger database: " << m_trig_database << Root::GEndl;
361
362 cq_lumi = new CoolQuery(m_lumi_database, triggerchain);
363 if (!cq_lumi->openDbConn()) {
364 delete cq_lumi;
365 return;
366 }
367
368 cq_trigger = new CoolQuery(m_trig_database, triggerchain);
369 if (!cq_trigger->openDbConn()) {
370 delete cq_trigger;
371 return;
372 }
373
374 if (m_uselar) {
375 m_logger << Root::kINFO << "LAr database: " << m_lar_database << Root::GEndl;
376
377 cq_lar = new CoolQuery(m_lar_database, triggerchain);
378 if (!cq_lar->openDbConn()) {
379 delete cq_lar;
380 return;
381 }
382 }
383
384 if (m_usebs) {
385 m_logger << Root::kINFO << "Onl beamspot database: " << m_bs_database << Root::GEndl;
386
387 cq_bs = new CoolQuery(m_bs_database, triggerchain);
388 if (!cq_bs->openDbConn()) {
389 delete cq_bs;
390 return;
391 }
392 }
393
394 // initialize
395 m_lbstarttime = 0.;
396 m_lbendtime = 0.;
397 m_totalDelL = 0.;
398 m_totalL = 0.;
399 m_totalLRun = 0.;
400 m_totaltime = 0.;
401 m_instLumi = 0.;
402 m_delLumi = 0.;
403 m_intLumi = 0.;
404 m_deltaT = 0.;
405
406 m_totalPrescaleWLiveTime = 0.;
407 m_totalPrescale = 0.;
408 m_lumiWOPrescale = 0.;
409 m_lumiLAr = 0.;
410 m_TotaldeltaT = 0.;
411 m_livefrac = 0.;
412 m_livetime = 0.;
413 m_lartime = 0.;
414 m_larfrac = 0.;
415 m_bsvalid = 0.;
416 m_livetime_l1acc = 0;
417 m_l1acc = 0;
418 m_l2acc = 0;
419 m_l3acc = 0;
420 m_totall1acc = 0;
421 m_livtrig_totall1acc = 0;
422 m_totall2acc = 0;
423 m_totall3acc = 0;
424 m_l1prescale = 1.;
425 m_l2prescale = 1.;
426 m_l3prescale = 1.;
427 m_livetime_beforeprescale = 0;
428 m_livetime_afterprescale = 0;
429 m_afterprescale = 0;
430 m_beforeprescale = 0;
431 m_totall1befpresc = 0;
432 m_runnbr = 0;
433 m_lbstart = 0;
434 m_lbstop = 0;
435 m_lbstart_prev = 0;
436 m_lbstop_prev = 0;
437 m_totalgoodblock = 0;
438 m_totalbadblock = 0;
439 m_clumiblocknbr = 0;
440 m_clumiblocknbrend =0;
441 m_triglevel = 0;
442 m_lbcollectionname = "LumiBlocks";
443 m_triggerchain = triggerchain;
444 m_l1rate = 0.;
445 m_l2rate = 0.;
446 m_l3rate = 0.;
447 m_l1ratediveffxsec = 0.;
448 m_total_l1ratediveffxsec = 0.;
449 m_total_l1ratediveffxsecRun = 0.;
450 m_l1ratediveffxsec_recorded = 0.;
451 m_total_l1ratediveffxsec_recorded = 0.;
452 m_total_l1ratediveffxsecRun_recorded = 0.;
453 m_runnbr_prev = 0;
454
455 // Set to extreme values, to be updated below
456 m_minrun = 99999999;
457 m_maxrun = 0;
458
459 bool firstL1Missing = true;
460 bool firstHLTMissing = true;
461
462 std::set<cool::ValidityKey> lbrunset;
463 lbrunset.clear();
464
465 // collisions xml file settings
466 if (m_collsgrl!=0) { delete m_collsgrl; m_collsgrl=0; }
467 if (m_makecollList) {
468 m_collsgrl = new Root::TGoodRunsList();
469 m_collsgrl->SetVersion("30"); // lumicalc signature
470 m_collsgrl->AddMetaData("Query","Generated by LumiCalculator");
471 }
472
473 // Figure out trigger level
474 m_triglevel = cq_lumi->getTriggerLevel(triggerchain);
475 if (m_triglevel == 0){
476 // Check if on purpose, otherwise, give a warning
477 if (triggerchain == "None") {
478 m_logger << Root::kINFO << "No trigger specified!" << Root::GEndl;
479 } else {
480 m_logger << Root::kWARNING << "Invalid trigger: [" << triggerchain << "] - will proceed with no trigger defined!" << Root::GEndl;
481 }
482 }
483
484 // Look for b-jet triggers
485 regex_t regex;
486 int reti = regcomp(&regex, "_[[:digit:]]?b[[:digit:]]+_", REG_EXTENDED);
487 if (reti) m_logger << Root::kWARNING << "Could not compile regex!" << Root::GEndl;
488
489 reti = regexec(&regex, triggerchain.c_str(), 0, NULL, 0);
490 if ( !reti && !m_usebs) {
491 m_logger << Root::kWARNING << "Trigger: [" << triggerchain << "] appears to be a b-jet trigger, but online beamspot validity not included in livefraction!" << Root::GEndl;
492 m_logger << Root::kWARNING << "Probably need to specify --beamspot to get accurate luminosity!" << Root::GEndl;
493 }
494 regfree(&regex);
495
496 // get lumi channel id (this never changes)
497 if (m_lumimethod != "") {
498 m_Lumiid = cq_lumi->getLumiChannelId(m_lumimethod, m_parlumiestfolder);
499 } else {
500 m_Lumiid = m_lumichannel;
501 }
502
503
504 // do main LumiBlock loop
505 //==============================
506
507 IOVData<cool::Int32> L1preObj;
508 IOVData<cool::Int32> L1preOther;
509 IOVData<cool::Float> L2preObj;
510 IOVData<cool::Float> L3preObj;
513
514 std::map<cool::ValidityKey, CoolQuery::LumiFolderData> LumiDataMap;
515
516 // UTC nanoseconds since 1970
517 std::map<cool::ValidityKey, cool::UInt63> L1starttime_map;
518 std::map<cool::ValidityKey, cool::UInt63> L1endtime_map;
519
520 // Livetime maps
521 std::map<cool::ValidityKey, CoolQuery::L1CountFolderData> Livetime_map;
522 std::map<cool::ValidityKey, CoolQuery::L1CountFolderData> L1accept_map;
523
524 for(xAOD::LumiBlockRangeContainer::const_iterator it = iovc->begin(); it != iovc->end(); ++it){
525 IOVRange iovr(IOVTime((*it)->startRunNumber(),(*it)->startLumiBlockNumber()),
526 IOVTime((*it)->stopRunNumber(),(*it)->stopLumiBlockNumber()));
527
528 // Bookkeeping temporary results
529 m_t_totalDelL = 0.;
530 m_t_totalL = 0.;
531 m_t_totalLRun = 0.;
532 m_t_totaltime = 0.;
533 m_t_deltaT = 0.;
534 m_t_l1acc = 0;
535 m_t_l2acc = 0;
536 m_t_l3acc = 0;
537 m_t_totalgoodblock = 0;
538 m_t_totalbadblock = 0;
539 m_t_totall1befpresc = 0;
540 m_t_totalPrescaleWLiveTime = 0.;
541 m_t_totalPrescale = 0.;
542 m_t_lumiWOPrescale = 0. ;
543 m_t_lumiLAr = 0.;
544
545 m_runnbr = iovr.start().run();
546 m_lbstart = iovr.start().event();
547 m_lbstop = iovr.stop().event();
548
549 // Look for duplicate run/LB
550 if(m_lbstart_prev == m_lbstart && m_lbstop_prev == m_lbstop && m_runnbr_prev == m_runnbr){
551 m_logger << Root::kWARNING << "Skipping multiply stored IOVRange: [" << m_lbstart << "-" << m_lbstop << "]" << Root::GEndl;
552 continue;
553 }
554
555 // new run, reset its lumi and reload the trigger channels
556 if ( m_runnbr!=m_runnbr_prev ) {
557
558 if (m_runnbr < m_minrun) m_minrun = m_runnbr;
559 if (m_runnbr > m_maxrun) m_maxrun = m_runnbr;
560
561 m_totalLRun=0.;
562 m_total_l1ratediveffxsecRun=0.;
563 m_total_l1ratediveffxsecRun_recorded=0.;
564
565 // Set IOV range for full run
566 cq_trigger->setIOVForRun (m_runnbr);
567 cq_lumi->setIOVForRun(m_runnbr);
568
569 // get trigger channel ids
570 m_L1Valid = false;
571 m_L2Valid = false;
572 m_L3Valid = false;
573 m_LiveValid = false;
574 m_triggerlowerchains.clear();
575 // These are used to resolve L1 ORs like [L1_MU20,L1_MU21]
576 m_L1triggerchains.clear();
577 m_L1idList.clear();
578
579 std::string lowerch = "";
580
581 // This is really inefficient...
582 if(m_triglevel == 3){
583 m_L3id = cq_trigger->getHLTChannelId(triggerchain, m_parhltmenufolder);
584 m_L3Valid = cq_trigger->channelIdValid();
585
586 lowerch = cq_trigger->getHLTLowerChainName(triggerchain, m_parhltmenufolder);
587 m_triggerlowerchains.push_back(lowerch);
588 //
589 m_L2id = cq_trigger->getHLTChannelId(lowerch, m_parhltmenufolder);
590 if (m_L2id == UINT_MAX) { //nonsense value for getHLTChannelId
591 throw std::runtime_error("LumiCalculator::IntegrateLumi : getHLTChannelId returned invalid value.");
592 }
593 m_L2Valid = cq_trigger->channelIdValid();
594
595 lowerch = cq_trigger->getHLTLowerChainName(lowerch, m_parhltmenufolder);
596 m_triggerlowerchains.push_back(lowerch);
597 //
598 m_L1id = cq_trigger->getL1ChannelId(lowerch, m_parlvl1menufolder );
599 m_L1Valid = cq_trigger->channelIdValid();
600
601 }else if(m_triglevel == 2){ // Run2 or starting at L2 in Run1
602 lowerch = cq_trigger->getHLTLowerChainName(triggerchain, m_parhltmenufolder);
603 m_triggerlowerchains.push_back(lowerch);
604 m_L2id = cq_trigger->getHLTChannelId(triggerchain, m_parhltmenufolder);
605 m_L2Valid = cq_trigger->channelIdValid();
606
607 ParseL1Trigger(lowerch, cq_trigger);
608
609 }else if(m_triglevel == 1){
610
611 ParseL1Trigger(triggerchain, cq_trigger);
612 }
613
614 if (firstL1Missing && (!m_L1Valid && m_triglevel > 0)) {
615 firstL1Missing = false;
616 cq_trigger->printL1Triggers(m_parlvl1menufolder);
617 }
618
619 if (firstHLTMissing && ((!m_L2Valid && m_triglevel > 1) || (!m_L3Valid && m_triglevel > 2)) ) {
620 firstHLTMissing = false;
621 cq_trigger->printHLTTriggers(m_parhltmenufolder);
622 }
623
624 // Should we check here for btag trigger and onl beamspot?
625
626 // Do we use dedicated livetime trigger?
627 if(m_uselivetrigger){
628 m_LiveL1id = cq_trigger->getL1ChannelId(m_livetrigger, m_parlvl1menufolder);
629 m_LiveValid = cq_trigger->channelIdValid();
630 }else{
631 // then fall back to the original trigger
632 m_LiveL1id = m_L1id;
633 m_LiveValid = m_L1Valid;
634 }
635
636 if (!m_LiveValid) {
637 if (m_uselivetrigger)
638 m_logger << Root::kWARNING << "Runnumber: [" << m_runnbr << "] can't find livefraction trigger [" << m_livetrigger << "]! Livefraction won't be calculated!" << Root::GEndl;
639 else
640 m_logger << Root::kWARNING << "Runnumber: [" << m_runnbr << "] can't find trigger [" << triggerchain << "]! Livefraction won't be calculated!" << Root::GEndl;
641 }
642
643 // Also load all deadtime counters here
644
645 // Need to fix this for non-MC. HLT counters now available, but not used here...
646
647 //-------------------------
648 // Load livetime information
649 Livetime_map.clear();
650 L1accept_map.clear();
651
652 if (m_LiveValid)
653 Livetime_map = cq_trigger->getL1CountFolderData(m_parlumilvl1folder, m_LiveL1id);
654
655
656 if (m_L1Valid && m_triglevel >= 1)
657 L1accept_map = cq_trigger->getL1CountFolderData(m_parlumilvl1folder, m_L1id);
658
659 // UTC nanoseconds since 1970
660 L1starttime_map.clear();
661 L1endtime_map.clear();
662 L1starttime_map = cq_trigger->getObjMapFromFolderAtChan<cool::UInt63>("StartTime", m_parlvl1lblbfolder, 0);
663 L1endtime_map = cq_trigger->getObjMapFromFolderAtChan<cool::UInt63>("EndTime", m_parlvl1lblbfolder, 0);
664
665 //---------------------------
666 // Load LAr defects
667 LArObj.clear();
668
669 if (m_uselar) {
670 cool::ValidityKey runstarttime = L1starttime_map.begin()->second;
671 cool::ValidityKey runendtime = L1endtime_map.rbegin()->second;
672
673 cq_lar->setIOV(runstarttime, runendtime);
674 LArObj = cq_lar->getIOVData<cool::UInt32>("EventVeto", m_parlareventvetofolder, 0, m_lartag);
675 }
676
677 //---------------------------
678 // Load Onl beamspot
679 BSObj.clear();
680
681 if (m_usebs) {
682 cq_bs->setIOVForRun(m_runnbr);
683 BSObj = cq_bs->getIOVData<cool::Int32>("status", m_paronlbeamspotfolder, 0, m_bstag);
684 }
685
686
687 //-----------------------------
688 // Load luminosity for this run
689 LumiDataMap.clear();
690 LumiDataMap = cq_lumi->getLumiFolderData(m_parlumiestfolder, m_lumitag, m_Lumiid);
691
692
693 } // End of new run stuff
694
695 m_lbstart_prev = m_lbstart;
696 m_lbstop_prev = m_lbstop;
697 m_runnbr_prev = m_runnbr;
698
699 // Update DB for this specific IOV range
700 cq_trigger->setIOV(iovr.start().re_time(), iovr.stop().re_time());
701
702 // Print this here (will be output for each contiguous LB range in XML file)
703 m_logger << Root::kINFO << std::left << "-----------------------------------" << Root::GEndl;
704 m_logger << Root::kINFO << "Beginning calculation for ";
705 if (m_triglevel > 0) {
706 m_logger << Root::kINFO << "Trigger " << m_triggerchain;
707 if(m_triggerlowerchains.size() > 0){
708 m_logger << Root::kINFO << "[";
709 for(std::vector<std::string>::iterator sit = m_triggerlowerchains.begin(); sit != m_triggerlowerchains.end(); ++sit){
710 m_logger << Root::kINFO << "==>" << *sit;
711 }
712 m_logger << Root::kINFO << "], ";
713 }
714 }
715 m_logger << Root::kINFO << "Run " << m_runnbr << " LB [" << m_lbstart << "-" << m_lbstop << "]" << Root::GEndl;
716 m_logger << Root::kINFO << std::left << "-----------------------------------" << Root::GEndl;
717
718 //----------------------------
719 // Load prescales for this IOV
720 L1preObj.clear();
721 L2preObj.clear();
722 L3preObj.clear();
723 L1preOther.clear(); // Used for mutliple L1 triggers
724
725 if(m_L1Valid) {
726 // Normal case
727 if (m_L1idList.size() <= 1) {
728 L1preObj = cq_trigger->getIOVData<cool::Int32>("Lvl1Prescale", m_parlvl1prescalesfolder, m_L1id);
729
730 // Multiple L1 items defined
731 } else {
732
733 m_logger << Root::kINFO << "Resolving multiple L1 prescales:" << Root::GEndl;
734
735 // Load first one
736 L1preObj = cq_trigger->getIOVData<cool::Int32>("Lvl1Prescale", m_parlvl1prescalesfolder, m_L1idList[0]);
737
738 // Dump prescales
739 std::list< std::pair<IOVRange, cool::Int32> >::iterator it;
740
741 m_logger << Root::kINFO << std::setw(10) << std::left << m_L1triggerchains[0];
742 for(it = L1preObj.data.begin(); it != L1preObj.data.end(); ++it) {
743 m_logger << Root::kINFO << std::setw(1) << std::left << "[" << it->first.start().event() << "," << it->first.stop().event()-1 << "]:" ;
744
745 if (it->second < 0)
746 m_logger << it->second << ", ";
747 else
748 m_logger << (0xFFFFFF / float(0x1000000 - it->second)) << ", ";
749 }
750 m_logger << Root::kINFO << Root::GEndl;
751
752 // Now loop over the remaining and replace
753 for (unsigned int iid=1; iid < m_L1idList.size(); iid++) {
754 L1preOther = cq_trigger->getIOVData<cool::Int32>("Lvl1Prescale", m_parlvl1prescalesfolder, m_L1idList[iid]);
755
756 // Dump prescales
757 m_logger << Root::kINFO << std::setw(10) << std::left << m_L1triggerchains[iid];
758 for(it = L1preOther.data.begin(); it != L1preOther.data.end(); ++it) {
759 m_logger << Root::kINFO << std::setw(1) << std::left << "[" << it->first.start().event() << "," << it->first.stop().event()-1 << "]:" ;
760
761 if (it->second < 0)
762 m_logger << it->second << ", ";
763 else
764 m_logger << (0xFFFFFF / float(0x1000000 - it->second)) << ", ";
765 }
766 m_logger << Root::kINFO << Root::GEndl;
767
768 // Iterate through both lists and keep lowest non-negative prescale
769 std::list< std::pair<IOVRange, cool::Int32> >::iterator it1;
770 std::list< std::pair<IOVRange, cool::Int32> >::iterator it2;
771 for(it1 = L1preObj.data.begin(), it2 = L1preOther.data.begin(); it1 != L1preObj.data.end(); ++it1, ++it2) {
772
773 // -1 is disabled, otherwise this is the event count to prescale
774 if ((it2->second > 0) && (it1->second > it2->second)) {
775 // Check if they are both prescaled (this is bad)
776 if ((it1->second > 1) && (it2->second > 1)) {
777 m_logger << Root::kWARNING << "L1 Prescales combined with both triggers prescaled for Run " << m_runnbr << "!" << Root::GEndl;
778 }
779 it1->second = it2->second;// See if this works...
780 }
781 }
782 }
783
784 // OK done, lets check the result
785 // Dump prescales
786 m_logger << Root::kINFO << std::setw(10) << std::left << "L1 Pre:";
787 for(it = L1preObj.data.begin(); it != L1preObj.data.end(); ++it) {
788 m_logger << Root::kINFO << std::setw(1) << std::left << "[" << it->first.start().event() << "," << it->first.stop().event()-1 << "]:" ;
789
790 if (it->second < 0)
791 m_logger << it->second << ", ";
792 else
793 m_logger << (0xFFFFFF / float(0x1000000 - it->second)) << ", ";
794 }
795 m_logger << Root::kINFO << Root::GEndl;
796
797 }
798
799 }
800
801 if(m_L2Valid && m_L2id != UINT_MAX) {
802 if (isrun2) {
803 L2preObj = cq_trigger->getIOVData<cool::Float>("Prescale", m_parhltprescalesfolder, 20000 + m_L2id);
804 }
805 else {
806 L2preObj = cq_trigger->getIOVData<cool::Float>("Prescale", m_parhltprescalesfolder, 2*m_L2id);
807 }
808 }
809
810 if(m_L3Valid) {
811 L3preObj = cq_trigger->getIOVData<cool::Float>("Prescale", m_parhltprescalesfolder, 2*m_L3id+1);
812 }
813
814 // Reload the time map to get the ATLAS range
815 L1starttime_map.clear();
816 L1endtime_map.clear();
817 L1starttime_map = cq_trigger->getObjMapFromFolderAtChan<cool::UInt63>("StartTime", m_parlvl1lblbfolder, 0);
818 L1endtime_map = cq_trigger->getObjMapFromFolderAtChan<cool::UInt63>("EndTime", m_parlvl1lblbfolder, 0);
819
820 // Restrict lb range if necessary based on actual ATLAS run/lb values
821 if (L1starttime_map.begin()->first > iovr.start().re_time() || L1starttime_map.rbegin()->first < iovr.stop().re_time()) {
822 m_lbstart = (L1starttime_map.begin()->first & 0xFFFFFFFF);
823 m_lbstop = (L1starttime_map.rbegin()->first & 0xFFFFFFFF);
824 m_logger << Root::kINFO << "Restricting to valid ATLAS lumi block range [" << m_lbstart << "-" << m_lbstop << "]" << Root::GEndl;
825 }
826
827
828 //
829 // Dont assume that all lumi blocks in this range are present in the lumi DB.
830 // Loop over all lumi blocks specifically and flag any missing lumi entries as bad lumi blocks.
831 //
832
833 // Counters to check for missing LB ranges
834 int firstMissing = -1;
835 int lastMissing = -1;
836
837 for (cool::ValidityKey currentVK = L1starttime_map.begin()->first; currentVK <= L1starttime_map.rbegin()->first; currentVK++) {
838
839 // Current IOVTime
840 IOVTime curIOV;
841 curIOV.setRETime(currentVK);
842
843 m_clumiblocknbr = curIOV.event();
844 m_clumiblocknbrend = curIOV.event()+1;
845
846 // Check here for duplicate lumi blocks by explicitly keeping track of every run/lb seen
847 if (lbrunset.count(curIOV.re_time()) != 0) {
848 m_logger << Root::kWARNING << "Skipping duplicate [run,lumiblock]: " << curIOV << " !" << Root::GEndl;
849 continue;
850 } else {
851 lbrunset.insert(curIOV.re_time());
852 }
853
854 // Not really good, just all seen
855 m_totalgoodblock += 1;
856 m_t_totalgoodblock += 1;
857
858 // Find luminosity record
859 std::map<cool::ValidityKey, CoolQuery::LumiFolderData>::iterator itOL = LumiDataMap.find(currentVK);
860
861 // Check if we didn't find anything
862 if (itOL == LumiDataMap.end()) {
863
864 m_t_totalbadblock++;
865 m_totalbadblock++;
866
867 if (firstMissing < 0) {
868 // Is this the first missing one? If so, make a note of it and go on
869 firstMissing = curIOV.event();
870 lastMissing = firstMissing;
871 } else if (int(curIOV.event()) == (lastMissing+1)) {
872 // Is this a contiguous missing lumi block?
873 lastMissing = curIOV.event();
874 } else {
875 // Not contiguous, print previous
876 if (firstMissing == lastMissing) {
877 m_logger << Root::kWARNING << "Luminosity info not found for Run " << m_runnbr << " LB " << firstMissing << " !" << Root::GEndl;
878 } else {
879 m_logger << Root::kWARNING << "Luminosity info not found for Run " << m_runnbr << " LB [" << firstMissing << "-" << lastMissing << "] !" << Root::GEndl;
880 }
881 firstMissing = curIOV.event();
882 lastMissing = firstMissing;
883 }
884
885 // If last time through loop, print this also
886 if (currentVK == L1starttime_map.rbegin()->first) {
887 if (firstMissing == lastMissing) {
888 m_logger << Root::kWARNING << "Luminosity info not found for Run " << m_runnbr << " LB " << firstMissing << " !" << Root::GEndl;
889 } else {
890 m_logger << Root::kWARNING << "Luminosity info not found for Run " << m_runnbr << " LB [" << firstMissing << "-" << lastMissing << "] !" << Root::GEndl;
891 }
892 firstMissing = -1;
893 lastMissing = -1;
894 }
895
896 // Skip rest of this LB
897 continue;
898
899 } else {
900
901 // Check if we had previous missing block
902 if (firstMissing >= 0) {
903 if (firstMissing == lastMissing) {
904 m_logger << Root::kWARNING << "Luminosity info not found for Run " << m_runnbr << " LB " << firstMissing << " !" << Root::GEndl;
905 } else {
906 m_logger << Root::kWARNING << "Luminosity info not found for Run " << m_runnbr << " LB [" << firstMissing << "-" << lastMissing << "] !" << Root::GEndl;
907 }
908
909 firstMissing = -1;
910 lastMissing = -1;
911 }
912
913 }
914
915
916 // given in units of ub^{-1} = 10^{30} cm^{-2}
917 m_instLumi = (itOL->second).LBAvInstLumi;
918 m_AvEvtsPerBX = (itOL->second).LBAvEvtsPerBX;
919 m_Valid = (itOL->second).Valid;
920
921
922 // Clear values in case trigger isn't defined
923 m_l1acc = 0.;
924 m_beforeprescale = 0.;
925 m_afterprescale = 0.;
926 m_l2acc = 0.;
927 m_l3acc = 0.;
928
929 // Store dummy prescale values at start
930 m_l1prescale = -1.;
931 m_l2prescale = -1.;
932 m_l3prescale = -1.;
933
934 // Some trigger is defined. Get prescales and values here
935 if(m_L1Valid && m_triglevel > 0) {
936
937 // Get L1 prescale
938 if (isrun2) {
939 if (L1preObj.getValue(curIOV) < 0)
940 m_l1prescale = L1preObj.getValue(curIOV);
941 else
942 m_l1prescale = 0xFFFFFF / float(0x1000000 - L1preObj.getValue(curIOV));
943 }
944 else {
945 m_l1prescale = L1preObj.getValue(curIOV);
946 }
947
948 if (m_triglevel >=2) {
949 if(m_L2Valid) {
950
951 // Get L2 prescale
952 m_l2prescale = L2preObj.getValue(curIOV);
953
954 }
955 // Else, prescale stays at -1.
956
957 } else {
958 // Force L2 'passthrough'
959 m_l2prescale = 1.;
960 }
961
962 if(m_triglevel == 3){
963 if (m_L3Valid) {
964
965 // Get L3 prescale
966 m_l3prescale = L3preObj.getValue(curIOV);
967
968 }
969 // Else, L3 prescale stays at -1.
970 } else {
971 // Force L3 'passthrough'
972 m_l3prescale = 1.;
973 }
974 }
975
976 //-------------------------------
977 // Calculate livetime from a dedicated not rare trigger if user requested
978 auto livetime_it = Livetime_map.find(currentVK);
979 if (livetime_it == Livetime_map.end()) {
980 throw std::runtime_error("LumiCalculator::IntegrateLumi Start or End times not found in map.");
981 }
982 CoolQuery::L1CountFolderData l1count = livetime_it->second;
983
984 m_livetime_beforeprescale = l1count.BeforePrescale;
985 m_livetime_afterprescale = l1count.AfterPrescale;
986 m_livetime_l1acc = l1count.L1Accept;
987 if(m_livetime_afterprescale > 0.){
988 m_livefrac = m_livetime_l1acc/(float)m_livetime_afterprescale;
989 }else{
990 m_livefrac = 0.0;
991 }
992
993 // Hack for bad SFO in particular run
994 if (m_runnbr == 286367) {
995 m_livefrac *= 5./6.;
996 }
997
998 if (m_runnbr == 281385) {
999 if (m_clumiblocknbr <= 196) {
1000 m_livefrac *= 4./6.;
1001 } else if (m_clumiblocknbr <= 374) {
1002 m_livefrac *= 5./6.;
1003 }
1004 }
1005
1006 // Check for low statistics in afterprescale counts
1007 if(m_livetime_beforeprescale > 0 && m_livetime_afterprescale <= 0 ){
1008 std::string ttrig = "";
1009 ttrig = triggerchain;
1010 if(m_uselivetrigger)ttrig = m_livetrigger;
1011 m_logger << Root::kWARNING << "L1 counts after prescale (before veto) are 0.0 for trigger " << std::move(ttrig) << "! Livefraction set to zero!" << Root::GEndl;
1012 m_logger << Root::kWARNING << "Try using a high rate L1 trigger for livetime calculation: --livetrigger=<high rate L1 trigger> " << Root::GEndl;
1013 m_logger << Root::kINFO << m_runnbr << "[" << m_clumiblocknbr << "]: L1Acc: " << m_l1acc << ", AfterPrescale: " << m_afterprescale << ", L1Presc: " << m_l1prescale << Root::GEndl;
1014 }
1015
1016 //------------------------
1017 // Calculate LAr veto time
1018 auto itStartTime = L1starttime_map.find(currentVK);
1019 auto itEndTime = L1endtime_map.find(currentVK);
1020 if (itStartTime == L1starttime_map.end() or itEndTime == L1endtime_map.end()){
1021 throw std::runtime_error("LumiCalculator::IntegrateLumi Start or End times not found in map.");
1022 }
1023 cool::ValidityKey lbstarttime = L1starttime_map.find(currentVK)->second;
1024 cool::ValidityKey lbendtime = L1endtime_map.find(currentVK)->second;
1025
1026 m_lartime = 0.; // Time to exclude in seconds
1027 if (m_uselar) {
1028
1030 lbstart.setRETime(lbstarttime);
1031 lbend.setRETime(lbendtime);
1032 IOVRange range(lbstart, lbend);
1033
1034 std::list<std::pair<IOVRange, cool::UInt32> > larlist;
1035 larlist = LArObj.getOverlap(range);
1036
1037 for (std::list<std::pair<IOVRange, cool::UInt32> >::iterator it = larlist.begin(); it != larlist.end(); ++it) {
1038 if (it->second == 0) continue;
1039 float dtime = (it->first.stop().re_time() - it->first.start().re_time())/1.E9;
1040 m_lartime += dtime;
1041 if (m_verbose == true) {
1042 m_logger << Root::kINFO << "Found LAr veto from " << it->first.start().re_time() << " to " << it->first.stop().re_time() << " = " << dtime << " seconds" << Root::GEndl;
1043 }
1044 }
1045 }
1046
1047 //------------------------
1048 // Check for online BS
1049 m_bsvalid = 1.;
1050 if (m_usebs) {
1051
1052 // Read beamspot validity
1053 bool valid = true;
1054 int status = BSObj.getValue(curIOV);
1055 //m_logger << Root::kINFO << "Found online beamspot status = " << status << Root::GEndl;
1056
1057 if (status != 7) valid = false;
1058
1059 if (!valid) {
1060 m_bsvalid = 0.0;
1061 m_livefrac = 0.0;
1062 if(m_verbose == true){
1063 m_logger << Root::kINFO << m_runnbr << "[" << m_clumiblocknbr << "]: Online beamspot invalid with Lumi=" << m_instLumi << " 10^30 cm-2 s-1" << Root::GEndl;
1064 }
1065 }
1066 }
1067
1068 //-------------------------------
1069
1070 {
1071 auto it = L1accept_map.find(currentVK);
1072 if (it == L1accept_map.end()) {
1073 throw std::runtime_error("LumiCalculator::IntegrateLumi Start or End times not found in map.");
1074 }
1075 l1count = it->second;
1076 }
1077 m_beforeprescale = l1count.BeforePrescale;
1078 m_beforeprescaleof = false;
1079 m_afterprescale = l1count.AfterPrescale;
1080 m_afterprescaleof = false;
1081 m_l1acc = l1count.L1Accept;
1082 m_l1accof = false;
1083
1084 m_deltaT = (lbendtime-lbstarttime)/1.E9;
1085 m_lbstarttime = lbstarttime/1.E9;
1086 m_lbendtime = lbendtime/1.E9;
1087
1088 if (m_deltaT > 0.) m_larfrac = 1.-m_lartime/m_deltaT;
1089
1090 // For online lumi database case one needs to make some changes:
1091 if(m_onlinelumi == true){
1092 // In Valid UInt32 type value some information is encoded:
1093 // see: https://twiki.cern.ch/twiki/bin/view/Atlas/CoolOnlineData#Folder_TRIGGER_LUMI_LBLESTONL
1094 cool::UInt32 tempValid = (m_Valid & 1023);
1095
1096 if(tempValid == 10){
1097 // in this bits, Value 301 and 302 means MBTS and LUCID lumi for which some care is needed
1098 if((m_Valid >> 22) == 301
1099 || (m_Valid >> 22) == 302
1100 || (m_Valid >> 22) == 101
1101 || (m_Valid >> 22) == 102
1102 || (m_Valid >> 22) == 103
1103 || (m_Valid >> 22) == 104
1104 ){
1105 m_Valid = 0;
1106 }else{
1107 m_Valid = tempValid;
1108 }
1109 }else{
1110 m_Valid = tempValid;
1111 }
1112
1113 } else {
1114
1115 // For offline, we also need to strip out the preferred channel value from the validity word
1116 m_Valid &= 0xFFFF;
1117 }
1118
1119 // Dump out debugging information
1120 if(m_verbose == true){
1121 m_logger << Root::kINFO << m_runnbr << "[" << m_clumiblocknbr << "]: L1Acc: " << m_l1acc;
1122 if(m_uselivetrigger) m_logger << ", Livetime trigger L1Acc: " << m_livetime_l1acc;
1123 m_logger << ", InstLumi: " << m_instLumi << ", deltaT: " << m_deltaT << ", AvEvtsPerBX: " << m_AvEvtsPerBX << ", BeforePrescale: " << m_beforeprescale << ", AfterPrescale: " << m_afterprescale;
1124 if (m_uselivetrigger) m_logger << ", Livetime trigger BeforePrescale: " << m_livetime_beforeprescale << " Livetime trigger AfterPrescale: " << m_livetime_afterprescale;
1125 if (isrun2) {
1126 m_logger << ", Livefrac: " << m_livefrac << ", L1Presc: " << m_l1prescale << ", HLTPresc: " << m_l2prescale << ", Valid: " << m_Valid;
1127 } else {
1128 m_logger << ", Livefrac: " << m_livefrac << ", L1Presc: " << m_l1prescale << ", L2Presc: " << m_l2prescale << ", L3Presc: " << m_l3prescale << ", Valid: " << m_Valid;
1129 }
1130 if (m_uselar) m_logger << ", LAr ready fraction: " << m_larfrac;
1131 m_logger << Root::GEndl;
1132 }
1133
1134 // Check if we have valid luminosity
1135 // Just need to check lowest digit. 10 -> BCID blind, 100 -> mu not valid
1136 if(m_Valid%10 != 0){
1137
1138 // Invalid luminosity entry, call it bad
1139 m_instLumi = 0.0;
1140 m_totalbadblock += 1;
1141 m_t_totalbadblock += 1;
1142 m_logger << Root::kWARNING << "Skipping lumiblock " << m_runnbr << "[" << m_clumiblocknbr << "] with invalid inst. lumi. (valid=" << m_Valid << ")!" << Root::GEndl;
1143
1144 } else if ((m_triglevel > 0) && (m_l1prescale < 0. || m_l2prescale < 0. || m_l3prescale < 0.)) {
1145
1146 // Disabled trigger, call bad but still record delivered luminosity
1147 m_totalbadblock += 1;
1148 m_t_totalbadblock += 1;
1149 m_logger << Root::kWARNING << "Lumiblock " << m_runnbr << "[" << m_clumiblocknbr << "] has a disabled or incorrectly specified trigger.! " << Root::GEndl;
1150
1151 }
1152
1153 //========================================================//
1154 // L U M I N O S I T Y C A L C U L A T I O N H E R E //
1155 //========================================================//
1156 // The actual calculation of integrated luminosity and
1157 // accumulation of some variables
1158 m_totalDelL += (m_deltaT*m_instLumi); // delivered lumi
1159 m_t_totalDelL += (m_deltaT*m_instLumi); // delivered lumi
1160 m_t_deltaT += m_deltaT;
1161
1162 // Count up everything
1163 m_totall1acc += m_l1acc;
1164 m_livtrig_totall1acc += m_livetime_l1acc;
1165 m_t_l1acc += m_l1acc;
1166 m_totall1befpresc += m_beforeprescale;
1167 m_t_totall1befpresc += m_beforeprescale;
1168 m_totall2acc += m_l2acc;
1169 m_t_l2acc += m_l2acc;
1170 m_totall3acc += m_l3acc;
1171 m_t_l3acc += m_l3acc;
1172 m_livetime = m_livefrac*m_deltaT;
1173 m_totaltime += m_livetime;
1174 m_t_totaltime += m_livetime;
1175
1176 double totalPrescale = m_l1prescale * m_l2prescale * m_l3prescale;
1177
1178 // Check for disabled triggers
1179 if (m_l1prescale < 0. ) totalPrescale = 0.;
1180 if (m_l2prescale < 0. ) totalPrescale = 0.;
1181 if (m_l3prescale < 0. ) totalPrescale = 0.;
1182
1183 // Check also for no trigger
1184 if (m_triglevel == 0) totalPrescale = 1.;
1185
1186 m_intLumi = 0.;
1187
1188 m_lumiWOPrescale += m_livetime*m_instLumi ;
1189 m_t_lumiWOPrescale += m_livetime*m_instLumi ;
1190
1191 m_lumiLAr += m_livetime*m_larfrac*m_instLumi;
1192 m_t_lumiLAr += m_livetime*m_larfrac*m_instLumi;
1193
1194 if (totalPrescale > 0.) {
1195 m_totalPrescaleWLiveTime += m_livetime/totalPrescale;
1196 m_t_totalPrescaleWLiveTime += m_livetime/totalPrescale;
1197 m_totalPrescale += 1./totalPrescale;
1198 m_t_totalPrescale += 1./totalPrescale;
1199 m_intLumi = m_larfrac * m_livetime * m_instLumi/totalPrescale; // <<<--- T H E F O R M U L A
1200 }
1201
1202
1203 m_totalL += m_intLumi;
1204 m_t_totalL += m_intLumi;
1205 m_totalLRun += m_intLumi;
1206
1207
1208 // MB: trigger rates, note that livefrac drops out of ratio
1209 m_l1rate = ( m_livetime>0 ? m_l1acc / m_livetime : 0. );
1210 m_l2rate = ( m_livetime>0 ? m_l2acc / m_livetime : 0. );
1211 m_l3rate = ( m_livetime>0 ? m_l3acc / m_livetime : 0. );
1212
1213 // MB: delivered lumi
1214 m_l1ratediveffxsec = (float)m_afterprescale/( m_deltaT*m_effxsec );
1215 m_total_l1ratediveffxsec += (float)m_afterprescale / m_effxsec ;
1216 m_total_l1ratediveffxsecRun += (float)m_afterprescale / m_effxsec ;
1217
1218 // MB: recorded lumi
1219 m_l1ratediveffxsec_recorded = (float)m_l1acc /( m_deltaT*m_effxsec );
1220 m_total_l1ratediveffxsec_recorded += (float)m_l1acc / m_effxsec ;
1221 m_total_l1ratediveffxsecRun_recorded += (float)m_l1acc / m_effxsec ;
1222
1223 if (m_collsgrl!=0) {
1224 if ( m_l1rate>=m_mintrigrate )
1225 m_collsgrl->AddRunLumiBlock(m_runnbr,m_clumiblocknbr);
1226 }
1227
1228 if (m_verbose) {
1229 if (m_effxsec!=1.0)
1230 m_logger << Root::kINFO << "L1rate a/ prescale: " << m_afterprescale << ", Delivered LumiFromL1rate (/ub): " << m_l1ratediveffxsec << ", Delivered TotalLumiFromL1rate (/ub): " << m_total_l1ratediveffxsec
1231 << Root::GEndl;
1232 }
1233
1234 if(m_LumiTree != 0)m_LumiTree->Fill();
1235 } // End of loop over lumi blocks
1236
1237 // Print IOV summary
1238 m_logger << Root::kINFO<< std::setw(10) << std::left << ">== Trigger : " << triggerchain << Root::GEndl;
1239 m_logger << Root::kINFO<< std::setw(10) << std::right << "Run" << std::setw(10) << std::right << "L1-Acc" << std::setw(10) << std::right << "L2-Acc" << std::setw(10) << std::right << "L3-Acc" << std::setw(10) << std::right << "LiveTime" << std::setw(18) << std::right << "IntL rec.(ub^-1)" << std::setw(18) << std::right << "IntL del.(ub^-1)" << Root::GEndl;
1240 m_logger << Root::kINFO<< std::setw(10) << std::right << m_runnbr << std::setw(10) << std::right << m_t_l1acc << std::setw(10) << std::right << m_t_l2acc << std::setw(10) << std::right << m_t_l3acc << std::setw(10) << std::right << m_t_totaltime << std::setw(18) << std::right << m_t_totalL << std::setw(18) << std::right << m_t_totalDelL << Root::GEndl;
1241 // m_logger << Root::kINFO << std::setw(10) << std::right << "BeforePrescale" << std::setw(10) << std::right << m_totall1befpresc << std::setw(10) << std::right << "" << std::setw(10) << std::right << "" << std::setw(10) << std::right << "" << std::setw(10) << std::right << "" << std::setw(10) << std::right << "" << std::setw(10) << m_TotaldeltaT << std::setw(14) << std::right << m_totalL << Root::GEndl;
1242
1243 m_logger << Root::kINFO<< std::setw(10) << std::left << "L1/2/3 accept: " << std::setw(10) << std::left << m_t_l1acc << std::setw(10) << std::left << m_t_l2acc << std::setw(10) << std::left << m_t_l3acc << Root::GEndl;
1244 m_logger << Root::kINFO << std::setw(10) << std::left << "L1BeforePresc: " << std::setw(10) << std::left << m_t_totall1befpresc << Root::GEndl;
1245 m_logger << Root::kINFO << std::setw(10) << std::left << "Livetime : " << m_t_totaltime << Root::GEndl;
1246 m_logger << Root::kINFO << std::setw(10) << std::left << "Prescale Weighted Livetime: " << m_t_totalPrescaleWLiveTime << Root::GEndl;
1247 m_logger << Root::kINFO<< std::setw(10) << std::left << "Good LBs : " << m_t_totalgoodblock - m_t_totalbadblock << Root::GEndl;
1248 m_logger << Root::kINFO<< std::setw(10) << std::left << "Bad LBs : " << m_t_totalbadblock << Root::GEndl;
1249
1250 if ( m_effxsec==1.0 ) {
1251 m_logger << Root::kINFO << std::setw(10) << std::left << "IntL delivered (ub^-1) : " << m_t_totalDelL << Root::GEndl;
1252 m_logger << Root::kINFO << std::setw(10) << std::left << "IntL after livefraction (ub^-1): " << m_t_lumiWOPrescale << Root::GEndl;
1253 if (m_uselar)
1254 m_logger << Root::kINFO << std::setw(10) << std::left << "IntL after LAr fraction (ub^-1): " << m_t_lumiLAr << Root::GEndl;
1255 m_logger << Root::kINFO << std::setw(10) << std::left << "IntL recorded after prescale (ub^-1) : " << m_t_totalL << Root::GEndl;
1256 } else {
1257 m_logger << Root::kINFO << std::setw(10) << std::left << "IntL delived (ub^-1) : " << m_total_l1ratediveffxsec << Root::GEndl;
1258 m_logger << Root::kINFO << std::setw(10) << std::left << "IntL recorded (ub^-1) : " << m_total_l1ratediveffxsec_recorded << Root::GEndl;
1259 }
1260
1261 // Print prescales as range of actual lumi blocks these apply to
1262 if(m_triglevel >= 1){
1263 // Print L1 Prescale values:
1264 m_logger << Root::kINFO << std::setw(10) << std::left << "L1 Prescales: ";
1265
1266 std::list< std::pair<IOVRange, cool::Int32> >::iterator it;
1267 for(it = L1preObj.data.begin(); it != L1preObj.data.end(); ++it) {
1268 m_logger << Root::kINFO << std::setw(1) << std::left << "[" << it->first.start().event() << "," << it->first.stop().event()-1 << "]:" ;
1269 if (isrun2) {
1270 if (it->second < 0)
1271 m_logger << it->second << ", ";
1272 else
1273 m_logger << (0xFFFFFF / float(0x1000000 - it->second)) << ", ";
1274 } else {
1275 m_logger << it->second << ", ";
1276 }
1277 }
1278 m_logger << Root::kINFO << Root::GEndl;
1279 }
1280
1281 if(m_triglevel >= 2){
1282 // Print L2 Prescale values:
1283 if (isrun2) {
1284 m_logger << Root::kINFO << std::setw(10) << std::left << "HLT Prescales: ";
1285 }
1286 else {
1287 m_logger << Root::kINFO << std::setw(10) << std::left << "L2 Prescales: ";
1288 }
1289
1290 std::list< std::pair<IOVRange, cool::Float> >::iterator it;
1291 for(it = L2preObj.data.begin(); it != L2preObj.data.end(); ++it) {
1292 m_logger << Root::kINFO << std::setw(1) << std::left << "[" << it->first.start().event() << "," << it->first.stop().event()-1 << "]:" << it->second << ", ";
1293 }
1294 m_logger << Root::kINFO << Root::GEndl;
1295 }
1296
1297 if(m_triglevel == 3){
1298 // Print L3 Prescale values:
1299 m_logger << Root::kINFO << std::setw(10) << std::left << "L3 Prescales: ";
1300
1301 std::list< std::pair<IOVRange, cool::Float> >::iterator it;
1302 for(it = L3preObj.data.begin(); it != L3preObj.data.end(); ++it) {
1303 m_logger << Root::kINFO << std::setw(1) << std::left << "[" << it->first.start().event() << "," << it->first.stop().event()-1 << "]:" << it->second << ", ";
1304 }
1305 m_logger << Root::kINFO << Root::GEndl;
1306 }
1307
1308 } // end lb collection loop
1309
1310
1311 // ------------------------------------------------------------------------------------------------
1312 // MB : Print total only at end of LB loop:
1313 m_logger << Root::kINFO << std::left << "-----------------------------------" << Root::GEndl;
1314 m_logger << Root::kINFO << std::left << " LumiCalculator summary" << Root::GEndl;
1315 m_logger << Root::kINFO << std::left << "-----------------------------------" << Root::GEndl;
1316 m_logger << Root::kINFO<< std::setw(10) << std::right << "Total" << std::setw(10) << std::right << "L1-Acc" << std::setw(10) << std::right << "L2-Acc" << std::setw(10) << std::right << "L3-Acc" <<
1317 std::setw(10) << std::right << "LiveTime" << std::setw(18) << std::right << "IntL rec.(ub^-1)" << std::setw(18) << std::right << "IntL del.(ub^-1)" << Root::GEndl;
1318 m_logger << Root::kINFO<< std::setw(10) << std::right << "" << std::setw(10) << std::right << m_totall1acc << std::setw(10) << std::right << m_totall2acc << std::setw(10) << std::right
1319 << m_totall3acc << std::setw(10) << std::right << m_totaltime << std::setw(18) << std::right << m_totalL << std::setw(18) << std::right << m_totalDelL << Root::GEndl;
1320 m_logger << Root::kINFO<< std::setw(10) << std::left << "Total L1/2/3 accept: " << std::setw(10) << std::left << m_totall1acc << std::setw(10) << std::left << m_totall2acc << std::setw(10)
1321 << std::left << m_totall3acc << Root::GEndl;
1322 if(m_uselivetrigger)m_logger << Root::kINFO<< std::setw(10) << std::left << "Total L1 livetime trigger accept: " << std::setw(10) << std::left << m_livtrig_totall1acc << Root::GEndl;
1323
1324 m_logger << Root::kINFO << std::setw(10) << std::left << "First Run: " << std::setw(10) << std::left << m_minrun << Root::GEndl;
1325 m_logger << Root::kINFO << std::setw(10) << std::left << "Last Run: " << std::setw(10) << std::left << m_maxrun << Root::GEndl;
1326 m_logger << Root::kINFO << std::setw(10) << std::left << "Total L1BeforePresc: " << std::setw(10) << std::left << m_totall1befpresc << Root::GEndl;
1327 m_logger << Root::kINFO<< std::setw(10) << std::left << "Total Livetime : " << m_totaltime << Root::GEndl;
1328 m_logger << Root::kINFO << std::setw(10) << std::left << "Total prescale weighted Livetime: " << m_totalPrescaleWLiveTime << Root::GEndl;
1329 m_logger << Root::kINFO<< std::setw(10) << std::left << "Total Good LBs : " << m_totalgoodblock - m_totalbadblock << Root::GEndl;
1330 m_logger << Root::kINFO<< std::setw(10) << std::left << "Total Bad LBs : " << m_totalbadblock << Root::GEndl;
1331 m_logger << Root::kINFO << std::setw(10) << std::left << "Total IntL delivered (ub^-1) : " << m_totalDelL << Root::GEndl;
1332 m_logger << Root::kINFO << std::setw(10) << std::left << "Total IntL after livefraction (ub^-1): " << m_lumiWOPrescale << Root::GEndl;
1333 if (m_uselar)
1334 m_logger << Root::kINFO << std::setw(10) << std::left << "Total IntL after LAr fraction (ub^-1): " << m_lumiLAr << Root::GEndl;
1335 m_logger << Root::kINFO << std::setw(10) << std::left << "Total IntL recorded (ub^-1) : " << m_totalL << Root::GEndl;
1336
1337 // ------------------------------------------------------------------------------------------------
1338
1339 if(m_makecollList == true){
1340 // store collisions xml file on demand
1341 if (m_collsgrl!=0) {
1342 TString collisionsxml = "collisions_" + m_collsgrl->GetSuggestedName() + ".xml";
1344 writer.SetGoodRunsList( *m_collsgrl );
1345 writer.SetFilename( collisionsxml.Data() );
1346 writer.WriteXMLFile();
1347 // can now delete grl
1348 delete m_collsgrl; m_collsgrl=0;
1349 }
1350 }
1351
1352
1353 // Creating monitoring plots on demand
1354 if(m_makePlots == true) { this->MakePlots(triggerchain); }
1355
1356 delete cq_trigger;
1357 delete cq_lumi;
1358 delete cq_lar;
1359 delete cq_bs;
1360}
1361
1362// Deal with composite L1 trigger
1363void
1364LumiCalculator::ParseL1Trigger(const std::string& lowerch, CoolQuery * cq_trigger) {
1365
1366 //
1367 // Check if we have multiple entries
1368 size_t last = 0;
1369 size_t next = 0;
1370 cool::ChannelId id;
1371 bool valid;
1372
1373 m_L1id = 0;
1374 m_L1Valid = false;
1375
1376 if (lowerch.find(',', last) == std::string::npos) {
1377 // Normal case
1378 m_L1id = cq_trigger->getL1ChannelId(lowerch, m_parlvl1menufolder);
1379 m_L1Valid = cq_trigger->channelIdValid();
1380
1381 } else {
1382 m_logger << Root::kINFO << "L1 item is composite: " << lowerch << Root::GEndl;
1383 do {
1384 next = lowerch.find(',', last);
1385 // Check if these are valid before using them
1386 id = cq_trigger->getL1ChannelId(lowerch.substr(last, next-last), m_parlvl1menufolder);
1387 valid = cq_trigger->channelIdValid();
1388 if (valid) {
1389 m_L1triggerchains.push_back(lowerch.substr(last, next-last));
1390 m_L1idList.push_back(m_L1id);
1391 m_L1id = id;
1392 m_L1Valid = true;
1393 } else {
1394 m_logger << Root::kINFO << lowerch.substr(last, next-last) << " Invalid" << Root::GEndl;
1395 }
1396
1397 last = next + 1;
1398 } while (next != std::string::npos);
1399
1400
1401 }
1402}
1403
1404// ---------------------------------------------------------------------------------
1405// Utility to print lumicalc summary results
1406void
1408
1409 os << std::left << "-----------------------------------" << std::endl;
1410 os << std::left << " LumiCalculator summary" << std::endl;
1411 os << std::left << "-----------------------------------" << std::endl;
1412 os << std::setw(10) << std::left << "Trigger: " << std::setw(10) << std::left << m_triggerchain << std::endl;
1413 os << std::setw(10) << std::left << "First Run: " << std::setw(10) << std::left << m_minrun << std::endl;
1414 os << std::setw(10) << std::left << "Last Run: " << std::setw(10) << std::left << m_maxrun << std::endl;
1415 os << std::setw(10) << std::left << "Total L1BeforePresc: " << std::setw(10) << std::left << m_totall1befpresc << std::endl;
1416 os << std::setw(10) << std::left << "Total Livetime : " << m_totaltime << std::endl;
1417 os << std::setw(10) << std::left << "Total prescale weighted Livetime: " << m_totalPrescaleWLiveTime << std::endl;
1418 os << std::setw(10) << std::left << "Total Good LBs : " << m_totalgoodblock - m_totalbadblock << std::endl;
1419 os << std::setw(10) << std::left << "Total Bad LBs : " << m_totalbadblock << std::endl;
1420 os << std::setw(10) << std::left << "Total IntL delivered (ub^-1) : " << m_totalDelL << std::endl;
1421 os << std::setw(10) << std::left << "Total IntL after livefraction (ub^-1): " << m_lumiWOPrescale << std::endl;
1422 if (m_uselar)
1423 os << std::setw(10) << std::left << "Total IntL after LAr fraction (ub^-1): " << m_lumiLAr << std::endl;
1424 os << std::setw(10) << std::left << "Total IntL recorded (ub^-1) : " << m_totalL << std::endl;
1425
1426}
1427
1428void
1429LumiCalculator::DoHistogramAdmin(const uint32_t& runnbr, const TString& trigName, const float& effxsec)
1430{
1431 // rebin the histograms once number of LBs is known
1432 int maxlb = 5000;
1433
1434 m_ntrigplbVec.push_back( new TH1F(Form("run%d_ntrigplb",runnbr), Form("Run %d",runnbr) , maxlb, 0., float(maxlb)) );
1435 m_trigrateplbVec.push_back( new TH1F(Form("run%d_trigrateplb",runnbr), Form("Run %d",runnbr) , maxlb, 0., float(maxlb)) );
1436 m_lumiplbVec.push_back( new TH1F(Form("run%d_peaklumiplb",runnbr), Form("Run %d",runnbr) , maxlb, 0., float(maxlb)) );
1437 m_lumitrigrateplbVec.push_back( new TH1F(Form("run%d_peaklumitrigrateplb",runnbr), Form("Run %d",runnbr) , maxlb, 0., float(maxlb)) );
1438 m_intlumiVec.push_back( new TH1F(Form("run%d_intlumi",runnbr), Form("Run %d",runnbr) , maxlb, 0., float(maxlb)) );
1439 m_intlumitrigrateVec.push_back( new TH1F(Form("run%d_intlumitrigrate",runnbr), Form("Run %d",runnbr) , maxlb, 0., float(maxlb)) );
1440 m_lumitrigrateplb_recordedVec.push_back( new TH1F(Form("run%d_peakrecordedlumitrigrateplb",runnbr), Form("Run %d",runnbr) , maxlb, 0., float(maxlb)) );
1441 m_intlumitrigrate_recordedVec.push_back( new TH1F(Form("run%d_intrecordedlumitrigrate",runnbr), Form("Run %d",runnbr) , maxlb, 0., float(maxlb)) );
1442
1443 m_ntrigplb = *(m_ntrigplbVec.rbegin());
1444 m_trigrateplb = *(m_trigrateplbVec.rbegin());
1445 m_lumiplb = *(m_lumiplbVec.rbegin());
1447 m_intlumi = *(m_intlumiVec.rbegin());
1451
1452 this->SetHistogramStyle(m_ntrigplb, Form("Run = %d",runnbr), "Luminosity block number", Form("# %s triggers / LB",trigName.Data()));
1453 this->SetHistogramStyle(m_trigrateplb, Form("Run = %d",runnbr), "Luminosity block number", Form("%s trigger rate / LB",trigName.Data()));
1454 this->SetHistogramStyle(m_lumiplb, Form("Run = %d",runnbr), "Luminosity block number", "Delivered luminosity (#mub^{-1}/s)");
1455 this->SetHistogramStyle(m_lumitrigrateplb, Form("Efficiency * x-sec = %.1f #mub, Run = %d",effxsec,runnbr), "Luminosity block number", Form("%s luminosity (#mub^{-1}/s)",trigName.Data()));
1456 this->SetHistogramStyle(m_intlumi, Form("Run = %d",runnbr), "Luminosity block number", "Integrated delivered luminosity (#mub^{-1})");
1457 this->SetHistogramStyle(m_intlumitrigrate, Form("Efficiency * x-sec = %.1f #mub, Run = %d",effxsec,runnbr), "Luminosity block number", Form("%s Integrated luminosity (#mub^{-1})",trigName.Data()));
1458 this->SetHistogramStyle(m_lumitrigrateplb_recorded, Form("Efficiency * x-sec = %.1f #mub, Run = %d",effxsec,runnbr), "Luminosity block number", Form("%s Recorded luminosity (#mub^{-1}/s)",trigName.Data()));
1459 this->SetHistogramStyle(m_intlumitrigrate_recorded, Form("Efficiency * x-sec = %.1f #mub, Run = %d",effxsec,runnbr), "Luminosity block number", Form("%s Integrated recorded luminosity (#mub^{-1})",trigName.Data()));
1460}
1461
1462
1463void
1464LumiCalculator::RebinHistograms(const int& nbins, const double& start, const double& end)
1465{
1466 m_ntrigplb->SetBins(nbins,start,end);
1467 m_trigrateplb->SetBins(nbins,start,end);
1468 m_lumiplb->SetBins(nbins,start,end);
1469 m_lumitrigrateplb->SetBins(nbins,start,end);
1470 m_intlumi->SetBins(nbins,start,end);
1471 m_intlumitrigrate->SetBins(nbins,start,end);
1472 m_lumitrigrateplb_recorded->SetBins(nbins,start,end);
1473 m_intlumitrigrate_recorded->SetBins(nbins,start,end);
1474}
1475
1476
1477void
1478LumiCalculator::SetHistogramStyle(TH1F* hist, const char* title, const char* xaxis, const char* yaxis)
1479{
1480 hist->SetFillColor(33); // light grey, blueish
1481 if (title!=0) hist->SetTitle(title);
1482 if (xaxis!=0) hist->GetXaxis()->SetTitle(xaxis);
1483 if (yaxis!=0) hist->GetYaxis()->SetTitle(yaxis);
1484
1485 hist->GetXaxis()->SetLabelFont(52);
1486 hist->GetXaxis()->SetLabelSize(0.04);
1487 hist->GetXaxis()->SetTitleSize(0.05);
1488 hist->GetXaxis()->SetTitleOffset(1.28);
1489 hist->GetXaxis()->SetTitleFont(42);
1490
1491 hist->GetYaxis()->SetLabelFont(52);
1492 hist->GetYaxis()->SetLabelSize(0.04);
1493 hist->GetYaxis()->SetTitleSize(0.05);
1494 hist->GetYaxis()->SetTitleOffset(1.25);
1495 hist->GetYaxis()->SetTitleFont(42);
1496
1497 hist->SetStats(false);
1498
1499 hist->SetLineWidth(2);
1500}
1501
1502void
1504{
1505 if(m_LumiTree != 0){
1506 // rebin and fill histograms:
1507 m_LumiTree->SetBranchAddress("LBStart", &m_clumiblocknbr);
1508 m_LumiTree->SetBranchAddress("L1AfterPrescale", &m_afterprescale);
1509 m_LumiTree->SetBranchAddress("L1Rate", &m_l1rate);
1510 m_LumiTree->SetBranchAddress("IntLumi",&m_intLumi);
1511 m_LumiTree->SetBranchAddress("L1Ratediveffxsec",&m_l1ratediveffxsec);
1512 m_LumiTree->SetBranchAddress("TotalLumi",&m_totalL);
1513 m_LumiTree->SetBranchAddress("Total_L1Ratediveffxsec",&m_total_l1ratediveffxsec);
1514 m_LumiTree->SetBranchAddress("TotalLumiRun",&m_totalLRun);
1515 m_LumiTree->SetBranchAddress("Total_L1RatediveffxsecRun",&m_total_l1ratediveffxsecRun);
1516 m_LumiTree->SetBranchAddress("RunNbr", &m_runnbr);
1517 m_LumiTree->SetBranchAddress("L1RatediveffxsecRecorded",&m_l1ratediveffxsec_recorded);
1518 m_LumiTree->SetBranchAddress("Total_L1RatediveffxsecRecorded",&m_total_l1ratediveffxsec_recorded);
1519 m_LumiTree->SetBranchAddress("Total_L1RatediveffxsecRunRecorded",&m_total_l1ratediveffxsecRun_recorded);
1520 m_LumiTree->SetBranchAddress("AvergeInteractionPerXing",&m_AvEvtsPerBX);
1521
1522 // get first and last run number
1523 m_LumiTree->GetEntry(0);
1524 int runnbrstart = m_runnbr;
1525 m_LumiTree->GetEntry(m_LumiTree->GetEntries()-1);
1526 int runnbrend = m_runnbr;
1527
1528 // makeup of integrated lumi histograms
1529 m_intlumiruns = new TH1F("intlumiruns","Luminosity",1,0.,1.);
1530 m_intlumitrigrateruns = new TH1F("intlumitrigrateruns","Delivered luminosity",1,0.,1.);
1531 m_intlumitrigrateruns_recorded = new TH1F("intlumitrigrateruns_recorded","Recorded luminosity",1,0.,1.);
1532 m_intlumiruns->SetBins(runnbrend-runnbrstart+10,float(runnbrstart),float(runnbrend+10));
1533 m_intlumitrigrateruns->SetBins(runnbrend-runnbrstart+10,float(runnbrstart),float(runnbrend+10));
1534 m_intlumitrigrateruns_recorded->SetBins(runnbrend-runnbrstart+10,float(runnbrstart),float(runnbrend+10));
1535
1536 // Lumi-weighted average interactions per crossing
1537 m_avgintperbx = new TH1F("avgintperbx", "Avg Int/BX", 1000, 0., 100.);
1538 m_avgintperbx->SetTitle("Lumi-weighted Interactions per BX");
1539 this->SetHistogramStyle(m_avgintperbx, "Lumi-weighted Average Interactions per BX", "Average Interactions per BX", "Recorded Luminosity (mb-1)");
1540
1541 // loop over the lumi tree
1542 m_runnbr_prev = 0;
1543 int nlbs(0);
1544 float totalL(0.), total_l1ratediveffxsec(0.), total_l1ratediveffxsec_recorded(0.), totalLRun(0.), total_l1ratediveffxsecRun(0.), total_l1ratediveffxsecRun_recorded(0.);
1545 for(int i=0; i < m_LumiTree->GetEntries(); i++){
1546 m_LumiTree->GetEntry(i);
1547
1548 // do histogram admin first
1549 if ( m_runnbr!=m_runnbr_prev ) {
1550 // first rebin prev histograms
1551 if (m_runnbr_prev>0) {
1552 m_intlumi->SetTitle(Form("Delivered luminosity = %.1f /#mu b, Run = %d",totalLRun,m_runnbr_prev));
1553 m_intlumitrigrate->SetTitle(Form("Delivered luminosity = %.1f /#mu b, Efficiency * x-sec = %.1f #mu b, Run = %d",total_l1ratediveffxsecRun,m_effxsec,m_runnbr_prev));
1554 m_intlumitrigrate_recorded->SetTitle(Form("Recorded luminosity = %.1f /#mu b, Efficiency * x-sec = %.1f #mu b, Run = %d",total_l1ratediveffxsecRun_recorded,m_effxsec,m_runnbr_prev));
1555 this->RebinHistograms(nlbs+10,0,double(nlbs+10));
1556 }
1557 // create new histograms
1559 // fill cullumative luminosity
1560 if (m_runnbr_prev>0) {
1561 for (uint32_t j=m_runnbr_prev; j<m_runnbr; ++j) {
1562 m_intlumiruns->Fill(j,totalL);
1563 m_intlumitrigrateruns->Fill(j,total_l1ratediveffxsec);
1564 m_intlumitrigrateruns_recorded->Fill(j,total_l1ratediveffxsec_recorded);
1565 }
1566 }
1568 }
1569 totalL = m_totalL;
1570 totalLRun = m_totalLRun;
1571 total_l1ratediveffxsec = m_total_l1ratediveffxsec;
1572 total_l1ratediveffxsecRun = m_total_l1ratediveffxsecRun;
1573 total_l1ratediveffxsec_recorded = m_total_l1ratediveffxsec_recorded;
1574 total_l1ratediveffxsecRun_recorded = m_total_l1ratediveffxsecRun_recorded;
1575 nlbs = m_clumiblocknbr;
1576
1585 m_avgintperbx->Fill(m_AvEvtsPerBX, m_intLumi); // Lumi-weighted mu
1586
1587 } // end tree loop
1588 m_intlumiruns->Fill(runnbrend,totalL);
1589 m_intlumitrigrateruns->Fill(runnbrend,total_l1ratediveffxsec);
1590 m_intlumitrigrateruns_recorded->Fill(runnbrend,total_l1ratediveffxsec_recorded);
1591
1592 // finish histograms make-up
1593 m_intlumiruns->SetMinimum(0.);
1594 m_intlumitrigrateruns->SetMinimum(0.);
1595 this->RebinHistograms(nlbs+10,0,double(nlbs+10));
1596 this->SetHistogramStyle(m_intlumiruns, Form("Delivered luminosity = %.1f /ub",totalL), "Run number", "Luminosity (#mu b^{-1})");
1597 this->SetHistogramStyle(m_intlumitrigrateruns, Form("Delivered luminosity = %.1f /#mu b, Recorded luminosity = %.1f /#mu b", //, Efficiency * x-sec = %.1f mb",
1598 total_l1ratediveffxsec,total_l1ratediveffxsec_recorded/*,m_effxsec*/), "Run number", Form("%s Luminosity (#mu b^{-1})",triggerchain.c_str()));
1599 this->SetHistogramStyle(m_intlumitrigrateruns_recorded, Form("Delivered luminosity = %.1f /#mub, Recorded luminosity = %.1f /#mub", //, Efficiency * x-sec = %.1f #mub",
1600 total_l1ratediveffxsec,total_l1ratediveffxsec_recorded/*,m_effxsec*/), "Run number", Form("%s Luminosity (#mu b^{-1})",triggerchain.c_str()));
1601
1602 std::vector<TH1F*>::iterator itr;
1603
1604 // and store the histograms
1605 TString histFileName = TString("ilumicalc_histograms_") + TString(triggerchain) + ( runnbrstart==runnbrend ? Form("_%d_",runnbrstart) : Form("_%d-%d_",runnbrstart,runnbrend)) + TString(m_lumitag) + TString(".root");
1606 TFile *ff = new TFile(histFileName.Data(),"recreate");
1607 m_avgintperbx->Write();
1608 if (m_effxsec==1.0) {
1609 m_intlumiruns->Write();
1610 for (itr=m_lumiplbVec.begin(); itr!=m_lumiplbVec.end(); ++itr) { (*itr)->Write(); }
1611 for (itr=m_intlumiVec.begin(); itr!=m_intlumiVec.end(); ++itr) { (*itr)->Write(); }
1612 }
1613 for (itr=m_ntrigplbVec.begin(); itr!=m_ntrigplbVec.end(); ++itr) { (*itr)->Write(); }
1614 for (itr=m_trigrateplbVec.begin(); itr!=m_trigrateplbVec.end(); ++itr) { (*itr)->Write(); }
1615 if (m_effxsec!=1.0) { // results only make sense when proper xsec is provided externally
1616 m_intlumitrigrateruns->Write();
1618 for (itr=m_lumitrigrateplbVec.begin(); itr!=m_lumitrigrateplbVec.end(); ++itr) { (*itr)->Write(); }
1619 for (itr=m_intlumitrigrateVec.begin(); itr!=m_intlumitrigrateVec.end(); ++itr) { (*itr)->Write(); }
1620 for (itr=m_lumitrigrateplb_recordedVec.begin(); itr!=m_lumitrigrateplb_recordedVec.end(); ++itr) { (*itr)->Write(); }
1621 for (itr=m_intlumitrigrate_recordedVec.begin(); itr!=m_intlumitrigrate_recordedVec.end(); ++itr) { (*itr)->Write(); }
1622 }
1623 m_LumiTree->Write();
1624
1625 // And write out the lumi tag information
1626 TObjString lumiTag(m_lumitag.c_str());
1627 lumiTag.Write("lumiTag");
1628
1629 TObjString larTag(m_lartag.c_str());
1630 larTag.Write("larTag");
1631
1632
1633 ff->Close();
1634 delete ff;
1635 m_logger << Root::kINFO << "Histograms stored as : " << histFileName << Root::GEndl;
1636
1637 }else{
1638 m_logger << Root::kWARNING << "LumiTree pointer does not exist! : " << Root::GEndl;
1639 }
1640
1641
1642}
1643
1644
std::shared_ptr< HepMC3::Writer > writer
#define ATLAS_NOT_THREAD_SAFE
getNoisyStrip() Find noisy strips from hitmaps and write out into xml/db formats
cool::ChannelId getHLTChannelId(const std::string &trigger, const std::string &folder_name)
void setIOV(const cool::ValidityKey start, const cool::ValidityKey stop)
Definition CoolQuery.cxx:96
IOVData< T > getIOVData(const std::string &name, const std::string &folder_name, const cool::ChannelId &id, const std::string &tag="")
Definition CoolQuery.h:426
std::map< cool::ValidityKey, L1CountFolderData > getL1CountFolderData(const std::string &folder_name, const cool::ChannelId &id)
unsigned int getTriggerLevel(const std::string &triggername)
Definition CoolQuery.cxx:79
void setIOVForRun(unsigned int runnum)
void printHLTTriggers(const std::string &folder_name)
std::map< cool::ValidityKey, LumiFolderData > getLumiFolderData(const std::string &folder_name, const std::string &tag, const cool::ChannelId &id)
void printL1Triggers(const std::string &folder_name)
bool openDbConn()
Definition CoolQuery.cxx:30
cool::ChannelId getL1ChannelId(const std::string &trigger, const std::string &folder_name)
cool::ChannelId getLumiChannelId(const std::string &lumimethod, const std::string &folder_name)
bool channelIdValid()
std::string getHLTLowerChainName(const std::string &trigger, const std::string &folder_name)
std::map< cool::ValidityKey, T > getObjMapFromFolderAtChan(const std::string &obj_name, const std::string &folder_name, const cool::ChannelId &id)
Definition CoolQuery.h:380
DataModel_detail::const_iterator< DataVector > const_iterator
Definition DataVector.h:838
const_iterator end() const noexcept
Return a const_iterator pointing past the end of the collection.
const_iterator begin() const noexcept
Return a const_iterator pointing at the beginning of the collection.
void clear()
Definition CoolQuery.h:70
T getValue(IOVTime time)
Definition CoolQuery.h:104
std::list< std::pair< IOVRange, T > > getOverlap(const IOVRange &range)
Definition CoolQuery.h:178
std::list< std::pair< IOVRange, T > > data
Definition CoolQuery.h:86
Validity Range object.
Definition IOVRange.h:30
const IOVTime & stop() const
Definition IOVRange.h:39
const IOVTime & start() const
Definition IOVRange.h:38
Basic time unit for IOVSvc.
Definition IOVTime.h:33
void setRETime(uint64_t time) noexcept
Definition IOVTime.cxx:84
uint32_t event() const noexcept
Definition IOVTime.h:106
uint32_t run() const noexcept
Definition IOVTime.h:105
uint64_t re_time() const noexcept
Definition IOVTime.h:107
unsigned int m_t_totalbadblock
std::vector< TH1F * > m_intlumitrigrate_recordedVec
cool::ChannelId m_L3id
unsigned int m_t_l3acc
std::string m_bs_database
unsigned int m_l2acc
cool::ChannelId m_Lumiid
unsigned int m_t_l1acc
Root::TMsgLogger m_logger
std::string m_bsonl
std::string m_lumionl
std::string m_parlumilvl1folder
float m_total_l1ratediveffxsecRun_recorded
unsigned int m_livtrig_totall1acc
std::string m_parlvl1menufolder
void UseOnlineLumi(bool online)
std::vector< TH1F * > m_intlumiVec
float m_l1ratediveffxsec_recorded
std::string m_lartag
unsigned int m_l1acc
void UseLumiChannel(int chan)
void UseLumiMethod(const std::string &method)
void SetHistogramStyle(TH1F *hist, const char *title=0, const char *xaxis=0, const char *yaxis=0)
ULong64_t m_livetime_afterprescale
cool::ChannelId m_L1id
std::vector< std::string > m_L1triggerchains
void RebinHistograms(const int &nbins, const double &start, const double &end)
std::string m_livetrigger
std::string m_lumimethod
void setTree(TTree *tree=0)
float m_total_l1ratediveffxsec
unsigned int m_t_totalgoodblock
float m_total_l1ratediveffxsecRun
std::string m_lbcollname
float m_total_l1ratediveffxsec_recorded
void UseLArNoiseDB(bool lar, const std::string &lardb)
std::string m_trigger
unsigned int m_t_totall1befpresc
std::string m_parlvl1prescalesfolder
float m_t_totalPrescaleWLiveTime
std::vector< TH1F * > m_trigrateplbVec
uint32_t m_lbstop_prev
unsigned int m_totalbadblock
std::string m_parlareventvetofolder
ULong64_t m_afterprescale
unsigned int m_l3acc
std::string m_lar_database
uint32_t m_lbstart_prev
unsigned int m_livetime_l1acc
std::vector< TH1F * > m_lumitrigrateplbVec
std::string m_parlvl1lblbfolder
unsigned int m_totall1befpresc
void UseMC(bool mc=true)
void SetCollName(const std::string &lbcollname)
std::string m_parhltmenufolder
std::vector< cool::ChannelId > m_L1idList
void UseLiveTrigger(bool live, std::string &livetrigger)
Root::TGoodRunsList * m_collsgrl
void DoHistogramAdmin(const uint32_t &runnbr, const TString &trigName, const float &effxsec)
unsigned int m_totall1acc
std::string m_trig_database
TH1F * m_intlumitrigrateruns
unsigned int m_totall3acc
unsigned int m_t_l2acc
uint32_t m_runnbr_prev
std::string m_lumi_database
void printSummary(std::ostream &os)
void ParseL1Trigger(const std::string &lowerch, CoolQuery *cq_trigger)
std::string m_parlumihltfolder
std::string m_bstag
std::vector< TH1F * > m_intlumitrigrateVec
std::vector< TH1F * > m_ntrigplbVec
std::string m_lumioff
ULong64_t m_livetime_beforeprescale
uint32_t m_clumiblocknbrend
void Verbose(bool verbose=true)
void UseBeamspot(bool bs, const std::string &bstag)
cool::ChannelId m_LiveL1id
TH1F * m_intlumitrigrateruns_recorded
std::string m_lumitag
ULong64_t m_beforeprescale
std::vector< TH1F * > m_lumiplbVec
std::string m_parhltprescalesfolder
void MakePlots(bool plots)
unsigned int m_totall2acc
std::string m_paronlbeamspotfolder
std::vector< TH1F * > m_lumitrigrateplb_recordedVec
unsigned int m_maxrun
std::string m_data_db
TH1F * m_intlumitrigrate_recorded
unsigned int m_triglevel
unsigned int m_totalgoodblock
float m_totalPrescaleWLiveTime
std::string m_laroff
TH1F * m_lumitrigrateplb_recorded
void UseLumiTag(const std::string &tag)
std::string m_triggerchain
cool::ChannelId m_L2id
uint32_t m_clumiblocknbr
unsigned int m_minrun
bool verbose
Definition hcg.cxx:73
static std::vector< std::string > triggerchain
Definition iLumiCalc.h:34
static std::vector< uint32_t > lbstart
Definition iLumiCalc.h:38
static std::string livetrigger
Definition iLumiCalc.h:35
static std::vector< uint32_t > lbend
Definition iLumiCalc.h:39
@ kWARNING
Definition TMsgLogger.h:41
@ kINFO
Definition TMsgLogger.h:40
LumiBlockRangeContainer_v1 LumiBlockRangeContainer
Declare the latest version of the container.
TChain * tree