JHUGen MELA  JHUGen v7.5.6, MELA v2.4.2
Matrix element calculations as used in JHUGen.
mod_JHUGen.F90
Go to the documentation of this file.
1 MODULE modjhugen
2  use modparameters
3  use modjhugenmela
4  use modkinematics
5 #if compiler==1
6  use ifport
7 #endif
8  implicit none
9 
10  public :: initfirsttime,resetpdfs
11 
12  CONTAINS
13 
14 
15 
16 SUBROUTINE initfirsttime(pdftable,pdfstrlength,pdfmember,collider_sqrts)
18  use modkinematics
19  use modmisc
20  implicit none
21  integer :: pdfstrlength
22  character(len=pdfstrlength) pdftable
23  integer :: pdfmember
24  real(8), intent(in) :: collider_sqrts
25 
26  call setjhugendefaults()
27 
28  includeinterference=.true.
29  includegammastar=.true.
30  mphotoncutoff=4d0*gev
31  widthscheme=0
32  pdfset=3 ! 1: CTEQ6L1 2: MRSW with best fit, 2xx: MSTW with eigenvector set xx=01..40
33  lhapdfstring = pdftable
34  lhapdfmember = pdfmember
35  lenlhapdfstring = pdfstrlength
36 
37 !---------------------------
38  call printlogo(io_stdout, "JHUGen MELA")
39 !---------------------------
40  call resetmubarhgabarh()
41 !---------------------------
42  call computeewvariables()
44  print *, "JHUGen CKM initialization"
45  print *, "Vud = ",vckm_ud
46  print *, "Vus = ",vckm_us
47  print *, "Vub = ",vckm_ub
48  print *, "Vcd = ",vckm_cd
49  print *, "Vcs = ",vckm_cs
50  print *, "Vcb = ",vckm_cb
51  print *, "Vtd = ",vckm_td
52  print *, "Vts = ",vckm_ts
53  print *, "Vtb = ",vckm_tb
54 !---------------------------
55 
56 #if useLHAPDF==1
57  if( lhapdfstring.eq."" ) then
58  print *, "Need to specify pdf file for LHAPDF"
59  stop
60  endif
61 #else
62  if( lhapdfstring.eq."" .and. pdfset.eq.3) then
63  print *, "Need to specify pdf file for NNPDF"
64  stop
65  endif
66 #endif
67 
68  if( isazdecay(decaymode1) ) then
69  m_v = m_z
70  ga_v= ga_z
71  elseif( isawdecay(decaymode1) ) then
72  m_v = m_w
73  ga_v= ga_w
74  elseif( isaphoton(decaymode1) ) then
75  m_v = 0d0
76  ga_v= 0d0
77  endif
78 
79  call initpdfs()
80  IF( collider.EQ.1) THEN
81  collider_energy = collider_sqrts
82  ELSE
83  print *,"Collider not implemented."
84  stop
85  ENDIF
86 
87 #if useCollier==1
88  collier_maxnloopprops = -1
89  collier_maxrank = -1
90 #endif
91  call initcollier(4,3) ! Arguments for ggZH
92 
93  return
94 
95 END SUBROUTINE
96 
97 SUBROUTINE setjhugendefaults()
100  implicit none
101 
102  call setdefaultckm()
103 
104  collider=1
105  vegasit1=-1
106  vegasnc0=-1
107  vegasnc1=-1
108  vegasnc2=-1
109  process = 0 ! select 0, 1 or 2 to represent the spin of the resonance
110  pchannel=2
111  decaymode1=0 ! Z/W+
112  decaymode2=0 ! Z/W-
113 ! ! DecayMode=0: Z --> l+ l- (l=e,mu)
114 ! ! DecayMode=1: Z --> q qbar (q=u,d,c,s,b)
115 ! ! DecayMode=2: Z --> tau+ tau-
116 ! ! DecayMode=3: Z --> nu nubar (nu=nu_e,nu_mu,nu_tau)
117 ! ! DecayMode=4: W --> l nu_l (l=e,mu)
118 ! ! DecayMode=5: W --> q qbar' (q=u,c, qbar'=d,s)
119 ! ! DecayMode=6: W --> tau nu_tau
120 ! ! DecayMode=7: photon
121 ! ! DecayMode=8: Z --> l+ l- (l=e,mu,tau)
122 ! ! DecayMode=9: Z --> anything
123 ! ! DecayMode=10: W --> l nu_l (l=e,mu,tau)
124 ! ! DecayMode=11: W --> anything
125  topdecays=-1
126  taudecays=-1
127  h_dk = .false.
128  unweighted =.true.
129  mufacmultiplier = 1d0
130  murenmultiplier = 1d0
133  offshellreson=.true.
134  offshellv1=.true.
135  offshellv2=.true.
136 
137  lheprodfile=""
138  readlhefile=.false.
139  convertlhefile=.false.
140  readcsmax=.false.
141  readpmzz = .false.
142  pmzzfile="PMZZdistribution.out"
143  pmzzevals=-1
144  doprintpmzz = .false.
145  printpmzzintervals = 20
146  generateevents=.false.
147  requestnleptons = -1
148  requestos=-1
149  requestossf=-1
150  counttauasany = .true.
152 
153  call setuphashes()
154 
155 END SUBROUTINE setjhugendefaults
156 
157 SUBROUTINE initpdfvalues()
159  use modkinematics
160  implicit none
161 
162 #if useLHAPDF==0
163  IF (alphas_mz .LE. 0d0) THEN
164  WRITE(6,*) .le.'alphas_mz 0:',alphas_mz
165  WRITE(6,*) 'continuing with alphas_mz=0.118'
166  alphas_mz=0.118d0
167  ENDIF
168 #endif
169 
170  mu_fact = m_reso ! Set pdf scale to resonance mass by default, later changed as necessary in the EvalWeighted/EvalUnweighted subroutines
171  mu_ren = m_reso ! Set renorm. scale to resonance mass by default, later changed as necessary in the EvalWeighted/EvalUnweighted subroutines
172  call evalalphas() ! Set alphas at default Mu_Ren. Notice ModParameters::ComputeQCDVariables is automatically called!
173  return
174 END SUBROUTINE
175 
176 SUBROUTINE initpdfs()
177 
178  use modparameters
179  use modkinematics
180  implicit none
181 
182 #if useLHAPDF==1
183 
184  implicit none
185  DOUBLE PRECISION alphasPDF
186 
187  call initpdfsetbyname(trim(lhapdfstring,lenlhapdfstring)) ! Let LHAPDF handle everything
188  call initpdf(lhapdfmember)
189 
190  alphas_mz=alphaspdf(zmass_pdf)
191  ! Dummy initialization, just in case. These values are not used.
192  !nloops_pdf = 1
193  zmass_pdf = m_z
194 
195 #else
196 
197  zmass_pdf = m_z ! Take zmass_pdf=M_Z in pdfs that do not specify this value
198 
199  if( pdfset.eq.1 ) then ! CTEQ6L1
200  call setctq6(4) ! 4 CTEQ6L1 Leading Order cteq6l1.tbl
201 
202  alphas_mz=0.130d0
203  !nloops_pdf=1
204  elseif( pdfset.eq.3 ) then ! NNPDF 3.0 LO with a_s=0.13
206  call nninitpdf(lhapdfmember)
207 
208  alphas_mz=0.130d0
209  !nloops_pdf=1
210  zmass_pdf=91.199996948242188d0*gev
211  elseif( (pdfset.eq.2) .or. (pdfset.ge.201 .and. pdfset.le.240) ) then ! MSTW2008 and variations
212  alphas_mz=0.13939d0
213  !nloops_pdf=1
214  else ! Everything else
215  write(6,*) "mod_JHUGen.F90::InitPDFs: PDFSet",pdfset,"QCD parameters are unknown. Please double-check! Stopping JHUGen..."
216  stop
217  ! Could also have used these instead of the stop statement, but why introduce arbitrary number?
218  !alphas_mz = 0.13229060d0
219  !nloops_pdf = 1
220  endif
221 
222 #endif
223 
224  call initpdfvalues() ! Call this only once
225  return
226 END SUBROUTINE
227 
228 subroutine resetpdfs(pdftable,pdfstrlength,pdfmember,pdfid)
229 implicit none
230 integer :: pdfstrlength
231 character(len=pdfstrlength) pdftable
232 integer :: pdfmember
233 integer, intent(in) :: pdfid
234  pdfset=pdfid ! 1: CTEQ6L1 2: MRSW with best fit, 2xx: MSTW with eigenvector set xx=01..40 3: NNPDF3.0
235  lhapdfstring = pdftable
236  lhapdfmember = pdfmember
237  lenlhapdfstring = pdfstrlength
238 #if useLHAPDF==1
239  if( lhapdfstring.eq."" ) then
240  print *, "Need to specify pdf file for LHAPDF"
241  stop
242  endif
243 #else
244  if( lhapdfstring.eq."" .and. pdfset.eq.3) then
245  print *, "Need to specify pdf file for NNPDF"
246  stop
247  endif
248 #endif
249  call initpdfs()
250 end subroutine
251 
252 
253 END MODULE modjhugen
modparameters::alphas_mz
real(dp), public alphas_mz
Definition: mod_Parameters.F90:270
modparameters::includegammastar
logical, public includegammastar
Definition: mod_Parameters.F90:213
modparameters::vegasnc2
integer, public vegasnc2
Definition: mod_Parameters.F90:18
modhashcollection::setuphashes
subroutine setuphashes()
Definition: mod_HashCollection.F90:1042
modhashcollection
Definition: mod_HashCollection.F90:1
modparameters::pdfset
integer, public pdfset
Definition: mod_Parameters.F90:69
modparameters::mphotoncutoff
real(8), public mphotoncutoff
Definition: mod_Parameters.F90:215
modparameters::lhapdfmember
integer, public lhapdfmember
Definition: mod_Parameters.F90:68
modmisc::printlogo
subroutine printlogo(TheUnit, title)
Definition: mod_Misc.F90:1531
modparameters::vckm_cb
real(8), public vckm_cb
Definition: mod_Parameters.F90:282
modparameters::pmzzevals
integer, public pmzzevals
Definition: mod_Parameters.F90:18
modparameters::doprintpmzz
logical, public doprintpmzz
Definition: mod_Parameters.F90:28
modparameters::pmzzfile
character(len=500), public pmzzfile
Definition: mod_Parameters.F90:137
modparameters::vckm_cs
real(8), public vckm_cs
Definition: mod_Parameters.F90:278
modparameters::vegasit1
integer, public vegasit1
Definition: mod_Parameters.F90:18
modparameters::mufacmultiplier
real(8), public mufacmultiplier
Definition: mod_Parameters.F90:21
modparameters::offshellreson
logical, public offshellreson
Definition: mod_Parameters.F90:28
modparameters::counttauasany
logical, public counttauasany
Definition: mod_Parameters.F90:29
modparameters::mu_fact
real(8), public mu_fact
Definition: mod_Parameters.F90:267
modparameters::murenmultiplier
real(8), public murenmultiplier
Definition: mod_Parameters.F90:21
modjhugen
Definition: mod_JHUGen.F90:1
modparameters::ga_w
real(8), public ga_w
Definition: mod_Parameters.F90:229
modjhugenmela::resetmubarhgabarh
subroutine, public resetmubarhgabarh()
Definition: mod_JHUGenMELA.F90:112
modparameters::initcollier
subroutine initcollier(Nmax, Rmax)
Definition: mod_Parameters.F90:3397
modparameters::topdecays
integer, public topdecays
Definition: mod_Parameters.F90:17
modparameters::krenfacscheme_default
integer, parameter, public krenfacscheme_default
Definition: mod_Parameters.F90:32
modkinematics
Definition: mod_Kinematics.F90:1
modparameters::facscheme
integer, public facscheme
Definition: mod_Parameters.F90:20
modparameters::lheprodfile
character(len=500) lheprodfile
Definition: mod_Parameters.F90:64
modparameters::h_dk
logical, public h_dk
Definition: mod_Parameters.F90:180
modparameters::isaphoton
logical function isaphoton(DKMode)
Definition: mod_Parameters.F90:2313
modparameters::convertlhefile
logical, public convertlhefile
Definition: mod_Parameters.F90:28
modjhugen::initfirsttime
subroutine, public initfirsttime(pdftable, pdfstrlength, pdfmember, collider_sqrts)
Definition: mod_JHUGen.F90:17
modparameters::vckm_tb
real(8), public vckm_tb
Definition: mod_Parameters.F90:280
modparameters::writefailedevents
integer, public writefailedevents
Definition: mod_Parameters.F90:30
modparameters::computeewvariables
subroutine computeewvariables()
Definition: mod_Parameters.F90:2785
modparameters::vegasnc1
integer, public vegasnc1
Definition: mod_Parameters.F90:18
modparameters::unweighted
logical, public unweighted
Definition: mod_Parameters.F90:28
modparameters::offshellv1
logical, public offshellv1
Definition: mod_Parameters.F90:28
modparameters::vckm_ts
real(8), public vckm_ts
Definition: mod_Parameters.F90:281
modparameters::pchannel
integer, public pchannel
Definition: mod_Parameters.F90:17
modparameters::includeinterference
logical, public includeinterference
Definition: mod_Parameters.F90:77
modparameters::requestnleptons
integer, dimension(1:2), public requestnleptons
Definition: mod_Parameters.F90:24
modparameters::printpmzzintervals
integer, public printpmzzintervals
Definition: mod_Parameters.F90:143
nninitpdf
subroutine nninitpdf(irep)
Definition: NNPDFDriver.f:50
modparameters::setdefaultckm
subroutine setdefaultckm()
Definition: mod_Parameters.F90:2774
modparameters::isawdecay
logical function isawdecay(DKMode)
Definition: mod_Parameters.F90:2278
modparameters::gev
real(8), parameter, public gev
Definition: mod_Parameters.F90:93
modparameters::collider
integer, public collider
Definition: mod_Parameters.F90:17
modjhugen::resetpdfs
subroutine, public resetpdfs(pdftable, pdfstrlength, pdfmember, pdfid)
Definition: mod_JHUGen.F90:229
modparameters::renscheme
integer, public renscheme
Definition: mod_Parameters.F90:20
setctq6
subroutine setctq6(Iset)
Definition: Cteq61Pdf.f:138
modparameters::lhapdfstring
character(len=100) lhapdfstring
Definition: mod_Parameters.F90:66
modparameters::readlhefile
logical, public readlhefile
Definition: mod_Parameters.F90:28
modparameters::ga_z
real(8), public ga_z
Definition: mod_Parameters.F90:227
modparameters::m_w
real(8), public m_w
Definition: mod_Parameters.F90:228
modparameters::zmass_pdf
real(8), public zmass_pdf
Definition: mod_Parameters.F90:266
modparameters::m_v
real(8), public m_v
Definition: mod_Parameters.F90:78
modparameters::requestossf
integer, dimension(1:2), public requestossf
Definition: mod_Parameters.F90:27
modparameters::readcsmax
logical, public readcsmax
Definition: mod_Parameters.F90:29
modparameters::computeckmelements
subroutine computeckmelements(inVCKM_ud, inVCKM_us, inVCKM_cd, inVCKM_cs, inVCKM_ts, inVCKM_tb, inVCKM_ub, inVCKM_cb, inVCKM_td)
Definition: mod_Parameters.F90:2708
modparameters::ga_v
real(8), public ga_v
Definition: mod_Parameters.F90:78
modparameters
Definition: mod_Parameters.F90:1
modparameters::vckm_ud
real(8), public vckm_ud
Definition: mod_Parameters.F90:276
modmisc
Definition: mod_Misc.F90:1
modparameters::vckm_cd
real(8), public vckm_cd
Definition: mod_Parameters.F90:279
modparameters::vegasnc0
integer, public vegasnc0
Definition: mod_Parameters.F90:18
modparameters::m_z
real(8), public m_z
Definition: mod_Parameters.F90:226
modparameters::vckm_ub
real(8), public vckm_ub
Definition: mod_Parameters.F90:283
modparameters::mu_ren
real(8), public mu_ren
Definition: mod_Parameters.F90:268
modjhugen::initpdfs
subroutine initpdfs()
Definition: mod_JHUGen.F90:177
modparameters::m_reso
real(8), public m_reso
Definition: mod_Parameters.F90:230
modparameters::lenlhapdfstring
integer, public lenlhapdfstring
Definition: mod_Parameters.F90:68
modparameters::requestos
integer, dimension(1:2), public requestos
Definition: mod_Parameters.F90:26
modjhugen::initpdfvalues
subroutine initpdfvalues()
Definition: mod_JHUGen.F90:158
modjhugen::setjhugendefaults
subroutine setjhugendefaults()
Definition: mod_JHUGen.F90:98
modparameters::vckm_us
real(8), public vckm_us
Definition: mod_Parameters.F90:277
modparameters::io_stdout
integer, parameter io_stdout
Definition: mod_Parameters.F90:1138
modparameters::isazdecay
logical function isazdecay(DKMode)
Definition: mod_Parameters.F90:2243
modparameters::decaymode1
integer, public decaymode1
Definition: mod_Parameters.F90:17
modparameters::offshellv2
logical, public offshellv2
Definition: mod_Parameters.F90:28
nnpdfdriver
subroutine nnpdfdriver(gridfilename, lenfilename)
Definition: NNPDFDriver.f:19
modparameters::collider_energy
real(8), public collider_energy
Definition: mod_Parameters.F90:19
modparameters::evalalphas
subroutine evalalphas()
Definition: mod_Parameters.F90:2850
modparameters::decaymode2
integer, public decaymode2
Definition: mod_Parameters.F90:17
modparameters::generateevents
logical, public generateevents
Definition: mod_Parameters.F90:29
modparameters::process
integer, public process
Definition: mod_Parameters.F90:17
modjhugenmela
Definition: mod_JHUGenMELA.F90:1
modparameters::taudecays
integer, public taudecays
Definition: mod_Parameters.F90:17
modparameters::widthscheme
integer, public widthscheme
Definition: mod_Parameters.F90:131
modparameters::readpmzz
logical, public readpmzz
Definition: mod_Parameters.F90:136
modparameters::vckm_td
real(8), public vckm_td
Definition: mod_Parameters.F90:284