CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SeedGeneratorFromRegionHitsEDProducer.cc
Go to the documentation of this file.
2 
5 
7 
11 
14 
17 
20 
23 
24 
26  const edm::ParameterSet& cfg)
27  : theConfig(cfg), theGenerator(0), theRegionProducer(0),
28  theClusterCheck(cfg.getParameter<edm::ParameterSet>("ClusterCheckPSet"),consumesCollector()),
29  theMerger_(0)
30 {
31  theSilentOnClusterCheck = cfg.getParameter<edm::ParameterSet>("ClusterCheckPSet").getUntrackedParameter<bool>("silentClusterCheck",false);
32 
33  moduleName = cfg.getParameter<std::string>("@module_label");
34 
35  // seed merger & its settings
36  if ( cfg.exists("SeedMergerPSet")) {
37  edm::ParameterSet mergerPSet = theConfig.getParameter<edm::ParameterSet>( "SeedMergerPSet" );
39  theMerger_->setTTRHBuilderLabel( mergerPSet.getParameter<std::string>( "ttrhBuilderLabel" ) );
40  theMerger_->setMergeTriplets( mergerPSet.getParameter<bool>( "mergeTriplets" ) );
41  theMerger_->setAddRemainingTriplets( mergerPSet.getParameter<bool>( "addRemainingTriplets" ) );
42  theMerger_->setLayerListName( mergerPSet.getParameter<std::string>( "layerListName" ) );
43  }
44 
45  edm::ParameterSet regfactoryPSet =
46  theConfig.getParameter<edm::ParameterSet>("RegionFactoryPSet");
47  std::string regfactoryName = regfactoryPSet.getParameter<std::string>("ComponentName");
48  theRegionProducer = TrackingRegionProducerFactory::get()->create(regfactoryName,regfactoryPSet, consumesCollector());
49 
50  produces<TrajectorySeedCollection>();
51 }
52 
54 {
55  delete theRegionProducer; theRegionProducer = nullptr;
56 }
57 
59  delete theGenerator;
60 }
61 
63 {
64  edm::ParameterSet hitsfactoryPSet =
65  theConfig.getParameter<edm::ParameterSet>("OrderedHitsFactoryPSet");
66  std::string hitsfactoryName = hitsfactoryPSet.getParameter<std::string>("ComponentName");
67  OrderedHitsGenerator* hitsGenerator =
68  OrderedHitsGeneratorFactory::get()->create( hitsfactoryName, hitsfactoryPSet);
69 
70  edm::ParameterSet comparitorPSet =
71  theConfig.getParameter<edm::ParameterSet>("SeedComparitorPSet");
72  std::string comparitorName = comparitorPSet.getParameter<std::string>("ComponentName");
73  SeedComparitor * aComparitor = (comparitorName == "none") ?
74  0 : SeedComparitorFactory::get()->create( comparitorName, comparitorPSet);
75 
76  edm::ParameterSet creatorPSet =
77  theConfig.getParameter<edm::ParameterSet>("SeedCreatorPSet");
78  std::string creatorName = creatorPSet.getParameter<std::string>("ComponentName");
79  SeedCreator * aCreator = SeedCreatorFactory::get()->create( creatorName, creatorPSet);
80 
81  theGenerator = new SeedGeneratorFromRegionHits(hitsGenerator, aComparitor, aCreator);
82 }
83 
85 {
86  std::auto_ptr<TrajectorySeedCollection> triplets(new TrajectorySeedCollection());
87  std::auto_ptr<TrajectorySeedCollection> quadruplets( new TrajectorySeedCollection() );
88 
89  //protection for big ass events...
90  size_t clustsOrZero = theClusterCheck.tooManyClusters(ev);
91  if (clustsOrZero){
93  edm::LogError("TooManyClusters") << "Found too many clusters (" << clustsOrZero << "), bailing out.\n";
94  ev.put(triplets);
95  return ;
96  }
97 
98  typedef std::vector<TrackingRegion* > Regions;
99  typedef Regions::const_iterator IR;
100  Regions regions = theRegionProducer->regions(ev,es);
101  if (theMerger_!=0)
102  theMerger_->update(es);
103 
104  for (IR ir=regions.begin(), irEnd=regions.end(); ir < irEnd; ++ir) {
105  const TrackingRegion & region = **ir;
106 
107  // make job
108  theGenerator->run(*triplets, region, ev,es);
109  // std::cout << "created seeds for " << moduleName << " " << triplets->size() << std::endl;
110 
111 
112  // make quadruplets
113  // (TODO: can partly be propagated to the merger)
114  if ( theMerger_ !=0 ) {
115  TrajectorySeedCollection const& tempQuads = theMerger_->mergeTriplets( *triplets, region, es, theConfig ); //@@
116  for( TrajectorySeedCollection::const_iterator qIt = tempQuads.begin(); qIt < tempQuads.end(); ++qIt ) {
117  quadruplets->push_back( *qIt );
118  }
119  }
120  }
121 
122  // clear memory
123  for (IR ir=regions.begin(), irEnd=regions.end(); ir < irEnd; ++ir) delete (*ir);
124 
125  // put to event
126  if ( theMerger_!=0)
127  ev.put(quadruplets);
128  else
129  ev.put(triplets);
130 }
T getParameter(std::string const &) const
bool exists(std::string const &parameterName) const
checks if a parameter exists
virtual void produce(edm::Event &ev, const edm::EventSetup &es) override
std::vector< TrajectorySeed > TrajectorySeedCollection
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.
size_t tooManyClusters(const edm::Event &e) const
void update(const edm::EventSetup &)
const OrderedSeedingHits & mergeTriplets(const OrderedSeedingHits &, const edm::EventSetup &)
void run(TrajectorySeedCollection &seedCollection, const TrackingRegion &region, const edm::Event &ev, const edm::EventSetup &es)
return(e1-e2)*(e1-e2)+dp *dp
virtual void beginRun(edm::Run const &run, const edm::EventSetup &es) override
virtual void endRun(edm::Run const &run, const edm::EventSetup &es) override
T get(const Candidate &c)
Definition: component.h:55
Definition: Run.h:41
virtual std::vector< TrackingRegion * > regions(const edm::Event &ev, const edm::EventSetup &es) const =0