CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
RochesterCorrMuonProducerT.h
Go to the documentation of this file.
1 #ifndef TauAnalysis_MCEmbeddingTools_RochesterCorrMuonProducerT_h
2 #define TauAnalysis_MCEmbeddingTools_RochesterCorrMuonProducerT_h
3 
26 
28 
30 
31 #include <TLorentzVector.h>
32 
33 #include <vector>
34 
35 template<typename T>
37 {
38  typedef std::vector<T> MuonCollection;
39 
40  public:
41 
43  : moduleLabel_(cfg.getParameter<std::string>("@module_label"))
44  {
45  src_ = cfg.getParameter<edm::InputTag>("src");
46  isMC_ = cfg.getParameter<bool>("isMC");
47 
48  produces<MuonCollection>("");
49  }
51 
52  private:
53 
54  void produce(edm::Event& evt, const edm::EventSetup& es)
55  {
56  //std::cout << "<RochesterCorrMuonProducer>:" << std::endl;
57 
58  std::auto_ptr<MuonCollection> correctedMuonCollection(new MuonCollection);
59 
60  edm::Handle<MuonCollection> uncorrectedMuonCollection;
61  evt.getByLabel(src_, uncorrectedMuonCollection);
62 
63  int idx = 0;
64  for ( typename MuonCollection::const_iterator uncorrectedMuon = uncorrectedMuonCollection->begin();
65  uncorrectedMuon != uncorrectedMuonCollection->end(); ++uncorrectedMuon ) {
66  //std::cout << " uncorrected Muon #" << idx << ": Pt = " << uncorrectedMuon->pt() << ","
67  // << " eta = " << uncorrectedMuon->eta() << ", phi = " << uncorrectedMuon->phi() << ", mass = " << uncorrectedMuon->mass() << std::endl;
68 
69  TLorentzVector muonP4(uncorrectedMuon->px(), uncorrectedMuon->py(), uncorrectedMuon->pz(), uncorrectedMuon->energy());
70 
71  float error;
72  if ( isMC_ ) algorithm_.momcor_mc(muonP4, uncorrectedMuon->charge(), 0, 0, error);
73  else algorithm_.momcor_data(muonP4, uncorrectedMuon->charge(), 0, 0, error);
74 
75  T correctedMuon(*uncorrectedMuon);
76  // CV: use PolarLorentzVector to preserve muon mass (Rochester muon corrections change muon mass)
77  correctedMuon.setP4(reco::Candidate::PolarLorentzVector(muonP4.Pt(), muonP4.Eta(), muonP4.Phi(), uncorrectedMuon->mass()));
78  //std::cout << " corrected Muon #" << idx << ": Pt = " << correctedMuon.pt() << ","
79  // << " eta = " << correctedMuon.eta() << ", phi = " << correctedMuon.phi() << ", mass = " << correctedMuon.mass() << std::endl;
80 
81  correctedMuonCollection->push_back(correctedMuon);
82 
83  ++idx;
84  }
85 
86  evt.put(correctedMuonCollection);
87  }
88 
90 
92  bool isMC_;
93 
95 };
96 
97 #endif
98 
99 
T getParameter(std::string const &) const
RochesterCorrMuonProducerT(const edm::ParameterSet &cfg)
tuple cfg
Definition: looper.py:293
void produce(edm::Event &evt, const edm::EventSetup &es)
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:120
void momcor_mc(TLorentzVector &, float, float, int, float &)
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:420
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
void momcor_data(TLorentzVector &, float, float, int, float &)
long double T
math::PtEtaPhiMLorentzVector PolarLorentzVector
Lorentz vector.
Definition: Candidate.h:39