herwig
is hosted by
Hepforge
,
IPPP Durham
Herwig
7.2.1
MatrixElement
Matchbox
Cuts
InvariantMassCut.h
1
// -*- C++ -*-
2
//
3
// InvariantMassCut.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_InvariantMassCut_H
10
#define Herwig_InvariantMassCut_H
11
//
12
// This is the declaration of the InvariantMassCut class.
13
//
14
15
#include "ThePEG/Cuts/TwoCutBase.h"
16
#include "ThePEG/PDT/MatcherBase.h"
17
18
namespace
Herwig
{
19
20
using namespace
ThePEG
;
21
31
class
InvariantMassCut
:
public
TwoCutBase
{
32
33
public
:
34
40
InvariantMassCut
() :
41
theMinMass(0*GeV), theMaxMass(Constants::
MaxEnergy
),
42
theSameFlavourOnly(false), theOppositeSignOnly(false) {}
44
45
public
:
46
55
virtual
bool
passCuts(tcCutsPtr parent,
tcPDPtr
pitype,
tcPDPtr
pjtype,
56
LorentzMomentum
pi
, LorentzMomentum pj,
57
bool
inci =
false
,
bool
incj =
false
)
const
;
58
63
virtual
Energy2
minSij
(
tcPDPtr
,
tcPDPtr
)
const
{
return
ZERO
; }
64
70
virtual
Energy2
minTij
(
tcPDPtr
,
tcPDPtr
)
const
{
return
ZERO
; }
71
77
virtual
double
minDeltaR
(
tcPDPtr
,
tcPDPtr
)
const
{
return
ZERO
; }
78
89
virtual
Energy
minKTClus
(
tcPDPtr
,
tcPDPtr
)
const
{
return
ZERO
; }
90
97
virtual
double
minDurham
(
tcPDPtr
,
tcPDPtr
)
const
{
return
ZERO
; }
98
100
104
virtual
void
describe()
const
;
105
106
public
:
107
111
Energy
minMass
()
const
{
return
theMinMass; }
112
116
Energy
maxMass
()
const
{
return
theMaxMass; }
117
121
bool
sameFlavourOnly
()
const
{
return
theSameFlavourOnly; }
122
126
bool
oppositeSignOnly
()
const
{
return
theOppositeSignOnly; }
127
132
Ptr<MatcherBase>::tptr
firstMatcher
()
const
{
return
theFirstMatcher; }
133
Ptr<MatcherBase>::tptr
secondMatcher()
const
{
return
theSecondMatcher; }
134
135
protected
:
136
140
int
family(
long
id
)
const
;
141
142
public
:
143
150
void
persistentOutput(
PersistentOStream
& os)
const
;
151
157
void
persistentInput(
PersistentIStream
& is,
int
version);
159
166
static
void
Init();
167
168
protected
:
169
176
virtual
IBPtr
clone()
const
;
177
182
virtual
IBPtr
fullclone()
const
;
184
185
private
:
186
190
Energy
theMinMass
;
191
195
Energy
theMaxMass
;
196
201
bool
theSameFlavourOnly
;
202
207
bool
theOppositeSignOnly
;
208
213
Ptr<MatcherBase>::ptr
theFirstMatcher
;
214
Ptr<MatcherBase>::ptr
theSecondMatcher;
215
216
private
:
217
222
InvariantMassCut
& operator=(
const
InvariantMassCut
&) =
delete
;
223
224
};
225
226
}
227
228
#endif
/* Herwig_InvariantMassCut_H */
Herwig::InvariantMassCut::theMaxMass
Energy theMaxMass
The maximal allowed mass cut value.
Definition:
InvariantMassCut.h:195
ThePEG::PersistentIStream
Herwig::InvariantMassCut::minSij
virtual Energy2 minSij(tcPDPtr, tcPDPtr) const
Return the minimum allowed squared invariant mass of two outgoing partons of type pi and pj...
Definition:
InvariantMassCut.h:63
ThePEG::Constants::MaxEnergy
constexpr Energy MaxEnergy
Ptr< MatcherBase >::tptr
transient_pointer tptr
Herwig::InvariantMassCut::oppositeSignOnly
bool oppositeSignOnly() const
Return whether cut acts on opposite-sign fermions only.
Definition:
InvariantMassCut.h:126
ThePEG::PersistentOStream
ThePEG::IBPtr
ThePEG::Ptr< InterfacedBase >::pointer IBPtr
Herwig::InvariantMassCut::minDurham
virtual double minDurham(tcPDPtr, tcPDPtr) const
Return the minimum allowed value of the Durham -algorithms distance measure.
Definition:
InvariantMassCut.h:97
ThePEG::Constants::pi
constexpr double pi
ThePEG
Herwig::InvariantMassCut::theFirstMatcher
Ptr< MatcherBase >::ptr theFirstMatcher
Matchers for a pair of particles to cut on.
Definition:
InvariantMassCut.h:213
Herwig::InvariantMassCut::theOppositeSignOnly
bool theOppositeSignOnly
Whether the cut is active on opposite-sign fermions only (ignored for pairs not consisting of two fer...
Definition:
InvariantMassCut.h:207
Herwig::InvariantMassCut::minMass
Energy minMass() const
Return the minimal allowed invariant mass.
Definition:
InvariantMassCut.h:111
Herwig::InvariantMassCut
This class implements a cut on the invariant mass.
Definition:
InvariantMassCut.h:31
Herwig::InvariantMassCut::maxMass
Energy maxMass() const
Return the maximal allowed invariant mass.
Definition:
InvariantMassCut.h:116
Herwig::InvariantMassCut::minDeltaR
virtual double minDeltaR(tcPDPtr, tcPDPtr) const
Return the minimum allowed value of of two outgoing partons of type pi and pj.
Definition:
InvariantMassCut.h:77
Herwig::InvariantMassCut::minTij
virtual Energy2 minTij(tcPDPtr, tcPDPtr) const
Return the minimum allowed value of the negative of the squared invariant mass of an incoming parton ...
Definition:
InvariantMassCut.h:70
Herwig::InvariantMassCut::theMinMass
Energy theMinMass
The minimal allowed mass cut value.
Definition:
InvariantMassCut.h:190
ThePEG::tcPDPtr
ThePEG::Ptr< ParticleData >::transient_const_pointer tcPDPtr
ThePEG::TwoCutBase
Herwig::InvariantMassCut::InvariantMassCut
InvariantMassCut()
The default constructor.
Definition:
InvariantMassCut.h:40
Herwig::InvariantMassCut::minKTClus
virtual Energy minKTClus(tcPDPtr, tcPDPtr) const
Return the minimum allowed value of the longitudinally invariant -algorithms distance measure...
Definition:
InvariantMassCut.h:89
Herwig::InvariantMassCut::theSameFlavourOnly
bool theSameFlavourOnly
Whether the cut is active on same-flavour fermions only (ignored for pairs not consisting of two ferm...
Definition:
InvariantMassCut.h:201
Herwig
-*- C++ -*-
Definition:
BasicConsistency.h:17
Herwig::InvariantMassCut::firstMatcher
Ptr< MatcherBase >::tptr firstMatcher() const
Return the matchers for a pair of particles to cut on.
Definition:
InvariantMassCut.h:132
Herwig::InvariantMassCut::sameFlavourOnly
bool sameFlavourOnly() const
Return whether cut acts on same-flavour fermions only.
Definition:
InvariantMassCut.h:121
ZERO
constexpr ZeroUnit ZERO
Ptr< MatcherBase >::ptr
pointer ptr
Generated on Sat Apr 11 2020 14:50:30 for Herwig by
1.8.13