CMS 3D CMS Logo

List of all members | Public Types | Public Member Functions | Protected Types | Protected Member Functions | Protected Attributes
Py8toJetInput Class Reference

#include <Py8toJetInput.h>

Inheritance diagram for Py8toJetInput:
Py8toJetInputHEPEVT

Public Types

typedef Pythia8::Event Event
 
typedef Pythia8::Particle Particle
 

Public Member Functions

virtual const std::vector< fastjet::PseudoJet > fillJetAlgoInput (const Event &, const Event &, const lhef::LHEEvent *lhee=0, const std::vector< int > *partonList=0)
 
 Py8toJetInput ()
 
void setJetEtaMax (double max)
 
virtual ~Py8toJetInput ()
 

Protected Types

enum  partonTypes { ID_TOP =6, ID_GLUON =21, ID_PHOTON =22 }
 

Protected Member Functions

int getAncestor (int, const Event &, const Event &)
 

Protected Attributes

double fJetEtaMax
 
std::vector< fastjet::PseudoJet > fJetInput
 

Detailed Description

Definition at line 13 of file Py8toJetInput.h.

Member Typedef Documentation

typedef Pythia8::Event Py8toJetInput::Event

Definition at line 16 of file Py8toJetInput.h.

typedef Pythia8::Particle Py8toJetInput::Particle

Definition at line 17 of file Py8toJetInput.h.

Member Enumeration Documentation

Enumerator
ID_TOP 
ID_GLUON 
ID_PHOTON 

Definition at line 29 of file Py8toJetInput.h.

Constructor & Destructor Documentation

Py8toJetInput::Py8toJetInput ( )
inline

Definition at line 19 of file Py8toJetInput.h.

19 : fJetEtaMax(10.) {}
double fJetEtaMax
Definition: Py8toJetInput.h:30
virtual Py8toJetInput::~Py8toJetInput ( )
inlinevirtual

Definition at line 20 of file Py8toJetInput.h.

20 {}

Member Function Documentation

const std::vector< fastjet::PseudoJet > Py8toJetInput::fillJetAlgoInput ( const Event event,
const Event workEvent,
const lhef::LHEEvent lhee = 0,
const std::vector< int > *  partonList = 0 
)
virtual

Reimplemented in Py8toJetInputHEPEVT.

Definition at line 8 of file Py8toJetInput.cc.

References end, PVValHelper::eta, spr::find(), fJetEtaMax, fJetInput, lhef::LHEEvent::getHEPEUP(), mps_fire::i, ID_PHOTON, ID_TOP, training_settings::idx, lhef::HEPEUP::MOTHUP, edm::shift, and mps_update::status.

Referenced by JetMatchingHook::doVetoPartonLevel().

11 {
12 
13  fJetInput.clear();
14 
15  Event workEventJet = workEvent;
16 
17  const lhef::HEPEUP& hepeup = *lhee->getHEPEUP();
18 
19  std::set< int > typeSet[3];
20 
21  // FIXME !!!
22  // This is not safe to assume it's 3 because we're passing in a pointer
23  // and we do NOT know what the actuial size is. I'll have to improve it.
24  //
25  for ( size_t i=0; i<3; i++ )
26  {
27  typeSet[i].clear();
28  for ( size_t j=0; j<typeIdx[i].size(); j++ )
29  {
30 
31  // HEPEUP and Py8 event record are shifted by 3
32  // (system particle + 2 beam particles in Py8 event)
33  // ... EXCEPT FOR THE DECAY PRODUCTS IF DONE AT THE ME LEVEL !!!!!
34  // ... in such case, the decay productes get gutted out of the sequence
35  // and get placed in Py8 Event later on...
36  // so we need to figure out the shift
37  int pos = typeIdx[i][j];
38  int shift = 3;
39  // skip the first 2 entrirs in HEPEUP - they're incoming partons
40  for ( int ip=2; ip<pos; ip++ )
41  {
42  // alternative can be: moth1 != 1 && moth2 !=2...
43  // but moth1 == moth2 means pointer to the same mother t
44  // that can only be a resonance, unless moth1==moth2==0
45  //
46  if ( hepeup.MOTHUP[ip].first == hepeup.MOTHUP[ip].second )
47  {
48  shift--;
49  }
50  }
51  pos += shift;
52  // typeSet[i].insert( event[pos].daughter1() );
53  typeSet[i].insert( pos );
54  }
55  }
56 
57  // int iTypeEnd = (typeIdx[1].empty()) ? 1 : 2;
58 
59 // --> FIXME !!!
60  int iType = 0; // only LIGHT jets for now
61  int jetAllow = 0; // hardcoded for now for the same value as is set in Steve's example
62  // at present, not even in use...
63  // int jetMatch = 0; // hardcoded for now for the same value as is set in Steve's example
64 
65  // Loop over particles and decide what to pass to the jet algorithm
66  for (int i = 0; i < workEventJet.size(); ++i)
67  {
68 
69  if (!workEventJet[i].isFinal()) continue;
70 
71  // jetAllow option to disallow certain particle types
72  if (jetAllow == 1)
73  {
74 
75  // Original AG+Py6 algorithm explicitly excludes tops,
76  // leptons and photons.
77  int id = workEventJet[i].idAbs();
78  if ((id >= 11 && id <= 16) || id == ID_TOP || id == ID_PHOTON)
79  {
80  workEventJet[i].statusNeg();
81  continue;
82  }
83  }
84 
85  // Get the index of this particle in original event
86  int idx = workEventJet[i].daughter1();
87 
88  // Start with particle idx, and afterwards track mothers
89  while (true)
90  {
91 
92  // Light jets
93  if (iType == 0)
94  {
95 
96  // Do not include if originates from heavy jet or 'other'
97  if (typeSet[1].find(idx) != typeSet[1].end() ||
98  typeSet[2].find(idx) != typeSet[2].end())
99  {
100  workEventJet[i].statusNeg();
101  break;
102  }
103 
104  // Made it to start of event record so done
105  if (idx == 0)
106  {
107  break;
108  }
109  // Otherwise next mother and continue
110  idx = event[idx].mother1();
111 
112  // Heavy jets
113  }
114  else if (iType == 1)
115  {
116 
117  // Only include if originates from heavy jet
118  if (typeSet[1].find(idx) != typeSet[1].end()) break;
119 
120  // Made it to start of event record with no heavy jet mother,
121  // so DO NOT include particle
122  if (idx == 0)
123  {
124  workEventJet[i].statusNeg();
125  break;
126  }
127 
128  // Otherwise next mother and continue
129  idx = event[idx].mother1();
130 
131  } // if (iType)
132  } // while (true)
133  } // for (i)
134 
135  // For jetMatch = 2, insert ghost particles corresponding to
136  // each hard parton in the original process
137 /*
138  if (jetMatch > 0)
139  {
140 
141  for (int i = 0; i < int(typeIdx[iType].size()); i++)
142  {
143  // Get y/phi of the parton
144  Vec4 pIn = eventProcess[typeIdx[iType][i]].p();
145  double y = Vec4y(pIn);
146  double phi = pIn.phi();
147 
148  // Create a ghost particle and add to the workEventJet
149  double e = MG5_GHOSTENERGY;
150  double e2y = exp(2. * y);
151  double pz = e * (e2y - 1.) / (e2y + 1.);
152  double pt = sqrt(e*e - pz*pz);
153  double px = pt * cos(phi);
154  double py = pt * sin(phi);
155  workEventJet.append(Particle(ID_GLUON, 99, 0, 0, 0, 0, 0, 0,
156  px, py, pz, e, 0., 0, 9.));
157 
158  } // for (i)
159  } // if (jetMatch == 2)
160 */
161 
162  for ( int i=0; i<workEventJet.size(); i++ )
163  {
164 
165  // fisrt, weed out all entries marked with statusNeg();
166  //
167  if ( workEventJet[i].status() < 0 ) continue;
168 
169 
170  // now weed out all entries above etaMax
171  // in principle, we can use etaJetMaxAlgo because it's set equal to etaJetMax
172  // as for etaJetMax, it gets set to memain_.etaclmax
173  //
174  if ( fabs(workEventJet[i].eta()) > fJetEtaMax ) continue ;
175 
176  // need to double check if native FastJet understands Py8 Event format
177  // in general, PseudoJet gets formed from (px,py,pz,E)
178  //
179  fastjet::PseudoJet partTmp = workEventJet[i];
180  fJetInput.push_back( partTmp );
181  }
182 
183  return fJetInput;
184 
185 }
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
const HEPEUP * getHEPEUP() const
Definition: LHEEvent.h:43
std::vector< std::pair< int, int > > MOTHUP
Definition: LesHouches.h:236
#define end
Definition: vmac.h:37
double fJetEtaMax
Definition: Py8toJetInput.h:30
static unsigned int const shift
std::vector< fastjet::PseudoJet > fJetInput
Definition: Py8toJetInput.h:34
int Py8toJetInput::getAncestor ( int  pos,
const Event fullEvent,
const Event workEvent 
)
protected

Definition at line 188 of file Py8toJetInput.cc.

References funct::abs(), counter, and mps_update::status.

Referenced by Py8toJetInputHEPEVT::fillJetAlgoInput().

189 {
190 
191  int parentId = fullEvent[pos].mother1();
192  int parentPrevId = 0;
193  int counter = pos;
194 
195  while ( parentId > 0 )
196  {
197  if ( parentId == fullEvent[counter].mother2() ) // carbon copy, keep walking up
198  {
199  parentPrevId = parentId;
200  counter = parentId;
201  parentId = fullEvent[parentPrevId].mother1();
202  continue;
203  }
204 
205  // we get here if not a carbon copy
206 
207  // let's check if it's a normal process, etc.
208  //
209  if ( (parentId < parentPrevId) || parentId < fullEvent[counter].mother2() ) // normal process
210  {
211 
212  // first of all, check if hard block
213  if ( abs(fullEvent[counter].status()) == 22 || abs(fullEvent[counter].status()) == 23 )
214  {
215  // yes, it's the hard block
216  // we got what we want, and can exit now !
217  parentId = counter;
218  break;
219  }
220  else
221  {
222  parentPrevId = parentId;
223  parentId = fullEvent[parentPrevId].mother1();
224  }
225  }
226  else if ( parentId > parentPrevId || parentId > pos ) // "circular"/"forward-pointing parent" - intermediate process
227  {
228  parentId = -1;
229  break;
230  }
231 
232  // additional checks... although we shouldn't be geeting here all that much...
233  //
234  if ( abs(fullEvent[parentId].status()) == 22 || abs(fullEvent[parentId].status())== 23 ) // hard block
235  {
236  break;
237  }
238  if ( abs(fullEvent[parentId].status()) < 22 ) // incoming
239  {
240  parentId = -1;
241  break;
242  }
243  }
244 
245  return parentId;
246 
247 }
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
static std::atomic< unsigned int > counter
void Py8toJetInput::setJetEtaMax ( double  max)
inline

Definition at line 25 of file Py8toJetInput.h.

References hpstanc_transforms::max.

Referenced by JetMatchingHook::init().

Member Data Documentation

double Py8toJetInput::fJetEtaMax
protected

Definition at line 30 of file Py8toJetInput.h.

Referenced by fillJetAlgoInput().

std::vector<fastjet::PseudoJet> Py8toJetInput::fJetInput
protected

Definition at line 34 of file Py8toJetInput.h.

Referenced by fillJetAlgoInput(), and Py8toJetInputHEPEVT::fillJetAlgoInput().