Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 11:58:41

0001 // -*- C++ -*-
0002 
0003 // user include files
0004 #include "FWCore/Framework/interface/MakerMacros.h"
0005 #include "Calibration/EcalCalibAlgos/interface/ElectronRecalibSuperClusterAssociator.h"
0006 
0007 #include "DataFormats/GsfTrackReco/interface/GsfTrack.h"
0008 #include <iostream>
0009 #include "DataFormats/Math/interface/LorentzVector.h"
0010 #include "DataFormats/Math/interface/deltaR.h"
0011 
0012 using namespace reco;
0013 using namespace edm;
0014 
0015 ElectronRecalibSuperClusterAssociator::ElectronRecalibSuperClusterAssociator(const edm::ParameterSet& iConfig) {
0016 #ifdef DEBUG
0017   std::cout << "ElectronRecalibSuperClusterAssociator::ElectronRecalibSuperClusterAssociator" << std::endl;
0018 #endif
0019 
0020   //register your products
0021   produces<GsfElectronCollection>();
0022   produces<GsfElectronCoreCollection>();
0023 
0024   superClusterCollectionEB_ = iConfig.getParameter<edm::InputTag>("superClusterCollectionEB");
0025   superClusterCollectionEE_ = iConfig.getParameter<edm::InputTag>("superClusterCollectionEE");
0026 
0027   outputLabel_ = iConfig.getParameter<std::string>("outputLabel");
0028   electronSrc_ = iConfig.getParameter<edm::InputTag>("electronSrc");
0029 
0030   electronToken_ = consumes<reco::GsfElectronCollection>(electronSrc_);
0031   ebScToken_ = consumes<reco::SuperClusterCollection>(superClusterCollectionEB_);
0032   eeScToken_ = consumes<reco::SuperClusterCollection>(superClusterCollectionEE_);
0033 
0034 #ifdef DEBUG
0035   std::cout << "ElectronRecalibSuperClusterAssociator::ElectronRecalibSuperClusterAssociator::end" << std::endl;
0036 #endif
0037 }
0038 
0039 ElectronRecalibSuperClusterAssociator::~ElectronRecalibSuperClusterAssociator() {}
0040 
0041 // ------------ method called to produce the data  ------------
0042 void ElectronRecalibSuperClusterAssociator::produce(edm::Event& e, const edm::EventSetup& iSetup) {
0043 #ifdef DEBUG
0044   std::cout << "GEDElectronRecalibSuperClusterAssociator::produce" << std::endl;
0045 #endif
0046 
0047   // Create the output collections
0048   auto pOutEle = std::make_unique<GsfElectronCollection>();
0049   auto pOutEleCore = std::make_unique<GsfElectronCoreCollection>();
0050 
0051   // Get SuperClusters in EB
0052   Handle<reco::SuperClusterCollection> superClusterEBHandle;
0053   e.getByToken(ebScToken_, superClusterEBHandle);
0054 
0055 #ifdef DEBUG
0056   std::cout << "EB scCollection->size()" << superClusterEBHandle->size() << std::endl;
0057 #endif
0058 
0059   // Get SuperClusters in EE
0060   Handle<reco::SuperClusterCollection> superClusterEEHandle;
0061   e.getByToken(eeScToken_, superClusterEEHandle);
0062 
0063 #ifdef DEBUG
0064   std::cout << "EE scCollection->size()" << superClusterEEHandle->size() << std::endl;
0065 #endif
0066 
0067   // Get Electrons
0068   edm::Handle<reco::GsfElectronCollection> eleHandle;
0069   e.getByToken(electronToken_, eleHandle);
0070 
0071   GsfElectronCoreRefProd rEleCore = e.getRefBeforePut<GsfElectronCoreCollection>();
0072   edm::Ref<GsfElectronCoreCollection>::key_type idxEleCore = 0;
0073 
0074   for (reco::GsfElectronCollection::const_iterator eleIt = eleHandle->begin(); eleIt != eleHandle->end(); eleIt++) {
0075     float DeltaRMineleSCbarrel(0.15);  //initial minDeltaR
0076     float DeltaRMineleSCendcap(0.15);
0077     const reco::SuperCluster* nearestSCbarrel = nullptr;
0078     const reco::SuperCluster* nearestSCendcap = nullptr;
0079     int iscRef = -1, iscRefendcap = -1;
0080     int iSC = 0;
0081 
0082     if (eleIt->trackerDrivenSeed()) {
0083       edm::LogError("trackerDriven") << "skipping trackerDriven electrons";
0084       continue;
0085     }
0086     // first loop is on EB superClusters
0087     iSC = 0;
0088     for (reco::SuperClusterCollection::const_iterator scIt = superClusterEBHandle->begin();
0089          scIt != superClusterEBHandle->end();
0090          scIt++, iSC++) {
0091       double DeltaReleSC = sqrt(reco::deltaR2(eleIt->eta(), eleIt->phi(), scIt->eta(), scIt->phi()));
0092 
0093       if (DeltaReleSC < DeltaRMineleSCbarrel)  //save the nearest SC
0094       {
0095         DeltaRMineleSCbarrel = DeltaReleSC;
0096         nearestSCbarrel = &*scIt;
0097         iscRef = iSC;
0098       }
0099 #ifdef DEBUG
0100       std::cout << "EB: " << scIt - superClusterEBHandle->begin() << " " << iSC << " " << iscRef << "\t"
0101                 << std::setprecision(4) << scIt->energy() << " " << scIt->eta() << " " << scIt->phi() << "\t--\t"
0102                 << eleIt->energy() << " " << eleIt->eta() << " " << eleIt->phi() << "\t" << DeltaRMineleSCbarrel
0103                 << std::endl;
0104 #endif
0105     }
0106 
0107     // second loop is on EE superClusters
0108     iSC = 0;
0109     for (reco::SuperClusterCollection::const_iterator scIt = superClusterEEHandle->begin();
0110          scIt != superClusterEEHandle->end();
0111          scIt++, iSC++) {
0112 #ifdef DEBUG
0113       std::cout << "EE: " << scIt - superClusterEEHandle->begin() << " " << iSC << " " << iscRef << "\t"
0114                 << std::setprecision(4) << scIt->energy() << " " << scIt->eta() << " " << scIt->phi() << "\t--\t "
0115                 << eleIt->energy() << " " << eleIt->eta() << " " << eleIt->phi() << "\t" << DeltaRMineleSCendcap
0116                 << std::endl;
0117 #endif
0118 
0119       double DeltaReleSC = sqrt(reco::deltaR2(eleIt->eta(), eleIt->phi(), scIt->eta(), scIt->phi()));
0120 
0121       if (DeltaReleSC < DeltaRMineleSCendcap) {
0122         DeltaRMineleSCendcap = DeltaReleSC;
0123         nearestSCendcap = &*scIt;
0124         iscRefendcap = iSC;
0125       }
0126     }
0127     if (eleIt->isEB() && DeltaRMineleSCbarrel > DeltaRMineleSCendcap) {
0128       edm::LogError("ElectronRecalibAssociator") << "EB electron, but nearest SC is in EE";
0129       ;
0130       continue;
0131     }
0132 
0133     if (eleIt->isEB() && nearestSCbarrel) {
0134       pOutEleCore->push_back(*eleIt->core());  // clone the old core and add to the collection of new cores
0135       reco::GsfElectronCoreRef newEleCoreRef(rEleCore,
0136                                              idxEleCore++);  // reference to the new electron core in the new collection
0137       reco::GsfElectronCore& newEleCore = pOutEleCore->back();                           // pick the clone
0138       reco::SuperClusterRef scRef(reco::SuperClusterRef(superClusterEBHandle, iscRef));  // Reference to the new SC
0139 #ifndef CMSSW_5_3_X
0140       newEleCore.setParentSuperCluster(scRef);  // mustache
0141 #endif
0142       newEleCore.setSuperCluster(scRef);  // let's check this! if it is possible to recreate the pfSC
0143 
0144       pOutEle->push_back(reco::GsfElectron(*eleIt, newEleCoreRef));
0145       reco::GsfElectron& newEle = pOutEle->back();
0146 
0147       //-- first possibility: set the new p4SC using refined SC
0148       newEle.setP4(reco::GsfElectron::P4_FROM_SUPER_CLUSTER,
0149                    eleIt->p4(reco::GsfElectron::P4_FROM_SUPER_CLUSTER),
0150                    eleIt->p4Error(reco::GsfElectron::P4_FROM_SUPER_CLUSTER),
0151                    false);  //*newEle.superCluster()->energy()/eleIt->superCluster()->energy());
0152 
0153       //-- second possibility: set the new p4SC using mustache SC
0154       //newEle.setP4(reco::GsfElectron::P4_FROM_SUPER_CLUSTER, eleIt->p4(reco::GsfElectron::P4_FROM_SUPER_CLUSTER)*newEle.parentSuperCluster()->energy()/eleIt->parentSuperCluster()->energy(), eleIt->p4Error(reco::GsfElectron::P4_FROM_SUPER_CLUSTER), false);
0155 
0156       //-- update the correctedEcalEnergy
0157       newEle.setCorrectedEcalEnergy(eleIt->ecalEnergy() *
0158                                     (scRef->energy() / eleIt->p4(reco::GsfElectron::P4_FROM_SUPER_CLUSTER).energy()));
0159       newEle.setCorrectedEcalEnergyError(eleIt->ecalEnergyError() * (scRef->energy() / eleIt->ecalEnergy()));
0160 
0161     } else if (!(eleIt->isEB()) && nearestSCendcap) {
0162       pOutEleCore->push_back(*eleIt->core());  // clone the old core and add to the collection of new cores
0163       reco::GsfElectronCoreRef newEleCoreRef(rEleCore,
0164                                              idxEleCore++);  // reference to the new electron core in the new collection
0165       reco::GsfElectronCore& newEleCore = pOutEleCore->back();  // pick the clone
0166       reco::SuperClusterRef scRef(
0167           reco::SuperClusterRef(superClusterEEHandle, iscRefendcap));  // Reference to the new SC
0168 #ifndef CMSSW_5_3_X
0169       newEleCore.setParentSuperCluster(scRef);  // mustache
0170 #endif
0171       newEleCore.setSuperCluster(scRef);  // let's check this! if it is possible to recreate the pfSC
0172 
0173       pOutEle->push_back(reco::GsfElectron(*eleIt, newEleCoreRef));
0174       reco::GsfElectron& newEle = pOutEle->back();
0175 
0176       //-- first possibility: set the new p4SC using refined SC
0177       newEle.setP4(reco::GsfElectron::P4_FROM_SUPER_CLUSTER,
0178                    eleIt->p4(reco::GsfElectron::P4_FROM_SUPER_CLUSTER),
0179                    eleIt->p4Error(reco::GsfElectron::P4_FROM_SUPER_CLUSTER),
0180                    false);  //*newEle.superCluster()->energy()/eleIt->superCluster()->energy());
0181 
0182       //-- second possibility: set the new p4SC using mustache SC
0183       //newEle.setP4(reco::GsfElectron::P4_FROM_SUPER_CLUSTER, eleIt->p4(reco::GsfElectron::P4_FROM_SUPER_CLUSTER)*newEle.parentSuperCluster()->energy()/eleIt->parentSuperCluster()->energy(), eleIt->p4Error(reco::GsfElectron::P4_FROM_SUPER_CLUSTER), false);
0184 
0185       //-- update the correctedEcalEnergy
0186       newEle.setCorrectedEcalEnergy(eleIt->ecalEnergy() *
0187                                     (scRef->energy() / eleIt->p4(reco::GsfElectron::P4_FROM_SUPER_CLUSTER).energy()));
0188       newEle.setCorrectedEcalEnergyError(eleIt->ecalEnergyError() * (scRef->energy() / eleIt->ecalEnergy()));
0189     } else {
0190       edm::LogError("Failed SC association") << "No SC to be associated to the electron";
0191     }
0192   }
0193 
0194 #ifdef DEBUG
0195   std::cout << "Filled new electrons  " << pOutEle->size() << std::endl;
0196   std::cout << "Filled new electronsCore  " << pOutEleCore->size() << std::endl;
0197 #endif
0198 
0199   // put result into the Event
0200 
0201   e.put(std::move(pOutEle));
0202   e.put(std::move(pOutEleCore));
0203 }
0204 
0205 DEFINE_FWK_MODULE(ElectronRecalibSuperClusterAssociator);