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 
39  genTag_ = pset.getUntrackedParameter<edm::InputTag> ("GenTag", edm::InputTag("genParticles"));
40 
41  produces<double>();
42 }
43 
46 
49 
52 
55 
56  if (iEvent.isRealData()) return;
57 
59  iEvent.getByLabel(genTag_, genParticles);
60  unsigned int gensize = genParticles->size();
61 
62  std::auto_ptr<double> weight (new double);
63 
64  // Set a default weight to start with
65  (*weight) = 1.;
66 
67  // Find the boson at the hard scattering level
68  const reco::GenParticle* boson = 0;
69  int parton1Key = -1;
70  int parton2Key = -1;
71  for (unsigned int i = 0; i<gensize; ++i) {
72  const reco::GenParticle& part = (*genParticles)[i];
73  int status = abs(part.status());
74  if (status!=3) continue;
75  if (part.numberOfMothers()!=2) continue;
76  int partId = abs(part.pdgId());
77  if (status==3 && (partId==23||abs(partId)==24)) {
78  boson = &(*genParticles)[i];
79  parton1Key = part.motherRef(0).key();
80  parton2Key = part.motherRef(1).key();
81  break;
82  }
83  }
84 
85  // Consider only photons near the hard-scattering process
86  const reco::GenParticle* photon = 0;
87  if (boson) {
88  for (unsigned int i = 0; i<gensize; ++i) {
89  photon = 0;
90  const reco::GenParticle& part = (*genParticles)[i];
91  int status = abs(part.status());
92  if (status!=1) continue;
93  int partId = abs(part.pdgId());
94  if (partId!=22) continue;
95  if (part.numberOfMothers()!=1) continue;
96  int keyM = part.motherRef(0).key();
97  const reco::GenParticle* mother = &(*genParticles)[keyM];
98  if (mother->status()!=3) continue;
99  int mId = mother->pdgId();
100  if (abs(mId)>6 && mId!=2212) continue;
101  for (unsigned int j=0; j<mother->numberOfDaughters(); ++j){
102  int keyD = mother->daughterRef(j).key();
103  if (keyD==parton1Key || keyD==parton2Key) {
104  photon = &part;
105  break;
106  }
107  }
108  if (photon) break;
109  }
110  }
111 
112  if (boson && photon) {
113  math::XYZTLorentzVector smom = boson->p4() + photon->p4();
114  double s = smom.M2();
115  double sqrts = smom.M();
116 
117  // Go to CM using the boost direction of the boson+photon system
118  ROOT::Math::Boost cmboost(smom.BoostToCM());
119  math::XYZTLorentzVector photonCM(cmboost(photon->p4()));
120  double pcostheta = ( smom.x()*photonCM.x()
121  + smom.y()*photonCM.y()
122  + smom.z()*photonCM.z() ) / smom.P();
123 
124  // Determine kinematic invariants
125  double t = - sqrts * (photonCM.t()-pcostheta);
126  double MV = boson->mass();
127  double u = MV*MV - s - t;
128  (*weight) = 1. - 2*t*u/(s*s+MV*MV*MV*MV);
129  //printf(">>>>>>>>> s %f t %f u %f, MV %f, weight = %f\n", s, t, u, MV, (*weight));
130  }
131 
132  iEvent.put(weight);
133 }
134 
virtual void endJob() override
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
daughters::value_type daughterRef(size_type i) const
reference to daughter at given position
virtual const LorentzVector & p4() const GCC11_FINAL
four-momentum Lorentz vector
virtual int pdgId() const GCC11_FINAL
PDG identifier.
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
virtual void beginJob() override
bool isRealData() const
Definition: EventBase.h:60
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
int iEvent
Definition: GenABIO.cc:243
virtual int status() const GCC11_FINAL
status word
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
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
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:390
virtual float mass() const GCC11_FINAL
mass
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