Main Page | Namespace List | Class Hierarchy | Class List | Directories | File List | Class Members | File Members

TauAnalysis.cc

Go to the documentation of this file.
00001 #include <iostream>
00002 #include <sstream>
00003 #include <stdio.h>
00004 #include "HepMC/GenEvent.h"
00005 #include "HepMC/IO_GenEvent.h"
00006 #include "HepMC/GenParticle.h"
00007 #include "HepMC/GenVertex.h"
00008 #include "HepMC/IO_AsciiParticles.h"
00009 #include "HepMC/SimpleVector.h"
00010 #include "CLHEP/Vector/LorentzVector.h"
00011 
00012 using namespace std;
00013 
00014 // ROOT headers
00015 #include "TH1.h"
00016 #include "TH2.h"
00017 #include "TFile.h"
00018 #include "TMath.h"
00019 #include "TLorentzVector.h"
00020 
00021 //Tau analysis header
00022 #include "../include/TauAnalysis.h"
00023 
00024 //**********************
00025 
00026 /// empty default constructor
00027 TauAnalysis::TauAnalysis()
00028 {
00029 }
00030 
00031 /// empty default destructor
00032 TauAnalysis::~TauAnalysis()
00033 {
00034 }
00035 
00036 int TauAnalysis::Init(double tr_max_eta, double tr_min_pt)
00037 {
00038   // specify default eta cut                                                                                                           
00039   m_max_eta=tr_max_eta;
00040 
00041   // specify default pt cut
00042   m_min_pt=tr_min_pt;
00043     
00044   // default Output file name
00045   m_outputFileName="TauAnalysis.root";
00046   m_outputRootDir="Tau";
00047 
00048   //declaration of histograms
00049   m_evtnr=initHist(string("evtnr"),string("Event number"),string("Eventnumber"),100, 0., 99.);
00050   m_taupt=initHist(string("taupt"),string("transversal momentum of tau+ and tau-"),string("p_{T} [GeV]"),100, 0., 200.);
00051   m_taupt_log=initHist(string("taupt_logscale"),string("transversal momentum of tau+ and tau- - logscale -"),
00052                        string("p_{T} [GeV]"),100, 0., 200.);
00053   m_taueta=initHist(string("taueta"),string("eta of tau+ and tau-"),string("#eta"),60, -6., 6.);
00054   m_tauphi=initHist(string("tauphi"),string("phi of tau+ and tau-"),string("#phi"),48, -3.15, 3.15);
00055   m_ptstable=initHist(string("ptstable"),string("transveral momentum of stable particles (without neutrinos)"),
00056                       string("p_{T} [GeV]"),60, 0., 15.);
00057   m_ptstable_log=initHist(string("ptstable_logscale"),string("transveral momentum of stable particles (without neutrinos) - logscale"),
00058                           string("p_{T} [GeV]"),60, 0., 15.);
00059   m_etastable=initHist(string("etastable"),string("eta of stable particles (without neutrinos)"),string("#eta"),100, -6., 6.);
00060   m_phistable=initHist(string("phistable"),string("phi of stable particles (without neutrinos)"),string("#phi"),100, -3.15, 3.15);
00061   m_E_stable=initHist(string("E_stable"),string("energy of stable particles (without neutrinos)"),string("E [GeV]"),50, 0., 100.);
00062   m_E_stable_log=initHist(string("E_stable_logscale"),string("energy of stable particles (without neutrinos) - logscale -"),
00063                           string("E [GeV]"),50, 0., 100.);
00064   m_ptstable_charged=initHist(string("ptstable_charged"),string("transveral momentum of charged stable particles"),
00065                               string("p_{T} [GeV]"),200, 0., 50.);
00066   m_ptstable_charged_log=initHist(string("ptstable_charged_logscale"),
00067                                   string("transveral momentum of charged stable particles - logscale -"),
00068                                   string("p_{T} [GeV]"),200, 0., 50.);
00069   m_etastable_charged=initHist(string("etastable_charged"),string("eta of charged stable particles"),string("#eta"),100, -6., 6.);
00070   m_phistable_charged=initHist(string("phistable_charged"),string("phi of charged stable particles"),string("#phi"),100, -3.15, 3.15);
00071   m_E_stable_charged=initHist(string("E_stable_without_neutrinos"),string("energy of charged stable particles"),
00072                               string("E [GeV]"),100, 0., 200.);
00073   m_E_stable_charged_log=initHist(string("E_stable_without_neutrinos_logscale"),
00074                                   string("energy of charged stable particles - logscale -"),string("E [GeV]"),100, 0., 200.);
00075   m_cmultpart=initHist(string("cmultpart"),string("number of charged particles in the event"),string("charged particle multiplicity"),
00076                        120, 0., 600.);
00077   m_nTrack_tau=initHist(string("nTrack_tau"),string("number of tracks of tau decay"),string("number of tracks"),16, -0.5, 15.5);
00078   m_pt_hightrack=initHist(string("pt_hightrack"),string("transversal momentum of highest p_{T}-track"),string("p_{T} [GeV]"),
00079                           200, 0., 100.);
00080   m_pt_hightrack_log=initHist(string("pt_hightrack_logscale"),string("transversal momentum of highest p_{T}-track - logscale -"),
00081                               string("p_{T} [GeV]"),200, 0., 100.);
00082   
00083   m_jet_count=baseAnalysis::initHist(string("jet_count"), string("Nr of Jets"), string("Nr Jets"),16, -0.5, 15.5);
00084   m_jet_pt=baseAnalysis::initHist(string("leadingjet_pt_high"),string("Leading Jet Transverse Momentum"),string("Jet p_{T} [GeV]"), 200, 0.,1000.);
00085   m_jet_pt_log=baseAnalysis::initHist(string("leadingjet_pt_high_logscale"),string("Leading Jet Transverse Momentum -logscale-"),string("Jet p_{T} [GeV]"), 200, 0.,1000.);
00086   
00087   return true;
00088 }
00089 
00090 int TauAnalysis::Process(HepMC::GenEvent* hepmcevt) 
00091 {
00092   
00093   int cmult = 0; //charged particle multiplicity
00094   int tau_barcode = 0;
00095   int taubar_barcode = 0;
00096   
00097   //variable for status of particle
00098   int properStatus = 0;
00099   
00100   m_evtnr->Fill(hepmcevt->event_number());
00101 
00102   // loop over the particles and select pdgid and pt
00103   for ( HepMC::GenEvent::particle_const_iterator p = 
00104           hepmcevt->particles_begin(); p != hepmcevt->particles_end(); ++p )
00105     {
00106       // use the pdg id to identify the particles
00107       int pid = (*p)->pdg_id();
00108       
00109       if (pid == 23)
00110         {
00111           if (!((*p)->end_vertex()))
00112             {
00113               std::cout<<"CRASH!!!!!"<<std::endl;
00114               std::cout<<"event number: "<<hepmcevt->event_number()<<std::endl;
00115               std::cout<<"pdgId (of Z): "<<pid<<std::endl;
00116               std::cout<<"particle number (Z barcode): "<<(*p)->barcode()<<std::endl;
00117               std::cout<<"(*p)->status(): "<<(*p)->status()<<std::endl;
00118               std::cout<<"(*p)->end_vertex(): "<<(*p)->end_vertex()<<std::endl;
00119               
00120               std::cout<<"skip event"<<std::endl;
00121               return false;
00122             }
00123       
00124           // find the last taus in the event and 
00125           // store their barcodes  to find them back in the event
00126           for(HepMC::GenVertex::particle_iterator iChild=(*p)->end_vertex()->particles_begin(HepMC::descendants); 
00127               iChild!=(*p)->end_vertex()->particles_end(HepMC::descendants); iChild++)
00128             {
00129               if( (*iChild)->pdg_id() != pid )
00130                 {
00131                   if ((*iChild)->pdg_id() == 15) tau_barcode = (*iChild)->barcode();
00132                   if ((*iChild)->pdg_id() == -15)  taubar_barcode = (*iChild)->barcode();
00133                 }
00134             }
00135         }
00136       
00137       // take all stable particles of this event
00138      if(!IsFinalStateParticle((*p))) continue;
00139       
00140       if(! ( abs(pid)==12 || abs(pid)==14 || abs(pid)==16 ) ) //energy stable particles without neutrinos
00141         {
00142           m_ptstable->Fill( (*p)->momentum().perp() );
00143           m_ptstable_log->Fill( (*p)->momentum().perp() );
00144           m_etastable->Fill( (*p)->momentum().eta() );
00145           m_phistable->Fill( (*p)->momentum().phi() );
00146           m_E_stable->Fill( (*p)->momentum().e() );
00147           m_E_stable_log->Fill( (*p)->momentum().e() );
00148         }
00149       
00150       // cut out the neutral particles  (charge is not stored in HepMC)
00151       // just remove the neutrinos, gammas and neutrons from 
00152       // the final state
00153       if(!chargedParticle(pid)) continue;
00154       
00155       m_ptstable_charged->Fill( (*p)->momentum().perp() );
00156       m_ptstable_charged_log->Fill( (*p)->momentum().perp() );
00157       m_etastable_charged->Fill( (*p)->momentum().eta() );
00158       m_phistable_charged->Fill( (*p)->momentum().phi() );
00159       m_E_stable_charged->Fill( (*p)->momentum().e() );
00160       m_E_stable_charged_log->Fill( (*p)->momentum().e() );
00161       
00162       cmult++;
00163     }  // end of loop over particles
00164   
00165   // storage of input particles for jet
00166   CLHEP::HepLorentzVector lv_leadingJet(0,0,0,0);
00167   
00168   // check if jets were found
00169   //if(!m_inclusive_jets.size()) std::cout<<"TauAnalysis::no jets found"<<std::endl;
00170 
00171   //Fill the histograms
00172   m_jet_count->Fill(m_inclusive_jets.size());
00173 
00174   if(m_inclusive_jets.size()){
00175     m_jet_pt->Fill(m_inclusive_jets[0].perp());
00176     m_jet_pt_log->Fill(m_inclusive_jets[0].perp());
00177   }
00178   
00179   // now we have the taupair - let's histogramm pt, eta and the pt sum
00180   HepMC::GenParticle *tau = hepmcevt->barcode_to_particle(tau_barcode);
00181   HepMC::GenParticle *taubar = hepmcevt->barcode_to_particle(taubar_barcode);
00182   if(tau && taubar) 
00183     {
00184       m_taupt->Fill( tau->momentum().perp() );
00185       m_taupt->Fill( taubar->momentum().perp() );
00186       m_taupt_log->Fill( tau->momentum().perp() );
00187       m_taupt_log->Fill( taubar->momentum().perp() );
00188       m_taueta->Fill( tau->momentum().eta() );
00189       m_taueta->Fill( taubar->momentum().eta() );
00190       m_tauphi->Fill( tau->momentum().phi() );
00191       m_tauphi->Fill( taubar->momentum().phi() );
00192       
00193       m_cmultpart->Fill(cmult);
00194     }
00195   
00196   //checking for children of tau, needed for nTrack and the highest pT track
00197   
00198   //tau-
00199   properStatus = 2;
00200   
00201   if (tau && ( tau->status() == properStatus ) )
00202     {
00203       int nTrack = 0; //number of tracks in tau decay
00204       double highpT_Track = 0.; //track with hightest pT coming from tau decay
00205       
00206       //important for the tracks (in the following the commented lines apply to the first generation of tau decay:
00207       //if you look only for the first generation as descendant (child), than the child can be stable or not stable,
00208       //only if you look for more than one generation as descendants, you need the proper status on stable;
00209       //this means that the final decay product of the tau has to be stable
00210       
00211       //only first generation
00212       //for(HepMC::GenVertex::particle_iterator iChild=tau->end_vertex()->particles_begin(HepMC::children); 
00213       //  iChild!=tau->end_vertex()->particles_end(HepMC::children); iChild++)
00214       
00215       //all generations as descendants
00216       for(HepMC::GenVertex::particle_iterator iChild=tau->end_vertex()->particles_begin(HepMC::descendants); 
00217           iChild!=tau->end_vertex()->particles_end(HepMC::descendants); iChild++)
00218         {
00219           if (!tau->end_vertex())
00220             {
00221               std::cout<<"CRASH!!!!!"<<std::endl;
00222               std::cout<<"event number: "<<hepmcevt->event_number()<<std::endl;
00223               std::cout<<"particle number (tau barcode): "<<tau->barcode()<<std::endl;
00224               std::cout<<"tau->status(): "<<tau->status()<<std::endl;
00225               std::cout<<"tau->end_vertex(): "<<tau->end_vertex()<<std::endl;
00226               
00227               std::cout<<"skip event"<<std::endl;
00228               return false;
00229             }
00230           
00231           if( (*iChild)->pdg_id() != tau->pdg_id() )
00232             {
00233               int iChild_pid = (*iChild)->pdg_id();
00234               properStatus = 1;
00235               
00236               if( (*iChild)->status() == properStatus && chargedParticle(iChild_pid)==true)
00237                 {
00238                   nTrack = nTrack + 1;
00239                   if ( (*iChild)->momentum().perp() > highpT_Track )
00240                     {
00241                       highpT_Track = (*iChild)->momentum().perp();
00242                     }
00243                 }
00244             }
00245         }
00246 
00247       //only fill nTrack and pt of highest track if a tau decays
00248       if (nTrack > 0) 
00249         {
00250           m_nTrack_tau->Fill(nTrack);
00251           m_pt_hightrack->Fill(highpT_Track);
00252           m_pt_hightrack_log->Fill(highpT_Track);
00253         }
00254     }
00255   
00256   //tau+ (tbar)
00257   properStatus = 2;
00258   
00259   if (taubar && ( taubar->status() == properStatus ) )
00260     {
00261       int nTrack = 0; //number of tracks in tau decay
00262       double highpT_Track = 0.; //track with hightest pT coming from tau decay
00263       
00264       //important for the tracks (in the following the commented lines apply to the first generation of tau decay:
00265       //if you look only for the first generation as descendant (child), than the child can be stable or not stable,
00266       //only if you look for more than one generation as descendants, you need the proper status on stable;
00267       //this means that the final decay product of the tau has to be stable
00268       
00269       //only first generation
00270       //for(HepMC::GenVertex::particle_iterator iChild=taubar->end_vertex()->particles_begin(HepMC::children); 
00271       //  iChild!=taubar->end_vertex()->particles_end(HepMC::children); iChild++)
00272       
00273       if (!taubar->end_vertex())
00274         {
00275           std::cout<<"CRASH!!!!!"<<std::endl;
00276           std::cout<<"event number: "<<hepmcevt->event_number()<<std::endl;
00277           std::cout<<"particle number (tau barcode): "<<taubar->barcode()<<std::endl;
00278           std::cout<<"taubar->status(): "<<taubar->status()<<std::endl;
00279           std::cout<<"taubar->end_vertex(): "<<taubar->end_vertex()<<std::endl;
00280           
00281           std::cout<<"skip event"<<std::endl;
00282           return false;
00283         }
00284 
00285       //all generations as descendants
00286       for(HepMC::GenVertex::particle_iterator iChild=taubar->end_vertex()->particles_begin(HepMC::descendants); 
00287           iChild!=taubar->end_vertex()->particles_end(HepMC::descendants); iChild++)
00288         {
00289           if( (*iChild)->pdg_id() != taubar->pdg_id() )
00290             {
00291               int iChild_pid = (*iChild)->pdg_id();
00292               properStatus = 1;
00293               
00294               if( (*iChild)->status() == properStatus && chargedParticle(iChild_pid)==true)
00295                 {
00296                   nTrack = nTrack + 1;
00297                   if ( (*iChild)->momentum().perp() > highpT_Track )
00298                     {
00299                       highpT_Track = (*iChild)->momentum().perp();
00300                     }
00301                 }
00302             }
00303           
00304         }
00305       //only fill nTrack and pt of highest track if a tau decays
00306       if (nTrack > 0) 
00307         {
00308           m_nTrack_tau->Fill(nTrack);
00309           m_pt_hightrack->Fill(highpT_Track);
00310           m_pt_hightrack_log->Fill(highpT_Track);
00311         }
00312     }
00313   return true;
00314 } //end of Process

Generated on Thu Jul 23 14:57:36 2009 for HepMCAnalysis by  doxygen 1.3.9.1