CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_5/src/Validation/EventGenerator/plugins/MBUEandQCDValidation.cc

Go to the documentation of this file.
00001 /*Class MBUEandQCDValidation
00002  *  
00003  *  Class to fill dqm monitor elements from existing EDM file
00004  *
00005  *  $Date: 2011/12/29 10:53:11 $
00006  *  $Revision: 1.4 $
00007  */
00008  
00009 #include "Validation/EventGenerator/interface/MBUEandQCDValidation.h"
00010 #include "Validation/EventGenerator/interface/HepMCValidationHelper.h"
00011 
00012 #include "CLHEP/Units/defs.h"
00013 #include "CLHEP/Units/PhysicalConstants.h"
00014 #include "CLHEP/Units/SystemOfUnits.h"
00015 
00016 #include "fastjet/JetDefinition.hh"
00017 #include "fastjet/ClusterSequence.hh"
00018 
00019 using namespace edm;
00020 
00021 MBUEandQCDValidation::MBUEandQCDValidation(const edm::ParameterSet& iPSet): 
00022   _wmanager(iPSet),
00023   hepmcCollection_(iPSet.getParameter<edm::InputTag>("hepmcCollection")),
00024   genchjetCollection_(iPSet.getParameter<edm::InputTag>("genChjetsCollection")),
00025   genjetCollection_(iPSet.getParameter<edm::InputTag>("genjetsCollection")),
00026   verbosity_(iPSet.getUntrackedParameter<unsigned int>("verbosity",0))
00027 {    
00028   dbe = 0;
00029   dbe = edm::Service<DQMStore>().operator->();
00030 
00031   hepmcGPCollection.reserve(initSize);
00032   hepmcCharge.reserve(initSize);
00033 
00034   theCalo= new CaloCellManager(verbosity_);
00035 
00036   eneInCell.resize(CaloCellManager::nCaloCell);
00037 
00038 }
00039 
00040 MBUEandQCDValidation::~MBUEandQCDValidation() {
00041 
00042   delete theCalo;
00043   
00044 }
00045 
00046 void MBUEandQCDValidation::beginJob()
00047 {
00048   if(dbe){
00050         dbe->setCurrentFolder("Generator/MBUEandQCD");
00051         
00053     
00054     // Number of analyzed events
00055     nEvt = dbe->book1D("nEvt", "n analyzed Events", 1, 0., 1.);
00056 
00057     // Number of events with no forward trigger
00058     nNoFwdTrig = dbe->book1D("nNoFwdTrig", "n Events no forward trigger", 1, 0., 1.);
00059         
00060     // Number of events with a single arm forward trigger
00061     nSaFwdTrig = dbe->book1D("nSaFwdTrig", "n Events single arm forward trigger", 1, 0., 1.);
00062 
00063     // Number of events with b quark
00064     nbquark = dbe->book1D("nbquark", "n Events with b quark", 1, 0., 1.);
00065     
00066     // Number of events with c and b quark
00067     ncandbquark = dbe->book1D("ncandbquark", "n Events with c and b quark", 1, 0., 1.);
00068     
00069     // Number of events with c and no b quark
00070     ncnobquark = dbe->book1D("ncnobquark", "n Events with c and no b quark", 1, 0., 1.);
00071     
00072 
00073     // Number of selected events for QCD-09-010
00074     nEvt1 = dbe->book1D("nEvt1", "n Events QCD-09-010", 1, 0., 1.);
00075     // dNchdpt QCD-09-010
00076         dNchdpt1 = dbe->book1D("dNchdpt1", "dNchdpt QCD-09-010", 30, 0., 6.); 
00077     // dNchdeta QCD-09-010
00078     dNchdeta1 = dbe->book1D("dNchdeta1", "dNchdeta QCD-09-010", 10, -2.5, 2.5);
00079     // Number of selected events for QCD-10-001
00080 
00081     nEvt2 = dbe->book1D("nEvt2", "n Events QCD-10-001", 1, 0., 1.);
00082     // Leading track pt QCD-10-001
00083     leadTrackpt = dbe->book1D("leadTrackpt", "leading track pt QCD-10-001", 200, 0., 100.);
00084     // Leading track eta QCD-10-001
00085     leadTracketa = dbe->book1D("leadTracketa", "leading track eta QCD-10-001", 50., -2.5,2.5);
00086     // transverse charged particle density vs leading track pt
00087     nChaDenLpt = dbe->bookProfile("nChaDenLpt", "charged density vs leading pt", 200, 0., 100., 0., 100., " ");
00088     // transverse charged particle density vs leading track pt
00089     sptDenLpt = dbe->bookProfile("sptDenLpt", "sum pt density vs leading pt", 200, 0., 100., 0., 300., " ");
00090     // dNchdpt QCD-10-001 transverse
00091         dNchdpt2 = dbe->book1D("dNchdpt2", "dNchdpt QCD-10-001", 200, 0., 100.); 
00092     // dNchdeta QCD-10-001 transverse
00093     dNchdeta2 = dbe->book1D("dNchdeta2", "dNchdeta QCD-10-001", 50, -2.5, 2.5);
00094     // nCha QCD-10-001 transverse
00095     nCha = dbe->book1D("nCha", "n charged QCD-10-001", 100, 0., 100.);
00096     // dNchdSpt transverse
00097     dNchdSpt = dbe->book1D("dNchdSpt", "dNchdSpt QCD-10-001", 300, 0., 300.);
00098     // dNchdphi
00099     dNchdphi = dbe->bookProfile("dNchdphi", "dNchdphi QCD-10-001", nphiBin, -180., 180., 0., 30., " ");
00100     // dSptdphi
00101     dSptdphi = dbe->bookProfile("dSptdphi", "dSptdphi QCD-10-001", nphiBin, -180., 180., 0., 30., " ");
00102 
00103     // number of charged jets QCD-10-001
00104     nChj = dbe->book1D("nChj", "n charged jets QCD-10-001", 30, 0, 30.);
00105     // dNchjdeta QCD-10-001
00106     dNchjdeta = dbe->book1D("dNchjdeta", "dNchjdeta QCD-10-001", 50, -2.5, 2.5);
00107     // dNchjdpt QCD-10-001
00108     dNchjdpt = dbe->book1D("dNchjdpt", "dNchjdpt QCD-10-001", 100, 0., 100.);
00109     // leading charged jet pt QCD-10-001
00110     leadChjpt = dbe->book1D("leadChjpt", "leadChjpt QCD-10-001", 100, 0., 100.);
00111     // leading charged jet eta QCD-10-001
00112     leadChjeta = dbe->book1D("leadChjeta", "leadChjeta QCD-10-001", 50, -2.5, 2.5);
00113     // (pt1+pt2)/ptot
00114     pt1pt2optotch = dbe->bookProfile("pt1pt2optotch", "sum 2 leading jets over ptot", 50, 0., 100., 0., 1., " ");
00115 
00116     // particle rates in tracker acceptance
00117     nPPbar = dbe->book1D("nPPbar", "nPPbar QCD-10-001", 30, 0., 30.);
00118     nKpm = dbe->book1D("nKpm", "nKpm QCD-10-001", 30, 0., 30.);
00119     nK0s = dbe->book1D("nK0s", "nK0s QCD-10-001", 30, 0., 30.);
00120     nL0 = dbe->book1D("nL0", "nL0 QCD-10-001", 30, 0., 30.);
00121     nXim = dbe->book1D("nXim", "nXim QCD-10-001", 30, 0., 30.);
00122     nOmega = dbe->book1D("nOmega", "nOmega QCD-10-001", 30, 0., 30.);
00123 
00124     pPPbar = dbe->book1D("pPPbar", "Log10(pt) PPbar QCD-10-001", 25, -2., 3.);
00125     pKpm = dbe->book1D("pKpm", "Log10(pt) Kpm QCD-10-001", 25, -2., 3.);
00126     pK0s = dbe->book1D("pK0s", "Log10(pt) K0s QCD-10-001", 25, -2., 3.);
00127     pL0 = dbe->book1D("pL0", "Log10(pt) L0 QCD-10-001", 25, -2., 3.);
00128     pXim = dbe->book1D("pXim", "Log10(pt) Xim QCD-10-001", 25, -2., 3.);
00129     pOmega = dbe->book1D("pOmega", "Log10(pt) Omega QCD-10-001", 25, -2., 3.);
00130 
00131     // neutral rate in the barrel + HF acceptance
00132     nNNbar = dbe->book1D("nNNbar", "nNNbar QCD-10-001", 30, 0., 30.);
00133     nGamma = dbe->book1D("nGamma", "nGamma QCD-10-001", 50, 0., 200.);
00134 
00135     pNNbar = dbe->book1D("pNNbar", "Log10(pt) NNbar QCD-10-001", 25, -2., 3.);
00136     pGamma = dbe->book1D("pGamma", "Log10(pt) Gamma QCD-10-001", 25, -2., 3.);
00137 
00138     // highest pt electron spectrum
00139     elePt = dbe->book1D("elePt", "highest pt electron Log10(pt)", 30, -2., 4.);
00140 
00141     // highest pt muon spectrum
00142     muoPt = dbe->book1D("muoPt", "highest pt muon Log10(pt)", 30, -2., 4.);
00143 
00144 
00145     // number of selected di-jet events
00146     nDijet = dbe->book1D("nDijet", "n Dijet Events", 1, 0., 1.);
00147     // number of jets 
00148     nj = dbe->book1D("nj", "n jets ", 30, 0, 30.);
00149     // dNjdeta 
00150     dNjdeta = dbe->book1D("dNjdeta", "dNjdeta ", 50, -5., 5.);
00151     // dNjdpt 
00152     dNjdpt = dbe->book1D("dNjdpt", "dNjdpt ", 60, 0., 300.);
00153     // (pt1+pt2)/ptot
00154     pt1pt2optot = dbe->bookProfile("pt1pt2optot", "sum 2 leading jets over Et tot ", 60, 0., 300., 0., 1., " ");
00155     // pt1-pt2
00156     pt1pt2balance = dbe->book1D("pt1pt2balance", "2 leading jets pt difference ", 10, 0., 1.);
00157     // pt1 pt2 Delta phi
00158     pt1pt2Dphi = dbe->book1D("pt1pt2Dphi", "pt1 pt2 delta phi ", nphiBin, 0., 180.);
00159     // pt1 pt2 invariant mass
00160     pt1pt2InvM = dbe->book1D("pt1pt2InvM", "pt1 pt2 invariant mass ", 60, 0., 600.);
00161     // pt3 fraction
00162     pt3Frac = dbe->book1D("pt3Frac", "2 pt3 over pt1+pt2 ", 30, 0., 1.);
00163     // sum of jets Et
00164     sumJEt = dbe->book1D("sumJEt", "sum Jet Et ", 60, 0., 300.);
00165     // fraction of missing Et over sum of jets Et
00166     missEtosumJEt = dbe->book1D("missEtosumJEt", "missing Et over sumJet Et ", 30, 0., 1.);
00167     // sum of final state particle Pt
00168     sumPt = dbe->book1D("sumPt", "sum particle Pt ", 60, 0., 600.);
00169     // sum of final state charged particle Pt
00170     sumChPt = dbe->book1D("sumChPt", "sum charged particle Pt ", 60, 0., 300.);
00171 
00172     //Number of selected events for the HF energy flux analysis
00173     nHFflow = dbe->book1D("nHFflow", "n HF flow events", 1, 0., 1.);
00174     //Forward energy flow for MinBias BSC selection
00175     dEdetaHFmb = dbe->bookProfile("dEdetaHFmb", "dEdeta HF MinBias", (int)CaloCellManager::nForwardEta, 0, (double)CaloCellManager::nForwardEta, 0., 300., " ");
00176     //Forward energy flow for QCD dijet selection
00177     dEdetaHFdj = dbe->bookProfile("dEdetaHFdj", "dEdeta HF QCD dijet", (int)CaloCellManager::nForwardEta, 0, (double)CaloCellManager::nForwardEta, 0., 300., " ");
00178 
00179     // FWD-10-001 like diffraction analysis
00180     nHFSD = dbe->book1D("nHFSD","n single diffraction in HF", 1, 0., 1.);
00181     // E-pz HF-
00182     EmpzHFm = dbe->book1D("EmpzHFm", "E-pz HF- SD", 40, 0., 200.);
00183     // Number of cells above threshold
00184     ntHFm = dbe->book1D("ntHFm", "number of HF- tower SD", 20, 0., 20.);
00185     // Energy in HF-
00186     eneHFmSel = dbe->book1D("eneHFmSel", "energy in HF-", 40, 0., 200.);
00187 
00188     // number of jets accepted in the 'Jet-Multiplicity' analysis
00189     _JM25njets = dbe->book1D("JM25njets", "n jets", 15, 0, 15.);    
00190     _JM25ht = dbe->book1D("JM25ht", "HT", 80, 0, 800.);    
00191     _JM25pt1 = dbe->book1D("JM25pt1", "pt", 40, 0, 200.);    
00192     _JM25pt2 = dbe->book1D("JM25pt2", "pt", 40, 0, 200.);    
00193     _JM25pt3 = dbe->book1D("JM25pt3", "pt", 40, 0, 200.);    
00194     _JM25pt4 = dbe->book1D("JM25pt4", "pt", 40, 0, 200.);    
00195 
00196     _JM80njets = dbe->book1D("JM80njets", "n jets", 15, 0, 15.);    
00197     _JM80ht = dbe->book1D("JM80ht", "HT", 80, 300, 1100.);    
00198     _JM80pt1 = dbe->book1D("JM80pt1", "pt", 40, 60, 260.);    
00199     _JM80pt2 = dbe->book1D("JM80pt2", "pt", 40, 60, 260.);    
00200     _JM80pt3 = dbe->book1D("JM80pt3", "pt", 40, 60, 260.);    
00201     _JM80pt4 = dbe->book1D("JM80pt4", "pt", 40, 60, 260.);    
00202 
00203 
00204     // differential jet rates
00205     djr10 = dbe->book1D("djr10", "Differential Jet Rate 1#rightarrow0", 60, -1., 5.);
00206     djr21 = dbe->book1D("djr21", "Differential Jet Rate 2#rightarrow1", 60, -1., 5.);
00207     djr32 = dbe->book1D("djr32", "Differential Jet Rate 3#rightarrow2", 60, -1., 5.);
00208     djr43 = dbe->book1D("djr43", "Differential Jet Rate 4#rightarrow3", 60, -1., 5.);
00209 
00210     // sumET analysis
00211     _sumEt = dbe->book1D("sumET", "Sum of stable particles Et", 150, 0, 600.);
00212     _sumEt1 = dbe->book1D("sumET1", "Sum of stable particles Et (eta<0.5)", 150, 0, 200.);
00213     _sumEt2 = dbe->book1D("sumET2", "Sum of stable particles Et (0.5<eta<1.0)", 150, 0, 200.);
00214     _sumEt3 = dbe->book1D("sumET3", "Sum of stable particles Et (1.0<eta<1.5)", 150, 0, 200.);
00215     _sumEt4 = dbe->book1D("sumET4", "Sum of stable particles Et (1.5<eta<2.0)", 150, 0, 200.);
00216     _sumEt5 = dbe->book1D("sumET5", "Sum of stable particles Et (2.0<eta<5.0)", 150, 0, 200.);
00217     
00218 
00219   }
00220   return;
00221 }
00222 
00223 void MBUEandQCDValidation::endJob(){return;}
00224 void MBUEandQCDValidation::beginRun(const edm::Run& iRun,const edm::EventSetup& iSetup){ 
00226   iSetup.getData( fPDGTable );
00227   return;
00228 }
00229 void MBUEandQCDValidation::endRun(const edm::Run& iRun,const edm::EventSetup& iSetup){return;}
00230 void MBUEandQCDValidation::analyze(const edm::Event& iEvent,const edm::EventSetup& iSetup)
00231 { 
00232 
00234   edm::Handle<HepMCProduct> evt;
00235   iEvent.getByLabel(hepmcCollection_, evt);
00236 
00237   //Get HepMC EVENT
00238   HepMC::GenEvent *myGenEvent = new HepMC::GenEvent(*(evt->GetEvent()));
00239 
00240   double weight = _wmanager.weight(iEvent);
00241 
00242 
00243   if ( verbosity_ > 0 ) { myGenEvent->print(); }
00244 
00245   double binW = 1.;
00246   
00247   hepmcGPCollection.clear();
00248   hepmcCharge.clear();
00249   for (unsigned int i = 0; i < eneInCell.size(); i++) { eneInCell[i] = 0.; }
00250 
00251   nEvt->Fill(0.5,weight);
00252   
00253   //Looping through HepMC::GenParticle collection to search for status 1 particles
00254   double charge = 0.;
00255   unsigned int nb = 0;
00256   unsigned int nc = 0;
00257   for (HepMC::GenEvent::particle_const_iterator iter = myGenEvent->particles_begin(); iter != myGenEvent->particles_end(); ++iter){
00258     if ( std::fabs((*iter)->pdg_id()) == 4 ) { nc++; }
00259     if ( std::fabs((*iter)->pdg_id()) == 5 ) { nb++; }
00260     if ( (*iter)->status() == 1) {
00261       hepmcGPCollection.push_back(*iter);
00262       const HepPDT::ParticleData* PData = fPDGTable->particle(HepPDT::ParticleID((*iter)->pdg_id()));
00263       if(PData==0) { charge = -999.; }
00264       else
00265         charge = PData->charge();
00266 
00267       hepmcCharge.push_back(charge);
00268       
00269       if ( verbosity_ > 0 ) {
00270         std::cout << "HepMC " << std::setw(14) << std::fixed << (*iter)->barcode() 
00271                   << std::setw(14) << std::fixed << (*iter)->pdg_id() 
00272                   << std::setw(14) << std::fixed << (*iter)->momentum().perp() 
00273                   << std::setw(14) << std::fixed << (*iter)->momentum().eta() 
00274                   << std::setw(14) << std::fixed << (*iter)->momentum().phi() << std::endl;
00275       }
00276 
00277     }
00278   }
00279   
00280   int nBSCp = 0; int nBSCm = 0; double eneHFp = 0.; double eneHFm = 0.; int nChapt05 = 0; int nChaVtx = 0;
00281   for (unsigned int i = 0; i < hepmcGPCollection.size(); i++ ){
00282     if ( !isNeutrino(i) ) {
00283 
00284       // BSC trigger
00285 
00286       if ( hepmcGPCollection[i]->momentum().eta() > 3.23 && hepmcGPCollection[i]->momentum().eta() < 4.65 ) { nBSCp++; }
00287       if ( hepmcGPCollection[i]->momentum().eta() < -3.23 && hepmcGPCollection[i]->momentum().eta() > -4.65 ) { nBSCm++; }
00288 
00289       // number of charged particles in different selections
00290 
00291       if ( std::fabs(hepmcGPCollection[i]->momentum().eta()) < 2.5 && hepmcGPCollection[i]->momentum().perp() > 0.5 && isCharged(i) ) { nChapt05++; }
00292       if ( std::fabs(hepmcGPCollection[i]->momentum().eta()) < 2.5 && hepmcGPCollection[i]->momentum().perp() > 0.1 && isCharged(i) ) { nChaVtx++; }
00293       unsigned int theIndex = theCalo->getCellIndexFromAngle(hepmcGPCollection[i]->momentum().eta(),hepmcGPCollection[i]->momentum().phi());
00294       if ( theIndex < CaloCellManager::nCaloCell ) eneInCell[theIndex] += hepmcGPCollection[i]->momentum().rho();
00295     }
00296   }
00297 
00298   // Forward calorimeters energy
00299 
00300   for (unsigned int icell = CaloCellManager::nBarrelCell+CaloCellManager::nEndcapCell; icell < CaloCellManager::nCaloCell; icell++ ) {
00301     if ( theCalo->getCellFromIndex(icell)->getEtaMin() < 0. ) { eneHFm += eneInCell[icell]; }
00302     else { eneHFp += eneInCell[icell]; }
00303   }
00304 
00305   // QCD-09-010 selection
00306   bool sel1 = false;
00307   if ( (nBSCp > 0 || nBSCm > 0) && eneHFp >= 3. && eneHFm >= 3. ) { sel1 = true; }
00308 
00309   // QCD-10-001 selection
00310   bool sel2 = false;
00311   if ( (nBSCp >0 || nBSCm > 0) && nChaVtx >= 3 && nChapt05 > 1 ) { sel2 = true; }
00312   
00313   // no forward trigger selection
00314   bool sel3 = false;
00315   if ( nBSCp == 0 && nBSCm == 0  ) { sel3 = true; }
00316   
00317   // single arm forward trigger selection
00318   bool sel4 = false;
00319   if ( ( nBSCp>0 && nBSCm == 0 ) || ( nBSCm>0 && nBSCp == 0 ) ) { sel4 = true; }
00320   
00321   // BSC selection
00322   bool sel5 = false;
00323   if ( nBSCp > 0 && nBSCm > 0 ) { sel5 = true; }
00324   
00325   // basic JME-10-001, FWD-10-002 and Jet-Multiplicity selection
00326   bool sel6 = false;
00327   if ( sel5 && nChaVtx > 3 ) { sel6 = true; }
00328 
00329   // FWD-10-001 selection
00330   bool sel7 = false;
00331   if ( nChaVtx >= 3 && nBSCm > 0 && eneHFp < 8. ) { sel7 = true; }
00332 
00333   // Fill selection histograms
00334   if ( sel1 ) nEvt1->Fill(0.5,weight);
00335   if ( sel2 ) nEvt2->Fill(0.5,weight);
00336   if ( sel3 ) nNoFwdTrig->Fill(0.5,weight);
00337   if ( sel4 ) nSaFwdTrig->Fill(0.5,weight);
00338   if ( sel6 ) nHFflow->Fill(0.5,weight);
00339   if ( sel7 ) nHFSD->Fill(0.5,weight);
00340   
00341   if ( nb > 0 ) nbquark->Fill(0.5,weight);
00342   if ( nb > 0 && nc > 0 ) ncandbquark->Fill(0.5,weight);
00343   if ( nb == 0 && nc > 0 ) ncnobquark->Fill(0.5,weight);
00344 
00345   // track analyses 
00346   double ptMax = 0.;
00347   unsigned int iMax = 0;
00348   double ptot = 0.;
00349   unsigned int ppbar = 0; unsigned int nnbar = 0; unsigned int kpm = 0; unsigned int k0s = 0; unsigned int l0 = 0; unsigned int gamma = 0; 
00350   unsigned int xim = 0; unsigned int omega = 0;
00351   unsigned int ele = 0; unsigned int muo = 0;
00352   unsigned int eleMax = 0;
00353   unsigned int muoMax = 0;
00354 
00355   std::vector<double> hfMB (CaloCellManager::nForwardEta,0);
00356   std::vector<double> hfDJ (CaloCellManager::nForwardEta,0);
00357 
00358   for (unsigned int i = 0; i < hepmcGPCollection.size(); i++ ){
00359     double eta = hepmcGPCollection[i]->momentum().eta();
00360     double pt = hepmcGPCollection[i]->momentum().perp();
00361     int pdgId = hepmcGPCollection[i]->pdg_id();
00362     if ( isCharged(i) && std::fabs(eta) < 2.5 ) {
00363       if ( sel1 ) {
00364         // QCD-09-010 
00365         binW = dNchdpt1->getTH1()->GetBinWidth(1);
00366         dNchdpt1->Fill(pt,1./binW); // weight to account for the pt bin width
00367         binW = dNchdeta1->getTH1()->GetBinWidth(1);
00368         dNchdeta1->Fill(eta,1./binW); // weight to account for the eta bin width
00369       }
00370       // search for the leading track QCD-10-001
00371       if ( sel2 ) {
00372         if ( pt > ptMax ) { ptMax = pt; iMax = i; }
00373         ptot += pt;
00374         
00375         // identified charged particle
00376         if (std::abs(pdgId) == 2212) {
00377           ppbar++;
00378           pPPbar->Fill(std::log10(pt),weight);
00379         }
00380         else if (std::abs(pdgId) == 321) {
00381           kpm++;
00382           pKpm->Fill(std::log10(pt),weight);
00383         }
00384         else if (std::abs(pdgId) == 3312) {
00385           xim++;
00386           pXim->Fill(std::log10(pt),weight);
00387         }
00388         else if (std::abs(pdgId) == 3334) {
00389           omega++;
00390           pOmega->Fill(std::log10(pt),weight);
00391         }
00392         else if (std::abs(pdgId) == 11) {
00393           ele++;
00394           eleMax = i;
00395         }
00396         else if (std::abs(pdgId) == 13) {
00397           muo++;
00398           muoMax = i;
00399         }
00400       }
00401     }
00402     else if ( sel2 && isNeutral(i) && std::fabs(eta) < 2.5 ) {
00403       if (std::abs(pdgId) == 310) {
00404         k0s++;
00405         pK0s->Fill(std::log10(pt),weight);
00406       }
00407       else if (std::abs(pdgId) == 3122) {
00408         l0++;
00409         pL0->Fill(std::log10(pt),weight);
00410       }
00411     }
00412     else if ( sel2 && isNeutral(i) && std::fabs(eta) < 5.19 ) {
00413       if (std::abs(pdgId) == 2112) {
00414         nnbar++;
00415         pNNbar->Fill(std::log10(pt),weight);
00416       }
00417       else if (std::abs(pdgId) == 22) {
00418         gamma++;
00419         pGamma->Fill(std::log10(pt),weight);
00420       }
00421     }
00422     unsigned int iBin = getHFbin(eta);
00423     if ( sel6 && !isNeutrino(i) &&  iBin < CaloCellManager::nForwardEta ) {
00424       hfMB[iBin] += hepmcGPCollection[i]->momentum().rho();
00425     }
00426   }
00427   nPPbar->Fill(ppbar,weight);
00428   nNNbar->Fill(nnbar,weight);
00429   nKpm->Fill(kpm,weight);
00430   nK0s->Fill(k0s,weight);
00431   nL0->Fill(l0,weight);
00432   nXim->Fill(xim,weight);
00433   nOmega->Fill(omega,weight);
00434   nGamma->Fill(gamma,weight);
00435 
00436   if ( ele > 0 ) elePt->Fill(std::log10(hepmcGPCollection[eleMax]->momentum().perp()),weight);
00437   if ( muo > 0 ) muoPt->Fill(std::log10(hepmcGPCollection[muoMax]->momentum().perp()),weight);
00438 
00439   leadTrackpt->Fill(hepmcGPCollection[iMax]->momentum().perp(),weight); 
00440   leadTracketa->Fill(hepmcGPCollection[iMax]->momentum().eta(),weight); 
00441 
00442   std::vector<double> theEtaRanges(theCalo->getEtaRanges());
00443 
00444   for (unsigned int i = 0; i < CaloCellManager::nForwardEta; i++ ) {
00445     binW = theEtaRanges[CaloCellManager::nBarrelEta+CaloCellManager::nEndcapEta+i+1]-theEtaRanges[CaloCellManager::nBarrelEta+CaloCellManager::nEndcapEta+i];
00446     dEdetaHFmb->Fill(i+0.5,hfMB[i]/binW);
00447   }
00448 
00449   // FWD-10-001
00450 
00451   if ( sel7 ) {
00452 
00453     double empz = 0.;
00454     unsigned int nCellOvTh = 0;
00455     double threshold = 0.;
00456 
00457     for (unsigned int icell = 0; icell < eneInCell.size(); icell++ ) {
00458 
00459       if ( theCalo->getCellFromIndex(icell)->getSubSys() != CaloCellId::Forward ) { threshold = 3.; }
00460       else { threshold = 4.; }
00461 
00462       if ( eneInCell[icell] > threshold ) {
00463         if ( theCalo->getCellFromIndex(icell)->getSubSys() == CaloCellId::Forward ) { nCellOvTh++; } 
00464         empz += eneInCell[icell]*(1.-std::cos(theCalo->getCellFromIndex(icell)->getThetaCell()));
00465       }
00466 
00467     }
00468     
00469     EmpzHFm->Fill(empz,weight);
00470     ntHFm->Fill(nCellOvTh,weight);
00471     eneHFmSel->Fill(eneHFm,weight);
00472 
00473   }
00474   
00475   // QCD-10-001
00476   double phiMax = hepmcGPCollection[iMax]->momentum().phi();
00477   std::vector<unsigned int> nchvsphi (nphiBin,0);
00478   std::vector<double> sptvsphi (nphiBin,0.);
00479   unsigned int nChaTra = 0;
00480   double sptTra = 0.;
00481   
00482   double binPhiW = 360./nphiBin;
00483   if ( sel2 ) {
00484     for (unsigned int i = 0; i < hepmcGPCollection.size(); i++ ){
00485       if ( isCharged(i) && std::fabs(hepmcGPCollection[i]->momentum().eta()) < 2. ) {
00486         double thePhi = (hepmcGPCollection[i]->momentum().phi()-phiMax)/CLHEP::degree;
00487         if ( thePhi < -180. ) { thePhi += 360.; }
00488         else if ( thePhi > 180. ) { thePhi -= 360.; }
00489         unsigned int thePhiBin = (int)((thePhi+180.)/binPhiW);
00490         if ( thePhiBin == nphiBin ) { thePhiBin -= 1; }
00491         nchvsphi[thePhiBin]++;
00492         sptvsphi[thePhiBin] += hepmcGPCollection[i]->momentum().perp();
00493         // analysis in the transverse region
00494         if ( std::fabs(thePhi) > 60. && std::fabs(thePhi) < 120. ) {
00495           nChaTra++;
00496           sptTra += hepmcGPCollection[i]->momentum().perp();
00497           binW = dNchdpt2->getTH1()->GetBinWidth(1);
00498           dNchdpt2->Fill(hepmcGPCollection[i]->momentum().perp(),1./binW); // weight to account for the pt bin width
00499           binW = dNchdeta2->getTH1()->GetBinWidth(1);
00500           dNchdeta2->Fill(hepmcGPCollection[i]->momentum().eta(),1./binW);  // weight to account for the eta bin width
00501         }
00502       }
00503     }
00504     nCha->Fill(nChaTra,weight);
00505     binW = dNchdSpt->getTH1()->GetBinWidth(1);
00506     dNchdSpt->Fill(sptTra,1.);
00507     //how do one apply weights to a profile? MonitorElement doesn't allow to 
00508     nChaDenLpt->Fill(hepmcGPCollection[iMax]->momentum().perp(),nChaTra/4./CLHEP::twopi);
00509     sptDenLpt->Fill(hepmcGPCollection[iMax]->momentum().perp(),sptTra/4./CLHEP::twopi);
00510     for ( unsigned int i = 0; i < nphiBin; i++ ) {
00511       double thisPhi = -180.+(i+0.5)*binPhiW;
00512       dNchdphi->Fill(thisPhi,nchvsphi[i]/binPhiW/4.); // density in phi and eta
00513       dSptdphi->Fill(thisPhi,sptvsphi[i]/binPhiW/4.); // density in phi and eta
00514     }
00515   }
00516   
00517   // Gather information in the charged GenJet collection
00518   edm::Handle<reco::GenJetCollection> genChJets;
00519   iEvent.getByLabel(genchjetCollection_, genChJets );
00520   
00521   unsigned int nJets = 0;
00522   double pt1 = 0.; 
00523   double pt2 = 0.;
00524   reco::GenJetCollection::const_iterator ij1 = genChJets->begin();
00525   reco::GenJetCollection::const_iterator ij2 = genChJets->begin();
00526   if ( sel2 ) {
00527     for (reco::GenJetCollection::const_iterator iter=genChJets->begin();iter!=genChJets->end();++iter){
00528       double eta = (*iter).eta();
00529       double pt = (*iter).pt();
00530       if ( verbosity_ > 0 ) { 
00531         std::cout << "GenJet " << std::setw(14) << std::fixed << (*iter).pt() 
00532                   << std::setw(14) << std::fixed << (*iter).eta() 
00533                   << std::setw(14) << std::fixed << (*iter).phi() << std::endl;
00534       }
00535       if ( std::fabs(eta) < 2. ) {
00536         nJets++;
00537         binW = dNchjdeta->getTH1()->GetBinWidth(1);
00538         dNchjdeta->Fill(eta,1./binW);
00539         binW = dNchjdpt->getTH1()->GetBinWidth(1);
00540         dNchjdpt->Fill(pt,1./binW);
00541         if ( pt >= pt1 ) { pt1 = pt; ij1 = iter; }
00542         if ( pt < pt1 && pt >= pt2 ) { pt2 = pt; ij2 = iter; }
00543       }
00544     }
00545     
00546     nChj->Fill(nJets,weight);
00547     if ( nJets > 0 && ij1 != genChJets->end() ) {
00548       leadChjpt->Fill(pt1,weight);
00549       leadChjeta->Fill((*ij1).eta(),weight);
00550       if ( nJets > 1 && ij2 != genChJets->end() ) {
00551         pt1pt2optotch->Fill(pt1+pt2,(pt1+pt2)/ptot);
00552       }
00553     }
00554   }
00555     
00556   
00557   // Gather information in the GenJet collection
00558   edm::Handle<reco::GenJetCollection> genJets;
00559   iEvent.getByLabel(genjetCollection_, genJets );
00560   
00561   nJets = 0;
00562   pt1 = 0.; 
00563   pt2 = 0.;
00564   double pt3 = 0.;
00565 
00566 
00567   // needed for Jet-Multiplicity Analysis
00568   int jm25njets  = 0;
00569   double jm25HT  = 0.;
00570   double jm25pt1 = 0.;
00571   double jm25pt2 = 0.;
00572   double jm25pt3 = 0.;
00573   double jm25pt4 = 0.;
00574 
00575   int jm80njets  = 0;
00576   double jm80HT  = 0.;
00577   double jm80pt1 = 0.;
00578   double jm80pt2 = 0.;
00579   double jm80pt3 = 0.;
00580   double jm80pt4 = 0.;
00581 
00582 
00583 
00584   reco::GenJetCollection::const_iterator ij3 = genJets->begin();
00585   if ( sel6 ) {
00586     for (reco::GenJetCollection::const_iterator iter=genJets->begin();iter!=genJets->end();++iter){
00587       double eta = (*iter).eta();
00588       double pt = (*iter).pt();
00589       if ( verbosity_ > 0 ) {
00590         std::cout << "GenJet " << std::setw(14) << std::fixed << (*iter).pt() 
00591                   << std::setw(14) << std::fixed << (*iter).eta() 
00592                   << std::setw(14) << std::fixed << (*iter).phi() << std::endl;
00593       }
00594       if ( std::fabs(eta) < 5. ) {
00595         nJets++;
00596         if ( pt >= pt1 ) { pt1 = pt; ij1 = iter; }
00597         if ( pt < pt1 && pt >= pt2 ) { pt2 = pt; ij2 = iter; }
00598         if ( pt < pt2 && pt >= pt3 ) { pt3 = pt; ij3 = iter; }
00599       }
00600 
00601       // find variables for Jet-Multiplicity Analysis
00602       if(fabs(iter->eta()) < 3. && iter->pt()>25.) {
00603         jm25njets++;
00604         jm25HT += iter->pt();
00605         if(iter->pt()>jm25pt1) {
00606           jm25pt4 = jm25pt3;
00607           jm25pt3 = jm25pt2;
00608           jm25pt2 = jm25pt1;
00609           jm25pt1 = iter->pt();
00610         } else if(iter->pt()>jm25pt2) {
00611           jm25pt4 = jm25pt3;
00612           jm25pt3 = jm25pt2;
00613           jm25pt2 = iter->pt();
00614         } else if(iter->pt()>jm25pt3) {
00615           jm25pt4 = jm25pt3;
00616           jm25pt3 = iter->pt();
00617         } else if(iter->pt()>jm25pt4) {
00618           jm25pt4 = iter->pt();
00619         }
00620         // even harder jets...
00621         if(iter->pt()>80.) {
00622           jm80njets++;
00623           jm80HT += iter->pt();
00624           if(iter->pt()>jm80pt1) {
00625             jm80pt4 = jm80pt3;
00626             jm80pt3 = jm80pt2;
00627             jm80pt2 = jm80pt1;
00628             jm80pt1 = iter->pt();
00629           } else if(iter->pt()>jm80pt2) {
00630             jm80pt4 = jm80pt3;
00631             jm80pt3 = jm80pt2;
00632             jm80pt2 = iter->pt();
00633           } else if(iter->pt()>jm80pt3) {
00634             jm80pt4 = jm80pt3;
00635             jm80pt3 = iter->pt();
00636           } else if(iter->pt()>jm80pt4) {
00637             jm80pt4 = iter->pt();
00638           }
00639         }
00640         
00641       }
00642 
00643       if(jm25njets>3) {
00644         _JM25njets ->Fill(jm25njets,weight);
00645         _JM25ht    ->Fill(jm25HT,weight);
00646         _JM25pt1   ->Fill(jm25pt1,weight);
00647         _JM25pt2   ->Fill(jm25pt2,weight);
00648         _JM25pt3   ->Fill(jm25pt3,weight);
00649         _JM25pt4   ->Fill(jm25pt4,weight);
00650       }
00651       if(jm80njets>3) {
00652         _JM80njets ->Fill(jm80njets,weight);
00653         _JM80ht    ->Fill(jm80HT,weight);
00654         _JM80pt1   ->Fill(jm80pt1,weight);
00655         _JM80pt2   ->Fill(jm80pt2,weight);
00656         _JM80pt3   ->Fill(jm80pt3,weight);
00657         _JM80pt4   ->Fill(jm80pt4,weight);
00658       }
00659     }
00660     
00661     // select a di-jet event JME-10-001 variant
00662     double sumJetEt = 0; double sumPartPt = 0.; double sumChPartPt = 0.;
00663     double jpx = 0; double jpy = 0;
00664     if ( nJets >= 2 && ij1 != genJets->end() && ij2 != genJets->end() ) {
00665       if ( (*ij1).pt() > 25. && (*ij1).pt() > 25. ) {
00666         double deltaPhi = std::fabs((*ij1).phi()-(*ij2).phi())/CLHEP::degree;
00667         if ( deltaPhi > 180. ) deltaPhi = 360.-deltaPhi;
00668         pt1pt2Dphi->Fill(deltaPhi,weight);
00669         if ( std::fabs(deltaPhi) > 2.5*CLHEP::degree ) {
00670 
00671           nDijet->Fill(0.5,weight);
00672 
00673           for (unsigned int i = 0; i < hepmcGPCollection.size(); i++ ){
00674             double eta = hepmcGPCollection[i]->momentum().eta();
00675             unsigned int iBin = getHFbin(eta);
00676             if ( !isNeutrino(i) &&  iBin < CaloCellManager::nForwardEta ) {
00677               hfDJ[iBin] += hepmcGPCollection[i]->momentum().rho();
00678             }
00679             if ( !isNeutrino(i) && std::fabs(eta) < 5. ) {
00680               sumPartPt += hepmcGPCollection[i]->momentum().perp();
00681               if ( isCharged(i) ) {
00682                 sumChPartPt += hepmcGPCollection[i]->momentum().perp();
00683               }
00684             }
00685           }
00686           for (unsigned int i = 0; i < CaloCellManager::nForwardEta; i++ ) {
00687             binW = theEtaRanges[CaloCellManager::nBarrelEta+CaloCellManager::nEndcapEta+i+1]-theEtaRanges[CaloCellManager::nBarrelEta+CaloCellManager::nEndcapEta+i];
00688             dEdetaHFdj->Fill(i+0.5,hfDJ[i]/binW);
00689           }
00690 
00691           double invMass = (*ij1).energy()*(*ij2).energy()-(*ij1).px()*(*ij2).px()-(*ij1).py()*(*ij2).py()-(*ij1).pz()*(*ij2).pz();
00692           invMass = std::sqrt(invMass);
00693           pt1pt2InvM->Fill(invMass,weight);
00694 
00695           sumPt->Fill(sumPartPt,weight);
00696           sumChPt->Fill(sumChPartPt,weight);
00697 
00698           unsigned int nSelJets = 0;
00699           for (reco::GenJetCollection::const_iterator iter=genJets->begin();iter!=genJets->end();++iter){
00700             double pt = (*iter).pt();
00701             double eta = (*iter).eta();
00702             if ( std::fabs(eta) < 5. ) { 
00703               nSelJets++; 
00704               binW = dNjdeta->getTH1()->GetBinWidth(1);
00705               dNjdeta->Fill(eta,1./binW*weight);
00706               binW = dNjdpt->getTH1()->GetBinWidth(1);
00707               dNjdpt->Fill(pt,1./binW*weight);
00708               sumJetEt += (*iter).pt();
00709               jpx += (*iter).px();
00710               jpy += (*iter).py();
00711             }
00712           }
00713 
00714           nj->Fill(nSelJets,weight);
00715           double mEt = std::sqrt(jpx*jpx+jpy*jpy);
00716           sumJEt->Fill(sumJetEt,weight);
00717           missEtosumJEt->Fill(mEt/sumJetEt,weight);
00718 
00719           if ( nSelJets >= 3 ) { pt3Frac->Fill((*ij3).pt()/(pt1+pt2),weight); }
00720 
00721           pt1pt2optot->Fill(pt1+pt2,(pt1+pt2)/sumJetEt);
00722           pt1pt2balance->Fill((pt1-pt2)/(pt1+pt2),weight);
00723         }
00724       }      
00725     }
00726   }
00727     
00728  
00729   //compute differential jet rates
00730   std::vector<const HepMC::GenParticle*> qcdActivity;
00731   HepMCValidationHelper::removeIsolatedLeptons(myGenEvent, 0.2, 3., qcdActivity);
00732   //HepMCValidationHelper::allStatus1(myGenEvent, qcdActivity);
00733   //fill PseudoJets to use fastjet
00734   std::vector<fastjet::PseudoJet> vecs;
00735   int counterUser = 1;
00736   std::vector<const HepMC::GenParticle*>::const_iterator iqcdact;
00737   for (iqcdact = qcdActivity.begin(); iqcdact != qcdActivity.end(); ++iqcdact){
00738     const HepMC::FourVector& fmom = (*iqcdact)->momentum(); 
00739     fastjet::PseudoJet pseudoJet(fmom.px(), fmom.py(), fmom.pz(), fmom.e());
00740     pseudoJet.set_user_index(counterUser);
00741     vecs.push_back(pseudoJet);
00742     ++counterUser;
00743   }
00744   //compute jets
00745   fastjet::ClusterSequence cseq(vecs, fastjet::JetDefinition(fastjet::kt_algorithm, 1., fastjet::E_scheme)); 
00746   //access the cluster sequence and get the relevant info
00747   djr10->Fill(std::log10(sqrt(cseq.exclusive_dmerge(0))),weight);
00748   djr21->Fill(std::log10(sqrt(cseq.exclusive_dmerge(1))),weight);
00749   djr32->Fill(std::log10(sqrt(cseq.exclusive_dmerge(2))),weight);
00750   djr43->Fill(std::log10(sqrt(cseq.exclusive_dmerge(3))),weight);
00751   
00752 
00753   // compute sumEt for all stable particles
00754   std::vector<const HepMC::GenParticle*> allStable;
00755   HepMCValidationHelper::allStatus1(myGenEvent, allStable);
00756   
00757   double sumEt  = 0.;
00758   double sumEt1 = 0.;
00759   double sumEt2 = 0.;
00760   double sumEt3 = 0.;
00761   double sumEt4 = 0.;
00762   double sumEt5 = 0.;
00763 
00764   for(std::vector<const HepMC::GenParticle*>::const_iterator iter=allStable.begin();
00765       iter != allStable.end(); ++iter) {
00766 
00767     double thisEta=fabs((*iter)->momentum().eta());
00768     
00769     if(thisEta < 5.) {
00770       const HepMC::FourVector mom=(*iter)->momentum();
00771       double px=mom.px();
00772       double py=mom.py();
00773       double pz=mom.pz();
00774       double E=mom.e();
00775       double thisSumEt = (
00776                           sqrt(px*px + py*py)*E /
00777                           sqrt(px*px + py*py + pz*pz)
00778                           );
00779       sumEt += thisSumEt;
00780       if(thisEta<1.0) sumEt1 += thisSumEt;
00781       else if(thisEta<2.0) sumEt2 += thisSumEt;
00782       else if(thisEta<3.0) sumEt3 += thisSumEt;
00783       else if(thisEta<4.0) sumEt4 += thisSumEt;
00784       else sumEt5 += thisSumEt;
00785       
00786     }
00787   }
00788   
00789   if(sumEt>0.)
00790     _sumEt->Fill(sumEt,weight);
00791   if(sumEt1>0.)
00792     _sumEt1->Fill(sumEt1,weight);
00793   if(sumEt2>0.)
00794     _sumEt2->Fill(sumEt2,weight);
00795   if(sumEt3>0.)
00796     _sumEt3->Fill(sumEt3,weight);
00797   if(sumEt4>0.)
00798     _sumEt4->Fill(sumEt4,weight);
00799   if(sumEt5>0.)
00800     _sumEt5->Fill(sumEt5,weight);
00801   
00802   delete myGenEvent;
00803 }//analyze
00804 
00805 bool MBUEandQCDValidation::isCharged(unsigned int i){
00806   
00807   bool status = false;
00808   if ( hepmcGPCollection.size() < i+1 ) { return status; }
00809   else { status = (hepmcCharge[i] != 0. && hepmcCharge[i] != -999.); }
00810   return status;
00811 
00812 }
00813 
00814 bool MBUEandQCDValidation::isNeutral(unsigned int i){
00815   
00816   bool status = false;
00817   int pdgId = std::abs(hepmcGPCollection[i]->pdg_id());
00818   if ( hepmcGPCollection.size() < i+1 ) { return status; }
00819   else { status = (hepmcCharge[i] == 0. && pdgId != 12 && pdgId != 14 && pdgId != 16) ; }
00820   return status;
00821 
00822 }
00823 
00824 bool MBUEandQCDValidation::isNeutrino(unsigned int i){
00825   
00826   bool status = false;
00827   int pdgId = std::abs(hepmcGPCollection[i]->pdg_id());
00828   if ( hepmcGPCollection.size() < i+1 ) { return status; }
00829   else { status = (pdgId == 12 || pdgId == 14 || pdgId == 16) ; }
00830   return status;
00831 
00832 }
00833 
00834 unsigned int MBUEandQCDValidation::getHFbin(double eta) {
00835 
00836   unsigned int iBin = 999;
00837 
00838   std::vector<double> theEtaRanges(theCalo->getEtaRanges());
00839 
00840   for (unsigned int i = CaloCellManager::nBarrelEta+CaloCellManager::nEndcapEta; 
00841        i < CaloCellManager::nBarrelEta+CaloCellManager::nEndcapEta+CaloCellManager::nForwardEta; i++ ){
00842     if ( std::fabs(eta) >= theEtaRanges[i] && std::fabs(eta) < theEtaRanges[i+1] ) 
00843       { iBin = i-CaloCellManager::nBarrelEta-CaloCellManager::nEndcapEta; }
00844   }
00845 
00846   return iBin;
00847 
00848 }