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 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217
/*
 * \file DTLocalTriggerSynchTask.cc
 *
 * \author C. Battilana - CIEMAT
 *
*/

#include "DQM/DTMonitorModule/src/DTLocalTriggerSynchTask.h"

// Framework
#include "FWCore/Framework/interface/EventSetup.h"
#include "FWCore/Framework/interface/ConsumesCollector.h"
#include "FWCore/Utilities/interface/InputTag.h"

// Geometry
#include "DataFormats/GeometryVector/interface/Pi.h"
#include "Geometry/DTGeometry/interface/DTGeometry.h"
#include "Geometry/DTGeometry/interface/DTLayer.h"
#include "Geometry/DTGeometry/interface/DTSuperLayer.h"
#include "Geometry/DTGeometry/interface/DTTopology.h"

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

// DT Digi
#include "DataFormats/DTDigi/interface/DTDigi.h"
#include "DataFormats/DTDigi/interface/DTDigiCollection.h"

//Root
#include "TH1.h"
#include "TAxis.h"

#include <sstream>
#include <iostream>
#include <fstream>

using namespace edm;
using namespace std;

DTLocalTriggerSynchTask::DTLocalTriggerSynchTask(const edm::ParameterSet& ps)
    : nevents(0),
      tTrigSync{DTTTrigSyncFactory::get()->create(ps.getParameter<std::string>("tTrigMode"),
                                                  ps.getParameter<edm::ParameterSet>("tTrigModeConfig"),
                                                  consumesCollector())},
      muonGeomToken_(esConsumes<edm::Transition::BeginRun>()) {
  edm::LogVerbatim("DTLocalTriggerSynchTask") << "[DTLocalTriggerSynchTask]: Constructor" << endl;
  tm_Token_ = consumes<L1MuDTChambPhContainer>(ps.getUntrackedParameter<edm::InputTag>("TMInputTag"));
  seg_Token_ = consumes<DTRecSegment4DCollection>(ps.getUntrackedParameter<edm::InputTag>("SEGInputTag"));

  bxTime = ps.getParameter<double>("bxTimeInterval");  // CB move this to static const or DB
  rangeInBX = ps.getParameter<bool>("rangeWithinBX");
  nBXLow = ps.getParameter<int>("nBXLow");
  nBXHigh = ps.getParameter<int>("nBXHigh");
  angleRange = ps.getParameter<double>("angleRange");
  minHitsPhi = ps.getParameter<int>("minHitsPhi");
  baseDirectory = ps.getParameter<string>("baseDir");
}

DTLocalTriggerSynchTask::~DTLocalTriggerSynchTask() {
  edm::LogVerbatim("DTLocalTriggerSynchTask") << "[DTLocalTriggerSynchTask]: analyzed " << nevents << " events" << endl;
}

void DTLocalTriggerSynchTask::bookHistograms(DQMStore::IBooker& ibooker,
                                             edm::Run const& iRun,
                                             edm::EventSetup const& context) {
  edm::LogVerbatim("DTLocalTriggerSynchTask") << "[DTLocalTriggerSynchTask]: Book Histograms" << endl;

  ibooker.setCurrentFolder(baseDir());
  ibooker.bookFloat("BXTimeSpacing")->Fill(bxTime);

  tTrigSync->setES(context);

  std::vector<const DTChamber*>::const_iterator chambIt = muonGeom->chambers().begin();
  std::vector<const DTChamber*>::const_iterator chambEnd = muonGeom->chambers().end();

  for (; chambIt != chambEnd; ++chambIt) {
    bookHistos(ibooker, (*chambIt)->id());
    triggerHistos[(*chambIt)->id().rawId()]["tTrig_SL1"]->Fill(tTrigSync->offset(DTWireId((*chambIt)->id(), 1, 1, 2)));
    triggerHistos[(*chambIt)->id().rawId()]["tTrig_SL3"]->Fill(tTrigSync->offset(DTWireId((*chambIt)->id(), 3, 1, 2)));
  }
}

void DTLocalTriggerSynchTask::dqmBeginRun(const Run& run, const EventSetup& context) {
  muonGeom = &context.getData(muonGeomToken_);
}

void DTLocalTriggerSynchTask::analyze(const edm::Event& event, const edm::EventSetup& context) {
  nevents++;

  for (int i = 0; i < 5; ++i) {
    for (int j = 0; j < 6; ++j) {
      for (int k = 0; k < 13; ++k) {
        phCodeBestTM[j][i][k] = -1;
      }
    }
  }

  // Get best TM triggers
  edm::Handle<L1MuDTChambPhContainer> l1DTTPGPh;
  event.getByToken(tm_Token_, l1DTTPGPh);
  vector<L1MuDTChambPhDigi> const* phTrigs = l1DTTPGPh->getContainer();

  vector<L1MuDTChambPhDigi>::const_iterator iph = phTrigs->begin();
  vector<L1MuDTChambPhDigi>::const_iterator iphe = phTrigs->end();
  for (; iph != iphe; ++iph) {
    int phwheel = iph->whNum();
    int phsec = iph->scNum() + 1;  // DTTF[0-11] -> DT[1-12] Sector Numbering
    int phst = iph->stNum();
    int phcode = iph->code();

    if (phcode > phCodeBestTM[phwheel + 3][phst][phsec] && phcode < 7) {
      phCodeBestTM[phwheel + 3][phst][phsec] = phcode;
    }
  }

  //Get best segments (highest number of phi hits)
  vector<const DTRecSegment4D*> bestSegments4D;
  Handle<DTRecSegment4DCollection> segments4D;
  event.getByToken(seg_Token_, segments4D);
  DTRecSegment4DCollection::const_iterator track;
  DTRecSegment4DCollection::id_iterator chambIdIt;

  for (chambIdIt = segments4D->id_begin(); chambIdIt != segments4D->id_end(); ++chambIdIt) {
    DTRecSegment4DCollection::range range = segments4D->get(*chambIdIt);
    const DTRecSegment4D* best = nullptr;
    int hitsBest = 0;
    int hits = 0;

    for (track = range.first; track != range.second; ++track) {
      if ((*track).hasPhi()) {
        hits = (*track).phiSegment()->degreesOfFreedom() + 2;
        if (hits > hitsBest) {
          best = &(*track);
          hitsBest = hits;
        }
      }
    }
    if (best) {
      bestSegments4D.push_back(best);
    }
  }

  // Filling histos
  vector<const DTRecSegment4D*>::const_iterator bestSegIt = bestSegments4D.begin();
  vector<const DTRecSegment4D*>::const_iterator bestSegEnd = bestSegments4D.end();
  for (; bestSegIt != bestSegEnd; ++bestSegIt) {
    float dir = atan((*bestSegIt)->localDirection().x() / (*bestSegIt)->localDirection().z()) * 180 /
                Geom::pi();  // CB cerca un modo migliore x farlo
    const DTRecSegment2D* seg2D = (*bestSegIt)->phiSegment();
    int nHitsPhi = seg2D->degreesOfFreedom() + 2;
    DTChamberId chambId = (*bestSegIt)->chamberId();
    map<string, MonitorElement*>& innerME = triggerHistos[chambId.rawId()];

    if (fabs(dir) < angleRange && nHitsPhi >= minHitsPhi && seg2D->ist0Valid()) {
      float t0seg = (*bestSegIt)->phiSegment()->t0();
      float tTrig = (tTrigSync->offset(DTWireId(chambId, 1, 1, 2)) + tTrigSync->offset(DTWireId(chambId, 3, 1, 2))) / 2;
      float time = tTrig + t0seg;
      float htime = rangeInBX ? time - int(time / bxTime) * bxTime : time - int(tTrig / bxTime) * bxTime;

      int wheel = chambId.wheel();
      int sector = chambId.sector();
      int station = chambId.station();
      int scsector = sector > 12 ? sector == 13 ? 4 : 10 : sector;

      int qualTM = phCodeBestTM[wheel + 3][station][scsector];

      if (fabs(t0seg) > 0.01) {
        innerME.find("SEG_TrackCrossingTime")->second->Fill(htime);
        if (qualTM >= 0)
          innerME.find("TM_TrackCrossingTimeAll")->second->Fill(htime);
        if (qualTM == 6)
          innerME.find("TM_TrackCrossingTimeHH")->second->Fill(htime);
      }
    }
  }
}

void DTLocalTriggerSynchTask::bookHistos(DQMStore::IBooker& ibooker, const DTChamberId& dtChId) {
  stringstream wheel;
  wheel << dtChId.wheel();
  stringstream station;
  station << dtChId.station();
  stringstream sector;
  sector << dtChId.sector();
  uint32_t chRawId = dtChId.rawId();

  ibooker.setCurrentFolder(baseDir() + "/Wheel" + wheel.str() + "/Sector" + sector.str() + "/Station" + station.str());

  std::vector<string> histoTags = {"SEG_TrackCrossingTime", "TM_TrackCrossingTimeAll", "TM_TrackCrossingTimeHH"};

  float min = rangeInBX ? 0 : nBXLow * bxTime;
  float max = rangeInBX ? bxTime : nBXHigh * bxTime;
  int nbins = static_cast<int>(ceil(rangeInBX ? bxTime : (nBXHigh - nBXLow) * bxTime));

  for (const auto& histoTag : histoTags) {
    string histoName =
        histoTag + (rangeInBX ? "InBX" : "") + "_W" + wheel.str() + "_Sec" + sector.str() + "_St" + station.str();
    edm::LogVerbatim("DTLocalTriggerSynchTask")
        << "[DTLocalTriggerSynchTask]: booking " << baseDir() + "/Wheel" << wheel.str() << "/Sector" << sector.str()
        << "/Station" << station.str() << "/" << histoName << endl;

    triggerHistos[chRawId][histoTag] = ibooker.book1D(histoName.c_str(), "Track time distribution", nbins, min, max);
  }

  string floatTag[2] = {"tTrig_SL1", "tTrig_SL3"};

  for (int iFloat = 0; iFloat < 2; ++iFloat) {
    string floatName = floatTag[iFloat] + "_W" + wheel.str() + "_Sec" + sector.str() + "_St" + station.str();
    triggerHistos[chRawId][floatTag[iFloat]] = ibooker.bookFloat(floatName);
  }
}

// Local Variables:
// show-trailing-whitespace: t
// truncate-lines: t
// End: