DTTPAnalyzer

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
/** \class DTTPAnalyzer
 *
 *  \author A. Vilela Pereira
 */

#include "FWCore/Framework/interface/one/EDAnalyzer.h"
#include "FWCore/Framework/interface/ConsumesCollector.h"
#include "FWCore/Utilities/interface/InputTag.h"
#include "FWCore/Framework/interface/ESHandle.h"
#include "DataFormats/DTDigi/interface/DTDigiCollection.h"
#include "Geometry/Records/interface/MuonGeometryRecord.h"

#include <string>
#include <map>

//class DTT0;
class DTLayerId;
class DTWireId;
class DTGeometry;
class DTTTrigBaseSync;
class TFile;

class DTTPAnalyzer : public edm::one::EDAnalyzer<> {
public:
  DTTPAnalyzer(const edm::ParameterSet&);
  ~DTTPAnalyzer() override;

  void analyze(const edm::Event&, const edm::EventSetup&) override;
  void endJob() override;

private:
  std::string getHistoName(const DTLayerId&);

  const bool subtractT0_;
  const edm::EDGetTokenT<DTDigiCollection> digiToken_;

  TFile* rootFile_;
  edm::ESHandle<DTGeometry> dtGeom_;
  const edm::ESGetToken<DTGeometry, MuonGeometryRecord> dtGeomToken_;
  std::unique_ptr<DTTTrigBaseSync> tTrigSync_;

  // Map of the t0 and sigma histos by layer
  std::map<DTWireId, int> nDigisPerWire_;
  std::map<DTWireId, double> sumWPerWire_;
  std::map<DTWireId, double> sumW2PerWire_;
};

#include "FWCore/Framework/interface/Event.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "FWCore/Framework/interface/MakerMacros.h"

#include "CalibMuon/DTDigiSync/interface/DTTTrigSyncFactory.h"
#include "CalibMuon/DTDigiSync/interface/DTTTrigBaseSync.h"

#include "Geometry/Records/interface/MuonNumberingRecord.h"
#include "Geometry/DTGeometry/interface/DTGeometry.h"

#include "DataFormats/GeometryVector/interface/GlobalPoint.h"

#include "CondFormats/DTObjects/interface/DTT0.h"

#include "TH1F.h"
#include "TFile.h"

DTTPAnalyzer::DTTPAnalyzer(const edm::ParameterSet& pset)
    : subtractT0_(pset.getParameter<bool>("subtractT0")),
      digiToken_(consumes<DTDigiCollection>(pset.getParameter<edm::InputTag>("digiLabel"))),
      dtGeomToken_(esConsumes()) {
  std::string rootFileName = pset.getUntrackedParameter<std::string>("rootFileName");
  rootFile_ = new TFile(rootFileName.c_str(), "RECREATE");
  rootFile_->cd();

  if (subtractT0_)
    tTrigSync_ = DTTTrigSyncFactory::get()->create(pset.getParameter<std::string>("tTrigMode"),
                                                   pset.getParameter<edm::ParameterSet>("tTrigModeConfig"),
                                                   consumesCollector());
}

DTTPAnalyzer::~DTTPAnalyzer() { rootFile_->Close(); }

void DTTPAnalyzer::analyze(const edm::Event& event, const edm::EventSetup& setup) {
  // Get the digis from the event
  const edm::Handle<DTDigiCollection>& digis = event.getHandle(digiToken_);

  if (subtractT0_) {
    tTrigSync_->setES(setup);
  }
  // Get the DT Geometry
  dtGeom_ = setup.getHandle(dtGeomToken_);

  // Iterate through all digi collections ordered by LayerId
  DTDigiCollection::DigiRangeIterator dtLayerIt;
  for (dtLayerIt = digis->begin(); dtLayerIt != digis->end(); ++dtLayerIt) {
    // Get the iterators over the digis associated with this LayerId
    const DTDigiCollection::Range& digiRange = (*dtLayerIt).second;

    // Get the layerId
    const DTLayerId layerId = (*dtLayerIt).first;  //FIXME: check to be in the right sector

    // Loop over all digis in the given layer
    for (DTDigiCollection::const_iterator digi = digiRange.first; digi != digiRange.second; ++digi) {
      const DTWireId wireId(layerId, (*digi).wire());

      double t0 = (*digi).countsTDC();

      //FIXME: Reject digis not coming from TP

      if (subtractT0_) {
        const DTLayer* layer = nullptr;  //fake
        const GlobalPoint glPt;          //fake
        double offset = tTrigSync_->offset(layer, wireId, glPt);
        t0 -= offset;
      }

      if (nDigisPerWire_.find(wireId) == nDigisPerWire_.end()) {
        nDigisPerWire_[wireId] = 0;
        sumWPerWire_[wireId] = 0.;
        sumW2PerWire_[wireId] = 0.;
      }

      ++nDigisPerWire_[wireId];
      sumWPerWire_[wireId] += t0;
      sumW2PerWire_[wireId] += t0 * t0;
    }
  }
}

void DTTPAnalyzer::endJob() {
  rootFile_->cd();
  std::map<DTLayerId, TH1F*> meanHistoMap;
  std::map<DTLayerId, TH1F*> sigmaHistoMap;
  for (std::map<DTWireId, int>::const_iterator wireIdIt = nDigisPerWire_.begin(); wireIdIt != nDigisPerWire_.end();
       ++wireIdIt) {
    DTWireId wireId((*wireIdIt).first);

    int nDigis = nDigisPerWire_[wireId];
    double sumW = sumWPerWire_[wireId];
    double sumW2 = sumW2PerWire_[wireId];

    double mean = sumW / nDigis;
    double rms = sumW2 / nDigis - mean * mean;
    rms = sqrt(rms);

    DTLayerId layerId = wireId.layerId();
    if (meanHistoMap.find(layerId) == meanHistoMap.end()) {
      std::string histoName = getHistoName(layerId);
      const int firstChannel = dtGeom_->layer(layerId)->specificTopology().firstChannel();
      const int nWires = dtGeom_->layer(layerId)->specificTopology().channels();
      TH1F* meanHistoTP = new TH1F((histoName + "_tpMean").c_str(),
                                   "mean from test pulses by channel",
                                   nWires,
                                   firstChannel,
                                   (firstChannel + nWires));
      TH1F* sigmaHistoTP = new TH1F((histoName + "_tpSigma").c_str(),
                                    "sigma from test pulses by channel",
                                    nWires,
                                    firstChannel,
                                    (firstChannel + nWires));
      meanHistoMap[layerId] = meanHistoTP;
      sigmaHistoMap[layerId] = sigmaHistoTP;
    }
    // Fill the histograms
    int nBin = meanHistoMap[layerId]->GetXaxis()->FindFixBin(wireId.wire());
    meanHistoMap[layerId]->SetBinContent(nBin, mean);
    sigmaHistoMap[layerId]->SetBinContent(nBin, rms);
  }

  for (std::map<DTLayerId, TH1F*>::const_iterator key = meanHistoMap.begin(); key != meanHistoMap.end(); ++key) {
    meanHistoMap[(*key).first]->Write();
    sigmaHistoMap[(*key).first]->Write();
  }
}

std::string DTTPAnalyzer::getHistoName(const DTLayerId& lId) {
  std::string histoName;
  std::stringstream theStream;
  theStream << "Ch_" << lId.wheel() << "_" << lId.station() << "_" << lId.sector() << "_SL" << lId.superlayer() << "_L"
            << lId.layer();
  theStream >> histoName;
  return histoName;
}

DEFINE_FWK_MODULE(DTTPAnalyzer);