test
CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PixelUnpackingRegions.cc
Go to the documentation of this file.
1 //
2 //
4 
9 
13 
19 
20 #include <algorithm>
21 #include <iterator>
22 
23 // local convenience functions
24 namespace
25 {
26 bool isBPIXFED(unsigned int fed) {return fed< 32;}
27 bool isFPIXFED(unsigned int fed) {return fed>=32;}
28 bool isBPIXModule(unsigned int id) {return DetId(id).subdetId() == PixelSubdetector::PixelBarrel;}
29 bool isFPIXModule(unsigned int id) {return DetId(id).subdetId() == PixelSubdetector::PixelEndcap;}
30 
31 inline std::ostream& operator<<(std::ostream& s, const PixelUnpackingRegions::Module& m)
32 {
33  s<< (isBPIXModule(m.id) ? "BPIX " : "FPIX ") <<m.id<<" "<<m.fed<<" "<<m.phi<<" "<<m.x<<" "<<m.y<<" "<<m.z<<" "<<sqrt(std::pow(m.x,2)+std::pow(m.y,2));
34  return s;
35 }
36 }
37 
38 
39 
41 {
42  edm::ParameterSet regPSet = conf.getParameter<edm::ParameterSet>("Regions");
43  beamSpotTag_ = regPSet.getParameter<edm::InputTag>("beamSpot");
44  inputs_ = regPSet.getParameter<std::vector<edm::InputTag> >("inputs");
45  dPhi_ = regPSet.getParameter<std::vector<double> >("deltaPhi");
46  maxZ_ = regPSet.getParameter<std::vector<double> >("maxZ");
47 
48  tBeamSpot = iC.consumes<reco::BeamSpot>(beamSpotTag_);
49  for (unsigned int t=0; t<inputs_.size(); t++ ) tCandidateView.push_back(iC.consumes< reco::CandidateView >(inputs_[t]));
50 
51  if (inputs_.size() != dPhi_.size() || dPhi_.size() != maxZ_.size() )
52  {
53  edm::LogError("PixelUnpackingRegions")<<"Not the same size of config parameters vectors!\n"
54  <<" inputs "<<inputs_.size()<<" deltaPhi "<<dPhi_.size() <<" maxZ "<< maxZ_.size();
55  }
56 }
57 
58 
60 {
61  feds_.clear();
62  modules_.clear();
63  nreg_ = 0;
64 
65  initialize(es);
66 
68  e.getByToken(tBeamSpot, beamSpot);
69  beamSpot_ = beamSpot->position();
70  //beamSpot_ = math::XYZPoint(0.,0.,0.);
71 
72  size_t ninputs = inputs_.size();
73  for(size_t input = 0; input < ninputs; ++input)
74  {
77 
78  size_t n = h->size();
79  for(size_t i = 0; i < n; ++i )
80  {
81  const reco::Candidate & c = (*h)[i];
82 
83  // different input collections can have different dPhi and maxZ
85  addRegion(r);
86  }
87  }
88 }
89 
90 
92 {
93  // initialize cabling map or update it if necessary
94  // and re-cache modules information
96  {
98  es.get<SiPixelFedCablingMapRcd>().get( cablingMap );
99  cabling_ = cablingMap->cablingTree();
100 
102  // get the TrackerGeom
103  es.get<TrackerDigiGeometryRecord>().get( geom );
104 
105  phiBPIX_.clear();
106  phiFPIXp_.clear();
107  phiFPIXm_.clear();
108 
109  phiBPIX_.reserve(1024);
110  phiFPIXp_.reserve(512);
111  phiFPIXm_.reserve(512);
112 
113  auto it = geom->dets().begin();
114  for ( ; it != geom->dets().end(); ++it)
115  {
116  int subdet = (*it)->geographicalId().subdetId();
117  if (! (subdet == PixelSubdetector::PixelBarrel ||
118  subdet == PixelSubdetector::PixelEndcap) ) continue;
119 
120  Module m;
121 
122  m.x = (*it)->position().x();
123  m.y = (*it)->position().y();
124  m.z = (*it)->position().z();
125 
126  m.phi = (*it)->position().phi();
127 
128  m.id = (*it)->geographicalId().rawId();
129  const std::vector<sipixelobjects::CablingPathToDetUnit> path2det = cabling_->pathToDetUnit(m.id);
130 
131  m.fed = path2det[0].fed;
132  assert(m.fed<40);
133 
134  if (subdet == PixelSubdetector::PixelBarrel)
135  {
136  phiBPIX_.push_back(m);
137  }
138  else if (subdet == PixelSubdetector::PixelEndcap)
139  {
140  if (m.z > 0.) phiFPIXp_.push_back(m);
141  else phiFPIXm_.push_back(m);
142  }
143  }
144 
145  // pre-sort by phi
146  std::sort(phiBPIX_.begin(), phiBPIX_.end());
147  std::sort(phiFPIXp_.begin(), phiFPIXp_.end());
148  std::sort(phiFPIXm_.begin(), phiFPIXm_.end());
149  }
150 }
151 
152 
154 {
155  ++nreg_;
156 
157  float phi = r.v.phi();
158 
159  Module lo(phi - r.dPhi);
160  Module hi(phi + r.dPhi);
161 
162  addRegionLocal(r, phiBPIX_, lo, hi);
163  if (r.v.eta() > 1.)
164  {
165  addRegionLocal(r, phiFPIXp_, lo, hi);
166  }
167  if (r.v.eta() < -1.)
168  {
169  addRegionLocal(r, phiFPIXm_, lo, hi);
170  }
171 }
172 
173 
174 void PixelUnpackingRegions::addRegionLocal(Region &r, std::vector<Module> &container,const Module& _lo,const Module& _hi)
175 {
176  Module lo = _lo;
177  Module hi = _hi;
178  Module pi_m(-M_PI);
179  Module pi_p( M_PI);
180 
181  std::vector<Module>::const_iterator a, b;
182 
183  if (lo.phi >= -M_PI && hi.phi <= M_PI) // interval doesn't cross the +-pi overlap
184  {
185  a = lower_bound(container.begin(), container.end(), lo);
186  b = upper_bound(container.begin(), container.end(), hi);
187  gatherFromRange(r, a, b);
188  }
189  else // interval is torn by the +-pi overlap
190  {
191  if (hi.phi > M_PI) hi.phi -= 2.*M_PI;
192  a = lower_bound(container.begin(), container.end(), pi_m);
193  b = upper_bound(container.begin(), container.end(), hi);
194  gatherFromRange(r, a, b);
195 
196  if (lo.phi < -M_PI) lo.phi += 2.*M_PI;
197  a = lower_bound(container.begin(), container.end(), lo);
198  b = upper_bound(container.begin(), container.end(), pi_p);
199  gatherFromRange(r, a, b);
200  }
201 }
202 
203 
204 void PixelUnpackingRegions::gatherFromRange(Region &r, std::vector<Module>::const_iterator a, std::vector<Module>::const_iterator b)
205 {
206  for(; a != b; ++a)
207  {
208  // projection in r's direction onto beam's z
209  float zmodule = a->z - ( (a->x - beamSpot_.x())*r.cosphi + (a->y - beamSpot_.y())*r.sinphi ) * r.atantheta;
210 
211  // do not include modules that project too far in z
212  if ( std::abs(zmodule) > r.maxZ ) continue;
213 
214  feds_.insert(a->fed);
215  modules_.insert(a->id);
216  }
217 }
218 
219 
220 bool PixelUnpackingRegions::mayUnpackFED(unsigned int fed_n) const
221 {
222  if (feds_.count(fed_n)) return true;
223  return false;
224 }
225 
227 {
228  return std::count_if(feds_.begin(), feds_.end(), isBPIXFED );
229 }
230 
232 {
233  return std::count_if(feds_.begin(), feds_.end(), isFPIXFED );
234 }
235 
236 
237 bool PixelUnpackingRegions::mayUnpackModule(unsigned int id) const
238 {
239  if (modules_.count(id)) return true;
240  return false;
241 }
242 
244 {
245  return std::count_if(modules_.begin(), modules_.end(), isBPIXModule );
246 }
247 
249 {
250  return std::count_if(modules_.begin(), modules_.end(), isFPIXModule );
251 }
T getParameter(std::string const &) const
bool mayUnpackModule(unsigned int id) const
check whether a module has to be unpacked
int i
Definition: DBlmapReader.cc:9
void addRegionLocal(Region &r, std::vector< Module > &container, const Module &lo, const Module &hi)
FWCore Framework interface EventSetupRecordImplementation h
Helper function to determine trigger accepts.
unsigned int nForwardFEDs() const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
std::vector< edm::EDGetTokenT< reco::CandidateView > > tCandidateView
assert(m_qm.get())
void gatherFromRange(Region &r, std::vector< Module >::const_iterator, std::vector< Module >::const_iterator)
std::ostream & operator<<(std::ostream &out, const ALILine &li)
Definition: ALILine.cc:188
PixelUnpackingRegions(const edm::ParameterSet &, edm::ConsumesCollector &&iC)
static std::string const input
Definition: EdmProvDump.cc:44
std::vector< Module > phiFPIXm_
std::unique_ptr< SiPixelFedCablingTree > cabling_
virtual Vector momentum() const =0
spatial momentum vector
edm::ESWatcher< SiPixelFedCablingMapRcd > watcherSiPixelFedCablingMap_
T sqrt(T t)
Definition: SSEVec.h:18
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
unsigned int nBarrelModules() const
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
#define M_PI
void run(const edm::Event &e, const edm::EventSetup &es)
has to be run during each event
Definition: DetId.h:18
std::vector< edm::InputTag > inputs_
std::vector< Module > phiFPIXp_
const T & get() const
Definition: EventSetup.h:56
double b
Definition: hdecay.h:120
bool check(const edm::EventSetup &iSetup)
Definition: ESWatcher.h:57
unsigned int nForwardModules() const
double a
Definition: hdecay.h:121
unsigned int nBarrelFEDs() const
bool mayUnpackFED(unsigned int fed_n) const
check whether a FED has to be unpacked
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_
std::vector< Module > phiBPIX_
edm::EDGetTokenT< reco::BeamSpot > tBeamSpot
std::set< unsigned int > feds_
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
std::set< unsigned int > modules_