CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
List of all members | Public Member Functions
hitTripletEDProducerT::Impl< T_HitTripletGenerator, T_SeedingHitSets, T_IntermediateHitTriplets > Class Template Reference

#include <HitTripletEDProducerT.h>

Inheritance diagram for hitTripletEDProducerT::Impl< T_HitTripletGenerator, T_SeedingHitSets, T_IntermediateHitTriplets >:
hitTripletEDProducerT::ImplGeneratorBase< T_HitTripletGenerator > hitTripletEDProducerT::ImplBase

Public Member Functions

 Impl (const edm::ParameterSet &iConfig, edm::ConsumesCollector &iC)
 
void produce (const IntermediateHitDoublets &regionDoublets, edm::Event &iEvent, const edm::EventSetup &iSetup) override
 
void produces (edm::ProducesCollector producesCollector) const override
 
 ~Impl () override=default
 
- Public Member Functions inherited from hitTripletEDProducerT::ImplGeneratorBase< T_HitTripletGenerator >
 ImplGeneratorBase (const edm::ParameterSet &iConfig, edm::ConsumesCollector &iC)
 
 ~ImplGeneratorBase () override=default
 
- Public Member Functions inherited from hitTripletEDProducerT::ImplBase
 ImplBase ()=default
 
virtual ~ImplBase ()=default
 

Additional Inherited Members

- Protected Attributes inherited from hitTripletEDProducerT::ImplGeneratorBase< T_HitTripletGenerator >
T_HitTripletGenerator generator_
 
- Protected Attributes inherited from hitTripletEDProducerT::ImplBase
edm::RunningAverage localRA_
 

Detailed Description

template<typename T_HitTripletGenerator, typename T_SeedingHitSets, typename T_IntermediateHitTriplets>
class hitTripletEDProducerT::Impl< T_HitTripletGenerator, T_SeedingHitSets, T_IntermediateHitTriplets >

Definition at line 71 of file HitTripletEDProducerT.h.

Constructor & Destructor Documentation

template<typename T_HitTripletGenerator , typename T_SeedingHitSets , typename T_IntermediateHitTriplets >
hitTripletEDProducerT::Impl< T_HitTripletGenerator, T_SeedingHitSets, T_IntermediateHitTriplets >::Impl ( const edm::ParameterSet iConfig,
edm::ConsumesCollector iC 
)
inline
template<typename T_HitTripletGenerator , typename T_SeedingHitSets , typename T_IntermediateHitTriplets >
hitTripletEDProducerT::Impl< T_HitTripletGenerator, T_SeedingHitSets, T_IntermediateHitTriplets >::~Impl ( )
overridedefault

Member Function Documentation

template<typename T_HitTripletGenerator , typename T_SeedingHitSets , typename T_IntermediateHitTriplets >
void hitTripletEDProducerT::Impl< T_HitTripletGenerator, T_SeedingHitSets, T_IntermediateHitTriplets >::produce ( const IntermediateHitDoublets regionDoublets,
edm::Event iEvent,
const edm::EventSetup iSetup 
)
inlineoverridevirtual

Implements hitTripletEDProducerT::ImplBase.

Definition at line 82 of file HitTripletEDProducerT.h.

References a, b, IntermediateHitDoublets::empty(), Exception, funct::exp(), LayerHitMapCache::extend(), newFWLiteAna::found, hitTripletEDProducerT::ImplGeneratorBase< T_HitTripletGenerator >::generator_, cmsLHEtoEOSManager::l, IntermediateHitDoublets::layerPairsSize(), LayerTriplets::layers(), hitTripletEDProducerT::ImplBase::localRA_, LogDebug, LogTrace, HLT_FULL_cff::region, IntermediateHitDoublets::regionSize(), IntermediateHitDoublets::seedingLayerHits(), OrderedHitTriplets::size(), edm::RunningAverage::update(), and edm::RunningAverage::upper().

84  {
85  const SeedingLayerSetsHits& seedingLayerHits = regionDoublets.seedingLayerHits();
86 
87  auto seedingHitSetsProducer = T_SeedingHitSets();
88  auto intermediateHitTripletsProducer = T_IntermediateHitTriplets(&seedingLayerHits);
89 
90  if (regionDoublets.empty()) {
91  seedingHitSetsProducer.putEmpty(iEvent);
92  intermediateHitTripletsProducer.putEmpty(iEvent);
93  return;
94  }
95 
96  seedingHitSetsProducer.reserve(regionDoublets.regionSize(), this->localRA_.upper());
97  intermediateHitTripletsProducer.reserve(regionDoublets.regionSize(), this->localRA_.upper());
98 
99  // match-making of pair and triplet layers
100  std::vector<LayerTriplets::LayerSetAndLayers> trilayers = LayerTriplets::layers(seedingLayerHits);
101 
102  OrderedHitTriplets triplets;
103  triplets.reserve(this->localRA_.upper());
104  size_t triplets_total = 0;
105 
106  LogDebug("HitTripletEDProducer") << "Creating triplets for " << regionDoublets.regionSize() << " regions, and "
107  << trilayers.size() << " pair+3rd layers from "
108  << regionDoublets.layerPairsSize() << " layer pairs";
109 
110  for (const auto& regionLayerPairs : regionDoublets) {
111  const TrackingRegion& region = regionLayerPairs.region();
112 
113  auto hitCachePtr_filler_shs = seedingHitSetsProducer.beginRegion(&region, nullptr);
114  auto hitCachePtr_filler_iht =
115  intermediateHitTripletsProducer.beginRegion(&region, std::get<0>(hitCachePtr_filler_shs));
116  auto hitCachePtr = std::get<0>(hitCachePtr_filler_iht);
117 
118  LayerHitMapCache& hitCache = *hitCachePtr;
119  hitCache.extend(regionLayerPairs.layerHitMapCache());
120 
121  LogTrace("HitTripletEDProducer") << " starting region";
122 
123  for (const auto& layerPair : regionLayerPairs) {
124  LogTrace("HitTripletEDProducer")
125  << " starting layer pair " << layerPair.innerLayerIndex() << "," << layerPair.outerLayerIndex();
126 
127  auto found = std::find_if(trilayers.begin(), trilayers.end(), [&](const LayerTriplets::LayerSetAndLayers& a) {
128  return a.first[0].index() == layerPair.innerLayerIndex() &&
129  a.first[1].index() == layerPair.outerLayerIndex();
130  });
131  if (found == trilayers.end()) {
132  auto exp = cms::Exception("LogicError") << "Did not find the layer pair from vector<pair+third layers>. "
133  "This is a sign of some internal inconsistency\n";
134  exp << "I was looking for layer pair " << layerPair.innerLayerIndex() << "," << layerPair.outerLayerIndex()
135  << ". Triplets have the following pairs:\n";
136  for (const auto& a : trilayers) {
137  exp << " " << a.first[0].index() << "," << a.first[1].index() << ": 3rd layers";
138  for (const auto& b : a.second) {
139  exp << " " << b.index();
140  }
141  exp << "\n";
142  }
143  throw exp;
144  }
145  const auto& thirdLayers = found->second;
146 
147  this->generator_.hitTriplets(region,
148  triplets,
149  iEvent,
150  iSetup,
151  layerPair.doublets(),
152  thirdLayers,
153  intermediateHitTripletsProducer.tripletLastLayerIndexVector(),
154  hitCache);
155 
156 #ifdef EDM_ML_DEBUG
157  LogTrace("HitTripletEDProducer")
158  << " created " << triplets.size() << " triplets for layer pair " << layerPair.innerLayerIndex() << ","
159  << layerPair.outerLayerIndex() << " and 3rd layers";
160  for (const auto& l : thirdLayers) {
161  LogTrace("HitTripletEDProducer") << " " << l.index();
162  }
163 #endif
164 
165  triplets_total += triplets.size();
166  seedingHitSetsProducer.fill(std::get<1>(hitCachePtr_filler_shs), triplets);
167  intermediateHitTripletsProducer.fill(
168  std::get<1>(hitCachePtr_filler_iht), layerPair.layerPair(), thirdLayers, triplets);
169 
170  triplets.clear();
171  }
172  }
173  this->localRA_.update(triplets_total);
174 
175  seedingHitSetsProducer.put(iEvent);
176  intermediateHitTripletsProducer.put(iEvent);
177  }
std::vector< LayerSetAndLayers > layers(const SeedingLayerSetsHits &sets)
Definition: LayerTriplets.cc:4
Exp< T >::type exp(const T &t)
Definition: Exp.h:22
#define LogTrace(id)
int upper() const
std::pair< LayerSet, std::vector< Layer > > LayerSetAndLayers
Definition: LayerTriplets.h:16
double b
Definition: hdecay.h:118
void extend(const LayerHitMapCache &other)
double a
Definition: hdecay.h:119
unsigned int size() const override
void update(unsigned int q)
const SeedingLayerSetsHits & seedingLayerHits() const
#define LogDebug(id)
template<typename T_HitTripletGenerator , typename T_SeedingHitSets , typename T_IntermediateHitTriplets >
void hitTripletEDProducerT::Impl< T_HitTripletGenerator, T_SeedingHitSets, T_IntermediateHitTriplets >::produces ( edm::ProducesCollector  producesCollector) const
inlineoverridevirtual

Implements hitTripletEDProducerT::ImplBase.

Definition at line 77 of file HitTripletEDProducerT.h.

77  {
78  T_SeedingHitSets::produces(producesCollector);
79  T_IntermediateHitTriplets::produces(producesCollector);
80  };