CMS 3D CMS Logo

Public Member Functions | Private Member Functions | Private Attributes

CaloTriggerAnalyzer2 Class Reference

#include <CaloTriggerAnalyzer2.h>

Inheritance diagram for CaloTriggerAnalyzer2:
edm::EDAnalyzer

List of all members.

Public Member Functions

 CaloTriggerAnalyzer2 (const edm::ParameterSet &)
 ~CaloTriggerAnalyzer2 ()

Private Member Functions

virtual void analyze (const edm::Event &, const edm::EventSetup &)
virtual void clearVectors ()
virtual void fillLHCBranches (const reco::Candidate *L1object, const reco::Candidate *genParticle)
virtual void fillSLHCBranches (const reco::Candidate *L1object, const reco::Candidate *genParticle)

Private Attributes

TH1F * absEta
TTree * CandTree
TH1F * dEta
TH1F * dPhi
TH1F * dPt
double DR_
TH1F * eta
TH1F * etaDenom
TH1F * etaNum
TTree * EventTree
std::vector< double > gEta
float ghighMPt
float ghighPt
std::vector< double > gPhi
std::vector< double > gPt
float gsecondPt
float highestGenPt
TH1F * highestPt
TH1F * highestPtGen
float highPt
const reco::CandidateL1closest
const reco::CandidateL1object
TH1F * LHCalldR
std::vector< double > LHCdEta
std::vector< double > LHCdEtac
bool LHCdoubleTrigger
std::vector< double > LHCdPhi
std::vector< double > LHCdPhic
std::vector< double > LHCdPt
std::vector< double > LHCdPtc
std::vector< double > LHCdR
std::vector< double > LHCdRc
float LHChighPt
std::vector< double > LHCL1Eta
std::vector< double > LHCL1Etac
std::vector< double > LHCL1Phi
std::vector< double > LHCL1Phic
std::vector< double > LHCL1Pt
std::vector< double > LHCL1Ptc
std::vector< bool > LHCpassDoubleThresh
std::vector< bool > LHCpassSingleThresh
std::vector< double > LHCRPt
std::vector< double > LHCRPtc
float LHCsecondPt
bool LHCsingleTrigger
edm::InputTag LHCsrc_
double maxEta_
int nEvents
TH1F * pt
TH1F * ptDenom
TH1F * ptNum
edm::InputTag ref_
TH1F * RPt
TProfile * RPtEta
TProfile * RPtEtaFull
float secondGenPt
TH1F * secondPt
float secondPtf
TH1F * secondPtGen
TH1F * SLHCalldR
std::vector< double > SLHCdEta
std::vector< double > SLHCdEtac
bool SLHCdoubleTrigger
std::vector< double > SLHCdPhi
std::vector< double > SLHCdPhic
std::vector< double > SLHCdPt
std::vector< double > SLHCdPtc
std::vector< double > SLHCdR
std::vector< double > SLHCdRc
float SLHChighPt
std::vector< double > SLHCL1Eta
std::vector< double > SLHCL1Etac
std::vector< double > SLHCL1Phi
std::vector< double > SLHCL1Phic
std::vector< double > SLHCL1Pt
std::vector< double > SLHCL1Ptc
std::vector< bool > SLHCpassDoubleThresh
std::vector< bool > SLHCpassSingleThresh
std::vector< double > SLHCRPt
std::vector< double > SLHCRPtc
float SLHCsecondPt
bool SLHCsingleTrigger
edm::InputTag SLHCsrc_
double threshold_

Detailed Description

Definition at line 27 of file CaloTriggerAnalyzer2.h.


Constructor & Destructor Documentation

CaloTriggerAnalyzer2::CaloTriggerAnalyzer2 ( const edm::ParameterSet iConfig) [explicit]

Definition at line 7 of file CaloTriggerAnalyzer2.cc.

References absEta, CandTree, dEta, dPhi, dPt, etaDenom, etaNum, EventTree, gEta, ghighMPt, ghighPt, gPhi, gPt, gsecondPt, highestPt, highestPtGen, LHCalldR, LHCdEta, LHCdoubleTrigger, LHCdPhi, LHCdPt, LHCdR, LHChighPt, LHCL1Eta, LHCL1Phi, LHCL1Pt, LHCpassSingleThresh, LHCRPt, LHCsecondPt, LHCsingleTrigger, nEvents, pt, ptDenom, ptNum, RPt, RPtEta, RPtEtaFull, secondPt, secondPtGen, SLHCalldR, SLHCdEta, SLHCdoubleTrigger, SLHCdPhi, SLHCdPt, SLHCdR, SLHChighPt, SLHCL1Eta, SLHCL1Phi, SLHCL1Pt, SLHCpassSingleThresh, SLHCRPt, SLHCsecondPt, and SLHCsingleTrigger.

                                                                        :
  ref_(iConfig.getParameter<edm::InputTag>("ref")),
  SLHCsrc_(iConfig.getParameter<edm::InputTag>("SLHCsrc")),
  LHCsrc_(iConfig.getParameter<edm::InputTag>("LHCsrc")),
  DR_(iConfig.getParameter<double>("deltaR")),
  threshold_(iConfig.getParameter<double>("threshold")),
  maxEta_(iConfig.getUntrackedParameter<double>("maxEta",2.5))
{
   //now do what ever initialization is needed

  edm::Service<TFileService> fs;

  EventTree = fs->make<TTree>("EventTree","Tree containing event-by-event info");
  CandTree = fs->make<TTree>("CandTree","Tree containing candidate info");

  CandTree->Branch("gPt",&gPt);
  CandTree->Branch("gEta",&gEta);
  CandTree->Branch("gPhi",&gPhi);
 
  CandTree->Branch("SLHCpassSingleThresh",&SLHCpassSingleThresh);
  CandTree->Branch("SLHCL1Pt",&SLHCL1Pt);
  CandTree->Branch("SLHCL1Eta",&SLHCL1Eta);
  CandTree->Branch("SLHCL1Phi",&SLHCL1Phi);
  CandTree->Branch("SLHCdR",&SLHCdR);
  CandTree->Branch("SLHCdPt",&SLHCdPt);
  CandTree->Branch("SLHCdEta",&SLHCdEta);
  CandTree->Branch("SLHCdPhi",&SLHCdPhi);
  CandTree->Branch("SLHCRPt",&SLHCRPt);
  
  CandTree->Branch("LHCpassSingleThresh",&LHCpassSingleThresh);
  CandTree->Branch("LHCL1Pt",&LHCL1Pt);
  CandTree->Branch("LHCL1Eta",&LHCL1Eta);
  CandTree->Branch("LHCL1Phi",&LHCL1Phi);
  CandTree->Branch("LHCdR",&LHCdR);
  CandTree->Branch("LHCdPt",&LHCdPt);
  CandTree->Branch("LHCdEta",&LHCdEta);
  CandTree->Branch("LHCdPhi",&LHCdPhi);
  CandTree->Branch("LHCRPt",&LHCRPt);
    
  EventTree->Branch("ghighPt",&ghighPt);
  EventTree->Branch("gsecondPt",&gsecondPt);
  EventTree->Branch("ghighMPt",&ghighMPt);
  EventTree->Branch("SLHChighPt",&SLHChighPt);
  EventTree->Branch("SLHCsecondPt",&SLHCsecondPt);
  EventTree->Branch("LHChighPt",&LHChighPt);
  EventTree->Branch("LHCsecondPt",&LHCsecondPt);
  EventTree->Branch("nEvents",&nEvents);
  EventTree->Branch("SLHCsingleTrigger",&SLHCsingleTrigger);
  EventTree->Branch("SLHCdoubleTrigger",&SLHCdoubleTrigger);
  EventTree->Branch("LHCsingleTrigger",&LHCsingleTrigger);
  EventTree->Branch("LHCdoubleTrigger",&LHCdoubleTrigger);

  ptNum    = fs->make<TH1F>( "ptNum"   , "ptNum"   , 20  ,  0., 100. );
  ptDenom  = fs->make<TH1F>( "ptDenom" , "ptDenom" , 20  ,  0., 100. );
  etaNum   = fs->make<TH1F>( "etaNum"  , "etaNum"  , 20  ,  -2.5, 2.5 );
  etaDenom = fs->make<TH1F>( "etaDenom", "etaDenom", 20  ,  -2.5, 2.5 );
  pt       = fs->make<TH1F>( "pt"      , "pt", 20  ,  0. , 100. );
  highestPt= fs->make<TH1F>( "highestPt"      , "highestPt", 50  ,  0. , 100. );
  secondPt = fs->make<TH1F>( "secondHighestPt", "secondHighestPt", 50  ,  0. , 100. );
  highestPtGen=fs->make<TH1F>("highestPtGen","highestPtGen",100, 0., 200.);
  secondPtGen=fs->make<TH1F>("secondPtGen","secondPtGen",100, 0., 200.);
  dPt      = fs->make<TH1F>( "dPt"      , "dPt", 50  , -1  , 1 );
  dEta     = fs->make<TH1F>( "dEta"      , "dEta", 50  , -0.5  , 0.5 );
  dPhi     = fs->make<TH1F>( "dPhi"      , "dPhi", 50  , -0.5  , 0.5 );
  RPt = fs->make<TH1F>("RPt", "Pt_Ratio", 10, 0, 2.);
  RPtEta = fs->make<TProfile>("RPtEta","Pt_Ratio as fcn of abs(eta)",13,0.0,2.6,0,2);
  RPtEtaFull = fs->make<TProfile>("RPtEtaFull","Pt_Ratio as fcn of eta",26,-2.6,2.6,0,2);
  absEta = fs->make<TH1F>("abseta","abs(eta)",13,0.,2.6); 
  SLHCalldR = fs->make<TH1F>("SLHCalldR","delta(R) (SLHC)",50,0,2);
  LHCalldR = fs->make<TH1F>("LHCalldR","delta(R) (LHC)",50,0,2);

}
CaloTriggerAnalyzer2::~CaloTriggerAnalyzer2 ( )

Definition at line 81 of file CaloTriggerAnalyzer2.cc.

{

   // do anything here that needs to be done at desctruction time
   // (e.g. close files, deallocate resources etc.)

}

Member Function Documentation

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

Implements edm::EDAnalyzer.

Definition at line 96 of file CaloTriggerAnalyzer2.cc.

References absEta, CandTree, clearVectors(), dEta, dPhi, dPt, DR_, etaDenom, etaNum, EventTree, fillLHCBranches(), fillSLHCBranches(), gEta, edm::Event::getByLabel(), ghighMPt, ghighPt, gPhi, gPt, gsecondPt, highestGenPt, highestPt, highestPtGen, highPt, i, j, L1closest, L1object, LHCalldR, LHCdoubleTrigger, LHChighPt, LHCpassDoubleThresh, LHCpassSingleThresh, LHCsecondPt, LHCsingleTrigger, LHCsrc_, maxEta_, nEvents, NULL, pt, ptDenom, ptNum, ref_, RPt, RPtEta, RPtEtaFull, secondGenPt, secondPt, secondPtf, secondPtGen, SLHCalldR, SLHCdoubleTrigger, SLHChighPt, SLHCpassDoubleThresh, SLHCpassSingleThresh, SLHCsecondPt, SLHCsingleTrigger, SLHCsrc_, alcazmumu_cfi::src, and threshold_.

{
  //  std::cout << src_ << std::endl;
  edm::Handle<edm::View<reco::Candidate> > SLHCsrc;
  edm::Handle<edm::View<reco::Candidate> > LHCsrc;
  edm::Handle<edm::View<reco::Candidate> > ref;
  edm::Handle<edm::View<reco::Candidate> > src;

  bool gotSLHCsrc = iEvent.getByLabel(SLHCsrc_,SLHCsrc);
  bool gotLHCsrc = iEvent.getByLabel(LHCsrc_,LHCsrc);
  bool gotRef = iEvent.getByLabel(ref_,ref);
  bool gotSrc =iEvent.getByLabel(SLHCsrc_,src);

  //clear variables...
  clearVectors();

  // const reco::Candidate * L1object;

  if(gotSrc) {
    highPt=0;
    secondPtf=0;
    for(edm::View<reco::Candidate>::const_iterator i = src->begin(); i!= src->end();++i)
      {
        pt->Fill(i->pt());
        if (i->pt()>highPt){
          secondPtf=highPt;
          highPt=i->pt();
        } else if (i->pt()>secondPtf){
          secondPtf=i->pt();
        }
      }

    if(src->size()>0)
      highestPt->Fill(highPt);
    else
      highestPt->Fill(0.0);


    if(src->size()>1)
      secondPt->Fill(secondPtf);
    else
      secondPt->Fill(0.0);
  }

  if(gotRef) {
    //get highest Pt gen object--loop over all to make sure it's the highest!
    highestGenPt=-1.0;
    secondGenPt=-1.0;
    ghighMPt=0;
    for(edm::View<reco::Candidate>::const_iterator i = ref->begin(); i!= ref->end();++i){

      if (i->pt()>highestGenPt){
        highestGenPt=i->pt();
        highestPtGen->Fill(i->pt());
      } else if (i->pt()>secondGenPt){
        secondGenPt=i->pt();
        secondPtGen->Fill(i->pt());
      }
      if(fabs(i->eta())<maxEta_&&i->pt()>threshold_/2.)
        {
          gPt.push_back(i->pt());
          gEta.push_back(i->eta());
          gPhi.push_back(i->phi());

          ptDenom->Fill(i->pt());
          etaDenom->Fill(i->eta());
          
          //printf("ref pt = %f  eta = %f phi = %f\n",i->pt(),i->eta(),i->phi());
          
          if(gotSLHCsrc) {
            bool matched=false;
            bool halfmatched=false;
            double mindR=900;
            L1object=NULL;
            L1closest=NULL;
            math::XYZTLorentzVector highestV(0.0001,0.,0.,0.0002);
            for(edm::View<reco::Candidate>::const_iterator j = SLHCsrc->begin(); j!= SLHCsrc->end();++j)
              {
                SLHCalldR->Fill(ROOT::Math::VectorUtil::DeltaR(i->p4(),j->p4()));
                
                if (ROOT::Math::VectorUtil::DeltaR(i->p4(),j->p4())<mindR){
                  mindR=ROOT::Math::VectorUtil::DeltaR(i->p4(),j->p4());
                  L1closest=&(*j);
                }

                if(ROOT::Math::VectorUtil::DeltaR(i->p4(),j->p4())<DR_) {
                  if(j->pt()>threshold_/2){
                    if(j->pt()>threshold_){
                      //printf("matched pt = %f  eta = %f phi = %f\n",j->pt(),j->eta(),j->phi());
                      matched=true;
                    }
                    halfmatched=true;
                  }
                  if(j->pt()>highestV.pt()){
                    highestV = j->p4();
                    L1object=&(*j);
                  }
                }
              } //end for loop over L1 objects
            
            //Fill L1 object tree
            fillSLHCBranches(L1object, &(*i));

            if(halfmatched) {
              SLHCpassDoubleThresh.push_back(true);
            } else{
              SLHCpassDoubleThresh.push_back(false);
            }
            if(matched) {
              SLHCpassSingleThresh.push_back(true);

              if (i->pt()>ghighMPt)
                ghighMPt=i->pt();
              
              RPt->Fill(i->pt()/highestV.pt());
              RPtEtaFull->Fill(highestV.eta(), i->pt()/highestV.pt());
              RPtEta->Fill(fabs(highestV.eta()), i->pt()/highestV.pt());
              
              //                printf("matched abs(eta) = %f \n", fabs(highestV.eta()) );
              absEta->Fill(fabs(highestV.eta()));
              dPt->Fill((highestV.pt()-i->pt())/i->pt());
              dEta->Fill(highestV.eta()-i->eta());
              dPhi->Fill(highestV.phi()-i->phi());
              
              ptNum->Fill(i->pt());
              etaNum->Fill(i->eta());
            } else {
              SLHCpassSingleThresh.push_back(false);
            }
          } //end (if gotSLHCsrc)

          if(gotLHCsrc) {
            bool matched=false;
            bool halfmatched=false;
            double mindR=900;
            L1object=NULL;
            printf("LHC src is %i in size \n",LHCsrc->size());
            math::XYZTLorentzVector highestV(0.0001,0.,0.,0.0002);
            for(edm::View<reco::Candidate>::const_iterator j = LHCsrc->begin(); j!= LHCsrc->end();++j)
              {
                LHCalldR->Fill(ROOT::Math::VectorUtil::DeltaR(i->p4(),j->p4()));

                if (ROOT::Math::VectorUtil::DeltaR(i->p4(),j->p4())<mindR){
                  mindR=ROOT::Math::VectorUtil::DeltaR(i->p4(),j->p4());
                  L1closest=&(*j);
                }
                
                if(ROOT::Math::VectorUtil::DeltaR(i->p4(),j->p4())<DR_) {
                  if(j->pt()>threshold_/2){
                    if(j->pt()>threshold_){
                      //printf("matched pt = %f  eta = %f phi = %f\n",j->pt(),j->eta(),j->phi());
                      matched=true;
                    }
                    halfmatched=true;
                  }
                  if(j->pt()>highestV.pt()){
                    highestV = j->p4();
                    L1object=&(*j); //Match to highest Pt object in cone.  If none is found, the written object is the closest in dR.
                  }
                }
              } //end for loop over L1 objects                                                                                                                

            //Fill L1 object tree   
            fillLHCBranches(L1object, &(*i));
            
            if(halfmatched) {
              LHCpassDoubleThresh.push_back(true);
            } else{
              LHCpassDoubleThresh.push_back(false);
            }
          
            if(matched) {
              LHCpassSingleThresh.push_back(true);
              RPt->Fill(i->pt()/highestV.pt());
              RPtEtaFull->Fill(highestV.eta(), i->pt()/highestV.pt());
              RPtEta->Fill(fabs(highestV.eta()), i->pt()/highestV.pt());

              //                printf("matched abs(eta) = %f \n", fabs(highestV.eta()) );                                                                    
              absEta->Fill(fabs(highestV.eta()));
              dPt->Fill((highestV.pt()-i->pt())/i->pt());
              dEta->Fill(highestV.eta()-i->eta());
              dPhi->Fill(highestV.phi()-i->phi());

              ptNum->Fill(i->pt());
              etaNum->Fill(i->eta());
            } else {
              LHCpassSingleThresh.push_back(false);
            }
          } //end (if gotLHCsrc) 
          
        } //end (if GEN in acceptance)
    } //end GEN loop 
  }  //end (if gotRef)
  CandTree->Fill();
  
  //-------------------Event-by-event work-----------------------
  nEvents=0; //used for filling number of events.

  if(gotRef) {
    ghighPt=0.0;
    gsecondPt=0.0;
    for(edm::View<reco::Candidate>::const_iterator i = ref->begin(); i!= ref->end();++i)
      {
        if (i->pt()>ghighPt && fabs(i->eta())<maxEta_){
          gsecondPt=ghighPt;
          ghighPt=i->pt();
        } else if (i->pt()>gsecondPt && fabs(i->eta())<maxEta_){
          gsecondPt=i->pt();
        }
      }
  }

  if(gotSLHCsrc) {
    SLHChighPt=0.0;
    SLHCsecondPt=0.0;
    for(edm::View<reco::Candidate>::const_iterator i = SLHCsrc->begin(); i!= SLHCsrc->end();++i)
      {
        if (i->pt()>SLHChighPt){
          SLHCsecondPt=SLHChighPt;
          SLHChighPt=i->pt();
        } else if (i->pt()>SLHCsecondPt){
          SLHCsecondPt=i->pt();
        }
      }
  }

  if(gotLHCsrc) {
    LHChighPt=0.0;
    LHCsecondPt=0.0;
    for(edm::View<reco::Candidate>::const_iterator i = LHCsrc->begin(); i!= LHCsrc->end();++i)
      {
        if (i->pt()>LHChighPt){
          LHCsecondPt=LHChighPt;
          LHChighPt=i->pt();
        } else if (i->pt()>LHCsecondPt){
          LHCsecondPt=i->pt();
        }
      }
  }

  SLHCsingleTrigger=false;
  for(unsigned int i=0; i<SLHCpassSingleThresh.size(); i++){
    if(SLHCpassSingleThresh.at(i))
      SLHCsingleTrigger=true;
  }
  SLHCdoubleTrigger=false;
  unsigned int nSLHCdouble=0;
  for(unsigned int i=0; i<SLHCpassDoubleThresh.size(); i++){
    if(SLHCpassDoubleThresh.at(i))
      nSLHCdouble++;
  }
  if(nSLHCdouble>=2)
    SLHCdoubleTrigger=true;

  LHCsingleTrigger=false;
  for(unsigned int i=0; i<LHCpassSingleThresh.size(); i++){
    if(LHCpassSingleThresh.at(i))
      LHCsingleTrigger=true;
  }
  LHCdoubleTrigger=false;
  unsigned int nLHCdouble=0;
  for(unsigned int i=0; i<LHCpassDoubleThresh.size(); i++){
    if(LHCpassDoubleThresh.at(i))
      nLHCdouble++;
  }
  if(nLHCdouble>=2)
    LHCdoubleTrigger=true;

  // printf("%s \n", (SLHCsingleTrigger && ghighPt>37.0)?"true":"false");

  EventTree->Fill();
}
void CaloTriggerAnalyzer2::clearVectors ( ) [private, virtual]

Definition at line 369 of file CaloTriggerAnalyzer2.cc.

References gEta, gPhi, gPt, LHCdEta, LHCdPhi, LHCdPt, LHCdR, LHCL1Eta, LHCL1Phi, LHCL1Pt, LHCpassDoubleThresh, LHCpassSingleThresh, LHCRPt, SLHCdEta, SLHCdPhi, SLHCdPt, SLHCdR, SLHCL1Eta, SLHCL1Phi, SLHCL1Pt, SLHCpassDoubleThresh, SLHCpassSingleThresh, and SLHCRPt.

Referenced by analyze().

                                        {
  gPt.clear();
  gEta.clear();
  gPhi.clear();

  SLHCL1Pt.clear();
  SLHCL1Eta.clear();
  SLHCL1Phi.clear();
  SLHCdR.clear();
  SLHCdPt.clear();
  SLHCdEta.clear();
  SLHCdPhi.clear();
  SLHCRPt.clear();
  SLHCpassSingleThresh.clear();
  SLHCpassDoubleThresh.clear();
 
  LHCL1Pt.clear();
  LHCL1Eta.clear();
  LHCL1Phi.clear();
  LHCdR.clear();
  LHCdPt.clear();
  LHCdEta.clear();
  LHCdPhi.clear();
  LHCRPt.clear();
  LHCpassSingleThresh.clear();
  LHCpassDoubleThresh.clear();
}
void CaloTriggerAnalyzer2::fillLHCBranches ( const reco::Candidate L1object,
const reco::Candidate genParticle 
) [private, virtual]

Definition at line 420 of file CaloTriggerAnalyzer2.cc.

References reco::Candidate::eta(), LHCdEta, LHCdPhi, LHCdPt, LHCdR, LHCL1Eta, LHCL1Phi, LHCL1Pt, LHCRPt, NULL, reco::Candidate::p4(), reco::Candidate::phi(), and reco::Candidate::pt().

Referenced by analyze().

                                                                                                        {
  if(cand==NULL){
    //    printf("No match found for LHC!\n");
    LHCL1Pt.push_back(-1);
    LHCL1Eta.push_back(-137);
    LHCL1Phi.push_back(-137);
    LHCdR.push_back(-137);
    LHCdPt.push_back(-137);
    LHCdEta.push_back(-137);
    LHCdPhi.push_back(-137);
    LHCRPt.push_back(-137);
  } else{
    printf("LHC L1Pt=%f",cand->pt());
    LHCL1Pt.push_back(cand->pt());
    LHCL1Eta.push_back(cand->eta());
    LHCL1Phi.push_back(cand->phi());
    LHCdR.push_back(ROOT::Math::VectorUtil::DeltaR(cand->p4(),genParticle->p4()));
    LHCdPt.push_back((cand->pt()-genParticle->pt())/genParticle->pt());
    LHCdEta.push_back(cand->eta()-genParticle->eta());
    LHCdPhi.push_back(cand->phi()-genParticle->phi());
    LHCRPt.push_back(genParticle->pt()/cand->pt());
  }
}
void CaloTriggerAnalyzer2::fillSLHCBranches ( const reco::Candidate L1object,
const reco::Candidate genParticle 
) [private, virtual]

Definition at line 397 of file CaloTriggerAnalyzer2.cc.

References reco::Candidate::eta(), NULL, reco::Candidate::p4(), reco::Candidate::phi(), reco::Candidate::pt(), SLHCdEta, SLHCdPhi, SLHCdPt, SLHCdR, SLHCL1Eta, SLHCL1Phi, SLHCL1Pt, and SLHCRPt.

Referenced by analyze().

                                                                                                         {
  if(cand==NULL){
    //printf("No match found for SLHC!\n");
    SLHCL1Pt.push_back(-1);
    SLHCL1Eta.push_back(-137);
    SLHCL1Phi.push_back(-137);
    SLHCdR.push_back(-137);
    SLHCdPt.push_back(-137);
    SLHCdEta.push_back(-137);
    SLHCdPhi.push_back(-137);
    SLHCRPt.push_back(-137);
  } else {
    SLHCL1Pt.push_back(cand->pt());
    SLHCL1Eta.push_back(cand->eta());
    SLHCL1Phi.push_back(cand->phi());
    SLHCdR.push_back(ROOT::Math::VectorUtil::DeltaR(cand->p4(),genParticle->p4()));
    SLHCdPt.push_back((cand->pt()-genParticle->pt())/genParticle->pt());
    SLHCdEta.push_back(cand->eta()-genParticle->eta());
    SLHCdPhi.push_back(cand->phi()-genParticle->phi());
    SLHCRPt.push_back(genParticle->pt()/cand->pt());
  }
}

Member Data Documentation

Definition at line 129 of file CaloTriggerAnalyzer2.h.

Referenced by analyze(), and CaloTriggerAnalyzer2().

Definition at line 47 of file CaloTriggerAnalyzer2.h.

Referenced by analyze(), and CaloTriggerAnalyzer2().

TH1F* CaloTriggerAnalyzer2::dEta [private]

Definition at line 122 of file CaloTriggerAnalyzer2.h.

Referenced by analyze(), and CaloTriggerAnalyzer2().

TH1F* CaloTriggerAnalyzer2::dPhi [private]

Definition at line 123 of file CaloTriggerAnalyzer2.h.

Referenced by analyze(), and CaloTriggerAnalyzer2().

TH1F* CaloTriggerAnalyzer2::dPt [private]

Definition at line 121 of file CaloTriggerAnalyzer2.h.

Referenced by analyze(), and CaloTriggerAnalyzer2().

double CaloTriggerAnalyzer2::DR_ [private]

Definition at line 42 of file CaloTriggerAnalyzer2.h.

Referenced by analyze().

TH1F* CaloTriggerAnalyzer2::eta [private]

Definition at line 114 of file CaloTriggerAnalyzer2.h.

Definition at line 119 of file CaloTriggerAnalyzer2.h.

Referenced by analyze(), and CaloTriggerAnalyzer2().

Definition at line 118 of file CaloTriggerAnalyzer2.h.

Referenced by analyze(), and CaloTriggerAnalyzer2().

Definition at line 46 of file CaloTriggerAnalyzer2.h.

Referenced by analyze(), and CaloTriggerAnalyzer2().

std::vector<double> CaloTriggerAnalyzer2::gEta [private]

Definition at line 50 of file CaloTriggerAnalyzer2.h.

Referenced by analyze(), CaloTriggerAnalyzer2(), and clearVectors().

Definition at line 93 of file CaloTriggerAnalyzer2.h.

Referenced by analyze(), and CaloTriggerAnalyzer2().

Definition at line 91 of file CaloTriggerAnalyzer2.h.

Referenced by analyze(), and CaloTriggerAnalyzer2().

std::vector<double> CaloTriggerAnalyzer2::gPhi [private]

Definition at line 51 of file CaloTriggerAnalyzer2.h.

Referenced by analyze(), CaloTriggerAnalyzer2(), and clearVectors().

std::vector<double> CaloTriggerAnalyzer2::gPt [private]

Definition at line 49 of file CaloTriggerAnalyzer2.h.

Referenced by analyze(), CaloTriggerAnalyzer2(), and clearVectors().

Definition at line 92 of file CaloTriggerAnalyzer2.h.

Referenced by analyze(), and CaloTriggerAnalyzer2().

Definition at line 106 of file CaloTriggerAnalyzer2.h.

Referenced by analyze().

Definition at line 124 of file CaloTriggerAnalyzer2.h.

Referenced by analyze(), and CaloTriggerAnalyzer2().

Definition at line 126 of file CaloTriggerAnalyzer2.h.

Referenced by analyze(), and CaloTriggerAnalyzer2().

Definition at line 108 of file CaloTriggerAnalyzer2.h.

Referenced by analyze().

Definition at line 112 of file CaloTriggerAnalyzer2.h.

Referenced by analyze().

Definition at line 111 of file CaloTriggerAnalyzer2.h.

Referenced by analyze().

Definition at line 131 of file CaloTriggerAnalyzer2.h.

Referenced by analyze(), and CaloTriggerAnalyzer2().

std::vector<double> CaloTriggerAnalyzer2::LHCdEta [private]

Definition at line 77 of file CaloTriggerAnalyzer2.h.

Referenced by CaloTriggerAnalyzer2(), clearVectors(), and fillLHCBranches().

std::vector<double> CaloTriggerAnalyzer2::LHCdEtac [private]

Definition at line 85 of file CaloTriggerAnalyzer2.h.

Definition at line 104 of file CaloTriggerAnalyzer2.h.

Referenced by analyze(), and CaloTriggerAnalyzer2().

std::vector<double> CaloTriggerAnalyzer2::LHCdPhi [private]

Definition at line 78 of file CaloTriggerAnalyzer2.h.

Referenced by CaloTriggerAnalyzer2(), clearVectors(), and fillLHCBranches().

std::vector<double> CaloTriggerAnalyzer2::LHCdPhic [private]

Definition at line 86 of file CaloTriggerAnalyzer2.h.

std::vector<double> CaloTriggerAnalyzer2::LHCdPt [private]

Definition at line 76 of file CaloTriggerAnalyzer2.h.

Referenced by CaloTriggerAnalyzer2(), clearVectors(), and fillLHCBranches().

std::vector<double> CaloTriggerAnalyzer2::LHCdPtc [private]

Definition at line 84 of file CaloTriggerAnalyzer2.h.

std::vector<double> CaloTriggerAnalyzer2::LHCdR [private]

Definition at line 75 of file CaloTriggerAnalyzer2.h.

Referenced by CaloTriggerAnalyzer2(), clearVectors(), and fillLHCBranches().

std::vector<double> CaloTriggerAnalyzer2::LHCdRc [private]

Definition at line 83 of file CaloTriggerAnalyzer2.h.

Definition at line 96 of file CaloTriggerAnalyzer2.h.

Referenced by analyze(), and CaloTriggerAnalyzer2().

std::vector<double> CaloTriggerAnalyzer2::LHCL1Eta [private]

Definition at line 73 of file CaloTriggerAnalyzer2.h.

Referenced by CaloTriggerAnalyzer2(), clearVectors(), and fillLHCBranches().

std::vector<double> CaloTriggerAnalyzer2::LHCL1Etac [private]

Definition at line 81 of file CaloTriggerAnalyzer2.h.

std::vector<double> CaloTriggerAnalyzer2::LHCL1Phi [private]

Definition at line 74 of file CaloTriggerAnalyzer2.h.

Referenced by CaloTriggerAnalyzer2(), clearVectors(), and fillLHCBranches().

std::vector<double> CaloTriggerAnalyzer2::LHCL1Phic [private]

Definition at line 82 of file CaloTriggerAnalyzer2.h.

std::vector<double> CaloTriggerAnalyzer2::LHCL1Pt [private]

Definition at line 72 of file CaloTriggerAnalyzer2.h.

Referenced by CaloTriggerAnalyzer2(), clearVectors(), and fillLHCBranches().

std::vector<double> CaloTriggerAnalyzer2::LHCL1Ptc [private]

Definition at line 80 of file CaloTriggerAnalyzer2.h.

std::vector<bool> CaloTriggerAnalyzer2::LHCpassDoubleThresh [private]

Definition at line 89 of file CaloTriggerAnalyzer2.h.

Referenced by analyze(), and clearVectors().

std::vector<bool> CaloTriggerAnalyzer2::LHCpassSingleThresh [private]

Definition at line 88 of file CaloTriggerAnalyzer2.h.

Referenced by analyze(), CaloTriggerAnalyzer2(), and clearVectors().

std::vector<double> CaloTriggerAnalyzer2::LHCRPt [private]

Definition at line 79 of file CaloTriggerAnalyzer2.h.

Referenced by CaloTriggerAnalyzer2(), clearVectors(), and fillLHCBranches().

std::vector<double> CaloTriggerAnalyzer2::LHCRPtc [private]

Definition at line 87 of file CaloTriggerAnalyzer2.h.

Definition at line 97 of file CaloTriggerAnalyzer2.h.

Referenced by analyze(), and CaloTriggerAnalyzer2().

Definition at line 103 of file CaloTriggerAnalyzer2.h.

Referenced by analyze(), and CaloTriggerAnalyzer2().

Definition at line 41 of file CaloTriggerAnalyzer2.h.

Referenced by analyze().

Definition at line 44 of file CaloTriggerAnalyzer2.h.

Referenced by analyze().

Definition at line 99 of file CaloTriggerAnalyzer2.h.

Referenced by analyze(), and CaloTriggerAnalyzer2().

TH1F* CaloTriggerAnalyzer2::pt [private]

Definition at line 120 of file CaloTriggerAnalyzer2.h.

Referenced by analyze(), and CaloTriggerAnalyzer2().

Definition at line 117 of file CaloTriggerAnalyzer2.h.

Referenced by analyze(), and CaloTriggerAnalyzer2().

TH1F* CaloTriggerAnalyzer2::ptNum [private]

Definition at line 116 of file CaloTriggerAnalyzer2.h.

Referenced by analyze(), and CaloTriggerAnalyzer2().

Definition at line 39 of file CaloTriggerAnalyzer2.h.

Referenced by analyze().

TH1F* CaloTriggerAnalyzer2::RPt [private]

Definition at line 128 of file CaloTriggerAnalyzer2.h.

Referenced by analyze(), and CaloTriggerAnalyzer2().

TProfile* CaloTriggerAnalyzer2::RPtEta [private]

Definition at line 132 of file CaloTriggerAnalyzer2.h.

Referenced by analyze(), and CaloTriggerAnalyzer2().

TProfile* CaloTriggerAnalyzer2::RPtEtaFull [private]

Definition at line 133 of file CaloTriggerAnalyzer2.h.

Referenced by analyze(), and CaloTriggerAnalyzer2().

Definition at line 107 of file CaloTriggerAnalyzer2.h.

Referenced by analyze().

Definition at line 125 of file CaloTriggerAnalyzer2.h.

Referenced by analyze(), and CaloTriggerAnalyzer2().

Definition at line 109 of file CaloTriggerAnalyzer2.h.

Referenced by analyze().

Definition at line 127 of file CaloTriggerAnalyzer2.h.

Referenced by analyze(), and CaloTriggerAnalyzer2().

Definition at line 130 of file CaloTriggerAnalyzer2.h.

Referenced by analyze(), and CaloTriggerAnalyzer2().

std::vector<double> CaloTriggerAnalyzer2::SLHCdEta [private]

Definition at line 58 of file CaloTriggerAnalyzer2.h.

Referenced by CaloTriggerAnalyzer2(), clearVectors(), and fillSLHCBranches().

std::vector<double> CaloTriggerAnalyzer2::SLHCdEtac [private]

Definition at line 66 of file CaloTriggerAnalyzer2.h.

Definition at line 102 of file CaloTriggerAnalyzer2.h.

Referenced by analyze(), and CaloTriggerAnalyzer2().

std::vector<double> CaloTriggerAnalyzer2::SLHCdPhi [private]

Definition at line 59 of file CaloTriggerAnalyzer2.h.

Referenced by CaloTriggerAnalyzer2(), clearVectors(), and fillSLHCBranches().

std::vector<double> CaloTriggerAnalyzer2::SLHCdPhic [private]

Definition at line 67 of file CaloTriggerAnalyzer2.h.

std::vector<double> CaloTriggerAnalyzer2::SLHCdPt [private]

Definition at line 57 of file CaloTriggerAnalyzer2.h.

Referenced by CaloTriggerAnalyzer2(), clearVectors(), and fillSLHCBranches().

std::vector<double> CaloTriggerAnalyzer2::SLHCdPtc [private]

Definition at line 65 of file CaloTriggerAnalyzer2.h.

std::vector<double> CaloTriggerAnalyzer2::SLHCdR [private]

Definition at line 56 of file CaloTriggerAnalyzer2.h.

Referenced by CaloTriggerAnalyzer2(), clearVectors(), and fillSLHCBranches().

std::vector<double> CaloTriggerAnalyzer2::SLHCdRc [private]

Definition at line 64 of file CaloTriggerAnalyzer2.h.

Definition at line 94 of file CaloTriggerAnalyzer2.h.

Referenced by analyze(), and CaloTriggerAnalyzer2().

std::vector<double> CaloTriggerAnalyzer2::SLHCL1Eta [private]

Definition at line 54 of file CaloTriggerAnalyzer2.h.

Referenced by CaloTriggerAnalyzer2(), clearVectors(), and fillSLHCBranches().

std::vector<double> CaloTriggerAnalyzer2::SLHCL1Etac [private]

Definition at line 62 of file CaloTriggerAnalyzer2.h.

std::vector<double> CaloTriggerAnalyzer2::SLHCL1Phi [private]

Definition at line 55 of file CaloTriggerAnalyzer2.h.

Referenced by CaloTriggerAnalyzer2(), clearVectors(), and fillSLHCBranches().

std::vector<double> CaloTriggerAnalyzer2::SLHCL1Phic [private]

Definition at line 63 of file CaloTriggerAnalyzer2.h.

std::vector<double> CaloTriggerAnalyzer2::SLHCL1Pt [private]

Definition at line 53 of file CaloTriggerAnalyzer2.h.

Referenced by CaloTriggerAnalyzer2(), clearVectors(), and fillSLHCBranches().

std::vector<double> CaloTriggerAnalyzer2::SLHCL1Ptc [private]

Definition at line 61 of file CaloTriggerAnalyzer2.h.

std::vector<bool> CaloTriggerAnalyzer2::SLHCpassDoubleThresh [private]

Definition at line 70 of file CaloTriggerAnalyzer2.h.

Referenced by analyze(), and clearVectors().

std::vector<bool> CaloTriggerAnalyzer2::SLHCpassSingleThresh [private]

Definition at line 69 of file CaloTriggerAnalyzer2.h.

Referenced by analyze(), CaloTriggerAnalyzer2(), and clearVectors().

std::vector<double> CaloTriggerAnalyzer2::SLHCRPt [private]

Definition at line 60 of file CaloTriggerAnalyzer2.h.

Referenced by CaloTriggerAnalyzer2(), clearVectors(), and fillSLHCBranches().

std::vector<double> CaloTriggerAnalyzer2::SLHCRPtc [private]

Definition at line 68 of file CaloTriggerAnalyzer2.h.

Definition at line 95 of file CaloTriggerAnalyzer2.h.

Referenced by analyze(), and CaloTriggerAnalyzer2().

Definition at line 101 of file CaloTriggerAnalyzer2.h.

Referenced by analyze(), and CaloTriggerAnalyzer2().

Definition at line 40 of file CaloTriggerAnalyzer2.h.

Referenced by analyze().

Definition at line 43 of file CaloTriggerAnalyzer2.h.

Referenced by analyze().