Line Code
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207
#include "CommonTools/RecoAlgos/interface/ClusterStorer.h"

#include "FWCore/Utilities/interface/Exception.h"

#include "DataFormats/TrackingRecHit/interface/TrackingRecHit.h"
#include "DataFormats/TrackerRecHit2D/interface/SiPixelRecHit.h"
#include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit1D.h"
#include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit2D.h"
#include "DataFormats/TrackerRecHit2D/interface/ProjectedSiStripRecHit2D.h"
#include "DataFormats/TrackerRecHit2D/interface/SiStripMatchedRecHit2D.h"
#include "DataFormats/TrackerRecHit2D/interface/Phase2TrackerRecHit1D.h"
#include "DataFormats/TrackerRecHit2D/interface/VectorHit.h"
// FastSim hits:
#include "DataFormats/TrackerRecHit2D/interface/FastTrackerRecHit.h"
#include "DataFormats/TrackerRecHit2D/interface/FastProjectedTrackerRecHit.h"
#include "DataFormats/TrackerRecHit2D/interface/FastMatchedTrackerRecHit.h"

#include "DataFormats/SiStripDetId/interface/SiStripDetId.h"

namespace helper {

  // -------------------------------------------------------------
  //FIXME (push cluster pointers...)
  void ClusterStorer::addCluster(TrackingRecHitCollection &hits, size_t index) {
    TrackingRecHit &newHit = hits[index];
    const std::type_info &hit_type = typeid(newHit);
    if (hit_type == typeid(SiPixelRecHit)) {
      //std::cout << "|  It is a Pixel hit !!" << std::endl;
      pixelClusterRecords_.push_back(PixelClusterHitRecord(static_cast<SiPixelRecHit &>(newHit), hits, index));
    } else if (hit_type == typeid(SiStripRecHit1D)) {
      //std::cout << "|   It is a SiStripRecHit1D hit !!" << std::endl;
      stripClusterRecords_.push_back(StripClusterHitRecord(static_cast<SiStripRecHit1D &>(newHit), hits, index));
    } else if (hit_type == typeid(SiStripRecHit2D)) {
      //std::cout << "|   It is a SiStripRecHit2D hit !!" << std::endl;
      stripClusterRecords_.push_back(StripClusterHitRecord(static_cast<SiStripRecHit2D &>(newHit), hits, index));
    } else if (hit_type == typeid(SiStripMatchedRecHit2D)) {
      //std::cout << "|   It is a SiStripMatchedRecHit2D hit !!" << std::endl;
      SiStripMatchedRecHit2D &mhit = static_cast<SiStripMatchedRecHit2D &>(newHit);
      stripClusterRecords_.push_back(StripClusterHitRecord(mhit.monoHit(), hits, index));
      stripClusterRecords_.push_back(StripClusterHitRecord(mhit.stereoHit(), hits, index));
    } else if (hit_type == typeid(ProjectedSiStripRecHit2D)) {
      //std::cout << "|   It is a ProjectedSiStripRecHit2D hit !!" << std::endl;
      ProjectedSiStripRecHit2D &phit = static_cast<ProjectedSiStripRecHit2D &>(newHit);
      stripClusterRecords_.push_back(StripClusterHitRecord(phit.originalHit(), hits, index));
    } else if (hit_type == typeid(Phase2TrackerRecHit1D)) {
      //FIXME:: this is just temporary solution for phase2,
      //it is not really running in the phase2 tracking wf - yet...
      phase2OTClusterRecords_.push_back(
          Phase2OTClusterHitRecord(static_cast<Phase2TrackerRecHit1D &>(newHit), hits, index));
    } else if (hit_type == typeid(VectorHit)) {
      //FIXME:: this is just temporary solution for phase2,
      //the VectorHit has 2 clusters but just a hit!
      phase2OTClusterRecords_.push_back(Phase2OTClusterHitRecord(static_cast<VectorHit &>(newHit), hits, index));
    } else {
      if (hit_type == typeid(FastTrackerRecHit) || hit_type == typeid(FastProjectedTrackerRecHit) ||
          hit_type == typeid(FastMatchedTrackerRecHit)) {
        //std::cout << "|   It is a " << hit_type.name() << " hit !!" << std::endl;
        // FastSim hits: Do nothing instead of caring about FastSim clusters,
        //               not even sure whether these really exist.
        //               At least predecessor code in TrackSelector and MuonSelector
        //               did not treat them.
      } else {
        // through for unknown types
        throw cms::Exception("UnknownHitType") << "helper::ClusterStorer::addCluster: "
                                               << "Unknown hit type " << hit_type.name() << ".\n";
      }
    }  // end 'switch' on hit type
  }

  // -------------------------------------------------------------
  void ClusterStorer::clear() {
    pixelClusterRecords_.clear();
    stripClusterRecords_.clear();
    phase2OTClusterRecords_.clear();
  }

  // -------------------------------------------------------------
  void ClusterStorer::processAllClusters(
      edmNew::DetSetVector<SiPixelCluster> &pixelDsvToFill,
      edm::RefProd<edmNew::DetSetVector<SiPixelCluster> > refPixelClusters,
      edmNew::DetSetVector<SiStripCluster> &stripDsvToFill,
      edm::RefProd<edmNew::DetSetVector<SiStripCluster> > refStripClusters,
      edmNew::DetSetVector<Phase2TrackerCluster1D> &phase2OTDsvToFill,
      edm::RefProd<edmNew::DetSetVector<Phase2TrackerCluster1D> > refPhase2OTClusters) {
    if (!pixelClusterRecords_.empty()) {
      this->processClusters<SiPixelRecHit, SiPixelCluster>(pixelClusterRecords_, pixelDsvToFill, refPixelClusters);
    }
    if (!stripClusterRecords_.empty()) {
      // All we need from the HitType 'SiStripRecHit2D' is the
      // typedef for 'SiStripRecHit2D::ClusterRef'.
      // The fact that rekey<SiStripRecHit2D> is called is irrelevant since
      // ClusterHitRecord<typename SiStripRecHit2D::ClusterRef>::rekey<RecHitType>
      // is specialised such that 'RecHitType' is not used...
      this->processClusters<SiStripRecHit2D, SiStripCluster>(stripClusterRecords_, stripDsvToFill, refStripClusters);
    }
    if (!phase2OTClusterRecords_.empty()) {
      this->processClusters<Phase2TrackerRecHit1D, Phase2TrackerCluster1D>(
          phase2OTClusterRecords_, phase2OTDsvToFill, refPhase2OTClusters);
    }
  }

  //-------------------------------------------------------------
  template <typename HitType, typename ClusterType>
  void ClusterStorer::processClusters(std::vector<ClusterHitRecord<typename HitType::ClusterRef> > &clusterRecords,
                                      edmNew::DetSetVector<ClusterType> &dsvToFill,
                                      edm::RefProd<edmNew::DetSetVector<ClusterType> > &refprod) {
    std::sort(clusterRecords.begin(), clusterRecords.end());  // this sorts them by detid
    typedef typename std::vector<ClusterHitRecord<typename HitType::ClusterRef> >::iterator RIT;
    RIT it = clusterRecords.begin(), end = clusterRecords.end();
    size_t clusters = 0;
    while (it != end) {
      RIT it2 = it;
      uint32_t detid = it->detid();

      // first isolate all clusters on the same detid
      while ((it2 != end) && (it2->detid() == detid)) {
        ++it2;
      }
      // now [it, it2] bracket one detid

      // then prepare to copy the clusters
      typename edmNew::DetSetVector<ClusterType>::FastFiller filler(dsvToFill, detid);
      typename HitType::ClusterRef lastRef, newRef;
      for (; it != it2; ++it) {  // loop on the detid
        // first check if we need to clone the hit
        if (it->clusterRef() != lastRef) {
          lastRef = it->clusterRef();
          // clone cluster
          filler.push_back(*lastRef);
          // make new ref
          newRef = typename HitType::ClusterRef(refprod, clusters++);
        }
        it->template rekey<HitType>(newRef);
      }  // end of the loop on a single detid

    }  // end of the loop on all clusters
  }

  //-------------------------------------------------------------
  // helper classes
  //-------------------------------------------------------------
  // FIXME (migrate to new RTTI and interface)
  // generic rekey (in practise for pixel only...)
  template <typename ClusterRefType>  // template for class
  template <typename RecHitType>      // template for member function
  void ClusterStorer::ClusterHitRecord<ClusterRefType>::rekey(const ClusterRefType &newRef) {
    TrackingRecHit &genericHit = (*hits_)[index_];
    RecHitType *hit = nullptr;
    if (genericHit.geographicalId().rawId() == detid_) {  // a hit on this det, so it's simple
      hit = dynamic_cast<RecHitType *>(&genericHit);      //static_cast<RecHitType *>(&genericHit);
    }
    assert(hit != nullptr);
    assert(hit->cluster() == ref_);  // otherwise something went wrong
    hit->setClusterRef(newRef);
  }

  // -------------------------------------------------------------
  // Specific rekey for class template ClusterRefType = SiStripRecHit2D::ClusterRef,
  // RecHitType is not used.
  template <>
  template <typename RecHitType>  // or template<> to specialise also here?
  void ClusterStorer::ClusterHitRecord<SiStripRecHit2D::ClusterRef>::rekey(const SiStripRecHit2D::ClusterRef &newRef) {
    TrackingRecHit &genericHit = (*hits_)[index_];
    const std::type_info &hit_type = typeid(genericHit);

    OmniClusterRef *cluRef = nullptr;
    if (typeid(SiStripRecHit1D) == hit_type) {
      cluRef = &static_cast<SiStripRecHit1D &>(genericHit).omniCluster();
    } else if (typeid(SiStripRecHit2D) == hit_type) {
      cluRef = &static_cast<SiStripRecHit2D &>(genericHit).omniCluster();
    } else if (typeid(SiStripMatchedRecHit2D) == hit_type) {
      SiStripMatchedRecHit2D &mhit = static_cast<SiStripMatchedRecHit2D &>(genericHit);
      cluRef = (SiStripDetId(detid_).stereo() ? &mhit.stereoClusterRef() : &mhit.monoClusterRef());
    } else if (typeid(ProjectedSiStripRecHit2D) == hit_type) {
      cluRef = &static_cast<ProjectedSiStripRecHit2D &>(genericHit).omniCluster();
    }

    assert(cluRef != nullptr);            // to catch missing RecHit types
    assert(cluRef->key() == ref_.key());  // otherwise something went wrong
    (*cluRef) = OmniClusterRef(newRef);
  }

  // -------------------------------------------------------------
  // Specific rekey for class template ClusterRefType = Phase2TrackerRecHit1D::ClusterRef
  template <>
  template <typename RecHitType>
  void ClusterStorer::ClusterHitRecord<Phase2TrackerRecHit1D::ClusterRef>::rekey(
      const Phase2TrackerRecHit1D::ClusterRef &newRef) {
    TrackingRecHit &genericHit = (*hits_)[index_];
    const std::type_info &hit_type = typeid(genericHit);

    OmniClusterRef *cluRef = nullptr;
    if (typeid(Phase2TrackerRecHit1D) == hit_type) {
      cluRef = &static_cast<Phase2TrackerRecHit1D &>(genericHit).omniCluster();
    } else if (typeid(VectorHit) == hit_type) {
      VectorHit &vHit = static_cast<VectorHit &>(genericHit);
      // FIXME: this essentially uses a hack
      // https://github.com/cms-sw/cmssw/blob/master/DataFormats/TrackerCommon/interface/TrackerTopology.h#L291
      cluRef = (SiStripDetId(detid_).stereo() ? &vHit.upperClusterRef() : &vHit.lowerClusterRef());
    }

    assert(cluRef != nullptr);            // to catch missing RecHit types
    assert(cluRef->key() == ref_.key());  // otherwise something went wrong
    (*cluRef) = OmniClusterRef(newRef);
  }

}  // namespace helper