Qrack  10.0
General classical-emulating-quantum development framework
qengine_cpu.hpp
Go to the documentation of this file.
1 //
3 // (C) Daniel Strano and the Qrack contributors 2017-2023. All rights reserved.
4 //
5 // This is a multithreaded, universal quantum register simulation, allowing
6 // (nonphysical) register cloning and direct measurement of probability and
7 // phase, to leverage what advantages classical emulation of qubits can have.
8 //
9 // Licensed under the GNU Lesser General Public License V3.
10 // See LICENSE.md in the project root or https://www.gnu.org/licenses/lgpl-3.0.en.html
11 // for details.
12 
13 #pragma once
14 
15 #include "qengine.hpp"
16 #include "statevector.hpp"
17 
18 #if ENABLE_QUNIT_CPU_PARALLEL && ENABLE_PTHREAD
19 #include "common/dispatchqueue.hpp"
20 #endif
21 
22 namespace Qrack {
23 
24 class QEngineCPU;
25 typedef std::shared_ptr<QEngineCPU> QEngineCPUPtr;
26 
27 template <class BidirectionalIterator>
28 void reverse(BidirectionalIterator first, BidirectionalIterator last, const bitCapInt& stride);
29 template <class BidirectionalIterator>
30 void rotate(
31  BidirectionalIterator first, BidirectionalIterator middle, BidirectionalIterator last, const bitCapInt& stride);
32 
36 class QEngineCPU : public QEngine {
37 protected:
39 #if ENABLE_QUNIT_CPU_PARALLEL && ENABLE_PTHREAD
40  DispatchQueue dispatchQueue;
41 #endif
42  double logFidelity;
43  bool isSparse;
44 
45  using QEngine::Copy;
46  void Copy(QInterfacePtr orig) { Copy(std::dynamic_pointer_cast<QEngineCPU>(orig)); }
47  void Copy(QEngineCPUPtr orig)
48  {
49  QEngine::Copy(std::dynamic_pointer_cast<QEngine>(orig));
50  stateVec = orig->stateVec;
51  logFidelity = orig->logFidelity;
52  isSparse = orig->isSparse;
53  }
54 
55  StateVectorSparsePtr CastStateVecSparse() { return std::dynamic_pointer_cast<StateVectorSparse>(stateVec); }
56 
57 public:
58  QEngineCPU(bitLenInt qBitCount, const bitCapInt& initState, qrack_rand_gen_ptr rgp = nullptr,
59  const complex& phaseFac = CMPLX_DEFAULT_ARG, bool doNorm = false, bool randomGlobalPhase = true,
60  bool ignored = false, int64_t ignored2 = -1, bool useHardwareRNG = true, bool useSparseStateVec = false,
61  real1_f norm_thresh = REAL1_EPSILON, std::vector<int64_t> ignored3 = {}, bitLenInt ignored4 = 0U,
62  real1_f ignored5 = _qrack_qunit_sep_thresh);
63 
64  ~QEngineCPU() { Dump(); }
65 
66  void Finish()
67  {
68 #if ENABLE_QUNIT_CPU_PARALLEL && ENABLE_PTHREAD
69  dispatchQueue.finish();
70 #endif
71  };
72 
73  bool isFinished()
74  {
75 #if ENABLE_QUNIT_CPU_PARALLEL && ENABLE_PTHREAD
76  return dispatchQueue.isFinished();
77 #else
78  return true;
79 #endif
80  }
81 
82  void Dump()
83  {
84 #if ENABLE_QUNIT_CPU_PARALLEL && ENABLE_PTHREAD
85  dispatchQueue.dump();
86 #endif
87  }
88 
90  {
91  if (!stateVec) {
92  return ZERO_R1_F;
93  }
94 
96  }
97 
99  {
100  Dump();
101  FreeStateVec();
103  logFidelity = 0.0;
104  }
105 
107  {
108  if (isSparse) {
109  logFidelity += (double)log(CastStateVecSparse()->truncate_to_size(QRACK_SPARSE_MAX_KEYS));
110  }
111  }
112 
113  void SparseRenorm(real1_f totFidelityLoss)
114  {
116  if ((totFidelityLoss > REAL1_EPSILON) && (sv->size() <= QRACK_SPARSE_MAX_KEYS)) {
117  const double f = 1.0 - (double)totFidelityLoss;
118  logFidelity += (double)log(f);
119  sv->mult((real1_f)std::sqrt(1.0 / f));
120  } else {
121  logFidelity += (double)log(sv->truncate_to_size(QRACK_SPARSE_MAX_KEYS));
122  }
123  }
124 
126  {
127  if (isSparse) {
128  return CastStateVecSparse()->size();
129  }
130 
131  return maxQPower;
132  }
133 
134  double GetUnitaryFidelity() { return exp(logFidelity); }
136 
137  bool IsZeroAmplitude() { return !stateVec; }
138  void GetAmplitudePage(complex* pagePtr, bitCapIntOcl offset, bitCapIntOcl length);
139  void SetAmplitudePage(const complex* pagePtr, bitCapIntOcl offset, bitCapIntOcl length);
140  void SetAmplitudePage(
141  QEnginePtr pageEnginePtr, bitCapIntOcl srcOffset, bitCapIntOcl dstOffset, bitCapIntOcl length);
142  void ShuffleBuffers(QEnginePtr engine);
143  void CopyStateVec(QEnginePtr src);
144 
146 
147  void QueueSetDoNormalize(bool doNorm)
148  {
149  Dispatch(1U, [this, doNorm] { doNormalize = doNorm; });
150  }
151  void QueueSetRunningNorm(real1_f runningNrm)
152  {
153  Dispatch(1U, [this, runningNrm] { runningNorm = runningNrm; });
154  }
155 
156  void SetQuantumState(const complex* inputState);
157  void GetQuantumState(complex* outputState);
158  void GetProbs(real1* outputProbs);
159  complex GetAmplitude(const bitCapInt& perm);
160  void SetAmplitude(const bitCapInt& perm, const complex& amp);
161 
162  using QEngine::Compose;
164  bitLenInt Compose(QInterfacePtr toCopy) { return Compose(std::dynamic_pointer_cast<QEngineCPU>(toCopy)); }
165  std::map<QInterfacePtr, bitLenInt> Compose(std::vector<QInterfacePtr> toCopy);
166  bitLenInt Compose(QEngineCPUPtr toCopy, bitLenInt start);
168  {
169  return Compose(std::dynamic_pointer_cast<QEngineCPU>(toCopy), start);
170  }
171 
172  using QEngine::Decompose;
173  void Decompose(bitLenInt start, QInterfacePtr dest);
174 
175  void Dispose(bitLenInt start, bitLenInt length);
176  void Dispose(bitLenInt start, bitLenInt length, const bitCapInt& disposedPerm);
177 
178  using QEngine::Allocate;
179  bitLenInt Allocate(bitLenInt start, bitLenInt length);
180 
183  void XMask(const bitCapInt& mask);
184  void PhaseParity(real1_f radians, const bitCapInt& mask);
185  void PhaseRootNMask(bitLenInt n, const bitCapInt& mask);
186 
193  void ROL(bitLenInt shift, bitLenInt start, bitLenInt length);
194 #if ENABLE_ALU
195  void INC(const bitCapInt& toAdd, bitLenInt start, bitLenInt length);
196  void CINC(const bitCapInt& toAdd, bitLenInt inOutStart, bitLenInt length, const std::vector<bitLenInt>& controls);
197  void INCS(const bitCapInt& toAdd, bitLenInt start, bitLenInt length, bitLenInt overflowIndex);
198 #if ENABLE_BCD
199  void INCBCD(const bitCapInt& toAdd, bitLenInt start, bitLenInt length);
200 #endif
201  void MUL(const bitCapInt& toMul, bitLenInt inOutStart, bitLenInt carryStart, bitLenInt length);
202  void DIV(const bitCapInt& toDiv, bitLenInt inOutStart, bitLenInt carryStart, bitLenInt length);
203  void MULModNOut(
204  const bitCapInt& toMul, const bitCapInt& modN, bitLenInt inStart, bitLenInt outStart, bitLenInt length);
205  void IMULModNOut(
206  const bitCapInt& toMul, const bitCapInt& modN, bitLenInt inStart, bitLenInt outStart, bitLenInt length);
207  void POWModNOut(
208  const bitCapInt& base, const bitCapInt& modN, bitLenInt inStart, bitLenInt outStart, bitLenInt length);
209  void CMUL(const bitCapInt& toMul, bitLenInt inOutStart, bitLenInt carryStart, bitLenInt length,
210  const std::vector<bitLenInt>& controls);
211  void CDIV(const bitCapInt& toDiv, bitLenInt inOutStart, bitLenInt carryStart, bitLenInt length,
212  const std::vector<bitLenInt>& controls);
213  void CMULModNOut(const bitCapInt& toMul, const bitCapInt& modN, bitLenInt inStart, bitLenInt outStart,
214  bitLenInt length, const std::vector<bitLenInt>& controls);
215  void CIMULModNOut(const bitCapInt& toMul, const bitCapInt& modN, bitLenInt inStart, bitLenInt outStart,
216  bitLenInt length, const std::vector<bitLenInt>& controls);
217  void CPOWModNOut(const bitCapInt& base, const bitCapInt& modN, bitLenInt inStart, bitLenInt outStart,
218  bitLenInt length, const std::vector<bitLenInt>& controls);
219  void FullAdd(bitLenInt inputBit1, bitLenInt inputBit2, bitLenInt carryInSumOut, bitLenInt carryOut);
220  void IFullAdd(bitLenInt inputBit1, bitLenInt inputBit2, bitLenInt carryInSumOut, bitLenInt carryOut);
221  bitCapInt IndexedLDA(bitLenInt indexStart, bitLenInt indexLength, bitLenInt valueStart, bitLenInt valueLength,
222  const unsigned char* values, bool resetValue = true);
223  bitCapInt IndexedADC(bitLenInt indexStart, bitLenInt indexLength, bitLenInt valueStart, bitLenInt valueLength,
224  bitLenInt carryIndex, const unsigned char* values);
225  bitCapInt IndexedSBC(bitLenInt indexStart, bitLenInt indexLength, bitLenInt valueStart, bitLenInt valueLength,
226  bitLenInt carryIndex, const unsigned char* values);
227  void Hash(bitLenInt start, bitLenInt length, const unsigned char* values);
228  void CPhaseFlipIfLess(const bitCapInt& greaterPerm, bitLenInt start, bitLenInt length, bitLenInt flagIndex);
229  void PhaseFlipIfLess(const bitCapInt& greaterPerm, bitLenInt start, bitLenInt length);
230 
232 #endif
233 
240  void SetPermutation(const bitCapInt& perm, const complex& phaseFac = CMPLX_DEFAULT_ARG);
242  void UniformlyControlledSingleBit(const std::vector<bitLenInt>& controls, bitLenInt qubitIndex,
243  const complex* mtrxs, const std::vector<bitCapInt>& mtrxSkipPowers, const bitCapInt& mtrxSkipValueMask);
244  void UniformParityRZ(const bitCapInt& mask, real1_f angle);
245  void CUniformParityRZ(const std::vector<bitLenInt>& controls, const bitCapInt& mask, real1_f angle);
246 
255  real1_f Prob(bitLenInt qubitIndex);
256  real1_f CtrlOrAntiProb(bool controlState, bitLenInt control, bitLenInt target);
257  real1_f ProbReg(bitLenInt start, bitLenInt length, const bitCapInt& permutation);
258  real1_f ProbMask(const bitCapInt& mask, const bitCapInt& permutation);
259  real1_f ProbParity(const bitCapInt& mask);
261  bitCapInt MAll();
262  bool ForceMParity(const bitCapInt& mask, bool result, bool doForce = true);
263  void NormalizeState(
264  real1_f nrm = REAL1_DEFAULT_ARG, real1_f norm_thresh = REAL1_DEFAULT_ARG, real1_f phaseArg = ZERO_R1_F);
265  real1_f SumSqrDiff(QInterfacePtr toCompare) { return SumSqrDiff(std::dynamic_pointer_cast<QEngineCPU>(toCompare)); }
266  real1_f SumSqrDiff(QEngineCPUPtr toCompare);
269 
272 protected:
273  real1_f GetExpectation(bitLenInt valueStart, bitLenInt valueLength);
274 
276  {
277  if (isSparse) {
278  return std::make_shared<StateVectorSparse>(elemCount);
279  } else {
280  return std::make_shared<StateVectorArray>(elemCount);
281  }
282  }
284  void FreeStateVec(complex* sv = nullptr) { stateVec = nullptr; }
285 
286  void Dispatch(bitCapIntOcl workItemCount, DispatchFn fn)
287  {
288 #if ENABLE_QUNIT_CPU_PARALLEL && ENABLE_PTHREAD
289  if ((workItemCount >= pow2Ocl(GetPreferredConcurrencyPower())) && (workItemCount < GetStride())) {
290  dispatchQueue.dispatch(fn);
291  } else {
292  Finish();
293  fn();
294  }
295 #else
296  fn();
297 #endif
298  }
299 
300  void DecomposeDispose(bitLenInt start, bitLenInt length, QEngineCPUPtr dest);
301  void Apply2x2(bitCapIntOcl offset1, bitCapIntOcl offset2, const complex* mtrx, bitLenInt bitCount,
302  const bitCapIntOcl* qPowersSorted, bool doCalcNorm, real1_f norm_thresh = REAL1_DEFAULT_ARG);
303  void UpdateRunningNorm(real1_f norm_thresh = REAL1_DEFAULT_ARG);
304  using QEngine::ApplyM;
305  void ApplyM(const bitCapInt& mask, const bitCapInt& result, const complex& nrm);
306 
307 #if ENABLE_ALU
308  void INCDECC(const bitCapInt& toMod, bitLenInt inOutStart, bitLenInt length, bitLenInt carryIndex);
309  void INCDECSC(const bitCapInt& toMod, bitLenInt inOutStart, bitLenInt length, bitLenInt carryIndex);
310  void INCDECSC(
311  const bitCapInt& toMod, bitLenInt inOutStart, bitLenInt length, bitLenInt overflowIndex, bitLenInt carryIndex);
312 #if ENABLE_BCD
313  void INCDECBCDC(const bitCapInt& toMod, bitLenInt inOutStart, bitLenInt length, bitLenInt carryIndex);
314 #endif
315 
316  typedef std::function<bitCapIntOcl(const bitCapIntOcl&, const bitCapIntOcl&)> IOFn;
317  void MULDIV(const IOFn& inFn, const IOFn& outFn, const bitCapInt& toMul, const bitLenInt& inOutStart,
318  const bitLenInt& carryStart, const bitLenInt& length);
319  void CMULDIV(const IOFn& inFn, const IOFn& outFn, const bitCapInt& toMul, const bitLenInt& inOutStart,
320  const bitLenInt& carryStart, const bitLenInt& length, const std::vector<bitLenInt>& controls);
321 
322  typedef std::function<bitCapIntOcl(const bitCapIntOcl&)> MFn;
323  void ModNOut(const MFn& kernelFn, const bitCapInt& modN, const bitLenInt& inStart, const bitLenInt& outStart,
324  const bitLenInt& length, const bool& inverse = false);
325  void CModNOut(const MFn& kernelFn, const bitCapInt& modN, const bitLenInt& inStart, const bitLenInt& outStart,
326  const bitLenInt& length, const std::vector<bitLenInt>& controls, const bool& inverse = false);
327 #endif
328 };
329 } // namespace Qrack
Definition: dispatchqueue.hpp:33
void finish()
Definition: dispatchqueue.cpp:40
void dispatch(const DispatchFn &op)
Definition: dispatchqueue.cpp:66
void dump()
Definition: dispatchqueue.cpp:51
bool isFinished()
Definition: dispatchqueue.hpp:51
bitLenInt GetPreferredConcurrencyPower()
Definition: parallel_for.hpp:43
bitCapIntOcl GetStride()
Definition: parallel_for.hpp:42
General purpose QEngineCPU implementation.
Definition: qengine_cpu.hpp:36
void PhaseRootNMask(bitLenInt n, const bitCapInt &mask)
Masked PhaseRootN gate.
Definition: state.cpp:764
bool isFinished()
Returns "false" if asynchronous work is still running, and "true" if all previously dispatched asynch...
Definition: qengine_cpu.hpp:73
bitLenInt Compose(QEngineCPUPtr toCopy)
Combine (a copy of) another QEngineCPU with this one, after the last bit index of this one.
Definition: state.cpp:977
void DecomposeDispose(bitLenInt start, bitLenInt length, QEngineCPUPtr dest)
Minimally decompose a set of contigious bits from the separable unit.
Definition: state.cpp:1263
void Dispatch(bitCapIntOcl workItemCount, DispatchFn fn)
Definition: qengine_cpu.hpp:286
void QueueSetRunningNorm(real1_f runningNrm)
Add an operation to the (OpenCL) queue, to set the value of runningNorm, which is the normalization c...
Definition: qengine_cpu.hpp:151
complex GetAmplitude(const bitCapInt &perm)
Get the representational amplitude of a full permutation.
Definition: state.cpp:184
StateVectorPtr stateVec
Definition: qengine_cpu.hpp:38
void SparseRenorm(real1_f totFidelityLoss)
Definition: qengine_cpu.hpp:113
void Copy(QInterfacePtr orig)
Definition: qengine_cpu.hpp:46
void ZeroAmplitudes()
Set all amplitudes to 0, and optionally temporarily deallocate state vector RAM.
Definition: qengine_cpu.hpp:98
void SetAmplitude(const bitCapInt &perm, const complex &amp)
Sets the representational amplitude of a full permutation.
Definition: state.cpp:200
void TruncateBySize()
Definition: qengine_cpu.hpp:106
QEnginePtr CloneEmpty()
Clone this QEngine's settings, with a zeroed state vector.
Definition: utility.cpp:35
virtual void ApplyM(const bitCapInt &qPower, bool result, const complex &nrm)
Definition: qengine.hpp:167
bool isSparse
Definition: qengine_cpu.hpp:43
void CopyStateVec(QEnginePtr src)
Exactly copy the state vector of a different QEngine instance.
Definition: state.cpp:162
void FreeStateVec(complex *sv=nullptr)
Definition: qengine_cpu.hpp:284
void UpdateRunningNorm(real1_f norm_thresh=REAL1_DEFAULT_ARG)
Force a calculation of the norm of the state vector, in order to make it unit length before the next ...
Definition: state.cpp:1949
void Dispose(bitLenInt start, bitLenInt length)
Minimally decompose a set of contiguous bits from the separably composed unit, and discard the separa...
Definition: state.cpp:1415
bitLenInt Compose(QInterfacePtr toCopy, bitLenInt start)
Compose() a QInterface peer, inserting its qubit into index order at start index.
Definition: qengine_cpu.hpp:167
void GetQuantumState(complex *outputState)
Get pure quantum state, in unsigned int permutation basis.
Definition: state.cpp:268
bitLenInt Compose(QInterfacePtr toCopy)
Combine another QInterface with this one, after the last bit index of this one.
Definition: qengine_cpu.hpp:164
void XMask(const bitCapInt &mask)
Masked X gate.
Definition: state.cpp:673
std::function< bitCapIntOcl(const bitCapIntOcl &, const bitCapIntOcl &)> IOFn
Definition: qengine_cpu.hpp:316
void Apply2x2(bitCapIntOcl offset1, bitCapIntOcl offset2, const complex *mtrx, bitLenInt bitCount, const bitCapIntOcl *qPowersSorted, bool doCalcNorm, real1_f norm_thresh=REAL1_DEFAULT_ARG)
Definition: state.cpp:530
StateVectorPtr AllocStateVec(bitCapIntOcl elemCount)
Definition: qengine_cpu.hpp:275
void GetProbs(real1 *outputProbs)
Get all probabilities, in unsigned int permutation basis.
Definition: state.cpp:284
QEngineCPU(bitLenInt qBitCount, const bitCapInt &initState, qrack_rand_gen_ptr rgp=nullptr, const complex &phaseFac=CMPLX_DEFAULT_ARG, bool doNorm=false, bool randomGlobalPhase=true, bool ignored=false, int64_t ignored2=-1, bool useHardwareRNG=true, bool useSparseStateVec=false, real1_f norm_thresh=REAL1_EPSILON, std::vector< int64_t > ignored3={}, bitLenInt ignored4=0U, real1_f ignored5=_qrack_qunit_sep_thresh)
Initialize a coherent unit with qBitCount number of bits, to initState unsigned integer permutation s...
Definition: state.cpp:35
void ResetStateVec(StateVectorPtr sv)
Definition: qengine_cpu.hpp:283
void ShuffleBuffers(QEnginePtr engine)
Swap the high half of this engine with the low half of another.
Definition: state.cpp:131
real1_f GetExpectation(bitLenInt valueStart, bitLenInt valueLength)
Definition: utility.cpp:70
void CModNOut(const MFn &kernelFn, const bitCapInt &modN, const bitLenInt &inStart, const bitLenInt &outStart, const bitLenInt &length, const std::vector< bitLenInt > &controls, const bool &inverse=false)
Definition: arithmetic.cpp:632
std::function< bitCapIntOcl(const bitCapIntOcl &)> MFn
Definition: qengine_cpu.hpp:322
void Copy(QEngineCPUPtr orig)
Definition: qengine_cpu.hpp:47
~QEngineCPU()
Definition: qengine_cpu.hpp:64
void PhaseParity(real1_f radians, const bitCapInt &mask)
Parity phase gate.
Definition: state.cpp:717
void SetAmplitudePage(const complex *pagePtr, bitCapIntOcl offset, bitCapIntOcl length)
Copy a "page" of amplitudes from pagePtr into this QEngine's internal state.
Definition: state.cpp:77
void MULDIV(const IOFn &inFn, const IOFn &outFn, const bitCapInt &toMul, const bitLenInt &inOutStart, const bitLenInt &carryStart, const bitLenInt &length)
Definition: arithmetic.cpp:384
virtual QInterfacePtr Decompose(bitLenInt start, bitLenInt length)
Definition: qengine.hpp:293
void INCDECC(const bitCapInt &toMod, bitLenInt inOutStart, bitLenInt length, bitLenInt carryIndex)
Add integer (without sign, with carry)
Definition: arithmetic.cpp:165
void CMULDIV(const IOFn &inFn, const IOFn &outFn, const bitCapInt &toMul, const bitLenInt &inOutStart, const bitLenInt &carryStart, const bitLenInt &length, const std::vector< bitLenInt > &controls)
Definition: arithmetic.cpp:450
bool IsZeroAmplitude()
Returns "true" only if amplitudes are all totally 0.
Definition: qengine_cpu.hpp:137
void QueueSetDoNormalize(bool doNorm)
Add an operation to the (OpenCL) queue, to set the value of doNormalize, which controls whether to au...
Definition: qengine_cpu.hpp:147
double logFidelity
Definition: qengine_cpu.hpp:42
void GetAmplitudePage(complex *pagePtr, bitCapIntOcl offset, bitCapIntOcl length)
Copy a "page" of amplitudes from this QEngine's internal state, into pagePtr.
Definition: state.cpp:63
void SetQuantumState(const complex *inputState)
Set arbitrary pure quantum state, in unsigned int permutation basis.
Definition: state.cpp:254
bitLenInt Allocate(bitLenInt start, bitLenInt length)
Allocate new "length" count of |0> state qubits at specified qubit index start position.
Definition: utility.cpp:54
real1_f FirstNonzeroPhase()
Get phase of lowest permutation nonzero amplitude.
Definition: qengine_cpu.hpp:89
StateVectorSparsePtr CastStateVecSparse()
Definition: qengine_cpu.hpp:55
void Dump()
If asynchronous work is still running, let the simulator know that it can be aborted.
Definition: qengine_cpu.hpp:82
double GetUnitaryFidelity()
When "Schmidt-decomposition rounding parameter" ("SDRP") is being used, starting from initial 1....
Definition: qengine_cpu.hpp:134
void INCDECBCDC(const bitCapInt &toMod, bitLenInt inOutStart, bitLenInt length, bitLenInt carryIndex)
Add BCD integer (without sign, with carry)
Definition: arithmetic.cpp:815
void ResetUnitaryFidelity()
Reset the internal fidelity calculation tracker to 1.0.
Definition: qengine_cpu.hpp:135
void INCDECSC(const bitCapInt &toMod, bitLenInt inOutStart, bitLenInt length, bitLenInt carryIndex)
Common driver method behind INCSC and DECSC (without overflow flag)
Definition: arithmetic.cpp:274
bitCapInt GetAmplitudeCount()
Count of amplitudes, which might be less than state vector if sparse or compressed.
Definition: qengine_cpu.hpp:125
void ModNOut(const MFn &kernelFn, const bitCapInt &modN, const bitLenInt &inStart, const bitLenInt &outStart, const bitLenInt &length, const bool &inverse=false)
Definition: arithmetic.cpp:557
void Finish()
If asynchronous work is still running, block until it finishes.
Definition: qengine_cpu.hpp:66
Abstract QEngine implementation, for all "Schroedinger method" engines.
Definition: qengine.hpp:31
virtual void Copy(QInterfacePtr orig)
Copy this QInterface.
Definition: qinterface.hpp:222
virtual void ApplyM(const bitCapInt &qPower, bool result, const complex &nrm)
Definition: qengine.hpp:167
real1 runningNorm
The value stored in runningNorm should always be the total probability implied by the norm of all amp...
Definition: qengine.hpp:39
virtual void Decompose(bitLenInt start, QInterfacePtr dest)=0
Minimally decompose a set of contiguous bits from the separably composed unit, into "destination".
bitCapInt maxQPower
Definition: qinterface.hpp:149
virtual bitLenInt Allocate(bitLenInt length)
Allocate new "length" count of |0> state qubits at end of qubit index position.
Definition: qinterface.hpp:477
virtual bitLenInt Compose(QInterfacePtr toCopy)
Combine another QInterface with this one, after the last bit index of this one.
Definition: qinterface.hpp:371
bool doNormalize
Definition: qinterface.hpp:143
Half-precision floating-point type.
Definition: half.hpp:2222
void IMULModNOut(const bitCapInt &toMul, const bitCapInt &modN, bitLenInt inStart, bitLenInt outStart, bitLenInt length)
Inverse of multiplication modulo N by integer, (out of place)
Definition: arithmetic.cpp:609
void DIV(const bitCapInt &toDiv, bitLenInt inOutStart, bitLenInt carryStart, bitLenInt length)
Divide by integer.
Definition: arithmetic.cpp:436
void IFullAdd(bitLenInt inputBit1, bitLenInt inputBit2, bitLenInt carryInSumOut, bitLenInt carryOut)
Inverse of FullAdd.
Definition: arithmetic.cpp:1358
void Hash(bitLenInt start, bitLenInt length, const unsigned char *values)
Transform a length of qubit register via lookup through a hash table.
Definition: arithmetic.cpp:1230
void CINC(const bitCapInt &toAdd, bitLenInt inOutStart, bitLenInt length, const std::vector< bitLenInt > &controls)
Add integer (without sign, with controls)
Definition: arithmetic.cpp:111
bitCapInt IndexedADC(bitLenInt indexStart, bitLenInt indexLength, bitLenInt valueStart, bitLenInt valueLength, bitLenInt carryIndex, const unsigned char *values)
Add based on an indexed load from classical memory.
Definition: arithmetic.cpp:984
void INC(const bitCapInt &toAdd, bitLenInt start, bitLenInt length)
Add integer (without sign)
Definition: arithmetic.cpp:68
void CPhaseFlipIfLess(const bitCapInt &greaterPerm, bitLenInt start, bitLenInt length, bitLenInt flagIndex)
The 6502 uses its carry flag also as a greater-than/less-than flag, for the CMP operation.
Definition: arithmetic.cpp:1443
void CMUL(const bitCapInt &toMul, bitLenInt inOutStart, bitLenInt carryStart, bitLenInt length, const std::vector< bitLenInt > &controls)
Controlled multiplication by integer.
Definition: arithmetic.cpp:515
void CDIV(const bitCapInt &toDiv, bitLenInt inOutStart, bitLenInt carryStart, bitLenInt length, const std::vector< bitLenInt > &controls)
Controlled division by power of integer.
Definition: arithmetic.cpp:537
void CIMULModNOut(const bitCapInt &toMul, const bitCapInt &modN, bitLenInt inStart, bitLenInt outStart, bitLenInt length, const std::vector< bitLenInt > &controls)
Inverse of controlled multiplication modulo N by integer, (out of place)
Definition: arithmetic.cpp:713
void CMULModNOut(const bitCapInt &toMul, const bitCapInt &modN, bitLenInt inStart, bitLenInt outStart, bitLenInt length, const std::vector< bitLenInt > &controls)
Controlled multiplication modulo N by integer, (out of place)
Definition: arithmetic.cpp:699
bitCapInt IndexedLDA(bitLenInt indexStart, bitLenInt indexLength, bitLenInt valueStart, bitLenInt valueLength, const unsigned char *values, bool resetValue=true)
Set 8 bit register bits based on read from classical memory.
Definition: arithmetic.cpp:910
void INCBCD(const bitCapInt &toAdd, bitLenInt start, bitLenInt length)
Add BCD integer (without sign)
Definition: arithmetic.cpp:739
void CPOWModNOut(const bitCapInt &base, const bitCapInt &modN, bitLenInt inStart, bitLenInt outStart, bitLenInt length, const std::vector< bitLenInt > &controls)
Controlled, raise a classical base to a quantum power, modulo N, (out of place)
Definition: arithmetic.cpp:725
void POWModNOut(const bitCapInt &base, const bitCapInt &modN, bitLenInt inStart, bitLenInt outStart, bitLenInt length)
Raise a classical base to a quantum power, modulo N, (out of place)
Definition: arithmetic.cpp:620
void PhaseFlipIfLess(const bitCapInt &greaterPerm, bitLenInt start, bitLenInt length)
This is an expedient for an adaptive Grover's search for a function's global minimum.
Definition: arithmetic.cpp:1468
void MUL(const bitCapInt &toMul, bitLenInt inOutStart, bitLenInt carryStart, bitLenInt length)
Multiply by integer.
Definition: arithmetic.cpp:420
void FullAdd(bitLenInt inputBit1, bitLenInt inputBit2, bitLenInt carryInSumOut, bitLenInt carryOut)
Quantum analog of classical "Full Adder" gate.
Definition: arithmetic.cpp:1274
void INCS(const bitCapInt &toAdd, bitLenInt start, bitLenInt length, bitLenInt overflowIndex)
Add an integer to the register, with sign and without carry.
Definition: arithmetic.cpp:217
void MULModNOut(const bitCapInt &toMul, const bitCapInt &modN, bitLenInt inStart, bitLenInt outStart, bitLenInt length)
Multiplication modulo N by integer, (out of place)
Definition: arithmetic.cpp:596
bitCapInt IndexedSBC(bitLenInt indexStart, bitLenInt indexLength, bitLenInt valueStart, bitLenInt valueLength, bitLenInt carryIndex, const unsigned char *values)
Subtract based on an indexed load from classical memory.
Definition: arithmetic.cpp:1105
void ROL(bitLenInt shift, bitLenInt start, bitLenInt length)
"Circular shift left" - shift bits left, and carry last bits.
Definition: arithmetic.cpp:23
virtual void UniformlyControlledSingleBit(const std::vector< bitLenInt > &controls, bitLenInt qubit, const complex *mtrxs)
Apply a "uniformly controlled" arbitrary single bit unitary transformation.
Definition: qinterface.hpp:634
virtual void U(bitLenInt target, real1_f theta, real1_f phi, real1_f lambda)
General unitary gate.
Definition: rotational.cpp:18
void UniformlyControlledSingleBit(const std::vector< bitLenInt > &controls, bitLenInt qubitIndex, const complex *mtrxs, const std::vector< bitCapInt > &mtrxSkipPowers, const bitCapInt &mtrxSkipValueMask)
Definition: state.cpp:802
void CUniformParityRZ(const std::vector< bitLenInt > &controls, const bitCapInt &mask, real1_f angle)
If the controls are set and the target qubit set parity is odd, this applies a phase factor of .
Definition: state.cpp:934
void SetPermutation(const bitCapInt &perm, const complex &phaseFac=CMPLX_DEFAULT_ARG)
Set to a specific permutation of all qubits.
Definition: state.cpp:225
void UniformParityRZ(const bitCapInt &mask, real1_f angle)
If the target qubit set parity is odd, this applies a phase factor of .
Definition: state.cpp:908
real1_f ProbParity(const bitCapInt &mask)
Overall probability of any odd permutation of the masked set of bits.
Definition: state.cpp:1664
virtual real1_f FirstNonzeroPhase()
Get phase of lowest permutation nonzero amplitude.
Definition: qinterface.hpp:3021
real1_f CtrlOrAntiProb(bool controlState, bitLenInt control, bitLenInt target)
PSEUDO-QUANTUM Direct measure of bit probability to be in |1> state, if control is in |0>/|1>,...
Definition: state.cpp:1535
real1_f Prob(bitLenInt qubitIndex)
PSEUDO-QUANTUM Direct measure of bit probability to be in |1> state.
Definition: state.cpp:1463
bool ForceMParity(const bitCapInt &mask, bool result, bool doForce=true)
Act as if is a measurement of parity of the masked set of qubits was applied, except force the (usual...
Definition: state.cpp:1763
QInterfacePtr Clone()
Clone this QInterface.
Definition: utility.cpp:17
void NormalizeState(real1_f nrm=REAL1_DEFAULT_ARG, real1_f norm_thresh=REAL1_DEFAULT_ARG, real1_f phaseArg=ZERO_R1_F)
Apply the normalization factor found by UpdateRunningNorm() or on the fly by a single bit gate.
Definition: state.cpp:1897
QInterfacePtr Copy()
Copy this QInterface.
Definition: utility.cpp:45
bitCapInt MAll()
Measure permutation state of all coherent bits.
Definition: state.cpp:1737
bitCapInt HighestProbAll()
Get highest probability permutation.
Definition: state.cpp:1706
real1_f ProbReg(bitLenInt start, bitLenInt length, const bitCapInt &permutation)
Direct measure of register permutation probability.
Definition: state.cpp:1589
real1_f ProbMask(const bitCapInt &mask, const bitCapInt &permutation)
Direct measure of masked permutation probability.
Definition: state.cpp:1625
real1_f SumSqrDiff(QInterfacePtr toCompare)
Calculates (1 - <\psi_e|\psi_c>) between states |\psi_c> and |\psi_e>.
Definition: qengine_cpu.hpp:265
GLOSSARY: bitLenInt - "bit-length integer" - unsigned integer ID of qubit position in register bitCap...
Definition: complex16x2simd.hpp:25
std::shared_ptr< QEngine > QEnginePtr
Definition: qrack_types.hpp:159
std::shared_ptr< QInterface > QInterfacePtr
Definition: qinterface.hpp:29
const real1_f _qrack_qunit_sep_thresh
Definition: qrack_functions.hpp:249
void rotate(BidirectionalIterator first, BidirectionalIterator middle, BidirectionalIterator last, const bitCapInt &stride)
const size_t QRACK_SPARSE_MAX_KEYS
Definition: qrack_functions.hpp:262
std::function< void(void)> DispatchFn
Definition: dispatchqueue.hpp:31
std::complex< real1 > complex
Definition: qrack_types.hpp:136
std::shared_ptr< StateVectorSparse > StateVectorSparsePtr
Definition: qrack_types.hpp:155
std::shared_ptr< QEngineCPU > QEngineCPUPtr
Definition: qengine_cpu.hpp:24
QRACK_CONST real1 REAL1_EPSILON
Definition: qrack_types.hpp:208
QRACK_CONST real1 ZERO_R1
Definition: qrack_types.hpp:191
float real1_f
Definition: qrack_types.hpp:103
QRACK_CONST complex CMPLX_DEFAULT_ARG
Definition: qrack_types.hpp:267
std::shared_ptr< StateVector > StateVectorPtr
Definition: qrack_types.hpp:151
void reverse(BidirectionalIterator first, BidirectionalIterator last, const bitCapInt &stride)
bitCapIntOcl pow2Ocl(const bitLenInt &p)
Definition: qrack_functions.hpp:144
uint32 sqrt(uint32 &r, int &exp)
Fixed point square root.
Definition: half.hpp:1654
half log(half arg)
Natural logarithm.
Definition: half.hpp:3318
half exp(half arg)
Exponential function.
Definition: half.hpp:3201
#define REAL1_DEFAULT_ARG
Definition: qrack_types.hpp:185
#define bitLenInt
Definition: qrack_types.hpp:42
#define ZERO_R1_F
Definition: qrack_types.hpp:168
#define qrack_rand_gen_ptr
Definition: qrack_types.hpp:164
#define bitCapInt
Definition: qrack_types.hpp:66
#define bitCapIntOcl
Definition: qrack_types.hpp:54