ffsim.qiskit¶
Code that uses Qiskit, e.g. for constructing quantum circuits.
- class ffsim.qiskit.DiagCoulombEvolutionJW(norb, mat, time, *, z_representation=False, label=None)[source]¶
Bases:
Gate
Diagonal Coulomb evolution under the Jordan-Wigner transformation.
The diagonal Coulomb evolution gate has the unitary
\[\exp\left(-i t \sum_{\sigma, \tau, i, j} Z^{(\sigma \tau)}_{ij} n_{\sigma, i} n_{\tau, j} / 2\right)\]where \(n_{\sigma, i}\) denotes the number operator on orbital \(i\) with spin \(\sigma\), \(Z^{(\sigma \tau)}\) is a real-valued matrix.
This gate assumes that qubits are ordered such that the first norb qubits correspond to the alpha orbitals and the last norb qubits correspond to the beta orbitals.
- __init__(norb, mat, time, *, z_representation=False, label=None)[source]¶
Create new diagonal Coulomb evolution gate.
- Parameters:
norb (
int
) – The number of spatial orbitals.mat (
ndarray
|tuple
[ndarray
|None
,ndarray
|None
,ndarray
|None
]) – The diagonal Coulomb matrix \(Z\). You can pass either a single Numpy array specifying the coefficients to use for all spin interactions, or you can pass a tuple of three Numpy arrays specifying independent coefficients for alpha-alpha, alpha-beta, and beta-beta interactions (in that order). If passing a tuple, you can set a tuple element toNone
to indicate the absence of interactions of that type. The alpha-alpha and beta-beta matrices are assumed to be symmetric, and only their upper triangular entries are used.time (
float
) – The evolution time.z_representation (
bool
) – Whether the input matrices are in the “Z” representation.label (
str
|None
) – The label of the gate.
- class ffsim.qiskit.DiagCoulombEvolutionSpinlessJW(norb, mat, time, *, label=None)[source]¶
Bases:
Gate
Spinless diagonal Coulomb evolution under the Jordan-Wigner transformation.
The spinless diagonal Coulomb evolution gate has the unitary
\[\exp\left(-i t \sum_{i, j} Z^{ij} n_i n_j / 2\right)\]where \(n_i\) denotes the number operator on orbital \(i\) and \(Z\) is a real symmetric matrix.
- __init__(norb, mat, time, *, label=None)[source]¶
Create new diagonal Coulomb evolution gate.
- Parameters:
norb (
int
) – The number of spatial orbitals.mat (
ndarray
) – The diagonal Coulomb matrix \(Z\). It is assumed to be symmetric, and only its upper triangular entries are used.time (
float
) – The evolution time.label (
str
|None
) – The label of the gate.
- class ffsim.qiskit.DropNegligible(atol=1e-08)[source]¶
Bases:
TransformationPass
Drop gates with negligible effects.
- class ffsim.qiskit.FfsimSampler(*, default_shots=1024, global_depolarizing=0.0, seed=None)[source]¶
Bases:
BaseSamplerV2
Implementation of the Qiskit Sampler primitive backed by ffsim.
- __init__(*, default_shots=1024, global_depolarizing=0.0, seed=None)[source]¶
Initialize the ffsim sampler.
- Parameters:
default_shots (
int
) – The default shots to use if not specified during run.global_depolarizing (
float
) – Depolarizing probability for a noisy simulation. Specifies the probability of sampling from the uniform distribution instead of the state vector.seed (
Generator
|int
|None
) – A seed to initialize the pseudorandom number generator. Should be a valid input tonp.random.default_rng
.
- run(pubs, *, shots=None)[source]¶
Run and collect samples from each pub.
- Parameters:
pubs (
Iterable
[Union
[QuantumCircuit
,Tuple
[QuantumCircuit
],Tuple
[QuantumCircuit
,Mapping
[Union
[Parameter
,str
,Tuple
[Union
[Parameter
,str
],...
]],Union
[Buffer
,_SupportsArray
[dtype
[Any
]],_NestedSequence
[_SupportsArray
[dtype
[Any
]]],bool
,int
,float
,complex
,str
,bytes
,_NestedSequence
[bool
|int
|float
|complex
|str
|bytes
]]]],Tuple
[QuantumCircuit
,Mapping
[Union
[Parameter
,str
,Tuple
[Union
[Parameter
,str
],...
]],Union
[Buffer
,_SupportsArray
[dtype
[Any
]],_NestedSequence
[_SupportsArray
[dtype
[Any
]]],bool
,int
,float
,complex
,str
,bytes
,_NestedSequence
[bool
|int
|float
|complex
|str
|bytes
]]],Optional
[Integral
]]]]) – An iterable of pub-like objects. For example, a list of circuits or tuples(circuit, parameter_values)
.shots (
int
|None
) – The total number of shots to sample for each sampler pub that does not specify its own shots. IfNone
, the primitive’s default shots value will be used, which can vary by implementation.
- Return type:
PrimitiveJob
[PrimitiveResult
[SamplerPubResult
]]- Returns:
The job object of Sampler’s result.
- class ffsim.qiskit.GivensAnsatzOpJW(givens_ansatz_op, *, label=None)[source]¶
Bases:
Gate
Givens rotation ansatz operator under the Jordan-Wigner transformation.
See
ffsim.GivensAnsatzOp
for a description of this gate’s unitary.- __init__(givens_ansatz_op, *, label=None)[source]¶
Create a new Givens ansatz operator gate.
- Parameters:
givens_ansatz_op (
GivensAnsatzOp
) – The Givens rotation ansatz operator.label (
str
|None
) – The label of the gate.
- class ffsim.qiskit.GivensAnsatzOpSpinlessJW(givens_ansatz_op, *, label=None)[source]¶
Bases:
Gate
Spinless Givens rotation ansatz operator under the Jordan-Wigner transformation.
Like
GivensAnsatzOpJW
but only acts on a single spin species.- __init__(givens_ansatz_op, *, label=None)[source]¶
Create a new Givens ansatz operator gate.
- Parameters:
givens_ansatz_op (
GivensAnsatzOp
) – The Givens rotation ansatz operator.label (
str
|None
) – The label of the gate.
- class ffsim.qiskit.MergeOrbitalRotations(*args, **kwargs)[source]¶
Bases:
TransformationPass
Merge consecutive orbital rotation gates.
- class ffsim.qiskit.NumNumAnsatzOpSpinBalancedJW(num_num_ansatz_op, *, label=None)[source]¶
Bases:
Gate
Spin-balanced number-number ansatz under the Jordan-Wigner transformation.
See
NumNumAnsatzOpSpinBalanced
for a description of this gate’s unitary.This gate assumes that qubits are ordered such that the first norb qubits correspond to the alpha orbitals and the last norb qubits correspond to the beta orbitals.
- __init__(num_num_ansatz_op, *, label=None)[source]¶
Create a new number-number ansatz operator gate.
- Parameters:
num_num_ansatz_op (
NumNumAnsatzOpSpinBalanced
) – The number-number ansatz operator.label (
str
|None
) – The label of the gate.
- class ffsim.qiskit.NumOpSumEvolutionJW(norb, coeffs, time, *, label=None)[source]¶
Bases:
Gate
Number operator sum evolution under the Jordan-Wigner transformation.
The number operator sum evolution gate has the unitary
\[\exp\left(-i t \sum_{\sigma, i} \lambda^{(\sigma)}_i n_{\sigma, i}\right)\]where \(n_{\sigma, i}\) denotes the number operator on orbital \(i\) with spin \(\sigma\) and the \(\lambda_i\) are real numbers.
This gate assumes that qubits are ordered such that the first norb qubits correspond to the alpha orbitals and the last norb qubits correspond to the beta orbitals.
- __init__(norb, coeffs, time, *, label=None)[source]¶
Create new number operator sum evolution gate.
- Parameters:
norb (
int
) – The number of spatial orbitals.coeffs (
ndarray
|tuple
[ndarray
|None
,ndarray
|None
]) – The coefficients of the linear combination. You can pass either a single Numpy array specifying the coefficients to apply to both spin sectors, or you can pass a pair of Numpy arrays specifying independent coefficients for spin alpha and spin beta. If passing a pair, you can useNone
for one of the values in the pair to indicate that no operation should be applied to that spin sector.time (
float
) – The evolution time.label (
str
|None
) – The label of the gate.
- class ffsim.qiskit.NumOpSumEvolutionSpinlessJW(norb, coeffs, time, *, label=None)[source]¶
Bases:
Gate
Spinless number operator sum evolution under the Jordan-Wigner transformation.
The spinless number operator sum evolution gate has the unitary
\[\exp\left(-i t \sum_{i} \lambda_i n_{i}\right)\]where \(n_i\) denotes the number operator on orbital \(i\) and the \(\lambda_i\) are real numbers.
- class ffsim.qiskit.OrbitalRotationJW(norb, orbital_rotation, *, label=None, validate=True, rtol=1e-05, atol=1e-08)[source]¶
Bases:
Gate
Orbital rotation under the Jordan-Wigner transformation.
An orbital rotation maps creation operators as
\[a^\dagger_{\sigma, i} \mapsto \sum_{j} U_{ji} a^\dagger_{\sigma, j}\]where \(U\) is a unitary matrix. This is equivalent to applying the transformation given by
\[\prod_{\sigma} \exp\left(\sum_{ij} \log(U)_{ij} a^\dagger_{\sigma, i} a_{\sigma, j}\right)\]This gate assumes that qubits are ordered such that the first norb qubits correspond to the alpha orbitals and the last norb qubits correspond to the beta orbitals.
- __init__(norb, orbital_rotation, *, label=None, validate=True, rtol=1e-05, atol=1e-08)[source]¶
Create new orbital rotation gate.
- Parameters:
norb (
int
) – The number of spatial orbitals.orbital_rotation (
ndarray
|tuple
[ndarray
|None
,ndarray
|None
]) – The orbital rotation. You can pass either a single Numpy array specifying the orbital rotation to apply to both spin sectors, or you can pass a pair of Numpy arrays specifying independent orbital rotations for spin alpha and spin beta. If passing a pair, you can useNone
for one of the values in the pair to indicate that no operation should be applied to that spin sector.label (
str
|None
) – The label of the gate.validate (
bool
) – Whether to check that the input orbital rotation(s) is unitary and raise an error if it isn’t.rtol (
float
) – Relative numerical tolerance for input validation.atol (
float
) – Absolute numerical tolerance for input validation.
- Raises:
ValueError – The input matrix is not unitary.
- class ffsim.qiskit.OrbitalRotationSpinlessJW(norb, orbital_rotation, *, label=None, validate=True, rtol=1e-05, atol=1e-08)[source]¶
Bases:
Gate
Orbital rotation under the Jordan-Wigner transformation, spinless version.
Like
OrbitalRotationJW
but only acts on a single spin species.- __init__(norb, orbital_rotation, *, label=None, validate=True, rtol=1e-05, atol=1e-08)[source]¶
Create new orbital rotation gate.
- Parameters:
norb (
int
) – The number of spatial orbitals.orbital_rotation (
ndarray
) – The orbital rotation.label (
str
|None
) – The label of the gate.validate (
bool
) – Whether to check that the input orbital rotation(s) is unitary and raise an error if it isn’t.rtol (
float
) – Relative numerical tolerance for input validation.atol (
float
) – Absolute numerical tolerance for input validation.
- Raises:
ValueError – The input matrix is not unitary.
- ffsim.qiskit.PRE_INIT = <qiskit.transpiler.passmanager.PassManager object>¶
Pass manager recommended for the Qiskit transpiler
pre_init
stage.See
pre_init_passes()
for a description of the transpiler passes included in this pass manager.
- class ffsim.qiskit.PrepareHartreeFockJW(norb, nelec, label=None)[source]¶
Bases:
Gate
Gate that prepares the Hartree-Fock state (under JWT) from the all zeros state.
This gate assumes the Jordan-Wigner transformation (JWT).
This gate is meant to be applied to the all zeros state. It decomposes simply as a sequence of X gates that prepares the Hartree-Fock electronic configuration.
This gate assumes that qubits are ordered such that the first norb qubits correspond to the alpha orbitals and the last norb qubits correspond to the beta orbitals.
- class ffsim.qiskit.PrepareHartreeFockSpinlessJW(norb, nelec, label=None)[source]¶
Bases:
Gate
Prepare the Hartree-Fock state (under JWT) from the zero state, spinless.
Like
PrepareHartreeFockJW
but only acts on a single spin species.
- class ffsim.qiskit.PrepareSlaterDeterminantJW(norb, occupied_orbitals, orbital_rotation=None, *, label=None, validate=True, rtol=1e-05, atol=1e-08)[source]¶
Bases:
Gate
Gate that prepares a Slater determinant (under JWT) from the all zeros state.
This gate assumes the Jordan-Wigner transformation (JWT).
A Slater determinant is a state of the form
\[\mathcal{U} \lvert x \rangle,\]where \(\mathcal{U}\) is an orbital rotation and \(\lvert x \rangle\) is an electronic configuration (computational basis state). The reason this gate exists (when
OrbitalRotationJW
already exists) is that the preparation of a Slater determinant has a more optimized circuit than a generic orbital rotation.This gate is meant to be applied to the all zeros state. Its behavior when applied to any other state is not guaranteed. The global phase of the prepared state may be arbitrary.
This gate assumes that qubits are ordered such that the first norb qubits correspond to the alpha orbitals and the last norb qubits correspond to the beta orbitals.
Reference: arXiv:1711.05395
- __init__(norb, occupied_orbitals, orbital_rotation=None, *, label=None, validate=True, rtol=1e-05, atol=1e-08)[source]¶
Create new Slater determinant preparation gate.
- Parameters:
norb (
int
) – The number of spatial orbitals.occupied_orbitals (
tuple
[Sequence
[int
],Sequence
[int
]]) – The occupied orbitals in the electonic configuration. This is a pair of lists of integers, where the first list specifies the spin alpha orbitals and the second list specifies the spin beta orbitals.orbital_rotation (
ndarray
|tuple
[ndarray
|None
,ndarray
|None
] |None
) – The optional orbital rotation. You can pass either a single Numpy array specifying the orbital rotation to apply to both spin sectors, or you can pass a pair of Numpy arrays specifying independent orbital rotations for spin alpha and spin beta. If passing a pair, you can useNone
for one of the values in the pair to indicate that no operation should be applied to that spin sector.label (
str
|None
) – The label of the gate.validate (
bool
) – Whether to check that the input orbital rotation(s) is unitary and raise an error if it isn’t.rtol (
float
) – Relative numerical tolerance for input validation.atol (
float
) – Absolute numerical tolerance for input validation.
- Raises:
ValueError – The input orbital rotation matrix is not unitary.
- class ffsim.qiskit.PrepareSlaterDeterminantSpinlessJW(norb, occupied_orbitals, orbital_rotation=None, *, label=None, validate=True, rtol=1e-05, atol=1e-08)[source]¶
Bases:
Gate
Prepare a Slater determinant (under JWT) from the zero state, spinless version.
Like
PrepareSlaterDeterminantJW
but only acts on a single spin species.- __init__(norb, occupied_orbitals, orbital_rotation=None, *, label=None, validate=True, rtol=1e-05, atol=1e-08)[source]¶
Create new Slater determinant preparation gate.
- Parameters:
norb (
int
) – The number of spatial orbitals.occupied_orbitals (
Sequence
[int
]) – The occupied orbitals in the electonic configuration.orbital_rotation (
ndarray
|None
) – The optional orbital rotation.label (
str
|None
) – The label of the gate.validate (
bool
) – Whether to check that the input orbital rotation(s) is unitary and raise an error if it isn’t.rtol (
float
) – Relative numerical tolerance for input validation.atol (
float
) – Absolute numerical tolerance for input validation.
- Raises:
ValueError – The input orbital rotation matrix is not unitary.
- class ffsim.qiskit.SimulateTrotterDiagCoulombSplitOpJW(hamiltonian, time, *, n_steps=1, order=0, label=None)[source]¶
Bases:
Gate
Split operator Trotter evolution of diagonal Coulomb Hamiltonian, Jordan-Wigner.
This gate assumes that qubits are ordered such that the first norb qubits correspond to the alpha orbitals and the last norb qubits correspond to the beta orbitals.
- __init__(hamiltonian, time, *, n_steps=1, order=0, label=None)[source]¶
Create diagonal Coulomb split-operator Trotter evolution gate.
- Parameters:
norb – The number of spatial orbitals.
hamiltonian (
DiagonalCoulombHamiltonian
) – The Hamiltonian.time (
float
) – The evolution time.n_steps (
int
) – The number of Trotter steps.order (
int
) – The order of the Trotter decomposition.label (
str
|None
) – The label of the gate.
- class ffsim.qiskit.SimulateTrotterDoubleFactorizedJW(hamiltonian, time, *, n_steps=1, order=0, label=None)[source]¶
Bases:
Gate
Trotter time evolution of double-factorized Hamiltonian, Jordan-Wigner.
This gate assumes that qubits are ordered such that the first norb qubits correspond to the alpha orbitals and the last norb qubits correspond to the beta orbitals.
- __init__(hamiltonian, time, *, n_steps=1, order=0, label=None)[source]¶
Create double-factorized Trotter evolution gate.
- Parameters:
norb – The number of spatial orbitals.
hamiltonian (
DoubleFactorizedHamiltonian
) – The Hamiltonian.time (
float
) – The evolution time.n_steps (
int
) – The number of Trotter steps.order (
int
) – The order of the Trotter decomposition.label (
str
|None
) – The label of the gate.
- class ffsim.qiskit.UCJOpSpinBalancedJW(ucj_op, *, label=None)[source]¶
Bases:
Gate
Spin-balanced UCJ operator under the Jordan-Wigner transformation.
See
ffsim.UCJOpSpinBalanced
for a description of this gate’s unitary.This gate assumes that qubits are ordered such that the first norb qubits correspond to the alpha orbitals and the last norb qubits correspond to the beta orbitals.
- __init__(ucj_op, *, label=None)[source]¶
Create a new spin-balanced unitary cluster Jastrow (UCJ) gate.
- Parameters:
ucj_op (
UCJOpSpinBalanced
) – The UCJ operator.label (
str
|None
) – The label of the gate.
- class ffsim.qiskit.UCJOpSpinUnbalancedJW(ucj_op, *, label=None)[source]¶
Bases:
Gate
Spin-unbalanced UCJ operator under the Jordan-Wigner transformation.
See
ffsim.UCJOpSpinUnbalanced
for a description of this gate’s unitary.This gate assumes that qubits are ordered such that the first norb qubits correspond to the alpha orbitals and the last norb qubits correspond to the beta orbitals.
- __init__(ucj_op, *, label=None)[source]¶
Create a new spin-unbalanced unitary cluster Jastrow (UCJ) gate.
- Parameters:
ucj_op (
UCJOpSpinUnbalanced
) – The UCJ operator.label (
str
|None
) – The label of the gate.
- class ffsim.qiskit.UCJOpSpinlessJW(ucj_op, *, label=None)[source]¶
Bases:
Gate
Spinless UCJ operator under the Jordan-Wigner transformation.
See
ffsim.UCJOpSpinless
for a description of this gate’s unitary.- __init__(ucj_op, *, label=None)[source]¶
Create a new spinless unitary cluster Jastrow (UCJ) gate.
- Parameters:
ucj_op (
UCJOpSpinless
) – The UCJ operator.label (
str
|None
) – The label of the gate.
- ffsim.qiskit.ffsim_vec_to_qiskit_vec(vec, norb, nelec)[source]¶
Convert an ffsim state vector to a Qiskit state vector.
- Parameters:
vec (
ndarray
) – A state vector in ffsim/PySCF format. It should be a one-dimensional vector of lengthcomb(norb, n_alpha) * comb(norb, n_beta)
.norb (
int
) – The number of spatial orbitals.nelec (
int
|tuple
[int
,int
]) – Either a single integer representing the number of fermions for a spinless system, or a pair of integers storing the numbers of spin alpha and spin beta fermions.
- Return type:
ndarray
- ffsim.qiskit.final_state_vector(circuit)[source]¶
Return the final state vector of a fermionic quantum circuit.
- Parameters:
circuit (
QuantumCircuit
) – The circuit composed of fermionic gates.- Return type:
- Returns:
The final state vector that results from applying the circuit to the vacuum state.
- ffsim.qiskit.jordan_wigner(op, n_qubits=None)[source]¶
Jordan-Wigner transformation.
Transform a fermion operator to a qubit operator using the Jordan-Wigner transformation.
- Parameters:
op (
FermionOperator
) – The fermion operator to transform.n_qubits (
int
|None
) – The number of qubits to include in the output qubit operator. If not specified, the minimum number of qubits needed to accommodate the fermion operator will be used. Must be non-negative.
- Return type:
SparsePauliOp
- Returns:
The qubit operator as a Qiskit SparsePauliOp.
- Raises:
ValueError – Number of qubits was negative.
ValueError – Number of qubits was not enough to accommodate the fermion operator.
- ffsim.qiskit.pre_init_passes()[source]¶
Yield transpiler passes recommended for the Qiskit transpiler
pre_init
stage.The following transpiler passes are yielded: :rtype:
Iterator
[BasePass
]Decompose pass that decomposes
PrepareHartreeFockJW
andUCJOperatorJW
gates to expose the underlyingPrepareSlaterDeterminantJW
andOrbitalRotationJW
gates.MergeOrbitalRotations
pass to merge the Slater determinant preparation and orbital rotation gates.
- ffsim.qiskit.qiskit_vec_to_ffsim_vec(vec, norb, nelec)[source]¶
Convert a Qiskit state vector to an ffsim state vector.
- Parameters:
vec (
ndarray
) – A state vector in Qiskit format. It should be a one-dimensional vector of length2 ** (2 * norb)
.norb (
int
) – The number of spatial orbitals.nelec (
int
|tuple
[int
,int
]) – Either a single integer representing the number of fermions for a spinless system, or a pair of integers storing the numbers of spin alpha and spin beta fermions.
- Return type:
ndarray