CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
PixelUnpackingRegions.h
Go to the documentation of this file.
1 #ifndef EventFilter_SiPixelRawToDigi_interface_PixelUnpackingRegions_h
2 #define EventFilter_SiPixelRawToDigi_interface_PixelUnpackingRegions_h
3 
12 
19 
20 #include <cmath>
21 #include <vector>
22 #include <set>
23 
33 public:
38  struct Region {
39  Region(const math::XYZVector& dir, float dphi = 0.5f, float maxz = 24.f) : v(dir), dPhi(dphi), maxZ(maxz) {
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 
50 
52 
54  void run(const edm::Event& e, const edm::EventSetup& es);
55 
57  bool mayUnpackFED(unsigned int fed_n) const;
58 
60  bool mayUnpackModule(unsigned int id) const;
61 
63  const std::set<unsigned int>* modulesToUnpack() const { return &modules_; }
64 
66  unsigned int nFEDs() const { return feds_.size(); }
67  unsigned int nModules() const { return modules_.size(); }
68  unsigned int nBarrelModules() const;
69  unsigned int nForwardModules() const;
70  unsigned int nRegions() const { return nreg_; }
71 
72  struct Module {
73  float phi;
74  float x, y, z;
75  unsigned int id;
76  unsigned int fed;
77 
78  Module() {}
79  Module(float ph) : phi(ph), x(0.f), y(0.f), z(0.f), id(0), fed(0) {}
80 
81  bool operator<(const Module& m) const {
82  if (phi < m.phi)
83  return true;
84  if (phi == m.phi && id < m.id)
85  return true;
86  return false;
87  }
88  };
89 
90 private:
91  // input parameters
92  std::vector<edm::InputTag> inputs_;
93  std::vector<double> dPhi_;
94  std::vector<double> maxZ_;
96 
98  std::vector<edm::EDGetTokenT<reco::CandidateView>> tCandidateView;
99 
100  std::set<unsigned int> feds_;
101  std::set<unsigned int> modules_;
102  unsigned int nreg_;
103 
105  void initialize(const edm::EventSetup& es);
106 
107  // add a new direction of a region of interest
108  void addRegion(Region& r);
109 
110  // gather info into feds_ and modules_ from a range of a Module vector
111  void gatherFromRange(Region& r, std::vector<Module>::const_iterator, std::vector<Module>::const_iterator);
112 
113  // addRegion for a local (BPIX or +-FPIX) container
114  void addRegionLocal(Region& r, std::vector<Module>& container, const Module& lo, const Module& hi);
115 
116  // local containers of barrel and endcaps Modules sorted by phi
117  std::vector<Module> phiBPIX_;
118  std::vector<Module> phiFPIXp_;
119  std::vector<Module> phiFPIXm_;
120 
121  std::unique_ptr<SiPixelFedCablingTree> cabling_;
123 
127 };
128 
129 #endif // EventFilter_SiPixelRawToDigi_interface_PixelUnpackingRegions_h
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 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
PixelUnpackingRegions(const edm::ParameterSet &, edm::ConsumesCollector &&iC)
std::vector< Module > phiFPIXm_
std::unique_ptr< SiPixelFedCablingTree > cabling_
edm::ESWatcher< SiPixelFedCablingMapRcd > watcherSiPixelFedCablingMap_
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
edm::ESGetToken< SiPixelFedCablingMap, SiPixelFedCablingMapRcd > cablingMapToken_
std::vector< edm::InputTag > inputs_
__constant__ float const maxz[nPairs]
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:12
unsigned int nForwardModules() const
unsigned int nFEDs() const
various informational accessors:
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)
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > trackerGeomToken_
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_