11
11
#***************************************************************************************************
12
12
import functools as _functools
13
13
import itertools as _itertools
14
-
14
+ import warnings
15
15
import numpy as _np
16
+ import scipy .linalg as _spl
17
+ import scipy .optimize as _spo
16
18
17
19
from .complementeffect import ComplementPOVMEffect
18
20
from .composedeffect import ComposedPOVMEffect
34
36
from pygsti .baseobjs import statespace as _statespace
35
37
from pygsti .tools import basistools as _bt
36
38
from pygsti .tools import optools as _ot
39
+ from pygsti .tools import sum_of_negative_choi_eigenvalues_gate
40
+ from pygsti .baseobjs import Basis
37
41
38
42
# Avoid circular import
39
43
import pygsti .modelmembers as _mm
@@ -418,7 +422,7 @@ def povm_type_from_op_type(op_type):
418
422
return povm_type_preferences
419
423
420
424
421
- def convert (povm , to_type , basis , ideal_povm = None , flatten_structure = False ):
425
+ def convert (povm , to_type , basis , ideal_povm = None , flatten_structure = False , cp_penalty = 1e-7 ):
422
426
"""
423
427
TODO: update docstring
424
428
Convert a POVM to a new type of parameterization.
@@ -450,13 +454,22 @@ def convert(povm, to_type, basis, ideal_povm=None, flatten_structure=False):
450
454
are separately converted, leaving the original POVM's structure
451
455
unchanged. When `True`, composed and embedded operations are "flattened"
452
456
into a single POVM of the requested `to_type`.
457
+
458
+ cp_penalty : float, optional (default 1e-7)
459
+ Converting SPAM operations to an error generator representation may
460
+ introduce trivial gauge degrees of freedom. These gauge degrees of freedom
461
+ are called trivial because they quite literally do not change the dense representation
462
+ (i.e. Hilbert-Schmidt vectors) at all. Despite being trivial, error generators along
463
+ this trivial gauge orbit may be non-CP, so this cptp penalty is used to favor channels
464
+ within this gauge orbit which are CPTP.
453
465
454
466
Returns
455
467
-------
456
468
POVM
457
469
The converted POVM vector, usually a distinct
458
470
object from the object passed as input.
459
471
"""
472
+
460
473
to_types = to_type if isinstance (to_type , (tuple , list )) else (to_type ,) # HACK to support multiple to_type values
461
474
error_msgs = {}
462
475
@@ -487,6 +500,7 @@ def convert(povm, to_type, basis, ideal_povm=None, flatten_structure=False):
487
500
from ..operations import LindbladParameterization as _LindbladParameterization
488
501
lndtype = _LindbladParameterization .cast (to_type )
489
502
503
+
490
504
#Construct a static "base" POVM
491
505
if isinstance (povm , ComputationalBasisPOVM ): # special easy case
492
506
base_povm = ComputationalBasisPOVM (povm .state_space .num_qubits , povm .evotype ) # just copy it?
@@ -498,14 +512,105 @@ def convert(povm, to_type, basis, ideal_povm=None, flatten_structure=False):
498
512
for lbl , vec in povm .items ()]
499
513
else :
500
514
raise RuntimeError ('Evotype must be compatible with Hilbert ops to use pure effects' )
501
- except Exception : # try static mixed states next:
502
- base_items = [(lbl , convert_effect (vec , 'static' , basis , idl .get (lbl , None ), flatten_structure ))
503
- for lbl , vec in povm .items ()]
515
+ except RuntimeError : # try static mixed states next:
516
+ #if idl.get(lbl,None) is not None:
517
+
518
+ base_items = []
519
+ for lbl , vec in povm .items ():
520
+ ideal_effect = idl .get (lbl ,None )
521
+ if ideal_effect is not None :
522
+ base_items .append ((lbl , convert_effect (ideal_effect , 'static' , basis , ideal_effect , flatten_structure )))
523
+ else :
524
+ base_items .append ((lbl , convert_effect (vec , 'static' , basis , idl .get (lbl , None ), flatten_structure )))
504
525
base_povm = UnconstrainedPOVM (base_items , povm .evotype , povm .state_space )
505
526
506
527
proj_basis = 'PP' if povm .state_space .is_entirely_qubits else basis
507
528
errorgen = _LindbladErrorgen .from_error_generator (povm .state_space .dim , lndtype , proj_basis ,
508
529
basis , truncate = True , evotype = povm .evotype )
530
+
531
+ #Collect all ideal effects
532
+ base_dense_effects = []
533
+ for item in base_items :
534
+ dense_effect = item [1 ].to_dense ()
535
+ base_dense_effects .append (dense_effect .reshape ((1 ,len (dense_effect ))))
536
+
537
+ dense_ideal_povm = _np .concatenate (base_dense_effects , axis = 0 )
538
+
539
+ #Collect all noisy effects
540
+ dense_effects = []
541
+ for effect in povm .values ():
542
+ dense_effect = effect .to_dense ()
543
+ dense_effects .append (dense_effect .reshape ((1 ,len (dense_effect ))))
544
+
545
+ dense_povm = _np .concatenate (dense_effects , axis = 0 )
546
+
547
+ #It is often the case that there are more error generators than physical degrees of freedom in the POVM
548
+ #We define a function which finds linear comb. of errgens that span these degrees of freedom.
549
+ #This has been called "the trivial gauge", and this function is meant to avoid it
550
+ def calc_physical_subspace (dense_ideal_povm , epsilon = 1e-4 ):
551
+
552
+ degrees_of_freedom = (dense_ideal_povm .shape [0 ] - 1 ) * dense_ideal_povm .shape [1 ]
553
+ errgen = _LindbladErrorgen .from_error_generator (povm .state_space .dim , parameterization = to_type )
554
+
555
+ if degrees_of_freedom > errgen .num_params :
556
+ warnings .warn ("POVM has more degrees of freedom than the available number of parameters, representation in this parameterization is not guaranteed" )
557
+ exp_errgen = _ExpErrorgenOp (errgen )
558
+
559
+ num_errgens = errgen .num_params
560
+ #TODO: Maybe we can use the num of params instead of number of matrix entries, as some of them are linearly dependent.
561
+ #i.e E0 completely determines E1 if those are the only two povm elements (E0 + E1 = Identity)
562
+ num_entries = dense_ideal_povm .size
563
+
564
+ #Compute the jacobian with respect to the error generators. This will allow us to see which
565
+ #error generators change the POVM entries
566
+ J = _np .zeros ((num_entries ,num_errgens ))
567
+ new_vec = _np .zeros (num_errgens )
568
+ for i in range (num_errgens ):
569
+
570
+ new_vec [i ] = epsilon
571
+ exp_errgen .from_vector (new_vec )
572
+ new_vec [i ] = 0
573
+ vectorized_povm = _np .zeros (num_entries )
574
+ perturbed_povm = (dense_ideal_povm @ exp_errgen .to_dense () - dense_ideal_povm )/ epsilon
575
+
576
+ vectorized_povm = perturbed_povm .flatten (order = 'F' )
577
+
578
+ J [:,i ] = vectorized_povm
579
+
580
+ _ ,S ,Vt = _np .linalg .svd (J , full_matrices = False )
581
+
582
+ #Only return nontrivial singular vectors
583
+ Vt = Vt [S > 1e-13 , :].reshape ((- 1 , Vt .shape [1 ]))
584
+ return Vt
585
+
586
+
587
+ phys_directions = calc_physical_subspace (dense_ideal_povm )
588
+
589
+ #We use optimization to find the best error generator representation
590
+ #we only vary physical directions, not independent error generators
591
+ def _objfn (v ):
592
+
593
+ #For some reason adding the sum_of_negative_choi_eigenvalues_gate term
594
+ #resulted in minimize() sometimes choosing NaN values for v. There are
595
+ #two stack exchange issues showing this problem with no solution.
596
+ if _np .isnan (v ).any ():
597
+ v = _np .zeros (len (v ))
598
+
599
+ L_vec = _np .zeros (len (phys_directions [0 ]))
600
+ for coeff , phys_direction in zip (v ,phys_directions ):
601
+ L_vec += coeff * phys_direction
602
+ errorgen .from_vector (L_vec )
603
+ proc_matrix = _spl .expm (errorgen .to_dense ())
604
+
605
+ return _np .linalg .norm (dense_povm - dense_ideal_povm @ proc_matrix ) + cp_penalty * sum_of_negative_choi_eigenvalues_gate (proc_matrix , basis )
606
+
607
+ soln = _spo .minimize (_objfn , _np .zeros (len (phys_directions ), 'd' ), method = "Nelder-Mead" , options = {},
608
+ tol = 1e-13 )
609
+ if not soln .success and soln .fun > 1e-6 : # not "or" because success is often not set correctly
610
+ raise ValueError ("Failed to find an errorgen such that <ideal|exp(errorgen) = <effect|" )
611
+ errgen_vec = _np .linalg .lstsq (phys_directions , soln .x )[0 ]
612
+ errorgen .from_vector (errgen_vec )
613
+
509
614
EffectiveExpErrorgen = _IdentityPlusErrorgenOp if lndtype .meta == '1+' else _ExpErrorgenOp
510
615
return ComposedPOVM (EffectiveExpErrorgen (errorgen ), base_povm , mx_basis = basis )
511
616
0 commit comments