CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PFCandidateMixer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: PFCandidateMixer
4 // Class: PFCandidateMixer
5 //
13 //
14 // Original Author: Tomasz Maciej Frueboes
15 // Created: Wed Dec 9 16:14:56 CET 2009
16 // $Id: PFCandidateMixer.cc,v 1.5 2011/06/24 13:52:42 aburgmei Exp $
17 //
18 //
19 
20 
21 // system include files
22 #include <memory>
23 
24 // user include files
27 
30 
35 
40 
41 //
42 // class decleration
43 //
44 
46  public:
47  explicit PFCandidateMixer(const edm::ParameterSet&);
49 
50  private:
51  virtual void beginJob() ;
52  virtual void produce(edm::Event&, const edm::EventSetup&);
53  virtual void endJob() ;
54 
58  // ----------member data ---------------------------
59 };
60 
61 //
62 // constants, enums and typedefs
63 //
64 
65 
66 //
67 // static data member definitions
68 //
69 
70 //
71 // constructors and destructor
72 //
74  _col1(iConfig.getUntrackedParameter<edm::InputTag>("col1") ),
75  _col2(iConfig.getUntrackedParameter<edm::InputTag>("col2") ),
76  _trackCol(iConfig.getUntrackedParameter<edm::InputTag>("trackCol") )
77 {
78 
79  produces< std::vector< reco::PFCandidate > >();
80 
81 }
82 
84 {
85 
86  // do anything here that needs to be done at desctruction time
87  // (e.g. close files, deallocate resources etc.)
88 
89 }
90 
91 
92 //
93 // member functions
94 //
95 
96 // ------------ method called to produce the data ------------
97 void
99 {
100  using namespace edm;
101  using namespace reco;
102 
104  iEvent.getByLabel( _trackCol, trackCol);
105 
106 
107  std::vector< Handle<PFCandidateCollection> > colVec;
110  iEvent.getByLabel(_col1,pfIn1);
111  iEvent.getByLabel(_col2,pfIn2);
112 
113  colVec.push_back(pfIn1);
114  colVec.push_back(pfIn2);
115 
116  std::auto_ptr<std::vector< reco::PFCandidate > > pOut(new std::vector< reco::PFCandidate > );
117 
118  std::vector< Handle<PFCandidateCollection> >::iterator itCol= colVec.begin();
119  std::vector< Handle<PFCandidateCollection> >::iterator itColE= colVec.end();
120 
121  int iCol = 0;
122  for (;itCol!=itColE; ++itCol){
123  if (!itCol->isValid()) {
124  std::cout << "Whoops!" << std::endl;
125  }
126  PFCandidateConstIterator it = (*itCol)->begin();
127  PFCandidateConstIterator itE = (*itCol)->end();
128  for (;it!=itE;++it) {
129  PFCandidate cand(*it);
130  size_t i = 0;
131  bool found = false;
132  double minDR = 9999.;
133  int iMinDr = -1;
134  if (it->trackRef().isNonnull()) {
135  for ( i = 0 ; i < trackCol->size(); ++i){
136  if ( reco::deltaR( *(it->trackRef()), trackCol->at(i) )<0.001 ) {
137  found = true;
138  break;
139  }
140  double dr = reco::deltaR( *(it->trackRef()), trackCol->at(i) );
141  if ( dr < minDR) {
142  iMinDr = i;
143  minDR = dr;
144  }
145  }
146  }
147  if ( found ){ // ref was found, overwrite in PFCand
148  reco::TrackRef trref(trackCol,i);
149  cand.setTrackRef(trref);
150  //std::cout << " YY track ok"<<std::endl;
151 
152  } else { // keep orginall ref
153  if (it->trackRef().isNonnull()) {
154  std::cout << " XXXXXXXXXXX track not found "
155  << " col " << iCol
156  << " ch " << it->charge()
157  << " id " << it->pdgId()
158  << " pt " << it->pt()
159  << " track: eta " << it->trackRef()->eta()
160  << " pt: " << it->trackRef()->pt()
161  << " charge: " << it->trackRef()->charge()
162  << std::endl;
163  std::cout << " minDR=" << minDR << std::endl;
164  if ( iMinDr >= 0 ) {
165  std::cout
166  << " closest track pt=" << trackCol->at(iMinDr).pt()
167  << " ch=" << trackCol->at(iMinDr).charge()
168  << std::endl;
169  }
170  edm::Provenance prov=iEvent.getProvenance(it->trackRef().id());
172  std::cout << " trackref in PFCand came from: " << tag.encode() << std::endl;
173 
174 
175  }
176  }
177  pOut->push_back(cand);
178  }
179  ++iCol;
180  }
181 
182 
183  iEvent.put(pOut);
184 
185 
186 }
187 
188 // ------------ method called once each job just before starting event loop ------------
189 void
191 {
192 }
193 
194 // ------------ method called once each job just after ending the event loop ------------
195 void
197 }
198 
199 //define this as a plug-in
int i
Definition: DBlmapReader.cc:9
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
virtual void beginJob()
virtual void produce(edm::Event &, const edm::EventSetup &)
PFCandidateCollection::const_iterator PFCandidateConstIterator
iterator
std::string const & processName() const
Definition: Provenance.h:63
edm::InputTag _col2
double deltaR(double eta1, double phi1, double eta2, double phi2)
Definition: deltaR.h:19
int iEvent
Definition: GenABIO.cc:243
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:85
edm::InputTag _col1
virtual void endJob()
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
std::string const & moduleLabel() const
Definition: Provenance.h:62
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:33
edm::InputTag _trackCol
PFCandidateMixer(const edm::ParameterSet &)
tuple cout
Definition: gather_cfg.py:121
std::string const & productInstanceName() const
Definition: Provenance.h:64
void setTrackRef(const reco::TrackRef &ref)
set track reference
Definition: PFCandidate.cc:322
Provenance getProvenance(BranchID const &theID) const
Definition: Event.cc:60