The Quantum Exact Simulation Toolkit v4.0.0
Loading...
Searching...
No Matches
densitymatrix.cpp
1/** @file
2 * Integration tests which combine many QuEST API functions in
3 * a single test using Quregs which are too large to validate via
4 * serial replication (like the unit tests perform), which instead
5 * utilise known scale-invariant analytic properties.
6 *
7 * @author Tyson Jones
8 */
9
10#include "quest/include/quest.h"
11
12#include <catch2/catch_test_macros.hpp>
13#include <catch2/matchers/catch_matchers_floating_point.hpp>
14
15#include "tests/utils/macros.hpp"
16#include "tests/utils/cache.hpp"
17#include "tests/utils/compare.hpp"
18#include "tests/utils/random.hpp"
19
20#include <algorithm>
21#include <string>
22#include <vector>
23#include <tuple>
24
25using std::vector;
26using std::tuple;
27using namespace Catch::Matchers;
28
29
30#define TEST_TAG \
31 LABEL_INTEGRATION_TAG
32
33
34
35void testDensityMatrixEvolution(Qureg psi, Qureg rho) {
36 DEMAND( psi.numQubits == rho.numQubits );
37
38 // set ||rho>> = |psi><psi|
40 initPureState(rho, psi);
41
42 // we will check all alculations produced within 'eps' of expected
43 qreal eps = std::max({(qreal) 1E-5, getTestAbsoluteEpsilon()});
44 REQUIRE_THAT( calcPurity(rho), WithinAbs(1, eps) );
45 REQUIRE_THAT( calcPurity(psi), WithinAbs(1, eps) );
46 REQUIRE_THAT( calcTotalProb(rho), WithinAbs(1, eps) );
47 REQUIRE_THAT( calcTotalProb(psi), WithinAbs(1, eps) );
48 REQUIRE_THAT( calcFidelity(rho, psi), WithinAbs(1, eps) );
49
50 // maximum size of tested any-target operators
51 int maxNumCompMatrTargs = std::min({6, (int) psi.logNumAmpsPerNode, (int) rho.logNumAmpsPerNode});
52 int maxNumDiagMatrTargs = std::min({8, psi.numQubits});
53 int maxNumPauliStrTargs = psi.numQubits;
54 int maxNumPauliGadTargs = psi.numQubits;
55 int maxNumPhaseGadTargs = psi.numQubits;
56
57 // below, we will apply every operation which invokes a unique backend
58 // function, with random parameters which (if they were deterministic
59 // and exhaustive) should access every edge-case including all templated
60 // optimisations. Each operation (except a few anomalously demanding ones)
61 // will be repeated with re-randomised parameters a total of 'numReps' times
62 int numReps = 10 * psi.numQubits;
63
64 vector<int> ctrls;
65 vector<int> targs;
66 vector<int> states;
67
68 // apply CompMatr1
69 for (int r=0; r<numReps; r++) {
70 auto [ctrls,states,targs] = getRandomVariNumCtrlsStatesTargs(psi.numQubits, 1,1);
72 applyMultiStateControlledCompMatr1(psi, ctrls.data(), states.data(), ctrls.size(), targs[0], matr);
73 applyMultiStateControlledCompMatr1(rho, ctrls.data(), states.data(), ctrls.size(), targs[0], matr);
74 }
75
76 // apply CompMatr2
77 for (int r=0; r<numReps; r++) {
78 auto [ctrls,states,targs] = getRandomVariNumCtrlsStatesTargs(psi.numQubits, 2,2);
80 applyMultiStateControlledCompMatr2(psi, ctrls.data(), states.data(), ctrls.size(), targs[0], targs[1], matr);
81 applyMultiStateControlledCompMatr2(rho, ctrls.data(), states.data(), ctrls.size(), targs[0], targs[1], matr);
82 }
83
84 // apply CompMatr
85 for (int r=0; r<numReps; r++) {
86 auto [ctrls,states,targs] = getRandomVariNumCtrlsStatesTargs(psi.numQubits, 1,maxNumCompMatrTargs);
87 CompMatr matr = createCompMatr(targs.size());
88 setCompMatr(matr, getRandomUnitary(targs.size()));
89 applyMultiStateControlledCompMatr(psi, ctrls.data(), states.data(), ctrls.size(), targs.data(), targs.size(), matr);
90 applyMultiStateControlledCompMatr(rho, ctrls.data(), states.data(), ctrls.size(), targs.data(), targs.size(), matr);
91 destroyCompMatr(matr);
92 }
93
94 // apply DiagMatr1
95 for (int r=0; r<numReps; r++) {
96 auto [ctrls,states,targs] = getRandomVariNumCtrlsStatesTargs(psi.numQubits, 1,1);
97 DiagMatr1 matr = getDiagMatr1(getDiagonals(getRandomDiagonalUnitary(1)));
98 applyMultiStateControlledDiagMatr1(psi, ctrls.data(), states.data(), ctrls.size(), targs[0], matr);
99 applyMultiStateControlledDiagMatr1(rho, ctrls.data(), states.data(), ctrls.size(), targs[0], matr);
100 }
101
102 // apply DiagMatr2
103 for (int r=0; r<numReps; r++) {
104 auto [ctrls,states,targs] = getRandomVariNumCtrlsStatesTargs(psi.numQubits, 2,2);
105 DiagMatr2 matr = getDiagMatr2(getDiagonals(getRandomDiagonalUnitary(2)));
106 applyMultiStateControlledDiagMatr2(psi, ctrls.data(), states.data(), ctrls.size(), targs[0], targs[1], matr);
107 applyMultiStateControlledDiagMatr2(rho, ctrls.data(), states.data(), ctrls.size(), targs[0], targs[1], matr);
108 }
109
110 // apply DiagMatr
111 for (int r=0; r<numReps; r++) {
112 auto [ctrls,states,targs] = getRandomVariNumCtrlsStatesTargs(psi.numQubits, 1,maxNumDiagMatrTargs);
113 DiagMatr matr = createDiagMatr(targs.size());
114 setDiagMatr(matr, getDiagonals(getRandomDiagonalUnitary(targs.size())));
115 applyMultiStateControlledDiagMatr(psi, ctrls.data(), states.data(), ctrls.size(), targs.data(), targs.size(), matr);
116 applyMultiStateControlledDiagMatr(rho, ctrls.data(), states.data(), ctrls.size(), targs.data(), targs.size(), matr);
117 destroyDiagMatr(matr);
118 }
119
120 // apply DiagMatr raised to power (real to retain unitarity)
121 for (int r=0; r<numReps; r++) {
122 auto [ctrls,states,targs] = getRandomVariNumCtrlsStatesTargs(psi.numQubits, 1,maxNumDiagMatrTargs);
123 DiagMatr matr = createDiagMatr(targs.size());
124 qreal expo = getRandomReal(0, 10);
125 setDiagMatr(matr, getDiagonals(getRandomDiagonalUnitary(targs.size())));
126 applyMultiStateControlledDiagMatrPower(psi, ctrls.data(), states.data(), ctrls.size(), targs.data(), targs.size(), matr, expo);
127 applyMultiStateControlledDiagMatrPower(rho, ctrls.data(), states.data(), ctrls.size(), targs.data(), targs.size(), matr, expo);
128 destroyDiagMatr(matr);
129 }
130
131 // FullStateDiagMatr and FullStateDiagMatrPower omitted for now
132
133 // apply SWAP
134 for (int r=0; r<numReps; r++) {
135 auto [ctrls,states,targs] = getRandomVariNumCtrlsStatesTargs(psi.numQubits, 2,2);
136 applyMultiStateControlledSwap(psi, ctrls.data(), states.data(), ctrls.size(), targs[0], targs[1]);
137 applyMultiStateControlledSwap(rho, ctrls.data(), states.data(), ctrls.size(), targs[0], targs[1]);
138 }
139
140 // apply PauliStr
141 for (int r=0; r<numReps; r++) {
142 auto [ctrls,states,targs] = getRandomVariNumCtrlsStatesTargs(psi.numQubits, 1,maxNumPauliStrTargs);
143 PauliStr str = getRandomPauliStr(targs);
144 applyMultiStateControlledPauliStr(psi, ctrls.data(), states.data(), ctrls.size(), str);
145 applyMultiStateControlledPauliStr(rho, ctrls.data(), states.data(), ctrls.size(), str);
146 }
147
148 // apply Pauli gadget
149 for (int r=0; r<numReps; r++) {
150 auto [ctrls,states,targs] = getRandomVariNumCtrlsStatesTargs(psi.numQubits, 1,maxNumPauliGadTargs);
151 PauliStr str = getRandomPauliStr(targs);
152 qreal phi = getRandomReal(-2 * 3.14, 2 * 3.14);
153 applyMultiStateControlledPauliGadget(psi, ctrls.data(), states.data(), ctrls.size(), str, phi);
154 applyMultiStateControlledPauliGadget(rho, ctrls.data(), states.data(), ctrls.size(), str, phi);
155 }
156
157 // apply phase gadet
158 for (int r=0; r<numReps; r++) {
159 auto [ctrls,states,targs] = getRandomVariNumCtrlsStatesTargs(psi.numQubits, 1,maxNumPhaseGadTargs);
160 qreal phi = getRandomReal(-2 * 3.14, 2 * 3.14);
161 applyMultiStateControlledPhaseGadget(psi, ctrls.data(), states.data(), ctrls.size(), targs.data(), targs.size(), phi);
162 applyMultiStateControlledPhaseGadget(rho, ctrls.data(), states.data(), ctrls.size(), targs.data(), targs.size(), phi);
163 }
164
165 // confirm purity and normalisation was maintained
166 REQUIRE_THAT( calcPurity(rho), WithinAbs(1, eps) );
167 REQUIRE_THAT( calcPurity(psi), WithinAbs(1, eps) );
168 REQUIRE_THAT( calcTotalProb(rho), WithinAbs(1, eps) );
169 REQUIRE_THAT( calcTotalProb(psi), WithinAbs(1, eps) );
170
171 // confirm states agree
172 REQUIRE_THAT( calcFidelity(rho, psi), WithinAbs(1, eps) );
173
174 // confirm expectation values agree
175 for (int r=0; r<numReps; r++) {
176 PauliStr str = getRandomPauliStr(rho.numQubits);
177 qreal psiExpec = calcExpecPauliStr(psi, str);
178 qreal rhoExpec = calcExpecPauliStr(rho, str);
179 REQUIRE_THAT( psiExpec, WithinAbs(rhoExpec, eps) );
180 }
181
182 // confirm all one-qubit probabilities agree
183 for (int q=0; q<psi.numQubits; q++) {
184 qreal psiProb = calcProbOfQubitOutcome(psi, q, 0);
185 qreal rhoProb = calcProbOfQubitOutcome(rho, q, 0);
186 REQUIRE_THAT( psiProb, WithinAbs(rhoProb, eps) );
187 }
188
189 // confirm some multi-qubit (1-6) probabilities agree
190 for (int r=0; r<numReps; r++) {
191 int numTargets = getRandomInt(1, 6+1);
192 vector<int> targets = getRandomSubRange(0, psi.numQubits, numTargets);
193 vector<int> outcomes = getRandomOutcomes(targets.size());
194 qreal psiProb = calcProbOfMultiQubitOutcome(psi, targets.data(), outcomes.data(), numTargets);
195 qreal rhoProb = calcProbOfMultiQubitOutcome(rho, targets.data(), outcomes.data(), numTargets);
196 REQUIRE_THAT( psiProb, WithinAbs(rhoProb, eps) );
197 }
198
199 // confirm some basis state probabilities agree
200 for (int r=0; r<numReps; r++) {
201 int ind = getRandomInt(0, psi.numAmps);
202 qreal psiProb = calcProbOfBasisState(psi, ind);
203 qreal rhoProb = calcProbOfBasisState(rho, ind);
204 REQUIRE_THAT( psiProb, WithinAbs(rhoProb, eps) );
205 }
206
207 // mix all canonical channels...
208 for (int r=0; r<numReps; r++) {
209 vector<int> targs = getRandomSubRange(0, rho.numQubits, 2);
210 mixDephasing(rho, targs[0], getRandomReal(0,1/2.));
211 mixTwoQubitDephasing(rho, targs[0], targs[1], getRandomReal(0,3/4.));
212 mixDepolarising(rho, targs[0], getRandomReal(0,3/4.));
213 mixTwoQubitDepolarising(rho, targs[0], targs[1], getRandomReal(0,15/16.));
214 mixDamping(rho, targs[0], getRandomReal(0,1));
215 mixPaulis(rho, targs[0], getRandomReal(0,1/10.), getRandomReal(0,1/20.), getRandomReal(0,1/30.));
216 }
217
218 // to confirm total probability unaffected
219 REQUIRE_THAT( calcTotalProb(rho), WithinAbs(1, eps) );
220
221 // confirm damping restores outcome probability
222 int qubitInd = 0;
223 int qubitOutcome = 0;
224 qreal dampProb = 1;
225 qreal outcomeProb = 1;
226 mixDamping(rho, qubitInd, dampProb);
227 REQUIRE_THAT( calcProbOfQubitOutcome(rho, qubitInd, qubitOutcome), WithinAbs(outcomeProb, eps) );
228}
229
230
231/// @ingroup integrationtests
232TEST_CASE( "density evolution", TEST_TAG ) {
233
234 auto deployments = getSupportedDeployments();
235
236 // try all combination of statevec and density-matrix deploments
237 for (auto [rhoDeploy, rhoMPI, rhoGPU, rhoOMP] : deployments) {
238 for (auto [psiDeploy, psiMPI, psiGPU, psiOMP] : deployments) {
239
240 // some combinations are illegal
241 if (psiMPI && !rhoMPI)
242 continue;
243
244 // Qureg size determined by slowest deployment
245 int numQubits = 6;
246 if (rhoMPI && rhoMPI) numQubits = 12;
247 if (rhoOMP && rhoOMP) numQubits = 12;
248 if (rhoGPU && psiGPU) numQubits = 14;
249
250 auto label = (
251 "rho = " + rhoDeploy + ", " +
252 "psi = " + psiDeploy + ", " +
253 "qubits = " + std::to_string(numQubits));
254
255 DYNAMIC_SECTION( label ) {
256
257 Qureg psi = createCustomQureg(numQubits, 0, psiMPI,psiGPU,psiOMP);
258 Qureg rho = createCustomQureg(numQubits, 1, rhoMPI,rhoGPU,rhoOMP);
259
260 testDensityMatrixEvolution(psi, rho);
261
262 destroyQureg(psi);
263 destroyQureg(rho);
264 }
265 }
266 }
267}
qreal calcFidelity(Qureg qureg, Qureg other)
qreal calcExpecPauliStr(Qureg qureg, PauliStr str)
qreal calcProbOfQubitOutcome(Qureg qureg, int qubit, int outcome)
qreal calcProbOfBasisState(Qureg qureg, qindex index)
qreal calcProbOfMultiQubitOutcome(Qureg qureg, int *qubits, int *outcomes, int numQubits)
qreal calcPurity(Qureg qureg)
qreal calcTotalProb(Qureg qureg)
void mixDepolarising(Qureg qureg, int qubit, qreal prob)
void mixTwoQubitDephasing(Qureg qureg, int qubit1, int qubit2, qreal prob)
void mixDamping(Qureg qureg, int qubit, qreal prob)
void mixPaulis(Qureg qureg, int qubit, qreal probX, qreal probY, qreal probZ)
void mixTwoQubitDepolarising(Qureg qureg, int qubit1, int qubit2, qreal prob)
void mixDephasing(Qureg qureg, int qubit, qreal prob)
void initRandomPureState(Qureg qureg)
void initPureState(Qureg qureg, Qureg pure)
CompMatr createCompMatr(int numQubits)
Definition matrices.cpp:211
DiagMatr createDiagMatr(int numQubits)
Definition matrices.cpp:246
void destroyDiagMatr(DiagMatr matrix)
Definition matrices.cpp:403
void destroyCompMatr(CompMatr matrix)
Definition matrices.cpp:402
static CompMatr2 getCompMatr2(qcomp **in)
Definition matrices.h:327
static CompMatr1 getCompMatr1(qcomp **in)
Definition matrices.h:311
static DiagMatr2 getDiagMatr2(qcomp *in)
Definition matrices.h:359
static DiagMatr1 getDiagMatr1(qcomp *in)
Definition matrices.h:345
void setDiagMatr(DiagMatr out, qcomp *in)
Definition matrices.cpp:435
void setCompMatr(CompMatr matr, qcomp **vals)
Definition matrices.cpp:427
void applyMultiStateControlledCompMatr1(Qureg qureg, int *controls, int *states, int numControls, int target, CompMatr1 matrix)
void applyMultiStateControlledCompMatr2(Qureg qureg, int *controls, int *states, int numControls, int target1, int target2, CompMatr2 matr)
void applyMultiStateControlledCompMatr(Qureg qureg, int *controls, int *states, int numControls, int *targets, int numTargets, CompMatr matr)
void applyMultiStateControlledDiagMatr1(Qureg qureg, int *controls, int *states, int numControls, int target, DiagMatr1 matr)
void applyMultiStateControlledDiagMatr2(Qureg qureg, int *controls, int *states, int numControls, int target1, int target2, DiagMatr2 matr)
void applyMultiStateControlledDiagMatrPower(Qureg qureg, int *controls, int *states, int numControls, int *targets, int numTargets, DiagMatr matrix, qcomp exponent)
void applyMultiStateControlledDiagMatr(Qureg qureg, int *controls, int *states, int numControls, int *targets, int numTargets, DiagMatr matrix)
void applyMultiStateControlledPauliGadget(Qureg qureg, int *controls, int *states, int numControls, PauliStr str, qreal angle)
void applyMultiStateControlledPauliStr(Qureg qureg, int *controls, int *states, int numControls, PauliStr str)
void applyMultiStateControlledPhaseGadget(Qureg qureg, int *controls, int *states, int numControls, int *targets, int numTargets, qreal angle)
void applyMultiStateControlledSwap(Qureg qureg, int *controls, int *states, int numControls, int qubit1, int qubit2)
Qureg createCustomQureg(int numQubits, int isDensMatr, int useDistrib, int useGpuAccel, int useMultithread)
Definition qureg.cpp:273
void destroyQureg(Qureg qureg)
Definition qureg.cpp:330
qmatrix getRandomUnitary(int numQb)
Definition random.cpp:348
qreal getRandomReal(qreal min, qreal maxExcl)
Definition random.cpp:63
int getRandomInt(int min, int maxExcl)
Definition random.cpp:90
TEST_CASE("calcExpecPauliStr", TEST_CATEGORY)
Definition qureg.h:49