CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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 
31 public:
32  explicit Phase2TrackerClusterizer(const edm::ParameterSet& conf);
33  ~Phase2TrackerClusterizer() override;
34  void produce(edm::Event& event, const edm::EventSetup& eventSetup) override;
35 
36 private:
37 #ifdef VERIFY_PH2_TK_CLUS
38  std::unique_ptr<Phase2TrackerClusterizerAlgorithm> clusterizer_;
39 #endif
41 };
42 
43 /*
44  * Initialise the producer
45  */
46 
48  :
49 #ifdef VERIFY_PH2_TK_CLUS
50  clusterizer_(new Phase2TrackerClusterizerAlgorithm(conf.getParameter<unsigned int>("maxClusterSize"),
51  conf.getParameter<unsigned int>("maxNumberClusters"))),
52 #endif
53  token_(consumes<edm::DetSetVector<Phase2TrackerDigi> >(conf.getParameter<edm::InputTag>("src"))) {
54  produces<Phase2TrackerCluster1DCollectionNew>();
55 }
56 
58 
59 /*
60  * Clusterize the events
61  */
62 
64  // Get the Digis
66  event.getByToken(token_, digis);
67 
68 #ifdef VERIFY_PH2_TK_CLUS
69  // Get the geometry
71  eventSetup.get<TrackerDigiGeometryRecord>().get(geomHandle);
72  const TrackerGeometry* tkGeom(&(*geomHandle));
73  // Global container for the clusters of each modules
74  auto outputClustersOld = std::make_unique<Phase2TrackerCluster1DCollectionNew>();
75 #endif
76  auto outputClusters = std::make_unique<Phase2TrackerCluster1DCollectionNew>();
77 
78  // Go over all the modules
79  for (const auto& DSViter : *digis) {
80  DetId detId(DSViter.detId());
81 
82  Phase2TrackerCluster1DCollectionNew::FastFiller clusters(*outputClusters, DSViter.detId());
84  algo.clusterizeDetUnit(DSViter, clusters);
85  if (clusters.empty())
86  clusters.abort();
87 
88 #ifdef VERIFY_PH2_TK_CLUS
89  if (!clusters.empty()) {
90  auto cp = clusters[0].column();
91  auto sp = clusters[0].firstStrip();
92  for (auto const& cl : clusters) {
93  if (cl.column() < cp)
94  std::cout << "column not in order! " << std::endl;
95  if (cl.column() == cp && cl.firstStrip() < sp)
96  std::cout << "strip not in order! " << std::endl;
97  cp = cl.column();
98  sp = cl.firstStrip();
99  }
100  }
101 #endif
102 
103 #ifdef VERIFY_PH2_TK_CLUS
104  // Geometry
105  const GeomDetUnit* geomDetUnit(tkGeom->idToDetUnit(detId));
106  const PixelGeomDetUnit* pixDet = dynamic_cast<const PixelGeomDetUnit*>(geomDetUnit);
107  if (!pixDet)
108  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())
120  clustersOld.abort();
121 
122  if (clusters.size() != clustersOld.size()) {
123  std::cout << "SIZEs " << int(detId) << ' ' << clusters.size() << ' ' << clustersOld.size() << std::endl;
124  for (auto const& cl : clusters)
125  std::cout << cl.size() << ' ' << cl.threshold() << ' ' << cl.firstRow() << ' ' << cl.column() << std::endl;
126  std::cout << "Old " << std::endl;
127  for (auto const& cl : clustersOld)
128  std::cout << cl.size() << ' ' << cl.threshold() << ' ' << cl.firstRow() << ' ' << cl.column() << std::endl;
129  }
130 #endif
131  }
132 
133 #ifdef VERIFY_PH2_TK_CLUS
134  // std::cout << "SIZEs " << outputClusters->dataSize() << ' ' << outputClustersOld->dataSize() << std::endl;
135  assert(outputClusters->dataSize() == outputClustersOld->dataSize());
136  for (auto i = 0U; i < outputClusters->dataSize(); ++i) {
137  assert(outputClusters->data()[i].size() == outputClustersOld->data()[i].size());
138  assert(outputClusters->data()[i].threshold() == outputClustersOld->data()[i].threshold());
139  assert(outputClusters->data()[i].firstRow() == outputClustersOld->data()[i].firstRow());
140  assert(outputClusters->data()[i].column() == outputClustersOld->data()[i].column());
141  }
142 #endif
143 
144  // Add the data to the output
145  outputClusters->shrink_to_fit();
146  event.put(std::move(outputClusters));
147 }
148 
const TrackerGeomDet * idToDetUnit(DetId) const override
Return the pointer to the GeomDetUnit corresponding to a given DetId.
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
assert(be >=bs)
Phase2TrackerClusterizer(const edm::ParameterSet &conf)
tuple cl
Definition: haddnano.py:49
def move
Definition: eostools.py:511
void produce(edm::Event &event, const edm::EventSetup &eventSetup) override
Definition: DetId.h:17
T get() const
Definition: EventSetup.h:82
tuple cout
Definition: gather_cfg.py:144
edm::EDGetTokenT< edm::DetSetVector< Phase2TrackerDigi > > token_