Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 /*
0002  * \file DTLocalLutTriggerTask.cc
0003  *
0004  * \author D. Fasanella - INFN Bologna
0005  *
0006  */
0007 
0008 #include "DQM/DTMonitorModule/src/DTLocalTriggerLutTask.h"
0009 
0010 // Framework
0011 #include "FWCore/Framework/interface/EventSetup.h"
0012 
0013 // DT trigger
0014 #include "DQM/DTMonitorModule/interface/DTTrigGeomUtils.h"
0015 
0016 // Geometry
0017 #include "DataFormats/GeometryVector/interface/Pi.h"
0018 #include "Geometry/DTGeometry/interface/DTGeometry.h"
0019 #include "Geometry/DTGeometry/interface/DTLayer.h"
0020 #include "Geometry/DTGeometry/interface/DTTopology.h"
0021 
0022 //Root
0023 #include "TH1.h"
0024 #include "TAxis.h"
0025 
0026 #include <sstream>
0027 #include <iostream>
0028 #include <cmath>
0029 
0030 using namespace edm;
0031 using namespace std;
0032 
0033 DTLocalTriggerLutTask::DTLocalTriggerLutTask(const edm::ParameterSet& ps)
0034     : muonGeomToken_(
0035           esConsumes<edm::Transition::BeginRun>(edm::ESInputTag("", ps.getUntrackedParameter<string>("geomLabel")))),
0036       trigGeomUtils(nullptr) {
0037   LogTrace("DTDQM|DTMonitorModule|DTLocalTriggerLutTask") << "[DTLocalTriggerLutTask]: Constructor" << endl;
0038 
0039   tm_TokenIn_ = consumes<L1MuDTChambPhContainer>(ps.getUntrackedParameter<InputTag>("inputTagTMin"));
0040   tm_TokenOut_ = consumes<L1MuDTChambPhContainer>(ps.getUntrackedParameter<InputTag>("inputTagTMout"));
0041   seg_Token_ = consumes<DTRecSegment4DCollection>(ps.getUntrackedParameter<InputTag>("inputTagSEG"));
0042 
0043   overUnderIn = ps.getUntrackedParameter<bool>("rebinOutFlowsInGraph");
0044   detailedAnalysis = ps.getUntrackedParameter<bool>("detailedAnalysis");
0045 
0046   if (detailedAnalysis) {
0047     nPhiBins = 401;
0048     rangePhi = 10.025;
0049     nPhibBins = 401;
0050     rangePhiB = 10.025;
0051   } else {
0052     nPhiBins = 51;
0053     rangePhi = 5.1;
0054     nPhibBins = 51;
0055     rangePhiB = 10.2;
0056   }
0057 
0058   baseFolder = "DT/03-LocalTrigger-TM/";
0059   parameters = ps;
0060 
0061   nEvents = 0;
0062   nLumis = 0;
0063 }
0064 
0065 DTLocalTriggerLutTask::~DTLocalTriggerLutTask() {
0066   LogTrace("DTDQM|DTMonitorModule|DTLocalTriggerLutTask")
0067       << "[DTLocalTriggerLutTask]: analyzed " << nEvents << " events" << endl;
0068   if (trigGeomUtils) {
0069     delete trigGeomUtils;
0070   }
0071 }
0072 
0073 void DTLocalTriggerLutTask::dqmBeginRun(const edm::Run& run, const edm::EventSetup& context) {
0074   // Get the geometry
0075   muonGeom = &context.getData(muonGeomToken_);
0076   trigGeomUtils = new DTTrigGeomUtils(muonGeom);
0077 }
0078 
0079 void DTLocalTriggerLutTask::bookHistos(DQMStore::IBooker& ibooker, DTChamberId chId) {
0080   stringstream wheel;
0081   wheel << chId.wheel();
0082   stringstream sector;
0083   sector << chId.sector();
0084   stringstream station;
0085   station << chId.station();
0086 
0087   ibooker.setCurrentFolder(topFolder() + "Wheel" + wheel.str() + "/Sector" + sector.str() + "/Station" + station.str() +
0088                            "/Segment");
0089 
0090   string chTag = "_W" + wheel.str() + "_Sec" + sector.str() + "_St" + station.str();
0091   std::map<std::string, MonitorElement*>& chambMap = chHistos[chId.rawId()];
0092 
0093   string hName = "TM_PhiResidualIn";
0094 
0095   chambMap[hName] = ibooker.book1D(hName + chTag,
0096                                    "Trigger local position In - Segment local position (correlated triggers)",
0097                                    nPhiBins,
0098                                    -rangePhi,
0099                                    rangePhi);
0100 
0101   hName = "TM_PhiResidualOut";
0102 
0103   chambMap[hName] = ibooker.book1D(hName + chTag,
0104                                    "Trigger local position Out - Segment local position (correlated triggers)",
0105                                    nPhiBins,
0106                                    -rangePhi,
0107                                    rangePhi);
0108 
0109   hName = "TM_PhibResidualIn";
0110 
0111   chambMap[hName] = ibooker.book1D(hName + chTag,
0112                                    "Trigger local direction In - Segment local direction (correlated triggers)",
0113                                    nPhibBins,
0114                                    -rangePhiB,
0115                                    rangePhiB);
0116 
0117   hName = "TM_PhibResidualOut";
0118 
0119   chambMap[hName] = ibooker.book1D(hName + chTag,
0120                                    "Trigger local direction Out - Segment local direction (correlated triggers)",
0121                                    nPhibBins,
0122                                    -rangePhiB,
0123                                    rangePhiB);
0124 
0125   if (detailedAnalysis) {
0126     hName = "TM_PhitkvsPhitrig";
0127 
0128     chambMap[hName] =
0129         ibooker.book2D(hName + chTag, "Local position: segment vs trigger", 100, -500., 500., 100, -500., 500.);
0130     hName = "TM_PhibtkvsPhibtrig";
0131 
0132     chambMap[hName] =
0133         ibooker.book2D(hName + chTag, "Local direction : segment vs trigger", 200, -40., 40., 200, -40., 40.);
0134     hName = "TM_PhibResidualvsTkPos";
0135 
0136     chambMap[hName] =
0137         ibooker.book2D(hName + chTag, "Local direction residual vs Segment Position", 100, -500., 500., 200, -10., 10.);
0138     hName = "TM_PhiResidualvsTkPos";
0139 
0140     chambMap[hName] =
0141         ibooker.book2D(hName + chTag, "Local Position residual vs Segment Position", 100, -500., 500., 200, -10., 10.);
0142   }
0143 }
0144 
0145 void DTLocalTriggerLutTask::bookHistograms(DQMStore::IBooker& ibooker,
0146                                            edm::Run const& run,
0147                                            edm::EventSetup const& context) {
0148   LogTrace("DTDQM|DTMonitorModule|DTLocalTriggerLutTask") << "[DTLocalTriggerLutTask]: bookHistograms" << endl;
0149 
0150   std::vector<const DTChamber*>::const_iterator chambIt = muonGeom->chambers().begin();
0151   std::vector<const DTChamber*>::const_iterator chambEnd = muonGeom->chambers().end();
0152 
0153   for (; chambIt != chambEnd; ++chambIt)
0154     bookHistos(ibooker, (*chambIt)->id());
0155 }
0156 
0157 void DTLocalTriggerLutTask::beginLuminosityBlock(const LuminosityBlock& lumiSeg, const EventSetup& context) {
0158   nLumis++;
0159   LogTrace("DTDQM|DTMonitorModule|DTLocalTriggerLutTask") << "[DTLocalTriggerLutTask]: Begin of LS transition" << endl;
0160 
0161   if (nLumis % parameters.getUntrackedParameter<int>("ResetCycle") == 0) {
0162     LogTrace("DTDQM|DTMonitorModule|DTLocalTriggerLutTask") << "[DTLocalTriggerLutTask]: Cleaning histos" << endl;
0163     map<uint32_t, map<string, MonitorElement*> >::const_iterator chambIt = chHistos.begin();
0164     map<uint32_t, map<string, MonitorElement*> >::const_iterator chambEnd = chHistos.end();
0165 
0166     for (; chambIt != chambEnd; ++chambIt) {
0167       map<string, MonitorElement*>::const_iterator histoIt = chambIt->second.begin();
0168       map<string, MonitorElement*>::const_iterator histoEnd = chambIt->second.end();
0169       for (; histoIt != histoEnd; ++histoIt) {
0170         histoIt->second->Reset();
0171       }
0172     }
0173   }
0174 }
0175 
0176 void DTLocalTriggerLutTask::analyze(const edm::Event& e, const edm::EventSetup& c) {
0177   nEvents++;
0178 
0179   edm::Handle<L1MuDTChambPhContainer> trigHandleIn;
0180   e.getByToken(tm_TokenIn_, trigHandleIn);
0181   edm::Handle<L1MuDTChambPhContainer> trigHandleOut;
0182   e.getByToken(tm_TokenOut_, trigHandleOut);
0183 
0184   vector<L1MuDTChambPhDigi> const* trigsIn = trigHandleIn->getContainer();
0185   searchTMBestIn(trigsIn);
0186   vector<L1MuDTChambPhDigi> const* trigsOut = trigHandleOut->getContainer();
0187   searchTMBestOut(trigsOut);
0188 
0189   Handle<DTRecSegment4DCollection> segments4D;
0190   e.getByToken(seg_Token_, segments4D);
0191   DTRecSegment4DCollection::id_iterator chamberId;
0192 
0193   // Preliminary loop finds best 4D Segment and high quality ones
0194   vector<const DTRecSegment4D*> best4DSegments;
0195 
0196   for (chamberId = segments4D->id_begin(); chamberId != segments4D->id_end(); ++chamberId) {
0197     DTRecSegment4DCollection::range rangeInCh = segments4D->get(*chamberId);
0198     DTRecSegment4DCollection::const_iterator trackIt = rangeInCh.first;
0199     DTRecSegment4DCollection::const_iterator trackEnd = rangeInCh.second;
0200 
0201     const DTRecSegment4D* tmpBest = nullptr;
0202     int tmpdof = 0;
0203     int dof = 0;
0204 
0205     for (; trackIt != trackEnd; ++trackIt) {
0206       if (trackIt->hasPhi()) {
0207         dof = trackIt->phiSegment()->degreesOfFreedom();
0208         if (dof > tmpdof) {
0209           tmpBest = &(*trackIt);
0210           tmpdof = dof;
0211         }
0212       }
0213     }
0214 
0215     if (tmpBest)
0216       best4DSegments.push_back(tmpBest);
0217   }
0218 
0219   vector<const DTRecSegment4D*>::const_iterator bestTrackIt = best4DSegments.begin();
0220   vector<const DTRecSegment4D*>::const_iterator bestTrackEnd = best4DSegments.end();
0221 
0222   for (; bestTrackIt != bestTrackEnd; ++bestTrackIt) {
0223     if ((*bestTrackIt)->hasPhi()) {
0224       DTChamberId chId = (*bestTrackIt)->chamberId();
0225       int nHitsPhi = (*bestTrackIt)->phiSegment()->degreesOfFreedom() + 2;
0226 
0227       int wheel = chId.wheel();
0228       int station = chId.station();
0229       int scsector = 0;
0230       float trackPosPhi, trackPosEta, trackDirPhi, trackDirEta;
0231       trigGeomUtils->computeSCCoordinates((*bestTrackIt), scsector, trackPosPhi, trackDirPhi, trackPosEta, trackDirEta);
0232 
0233       map<string, MonitorElement*>& chMap = chHistos[chId.rawId()];
0234       // In part
0235       if (trigQualBestIn[wheel + 3][station][scsector] > 3 &&  // residuals only for correlate triggers
0236           trigQualBestIn[wheel + 3][station][scsector] < 7 && nHitsPhi >= 7) {
0237         float trigPos = trigGeomUtils->trigPos(trigBestIn[wheel + 3][station][scsector]);
0238         float trigDir = trigGeomUtils->trigDir(trigBestIn[wheel + 3][station][scsector]);
0239         trigGeomUtils->trigToSeg(station, trigPos, trackDirPhi);
0240 
0241         double deltaPos = trigPos - trackPosPhi;
0242         deltaPos = overUnderIn ? max(min(deltaPos, rangePhi - 0.01), -rangePhi + 0.01) : deltaPos;
0243         double deltaDir = trigDir - trackDirPhi;
0244         deltaDir = overUnderIn ? max(min(deltaDir, rangePhiB - 0.01), -rangePhiB + 0.01) : deltaDir;
0245         chMap.find("TM_PhiResidualIn")->second->Fill(deltaPos);
0246         chMap.find("TM_PhibResidualIn")->second->Fill(deltaDir);
0247 
0248         if (detailedAnalysis) {
0249           chMap.find("TM_PhitkvsPhitrig")->second->Fill(trigPos, trackPosPhi);
0250           chMap.find("TM_PhibtkvsPhibtrig")->second->Fill(trigDir, trackDirPhi);
0251           chMap.find("TM_PhibResidualvsTkPos")->second->Fill(trackPosPhi, trigDir - trackDirPhi);
0252           chMap.find("TM_PhiResidualvsTkPos")->second->Fill(trackPosPhi, trigPos - trackPosPhi);
0253         }
0254       }
0255 
0256       // Out part
0257       if (trigQualBestOut[wheel + 3][station][scsector] > 3 &&  // residuals only for correlate triggers
0258           trigQualBestOut[wheel + 3][station][scsector] < 7 && nHitsPhi >= 7) {
0259         float trigPos = trigGeomUtils->trigPos(trigBestOut[wheel + 3][station][scsector]);
0260         float trigDir = trigGeomUtils->trigDir(trigBestOut[wheel + 3][station][scsector]);
0261         trigGeomUtils->trigToSeg(station, trigPos, trackDirPhi);
0262 
0263         double deltaPos = trigPos - trackPosPhi;
0264         deltaPos = overUnderIn ? max(min(deltaPos, rangePhi - 0.01), -rangePhi + 0.01) : deltaPos;
0265         double deltaDir = trigDir - trackDirPhi;
0266         deltaDir = overUnderIn ? max(min(deltaDir, rangePhiB - 0.01), -rangePhiB + 0.01) : deltaDir;
0267         chMap.find("TM_PhiResidualOut")->second->Fill(deltaPos);
0268         chMap.find("TM_PhibResidualOut")->second->Fill(deltaDir);
0269       }
0270     }
0271   }
0272 }
0273 
0274 void DTLocalTriggerLutTask::searchTMBestIn(std::vector<L1MuDTChambPhDigi> const* trigs) {
0275   string histoType;
0276   string histoTag;
0277 
0278   // define best quality trigger segment
0279   // start from 1 and zero is kept empty
0280   for (int st = 0; st <= 4; ++st)
0281     for (int wh = 0; wh <= 5; ++wh)
0282       for (int sec = 0; sec <= 12; ++sec)
0283         trigQualBestIn[wh][st][sec] = -1;
0284 
0285   vector<L1MuDTChambPhDigi>::const_iterator trigIt = trigs->begin();
0286   vector<L1MuDTChambPhDigi>::const_iterator trigEnd = trigs->end();
0287   for (; trigIt != trigEnd; ++trigIt) {
0288     int wh = trigIt->whNum();
0289     int sec = trigIt->scNum() + 1;  // DTTF -> DT sector range transform
0290     int st = trigIt->stNum();
0291     int qual = trigIt->code();
0292 
0293     if (qual > trigQualBestIn[wh + wheelArrayShift][st][sec] && qual < 7) {
0294       trigQualBestIn[wh + wheelArrayShift][st][sec] = qual;
0295       trigBestIn[wh + wheelArrayShift][st][sec] = &(*trigIt);
0296     }
0297   }
0298 }
0299 
0300 void DTLocalTriggerLutTask::searchTMBestOut(std::vector<L1MuDTChambPhDigi> const* trigs) {
0301   string histoType;
0302   string histoTag;
0303 
0304   // define best quality trigger segment
0305   //   // start from 1 and zero is kept empty
0306   for (int st = 0; st <= 4; ++st)
0307     for (int wh = 0; wh <= 5; ++wh)
0308       for (int sec = 0; sec <= 12; ++sec)
0309         trigQualBestOut[wh][st][sec] = -1;
0310 
0311   vector<L1MuDTChambPhDigi>::const_iterator trigIt = trigs->begin();
0312   vector<L1MuDTChambPhDigi>::const_iterator trigEnd = trigs->end();
0313   for (; trigIt != trigEnd; ++trigIt) {
0314     int wh = trigIt->whNum();
0315     int sec = trigIt->scNum() + 1;  // DTTF -> DT sector range transform
0316     int st = trigIt->stNum();
0317     int qual = trigIt->code();
0318 
0319     if (qual > trigQualBestOut[wh + wheelArrayShift][st][sec] && qual < 7) {
0320       trigQualBestOut[wh + wheelArrayShift][st][sec] = qual;
0321       trigBestOut[wh + wheelArrayShift][st][sec] = &(*trigIt);
0322     }
0323   }
0324 }
0325 
0326 // Local Variables:
0327 // show-trailing-whitespace: t
0328 // truncate-lines: t
0329 // End: