DESY Hbb Analysis Framework
AnalysisRecoMuonsExample.cc
Go to the documentation of this file.
1 #include <string>
2 #include <iostream>
3 #include <vector>
4 
5 #include "TFile.h"
6 #include "TFileCollection.h"
7 #include "TChain.h"
8 #include "TH1.h"
9 
11 
12 using namespace std;
13 using namespace analysis;
14 using namespace analysis::tools;
15 
16 
17 // =============================================================================================
18 int main(int argc, char * argv[])
19 {
20 // bool isMC = false;
21  TH1::SetDefaultSumw2(); // proper treatment of errors when scaling histograms
22 
23  // Input files list
24  std::string inputList = "rootFileList.txt";
25  Analysis analysis(inputList,"RecoNtuple/Events/EventInfo");
26 
27  // Physics Objects Collections
28  analysis.addTree<RecoMuon> ("Muons","RecoNtuple/Events/muons1Leg");
29 
30 
31  TFile hout("histograms_muons.root","recreate");
32 
33  std::map<std::string, TH1F*> h1;
34  h1["n"] = new TH1F("n" , "" , 30, 0, 30);
35  h1["pt"] = new TH1F("pt" , "" , 100, 0, 770);
36  h1["eta"] = new TH1F("eta" , "" , 100, -5.30, 3.30);
37  h1["phi"] = new TH1F("phi" , "" , 100, -3.7, 3.7);
38  h1["q"] = new TH1F("q", "", 4, -2, 2);
39  h1["e"] = new TH1F("e", "", 100, 0, 10);
40 
41  // Analysis of events
42  std::cout << "This analysis has " << analysis.size() << " events" << std::endl;
43  int nevents = analysis.size();
44 
45 // int nFired = 0;
46  for ( int i = 0 ; i < nevents ; ++i )
47  {
48  if ( i > 0 && i%100000 == 0 ) std::cout << i << " events processed..." << std::endl;
49 
50  analysis.event(i);
51 
52  // Muons
53  auto muons = analysis.collection<RecoMuon>("Muons");
54  int nmuons = 0;
55  for ( int m = 0 ; m < muons->size() ; ++m )
56  {
57  RecoMuon muon = muons->at(m);
58  h1["pt"] -> Fill(muon.pt());
59  h1["eta"] -> Fill(muon.eta());
60  h1["phi"] -> Fill(muon.phi());
61  h1["q"] -> Fill(muon.q());
62  h1["e"] -> Fill(muon.e());
63  ++nmuons;
64 
65  }
66  h1["n"] -> Fill(nmuons);
67 
68  }
69 
70  for (auto & ih1 : h1)
71  {
72  ih1.second -> Write();
73  }
74 
75 //
76 }
77 
std::shared_ptr< Collection< Object > > collection(const std::string &unique_name)
Definition: Analysis.h:336
float eta() const
returns the pseudorapidity
Definition: Candidate.cc:134
float e() const
returns the energy
Definition: Candidate.cc:136
float phi() const
returns the azimuthal angle
Definition: Candidate.cc:135
float pt() const
returns the transverse momentum
Definition: Candidate.cc:133
void event(const int &event, const bool &addCollections=true)
Definition: Analysis.cc:94
std::shared_ptr< PhysicsObjectTree< Object > > addTree(const std::string &unique_name, const std::string &path)
Definition: Analysis.h:273
int main(int argc, char *argv[])
int q() const
returns the charge
Definition: Candidate.cc:139