test
CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SeedCombiner.cc
Go to the documentation of this file.
1 #include "SeedCombiner.h"
2 
7 
9 
13 
16 
20 
22 
25 
26 using namespace edm;
27 //using namespace reco;
28 
30 {
31  inputCollections_ = edm::vector_transform(
32  cfg.getParameter<std::vector<edm::InputTag> >("seedCollections"),
33  [this](edm::InputTag const & tag) { return consumes<TrajectorySeedCollection>(tag); } );
34  produces<TrajectorySeedCollection>();
35  reKeing_=false;
36  if (cfg.exists("clusterRemovalInfos")){
37  clusterRemovalInfos_=cfg.getParameter<std::vector<edm::InputTag> >("clusterRemovalInfos");
38  clusterRemovalTokens_.resize(clusterRemovalInfos_.size());
39  for (unsigned int i=0;i<clusterRemovalInfos_.size();++i)
40  if (!(clusterRemovalInfos_[i]==edm::InputTag("")))
41  clusterRemovalTokens_[i] = consumes<reco::ClusterRemovalInfo>(clusterRemovalInfos_[i]);
42  if (clusterRemovalInfos_.size()!=0 && clusterRemovalInfos_.size()==inputCollections_.size()) reKeing_=true;
43  }
44 }
45 
46 
48 {
49 }
50 
51 
53 {
54  // Read inputs, and count total seeds
55  size_t ninputs = inputCollections_.size();
56  size_t nseeds = 0;
57  std::vector<Handle<TrajectorySeedCollection > > seedCollections(ninputs);
58  for (size_t i = 0; i < ninputs; ++i) {
59  ev.getByToken(inputCollections_[i], seedCollections[i]);
60  nseeds += seedCollections[i]->size();
61  }
62 
63  // Prepare output collections, with the correct capacity
64  std::auto_ptr<TrajectorySeedCollection> result(new TrajectorySeedCollection());
65  result->reserve( nseeds );
66 
67  // Write into output collection
68  unsigned int iSC=0,iSC_max=seedCollections.size();
69  for (;iSC!=iSC_max;++iSC){
70  Handle<TrajectorySeedCollection> & collection=seedCollections[iSC];
71  if (reKeing_ && !(clusterRemovalInfos_[iSC]==edm::InputTag(""))){
73 
74  for (TrajectorySeedCollection::const_iterator iS=collection->begin();
75  iS!=collection->end();++iS){
76  TrajectorySeed::recHitContainer newRecHitContainer;
77  newRecHitContainer.reserve(iS->nHits());
78  TrajectorySeed::const_iterator iH=iS->recHits().first;
79  TrajectorySeed::const_iterator iH_end=iS->recHits().second;
80  //loop seed rechits, copy over and rekey.
81  for (;iH!=iH_end;++iH){
82  newRecHitContainer.push_back(*iH);
83  refSetter.reKey(&newRecHitContainer.back());
84  }
85  result->push_back(TrajectorySeed(iS->startingState(),
86  std::move(newRecHitContainer),
87  iS->direction()));
88  }
89  }else{
90  //just insert the new seeds as they are
91  result->insert(result->end(), collection->begin(), collection->end());
92  }
93  }
94 
95  // Save result into the event
96  ev.put(result);
97 }
SeedCombiner(const edm::ParameterSet &cfg)
Definition: SeedCombiner.cc:29
T getParameter(std::string const &) const
std::vector< edm::InputTag > clusterRemovalInfos_
Definition: SeedCombiner.h:24
int i
Definition: DBlmapReader.cc:9
reference back()
Definition: OwnVector.h:321
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
void reKey(TrackingRecHit *hit) const
bool ev
auto vector_transform(std::vector< InputType > const &input, Function predicate) -> std::vector< typename std::remove_cv< typename std::remove_reference< decltype(predicate(input.front()))>::type >::type >
Definition: transform.h:11
void push_back(D *&d)
Definition: OwnVector.h:274
std::vector< TrajectorySeed > TrajectorySeedCollection
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
recHitContainer::const_iterator const_iterator
tuple result
Definition: query.py:137
std::vector< edm::EDGetTokenT< reco::ClusterRemovalInfo > > clusterRemovalTokens_
Definition: SeedCombiner.h:25
std::vector< edm::EDGetTokenT< TrajectorySeedCollection > > inputCollections_
Definition: SeedCombiner.h:22
virtual void produce(edm::Event &ev, const edm::EventSetup &es) override
Definition: SeedCombiner.cc:52
void reserve(size_t)
Definition: OwnVector.h:268