seemps.analysis.integration.integrate_mps#

seemps.analysis.integration.integrate_mps(mps, domain, mps_order='A')[source]#

Returns the integral of a multivariate function represented as a MPS. Uses the ‘best possible’ quadrature rule according to the intervals that compose the mesh. Intervals of type RegularInterval employ high-order Newton-Côtes rules, while those of type ChebyshevInterval employ Clenshaw-Curtis rules.

Parameters:
mpsMPS

The MPS representation of the multivariate function to be integrated.

domainUnion[Interval, Mesh]

An object defining the discretization domain of the function. Can be either an Interval or a Mesh given by a collection of intervals. The quadrature rules are selected based on the properties of these intervals.

mps_orderstr, default=’A’

Specifies the ordering of the qubits in the quadrature. Possible values:. - ‘A’: Qubits are serially ordered (by variable). - ‘B’: Qubits are interleaved (by significance).

Returns:
Union[complex, float]

The integral of the MPS representation of the function discretized in the given Mesh.

Notes

  • This algorithm assumes that all function variables are in the standard MPS form, i.e.

quantized in base 2, and are discretized either on a RegularInterval or `ChebyshevInterval.

  • For more general structures, the quadrature MPS can be constructed

using the univariate quadrature rules and the mps_tensor_product routine, which can be subsequently contracted using the scprod routine.

Examples

# Integrate a given bivariate function using the Clenshaw-Curtis quadrature.
# Assuming that the MPS is already loaded (for example, using TT-Cross or Chebyshev).
mps_function_2d = ...

# Define a domain that matches the MPS to integrate.
start, stop = -1, 1
n_qubits = 10
interval = ChebyshevInterval(-1, 1, 2**n_qubits, endpoints=True)
mesh = Mesh([interval, interval])

# Integrate the MPS on the given discretization domain.
integral = integrate_mps(mps_function_2d, mesh)