Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 /*
0002  * \file DTLocalTriggerSynchTask.cc
0003  *
0004  * \author C. Battilana - CIEMAT
0005  *
0006 */
0007 
0008 #include "DQM/DTMonitorModule/src/DTLocalTriggerSynchTask.h"
0009 
0010 // Framework
0011 #include "FWCore/Framework/interface/EventSetup.h"
0012 #include "FWCore/Framework/interface/ConsumesCollector.h"
0013 #include "FWCore/Utilities/interface/InputTag.h"
0014 
0015 // Geometry
0016 #include "DataFormats/GeometryVector/interface/Pi.h"
0017 #include "Geometry/DTGeometry/interface/DTGeometry.h"
0018 #include "Geometry/DTGeometry/interface/DTLayer.h"
0019 #include "Geometry/DTGeometry/interface/DTSuperLayer.h"
0020 #include "Geometry/DTGeometry/interface/DTTopology.h"
0021 
0022 // tTrigs
0023 #include "CalibMuon/DTDigiSync/interface/DTTTrigSyncFactory.h"
0024 #include "CalibMuon/DTDigiSync/interface/DTTTrigBaseSync.h"
0025 
0026 // DT Digi
0027 #include "DataFormats/DTDigi/interface/DTDigi.h"
0028 #include "DataFormats/DTDigi/interface/DTDigiCollection.h"
0029 
0030 //Root
0031 #include "TH1.h"
0032 #include "TAxis.h"
0033 
0034 #include <sstream>
0035 #include <iostream>
0036 #include <fstream>
0037 
0038 using namespace edm;
0039 using namespace std;
0040 
0041 DTLocalTriggerSynchTask::DTLocalTriggerSynchTask(const edm::ParameterSet& ps)
0042     : nevents(0),
0043       tTrigSync{DTTTrigSyncFactory::get()->create(ps.getParameter<std::string>("tTrigMode"),
0044                                                   ps.getParameter<edm::ParameterSet>("tTrigModeConfig"),
0045                                                   consumesCollector())},
0046       muonGeomToken_(esConsumes<edm::Transition::BeginRun>()) {
0047   edm::LogVerbatim("DTLocalTriggerSynchTask") << "[DTLocalTriggerSynchTask]: Constructor" << endl;
0048   tm_Token_ = consumes<L1MuDTChambPhContainer>(ps.getUntrackedParameter<edm::InputTag>("TMInputTag"));
0049   seg_Token_ = consumes<DTRecSegment4DCollection>(ps.getUntrackedParameter<edm::InputTag>("SEGInputTag"));
0050 
0051   bxTime = ps.getParameter<double>("bxTimeInterval");  // CB move this to static const or DB
0052   rangeInBX = ps.getParameter<bool>("rangeWithinBX");
0053   nBXLow = ps.getParameter<int>("nBXLow");
0054   nBXHigh = ps.getParameter<int>("nBXHigh");
0055   angleRange = ps.getParameter<double>("angleRange");
0056   minHitsPhi = ps.getParameter<int>("minHitsPhi");
0057   baseDirectory = ps.getParameter<string>("baseDir");
0058 }
0059 
0060 DTLocalTriggerSynchTask::~DTLocalTriggerSynchTask() {
0061   edm::LogVerbatim("DTLocalTriggerSynchTask") << "[DTLocalTriggerSynchTask]: analyzed " << nevents << " events" << endl;
0062 }
0063 
0064 void DTLocalTriggerSynchTask::bookHistograms(DQMStore::IBooker& ibooker,
0065                                              edm::Run const& iRun,
0066                                              edm::EventSetup const& context) {
0067   edm::LogVerbatim("DTLocalTriggerSynchTask") << "[DTLocalTriggerSynchTask]: Book Histograms" << endl;
0068 
0069   ibooker.setCurrentFolder(baseDir());
0070   ibooker.bookFloat("BXTimeSpacing")->Fill(bxTime);
0071 
0072   tTrigSync->setES(context);
0073 
0074   std::vector<const DTChamber*>::const_iterator chambIt = muonGeom->chambers().begin();
0075   std::vector<const DTChamber*>::const_iterator chambEnd = muonGeom->chambers().end();
0076 
0077   for (; chambIt != chambEnd; ++chambIt) {
0078     bookHistos(ibooker, (*chambIt)->id());
0079     triggerHistos[(*chambIt)->id().rawId()]["tTrig_SL1"]->Fill(tTrigSync->offset(DTWireId((*chambIt)->id(), 1, 1, 2)));
0080     triggerHistos[(*chambIt)->id().rawId()]["tTrig_SL3"]->Fill(tTrigSync->offset(DTWireId((*chambIt)->id(), 3, 1, 2)));
0081   }
0082 }
0083 
0084 void DTLocalTriggerSynchTask::dqmBeginRun(const Run& run, const EventSetup& context) {
0085   muonGeom = &context.getData(muonGeomToken_);
0086 }
0087 
0088 void DTLocalTriggerSynchTask::analyze(const edm::Event& event, const edm::EventSetup& context) {
0089   nevents++;
0090 
0091   for (int i = 0; i < 5; ++i) {
0092     for (int j = 0; j < 6; ++j) {
0093       for (int k = 0; k < 13; ++k) {
0094         phCodeBestTM[j][i][k] = -1;
0095       }
0096     }
0097   }
0098 
0099   // Get best TM triggers
0100   edm::Handle<L1MuDTChambPhContainer> l1DTTPGPh;
0101   event.getByToken(tm_Token_, l1DTTPGPh);
0102   vector<L1MuDTChambPhDigi> const* phTrigs = l1DTTPGPh->getContainer();
0103 
0104   vector<L1MuDTChambPhDigi>::const_iterator iph = phTrigs->begin();
0105   vector<L1MuDTChambPhDigi>::const_iterator iphe = phTrigs->end();
0106   for (; iph != iphe; ++iph) {
0107     int phwheel = iph->whNum();
0108     int phsec = iph->scNum() + 1;  // DTTF[0-11] -> DT[1-12] Sector Numbering
0109     int phst = iph->stNum();
0110     int phcode = iph->code();
0111 
0112     if (phcode > phCodeBestTM[phwheel + 3][phst][phsec] && phcode < 7) {
0113       phCodeBestTM[phwheel + 3][phst][phsec] = phcode;
0114     }
0115   }
0116 
0117   //Get best segments (highest number of phi hits)
0118   vector<const DTRecSegment4D*> bestSegments4D;
0119   Handle<DTRecSegment4DCollection> segments4D;
0120   event.getByToken(seg_Token_, segments4D);
0121   DTRecSegment4DCollection::const_iterator track;
0122   DTRecSegment4DCollection::id_iterator chambIdIt;
0123 
0124   for (chambIdIt = segments4D->id_begin(); chambIdIt != segments4D->id_end(); ++chambIdIt) {
0125     DTRecSegment4DCollection::range range = segments4D->get(*chambIdIt);
0126     const DTRecSegment4D* best = nullptr;
0127     int hitsBest = 0;
0128     int hits = 0;
0129 
0130     for (track = range.first; track != range.second; ++track) {
0131       if ((*track).hasPhi()) {
0132         hits = (*track).phiSegment()->degreesOfFreedom() + 2;
0133         if (hits > hitsBest) {
0134           best = &(*track);
0135           hitsBest = hits;
0136         }
0137       }
0138     }
0139     if (best) {
0140       bestSegments4D.push_back(best);
0141     }
0142   }
0143 
0144   // Filling histos
0145   vector<const DTRecSegment4D*>::const_iterator bestSegIt = bestSegments4D.begin();
0146   vector<const DTRecSegment4D*>::const_iterator bestSegEnd = bestSegments4D.end();
0147   for (; bestSegIt != bestSegEnd; ++bestSegIt) {
0148     float dir = atan((*bestSegIt)->localDirection().x() / (*bestSegIt)->localDirection().z()) * 180 /
0149                 Geom::pi();  // CB cerca un modo migliore x farlo
0150     const DTRecSegment2D* seg2D = (*bestSegIt)->phiSegment();
0151     int nHitsPhi = seg2D->degreesOfFreedom() + 2;
0152     DTChamberId chambId = (*bestSegIt)->chamberId();
0153     map<string, MonitorElement*>& innerME = triggerHistos[chambId.rawId()];
0154 
0155     if (fabs(dir) < angleRange && nHitsPhi >= minHitsPhi && seg2D->ist0Valid()) {
0156       float t0seg = (*bestSegIt)->phiSegment()->t0();
0157       float tTrig = (tTrigSync->offset(DTWireId(chambId, 1, 1, 2)) + tTrigSync->offset(DTWireId(chambId, 3, 1, 2))) / 2;
0158       float time = tTrig + t0seg;
0159       float htime = rangeInBX ? time - int(time / bxTime) * bxTime : time - int(tTrig / bxTime) * bxTime;
0160 
0161       int wheel = chambId.wheel();
0162       int sector = chambId.sector();
0163       int station = chambId.station();
0164       int scsector = sector > 12 ? sector == 13 ? 4 : 10 : sector;
0165 
0166       int qualTM = phCodeBestTM[wheel + 3][station][scsector];
0167 
0168       if (fabs(t0seg) > 0.01) {
0169         innerME.find("SEG_TrackCrossingTime")->second->Fill(htime);
0170         if (qualTM >= 0)
0171           innerME.find("TM_TrackCrossingTimeAll")->second->Fill(htime);
0172         if (qualTM == 6)
0173           innerME.find("TM_TrackCrossingTimeHH")->second->Fill(htime);
0174       }
0175     }
0176   }
0177 }
0178 
0179 void DTLocalTriggerSynchTask::bookHistos(DQMStore::IBooker& ibooker, const DTChamberId& dtChId) {
0180   stringstream wheel;
0181   wheel << dtChId.wheel();
0182   stringstream station;
0183   station << dtChId.station();
0184   stringstream sector;
0185   sector << dtChId.sector();
0186   uint32_t chRawId = dtChId.rawId();
0187 
0188   ibooker.setCurrentFolder(baseDir() + "/Wheel" + wheel.str() + "/Sector" + sector.str() + "/Station" + station.str());
0189 
0190   std::vector<string> histoTags = {"SEG_TrackCrossingTime", "TM_TrackCrossingTimeAll", "TM_TrackCrossingTimeHH"};
0191 
0192   float min = rangeInBX ? 0 : nBXLow * bxTime;
0193   float max = rangeInBX ? bxTime : nBXHigh * bxTime;
0194   int nbins = static_cast<int>(ceil(rangeInBX ? bxTime : (nBXHigh - nBXLow) * bxTime));
0195 
0196   for (const auto& histoTag : histoTags) {
0197     string histoName =
0198         histoTag + (rangeInBX ? "InBX" : "") + "_W" + wheel.str() + "_Sec" + sector.str() + "_St" + station.str();
0199     edm::LogVerbatim("DTLocalTriggerSynchTask")
0200         << "[DTLocalTriggerSynchTask]: booking " << baseDir() + "/Wheel" << wheel.str() << "/Sector" << sector.str()
0201         << "/Station" << station.str() << "/" << histoName << endl;
0202 
0203     triggerHistos[chRawId][histoTag] = ibooker.book1D(histoName.c_str(), "Track time distribution", nbins, min, max);
0204   }
0205 
0206   string floatTag[2] = {"tTrig_SL1", "tTrig_SL3"};
0207 
0208   for (int iFloat = 0; iFloat < 2; ++iFloat) {
0209     string floatName = floatTag[iFloat] + "_W" + wheel.str() + "_Sec" + sector.str() + "_St" + station.str();
0210     triggerHistos[chRawId][floatTag[iFloat]] = ibooker.bookFloat(floatName);
0211   }
0212 }
0213 
0214 // Local Variables:
0215 // show-trailing-whitespace: t
0216 // truncate-lines: t
0217 // End: