CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_1/src/Validation/EventGenerator/plugins/TauValidation.cc

Go to the documentation of this file.
00001 /*class TauValidation
00002  *  
00003  *  Class to fill dqm monitor elements from existing EDM file
00004  *
00005  *  $Date: 2011/12/29 10:53:11 $
00006  *  $Revision: 1.14 $
00007  */
00008  
00009 #include "Validation/EventGenerator/interface/TauValidation.h"
00010 
00011 #include "CLHEP/Units/defs.h"
00012 #include "CLHEP/Units/PhysicalConstants.h"
00013 
00014 #include "DataFormats/Math/interface/LorentzVector.h"
00015 
00016 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
00017 
00018 using namespace edm;
00019 
00020 TauValidation::TauValidation(const edm::ParameterSet& iPSet): 
00021   _wmanager(iPSet),
00022   hepmcCollection_(iPSet.getParameter<edm::InputTag>("hepmcCollection")),
00023   tauEtCut(iPSet.getParameter<double>("tauEtCutForRtau"))
00024 {    
00025   dbe = 0;
00026   dbe = edm::Service<DQMStore>().operator->();
00027 }
00028 
00029 TauValidation::~TauValidation() {}
00030 
00031 void TauValidation::beginJob()
00032 {
00033   if(dbe){
00035     dbe->setCurrentFolder("Generator/Tau");
00036     
00037     // Number of analyzed events
00038     nEvt = dbe->book1D("nEvt", "n analyzed Events", 1, 0., 1.);
00039 
00040     //Kinematics
00041     TauPt            = dbe->book1D("TauPt","Tau pT", 100 ,0,100);
00042     TauEta           = dbe->book1D("TauEta","Tau eta", 100 ,-2.5,2.5);
00043     TauPhi           = dbe->book1D("TauPhi","Tau phi", 100 ,-3.14,3.14);
00044     TauProngs        = dbe->book1D("TauProngs","Tau n prongs", 7 ,0,7);
00045     TauDecayChannels = dbe->book1D("TauDecayChannels","Tau decay channels", 12 ,0,12);
00046         TauDecayChannels->setBinLabel(1+undetermined,"?");
00047         TauDecayChannels->setBinLabel(1+electron,"e");
00048         TauDecayChannels->setBinLabel(1+muon,"mu");
00049         TauDecayChannels->setBinLabel(1+pi,"#pi^{#pm}");
00050         TauDecayChannels->setBinLabel(1+rho,"#rho^{#pm}");
00051         TauDecayChannels->setBinLabel(1+a1,"a_{1}^{#pm}");
00052         TauDecayChannels->setBinLabel(1+pi1pi0,"#pi^{#pm}#pi^{0}");
00053         TauDecayChannels->setBinLabel(1+pinpi0,"#pi^{#pm}n#pi^{0}");
00054         TauDecayChannels->setBinLabel(1+tripi,"3#pi^{#pm}");
00055 //      TauDecayChannels->setBinLabel(1+tripinpi0,"3#pi^{#pm}n#pi^{0}");
00056         TauDecayChannels->setBinLabel(1+K,"K");
00057         TauDecayChannels->setBinLabel(1+Kstar,"K^{*}");
00058         TauDecayChannels->setBinLabel(1+stable,"Stable");
00059 
00060     TauMothers        = dbe->book1D("TauMothers","Tau mother particles", 10 ,0,10);
00061         TauMothers->setBinLabel(1+other,"?");
00062         TauMothers->setBinLabel(1+gamma,"#gamma");
00063         TauMothers->setBinLabel(1+Z,"Z");
00064         TauMothers->setBinLabel(1+W,"W");
00065         TauMothers->setBinLabel(1+HSM,"H_{SM}/h^{0}");
00066         TauMothers->setBinLabel(1+H0,"H^{0}");
00067         TauMothers->setBinLabel(1+A0,"A^{0}");
00068         TauMothers->setBinLabel(1+Hpm,"H^{#pm}");
00069 
00070     TauRtauW          = dbe->book1D("TauRtauW","W->Tau p(leading track)/E(visible tau)", 50 ,0,1);
00071         TauRtauW->setAxisTitle("rtau");
00072     TauRtauHpm        = dbe->book1D("TauRtauHpm","Hpm->Tau p(leading track)/E(visible tau)", 50 ,0,1);
00073         TauRtauHpm->setAxisTitle("rtau");
00074     TauSpinEffectsW   = dbe->book1D("TauSpinEffectsW","Pion energy in W rest frame", 50 ,0,1);
00075         TauSpinEffectsW->setAxisTitle("Energy");
00076     TauSpinEffectsHpm = dbe->book1D("TauSpinEffectsHpm","Pion energy in Hpm rest frame", 50 ,0,1);
00077         TauSpinEffectsHpm->setAxisTitle("Energy");
00078     TauSpinEffectsZ   = dbe->book1D("TauSpinEffectsZ","Mass of pi+ pi-", 22 ,0,1.1);
00079         TauSpinEffectsZ->setAxisTitle("M_{#pi^{+}#pi^{-}}");
00080 
00081     TauPhotonsN        = dbe->book1D("TauPhotonsN","Photons radiating from tau", 2 ,0,2);
00082         TauPhotonsN->setBinLabel(1,"Number of taus");
00083         TauPhotonsN->setBinLabel(2,"Number of taus radiating photons");
00084     TauPhotonsPt       = dbe->book1D("TauPhotonsPt","Photon pt radiating from tau", 2 ,0,2);
00085         TauPhotonsPt->setBinLabel(1,"Sum of tau pt");
00086         TauPhotonsPt->setBinLabel(2,"Sum of tau pt radiated by photons");
00087   }
00088 
00089   return;
00090 }
00091 
00092 void TauValidation::endJob(){
00093   return;
00094 }
00095 
00096 void TauValidation::beginRun(const edm::Run& iRun,const edm::EventSetup& iSetup)
00097 {
00099   iSetup.getData( fPDGTable );
00100   return;
00101 }
00102 void TauValidation::endRun(const edm::Run& iRun,const edm::EventSetup& iSetup){return;}
00103 void TauValidation::analyze(const edm::Event& iEvent,const edm::EventSetup& iSetup)
00104 { 
00105   
00107   edm::Handle<HepMCProduct> evt;
00108   iEvent.getByLabel(hepmcCollection_, evt);
00109 
00110   //Get EVENT
00111   HepMC::GenEvent *myGenEvent = new HepMC::GenEvent(*(evt->GetEvent()));
00112 
00113   double weight = _wmanager.weight(iEvent);
00114 
00115 
00116   nEvt->Fill(0.5,weight);
00117 
00118   // find taus
00119   for(HepMC::GenEvent::particle_const_iterator iter = myGenEvent->particles_begin(); iter != myGenEvent->particles_end(); ++iter) {
00120     if ((*iter)->status()==3)
00121     {
00122       if(abs((*iter)->pdg_id())==15){
00123         TauPt->Fill((*iter)->momentum().perp(),weight);
00124         TauEta->Fill((*iter)->momentum().eta(),weight);
00125         TauPhi->Fill((*iter)->momentum().phi(),weight);
00126         int mother  = tauMother(*iter, weight);
00127         int decaychannel = tauDecayChannel(*iter, weight);
00128         tauProngs(*iter, weight);
00129         rtau(*iter,mother,decaychannel, weight);
00130         spinEffects(*iter,mother,decaychannel, weight);
00131         photons(*iter, weight);
00132       }
00133       if(abs((*iter)->pdg_id())==23){
00134         spinEffectsZ(*iter, weight);
00135       }
00136     }
00137   }
00138 
00139   delete myGenEvent;
00140 }//analyze
00141 
00142 int TauValidation::findMother(const HepMC::GenParticle* tau){
00143         int mother_pid = 0;
00144 
00145         if ( tau->production_vertex() ) {
00146                 HepMC::GenVertex::particle_iterator mother;
00147                 for (mother = tau->production_vertex()->particles_begin(HepMC::parents);
00148                      mother!= tau->production_vertex()->particles_end(HepMC::parents); ++mother ) {
00149                         mother_pid = (*mother)->pdg_id();
00150                         if(mother_pid == tau->pdg_id()) mother_pid = -1;//findMother(*mother);
00151                         //std::cout << " parent " << mother_pid << std::endl;
00152                 }
00153         }
00154         return mother_pid;
00155 }
00156 
00157 int TauValidation::tauMother(const HepMC::GenParticle* tau, double weight){
00158 
00159         if(abs(tau->pdg_id()) != 15 ) return -1;
00160 
00161         int mother_pid = findMother(tau);
00162         if(mother_pid == -1) return -1;
00163 
00164         int label = other;
00165         if(abs(mother_pid) == 24) label = W;
00166         if(abs(mother_pid) == 23) label = Z;
00167         if(abs(mother_pid) == 22) label = gamma;
00168         if(abs(mother_pid) == 25) label = HSM;
00169         if(abs(mother_pid) == 35) label = H0;
00170         if(abs(mother_pid) == 36) label = A0;
00171         if(abs(mother_pid) == 37) label = Hpm;
00172 
00173         TauMothers->Fill(label,weight);
00174 
00175         return mother_pid;
00176 }
00177 
00178 int TauValidation::tauProngs(const HepMC::GenParticle* tau, double weight){
00179 
00180         int nProngs = 0;
00181         if ( tau->end_vertex() ) {
00182                 HepMC::GenVertex::particle_iterator des;
00183                 for(des = tau->end_vertex()->particles_begin(HepMC::descendants);
00184                     des!= tau->end_vertex()->particles_end(HepMC::descendants);++des ) {
00185                         int pid = (*des)->pdg_id();
00186                         if(abs(pid) == 15) return tauProngs(*des, weight);
00187                         if((*des)->status() != 1) continue; // dont count unstable particles
00188 
00189                         const HepPDT::ParticleData*  pd = fPDGTable->particle((*des)->pdg_id ());
00190                         int charge = (int) pd->charge();
00191                         if(charge == 0) continue;
00192                         //std::cout << "TauValidation::tauProngs barcode=" << (*des)->barcode() << " pid=" 
00193                         //          << pid << " mom=" << tauMother(*des) << " status=" 
00194                         //          << (*des)->status() << " charge=" << charge << std::endl;
00195                         nProngs++;
00196                 }
00197         }
00198         TauProngs->Fill(nProngs,weight);
00199         return nProngs;
00200 }
00201 
00202 int TauValidation::findTauDecayChannel(const HepMC::GenParticle* tau){
00203 
00204         int channel = undetermined;
00205         if(tau->status() == 1) channel = stable;
00206 
00207         int allCount   = 0,
00208             eCount     = 0,
00209             muCount    = 0,
00210             pi0Count   = 0,
00211             piCount    = 0,
00212             rhoCount   = 0,
00213             a1Count    = 0,
00214             KCount     = 0,
00215             KstarCount = 0;
00216 
00217         if ( tau->end_vertex() ) {
00218               HepMC::GenVertex::particle_iterator des;
00219               for(des = tau->end_vertex()->particles_begin(HepMC::descendants);
00220                   des!= tau->end_vertex()->particles_end(HepMC::descendants);++des ) {
00221 
00222                         //if(abs(tauMother(*des)) != 15) continue;
00223                         int pid = (*des)->pdg_id();
00224                         //std::cout << " barcode=" << (*des)->barcode() << " pid="
00225                         //          << pid << " mom=" << tauMother(*des) << " status="
00226                         //          << (*des)->status() << std::endl;
00227 
00228                         if(abs(pid) == 15) return findTauDecayChannel(*des);
00229 
00230                         allCount++;
00231                         if(abs(pid) == 11)    eCount++;
00232                         if(abs(pid) == 13)    muCount++;
00233                         if(abs(pid) == 111)   pi0Count++;
00234                         if(abs(pid) == 211)   piCount++;
00235                         if(abs(pid) == 213)   rhoCount++;
00236                         if(abs(pid) == 20213) a1Count++;
00237                         if(abs(pid) == 321)   KCount++;
00238                         if(abs(pid) == 323)   KstarCount++;
00239 
00240                 }
00241         }
00242 
00243         if(KCount == 1 && allCount == 2)  channel = K;
00244         if(KstarCount == 1 && allCount == 2)  channel = Kstar;
00245         if(a1Count == 1 && allCount == 2)  channel = a1;
00246         if(rhoCount == 1 && allCount == 2)  channel = rho;
00247 
00248         if(piCount == 1 && pi0Count == 0) channel = pi;
00249         if(piCount == 1 && pi0Count == 1) channel = pi1pi0;
00250         if(piCount == 1 && pi0Count > 1)  channel = pinpi0;
00251 
00252         if(piCount == 3 && pi0Count == 0) channel = tripi;
00253 //        if(piCount == 3 && pi0Count > 0)  channel = tripinpi0;
00254 
00255         if(eCount == 1)                   channel = electron;
00256         if(muCount == 1)                  channel = muon;
00257 
00258         return channel;
00259 }
00260 
00261 int TauValidation::tauDecayChannel(const HepMC::GenParticle* tau, double weight){
00262         int channel = findTauDecayChannel(tau);
00263         TauDecayChannels->Fill(channel,weight);
00264         return channel;
00265 }
00266 
00267 void TauValidation::rtau(const HepMC::GenParticle* tau,int mother, int decay, double weight){
00268 
00269         if(decay != pi1pi0) return; // polarization only for 1-prong hadronic taus with one neutral pion to make a clean case
00270 
00271         if(tau->momentum().perp() < tauEtCut) return; // rtau visible only for boosted taus
00272         
00273         double rTau = 0;
00274         double ltrack = leadingPionMomentum(tau, weight);
00275         double visibleTauE = visibleTauEnergy(tau);
00276 
00277         if(visibleTauE != 0) rTau = ltrack/visibleTauE;
00278 
00279         if(abs(mother) == 24) TauRtauW->Fill(rTau,weight);
00280         if(abs(mother) == 37) TauRtauHpm->Fill(rTau,weight); 
00281 }
00282 
00283 void TauValidation::spinEffects(const HepMC::GenParticle* tau,int mother, int decay, double weight){
00284 
00285         if(decay != pi) return; // polarization only for 1-prong hadronic taus with no neutral pions
00286 
00287         TLorentzVector momP4 = motherP4(tau);
00288         TLorentzVector pionP4 = leadingPionP4(tau);
00289 
00290         pionP4.Boost(-1*momP4.BoostVector());
00291 
00292         double energy = pionP4.E()/(momP4.M()/2);
00293 
00294         if(abs(mother) == 24) TauSpinEffectsW->Fill(energy,weight);     
00295         if(abs(mother) == 37) TauSpinEffectsHpm->Fill(energy,weight);
00296 }
00297 
00298 void TauValidation::spinEffectsZ(const HepMC::GenParticle* boson, double weight){
00299 
00300         TLorentzVector tautau(0,0,0,0);
00301         TLorentzVector pipi(0,0,0,0);
00302 
00303         int nSinglePionDecays = 0;
00304         if ( boson->end_vertex() ) {
00305               HepMC::GenVertex::particle_iterator des;
00306               for(des = boson->end_vertex()->particles_begin(HepMC::descendants);
00307                   des!= boson->end_vertex()->particles_end(HepMC::descendants);++des ) {
00308 
00309                         int pid = (*des)->pdg_id();
00310                         /*std::cout << " barcode=" << (*des)->barcode() << " pid="
00311                                   << pid << " mom=" << findMother(*des) << " status="
00312                                   << (*des)->status() << " px="
00313                                   << (*des)->momentum().px() << " decay=" 
00314                                   << findTauDecayChannel(*des) << std::endl;
00315                         */
00316                         if(abs(findMother(*des)) != 15 &&
00317                            abs(pid) == 15 && findTauDecayChannel(*des) == pi){
00318                           nSinglePionDecays++;
00319                           tautau += TLorentzVector((*des)->momentum().px(),
00320                                                    (*des)->momentum().py(),
00321                                                    (*des)->momentum().pz(),
00322                                                    (*des)->momentum().e());
00323                           pipi += leadingPionP4(*des);
00324                         }
00325                 }
00326         }
00327         if(nSinglePionDecays == 2 && tautau.M() != 0) {
00328           TauSpinEffectsZ->Fill(pipi.M()/tautau.M(),weight);
00329         }
00330 }
00331 
00332 double TauValidation::leadingPionMomentum(const HepMC::GenParticle* tau, double weight){
00333         return leadingPionP4(tau).P();
00334 }
00335 
00336 TLorentzVector TauValidation::leadingPionP4(const HepMC::GenParticle* tau){
00337 
00338         TLorentzVector p4(0,0,0,0);
00339 
00340         if ( tau->end_vertex() ) {
00341               HepMC::GenVertex::particle_iterator des;
00342               for(des = tau->end_vertex()->particles_begin(HepMC::descendants);
00343                   des!= tau->end_vertex()->particles_end(HepMC::descendants);++des ) {
00344                         int pid = (*des)->pdg_id();
00345 
00346                         if(abs(pid) == 15) return leadingPionP4(*des);
00347 
00348                         if(abs(pid) != 211) continue;
00349  
00350                         if((*des)->momentum().rho() > p4.P()) {
00351                                 p4 = TLorentzVector((*des)->momentum().px(),
00352                                                     (*des)->momentum().py(),
00353                                                     (*des)->momentum().pz(),
00354                                                     (*des)->momentum().e());
00355                         } 
00356                 }
00357         }
00358 
00359         return p4;
00360 }
00361 
00362 TLorentzVector TauValidation::motherP4(const HepMC::GenParticle* tau){
00363 
00364         TLorentzVector p4(0,0,0,0);
00365 
00366         if ( tau->production_vertex() ) {
00367                 HepMC::GenVertex::particle_iterator mother;
00368                 for (mother = tau->production_vertex()->particles_begin(HepMC::parents);
00369                      mother!= tau->production_vertex()->particles_end(HepMC::parents); ++mother ) {
00370                         //mother_pid = (*mother)->pdg_id();
00371                         //std::cout << " parent " << mother_pid << std::endl;
00372                         p4 = TLorentzVector((*mother)->momentum().px(),
00373                                             (*mother)->momentum().py(),
00374                                             (*mother)->momentum().pz(),
00375                                             (*mother)->momentum().e());
00376                 }
00377         }
00378 
00379         return p4;
00380 }
00381 
00382 double TauValidation::visibleTauEnergy(const HepMC::GenParticle* tau){
00383         TLorentzVector p4(tau->momentum().px(),
00384                           tau->momentum().py(),
00385                           tau->momentum().pz(),
00386                           tau->momentum().e());
00387 
00388         if ( tau->end_vertex() ) {
00389               HepMC::GenVertex::particle_iterator des;
00390               for(des = tau->end_vertex()->particles_begin(HepMC::descendants);
00391                   des!= tau->end_vertex()->particles_end(HepMC::descendants);++des ) {
00392                         int pid = (*des)->pdg_id();
00393 
00394                         if(abs(pid) == 15) return visibleTauEnergy(*des);
00395 
00396                         if(abs(pid) == 12 || abs(pid) == 14 || abs(pid) == 16) {
00397                                 p4 -= TLorentzVector((*des)->momentum().px(),
00398                                                      (*des)->momentum().py(),
00399                                                      (*des)->momentum().pz(),
00400                                                      (*des)->momentum().e());
00401                         }
00402                 }
00403         }
00404 
00405         return p4.E();
00406 }
00407 
00408 void TauValidation::photons(const HepMC::GenParticle* tau, double weight){
00409 
00410         if ( tau->end_vertex() ) {
00411               double photonFromTauPtSum = 0;
00412               bool photonFromTau = false;
00413               HepMC::GenVertex::particle_iterator des;
00414               for(des = tau->end_vertex()->particles_begin(HepMC::descendants);
00415                   des!= tau->end_vertex()->particles_end(HepMC::descendants);++des ) {
00416                         int pid = (*des)->pdg_id();
00417                         if(pid == 22) {
00418                                 photonFromTauPtSum += (*des)->momentum().perp();
00419                                 photonFromTau = true;
00420                         } 
00421               }
00422               
00423               TauPhotonsN->Fill(0.5,weight);
00424               //doesn't seems like it makes sense to use a weight below  
00425               TauPhotonsPt->Fill(0.5,tau->momentum().perp());
00426               if(photonFromTau) {
00427                 TauPhotonsN->Fill(1.5,weight);
00428                 //doesn't seems like it makes sense to use a weight below  
00429                 TauPhotonsPt->Fill(1.5,photonFromTauPtSum);
00430               }
00431         }
00432 }
00433