Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-06-22 02:24:06

0001 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0002 #include "FWCore/Framework/interface/Event.h"
0003 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0004 #include "CommonTools/ParticleFlow/interface/PFClusterWidthAlgo.h"
0005 #include "RecoEcal/EgammaCoreTools/interface/Mustache.h"
0006 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
0007 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidateElectronExtra.h"
0008 #include "DataFormats/EgammaReco/interface/PreshowerCluster.h"
0009 #include "DataFormats/EgammaReco/interface/SuperCluster.h"
0010 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
0011 #include "DataFormats/EgammaCandidates/interface/GsfElectronCore.h"
0012 #include "DataFormats/ParticleFlowReco/interface/PFBlockElement.h"
0013 #include "DataFormats/ParticleFlowReco/interface/PFBlockElementGsfTrack.h"
0014 #include "FWCore/Framework/interface/stream/EDProducer.h"
0015 #include "FWCore/Framework/interface/MakerMacros.h"
0016 #include "DataFormats/Common/interface/ValueMap.h"
0017 #include "DataFormats/ParticleFlowReco/interface/PFCluster.h"
0018 #include "DataFormats/ParticleFlowReco/interface/PFBlock.h"
0019 #include "DataFormats/EgammaReco/interface/BasicCluster.h"
0020 #include "CondFormats/EcalObjects/interface/EcalMustacheSCParameters.h"
0021 #include "CondFormats/DataRecord/interface/EcalMustacheSCParametersRcd.h"
0022 
0023 class DetId;
0024 namespace edm {
0025   class EventSetup;
0026 }  // namespace edm
0027 
0028 class PFElectronTranslator : public edm::stream::EDProducer<> {
0029 public:
0030   explicit PFElectronTranslator(const edm::ParameterSet&);
0031   ~PFElectronTranslator() override;
0032 
0033   void produce(edm::Event&, const edm::EventSetup&) override;
0034   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0035 
0036   typedef std::vector<edm::Handle<edm::ValueMap<double>>> IsolationValueMaps;
0037 
0038 private:
0039   // to retrieve the collection from the event
0040   bool fetchCandidateCollection(edm::Handle<reco::PFCandidateCollection>& c,
0041                                 const edm::EDGetTokenT<reco::PFCandidateCollection>& token,
0042                                 const edm::Event& iEvent) const;
0043   // to retrieve the collection from the event
0044   void fetchGsfCollection(edm::Handle<reco::GsfTrackCollection>& c,
0045                           const edm::EDGetTokenT<reco::GsfTrackCollection>& token,
0046                           const edm::Event& iEvent) const;
0047 
0048   // makes a basic cluster from PFBlockElement and add it to the collection ; the corrected energy is taken
0049   // from the PFCandidate
0050   void createBasicCluster(const reco::PFBlockElement&,
0051                           reco::BasicClusterCollection& basicClusters,
0052                           std::vector<const reco::PFCluster*>&,
0053                           const reco::PFCandidate& coCandidate) const;
0054   // makes a preshower cluster from of PFBlockElement and add it to the collection
0055   void createPreshowerCluster(const reco::PFBlockElement& PFBE,
0056                               reco::PreshowerClusterCollection& preshowerClusters,
0057                               unsigned plane) const;
0058 
0059   // make a super cluster from its ingredients and add it to the collection
0060   void createSuperClusters(const reco::PFCandidateCollection&, reco::SuperClusterCollection& superClusters) const;
0061 
0062   // make GsfElectronCores from ingredients
0063   void createGsfElectronCores(reco::GsfElectronCoreCollection&) const;
0064 
0065   // create the basic cluster Ptr
0066   void createBasicClusterPtrs(const edm::OrphanHandle<reco::BasicClusterCollection>& basicClustersHandle);
0067 
0068   // create the preshower cluster Refs
0069   void createPreshowerClusterPtrs(const edm::OrphanHandle<reco::PreshowerClusterCollection>& preshowerClustersHandle);
0070 
0071   // create the super cluster Refs
0072   void createSuperClusterGsfMapRefs(const edm::OrphanHandle<reco::SuperClusterCollection>& superClustersHandle);
0073 
0074   // create the GsfElectronCore Refs
0075   void createGsfElectronCoreRefs(const edm::OrphanHandle<reco::GsfElectronCoreCollection>& gsfElectronCoreHandle);
0076 
0077   // create the GsfElectrons
0078   void createGsfElectrons(const reco::PFCandidateCollection&,
0079                           const IsolationValueMaps& isolationValues,
0080                           reco::GsfElectronCollection&);
0081 
0082   // The following methods are used to fill the value maps
0083   void fillMVAValueMap(edm::Event& iEvent, edm::ValueMap<float>::Filler& filler);
0084   void fillValueMap(edm::Event& iEvent, edm::ValueMap<float>::Filler& filler) const;
0085   void fillSCRefValueMap(edm::Event& iEvent, edm::ValueMap<reco::SuperClusterRef>::Filler& filler) const;
0086   void getAmbiguousGsfTracks(const reco::PFBlockElement& PFBE, std::vector<reco::GsfTrackRef>&) const;
0087 
0088   const reco::PFCandidate& correspondingDaughterCandidate(const reco::PFCandidate& cand,
0089                                                           const reco::PFBlockElement& pfbe) const;
0090 
0091 private:
0092   edm::EDGetTokenT<reco::PFCandidateCollection> inputTokenPFCandidates_;
0093   edm::EDGetTokenT<reco::PFCandidateCollection> inputTokenPFCandidateElectrons_;
0094   edm::EDGetTokenT<reco::GsfTrackCollection> inputTokenGSFTracks_;
0095   std::vector<edm::EDGetTokenT<edm::ValueMap<double>>> inputTokenIsoVals_;
0096   std::string PFBasicClusterCollection_;
0097   std::string PFPreshowerClusterCollection_;
0098   std::string PFSuperClusterCollection_;
0099   std::string PFMVAValueMap_;
0100   std::string PFSCValueMap_;
0101   std::string GsfElectronCoreCollection_;
0102   std::string GsfElectronCollection_;
0103   double MVACut_;
0104   bool checkStatusFlag_;
0105 
0106   // The following vectors correspond to a GSF track, but the order is not
0107   // the order of the tracks in the GSF track collection
0108   std::vector<reco::GsfTrackRef> GsfTrackRef_;
0109   // the list of candidatePtr
0110   std::vector<reco::CandidatePtr> CandidatePtr_;
0111   //the list of KfTrackRef
0112   std::vector<reco::TrackRef> kfTrackRef_;
0113   // the list of ambiguous tracks
0114   std::vector<std::vector<reco::GsfTrackRef>> ambiguousGsfTracks_;
0115   // the collection of basic clusters associated to a GSF track
0116   std::vector<reco::BasicClusterCollection> basicClusters_;
0117   // the correcsponding PFCluster ref
0118   std::vector<std::vector<const reco::PFCluster*>> pfClusters_;
0119   // the collection of preshower clusters associated to a GSF track
0120   std::vector<reco::PreshowerClusterCollection> preshowerClusters_;
0121   // the super cluster collection (actually only one) associated to a GSF track
0122   std::vector<reco::SuperClusterCollection> superClusters_;
0123   // the references to the basic clusters associated to a GSF track
0124   std::vector<reco::CaloClusterPtrVector> basicClusterPtr_;
0125   // the references to the basic clusters associated to a GSF track
0126   std::vector<reco::CaloClusterPtrVector> preshowerClusterPtr_;
0127   // the references to the GsfElectonCore associated to a GSF track
0128   std::vector<reco::GsfElectronCoreRef> gsfElectronCoreRefs_;
0129   // keep track of the index of the PF Candidate
0130   std::vector<int> gsfPFCandidateIndex_;
0131   // maps to ease the creation of the Value Maps
0132   std::map<reco::GsfTrackRef, reco::SuperClusterRef> scMap_;
0133   std::map<reco::GsfTrackRef, float> gsfMvaMap_;
0134 
0135   // Mustache SC parameters
0136   edm::ESGetToken<EcalMustacheSCParameters, EcalMustacheSCParametersRcd> ecalMustacheSCParametersToken_;
0137   const EcalMustacheSCParameters* mustacheSCParams_;
0138 
0139   bool emptyIsOk_;
0140 };
0141 
0142 DEFINE_FWK_MODULE(PFElectronTranslator);
0143 
0144 void PFElectronTranslator::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0145   edm::ParameterSetDescription desc;
0146   desc.add<edm::InputTag>("PFCandidate", {"pfSelectedElectrons"});
0147   desc.add<edm::InputTag>("PFCandidateElectron", {"particleFlowTmp:electrons"});
0148   desc.add<edm::InputTag>("GSFTracks", {"electronGsfTracks"});
0149   desc.add<std::string>("PFBasicClusters", "pf");
0150   desc.add<std::string>("PFPreshowerClusters", "pf");
0151   desc.add<std::string>("PFSuperClusters", "pf");
0152   desc.add<std::string>("ElectronMVA", "pf");
0153   desc.add<std::string>("ElectronSC", "pf");
0154   desc.add<std::string>("PFGsfElectron", "pf");
0155   desc.add<std::string>("PFGsfElectronCore", "pf");
0156   desc.add<edm::ParameterSetDescription>("MVACutBlock", {});
0157   desc.add<bool>("CheckStatusFlag", true);
0158   desc.add<bool>("useIsolationValues", false);
0159   {
0160     edm::ParameterSetDescription pset;
0161     pset.add<edm::InputTag>("pfSumChargedHadronPt", {"elPFIsoValueCharged04PFId"});
0162     pset.add<edm::InputTag>("pfSumPhotonEt", {"elPFIsoValueGamma04PFId"});
0163     pset.add<edm::InputTag>("pfSumNeutralHadronEt", {"elPFIsoValueNeutral04PFId"});
0164     pset.add<edm::InputTag>("pfSumPUPt", {"elPFIsoValuePU04PFId"});
0165     desc.add<edm::ParameterSetDescription>("isolationValues", pset);
0166   }
0167   desc.add<bool>("emptyIsOk", false);
0168   descriptions.addWithDefaultLabel(desc);
0169 }
0170 
0171 PFElectronTranslator::PFElectronTranslator(const edm::ParameterSet& iConfig) {
0172   inputTokenPFCandidates_ = consumes<reco::PFCandidateCollection>(iConfig.getParameter<edm::InputTag>("PFCandidate"));
0173   inputTokenPFCandidateElectrons_ =
0174       consumes<reco::PFCandidateCollection>(iConfig.getParameter<edm::InputTag>("PFCandidateElectron"));
0175   inputTokenGSFTracks_ = consumes<reco::GsfTrackCollection>(iConfig.getParameter<edm::InputTag>("GSFTracks"));
0176 
0177   bool useIsolationValues = iConfig.getParameter<bool>("useIsolationValues");
0178   if (useIsolationValues) {
0179     const auto& isoVals = iConfig.getParameter<edm::ParameterSet>("isolationValues");
0180     if (isoVals.empty())
0181       throw cms::Exception("PFElectronTranslator|InternalError") << "Missing ParameterSet isolationValues";
0182     else {
0183       inputTokenIsoVals_.push_back(
0184           consumes<edm::ValueMap<double>>(isoVals.getParameter<edm::InputTag>("pfSumChargedHadronPt")));
0185       inputTokenIsoVals_.push_back(
0186           consumes<edm::ValueMap<double>>(isoVals.getParameter<edm::InputTag>("pfSumPhotonEt")));
0187       inputTokenIsoVals_.push_back(
0188           consumes<edm::ValueMap<double>>(isoVals.getParameter<edm::InputTag>("pfSumNeutralHadronEt")));
0189       inputTokenIsoVals_.push_back(consumes<edm::ValueMap<double>>(isoVals.getParameter<edm::InputTag>("pfSumPUPt")));
0190     }
0191   }
0192 
0193   PFBasicClusterCollection_ = iConfig.getParameter<std::string>("PFBasicClusters");
0194   PFPreshowerClusterCollection_ = iConfig.getParameter<std::string>("PFPreshowerClusters");
0195   PFSuperClusterCollection_ = iConfig.getParameter<std::string>("PFSuperClusters");
0196   GsfElectronCoreCollection_ = iConfig.getParameter<std::string>("PFGsfElectronCore");
0197   GsfElectronCollection_ = iConfig.getParameter<std::string>("PFGsfElectron");
0198 
0199   PFMVAValueMap_ = iConfig.getParameter<std::string>("ElectronMVA");
0200   PFSCValueMap_ = iConfig.getParameter<std::string>("ElectronSC");
0201   MVACut_ = (iConfig.getParameter<edm::ParameterSet>("MVACutBlock")).getParameter<double>("MVACut");
0202   checkStatusFlag_ = iConfig.getParameter<bool>("CheckStatusFlag");
0203   emptyIsOk_ = iConfig.getParameter<bool>("emptyIsOk");
0204 
0205   ecalMustacheSCParametersToken_ = esConsumes<EcalMustacheSCParameters, EcalMustacheSCParametersRcd>();
0206 
0207   produces<reco::BasicClusterCollection>(PFBasicClusterCollection_);
0208   produces<reco::PreshowerClusterCollection>(PFPreshowerClusterCollection_);
0209   produces<reco::SuperClusterCollection>(PFSuperClusterCollection_);
0210   produces<reco::GsfElectronCoreCollection>(GsfElectronCoreCollection_);
0211   produces<reco::GsfElectronCollection>(GsfElectronCollection_);
0212   produces<edm::ValueMap<float>>(PFMVAValueMap_);
0213   produces<edm::ValueMap<reco::SuperClusterRef>>(PFSCValueMap_);
0214 }
0215 
0216 PFElectronTranslator::~PFElectronTranslator() {}
0217 
0218 void PFElectronTranslator::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0219   mustacheSCParams_ = &iSetup.getData(ecalMustacheSCParametersToken_);
0220 
0221   auto gsfElectronCores_p = std::make_unique<reco::GsfElectronCoreCollection>();
0222 
0223   auto gsfElectrons_p = std::make_unique<reco::GsfElectronCollection>();
0224 
0225   auto superClusters_p = std::make_unique<reco::SuperClusterCollection>();
0226 
0227   auto basicClusters_p = std::make_unique<reco::BasicClusterCollection>();
0228 
0229   auto psClusters_p = std::make_unique<reco::PreshowerClusterCollection>();
0230 
0231   auto mvaMap_p = std::make_unique<edm::ValueMap<float>>();
0232   edm::ValueMap<float>::Filler mvaFiller(*mvaMap_p);
0233 
0234   auto scMap_p = std::make_unique<edm::ValueMap<reco::SuperClusterRef>>();
0235   edm::ValueMap<reco::SuperClusterRef>::Filler scRefFiller(*scMap_p);
0236 
0237   edm::Handle<reco::PFCandidateCollection> pfCandidates;
0238   bool status = fetchCandidateCollection(pfCandidates, inputTokenPFCandidates_, iEvent);
0239 
0240   IsolationValueMaps isolationValues(inputTokenIsoVals_.size());
0241   for (size_t j = 0; j < inputTokenIsoVals_.size(); ++j) {
0242     isolationValues[j] = iEvent.getHandle(inputTokenIsoVals_[j]);
0243   }
0244 
0245   // clear the vectors
0246   GsfTrackRef_.clear();
0247   CandidatePtr_.clear();
0248   ambiguousGsfTracks_.clear();
0249   kfTrackRef_.clear();
0250   basicClusters_.clear();
0251   pfClusters_.clear();
0252   preshowerClusters_.clear();
0253   superClusters_.clear();
0254   basicClusterPtr_.clear();
0255   preshowerClusterPtr_.clear();
0256   gsfPFCandidateIndex_.clear();
0257   gsfElectronCoreRefs_.clear();
0258   scMap_.clear();
0259 
0260   // loop on the candidates
0261   //CC@@
0262   // we need first to create AND put the SuperCluster,
0263   // basic clusters and presh clusters collection
0264   // in order to get a working Handle
0265   unsigned ncand = (status) ? pfCandidates->size() : 0;
0266   unsigned iGSF = 0;
0267   for (unsigned i = 0; i < ncand; ++i) {
0268     const reco::PFCandidate& cand = (*pfCandidates)[i];
0269     if (cand.particleId() != reco::PFCandidate::e)
0270       continue;
0271     if (cand.gsfTrackRef().isNull())
0272       continue;
0273     // Note that -1 will still cut some total garbage candidates
0274     // Fill the MVA map
0275     if (cand.mva_e_pi() < MVACut_)
0276       continue;
0277 
0278     // Check the status flag
0279     if (checkStatusFlag_ && !cand.electronExtraRef()->electronStatus(reco::PFCandidateElectronExtra::Selected)) {
0280       continue;
0281     }
0282 
0283     GsfTrackRef_.push_back(cand.gsfTrackRef());
0284     kfTrackRef_.push_back(cand.trackRef());
0285     gsfPFCandidateIndex_.push_back(i);
0286 
0287     reco::PFCandidatePtr ptrToPFElectron(pfCandidates, i);
0288     //CandidatePtr_.push_back(ptrToPFElectron->sourceCandidatePtr(0));
0289     CandidatePtr_.push_back(ptrToPFElectron);
0290 
0291     basicClusters_.push_back(reco::BasicClusterCollection());
0292     pfClusters_.push_back(std::vector<const reco::PFCluster*>());
0293     preshowerClusters_.push_back(reco::PreshowerClusterCollection());
0294     ambiguousGsfTracks_.push_back(std::vector<reco::GsfTrackRef>());
0295 
0296     for (unsigned iele = 0; iele < cand.elementsInBlocks().size(); ++iele) {
0297       // first get the block
0298       reco::PFBlockRef blockRef = cand.elementsInBlocks()[iele].first;
0299       //
0300       unsigned elementIndex = cand.elementsInBlocks()[iele].second;
0301       // check it actually exists
0302       if (blockRef.isNull())
0303         continue;
0304 
0305       // then get the elements of the block
0306       const edm::OwnVector<reco::PFBlockElement>& elements = (*blockRef).elements();
0307 
0308       const reco::PFBlockElement& pfbe(elements[elementIndex]);
0309       // The first ECAL element should be the cluster associated to the GSF; defined as the seed
0310       if (pfbe.type() == reco::PFBlockElement::ECAL) {
0311         //    const reco::PFCandidate * coCandidate = &cand;
0312         // the Brem photons are saved as daughter PFCandidate; this
0313         // is convenient to access the corrected energy
0314         //    std::cout << " Found candidate "  << correspondingDaughterCandidate(coCandidate,pfbe) << " " << coCandidate << std::endl;
0315         createBasicCluster(pfbe, basicClusters_[iGSF], pfClusters_[iGSF], correspondingDaughterCandidate(cand, pfbe));
0316       }
0317       if (pfbe.type() == reco::PFBlockElement::PS1) {
0318         createPreshowerCluster(pfbe, preshowerClusters_[iGSF], 1);
0319       }
0320       if (pfbe.type() == reco::PFBlockElement::PS2) {
0321         createPreshowerCluster(pfbe, preshowerClusters_[iGSF], 2);
0322       }
0323       if (pfbe.type() == reco::PFBlockElement::GSF) {
0324         getAmbiguousGsfTracks(pfbe, ambiguousGsfTracks_[iGSF]);
0325       }
0326 
0327     }  // loop on the elements
0328 
0329     // save the basic clusters
0330     basicClusters_p->insert(basicClusters_p->end(), basicClusters_[iGSF].begin(), basicClusters_[iGSF].end());
0331     // save the preshower clusters
0332     psClusters_p->insert(psClusters_p->end(), preshowerClusters_[iGSF].begin(), preshowerClusters_[iGSF].end());
0333 
0334     ++iGSF;
0335   }  // loop on PFCandidates
0336 
0337   //Save the basic clusters and get an handle as to be able to create valid Refs (thanks to Claude)
0338   //  std::cout << " Number of basic clusters " << basicClusters_p->size() << std::endl;
0339   const edm::OrphanHandle<reco::BasicClusterCollection> bcRefProd =
0340       iEvent.put(std::move(basicClusters_p), PFBasicClusterCollection_);
0341 
0342   //preshower clusters
0343   const edm::OrphanHandle<reco::PreshowerClusterCollection> psRefProd =
0344       iEvent.put(std::move(psClusters_p), PFPreshowerClusterCollection_);
0345 
0346   // now that the Basic clusters are in the event, the Ref can be created
0347   createBasicClusterPtrs(bcRefProd);
0348   // now that the preshower clusters are in the event, the Ref can be created
0349   createPreshowerClusterPtrs(psRefProd);
0350 
0351   // and now the Super cluster can be created with valid references
0352   if (status)
0353     createSuperClusters(*pfCandidates, *superClusters_p);
0354 
0355   // Let's put the super clusters in the event
0356   const edm::OrphanHandle<reco::SuperClusterCollection> scRefProd =
0357       iEvent.put(std::move(superClusters_p), PFSuperClusterCollection_);
0358   // create the super cluster Ref
0359   createSuperClusterGsfMapRefs(scRefProd);
0360 
0361   // Now create the GsfElectronCoers
0362   createGsfElectronCores(*gsfElectronCores_p);
0363   // Put them in the as to get to be able to build a Ref
0364   const edm::OrphanHandle<reco::GsfElectronCoreCollection> gsfElectronCoreRefProd =
0365       iEvent.put(std::move(gsfElectronCores_p), GsfElectronCoreCollection_);
0366 
0367   // now create the Refs
0368   createGsfElectronCoreRefs(gsfElectronCoreRefProd);
0369 
0370   // now make the GsfElectron
0371   createGsfElectrons(*pfCandidates, isolationValues, *gsfElectrons_p);
0372   iEvent.put(std::move(gsfElectrons_p), GsfElectronCollection_);
0373 
0374   fillMVAValueMap(iEvent, mvaFiller);
0375   mvaFiller.fill();
0376 
0377   fillSCRefValueMap(iEvent, scRefFiller);
0378   scRefFiller.fill();
0379 
0380   // MVA map
0381   iEvent.put(std::move(mvaMap_p), PFMVAValueMap_);
0382   // Gsf-SC map
0383   iEvent.put(std::move(scMap_p), PFSCValueMap_);
0384 }
0385 
0386 bool PFElectronTranslator::fetchCandidateCollection(edm::Handle<reco::PFCandidateCollection>& c,
0387                                                     const edm::EDGetTokenT<reco::PFCandidateCollection>& token,
0388                                                     const edm::Event& iEvent) const {
0389   c = iEvent.getHandle(token);
0390 
0391   if (!c.isValid() && !emptyIsOk_) {
0392     std::ostringstream err;
0393     err << " cannot get PFCandidates " << std::endl;
0394     edm::LogError("PFElectronTranslator") << err.str();
0395   }
0396   return c.isValid();
0397 }
0398 
0399 void PFElectronTranslator::fetchGsfCollection(edm::Handle<reco::GsfTrackCollection>& c,
0400                                               const edm::EDGetTokenT<reco::GsfTrackCollection>& token,
0401                                               const edm::Event& iEvent) const {
0402   c = iEvent.getHandle(token);
0403 
0404   if (!c.isValid()) {
0405     std::ostringstream err;
0406     err << " cannot get GSFTracks " << std::endl;
0407     edm::LogError("PFElectronTranslator") << err.str();
0408     throw cms::Exception("MissingProduct", err.str());
0409   }
0410 }
0411 
0412 // The basic cluster is a copy of the PFCluster -> the energy is not corrected
0413 // It should be possible to get the corrected energy (including the associated PS energy)
0414 // from the PFCandidate daugthers ; Needs some work
0415 void PFElectronTranslator::createBasicCluster(const reco::PFBlockElement& PFBE,
0416                                               reco::BasicClusterCollection& basicClusters,
0417                                               std::vector<const reco::PFCluster*>& pfClusters,
0418                                               const reco::PFCandidate& coCandidate) const {
0419   const reco::PFClusterRef& myPFClusterRef = PFBE.clusterRef();
0420   if (myPFClusterRef.isNull())
0421     return;
0422 
0423   const reco::PFCluster& myPFCluster(*myPFClusterRef);
0424   pfClusters.push_back(&myPFCluster);
0425   //  std::cout << " Creating BC " << myPFCluster.energy() << " " << coCandidate.ecalEnergy() <<" "<<  coCandidate.rawEcalEnergy() <<std::endl;
0426   //  std::cout << " # hits " << myPFCluster.hitsAndFractions().size() << std::endl;
0427 
0428   //  basicClusters.push_back(reco::CaloCluster(myPFCluster.energy(),
0429   basicClusters.push_back(reco::CaloCluster(
0430       //        myPFCluster.energy(),
0431       coCandidate.rawEcalEnergy(),
0432       myPFCluster.position(),
0433       myPFCluster.caloID(),
0434       myPFCluster.hitsAndFractions(),
0435       myPFCluster.algo(),
0436       myPFCluster.seed()));
0437 }
0438 
0439 void PFElectronTranslator::createPreshowerCluster(const reco::PFBlockElement& PFBE,
0440                                                   reco::PreshowerClusterCollection& preshowerClusters,
0441                                                   unsigned plane) const {
0442   const reco::PFClusterRef& myPFClusterRef = PFBE.clusterRef();
0443   preshowerClusters.push_back(reco::PreshowerCluster(
0444       myPFClusterRef->energy(), myPFClusterRef->position(), myPFClusterRef->hitsAndFractions(), plane));
0445 }
0446 
0447 void PFElectronTranslator::createBasicClusterPtrs(
0448     const edm::OrphanHandle<reco::BasicClusterCollection>& basicClustersHandle) {
0449   unsigned size = GsfTrackRef_.size();
0450   unsigned basicClusterCounter = 0;
0451   basicClusterPtr_.resize(size);
0452 
0453   for (unsigned iGSF = 0; iGSF < size; ++iGSF)  // loop on tracks
0454   {
0455     unsigned nbc = basicClusters_[iGSF].size();
0456     for (unsigned ibc = 0; ibc < nbc; ++ibc)  // loop on basic clusters
0457     {
0458       //      std::cout <<  "Track "<< iGSF << " ref " << basicClusterCounter << std::endl;
0459       reco::CaloClusterPtr bcPtr(basicClustersHandle, basicClusterCounter);
0460       basicClusterPtr_[iGSF].push_back(bcPtr);
0461       ++basicClusterCounter;
0462     }
0463   }
0464 }
0465 
0466 void PFElectronTranslator::createPreshowerClusterPtrs(
0467     const edm::OrphanHandle<reco::PreshowerClusterCollection>& preshowerClustersHandle) {
0468   unsigned size = GsfTrackRef_.size();
0469   unsigned psClusterCounter = 0;
0470   preshowerClusterPtr_.resize(size);
0471 
0472   for (unsigned iGSF = 0; iGSF < size; ++iGSF)  // loop on tracks
0473   {
0474     unsigned nbc = preshowerClusters_[iGSF].size();
0475     for (unsigned ibc = 0; ibc < nbc; ++ibc)  // loop on basic clusters
0476     {
0477       //      std::cout <<  "Track "<< iGSF << " ref " << basicClusterCounter << std::endl;
0478       reco::CaloClusterPtr psPtr(preshowerClustersHandle, psClusterCounter);
0479       preshowerClusterPtr_[iGSF].push_back(psPtr);
0480       ++psClusterCounter;
0481     }
0482   }
0483 }
0484 
0485 void PFElectronTranslator::createSuperClusterGsfMapRefs(
0486     const edm::OrphanHandle<reco::SuperClusterCollection>& superClustersHandle) {
0487   unsigned size = GsfTrackRef_.size();
0488 
0489   for (unsigned iGSF = 0; iGSF < size; ++iGSF)  // loop on tracks
0490   {
0491     edm::Ref<reco::SuperClusterCollection> scRef(superClustersHandle, iGSF);
0492     scMap_[GsfTrackRef_[iGSF]] = scRef;
0493   }
0494 }
0495 
0496 void PFElectronTranslator::fillMVAValueMap(edm::Event& iEvent, edm::ValueMap<float>::Filler& filler) {
0497   gsfMvaMap_.clear();
0498   edm::Handle<reco::PFCandidateCollection> pfCandidates;
0499   bool status = fetchCandidateCollection(pfCandidates, inputTokenPFCandidateElectrons_, iEvent);
0500 
0501   unsigned ncand = (status) ? pfCandidates->size() : 0;
0502   for (unsigned i = 0; i < ncand; ++i) {
0503     const reco::PFCandidate& cand = (*pfCandidates)[i];
0504     if (cand.particleId() != reco::PFCandidate::e)
0505       continue;
0506     if (cand.gsfTrackRef().isNull())
0507       continue;
0508     // Fill the MVA map
0509     gsfMvaMap_[cand.gsfTrackRef()] = cand.mva_e_pi();
0510   }
0511 
0512   edm::Handle<reco::GsfTrackCollection> gsfTracks;
0513   fetchGsfCollection(gsfTracks, inputTokenGSFTracks_, iEvent);
0514   unsigned ngsf = gsfTracks->size();
0515   std::vector<float> values;
0516   for (unsigned igsf = 0; igsf < ngsf; ++igsf) {
0517     reco::GsfTrackRef theTrackRef(gsfTracks, igsf);
0518     std::map<reco::GsfTrackRef, float>::const_iterator itcheck = gsfMvaMap_.find(theTrackRef);
0519     if (itcheck == gsfMvaMap_.end()) {
0520       //      edm::LogWarning("PFElectronTranslator") << "MVA Map, missing GSF track ref " << std::endl;
0521       values.push_back(-99.);
0522       //      std::cout << " Push_back -99. " << std::endl;
0523     } else {
0524       //      std::cout <<  " Value " << itcheck->second << std::endl;
0525       values.push_back(itcheck->second);
0526     }
0527   }
0528   filler.insert(gsfTracks, values.begin(), values.end());
0529 }
0530 
0531 void PFElectronTranslator::fillSCRefValueMap(edm::Event& iEvent,
0532                                              edm::ValueMap<reco::SuperClusterRef>::Filler& filler) const {
0533   edm::Handle<reco::GsfTrackCollection> gsfTracks;
0534   fetchGsfCollection(gsfTracks, inputTokenGSFTracks_, iEvent);
0535   unsigned ngsf = gsfTracks->size();
0536   std::vector<reco::SuperClusterRef> values;
0537   for (unsigned igsf = 0; igsf < ngsf; ++igsf) {
0538     reco::GsfTrackRef theTrackRef(gsfTracks, igsf);
0539     std::map<reco::GsfTrackRef, reco::SuperClusterRef>::const_iterator itcheck = scMap_.find(theTrackRef);
0540     if (itcheck == scMap_.end()) {
0541       //      edm::LogWarning("PFElectronTranslator") << "SCRef Map, missing GSF track ref" << std::endl;
0542       values.push_back(reco::SuperClusterRef());
0543     } else {
0544       values.push_back(itcheck->second);
0545     }
0546   }
0547   filler.insert(gsfTracks, values.begin(), values.end());
0548 }
0549 
0550 void PFElectronTranslator::createSuperClusters(const reco::PFCandidateCollection& pfCand,
0551                                                reco::SuperClusterCollection& superClusters) const {
0552   unsigned nGSF = GsfTrackRef_.size();
0553   for (unsigned iGSF = 0; iGSF < nGSF; ++iGSF) {
0554     // Computes energy position a la e/gamma
0555     double sclusterE = 0;
0556     double posX = 0.;
0557     double posY = 0.;
0558     double posZ = 0.;
0559 
0560     unsigned nbasics = basicClusters_[iGSF].size();
0561     for (unsigned ibc = 0; ibc < nbasics; ++ibc) {
0562       double e = basicClusters_[iGSF][ibc].energy();
0563       sclusterE += e;
0564       posX += e * basicClusters_[iGSF][ibc].position().X();
0565       posY += e * basicClusters_[iGSF][ibc].position().Y();
0566       posZ += e * basicClusters_[iGSF][ibc].position().Z();
0567     }
0568     posX /= sclusterE;
0569     posY /= sclusterE;
0570     posZ /= sclusterE;
0571 
0572     if (pfCand[gsfPFCandidateIndex_[iGSF]].gsfTrackRef() != GsfTrackRef_[iGSF]) {
0573       edm::LogError("PFElectronTranslator") << " Major problem in PFElectron Translator" << std::endl;
0574     }
0575 
0576     // compute the width
0577     PFClusterWidthAlgo pfwidth(pfClusters_[iGSF]);
0578 
0579     double correctedEnergy = pfCand[gsfPFCandidateIndex_[iGSF]].ecalEnergy();
0580     reco::SuperCluster mySuperCluster(correctedEnergy, math::XYZPoint(posX, posY, posZ));
0581     // protection against empty basic cluster collection ; the value is -2 in this case
0582     if (nbasics) {
0583       //      std::cout << "SuperCluster creation; energy " << pfCand[gsfPFCandidateIndex_[iGSF]].ecalEnergy();
0584       //      std::cout << " " <<   pfCand[gsfPFCandidateIndex_[iGSF]].rawEcalEnergy() << std::endl;
0585       //      std::cout << "Seed energy from basic " << basicClusters_[iGSF][0].energy() << std::endl;
0586       mySuperCluster.setSeed(basicClusterPtr_[iGSF][0]);
0587     } else {
0588       //      std::cout << "SuperCluster creation ; seed energy " << 0 << std::endl;
0589       //      std::cout << "SuperCluster creation ; energy " << pfCand[gsfPFCandidateIndex_[iGSF]].ecalEnergy();
0590       //      std::cout << " " <<   pfCand[gsfPFCandidateIndex_[iGSF]].rawEcalEnergy() << std::endl;
0591       //      std::cout << " No seed found " << 0 << std::endl;
0592       //      std::cout << " MVA " << pfCand[gsfPFCandidateIndex_[iGSF]].mva_e_pi() << std::endl;
0593       mySuperCluster.setSeed(reco::CaloClusterPtr());
0594     }
0595     // the seed should be the first basic cluster
0596 
0597     for (unsigned ibc = 0; ibc < nbasics; ++ibc) {
0598       mySuperCluster.addCluster(basicClusterPtr_[iGSF][ibc]);
0599       //      std::cout <<"Adding Ref to SC " << basicClusterPtr_[iGSF][ibc].index() << std::endl;
0600       const std::vector<std::pair<DetId, float>>& v1 = basicClusters_[iGSF][ibc].hitsAndFractions();
0601       //      std::cout << " Number of cells " << v1.size() << std::endl;
0602       for (std::vector<std::pair<DetId, float>>::const_iterator diIt = v1.begin(); diIt != v1.end(); ++diIt) {
0603         //      std::cout << " Adding DetId " << (diIt->first).rawId() << " " << diIt->second << std::endl;
0604         mySuperCluster.addHitAndFraction(diIt->first, diIt->second);
0605       }  // loop over rechits
0606     }
0607 
0608     unsigned nps = preshowerClusterPtr_[iGSF].size();
0609     for (unsigned ips = 0; ips < nps; ++ips) {
0610       mySuperCluster.addPreshowerCluster(preshowerClusterPtr_[iGSF][ips]);
0611     }
0612 
0613     // Set the preshower energy
0614     mySuperCluster.setPreshowerEnergy(pfCand[gsfPFCandidateIndex_[iGSF]].pS1Energy() +
0615                                       pfCand[gsfPFCandidateIndex_[iGSF]].pS2Energy());
0616 
0617     // Set the cluster width
0618     mySuperCluster.setEtaWidth(pfwidth.pflowEtaWidth());
0619     mySuperCluster.setPhiWidth(pfwidth.pflowPhiWidth());
0620     // Force the computation of rawEnergy_ of the reco::SuperCluster
0621     mySuperCluster.rawEnergy();
0622     superClusters.push_back(mySuperCluster);
0623   }
0624 }
0625 
0626 const reco::PFCandidate& PFElectronTranslator::correspondingDaughterCandidate(const reco::PFCandidate& cand,
0627                                                                               const reco::PFBlockElement& pfbe) const {
0628   unsigned refindex = pfbe.index();
0629   //  std::cout << " N daughters " << cand.numberOfDaughters() << std::endl;
0630   reco::PFCandidate::const_iterator myDaughterCandidate = cand.begin();
0631   reco::PFCandidate::const_iterator itend = cand.end();
0632 
0633   for (; myDaughterCandidate != itend; ++myDaughterCandidate) {
0634     const reco::PFCandidate* myPFCandidate = (const reco::PFCandidate*)&*myDaughterCandidate;
0635     if (myPFCandidate->elementsInBlocks().size() != 1) {
0636       //      std::cout << " Daughter with " << myPFCandidate.elementsInBlocks().size()<< " element in block " << std::endl;
0637       return cand;
0638     }
0639     if (myPFCandidate->elementsInBlocks()[0].second == refindex) {
0640       //      std::cout << " Found it " << cand << std::endl;
0641       return *myPFCandidate;
0642     }
0643   }
0644   return cand;
0645 }
0646 
0647 void PFElectronTranslator::createGsfElectronCores(reco::GsfElectronCoreCollection& gsfElectronCores) const {
0648   unsigned nGSF = GsfTrackRef_.size();
0649   for (unsigned iGSF = 0; iGSF < nGSF; ++iGSF) {
0650     reco::GsfElectronCore myElectronCore(GsfTrackRef_[iGSF]);
0651     myElectronCore.setCtfTrack(kfTrackRef_[iGSF], -1.);
0652     std::map<reco::GsfTrackRef, reco::SuperClusterRef>::const_iterator itcheck = scMap_.find(GsfTrackRef_[iGSF]);
0653     if (itcheck != scMap_.end())
0654       myElectronCore.setParentSuperCluster(itcheck->second);
0655     gsfElectronCores.push_back(myElectronCore);
0656   }
0657 }
0658 
0659 void PFElectronTranslator::createGsfElectronCoreRefs(
0660     const edm::OrphanHandle<reco::GsfElectronCoreCollection>& gsfElectronCoreHandle) {
0661   unsigned size = GsfTrackRef_.size();
0662 
0663   for (unsigned iGSF = 0; iGSF < size; ++iGSF)  // loop on tracks
0664   {
0665     edm::Ref<reco::GsfElectronCoreCollection> elecCoreRef(gsfElectronCoreHandle, iGSF);
0666     gsfElectronCoreRefs_.push_back(elecCoreRef);
0667   }
0668 }
0669 
0670 void PFElectronTranslator::getAmbiguousGsfTracks(const reco::PFBlockElement& PFBE,
0671                                                  std::vector<reco::GsfTrackRef>& tracks) const {
0672   const reco::PFBlockElementGsfTrack* GsfEl = dynamic_cast<const reco::PFBlockElementGsfTrack*>(&PFBE);
0673   if (GsfEl == nullptr)
0674     return;
0675   const std::vector<reco::GsfPFRecTrackRef>& ambPFRecTracks(GsfEl->GsftrackRefPF()->convBremGsfPFRecTrackRef());
0676   unsigned ntracks = ambPFRecTracks.size();
0677   for (unsigned it = 0; it < ntracks; ++it) {
0678     tracks.push_back(ambPFRecTracks[it]->gsfTrackRef());
0679   }
0680 }
0681 
0682 void PFElectronTranslator::createGsfElectrons(const reco::PFCandidateCollection& pfcand,
0683                                               const IsolationValueMaps& isolationValues,
0684                                               reco::GsfElectronCollection& gsfelectrons) {
0685   unsigned size = GsfTrackRef_.size();
0686 
0687   for (unsigned iGSF = 0; iGSF < size; ++iGSF)  // loop on tracks
0688   {
0689     const reco::PFCandidate& pfCandidate(pfcand[gsfPFCandidateIndex_[iGSF]]);
0690     // Electron
0691     reco::GsfElectron myElectron(gsfElectronCoreRefs_[iGSF]);
0692     // Warning set p4 error !
0693     myElectron.setP4(reco::GsfElectron::P4_PFLOW_COMBINATION, pfCandidate.p4(), pfCandidate.deltaP(), true);
0694 
0695     // MVA inputs
0696     reco::GsfElectron::MvaInput myMvaInput;
0697     myMvaInput.earlyBrem = pfCandidate.electronExtraRef()->mvaVariable(reco::PFCandidateElectronExtra::MVA_FirstBrem);
0698     myMvaInput.lateBrem = pfCandidate.electronExtraRef()->mvaVariable(reco::PFCandidateElectronExtra::MVA_LateBrem);
0699     myMvaInput.deltaEta =
0700         pfCandidate.electronExtraRef()->mvaVariable(reco::PFCandidateElectronExtra::MVA_DeltaEtaTrackCluster);
0701     myMvaInput.sigmaEtaEta = pfCandidate.electronExtraRef()->sigmaEtaEta();
0702     myMvaInput.hadEnergy = pfCandidate.electronExtraRef()->hadEnergy();
0703 
0704     // Mustache
0705     reco::Mustache myMustache(mustacheSCParams_);
0706     myMustache.MustacheID(
0707         *(myElectron.parentSuperCluster()), myMvaInput.nClusterOutsideMustache, myMvaInput.etOutsideMustache);
0708 
0709     myElectron.setMvaInput(myMvaInput);
0710 
0711     // MVA output
0712     reco::GsfElectron::MvaOutput myMvaOutput;
0713     myMvaOutput.status = pfCandidate.electronExtraRef()->electronStatus();
0714     myMvaOutput.mva_e_pi = pfCandidate.mva_e_pi();
0715     myElectron.setMvaOutput(myMvaOutput);
0716 
0717     // ambiguous tracks
0718     unsigned ntracks = ambiguousGsfTracks_[iGSF].size();
0719     for (unsigned it = 0; it < ntracks; ++it) {
0720       myElectron.addAmbiguousGsfTrack(ambiguousGsfTracks_[iGSF][it]);
0721     }
0722 
0723     // isolation
0724     if (!isolationValues.empty()) {
0725       reco::GsfElectron::PflowIsolationVariables myPFIso;
0726       myPFIso.sumChargedHadronPt = (*isolationValues[0])[CandidatePtr_[iGSF]];
0727       myPFIso.sumPhotonEt = (*isolationValues[1])[CandidatePtr_[iGSF]];
0728       myPFIso.sumNeutralHadronEt = (*isolationValues[2])[CandidatePtr_[iGSF]];
0729       myPFIso.sumPUPt = (*isolationValues[3])[CandidatePtr_[iGSF]];
0730       myElectron.setPfIsolationVariables(myPFIso);
0731     }
0732 
0733     gsfelectrons.push_back(myElectron);
0734   }
0735 }