00001
00002
00003
00004
00005
00006
00007
00008
00009 #include "Validation/EventGenerator/interface/BasicGenParticleValidation.h"
00010
00011 #include "CLHEP/Units/defs.h"
00012 #include "CLHEP/Units/PhysicalConstants.h"
00013
00014 using namespace edm;
00015
00016 BasicGenParticleValidation::BasicGenParticleValidation(const edm::ParameterSet& iPSet):
00017 hepmcCollection_(iPSet.getParameter<edm::InputTag>("hepmcCollection")),
00018 genparticleCollection_(iPSet.getParameter<edm::InputTag>("genparticleCollection")),
00019 genjetCollection_(iPSet.getParameter<edm::InputTag>("genjetsCollection")),
00020 matchPr_(iPSet.getParameter<double>("matchingPrecision")),
00021 verbosity_(iPSet.getUntrackedParameter<unsigned int>("verbosity",0))
00022 {
00023 dbe = 0;
00024 dbe = edm::Service<DQMStore>().operator->();
00025 }
00026
00027 BasicGenParticleValidation::~BasicGenParticleValidation() {}
00028
00029 void BasicGenParticleValidation::beginJob()
00030 {
00031 if(dbe){
00033 dbe->setCurrentFolder("Generator/GenParticles");
00034
00036
00037
00038 nEvt = dbe->book1D("nEvt", "n analyzed Events", 1, 0., 1.);
00039
00041 genPMultiplicity = dbe->book1D("genPMultiplicty", "Log(No. all GenParticles)", 50, -1, 5);
00042
00043 genMatched = dbe->book1D("genMatched", "Difference reco - matched", 50, -25, 25);
00044
00045 multipleMatching = dbe->book1D("multipleMatching", "multiple reco HepMC matching", 50, 0, 50);
00046
00047 matchedResolution = dbe->book1D("matchedResolution", "log10(momentum difference of matched particles)", 70, -10., -3.);
00048
00049
00050 genJetMult = dbe->book1D("genJetMult", "GenJet multiplicity", 50, 0, 50);
00051 genJetEnergy = dbe->book1D("genJetEnergy", "Log10(GenJet energy)", 60, -1, 5);
00052 genJetPt = dbe->book1D("genJetPt", "Log10(GenJet pt)", 60, -1, 5);
00053 genJetEta = dbe->book1D("genJetEta", "GenJet eta", 220, -11, 11);
00054 genJetPhi = dbe->book1D("genJetPhi", "GenJet phi", 360, -180, 180);
00055 genJetDeltaEtaMin = dbe->book1D("genJetDeltaEtaMin", "GenJet minimum rapidity gap", 30, 0, 30);
00056
00057 genJetPto1 = dbe->book1D("genJetPto1", "GenJet multiplicity above 1 GeV", 50, 0, 50);
00058 genJetPto10 = dbe->book1D("genJetPto10", "GenJet multiplicity above 10 GeV", 50, 0, 50);
00059 genJetPto100 = dbe->book1D("genJetPto100", "GenJet multiplicity above 100 GeV", 50, 0, 50);
00060 genJetCentral = dbe->book1D("genJetCentral", "GenJet multiplicity |eta|.lt.2.5", 50, 0, 50);
00061
00062 genJetTotPt = dbe->book1D("genJetTotPt", "Log10(GenJet total pt)", 100, -5, 5);
00063
00064 }
00065 return;
00066 }
00067
00068 void BasicGenParticleValidation::endJob(){return;}
00069 void BasicGenParticleValidation::beginRun(const edm::Run& iRun,const edm::EventSetup& iSetup)
00070 {
00072 iSetup.getData( fPDGTable );
00073 return;
00074 }
00075 void BasicGenParticleValidation::endRun(const edm::Run& iRun,const edm::EventSetup& iSetup){return;}
00076 void BasicGenParticleValidation::analyze(const edm::Event& iEvent,const edm::EventSetup& iSetup)
00077 {
00078
00079 unsigned int initSize = 1000;
00080
00082 edm::Handle<HepMCProduct> evt;
00083 iEvent.getByLabel(hepmcCollection_, evt);
00084
00085
00086 HepMC::GenEvent *myGenEvent = new HepMC::GenEvent(*(evt->GetEvent()));
00087
00088 nEvt->Fill(0.5);
00089
00090 std::vector<const HepMC::GenParticle*> hepmcGPCollection;
00091 std::vector<int> barcodeList;
00092 hepmcGPCollection.reserve(initSize);
00093 barcodeList.reserve(initSize);
00094
00095
00096 for (HepMC::GenEvent::particle_const_iterator iter = myGenEvent->particles_begin(); iter != myGenEvent->particles_end(); ++iter){
00097 if ( (*iter)->status() == 1) {
00098 hepmcGPCollection.push_back(*iter);
00099 barcodeList.push_back((*iter)->barcode());
00100 if ( verbosity_ > 0 ) {
00101 std::cout << "HepMC " << std::setw(14) << std::fixed << (*iter)->pdg_id() << std::setw(14) << std::fixed << (*iter)->momentum().px() << std::setw(14) << std::fixed
00102 << (*iter)->momentum().py() << std::setw(14) << std::fixed << (*iter)->momentum().pz() << std::endl;
00103 }
00104 }
00105 }
00106
00107
00108
00109 edm::Handle<reco::GenParticleCollection> genParticles;
00110 iEvent.getByLabel(genparticleCollection_, genParticles );
00111
00112 std::vector<const reco::GenParticle*> particles;
00113 particles.reserve(initSize);
00114 for (reco::GenParticleCollection::const_iterator iter=genParticles->begin();iter!=genParticles->end();++iter){
00115 if ( (*iter).status() == 1) {
00116 particles.push_back(&*iter);
00117 if ( verbosity_ > 0 ) {
00118 std::cout << "reco " << std::setw(14) << std::fixed << (*iter).pdgId() << std::setw(14) << std::fixed << (*iter).px()
00119 << std::setw(14) << std::fixed << (*iter).py() << std::setw(14) << std::fixed << (*iter).pz() << std::endl;
00120 }
00121 }
00122 }
00123
00124 unsigned int nReco = particles.size();
00125 unsigned int nHepMC = hepmcGPCollection.size();
00126
00127 genPMultiplicity->Fill(std::log10(nReco));
00128
00129
00130 std::vector<int> hepmcMatchIndex;
00131 hepmcMatchIndex.reserve(initSize);
00132
00133
00134
00135
00136
00137 if ( nReco != nHepMC ) {
00138 edm::LogWarning("CollectionSizeInconsistency") << "Collection size inconsistency: HepMC::GenParticle = " << nHepMC << " reco::GenParticle = " << nReco;
00139 }
00140
00141
00142
00143 for ( unsigned int i = 0; i < nReco; ++i ){
00144 for ( unsigned int j = 0; j < nHepMC; ++j ){
00145 if ( matchParticles( hepmcGPCollection[j], particles[i] ) ) {
00146 hepmcMatchIndex.push_back((int)j);
00147 if ( hepmcGPCollection[j]->momentum().rho() != 0. ) {
00148 double reso = 1.-particles[i]->p()/hepmcGPCollection[j]->momentum().rho();
00149 if ( verbosity_ > 0 ) {
00150 std::cout << "Matching momentum: reco = " << particles[i]->p() << " HepMC = "
00151 << hepmcGPCollection[j]->momentum().rho() << " resoultion = " << reso << std::endl;
00152 }
00153 matchedResolution->Fill(std::log10(std::fabs(reso))); }
00154 continue;
00155 }
00156 }
00157 }
00158
00159
00160
00161 unsigned int nMatched = hepmcMatchIndex.size();
00162
00163 if ( nMatched != nReco ) {
00164 edm::LogWarning("IncorrectMatching") << "Incorrect number of matched indexes: GenParticle = " << nReco << " matched indexes = " << nMatched;
00165 }
00166 genMatched->Fill(int(nReco-nMatched));
00167
00168 unsigned int nWrMatch = 0;
00169
00170 for ( unsigned int i = 0; i < nMatched; ++i ){
00171 for (unsigned int j = i+1; j < nMatched; ++j ){
00172 if ( hepmcMatchIndex[i] == hepmcMatchIndex[j] ) {
00173 int theIndex = hepmcMatchIndex[i];
00174 edm::LogWarning("DuplicatedMatching") << "Multiple matching occurencies for GenParticle barcode = " << barcodeList[theIndex];
00175 nWrMatch++;
00176 }
00177 }
00178 }
00179 multipleMatching->Fill(int(nWrMatch));
00180
00181
00182 edm::Handle<reco::GenJetCollection> genJets;
00183 iEvent.getByLabel(genjetCollection_, genJets );
00184
00185 int nJets = 0;
00186 int nJetso1 = 0;
00187 int nJetso10 = 0;
00188 int nJetso100 = 0;
00189 int nJetsCentral = 0;
00190 double totPt = 0.;
00191
00192 std::vector<double> jetEta;
00193 jetEta.reserve(initSize);
00194
00195 for (reco::GenJetCollection::const_iterator iter=genJets->begin();iter!=genJets->end();++iter){
00196 nJets++;
00197 double pt = (*iter).pt();
00198 totPt += pt;
00199 if (pt > 1.) nJetso1++;
00200 if (pt > 10.) nJetso10++;
00201 if (pt > 100.) nJetso100++;
00202 double eta = (*iter).eta();
00203 if ( std::fabs(eta) < 2.5 ) nJetsCentral++;
00204 jetEta.push_back(eta);
00205
00206 genJetEnergy->Fill(std::log10((*iter).energy()));
00207 genJetPt->Fill(std::log10(pt));
00208 genJetEta->Fill(eta);
00209 genJetPhi->Fill((*iter).phi()/CLHEP::degree);
00210 }
00211
00212 genJetMult->Fill(nJets);
00213 genJetPto1->Fill(nJetso1);
00214 genJetPto10->Fill(nJetso10);
00215 genJetPto100->Fill(nJetso100);
00216 genJetCentral->Fill(nJetsCentral);
00217
00218 genJetTotPt->Fill(std::log10(totPt));
00219
00220 double deltaEta = 999.;
00221 if ( jetEta.size() > 1 ) {
00222 for (unsigned int i = 0; i < jetEta.size(); i++){
00223 for (unsigned int j = i+1; j < jetEta.size(); j++){
00224 deltaEta = std::min(deltaEta,std::fabs(jetEta[i]-jetEta[j]));
00225 }
00226 }
00227 }
00228
00229 genJetDeltaEtaMin->Fill(deltaEta);
00230
00231 delete myGenEvent;
00232 }
00233
00234 bool BasicGenParticleValidation::matchParticles(const HepMC::GenParticle*& hepmcP, const reco::GenParticle*& recoP){
00235
00236 bool state = false;
00237
00238 if ( hepmcP->pdg_id() != recoP->pdgId() ) return state;
00239 if ( std::fabs(hepmcP->momentum().px()-recoP->px()) < std::fabs(matchPr_*hepmcP->momentum().px()) &&
00240 std::fabs(hepmcP->momentum().py()-recoP->py()) < std::fabs(matchPr_*hepmcP->momentum().py()) &&
00241 std::fabs(hepmcP->momentum().pz()-recoP->pz()) < std::fabs(matchPr_*hepmcP->momentum().pz())) {
00242 state = true; }
00243
00244 return state;
00245
00246 }