StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
History.h
1 // History.h is a part of the PYTHIA event generator.
2 // Copyright (C) 2020 Torbjorn Sjostrand.
3 // PYTHIA is licenced under the GNU GPL v2 or later, see COPYING for details.
4 // Please respect the MCnet Guidelines, see GUIDELINES for details.
5 
6 // This file is written by Stefan Prestel.
7 // It contains the main class for matrix element merging.
8 // Header file for the Clustering and History classes.
9 
10 #ifndef Pythia8_History_H
11 #define Pythia8_History_H
12 
13 #include "Pythia8/Basics.h"
14 #include "Pythia8/BeamParticle.h"
15 #include "Pythia8/Event.h"
16 #include "Pythia8/Info.h"
17 #include "Pythia8/ParticleData.h"
18 #include "Pythia8/PartonLevel.h"
19 #include "Pythia8/PythiaStdlib.h"
20 #include "Pythia8/Settings.h"
21 #include "Pythia8/StandardModel.h"
22 #include "Pythia8/SimpleWeakShowerMEs.h"
23 
24 namespace Pythia8 {
25 
26 //==========================================================================
27 
28 // Declaration of Clustering class.
29 // This class holds information about one radiator, recoiler, emitted system.
30 // This class is a container class for History class use.
31 
32 class Clustering {
33 
34 public:
35 
36  // The emitted parton location.
37  int emitted;
38  // The emittor parton
39  int emittor;
40  // The recoiler parton.
41  int recoiler;
42  // The colour connected recoiler (Can be different for ISR)
43  int partner;
44  // The scale associated with this clustering.
45  double pTscale;
46  // The flavour of the radiator prior to the emission.
47  int flavRadBef;
48  // Spin of the radiator (-1 is left handed, +1 is right handed).
49  int spinRad;
50  // Spin of the emitted (-1 is left handed, +1 is right handed).
51  int spinEmt;
52  // Spin of the recoiler (-1 is left handed, +1 is right handed).
53  int spinRec;
54  // Spin of the radiator before the splitting.
55  int spinRadBef;
56  // The radiator before the splitting.
57  int radBef;
58  // The recoiler before the splitting.
59  int recBef;
60 
61  // Default constructor
62  Clustering() : emitted(0), emittor(0), recoiler(0), partner(0), pTscale(),
63  flavRadBef(0), spinRad(9), spinEmt(9), spinRec(9), spinRadBef(9),
64  radBef(0), recBef(0) {}
65 
66  // Default destructor
67  ~Clustering(){}
68 
69  // Copy constructor
70  Clustering( const Clustering& inSystem ) :
71  emitted(inSystem.emitted), emittor(inSystem.emittor),
72  recoiler(inSystem.recoiler), partner(inSystem.partner),
73  pTscale(inSystem.pTscale), flavRadBef(inSystem.flavRadBef),
74  spinRad(inSystem.spinRad), spinEmt(inSystem.spinEmt),
75  spinRec(inSystem.spinRec), spinRadBef(inSystem.spinRad),
76  radBef(inSystem.radBef), recBef(inSystem.recBef) {}
77 
78  // Constructor with input
79  Clustering( int emtIn, int radIn, int recIn, int partnerIn,
80  double pTscaleIn, int flavRadBefIn = 0, int spinRadIn = 9,
81  int spinEmtIn = 9, int spinRecIn = 9, int spinRadBefIn = 9,
82  int radBefIn = 0, int recBefIn = 0)
83  : emitted(emtIn), emittor(radIn), recoiler(recIn), partner(partnerIn),
84  pTscale(pTscaleIn), flavRadBef(flavRadBefIn), spinRad(spinRadIn),
85  spinEmt(spinEmtIn), spinRec(spinRecIn), spinRadBef(spinRadBefIn),
86  radBef(radBefIn), recBef(recBefIn) {}
87 
88  // Function to return pythia pT scale of clustering
89  double pT() const { return pTscale; }
90 
91  // print for debug
92  void list() const;
93 
94 };
95 
96 //==========================================================================
97 
98 // Declaration of History class
99 //
100 // A History object represents an event in a given step in the CKKW-L
101 // clustering procedure. It defines a tree-like recursive structure,
102 // where the root node represents the state with n jets as given by
103 // the matrix element generator, and is characterized by the member
104 // variable mother being null. The leaves on the tree corresponds to a
105 // fully clustered paths where the original n-jets has been clustered
106 // down to the Born-level state. Also states which cannot be clustered
107 // down to the Born-level are possible - these will be called
108 // incomplete. The leaves are characterized by the vector of children
109 // being empty.
110 
111 class History {
112 
113 public:
114 
115  // The only constructor. Default arguments are used when creating
116  // the initial history node. The \a depth is the maximum number of
117  // clusterings requested. \a scalein is the scale at which the \a
118  // statein was clustered (should be set to the merging scale for the
119  // initial history node. \a beamAIn and beamBIn are needed to
120  // calcutate PDF ratios, \a particleDataIn to have access to the
121  // correct masses of particles. If \a isOrdered is true, the previous
122  // clusterings has been ordered. \a is the PDF ratio for this
123  // clustering (=1 for FSR clusterings). \a probin is the accumulated
124  // probabilities for the previous clusterings, and \ mothin is the
125  // previous history node (null for the initial node).
126  History( int depthIn,
127  double scalein,
128  Event statein,
129  Clustering c,
130  MergingHooksPtr mergingHooksPtrIn,
131  BeamParticle beamAIn,
132  BeamParticle beamBIn,
133  ParticleData* particleDataPtrIn,
134  Info* infoPtrIn,
135  PartonLevel* showersIn,
136  CoupSM* coupSMPtrIn,
137  bool isOrdered,
138  bool isStronglyOrdered,
139  bool isAllowed,
140  bool isNextInInput,
141  double probin,
142  History * mothin);
143 
144  // The destructor deletes each child.
145  ~History() {
146  for ( int i = 0, N = children.size(); i < N; ++i ) delete children[i];
147  }
148 
149  // Function to project paths onto desired paths.
150  bool projectOntoDesiredHistories();
151 
152  // For CKKW-L, NL3 and UMEPS:
153  // In the initial history node, select one of the paths according to
154  // the probabilities. This function should be called for the initial
155  // history node.
156  // IN trialShower* : Previously initialised trialShower object,
157  // to perform trial showering and as
158  // repository of pointers to initialise alphaS
159  // PartonSystems* : PartonSystems object needed to initialise
160  // shower objects
161  // OUT vector<double> : (Sukadov) , (alpha_S ratios) , (PDF ratios)
162  vector<double> weightCKKWL(PartonLevel* trial, AlphaStrong * asFSR,
163  AlphaStrong * asISR, AlphaEM * aemFSR, AlphaEM * aemISR, double RN);
164 
165 
166  // For default NL3:
167  // Return weight of virtual correction and subtractive for NL3 merging
168  vector<double> weightNL3Loop(PartonLevel* trial, double RN);
169  // Return O(\alpha_s)-term of CKKWL-weight for NL3 merging
170  vector<double> weightNL3First(PartonLevel* trial, AlphaStrong* asFSR,
171  AlphaStrong * asISR, AlphaEM * aemFSR, AlphaEM * aemISR, double RN,
172  Rndm* rndmPtr);
173  vector<double> weightNL3Tree(PartonLevel* trial, AlphaStrong * asFSR,
174  AlphaStrong * asISR, AlphaEM * aemFSR, AlphaEM * aemISR, double RN);
175 
176 
177  // For UMEPS:
178  vector<double> weightUMEPSTree(PartonLevel* trial, AlphaStrong * asFSR,
179  AlphaStrong * asISR, AlphaEM * aemFSR, AlphaEM * aemISR, double RN);
180  vector<double> weightUMEPSSubt(PartonLevel* trial, AlphaStrong * asFSR,
181  AlphaStrong * asISR, AlphaEM * aemFSR, AlphaEM * aemISR, double RN);
182 
183 
184  // For unitary NL3:
185  vector<double> weightUNLOPSTree(PartonLevel* trial, AlphaStrong * asFSR,
186  AlphaStrong * asISR, AlphaEM * aemFSR, AlphaEM * aemISR, double RN,
187  int depthIn = -1);
188  vector<double> weightUNLOPSSubt(PartonLevel* trial, AlphaStrong * asFSR,
189  AlphaStrong * asISR, AlphaEM * aemFSR, AlphaEM * aemISR, double RN,
190  int depthIn = -1);
191  vector<double> weightUNLOPSLoop(PartonLevel* trial, AlphaStrong * asFSR,
192  AlphaStrong * asISR, AlphaEM * aemFSR, AlphaEM * aemISR, double RN,
193  int depthIn = -1);
194  vector<double> weightUNLOPSSubtNLO(PartonLevel* trial, AlphaStrong * asFSR,
195  AlphaStrong * asISR, AlphaEM * aemFSR, AlphaEM * aemISR, double RN,
196  int depthIn = -1);
197  vector<double> weightUNLOPSFirst( int order, PartonLevel* trial,
198  AlphaStrong* asFSR, AlphaStrong * asISR, AlphaEM * aemFSR,
199  AlphaEM * aemISR, double RN, Rndm* rndmPtr );
200 
201  // Function to check if any allowed histories were found
202  bool foundAllowedHistories() {
203  return (children.size() > 0 && foundAllowedPath); }
204  // Function to check if any ordered histories were found
205  bool foundOrderedHistories() {
206  return (children.size() > 0 && foundOrderedPath); }
207  // Function to check if any ordered histories were found
208  bool foundCompleteHistories() {
209  return (children.size() > 0 && foundCompletePath); }
210 
211  // Function to set the state with complete scales for evolution
212  void getStartingConditions( const double RN, Event& outState );
213  // Function to get the state with complete scales for evolution
214  bool getClusteredEvent( const double RN, int nSteps, Event& outState);
215  // Function to get the first reclustered state above the merging scale.
216  bool getFirstClusteredEventAboveTMS( const double RN, int nDesired,
217  Event& process, int & nPerformed, bool updateProcess = true );
218  // Function to return the depth of the history (i.e. the number of
219  // reclustered splittings)
220  int nClusterings();
221 
222  // Function to get the lowest multiplicity reclustered event
223  Event lowestMultProc( const double RN) {
224  // Return lowest multiplicity state
225  return (select(RN)->state);
226  }
227 
228  // Calculate and return pdf ratio
229  double getPDFratio( int side, bool forSudakov, bool useHardPDF,
230  int flavNum, double xNum, double muNum,
231  int flavDen, double xDen, double muDen);
232 
233  // Envelope function that calls the recursive getWeakProb.
234  double getWeakProb();
235 
236  // Recursive function that returns the weak probability for the given path.
237  // Mode refers to which ME correction to use, 1 = sChannel, 2 = gluon
238  // channel, 3 = double quark t-channel, 4 is double quark u-channel.
239  double getWeakProb(vector<int>& mode, vector<Vec4>& mom,
240  vector<int> fermionLines);
241 
242  // return the weak probability of a single reclustering.
243  // Mode refers to which ME correction to use, 1 = sChannel, 2 = gluon
244  // channel, 3 = double quark t-channel, 4 is double quark u-channel.
245  double getSingleWeakProb(vector<int> &mode, vector<Vec4> &mom,
246  vector<int> fermionLines);
247 
248  // Find map between indecies in the current state and the state after
249  // the splitting.
250  // NOT IMPLEMENTED FOR MULTIPLE W/Z/GAMMA (NEED TO HAVE A WAY TO
251  // IDENTIFY THEM).
252  void findStateTransfer(map<int,int> &transfer);
253 
254  // Function to print the history that would be chosen from the random number
255  // RN. Mainly for debugging.
256  void printHistory( const double RN );
257  // Function to print the states in a history, starting from the hard process.
258  // Mainly for debugging.
259  void printStates();
260 
261  // Make Pythia class friend
262  friend class Pythia;
263  // Make Merging class friend
264  friend class Merging;
265 
266 private:
267 
268  // Number of trial emission to use for calculating the average number of
269  // emissions
270  static const int NTRIAL;
271 
272  // Function to set all scales in the sequence of states. This is a
273  // wrapper routine for setScales and setEventScales methods
274  void setScalesInHistory();
275 
276  // Function to find the index (in the mother histories) of the
277  // child history, thus providing a way access the path from both
278  // initial history (mother == 0) and final history (all children == 0)
279  // IN vector<int> : The index of each child in the children vector
280  // of the current history node will be saved in
281  // this vector
282  // NO OUTPUT
283  void findPath(vector<int>& out);
284 
285  // Functions to set the parton production scales and enforce
286  // ordering on the scales of the respective clusterings stored in
287  // the History node:
288  // Method will start from lowest multiplicity state and move to
289  // higher states, setting the production scales the shower would
290  // have used.
291  // When arriving at the highest multiplicity, the method will switch
292  // and go back in direction of lower states to check and enforce
293  // ordering for unordered histories.
294  // IN vector<int> : Vector of positions of the chosen child
295  // in the mother history to allow to move
296  // in direction initial->final along path
297  // bool : True: Move in direction low->high
298  // multiplicity and set production scales
299  // False: Move in direction high->low
300  // multiplicity and check and enforce
301  // ordering
302  // NO OUTPUT
303  void setScales( vector<int> index, bool forward);
304 
305  // Function to find a particle in all higher multiplicity events
306  // along the history path and set its production scale to the input
307  // scale
308  // IN int iPart : Parton in refEvent to be checked / rescaled
309  // Event& refEvent : Reference event for iPart
310  // double scale : Scale to be set as production scale for
311  // unchanged copies of iPart in subsequent steps
312  void scaleCopies(int iPart, const Event& refEvent, double rho);
313 
314  // Function to set the OVERALL EVENT SCALES [=state.scale()] to
315  // the scale of the last clustering
316  // NO INPUT
317  // NO OUTPUT
318  void setEventScales();
319 
320  // Function to print information on the reconstructed scales in one path.
321  // For debug only.
322  void printScales() { if ( mother ) mother->printScales();
323  cout << " size " << state.size() << " scale " << scale << " clusterIn "
324  << clusterIn.pT() << " state.scale() " << state.scale() << endl; }
325 
326  // Function to project paths onto desired paths.
327  bool trimHistories();
328  // Function to tag history for removal.
329  void remove(){ doInclude = false; }
330  // Function to return flag of allowed histories to choose from.
331  bool keep(){ return doInclude; }
332  // Function implementing checks on a paths, for deciding if the path should
333  // be considered valid or not.
334  bool keepHistory();
335  // Function to check if a path is ordered in evolution pT.
336  bool isOrderedPath( double maxscale );
337 
338  bool followsInputPath();
339 
340  // Function to check if all reconstucted states in a path pass the merging
341  // scale cut.
342  bool allIntermediateAboveRhoMS( double rhoms, bool good = true );
343  // Function to check if any ordered paths were found (and kept).
344  bool foundAnyOrderedPaths();
345 
346  // Functions to return the z value of the last ISR splitting
347  // NO INPUT
348  // OUTPUT double : z value of last ISR splitting in history
349  double zISR();
350 
351  // Functions to return the z value of the last FSR splitting
352  // NO INPUT
353  // OUTPUT double : z value of last FSR splitting in history
354  double zFSR();
355 
356  // Functions to return the pT scale of the last ISR splitting
357  // NO INPUT
358  // OUTPUT double : pT scale of last ISR splitting in history
359  double pTISR();
360 
361  // Functions to return the pT scale of the last FSR splitting
362  // NO INPUT
363  // OUTPUT double : pT scale of last FSR splitting in history
364  double pTFSR();
365 
366  // Functions to return the event with nSteps additional partons
367  // INPUT int : Number of splittings in the event,
368  // as counted from core 2->2 process
369  // OUTPUT Event : event with nSteps additional partons
370  Event clusteredState( int nSteps);
371 
372  // Function to choose a path from all paths in the tree
373  // according to their splitting probabilities
374  // IN double : Random number
375  // OUT History* : Leaf of history path chosen
376  History * select(double rnd);
377 
378  // For a full path, find the weight calculated from the ratio of
379  // couplings, the no-emission probabilities, and possible PDF
380  // ratios. This function should only be called for the last history
381  // node of a full path.
382  // IN TimeShower : Already initialised shower object to be used as
383  // trial shower
384  // double : alpha_s value used in ME calculation
385  // double : Maximal mass scale of the problem (e.g. E_CM)
386  // AlphaStrong: Initialised shower alpha_s object for FSR alpha_s
387  // ratio calculation
388  // AlphaStrong: Initialised shower alpha_s object for ISR alpha_s
389  // ratio calculation (can be different from previous)
390  vector<double> weightTree(PartonLevel* trial, double as0, double aem0,
391  double maxscale, double pdfScale, AlphaStrong * asFSR, AlphaStrong * asISR,
392  AlphaEM * aemFSR, AlphaEM * aemISR, vector<double>& asWeight,
393  vector<double>& aemWeight, vector<double>& pdfWeight);
394 
395  // Function to return the \alpha_s-ratio part of the CKKWL weight.
396  vector<double> weightTreeAlphaS( double as0, AlphaStrong * asFSR,
397  AlphaStrong * asISR, int njetMax = -1, bool asVarInME = false );
398  // Function to return the \alpha_em-ratio part of the CKKWL weight.
399  vector<double> weightTreeAlphaEM( double aem0, AlphaEM * aemFSR,
400  AlphaEM * aemISR, int njetMax = -1 );
401  // Function to return the PDF-ratio part of the CKKWL weight.
402  vector<double> weightTreePDFs( double maxscale, double pdfScale,
403  int njetMax = -1 );
404  // Function to return the no-emission probability part of the CKKWL weight.
405  vector<double> weightTreeEmissions( PartonLevel* trial, int type,
406  int njetMin, int njetMax, double maxscale );
407 
408  // Function to generate the O(\alpha_s)-term of the CKKWL-weight
409  double weightFirst(PartonLevel* trial, double as0, double muR,
410  double maxscale, AlphaStrong * asFSR, AlphaStrong * asISR, Rndm* rndmPtr );
411 
412  // Function to generate the O(\alpha_s)-term of the \alpha_s-ratios
413  // appearing in the CKKWL-weight.
414  double weightFirstAlphaS( double as0, double muR, AlphaStrong * asFSR,
415  AlphaStrong * asISR);
416  // Function to generate the O(\alpha_em)-term of the \alpha_em-ratios
417  // appearing in the CKKWL-weight.
418  double weightFirstAlphaEM( double aem0, double muR, AlphaEM * aemFSR,
419  AlphaEM * aemISR);
420  // Function to generate the O(\alpha_s)-term of the PDF-ratios
421  // appearing in the CKKWL-weight.
422  double weightFirstPDFs( double as0, double maxscale, double pdfScale,
423  Rndm* rndmPtr );
424  // Function to generate the O(\alpha_s)-term of the no-emission
425  // probabilities appearing in the CKKWL-weight.
426  double weightFirstEmissions(PartonLevel* trial, double as0, double maxscale,
427  AlphaStrong * asFSR, AlphaStrong * asISR, bool fixpdf, bool fixas );
428 
429  // Function to return the default factorisation scale of the hard process.
430  double hardFacScale(const Event& event);
431  // Function to return the default renormalisation scale of the hard process.
432  double hardRenScale(const Event& event);
433 
434  // Perform a trial shower using the \a pythia object between
435  // maxscale down to this scale and return the corresponding Sudakov
436  // form factor.
437  // IN trialShower : Shower object used as trial shower
438  // double : Maximum scale for trial shower branching
439  // OUT 0.0 : trial shower emission outside allowed pT range
440  // 1.0 : trial shower successful (any emission was below
441  // the minimal scale )
442  vector<double> doTrialShower(PartonLevel* trial, int type, double maxscale,
443  double minscale = 0.);
444 
445  // Function to bookkeep the indices of weights generated in countEmissions
446  bool updateind(vector<int> & ind, int i, int N);
447 
448  // Function to count number of emissions between two scales for NLO merging
449  vector<double> countEmissions(PartonLevel* trial, double maxscale,
450  double minscale, int showerType, double as0, AlphaStrong * asFSR,
451  AlphaStrong * asISR, int N, bool fixpdf, bool fixas);
452 
453  // Function to integrate PDF ratios between two scales over x and t,
454  // where the PDFs are always evaluated at the lower t-integration limit
455  double monteCarloPDFratios(int flav, double x, double maxScale,
456  double minScale, double pdfScale, double asME, Rndm* rndmPtr);
457 
458  // Default: Check if a ordered (and complete) path has been found in
459  // the initial node, in which case we will no longer be interested in
460  // any unordered paths.
461  bool onlyOrderedPaths();
462 
463  // Check if a strongly ordered (and complete) path has been found in the
464  // initial node, in which case we will no longer be interested in
465  // any unordered paths.
466  bool onlyStronglyOrderedPaths();
467 
468  // Check if an allowed (according to user-criterion) path has been found in
469  // the initial node, in which case we will no longer be interested in
470  // any forbidden paths.
471  bool onlyAllowedPaths();
472 
473  // When a full path has been found, register it with the initial
474  // history node.
475  // IN History : History to be registered as path
476  // bool : Specifying if clusterings so far were ordered
477  // bool : Specifying if path is complete down to 2->2 process
478  // OUT true if History object forms a plausible path (eg prob>0 ...)
479  bool registerPath(History & l, bool isOrdered, bool isStronglyOrdered,
480  bool isAllowed, bool isComplete);
481 
482  // For the history-defining state (and if necessary interfering
483  // states), find all possible clusterings.
484  // NO INPUT
485  // OUT vector of all (rad,rec,emt) systems
486  vector<Clustering> getAllQCDClusterings();
487 
488  // For one given state, find all possible clusterings.
489  // IN Event : state to be investigated
490  // OUT vector of all (rad,rec,emt) systems in the state
491  vector<Clustering> getQCDClusterings( const Event& event);
492 
493  // Function to construct (rad,rec,emt) triples from the event
494  // IN int : Position of Emitted in event record for which
495  // dipoles should be constructed
496  // int : Colour topogy to be tested
497  // 1= g -> qqbar, causing 2 -> 2 dipole splitting
498  // 2= q(bar) -> q(bar) g && g -> gg,
499  // causing a 2 -> 3 dipole splitting
500  // Event : event record to be checked for ptential partners
501  // OUT vector of all allowed radiator+recoiler+emitted triples
502  vector<Clustering> findQCDTriple (int emtTagIn, int colTopIn,
503  const Event& event, vector<int> posFinalPartn,
504  vector <int> posInitPartn );
505 
506  vector<Clustering> getAllEWClusterings();
507  vector<Clustering> getEWClusterings( const Event& event);
508  vector<Clustering> findEWTripleW( int emtTagIn, const Event& event,
509  vector<int> posFinalPartn, vector<int> posInitPartn );
510  vector<Clustering> findEWTripleZ( int emtTagIn, const Event& event,
511  vector<int> posFinalPartn, vector<int> posInitPartn );
512 
513  vector<Clustering> getAllSQCDClusterings();
514  vector<Clustering> getSQCDClusterings( const Event& event);
515  vector<Clustering> findSQCDTriple (int emtTagIn, int colTopIn,
516  const Event& event, vector<int> posFinalPartn,
517  vector <int> posInitPartn );
518 
519  // Function to attach (spin-dependent duplicates of) a clustering.
520  void attachClusterings (vector<Clustering>& clus, int iEmt, int iRad,
521  int iRec, int iPartner, double pT, const Event& event);
522 
523  // Calculate and return the probability of a clustering.
524  // IN Clustering : rad,rec,emt - System for which the splitting
525  // probability should be calcuated
526  // OUT splitting probability
527  double getProb(const Clustering & SystemIn);
528 
529  // Set up the beams (fill the beam particles with the correct
530  // current incoming particles) to allow calculation of splitting
531  // probability.
532  // For interleaved evolution, set assignments dividing PDFs into
533  // sea and valence content. This assignment is, until a history path
534  // is chosen and a first trial shower performed, not fully correct
535  // (since content is chosen form too high x and too low scale). The
536  // assignment used for reweighting will be corrected after trial
537  // showering
538  void setupBeams();
539 
540  // Calculate the PDF ratio used in the argument of the no-emission
541  // probability.
542  double pdfForSudakov();
543 
544  // Calculate the hard process matrix element to include in the selection
545  // probability.
546  double hardProcessME( const Event& event);
547 
548  // Perform the clustering of the current state and return the
549  // clustered state.
550  // IN Clustering : rad,rec,emt triple to be clustered to two partons
551  // OUT clustered state
552  Event cluster( Clustering & inSystem);
553 
554  // Function to get the flavour of the radiator before the splitting
555  // for clustering
556  // IN int : Position of the radiator after the splitting, in the event
557  // int : Position of the emitted after the splitting, in the event
558  // Event : Reference event
559  // OUT int : Flavour of the radiator before the splitting
560  int getRadBeforeFlav(const int radAfter, const int emtAfter,
561  const Event& event);
562 
563  // Function to get the spin of the radiator before the splitting
564  // IN int : Spin of the radiator after the splitting
565  // int : Spin of the emitted after the splitting
566  // OUT int : Spin of the radiator before the splitting
567  int getRadBeforeSpin(const int radAfter, const int emtAfter,
568  const int spinRadAfter, const int spinEmtAfter,
569  const Event& event);
570 
571  // Function to get the colour of the radiator before the splitting
572  // for clustering
573  // IN int : Position of the radiator after the splitting, in the event
574  // int : Position of the emitted after the splitting, in the event
575  // Event : Reference event
576  // OUT int : Colour of the radiator before the splitting
577  int getRadBeforeCol(const int radAfter, const int emtAfter,
578  const Event& event);
579 
580  // Function to get the anticolour of the radiator before the splitting
581  // for clustering
582  // IN int : Position of the radiator after the splitting, in the event
583  // int : Position of the emitted after the splitting, in the event
584  // Event : Reference event
585  // OUT int : Anticolour of the radiator before the splitting
586  int getRadBeforeAcol(const int radAfter, const int emtAfter,
587  const Event& event);
588 
589  // Function to get the parton connected to in by a colour line
590  // IN int : Position of parton for which partner should be found
591  // Event : Reference event
592  // OUT int : If a colour line connects the "in" parton with another
593  // parton, return the Position of the partner, else return 0
594  int getColPartner(const int in, const Event& event);
595 
596  // Function to get the parton connected to in by an anticolour line
597  // IN int : Position of parton for which partner should be found
598  // Event : Reference event
599  // OUT int : If an anticolour line connects the "in" parton with another
600  // parton, return the Position of the partner, else return 0
601  int getAcolPartner(const int in, const Event& event);
602 
603  // Function to get the list of partons connected to the particle
604  // formed by reclusterinf emt and rad by colour and anticolour lines
605  // IN int : Position of radiator in the clustering
606  // IN int : Position of emitted in the clustering
607  // Event : Reference event
608  // OUT vector<int> : List of positions of all partons that are connected
609  // to the parton that will be formed
610  // by clustering emt and rad.
611  vector<int> getReclusteredPartners(const int rad, const int emt,
612  const Event& event);
613 
614  // Function to extract a chain of colour-connected partons in
615  // the event
616  // IN int : Type of parton from which to start extracting a
617  // parton chain. If the starting point is a quark
618  // i.e. flavType = 1, a chain of partons that are
619  // consecutively connected by colour-lines will be
620  // extracted. If the starting point is an antiquark
621  // i.e. flavType =-1, a chain of partons that are
622  // consecutively connected by anticolour-lines
623  // will be extracted.
624  // IN int : Position of the parton from which a
625  // colour-connected chain should be derived
626  // IN Event : Refernence event
627  // IN/OUT vector<int> : Partons that should be excluded from the search.
628  // OUT vector<int> : Positions of partons along the chain
629  // OUT bool : Found singlet / did not find singlet
630  bool getColSinglet(const int flavType, const int iParton,
631  const Event& event, vector<int>& exclude,
632  vector<int>& colSinglet);
633 
634  // Function to check that a set of partons forms a colour singlet
635  // IN Event : Reference event
636  // IN vector<int> : Positions of the partons in the set
637  // OUT bool : Is a colour singlet / is not
638  bool isColSinglet( const Event& event, vector<int> system);
639  // Function to check that a set of partons forms a flavour singlet
640  // IN Event : Reference event
641  // IN vector<int> : Positions of the partons in the set
642  // IN int : Flavour of all the quarks in the set, if
643  // all quarks in a set should have a fixed flavour
644  // OUT bool : Is a flavour singlet / is not
645  bool isFlavSinglet( const Event& event,
646  vector<int> system, int flav=0);
647 
648  // Function to properly colour-connect the radiator to the rest of
649  // the event, as needed during clustering
650  // IN Particle& : Particle to be connected
651  // Particle : Recoiler forming a dipole with Radiator
652  // Event : event to which Radiator shall be appended
653  // OUT true : Radiator could be connected to the event
654  // false : Radiator could not be connected to the
655  // event or the resulting event was
656  // non-valid
657  bool connectRadiator( Particle& radiator, const int radType,
658  const Particle& recoiler, const int recType,
659  const Event& event );
660 
661  // Function to find a colour (anticolour) index in the input event
662  // IN int col : Colour tag to be investigated
663  // int iExclude1 : Identifier of first particle to be excluded
664  // from search
665  // int iExclude2 : Identifier of second particle to be excluded
666  // from search
667  // Event event : event to be searched for colour tag
668  // int type : Tag to define if col should be counted as
669  // colour (type = 1) [->find anti-colour index
670  // contracted with col]
671  // anticolour (type = 2) [->find colour index
672  // contracted with col]
673  // OUT int : Position of particle in event record
674  // contraced with col [0 if col is free tag]
675  int FindCol(int col, int iExclude1, int iExclude2,
676  const Event& event, int type, bool isHardIn);
677 
678  // Function to in the input event find a particle with quantum
679  // numbers matching those of the input particle
680  // IN Particle : Particle to be searched for
681  // Event : Event to be searched in
682  // OUT int : > 0 : Position of matching particle in event
683  // < 0 : No match in event
684  int FindParticle( const Particle& particle, const Event& event,
685  bool checkStatus = true );
686 
687  // Function to check if rad,emt,rec triple is allowed for clustering
688  // IN int rad,emt,rec,partner : Positions (in event record) of the three
689  // particles considered for clustering, and the
690  // correct colour-connected recoiler (=partner)
691  // Event event : Reference event
692  bool allowedClustering( int rad, int emt, int rec, int partner,
693  const Event& event );
694 
695  // Function to check if rad,emt,rec triple is results in
696  // colour singlet radBefore+recBefore
697  // IN int rad,emt,rec : Positions (in event record) of the three
698  // particles considered for clustering
699  // Event event : Reference event
700  bool isSinglett( int rad, int emt, int rec, const Event& event );
701 
702  // Function to check if event is sensibly constructed: Meaning
703  // that all colour indices are contracted and that the charge in
704  // initial and final states matches
705  // IN event : event to be checked
706  // OUT TRUE : event is properly construced
707  // FALSE : event not valid
708  bool validEvent( const Event& event );
709 
710  // Function to check whether two clusterings are identical, used
711  // for finding the history path in the mother -> children direction
712  bool equalClustering( Clustering clus1 , Clustering clus2 );
713 
714  // Chose dummy scale for event construction. By default, choose
715  // sHat for 2->Boson(->2)+ n partons processes and
716  // M_Boson for 2->Boson(->) processes
717  double choseHardScale( const Event& event ) const;
718 
719  // If the state has an incoming hadron return the flavour of the
720  // parton entering the hard interaction. Otherwise return 0
721  int getCurrentFlav(const int) const;
722 
723  // If the state has an incoming hadron return the x-value for the
724  // parton entering the hard interaction. Otherwise return 0.
725  double getCurrentX(const int) const;
726 
727  double getCurrentZ(const int rad, const int rec, const int emt,
728  int idRadBef = 0) const;
729 
730  // Function to compute "pythia pT separation" from Particle input
731  double pTLund(const Event& event, int radAfterBranch, int emtAfterBranch,
732  int recAfterBranch, int showerType, int idRadBef = 0);
733 
734  // Function to return the position of the initial line before (or after)
735  // a single (!) splitting.
736  int posChangedIncoming(const Event& event, bool before);
737 
738  // Function to give back the ratio of PDFs and PDF * splitting kernels
739  // needed to convert a splitting at scale pdfScale, chosen with running
740  // PDFs, to a splitting chosen with PDFs at a fixed scale mu. As needed to
741  // properly count emissions.
742  double pdfFactor( const Event& event, const int type, double pdfScale,
743  double mu );
744 
745  // Function giving the product of splitting kernels and PDFs so that the
746  // resulting flavour is given by flav. This is used as a helper routine
747  // to dgauss
748  double integrand(int flav, double x, double scaleInt, double z);
749 
750  // Function providing a list of possible new flavours after a w emssion
751  // from the input flavour.
752  vector<int> posFlavCKM(int flav);
753 
754  // Check if the new flavour structure is possible.
755  // If clusType is 1 final clustering is assumed, otherwise initial clustering
756  // is assumed.
757  bool checkFlavour(vector<int>& flavCounts, int flavRad, int flavRadBef,
758  int clusType);
759 
760  // Check if the weak recoil structure is allowed.
761  bool checkWeakRecoils(map<int,int> &allowedRecoils, bool isFirst = false);
762 
763  // Find the recoiler for an ISR scattered weak particle.
764 
765  int findISRRecoiler();
766 
767  // Reverse the boost carried out by the ISR.
768  // The three first momenta are given by the ME,
769  // the last two are filled in by this function.
770  void reverseBoostISR(Vec4& pMother, Vec4& pSister, Vec4& pPartner,
771  Vec4& pDaughter, Vec4& pRecoiler, int sign, double eCM, double& phi);
772 
773  // Check if an event reclustered into a 2 -> 2 dijet.
774  // (Only enabled if W reclustering is used).
775  bool isQCD2to2(const Event& event);
776 
777  // Check if an event reclustered into a 2 -> 1 Drell-Yan.
778  // (Only enabled if W reclustering is used).
779  bool isEW2to1(const Event& event);
780 
781  // Set selected child indices.
782  void setSelectedChild();
783 
784  // Setup the weak dipole showers.
785  void setupSimpleWeakShower(int nSteps);
786 
787  // Update weak dipoles after an emission.
788  void transferSimpleWeakShower(vector<int> &mode, vector<Vec4> &mom,
789  vector<int> fermionLines, vector<pair<int,int> > &dipoles, int nSteps);
790 
791  // Update the weak modes after an emission.
792  vector<int> updateWeakModes(vector<int>& weakModes,
793  map<int,int>& stateTransfer);
794 
795  // Update the fermion lines for the 2 -> 2 process. This is needed for
796  // the weak probabilities.
797  vector<int> updateWeakFermionLines(vector<int> fermionLines,
798  map<int,int>& stateTransfer);
799 
800  // Update the list of weak dipoles. This is needed to setup the PS correctly.
801  vector<pair<int,int> > updateWeakDipoles(vector<pair<int,int> > &dipoles,
802  map<int,int>& stateTransfer);
803 
804  // Setup the hard process information needed for calculating weak
805  // probabilities and setting up the shower.
806  void setupWeakHard(vector<int>& mode, vector<int>& fermionLines,
807  vector<Vec4>& mom);
808 
809  // Functions to retrieve scale information from external showers.
810  double getShowerPluginScale(const Event& event, int rad, int emt, int rec,
811  string key, double scalePythia);
812 
813  //----------------------------------------------------------------------//
814  // Class members.
815  //----------------------------------------------------------------------//
816 
817  // The state of the event correponding to this step in the
818  // reconstruction.
819  Event state;
820 
821  // The previous step from which this step has been clustered. If
822  // null, this is the initial step with the n-jet state generated by
823  // the matrix element.
824  History * mother;
825 
826  // The different steps corresponding to possible clusterings of this
827  // state.
828  vector<History *> children;
829 
830  // After a path is selected, store the child index.
831  int selectedChild;
832 
833  // The different paths which have been reconstructed indexed with
834  // the (incremental) corresponding probability. This map is empty
835  // unless this is the initial step (mother == 0).
836  map<double,History *> paths;
837 
838  // The sum of the probabilities of the full paths. This is zero
839  // unless this is the initial step (mother == 0).
840  double sumpath;
841 
842  // The different allowed paths after projection, indexed with
843  // the (incremental) corresponding probability. This map is empty
844  // unless this is the initial step (mother == 0).
845  map<double,History *> goodBranches, badBranches;
846  // The sum of the probabilities of allowed paths after projection. This is
847  // zero unless this is the initial step (mother == 0).
848  double sumGoodBranches, sumBadBranches;
849 
850  // This is set true if an ordered (and complete) path has been found
851  // and inserted in paths.
852  bool foundOrderedPath;
853 
854  // This is set true if a strongly ordered (and complete) path has been found
855  // and inserted in paths.
856  bool foundStronglyOrderedPath;
857 
858  // This is set true if an allowed (according to a user criterion) path has
859  // been found and inserted in paths.
860  bool foundAllowedPath;
861 
862  // This is set true if a complete (with the required number of
863  // clusterings) path has been found and inserted in paths.
864  bool foundCompletePath;
865 
866  // The scale of this step, corresponding to clustering which
867  // constructed the corresponding state (or the merging scale in case
868  // mother == 0).
869  double scale;
870 
871  // Flag indicating if a clustering in the construction of all histories is
872  // the next clustering demanded by inout clusterings in LesHouches 2.0
873  // accord.
874  bool nextInInput;
875 
876  // The probability associated with this step and the previous steps.
877  double prob;
878 
879  // The partons and scale of the last clustering.
880  Clustering clusterIn;
881  int iReclusteredOld, iReclusteredNew;
882 
883  // Flag to include the path amongst allowed paths.
884  bool doInclude;
885 
886  // Pointer to MergingHooks object to get all the settings.
887  MergingHooksPtr mergingHooksPtr;
888 
889  // The default constructor is private.
890  History() : mother(), selectedChild(), sumpath(), sumGoodBranches(),
891  sumBadBranches(), foundOrderedPath(), foundStronglyOrderedPath(),
892  foundAllowedPath(), foundCompletePath(), scale(), nextInInput(), prob(),
893  iReclusteredOld(), iReclusteredNew(), doInclude(), mergingHooksPtr(),
894  particleDataPtr(), infoPtr(), showers(), coupSMPtr(), sumScalarPT(),
895  probMaxSave(), depth(), minDepthSave() {}
896 
897  // The copy-constructor is private.
898  History(const History &) {}
899 
900  // The assignment operator is private.
901  History & operator=(const History &) {
902  return *this;
903  }
904 
905  // BeamParticle to get access to PDFs
906  BeamParticle beamA;
907  BeamParticle beamB;
908  // ParticleData needed to initialise the shower AND to get the
909  // correct masses of partons needed in calculation of probability
910  ParticleData* particleDataPtr;
911 
912  // Info object to have access to all information read from LHE file
913  Info* infoPtr;
914 
915  // Class for calculation weak shower ME.
916  SimpleWeakShowerMEs simpleWeakShowerMEs;
917 
918  // Pointer to showers, to simplify external clusterings.
919  PartonLevel* showers;
920 
921  // Pointer to standard model couplings.
922  CoupSM* coupSMPtr;
923 
924  // Minimal scalar sum of pT used in Herwig to choose history
925  double sumScalarPT;
926 
927  double probMaxSave;
928  double probMax() {
929  if (mother) return mother->probMax();
930  return probMaxSave;
931  }
932  void updateProbMax(double probIn, bool isComplete = false) {
933  if (mother) mother->updateProbMax(probIn, isComplete);
934  if (!isComplete && !foundCompletePath) return;
935  if (abs(probIn) > probMaxSave) probMaxSave = probIn;
936  }
937 
938  int depth, minDepthSave;
939  int minDepth() {
940  if ( mother ) return mother->minDepth();
941  return minDepthSave;
942  }
943  void updateMinDepth(int depthIn) {
944  if ( mother ) return mother->updateMinDepth(depthIn);
945  minDepthSave = (minDepthSave>0) ? min(minDepthSave,depthIn) : depthIn;
946  }
947 
948 };
949 
950 //==========================================================================
951 
952 } // end namespace Pythia8
953 
954 #endif // end Pythia8_History_H
Definition: AgUStep.h:26