CMS 3D CMS Logo

ClusterStorer.cc
Go to the documentation of this file.
2 
4 
12 // FastSim hits:
16 
17 
19 
20 namespace helper {
21 
22  // -------------------------------------------------------------
23  //FIXME (push cluster pointers...)
25  {
26  TrackingRecHit &newHit = hits[index];
27  const std::type_info &hit_type = typeid(newHit);
28  if (hit_type == typeid(SiPixelRecHit)) {
29  //std::cout << "| It is a Pixel hit !!" << std::endl;
30  pixelClusterRecords_.push_back(PixelClusterHitRecord(static_cast<SiPixelRecHit&>(newHit),
31  hits, index));
32  } else if (hit_type == typeid(SiStripRecHit1D)) {
33  //std::cout << "| It is a SiStripRecHit1D hit !!" << std::endl;
34  stripClusterRecords_.push_back(StripClusterHitRecord(static_cast<SiStripRecHit1D&>(newHit),
35  hits, index));
36  } else if (hit_type == typeid(SiStripRecHit2D)) {
37  //std::cout << "| It is a SiStripRecHit2D hit !!" << std::endl;
38  stripClusterRecords_.push_back(StripClusterHitRecord(static_cast<SiStripRecHit2D&>(newHit),
39  hits, index));
40  } else if (hit_type == typeid(SiStripMatchedRecHit2D)) {
41  //std::cout << "| It is a SiStripMatchedRecHit2D hit !!" << std::endl;
42  SiStripMatchedRecHit2D &mhit = static_cast<SiStripMatchedRecHit2D&>(newHit);
45  } else if (hit_type == typeid(ProjectedSiStripRecHit2D)) {
46  //std::cout << "| It is a ProjectedSiStripRecHit2D hit !!" << std::endl;
47  ProjectedSiStripRecHit2D &phit = static_cast<ProjectedSiStripRecHit2D&>(newHit);
49  } else if (hit_type == typeid(Phase2TrackerRecHit1D)) {
50  //FIXME:: this is just temporary solution for phase2,
51  //it is not really running in the phase2 tracking wf - yet...
52  //std::cout << "| It is a Phase2TrackerRecHit1D hit !!" << std::endl;
53  phase2OTClusterRecords_.push_back(Phase2OTClusterHitRecord(static_cast<Phase2TrackerRecHit1D&>(newHit), hits, index));
54  } else {
55  if (hit_type == typeid(FastTrackerRecHit)
56  || hit_type == typeid(FastProjectedTrackerRecHit)
57  || hit_type == typeid(FastMatchedTrackerRecHit)) {
58  //std::cout << "| It is a " << hit_type.name() << " hit !!" << std::endl;
59  // FastSim hits: Do nothing instead of caring about FastSim clusters,
60  // not even sure whether these really exist.
61  // At least predecessor code in TrackSelector and MuonSelector
62  // did not treat them.
63  } else {
64  // through for unknown types
65  throw cms::Exception("UnknownHitType") << "helper::ClusterStorer::addCluster: "
66  << "Unknown hit type " << hit_type.name()
67  << ".\n";
68  }
69  } // end 'switch' on hit type
70 
71  }
72 
73  // -------------------------------------------------------------
75  {
76  pixelClusterRecords_.clear();
77  stripClusterRecords_.clear();
78  }
79 
80  // -------------------------------------------------------------
81  void ClusterStorer::
86  {
87  if (!pixelClusterRecords_.empty()) {
88  this->processClusters<SiPixelRecHit, SiPixelCluster>
89  (pixelClusterRecords_, pixelDsvToFill, refPixelClusters);
90  }
91  if (!stripClusterRecords_.empty()) {
92  // All we need from the HitType 'SiStripRecHit2D' is the
93  // typedef for 'SiStripRecHit2D::ClusterRef'.
94  // The fact that rekey<SiStripRecHit2D> is called is irrelevant since
95  // ClusterHitRecord<typename SiStripRecHit2D::ClusterRef>::rekey<RecHitType>
96  // is specialised such that 'RecHitType' is not used...
97  this->processClusters<SiStripRecHit2D, SiStripCluster>
98  (stripClusterRecords_, stripDsvToFill, refStripClusters);
99  }
100  }
101 
102  //-------------------------------------------------------------
103  template<typename HitType, typename ClusterType>
104  void ClusterStorer::
108  {
109  std::sort(clusterRecords.begin(), clusterRecords.end()); // this sorts them by detid
110  typedef
111  typename std::vector<ClusterHitRecord<typename HitType::ClusterRef> >::const_iterator
112  RIT;
113  RIT it = clusterRecords.begin(), end = clusterRecords.end();
114  size_t clusters = 0;
115  while (it != end) {
116  RIT it2 = it;
117  uint32_t detid = it->detid();
118 
119  // first isolate all clusters on the same detid
120  while ( (it2 != end) && (it2->detid() == detid)) { ++it2; }
121  // now [it, it2] bracket one detid
122 
123  // then prepare to copy the clusters
124  typename edmNew::DetSetVector<ClusterType>::FastFiller filler(dsvToFill, detid);
125  typename HitType::ClusterRef lastRef, newRef;
126  for ( ; it != it2; ++it) { // loop on the detid
127  // first check if we need to clone the hit
128  if (it->clusterRef() != lastRef) {
129  lastRef = it->clusterRef();
130  // clone cluster
131  filler.push_back( *lastRef );
132  // make new ref
133  newRef = typename HitType::ClusterRef( refprod, clusters++ );
134  }
135  it->template rekey<HitType>(newRef);
136  } // end of the loop on a single detid
137 
138  } // end of the loop on all clusters
139  }
140 
141  //-------------------------------------------------------------
142  // helper classes
143  //-------------------------------------------------------------
144  // FIXME (migrate to new RTTI and interface)
145  // generic rekey (in practise for pixel only...)
146  template<typename ClusterRefType> // template for class
147  template<typename RecHitType> // template for member function
149  rekey(const ClusterRefType &newRef) const
150  {
151  TrackingRecHit & genericHit = (*hits_)[index_];
152  RecHitType *hit = 0;
153  if (genericHit.geographicalId().rawId() == detid_) { // a hit on this det, so it's simple
154  hit = dynamic_cast<RecHitType *>(&genericHit); //static_cast<RecHitType *>(&genericHit);
155  }
156  assert (hit != 0);
157  assert (hit->cluster() == ref_); // otherwise something went wrong
158  hit->setClusterRef(newRef);
159  }
160 
161  // -------------------------------------------------------------
162  // Specific rekey for class template ClusterRefType = SiStripRecHit2D::ClusterRef,
163  // RecHitType is not used.
164  template<>
165  template<typename RecHitType> // or template<> to specialise also here?
167  // rekey<SiStripRecHit2D>(const SiStripRecHit2D::ClusterRef &newRef) const
168  rekey(const SiStripRecHit2D::ClusterRef &newRef) const
169  {
170  TrackingRecHit &genericHit = (*hits_)[index_];
171  const std::type_info &hit_type = typeid(genericHit);
172 
173  OmniClusterRef * cluRef=0;
174  if (typeid(SiStripRecHit1D) == hit_type) {
175  cluRef = &static_cast<SiStripRecHit1D&>(genericHit).omniCluster();
176  } else if (typeid(SiStripRecHit2D) == hit_type) {
177  cluRef = &static_cast<SiStripRecHit2D&>(genericHit).omniCluster();
178  } else if (typeid(SiStripMatchedRecHit2D) == hit_type) {
179  SiStripMatchedRecHit2D &mhit = static_cast<SiStripMatchedRecHit2D&>(genericHit);
180  cluRef = (SiStripDetId(detid_).stereo() ? &mhit.stereoClusterRef() : &mhit.monoClusterRef());
181  } else if (typeid(ProjectedSiStripRecHit2D) == hit_type) {
182  cluRef = &static_cast<ProjectedSiStripRecHit2D&>(genericHit).originalHit().omniCluster();
183  }
184 
185  assert(cluRef != 0); // to catch missing RecHit types
186  assert(cluRef->key() == ref_.key()); // otherwise something went wrong
187  (*cluRef) = OmniClusterRef(newRef);
188  }
189 
190 } // end namespace 'helper'
A struct for clusters associated to hits.
Definition: ClusterStorer.h:47
void push_back(data_type const &d)
ClusterHitRecord< SiStripRecHit2D::ClusterRef > StripClusterHitRecord
Assuming that the ClusterRef is the same for all SiStripRecHit*:
Definition: ClusterStorer.h:77
Definition: helper.py:1
uint32_t stereo() const
Definition: SiStripDetId.h:160
OmniClusterRef const & stereoClusterRef() const
void rekey(const ClusterRefType &newRef) const
void addCluster(TrackingRecHitCollection &hits, size_t index)
add cluster of newHit to list (throws if hit is of unknown type)
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
ClusterHitRecord< Phase2TrackerRecHit1D::CluRef > Phase2OTClusterHitRecord
Definition: ClusterStorer.h:80
OmniClusterRef const & omniCluster() const
std::vector< Phase2OTClusterHitRecord > phase2OTClusterRecords_
Definition: ClusterStorer.h:95
#define end
Definition: vmac.h:37
ClusterHitRecord< SiPixelRecHit::ClusterRef > PixelClusterHitRecord
Definition: ClusterStorer.h:75
OmniClusterRef const & monoClusterRef() const
std::vector< PixelClusterHitRecord > pixelClusterRecords_
Definition: ClusterStorer.h:93
SiStripRecHit2D originalHit() const
Detector identifier class for the strip tracker.
Definition: SiStripDetId.h:17
SiStripRecHit2D stereoHit() const
void clear()
clear records
SiStripRecHit2D monoHit() const
void processClusters(std::vector< ClusterHitRecord< typename HitType::ClusterRef > > &clusterRecords, edmNew::DetSetVector< ClusterType > &dsvToFill, edm::RefProd< edmNew::DetSetVector< ClusterType > > &refprod)
void processAllClusters(edmNew::DetSetVector< SiPixelCluster > &pixelDsvToFill, edm::RefProd< edmNew::DetSetVector< SiPixelCluster > > refPixelClusters, edmNew::DetSetVector< SiStripCluster > &stripDsvToFill, edm::RefProd< edmNew::DetSetVector< SiStripCluster > > refStripClusters)
DetId geographicalId() const
unsigned int key() const
Our base class.
Definition: SiPixelRecHit.h:23
std::vector< StripClusterHitRecord > stripClusterRecords_
Definition: ClusterStorer.h:94