Go to the documentation of this file.00001
00002
00003
00004
00005
00006
00007
00008
00009
00018
00019
00020 #include "EgammaAnalysis/ElectronTools/plugins/CalibratedPatElectronProducer.h"
00021
00022 #include "FWCore/Framework/interface/Frameworkfwd.h"
00023 #include "FWCore/Framework/interface/MakerMacros.h"
00024 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00025 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00026 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
00027 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
00028 #include "FWCore/ParameterSet/interface/FileInPath.h"
00029
00030 #include "DataFormats/EgammaCandidates/interface/GsfElectronFwd.h"
00031 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
00032 #include "DataFormats/EgammaReco/interface/ElectronSeed.h"
00033 #include "DataFormats/GsfTrackReco/interface/GsfTrack.h"
00034 #include "DataFormats/TrackCandidate/interface/TrackCandidateCollection.h"
00035 #include "DataFormats/TrackingRecHit/interface/TrackingRecHitFwd.h"
00036 #include "DataFormats/EgammaCandidates/interface/GsfElectronFwd.h"
00037 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
00038 #include "DataFormats/PatCandidates/interface/Electron.h"
00039
00040 #include "EgammaAnalysis/ElectronTools/interface/ElectronEnergyCalibrator.h"
00041
00042 #include <iostream>
00043
00044 using namespace edm ;
00045 using namespace std ;
00046 using namespace reco ;
00047 using namespace pat ;
00048
00049 CalibratedPatElectronProducer::CalibratedPatElectronProducer( const edm::ParameterSet & cfg )
00050
00051 {
00052 produces<ElectronCollection>();
00053
00054 inputPatElectrons = cfg.getParameter<edm::InputTag>("inputPatElectronsTag");
00055 dataset = cfg.getParameter<std::string>("inputDataset");
00056 isMC = cfg.getParameter<bool>("isMC");
00057 updateEnergyError = cfg.getParameter<bool>("updateEnergyError");
00058 lumiRatio = cfg.getParameter<double>("lumiRatio");
00059 correctionsType = cfg.getParameter<int>("correctionsType");
00060 applyLinearityCorrection = cfg.getParameter<bool>("applyLinearityCorrection");
00061 combinationType = cfg.getParameter<int>("combinationType");
00062 verbose = cfg.getParameter<bool>("verbose");
00063 synchronization = cfg.getParameter<bool>("synchronization");
00064 combinationRegressionInputPath = cfg.getParameter<std::string>("combinationRegressionInputPath");
00065 scaleCorrectionsInputPath = cfg.getParameter<std::string>("scaleCorrectionsInputPath");
00066 linCorrectionsInputPath = cfg.getParameter<std::string>("linearityCorrectionsInputPath");
00067
00068
00069 if ( isMC && ( dataset != "Summer11" && dataset != "Fall11"
00070 && dataset != "Summer12" && dataset != "Summer12_DR53X_HCP2012"
00071 && dataset != "Summer12_LegacyPaper" ) )
00072 {
00073 throw cms::Exception("CalibratedPATElectronProducer|ConfigError") << "Unknown MC dataset";
00074 }
00075 if ( !isMC && ( dataset != "Prompt" && dataset != "ReReco"
00076 && dataset != "Jan16ReReco" && dataset != "ICHEP2012"
00077 && dataset != "Moriond2013" && dataset != "22Jan2013ReReco" ) )
00078 {
00079 throw cms::Exception("CalibratedPATElectronProducer|ConfigError") << "Unknown Data dataset";
00080 }
00081
00082
00083 if(combinationType!=3 && applyLinearityCorrection)
00084 {
00085 std::cout << "[CalibratedElectronProducer] "
00086 << "Warning: you chose combinationType!=3 and applyLinearityCorrection=True. Linearity corrections are only applied on top of combination 3." << std::endl;
00087 }
00088
00089
00090 std::cout << "[CalibratedPATElectronProducer] Correcting scale for dataset " << dataset << std::endl;
00091
00092
00093 std::string pathToDataCorr;
00094 switch ( correctionsType )
00095 {
00096 case 0:
00097 break;
00098 case 1:
00099 if ( verbose )
00100 {
00101 std::cout << "You choose regression 1 scale corrections" << std::endl;
00102 }
00103 break;
00104 case 2:
00105 if ( verbose )
00106 {
00107 std::cout << "You choose regression 2 scale corrections." << std::endl;
00108 }
00109 break;
00110 case 3:
00111 throw cms::Exception("CalibratedPATElectronProducer|ConfigError")
00112 << "You choose standard non-regression ecal energy scale corrections. They are not implemented yet.";
00113 break;
00114 default:
00115 throw cms::Exception("CalibratedPATElectronProducer|ConfigError")
00116 << "Unknown correctionsType !!!";
00117 }
00118
00119 theEnCorrector = new ElectronEnergyCalibrator
00120 (
00121 edm::FileInPath(scaleCorrectionsInputPath.c_str()).fullPath().c_str(),
00122 edm::FileInPath(linCorrectionsInputPath.c_str()).fullPath().c_str(),
00123 dataset,
00124 correctionsType,
00125 applyLinearityCorrection,
00126 lumiRatio,
00127 isMC,
00128 updateEnergyError,
00129 verbose,
00130 synchronization
00131 );
00132
00133 if ( verbose )
00134 {
00135 std::cout << "[CalibratedPATElectronProducer] "
00136 << "ElectronEnergyCalibrator object is created " << std::endl;
00137 }
00138
00139 myEpCombinationTool = new EpCombinationTool();
00140 myEpCombinationTool->init
00141 (
00142 edm::FileInPath(combinationRegressionInputPath.c_str()).fullPath().c_str(),
00143 "CombinationWeight"
00144 );
00145
00146 myCombinator = new ElectronEPcombinator();
00147
00148 if ( verbose )
00149 {
00150 std::cout << "[CalibratedPATElectronProducer] "
00151 << "Combination tools are created and initialized " << std::endl;
00152 }
00153 }
00154
00155
00156
00157 CalibratedPatElectronProducer::~CalibratedPatElectronProducer()
00158 {}
00159
00160 void CalibratedPatElectronProducer::produce( edm::Event & event, const edm::EventSetup & setup )
00161 {
00162
00163 edm::Handle<edm::View<reco::Candidate> > oldElectrons ;
00164 event.getByLabel(inputPatElectrons,oldElectrons) ;
00165 std::auto_ptr<ElectronCollection> electrons( new ElectronCollection ) ;
00166 ElectronCollection::const_iterator electron ;
00167 ElectronCollection::iterator ele ;
00168
00169 for
00170 (
00171 edm::View<reco::Candidate>::const_iterator ele=oldElectrons->begin();
00172 ele!=oldElectrons->end();
00173 ++ele
00174 )
00175 {
00176 const pat::ElectronRef elecsRef = edm::RefToBase<reco::Candidate>(oldElectrons,ele-oldElectrons->begin()).castTo<pat::ElectronRef>();
00177 pat::Electron clone = *edm::RefToBase<reco::Candidate>(oldElectrons,ele-oldElectrons->begin()).castTo<pat::ElectronRef>();
00178 electrons->push_back(clone);
00179 }
00180
00181 if (correctionsType != 0 )
00182 {
00183 for
00184 (
00185 ele = electrons->begin();
00186 ele != electrons->end() ;
00187 ++ele
00188 )
00189 {
00190 int elClass = -1;
00191 int run = event.run();
00192
00193 float r9 = ele->r9();
00194 double correctedEcalEnergy = ele->correctedEcalEnergy();
00195 double correctedEcalEnergyError = ele->correctedEcalEnergyError();
00196 double trackMomentum = ele->trackMomentumAtVtx().R();
00197 double trackMomentumError = ele->trackMomentumError();
00198 double combinedMomentum = ele->p();
00199 double combinedMomentumError = ele->p4Error(ele->candidateP4Kind());
00200
00201
00202
00203 if( !ele->ecalDrivenSeed() )
00204 {
00205 double error = 999. ;
00206 double momentum = (combinedMomentum<15. ? 15. : combinedMomentum);
00207 if ( ele->isEB() )
00208 {
00209 float parEB[3] = { 5.24e-02, 2.01e-01, 1.00e-02} ;
00210 error = momentum * sqrt( pow(parEB[0]/sqrt(momentum),2) + pow(parEB[1]/momentum,2) + pow(parEB[2],2) );
00211 }
00212 else if ( ele->isEE() )
00213 {
00214 float parEE[3] = { 1.46e-01, 9.21e-01, 1.94e-03} ;
00215 error = momentum * sqrt( pow(parEE[0]/sqrt(momentum),2) + pow(parEE[1]/momentum,2) + pow(parEE[2],2) );
00216 }
00217 combinedMomentumError = error;
00218 }
00219
00220 if (ele->classification() == reco::GsfElectron::GOLDEN) {elClass = 0;}
00221 if (ele->classification() == reco::GsfElectron::BIGBREM) {elClass = 1;}
00222 if (ele->classification() == reco::GsfElectron::BADTRACK) {elClass = 2;}
00223 if (ele->classification() == reco::GsfElectron::SHOWERING) {elClass = 3;}
00224 if (ele->classification() == reco::GsfElectron::GAP) {elClass = 4;}
00225
00226 SimpleElectron mySimpleElectron
00227 (
00228 run,
00229 elClass,
00230 r9,
00231 correctedEcalEnergy,
00232 correctedEcalEnergyError,
00233 trackMomentum,
00234 trackMomentumError,
00235 ele->ecalRegressionEnergy(),
00236 ele->ecalRegressionError(),
00237 combinedMomentum,
00238 combinedMomentumError,
00239 ele->superCluster()->eta(),
00240 ele->isEB(),
00241 isMC,
00242 ele->ecalDriven(),
00243 ele->trackerDrivenSeed()
00244 );
00245
00246
00247 if ( ele->core()->ecalDrivenSeed() || correctionsType==2 || combinationType==3 )
00248 {
00249 theEnCorrector->calibrate(mySimpleElectron);
00250
00251
00252
00253 switch ( combinationType )
00254 {
00255 case 0:
00256 if ( verbose )
00257 {
00258 std::cout << "[CalibratedPATElectronProducer] "
00259 << "You choose not to combine." << std::endl;
00260 }
00261 break;
00262 case 1:
00263 if ( verbose )
00264 {
00265 std::cout << "[CalibratedPATElectronProducer] "
00266 << "You choose corrected regression energy for standard combination" << std::endl;
00267 }
00268 myCombinator->setCombinationMode(1);
00269 myCombinator->combine(mySimpleElectron);
00270 break;
00271 case 2:
00272 if ( verbose )
00273 {
00274 std::cout << "[CalibratedPATElectronProducer] "
00275 << "You choose uncorrected regression energy for standard combination" << std::endl;
00276 }
00277 myCombinator->setCombinationMode(2);
00278 myCombinator->combine(mySimpleElectron);
00279 break;
00280 case 3:
00281 if ( verbose )
00282 {
00283 std::cout << "[CalibratedPATElectronProducer] "
00284 << "You choose regression combination." << std::endl;
00285 }
00286 myEpCombinationTool->combine(mySimpleElectron);
00287 theEnCorrector->correctLinearity(mySimpleElectron);
00288 break;
00289 default:
00290 throw cms::Exception("CalibratedPATElectronProducer|ConfigError")
00291 << "Unknown combination Type !!!" ;
00292 }
00293
00294 math::XYZTLorentzVector oldMomentum = ele->p4() ;
00295 math::XYZTLorentzVector newMomentum_ ;
00296 newMomentum_ = math::XYZTLorentzVector
00297 ( oldMomentum.x()*mySimpleElectron.getCombinedMomentum()/oldMomentum.t(),
00298 oldMomentum.y()*mySimpleElectron.getCombinedMomentum()/oldMomentum.t(),
00299 oldMomentum.z()*mySimpleElectron.getCombinedMomentum()/oldMomentum.t(),
00300 mySimpleElectron.getCombinedMomentum() ) ;
00301
00302 ele->correctMomentum
00303 (
00304 newMomentum_,
00305 mySimpleElectron.getTrackerMomentumError(),
00306 mySimpleElectron.getCombinedMomentumError()
00307 );
00308
00309 if ( verbose )
00310 {
00311 std::cout << "[CalibratedPATElectronProducer] Combined momentum after saving "
00312 << ele->p4().t() << std::endl;
00313 }
00314 }
00315 }
00316 } else
00317 {
00318 if ( verbose )
00319 {
00320 std::cout << "[CalibratedPATElectronProducer] "
00321 << "You choose not to correct. Uncorrected Regression Energy is taken." << std::endl;
00322 }
00323 }
00324
00325 event.put(electrons) ;
00326 }
00327
00328 #include "FWCore/Framework/interface/MakerMacros.h"
00329 #include "FWCore/Framework/interface/ESProducer.h"
00330 #include "FWCore/Framework/interface/ModuleFactory.h"
00331 DEFINE_FWK_MODULE(CalibratedPatElectronProducer);