Go to the documentation of this file.00001
00002 #include <memory>
00003
00004
00005 #include "FWCore/Framework/interface/Frameworkfwd.h"
00006
00007 #include "DataFormats/Common/interface/Handle.h"
00008 #include "DataFormats/Common/interface/EDProduct.h"
00009 #include "FWCore/Framework/interface/EDProducer.h"
00010 #include "FWCore/Framework/interface/Event.h"
00011 #include "FWCore/Framework/interface/EventSetup.h"
00012 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00013
00014 #include "FWCore/Framework/interface/MakerMacros.h"
00015 #include "FWCore/Framework/interface/EDFilter.h"
00016
00017 #include "FWCore/Framework/interface/ESHandle.h"
00018 #include "Geometry/CaloEventSetup/interface/CaloTopologyRecord.h"
00019 #include "Geometry/Records/interface/CaloGeometryRecord.h"
00020 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
00021 #include "DataFormats/EgammaCandidates/interface/GsfElectronFwd.h"
00022 #include "DataFormats/PatCandidates/interface/Electron.h"
00023 #include "DataFormats/VertexReco/interface/Vertex.h"
00024 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00025 #include "EgammaAnalysis/ElectronTools/interface/ElectronEnergyRegressionEvaluate.h"
00026
00027
00028
00029
00030
00031 using namespace std;
00032 using namespace reco;
00033 using namespace edm;
00034
00035 class ElectronRegressionEnergyProducer : public edm::EDFilter {
00036 public:
00037 explicit ElectronRegressionEnergyProducer(const edm::ParameterSet&);
00038 ~ElectronRegressionEnergyProducer();
00039 private:
00040 virtual bool filter(edm::Event&, const edm::EventSetup&);
00041
00042
00043 bool printDebug_;
00044 edm::InputTag electronTag_;
00045
00046 std::string regressionInputFile_;
00047 uint32_t energyRegressionType_;
00048
00049 std::string nameEnergyReg_;
00050 std::string nameEnergyErrorReg_;
00051
00052 bool geomInitialized_;
00053
00054 const CaloTopology* ecalTopology_;
00055 const CaloGeometry* caloGeometry_;
00056
00057 edm::InputTag recHitCollectionEB_;
00058 edm::InputTag recHitCollectionEE_;
00059
00060 ElectronEnergyRegressionEvaluate *regressionEvaluator;
00061
00062 };
00063
00064
00065 ElectronRegressionEnergyProducer::ElectronRegressionEnergyProducer(const edm::ParameterSet& iConfig) {
00066 printDebug_ = iConfig.getUntrackedParameter<bool>("printDebug", false);
00067 electronTag_ = iConfig.getParameter<edm::InputTag>("electronTag");
00068
00069 regressionInputFile_ = iConfig.getParameter<std::string>("regressionInputFile");
00070 energyRegressionType_ = iConfig.getParameter<uint32_t>("energyRegressionType");
00071
00072 nameEnergyReg_ = iConfig.getParameter<std::string>("nameEnergyReg");
00073 nameEnergyErrorReg_ = iConfig.getParameter<std::string>("nameEnergyErrorReg");
00074
00075 recHitCollectionEB_ = iConfig.getParameter<edm::InputTag>("recHitCollectionEB");
00076 recHitCollectionEE_ = iConfig.getParameter<edm::InputTag>("recHitCollectionEE");
00077
00078 produces<edm::ValueMap<double> >(nameEnergyReg_);
00079 produces<edm::ValueMap<double> >(nameEnergyErrorReg_);
00080
00081 regressionEvaluator = new ElectronEnergyRegressionEvaluate();
00082
00083
00084 ElectronEnergyRegressionEvaluate::ElectronEnergyRegressionType type = ElectronEnergyRegressionEvaluate::kNoTrkVar;
00085 if (energyRegressionType_ == 1) type = ElectronEnergyRegressionEvaluate::kNoTrkVar;
00086 else if (energyRegressionType_ == 2) type = ElectronEnergyRegressionEvaluate::kWithSubCluVar;
00087 else if (energyRegressionType_ == 3) type = ElectronEnergyRegressionEvaluate::kWithTrkVarV1;
00088 else if (energyRegressionType_ == 4) type = ElectronEnergyRegressionEvaluate::kWithTrkVarV2;
00089
00090
00091 regressionEvaluator->initialize(regressionInputFile_.c_str(),type);
00092
00093 geomInitialized_ = false;
00094
00095 }
00096
00097
00098 ElectronRegressionEnergyProducer::~ElectronRegressionEnergyProducer()
00099 {
00100 delete regressionEvaluator;
00101 }
00102
00103
00104 bool ElectronRegressionEnergyProducer::filter(edm::Event& iEvent, const edm::EventSetup& iSetup) {
00105
00106 assert(regressionEvaluator->isInitialized());
00107
00108 if (!geomInitialized_) {
00109 edm::ESHandle<CaloTopology> theCaloTopology;
00110 iSetup.get<CaloTopologyRecord>().get(theCaloTopology);
00111 ecalTopology_ = & (*theCaloTopology);
00112
00113 edm::ESHandle<CaloGeometry> theCaloGeometry;
00114 iSetup.get<CaloGeometryRecord>().get(theCaloGeometry);
00115 caloGeometry_ = & (*theCaloGeometry);
00116 geomInitialized_ = true;
00117 }
00118
00119 std::auto_ptr<edm::ValueMap<double> > regrEnergyMap(new edm::ValueMap<double>() );
00120 edm::ValueMap<double>::Filler energyFiller(*regrEnergyMap);
00121
00122 std::auto_ptr<edm::ValueMap<double> > regrEnergyErrorMap(new edm::ValueMap<double>() );
00123 edm::ValueMap<double>::Filler energyErrorFiller(*regrEnergyErrorMap);
00124
00125 Handle<reco::GsfElectronCollection> egCollection;
00126 iEvent.getByLabel(electronTag_,egCollection);
00127 const reco::GsfElectronCollection egCandidates = (*egCollection.product());
00128
00129 std::vector<double> energyValues;
00130 std::vector<double> energyErrorValues;
00131 energyValues.reserve(egCollection->size());
00132 energyErrorValues.reserve(egCollection->size());
00133
00134
00135
00136
00137 edm::Handle< EcalRecHitCollection > pEBRecHits;
00138 edm::Handle< EcalRecHitCollection > pEERecHits;
00139 iEvent.getByLabel( recHitCollectionEB_, pEBRecHits );
00140 iEvent.getByLabel( recHitCollectionEE_, pEERecHits );
00141
00142
00143
00144
00145 Handle<reco::VertexCollection> hVertexProduct;
00146 iEvent.getByLabel(edm::InputTag("offlinePrimaryVertices"),hVertexProduct);
00147 const reco::VertexCollection inVertices = *(hVertexProduct.product());
00148
00149
00150 Int_t nvertices = 0;
00151 for (reco::VertexCollection::const_iterator inV = inVertices.begin();
00152 inV != inVertices.end(); ++inV) {
00153
00154
00155 if (inV->ndof() >= 4
00156 && inV->position().Rho() <= 2.0
00157 && fabs(inV->z()) <= 24.0
00158 ) {
00159 nvertices++;
00160 }
00161 }
00162
00163
00164
00165
00166 double rho = 0;
00167 Handle<double> hRhoKt6PFJets;
00168 iEvent.getByLabel(edm::InputTag("kt6PFJets","rho"), hRhoKt6PFJets);
00169 rho = (*hRhoKt6PFJets);
00170
00171
00172 for ( reco::GsfElectronCollection::const_iterator egIter = egCandidates.begin();
00173 egIter != egCandidates.end(); ++egIter) {
00174
00175 const EcalRecHitCollection * recHits=0;
00176 if(egIter->isEB())
00177 recHits = pEBRecHits.product();
00178 else
00179 recHits = pEERecHits.product();
00180
00181 SuperClusterHelper mySCHelper(&(*egIter),recHits,ecalTopology_,caloGeometry_);
00182
00183 double energy=regressionEvaluator->calculateRegressionEnergy(&(*egIter),
00184 mySCHelper,
00185 rho,nvertices,
00186 printDebug_);
00187
00188 double error=regressionEvaluator->calculateRegressionEnergyUncertainty(&(*egIter),
00189 mySCHelper,
00190 rho,nvertices,
00191 printDebug_);
00192
00193 energyValues.push_back(energy);
00194 energyErrorValues.push_back(error);
00195
00196 }
00197
00198 energyFiller.insert( egCollection, energyValues.begin(), energyValues.end() );
00199 energyFiller.fill();
00200
00201 energyErrorFiller.insert( egCollection, energyErrorValues.begin(), energyErrorValues.end() );
00202 energyErrorFiller.fill();
00203
00204 iEvent.put(regrEnergyMap,nameEnergyReg_);
00205 iEvent.put(regrEnergyErrorMap,nameEnergyErrorReg_);
00206
00207 return true;
00208
00209 }
00210
00211
00212
00213 DEFINE_FWK_MODULE(ElectronRegressionEnergyProducer);
00214
00215
00216