CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
Type0PFMETcorrInputProducer Class Reference

#include <Type0PFMETcorrInputProducer.h>

Inheritance diagram for Type0PFMETcorrInputProducer:
edm::EDProducer edm::ProducerBase edm::ProductRegistryHelper

Public Member Functions

 Type0PFMETcorrInputProducer (const edm::ParameterSet &)
 
 ~Type0PFMETcorrInputProducer ()
 
- Public Member Functions inherited from edm::EDProducer
 EDProducer ()
 
virtual ~EDProducer ()
 
- Public Member Functions inherited from edm::ProducerBase
 ProducerBase ()
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
boost::function< void(const
BranchDescription &)> 
registrationCallback () const
 used by the fwk to register list of products More...
 
virtual ~ProducerBase ()
 

Private Member Functions

void produce (edm::Event &, const edm::EventSetup &)
 

Private Attributes

TFormula * correction_
 
double minDz_
 
std::string moduleLabel_
 
edm::InputTag srcHardScatterVertex_
 
edm::InputTag srcPFCandidateToVertexAssociations_
 

Additional Inherited Members

- Public Types inherited from edm::EDProducer
typedef EDProducer ModuleType
 
typedef WorkerT< EDProducerWorkerType
 
- Public Types inherited from edm::ProducerBase
typedef
ProductRegistryHelper::TypeLabelList 
TypeLabelList
 
- Static Public Member Functions inherited from edm::EDProducer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::EDProducer
CurrentProcessingContext const * currentContext () const
 
- Protected Member Functions inherited from edm::ProducerBase
template<class TProducer , class TMethod >
void callWhenNewProductsRegistered (TProducer *iProd, TMethod iMethod)
 

Detailed Description

Compute Type 0 (PF)MET corrections ( https://indico.cern.ch/getFile.py/access?contribId=4&resId=1&materialId=slides&confId=161159 )

Author
Christian Veelken, LLR
Version
Revision:
1.1
Id:
Type0PFMETcorrInputProducer.h,v 1.1 2012/04/19 17:59:40 veelken Exp

Definition at line 27 of file Type0PFMETcorrInputProducer.h.

Constructor & Destructor Documentation

Type0PFMETcorrInputProducer::Type0PFMETcorrInputProducer ( const edm::ParameterSet cfg)
explicit

Definition at line 16 of file Type0PFMETcorrInputProducer.cc.

References correction_, edm::ParameterSet::getParameter(), minDz_, moduleLabel_, srcHardScatterVertex_, and srcPFCandidateToVertexAssociations_.

17  : moduleLabel_(cfg.getParameter<std::string>("@module_label")),
18  correction_(0)
19 {
20  srcPFCandidateToVertexAssociations_ = cfg.getParameter<edm::InputTag>("srcPFCandidateToVertexAssociations");
21  srcHardScatterVertex_ = cfg.getParameter<edm::InputTag>("srcHardScatterVertex");
22 
23  edm::ParameterSet cfgCorrection_function = cfg.getParameter<edm::ParameterSet>("correction");
24  std::string corrFunctionName = std::string(moduleLabel_).append("correction");
25  std::string corrFunctionFormula = cfgCorrection_function.getParameter<std::string>("formula");
26  correction_ = new TFormula(corrFunctionName.data(), corrFunctionFormula.data());
27  int numParameter = correction_->GetNpar();
28  for ( int iParameter = 0; iParameter < numParameter; ++iParameter ) {
29  std::string parName = Form("par%i", iParameter);
30  double parValue = cfgCorrection_function.getParameter<double>(parName);
31  correction_->SetParameter(iParameter, parValue);
32  }
33 
34  minDz_ = cfg.getParameter<double>("minDz");
35 
36  produces<CorrMETData>();
37 }
T getParameter(std::string const &) const
Type0PFMETcorrInputProducer::~Type0PFMETcorrInputProducer ( )

Definition at line 39 of file Type0PFMETcorrInputProducer.cc.

References correction_.

40 {
41  delete correction_;
42 }

Member Function Documentation

void Type0PFMETcorrInputProducer::produce ( edm::Event evt,
const edm::EventSetup es 
)
privatevirtual

Implements edm::EDProducer.

Definition at line 44 of file Type0PFMETcorrInputProducer.cc.

References correction_, reco::PFCandidate::e, edm::Event::getByLabel(), reco::PFCandidate::h, minDz_, reco::PFCandidate::mu, reco::LeafCandidate::p4(), reco::PFCandidate::particleId(), pfMETCorrectionType0_cfi::pfCandidateToVertexAssociation, phi, edm::Event::put(), srcHardScatterVertex_, and srcPFCandidateToVertexAssociations_.

45 {
46  edm::Handle<reco::VertexCollection> hardScatterVertex;
47  evt.getByLabel(srcHardScatterVertex_, hardScatterVertex);
48 
49  edm::Handle<PFCandToVertexAssMap> pfCandidateToVertexAssociations;
50  evt.getByLabel(srcPFCandidateToVertexAssociations_, pfCandidateToVertexAssociations);
51 
52  std::auto_ptr<CorrMETData> pfMEtCorrection(new CorrMETData());
53 
54  for ( PFCandToVertexAssMap::const_iterator pfCandidateToVertexAssociation = pfCandidateToVertexAssociations->begin();
55  pfCandidateToVertexAssociation != pfCandidateToVertexAssociations->end(); ++pfCandidateToVertexAssociation ) {
57  const PFCandQualityPairVector& pfCandidates_vertex = pfCandidateToVertexAssociation->val;
58 
59  bool isHardScatterVertex = false;
60  for ( reco::VertexCollection::const_iterator hardScatterVertex_i = hardScatterVertex->begin();
61  hardScatterVertex_i != hardScatterVertex->end(); ++hardScatterVertex_i ) {
62  if ( TMath::Abs(vertex->position().z() - hardScatterVertex_i->position().z()) < minDz_ ) {
63  isHardScatterVertex = true;
64  break;
65  }
66  }
67 
68  if ( !isHardScatterVertex ) {
69  reco::Candidate::LorentzVector sumChargedPFCandP4_vertex;
70  for ( PFCandQualityPairVector::const_iterator pfCandidate_vertex = pfCandidates_vertex.begin();
71  pfCandidate_vertex != pfCandidates_vertex.end(); ++pfCandidate_vertex ) {
72  const reco::PFCandidate& pfCandidate = (*pfCandidate_vertex->first);
73  if ( pfCandidate.particleId() == reco::PFCandidate::h ||
74  pfCandidate.particleId() == reco::PFCandidate::e ||
75  pfCandidate.particleId() == reco::PFCandidate::mu ) {
76  sumChargedPFCandP4_vertex += pfCandidate.p4();
77  }
78  }
79 
80  double pt = sumChargedPFCandP4_vertex.pt();
81  double phi = sumChargedPFCandP4_vertex.phi();
82  double ptCorr = correction_->Eval(pt);
83  double pxCorr = TMath::Cos(phi)*ptCorr;
84  double pyCorr = TMath::Sin(phi)*ptCorr;
85 
86  pfMEtCorrection->mex += pxCorr;
87  pfMEtCorrection->mey += pyCorr;
88  }
89  }
90 
91  evt.put(pfMEtCorrection);
92 }
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:85
std::vector< PFCandQualityPair > PFCandQualityPairVector
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:38
a MET correction term
Definition: CorrMETData.h:14
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:33
virtual ParticleType particleId() const
Definition: PFCandidate.h:324
virtual const LorentzVector & p4() const
four-momentum Lorentz vector
Definition: DDAxes.h:10

Member Data Documentation

TFormula* Type0PFMETcorrInputProducer::correction_
private
double Type0PFMETcorrInputProducer::minDz_
private

Definition at line 45 of file Type0PFMETcorrInputProducer.h.

Referenced by produce(), and Type0PFMETcorrInputProducer().

std::string Type0PFMETcorrInputProducer::moduleLabel_
private
edm::InputTag Type0PFMETcorrInputProducer::srcHardScatterVertex_
private

Definition at line 41 of file Type0PFMETcorrInputProducer.h.

Referenced by produce(), and Type0PFMETcorrInputProducer().

edm::InputTag Type0PFMETcorrInputProducer::srcPFCandidateToVertexAssociations_
private

Definition at line 40 of file Type0PFMETcorrInputProducer.h.

Referenced by produce(), and Type0PFMETcorrInputProducer().