CMS 3D CMS Logo

PixelUnpackingRegions.cc
Go to the documentation of this file.
1 //
2 //
4 
9 
11 
18 
19 #include <algorithm>
20 #include <iterator>
21 
22 // local convenience functions
23 namespace {
24  bool isBPIXModule(unsigned int id) { return DetId(id).subdetId() == PixelSubdetector::PixelBarrel; }
25  bool isFPIXModule(unsigned int id) { return DetId(id).subdetId() == PixelSubdetector::PixelEndcap; }
26 } // namespace
27 
29  edm::ParameterSet regPSet = conf.getParameter<edm::ParameterSet>("Regions");
30  beamSpotTag_ = regPSet.getParameter<edm::InputTag>("beamSpot");
31  inputs_ = regPSet.getParameter<std::vector<edm::InputTag> >("inputs");
32  dPhi_ = regPSet.getParameter<std::vector<double> >("deltaPhi");
33  maxZ_ = regPSet.getParameter<std::vector<double> >("maxZ");
36 
37  tBeamSpot = iC.consumes<reco::BeamSpot>(beamSpotTag_);
38  for (unsigned int t = 0; t < inputs_.size(); t++)
39  tCandidateView.push_back(iC.consumes<reco::CandidateView>(inputs_[t]));
40 
41  if (inputs_.size() != dPhi_.size() || dPhi_.size() != maxZ_.size()) {
42  edm::LogError("PixelUnpackingRegions")
43  << "Not the same size of config parameters vectors!\n"
44  << " inputs " << inputs_.size() << " deltaPhi " << dPhi_.size() << " maxZ " << maxZ_.size();
45  }
46 }
47 
49  feds_.clear();
50  modules_.clear();
51  nreg_ = 0;
52 
53  initialize(es);
54 
56  e.getByToken(tBeamSpot, beamSpot);
57  beamSpot_ = beamSpot->position();
58  //beamSpot_ = math::XYZPoint(0.,0.,0.);
59 
60  size_t ninputs = inputs_.size();
61  for (size_t input = 0; input < ninputs; ++input) {
63  e.getByToken(tCandidateView[input], h);
64 
65  size_t n = h->size();
66  for (size_t i = 0; i < n; ++i) {
67  const reco::Candidate& c = (*h)[i];
68 
69  // different input collections can have different dPhi and maxZ
70  Region r(c.momentum(), dPhi_[input], maxZ_[input]);
71  addRegion(r);
72  }
73  }
74 }
75 
77  // initialize cabling map or update it if necessary
78  // and re-cache modules information
81  cabling_ = cablingMap->cablingTree();
82 
83  // get the TrackerGeom
85 
86  // switch on the phase1
87  unsigned int fedMin = FEDNumbering::MINSiPixelFEDID; // phase0
88  unsigned int fedMax = FEDNumbering::MAXSiPixelFEDID;
89  if ((geom->isThere(GeomDetEnumerators::P1PXB)) && (geom->isThere(GeomDetEnumerators::P1PXEC))) {
90  fedMin = FEDNumbering::MINSiPixeluTCAFEDID; // phase1
92  }
93 
94  phiBPIX_.clear();
95  phiFPIXp_.clear();
96  phiFPIXm_.clear();
97 
98  phiBPIX_.reserve(1024);
99  phiFPIXp_.reserve(512);
100  phiFPIXm_.reserve(512);
101 
102  auto it = geom->dets().begin();
103  for (; it != geom->dets().end(); ++it) {
104  int subdet = (*it)->geographicalId().subdetId();
105  if (!(subdet == PixelSubdetector::PixelBarrel || subdet == PixelSubdetector::PixelEndcap))
106  continue;
107 
108  Module m;
109 
110  m.x = (*it)->position().x();
111  m.y = (*it)->position().y();
112  m.z = (*it)->position().z();
113 
114  m.phi = (*it)->position().phi();
115 
116  m.id = (*it)->geographicalId().rawId();
117  const std::vector<sipixelobjects::CablingPathToDetUnit> path2det = cabling_->pathToDetUnit(m.id);
118 
119  m.fed = path2det[0].fed;
120  assert((m.fed <= fedMax) && (m.fed >= fedMin));
121 
122  if (subdet == PixelSubdetector::PixelBarrel) {
123  phiBPIX_.push_back(m);
124  } else if (subdet == PixelSubdetector::PixelEndcap) {
125  if (m.z > 0.)
126  phiFPIXp_.push_back(m);
127  else
128  phiFPIXm_.push_back(m);
129  }
130  }
131 
132  // pre-sort by phi
133  std::sort(phiBPIX_.begin(), phiBPIX_.end());
134  std::sort(phiFPIXp_.begin(), phiFPIXp_.end());
135  std::sort(phiFPIXm_.begin(), phiFPIXm_.end());
136  }
137 }
138 
140  ++nreg_;
141 
142  float phi = r.v.phi();
143 
144  Module lo(phi - r.dPhi);
145  Module hi(phi + r.dPhi);
146 
147  addRegionLocal(r, phiBPIX_, lo, hi);
148  if (r.v.eta() > 1.) {
149  addRegionLocal(r, phiFPIXp_, lo, hi);
150  }
151  if (r.v.eta() < -1.) {
152  addRegionLocal(r, phiFPIXm_, lo, hi);
153  }
154 }
155 
157  std::vector<Module>& container,
158  const Module& _lo,
159  const Module& _hi) {
160  Module lo = _lo;
161  Module hi = _hi;
162  Module pi_m(-M_PI);
163  Module pi_p(M_PI);
164 
165  std::vector<Module>::const_iterator a, b;
166 
167  if (lo.phi >= -M_PI && hi.phi <= M_PI) // interval doesn't cross the +-pi overlap
168  {
169  a = lower_bound(container.begin(), container.end(), lo);
170  b = upper_bound(container.begin(), container.end(), hi);
171  gatherFromRange(r, a, b);
172  } else // interval is torn by the +-pi overlap
173  {
174  if (hi.phi > M_PI)
175  hi.phi -= 2. * M_PI;
176  a = lower_bound(container.begin(), container.end(), pi_m);
177  b = upper_bound(container.begin(), container.end(), hi);
178  gatherFromRange(r, a, b);
179 
180  if (lo.phi < -M_PI)
181  lo.phi += 2. * M_PI;
182  a = lower_bound(container.begin(), container.end(), lo);
183  b = upper_bound(container.begin(), container.end(), pi_p);
184  gatherFromRange(r, a, b);
185  }
186 }
187 
189  std::vector<Module>::const_iterator a,
190  std::vector<Module>::const_iterator b) {
191  for (; a != b; ++a) {
192  // projection in r's direction onto beam's z
193  float zmodule = a->z - ((a->x - beamSpot_.x()) * r.cosphi + (a->y - beamSpot_.y()) * r.sinphi) * r.atantheta;
194 
195  // do not include modules that project too far in z
196  if (std::abs(zmodule) > r.maxZ)
197  continue;
198 
199  feds_.insert(a->fed);
200  modules_.insert(a->id);
201  }
202 }
203 
204 bool PixelUnpackingRegions::mayUnpackFED(unsigned int fed_n) const {
205  if (feds_.count(fed_n))
206  return true;
207  return false;
208 }
209 
210 bool PixelUnpackingRegions::mayUnpackModule(unsigned int id) const {
211  if (modules_.count(id))
212  return true;
213  return false;
214 }
215 
217  return std::count_if(modules_.begin(), modules_.end(), isBPIXModule);
218 }
219 
221  return std::count_if(modules_.begin(), modules_.end(), isFPIXModule);
222 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
std::unique_ptr< SiPixelFedCablingTree > cablingTree() const
void addRegionLocal(Region &r, std::vector< Module > &container, const Module &lo, const Module &hi)
std::vector< edm::EDGetTokenT< reco::CandidateView > > tCandidateView
void gatherFromRange(Region &r, std::vector< Module >::const_iterator, std::vector< Module >::const_iterator)
Log< level::Error, false > LogError
unsigned int nForwardModules() const
assert(be >=bs)
PixelUnpackingRegions(const edm::ParameterSet &, edm::ConsumesCollector &&iC)
static std::string const input
Definition: EdmProvDump.cc:50
std::vector< Module > phiFPIXm_
std::unique_ptr< SiPixelFedCablingTree > cabling_
edm::ESWatcher< SiPixelFedCablingMapRcd > watcherSiPixelFedCablingMap_
Definition: EPCuts.h:4
bool mayUnpackModule(unsigned int id) const
check whether a module has to be unpacked
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:130
#define M_PI
void run(const edm::Event &e, const edm::EventSetup &es)
has to be run during each event
Definition: DetId.h:17
edm::ESGetToken< SiPixelFedCablingMap, SiPixelFedCablingMapRcd > cablingMapToken_
std::vector< edm::InputTag > inputs_
std::vector< Module > phiFPIXp_
double b
Definition: hdecay.h:120
unsigned int nBarrelModules() const
bool check(const edm::EventSetup &iSetup)
Definition: ESWatcher.h:57
double a
Definition: hdecay.h:121
std::vector< double > dPhi_
void initialize(const edm::EventSetup &es)
run by the run method: (re)initialize the cabling data when it&#39;s necessary
std::vector< double > maxZ_
bool mayUnpackFED(unsigned int fed_n) const
check whether a FED has to be unpacked
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > trackerGeomToken_
std::vector< Module > phiBPIX_
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
edm::EDGetTokenT< reco::BeamSpot > tBeamSpot
std::set< unsigned int > feds_
std::set< unsigned int > modules_