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