CMS 3D CMS Logo

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