CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ZmumuPFEmbedder.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: ZmumuPFEmbedder
4 // Class: ZmumuPFEmbedder
5 //
13 //
14 // Original Author: Tomasz Maciej Frueboes
15 // Created: Wed Dec 9 16:14:56 CET 2009
16 // $Id: ZmumuPFEmbedder.cc,v 1.12 2012/01/27 13:17:12 aburgmei Exp $
17 //
18 //
19 
20 
21 // system include files
22 #include <memory>
23 
24 // user include files
27 
30 
36 
39 
43 
46 //
47 // class decleration
48 //
49 
51  public:
52  explicit ZmumuPFEmbedder(const edm::ParameterSet&);
54 
55  private:
56  virtual void beginJob() ;
57  virtual void produce(edm::Event&, const edm::EventSetup&);
58  void producePFCandColl(edm::Event&, const std::vector< reco::Particle::LorentzVector > * toBeAdded );
59  void produceTrackColl(edm::Event&, const std::vector< reco::Particle::LorentzVector > * toBeAdded );
60  virtual void endJob() ;
61 
65 
66  // ----------member data ---------------------------
67 };
68 
69 //
70 // constants, enums and typedefs
71 //
72 
73 
74 //
75 // static data member definitions
76 //
77 
78 //
79 // constructors and destructor
80 //
82  : _tracks(iConfig.getParameter<edm::InputTag>("tracks")),
83  _selectedMuons(iConfig.getParameter<edm::InputTag>("selectedMuons")),
84  _useCombinedCandidate(iConfig.getUntrackedParameter<bool>("useCombinedCandidate", false))
85 {
86 
87  //register your products
88  // produces< std::vector< reco::Muon > >("zMusExtracted"); //
89  produces<reco::TrackCollection>("tracks");
90  produces< std::vector< reco::PFCandidate > >("pfCands");
91 
92 
93 }
94 
96 {
97 
98  // do anything here that needs to be done at desctruction time
99  // (e.g. close files, deallocate resources etc.)
100 
101 }
102 
103 
104 //
105 // member functions
106 //
107 
108 // ------------ method called to produce the data ------------
109 void
111 {
112  std::vector< reco::Particle::LorentzVector > toBeAdded;
113 
115  {
117  if (iEvent.getByLabel(_selectedMuons, combCandidatesHandle) && combCandidatesHandle->size()>0)
118  for (size_t idx = 0; idx < combCandidatesHandle->at(0).numberOfDaughters(); ++idx) // use only the first combined candidate
119  toBeAdded.push_back(combCandidatesHandle->at(0).daughter(idx)->p4());
120  }
121  else
122  {
123  edm::Handle< edm::View< reco::Muon > > selectedZMuonsHandle;
124  if (iEvent.getByLabel(_selectedMuons, selectedZMuonsHandle))
125  for (size_t idx = 0; idx < selectedZMuonsHandle->size(); ++idx)
126  toBeAdded.push_back(selectedZMuonsHandle->at(idx).p4());
127  }
128 
129  if (toBeAdded.size() == 0)
130  return;
131 
132  producePFCandColl(iEvent, &toBeAdded);
133  produceTrackColl(iEvent, &toBeAdded);
134 }
135 
136 // produces clean PFCandidate collection w/o muon candidates.
137 void ZmumuPFEmbedder::producePFCandColl(edm::Event & iEvent, const std::vector< reco::Particle::LorentzVector > * toBeAdded )
138 {
140  iEvent.getByLabel("particleFlow",pfIn);
141 
142  //std::vector< reco::Muon > toBeAdded;
143  // TODO - check col size
144  //reco::Muon l1 = *(toBeAdded->begin());
145  //reco::Muon l2 = *(toBeAdded->rbegin());
146 
147  std::auto_ptr<std::vector< reco::PFCandidate > > newCol(new std::vector< reco::PFCandidate > );
148 
149  //get selected muons
150  // iterate over pfcandidates, make copy if its not a selected muon
151  reco::PFCandidateConstIterator it = pfIn->begin();
152  reco::PFCandidateConstIterator itE = pfIn->end();
153 
154  for (;it!=itE;++it) {
155  int pdg = std::abs( it->pdgId() );
156  double minDR = 10;
157  /*
158  double dr1 = reco::deltaR( *it, l1);
159  double dr2 = reco::deltaR( *it, l2);
160  */
161  std::vector< reco::Particle::LorentzVector >::const_iterator itSelectedMu = toBeAdded->begin();
162  std::vector< reco::Particle::LorentzVector >::const_iterator itSelectedMuE = toBeAdded->end();
163  for (; itSelectedMu != itSelectedMuE; ++itSelectedMu ){
164  double dr = reco::deltaR( *it, *itSelectedMu);
165  if (dr < minDR) minDR = dr;
166  }
167 
168  //if ( pdg == 13 && (dr1 < 0.001 || dr2 < 0.002 ) ) { // it is a selected muon, do not copy
169  if ( pdg == 13 && (minDR < 0.001 ) ) { // it is a selected muon, do not copy
170 
171 
172  } else {
173  newCol->push_back(*it);
174  }
175  }
176 
177  iEvent.put(newCol, "pfCands");
178 }
179 
180 // produces clean track collection w/o muon tracks.
181 void ZmumuPFEmbedder::produceTrackColl(edm::Event & iEvent, const std::vector< reco::Particle::LorentzVector > * toBeAdded )
182 {
184  iEvent.getByLabel( _tracks, tks);
185 
186  std::auto_ptr< reco::TrackCollection > newCol(new reco::TrackCollection );
187 
188  double epsilon = 0.00001;
189  unsigned int nMatched = 0;
190 
191  for ( reco::TrackCollection::const_iterator it = tks->begin() ; it != tks->end(); ++it)
192  {
193  bool ok = true;
194  for ( std::vector< reco::Particle::LorentzVector >::const_iterator itTBA = toBeAdded->begin();
195  itTBA != toBeAdded->end();
196  ++itTBA)
197  {
198  double dr = reco::deltaR( *it, *itTBA);
199  //std::cout << "TTTT " << dr << std::endl;
200  if (dr < epsilon) {
201  ++ nMatched;
202  ok = false;
203  }
204  }
205  if (ok) newCol->push_back(*it);
206  }
207  if (nMatched!=toBeAdded->size() ) std::cout << "TTT ARGGGHGH " << nMatched << std::endl;
208 
209  iEvent.put(newCol, "tracks");
210 }
211 // ------------ method called once each job just before starting event loop ------------
212 void
214 {
215 }
216 
217 // ------------ method called once each job just after ending the event loop ------------
218 void
220 }
221 
222 //define this as a plug-in
virtual void endJob()
virtual void beginJob()
void produceTrackColl(edm::Event &, const std::vector< reco::Particle::LorentzVector > *toBeAdded)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:10
#define abs(x)
Definition: mlp_lapack.h:159
PFCandidateCollection::const_iterator PFCandidateConstIterator
iterator
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
ZmumuPFEmbedder(const edm::ParameterSet &)
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
edm::InputTag _tracks
tuple cout
Definition: gather_cfg.py:121
void producePFCandColl(edm::Event &, const std::vector< reco::Particle::LorentzVector > *toBeAdded)
const double epsilon
virtual void produce(edm::Event &, const edm::EventSetup &)
edm::InputTag _selectedMuons