herwig
is hosted by
Hepforge
,
IPPP Durham
Herwig
7.3.0
MatrixElement
Matchbox
Cuts
MatchboxDeltaRCut.h
1
// -*- C++ -*-
2
//
3
// MatchboxDeltaRCut.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_MatchboxDeltaRCut_H
10
#define Herwig_MatchboxDeltaRCut_H
11
//
12
// This is the declaration of the MatchboxDeltaRCut 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
MatchboxDeltaRCut
:
public
TwoCutBase
{
32
33
public
:
34
38
MatchboxDeltaRCut
();
39
40
public
:
41
54
virtual
Energy
minDeltaMeasureCuts
(
tcPDPtr
,
tcPDPtr
)
const
{
return
ZERO
; }
55
// virtual Energy minDeltaMeasureCuts(tcPDPtr pi, tcPDPtr pj) const { return ZERO; }
56
61
virtual
Energy
minKTClus
(
tcPDPtr
,
tcPDPtr
)
const
{
return
ZERO
; }
62
// virtual Energy minKTClus(tcPDPtr pi, tcPDPtr pj) const { return ZERO; }
63
68
virtual
Energy2
minSij
(
tcPDPtr
,
tcPDPtr
)
const
{
return
ZERO
; }
69
// virtual Energy2 minSij(tcPDPtr pi, tcPDPtr pj) const { return ZERO; }
70
76
virtual
Energy2
minTij
(
tcPDPtr
,
tcPDPtr
)
const
{
return
ZERO
; }
77
// virtual Energy2 minTij(tcPDPtr pi, tcPDPtr po) const { return ZERO; }
78
84
virtual
double
minDeltaR
(
tcPDPtr
,
tcPDPtr
)
const
{
return
ZERO
; }
85
// virtual double minDeltaR(tcPDPtr pi, tcPDPtr pj) const { return ZERO; }
86
93
virtual
double
minDurham
(
tcPDPtr
,
tcPDPtr
)
const
{
return
0.0; }
94
// virtual double minDurham(tcPDPtr pi, tcPDPtr pj) const { return 0.0; }
95
102
virtual
bool
passCuts
(tcCutsPtr parent,
tcPDPtr
pitype,
tcPDPtr
pjtype,
103
LorentzMomentum pi, LorentzMomentum pj,
104
bool
inci =
false
,
bool
incj =
false
)
const
;
105
109
virtual
void
describe
()
const
;
111
112
public
:
113
117
double
deltaRMin
()
const
{
return
theDeltaRMin
; }
118
double
deltaRMax()
const
{
return
theDeltaRMax; }
119
123
double
deltaYMin
()
const
{
return
theDeltaYMin
; }
124
double
deltaYMax()
const
{
return
theDeltaYMax; }
125
129
double
deltaPhiMin
()
const
{
return
theDeltaPhiMin
; }
130
double
deltaPhiMax()
const
{
return
theDeltaPhiMax; }
131
136
Ptr<MatcherBase>::tptr
firstMatcher
()
const
{
return
theFirstMatcher
; }
137
Ptr<MatcherBase>::tptr secondMatcher()
const
{
return
theSecondMatcher; }
138
139
public
:
140
147
void
persistentOutput
(
PersistentOStream
& os)
const
;
148
154
void
persistentInput
(
PersistentIStream
& is,
int
version);
156
163
static
void
Init
();
164
165
protected
:
166
173
virtual
IBPtr
clone
()
const
;
174
179
virtual
IBPtr
fullclone
()
const
;
181
182
183
// If needed, insert declarations of virtual function defined in the
184
// InterfacedBase class here (using ThePEG-interfaced-decl in Emacs).
185
186
private
:
187
191
double
theDeltaRMin
;
192
double
theDeltaRMax;
193
197
double
theDeltaYMin
;
198
double
theDeltaYMax;
199
203
double
theDeltaPhiMin
;
204
double
theDeltaPhiMax;
205
210
Ptr<MatcherBase>::ptr
theFirstMatcher
;
211
Ptr<MatcherBase>::ptr theSecondMatcher;
212
213
private
:
214
219
MatchboxDeltaRCut
&
operator=
(
const
MatchboxDeltaRCut
&) =
delete
;
220
221
};
222
223
}
224
225
#endif
/* Herwig_MatchboxDeltaRCut_H */
Herwig::MatchboxDeltaRCut
MatchboxDeltaRCut implements cuts related to the separation in the legoplot plane.
Definition:
MatchboxDeltaRCut.h:31
Herwig::MatchboxDeltaRCut::theDeltaYMin
double theDeltaYMin
The minimum and maximum allowed rapidity separation.
Definition:
MatchboxDeltaRCut.h:197
Herwig::MatchboxDeltaRCut::theFirstMatcher
Ptr< MatcherBase >::ptr theFirstMatcher
Matchers for a pair of particles to cut on.
Definition:
MatchboxDeltaRCut.h:210
Herwig::MatchboxDeltaRCut::fullclone
virtual IBPtr fullclone() const
Make a clone of this object, possibly modifying the cloned object to make it sane.
Herwig::MatchboxDeltaRCut::theDeltaRMin
double theDeltaRMin
The minimum and maximum allowed legoplot separation.
Definition:
MatchboxDeltaRCut.h:191
Herwig::MatchboxDeltaRCut::minKTClus
virtual Energy minKTClus(tcPDPtr, tcPDPtr) const
Return the minimum allowed value of the longitudinally invariant -algorithms distance measure.
Definition:
MatchboxDeltaRCut.h:61
Herwig::MatchboxDeltaRCut::theDeltaPhiMin
double theDeltaPhiMin
The minimum and maximum allowed azimuthal separation.
Definition:
MatchboxDeltaRCut.h:203
Herwig::MatchboxDeltaRCut::persistentOutput
void persistentOutput(PersistentOStream &os) const
Function used to write out object persistently.
Herwig::MatchboxDeltaRCut::deltaPhiMin
double deltaPhiMin() const
Return the minimum and maximum allowed azimuthal separation.
Definition:
MatchboxDeltaRCut.h:129
Herwig::MatchboxDeltaRCut::MatchboxDeltaRCut
MatchboxDeltaRCut()
The default constructor.
Herwig::MatchboxDeltaRCut::operator=
MatchboxDeltaRCut & operator=(const MatchboxDeltaRCut &)=delete
The assignment operator is private and must never be called.
Herwig::MatchboxDeltaRCut::minDeltaR
virtual double minDeltaR(tcPDPtr, tcPDPtr) const
Return the minimum allowed value of of two outgoing partons of type pi and pj.
Definition:
MatchboxDeltaRCut.h:84
Herwig::MatchboxDeltaRCut::Init
static void Init()
The standard Init function used to initialize the interfaces.
Herwig::MatchboxDeltaRCut::minSij
virtual Energy2 minSij(tcPDPtr, tcPDPtr) const
Return the minimum allowed squared invariant mass of two outgoing partons of type pi and pj.
Definition:
MatchboxDeltaRCut.h:68
Herwig::MatchboxDeltaRCut::persistentInput
void persistentInput(PersistentIStream &is, int version)
Function used to read in object persistently.
Herwig::MatchboxDeltaRCut::deltaRMin
double deltaRMin() const
Return the minimum and maximum allowed legoplot separation.
Definition:
MatchboxDeltaRCut.h:117
Herwig::MatchboxDeltaRCut::deltaYMin
double deltaYMin() const
Return the minimum and maximum allowed rapidity separation.
Definition:
MatchboxDeltaRCut.h:123
Herwig::MatchboxDeltaRCut::minDurham
virtual double minDurham(tcPDPtr, tcPDPtr) const
Return the minimum allowed value of the Durham -algorithms distance measure.
Definition:
MatchboxDeltaRCut.h:93
Herwig::MatchboxDeltaRCut::describe
virtual void describe() const
Describe the currently active cuts in the log file.
Herwig::MatchboxDeltaRCut::minDeltaMeasureCuts
virtual Energy minDeltaMeasureCuts(tcPDPtr, tcPDPtr) const
Return the minimum allowed value of the longitudinally invariant -algorithms distance measure.
Definition:
MatchboxDeltaRCut.h:54
Herwig::MatchboxDeltaRCut::passCuts
virtual bool passCuts(tcCutsPtr parent, tcPDPtr pitype, tcPDPtr pjtype, LorentzMomentum pi, LorentzMomentum pj, bool inci=false, bool incj=false) const
Return true if a pair of particles with type pitype and pjtype and momenta pi and pj respectively pas...
Herwig::MatchboxDeltaRCut::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:
MatchboxDeltaRCut.h:76
Herwig::MatchboxDeltaRCut::firstMatcher
Ptr< MatcherBase >::tptr firstMatcher() const
Return the matchers for a pair of particles to cut on.
Definition:
MatchboxDeltaRCut.h:136
Herwig::MatchboxDeltaRCut::clone
virtual IBPtr clone() const
Make a simple clone of this object.
ThePEG::PersistentIStream
ThePEG::PersistentOStream
ThePEG::TwoCutBase
Herwig
-*- C++ -*-
Definition:
BasicConsistency.h:17
ThePEG
ThePEG::IBPtr
ThePEG::Ptr< InterfacedBase >::pointer IBPtr
ThePEG::ZERO
constexpr ZeroUnit ZERO
ThePEG::tcPDPtr
ThePEG::Ptr< ParticleData >::transient_const_pointer tcPDPtr
Generated on Thu Jun 20 2024 17:50:52 for Herwig by
1.9.6