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
//
//
#include "EventFilter/SiPixelRawToDigi/interface/PixelUnpackingRegions.h"

#include "FWCore/Framework/interface/ESHandle.h"
#include "FWCore/Framework/interface/ESTransientHandle.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "FWCore/MessageLogger/interface/MessageLogger.h"

#include "Geometry/CommonDetUnit/interface/PixelGeomDetUnit.h"

#include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
#include "DataFormats/BeamSpot/interface/BeamSpot.h"
#include "DataFormats/FEDRawData/interface/FEDNumbering.h"
#include "DataFormats/Candidate/interface/LeafCandidate.h"
#include "FWCore/Framework/interface/ConsumesCollector.h"
#include "DataFormats/FEDRawData/interface/FEDNumbering.h"

#include <algorithm>
#include <iterator>

// local convenience functions
namespace {
  bool isBPIXModule(unsigned int id) { return DetId(id).subdetId() == PixelSubdetector::PixelBarrel; }
  bool isFPIXModule(unsigned int id) { return DetId(id).subdetId() == PixelSubdetector::PixelEndcap; }
}  // namespace

PixelUnpackingRegions::PixelUnpackingRegions(const edm::ParameterSet& conf, edm::ConsumesCollector&& iC) {
  edm::ParameterSet regPSet = conf.getParameter<edm::ParameterSet>("Regions");
  beamSpotTag_ = regPSet.getParameter<edm::InputTag>("beamSpot");
  inputs_ = regPSet.getParameter<std::vector<edm::InputTag> >("inputs");
  dPhi_ = regPSet.getParameter<std::vector<double> >("deltaPhi");
  maxZ_ = regPSet.getParameter<std::vector<double> >("maxZ");
  trackerGeomToken_ = iC.esConsumes<TrackerGeometry, TrackerDigiGeometryRecord>();
  cablingMapToken_ = iC.esConsumes<SiPixelFedCablingMap, SiPixelFedCablingMapRcd>();

  tBeamSpot = iC.consumes<reco::BeamSpot>(beamSpotTag_);
  for (unsigned int t = 0; t < inputs_.size(); t++)
    tCandidateView.push_back(iC.consumes<reco::CandidateView>(inputs_[t]));

  if (inputs_.size() != dPhi_.size() || dPhi_.size() != maxZ_.size()) {
    edm::LogError("PixelUnpackingRegions")
        << "Not the same size of config parameters vectors!\n"
        << "   inputs " << inputs_.size() << "  deltaPhi " << dPhi_.size() << "  maxZ " << maxZ_.size();
  }
}

void PixelUnpackingRegions::run(const edm::Event& e, const edm::EventSetup& es) {
  feds_.clear();
  modules_.clear();
  nreg_ = 0;

  initialize(es);

  edm::Handle<reco::BeamSpot> beamSpot;
  e.getByToken(tBeamSpot, beamSpot);
  beamSpot_ = beamSpot->position();
  //beamSpot_ = math::XYZPoint(0.,0.,0.);

  size_t ninputs = inputs_.size();
  for (size_t input = 0; input < ninputs; ++input) {
    edm::Handle<reco::CandidateView> h;
    e.getByToken(tCandidateView[input], h);

    size_t n = h->size();
    for (size_t i = 0; i < n; ++i) {
      const reco::Candidate& c = (*h)[i];

      // different input collections can have different dPhi and maxZ
      Region r(c.momentum(), dPhi_[input], maxZ_[input]);
      addRegion(r);
    }
  }
}

void PixelUnpackingRegions::initialize(const edm::EventSetup& es) {
  // initialize cabling map or update it if necessary
  // and re-cache modules information
  if (watcherSiPixelFedCablingMap_.check(es)) {
    edm::ESHandle<SiPixelFedCablingMap> cablingMap = es.getHandle(cablingMapToken_);
    cabling_ = cablingMap->cablingTree();

    // get the TrackerGeom
    edm::ESHandle<TrackerGeometry> geom = es.getHandle(trackerGeomToken_);

    // switch on the phase1
    unsigned int fedMin = FEDNumbering::MINSiPixelFEDID;  // phase0
    unsigned int fedMax = FEDNumbering::MAXSiPixelFEDID;
    if ((geom->isThere(GeomDetEnumerators::P1PXB)) && (geom->isThere(GeomDetEnumerators::P1PXEC))) {
      fedMin = FEDNumbering::MINSiPixeluTCAFEDID;  // phase1
      fedMax = FEDNumbering::MAXSiPixeluTCAFEDID;
    }

    phiBPIX_.clear();
    phiFPIXp_.clear();
    phiFPIXm_.clear();

    phiBPIX_.reserve(1024);
    phiFPIXp_.reserve(512);
    phiFPIXm_.reserve(512);

    auto it = geom->dets().begin();
    for (; it != geom->dets().end(); ++it) {
      int subdet = (*it)->geographicalId().subdetId();
      if (!(subdet == PixelSubdetector::PixelBarrel || subdet == PixelSubdetector::PixelEndcap))
        continue;

      Module m;

      m.x = (*it)->position().x();
      m.y = (*it)->position().y();
      m.z = (*it)->position().z();

      m.phi = (*it)->position().phi();

      m.id = (*it)->geographicalId().rawId();
      const std::vector<sipixelobjects::CablingPathToDetUnit> path2det = cabling_->pathToDetUnit(m.id);

      m.fed = path2det[0].fed;
      assert((m.fed <= fedMax) && (m.fed >= fedMin));

      if (subdet == PixelSubdetector::PixelBarrel) {
        phiBPIX_.push_back(m);
      } else if (subdet == PixelSubdetector::PixelEndcap) {
        if (m.z > 0.)
          phiFPIXp_.push_back(m);
        else
          phiFPIXm_.push_back(m);
      }
    }

    // pre-sort by phi
    std::sort(phiBPIX_.begin(), phiBPIX_.end());
    std::sort(phiFPIXp_.begin(), phiFPIXp_.end());
    std::sort(phiFPIXm_.begin(), phiFPIXm_.end());
  }
}

void PixelUnpackingRegions::addRegion(Region& r) {
  ++nreg_;

  float phi = r.v.phi();

  Module lo(phi - r.dPhi);
  Module hi(phi + r.dPhi);

  addRegionLocal(r, phiBPIX_, lo, hi);
  if (r.v.eta() > 1.) {
    addRegionLocal(r, phiFPIXp_, lo, hi);
  }
  if (r.v.eta() < -1.) {
    addRegionLocal(r, phiFPIXm_, lo, hi);
  }
}

void PixelUnpackingRegions::addRegionLocal(Region& r,
                                           std::vector<Module>& container,
                                           const Module& _lo,
                                           const Module& _hi) {
  Module lo = _lo;
  Module hi = _hi;
  Module pi_m(-M_PI);
  Module pi_p(M_PI);

  std::vector<Module>::const_iterator a, b;

  if (lo.phi >= -M_PI && hi.phi <= M_PI)  // interval doesn't cross the +-pi overlap
  {
    a = lower_bound(container.begin(), container.end(), lo);
    b = upper_bound(container.begin(), container.end(), hi);
    gatherFromRange(r, a, b);
  } else  // interval is torn by the +-pi overlap
  {
    if (hi.phi > M_PI)
      hi.phi -= 2. * M_PI;
    a = lower_bound(container.begin(), container.end(), pi_m);
    b = upper_bound(container.begin(), container.end(), hi);
    gatherFromRange(r, a, b);

    if (lo.phi < -M_PI)
      lo.phi += 2. * M_PI;
    a = lower_bound(container.begin(), container.end(), lo);
    b = upper_bound(container.begin(), container.end(), pi_p);
    gatherFromRange(r, a, b);
  }
}

void PixelUnpackingRegions::gatherFromRange(Region& r,
                                            std::vector<Module>::const_iterator a,
                                            std::vector<Module>::const_iterator b) {
  for (; a != b; ++a) {
    // projection in r's direction onto beam's z
    float zmodule = a->z - ((a->x - beamSpot_.x()) * r.cosphi + (a->y - beamSpot_.y()) * r.sinphi) * r.atantheta;

    // do not include modules that project too far in z
    if (std::abs(zmodule) > r.maxZ)
      continue;

    feds_.insert(a->fed);
    modules_.insert(a->id);
  }
}

bool PixelUnpackingRegions::mayUnpackFED(unsigned int fed_n) const {
  if (feds_.count(fed_n))
    return true;
  return false;
}

bool PixelUnpackingRegions::mayUnpackModule(unsigned int id) const {
  if (modules_.count(id))
    return true;
  return false;
}

unsigned int PixelUnpackingRegions::nBarrelModules() const {
  return std::count_if(modules_.begin(), modules_.end(), isBPIXModule);
}

unsigned int PixelUnpackingRegions::nForwardModules() const {
  return std::count_if(modules_.begin(), modules_.end(), isFPIXModule);
}