SiPixelPhase1TrackEfficiency

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 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562
// -*- C++ -*-
//
// Package:     SiPixelPhase1TrackEfficiency
// Class:       SiPixelPhase1TrackEfficiency
//

// Original Author: Marcel Schneider

#include "DQM/SiPixelPhase1Common/interface/SiPixelPhase1Base.h"
#include "DataFormats/SiPixelCluster/interface/SiPixelCluster.h"
#include "DataFormats/TrackReco/interface/TrackFwd.h"
#include "DataFormats/VertexReco/interface/VertexFwd.h"
#include "FWCore/Framework/interface/MakerMacros.h"

#include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
#include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
#include "Geometry/CommonDetUnit/interface/PixelGeomDetUnit.h"
#include "Geometry/CommonTopologies/interface/PixelTopology.h"
#include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
#include "DataFormats/SiStripDetId/interface/StripSubdetector.h"
#include "DataFormats/TrackReco/interface/Track.h"
#include "DataFormats/TrackerRecHit2D/interface/SiPixelRecHit.h"
#include "DataFormats/VertexReco/interface/Vertex.h"
#include "DataFormats/VertexReco/interface/VertexFwd.h"
#include "TrackingTools/PatternTools/interface/TrajTrackAssociation.h"
#include "TrackingTools/TrajectoryState/interface/FreeTrajectoryState.h"
#include "RecoTracker/MeasurementDet/interface/MeasurementTrackerEvent.h"
#include "TrackingTools/KalmanUpdators/interface/Chi2MeasurementEstimatorBase.h"
#include "TrackingTools/Records/interface/TrackingComponentsRecord.h"

#include "TrackingTools/KalmanUpdators/interface/Chi2MeasurementEstimator.h"
#include "TrackingTools/GeomPropagators/interface/Propagator.h"
#include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
#include "RecoTracker/MeasurementDet/interface/MeasurementTracker.h"
#include "RecoTracker/Record/interface/CkfComponentsRecord.h"
#include "TrackingTools/MeasurementDet/interface/LayerMeasurements.h"
#include "RecoLocalTracker/ClusterParameterEstimator/interface/PixelClusterParameterEstimator.h"

///commnet

namespace {

  class SiPixelPhase1TrackEfficiency final : public SiPixelPhase1Base {
    enum { VALID, MISSING, INACTIVE, EFFICIENCY, VERTICES };

  public:
    explicit SiPixelPhase1TrackEfficiency(const edm::ParameterSet& conf);
    void analyze(const edm::Event&, const edm::EventSetup&) override;

  private:
    edm::EDGetTokenT<edmNew::DetSetVector<SiPixelCluster>> clustersToken_;
    edm::EDGetTokenT<reco::TrackCollection> tracksToken_;
    edm::EDGetTokenT<reco::VertexCollection> vtxToken_;
    edm::EDGetTokenT<TrajTrackAssociationCollection> trajTrackCollectionToken_;
    edm::EDGetTokenT<MeasurementTrackerEvent> tracker_;  //new
    bool applyVertexCut_;

    edm::ESGetToken<TrackerTopology, TrackerTopologyRcd> trackerTopoToken_;
    edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> trackerGeomToken_;
    edm::ESGetToken<Propagator, TrackingComponentsRecord> propagatorToken_;
    edm::ESGetToken<Chi2MeasurementEstimatorBase, TrackingComponentsRecord> chi2MeasurementEstimatorBaseToken_;
    edm::ESGetToken<MeasurementTracker, CkfComponentsRecord> measurementTrackerToken_;
    edm::ESGetToken<PixelClusterParameterEstimator, TkPixelCPERecord> pixelClusterParameterEstimatorToken_;

    const TrackerTopology* trackerTopology_;
    const Propagator* trackerPropagator_;
    const MeasurementEstimator* chi2MeasurementEstimator_;
  };

  SiPixelPhase1TrackEfficiency::SiPixelPhase1TrackEfficiency(const edm::ParameterSet& iConfig)
      : SiPixelPhase1Base(iConfig)  //,
  {
    tracker_ = consumes<MeasurementTrackerEvent>(iConfig.getParameter<edm::InputTag>("tracker"));
    tracksToken_ = consumes<reco::TrackCollection>(iConfig.getParameter<edm::InputTag>("tracks"));
    vtxToken_ = consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("primaryvertices"));
    applyVertexCut_ = iConfig.getUntrackedParameter<bool>("VertexCut", true);
    trajTrackCollectionToken_ =
        consumes<TrajTrackAssociationCollection>(iConfig.getParameter<edm::InputTag>("trajectoryInput"));
    clustersToken_ = consumes<edmNew::DetSetVector<SiPixelCluster>>(iConfig.getParameter<edm::InputTag>("clusters"));

    trackerTopoToken_ = esConsumes<TrackerTopology, TrackerTopologyRcd>();
    trackerGeomToken_ = esConsumes<TrackerGeometry, TrackerDigiGeometryRecord>();
    propagatorToken_ = esConsumes<Propagator, TrackingComponentsRecord>(edm::ESInputTag("", "PropagatorWithMaterial"));
    chi2MeasurementEstimatorBaseToken_ =
        esConsumes<Chi2MeasurementEstimatorBase, TrackingComponentsRecord>(edm::ESInputTag("", "Chi2"));
    measurementTrackerToken_ = esConsumes<MeasurementTracker, CkfComponentsRecord>();
    pixelClusterParameterEstimatorToken_ =
        esConsumes<PixelClusterParameterEstimator, TkPixelCPERecord>(edm::ESInputTag("", "PixelCPEGeneric"));
  }

  void SiPixelPhase1TrackEfficiency::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
    if (!checktrigger(iEvent, iSetup, DCS))
      return;

    // get geometry
    edm::ESHandle<TrackerGeometry> tracker = iSetup.getHandle(trackerGeomToken_);
    assert(tracker.isValid());

    // get primary vertex
    edm::Handle<reco::VertexCollection> vertices;
    iEvent.getByToken(vtxToken_, vertices);

    // TrackerTopology for module informations
    edm::ESHandle<TrackerTopology> trackerTopologyHandle = iSetup.getHandle(trackerTopoToken_);
    trackerTopology_ = trackerTopologyHandle.product();

    // Tracker propagator for propagating tracks to other layers
    edm::ESHandle<Propagator> propagatorHandle = iSetup.getHandle(propagatorToken_);
    std::unique_ptr<Propagator> propagatorUniquePtr(propagatorHandle.product()->clone());
    trackerPropagator_ = propagatorUniquePtr.get();
    const_cast<Propagator*>(trackerPropagator_)->setPropagationDirection(oppositeToMomentum);

    // Measurement estimator
    edm::ESHandle<Chi2MeasurementEstimatorBase> chi2MeasurementEstimatorHandle =
        iSetup.getHandle(chi2MeasurementEstimatorBaseToken_);
    chi2MeasurementEstimator_ = chi2MeasurementEstimatorHandle.product();

    //Tracker
    edm::Handle<MeasurementTrackerEvent> trackerMeas;
    iEvent.getByToken(tracker_, trackerMeas);

    edm::ESHandle<MeasurementTracker> measurementTrackerHandle = iSetup.getHandle(measurementTrackerToken_);

    //vertices
    if (!vertices.isValid())
      return;
    histo[VERTICES].fill(vertices->size(), DetId(0), &iEvent);
    if (applyVertexCut_ && vertices->empty())
      return;

    // should be used for weird cuts
    //const auto primaryVertex = vertices->at(0);

    // get the map
    edm::Handle<reco::TrackCollection> tracks;
    iEvent.getByToken(tracksToken_, tracks);
    if (!tracks.isValid())
      return;

    //new
    edm::Handle<TrajTrackAssociationCollection> trajTrackCollectionHandle;
    iEvent.getByToken(trajTrackCollectionToken_, trajTrackCollectionHandle);
    if (!trajTrackCollectionHandle.isValid())
      return;

    //Access Pixel Clusters
    edm::Handle<edmNew::DetSetVector<SiPixelCluster>> siPixelClusters;
    iEvent.getByToken(clustersToken_, siPixelClusters);
    if (!siPixelClusters.isValid())
      return;
    //

    edm::ESHandle<PixelClusterParameterEstimator> cpEstimator = iSetup.getHandle(pixelClusterParameterEstimatorToken_);
    if (!cpEstimator.isValid())
      return;

    const PixelClusterParameterEstimator& cpe(*cpEstimator);
    const TrackerGeometry* tkgeom = &(*tracker);

    //////////////////////////////////////////////////////////////////////////////////////////

    // Hp cut
    static constexpr int TRACK_QUALITY_HIGH_PURITY_BIT = 2;
    static constexpr int TRACK_QUALITY_HIGH_PURITY_MASK = 1 << TRACK_QUALITY_HIGH_PURITY_BIT;

    // Pt cut
    static constexpr float TRACK_PT_CUT_VAL = 1.0f;

    // Nstrip cut
    static constexpr int TRACK_NSTRIP_CUT_VAL = 10;

    //D0
    static constexpr std::array<float, 4> TRACK_D0_CUT_BARREL_VAL = {{0.01f, 0.02f, 0.02f, 0.02f}};
    static constexpr float TRACK_D0_CUT_FORWARD_VAL = 0.05f;

    //Dz
    static constexpr float TRACK_DZ_CUT_BARREL_VAL = 0.01f;
    static constexpr float TRACK_DZ_CUT_FORWARD_VAL = 0.5f;

    TrajectoryStateOnSurface tsosPXB2;
    bool valid_layerFrom = false;

    const GeometricSearchTracker* gst_ = trackerMeas->geometricSearchTracker();
    const auto* pxbLayer1_ = gst_->pixelBarrelLayers().front();
    const LayerMeasurements* theLayerMeasurements_ = new LayerMeasurements(*measurementTrackerHandle, *trackerMeas);

    std::vector<TrajectoryMeasurement> expTrajMeasurements;
    std::vector<std::pair<int, bool[3]>> eff_pxb1_vector;

    for (const auto& pair : *trajTrackCollectionHandle) {
      const edm::Ref<std::vector<Trajectory>> traj = pair.key;
      const reco::TrackRef track = pair.val;

      expTrajMeasurements.clear();
      eff_pxb1_vector.clear();
      //this cut is needed to be consisten with residuals calculation
      if (applyVertexCut_ &&
          (track->pt() < 0.75 || std::abs(track->dxy(vertices->at(0).position())) > 5 * track->dxyError()))
        continue;

      bool isBpixtrack = false;
      bool isFpixtrack = false;
      int nStripHits = 0;
      int nBpixL1Hits = 0;
      int nBpixL2Hits = 0;
      int nBpixL3Hits = 0;
      int nBpixL4Hits = 0;
      int nFpixD1Hits = 0;
      int nFpixD2Hits = 0;
      int nFpixD3Hits = 0;
      bool passcuts = true;

      // first, look at the full track to see whether it is good
      // auto const & trajParams = track.extra()->trajParams();

      auto hb = track->recHitsBegin();
      for (unsigned int h = 0; h < track->recHitsSize(); h++) {
        auto hit = *(hb + h);
        if (!hit->isValid())
          continue;

        DetId id = hit->geographicalId();
        uint32_t subdetid = (id.subdetId());

        //Check the location of valid hit
        if (subdetid == PixelSubdetector::PixelBarrel && hit->isValid()) {
          isBpixtrack = true;
          if (trackerTopology_->pxbLayer(id) == 1)
            nBpixL1Hits++;
          else if (trackerTopology_->pxbLayer(id) == 2)
            nBpixL2Hits++;
          else if (trackerTopology_->pxbLayer(id) == 3)
            nBpixL3Hits++;
          else if (trackerTopology_->pxbLayer(id) == 4)
            nBpixL4Hits++;
        } else if (subdetid == PixelSubdetector::PixelEndcap && hit->isValid()) {
          isFpixtrack = true;
          if (trackerTopology_->pxfDisk(id) == 1)
            nFpixD1Hits++;
          else if (trackerTopology_->pxfDisk(id) == 2)
            nFpixD2Hits++;
          else if (trackerTopology_->pxfDisk(id) == 3)
            nFpixD3Hits++;
        }

        // count strip hits
        if (subdetid == StripSubdetector::TIB || subdetid == StripSubdetector::TOB ||
            subdetid == StripSubdetector::TID || subdetid == StripSubdetector::TEC)
          nStripHits++;
      }

      if (!isBpixtrack && !isFpixtrack)
        continue;

      // Hp cut
      if (!((track->qualityMask() & TRACK_QUALITY_HIGH_PURITY_MASK) >> TRACK_QUALITY_HIGH_PURITY_BIT))
        passcuts = false;

      // Pt cut
      if (!(TRACK_PT_CUT_VAL < track->pt()))
        passcuts = false;

      // Nstrip cut
      if (!(TRACK_NSTRIP_CUT_VAL < nStripHits))
        passcuts = false;

      // then, look at each hit
      for (unsigned int h = 0; h < track->recHitsSize(); h++) {
        bool passcuts_hit = true;
        auto hit = *(hb + h);

        DetId id = hit->geographicalId();
        uint32_t subdetid = (id.subdetId());
        if (subdetid != PixelSubdetector::PixelBarrel && subdetid != PixelSubdetector::PixelEndcap)
          continue;

        bool isHitValid = hit->getType() == TrackingRecHit::valid;
        bool isHitMissing = hit->getType() == TrackingRecHit::missing;
        bool isHitInactive = hit->getType() == TrackingRecHit::inactive;

        //D0
        if (subdetid == PixelSubdetector::PixelBarrel) {
          if (!((std::abs(track->dxy(vertices->at(0).position())) * -1.0) <
                TRACK_D0_CUT_BARREL_VAL[trackerTopology_->pxbLayer(id) - 1]))
            passcuts_hit = false;
        } else if (subdetid == PixelSubdetector::PixelEndcap) {
          if (!((std::abs(track->dxy(vertices->at(0).position())) * -1.0) < TRACK_D0_CUT_FORWARD_VAL))
            passcuts_hit = false;
        }

        //Dz
        if (subdetid == PixelSubdetector::PixelBarrel) {
          if (!(std::abs(track->dz(vertices->at(0).position())) < TRACK_DZ_CUT_BARREL_VAL))
            passcuts_hit = false;
        } else if (subdetid == PixelSubdetector::PixelEndcap) {
          if (!(std::abs(track->dz(vertices->at(0).position())) < TRACK_DZ_CUT_FORWARD_VAL))
            passcuts_hit = false;
        }

        // Pixhit cut
        if (subdetid == PixelSubdetector::PixelBarrel) {
          if (trackerTopology_->pxbLayer(id) == 1) {
            if (!((nBpixL2Hits > 0 && nBpixL3Hits > 0 && nBpixL4Hits > 0) ||
                  (nBpixL2Hits > 0 && nBpixL3Hits > 0 && nFpixD1Hits > 0) ||
                  (nBpixL2Hits > 0 && nFpixD1Hits > 0 && nFpixD2Hits > 0) ||
                  (nFpixD1Hits > 0 && nFpixD2Hits > 0 && nFpixD3Hits > 0)))
              passcuts_hit = false;
          } else if (trackerTopology_->pxbLayer(id) == 2) {
            if (!((nBpixL1Hits > 0 && nBpixL3Hits > 0 && nBpixL4Hits > 0) ||
                  (nBpixL1Hits > 0 && nBpixL3Hits > 0 && nFpixD1Hits > 0) ||
                  (nBpixL1Hits > 0 && nFpixD1Hits > 0 && nFpixD2Hits > 0)))
              passcuts_hit = false;
          } else if (trackerTopology_->pxbLayer(id) == 3) {
            if (!((nBpixL1Hits > 0 && nBpixL2Hits > 0 && nBpixL4Hits > 0) ||
                  (nBpixL1Hits > 0 && nBpixL2Hits > 0 && nFpixD1Hits > 0)))
              passcuts_hit = false;
          } else if (trackerTopology_->pxbLayer(id) == 4)
            if (!((nBpixL1Hits > 0 && nBpixL2Hits > 0 && nBpixL3Hits > 0)))
              passcuts_hit = false;
        } else if (subdetid == PixelSubdetector::PixelEndcap) {
          if (trackerTopology_->pxfDisk(id) == 1) {
            if (!((nBpixL1Hits > 0 && nBpixL2Hits > 0 && nBpixL3Hits > 0) ||
                  (nBpixL1Hits > 0 && nBpixL2Hits > 0 && nFpixD2Hits > 0) ||
                  (nBpixL1Hits > 0 && nFpixD2Hits > 0 && nFpixD3Hits > 0)))
              passcuts_hit = false;
          } else if (trackerTopology_->pxfDisk(id) == 2) {
            if (!((nBpixL1Hits > 0 && nBpixL2Hits > 0 && nFpixD1Hits > 0) ||
                  (nBpixL1Hits > 0 && nFpixD1Hits > 0 && nFpixD3Hits > 0)))
              passcuts_hit = false;
          } else if (trackerTopology_->pxfDisk(id) == 3) {
            if (!((nBpixL1Hits > 0 && nFpixD1Hits > 0 && nFpixD2Hits > 0)))
              passcuts_hit = false;
          }
        }
        /* 
      //Fiducial Cut - will work on it later
      const SiPixelRecHit* pixhit = dynamic_cast<const SiPixelRecHit*>(hit);
      const PixelGeomDetUnit* geomdetunit = dynamic_cast<const PixelGeomDetUnit*> ( tracker->idToDet(id) );
      const PixelTopology& topol = geomdetunit->specificTopology();
      
      LocalPoint lp;
      if (pixhit) {
	lp = pixhit->localPosition();
      }
      
      MeasurementPoint mp = topol.measurementPosition(lp);
      int row = (int) mp.x() % 80;
      int col = (int) mp.y() % 52;
      
      int centerrow = 40;
      int centercol = 26;
      
      if (!((col < (centercol + 10)) && (col > (centercol - 10)) && (row < (centerrow + 10)) && (row > (centerrow -10 )))) passcuts_hit = false;
      */

        if (passcuts_hit && passcuts) {
          if (!(subdetid == PixelSubdetector::PixelBarrel && trackerTopology_->pxbLayer(id) == 1)) {
            if (isHitValid) {
              histo[VALID].fill(id, &iEvent);
              histo[EFFICIENCY].fill(1, id, &iEvent);
            }
            if (isHitMissing) {
              histo[MISSING].fill(id, &iEvent);
              histo[EFFICIENCY].fill(0, id, &iEvent);
            }
            if (isHitInactive) {
              histo[INACTIVE].fill(id, &iEvent);
            }
          }
        }
      }

      ///////////////////////////////////////////////layer 1 specific here/////////////////////////////////////////////////////////////////////
      valid_layerFrom = false;

      //propagation only from PXB2 and PXD1, more cuts later
      for (const auto& tm : traj->measurements()) {
        if (tm.recHit().get() && tm.recHitR().isValid()) {
          DetId where = tm.recHitR().geographicalId();
          int source_det = where.subdetId();

          if (source_det == PixelSubdetector::SubDetector::PixelBarrel) {
            int source_layer = trackerTopology_->pxbLayer(where);
            if (source_layer == 2) {
              if (tm.updatedState().isValid()) {
                tsosPXB2 = tm.updatedState();
                valid_layerFrom = true;
              }
            }
          }

          if (source_det == PixelSubdetector::SubDetector::PixelEndcap) {
            int source_layer = trackerTopology_->pxfDisk(where);
            if (source_layer == 1) {
              if (tm.updatedState().isValid()) {
                tsosPXB2 = tm.updatedState();
                valid_layerFrom = true;
              }
            }
          }
        }
      }  //uodated tsosPXB2 here

      if (!valid_layerFrom)
        continue;
      if (!tsosPXB2.isValid())
        continue;

      //propagation A: Calculate the efficiency by the distance to the closest cluster
      expTrajMeasurements =
          theLayerMeasurements_->measurements(*pxbLayer1_, tsosPXB2, *trackerPropagator_, *chi2MeasurementEstimator_);
      auto compDets = pxbLayer1_->compatibleDets(tsosPXB2, *trackerPropagator_, *chi2MeasurementEstimator_);
      std::pair<int, bool[3]> eff_map;

      //Fiducial Cut, only calculate the efficiency of the central pixels
      for (uint p = 0; p < expTrajMeasurements.size(); p++) {
        bool valid = false;
        bool missing = false;
        bool passcuts_hit = true;
        TrajectoryMeasurement pxb1TM(expTrajMeasurements[p]);
        const auto& pxb1Hit = pxb1TM.recHit();
        bool inactive = (pxb1Hit->getType() == TrackingRecHit::inactive);
        int detidHit = pxb1Hit->geographicalId();
        if (detidHit == 0)
          continue;

        const SiPixelRecHit* pixhit = dynamic_cast<const SiPixelRecHit*>(pxb1Hit->hit());
        const PixelGeomDetUnit* geomdetunit = dynamic_cast<const PixelGeomDetUnit*>(tracker->idToDet(detidHit));
        const PixelTopology& topol = geomdetunit->specificTopology();

        if (!pixhit)
          continue;

        LocalPoint lp = pixhit->localPosition();
        MeasurementPoint mp = topol.measurementPosition(lp);
        const int nRows = topol.rowsperroc();
        const int nColumns = topol.colsperroc();
        int row = (int)mp.x() % nRows;
        int col = (int)mp.y() % nColumns;

        int centerrow = nRows / 2;
        int centercol = nColumns / 2;

        if (!((col < (centercol + 10)) && (col > (centercol - 10)) && (row < (centerrow + 10)) &&
              (row > (centerrow - 10))))
          passcuts_hit = false;

        //Access the distance to the closest cluster
        for (const auto& detAndState : compDets) {
          const auto& pXb1_lpos = detAndState.second.localPosition();
          if (pxb1Hit->geographicalId().rawId() != detAndState.first->geographicalId().rawId())
            continue;
          int detid = detAndState.first->geographicalId().rawId();

          for (edmNew::DetSetVector<SiPixelCluster>::const_iterator iter_cl = siPixelClusters->begin();
               iter_cl != siPixelClusters->end();
               iter_cl++) {
            DetId detId(iter_cl->id());
            float minD[2], minDist = 10000;
            minD[0] = minD[1] = 10000.;
            if (detId.rawId() != detAndState.first->geographicalId().rawId())
              continue;

            const PixelGeomDetUnit* pixdet = (const PixelGeomDetUnit*)tkgeom->idToDetUnit(detId);
            edmNew::DetSet<SiPixelCluster>::const_iterator itCluster = iter_cl->begin();
            if (passcuts_hit) {
              for (; itCluster != iter_cl->end(); ++itCluster) {
                LocalPoint lp(itCluster->x(), itCluster->y(), 0.);
                PixelClusterParameterEstimator::ReturnType params = cpe.getParameters(*itCluster, *pixdet);
                lp = std::get<0>(params);

                float Xdist = abs(lp.x() - pXb1_lpos.x());
                float Ydist = abs(lp.y() - pXb1_lpos.y());
                float dist = sqrt(Xdist * Xdist + Ydist * Ydist);
                if (dist < minDist) {
                  minDist = dist;
                  minD[0] = Xdist;
                  minD[1] = Ydist;
                }
              }

              if ((minD[0] < 0.02) && (minD[1] < 0.02)) {
                valid = true;
                missing = false;
                inactive = false;
              } else if (inactive) {
                valid = false;
                missing = false;
              } else {
                missing = true;
                valid = false;
              }
            }
          }

          //cuts: exactly the same as for other hits but assuming PXB1

          //D0
          if (!((std::abs(track->dxy(vertices->at(0).position())) * -1.0) <
                TRACK_D0_CUT_BARREL_VAL[trackerTopology_->pxbLayer(detid) - 1]))
            passcuts_hit = false;
          //Dz
          if (!(std::abs(track->dz(vertices->at(0).position())) < TRACK_DZ_CUT_BARREL_VAL))
            passcuts_hit = false;
          // Pixhit cut
          if (!((nBpixL2Hits > 0 && nBpixL3Hits > 0 && nBpixL4Hits > 0) ||
                (nBpixL2Hits > 0 && nBpixL3Hits > 0 && nFpixD1Hits > 0) ||
                (nBpixL2Hits > 0 && nFpixD1Hits > 0 && nFpixD2Hits > 0) ||
                (nFpixD1Hits > 0 && nFpixD2Hits > 0 && nFpixD3Hits > 0)))
            passcuts_hit = false;
          bool found_det = false;

          if (passcuts && passcuts_hit) {
            for (unsigned int i_eff = 0; i_eff < eff_pxb1_vector.size(); i_eff++) {
              //in case found hit in the same det, take only the valid hit
              if (eff_pxb1_vector[i_eff].first == detid) {
                found_det = true;
                if (eff_pxb1_vector[i_eff].second[0] == false && valid == true) {
                  eff_pxb1_vector[i_eff].second[0] = valid;
                  eff_pxb1_vector[i_eff].second[1] = missing;
                  eff_pxb1_vector[i_eff].second[2] = inactive;
                }
              }
            }
            if (!found_det) {
              eff_map.first = detid;
              eff_map.second[0] = valid;
              eff_map.second[1] = missing;
              eff_map.second[2] = inactive;
              eff_pxb1_vector.push_back(eff_map);
            }
          }
        }
      }

      if (eff_pxb1_vector.size() == 1) {
        //eff map is filled -> decide what to do for double hits, ie eff_pxb1_vector.size>1 ... if 1 just use MISSING and VALID as usual

        if (eff_pxb1_vector[0].second[0]) {
          histo[VALID].fill(eff_pxb1_vector[0].first, &iEvent);
          histo[EFFICIENCY].fill(1, eff_pxb1_vector[0].first, &iEvent);
        }
        if (eff_pxb1_vector[0].second[1]) {
          histo[MISSING].fill(eff_pxb1_vector[0].first, &iEvent);
          histo[EFFICIENCY].fill(0, eff_pxb1_vector[0].first, &iEvent);
        }

        if (eff_pxb1_vector[0].second[2]) {
          histo[INACTIVE].fill(eff_pxb1_vector[0].first, &iEvent);
        }
      }

    }  //trajTrackHandle

    histo[VALID].executePerEventHarvesting(&iEvent);
    histo[MISSING].executePerEventHarvesting(&iEvent);
    histo[INACTIVE].executePerEventHarvesting(&iEvent);
  }

}  // namespace

DEFINE_FWK_MODULE(SiPixelPhase1TrackEfficiency);