ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/MitProd/TreeFiller/src/FillerPFMet.cc
Revision: 1.7
Committed: Mon Oct 3 16:15:50 2011 UTC (13 years, 7 months ago) by ksung
Content type: text/plain
Branch: MAIN
CVS Tags: Mit_025b, Mit_025a, Mit_025
Changes since 1.6: +13 -5 lines
Log Message:
changes to handle clustered PFMET variables

File Contents

# User Rev Content
1 ksung 1.7 // $Id: FillerPFMet.cc,v 1.6 2011/09/12 13:48:16 ksung Exp $
2 bendavid 1.1
3     #include "MitProd/TreeFiller/interface/FillerPFMet.h"
4     #include "DataFormats/METReco/interface/PFMET.h"
5 loizides 1.3 #include "DataFormats/METReco/interface/PFMETCollection.h"
6 bendavid 1.1 #include "MitAna/DataTree/interface/Names.h"
7 loizides 1.3 #include "MitAna/DataTree/interface/PFMetCol.h"
8     #include "MitProd/ObjectService/interface/ObjectService.h"
9 bendavid 1.1
10     using namespace std;
11     using namespace edm;
12     using namespace mithep;
13    
14     //--------------------------------------------------------------------------------------------------
15     FillerPFMet::FillerPFMet(const ParameterSet &cfg, const char *name, bool active) :
16     BaseFiller(cfg,name,active),
17 ksung 1.7 edmName_(Conf().getUntrackedParameter<edm::InputTag>("edmName")),
18 bendavid 1.1 mitName_(Conf().getUntrackedParameter<string>("mitName",Names::gkCaloMetBrn)),
19     pfMets_(new mithep::PFMetArr)
20     {
21     // Constructor.
22     }
23    
24     //--------------------------------------------------------------------------------------------------
25     FillerPFMet::~FillerPFMet()
26     {
27     // Destructor.
28    
29     delete pfMets_;
30     }
31    
32     //--------------------------------------------------------------------------------------------------
33 bendavid 1.5 void FillerPFMet::BookDataBlock(TreeWriter &tws)
34 bendavid 1.1 {
35     // Add mets branch to tree.
36    
37 loizides 1.2 tws.AddBranch(mitName_,&pfMets_);
38     OS()->add<mithep::PFMetArr>(pfMets_,mitName_);
39 bendavid 1.1 }
40    
41     //--------------------------------------------------------------------------------------------------
42     void FillerPFMet::FillDataBlock(const edm::Event &event,
43 loizides 1.2 const edm::EventSetup &setup)
44 bendavid 1.1 {
45     // Fill missing energy from edm collection into our collection.
46    
47     pfMets_->Delete();
48    
49 ksung 1.7 reco::PFMETCollection inPFMets;
50    
51 bendavid 1.1 Handle<reco::PFMETCollection> hCaloMetProduct;
52 ksung 1.7 //GetProduct(edmName_, hCaloMetProduct, event);
53     event.getByLabel(edmName_,hCaloMetProduct);
54     if(hCaloMetProduct.isValid()) {
55     inPFMets = *(hCaloMetProduct.product());
56     } else {
57     Handle<reco::PFMET> hSingleMetProduct;
58     event.getByLabel(edmName_, hSingleMetProduct);
59     inPFMets.push_back(*hSingleMetProduct.product());
60     }
61 bendavid 1.1
62     // loop through all mets
63     for (reco::PFMETCollection::const_iterator inPFMet = inPFMets.begin();
64     inPFMet != inPFMets.end(); ++inPFMet) {
65    
66     mithep::PFMet *pfMet = pfMets_->Allocate();
67     new (pfMet) mithep::PFMet(inPFMet->px(), inPFMet->py());
68    
69 loizides 1.2 // fill Met base class data
70 bendavid 1.1 pfMet->SetSumEt(inPFMet->sumEt());
71     pfMet->SetElongitudinal(inPFMet->e_longitudinal());
72     for(unsigned i=0; i<inPFMet->mEtCorr().size(); i++) {
73     pfMet->PushCorrectionX(inPFMet->mEtCorr()[i].mex);
74     pfMet->PushCorrectionY(inPFMet->mEtCorr()[i].mey);
75     pfMet->PushCorrectionSumEt(inPFMet->mEtCorr()[i].sumet);
76     }
77    
78 loizides 1.2 // fill PFMet class data
79 ksung 1.6 pfMet->SetPFMetSig(inPFMet->significance());
80 bendavid 1.1 pfMet->SetNeutralEMFraction(inPFMet->NeutralEMFraction());
81     pfMet->SetNeutralHadFraction(inPFMet->NeutralHadFraction());
82     pfMet->SetChargedEMFraction(inPFMet->ChargedEMFraction());
83     pfMet->SetChargedHadFraction(inPFMet->ChargedHadFraction());
84     pfMet->SetMuonFraction(inPFMet->MuonFraction());
85    
86     }
87     pfMets_->Trim();
88     }