Qrack  7.0
General classical-emulating-quantum development framework
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
Qrack::QPager Class Reference

A "Qrack::QPager" splits a "Qrack::QEngine" implementation into equal-length "pages." This helps both optimization and distribution of a single coherent quantum register across multiple devices. More...

#include <qpager.hpp>

Inheritance diagram for Qrack::QPager:
Inheritance graph
[legend]
Collaboration diagram for Qrack::QPager:
Collaboration graph
[legend]

Public Member Functions

 QPager (std::vector< QInterfaceEngine > eng, bitLenInt qBitCount, bitCapInt initState=0, qrack_rand_gen_ptr rgp=nullptr, complex phaseFac=CMPLX_DEFAULT_ARG, bool doNorm=false, bool ignored=false, bool useHostMem=false, int deviceId=-1, bool useHardwareRNG=true, bool useSparseStateVec=false, real1_f norm_thresh=REAL1_EPSILON, std::vector< int > devList={}, bitLenInt qubitThreshold=0, real1_f separation_thresh=FP_NORM_EPSILON)
 
 QPager (bitLenInt qBitCount, bitCapInt initState=0, qrack_rand_gen_ptr rgp=nullptr, complex phaseFac=CMPLX_DEFAULT_ARG, bool doNorm=false, bool ignored=false, bool useHostMem=false, int deviceId=-1, bool useHardwareRNG=true, bool useSparseStateVec=false, real1_f norm_thresh=REAL1_EPSILON, std::vector< int > devList={}, bitLenInt qubitThreshold=0, real1_f separation_thresh=FP_NORM_EPSILON)
 
 QPager (QEnginePtr enginePtr, std::vector< QInterfaceEngine > eng, bitLenInt qBitCount, bitCapInt ignored=0, qrack_rand_gen_ptr rgp=nullptr, complex phaseFac=CMPLX_DEFAULT_ARG, bool doNorm=false, bool ignored2=false, bool useHostMem=false, int deviceId=-1, bool useHardwareRNG=true, bool useSparseStateVec=false, real1_f norm_thresh=REAL1_EPSILON, std::vector< int > devList={}, bitLenInt qubitThreshold=0, real1_f separation_thresh=FP_NORM_EPSILON)
 
virtual void SetConcurrency (uint32_t threadsPerEngine)
 Set the number of threads in parallel for loops, per component QEngine. More...
 
virtual QEnginePtr ReleaseEngine ()
 
virtual void LockEngine (QEnginePtr eng)
 
virtual void SetQuantumState (const complex *inputState)
 Set an arbitrary pure quantum state representation. More...
 
virtual void GetQuantumState (complex *outputState)
 Get the pure quantum state representation. More...
 
virtual void GetProbs (real1 *outputProbs)
 Get the pure quantum state representation. More...
 
virtual complex GetAmplitude (bitCapInt perm)
 Get the representational amplitude of a full permutation. More...
 
virtual void SetAmplitude (bitCapInt perm, complex amp)
 Sets the representational amplitude of a full permutation. More...
 
real1_f ProbAll (bitCapInt fullRegister)
 Direct measure of full permutation probability. More...
 
virtual void SetPermutation (bitCapInt perm, complex phaseFac=CMPLX_DEFAULT_ARG)
 Set to a specific permutation of all qubits. More...
 
virtual bitLenInt Compose (QPagerPtr toCopy)
 
virtual bitLenInt Compose (QInterfacePtr toCopy)
 Combine another QInterface with this one, after the last bit index of this one. More...
 
virtual void Decompose (bitLenInt start, QInterfacePtr dest)
 Minimally decompose a set of contiguous bits from the separably composed unit, into "destination". More...
 
virtual void Decompose (bitLenInt start, QPagerPtr dest)
 
virtual QInterfacePtr Decompose (bitLenInt start, bitLenInt length)
 Schmidt decompose a length of qubits. More...
 
virtual void Dispose (bitLenInt start, bitLenInt length)
 Minimally decompose a set of contiguous bits from the separably composed unit, and discard the separable bits from index "start" for "length.". More...
 
virtual void Dispose (bitLenInt start, bitLenInt length, bitCapInt disposedPerm)
 Dispose a a contiguous set of qubits that are already in a permutation eigenstate. More...
 
virtual void Mtrx (const complex *mtrx, bitLenInt qubitIndex)
 Apply an arbitrary single bit unitary transformation. More...
 
virtual void Phase (complex topLeft, complex bottomRight, bitLenInt qubitIndex)
 Apply a single bit transformation that only effects phase. More...
 
virtual void Invert (complex topRight, complex bottomLeft, bitLenInt qubitIndex)
 Apply a single bit transformation that reverses bit probability and might effect phase. More...
 
virtual void MCMtrx (const bitLenInt *controls, bitLenInt controlLen, const complex *mtrx, bitLenInt target)
 Apply an arbitrary single bit unitary transformation, with arbitrary control bits. More...
 
virtual void MACMtrx (const bitLenInt *controls, bitLenInt controlLen, const complex *mtrx, bitLenInt target)
 Apply an arbitrary single bit unitary transformation, with arbitrary (anti-)control bits. More...
 
virtual void UniformParityRZ (bitCapInt mask, real1_f angle)
 If the target qubit set parity is odd, this applies a phase factor of \(e^{i angle}\). More...
 
virtual void CUniformParityRZ (const bitLenInt *controls, bitLenInt controlLen, bitCapInt mask, real1_f angle)
 If the controls are set and the target qubit set parity is odd, this applies a phase factor of \(e^{i angle}\). More...
 
virtual void XMask (bitCapInt mask)
 Masked X gate. More...
 
virtual void ZMask (bitCapInt mask)
 Masked Z gate. More...
 
virtual void PhaseParity (real1_f radians, bitCapInt mask)
 Parity phase gate. More...
 
virtual bool ForceM (bitLenInt qubit, bool result, bool doForce=true, bool doApply=true)
 Act as if is a measurement was applied, except force the (usually random) result. More...
 
virtual bool M (bitLenInt q)
 
virtual void X (bitLenInt q)
 
virtual void INC (bitCapInt toAdd, bitLenInt start, bitLenInt length)
 Add integer (without sign) More...
 
virtual void DEC (bitCapInt toSub, bitLenInt start, bitLenInt length)
 Add integer (without sign) More...
 
virtual void INCC (bitCapInt toAdd, bitLenInt start, bitLenInt length, bitLenInt carryIndex)
 Add integer (without sign, with carry) More...
 
virtual void DECC (bitCapInt toSub, bitLenInt start, bitLenInt length, bitLenInt carryIndex)
 Subtract classical integer (without sign, with carry) More...
 
virtual void INCS (bitCapInt toAdd, bitLenInt start, bitLenInt length, bitLenInt overflowIndex)
 Add a classical integer to the register, with sign and without carry. More...
 
virtual void DECS (bitCapInt toSub, bitLenInt start, bitLenInt length, bitLenInt overflowIndex)
 Add a classical integer to the register, with sign and without carry. More...
 
virtual void CINC (bitCapInt toAdd, bitLenInt inOutStart, bitLenInt length, const bitLenInt *controls, bitLenInt controlLen)
 Add integer (without sign, with controls) More...
 
virtual void CDEC (bitCapInt toSub, bitLenInt inOutStart, bitLenInt length, const bitLenInt *controls, bitLenInt controlLen)
 Subtract integer (without sign, with controls) More...
 
virtual void INCDECC (bitCapInt toAdd, bitLenInt start, bitLenInt length, bitLenInt carryIndex)
 Common driver method behind INCC and DECC (without sign, with carry) More...
 
virtual void INCDECSC (bitCapInt toAdd, bitLenInt start, bitLenInt length, bitLenInt overflowIndex, bitLenInt carryIndex)
 Common driver method behind INCSC and DECSC (with overflow flag) More...
 
virtual void INCDECSC (bitCapInt toAdd, bitLenInt start, bitLenInt length, bitLenInt carryIndex)
 Common driver method behind INCSC and DECSC (without overflow flag) More...
 
virtual void INCBCD (bitCapInt toAdd, bitLenInt start, bitLenInt length)
 Add classical BCD integer (without sign) More...
 
virtual void INCDECBCDC (bitCapInt toAdd, bitLenInt start, bitLenInt length, bitLenInt carryIndex)
 Common driver method behind INCSC and DECSC (without overflow flag) More...
 
virtual void MUL (bitCapInt toMul, bitLenInt inOutStart, bitLenInt carryStart, bitLenInt length)
 Multiply by integer. More...
 
virtual void DIV (bitCapInt toDiv, bitLenInt inOutStart, bitLenInt carryStart, bitLenInt length)
 Divide by integer. More...
 
virtual void MULModNOut (bitCapInt toMul, bitCapInt modN, bitLenInt inStart, bitLenInt outStart, bitLenInt length)
 Multiplication modulo N by integer, (out of place) More...
 
virtual void IMULModNOut (bitCapInt toMul, bitCapInt modN, bitLenInt inStart, bitLenInt outStart, bitLenInt length)
 Inverse of multiplication modulo N by integer, (out of place) More...
 
virtual void POWModNOut (bitCapInt base, bitCapInt modN, bitLenInt inStart, bitLenInt outStart, bitLenInt length)
 Raise a classical base to a quantum power, modulo N, (out of place) More...
 
virtual void CMUL (bitCapInt toMul, bitLenInt inOutStart, bitLenInt carryStart, bitLenInt length, const bitLenInt *controls, bitLenInt controlLen)
 Controlled multiplication by integer. More...
 
virtual void CDIV (bitCapInt toDiv, bitLenInt inOutStart, bitLenInt carryStart, bitLenInt length, const bitLenInt *controls, bitLenInt controlLen)
 Controlled division by power of integer. More...
 
virtual void CMULModNOut (bitCapInt toMul, bitCapInt modN, bitLenInt inStart, bitLenInt outStart, bitLenInt length, const bitLenInt *controls, bitLenInt controlLen)
 Controlled multiplication modulo N by integer, (out of place) More...
 
virtual void CIMULModNOut (bitCapInt toMul, bitCapInt modN, bitLenInt inStart, bitLenInt outStart, bitLenInt length, const bitLenInt *controls, bitLenInt controlLen)
 Inverse of controlled multiplication modulo N by integer, (out of place) More...
 
virtual void CPOWModNOut (bitCapInt base, bitCapInt modN, bitLenInt inStart, bitLenInt outStart, bitLenInt length, const bitLenInt *controls, bitLenInt controlLen)
 Controlled, raise a classical base to a quantum power, modulo N, (out of place) More...
 
virtual bitCapInt IndexedLDA (bitLenInt indexStart, bitLenInt indexLength, bitLenInt valueStart, bitLenInt valueLength, const unsigned char *values, bool resetValue=true)
 Set 8 bit register bits by a superposed index-offset-based read from classical memory. More...
 
virtual bitCapInt IndexedADC (bitLenInt indexStart, bitLenInt indexLength, bitLenInt valueStart, bitLenInt valueLength, bitLenInt carryIndex, const unsigned char *values)
 Add to entangled 8 bit register state with a superposed index-offset-based read from classical memory. More...
 
virtual bitCapInt IndexedSBC (bitLenInt indexStart, bitLenInt indexLength, bitLenInt valueStart, bitLenInt valueLength, bitLenInt carryIndex, const unsigned char *values)
 Subtract from an entangled 8 bit register state with a superposed index-offset-based read from classical memory. More...
 
virtual void Hash (bitLenInt start, bitLenInt length, const unsigned char *values)
 Transform a length of qubit register via lookup through a hash table. More...
 
virtual void CPhaseFlipIfLess (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. More...
 
virtual void PhaseFlipIfLess (bitCapInt greaterPerm, bitLenInt start, bitLenInt length)
 This is an expedient for an adaptive Grover's search for a function's global minimum. More...
 
virtual void Swap (bitLenInt qubitIndex1, bitLenInt qubitIndex2)
 Swap values of two bits in register. More...
 
virtual void ISwap (bitLenInt qubitIndex1, bitLenInt qubitIndex2)
 Swap values of two bits in register, and apply phase factor of i if bits are different. More...
 
virtual void SqrtSwap (bitLenInt qubitIndex1, bitLenInt qubitIndex2)
 Square root of Swap gate. More...
 
virtual void ISqrtSwap (bitLenInt qubitIndex1, bitLenInt qubitIndex2)
 Inverse square root of Swap gate. More...
 
virtual void FSim (real1_f theta, real1_f phi, bitLenInt qubitIndex1, bitLenInt qubitIndex2)
 The 2-qubit "fSim" gate, (useful in the simulation of particles with fermionic statistics) More...
 
virtual real1_f Prob (bitLenInt qubitIndex)
 Direct measure of bit probability to be in |1> state. More...
 
virtual real1_f ProbMask (bitCapInt mask, bitCapInt permutation)
 Direct measure of masked permutation probability. More...
 
virtual real1_f ProbParity (bitCapInt mask)
 Overall probability of any odd permutation of the masked set of bits. More...
 
virtual bool ForceMParity (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 (usually random) result. More...
 
virtual real1_f ExpectationBitsAll (const bitLenInt *bits, bitLenInt length, bitCapInt offset=0)
 Get permutation expectation value of bits. More...
 
virtual 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 probability or measurement operation. More...
 
virtual void NormalizeState (real1_f nrm=REAL1_DEFAULT_ARG, real1_f norm_thresh=REAL1_DEFAULT_ARG, real1_f phaseArg=ZERO_R1)
 Apply the normalization factor found by UpdateRunningNorm() or on the fly by a single bit gate. More...
 
virtual void Finish ()
 If asynchronous work is still running, block until it finishes. More...
 
virtual bool isFinished ()
 Returns "false" if asynchronous work is still running, and "true" if all previously dispatched asynchronous work is done. More...
 
virtual void Dump ()
 If asynchronous work is still running, let the simulator know that it can be aborted. More...
 
virtual QInterfacePtr Clone ()
 Clone this QInterface. More...
 
virtual void SetDevice (int dID, bool forceReInit=false)
 Set the device index, if more than one device is available. More...
 
virtual int64_t GetDevice ()
 Get the device index. More...
 
bitCapIntOcl GetMaxSize ()
 
virtual real1_f SumSqrDiff (QInterfacePtr toCompare)
 
virtual real1_f SumSqrDiff (QPagerPtr toCompare)
 
- Public Member Functions inherited from Qrack::QAlu
virtual void INCSC (bitCapInt toAdd, bitLenInt start, bitLenInt length, bitLenInt overflowIndex, bitLenInt carryIndex)
 Add a classical integer to the register, with sign and with carry. More...
 
virtual void INCSC (bitCapInt toAdd, bitLenInt start, bitLenInt length, bitLenInt carryIndex)
 Add a classical integer to the register, with sign and with (phase-based) carry. More...
 
virtual void DECSC (bitCapInt toSub, bitLenInt start, bitLenInt length, bitLenInt overflowIndex, bitLenInt carryIndex)
 Subtract a classical integer from the register, with sign and with carry. More...
 
virtual void DECSC (bitCapInt toSub, bitLenInt start, bitLenInt length, bitLenInt carryIndex)
 Subtract a classical integer from the register, with sign and with carry. More...
 
virtual void DECBCD (bitCapInt toSub, bitLenInt start, bitLenInt length)
 Subtract classical BCD integer (without sign) More...
 
virtual void INCBCDC (bitCapInt toAdd, bitLenInt start, bitLenInt length, bitLenInt carryIndex)
 Add classical BCD integer (without sign, with carry) More...
 
virtual void DECBCDC (bitCapInt toSub, bitLenInt start, bitLenInt length, bitLenInt carryIndex)
 Subtract BCD integer (without sign, with carry) More...
 
- Public Member Functions inherited from Qrack::QParity
virtual bool MParity (bitCapInt mask)
 Measure (and collapse) parity of the masked set of qubits. More...
 
- Public Member Functions inherited from Qrack::QInterface
 QInterface (bitLenInt n, qrack_rand_gen_ptr rgp=nullptr, bool doNorm=false, bool useHardwareRNG=true, bool randomGlobalPhase=true, real1_f norm_thresh=REAL1_EPSILON)
 
 QInterface ()
 Default constructor, primarily for protected internal use. More...
 
virtual ~QInterface ()
 
void SetRandomSeed (uint32_t seed)
 
virtual bitLenInt GetQubitCount ()
 Get the count of bits in this register. More...
 
virtual bitCapInt GetMaxQPower ()
 Get the maximum number of basis states, namely \( 2^n \) for \( n \) qubits. More...
 
virtual bool GetIsArbitraryGlobalPhase ()
 
real1_f Rand ()
 Generate a random real number between 0 and 1. More...
 
virtual std::map< QInterfacePtr, bitLenIntCompose (std::vector< QInterfacePtr > toCopy)
 
virtual bitLenInt Compose (QInterfacePtr toCopy, bitLenInt start)
 
virtual void MCPhase (const bitLenInt *controls, bitLenInt controlLen, complex topLeft, complex bottomRight, bitLenInt target)
 Apply a single bit transformation that only effects phase, with arbitrary control bits. More...
 
virtual void MCInvert (const bitLenInt *controls, bitLenInt controlLen, complex topRight, complex bottomLeft, bitLenInt target)
 Apply a single bit transformation that reverses bit probability and might effect phase, with arbitrary control bits. More...
 
virtual void MACPhase (const bitLenInt *controls, bitLenInt controlLen, complex topLeft, complex bottomRight, bitLenInt target)
 Apply a single bit transformation that only effects phase, with arbitrary (anti-)control bits. More...
 
virtual void MACInvert (const bitLenInt *controls, bitLenInt controlLen, complex topRight, complex bottomLeft, bitLenInt target)
 Apply a single bit transformation that reverses bit probability and might effect phase, with arbitrary (anti-)control bits. More...
 
virtual void UniformlyControlledSingleBit (const bitLenInt *controls, bitLenInt controlLen, bitLenInt qubitIndex, const complex *mtrxs)
 Apply a "uniformly controlled" arbitrary single bit unitary transformation. More...
 
virtual void UniformlyControlledSingleBit (const bitLenInt *controls, bitLenInt controlLen, bitLenInt qubitIndex, const complex *mtrxs, const bitCapInt *mtrxSkipPowers, bitLenInt mtrxSkipLen, bitCapInt mtrxSkipValueMask)
 
virtual void TimeEvolve (Hamiltonian h, real1_f timeDiff)
 To define a Hamiltonian, give a vector of controlled single bit gates ("HamiltonianOp" instances) that are applied by left-multiplication in low-to-high vector index order on the state vector. More...
 
virtual void CSwap (const bitLenInt *controls, bitLenInt controlLen, bitLenInt qubit1, bitLenInt qubit2)
 Apply a swap with arbitrary control bits. More...
 
virtual void AntiCSwap (const bitLenInt *controls, bitLenInt controlLen, bitLenInt qubit1, bitLenInt qubit2)
 Apply a swap with arbitrary (anti) control bits. More...
 
virtual void CSqrtSwap (const bitLenInt *controls, bitLenInt controlLen, bitLenInt qubit1, bitLenInt qubit2)
 Apply a square root of swap with arbitrary control bits. More...
 
virtual void AntiCSqrtSwap (const bitLenInt *controls, bitLenInt controlLen, bitLenInt qubit1, bitLenInt qubit2)
 Apply a square root of swap with arbitrary (anti) control bits. More...
 
virtual void CISqrtSwap (const bitLenInt *controls, bitLenInt controlLen, bitLenInt qubit1, bitLenInt qubit2)
 Apply an inverse square root of swap with arbitrary control bits. More...
 
virtual void AntiCISqrtSwap (const bitLenInt *controls, bitLenInt controlLen, bitLenInt qubit1, bitLenInt qubit2)
 Apply an inverse square root of swap with arbitrary (anti) control bits. More...
 
virtual void CCNOT (bitLenInt control1, bitLenInt control2, bitLenInt target)
 Doubly-controlled NOT gate. More...
 
virtual void AntiCCNOT (bitLenInt control1, bitLenInt control2, bitLenInt target)
 Anti doubly-controlled NOT gate. More...
 
virtual void CNOT (bitLenInt control, bitLenInt target)
 Controlled NOT gate. More...
 
virtual void AntiCNOT (bitLenInt control, bitLenInt target)
 Anti controlled NOT gate. More...
 
virtual void CY (bitLenInt control, bitLenInt target)
 Controlled Y gate. More...
 
virtual void AntiCY (bitLenInt control, bitLenInt target)
 Anti controlled Y gate. More...
 
virtual void CCY (bitLenInt control1, bitLenInt control2, bitLenInt target)
 Doubly-Controlled Y gate. More...
 
virtual void AntiCCY (bitLenInt control1, bitLenInt control2, bitLenInt target)
 Anti doubly-controlled Y gate. More...
 
virtual void CZ (bitLenInt control, bitLenInt target)
 Controlled Z gate. More...
 
virtual void AntiCZ (bitLenInt control, bitLenInt target)
 Anti controlled Z gate. More...
 
virtual void CCZ (bitLenInt control1, bitLenInt control2, bitLenInt target)
 Doubly-Controlled Z gate. More...
 
virtual void AntiCCZ (bitLenInt control1, bitLenInt control2, bitLenInt target)
 Anti doubly-controlled Z gate. More...
 
virtual void U (bitLenInt target, real1_f theta, real1_f phi, real1_f lambda)
 General unitary gate. More...
 
virtual void U2 (bitLenInt target, real1_f phi, real1_f lambda)
 2-parameter unitary gate More...
 
virtual void IU2 (bitLenInt target, real1_f phi, real1_f lambda)
 Inverse 2-parameter unitary gate. More...
 
virtual void AI (bitLenInt target, real1_f azimuth, real1_f inclination)
 "Azimuth, Inclination" (RY-RZ) More...
 
virtual void IAI (bitLenInt target, real1_f azimuth, real1_f inclination)
 Invert "Azimuth, Inclination" (RY-RZ) More...
 
virtual void CU (const bitLenInt *controls, bitLenInt controlLen, bitLenInt target, real1_f theta, real1_f phi, real1_f lambda)
 Controlled general unitary gate. More...
 
virtual void AntiCU (const bitLenInt *controls, bitLenInt controlLen, bitLenInt target, real1_f theta, real1_f phi, real1_f lambda)
 (Anti-)Controlled general unitary gate More...
 
virtual void H (bitLenInt qubitIndex)
 Hadamard gate. More...
 
virtual void SqrtH (bitLenInt qubitIndex)
 Square root of Hadamard gate. More...
 
virtual void SH (bitLenInt qubitIndex)
 Y-basis transformation gate. More...
 
virtual void HIS (bitLenInt qubitIndex)
 Y-basis (inverse) transformation gate. More...
 
virtual void S (bitLenInt qubitIndex)
 S gate. More...
 
virtual void IS (bitLenInt qubitIndex)
 Inverse S gate. More...
 
virtual void T (bitLenInt qubitIndex)
 T gate. More...
 
virtual void IT (bitLenInt qubitIndex)
 Inverse T gate. More...
 
virtual void PhaseRootN (bitLenInt n, bitLenInt qubitIndex)
 "PhaseRootN" gate More...
 
virtual void IPhaseRootN (bitLenInt n, bitLenInt qubitIndex)
 Inverse "PhaseRootN" gate. More...
 
virtual void Y (bitLenInt qubitIndex)
 Y gate. More...
 
virtual void YMask (bitCapInt mask)
 Masked Y gate. More...
 
virtual void Z (bitLenInt qubitIndex)
 Z gate. More...
 
virtual void SqrtX (bitLenInt qubitIndex)
 Square root of X gate. More...
 
virtual void ISqrtX (bitLenInt qubitIndex)
 Inverse square root of X gate. More...
 
virtual void SqrtXConjT (bitLenInt qubitIndex)
 Phased square root of X gate. More...
 
virtual void ISqrtXConjT (bitLenInt qubitIndex)
 Inverse phased square root of X gate. More...
 
virtual void SqrtY (bitLenInt qubitIndex)
 Square root of Y gate. More...
 
virtual void ISqrtY (bitLenInt qubitIndex)
 Square root of Y gate. More...
 
virtual void CH (bitLenInt control, bitLenInt target)
 Controlled H gate. More...
 
virtual void AntiCH (bitLenInt control, bitLenInt target)
 (Anti-)controlled H gate More...
 
virtual void CS (bitLenInt control, bitLenInt target)
 Controlled S gate. More...
 
virtual void AntiCS (bitLenInt control, bitLenInt target)
 (Anti-)controlled S gate More...
 
virtual void CIS (bitLenInt control, bitLenInt target)
 Controlled inverse S gate. More...
 
virtual void AntiCIS (bitLenInt control, bitLenInt target)
 (Anti-)controlled inverse S gate More...
 
virtual void CT (bitLenInt control, bitLenInt target)
 Controlled T gate. More...
 
virtual void CIT (bitLenInt control, bitLenInt target)
 Controlled inverse T gate. More...
 
virtual void CPhaseRootN (bitLenInt n, bitLenInt control, bitLenInt target)
 Controlled "PhaseRootN" gate. More...
 
virtual void AntiCPhaseRootN (bitLenInt n, bitLenInt control, bitLenInt target)
 (Anti-)controlled "PhaseRootN" gate More...
 
virtual void CIPhaseRootN (bitLenInt n, bitLenInt control, bitLenInt target)
 Controlled inverse "PhaseRootN" gate. More...
 
virtual void AntiCIPhaseRootN (bitLenInt n, bitLenInt control, bitLenInt target)
 (Anti-)controlled inverse "PhaseRootN" gate More...
 
virtual void AND (bitLenInt inputBit1, bitLenInt inputBit2, bitLenInt outputBit)
 Quantum analog of classical "AND" gate. More...
 
virtual void OR (bitLenInt inputBit1, bitLenInt inputBit2, bitLenInt outputBit)
 Quantum analog of classical "OR" gate. More...
 
virtual void XOR (bitLenInt inputBit1, bitLenInt inputBit2, bitLenInt outputBit)
 Quantum analog of classical "XOR" gate. More...
 
virtual void CLAND (bitLenInt inputQBit, bool inputClassicalBit, bitLenInt outputBit)
 Quantum analog of classical "AND" gate. More...
 
virtual void CLOR (bitLenInt inputQBit, bool inputClassicalBit, bitLenInt outputBit)
 Quantum analog of classical "OR" gate. More...
 
virtual void CLXOR (bitLenInt inputQBit, bool inputClassicalBit, bitLenInt outputBit)
 Quantum analog of classical "XOR" gate. More...
 
virtual void NAND (bitLenInt inputBit1, bitLenInt inputBit2, bitLenInt outputBit)
 Quantum analog of classical "NAND" gate. More...
 
virtual void NOR (bitLenInt inputBit1, bitLenInt inputBit2, bitLenInt outputBit)
 Quantum analog of classical "NOR" gate. More...
 
virtual void XNOR (bitLenInt inputBit1, bitLenInt inputBit2, bitLenInt outputBit)
 Quantum analog of classical "XNOR" gate. More...
 
virtual void CLNAND (bitLenInt inputQBit, bool inputClassicalBit, bitLenInt outputBit)
 Quantum analog of classical "NAND" gate. More...
 
virtual void CLNOR (bitLenInt inputQBit, bool inputClassicalBit, bitLenInt outputBit)
 Quantum analog of classical "NOR" gate. More...
 
virtual void CLXNOR (bitLenInt inputQBit, bool inputClassicalBit, bitLenInt outputBit)
 Quantum analog of classical "XNOR" gate. More...
 
virtual void UniformlyControlledRY (const bitLenInt *controls, bitLenInt controlLen, bitLenInt qubitIndex, const real1 *angles)
 Apply a "uniformly controlled" rotation of a bit around the Pauli Y axis. More...
 
virtual void UniformlyControlledRZ (const bitLenInt *controls, bitLenInt controlLen, bitLenInt qubitIndex, const real1 *angles)
 Apply a "uniformly controlled" rotation of a bit around the Pauli Z axis. More...
 
virtual void RT (real1_f radians, bitLenInt qubitIndex)
 Phase shift gate. More...
 
virtual void RX (real1_f radians, bitLenInt qubitIndex)
 X axis rotation gate. More...
 
virtual void RY (real1_f radians, bitLenInt qubitIndex)
 Y axis rotation gate. More...
 
virtual void RZ (real1_f radians, bitLenInt qubitIndex)
 Z axis rotation gate. More...
 
virtual void CRZ (real1_f radians, bitLenInt control, bitLenInt target)
 Controlled Z axis rotation gate. More...
 
virtual void RTDyad (int numerator, int denomPower, bitLenInt qubitIndex)
 Dyadic fraction phase shift gate. More...
 
virtual void RXDyad (int numerator, int denomPower, bitLenInt qubitIndex)
 Dyadic fraction X axis rotation gate. More...
 
virtual void Exp (real1_f radians, bitLenInt qubitIndex)
 (Identity) Exponentiation gate More...
 
virtual void Exp (const bitLenInt *controls, bitLenInt controlLen, bitLenInt qubit, const complex *matrix2x2, bool antiCtrled=false)
 Imaginary exponentiation of arbitrary 2x2 gate. More...
 
virtual void ExpDyad (int numerator, int denomPower, bitLenInt qubitIndex)
 Dyadic fraction (identity) exponentiation gate. More...
 
virtual void ExpX (real1_f radians, bitLenInt qubitIndex)
 Pauli X exponentiation gate. More...
 
virtual void ExpXDyad (int numerator, int denomPower, bitLenInt qubitIndex)
 Dyadic fraction Pauli X exponentiation gate. More...
 
virtual void ExpY (real1_f radians, bitLenInt qubitIndex)
 Pauli Y exponentiation gate. More...
 
virtual void ExpYDyad (int numerator, int denomPower, bitLenInt qubitIndex)
 Dyadic fraction Pauli Y exponentiation gate. More...
 
virtual void ExpZ (real1_f radians, bitLenInt qubitIndex)
 Pauli Z exponentiation gate. More...
 
virtual void ExpZDyad (int numerator, int denomPower, bitLenInt qubitIndex)
 Dyadic fraction Pauli Z exponentiation gate. More...
 
virtual void CRX (real1_f radians, bitLenInt control, bitLenInt target)
 Controlled X axis rotation gate. More...
 
virtual void CRXDyad (int numerator, int denomPower, bitLenInt control, bitLenInt target)
 Controlled dyadic fraction X axis rotation gate. More...
 
virtual void RYDyad (int numerator, int denomPower, bitLenInt qubitIndex)
 Dyadic fraction Y axis rotation gate. More...
 
virtual void CRY (real1_f radians, bitLenInt control, bitLenInt target)
 Controlled Y axis rotation gate. More...
 
virtual void CRYDyad (int numerator, int denomPower, bitLenInt control, bitLenInt target)
 Controlled dyadic fraction y axis rotation gate. More...
 
virtual void RZDyad (int numerator, int denomPower, bitLenInt qubitIndex)
 Dyadic fraction Z axis rotation gate. More...
 
virtual void CRZDyad (int numerator, int denomPower, bitLenInt control, bitLenInt target)
 Controlled dyadic fraction Z axis rotation gate. More...
 
virtual void CRT (real1_f radians, bitLenInt control, bitLenInt target)
 Controlled "phase shift gate". More...
 
virtual void CRTDyad (int numerator, int denomPower, bitLenInt control, bitLenInt target)
 Controlled dyadic fraction "phase shift gate". More...
 
virtual void H (bitLenInt start, bitLenInt length)
 Bitwise Hadamard. More...
 
virtual void X (bitLenInt start, bitLenInt length)
 Bitwise Pauli X (or logical "NOT") operator. More...
 
virtual void ROL (bitLenInt shift, bitLenInt start, bitLenInt length)
 Circular shift left - shift bits left, and carry last bits. More...
 
virtual void ROR (bitLenInt shift, bitLenInt start, bitLenInt length)
 Circular shift right - shift bits right, and carry first bits. More...
 
virtual void ASL (bitLenInt shift, bitLenInt start, bitLenInt length)
 Arithmetic shift left, with last 2 bits as sign and carry. More...
 
virtual void ASR (bitLenInt shift, bitLenInt start, bitLenInt length)
 Arithmetic shift right, with last 2 bits as sign and carry. More...
 
virtual void LSL (bitLenInt shift, bitLenInt start, bitLenInt length)
 Logical shift left, filling the extra bits with |0> More...
 
virtual void LSR (bitLenInt shift, bitLenInt start, bitLenInt length)
 Logical shift right, filling the extra bits with |0> More...
 
virtual void FullAdd (bitLenInt inputBit1, bitLenInt inputBit2, bitLenInt carryInSumOut, bitLenInt carryOut)
 Quantum analog of classical "Full Adder" gate. More...
 
virtual void IFullAdd (bitLenInt inputBit1, bitLenInt inputBit2, bitLenInt carryInSumOut, bitLenInt carryOut)
 Inverse of FullAdd. More...
 
virtual void CFullAdd (const bitLenInt *controls, bitLenInt controlLen, bitLenInt inputBit1, bitLenInt inputBit2, bitLenInt carryInSumOut, bitLenInt carryOut)
 Controlled quantum analog of classical "Full Adder" gate. More...
 
virtual void CIFullAdd (const bitLenInt *controls, bitLenInt controlLen, bitLenInt inputBit1, bitLenInt inputBit2, bitLenInt carryInSumOut, bitLenInt carryOut)
 Inverse of CFullAdd. More...
 
virtual void ADC (bitLenInt input1, bitLenInt input2, bitLenInt output, bitLenInt length, bitLenInt carry)
 Add a quantum integer to a quantum integer, with carry. More...
 
virtual void IADC (bitLenInt input1, bitLenInt input2, bitLenInt output, bitLenInt length, bitLenInt carry)
 Inverse of ADC. More...
 
virtual void CADC (const bitLenInt *controls, bitLenInt controlLen, bitLenInt input1, bitLenInt input2, bitLenInt output, bitLenInt length, bitLenInt carry)
 Add a quantum integer to a quantum integer, with carry and with controls. More...
 
virtual void CIADC (const bitLenInt *controls, bitLenInt controlLen, bitLenInt input1, bitLenInt input2, bitLenInt output, bitLenInt length, bitLenInt carry)
 Inverse of CADC. More...
 
virtual void QFT (bitLenInt start, bitLenInt length, bool trySeparate=false)
 Quantum Fourier Transform - Apply the quantum Fourier transform to the register. More...
 
virtual void QFTR (const bitLenInt *qubits, bitLenInt length, bool trySeparate=false)
 Quantum Fourier Transform (random access) - Apply the quantum Fourier transform to the register. More...
 
virtual void IQFT (bitLenInt start, bitLenInt length, bool trySeparate=false)
 Inverse Quantum Fourier Transform - Apply the inverse quantum Fourier transform to the register. More...
 
virtual void IQFTR (const bitLenInt *qubits, bitLenInt length, bool trySeparate=false)
 Inverse Quantum Fourier Transform (random access) - Apply the inverse quantum Fourier transform to the register. More...
 
virtual void ZeroPhaseFlip (bitLenInt start, bitLenInt length)
 Reverse the phase of the state where the register equals zero. More...
 
virtual void PhaseFlip ()
 Phase flip always - equivalent to Z X Z X on any bit in the QInterface. More...
 
virtual void SetReg (bitLenInt start, bitLenInt length, bitCapInt value)
 Set register bits to given permutation. More...
 
virtual bitCapInt MReg (bitLenInt start, bitLenInt length)
 Measure permutation state of a register. More...
 
virtual bitCapInt MAll ()
 Measure permutation state of all coherent bits. More...
 
virtual bitCapInt ForceMReg (bitLenInt start, bitLenInt length, bitCapInt result, bool doForce=true, bool doApply=true)
 Act as if is a measurement was applied, except force the (usually random) result. More...
 
virtual bitCapInt M (const bitLenInt *bits, bitLenInt length)
 Measure bits with indices in array, and return a mask of the results. More...
 
virtual bitCapInt ForceM (const bitLenInt *bits, bitLenInt length, const bool *values, bool doApply=true)
 Measure bits with indices in array, and return a mask of the results. More...
 
virtual void Reverse (bitLenInt first, bitLenInt last)
 Reverse all of the bits in a sequence. More...
 
virtual real1_f ProbReg (bitLenInt start, bitLenInt length, bitCapInt permutation)
 Direct measure of register permutation probability. More...
 
virtual void ProbMaskAll (bitCapInt mask, real1 *probsArray)
 Direct measure of masked permutation probability. More...
 
virtual void ProbBitsAll (const bitLenInt *bits, bitLenInt length, real1 *probsArray)
 Direct measure of listed permutation probability. More...
 
virtual std::map< bitCapInt, int > MultiShotMeasureMask (const bitCapInt *qPowers, bitLenInt qPowerCount, unsigned shots)
 Statistical measure of masked permutation probability. More...
 
virtual void MultiShotMeasureMask (const bitCapInt *qPowers, bitLenInt qPowerCount, unsigned shots, unsigned *shotsArray)
 Statistical measure of masked permutation probability (returned as array) More...
 
virtual void SetBit (bitLenInt qubitIndex1, bool value)
 Set individual bit to pure |0> (false) or |1> (true) state. More...
 
virtual bool ApproxCompare (QInterfacePtr toCompare, real1_f error_tol=TRYDECOMPOSE_EPSILON)
 Compare state vectors approximately, component by component, to determine whether this state vector is the same as the target. More...
 
virtual bool TryDecompose (bitLenInt start, QInterfacePtr dest, real1_f error_tol=TRYDECOMPOSE_EPSILON)
 
virtual bool isBinaryDecisionTree ()
 Returns "true" if current state representation is definitely a binary decision tree, "false" if it is definitely not, or "true" if it cannot be determined. More...
 
virtual bool isClifford ()
 Returns "true" if current state is identifiably within the Clifford set, or "false" if it is not or cannot be determined. More...
 
virtual bool isClifford (bitLenInt qubit)
 Returns "true" if current qubit state is identifiably within the Clifford set, or "false" if it is not or cannot be determined. More...
 
virtual bool TrySeparate (const bitLenInt *qubits, bitLenInt length, real1_f error_tol)
 Qrack::QUnit types maintain explicit separation of representations of qubits, which reduces memory usage and increases gate speed. More...
 
virtual bool TrySeparate (bitLenInt qubit)
 Single-qubit TrySeparate() More...
 
virtual bool TrySeparate (bitLenInt qubit1, bitLenInt qubit2)
 Two-qubit TrySeparate() More...
 
virtual void SetReactiveSeparate (bool isAggSep)
 Set reactive separation option (on by default if available) More...
 
virtual bool GetReactiveSeparate ()
 Get reactive separation option. More...
 
bitCapIntOcl GetMaxSize ()
 Get maximum number of amplitudes that can be allocated on current device. More...
 
virtual real1_f FirstNonzeroPhase ()
 Get phase of lowest permutation nonzero amplitude. More...
 
- Public Member Functions inherited from Qrack::ParallelFor
 ParallelFor ()
 
virtual ~ParallelFor ()
 
void SetConcurrencyLevel (unsigned num)
 
unsigned GetConcurrencyLevel ()
 
bitCapIntOcl GetStride ()
 
void par_for_inc (const bitCapIntOcl begin, const bitCapIntOcl itemCount, IncrementFunc, ParallelFunc fn)
 Iterate through the permutations a maximum of end-begin times, allowing the caller to control the incrementation offset through 'inc'. More...
 
void par_for (const bitCapIntOcl begin, const bitCapIntOcl end, ParallelFunc fn)
 Call fn once for every numerical value between begin and end. More...
 
void par_for_skip (const bitCapIntOcl begin, const bitCapIntOcl end, const bitCapIntOcl skipPower, const bitLenInt skipBitCount, ParallelFunc fn)
 Skip over the skipPower bits. More...
 
void par_for_mask (const bitCapIntOcl, const bitCapIntOcl, const bitCapIntOcl *maskArray, const bitLenInt maskLen, ParallelFunc fn)
 Skip over the bits listed in maskArray in the same fashion as par_for_skip. More...
 
void par_for_set (const std::set< bitCapIntOcl > &sparseSet, ParallelFunc fn)
 Iterate over a sparse state vector. More...
 
void par_for_set (const std::vector< bitCapIntOcl > &sparseSet, ParallelFunc fn)
 Iterate over a sparse state vector. More...
 
void par_for_sparse_compose (const std::vector< bitCapIntOcl > &lowSet, const std::vector< bitCapIntOcl > &highSet, const bitLenInt &highStart, ParallelFunc fn)
 Iterate over the power set of 2 sparse state vectors. More...
 
void par_for_qbdt (const bitCapInt begin, const bitCapInt end, BdtFunc fn)
 Iterate over a QBDT tree. More...
 
real1_f par_norm (const bitCapIntOcl maxQPower, const StateVectorPtr stateArray, real1_f norm_thresh=ZERO_R1)
 Calculate the normal for the array, (with flooring). More...
 
real1_f par_norm_exact (const bitCapIntOcl maxQPower, const StateVectorPtr stateArray)
 Calculate the normal for the array, (without flooring.) More...
 

Protected Member Functions

QEnginePtr MakeEngine (bitLenInt length, bitCapInt perm, int deviceId)
 
virtual void SetQubitCount (bitLenInt qb)
 
bitCapInt pageMaxQPower ()
 
bitLenInt pagedQubitCount ()
 
bitLenInt qubitsPerPage ()
 
void CombineEngines (bitLenInt thresholdBits)
 
void CombineEngines ()
 
void SeparateEngines (bitLenInt thresholdBits, bool noBaseFloor=false)
 
void SeparateEngines ()
 
template<typename Qubit1Fn >
void SingleBitGate (bitLenInt target, Qubit1Fn fn, bool isSqiCtrl=false, bool isAnti=false)
 
template<typename Qubit1Fn >
void MetaControlled (bool anti, const std::vector< bitLenInt > &controls, bitLenInt target, Qubit1Fn fn, const complex *mtrx, bool isSqiCtrl=false, bool isIntraCtrled=false)
 
template<typename Qubit1Fn >
void SemiMetaControlled (bool anti, std::vector< bitLenInt > controls, bitLenInt target, Qubit1Fn fn)
 
void MetaSwap (bitLenInt qubit1, bitLenInt qubit2, bool isIPhaseFac)
 
void SemiMetaSwap (bitLenInt qubit1, bitLenInt qubit2, bool isIPhaseFac)
 
template<typename F >
void CombineAndOp (F fn, std::vector< bitLenInt > bits)
 
template<typename F >
void CombineAndOpControlled (F fn, std::vector< bitLenInt > bits, const bitLenInt *controls, bitLenInt controlLen)
 
void ApplySingleEither (bool isInvert, complex top, complex bottom, bitLenInt target)
 
void ApplyEitherControlledSingleBit (bool anti, const bitLenInt *controls, bitLenInt controlLen, bitLenInt target, const complex *mtrx)
 
void Init ()
 
- Protected Member Functions inherited from Qrack::QInterface
template<typename GateFunc >
void ControlledLoopFixture (bitLenInt length, GateFunc gate)
 
void FreeAligned (void *toFree)
 
complex GetNonunitaryPhase ()
 
template<typename Fn >
void MACWrapper (const bitLenInt *controls, bitLenInt controlLen, Fn fn)
 

Protected Attributes

std::vector< QInterfaceEngineengines
 
QInterfaceEngine rootEngine
 
int devID
 
complex phaseFactor
 
bool useHostRam
 
bool useRDRAND
 
bool isSparse
 
std::vector< QEnginePtrqPages
 
std::vector< int > deviceIDs
 
bool useHardwareThreshold
 
bool useGpuThreshold
 
bitLenInt segmentGlobalQb
 
bitLenInt minPageQubits
 
bitLenInt maxPageQubits
 
bitLenInt thresholdQubitsPerPage
 
bitLenInt baseQubitsPerPage
 
bitCapInt basePageCount
 
bitCapIntOcl basePageMaxQPower
 
bitLenInt maxQubits
 
- Protected Attributes inherited from Qrack::QInterface
bitLenInt qubitCount
 
bitCapInt maxQPower
 
uint32_t randomSeed
 
qrack_rand_gen_ptr rand_generator
 
std::uniform_real_distribution< real1_frand_distribution
 
std::shared_ptr< RdRandomhardware_rand_generator
 
bool doNormalize
 
bool randGlobalPhase
 
bool useRDRAND
 
real1 amplitudeFloor
 

Additional Inherited Members

- Static Protected Member Functions inherited from Qrack::QInterface
static real1_f normHelper (complex c)
 
static real1_f clampProb (real1_f toClamp)
 

Detailed Description

A "Qrack::QPager" splits a "Qrack::QEngine" implementation into equal-length "pages." This helps both optimization and distribution of a single coherent quantum register across multiple devices.

Constructor & Destructor Documentation

◆ QPager() [1/3]

Qrack::QPager::QPager ( std::vector< QInterfaceEngine eng,
bitLenInt  qBitCount,
bitCapInt  initState = 0,
qrack_rand_gen_ptr  rgp = nullptr,
complex  phaseFac = CMPLX_DEFAULT_ARG,
bool  doNorm = false,
bool  ignored = false,
bool  useHostMem = false,
int  deviceId = -1,
bool  useHardwareRNG = true,
bool  useSparseStateVec = false,
real1_f  norm_thresh = REAL1_EPSILON,
std::vector< int >  devList = {},
bitLenInt  qubitThreshold = 0,
real1_f  separation_thresh = FP_NORM_EPSILON 
)

◆ QPager() [2/3]

Qrack::QPager::QPager ( bitLenInt  qBitCount,
bitCapInt  initState = 0,
qrack_rand_gen_ptr  rgp = nullptr,
complex  phaseFac = CMPLX_DEFAULT_ARG,
bool  doNorm = false,
bool  ignored = false,
bool  useHostMem = false,
int  deviceId = -1,
bool  useHardwareRNG = true,
bool  useSparseStateVec = false,
real1_f  norm_thresh = REAL1_EPSILON,
std::vector< int >  devList = {},
bitLenInt  qubitThreshold = 0,
real1_f  separation_thresh = FP_NORM_EPSILON 
)
inline

◆ QPager() [3/3]

Qrack::QPager::QPager ( QEnginePtr  enginePtr,
std::vector< QInterfaceEngine eng,
bitLenInt  qBitCount,
bitCapInt  ignored = 0,
qrack_rand_gen_ptr  rgp = nullptr,
complex  phaseFac = CMPLX_DEFAULT_ARG,
bool  doNorm = false,
bool  ignored2 = false,
bool  useHostMem = false,
int  deviceId = -1,
bool  useHardwareRNG = true,
bool  useSparseStateVec = false,
real1_f  norm_thresh = REAL1_EPSILON,
std::vector< int >  devList = {},
bitLenInt  qubitThreshold = 0,
real1_f  separation_thresh = FP_NORM_EPSILON 
)

Member Function Documentation

◆ ApplyEitherControlledSingleBit()

void Qrack::QPager::ApplyEitherControlledSingleBit ( bool  anti,
const bitLenInt controls,
bitLenInt  controlLen,
bitLenInt  target,
const complex mtrx 
)
protected

◆ ApplySingleEither()

void Qrack::QPager::ApplySingleEither ( bool  isInvert,
complex  top,
complex  bottom,
bitLenInt  target 
)
protected

◆ CDEC()

virtual void Qrack::QPager::CDEC ( bitCapInt  toSub,
bitLenInt  start,
bitLenInt  length,
const bitLenInt controls,
bitLenInt  controlLen 
)
inlinevirtual

Subtract integer (without sign, with controls)

Reimplemented from Qrack::QAlu.

◆ CDIV()

void Qrack::QPager::CDIV ( bitCapInt  toDiv,
bitLenInt  start,
bitLenInt  carryStart,
bitLenInt  length,
const bitLenInt controls,
bitLenInt  controlLen 
)
virtual

Controlled division by power of integer.

Implements Qrack::QAlu.

◆ CIMULModNOut()

void Qrack::QPager::CIMULModNOut ( bitCapInt  toMul,
bitCapInt  modN,
bitLenInt  inStart,
bitLenInt  outStart,
bitLenInt  length,
const bitLenInt controls,
bitLenInt  controlLen 
)
virtual

Inverse of controlled multiplication modulo N by integer, (out of place)

Implements Qrack::QAlu.

◆ CINC()

virtual void Qrack::QPager::CINC ( bitCapInt  toAdd,
bitLenInt  start,
bitLenInt  length,
const bitLenInt controls,
bitLenInt  controlLen 
)
inlinevirtual

Add integer (without sign, with controls)

Implements Qrack::QAlu.

◆ Clone()

QInterfacePtr Qrack::QPager::Clone ( )
virtual

Clone this QInterface.

Implements Qrack::QInterface.

◆ CMUL()

void Qrack::QPager::CMUL ( bitCapInt  toMul,
bitLenInt  start,
bitLenInt  carryStart,
bitLenInt  length,
const bitLenInt controls,
bitLenInt  controlLen 
)
virtual

Controlled multiplication by integer.

Implements Qrack::QAlu.

◆ CMULModNOut()

void Qrack::QPager::CMULModNOut ( bitCapInt  toMul,
bitCapInt  modN,
bitLenInt  inStart,
bitLenInt  outStart,
bitLenInt  length,
const bitLenInt controls,
bitLenInt  controlLen 
)
virtual

Controlled multiplication modulo N by integer, (out of place)

Implements Qrack::QAlu.

◆ CombineAndOp()

template<typename F >
void Qrack::QPager::CombineAndOp ( fn,
std::vector< bitLenInt bits 
)
protected

◆ CombineAndOpControlled()

template<typename F >
void Qrack::QPager::CombineAndOpControlled ( fn,
std::vector< bitLenInt bits,
const bitLenInt controls,
bitLenInt  controlLen 
)
protected

◆ CombineEngines() [1/2]

void Qrack::QPager::CombineEngines ( bitLenInt  thresholdBits)
protected

◆ CombineEngines() [2/2]

void Qrack::QPager::CombineEngines ( )
inlineprotected

◆ Compose() [1/2]

bitLenInt Qrack::QPager::Compose ( QPagerPtr  toCopy)
virtual

◆ Compose() [2/2]

virtual bitLenInt Qrack::QPager::Compose ( QInterfacePtr  toCopy)
inlinevirtual

Combine another QInterface with this one, after the last bit index of this one.

"Compose" combines the quantum description of state of two independent QInterface objects into one object, containing the full permutation basis of the full object. The "inputState" bits are added after the last qubit index of the QInterface to which we "Compose." Informally, "Compose" is equivalent to "just setting another group of qubits down next to the first" without interacting them. Schroedinger's equation can form a description of state for two independent subsystems at once or "separable quantum subsystems" without interacting them. Once the description of state of the independent systems is combined, we can interact them, and we can describe their entanglements to each other, in which case they are no longer independent. A full entangled description of quantum state is not possible for two independent quantum subsystems until we "Compose" them.

"Compose" multiplies the probabilities of the indepedent permutation states of the two subsystems to find the probabilites of the entire set of combined permutations, by simple combinatorial reasoning. If the probablity of the "left-hand" subsystem being in |00> is 1/4, and the probablity of the "right-hand" subsystem being in |101> is 1/8, than the probability of the combined |00101> permutation state is 1/32, and so on for all permutations of the new combined state.

If the programmer doesn't want to "cheat" quantum mechanically, then the original copy of the state which is duplicated into the larger QInterface should be "thrown away" to satisfy "no clone theorem." This is not semantically enforced in Qrack, because optimization of an emulator might be acheived by "cloning" "under-the-hood" while only exposing a quantum mechanically consistent API or instruction set.

Returns the quantum bit offset that the QInterface was appended at, such that bit 5 in toCopy is equal to offset+5 in this object.

Reimplemented from Qrack::QInterface.

◆ CPhaseFlipIfLess()

void Qrack::QPager::CPhaseFlipIfLess ( bitCapInt  greaterPerm,
bitLenInt  start,
bitLenInt  length,
bitLenInt  flagIndex 
)
virtual

The 6502 uses its carry flag also as a greater-than/less-than flag, for the CMP operation.

Implements Qrack::QAlu.

◆ CPOWModNOut()

void Qrack::QPager::CPOWModNOut ( bitCapInt  base,
bitCapInt  modN,
bitLenInt  inStart,
bitLenInt  outStart,
bitLenInt  length,
const bitLenInt controls,
bitLenInt  controlLen 
)
virtual

Controlled, raise a classical base to a quantum power, modulo N, (out of place)

Implements Qrack::QAlu.

◆ CUniformParityRZ()

void Qrack::QPager::CUniformParityRZ ( const bitLenInt controls,
bitLenInt  controlLen,
bitCapInt  mask,
real1_f  angle 
)
virtual

If the controls are set and the target qubit set parity is odd, this applies a phase factor of \(e^{i angle}\).

If the controls are set and the target qubit set parity is even, this applies the conjugate, \(e^{-i angle}\). Otherwise, do nothing if any control is not set.

Implements Qrack::QParity.

◆ DEC()

virtual void Qrack::QPager::DEC ( bitCapInt  toSub,
bitLenInt  start,
bitLenInt  length 
)
inlinevirtual

Add integer (without sign)

Subtract integer (without sign)

Implements Qrack::QAlu.

◆ DECC()

virtual void Qrack::QPager::DECC ( bitCapInt  toSub,
bitLenInt  start,
bitLenInt  length,
bitLenInt  carryIndex 
)
inlinevirtual

Subtract classical integer (without sign, with carry)

Subtract integer (without sign, with carry)

Reimplemented from Qrack::QAlu.

◆ Decompose() [1/3]

virtual void Qrack::QPager::Decompose ( bitLenInt  start,
QInterfacePtr  dest 
)
inlinevirtual

Minimally decompose a set of contiguous bits from the separably composed unit, into "destination".

Minimally decompose a set of contigious bits from the separably composed unit. The length of this separable unit is reduced by the length of bits decomposed, and the bits removed are output in the destination QInterface pointer. The destination object must be initialized to the correct number of bits, in 0 permutation state. For quantum mechanical accuracy, the bit set removed and the bit set left behind should be quantum mechanically "separable."

Like how "Compose" is like "just setting another group of qubits down next to the first," if two sets of qubits are not entangled, then "Decompose" is like "just moving a few qubits away from the rest." Schroedinger's equation does not require bits to be explicitly interacted in order to describe their permutation basis, and the descriptions of state of separable subsystems, those which are not entangled with other subsystems, are just as easily removed from the description of state. (This is equivalent to a "Schmidt decomposition.")

If we have for example 5 qubits, and we wish to separate into "left" and "right" subsystems of 3 and 2 qubits, we sum probabilities of one permutation of the "left" three over ALL permutations of the "right" two, for all permutations, and vice versa, like so:

\( P(|1000>|xy>) = P(|1000 00>) + P(|1000 10>) + P(|1000 01>) + P(|1000 11>). \)

If the subsystems are not "separable," i.e. if they are entangled, this operation is not well-motivated, and its output is not necessarily defined. (The summing of probabilities over permutations of subsytems will be performed as described above, but this is not quantum mechanically meaningful.) To ensure that the subsystem is "separable," i.e. that it has no entanglements to other subsystems in the QInterface, it can be measured with M(), or else all qubits other than the subsystem can be measured.

Implements Qrack::QInterface.

◆ Decompose() [2/3]

void Qrack::QPager::Decompose ( bitLenInt  start,
QPagerPtr  dest 
)
virtual

◆ Decompose() [3/3]

QInterfacePtr Qrack::QPager::Decompose ( bitLenInt  start,
bitLenInt  length 
)
virtual

Schmidt decompose a length of qubits.

Implements Qrack::QInterface.

◆ DECS()

virtual void Qrack::QPager::DECS ( bitCapInt  toSub,
bitLenInt  start,
bitLenInt  length,
bitLenInt  overflowIndex 
)
inlinevirtual

Add a classical integer to the register, with sign and without carry.

Subtract an integer from the register, with sign and without carry.

Because the register length is an arbitrary number of bits, the sign bit position on the integer to add is variable. Hence, the integer to add is specified as cast to an unsigned format, with the sign bit assumed to be set at the appropriate position before the cast.

Implements Qrack::QAlu.

◆ Dispose() [1/2]

void Qrack::QPager::Dispose ( bitLenInt  start,
bitLenInt  length 
)
virtual

Minimally decompose a set of contiguous bits from the separably composed unit, and discard the separable bits from index "start" for "length.".

Minimally decompose a set of contigious bits from the separably composed unit. The length of this separable unit is reduced by the length of bits decomposed, and the bits removed are output in the destination QInterface pointer. The destination object must be initialized to the correct number of bits, in 0 permutation state. For quantum mechanical accuracy, the bit set removed and the bit set left behind should be quantum mechanically "separable."

Like how "Compose" is like "just setting another group of qubits down next to the first," if two sets of qubits are not entangled, then "Decompose" is like "just moving a few qubits away from the rest." Schroedinger's equation does not require bits to be explicitly interacted in order to describe their permutation basis, and the descriptions of state of separable subsystems, those which are not entangled with other subsystems, are just as easily removed from the description of state. (This is equivalent to a "Schmidt decomposition.")

If we have for example 5 qubits, and we wish to separate into "left" and "right" subsystems of 3 and 2 qubits, we sum probabilities of one permutation of the "left" three over ALL permutations of the "right" two, for all permutations, and vice versa, like so:

\( P(|1000>|xy>) = P(|1000 00>) + P(|1000 10>) + P(|1000 01>) + P(|1000 11>). \)

If the subsystems are not "separable," i.e. if they are entangled, this operation is not well-motivated, and its output is not necessarily defined. (The summing of probabilities over permutations of subsytems will be performed as described above, but this is not quantum mechanically meaningful.) To ensure that the subsystem is "separable," i.e. that it has no entanglements to other subsystems in the QInterface, it can be measured with M(), or else all qubits other than the subsystem can be measured.

Implements Qrack::QInterface.

◆ Dispose() [2/2]

void Qrack::QPager::Dispose ( bitLenInt  start,
bitLenInt  length,
bitCapInt  disposedPerm 
)
virtual

Dispose a a contiguous set of qubits that are already in a permutation eigenstate.

Implements Qrack::QInterface.

◆ DIV()

void Qrack::QPager::DIV ( bitCapInt  toDiv,
bitLenInt  start,
bitLenInt  carryStart,
bitLenInt  length 
)
virtual

Divide by integer.

Implements Qrack::QAlu.

◆ Dump()

virtual void Qrack::QPager::Dump ( )
inlinevirtual

If asynchronous work is still running, let the simulator know that it can be aborted.

Note that this method is typically used internally where appropriate, such that user code typically does not call Dump().

Reimplemented from Qrack::QInterface.

◆ ExpectationBitsAll()

real1_f Qrack::QPager::ExpectationBitsAll ( const bitLenInt bits,
bitLenInt  length,
bitCapInt  offset = 0 
)
virtual

Get permutation expectation value of bits.

The permutation expectation value of all included bits is returned, with bits valued from low to high as the order of the "bits" array parameter argument.

Warning
PSEUDO-QUANTUM

Reimplemented from Qrack::QInterface.

◆ Finish()

virtual void Qrack::QPager::Finish ( )
inlinevirtual

If asynchronous work is still running, block until it finishes.

Note that this is never necessary to get correct, timely return values. QEngines and other layers will always internally "Finish" when necessary for correct return values. This is primarily for debugging and benchmarking.

Reimplemented from Qrack::QInterface.

◆ ForceM()

bool Qrack::QPager::ForceM ( bitLenInt  qubit,
bool  result,
bool  doForce = true,
bool  doApply = true 
)
virtual

Act as if is a measurement was applied, except force the (usually random) result.

Warning
PSEUDO-QUANTUM

Implements Qrack::QInterface.

◆ ForceMParity()

virtual bool Qrack::QPager::ForceMParity ( bitCapInt  mask,
bool  result,
bool  doForce = true 
)
inlinevirtual

Act as if is a measurement of parity of the masked set of qubits was applied, except force the (usually random) result.

Warning
PSEUDO-QUANTUM

Implements Qrack::QParity.

◆ FSim()

void Qrack::QPager::FSim ( real1_f  theta,
real1_f  phi,
bitLenInt  qubitIndex1,
bitLenInt  qubitIndex2 
)
virtual

The 2-qubit "fSim" gate, (useful in the simulation of particles with fermionic statistics)

Implements Qrack::QInterface.

◆ GetAmplitude()

virtual complex Qrack::QPager::GetAmplitude ( bitCapInt  perm)
inlinevirtual

Get the representational amplitude of a full permutation.

Warning
PSEUDO-QUANTUM

Implements Qrack::QInterface.

◆ GetDevice()

virtual int64_t Qrack::QPager::GetDevice ( )
inlinevirtual

Get the device index.

("-1" is default).

Reimplemented from Qrack::QInterface.

◆ GetMaxSize()

bitCapIntOcl Qrack::QPager::GetMaxSize ( )
inline

◆ GetProbs()

void Qrack::QPager::GetProbs ( real1 outputProbs)
virtual

Get the pure quantum state representation.

Warning
PSEUDO-QUANTUM

Implements Qrack::QInterface.

◆ GetQuantumState()

void Qrack::QPager::GetQuantumState ( complex outputState)
virtual

Get the pure quantum state representation.

Warning
PSEUDO-QUANTUM

Implements Qrack::QInterface.

◆ Hash()

void Qrack::QPager::Hash ( bitLenInt  start,
bitLenInt  length,
const unsigned char *  values 
)
virtual

Transform a length of qubit register via lookup through a hash table.

The hash table must be a one-to-one function, otherwise the behavior of this method is undefined. The value array definition convention is the same as IndexedLDA(). Essentially, this is an IndexedLDA() operation that replaces the index register with the value register, but the lookup table must therefore be one-to-one, for this operation to be unitary, as required.

Implements Qrack::QAlu.

◆ IMULModNOut()

void Qrack::QPager::IMULModNOut ( bitCapInt  toMul,
bitCapInt  modN,
bitLenInt  inStart,
bitLenInt  outStart,
bitLenInt  length 
)
virtual

Inverse of multiplication modulo N by integer, (out of place)

Implements Qrack::QAlu.

◆ INC()

virtual void Qrack::QPager::INC ( bitCapInt  toAdd,
bitLenInt  start,
bitLenInt  length 
)
inlinevirtual

Add integer (without sign)

Implements Qrack::QAlu.

◆ INCBCD()

void Qrack::QPager::INCBCD ( bitCapInt  toAdd,
bitLenInt  start,
bitLenInt  length 
)
virtual

Add classical BCD integer (without sign)

Implements Qrack::QAlu.

◆ INCC()

virtual void Qrack::QPager::INCC ( bitCapInt  toAdd,
bitLenInt  start,
bitLenInt  length,
bitLenInt  carryIndex 
)
inlinevirtual

Add integer (without sign, with carry)

Reimplemented from Qrack::QAlu.

◆ INCDECBCDC()

void Qrack::QPager::INCDECBCDC ( bitCapInt  toMod,
bitLenInt  start,
bitLenInt  length,
bitLenInt  carryIndex 
)
virtual

Common driver method behind INCSC and DECSC (without overflow flag)

Implements Qrack::QAlu.

◆ INCDECC()

virtual void Qrack::QPager::INCDECC ( bitCapInt  toMod,
bitLenInt  start,
bitLenInt  length,
bitLenInt  carryIndex 
)
inlinevirtual

Common driver method behind INCC and DECC (without sign, with carry)

Implements Qrack::QAlu.

◆ INCDECSC() [1/2]

void Qrack::QPager::INCDECSC ( bitCapInt  toMod,
bitLenInt  start,
bitLenInt  length,
bitLenInt  overflowIndex,
bitLenInt  carryIndex 
)
virtual

Common driver method behind INCSC and DECSC (with overflow flag)

Implements Qrack::QAlu.

◆ INCDECSC() [2/2]

void Qrack::QPager::INCDECSC ( bitCapInt  toMod,
bitLenInt  start,
bitLenInt  length,
bitLenInt  carryIndex 
)
virtual

Common driver method behind INCSC and DECSC (without overflow flag)

Implements Qrack::QAlu.

◆ INCS()

virtual void Qrack::QPager::INCS ( bitCapInt  toAdd,
bitLenInt  start,
bitLenInt  length,
bitLenInt  overflowIndex 
)
inlinevirtual

Add a classical integer to the register, with sign and without carry.

Implements Qrack::QAlu.

◆ IndexedADC()

bitCapInt Qrack::QPager::IndexedADC ( bitLenInt  indexStart,
bitLenInt  indexLength,
bitLenInt  valueStart,
bitLenInt  valueLength,
bitLenInt  carryIndex,
const unsigned char *  values 
)
virtual

Add to entangled 8 bit register state with a superposed index-offset-based read from classical memory.

"inputStart" is the start index of 8 qubits that act as an index into the 256 byte "values" array. The "outputStart" bits would usually already be entangled with the "inputStart" bits via a IndexedLDA() operation. With the "inputStart" bits being a "key" and the "outputStart" bits being a value, the permutation state |key, value> is mapped to |key, value + values[key]>. This is similar to classical parallel addition of two arrays. However, when either of the registers are measured, both registers will collapse into one random VALID key-value pair, with any addition or subtraction done to the "value." See IndexedLDA() for context.

FOR BEST EFFICIENCY, the "values" array should be allocated aligned to a 64-byte boundary. (See the unit tests suite code for an example of how to align the allocation.)

While a QInterface represents an interacting set of qubit-based registers, or a virtual quantum chip, the registers need to interact in some way with (classical or quantum) RAM. IndexedLDA is a RAM access method similar to the X addressing mode of the MOS 6502 chip, if the X register can be in a state of coherent superposition when it loads from RAM. "IndexedADC" and "IndexedSBC" perform add and subtract (with carry) operations on a state usually initially prepared with IndexedLDA().

Implements Qrack::QAlu.

◆ IndexedLDA()

bitCapInt Qrack::QPager::IndexedLDA ( bitLenInt  indexStart,
bitLenInt  indexLength,
bitLenInt  valueStart,
bitLenInt  valueLength,
const unsigned char *  values,
bool  resetValue = true 
)
virtual

Set 8 bit register bits by a superposed index-offset-based read from classical memory.

"inputStart" is the start index of 8 qubits that act as an index into the 256 byte "values" array. The "outputStart" bits are first cleared, then the separable |input, 00000000> permutation state is mapped to |input, values[input]>, with "values[input]" placed in the "outputStart" register. FOR BEST EFFICIENCY, the "values" array should be allocated aligned to a 64-byte boundary. (See the unit tests suite code for an example of how to align the allocation.)

While a QInterface represents an interacting set of qubit-based registers, or a virtual quantum chip, the registers need to interact in some way with (classical or quantum) RAM. IndexedLDA is a RAM access method similar to the X addressing mode of the MOS 6502 chip, if the X register can be in a state of coherent superposition when it loads from RAM.

The physical motivation for this addressing mode can be explained as follows: say that we have a superconducting quantum interface device (SQUID) based chip. SQUIDs have already been demonstrated passing coherently superposed electrical currents. In a sufficiently quantum-mechanically isolated qubit chip with a classical cache, with both classical RAM and registers likely cryogenically isolated from the environment, SQUIDs could (hopefully) pass coherently superposed electrical currents into the classical RAM cache to load values into a qubit register. The state loaded would be a superposition of the values of all RAM to which coherently superposed electrical currents were passed.

In qubit system similar to the MOS 6502, say we have qubit-based "accumulator" and "X index" registers, and say that we start with a superposed X index register. In (classical) X addressing mode, the X index register value acts an offset into RAM from a specified starting address. The X addressing mode of a LoaD Accumulator (LDA) instruction, by the physical mechanism described above, should load the accumulator in quantum parallel with the values of every different address of RAM pointed to in superposition by the X index register. The superposed values in the accumulator are entangled with those in the X index register, by way of whatever values the classical RAM pointed to by X held at the time of the load. (If the RAM at index "36" held an unsigned char value of "27," then the value "36" in the X index register becomes entangled with the value "27" in the accumulator, and so on in quantum parallel for all superposed values of the X index register, at once.) If the X index register or accumulator are then measured, the two registers will both always collapse into a random but valid key-value pair of X index offset and value at that classical RAM address.

Note that a "superposed store operation in classical RAM" is not possible by analagous reasoning. Classical RAM would become entangled with both the accumulator and the X register. When the state of the registers was collapsed, we would find that only one "store" operation to a single memory address had actually been carried out, consistent with the address offset in the collapsed X register and the byte value in the collapsed accumulator. It would not be possible by this model to write in quantum parallel to more than one address of classical memory at a time.

Implements Qrack::QAlu.

◆ IndexedSBC()

bitCapInt Qrack::QPager::IndexedSBC ( bitLenInt  indexStart,
bitLenInt  indexLength,
bitLenInt  valueStart,
bitLenInt  valueLength,
bitLenInt  carryIndex,
const unsigned char *  values 
)
virtual

Subtract from an entangled 8 bit register state with a superposed index-offset-based read from classical memory.

"inputStart" is the start index of 8 qubits that act as an index into the 256 byte "values" array. The "outputStart" bits would usually already be entangled with the "inputStart" bits via a IndexedLDA() operation. With the "inputStart" bits being a "key" and the "outputStart" bits being a value, the permutation state |key, value> is mapped to |key, value - values[key]>. This is similar to classical parallel addition of two arrays. However, when either of the registers are measured, both registers will collapse into one random VALID key-value pair, with any addition or subtraction done to the "value." See QInterface::IndexedLDA for context.

FOR BEST EFFICIENCY, the "values" array should be allocated aligned to a 64-byte boundary. (See the unit tests suite code for an example of how to align the allocation.)

While a QInterface represents an interacting set of qubit-based registers, or a virtual quantum chip, the registers need to interact in some way with (classical or quantum) RAM. IndexedLDA is a RAM access method similar to the X addressing mode of the MOS 6502 chip, if the X register can be in a state of coherent superposition when it loads from RAM. "IndexedADC" and "IndexedSBC" perform add and subtract (with carry) operations on a state usually initially prepared with IndexedLDA().

Implements Qrack::QAlu.

◆ Init()

void Qrack::QPager::Init ( )
protected

◆ Invert()

virtual void Qrack::QPager::Invert ( complex  topRight,
complex  bottomLeft,
bitLenInt  qubitIndex 
)
inlinevirtual

Apply a single bit transformation that reverses bit probability and might effect phase.

Reimplemented from Qrack::QInterface.

◆ isFinished()

virtual bool Qrack::QPager::isFinished ( )
inlinevirtual

Returns "false" if asynchronous work is still running, and "true" if all previously dispatched asynchronous work is done.

Reimplemented from Qrack::QInterface.

◆ ISqrtSwap()

void Qrack::QPager::ISqrtSwap ( bitLenInt  qubitIndex1,
bitLenInt  qubitIndex2 
)
virtual

Inverse square root of Swap gate.

Reimplemented from Qrack::QInterface.

◆ ISwap()

void Qrack::QPager::ISwap ( bitLenInt  qubitIndex1,
bitLenInt  qubitIndex2 
)
virtual

Swap values of two bits in register, and apply phase factor of i if bits are different.

Reimplemented from Qrack::QInterface.

◆ LockEngine()

virtual void Qrack::QPager::LockEngine ( QEnginePtr  eng)
inlinevirtual

◆ M()

virtual bool Qrack::QPager::M ( bitLenInt  q)
inlinevirtual

Implements Qrack::QAlu.

◆ MACMtrx()

virtual void Qrack::QPager::MACMtrx ( const bitLenInt controls,
bitLenInt  controlLen,
const complex mtrx,
bitLenInt  target 
)
inlinevirtual

Apply an arbitrary single bit unitary transformation, with arbitrary (anti-)control bits.

Reimplemented from Qrack::QInterface.

◆ MakeEngine()

QEnginePtr Qrack::QPager::MakeEngine ( bitLenInt  length,
bitCapInt  perm,
int  deviceId 
)
protected

◆ MCMtrx()

virtual void Qrack::QPager::MCMtrx ( const bitLenInt controls,
bitLenInt  controlLen,
const complex mtrx,
bitLenInt  target 
)
inlinevirtual

Apply an arbitrary single bit unitary transformation, with arbitrary control bits.

Implements Qrack::QInterface.

◆ MetaControlled()

template<typename Qubit1Fn >
void Qrack::QPager::MetaControlled ( bool  anti,
const std::vector< bitLenInt > &  controls,
bitLenInt  target,
Qubit1Fn  fn,
const complex mtrx,
bool  isSqiCtrl = false,
bool  isIntraCtrled = false 
)
protected

◆ MetaSwap()

void Qrack::QPager::MetaSwap ( bitLenInt  qubit1,
bitLenInt  qubit2,
bool  isIPhaseFac 
)
protected

◆ Mtrx()

void Qrack::QPager::Mtrx ( const complex mtrx,
bitLenInt  qubitIndex 
)
virtual

Apply an arbitrary single bit unitary transformation.

Implements Qrack::QInterface.

◆ MUL()

void Qrack::QPager::MUL ( bitCapInt  toMul,
bitLenInt  start,
bitLenInt  carryStart,
bitLenInt  length 
)
virtual

Multiply by integer.

Implements Qrack::QAlu.

◆ MULModNOut()

void Qrack::QPager::MULModNOut ( bitCapInt  toMul,
bitCapInt  modN,
bitLenInt  inStart,
bitLenInt  outStart,
bitLenInt  length 
)
virtual

Multiplication modulo N by integer, (out of place)

Implements Qrack::QAlu.

◆ NormalizeState()

void Qrack::QPager::NormalizeState ( real1_f  nrm = REAL1_DEFAULT_ARG,
real1_f  norm_thresh = REAL1_DEFAULT_ARG,
real1_f  phaseArg = ZERO_R1 
)
virtual

Apply the normalization factor found by UpdateRunningNorm() or on the fly by a single bit gate.

(On an actual quantum computer, the state should never require manual normalization.)

Warning
PSEUDO-QUANTUM

Implements Qrack::QInterface.

◆ pagedQubitCount()

bitLenInt Qrack::QPager::pagedQubitCount ( )
inlineprotected

◆ pageMaxQPower()

bitCapInt Qrack::QPager::pageMaxQPower ( )
inlineprotected

◆ Phase()

virtual void Qrack::QPager::Phase ( complex  topLeft,
complex  bottomRight,
bitLenInt  qubitIndex 
)
inlinevirtual

Apply a single bit transformation that only effects phase.

Reimplemented from Qrack::QInterface.

◆ PhaseFlipIfLess()

void Qrack::QPager::PhaseFlipIfLess ( bitCapInt  greaterPerm,
bitLenInt  start,
bitLenInt  length 
)
virtual

This is an expedient for an adaptive Grover's search for a function's global minimum.

Implements Qrack::QAlu.

◆ PhaseParity()

void Qrack::QPager::PhaseParity ( real1_f  radians,
bitCapInt  mask 
)
virtual

Parity phase gate.

Applies e^(i*angle) phase factor to all combinations of bits with odd parity, based upon permutations of qubits.

Reimplemented from Qrack::QInterface.

◆ POWModNOut()

void Qrack::QPager::POWModNOut ( bitCapInt  base,
bitCapInt  modN,
bitLenInt  inStart,
bitLenInt  outStart,
bitLenInt  length 
)
virtual

Raise a classical base to a quantum power, modulo N, (out of place)

Implements Qrack::QAlu.

◆ Prob()

real1_f Qrack::QPager::Prob ( bitLenInt  qubitIndex)
virtual

Direct measure of bit probability to be in |1> state.

Warning
PSEUDO-QUANTUM

Implements Qrack::QInterface.

◆ ProbAll()

real1_f Qrack::QPager::ProbAll ( bitCapInt  fullRegister)
inlinevirtual

Direct measure of full permutation probability.

Warning
PSEUDO-QUANTUM

Reimplemented from Qrack::QInterface.

◆ ProbMask()

real1_f Qrack::QPager::ProbMask ( bitCapInt  mask,
bitCapInt  permutation 
)
virtual

Direct measure of masked permutation probability.

Returns probability of permutation of the mask.

"mask" masks the bits to check the probability of. "permutation" sets the 0 or 1 value for each bit in the mask. Bits which are set in the mask can be set to 0 or 1 in the permutation, while reset bits in the mask should be 0 in the permutation.

Warning
PSEUDO-QUANTUM

Reimplemented from Qrack::QInterface.

◆ ProbParity()

virtual real1_f Qrack::QPager::ProbParity ( bitCapInt  mask)
inlinevirtual

Overall probability of any odd permutation of the masked set of bits.

Implements Qrack::QParity.

◆ qubitsPerPage()

bitLenInt Qrack::QPager::qubitsPerPage ( )
inlineprotected

◆ ReleaseEngine()

virtual QEnginePtr Qrack::QPager::ReleaseEngine ( )
inlinevirtual

◆ SemiMetaControlled()

template<typename Qubit1Fn >
void Qrack::QPager::SemiMetaControlled ( bool  anti,
std::vector< bitLenInt controls,
bitLenInt  target,
Qubit1Fn  fn 
)
protected

◆ SemiMetaSwap()

void Qrack::QPager::SemiMetaSwap ( bitLenInt  qubit1,
bitLenInt  qubit2,
bool  isIPhaseFac 
)
protected

◆ SeparateEngines() [1/2]

void Qrack::QPager::SeparateEngines ( bitLenInt  thresholdBits,
bool  noBaseFloor = false 
)
protected

◆ SeparateEngines() [2/2]

void Qrack::QPager::SeparateEngines ( )
inlineprotected

◆ SetAmplitude()

virtual void Qrack::QPager::SetAmplitude ( bitCapInt  perm,
complex  amp 
)
inlinevirtual

Sets the representational amplitude of a full permutation.

Warning
PSEUDO-QUANTUM

Implements Qrack::QInterface.

◆ SetConcurrency()

virtual void Qrack::QPager::SetConcurrency ( uint32_t  threadsPerEngine)
inlinevirtual

Set the number of threads in parallel for loops, per component QEngine.

Reimplemented from Qrack::QInterface.

◆ SetDevice()

virtual void Qrack::QPager::SetDevice ( int  dID,
bool  forceReInit = false 
)
inlinevirtual

Set the device index, if more than one device is available.

Reimplemented from Qrack::QInterface.

◆ SetPermutation()

void Qrack::QPager::SetPermutation ( bitCapInt  perm,
complex  phaseFac = CMPLX_DEFAULT_ARG 
)
virtual

Set to a specific permutation of all qubits.

Reimplemented from Qrack::QInterface.

◆ SetQuantumState()

void Qrack::QPager::SetQuantumState ( const complex inputState)
virtual

Set an arbitrary pure quantum state representation.

Warning
PSEUDO-QUANTUM

Implements Qrack::QInterface.

◆ SetQubitCount()

virtual void Qrack::QPager::SetQubitCount ( bitLenInt  qb)
inlineprotectedvirtual

Reimplemented from Qrack::QInterface.

◆ SingleBitGate()

template<typename Qubit1Fn >
void Qrack::QPager::SingleBitGate ( bitLenInt  target,
Qubit1Fn  fn,
bool  isSqiCtrl = false,
bool  isAnti = false 
)
protected

◆ SqrtSwap()

void Qrack::QPager::SqrtSwap ( bitLenInt  qubitIndex1,
bitLenInt  qubitIndex2 
)
virtual

Square root of Swap gate.

Reimplemented from Qrack::QInterface.

◆ SumSqrDiff() [1/2]

virtual real1_f Qrack::QPager::SumSqrDiff ( QInterfacePtr  toCompare)
inlinevirtual

Implements Qrack::QInterface.

◆ SumSqrDiff() [2/2]

real1_f Qrack::QPager::SumSqrDiff ( QPagerPtr  toCompare)
virtual

◆ Swap()

void Qrack::QPager::Swap ( bitLenInt  qubitIndex1,
bitLenInt  qubitIndex2 
)
virtual

Swap values of two bits in register.

Reimplemented from Qrack::QInterface.

◆ UniformParityRZ()

void Qrack::QPager::UniformParityRZ ( bitCapInt  mask,
real1_f  angle 
)
virtual

If the target qubit set parity is odd, this applies a phase factor of \(e^{i angle}\).

If the target qubit set parity is even, this applies the conjugate, e^{-i angle}.

Reimplemented from Qrack::QParity.

◆ UpdateRunningNorm()

void Qrack::QPager::UpdateRunningNorm ( real1_f  norm_thresh = REAL1_DEFAULT_ARG)
virtual

Force a calculation of the norm of the state vector, in order to make it unit length before the next probability or measurement operation.

(On an actual quantum computer, the state should never require manual normalization.)

Warning
PSEUDO-QUANTUM

Implements Qrack::QInterface.

◆ X()

virtual void Qrack::QPager::X ( bitLenInt  q)
inlinevirtual

Implements Qrack::QAlu.

◆ XMask()

void Qrack::QPager::XMask ( bitCapInt  mask)
virtual

Masked X gate.

Applies the Pauli "X" operator to all qubits in the mask. A qubit index "n" is in the mask if (((1 << n) & mask)

0). The Pauli "X" operator is equivalent to a logical "NOT."

Reimplemented from Qrack::QInterface.

◆ ZMask()

virtual void Qrack::QPager::ZMask ( bitCapInt  mask)
inlinevirtual

Masked Z gate.

Applies the Pauli "Z" operator to all qubits in the mask. A qubit index "n" is in the mask if (((1 << n) & mask)

0). The Pauli "Z" operator reverses the phase of |1> and leaves |0> unchanged.

Reimplemented from Qrack::QInterface.

Member Data Documentation

◆ basePageCount

bitCapInt Qrack::QPager::basePageCount
protected

◆ basePageMaxQPower

bitCapIntOcl Qrack::QPager::basePageMaxQPower
protected

◆ baseQubitsPerPage

bitLenInt Qrack::QPager::baseQubitsPerPage
protected

◆ deviceIDs

std::vector<int> Qrack::QPager::deviceIDs
protected

◆ devID

int Qrack::QPager::devID
protected

◆ engines

std::vector<QInterfaceEngine> Qrack::QPager::engines
protected

◆ isSparse

bool Qrack::QPager::isSparse
protected

◆ maxPageQubits

bitLenInt Qrack::QPager::maxPageQubits
protected

◆ maxQubits

bitLenInt Qrack::QPager::maxQubits
protected

◆ minPageQubits

bitLenInt Qrack::QPager::minPageQubits
protected

◆ phaseFactor

complex Qrack::QPager::phaseFactor
protected

◆ qPages

std::vector<QEnginePtr> Qrack::QPager::qPages
protected

◆ rootEngine

QInterfaceEngine Qrack::QPager::rootEngine
protected

◆ segmentGlobalQb

bitLenInt Qrack::QPager::segmentGlobalQb
protected

◆ thresholdQubitsPerPage

bitLenInt Qrack::QPager::thresholdQubitsPerPage
protected

◆ useGpuThreshold

bool Qrack::QPager::useGpuThreshold
protected

◆ useHardwareThreshold

bool Qrack::QPager::useHardwareThreshold
protected

◆ useHostRam

bool Qrack::QPager::useHostRam
protected

◆ useRDRAND

bool Qrack::QPager::useRDRAND
protected

The documentation for this class was generated from the following files: