15 #include "ThePEG/Interface/Interfaced.h" 17 #include "OneCutBase.h" 18 #include "TwoCutBase.h" 19 #include "MultiCutBase.h" 20 #include "JetFinder.h" 124 initSubProcess(Energy2 shat,
double yhat,
bool mirror =
false)
const;
295 template <
typename T>
296 vector<typename Ptr<T>::transient_const_pointer>
303 template <
typename T>
304 vector<typename Ptr<T>::transient_const_pointer>
311 template <
typename T>
312 vector<typename Ptr<T>::transient_const_pointer>
430 bool yHat(
double y)
const;
436 double x1Min()
const;
442 double x1Max()
const;
447 bool x1(
double x)
const;
453 double x2Min()
const;
459 double x2Max()
const;
464 bool x2(
double x)
const;
736 #include "ThePEG/Utilities/ClassTraits.h" 756 static string className() {
return "ThePEG::Cuts"; }
Energy maxMHatMin() const
Helper function used by the interface.
Energy2 currentSHat() const
The invariant mass squared of the hard sub-process of the event being considered. ...
Energy2 minSij(tcPDPtr pi, tcPDPtr pj) const
Return the minimum allowed squared invariant mass of two outgoing partons of type pi and pj...
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.
bool yHat(double y) const
Check if the given is within the cuts.
double currentYHat() const
The total rapidity of hard sub-process (wrt.
double minYHatMax() const
Helper function used by the interface.
double theY
The total rapidity of the colliding particles corresponding to the maximum invariant mass squared...
double maxEta(tcPDPtr p) const
Return the maximum allowed pseudo-rapidity of an outgoing parton of the given type.
double Y() const
The total rapidity of the colliding particles corresponding to the maximum invariant mass squared...
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...
Energy minKT(tcPDPtr p) const
Return the minimum allowed value of the transverse momentum of an outgoing parton.
double yHatMin() const
The minimum value of the rapidity of the hard sub-process (wrt.
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 theX2Min
The minimum value of the negative light-cone fraction of the hard sub-process.
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.
double theX1Min
The minimum value of the positive light-cone fraction of the hard sub-process.
Energy2 maxScaleMin() const
Helper function used by the interface.
Ptr< JetFinder >::ptr theJetFinder
An optional jet finder used to define cuts on the level of reconstructed jets.
Energy mHatMin() const
The minimum allowed value of .
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.
vector< tcPDPtr > tcPDVector
A vector of transient pointers to const ParticleData objects.
void add(tOneCutPtr c)
Add a OneCutBase object.
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...
TransientConstRCPtr is a simple wrapper around a bare const pointer which can be assigned to and from...
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 minYStar(tcPDPtr p) const
Return the minimum allowed rapidity of an outgoing parton of the given type in the center-of-mass sys...
TransientRCPtr is a simple wrapper around a bare pointer which can be assigned to and from an RCPtr a...
double theCurrentYHat
The total rapidity of hard sub-process (wrt.
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 x1Max() const
The maximum value of the positive light-cone fraction of the hard sub-process.
const TwoCutVector & twoCuts() const
Return the objects defining cuts on pairs of particles in the hard sub-process.
This is the main namespace within which all identifiers in ThePEG are declared.
virtual void doinitrun()
Initialize this object.
Energy minKTClus(tcPDPtr pi, tcPDPtr pj) const
Return the minimum allowed value of the longitudinally invariant -algorithms distance measure...
Energy2 sHatMin() const
The minimum allowed value of .
Energy2 theScaleMax
The maximum allowed value of the scale to be used in PDF's and coupling constants.
Energy theMHatMin
The minimum allowed value of .
bool theSubMirror
Set to true if a matrix element is should be using this cut and is mirrored along the z-axis ...
Energy2 theSMax
The maximum allowed total invariant mass squared allowed for events to be considered.
Cuts & operator=(const Cuts &)
The assignment operator is private and must never be called.
A SubProcess object represents a hard sub-process in a collision.
bool x2(double x) const
Check if the given is within the cuts.
QTY< 0, 1, 0 >::Type Energy
Energy.
static ClassDescription< Cuts > initCuts
The static object used to initialize the description of this class.
bool subMirror() const
Set true if a matrix element is should be using this cut and is mirrored along the z-axis ...
double minRapidityMax(tcPDPtr p) const
Return the minimum allowed rapidity of an outgoing parton of the given type.
virtual IBPtr clone() const
Make a simple clone of this object.
Energy2 scaleMin() const
The minimum allowed value of the scale to be used in PDF's and coupling constants.
vector< OneCutPtr > OneCutVector
A vector of OneCutBase pointers.
bool sHat(Energy2 sh) const
Check if the given is within the cuts.
virtual void describe() const
Describe the currently active cuts in the log file.
vector< tcPPtr > tcPVector
A vector of transient pointers to const Particle objects.
Energy2 scaleMax() const
The maximum allowed value of the scale to be used in PDF's and coupling constants.
vector< typename Ptr< T >::transient_const_pointer > twoCutObjects() const
Return a vector of pointers to objects of the given class (with base class TwoCutBase).
const OneCutVector & oneCuts() const
Return the objects defining cuts on single outgoing partons from the hard sub-process.
const MultiCutVector & multiCuts() const
Return the objects defining cuts on sets of outgoing particles from the hard sub-process.
double minDeltaR(tcPDPtr pi, tcPDPtr pj) const
Return the minimum allowed value of of two outgoing partons of type pi and pj.
double yHatMax() const
The maximum value of the rapidity of the hard sub-process (wrt.
Cuts(Energy MhatMin=2 *GeV)
The default constructor.
MultiCutVector theMultiCuts
The objects defining cuts on sets of outgoing particles from the hard sub-process.
Energy2 SMax() const
The maximum allowed total invariant mass squared allowed for events to be considered.
double x2Max() const
The maximum value of the negative light-cone fraction of the hard sub-process.
Energy2 minScaleMax() const
Helper function used by the interface.
Energy mHatMax() const
The maximum allowed value of .
double minX1Max() const
Helper function used by the interface.
ThePEG::Ptr< ParticleData >::transient_const_pointer tcPDPtr
Alias for a transient pointer to a const ParticleData .
void add(tMultiCutPtr c)
Add a MultiCutBase object.
double maxX1Min() const
Helper function used by the interface.
double minX2Max() const
Helper function used by the interface.
double x1Min() const
The minimum value of the positive light-cone fraction of the hard sub-process.
RCPtr is a reference counted (smart) pointer.
double minEta(tcPDPtr p) const
Return the minimum allowed pseudo-rapidity of an outgoing parton of the given type.
TwoCutVector theTwoCuts
The objects defining cuts on pairs of particles in the hard sub-process.
Energy minMHatMax() const
Helper function used by the interface.
static void Init()
The standard Init function used to initialize the interfaces.
This is the decalaration of the Collision class.
The Interfaced class is derived from the InterfacedBase class adding a couple of things particular to...
Energy theMHatMax
The maximum allowed value of .
vector< TwoCutPtr > TwoCutVector
A vector of TwoCutBase pointers.
bool x1(double x) const
Check if the given is within the cuts.
double theYHatMax
The maximum value of the rapidity of the hard sub-process (wrt.
double theX2Max
The maximum value of the negative light-cone fraction of the hard sub-process.
bool scale(Energy2 Q2) const
Check if the given scale is within the cuts.
void persistentOutput(PersistentOStream &os) const
Function used to write out object persistently.
The default concrete implementation of ClassTraitsBase.
double maxX2Min() const
Helper function used by the interface.
const double epsilon
The smallest non-zero double.
Energy2 theCurrentSHat
The invariant mass squared of the hard sub-process of the event being considered. ...
OneCutVector theOneCuts
The objects defining cuts on single outgoing partons from the hard sub-process.
virtual ~Cuts()
The destructor.
Energy2 theScaleMin
The minimum allowed value of the scale to be used in PDF's and coupling constants.
virtual void initialize(Energy2 smax, double Y)
Initialize this object specifying the maximum total invariant mass squared, smax, and the total rapid...
double minDurham(tcPDPtr pi, tcPDPtr pj) const
Return the minimum allowed value of the Durham -algorithms distance measure.
double maxYHatMin() const
Helper function used by the interface.
Ptr< JetFinder >::tptr jetFinder() const
Return the jet finder.
BaseClassTraits describes the base classes of the templated class.
double maxRapidityMin(tcPDPtr p) const
Return the maximum allowed rapidity of an outgoing parton of the given type.
vector< MultiCutPtr > MultiCutVector
A vector of MultiCutBase pointers.
virtual IBPtr fullclone() const
Make a clone of this object, possibly modifying the cloned object to make it sane.
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...
double x2Min() const
The minimum value of the negative light-cone fraction of the hard sub-process.
Energy2 sHatMax() const
The maximum allowed value of .
vector< typename Ptr< T >::transient_const_pointer > oneCutObjects() const
Return a vector of pointers to objects of the given class (with base class OneCutBase).
Cuts is a class for implementing kinematical cuts in ThePEG.
void add(tTwoCutPtr c)
Add a TwoCutBase object.
The templated ClassTraitsBase class defines a set of default information about classes used by ThePEG...
double maxYStar(tcPDPtr p) const
Return the minimum allowed rapidity of an outgoing parton of the given type in the center-of-mass sys...
virtual void initEvent()
Initialize this object for a new event.