Diagonal

Propagators for diagonal operators in the computational basis. These are commonly used for phase/cost operators.

Propagator Class

class quop_mpi.propagator.diagonal.unitary(*args, **kwargs)

Bases: Unitary

Compute the action of a mixing unitary with a phase_shift operator or a sequence of mixing-unitaries with phase_shift operators (see the unitary_n_params attribute below).

Inheritance Diagram:

digraph "sphinx-ext-graphviz" {
    rankdir="LR";
    node [fontsize="10"];
    Unitary[label="quop_mpi.Unitary", shape="rectangle"];
    unitary[label="quop_mpi.propagator.phase_shift.unitary", shape="rectangle"];

    Unitary -> unitary;
}

See quop_mpi.Unitary.

Attributes:
unitary_type

'phase_shift'

planner

false

unitary_n_params

Set on initialisation to 1 or more. If unitary_n_parameters > 1, the Operator Function must return a list[csr_matrix] of length unitary_n_parameters containing csr_matrix partitions of of local_i rows.

propagate(gammas)

Simulation of the action of a :term`unitary`.

When implemented, propagate contains a call to a method (typically a contained in a complied Python extension module) that takes the class attributes initial_state, final_state and MPI_COMM, together with attributes describing the parallel partitioning scheme and variational parameters x, as input. The action of the unitary is computed in MPI parallel, with the computed result written to final_state.

Warning

Not implemented by the base Unitary class.

Parameters:
xndarray[float64]

a 1-D real array of n_params variational parameters

Examples

def propagate(self, x):

    external_propagator(
        x, self.partition_table, self.initial_state,
        self.final_state, self.MPI_COMM )

Operators

Pre-defined operator functions for diagonal propagators.

Predefined Operator Functions and related utility for quop_mpi.propagator.diagonal.unitary.

Return Format

An Operator Function for 'diagonal' unitary instances returns:

diagonalndarray[float64] or list[ndarray[float64]]

A 1-D real array of size local_i containing the diagonal elements of the operator for global indices local_i_offset to local_i_offset + local_i - 1.

Alternatively, a list of such arrays for multiple diagonal operators (requires unitary_n_params to match the list length).

If the Operator function returns list[ndarray[float64]], the unitary instance must be initialised with unitary_n_parameters equal to the length of returned list. The resulting unitary is then equivalent to a sequence of phase-shift unitaries with independently parameterised unitary parameters.

Propagation Method

The diagonal propagator computes \(e^{-itH}|\psi\rangle\) via direct element-wise multiplication:

\[|\psi'\rangle_k = e^{-it H_{kk}} |\psi\rangle_k\]

where \(H_{kk}\) are the diagonal elements. This is the most efficient propagation method as it requires only \(O(N)\) operations with no communication between MPI ranks.

quop_mpi.propagator.diagonal.operator.array(partition_table: list[int], MPI_COMM: Intracomm, array: np.ndarray[np.float64]) np.ndarray[np.float64]

Define the diagonal of the phase-shift unitary using a Numpy array.

An Operator Function for the quop_mpi.propagate.diagonal.unitary class.

Note

For memory efficiency, array can be present as an ndarray[float64] at MPI_COMM.rank == 0 only and None at all other ranks in MPI_COMM.

Parameters:
partition_tablelist[int]

describes the parallel partitioning scheme, quop_mpi.Ansatz attribute

MPI_COMMIntracomm

MPI communicator, quop_mpi.Ansatz attribute

arraynp.ndarray[np.float64]

a 1-D real array of size system size

Returns:
ndarray[float64]

a 1-D real array containing local_i elements of the operator diagonal with global index offset local_i_offset

quop_mpi.propagator.diagonal.operator.cartesian(system_size: int, local_i: int, local_i_offset: int, Ns: list[int], deltas: list[float], mins: list[float], function: Callable, *args, **kwargs) np.ndarray[np.float64]

TODO:UPDATE Generate the diagonal of a phase-shift unitary operator using a Python function defined in discrete Cartesian coordinates.

An Observables Function. Depending on wether QVA simulation is defined using the Ansatz class directly or with a predefined Ansatz subclass from the algorithm submodule, the following arguments must be defined in a corresponding FunctionDict on initialisation of the unitary instance:

Additional positional and keyword arguments in the FunctionDict are passed to function.

The function argument must conform to the signature,

def function(x: ndarray[float64], *args, *kwargs) -> float

where x is a 1-D array containing a len(Ns) -dimensional grid point.

Parameters:
system_sizeint

size of the simulated QVA

local_iint

size of the local system state partition, quop_mpi.Unitary attribute

local_i_offsetint

global index offset of the local system state partition, quop_mpi.Unitary attribute

Nslist[int]

the number of qubits assigned to each dimension of the cartesian grid such that there is 2 ** Ns[d] grid points per dimension d

deltaslist[float]

step size in each Cartesian coordinate

minslist[float]

lower bound of each Cartesian coordinate

functionCallable

a Python function that takes a list of len(Ns) real coordinate values and returns a float

Returns:
ndarray[float64]

a 1-D real array containing local_i elements of the operator diagonal with global index offset local_i_offset

See also

setup_cartesian

compute deltas and mins

cartesian_scaled

alternative to cartesian, scales function between 0 and an upper bound.

quop_mpi.propagator.diagonal.operator.cartesian_scaled(system_size: int, local_i: int, local_i_offset: int, MPI_COMM: Intracomm, Ns: list[int], deltas: list[float], mins: list[float], function: Callable, coeff: float, *args, **kwargs) np.ndarray[np.float64]

Generate the diagonal of a phase-shift unitary operator using a Python function defined in discrete Cartesian coordinates with the function scaled between 0 and coeff.

An Observables Function. Depending on wether QVA simulation is defined using the Ansatz class directly or with a predefined Ansatz subclass from the algorithm submodule, the following arguments must be defined in a corresponding FunctionDict on initialisation of the unitary instance:

Additional positional and keyword arguments in the FunctionDict are passed to function.

The function argument must conform to the signature,

def function(x: ndarray[float64], *args, *kwargs) -> float

where x is a 1-D array containing a len(Ns) -dimensional grid point.

Parameters:
system_sizeint

size of the simulated QVA

local_iint

size of the local system state partition, quop_mpi.Unitary attribute

local_i_offsetint

global index offset of the local system state partition, quop_mpi.Unitary attribute

Nslist[int]

the number of qubits assigned to each dimension of the cartesian grid such that there is 2 ** Ns[d] grid points per dimension d

deltaslist[float]

step size in each Cartesian coordinate

minslist[float]

lower bound of each Cartesian coordinate

functionCallable

a Python function that takes a list of len(Ns) real coordinate values and returns a float

coefffloat

a positive real number, the upper bound of the scaling range

Returns:
ndarray[float64]

a 1-D real array containing local_i elements of the operator diagonal with global index offset local_i_offset

See also

setup_cartesian

compute deltas and mins

cartesian_scaled

alternative to cartesian_scaled, does not scale function

quop_mpi.propagator.diagonal.operator.csv(partition_table: list[int], MPI_COMM: Intracomm, filename: Callable, *args, **kwargs) np.ndarray[np.float64]

Load the diagonal of a phase-shift unitary using pandas.

An Operator Function for the quop_mpi.propagate.diagonal.unitary class. The filename argument must be defined in a corresponding FunctionDict on initialisation of the unitary instance. Additional keyword arguments in the FunctionDict are passed to the pandas.read_csv method.

Parameters:
partition_tablelist[int]

describes the parallel partitioning scheme, quop_mpi.Ansatz attribute

MPI_COMMIntracomm

MPI communicator, quop_mpi.Ansatz attribute

filenameCallable

path to a *.csv file

Returns:
ndarray[float64]

a 1-D real array containing local_i elements of the operator diagonal with global index offset local_i_offset

quop_mpi.propagator.diagonal.operator.hdf5(partition_table: int, MPI_COMM: Intracomm, filename: str, dataset_name: str) np.ndarray[np.float64]

Load the diagonal of a phase-shift unitary using HDF5 for Python.

An Operator Function for the quop_mpi.propagate.diagonal.unitary class. The filename and dataset_name arguments must be defined in a corresponding FunctionDict on initialisation of the unitary instance. Additional positional and keyword arguments in the FunctionDict are passed to the h5py.File method.

Parameters:
partition_tableint

describes the parallel partitioning scheme, quop_mpi.Ansatz attribute

MPI_COMMIntracomm

MPI communicator, quop_mpi.Ansatz attribute

filenamestr

path to a *.h5 file

dataset_namestr

path to the dataset containing a ndarray[float64] of size system size

Returns:
np.ndarray[np.float64]

a 1-D real array containing local_i elements of the operator diagonal with global index offset local_i_offset

quop_mpi.propagator.diagonal.operator.serial(partition_table: list[int], MPI_COMM: Intracomm, variational_parameters: np.ndarray[np.float64], function: Callable, *args, **kwargs) np.ndarray[np.float64] | list[np.ndarray[np.float64]]

Generate the diagonal of the operator for one or more sequential phase-shift unitaries using a serial Python function.

An Operator Function for the quop_mpi.propagator.diagonal.unitary class. The function argument must be defined in a corresponding FunctionDict on initialisation of the unitary instance. Additional positional and keyword arguments contained in the FunctionDict are passed to function.

The function argument must conform to the signature,

def function(*args, *kwargs) -> (ndarray[float64] | list[ndarray[float64]])

where the output is a 1-D real array of type ndarray[float64] and length system size, or list containing one or more 1-D real arrays of type ndarray[float64] and length system_size.

Parameters:
partition_tablelist[int]

describes the parallel partitioning scheme, quop_mpi.Ansatz attribute

MPI_COMMIntracomm

MPI communicator, quop_mpi.Ansatz attribute

variational_parametersndarray[float64]

operator parameters, passed to function if unitary.operator_n_params > 0, quop_mpi.Unitary attribute

functionCallable

returns one or more ndarray[float64] of size system size corresponding to the diagonal of the operator of one or more phase-shift unitaries

Returns:
ndarray[float64] or list[ndarray[float64]]

a 1-d real array or list of 1-D real arrays containing a local_i elements of the operator diagonal with global index offset local_i_offset

quop_mpi.propagator.diagonal.operator.setup_cartesian(Ns: list[int], bounds: list[list[float]]) list[list[float]]

Compute the step-size and minimum coordinate values in each dimension of a Cartesian grid.

Parameters:
Nslist[int]

the number of qubits assigned to each dimension of the Cartesian grid such that there is 2 ** Ns[d] grid points per dimension d

boundslist[list[float]]

the lower and upper bounds of each dimension where len(Ns) == len(bounds)

Returns:
list[list[float]]

the step-size, deltas, and the minimum, mins, in each Cartesian coordinate