Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:33:38

0001 // -*- C++ -*-
0002 //
0003 // Package:     SiPixelPhase1RecHitsV
0004 // Class:       SiPixelPhase1RecHitsV
0005 //
0006 
0007 // Original Author: Marcel Schneider
0008 // Additional Authors: Alexander Morton - modifying code for validation use
0009 
0010 #include "Validation/SiPixelPhase1RecHitsV/interface/SiPixelPhase1RecHitsV.h"
0011 #include "FWCore/Framework/interface/MakerMacros.h"
0012 
0013 #include "FWCore/Framework/interface/ESHandle.h"
0014 #include "DataFormats/GeometryVector/interface/LocalPoint.h"
0015 
0016 #include "SimDataFormats/TrackingHit/interface/PSimHit.h"
0017 #include "SimDataFormats/TrackingHit/interface/PSimHitContainer.h"
0018 
0019 SiPixelPhase1RecHitsV::SiPixelPhase1RecHitsV(const edm::ParameterSet& iConfig)
0020     : SiPixelPhase1Base(iConfig),
0021       trackerHitAssociatorConfig_(iConfig, consumesCollector()),
0022       srcToken_(consumes<SiPixelRecHitCollection>(iConfig.getParameter<edm::InputTag>("src"))) {}
0023 
0024 void SiPixelPhase1RecHitsV::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0025   edm::Handle<SiPixelRecHitCollection> input;
0026   iEvent.getByToken(srcToken_, input);
0027   if (!input.isValid())
0028     return;
0029 
0030   TrackerHitAssociator associate(iEvent, trackerHitAssociatorConfig_);
0031 
0032   SiPixelRecHitCollection::const_iterator it;
0033   for (it = input->begin(); it != input->end(); ++it) {
0034     auto id = DetId(it->detId());
0035 
0036     for (SiPixelRecHit const& rechit : *it) {
0037       std::vector<PSimHit> associateSimHit;
0038       associateSimHit = associate.associateHit(rechit);
0039       std::vector<PSimHit>::const_iterator closestIt = associateSimHit.begin();
0040 
0041       LocalPoint lp = rechit.localPosition();
0042       float rechit_x = lp.x();
0043       float rechit_y = lp.y();
0044 
0045       LocalError lerr = rechit.localPositionError();
0046       float lerr_x = sqrt(lerr.xx());
0047       float lerr_y = sqrt(lerr.yy());
0048 
0049       // loop over associated sim hits and find the closest
0050       if (!associateSimHit.empty()) {
0051         float closestSimHit = 9999.9;
0052 
0053         for (std::vector<PSimHit>::const_iterator m = associateSimHit.begin(); m < associateSimHit.end(); m++) {
0054           float sim_x1((*m).entryPoint().x()), sim_x2((*m).exitPoint().x()), sim_xpos(0.5 * (sim_x1 + sim_x2));
0055           float sim_y1((*m).entryPoint().y()), sim_y2((*m).exitPoint().y()), sim_ypos(0.5 * (sim_y1 + sim_y2));
0056 
0057           float xres(std::abs(sim_xpos - rechit_x)), yres(std::abs(sim_ypos - rechit_y));
0058           float dist = std::sqrt(xres * xres + yres * yres);
0059 
0060           if (dist < closestSimHit) {
0061             closestSimHit = dist;
0062             closestIt = m;
0063           }
0064         }
0065       }
0066 
0067       // Sim Hit stuff
0068 
0069       if (!associateSimHit.empty()) {
0070         const PSimHit& simHit = *closestIt;
0071 
0072         int bunch = simHit.eventId().bunchCrossing();
0073         int event = simHit.eventId().event();
0074 
0075         float sim_x1(simHit.entryPoint().x()), sim_x2(simHit.exitPoint().x()), sim_xpos(0.5 * (sim_x1 + sim_x2));
0076         float sim_y1(simHit.entryPoint().y()), sim_y2(simHit.exitPoint().y()), sim_ypos(0.5 * (sim_y1 + sim_y2));
0077 
0078         float res_x = (rechit_x - sim_xpos) * 10000.0;
0079         float res_y = (rechit_y - sim_ypos) * 10000.0;
0080 
0081         float pull_x = (rechit_x - sim_xpos) / lerr_x;
0082         float pull_y = (rechit_y - sim_ypos) / lerr_y;
0083 
0084         // Now Plotting stuff
0085         if (bunch == 0)
0086           histo[IN_TIME_BUNCH].fill(bunch, id, &iEvent);
0087         if (bunch != 0)
0088           histo[OUT_TIME_BUNCH].fill(bunch, id, &iEvent);
0089 
0090         histo[NSIMHITS].fill(event, id, &iEvent);
0091 
0092         histo[RECHIT_X].fill(rechit_x, id, &iEvent);
0093         histo[RECHIT_Y].fill(rechit_y, id, &iEvent);
0094 
0095         histo[RES_X].fill(res_x, id, &iEvent);
0096         histo[RES_Y].fill(res_y, id, &iEvent);
0097 
0098         histo[ERROR_X].fill(lerr_x, id, &iEvent);
0099         histo[ERROR_Y].fill(lerr_y, id, &iEvent);
0100 
0101         histo[PULL_X].fill(pull_x, id, &iEvent);
0102         histo[PULL_Y].fill(pull_y, id, &iEvent);
0103       }
0104     }
0105   }
0106 }
0107 
0108 DEFINE_FWK_MODULE(SiPixelPhase1RecHitsV);