CMS 3D CMS Logo

Public Member Functions | Private Attributes

AnalysisTasksAnalyzerJEC Class Reference

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

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

Inheritance diagram for AnalysisTasksAnalyzerJEC:
edm::BasicAnalyzer

List of all members.

Public Member Functions

 AnalysisTasksAnalyzerJEC (const edm::ParameterSet &cfg, TFileDirectory &fs)
 default constructor
void analyze (const edm::EventBase &event)
 everything that needs to be done during the event loop
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 ~AnalysisTasksAnalyzerJEC ()
 default destructor

Private Attributes

bool help_
std::map< std::string, TH2 * > hists_
 histograms
std::string jecLevel_
std::string jecSetLabel_
unsigned int jetInEvents_
edm::InputTag Jets_
 input tag for mouns
std::string outputFileName_
std::string patJetCorrFactors_

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 AnalysisTasksAnalyzerJEC.h.


Constructor & Destructor Documentation

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

default constructor

Definition at line 9 of file AnalysisTasksAnalyzerJEC.cc.

References hists_, jecLevel_, and TFileDirectory::make().

                                                                                                : 
  edm::BasicAnalyzer::BasicAnalyzer(cfg, fs),
  Jets_(cfg.getParameter<edm::InputTag>("Jets")),
  jecLevel_(cfg.getParameter<std::string>("jecLevel")),
  jecSetLabel_(cfg.getParameter<std::string>("jecSetLabel"))
{
  hists_["Response"  ] = fs.make<TH2F>("Response" , "response "+TString(jecLevel_)+"; #eta;p_{T}(reco)/p_{T}(gen)"  ,  5,  -3.3, 3.3, 100, 0.4, 1.6);
}
AnalysisTasksAnalyzerJEC::~AnalysisTasksAnalyzerJEC ( ) [virtual]

default destructor

deconstructor

Definition at line 18 of file AnalysisTasksAnalyzerJEC.cc.

{
}

Member Function Documentation

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

everything that needs to be done during the event loop

Implements edm::BasicAnalyzer.

Definition at line 23 of file AnalysisTasksAnalyzerJEC.cc.

References hists_, jecLevel_, jecSetLabel_, METSkim_cff::Jets, Jets_, and gen::k.

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

  // Handle to the Jet collection
  edm::Handle<std::vector<Jet> > Jets;
  event.getByLabel(Jets_, Jets);

  // loop Jet collection and fill histograms
  for(std::vector<Jet>::const_iterator jet_it=Jets->begin(); jet_it!=Jets->end(); ++jet_it){  


      for(unsigned int k=0; k< jet_it->availableJECSets().size(); ++k){
          edm::LogInfo ("hint1") <<"\n \n *************** HINT 1 *************** \n \n Label of available JEC Set: "<< jet_it->availableJECSets()[k] 
                                 <<"\n \n You found out which JEC Sets were created within the PAT tuple creation! The wrong label caused the segmentation fault in the next for-loop where you ask for JEC levels of a JEC Set that does not exist. Correct the JEC label in your config file and eliminate the segmentation fault. \n  *********************************************** \n"; 
      }                 
      for(unsigned int k=0; k< jet_it->availableJECLevels().size(); ++k){
          edm::LogInfo("hint2")<<" \n \n  Label of available JEC level: "<< jet_it->availableJECLevels(jecSetLabel_)[k] << "\n \n "; 
      }
      edm::LogInfo("hint2")<<"\n \n  *************** HINT 2 ************** \n You did it correctly congratulations!!!!  And you found out above which JEC levels are saved within the jets. We want to investigate these in the following with the response function. At the moment you are trying to access the JEC Level: "
                           << jecLevel_ << ". Does it exist? This causes the error below (see 'Exception Message')! Correct the label of the correction level to an existing one that you are interested in. \n ******************************************* \n " <<std::endl;                   
      if(jet_it->genParticlesSize()>0){
          hists_["Response" ]->Fill( jet_it->correctedJet(jecLevel_).eta(), jet_it->correctedJet(jecLevel_).pt()/ jet_it->genParticle(0)->pt());
      }
  }
}
void AnalysisTasksAnalyzerJEC::beginJob ( void  ) [inline, virtual]

everything that needs to be done before the event loop

Implements edm::BasicAnalyzer.

Definition at line 27 of file AnalysisTasksAnalyzerJEC.h.

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

everything that needs to be done after the event loop

Implements edm::BasicAnalyzer.

Definition at line 29 of file AnalysisTasksAnalyzerJEC.h.

{};

Member Data Documentation

Definition at line 39 of file AnalysisTasksAnalyzerJEC.h.

std::map<std::string, TH2*> AnalysisTasksAnalyzerJEC::hists_ [private]

histograms

Definition at line 43 of file AnalysisTasksAnalyzerJEC.h.

Referenced by AnalysisTasksAnalyzerJEC(), and analyze().

std::string AnalysisTasksAnalyzerJEC::jecLevel_ [private]

Definition at line 36 of file AnalysisTasksAnalyzerJEC.h.

Referenced by AnalysisTasksAnalyzerJEC(), and analyze().

Definition at line 37 of file AnalysisTasksAnalyzerJEC.h.

Referenced by analyze().

Definition at line 41 of file AnalysisTasksAnalyzerJEC.h.

input tag for mouns

Definition at line 35 of file AnalysisTasksAnalyzerJEC.h.

Referenced by analyze().

Definition at line 40 of file AnalysisTasksAnalyzerJEC.h.

Definition at line 38 of file AnalysisTasksAnalyzerJEC.h.