herwig is hosted by Hepforge, IPPP Durham
Herwig  7.2.1
IdentifiedParticleCut.h
1 // -*- C++ -*-
2 //
3 // IdentifiedParticleCut.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_IdentifiedParticleCut_H
10 #define Herwig_IdentifiedParticleCut_H
11 //
12 // This is the declaration of the IdentifiedParticleCut class.
13 //
14 
15 #include "ThePEG/Cuts/OneCutBase.h"
16 #include "ThePEG/PDT/MatcherBase.h"
17 
18 namespace Herwig {
19 
20 using namespace ThePEG;
21 
32 
33 public:
34 
41 
45  virtual ~IdentifiedParticleCut();
47 
48 public:
49 
56  virtual Energy minKT(tcPDPtr) const { return ZERO; }
57 
63  virtual double minEta(tcPDPtr) const { return -Constants::MaxRapidity; }
64 
70  virtual double maxEta(tcPDPtr) const { return Constants::MaxRapidity; }
71 
77  virtual bool passCuts(tcCutsPtr parent,
78  tcPDPtr ptype, LorentzMomentum p) const;
79 
83  virtual void describe() const;
85 
86 public:
87 
91  Energy ptMin() const { return thePtMin; }
92 
96  Energy ptMax() const { return thePtMax; }
97 
101  const vector<pair<double,double> >& yRanges() const { return theYRanges; }
102 
106  Ptr<MatcherBase>::tptr matcher() const { return theMatcher; }
107 
108 public:
109 
116  void persistentOutput(PersistentOStream & os) const;
117 
123  void persistentInput(PersistentIStream & is, int version);
125 
132  static void Init();
133 
134 protected:
135 
142  virtual IBPtr clone() const;
143 
148  virtual IBPtr fullclone() const;
150 
151 
152 // If needed, insert declarations of virtual function defined in the
153 // InterfacedBase class here (using ThePEG-interfaced-decl in Emacs).
154 
155 private:
156 
160  string doYRange(string);
161 
165  Energy thePtMin;
166 
170  Energy thePtMax;
171 
175  vector<pair<double,double> > theYRanges;
176 
181 
182 private:
183 
188  IdentifiedParticleCut & operator=(const IdentifiedParticleCut &) = delete;
189 
190 };
191 
192 }
193 
194 #endif /* Herwig_IdentifiedParticleCut_H */
transient_pointer tptr
const vector< pair< double, double > > & yRanges() const
Return the rapidity ranges.
Ptr< MatcherBase >::ptr theMatcher
A matcher for particles to cut on.
ThePEG::Ptr< InterfacedBase >::pointer IBPtr
Energy ptMin() const
Return the minimum pt.
vector< pair< double, double > > theYRanges
The rapidity ranges.
virtual Energy minKT(tcPDPtr) const
Return the minimum allowed value of the transverse momentum of an outgoing parton.
IdentifiedParticleCut implements cuts on single momenta.
ThePEG::Ptr< ParticleData >::transient_const_pointer tcPDPtr
Energy ptMax() const
Return the maximum pt.
Ptr< MatcherBase >::tptr matcher() const
Return the matcher for particles to cut on.
-*- C++ -*-
virtual double minEta(tcPDPtr) const
Return the minimum allowed pseudo-rapidity of an outgoing parton of the given type.
virtual double maxEta(tcPDPtr) const
Return the maximum allowed pseudo-rapidity of an outgoing parton of the given type.
constexpr ZeroUnit ZERO