CMS 3D CMS Logo

Public Member Functions | Private Attributes

WValidation Class Reference

#include <WValidation.h>

Inheritance diagram for WValidation:
edm::EDAnalyzer

List of all members.

Public Member Functions

virtual void analyze (const edm::Event &, const edm::EventSetup &)
virtual void beginJob ()
virtual void beginRun (const edm::Run &, const edm::EventSetup &)
virtual void endJob ()
virtual void endRun (const edm::Run &, const edm::EventSetup &)
 WValidation (const edm::ParameterSet &)
virtual ~WValidation ()

Private Attributes

int _flavor
 decay flavor
std::string _name
 decay flavor name
MonitorElementcos_theta_gamma_lepton
DQMStoredbe
 ME's "container".
edm::ESHandle
< HepPDT::ParticleDataTable
fPDGTable
 PDT table.
MonitorElementgamma_energy
edm::InputTag hepmcCollection_
MonitorElementlepeta
MonitorElementlepmet_mT
MonitorElementlepmet_mTPeak
MonitorElementlepmet_pt
MonitorElementlepmet_ptLog
MonitorElementlepmet_rap
MonitorElementleppt
MonitorElementmet
MonitorElementnEvt
MonitorElementWdaughters
MonitorElementWmass
MonitorElementWmassPeak
MonitorElementWpt
MonitorElementWptLog
MonitorElementWrap

Detailed Description

Definition at line 34 of file WValidation.h.


Constructor & Destructor Documentation

WValidation::WValidation ( const edm::ParameterSet iPSet) [explicit]

Definition at line 20 of file WValidation.cc.

References dbe, and cmsCodeRules::cppFunctionSkipper::operator.

                                                    :  
  hepmcCollection_(iPSet.getParameter<edm::InputTag>("hepmcCollection")),
  _flavor(iPSet.getParameter<int>("decaysTo")),
  _name(iPSet.getParameter<std::string>("name")) 
{    
  dbe = 0;
  dbe = edm::Service<DQMStore>().operator->();
}
WValidation::~WValidation ( ) [virtual]

Definition at line 29 of file WValidation.cc.

{}

Member Function Documentation

void WValidation::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup 
) [virtual]

Gathering the HepMCProduct information

Implements edm::EDAnalyzer.

Definition at line 74 of file WValidation.cc.

References _flavor, abs, funct::cos(), cos_theta_gamma_lepton, MonitorElement::Fill(), HepMCValidationHelper::findFSRPhotons(), gamma_energy, HepMCValidationHelper::genMet(), edm::Event::getByLabel(), hepmcCollection_, lepeta, lepmet_mT, lepmet_mTPeak, lepmet_pt, lepmet_ptLog, leppt, met, nEvt, edm::es::products(), python::multivaluedict::sort(), HepMCValidationHelper::sortByPt(), matplotRender::t, Wdaughters, Wmass, WmassPeak, Wpt, WptLog, Wrap, x, detailsBasic3DVector::y, and z.

{ 
  
  // we *DO NOT* rely on a Z entry in the particle listings!

  edm::Handle<HepMCProduct> evt;
  iEvent.getByLabel(hepmcCollection_, evt);

  //Get EVENT
  const HepMC::GenEvent *myGenEvent = evt->GetEvent();

  nEvt->Fill(0.5);

  std::vector<const HepMC::GenParticle*> allleptons; 
  std::vector<const HepMC::GenParticle*> allneutrinos; 

  //requires status 1 for leptons and neutrinos (except tau)
  int requiredstatus = (abs(_flavor) == 11 || abs(_flavor) ==13 ) ? 1 : 3;

  bool vetotau = true; //(abs(_flavor) == 11 || abs(_flavor) == 12 || abs(_flavor) ==13 || abs(_flavor) ==14 || abs(_flavor) ==16) ? true : false;  

  for(HepMC::GenEvent::particle_const_iterator iter = myGenEvent->particles_begin(); iter != myGenEvent->particles_end(); ++iter) {
    if (vetotau) {
      if ((*iter)->status()==3 && abs((*iter)->pdg_id() == 15) ) return;
    }
    if((*iter)->status()==requiredstatus) {
      //@todo: improve this selection   
      if((*iter)->pdg_id()==_flavor)
        allleptons.push_back(*iter);
      else if (abs((*iter)->pdg_id()) == abs(_flavor)+1)
        allneutrinos.push_back(*iter);  
    }
  }
 
  //nothing to do if we don't have 2 particles
  if (allleptons.empty() || allneutrinos.empty()) return;

  //sort them in pt
  std::sort(allleptons.begin(), allleptons.end(), HepMCValidationHelper::sortByPt); 
  std::sort(allneutrinos.begin(), allneutrinos.end(), HepMCValidationHelper::sortByPt); 

  //get the first lepton and the first neutrino, and check that one is particle one is antiparticle (product of pdgids < 0) 
  std::vector<const HepMC::GenParticle*> products;
  if (allleptons.front()->pdg_id() * allneutrinos.front()->pdg_id() > 0) return;        

  //require at least 20 GeV on the lepton
  if (allleptons.front()->momentum().perp() < 20. || allneutrinos.front()->momentum().perp() < 20. ) return;

  //find possible qed fsr photons
  std::vector<const HepMC::GenParticle*> selectedLepton;
  selectedLepton.push_back(allleptons.front()); 
  std::vector<const HepMC::GenParticle*> fsrphotons;
  HepMCValidationHelper::findFSRPhotons(selectedLepton, myGenEvent, 0.1, fsrphotons);

  Wdaughters->Fill(allleptons.front()->pdg_id()); 
  Wdaughters->Fill(allneutrinos.front()->pdg_id()); 
 
  //assemble FourMomenta
  TLorentzVector lep1(allleptons[0]->momentum().x(), allleptons[0]->momentum().y(), allleptons[0]->momentum().z(), allleptons[0]->momentum().t()); 
  TLorentzVector lep2(allneutrinos[0]->momentum().x(), allneutrinos[0]->momentum().y(), allneutrinos[0]->momentum().z(), allneutrinos[0]->momentum().t()); 
  TLorentzVector dilepton_mom = lep1 + lep2;
  TLorentzVector dilepton_andphoton_mom = dilepton_mom;
  std::vector<TLorentzVector> gammasMomenta;
  for (unsigned int ipho = 0; ipho < fsrphotons.size(); ++ipho){
    TLorentzVector phomom(fsrphotons[ipho]->momentum().x(), fsrphotons[ipho]->momentum().y(), fsrphotons[ipho]->momentum().z(), fsrphotons[ipho]->momentum().t()); 
    dilepton_andphoton_mom += phomom;
    Wdaughters->Fill(fsrphotons[ipho]->pdg_id());
    gammasMomenta.push_back(phomom);
  }  
  //Fill "true" W histograms
  Wmass->Fill(dilepton_andphoton_mom.M());
  WmassPeak->Fill(dilepton_andphoton_mom.M());
  Wpt->Fill(dilepton_andphoton_mom.Pt());
  WptLog->Fill(log10(dilepton_andphoton_mom.Pt())); 
  Wrap->Fill(dilepton_andphoton_mom.Rapidity());

  TLorentzVector met_mom = HepMCValidationHelper::genMet(myGenEvent, -3., 3.);
  TLorentzVector lep1T(lep1.Px(), lep1.Py(), 0., lep1.Et());
  TLorentzVector lepmet_mom = lep1T + met_mom;
  //Fill lepmet histograms
  lepmet_mT->Fill(lepmet_mom.M());
  lepmet_mTPeak->Fill(lepmet_mom.M());
  lepmet_pt->Fill(lepmet_mom.Pt());
  lepmet_ptLog->Fill(log10(lepmet_mom.Pt()));

  //Fill lepton histograms 
  leppt->Fill(lep1.Pt());
  lepeta->Fill(lep1.Eta());
  met->Fill(met_mom.Pt());      

  //boost everything in the W frame
  TVector3 boost = dilepton_andphoton_mom.BoostVector();
  boost*=-1.;
  lep1.Boost(boost);
  lep2.Boost(boost);
  for (unsigned int ipho = 0; ipho < gammasMomenta.size(); ++ipho){
    gammasMomenta[ipho].Boost(boost);
  }
  std::sort(gammasMomenta.begin(), gammasMomenta.end(), HepMCValidationHelper::GreaterByE<TLorentzVector>);

  //fill gamma histograms
  if (gammasMomenta.size() != 0 && dilepton_andphoton_mom.M() > 50.) {
    gamma_energy->Fill(gammasMomenta.front().E());
    double dphi = lep1.DeltaR(gammasMomenta.front());
    cos_theta_gamma_lepton->Fill(cos(dphi));
  } 


}//analyze
void WValidation::beginJob ( void  ) [virtual]

Setting the DQM top directories

Reimplemented from edm::EDAnalyzer.

Definition at line 31 of file WValidation.cc.

References _name, DQMStore::book1D(), cos_theta_gamma_lepton, dbe, gamma_energy, lepeta, lepmet_mT, lepmet_mTPeak, lepmet_pt, lepmet_ptLog, leppt, met, nEvt, DQMStore::setCurrentFolder(), Wdaughters, Wmass, WmassPeak, Wpt, WptLog, and Wrap.

{
  if(dbe){
    std::string folderName = "Generator/W";
    folderName+=_name;
    dbe->setCurrentFolder(folderName.c_str());
    
    // Number of analyzed events
    nEvt = dbe->book1D("nEvt", "n analyzed Events", 1, 0., 1.);
    
    //Kinematics
    Wmass = dbe->book1D("Wmass","inv. Mass W", 70 ,0,140);
    WmassPeak = dbe->book1D("WmassPeak","inv. Mass W", 80 ,80 ,100);
    Wpt = dbe->book1D("Wpt","W pt",100,0,200);
    WptLog = dbe->book1D("WptLog","log(W pt)",100,0.,5.);
    Wrap = dbe->book1D("Wrap", "W y", 100, -5, 5);
    Wdaughters = dbe->book1D("Wdaughters", "W daughters", 60, -30, 30);

    lepmet_mT = dbe->book1D("lepmet_mT","lepton-met transverse mass", 70 ,0,140);
    lepmet_mTPeak = dbe->book1D("lepmet_mTPeak","lepton-met transverse mass", 80 ,80 ,100);
    lepmet_pt = dbe->book1D("lepmet_pt","lepton-met",100,0,200);
    lepmet_ptLog = dbe->book1D("lepmet_ptLog","log(lepton-met pt)",100,0.,5.);

    gamma_energy = dbe->book1D("gamma_energy", "photon energy in W rest frame", 200, 0., 100.);
    cos_theta_gamma_lepton = dbe->book1D("cos_theta_gamma_lepton",      "cos_theta_gamma_lepton in W rest frame",      200, -1, 1);

    leppt = dbe->book1D("leadpt","lepton pt", 200, 0., 200.);    
    met   = dbe->book1D("met","met", 200, 0., 200.);    
    lepeta = dbe->book1D("leadeta","leading lepton eta", 100, -5., 5.);

  }
  return;
}
void WValidation::beginRun ( const edm::Run iRun,
const edm::EventSetup iSetup 
) [virtual]

Get PDT Table

Reimplemented from edm::EDAnalyzer.

Definition at line 67 of file WValidation.cc.

References fPDGTable, and edm::EventSetup::getData().

{
  iSetup.getData( fPDGTable );
  return;
}
void WValidation::endJob ( void  ) [virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 66 of file WValidation.cc.

{return;}
void WValidation::endRun ( const edm::Run iRun,
const edm::EventSetup iSetup 
) [virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 73 of file WValidation.cc.

{return;}

Member Data Documentation

int WValidation::_flavor [private]

decay flavor

Definition at line 62 of file WValidation.h.

Referenced by analyze().

std::string WValidation::_name [private]

decay flavor name

Definition at line 64 of file WValidation.h.

Referenced by beginJob().

Definition at line 59 of file WValidation.h.

Referenced by analyze(), and beginJob().

ME's "container".

Definition at line 53 of file WValidation.h.

Referenced by beginJob(), and WValidation().

PDT table.

Definition at line 50 of file WValidation.h.

Referenced by beginRun().

Definition at line 59 of file WValidation.h.

Referenced by analyze(), and beginJob().

Definition at line 47 of file WValidation.h.

Referenced by analyze().

Definition at line 58 of file WValidation.h.

Referenced by analyze(), and beginJob().

Definition at line 57 of file WValidation.h.

Referenced by analyze(), and beginJob().

Definition at line 57 of file WValidation.h.

Referenced by analyze(), and beginJob().

Definition at line 57 of file WValidation.h.

Referenced by analyze(), and beginJob().

Definition at line 57 of file WValidation.h.

Referenced by analyze(), and beginJob().

Definition at line 57 of file WValidation.h.

Definition at line 58 of file WValidation.h.

Referenced by analyze(), and beginJob().

Definition at line 58 of file WValidation.h.

Referenced by analyze(), and beginJob().

Definition at line 55 of file WValidation.h.

Referenced by analyze(), and beginJob().

Definition at line 56 of file WValidation.h.

Referenced by analyze(), and beginJob().

Definition at line 56 of file WValidation.h.

Referenced by analyze(), and beginJob().

Definition at line 56 of file WValidation.h.

Referenced by analyze(), and beginJob().

Definition at line 56 of file WValidation.h.

Referenced by analyze(), and beginJob().

Definition at line 56 of file WValidation.h.

Referenced by analyze(), and beginJob().

Definition at line 56 of file WValidation.h.

Referenced by analyze(), and beginJob().