How to implement new cuts
- Keywords:
- cuts madevent LO NLO
- Last updated by:
- Olivier Mattelaer
For version before 3.5.0, see https:/
Since 3.5.0:
1: You can prepare a file somewhere in your filesystem containing the fortran definition of your cuts (see below for LO and NLO template)
2: You can pass the path to that file within the run_card.dat within the custom_fcts entry of the run_card (this entry expects a list of path), i.e. with the following line:
['/home/
3: If your function needs extra parameter to be controlled directly via the run_card by adding a line in the run_card like "10.0 = min_newvar". Your fortran function could then use the parameter "min_newvar" as long as you keep the "include run.inc"
Remember that madgraph works in the center of mass of the parton and that he is using symmetry to speed-up the code.
(including flavor, beam symmetry, identical particle,...) Therefore, not all cuts are possible to implement within the code.
==================
Template for LO cuts:
==================
logical FUNCTION dummy_cuts(P)
C******
C INPUT:
C P(0:3,1) MOMENTUM OF INCOMING PARTON (partonic center of mass)
C P(0:3,2) MOMENTUM OF INCOMING PARTON
C P(0:3,3) MOMENTUM OF ...
C ALL MOMENTA ARE IN THE REST FRAME of the parton
C OUTPUT:
C TRUE IF EVENTS PASSES ALL CUTS LISTED
C******
C For the demo: implement pt cut on particle 3 and 4 via custom entry of the run_card
C my_param and my_param2
IMPLICIT NONE
c
c Constants
c
include 'genps.inc'
include 'nexternal.inc'
include 'run.inc'
C
C ARGUMENTS
C
REAL*8 P(0:3,nexternal)
C
C PARAMETERS
C
real*8 PI
parameter( PI = 3.1415926535897
c
c particle identification
c
LOGICAL IS_A_J(
LOGICAL IS_A_B(
LOGICAL IS_A_NU(
logical do_cuts(nexternal)
COMMON /TO_SPECISA/
. IS_A_ONIUM, do_cuts
double precision pt
dummy_
if (pt(P(0,
return
endif
if (pt(P(0,
return
endif
==================
Template for NLO cuts:
==================
c NOTE THAT ONLY IRC-SAFE CUTS CAN BE APPLIED OTHERWISE THE INTEGRATION
c MIGHT NOT CONVERGE
c
logical function dummy_cuts(
implicit none
c This includes the 'nexternal' parameter that labels the number of
c particles in the (n+1)-body process
include 'nexternal.inc'
c This include file contains common blocks filled with the cuts defined
c in the run_card.dat (including custom set entry)
include 'run.inc'
include 'cuts.inc'
c
c This is an array which is '-1' for initial state and '1' for final
c state particles
integer istatus(nexternal)
c This is an array with (simplified) PDG codes for the particles. Note
c that channels that are combined (i.e. they have the same matrix
c elements) are given only 1 set of PDG codes. This means, e.g., that
c when using a 5-flavour scheme calculation (massless b quark), no
c b-tagging can be applied.
integer iPDG(nexternal)
c The array of the momenta and masses of the initial and final state
c particles in the lab frame. The format is "E, px, py, pz, mass", while
c the second dimension loops over the particles in the process. Note
c that these are the (n+1)-body particles; for the n-body there is one
c momenta equal to all zero's (this is not necessarily the last particle
c in the list). If one uses IR-safe obserables only, there should be no
c difficulty in using this.
double precision p(0:4,nexternal)
c
C external functions that can be used. Some are defined in cuts.f
C others are in ./Source/
REAL*8 R2_04,invm2_
external R2_04,invm2_
dummy_cuts = .true. ! event is okay; otherwise it is changed
c$$$C EXAMPLE: cut on top quark pT
c$$$C Note that PDG specific cut are more optimised than simple user cut
c$$$ do i=1,nexternal ! loop over all external particles
c$$$ if (istatus(i).eq.1 ! final state particle
c$$$ & .and. abs(ipdg(i)).eq.6) then ! top quark
c$$$C apply the pT cut (pT should be large than 200 GeV for the event to
c$$$C pass cuts)
c$$$ if ( p(1,i)**2+p(2,i)**2 .lt. 200d0**2 ) then
c$$$C momenta do not pass cuts. Set passcuts_user to false and return
c$$$ dummy_cuts=.false.
c$$$ return
c$$$ endif
c$$$ endif
c$$$ enddo
c
return
end