CMS 3D CMS Logo

HGCalTriggerGeometryGenericMapping.h
Go to the documentation of this file.
1 #ifndef __L1Trigger_L1THGCal_HGCalTriggerGeometryGenericMapping_h__
2 #define __L1Trigger_L1THGCal_HGCalTriggerGeometryGenericMapping_h__
3 
5 
6 /*******
7  *
8  * class: HGCalTriggerGeometryGenericMapping
9  * author: L.Gray (FNAL)
10  * date: 26 July, 2015
11  *
12  * This class is the base class for all HGCal trigger geometry definitions.
13  * Typically this is just a ganging of cells and nearest-neighbour relationships.
14  *
15  * Classes for TriggerCells and Modules are defined. They are containers for
16  * maps relating modules to modules and trigger cells to trigger cells and
17  * raw HGC cells.
18  *
19  * It is up to the user of the class to define what these mappings are through the
20  * initialize() function. This function takes the full HGC geometry as input and
21  * creates all the necessary maps for navigating the trigger system. All of the
22  * containers for parts of the map are available through protected members of the base
23  * class.
24  *
25  * All unsigned ints used here are DetIds, or will become them once we sort out the
26  * full nature of the trigger detids.
27  *******/
28 
30 
32  class TriggerCell {
33  public:
34  typedef std::unordered_set<unsigned> list_type;
35 
36  TriggerCell(unsigned tc_id, unsigned mod_id, const GlobalPoint& pos,
37  const list_type& neighbs, const list_type& comps) :
38  trigger_cell_id_(tc_id),
39  module_id_(mod_id),
40  position_(pos),
41  neighbours_(neighbs),
42  components_(comps)
43  {}
45 
46 
47  unsigned triggerCellId() const { return trigger_cell_id_; }
48  unsigned moduleId() const { return module_id_; }
49 
50  bool containsCell(const unsigned cell) const {
51  return ( components_.find(cell) != components_.end() );
52  }
53 
54  const GlobalPoint& position() const { return position_; }
55 
56  const std::unordered_set<unsigned>& neighbours() const { return neighbours_; }
57  const std::unordered_set<unsigned>& components() const { return components_; }
58 
59  private:
60  unsigned trigger_cell_id_; // the ID of this trigger cell
61  unsigned module_id_; // module this TC belongs to
63  list_type neighbours_; // neighbouring trigger cells
64  list_type components_; // contained HGC cells
65  };
66 
67  class Module {
68  public:
69  typedef std::unordered_set<unsigned> list_type;
70  typedef std::unordered_multimap<unsigned,unsigned> tc_map_type;
71 
72  Module(unsigned mod_id, const GlobalPoint& pos,
73  const list_type& neighbs, const list_type& comps,
74  const tc_map_type& tc_comps):
75  module_id_(mod_id),
76  position_(pos),
77  neighbours_(neighbs),
78  components_(comps),
79  tc_components_(tc_comps)
80  {}
81  ~Module() {}
82 
83  unsigned moduleId() const { return module_id_; }
84 
85  bool containsTriggerCell(const unsigned trig_cell) const {
86  return ( components_.find(trig_cell) != components_.end() );
87  }
88 
89  bool containsCell(const unsigned cell) const {
90  for( const auto& value : tc_components_ ) {
91  if( value.second == cell ) return true;
92  }
93  return false;
94  }
95 
96  const GlobalPoint& position() const { return position_; }
97 
98  const list_type& neighbours() const { return neighbours_; }
99  const list_type& components() const { return components_; }
100 
101  const tc_map_type& triggerCellComponents() const { return tc_components_; }
102 
103  private:
104  unsigned module_id_; // module this TC belongs to
106  list_type neighbours_; // neighbouring Modules
107  list_type components_; // contained HGC trigger cells
108  tc_map_type tc_components_; // cells contained by trigger cells
109  };
110 }
111 
113  public:
114 
115  typedef std::unordered_map<unsigned,std::unique_ptr<const HGCalTriggerGeometry::Module> > module_map;
116  typedef std::unordered_map<unsigned,std::unique_ptr<const HGCalTriggerGeometry::TriggerCell> > trigger_cell_map;
117 
120 
121  // non-const access to the geometry class
122  virtual void reset() override final;
123 
124  virtual unsigned getTriggerCellFromCell( const unsigned cell_det_id ) const override final;
125  virtual unsigned getModuleFromCell( const unsigned cell_det_id ) const override final;
126  virtual unsigned getModuleFromTriggerCell( const unsigned trigger_cell_det_id ) const override final;
127 
128  virtual geom_set getCellsFromTriggerCell( const unsigned cell_det_id ) const override final;
129  virtual geom_set getCellsFromModule( const unsigned cell_det_id ) const override final;
130  virtual geom_set getTriggerCellsFromModule( const unsigned trigger_cell_det_id ) const override final;
131 
132  virtual geom_ordered_set getOrderedCellsFromModule( const unsigned cell_det_id ) const override final;
133  virtual geom_ordered_set getOrderedTriggerCellsFromModule( const unsigned trigger_cell_det_id ) const override final;
134 
135  virtual geom_set getNeighborsFromTriggerCell( const unsigned trigger_cell_det_id ) const override final;
136 
137  virtual GlobalPoint getTriggerCellPosition(const unsigned trigger_cell_det_id) const override final;
138  virtual GlobalPoint getModulePosition(const unsigned module_det_id) const override final;
139 
140  bool validTriggerCell( const unsigned trigger_cell_det_id ) const final;
141  bool disconnectedModule(const unsigned module_id) const final;
142  unsigned triggerLayer(const unsigned id) const final;
143 
144  protected:
147 
148  module_map modules_;
149  trigger_cell_map trigger_cells_;
150 
151 };
152 
153 
154 #endif
std::unordered_set< unsigned > list_type
const tc_map_type & triggerCellComponents() const
Module(unsigned mod_id, const GlobalPoint &pos, const list_type &neighbs, const list_type &comps, const tc_map_type &tc_comps)
std::unordered_multimap< unsigned, unsigned > tc_map_type
Definition: value.py:1
TriggerCell(unsigned tc_id, unsigned mod_id, const GlobalPoint &pos, const list_type &neighbs, const list_type &comps)
const std::unordered_set< unsigned > & neighbours() const
std::unordered_map< unsigned, std::unique_ptr< const HGCalTriggerGeometry::TriggerCell > > trigger_cell_map
std::unordered_map< unsigned, std::unique_ptr< const HGCalTriggerGeometry::Module > > module_map
bool containsTriggerCell(const unsigned trig_cell) const
std::set< unsigned > geom_ordered_set
std::unordered_set< unsigned > geom_set
bool containsCell(const unsigned cell) const
void reset(double vett[256])
Definition: TPedValues.cc:11
const std::unordered_set< unsigned > & components() const