@@ -845,6 +845,76 @@ def frame_potential(
845
845
return potential / samples ** 2
846
846
847
847
848
+ def quantum_fisher_information_matrix (
849
+ circuit ,
850
+ parameters = None ,
851
+ initial_state = None ,
852
+ return_complex : bool = True ,
853
+ backend = None ,
854
+ ):
855
+ """Calculate the Quantum Fisher Information Matrix (QFIM) of a parametrized ``circuit``.
856
+
857
+ Given a set of ``parameters`` :math:`\\ theta = \\ {\\ theta_{k}\\ }_{k\\ in[M]}` and a
858
+ parameterized unitary ``circuit`` :math:`U(\\ theta)` acting on an ``initial_state``
859
+ :math:`\\ ket{\\ phi}`, the QFIM is such that its elements can be calculated as
860
+
861
+ .. math::
862
+ \\ mathbf{F}_{jk} = 4 \\ , \\ text{Re}\\ left\\ { \\ braket{\\ partial_{j} \\ psi | \\ partial_{k}
863
+ \\ psi} - \\ braket{\\ partial_{j} \\ psi | \\ psi}\\ !\\ braket{\\ psi | \\ partial_{k} \\ psi}
864
+ \\ right\\ } \\ , ,
865
+
866
+ where we have used the short notations :math:`\\ ket{\\ psi} \\ equiv \\ ket{\\ psi(\\ theta)}
867
+ = U(\\ theta) \\ ket{\\ phi}`, and :math:`\\ ket{\\ partial_{k} \\ psi} \\ equiv \\ frac{\\ partial}
868
+ {\\ partial\\ theta_{k}} \\ ket{\\ psi(\\ theta)}`.
869
+ If the ``initial_state`` :math:`\\ ket{\\ phi}` is not specified, it defaults to
870
+ :math:`\\ ket{0}^{\\ otimes n}`.
871
+
872
+ Args:
873
+ circuit (:class:`qibo.models.circuit.Circuit`): parametrized circuit :math:`U(\\ theta)`.
874
+ parameters (ndarray, optional): parameters whose QFIM to calculate.
875
+ If ``None``, QFIM is calculated with the paremeters from ``circuit``, i.e.
876
+ ``parameters = circuit.get_parameters()``. Defaults to ``None``.
877
+ initial_state (ndarray, optional): Initial configuration. It can be specified
878
+ by the setting the state vector using an array or a circuit. If ``None``,
879
+ the initial state is :math:`\\ ket{0}^{\\ otimes n}`. Defaults to ``None``.
880
+ return_complex (bool, optional): If ``True``, calculates the Jacobian matrix
881
+ of real and imaginary parts of :math:`\\ ket{\\ psi(\\ theta)}`. If ``False``,
882
+ calculates only the Jacobian matrix of the real part. Defaults to ``True``.
883
+ backend (:class:`qibo.backends.abstract.Backend`, optional): backend to be used
884
+ in the execution. If ``None``, it uses :class:`qibo.backends.GlobalBackend`.
885
+ Defaults to ``None``.
886
+
887
+ Returns:
888
+ ndarray: Quantum Fisher Information :math:`\\ mathbf{F}`.
889
+ """
890
+ backend = _check_backend (backend )
891
+
892
+ if parameters is None :
893
+ parameters = circuit .get_parameters ()
894
+ parameters = backend .cast (parameters , dtype = float ).flatten ()
895
+
896
+ jacobian = backend .calculate_jacobian_matrix (
897
+ circuit , parameters , initial_state , return_complex
898
+ )
899
+
900
+ if return_complex :
901
+ jacobian = jacobian [0 ] + 1j * jacobian [1 ]
902
+
903
+ jacobian = backend .cast (jacobian , dtype = np .complex128 )
904
+
905
+ copied = circuit .copy (deep = True )
906
+ copied .set_parameters (parameters )
907
+
908
+ state = backend .execute_circuit (copied , initial_state = initial_state ).state ()
909
+
910
+ overlaps = jacobian .T @ state
911
+
912
+ qfim = jacobian .T @ jacobian
913
+ qfim = qfim - backend .np .outer (overlaps , backend .np .conj (overlaps .T ))
914
+
915
+ return 4 * backend .np .real (qfim )
916
+
917
+
848
918
def _check_hermitian (matrix , backend = None ):
849
919
"""Checks if a given matrix is Hermitian.
850
920
0 commit comments