Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Merge attempt #166

Merged
merged 3 commits into from
Jul 23, 2013
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
94 changes: 46 additions & 48 deletions Validation/EventGenerator/plugins/TauValidation.cc
Original file line number Diff line number Diff line change
Expand Up @@ -2,10 +2,7 @@
*
* Class to fill dqm monitor elements from existing EDM file
*
* $Date: 2013/03/06 01:51:11 $
* $Revision: 1.28 $
*/

#include "Validation/EventGenerator/interface/TauValidation.h"

#include "CLHEP/Units/defs.h"
Expand Down Expand Up @@ -159,8 +156,7 @@ void TauValidation::beginRun(const edm::Run& iRun,const edm::EventSetup& iSetup)
void TauValidation::endRun(const edm::Run& iRun,const edm::EventSetup& iSetup){return;}
void TauValidation::analyze(const edm::Event& iEvent,const edm::EventSetup& iSetup)
{

///Gathering the HepMCProduct information
///Gathering the HepMCProduct information
edm::Handle<HepMCProduct> evt;
iEvent.getByLabel(hepmcCollection_, evt);

Expand All @@ -176,11 +172,9 @@ void TauValidation::analyze(const edm::Event& iEvent,const edm::EventSetup& iSet
weight = 1.0;
if(*(WT.product())>1e-3 && *(WT.product())<=10.0) weight=(*(WT.product()));
else {weight=1.0;}
std::cout << "WT " << (*(WT.product())) << std::endl;
*/
///////////////////////////////////////////////


// find taus
for(HepMC::GenEvent::particle_const_iterator iter = myGenEvent->particles_begin(); iter != myGenEvent->particles_end(); iter++) {
if(abs((*iter)->pdg_id())==23){
Expand All @@ -205,55 +199,59 @@ void TauValidation::analyze(const edm::Event& iEvent,const edm::EventSetup& iSet
//
TauDecay_CMSSW TD;
unsigned int jak_id, TauBitMask;
TD.AnalyzeTau((*iter),jak_id,TauBitMask,false,false);
JAKID->Fill(jak_id,weight);
if(jak_id<=NJAKID){
int tcharge=(*iter)->pdg_id()/abs((*iter)->pdg_id());
std::vector<HepMC::GenParticle*> part=TD.Get_TauDecayProducts();
spinEffects(*iter,mother,jak_id,part,weight);
TLorentzVector LVQ(0,0,0,0);
TLorentzVector LVS12(0,0,0,0);
TLorentzVector LVS13(0,0,0,0);
TLorentzVector LVS23(0,0,0,0);
bool haspart1=false;
for(unsigned int i=0;i<part.size();i++){
if(TD.isTauFinalStateParticle(part.at(i)->pdg_id()) &&
abs(part.at(i)->pdg_id())!=PdtPdgMini::nu_e &&
abs(part.at(i)->pdg_id())!=PdtPdgMini::nu_mu &&
abs(part.at(i)->pdg_id())!=PdtPdgMini::nu_tau ){
TLorentzVector LV(part.at(i)->momentum().px(),part.at(i)->momentum().py(),part.at(i)->momentum().pz(),part.at(i)->momentum().e());
LVQ+=LV;
if(jak_id==TauDecay::JAK_A1_3PI ||
jak_id==TauDecay::JAK_KPIK ||
jak_id==TauDecay::JAK_KPIPI
){
if((tcharge==part.at(i)->pdg_id()/abs(part.at(i)->pdg_id()) && TD.nProng(TauBitMask)==3) || (jak_id==TauDecay::JAK_A1_3PI && TD.nProng(TauBitMask)==1 && abs(part.at(i)->pdg_id())==PdtPdgMini::pi_plus) ){
LVS13+=LV;
LVS23+=LV;
}
else{
LVS12+=LV;
if(!haspart1 && ((jak_id==TauDecay::JAK_A1_3PI) || (jak_id!=TauDecay::JAK_A1_3PI && abs(part.at(i)->pdg_id())==PdtPdgMini::K_plus) )){
if(TD.AnalyzeTau((*iter),jak_id,TauBitMask,false,false)){
JAKID->Fill(jak_id,weight);
if(jak_id<=NJAKID){
int tcharge=(*iter)->pdg_id()/abs((*iter)->pdg_id());
std::vector<HepMC::GenParticle*> part=TD.Get_TauDecayProducts();
spinEffects(*iter,mother,jak_id,part,weight);
TLorentzVector LVQ(0,0,0,0);
TLorentzVector LVS12(0,0,0,0);
TLorentzVector LVS13(0,0,0,0);
TLorentzVector LVS23(0,0,0,0);
bool haspart1=false;
for(unsigned int i=0;i<part.size();i++){
if(TD.isTauFinalStateParticle(part.at(i)->pdg_id()) &&
abs(part.at(i)->pdg_id())!=PdtPdgMini::nu_e &&
abs(part.at(i)->pdg_id())!=PdtPdgMini::nu_mu &&
abs(part.at(i)->pdg_id())!=PdtPdgMini::nu_tau ){
TLorentzVector LV(part.at(i)->momentum().px(),part.at(i)->momentum().py(),part.at(i)->momentum().pz(),part.at(i)->momentum().e());
LVQ+=LV;
if(jak_id==TauDecay::JAK_A1_3PI ||
jak_id==TauDecay::JAK_KPIK ||
jak_id==TauDecay::JAK_KPIPI
){
if((tcharge==part.at(i)->pdg_id()/abs(part.at(i)->pdg_id()) && TD.nProng(TauBitMask)==3) || (jak_id==TauDecay::JAK_A1_3PI && TD.nProng(TauBitMask)==1 && abs(part.at(i)->pdg_id())==PdtPdgMini::pi_plus) ){
LVS13+=LV;
haspart1=true;
LVS23+=LV;
}
else{
LVS23+=LV;
LVS12+=LV;
if(!haspart1 && ((jak_id==TauDecay::JAK_A1_3PI) || (jak_id!=TauDecay::JAK_A1_3PI && abs(part.at(i)->pdg_id())==PdtPdgMini::K_plus) )){
LVS13+=LV;
haspart1=true;
}
else{
LVS23+=LV;
}
}
}
}
}
part.clear();
JAKInvMass.at(jak_id).at(0)->Fill(LVQ.M(),weight);
if(jak_id==TauDecay::JAK_A1_3PI ||
jak_id==TauDecay::JAK_KPIK ||
jak_id==TauDecay::JAK_KPIPI
){
JAKInvMass.at(jak_id).at(1)->Fill(LVS13.M(),weight);
JAKInvMass.at(jak_id).at(2)->Fill(LVS23.M(),weight);
JAKInvMass.at(jak_id).at(3)->Fill(LVS12.M(),weight);
}
}
part.clear();
JAKInvMass.at(jak_id).at(0)->Fill(LVQ.M(),weight);
if(jak_id==TauDecay::JAK_A1_3PI ||
jak_id==TauDecay::JAK_KPIK ||
jak_id==TauDecay::JAK_KPIPI
){
JAKInvMass.at(jak_id).at(1)->Fill(LVS13.M(),weight);
JAKInvMass.at(jak_id).at(2)->Fill(LVS23.M(),weight);
JAKInvMass.at(jak_id).at(3)->Fill(LVS12.M(),weight);
}
}
else{
JAKID->Fill(jak_id,weight);
}
}
}
Expand Down
12 changes: 7 additions & 5 deletions Validation/EventGenerator/src/TauDecay_CMSSW.cc
Original file line number Diff line number Diff line change
Expand Up @@ -21,12 +21,14 @@ bool TauDecay_CMSSW::AnalyzeTau(HepMC::GenParticle *Tau,unsigned int &JAK_ID,uns
if(abs(Tau->pdg_id())==PdtPdgMini::tau_minus){ // check that it is a tau
unsigned int Tauidx=TauDecayProducts.size();
HepMC::GenVertex::particle_iterator des;
for(des = Tau->end_vertex()->particles_begin(HepMC::children);
des!= Tau->end_vertex()->particles_end(HepMC::children);++des ) {
Analyze((*des),Tauidx,dores,dopi0);
if( Tau->end_vertex()){
for(des = Tau->end_vertex()->particles_begin(HepMC::children);
des!= Tau->end_vertex()->particles_end(HepMC::children);++des ) {
Analyze((*des),Tauidx,dores,dopi0);
}
ClassifyDecayMode(JAK_ID,TauBitMask);
return true;
}
ClassifyDecayMode(JAK_ID,TauBitMask);
return true;
}
return false;
}
Expand Down