thepeg is hosted by Hepforge, IPPP Durham
ThePEG 2.3.0
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
23namespace ThePEG {
24
53class Cuts: public Interfaced {
54
55public:
56
60 typedef vector<OneCutPtr> OneCutVector;
61
65 typedef vector<TwoCutPtr> TwoCutVector;
66
70 typedef vector<MultiCutPtr> MultiCutVector;
71
72public:
73
79 Cuts(Energy MhatMin=2*GeV);
81
82public:
83
94 virtual void initialize(Energy2 smax, double Y);
95
102 virtual void initEvent();
103
119 virtual bool
120 initSubProcess(Energy2 shat, double yhat, bool mirror = false) const;
122
132 virtual bool passCuts(const tcPDVector & ptype, const vector<LorentzMomentum> & p,
133 tcPDPtr t1 = tcPDPtr(), tcPDPtr t2 = tcPDPtr()) const;
134
142 bool passCuts(const tcPVector & p,
143 tcPDPtr t1 = tcPDPtr(), tcPDPtr t2 = tcPDPtr()) const;
144
150 bool passCuts(const SubProcess & sub) const;
151
156 bool passCuts(const Collision & coll) const;
158
169
179
187 double minDeltaR(tcPDPtr pi, tcPDPtr pj) const;
188
202
210 double minDurham(tcPDPtr pi, tcPDPtr pj) const;
211
219
226 double minEta(tcPDPtr p) const;
227
234 double maxEta(tcPDPtr p) const;
235
242 double minRapidityMax(tcPDPtr p) const;
243
250 double maxRapidityMin(tcPDPtr p) const;
251
257 double minYStar(tcPDPtr p) const;
258
264 double maxYStar(tcPDPtr p) const;
265
273 Energy2 minS(const tcPDVector & pv) const;
274
282 Energy2 maxS(const tcPDVector & pv) const;
284
291 template <typename T>
292 vector<typename Ptr<T>::transient_const_pointer>
294
299 template <typename T>
300 vector<typename Ptr<T>::transient_const_pointer>
302
307 template <typename T>
308 vector<typename Ptr<T>::transient_const_pointer>
310
315 const OneCutVector& oneCuts() const { return theOneCuts; }
316
321 const TwoCutVector& twoCuts() const { return theTwoCuts; }
322
327 const MultiCutVector& multiCuts() const { return theMultiCuts; }
328
333
337 void add(tOneCutPtr c) { theOneCuts.push_back(c); }
338
342 void add(tTwoCutPtr c) { theTwoCuts.push_back(c); }
343
347 void add(tMultiCutPtr c) { theMultiCuts.push_back(c); }
349
350public:
351
358 Energy2 SMax() const { return theSMax; }
359
360
365 double Y() const { return theY; }
366
372
378 double currentYHat() const { return theCurrentYHat; }
379
381
387 Energy2 sHatMin() const { return max(sqr(theMHatMin), theX1Min*theX2Min*SMax()); }
388
392 Energy2 sHatMax() const { return min(sqr(theMHatMax), theX1Max*theX2Max*SMax()); }
393
397 bool sHat(Energy2 sh) const {
398 return sh > sHatMin() && sh <= sHatMax()*(1.0 + 1000.0*Constants::epsilon);
399 }
400
404 Energy mHatMin() const { return max(theMHatMin, sqrt(theX1Min*theX2Min*SMax())); }
405
409 Energy mHatMax() const { return min(theMHatMax, sqrt(theX1Max*theX2Max*SMax())); }
410
415 double yHatMin() const;
416
421 double yHatMax() const;
422
426 bool yHat(double y) const;
427
432 double x1Min() const;
433
438 double x1Max() const;
439
443 bool x1(double x) const;
444
449 double x2Min() const;
450
455 double x2Max() const;
456
460 bool x2(double x) const;
461
466 Energy2 scaleMin() const { return theScaleMin; }
467
472 Energy2 scaleMax() const { return theScaleMax; }
473
477 bool scale(Energy2 Q2) const { return Q2 > scaleMin() && Q2 < scaleMax(); }
478
483 bool subMirror() const { return theSubMirror; }
484
488 double cutWeight() const { return theCutWeight; }
489
494 void lastCutWeight(double w) const { theLastCutWeight = w; }
495
500
504 template<class CutType, class Value>
505 bool isInside(const Value& v, const Value& lower, const Value& upper, double& weight) const {
506 if ( !fuzzy() ) {
507 if ( v >= lower && v <= upper )
508 return true;
509 weight = 0.0;
510 return false;
511 }
512 return fuzzy()->isInside<CutType>(v,lower,upper,weight);
513 }
514
518 template<class CutType, class Value>
519 bool isLessThan(const Value& v, const Value& upper, double& weight) const {
520 if ( !fuzzy() ) {
521 if ( v <= upper )
522 return true;
523 weight = 0.0;
524 return false;
525 }
526 return fuzzy()->isLessThan<CutType>(v,upper,weight);
527 }
528
532 template<class CutType, class Value>
533 bool isLargerThan(const Value& v, const Value& lower, double& weight) const {
534 if ( !fuzzy() ) {
535 if ( v >= lower )
536 return true;
537 weight = 0.0;
538 return false;
539 }
540 return fuzzy()->isLargerThan<CutType>(v,lower,weight);
541 }
543
544public:
545
549 virtual void describe() const;
550
551protected:
552
559 virtual void doinitrun();
561
562public:
563
571
577 void persistentInput(PersistentIStream & is, int version);
579
586 static void Init();
587
588protected:
589
596 virtual IBPtr clone() const;
597
602 virtual IBPtr fullclone() const;
604
605private:
606
611
616
620 double maxYHatMin() const;
621
625 double minYHatMax() const;
626
630 double maxX1Min() const;
631
635 double minX1Max() const;
636
640 double maxX2Min() const;
641
645 double minX2Max() const;
646
651
656
657private:
658
664
669 double theY;
670
676
682 mutable double theCurrentYHat;
683
688
693
699
705
710 double theX1Min;
711
716 double theX1Max;
717
722 double theX2Min;
723
728 double theX2Max;
729
735
741
747
753
759
765
770 mutable bool theSubMirror;
771
775 mutable double theCutWeight;
776
781 mutable double theLastCutWeight;
782
787
788private:
789
795
800 Cuts & operator=(const Cuts &) = delete;
801
802};
803
804}
805
806#include "ThePEG/Utilities/ClassTraits.h"
807
808namespace ThePEG {
809
814template <>
815struct BaseClassTrait<Cuts,1> {
817 typedef Interfaced NthBase;
818};
819
822template <>
823struct ClassTraits<Cuts>
824 : public ClassTraitsBase<Cuts> {
826 static string className() { return "ThePEG::Cuts"; }
827};
828
831}
832
833#endif /* THEPEG_Cuts_H */
A concreate implementation of ClassDescriptionBase describing a concrete class with persistent data.
This is the decalaration of the Collision class.
Definition: Collision.h:34
Cuts is a class for implementing kinematical cuts in ThePEG.
Definition: Cuts.h:53
const TwoCutVector & twoCuts() const
Return the objects defining cuts on pairs of particles in the hard sub-process.
Definition: Cuts.h:321
double minDeltaR(tcPDPtr pi, tcPDPtr pj) const
Return the minimum allowed value of of two outgoing partons of type pi and pj.
bool passCuts(const SubProcess &sub) const
Check if the incoming and outgoing particles in the given sub-process passes the cuts.
Energy mHatMax() const
The maximum allowed value of .
Definition: Cuts.h:409
virtual void initialize(Energy2 smax, double Y)
Initialize this object specifying the maximum total invariant mass squared, smax, and the total rapid...
const OneCutVector & oneCuts() const
Return the objects defining cuts on single outgoing partons from the hard sub-process.
Definition: Cuts.h:315
bool passCuts(const Collision &coll) const
Check if the given collision passes the cuts.
Ptr< FuzzyTheta >::ptr theFuzzyTheta
The fuzziness object.
Definition: Cuts.h:786
double minX2Max() const
Helper function used by the interface.
double minEta(tcPDPtr p) const
Return the minimum allowed pseudo-rapidity of an outgoing parton of the given type.
bool x1(double x) const
Check if the given is within the cuts.
double cutWeight() const
Return the overall cut weight.
Definition: Cuts.h:488
Energy theMHatMin
The minimum allowed value of .
Definition: Cuts.h:687
Energy2 sHatMin() const
The minimum allowed value of .
Definition: Cuts.h:387
void add(tOneCutPtr c)
Add a OneCutBase object.
Definition: Cuts.h:337
Energy minMHatMax() const
Helper function used by the interface.
double minX1Max() const
Helper function used by the interface.
void lastCutWeight(double w) const
Set the cut weight as appropriate from the call to the last n-cut object.
Definition: Cuts.h:494
double theYHatMin
The minimum value of the rapidity of the hard sub-process (wrt.
Definition: Cuts.h:698
TwoCutVector theTwoCuts
The objects defining cuts on pairs of particles in the hard sub-process.
Definition: Cuts.h:752
vector< typename Ptr< T >::transient_const_pointer > multiCutObjects() const
Return a vector of pointers to objects of the given class (with base class MultiCutBase).
double x2Max() const
The maximum value of the negative light-cone fraction of the hard sub-process.
virtual IBPtr fullclone() const
Make a clone of this object, possibly modifying the cloned object to make it sane.
double theX2Min
The minimum value of the negative light-cone fraction of the hard sub-process.
Definition: Cuts.h:722
virtual void doinitrun()
Initialize this object.
double maxEta(tcPDPtr p) const
Return the maximum allowed pseudo-rapidity of an outgoing parton of the given type.
Energy2 minScaleMax() const
Helper function used by the interface.
void add(tTwoCutPtr c)
Add a TwoCutBase object.
Definition: Cuts.h:342
double x2Min() const
The minimum value of the negative light-cone fraction of the hard sub-process.
double theYHatMax
The maximum value of the rapidity of the hard sub-process (wrt.
Definition: Cuts.h:704
Cuts & operator=(const Cuts &)=delete
The assignment operator is private and must never be called.
Energy2 minTij(tcPDPtr pi, tcPDPtr po) const
Return the minimum allowed value of the negative of the squared invariant mass of an incoming parton ...
double maxX2Min() const
Helper function used by the interface.
double theCurrentYHat
The total rapidity of hard sub-process (wrt.
Definition: Cuts.h:682
const MultiCutVector & multiCuts() const
Return the objects defining cuts on sets of outgoing particles from the hard sub-process.
Definition: Cuts.h:327
Ptr< JetFinder >::ptr theJetFinder
An optional jet finder used to define cuts on the level of reconstructed jets.
Definition: Cuts.h:764
Ptr< FuzzyTheta >::tcptr fuzzy() const
Return the fuzziness object.
Definition: Cuts.h:499
double currentYHat() const
The total rapidity of hard sub-process (wrt.
Definition: Cuts.h:378
double maxYHatMin() const
Helper function used by the interface.
void persistentOutput(PersistentOStream &os) const
Function used to write out object persistently.
Energy2 currentSHat() const
The invariant mass squared of the hard sub-process of the event being considered.
Definition: Cuts.h:371
OneCutVector theOneCuts
The objects defining cuts on single outgoing partons from the hard sub-process.
Definition: Cuts.h:746
Energy2 theSMax
The maximum allowed total invariant mass squared allowed for events to be considered.
Definition: Cuts.h:663
double maxYStar(tcPDPtr p) const
Return the minimum allowed rapidity of an outgoing parton of the given type in the center-of-mass sys...
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:519
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.
double maxRapidityMin(tcPDPtr p) const
Return the maximum allowed rapidity of an outgoing parton of the given type.
MultiCutVector theMultiCuts
The objects defining cuts on sets of outgoing particles from the hard sub-process.
Definition: Cuts.h:758
vector< TwoCutPtr > TwoCutVector
A vector of TwoCutBase pointers.
Definition: Cuts.h:65
double minYHatMax() const
Helper function used by the interface.
double yHatMax() const
The maximum value of the rapidity of the hard sub-process (wrt.
Energy2 scaleMin() const
The minimum allowed value of the scale to be used in PDF's and coupling constants.
Definition: Cuts.h:466
double x1Max() const
The maximum value of the positive light-cone fraction of the hard sub-process.
double x1Min() const
The minimum value of the positive light-cone fraction of the hard sub-process.
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...
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...
double theX2Max
The maximum value of the negative light-cone fraction of the hard sub-process.
Definition: Cuts.h:728
Energy maxMHatMin() const
Helper function used by the interface.
Cuts(Energy MhatMin=2 *GeV)
The default constructor.
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:505
Energy minKTClus(tcPDPtr pi, tcPDPtr pj) const
Return the minimum allowed value of the longitudinally invariant -algorithms distance measure.
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:533
Energy2 sHatMax() const
The maximum allowed value of .
Definition: Cuts.h:392
double maxX1Min() const
Helper function used by the interface.
bool scale(Energy2 Q2) const
Check if the given scale is within the cuts.
Definition: Cuts.h:477
vector< OneCutPtr > OneCutVector
A vector of OneCutBase pointers.
Definition: Cuts.h:60
double Y() const
The total rapidity of the colliding particles corresponding to the maximum invariant mass squared,...
Definition: Cuts.h:365
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 minRapidityMax(tcPDPtr p) const
Return the minimum allowed rapidity of an outgoing parton of the given type.
void add(tMultiCutPtr c)
Add a MultiCutBase object.
Definition: Cuts.h:347
bool sHat(Energy2 sh) const
Check if the given is within the cuts.
Definition: Cuts.h:397
virtual void describe() const
Describe the currently active cuts in the log file.
Energy mHatMin() const
The minimum allowed value of .
Definition: Cuts.h:404
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...
virtual IBPtr clone() const
Make a simple clone of this object.
Energy2 theScaleMin
The minimum allowed value of the scale to be used in PDF's and coupling constants.
Definition: Cuts.h:734
double theX1Min
The minimum value of the positive light-cone fraction of the hard sub-process.
Definition: Cuts.h:710
Energy2 theScaleMax
The maximum allowed value of the scale to be used in PDF's and coupling constants.
Definition: Cuts.h:740
bool passCuts(const tcPVector &p, tcPDPtr t1=tcPDPtr(), tcPDPtr t2=tcPDPtr()) const
Check if the outgoing particles from a sub-process passes the cuts.
void persistentInput(PersistentIStream &is, int version)
Function used to read in object persistently.
Energy minKT(tcPDPtr p) const
Return the minimum allowed value of the transverse momentum of an outgoing parton.
double theCutWeight
The overall cut weight.
Definition: Cuts.h:775
Energy2 SMax() const
The maximum allowed total invariant mass squared allowed for events to be considered.
Definition: Cuts.h:358
double theX1Max
The maximum value of the positive light-cone fraction of the hard sub-process.
Definition: Cuts.h:716
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:770
virtual void initEvent()
Initialize this object for a new event.
Energy2 maxScaleMin() const
Helper function used by the interface.
double minDurham(tcPDPtr pi, tcPDPtr pj) const
Return the minimum allowed value of the Durham -algorithms distance measure.
Ptr< JetFinder >::tptr jetFinder() const
Return the jet finder.
Definition: Cuts.h:332
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:483
static ClassDescription< Cuts > initCuts
The static object used to initialize the description of this class.
Definition: Cuts.h:794
double yHatMin() const
The minimum value of the rapidity of the hard sub-process (wrt.
Energy2 theCurrentSHat
The invariant mass squared of the hard sub-process of the event being considered.
Definition: Cuts.h:675
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 theY
The total rapidity of the colliding particles corresponding to the maximum invariant mass squared,...
Definition: Cuts.h:669
Energy2 minSij(tcPDPtr pi, tcPDPtr pj) const
Return the minimum allowed squared invariant mass of two outgoing partons of type pi and pj.
vector< typename Ptr< T >::transient_const_pointer > oneCutObjects() const
Return a vector of pointers to objects of the given class (with base class OneCutBase).
Energy2 scaleMax() const
The maximum allowed value of the scale to be used in PDF's and coupling constants.
Definition: Cuts.h:472
static void Init()
The standard Init function used to initialize the interfaces.
bool yHat(double y) const
Check if the given is within the cuts.
vector< MultiCutPtr > MultiCutVector
A vector of MultiCutBase pointers.
Definition: Cuts.h:70
bool x2(double x) const
Check if the given is within the cuts.
Energy theMHatMax
The maximum allowed value of .
Definition: Cuts.h:692
double theLastCutWeight
The cut weight as appropriate from the call to the last n-cut object.
Definition: Cuts.h:781
The Interfaced class is derived from the InterfacedBase class adding a couple of things particular to...
Definition: Interfaced.h:38
PersistentIStream is used to read persistent objects from a stream where they were previously written...
PersistentOStream is used to write objects persistently to a stream from which they can be read in ag...
RCPtr is a reference counted (smart) pointer.
Definition: RCPtr.h:60
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
A SubProcess object represents a hard sub-process in a collision.
Definition: SubProcess.h:33
constexpr double epsilon
The smallest non-zero double.
Definition: Constants.h:63
This is the main namespace within which all identifiers in ThePEG are declared.
Definition: FactoryBase.h:28
ThePEG::Ptr< ParticleData >::transient_const_pointer tcPDPtr
Alias for a transient pointer to a const ParticleData .
Definition: Pointers.h:64
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
vector< tcPPtr > tcPVector
A vector of transient pointers to const Particle objects.
Definition: Containers.h:85
vector< tcPDPtr > tcPDVector
A vector of transient pointers to const ParticleData objects.
Definition: Containers.h:42
int NthBase
The type of the BaseN'th base class (int means there are no further base classes).
Definition: ClassTraits.h:161
static string className()
Return the name of class T.
Definition: ClassTraits.h:66