JHUGen MELA  v2.4.1
Matrix element calculations as used in JHUGen. MELA is an important tool that was used for the Higgs boson discovery and for precise measurements of its structure and interactions. Please see the website https://spin.pha.jhu.edu/ and papers cited there for more details, and kindly cite those papers when using this code.
Public Member Functions | List of all members
collier_tensors::aten_cll Interface Reference

Public Member Functions

subroutine aten_main_cll (TA, TAuv, masses2, rmax, TAerr)
 
subroutine aten_list_cll (TA, TAuv, masses2, rmax, TAerr)
 
subroutine aten_args_cll (TA, TAuv, m02, rmax, TAerr)
 
subroutine aten_args_list_cll (TA, TAuv, m02, rmax, TAerr)
 

Detailed Description

Definition at line 44 of file collier_tensors.F90.

Member Function/Subroutine Documentation

◆ aten_args_cll()

subroutine collier_tensors::aten_cll::aten_args_cll ( double complex, dimension(0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  TA,
double complex, dimension(0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  TAuv,
double complex, intent(in)  m02,
integer, intent(in)  rmax,
double precision, dimension(0:rmax), intent(out), optional  TAerr 
)

Definition at line 446 of file collier_tensors.F90.

446 
447  integer, intent(in) :: rmax
448  double complex,intent(in) :: m02
449  double complex, intent(out) :: TA(0:rmax,0:rmax,0:rmax,0:rmax)
450  double complex, intent(out) :: TAuv(0:rmax,0:rmax,0:rmax,0:rmax)
451  double precision, intent(out), optional :: TAerr(0:rmax)
452  double complex :: TA2(0:rmax,0:rmax,0:rmax,0:rmax), TAuv2(0:rmax,0:rmax,0:rmax,0:rmax)
453  double complex :: CA(0:rmax/2), CAuv(0:rmax/2)
454  double precision :: CAerr(0:rmax),TAerr_aux(0:rmax),TAerr_aux2(0:rmax)
455  double complex :: args(1),masses2(0:0)
456  double precision :: TAdiff(0:rmax),norm(0:rmax),norm_coli,norm_dd,TAacc(0:rmax)
457  integer :: r,n0,n1,n2,n3
458  logical :: eflag
459 
460  if (1.gt.nmax_cll) then
461  call seterrflag_cll(-10)
462  call errout_cll('Aten_cll','Nmax_cll smaller 1',eflag,.true.)
463  write(nerrout_cll,*) 'Nmax_cll =',nmax_cll
464  write(nerrout_cll,*) 'Reinitialize COLLIER with Nmax_cll >= 1'
465  call propagateerrflag_cll
466  return
467  end if
468  if (rmax.gt.rmax_cll) then
469  call seterrflag_cll(-10)
470  call errout_cll('Aten_cll','argument rmax larger than rmax_cll',eflag,.true.)
471  write(nerrout_cll,*) 'rmax =',rmax
472  write(nerrout_cll,*) 'rmax_cll =',rmax_cll
473  write(nerrout_cll,*) 'Reinitialize COLLIER with rmax_cll >= ',rmax
474  call propagateerrflag_cll
475  return
476  end if
477 
478  args(1) = m02
479  masses2(0) = m02
480  call setmasterfname_cll('Aten_cll')
481  call setmastern_cll(1)
482  call setmasterr_cll(rmax)
483  call setmasterargs_cll(1,args)
484 
485  call settencache_cll(tenred_cll-1)
486 
487  if (mode_cll.eq.3) then
488  ! calculate tensor with coefficients from COLI
489  mode_cll = 1
490  call a_cll(ca,cauv,m02,rmax,caerr,0)
491  call calctensora(ta,tauv,taerr_aux,ca,cauv,caerr,rmax)
492 
493  ! calculate tensor with coefficients from DD
494  mode_cll = 2
495  call a_cll(ca,cauv,m02,rmax,caerr,0)
496  call calctensora(ta2,tauv2,taerr_aux2,ca,cauv,caerr,rmax)
497 
498  ! comparison --> take better result
499  mode_cll = 3
500  do r=0,rmax
501  norm_coli=0d0
502  norm_dd=0d0
503  do n0=0,r
504  do n1=0,r-n0
505  do n2=0,r-n0-n1
506  n3=r-n0-n1-n2
507  norm_coli = max(norm_coli,abs(ta(n0,n1,n2,n3)))
508  norm_dd = max(norm_dd,abs(ta2(n0,n1,n2,n3)))
509  end do
510  end do
511  end do
512  if (norm_coli.eq.0d0) then
513  norm_coli = abs(masses2(0))
514  if(norm_coli.ne.0d0) then
515  norm_coli=norm_coli**(1+real(r)/2)
516  else
517  norm_coli=muuv2_cll**(1+real(r)/2)
518  end if
519  end if
520  if (norm_dd.eq.0d0) then
521  norm_dd = abs(masses2(0))
522  if(norm_dd.ne.0d0) then
523  norm_dd=norm_dd**(1+real(r)/2)
524  else
525  norm_dd=muuv2_cll**(1+real(r)/2)
526  end if
527  end if
528  norm(r) = min(norm_coli,norm_dd)
529  end do
530 
531  call checktena_cll(ta,ta2,masses2,norm,rmax,tadiff)
532 
533  if (taerr_aux(rmax).lt.taerr_aux2(rmax)) then
534  if (present(taerr)) taerr = max(taerr_aux,tadiff*norm)
535  do r=0,rmax
536  taacc(r) = max(taerr_aux(r)/norm(r),tadiff(r))
537  end do
538  if (monitoring) pointscntaten_coli = pointscntaten_coli + 1
539  else
540  ta = ta2
541  tauv = tauv2
542  if (present(taerr)) taerr = max(taerr_aux2,tadiff*norm)
543  do r=0,rmax
544  taacc(r) = max(taerr_aux2(r)/norm(r),tadiff(r))
545  end do
546  if (monitoring) pointscntaten_dd = pointscntaten_dd + 1
547  end if
548 
549  else
550  call a_cll(ca,cauv,m02,rmax,caerr,0)
551  call calctensora(ta,tauv,taerr_aux,ca,cauv,caerr,rmax)
552  if (present(taerr)) taerr = taerr_aux
553  do r=0,rmax
554  norm(r)=0d0
555  do n0=0,r
556  do n1=0,r-n0
557  do n2=0,r-n0-n1
558  n3=r-n0-n1-n2
559  norm(r) = max(norm(r),abs(ta(n0,n1,n2,n3)))
560  end do
561  end do
562  end do
563  if (norm(r).eq.0d0) then
564  norm(r) = abs(masses2(0))
565  if(norm(r).ne.0d0) then
566  norm(r)=norm(r)**(1+real(r)/2)
567  else
568  norm(r)=muuv2_cll**(1+real(r)/2)
569  end if
570  end if
571  end do
572  do r=0,rmax
573  taacc(r) = taerr_aux(r)/norm(r)
574  end do
575 
576  end if
577 
578  call propagateaccflag_cll(taacc,rmax)
579  call propagateerrflag_cll
580 
581  if (monitoring) then
582  pointscntaten_cll = pointscntaten_cll + 1
583 
584  if(maxval(taacc).gt.reqacc_cll) accpointscntaten_cll = accpointscntaten_cll + 1
585 
586  if(maxval(taacc).gt.critacc_cll) then
587  critpointscntaten_cll = critpointscntaten_cll + 1
588  if ( critpointscntaten_cll.le.noutcritpointsmax_cll(1) ) then
589  call critpointsout_cll('TAten_cll',0,maxval(taacc),critpointscntaten_cll)
590  if( critpointscntaten_cll.eq.noutcritpointsmax_cll(1)) then
591  write(ncpout_cll,*) ' Further output of Critical Points for TAten_cll suppressed'
592  write(ncpout_cll,*)
593  endif
594 #ifdef CritPoints2
595  call critpointsout2_cll('TAten_cll',0,maxval(taacc),critpointscntaten_cll)
596  if( critpointscntaten_cll.eq.noutcritpointsmax_cll(1)) then
597  write(ncpout2_cll,*) ' Further output of Critical Points for TAten_cll suppressed'
598  write(ncpout2_cll,*)
599  endif
600 #endif
601  end if
602  end if
603  end if
604 

◆ aten_args_list_cll()

subroutine collier_tensors::aten_cll::aten_args_list_cll ( double complex, dimension(:), intent(out)  TA,
double complex, dimension(:), intent(out)  TAuv,
double complex, intent(in)  m02,
integer, intent(in)  rmax,
double precision, dimension(0:), intent(out), optional  TAerr 
)

Definition at line 617 of file collier_tensors.F90.

617 
618  integer, intent(in) :: rmax
619  double complex,intent(in) :: m02
620  double complex, intent(out) :: TA(:),TAuv(:)
621  double precision, intent(out), optional :: TAerr(0:)
622  integer :: r,i
623  logical :: eflag
624 
625  if (1.gt.nmax_cll) then
626  call seterrflag_cll(-10)
627  call errout_cll('Aten_cll','Nmax_cll smaller 1',eflag,.true.)
628  write(nerrout_cll,*) 'Nmax_cll =',nmax_cll
629  write(nerrout_cll,*) 'Reinitialize COLLIER with Nmax_cll >= 1'
630  call propagateerrflag_cll
631  return
632  end if
633  if (rmax.gt.rmax_cll) then
634  call seterrflag_cll(-10)
635  call errout_cll('Aten_cll','argument rmax larger than rmax_cll',eflag,.true.)
636  write(nerrout_cll,*) 'rmax =',rmax
637  write(nerrout_cll,*) 'rmax_cll =',rmax_cll
638  write(nerrout_cll,*) 'Reinitialize COLLIER with rmax_cll >= ',rmax
639  call propagateerrflag_cll
640  return
641  end if
642 
643  call aten_args_list_checked_cll(ta,tauv,m02,rmax,taerr)
644 

◆ aten_list_cll()

subroutine collier_tensors::aten_cll::aten_list_cll ( double complex, dimension(:), intent(out)  TA,
double complex, dimension(:), intent(out)  TAuv,
double complex, dimension(0:0), intent(in)  masses2,
integer, intent(in)  rmax,
double precision, dimension(0:rmax), intent(out), optional  TAerr 
)

Definition at line 274 of file collier_tensors.F90.

274 
275  integer, intent(in) :: rmax
276  double complex,intent(in) :: masses2(0:0)
277  double complex, intent(out) :: TA(:),TAuv(:)
278  double precision, intent(out), optional :: TAerr(0:rmax)
279  integer :: r,i
280  logical :: eflag
281 
282  if (1.gt.nmax_cll) then
283  call seterrflag_cll(-10)
284  call errout_cll('Aten_cll','Nmax_cll smaller 1',eflag,.true.)
285  write(nerrout_cll,*) 'Nmax_cll =',nmax_cll
286  write(nerrout_cll,*) 'Reinitialize COLLIER with Nmax_cll >= 1'
287  call propagateerrflag_cll
288  return
289  end if
290  if (rmax.gt.rmax_cll) then
291  call seterrflag_cll(-10)
292  call errout_cll('Aten_cll','argument rmax larger than rmax_cll',eflag,.true.)
293  write(nerrout_cll,*) 'rmax =',rmax
294  write(nerrout_cll,*) 'rmax_cll =',rmax_cll
295  write(nerrout_cll,*) 'Reinitialize COLLIER with rmax_cll >= ',rmax
296  call propagateerrflag_cll
297  return
298  end if
299 
300  call aten_list_checked_cll(ta,tauv,masses2,rmax,taerr)
301 

◆ aten_main_cll()

subroutine collier_tensors::aten_cll::aten_main_cll ( double complex, dimension(0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  TA,
double complex, dimension(0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  TAuv,
double complex, dimension(0:0), intent(in)  masses2,
integer, intent(in)  rmax,
double precision, dimension(0:rmax), intent(out), optional  TAerr 
)

Definition at line 104 of file collier_tensors.F90.

104 
105  integer, intent(in) :: rmax
106  double complex,intent(in) :: masses2(0:0)
107  double complex, intent(out) :: TA(0:rmax,0:rmax,0:rmax,0:rmax)
108  double complex, intent(out) :: TAuv(0:rmax,0:rmax,0:rmax,0:rmax)
109  double precision, intent(out), optional :: TAerr(0:rmax)
110  double complex :: TA2(0:rmax,0:rmax,0:rmax,0:rmax), TAuv2(0:rmax,0:rmax,0:rmax,0:rmax)
111  double complex :: CA(0:rmax/2), CAuv(0:rmax/2)
112  double precision :: CAerr(0:rmax),TAerr_aux(0:rmax),TAerr_aux2(0:rmax)
113  double complex :: args(1)
114  double precision :: TAdiff(0:rmax),norm(0:rmax),norm_coli,norm_dd,TAacc(0:rmax)
115  integer :: r,n0,n1,n2,n3
116  logical :: eflag
117 
118  if (1.gt.nmax_cll) then
119  call seterrflag_cll(-10)
120  call errout_cll('Aten_cll','Nmax_cll smaller 1',eflag,.true.)
121  write(nerrout_cll,*) 'Nmax_cll =',nmax_cll
122  write(nerrout_cll,*) 'Reinitialize COLLIER with Nmax_cll >= 1'
123  call propagateerrflag_cll
124  return
125  end if
126  if (rmax.gt.rmax_cll) then
127  call seterrflag_cll(-10)
128  call errout_cll('Aten_cll','argument rmax larger than rmax_cll',eflag,.true.)
129  write(nerrout_cll,*) 'rmax =',rmax
130  write(nerrout_cll,*) 'rmax_cll =',rmax_cll
131  write(nerrout_cll,*) 'Reinitialize COLLIER with rmax_cll >= ',rmax
132  call propagateerrflag_cll
133  return
134  end if
135 
136  args(1) = masses2(0)
137  call setmasterfname_cll('Aten_cll')
138  call setmastern_cll(1)
139  call setmasterr_cll(rmax)
140  call setmasterargs_cll(1,args)
141 
142  call settencache_cll(tenred_cll-1)
143 
144  if (mode_cll.eq.3) then
145  ! calculate tensor with coefficients from COLI
146  mode_cll = 1
147  call a_cll(ca,cauv,masses2(0),rmax,caerr,0)
148  call calctensora(ta,tauv,taerr_aux,ca,cauv,caerr,rmax)
149 
150  ! calculate tensor with coefficients from DD
151  mode_cll = 2
152  call a_cll(ca,cauv,masses2(0),rmax,caerr,0)
153  call calctensora(ta2,tauv2,taerr_aux2,ca,cauv,caerr,rmax)
154 
155  ! comparison --> take better result
156  mode_cll = 3
157  do r=0,rmax
158  norm_coli=0d0
159  norm_dd=0d0
160  do n0=0,r
161  do n1=0,r-n0
162  do n2=0,r-n0-n1
163  n3=r-n0-n1-n2
164  norm_coli = max(norm_coli,abs(ta(n0,n1,n2,n3)))
165  norm_dd = max(norm_dd,abs(ta2(n0,n1,n2,n3)))
166  end do
167  end do
168  end do
169  if (norm_coli.eq.0d0) then
170  norm_coli = abs(masses2(0))
171  if(norm_coli.ne.0d0) then
172  norm_coli=norm_coli**(1+real(r)/2)
173  else
174  norm_coli=muuv2_cll**(1+real(r)/2)
175  end if
176  end if
177  if (norm_dd.eq.0d0) then
178  norm_dd = abs(masses2(0))
179  if(norm_dd.ne.0d0) then
180  norm_dd=norm_dd**(1+real(r)/2)
181  else
182  norm_dd=muuv2_cll**(1+real(r)/2)
183  end if
184  end if
185  norm(r) = min(norm_coli,norm_dd)
186  end do
187 
188  call checktena_cll(ta,ta2,masses2,norm,rmax,tadiff)
189 
190  if (taerr_aux(rmax).lt.taerr_aux2(rmax)) then
191  if (present(taerr)) taerr = max(taerr_aux,tadiff*norm)
192  do r=0,rmax
193  taacc(r) = max(taerr_aux(r)/norm(r),tadiff(r))
194  end do
195  if (monitoring) pointscntaten_coli = pointscntaten_coli + 1
196  else
197  ta = ta2
198  tauv = tauv2
199  if (present(taerr)) taerr = max(taerr_aux2,tadiff*norm)
200  do r=0,rmax
201  taacc(r) = max(taerr_aux2(r)/norm(r),tadiff(r))
202  end do
203  if (monitoring) pointscntaten_dd = pointscntaten_dd + 1
204  end if
205 
206  else
207  call a_cll(ca,cauv,masses2(0),rmax,caerr,0)
208  call calctensora(ta,tauv,taerr_aux,ca,cauv,caerr,rmax)
209  if (present(taerr)) taerr = taerr_aux
210  do r=0,rmax
211  norm(r)=0d0
212  do n0=0,r
213  do n1=0,r-n0
214  do n2=0,r-n0-n1
215  n3=r-n0-n1-n2
216  norm(r) = max(norm(r),abs(ta(n0,n1,n2,n3)))
217  end do
218  end do
219  end do
220  if (norm(r).eq.0d0) then
221  norm(r) = abs(masses2(0))
222  if(norm(r).ne.0d0) then
223  norm(r)=norm(r)**(1+real(r)/2)
224  else
225  norm(r)=muuv2_cll**(1+real(r)/2)
226  end if
227  end if
228  end do
229  do r=0,rmax
230  taacc(r) = taerr_aux(r)/norm(r)
231  end do
232 
233  end if
234 
235  call propagateaccflag_cll(taacc,rmax)
236  call propagateerrflag_cll
237 
238  if (monitoring) then
239  pointscntaten_cll = pointscntaten_cll + 1
240 
241  if(maxval(taacc).gt.reqacc_cll) accpointscntaten_cll = accpointscntaten_cll + 1
242 
243  if(maxval(taacc).gt.critacc_cll) then
244  critpointscntaten_cll = critpointscntaten_cll + 1
245  if ( critpointscntaten_cll.le.noutcritpointsmax_cll(1) ) then
246  call critpointsout_cll('TAten_cll',0,maxval(taacc),critpointscntaten_cll)
247  if( critpointscntaten_cll.eq.noutcritpointsmax_cll(1)) then
248  write(ncpout_cll,*) ' Further output of Critical Points for TAten_cll suppressed'
249  write(ncpout_cll,*)
250  endif
251 #ifdef CritPoints2
252  call critpointsout2_cll('TAten_cll',0,maxval(taacc),critpointscntaten_cll)
253  if( critpointscntaten_cll.eq.noutcritpointsmax_cll(1)) then
254  write(ncpout2_cll,*) ' Further output of Critical Points for TAten_cll suppressed'
255  write(ncpout2_cll,*)
256  endif
257 #endif
258  end if
259  end if
260  end if
261 

The documentation for this interface was generated from the following file:
endif
O0 g endif() string(TOLOWER "$
Definition: CMakeLists.txt:143