5. Reference Guide

This section describes all public functions and classes in pySecDec.

5.1. Algebra

Implementation of a simple computer algebra system.

class pySecDec.algebra.ExponentiatedPolynomial(expolist, coeffs, exponent=1, polysymbols='x', copy=True)

Like Polynomial, but with a global exponent. polynomial^{exponent}

Parameters:
  • expolist – iterable of iterables; The variable’s powers for each term.
  • coeffs – iterable; The coefficients of the polynomial.
  • exponent – object, optional; The global exponent.
  • polysymbols

    iterable or string, optional; The symbols to be used for the polynomial variables when converted to string. If a string is passed, the variables will be consecutively numbered.

    For example: expolist=[[2,0],[1,1]] coeffs=[“A”,”B”]
    • polysymbols=’x’ (default) <-> “A*x0**2 + B*x0*x1”
    • polysymbols=[‘x’,’y’] <-> “A*x**2 + B*x*y”
  • copy

    bool; Whether or not to copy the expolist, the coeffs, and the exponent.

    Note

    If copy is False, it is assumed that the expolist, the coeffs and the exponent have the correct type.

copy()

Return a copy of a Polynomial or a subclass.

derive(index)

Generate the derivative by the parameter indexed index.

Parameters:index – integer; The index of the paramater to derive by.
simplify()

Apply the identity <something>**0 = 1 or <something>**1 = <something> or 1**<something> = 1 if possible, otherwise call the simplify method of the base class. Convert exponent to symbol if possible.

pySecDec.algebra.Expression(expression, polysymbols, follow_functions=False)

Convert a sympy expression to an expression in terms of this module.

Parameters:
  • expression – string or sympy expression; The expression to be converted
  • polysymbols – iterable of strings or sympy symbols; The symbols to be stored as expolists (see Polynomial) where possible.
  • follow_functions – bool, optional (default = False); If true, return the converted expression and a list of Function that occur in the expression.
class pySecDec.algebra.Function(symbol, *arguments, **kwargs)

Symbolic function that can take care of parameter transformations. It keeps track of all taken derivatives: When derive() is called, save the multiindex of the taken derivative.

The derivative multiindices are the keys in the dictionary self.derivative_tracks. The values are lists with two elements: Its first element is the index to derive the derivative indicated by the multiindex in the second element by, in order to abtain the derivative indicated by the key:

>>> from pySecDec.algebra import Polynomial, Function
>>> x = Polynomial.from_expression('x', ['x','y'])
>>> y = Polynomial.from_expression('y', ['x','y'])
>>> poly = x**2*y + y**2
>>> func = Function('f', x, y)
>>> ddfuncd0d1 = func.derive(0).derive(1)
>>> func
Function(f( + (1)*x, + (1)*y), derivative_tracks = {(1, 0): [0, (0, 0)], (1, 1): [1, (1, 0)]})
>>> func.derivative_tracks
{(1, 0): [0, (0, 0)], (1, 1): [1, (1, 0)]}
>>> func.compute_derivatives(poly)
{(1, 0):  + (2)*x*y, (1, 1):  + (2)*x}
Parameters:
  • symbol – string; The symbol to be used to represent the Function.
  • arguments – arbitrarily many _Expression; The arguments of the Function.
  • copy – bool; Whether or not to copy the arguments.
compute_derivatives(expression=None)

Compute all derivatives of expression that are mentioned in self.derivative_tracks. The purpose of this function is to avoid computing the same derivatives multiple times.

Parameters:expression_Expression, optional; The expression to compute the derivatives of. If not provided, the derivatives are shown as in terms of the function‘s derivatives dfd<index>.
copy()

Return a copy of a Function.

derive(index)

Generate the derivative by the parameter indexed index. The derivative of a function with symbol f by some index is denoted as dfd<index>.

Parameters:index – integer; The index of the paramater to derive by.
replace(expression, index, value, remove=False)

Replace a variable in an expression by a number or a symbol. The entries in all expolist of the underlying Polynomial are set to zero. The coefficients are modified according to value and the powers indicated in the expolist.

Parameters:
  • expression_Expression; The expression to replace the variable.
  • index – integer; The index of the variable to be replaced.
  • value – number or sympy expression; The value to insert for the chosen variable.
  • remove – bool; Whether or not to remove the replaced parameter from the parameters in the expression.
simplify()

Simplify the arguments.

class pySecDec.algebra.Log(arg, copy=True)

The (natural) logarithm to base e (2.718281828459..). Store the expressions log(arg).

Parameters:
  • arg_Expression; The argument of the logarithm.
  • copy – bool; Whether or not to copy the arg.
copy()

Return a copy of a Log.

derive(index)

Generate the derivative by the parameter indexed index.

Parameters:index – integer; The index of the paramater to derive by.
replace(expression, index, value, remove=False)

Replace a variable in an expression by a number or a symbol. The entries in all expolist of the underlying Polynomial are set to zero. The coefficients are modified according to value and the powers indicated in the expolist.

Parameters:
  • expression_Expression; The expression to replace the variable.
  • index – integer; The index of the variable to be replaced.
  • value – number or sympy expression; The value to insert for the chosen variable.
  • remove – bool; Whether or not to remove the replaced parameter from the parameters in the expression.
simplify()

Apply log(1) = 0.

class pySecDec.algebra.LogOfPolynomial(expolist, coeffs, polysymbols='x', copy=True)

The natural logarithm of a Polynomial.

Parameters:
  • expolist – iterable of iterables; The variable’s powers for each term.
  • coeffs – iterable; The coefficients of the polynomial.
  • exponent – object, optional; The global exponent.
  • polysymbols

    iterable or string, optional; The symbols to be used for the polynomial variables when converted to string. If a string is passed, the variables will be consecutively numbered.

    For example: expolist=[[2,0],[1,1]] coeffs=[“A”,”B”]
    • polysymbols=’x’ (default) <-> “A*x0**2 + B*x0*x1”
    • polysymbols=[‘x’,’y’] <-> “A*x**2 + B*x*y”
derive(index)

Generate the derivative by the parameter indexed index.

Parameters:index – integer; The index of the paramater to derive by.
static from_expression(expression, polysymbols)

Alternative constructor. Construct the LogOfPolynomial from an algebraic expression.

Parameters:
  • expression – string or sympy expression; The algebraic representation of the polynomial, e.g. “5*x1**2 + x1*x2”
  • polysymbols – iterable of strings or sympy symbols; The symbols to be interpreted as the polynomial variables, e.g. “[‘x1’,’x2’]”.
simplify()

Apply the identity log(1) = 0, otherwise call the simplify method of the base class.

class pySecDec.algebra.Polynomial(expolist, coeffs, polysymbols='x', copy=True)

Container class for polynomials. Store a polynomial as list of lists counting the powers of the variables. For example the polynomial “x1**2 + x1*x2” is stored as [[2,0],[1,1]].

Coefficients are stored in a separate list of strings, e.g. “A*x0**2 + B*x0*x1” <-> [[2,0],[1,1]] and [“A”,”B”].

Parameters:
  • expolist

    iterable of iterables; The variable’s powers for each term.

    Hint

    Negative powers are allowed.

  • coeffs – 1d array-like with numerical or sympy-symbolic (see http://www.sympy.org/) content, e.g. [x,1,2] where x is a sympy symbol; The coefficients of the polynomial.
  • polysymbols

    iterable or string, optional; The symbols to be used for the polynomial variables when converted to string. If a string is passed, the variables will be consecutively numbered.

    For example: expolist=[[2,0],[1,1]] coeffs=[“A”,”B”]
    • polysymbols=’x’ (default) <-> “A*x0**2 + B*x0*x1”
    • polysymbols=[‘x’,’y’] <-> “A*x**2 + B*x*y”
  • copy

    bool; Whether or not to copy the expolist and the coeffs.

    Note

    If copy is False, it is assumed that the expolist and the coeffs have the correct type.

becomes_zero_for(zero_params)

Return True if the polynomial becomes zero if the parameters passed in zero_params are set to zero. Otherwise, return False.

Parameters:zero_params – iterable of integers; The indices of the parameters to be checked.
copy()

Return a copy of a Polynomial or a subclass.

derive(index)

Generate the derivative by the parameter indexed index.

Parameters:index – integer; The index of the paramater to derive by.
static from_expression(expression, polysymbols)

Alternative constructor. Construct the polynomial from an algebraic expression.

Parameters:
  • expression – string or sympy expression; The algebraic representation of the polynomial, e.g. “5*x1**2 + x1*x2”
  • polysymbols – iterable of strings or sympy symbols; The symbols to be interpreted as the polynomial variables, e.g. “[‘x1’,’x2’]”.
has_constant_term(indices=None)

Return True if the polynomial can be written as:

const + ...

Otherwise, return False.

Parameters:indices – list of integers or None; The indices of the polysymbols to consider. If None (default) all indices are taken into account.
replace(expression, index, value, remove=False)

Replace a variable in an expression by a number or a symbol. The entries in all expolist of the underlying Polynomial are set to zero. The coefficients are modified according to value and the powers indicated in the expolist.

Parameters:
  • expression_Expression; The expression to replace the variable.
  • index – integer; The index of the variable to be replaced.
  • value – number or sympy expression; The value to insert for the chosen variable.
  • remove – bool; Whether or not to remove the replaced parameter from the parameters in the expression.
simplify(deep=True)

Combine terms that have the same exponents of the variables.

Parameters:deep – bool; If True (default) call the simplify method of the coefficients if they are of type _Expression.
class pySecDec.algebra.Pow(base, exponent, copy=True)

Exponential. Store two expressions A and B to be interpreted as the exponential A**B.

Parameters:
  • base_Expression; The base A of the exponential.
  • exponent_Expression; The exponent B.
  • copy – bool; Whether or not to copy base and exponent.
copy()

Return a copy of a Pow.

derive(index)

Generate the derivative by the parameter indexed index.

Parameters:index – integer; The index of the paramater to derive by.
replace(expression, index, value, remove=False)

Replace a variable in an expression by a number or a symbol. The entries in all expolist of the underlying Polynomial are set to zero. The coefficients are modified according to value and the powers indicated in the expolist.

Parameters:
  • expression_Expression; The expression to replace the variable.
  • index – integer; The index of the variable to be replaced.
  • value – number or sympy expression; The value to insert for the chosen variable.
  • remove – bool; Whether or not to remove the replaced parameter from the parameters in the expression.
simplify()

Apply the identity <something>**0 = 1 or <something>**1 = <something> or 1**<something> = 1 if possible. Convert to ExponentiatedPolynomial or Polynomial if possible.

class pySecDec.algebra.Product(*factors, **kwargs)

Product of polynomials. Store one or polynomials p_i to be interpreted as product \prod_i p_i.

Parameters:
  • factors – arbitrarily many instances of Polynomial; The factors p_i.
  • copy – bool; Whether or not to copy the factors.

p_i can be accessed with self.factors[i].

Example:

p = Product(p0, p1)
p0 = p.factors[0]
p1 = p.factors[1]
copy()

Return a copy of a Product.

derive(index)

Generate the derivative by the parameter indexed index. Return an instance of the optimized ProductRule.

Parameters:index – integer; The index of the paramater to derive by.
replace(expression, index, value, remove=False)

Replace a variable in an expression by a number or a symbol. The entries in all expolist of the underlying Polynomial are set to zero. The coefficients are modified according to value and the powers indicated in the expolist.

Parameters:
  • expression_Expression; The expression to replace the variable.
  • index – integer; The index of the variable to be replaced.
  • value – number or sympy expression; The value to insert for the chosen variable.
  • remove – bool; Whether or not to remove the replaced parameter from the parameters in the expression.
simplify()

If one or more of self.factors is a Product, replace it by its factors. If only one factor is present, return that factor. Remove factors of one and zero.

class pySecDec.algebra.ProductRule(*expressions, **kwargs)

Store an expression of the form

\sum_i c_i \prod_j \prod_k
\left( \frac{d}{dx_k} \right)^{n_{ijk}}
f_j \left( \lbrace x_k \rbrace \right)

The main reason for introducing this class is a speedup when calculating derivatives. In particular, this class implements simplifications such that the number of terms grows less than exponentially (scaling of the naive implementation of the product rule) with the number of derivatives.

Parameters:expressions – arbitrarily many expressions; The expressions f_j.
copy()

Return a copy of a ProductRule.

derive(index)

Generate the derivative by the parameter indexed index. Note that this class is particularly designed to hold derivatives of a product.

Parameters:index – integer; The index of the paramater to derive by.
replace(index, value, remove=False)

Replace a variable in an expression by a number or a symbol. The entries in all expolist of the underlying Polynomial are set to zero. The coefficients are modified according to value and the powers indicated in the expolist.

Parameters:
  • expression_Expression; The expression to replace the variable.
  • index – integer; The index of the variable to be replaced.
  • value – number or sympy expression; The value to insert for the chosen variable.
  • remove – bool; Whether or not to remove the replaced parameter from the parameters in the expression.
simplify()

Combine terms that have the same derivatives of the expressions.

to_sum()

Convert the ProductRule to Sum

class pySecDec.algebra.Sum(*summands, **kwargs)

Sum of polynomials. Store one or polynomials p_i to be interpreted as product \sum_i p_i.

Parameters:
  • summands – arbitrarily many instances of Polynomial; The summands p_i.
  • copy – bool; Whether or not to copy the summands.

p_i can be accessed with self.summands[i].

Example:

p = Sum(p0, p1)
p0 = p.summands[0]
p1 = p.summands[1]
copy()

Return a copy of a Sum.

derive(index)

Generate the derivative by the parameter indexed index.

Parameters:index – integer; The index of the paramater to derive by.
replace(expression, index, value, remove=False)

Replace a variable in an expression by a number or a symbol. The entries in all expolist of the underlying Polynomial are set to zero. The coefficients are modified according to value and the powers indicated in the expolist.

Parameters:
  • expression_Expression; The expression to replace the variable.
  • index – integer; The index of the variable to be replaced.
  • value – number or sympy expression; The value to insert for the chosen variable.
  • remove – bool; Whether or not to remove the replaced parameter from the parameters in the expression.
simplify()

If one or more of self.summands is a Sum, replace it by its summands. If only one summand is present, return that summand. Remove zero from sums.

5.2. Loop Integral

This module defines routines to Feynman parametrize a loop integral and build a c++ package that numerically integrates over the sector decomposed integrand.

5.2.1. Feynman Parametrization

Routines to Feynman parametrize a loop integral.

class pySecDec.loop_integral.LoopIntegral(*args, **kwargs)

Container class for loop integrals. The main purpose of this class is to convert a loop integral from the momentum representation to the Feynman parameter representation.

It is possible to provide either the graph of the loop integrals as adjacency list, or the propagators.

The Feynman parametrized integral is a product of the following expressions that are accessible as member properties:

  • self.regulator ** self.regulator_power
  • self.Gamma_factor
  • self.exponentiated_U
  • self.exponentiated_F
  • self.numerator
  • self.measure,

where self is an instance of either LoopIntegralFromGraph or LoopIntegralFromPropagators.

When inverse propagators or nonnumerical propagator powers are present (see powerlist), some Feynman_parameters drop out of the integral. The variables to integrate over can be accessed as self.integration_variables.

While self.numerator describes the numerator polynomial generated by tensor numerators or inverse propagators, self.measure contains the monomial associated with the integration measure in the case of propagator powers \neq 1. The Gamma functions in the denominator belonging to the measure, however, are multiplied to the overall Gamma factor given by self.Gamma_factor. The overall sign (-1)^{N_\nu} is included in self.numerator.

See also

class pySecDec.loop_integral.LoopIntegralFromGraph(internal_lines, external_lines, replacement_rules=[], Feynman_parameters='x', regulator='eps', regulator_power=0, dimensionality='4-2*eps', powerlist=[])

Construct the Feynman parametrization of a loop integral from the graph using the cut construction method.

Example:

>>> from pySecDec.loop_integral import *
>>> internal_lines = [['0',[1,2]], ['m',[2,3]], ['m',[3,1]]]
>>> external_lines = [['p1',1],['p2',2],['-p12',3]]
>>> li = LoopIntegralFromGraph(internal_lines, external_lines)
>>> li.exponentiated_U
( + (1)*x0 + (1)*x1 + (1)*x2)**(2*eps - 1)
>>> li.exponentiated_F
( + (m**2)*x2**2 + (2*m**2 - p12**2)*x1*x2 + (m**2)*x1**2 + (m**2 - p1**2)*x0*x2 + (m**2 - p2**2)*x0*x1)**(-eps - 1)
Parameters:
  • internal_lines – iterable of internal line specification, consisting of string or sympy expression for mass and a pair of strings or numbers for the vertices, e.g. [[‘m’, [1,2]], [‘0’, [2,1]]].
  • external_lines – iterable of external line specification, consisting of string or sympy expression for external momentum and a strings or number for the vertex, e.g. [[‘p1’, 1], [‘p2’, 2]].
  • replacement_rules – iterable of iterables with two strings or sympy expressions, optional; Symbolic replacements to be made for the external momenta, e.g. definition of Mandelstam variables. Example: [(‘p1*p2’, ‘s’), (‘p1**2’, 0)] where p1 and p2 are external momenta. It is also possible to specify vector replacements, for example [(‘p4’, ‘-(p1+p2+p3)’)].
  • Feynman_parameters – iterable or string, optional; The symbols to be used for the Feynman parameters. If a string is passed, the Feynman parameter variables will be consecutively numbered starting from zero.
  • regulator

    string or sympy symbol, optional; The symbol to be used for the dimensional regulator (typically \epsilon or \epsilon_D)

    Note

    If you change this symbol, you have to adapt the dimensionality accordingly.

  • regulator_power

    integer; An additional factor to the numerator.

    See also

    LoopIntegral

  • dimensionality – string or sympy expression, optional; The dimensionality; typically 4-2\epsilon, which is the default value.
  • powerlist – iterable, optional; The powers of the propergators, possibly dependent on the regulator. In case of negative powers, the numerator is constructed by taking derivatives with respect to the corresponding Feynman parameters as explained in Section 3.2.4 of Ref. [BHJ+15]. If negative powers are combined with a tensor numerator, the derivatives act on the Feynman-parametrized tensor numerator as well, which leads to a consistent result.
class pySecDec.loop_integral.LoopIntegralFromPropagators(propagators, loop_momenta, external_momenta=[], Lorentz_indices=[], numerator=1, metric_tensor='g', replacement_rules=[], Feynman_parameters='x', regulator='eps', regulator_power=0, dimensionality='4-2*eps', powerlist=[])

Construct the Feynman parametrization of a loop integral from the algebraic momentum representation.

See also

[Hei08], [GKR+11]

Example:

>>> from pySecDec.loop_integral import *
>>> propagators = ['k**2', '(k - p)**2']
>>> loop_momenta = ['k']
>>> li = LoopIntegralFromPropagators(propagators, loop_momenta)
>>> li.exponentiated_U
( + (1)*x0 + (1)*x1)**(2*eps - 2)
>>> li.exponentiated_F
( + (-p**2)*x0*x1)**(-eps)

The 1st (U) and 2nd (F) Symanzik polynomials and their exponents can also be accessed independently:

>>> li.U
 + (1)*x0 + (1)*x1
>>> li.F
 + (-p**2)*x0*x1
>>>
>>> li.exponent_U
2*eps - 2
>>> li.exponent_F
-eps
Parameters:
  • propagators – iterable of strings or sympy expressions; The propagators, e.g. [‘k1**2’, ‘(k1-k2)**2 - m1**2’].
  • loop_momenta – iterable of strings or sympy expressions; The loop momenta, e.g. [‘k1’,’k2’].
  • external_momenta

    iterable of strings or sympy expressions, optional; The external momenta, e.g. [‘p1’,’p2’]. Specifying the external_momenta is only required when a numerator is to be constructed.

    See also

    parameter numerator

  • Lorentz_indices

    iterable of strings or sympy expressions, optional; Symbols to be used as Lorentz indices in the numerator.

    See also

    parameter numerator

  • numerator

    string or sympy expression, optional; The numerator of the loop integral. Scalar products must be passed in index notation e.g. “k1(mu)*k2(mu)”. The numerator should be a sum of products of exclusively: * numbers * scalar products (e.g. “p1(mu)*k1(mu)*p1(nu)*k2(nu)”) * symbols (e.g. “m”)

    Examples:
    • p1(mu)*k1(mu)*p1(nu)*k2(nu) + 4*s*eps*k1(mu)*k1(mu)
    • p1(mu)*(k1(mu) + k2(mu))*p1(nu)*k2(nu)
    • p1(mu)*k1(mu)*my_function(eps)

    Warning

    All Lorentz indices (including the contracted ones and also including the numbers that have been used) must be explicitly defined using the parameter Lorentz_indices.

    Warning

    It is assumed that the numerator is and all its derivatives by the regulator are finite and defined if \epsilon=0 is inserted explicitly. In particular, if user defined functions (like in the example p1(mu)*k1(mu)*my_function(eps)) appear, make sure that my_function(0) is finite.

    Hint

    In order to mimic a singular user defined function, use the parameter regulator_power. For example, instead of numerator = gamma(eps) you could enter numerator = eps_times_gamma(eps) in conjunction with regulator_power = -1.

    Hint

    It is possible to use numbers as indices, for example p1(mu)*p2(mu)*k1(nu)*k2(nu) = p1(1)*p2(1)*k1(2)*k2(2).

    Hint

    The numerator may have uncontracted indices, e.g. k1(mu)*k2(nu).

  • metric_tensor – string or sympy symbol, optional; The symbol to be used for the (Minkowski) metric tensor g^{\mu\nu}.
  • replacement_rules – iterable of iterables with two strings or sympy expressions, optional; Symbolic replacements to be made for the external momenta, e.g. definition of Mandelstam variables. Example: [(‘p1*p2’, ‘s’), (‘p1**2’, 0)] where p1 and p2 are external momenta. It is also possible to specify vector replacements, for example [(‘p4’, ‘-(p1+p2+p3)’)].
  • Feynman_parameters – iterable or string, optional; The symbols to be used for the Feynman parameters. If a string is passed, the Feynman parameter variables will be consecutively numbered starting from zero.
  • regulator

    string or sympy symbol, optional; The symbol to be used for the dimensional regulator (typically \epsilon or \epsilon_D)

    Note

    If you change this symbol, you have to adapt the dimensionality accordingly.

  • regulator_power

    integer; An additional factor to the numerator.

    See also

    LoopIntegral

  • dimensionality – string or sympy expression, optional; The dimensionality; typically 4-2\epsilon, which is the default value.
  • powerlist – iterable, optional; The powers of the propergators, possibly dependent on the regulator. In case of negative powers, the numerator is constructed by taking derivatives with respect to the corresponding Feynman parameters as explained in Section 3.2.4 of Ref. [BHJ+15]. If negative powers are combined with a tensor numerator, the derivatives act on the Feynman-parametrized tensor numerator as well, which leads to a consistent result.

5.2.2. Loop Package

This module contains the function that generates a c++ package.

pySecDec.loop_integral.loop_package(name, loop_integral, requested_order, real_parameters=[], complex_parameters=[], contour_deformation=True, additional_prefactor=1, form_optimization_level=2, form_work_space='500M', decomposition_method='iterative', normaliz_executable='normaliz', enforce_complex=False, split=False, ibp_power_goal=-1, use_dreadnaut=True)

Decompose, subtract and expand a Feynman parametrized loop integral. Return it as c++ package.

See also

This function is a wrapper around pySecDec.code_writer.make_package().

See also

The generated library is described in Generated C++ Libraries.

Parameters:
  • name – string; The name of the c++ namespace and the output directory.
  • loop_integralpySecDec.loop_integral.LoopIntegral; The loop integral to be computed.
  • requested_orders – integer; Compute the expansion in the regulator to this order.
  • real_parameters – iterable of strings or sympy symbols, optional; Parameters to be interpreted as real numbers, e.g. Mandelstam invariants and masses.
  • complex_parameters – iterable of strings or sympy symbols, optional; Parameters to be interpreted as complex numbers. To use the complex mass scheme, define the masses as complex parameters.
  • contour_deformation – bool, optional; Whether or not to produce code for contour deformation. Default: True.
  • additional_prefactor – string or sympy expression, optional; An additional factor to be multiplied to the loop integral. It may depend on the regulator, the real_parameters, and the complex_parameters.
  • form_optimization_level – integer out of the interval [0,3], optional; The optimization level to be used in FORM. Default: 2.
  • form_work_space – string, optional; The FORM WorkSpace. Default: '500M'.
  • decomposition_method

    string, optional; The strategy for decomposing the polynomials. The following strategies are available:

    • ‘iterative’ (default)
    • ‘geometric’
    • ‘geometric_ku’

    Note

    For ‘geometric’ and ‘geometric_ku’, the third-party program “normaliz” is needed. See The Geomethod and Normaliz.

  • normaliz_executable – string, optional; The command to run normaliz. normaliz is only required if decomposition_method is set to ‘geometric’ or ‘geometric_ku’. Default: ‘normaliz’
  • enforce_complex – bool, optional; Whether or not the generated integrand functions should have a complex return type even though they might be purely real. The return type of the integrands is automatically complex if contour_deformation is True or if there are complex_parameters. In other cases, the calculation can typically be kept purely real. Most commonly, this flag is needed if log(<negative real>) occurs in one of the integrand functions. However, pySecDec will suggest setting this flag to True in that case. Default: False
  • split – bool, optional; Whether or not to split the integration domain in order to map singularities from 1 to 0. Set this option to True if you have singularties when one or more integration variables are one. Default: False
  • ibp_power_goal

    integer, optional; The power_goal that is forwarded to integrate_by_parts().

    This option controls how the subtraction terms are generated. Setting it to -numpy.inf disables integrate_by_parts(), while 0 disables integrate_pole_part().

    See also

    To generate the subtraction terms, this function first calls integrate_by_parts() for each integration variable with the give ibp_power_goal. Then integrate_pole_part() is called.

    Default: -1

  • use_dreadnaut – bool or string, optional; Whether or not to use squash_symmetry_redundant_sectors_dreadnaut() to find sector symmetries. If given a string, interpret that string as the command line executable dreadnaut. If True, try $SECDEC_CONTRIB/bin/dreadnaut and, if the environment variable $SECDEC_CONTRIB is not set, dreadnaut. Default: True

5.2.3. Drawing Feynman Diagrams

Use the following function to draw Feynman diagrams.

pySecDec.loop_integral.draw.plot_diagram(internal_lines, external_lines, filename, powerlist=None, neato='neato', extension='pdf', Gstart=0)

Draw a Feynman diagram using Graphviz (neato).

Thanks to Viktor Papara <papara@mpp.mpg.de> for his major contributions to this function.

Note

This function requires the command line tool neato. See also Drawing Feynman Diagrams with neato.

Warning

The target is overwritten without prompt if it exists already.

Parameters:
  • internal_lines – list; Adjacency list of internal lines, e.g. [['m',['a',4]],['m',[4,5]], ['m',['a',5]],[0,[1,2]],[0,[4,1]],[0,[2,5]]]
  • external_lines – list; Adjacency list of external lines, e.g. [[‘p1’,1],[‘p2’,2],[‘p3’,’a’]]
  • filename – string; The name of the output file. The generated file gets this name plus the file extension.
  • powerlist – list, optional; The powers of the propagators defined by the internal_lines.
  • neato – string, default: “neato”; The shell command to call “neato”.
  • extension – string, default: “pdf”; The file extension. This also defines the output format.
  • Gstart – nonnegative int; The is value is passed to “neato” with the “-Gstart” option. Try changing this value if the visualization looks bad.

5.3. Decomposition

The core of sector decomposition. This module implements the actual decomposition routines.

5.3.1. Common

This module collects routines that are used by multiple decompition modules.

class pySecDec.decomposition.Sector(cast, other=[], Jacobian=None)

Container class for sectors that arise during the sector decomposition.

Parameters:
  • cast – iterable of algebra.Product or of algebra.Polynomial; The polynomials to be cast to the form <monomial> * (const + ...)
  • other – iterable of algebra.Polynomial, optional; All variable transformations are applied to these polynomials but it is not attempted to achieve the form <monomial> * (const + ...)
  • Jacobianalgebra.Polynomial with one term, optional; The Jacobian determinant of this sector. If not provided, the according unit monomial (1*x0^0*x1^0...) is assumed.
pySecDec.decomposition.squash_symmetry_redundant_sectors_sort(sectors, sort_function)

Reduce a list of sectors by squashing duplicates with equal integral.

If two sectors only differ by a permutation of the polysymbols (to be interpreted as integration variables over some inteval), then the two sectors integrate to the same value. Thus we can drop one of them and count the other twice. The multiple counting of a sector is accounted for by increasing the coefficient of the Jacobian by one.

Equivalence up to permutation is established by applying the sort_function to each sector, this brings them into a canonical form. Sectors with identical canonical forms differ only by a permutation.

Note: whether all symmetries are found depends on the choice of sort_function. Neither pySecDec.matrix_sort.iterative_sort() nor pySecDec.matrix_sort.Pak_sort() find all symmetries.

See also: squash_symmetry_redundant_sectors_dreadnaut()

Example:

>>> from pySecDec.algebra import Polynomial
>>> from pySecDec.decomposition import Sector
>>> from pySecDec.decomposition import squash_symmetry_redundant_sectors_sort
>>> from pySecDec.matrix_sort import Pak_sort
>>>
>>> poly = Polynomial([(0,1),(1,0)], ['a','b'])
>>> swap = Polynomial([(1,0),(0,1)], ['a','b'])
>>> Jacobian_poly = Polynomial([(1,0)], [3]) # three
>>> Jacobian_swap = Polynomial([(0,1)], [5]) # five
>>> sectors = (
...               Sector([poly],Jacobian=Jacobian_poly),
...               Sector([swap],Jacobian=Jacobian_swap)
...           )
>>>
>>> reduced_sectors = squash_symmetry_redundant_sectors_sort(sectors,
...                   Pak_sort)
>>> len(reduced_sectors) # symmetry x0 <--> x1
1
>>> # The Jacobians are added together to account
>>> # for the double counting of the sector.
>>> reduced_sectors[0].Jacobian
 + (8)*x0
Parameters:
pySecDec.decomposition.squash_symmetry_redundant_sectors_dreadnaut(sectors, dreadnaut='dreadnaut', workdir='dreadnaut_tmp', keep_workdir=False)

Reduce a list of sectors by squashing duplicates with equal integral.

Each Sector is converted to a Polynomial which is represented as a graph following the example of [BKAP] (v2.6 Figure 7, Isotopy of matrices).

We first multiply each polynomial in the sector by a unique tag then sum the polynomials of the sector, this converts a sector to a polynomial. Next, we convert the expolist of the resulting polynomial to a graph where each unique exponent in the expolist is considered to be a different symbol. Each unique coefficient in the polynomial`s coeffs is assigned a vertex and connected to the row vertex of any term it multiplies. The external program dreadnaut is then used to bring the graph into a canonical form and provide a hash. Sectors with equivalent hashes may be identical, their canonical graphs are compared and if they are identical the sectors are combined.

Note

This function calls the command line executable of dreadnaut [BKAP]. It has been tested with dreadnaut version nauty26r7.

See also: squash_symmetry_redundant_sectors_sort()

Parameters:
  • sectors – iterable of Sector; the sectors to be reduced.
  • dreadnaut – string; The shell command to run dreadnaut.
  • workdir

    string; The directory for the communication with dreadnaut. A directory with the specified name will be created in the current working directory. If the specified directory name already exists, an OSError is raised.

    Note

    The communication with dreadnaut is done via files.

  • keep_workdir – bool; Whether or not to delete the workdir after execution.

5.3.2. Iterative

The iterative sector decomposition routines.

exception pySecDec.decomposition.iterative.EndOfDecomposition

This exception is raised if the function iteration_step() is called although the sector is already in standard form.

pySecDec.decomposition.iterative.find_singular_set(sector, indices=None)

Function within the iterative sector decomposition procedure which heuristically chooses an optimal decomposition set. The strategy was introduced in arXiv:hep-ph/0004013 [BH00] and is described in 4.2.2 of arXiv:1410.7939 [Bor14]. Return a list of indices.

Parameters:
  • sectorSector; The sector to be decomposed.
  • indices – iterable of integers or None; The indices of the parameters to be considered as integration variables. By default (indices=None), all parameters are considered as integration variables.
pySecDec.decomposition.iterative.iteration_step(sector, indices=None)

Run a single step of the iterative sector decomposition as described in chapter 3.2 (part II) of arXiv:0803.4177v2 [Hei08]. Return an iterator of Sector - the arising subsectors.

Parameters:
  • sectorSector; The sector to be decomposed.
  • indices – iterable of integers or None; The indices of the parameters to be considered as integration variables. By default (indices=None), all parameters are considered as integration variables.
pySecDec.decomposition.iterative.iterative_decomposition(sector, indices=None)

Run the iterative sector decomposition as described in chapter 3.2 (part II) of arXiv:0803.4177v2 [Hei08]. Return an iterator of Sector - the arising subsectors.

Parameters:
  • sectorSector; The sector to be decomposed.
  • indices – iterable of integers or None; The indices of the parameters to be considered as integration variables. By default (indices=None), all parameters are considered as integration variables.
pySecDec.decomposition.iterative.primary_decomposition(sector, indices=None)

Perform the primary decomposition as described in chapter 3.2 (part I) of arXiv:0803.4177v2 [Hei08]. Return a list of Sector - the primary sectors. For N Feynman parameters, there are N primary sectors where the i-th Feynman parameter is set to 1 in sector i.

Parameters:
  • sectorSector; The container holding the polynomials (typically U and F) to eliminate the Dirac delta from.
  • indices – iterable of integers or None; The indices of the parameters to be considered as integration variables. By default (indices=None), all parameters are considered as integration variables.
pySecDec.decomposition.iterative.primary_decomposition_polynomial(polynomial, indices=None)

Perform the primary decomposition on a single polynomial.

Parameters:
  • polynomialalgebra.Polynomial; The polynomial to eliminate the Dirac delta from.
  • indices – iterable of integers or None; The indices of the parameters to be considered as integration variables. By default (indices=None), all parameters are considered as integration variables.
pySecDec.decomposition.iterative.remap_parameters(singular_parameters, Jacobian, *polynomials)

Remap the Feynman parameters according to eq. (16) of arXiv:0803.4177v2 [Hei08]. The parameter whose index comes first in singular_parameters is kept fix.

The remapping is done in place; i.e. the polynomials are NOT copied.

Parameters:
  • singular_parameters – list of integers; The indices \alpha_r such that at least one of polynomials becomes zero if all t_{\alpha_r} \rightarrow 0.
  • JacobianPolynomial; The Jacobian determinant is multiplied to this polynomial.
  • polynomials – abritrarily many instances of algebra.Polynomial where all of these have an equal number of variables; The polynomials of Feynman parameters to be remapped. These are typically F and U.

Example:

remap_parameters([1,2], Jacobian, F, U)

5.3.3. Geometric

The geometric sector decomposition routines.

pySecDec.decomposition.geometric.Cheng_Wu(sector, index=-1)

Replace one Feynman parameter by one. This means integrating out the Dirac delta according to the Cheng-Wu theorem.

Parameters:
  • sectorSector; The container holding the polynomials (typically U and F) to eliminate the Dirac delta from.
  • index – integer, optional; The index of the Feynman parameter to eliminate. Default: -1 (the last Feynman parameter)
class pySecDec.decomposition.geometric.Polytope(vertices=None, facets=None)

Representation of a polytope defined by either its vertices or its facets. Call complete_representation() to translate from one to the other representation.

Parameters:
  • vertices – two dimensional array; The polytope in vertex representation. Each row is interpreted as one vertex.
  • facets

    two dimensional array; The polytope in facet representation. Each row represents one facet F. A row in facets is interpreted as one normal vector n_F with additionally the constant a_F in the last column. The points v of the polytope obey

    \bigcap_F \left( {\langle n_F, v \rangle} + a_F \right) \ge 0

complete_representation(normaliz='normaliz', workdir='normaliz_tmp', keep_workdir=False)

Transform the vertex representation of a polytope to the facet representation or the other way round. Remove surplus entries in self.facets or self.vertices.

Note

This function calls the command line executable of normaliz [BIR]. It has been tested with normaliz versions 3.0.0, 3.1.0, and 3.1.1.

Parameters:
  • normaliz – string; The shell command to run normaliz.
  • workdir

    string; The directory for the communication with normaliz. A directory with the specified name will be created in the current working directory. If the specified directory name already exists, an OSError is raised.

    Note

    The communication with normaliz is done via files.

  • keep_workdir – bool; Whether or not to delete the workdir after execution.
vertex_incidence_lists()

Return for each vertex the list of facets it lies in (as dictonary). The keys of the output dictonary are the vertices while the values are the indices of the facets in self.facets.

pySecDec.decomposition.geometric.convex_hull(*polynomials)

Calculate the convex hull of the Minkowski sum of all polynomials in the input. The algorithm sets all coefficients to one first and then only keeps terms of the polynomial product that have coefficient 1. Return the list of these entries in the expolist of the product of all input polynomials.

Parameters:polynomials – abritrarily many instances of Polynomial where all of these have an equal number of variables; The polynomials to calculate the convex hull for.
pySecDec.decomposition.geometric.generate_fan(*polynomials)

Calculate the fan of the polynomials in the input. The rays of a cone are given by the exponent vectors after factoring out a monomial together with the standard basis vectors. Each choice of factored out monomials gives a different cone. Only full (N-) dimensional cones in R^N_{\geq 0} need to be considered.

Parameters:polynomials – abritrarily many instances of Polynomial where all of these have an equal number of variables; The polynomials to calculate the fan for.
pySecDec.decomposition.geometric.geometric_decomposition(sector, indices=None, normaliz='normaliz', workdir='normaliz_tmp')

Run the sector decomposition using the geomethod as described in [BHJ+15].

Note

This function calls the command line executable of normaliz [BIR]. It has been tested with normaliz versions 3.0.0, 3.1.0, and 3.1.1.

Parameters:
  • sectorSector; The sector to be decomposed.
  • indices – list of integers or None; The indices of the parameters to be considered as integration variables. By default (indices=None), all parameters are considered as integration variables.
  • normaliz – string; The shell command to run normaliz.
  • workdir

    string; The directory for the communication with normaliz. A directory with the specified name will be created in the current working directory. If the specified directory name already exists, an OSError is raised.

    Note

    The communication with normaliz is done via files.

pySecDec.decomposition.geometric.geometric_decomposition_ku(sector, indices=None, normaliz='normaliz', workdir='normaliz_tmp')

Run the sector decomposition using the original geometric decomposition strategy by Kaneko and Ueda as described in [KU10].

Note

This function calls the command line executable of normaliz [BIR]. It has been tested with normaliz versions 3.0.0, 3.1.0, and 3.1.1.

Parameters:
  • sectorSector; The sector to be decomposed.
  • indices – list of integers or None; The indices of the parameters to be considered as integration variables. By default (indices=None), all parameters are considered as integration variables.
  • normaliz – string; The shell command to run normaliz.
  • workdir

    string; The directory for the communication with normaliz. A directory with the specified name will be created in the current working directory. If the specified directory name already exists, an OSError is raised.

    Note

    The communication with normaliz is done via files.

pySecDec.decomposition.geometric.transform_variables(polynomial, transformation, polysymbols='y')

Transform the parameters x_i of a pySecDec.algebra.Polynomial,

x_i \rightarrow \prod_j x_j^{T_{ij}}

, where T_{ij} is the transformation matrix.

Parameters:
  • polynomialpySecDec.algebra.Polynomial; The polynomial to transform the variables in.
  • transformation – two dimensional array; The transformation matrix T_{ij}.
  • polysymbols – string or iterable of strings; The symbols for the new variables. This argument is passed to the default constructor of pySecDec.algebra.Polynomial. Refer to the documentation of pySecDec.algebra.Polynomial for further details.
pySecDec.decomposition.geometric.triangulate(cone, normaliz='normaliz', workdir='normaliz_tmp', keep_workdir=False, switch_representation=False)

Split a cone into simplicial cones; i.e. cones defined by exactly D rays where D is the dimensionality.

Note

This function calls the command line executable of normaliz [BIR]. It has been tested with normaliz versions 3.0.0, 3.1.0, and 3.1.1.

Parameters:
  • cone – two dimensional array; The defining rays of the cone.
  • normaliz – string; The shell command to run normaliz.
  • workdir

    string; The directory for the communication with normaliz. A directory with the specified name will be created in the current working directory. If the specified directory name already exists, an OSError is raised.

    Note

    The communication with normaliz is done via files.

  • keep_workdir – bool; Whether or not to delete the workdir after execution.
  • switch_representation – bool; Whether or not to switch between facet and vertex/ray representation.

5.3.4. Splitting

Routines to split the integration between 0 and 1. This maps singularities from 1 to 0.

pySecDec.decomposition.splitting.find_singular_sets_at_one(polynomial)

Find all possible sets of parameters such that the polynomial‘s constant term vanishes if these parameters are set to one.

Example:

>>> from pySecDec.algebra import Polynomial
>>> from pySecDec.decomposition.splitting import find_singular_sets_at_one
>>> polysymbols = ['x0', 'x1']
>>> poly = Polynomial.from_expression('1 - 10*x0 - x1', polysymbols)
>>> find_singular_sets_at_one(poly)
[(1,)]
Parameters:polynomialPolynomial; The polynomial to search in.
pySecDec.decomposition.splitting.remap_one_to_zero(polynomial, *indices)

Apply the transformation x \rightarrow 1 - x to polynomial for the parameters of the given indices.

Parameters:
  • polynomialPolynomial; The polynomial to apply the transformation to.
  • indices – arbitrarily many int; The indices of the polynomial.polysymbols to apply the transformation to.

Example:

>>> from pySecDec.algebra import Polynomial
>>> from pySecDec.decomposition.splitting import remap_one_to_zero
>>> polysymbols = ['x0']
>>> polynomial = Polynomial.from_expression('x0', polysymbols)
>>> remap_one_to_zero(polynomial, 0)
 + (1) + (-1)*x0
pySecDec.decomposition.splitting.split(sector, seed, *indices)

Split the integration interval [0,1] for the parameters given by indices. The splitting point is fixed using numpy’s random number generator.

Return an iterator of Sector - the arising subsectors.

Parameters:sectorSector; The sector to be split.
:param seed;
integer; The seed for the random number generator that is used to fix the splitting point.
Parameters:indices – arbitrarily many integers; The indices of the variables to be split.
pySecDec.decomposition.splitting.split_singular(sector, seed, indices=[])

Split the integration interval [0,1] for the parameters that can lead to singularities at one for the polynomials in sector.cast.

Return an iterator of Sector - the arising subsectors.

Parameters:
  • sectorSector; The sector to be split.
  • seed – integer; The seed for the random number generator that is used to fix the splitting point.
  • indices – iterables of integers; The indices of the variables to be split if required. An empty iterator means that all variables may potentially be split.

5.4. Matrix Sort

Algorithms to sort a matrix when column and row permutations are allowed.

pySecDec.matrix_sort.Pak_sort(matrix)

Inplace modify the matrix to some ordering, when permutations of rows and columns (excluding the first) are allowed. The implementation of this function is describedin chapter 2 of [Pak11].

Note

This function may result in different orderings depending on the initial ordering.

See also

iterative_sort()

Parameters:matrix – 2D array-like; The matrix to be canonicalized.
pySecDec.matrix_sort.iterative_sort(matrix)

Inplace modify the matrix to some ordering, when permutations of rows and columns (excluding the first) are allowed.

Note

This function may result in different orderings depending on the initial ordering.

See also

Pak_sort()

Parameters:matrix – 2D array-like; The matrix to be canonicalized.

5.5. Subtraction

Routines to isolate the divergencies in an \epsilon expansion.

pySecDec.subtraction.integrate_by_parts(polyprod, power_goal, *indices)

Repeatedly apply integration by parts,

\int_0^1
{
    dt_j t_j^{(a - b \epsilon_1 - c \epsilon_2 + ...)}
    \mathcal{I} (t_j,\{t_{i \neq j}\}, \epsilon_1, \epsilon_2, ...)
}
=
\frac{1}{a + 1 - b \epsilon_1 - c \epsilon_2 - ...}
\left(
\mathcal{I} (1,\{t_{i \neq j}\}, \epsilon_1, \epsilon_2, ...)
-
\int_0^1
{
    dt_j t_j^{(a + 1 - b \epsilon_1 - c \epsilon_2 + ...)}
    \mathcal{I}^{\prime} (t_j,\{t_{i \neq j}\}, \epsilon_1, \epsilon_2, ...)
}
\right)

, where \mathcal{I}^{\prime} denotes the derivative of \mathcal{I} with respect to t_j. The iteration stops, when a>= power_goal.

See also

This function provides an alternative to integrate_pole_part().

Parameters:
  • polyprodalgebra.Product of the form <product of <monomial>**(a_j + ...)> * <regulator poles of cal_I> * <cal_I>; The input product as decribed above. The <product of <monomial>**(a_j + ...)> should be a pySecDec.algebra.Product of <monomial>**(a_j + ...). as described below. The <monomial>**(a_j + ...) should be an pySecDec.algebra.ExponentiatedPolynomial with exponent being a Polynomial of the regulators \epsilon_1, \epsilon_2, .... Although no dependence on the Feynman parameters is expected in the exponent, the polynomial variables should be the Feynman parameters and the regulators. The constant term of the exponent should be numerical. The polynomial variables of monomial and the other factors (interpreted as \mathcal{I}) are interpreted as the Feynman parameters and the epsilon regulators. Make sure that the last factor (<cal_I>) is defined and finite for \epsilon = 0. All poles for \epsilon \rightarrow 0 should be made explicit by putting them into <regulator poles of cal_I> as pySecDec.algebra.Pow with exponent = -1 and the base of type pySecDec.algebra.Polynomial.
  • power_goal – number, e.g. float, integer, ...; The stopping criterion for the iteration.
  • indices – arbitrarily many integers; The index/indices of the parameter(s) to partially integrate. j in the formulae above.

Return the pole part and the numerically integrable remainder as a list. Each returned list element has the same structure as the input polyprod.

pySecDec.subtraction.integrate_pole_part(polyprod, *indices)

Transform an integral of the form

\int_0^1
{
    dt_j t_j^{(a - b \epsilon_1 - c \epsilon_2 + ...)}
    \mathcal{I} (t_j,\{t_{i \neq j}\}, \epsilon_1, \epsilon_2, ...)
}

into the form

\sum_{p=0}^{|a|-1}
{
    \frac{1}{a + p + 1 - b \epsilon_1 - c \epsilon_2 - ...}
    \frac{\mathcal{I}^{(p)} (0,\{t_{i \neq j}\}, \epsilon_1, \epsilon_2, ...)}{p!} +
    \int_0^1
    {
        dt_j t_j^{(a - b \epsilon_1 - c \epsilon_2 + ...)}
        R(t_j,\{t_{i \neq j}\}, \epsilon_1, \epsilon_2, ...)
    }
}

, where \mathcal{I}^{(p)} denotes the p-th derivative of \mathcal{I} with respect to t_j. The equations above are to be understood schematically.

See also

This function implements the transformation from equation (19) to (21) as described in arXiv:0803.4177v2 [Hei08].

Parameters:
  • polyprodalgebra.Product of the form <product of <monomial>**(a_j + ...)> * <regulator poles of cal_I> * <cal_I>; The input product as decribed above. The <product of <monomial>**(a_j + ...)> should be a pySecDec.algebra.Product of <monomial>**(a_j + ...). as described below. The <monomial>**(a_j + ...) should be an pySecDec.algebra.ExponentiatedPolynomial with exponent being a Polynomial of the regulators \epsilon_1, \epsilon_2, .... Although no dependence on the Feynman parameters is expected in the exponent, the polynomial variables should be the Feynman parameters and the regulators. The constant term of the exponent should be numerical. The polynomial variables of monomial and the other factors (interpreted as \mathcal{I}) are interpreted as the Feynman parameters and the epsilon regulators. Make sure that the last factor (<cal_I>) is defined and finite for \epsilon = 0. All poles for \epsilon \rightarrow 0 should be made explicit by putting them into <regulator poles of cal_I> as pySecDec.algebra.Pow with exponent = -1 and the base of type pySecDec.algebra.Polynomial.
  • indices – arbitrarily many integers; The index/indices of the parameter(s) to partially integrate. j in the formulae above.

Return the pole part and the numerically integrable remainder as a list. That is the sum and the integrand of equation (21) in arXiv:0803.4177v2 [Hei08]. Each returned list element has the same structure as the input polyprod.

pySecDec.subtraction.pole_structure(monomial_product, *indices)

Return a list of the unregulated exponents of the parameters specified by indices in monomial_product.

Parameters:
  • monomial_productpySecDec.algebra.ExponentiatedPolynomial with exponent being a Polynomial; The monomials of the subtraction to extract the pole structure from.
  • indices – arbitrarily many integers; The index/indices of the parameter(s) to partially investigate.

5.6. Expansion

Routines to series expand singular and nonsingular expressions.

exception pySecDec.expansion.OrderError

This exception is raised if an expansion to a lower than the lowest order of an expression is requested.

pySecDec.expansion.expand_Taylor(expression, indices, orders)

Series/Taylor expand a nonsingular expression around zero.

Return a algebra.Polynomial - the series expansion.

Parameters:
  • expression – an expression composed of the types defined in the module algebra; The expression to be series expanded.
  • indices – integer or iterable of integers; The indices of the parameters to expand. The ordering of the indices defines the ordering of the expansion.
  • order – integer or iterable of integers; The order to which the expansion is to be calculated.
pySecDec.expansion.expand_singular(product, indices, orders)

Series expand a potentially singular expression of the form

{
    \frac{a_N \epsilon_0 + b_N \epsilon_1 + ...}
         {a_D \epsilon_0 + b_D \epsilon_1 + ...}
}

Return a algebra.Polynomial - the series expansion.

See also

To expand more general expressions use expand_sympy().

Parameters:
  • productalgebra.Product with factors of the form <polynomial> and <polynomial> ** -1; The expression to be series expanded.
  • indices – integer or iterable of integers; The indices of the parameters to expand. The ordering of the indices defines the ordering of the expansion.
  • order – integer or iterable of integers; The order to which the expansion is to be calculated.
pySecDec.expansion.expand_sympy(expression, variables, orders)

Expand a sympy expression in the variables to given orders. Return the expansion as nested pySecDec.algebra.Polynomial.

See also

This function is a generalization of expand_singular().

Parameters:
  • expression – string or sympy expression; The expression to be expanded
  • variables – iterable of strings or sympy symbols; The variables to expand the expression in.
  • orders – iterable of integers; The orders to expand to.

5.7. Code Writer

This module collects routines to create a c++ library.

5.7.1. Make Package

This is the main function of pySecDec.

pySecDec.code_writer.make_package(name, integration_variables, regulators, requested_orders, polynomials_to_decompose, polynomial_names=[], other_polynomials=[], prefactor=1, remainder_expression=1, functions=[], real_parameters=[], complex_parameters=[], form_optimization_level=2, form_work_space='500M', form_insertion_depth=5, contour_deformation_polynomial=None, positive_polynomials=[], decomposition_method='iterative_no_primary', normaliz_executable='normaliz', enforce_complex=False, split=False, ibp_power_goal=-1, use_dreadnaut=True)

Decompose, subtract and expand an expression. Return it as c++ package.

See also

In order to decompose a loop integral, use the function pySecDec.loop_integral.loop_package().

See also

The generated library is described in Generated C++ Libraries.

Parameters:
  • name – string; The name of the c++ namepace and the output directory.
  • integration_variables – iterable of strings or sympy symbols; The variables that are to be integrated from 0 to 1.
  • regulators – iterable of strings or sympy symbols; The UV/IR regulators of the integral.
  • requested_orders – iterable of integers; Compute the expansion in the regulators to these orders.
  • polynomials_to_decompose – iterable of strings or sympy expressions or pySecDec.algebra.ExponentiatedPolynomial or pySecDec.algebra.Polynomial; The polynomials to be decomposed.
  • polynomial_names – iterable of strings; Assign symbols for the polynomials_to_decompose. These can be referenced in the other_polynomials; see other_polynomials for details.
  • other_polynomials

    iterable of strings or sympy expressions or pySecDec.algebra.ExponentiatedPolynomial or pySecDec.algebra.Polynomial; Additional polynomials where no decomposition is attempted. The symbols defined in polynomial_names can be used to reference the polynomials_to_decompose. This is particularly useful when computing loop integrals where the “numerator” can depend on the first and second Symanzik polynomials.

    Example (1-loop bubble with numerator):

    >>> polynomials_to_decompose = ["(x0 + x1)**(2*eps - 4)",
    ...                             "(-p**2*x0*x1)**(-eps))"]
    >>> polynomial_names = ["U", "F"]
    >>> other_polynomials = ["""   (eps - 1)*s*U**2
    ...                          + (eps - 2)*F
    ...                          - (eps - 1)*2*s*x0*U
    ...                          + (eps - 1)*s*x0**2"""]
    

    Note that the polynomial_names refer to the polynomials_to_decompose without their exponents.

  • prefactor – string or sympy expression, optional; A factor that does not depend on the integration variables.
  • remainder_expression

    string or sympy expression or pySecDec.algebra._Expression, optional; An additional factor.

    Dummy function must be provided with all arguments, e.g. remainder_expression='exp(eps)*f(x0,x1)'. In addition, all dummy function must be listed in functions.

  • functions

    iterable of strings or sympy symbols, optional; Function symbols occuring in remainder_expression, e.g.``[‘f’]``.

    Note

    Only user-defined functions that are provided as c++-callable code should be mentioned here. Listing basic mathematical functions (e.g. log, pow, exp, sqrt, ...) is not required and considered an error to avoid name conflicts.

    Note

    The power function pow and the logarithm log use the nonstandard continuation with an infinitesimal negative imaginary part on the negative real axis (e.g. log(-1) = -i*pi).

  • real_parameters – iterable of strings or sympy symbols, optional; Symbols to be interpreted as real variables.
  • complex_parameters – iterable of strings or sympy symbols, optional; Symbols to be interpreted as complex variables.
  • form_optimization_level – integer out of the interval [0,3], optional; The optimization level to be used in FORM. Default: 2.
  • form_work_space – string, optional; The FORM WorkSpace. Default: '500M'.
  • form_insertion_depth – nonnegative integer, optional; How deep FORM should try to resolve nested function calls. Default: 5.
  • contour_deformation_polynomial – string or sympy symbol, optional; The name of the polynomial in polynomial_names that is to be continued to the complex plane according to a - i\delta prescription. For loop integrals, this is the second Symanzik polynomial F. If not provided, no code for contour deformation is created.
  • positive_polynomials – iterable of strings or sympy symbols, optional; The names of the polynomials in polynomial_names that should always have a positive real part. For loop integrals, this applies to the first Symanzik polynomial U. If not provided, no polynomial is checked for positiveness. If contour_deformation_polynomial is None, this parameter is ignored.
  • decomposition_method

    string, optional; The strategy to decompose the polynomials. The following strategies are available:

    • ‘iterative_no_primary’ (default)
    • ‘geometric_no_primary’
    • ‘iterative’
    • ‘geometric’
    • ‘geometric_ku’

    ‘iterative’, ‘geometric’, and ‘geometric_ku’ are only valid for loop integrals. An end user should always use ‘iterative_no_primary’ or ‘geometric_no_primary’ here. In order to compute loop integrals, please use the function pySecDec.loop_integral.loop_package().

  • normaliz_executable – string, optional; The command to run normaliz. normaliz is only required if decomposition_method starts with ‘geometric’. Default: ‘normaliz’
  • enforce_complex – bool, optional; Whether or not the generated integrand functions should have a complex return type even though they might be purely real. The return type of the integrands is automatically complex if contour_deformation is True or if there are complex_parameters. In other cases, the calculation can typically be kept purely real. Most commonly, this flag is needed if log(<negative real>) occurs in one of the integrand functions. However, pySecDec will suggest setting this flag to True in that case. Default: False
  • split – bool or integer, optional; Whether or not to split the integration domain in order to map singularities from 1 to 0. Set this option to True if you have singularties when one or more integration variables are one. If an integer is passed, that integer is used as seed to generate the splitting point. Default: False
  • ibp_power_goal

    integer, optional; The power_goal that is forwarded to integrate_by_parts().

    This option controls how the subtraction terms are generated. Setting it to -numpy.inf disables integrate_by_parts(), while 0 disables integrate_pole_part().

    See also

    To generate the subtraction terms, this function first calls integrate_by_parts() for each integration variable with the give ibp_power_goal. Then integrate_pole_part() is called.

    Default: -1

  • use_dreadnaut – bool or string, optional; Whether or not to use squash_symmetry_redundant_sectors_dreadnaut() to find sector symmetries. If given a string, interpret that string as the command line executable dreadnaut. If True, try $SECDEC_CONTRIB/bin/dreadnaut and, if the environment variable $SECDEC_CONTRIB is not set, dreadnaut. Default: True

5.7.2. Template Parser

Functions to generate c++ sources from template files.

pySecDec.code_writer.template_parser.parse_template_file(src, dest, replacements={})

Copy a file from src to dest replacing %(...) instructions in the standard python way.

Warning

If the file specified in dest exists, it is overwritten without prompt.

Parameters:
  • src – str; The path to the template file.
  • dest – str; The path to the destination file.
  • replacements

    dict; The replacements to be performed. The standard python replacement rules apply:

    >>> '%(var)s = %(value)i' % dict(
    ...     var = 'my_variable',
    ...     value = 5)
    'my_variable = 5'
    
pySecDec.code_writer.template_parser.parse_template_tree(src, dest, replacements_in_files={}, filesystem_replacements={})

Copy a directory tree from src to dest using parse_template_file() for each file and replacing the filenames according to filesystem_replacements.

Parameters:
  • src – str; The path to the template directory.
  • dest – str; The path to the destination directory.
  • replacements_in_files

    dict; The replacements to be performed in the files. The standard python replacement rules apply:

    >>> '%(var)s = %(value)i' % dict(
    ...     var = 'my_variable',
    ...     value = 5)
    'my_variable = 5'
    
  • filesystem_replacements – dict; Renaming rules for the destination files. and directories. If a file or directory name in the source tree src matches a key in this dictionary, it is renamed to the corresponding value. If the value is None, the corresponding file is ignored.

5.8. Generated C++ Libraries

A C++ Library to numerically compute a given integral (loop integral) can be generated by the make_package() (loop_package()) functions. The name passed to the make_package() or loop_package() function will be used as the C++ namespace of the generated library. A program demonstrating the use of the C++ library is generated for each integral and written to name/integrate_name.cpp. Here we document the C++ library API.

See also

C++ Interface

typedef double real_t

The real type used by the library.

typedef std::complex<real_t> complex_t

The complex type used by the library.

type integrand_return_t

The return type of the integrand function. If the integral has complex parameters or uses contour deformation or if enforce_complex is set to True in the call to make_package() or loop_package() then integrand_return_t is complex_t. Otherwise integrand_return_t is real_t.

template<typename T>
using nested_series_t = secdecutil::Series<secdecutil::Series<...<T>>>

A potentially nested secdecutil::Series representing the series expansion in each of the regulators. If the integral depends on only one regulator (for example, a loop integral generated with loop_package()) this type will be a secdecutil::Series. For integrals that depend on multiple regulators then this will be a series of series representing the multivariate series. This type can be used to write code that can handle integrals depending on arbitrarily many regulators.

typedef secdecutil::IntegrandContainer<integrand_return_t, real_t const *const > integrand_t

The type of the integrand. Within the generated C++ library integrands are stored in a container along with the number of integration variables upon which they depend. These containers can be passed to an integrator for numerical integration.

const unsigned int number_of_sectors

The number of sectors generated by the sector decomposition.

const unsigned int number_of_regulators

The number of regulators on which the integral depends.

const unsigned int number_of_real_parameters

The number of real parameters on which the integral depends.

const std::vector<std::string> names_of_real_parameters

An ordered vector of string representations of the names of the real parameters.

const unsigned int number_of_complex_parameters

The number of complex parameters on which the integral depends.

const std::vector<std::string> names_of_complex_parameters

An ordered vector of string representations of the names of the complex parameters.

const std::vector<int> lowest_orders

A vector of the lowest order of each regulator which appears in the integral, not including the prefactor.

const std::vector<int> highest_orders

A vector of the highest order of each regulator which appears in the integral, not including the prefactor. This depends on the requested_orders and prefactor/additional_prefactor parameter passed to make_package() or loop_package(). In the case of loop_package() it also depends on the \Gamma-function prefactor of the integral which appears upon Feynman parametrization.

const std::vector<int> lowest_prefactor_orders

A vector of the lowest order of each regulator which appears in the prefactor of the integral.

const std::vector<int> highest_prefactor_orders

A vector of the highest order of each regulator which appears in the prefactor of the integral.

const std::vector<int> requested_orders

A vector of the requested orders of each regulator used to generate the C++ library, i.e. the requested_orders parameter passed to make_package() or loop_package().

const std::vector<nested_series_t<sector_container_t>> sectors

A low level interface for obtaining the underlying integrand C++ functions.

Warning

The precise definition and usage of sectors is likely to change in future versions of pySecDec.

nested_series_t<integrand_return_t> prefactor(const std::vector<real_t> &real_parameters, const std::vector<complex_t> &complex_parameters)

The series expansion of the integral prefactor evaluated with the given parameters. If the library was generated using make_package() it will be equal to the prefactor passed to make_package(). If the library was generated with loop_package() it will be the product of the additional_prefactor passed to loop_package() and the \Gamma-function prefactor of the integral which appears upon Feynman parametrization.

const std::vector<std::vector<real_t>> pole_structures

A vector of the powers of the monomials that can be factored out of each sector of the polynomial during the decomposition.

Example: an integral depending on variables x and y may have two sectors, the first may have a monomial x^{-1} y^{-2} factored out and the second may have a monomial x^{-1} factored out during the decomposition. The resulting pole_structures would read { {-1,-2}, {-1,0} }. Poles of type x^{-1} are known as logarithmic poles, poles of type x^{-2} are known as linear poles.

std::vector<nested_series_t<integrand_t>> make_integrands(const std::vector<real_t> &real_parameters, const std::vector<complex_t> &complex_parameters)

(without contour deformation)

std::vector<nested_series_t<integrand_t>> make_integrands(const std::vector<real_t> &real_parameters, const std::vector<complex_t> &complex_parameters, unsigned number_of_presamples = 100000, real_t deformation_parameters_maximum = 1., real_t deformation_parameters_minimum = 1.e-5, real_t deformation_parameters_decrease_factor = 0.9)

(with contour deformation)

Gives a vector containing the series expansions of individual sectors of the integrand after sector decomposition with the specified real_paraemters and complex_parameters bound. Each element of the vector contains the series expansion of an individual sector. The series consists of instances of secdecutil::IntegrandContainer which contain the integrand functions and the number of integration variables upon which they depend. The real and complex parameters are bound to the values passed in real_parameters and complex_parameters. If enabled, contour deformation is controlled by the parameters number_of_presamples, deformation_parameters_maximum, deformation_parameters_minimum, deformation_parameters_decrease_factor which are documented in pySecDec.integral_interface.IntegralLibrary.

Passing the integrand_t to the secdecutil::Integrator::integrate() function of an instance of a particular secdecutil::Integrator will return the numerically evaluated integral. To integrate all orders of all sectors secdecutil::deep_apply() can be used.

Note

This is the recommended way to access the integrand functions.

5.9. Integral Interface

An interface to libraries generated by pySecDec.code_writer.make_package() or pySecDec.loop_integral.loop_package().

class pySecDec.integral_interface.CPPIntegrator

Abstract base class for integrators to be used with an IntegralLibrary. This class holds a pointer to the c++ integrator and defines the destructor.

class pySecDec.integral_interface.Cuhre(integral_library, epsrel=0.01, epsabs=1e-07, flags=0, mineval=0, maxeval=1000000, key=0, real_complex_together=False)

Wrapper for the Cuhre integrator defined in the cuba library.

Parameters:integral_libraryIntegralLibrary; The integral to be computed with this integrator.

The other options are defined in the cuba manual.

class pySecDec.integral_interface.Divonne(integral_library, epsrel=0.01, epsabs=1e-07, flags=0, seed=0, mineval=0, maxeval=1000000, key1=2000, key2=1, key3=1, maxpass=4, border=0.0, maxchisq=1.0, mindeviation=0.15, real_complex_together=False)

Wrapper for the Divonne integrator defined in the cuba library.

Parameters:integral_libraryIntegralLibrary; The integral to be computed with this integrator.

The other options are defined in the cuba manual.

class pySecDec.integral_interface.IntegralLibrary(shared_object_path)

Interface to a c++ library produced by make_package() or loop_package().

Parameters:shared_object_path

str; The path to the file “<name>_pylink.so” that can be built by the command

$ make pylink

in the root directory of the c++ library.

Instances of this class can be called with the following arguments:

Parameters:
  • real_parameters – iterable of float; The real_parameters of the library.
  • complex_parameters – iterable of complex; The complex parameters of the library.
  • together – bool; Whether to integrate the sum of all sectors or to integrate the sectors separately. Default: True.
  • number_of_presamples – unsigned int, optional; The number of samples used for the comtour optimization. This option is ignored if the integral library was created without deformation. Default: 100000.
  • deformation_parameters_maximum – float, optional; The maximal value the deformation parameters \lambda_i can obtain. If number_of_presamples=0, all \lambda_i are set to this value. This option is ignored if the integral library was created without deformation. Default: 1.0.
  • deformation_parameters_minimum – float, optional; This option is ignored if the integral library was created without deformation. Default: 1e-5.
  • deformation_parameters_decrease_factor – float, optional; If the sign check with the optimized \lambda_i fails, all \lambda_i are multiplied by this value until the sign check passes. This option is ignored if the integral library was created without deformation. Default: 0.9.

The call operator returns three strings: * The integral without its prefactor * The prefactor * The integral multiplied by the prefactor

The integrator cen be configured by calling the member methods use_Vegas(), use_Suave(), use_Divonne(), and use_Cuhre(). The available options are listed in the documentation of Vegas, Suave, Divonne, and Cuhre, respectively. If not specified otherwise, Vegas is used with its default arguments. For details about the options, refer to the cuba manual.

Further information about the library is stored in the member variable info of type dict.

class pySecDec.integral_interface.Suave(integral_library, epsrel=0.01, epsabs=1e-07, flags=0, seed=0, mineval=0, maxeval=1000000, nnew=1000, nmin=10, flatness=25.0, real_complex_together=False)

Wrapper for the Suave integrator defined in the cuba library.

Parameters:integral_libraryIntegralLibrary; The integral to be computed with this integrator.

The other options are defined in the cuba manual.

class pySecDec.integral_interface.Vegas(integral_library, epsrel=0.01, epsabs=1e-07, flags=0, seed=0, mineval=0, maxeval=1000000, nstart=1000, nincrease=500, nbatch=1000, real_complex_together=False)

Wrapper for the Vegas integrator defined in the cuba library.

Parameters:integral_libraryIntegralLibrary; The integral to be computed with this integrator.

The other options are defined in the cuba manual.

5.10. Miscellaneous

Collection of general-purpose helper functions.

pySecDec.misc.adjugate(M)

Calculate the adjugate of a matrix.

Parameters:M – a square-matrix-like array;
pySecDec.misc.all_pairs(iterable)

Return all possible pairs of a given set. all_pairs([1,2,3,4]) --> [(1,2),(3,4)] [(1,3),(2,4)] [(1,4),(2,3)]

Parameters:iterable – iterable; The set to be split into all possible pairs.
pySecDec.misc.argsort_2D_array(array)

Sort a 2D array according to its row entries. The idea is to bring identical rows together.

See also

If your array is not two dimesional use argsort_ND_array().

Example:
input   sorted
1 2 3   1 2 3
2 3 4   1 2 3
1 2 3   2 3 4

Return the indices like numpy’s argsort() would.

Parameters:array – 2D array; The array to be argsorted.
pySecDec.misc.argsort_ND_array(array)

Like argsort_2D_array(), this function groups identical entries in an array with any dimensionality greater than (or equal to) two together.

Return the indices like numpy’s argsort() would.

Parameters:array – ND array, N>=2; The array to be argsorted.
pySecDec.misc.assert_degree_at_most_max_degree(expression, variables, max_degree, error_message)

Assert that expression is a polynomial of degree less or equal max_degree in the variables.

pySecDec.misc.cached_property(method)

Like the builtin property to be used as decorator but the method is only called once per instance.

Example:

class C(object):
    'Sum up the numbers from one to `N`.'
    def __init__(self, N):
        self.N = N
    @cached_property
    def sum(self):
        result = 0
        for i in range(1, self.N + 1):
            result += i
        return result
pySecDec.misc.det(M)

Calculate the determinant of a matrix.

Parameters:M – a square-matrix-like array;
pySecDec.misc.doc(docstring)

Decorator that replaces a function’s docstring with docstring.

Example:

@doc('documentation of `some_funcion`')
def some_function(*args, **kwargs):
    pass
pySecDec.misc.lowest_order(expression, variable)

Find the lowest order of expression‘s series expansion in variable.

Example:

>>> from pySecDec.misc import lowest_order
>>> lowest_order('exp(eps)', 'eps')
0
>>> lowest_order('gamma(eps)', 'eps')
-1
Parameters:
  • expression – string or sympy expression; The expression to compute the lowest expansion order of.
  • variable – string or sympy expression; The variable in which to expand.
pySecDec.misc.missing(full, part)

Return the elements in full that are not contained in part. Raise ValueError if an element is in part but not in full. missing([1,2,3], [1]) --> [2,3] missing([1,2,3,1], [1,2]) --> [3,1] missing([1,2,3], [1,'a']) --> ValueError

Parameters:
  • full – iterable; The set of elements to complete part with.
  • part – iterable; The set to be completed to a superset of full.
pySecDec.misc.powerset(iterable, min_length=0, stride=1)

Return an iterator over the powerset of a given set. powerset([1,2,3]) --> () (1,) (2,) (3,) (1,2) (1,3) (2,3) (1,2,3)

Parameters:
  • iterable – iterable; The set to generate the powerset for.
  • min_length – integer, optional; Only generate sets with minimal given length. Default: 0.
  • stride – integer; Only generate sets that have a multiple of stride elements. powerset([1,2,3], stride=2) --> () (1,2) (1,3) (2,3)
pySecDec.misc.rangecomb(low, high)

Return an iterator over the occuring orders in a multivariate series expansion between low and high.

Parameters:
  • low – vector-like array; The lowest orders.
  • high – vector-like array; The highest orders.

Example:

>>> from pySecDec.misc import rangecomb
>>> all_orders = rangecomb([-1,-2], [0,0])
>>> list(all_orders)
[(-1, -2), (-1, -1), (-1, 0), (0, -2), (0, -1), (0, 0)]
pySecDec.misc.sympify_symbols(iterable, error_message, allow_number=False)

sympify each item in iterable and assert that it is a symbol.