Go to the documentation of this file.00001
00002
00003
00004
00005
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022 #include <memory>
00023
00024
00025 #include "FWCore/Framework/interface/Frameworkfwd.h"
00026 #include "FWCore/Framework/interface/EDProducer.h"
00027
00028 #include "FWCore/Framework/interface/Event.h"
00029 #include "FWCore/Framework/interface/MakerMacros.h"
00030
00031 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00032 #include <DataFormats/RecoCandidate/interface/RecoCandidate.h>
00033 #include <DataFormats/Candidate/interface/CompositeRefCandidate.h>
00034 #include <DataFormats/MuonReco/interface/Muon.h>
00035
00036 #include <DataFormats/ParticleFlowCandidate/interface/PFCandidate.h>
00037 #include <DataFormats/ParticleFlowCandidate/interface/PFCandidateFwd.h>
00038 #include "DataFormats/TrackReco/interface/Track.h"
00039 #include "DataFormats/Math/interface/deltaR.h"
00040
00041
00042
00043
00044
00045 class PFCandidateMixer : public edm::EDProducer {
00046 public:
00047 explicit PFCandidateMixer(const edm::ParameterSet&);
00048 ~PFCandidateMixer();
00049
00050 private:
00051 virtual void beginJob() ;
00052 virtual void produce(edm::Event&, const edm::EventSetup&);
00053 virtual void endJob() ;
00054
00055 edm::InputTag _col1;
00056 edm::InputTag _col2;
00057 edm::InputTag _trackCol;
00058
00059 };
00060
00061
00062
00063
00064
00065
00066
00067
00068
00069
00070
00071
00072
00073 PFCandidateMixer::PFCandidateMixer(const edm::ParameterSet& iConfig):
00074 _col1(iConfig.getUntrackedParameter<edm::InputTag>("col1") ),
00075 _col2(iConfig.getUntrackedParameter<edm::InputTag>("col2") ),
00076 _trackCol(iConfig.getUntrackedParameter<edm::InputTag>("trackCol") )
00077 {
00078
00079 produces< std::vector< reco::PFCandidate > >();
00080
00081 }
00082
00083 PFCandidateMixer::~PFCandidateMixer()
00084 {
00085
00086
00087
00088
00089 }
00090
00091
00092
00093
00094
00095
00096
00097 void
00098 PFCandidateMixer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00099 {
00100 using namespace edm;
00101 using namespace reco;
00102
00103 Handle< std::vector<reco::Track> > trackCol;
00104 iEvent.getByLabel( _trackCol, trackCol);
00105
00106
00107 std::vector< Handle<PFCandidateCollection> > colVec;
00108 Handle<PFCandidateCollection> pfIn1;
00109 Handle<PFCandidateCollection> pfIn2;
00110 iEvent.getByLabel(_col1,pfIn1);
00111 iEvent.getByLabel(_col2,pfIn2);
00112
00113 colVec.push_back(pfIn1);
00114 colVec.push_back(pfIn2);
00115
00116 std::auto_ptr<std::vector< reco::PFCandidate > > pOut(new std::vector< reco::PFCandidate > );
00117
00118 std::vector< Handle<PFCandidateCollection> >::iterator itCol= colVec.begin();
00119 std::vector< Handle<PFCandidateCollection> >::iterator itColE= colVec.end();
00120
00121 int iCol = 0;
00122 for (;itCol!=itColE; ++itCol){
00123 if (!itCol->isValid()) {
00124 std::cout << "Whoops!" << std::endl;
00125 }
00126 PFCandidateConstIterator it = (*itCol)->begin();
00127 PFCandidateConstIterator itE = (*itCol)->end();
00128 for (;it!=itE;++it) {
00129 PFCandidate cand(*it);
00130 size_t i = 0;
00131 bool found = false;
00132 double minDR = 9999.;
00133 int iMinDr = -1;
00134 if (it->trackRef().isNonnull()) {
00135 for ( i = 0 ; i < trackCol->size(); ++i){
00136 if ( reco::deltaR( *(it->trackRef()), trackCol->at(i) )<0.001 ) {
00137 found = true;
00138 break;
00139 }
00140 double dr = reco::deltaR( *(it->trackRef()), trackCol->at(i) );
00141 if ( dr < minDR) {
00142 iMinDr = i;
00143 minDR = dr;
00144 }
00145 }
00146 }
00147 if ( found ){
00148 reco::TrackRef trref(trackCol,i);
00149 cand.setTrackRef(trref);
00150
00151
00152 } else {
00153 if (it->trackRef().isNonnull()) {
00154 std::cout << " XXXXXXXXXXX track not found "
00155 << " col " << iCol
00156 << " ch " << it->charge()
00157 << " id " << it->pdgId()
00158 << " pt " << it->pt()
00159 << " track: eta " << it->trackRef()->eta()
00160 << " pt: " << it->trackRef()->pt()
00161 << " charge: " << it->trackRef()->charge()
00162 << std::endl;
00163 std::cout << " minDR=" << minDR << std::endl;
00164 if ( iMinDr >= 0 ) {
00165 std::cout
00166 << " closest track pt=" << trackCol->at(iMinDr).pt()
00167 << " ch=" << trackCol->at(iMinDr).charge()
00168 << std::endl;
00169 }
00170 edm::Provenance prov=iEvent.getProvenance(it->trackRef().id());
00171 edm::InputTag tag(prov.moduleLabel(), prov.productInstanceName(), prov.processName());
00172 std::cout << " trackref in PFCand came from: " << tag.encode() << std::endl;
00173
00174
00175 }
00176 }
00177 pOut->push_back(cand);
00178 }
00179 ++iCol;
00180 }
00181
00182
00183 iEvent.put(pOut);
00184
00185
00186 }
00187
00188
00189 void
00190 PFCandidateMixer::beginJob()
00191 {
00192 }
00193
00194
00195 void
00196 PFCandidateMixer::endJob() {
00197 }
00198
00199
00200 DEFINE_FWK_MODULE(PFCandidateMixer);