CMS 3D CMS Logo

Phase2TrackerClusterizer.cc
Go to the documentation of this file.
11 
12 #ifdef VERIFY_PH2_TK_CLUS
14 #endif
16 
21 
26 
27 #include <vector>
28 #include <memory>
29 
30 
31 
33 
34  public:
35  explicit Phase2TrackerClusterizer(const edm::ParameterSet& conf);
36  ~Phase2TrackerClusterizer() override;
37  void produce(edm::Event& event, const edm::EventSetup& eventSetup) override;
38 
39  private:
40 #ifdef VERIFY_PH2_TK_CLUS
41  std::unique_ptr< Phase2TrackerClusterizerAlgorithm > clusterizer_;
42 #endif
44 
45 };
46 
47 
48  /*
49  * Initialise the producer
50  */
51 
53 #ifdef VERIFY_PH2_TK_CLUS
54  clusterizer_(new Phase2TrackerClusterizerAlgorithm(conf.getParameter< unsigned int >("maxClusterSize"), conf.getParameter< unsigned int >("maxNumberClusters"))),
55 #endif
56  token_(consumes< edm::DetSetVector< Phase2TrackerDigi > >(conf.getParameter<edm::InputTag>("src"))) {
57  produces< Phase2TrackerCluster1DCollectionNew >();
58  }
59 
61 
62  /*
63  * Clusterize the events
64  */
65 
67 
68  // Get the Digis
70  event.getByToken(token_, digis);
71 
72 #ifdef VERIFY_PH2_TK_CLUS
73  // Get the geometry
75  eventSetup.get< TrackerDigiGeometryRecord >().get(geomHandle);
76  const TrackerGeometry* tkGeom(&(*geomHandle));
77  // Global container for the clusters of each modules
78  auto outputClustersOld = std::make_unique<Phase2TrackerCluster1DCollectionNew>();
79 #endif
80  auto outputClusters = std::make_unique<Phase2TrackerCluster1DCollectionNew>();
81 
82  // Go over all the modules
83  for (auto DSViter : *digis) {
84  DetId detId(DSViter.detId());
85 
86  Phase2TrackerCluster1DCollectionNew::FastFiller clusters(*outputClusters, DSViter.detId());
88  algo.clusterizeDetUnit(DSViter, clusters);
89  if (clusters.empty()) clusters.abort();
90 
91 #ifdef VERIFY_PH2_TK_CLUS
92  if (!clusters.empty()) {
93  auto cp = clusters[0].column();
94  auto sp = clusters[0].firstStrip();
95  for (auto const & cl : clusters) {
96  if (cl.column()<cp) std::cout << "column not in order! " << std::endl;
97  if (cl.column()==cp && cl.firstStrip()<sp) std::cout << "strip not in order! " << std::endl;
98  cp = cl.column();
99  sp = cl.firstStrip();
100  }
101  }
102 #endif
103 
104 #ifdef VERIFY_PH2_TK_CLUS
105  // Geometry
106  const GeomDetUnit* geomDetUnit(tkGeom->idToDetUnit(detId));
107  const PixelGeomDetUnit* pixDet = dynamic_cast< const PixelGeomDetUnit* >(geomDetUnit);
108  if (!pixDet) assert(0);
109 
110  // Container for the clusters that will be produced for this modules
111  Phase2TrackerCluster1DCollectionNew::FastFiller clustersOld(*outputClustersOld, DSViter.detId());
112 
113  // Setup the clusterizer algorithm for this detector (see ClusterizerAlgorithm for more details)
114  clusterizer_->setup(pixDet);
115 
116  // Pass the list of Digis to the main algorithm
117  // This function will store the clusters in the previously created container
118  clusterizer_->clusterizeDetUnit(DSViter, clustersOld);
119  if (clustersOld.empty()) clustersOld.abort();
120 
121  if (clusters.size() != clustersOld.size()) {
122  std::cout << "SIZEs " << int(detId) << ' ' << clusters.size() << ' ' << clustersOld.size() << std::endl;
123  for (auto const & cl : clusters) std::cout << cl.size() << ' ' << cl.threshold() << ' ' << cl.firstRow() << ' ' << cl.column() << std::endl;
124  std::cout << "Old " << std::endl;
125  for (auto const & cl : clustersOld) std::cout << cl.size() << ' ' << cl.threshold() << ' ' << cl.firstRow() << ' ' << cl.column() << std::endl;
126  }
127 #endif
128 
129  }
130 
131 #ifdef VERIFY_PH2_TK_CLUS
132  // std::cout << "SIZEs " << outputClusters->dataSize() << ' ' << outputClustersOld->dataSize() << std::endl;
133  assert(outputClusters->dataSize()==outputClustersOld->dataSize());
134  for (auto i=0U;i<outputClusters->dataSize(); ++i) {
135  assert(outputClusters->data()[i].size() == outputClustersOld->data()[i].size());
136  assert(outputClusters->data()[i].threshold() == outputClustersOld->data()[i].threshold());
137  assert(outputClusters->data()[i].firstRow() == outputClustersOld->data()[i].firstRow());
138  assert(outputClusters->data()[i].column() == outputClustersOld->data()[i].column());
139  }
140 #endif
141 
142  // Add the data to the output
143  outputClusters->shrink_to_fit();
144  event.put(std::move(outputClusters));
145  }
146 
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
edm::EDGetTokenT< edm::DetSetVector< Phase2TrackerDigi > > token_
Phase2TrackerClusterizer(const edm::ParameterSet &conf)
const TrackerGeomDet * idToDetUnit(DetId) const override
Return the pointer to the GeomDetUnit corresponding to a given DetId.
void produce(edm::Event &event, const edm::EventSetup &eventSetup) override
Definition: DetId.h:18
HLT enums.
T get() const
Definition: EventSetup.h:63
def move(src, dest)
Definition: eostools.py:510
Definition: event.py:1