CMS 3D CMS Logo

PixelUnpackingRegions.h
Go to the documentation of this file.
1 #ifndef PixelUnpackingRegions_H
2 #define PixelUnpackingRegions_H
3 
11 
16 
17 #include <cmath>
18 #include <vector>
19 #include <set>
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 nModules() const { return modules_.size(); }
72  unsigned int nBarrelModules() const;
73  unsigned int nForwardModules() const;
74  unsigned int nRegions() const { return nreg_; }
75 
76  struct Module
77  {
78  float phi;
79  float x, y, z;
80  unsigned int id;
81  unsigned int fed;
82 
83  Module() {}
84  Module(float ph) : phi(ph), x(0.f), y(0.f), z(0.f), id(0), fed(0) {}
85 
86  bool operator < (const Module& m) const
87  {
88  if(phi < m.phi) return true;
89  if(phi == m.phi && id < m.id) return true;
90  return false;
91  }
92  };
93 
94 private:
95 
96  // input parameters
97  std::vector<edm::InputTag> inputs_;
98  std::vector<double> dPhi_;
99  std::vector<double> maxZ_;
101 
103  std::vector<edm::EDGetTokenT<reco::CandidateView>> tCandidateView;
104 
105  std::set<unsigned int> feds_;
106  std::set<unsigned int> modules_;
107  unsigned int nreg_;
108 
110  void initialize(const edm::EventSetup& es);
111 
112  // add a new direction of a region of interest
113  void addRegion(Region &r);
114 
115  // gather info into feds_ and modules_ from a range of a Module vector
116  void gatherFromRange(Region &r, std::vector<Module>::const_iterator, std::vector<Module>::const_iterator);
117 
118  // addRegion for a local (BPIX or +-FPIX) container
119  void addRegionLocal(Region &r, std::vector<Module> &container, const Module& lo,const Module& hi);
120 
121  // local containers of barrel and endcaps Modules sorted by phi
122  std::vector<Module> phiBPIX_;
123  std::vector<Module> phiFPIXp_;
124  std::vector<Module> phiFPIXm_;
125 
126  std::unique_ptr<SiPixelFedCablingTree> cabling_;
128 
130 };
131 
132 #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 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_
bool operator<(const FedChannelConnection &, const FedChannelConnection &)
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:
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_
std::set< unsigned int > modules_