Line Code
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193

/*
 *  See header file for a description of this class.
 *
 *  \author S. Bolognesi - INFN Torino
 */

#include "DTTTrigAnalyzer.h"

#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "FWCore/Framework/interface/EventSetup.h"
#include "CondFormats/DTObjects/interface/DTTtrig.h"
#include "CondFormats/DataRecord/interface/DTTtrigRcd.h"
#include <iostream>
#include "TFile.h"
#include "TH1D.h"
#include "TString.h"

using namespace edm;
using namespace std;

DTTTrigAnalyzer::DTTTrigAnalyzer(const ParameterSet &pset) {
  // The root file which will contain the histos
  string rootFileName = pset.getUntrackedParameter<string>("rootFileName");
  theFile = new TFile(rootFileName.c_str(), "RECREATE");
  theFile->cd();
  //The k factor to compute ttrig
  //kfactor = pset.getUntrackedParameter<double>("kfactor",0);
  ttrigToken_ =
      esConsumes<edm::Transition::BeginRun>(edm::ESInputTag("", pset.getUntrackedParameter<string>("dbLabel")));
}

DTTTrigAnalyzer::~DTTTrigAnalyzer() { theFile->Close(); }

void DTTTrigAnalyzer::beginRun(const edm::Run &, const edm::EventSetup &eventSetup) {
  ESHandle<DTTtrig> tTrig = eventSetup.getHandle(ttrigToken_);
  tTrigMap = &*tTrig;
  cout << "[DTTTrigAnalyzer] TTrig version: " << tTrig->version() << endl;
}

void DTTTrigAnalyzer::endJob() {
  //   static const double convToNs = 25./32.;
  // Loop over DB entries
  for (DTTtrig::const_iterator ttrig = tTrigMap->begin(); ttrig != tTrigMap->end(); ++ttrig) {
    DTWireId wireId(
        (*ttrig).first.wheelId, (*ttrig).first.stationId, (*ttrig).first.sectorId, (*ttrig).first.slId, 0, 0);
    float tmean;
    float sigma;
    float kFactor;

    DetId detId(wireId.rawId());
    // ttrig and rms are ns
    tTrigMap->get(detId, tmean, sigma, kFactor, DTTimeUnits::ns);
    float ttcor = tmean + kFactor * sigma;
    cout << "Wire: " << wireId << endl
         << " Ttrig (ns): " << ttcor << endl
         << " tmean (ns): " << tmean << endl
         << " sigma (ns): " << sigma << endl
         << " kfactor: " << kFactor << endl;

    //Define an histo for each wheel and each superlayer type
    TH1D *hTTrigHisto = theTTrigHistoMap[make_pair(wireId.wheel(), wireId.superlayer())];
    TH1D *hTMeanHisto = theTMeanHistoMap[make_pair(wireId.wheel(), wireId.superlayer())];
    TH1D *hSigmaHisto = theSigmaHistoMap[make_pair(wireId.wheel(), wireId.superlayer())];
    TH1D *hKFactorHisto = theKFactorHistoMap[make_pair(wireId.wheel(), wireId.superlayer())];
    if (hTTrigHisto == 0) {
      theFile->cd();
      TString name = getHistoName(wireId).c_str();
      if (wireId.superlayer() != 2) {
        hTTrigHisto = new TH1D(name + "_TTrig", "TTrig calibrated from TB per superlayer", 50, 0, 50);
        hTMeanHisto = new TH1D(name + "_TMean", "TMean calibrated from TB per superlayer", 50, 0, 50);
        hSigmaHisto = new TH1D(name + "_Sigma", "Sigma calibrated from TB per superlayer", 50, 0, 50);
        hKFactorHisto = new TH1D(name + "_KFactor", "KFactor calibrated from TB per superlayer", 50, 0, 50);
      } else {
        hTTrigHisto = new TH1D(name + "_TTrig", "TTrig calibrated from TB per superlayer", 36, 0, 36);
        hTMeanHisto = new TH1D(name + "_TMean", "TMean calibrated from TB per superlayer", 36, 0, 36);
        hSigmaHisto = new TH1D(name + "_Sigma", "Sigma calibrated from TB per superlayer", 36, 0, 36);
        hKFactorHisto = new TH1D(name + "_KFactor", "KFactor calibrated from TB per superlayer", 36, 0, 36);
      }
      theTTrigHistoMap[make_pair(wireId.wheel(), wireId.superlayer())] = hTTrigHisto;
      theTMeanHistoMap[make_pair(wireId.wheel(), wireId.superlayer())] = hTMeanHisto;
      theSigmaHistoMap[make_pair(wireId.wheel(), wireId.superlayer())] = hSigmaHisto;
      theKFactorHistoMap[make_pair(wireId.wheel(), wireId.superlayer())] = hKFactorHisto;
    }

    //Fill the histos and set the bin label
    int binNumber = wireId.sector() + 12 * (wireId.station() - 1);
    //    hTTrigHisto->SetBinContent(binNumber,ttrig);
    hTTrigHisto->SetBinContent(binNumber, ttcor);
    hTMeanHisto->SetBinContent(binNumber, tmean);
    hSigmaHisto->SetBinContent(binNumber, sigma);
    hKFactorHisto->SetBinContent(binNumber, kFactor);
    string labelName;
    stringstream theStream;
    if (wireId.sector() == 1)
      theStream << "MB" << wireId.station() << "_Sec" << wireId.sector();
    else
      theStream << "Sec" << wireId.sector();
    theStream >> labelName;
    hTTrigHisto->GetXaxis()->SetBinLabel(binNumber, labelName.c_str());
    hTMeanHisto->GetXaxis()->SetBinLabel(binNumber, labelName.c_str());
    hSigmaHisto->GetXaxis()->SetBinLabel(binNumber, labelName.c_str());
    hKFactorHisto->GetXaxis()->SetBinLabel(binNumber, labelName.c_str());

    //Define a distribution for each wheel,station and each superlayer type
    vector<int> Wh_St_SL;
    Wh_St_SL.push_back(wireId.wheel());
    Wh_St_SL.push_back(wireId.station());
    Wh_St_SL.push_back(wireId.superlayer());
    TH1D *hTTrigDistrib = theTTrigDistribMap[Wh_St_SL];
    TH1D *hTMeanDistrib = theTMeanDistribMap[Wh_St_SL];
    TH1D *hSigmaDistrib = theSigmaDistribMap[Wh_St_SL];
    TH1D *hKFactorDistrib = theKFactorDistribMap[Wh_St_SL];
    if (hTTrigDistrib == 0) {
      theFile->cd();
      TString name = getDistribName(wireId).c_str();
      hTTrigDistrib = new TH1D(name + "_TTrig", "TTrig calibrated from TB per superlayer", 500, 100, 600);
      hTMeanDistrib = new TH1D(name + "_TMean", "TMean calibrated from TB per superlayer", 500, 100, 600);
      hSigmaDistrib = new TH1D(name + "_Sigma", "Sigma calibrated from TB per superlayer", 50, 0, 50);
      hKFactorDistrib = new TH1D(name + "_KFactor", "KFactor calibrated from TB per superlayer", 200, -10.0, 10.0);
      theTTrigDistribMap[Wh_St_SL] = hTTrigDistrib;
      theTMeanDistribMap[Wh_St_SL] = hTMeanDistrib;
      theSigmaDistribMap[Wh_St_SL] = hSigmaDistrib;
      theKFactorDistribMap[Wh_St_SL] = hKFactorDistrib;
    }
    //Fill the distributions
    //    hTTrigDistrib->Fill(ttrig);
    hTTrigDistrib->Fill(ttcor);
    hTMeanDistrib->Fill(tmean);
    hSigmaDistrib->Fill(sigma);
    hKFactorDistrib->Fill(kFactor);
  }

  //Write histos in a .root file
  theFile->cd();
  for (map<pair<int, int>, TH1D *>::const_iterator lHisto = theTTrigHistoMap.begin(); lHisto != theTTrigHistoMap.end();
       ++lHisto) {
    (*lHisto).second->GetXaxis()->LabelsOption("v");
    (*lHisto).second->Write();
  }
  for (map<pair<int, int>, TH1D *>::const_iterator lHisto = theTMeanHistoMap.begin(); lHisto != theTMeanHistoMap.end();
       ++lHisto) {
    (*lHisto).second->GetXaxis()->LabelsOption("v");
    (*lHisto).second->Write();
  }
  for (map<pair<int, int>, TH1D *>::const_iterator lHisto = theSigmaHistoMap.begin(); lHisto != theSigmaHistoMap.end();
       ++lHisto) {
    (*lHisto).second->GetXaxis()->LabelsOption("v");
    (*lHisto).second->Write();
  }
  for (map<pair<int, int>, TH1D *>::const_iterator lHisto = theKFactorHistoMap.begin();
       lHisto != theKFactorHistoMap.end();
       ++lHisto) {
    (*lHisto).second->GetXaxis()->LabelsOption("v");
    (*lHisto).second->Write();
  }
  for (map<vector<int>, TH1D *>::const_iterator lDistrib = theTTrigDistribMap.begin();
       lDistrib != theTTrigDistribMap.end();
       ++lDistrib) {
    (*lDistrib).second->Write();
  }
  for (map<vector<int>, TH1D *>::const_iterator lDistrib = theTMeanDistribMap.begin();
       lDistrib != theTMeanDistribMap.end();
       ++lDistrib) {
    (*lDistrib).second->Write();
  }
  for (map<vector<int>, TH1D *>::const_iterator lDistrib = theSigmaDistribMap.begin();
       lDistrib != theSigmaDistribMap.end();
       ++lDistrib) {
    (*lDistrib).second->Write();
  }
  for (map<vector<int>, TH1D *>::const_iterator lDistrib = theKFactorDistribMap.begin();
       lDistrib != theKFactorDistribMap.end();
       ++lDistrib) {
    (*lDistrib).second->Write();
  }
}

string DTTTrigAnalyzer::getHistoName(const DTWireId &wId) const {
  string histoName;
  stringstream theStream;
  theStream << "Wheel" << wId.wheel() << "_SL" << wId.superlayer();
  theStream >> histoName;
  return histoName;
}

string DTTTrigAnalyzer::getDistribName(const DTWireId &wId) const {
  string histoName;
  stringstream theStream;
  theStream << "Wheel" << wId.wheel() << "_Station" << wId.station() << "_SL" << wId.superlayer();
  theStream >> histoName;
  return histoName;
}