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_coefs::f_cll Interface Reference

Public Member Functions

subroutine f_main_cll (F, Fuv, p10, p21, p32, p43, p54, p50, p20, p31, p42, p53, p40, p51, p30, p41, p52, m02, m12, m22, m32, m42, m52, rmax, Ferr, id_in, Ferr2)
 
subroutine f_arrays_cll (F, Fuv, MomInv, masses2, rmax, Ferr, Ferr2)
 
subroutine f_list_cll (F, Fuv, p10, p21, p32, p43, p54, p50, p20, p31, p42, p53, p40, p51, p30, p41, p52, m02, m12, m22, m32, m42, m52, rmax, Ferr, Ferr2)
 
subroutine f_arrays_list_cll (F, Fuv, MomInv, masses2, rmax, Ferr, Ferr2)
 

Detailed Description

Definition at line 64 of file collier_coefs.F90.

Member Function/Subroutine Documentation

◆ f_arrays_cll()

subroutine collier_coefs::f_cll::f_arrays_cll ( double complex, dimension(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  F,
double complex, dimension(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  Fuv,
double complex, dimension(15), intent(in)  MomInv,
double complex, dimension(0:5), intent(in)  masses2,
integer, intent(in)  rmax,
double precision, dimension(0:rmax), intent(out), optional  Ferr,
double precision, dimension(0:rmax), intent(out), optional  Ferr2 
)

Definition at line 2755 of file collier_coefs.F90.

2755 
2756  integer, intent(in) :: rmax
2757  double complex, intent(in) :: MomInv(15), masses2(0:5)
2758  double complex, intent(out) :: F(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax)
2759  double complex, intent(out) :: Fuv(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax)
2760  double precision, optional, intent(out) ::Ferr(0:rmax),Ferr2(0:rmax)
2761  double precision :: Ferraux(0:rmax),Ferr2aux(0:rmax)
2762 
2763  if (present(ferr)) then
2764  if (present(ferr2)) then
2765  call f_main_cll(f,fuv,mominv(1),mominv(2),mominv(3),mominv(4),mominv(5),mominv(6), &
2766  mominv(7),mominv(8),mominv(9),mominv(10),mominv(11),mominv(12), &
2767  mominv(13),mominv(14),mominv(15),masses2(0),masses2(1), &
2768  masses2(2),masses2(3),masses2(4),masses2(5),rmax,ferr,ferr2=ferr2)
2769  else
2770  call f_main_cll(f,fuv,mominv(1),mominv(2),mominv(3),mominv(4),mominv(5),mominv(6), &
2771  mominv(7),mominv(8),mominv(9),mominv(10),mominv(11),mominv(12), &
2772  mominv(13),mominv(14),mominv(15),masses2(0),masses2(1), &
2773  masses2(2),masses2(3),masses2(4),masses2(5),rmax,ferr)
2774  end if
2775  else
2776  if (present(ferr2)) then
2777  call f_main_cll(f,fuv,mominv(1),mominv(2),mominv(3),mominv(4),mominv(5),mominv(6), &
2778  mominv(7),mominv(8),mominv(9),mominv(10),mominv(11),mominv(12), &
2779  mominv(13),mominv(14),mominv(15),masses2(0),masses2(1), &
2780  masses2(2),masses2(3),masses2(4),masses2(5),rmax,ferraux,ferr2=ferr2)
2781  else
2782  call f_main_cll(f,fuv,mominv(1),mominv(2),mominv(3),mominv(4),mominv(5),mominv(6), &
2783  mominv(7),mominv(8),mominv(9),mominv(10),mominv(11),mominv(12), &
2784  mominv(13),mominv(14),mominv(15),masses2(0),masses2(1), &
2785  masses2(2),masses2(3),masses2(4),masses2(5),rmax,ferraux)
2786  end if
2787  end if
2788 

◆ f_arrays_list_cll()

subroutine collier_coefs::f_cll::f_arrays_list_cll ( double complex, dimension(:), intent(out)  F,
double complex, dimension(:), intent(out)  Fuv,
double complex, dimension(15), intent(in)  MomInv,
double complex, dimension(0:5), intent(in)  masses2,
integer, intent(in)  rmax,
double precision, dimension(0:rmax), intent(out), optional  Ferr,
double precision, dimension(0:rmax), intent(out), optional  Ferr2 
)

Definition at line 2897 of file collier_coefs.F90.

2897 
2898  integer, intent(in) :: rmax
2899  double complex, intent(in) :: MomInv(15), masses2(0:5)
2900  double complex, intent(out) :: F(:),Fuv(:)
2901  double precision, optional, intent(out) ::Ferr(0:rmax),Ferr2(0:rmax)
2902  logical :: eflag
2903 
2904  if (6.gt.nmax_cll) then
2905  call seterrflag_cll(-10)
2906  call errout_cll('F_cll','Nmax_cll smaller 6',eflag,.true.)
2907  write(nerrout_cll,*) 'Nmax_cll =',nmax_cll
2908  write(nerrout_cll,*) 'Reinitialize COLLIER with Nmax_cll >= 6'
2909  call propagateerrflag_cll
2910  return
2911  end if
2912  if (rmax.gt.rmax_cll) then
2913  call seterrflag_cll(-10)
2914  call errout_cll('F_cll','argument rmax larger than rmax_cll',eflag,.true.)
2915  write(nerrout_cll,*) 'rmax =',rmax
2916  write(nerrout_cll,*) 'rmax_cll =',rmax_cll
2917  write(nerrout_cll,*) 'Reinitialize COLLIER with rmax_cll >= ',rmax
2918  call propagateerrflag_cll
2919  return
2920  end if
2921 
2922  call f_arrays_list_checked_cll(f,fuv,mominv,masses2,rmax,ferr,ferr2)
2923 

◆ f_list_cll()

subroutine collier_coefs::f_cll::f_list_cll ( double complex, dimension(:), intent(out)  F,
double complex, dimension(:), intent(out)  Fuv,
double complex, intent(in)  p10,
double complex, intent(in)  p21,
double complex, intent(in)  p32,
double complex, intent(in)  p43,
double complex, intent(in)  p54,
double complex, intent(in)  p50,
double complex, intent(in)  p20,
double complex, intent(in)  p31,
double complex, intent(in)  p42,
double complex, intent(in)  p53,
double complex, intent(in)  p40,
double complex, intent(in)  p51,
double complex, intent(in)  p30,
double complex, intent(in)  p41,
double complex, intent(in)  p52,
double complex, intent(in)  m02,
double complex, intent(in)  m12,
double complex, intent(in)  m22,
double complex, intent(in)  m32,
double complex, intent(in)  m42,
double complex, intent(in)  m52,
integer, intent(in)  rmax,
double precision, dimension(0:rmax), intent(out), optional  Ferr,
double precision, dimension(0:rmax), intent(out), optional  Ferr2 
)

Definition at line 2803 of file collier_coefs.F90.

2803 
2804  integer, intent(in) :: rmax
2805  double complex, intent(in) :: p10,p21,p32,p43,p54,p50,p20,p31,p42,p53,p40
2806  double complex, intent(in) :: p51,p30,p41,p52,m02,m12,m22,m32,m42,m52
2807  double complex, intent(out) :: F(:),Fuv(:)
2808  double precision, optional, intent(out) ::Ferr(0:rmax),Ferr2(0:rmax)
2809  logical :: eflag
2810 
2811  if (6.gt.nmax_cll) then
2812  call seterrflag_cll(-10)
2813  call errout_cll('F_cll','Nmax_cll smaller 6',eflag,.true.)
2814  write(nerrout_cll,*) 'Nmax_cll =',nmax_cll
2815  write(nerrout_cll,*) 'Reinitialize COLLIER with Nmax_cll >= 6'
2816  call propagateerrflag_cll
2817  return
2818  end if
2819  if (rmax.gt.rmax_cll) then
2820  call seterrflag_cll(-10)
2821  call errout_cll('F_cll','argument rmax larger than rmax_cll',eflag,.true.)
2822  write(nerrout_cll,*) 'rmax =',rmax
2823  write(nerrout_cll,*) 'rmax_cll =',rmax_cll
2824  write(nerrout_cll,*) 'Reinitialize COLLIER with rmax_cll >= ',rmax
2825  call propagateerrflag_cll
2826  return
2827  end if
2828 
2829  call f_list_checked_cll(f,fuv,p10,p21,p32,p43,p54,p50,p20,p31,p42,p53,p40, &
2830  p51,p30,p41,p52,m02,m12,m22,m32,m42,m52,rmax,ferr,ferr2)
2831 

◆ f_main_cll()

subroutine collier_coefs::f_cll::f_main_cll ( double complex, dimension(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  F,
double complex, dimension(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  Fuv,
double complex, intent(in)  p10,
double complex, intent(in)  p21,
double complex, intent(in)  p32,
double complex, intent(in)  p43,
double complex, intent(in)  p54,
double complex, intent(in)  p50,
double complex, intent(in)  p20,
double complex, intent(in)  p31,
double complex, intent(in)  p42,
double complex, intent(in)  p53,
double complex, intent(in)  p40,
double complex, intent(in)  p51,
double complex, intent(in)  p30,
double complex, intent(in)  p41,
double complex, intent(in)  p52,
double complex, intent(in)  m02,
double complex, intent(in)  m12,
double complex, intent(in)  m22,
double complex, intent(in)  m32,
double complex, intent(in)  m42,
double complex, intent(in)  m52,
integer, intent(in)  rmax,
double precision, dimension(0:rmax), intent(out), optional  Ferr,
integer, intent(in), optional  id_in,
double precision, dimension(0:rmax), intent(out), optional  Ferr2 
)

Definition at line 2371 of file collier_coefs.F90.

2371 
2372  integer, intent(in) :: rmax
2373  double complex, intent(in) :: p10,p21,p32,p43,p54,p50,p20,p31,p42,p53,p40
2374  double complex, intent(in) :: p51,p30,p41,p52,m02,m12,m22,m32,m42,m52
2375  double complex, intent(out) :: F(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax)
2376  double complex, intent(out) :: Fuv(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax)
2377  double precision, optional, intent(out) ::Ferr(0:rmax),Ferr2(0:rmax)
2378  double precision :: q10,q21,q32,q43,q54,q50,q20,q31,q42,q53,q40
2379  double precision :: q51,q30,q41,q52
2380  double complex :: mm02,mm12,mm22,mm32,mm42,mm52
2381  integer, optional, intent(in) :: id_in
2382  double complex :: F2uv(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax)
2383  double complex :: F2(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax)
2384  double complex :: Fdd(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax)
2385  double precision :: Ferraux(0:rmax),Ferr2aux(0:rmax),Fdiff(0:rmax)
2386  double complex :: elimcminf2
2387  double complex :: args(21)
2388  integer :: n0,rank,errflag,id
2389  double precision :: accrelDD(0:rmax_DD),accabsDD(0:rmax_DD),Facc(0:rmax),norm,norm_coli,norm_dd,Facc2(0:rmax)
2390  double precision :: accrel2DD(0:rmax_DD),accabs2DD(0:rmax_DD)
2391  integer :: accflagDD,errflagDD,NDD,rankDD
2392  logical :: mflag,eflag
2393  integer :: r,n1,n2,n3,n4,n5
2394 
2395  if (6.gt.nmax_cll) then
2396  call seterrflag_cll(-10)
2397  call errout_cll('F_cll','Nmax_cll smaller 6',eflag,.true.)
2398  write(nerrout_cll,*) 'Nmax_cll =',nmax_cll
2399  write(nerrout_cll,*) 'Reinitialize COLLIER with Nmax_cll >= 6'
2400  call propagateerrflag_cll
2401  return
2402  end if
2403  if (rmax.gt.rmax_cll) then
2404  call seterrflag_cll(-10)
2405  call errout_cll('F_cll','argument rmax larger than rmax_cll',eflag,.true.)
2406  write(nerrout_cll,*) 'rmax =',rmax
2407  write(nerrout_cll,*) 'rmax_cll =',rmax_cll
2408  write(nerrout_cll,*) 'Reinitialize COLLIER with rmax_cll >= ',rmax
2409  call propagateerrflag_cll
2410  return
2411  end if
2412 
2413  mflag=.true.
2414  if (present(id_in)) then
2415  mflag=.false.
2416  id = id_in
2417  else
2418  id = 0
2419  end if
2420 
2421  if (mflag) then
2422  ! set ID of master call
2423  args(1) = p10
2424  args(2) = p21
2425  args(3) = p32
2426  args(4) = p43
2427  args(5) = p54
2428  args(6) = p50
2429  args(7) = p20
2430  args(8) = p31
2431  args(9) = p42
2432  args(10) = p53
2433  args(11) = p40
2434  args(12) = p51
2435  args(13) = p30
2436  args(14) = p41
2437  args(15) = p52
2438  args(16) = m02
2439  args(17) = m12
2440  args(18) = m22
2441  args(19) = m32
2442  args(20) = m42
2443  args(21) = m52
2444  call setmasterfname_cll('F_cll')
2445  call setmastern_cll(6)
2446  call setmasterr_cll(rmax)
2447  call setmasterargs_cll(21,args)
2448 
2449  call settencache_cll(never_tenred_cll)
2450  end if
2451 
2452 
2453  select case (mode_cll)
2454 
2455  case (1)
2456  ! calculate loop integral using
2457  ! COLI implementation by AD/LH
2458 
2459  call calcf(f,fuv,p10,p21,p32,p43,p54,p50,p20,p31,p42,p53,p40, &
2460  p51,p30,p41,p52,m02,m12,m22,m32,m42,m52,rmax,id,ferraux,ferr2aux)
2461 
2462  norm = abs(f(0,0,0,0,0,0))
2463  do r=1,rmax
2464  do n1=0,r
2465  do n2=0,r-n1
2466  do n3=0,r-n1-n2
2467  do n4=0,r-n1-n2-n3
2468  n5=r-n1-n2-n3-n4
2469  norm = max(norm,abs(f(0,n1,n2,n3,n4,n5)))
2470  end do
2471  end do
2472  end do
2473  end do
2474  end do
2475  if (norm.eq.0d0) then
2476  norm = max(abs(p10),abs(p21),abs(p32),abs(p43),abs(p54), &
2477  abs(p50),abs(p20),abs(p31),abs(p42),abs(p53), &
2478  abs(p40),abs(p51),abs(p30),abs(p41),abs(p52), &
2479  abs(m02),abs(m12),abs(m22),abs(m32),abs(m42),abs(m52))
2480  if(norm.ne.0d0) then
2481  norm=1d0/norm**4
2482  else
2483  norm=1d0/muir2_cll**4
2484  end if
2485  end if
2486  if (norm.ne.0d0) then
2487  facc = ferraux/norm
2488  facc2 = ferr2aux/norm
2489  else
2490  facc = 0d0
2491  facc2 = 0d0
2492  end if
2493 
2494  if (present(ferr)) ferr = ferraux
2495  if (present(ferr2)) ferr2 = ferr2aux
2496 
2497  if (mflag) call propagateaccflag_cll(facc,rmax)
2498 
2499 
2500  case (2)
2501  ! calculate loop integral using
2502  ! DD implementation by SD
2503 
2504  id=0
2505  if (rmax.gt.6) then
2506  call seterrflag_cll(-10)
2507  call errout_cll('F_cll','rank higher than maximum rank implemented in DD library',eflag)
2508  if(eflag) then
2509  write(nerrout_cll,*) 'F_cll: 6-point function of rank>6 not implemented in DD library'
2510  end if
2511  end if
2512 
2513 
2514  ! replace small masses by DD-identifiers
2515  q10 = dreal(getminf2dd_cll(p10))
2516  q21 = dreal(getminf2dd_cll(p21))
2517  q32 = dreal(getminf2dd_cll(p32))
2518  q43 = dreal(getminf2dd_cll(p43))
2519  q54 = dreal(getminf2dd_cll(p54))
2520  q50 = dreal(getminf2dd_cll(p50))
2521  q20 = dreal(getminf2dd_cll(p20))
2522  q31 = dreal(getminf2dd_cll(p31))
2523  q42 = dreal(getminf2dd_cll(p42))
2524  q53 = dreal(getminf2dd_cll(p53))
2525  q40 = dreal(getminf2dd_cll(p40))
2526  q51 = dreal(getminf2dd_cll(p51))
2527  q30 = dreal(getminf2dd_cll(p30))
2528  q41 = dreal(getminf2dd_cll(p41))
2529  q52 = dreal(getminf2dd_cll(p52))
2530  mm02 = getminf2dd_cll(m02)
2531  mm12 = getminf2dd_cll(m12)
2532  mm22 = getminf2dd_cll(m22)
2533  mm32 = getminf2dd_cll(m32)
2534  mm42 = getminf2dd_cll(m42)
2535  mm52 = getminf2dd_cll(m52)
2536 
2537  rank = rmax
2538  call f_dd(fdd,q10,q21,q32,q43,q54,q50,q20,q31,q42,q53,q40, &
2539  q51,q30,q41,q52,mm02,mm12,mm22,mm32,mm42,mm52,rank,id)
2540  f(0:rank/2,0:rank,0:rank,0:rank,0:rank,0:rank) = fdd(0:rank/2,0:rank,0:rank,0:rank,0:rank,0:rank)
2541  fuv = 0d0
2542 
2543  call ddgetacc(accreldd,accabsdd,accrel2dd,accabs2dd,ndd,rankdd,accflagdd,errflagdd,id)
2544  if (present(ferr)) ferr(0:rmax) = accabsdd(0:rmax)
2545  if (present(ferr2)) ferr2(0:rmax) = accabs2dd(0:rmax)
2546 
2547  norm = abs(f(0,0,0,0,0,0))
2548  do r=1,rmax
2549  do n1=0,r
2550  do n2=0,r-n1
2551  do n3=0,r-n1-n2
2552  do n4=0,r-n1-n2-n3
2553  n5=r-n1-n2-n3-n4
2554  norm = max(norm,abs(f(0,n1,n2,n3,n4,n5)))
2555  end do
2556  end do
2557  end do
2558  end do
2559  end do
2560  if (norm.eq.0d0) then
2561  norm = max(abs(p10),abs(p21),abs(p32),abs(p43),abs(p54), &
2562  abs(p50),abs(p20),abs(p31),abs(p42),abs(p53), &
2563  abs(p40),abs(p51),abs(p30),abs(p41),abs(p52), &
2564  abs(m02),abs(m12),abs(m22),abs(m32),abs(m42),abs(m52))
2565  if(norm.ne.0d0) then
2566  norm=1d0/norm**4
2567  else
2568  norm=1d0/muir2_cll**4
2569  end if
2570  end if
2571  if (norm.ne.0d0) then
2572  facc = accabsdd(0:rmax)/norm
2573  facc2 = accabs2dd(0:rmax)/norm
2574  else
2575  facc = 0d0
2576  facc2 = 0d0
2577  end if
2578 
2579  if (mflag) call propagateaccflag_cll(facc,rmax)
2580 
2581 
2582  case (3)
2583  ! cross-check mode
2584  ! compare results for loop integral
2585  ! from COLI implementation by AD/LH and
2586  ! from DD implementation by SD
2587 
2588  ! calculate loop integral
2589  call calcf(f,fuv,p10,p21,p32,p43,p54,p50,p20,p31,p42,p53,p40, &
2590  p51,p30,p41,p52,m02,m12,m22,m32,m42,m52,rmax,id,ferraux,ferr2aux)
2591 
2592 
2593  if (rmax.gt.6) then
2594  call seterrflag_cll(-10)
2595  call errout_cll('F_cll','rank higher than maximum rank implemented in DD library',eflag)
2596  if(eflag) then
2597  write(nerrout_cll,*) 'F_cll: 6-point function of rank>6 not implemented in DD library'
2598  end if
2599  end if
2600 
2601 
2602  ! replace small masses by DD-identifiers
2603  q10 = dreal(getminf2dd_cll(p10))
2604  q21 = dreal(getminf2dd_cll(p21))
2605  q32 = dreal(getminf2dd_cll(p32))
2606  q43 = dreal(getminf2dd_cll(p43))
2607  q54 = dreal(getminf2dd_cll(p54))
2608  q50 = dreal(getminf2dd_cll(p50))
2609  q20 = dreal(getminf2dd_cll(p20))
2610  q31 = dreal(getminf2dd_cll(p31))
2611  q42 = dreal(getminf2dd_cll(p42))
2612  q53 = dreal(getminf2dd_cll(p53))
2613  q40 = dreal(getminf2dd_cll(p40))
2614  q51 = dreal(getminf2dd_cll(p51))
2615  q30 = dreal(getminf2dd_cll(p30))
2616  q41 = dreal(getminf2dd_cll(p41))
2617  q52 = dreal(getminf2dd_cll(p52))
2618  mm02 = getminf2dd_cll(m02)
2619  mm12 = getminf2dd_cll(m12)
2620  mm22 = getminf2dd_cll(m22)
2621  mm32 = getminf2dd_cll(m32)
2622  mm42 = getminf2dd_cll(m42)
2623  mm52 = getminf2dd_cll(m52)
2624 
2625  id=0
2626  rank = rmax
2627  call f_dd(fdd,q10,q21,q32,q43,q54,q50,q20,q31,q42,q53,q40, &
2628  q51,q30,q41,q52,mm02,mm12,mm22,mm32,mm42,mm52,rank,id)
2629  f2(0:rank/2,0:rank,0:rank,0:rank,0:rank,0:rank) = fdd(0:rank/2,0:rank,0:rank,0:rank,0:rank,0:rank)
2630  f2uv = 0d0
2631 
2632  call ddgetacc(accreldd,accabsdd,accrel2dd,accabs2dd,ndd,rankdd,accflagdd,errflagdd,id)
2633 
2634  norm_coli = abs(f(0,0,0,0,0,0))
2635  norm_dd = abs(f2(0,0,0,0,0,0))
2636  do r=1,rmax
2637  do n1=0,r
2638  do n2=0,r-n1
2639  do n3=0,r-n1-n2
2640  do n4=0,r-n1-n2-n3
2641  n5=r-n1-n2-n3-n4
2642  norm_coli = max(norm_coli,abs(f(0,n1,n2,n3,n4,n5)))
2643  norm_dd = max(norm_dd,abs(f2(0,n1,n2,n3,n4,n5)))
2644  end do
2645  end do
2646  end do
2647  end do
2648  end do
2649  if (norm_coli.eq.0d0) then
2650  norm_coli = max(abs(p10),abs(p21),abs(p32),abs(p43),abs(p54), &
2651  abs(p50),abs(p20),abs(p31),abs(p42),abs(p53), &
2652  abs(p40),abs(p51),abs(p30),abs(p41),abs(p52), &
2653  abs(m02),abs(m12),abs(m22),abs(m32),abs(m42),abs(m52))
2654  if(norm_coli.ne.0d0) then
2655  norm_coli=1d0/norm_coli**4
2656  else
2657  norm_coli=1d0/muir2_cll**4
2658  end if
2659  end if
2660  if (norm_dd.eq.0d0) then
2661  norm_dd = max(abs(p10),abs(p21),abs(p32),abs(p43),abs(p54), &
2662  abs(p50),abs(p20),abs(p31),abs(p42),abs(p53), &
2663  abs(p40),abs(p51),abs(p30),abs(p41),abs(p52), &
2664  abs(m02),abs(m12),abs(m22),abs(m32),abs(m42),abs(m52))
2665  if(norm_dd.ne.0d0) then
2666  norm_dd=1d0/norm_dd**4
2667  else
2668  norm_dd=1d0/muir2_cll**4
2669  end if
2670  end if
2671  norm = min(norm_coli,norm_dd)
2672 
2673  ! cross-check
2674  call checkcoefsf_cll(f,f2,p10,p21,p32,p43,p54,p50,p20,p31,p42,p53,p40, &
2675  p51,p30,p41,p52,m02,m12,m22,m32,m42,m52,rmax,norm,fdiff)
2676 
2677 
2678  if (ferraux(rmax).lt.accabsdd(rmax)) then
2679  if (present(ferr)) ferr = max(ferraux,fdiff)
2680  if (present(ferr2)) ferr2 = ferr2aux
2681  if (norm.ne.0d0) then
2682  facc = max(ferraux/norm_coli,fdiff/norm)
2683  facc2 = ferr2aux/norm_coli
2684  else
2685  facc = fdiff
2686  facc2 = 0d0
2687  end if
2688  if (monitoring) pointscntf_coli = pointscntf_coli + 1
2689  else
2690  f = f2
2691  fuv = f2uv
2692  if (present(ferr)) ferr = max(accabsdd(0:rmax),fdiff)
2693  if (present(ferr2)) ferr2 = accabs2dd(0:rmax)
2694  if (norm.ne.0d0) then
2695  facc = max(accabsdd(0:rmax)/norm_dd,fdiff/norm)
2696  facc2 = accabs2dd(0:rmax)/norm_dd
2697  else
2698  facc = fdiff
2699  facc2 = 0d0
2700  end if
2701  if (monitoring) pointscntf_dd = pointscntf_dd + 1
2702  end if
2703 
2704  if (mflag) call propagateaccflag_cll(facc,rmax)
2705 
2706  end select
2707 
2708  if (mflag) call propagateerrflag_cll
2709 
2710  if (monitoring) then
2711  pointscntf_cll = pointscntf_cll + 1
2712 
2713  if(maxval(facc).gt.reqacc_cll) accpointscntf_cll = accpointscntf_cll + 1
2714 
2715  if(maxval(facc).gt.critacc_cll) then
2716  critpointscntf_cll = critpointscntf_cll + 1
2717  if ( critpointscntf_cll.le.noutcritpointsmax_cll(6) ) then
2718  call critpointsout_cll('F_cll',0,maxval(facc), critpointscntf_cll)
2719  if( critpointscntf_cll.eq.noutcritpointsmax_cll(6)) then
2720  write(ncpout_cll,*) ' Further output of Critical Points for F_cll suppressed '
2721  write(ncpout_cll,*)
2722  endif
2723  end if
2724  end if
2725 
2726 #ifdef CritPoints2
2727  if(maxval(facc2).gt.reqacc_cll) accpointscntf2_cll = accpointscntf2_cll + 1
2728 
2729  if(maxval(facc2).gt.critacc_cll) then
2730  critpointscntf2_cll = critpointscntf2_cll + 1
2731  if ( critpointscntf2_cll.le.noutcritpointsmax_cll(6) ) then
2732  call critpointsout2_cll('F_cll',0,maxval(facc2), critpointscntf2_cll)
2733  if( critpointscntf2_cll.eq.noutcritpointsmax_cll(6)) then
2734  write(ncpout2_cll,*) ' Further output of Critical Points for F_cll suppressed '
2735  write(ncpout2_cll,*)
2736  endif
2737  end if
2738  end if
2739 #endif
2740 
2741  end if
2742 

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