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 
9 
12 
13 #include <cmath>
14 #include <vector>
15 #include <set>
16 #include <boost/scoped_ptr.hpp>
17 
18 
28 {
29 public:
30 
35  struct Region
36  {
37  Region(const math::XYZVector &dir, float dphi = 0.5f, float maxz = 24.f):
38  v(dir), dPhi(dphi), maxZ(maxz)
39  {
40  cosphi = v.x()/v.rho();
41  sinphi = v.y()/v.rho();
42  atantheta = v.z()/v.rho();
43  }
45  float dPhi, maxZ;
47  };
48 
49 
51 
53 
55  void run(const edm::Event& e, const edm::EventSetup& es);
56 
58  bool mayUnpackFED(unsigned int fed_n) const;
59 
61  bool mayUnpackModule(unsigned int id) const;
62 
64  const std::set<unsigned int> * modulesToUnpack() const {return &modules_;}
65 
67  unsigned int nFEDs() const { return feds_.size(); }
68  unsigned int nBarrelFEDs() const;
69  unsigned int nForwardFEDs() const;
70  unsigned int nModules() const { return modules_.size(); }
71  unsigned int nBarrelModules() const;
72  unsigned int nForwardModules() const;
73  unsigned int nRegions() const { return nreg_; }
74 
75  struct Module
76  {
77  float phi;
78  float x, y, z;
79  unsigned int id;
80  unsigned int fed;
81 
82  Module() {}
83  Module(float ph) : phi(ph), x(0.f), y(0.f), z(0.f), id(0), fed(0) {}
84 
85  bool operator < (const Module& m) const
86  {
87  if(phi < m.phi) return true;
88  if(phi == m.phi && id < m.id) return true;
89  return false;
90  }
91  };
92 
93 private:
94 
95  // input parameters
96  std::vector<edm::InputTag> inputs_;
97  std::vector<double> dPhi_;
98  std::vector<double> maxZ_;
100 
101  std::set<unsigned int> feds_;
102  std::set<unsigned int> modules_;
103  unsigned int nreg_;
104 
106  void initialize(const edm::EventSetup& es);
107 
108  // add a new direction of a region of interest
109  void addRegion(Region &r);
110 
111  // gather info into feds_ and modules_ from a range of a Module vector
112  void gatherFromRange(Region &r, std::vector<Module>::const_iterator, std::vector<Module>::const_iterator);
113 
114  // addRegion for a local (BPIX or +-FPIX) container
115  void addRegionLocal(Region &r, std::vector<Module> &container, Module lo, Module hi);
116 
117  // local containers of barrel and endcaps Modules sorted by phi
118  std::vector<Module> phiBPIX_;
119  std::vector<Module> phiFPIXp_;
120  std::vector<Module> phiFPIXm_;
121 
122  boost::scoped_ptr<SiPixelFedCabling> cabling_;
124 
126 };
127 
128 #endif
bool mayUnpackModule(unsigned int id) const
check whether a module has to be unpacked
unsigned int nForwardFEDs() const
unsigned int nRegions() const
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_
std::vector< Module > phiFPIXm_
PixelUnpackingRegions(const edm::ParameterSet &)
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
void addRegionLocal(Region &r, std::vector< Module > &container, Module lo, Module hi)
std::vector< edm::InputTag > inputs_
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
std::vector< Module > phiFPIXp_
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:13
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_
std::set< unsigned int > feds_
bool operator<(const Module &m) const
std::set< unsigned int > modules_