Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:07:59

0001 #include "DQM/L1TMonitorClient/interface/L1TDTTPGClient.h"
0002 
0003 #include "FWCore/ServiceRegistry/interface/Service.h"
0004 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0005 #include "FWCore/Framework/interface/ESHandle.h"
0006 #include "FWCore/Framework/interface/EventSetup.h"
0007 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0008 #include "DQMServices/Core/interface/DQMStore.h"
0009 #include "TRandom.h"
0010 
0011 #include <TF1.h>
0012 #include <cstdio>
0013 #include <sstream>
0014 #include <cmath>
0015 #include <TProfile.h>
0016 #include <TProfile2D.h>
0017 
0018 using namespace edm;
0019 using namespace std;
0020 
0021 L1TDTTPGClient::L1TDTTPGClient(const edm::ParameterSet &ps) {
0022   parameters_ = ps;
0023   initialize();
0024 }
0025 
0026 L1TDTTPGClient::~L1TDTTPGClient() { LogInfo("TriggerDQM") << "[TriggerDQM]: ending... "; }
0027 
0028 //--------------------------------------------------------
0029 void L1TDTTPGClient::initialize() {
0030   counterLS_ = 0;
0031   counterEvt_ = 0;
0032 
0033   // get back-end interface
0034   //dbe_ = Service<DQMStore>().operator->();
0035 
0036   // base folder for the contents of this job
0037   monitorName_ = parameters_.getUntrackedParameter<string>("monitorName", "");
0038   //  cout << "Monitor name = " << monitorName_ << endl;
0039   prescaleLS_ = parameters_.getUntrackedParameter<int>("prescaleLS", -1);
0040   //  cout << "DQM lumi section prescale = " << prescaleLS_ << " lumi section(s)"<< endl;
0041   prescaleEvt_ = parameters_.getUntrackedParameter<int>("prescaleEvt", -1);
0042   //  cout << "DQM event prescale = " << prescaleEvt_ << " events(s)"<< endl;
0043   output_dir_ = parameters_.getUntrackedParameter<string>("output_dir", "");
0044   //  cout << "DQM output dir = " << output_dir_ << endl;
0045   input_dir_ = parameters_.getUntrackedParameter<string>("input_dir", "");
0046   //  cout << "DQM input dir = " << input_dir_ << endl;
0047 
0048   LogInfo("TriggerDQM");
0049 }
0050 
0051 //--------------------------------------------------------
0052 void L1TDTTPGClient::dqmEndJob(DQMStore::IBooker &ibooker, DQMStore::IGetter &igetter) {
0053   ibooker.setCurrentFolder(output_dir_);
0054 
0055   // booking
0056 
0057   dttpgphmapcorrf = ibooker.book2D(
0058       "DT_TPG_phi_map_corr_frac", "Fraction of correlated best triggers per station", 20, 1, 21, 12, 0, 12);
0059   dttpgphmap2ndf =
0060       ibooker.book2D("DT_TPG_phi_map_2nd_frac", "Fraction of second tracks per station", 20, 1, 21, 12, 0, 12);
0061   dttpgphmapbxf[0] =
0062       ibooker.book2D("DT_TPG_phi_map_bx-1_frac", "Fraction of triggers per station (BX=-1)", 20, 1, 21, 12, 0, 12);
0063   dttpgphmapbxf[1] =
0064       ibooker.book2D("DT_TPG_phi_map_bx0_frac", "Fraction of triggers per station (BX=0)", 20, 1, 21, 12, 0, 12);
0065   dttpgphmapbxf[2] =
0066       ibooker.book2D("DT_TPG_phi_map_bx+1_frac", "Fraction of triggers per station (BX=1)", 20, 1, 21, 12, 0, 12);
0067   setMapPhLabel(dttpgphmapcorrf);
0068   setMapPhLabel(dttpgphmap2ndf);
0069   setMapPhLabel(dttpgphmapbxf[0]);
0070   setMapPhLabel(dttpgphmapbxf[1]);
0071   setMapPhLabel(dttpgphmapbxf[2]);
0072 
0073   dttpgthmaphf = ibooker.book2D(
0074       "DT_TPG_theta_map_corr_frac", "Fraction of H quality best triggers per station", 15, 1, 16, 12, 0, 12);
0075   dttpgthmapbxf[0] =
0076       ibooker.book2D("DT_TPG_theta_map_bx-1_frac", "Fraction of triggers per station (BX=-1)", 15, 1, 16, 12, 0, 12);
0077   dttpgthmapbxf[1] =
0078       ibooker.book2D("DT_TPG_theta_map_bx0_frac", "Fraction of triggers per station (BX=0)", 15, 1, 16, 12, 0, 12);
0079   dttpgthmapbxf[2] =
0080       ibooker.book2D("DT_TPG_theta_map_bx+1_frac", "Fraction of triggers per station (BX=1)", 15, 1, 16, 12, 0, 12);
0081   setMapThLabel(dttpgthmaphf);
0082   setMapThLabel(dttpgthmapbxf[0]);
0083   setMapThLabel(dttpgthmapbxf[1]);
0084   setMapThLabel(dttpgthmapbxf[2]);
0085 
0086   //   cout << "L1TDTTPGClient::analyze" << endl;
0087   counterEvt_++;
0088   if (prescaleEvt_ < 1)
0089     return;
0090   if (prescaleEvt_ > 0 && counterEvt_ % prescaleEvt_ != 0)
0091     return;
0092 
0093   string nName = "DT_TPG_phi_best_map_corr";
0094   string dName = "DT_TPG_phi_best_map";
0095   makeRatioHisto(igetter, dttpgphmapcorrf, nName, dName);
0096   dName = "DT_TPG_phi_map";
0097   nName = "DT_TPG_phi_map_2nd";
0098   makeRatioHisto(igetter, dttpgphmap2ndf, nName, dName);
0099   nName = "DT_TPG_phi_map_bx-1";
0100   makeRatioHisto(igetter, dttpgphmapbxf[0], nName, dName);
0101   nName = "DT_TPG_phi_map_bx0";
0102   makeRatioHisto(igetter, dttpgphmapbxf[1], nName, dName);
0103   nName = "DT_TPG_phi_map_bx+1";
0104   makeRatioHisto(igetter, dttpgphmapbxf[2], nName, dName);
0105 
0106   nName = "DT_TPG_theta_best_map_h";
0107   dName = "DT_TPG_theta_best_map";
0108   makeRatioHisto(igetter, dttpgthmaphf, nName, dName);
0109   dName = "DT_TPG_theta_map";
0110   nName = "DT_TPG_theta_map_bx-1";
0111   makeRatioHisto(igetter, dttpgthmapbxf[0], nName, dName);
0112   nName = "DT_TPG_theta_map_bx0";
0113   makeRatioHisto(igetter, dttpgthmapbxf[1], nName, dName);
0114   nName = "DT_TPG_theta_map_bx+1";
0115   makeRatioHisto(igetter, dttpgthmapbxf[2], nName, dName);
0116 }
0117 
0118 void L1TDTTPGClient::makeRatioHisto(DQMStore::IGetter &igetter, MonitorElement *ratioME, string &nName, string &dName) {
0119   TH2F *numerator;
0120   TH2F *denominator;
0121 
0122   denominator = this->get2DHisto(input_dir_ + "/" + dName, igetter);
0123   numerator = this->get2DHisto(input_dir_ + "/" + nName, igetter);
0124 
0125   if (numerator && denominator) {
0126     TH2F *ratio = ratioME->getTH2F();
0127     if (ratio) {
0128       ratio->Divide(numerator, denominator);
0129     } else {
0130       LogInfo("TriggerDQM") << "[TriggerDQM]: ratio histo named \"" << ratioME->getName() << "\" not found!" << endl;
0131     }
0132   } else {
0133     if (!numerator)
0134       LogInfo("TriggerDQM") << "[TriggerDQM]: numerator histo \"" << nName << "\" not found!" << endl;
0135     if (!denominator)
0136       LogInfo("TriggerDQM") << "[TriggerDQM]: denominator histo \"" << dName << "\" not found!" << endl;
0137   }
0138 }
0139 
0140 TH1F *L1TDTTPGClient::get1DHisto(string meName, DQMStore::IGetter &igetter) {
0141   MonitorElement *me_ = igetter.get(meName);
0142 
0143   if (!me_) {
0144     LogInfo("TriggerDQM") << "ME NOT FOUND.";
0145     return nullptr;
0146   }
0147 
0148   return me_->getTH1F();
0149 }
0150 
0151 TH2F *L1TDTTPGClient::get2DHisto(string meName, DQMStore::IGetter &igetter) {
0152   MonitorElement *me_ = igetter.get(meName);
0153 
0154   if (!me_) {
0155     LogInfo("TriggerDQM") << "ME NOT FOUND.";
0156     return nullptr;
0157   }
0158 
0159   return me_->getTH2F();
0160 }
0161 
0162 TProfile2D *L1TDTTPGClient::get2DProfile(string meName, DQMStore::IGetter &igetter) {
0163   MonitorElement *me_ = igetter.get(meName);
0164 
0165   if (!me_) {
0166     LogInfo("TriggerDQM") << "ME NOT FOUND.";
0167     return nullptr;
0168   }
0169 
0170   return me_->getTProfile2D();
0171 }
0172 
0173 TProfile *L1TDTTPGClient::get1DProfile(string meName, DQMStore::IGetter &igetter) {
0174   MonitorElement *me_ = igetter.get(meName);
0175 
0176   if (!me_) {
0177     LogInfo("TriggerDQM") << "ME NOT FOUND.";
0178     return nullptr;
0179   }
0180 
0181   return me_->getTProfile();
0182 }
0183 
0184 void L1TDTTPGClient::setMapPhLabel(MonitorElement *me) {
0185   me->setAxisTitle("DTTF Sector", 2);
0186   for (int i = 0; i < 5; i++) {
0187     ostringstream wheel;
0188     wheel << i - 2;
0189     me->setBinLabel(1 + i * 4, "Wheel " + wheel.str(), 1);
0190   }
0191 }
0192 
0193 void L1TDTTPGClient::setMapThLabel(MonitorElement *me) {
0194   me->setAxisTitle("DTTF Sector", 2);
0195   for (int i = 0; i < 5; i++) {
0196     ostringstream wheel;
0197     wheel << i - 2;
0198     me->setBinLabel(1 + i * 3, "Wheel " + wheel.str(), 1);
0199   }
0200 }