CMS 3D CMS Logo

ZMuMuOneTrackUserData.cc
Go to the documentation of this file.
10 
16 
20 
21 #include <vector>
22 
23 using namespace edm;
24 using namespace std;
25 using namespace reco;
26 using namespace isodeposit;
27 //using namespace pat;
28 
30 public:
33 
34 private:
35  void produce(edm::Event &, const edm::EventSetup &) override;
36 
41  double alpha_, beta_;
42  string hltPath_;
43  int counter;
44 };
45 
47  : srcToken_(consumes<std::vector<reco::CompositeCandidate> >(cfg.getParameter<InputTag>("src"))),
48  beamSpotToken_(consumes<BeamSpot>(cfg.getParameter<InputTag>("beamSpot"))),
49  primaryVerticesToken_(consumes<VertexCollection>(cfg.getParameter<InputTag>("primaryVertices"))),
50  zGenParticlesMatchToken_(consumes<GenParticleMatch>(cfg.getParameter<InputTag>("zGenParticlesMatch"))),
51  alpha_(cfg.getParameter<double>("alpha")),
52  beta_(cfg.getParameter<double>("beta")),
53  hltPath_(cfg.getParameter<std::string>("hltPath")) {
54  produces<vector<pat::CompositeCandidate> >();
55 }
56 
59  evt.getByToken(srcToken_, dimuons);
60 
61  Handle<BeamSpot> beamSpotHandle;
62  if (!evt.getByToken(beamSpotToken_, beamSpotHandle)) {
63  std::cout << ">>> No beam spot found !!!" << std::endl;
64  }
65 
66  Handle<VertexCollection> primaryVertices; // Collection of primary Vertices
67  if (!evt.getByToken(primaryVerticesToken_, primaryVertices)) {
68  std::cout << ">>> No primary vertices found !!!" << std::endl;
69  }
70 
71  bool isMCMatchTrue = false;
72 
74  if (evt.getByToken(zGenParticlesMatchToken_, zGenParticlesMatch)) {
75  isMCMatchTrue = true;
76  }
77 
78  //cout<<"isMCMatchTrue"<<isMCMatchTrue <<endl;
79  unique_ptr<vector<pat::CompositeCandidate> > dimuonColl(new vector<pat::CompositeCandidate>());
80 
81  for (unsigned int i = 0; i < dimuons->size(); ++i) {
82  const CompositeCandidate &z = (*dimuons)[i];
83  //CandidateBaseRef zRef = dimuons ->refAt(i);
85  pat::CompositeCandidate dimuon(z);
86 
87  float trueMass, truePt, trueEta, truePhi, trueY;
88  if (isMCMatchTrue) {
89  GenParticleRef trueZRef = (*zGenParticlesMatch)[zRef];
90  //CandidateRef trueZRef = trueZIter->val;
91  if (trueZRef.isNonnull()) {
92  const Candidate &z = *trueZRef;
93  trueMass = z.mass();
94  truePt = z.pt();
95  trueEta = z.eta();
96  truePhi = z.phi();
97  trueY = z.rapidity();
98  } else {
99  trueMass = -100;
100  truePt = -100;
101  trueEta = -100;
102  truePhi = -100;
103  trueY = -100;
104  }
105 
106  dimuon.addUserFloat("TrueMass", trueMass);
107  dimuon.addUserFloat("TruePt", truePt);
108  dimuon.addUserFloat("TrueEta", trueEta);
109  dimuon.addUserFloat("TruePhi", truePhi);
110  dimuon.addUserFloat("TrueY", trueY);
111  }
112 
113  dimuonColl->push_back(dimuon);
114  }
115 
116  evt.put(std::move(dimuonColl));
117 }
118 
120 
Analysis-level particle class.
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:131
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:238
EDGetTokenT< GenParticleMatch > zGenParticlesMatchToken_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
ZMuMuUserDataOneTrack(const edm::ParameterSet &)
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
void addUserFloat(const std::string &label, float data, const bool overwrite=false)
Set user-defined float.
Definition: PATObject.h:897
EDGetTokenT< VertexCollection > primaryVerticesToken_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
EDGetTokenT< BeamSpot > beamSpotToken_
EDGetTokenT< std::vector< reco::CompositeCandidate > > srcToken_
virtual double eta() const =0
momentum pseudorapidity
virtual double pt() const =0
transverse momentum
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
virtual double mass() const =0
mass
virtual double rapidity() const =0
rapidity
fixed size matrix
HLT enums.
void produce(edm::Event &, const edm::EventSetup &) override
virtual double phi() const =0
momentum azimuthal angle
def move(src, dest)
Definition: eostools.py:511