Qrack  9.13
General classical-emulating-quantum development framework
wasm_api.hpp
Go to the documentation of this file.
1 // Copyright (c) Microsoft Corporation. All rights reserved.
2 // Licensed under the MIT License.
3 //
4 // (Extensively modified and adapted by Daniel Strano in unitaryfund/qrack)
5 
6 #pragma once
7 
8 #include "common/pauli.hpp"
10 #include "common/qrack_types.hpp"
11 
12 #include <set>
13 #include <string>
14 #include <vector>
15 
25 namespace Qrack {
26 
27 typedef uint64_t quid;
28 
31  bool val;
33  : qid(q)
34  , val(v)
35  {
36  // Intentionally left blank
37  }
38 };
39 
44  : qid(q)
45  , val(v)
46  {
47  // Intentionally left blank
48  }
49 };
50 
55  : qid(q)
56  , val(v)
57  {
58  // Intentionally left blank
59  }
60 };
61 
66  : qid(q)
67  , b(basis)
68  {
69  // Intentionally left blank
70  }
71 };
72 
73 struct QubitU3Basis {
75  real1_f b[3U];
76  QubitU3Basis(bitLenInt q, std::vector<real1_f> basis)
77  : qid(q)
78  {
79  if (basis.size() != 3U) {
80  throw std::invalid_argument("QubitU3Basis must have 3 basis angles!");
81  }
82  for (int i = 0; i < 3; ++i) {
83  b[i] = basis[i];
84  }
85  }
86 };
87 
90  complex b[4U];
91  QubitMatrixBasis(bitLenInt q, std::vector<complex> basis)
92  : qid(q)
93  {
94  if (basis.size() != 4U) {
95  throw std::invalid_argument("QubitMatrixBasis must have 4 matrix components for basis!");
96  }
97  for (int i = 0; i < 4; ++i) {
98  b[i] = basis[i];
99  }
100  }
101 };
102 
105  real1_f b[3U];
106  real1_f e[2U];
107  QubitU3BasisEigenVal(bitLenInt q, std::vector<real1_f> basis, std::vector<real1_f> ex)
108  : qid(q)
109  {
110  if (basis.size() != 3U) {
111  throw std::invalid_argument("QubitU3BasisEigenVal must have 3 basis angles!");
112  }
113  if (ex.size() != 2U) {
114  throw std::invalid_argument("QubitU3BasisEigenVal must have 2 eigenvalues!");
115  }
116  for (int i = 0; i < 3; ++i) {
117  b[i] = basis[i];
118  }
119  e[0U] = ex[0U];
120  e[1U] = ex[1U];
121  }
122 };
123 
126  complex b[4U];
127  real1_f e[2U];
128  QubitMatrixBasisEigenVal(bitLenInt q, std::vector<complex> basis, std::vector<real1_f> ex)
129  : qid(q)
130  {
131  if (basis.size() != 4U) {
132  throw std::invalid_argument("QubitMatrixBasisEigenVal must have 4 matrix components for basis!");
133  }
134  if (ex.size() != 2U) {
135  throw std::invalid_argument("QubitMatrixBasisEigenVal must have 2 eigenvalues!");
136  }
137  for (int i = 0; i < 4; ++i) {
138  b[i] = basis[i];
139  }
140  e[0U] = ex[0U];
141  e[1U] = ex[1U];
142  }
143 };
144 
159  bitLenInt q, bool tn, bool md, bool sd, bool sh, bool bdt, bool pg, bool nw, bool hy, bool oc, bool hp);
160 
161 // Utility
162 
166 quid init_count(bitLenInt q, bool dm);
167 
171 quid init();
172 
176 quid init_clone(quid sid);
177 
182 
186 void destroy(quid sid);
187 
191 void seed(quid sid, unsigned s);
192 
196 void set_concurrency(quid sid, unsigned p);
197 
201 void allocateQubit(quid sid, bitLenInt qid);
205 bool release(quid sid, bitLenInt q);
213 void SetPermutation(quid sid, bitCapInt p);
214 
218 void qstabilizer_out_to_file(quid sid, std::string f);
223 void qstabilizer_in_from_file(quid sid, std::string f);
224 
228 std::vector<real1> ProbAll(quid sid, std::vector<bitLenInt> q);
232 real1_f Prob(quid sid, bitLenInt q);
236 real1_f ProbRdm(quid sid, bitLenInt q);
240 real1_f PermutationProb(quid sid, std::vector<QubitIndexState> q);
244 real1_f PermutationProbRdm(quid sid, std::vector<QubitIndexState> q, bool r);
248 real1_f PermutationExpectation(quid sid, std::vector<bitLenInt> q);
252 real1_f PermutationExpectationRdm(quid sid, std::vector<bitLenInt> q, bool r);
256 real1_f FactorizedExpectation(quid sid, std::vector<QubitIntegerExpVar> q);
261 real1_f FactorizedExpectationRdm(quid sid, std::vector<QubitIntegerExpVar> q, bool r);
265 real1_f FactorizedExpectationFp(quid sid, std::vector<QubitRealExpVar> q);
270 real1_f FactorizedExpectationFpRdm(quid sid, std::vector<QubitRealExpVar> q, bool r);
274 real1_f UnitaryExpectation(quid sid, std::vector<QubitU3Basis> q);
278 real1_f MatrixExpectation(quid sid, std::vector<QubitMatrixBasis> q);
282 real1_f UnitaryExpectationEigenVal(quid sid, std::vector<QubitU3BasisEigenVal> q);
286 real1_f MatrixExpectationEigenVal(quid sid, std::vector<QubitMatrixBasisEigenVal> q);
290 real1_f PauliExpectation(quid sid, std::vector<QubitPauliBasis> q);
294 real1_f Variance(quid sid, std::vector<bitLenInt> q);
298 real1_f VarianceRdm(quid sid, std::vector<bitLenInt> q, bool r);
302 real1_f FactorizedVariance(quid sid, std::vector<QubitIntegerExpVar> q);
307 real1_f FactorizedVarianceRdm(quid sid, std::vector<QubitIntegerExpVar> q, bool r);
311 real1_f FactorizedVarianceFp(quid sid, std::vector<QubitRealExpVar> q);
316 real1_f FactorizedVarianceFpRdm(quid sid, std::vector<QubitRealExpVar> q, bool r);
320 real1_f UnitaryVariance(quid sid, std::vector<QubitU3Basis> q);
324 real1_f MatrixVariance(quid sid, std::vector<QubitMatrixBasis> q);
328 real1_f UnitaryVarianceEigenVal(quid sid, std::vector<QubitU3BasisEigenVal> q);
332 real1_f MatrixVarianceEigenVal(quid sid, std::vector<QubitMatrixBasisEigenVal> q);
336 real1_f PauliVariance(quid sid, std::vector<QubitPauliBasis> q);
337 
341 size_t random_choice(quid sid, std::vector<real1> p);
342 
346 void PhaseParity(quid sid, real1_f lambda, std::vector<bitLenInt> q);
347 
351 void PhaseRootN(quid sid, bitLenInt p, std::vector<bitLenInt> q);
352 
356 real1_f JointEnsembleProbability(quid sid, std::vector<QubitPauliBasis> q);
357 
358 // SPAM and non-unitary
359 
363 bool M(quid sid, bitLenInt q);
367 bool ForceM(quid sid, bitLenInt q, bool r);
371 bool Measure(quid sid, std::vector<QubitPauliBasis> q);
375 bitCapInt MAll(quid sid);
379 std::vector<long long unsigned int> MeasureShots(quid sid, std::vector<bitLenInt> q, unsigned s);
383 void ResetAll(quid sid);
384 
385 // single-qubit gates
386 void X(quid sid, bitLenInt q);
387 void Y(quid sid, bitLenInt q);
388 void Z(quid sid, bitLenInt q);
389 void H(quid sid, bitLenInt q);
390 void S(quid sid, bitLenInt q);
391 void SX(quid sid, bitLenInt q);
392 void SY(quid sid, bitLenInt q);
393 void T(quid sid, bitLenInt q);
394 void AdjS(quid sid, bitLenInt q);
395 void AdjSX(quid sid, bitLenInt q);
396 void AdjSY(quid sid, bitLenInt q);
397 void AdjT(quid sid, bitLenInt q);
398 void U(quid sid, bitLenInt q, real1_f theta, real1_f phi, real1_f lambda);
399 void Mtrx(quid sid, std::vector<complex> m, bitLenInt q);
400 
401 // multi-controlled single-qubit gates
402 void MCX(quid sid, std::vector<bitLenInt> c, bitLenInt q);
403 void MCY(quid sid, std::vector<bitLenInt> c, bitLenInt q);
404 void MCZ(quid sid, std::vector<bitLenInt> c, bitLenInt q);
405 void MCH(quid sid, std::vector<bitLenInt> c, bitLenInt q);
406 void MCS(quid sid, std::vector<bitLenInt> c, bitLenInt q);
407 void MCT(quid sid, std::vector<bitLenInt> c, bitLenInt q);
408 void MCAdjS(quid sid, std::vector<bitLenInt> c, bitLenInt q);
409 void MCAdjT(quid sid, std::vector<bitLenInt> c, bitLenInt q);
410 void MCU(quid sid, std::vector<bitLenInt> c, bitLenInt q, real1_f theta, real1_f phi, real1_f lambda);
411 void MCMtrx(quid sid, std::vector<bitLenInt> c, std::vector<complex> m, bitLenInt q);
412 // multi-("anti"-) controlled single-qubits gates (that activate when all controls are |0>)
413 void MACX(quid sid, std::vector<bitLenInt> c, bitLenInt q);
414 void MACY(quid sid, std::vector<bitLenInt> c, bitLenInt q);
415 void MACZ(quid sid, std::vector<bitLenInt> c, bitLenInt q);
416 void MACH(quid sid, std::vector<bitLenInt> c, bitLenInt q);
417 void MACS(quid sid, std::vector<bitLenInt> c, bitLenInt q);
418 void MACT(quid sid, std::vector<bitLenInt> c, bitLenInt q);
419 void MACAdjS(quid sid, std::vector<bitLenInt> c, bitLenInt q);
420 void MACAdjT(quid sid, std::vector<bitLenInt> c, bitLenInt q);
421 void MACU(quid sid, std::vector<bitLenInt> c, bitLenInt q, real1_f theta, real1_f phi, real1_f lambda);
422 void MACMtrx(quid sid, std::vector<bitLenInt> c, std::vector<complex> m, bitLenInt q);
423 
427 void UCMtrx(quid sid, std::vector<bitLenInt> c, std::vector<complex> m, bitLenInt q, bitCapIntOcl p);
431 void Multiplex1Mtrx(quid sid, std::vector<bitLenInt> c, bitLenInt q, std::vector<complex> m);
432 
433 // coalesced single-qubit gates
434 void MX(quid sid, std::vector<bitLenInt> q);
435 void MY(quid sid, std::vector<bitLenInt> q);
436 void MZ(quid sid, std::vector<bitLenInt> q);
437 
438 // single-qubit rotations
439 void R(quid sid, real1_f phi, QubitPauliBasis q);
440 // multi-controlled single-qubit rotations
441 void MCR(quid sid, real1_f phi, std::vector<bitLenInt> c, QubitPauliBasis q);
442 
443 // exponential of Pauli operators
444 void Exp(quid sid, real1_f phi, std::vector<QubitPauliBasis> q);
445 // multi-controlled exponential of Pauli operators
446 void MCExp(quid sid, real1_f phi, std::vector<bitLenInt> c, std::vector<QubitPauliBasis> q);
447 
448 // swap variants
449 void SWAP(quid sid, bitLenInt qi1, bitLenInt qi2);
450 void ISWAP(quid sid, bitLenInt qi1, bitLenInt qi2);
451 void AdjISWAP(quid sid, bitLenInt qi1, bitLenInt qi2);
452 void FSim(quid sid, real1_f theta, real1_f phi, bitLenInt qi1, bitLenInt qi2);
453 void CSWAP(quid sid, std::vector<bitLenInt> c, bitLenInt qi1, bitLenInt qi2);
454 void ACSWAP(quid sid, std::vector<bitLenInt> c, bitLenInt qi1, bitLenInt qi2);
455 
456 // Schmidt decomposition
457 void Compose(quid sid1, quid sid2, std::vector<bitLenInt> q);
458 quid Decompose(quid sid, std::vector<bitLenInt> q);
459 void Dispose(quid sid, std::vector<bitLenInt> q);
460 
461 // Quantum boolean (Toffoli) operations:
462 // Two qubits input
463 void AND(quid sid, bitLenInt qi1, bitLenInt qi2, bitLenInt qo);
464 void OR(quid sid, bitLenInt qi1, bitLenInt qi2, bitLenInt qo);
465 void XOR(quid sid, bitLenInt qi1, bitLenInt qi2, bitLenInt qo);
466 void NAND(quid sid, bitLenInt qi1, bitLenInt qi2, bitLenInt qo);
467 void NOR(quid sid, bitLenInt qi1, bitLenInt qi2, bitLenInt qo);
468 void XNOR(quid sid, bitLenInt qi1, bitLenInt qi2, bitLenInt qo);
469 // One qubit, one classical bit input
470 void CLAND(quid sid, bool ci, bitLenInt qi, bitLenInt qo);
471 void CLOR(quid sid, bool ci, bitLenInt qi, bitLenInt qo);
472 void CLXOR(quid sid, bool ci, bitLenInt qi, bitLenInt qo);
473 void CLNAND(quid sid, bool ci, bitLenInt qi, bitLenInt qo);
474 void CLNOR(quid sid, bool ci, bitLenInt qi, bitLenInt qo);
475 void CLXNOR(quid sid, bool ci, bitLenInt qi, bitLenInt qo);
476 
480 void QFT(quid sid, std::vector<bitLenInt> q);
484 void IQFT(quid sid, std::vector<bitLenInt> q);
485 
486 #if ENABLE_ALU
487 // Arithmetic logic unit:
488 // Two's complement
489 void ADD(quid sid, bitCapInt a, std::vector<bitLenInt> q);
490 void SUB(quid sid, bitCapInt a, std::vector<bitLenInt> q);
491 // Overflow
492 void ADDS(quid sid, bitCapInt a, bitLenInt s, std::vector<bitLenInt> q);
493 void SUBS(quid sid, bitCapInt a, bitLenInt s, std::vector<bitLenInt> q);
494 // Controlled
495 void MCADD(quid sid, bitCapInt a, std::vector<bitLenInt> c, std::vector<bitLenInt> q);
496 void MCSUB(quid sid, bitCapInt a, std::vector<bitLenInt> c, std::vector<bitLenInt> q);
497 // In-place
498 void MUL(quid sid, bitCapInt a, std::vector<bitLenInt> q, std::vector<bitLenInt> o);
499 void DIV(quid sid, bitCapInt a, std::vector<bitLenInt> q, std::vector<bitLenInt> o);
500 // Modulo, out-of-place
501 void MULN(quid sid, bitCapInt a, bitCapInt m, std::vector<bitLenInt> q, std::vector<bitLenInt> o);
502 void DIVN(quid sid, bitCapInt a, bitCapInt m, std::vector<bitLenInt> q, std::vector<bitLenInt> o);
503 void POWN(quid sid, bitCapInt a, bitCapInt m, std::vector<bitLenInt> q, std::vector<bitLenInt> o);
504 // Controlled in-place
505 void MCMUL(quid sid, bitCapInt a, std::vector<bitLenInt> c, std::vector<bitLenInt> q, std::vector<bitLenInt> o);
506 void MCDIV(quid sid, bitCapInt a, std::vector<bitLenInt> c, std::vector<bitLenInt> q, std::vector<bitLenInt> o);
507 // Controlled modulo, out-of-place
508 void MCMULN(
509  quid sid, bitCapInt a, std::vector<bitLenInt> c, bitCapInt m, std::vector<bitLenInt> q, std::vector<bitLenInt> o);
510 void MCDIVN(
511  quid sid, bitCapInt a, std::vector<bitLenInt> c, bitCapInt m, std::vector<bitLenInt> q, std::vector<bitLenInt> o);
512 void MCPOWN(
513  quid sid, bitCapInt a, std::vector<bitLenInt> c, bitCapInt m, std::vector<bitLenInt> q, std::vector<bitLenInt> o);
514 
515 #if 0
516 // Amplitude amplification
517 void LDA(quid sid, std::vector<bitLenInt> qi, std::vector<bitLenInt> qv, std::vector<unsigned char> t);
518 void ADC(quid sid, bitLenInt s, std::vector<bitLenInt> qi, std::vector<bitLenInt> qv, std::vector<unsigned char> t);
519 void SBC(quid sid, bitLenInt s, std::vector<bitLenInt> qi, std::vector<bitLenInt> qv, std::vector<unsigned char> t);
520 void Hash(quid sid, std::vector<bitLenInt> q, std::vector<unsigned char> t);
521 #endif
522 #endif
523 
524 // Utility functions
529 bool TrySeparate1Qb(quid sid, bitLenInt qi1);
534 bool TrySeparate2Qb(quid sid, bitLenInt qi1, bitLenInt qi2);
539 bool TrySeparateTol(quid sid, std::vector<bitLenInt> q, real1_f tol);
543 void Separate(quid sid, std::vector<bitLenInt> q);
547 double GetUnitaryFidelity(quid sid);
551 void ResetUnitaryFidelity(quid sid);
555 void SetSdrp(quid sid, double sdrp);
559 void SetNcrp(quid sid, double sdrp);
563 void SetReactiveSeparate(quid sid, bool irs);
567 void SetTInjection(quid sid, bool iti);
571 void SetNoiseParameter(quid sid, double np);
575 void Normalize(quid sid);
576 
581 quid init_qneuron(quid sid, std::vector<bitLenInt> c, bitLenInt q, QNeuronActivationFn f, real1_f a, real1_f tol);
585 quid clone_qneuron(quid nid);
589 void destroy_qneuron(quid nid);
593 void set_qneuron_angles(quid nid, std::vector<real1> angles);
597 std::vector<real1> get_qneuron_angles(quid nid);
601 void set_qneuron_alpha(quid nid, real1_f alpha);
617 real1_f qneuron_predict(quid nid, bool e, bool r);
621 real1_f qneuron_unpredict(quid nid, bool e);
625 real1_f qneuron_learn_cycle(quid nid, bool e);
629 void qneuron_learn(quid nid, real1_f eta, bool e, bool r);
633 void qneuron_learn_permutation(quid nid, real1_f eta, bool e, bool r);
634 
635 // Quantum circuit objects
636 quid init_qcircuit(bool collapse, bool clifford);
639 quid qcircuit_past_light_cone(quid cid, std::set<bitLenInt> q);
640 void destroy_qcircuit(quid cid);
642 void qcircuit_swap(quid cid, bitLenInt q1, bitLenInt q2);
643 void qcircuit_append_1qb(quid cid, std::vector<complex> m, bitLenInt q);
644 void qcircuit_append_mc(quid cid, std::vector<complex> m, std::vector<bitLenInt> c, bitLenInt q, bitCapInt p);
645 void qcircuit_run(quid cid, quid sid);
646 void qcircuit_out_to_file(quid cid, std::string f);
647 void qcircuit_in_from_file(quid cid, std::string f);
648 std::string qcircuit_out_to_string(quid cid);
649 } // namespace Qrack
GLOSSARY: bitLenInt - "bit-length integer" - unsigned integer ID of qubit position in register bitCap...
Definition: complex16x2simd.hpp:25
real1_f qneuron_learn_cycle(quid nid, bool e)
Train a quantum neuron for one epoch, and also uncompute the intermediate side-effects.
Definition: wasm_api.cpp:2630
bool TrySeparate1Qb(quid sid, bitLenInt qi1)
Try to factorize a single-qubit subsystem out of "bulk" simulator state.
Definition: wasm_api.cpp:2367
real1_f ProbRdm(quid sid, bitLenInt q)
"Reduced density matrix" Z-basis expectation value of qubit
Definition: wasm_api.cpp:1802
void set_qneuron_alpha(quid nid, real1_f alpha)
Set the "leakage" parameter for "leaky" quantum neuron activation functions.
Definition: wasm_api.cpp:2582
void MCT(quid sid, std::vector< bitLenInt > c, bitLenInt q)
(External API) Controlled "T" Gate
Definition: wasm_api.cpp:1221
void CLNOR(quid sid, bool ci, bitLenInt qi, bitLenInt qo)
Definition: wasm_api.cpp:1738
void X(quid sid, bitLenInt q)
(External API) "X" Gate
Definition: wasm_api.cpp:1035
quid qcircuit_past_light_cone(quid cid, std::set< bitLenInt > q)
Definition: wasm_api.cpp:2728
void MACS(quid sid, std::vector< bitLenInt > c, bitLenInt q)
(External API) "Anti-"Controlled "S" Gate
Definition: wasm_api.cpp:1311
void Y(quid sid, bitLenInt q)
(External API) "Y" Gate
Definition: wasm_api.cpp:1044
void ADD(quid sid, bitCapInt a, std::vector< bitLenInt > q)
Definition: wasm_api.cpp:2178
void MCMULN(quid sid, bitCapInt a, std::vector< bitLenInt > c, bitCapInt m, std::vector< bitLenInt > q, std::vector< bitLenInt > o)
Definition: wasm_api.cpp:2297
void MY(quid sid, std::vector< bitLenInt > q)
(External API) Multiple "Y" Gate
Definition: wasm_api.cpp:1411
void XOR(quid sid, bitLenInt qi1, bitLenInt qi2, bitLenInt qo)
Definition: wasm_api.cpp:1690
void SetPermutation(quid sid, bitCapInt p)
Set bit string permutation eigenstate of simulator instance.
Definition: wasm_api.cpp:1026
void MCADD(quid sid, bitCapInt a, std::vector< bitLenInt > c, std::vector< bitLenInt > q)
Definition: wasm_api.cpp:2199
real1_f JointEnsembleProbability(quid sid, std::vector< QubitPauliBasis > q)
Overall probability of any odd permutation of the masked set of bits.
Definition: wasm_api.cpp:933
void ISWAP(quid sid, bitLenInt qi1, bitLenInt qi2)
Definition: wasm_api.cpp:1570
quid init()
"Quasi-default constructor" (for an empty simulator)
Definition: wasm_api.cpp:640
void Compose(quid sid1, quid sid2, std::vector< bitLenInt > q)
Definition: wasm_api.cpp:1600
void allocateQubit(quid sid, bitLenInt qid)
Allocate new qubit with ID.
Definition: wasm_api.cpp:958
void POWN(quid sid, bitCapInt a, bitCapInt m, std::vector< bitLenInt > q, std::vector< bitLenInt > o)
Definition: wasm_api.cpp:2260
void MCPOWN(quid sid, bitCapInt a, std::vector< bitLenInt > c, bitCapInt m, std::vector< bitLenInt > q, std::vector< bitLenInt > o)
Definition: wasm_api.cpp:2325
real1_f MatrixVariance(quid sid, std::vector< QubitMatrixBasis > q)
Get the single-qubit (2x2) operator variance for the array of qubits and bases.
Definition: wasm_api.cpp:2043
void set_qneuron_angles(quid nid, std::vector< real1 > angles)
Set the (RY-rotation) angle parameters for each permutation of quantum neuron input qubits.
Definition: wasm_api.cpp:2550
void qcircuit_append_1qb(quid cid, std::vector< complex > m, bitLenInt q)
void SetNoiseParameter(quid sid, double np)
Set noise parameter (for QInterfaceNoisy)
Definition: wasm_api.cpp:2447
quid init_qcircuit_clone(quid cid)
Definition: wasm_api.cpp:2724
void SetReactiveSeparate(quid sid, bool irs)
Turn off/on "reactive separation" feature (for less/more aggressive automatic state factorization)
Definition: wasm_api.cpp:2435
bitLenInt get_qcircuit_qubit_count(quid cid)
Definition: wasm_api.cpp:2739
void MCH(quid sid, std::vector< bitLenInt > c, bitLenInt q)
(External API) Controlled "H" Gate
Definition: wasm_api.cpp:1200
void SetSdrp(quid sid, double sdrp)
Set "Schmidt decomposition rounding parameter" (SDRP) value (see https://arxiv.org/abs/2304....
Definition: wasm_api.cpp:2423
void qcircuit_out_to_file(quid cid, std::string f)
Definition: wasm_api.cpp:2793
void qstabilizer_out_to_file(quid sid, std::string f)
Output stabilizer simulation tableau to file (or raise exception for "get_error()" if incompatible si...
Definition: wasm_api.cpp:848
void SY(quid sid, bitLenInt q)
(External API) Square root of Y gate
Definition: wasm_api.cpp:1089
void MUL(quid sid, bitCapInt a, std::vector< bitLenInt > q, std::vector< bitLenInt > o)
Definition: wasm_api.cpp:2220
bool TrySeparateTol(quid sid, std::vector< bitLenInt > q, real1_f tol)
Try to factorize a qubit subsystem out of "bulk" simulator state.
Definition: wasm_api.cpp:2379
void MCX(quid sid, std::vector< bitLenInt > c, bitLenInt q)
(External API) Controlled "X" Gate
Definition: wasm_api.cpp:1173
void XNOR(quid sid, bitLenInt qi1, bitLenInt qi2, bitLenInt qo)
Definition: wasm_api.cpp:1708
void MCY(quid sid, std::vector< bitLenInt > c, bitLenInt q)
(External API) Controlled "Y" Gate
Definition: wasm_api.cpp:1182
void T(quid sid, bitLenInt q)
(External API) "T" Gate
Definition: wasm_api.cpp:1098
void IQFT(quid sid, std::vector< bitLenInt > q)
(Inverse) Quantum Fourier Transform
Definition: wasm_api.cpp:2165
void MACT(quid sid, std::vector< bitLenInt > c, bitLenInt q)
(External API) "Anti-"Controlled "T" Gate
Definition: wasm_api.cpp:1320
void MACU(quid sid, std::vector< bitLenInt > c, bitLenInt q, real1_f theta, real1_f phi, real1_f lambda)
(External API) Controlled 3-parameter unitary gate
Definition: wasm_api.cpp:1347
void MCAdjS(quid sid, std::vector< bitLenInt > c, bitLenInt q)
(External API) Controlled Inverse "S" Gate
Definition: wasm_api.cpp:1230
void Mtrx(quid sid, std::vector< complex > m, bitLenInt q)
(External API) 2x2 complex matrix unitary gate
Definition: wasm_api.cpp:1152
void MCDIV(quid sid, bitCapInt a, std::vector< bitLenInt > c, std::vector< bitLenInt > q, std::vector< bitLenInt > o)
Definition: wasm_api.cpp:2284
void CSWAP(quid sid, std::vector< bitLenInt > c, bitLenInt qi1, bitLenInt qi2)
Definition: wasm_api.cpp:1588
bool Measure(quid sid, std::vector< QubitPauliBasis > q)
Each in its specified Pauli basis, collapse an ensemble of qubits jointly via measurement.
Definition: wasm_api.cpp:1527
real1_f UnitaryVariance(quid sid, std::vector< QubitU3Basis > q)
Get the single-qubit (3-parameter) operator variance for the array of qubits and bases.
Definition: wasm_api.cpp:2011
quid clone_qneuron(quid nid)
"Clone" a quantum neuron (which is a classical state)
Definition: wasm_api.cpp:2503
void AdjISWAP(quid sid, bitLenInt qi1, bitLenInt qi2)
Definition: wasm_api.cpp:1576
void UCMtrx(quid sid, std::vector< bitLenInt > c, std::vector< complex > m, bitLenInt q, bitCapIntOcl p)
Multi-controlled gate that activates only for the specified permutation of controls,...
Definition: wasm_api.cpp:1371
real1_f FactorizedExpectation(quid sid, std::vector< QubitIntegerExpVar > q)
Expectation value for bit-string integer from group of qubits with per-qubit integer expectation valu...
Definition: wasm_api.cpp:1898
real1_f UnitaryExpectationEigenVal(quid sid, std::vector< QubitU3BasisEigenVal > q)
Get the single-qubit (3-parameter) operator expectation value for the array of qubits and bases.
Definition: wasm_api.cpp:2073
void MCMUL(quid sid, bitCapInt a, std::vector< bitLenInt > c, std::vector< bitLenInt > q, std::vector< bitLenInt > o)
Definition: wasm_api.cpp:2271
void seed(quid sid, unsigned s)
"Seed" random number generator (if pseudo-random Mersenne twister is in use)
Definition: wasm_api.cpp:833
real1_f FactorizedVariance(quid sid, std::vector< QubitIntegerExpVar > q)
Variance for bit-string integer from group of qubits with per-qubit integer variance.
Definition: wasm_api.cpp:1906
void ACSWAP(quid sid, std::vector< bitLenInt > c, bitLenInt qi1, bitLenInt qi2)
Definition: wasm_api.cpp:1594
void MACZ(quid sid, std::vector< bitLenInt > c, bitLenInt q)
(External API) "Anti-"Controlled "Z" Gate
Definition: wasm_api.cpp:1290
void MCAdjT(quid sid, std::vector< bitLenInt > c, bitLenInt q)
(External API) Controlled Inverse "T" Gate
Definition: wasm_api.cpp:1239
void U(quid sid, bitLenInt q, real1_f theta, real1_f phi, real1_f lambda)
(External API) 3-parameter unitary gate
Definition: wasm_api.cpp:1143
bool M(quid sid, bitLenInt q)
Measure single qubit (according to Born rules) and return the result.
Definition: wasm_api.cpp:1499
size_t random_choice(quid sid, std::vector< real1 > p)
Select from a distribution of "p.size()" count of elements according to the discrete probabilities in...
Definition: wasm_api.cpp:882
void AdjSY(quid sid, bitLenInt q)
(External API) Inverse square root of Y gate
Definition: wasm_api.cpp:1125
void CLNAND(quid sid, bool ci, bitLenInt qi, bitLenInt qo)
Definition: wasm_api.cpp:1732
std::complex< real1 > complex
Definition: qrack_types.hpp:128
void FSim(quid sid, real1_f theta, real1_f phi, bitLenInt qi1, bitLenInt qi2)
Definition: wasm_api.cpp:1582
real1_f VarianceRdm(quid sid, std::vector< bitLenInt > q, bool r)
"Reduced density matrix" variance for bit-string integer equivalent of specified arbitrary group of q...
Definition: wasm_api.cpp:1874
void MCMtrx(quid sid, std::vector< bitLenInt > c, std::vector< complex > m, bitLenInt q)
(External API) Controlled 2x2 complex matrix unitary gate
Definition: wasm_api.cpp:1257
void PhaseParity(quid sid, real1_f lambda, std::vector< bitLenInt > q)
Applies e^(i*angle) phase factor to all combinations of bits with odd parity, based upon permutations...
Definition: wasm_api.cpp:904
real1_f FactorizedVarianceFp(quid sid, std::vector< QubitRealExpVar > q)
Variance for bit-string integer from group of qubits with per-qubit real1 variance.
Definition: wasm_api.cpp:1959
real1_f FactorizedExpectationFp(quid sid, std::vector< QubitRealExpVar > q)
Expectation value for bit-string integer from group of qubits with per-qubit real1 expectation value.
Definition: wasm_api.cpp:1951
real1_f Variance(quid sid, std::vector< bitLenInt > q)
Variance for bit-string integer equivalent of specified arbitrary group of qubits.
Definition: wasm_api.cpp:1868
void qneuron_learn_permutation(quid nid, real1_f eta, bool e, bool r)
Train a quantum neuron for one epoch, assuming that the input state is a Z-basis eigenstate.
Definition: wasm_api.cpp:2653
void AdjT(quid sid, bitLenInt q)
(External API) Inverse "T" Gate
Definition: wasm_api.cpp:1134
bitCapInt MAll(quid sid)
Measure all qubits (according to Born rules) and return the result as a bit string (integer).
Definition: wasm_api.cpp:1518
void Exp(quid sid, real1_f phi, std::vector< QubitPauliBasis > q)
(External API) Exponentiation of Pauli operators
Definition: wasm_api.cpp:1447
real1_f PermutationProb(quid sid, std::vector< QubitIndexState > q)
Probability of specified (single) permutation of any arbitrary group of qubits.
Definition: wasm_api.cpp:1824
void S(quid sid, bitLenInt q)
(External API) "S" Gate
Definition: wasm_api.cpp:1071
void MACX(quid sid, std::vector< bitLenInt > c, bitLenInt q)
(External API) "Anti-"Controlled "X" Gate
Definition: wasm_api.cpp:1272
void destroy(quid sid)
"Destroy" or release simulator allocation
Definition: wasm_api.cpp:820
QNeuronActivationFn get_qneuron_activation_fn(quid nid)
Get the activation function for a quantum neuron.
Definition: wasm_api.cpp:2600
quid init_qneuron(quid sid, std::vector< bitLenInt > c, bitLenInt q, QNeuronActivationFn f, real1_f a, real1_f tol)
Initialize a "quantum neuron" that takes a list of qubit "controls" for input and acts on a single "t...
Definition: wasm_api.cpp:2460
void destroy_qcircuit(quid cid)
Definition: wasm_api.cpp:2730
void MACH(quid sid, std::vector< bitLenInt > c, bitLenInt q)
(External API) "Anti-"Controlled "H" Gate
Definition: wasm_api.cpp:1299
void set_qneuron_activation_fn(quid nid, QNeuronActivationFn f)
Set the activation function for a quantum neuron.
Definition: wasm_api.cpp:2594
real1_f MatrixExpectation(quid sid, std::vector< QubitMatrixBasis > q)
Get the single-qubit (2x2) operator expectation value for the array of qubits and bases.
Definition: wasm_api.cpp:2038
uint64_t quid
Definition: wasm_api.hpp:27
void ADDS(quid sid, bitCapInt a, bitLenInt s, std::vector< bitLenInt > q)
Definition: wasm_api.cpp:2188
void Separate(quid sid, std::vector< bitLenInt > q)
Force (inexact) factorization of qubit subsystem out of "bulk" simulator state.
Definition: wasm_api.cpp:2385
void H(quid sid, bitLenInt q)
(External API) Walsh-Hadamard transform applied for simulator ID and qubit ID
Definition: wasm_api.cpp:1062
void QFT(quid sid, std::vector< bitLenInt > q)
Quantum Fourier Transform.
Definition: wasm_api.cpp:2154
void qcircuit_in_from_file(quid cid, std::string f)
Definition: wasm_api.cpp:2803
void PhaseRootN(quid sid, bitLenInt p, std::vector< bitLenInt > q)
Applies a -2 * PI_R1 / (2^N) phase rotation to each qubit.
Definition: wasm_api.cpp:906
void OR(quid sid, bitLenInt qi1, bitLenInt qi2, bitLenInt qo)
Definition: wasm_api.cpp:1684
real1_f PauliExpectation(quid sid, std::vector< QubitPauliBasis > q)
Pauli operator expectation value for the array of qubits and bases.
Definition: wasm_api.cpp:2147
bool ForceM(quid sid, bitLenInt q, bool r)
(PSEUDO-QUANTUM:) Force measurement result of single qubit (and return the result)
Definition: wasm_api.cpp:1508
void SX(quid sid, bitLenInt q)
(External API) Square root of X gate
Definition: wasm_api.cpp:1080
void R(quid sid, real1_f phi, QubitPauliBasis q)
(External API) Rotation around Pauli axes
Definition: wasm_api.cpp:1429
void DIV(quid sid, bitCapInt a, std::vector< bitLenInt > q, std::vector< bitLenInt > o)
Definition: wasm_api.cpp:2230
void MACAdjT(quid sid, std::vector< bitLenInt > c, bitLenInt q)
(External API) "Anti-"Controlled Inverse "T" Gate
Definition: wasm_api.cpp:1338
void MCDIVN(quid sid, bitCapInt a, std::vector< bitLenInt > c, bitCapInt m, std::vector< bitLenInt > q, std::vector< bitLenInt > o)
Definition: wasm_api.cpp:2311
void qcircuit_swap(quid cid, bitLenInt q1, bitLenInt q2)
Definition: wasm_api.cpp:2744
real1_f UnitaryVarianceEigenVal(quid sid, std::vector< QubitU3BasisEigenVal > q)
Get the single-qubit (3-parameter) operator variance for the array of qubits and bases.
Definition: wasm_api.cpp:2081
void CLOR(quid sid, bool ci, bitLenInt qi, bitLenInt qo)
Definition: wasm_api.cpp:1720
real1_f PermutationExpectationRdm(quid sid, std::vector< bitLenInt > q, bool r)
"Reduced density matrix" expectation value for bit-string integer equivalent of specified arbitrary g...
Definition: wasm_api.cpp:1860
bitLenInt num_qubits(quid sid)
Total count of qubits in simulator instance.
Definition: wasm_api.cpp:1020
void CLXOR(quid sid, bool ci, bitLenInt qi, bitLenInt qo)
Definition: wasm_api.cpp:1726
real1_f PermutationProbRdm(quid sid, std::vector< QubitIndexState > q, bool r)
"Reduced density matrix" probability of specified (single) permutation of any arbitrary group of qubi...
Definition: wasm_api.cpp:1830
void MCSUB(quid sid, bitCapInt a, std::vector< bitLenInt > c, std::vector< bitLenInt > q)
Definition: wasm_api.cpp:2209
void qstabilizer_in_from_file(quid sid, std::string f)
Initialize stabilizer simulation from a tableau file (or raise exception for "get_error()" if incompa...
Definition: wasm_api.cpp:862
void AdjSX(quid sid, bitLenInt q)
(External API) Inverse square root of X gate
Definition: wasm_api.cpp:1116
void MZ(quid sid, std::vector< bitLenInt > q)
(External API) Multiple "Z" Gate
Definition: wasm_api.cpp:1420
void MULN(quid sid, bitCapInt a, bitCapInt m, std::vector< bitLenInt > q, std::vector< bitLenInt > o)
Definition: wasm_api.cpp:2240
real1_f Prob(quid sid, bitLenInt q)
Z-basis expectation value of qubit.
Definition: wasm_api.cpp:1796
void destroy_qneuron(quid nid)
"Destroy" or release simulator allocation
Definition: wasm_api.cpp:2540
real1_f FactorizedExpectationRdm(quid sid, std::vector< QubitIntegerExpVar > q, bool r)
"Reduced density matrix" Expectation value for bit-string integer from group of qubits with per-qubit...
Definition: wasm_api.cpp:1915
real1_f get_qneuron_alpha(quid nid)
Get the "leakage" parameter for "leaky" quantum neuron activation functions.
Definition: wasm_api.cpp:2588
void SetTInjection(quid sid, bool iti)
Turn off/on "T-injection" feature (for "near-Clifford" simulation with RZ gates)
Definition: wasm_api.cpp:2441
float real1_f
Definition: qrack_types.hpp:95
void qneuron_learn(quid nid, real1_f eta, bool e, bool r)
Train a quantum neuron for one epoh (without uncomputing any intermediate side-effects)
Definition: wasm_api.cpp:2642
quid Decompose(quid sid, std::vector< bitLenInt > q)
Definition: wasm_api.cpp:1630
void Dispose(quid sid, std::vector< bitLenInt > q)
Definition: wasm_api.cpp:1656
void CLAND(quid sid, bool ci, bitLenInt qi, bitLenInt qo)
Definition: wasm_api.cpp:1714
std::string qcircuit_out_to_string(quid cid)
Definition: wasm_api.cpp:2813
quid init_qcircuit(bool collapse, bool clifford)
Definition: wasm_api.cpp:2664
quid qcircuit_inverse(quid cid)
Definition: wasm_api.cpp:2726
void qcircuit_append_mc(quid cid, std::vector< complex > m, std::vector< bitLenInt > c, bitLenInt q, bitCapInt p)
void NAND(quid sid, bitLenInt qi1, bitLenInt qi2, bitLenInt qo)
Definition: wasm_api.cpp:1696
real1_f qneuron_unpredict(quid nid, bool e)
Perform the inverse of quantum neuron inference (for "uncomputation")
Definition: wasm_api.cpp:2618
void SetNcrp(quid sid, double sdrp)
Set "Near-Clifford rounding parameter".
Definition: wasm_api.cpp:2429
void set_concurrency(quid sid, unsigned p)
Set CPU concurrency (if build isn't serial)
Definition: wasm_api.cpp:842
bool TrySeparate2Qb(quid sid, bitLenInt qi1, bitLenInt qi2)
Try to factorize a two-qubit subsystem out of "bulk" simulator state.
Definition: wasm_api.cpp:2373
double GetUnitaryFidelity(quid sid)
Report fidelity for "Schmidt decomposition rounding parameter" (SDRP) and "near-Clifford rounding".
Definition: wasm_api.cpp:2411
void ResetUnitaryFidelity(quid sid)
Reset fidelity to 1 for "Schmidt decomposition rounding parameter" (SDRP) and "near-Clifford rounding...
Definition: wasm_api.cpp:2417
void ResetAll(quid sid)
Set simulator to |0> permutation state.
Definition: wasm_api.cpp:949
void MACAdjS(quid sid, std::vector< bitLenInt > c, bitLenInt q)
(External API) "Anti-"Controlled Inverse "S" Gate
Definition: wasm_api.cpp:1329
void DIVN(quid sid, bitCapInt a, bitCapInt m, std::vector< bitLenInt > q, std::vector< bitLenInt > o)
Definition: wasm_api.cpp:2250
std::vector< real1 > get_qneuron_angles(quid nid)
Get the (RY-rotation) angle parameters for each permutation of quantum neuron input qubits.
Definition: wasm_api.cpp:2562
real1_f qneuron_predict(quid nid, bool e, bool r)
Infer quantum neuron output from inputs (after training)
Definition: wasm_api.cpp:2606
void CLXNOR(quid sid, bool ci, bitLenInt qi, bitLenInt qo)
Definition: wasm_api.cpp:1744
void SWAP(quid sid, bitLenInt qi1, bitLenInt qi2)
Definition: wasm_api.cpp:1564
void NOR(quid sid, bitLenInt qi1, bitLenInt qi2, bitLenInt qo)
Definition: wasm_api.cpp:1702
void SUBS(quid sid, bitCapInt a, bitLenInt s, std::vector< bitLenInt > q)
Definition: wasm_api.cpp:2193
std::vector< long long unsigned int > MeasureShots(quid sid, std::vector< bitLenInt > q, unsigned s)
Repeat (Z-basis) measurement of a set of qubits for a count of "shots" (without collapsing the simula...
Definition: wasm_api.cpp:1538
quid init_count(bitLenInt q, bool dm)
"Default optimal" (BQP-complete-targeted) simulator type initialization (with "direct memory" option)
Definition: wasm_api.cpp:645
void AND(quid sid, bitLenInt qi1, bitLenInt qi2, bitLenInt qo)
Definition: wasm_api.cpp:1678
real1_f PauliVariance(quid sid, std::vector< QubitPauliBasis > q)
Pauli operator variance for the array of qubits and bases.
Definition: wasm_api.cpp:2152
real1_f PermutationExpectation(quid sid, std::vector< bitLenInt > q)
Expectation value for bit-string integer equivalent of specified arbitrary group of qubits.
Definition: wasm_api.cpp:1851
void MCU(quid sid, std::vector< bitLenInt > c, bitLenInt q, real1_f theta, real1_f phi, real1_f lambda)
(External API) Controlled 3-parameter unitary gate
Definition: wasm_api.cpp:1248
void Normalize(quid sid)
Normalize the state (which should never be necessary unless Decompose() is "abused")
Definition: wasm_api.cpp:2453
void AdjS(quid sid, bitLenInt q)
(External API) Inverse "S" Gate
Definition: wasm_api.cpp:1107
void MACY(quid sid, std::vector< bitLenInt > c, bitLenInt q)
(External API) "Anti-"Controlled "Y" Gate
Definition: wasm_api.cpp:1281
Pauli
Enumerated list of Pauli bases.
Definition: pauli.hpp:19
std::vector< real1 > ProbAll(quid sid, std::vector< bitLenInt > q)
Get the probabilities of all permutations of the requested subset of qubits.
Definition: wasm_api.cpp:1750
real1_f FactorizedVarianceFpRdm(quid sid, std::vector< QubitRealExpVar > q, bool r)
"Reduced density matrix" variance for bit-string integer from group of qubits with per-qubit real1 va...
Definition: wasm_api.cpp:1977
void SUB(quid sid, bitCapInt a, std::vector< bitLenInt > q)
Definition: wasm_api.cpp:2183
real1_f FactorizedExpectationFpRdm(quid sid, std::vector< QubitRealExpVar > q, bool r)
"Reduced density matrix" Expectation value for bit-string integer from group of qubits with per-qubit...
Definition: wasm_api.cpp:1968
void qcircuit_run(quid cid, quid sid)
Definition: wasm_api.cpp:2787
void MCR(quid sid, real1_f phi, std::vector< bitLenInt > c, QubitPauliBasis q)
(External API) Controlled rotation around Pauli axes
Definition: wasm_api.cpp:1438
void MCExp(quid sid, real1_f phi, std::vector< bitLenInt > c, std::vector< QubitPauliBasis > q)
(External API) Controlled exponentiation of Pauli operators
Definition: wasm_api.cpp:1473
void Z(quid sid, bitLenInt q)
(External API) "Z" Gate
Definition: wasm_api.cpp:1053
real1_f UnitaryExpectation(quid sid, std::vector< QubitU3Basis > q)
Get the single-qubit (3-parameter) operator expectation value for the array of qubits and bases.
Definition: wasm_api.cpp:2006
QNeuronActivationFn
Enumerated list of activation functions.
Definition: qneuron_activation_function.hpp:19
real1_f FactorizedVarianceRdm(quid sid, std::vector< QubitIntegerExpVar > q, bool r)
"Reduced density matrix" variance for bit-string integer from group of qubits with per-qubit integer ...
Definition: wasm_api.cpp:1924
void Multiplex1Mtrx(quid sid, std::vector< bitLenInt > c, bitLenInt q, std::vector< complex > m)
Multi-controlled, single-target multiplexer gate.
Definition: wasm_api.cpp:1383
quid init_count_type(bitLenInt q, bool tn, bool md, bool sd, bool sh, bool bdt, bool pg, bool nw, bool hy, bool oc, bool hp)
Options for simulator type in initialization (any set of options theoretically functions together): t...
Definition: wasm_api.cpp:555
quid init_qbdd_count(bitLenInt q)
"Default optimal" (BQP-complete-targeted) simulator type initialization (with "direct memory" option)
Definition: wasm_api.cpp:774
void MX(quid sid, std::vector< bitLenInt > q)
(External API) Multiple "X" Gate
Definition: wasm_api.cpp:1402
void MACMtrx(quid sid, std::vector< bitLenInt > c, std::vector< complex > m, bitLenInt q)
(External API) Controlled 2x2 complex matrix unitary gate
Definition: wasm_api.cpp:1356
quid init_clone(quid sid)
"Clone" simulator (no-clone theorem does not apply to classical simulation)
Definition: wasm_api.cpp:690
void MCS(quid sid, std::vector< bitLenInt > c, bitLenInt q)
(External API) Controlled "S" Gate
Definition: wasm_api.cpp:1212
real1_f MatrixExpectationEigenVal(quid sid, std::vector< QubitMatrixBasisEigenVal > q)
Get the single-qubit (2x2) operator expectation value for the array of qubits and bases.
Definition: wasm_api.cpp:2115
void MCZ(quid sid, std::vector< bitLenInt > c, bitLenInt q)
(External API) Controlled "Z" Gate
Definition: wasm_api.cpp:1191
real1_f MatrixVarianceEigenVal(quid sid, std::vector< QubitMatrixBasisEigenVal > q)
Get the single-qubit (2x2) operator variance for the array of qubits and bases.
Definition: wasm_api.cpp:2123
bool release(quid sid, bitLenInt q)
Release qubit ID.
Definition: wasm_api.cpp:988
MICROSOFT_QUANTUM_DECL void SBC(_In_ uintq sid, uintq s, _In_ uintq ni, _In_reads_(ni) uintq *qi, _In_ uintq nv, _In_reads_(nv) uintq *qv, unsigned char *t)
Definition: pinvoke_api.cpp:3322
MICROSOFT_QUANTUM_DECL void Hash(_In_ uintq sid, _In_ uintq n, _In_reads_(n) uintq *q, unsigned char *t)
Definition: pinvoke_api.cpp:3335
MICROSOFT_QUANTUM_DECL void ADC(_In_ uintq sid, uintq s, _In_ uintq ni, _In_reads_(ni) uintq *qi, _In_ uintq nv, _In_reads_(nv) uintq *qv, unsigned char *t)
Definition: pinvoke_api.cpp:3309
MICROSOFT_QUANTUM_DECL void LDA(_In_ uintq sid, _In_ uintq ni, _In_reads_(ni) uintq *qi, _In_ uintq nv, _In_reads_(nv) uintq *qv, unsigned char *t)
Definition: pinvoke_api.cpp:3296
#define bitLenInt
Definition: qrack_types.hpp:38
#define bitCapInt
Definition: qrack_types.hpp:62
#define bitCapIntOcl
Definition: qrack_types.hpp:50
Definition: wasm_api.hpp:29
QubitIndexState(bitLenInt q, bool v)
Definition: wasm_api.hpp:32
bool val
Definition: wasm_api.hpp:31
bitLenInt qid
Definition: wasm_api.hpp:30
Definition: wasm_api.hpp:40
QubitIntegerExpVar(bitLenInt q, bitCapInt v)
Definition: wasm_api.hpp:43
bitLenInt qid
Definition: wasm_api.hpp:41
bitCapInt val
Definition: wasm_api.hpp:42
Definition: wasm_api.hpp:124
complex b[4U]
Definition: wasm_api.hpp:126
real1_f e[2U]
Definition: wasm_api.hpp:127
bitLenInt qid
Definition: wasm_api.hpp:125
QubitMatrixBasisEigenVal(bitLenInt q, std::vector< complex > basis, std::vector< real1_f > ex)
Definition: wasm_api.hpp:128
Definition: wasm_api.hpp:88
complex b[4U]
Definition: wasm_api.hpp:90
QubitMatrixBasis(bitLenInt q, std::vector< complex > basis)
Definition: wasm_api.hpp:91
bitLenInt qid
Definition: wasm_api.hpp:89
Definition: wasm_api.hpp:62
Pauli b
Definition: wasm_api.hpp:64
QubitPauliBasis(bitLenInt q, Pauli basis)
Definition: wasm_api.hpp:65
bitLenInt qid
Definition: wasm_api.hpp:63
Definition: wasm_api.hpp:51
QubitRealExpVar(bitLenInt q, real1_f v)
Definition: wasm_api.hpp:54
bitLenInt qid
Definition: wasm_api.hpp:52
real1_f val
Definition: wasm_api.hpp:53
Definition: wasm_api.hpp:103
bitLenInt qid
Definition: wasm_api.hpp:104
QubitU3BasisEigenVal(bitLenInt q, std::vector< real1_f > basis, std::vector< real1_f > ex)
Definition: wasm_api.hpp:107
real1_f e[2U]
Definition: wasm_api.hpp:106
real1_f b[3U]
Definition: wasm_api.hpp:105
Definition: wasm_api.hpp:73
bitLenInt qid
Definition: wasm_api.hpp:74
real1_f b[3U]
Definition: wasm_api.hpp:75
QubitU3Basis(bitLenInt q, std::vector< real1_f > basis)
Definition: wasm_api.hpp:76