.. role:: bash(code) :language: bash *************** Getting Started *************** In this section, we show how to implement quantum computing simulation using cuStateVec. Firstly, We describe how to install the library and how to compile it. Then, we present an example code to perform common steps in cuStateVec. ============================ Installation and Compilation ============================ Download the cuQuantum package (which cuStateVec is part of) from https://developer.nvidia.com/cuQuantum-downloads. ----- Linux ----- Assuming cuQuantum has been extracted in ``CUQUANTUM_ROOT``, we update the library path accordingly: .. code-block:: bash export LD_LIBRARY_PATH=${CUQUANTUM_ROOT}/lib:${LD_LIBRARY_PATH} We can compile the sample code we will discuss below (``statevec_example.cu``) via the following command: .. code-block:: bash nvcc statevec_example.cu -I${CUQUANTUM_ROOT}/include -L${CUQUANTUM_ROOT}/lib -lcustatevec -o statevec_example .. note:: Depending on the source of the cuQuantum package, you may need to replace ``lib`` above by ``lib64``. .. _cuStateVec C example: ============ Code Example ============ The following code example shows the common steps to use cuStateVec. Here we apply a Toffoli gate, which inverts the third bit when the first two bits are both 1. .. NOTE: the C example below is translated to Python in sphinx/python/overview.rst, so make sure both codes are in sync! .. figure:: ./figures/toffoli.png :width: 200px :align: center .. code-block:: cpp #include // cudaMalloc, cudaMemcpy, etc. #include // cuDoubleComplex #include // custatevecApplyMatrix #include // printf #include // EXIT_FAILURE int main(void) { const int nIndexBits = 3; const int nSvSize = (1 << nIndexBits); const int nTargets = 1; const int nControls = 2; const int adjoint = 0; int targets[] = {2}; int controls[] = {0, 1}; cuDoubleComplex h_sv[] = {{ 0.0, 0.0}, { 0.0, 0.1}, { 0.1, 0.1}, { 0.1, 0.2}, { 0.2, 0.2}, { 0.3, 0.3}, { 0.3, 0.4}, { 0.4, 0.5}}; cuDoubleComplex h_sv_result[] = {{ 0.0, 0.0}, { 0.0, 0.1}, { 0.1, 0.1}, { 0.4, 0.5}, { 0.2, 0.2}, { 0.3, 0.3}, { 0.3, 0.4}, { 0.1, 0.2}}; cuDoubleComplex matrix[] = {{0.0, 0.0}, {1.0, 0.0}, {1.0, 0.0}, {0.0, 0.0}}; cuDoubleComplex *d_sv; cudaMalloc((void**)&d_sv, nSvSize * sizeof(cuDoubleComplex)); cudaMemcpy(d_sv, h_sv, nSvSize * sizeof(cuDoubleComplex), cudaMemcpyHostToDevice); //-------------------------------------------------------------------------- // custatevec handle initialization custatevecHandle_t handle; custatevecCreate(&handle); void* extraWorkspace = nullptr; size_t extraWorkspaceSizeInBytes = 0; // check the size of external workspace custatevecApplyMatrixGetWorkspaceSize( handle, CUDA_C_64F, nIndexBits, matrix, CUDA_C_64F, CUSTATEVEC_MATRIX_LAYOUT_ROW, adjoint, nTargets, nControls, CUSTATEVEC_COMPUTE_64F, &extraWorkspaceSizeInBytes); // allocate external workspace if necessary if (extraWorkspaceSizeInBytes > 0) cudaMalloc(&extraWorkspace, extraWorkspaceSizeInBytes); // apply gate custatevecApplyMatrix( handle, d_sv, CUDA_C_64F, nIndexBits, matrix, CUDA_C_64F, CUSTATEVEC_MATRIX_LAYOUT_ROW, adjoint, targets, nTargets, controls, nullptr, nControls, CUSTATEVEC_COMPUTE_64F, extraWorkspace, extraWorkspaceSizeInBytes); // destroy handle custatevecDestroy(handle); //-------------------------------------------------------------------------- cudaMemcpy(h_sv, d_sv, nSvSize * sizeof(cuDoubleComplex), cudaMemcpyDeviceToHost); bool correct = true; for (int i = 0; i < nSvSize; i++) { if ((h_sv[i].x != h_sv_result[i].x) || (h_sv[i].y != h_sv_result[i].y)) { correct = false; break; } } if (correct) printf("example PASSED\n"); else printf("example FAILED: wrong result\n"); cudaFree(d_sv); if (extraWorkspaceSizeInBytes) cudaFree(extraWorkspace); return EXIT_SUCCESS; } More samples can be found in the `NVIDIA/cuQuantum `_ repository. =========== Useful tips =========== * For debugging, the environment variable :bash:`CUSTATEVEC_LOG_LEVEL=n` can be set. The level `n` = 0, 1, ..., 5 corresponds to the logger level as described and used in `custatevecLoggerSetLevel`. The environment variable :bash:`CUSTATEVEC_LOG_FILE=` can be used to direct the log output to a custom file at ```` instead of stdout.