File indexing completed on 2023-10-25 10:01:50
0001 #include "FWCore/Framework/interface/Event.h"
0002 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0003 #include "FWCore/Framework/interface/MakerMacros.h"
0004 #include "FWCore/Framework/interface/ESHandle.h"
0005 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0006 #include "RecoParticleFlow/PFProducer/test/PFSuperClusterReader.h"
0007
0008 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
0009 #include "DataFormats/EgammaReco/interface/BasicCluster.h"
0010 #include "DataFormats/EgammaReco/interface/PreshowerCluster.h"
0011 #include "DataFormats/EgammaReco/interface/SuperCluster.h"
0012 #include "DataFormats/GsfTrackReco/interface/GsfTrack.h"
0013 #include "DataFormats/ParticleFlowReco/interface/PFBlockElement.h"
0014 #include "DataFormats/ParticleFlowReco/interface/PFBlockFwd.h"
0015 #include "DataFormats/ParticleFlowReco/interface/PFBlock.h"
0016 #include "DataFormats/ParticleFlowReco/interface/PFCluster.h"
0017
0018 PFSuperClusterReader::PFSuperClusterReader(const edm::ParameterSet& iConfig)
0019 : inputTagGSFTracks_(iConfig.getParameter<edm::InputTag>("GSFTracks")),
0020 inputTagValueMapSC_(iConfig.getParameter<edm::InputTag>("SuperClusterRefMap")),
0021 inputTagValueMapMVA_(iConfig.getParameter<edm::InputTag>("MVAMap")),
0022 inputTagPFCandidates_(iConfig.getParameter<edm::InputTag>("PFCandidate")),
0023 trackToken_(consumes<reco::GsfTrackCollection>(inputTagGSFTracks_)),
0024 pfCandToken_(consumes<reco::PFCandidateCollection>(inputTagPFCandidates_)),
0025 pfClusToken_(consumes<edm::ValueMap<reco::SuperClusterRef> >(inputTagValueMapSC_)),
0026 pfMapToken_(consumes<edm::ValueMap<float> >(inputTagValueMapMVA_)) {}
0027
0028 void PFSuperClusterReader::analyze(const edm::Event& iEvent, const edm::EventSetup& c) {
0029 const edm::Handle<reco::GsfTrackCollection>& gsfTracksH = iEvent.getHandle(trackToken_);
0030 if (!gsfTracksH.isValid()) {
0031 std::ostringstream err;
0032 err << " cannot get GsfTracks: " << inputTagGSFTracks_ << std::endl;
0033 edm::LogError("PFSuperClusterReader") << err.str();
0034 throw cms::Exception("MissingProduct", err.str());
0035 }
0036
0037 const edm::Handle<reco::PFCandidateCollection>& pfCandidatesH = iEvent.getHandle(pfCandToken_);
0038 if (!pfCandidatesH.isValid()) {
0039 std::ostringstream err;
0040 err << " cannot get PFCandidates: " << inputTagPFCandidates_ << std::endl;
0041 edm::LogError("PFSuperClusterReader") << err.str();
0042 throw cms::Exception("MissingProduct", err.str());
0043 }
0044
0045 const edm::Handle<edm::ValueMap<reco::SuperClusterRef> >& pfClusterTracksH = iEvent.getHandle(pfClusToken_);
0046 if (!pfClusterTracksH.isValid()) {
0047 std::ostringstream err;
0048 err << " cannot get SuperClusterRef Map: " << inputTagValueMapSC_ << std::endl;
0049 edm::LogError("PFSuperClusterReader") << err.str();
0050 throw cms::Exception("MissingProduct", err.str());
0051 }
0052
0053 const edm::Handle<edm::ValueMap<float> >& pfMVAH = iEvent.getHandle(pfMapToken_);
0054 if (!pfMVAH.isValid()) {
0055 std::ostringstream err;
0056 err << " cannot get MVA Map: " << inputTagValueMapSC_ << std::endl;
0057 edm::LogError("PFSuperClusterReader") << err.str();
0058 throw cms::Exception("MissingProduct", err.str());
0059 }
0060
0061 const edm::ValueMap<reco::SuperClusterRef>& mySCValueMap = *pfClusterTracksH;
0062 const edm::ValueMap<float>& myMVAValueMap = *pfMVAH;
0063
0064 unsigned ngsf = gsfTracksH->size();
0065 for (unsigned igsf = 0; igsf < ngsf; ++igsf) {
0066 reco::GsfTrackRef theTrackRef(gsfTracksH, igsf);
0067 if (mySCValueMap[theTrackRef].isNull()) {
0068 continue;
0069 }
0070 const reco::SuperCluster& mySuperCluster(*(mySCValueMap[theTrackRef]));
0071 float mva = myMVAValueMap[theTrackRef];
0072 float eta = mySuperCluster.position().eta();
0073 float et = mySuperCluster.energy() * sin(mySuperCluster.position().theta());
0074 edm::LogVerbatim("PFSuperClusterReader")
0075 << " Super Cluster energy, eta, Et , EtaWidth, PhiWidth " << mySuperCluster.energy() << " " << eta << " " << et
0076 << " " << mySuperCluster.etaWidth() << " " << mySuperCluster.phiWidth();
0077
0078 const reco::PFCandidate* myPFCandidate = findPFCandidate(pfCandidatesH.product(), theTrackRef);
0079
0080 if (mySuperCluster.seed().isNull()) {
0081 continue;
0082 }
0083 edm::LogVerbatim("PFSuperClusterReader") << " Super Cluster seed energy " << mySuperCluster.seed()->energy();
0084 edm::LogVerbatim("PFSuperClusterReader") << " Preshower contribution " << mySuperCluster.preshowerEnergy();
0085 edm::LogVerbatim("PFSuperClusterReader") << " MVA value " << mva;
0086 edm::LogVerbatim("PFSuperClusterReader") << " List of basic clusters ";
0087 reco::CaloCluster_iterator it = mySuperCluster.clustersBegin();
0088 reco::CaloCluster_iterator it_end = mySuperCluster.clustersEnd();
0089 float etotbasic = 0;
0090 for (; it != it_end; ++it) {
0091 edm::LogVerbatim("PFSuperClusterReader") << " Basic cluster " << (*it)->energy();
0092 etotbasic += (*it)->energy();
0093 }
0094 it = mySuperCluster.preshowerClustersBegin();
0095 it_end = mySuperCluster.preshowerClustersEnd();
0096 for (; it != it_end; ++it) {
0097 edm::LogVerbatim("PFSuperClusterReader") << " Preshower cluster " << (*it)->energy();
0098 }
0099
0100 edm::LogVerbatim("PFSuperClusterReader") << " Comparison with PFCandidate : Energy " << myPFCandidate->ecalEnergy()
0101 << " SC : " << mySuperCluster.energy();
0102 std::ostringstream st1;
0103 st1 << " Sum of Basic clusters :" << etotbasic;
0104 st1 << " Calibrated preshower energy : " << mySuperCluster.preshowerEnergy();
0105 etotbasic += mySuperCluster.preshowerEnergy();
0106 st1 << " Basic Clusters + PS :" << etotbasic;
0107 edm::LogVerbatim("PFSuperClusterReader") << st1.str();
0108 edm::LogVerbatim("PFSuperClusterReader") << " Summary " << mySuperCluster.energy() << " " << eta << " " << et << " "
0109 << mySuperCluster.preshowerEnergy() << " " << mva << std::endl;
0110 }
0111 }
0112
0113 const reco::PFCandidate* PFSuperClusterReader::findPFCandidate(const reco::PFCandidateCollection* coll,
0114 const reco::GsfTrackRef& ref) {
0115 const reco::PFCandidate* result = 0;
0116 unsigned ncand = coll->size();
0117 for (unsigned icand = 0; icand < ncand; ++icand) {
0118 if (!(*coll)[icand].gsfTrackRef().isNull() && (*coll)[icand].gsfTrackRef() == ref) {
0119 result = &((*coll)[icand]);
0120 return result;
0121 }
0122 }
0123 return result;
0124 }
0125
0126 DEFINE_FWK_MODULE(PFSuperClusterReader);