ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/Morgan/src/METAnalyzer.cc
Revision: 1.6
Committed: Wed Jun 10 11:17:06 2009 UTC (15 years, 10 months ago) by lethuill
Content type: text/plain
Branch: MAIN
CVS Tags: all_3_3_2_01, all_3_2_5_02, all_3_2_5_01, all_2_2_9_03, all_2_2_9_02, all_2_2_9_01, HEAD
Branch point for: CMSSW_2_2_X_br
Changes since 1.5: +164 -137 lines
Log Message:
Better protection against missing collection / Cleaning data format selection / Last iteration for migration to PAT of Photons

File Contents

# User Rev Content
1 lethuill 1.2 #include "../interface/METAnalyzer.h"
2 lethuill 1.1
3     using namespace std;
4     using namespace reco;
5     using namespace edm;
6    
7    
8 lethuill 1.3 METAnalyzer::METAnalyzer(const edm::ParameterSet& producersNames, const edm::ParameterSet& myConfig, int verbosity):verbosity_(verbosity)
9 lethuill 1.1 {
10 lethuill 1.6 dataType_ = producersNames.getUntrackedParameter<string>("dataType","unknown");
11     metProducer_ = producersNames.getParameter<edm::InputTag>("metProducer");
12     useMC_ = myConfig.getUntrackedParameter<bool>("doMETMC");
13     allowMissingCollection_ = producersNames.getUntrackedParameter<bool>("allowMissingCollection", false);
14 lethuill 1.1 }
15    
16     METAnalyzer::~METAnalyzer()
17     {
18     }
19    
20 lethuill 1.6 bool METAnalyzer::process(const edm::Event& iEvent, TClonesArray* rootMET)
21 lethuill 1.1 {
22 lethuill 1.6
23     unsigned int nMETs=0;
24    
25     edm::Handle < std::vector <reco::CaloMET> > recoMETs;
26     if( dataType_=="RECO" )
27     {
28     {
29     try
30     {
31     iEvent.getByLabel(metProducer_, recoMETs);
32     nMETs = recoMETs->size();
33     }
34     catch (cms::Exception& exception)
35     {
36     if ( !allowMissingCollection_ )
37     {
38     cout << " ##### ERROR IN METAnalyzer::process => reco::CaloMET collection is missing #####"<<endl;
39     throw exception;
40     }
41     if(verbosity_>1) cout << " ===> No reco::CaloMET collection, skip MET info" << endl;
42     return false;
43     }
44     }
45     }
46    
47     edm::Handle < std::vector <pat::MET> > patMETs;
48     if( dataType_=="PAT" )
49     {
50     {
51     try
52     {
53     iEvent.getByLabel(metProducer_, patMETs);
54     nMETs = patMETs->size();
55     }
56     catch (cms::Exception& exception)
57     {
58     if ( !allowMissingCollection_ )
59     {
60     cout << " ##### ERROR IN METAnalyzer::process => pat::MET collection is missing #####"<<endl;
61     throw exception;
62     }
63     if(verbosity_>1) cout << " ===> No pat::MET collection, skip MET info" << endl;
64     return false;
65     }
66     }
67     }
68    
69     if(verbosity_>1) std::cout << " Number of MET objects = " << nMETs << " Label: " << metProducer_.label() << " Instance: " << metProducer_.instance() << std::endl;
70    
71     for (unsigned int j=0; j<nMETs; j++)
72     {
73     const reco::Candidate* met = 0;
74     if( dataType_=="RECO" ) met = (const reco::Candidate*) ( & ((*recoMETs)[j]) );
75     if( dataType_=="PAT" ) met = (const reco::Candidate*) ( & ((*patMETs)[j]) );
76    
77     TRootMET localMET(
78     met->px()
79     ,met->py()
80     ,met->pz()
81     ,met->energy()
82     ,met->vx()
83     ,met->vy()
84     ,met->vz()
85     );
86    
87    
88     if( dataType_=="RECO" )
89     {
90     // Some specific methods to reco::MET
91     /*
92     const reco::CaloMET *recoMET = dynamic_cast<const reco::CaloMET*>(&*met);
93     localMET.setCaloMETFraction(
94     recoMET->maxEtInEmTowers()
95     ,recoMET->maxEtInHadTowers()
96     ,recoMET->hadEtInHO()
97     ,recoMET->hadEtInHB()
98     ,recoMET->hadEtInHF()
99     ,recoMET->hadEtInHE()
100     ,recoMET->emEtInEB()
101     ,recoMET->emEtInEE()
102     ,recoMET->emEtInHF()
103     ,recoMET->etFractionHadronic()
104     ,recoMET->emEtFraction()
105     ,recoMET->metSignificance()
106     ,recoMET->CaloMETInpHF()
107     ,recoMET->CaloMETInmHF()
108     ,recoMET->CaloSETInpHF()
109     ,recoMET->CaloSETInmHF()
110     ,recoMET->CaloMETPhiInpHF()
111     ,recoMET->CaloMETPhiInmHF()
112     );
113     */
114    
115     }
116    
117     if( dataType_=="PAT" )
118     {
119     // Some specific methods to pat::MET
120     const pat::MET *patMET = dynamic_cast<const pat::MET*>(&*met);
121    
122     localMET.setCaloMETFraction(
123     patMET->maxEtInEmTowers()
124     ,patMET->maxEtInHadTowers()
125     ,patMET->hadEtInHO()
126     ,patMET->hadEtInHB()
127     ,patMET->hadEtInHF()
128     ,patMET->hadEtInHE()
129     ,patMET->emEtInEB()
130     ,patMET->emEtInEE()
131     ,patMET->emEtInHF()
132     ,patMET->etFractionHadronic()
133     ,patMET->emEtFraction()
134     ,patMET->metSignificance()
135     ,patMET->CaloMETInpHF()
136     ,patMET->CaloMETInmHF()
137     ,patMET->CaloSETInpHF()
138     ,patMET->CaloSETInmHF()
139     ,patMET->CaloMETPhiInpHF()
140     ,patMET->CaloMETPhiInmHF()
141     );
142    
143     //pat::MET::UncorectionType ix;
144     //ix = pat::MET::uncorrALL;
145     localMET.setCorExALL(patMET->corEx(pat::MET::uncorrALL));
146     localMET.setCorExJES(patMET->corEx(pat::MET::uncorrJES));
147     localMET.setCorExMUON(patMET->corEx(pat::MET::uncorrMUON));
148    
149     localMET.setCorEyALL(patMET->corEy(pat::MET::uncorrALL));
150     localMET.setCorEyJES(patMET->corEy(pat::MET::uncorrJES));
151     localMET.setCorEyMUON(patMET->corEy(pat::MET::uncorrMUON));
152    
153     localMET.setCorSumEtALL(patMET->corSumEt(pat::MET::uncorrALL));
154     localMET.setCorSumEtJES(patMET->corSumEt(pat::MET::uncorrJES));
155     localMET.setCorSumEtMUON(patMET->corSumEt(pat::MET::uncorrMUON));
156    
157     localMET.setUncorrectedPhiALL(patMET->uncorrectedPhi(pat::MET::uncorrALL));
158     localMET.setUncorrectedPhiJES(patMET->uncorrectedPhi(pat::MET::uncorrJES));
159     localMET.setUncorrectedPhiMUON(patMET->uncorrectedPhi(pat::MET::uncorrMUON));
160    
161     localMET.setUncorrectedPtALL(patMET->uncorrectedPt(pat::MET::uncorrALL));
162     localMET.setUncorrectedPtJES(patMET->uncorrectedPt(pat::MET::uncorrJES));
163     localMET.setUncorrectedPtMUON(patMET->uncorrectedPt(pat::MET::uncorrMUON));
164    
165     if(useMC_)
166     {
167     // MC truth associator index
168     if ((patMET->genParticleRef()).isNonnull()) {
169     localMET.setGenParticleIndex((patMET->genParticleRef()).index());
170     } else {
171     localMET.setGenParticleIndex(-1);
172     }
173     }
174     }
175    
176     new( (*rootMET)[j] ) TRootMET(localMET);
177     if(verbosity_>2) cout << " ["<< setw(3) << j << "] " << localMET << endl;
178     }
179    
180     return true;
181 lethuill 1.1 }