thepeg is hosted by Hepforge, IPPP Durham
ThePEG  2.2.1
Cuts.h
1 // -*- C++ -*-
2 //
3 // Cuts.h is a part of ThePEG - Toolkit for HEP Event Generation
4 // Copyright (C) 1999-2019 Leif Lonnblad
5 //
6 // ThePEG 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 THEPEG_Cuts_H
10 #define THEPEG_Cuts_H
11 //
12 // This is the declaration of the Cuts class.
13 //
14 
15 #include "ThePEG/Interface/Interfaced.h"
16 #include "Cuts.fh"
17 #include "OneCutBase.h"
18 #include "TwoCutBase.h"
19 #include "MultiCutBase.h"
20 #include "JetFinder.h"
21 #include "FuzzyTheta.h"
22 
23 namespace ThePEG {
24 
53 class Cuts: public Interfaced {
54 
55 public:
56 
60  typedef vector<OneCutPtr> OneCutVector;
61 
65  typedef vector<TwoCutPtr> TwoCutVector;
66 
70  typedef vector<MultiCutPtr> MultiCutVector;
71 
72 public:
73 
79  Cuts(Energy MhatMin=2*GeV);
80 
84  virtual ~Cuts();
86 
87 public:
88 
99  virtual void initialize(Energy2 smax, double Y);
100 
107  virtual void initEvent();
108 
124  virtual bool
125  initSubProcess(Energy2 shat, double yhat, bool mirror = false) const;
127 
137  virtual bool passCuts(const tcPDVector & ptype, const vector<LorentzMomentum> & p,
138  tcPDPtr t1 = tcPDPtr(), tcPDPtr t2 = tcPDPtr()) const;
139 
147  bool passCuts(const tcPVector & p,
148  tcPDPtr t1 = tcPDPtr(), tcPDPtr t2 = tcPDPtr()) const;
149 
155  bool passCuts(const SubProcess & sub) const;
156 
161  bool passCuts(const Collision & coll) const;
163 
173  Energy2 minSij(tcPDPtr pi, tcPDPtr pj) const;
174 
183  Energy2 minTij(tcPDPtr pi, tcPDPtr po) const;
184 
192  double minDeltaR(tcPDPtr pi, tcPDPtr pj) const;
193 
206  Energy minKTClus(tcPDPtr pi, tcPDPtr pj) const;
207 
215  double minDurham(tcPDPtr pi, tcPDPtr pj) const;
216 
223  Energy minKT(tcPDPtr p) const;
224 
231  double minEta(tcPDPtr p) const;
232 
239  double maxEta(tcPDPtr p) const;
240 
247  double minRapidityMax(tcPDPtr p) const;
248 
255  double maxRapidityMin(tcPDPtr p) const;
256 
262  double minYStar(tcPDPtr p) const;
263 
269  double maxYStar(tcPDPtr p) const;
270 
278  Energy2 minS(const tcPDVector & pv) const;
279 
287  Energy2 maxS(const tcPDVector & pv) const;
289 
296  template <typename T>
297  vector<typename Ptr<T>::transient_const_pointer>
298  oneCutObjects() const;
299 
304  template <typename T>
305  vector<typename Ptr<T>::transient_const_pointer>
306  twoCutObjects() const;
307 
312  template <typename T>
313  vector<typename Ptr<T>::transient_const_pointer>
314  multiCutObjects() const;
315 
320  const OneCutVector& oneCuts() const { return theOneCuts; }
321 
326  const TwoCutVector& twoCuts() const { return theTwoCuts; }
327 
332  const MultiCutVector& multiCuts() const { return theMultiCuts; }
333 
338 
342  void add(tOneCutPtr c) { theOneCuts.push_back(c); }
343 
347  void add(tTwoCutPtr c) { theTwoCuts.push_back(c); }
348 
352  void add(tMultiCutPtr c) { theMultiCuts.push_back(c); }
354 
355 public:
356 
363  Energy2 SMax() const { return theSMax; }
364 
365 
370  double Y() const { return theY; }
371 
376  Energy2 currentSHat() const { return theCurrentSHat; }
377 
383  double currentYHat() const { return theCurrentYHat; }
384 
386 
392  Energy2 sHatMin() const { return max(sqr(theMHatMin), theX1Min*theX2Min*SMax()); }
393 
397  Energy2 sHatMax() const { return min(sqr(theMHatMax), theX1Max*theX2Max*SMax()); }
398 
402  bool sHat(Energy2 sh) const {
403  return sh > sHatMin() && sh <= sHatMax()*(1.0 + 1000.0*Constants::epsilon);
404  }
405 
409  Energy mHatMin() const { return max(theMHatMin, sqrt(theX1Min*theX2Min*SMax())); }
410 
414  Energy mHatMax() const { return min(theMHatMax, sqrt(theX1Max*theX2Max*SMax())); }
415 
420  double yHatMin() const;
421 
426  double yHatMax() const;
427 
431  bool yHat(double y) const;
432 
437  double x1Min() const;
438 
443  double x1Max() const;
444 
448  bool x1(double x) const;
449 
454  double x2Min() const;
455 
460  double x2Max() const;
461 
465  bool x2(double x) const;
466 
471  Energy2 scaleMin() const { return theScaleMin; }
472 
477  Energy2 scaleMax() const { return theScaleMax; }
478 
482  bool scale(Energy2 Q2) const { return Q2 > scaleMin() && Q2 < scaleMax(); }
483 
488  bool subMirror() const { return theSubMirror; }
489 
493  double cutWeight() const { return theCutWeight; }
494 
499  void lastCutWeight(double w) const { theLastCutWeight = w; }
500 
505 
509  template<class CutType, class Value>
510  bool isInside(const Value& v, const Value& lower, const Value& upper, double& weight) const {
511  if ( !fuzzy() ) {
512  if ( v >= lower && v <= upper )
513  return true;
514  weight = 0.0;
515  return false;
516  }
517  return fuzzy()->isInside<CutType>(v,lower,upper,weight);
518  }
519 
523  template<class CutType, class Value>
524  bool isLessThan(const Value& v, const Value& upper, double& weight) const {
525  if ( !fuzzy() ) {
526  if ( v <= upper )
527  return true;
528  weight = 0.0;
529  return false;
530  }
531  return fuzzy()->isLessThan<CutType>(v,upper,weight);
532  }
533 
537  template<class CutType, class Value>
538  bool isLargerThan(const Value& v, const Value& lower, double& weight) const {
539  if ( !fuzzy() ) {
540  if ( v >= lower )
541  return true;
542  weight = 0.0;
543  return false;
544  }
545  return fuzzy()->isLargerThan<CutType>(v,lower,weight);
546  }
548 
549 public:
550 
554  virtual void describe() const;
555 
556 protected:
557 
564  virtual void doinitrun();
566 
567 public:
568 
575  void persistentOutput(PersistentOStream & os) const;
576 
582  void persistentInput(PersistentIStream & is, int version);
584 
591  static void Init();
592 
593 protected:
594 
601  virtual IBPtr clone() const;
602 
607  virtual IBPtr fullclone() const;
609 
610 private:
611 
615  Energy maxMHatMin() const;
616 
620  Energy minMHatMax() const;
621 
625  double maxYHatMin() const;
626 
630  double minYHatMax() const;
631 
635  double maxX1Min() const;
636 
640  double minX1Max() const;
641 
645  double maxX2Min() const;
646 
650  double minX2Max() const;
651 
655  Energy2 maxScaleMin() const;
656 
660  Energy2 minScaleMax() const;
661 
662 private:
663 
669 
674  double theY;
675 
681 
687  mutable double theCurrentYHat;
688 
693 
698 
703  double theYHatMin;
704 
709  double theYHatMax;
710 
715  double theX1Min;
716 
721  double theX1Max;
722 
727  double theX2Min;
728 
733  double theX2Max;
734 
740 
746 
751  OneCutVector theOneCuts;
752 
757  TwoCutVector theTwoCuts;
758 
763  MultiCutVector theMultiCuts;
764 
770 
775  mutable bool theSubMirror;
776 
780  mutable double theCutWeight;
781 
786  mutable double theLastCutWeight;
787 
792 
793 private:
794 
800 
805  Cuts & operator=(const Cuts &) = delete;
806 
807 };
808 
809 }
810 
811 #include "ThePEG/Utilities/ClassTraits.h"
812 
813 namespace ThePEG {
814 
819 template <>
820 struct BaseClassTrait<Cuts,1> {
822  typedef Interfaced NthBase;
823 };
824 
827 template <>
828 struct ClassTraits<Cuts>
829  : public ClassTraitsBase<Cuts> {
831  static string className() { return "ThePEG::Cuts"; }
832 };
833 
836 }
837 
838 #endif /* THEPEG_Cuts_H */
double x2Max() const
The maximum value of the negative light-cone fraction of the hard sub-process.
double currentYHat() const
The total rapidity of hard sub-process (wrt.
Definition: Cuts.h:383
Energy maxMHatMin() const
Helper function used by the interface.
virtual IBPtr fullclone() const
Make a clone of this object, possibly modifying the cloned object to make it sane.
PersistentIStream is used to read persistent objects from a stream where they were previously written...
void persistentInput(PersistentIStream &is, int version)
Function used to read in object persistently.
Energy minKTClus(tcPDPtr pi, tcPDPtr pj) const
Return the minimum allowed value of the longitudinally invariant -algorithms distance measure...
double x2Min() const
The minimum value of the negative light-cone fraction of the hard sub-process.
double theY
The total rapidity of the colliding particles corresponding to the maximum invariant mass squared...
Definition: Cuts.h:674
Ptr< FuzzyTheta >::tcptr fuzzy() const
Return the fuzziness object.
Definition: Cuts.h:504
double theX2Min
The minimum value of the negative light-cone fraction of the hard sub-process.
Definition: Cuts.h:727
A concreate implementation of ClassDescriptionBase describing a concrete class with persistent data...
double theX1Max
The maximum value of the positive light-cone fraction of the hard sub-process.
Definition: Cuts.h:721
double theX1Min
The minimum value of the positive light-cone fraction of the hard sub-process.
Definition: Cuts.h:715
bool isInside(const Value &v, const Value &lower, const Value &upper, double &weight) const
Check for value inside the given bounds and update the weight.
Definition: Cuts.h:510
bool isLessThan(const Value &v, const Value &upper, double &weight) const
Check for value inside the given bounds and update the weight.
Definition: Cuts.h:524
Ptr< JetFinder >::ptr theJetFinder
An optional jet finder used to define cuts on the level of reconstructed jets.
Definition: Cuts.h:769
PersistentOStream is used to write objects persistently to a stream from which they can be read in ag...
double theYHatMin
The minimum value of the rapidity of the hard sub-process (wrt.
Definition: Cuts.h:703
vector< tcPDPtr > tcPDVector
A vector of transient pointers to const ParticleData objects.
Definition: Containers.h:42
void add(tOneCutPtr c)
Add a OneCutBase object.
Definition: Cuts.h:342
TransientConstRCPtr is a simple wrapper around a bare const pointer which can be assigned to and from...
Definition: RCPtr.h:696
TransientRCPtr is a simple wrapper around a bare pointer which can be assigned to and from an RCPtr a...
Definition: RCPtr.h:519
double minDurham(tcPDPtr pi, tcPDPtr pj) const
Return the minimum allowed value of the Durham -algorithms distance measure.
constexpr auto sqr(const T &x) -> decltype(x *x)
The square function should really have been included in the standard C++ library. ...
Definition: ThePEG.h:117
double maxYStar(tcPDPtr p) const
Return the minimum allowed rapidity of an outgoing parton of the given type in the center-of-mass sys...
Energy2 currentSHat() const
The invariant mass squared of the hard sub-process of the event being considered. ...
Definition: Cuts.h:376
double theCurrentYHat
The total rapidity of hard sub-process (wrt.
Definition: Cuts.h:687
bool yHat(double y) const
Check if the given is within the cuts.
Energy2 minS(const tcPDVector &pv) const
Return the minimum allowed value of the squared invariant mass of a set of outgoing partons of the gi...
vector< typename Ptr< T >::transient_const_pointer > twoCutObjects() const
Return a vector of pointers to objects of the given class (with base class TwoCutBase).
double minX2Max() const
Helper function used by the interface.
This is the main namespace within which all identifiers in ThePEG are declared.
Definition: FactoryBase.h:28
virtual void doinitrun()
Initialize this object.
bool x1(double x) const
Check if the given is within the cuts.
Energy2 theScaleMax
The maximum allowed value of the scale to be used in PDF&#39;s and coupling constants.
Definition: Cuts.h:745
Energy2 SMax() const
The maximum allowed total invariant mass squared allowed for events to be considered.
Definition: Cuts.h:363
double minYHatMax() const
Helper function used by the interface.
Energy theMHatMin
The minimum allowed value of .
Definition: Cuts.h:692
bool theSubMirror
Set to true if a matrix element is should be using this cut and is mirrored along the z-axis ...
Definition: Cuts.h:775
constexpr double epsilon
The smallest non-zero double.
Definition: Constants.h:63
Energy mHatMin() const
The minimum allowed value of .
Definition: Cuts.h:409
Energy2 theSMax
The maximum allowed total invariant mass squared allowed for events to be considered.
Definition: Cuts.h:668
double maxEta(tcPDPtr p) const
Return the maximum allowed pseudo-rapidity of an outgoing parton of the given type.
A SubProcess object represents a hard sub-process in a collision.
Definition: SubProcess.h:33
Ptr< FuzzyTheta >::ptr theFuzzyTheta
The fuzziness object.
Definition: Cuts.h:791
static ClassDescription< Cuts > initCuts
The static object used to initialize the description of this class.
Definition: Cuts.h:799
Energy2 scaleMax() const
The maximum allowed value of the scale to be used in PDF&#39;s and coupling constants.
Definition: Cuts.h:477
double yHatMax() const
The maximum value of the rapidity of the hard sub-process (wrt.
double Y() const
The total rapidity of the colliding particles corresponding to the maximum invariant mass squared...
Definition: Cuts.h:370
Cuts & operator=(const Cuts &)=delete
The assignment operator is private and must never be called.
double minX1Max() const
Helper function used by the interface.
double minRapidityMax(tcPDPtr p) const
Return the minimum allowed rapidity of an outgoing parton of the given type.
vector< OneCutPtr > OneCutVector
A vector of OneCutBase pointers.
Definition: Cuts.h:60
double x1Max() const
The maximum value of the positive light-cone fraction of the hard sub-process.
vector< tcPPtr > tcPVector
A vector of transient pointers to const Particle objects.
Definition: Containers.h:85
Energy minKT(tcPDPtr p) const
Return the minimum allowed value of the transverse momentum of an outgoing parton.
Ptr< JetFinder >::tptr jetFinder() const
Return the jet finder.
Definition: Cuts.h:337
double minDeltaR(tcPDPtr pi, tcPDPtr pj) const
Return the minimum allowed value of of two outgoing partons of type pi and pj.
bool subMirror() const
Set true if a matrix element is should be using this cut and is mirrored along the z-axis ...
Definition: Cuts.h:488
double x1Min() const
The minimum value of the positive light-cone fraction of the hard sub-process.
Cuts(Energy MhatMin=2 *GeV)
The default constructor.
MultiCutVector theMultiCuts
The objects defining cuts on sets of outgoing particles from the hard sub-process.
Definition: Cuts.h:763
ThePEG::Ptr< ParticleData >::transient_const_pointer tcPDPtr
Alias for a transient pointer to a const ParticleData .
Definition: Pointers.h:64
void add(tMultiCutPtr c)
Add a MultiCutBase object.
Definition: Cuts.h:352
Energy2 minSij(tcPDPtr pi, tcPDPtr pj) const
Return the minimum allowed squared invariant mass of two outgoing partons of type pi and pj...
Energy2 maxS(const tcPDVector &pv) const
Return the maximum allowed value of the squared invariant mass of a set of outgoing partons of the gi...
void persistentOutput(PersistentOStream &os) const
Function used to write out object persistently.
bool scale(Energy2 Q2) const
Check if the given scale is within the cuts.
Definition: Cuts.h:482
virtual bool initSubProcess(Energy2 shat, double yhat, bool mirror=false) const
Set information about the invariant mass squared, shat, and rapidity, yhat, of the hard sub-process...
RCPtr is a reference counted (smart) pointer.
Definition: RCPtr.h:60
const TwoCutVector & twoCuts() const
Return the objects defining cuts on pairs of particles in the hard sub-process.
Definition: Cuts.h:326
TwoCutVector theTwoCuts
The objects defining cuts on pairs of particles in the hard sub-process.
Definition: Cuts.h:757
static void Init()
The standard Init function used to initialize the interfaces.
This is the decalaration of the Collision class.
Definition: Collision.h:34
The Interfaced class is derived from the InterfacedBase class adding a couple of things particular to...
Definition: Interfaced.h:38
Energy2 scaleMin() const
The minimum allowed value of the scale to be used in PDF&#39;s and coupling constants.
Definition: Cuts.h:471
Energy2 sHatMin() const
The minimum allowed value of .
Definition: Cuts.h:392
Energy theMHatMax
The maximum allowed value of .
Definition: Cuts.h:697
vector< typename Ptr< T >::transient_const_pointer > oneCutObjects() const
Return a vector of pointers to objects of the given class (with base class OneCutBase).
vector< TwoCutPtr > TwoCutVector
A vector of TwoCutBase pointers.
Definition: Cuts.h:65
Energy2 sHatMax() const
The maximum allowed value of .
Definition: Cuts.h:397
vector< typename Ptr< T >::transient_const_pointer > multiCutObjects() const
Return a vector of pointers to objects of the given class (with base class MultiCutBase).
Energy2 maxScaleMin() const
Helper function used by the interface.
double theYHatMax
The maximum value of the rapidity of the hard sub-process (wrt.
Definition: Cuts.h:709
double minYStar(tcPDPtr p) const
Return the minimum allowed rapidity of an outgoing parton of the given type in the center-of-mass sys...
double cutWeight() const
Return the overall cut weight.
Definition: Cuts.h:493
double maxRapidityMin(tcPDPtr p) const
Return the maximum allowed rapidity of an outgoing parton of the given type.
double theX2Max
The maximum value of the negative light-cone fraction of the hard sub-process.
Definition: Cuts.h:733
The default concrete implementation of ClassTraitsBase.
Definition: ClassTraits.h:134
bool isLargerThan(const Value &v, const Value &lower, double &weight) const
Check for value inside the given bounds and update the weight.
Definition: Cuts.h:538
double maxX2Min() const
Helper function used by the interface.
double theCutWeight
The overall cut weight.
Definition: Cuts.h:780
Energy minMHatMax() const
Helper function used by the interface.
Energy2 theCurrentSHat
The invariant mass squared of the hard sub-process of the event being considered. ...
Definition: Cuts.h:680
double theLastCutWeight
The cut weight as appropriate from the call to the last n-cut object.
Definition: Cuts.h:786
virtual void describe() const
Describe the currently active cuts in the log file.
OneCutVector theOneCuts
The objects defining cuts on single outgoing partons from the hard sub-process.
Definition: Cuts.h:751
virtual ~Cuts()
The destructor.
void lastCutWeight(double w) const
Set the cut weight as appropriate from the call to the last n-cut object.
Definition: Cuts.h:499
virtual bool passCuts(const tcPDVector &ptype, const vector< LorentzMomentum > &p, tcPDPtr t1=tcPDPtr(), tcPDPtr t2=tcPDPtr()) const
Check if the outgoing particles, with the given types and momenta, from a sub-process passes the cuts...
Energy2 theScaleMin
The minimum allowed value of the scale to be used in PDF&#39;s and coupling constants.
Definition: Cuts.h:739
virtual void initialize(Energy2 smax, double Y)
Initialize this object specifying the maximum total invariant mass squared, smax, and the total rapid...
Energy2 minScaleMax() const
Helper function used by the interface.
double yHatMin() const
The minimum value of the rapidity of the hard sub-process (wrt.
const MultiCutVector & multiCuts() const
Return the objects defining cuts on sets of outgoing particles from the hard sub-process.
Definition: Cuts.h:332
bool sHat(Energy2 sh) const
Check if the given is within the cuts.
Definition: Cuts.h:402
BaseClassTraits describes the base classes of the templated class.
Definition: ClassTraits.h:156
double minEta(tcPDPtr p) const
Return the minimum allowed pseudo-rapidity of an outgoing parton of the given type.
double maxYHatMin() const
Helper function used by the interface.
vector< MultiCutPtr > MultiCutVector
A vector of MultiCutBase pointers.
Definition: Cuts.h:70
virtual IBPtr clone() const
Make a simple clone of this object.
Energy mHatMax() const
The maximum allowed value of .
Definition: Cuts.h:414
double maxX1Min() const
Helper function used by the interface.
const OneCutVector & oneCuts() const
Return the objects defining cuts on single outgoing partons from the hard sub-process.
Definition: Cuts.h:320
Cuts is a class for implementing kinematical cuts in ThePEG.
Definition: Cuts.h:53
void add(tTwoCutPtr c)
Add a TwoCutBase object.
Definition: Cuts.h:347
The templated ClassTraitsBase class defines a set of default information about classes used by ThePEG...
Definition: ClassTraits.h:52
bool x2(double x) const
Check if the given is within the cuts.
virtual void initEvent()
Initialize this object for a new event.
Energy2 minTij(tcPDPtr pi, tcPDPtr po) const
Return the minimum allowed value of the negative of the squared invariant mass of an incoming parton ...