CMS 3D CMS Logo

Public Member Functions | Private Attributes

BasicMuonAnalyzer Class Reference

Example class that can be used both within FWLite and within the full framework. More...

#include <PhysicsTools/UtilAlgos/interface/BasicMuonAnalyzer.h>

Inheritance diagram for BasicMuonAnalyzer:
edm::BasicAnalyzer

List of all members.

Public Member Functions

void analyze (const edm::EventBase &event)
 everything that needs to be done during the event loop
 BasicMuonAnalyzer (const edm::ParameterSet &cfg, TFileDirectory &fs)
 default constructor
void beginJob ()
 everything that needs to be done before the event loop
void endJob ()
 everything that needs to be done after the event loop
virtual ~BasicMuonAnalyzer ()
 default destructor

Private Attributes

std::map< std::string, TH1 * > hists_
 histograms
edm::InputTag muons_
 input tag for mouns

Detailed Description

Example class that can be used both within FWLite and within the full framework.

This is an example for keeping classes that can be used both within FWLite and within the full framework. The class is derived from the BasicAnalyzer base class, which is an interface for the two wrapper classes EDAnalyzerWrapper and FWLiteAnalyzerWrapper. The latter provides basic configuration file reading and event looping equivalent to the FWLiteHistograms executable of this package. You can see the FWLiteAnalyzerWrapper class at work in the FWLiteWithBasicAnalyzer executable of this package.

Definition at line 19 of file BasicMuonAnalyzer.h.


Constructor & Destructor Documentation

BasicMuonAnalyzer::BasicMuonAnalyzer ( const edm::ParameterSet cfg,
TFileDirectory fs 
)

default constructor

Definition at line 8 of file BasicMuonAnalyzer.cc.

References hists_, and TFileDirectory::make().

                                                                                  : 
  edm::BasicAnalyzer::BasicAnalyzer(cfg, fs),
  muons_(cfg.getParameter<edm::InputTag>("muons"))
{
  hists_["muonPt"  ] = fs.make<TH1F>("muonPt"  , "pt"  ,  100,  0., 300.);
  hists_["muonEta" ] = fs.make<TH1F>("muonEta" , "eta" ,  100, -3.,   3.);
  hists_["muonPhi" ] = fs.make<TH1F>("muonPhi" , "phi" ,  100, -5.,   5.); 
  hists_["mumuMass"] = fs.make<TH1F>("mumuMass", "mass",   90, 30., 120.);
}
virtual BasicMuonAnalyzer::~BasicMuonAnalyzer ( ) [inline, virtual]

default destructor

Definition at line 25 of file BasicMuonAnalyzer.h.

{};

Member Function Documentation

void BasicMuonAnalyzer::analyze ( const edm::EventBase event) [virtual]

everything that needs to be done during the event loop

Implements edm::BasicAnalyzer.

Definition at line 20 of file BasicMuonAnalyzer.cc.

References hists_, patZpeak::muons, and muons_.

{
  // define what muon you are using; this is necessary as FWLite is not 
  // capable of reading edm::Views
  using reco::Muon;

  // Handle to the muon collection
  edm::Handle<std::vector<Muon> > muons;
  event.getByLabel(muons_, muons);

  // loop muon collection and fill histograms
  for(std::vector<Muon>::const_iterator mu1=muons->begin(); mu1!=muons->end(); ++mu1){
    hists_["muonPt" ]->Fill( mu1->pt () );
    hists_["muonEta"]->Fill( mu1->eta() );
    hists_["muonPhi"]->Fill( mu1->phi() );
    if( mu1->pt()>20 && fabs(mu1->eta())<2.1 ){
      for(std::vector<Muon>::const_iterator mu2=muons->begin(); mu2!=muons->end(); ++mu2){
        if(mu2>mu1){ // prevent double conting
          if( mu1->charge()*mu2->charge()<0 ){ // check only muon pairs of unequal charge 
            if( mu2->pt()>20 && fabs(mu2->eta())<2.1 ){
              hists_["mumuMass"]->Fill( (mu1->p4()+mu2->p4()).mass() );
            }
          }
        }
      }
    }
  }
}
void BasicMuonAnalyzer::beginJob ( void  ) [inline, virtual]

everything that needs to be done before the event loop

Implements edm::BasicAnalyzer.

Definition at line 27 of file BasicMuonAnalyzer.h.

{};
void BasicMuonAnalyzer::endJob ( void  ) [inline, virtual]

everything that needs to be done after the event loop

Implements edm::BasicAnalyzer.

Definition at line 29 of file BasicMuonAnalyzer.h.

{};

Member Data Documentation

std::map<std::string, TH1*> BasicMuonAnalyzer::hists_ [private]

histograms

Definition at line 37 of file BasicMuonAnalyzer.h.

Referenced by analyze(), and BasicMuonAnalyzer().

input tag for mouns

Definition at line 35 of file BasicMuonAnalyzer.h.

Referenced by analyze().