CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/GeneratorInterface/GenFilters/src/HerwigMaxPtPartonFilter.cc

Go to the documentation of this file.
00001 /*
00002 Author: Brian L. Dorney
00003 Date: July 29th 2010
00004 Version: 2.2
00005 First Release In: CMSSW_3_8_X
00006 
00007 Modified From: PythiaFilter.cc
00008 
00009 Special Thanks to Filip Moortgat
00010 
00011 PURPOSE: This Filter is designed to run on Herwig Monte Carlo Event Files
00012 (Pythia status codes are assumed to NOT BE EMULATED!!!!)
00013 
00014 For a description of Herwig Status Codes, See: 
00015 https://webber.home.cern.ch/webber/hw65_manual.html#htoc96
00016 (Section 8.3.1)
00017 
00018 This Filter Finds all final state quarks (pdg_id=1,2,3,4 or 5, status=158 or 159) with Pt>1 GeV/c
00019 that occur before the first cluster (pdg_id=91) appears in the event cascade. This is done per event.
00020 
00021 Then a histogram (which is RESET EACH EVENT) 2D histogram is formed, the Final State quarks
00022 are then plotted in eta-phi space.  These histogram entries are weighted by the quark Pt.
00023 
00024 This forms a 2D eta-phi "jet" topology for each event, and acts as a very rudimentary jet algorithm
00025 
00026 The maximum bin entry (i.e. "jet") in this histogram is the highest Pt "Jet" in the event.
00027 
00028 This is then used for filtering.
00029 
00030 The size of each bin in this 2D histogram corresponds roughly to a cone radius of 0.5
00031 
00032 i.e. This Filter Checks:
00033 minptcut <= Highest Pt "Jet" < maxptcut
00034 
00035 If this is true, the event is accepted.
00036 */
00037 
00038 #include "GeneratorInterface/GenFilters/interface/HerwigMaxPtPartonFilter.h"
00039 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
00040 
00041 #include <math.h>
00042 
00043 using namespace edm;
00044 using namespace std;
00045 
00046 HerwigMaxPtPartonFilter::HerwigMaxPtPartonFilter(const edm::ParameterSet& iConfig) :
00047   label_(iConfig.getUntrackedParameter("moduleLabel",std::string("generator"))),
00048   minptcut(iConfig.getUntrackedParameter("MinPt", 0.)),
00049   maxptcut(iConfig.getUntrackedParameter("MaxPt", 10000.)),
00050   processID(iConfig.getUntrackedParameter("ProcessID", 0)){
00051     //now do what ever initialization is needed
00052 
00053   
00054   hFSPartons_JS_PtWgting = new TH2D("hFSPartons_JS_PtWgting","#phi-#eta Space of FS Partons (p_{T} wgt'ing)",20,-5.205,5.205,32,-M_PI,M_PI);
00055   
00056 }
00057  
00058  
00059 HerwigMaxPtPartonFilter::~HerwigMaxPtPartonFilter(){
00060   
00061   // do anything here that needs to be done at desctruction time
00062   // (e.g. close files, deallocate resources etc.)
00063 
00064   if(hFSPartons_JS_PtWgting) delete hFSPartons_JS_PtWgting;
00065   
00066 }
00067  
00068  
00069 //
00070 // member functions
00071 //
00072  
00073 // ------------ method called to produce the data  ------------
00074 bool HerwigMaxPtPartonFilter::filter(edm::Event& iEvent, const edm::EventSetup& iSetup){
00075   
00076   //Histogram, reset each event
00077   hFSPartons_JS_PtWgting->Reset();
00078   
00079   bool accepted = false; //The Accept/Reject Variable
00080   bool isFSQuark = false; //Keeps track of whether a particle is a Final State Quark
00081   double maxPartonPt=0.0; //Self Explanatory
00082 
00083   //int ChosenPartonId=0, ChosenPartonSt=0;
00084 
00085   int pos1stCluster=0; //keeps track of the position of the first herwig cluster in the event
00086 
00087   //This is when Hadronization within the event occurs.
00088   long counter = 0; //keeps track of the particle index in the event
00089   
00090   Handle<HepMCProduct> evt;
00091   iEvent.getByLabel(label_, evt);
00092   
00093   const HepMC::GenEvent * myGenEvent = evt->GetEvent();
00094   
00095   
00096   if(processID == 0 || processID == myGenEvent->signal_process_id()) { //processId if statement
00097     
00098     //Find the Position of the 1st Herwig Cluster
00099     for ( HepMC::GenEvent::particle_const_iterator p = myGenEvent->particles_begin();p != myGenEvent->particles_end(); ++p ) {     
00100       if(abs((*p)->pdg_id())==91){
00101         break;
00102       }
00103       pos1stCluster++; //Starts at Zero, like the Collection
00104     }
00105     
00106     //Loop through the all particles in the event
00107     for ( HepMC::GenEvent::particle_const_iterator p = myGenEvent->particles_begin();p != myGenEvent->particles_end(); ++p ) {     
00108       //"Garbage" Cut, 1 GeV/c Pt Cut on All Particles Considered
00109       if((*p)->momentum().perp()>1.0){
00110         //Final State Quark criterion
00111         if(abs((*p)->pdg_id())==1 || abs((*p)->pdg_id())==2 || abs((*p)->pdg_id())==3 || abs((*p)->pdg_id())==4 || abs((*p)->pdg_id())==5){
00112           if( counter<pos1stCluster && ((*p)->status()==158 || (*p)->status()==159) ){
00113             isFSQuark=true;
00114           }
00115         }//end if FS Quark criterion
00116       }//End "Garbage" Cut
00117       
00118       if(isFSQuark){
00119         hFSPartons_JS_PtWgting->Fill( (*p)->momentum().eta(), (*p)->momentum().phi(), (*p)->momentum().perp()); //weighted by Particle Pt
00120       }
00121       
00122       counter++;
00123       isFSQuark=false;
00124     } //end all particles loop
00125     
00126     maxPartonPt=hFSPartons_JS_PtWgting->GetMaximum();
00127     
00128     //The Actual Filtering Process
00129     if(maxPartonPt>=minptcut && maxPartonPt<maxptcut){
00130       accepted=true; //Accept the Event
00131 
00132     }//End Filtering
00133   }//end processId if statement
00134   
00135   else{ accepted = true; }
00136  
00137   if (accepted){return true; } 
00138   else {return false;}
00139 }