1c532df63SYohann // Copyright (c) 2017-2018, Lawrence Livermore National Security, LLC. 2c532df63SYohann // Produced at the Lawrence Livermore National Laboratory. LLNL-CODE-734707. 3c532df63SYohann // All Rights reserved. See files LICENSE and NOTICE for details. 4c532df63SYohann // 5c532df63SYohann // This file is part of CEED, a collection of benchmarks, miniapps, software 6c532df63SYohann // libraries and APIs for efficient high-order finite element and spectral 7c532df63SYohann // element discretizations for exascale applications. For more information and 8c532df63SYohann // source code availability see http://github.com/ceed. 9c532df63SYohann // 10c532df63SYohann // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC, 11c532df63SYohann // a collaborative effort of two U.S. Department of Energy organizations (Office 12c532df63SYohann // of Science and the National Nuclear Security Administration) responsible for 13c532df63SYohann // the planning and preparation of a capable exascale ecosystem, including 14c532df63SYohann // software, applications, hardware, advanced system engineering and early 15c532df63SYohann // testbed platforms, in support of the nation's exascale computing imperative. 16c532df63SYohann 17c532df63SYohann #include <ceed-backend.h> 18c532df63SYohann #include <ceed.h> 19c532df63SYohann #include "ceed-cuda-shared.h" 20c532df63SYohann #include "../cuda/ceed-cuda.h" 21c532df63SYohann 22c532df63SYohann //********************* 23c532df63SYohann // shared mem kernels 24c532df63SYohann static const char *kernelsShared = QUOTE( 25c532df63SYohann 26c532df63SYohann inline __device__ void add(CeedScalar *r_V, const CeedScalar *r_U) { 27c532df63SYohann for (int i = 0; i < Q1D; i++) 28c532df63SYohann r_V[i] += r_U[i]; 29c532df63SYohann } 30c532df63SYohann 31c532df63SYohann ////////// 32c532df63SYohann // 1D // 33c532df63SYohann ////////// 34c532df63SYohann 35c532df63SYohann inline __device__ void readDofs1d(const int elem, const int tidx, 36d94769d2SYohann Dudouit const int tidy, const int tidz,const int comp, 37*7f823360Sjeremylt const int nelem, const CeedScalar *d_U, 38*7f823360Sjeremylt CeedScalar *slice) { 39c532df63SYohann for (int i = 0; i < P1D; i++) 40d94769d2SYohann Dudouit slice[i+tidz*Q1D] = d_U[i + comp*P1D + elem*BASIS_NCOMP*P1D]; 41c532df63SYohann for (int i = P1D; i < Q1D; i++) 42d94769d2SYohann Dudouit slice[i+tidz*Q1D] = 0.0; 43c532df63SYohann } 44c532df63SYohann 45c532df63SYohann inline __device__ void writeDofs1d(const int elem, const int tidx, 46c532df63SYohann const int tidy, const int comp, 47288c0443SJeremy L Thompson const int nelem, const CeedScalar &r_V, 48288c0443SJeremy L Thompson CeedScalar *d_V) { 49c532df63SYohann if (tidx<P1D) { 50c532df63SYohann d_V[tidx + comp*P1D + elem*BASIS_NCOMP*P1D] = r_V; 51c532df63SYohann } 52c532df63SYohann } 53c532df63SYohann 54c532df63SYohann inline __device__ void readQuads1d(const int elem, const int tidx, 55d94769d2SYohann Dudouit const int tidy, const int tidz, const int comp, 56288c0443SJeremy L Thompson const int dim, const int nelem, 57288c0443SJeremy L Thompson const CeedScalar *d_U, CeedScalar *slice) { 58c532df63SYohann for (int i = 0; i < Q1D; i++) 594d537eeaSYohann slice[i+tidz*Q1D] = d_U[i + elem*Q1D + comp*Q1D*nelem + 604d537eeaSYohann dim*BASIS_NCOMP*nelem*Q1D]; 61c532df63SYohann } 62c532df63SYohann 63c532df63SYohann inline __device__ void writeQuads1d(const int elem, const int tidx, 64c532df63SYohann const int tidy, const int comp, 65288c0443SJeremy L Thompson const int dim, const int nelem, 66288c0443SJeremy L Thompson const CeedScalar &r_V, CeedScalar *d_V) { 67c532df63SYohann d_V[tidx + elem*Q1D + comp*Q1D*nelem + dim*BASIS_NCOMP*nelem*Q1D] = r_V; 68c532df63SYohann } 69c532df63SYohann 70c532df63SYohann inline __device__ void ContractX1d(CeedScalar *slice, const int tidx, 71d94769d2SYohann Dudouit const int tidy, const int tidz, 72288c0443SJeremy L Thompson const CeedScalar &U, const CeedScalar *B, 73288c0443SJeremy L Thompson CeedScalar &V) { 74c532df63SYohann V = 0.0; 75c532df63SYohann for (int i = 0; i < P1D; ++i) { 76d94769d2SYohann Dudouit V += B[i + tidx*P1D] * slice[i+tidz*Q1D];//contract x direction 77c532df63SYohann } 78c532df63SYohann } 79c532df63SYohann 80c532df63SYohann inline __device__ void ContractTransposeX1d(CeedScalar *slice, const int tidx, 81d94769d2SYohann Dudouit const int tidy, const int tidz, 82c532df63SYohann const CeedScalar &U, const CeedScalar *B, CeedScalar &V) { 83c532df63SYohann V = 0.0; 84c532df63SYohann for (int i = 0; i < Q1D; ++i) { 85d94769d2SYohann Dudouit V += B[tidx + i*P1D] * slice[i+tidz*Q1D];//contract x direction 86c532df63SYohann } 87c532df63SYohann } 88c532df63SYohann 89c532df63SYohann inline __device__ void interp1d(const CeedInt nelem, const int transpose, 90288c0443SJeremy L Thompson const CeedScalar *c_B, 91288c0443SJeremy L Thompson const CeedScalar *__restrict__ d_U, 92c532df63SYohann CeedScalar *__restrict__ d_V, 93c532df63SYohann CeedScalar *slice) { 94c532df63SYohann CeedScalar r_V; 95c532df63SYohann CeedScalar r_t; 96c532df63SYohann 97c532df63SYohann const int tidx = threadIdx.x; 98c532df63SYohann const int tidy = threadIdx.y; 99d94769d2SYohann Dudouit const int tidz = threadIdx.z; 100c532df63SYohann 101c532df63SYohann 102c532df63SYohann for (CeedInt elem = blockIdx.x*blockDim.z + threadIdx.z; elem < nelem; 103c532df63SYohann elem += gridDim.x*blockDim.z) { 104c532df63SYohann for(int comp=0; comp<BASIS_NCOMP; comp++) { 105c532df63SYohann if(!transpose) { 106d94769d2SYohann Dudouit readDofs1d(elem, tidx, tidy, tidz, comp, nelem, d_U, slice); 107d94769d2SYohann Dudouit ContractX1d(slice, tidx, tidy, tidz, r_t, c_B, r_V); 108c532df63SYohann writeQuads1d(elem, tidx, tidy, comp, 0, nelem, r_V, d_V); 109c532df63SYohann } else { 110d94769d2SYohann Dudouit readQuads1d(elem, tidx, tidy, tidz, comp, 0, nelem, d_U, slice); 111d94769d2SYohann Dudouit ContractTransposeX1d(slice, tidx, tidy, tidz, r_t, c_B, r_V); 112c532df63SYohann writeDofs1d(elem, tidx, tidy, comp, nelem, r_V, d_V); 113c532df63SYohann } 114c532df63SYohann } 115c532df63SYohann } 116c532df63SYohann } 117c532df63SYohann 118c532df63SYohann inline __device__ void grad1d(const CeedInt nelem, const int transpose, 119c532df63SYohann const CeedScalar *c_B, const CeedScalar *c_G, 120288c0443SJeremy L Thompson const CeedScalar *__restrict__ d_U, 121288c0443SJeremy L Thompson CeedScalar *__restrict__ d_V, 122c532df63SYohann CeedScalar *slice) { 123c532df63SYohann CeedScalar r_U; 124c532df63SYohann CeedScalar r_V; 125c532df63SYohann 126c532df63SYohann const int tidx = threadIdx.x; 127d94769d2SYohann Dudouit const int tidy = threadIdx.y; 128d94769d2SYohann Dudouit const int tidz = threadIdx.z; 129c532df63SYohann int dim; 130c532df63SYohann 131c532df63SYohann for (CeedInt elem = blockIdx.x*blockDim.z + threadIdx.z; elem < nelem; 132c532df63SYohann elem += gridDim.x*blockDim.z) { 133c532df63SYohann for(int comp=0; comp<BASIS_NCOMP; comp++) { 134c532df63SYohann if(!transpose) { 135d94769d2SYohann Dudouit readDofs1d(elem, tidx, tidy, tidz, comp, nelem, d_U, slice); 136d94769d2SYohann Dudouit ContractX1d(slice, tidx, tidy, tidz, r_U, c_G, r_V); 137c532df63SYohann dim = 0; 138c532df63SYohann writeQuads1d(elem, tidx, tidy, comp, dim, nelem, r_V, d_V); 139c532df63SYohann } else { 140c532df63SYohann dim = 0; 141d94769d2SYohann Dudouit readQuads1d(elem, tidx, tidy, tidz, comp, dim, nelem, d_U, slice); 142d94769d2SYohann Dudouit ContractTransposeX1d(slice, tidx, tidy, tidz, r_U, c_G, r_V); 143c532df63SYohann writeDofs1d(elem, tidx, tidy, comp, nelem, r_V, d_V); 144c532df63SYohann } 145c532df63SYohann } 146c532df63SYohann } 147c532df63SYohann } 148c532df63SYohann ////////// 149c532df63SYohann // 2D // 150c532df63SYohann ////////// 151c532df63SYohann 152c532df63SYohann inline __device__ void readDofs2d(const int elem, const int tidx, 153c532df63SYohann const int tidy, const int comp, 154288c0443SJeremy L Thompson const int nelem, const CeedScalar *d_U, 155288c0443SJeremy L Thompson CeedScalar &U) { 156c532df63SYohann U = (tidx<P1D 157*7f823360Sjeremylt && tidy<P1D) ? d_U[tidx + tidy*P1D + comp*P1D*P1D + 158*7f823360Sjeremylt elem*BASIS_NCOMP*P1D*P1D ] : 159c532df63SYohann 0.0; 160c532df63SYohann } 161c532df63SYohann 162c532df63SYohann inline __device__ void writeDofs2d(const int elem, const int tidx, 163c532df63SYohann const int tidy, const int comp, 164288c0443SJeremy L Thompson const int nelem, const CeedScalar &r_V, 165288c0443SJeremy L Thompson CeedScalar *d_V) { 166c532df63SYohann if (tidx<P1D && tidy<P1D) { 167c532df63SYohann d_V[tidx + tidy*P1D + comp*P1D*P1D + elem*BASIS_NCOMP*P1D*P1D ] = r_V; 168c532df63SYohann } 169c532df63SYohann } 170c532df63SYohann 171c532df63SYohann inline __device__ void readQuads2d(const int elem, const int tidx, 172c532df63SYohann const int tidy, const int comp, 173288c0443SJeremy L Thompson const int dim, const int nelem, 174288c0443SJeremy L Thompson const CeedScalar *d_U, CeedScalar &U ) { 175c532df63SYohann U = d_U[tidx + tidy*Q1D + elem*Q1D*Q1D + comp*Q1D*Q1D*nelem + 176c532df63SYohann dim*BASIS_NCOMP*nelem*Q1D*Q1D]; 177c532df63SYohann } 178c532df63SYohann 179c532df63SYohann inline __device__ void writeQuads2d(const int elem, const int tidx, 180c532df63SYohann const int tidy, const int comp, 181288c0443SJeremy L Thompson const int dim, const int nelem, 182288c0443SJeremy L Thompson const CeedScalar &r_V, CeedScalar *d_V) { 183c532df63SYohann d_V[tidx + tidy*Q1D + elem*Q1D*Q1D + comp*Q1D*Q1D*nelem + 184c532df63SYohann dim*BASIS_NCOMP*nelem*Q1D*Q1D ] = r_V; 185c532df63SYohann } 186c532df63SYohann 187c532df63SYohann inline __device__ void ContractX2d(CeedScalar *slice, const int tidx, 1884247ecf3SYohann Dudouit const int tidy, const int tidz, 189288c0443SJeremy L Thompson const CeedScalar &U, const CeedScalar *B, 190288c0443SJeremy L Thompson CeedScalar &V) { 1914247ecf3SYohann Dudouit slice[tidx+tidy*Q1D+tidz*Q1D*Q1D] = U; 192c532df63SYohann __syncthreads(); 193c532df63SYohann V = 0.0; 194c532df63SYohann for (int i = 0; i < P1D; ++i) { 1954247ecf3SYohann Dudouit V += B[i + tidx*P1D] * slice[i + tidy*Q1D + tidz*Q1D*Q1D];//contract x direction 196c532df63SYohann } 197c532df63SYohann __syncthreads(); 198c532df63SYohann } 199c532df63SYohann 200c532df63SYohann inline __device__ void ContractY2d(CeedScalar *slice, const int tidx, 2014247ecf3SYohann Dudouit const int tidy, const int tidz, 202288c0443SJeremy L Thompson const CeedScalar &U, const CeedScalar *B, 203288c0443SJeremy L Thompson CeedScalar &V) { 2044247ecf3SYohann Dudouit slice[tidx+tidy*Q1D+tidz*Q1D*Q1D] = U; 205c532df63SYohann __syncthreads(); 206c532df63SYohann V = 0.0; 207c532df63SYohann for (int i = 0; i < P1D; ++i) { 2084247ecf3SYohann Dudouit V += B[i + tidy*P1D] * slice[tidx + i*Q1D + tidz*Q1D*Q1D];//contract y direction 209c532df63SYohann } 210c532df63SYohann __syncthreads(); 211c532df63SYohann } 212c532df63SYohann 213c532df63SYohann inline __device__ void ContractTransposeY2d(CeedScalar *slice, const int tidx, 2144247ecf3SYohann Dudouit const int tidy, const int tidz, 215c532df63SYohann const CeedScalar &U, const CeedScalar *B, CeedScalar &V) { 2164247ecf3SYohann Dudouit slice[tidx+tidy*Q1D+tidz*Q1D*Q1D] = U; 217c532df63SYohann __syncthreads(); 218c532df63SYohann V = 0.0; 219c532df63SYohann if (tidy<P1D) { 220c532df63SYohann for (int i = 0; i < Q1D; ++i) { 2214247ecf3SYohann Dudouit V += B[tidy + i*P1D] * slice[tidx + i*Q1D + tidz*Q1D*Q1D];//contract y direction 222c532df63SYohann } 223c532df63SYohann } 224c532df63SYohann __syncthreads(); 225c532df63SYohann } 226c532df63SYohann 227c532df63SYohann inline __device__ void ContractTransposeX2d(CeedScalar *slice, const int tidx, 2284247ecf3SYohann Dudouit const int tidy, const int tidz, 229c532df63SYohann const CeedScalar &U, const CeedScalar *B, CeedScalar &V) { 2304247ecf3SYohann Dudouit slice[tidx+tidy*Q1D+tidz*Q1D*Q1D] = U; 231c532df63SYohann __syncthreads(); 232c532df63SYohann V = 0.0; 233c532df63SYohann if (tidx<P1D) { 234c532df63SYohann for (int i = 0; i < Q1D; ++i) { 2354247ecf3SYohann Dudouit V += B[tidx + i*P1D] * slice[i + tidy*Q1D + tidz*Q1D*Q1D];//contract x direction 236c532df63SYohann } 237c532df63SYohann } 238c532df63SYohann __syncthreads(); 239c532df63SYohann } 240c532df63SYohann 241c532df63SYohann inline __device__ void interp2d(const CeedInt nelem, const int transpose, 242288c0443SJeremy L Thompson const CeedScalar *c_B, 243288c0443SJeremy L Thompson const CeedScalar *__restrict__ d_U, 244c532df63SYohann CeedScalar *__restrict__ d_V, 245c532df63SYohann CeedScalar *slice) { 246c532df63SYohann CeedScalar r_V; 247c532df63SYohann CeedScalar r_t; 248c532df63SYohann 249c532df63SYohann const int tidx = threadIdx.x; 250c532df63SYohann const int tidy = threadIdx.y; 2514247ecf3SYohann Dudouit const int tidz = threadIdx.z; 2524247ecf3SYohann Dudouit const int blockElem = tidz/BASIS_NCOMP; 2534247ecf3SYohann Dudouit const int elemsPerBlock = blockDim.z/BASIS_NCOMP; 2544247ecf3SYohann Dudouit const int comp = tidz%BASIS_NCOMP; 255c532df63SYohann 2564247ecf3SYohann Dudouit for (CeedInt elem = blockIdx.x*elemsPerBlock + blockElem; elem < nelem; 2574247ecf3SYohann Dudouit elem += gridDim.x*elemsPerBlock) { 2584247ecf3SYohann Dudouit const int comp = tidz%BASIS_NCOMP; 259c532df63SYohann r_V = 0.0; 260c532df63SYohann r_t = 0.0; 261c532df63SYohann if(!transpose) { 262c532df63SYohann readDofs2d(elem, tidx, tidy, comp, nelem, d_U, r_V); 2634247ecf3SYohann Dudouit ContractX2d(slice, tidx, tidy, tidz, r_V, c_B, r_t); 2644247ecf3SYohann Dudouit ContractY2d(slice, tidx, tidy, tidz, r_t, c_B, r_V); 265c532df63SYohann writeQuads2d(elem, tidx, tidy, comp, 0, nelem, r_V, d_V); 266c532df63SYohann } else { 267c532df63SYohann readQuads2d(elem, tidx, tidy, comp, 0, nelem, d_U, r_V); 2684247ecf3SYohann Dudouit ContractTransposeY2d(slice, tidx, tidy, tidz, r_V, c_B, r_t); 2694247ecf3SYohann Dudouit ContractTransposeX2d(slice, tidx, tidy, tidz, r_t, c_B, r_V); 270c532df63SYohann writeDofs2d(elem, tidx, tidy, comp, nelem, r_V, d_V); 271c532df63SYohann } 272c532df63SYohann } 273c532df63SYohann } 274c532df63SYohann 275c532df63SYohann inline __device__ void grad2d(const CeedInt nelem, const int transpose, 276c532df63SYohann const CeedScalar *c_B, const CeedScalar *c_G, 277*7f823360Sjeremylt const CeedScalar *__restrict__ d_U, 278*7f823360Sjeremylt CeedScalar *__restrict__ d_V, CeedScalar *slice) { 279c532df63SYohann CeedScalar r_U; 280c532df63SYohann CeedScalar r_V; 281c532df63SYohann CeedScalar r_t; 282c532df63SYohann 283c532df63SYohann const int tidx = threadIdx.x; 284c532df63SYohann const int tidy = threadIdx.y; 2854247ecf3SYohann Dudouit const int tidz = threadIdx.z; 2864247ecf3SYohann Dudouit const int blockElem = tidz/BASIS_NCOMP; 2874247ecf3SYohann Dudouit const int elemsPerBlock = blockDim.z/BASIS_NCOMP; 2884247ecf3SYohann Dudouit const int comp = tidz%BASIS_NCOMP; 289c532df63SYohann int dim; 290c532df63SYohann 2914247ecf3SYohann Dudouit for (CeedInt elem = blockIdx.x*elemsPerBlock + blockElem; elem < nelem; 2924247ecf3SYohann Dudouit elem += gridDim.x*elemsPerBlock) { 293c532df63SYohann if(!transpose) { 294c532df63SYohann readDofs2d(elem, tidx, tidy, comp, nelem, d_U, r_U); 2954247ecf3SYohann Dudouit ContractX2d(slice, tidx, tidy, tidz, r_U, c_G, r_t); 2964247ecf3SYohann Dudouit ContractY2d(slice, tidx, tidy, tidz, r_t, c_B, r_V); 297c532df63SYohann dim = 0; 298c532df63SYohann writeQuads2d(elem, tidx, tidy, comp, dim, nelem, r_V, d_V); 2994247ecf3SYohann Dudouit ContractX2d(slice, tidx, tidy, tidz, r_U, c_B, r_t); 3004247ecf3SYohann Dudouit ContractY2d(slice, tidx, tidy, tidz, r_t, c_G, r_V); 301c532df63SYohann dim = 1; 302c532df63SYohann writeQuads2d(elem, tidx, tidy, comp, dim, nelem, r_V, d_V); 303c532df63SYohann } else { 304c532df63SYohann dim = 0; 305c532df63SYohann readQuads2d(elem, tidx, tidy, comp, dim, nelem, d_U, r_U); 3064247ecf3SYohann Dudouit ContractTransposeY2d(slice, tidx, tidy, tidz, r_U, c_B, r_t); 3074247ecf3SYohann Dudouit ContractTransposeX2d(slice, tidx, tidy, tidz, r_t, c_G, r_V); 308c532df63SYohann dim = 1; 309c532df63SYohann readQuads2d(elem, tidx, tidy, comp, dim, nelem, d_U, r_U); 3104247ecf3SYohann Dudouit ContractTransposeY2d(slice, tidx, tidy, tidz, r_U, c_G, r_t); 3114247ecf3SYohann Dudouit ContractTransposeX2d(slice, tidx, tidy, tidz, r_t, c_B, r_U); 312c532df63SYohann r_V+=r_U; 313c532df63SYohann writeDofs2d(elem, tidx, tidy, comp, nelem, r_V, d_V); 314c532df63SYohann } 315c532df63SYohann } 316c532df63SYohann } 317c532df63SYohann ////////// 318c532df63SYohann // 3D // 319c532df63SYohann ////////// 320c532df63SYohann 321c532df63SYohann inline __device__ void readDofs3d(const int elem, const int tidx, 322c532df63SYohann const int tidy, const int comp, 323*7f823360Sjeremylt const int nelem, const CeedScalar *d_U, 324*7f823360Sjeremylt CeedScalar *r_U) { 325c532df63SYohann for (int i = 0; i < P1D; i++) 326c532df63SYohann r_U[i] = (tidx<P1D 327c532df63SYohann && tidy<P1D) ? d_U[tidx + tidy*P1D + i*P1D*P1D + comp*P1D*P1D*P1D + 328c532df63SYohann elem*BASIS_NCOMP*P1D*P1D*P1D ] : 0.0; 329c532df63SYohann for (int i = P1D; i < Q1D; i++) 330c532df63SYohann r_U[i] = 0.0; 331c532df63SYohann } 332c532df63SYohann 333c532df63SYohann inline __device__ void readQuads3d(const int elem, const int tidx, 334c532df63SYohann const int tidy, const int comp, 335*7f823360Sjeremylt const int dim, const int nelem, 336*7f823360Sjeremylt const CeedScalar *d_U, CeedScalar *r_U) { 337c532df63SYohann for (int i = 0; i < Q1D; i++) 338c532df63SYohann r_U[i] = d_U[tidx + tidy*Q1D + i*Q1D*Q1D + elem*Q1D*Q1D*Q1D + 339c532df63SYohann comp*Q1D*Q1D*Q1D*nelem + dim*BASIS_NCOMP*nelem*Q1D*Q1D*Q1D]; 340c532df63SYohann } 341c532df63SYohann 342c532df63SYohann inline __device__ void writeDofs3d(const int elem, const int tidx, 343c532df63SYohann const int tidy, const int comp, 344*7f823360Sjeremylt const int nelem, const CeedScalar *r_V, 345*7f823360Sjeremylt CeedScalar *d_V) { 346c532df63SYohann if (tidx<P1D && tidy<P1D) { 347c532df63SYohann for (int i = 0; i < P1D; i++) 348c532df63SYohann d_V[tidx + tidy*P1D + i*P1D*P1D + comp*P1D*P1D*P1D + 349c532df63SYohann elem*BASIS_NCOMP*P1D*P1D*P1D ] = r_V[i]; 350c532df63SYohann } 351c532df63SYohann } 352c532df63SYohann 353c532df63SYohann inline __device__ void writeQuads3d(const int elem, const int tidx, 354c532df63SYohann const int tidy, const int comp, 355*7f823360Sjeremylt const int dim, const int nelem, 356*7f823360Sjeremylt const CeedScalar *r_V, CeedScalar *d_V) { 357c532df63SYohann for (int i = 0; i < Q1D; i++) 358c532df63SYohann d_V[tidx + tidy*Q1D + i*Q1D*Q1D + elem*Q1D*Q1D*Q1D + comp*Q1D*Q1D*Q1D*nelem + 359c532df63SYohann dim*BASIS_NCOMP*nelem*Q1D*Q1D*Q1D ] = r_V[i]; 360c532df63SYohann } 361c532df63SYohann 362c532df63SYohann inline __device__ void ContractX3d(CeedScalar *slice, const int tidx, 363698ebc35SYohann Dudouit const int tidy, const int tidz, 364*7f823360Sjeremylt const CeedScalar *U, const CeedScalar *B, 365*7f823360Sjeremylt CeedScalar *V) { 366c532df63SYohann for (int k = 0; k < P1D; ++k) { 367698ebc35SYohann Dudouit slice[tidx+tidy*Q1D+tidz*Q1D*Q1D] = U[k]; 368c532df63SYohann __syncthreads(); 369c532df63SYohann V[k] = 0.0; 370c532df63SYohann for (int i = 0; i < P1D; ++i) { 3714d537eeaSYohann V[k] += B[i + tidx*P1D] * slice[i + tidy*Q1D + 3724d537eeaSYohann tidz*Q1D*Q1D];//contract x direction 373c532df63SYohann } 374c532df63SYohann __syncthreads(); 375c532df63SYohann } 376c532df63SYohann } 377c532df63SYohann 378c532df63SYohann inline __device__ void ContractY3d(CeedScalar *slice, const int tidx, 379698ebc35SYohann Dudouit const int tidy, const int tidz, 380*7f823360Sjeremylt const CeedScalar *U, const CeedScalar *B, 381*7f823360Sjeremylt CeedScalar *V) { 382c532df63SYohann for (int k = 0; k < P1D; ++k) { 383698ebc35SYohann Dudouit slice[tidx+tidy*Q1D+tidz*Q1D*Q1D] = U[k]; 384c532df63SYohann __syncthreads(); 385c532df63SYohann V[k] = 0.0; 386c532df63SYohann for (int i = 0; i < P1D; ++i) { 3874d537eeaSYohann V[k] += B[i + tidy*P1D] * slice[tidx + i*Q1D + 3884d537eeaSYohann tidz*Q1D*Q1D];//contract y direction 389c532df63SYohann } 390c532df63SYohann __syncthreads(); 391c532df63SYohann } 392c532df63SYohann } 393c532df63SYohann 394c532df63SYohann inline __device__ void ContractZ3d(CeedScalar *slice, const int tidx, 395698ebc35SYohann Dudouit const int tidy, const int tidz, 396*7f823360Sjeremylt const CeedScalar *U, const CeedScalar *B, 397*7f823360Sjeremylt CeedScalar *V) { 398c532df63SYohann for (int k = 0; k < Q1D; ++k) { 399c532df63SYohann V[k] = 0.0; 400c532df63SYohann for (int i = 0; i < P1D; ++i) { 401c532df63SYohann V[k] += B[i + k*P1D] * U[i];//contract z direction 402c532df63SYohann } 403c532df63SYohann } 404c532df63SYohann } 405c532df63SYohann 406c532df63SYohann inline __device__ void ContractTransposeZ3d(CeedScalar *slice, const int tidx, 407698ebc35SYohann Dudouit const int tidy, const int tidz, 408c532df63SYohann const CeedScalar *U, const CeedScalar *B, CeedScalar *V) { 409c532df63SYohann for (int k = 0; k < Q1D; ++k) { 410c532df63SYohann V[k] = 0.0; 411c532df63SYohann if (k<P1D) { 412c532df63SYohann for (int i = 0; i < Q1D; ++i) { 413c532df63SYohann V[k] += B[k + i*P1D] * U[i];//contract z direction 414c532df63SYohann } 415c532df63SYohann } 416c532df63SYohann } 417c532df63SYohann } 418c532df63SYohann 419c532df63SYohann inline __device__ void ContractTransposeY3d(CeedScalar *slice, const int tidx, 420698ebc35SYohann Dudouit const int tidy, const int tidz, 421c532df63SYohann const CeedScalar *U, const CeedScalar *B, CeedScalar *V) { 422c532df63SYohann for (int k = 0; k < P1D; ++k) { 423698ebc35SYohann Dudouit slice[tidx+tidy*Q1D+tidz*Q1D*Q1D] = U[k]; 424c532df63SYohann __syncthreads(); 425c532df63SYohann V[k] = 0.0; 426c532df63SYohann if (tidy<P1D) { 427c532df63SYohann for (int i = 0; i < Q1D; ++i) { 4284d537eeaSYohann V[k] += B[tidy + i*P1D] * slice[tidx + i*Q1D + 4294d537eeaSYohann tidz*Q1D*Q1D];//contract y direction 430c532df63SYohann } 431c532df63SYohann } 432c532df63SYohann __syncthreads(); 433c532df63SYohann } 434c532df63SYohann } 435c532df63SYohann 436c532df63SYohann inline __device__ void ContractTransposeX3d(CeedScalar *slice, const int tidx, 437698ebc35SYohann Dudouit const int tidy, const int tidz, 438c532df63SYohann const CeedScalar *U, const CeedScalar *B, CeedScalar *V) { 439c532df63SYohann for (int k = 0; k < P1D; ++k) { 440698ebc35SYohann Dudouit slice[tidx+tidy*Q1D+tidz*Q1D*Q1D] = U[k]; 441c532df63SYohann __syncthreads(); 442c532df63SYohann V[k] = 0.0; 443c532df63SYohann if (tidx<P1D) { 444c532df63SYohann for (int i = 0; i < Q1D; ++i) { 4454d537eeaSYohann V[k] += B[tidx + i*P1D] * slice[i + tidy*Q1D + 4464d537eeaSYohann tidz*Q1D*Q1D];//contract x direction 447c532df63SYohann } 448c532df63SYohann } 449c532df63SYohann __syncthreads(); 450c532df63SYohann } 451c532df63SYohann } 452c532df63SYohann 453c532df63SYohann inline __device__ void interp3d(const CeedInt nelem, const int transpose, 454*7f823360Sjeremylt const CeedScalar *c_B, 455*7f823360Sjeremylt const CeedScalar *__restrict__ d_U, 456c532df63SYohann CeedScalar *__restrict__ d_V, 457c532df63SYohann CeedScalar *slice) { 458c532df63SYohann CeedScalar r_V[Q1D]; 459c532df63SYohann CeedScalar r_t[Q1D]; 460c532df63SYohann 461c532df63SYohann const int tidx = threadIdx.x; 462c532df63SYohann const int tidy = threadIdx.y; 463698ebc35SYohann Dudouit const int tidz = threadIdx.z; 464698ebc35SYohann Dudouit const int blockElem = tidz/BASIS_NCOMP; 465698ebc35SYohann Dudouit const int elemsPerBlock = blockDim.z/BASIS_NCOMP; 466698ebc35SYohann Dudouit const int comp = tidz%BASIS_NCOMP; 467c532df63SYohann 468698ebc35SYohann Dudouit for (CeedInt elem = blockIdx.x*elemsPerBlock + blockElem; elem < nelem; 469698ebc35SYohann Dudouit elem += gridDim.x*elemsPerBlock) { 470c532df63SYohann for (int i = 0; i < Q1D; ++i) { 471c532df63SYohann r_V[i] = 0.0; 472c532df63SYohann r_t[i] = 0.0; 473c532df63SYohann } 474c532df63SYohann if(!transpose) { 475c532df63SYohann readDofs3d(elem, tidx, tidy, comp, nelem, d_U, r_V); 476698ebc35SYohann Dudouit ContractX3d(slice, tidx, tidy, tidz, r_V, c_B, r_t); 477698ebc35SYohann Dudouit ContractY3d(slice, tidx, tidy, tidz, r_t, c_B, r_V); 478698ebc35SYohann Dudouit ContractZ3d(slice, tidx, tidy, tidz, r_V, c_B, r_t); 479c532df63SYohann writeQuads3d(elem, tidx, tidy, comp, 0, nelem, r_t, d_V); 480c532df63SYohann } else { 481c532df63SYohann readQuads3d(elem, tidx, tidy, comp, 0, nelem, d_U, r_V); 482698ebc35SYohann Dudouit ContractTransposeZ3d(slice, tidx, tidy, tidz, r_V, c_B, r_t); 483698ebc35SYohann Dudouit ContractTransposeY3d(slice, tidx, tidy, tidz, r_t, c_B, r_V); 484698ebc35SYohann Dudouit ContractTransposeX3d(slice, tidx, tidy, tidz, r_V, c_B, r_t); 485c532df63SYohann writeDofs3d(elem, tidx, tidy, comp, nelem, r_t, d_V); 486c532df63SYohann } 487c532df63SYohann } 488c532df63SYohann } 489c532df63SYohann 490c532df63SYohann inline __device__ void grad3d(const CeedInt nelem, const int transpose, 491c532df63SYohann const CeedScalar *c_B, const CeedScalar *c_G, 492*7f823360Sjeremylt const CeedScalar *__restrict__ d_U, 493*7f823360Sjeremylt CeedScalar *__restrict__ d_V, 494c532df63SYohann CeedScalar *slice) { 495c532df63SYohann //use P1D for one of these 496c532df63SYohann CeedScalar r_U[Q1D]; 497c532df63SYohann CeedScalar r_V[Q1D]; 498c532df63SYohann CeedScalar r_t[Q1D]; 499c532df63SYohann 500c532df63SYohann const int tidx = threadIdx.x; 501c532df63SYohann const int tidy = threadIdx.y; 502698ebc35SYohann Dudouit const int tidz = threadIdx.z; 503698ebc35SYohann Dudouit const int blockElem = tidz/BASIS_NCOMP; 504698ebc35SYohann Dudouit const int elemsPerBlock = blockDim.z/BASIS_NCOMP; 505698ebc35SYohann Dudouit const int comp = tidz%BASIS_NCOMP; 506c532df63SYohann int dim; 507c532df63SYohann 508698ebc35SYohann Dudouit for (CeedInt elem = blockIdx.x*elemsPerBlock + blockElem; elem < nelem; 509698ebc35SYohann Dudouit elem += gridDim.x*elemsPerBlock) { 510c532df63SYohann if(!transpose) { 511c532df63SYohann readDofs3d(elem, tidx, tidy, comp, nelem, d_U, r_U); 512698ebc35SYohann Dudouit ContractX3d(slice, tidx, tidy, tidz, r_U, c_G, r_V); 513698ebc35SYohann Dudouit ContractY3d(slice, tidx, tidy, tidz, r_V, c_B, r_t); 514698ebc35SYohann Dudouit ContractZ3d(slice, tidx, tidy, tidz, r_t, c_B, r_V); 515c532df63SYohann dim = 0; 516c532df63SYohann writeQuads3d(elem, tidx, tidy, comp, dim, nelem, r_V, d_V); 517698ebc35SYohann Dudouit ContractX3d(slice, tidx, tidy, tidz, r_U, c_B, r_V); 518698ebc35SYohann Dudouit ContractY3d(slice, tidx, tidy, tidz, r_V, c_G, r_t); 519698ebc35SYohann Dudouit ContractZ3d(slice, tidx, tidy, tidz, r_t, c_B, r_V); 520c532df63SYohann dim = 1; 521c532df63SYohann writeQuads3d(elem, tidx, tidy, comp, dim, nelem, r_V, d_V); 522698ebc35SYohann Dudouit ContractX3d(slice, tidx, tidy, tidz, r_U, c_B, r_V); 523698ebc35SYohann Dudouit ContractY3d(slice, tidx, tidy, tidz, r_V, c_B, r_t); 524698ebc35SYohann Dudouit ContractZ3d(slice, tidx, tidy, tidz, r_t, c_G, r_V); 525c532df63SYohann dim = 2; 526c532df63SYohann writeQuads3d(elem, tidx, tidy, comp, dim, nelem, r_V, d_V); 527c532df63SYohann } else { 528c532df63SYohann dim = 0; 529c532df63SYohann readQuads3d(elem, tidx, tidy, comp, dim, nelem, d_U, r_U); 530698ebc35SYohann Dudouit ContractTransposeZ3d(slice, tidx, tidy, tidz, r_U, c_B, r_t); 531698ebc35SYohann Dudouit ContractTransposeY3d(slice, tidx, tidy, tidz, r_t, c_B, r_U); 532698ebc35SYohann Dudouit ContractTransposeX3d(slice, tidx, tidy, tidz, r_U, c_G, r_V); 533c532df63SYohann dim = 1; 534c532df63SYohann readQuads3d(elem, tidx, tidy, comp, dim, nelem, d_U, r_U); 535698ebc35SYohann Dudouit ContractTransposeZ3d(slice, tidx, tidy, tidz, r_U, c_B, r_t); 536698ebc35SYohann Dudouit ContractTransposeY3d(slice, tidx, tidy, tidz, r_t, c_G, r_U); 537698ebc35SYohann Dudouit ContractTransposeX3d(slice, tidx, tidy, tidz, r_U, c_B, r_t); 538c532df63SYohann add(r_V, r_t); 539c532df63SYohann dim = 2; 540c532df63SYohann readQuads3d(elem, tidx, tidy, comp, dim, nelem, d_U, r_U); 541698ebc35SYohann Dudouit ContractTransposeZ3d(slice, tidx, tidy, tidz, r_U, c_G, r_t); 542698ebc35SYohann Dudouit ContractTransposeY3d(slice, tidx, tidy, tidz, r_t, c_B, r_U); 543698ebc35SYohann Dudouit ContractTransposeX3d(slice, tidx, tidy, tidz, r_U, c_B, r_t); 544c532df63SYohann add(r_V, r_t); 545c532df63SYohann writeDofs3d(elem, tidx, tidy, comp, nelem, r_V, d_V); 546c532df63SYohann } 547c532df63SYohann } 548c532df63SYohann } 549c532df63SYohann 550c532df63SYohann ///////////// 551c532df63SYohann // Kernels // 552c532df63SYohann ///////////// 553c532df63SYohann extern "C" __global__ void interp(const CeedInt nelem, const int transpose, 554*7f823360Sjeremylt const CeedScalar *c_B, 555*7f823360Sjeremylt const CeedScalar *__restrict__ d_U, 556c532df63SYohann CeedScalar *__restrict__ d_V) { 557074be161SYohann Dudouit extern __shared__ double slice[]; 558c532df63SYohann if (BASIS_DIM==1) { 559c532df63SYohann interp1d(nelem, transpose, c_B, d_U, d_V, slice); 560c532df63SYohann } else if (BASIS_DIM==2) { 561c532df63SYohann interp2d(nelem, transpose, c_B, d_U, d_V, slice); 562c532df63SYohann } else if (BASIS_DIM==3) { 563c532df63SYohann interp3d(nelem, transpose, c_B, d_U, d_V, slice); 564c532df63SYohann } 565c532df63SYohann } 566c532df63SYohann 567c532df63SYohann extern "C" __global__ void grad(const CeedInt nelem, const int transpose, 568c532df63SYohann const CeedScalar *c_B, const CeedScalar *c_G, 569*7f823360Sjeremylt const CeedScalar *__restrict__ d_U, 570*7f823360Sjeremylt CeedScalar *__restrict__ d_V) { 571074be161SYohann Dudouit extern __shared__ double slice[]; 572c532df63SYohann if (BASIS_DIM==1) { 573c532df63SYohann grad1d(nelem, transpose, c_B, c_G, d_U, d_V, slice); 574c532df63SYohann } else if (BASIS_DIM==2) { 575c532df63SYohann grad2d(nelem, transpose, c_B, c_G, d_U, d_V, slice); 576c532df63SYohann } else if (BASIS_DIM==3) { 577c532df63SYohann grad3d(nelem, transpose, c_B, c_G, d_U, d_V, slice); 578c532df63SYohann } 579c532df63SYohann } 580c532df63SYohann 581c532df63SYohann ///////////// 582c532df63SYohann // Weights // 583c532df63SYohann ///////////// 584c532df63SYohann __device__ void weight1d(const CeedInt nelem, const CeedScalar *qweight1d, 585c532df63SYohann CeedScalar *w) { 586074be161SYohann Dudouit const int tid = threadIdx.x; 587074be161SYohann Dudouit const CeedScalar weight = qweight1d[tid]; 588074be161SYohann Dudouit for (CeedInt elem = blockIdx.x*blockDim.y + threadIdx.y; elem < nelem; 589074be161SYohann Dudouit elem += gridDim.x*blockDim.y) { 590074be161SYohann Dudouit const int ind = elem*Q1D + tid; 591074be161SYohann Dudouit w[ind] = weight; 592c532df63SYohann } 593c532df63SYohann } 594c532df63SYohann 595c532df63SYohann __device__ void weight2d(const CeedInt nelem, const CeedScalar *qweight1d, 596c532df63SYohann CeedScalar *w) { 597074be161SYohann Dudouit const int i = threadIdx.x; 598074be161SYohann Dudouit const int j = threadIdx.y; 599074be161SYohann Dudouit const CeedScalar weight = qweight1d[i]*qweight1d[j]; 600074be161SYohann Dudouit for (CeedInt elem = blockIdx.x*blockDim.z + threadIdx.z; elem < nelem; 601074be161SYohann Dudouit elem += gridDim.x*blockDim.z) { 602074be161SYohann Dudouit const int ind = elem*Q1D*Q1D + i + j*Q1D; 603074be161SYohann Dudouit w[ind] = weight; 604c532df63SYohann } 605c532df63SYohann } 606c532df63SYohann 607c532df63SYohann __device__ void weight3d(const CeedInt nelem, const CeedScalar *qweight1d, 608c532df63SYohann CeedScalar *w) { 609074be161SYohann Dudouit const int i = threadIdx.x; 610074be161SYohann Dudouit const int j = threadIdx.y; 611074be161SYohann Dudouit const int k = threadIdx.z; 612074be161SYohann Dudouit const CeedScalar weight = qweight1d[i]*qweight1d[j]*qweight1d[k]; 613074be161SYohann Dudouit for (int e = blockIdx.x; e < nelem; e += gridDim.x) { 614074be161SYohann Dudouit const int ind = e*Q1D*Q1D*Q1D + i + j*Q1D + k*Q1D*Q1D; 615074be161SYohann Dudouit w[ind] = weight; 616c532df63SYohann } 617c532df63SYohann } 618c532df63SYohann 619c532df63SYohann extern "C" __global__ void weight(const CeedInt nelem, 620*7f823360Sjeremylt const CeedScalar *__restrict__ qweight1d, 621*7f823360Sjeremylt CeedScalar *__restrict__ v) { 622c532df63SYohann if (BASIS_DIM==1) { 623c532df63SYohann weight1d(nelem, qweight1d, v); 624c532df63SYohann } else if (BASIS_DIM==2) { 625c532df63SYohann weight2d(nelem, qweight1d, v); 626c532df63SYohann } else if (BASIS_DIM==3) { 627c532df63SYohann weight3d(nelem, qweight1d, v); 628c532df63SYohann } 629c532df63SYohann } 630c532df63SYohann 631c532df63SYohann ); 632c532df63SYohann 633c532df63SYohann int CeedCudaInitInterp(CeedScalar *d_B, CeedInt P1d, CeedInt Q1d, 634c532df63SYohann CeedScalar **c_B); 635c532df63SYohann int CeedCudaInitInterpGrad(CeedScalar *d_B, CeedScalar *d_G, CeedInt P1d, 636*7f823360Sjeremylt CeedInt Q1d, CeedScalar **c_B_ptr, 637*7f823360Sjeremylt CeedScalar **c_G_ptr); 638c532df63SYohann 639c532df63SYohann int CeedBasisApplyTensor_Cuda_shared(CeedBasis basis, const CeedInt nelem, 640c532df63SYohann CeedTransposeMode tmode, 641*7f823360Sjeremylt CeedEvalMode emode, CeedVector u, 642*7f823360Sjeremylt CeedVector v) { 643c532df63SYohann int ierr; 644c532df63SYohann Ceed ceed; 645c532df63SYohann ierr = CeedBasisGetCeed(basis, &ceed); CeedChk(ierr); 646c532df63SYohann Ceed_Cuda_shared *ceed_Cuda; 647c532df63SYohann CeedGetData(ceed, (void *) &ceed_Cuda); CeedChk(ierr); 648c532df63SYohann CeedBasis_Cuda_shared *data; 649c532df63SYohann CeedBasisGetData(basis, (void *)&data); CeedChk(ierr); 650c532df63SYohann const CeedInt transpose = tmode == CEED_TRANSPOSE; 6514247ecf3SYohann Dudouit CeedInt dim, ncomp; 652074be161SYohann Dudouit ierr = CeedBasisGetDimension(basis, &dim); CeedChk(ierr); 6534247ecf3SYohann Dudouit ierr = CeedBasisGetNumComponents(basis, &ncomp); CeedChk(ierr); 654c532df63SYohann 655c532df63SYohann const CeedScalar *d_u; 656c532df63SYohann CeedScalar *d_v; 657c532df63SYohann if(emode!=CEED_EVAL_WEIGHT) { 658c532df63SYohann ierr = CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u); CeedChk(ierr); 659c532df63SYohann } 660c532df63SYohann ierr = CeedVectorGetArray(v, CEED_MEM_DEVICE, &d_v); CeedChk(ierr); 661c532df63SYohann 662c532df63SYohann if (tmode == CEED_TRANSPOSE) { 663c532df63SYohann CeedInt length; 664c532df63SYohann ierr = CeedVectorGetLength(v, &length); CeedChk(ierr); 665c532df63SYohann ierr = cudaMemset(d_v, 0, length * sizeof(CeedScalar)); CeedChk(ierr); 666c532df63SYohann } 667c532df63SYohann if (emode == CEED_EVAL_INTERP) { 668c532df63SYohann CeedInt P1d, Q1d; 669c532df63SYohann ierr = CeedBasisGetNumNodes1D(basis, &P1d); CeedChk(ierr); 670c532df63SYohann ierr = CeedBasisGetNumQuadraturePoints1D(basis, &Q1d); CeedChk(ierr); 671c532df63SYohann ierr = CeedCudaInitInterp(data->d_interp1d, P1d, Q1d, &data->c_B); 672c532df63SYohann CeedChk(ierr); 673c532df63SYohann void *interpargs[] = {(void *) &nelem, (void *) &transpose, &data->c_B, &d_u, &d_v}; 6744d537eeaSYohann if (dim==1) { 675d94769d2SYohann Dudouit CeedInt elemsPerBlock = 32; 6764d537eeaSYohann CeedInt grid = nelem/elemsPerBlock + ( (nelem/elemsPerBlock*elemsPerBlock<nelem) 6774d537eeaSYohann ? 1 : 0 ); 678d94769d2SYohann Dudouit CeedInt sharedMem = elemsPerBlock*Q1d*sizeof(CeedScalar); 6794d537eeaSYohann ierr = CeedRunKernelDimSharedCuda(ceed, data->interp, grid, Q1d, 1, 6804d537eeaSYohann elemsPerBlock, sharedMem, 681c532df63SYohann interpargs); 682c532df63SYohann CeedChk(ierr); 683074be161SYohann Dudouit } else if (dim==2) { 6844247ecf3SYohann Dudouit const CeedInt optElems[7] = {0,32,8,6,4,2,8}; 6854247ecf3SYohann Dudouit CeedInt elemsPerBlock = Q1d < 7 ? optElems[Q1d]/ncomp : 1; 6864d537eeaSYohann CeedInt grid = nelem/elemsPerBlock + ( (nelem/elemsPerBlock*elemsPerBlock<nelem) 6874d537eeaSYohann ? 1 : 0 ); 6884247ecf3SYohann Dudouit CeedInt sharedMem = ncomp*elemsPerBlock*Q1d*Q1d*sizeof(CeedScalar); 6894d537eeaSYohann ierr = CeedRunKernelDimSharedCuda(ceed, data->interp, grid, Q1d, Q1d, 6904d537eeaSYohann ncomp*elemsPerBlock, sharedMem, 691074be161SYohann Dudouit interpargs); 692074be161SYohann Dudouit CeedChk(ierr); 693074be161SYohann Dudouit } else if (dim==3) { 6943f63d318SYohann Dudouit CeedInt elemsPerBlock = 1; 6954d537eeaSYohann CeedInt grid = nelem/elemsPerBlock + ( (nelem/elemsPerBlock*elemsPerBlock<nelem) 6964d537eeaSYohann ? 1 : 0 ); 697698ebc35SYohann Dudouit CeedInt sharedMem = ncomp*elemsPerBlock*Q1d*Q1d*sizeof(CeedScalar); 6984d537eeaSYohann ierr = CeedRunKernelDimSharedCuda(ceed, data->interp, grid, Q1d, Q1d, 6994d537eeaSYohann ncomp*elemsPerBlock, sharedMem, 700074be161SYohann Dudouit interpargs); 701074be161SYohann Dudouit CeedChk(ierr); 702074be161SYohann Dudouit } 703c532df63SYohann } else if (emode == CEED_EVAL_GRAD) { 704c532df63SYohann CeedInt P1d, Q1d; 705c532df63SYohann ierr = CeedBasisGetNumNodes1D(basis, &P1d); CeedChk(ierr); 706c532df63SYohann ierr = CeedBasisGetNumQuadraturePoints1D(basis, &Q1d); CeedChk(ierr); 707c532df63SYohann ierr = CeedCudaInitInterpGrad(data->d_interp1d, data->d_grad1d, P1d, 708c532df63SYohann Q1d, &data->c_B, &data->c_G); 709c532df63SYohann CeedChk(ierr); 710c532df63SYohann void *gradargs[] = {(void *) &nelem, (void *) &transpose, &data->c_B, &data->c_G, &d_u, &d_v}; 7114d537eeaSYohann if (dim==1) { 712d94769d2SYohann Dudouit CeedInt elemsPerBlock = 32; 7134d537eeaSYohann CeedInt grid = nelem/elemsPerBlock + ( (nelem/elemsPerBlock*elemsPerBlock<nelem) 7144d537eeaSYohann ? 1 : 0 ); 715d94769d2SYohann Dudouit CeedInt sharedMem = elemsPerBlock*Q1d*sizeof(CeedScalar); 7164d537eeaSYohann ierr = CeedRunKernelDimSharedCuda(ceed, data->grad, grid, Q1d, 1, elemsPerBlock, 7174d537eeaSYohann sharedMem, 718c532df63SYohann gradargs); 719c532df63SYohann CeedChk(ierr); 720074be161SYohann Dudouit } else if (dim==2) { 7214247ecf3SYohann Dudouit const CeedInt optElems[7] = {0,32,8,6,4,2,8}; 7224247ecf3SYohann Dudouit CeedInt elemsPerBlock = Q1d < 7 ? optElems[Q1d]/ncomp : 1; 7234d537eeaSYohann CeedInt grid = nelem/elemsPerBlock + ( (nelem/elemsPerBlock*elemsPerBlock<nelem) 7244d537eeaSYohann ? 1 : 0 ); 7254247ecf3SYohann Dudouit CeedInt sharedMem = ncomp*elemsPerBlock*Q1d*Q1d*sizeof(CeedScalar); 7264d537eeaSYohann ierr = CeedRunKernelDimSharedCuda(ceed, data->grad, grid, Q1d, Q1d, 7274d537eeaSYohann ncomp*elemsPerBlock, sharedMem, 728074be161SYohann Dudouit gradargs); 729074be161SYohann Dudouit CeedChk(ierr); 730074be161SYohann Dudouit } else if (dim==3) { 7313f63d318SYohann Dudouit CeedInt elemsPerBlock = 1; 7324d537eeaSYohann CeedInt grid = nelem/elemsPerBlock + ( (nelem/elemsPerBlock*elemsPerBlock<nelem) 7334d537eeaSYohann ? 1 : 0 ); 734698ebc35SYohann Dudouit CeedInt sharedMem = ncomp*elemsPerBlock*Q1d*Q1d*sizeof(CeedScalar); 7354d537eeaSYohann ierr = CeedRunKernelDimSharedCuda(ceed, data->grad, grid, Q1d, Q1d, 7364d537eeaSYohann ncomp*elemsPerBlock, sharedMem, 737074be161SYohann Dudouit gradargs); 738074be161SYohann Dudouit CeedChk(ierr); 739074be161SYohann Dudouit } 740c532df63SYohann } else if (emode == CEED_EVAL_WEIGHT) { 741074be161SYohann Dudouit CeedInt Q1d; 742074be161SYohann Dudouit ierr = CeedBasisGetNumQuadraturePoints1D(basis, &Q1d); CeedChk(ierr); 743c532df63SYohann void *weightargs[] = {(void *) &nelem, (void *) &data->d_qweight1d, &d_v}; 744074be161SYohann Dudouit if(dim==1) { 745074be161SYohann Dudouit const CeedInt elemsPerBlock = 32/Q1d; 7464d537eeaSYohann const CeedInt gridsize = nelem/elemsPerBlock + ( ( 7474d537eeaSYohann nelem/elemsPerBlock*elemsPerBlock<nelem)? 1 : 0 ); 748*7f823360Sjeremylt ierr = CeedRunKernelDimCuda(ceed, data->weight, gridsize, Q1d, 749*7f823360Sjeremylt elemsPerBlock, 1, weightargs); 7501226057fSYohann Dudouit CeedChk(ierr); 751074be161SYohann Dudouit } else if(dim==2) { 752717ff8a3SYohann Dudouit const CeedInt optElems = 32/(Q1d*Q1d); 753717ff8a3SYohann Dudouit const CeedInt elemsPerBlock = optElems>0?optElems:1; 7544d537eeaSYohann const CeedInt gridsize = nelem/elemsPerBlock + ( ( 7554d537eeaSYohann nelem/elemsPerBlock*elemsPerBlock<nelem)? 1 : 0 ); 7564d537eeaSYohann ierr = CeedRunKernelDimCuda(ceed, data->weight, gridsize, Q1d, Q1d, 7574d537eeaSYohann elemsPerBlock, weightargs); 7581226057fSYohann Dudouit CeedChk(ierr); 759074be161SYohann Dudouit } else if(dim==3) { 760074be161SYohann Dudouit const CeedInt gridsize = nelem; 7614d537eeaSYohann ierr = CeedRunKernelDimCuda(ceed, data->weight, gridsize, Q1d, Q1d, Q1d, 7624d537eeaSYohann weightargs); 7631226057fSYohann Dudouit CeedChk(ierr); 764074be161SYohann Dudouit } 765c532df63SYohann } 766c532df63SYohann 767c532df63SYohann if(emode!=CEED_EVAL_WEIGHT) { 768c532df63SYohann ierr = CeedVectorRestoreArrayRead(u, &d_u); CeedChk(ierr); 769c532df63SYohann } 770c532df63SYohann ierr = CeedVectorRestoreArray(v, &d_v); CeedChk(ierr); 771c532df63SYohann 772c532df63SYohann return 0; 773c532df63SYohann } 774c532df63SYohann 775c532df63SYohann static int CeedBasisDestroy_Cuda_shared(CeedBasis basis) { 776c532df63SYohann int ierr; 777c532df63SYohann Ceed ceed; 778c532df63SYohann ierr = CeedBasisGetCeed(basis, &ceed); CeedChk(ierr); 779c532df63SYohann 780c532df63SYohann CeedBasis_Cuda_shared *data; 781c532df63SYohann ierr = CeedBasisGetData(basis, (void *) &data); CeedChk(ierr); 782c532df63SYohann 783c532df63SYohann CeedChk_Cu(ceed, cuModuleUnload(data->module)); 784c532df63SYohann 785c532df63SYohann ierr = cudaFree(data->d_qweight1d); CeedChk_Cu(ceed, ierr); 786c532df63SYohann ierr = cudaFree(data->d_interp1d); CeedChk_Cu(ceed, ierr); 787c532df63SYohann ierr = cudaFree(data->d_grad1d); CeedChk_Cu(ceed, ierr); 788c532df63SYohann 789c532df63SYohann ierr = CeedFree(&data); CeedChk(ierr); 790c532df63SYohann 791c532df63SYohann return 0; 792c532df63SYohann } 793c532df63SYohann 794c532df63SYohann int CeedBasisCreateTensorH1_Cuda_shared(CeedInt dim, CeedInt P1d, CeedInt Q1d, 795c532df63SYohann const CeedScalar *interp1d, 796c532df63SYohann const CeedScalar *grad1d, 797c532df63SYohann const CeedScalar *qref1d, 798c532df63SYohann const CeedScalar *qweight1d, 799c532df63SYohann CeedBasis basis) { 800c532df63SYohann int ierr; 801c532df63SYohann Ceed ceed; 802c532df63SYohann ierr = CeedBasisGetCeed(basis, &ceed); CeedChk(ierr); 8034d537eeaSYohann if (Q1d<P1d) { 8041226057fSYohann Dudouit return CeedError(ceed, 1, "Backend does not implement underintegrated basis."); 8051226057fSYohann Dudouit } 806c532df63SYohann CeedBasis_Cuda_shared *data; 807c532df63SYohann ierr = CeedCalloc(1, &data); CeedChk(ierr); 808c532df63SYohann 809c532df63SYohann const CeedInt qBytes = Q1d * sizeof(CeedScalar); 810c532df63SYohann ierr = cudaMalloc((void **)&data->d_qweight1d, qBytes); CeedChk_Cu(ceed, ierr); 811c532df63SYohann ierr = cudaMemcpy(data->d_qweight1d, qweight1d, qBytes, 812c532df63SYohann cudaMemcpyHostToDevice); CeedChk_Cu(ceed, ierr); 813c532df63SYohann 814c532df63SYohann const CeedInt iBytes = qBytes * P1d; 815c532df63SYohann ierr = cudaMalloc((void **)&data->d_interp1d, iBytes); CeedChk_Cu(ceed, ierr); 816c532df63SYohann ierr = cudaMemcpy(data->d_interp1d, interp1d, iBytes, 817c532df63SYohann cudaMemcpyHostToDevice); CeedChk_Cu(ceed, ierr); 818c532df63SYohann 819c532df63SYohann ierr = cudaMalloc((void **)&data->d_grad1d, iBytes); CeedChk_Cu(ceed, ierr); 820c532df63SYohann ierr = cudaMemcpy(data->d_grad1d, grad1d, iBytes, 821c532df63SYohann cudaMemcpyHostToDevice); CeedChk_Cu(ceed, ierr); 822c532df63SYohann 823ac421f39SYohann data->d_collograd1d = NULL; 824ac421f39SYohann if (dim==3 && Q1d >= P1d) { 825ac421f39SYohann CeedScalar *collograd1d; 826ac421f39SYohann ierr = CeedMalloc(Q1d*Q1d, &collograd1d); CeedChk(ierr); 827ac421f39SYohann ierr = CeedBasisGetCollocatedGrad(basis, collograd1d); CeedChk(ierr); 828ac421f39SYohann ierr = cudaMalloc((void **)&data->d_collograd1d, qBytes * Q1d); 829ac421f39SYohann CeedChk_Cu(ceed, ierr); 830ac421f39SYohann ierr = cudaMemcpy(data->d_collograd1d, collograd1d, qBytes * Q1d, 831ac421f39SYohann cudaMemcpyHostToDevice); CeedChk_Cu(ceed, ierr); 832ac421f39SYohann } 833ac421f39SYohann 834c532df63SYohann CeedInt ncomp; 835c532df63SYohann ierr = CeedBasisGetNumComponents(basis, &ncomp); CeedChk(ierr); 8364a6d4bbdSYohann Dudouit ierr = CeedCompileCuda(ceed, kernelsShared, &data->module, 7, 837c532df63SYohann "Q1D", Q1d, 838c532df63SYohann "P1D", P1d, 839c532df63SYohann "BASIS_BUF_LEN", ncomp * CeedIntPow(Q1d > P1d ? 840c532df63SYohann Q1d : P1d, dim), 841c532df63SYohann "BASIS_DIM", dim, 842c532df63SYohann "BASIS_NCOMP", ncomp, 843c532df63SYohann "BASIS_ELEMSIZE", CeedIntPow(P1d, dim), 844c532df63SYohann "BASIS_NQPT", CeedIntPow(Q1d, dim) 845c532df63SYohann ); CeedChk(ierr); 8464a6d4bbdSYohann Dudouit ierr = CeedGetKernelCuda(ceed, data->module, "interp", &data->interp); 847c532df63SYohann CeedChk(ierr); 8484a6d4bbdSYohann Dudouit ierr = CeedGetKernelCuda(ceed, data->module, "grad", &data->grad); 849c532df63SYohann CeedChk(ierr); 8504a6d4bbdSYohann Dudouit ierr = CeedGetKernelCuda(ceed, data->module, "weight", &data->weight); 851c532df63SYohann CeedChk(ierr); 852c532df63SYohann 853c532df63SYohann ierr = CeedBasisSetData(basis, (void *)&data); 854c532df63SYohann CeedChk(ierr); 855c532df63SYohann ierr = CeedSetBackendFunction(ceed, "Basis", basis, "Apply", 856c532df63SYohann CeedBasisApplyTensor_Cuda_shared); 857c532df63SYohann CeedChk(ierr); 858c532df63SYohann ierr = CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", 859c532df63SYohann CeedBasisDestroy_Cuda_shared); 860c532df63SYohann CeedChk(ierr); 861c532df63SYohann return 0; 862c532df63SYohann } 863