CMS 3D CMS Logo

SiPixelRecHitConverter.cc
Go to the documentation of this file.
1 
12 //---------------------------------------------------------------------------
40 //---------------------------------------------------------------------------
41 
42 //--- Base class for CPEs:
43 
45 
46 //--- Geometry + DataFormats
51 
52 //--- Framework
57 
64 
65 // Geometry
68 
69 // Data Formats
73 
74 // STL
75 #include <vector>
76 #include <memory>
77 #include <string>
78 #include <iostream>
79 
80 // MessageLogger
82 
84 
85 using namespace std;
86 
87 namespace cms {
88 
90  public:
91  //--- Constructor, virtual destructor (just in case)
92  explicit SiPixelRecHitConverter(const edm::ParameterSet& conf);
93  ~SiPixelRecHitConverter() override;
94 
95  //--- Factory method to make CPE's depending on the ParameterSet
96  //--- Not sure if we need to make more than one CPE to run concurrently
97  //--- on different parts of the detector (e.g., one for the barrel and the
98  //--- one for the forward). The way the CPE's are written now, it's
99  //--- likely we can use one (and they will switch internally), or
100  //--- make two of the same but configure them differently. We need a more
101  //--- realistic use case...
102 
103  //--- The top-level event method.
104  void produce(edm::Event& e, const edm::EventSetup& c) override;
105 
106  //--- Execute the position estimator algorithm(s).
107  //--- New interface with DetSetVector
110  TrackerGeometry const& geom);
111 
114  TrackerGeometry const& geom);
115 
116  private:
117  // TO DO: maybe allow a map of pointers?
119  PixelCPEBase const* cpe_ = nullptr; // What we got (for now, one ptr to base class)
125  bool m_newCont; // save also in emdNew::DetSetVector
126  };
127 
128  //---------------------------------------------------------------------------
130  //---------------------------------------------------------------------------
131  SiPixelRecHitConverter::SiPixelRecHitConverter(edm::ParameterSet const& conf)
132  : src_(conf.getParameter<edm::InputTag>("src")),
133  tPixelCluster_(consumes<edmNew::DetSetVector<SiPixelCluster>>(src_)),
134  tPut_(produces<SiPixelRecHitCollection>()),
135  tTrackerGeom_(esConsumes<TrackerGeometry, TrackerDigiGeometryRecord>()),
137  edm::ESInputTag("", conf.getParameter<std::string>("CPE")))) {}
138 
139  // Destructor
141 
142  //---------------------------------------------------------------------------
144  //---------------------------------------------------------------------------
146  // Step A.1: get input data
148  e.getByToken(tPixelCluster_, input);
149 
150  // Step A.2: get event setup
151  auto const& geom = es.getData(tTrackerGeom_);
152 
153  // Step B: create empty output collection
155 
156  // Step B*: create CPE
157  cpe_ = dynamic_cast<const PixelCPEBase*>(&es.getData(tCPE_));
158 
159  // Step C: Iterate over DetIds and invoke the strip CPE algorithm
160  // on each DetUnit
161 
162  run(input, output, geom);
163 
164  output.shrink_to_fit();
165  e.emplace(tPut_, std::move(output));
166  }
167 
168  //---------------------------------------------------------------------------
172  //---------------------------------------------------------------------------
175  TrackerGeometry const& geom) {
176  if (!cpe_) {
177  edm::LogError("SiPixelRecHitConverter") << " at least one CPE is not ready -- can't run!";
178  // TO DO: throw an exception here? The user may want to know...
179  assert(0);
180  return; // clusterizer is invalid, bail out
181  }
182 
183  int numberOfDetUnits = 0;
184  int numberOfClusters = 0;
185 
186  const edmNew::DetSetVector<SiPixelCluster>& input = *inputhandle;
187 
189 
190  for (; DSViter != input.end(); DSViter++) {
191  numberOfDetUnits++;
192  unsigned int detid = DSViter->detId();
193  DetId detIdObject(detid);
194  const GeomDetUnit* genericDet = geom.idToDetUnit(detIdObject);
195  const PixelGeomDetUnit* pixDet = dynamic_cast<const PixelGeomDetUnit*>(genericDet);
196  assert(pixDet);
197  SiPixelRecHitCollectionNew::FastFiller recHitsOnDetUnit(output, detid);
198 
199  edmNew::DetSet<SiPixelCluster>::const_iterator clustIt = DSViter->begin(), clustEnd = DSViter->end();
200 
201  for (; clustIt != clustEnd; clustIt++) {
202  numberOfClusters++;
203  std::tuple<LocalPoint, LocalError, SiPixelRecHitQuality::QualWordType> tuple =
204  cpe_->getParameters(*clustIt, *genericDet);
205  LocalPoint lp(std::get<0>(tuple));
206  LocalError le(std::get<1>(tuple));
207  SiPixelRecHitQuality::QualWordType rqw(std::get<2>(tuple));
208  // Create a persistent edm::Ref to the cluster
210  edmNew::makeRefTo(inputhandle, clustIt);
211  // Make a RecHit and add it to the DetSet
212  // old : recHitsOnDetUnit.push_back( new SiPixelRecHit( lp, le, detIdObject, &*clustIt) );
213  SiPixelRecHit hit(lp, le, rqw, *genericDet, cluster);
214  //
215  // Now save it =================
216  recHitsOnDetUnit.push_back(hit);
217  // =============================
218 
219  // std::cout << "SiPixelRecHitConverterVI " << numberOfClusters << ' '<< lp << " " << le << std::endl;
220  } // <-- End loop on Clusters
221 
222  // LogDebug("SiPixelRecHitConverter")
223  //std::cout << "SiPixelRecHitConverterVI "
224  // << " Found " << recHitsOnDetUnit.size() << " RecHits on " << detid //;
225  // << std::endl;
226 
227  } // <-- End loop on DetUnits
228 
229  // LogDebug ("SiPixelRecHitConverter")
230  // std::cout << "SiPixelRecHitConverterVI "
231  // << cpeName_ << " converted " << numberOfClusters
232  // << " SiPixelClusters into SiPixelRecHits, in "
233  // << numberOfDetUnits << " DetUnits." //;
234  // << std::endl;
235  }
236 } // end of namespace cms
237 
239 
PixelClusterParameterEstimator
Definition: PixelClusterParameterEstimator.h:15
PixelCPEBase::getParameters
ReturnType getParameters(const SiPixelCluster &cl, const GeomDetUnit &det) const override
Definition: PixelCPEBase.h:143
Handle.h
input
static const std::string input
Definition: EdmProvDump.cc:48
MessageLogger.h
cms::SiPixelRecHitConverter::src_
const edm::InputTag src_
Definition: SiPixelRecHitConverter.cc:120
TrackerGeometry.h
cms::SiPixelRecHitConverter::cpe_
PixelCPEBase const * cpe_
const PixelClusterParameterEstimator * cpe_; // what we got (for now, one ptr to base class)
Definition: SiPixelRecHitConverter.cc:119
GeomDet
Definition: GeomDet.h:27
cms::SiPixelRecHitConverter::tPut_
const edm::EDPutTokenT< SiPixelRecHitCollection > tPut_
Definition: SiPixelRecHitConverter.cc:122
ESHandle.h
ESInputTag
convertSQLitetoXML_cfg.output
output
Definition: convertSQLitetoXML_cfg.py:72
edm::EDGetTokenT
Definition: EDGetToken.h:33
edm
HLT enums.
Definition: AlignableModifier.h:19
edm::EDPutTokenT
Definition: EDPutToken.h:33
SiPixelCluster.h
HLT_FULL_cff.InputTag
InputTag
Definition: HLT_FULL_cff.py:85964
edmNew::makeRefTo
edm::Ref< typename HandleT::element_type, typename HandleT::element_type::value_type::value_type > makeRefTo(const HandleT &iHandle, typename HandleT::element_type::value_type::const_iterator itIter)
Definition: DetSetVectorNew.h:689
edmNew::DetSetVector::const_iterator
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
Definition: DetSetVectorNew.h:197
TkPixelCPERecord
Definition: TkPixelCPERecord.h:18
cms::cuda::assert
assert(be >=bs)
EDProducer.h
PixelCPEBase.h
edmNew::DetSetVector::begin
const_iterator begin(bool update=false) const
Definition: DetSetVectorNew.h:530
SiPixelCluster
Pixel cluster – collection of neighboring pixels above threshold.
Definition: SiPixelCluster.h:27
edmNew::DetSetVector::FastFiller::push_back
void push_back(data_type const &d)
Definition: DetSetVectorNew.h:274
edm::Handle
Definition: AssociativeIterator.h:50
ESGetToken.h
SiPixelRecHit
Our base class.
Definition: SiPixelRecHit.h:23
edmNew
Definition: DetSet2RangeMap.h:11
edm::Ref
Definition: AssociativeIterator.h:58
DetId
Definition: DetId.h:17
MakerMacros.h
cms::SiPixelRecHitConverter::tCPE_
const edm::ESGetToken< PixelClusterParameterEstimator, TkPixelCPERecord > tCPE_
Definition: SiPixelRecHitConverter.cc:124
cms::SiPixelRecHitConverter::~SiPixelRecHitConverter
~SiPixelRecHitConverter() override
Definition: SiPixelRecHitConverter.cc:140
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
SiPixelRecHitQuality::QualWordType
unsigned int QualWordType
Definition: SiPixelRecHitQuality.h:10
cms::SiPixelRecHitConverter::run
void run(const edmNew::DetSetVector< SiPixelCluster > &input, SiPixelRecHitCollectionNew &output, TrackerGeometry const &geom)
PixelGeomDetUnit
Definition: PixelGeomDetUnit.h:15
TrackerDigiGeometryRecord
Definition: TrackerDigiGeometryRecord.h:15
relativeConstraints.geom
geom
Definition: relativeConstraints.py:72
Point3DBase< float, LocalTag >
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
TrackerDigiGeometryRecord.h
SiPixelRecHitCollection.h
edm::ParameterSet
Definition: ParameterSet.h:47
Event.h
EDPutToken.h
LocalError
Definition: LocalError.h:12
edm::stream::EDProducer
Definition: EDProducer.h:38
edm::EventSetup
Definition: EventSetup.h:57
cms::SiPixelRecHitConverter
Definition: SiPixelRecHitConverter.cc:89
DetSetVector.h
HltBtagPostValidation_cff.c
c
Definition: HltBtagPostValidation_cff.py:31
edm::LogError
Log< level::Error, false > LogError
Definition: MessageLogger.h:123
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord >
InputTag.h
edm::EventSetup::getData
bool getData(T &iHolder) const
Definition: EventSetup.h:120
edmNew::DetSetVector
Definition: DetSetNew.h:13
cms::SiPixelRecHitConverter::produce
void produce(edm::Event &e, const edm::EventSetup &c) override
The "Event" entrypoint: gets called by framework for every event.
Definition: SiPixelRecHitConverter.cc:145
eostools.move
def move(src, dest)
Definition: eostools.py:511
std
Definition: JetResolutionObject.h:76
writedatasetfile.run
run
Definition: writedatasetfile.py:27
Ref.h
DetId.h
cms::SiPixelRecHitConverter::tPixelCluster_
const edm::EDGetTokenT< edmNew::DetSetVector< SiPixelCluster > > tPixelCluster_
Definition: SiPixelRecHitConverter.cc:121
edmNew::DetSetVector::end
const_iterator end(bool update=false) const
Definition: DetSetVectorNew.h:535
PixelCPEBase
Definition: PixelCPEBase.h:54
PixelGeomDetUnit.h
EventSetup.h
DetSet2RangeMap.h
TkPixelCPERecord.h
cms::SiPixelRecHitConverter::m_newCont
bool m_newCont
Definition: SiPixelRecHitConverter.cc:125
edmNew::DetSetVector::FastFiller
Definition: DetSetVectorNew.h:202
ParameterSet.h
edm::Event
Definition: Event.h:73
edm::InputTag
Definition: InputTag.h:15
hit
Definition: SiStripHitEffFromCalibTree.cc:88
cms::SiPixelRecHitConverter::tTrackerGeom_
const edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > tTrackerGeom_
Definition: SiPixelRecHitConverter.cc:123
cms
Namespace of DDCMS conversion namespace.
Definition: ProducerAnalyzer.cc:21
MillePedeFileConverter_cfg.e
e
Definition: MillePedeFileConverter_cfg.py:37
TrackerGeometry
Definition: TrackerGeometry.h:14
edmNew::DetSet::const_iterator
const data_type * const_iterator
Definition: DetSetNew.h:31