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 
6 
8 
12 
15 
18 
21 
24 
25 
27  const edm::ParameterSet& cfg)
28  : theConfig(cfg), theRegionProducer(0),
29  theClusterCheck(cfg.getParameter<edm::ParameterSet>("ClusterCheckPSet"),consumesCollector()),
30  theMerger_(0)
31 {
32  theSilentOnClusterCheck = cfg.getParameter<edm::ParameterSet>("ClusterCheckPSet").getUntrackedParameter<bool>("silentClusterCheck",false);
33 
34  moduleName = cfg.getParameter<std::string>("@module_label");
35 
36  // seed merger & its settings
38  if ( cfg.exists("SeedMergerPSet")) {
39  edm::ParameterSet mergerPSet = theConfig.getParameter<edm::ParameterSet>( "SeedMergerPSet" );
40  theMerger_=new QuadrupletSeedMerger(mergerPSet.getParameter<edm::ParameterSet>( "layerList" ), iC);
41  theMerger_->setTTRHBuilderLabel( mergerPSet.getParameter<std::string>( "ttrhBuilderLabel" ) );
42  theMerger_->setMergeTriplets( mergerPSet.getParameter<bool>( "mergeTriplets" ) );
43  theMerger_->setAddRemainingTriplets( mergerPSet.getParameter<bool>( "addRemainingTriplets" ) );
44  }
45 
46  edm::ParameterSet regfactoryPSet =
47  theConfig.getParameter<edm::ParameterSet>("RegionFactoryPSet");
48  std::string regfactoryName = regfactoryPSet.getParameter<std::string>("ComponentName");
49  theRegionProducer = TrackingRegionProducerFactory::get()->create(regfactoryName,regfactoryPSet, consumesCollector());
50 
51  edm::ParameterSet hitsfactoryPSet =
52  theConfig.getParameter<edm::ParameterSet>("OrderedHitsFactoryPSet");
53  std::string hitsfactoryName = hitsfactoryPSet.getParameter<std::string>("ComponentName");
54  OrderedHitsGenerator* hitsGenerator =
55  OrderedHitsGeneratorFactory::get()->create( hitsfactoryName, hitsfactoryPSet, iC);
56 
57  edm::ParameterSet comparitorPSet =
58  theConfig.getParameter<edm::ParameterSet>("SeedComparitorPSet");
59  std::string comparitorName = comparitorPSet.getParameter<std::string>("ComponentName");
60  SeedComparitor * aComparitor = (comparitorName == "none") ?
61  0 : SeedComparitorFactory::get()->create( comparitorName, comparitorPSet);
62 
63  edm::ParameterSet creatorPSet =
64  theConfig.getParameter<edm::ParameterSet>("SeedCreatorPSet");
65  std::string creatorName = creatorPSet.getParameter<std::string>("ComponentName");
66  SeedCreator * aCreator = SeedCreatorFactory::get()->create( creatorName, creatorPSet);
67 
68  theGenerator.reset(new SeedGeneratorFromRegionHits(hitsGenerator, aComparitor, aCreator));
69 
70  produces<TrajectorySeedCollection>();
71 }
72 
74 {
75  delete theRegionProducer; theRegionProducer = nullptr;
76 }
77 
79 }
80 
82 {
83 }
84 
86 {
87  std::auto_ptr<TrajectorySeedCollection> triplets(new TrajectorySeedCollection());
88  std::auto_ptr<TrajectorySeedCollection> quadruplets( new TrajectorySeedCollection() );
89 
90  //protection for big ass events...
91  size_t clustsOrZero = theClusterCheck.tooManyClusters(ev);
92  if (clustsOrZero){
94  edm::LogError("TooManyClusters") << "Found too many clusters (" << clustsOrZero << "), bailing out.\n";
95  ev.put(triplets);
96  return ;
97  }
98 
99  typedef std::vector<TrackingRegion* > Regions;
100  typedef Regions::const_iterator IR;
101  Regions regions = theRegionProducer->regions(ev,es);
102  if (theMerger_!=0)
103  theMerger_->update(es);
104 
105  for (IR ir=regions.begin(), irEnd=regions.end(); ir < irEnd; ++ir) {
106  const TrackingRegion & region = **ir;
107 
108  // make job
109  theGenerator->run(*triplets, region, ev,es);
110  // std::cout << "created seeds for " << moduleName << " " << triplets->size() << std::endl;
111 
112 
113  // make quadruplets
114  // (TODO: can partly be propagated to the merger)
115  if ( theMerger_ !=0 ) {
116  TrajectorySeedCollection const& tempQuads = theMerger_->mergeTriplets( *triplets, region, es, theConfig ); //@@
117  for( TrajectorySeedCollection::const_iterator qIt = tempQuads.begin(); qIt < tempQuads.end(); ++qIt ) {
118  quadruplets->push_back( *qIt );
119  }
120  }
121  }
122 
123  // clear memory
124  for (IR ir=regions.begin(), irEnd=regions.end(); ir < irEnd; ++ir) delete (*ir);
125 
126  // put to event
127  if ( theMerger_!=0)
128  ev.put(quadruplets);
129  else
130  ev.put(triplets);
131 }
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 &)
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
std::unique_ptr< SeedGeneratorFromRegionHits > theGenerator
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