herwig
is hosted by
Hepforge
,
IPPP Durham
Herwig
7.3.0
Shower
Dipole
Kernels
DipoleSplittingKernel.h
1
// -*- C++ -*-
2
//
3
// DipoleSplittingKernel.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_DipoleSplittingKernel_H
10
#define HERWIG_DipoleSplittingKernel_H
11
//
12
// This is the declaration of the DipoleSplittingKernel class.
13
//
14
15
#include "ThePEG/Handlers/HandlerBase.h"
16
#include "ThePEG/StandardModel/AlphaSBase.h"
17
#include "ThePEG/PDF/PDF.h"
18
19
#include "Herwig/Shower/Dipole/Utility/PDFRatio.h"
20
#include "Herwig/Shower/Dipole/Base/DipoleSplittingInfo.h"
21
#include "Herwig/Shower/Dipole/Kinematics/DipoleSplittingKinematics.h"
22
23
#include "ThePEG/EventRecord/RhoDMatrix.h"
24
#include "Herwig/Decay/DecayMatrixElement.h"
25
#include "Herwig/Decay/TwoBodyDecayMatrixElement.h"
26
27
namespace
Herwig
{
28
29
using namespace
ThePEG
;
30
41
class
DipoleSplittingKernel
:
public
HandlerBase
{
42
43
public
:
44
48
DipoleSplittingKernel
();
49
50
public
:
51
55
Ptr<AlphaSBase>::tptr
alphaS
()
const
{
return
theAlphaS
; }
56
60
void
alphaS
(Ptr<AlphaSBase>::tptr ap) {
theAlphaS
= ap; }
61
65
Ptr<DipoleSplittingKinematics>::tptr
splittingKinematics
()
const
{
66
return
theSplittingKinematics
;
67
}
68
72
Ptr<DipoleMCCheck>::ptr
mcCheck
()
const
{
return
theMCCheck
; }
73
77
void
splittingKinematics
(Ptr<DipoleSplittingKinematics>::tptr sp) {
78
theSplittingKinematics
= sp;
79
}
80
84
Ptr<PDFRatio>::tptr
pdfRatio
()
const
{
return
thePDFRatio
; }
85
89
void
pdfRatio
(Ptr<PDFRatio>::tptr sp) {
thePDFRatio
= sp; }
90
99
virtual
int
nDimAdditional
()
const
{
return
0; }
100
104
void
renormalizationScaleFreeze
(Energy s) {
theRenormalizationScaleFreeze
= s; }
105
109
void
factorizationScaleFreeze
(Energy s) {
theFactorizationScaleFreeze
= s; }
110
114
Energy
renormalizationScaleFreeze
()
const
{
return
theRenormalizationScaleFreeze
; }
115
119
Energy
factorizationScaleFreeze
()
const
{
return
theFactorizationScaleFreeze
; }
120
121
public
:
122
127
virtual
bool
canHandle
(
const
DipoleIndex
&)
const
= 0;
128
134
virtual
bool
canHandleEquivalent
(
const
DipoleIndex
& a,
135
const
DipoleSplittingKernel
& sk,
136
const
DipoleIndex
& b)
const
= 0;
137
142
virtual
tcPDPtr
emitter
(
const
DipoleIndex
&)
const
= 0;
143
148
virtual
tcPDPtr
emission
(
const
DipoleIndex
&)
const
= 0;
149
154
virtual
tcPDPtr
spectator
(
const
DipoleIndex
&)
const
= 0;
155
160
PDPtr
flavour
()
const
{
return
theFlavour
; }
161
166
bool
strictLargeN
()
const
{
return
theStrictLargeN
; }
167
168
public
:
169
174
virtual
void
startPresampling
(
const
DipoleIndex
&) {
175
presampling
=
true
;
176
}
177
182
virtual
void
stopPresampling
(
const
DipoleIndex
&) {
183
presampling
=
false
;
184
}
185
190
unsigned
long
presamplingPoints
()
const
{
return
thePresamplingPoints
; }
191
196
unsigned
long
maxtry
()
const
{
return
theMaxtry
; }
197
202
unsigned
long
freezeGrid
()
const
{
return
theFreezeGrid
; }
203
208
void
freezeGrid
(
unsigned
long
n) {
theFreezeGrid
= n; }
209
213
void
detuning
(
double
d) {
theDetuning
= d; }
214
218
double
detuning
()
const
{
return
theDetuning
; }
219
224
virtual
double
evaluate
(
const
DipoleSplittingInfo
&)
const
= 0;
225
230
virtual
vector< pair<int, Complex> >
generatePhi
(
const
DipoleSplittingInfo
& dInfo,
const
RhoDMatrix
& rho)
const
= 0;
231
235
virtual
DecayMEPtr
matrixElement
(
const
DipoleSplittingInfo
& dInfo)
const
= 0;
236
240
void
clearAlphaPDFCache
()
const
{
241
theAlphaSCache
.clear();
242
thePDFCache
.clear();
243
}
244
249
virtual
void
accept
(
const
DipoleSplittingInfo
&,
double
,
double
, map<string,double>&)
const
;
250
255
virtual
void
veto
(
const
DipoleSplittingInfo
&,
double
,
double
, map<string,double>&)
const
;
256
264
virtual
bool
haveOverestimate
(
const
DipoleSplittingInfo
&)
const
{
return
false
; }
265
270
virtual
double
overestimate
(
const
DipoleSplittingInfo
&)
const
{
return
-1.; }
271
278
virtual
double
invertOverestimateIntegral
(
const
DipoleSplittingInfo
&,
double
)
const
{
279
return
-1.;
280
}
281
285
bool
useThisKernel
()
const
{
286
return
theUseThisKernel
;
287
}
288
289
public
:
290
294
double
factorizationScaleFactor
()
const
{
return
theFactorizationScaleFactor
; }
295
299
void
factorizationScaleFactor
(
double
f) {
theFactorizationScaleFactor
= f; }
300
304
double
renormalizationScaleFactor
()
const
{
return
theRenormalizationScaleFactor
; }
305
309
void
renormalizationScaleFactor
(
double
f) {
theRenormalizationScaleFactor
= f; }
310
314
unsigned
int
cmwScheme
()
const
{
return
theCMWScheme
;}
315
316
protected
:
317
321
double
alphaPDF
(
const
DipoleSplittingInfo
&,
322
Energy optScale =
ZERO
,
323
double
rScaleFactor = 1.0,
324
double
fScaleFactor = 1.0)
const
;
325
330
bool
virtualitySplittingScale
()
const
{
return
theVirtualitySplittingScale
; }
331
332
public
:
333
340
void
persistentOutput
(
PersistentOStream
& os)
const
;
341
347
void
persistentInput
(
PersistentIStream
& is,
int
version);
349
356
static
void
Init
();
357
358
359
// If needed, insert declarations of virtual function defined in the
360
// InterfacedBase class here (using ThePEG-interfaced-decl in Emacs).
361
362
private
:
363
367
Ptr<AlphaSBase>::ptr
theAlphaS
;
368
373
Energy
theScreeningScale
;
374
378
Ptr<DipoleSplittingKinematics>::ptr
theSplittingKinematics
;
379
384
Ptr<PDFRatio>::ptr
thePDFRatio
;
385
390
unsigned
long
thePresamplingPoints
;
391
396
unsigned
long
theMaxtry
;
397
402
static
double
theMaxPDFRatio
;
403
408
unsigned
long
theFreezeGrid
;
409
413
double
theDetuning
;
414
419
PDPtr
theFlavour
;
420
424
Ptr<DipoleMCCheck>::ptr
theMCCheck
;
425
430
bool
theStrictLargeN
;
431
435
double
theFactorizationScaleFactor
;
436
440
double
theRenormalizationScaleFactor
;
441
445
Energy
theRenormalizationScaleFreeze
;
446
450
Energy
theFactorizationScaleFreeze
;
451
456
bool
theVirtualitySplittingScale
;
457
462
unsigned
int
theCMWScheme
=0;
463
467
mutable
map<double,double>
theAlphaSCache
;
468
472
mutable
map<double,double>
thePDFCache
;
473
477
bool
presampling
;
478
482
bool
theUseThisKernel
=
true
;
483
484
485
private
:
486
491
static
AbstractClassDescription<DipoleSplittingKernel>
initDipoleSplittingKernel
;
492
497
DipoleSplittingKernel
&
operator=
(
const
DipoleSplittingKernel
&) =
delete
;
498
499
};
500
501
}
502
503
#include "ThePEG/Utilities/ClassTraits.h"
504
505
namespace
ThePEG
{
506
511
template
<>
512
struct
BaseClassTrait
<
Herwig
::DipoleSplittingKernel,1> {
514
typedef
HandlerBase
NthBase
;
515
};
516
519
template
<>
520
struct
ClassTraits
<
Herwig
::DipoleSplittingKernel>
521
:
public
ClassTraitsBase
<Herwig::DipoleSplittingKernel> {
523
static
string
className
() {
return
"Herwig::DipoleSplittingKernel"
; }
531
static
string
library
() {
return
"HwDipoleShower.so"
; }
532
};
533
536
}
537
538
#endif
/* HERWIG_DipoleSplittingKernel_H */
Herwig::DipoleIndex
DipoleIndex is used to index splitting generators for a particular dipole.
Definition:
DipoleSplittingInfo.h:35
Herwig::DipoleSplittingInfo
DipoleSplittingInfo contains all parameters to generate a full dipole splitting.
Definition:
DipoleSplittingInfo.h:209
Herwig::DipoleSplittingKernel
DipoleSplittingKernel is the base class for all kernels used within the dipole shower.
Definition:
DipoleSplittingKernel.h:41
Herwig::DipoleSplittingKernel::strictLargeN
bool strictLargeN() const
Return true, if this splitting kernel is supposed to work in a strict large-N limit,...
Definition:
DipoleSplittingKernel.h:166
Herwig::DipoleSplittingKernel::persistentInput
void persistentInput(PersistentIStream &is, int version)
Function used to read in object persistently.
Herwig::DipoleSplittingKernel::renormalizationScaleFactor
void renormalizationScaleFactor(double f)
Set the renormalization scale factor.
Definition:
DipoleSplittingKernel.h:309
Herwig::DipoleSplittingKernel::overestimate
virtual double overestimate(const DipoleSplittingInfo &) const
Return the overestimate to this splitting kernel for the given dipole splitting.
Definition:
DipoleSplittingKernel.h:270
Herwig::DipoleSplittingKernel::alphaS
Ptr< AlphaSBase >::tptr alphaS() const
Return the alpha_s to be used.
Definition:
DipoleSplittingKernel.h:55
Herwig::DipoleSplittingKernel::factorizationScaleFreeze
Energy factorizationScaleFreeze() const
Get the freezing value for the factorization scale.
Definition:
DipoleSplittingKernel.h:119
Herwig::DipoleSplittingKernel::thePDFRatio
Ptr< PDFRatio >::ptr thePDFRatio
An optional PDF ratio object to be used when evaluating this kernel.
Definition:
DipoleSplittingKernel.h:384
Herwig::DipoleSplittingKernel::theRenormalizationScaleFactor
double theRenormalizationScaleFactor
The renormalization scale factor.
Definition:
DipoleSplittingKernel.h:440
Herwig::DipoleSplittingKernel::theScreeningScale
Energy theScreeningScale
An optional 'colour screening' scale for alternative intrinsic pt generation.
Definition:
DipoleSplittingKernel.h:373
Herwig::DipoleSplittingKernel::spectator
virtual tcPDPtr spectator(const DipoleIndex &) const =0
Return the spectator data after splitting, given a dipole index.
Herwig::DipoleSplittingKernel::theMCCheck
Ptr< DipoleMCCheck >::ptr theMCCheck
Pointer to a check histogram object.
Definition:
DipoleSplittingKernel.h:424
Herwig::DipoleSplittingKernel::presampling
bool presampling
True, if we are presampling.
Definition:
DipoleSplittingKernel.h:477
Herwig::DipoleSplittingKernel::theFlavour
PDPtr theFlavour
The flavour produced, if this cannot be determined from the dipole.
Definition:
DipoleSplittingKernel.h:419
Herwig::DipoleSplittingKernel::flavour
PDPtr flavour() const
Return the flavour produced, if this cannot be determined from the dipole.
Definition:
DipoleSplittingKernel.h:160
Herwig::DipoleSplittingKernel::thePDFCache
map< double, double > thePDFCache
Cache for PDF evaluations.
Definition:
DipoleSplittingKernel.h:472
Herwig::DipoleSplittingKernel::freezeGrid
unsigned long freezeGrid() const
Return the number of accepted points after which the grid should be frozen.
Definition:
DipoleSplittingKernel.h:202
Herwig::DipoleSplittingKernel::canHandleEquivalent
virtual bool canHandleEquivalent(const DipoleIndex &a, const DipoleSplittingKernel &sk, const DipoleIndex &b) const =0
Return true, if this splitting kernel is the same for the given index a, as the given splitting kerne...
Herwig::DipoleSplittingKernel::thePresamplingPoints
unsigned long thePresamplingPoints
The number of points to presample this splitting generator.
Definition:
DipoleSplittingKernel.h:390
Herwig::DipoleSplittingKernel::theFactorizationScaleFreeze
Energy theFactorizationScaleFreeze
A freezing value for the factorization scale.
Definition:
DipoleSplittingKernel.h:450
Herwig::DipoleSplittingKernel::evaluate
virtual double evaluate(const DipoleSplittingInfo &) const =0
Evaluate this splitting kernel for the given dipole splitting.
Herwig::DipoleSplittingKernel::splittingKinematics
Ptr< DipoleSplittingKinematics >::tptr splittingKinematics() const
Return the splitting kinematics object.
Definition:
DipoleSplittingKernel.h:65
Herwig::DipoleSplittingKernel::freezeGrid
void freezeGrid(unsigned long n)
Set the number of accepted points after which the grid should be frozen.
Definition:
DipoleSplittingKernel.h:208
Herwig::DipoleSplittingKernel::detuning
double detuning() const
Return the detuning factor applied to the sampling overestimate kernel.
Definition:
DipoleSplittingKernel.h:218
Herwig::DipoleSplittingKernel::renormalizationScaleFreeze
void renormalizationScaleFreeze(Energy s)
Set the freezing value for the renormalization scale.
Definition:
DipoleSplittingKernel.h:104
Herwig::DipoleSplittingKernel::detuning
void detuning(double d)
Set a detuning factor to be applied to the sampling overestimate kernel.
Definition:
DipoleSplittingKernel.h:213
Herwig::DipoleSplittingKernel::accept
virtual void accept(const DipoleSplittingInfo &, double, double, map< string, double > &) const
Update the variations vector at the given splitting using the indicated kernel and overestimate value...
Herwig::DipoleSplittingKernel::pdfRatio
void pdfRatio(Ptr< PDFRatio >::tptr sp)
Set the PDFRatio object.
Definition:
DipoleSplittingKernel.h:89
Herwig::DipoleSplittingKernel::theMaxtry
unsigned long theMaxtry
The maximum number of trials to generate a splitting.
Definition:
DipoleSplittingKernel.h:396
Herwig::DipoleSplittingKernel::theFactorizationScaleFactor
double theFactorizationScaleFactor
The factorization scale factor.
Definition:
DipoleSplittingKernel.h:435
Herwig::DipoleSplittingKernel::splittingKinematics
void splittingKinematics(Ptr< DipoleSplittingKinematics >::tptr sp)
Set the splitting kinematics object.
Definition:
DipoleSplittingKernel.h:77
Herwig::DipoleSplittingKernel::theSplittingKinematics
Ptr< DipoleSplittingKinematics >::ptr theSplittingKinematics
The splitting kinematics to be used.
Definition:
DipoleSplittingKernel.h:378
Herwig::DipoleSplittingKernel::maxtry
unsigned long maxtry() const
Return the maximum number of trials to generate a splitting.
Definition:
DipoleSplittingKernel.h:196
Herwig::DipoleSplittingKernel::theCMWScheme
unsigned int theCMWScheme
Implementing CMW in the kernels.
Definition:
DipoleSplittingKernel.h:462
Herwig::DipoleSplittingKernel::cmwScheme
unsigned int cmwScheme() const
Return the CMW sheme used by the kernel.
Definition:
DipoleSplittingKernel.h:314
Herwig::DipoleSplittingKernel::nDimAdditional
virtual int nDimAdditional() const
Return the number of additional parameter random variables needed to evaluate this kernel except the ...
Definition:
DipoleSplittingKernel.h:99
Herwig::DipoleSplittingKernel::renormalizationScaleFactor
double renormalizationScaleFactor() const
Get the renormalization scale factor.
Definition:
DipoleSplittingKernel.h:304
Herwig::DipoleSplittingKernel::alphaS
void alphaS(Ptr< AlphaSBase >::tptr ap)
Set the alpha_s to be used.
Definition:
DipoleSplittingKernel.h:60
Herwig::DipoleSplittingKernel::startPresampling
virtual void startPresampling(const DipoleIndex &)
Inform this splitting kernel, that it is being presampled until a call to stopPresampling.
Definition:
DipoleSplittingKernel.h:174
Herwig::DipoleSplittingKernel::renormalizationScaleFreeze
Energy renormalizationScaleFreeze() const
Get the freezing value for the renormalization scale.
Definition:
DipoleSplittingKernel.h:114
Herwig::DipoleSplittingKernel::generatePhi
virtual vector< pair< int, Complex > > generatePhi(const DipoleSplittingInfo &dInfo, const RhoDMatrix &rho) const =0
Evaluate rho_ii' V_ijk V*_i'jk / equivalent for initial-state splitting, required for generating spin...
Herwig::DipoleSplittingKernel::factorizationScaleFreeze
void factorizationScaleFreeze(Energy s)
Set the freezing value for the factorization scale.
Definition:
DipoleSplittingKernel.h:109
Herwig::DipoleSplittingKernel::initDipoleSplittingKernel
static AbstractClassDescription< DipoleSplittingKernel > initDipoleSplittingKernel
The static object used to initialize the description of this class.
Definition:
DipoleSplittingKernel.h:491
Herwig::DipoleSplittingKernel::stopPresampling
virtual void stopPresampling(const DipoleIndex &)
Inform this splitting kernel, that it is not being presampled until a call to startPresampling.
Definition:
DipoleSplittingKernel.h:182
Herwig::DipoleSplittingKernel::emitter
virtual tcPDPtr emitter(const DipoleIndex &) const =0
Return the emitter data after splitting, given a dipole index.
Herwig::DipoleSplittingKernel::persistentOutput
void persistentOutput(PersistentOStream &os) const
Function used to write out object persistently.
Herwig::DipoleSplittingKernel::theMaxPDFRatio
static double theMaxPDFRatio
The maximum value for any pdf ratio.
Definition:
DipoleSplittingKernel.h:402
Herwig::DipoleSplittingKernel::factorizationScaleFactor
void factorizationScaleFactor(double f)
Set the factorization scale factor.
Definition:
DipoleSplittingKernel.h:299
Herwig::DipoleSplittingKernel::haveOverestimate
virtual bool haveOverestimate(const DipoleSplittingInfo &) const
Return true, if this kernel is capable of delivering an overestimate to the kernel,...
Definition:
DipoleSplittingKernel.h:264
Herwig::DipoleSplittingKernel::theAlphaSCache
map< double, double > theAlphaSCache
Cache for alphas evaluations.
Definition:
DipoleSplittingKernel.h:467
Herwig::DipoleSplittingKernel::invertOverestimateIntegral
virtual double invertOverestimateIntegral(const DipoleSplittingInfo &, double) const
Invert the integral over the overestimate w.r.t.
Definition:
DipoleSplittingKernel.h:278
Herwig::DipoleSplittingKernel::canHandle
virtual bool canHandle(const DipoleIndex &) const =0
Return true, if this splitting kernel applies to the given dipole index.
Herwig::DipoleSplittingKernel::clearAlphaPDFCache
void clearAlphaPDFCache() const
Clear the alphaPDF cache.
Definition:
DipoleSplittingKernel.h:240
Herwig::DipoleSplittingKernel::theAlphaS
Ptr< AlphaSBase >::ptr theAlphaS
The alpha_s to be used.
Definition:
DipoleSplittingKernel.h:367
Herwig::DipoleSplittingKernel::operator=
DipoleSplittingKernel & operator=(const DipoleSplittingKernel &)=delete
The assignment operator is private and must never be called.
Herwig::DipoleSplittingKernel::DipoleSplittingKernel
DipoleSplittingKernel()
The default constructor.
Herwig::DipoleSplittingKernel::virtualitySplittingScale
bool virtualitySplittingScale() const
Return true, if the virtuality of the splitting should be used as the argument of alphas rather than ...
Definition:
DipoleSplittingKernel.h:330
Herwig::DipoleSplittingKernel::mcCheck
Ptr< DipoleMCCheck >::ptr mcCheck() const
Return the mc check object.
Definition:
DipoleSplittingKernel.h:72
Herwig::DipoleSplittingKernel::alphaPDF
double alphaPDF(const DipoleSplittingInfo &, Energy optScale=ZERO, double rScaleFactor=1.0, double fScaleFactor=1.0) const
Return the common factor of (alphas/2pi)*(pdf ratio)
Herwig::DipoleSplittingKernel::factorizationScaleFactor
double factorizationScaleFactor() const
Get the factorization scale factor.
Definition:
DipoleSplittingKernel.h:294
Herwig::DipoleSplittingKernel::emission
virtual tcPDPtr emission(const DipoleIndex &) const =0
Return the emission data after splitting, given a dipole index.
Herwig::DipoleSplittingKernel::presamplingPoints
unsigned long presamplingPoints() const
Return the number of points to presample this splitting generator.
Definition:
DipoleSplittingKernel.h:190
Herwig::DipoleSplittingKernel::useThisKernel
bool useThisKernel() const
Definition:
DipoleSplittingKernel.h:285
Herwig::DipoleSplittingKernel::pdfRatio
Ptr< PDFRatio >::tptr pdfRatio() const
Return the PDFRatio object.
Definition:
DipoleSplittingKernel.h:84
Herwig::DipoleSplittingKernel::theFreezeGrid
unsigned long theFreezeGrid
Return the number of accepted points after which the grid should be frozen.
Definition:
DipoleSplittingKernel.h:408
Herwig::DipoleSplittingKernel::veto
virtual void veto(const DipoleSplittingInfo &, double, double, map< string, double > &) const
Update the variations vector at the given splitting using the indicated kernel and overestimate value...
Herwig::DipoleSplittingKernel::theDetuning
double theDetuning
The detuning factor applied to the sampling overestimate kernel.
Definition:
DipoleSplittingKernel.h:413
Herwig::DipoleSplittingKernel::theUseThisKernel
bool theUseThisKernel
True, if the kernel should be used.
Definition:
DipoleSplittingKernel.h:482
Herwig::DipoleSplittingKernel::theVirtualitySplittingScale
bool theVirtualitySplittingScale
True, if the virtuality of the splitting should be used as the argument of alphas rather than the pt.
Definition:
DipoleSplittingKernel.h:456
Herwig::DipoleSplittingKernel::Init
static void Init()
The standard Init function used to initialize the interfaces.
Herwig::DipoleSplittingKernel::matrixElement
virtual DecayMEPtr matrixElement(const DipoleSplittingInfo &dInfo) const =0
Return the completely spin-unaveraged (i.e.
Herwig::DipoleSplittingKernel::theStrictLargeN
bool theStrictLargeN
True, if this splitting kernel is supposed to work in a strict large-N limit, i.e.
Definition:
DipoleSplittingKernel.h:430
Herwig::DipoleSplittingKernel::theRenormalizationScaleFreeze
Energy theRenormalizationScaleFreeze
A freezing value for the renormalization scale.
Definition:
DipoleSplittingKernel.h:445
ThePEG::AbstractClassDescription
ThePEG::HandlerBase
ThePEG::PersistentIStream
ThePEG::PersistentOStream
ThePEG::RhoDMatrix
Herwig
-*- C++ -*-
Definition:
BasicConsistency.h:17
ThePEG
ThePEG::ZERO
constexpr ZeroUnit ZERO
ThePEG::tcPDPtr
ThePEG::Ptr< ParticleData >::transient_const_pointer tcPDPtr
ThePEG::PDPtr
ThePEG::Ptr< ParticleData >::pointer PDPtr
ThePEG::BaseClassTrait
ThePEG::BaseClassTrait::NthBase
int NthBase
ThePEG::ClassTraitsBase
ThePEG::ClassTraitsBase::className
static string className()
ThePEG::ClassTraitsBase::library
static string library()
ThePEG::ClassTraits
Generated on Thu Jun 20 2024 17:50:53 for Herwig by
1.9.6