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 
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 
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 
edmNew::DetSetVector::FastFiller::abort
void abort()
Definition: DetSetVectorNew.h:236
mps_fire.i
i
Definition: mps_fire.py:428
MessageLogger.h
TrackerGeometry.h
GeomDet
Definition: GeomDet.h:27
Phase2TrackerClusterizer::Phase2TrackerClusterizer
Phase2TrackerClusterizer(const edm::ParameterSet &conf)
Definition: Phase2TrackerClusterizer.cc:47
ESHandle.h
Phase2TrackerClusterizer::~Phase2TrackerClusterizer
~Phase2TrackerClusterizer() override
Definition: Phase2TrackerClusterizer.cc:57
Phase2TrackerClusterizer::produce
void produce(edm::Event &event, const edm::EventSetup &eventSetup) override
Definition: Phase2TrackerClusterizer.cc:63
edm::EDGetTokenT
Definition: EDGetToken.h:33
edm
HLT enums.
Definition: AlignableModifier.h:19
gather_cfg.cout
cout
Definition: gather_cfg.py:144
HLT_FULL_cff.InputTag
InputTag
Definition: HLT_FULL_cff.py:89281
Phase2TrackerClusterizerSequentialAlgorithm.h
cms::cuda::assert
assert(be >=bs)
EDProducer.h
hgcal_conditions::parameters
Definition: HGCConditions.h:86
edm::Handle
Definition: AssociativeIterator.h:50
TrackerGeometry::idToDetUnit
const TrackerGeomDet * idToDetUnit(DetId) const override
Return the pointer to the GeomDetUnit corresponding to a given DetId.
Definition: TrackerGeometry.cc:183
GetRecoTauVFromDQM_MC_cff.cl
cl
Definition: GetRecoTauVFromDQM_MC_cff.py:38
DetId
Definition: DetId.h:17
cmsdt::algo
algo
Definition: constants.h:165
MakerMacros.h
TrackerTopology.h
edm::EventSetup::get
T get() const
Definition: EventSetup.h:87
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
PixelGeomDetUnit
Definition: PixelGeomDetUnit.h:15
TrackerDigiGeometryRecord
Definition: TrackerDigiGeometryRecord.h:15
edm::ESHandle< TrackerGeometry >
slimmedTrackExtras_cff.outputClusters
outputClusters
Definition: slimmedTrackExtras_cff.py:14
Phase2TrackerClusterizerSequentialAlgorithm
Definition: Phase2TrackerClusterizerSequentialAlgorithm.h:9
mitigatedMETSequence_cff.U
U
Definition: mitigatedMETSequence_cff.py:36
bsc_activity_cfg.clusters
clusters
Definition: bsc_activity_cfg.py:36
TrackerDigiGeometryRecord.h
Phase2TrackerDigi.h
edm::ParameterSet
Definition: ParameterSet.h:47
Event.h
Phase2TrackerClusterizerAlgorithm
Definition: Phase2TrackerClusterizerAlgorithm.h:13
ModuleDef.h
createfilelist.int
int
Definition: createfilelist.py:10
edm::stream::EDProducer
Definition: EDProducer.h:38
edm::EventSetup
Definition: EventSetup.h:58
DetSetVector.h
get
#define get
Phase2TrackerClusterizer::token_
edm::EDGetTokenT< edm::DetSetVector< Phase2TrackerDigi > > token_
Definition: Phase2TrackerClusterizer.cc:40
Phase2TrackerClusterizerAlgorithm.h
InputTag.h
Phase2TrackerClusterizer
Definition: Phase2TrackerClusterizer.cc:30
Phase2TrackerDigi
Definition: Phase2TrackerDigi.h:12
GeomDet.h
eostools.move
def move(src, dest)
Definition: eostools.py:511
DetId.h
PixelGeomDetUnit.h
EventSetup.h
edmNew::DetSetVector::FastFiller
Definition: DetSetVectorNew.h:202
ConsumesCollector.h
ParameterSet.h
event
Definition: event.py:1
edm::Event
Definition: Event.h:73
CommonMethods.cp
def cp(fromDir, toDir, listOfFiles, overwrite=False, smallList=False)
Definition: CommonMethods.py:192
TrackerGeometry
Definition: TrackerGeometry.h:14