CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PixelUnpackingRegions.h
Go to the documentation of this file.
1 #ifndef PixelUnpackingRegions_H
2 #define PixelUnpackingRegions_H
3 
11 
15 
16 #include <cmath>
17 #include <vector>
18 #include <set>
19 #include <boost/scoped_ptr.hpp>
20 
21 
31 {
32 public:
33 
38  struct Region
39  {
40  Region(const math::XYZVector &dir, float dphi = 0.5f, float maxz = 24.f):
41  v(dir), dPhi(dphi), maxZ(maxz)
42  {
43  cosphi = v.x()/v.rho();
44  sinphi = v.y()/v.rho();
45  atantheta = v.z()/v.rho();
46  }
48  float dPhi, maxZ;
50  };
51 
52 
54 
56 
58  void run(const edm::Event& e, const edm::EventSetup& es);
59 
61  bool mayUnpackFED(unsigned int fed_n) const;
62 
64  bool mayUnpackModule(unsigned int id) const;
65 
67  const std::set<unsigned int> * modulesToUnpack() const {return &modules_;}
68 
70  unsigned int nFEDs() const { return feds_.size(); }
71  unsigned int nBarrelFEDs() const;
72  unsigned int nForwardFEDs() const;
73  unsigned int nModules() const { return modules_.size(); }
74  unsigned int nBarrelModules() const;
75  unsigned int nForwardModules() const;
76  unsigned int nRegions() const { return nreg_; }
77 
78  struct Module
79  {
80  float phi;
81  float x, y, z;
82  unsigned int id;
83  unsigned int fed;
84 
85  Module() {}
86  Module(float ph) : phi(ph), x(0.f), y(0.f), z(0.f), id(0), fed(0) {}
87 
88  bool operator < (const Module& m) const
89  {
90  if(phi < m.phi) return true;
91  if(phi == m.phi && id < m.id) return true;
92  return false;
93  }
94  };
95 
96 private:
97 
98  // input parameters
99  std::vector<edm::InputTag> inputs_;
100  std::vector<double> dPhi_;
101  std::vector<double> maxZ_;
103 
105  std::vector<edm::EDGetTokenT<reco::CandidateView>> tCandidateView;
106 
107  std::set<unsigned int> feds_;
108  std::set<unsigned int> modules_;
109  unsigned int nreg_;
110 
112  void initialize(const edm::EventSetup& es);
113 
114  // add a new direction of a region of interest
115  void addRegion(Region &r);
116 
117  // gather info into feds_ and modules_ from a range of a Module vector
118  void gatherFromRange(Region &r, std::vector<Module>::const_iterator, std::vector<Module>::const_iterator);
119 
120  // addRegion for a local (BPIX or +-FPIX) container
121  void addRegionLocal(Region &r, std::vector<Module> &container, const Module& lo,const Module& hi);
122 
123  // local containers of barrel and endcaps Modules sorted by phi
124  std::vector<Module> phiBPIX_;
125  std::vector<Module> phiFPIXp_;
126  std::vector<Module> phiFPIXm_;
127 
128  boost::scoped_ptr<SiPixelFedCabling> cabling_;
130 
132 };
133 
134 #endif
bool mayUnpackModule(unsigned int id) const
check whether a module has to be unpacked
void addRegionLocal(Region &r, std::vector< Module > &container, const Module &lo, const Module &hi)
unsigned int nForwardFEDs() const
unsigned int nRegions() const
std::vector< edm::EDGetTokenT< reco::CandidateView > > tCandidateView
void gatherFromRange(Region &r, std::vector< Module >::const_iterator, std::vector< Module >::const_iterator)
const std::set< unsigned int > * modulesToUnpack() const
full set of module ids to unpack
boost::scoped_ptr< SiPixelFedCabling > cabling_
PixelUnpackingRegions(const edm::ParameterSet &, edm::ConsumesCollector &&iC)
std::vector< Module > phiFPIXm_
edm::ESWatcher< SiPixelFedCablingMapRcd > watcherSiPixelFedCablingMap_
double f[11][100]
unsigned int nBarrelModules() const
void run(const edm::Event &e, const edm::EventSetup &es)
has to be run during each event
unsigned int nModules() const
std::vector< edm::InputTag > inputs_
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
std::vector< Module > phiFPIXp_
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
unsigned int nForwardModules() const
unsigned int nFEDs() const
various informational accessors:
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_
Region(const math::XYZVector &dir, float dphi=0.5f, float maxz=24.f)
dbl *** dir
Definition: mlp_gen.cc:35
std::vector< Module > phiBPIX_
edm::EDGetTokenT< reco::BeamSpot > tBeamSpot
std::set< unsigned int > feds_
bool operator<(const Module &m) const
std::set< unsigned int > modules_