#include "QuEST.h"
#include "QuEST_internal.h"
#include "QuEST_precision.h"
#include "mt19937ar.h"
#include "QuEST_cpu_internal.h"
#include <math.h>
#include <stdio.h>
#include <stdlib.h>
#include <stdint.h>
#include <assert.h>
Go to the source code of this file.
Functions | |
DiagonalOp | agnostic_createDiagonalOp (int numQubits, QuESTEnv env) |
void | agnostic_destroyDiagonalOp (DiagonalOp op) |
void | agnostic_setDiagonalOpElems (DiagonalOp op, long long int startInd, qreal *real, qreal *imag, long long int numElems) |
void | agnostic_syncDiagonalOp (DiagonalOp op) |
void | alternateNormZeroingSomeAmpBlocks (Qureg qureg, qreal norm, int normFirst, long long int startAmpInd, long long int numAmps, long long int blockSize) |
void | copyStateFromGPU (Qureg qureg) |
In GPU mode, this copies the state-vector (or density matrix) from GPU memory (qureg.deviceStateVec) to RAM (qureg.stateVec), where it can be accessed/modified by the user. More... | |
void | copyStateToGPU (Qureg qureg) |
In GPU mode, this copies the state-vector (or density matrix) from RAM (qureg.stateVec) to VRAM / GPU-memory (qureg.deviceStateVec), which is the version operated upon by other calls to the API. More... | |
void | densmatr_applyDiagonalOpLocal (Qureg qureg, DiagonalOp op) |
Complex | densmatr_calcExpecDiagonalOpLocal (Qureg qureg, DiagonalOp op) |
qreal | densmatr_calcFidelityLocal (Qureg qureg, Qureg pureState) |
computes a few dens-columns-worth of (vec^*T) dens * vec More... | |
qreal | densmatr_calcHilbertSchmidtDistanceSquaredLocal (Qureg a, Qureg b) |
computes Tr((a-b) conjTrans(a-b)) = sum of abs values of (a-b) More... | |
qreal | densmatr_calcInnerProductLocal (Qureg a, Qureg b) |
computes Tr(conjTrans(a) b) = sum of (a_ij^* b_ij) More... | |
qreal | densmatr_calcPurityLocal (Qureg qureg) |
void | densmatr_collapseToKnownProbOutcome (Qureg qureg, int measureQubit, int outcome, qreal totalStateProb) |
Renorms (/prob) every | * outcome * >< * outcome * | state, setting all others to zero. More... | |
qreal | densmatr_findProbabilityOfZeroLocal (Qureg qureg, int measureQubit) |
void | densmatr_initClassicalState (Qureg qureg, long long int stateInd) |
void | densmatr_initPlusState (Qureg qureg) |
void | densmatr_initPureStateLocal (Qureg targetQureg, Qureg copyQureg) |
void | densmatr_mixDampingDistributed (Qureg qureg, int targetQubit, qreal damping) |
void | densmatr_mixDampingLocal (Qureg qureg, int targetQubit, qreal damping) |
void | densmatr_mixDensityMatrix (Qureg combineQureg, qreal otherProb, Qureg otherQureg) |
void | densmatr_mixDephasing (Qureg qureg, int targetQubit, qreal dephase) |
void | densmatr_mixDepolarisingDistributed (Qureg qureg, int targetQubit, qreal depolLevel) |
void | densmatr_mixDepolarisingLocal (Qureg qureg, int targetQubit, qreal depolLevel) |
void | densmatr_mixTwoQubitDephasing (Qureg qureg, int qubit1, int qubit2, qreal dephase) |
void | densmatr_mixTwoQubitDepolarisingDistributed (Qureg qureg, int targetQubit, int qubit2, qreal delta, qreal gamma) |
void | densmatr_mixTwoQubitDepolarisingLocal (Qureg qureg, int qubit1, int qubit2, qreal delta, qreal gamma) |
void | densmatr_mixTwoQubitDepolarisingLocalPart1 (Qureg qureg, int qubit1, int qubit2, qreal delta) |
void | densmatr_mixTwoQubitDepolarisingQ1LocalQ2DistributedPart3 (Qureg qureg, int targetQubit, int qubit2, qreal delta, qreal gamma) |
void | densmatr_oneQubitDegradeOffDiagonal (Qureg qureg, int targetQubit, qreal retain) |
int | getBitMaskParity (long long int mask) |
void | normaliseSomeAmps (Qureg qureg, qreal norm, long long int startInd, long long int numAmps) |
int | qsortComp (const void *a, const void *b) |
void | statevec_applyDiagonalOp (Qureg qureg, DiagonalOp op) |
Complex | statevec_calcExpecDiagonalOpLocal (Qureg qureg, DiagonalOp op) |
Complex | statevec_calcInnerProductLocal (Qureg bra, Qureg ket) |
void | statevec_cloneQureg (Qureg targetQureg, Qureg copyQureg) |
void | statevec_collapseToKnownProbOutcomeDistributedRenorm (Qureg qureg, int measureQubit, qreal totalProbability) |
Renormalise parts of the state vector where measureQubit=0 or 1, based on the total probability of that qubit being in state 0 or 1. More... | |
void | statevec_collapseToKnownProbOutcomeLocal (Qureg qureg, int measureQubit, int outcome, qreal totalProbability) |
Update the state vector to be consistent with measuring measureQubit=0 if outcome=0 and measureQubit=1 if outcome=1. More... | |
void | statevec_collapseToOutcomeDistributedSetZero (Qureg qureg) |
Set all amplitudes in one chunk to 0. More... | |
void | statevec_compactUnitaryDistributed (Qureg qureg, Complex rot1, Complex rot2, ComplexArray stateVecUp, ComplexArray stateVecLo, ComplexArray stateVecOut) |
Rotate a single qubit in the state vector of probability amplitudes, given two complex numbers alpha and beta, and a subset of the state vector with upper and lower block values stored seperately. More... | |
void | statevec_compactUnitaryLocal (Qureg qureg, int targetQubit, Complex alpha, Complex beta) |
int | statevec_compareStates (Qureg mq1, Qureg mq2, qreal precision) |
void | statevec_controlledCompactUnitaryDistributed (Qureg qureg, int controlQubit, Complex rot1, Complex rot2, ComplexArray stateVecUp, ComplexArray stateVecLo, ComplexArray stateVecOut) |
Rotate a single qubit in the state vector of probability amplitudes, given two complex numbers alpha and beta and a subset of the state vector with upper and lower block values stored seperately. More... | |
void | statevec_controlledCompactUnitaryLocal (Qureg qureg, int controlQubit, int targetQubit, Complex alpha, Complex beta) |
void | statevec_controlledNotDistributed (Qureg qureg, int controlQubit, ComplexArray stateVecIn, ComplexArray stateVecOut) |
Rotate a single qubit by {{0,1},{1,0}. More... | |
void | statevec_controlledNotLocal (Qureg qureg, int controlQubit, int targetQubit) |
void | statevec_controlledPauliYDistributed (Qureg qureg, int controlQubit, ComplexArray stateVecIn, ComplexArray stateVecOut, int conjFac) |
void | statevec_controlledPauliYLocal (Qureg qureg, int controlQubit, int targetQubit, int conjFac) |
void | statevec_controlledPhaseFlip (Qureg qureg, int idQubit1, int idQubit2) |
void | statevec_controlledPhaseShift (Qureg qureg, int idQubit1, int idQubit2, qreal angle) |
void | statevec_controlledUnitaryDistributed (Qureg qureg, int controlQubit, Complex rot1, Complex rot2, ComplexArray stateVecUp, ComplexArray stateVecLo, ComplexArray stateVecOut) |
Rotate a single qubit in the state vector of probability amplitudes, given two complex numbers alpha and beta and a subset of the state vector with upper and lower block values stored seperately. More... | |
void | statevec_controlledUnitaryLocal (Qureg qureg, int controlQubit, int targetQubit, ComplexMatrix2 u) |
void | statevec_createQureg (Qureg *qureg, int numQubits, QuESTEnv env) |
void | statevec_destroyQureg (Qureg qureg, QuESTEnv env) |
qreal | statevec_findProbabilityOfZeroDistributed (Qureg qureg) |
Measure the probability of a specified qubit being in the zero state across all amplitudes held in this chunk. More... | |
qreal | statevec_findProbabilityOfZeroLocal (Qureg qureg, int measureQubit) |
Measure the total probability of a specified qubit being in the zero state across all amplitudes in this chunk. More... | |
void | statevec_getEnvironmentString (QuESTEnv env, Qureg qureg, char str[200]) |
void | statevec_hadamardDistributed (Qureg qureg, ComplexArray stateVecUp, ComplexArray stateVecLo, ComplexArray stateVecOut, int updateUpper) |
Rotate a single qubit by {{1,1},{1,-1}}/sqrt2. More... | |
void | statevec_hadamardLocal (Qureg qureg, int targetQubit) |
void | statevec_initBlankState (Qureg qureg) |
void | statevec_initClassicalState (Qureg qureg, long long int stateInd) |
void | statevec_initDebugState (Qureg qureg) |
Initialise the state vector of probability amplitudes to an (unphysical) state with each component of each probability amplitude a unique floating point value. More... | |
void | statevec_initPlusState (Qureg qureg) |
int | statevec_initStateFromSingleFile (Qureg *qureg, char filename[200], QuESTEnv env) |
void | statevec_initStateOfSingleQubit (Qureg *qureg, int qubitId, int outcome) |
Initialise the state vector of probability amplitudes such that one qubit is set to 'outcome' and all other qubits are in an equal superposition of zero and one. More... | |
void | statevec_initZeroState (Qureg qureg) |
void | statevec_multiControlledMultiQubitUnitaryLocal (Qureg qureg, long long int ctrlMask, int *targs, int numTargs, ComplexMatrixN u) |
void | statevec_multiControlledPhaseFlip (Qureg qureg, int *controlQubits, int numControlQubits) |
void | statevec_multiControlledPhaseShift (Qureg qureg, int *controlQubits, int numControlQubits, qreal angle) |
void | statevec_multiControlledTwoQubitUnitaryLocal (Qureg qureg, long long int ctrlMask, int q1, int q2, ComplexMatrix4 u) |
void | statevec_multiControlledUnitaryDistributed (Qureg qureg, int targetQubit, long long int ctrlQubitsMask, long long int ctrlFlipMask, Complex rot1, Complex rot2, ComplexArray stateVecUp, ComplexArray stateVecLo, ComplexArray stateVecOut) |
Apply a unitary operation to a single qubit in the state vector of probability amplitudes, given a subset of the state vector with upper and lower block values stored seperately. More... | |
void | statevec_multiControlledUnitaryLocal (Qureg qureg, int targetQubit, long long int ctrlQubitsMask, long long int ctrlFlipMask, ComplexMatrix2 u) |
void | statevec_multiRotateZ (Qureg qureg, long long int mask, qreal angle) |
void | statevec_pauliXDistributed (Qureg qureg, ComplexArray stateVecIn, ComplexArray stateVecOut) |
Rotate a single qubit by {{0,1},{1,0}. More... | |
void | statevec_pauliXLocal (Qureg qureg, int targetQubit) |
void | statevec_pauliYDistributed (Qureg qureg, ComplexArray stateVecIn, ComplexArray stateVecOut, int updateUpper, int conjFac) |
Rotate a single qubit by +-{{0,-i},{i,0}. More... | |
void | statevec_pauliYLocal (Qureg qureg, int targetQubit, int conjFac) |
void | statevec_phaseShiftByTerm (Qureg qureg, int targetQubit, Complex term) |
void | statevec_reportStateToScreen (Qureg qureg, QuESTEnv env, int reportRank) |
void | statevec_setAmps (Qureg qureg, long long int startInd, qreal *reals, qreal *imags, long long int numAmps) |
void | statevec_setWeightedQureg (Complex fac1, Qureg qureg1, Complex fac2, Qureg qureg2, Complex facOut, Qureg out) |
void | statevec_swapQubitAmpsDistributed (Qureg qureg, int pairRank, int qb1, int qb2) |
qureg.pairStateVec contains the entire set of amplitudes of the paired node which includes the set of all amplitudes which need to be swapped between |..0..1..> and |..1..0..> More... | |
void | statevec_swapQubitAmpsLocal (Qureg qureg, int qb1, int qb2) |
It is ensured that all amplitudes needing to be swapped are on this node. More... | |
void | statevec_unitaryDistributed (Qureg qureg, Complex rot1, Complex rot2, ComplexArray stateVecUp, ComplexArray stateVecLo, ComplexArray stateVecOut) |
Apply a unitary operation to a single qubit given a subset of the state vector with upper and lower block values stored seperately. More... | |
void | statevec_unitaryLocal (Qureg qureg, int targetQubit, ComplexMatrix2 u) |
void | zeroSomeAmps (Qureg qureg, long long int startInd, long long int numAmps) |
Detailed Description
The core of the CPU backend functionality. The CPU/MPI implementations of the pure state functions in ../QuEST_ops_pure.h are in QuEST_cpu_local.c and QuEST_cpu_distributed.c which mostly wrap the core functions defined here. Some additional hardware-agnostic functions are defined here
Definition in file QuEST_cpu.c.
Function Documentation
◆ agnostic_createDiagonalOp()
DiagonalOp agnostic_createDiagonalOp | ( | int | numQubits, |
QuESTEnv | env | ||
) |
Definition at line 1335 of file QuEST_cpu.c.
References DiagonalOp::chunkId, DiagonalOp::imag, DiagonalOp::numChunks, DiagonalOp::numElemsPerChunk, DiagonalOp::numQubits, QuESTEnv::numRanks, qreal, QuESTEnv::rank, and DiagonalOp::real.
◆ agnostic_destroyDiagonalOp()
void agnostic_destroyDiagonalOp | ( | DiagonalOp | op | ) |
Definition at line 1357 of file QuEST_cpu.c.
References DiagonalOp::imag, and DiagonalOp::real.
◆ agnostic_setDiagonalOpElems()
void agnostic_setDiagonalOpElems | ( | DiagonalOp | op, |
long long int | startInd, | ||
qreal * | real, | ||
qreal * | imag, | ||
long long int | numElems | ||
) |
Definition at line 3842 of file QuEST_cpu.c.
References DiagonalOp::chunkId, DiagonalOp::imag, DiagonalOp::numElemsPerChunk, qreal, and DiagonalOp::real.
◆ agnostic_syncDiagonalOp()
void agnostic_syncDiagonalOp | ( | DiagonalOp | op | ) |
Definition at line 1362 of file QuEST_cpu.c.
◆ alternateNormZeroingSomeAmpBlocks()
void alternateNormZeroingSomeAmpBlocks | ( | Qureg | qureg, |
qreal | norm, | ||
int | normFirst, | ||
long long int | startAmpInd, | ||
long long int | numAmps, | ||
long long int | blockSize | ||
) |
Definition at line 754 of file QuEST_cpu.c.
References normaliseSomeAmps(), and zeroSomeAmps().
Referenced by densmatr_collapseToKnownProbOutcome().
◆ densmatr_applyDiagonalOpLocal()
void densmatr_applyDiagonalOpLocal | ( | Qureg | qureg, |
DiagonalOp | op | ||
) |
Definition at line 3696 of file QuEST_cpu.c.
References Qureg::numAmpsPerChunk, DiagonalOp::numQubits, Qureg::pairStateVec, qreal, and Qureg::stateVec.
Referenced by densmatr_applyDiagonalOp().
◆ densmatr_calcExpecDiagonalOpLocal()
Complex densmatr_calcExpecDiagonalOpLocal | ( | Qureg | qureg, |
DiagonalOp | op | ||
) |
Definition at line 3781 of file QuEST_cpu.c.
References Qureg::chunkId, Complex::imag, DiagonalOp::imag, Qureg::numAmpsPerChunk, Qureg::numQubitsRepresented, qreal, Complex::real, DiagonalOp::real, and Qureg::stateVec.
Referenced by densmatr_calcExpecDiagonalOp().
◆ densmatr_calcFidelityLocal()
computes a few dens-columns-worth of (vec^*T) dens * vec
Definition at line 990 of file QuEST_cpu.c.
References Qureg::chunkId, Qureg::numAmpsPerChunk, Qureg::numAmpsTotal, Qureg::pairStateVec, qreal, and Qureg::stateVec.
Referenced by densmatr_calcFidelity().
◆ densmatr_calcHilbertSchmidtDistanceSquaredLocal()
computes Tr((a-b) conjTrans(a-b)) = sum of abs values of (a-b)
Definition at line 923 of file QuEST_cpu.c.
References Qureg::numAmpsPerChunk, qreal, and Qureg::stateVec.
Referenced by densmatr_calcHilbertSchmidtDistance().
◆ densmatr_calcInnerProductLocal()
computes Tr(conjTrans(a) b) = sum of (a_ij^* b_ij)
Definition at line 958 of file QuEST_cpu.c.
References Qureg::numAmpsPerChunk, qreal, and Qureg::stateVec.
Referenced by densmatr_calcInnerProduct().
◆ densmatr_calcPurityLocal()
Definition at line 861 of file QuEST_cpu.c.
References Qureg::numAmpsPerChunk, qreal, and Qureg::stateVec.
Referenced by densmatr_calcPurity().
◆ densmatr_collapseToKnownProbOutcome()
void densmatr_collapseToKnownProbOutcome | ( | Qureg | qureg, |
int | measureQubit, | ||
int | outcome, | ||
qreal | totalStateProb | ||
) |
Renorms (/prob) every | * outcome * >< * outcome * | state, setting all others to zero.
Definition at line 785 of file QuEST_cpu.c.
References alternateNormZeroingSomeAmpBlocks(), Qureg::chunkId, extractBit(), normaliseSomeAmps(), Qureg::numAmpsPerChunk, Qureg::numQubitsRepresented, and zeroSomeAmps().
◆ densmatr_findProbabilityOfZeroLocal()
Definition at line 3151 of file QuEST_cpu.c.
References Qureg::chunkId, extractBit(), Qureg::numAmpsPerChunk, Qureg::numQubitsRepresented, qreal, and Qureg::stateVec.
Referenced by densmatr_calcProbOfOutcome().
◆ densmatr_initClassicalState()
void densmatr_initClassicalState | ( | Qureg | qureg, |
long long int | stateInd | ||
) |
Definition at line 1115 of file QuEST_cpu.c.
References Qureg::chunkId, Qureg::numAmpsPerChunk, Qureg::numQubitsRepresented, qreal, and Qureg::stateVec.
◆ densmatr_initPlusState()
void densmatr_initPlusState | ( | Qureg | qureg | ) |
Definition at line 1154 of file QuEST_cpu.c.
References Qureg::numAmpsPerChunk, Qureg::numQubitsRepresented, qreal, and Qureg::stateVec.
◆ densmatr_initPureStateLocal()
Definition at line 1184 of file QuEST_cpu.c.
References Qureg::chunkId, Qureg::numAmpsPerChunk, Qureg::numAmpsTotal, Qureg::pairStateVec, qreal, and Qureg::stateVec.
Referenced by densmatr_initPureState().
◆ densmatr_mixDampingDistributed()
Definition at line 300 of file QuEST_cpu.c.
References Qureg::chunkId, densmatr_oneQubitDegradeOffDiagonal(), extractBit(), Qureg::numAmpsPerChunk, Qureg::numQubitsRepresented, Qureg::pairStateVec, qreal, and Qureg::stateVec.
Referenced by densmatr_mixDamping().
◆ densmatr_mixDampingLocal()
Definition at line 174 of file QuEST_cpu.c.
References Qureg::chunkId, Qureg::numAmpsPerChunk, Qureg::numQubitsRepresented, qreal, and Qureg::stateVec.
Referenced by densmatr_mixDamping().
◆ densmatr_mixDensityMatrix()
Definition at line 890 of file QuEST_cpu.c.
References Qureg::numAmpsPerChunk, qreal, and Qureg::stateVec.
◆ densmatr_mixDephasing()
Definition at line 79 of file QuEST_cpu.c.
References densmatr_oneQubitDegradeOffDiagonal(), and qreal.
Referenced by densmatr_mixDepolarisingDistributed().
◆ densmatr_mixDepolarisingDistributed()
Definition at line 224 of file QuEST_cpu.c.
References Qureg::chunkId, densmatr_mixDephasing(), extractBit(), Qureg::numAmpsPerChunk, Qureg::numQubitsRepresented, Qureg::pairStateVec, and Qureg::stateVec.
Referenced by densmatr_mixDepolarising().
◆ densmatr_mixDepolarisingLocal()
Definition at line 125 of file QuEST_cpu.c.
References Qureg::chunkId, Qureg::numAmpsPerChunk, Qureg::numQubitsRepresented, qreal, and Qureg::stateVec.
Referenced by densmatr_mixDepolarising().
◆ densmatr_mixTwoQubitDephasing()
Definition at line 84 of file QuEST_cpu.c.
References Qureg::chunkId, Qureg::numAmpsPerChunk, Qureg::numQubitsRepresented, qreal, and Qureg::stateVec.
◆ densmatr_mixTwoQubitDepolarisingDistributed()
void densmatr_mixTwoQubitDepolarisingDistributed | ( | Qureg | qureg, |
int | targetQubit, | ||
int | qubit2, | ||
qreal | delta, | ||
qreal | gamma | ||
) |
Definition at line 541 of file QuEST_cpu.c.
References Qureg::chunkId, extractBit(), Qureg::numAmpsPerChunk, Qureg::numQubitsRepresented, Qureg::pairStateVec, and Qureg::stateVec.
Referenced by densmatr_mixTwoQubitDepolarising().
◆ densmatr_mixTwoQubitDepolarisingLocal()
void densmatr_mixTwoQubitDepolarisingLocal | ( | Qureg | qureg, |
int | qubit1, | ||
int | qubit2, | ||
qreal | delta, | ||
qreal | gamma | ||
) |
Definition at line 387 of file QuEST_cpu.c.
References Qureg::chunkId, Qureg::numAmpsPerChunk, Qureg::numQubitsRepresented, qreal, and Qureg::stateVec.
Referenced by densmatr_mixTwoQubitDepolarising().
◆ densmatr_mixTwoQubitDepolarisingLocalPart1()
void densmatr_mixTwoQubitDepolarisingLocalPart1 | ( | Qureg | qureg, |
int | qubit1, | ||
int | qubit2, | ||
qreal | delta | ||
) |
Definition at line 488 of file QuEST_cpu.c.
References Qureg::chunkId, Qureg::numAmpsPerChunk, Qureg::numQubitsRepresented, qreal, and Qureg::stateVec.
Referenced by densmatr_mixTwoQubitDepolarising().
◆ densmatr_mixTwoQubitDepolarisingQ1LocalQ2DistributedPart3()
void densmatr_mixTwoQubitDepolarisingQ1LocalQ2DistributedPart3 | ( | Qureg | qureg, |
int | targetQubit, | ||
int | qubit2, | ||
qreal | delta, | ||
qreal | gamma | ||
) |
Definition at line 632 of file QuEST_cpu.c.
References Qureg::chunkId, extractBit(), Qureg::numAmpsPerChunk, Qureg::numQubitsRepresented, Qureg::pairStateVec, and Qureg::stateVec.
Referenced by densmatr_mixTwoQubitDepolarising().
◆ densmatr_oneQubitDegradeOffDiagonal()
Definition at line 48 of file QuEST_cpu.c.
References Qureg::chunkId, Qureg::numAmpsPerChunk, Qureg::numQubitsRepresented, and Qureg::stateVec.
Referenced by densmatr_mixDampingDistributed(), and densmatr_mixDephasing().
◆ getBitMaskParity()
int getBitMaskParity | ( | long long int | mask | ) |
◆ normaliseSomeAmps()
Definition at line 744 of file QuEST_cpu.c.
References Qureg::stateVec.
Referenced by alternateNormZeroingSomeAmpBlocks(), and densmatr_collapseToKnownProbOutcome().
◆ qsortComp()
int qsortComp | ( | const void * | a, |
const void * | b | ||
) |
Definition at line 1842 of file QuEST_cpu.c.
Referenced by statevec_multiControlledMultiQubitUnitaryLocal().
◆ statevec_applyDiagonalOp()
void statevec_applyDiagonalOp | ( | Qureg | qureg, |
DiagonalOp | op | ||
) |
Definition at line 3661 of file QuEST_cpu.c.
References DiagonalOp::imag, Qureg::numAmpsPerChunk, qreal, DiagonalOp::real, and Qureg::stateVec.
◆ statevec_calcExpecDiagonalOpLocal()
Complex statevec_calcExpecDiagonalOpLocal | ( | Qureg | qureg, |
DiagonalOp | op | ||
) |
Definition at line 3738 of file QuEST_cpu.c.
References Complex::imag, DiagonalOp::imag, Qureg::numAmpsPerChunk, qreal, Complex::real, DiagonalOp::real, and Qureg::stateVec.
Referenced by statevec_calcExpecDiagonalOp().
◆ statevec_calcInnerProductLocal()
Definition at line 1071 of file QuEST_cpu.c.
References Complex::imag, Qureg::numAmpsPerChunk, qreal, Complex::real, and Qureg::stateVec.
Referenced by statevec_calcInnerProduct().
◆ statevec_cloneQureg()
Definition at line 1506 of file QuEST_cpu.c.
References Qureg::numAmpsPerChunk, qreal, and Qureg::stateVec.
◆ statevec_collapseToKnownProbOutcomeDistributedRenorm()
void statevec_collapseToKnownProbOutcomeDistributedRenorm | ( | Qureg | qureg, |
int | measureQubit, | ||
qreal | totalProbability | ||
) |
Renormalise parts of the state vector where measureQubit=0 or 1, based on the total probability of that qubit being in state 0 or 1.
Measure in Zero performs an irreversible change to the state vector: it updates the vector according to the event that the value 'outcome' has been measured on the qubit indicated by measureQubit (where this label starts from 0, of course). It achieves this by setting all inconsistent amplitudes to 0 and then renormalising based on the total probability of measuring measureQubit=0 if outcome=0 and measureQubit=1 if outcome=1. In the distributed version, one block (with measureQubit=0 in the first half of the block and measureQubit=1 in the second half of the block) is spread over multiple chunks, meaning that each chunks performs only renormalisation or only setting amplitudes to 0. This function handles the renormalisation.
- Parameters
-
[in,out] qureg object representing the set of qubits [in] measureQubit qubit to measure [in] totalProbability probability of qubit measureQubit being zero
Definition at line 3462 of file QuEST_cpu.c.
References Qureg::numAmpsPerChunk, qreal, and Qureg::stateVec.
Referenced by statevec_collapseToKnownProbOutcome().
◆ statevec_collapseToKnownProbOutcomeLocal()
void statevec_collapseToKnownProbOutcomeLocal | ( | Qureg | qureg, |
int | measureQubit, | ||
int | outcome, | ||
qreal | totalProbability | ||
) |
Update the state vector to be consistent with measuring measureQubit=0 if outcome=0 and measureQubit=1 if outcome=1.
Performs an irreversible change to the state vector: it updates the vector according to the event that an outcome have been measured on the qubit indicated by measureQubit (where this label starts from 0, of course). It achieves this by setting all inconsistent amplitudes to 0 and then renormalising based on the total probability of measuring measureQubit=0 or 1 according to the value of outcome. In the local version, one or more blocks (with measureQubit=0 in the first half of the block and measureQubit=1 in the second half of the block) fit entirely into one chunk.
- Parameters
-
[in,out] qureg object representing the set of qubits [in] measureQubit qubit to measure [in] totalProbability probability of qubit measureQubit being either zero or one [in] outcome to measure the probability of and set the state to – either zero or one
Definition at line 3380 of file QuEST_cpu.c.
References Qureg::numAmpsPerChunk, qreal, and Qureg::stateVec.
Referenced by statevec_collapseToKnownProbOutcome().
◆ statevec_collapseToOutcomeDistributedSetZero()
void statevec_collapseToOutcomeDistributedSetZero | ( | Qureg | qureg | ) |
Set all amplitudes in one chunk to 0.
Measure in Zero performs an irreversible change to the state vector: it updates the vector according to the event that a zero have been measured on the qubit indicated by measureQubit (where this label starts from 0, of course). It achieves this by setting all inconsistent amplitudes to 0 and then renormalising based on the total probability of measuring measureQubit=0 or 1. In the distributed version, one block (with measureQubit=0 in the first half of the block and measureQubit=1 in the second half of the block) is spread over multiple chunks, meaning that each chunks performs only renormalisation or only setting amplitudes to 0. This function handles setting amplitudes to 0.
- Parameters
-
[in,out] qureg object representing the set of qubits [in] measureQubit qubit to measure
Definition at line 3501 of file QuEST_cpu.c.
References Qureg::numAmpsPerChunk, qreal, and Qureg::stateVec.
Referenced by statevec_collapseToKnownProbOutcome().
◆ statevec_compactUnitaryDistributed()
void statevec_compactUnitaryDistributed | ( | Qureg | qureg, |
Complex | rot1, | ||
Complex | rot2, | ||
ComplexArray | stateVecUp, | ||
ComplexArray | stateVecLo, | ||
ComplexArray | stateVecOut | ||
) |
Rotate a single qubit in the state vector of probability amplitudes, given two complex numbers alpha and beta, and a subset of the state vector with upper and lower block values stored seperately.
- Parameters
-
[in,out] qureg object representing the set of qubits [in] rot1 rotation angle [in] rot2 rotation angle [in] stateVecUp probability amplitudes in upper half of a block [in] stateVecLo probability amplitudes in lower half of a block [out] stateVecOut array section to update (will correspond to either the lower or upper half of a block)
Definition at line 2001 of file QuEST_cpu.c.
References Complex::imag, Qureg::numAmpsPerChunk, qreal, and Complex::real.
Referenced by statevec_compactUnitary().
◆ statevec_compactUnitaryLocal()
Definition at line 1688 of file QuEST_cpu.c.
References Complex::imag, Qureg::numAmpsPerChunk, qreal, Complex::real, and Qureg::stateVec.
Referenced by statevec_compactUnitary().
◆ statevec_compareStates()
Definition at line 1675 of file QuEST_cpu.c.
References Qureg::numAmpsPerChunk, qreal, and Qureg::stateVec.
◆ statevec_controlledCompactUnitaryDistributed()
void statevec_controlledCompactUnitaryDistributed | ( | Qureg | qureg, |
int | controlQubit, | ||
Complex | rot1, | ||
Complex | rot2, | ||
ComplexArray | stateVecUp, | ||
ComplexArray | stateVecLo, | ||
ComplexArray | stateVecOut | ||
) |
Rotate a single qubit in the state vector of probability amplitudes, given two complex numbers alpha and beta and a subset of the state vector with upper and lower block values stored seperately.
Only perform the rotation where the control qubit is one.
- Parameters
-
[in,out] qureg object representing the set of qubits [in] controlQubit qubit to determine whether or not to perform a rotation [in] rot1 rotation angle [in] rot2 rotation angle [in] stateVecUp probability amplitudes in upper half of a block [in] stateVecLo probability amplitudes in lower half of a block [out] stateVecOut array section to update (will correspond to either the lower or upper half of a block)
Definition at line 2319 of file QuEST_cpu.c.
References Qureg::chunkId, extractBit(), Complex::imag, Qureg::numAmpsPerChunk, qreal, and Complex::real.
Referenced by statevec_controlledCompactUnitary().
◆ statevec_controlledCompactUnitaryLocal()
void statevec_controlledCompactUnitaryLocal | ( | Qureg | qureg, |
int | controlQubit, | ||
int | targetQubit, | ||
Complex | alpha, | ||
Complex | beta | ||
) |
Definition at line 2101 of file QuEST_cpu.c.
References Qureg::chunkId, extractBit(), Complex::imag, Qureg::numAmpsPerChunk, qreal, Complex::real, and Qureg::stateVec.
Referenced by statevec_controlledCompactUnitary().
◆ statevec_controlledNotDistributed()
void statevec_controlledNotDistributed | ( | Qureg | qureg, |
int | controlQubit, | ||
ComplexArray | stateVecIn, | ||
ComplexArray | stateVecOut | ||
) |
Rotate a single qubit by {{0,1},{1,0}.
Operate on a subset of the state vector with upper and lower block values stored seperately. This rotation is just swapping upper and lower values, and stateVecIn must already be the correct section for this chunk. Only perform the rotation for elements where controlQubit is one.
- Parameters
-
[in,out] qureg object representing the set of qubits [in] stateVecIn probability amplitudes in lower or upper half of a block depending on chunkId [out] stateVecOut array section to update (will correspond to either the lower or upper half of a block)
Definition at line 2646 of file QuEST_cpu.c.
References Qureg::chunkId, extractBit(), Qureg::numAmpsPerChunk, and qreal.
Referenced by statevec_controlledNot().
◆ statevec_controlledNotLocal()
void statevec_controlledNotLocal | ( | Qureg | qureg, |
int | controlQubit, | ||
int | targetQubit | ||
) |
Definition at line 2584 of file QuEST_cpu.c.
References Qureg::chunkId, extractBit(), Qureg::numAmpsPerChunk, qreal, and Qureg::stateVec.
Referenced by statevec_controlledNot().
◆ statevec_controlledPauliYDistributed()
void statevec_controlledPauliYDistributed | ( | Qureg | qureg, |
int | controlQubit, | ||
ComplexArray | stateVecIn, | ||
ComplexArray | stateVecOut, | ||
int | conjFac | ||
) |
Definition at line 2830 of file QuEST_cpu.c.
References Qureg::chunkId, extractBit(), Qureg::numAmpsPerChunk, and qreal.
Referenced by statevec_controlledPauliY(), and statevec_controlledPauliYConj().
◆ statevec_controlledPauliYLocal()
void statevec_controlledPauliYLocal | ( | Qureg | qureg, |
int | controlQubit, | ||
int | targetQubit, | ||
int | conjFac | ||
) |
Definition at line 2776 of file QuEST_cpu.c.
References Qureg::chunkId, extractBit(), Qureg::numAmpsPerChunk, qreal, and Qureg::stateVec.
Referenced by statevec_controlledPauliY(), and statevec_controlledPauliYConj().
◆ statevec_controlledPhaseFlip()
void statevec_controlledPhaseFlip | ( | Qureg | qureg, |
int | idQubit1, | ||
int | idQubit2 | ||
) |
Definition at line 3300 of file QuEST_cpu.c.
References Qureg::chunkId, extractBit(), Qureg::numAmpsPerChunk, qreal, and Qureg::stateVec.
◆ statevec_controlledPhaseShift()
Definition at line 3019 of file QuEST_cpu.c.
References Qureg::chunkId, extractBit(), Qureg::numAmpsPerChunk, qreal, and Qureg::stateVec.
◆ statevec_controlledUnitaryDistributed()
void statevec_controlledUnitaryDistributed | ( | Qureg | qureg, |
int | controlQubit, | ||
Complex | rot1, | ||
Complex | rot2, | ||
ComplexArray | stateVecUp, | ||
ComplexArray | stateVecLo, | ||
ComplexArray | stateVecOut | ||
) |
Rotate a single qubit in the state vector of probability amplitudes, given two complex numbers alpha and beta and a subset of the state vector with upper and lower block values stored seperately.
Only perform the rotation where the control qubit is one.
- Parameters
-
[in,out] qureg object representing the set of qubits [in] controlQubit qubit to determine whether or not to perform a rotation [in] rot1 rotation angle [in] rot2 rotation angle [in] stateVecUp probability amplitudes in upper half of a block [in] stateVecLo probability amplitudes in lower half of a block [out] stateVecOut array section to update (will correspond to either the lower or upper half of a block)
Definition at line 2381 of file QuEST_cpu.c.
References Qureg::chunkId, extractBit(), Complex::imag, Qureg::numAmpsPerChunk, qreal, and Complex::real.
Referenced by statevec_controlledUnitary().
◆ statevec_controlledUnitaryLocal()
void statevec_controlledUnitaryLocal | ( | Qureg | qureg, |
int | controlQubit, | ||
int | targetQubit, | ||
ComplexMatrix2 | u | ||
) |
Definition at line 2241 of file QuEST_cpu.c.
References Qureg::chunkId, extractBit(), ComplexMatrix2::imag, Qureg::numAmpsPerChunk, qreal, ComplexMatrix2::real, and Qureg::stateVec.
Referenced by statevec_controlledUnitary().
◆ statevec_createQureg()
Definition at line 1279 of file QuEST_cpu.c.
References Qureg::chunkId, Qureg::isDensityMatrix, Qureg::numAmpsPerChunk, Qureg::numAmpsTotal, Qureg::numChunks, Qureg::numQubitsInStateVec, QuESTEnv::numRanks, Qureg::pairStateVec, QuESTEnv::rank, and Qureg::stateVec.
◆ statevec_destroyQureg()
Definition at line 1317 of file QuEST_cpu.c.
References Qureg::numAmpsPerChunk, Qureg::numAmpsTotal, Qureg::numQubitsInStateVec, QuESTEnv::numRanks, Qureg::pairStateVec, and Qureg::stateVec.
◆ statevec_findProbabilityOfZeroDistributed()
Measure the probability of a specified qubit being in the zero state across all amplitudes held in this chunk.
Size of regions to skip is a multiple of chunkSize. The results are communicated and aggregated by the caller
- Parameters
-
[in] qureg object representing the set of qubits
- Returns
- probability of qubit measureQubit being zero
Definition at line 3262 of file QuEST_cpu.c.
References Qureg::numAmpsPerChunk, qreal, and Qureg::stateVec.
Referenced by statevec_calcProbOfOutcome().
◆ statevec_findProbabilityOfZeroLocal()
Measure the total probability of a specified qubit being in the zero state across all amplitudes in this chunk.
Size of regions to skip is less than the size of one chunk.
- Parameters
-
[in] qureg object representing the set of qubits [in] measureQubit qubit to measure
- Returns
- probability of qubit measureQubit being zero
Definition at line 3206 of file QuEST_cpu.c.
References Qureg::numAmpsPerChunk, qreal, and Qureg::stateVec.
Referenced by statevec_calcProbOfOutcome().
◆ statevec_getEnvironmentString()
Definition at line 1390 of file QuEST_cpu.c.
References Qureg::numQubitsInStateVec, and QuESTEnv::numRanks.
◆ statevec_hadamardDistributed()
void statevec_hadamardDistributed | ( | Qureg | qureg, |
ComplexArray | stateVecUp, | ||
ComplexArray | stateVecLo, | ||
ComplexArray | stateVecOut, | ||
int | updateUpper | ||
) |
Rotate a single qubit by {{1,1},{1,-1}}/sqrt2.
Operate on a subset of the state vector with upper and lower block values stored seperately. This rotation is just swapping upper and lower values, and stateVecIn must already be the correct section for this chunk
- Parameters
-
[in,out] qureg object representing the set of qubits [in] stateVecIn probability amplitudes in lower or upper half of a block depending on chunkId [in] updateUpper flag, 1: updating upper values, 0: updating lower values in block [out] stateVecOut array section to update (will correspond to either the lower or upper half of a block)
Definition at line 2932 of file QuEST_cpu.c.
References Qureg::numAmpsPerChunk, and qreal.
Referenced by statevec_hadamard().
◆ statevec_hadamardLocal()
void statevec_hadamardLocal | ( | Qureg | qureg, |
int | targetQubit | ||
) |
Definition at line 2872 of file QuEST_cpu.c.
References Qureg::numAmpsPerChunk, qreal, and Qureg::stateVec.
Referenced by statevec_hadamard().
◆ statevec_initBlankState()
void statevec_initBlankState | ( | Qureg | qureg | ) |
Definition at line 1398 of file QuEST_cpu.c.
References Qureg::numAmpsPerChunk, qreal, and Qureg::stateVec.
Referenced by statevec_initZeroState().
◆ statevec_initClassicalState()
void statevec_initClassicalState | ( | Qureg | qureg, |
long long int | stateInd | ||
) |
Definition at line 1470 of file QuEST_cpu.c.
References Qureg::chunkId, Qureg::numAmpsPerChunk, qreal, and Qureg::stateVec.
◆ statevec_initDebugState()
void statevec_initDebugState | ( | Qureg | qureg | ) |
Initialise the state vector of probability amplitudes to an (unphysical) state with each component of each probability amplitude a unique floating point value.
For debugging processes
- Parameters
-
[in,out] qureg object representing the set of qubits to be initialised
Definition at line 1591 of file QuEST_cpu.c.
References Qureg::chunkId, Qureg::numAmpsPerChunk, qreal, and Qureg::stateVec.
◆ statevec_initPlusState()
void statevec_initPlusState | ( | Qureg | qureg | ) |
Definition at line 1438 of file QuEST_cpu.c.
References Qureg::numAmpsPerChunk, Qureg::numChunks, qreal, and Qureg::stateVec.
◆ statevec_initStateFromSingleFile()
Definition at line 1625 of file QuEST_cpu.c.
References Qureg::chunkId, Qureg::numAmpsPerChunk, Qureg::numChunks, qreal, Qureg::stateVec, and syncQuESTEnv().
◆ statevec_initStateOfSingleQubit()
void statevec_initStateOfSingleQubit | ( | Qureg * | qureg, |
int | qubitId, | ||
int | outcome | ||
) |
Initialise the state vector of probability amplitudes such that one qubit is set to 'outcome' and all other qubits are in an equal superposition of zero and one.
- Parameters
-
[in,out] qureg object representing the set of qubits to be initialised [in] qubitId id of qubit to set to state 'outcome' [in] value of qubit 'qubitId'
Definition at line 1545 of file QuEST_cpu.c.
References Qureg::chunkId, extractBit(), Qureg::numAmpsPerChunk, Qureg::numChunks, qreal, and Qureg::stateVec.
◆ statevec_initZeroState()
void statevec_initZeroState | ( | Qureg | qureg | ) |
Definition at line 1428 of file QuEST_cpu.c.
References Qureg::chunkId, Qureg::stateVec, and statevec_initBlankState().
◆ statevec_multiControlledMultiQubitUnitaryLocal()
void statevec_multiControlledMultiQubitUnitaryLocal | ( | Qureg | qureg, |
long long int | ctrlMask, | ||
int * | targs, | ||
int | numTargs, | ||
ComplexMatrixN | u | ||
) |
Definition at line 1846 of file QuEST_cpu.c.
References Qureg::chunkId, extractBit(), flipBit(), ComplexMatrixN::imag, insertZeroBit(), Qureg::numAmpsPerChunk, ComplexMatrixN::numQubits, qreal, qsortComp(), ComplexMatrixN::real, and Qureg::stateVec.
Referenced by statevec_multiControlledMultiQubitUnitary().
◆ statevec_multiControlledPhaseFlip()
void statevec_multiControlledPhaseFlip | ( | Qureg | qureg, |
int * | controlQubits, | ||
int | numControlQubits | ||
) |
Definition at line 3331 of file QuEST_cpu.c.
References Qureg::chunkId, getQubitBitMask(), Qureg::numAmpsPerChunk, qreal, and Qureg::stateVec.
◆ statevec_multiControlledPhaseShift()
void statevec_multiControlledPhaseShift | ( | Qureg | qureg, |
int * | controlQubits, | ||
int | numControlQubits, | ||
qreal | angle | ||
) |
Definition at line 3059 of file QuEST_cpu.c.
References Qureg::chunkId, getQubitBitMask(), Qureg::numAmpsPerChunk, qreal, and Qureg::stateVec.
◆ statevec_multiControlledTwoQubitUnitaryLocal()
void statevec_multiControlledTwoQubitUnitaryLocal | ( | Qureg | qureg, |
long long int | ctrlMask, | ||
int | q1, | ||
int | q2, | ||
ComplexMatrix4 | u | ||
) |
Definition at line 1747 of file QuEST_cpu.c.
References Qureg::chunkId, flipBit(), ComplexMatrix4::imag, insertTwoZeroBits(), Qureg::numAmpsPerChunk, qreal, ComplexMatrix4::real, and Qureg::stateVec.
Referenced by statevec_multiControlledTwoQubitUnitary().
◆ statevec_multiControlledUnitaryDistributed()
void statevec_multiControlledUnitaryDistributed | ( | Qureg | qureg, |
int | targetQubit, | ||
long long int | ctrlQubitsMask, | ||
long long int | ctrlFlipMask, | ||
Complex | rot1, | ||
Complex | rot2, | ||
ComplexArray | stateVecUp, | ||
ComplexArray | stateVecLo, | ||
ComplexArray | stateVecOut | ||
) |
Apply a unitary operation to a single qubit in the state vector of probability amplitudes, given a subset of the state vector with upper and lower block values stored seperately.
Only perform the rotation where all the control qubits are 1.
- Parameters
-
[in,out] qureg object representing the set of qubits [in] targetQubit qubit to rotate [in] ctrlQubitsMask a bit mask indicating whether each qubit is a control (1) or not (0) [in] ctrlFlipMask a bit mask indicating whether each qubit (only controls are relevant) should be flipped when checking the control condition [in] rot1 rotation angle [in] rot2 rotation angle [in] stateVecUp probability amplitudes in upper half of a block [in] stateVecLo probability amplitudes in lower half of a block [out] stateVecOut array section to update (will correspond to either the lower or upper half of a block)
Definition at line 2447 of file QuEST_cpu.c.
References Qureg::chunkId, Complex::imag, Qureg::numAmpsPerChunk, qreal, and Complex::real.
Referenced by statevec_multiControlledUnitary().
◆ statevec_multiControlledUnitaryLocal()
void statevec_multiControlledUnitaryLocal | ( | Qureg | qureg, |
int | targetQubit, | ||
long long int | ctrlQubitsMask, | ||
long long int | ctrlFlipMask, | ||
ComplexMatrix2 | u | ||
) |
Definition at line 2173 of file QuEST_cpu.c.
References Qureg::chunkId, ComplexMatrix2::imag, Qureg::numAmpsPerChunk, qreal, ComplexMatrix2::real, and Qureg::stateVec.
Referenced by statevec_multiControlledUnitary().
◆ statevec_multiRotateZ()
Definition at line 3109 of file QuEST_cpu.c.
References Qureg::chunkId, getBitMaskParity(), Qureg::numAmpsPerChunk, qreal, and Qureg::stateVec.
◆ statevec_pauliXDistributed()
void statevec_pauliXDistributed | ( | Qureg | qureg, |
ComplexArray | stateVecIn, | ||
ComplexArray | stateVecOut | ||
) |
Rotate a single qubit by {{0,1},{1,0}.
Operate on a subset of the state vector with upper and lower block values stored seperately. This rotation is just swapping upper and lower values, and stateVecIn must already be the correct section for this chunk
- Remarks
- Qubits are zero-based and the
the first qubit is the rightmost
- Parameters
-
[in,out] qureg object representing the set of qubits [in] stateVecIn probability amplitudes in lower or upper half of a block depending on chunkId [out] stateVecOut array section to update (will correspond to either the lower or upper half of a block)
Definition at line 2556 of file QuEST_cpu.c.
References Qureg::numAmpsPerChunk, and qreal.
Referenced by statevec_pauliX().
◆ statevec_pauliXLocal()
void statevec_pauliXLocal | ( | Qureg | qureg, |
int | targetQubit | ||
) |
Definition at line 2498 of file QuEST_cpu.c.
References Qureg::numAmpsPerChunk, qreal, and Qureg::stateVec.
Referenced by statevec_pauliX().
◆ statevec_pauliYDistributed()
void statevec_pauliYDistributed | ( | Qureg | qureg, |
ComplexArray | stateVecIn, | ||
ComplexArray | stateVecOut, | ||
int | updateUpper, | ||
int | conjFac | ||
) |
Rotate a single qubit by +-{{0,-i},{i,0}.
Operate on a subset of the state vector with upper and lower block values stored seperately. This rotation is just swapping upper and lower values, and stateVecIn must already be the correct section for this chunk
- Remarks
- Qubits are zero-based and the
the first qubit is the rightmost
- Parameters
-
[in,out] qureg object representing the set of qubits [in] stateVecIn probability amplitudes in lower or upper half of a block depending on chunkId [in] updateUpper flag, 1: updating upper values, 0: updating lower values in block [out] stateVecOut array section to update (will correspond to either the lower or upper half of a block)
Definition at line 2739 of file QuEST_cpu.c.
References Qureg::numAmpsPerChunk, and qreal.
Referenced by statevec_pauliY(), and statevec_pauliYConj().
◆ statevec_pauliYLocal()
void statevec_pauliYLocal | ( | Qureg | qureg, |
int | targetQubit, | ||
int | conjFac | ||
) |
Definition at line 2682 of file QuEST_cpu.c.
References Qureg::numAmpsPerChunk, qreal, and Qureg::stateVec.
Referenced by statevec_pauliY(), and statevec_pauliYConj().
◆ statevec_phaseShiftByTerm()
Definition at line 2978 of file QuEST_cpu.c.
References Qureg::chunkId, extractBit(), Complex::imag, Qureg::numAmpsPerChunk, qreal, Complex::real, and Qureg::stateVec.
◆ statevec_reportStateToScreen()
Definition at line 1366 of file QuEST_cpu.c.
References Qureg::chunkId, Qureg::numAmpsPerChunk, Qureg::numChunks, Qureg::numQubitsInStateVec, Qureg::stateVec, and syncQuESTEnv().
◆ statevec_setAmps()
void statevec_setAmps | ( | Qureg | qureg, |
long long int | startInd, | ||
qreal * | reals, | ||
qreal * | imags, | ||
long long int | numAmps | ||
) |
Definition at line 1237 of file QuEST_cpu.c.
References Qureg::chunkId, Qureg::numAmpsPerChunk, qreal, and Qureg::stateVec.
◆ statevec_setWeightedQureg()
void statevec_setWeightedQureg | ( | Complex | fac1, |
Qureg | qureg1, | ||
Complex | fac2, | ||
Qureg | qureg2, | ||
Complex | facOut, | ||
Qureg | out | ||
) |
Definition at line 3619 of file QuEST_cpu.c.
References Complex::imag, Qureg::numAmpsPerChunk, qreal, Complex::real, and Qureg::stateVec.
◆ statevec_swapQubitAmpsDistributed()
void statevec_swapQubitAmpsDistributed | ( | Qureg | qureg, |
int | pairRank, | ||
int | qb1, | ||
int | qb2 | ||
) |
qureg.pairStateVec contains the entire set of amplitudes of the paired node which includes the set of all amplitudes which need to be swapped between |..0..1..> and |..1..0..>
Definition at line 3579 of file QuEST_cpu.c.
References Qureg::chunkId, flipBit(), isOddParity(), Qureg::numAmpsPerChunk, Qureg::pairStateVec, qreal, and Qureg::stateVec.
Referenced by statevec_swapQubitAmps().
◆ statevec_swapQubitAmpsLocal()
void statevec_swapQubitAmpsLocal | ( | Qureg | qureg, |
int | qb1, | ||
int | qb2 | ||
) |
It is ensured that all amplitudes needing to be swapped are on this node.
This means that amplitudes for |a 0..0..> to |a 1..1..> all exist on this node and each node has a different bit-string prefix "a". The prefix 'a' (and ergo, the chunkID) don't enter the calculations for the offset of |a 0..1..> and |a 1..0..> from |a 0..0..> and ergo are not featured below.
Definition at line 3536 of file QuEST_cpu.c.
References flipBit(), insertTwoZeroBits(), Qureg::numAmpsPerChunk, qreal, and Qureg::stateVec.
Referenced by statevec_swapQubitAmps().
◆ statevec_unitaryDistributed()
void statevec_unitaryDistributed | ( | Qureg | qureg, |
Complex | rot1, | ||
Complex | rot2, | ||
ComplexArray | stateVecUp, | ||
ComplexArray | stateVecLo, | ||
ComplexArray | stateVecOut | ||
) |
Apply a unitary operation to a single qubit given a subset of the state vector with upper and lower block values stored seperately.
- Remarks
- Qubits are zero-based and the first qubit is the rightmost
- Parameters
-
[in,out] qureg object representing the set of qubits [in] u unitary matrix to apply [in] stateVecUp probability amplitudes in upper half of a block [in] stateVecLo probability amplitudes in lower half of a block [out] stateVecOut array section to update (will correspond to either the lower or upper half of a block)
Definition at line 2056 of file QuEST_cpu.c.
References Complex::imag, Qureg::numAmpsPerChunk, qreal, and Complex::real.
Referenced by statevec_unitary().
◆ statevec_unitaryLocal()
void statevec_unitaryLocal | ( | Qureg | qureg, |
int | targetQubit, | ||
ComplexMatrix2 | u | ||
) |
Definition at line 1932 of file QuEST_cpu.c.
References ComplexMatrix2::imag, Qureg::numAmpsPerChunk, qreal, ComplexMatrix2::real, and Qureg::stateVec.
Referenced by statevec_unitary().
◆ zeroSomeAmps()
void zeroSomeAmps | ( | Qureg | qureg, |
long long int | startInd, | ||
long long int | numAmps | ||
) |
Definition at line 734 of file QuEST_cpu.c.
References Qureg::stateVec.
Referenced by alternateNormZeroingSomeAmpBlocks(), and densmatr_collapseToKnownProbOutcome().