JHUGen MELA  JHUGen v7.5.6, MELA v2.4.2
Matrix element calculations as used in JHUGen.
collier_coefs::b_cll Interface Reference

Public Member Functions

subroutine b_main_cll (B, Buv, p10, m02, m12, rmax, Berr, id_in)
 
subroutine b_arrays_cll (B, Buv, MomInv, masses2, rmax, Berr)
 
subroutine b_list_cll (B, Buv, p10, m02, m12, rmax, Berr)
 
subroutine b_arrays_list_cll (B, Buv, MomInv, masses2, rmax, Berr)
 

Detailed Description

Definition at line 40 of file collier_coefs.F90.

Member Function/Subroutine Documentation

◆ b_arrays_cll()

subroutine collier_coefs::b_cll::b_arrays_cll ( double complex, dimension(0:rmax/2,0:rmax), intent(out)  B,
double complex, dimension(0:rmax/2,0:rmax), intent(out)  Buv,
double complex, dimension(1), intent(in)  MomInv,
double complex, dimension(0:1), intent(in)  masses2,
integer, intent(in)  rmax,
double precision, dimension(0:rmax), intent(out), optional  Berr 
)

Definition at line 534 of file collier_coefs.F90.

534 
535  integer, intent(in) :: rmax
536  double complex, intent(in) :: MomInv(1), masses2(0:1)
537  double complex, intent(out) :: Buv(0:rmax/2,0:rmax)
538  double complex, intent(out) :: B(0:rmax/2,0:rmax)
539  double precision, optional, intent(out) :: Berr(0:rmax)
540  double precision :: Berraux(0:rmax)
541 
542  if (present(berr)) then
543  call b_main_cll(b,buv,mominv(1),masses2(0),masses2(1),rmax,berr)
544  else
545  call b_main_cll(b,buv,mominv(1),masses2(0),masses2(1),rmax,berraux)
546  end if
547 

◆ b_arrays_list_cll()

subroutine collier_coefs::b_cll::b_arrays_list_cll ( double complex, dimension(1:), intent(out)  B,
double complex, dimension(1:), intent(out)  Buv,
double complex, dimension(1), intent(in)  MomInv,
double complex, dimension(0:1), intent(in)  masses2,
integer, intent(in)  rmax,
double precision, dimension(0:rmax), intent(out), optional  Berr 
)

Definition at line 630 of file collier_coefs.F90.

630 
631  integer, intent(in) :: rmax
632  double complex, intent(in) :: MomInv(1),masses2(0:1)
633  double precision, optional, intent(out) :: Berr(0:rmax)
634  double complex, intent(out) :: Buv(1:),B(1:)
635  logical :: eflag
636 
637  if (2.gt.nmax_cll) then
638  call seterrflag_cll(-10)
639  call errout_cll('B_cll','Nmax_cll smaller 2',eflag,.true.)
640  write(nerrout_cll,*) 'Nmax_cll =',nmax_cll
641  write(nerrout_cll,*) 'Reinitialize COLLIER with Nmax_cll >= 2'
642  call propagateerrflag_cll
643  return
644  end if
645  if (rmax.gt.rmax_cll) then
646  call seterrflag_cll(-10)
647  call errout_cll('B_cll','argument rmax larger than rmax_cll',eflag,.true.)
648  write(nerrout_cll,*) 'rmax =',rmax
649  write(nerrout_cll,*) 'rmax_cll =',rmax_cll
650  write(nerrout_cll,*) 'Reinitialize COLLIER with rmax_cll >= ',rmax
651  call propagateerrflag_cll
652  return
653  end if
654 
655  call b_arrays_list_checked_cll(b,buv,mominv,masses2,rmax,berr)
656 

◆ b_list_cll()

subroutine collier_coefs::b_cll::b_list_cll ( double complex, dimension(1:), intent(out)  B,
double complex, dimension(1:), intent(out)  Buv,
double complex, intent(in)  p10,
double complex, intent(in)  m02,
double complex, intent(in)  m12,
integer, intent(in)  rmax,
double precision, dimension(0:rmax), intent(out), optional  Berr 
)

Definition at line 560 of file collier_coefs.F90.

560 
561  integer, intent(in) :: rmax
562  double complex, intent(in) :: p10,m02,m12
563  double complex, intent(out) :: Buv(1:),B(1:)
564  double precision, optional, intent(out) :: Berr(0:rmax)
565  double precision :: Berraux(0:rmax)
566  logical :: eflag
567 
568  if (2.gt.nmax_cll) then
569  call seterrflag_cll(-10)
570  call errout_cll('B_cll','Nmax_cll smaller 2',eflag,.true.)
571  write(nerrout_cll,*) 'Nmax_cll =',nmax_cll
572  write(nerrout_cll,*) 'Reinitialize COLLIER with Nmax_cll >= 2'
573  call propagateerrflag_cll
574  return
575  end if
576  if (rmax.gt.rmax_cll) then
577  call seterrflag_cll(-10)
578  call errout_cll('B_cll','argument rmax larger than rmax_cll',eflag,.true.)
579  write(nerrout_cll,*) 'rmax =',rmax
580  write(nerrout_cll,*) 'rmax_cll =',rmax_cll
581  write(nerrout_cll,*) 'Reinitialize COLLIER with rmax_cll >= ',rmax
582  call propagateerrflag_cll
583  return
584  end if
585 
586  call b_list_checked_cll(b,buv,p10,m02,m12,rmax,berr)
587 

◆ b_main_cll()

subroutine collier_coefs::b_cll::b_main_cll ( double complex, dimension(0:rmax/2,0:rmax), intent(out)  B,
double complex, dimension(0:rmax/2,0:rmax), intent(out)  Buv,
double complex, intent(in)  p10,
double complex, intent(in)  m02,
double complex, intent(in)  m12,
integer, intent(in)  rmax,
double precision, dimension(0:rmax), intent(out), optional  Berr,
integer, intent(in), optional  id_in 
)

Definition at line 327 of file collier_coefs.F90.

327 
328  integer, intent(in) :: rmax
329  double complex, intent(in) :: p10,m02,m12
330  double precision :: q10
331  double complex :: mm02,mm12
332  double complex, intent(out) :: Buv(0:rmax/2,0:rmax)
333  double complex, intent(out) :: B(0:rmax/2,0:rmax)
334  double precision, optional, intent(out) :: Berr(0:rmax)
335  integer, optional, intent(in) :: id_in
336  double complex :: B2uv(0:rmax/2,0:rmax), B2(0:rmax/2,0:rmax)
337  double complex :: Bcoliuv(0:rmax,0:rmax)
338  double complex :: Bcoli(0:rmax,0:rmax)
339  double complex :: Bdduv(0:rmax,0:rmax)
340  double complex :: Bdd(0:rmax,0:rmax)
341  double precision :: Berraux(0:rmax),Bdiff(0:rmax)
342  double complex :: args(3)
343  integer :: n0,rank,errflag,id,r
344  double precision :: accrelDD(0:rmax_DD),accabsDD(0:rmax_DD)
345  double precision :: accrel2DD(0:rmax_DD),accabs2DD(0:rmax_DD)
346  double precision :: Bacc(0:rmax),Bacc2(0:rmax),norm,norm_coli,norm_dd
347  integer :: accflagDD,errflagDD,NDD,rankDD
348  logical :: mflag,eflag
349 
350  if (2.gt.nmax_cll) then
351  call seterrflag_cll(-10)
352  call errout_cll('B_cll','Nmax_cll smaller 2',eflag,.true.)
353  write(nerrout_cll,*) 'Nmax_cll =',nmax_cll
354  write(nerrout_cll,*) 'Reinitialize COLLIER with Nmax_cll >= 2'
355  call propagateerrflag_cll
356  return
357  end if
358  if (rmax.gt.rmax_cll) then
359  call seterrflag_cll(-10)
360  call errout_cll('B_cll','argument rmax larger than rmax_cll',eflag,.true.)
361  write(nerrout_cll,*) 'rmax =',rmax
362  write(nerrout_cll,*) 'rmax_cll =',rmax_cll
363  write(nerrout_cll,*) 'Reinitialize COLLIER with rmax_cll >= ',rmax
364  call propagateerrflag_cll
365  return
366  end if
367 
368  mflag=.true.
369  if (present(id_in)) then
370  mflag=.false.
371  id = id_in
372  else
373  id = 0
374  end if
375 
376  if (mflag) then
377  ! set ID of master call
378  args(1) = p10
379  args(2) = m02
380  args(3) = m12
381  call setmasterfname_cll('B_cll')
382  call setmastern_cll(2)
383  call setmasterr_cll(rmax)
384  call setmasterargs_cll(3,args)
385 
386  call settencache_cll(never_tenred_cll)
387  end if
388 
389 
390  select case (mode_cll)
391 
392  case (1)
393  ! calculate loop integral using
394  ! COLI implementation by AD/LH
395 
396  call calcb(bcoli,bcoliuv,p10,m02,m12,rmax,id,berraux)
397 
398  norm = maxval(abs(bcoli(0,0:rmax)))
399  if (norm.ne.0d0) then
400  bacc = berraux/norm
401  else
402  bacc = berraux
403  end if
404 
405  if (present(berr)) then
406  berr = berraux
407  end if
408 
409  if (mflag) call propagateaccflag_cll(bacc,rmax)
410 
411  b(0:rmax/2,0:rmax) = bcoli(0:rmax/2,0:rmax)
412  buv(0:rmax/2,0:rmax) = bcoliuv(0:rmax/2,0:rmax)
413 
414  case (2)
415  ! calculate loop integral using
416  ! DD implementation by SD
417 
418  id=0
419 
420  ! replace small masses by DD-identifiers
421  q10 = dreal(getminf2dd_cll(p10))
422  mm02 = getminf2dd_cll(m02)
423  mm12 = getminf2dd_cll(m12)
424 
425  rank = rmax
426  call b_dd(bdd,bdduv,q10,mm02,mm12,rank,id)
427  do n0=0,rank/2
428  b(n0,0:rank) = bdd(n0,0:rank)
429  buv(n0,0:rank) = bdduv(n0,0:rank)
430  end do
431 
432  call ddgetacc(accreldd,accabsdd,accrel2dd,accabs2dd,ndd,rankdd,accflagdd,errflagdd,id)
433  if (present(berr)) then
434  berr(0:rmax) = accabsdd(0:rmax)
435  end if
436 
437  norm = maxval(abs(b(0,0:rmax)))
438  if (norm.ne.0d0) then
439  bacc = accabsdd(0:rmax)/norm
440  else
441  bacc = accabsdd(0:rmax)
442  end if
443  if (mflag) call propagateaccflag_cll(bacc,rmax)
444 
445 
446  case (3)
447  ! cross-check mode
448  ! compare results for loop integral
449  ! from COLI implementation by AD/LH and
450  ! from DD implementation by SD
451 
452  ! calculate loop integral using COLI
453  call calcb(bcoli,bcoliuv,p10,m02,m12,rmax,id,berraux)
454 
455  b(0:rmax/2,0:rmax) = bcoli(0:rmax/2,0:rmax)
456  buv(0:rmax/2,0:rmax) = bcoliuv(0:rmax/2,0:rmax)
457 
458 
459  ! calculate loop integral using DD
460 
461  id=0
462 
463  ! replace small masses by DD-identifiers
464  q10 = dreal(getminf2dd_cll(p10))
465  mm02 = getminf2dd_cll(m02)
466  mm12 = getminf2dd_cll(m12)
467 
468  rank = rmax
469  call b_dd(bdd,bdduv,q10,mm02,mm12,rank,0)
470  do n0=0,rank/2
471  b2(n0,0:rmax) = bdd(n0,0:rmax)
472  b2uv(n0,0:rmax) = bdduv(n0,0:rmax)
473  end do
474  call ddgetacc(accreldd,accabsdd,accrel2dd,accabs2dd,ndd,rankdd,accflagdd,errflagdd,id)
475 
476  norm_coli = maxval(abs(b(0,0:rmax)))
477  if (norm_coli.eq.0d0) norm_coli = 1d0
478  norm_dd = maxval(abs(b2(0,0:rmax)))
479  if (norm_dd.eq.0d0) norm_dd = 1d0
480  norm = min(norm_coli,norm_dd)
481 
482  ! cross-check
483  call checkcoefsb_cll(b,b2,p10,m02,m12,rmax,norm,bdiff)
484 
485  if (berraux(rmax).lt.accabsdd(rmax)) then
486  if (present(berr)) berr = max(berraux,bdiff)
487  bacc = max(berraux/norm_coli,bdiff/norm)
488  if (monitoring) pointscntb_coli = pointscntb_coli + 1
489  else
490  b = b2
491  buv = b2uv
492  if (present(berr)) berr = max(accabsdd(0:rmax),bdiff)
493  bacc = max(accabsdd(0:rmax)/norm_dd,bdiff/norm)
494  if (monitoring) pointscntb_dd = pointscntb_dd + 1
495  end if
496 
497  if (mflag) call propagateaccflag_cll(bacc,rmax)
498 
499  end select
500 
501  if (mflag) call propagateerrflag_cll
502 
503  if (monitoring) then
504  pointscntb_cll = pointscntb_cll + 1
505 
506  if(maxval(bacc).gt.reqacc_cll) accpointscntb_cll = accpointscntb_cll + 1
507 
508  if(maxval(bacc).gt.critacc_cll) then
509  critpointscntb_cll = critpointscntb_cll + 1
510  if ( critpointscntb_cll.le.noutcritpointsmax_cll(2) ) then
511  call critpointsout_cll('B_cll',0,maxval(bacc), critpointscntb_cll)
512  if( critpointscntb_cll.eq.noutcritpointsmax_cll(2)) then
513  write(ncpout_cll,*) ' Further output of Critical Points for B_cll suppressed '
514  write(ncpout_cll,*)
515  endif
516  end if
517  end if
518 
519  end if
520 
521 

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