CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ISRGammaWeightProducer.cc
Go to the documentation of this file.
1 #include <memory>
2 
5 
8 
10 
15 #include <Math/VectorUtil.h>
16 
17 //
18 // class declaration
19 //
21  public:
24 
25  virtual void beginJob() override ;
26  virtual void produce(edm::Event&, const edm::EventSetup&) override;
27  virtual void endJob() override ;
28 
29  private:
31  std::vector<double> isrBinEdges_;
32  std::vector<double> ptWeights_;
33 };
34 
35 
38  genToken_ = consumes<reco::GenParticleCollection>(pset.getUntrackedParameter<edm::InputTag> ("GenTag", edm::InputTag("genParticles")));
39 
40  produces<double>();
41 }
42 
45 
48 
51 
54 
55  if (iEvent.isRealData()) return;
56 
58  iEvent.getByToken(genToken_, genParticles);
59  unsigned int gensize = genParticles->size();
60 
61  std::auto_ptr<double> weight (new double);
62 
63  // Set a default weight to start with
64  (*weight) = 1.;
65 
66  // Find the boson at the hard scattering level
67  const reco::GenParticle* boson = 0;
68  int parton1Key = -1;
69  int parton2Key = -1;
70  for (unsigned int i = 0; i<gensize; ++i) {
71  const reco::GenParticle& part = (*genParticles)[i];
72  int status = abs(part.status());
73  if (status!=3) continue;
74  if (part.numberOfMothers()!=2) continue;
75  int partId = abs(part.pdgId());
76  if (status==3 && (partId==23||abs(partId)==24)) {
77  boson = &(*genParticles)[i];
78  parton1Key = part.motherRef(0).key();
79  parton2Key = part.motherRef(1).key();
80  break;
81  }
82  }
83 
84  // Consider only photons near the hard-scattering process
85  const reco::GenParticle* photon = 0;
86  if (boson) {
87  for (unsigned int i = 0; i<gensize; ++i) {
88  photon = 0;
89  const reco::GenParticle& part = (*genParticles)[i];
90  int status = abs(part.status());
91  if (status!=1) continue;
92  int partId = abs(part.pdgId());
93  if (partId!=22) continue;
94  if (part.numberOfMothers()!=1) continue;
95  int keyM = part.motherRef(0).key();
96  const reco::GenParticle* mother = &(*genParticles)[keyM];
97  if (mother->status()!=3) continue;
98  int mId = mother->pdgId();
99  if (abs(mId)>6 && mId!=2212) continue;
100  for (unsigned int j=0; j<mother->numberOfDaughters(); ++j){
101  int keyD = mother->daughterRef(j).key();
102  if (keyD==parton1Key || keyD==parton2Key) {
103  photon = &part;
104  break;
105  }
106  }
107  if (photon) break;
108  }
109  }
110 
111  if (boson && photon) {
112  math::XYZTLorentzVector smom = boson->p4() + photon->p4();
113  double s = smom.M2();
114  double sqrts = smom.M();
115 
116  // Go to CM using the boost direction of the boson+photon system
117  ROOT::Math::Boost cmboost(smom.BoostToCM());
118  math::XYZTLorentzVector photonCM(cmboost(photon->p4()));
119  double pcostheta = ( smom.x()*photonCM.x()
120  + smom.y()*photonCM.y()
121  + smom.z()*photonCM.z() ) / smom.P();
122 
123  // Determine kinematic invariants
124  double t = - sqrts * (photonCM.t()-pcostheta);
125  double MV = boson->mass();
126  double u = MV*MV - s - t;
127  (*weight) = 1. - 2*t*u/(s*s+MV*MV*MV*MV);
128  //printf(">>>>>>>>> s %f t %f u %f, MV %f, weight = %f\n", s, t, u, MV, (*weight));
129  }
130 
131  iEvent.put(weight);
132 }
133 
virtual void endJob() override
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
virtual int pdgId() const
PDG identifier.
edm::EDGetTokenT< reco::GenParticleCollection > genToken_
daughters::value_type daughterRef(size_type i) const
reference to daughter at given position
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:449
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
virtual int status() const
status word
virtual void beginJob() override
bool isRealData() const
Definition: EventBase.h:64
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
virtual double mass() const
mass
int iEvent
Definition: GenABIO.cc:230
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:113
virtual size_t numberOfMothers() const
number of mothers
virtual void produce(edm::Event &, const edm::EventSetup &) override
ISRGammaWeightProducer(const edm::ParameterSet &)
virtual size_t numberOfDaughters() const
number of daughters
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int j
Definition: DBlmapReader.cc:9
daughters::value_type motherRef(size_type i=0) const
reference to mother at given position
part
Definition: HCALResponse.h:20
std::vector< double > isrBinEdges_
std::vector< double > ptWeights_
int weight
Definition: histoStyle.py:50
tuple status
Definition: ntuplemaker.py:245
virtual const LorentzVector & p4() const
four-momentum Lorentz vector
Definition: LeafCandidate.h:99