herwig is hosted by Hepforge, IPPP Durham
Herwig 7.3.0
DipoleShowerHandler.h
1// -*- C++ -*-
2//
3// DipoleShowerHandler.h is a part of Herwig - A multi-purpose Monte Carlo event generator
4// Copyright (C) 2002-2019 The Herwig Collaboration
5//
6// Herwig is licenced under version 3 of the GPL, see COPYING for details.
7// Please respect the MCnet academic guidelines, see GUIDELINES for details.
8//
9#ifndef HERWIG_DipoleShowerHandler_H
10#define HERWIG_DipoleShowerHandler_H
11//
12// This is the declaration of the DipoleShowerHandler class.
13//
14
15#include "Herwig/Shower/ShowerHandler.h"
16
17#include "Herwig/Shower/Dipole/DipoleShowerHandler.fh"
18#include "Herwig/Shower/Dipole/Base/DipoleSplittingInfo.h"
19#include "Herwig/Shower/Dipole/Base/DipoleSplittingReweight.h"
20#include "Herwig/Shower/Dipole/Kernels/DipoleSplittingKernel.h"
21#include "Herwig/Shower/Dipole/Base/DipoleSplittingGenerator.h"
22#include "Herwig/Shower/Dipole/Base/DipoleEventRecord.h"
23#include "Herwig/Shower/Dipole/Base/DipoleEvolutionOrdering.h"
24#include "Herwig/Shower/Dipole/Base/DipoleEventReweight.h"
25#include "Herwig/Shower/Dipole/Utility/ConstituentReshuffler.h"
26#include "Herwig/Shower/Dipole/Utility/IntrinsicPtGenerator.h"
27#include "Herwig/MatrixElement/Matchbox/Base/MergerBase.h"
28#include "Herwig/MatrixElement/Matchbox/Matching/ShowerApproximation.h"
29
30#include "Herwig/Shower/Dipole/SpinCorrelations/DipoleVertexRecord.h"
31#include "Herwig/MatrixElement/Matchbox/Utility/DensityOperator.h"
32
33#include <tuple>
34
35namespace Herwig {
36
37using namespace ThePEG;
38
50
51
52 friend class Merger;
53
54public:
55
60
61public:
62
63
64 inline void colourPrint();
65
69 struct RedoShower {};
70
74 void addSplitting(Ptr<DipoleSplittingKernel>::ptr sp) {
75 kernels.push_back(sp);
76 }
77
81 void resetAlphaS(Ptr<AlphaSBase>::tptr);
82
83
84 virtual void cascade(tPVector);
85
89 void resetReweight(Ptr<DipoleSplittingReweight>::tptr);
90
95 virtual bool canHandleMatchboxTrunc() const { return false; }
96
101 virtual bool isReshuffling() const { return false; }
102
106 virtual Energy hardScale() const {
107 return muPt;
108 }
109
114 double as(Energy Q)const{return theGlobalAlphaS->value(sqr(Q));}
115
121 double Nf(Energy Q)const{return theGlobalAlphaS->Nf(sqr(Q));}
122
123 /*
124 * Access the vertex record
125 */
126 DipoleVertexRecord& vertexRecord() { return theVertexRecord; }
127
132
133
138 void setMerger(Ptr<MergerBase>::ptr mh){theMergingHelper=mh;}
139
140
141
142public:
143
148 const map<PPtr,size_t>& particleIndices() const { return eventRecord().particleIndices(); }
149
154 const cPDVector& particlesAfter() const { return eventRecord().particlesAfter(); }
155
160
165 const map<pair<vector<PDT::Colour>,pair<size_t,size_t> >,double>& correlatorMap() const {
167 }
168
172 Ptr<ColourBasis>::tptr colourBasis() {
174 }
175
180
181protected:
182
186 void addCandidates(PPair particles, list<DipoleSplittingInfo>& clist) const;
187
191 void getCandidates(list<DipoleSplittingInfo>& clist) const;
192
197
201 Energy nextSubleadingSplitting(Energy hardPt,
202 Energy optHardPt, Energy optCutoff,
203 const bool decay);
204
205protected:
206
207 typedef multimap<DipoleIndex,Ptr<DipoleSplittingGenerator>::ptr> GeneratorMap;
208
212 virtual tPPair cascade(tSubProPtr sub, XCombPtr xcomb) {
213 return cascade(sub,xcomb,ZERO,ZERO);
214 }
215
220 Energy optHardPt, Energy optCutoff);
221
227 Ptr<DipoleSplittingReweight>::tptr rw =
228 Ptr<DipoleSplittingReweight>::tptr());
229
233 void hardScales(Energy2 scale);
234
238 void hardScalesSubleading(list<DipoleSplittingInfo> candidates,Energy hardPt);
239
243 Ptr<DipoleEvolutionOrdering>::tptr evolutionOrdering() const { return theEvolutionOrdering; }
244
249
253 void decayConstituentReshuffle( PerturbativeProcessPtr decayProc);
254
258 GeneratorMap& generators() { return theGenerators; }
259
264
269
273 const vector<Ptr<DipoleSplittingKernel>::ptr>& splittingKernels() const {
274 return kernels;
275 }
276
280 const set<long>& offShellPartons() { return theColouredOffShellInShower; }
281
286 bool isMERegion(const Energy winnerScale,
287 const DipoleSplittingInfo & winner,
288 const list<Dipole>::iterator winnerDip);
289
294 bool realign();
295
299 virtual int showerPhaseSpaceOption() const {
300 return theZBoundaries;
301 }
302
303protected:
304
308 void doCascade(unsigned int& emDone,
309 Energy optHardPt = ZERO,
310 Energy optCutoff = ZERO,
311 const bool decay = false);
315 void setNEmissions(unsigned int n){nEmissions=n;}
316
317
323 const Dipole& dip,
324 pair<bool,bool> conf,
325 Energy optHardPt = ZERO,
326 Energy optCutoff = ZERO);
327
333 Energy optHardPt = ZERO,
334 Energy optCutoff = ZERO);
335
341 Energy optHardPt = ZERO,
342 Energy optCutoff = ZERO);
343
349 const DipoleIndex& index,
350 double emitterX, double spectatorX,
351 pair<bool,bool> conf,
352 tPPtr emitter, tPPtr spectator,
353 Energy startScale,
354 Energy optHardPt = ZERO,
355 Energy optCutoff = ZERO);
356
357public:
358
366
372 void persistentInput(PersistentIStream & is, int version);
374
381 static void Init();
382
383protected:
384
391 virtual IBPtr clone() const;
392
397 virtual IBPtr fullclone() const;
399
400
401// If needed, insert declarations of virtual function defined in the
402// InterfacedBase class here (using ThePEG-interfaced-decl in Emacs).
403
404
405protected:
406
414 virtual void doinit();
415
420 virtual void doinitrun();
421
426 virtual void dofinish();
428
429
430private:
431
435 vector<Ptr<DipoleSplittingKernel>::ptr> kernels;
436
440 Ptr<DipoleEvolutionOrdering>::ptr theEvolutionOrdering;
441
445 Ptr<ConstituentReshuffler>::ptr constituentReshuffler;
446
450 Ptr<IntrinsicPtGenerator>::ptr intrinsicPtGenerator;
451
455 Ptr<AlphaSBase>::ptr theGlobalAlphaS;
456
462
467 unsigned int nEmissions;
468
473
478
483
488 //bool theAnalyseSpinCorrelations;
489
494
499
504
510
520
526
531
536
541
546
551
552private:
553
564
569
570private:
571
576 GeneratorMap theGenerators;
577
582
587
591 unsigned int nTries;
592
597
602
608 vector<vector<std::tuple<Energy,double,bool> > > theWeightsVector;
609
613 size_t winnerIndex;//debug
614 size_t kernelIndex;//debug
615 size_t winningKernelIndex;//debug
616 vector<Energy> scales;//debug
617
618private:
619
624
629
633 Ptr<ShowerApproximation>::tptr theShowerApproximation;
634
639
644 unsigned long theFreezeGrid;
645
650
654 Ptr<DipoleEventReweight>::ptr theEventReweight;
655
659 Ptr<DipoleSplittingReweight>::ptr theSplittingReweight;
660
664 static bool firstWarn;
665
669 Energy maxPt;
670
674 Energy muPt;
675
676
683 Ptr<MergerBase>::ptr theMergingHelper;
684
685
691
697
702
703private:
704
710
716
717};
718
719}
720
721#include "ThePEG/Utilities/ClassTraits.h"
722
723namespace ThePEG {
724
729template <>
730struct BaseClassTrait<Herwig::DipoleShowerHandler,1> {
733};
734
737template <>
738struct ClassTraits<Herwig::DipoleShowerHandler>
739 : public ClassTraitsBase<Herwig::DipoleShowerHandler> {
741 static string className() { return "Herwig::DipoleShowerHandler"; }
749 static string library() { return "HwDipoleShower.so"; }
750};
751
754}
755
756#endif /* HERWIG_DipoleShowerHandler_H */
Here is the documentation of the DensityOperator class.
const map< pair< vector< PDT::Colour >, pair< size_t, size_t > >, double > & correlatorMap() const
Get the correlator map.
Ptr< ColourBasis >::tptr colourBasis()
Get the colour basis.
The DipoleEventRecord class is used internally by the dipole shower.
cPDVector & particlesAfter()
Return the particles after emission.
map< PPtr, size_t > & particleIndices()
Return the dictionary for the particles.
bool getContinueSubleadingNc() const
Get the continue subleading Nc flag.
DensityOperator & densityOperator()
Return the density operator.
DipoleIndex is used to index splitting generators for a particular dipole.
The DipoleShowerHandler class manages the showering using the dipole shower algorithm.
bool chainOrderVetoScales
Apply chain ordering to events from matrix element corrections.
void setNEmissions(unsigned int n)
Set the number of emissions.
double as(Energy Q) const
Calculate the alpha_s value the shower uses for Q.
void getCandidates(list< DipoleSplittingInfo > &clist) const
Get all possible radiating dipoles.
GeneratorMap theGenerators
The splitting generators indexed by the dipole indices they can work on.
Energy getWinner(SubleadingSplittingInfo &winner, Energy optHardPt=ZERO, Energy optCutoff=ZERO)
Get the winning splitting for the given dipole and configuration.
Energy maxPt
The shower starting scale for the last event encountered.
Ptr< DipoleEventReweight >::ptr theEventReweight
A pointer to the dipole event reweight object.
virtual IBPtr fullclone() const
Make a clone of this object, possibly modifying the cloned object to make it sane.
virtual void dofinish()
Finalize this object.
const map< PPtr, size_t > & particleIndices() const
Return the dictionary for the incoming particles and outgoing partons in the event,...
int realignmentScheme
Switch to record information required for the nearest neighbour analysis.
virtual int showerPhaseSpaceOption() const
The choice of z boundaries; 0 = restricted, 1 = open, 2 = mixed/other.
vector< vector< std::tuple< Energy, double, bool > > > theWeightsVector
Vector of candidate splittings containing a vector of the weights, scale and a bool for every step in...
void performSplitting(DipoleSplittingInfo &) const
Perform a splitting independent of any chains.
const DipoleVertexRecord & vertexRecord() const
Return the event record.
virtual tPPair cascade(tSubProPtr sub, XCombPtr xcomb)
The main method which manages the showering of a subprocess.
void persistentOutput(PersistentOStream &os) const
Function used to write out object persistently.
virtual Energy hardScale() const
Return the relevant hard scale to be used in the profile scales.
const vector< Ptr< DipoleSplittingKernel >::ptr > & splittingKernels() const
Return the splitting kernels.
Ptr< DipoleSplittingReweight >::ptr theSplittingReweight
A pointer to a global dipole splitting reweight.
bool isMERegion(const Energy winnerScale, const DipoleSplittingInfo &winner, const list< Dipole >::iterator winnerDip)
In a merging setup this function checks if the next shower configuration is part of the matrix elemen...
Ptr< ShowerApproximation >::tptr theShowerApproximation
The matching subtraction, if appropriate.
double negCMECScaling
Scaling factor for the negative colour matrix element corrections.
Energy2 densityOperatorCutoff
Cutoff scale for the invariants (e.g.
virtual bool canHandleMatchboxTrunc() const
Return true, if the shower handler can generate a truncated shower for POWHEG style events generated ...
bool doSubleadingNc
Switch on or off subleading Nc corrections.
bool continueSubleadingNc() const
Return the continue subleading Nc flag from the event record.
void resetAlphaS(Ptr< AlphaSBase >::tptr)
Reset the alpha_s for all splitting kernels.
void doCascade(unsigned int &emDone, Energy optHardPt=ZERO, Energy optCutoff=ZERO, const bool decay=false)
Perform the cascade.
virtual void doinit()
Initialize this object after the setup phase before saving an EventGenerator to disk.
unsigned int nTries
The number of shoer tries so far.
void persistentInput(PersistentIStream &is, int version)
Function used to read in object persistently.
GeneratorMap & generators()
Access the generator map.
double Nf(Energy Q) const
Return the number of scale dependent active flavours from the alpha_s object.
bool thePowhegDecayEmission
True if powheg style emissions are to be used in the decays.
static void Init()
The standard Init function used to initialize the interfaces.
void hardScales(Energy2 scale)
Setup the hard scales.
Ptr< DipoleEvolutionOrdering >::ptr theEvolutionOrdering
The evolution ordering considered.
bool realign()
Realign the event such as to have the incoming partons along thre beam axes.
void setMerger(Ptr< MergerBase >::ptr mh)
Set the pointer to the Merging Helper.
int verbosity
The verbosity level.
Energy theRenormalizationScaleFreeze
A freezing value for the renormalization scale.
Energy theFactorizationScaleFreeze
A freezing value for the factorization scale.
DipoleEventRecord & eventRecord()
Access the event record.
int theZBoundaries
The choice of z boundaries; 0 = restricted, 1 = open, 2 = mixed/other.
bool firstMCatNLOEmission
Perform the first MC@NLO emission only.
void decayConstituentReshuffle(PerturbativeProcessPtr decayProc)
Reshuffle to constituent mass shells.
tPPair cascade(tSubProPtr sub, XCombPtr xcomb, Energy optHardPt, Energy optCutoff)
The main method which manages the showering of a subprocess.
int currentReferenceWeight
Current reference weight used for partial unweighting of the subleading colour shower (updated each e...
void hardScalesSubleading(list< DipoleSplittingInfo > candidates, Energy hardPt)
Setup the hard scales of the dipoles after subleading emissions.
Ptr< AlphaSBase >::ptr theGlobalAlphaS
A global alpha_s to be used for all splitting kernels.
virtual void doinitrun()
Initialize this object.
bool didRealign
Whether or not we did realign the event.
set< long > theColouredOffShellInShower
PDG codes of the partons which can have an off-shell mass, this is fast storage for use during runnin...
DipoleShowerHandler()
The default constructor.
const cPDVector & particlesAfter() const
Return the particle data vector of the particles in the event, will be empty if subleading Nc is turn...
void resetReweight(Ptr< DipoleSplittingReweight >::tptr)
Reset the splitting reweight for all splitting kernels.
static ClassDescription< DipoleShowerHandler > initDipoleShowerHandler
The static object used to initialize the description of this class.
static bool firstWarn
True if no warnings have been issued yet.
const set< long > & offShellPartons()
Return the set of offshell parton ids.
Energy getWinner(DipoleSplittingInfo &winner, const Dipole &dip, pair< bool, bool > conf, Energy optHardPt=ZERO, Energy optCutoff=ZERO)
Get the winning splitting for the given dipole and configuration.
DipoleShowerHandler & operator=(const DipoleShowerHandler &)=delete
The assignment operator is private and must never be called.
DensityOperator & densityOperator()
Return the density operator from the event record.
const map< pair< vector< PDT::Colour >, pair< size_t, size_t > >, double > & correlatorMap() const
Return the colour matrix element correction map, will be empty if subleading Nc is turned off.
vector< Ptr< DipoleSplittingKernel >::ptr > kernels
The splitting kernels to be used.
Ptr< IntrinsicPtGenerator >::ptr intrinsicPtGenerator
The intrinsic pt generator to be used.
Ptr< ColourBasis >::tptr colourBasis()
Return the colour basis used in the density operator.
Energy getWinner(DipoleSplittingInfo &winner, Energy optHardPt=ZERO, Energy optCutoff=ZERO)
Get the winning splitting for the given dipole and configuration.
DipoleEventRecord theEventRecord
The evnt record used.
size_t winnerIndex
Winning candidate index.
Ptr< MergerBase >::ptr theMergingHelper
The merging helper takes care of merging multiple LO and NLO cross sections.
virtual IBPtr clone() const
Make a simple clone of this object.
void constituentReshuffle()
Reshuffle to constituent mass shells.
void addCandidates(PPair particles, list< DipoleSplittingInfo > &clist) const
Add the possible splitting candidates given a pair of emitting particles.
Energy nextSubleadingSplitting(Energy hardPt, Energy optHardPt, Energy optCutoff, const bool decay)
Generate the next subleading Nc improved splitting and return the scale.
virtual bool isReshuffling() const
Return true, if this cascade handler will perform reshuffling from hard process masses.
Ptr< ConstituentReshuffler >::ptr constituentReshuffler
The ConstituentReshuffler to be used.
double referenceWeight
Reference weight for the partial unweighting.
bool discardNoEmissions
Discard events which did not radiate.
Energy muPt
The shower hard scale for the last event encountered.
void getGenerators(const DipoleIndex &, Ptr< DipoleSplittingReweight >::tptr rw=Ptr< DipoleSplittingReweight >::tptr())
Build splitting generators for the given dipole index.
bool theDoCompensate
True, if sampler should apply compensation.
DipoleVertexRecord theVertexRecord
The vertex record.
void addSplitting(Ptr< DipoleSplittingKernel >::ptr sp)
Insert an additional splitting kernel.
vector< long > theInputColouredOffShellInShower
PDG codes of the partons which can have an off-shell mass, this is a vector that is interfaced so the...
size_t subleadingNcEmissionsLimit
Number of emissions to do subleading Nc corrections to.
double theDetuning
The detuning factor applied to the sampling overestimate kernel.
unsigned int nEmissions
Limit the number of emissions.
bool doPartialUnweighting
Switch on or off partial unweighting in the splitting generator.
unsigned long theFreezeGrid
Return the number of accepted points after which the grid should be frozen.
int densityOperatorEvolution
Integer used to set which method of evolving the density operator to use: 0 - Vijk is Eikonal but the...
Energy getWinner(DipoleSplittingInfo &winner, const DipoleIndex &index, double emitterX, double spectatorX, pair< bool, bool > conf, tPPtr emitter, tPPtr spectator, Energy startScale, Energy optHardPt=ZERO, Energy optCutoff=ZERO)
Get the winning splitting for the given dipole and configuration.
double cmecReweightFactor
Factor changing the acceptance probability for the veto algorithm.
bool doPartialUnweightingAtEmission
Switch on or off partial unweighting in after each subleading emission.
bool didRadiate
Whether or not we did radiate anything.
const DipoleEventRecord & eventRecord() const
Return the event record.
Ptr< DipoleEvolutionOrdering >::tptr evolutionOrdering() const
Return the evolution ordering.
DipoleSplittingInfo contains all parameters to generate a full dipole splitting.
Here is the documentation of the DipoleVertexRecord class.
The Dipole class is used by the dipole shower to represent a dipole of two coloured partons.
Definition: Dipole.h:30
This class is responsible for the handling of the merging process after the setup stage,...
Definition: Merger.h:95
This class is the main driver of the shower: it is responsible for the proper handling of all other s...
Definition: ShowerHandler.h:40
tDMPtr decay(PerturbativeProcessPtr, DecayProcessMap &decay, bool radPhotons=false) const
Decay a particle.
virtual void cascade()
The main method which manages the multiple interactions and starts the shower by calling cascade(sub,...
Generalized dipole splitting info to deal with subleading-N splittings.
-*- C++ -*-
ThePEG::Ptr< Particle >::transient_pointer tPPtr
pair< PPtr, PPtr > PPair
ThePEG::Ptr< InterfacedBase >::pointer IBPtr
ThePEG::Ptr< SubProcess >::transient_pointer tSubProPtr
vector< tPPtr > tPVector
pair< tPPtr, tPPtr > tPPair
constexpr ZeroUnit ZERO
constexpr auto sqr(const T &x) -> decltype(x *x)
vector< cPDPtr > cPDVector
ThePEG::Ptr< XComb >::pointer XCombPtr
Indicate a problem in the shower.
static string className()
static string library()