1 // Copyright (c) 2017-2018, Lawrence Livermore National Security, LLC. 2 // Produced at the Lawrence Livermore National Laboratory. LLNL-CODE-734707. 3 // All Rights reserved. See files LICENSE and NOTICE for details. 4 // 5 // This file is part of CEED, a collection of benchmarks, miniapps, software 6 // libraries and APIs for efficient high-order finite element and spectral 7 // element discretizations for exascale applications. For more information and 8 // source code availability see http://github.com/ceed. 9 // 10 // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC, 11 // a collaborative effort of two U.S. Department of Energy organizations (Office 12 // of Science and the National Nuclear Security Administration) responsible for 13 // the planning and preparation of a capable exascale ecosystem, including 14 // software, applications, hardware, advanced system engineering and early 15 // testbed platforms, in support of the nation's exascale computing imperative. 16 17 #define CEED_DEBUG_COLOR 12 18 19 #include <ceed.h> 20 #include <ceed-backend.h> 21 #include <iostream> 22 #include <string> 23 #include <sstream> 24 #include "ceed-hip-gen.h" 25 #include "../hip/ceed-hip.h" 26 #include "../hip-shared/ceed-hip-shared.h" 27 #include "../hip/ceed-hip-compile.h" 28 29 static const char *deviceFunctions = QUOTE( 30 31 //------------------------------------------------------------------------------ 32 // Typedefs 33 //------------------------------------------------------------------------------ 34 typedef struct { const CeedScalar* in[16]; CeedScalar* out[16]; } HipFields; 35 typedef struct { CeedInt* in[16]; CeedInt* out[16]; } HipFieldsInt; 36 37 typedef struct { 38 CeedInt tidx; 39 CeedInt tidy; 40 CeedInt tidz; 41 CeedInt tid; 42 CeedScalar* slice; 43 } BackendData; 44 45 //------------------------------------------------------------------------------ 46 // Load matrices for basis actions 47 //------------------------------------------------------------------------------ 48 template <int P, int Q> 49 inline __device__ void loadMatrix(BackendData& data, const CeedScalar* d_B, CeedScalar* B) { 50 for (CeedInt i = data.tid; i < P*Q; i += blockDim.x*blockDim.y*blockDim.z) 51 B[i] = d_B[i]; 52 } 53 54 //------------------------------------------------------------------------------ 55 // 1D 56 //------------------------------------------------------------------------------ 57 58 //------------------------------------------------------------------------------ 59 // L-vector -> E-vector, offsets provided 60 //------------------------------------------------------------------------------ 61 template <int NCOMP, int COMPSTRIDE, int P1d> 62 inline __device__ void readDofsOffset1d(BackendData& data, const CeedInt nnodes, const CeedInt elem, const CeedInt* indices, const CeedScalar* d_u, CeedScalar* r_u) { 63 if (data.tidx < P1d) { 64 const CeedInt node = data.tidx; 65 const CeedInt ind = indices[node + elem * P1d]; 66 for (CeedInt comp = 0; comp < NCOMP; ++comp) 67 r_u[comp] = d_u[ind + COMPSTRIDE * comp]; 68 } 69 } 70 71 //------------------------------------------------------------------------------ 72 // L-vector -> E-vector, strided 73 //------------------------------------------------------------------------------ 74 template <int NCOMP, int P1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM> 75 inline __device__ void readDofsStrided1d(BackendData& data, const CeedInt elem, const CeedScalar* d_u, CeedScalar* r_u) { 76 if (data.tidx < P1d) { 77 const CeedInt node = data.tidx; 78 const CeedInt ind = node * STRIDES_NODE + elem * STRIDES_ELEM; 79 for (CeedInt comp = 0; comp < NCOMP; ++comp) 80 r_u[comp] = d_u[ind + comp * STRIDES_COMP]; 81 } 82 } 83 84 //------------------------------------------------------------------------------ 85 // E-vector -> L-vector, offsets provided 86 //------------------------------------------------------------------------------ 87 template <int NCOMP, int COMPSTRIDE, int P1d> 88 inline __device__ void writeDofsOffset1d(BackendData& data, const CeedInt nnodes, const CeedInt elem, const CeedInt* indices, const CeedScalar* r_v, CeedScalar* d_v) { 89 if (data.tidx < P1d) { 90 const CeedInt node = data.tidx; 91 const CeedInt ind = indices[node + elem * P1d]; 92 for (CeedInt comp = 0; comp < NCOMP; ++comp) 93 atomicAdd(&d_v[ind + COMPSTRIDE * comp], r_v[comp]); 94 } 95 } 96 97 //------------------------------------------------------------------------------ 98 // E-vector -> L-vector, strided 99 //------------------------------------------------------------------------------ 100 template <int NCOMP, int P1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM> 101 inline __device__ void writeDofsStrided1d(BackendData& data, const CeedInt elem, const CeedScalar* r_v, CeedScalar* d_v) { 102 if (data.tidx < P1d) { 103 const CeedInt node = data.tidx; 104 const CeedInt ind = node * STRIDES_NODE + elem * STRIDES_ELEM; 105 for (CeedInt comp = 0; comp < NCOMP; ++comp) 106 d_v[ind + comp * STRIDES_COMP] += r_v[comp]; 107 } 108 } 109 110 //------------------------------------------------------------------------------ 111 // 1D tensor contraction x 112 //------------------------------------------------------------------------------ 113 template <int NCOMP, int P1d, int Q1d> 114 inline __device__ void ContractX1d(BackendData& data, const CeedScalar *U, const CeedScalar *B, CeedScalar *V) { 115 data.slice[data.tidx] = *U; 116 __syncthreads(); 117 *V = 0.0; 118 if (data.tidx < Q1d) 119 for (CeedInt i = 0; i < P1d; ++i) 120 *V += B[i + data.tidx*P1d] * data.slice[i]; // Contract x direction 121 __syncthreads(); 122 } 123 124 //------------------------------------------------------------------------------ 125 // 1D transpose tensor contraction x 126 //------------------------------------------------------------------------------ 127 template <int NCOMP, int P1d, int Q1d> 128 inline __device__ void ContractTransposeX1d(BackendData& data, const CeedScalar *U, const CeedScalar *B, CeedScalar *V) { 129 data.slice[data.tidx] = *U; 130 __syncthreads(); 131 *V = 0.0; 132 if (data.tidx < P1d) 133 for (CeedInt i = 0; i < Q1d; ++i) 134 *V += B[data.tidx + i*P1d] * data.slice[i]; // Contract x direction 135 __syncthreads(); 136 } 137 138 //------------------------------------------------------------------------------ 139 // 1D interpolate to quadrature points 140 //------------------------------------------------------------------------------ 141 template <int NCOMP, int P1d, int Q1d> 142 inline __device__ void interp1d(BackendData& data, const CeedScalar *__restrict__ r_U, const CeedScalar *c_B, CeedScalar *__restrict__ r_V) { 143 for (CeedInt comp = 0; comp < NCOMP; comp++) 144 ContractX1d<NCOMP, P1d, Q1d>(data, r_U + comp, c_B, r_V + comp); 145 } 146 147 //------------------------------------------------------------------------------ 148 // 1D interpolate transpose 149 //------------------------------------------------------------------------------ 150 template <int NCOMP, int P1d, int Q1d> 151 inline __device__ void interpTranspose1d(BackendData& data, const CeedScalar *__restrict__ r_U, const CeedScalar *c_B, CeedScalar *__restrict__ r_V) { 152 for (CeedInt comp=0; comp<NCOMP; comp++) 153 ContractTransposeX1d<NCOMP, P1d, Q1d>(data, r_U + comp, c_B, r_V + comp); 154 } 155 156 //------------------------------------------------------------------------------ 157 // 1D derivatives at quadrature points 158 //------------------------------------------------------------------------------ 159 template <int NCOMP, int P1d, int Q1d> 160 inline __device__ void grad1d(BackendData& data, const CeedScalar *__restrict__ r_U, const CeedScalar *c_B, const CeedScalar *c_G, CeedScalar *__restrict__ r_V) { 161 for (CeedInt comp = 0; comp < NCOMP; comp++) 162 ContractX1d<NCOMP, P1d, Q1d>(data, r_U + comp, c_G, r_V + comp); 163 } 164 165 //------------------------------------------------------------------------------ 166 // 1D derivatives transpose 167 //------------------------------------------------------------------------------ 168 template <int NCOMP, int P1d, int Q1d> 169 inline __device__ void gradTranspose1d(BackendData& data, const CeedScalar *__restrict__ r_U, const CeedScalar *c_B, const CeedScalar *c_G, CeedScalar *__restrict__ r_V) { 170 for (CeedInt comp = 0; comp < NCOMP; comp++) 171 ContractTransposeX1d<NCOMP, P1d, Q1d>(data, r_U + comp, c_G, r_V + comp); 172 } 173 174 //------------------------------------------------------------------------------ 175 // 2D 176 //------------------------------------------------------------------------------ 177 178 //------------------------------------------------------------------------------ 179 // L-vector -> E-vector, offsets provided 180 //------------------------------------------------------------------------------ 181 template <int NCOMP, int COMPSTRIDE, int P1d> 182 inline __device__ void readDofsOffset2d(BackendData& data, const CeedInt nnodes, const CeedInt elem, const CeedInt* indices, const CeedScalar* d_u, CeedScalar* r_u) { 183 if (data.tidx < P1d && data.tidy < P1d) { 184 const CeedInt node = data.tidx + data.tidy*P1d; 185 const CeedInt ind = indices[node + elem * P1d*P1d]; 186 for (CeedInt comp = 0; comp < NCOMP; ++comp) 187 r_u[comp] = d_u[ind + COMPSTRIDE * comp]; 188 } 189 } 190 191 //------------------------------------------------------------------------------ 192 // L-vector -> E-vector, strided 193 //------------------------------------------------------------------------------ 194 template <int NCOMP, int P1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM> 195 inline __device__ void readDofsStrided2d(BackendData& data, const CeedInt elem, const CeedScalar* d_u, CeedScalar* r_u) { 196 if (data.tidx < P1d && data.tidy < P1d) { 197 const CeedInt node = data.tidx + data.tidy*P1d; 198 const CeedInt ind = node * STRIDES_NODE + elem * STRIDES_ELEM; 199 for (CeedInt comp = 0; comp < NCOMP; ++comp) 200 r_u[comp] = d_u[ind + comp * STRIDES_COMP]; 201 } 202 } 203 204 //------------------------------------------------------------------------------ 205 // E-vector -> L-vector, offsets provided 206 //------------------------------------------------------------------------------ 207 template <int NCOMP, int COMPSTRIDE, int P1d> 208 inline __device__ void writeDofsOffset2d(BackendData& data, const CeedInt nnodes, const CeedInt elem, const CeedInt* indices, const CeedScalar* r_v, CeedScalar* d_v) { 209 if (data.tidx < P1d && data.tidy < P1d) { 210 const CeedInt node = data.tidx + data.tidy*P1d; 211 const CeedInt ind = indices[node + elem * P1d*P1d]; 212 for (CeedInt comp = 0; comp < NCOMP; ++comp) 213 atomicAdd(&d_v[ind + COMPSTRIDE * comp], r_v[comp]); 214 } 215 } 216 217 //------------------------------------------------------------------------------ 218 // E-vector -> L-vector, strided 219 //------------------------------------------------------------------------------ 220 template <int NCOMP, int P1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM> 221 inline __device__ void writeDofsStrided2d(BackendData& data, const CeedInt elem, const CeedScalar* r_v, CeedScalar* d_v) { 222 if (data.tidx < P1d && data.tidy < P1d) { 223 const CeedInt node = data.tidx + data.tidy*P1d; 224 const CeedInt ind = node * STRIDES_NODE + elem * STRIDES_ELEM; 225 for (CeedInt comp = 0; comp < NCOMP; ++comp) 226 d_v[ind + comp * STRIDES_COMP] += r_v[comp]; 227 } 228 } 229 230 //------------------------------------------------------------------------------ 231 // 2D tensor contraction x 232 //------------------------------------------------------------------------------ 233 template <int NCOMP, int P1d, int Q1d> 234 inline __device__ void ContractX2d(BackendData& data, const CeedScalar *U, const CeedScalar *B, CeedScalar *V) { 235 data.slice[data.tidx+data.tidy*T1d] = *U; 236 __syncthreads(); 237 *V = 0.0; 238 if (data.tidx < Q1d && data.tidy < P1d) 239 for (CeedInt i = 0; i < P1d; ++i) 240 *V += B[i + data.tidx*P1d] * data.slice[i + data.tidy*T1d]; // Contract x direction 241 __syncthreads(); 242 } 243 244 //------------------------------------------------------------------------------ 245 // 2D tensor contract y 246 //------------------------------------------------------------------------------ 247 template <int NCOMP, int P1d, int Q1d> 248 inline __device__ void ContractY2d(BackendData& data, const CeedScalar *U, const CeedScalar *B, CeedScalar *V) { 249 data.slice[data.tidx+data.tidy*T1d] = *U; 250 __syncthreads(); 251 *V = 0.0; 252 if (data.tidx < Q1d && data.tidy < Q1d) 253 for (CeedInt i = 0; i < P1d; ++i) 254 *V += B[i + data.tidy*P1d] * data.slice[data.tidx + i*T1d]; // Contract y direction 255 __syncthreads(); 256 } 257 258 //------------------------------------------------------------------------------ 259 // 2D transpose tensor contract y 260 //------------------------------------------------------------------------------ 261 template <int NCOMP, int P1d, int Q1d> 262 inline __device__ void ContractYTranspose2d(BackendData& data, const CeedScalar *U, const CeedScalar *B, CeedScalar *V) { 263 data.slice[data.tidx+data.tidy*T1d] = *U; 264 __syncthreads(); 265 *V = 0.0; 266 if (data.tidx < Q1d && data.tidy < P1d) 267 for (CeedInt i = 0; i < Q1d; ++i) 268 *V += B[data.tidy + i*P1d] * data.slice[data.tidx + i*T1d]; // Contract y direction 269 __syncthreads(); 270 } 271 272 //------------------------------------------------------------------------------ 273 // 2D transpose tensor contract x 274 //------------------------------------------------------------------------------ 275 template <int NCOMP, int P1d, int Q1d> 276 inline __device__ void ContractXTranspose2d(BackendData& data, const CeedScalar *U, const CeedScalar *B, CeedScalar *V) { 277 data.slice[data.tidx+data.tidy*T1d] = *U; 278 __syncthreads(); 279 *V = 0.0; 280 if (data.tidx < P1d && data.tidy < P1d) 281 for (CeedInt i = 0; i < Q1d; ++i) 282 *V += B[data.tidx + i*P1d] * data.slice[i + data.tidy*T1d]; // Contract x direction 283 __syncthreads(); 284 } 285 286 //------------------------------------------------------------------------------ 287 // 2D transpose tensor contract and add x 288 //------------------------------------------------------------------------------ 289 template <int NCOMP, int P1d, int Q1d> 290 inline __device__ void ContractXTransposeAdd2d(BackendData& data, const CeedScalar *U, const CeedScalar *B, CeedScalar *V) { 291 data.slice[data.tidx+data.tidy*T1d] = *U; 292 __syncthreads(); 293 if (data.tidx < P1d && data.tidy < P1d) 294 for (CeedInt i = 0; i < Q1d; ++i) 295 *V += B[data.tidx + i*P1d] * data.slice[i + data.tidy*T1d]; // Contract x direction 296 __syncthreads(); 297 } 298 299 //------------------------------------------------------------------------------ 300 // 2D interpolate to quadrature points 301 //------------------------------------------------------------------------------ 302 template <int NCOMP, int P1d, int Q1d> 303 inline __device__ void interp2d(BackendData& data, const CeedScalar *__restrict__ r_U, const CeedScalar *c_B, CeedScalar *__restrict__ r_V) { 304 CeedScalar r_t[1]; 305 for (CeedInt comp = 0; comp < NCOMP; comp++) { 306 ContractX2d<NCOMP, P1d, Q1d>(data, r_U + comp, c_B, r_t); 307 ContractY2d<NCOMP, P1d, Q1d>(data, r_t, c_B, r_V + comp); 308 } 309 } 310 311 //------------------------------------------------------------------------------ 312 // 2D interpolate transpose 313 //------------------------------------------------------------------------------ 314 template <int NCOMP, int P1d, int Q1d> 315 inline __device__ void interpTranspose2d(BackendData& data, const CeedScalar *__restrict__ r_U, const CeedScalar *c_B, CeedScalar *__restrict__ r_V) { 316 CeedScalar r_t[1]; 317 for (CeedInt comp = 0; comp < NCOMP; comp++) { 318 ContractYTranspose2d<NCOMP, P1d, Q1d>(data, r_U + comp, c_B, r_t); 319 ContractXTranspose2d<NCOMP, P1d, Q1d>(data, r_t, c_B, r_V + comp); 320 } 321 } 322 323 //------------------------------------------------------------------------------ 324 // 2D derivatives at quadrature points 325 //------------------------------------------------------------------------------ 326 template <int NCOMP, int P1d, int Q1d> 327 inline __device__ void grad2d(BackendData& data, const CeedScalar *__restrict__ r_U, const CeedScalar *c_B, const CeedScalar *c_G, CeedScalar *__restrict__ r_V) { 328 CeedScalar r_t[1]; 329 for (CeedInt comp = 0; comp < NCOMP; comp++) { 330 ContractX2d<NCOMP, P1d, Q1d>(data, r_U + comp, c_G, r_t); 331 ContractY2d<NCOMP, P1d, Q1d>(data, r_t, c_B, r_V + comp + 0*NCOMP); 332 ContractX2d<NCOMP, P1d, Q1d>(data, r_U + comp, c_B, r_t); 333 ContractY2d<NCOMP, P1d, Q1d>(data, r_t, c_G, r_V + comp + 1*NCOMP); 334 } 335 } 336 337 //------------------------------------------------------------------------------ 338 // 2D derivatives transpose 339 //------------------------------------------------------------------------------ 340 template <int NCOMP, int P1d, int Q1d> 341 inline __device__ void gradTranspose2d(BackendData& data, const CeedScalar *__restrict__ r_U, const CeedScalar *c_B, const CeedScalar *c_G, CeedScalar *__restrict__ r_V) { 342 CeedScalar r_t[1]; 343 for (CeedInt comp = 0; comp < NCOMP; comp++) { 344 ContractYTranspose2d<NCOMP, P1d, Q1d>(data, r_U + comp + 0*NCOMP, c_B, r_t); 345 ContractXTranspose2d<NCOMP, P1d, Q1d>(data, r_t, c_G, r_V + comp); 346 ContractYTranspose2d<NCOMP, P1d, Q1d>(data, r_U + comp + 1*NCOMP, c_G, r_t); 347 ContractXTransposeAdd2d<NCOMP, P1d, Q1d>(data, r_t, c_B, r_V + comp); 348 } 349 } 350 351 //------------------------------------------------------------------------------ 352 // 3D 353 //------------------------------------------------------------------------------ 354 355 //------------------------------------------------------------------------------ 356 // L-vector -> E-vector, offsets provided 357 //------------------------------------------------------------------------------ 358 template <int NCOMP, int COMPSTRIDE, int P1d> 359 inline __device__ void readDofsOffset3d(BackendData& data, const CeedInt nnodes, const CeedInt elem, const CeedInt* indices, const CeedScalar* d_u, CeedScalar* r_u) { 360 if (data.tidx < P1d && data.tidy < P1d) 361 for (CeedInt z = 0; z < P1d; ++z) { 362 const CeedInt node = data.tidx + data.tidy*P1d + z*P1d*P1d; 363 const CeedInt ind = indices[node + elem * P1d*P1d*P1d]; 364 for (CeedInt comp = 0; comp < NCOMP; ++comp) 365 r_u[z+comp*P1d] = d_u[ind + COMPSTRIDE * comp]; 366 } 367 } 368 369 //------------------------------------------------------------------------------ 370 // L-vector -> E-vector, strided 371 //------------------------------------------------------------------------------ 372 template <int NCOMP, int P1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM> 373 inline __device__ void readDofsStrided3d(BackendData& data, const CeedInt elem, const CeedScalar* d_u, CeedScalar* r_u) { 374 if (data.tidx < P1d && data.tidy < P1d) 375 for (CeedInt z = 0; z < P1d; ++z) { 376 const CeedInt node = data.tidx + data.tidy*P1d + z*P1d*P1d; 377 const CeedInt ind = node * STRIDES_NODE + elem * STRIDES_ELEM; 378 for (CeedInt comp = 0; comp < NCOMP; ++comp) 379 r_u[z+comp*P1d] = d_u[ind + comp * STRIDES_COMP]; 380 } 381 } 382 383 //------------------------------------------------------------------------------ 384 // E-vector -> Q-vector, offests provided 385 //------------------------------------------------------------------------------ 386 template <int NCOMP, int COMPSTRIDE, int Q1d> 387 inline __device__ void readSliceQuadsOffset3d(BackendData& data, const CeedInt nquads, const CeedInt elem, const CeedInt q, const CeedInt* indices, const CeedScalar* d_u, CeedScalar* r_u) { 388 if (data.tidx < Q1d && data.tidy < Q1d) { 389 const CeedInt node = data.tidx + data.tidy*Q1d + q*Q1d*Q1d; 390 const CeedInt ind = indices[node + elem * Q1d*Q1d*Q1d];; 391 for (CeedInt comp = 0; comp < NCOMP; ++comp) 392 r_u[comp] = d_u[ind + COMPSTRIDE * comp]; 393 } 394 } 395 396 //------------------------------------------------------------------------------ 397 // E-vector -> Q-vector, strided 398 //------------------------------------------------------------------------------ 399 template <int NCOMP, int Q1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM> 400 inline __device__ void readSliceQuadsStrided3d(BackendData& data, const CeedInt elem, const CeedInt q, const CeedScalar* d_u, CeedScalar* r_u) { 401 if (data.tidx < Q1d && data.tidy < Q1d) { 402 const CeedInt node = data.tidx + data.tidy*Q1d + q*Q1d*Q1d; 403 const CeedInt ind = node * STRIDES_NODE + elem * STRIDES_ELEM; 404 for (CeedInt comp = 0; comp < NCOMP; ++comp) 405 r_u[comp] = d_u[ind + comp * STRIDES_COMP]; 406 } 407 } 408 409 //------------------------------------------------------------------------------ 410 // E-vector -> L-vector, offsets provided 411 //------------------------------------------------------------------------------ 412 template <int NCOMP, int COMPSTRIDE, int P1d> 413 inline __device__ void writeDofsOffset3d(BackendData& data, const CeedInt nnodes, const CeedInt elem, const CeedInt* indices, const CeedScalar* r_v, CeedScalar* d_v) { 414 if (data.tidx < P1d && data.tidy < P1d) 415 for (CeedInt z = 0; z < P1d; ++z) { 416 const CeedInt node = data.tidx + data.tidy*P1d + z*P1d*P1d; 417 const CeedInt ind = indices[node + elem * P1d*P1d*P1d]; 418 for (CeedInt comp = 0; comp < NCOMP; ++comp) 419 atomicAdd(&d_v[ind + COMPSTRIDE * comp], r_v[z+comp*P1d]); 420 } 421 } 422 423 //------------------------------------------------------------------------------ 424 // E-vector -> L-vector, strided 425 //------------------------------------------------------------------------------ 426 template <int NCOMP, int P1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM> 427 inline __device__ void writeDofsStrided3d(BackendData& data, const CeedInt elem, const CeedScalar* r_v, CeedScalar* d_v) { 428 if (data.tidx < P1d && data.tidy < P1d) 429 for (CeedInt z = 0; z < P1d; ++z) { 430 const CeedInt node = data.tidx + data.tidy*P1d + z*P1d*P1d; 431 const CeedInt ind = node * STRIDES_NODE + elem * STRIDES_ELEM; 432 for (CeedInt comp = 0; comp < NCOMP; ++comp) 433 d_v[ind + comp * STRIDES_COMP] += r_v[z+comp*P1d]; 434 } 435 } 436 437 //------------------------------------------------------------------------------ 438 // 3D tensor contract x 439 //------------------------------------------------------------------------------ 440 template <int NCOMP, int P1d, int Q1d> 441 inline __device__ void ContractX3d(BackendData& data, const CeedScalar *U, const CeedScalar *B, CeedScalar *V) { 442 CeedScalar r_B[P1d]; 443 for (CeedInt i = 0; i < P1d; ++i) 444 r_B[i] = B[i + data.tidx*P1d]; 445 446 for (CeedInt k = 0; k < P1d; ++k) { 447 data.slice[data.tidx+data.tidy*T1d] = U[k]; 448 __syncthreads(); 449 V[k] = 0.0; 450 if (data.tidx < Q1d && data.tidy < P1d) 451 for (CeedInt i = 0; i < P1d; ++i) 452 V[k] += r_B[i] * data.slice[i + data.tidy*T1d]; // Contract x direction 453 __syncthreads(); 454 } 455 } 456 457 //------------------------------------------------------------------------------ 458 // 3D tensor contract y 459 //------------------------------------------------------------------------------ 460 template <int NCOMP, int P1d, int Q1d> 461 inline __device__ void ContractY3d(BackendData& data, const CeedScalar *U, const CeedScalar *B, CeedScalar *V) { 462 CeedScalar r_B[P1d]; 463 for (CeedInt i = 0; i < P1d; ++i) 464 r_B[i] = B[i + data.tidy*P1d]; 465 466 for (CeedInt k = 0; k < P1d; ++k) { 467 data.slice[data.tidx+data.tidy*T1d] = U[k]; 468 __syncthreads(); 469 V[k] = 0.0; 470 if (data.tidx < Q1d && data.tidy < Q1d) 471 for (CeedInt i = 0; i < P1d; ++i) 472 V[k] += r_B[i] * data.slice[data.tidx + i*T1d]; // Contract y direction 473 __syncthreads(); 474 } 475 } 476 477 //------------------------------------------------------------------------------ 478 // 3D tensor contract z 479 //------------------------------------------------------------------------------ 480 template <int NCOMP, int P1d, int Q1d> 481 inline __device__ void ContractZ3d(BackendData& data, const CeedScalar *U, const CeedScalar *B, CeedScalar *V) { 482 for (CeedInt k = 0; k < Q1d; ++k) { 483 V[k] = 0.0; 484 if (data.tidx < Q1d && data.tidy < Q1d) 485 for (CeedInt i = 0; i < P1d; ++i) 486 V[k] += B[i + k*P1d] * U[i]; // Contract z direction 487 } 488 } 489 490 //------------------------------------------------------------------------------ 491 // 3D transpose tensor contract z 492 //------------------------------------------------------------------------------ 493 template <int NCOMP, int P1d, int Q1d> 494 inline __device__ void ContractTransposeZ3d(BackendData& data, const CeedScalar *U, const CeedScalar *B, CeedScalar *V) { 495 for (CeedInt k = 0; k < P1d; ++k) { 496 V[k] = 0.0; 497 if (data.tidx < Q1d && data.tidy < Q1d) 498 for (CeedInt i = 0; i < Q1d; ++i) 499 V[k] += B[k + i*P1d] * U[i]; // Contract z direction 500 } 501 } 502 503 //------------------------------------------------------------------------------ 504 // 3D transpose tensor contract y 505 //------------------------------------------------------------------------------ 506 template <int NCOMP, int P1d, int Q1d> 507 inline __device__ void ContractTransposeY3d(BackendData& data, const CeedScalar *U, const CeedScalar *B, CeedScalar *V) { 508 CeedScalar r_B[Q1d]; 509 for (CeedInt i = 0; i < Q1d; ++i) 510 r_B[i] = B[data.tidy + i*P1d]; 511 512 for (CeedInt k = 0; k < P1d; ++k) { 513 data.slice[data.tidx+data.tidy*T1d] = U[k]; 514 __syncthreads(); 515 V[k] = 0.0; 516 if (data.tidx < Q1d && data.tidy < P1d) 517 for (CeedInt i = 0; i < Q1d; ++i) 518 V[k] += r_B[i] * data.slice[data.tidx + i*T1d]; // Contract y direction 519 __syncthreads(); 520 } 521 } 522 523 //------------------------------------------------------------------------------ 524 // 3D transpose tensor contract add y 525 //------------------------------------------------------------------------------ 526 template <int NCOMP, int P1d, int Q1d> 527 inline __device__ void ContractTransposeAddY3d(BackendData& data, const CeedScalar *U, const CeedScalar *B, CeedScalar *V) { 528 CeedScalar r_B[Q1d]; 529 for (CeedInt i = 0; i < Q1d; ++i) 530 r_B[i] = B[data.tidy + i*P1d]; 531 532 for (CeedInt k = 0; k < P1d; ++k) { 533 data.slice[data.tidx+data.tidy*T1d] = U[k]; 534 __syncthreads(); 535 if (data.tidx < Q1d && data.tidy < P1d) 536 for (CeedInt i = 0; i < Q1d; ++i) 537 V[k] += r_B[i] * data.slice[data.tidx + i*T1d]; // Contract y direction 538 __syncthreads(); 539 } 540 } 541 542 //------------------------------------------------------------------------------ 543 // 3D transpose tensor contract x 544 //------------------------------------------------------------------------------ 545 template <int NCOMP, int P1d, int Q1d> 546 inline __device__ void ContractTransposeX3d(BackendData& data, const CeedScalar *U, const CeedScalar *B, CeedScalar *V) { 547 CeedScalar r_B[Q1d]; 548 for (CeedInt i = 0; i < Q1d; ++i) 549 r_B[i] = B[data.tidx + i*P1d]; 550 551 for (CeedInt k = 0; k < P1d; ++k) { 552 data.slice[data.tidx+data.tidy*T1d] = U[k]; 553 __syncthreads(); 554 V[k] = 0.0; 555 if (data.tidx < P1d && data.tidy < P1d) 556 for (CeedInt i = 0; i < Q1d; ++i) 557 V[k] += r_B[i] * data.slice[i + data.tidy*T1d]; // Contract x direction 558 __syncthreads(); 559 } 560 } 561 562 //------------------------------------------------------------------------------ 563 // 3D transpose tensor contract add x 564 //------------------------------------------------------------------------------ 565 template <int NCOMP, int P1d, int Q1d> 566 inline __device__ void ContractTransposeAddX3d(BackendData& data, const CeedScalar *U, const CeedScalar *B, CeedScalar *V) { 567 CeedScalar r_B[Q1d]; 568 for (CeedInt i = 0; i < Q1d; ++i) 569 r_B[i] = B[data.tidx + i*P1d]; 570 571 for (CeedInt k = 0; k < P1d; ++k) { 572 data.slice[data.tidx+data.tidy*T1d] = U[k]; 573 __syncthreads(); 574 if (data.tidx < P1d && data.tidy < P1d) 575 for (CeedInt i = 0; i < Q1d; ++i) 576 V[k] += r_B[i] * data.slice[i + data.tidy*T1d]; // Contract x direction 577 __syncthreads(); 578 } 579 } 580 581 //------------------------------------------------------------------------------ 582 // 3D interpolate to quadrature points 583 //------------------------------------------------------------------------------ 584 template <int NCOMP, int P1d, int Q1d> 585 inline __device__ void interp3d(BackendData& data, const CeedScalar *__restrict__ r_U, const CeedScalar *c_B, CeedScalar *__restrict__ r_V) { 586 CeedScalar r_t1[T1d]; 587 CeedScalar r_t2[T1d]; 588 for (CeedInt comp = 0; comp < NCOMP; comp++) { 589 ContractX3d<NCOMP, P1d, Q1d>(data, r_U + comp*P1d, c_B, r_t1); 590 ContractY3d<NCOMP, P1d, Q1d>(data, r_t1, c_B, r_t2); 591 ContractZ3d<NCOMP, P1d, Q1d>(data, r_t2, c_B, r_V + comp*Q1d); 592 } 593 } 594 595 //------------------------------------------------------------------------------ 596 // 3D interpolate transpose 597 //------------------------------------------------------------------------------ 598 template <int NCOMP, int P1d, int Q1d> 599 inline __device__ void interpTranspose3d(BackendData& data, const CeedScalar *__restrict__ r_U, const CeedScalar *c_B, CeedScalar *__restrict__ r_V) { 600 CeedScalar r_t1[T1d]; 601 CeedScalar r_t2[T1d]; 602 for (CeedInt comp = 0; comp < NCOMP; comp++) { 603 ContractTransposeZ3d<NCOMP, P1d, Q1d>(data, r_U + comp*Q1d, c_B, r_t1); 604 ContractTransposeY3d<NCOMP, P1d, Q1d>(data, r_t1, c_B, r_t2); 605 ContractTransposeX3d<NCOMP, P1d, Q1d>(data, r_t2, c_B, r_V + comp*P1d); 606 } 607 } 608 609 //------------------------------------------------------------------------------ 610 // 3D derivatives at quadrature points 611 //------------------------------------------------------------------------------ 612 template <int NCOMP, int P1d, int Q1d> 613 inline __device__ void grad3d(BackendData& data, const CeedScalar *__restrict__ r_U, const CeedScalar *c_B, const CeedScalar *c_G, CeedScalar *__restrict__ r_V) { 614 CeedScalar r_t1[T1d]; 615 CeedScalar r_t2[T1d]; 616 for (CeedInt comp = 0; comp < NCOMP; comp++) { 617 ContractX3d<NCOMP, P1d, Q1d>(data, r_U + comp*P1d, c_G, r_t1); 618 ContractY3d<NCOMP, P1d, Q1d>(data, r_t1, c_B, r_t2); 619 ContractZ3d<NCOMP, P1d, Q1d>(data, r_t2, c_B, r_V + comp*Q1d + 0*NCOMP*Q1d); 620 ContractX3d<NCOMP, P1d, Q1d>(data, r_U + comp*P1d, c_B, r_t1); 621 ContractY3d<NCOMP, P1d, Q1d>(data, r_t1, c_G, r_t2); 622 ContractZ3d<NCOMP, P1d, Q1d>(data, r_t2, c_B, r_V + comp*Q1d + 1*NCOMP*Q1d); 623 ContractX3d<NCOMP, P1d, Q1d>(data, r_U + comp*P1d, c_B, r_t1); 624 ContractY3d<NCOMP, P1d, Q1d>(data, r_t1, c_B, r_t2); 625 ContractZ3d<NCOMP, P1d, Q1d>(data, r_t2, c_G, r_V + comp*Q1d + 2*NCOMP*Q1d); 626 } 627 } 628 629 //------------------------------------------------------------------------------ 630 // 3D derivatives transpose 631 //------------------------------------------------------------------------------ 632 template <int NCOMP, int P1d, int Q1d> 633 inline __device__ void gradTranspose3d(BackendData& data, const CeedScalar *__restrict__ r_U, const CeedScalar *c_B, const CeedScalar *c_G, CeedScalar *__restrict__ r_V) { 634 CeedScalar r_t1[T1d]; 635 CeedScalar r_t2[T1d]; 636 for (CeedInt comp = 0; comp < NCOMP; comp++) { 637 ContractTransposeZ3d<NCOMP, P1d, Q1d>(data, r_U + comp*Q1d + 0*NCOMP*Q1d, c_B, r_t1); 638 ContractTransposeY3d<NCOMP, P1d, Q1d>(data, r_t1, c_B, r_t2); 639 ContractTransposeX3d<NCOMP, P1d, Q1d>(data, r_t2, c_G, r_V + comp*P1d); 640 ContractTransposeZ3d<NCOMP, P1d, Q1d>(data, r_U + comp*Q1d + 1*NCOMP*Q1d, c_B, r_t1); 641 ContractTransposeY3d<NCOMP, P1d, Q1d>(data, r_t1, c_G, r_t2); 642 ContractTransposeAddX3d<NCOMP,P1d, Q1d>(data, r_t2, c_B, r_V + comp*P1d); 643 ContractTransposeZ3d<NCOMP, P1d, Q1d>(data, r_U + comp*Q1d + 2*NCOMP*Q1d, c_G, r_t1); 644 ContractTransposeY3d<NCOMP, P1d, Q1d>(data, r_t1, c_B, r_t2); 645 ContractTransposeAddX3d<NCOMP, P1d, Q1d>(data, r_t2, c_B, r_V + comp*P1d); 646 } 647 } 648 649 //------------------------------------------------------------------------------ 650 // 3D collocated derivatives computation 651 //------------------------------------------------------------------------------ 652 template <int NCOMP, int Q1d> 653 inline __device__ void gradCollo3d(BackendData& data, const CeedInt q, const CeedScalar *__restrict__ r_U, const CeedScalar *c_G, CeedScalar *__restrict__ r_V) { 654 if (data.tidx < Q1d && data.tidy < Q1d) { 655 for (CeedInt comp = 0; comp < NCOMP; ++comp) { 656 data.slice[data.tidx + data.tidy*T1d] = r_U[q + comp*Q1d]; 657 __syncthreads(); 658 // X derivative 659 r_V[comp+0*NCOMP] = 0.0; 660 for (CeedInt i = 0; i < Q1d; ++i) 661 r_V[comp+0*NCOMP] += c_G[i + data.tidx*Q1d] * data.slice[i + data.tidy*T1d]; // Contract x direction (X derivative) 662 // Y derivative 663 r_V[comp+1*NCOMP] = 0.0; 664 for (CeedInt i = 0; i < Q1d; ++i) 665 r_V[comp+1*NCOMP] += c_G[i + data.tidy*Q1d] * data.slice[data.tidx + i*T1d]; // Contract y direction (Y derivative) 666 // Z derivative 667 r_V[comp+2*NCOMP] = 0.0; 668 for (CeedInt i = 0; i < Q1d; ++i) 669 r_V[comp+2*NCOMP] += c_G[i + q*Q1d] * r_U[i + comp*Q1d]; // Contract z direction (Z derivative) 670 __syncthreads(); 671 } 672 } 673 } 674 675 //------------------------------------------------------------------------------ 676 // 3D collocated derivatives transpose 677 //------------------------------------------------------------------------------ 678 template <int NCOMP, int Q1d> 679 inline __device__ void gradColloTranspose3d(BackendData& data, const CeedInt q, const CeedScalar *__restrict__ r_U, const CeedScalar *c_G, CeedScalar *__restrict__ r_V) { 680 if (data.tidx < Q1d && data.tidy < Q1d) { 681 for (CeedInt comp = 0; comp < NCOMP; ++comp) { 682 // X derivative 683 data.slice[data.tidx + data.tidy*T1d] = r_U[comp + 0*NCOMP]; 684 __syncthreads(); 685 for (CeedInt i = 0; i < Q1d; ++i) 686 r_V[q+comp*Q1d] += c_G[data.tidx + i*Q1d] * data.slice[i + data.tidy*T1d]; // Contract x direction (X derivative) 687 __syncthreads(); 688 // Y derivative 689 data.slice[data.tidx + data.tidy*T1d] = r_U[comp + 1*NCOMP]; 690 __syncthreads(); 691 for (CeedInt i = 0; i < Q1d; ++i) 692 r_V[q+comp*Q1d] += c_G[data.tidy + i*Q1d] * data.slice[data.tidx + i*T1d]; // Contract y direction (Y derivative) 693 __syncthreads(); 694 // Z derivative 695 for (CeedInt i = 0; i < Q1d; ++i) 696 r_V[i+comp*Q1d] += c_G[i + q*Q1d] * r_U[comp + 2*NCOMP]; // PARTIAL contract z direction (Z derivative) 697 } 698 } 699 } 700 701 //------------------------------------------------------------------------------ 702 // 1D quadrature weights 703 //------------------------------------------------------------------------------ 704 template <int Q1d> 705 inline __device__ void weight1d(BackendData& data, const CeedScalar *qweight1d, CeedScalar *w) { 706 *w = (data.tidx < Q1d) ? qweight1d[data.tidx] : 0.0; 707 } 708 709 //------------------------------------------------------------------------------ 710 // 2D quadrature weights 711 //------------------------------------------------------------------------------ 712 template <int Q1d> 713 inline __device__ void weight2d(BackendData& data, const CeedScalar *qweight1d, CeedScalar *w) { 714 *w = (data.tidx < Q1d && data.tidy < Q1d) ? 715 qweight1d[data.tidx]*qweight1d[data.tidy] : 0.0; 716 } 717 718 //------------------------------------------------------------------------------ 719 // 3D quadrature weights 720 //------------------------------------------------------------------------------ 721 template <int Q1d> 722 inline __device__ void weight3d(BackendData& data, const CeedScalar *qweight1d, CeedScalar *w) { 723 const bool quad = (data.tidx < Q1d && data.tidy < Q1d); 724 const CeedScalar pw = quad ? qweight1d[data.tidx]*qweight1d[data.tidy] : 0.0; 725 for (CeedInt z = 0; z < Q1d; ++z) 726 w[z] = quad ? pw*qweight1d[z] : 0.0; 727 } 728 729 ); 730 //------------------------------------------------------------------------------ 731 // Build singe operator kernel 732 //------------------------------------------------------------------------------ 733 CEED_INTERN int CeedHipGenOperatorBuild(CeedOperator op) { 734 735 using std::ostringstream; 736 using std::string; 737 int ierr; 738 bool setupdone; 739 ierr = CeedOperatorIsSetupDone(op, &setupdone); CeedChk(ierr); 740 if (setupdone) return 0; 741 Ceed ceed; 742 ierr = CeedOperatorGetCeed(op, &ceed); CeedChk(ierr); 743 CeedOperator_Hip_gen *data; 744 ierr = CeedOperatorGetData(op, &data); CeedChk(ierr); 745 CeedQFunction qf; 746 CeedQFunction_Hip_gen *qf_data; 747 ierr = CeedOperatorGetQFunction(op, &qf); CeedChk(ierr); 748 ierr = CeedQFunctionGetData(qf, &qf_data); CeedChk(ierr); 749 CeedInt Q, P1d, Q1d = 0, numelements, elemsize, numinputfields, 750 numoutputfields, ncomp, dim = 0, lsize; 751 ierr = CeedOperatorGetNumQuadraturePoints(op, &Q); CeedChk(ierr); 752 ierr = CeedOperatorGetNumElements(op, &numelements); CeedChk(ierr); 753 ierr = CeedQFunctionGetNumArgs(qf, &numinputfields, &numoutputfields); 754 CeedChk(ierr); 755 CeedOperatorField *opinputfields, *opoutputfields; 756 ierr = CeedOperatorGetFields(op, &opinputfields, &opoutputfields); 757 CeedChk(ierr); 758 CeedQFunctionField *qfinputfields, *qfoutputfields; 759 ierr = CeedQFunctionGetFields(qf, &qfinputfields, &qfoutputfields); 760 CeedChk(ierr); 761 CeedEvalMode emode; 762 CeedBasis basis; 763 CeedBasis_Hip_shared *basis_data; 764 CeedElemRestriction Erestrict; 765 CeedElemRestriction_Hip *restr_data; 766 767 ostringstream code; 768 string devFunctions(deviceFunctions); 769 code << devFunctions; 770 771 string qFunction(qf_data->qFunctionSource); 772 string qFunctionName(qf_data->qFunctionName); 773 string oper; 774 oper = "CeedKernel_Hip_gen_" + qFunctionName; 775 776 code << "\n#define CEED_QFUNCTION(name) inline __device__ int name\n"; 777 code << "\n#define CeedPragmaSIMD\n"; 778 779 // Find dim and Q1d 780 bool useCollograd = true; 781 data->maxP1d = 0; 782 for (CeedInt i = 0; i < numinputfields; i++) { 783 ierr = CeedOperatorFieldGetBasis(opinputfields[i], &basis); CeedChk(ierr); 784 if (basis != CEED_BASIS_COLLOCATED) { 785 ierr = CeedBasisGetData(basis, &basis_data); CeedChk(ierr); 786 ierr = CeedQFunctionFieldGetEvalMode(qfinputfields[i], &emode); 787 CeedChk(ierr); 788 789 // Check for collocated gradient 790 useCollograd = useCollograd && basis_data->d_collograd1d; 791 792 // Collect dim and Q1d 793 ierr = CeedBasisGetDimension(basis, &dim); CeedChk(ierr); 794 bool isTensor; 795 ierr = CeedBasisIsTensor(basis, &isTensor); CeedChk(ierr); 796 if (isTensor) { 797 ierr = CeedBasisGetNumQuadraturePoints1D(basis, &Q1d); CeedChk(ierr); 798 ierr = CeedBasisGetNumNodes1D(basis, &P1d); CeedChk(ierr); 799 if (P1d>data->maxP1d) data->maxP1d = P1d; 800 } else { 801 // LCOV_EXCL_START 802 return CeedError(ceed, 1, "Backend does not implement operators with non-tensor basis"); 803 // LCOV_EXCL_STOP 804 } 805 } 806 } 807 // Check output bases for Q1d, dim as well 808 // The only imput basis might be CEED_BASIS_COLLOCATED 809 for (CeedInt i = 0; i < numoutputfields; i++) { 810 ierr = CeedOperatorFieldGetBasis(opoutputfields[i], &basis); CeedChk(ierr); 811 812 if (basis != CEED_BASIS_COLLOCATED) { 813 ierr = CeedBasisGetData(basis, &basis_data); CeedChk(ierr); 814 ierr = CeedQFunctionFieldGetEvalMode(qfoutputfields[i], &emode); 815 CeedChk(ierr); 816 817 // Collect dim and Q1d 818 ierr = CeedBasisGetDimension(basis, &dim); CeedChk(ierr); 819 bool isTensor; 820 ierr = CeedBasisIsTensor(basis, &isTensor); CeedChk(ierr); 821 if (isTensor) { 822 ierr = CeedBasisGetNumQuadraturePoints1D(basis, &Q1d); CeedChk(ierr); 823 } else { 824 // LCOV_EXCL_START 825 return CeedError(ceed, 1, "Backend does not implement operators with non-tensor basis"); 826 // LCOV_EXCL_STOP 827 } 828 829 // Check for collocated gradient 830 useCollograd = useCollograd && basis_data->d_collograd1d; 831 } 832 } 833 data->dim = dim; 834 data->Q1d = Q1d; 835 836 // Define CEED_Q_VLA 837 if (dim != 3 || useCollograd) { 838 code << "\n#define CEED_Q_VLA 1\n\n"; 839 } else { 840 code << "\n#define CEED_Q_VLA "<<Q1d<<"\n\n"; 841 } 842 843 code << qFunction; 844 845 // Setup 846 code << "\n// -----------------------------------------------------------------------------\n"; 847 code << "\nextern \"C\" __global__ void "<<oper<<"(CeedInt nelem, void* ctx, HipFieldsInt indices, HipFields fields, HipFields B, HipFields G, CeedScalar* W) {\n"; 848 for (CeedInt i = 0; i < numinputfields; i++) { 849 ierr = CeedQFunctionFieldGetEvalMode(qfinputfields[i], &emode); 850 CeedChk(ierr); 851 if (emode != CEED_EVAL_WEIGHT) { // Skip CEED_EVAL_WEIGHT 852 code << " const CeedScalar* d_u" <<i<<" = fields.in["<<i<<"];\n"; 853 } 854 } 855 856 for (CeedInt i = 0; i < numoutputfields; i++) { 857 code << " CeedScalar* d_v"<<i<<" = fields.out["<<i<<"];\n"; 858 } 859 860 code << " const CeedInt Dim = "<<dim<<";\n"; 861 code << " const CeedInt Q1d = "<<Q1d<<";\n"; 862 863 code << " HIP_DYNAMIC_SHARED( CeedScalar, slice)\n"; 864 code << " BackendData data;\n"; 865 code << " data.tidx = threadIdx.x;\n"; 866 code << " data.tidy = threadIdx.y;\n"; 867 code << " data.tidz = threadIdx.z;\n"; 868 code << " data.tid = threadIdx.x + threadIdx.y*blockDim.x + threadIdx.z*blockDim.y*blockDim.x;\n"; 869 code << " data.slice = slice+data.tidz*T1d"<<(dim>1?"*T1d":"")<<";\n"; 870 871 code << "\n // -- Input field constants and basis data --\n"; 872 //Initialize constants, and matrices B and G 873 for (CeedInt i = 0; i < numinputfields; i++) { 874 code << " // ---- Input field "<<i<<" ----\n"; 875 // Get elemsize, emode, ncomp 876 ierr = CeedOperatorFieldGetElemRestriction(opinputfields[i], &Erestrict); 877 CeedChk(ierr); 878 ierr = CeedElemRestrictionGetElementSize(Erestrict, &elemsize); 879 CeedChk(ierr); 880 ierr = CeedQFunctionFieldGetEvalMode(qfinputfields[i], &emode); 881 CeedChk(ierr); 882 ierr = CeedElemRestrictionGetNumComponents(Erestrict, &ncomp); 883 CeedChk(ierr); 884 885 // Set field constants 886 if (emode != CEED_EVAL_WEIGHT) { 887 ierr = CeedOperatorFieldGetBasis(opinputfields[i], &basis); CeedChk(ierr); 888 if (basis != CEED_BASIS_COLLOCATED) { 889 ierr = CeedBasisGetNumNodes1D(basis, &P1d); CeedChk(ierr); 890 code << " const CeedInt P_in_"<<i<<" = "<<P1d<<";\n"; 891 } else { 892 code << " const CeedInt P_in_"<<i<<" = "<<Q1d<<";\n"; 893 } 894 code << " const CeedInt ncomp_in_"<<i<<" = "<<ncomp<<";\n"; 895 } 896 897 // Load basis data 898 code << " // EvalMode: "<<CeedEvalModes[emode]<<"\n"; 899 switch (emode) { 900 case CEED_EVAL_NONE: 901 break; 902 case CEED_EVAL_INTERP: 903 ierr = CeedBasisGetData(basis, &basis_data); CeedChk(ierr); 904 data->B.in[i] = basis_data->d_interp1d; 905 code << " __shared__ double s_B_in_"<<i<<"["<<P1d*Q1d<<"];\n"; 906 code << " loadMatrix<P_in_"<<i<<",Q1d>(data, B.in["<<i<<"], s_B_in_"<<i<<");\n"; 907 break; 908 case CEED_EVAL_GRAD: 909 ierr = CeedBasisGetData(basis, &basis_data); CeedChk(ierr); 910 data->B.in[i] = basis_data->d_interp1d; 911 code << " __shared__ double s_B_in_"<<i<<"["<<P1d*Q1d<<"];\n"; 912 code << " loadMatrix<P_in_"<<i<<",Q1d>(data, B.in["<<i<<"], s_B_in_"<<i<<");\n"; 913 if (useCollograd) { 914 data->G.in[i] = basis_data->d_collograd1d; 915 code << " __shared__ double s_G_in_"<<i<<"["<<Q1d*Q1d<<"];\n"; 916 code << " loadMatrix<Q1d,Q1d>(data, G.in["<<i<<"], s_G_in_"<<i<<");\n"; 917 } else { 918 data->G.in[i] = basis_data->d_grad1d; 919 code << " __shared__ double s_G_in_"<<i<<"["<<P1d*Q1d<<"];\n"; 920 code << " loadMatrix<P_in_"<<i<<",Q1d>(data, G.in["<<i<<"], s_G_in_"<<i<<");\n"; 921 } 922 break; 923 case CEED_EVAL_WEIGHT: 924 break; // No action 925 case CEED_EVAL_DIV: 926 break; // TODO: Not implemented 927 case CEED_EVAL_CURL: 928 break; // TODO: Not implemented 929 } 930 } 931 932 code << "\n // -- Output field constants and basis data --\n"; 933 for (CeedInt i = 0; i < numoutputfields; i++) { 934 code << " // ---- Output field "<<i<<" ----\n"; 935 // Get elemsize, emode, ncomp 936 ierr = CeedOperatorFieldGetElemRestriction(opoutputfields[i], &Erestrict); 937 CeedChk(ierr); 938 ierr = CeedElemRestrictionGetElementSize(Erestrict, &elemsize); 939 CeedChk(ierr); 940 ierr = CeedQFunctionFieldGetEvalMode(qfoutputfields[i], &emode); 941 CeedChk(ierr); 942 ierr = CeedElemRestrictionGetNumComponents(Erestrict, &ncomp); 943 CeedChk(ierr); 944 945 // Set field constants 946 ierr = CeedOperatorFieldGetBasis(opoutputfields[i], &basis); CeedChk(ierr); 947 if (basis != CEED_BASIS_COLLOCATED) { 948 ierr = CeedBasisGetNumNodes1D(basis, &P1d); CeedChk(ierr); 949 code << " const CeedInt P_out_"<<i<<" = "<<P1d<<";\n"; 950 } else { 951 code << " const CeedInt P_out_"<<i<<" = "<<Q1d<<";\n"; 952 } 953 code << " const CeedInt ncomp_out_"<<i<<" = "<<ncomp<<";\n"; 954 955 // Load basis data 956 code << " // EvalMode: "<<CeedEvalModes[emode]<<"\n"; 957 switch (emode) { 958 case CEED_EVAL_NONE: 959 break; // No action 960 case CEED_EVAL_INTERP: 961 ierr = CeedBasisGetData(basis, &basis_data); CeedChk(ierr); 962 data->B.out[i] = basis_data->d_interp1d; 963 code << " __shared__ double s_B_out_"<<i<<"["<<P1d*Q1d<<"];\n"; 964 code << " loadMatrix<P_out_"<<i<<",Q1d>(data, B.out["<<i<<"], s_B_out_"<<i<<");\n"; 965 break; 966 case CEED_EVAL_GRAD: 967 ierr = CeedBasisGetData(basis, &basis_data); CeedChk(ierr); 968 data->B.out[i] = basis_data->d_interp1d; 969 code << " __shared__ double s_B_out_"<<i<<"["<<P1d*Q1d<<"];\n"; 970 code << " loadMatrix<P_out_"<<i<<",Q1d>(data, B.out["<<i<<"], s_B_out_"<<i<<");\n"; 971 if (useCollograd) { 972 data->G.out[i] = basis_data->d_collograd1d; 973 code << " __shared__ double s_G_out_"<<i<<"["<<Q1d*Q1d<<"];\n"; 974 code << " loadMatrix<Q1d,Q1d>(data, G.out["<<i<<"], s_G_out_"<<i<<");\n"; 975 } else { 976 data->G.out[i] = basis_data->d_grad1d; 977 code << " __shared__ double s_G_out_"<<i<<"["<<P1d*Q1d<<"];\n"; 978 code << " loadMatrix<P_out_"<<i<<",Q1d>(data, G.out["<<i<<"], s_G_out_"<<i<<");\n"; 979 } 980 break; 981 // LCOV_EXCL_START 982 case CEED_EVAL_WEIGHT: { 983 Ceed ceed; 984 ierr = CeedOperatorGetCeed(op, &ceed); CeedChk(ierr); 985 return CeedError(ceed, 1, 986 "CEED_EVAL_WEIGHT cannot be an output evaluation mode"); 987 break; // Should not occur 988 } 989 case CEED_EVAL_DIV: 990 break; // TODO: Not implemented 991 case CEED_EVAL_CURL: 992 break; // TODO: Not implemented 993 // LCOV_EXCL_STOP 994 } 995 } 996 code << "\n // -- Element loop --\n"; 997 code << " __syncthreads();\n"; 998 code << " for (CeedInt elem = blockIdx.x*blockDim.z + threadIdx.z; elem < nelem; elem += gridDim.x*blockDim.z) {\n"; 999 // Input basis apply if needed 1000 // Generate the correct eval mode code for each input 1001 code << " // -- Input field restrictions and basis actions --\n"; 1002 for (CeedInt i = 0; i < numinputfields; i++) { 1003 code << " // ---- Input field "<<i<<" ----\n"; 1004 // Get elemsize, emode, ncomp 1005 ierr = CeedOperatorFieldGetElemRestriction(opinputfields[i], &Erestrict); 1006 CeedChk(ierr); 1007 ierr = CeedElemRestrictionGetElementSize(Erestrict, &elemsize); 1008 CeedChk(ierr); 1009 ierr = CeedQFunctionFieldGetEvalMode(qfinputfields[i], &emode); 1010 CeedChk(ierr); 1011 ierr = CeedElemRestrictionGetNumComponents(Erestrict, &ncomp); 1012 CeedChk(ierr); 1013 1014 // Restriction 1015 if (emode != CEED_EVAL_WEIGHT && 1016 !((emode == CEED_EVAL_NONE) && useCollograd)) { 1017 code << " CeedScalar r_u"<<i<<"[ncomp_in_"<<i<<"*P_in_"<<i<<"];\n"; 1018 1019 bool isStrided; 1020 ierr = CeedElemRestrictionIsStrided(Erestrict, &isStrided); CeedChk(ierr); 1021 if (!isStrided) { 1022 ierr = CeedElemRestrictionGetLVectorSize(Erestrict, &lsize); 1023 CeedChk(ierr); 1024 code << " const CeedInt lsize_in_"<<i<<" = "<<lsize<<";\n"; 1025 CeedInt compstride; 1026 ierr = CeedElemRestrictionGetCompStride(Erestrict, &compstride); CeedChk(ierr); 1027 code << " // CompStride: "<<compstride<<"\n"; 1028 ierr = CeedElemRestrictionGetData(Erestrict, &restr_data); CeedChk(ierr); 1029 data->indices.in[i] = restr_data->d_ind; 1030 code << " readDofsOffset"<<dim<<"d<ncomp_in_"<<i<<", "<<compstride<<", P_in_"<<i<<">(data, lsize_in_"<<i<<", elem, indices.in["<<i<<"], d_u"<<i<<", r_u"<<i<<");\n"; 1031 } else { 1032 bool backendstrides; 1033 ierr = CeedElemRestrictionHasBackendStrides(Erestrict, &backendstrides); 1034 CeedChk(ierr); 1035 CeedInt nelem; 1036 ierr = CeedElemRestrictionGetNumElements(Erestrict, &nelem); 1037 CeedChk(ierr); 1038 CeedInt strides[3] = {1, elemsize*nelem, elemsize}; 1039 if (!backendstrides) { 1040 ierr = CeedElemRestrictionGetStrides(Erestrict, &strides); 1041 CeedChk(ierr); 1042 } 1043 code << " // Strides: {"<<strides[0]<<", "<<strides[1]<<", "<<strides[2]<<"}\n"; 1044 code << " readDofsStrided"<<dim<<"d<ncomp_in_"<<i<<",P_in_"<<i<<","<<strides[0]<<","<<strides[1]<<","<<strides[2]<<">(data, elem, d_u"<<i<<", r_u"<<i<<");\n"; 1045 } 1046 } 1047 1048 // Basis action 1049 code << " // EvalMode: "<<CeedEvalModes[emode]<<"\n"; 1050 switch (emode) { 1051 case CEED_EVAL_NONE: 1052 if (!useCollograd) { 1053 code << " CeedScalar* r_t"<<i<<" = r_u"<<i<<";\n"; 1054 } 1055 break; 1056 case CEED_EVAL_INTERP: 1057 code << " CeedScalar r_t"<<i<<"[ncomp_in_"<<i<<"*Q1d];\n"; 1058 code << " interp"<<dim<<"d<ncomp_in_"<<i<<",P_in_"<<i<<",Q1d>(data, r_u"<<i<<", s_B_in_"<<i<<", r_t"<<i<<");\n"; 1059 break; 1060 case CEED_EVAL_GRAD: 1061 if (useCollograd) { 1062 code << " CeedScalar r_t"<<i<<"[ncomp_in_"<<i<<"*Q1d];\n"; 1063 code << " interp"<<dim<<"d<ncomp_in_"<<i<<",P_in_"<<i<<",Q1d>(data, r_u"<<i<<", s_B_in_"<<i<<", r_t"<<i<<");\n"; 1064 } else { 1065 code << " CeedScalar r_t"<<i<<"[ncomp_in_"<<i<<"*Dim*Q1d];\n"; 1066 code << " grad"<<dim<<"d<ncomp_in_"<<i<<",P_in_"<<i<<",Q1d>(data, r_u"<<i<<", s_B_in_"<<i<<", s_G_in_"<<i<<", r_t"<<i<<");\n"; 1067 } 1068 break; 1069 case CEED_EVAL_WEIGHT: 1070 code << " CeedScalar r_t"<<i<<"[Q1d];\n"; 1071 ierr = CeedOperatorFieldGetBasis(opinputfields[i], &basis); CeedChk(ierr); 1072 ierr = CeedBasisGetData(basis, &basis_data); CeedChk(ierr); 1073 data->W = basis_data->d_qweight1d; 1074 code << " weight"<<dim<<"d<Q1d>(data, W, r_t"<<i<<");\n"; 1075 break; // No action 1076 case CEED_EVAL_DIV: 1077 break; // TODO: Not implemented 1078 case CEED_EVAL_CURL: 1079 break; // TODO: Not implemented 1080 } 1081 } 1082 1083 // Q function 1084 code << "\n // -- Output field setup --\n"; 1085 for (CeedInt i = 0; i < numoutputfields; i++) { 1086 code << "\n // ---- Output field "<<i<<" ----\n"; 1087 ierr = CeedQFunctionFieldGetEvalMode(qfoutputfields[i], &emode); 1088 CeedChk(ierr); 1089 if (emode==CEED_EVAL_GRAD) 1090 { 1091 if (useCollograd) { 1092 //Accumulator for gradient slices 1093 code << " CeedScalar r_tt"<<i<<"[ncomp_out_"<<i<<"*Q1d];\n"; 1094 code << " for (CeedInt i = 0; i < ncomp_out_"<<i<<"; ++i) {\n"; 1095 code << " for (CeedInt j = 0; j < Q1d; ++j) {\n"; 1096 code << " r_tt"<<i<<"[j + i*Q1d] = 0.0;\n"; 1097 code << " }\n"; 1098 code << " }\n"; 1099 } else { 1100 code << " CeedScalar r_tt"<<i<<"[ncomp_out_"<<i<<"*Dim*Q1d];\n"; 1101 } 1102 } 1103 if (emode==CEED_EVAL_NONE || emode==CEED_EVAL_INTERP) 1104 { 1105 code << " CeedScalar r_tt"<<i<<"[ncomp_out_"<<i<<"*Q1d];\n"; 1106 } 1107 } 1108 // We treat quadrature points per slice in 3d to save registers 1109 if (useCollograd) { 1110 code << "\n // Note: Collocated Gradient\n"; 1111 code << "#pragma unroll\n"; 1112 code << " for (CeedInt q=0; q<Q1d; q++) {\n"; 1113 code << " // -- Input fields --\n"; 1114 for (CeedInt i = 0; i < numinputfields; i++) { 1115 code << " // ---- Input field "<<i<<" ----\n"; 1116 // Get elemsize, emode, ncomp 1117 ierr = CeedQFunctionFieldGetEvalMode(qfinputfields[i], &emode); 1118 CeedChk(ierr); 1119 // Basis action 1120 code << " // EvalMode: "<<CeedEvalModes[emode]<<"\n"; 1121 switch (emode) { 1122 case CEED_EVAL_NONE: 1123 code << " CeedScalar r_q"<<i<<"[ncomp_in_"<<i<<"];\n"; 1124 1125 bool isStrided; 1126 ierr = CeedOperatorFieldGetElemRestriction(opinputfields[i], &Erestrict); CeedChk(ierr); 1127 ierr = CeedElemRestrictionGetElementSize(Erestrict, &elemsize); CeedChk(ierr); 1128 ierr = CeedElemRestrictionIsStrided(Erestrict, &isStrided); CeedChk(ierr); 1129 if (!isStrided) { 1130 ierr = CeedElemRestrictionGetLVectorSize(Erestrict, &lsize); 1131 CeedChk(ierr); 1132 code << " const CeedInt lsize_in_"<<i<<" = "<<lsize<<";\n"; 1133 CeedInt compstride; 1134 ierr = CeedElemRestrictionGetCompStride(Erestrict, &compstride); CeedChk(ierr); 1135 code << " // CompStride: "<<compstride<<"\n"; 1136 ierr = CeedElemRestrictionGetData(Erestrict, &restr_data); CeedChk(ierr); 1137 data->indices.in[i] = restr_data->d_ind; 1138 code << " readSliceQuadsOffset"<<"3d<ncomp_in_"<<i<<", "<<compstride<<", Q1d>(data, lsize_in_"<<i<<", elem, q, indices.in["<<i<<"], d_u"<<i<<", r_q"<<i<<");\n"; 1139 } else { 1140 bool backendstrides; 1141 ierr = CeedElemRestrictionHasBackendStrides(Erestrict, &backendstrides); 1142 CeedChk(ierr); 1143 CeedInt nelem; 1144 ierr = CeedElemRestrictionGetNumElements(Erestrict, &nelem); 1145 CeedChk(ierr); 1146 CeedInt strides[3] = {1, elemsize*nelem, elemsize}; 1147 if (!backendstrides) { 1148 ierr = CeedElemRestrictionGetStrides(Erestrict, &strides); 1149 CeedChk(ierr); 1150 } 1151 code << " // Strides: {"<<strides[0]<<", "<<strides[1]<<", "<<strides[2]<<"}\n"; 1152 code << " readSliceQuadsStrided"<<"3d<ncomp_in_"<<i<<",Q1d"","<<strides[0]<<","<<strides[1]<<","<<strides[2]<<">(data, elem, q, d_u"<<i<<", r_q"<<i<<");\n"; 1153 } 1154 break; 1155 case CEED_EVAL_INTERP: 1156 code << " CeedScalar r_q"<<i<<"[ncomp_in_"<<i<<"];\n"; 1157 code << " for (CeedInt j = 0; j < ncomp_in_"<<i<<" ; ++j) {\n"; 1158 code << " r_q"<<i<<"[j] = r_t"<<i<<"[q + j*Q1d];\n"; 1159 code << " }\n"; 1160 break; 1161 case CEED_EVAL_GRAD: 1162 code << " CeedScalar r_q"<<i<<"[ncomp_in_"<<i<<"*Dim];\n"; 1163 code << " gradCollo3d<ncomp_in_"<<i<<",Q1d>(data, q, r_t"<<i<<", s_G_in_"<<i<<", r_q"<<i<<");\n"; 1164 break; 1165 case CEED_EVAL_WEIGHT: 1166 code << " CeedScalar r_q"<<i<<"[1];\n"; 1167 code << " r_q"<<i<<"[0] = r_t"<<i<<"[q];\n"; 1168 break; // No action 1169 case CEED_EVAL_DIV: 1170 break; // TODO: Not implemented 1171 case CEED_EVAL_CURL: 1172 break; // TODO: Not implemented 1173 } 1174 } 1175 code << "\n // -- Output fields --\n"; 1176 for (CeedInt i = 0; i < numoutputfields; i++) { 1177 code << " // ---- Output field "<<i<<" ----\n"; 1178 ierr = CeedQFunctionFieldGetEvalMode(qfoutputfields[i], &emode); 1179 CeedChk(ierr); 1180 // Basis action 1181 switch (emode) { 1182 case CEED_EVAL_NONE: 1183 code << " CeedScalar r_qq"<<i<<"[ncomp_out_"<<i<<"];\n"; 1184 break; // No action 1185 case CEED_EVAL_INTERP: 1186 code << " CeedScalar r_qq"<<i<<"[ncomp_out_"<<i<<"];\n"; 1187 break; 1188 case CEED_EVAL_GRAD: 1189 code << " CeedScalar r_qq"<<i<<"[ncomp_out_"<<i<<"*Dim];\n"; 1190 break; 1191 case CEED_EVAL_WEIGHT: 1192 break; // Should not occur 1193 case CEED_EVAL_DIV: 1194 break; // TODO: Not implemented 1195 case CEED_EVAL_CURL: 1196 break; // TODO: Not implemented 1197 } 1198 } 1199 } else { 1200 code << "\n // Note: No Collocated Gradient\n"; 1201 code << " // -- Input fields --\n"; 1202 for (CeedInt i = 0; i < numinputfields; i++) { 1203 code << " // ---- Input field "<<i<<" ----\n"; 1204 code << " CeedScalar* r_q"<<i<<" = r_t"<<i<<";\n"; 1205 } 1206 code << " // -- Output fields --\n"; 1207 for (CeedInt i = 0; i < numoutputfields; i++) { 1208 code << " // ---- Output field "<<i<<" ----\n"; 1209 code << " CeedScalar* r_qq"<<i<<" = r_tt"<<i<<";\n"; 1210 } 1211 } 1212 code << "\n // -- QFunction Inputs and outputs --\n"; 1213 code << " CeedScalar* in["<<numinputfields<<"];\n"; 1214 for (CeedInt i = 0; i < numinputfields; i++) { 1215 code << " // ---- Input field "<<i<<" ----\n"; 1216 code << " in["<<i<<"] = r_q"<<i<<";\n"; 1217 } 1218 code << " CeedScalar* out["<<numoutputfields<<"];\n"; 1219 for (CeedInt i = 0; i < numoutputfields; i++) { 1220 code << " // ---- Output field "<<i<<" ----\n"; 1221 code << " out["<<i<<"] = r_qq"<<i<<";\n"; 1222 } 1223 code << "\n // -- Apply QFunction --\n"; 1224 code << " "<<qFunctionName<<"(ctx, "; 1225 if (dim != 3 || useCollograd) { 1226 code << "1"; 1227 } else { 1228 code << "Q1d"; 1229 } 1230 code << ", in, out);\n"; 1231 if (useCollograd) { 1232 code << "\n // Note: Collocated Gradient\n"; 1233 code << " // -- Output fields --\n"; 1234 for (CeedInt i = 0; i < numoutputfields; i++) { 1235 code << " // ---- Output field "<<i<<" ----\n"; 1236 ierr = CeedQFunctionFieldGetEvalMode(qfoutputfields[i], &emode); 1237 CeedChk(ierr); 1238 // Basis action 1239 code << " // EvalMode: "<<CeedEvalModes[emode]<<"\n"; 1240 switch (emode) { 1241 case CEED_EVAL_NONE: 1242 code << " for (CeedInt j = 0; j < ncomp_out_"<<i<<" ; ++j) {\n"; 1243 code << " r_tt"<<i<<"[q + j*Q1d] = r_qq"<<i<<"[j];\n"; 1244 code << " }\n"; 1245 break; // No action 1246 case CEED_EVAL_INTERP: 1247 code << " for (CeedInt j = 0; j < ncomp_out_"<<i<<" ; ++j) {\n"; 1248 code << " r_tt"<<i<<"[q + j*Q1d] = r_qq"<<i<<"[j];\n"; 1249 code << " }\n"; 1250 break; 1251 case CEED_EVAL_GRAD: 1252 code << " gradColloTranspose3d<ncomp_out_"<<i<<",Q1d>(data, q, r_qq"<<i<<", s_G_out_"<<i<<", r_tt"<<i<<");\n"; 1253 break; 1254 case CEED_EVAL_WEIGHT: 1255 break; // Should not occur 1256 case CEED_EVAL_DIV: 1257 break; // TODO: Not implemented 1258 case CEED_EVAL_CURL: 1259 break; // TODO: Not implemented 1260 } 1261 } 1262 code << " }\n"; 1263 } 1264 1265 // Output basis apply if needed 1266 // Generate the correct eval mode code for each output 1267 code << "\n // -- Output field basis action and restrictions --\n"; 1268 for (CeedInt i = 0; i < numoutputfields; i++) { 1269 code << " // ---- Output field "<<i<<" ----\n"; 1270 // Get elemsize, emode, ncomp 1271 ierr = CeedOperatorFieldGetElemRestriction(opoutputfields[i], &Erestrict); 1272 CeedChk(ierr); 1273 ierr = CeedElemRestrictionGetElementSize(Erestrict, &elemsize); 1274 CeedChk(ierr); 1275 ierr = CeedQFunctionFieldGetEvalMode(qfoutputfields[i], &emode); 1276 CeedChk(ierr); 1277 ierr = CeedElemRestrictionGetNumComponents(Erestrict, &ncomp); 1278 CeedChk(ierr); 1279 // Basis action 1280 code << " // EvalMode: "<<CeedEvalModes[emode]<<"\n"; 1281 switch (emode) { 1282 case CEED_EVAL_NONE: 1283 code << " CeedScalar* r_v"<<i<<" = r_tt"<<i<<";\n"; 1284 break; // No action 1285 case CEED_EVAL_INTERP: 1286 code << " CeedScalar r_v"<<i<<"[ncomp_out_"<<i<<"*P_out_"<<i<<"];\n"; 1287 code << " interpTranspose"<<dim<<"d<ncomp_out_"<<i<<",P_out_"<<i<<",Q1d>(data, r_tt"<<i<<", s_B_out_"<<i<<", r_v"<<i<<");\n"; 1288 break; 1289 case CEED_EVAL_GRAD: 1290 code << " CeedScalar r_v"<<i<<"[ncomp_out_"<<i<<"*P_out_"<<i<<"];\n"; 1291 if (useCollograd) { 1292 code << " interpTranspose"<<dim<<"d<ncomp_out_"<<i<<",P_out_"<<i<<",Q1d>(data, r_tt"<<i<<", s_B_out_"<<i<<", r_v"<<i<<");\n"; 1293 } else { 1294 code << " gradTranspose"<<dim<<"d<ncomp_out_"<<i<<",P_out_"<<i<<",Q1d>(data, r_tt"<<i<<", s_B_out_"<<i<<", s_G_out_"<<i<<", r_v"<<i<<");\n"; 1295 } 1296 break; 1297 // LCOV_EXCL_START 1298 case CEED_EVAL_WEIGHT: { 1299 Ceed ceed; 1300 ierr = CeedOperatorGetCeed(op, &ceed); CeedChk(ierr); 1301 return CeedError(ceed, 1, 1302 "CEED_EVAL_WEIGHT cannot be an output evaluation mode"); 1303 break; // Should not occur 1304 } 1305 case CEED_EVAL_DIV: 1306 break; // TODO: Not implemented 1307 case CEED_EVAL_CURL: 1308 break; // TODO: Not implemented 1309 // LCOV_EXCL_STOP 1310 } 1311 // Restriction 1312 bool isStrided; 1313 ierr = CeedElemRestrictionIsStrided(Erestrict, &isStrided); CeedChk(ierr); 1314 if (!isStrided) { 1315 ierr = CeedElemRestrictionGetLVectorSize(Erestrict, &lsize); 1316 CeedChk(ierr); 1317 code << " const CeedInt lsize_out_"<<i<<" = "<<lsize<<";\n"; 1318 CeedInt compstride; 1319 ierr = CeedElemRestrictionGetCompStride(Erestrict, &compstride); CeedChk(ierr); 1320 code << " // CompStride: "<<compstride<<"\n"; 1321 ierr = CeedElemRestrictionGetData(Erestrict, &restr_data); CeedChk(ierr); 1322 data->indices.out[i] = restr_data->d_ind; 1323 code << " writeDofsOffset"<<dim<<"d<ncomp_out_"<<i<<", "<<compstride<<", P_out_"<<i<<">(data, lsize_out_"<<i<<", elem, indices.out["<<i<<"], r_v"<<i<<", d_v"<<i<<");\n"; 1324 } else { 1325 bool backendstrides; 1326 ierr = CeedElemRestrictionHasBackendStrides(Erestrict, &backendstrides); 1327 CeedChk(ierr); 1328 CeedInt nelem; 1329 ierr = CeedElemRestrictionGetNumElements(Erestrict, &nelem); 1330 CeedChk(ierr); 1331 CeedInt strides[3] = {1, elemsize*nelem, elemsize}; 1332 if (!backendstrides) { 1333 ierr = CeedElemRestrictionGetStrides(Erestrict, &strides); 1334 CeedChk(ierr); 1335 } 1336 code << " // Strides: {"<<strides[0]<<", "<<strides[1]<<", "<<strides[2]<<"}\n"; 1337 code << " writeDofsStrided"<<dim<<"d<ncomp_out_"<<i<<",P_out_"<<i<<","<<strides[0]<<","<<strides[1]<<","<<strides[2]<<">(data, elem, r_v"<<i<<", d_v"<<i<<");\n"; 1338 } 1339 } 1340 1341 code << " }\n"; 1342 code << "}\n"; 1343 code << "// -----------------------------------------------------------------------------\n\n"; 1344 1345 // View kernel for debugging 1346 CeedDebug(code.str().c_str()); 1347 1348 ierr = CeedCompileHip(ceed, code.str().c_str(), &data->module, 1, 1349 "T1d", CeedIntMax(Q1d, data->maxP1d)); 1350 CeedChk(ierr); 1351 ierr = CeedGetKernelHip(ceed, data->module, oper.c_str(), &data->op); 1352 CeedChk(ierr); 1353 1354 ierr = CeedOperatorSetSetupDone(op); CeedChk(ierr); 1355 return 0; 1356 } 1357 //------------------------------------------------------------------------------ 1358