CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_5_3_14/src/EgammaAnalysis/ElectronTools/plugins/RegressionEnergyPatElectronProducer.cc

Go to the documentation of this file.
00001 #include "EgammaAnalysis/ElectronTools/plugins/RegressionEnergyPatElectronProducer.h"
00002 
00003 #include "FWCore/Framework/interface/Frameworkfwd.h"
00004 #include "FWCore/Framework/interface/MakerMacros.h"
00005 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00006 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00007 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
00008 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
00009 #include "FWCore/Framework/interface/ESHandle.h"
00010 #include "FWCore/Utilities/interface/EDMException.h"
00011 #include "Geometry/CaloEventSetup/interface/CaloTopologyRecord.h"
00012 #include "Geometry/Records/interface/CaloGeometryRecord.h"
00013 #include "DataFormats/PatCandidates/interface/Electron.h"
00014 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00015 #include "DataFormats/VertexReco/interface/Vertex.h"
00016 #include "EgammaAnalysis/ElectronTools/interface/SuperClusterHelper.h"
00017 
00018 #include <iostream>
00019 
00020 using namespace edm ;
00021 using namespace std ;
00022 using namespace reco ;
00023 using namespace pat ;
00024 
00025 
00026 RegressionEnergyPatElectronProducer::RegressionEnergyPatElectronProducer( const edm::ParameterSet & cfg )
00027 {
00028 
00029   inputElectrons_ = cfg.getParameter<edm::InputTag>("inputElectronsTag");
00030   inputCollectionType_ = cfg.getParameter<uint32_t>("inputCollectionType");
00031   rhoInputTag_ = cfg.getParameter<edm::InputTag>("rhoCollection");
00032   verticesInputTag_ = cfg.getParameter<edm::InputTag>("vertexCollection");
00033   energyRegressionType_ = cfg.getParameter<uint32_t>("energyRegressionType");
00034   regressionInputFile_ = cfg.getParameter<std::string>("regressionInputFile");
00035   recHitCollectionEB_ = cfg.getParameter<edm::InputTag>("recHitCollectionEB");
00036   recHitCollectionEE_ = cfg.getParameter<edm::InputTag>("recHitCollectionEE");
00037   nameEnergyReg_      = cfg.getParameter<std::string>("nameEnergyReg");
00038   nameEnergyErrorReg_ = cfg.getParameter<std::string>("nameEnergyErrorReg");
00039   debug_ = cfg.getUntrackedParameter<bool>("debug");
00040   useReducedRecHits_ = cfg.getParameter<bool>("useRecHitCollections");
00041   produceValueMaps_ = cfg.getParameter<bool>("produceValueMaps");
00042 
00043   // if Gsf Electrons; useReducedRecHits should be used)
00044   if(inputCollectionType_ == 0 && !useReducedRecHits_) {
00045     throw cms::Exception("InconsistentParameters") << " *** Inconsistent configuration : if you read GsfElectrons, you should set useRecHitCollections to true and provide the correcte values to recHitCollectionEB and recHitCollectionEE (most probably reducedEcalRecHitsEB and reducedEcalRecHitsEE )" << std::endl;
00046   }
00047 
00048   if(inputCollectionType_ == 0 && !produceValueMaps_) {
00049     std::cout << " You are running on GsfElectrons and the producer is not configured to produce ValueMaps with the results. In that case, it does not nothing !! " << std::endl;
00050   }
00051 
00052   if (inputCollectionType_ == 0) {
00053     // do nothing
00054   } else if (inputCollectionType_ == 1) {
00055     produces<ElectronCollection>();
00056   } else {
00057     throw cms::Exception("InconsistentParameters")  << " inputCollectionType should be either 0 (GsfElectrons) or 1 (pat::Electrons) " << std::endl;  
00058   }
00059 
00060 
00061   //set regression type
00062   ElectronEnergyRegressionEvaluate::ElectronEnergyRegressionType type = ElectronEnergyRegressionEvaluate::kNoTrkVar;
00063   if (energyRegressionType_ == 1) type = ElectronEnergyRegressionEvaluate::kNoTrkVar;
00064   else if (energyRegressionType_ == 2) type = ElectronEnergyRegressionEvaluate::kWithSubCluVar;
00065   else if (energyRegressionType_ == 3) type = ElectronEnergyRegressionEvaluate::kWithTrkVar;
00066 
00067   //load weights and initialize
00068   regressionEvaluator_ = new ElectronEnergyRegressionEvaluate();
00069   regressionEvaluator_->initialize(regressionInputFile_.c_str(),type);
00070 
00071   if(produceValueMaps_) {
00072     produces<edm::ValueMap<double> >(nameEnergyReg_);
00073     produces<edm::ValueMap<double> >(nameEnergyErrorReg_);
00074   }
00075 
00076 
00077   //****************************************************************************************
00078   //set up regression calculator
00079   //****************************************************************************************
00080 
00081   geomInitialized_ = false;
00082 
00083   std::cout << " Finished initialization " << std::endl;
00084 
00085 }
00086  
00087 RegressionEnergyPatElectronProducer::~RegressionEnergyPatElectronProducer()
00088 {
00089   delete regressionEvaluator_;
00090 
00091 }
00092 
00093 void RegressionEnergyPatElectronProducer::produce( edm::Event & event, const edm::EventSetup & setup )
00094 {
00095 
00096   assert(regressionEvaluator_->isInitialized());
00097   if (!geomInitialized_) {
00098     edm::ESHandle<CaloTopology> theCaloTopology;
00099     setup.get<CaloTopologyRecord>().get(theCaloTopology);
00100     ecalTopology_ = & (*theCaloTopology);
00101     
00102     edm::ESHandle<CaloGeometry> theCaloGeometry;
00103     setup.get<CaloGeometryRecord>().get(theCaloGeometry); 
00104     caloGeometry_ = & (*theCaloGeometry);
00105     geomInitialized_ = true;
00106   }
00107 
00108 
00109   //**************************************************************************
00110   //Get Number of Vertices
00111   //**************************************************************************
00112   Handle<reco::VertexCollection> hVertexProduct;
00113   event.getByLabel(verticesInputTag_,hVertexProduct);
00114   const reco::VertexCollection inVertices = *(hVertexProduct.product());  
00115 
00116   // loop through all vertices
00117   Int_t nvertices = 0;
00118   for (reco::VertexCollection::const_iterator inV = inVertices.begin(); 
00119        inV != inVertices.end(); ++inV) {
00120 
00121     //pass these vertex cuts
00122     if (inV->ndof() >= 4
00123         && inV->position().Rho() <= 2.0
00124         && fabs(inV->z()) <= 24.0
00125         ) {
00126       nvertices++;
00127     }
00128   }
00129 
00130   //**************************************************************************
00131   //Get Rho
00132   //**************************************************************************
00133   double rho = 0;
00134   Handle<double> hRhoKt6PFJets;
00135   event.getByLabel(rhoInputTag_, hRhoKt6PFJets);
00136   rho = (*hRhoKt6PFJets);
00137 
00138   //*************************************************************************
00139   // Get the RecHits
00140   //************************************************************************
00141 
00142   edm::Handle< EcalRecHitCollection > pEBRecHits;
00143   edm::Handle< EcalRecHitCollection > pEERecHits;
00144   if (useReducedRecHits_) {
00145     event.getByLabel( recHitCollectionEB_, pEBRecHits );
00146     event.getByLabel( recHitCollectionEE_, pEERecHits );
00147   }
00148 
00149   edm::Handle<GsfElectronCollection> gsfCollectionH ;
00150   edm::Handle<ElectronCollection> patCollectionH;
00151   if ( inputCollectionType_ == 0 ) {
00152     event.getByLabel ( inputElectrons_,gsfCollectionH )  ;
00153     nElectrons_ = gsfCollectionH->size();
00154   }
00155   if ( inputCollectionType_ == 1 ) {
00156     event.getByLabel ( inputElectrons_,patCollectionH )  ;
00157     nElectrons_ = patCollectionH->size();
00158   }
00159 
00160   // prepare the two even if only one is used 
00161   std::auto_ptr<ElectronCollection> patElectrons( new ElectronCollection ) ;
00162 
00163   // Fillers for ValueMaps:
00164   std::auto_ptr<edm::ValueMap<double> > regrEnergyMap(new edm::ValueMap<double>() );
00165   edm::ValueMap<double>::Filler energyFiller(*regrEnergyMap);
00166 
00167   std::auto_ptr<edm::ValueMap<double> > regrEnergyErrorMap(new edm::ValueMap<double>() );
00168   edm::ValueMap<double>::Filler energyErrorFiller(*regrEnergyErrorMap);
00169 
00170 
00171   // Reserve the vectors with outputs
00172   std::vector<double> energyValues;
00173   std::vector<double> energyErrorValues;
00174   energyValues.reserve(nElectrons_);
00175   energyErrorValues.reserve(nElectrons_);
00176  
00177 
00178   for(unsigned iele=0; iele < nElectrons_ ; ++iele) {
00179     
00180     const GsfElectron * ele = ( inputCollectionType_ == 0 ) ? &(*gsfCollectionH)[iele] : &(*patCollectionH)[iele] ;
00181     if (debug_) {
00182       std::cout << "***********************************************************************\n";
00183       std::cout << "Run Lumi Event: " << event.id().run() << " " << event.luminosityBlock() << " " << event.id().event() << "\n";
00184       std::cout << "Pat Electron : " << ele->pt() << " " << ele->eta() << " " << ele->phi() << "\n";
00185     }
00186     
00187     pat::Electron * myPatElectron = (inputCollectionType_ == 0 ) ?  0 : new pat::Electron((*patCollectionH)[iele]);
00188     // Get RecHit Collection 
00189     const EcalRecHitCollection * recHits=0;
00190     if (useReducedRecHits_) {
00191       if(ele->isEB()) {
00192         recHits = pEBRecHits.product();
00193       } else
00194         recHits = pEERecHits.product();
00195     } else {
00196       recHits = (*patCollectionH)[iele].recHits();
00197     } 
00198 
00199     SuperClusterHelper * mySCHelper = 0 ;
00200     if ( inputCollectionType_ == 0 ) { 
00201       mySCHelper = new SuperClusterHelper(&(*ele),recHits,ecalTopology_,caloGeometry_);
00202     } else if ( inputCollectionType_ == 1) {
00203       mySCHelper = new SuperClusterHelper( &(*patCollectionH)[iele], recHits,ecalTopology_,caloGeometry_);
00204     }
00205    
00206     // apply regression energy
00207     Double_t FinalMomentum = 0;
00208     Double_t FinalMomentumError = 0;
00209     Double_t RegressionMomentum = 0;
00210     Double_t RegressionMomentumError = 0;
00211 
00212     if (energyRegressionType_ == 1) {
00213 
00214         RegressionMomentum = regressionEvaluator_->regressionValueNoTrkVar( mySCHelper->rawEnergy(),
00215                                                                             mySCHelper->eta(),
00216                                                                             mySCHelper->phi(),
00217                                                                             mySCHelper->r9(),
00218                                                                             mySCHelper->etaWidth(),
00219                                                                             mySCHelper->phiWidth(),
00220                                                                             mySCHelper->clustersSize(),
00221                                                                             mySCHelper->hadronicOverEm(),
00222                                                                             rho, 
00223                                                                             nvertices, 
00224                                                                             mySCHelper->seedEta(),
00225                                                                             mySCHelper->seedPhi(),
00226                                                                             mySCHelper->seedEnergy(),
00227                                                                             mySCHelper->e3x3(),
00228                                                                             mySCHelper->e5x5(),
00229                                                                             mySCHelper->sigmaIetaIeta(),
00230                                                                             mySCHelper->spp(),
00231                                                                             mySCHelper->sep(),
00232                                                                             mySCHelper->eMax(),
00233                                                                             mySCHelper->e2nd(),
00234                                                                             mySCHelper->eTop(),
00235                                                                             mySCHelper->eBottom(),
00236                                                                             mySCHelper->eLeft(),
00237                                                                             mySCHelper->eRight(),
00238                                                                             mySCHelper->e2x5Max(),
00239                                                                             mySCHelper->e2x5Top(),
00240                                                                             mySCHelper->e2x5Bottom(),
00241                                                                             mySCHelper->e2x5Left(),
00242                                                                             mySCHelper->e2x5Right(),
00243                                                                             mySCHelper->ietaSeed(),
00244                                                                             mySCHelper->iphiSeed(),
00245                                                                             mySCHelper->etaCrySeed(),
00246                                                                             mySCHelper->phiCrySeed(),
00247                                                                             mySCHelper->preshowerEnergyOverRaw(),
00248                                                                             debug_);
00249         RegressionMomentumError = regressionEvaluator_->regressionUncertaintyNoTrkVar( 
00250                                                                                       mySCHelper->rawEnergy(),
00251                                                                                       mySCHelper->eta(),
00252                                                                                       mySCHelper->phi(),
00253                                                                                       mySCHelper->r9(),
00254                                                                                       mySCHelper->etaWidth(),
00255                                                                                       mySCHelper->phiWidth(),
00256                                                                                       mySCHelper->clustersSize(),
00257                                                                                       mySCHelper->hadronicOverEm(),
00258                                                                                       rho, 
00259                                                                                       nvertices, 
00260                                                                                       mySCHelper->seedEta(),
00261                                                                                       mySCHelper->seedPhi(),
00262                                                                                       mySCHelper->seedEnergy(),
00263                                                                                       mySCHelper->e3x3(),
00264                                                                                       mySCHelper->e5x5(),
00265                                                                                       mySCHelper->sigmaIetaIeta(),
00266                                                                                       mySCHelper->spp(),
00267                                                                                       mySCHelper->sep(),
00268                                                                                       mySCHelper->eMax(),
00269                                                                                       mySCHelper->e2nd(),
00270                                                                                       mySCHelper->eTop(),
00271                                                                                       mySCHelper->eBottom(),
00272                                                                                       mySCHelper->eLeft(),
00273                                                                                       mySCHelper->eRight(),
00274                                                                                       mySCHelper->e2x5Max(),
00275                                                                                       mySCHelper->e2x5Top(),
00276                                                                                       mySCHelper->e2x5Bottom(),
00277                                                                                       mySCHelper->e2x5Left(),
00278                                                                                       mySCHelper->e2x5Right(),
00279                                                                                       mySCHelper->ietaSeed(),
00280                                                                                       mySCHelper->iphiSeed(),
00281                                                                                       mySCHelper->etaCrySeed(),
00282                                                                                       mySCHelper->phiCrySeed(),
00283                                                                                       mySCHelper->preshowerEnergyOverRaw(),
00284                                                                                       debug_);
00285         
00286         // PAT method
00287         if(inputCollectionType_ == 1) {
00288           myPatElectron->setEcalRegressionEnergy(RegressionMomentum, RegressionMomentumError); 
00289         }
00290         energyValues.push_back(RegressionMomentum);
00291         energyErrorValues.push_back(RegressionMomentumError);   
00292 
00293 
00294     } else if (energyRegressionType_ == 2) {// ECAL regression with subcluster information
00295         RegressionMomentum = regressionEvaluator_->regressionValueWithSubClusters( 
00296                 mySCHelper->rawEnergy(),
00297                 mySCHelper->eta(),
00298                 mySCHelper->phi(),
00299                 mySCHelper->r9(),
00300                 mySCHelper->etaWidth(),
00301                 mySCHelper->phiWidth(),
00302                 mySCHelper->clustersSize(),
00303                 mySCHelper->hadronicOverEm(),
00304                 rho, 
00305                 nvertices, 
00306                 mySCHelper->seedEta(),
00307                 mySCHelper->seedPhi(),
00308                 mySCHelper->seedEnergy(),
00309                 mySCHelper->e3x3(),
00310                 mySCHelper->e5x5(),
00311                 mySCHelper->sigmaIetaIeta(),
00312                 mySCHelper->spp(),
00313                 mySCHelper->sep(),
00314                 mySCHelper->eMax(),
00315                 mySCHelper->e2nd(),
00316                 mySCHelper->eTop(),
00317                 mySCHelper->eBottom(),
00318                 mySCHelper->eLeft(),
00319                 mySCHelper->eRight(),
00320                 mySCHelper->e2x5Max(),
00321                 mySCHelper->e2x5Top(),
00322                 mySCHelper->e2x5Bottom(),
00323                 mySCHelper->e2x5Left(),
00324                 mySCHelper->e2x5Right(),
00325                 mySCHelper->ietaSeed(),
00326                 mySCHelper->iphiSeed(),
00327                 mySCHelper->etaCrySeed(),
00328                 mySCHelper->phiCrySeed(),
00329                 mySCHelper->preshowerEnergyOverRaw(),
00330                 ele->ecalDrivenSeed(),
00331                 ele->isEBEtaGap(),
00332                 ele->isEBPhiGap(),
00333                 ele->isEEDeeGap(),
00334                 mySCHelper->eSubClusters(),
00335                 mySCHelper->subClusterEnergy(1),
00336                 mySCHelper->subClusterEta(1),
00337                 mySCHelper->subClusterPhi(1),
00338                 mySCHelper->subClusterEmax(1),
00339                 mySCHelper->subClusterE3x3(1),
00340                 mySCHelper->subClusterEnergy(2),
00341                 mySCHelper->subClusterEta(2),
00342                 mySCHelper->subClusterPhi(2),
00343                 mySCHelper->subClusterEmax(2),
00344                 mySCHelper->subClusterE3x3(2),
00345                 mySCHelper->subClusterEnergy(3),
00346                 mySCHelper->subClusterEta(3),
00347                 mySCHelper->subClusterPhi(3),
00348                 mySCHelper->subClusterEmax(3),
00349                 mySCHelper->subClusterE3x3(3),
00350                         mySCHelper->nPreshowerClusters(),
00351                         mySCHelper->eESClusters(),
00352                         mySCHelper->esClusterEnergy(0),
00353                         mySCHelper->esClusterEta(0),
00354                         mySCHelper->esClusterPhi(0),
00355                         mySCHelper->esClusterEnergy(1),
00356                         mySCHelper->esClusterEta(1),
00357                         mySCHelper->esClusterPhi(1),
00358                         mySCHelper->esClusterEnergy(2),
00359                         mySCHelper->esClusterEta(2),
00360                         mySCHelper->esClusterPhi(2),
00361                 ele->isEB(),
00362                 debug_);
00363         RegressionMomentumError = regressionEvaluator_->regressionUncertaintyWithSubClusters( 
00364                 mySCHelper->rawEnergy(),
00365                 mySCHelper->eta(),
00366                 mySCHelper->phi(),
00367                 mySCHelper->r9(),
00368                 mySCHelper->etaWidth(),
00369                 mySCHelper->phiWidth(),
00370                 mySCHelper->clustersSize(),
00371                 mySCHelper->hadronicOverEm(),
00372                 rho, 
00373                 nvertices, 
00374                 mySCHelper->seedEta(),
00375                 mySCHelper->seedPhi(),
00376                 mySCHelper->seedEnergy(),
00377                 mySCHelper->e3x3(),
00378                 mySCHelper->e5x5(),
00379                 mySCHelper->sigmaIetaIeta(),
00380                 mySCHelper->spp(),
00381                 mySCHelper->sep(),
00382                 mySCHelper->eMax(),
00383                 mySCHelper->e2nd(),
00384                 mySCHelper->eTop(),
00385                 mySCHelper->eBottom(),
00386                 mySCHelper->eLeft(),
00387                 mySCHelper->eRight(),
00388                 mySCHelper->e2x5Max(),
00389                 mySCHelper->e2x5Top(),
00390                 mySCHelper->e2x5Bottom(),
00391                 mySCHelper->e2x5Left(),
00392                 mySCHelper->e2x5Right(),
00393                 mySCHelper->ietaSeed(),
00394                 mySCHelper->iphiSeed(),
00395                 mySCHelper->etaCrySeed(),
00396                 mySCHelper->phiCrySeed(),
00397                 mySCHelper->preshowerEnergyOverRaw(),
00398                 ele->ecalDrivenSeed(),
00399                 ele->isEBEtaGap(),
00400                 ele->isEBPhiGap(),
00401                 ele->isEEDeeGap(),
00402                 mySCHelper->eSubClusters(),
00403                 mySCHelper->subClusterEnergy(1),
00404                 mySCHelper->subClusterEta(1),
00405                 mySCHelper->subClusterPhi(1),
00406                 mySCHelper->subClusterEmax(1),
00407                 mySCHelper->subClusterE3x3(1),
00408                 mySCHelper->subClusterEnergy(2),
00409                 mySCHelper->subClusterEta(2),
00410                 mySCHelper->subClusterPhi(2),
00411                 mySCHelper->subClusterEmax(2),
00412                 mySCHelper->subClusterE3x3(2),
00413                 mySCHelper->subClusterEnergy(3),
00414                 mySCHelper->subClusterEta(3),
00415                 mySCHelper->subClusterPhi(3),
00416                 mySCHelper->subClusterEmax(3),
00417                 mySCHelper->subClusterE3x3(3),
00418                         mySCHelper->nPreshowerClusters(),
00419                         mySCHelper->eESClusters(),
00420                         mySCHelper->esClusterEnergy(0),
00421                         mySCHelper->esClusterEta(0),
00422                         mySCHelper->esClusterPhi(0),
00423                         mySCHelper->esClusterEnergy(1),
00424                         mySCHelper->esClusterEta(1),
00425                         mySCHelper->esClusterPhi(1),
00426                         mySCHelper->esClusterEnergy(2),
00427                         mySCHelper->esClusterEta(2),
00428                         mySCHelper->esClusterPhi(2),
00429                 ele->isEB(),
00430                 debug_);
00431 
00432         // PAT method
00433         if(inputCollectionType_ == 1) {
00434             myPatElectron->setEcalRegressionEnergy(RegressionMomentum, RegressionMomentumError); 
00435         }
00436         energyValues.push_back(RegressionMomentum);
00437         energyErrorValues.push_back(RegressionMomentumError);   
00438 
00439 
00440     }
00441       
00442       else if (energyRegressionType_ == 3) {
00443         RegressionMomentum = regressionEvaluator_->regressionValueWithTrkVar(ele->p(),
00444                                                                             mySCHelper->rawEnergy(),
00445                                                                             mySCHelper->eta(),
00446                                                                             mySCHelper->phi(),
00447                                                                             mySCHelper->etaWidth(),
00448                                                                             mySCHelper->phiWidth(),
00449                                                                             mySCHelper->clustersSize(),
00450                                                                             mySCHelper->hadronicOverEm(),
00451                                                                             mySCHelper->r9(),
00452                                                                             rho, 
00453                                                                             nvertices, 
00454                                                                             mySCHelper->seedEta(),
00455                                                                             mySCHelper->seedPhi(),
00456                                                                             mySCHelper->seedEnergy(),
00457                                                                             mySCHelper->e3x3(),
00458                                                                             mySCHelper->e5x5(),
00459                                                                             mySCHelper->sigmaIetaIeta(),
00460                                                                             mySCHelper->spp(),
00461                                                                             mySCHelper->sep(),
00462                                                                             mySCHelper->eMax(),
00463                                                                             mySCHelper->e2nd(),
00464                                                                             mySCHelper->eTop(),
00465                                                                             mySCHelper->eBottom(),
00466                                                                             mySCHelper->eLeft(),
00467                                                                             mySCHelper->eRight(),
00468                                                                             mySCHelper->e2x5Max(),
00469                                                                             mySCHelper->e2x5Top(),
00470                                                                             mySCHelper->e2x5Bottom(),
00471                                                                             mySCHelper->e2x5Left(),
00472                                                                             mySCHelper->e2x5Right(),
00473                                                                             ele->pt(),
00474                                                                             ele->trackMomentumAtVtx().R(),
00475                                                                             ele->fbrem(),
00476                                                                             ele->charge(),
00477                                                                             ele->eSuperClusterOverP(),
00478                                                                             mySCHelper->ietaSeed(),
00479                                                                             mySCHelper->iphiSeed(),
00480                                                                             mySCHelper->etaCrySeed(),
00481                                                                             mySCHelper->phiCrySeed(),
00482                                                                             mySCHelper->preshowerEnergy(),
00483                                                                             debug_);
00484         RegressionMomentumError = regressionEvaluator_->regressionUncertaintyWithTrkVar(
00485                                                                                        ele->p(),
00486                                                                                        mySCHelper->rawEnergy(),
00487                                                                                        mySCHelper->eta(),
00488                                                                                        mySCHelper->phi(),
00489                                                                                        mySCHelper->etaWidth(),
00490                                                                                        mySCHelper->phiWidth(),
00491                                                                                        mySCHelper->clustersSize(),
00492                                                                                        mySCHelper->hadronicOverEm(),
00493                                                                                        mySCHelper->r9(),
00494                                                                                        rho, 
00495                                                                                        nvertices, 
00496                                                                                        mySCHelper->seedEta(),
00497                                                                                        mySCHelper->seedPhi(),
00498                                                                                        mySCHelper->seedEnergy(),
00499                                                                                        mySCHelper->e3x3(),
00500                                                                                        mySCHelper->e5x5(),
00501                                                                                        mySCHelper->sigmaIetaIeta(),
00502                                                                                        mySCHelper->spp(),
00503                                                                                        mySCHelper->sep(),
00504                                                                                        mySCHelper->eMax(),
00505                                                                                        mySCHelper->e2nd(),
00506                                                                                        mySCHelper->eTop(),
00507                                                                                        mySCHelper->eBottom(),
00508                                                                                        mySCHelper->eLeft(),
00509                                                                                        mySCHelper->eRight(),
00510                                                                                        mySCHelper->e2x5Max(),
00511                                                                                        mySCHelper->e2x5Top(),
00512                                                                                        mySCHelper->e2x5Bottom(),
00513                                                                                        mySCHelper->e2x5Left(),
00514                                                                                        mySCHelper->e2x5Right(),
00515                                                                                        ele->pt(),
00516                                                                                        ele->trackMomentumAtVtx().R(),
00517                                                                                        ele->fbrem(),
00518                                                                                        ele->charge(),
00519                                                                                        ele->eSuperClusterOverP(),
00520                                                                                        mySCHelper->ietaSeed(),
00521                                                                                        mySCHelper->iphiSeed(),
00522                                                                                        mySCHelper->etaCrySeed(),
00523                                                                                        mySCHelper->phiCrySeed(),
00524                                                                                        mySCHelper->preshowerEnergy(),
00525                                                                                        debug_);
00526         FinalMomentum = RegressionMomentum;
00527         FinalMomentumError = RegressionMomentumError;
00528         math::XYZTLorentzVector oldMomentum = ele->p4();
00529         math::XYZTLorentzVector newMomentum = math::XYZTLorentzVector
00530           ( oldMomentum.x()*FinalMomentum/oldMomentum.t(),
00531             oldMomentum.y()*FinalMomentum/oldMomentum.t(),
00532             oldMomentum.z()*FinalMomentum/oldMomentum.t(),
00533             FinalMomentum ) ;
00534 
00535         myPatElectron->correctEcalEnergy(RegressionMomentum, RegressionMomentumError);
00536         myPatElectron->correctMomentum(newMomentum,ele->trackMomentumError(),FinalMomentumError);       
00537         
00538         energyValues.push_back(RegressionMomentum);
00539         energyErrorValues.push_back(RegressionMomentumError);   
00540       } else {
00541         cout << "Error: RegressionType = " << energyRegressionType_ << " is not supported.\n";
00542       }
00543 
00544       if(inputCollectionType_ == 1) {
00545         patElectrons->push_back(*myPatElectron);                                    
00546      }
00547       if (myPatElectron) delete myPatElectron;
00548       if (mySCHelper) delete mySCHelper;
00549   } // loop on electrons
00550 
00551   // Write the new collection in the event (AOD case)
00552   if(inputCollectionType_ == 1) {
00553     event.put(patElectrons) ;
00554   }
00555   
00556   // now AOD case: write ValueMaps
00557   if (produceValueMaps_) {
00558     
00559     if ( inputCollectionType_ ==0 ) {
00560       energyFiller.insert( gsfCollectionH, energyValues.begin(), energyValues.end() );
00561       energyErrorFiller.insert( gsfCollectionH, energyErrorValues.begin(), energyErrorValues.end() );
00562     } else if ( inputCollectionType_ ==1 ) {
00563       energyFiller.insert( patCollectionH, energyValues.begin(), energyValues.end() );
00564       energyErrorFiller.insert( patCollectionH, energyErrorValues.begin(), energyErrorValues.end() );
00565     }
00566      
00567     energyFiller.fill();
00568     energyErrorFiller.fill();
00569     event.put(regrEnergyMap,nameEnergyReg_);
00570     event.put(regrEnergyErrorMap,nameEnergyErrorReg_);
00571   }
00572   
00573 }
00574   
00575   
00576 #include "FWCore/Framework/interface/MakerMacros.h"
00577 #include "FWCore/Framework/interface/ESProducer.h"
00578 #include "FWCore/Framework/interface/ModuleFactory.h"
00579 DEFINE_FWK_MODULE(RegressionEnergyPatElectronProducer);