// Copyright (c) 2017-2018, Lawrence Livermore National Security, LLC. // Produced at the Lawrence Livermore National Laboratory. LLNL-CODE-734707. // All Rights reserved. See files LICENSE and NOTICE for details. // // This file is part of CEED, a collection of benchmarks, miniapps, software // libraries and APIs for efficient high-order finite element and spectral // element discretizations for exascale applications. For more information and // source code availability see http://github.com/ceed. // // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC, // a collaborative effort of two U.S. Department of Energy organizations (Office // of Science and the National Nuclear Security Administration) responsible for // the planning and preparation of a capable exascale ecosystem, including // software, applications, hardware, advanced system engineering and early // testbed platforms, in support of the nation's exascale computing imperative. #include #include #include "ceed-cuda-shared.h" #include "../cuda/ceed-cuda.h" //********************* // shared mem kernels static const char *kernelsShared = QUOTE( inline __device__ void add(CeedScalar *r_V, const CeedScalar *r_U) { for (int i = 0; i < Q1D; i++) r_V[i] += r_U[i]; } ////////// // 1D // ////////// inline __device__ void readDofs1d(const int elem, const int tidx, const int tidy, const int tidz,const int comp, const int nelem, const CeedScalar *d_U, CeedScalar *slice) { for (int i = 0; i < P1D; i++) slice[i+tidz*Q1D] = d_U[i + comp*P1D + elem*BASIS_NCOMP*P1D]; for (int i = P1D; i < Q1D; i++) slice[i+tidz*Q1D] = 0.0; } inline __device__ void writeDofs1d(const int elem, const int tidx, const int tidy, const int comp, const int nelem, const CeedScalar &r_V, CeedScalar *d_V) { if (tidxd_interp1d, P1d, Q1d, &data->c_B); CeedChk(ierr); void *interpargs[] = {(void *) &nelem, (void *) &transpose, &data->c_B, &d_u, &d_v}; if (dim==1) { CeedInt elemsPerBlock = 32; CeedInt grid = nelem/elemsPerBlock + ( (nelem/elemsPerBlock*elemsPerBlockinterp, grid, Q1d, 1, elemsPerBlock, sharedMem, interpargs); CeedChk(ierr); } else if (dim==2) { const CeedInt optElems[7] = {0,32,8,6,4,2,8}; CeedInt elemsPerBlock = Q1d < 7 ? optElems[Q1d]/ncomp : 1; CeedInt grid = nelem/elemsPerBlock + ( (nelem/elemsPerBlock*elemsPerBlockinterp, grid, Q1d, Q1d, ncomp*elemsPerBlock, sharedMem, interpargs); CeedChk(ierr); } else if (dim==3) { CeedInt elemsPerBlock = 1; CeedInt grid = nelem/elemsPerBlock + ( (nelem/elemsPerBlock*elemsPerBlockinterp, grid, Q1d, Q1d, ncomp*elemsPerBlock, sharedMem, interpargs); CeedChk(ierr); } } else if (emode == CEED_EVAL_GRAD) { CeedInt P1d, Q1d; ierr = CeedBasisGetNumNodes1D(basis, &P1d); CeedChk(ierr); ierr = CeedBasisGetNumQuadraturePoints1D(basis, &Q1d); CeedChk(ierr); ierr = CeedCudaInitInterpGrad(data->d_interp1d, data->d_grad1d, P1d, Q1d, &data->c_B, &data->c_G); CeedChk(ierr); void *gradargs[] = {(void *) &nelem, (void *) &transpose, &data->c_B, &data->c_G, &d_u, &d_v}; if (dim==1) { CeedInt elemsPerBlock = 32; CeedInt grid = nelem/elemsPerBlock + ( (nelem/elemsPerBlock*elemsPerBlockgrad, grid, Q1d, 1, elemsPerBlock, sharedMem, gradargs); CeedChk(ierr); } else if (dim==2) { const CeedInt optElems[7] = {0,32,8,6,4,2,8}; CeedInt elemsPerBlock = Q1d < 7 ? optElems[Q1d]/ncomp : 1; CeedInt grid = nelem/elemsPerBlock + ( (nelem/elemsPerBlock*elemsPerBlockgrad, grid, Q1d, Q1d, ncomp*elemsPerBlock, sharedMem, gradargs); CeedChk(ierr); } else if (dim==3) { CeedInt elemsPerBlock = 1; CeedInt grid = nelem/elemsPerBlock + ( (nelem/elemsPerBlock*elemsPerBlockgrad, grid, Q1d, Q1d, ncomp*elemsPerBlock, sharedMem, gradargs); CeedChk(ierr); } } else if (emode == CEED_EVAL_WEIGHT) { CeedInt Q1d; ierr = CeedBasisGetNumQuadraturePoints1D(basis, &Q1d); CeedChk(ierr); void *weightargs[] = {(void *) &nelem, (void *) &data->d_qweight1d, &d_v}; if(dim==1){ const CeedInt elemsPerBlock = 32/Q1d; const CeedInt gridsize = nelem/elemsPerBlock + ( (nelem/elemsPerBlock*elemsPerBlockweight, gridsize, Q1d, elemsPerBlock, 1, weightargs); } else if(dim==2) { const CeedInt optElems = 32/(Q1d*Q1d); const CeedInt elemsPerBlock = optElems>0?optElems:1; const CeedInt gridsize = nelem/elemsPerBlock + ( (nelem/elemsPerBlock*elemsPerBlockweight, gridsize, Q1d, Q1d, elemsPerBlock, weightargs); } else if(dim==3) { const CeedInt gridsize = nelem; ierr = run_kernel_dim(ceed, data->weight, gridsize, Q1d, Q1d, Q1d, weightargs); } } if(emode!=CEED_EVAL_WEIGHT) { ierr = CeedVectorRestoreArrayRead(u, &d_u); CeedChk(ierr); } ierr = CeedVectorRestoreArray(v, &d_v); CeedChk(ierr); return 0; } static int CeedBasisDestroy_Cuda_shared(CeedBasis basis) { int ierr; Ceed ceed; ierr = CeedBasisGetCeed(basis, &ceed); CeedChk(ierr); CeedBasis_Cuda_shared *data; ierr = CeedBasisGetData(basis, (void *) &data); CeedChk(ierr); CeedChk_Cu(ceed, cuModuleUnload(data->module)); ierr = cudaFree(data->d_qweight1d); CeedChk_Cu(ceed, ierr); ierr = cudaFree(data->d_interp1d); CeedChk_Cu(ceed, ierr); ierr = cudaFree(data->d_grad1d); CeedChk_Cu(ceed, ierr); ierr = CeedFree(&data); CeedChk(ierr); return 0; } int CeedBasisCreateTensorH1_Cuda_shared(CeedInt dim, CeedInt P1d, CeedInt Q1d, const CeedScalar *interp1d, const CeedScalar *grad1d, const CeedScalar *qref1d, const CeedScalar *qweight1d, CeedBasis basis) { int ierr; Ceed ceed; ierr = CeedBasisGetCeed(basis, &ceed); CeedChk(ierr); CeedBasis_Cuda_shared *data; ierr = CeedCalloc(1, &data); CeedChk(ierr); const CeedInt qBytes = Q1d * sizeof(CeedScalar); ierr = cudaMalloc((void **)&data->d_qweight1d, qBytes); CeedChk_Cu(ceed, ierr); ierr = cudaMemcpy(data->d_qweight1d, qweight1d, qBytes, cudaMemcpyHostToDevice); CeedChk_Cu(ceed, ierr); const CeedInt iBytes = qBytes * P1d; ierr = cudaMalloc((void **)&data->d_interp1d, iBytes); CeedChk_Cu(ceed, ierr); ierr = cudaMemcpy(data->d_interp1d, interp1d, iBytes, cudaMemcpyHostToDevice); CeedChk_Cu(ceed, ierr); ierr = cudaMalloc((void **)&data->d_grad1d, iBytes); CeedChk_Cu(ceed, ierr); ierr = cudaMemcpy(data->d_grad1d, grad1d, iBytes, cudaMemcpyHostToDevice); CeedChk_Cu(ceed, ierr); CeedInt ncomp; ierr = CeedBasisGetNumComponents(basis, &ncomp); CeedChk(ierr); ierr = compile(ceed, kernelsShared, &data->module, 7, "Q1D", Q1d, "P1D", P1d, "BASIS_BUF_LEN", ncomp * CeedIntPow(Q1d > P1d ? Q1d : P1d, dim), "BASIS_DIM", dim, "BASIS_NCOMP", ncomp, "BASIS_ELEMSIZE", CeedIntPow(P1d, dim), "BASIS_NQPT", CeedIntPow(Q1d, dim) ); CeedChk(ierr); ierr = get_kernel(ceed, data->module, "interp", &data->interp); CeedChk(ierr); ierr = get_kernel(ceed, data->module, "grad", &data->grad); CeedChk(ierr); ierr = get_kernel(ceed, data->module, "weight", &data->weight); CeedChk(ierr); ierr = CeedBasisSetData(basis, (void *)&data); CeedChk(ierr); ierr = CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyTensor_Cuda_shared); CeedChk(ierr); ierr = CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroy_Cuda_shared); CeedChk(ierr); return 0; } int CeedBasisCreateH1_Cuda_shared(CeedElemTopology topo, CeedInt dim, CeedInt ndof, CeedInt nqpts, const CeedScalar *interp, const CeedScalar *grad, const CeedScalar *qref, const CeedScalar *qweight, CeedBasis basis) { int ierr; Ceed ceed; ierr = CeedBasisGetCeed(basis, &ceed); CeedChk(ierr); return CeedError(ceed, 1, "Backend does not implement generic H1 basis"); }