CMS 3D CMS Logo

PixelQuadrupletMergerEDProducer.cc
Go to the documentation of this file.
10 
14 
15 // following are needed only to keep the same results
20 
22 public:
25 
26  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
27 
28  virtual void produce(edm::Event& iEvent, const edm::EventSetup& iSetup) override;
29 
30 private:
32 
34 
36 
37  // to keep old results
38  std::unique_ptr<SeedComparitor> comparitor_;
39  std::unique_ptr<SeedCreator> seedCreator_;
40 };
41 
43  tripletToken_(consumes<RegionsSeedingHitSets>(iConfig.getParameter<edm::InputTag>("triplets"))),
44  merger_(iConfig.getParameter<edm::ParameterSet>("layerList"), consumesCollector())
45 {
46  merger_.setTTRHBuilderLabel(iConfig.getParameter<std::string>("ttrhBuilderLabel"));
47  merger_.setMergeTriplets(iConfig.getParameter<bool>("mergeTriplets"));
48  merger_.setAddRemainingTriplets(iConfig.getParameter<bool>("addRemainingTriplets"));
49 
50  edm::ParameterSet comparitorPSet = iConfig.getParameter<edm::ParameterSet>("SeedComparitorPSet");
51  std::string comparitorName = comparitorPSet.getParameter<std::string>("ComponentName");
52  if(comparitorName != "none") {
53  auto iC = consumesCollector();
54  comparitor_.reset(SeedComparitorFactory::get()->create(comparitorName, comparitorPSet, iC));
55  }
56 
57  edm::ParameterSet creatorPSet = iConfig.getParameter<edm::ParameterSet>("SeedCreatorPSet");
58  std::string creatorName = creatorPSet.getParameter<std::string>("ComponentName");
59  if(creatorName != "none") // pixel tracking does not use seed creator
60  seedCreator_.reset(SeedCreatorFactory::get()->create( creatorName, creatorPSet));
61 
62  produces<RegionsSeedingHitSets>();
63  produces<TrajectorySeedCollection>(); // need to keep these in memory because TrajectorySeed owns its RecHits
64 }
65 
68 
69  desc.add<edm::InputTag>("triplets", edm::InputTag("hitTripletMergerEDProducer"));
70  desc.add<std::string>("ttrhBuilderLabel", "PixelTTRHBuilderWithoutAngle");
71  desc.add<bool>("mergeTriplets", true);
72  desc.add<bool>("addRemainingTriplets", false);
73  // This would be really on the responsibility of
74  // QuadrupletSeedMerger and SeedingLayerSetsBuilder. The former is
75  // almost obsolete by now (so I don't want to put effort there), and
76  // the latter is better dealt in the context of SeedingLayersEDProducer.
78  descLayers.setAllowAnything();
79  desc.add<edm::ParameterSetDescription>("layerList", descLayers);
80 
81  // to keep old results
82  edm::ParameterSetDescription descComparitor;
83  descComparitor.add<std::string>("ComponentName", "none");
84  descComparitor.setAllowAnything();
85  desc.add<edm::ParameterSetDescription>("SeedComparitorPSet", descComparitor);
86  edm::ParameterSetDescription descCreator;
87  descCreator.add<std::string>("ComponentName", "none");
88  descCreator.setAllowAnything();
89  desc.add<edm::ParameterSetDescription>("SeedCreatorPSet", descCreator);
90 
91  descriptions.add("pixelQuadrupletMergerEDProducer", desc);
92 }
93 
96  iEvent.getByToken(tripletToken_, htriplets);
97  const auto& regionTriplets = *htriplets;
98 
99  auto seedingHitSets = std::make_unique<RegionsSeedingHitSets>();
100  if(regionTriplets.empty()) {
101  iEvent.put(std::move(seedingHitSets));
102  return;
103  }
104  seedingHitSets->reserve(regionTriplets.regionSize(), localRA_.upper());
105 
106  // to keep old results
107  auto tmpSeedCollection = std::make_unique<TrajectorySeedCollection>();
108 
109  OrderedHitSeeds quadruplets;
110  quadruplets.reserve(localRA_.upper());
111 
112  OrderedHitSeeds tripletsPerRegion;
113  tripletsPerRegion.reserve(localRA_.upper());
114 
115  LogDebug("PixelQuadrupletMergerEDProducer") << "Creating quadruplets for " << regionTriplets.regionSize() << " regions from " << regionTriplets.size() << " triplets";
116  merger_.update(iSetup);
117 
118  // to keep old results
119  if(comparitor_) comparitor_->init(iEvent, iSetup);
120 
121  for(const auto& regionSeedingHitSets: regionTriplets) {
122  const TrackingRegion& region = regionSeedingHitSets.region();
123  auto seedingHitSetsFiller = seedingHitSets->beginRegion(&region);
124 
125 
126  // Keeping same resuls has been made really difficult...
127  // Especially when supporting both pixel tracking and seeding
128  // Following is from SeedGeneratorFromRegionHits
129  if(seedCreator_) {
130  seedCreator_->init(region, iSetup, comparitor_.get());
131  for(const auto& hits: regionSeedingHitSets) {
132  if(!comparitor_ || comparitor_->compatible(hits)) {
133  seedCreator_->makeSeed(*tmpSeedCollection, hits);
134  }
135  }
136 
137  // then convert seeds back to hits
138  // awful, but hopefully only temporary to preserve old results
139  for(const auto& seed: *tmpSeedCollection) {
140  auto hitRange = seed.recHits();
141  assert(std::distance(hitRange.first, hitRange.second) == 3);
142  tripletsPerRegion.emplace_back(static_cast<SeedingHitSet::ConstRecHitPointer>(&*(hitRange.first)),
143  static_cast<SeedingHitSet::ConstRecHitPointer>(&*(hitRange.first+1)),
144  static_cast<SeedingHitSet::ConstRecHitPointer>(&*(hitRange.first+2)));
145  }
146  }
147  else {
148  for(const auto& hits: regionSeedingHitSets) {
149  tripletsPerRegion.emplace_back(hits[0], hits[1], hits[2]);
150  }
151  }
152 
153  LogTrace("PixelQuadrupletEDProducer") << " starting region, number of triplets " << tripletsPerRegion.size();
154 
155  const auto& quadruplets = merger_.mergeTriplets(tripletsPerRegion, iSetup);
156 
157  LogTrace("PixelQuadrupletEDProducer") << " created " << quadruplets.size() << " quadruplets";
158 
159  for(size_t i=0; i!= quadruplets.size(); ++i) {
160  const auto& quad = quadruplets[i];
161  seedingHitSetsFiller.emplace_back(quad[0], quad[1], quad[2], quad[3]);
162  }
163 
164  tripletsPerRegion.clear();
165  }
166  localRA_.update(seedingHitSets->size());
167 
168  seedingHitSets->shrink_to_fit();
169  tmpSeedCollection->shrink_to_fit();
170  iEvent.put(std::move(seedingHitSets));
171  iEvent.put(std::move(tmpSeedCollection));
172 }
173 
#define LogDebug(id)
T getParameter(std::string const &) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:122
PixelQuadrupletMergerEDProducer(const edm::ParameterSet &iConfig)
void setTTRHBuilderLabel(std::string)
def create(alignables, pedeDump, additionalData, outputFile, config)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
void setAllowAnything()
allow any parameter label/value pairs
edm::EDGetTokenT< RegionsSeedingHitSets > tripletToken_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
virtual unsigned int size() const
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
int upper() const
int iEvent
Definition: GenABIO.cc:230
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.
std::unique_ptr< SeedCreator > seedCreator_
ParameterDescriptionBase * add(U const &iLabel, T const &value)
#define LogTrace(id)
void update(const edm::EventSetup &)
const OrderedSeedingHits & mergeTriplets(const OrderedSeedingHits &, const edm::EventSetup &)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
HLT enums.
virtual void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) override
std::unique_ptr< SeedComparitor > comparitor_
def move(src, dest)
Definition: eostools.py:510
T get(const Candidate &c)
Definition: component.h:55
void update(unsigned int q)