/
directx.py
661 lines (527 loc) · 28.2 KB
/
directx.py
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
""" Functions for generating Direct-(LGST, MC2GST, MLGST) models """
#***************************************************************************************************
# Copyright 2015, 2019 National Technology & Engineering Solutions of Sandia, LLC (NTESS).
# Under the terms of Contract DE-NA0003525 with NTESS, the U.S. Government retains certain rights
# in this software.
# Licensed under the Apache License, Version 2.0 (the "License"); you may not use this file except
# in compliance with the License. You may obtain a copy of the License at
# http://www.apache.org/licenses/LICENSE-2.0 or in the LICENSE file in the root pyGSTi directory.
#***************************************************************************************************
from .. import tools as _tools
from .. import construction as _construction
from .. import objects as _objs
from . import core as _core
def model_with_lgst_circuit_estimates(
circuits_to_estimate, dataset, prep_strs, effect_strs,
target_model, include_target_ops=True, op_label_aliases=None,
guess_model_for_gauge=None, circuit_labels=None, svd_truncate_to=None,
verbosity=0):
"""
Constructs a model that contains LGST estimates for circuits_to_estimate.
For each operation sequence s in circuits_to_estimate, the constructed model
contains the LGST estimate for s as separate gate, labeled either by
the corresponding element of circuit_labels or by the tuple of s itself.
Parameters
----------
circuits_to_estimate : list of Circuits or tuples
The operation sequences to estimate using LGST
dataset : DataSet
The data to use for LGST
prep_strs,effect_strs : list of Circuits
Fiducial Circuit lists used to construct a informationally complete
preparation and measurement.
target_model : Model
A model used by LGST to specify which operation labels should be estimated,
a guess for which gauge these estimates should be returned in, and
used to simplify operation sequences.
include_target_ops : bool, optional
If True, the operation labels in target_model will be included in the
returned model.
op_label_aliases : dictionary, optional
Dictionary whose keys are operation label "aliases" and whose values are tuples
corresponding to what that operation label should be expanded into before querying
the dataset. Defaults to the empty dictionary (no aliases defined)
e.g. opLabelAliases['Gx^3'] = ('Gx','Gx','Gx')
guess_model_for_gauge : Model, optional
A model used to compute a gauge transformation that is applied to
the LGST estimates. This gauge transformation is computed such that
if the estimated gates matched the model given, then the gate
matrices would match, i.e. the gauge would be the same as
the model supplied. Defaults to the target_model.
circuit_labels : list of strings, optional
A list of labels in one-to-one correspondence with the
operation sequence in circuits_to_estimate. These labels are
the keys to access the operation matrices in the returned
Model, i.e. op_matrix = returned_model[op_label]
svd_truncate_to : int, optional
The Hilbert space dimension to truncate the operation matrices to using
a SVD to keep only the largest svdToTruncateTo singular values of
the I_tildle LGST matrix. Zero means no truncation.
Defaults to dimension of `target_model`.
verbosity : int, optional
Verbosity value to send to do_lgst(...) call.
Returns
-------
Model
A model containing LGST estimates for all the requested
operation sequences and possibly the gates in target_model.
"""
opLabels = [] # list of operation labels for LGST to estimate
if op_label_aliases is None: aliases = {}
else: aliases = op_label_aliases.copy()
#Add operation sequences to estimate as aliases
if circuit_labels is not None:
assert(len(circuit_labels) == len(circuits_to_estimate))
for opLabel, opStr in zip(circuit_labels, circuits_to_estimate):
aliases[opLabel] = opStr.replace_layers_with_aliases(op_label_aliases)
opLabels.append(opLabel)
else:
for opStr in circuits_to_estimate:
newLabel = 'G' + '.'.join(map(str, tuple(opStr)))
aliases[newLabel] = opStr.replace_layers_with_aliases(op_label_aliases) # use circuit tuple as label
opLabels.append(newLabel)
#Add target model labels (not aliased) if requested
if include_target_ops and target_model is not None:
for targetOpLabel in target_model.operations:
if targetOpLabel not in opLabels: # very unlikely that this is false
opLabels.append(targetOpLabel)
return _core.do_lgst(dataset, prep_strs, effect_strs, target_model,
opLabels, aliases, guess_model_for_gauge,
svd_truncate_to, verbosity)
def direct_lgst_model(circuit_to_estimate, circuit_label, dataset,
prep_strs, effect_strs, target_model,
op_label_aliases=None, svd_truncate_to=None, verbosity=0):
"""
Constructs a model of LGST estimates for target gates and circuit_to_estimate.
Parameters
----------
circuit_to_estimate : Circuit or tuple
The single operation sequence to estimate using LGST
circuit_label : string
The label for the estimate of circuit_to_estimate.
i.e. op_matrix = returned_model[op_label]
dataset : DataSet
The data to use for LGST
prep_strs,effect_strs : list of Circuits
Fiducial Circuit lists used to construct a informationally complete
preparation and measurement.
target_model : Model
The target model used by LGST to extract operation labels and an initial gauge
op_label_aliases : dictionary, optional
Dictionary whose keys are operation label "aliases" and whose values are tuples
corresponding to what that operation label should be expanded into before querying
the dataset. Defaults to the empty dictionary (no aliases defined)
e.g. opLabelAliases['Gx^3'] = ('Gx','Gx','Gx')
svd_truncate_to : int, optional
The Hilbert space dimension to truncate the operation matrices to using
a SVD to keep only the largest svdToTruncateTo singular values of
the I_tildle LGST matrix. Zero means no truncation.
Defaults to dimension of `target_model`.
verbosity : int, optional
Verbosity value to send to do_lgst(...) call.
Returns
-------
Model
A model containing LGST estimates of circuit_to_estimate
and the gates of target_model.
"""
return model_with_lgst_circuit_estimates(
[circuit_to_estimate], dataset, prep_strs, effect_strs, target_model,
True, op_label_aliases, None, [circuit_label], svd_truncate_to,
verbosity)
def direct_lgst_models(circuits, dataset, prep_strs, effect_strs, target_model,
op_label_aliases=None, svd_truncate_to=None, verbosity=0):
"""
Constructs a dictionary with keys == operation sequences and values == Direct-LGST Models.
Parameters
----------
circuits : list of Circuit or tuple objects
The operation sequences to estimate using LGST. The elements of this list
are the keys of the returned dictionary.
dataset : DataSet
The data to use for all LGST estimates.
prep_strs,effect_strs : list of Circuits
Fiducial Circuit lists used to construct a informationally complete
preparation and measurement.
target_model : Model
The target model used by LGST to extract operation labels and an initial gauge
op_label_aliases : dictionary, optional
Dictionary whose keys are operation label "aliases" and whose values are tuples
corresponding to what that operation label should be expanded into before querying
the dataset. Defaults to the empty dictionary (no aliases defined)
e.g. opLabelAliases['Gx^3'] = ('Gx','Gx','Gx')
svd_truncate_to : int, optional
The Hilbert space dimension to truncate the operation matrices to using
a SVD to keep only the largest svdToTruncateTo singular values of
the I_tildle LGST matrix. Zero means no truncation.
Defaults to dimension of `target_model`.
verbosity : int, optional
Verbosity value to send to do_lgst(...) call.
Returns
-------
dict
A dictionary that relates each operation sequence of circuits to a
Model containing the LGST estimate of that operation sequence stored under
the operation label "GsigmaLbl", along with LGST estimates of the gates in
target_model.
"""
printer = _objs.VerbosityPrinter.build_printer(verbosity)
directLGSTmodels = {}
printer.log("--- Direct LGST precomputation ---")
with printer.progress_logging(1):
for i, sigma in enumerate(circuits):
printer.show_progress(i, len(circuits), prefix="--- Computing model for string -", suffix='---')
directLGSTmodels[sigma] = direct_lgst_model(
sigma, "GsigmaLbl", dataset, prep_strs, effect_strs, target_model,
op_label_aliases, svd_truncate_to, verbosity)
return directLGSTmodels
def direct_mc2gst_model(circuit_to_estimate, circuit_label, dataset,
prep_strs, effect_strs, target_model,
op_label_aliases=None, svd_truncate_to=None,
min_prob_clip_for_weighting=1e-4,
prob_clip_interval=(-1e6, 1e6), verbosity=0):
"""
Constructs a model of LSGST estimates for target gates and circuit_to_estimate.
Starting with a Direct-LGST estimate for circuit_to_estimate, runs LSGST
using the same strings that LGST would have used to estimate circuit_to_estimate
and each of the target gates. That is, LSGST is run with strings of the form:
1. prepStr
2. effectStr
3. prepStr + effectStr
4. prepStr + singleGate + effectStr
5. prepStr + circuit_to_estimate + effectStr
and the resulting Model estimate is returned.
Parameters
----------
circuit_to_estimate : Circuit
The single operation sequence to estimate using LSGST
circuit_label : string
The label for the estimate of circuit_to_estimate.
i.e. op_matrix = returned_mode[op_label]
dataset : DataSet
The data to use for LGST
prep_strs,effect_strs : list of Circuits
Fiducial Circuit lists used to construct a informationally complete
preparation and measurement.
target_model : Model
The target model used by LGST to extract operation labels and an initial gauge
op_label_aliases : dictionary, optional
Dictionary whose keys are operation label "aliases" and whose values are tuples
corresponding to what that operation label should be expanded into before querying
the dataset. Defaults to the empty dictionary (no aliases defined)
e.g. opLabelAliases['Gx^3'] = ('Gx','Gx','Gx')
svd_truncate_to : int, optional
The Hilbert space dimension to truncate the operation matrices to using
a SVD to keep only the largest svdToTruncateTo singular values of
the I_tildle LGST matrix. Zero means no truncation.
Defaults to dimension of `target_model`.
min_prob_clip_for_weighting : float, optional
defines the clipping interval for the statistical weight used
within the chi^2 function (see chi2fn).
prob_clip_interval : 2-tuple, optional
(min,max) to clip probabilities to within Model probability
computation routines (see Model.bulk_fill_probs)
verbosity : int, optional
Verbosity value to send to do_lgst(...) and do_mc2gst(...) calls.
Returns
-------
Model
A model containing LSGST estimates of circuit_to_estimate
and the gates of target_model.
"""
direct_lgst = model_with_lgst_circuit_estimates(
[circuit_to_estimate], dataset, prep_strs, effect_strs, target_model,
True, op_label_aliases, None, [circuit_label], svd_truncate_to, verbosity)
# LEXICOGRAPHICAL VS MATRIX ORDER
circuits = prep_strs + effect_strs + [prepStr + effectStr for prepStr in prep_strs for effectStr in effect_strs]
for opLabel in direct_lgst.operations:
circuits.extend([prepStr + _objs.Circuit((opLabel,)) + effectStr
for prepStr in prep_strs for effectStr in effect_strs])
aliases = {} if (op_label_aliases is None) else op_label_aliases.copy()
aliases[circuit_label] = circuit_to_estimate.replace_layers_with_aliases(op_label_aliases)
obuilder = _objs.Chi2Function.builder(regularization={'min_prob_clip_for_weighting': min_prob_clip_for_weighting},
penalties={'prob_clip_interval': prob_clip_interval})
bulk_circuits = _objs.BulkCircuitList(circuits, aliases)
_, direct_lsgst = _core.do_gst_fit(dataset, direct_lgst, bulk_circuits, optimizer=None,
objective_function_builder=obuilder, resource_alloc=None, cache=None,
verbosity=verbosity)
return direct_lsgst
def direct_mc2gst_models(circuits, dataset, prep_strs, effect_strs,
target_model, op_label_aliases=None,
svd_truncate_to=None, min_prob_clip_for_weighting=1e-4,
prob_clip_interval=(-1e6, 1e6), verbosity=0):
"""
Constructs a dictionary with keys == operation sequences and values == Direct-LSGST Models.
Parameters
----------
circuits : list of Circuit or tuple objects
The operation sequences to estimate using LSGST. The elements of this list
are the keys of the returned dictionary.
dataset : DataSet
The data to use for all LGST and LSGST estimates.
prep_strs,effect_strs : list of Circuits
Fiducial Circuit lists used to construct a informationally complete
preparation and measurement.
target_model : Model
The target model used by LGST to extract operation labels and an initial gauge
op_label_aliases : dictionary, optional
Dictionary whose keys are operation label "aliases" and whose values are tuples
corresponding to what that operation label should be expanded into before querying
the dataset. Defaults to the empty dictionary (no aliases defined)
e.g. opLabelAliases['Gx^3'] = ('Gx','Gx','Gx')
svd_truncate_to : int, optional
The Hilbert space dimension to truncate the operation matrices to using
a SVD to keep only the largest svdToTruncateTo singular values of
the I_tildle LGST matrix. Zero means no truncation.
Defaults to dimension of `target_model`.
min_prob_clip_for_weighting : float, optional
defines the clipping interval for the statistical weight used
within the chi^2 function (see chi2fn).
prob_clip_interval : 2-tuple, optional
(min,max) to clip probabilities to within Model probability
computation routines (see Model.bulk_fill_probs)
verbosity : int, optional
Verbosity value to send to do_lgst(...) and do_mc2gst(...) calls.
Returns
-------
dict
A dictionary that relates each operation sequence of circuits to a
Model containing the LSGST estimate of that operation sequence stored under
the operation label "GsigmaLbl", along with LSGST estimates of the gates in
target_model.
"""
printer = _objs.VerbosityPrinter.build_printer(verbosity)
directLSGSTmodels = {}
printer.log("--- Direct LSGST precomputation ---")
with printer.progress_logging(1):
for i, sigma in enumerate(circuits):
printer.show_progress(i, len(circuits), prefix="--- Computing model for string-", suffix='---')
directLSGSTmodels[sigma] = direct_mc2gst_model(
sigma, "GsigmaLbl", dataset, prep_strs, effect_strs, target_model,
op_label_aliases, svd_truncate_to, min_prob_clip_for_weighting,
prob_clip_interval, verbosity)
return directLSGSTmodels
def direct_mlgst_model(circuit_to_estimate, circuit_label, dataset,
prep_strs, effect_strs, target_model,
op_label_aliases=None, svd_truncate_to=None, min_prob_clip=1e-6,
prob_clip_interval=(-1e6, 1e6), verbosity=0):
"""
Constructs a model of MLEGST estimates for target gates and circuit_to_estimate.
Starting with a Direct-LGST estimate for circuit_to_estimate, runs MLEGST
using the same strings that LGST would have used to estimate circuit_to_estimate
and each of the target gates. That is, MLEGST is run with strings of the form:
1. prepStr
2. effectStr
3. prepStr + effectStr
4. prepStr + singleGate + effectStr
5. prepStr + circuit_to_estimate + effectStr
and the resulting Model estimate is returned.
Parameters
----------
circuit_to_estimate : Circuit or tuple
The single operation sequence to estimate using LSGST
circuit_label : string
The label for the estimate of circuit_to_estimate.
i.e. op_matrix = returned_model[op_label]
dataset : DataSet
The data to use for LGST
prep_strs,effect_strs : list of Circuits
Fiducial Circuit lists used to construct a informationally complete
preparation and measurement.
target_model : Model
The target model used by LGST to extract operation labels and an initial gauge
op_label_aliases : dictionary, optional
Dictionary whose keys are operation label "aliases" and whose values are tuples
corresponding to what that operation label should be expanded into before querying
the dataset. Defaults to the empty dictionary (no aliases defined)
e.g. opLabelAliases['Gx^3'] = ('Gx','Gx','Gx')
svd_truncate_to : int, optional
The Hilbert space dimension to truncate the operation matrices to using
a SVD to keep only the largest svdToTruncateTo singular values of
the I_tildle LGST matrix. Zero means no truncation.
Defaults to dimension of `target_model`.
min_prob_clip : float, optional
defines the minimum probability "patch point" used
within the logl function.
prob_clip_interval : 2-tuple, optional
(min,max) to clip probabilities to within Model probability
computation routines (see Model.bulk_fill_probs)
verbosity : int, optional
Verbosity value to send to do_lgst(...) and do_mlgst(...) calls.
Returns
-------
Model
A model containing MLEGST estimates of circuit_to_estimate
and the gates of target_model.
"""
direct_lgst = model_with_lgst_circuit_estimates(
[circuit_to_estimate], dataset, prep_strs, effect_strs, target_model,
True, op_label_aliases, None, [circuit_label], svd_truncate_to, verbosity)
# LEXICOGRAPHICAL VS MATRIX ORDER
circuits = prep_strs + effect_strs + [prepStr + effectStr for prepStr in prep_strs for effectStr in effect_strs]
for opLabel in direct_lgst.operations:
circuits.extend([prepStr + _objs.Circuit((opLabel,)) + effectStr
for prepStr in prep_strs for effectStr in effect_strs])
aliases = {} if (op_label_aliases is None) else op_label_aliases.copy()
aliases[circuit_label] = circuit_to_estimate.replace_layers_with_aliases(op_label_aliases)
obuilder = _objs.PoissonPicDeltaLogLFunction.builder(regularization={'min_prob_clip': min_prob_clip},
penalties={'prob_clip_interval': prob_clip_interval})
bulk_circuits = _objs.BulkCircuitList(circuits, aliases)
_, direct_mlegst = _core.do_gst_fit(dataset, direct_lgst, bulk_circuits, optimizer=None,
objective_function_builder=obuilder, resource_alloc=None, cache=None,
verbosity=verbosity)
return direct_mlegst
def direct_mlgst_models(circuits, dataset, prep_strs, effect_strs, target_model,
op_label_aliases=None, svd_truncate_to=None, min_prob_clip=1e-6,
prob_clip_interval=(-1e6, 1e6), verbosity=0):
"""
Constructs a dictionary with keys == operation sequences and values == Direct-MLEGST Models.
Parameters
----------
circuits : list of Circuit or tuple objects
The operation sequences to estimate using MLEGST. The elements of this list
are the keys of the returned dictionary.
dataset : DataSet
The data to use for all LGST and LSGST estimates.
prep_strs,effect_strs : list of Circuits
Fiducial Circuit lists used to construct a informationally complete
preparation and measurement.
target_model : Model
The target model used by LGST to extract operation labels and an initial gauge
op_label_aliases : dictionary, optional
Dictionary whose keys are operation label "aliases" and whose values are tuples
corresponding to what that operation label should be expanded into before querying
the dataset. Defaults to the empty dictionary (no aliases defined)
e.g. opLabelAliases['Gx^3'] = ('Gx','Gx','Gx')
svd_truncate_to : int, optional
The Hilbert space dimension to truncate the operation matrices to using
a SVD to keep only the largest svdToTruncateTo singular values of
the I_tildle LGST matrix. Zero means no truncation.
Defaults to dimension of `target_model`.
min_prob_clip : float, optional
defines the minimum probability "patch point" used
within the logl function.
prob_clip_interval : 2-tuple, optional
(min,max) to clip probabilities to within Model probability
computation routines (see Model.bulk_fill_probs)
verbosity : int, optional
Verbosity value to send to do_lgst(...) and do_mlgst(...) calls.
Returns
-------
dict
A dictionary that relates each operation sequence of circuits to a
Model containing the MLEGST estimate of that operation sequence stored under
the operation label "GsigmaLbl", along with MLEGST estimates of the gates in
target_model.
"""
printer = _objs.VerbosityPrinter.build_printer(verbosity)
directMLEGSTmodels = {}
printer.log("--- Direct MLEGST precomputation ---")
with printer.progress_logging(1):
for i, sigma in enumerate(circuits):
printer.show_progress(i, len(circuits), prefix="--- Computing model for string ", suffix="---")
directMLEGSTmodels[sigma] = direct_mlgst_model(
sigma, "GsigmaLbl", dataset, prep_strs, effect_strs, target_model,
op_label_aliases, svd_truncate_to, min_prob_clip,
prob_clip_interval, verbosity)
return directMLEGSTmodels
def focused_mc2gst_model(circuit_to_estimate, circuit_label, dataset,
prep_strs, effect_strs, start_model,
op_label_aliases=None, min_prob_clip_for_weighting=1e-4,
prob_clip_interval=(-1e6, 1e6), verbosity=0):
"""
Constructs a model containing a single LSGST estimate of circuit_to_estimate.
Starting with start_model, run LSGST with the same operation sequences that LGST
would use to estimate circuit_to_estimate. That is, LSGST is run with
strings of the form: prepStr + circuit_to_estimate + effectStr
and return the resulting Model.
Parameters
----------
circuit_to_estimate : Circuit or tuple
The single operation sequence to estimate using LSGST
circuit_label : string
The label for the estimate of circuit_to_estimate.
i.e. op_matrix = returned_model[op_label]
dataset : DataSet
The data to use for LGST
prep_strs,effect_strs : list of Circuits
Fiducial Circuit lists used to construct a informationally complete
preparation and measurement.
start_model : Model
The model to seed LSGST with. Often times obtained via LGST.
op_label_aliases : dictionary, optional
Dictionary whose keys are operation label "aliases" and whose values are tuples
corresponding to what that operation label should be expanded into before querying
the dataset. Defaults to the empty dictionary (no aliases defined)
e.g. opLabelAliases['Gx^3'] = ('Gx','Gx','Gx')
min_prob_clip_for_weighting : float, optional
defines the clipping interval for the statistical weight used
within the chi^2 function (see chi2fn).
prob_clip_interval : 2-tuple, optional
(min,max) to clip probabilities to within Model probability
computation routines (see Model.bulk_fill_probs)
verbosity : int, optional
Verbosity value to send do_mc2gst(...) call.
Returns
-------
Model
A model containing LSGST estimate of circuit_to_estimate.
"""
circuits = [prepStr + circuit_to_estimate + effectStr for prepStr in prep_strs for effectStr in effect_strs]
obuilder = _objs.Chi2Function.builder(regularization={'min_prob_clip_for_weighting': min_prob_clip_for_weighting},
penalties={'prob_clip_interval': prob_clip_interval})
bulk_circuits = _objs.BulkCircuitList(circuits, op_label_aliases)
_, focused_lsgst = _core.do_gst_fit(dataset, start_model, bulk_circuits, optimizer=None,
objective_function_builder=obuilder, resource_alloc=None, cache=None,
verbosity=verbosity)
focused_lsgst.operations[circuit_label] = _objs.FullDenseOp(
focused_lsgst.product(circuit_to_estimate)) # add desired string as a separate labeled gate
return focused_lsgst
def focused_mc2gst_models(circuits, dataset, prep_strs, effect_strs,
start_model, op_label_aliases=None,
min_prob_clip_for_weighting=1e-4,
prob_clip_interval=(-1e6, 1e6), verbosity=0):
"""
Constructs a dictionary with keys == operation sequences and values == Focused-LSGST Models.
Parameters
----------
circuits : list of Circuit or tuple objects
The operation sequences to estimate using LSGST. The elements of this list
are the keys of the returned dictionary.
dataset : DataSet
The data to use for all LGST and LSGST estimates.
prep_strs,effect_strs : list of Circuits
Fiducial Circuit lists used to construct a informationally complete
preparation and measurement.
start_model : Model
The model to seed LSGST with. Often times obtained via LGST.
op_label_aliases : dictionary, optional
Dictionary whose keys are operation label "aliases" and whose values are tuples
corresponding to what that operation label should be expanded into before querying
the dataset. Defaults to the empty dictionary (no aliases defined)
e.g. opLabelAliases['Gx^3'] = ('Gx','Gx','Gx')
min_prob_clip_for_weighting : float, optional
defines the clipping interval for the statistical weight used
within the chi^2 function (see chi2fn).
prob_clip_interval : 2-tuple, optional
(min,max) to clip probabilities to within Model probability
computation routines (see Model.bulk_fill_probs)
verbosity : int, optional
Verbosity value to send to do_mc2gst(...) call.
Returns
-------
dict
A dictionary that relates each operation sequence of circuits to a
Model containing the LSGST estimate of that operation sequence stored under
the operation label "GsigmaLbl".
"""
printer = _objs.VerbosityPrinter.build_printer(verbosity)
focusedLSGSTmodels = {}
printer.log("--- Focused LSGST precomputation ---")
with printer.progress_logging(1):
for i, sigma in enumerate(circuits):
printer.show_progress(i, len(circuits), prefix="--- Computing model for string", suffix='---')
focusedLSGSTmodels[sigma] = focused_mc2gst_model(
sigma, "GsigmaLbl", dataset, prep_strs, effect_strs, start_model,
op_label_aliases, min_prob_clip_for_weighting, prob_clip_interval, verbosity)
return focusedLSGSTmodels