CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
EgammaHLTRemoveDuplicatedSC.cc
Go to the documentation of this file.
7 
11 
12 #include <string>
13 
15 {
16 
17 
18  // the input producers
19 
20  sCInputProducer_ = ps.getParameter<edm::InputTag>("L1NonIsoUskimmedSC");
22  // set the producer parameters
23  outputCollection_ = ps.getParameter<std::string>("L1NonIsoSkimmedCollection");
24  produces<reco::SuperClusterCollection>(outputCollection_);
25 
26 
27 }
28 
30 {
31  ;
32 }
33 
34 void
36 {
37  using namespace edm;
38 
39 
40  // Get raw SuperClusters from the event
41  Handle<reco::SuperClusterCollection> UnskimmedSuperCluster;
42  evt.getByLabel(sCInputProducer_, UnskimmedSuperCluster);
43  Handle<reco::SuperClusterCollection> L1IsoSuperCluster;
44  evt.getByLabel(alreadyExistingSC_, L1IsoSuperCluster);
45  /*
46  edm::LogError("EgammaHLTRemoveDuplicatedSCError")
47  << "Error! can't get the rawSuperClusters "
48  << sCInputProducer_.label() ;
49  */
50 
51 
52  // Create a pointer to the existing SuperClusters
53  const reco::SuperClusterCollection *UnskimmedL1NonIsoSC = UnskimmedSuperCluster.product();
54  const reco::SuperClusterCollection *L1IsoSC = L1IsoSuperCluster.product();
55 
56  /*
57  for(reco::SuperClusterCollection::const_iterator it = L1IsoSC->begin(); it != L1IsoSC->end(); it++){
58  std::cout<<"L1 iso E, eta, phi: "<<it->energy()<<" "<<it->eta()<<" "<<it->phi()<<std::endl;
59  }
60  std::cout<<std::endl;
61  for(reco::SuperClusterCollection::const_iterator it = UnskimmedL1NonIsoSC->begin(); it != UnskimmedL1NonIsoSC->end(); it++){
62  std::cout<<"L1 Non iso (not skimmed) E, eta, phi: "<<it->energy()<<" "<<it->eta()<<" "<<it->phi()<<std::endl;
63  }
64  std::cout<<std::endl;
65  */
66 
67  // Define a collection of corrected SuperClusters to put back into the event
68  std::auto_ptr<reco::SuperClusterCollection> corrClusters(new reco::SuperClusterCollection);
69 
70  // Loop over raw clusters and make corrected ones
71  reco::SuperClusterCollection::const_iterator aClus;
72  reco::SuperClusterCollection::const_iterator cit;
73  for(aClus = UnskimmedL1NonIsoSC->begin(); aClus != UnskimmedL1NonIsoSC->end(); aClus++)
74  {
75  bool AlreadyThere = false;
76  //reco::SuperCluster newClus;
77  for(cit = L1IsoSC->begin(); cit != L1IsoSC->end(); cit++){
78  if( fabs(aClus->energy()- cit->energy()) < 0.5 && fabs(aClus->eta()- cit->eta()) < 0.0175 ){
79  float deltaphi=fabs( aClus->phi() - cit->phi() );
80  if(deltaphi>6.283185308) deltaphi -= 6.283185308;
81  if(deltaphi>3.141592654) deltaphi = 6.283185308-deltaphi;
82 
83  if( deltaphi < 0.035 ){AlreadyThere = true; break;}
84  }
85  }
86  // if(AlreadyThere){std::cout<<"AAAA: SC removed: "<<aClus->energy()<<" "<<aClus->eta()<<" "<<aClus->phi()<<std::endl;}
87  if(!AlreadyThere){corrClusters->push_back(*aClus);}
88  }
89 
90  // Put collection of corrected SuperClusters into the event
91  evt.put(corrClusters, outputCollection_);
92 
93 }
94 
95 
96 
T getParameter(std::string const &) const
virtual void produce(edm::Event &, const edm::EventSetup &)
EgammaHLTRemoveDuplicatedSC(const edm::ParameterSet &)
std::vector< SuperCluster > SuperClusterCollection
collection of SuperCluser objectr
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:85
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356