######################################################################## Here follows a list of the routines that become available if the module "avh_olo" is USEd. subroutine olo_precision( ndecimals ) subroutine olo_unit( iunit ,message ) subroutine olo_scale( muscale ) subroutine olo_onshell( thrs ) subroutine olo_a0( rslt ,mm ,rmu ) subroutine olo_an( rslt ,rank ,mm ,rmu ) subroutine olo_b0( rslt ,pp ,m1,m2 ,rmu ) subroutine olo_b11( b11,b00,b1,b0 ,pp,m1,m2 ,rmu ) subroutine olo_bn( rslt ,rank ,pp,m1,m2 ,rmu ) subroutine olo_c0( rslt ,p1,p2,p3 ,m1,m2,m3 ,rmu ) subroutine olo_d0( rslt ,p1,p2,p3,p4,p12,p23 ,m1,m2,m3,m4 ,rmu ) integer function olo_dp_precision() integer function olo_qp_precision() integer function olo_dd_precision() integer function olo_qd_precision() integer function olo_mp_precision() They are described one by one below. The routines for the evaluation of loop integrals are generic, ie they can be called with real and with complex arguments, and with or without the last argument. If the program is compiled with several levels of precision, the same routines are called with different kinds and types of input. The output is complex of the same kind/type as the input. The extensions "_a0", "_b0" etc. may be omitted, you can also just call olo( rslt ,mm ,rmu ) call olo( rslt ,rank ,mm ,rmu ) call olo( rslt ,pp ,m1,m2 ,rmu ) call olo( b11,b00,b1,b0 ,pp,m1,m2 ,rmu ) call olo( rslt ,rank ,pp,m1,m2 ,rmu ) call olo( rslt ,p1,p2,p3 ,m1,m2,m3 ,rmu ) call olo( rslt ,p1,p2,p3,p4,p12,p23 ,m1,m2,m3,m4 ,rmu ) If the names of routines in "avh_olo" lead to name clashes, realize that these can be cured in the USE statement, eg use avh_olo, my_olo=>olo and then call my_olo( rslt ,pp ,m1,m2 ,rmu ) etc. ######################################################################## subroutine olo_precision( ndecimals ) integer ,intent(in) :: ndecimals If OneLOop is operating at arbitrary precision, the number of decimals can be set with this routine. It is the responsibility of the user that this number corresponds to the actual precision of the variables. The precision can be changed during runtime. This has no effect on routine calls with Fortran intrinsic types, or types with fixed precision. ######################################################################## subroutine olo_unit( iunit ,message ) integer ,intent(in) :: iunit character(*),intent(in),optional :: message Set the units to which messages are send. The argument "message" specifies the type of message for which you want to set the unit. At the moment, the options and their defaults are 'error'--> 6, 'warning'--> 6, 'message'--> 6, 'printall'--> 0. -If the unit for 'printall' is positive, then all input and output for all calls to the loop functions will be printed to that unit. -If iunit is not positive, then nothing will be printed. -If the routine is called without the second argument, then all units will be set to iunit, except the unit for 'printall', which will be set to 0. ######################################################################## subroutine olo_scale( muscale ) real(kind(1d0)) ,intent(in) :: muscale The type/kind of the input is hard-wired strictly to real(kind(1d0)). Set the default renormalization scale. The input has the units of mass, not mass^2. If this routine is not called, the scale is set to 1. ######################################################################## subroutine olo_onshell( thrs ) real(kind(1d0)) ,intent(in) :: thrs The type/kind of the input is hard-wired strictly to real(kind(1d0)). Set threshold to consider internal masses identical zero and external squared momenta identical zero or equal to internal masses, if this leads to an IR divergent case. For example, if |p1-m1|