Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 11:58:30

0001 
0002 /*
0003  *  See header file for a description of this class.
0004  *
0005  *  \author S. Bolognesi - INFN Torino
0006  */
0007 
0008 #include "DTTTrigAnalyzer.h"
0009 
0010 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0011 #include "FWCore/Framework/interface/EventSetup.h"
0012 #include "CondFormats/DTObjects/interface/DTTtrig.h"
0013 #include "CondFormats/DataRecord/interface/DTTtrigRcd.h"
0014 #include <iostream>
0015 #include "TFile.h"
0016 #include "TH1D.h"
0017 #include "TString.h"
0018 
0019 using namespace edm;
0020 using namespace std;
0021 
0022 DTTTrigAnalyzer::DTTTrigAnalyzer(const ParameterSet &pset) {
0023   // The root file which will contain the histos
0024   string rootFileName = pset.getUntrackedParameter<string>("rootFileName");
0025   theFile = new TFile(rootFileName.c_str(), "RECREATE");
0026   theFile->cd();
0027   //The k factor to compute ttrig
0028   //kfactor = pset.getUntrackedParameter<double>("kfactor",0);
0029   ttrigToken_ =
0030       esConsumes<edm::Transition::BeginRun>(edm::ESInputTag("", pset.getUntrackedParameter<string>("dbLabel")));
0031 }
0032 
0033 DTTTrigAnalyzer::~DTTTrigAnalyzer() { theFile->Close(); }
0034 
0035 void DTTTrigAnalyzer::beginRun(const edm::Run &, const edm::EventSetup &eventSetup) {
0036   ESHandle<DTTtrig> tTrig = eventSetup.getHandle(ttrigToken_);
0037   tTrigMap = &*tTrig;
0038   cout << "[DTTTrigAnalyzer] TTrig version: " << tTrig->version() << endl;
0039 }
0040 
0041 void DTTTrigAnalyzer::endJob() {
0042   //   static const double convToNs = 25./32.;
0043   // Loop over DB entries
0044   for (DTTtrig::const_iterator ttrig = tTrigMap->begin(); ttrig != tTrigMap->end(); ++ttrig) {
0045     DTWireId wireId(
0046         (*ttrig).first.wheelId, (*ttrig).first.stationId, (*ttrig).first.sectorId, (*ttrig).first.slId, 0, 0);
0047     float tmean;
0048     float sigma;
0049     float kFactor;
0050 
0051     DetId detId(wireId.rawId());
0052     // ttrig and rms are ns
0053     tTrigMap->get(detId, tmean, sigma, kFactor, DTTimeUnits::ns);
0054     float ttcor = tmean + kFactor * sigma;
0055     cout << "Wire: " << wireId << endl
0056          << " Ttrig (ns): " << ttcor << endl
0057          << " tmean (ns): " << tmean << endl
0058          << " sigma (ns): " << sigma << endl
0059          << " kfactor: " << kFactor << endl;
0060 
0061     //Define an histo for each wheel and each superlayer type
0062     TH1D *hTTrigHisto = theTTrigHistoMap[make_pair(wireId.wheel(), wireId.superlayer())];
0063     TH1D *hTMeanHisto = theTMeanHistoMap[make_pair(wireId.wheel(), wireId.superlayer())];
0064     TH1D *hSigmaHisto = theSigmaHistoMap[make_pair(wireId.wheel(), wireId.superlayer())];
0065     TH1D *hKFactorHisto = theKFactorHistoMap[make_pair(wireId.wheel(), wireId.superlayer())];
0066     if (hTTrigHisto == 0) {
0067       theFile->cd();
0068       TString name = getHistoName(wireId).c_str();
0069       if (wireId.superlayer() != 2) {
0070         hTTrigHisto = new TH1D(name + "_TTrig", "TTrig calibrated from TB per superlayer", 50, 0, 50);
0071         hTMeanHisto = new TH1D(name + "_TMean", "TMean calibrated from TB per superlayer", 50, 0, 50);
0072         hSigmaHisto = new TH1D(name + "_Sigma", "Sigma calibrated from TB per superlayer", 50, 0, 50);
0073         hKFactorHisto = new TH1D(name + "_KFactor", "KFactor calibrated from TB per superlayer", 50, 0, 50);
0074       } else {
0075         hTTrigHisto = new TH1D(name + "_TTrig", "TTrig calibrated from TB per superlayer", 36, 0, 36);
0076         hTMeanHisto = new TH1D(name + "_TMean", "TMean calibrated from TB per superlayer", 36, 0, 36);
0077         hSigmaHisto = new TH1D(name + "_Sigma", "Sigma calibrated from TB per superlayer", 36, 0, 36);
0078         hKFactorHisto = new TH1D(name + "_KFactor", "KFactor calibrated from TB per superlayer", 36, 0, 36);
0079       }
0080       theTTrigHistoMap[make_pair(wireId.wheel(), wireId.superlayer())] = hTTrigHisto;
0081       theTMeanHistoMap[make_pair(wireId.wheel(), wireId.superlayer())] = hTMeanHisto;
0082       theSigmaHistoMap[make_pair(wireId.wheel(), wireId.superlayer())] = hSigmaHisto;
0083       theKFactorHistoMap[make_pair(wireId.wheel(), wireId.superlayer())] = hKFactorHisto;
0084     }
0085 
0086     //Fill the histos and set the bin label
0087     int binNumber = wireId.sector() + 12 * (wireId.station() - 1);
0088     //    hTTrigHisto->SetBinContent(binNumber,ttrig);
0089     hTTrigHisto->SetBinContent(binNumber, ttcor);
0090     hTMeanHisto->SetBinContent(binNumber, tmean);
0091     hSigmaHisto->SetBinContent(binNumber, sigma);
0092     hKFactorHisto->SetBinContent(binNumber, kFactor);
0093     string labelName;
0094     stringstream theStream;
0095     if (wireId.sector() == 1)
0096       theStream << "MB" << wireId.station() << "_Sec" << wireId.sector();
0097     else
0098       theStream << "Sec" << wireId.sector();
0099     theStream >> labelName;
0100     hTTrigHisto->GetXaxis()->SetBinLabel(binNumber, labelName.c_str());
0101     hTMeanHisto->GetXaxis()->SetBinLabel(binNumber, labelName.c_str());
0102     hSigmaHisto->GetXaxis()->SetBinLabel(binNumber, labelName.c_str());
0103     hKFactorHisto->GetXaxis()->SetBinLabel(binNumber, labelName.c_str());
0104 
0105     //Define a distribution for each wheel,station and each superlayer type
0106     vector<int> Wh_St_SL;
0107     Wh_St_SL.push_back(wireId.wheel());
0108     Wh_St_SL.push_back(wireId.station());
0109     Wh_St_SL.push_back(wireId.superlayer());
0110     TH1D *hTTrigDistrib = theTTrigDistribMap[Wh_St_SL];
0111     TH1D *hTMeanDistrib = theTMeanDistribMap[Wh_St_SL];
0112     TH1D *hSigmaDistrib = theSigmaDistribMap[Wh_St_SL];
0113     TH1D *hKFactorDistrib = theKFactorDistribMap[Wh_St_SL];
0114     if (hTTrigDistrib == 0) {
0115       theFile->cd();
0116       TString name = getDistribName(wireId).c_str();
0117       hTTrigDistrib = new TH1D(name + "_TTrig", "TTrig calibrated from TB per superlayer", 500, 100, 600);
0118       hTMeanDistrib = new TH1D(name + "_TMean", "TMean calibrated from TB per superlayer", 500, 100, 600);
0119       hSigmaDistrib = new TH1D(name + "_Sigma", "Sigma calibrated from TB per superlayer", 50, 0, 50);
0120       hKFactorDistrib = new TH1D(name + "_KFactor", "KFactor calibrated from TB per superlayer", 200, -10.0, 10.0);
0121       theTTrigDistribMap[Wh_St_SL] = hTTrigDistrib;
0122       theTMeanDistribMap[Wh_St_SL] = hTMeanDistrib;
0123       theSigmaDistribMap[Wh_St_SL] = hSigmaDistrib;
0124       theKFactorDistribMap[Wh_St_SL] = hKFactorDistrib;
0125     }
0126     //Fill the distributions
0127     //    hTTrigDistrib->Fill(ttrig);
0128     hTTrigDistrib->Fill(ttcor);
0129     hTMeanDistrib->Fill(tmean);
0130     hSigmaDistrib->Fill(sigma);
0131     hKFactorDistrib->Fill(kFactor);
0132   }
0133 
0134   //Write histos in a .root file
0135   theFile->cd();
0136   for (map<pair<int, int>, TH1D *>::const_iterator lHisto = theTTrigHistoMap.begin(); lHisto != theTTrigHistoMap.end();
0137        ++lHisto) {
0138     (*lHisto).second->GetXaxis()->LabelsOption("v");
0139     (*lHisto).second->Write();
0140   }
0141   for (map<pair<int, int>, TH1D *>::const_iterator lHisto = theTMeanHistoMap.begin(); lHisto != theTMeanHistoMap.end();
0142        ++lHisto) {
0143     (*lHisto).second->GetXaxis()->LabelsOption("v");
0144     (*lHisto).second->Write();
0145   }
0146   for (map<pair<int, int>, TH1D *>::const_iterator lHisto = theSigmaHistoMap.begin(); lHisto != theSigmaHistoMap.end();
0147        ++lHisto) {
0148     (*lHisto).second->GetXaxis()->LabelsOption("v");
0149     (*lHisto).second->Write();
0150   }
0151   for (map<pair<int, int>, TH1D *>::const_iterator lHisto = theKFactorHistoMap.begin();
0152        lHisto != theKFactorHistoMap.end();
0153        ++lHisto) {
0154     (*lHisto).second->GetXaxis()->LabelsOption("v");
0155     (*lHisto).second->Write();
0156   }
0157   for (map<vector<int>, TH1D *>::const_iterator lDistrib = theTTrigDistribMap.begin();
0158        lDistrib != theTTrigDistribMap.end();
0159        ++lDistrib) {
0160     (*lDistrib).second->Write();
0161   }
0162   for (map<vector<int>, TH1D *>::const_iterator lDistrib = theTMeanDistribMap.begin();
0163        lDistrib != theTMeanDistribMap.end();
0164        ++lDistrib) {
0165     (*lDistrib).second->Write();
0166   }
0167   for (map<vector<int>, TH1D *>::const_iterator lDistrib = theSigmaDistribMap.begin();
0168        lDistrib != theSigmaDistribMap.end();
0169        ++lDistrib) {
0170     (*lDistrib).second->Write();
0171   }
0172   for (map<vector<int>, TH1D *>::const_iterator lDistrib = theKFactorDistribMap.begin();
0173        lDistrib != theKFactorDistribMap.end();
0174        ++lDistrib) {
0175     (*lDistrib).second->Write();
0176   }
0177 }
0178 
0179 string DTTTrigAnalyzer::getHistoName(const DTWireId &wId) const {
0180   string histoName;
0181   stringstream theStream;
0182   theStream << "Wheel" << wId.wheel() << "_SL" << wId.superlayer();
0183   theStream >> histoName;
0184   return histoName;
0185 }
0186 
0187 string DTTTrigAnalyzer::getDistribName(const DTWireId &wId) const {
0188   string histoName;
0189   stringstream theStream;
0190   theStream << "Wheel" << wId.wheel() << "_Station" << wId.station() << "_SL" << wId.superlayer();
0191   theStream >> histoName;
0192   return histoName;
0193 }