Page MenuHomeHEPForge

driver.f
No OneTemporary

driver.f

!-----------------------
! madgraph - a Feynman Diagram package by Tim Stelzer and Bill Long
! (c) 1993
!
! Filename: driver.f
!-----------------------
program MADGRAPH
!*************************************************************************
! Main program for MadGraph
!*************************************************************************
implicit none
! Constants
include 'params.inc'
! Local Variables
integer tops(0:4,0:maxnodes,0:maxtops)
integer graphs(0:maxgraphs,0:maxlines) !graphs, topology internal lines
integer next,i,ntopold,ii,j
double precision sflows(maxflows,maxflows)
logical more
integer fortran
integer, parameter :: unit_par_def=80, unit_par_set=81
! External Functions
! Global Variables
integer iline(-maxlines:maxlines),idir(-maxlines:maxlines)
integer this_coup(max_coup) ,goal_coup(max_coup)
common/to_proc/iline,idir,this_coup,goal_coup
character*15 name
integer iname
common/to_name/iname,name
integer alen(maxgraphs)
integer wnum(-maxlines:maxlines),wcount,wlen(maxgraphs*maxlines)
common/toopt/wcount,wnum,wlen,alen
integer flows(0:2*maxlines,0:maxfactors,0:24,0:maxflows)
integer graphcolor(0:2,0:maxflows,0:maxgraphs)
common/to_color/graphcolor,flows
character*(4*max_particles) particle(4)
integer charge_c(max_particles)
integer iparticle(0:max_particles,0:4),inverse(max_particles)
common/to_model/iparticle, particle, inverse, charge_c
logical this_block(0:max_blocks), need_block(max_blocks)
common/to_blocks/this_block , need_block
logical sumgluons
common/to_cflow/sumgluons
! Data Table
! If you want the output function written using real*8 and complex*16
! style declarations, remove the comment character from the following
! line.
! data fortran /77/
! If you want the output function written using Fortran 90 style
! type declarations, remove the comment character from the following
! line.
data fortran /90/
!-----------
! Begin Code
!-----------
ntopold=0
call resetV
call Load_Particles
call Load_Interactions
open (unit=unit_par_def, status='scratch')
open (unit=unit_par_set, status='scratch')
call Load_Parameters (unit_par_def, unit_par_set) ! added by WK
print*,' '
c call writelogo
print*,' '
call readproc(more)
CC open(unit=isubfile,file='subproc.txt',status='unknown')
CC write(isubfile,'($a)') 'PROCESS = '
CC open(unit=22,status='scratch')
do while (more)
open(unit=91,status='scratch')
flows(0,0,0,0)=0
sumgluons=.true.
sumgluons=.false.
wcount=0
print*,' '
next=iline(0)
c if (next .ne. ntopold) then
c call gentops(next,tops)
c ntopold=next
c endif
c call SortVertex(tops)
call InsertFields(tops,graphs)
if (graphs(0,0) .eq. 0) then
print *, "There are no graphs for this process."
else
print*,'There are ',graphs(0,0),' graphs.'
graphcolor(0,0,0)=graphs(0,0) !Number of graphs
print*,' '
call saveproc(i,graphs(0,0))
call gen_ps(graphs,tops,iline(0))
if (abs(i) .gt. 0) then !New Process
call square(flows,sflows)
! open(unit=99,status='scratch')
call matrixcolor (graphcolor, sflows, fortran)
! rewind 99
rewind 91
c if (fortran .eq. 90) then
call writeamp_f90 (graphs(0,0),next,i,sflows,
& unit_par_def, unit_par_set)
c else
c call writesub(graphs(0,0),next,
c $ goal_coup(2),goal_coup(1),i)
c call writeamp(graphs(0,0),next,
c $ goal_coup(2),goal_coup(1),i,flows)
c end if
! close(99,status='delete')
endif
close(91,status='delete')
endif
c call resetV
write(*,*)
write(*,*)
call readproc(more)
enddo
close (unit_par_def)
close (unit_par_set)
CC rewind 22
CC call writedsig(iline(0),22)
CC close(22,status='delete')
CC close(unit=isubfile)
write(*,*)
write(*,*)
write(*,*) 'Thank you for using MadGraph'
write(*,*)
!************* These are for color flow information only ******************
! call writecolor(graphcolor,sflows,goodflow)
! call writeflows(flows,goodflow)
end
subroutine square(flows,sflows)
!**************************************************************************
! Square terms and simplify to get color factors for flows
!**************************************************************************
implicit none
! Constants
include 'params.inc'
c integer maxlines, maxfactors, maxterms, maxflows
c parameter (maxlines=8,maxfactors=9,maxterms=250,maxflows=200)
! Arguments
integer flows(0:2*maxlines,0:maxfactors,0:24,0:maxflows)
double precision sflows(maxflows,maxflows)
! Local Variables
integer color(0:2*maxlines,0:maxfactors,0:maxterms)
integer nflows,nfactors,nelements
integer iflow,ifact,ielement,iterm,cterm,fterm
integer cflow,cfact,cfactors
integer i,j
double precision sdiag
! Global Variables
integer isflows(2,maxflows,maxflows)
common/to_isflow/isflows
logical sumgluons
common/to_cflow/sumgluons
!-----------
! Begin Code
!-----------
nflows=flows(0,0,0,0)
sumgluons=.true.
do iflow=1,nflows
c write(*,*) 'Flow number ',iflow
c call PrintColors(flows(0,0,0,iflow))
call SimpColors(flows(0,0,0,iflow),-1) !Sum over internal gluon
c call PrintColors(flows(0,0,0,iflow))
enddo
do iflow=1,nflows
do cflow=iflow,nflows
do iterm=1,flows(0,0,0,iflow)
nfactors=flows(0,0,iterm,iflow)
do ifact=0,nfactors
nelements=abs(flows(0,ifact,iterm,iflow))
do ielement=0,nelements
color(ielement,ifact,iterm)=
& flows(ielement,ifact,iterm,iflow)
enddo
enddo
enddo !have all terms copied
color(0,0,0)=flows(0,0,0,iflow) !This is number of terms
do cterm=2,flows(0,0,0,cflow)
do iterm=1,flows(0,0,0,iflow)
call copyterm(color,iterm,-999,-999)
enddo
enddo !now we have enough set for all square terms
do cterm=1,flows(0,0,0,cflow)
cfactors=flows(0,0,cterm,cflow)
do iterm=1,flows(0,0,0,iflow)
nfactors=flows(0,0,iterm,iflow)
fterm = flows(0,0,0,iflow)*(cterm-1)+iterm
color(0,0,fterm) = nfactors+cfactors-1
color(0,1,fterm) = 2 !num and denom
color(1,1,fterm) = color(1,1,fterm)*
& flows(1,1,cterm,cflow)
color(2,1,fterm) = color(2,1,fterm)*
& flows(2,1,cterm,cflow)
do cfact=2,cfactors
nelements=abs(flows(0,cfact,cterm,cflow))
color(0,nfactors+cfact-1,fterm)=
& flows(0,cfact,cterm,cflow)
if (flows(0,cfact,cterm,cflow) .lt. 0) then !T Matrix
color(1,cfact+nfactors-1,fterm) =
& flows(2,cfact,cterm,cflow)
color(2,cfact+nfactors-1,fterm) =
& flows(1,cfact,cterm,cflow)
do ielement=3,nelements
color(ielement,cfact+nfactors-1,fterm)=
& flows(nelements-ielement+3,cfact,
& cterm,cflow)
enddo
else
do ielement=1,nelements
color(ielement,cfact+nfactors-1,fterm)=
& flows(nelements-ielement+1,cfact,
& cterm,cflow)
enddo
endif
enddo
enddo
enddo
call simpcolors(color,1)
call addterms(color)
if (color(0,0,0) .ge. 1) then
isflows(1,iflow,cflow)=color(1,1,1)
isflows(2,iflow,cflow)=color(2,1,1)
isflows(1,cflow,iflow)=color(1,1,1)
isflows(2,cflow,iflow)=color(2,1,1)
sflows(iflow,cflow)=dble(color(1,1,1))/
& dble(color(2,1,1))
sflows(cflow,iflow)=dble(color(1,1,1))/
& dble(color(2,1,1))
if (color(0,0,0) .gt. 1) then
print*,'Error More than one term',color(0,0,1)
call printcolors(color)
elseif (color(0,0,1) .gt. 1) then !one factor
print*,'One term but many factors',color(0,0,1),
& iflow,cflow,color(1,1,1)
call printcolors(color)
endif
else
sflows(iflow,cflow)=0.d0
sflows(cflow,iflow)=0.d0
endif
c if (iflow .eq. cflow .and. sflows(iflow,cflow) .lt. 0) then
c write(*,'(2i4,f21.15)') iflow,cflow,sflows(iflow,cflow)
c write(*,'(4i4)') color(1,1,1),color(2,1,1)
c endif
enddo
enddo
c
c Surprisingly diagonal elements aren't always biggest because
c of extra terms which subtract. So don't use this check anymore
c
c do i=1,nflows
c sdiag = abs(sflows(i,i))
c do j=1,nflows
c if (abs(sflows(j,i)) .gt. sdiag) then
c write(*,'(a,2i4,2f20.5)') 'Warning diag not largest',j,i,
c & sflows(j,i)/sdiag, sdiag
c write(*,'(a,2i4,1f10.5,2i20)') 'Warning',j,i,
c & sflows(j,i)/sdiag, isflows(1,j,i),isflows(2,j,i)
c
c endif
c enddo
c enddo
end

File Metadata

Mime Type
text/x-fortran
Expires
Fri, Apr 4, 9:37 PM (15 h, 58 m)
Storage Engine
blob
Storage Format
Raw Data
Storage Handle
4737423
Default Alt Text
driver.f (10 KB)

Event Timeline