Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #include "CommonTools/TriggerUtils/interface/GenericTriggerEventFlag.h"
0002 #include "DQM/SiPixelCommon/interface/SiPixelFolderOrganizer.h"
0003 #include "DQM/SiStripCommon/interface/SiStripFolderOrganizer.h"
0004 #include "DQM/SiStripCommon/interface/SiStripHistoId.h"
0005 #include "DQM/TrackerMonitorTrack/interface/MonitorTrackResiduals.h"
0006 #include "DQMServices/Core/interface/DQMStore.h"
0007 #include "DataFormats/DetId/interface/DetId.h"
0008 #include "DataFormats/TrackReco/interface/Track.h"
0009 #include "DataFormats/VertexReco/interface/Vertex.h"
0010 #include "FWCore/Framework/interface/Event.h"
0011 #include "FWCore/ServiceRegistry/interface/Service.h"
0012 #include "FWCore/Utilities/interface/Transition.h"
0013 
0014 template <TrackerType pixel_or_strip>
0015 MonitorTrackResidualsBase<pixel_or_strip>::MonitorTrackResidualsBase(const edm::ParameterSet &iConfig)
0016     : conf_(iConfig),
0017       tkDetMapToken_{esConsumes<TkDetMap, TrackerTopologyRcd, edm::Transition::BeginRun>()},
0018       trackerTopologyRunToken_{esConsumes<TrackerTopology, TrackerTopologyRcd, edm::Transition::BeginRun>()},
0019       trackerGeometryToken_{esConsumes<TrackerGeometry, TrackerDigiGeometryRecord, edm::Transition::BeginRun>()},
0020       trackerTopologyEventToken_{esConsumes<TrackerTopology, TrackerTopologyRcd>()},
0021       m_cacheID_(0),
0022       genTriggerEventFlag_(new GenericTriggerEventFlag(
0023           iConfig.getParameter<edm::ParameterSet>("genericTriggerEventPSet"), consumesCollector(), *this)),
0024       avalidator_(iConfig, consumesCollector()) {
0025   applyVertexCut_ = conf_.getUntrackedParameter<bool>("VertexCut", true);
0026   ModOn = conf_.getParameter<bool>("Mod_On");
0027   offlinePrimaryVerticesToken_ = consumes<reco::VertexCollection>(std::string("offlinePrimaryVertices"));
0028 }
0029 
0030 template <TrackerType pixel_or_strip>
0031 MonitorTrackResidualsBase<pixel_or_strip>::~MonitorTrackResidualsBase() {
0032   if (genTriggerEventFlag_)
0033     delete genTriggerEventFlag_;
0034 }
0035 
0036 template <TrackerType pixel_or_strip>
0037 void MonitorTrackResidualsBase<pixel_or_strip>::bookHistograms(DQMStore::IBooker &ibooker,
0038                                                                const edm::Run &run,
0039                                                                const edm::EventSetup &iSetup) {
0040   unsigned long long cacheID = iSetup.get<TrackerDigiGeometryRecord>().cacheIdentifier();
0041   if (m_cacheID_ != cacheID) {
0042     m_cacheID_ = cacheID;
0043     this->createMEs(ibooker, iSetup);
0044   }
0045   std::string topFolderName_ = "SiStrip";
0046   SiStripFolderOrganizer folder_organizer;
0047   folder_organizer.setSiStripFolderName(topFolderName_);
0048   const TkDetMap *tkDetMap = &iSetup.getData(tkDetMapToken_);
0049   tkhisto_ResidualsMean =
0050       std::make_unique<TkHistoMap>(tkDetMap, ibooker, topFolderName_, "TkHMap_ResidualsMean", 0.0, true);
0051 }
0052 
0053 template <TrackerType pixel_or_strip>
0054 void MonitorTrackResidualsBase<pixel_or_strip>::dqmBeginRun(edm::Run const &run, edm::EventSetup const &iSetup) {
0055   // Initialize the GenericTriggerEventFlag
0056   if (genTriggerEventFlag_->on())
0057     genTriggerEventFlag_->initRun(run, iSetup);
0058 }
0059 
0060 template <TrackerType pixel_or_strip>
0061 std::pair<std::string, int32_t> MonitorTrackResidualsBase<pixel_or_strip>::findSubdetAndLayer(
0062     uint32_t ModuleID, const TrackerTopology *tTopo) {
0063   std::string subdet = "";
0064   int32_t layer = 0;
0065   auto id = DetId(ModuleID);
0066   switch (id.subdetId()) {
0067     // Pixel Barrel, Endcap
0068     case 1:
0069       subdet = "BPIX";
0070       layer = tTopo->pxbLayer(id);
0071       break;
0072     case 2:
0073       subdet = "FPIX";
0074       layer = tTopo->pxfDisk(id) * (tTopo->pxfSide(ModuleID) == 1 ? -1 : +1);
0075       break;
0076     // Strip TIB, TID, TOB, TEC
0077     case 3:
0078       subdet = "TIB";
0079       layer = tTopo->tibLayer(id);
0080       break;
0081     case 4:
0082       subdet = "TID";
0083       layer = tTopo->tidWheel(id) * (tTopo->tidSide(ModuleID) == 1 ? -1 : +1);
0084       break;
0085     case 5:
0086       subdet = "TOB";
0087       layer = tTopo->tobLayer(id);
0088       break;
0089     case 6:
0090       subdet = "TEC";
0091       layer = tTopo->tecWheel(id) * (tTopo->tecSide(ModuleID) == 1 ? -1 : +1);
0092       break;
0093     default:
0094       // TODO: Fail loudly.
0095       subdet = "UNKNOWN";
0096       layer = 0;
0097   }
0098   return std::make_pair(subdet, layer);
0099 }
0100 
0101 template <TrackerType pixel_or_strip>
0102 void MonitorTrackResidualsBase<pixel_or_strip>::createMEs(DQMStore::IBooker &ibooker, const edm::EventSetup &iSetup) {
0103   // Retrieve tracker topology and geometry
0104   const TrackerTopology *const tTopo = &iSetup.getData(trackerTopologyRunToken_);
0105   const TrackerGeometry *TG = &iSetup.getData(trackerGeometryToken_);
0106 
0107   Parameters = conf_.getParameter<edm::ParameterSet>("TH1ResModules");
0108   int32_t i_residuals_Nbins = Parameters.getParameter<int32_t>("Nbinx");
0109   double d_residual_xmin = Parameters.getParameter<double>("xmin");
0110   double d_residual_xmax = Parameters.getParameter<double>("xmax");
0111   Parameters = conf_.getParameter<edm::ParameterSet>("TH1NormResModules");
0112   int32_t i_normres_Nbins = Parameters.getParameter<int32_t>("Nbinx");
0113   double d_normres_xmin = Parameters.getParameter<double>("xmin");
0114   double d_normres_xmax = Parameters.getParameter<double>("xmax");
0115 
0116   // use SistripHistoId for producing histogram id (and title)
0117   SiStripHistoId hidmanager;
0118 
0119   SiStripFolderOrganizer strip_organizer;
0120   auto pixel_organizer = SiPixelFolderOrganizer(false);
0121 
0122   // Collect list of modules from Tracker Geometry
0123   // book histo per each detector module
0124   auto ids = TG->detIds();  // or detUnitIds?
0125   for (DetId id : ids) {
0126     auto ModuleID = id.rawId();
0127     auto isPixel = id.subdetId() == 1 || id.subdetId() == 2;
0128     if (isPixel != (pixel_or_strip == TRACKERTYPE_PIXEL))
0129       continue;
0130 
0131     // Book module histogramms?
0132     if (ModOn) {
0133       switch (id.subdetId()) {
0134         case 1:
0135           pixel_organizer.setModuleFolder(ibooker, ModuleID, 0);
0136           break;
0137         case 2:
0138           pixel_organizer.setModuleFolder(ibooker, ModuleID, 0);
0139           break;
0140         default:
0141           strip_organizer.setDetectorFolder(ModuleID, tTopo);
0142       }
0143       {
0144         // this sounds strip specific but also works for pixel
0145         std::string hid = hidmanager.createHistoId("HitResidualsX", "det", ModuleID);
0146         std::string normhid = hidmanager.createHistoId("NormalizedHitResidualsX", "det", ModuleID);
0147         auto &histos = m_ModuleResiduals[std::make_pair("", ModuleID)];
0148         histos.x.base = ibooker.book1D(hid, hid, i_residuals_Nbins, d_residual_xmin, d_residual_xmax);
0149         histos.x.base->setAxisTitle("(x_{pred} - x_{rec})' [cm]");
0150         histos.x.normed = ibooker.book1D(normhid, normhid, i_normres_Nbins, d_normres_xmin, d_normres_xmax);
0151         histos.x.normed->setAxisTitle("(x_{pred} - x_{rec})'/#sigma");
0152       }
0153       {
0154         std::string hid = hidmanager.createHistoId("HitResidualsY", "det", ModuleID);
0155         std::string normhid = hidmanager.createHistoId("NormalizedHitResidualsY", "det", ModuleID);
0156         auto &histos = m_ModuleResiduals[std::make_pair("", ModuleID)];
0157         histos.y.base = ibooker.book1D(hid, hid, i_residuals_Nbins, d_residual_xmin, d_residual_xmax);
0158         histos.y.base->setAxisTitle("(y_{pred} - y_{rec})' [cm]");
0159         histos.y.normed = ibooker.book1D(normhid, normhid, i_normres_Nbins, d_normres_xmin, d_normres_xmax);
0160         histos.y.normed->setAxisTitle("(y_{pred} - y_{rec})'/#sigma");
0161       }
0162     }
0163 
0164     auto subdetandlayer = findSubdetAndLayer(ModuleID, tTopo);
0165     if (m_SubdetLayerResiduals.find(subdetandlayer) == m_SubdetLayerResiduals.end()) {
0166       // add new histograms
0167       auto &histos = m_SubdetLayerResiduals[subdetandlayer];
0168       switch (id.subdetId()) {
0169         // Pixel Barrel, Endcap
0170         // We can't use the folder organizer here (SiPixelActionExecutor.cc#1638
0171         // does the same)
0172         case 1:
0173           ibooker.setCurrentFolder("Pixel/Barrel");
0174           break;
0175         case 2:
0176           ibooker.setCurrentFolder("Pixel/Endcap");
0177           break;
0178         // All strip
0179         default:
0180           strip_organizer.setLayerFolder(ModuleID, tTopo, subdetandlayer.second);
0181       }
0182 
0183       auto isBarrel = subdetandlayer.first.find("B") != std::string::npos;
0184 
0185       auto xy = std::vector<std::pair<HistoPair &, const char *>>{std::make_pair(std::ref(histos.x), "X"),
0186                                                                   std::make_pair(std::ref(histos.y), "Y")};
0187       for (auto &histopair : xy) {
0188         // book histogramms on layer level, check for barrel/pixel only for
0189         // correct labeling
0190 
0191         // Skip the Y plots for strips.
0192         if (!isPixel && histopair.second[0] == 'Y')
0193           continue;
0194 
0195         std::string histoname = isPixel ? (  // Pixel name
0196                                               Form("HitResiduals%s_%s%d",
0197                                                    histopair.second,
0198                                                    isBarrel ? "L" : (subdetandlayer.second > 0 ? "Dp" : "Dm"),
0199                                                    std::abs(subdetandlayer.second)))
0200                                         : (Form("HitResiduals_%s__%s__%d",  // Strip TODO: We use a
0201                                                                             // legacy name.
0202                                                 subdetandlayer.first.c_str(),
0203                                                 isBarrel ? "Layer" : "wheel",
0204                                                 std::abs(subdetandlayer.second)));
0205 
0206         std::string histotitle = Form("HitResiduals %s on %s%s full %s %d",
0207                                       histopair.second,
0208                                       subdetandlayer.first.c_str(),
0209                                       isBarrel ? "" : (subdetandlayer.second > 0 ? "+" : "-"),
0210                                       isBarrel ? "Layer" : (isPixel ? "Disk" : "Wheel"),
0211                                       std::abs(subdetandlayer.second));
0212 
0213         std::string normhistoname = Form("Normalized%s", histoname.c_str());
0214         std::string normhistotitle = Form("Normalized%s", histotitle.c_str());
0215 
0216         // std::cout << "##### Booking: " << ibooker.pwd() << " title " <<
0217         // histoname << std::endl;
0218 
0219         histopair.first.base =
0220             ibooker.book1D(histoname.c_str(), histotitle.c_str(), i_residuals_Nbins, d_residual_xmin, d_residual_xmax);
0221         histopair.first.base->setAxisTitle("(x_{pred} - x_{rec})' [cm]");
0222 
0223         histopair.first.normed = ibooker.book1D(
0224             normhistoname.c_str(), normhistotitle.c_str(), i_normres_Nbins, d_normres_xmin, d_normres_xmax);
0225         histopair.first.normed->setAxisTitle("(x_{pred} - x_{rec})'/#sigma");
0226       }
0227     }
0228   }  // end loop over activeDets
0229 }
0230 
0231 template <TrackerType pixel_or_strip>
0232 void MonitorTrackResidualsBase<pixel_or_strip>::analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup) {
0233   auto vtracks = std::vector<TrackerValidationVariables::AVTrackStruct>();
0234   // Filter out events if Trigger Filtering is requested
0235   if (genTriggerEventFlag_->on() && !genTriggerEventFlag_->accept(iEvent, iSetup))
0236     return;
0237 
0238   edm::Handle<reco::VertexCollection> vertices;
0239   if (applyVertexCut_) {
0240     iEvent.getByToken(offlinePrimaryVerticesToken_, vertices);
0241     if (!vertices.isValid() || vertices->empty())
0242       return;
0243   }
0244 
0245   // Retrieve tracker topology from geometry
0246   const TrackerTopology *const tTopo = &iSetup.getData(trackerTopologyEventToken_);
0247 
0248   avalidator_.fillTrackQuantities(
0249       iEvent,
0250       iSetup,
0251       // tell the validator to only look at good tracks
0252       [&](const reco::Track &track) -> bool {
0253         return (!applyVertexCut_ ||
0254                 (track.pt() > 0.75 && abs(track.dxy(vertices->at(0).position())) < 5 * track.dxyError()));
0255       },
0256       vtracks);
0257 
0258   for (auto &track : vtracks) {
0259     for (auto &it : track.hits) {
0260       uint RawId = it.rawDetId;
0261 
0262       auto id = DetId(RawId);
0263       auto isPixel = id.subdetId() == 1 || id.subdetId() == 2;
0264       if (isPixel != (pixel_or_strip == TRACKERTYPE_PIXEL))
0265         continue;
0266 
0267       if (ModOn) {
0268         auto &mod_histos = m_ModuleResiduals[std::make_pair("", RawId)];
0269         mod_histos.x.base->Fill(it.resXprime);
0270         mod_histos.x.normed->Fill(it.resXprime / it.resXprimeErr);
0271         mod_histos.y.base->Fill(it.resYprime);
0272         mod_histos.y.normed->Fill(it.resYprime / it.resYprimeErr);
0273       }
0274 
0275       auto subdetandlayer = findSubdetAndLayer(RawId, tTopo);
0276       auto histos = m_SubdetLayerResiduals[subdetandlayer];
0277       // fill if its error is not zero
0278       if (it.resXprimeErr != 0 && histos.x.base) {
0279         histos.x.base->Fill(it.resXprime);
0280         histos.x.normed->Fill(it.resXprime / it.resXprimeErr);
0281         if (!isPixel)
0282           tkhisto_ResidualsMean->fill(RawId, it.resXprime);
0283       }
0284       if (it.resYprimeErr != 0 && histos.y.base) {
0285         histos.y.base->Fill(it.resYprime);
0286         histos.y.normed->Fill(it.resYprime / it.resYprimeErr);
0287       }
0288     }
0289   }
0290 }
0291 
0292 DEFINE_FWK_MODULE(MonitorTrackResiduals);
0293 DEFINE_FWK_MODULE(SiPixelMonitorTrackResiduals);