Goto Chapter: Top 1 Ind
 [Top of Book]  [Contents]   [Previous Chapter]   [Next Chapter] 

1 Loop integrals
 1.1 Introduction
 1.2 Examples
 1.3 Global variables
 1.4 Properties
 1.5 Attributes

  1.5-1 Dimension

  1.5-2 UnderlyingMatrix

  1.5-3 LoopDiagram

  1.5-4 Components

  1.5-5 LoopMomenta

  1.5-6 ExternalMomenta

  1.5-7 DimensionOfCoefficientsVector

  1.5-8 UnderlyingRing

  1.5-9 RelationsOfExternalMomenta

  1.5-10 RelationsMatrixOfExternalMomenta

  1.5-11 IndependentLorentzInvariants

  1.5-12 Propagators

  1.5-13 Numerators

  1.5-14 ExtraLorentzInvariants

  1.5-15 OriginalIBPGeneratingMatrixOfPropagators

  1.5-16 OriginalIBPGeneratingMatrixOfNumerators

  1.5-17 OriginalTaylorOfPropagators

  1.5-18 MatrixOfMomenta

  1.5-19 IBPGeneratingMatrixOfPropagators

  1.5-20 IBPGeneratingMatrixOfNumerators

  1.5-21 PairOfMatricesOfLoopDiagram

  1.5-22 RingOfExtraLorentzInvariants

  1.5-23 ReductionMatrixOfExtraLorentzInvariants

  1.5-24 RingOfIndependentLorentzInvariants

  1.5-25 ReductionMatrixOfIndependentLorentzInvariants

  1.5-26 RingOfPropagatorsAndNumeratorsAndExtraLorentzInvariants

  1.5-27 RingOfLoopDiagram

  1.5-28 ReductionMatrixOfPropagatorsAndNumeratorsAndExtraLorentzInvariants

  1.5-29 IBPGeneratingMatrixOfPropagatorsInIndependentLorentzInvariants

  1.5-30 IBPGeneratingMatrixOfPropagatorsInPropagators

  1.5-31 IBPGeneratingMatrixOfNumeratorsInPropagators

  1.5-32 PairOfMatricesOfLoopDiagramInIndependentLorentzInvariants

  1.5-33 PairOfMatricesOfLoopDiagramInPropagators

  1.5-34 IBPGeneratingMatrixOfLoopDiagramInPropagators

  1.5-35 MatrixOfIBPRelations

  1.5-36 BasisOfIBPRelations

  1.5-37 MatrixOfSpecialIBPRelations

  1.5-38 BasisOfSpecialIBPRelations

  1.5-39 MatrixOfIBPRelationsInWeylAlgebra

  1.5-40 BasisOfIBPRelationsInWeylAlgebra

  1.5-41 MatrixOfSpecialIBPRelationsInWeylAlgebra

  1.5-42 BasisOfSpecialIBPRelationsInWeylAlgebra

  1.5-43 FieldOfCoefficientsOfLoopDiagramInSingular

  1.5-44 FieldOfCoefficientsOfLoopDiagramInMaple

  1.5-45 FieldOfCoefficientsOfLoopDiagramInHecke

  1.5-46 MatrixOfCoefficientsOfParametricIBPs

  1.5-47 SymanzikPolynomials

  1.5-48 AssociatedWeylAlgebra

  1.5-49 GeneratorsOfScalelessSectors

  1.5-50 GeneratorsOfScalelessSectorsInWeylAlgebra
 1.6 Constructors
 1.7 Operations

1 Loop integrals

1.1 Introduction

This package is an ongoing work to compute master integrals using commutative and noncommutative methods from computational algebraic geometry.

1.2 Examples

1.2-1 Lorentz vectors
gap> LD := LoopDiagram( "l1..2", "p1..2", 4 );
<A loop diagram with loop momenta [ l1, l2 ] & external momenta [ p1, p2 ]>
gap> l1;
l1
gap> Display( l1 );
l1_0,l1_1,l1_2,l1_3
gap> -l1;
<A 4-vector>
gap> Display( -l1 );
-l1_0,-l1_1,-l1_2,-l1_3
gap> l1^0;
1
gap> l1^2;
l1_0^2-l1_1^2-l1_2^2-l1_3^2
gap> l1*l1*p1;
<A 4-vector>
gap> EntriesOfHomalgMatrix( UnderlyingMatrix( l1*l1*p1 ) );
[ l1_0^2*p1_0-l1_1^2*p1_0-l1_2^2*p1_0-l1_3^2*p1_0,
  l1_0^2*p1_1-l1_1^2*p1_1-l1_2^2*p1_1-l1_3^2*p1_1,
  l1_0^2*p1_2-l1_1^2*p1_2-l1_2^2*p1_2-l1_3^2*p1_2,
  l1_0^2*p1_3-l1_1^2*p1_3-l1_2^2*p1_3-l1_3^2*p1_3 ]
gap> l1*l1*p1*p1 = (l1^2)*(p1^2);
true
gap> l1*(l1*p1)*p1 = (l1*(l1*p1))*p1;
true

1.2-2 1 Loop tadpole
gap> LD := LoopDiagram( "k", "" : masses := "m" );
<A loop diagram with loop momenta [ k ] & external momenta [ ] &
 masses [ m ]>
gap> SetRelationsOfExternalMomenta( LD, [ ] );
gap> SetIndependentLorentzInvariants( LD, [ k^2 ] );
gap> SetPropagators( LD, [ k^2 - m^2 ] );
gap> SetNumerators( LD, [ ] );
gap> SetExtraLorentzInvariants( LD, [ ] );
gap> R := RingOfLoopDiagram( LD );
Q[m,d][D1]
gap> ibps := MatrixOfIBPRelations( LD );
<A 1 x 1 matrix over a residue class ring>
gap> ibp := ibps[1,1];
|[ -2*m^2*a1*D1_+d-2*a1 ]|
gap> ViewList( DecomposeInMonomials( ibp ) );
[ [ |[ -2*m^2*a1 ]|, |[ D1_ ]| ],
  [ |[ d-2*a1 ]|, |[ 1 ]| ] ]
gap> Ypol := HomalgRing( ibp );
Q[m,d][a1]<D1,D1_>/( D1*D1_-1 )
gap> ibpws := MatrixOfIBPRelationsInWeylAlgebra( LD );
<A 1 x 1 matrix over an external ring>
gap> ibpw := MatElm( ibpws, 1, 1 );
-2*m^2*A1+d-2*D1*A1-2
gap> W := HomalgRing( ibpw );
Q[m,d][D1]<A1>
gap> E12 := PairOfMatricesOfLoopDiagramInPropagators( LD );
[ <A 1 x 1 matrix over an external ring>,
  <A 1 x 1 matrix over an external ring> ]
gap> Display( E12[1] );
2*m^2+2*D1
gap> Display( E12[2] );
D1
gap> S := SyzygiesOfColumns( E12 );
<A non-zero 1 x 1 matrix over an external ring>
gap> Display( S );
D1
gap> Sibp := IBPRelation( CertainColumns( S, [ 1 ] ), LD );
|[ -2*m^2*a1+2*m^2+d*D1-2*a1*D1+2*D1 ]|
gap> ViewList( DecomposeInMonomials( Sibp ) );
[ [ |[ d-2*a1+2 ]|, |[ D1 ]| ],
  [ |[ -2*m^2*a1+2*m^2 ]|, |[ 1 ]| ] ]
gap> Sibps := MatrixOfSpecialIBPRelations( LD );
<A 1 x 1 matrix over a residue class ring>
gap> x := RightDivide( Sibps, ibps );
<An unevaluated 1 x 1 matrix over a residue class ring>
gap> x * ibps = Sibps;
true
gap> Display( x );
D1

modulo [ D1*D1_-1 ]
gap> y := RightDivide( ibps, Sibps );
<An unevaluated 1 x 1 matrix over a residue class ring>
gap> y * Sibps = ibps;
true
gap> Display( y );
D1_

modulo [ D1*D1_-1 ]
gap> bas := BasisOfIBPRelations( LD );
<A non-zero 1 x 1 matrix over a residue class ring>
gap> Sbas := BasisOfSpecialIBPRelations( LD );
<A non-zero 1 x 1 matrix over a residue class ring>
gap> Sbas = bas;
true
gap> ExportVariables( Ypol );
[ |[ m ]|, |[ d ]|, |[ a1 ]|, |[ D1 ]|, |[ D1_ ]| ]
gap> r := 2 * m^2 * ( a1 + 4 - 2 ) * D1_^(4-1);
|[ 2*m^2*a1*D1_^3+4*m^2*D1_^3 ]|
gap> n := DecideZero( r, bas );
|[ d*D1_^2-2*a1*D1_^2-4*D1_^2 ]|
gap> DecideZero( n - r, bas );
|[ 0 ]|
gap> s := D1_^(4 - 2) * ( d - 2 * a1 );
|[ d*D1_^2-2*a1*D1_^2-4*D1_^2 ]|
gap> s = n;
true
gap> t := ( d - 2 * a1 - 2 * ( 4 - 2 ) ) * D1_^(4-2);
|[ d*D1_^2-2*a1*D1_^2-4*D1_^2 ]|
gap> t = n;
true
gap> w := ( d - 2 * 1 - 2 * ( 4 - 2 ) ) * D1_^(4-2);
|[ d*D1_^2-6*D1_^2 ]|
gap> DecideZero( w, bas );
|[ d*D1_^2-6*D1_^2 ]|

1.2-3 1 Loop bubble
gap> LD := LoopDiagram( "l1", "k1" );
<A loop diagram with loop momenta [ l1 ] & external momenta [ k1 ]>
gap> s := k1*k1;
k1^2
gap> SetAbbreviation( s, "s" );
gap> rel := [ ];;
gap> SetRelationsOfExternalMomenta( LD, rel );
gap> SetIndependentLorentzInvariants( LD, [ ] );
gap> SetPropagators( LD, -[ l1^2, ( l1 + k1 )^2 ] );
gap> SetNumerators( LD, -[ ] );
gap> SetExtraLorentzInvariants( LD, [ s ] );
gap> R := RingOfLoopDiagram( LD );
Q[d,s][D1,D2]
gap> ibps := MatrixOfIBPRelations( LD );
<A 2 x 1 matrix over a residue class ring>
gap> ibp1 := ibps[1,1];
|[ -s*a2*D2_-a2*D1*D2_+d-2*a1-a2 ]|
gap> ViewList( DecomposeInMonomials( ibp1 ) );
[ [ |[ -a2 ]|, |[ D1*D2_ ]| ],
  [ |[ -s*a2 ]|, |[ D2_ ]| ],
  [ |[ d-2*a1-a2 ]|, |[ 1 ]| ] ]
gap> Ypol := HomalgRing( ibp1 );
Q[d,s][a1,a2]<D1,D1_,D2,D2_>/( D2*D2_-1, D1*D1_-1 )
gap> bas := BasisOfRows( ibps );
<A non-zero 4 x 1 matrix over a residue class ring>
gap> gen := GeneratorsOfScalelessSectors( LD );
<A 1 x 2 matrix over an external ring>
gap> Display( gen );
D2,D1
gap> gen2 := GeneratorsOfScalelessSectors( LD, [ 2, 2 ] );
<An unevaluated 1 x 2 matrix over an external ring>
gap> Display( gen2 );
D1*D2^2,D1^2*D2
gap> prel2 := ColumnReversedMatrixOfCoefficientsOfParametricIBPs( LD, 2 );
[ <A non-zero 5 x 6 matrix over an external ring>,
  [ D1_^2, D1_*D2_, D2_^2, D1_, D2_, 1 ] ]
gap> EntriesOfHomalgMatrix( prel2[1][4] );
[ 0, 0, 0,
  (d*s*a1-2*s*a1^2-2*s*a1), 0,
  (-d^2+3*d*a1+3*d*a2+d-2*a1^2-4*a1*a2-2*a1-2*a2^2-2*a2) ]
gap> EntriesOfHomalgMatrix( prel2[1][5] );
[ 0, 0, 0,
  0, (d*s*a2-2*s*a2^2-2*s*a2),
  (-d^2+3*d*a1+3*d*a2+d-2*a1^2-4*a1*a2-2*a1-2*a2^2-2*a2) ]

1.2-4 1 Loop massive on-shell
gap> LD := LoopDiagram( "k", "q" : masses := "m" );
<A loop diagram with loop momenta [ k ] & external momenta [ q ] &
 masses [ m ]>
gap> rel := [ q^2 - m^2 ];
[ q^2+(-m^2) ]
gap> SetRelationsOfExternalMomenta( LD, rel );
gap> SetIndependentLorentzInvariants( LD, [  ] );
gap> SetPropagators( LD, -[ k^2, k^2 + 2*q*k ] );
gap> SetNumerators( LD, -[ ] );
gap> SetExtraLorentzInvariants( LD, [ ] );
gap> R := RingOfLoopDiagram( LD );
Q[m,d][D1,D2]
gap> ibps := MatrixOfIBPRelations( LD );
<A 2 x 1 matrix over a residue class ring>
gap> Display( ibps );
-a2*D1*D2_+d-2*a1-a2,
2*m^2*a2*D2_-a1*D1_*D2+a2*D1*D2_+a1-a2

modulo [ D2*D2_-1, D1*D1_-1 ]
gap> ibp1 := ibps[1,1];
|[ -a2*D1*D2_+d-2*a1-a2 ]|
gap> ViewList( DecomposeInMonomials( ibp1 ) );
[ [ |[ -a2 ]|, |[ D1*D2_ ]| ],
  [ |[ d-2*a1-a2 ]|, |[ 1 ]| ] ]
gap> ibp2 := ibps[2,1];
|[ 2*m^2*a2*D2_-a1*D1_*D2+a2*D1*D2_+a1-a2 ]|
gap> ViewList( DecomposeInMonomials( ibp2 ) );
[ [ |[ -a1 ]|, |[ D1_*D2 ]| ],
  [ |[ a2 ]|, |[ D1*D2_ ]| ],
  [ |[ 2*m^2*a2 ]|, |[ D2_ ]| ],
  [ |[ a1-a2 ]|, |[ 1 ]| ] ]
gap> ibp2_red := DecideZero( ibp2, ibps{[ 1 ]} );
|[ 2*m^2*a2*D2_-a1*D1_*D2+d-a1-2*a2 ]|
gap> ViewList( DecomposeInMonomials( ibp2_red ) );
[ [ |[ -a1 ]|, |[ D1_*D2 ]| ],
  [ |[ 2*m^2*a2 ]|, |[ D2_ ]| ],
  [ |[ d-a1-2*a2 ]|, |[ 1 ]| ] ]
gap> Ypol := HomalgRing( ibp1 );
Q[m,d][a1,a2]<D1,D1_,D2,D2_>/( D2*D2_-1, D1*D1_-1 )
gap> gen := GeneratorsOfScalelessSectors( LD );
<A 1 x 1 matrix over an external ring>
gap> Display( gen );
D2

1.2-5 1 Loop massive bubble
gap> LD := LoopDiagram( "l1", "k1" : masses := "m" );
<A loop diagram with loop momenta [ l1 ] & external momenta [ k1 ] &
 masses [ m ]>
gap> s := k1*k1;
k1^2
gap> SetAbbreviation( s, "s" );
gap> rel := [ ];;
gap> SetRelationsOfExternalMomenta( LD, rel );
gap> SetIndependentLorentzInvariants( LD, [ ] );
gap> SetPropagators( LD, -[ l1^2 - m^2, ( l1 + k1 )^2 ] );
gap> SetNumerators( LD, -[ ] );
gap> SetExtraLorentzInvariants( LD, [ s ] );
gap> R := RingOfLoopDiagram( LD );
Q[m,d,s][D1,D2]
gap> ibps := MatrixOfIBPRelations( LD );
<A 2 x 1 matrix over a residue class ring>
gap> ibp1 := ibps[1,1];
|[ 2*m^2*a1*D1_+m^2*a2*D2_-s*a2*D2_-a2*D1*D2_+d-2*a1-a2 ]|
gap> ViewList( DecomposeInMonomials( ibp1 ) );
[ [ |[ -a2 ]|, |[ D1*D2_ ]| ],
  [ |[ 2*m^2*a1 ]|, |[ D1_ ]| ],
  [ |[ m^2*a2-s*a2 ]|, |[ D2_ ]| ],
  [ |[ d-2*a1-a2 ]|, |[ 1 ]| ] ]
gap> Ypol := HomalgRing( ibp1 );
Q[m,d,s][a1,a2]<D1,D1_,D2,D2_>/( D2*D2_-1, D1*D1_-1 )
gap> bas := BasisOfRows( ibps );
<A non-zero 7 x 1 matrix over a residue class ring>
gap> gen := GeneratorsOfScalelessSectors( LD );
<A 1 x 1 matrix over an external ring>
gap> Display( gen );
D1
gap> gen2 := GeneratorsOfScalelessSectors( LD, [ 2, 2 ] );
<An unevaluated 1 x 1 matrix over an external ring>
gap> Display( gen2 );
D1^2*D2

1.2-6 1 Loop double massive bubble
gap> LD := LoopDiagram( "l1", "k1" : masses := "m" );
<A loop diagram with loop momenta [ l1 ] & external momenta [ k1 ] &
 masses [ m ]>
gap> s := k1*k1;
k1^2
gap> SetAbbreviation( s, "s" );
gap> rel := [ ];;
gap> SetRelationsOfExternalMomenta( LD, rel );
gap> SetIndependentLorentzInvariants( LD, [ ] );
gap> SetPropagators( LD, -[ l1^2 + m^2, ( l1 + k1 )^2 + m^2 ] );
gap> SetNumerators( LD, -[ ] );
gap> SetExtraLorentzInvariants( LD, [ s ] );
gap> R := RingOfLoopDiagram( LD );
Q[m,d,s][D1,D2]
gap> ibps := MatrixOfIBPRelations( LD );
<A 2 x 1 matrix over a residue class ring>
gap> ibp1 := ibps[1,1];
|[ -2*m^2*a1*D1_-2*m^2*a2*D2_-s*a2*D2_-a2*D1*D2_+d-2*a1-a2 ]|
gap> ViewList( DecomposeInMonomials( ibp1 ) );
[ [ |[ -a2 ]|, |[ D1*D2_ ]| ],
  [ |[ -2*m^2*a1 ]|, |[ D1_ ]| ],
  [ |[ -2*m^2*a2-s*a2 ]|, |[ D2_ ]| ],
  [ |[ d-2*a1-a2 ]|, |[ 1 ]| ] ]
gap> Ypol := HomalgRing( ibp1 );
Q[m,d,s][a1,a2]<D1,D1_,D2,D2_>/( D2*D2_-1, D1*D1_-1 )
gap> bas := BasisOfRows( ibps );
<A non-zero 6 x 1 matrix over a residue class ring>
gap> gen := GeneratorsOfScalelessSectors( LD );
<A 1 x 1 matrix over an external ring>
gap> Display( gen );
D1*D2
gap> gen2 := GeneratorsOfScalelessSectors( LD, [ 2, 2 ] );
<An unevaluated 1 x 1 matrix over an external ring>
gap> Display( gen2 );
D1^2*D2^2

1.2-7 1 Loop box
gap> LD := LoopDiagram( "l1", "k1..2,k4" );
<A loop diagram with loop momenta [ l1 ] & external momenta [ k1, k2, k4 ]>
gap> s12 := 2*k1*k2;
2*k1*k2
gap> SetAbbreviation( s12, "s12" );
gap> s14 := 2*k1*k4;
2*k1*k4
gap> SetAbbreviation( s14, "s14" );
gap> rel1 := List( ExternalMomenta( LD ), k -> k^2 );
[ k1^2, k2^2, k4^2 ]
gap> rel2 := [ (k1+k2+k4)^2 ];;
gap> SetRelationsOfExternalMomenta( LD, Concatenation( rel1, rel2 ) );
gap> SetIndependentLorentzInvariants( LD,
>         [ l1^2, l1*k1, l1*k2, l1*k4, s12, s14 ] );
gap> SetPropagators( LD, -[ l1^2, (l1-k1)^2, (l1-k1-k2)^2, (l1+k4)^2 ] );
gap> SetNumerators( LD, -[ ] );
gap> SetExtraLorentzInvariants( LD, [ s12, s14 ] );
gap> R := RingOfLoopDiagram( LD );
Q[d,s12,s14][D1,D2,D3,D4]
gap> ibps := MatrixOfIBPRelations( LD );
<A 4 x 1 matrix over a residue class ring>
gap> ibp1 := ibps[1,1];
|[ -a2*D1*D2_-s12*a3*D3_-a3*D1*D3_-a4*D1*D4_+d-2*a1-a2-a3-a4 ]|
gap> ViewList( DecomposeInMonomials( ibp1 ) );
[ [ |[ -a2 ]|, |[ D1*D2_ ]| ],
  [ |[ -a3 ]|, |[ D1*D3_ ]| ],
  [ |[ -a4 ]|, |[ D1*D4_ ]| ],
  [ |[ -s12*a3 ]|, |[ D3_ ]| ],
  [ |[ d-2*a1-a2-a3-a4 ]|, |[ 1 ]| ] ]
gap> Ypol := HomalgRing( ibp1 );
Q[d,s12,s14][a1,a2,a3,a4]<D1,D1_,D2,D2_,D3,D3_,D4,D4_>/( D4*D4_-1, D3*D3_-1,\
  D2*D2_-1, D1*D1_-1 )
gap> E12 := PairOfMatricesOfLoopDiagramInPropagators( LD );
[ <A 4 x 4 matrix over an external ring>,
  <A 4 x 4 matrix over an external ring> ]
gap> Display( E12[1] );
2*D1,     D1-D2,     -s12+D2-D3,-D1+D4,
D1+D2,    D1-D2,     D2-D3,     s14-D1+D4,
s12+D1+D3,s12+D1-D2, D2-D3,     -s12-D1+D4,
D1+D4,    -s14+D1-D2,s14+D2-D3, -D1+D4
gap> S := SyzygiesOfColumns( E12 );
<A non-zero 4 x 8 matrix over an external ring>
gap> Sred := ReducedBasisOfColumnModule( BasisOfColumnModule( S ) );
<A non-zero 4 x 6 matrix over an external ring>
gap> Display( Sred );
D2-D4,D1-D3,s12*D4+2*D3*D4-2*D4^2,     s14*D3-2*D3^2+2*D3*D4,-D3*D4^2,D3^2*D4,
D4,   -D1,  -s12*D4+D2*D4-D3*D4+2*D4^2,-D1*D3-D3*D4,         D3*D4^2, 0,      
0,    -D1,  D1*D4+D2*D4,               -2*D1*D3,             0,       D1*D3*D4,
D2,   0,    2*D2*D4,                   -D1*D3-D2*D3,         D2*D3*D4,0       
gap> Display( EntriesOfHomalgMatrixAsListList( CertainColumns( Sred, [1 .. 3] ) ) );
[ [ D2-D4, D1-D3, s12*D4+2*D3*D4-2*D4^2 ],
  [ D4, -D1, -s12*D4+D2*D4-D3*D4+2*D4^2 ],
  [ 0, -D1, D1*D4+D2*D4 ],
  [ D2, 0, 2*D2*D4 ] ]
gap> Sibp1 := IBPRelation( CertainColumns( S, [ 1 ] ), LD );
|[ -s14*a2+s14*a4+d*D2-a1*D2-a2*D2-a3*D2-a4*D2-d*D4+a1*D4+a2*D4+a3*D4+a4*D4 ]|
gap> ViewList( DecomposeInMonomials( Sibp1 ) );
[ [ |[ d-a1-a2-a3-a4 ]|, |[ D2 ]| ],
  [ |[ -d+a1+a2+a3+a4 ]|, |[ D4 ]| ],
  [ |[ -s14*a2+s14*a4 ]|, |[ 1 ]| ] ]
gap> sibp1 := IBPRelation( CertainColumns( S, [ 1 ] ), LD, [ 1, 1, 1, 1 ] );
|[ d*D2-d*D4-4*D2+4*D4 ]|
gap> ViewList( DecomposeInMonomials( sibp1 ) );
[ [ |[ d-4 ]|, |[ D2 ]| ],
  [ |[ -d+4 ]|, |[ D4 ]| ] ]
gap> Sibp2 := IBPRelation( CertainColumns( S, [ 2 ] ), LD );
|[ -s12*a1+s12*a3+d*D1-a1*D1-a2*D1-a3*D1-a4*D1-d*D3+a1*D3+a2*D3+a3*D3+a4*D3 ]|
gap> ViewList( DecomposeInMonomials( Sibp2 ) );
[ [ |[ d-a1-a2-a3-a4 ]|, |[ D1 ]| ],
  [ |[ -d+a1+a2+a3+a4 ]|, |[ D3 ]| ],
  [ |[ -s12*a1+s12*a3 ]|, |[ 1 ]| ] ]
gap> sibp2 := IBPRelation( CertainColumns( S, [ 2 ] ), LD, [ 1, 1, 1, 1 ] );
|[ d*D1-d*D3-4*D1+4*D3 ]|
gap> ViewList( DecomposeInMonomials( sibp2 ) );
[ [ |[ d-4 ]|, |[ D1 ]| ],
  [ |[ -d+4 ]|, |[ D3 ]| ] ]
gap> Sibp3 := IBPRelation( CertainColumns( S, [ 3 ] ), LD );;
gap> ViewList( DecomposeInMonomials( Sibp3 ) );
[ [ |[ 2*d-2*a1-2*a2-2*a3-2*a4+2 ]|, |[ D3*D4 ]| ],
  [ |[ -2*d+2*a1+2*a2+2*a3+2*a4-2 ]|, |[ D4^2 ]| ],
  [ |[ -s14*a4+s14 ]|, |[ D1 ]| ],
  [ |[ -s12*a4+s12 ]|, |[ D2 ]| ],
  [ |[ -s14*a4+s14 ]|, |[ D3 ]| ],
  [ |[ d*s12-2*s12*a2-2*s14*a2-2*s12*a3-s12*a4+2*s14*a4+s12-2*s14 ]|, |[ D4 ]| ],
  [ |[ -s12*s14*a4+s12*s14 ]|, |[ 1 ]| ] ]
gap> sibp3 := IBPRelation( CertainColumns( S, [ 3 ] ), LD, [ 1, 1, 1, 1 ] );
|[ d*s12*D4+2*d*D3*D4-2*d*D4^2-4*s12*D4-2*s14*D4-6*D3*D4+6*D4^2 ]|
gap> ViewList( DecomposeInMonomials( sibp3 ) );
[ [ |[ 2*d-6 ]|, |[ D3*D4 ]| ],
  [ |[ -2*d+6 ]|, |[ D4^2 ]| ],
  [ |[ d*s12-4*s12-2*s14 ]|, |[ D4 ]| ] ]
gap> bas := BasisOfIBPRelations( LD );
<A non-zero 28 x 1 matrix over a residue class ring>
gap> Sbas := BasisOfSpecialIBPRelations( LD );
<A non-zero 28 x 1 matrix over a residue class ring>
gap> bas = Sbas;
true
gap> SymanzikPolynomials( LD );
[ z1+z2+z3+z4, -s12*z1*z3-s14*z2*z4 ]
gap> SymanzikPolynomials( LD, [ 1, 2, 3, 4 ] );
[ z1+z2+z3+z4, -s12*z1*z3-s14*z2*z4 ]
gap> SymanzikPolynomials( LD, [ 1, 2, 3 ] );
[ z1+z2+z3, -s12*z1*z3 ]
gap> SymanzikPolynomials( LD, [ 1, 2 ] );
[ z1+z2, 0 ]
gap> SymanzikPolynomials( LD, [ 1 ] );
[ z1, 0 ]
gap> SymanzikPolynomials( LD, [ ] );
[ 0, 0 ]
gap> gen := GeneratorsOfScalelessSectors( LD );
<A 1 x 4 matrix over an external ring>
gap> Display( gen );
D3*D4,D1*D4,D2*D3,D1*D2
gap> gen2 := GeneratorsOfScalelessSectors( LD, [ 2, 2, 2, 2 ] );
<An unevaluated 1 x 4 matrix over an external ring>
gap> Display( gen2 );
D1*D2*D3^2*D4^2,D1^2*D2*D3*D4^2,D1*D2^2*D3^2*D4,D1^2*D2^2*D3*D4
gap> prel2 := ColumnReversedMatrixOfCoefficientsOfParametricIBPs( LD, 2 );
[ <A non-zero 7 x 11 matrix over an external ring>,
  [ D1_*D3_, D1_*D4_, D2_*D4_, D3_*D4_, D1_, D2_, D3_, D4_, 1, D1, D2 ] ]

1.2-8 Kite
gap> LD := LoopDiagram( "l1..2", "p", 4 );
<A loop diagram with loop momenta [ l1, l2 ] & external momenta [ p ]>
gap> s := p^2;;
gap> SetAbbreviation( s, "s" );
gap> SetRelationsOfExternalMomenta( LD, [ ] );
gap> SetIndependentLorentzInvariants( LD, [ l1^2, l2^2, l1*l2, l1*p, l2*p, s ] );
gap> SetPropagators( LD, -[ l1^2, (l1+p)^2, l2^2, (l2-p)^2, (l1+l2)^2 ] );
gap> SetNumerators( LD, -[ ] );
gap> SetExtraLorentzInvariants( LD, [ s ] );
gap> e12 := PairOfMatricesOfLoopDiagramInIndependentLorentzInvariants( LD );
[ <An unevaluated non-zero 5 x 6 matrix over an external ring>,
  <An unevaluated non-zero 5 x 5 matrix over an external ring> ]
gap> Display( e12[1] );
-2*l1l1,-2*l1l2,-2*l1p,      0,      0,      0,
_[2,1], _[2,2], -2*l1p-2*s,  0,      0,      0,
0,      0,      0,           -2*l1l2,-2*l2l2,-2*l2p,
0,      0,      0,           _[4,4], _[4,5], -2*l2p+2*s,
_[5,1], _[5,2], -2*l1p-2*l2p,_[5,4], _[5,5], -2*l1p-2*l2p
gap> Display( e12[2] );
-l1l1,0,            0,    0,            0,
0,    -l1l1-2*l1p-s,0,    0,            0,
0,    0,            -l2l2,0,            0,
0,    0,            0,    -l2l2+2*l2p-s,0,
0,    0,            0,    0,            -l1l1-l2l2-2*l1l2
gap> E12 := PairOfMatricesOfLoopDiagramInPropagators( LD );
[ <A 5 x 6 matrix over an external ring>,
  <A 5 x 5 matrix over an external ring> ]
gap> Display( E12[1] );
2*D1,    -D1-D3+D5,  s-D1+D2,     0,          0,        0,
s+D1+D2, -s-D1-D4+D5,-s-D1+D2,    0,          0,        0,
0,       0,          0,           -D1-D3+D5,  2*D3,     -s+D3-D4,
0,       0,          0,           -s-D2-D3+D5,s+D3+D4,  s+D3-D4,
D1-D3+D5,-D1+D3+D5,  -D1+D2+D3-D4,D1-D3+D5,   -D1+D3+D5,-D1+D2+D3-D4
gap> Display( E12[2] );
D1,0, 0, 0, 0,
0, D2,0, 0, 0,
0, 0, D3,0, 0,
0, 0, 0, D4,0,
0, 0, 0, 0, D5
gap> S := SyzygiesOfColumns( E12 );
<A non-zero 6 x 11 matrix over an external ring>
gap> Display( EntriesOfHomalgMatrix( CertainColumns( S, [ 1 ] ) ) );
[ D1-D2, 0, D1, 0, D3-D4, -D3 ]
gap> Display( EntriesOfHomalgMatrix( CertainColumns( S, [ 2 ] ) ) );
[ 2*D2, 0, 0, s+D3+D4, s+D1+D2+2*D4-2*D5, -D1+D3+D5 ]
gap> Display( EntriesOfHomalgMatrix( CertainColumns( S, [ 3 ] ) ) );
[ s+2*D2+D3+D4-2*D5, s+D1+D2, -D1+D3-D5, 0, 2*D4, 0 ]
gap> Display( EntriesOfHomalgMatrix( CertainColumns( S, [ 4 ] ) ) );
[ 0, 0, 0, 0, D4*D5, 0 ]
gap> Sibp1 := IBPRelation( CertainColumns( S, [ 1 ] ), LD );
|[ -s*a1+s*a2-s*a3+s*a4+d*D1-a1*D1-a2*D1-a5*D1-d*D2+a1*D2+a2*D2+a5*D2\
   +d*D3-a3*D3-a4*D3-a5*D3-d*D4+a3*D4+a4*D4+a5*D4 ]|
gap> ViewList( DecomposeInMonomials( Sibp1 ) );
[ [ |[ d-a1-a2-a5 ]|, |[ D1 ]| ],
  [ |[ -d+a1+a2+a5 ]|, |[ D2 ]| ],
  [ |[ d-a3-a4-a5 ]|, |[ D3 ]| ],
  [ |[ -d+a3+a4+a5 ]|, |[ D4 ]| ],
  [ |[ -s*a1+s*a2-s*a3+s*a4 ]|, |[ 1 ]| ] ]
gap> sibp1 := IBPRelation( CertainColumns( S, [ 1 ] ), LD, [ 1, 1, 1, 1, 1 ] );
|[ d*D1-d*D2+d*D3-d*D4-3*D1+3*D2-3*D3+3*D4 ]|
gap> ViewList( DecomposeInMonomials( sibp1 ) );
[ [ |[ d-3 ]|, |[ D1 ]| ],
  [ |[ -d+3 ]|, |[ D2 ]| ],
  [ |[ d-3 ]|, |[ D3 ]| ],
  [ |[ -d+3 ]|, |[ D4 ]| ] ]
gap> Sibp2 := IBPRelation( CertainColumns( S, [ 2 ] ), LD );
|[ d*s-2*s*a2-2*s*a4-2*s*a5+d*D1-2*a2*D1-2*a4*D1-2*a5*D1\
   +3*d*D2-4*a1*D2-2*a2*D2-2*a3*D2-4*a5*D2+2*d*D4-2*a3*D4-2*a4*D4-2*a5*D4\
   -2*d*D5+2*a3*D5+2*a4*D5+2*a5*D5+2*s+2*D1+2*D2 ]|
gap> ViewList( DecomposeInMonomials( Sibp2 ) );
[ [ |[ d-2*a2-2*a4-2*a5+2 ]|, |[ D1 ]| ],
  [ |[ 3*d-4*a1-2*a2-2*a3-4*a5+2 ]|, |[ D2 ]| ],
  [ |[ 2*d-2*a3-2*a4-2*a5 ]|, |[ D4 ]| ],
  [ |[ -2*d+2*a3+2*a4+2*a5 ]|, |[ D5 ]| ],
  [ |[ d*s-2*s*a2-2*s*a4-2*s*a5+2*s ]|, |[ 1 ]| ] ]
gap> sibp2 := IBPRelation( CertainColumns( S, [ 2 ] ), LD, [ 1, 1, 1, 1, 1 ] );
|[ d*s+d*D1+3*d*D2+2*d*D4-2*d*D5-4*s-4*D1-10*D2-6*D4+6*D5 ]|
gap> ViewList( DecomposeInMonomials( sibp2 ) );
[ [ |[ d-4 ]|, |[ D1 ]| ],
  [ |[ 3*d-10 ]|, |[ D2 ]| ],
  [ |[ 2*d-6 ]|, |[ D4 ]| ],
  [ |[ -2*d+6 ]|, |[ D5 ]| ],
  [ |[ d*s-4*s ]|, |[ 1 ]| ] ]
gap> Sibp3 := IBPRelation( CertainColumns( S, [ 3 ] ), LD );
|[ d*s-2*s*a2-2*s*a4-2*s*a5+2*d*D2-2*a1*D2-2*a2*D2-2*a5*D2\
   +d*D3-2*a2*D3-2*a4*D3-2*a5*D3+3*d*D4-2*a1*D4-4*a3*D4-2*a4*D4-4*a5*D4\
   -2*d*D5+2*a1*D5+2*a2*D5+2*a5*D5+2*s+2*D3+2*D4 ]|
gap> ViewList( DecomposeInMonomials( Sibp3 ) );
[ [ |[ 2*d-2*a1-2*a2-2*a5 ]|, |[ D2 ]| ],
  [ |[ d-2*a2-2*a4-2*a5+2 ]|, |[ D3 ]| ],
  [ |[ 3*d-2*a1-4*a3-2*a4-4*a5+2 ]|, |[ D4 ]| ],
  [ |[ -2*d+2*a1+2*a2+2*a5 ]|, |[ D5 ]| ],
  [ |[ d*s-2*s*a2-2*s*a4-2*s*a5+2*s ]|, |[ 1 ]| ] ]
gap> sibp3 := IBPRelation( CertainColumns( S, [ 3 ] ), LD, [ 1, 1, 1, 1, 1 ] );
|[ d*s+2*d*D2+d*D3+3*d*D4-2*d*D5-4*s-6*D2-4*D3-10*D4+6*D5 ]|
gap> ViewList( DecomposeInMonomials( sibp3 ) );
[ [ |[ 2*d-6 ]|, |[ D2 ]| ],
  [ |[ d-4 ]|, |[ D3 ]| ],
  [ |[ 3*d-10 ]|, |[ D4 ]| ],
  [ |[ -2*d+6 ]|, |[ D5 ]| ],
  [ |[ d*s-4*s ]|, |[ 1 ]| ] ]
gap> Sibp4 := IBPRelation( CertainColumns( S, [ 4 ] ), LD );
|[ a5*D1*D4-a5*D3*D4-s*a4*D5-a4*D3*D5+d*D4*D5-2*a3*D4*D5-a4*D4*D5-a5*D4*D5\
   -D1*D4+D3*D4+s*D5+D3*D5+2*D4*D5 ]|
gap> ViewList( DecomposeInMonomials( Sibp4 ) );
[ [ |[ a5-1 ]|, |[ D1*D4 ]| ],
  [ |[ -a5+1 ]|, |[ D3*D4 ]| ],
  [ |[ -a4+1 ]|, |[ D3*D5 ]| ],
  [ |[ d-2*a3-a4-a5+2 ]|, |[ D4*D5 ]| ],
  [ |[ -s*a4+s ]|, |[ D5 ]| ] ]
gap> sibp4 := IBPRelation( CertainColumns( S, [ 4 ] ), LD, [ 1, 1, 1, 1, 1 ] );
|[ d*D4*D5-2*D4*D5 ]|
gap> ViewList( DecomposeInMonomials( sibp4 ) );
[ [ |[ d-2 ]|, |[ D4*D5 ]| ] ]
gap> LD := LoopDiagram( "l1..2", "p", 4 : abbreviation := false );
<A loop diagram with loop momenta [ l1, l2 ] & external momenta [ p ]>
gap> s := p^2;;
gap> SetAbbreviation( s, "s" );
gap> SetRelationsOfExternalMomenta( LD, [ ] );
gap> SetIndependentLorentzInvariants( LD, [ l1^2, l2^2, l1*l2, l1*p, l2*p, s ] );
gap> SetPropagators( LD, -[ l1^2, (l1+p)^2, l2^2, (l2-p)^2, (l1+l2)^2 ] );
gap> SetNumerators( LD, -[ ] );
gap> SetExtraLorentzInvariants( LD, [ s ] );
gap> e12 := PairOfMatricesOfLoopDiagramInIndependentLorentzInvariants( LD );
[ <An unevaluated non-zero 5 x 6 matrix over an external ring>,
  <An unevaluated non-zero 5 x 5 matrix over an external ring> ]
gap> Display( e12[1] );
-2*x1,     -2*x3,     -2*x4,     0,         0,         0,
-2*x1-2*x4,-2*x3-2*x5,-2*x4-2*x6,0,         0,         0,
0,         0,         0,         -2*x3,     -2*x2,     -2*x5,
0,         0,         0,         -2*x3+2*x4,-2*x2+2*x5,-2*x5+2*x6,
-2*x1-2*x3,-2*x2-2*x3,-2*x4-2*x5,-2*x1-2*x3,-2*x2-2*x3,-2*x4-2*x5
gap> Display( e12[2] );
-x1,0,          0,  0,          0,
0,  -x1-2*x4-x6,0,  0,          0,
0,  0,          -x2,0,          0,
0,  0,          0,  -x2+2*x5-x6,0,
0,  0,          0,  0,          -x1-x2-2*x3

1.2-9 2 Loop nonplanar ladder
gap> LD := LoopDiagram( "l1..2", "p1..2" );
<A loop diagram with loop momenta [ l1, l2 ] & external momenta [ p1, p2 ]>
gap> s := 2*p1*p2;;
gap> SetAbbreviation( s, "s" );
gap> rel := List( ExternalMomenta( LD ), p -> p^2 );;
gap> SetRelationsOfExternalMomenta( LD, rel );
gap> SetIndependentLorentzInvariants( LD,
>         [ l1^2, l2^2, l1*l2,
>           l1*p1, l2*p1, l1*p2, l2*p2,
>           s ] );
gap> SetPropagators( LD,
>         -[ (l1-p1)^2, (l1+p2)^2, l2^2, (l2-p1)^2,
>            (l1-l2)^2, (l1-l2+p2)^2 ] );
gap> SetNumerators( LD, -[ l1^2 ] );
gap> SetExtraLorentzInvariants( LD, [ s ] );
gap> E12 := PairOfMatricesOfLoopDiagramInPropagators( LD );
[ <A 6 x 8 matrix over an external ring>,
  <A 6 x 6 matrix over an external ring> ]
gap> Display( E12[1] );
D1+N7,    D4-D5+N7, -D1+N7,  s+D2-N7,0,        0,        0,     0,
D2+N7,    D2+D3-D6, -s-D1+N7,D2-N7,  0,        0,        0,     0,
0,        0,        0,       0,      D3-D5+N7, 2*D3,     D3-D4, _[3,8],
0,        0,        0,       0,      D1+D3-D5, D3+D4,    D3-D4, _[4,8],
-D3+D5+N7,-D3-D5+N7,_[5,3],  -D5+D6, D3-D5-N7, D3+D5-N7, _[5,7],D5-D6,
D2-D3+D5, D2-D3-D6, _[6,3],  -D5+D6, -D2+D3-D5,-D2+D3+D6,_[6,7],D5-D6
gap> S := SyzygiesOfColumns( E12 );
<A non-zero 8 x 48 matrix over an external ring>
gap> Display( EntriesOfHomalgMatrix( CertainColumns( S, [ 1 ] ) ) );
[ D1-D2, 0, D2, D1,
 -D3+D4, D1-D2+D3-D4-D5+D6, D2-D6, D4 ]
gap> Display( EntriesOfHomalgMatrix( CertainColumns( S, [ 2 ] ) ) );
[ s+2*D2-2*N7, 0, -D2+N7, -D1-N7,
  0, s+2*D2+2*D5-2*D6-2*N7, -D2-D5+D6+N7, -D3-D4 ]
gap> Display( EntriesOfHomalgMatrix( CertainColumns( S, [ 3 ] ) ) );
[ 0, 0, 0, 0,
  -D3*D6+D4*D6, D1*D6+2*D3*D6-2*D4*D6-D6*N7, -D3*D6-D5*D6+D6*N7, 0 ]

1.2-10 2 Loop cusp anomalous dimensions
gap> LD := LoopDiagram( "k1, k2", "v1, v2" );
<A loop diagram with loop momenta [ k1, k2 ] & external momenta [ v1, v2 ]>
gap> cos:= v1*v2;;
gap> SetAbbreviation( cos, "cos" );
gap> rel := [ v1^2 - 1, v2^2 - 1 ];;
gap> SetRelationsOfExternalMomenta( LD, rel );
gap> SetIndependentLorentzInvariants( LD,
>         [ k1^2, k1*k2, k1*v1, k1*v2, k2^2, k2*v1, k2*v2, cos ] );
gap> SetPropagators( LD,
>         -[ 2*k2*v1 - 1, 2*k2*v2 - 1, (k1 - k2)^2,
>            2*k1*v1 - 1, 2*k1*v2 - 1, k1^2 ] );
gap> SetNumerators( LD, -[ k2^2 ] );
gap> SetExtraLorentzInvariants( LD, [ cos ] );
gap> E12 := PairOfMatricesOfLoopDiagramInPropagators( LD );
[ <A 6 x 8 matrix over an external ring>,
  <A 6 x 6 matrix over an external ring> ]
gap> Display( E12[1] );
0,       0,        0,     0,     D4-1,     D1-1,    -2,    -2*cos,
0,       0,        0,     0,     D5-1,     D2-1,    -2*cos,-2,
D3+D6-N7,-D3+D6-N7,-D1+D4,-D2+D5,-D3-D6+N7,D3-D6+N7,D1-D4, D2-D5,
D4-1,    D1-1,     -2,    -2*cos,0,        0,       0,     0,
D5-1,    D2-1,     -2*cos,-2,    0,        0,       0,     0,
2*D6,    -D3+D6+N7,D4-1,  D5-1,  0,        0,       0,     0
gap> S := SyzygiesOfColumns( E12 );
<A non-zero 8 x 65 matrix over an external ring>
gap> EntriesOfHomalgMatrix( CertainColumns( S, [ 1 ] ) );
[ -2*D4*D6-2*D5*D6, 0, D5*D6, D4*D6,
  2*cos*D6-2*D6, -2*cos*D6-2*D4*D6-2*D5*D6+2*D6, D5*D6, D4*D6 ]
gap> EntriesOfHomalgMatrix( CertainColumns( S, [ 2 ] ) );
[ 2*cos*D6+2*D6, 0, -D6, -D6, 0, 2*cos*D6+2*D6, -D6, -D6 ]
gap> EntriesOfHomalgMatrix( CertainColumns( S, [ 3 ] ) );
[ D2*D4+2*D3*D4+D1*D5+2*D3*D5+2*D4*D5-2*D4*N7-2*D5*N7-D4-D5, -2*D4*D5+D4+D5,
  -D3*D5+D5*N7, -D3*D4+D4*N7, 2*cos*N7+D1+D2-2*N7,
  D2*D4+D1*D5-2*cos*N7-2*D4*N7-2*D5*N7-D1-D2+2*N7, D5*N7, D4*N7 ]
gap> LD := LoopDiagram( "k1, k2", "v1, v2" : rational := true );
<A loop diagram with loop momenta [ k1, k2 ] & external momenta [ v1, v2 ]>
gap> cos:= v1*v2;;
gap> SetAbbreviation( cos, "cos" );
gap> rel := [ v1^2 - 1, v2^2 - 1 ];;
gap> SetRelationsOfExternalMomenta( LD, rel );
gap> SetIndependentLorentzInvariants( LD,
>         [ k1^2, k1*k2, k1*v1, k1*v2, k2^2, k2*v1, k2*v2, cos ] );
gap> SetPropagators( LD,
>         -[ 2*k2*v1 - 1, 2*k2*v2 - 1, (k1 - k2)^2,
>            2*k1*v1 - 1, 2*k1*v2 - 1, k1^2 ] );
gap> SetNumerators( LD, -[ k2^2 ] );
gap> SetExtraLorentzInvariants( LD, [ cos ] );
gap> E12 := PairOfMatricesOfLoopDiagramInPropagators( LD );
[ <A 6 x 8 matrix over an external ring>,
  <A 6 x 6 matrix over an external ring> ]
gap> Display( E12[1] );
0,       0,        0,       0,       D4-1,     D1-1,    -2,      (-2*cos),
0,       0,        0,       0,       D5-1,     D2-1,    (-2*cos),-2,
D3+D6-N7,-D3+D6-N7,-D1+D4,  -D2+D5,  -D3-D6+N7,D3-D6+N7,D1-D4,   D2-D5,
D4-1,    D1-1,     -2,      (-2*cos),0,        0,       0,       0,
D5-1,    D2-1,     (-2*cos),-2,      0,        0,       0,       0,
2*D6,    -D3+D6+N7,D4-1,    D5-1,    0,        0,       0,       0
gap> s := SyzygiesOfColumns( E12 );
<A non-zero 8 x 35 matrix over an external ring>

1.2-11 2 Loop box
gap> LD := LoopDiagram( "l1..2", "k1..2,k4" );
<A loop diagram with loop momenta [ l1, l2 ] &
 external momenta [ k1, k2, k4 ]>
gap> s12 := 2*k1*k2;;
gap> SetAbbreviation( s12, "s12" );
gap> s14 := 2*k1*k4;;
gap> SetAbbreviation( s14, "s14" );
gap> rel1 := List( ExternalMomenta( LD ), k -> k^2 );;
gap> rel2 := [ (k1+k2+k4)^2 ];;
gap> SetRelationsOfExternalMomenta( LD, Concatenation( rel1, rel2 ) );
gap> SetIndependentLorentzInvariants( LD,
>         [ l1^2, l1*l2, l2^2, l1*k1, l1*k2, l1*k4,
>           l2*k1, l2*k2, l2*k4, s12, s14 ] );
gap> SetPropagators( LD,
>         -[ l1^2, (l1-k1)^2, (l1-k1-k2)^2,
>            l2^2, (l2-k4)^2, (l2+k1+k2)^2,
>            (l1+l2)^2 ] );
gap> SetNumerators( LD, -[ l1*k4, l2*k1 ] );
gap> SetExtraLorentzInvariants( LD, [ s12, s14 ] );
gap> E12 := PairOfMatricesOfLoopDiagramInPropagators( LD );
[ <A 7 x 10 matrix over an external ring>,
  <A 7 x 7 matrix over an external ring> ]
gap> S := SyzygiesOfColumns( E12 );
<A non-zero 10 x 95 matrix over an external ring>
gap> Display( EntriesOfHomalgMatrix( CertainColumns( S, [ 1 ] ) ) );
[ -D2+D3+2*N8, 0, -2*N8, D1, -D2, 0, -2*D4+D5+D6-2*N9, -D4+D5, -D4, D4+2*N9 ]
gap> Display( EntriesOfHomalgMatrix( CertainColumns( S, [ 2 ] ) ) );
[ D1-D2+2*N8, 0, -D1-2*N8, 0, -D2, 0, -D4+D5-2*N9, D5, 0, D4+2*N9 ]
gap> Display( EntriesOfHomalgMatrix( CertainColumns( S, [ 3 ] ) ) );
[ 0, 0, 0, 0, 0, 0, D5*D7-D6*D7-2*D7*N9, D4*D7+D5*D7, D4*D7, D4*D7+2*D7*N9 ]
gap> Sibp1 := IBPRelation( CertainColumns( S, [ 1 ] ), LD );;
gap> ViewList( DecomposeInMonomials( Sibp1 ) );
[ [ |[ -1/2*s12*a8-1/2*s14*a8 ]|, |[ D1*N8_ ]| ],
  [ |[ -1/2*s12*a9+1/2*s14*a9 ]|, |[ D4*N9_ ]| ],
  [ |[ -d+a1+a2+a3+a7+a8 ]|, |[ D2 ]| ],
  [ |[ d-a1-a2-a3-a7-a8 ]|, |[ D3 ]| ],
  [ |[ -2*d+2*a4+2*a5+2*a6+2*a7+2*a9 ]|, |[ D4 ]| ],
  [ |[ d-a4-a5-a6-a7-a9 ]|, |[ D5 ]| ],
  [ |[ d-a4-a5-a6-a7-a9 ]|, |[ D6 ]| ],
  [ |[ 2*d-2*a1-2*a2-2*a3-2*a7-2*a8 ]|, |[ N8 ]| ],
  [ |[ -2*d+2*a4+2*a5+2*a6+2*a7+2*a9 ]|, |[ N9 ]| ],
  [ |[ s12*a1+s14*a2-s12*a3+s12*a4-s14*a5-s12*a6-s14*a8+s14*a9 ]|, |[ 1 ]| ] ]
gap> sibp1 := IBPRelation( CertainColumns( S, [ 1 ] ), LD, [ 1, 1, 1, 1, 1, 1, 1, 1, 1 ] );
|[ -1/2*s12*D1*N8_-1/2*s14*D1*N8_-1/2*s12*D4*N9_+1/2*s14*D4*N9_-d*D2+d*D3\
   -2*d*D4+d*D5+d*D6+2*d*N8-2*d*N9+5*D2-5*D3+10*D4-5*D5-5*D6-10*N8+10*N9 ]|
gap> ViewList( DecomposeInMonomials( sibp1 ) );
[ [ |[ -1/2*s12-1/2*s14 ]|, |[ D1*N8_ ]| ],
  [ |[ -1/2*s12+1/2*s14 ]|, |[ D4*N9_ ]| ],
  [ |[ -d+5 ]|, |[ D2 ]| ],
  [ |[ d-5 ]|, |[ D3 ]| ],
  [ |[ -2*d+10 ]|, |[ D4 ]| ],
  [ |[ d-5 ]|, |[ D5 ]| ],
  [ |[ d-5 ]|, |[ D6 ]| ],
  [ |[ 2*d-10 ]|, |[ N8 ]| ],
  [ |[ -2*d+10 ]|, |[ N9 ]| ] ]
gap> Sibp2 := IBPRelation( CertainColumns( S, [ 2 ] ), LD );;
gap> ViewList( DecomposeInMonomials( Sibp2 ) );
[ [ |[ -1/2*s14*a8 ]|, |[ D1*N8_ ]| ],
  [ |[ 1/2*s14*a9 ]|, |[ D4*N9_ ]| ],
  [ |[ d-a1-a2-a3-a7-a8 ]|, |[ D1 ]| ],
  [ |[ -d+a1+a2+a3+a7+a8 ]|, |[ D2 ]| ],
  [ |[ -d+a4+a5+a6+a7+a9 ]|, |[ D4 ]| ],
  [ |[ d-a4-a5-a6-a7-a9 ]|, |[ D5 ]| ],
  [ |[ 2*d-2*a1-2*a2-2*a3-2*a7-2*a8 ]|, |[ N8 ]| ],
  [ |[ -2*d+2*a4+2*a5+2*a6+2*a7+2*a9 ]|, |[ N9 ]| ],
  [ |[ s14*a2-s14*a5-s14*a8+s14*a9 ]|, |[ 1 ]| ] ]
gap> sibp2 := IBPRelation( CertainColumns( S, [ 2 ] ), LD, [ 1, 1, 1, 1, 1, 1, 1, 1, 1 ] );
|[ -1/2*s14*D1*N8_+1/2*s14*D4*N9_+d*D1-d*D2-d*D4+d*D5\
   +2*d*N8-2*d*N9-5*D1+5*D2+5*D4-5*D5-10*N8+10*N9 ]|
gap> ViewList( DecomposeInMonomials( sibp2 ) );
[ [ |[ -1/2*s14 ]|, |[ D1*N8_ ]| ],
  [ |[ 1/2*s14 ]|, |[ D4*N9_ ]| ],
  [ |[ d-5 ]|, |[ D1 ]| ],
  [ |[ -d+5 ]|, |[ D2 ]| ],
  [ |[ -d+5 ]|, |[ D4 ]| ],
  [ |[ d-5 ]|, |[ D5 ]| ],
  [ |[ 2*d-10 ]|, |[ N8 ]| ],
  [ |[ -2*d+10 ]|, |[ N9 ]| ] ]
gap> Sibp3 := IBPRelation( CertainColumns( S, [ 3 ] ), LD );;
gap> ViewList( DecomposeInMonomials( Sibp3 ) );
[ [ |[ 1/2*s12*a9+1/2*s14*a9 ]|, |[ D4*D7*N9_ ]| ],
  [ |[ -a7+1 ]|, |[ D1*D4 ]| ],
  [ |[ a7-1 ]|, |[ D3*D4 ]| ],
  [ |[ a7-1 ]|, |[ D2*D5 ]| ],
  [ |[ -a7+1 ]|, |[ D1*D6 ]| ],
  [ |[ d-a4-a5-a6-a7-a9+1 ]|, |[ D5*D7 ]| ],
  [ |[ -d+a4+a5+a6+a7+a9-1 ]|, |[ D6*D7 ]| ],
  [ |[ -2*a7+2 ]|, |[ D4*N8 ]| ],
  [ |[ -2*a7+2 ]|, |[ D1*N9 ]| ],
  [ |[ -2*d+2*a4+2*a5+2*a6+2*a7+2*a9-2 ]|, |[ D7*N9 ]| ],
  [ |[ -4*a7+4 ]|, |[ N8*N9 ]| ],
  [ |[ -s12*a4-s14*a5+s12*a6+s14*a9 ]|, |[ D7 ]| ] ]
gap> sibp3 := IBPRelation( CertainColumns( S, [ 3 ] ), LD, [ 1, 1, 1, 1, 1, 1, 1, 1, 1 ] );
|[ 1/2*s12*D4*D7*N9_+1/2*s14*D4*D7*N9_+d*D5*D7-d*D6*D7\
   -2*d*D7*N9-4*D5*D7+4*D6*D7+8*D7*N9 ]|
gap> ViewList( DecomposeInMonomials( sibp3 ) );
[ [ |[ 1/2*s12+1/2*s14 ]|, |[ D4*D7*N9_ ]| ],
  [ |[ d-4 ]|, |[ D5*D7 ]| ],
  [ |[ -d+4 ]|, |[ D6*D7 ]| ],
  [ |[ -2*d+8 ]|, |[ D7*N9 ]| ] ]
gap> gen := GeneratorsOfScalelessSectors( LD );
<A 1 x 10 matrix over an external ring>
gap> EntriesOfHomalgMatrix( gen );
[ D6*D7*N8*N9, D4*D7*N8*N9, D3*D7*N8*N9, D1*D7*N8*N9, D4*D5*D6*N8*N9,
  D3*D5*D6*N8*N9, D2*D3*D6*N8*N9, D1*D4*D5*N8*N9, D1*D2*D4*N8*N9,
  D1*D2*D3*N8*N9 ]

1.2-12 Pentagon with off shell leg
gap> LD := LoopDiagram( "k1", "p1, p2, p3, p4" );
<A loop diagram with loop momenta [ k1 ] &
 external momenta [ p1, p2, p3, p4 ]>
gap> pp5:= 2*p1*p2 + 2*p1*p3 + 2*p2*p3 + 2*p1*p4 + 2*p2*p4 + 2*p3*p4;;
gap> SetAbbreviation( pp5, "pp5" );
gap> s12:= 2*p1*p2;;
gap> SetAbbreviation( s12, "s12" );
gap> s23:= 2*p2*p3;;
gap> SetAbbreviation( s23, "s23" );
gap> s34:= 2*p3*p4;;
gap> SetAbbreviation( s34, "s34" );
gap> s45:= 2*p1*p2 + 2*p1*p3 + 2*p2*p3;;
gap> SetAbbreviation( s45, "s45" );
gap> s51:= 2*p2*p3 + 2*p2*p4 + 2*p3*p4;;
gap> SetAbbreviation( s51, "s51" );
gap> rel := [ p1^2, p2^2, p3^2, p4^2 ];;
gap> SetRelationsOfExternalMomenta( LD, rel );
gap> SetIndependentLorentzInvariants( LD,
>         [ k1^2, k1*p1, k1*p2, k1*p3, k1*p4, pp5, s12, s23, s34, s45, s51 ] );
gap> SetPropagators( LD,
>         -[ k1^2, (k1 + p1)^2, (k1 + p1 + p2)^2, (k1 + p1 + p2 + p3)^2,
>            (k1 + p1 + p2 + p3 + p4)^2 ] );
gap> SetNumerators( LD, -[  ] );
gap> SetExtraLorentzInvariants( LD, [ pp5, s12, s23, s34, s45, s51 ] );
gap> E12 := PairOfMatricesOfLoopDiagramInPropagators( LD );
[ <A 5 x 5 matrix over an external ring>,
  <A 5 x 5 matrix over an external ring> ]
gap> Display( E12[1] );
2*D1,     -D1+D2,        s12-D2+D3,    -s12+s45-D3+D4,pp5-s45-D4+D5,
D1+D2,    -D1+D2,        -D2+D3,       s23-D3+D4,     -s23+s51-D4+D5,
s12+D1+D3,-s12-D1+D2,    -D2+D3,       -D3+D4,        s34-D4+D5,
s45+D1+D4,s23-s45-D1+D2, -s23-D2+D3,   -D3+D4,        -D4+D5,
pp5+D1+D5,-pp5+s51-D1+D2,s34-s51-D2+D3,-s34-D3+D4,    -D4+D5
gap> S := SyzygiesOfColumns( E12 );
<A non-zero 5 x 76 matrix over an external ring>
gap> Display( EntriesOfHomalgMatrix( CertainColumns( S, [ 1 ] ) ) );
[ D1*D4*D5-D3*D4*D5, D1*D4*D5, D1*D4*D5, 0, 0 ]

1.3 Global variables

1.3-1 LOOP_INTEGRALS
‣ LOOP_INTEGRALS( global variable )

1.4 Properties

1.4-1 IsLoopMomentum
‣ IsLoopMomentum( arg )( property )

Returns: true or false

1.4-2 IsExternalMomentum
‣ IsExternalMomentum( arg )( property )

Returns: true or false

1.5 Attributes

1.5-1 Dimension
‣ Dimension( arg )( attribute )

1.5-2 UnderlyingMatrix
‣ UnderlyingMatrix( arg )( attribute )

1.5-3 LoopDiagram
‣ LoopDiagram( arg )( attribute )

1.5-4 Components
‣ Components( arg )( attribute )

1.5-5 LoopMomenta
‣ LoopMomenta( LD )( attribute )

1.5-6 ExternalMomenta
‣ ExternalMomenta( LD )( attribute )

1.5-7 DimensionOfCoefficientsVector
‣ DimensionOfCoefficientsVector( LD )( attribute )

1.5-8 UnderlyingRing
‣ UnderlyingRing( LD )( attribute )

1.5-9 RelationsOfExternalMomenta
‣ RelationsOfExternalMomenta( LD )( attribute )

1.5-10 RelationsMatrixOfExternalMomenta
‣ RelationsMatrixOfExternalMomenta( LD )( attribute )

1.5-11 IndependentLorentzInvariants
‣ IndependentLorentzInvariants( LD )( attribute )

1.5-12 Propagators
‣ Propagators( LD )( attribute )

1.5-13 Numerators
‣ Numerators( LD )( attribute )

1.5-14 ExtraLorentzInvariants
‣ ExtraLorentzInvariants( LD )( attribute )

1.5-15 OriginalIBPGeneratingMatrixOfPropagators
‣ OriginalIBPGeneratingMatrixOfPropagators( LD )( attribute )

1.5-16 OriginalIBPGeneratingMatrixOfNumerators
‣ OriginalIBPGeneratingMatrixOfNumerators( LD )( attribute )

1.5-17 OriginalTaylorOfPropagators
‣ OriginalTaylorOfPropagators( LD )( attribute )

1.5-18 MatrixOfMomenta
‣ MatrixOfMomenta( LD )( attribute )

1.5-19 IBPGeneratingMatrixOfPropagators
‣ IBPGeneratingMatrixOfPropagators( LD )( attribute )

1.5-20 IBPGeneratingMatrixOfNumerators
‣ IBPGeneratingMatrixOfNumerators( LD )( attribute )

1.5-21 PairOfMatricesOfLoopDiagram
‣ PairOfMatricesOfLoopDiagram( LD )( attribute )

1.5-22 RingOfExtraLorentzInvariants
‣ RingOfExtraLorentzInvariants( LD )( attribute )

1.5-23 ReductionMatrixOfExtraLorentzInvariants
‣ ReductionMatrixOfExtraLorentzInvariants( LD )( attribute )

1.5-24 RingOfIndependentLorentzInvariants
‣ RingOfIndependentLorentzInvariants( LD )( attribute )

1.5-25 ReductionMatrixOfIndependentLorentzInvariants
‣ ReductionMatrixOfIndependentLorentzInvariants( LD )( attribute )

1.5-26 RingOfPropagatorsAndNumeratorsAndExtraLorentzInvariants
‣ RingOfPropagatorsAndNumeratorsAndExtraLorentzInvariants( LD )( attribute )

1.5-27 RingOfLoopDiagram
‣ RingOfLoopDiagram( LD )( attribute )

1.5-28 ReductionMatrixOfPropagatorsAndNumeratorsAndExtraLorentzInvariants
‣ ReductionMatrixOfPropagatorsAndNumeratorsAndExtraLorentzInvariants( LD )( attribute )

1.5-29 IBPGeneratingMatrixOfPropagatorsInIndependentLorentzInvariants
‣ IBPGeneratingMatrixOfPropagatorsInIndependentLorentzInvariants( LD )( attribute )

1.5-30 IBPGeneratingMatrixOfPropagatorsInPropagators
‣ IBPGeneratingMatrixOfPropagatorsInPropagators( LD )( attribute )

1.5-31 IBPGeneratingMatrixOfNumeratorsInPropagators
‣ IBPGeneratingMatrixOfNumeratorsInPropagators( LD )( attribute )

1.5-32 PairOfMatricesOfLoopDiagramInIndependentLorentzInvariants
‣ PairOfMatricesOfLoopDiagramInIndependentLorentzInvariants( LD )( attribute )

1.5-33 PairOfMatricesOfLoopDiagramInPropagators
‣ PairOfMatricesOfLoopDiagramInPropagators( LD )( attribute )

1.5-34 IBPGeneratingMatrixOfLoopDiagramInPropagators
‣ IBPGeneratingMatrixOfLoopDiagramInPropagators( LD )( attribute )

1.5-35 MatrixOfIBPRelations
‣ MatrixOfIBPRelations( LD )( attribute )

1.5-36 BasisOfIBPRelations
‣ BasisOfIBPRelations( LD )( attribute )

1.5-37 MatrixOfSpecialIBPRelations
‣ MatrixOfSpecialIBPRelations( LD )( operation )

1.5-38 BasisOfSpecialIBPRelations
‣ BasisOfSpecialIBPRelations( LD )( operation )

1.5-39 MatrixOfIBPRelationsInWeylAlgebra
‣ MatrixOfIBPRelationsInWeylAlgebra( LD )( attribute )

1.5-40 BasisOfIBPRelationsInWeylAlgebra
‣ BasisOfIBPRelationsInWeylAlgebra( LD )( attribute )

1.5-41 MatrixOfSpecialIBPRelationsInWeylAlgebra
‣ MatrixOfSpecialIBPRelationsInWeylAlgebra( LD )( attribute )

1.5-42 BasisOfSpecialIBPRelationsInWeylAlgebra
‣ BasisOfSpecialIBPRelationsInWeylAlgebra( LD )( attribute )

1.5-43 FieldOfCoefficientsOfLoopDiagramInSingular
‣ FieldOfCoefficientsOfLoopDiagramInSingular( LD )( attribute )

1.5-44 FieldOfCoefficientsOfLoopDiagramInMaple
‣ FieldOfCoefficientsOfLoopDiagramInMaple( LD )( attribute )

1.5-45 FieldOfCoefficientsOfLoopDiagramInHecke
‣ FieldOfCoefficientsOfLoopDiagramInHecke( LD )( attribute )

1.5-46 MatrixOfCoefficientsOfParametricIBPs
‣ MatrixOfCoefficientsOfParametricIBPs( LD )( attribute )

1.5-47 SymanzikPolynomials
‣ SymanzikPolynomials( LD )( attribute )

1.5-48 AssociatedWeylAlgebra
‣ AssociatedWeylAlgebra( LD )( attribute )

1.5-49 GeneratorsOfScalelessSectors
‣ GeneratorsOfScalelessSectors( LD )( attribute )

1.5-50 GeneratorsOfScalelessSectorsInWeylAlgebra
‣ GeneratorsOfScalelessSectorsInWeylAlgebra( LD )( attribute )

1.6 Constructors

1.6-1 LorentzVector
‣ LorentzVector( arg1, arg2 )( operation )

1.6-2 LorentzVector
‣ LorentzVector( arg )( operation )

1.6-3 LoopDiagram
‣ LoopDiagram( arg1, arg2, arg3 )( operation )
‣ LoopDiagram( arg1, arg2 )( operation )

1.7 Operations

1.7-1 SetAbbreviation
‣ SetAbbreviation( arg1, arg2 )( operation )

1.7-2 +
‣ +( arg1, arg2 )( operation )

1.7-3 -
‣ -( arg1, arg2 )( operation )

1.7-4 AdditiveInverse
‣ AdditiveInverse( arg )( attribute )

1.7-5 []
‣ []( arg1, arg2 )( operation )

1.7-6 *
‣ *( arg1, arg2 )( operation )

1.7-7 *
‣ *( arg1, arg2 )( operation )

1.7-8 *
‣ *( arg1, arg2 )( operation )

1.7-9 ExpressInExtraLorentzInvariants
‣ ExpressInExtraLorentzInvariants( arg1, arg2 )( operation )

1.7-10 ExpressInIndependentLorentzInvariants
‣ ExpressInIndependentLorentzInvariants( arg1, arg2 )( operation )

1.7-11 ExpressInPropagatorsAndNumeratorsAndExtraLorentzInvariants
‣ ExpressInPropagatorsAndNumeratorsAndExtraLorentzInvariants( arg1, arg2 )( operation )

1.7-12 JacobianOfCoefficientsVectorInPropagators
‣ JacobianOfCoefficientsVectorInPropagators( vec, LD )( operation )

1.7-13 DivergenceOfCoefficientsVectorOfLoopDiagram
‣ DivergenceOfCoefficientsVectorOfLoopDiagram( vec, LD )( operation )

1.7-14 DoubleShiftAlgebra
‣ DoubleShiftAlgebra( R )( operation )

1.7-15 RationalDoubleShiftAlgebra
‣ RationalDoubleShiftAlgebra( R )( operation )

1.7-16 IBPRelation
‣ IBPRelation( vec, LD )( operation )

1.7-17 IBPRelation
‣ IBPRelation( vec, LD, exponents )( operation )

1.7-18 MatrixOfIBPRelations
‣ MatrixOfIBPRelations( mat, LD )( operation )

1.7-19 MatrixOfIBPRelations
‣ MatrixOfIBPRelations( mat, LD, exponents )( operation )

1.7-20 MatrixOfIBPRelations
‣ MatrixOfIBPRelations( LD, exponents )( operation )

1.7-21 MatrixOfSpecialIBPRelations
‣ MatrixOfSpecialIBPRelations( LD, exponents )( operation )

1.7-22 IBPRelationInWeylAlgebra
‣ IBPRelationInWeylAlgebra( vec, LD )( operation )

1.7-23 MatrixOfIBPRelationsInWeylAlgebra
‣ MatrixOfIBPRelationsInWeylAlgebra( mat, LD )( operation )

1.7-24 MatrixOfCoefficientsOfIBPs
‣ MatrixOfCoefficientsOfIBPs( LD )( attribute )

1.7-25 MatrixOfCoefficientsOfIBPs
‣ MatrixOfCoefficientsOfIBPs( LD )( operation )

1.7-26 MatrixOfCoefficientsOfParametricIBPs
‣ MatrixOfCoefficientsOfParametricIBPs( LD, degree, ring )( operation )

1.7-27 MatrixOfCoefficientsOfParametricIBPs
‣ MatrixOfCoefficientsOfParametricIBPs( LD, degree )( operation )

1.7-28 ColumnReversedMatrixOfCoefficientsOfParametricIBPs
‣ ColumnReversedMatrixOfCoefficientsOfParametricIBPs( LD, degree, ring )( operation )

1.7-29 ColumnReversedMatrixOfCoefficientsOfParametricIBPs
‣ ColumnReversedMatrixOfCoefficientsOfParametricIBPs( LD, degree )( operation )

1.7-30 SymanzikPolynomials
‣ SymanzikPolynomials( LD, list_of_ones )( operation )

1.7-31 DegreesOfMonomialsOfProductOfSymanzikPolynomials
‣ DegreesOfMonomialsOfProductOfSymanzikPolynomials( LD, list_of_ones )( operation )

1.7-32 IsScalelessLoopIntegral
‣ IsScalelessLoopIntegral( LD, list_of_ones )( operation )

1.7-33 GeneratorsOfScalelessSectors
‣ GeneratorsOfScalelessSectors( LD, exponents )( operation )

1.7-34 NormalForm
‣ NormalForm( operator, G )( operation )

1.7-35 NormalForm
‣ NormalForm( operator, G, initial_integral )( operation )

1.7-36 NormalFormWrtInitialIntegral
‣ NormalFormWrtInitialIntegral( operator, G )( operation )
 [Top of Book]  [Contents]   [Previous Chapter]   [Next Chapter] 
Goto Chapter: Top 1 Ind

generated by GAPDoc2HTML