Go to the documentation of this file.00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015 #include <memory>
00016 #include "RecoMET/METProducers/interface/METProducer.h"
00017 #include "RecoMET/METAlgorithms/interface/SignCaloSpecificAlgo.h"
00018 #include "RecoMET/METAlgorithms/interface/SignAlgoResolutions.h"
00019 #include "RecoMET/METAlgorithms/interface/CaloSpecificAlgo.h"
00020 #include "RecoMET/METAlgorithms/interface/GenSpecificAlgo.h"
00021 #include "RecoMET/METAlgorithms/interface/TCMETAlgo.h"
00022 #include "RecoMET/METAlgorithms/interface/PFSpecificAlgo.h"
00023 #include "DataFormats/JetReco/interface/CaloJet.h"
00024 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
00025 #include "DataFormats/METReco/interface/CommonMETData.h"
00026 #include "DataFormats/Candidate/interface/LeafCandidate.h"
00027 #include "DataFormats/Candidate/interface/Particle.h"
00028 #include "DataFormats/Candidate/interface/Candidate.h"
00029 #include "DataFormats/Candidate/interface/CandidateFwd.h"
00030 #include "DataFormats/Common/interface/View.h"
00031
00032 using namespace edm;
00033 using namespace std;
00034 using namespace reco;
00035
00036 namespace cms
00037 {
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049
00050
00051 METProducer::METProducer(const edm::ParameterSet& iConfig) : conf_(iConfig), alg_() , resolutions_(0)
00052 {
00053 inputLabel = iConfig.getParameter<edm::InputTag>("src");
00054 inputType = iConfig.getParameter<std::string>("InputType");
00055 METtype = iConfig.getParameter<std::string>("METType");
00056 alias = iConfig.getParameter<std::string>("alias");
00057 globalThreshold = iConfig.getParameter<double>("globalThreshold");
00058 calculateSignificance_ = false ;
00059
00060 if( METtype == "CaloMET" )
00061 {
00062 noHF = iConfig.getParameter<bool>("noHF");
00063 produces<CaloMETCollection>().setBranchAlias(alias.c_str());
00064 calculateSignificance_ = iConfig.getParameter<bool>("calculateSignificance");
00065 }
00066 else if( METtype == "GenMET" )
00067 {
00068 onlyFiducial = iConfig.getParameter<bool>("onlyFiducialParticles");
00069 usePt = iConfig.getUntrackedParameter<bool>("usePt",false);
00070 produces<GenMETCollection>().setBranchAlias(alias.c_str());
00071 }
00072 else if( METtype == "PFMET" )
00073 {
00074 produces<PFMETCollection>().setBranchAlias(alias.c_str());
00075
00076 calculateSignificance_ = iConfig.getParameter<bool>("calculateSignificance");
00077
00078 if(calculateSignificance_){
00079 jetsLabel_ = iConfig.getParameter<edm::InputTag>("jets");
00080 }
00081
00082 }
00083 else if (METtype == "TCMET" )
00084 {
00085 produces<METCollection>().setBranchAlias(alias.c_str());
00086 TCMETAlgo ALGO;
00087
00088 int rfType_ = iConfig.getParameter<int>("rf_type");
00089 bool correctShowerTracks_ = iConfig.getParameter<bool>("correctShowerTracks");
00090
00091 if(correctShowerTracks_){
00092 showerRF_ = (*ALGO.getResponseFunction_shower());
00093 responseFunction_ = (*ALGO.getResponseFunction_noshower());
00094 }else{
00095
00096 if( rfType_ == 1 )
00097 responseFunction_ = (*ALGO.getResponseFunction_fit());
00098 else if( rfType_ == 2 )
00099 responseFunction_ = (*ALGO.getResponseFunction_mode());
00100
00101 }
00102 }
00103 else
00104 produces<METCollection>().setBranchAlias(alias.c_str());
00105
00106 if (calculateSignificance_ && ( METtype == "CaloMET" || METtype == "PFMET")){
00107 resolutions_ = new metsig::SignAlgoResolutions(iConfig);
00108
00109 }
00110 }
00111
00112
00113
00114
00115
00116 METProducer::METProducer() : alg_()
00117 {
00118 produces<METCollection>();
00119 }
00120
00121
00122
00123
00124
00125 METProducer::~METProducer() {}
00126
00127
00128
00129
00130
00131 void METProducer::produce(Event& event, const EventSetup& setup)
00132 {
00133
00134
00135
00136 edm::Handle<edm::View<Candidate> > input;
00137 event.getByLabel(inputLabel,input);
00138
00139
00140 CommonMETData output;
00141
00142
00143
00144
00145
00146
00147
00148
00149
00150
00151
00152
00153
00154
00155
00156
00157
00158 if( METtype == "CaloMET" )
00159 {
00160
00161 alg_.run(input, &output, globalThreshold);
00162
00163
00164 CaloSpecificAlgo calospecalgo;
00165 CaloMET calomet = calospecalgo.addInfo(input,output,noHF, globalThreshold);
00166
00167
00168 if( calculateSignificance_ )
00169 {
00170 SignCaloSpecificAlgo signcalospecalgo;
00171
00172
00173 signcalospecalgo.calculateBaseCaloMET(input,output,*resolutions_,noHF,globalThreshold);
00174 calomet.SetMetSignificance( signcalospecalgo.getSignificance() );
00175 calomet.setSignificanceMatrix(signcalospecalgo.getSignificanceMatrix());
00176 }
00177
00178 std::auto_ptr<CaloMETCollection> calometcoll;
00179 calometcoll.reset(new CaloMETCollection);
00180 calometcoll->push_back( calomet ) ;
00181 event.put( calometcoll );
00182
00183 }
00184
00185 else if( METtype == "TCMET" )
00186 {
00187 TCMETAlgo tcmetalgorithm;
00188 std::auto_ptr<METCollection> tcmetcoll;
00189 tcmetcoll.reset(new METCollection);
00190 tcmetcoll->push_back( tcmetalgorithm.CalculateTCMET(event, setup, conf_, &responseFunction_, &showerRF_) ) ;
00191 event.put( tcmetcoll );
00192 }
00193
00194 else if( METtype == "PFMET" )
00195 {
00196 alg_.run(input, &output, globalThreshold);
00197 PFSpecificAlgo pf;
00198 std::auto_ptr<PFMETCollection> pfmetcoll;
00199 pfmetcoll.reset (new PFMETCollection);
00200
00201
00202 if( calculateSignificance_ )
00203 {
00204
00205 edm::Handle<edm::View<reco::PFJet> > jets;
00206 event.getByLabel(jetsLabel_,jets);
00207 pf.runSignificance(*resolutions_, jets);
00208 }
00209 pfmetcoll->push_back( pf.addInfo(input, output) );
00210 event.put( pfmetcoll );
00211 }
00212
00213 else if( METtype == "GenMET" )
00214 {
00215 GenSpecificAlgo gen;
00216 std::auto_ptr<GenMETCollection> genmetcoll;
00217 genmetcoll.reset (new GenMETCollection);
00218 genmetcoll->push_back( gen.addInfo(input, &output, globalThreshold, onlyFiducial, usePt) );
00219 event.put( genmetcoll );
00220 }
00221 else
00222 {
00223 alg_.run(input, &output, globalThreshold);
00224 LorentzVector p4( output.mex, output.mey, 0.0, output.met);
00225 Point vtx(0,0,0);
00226 MET met( output.sumet, p4, vtx );
00227 std::auto_ptr<METCollection> metcoll;
00228 metcoll.reset(new METCollection);
00229 metcoll->push_back( met );
00230 event.put( metcoll );
00231 }
00232
00233 }
00234
00235 }