CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ZMuMuSaMassHistogram.cc
Go to the documentation of this file.
23 #include "TH1.h"
24 #include "TH2.h"
25 #include <vector>
26 #include <string>
27 #include <iostream>
28 #include <sstream>
29 
30 using namespace edm;
31 using namespace std;
32 using namespace reco;
33 using namespace isodeposit;
34 
35 
36 
38 public:
41 private:
42  virtual void analyze(const edm::Event& event, const edm::EventSetup& setup) override;
43  virtual void endJob() override;
45  int counter;
46  double min, max;
47  int Nbin;
48  TH1F * ZMassSa;
49  void histo(const TH1F* hist, char* cx, char* cy) const;
50 };
51 
52 void ZMuMuSaMassHistogram::histo(const TH1F* hist,char* cx, char*cy) const{
53  hist->GetXaxis()->SetTitle(cx);
54  hist->GetYaxis()->SetTitle(cy);
55  hist->GetXaxis()->SetTitleOffset(1);
56  hist->GetYaxis()->SetTitleOffset(1.2);
57  hist->GetXaxis()->SetTitleSize(0.04);
58  hist->GetYaxis()->SetTitleSize(0.04);
59  hist->GetXaxis()->SetLabelSize(0.03);
60  hist->GetYaxis()->SetLabelSize(0.03);
61 }
62 
63 
65  srcToken_(consumes<CandidateView>(pset.getParameter<InputTag>("src_m"))),
66  counter(0),
67  min(pset.getUntrackedParameter<double>("min")),
68  max(pset.getUntrackedParameter<double>("max")),
69  Nbin(pset.getUntrackedParameter<int>("nbin")) {
71  ZMassSa = fs->make<TH1F>("zMass","ZMass OneStandAlone (GeV/c^{2})",Nbin,min,max);
72 
73 }
74 
75 
78  event.getByToken(srcToken_, dimuons);
79  for(unsigned int i=0; i< dimuons->size(); ++ i ) {
80  const Candidate & zmm = (* dimuons)[i];
81  const Candidate * dau0 = zmm.daughter(0);
82  const Candidate * dau1 = zmm.daughter(1);
83  TrackRef stAloneTrack;
85  double mu_mass;
86  if(counter % 2 == 0) {
87  stAloneTrack = dau0->get<TrackRef,reco::StandAloneMuonTag>();
88  p4_0 = dau1->polarP4();
89  mu_mass = dau0->mass();
90  }
91  else{
92  stAloneTrack = dau1->get<TrackRef,reco::StandAloneMuonTag>();
93  p4_0= dau0->polarP4();
94  mu_mass = dau1->mass();
95  }
96 
97  Vector momentum = stAloneTrack->momentum();
98  Candidate::PolarLorentzVector p4_1(momentum.rho(), momentum.eta(),momentum.phi(), mu_mass);
99  double mass = (p4_0+p4_1).mass();
100  ZMassSa->Fill(mass);
101  ++counter;
102 
103  }
104 
105 
106 }
107 
108 
110 }
111 
113 
int i
Definition: DBlmapReader.cc:9
virtual const Candidate * daughter(size_type i) const =0
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode) ...
virtual float mass() const =0
mass
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
void histo(const TH1F *hist, char *cx, char *cy) const
stand alone muon component tag
Definition: RecoCandidate.h:72
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
virtual const PolarLorentzVector & polarP4() const =0
four-momentum Lorentz vector
EDGetTokenT< CandidateView > srcToken_
virtual void endJob() override
EventID const & min(EventID const &lh, EventID const &rh)
Definition: EventID.h:132
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
ZMuMuSaMassHistogram(const edm::ParameterSet &pset)
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
edm::Ref< TrackCollection > TrackRef
persistent reference to a Track
Definition: TrackFwd.h:14
virtual void analyze(const edm::Event &event, const edm::EventSetup &setup) override
static std::atomic< unsigned int > counter
T get() const
get a component
Definition: Candidate.h:219
void setup(std::vector< TH2F > &depth, std::string name, std::string units="")
EventID const & max(EventID const &lh, EventID const &rh)
Definition: EventID.h:137
math::PtEtaPhiMLorentzVector PolarLorentzVector
Lorentz vector.
Definition: Candidate.h:43