CMS 3D CMS Logo

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