xref: /libCEED/backends/hip-gen/ceed-hip-gen-operator-build.cpp (revision 7eaaafceb80f5cb5ccd1fd80420c1c1d8166ddd6)
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/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); CeedChkBackend(ierr);
740   if (setupdone) return CEED_ERROR_SUCCESS;
741   Ceed ceed;
742   ierr = CeedOperatorGetCeed(op, &ceed); CeedChkBackend(ierr);
743   CeedOperator_Hip_gen *data;
744   ierr = CeedOperatorGetData(op, &data); CeedChkBackend(ierr);
745   CeedQFunction qf;
746   CeedQFunction_Hip_gen *qf_data;
747   ierr = CeedOperatorGetQFunction(op, &qf); CeedChkBackend(ierr);
748   ierr = CeedQFunctionGetData(qf, &qf_data); CeedChkBackend(ierr);
749   CeedInt Q, P1d = 0, Q1d = 0, numelements, elemsize, numinputfields,
750           numoutputfields, ncomp, dim = 0, lsize;
751   ierr = CeedOperatorGetNumQuadraturePoints(op, &Q); CeedChkBackend(ierr);
752   ierr = CeedOperatorGetNumElements(op, &numelements); CeedChkBackend(ierr);
753   ierr = CeedQFunctionGetNumArgs(qf, &numinputfields, &numoutputfields);
754   CeedChkBackend(ierr);
755   CeedOperatorField *opinputfields, *opoutputfields;
756   ierr = CeedOperatorGetFields(op, &opinputfields, &opoutputfields);
757   CeedChkBackend(ierr);
758   CeedQFunctionField *qfinputfields, *qfoutputfields;
759   ierr = CeedQFunctionGetFields(qf, &qfinputfields, &qfoutputfields);
760   CeedChkBackend(ierr);
761   CeedEvalMode emode;
762   CeedBasis basis;
763   CeedBasis_Hip_shared *basis_data;
764   CeedElemRestriction Erestrict;
765   CeedElemRestriction_Hip *restr_data;
766 
767   // Check for restriction only identity operator
768   bool is_identity_qf;
769   ierr = CeedQFunctionIsIdentity(qf, &is_identity_qf); CeedChkBackend(ierr);
770   if (is_identity_qf) {
771     CeedEvalMode emodein, emodeout;
772     ierr = CeedQFunctionFieldGetEvalMode(qfinputfields[0], &emodein);  CeedChkBackend(ierr);
773     ierr = CeedQFunctionFieldGetEvalMode(qfoutputfields[0], &emodeout);  CeedChkBackend(ierr);
774     if (emodein == CEED_EVAL_NONE && emodeout == CEED_EVAL_NONE)
775       // LCOV_EXCL_START
776       return CeedError(ceed, CEED_ERROR_BACKEND,
777                        "Backend does not implement restriction only identity operators");
778     // LCOV_EXCL_STOP
779   }
780 
781   ostringstream code;
782   string devFunctions(deviceFunctions);
783   code << devFunctions;
784 
785   string qFunction(qf_data->qFunctionSource);
786   string qFunctionName(qf_data->qFunctionName);
787   string oper;
788   oper = "CeedKernel_Hip_gen_" + qFunctionName;
789 
790   code << "\n#define CEED_QFUNCTION(name) inline __device__ int name\n";
791   code << "#define CeedPragmaSIMD\n";
792   code << "#define CEED_ERROR_SUCCESS 0\n\n";
793 
794   // Find dim and Q1d
795   bool useCollograd = true;
796   data->maxP1d = 0;
797   for (CeedInt i = 0; i < numinputfields; i++) {
798     ierr = CeedOperatorFieldGetBasis(opinputfields[i], &basis); CeedChkBackend(ierr);
799     if (basis != CEED_BASIS_COLLOCATED) {
800       ierr = CeedBasisGetData(basis, &basis_data); CeedChkBackend(ierr);
801       ierr = CeedQFunctionFieldGetEvalMode(qfinputfields[i], &emode);
802       CeedChkBackend(ierr);
803 
804       // Check for collocated gradient
805       useCollograd = useCollograd && basis_data->d_collograd1d;
806 
807       // Collect dim and Q1d
808       ierr = CeedBasisGetDimension(basis, &dim); CeedChkBackend(ierr);
809       bool isTensor;
810       ierr = CeedBasisIsTensor(basis, &isTensor); CeedChkBackend(ierr);
811       if (isTensor) {
812         ierr = CeedBasisGetNumQuadraturePoints1D(basis, &Q1d); CeedChkBackend(ierr);
813         ierr = CeedBasisGetNumNodes1D(basis, &P1d); CeedChkBackend(ierr);
814         if (P1d>data->maxP1d) data->maxP1d = P1d;
815       } else {
816         // LCOV_EXCL_START
817         return CeedError(ceed, CEED_ERROR_BACKEND, "Backend does not implement operators with non-tensor basis");
818         // LCOV_EXCL_STOP
819         }
820     }
821   }
822   // Check output bases for Q1d, dim as well
823   //   The only imput basis might be CEED_BASIS_COLLOCATED
824   for (CeedInt i = 0; i < numoutputfields; i++) {
825     ierr = CeedOperatorFieldGetBasis(opoutputfields[i], &basis); CeedChkBackend(ierr);
826 
827     if (basis != CEED_BASIS_COLLOCATED) {
828       ierr = CeedBasisGetData(basis, &basis_data); CeedChkBackend(ierr);
829       ierr = CeedQFunctionFieldGetEvalMode(qfoutputfields[i], &emode);
830       CeedChkBackend(ierr);
831 
832       // Collect dim and Q1d
833       ierr = CeedBasisGetDimension(basis, &dim); CeedChkBackend(ierr);
834       bool isTensor;
835       ierr = CeedBasisIsTensor(basis, &isTensor); CeedChkBackend(ierr);
836       if (isTensor) {
837         ierr = CeedBasisGetNumQuadraturePoints1D(basis, &Q1d); CeedChkBackend(ierr);
838       } else {
839         // LCOV_EXCL_START
840         return CeedError(ceed, CEED_ERROR_BACKEND, "Backend does not implement operators with non-tensor basis");
841         // LCOV_EXCL_STOP
842         }
843 
844       // Check for collocated gradient
845       useCollograd = useCollograd && basis_data->d_collograd1d;
846     }
847   }
848   data->dim = dim;
849   data->Q1d = Q1d;
850 
851   // Define CEED_Q_VLA
852   if (dim != 3 || useCollograd) {
853     code << "\n#define CEED_Q_VLA 1\n\n";
854   } else {
855     code << "\n#define CEED_Q_VLA "<<Q1d<<"\n\n";
856   }
857 
858   code << qFunction;
859 
860   // Setup
861   code << "\n// -----------------------------------------------------------------------------\n";
862   code << "\nextern \"C\" __global__ void "<<oper<<"(CeedInt nelem, void* ctx, HipFieldsInt indices, HipFields fields, HipFields B, HipFields G, CeedScalar* W) {\n";
863   for (CeedInt i = 0; i < numinputfields; i++) {
864     ierr = CeedQFunctionFieldGetEvalMode(qfinputfields[i], &emode);
865     CeedChkBackend(ierr);
866     if (emode != CEED_EVAL_WEIGHT) { // Skip CEED_EVAL_WEIGHT
867       code << "  const CeedScalar* d_u" <<i<<" = fields.in["<<i<<"];\n";
868     }
869   }
870 
871   for (CeedInt i = 0; i < numoutputfields; i++) {
872     code << "  CeedScalar* d_v"<<i<<" = fields.out["<<i<<"];\n";
873   }
874 
875   code << "  const CeedInt Dim = "<<dim<<";\n";
876   code << "  const CeedInt Q1d = "<<Q1d<<";\n";
877 
878   code << "  HIP_DYNAMIC_SHARED( CeedScalar, slice)\n";
879   code << "  BackendData data;\n";
880   code << "  data.tidx = threadIdx.x;\n";
881   code << "  data.tidy = threadIdx.y;\n";
882   code << "  data.tidz = threadIdx.z;\n";
883   code << "  data.tid  = threadIdx.x + threadIdx.y*blockDim.x + threadIdx.z*blockDim.y*blockDim.x;\n";
884   code << "  data.slice = slice+data.tidz*T1d"<<(dim>1?"*T1d":"")<<";\n";
885 
886   code << "\n  // -- Input field constants and basis data --\n";
887   //Initialize constants, and matrices B and G
888   for (CeedInt i = 0; i < numinputfields; i++) {
889     code << "  // ---- Input field "<<i<<" ----\n";
890     // Get elemsize, emode, ncomp
891     ierr = CeedOperatorFieldGetElemRestriction(opinputfields[i], &Erestrict);
892     CeedChkBackend(ierr);
893     ierr = CeedElemRestrictionGetElementSize(Erestrict, &elemsize);
894     CeedChkBackend(ierr);
895     ierr = CeedQFunctionFieldGetEvalMode(qfinputfields[i], &emode);
896     CeedChkBackend(ierr);
897     ierr = CeedElemRestrictionGetNumComponents(Erestrict, &ncomp);
898     CeedChkBackend(ierr);
899 
900     // Set field constants
901     if (emode != CEED_EVAL_WEIGHT) {
902       ierr = CeedOperatorFieldGetBasis(opinputfields[i], &basis); CeedChkBackend(ierr);
903       if (basis != CEED_BASIS_COLLOCATED) {
904         ierr = CeedBasisGetNumNodes1D(basis, &P1d); CeedChkBackend(ierr);
905         code << "  const CeedInt P_in_"<<i<<" = "<<P1d<<";\n";
906       } else {
907         code << "  const CeedInt P_in_"<<i<<" = "<<Q1d<<";\n";
908       }
909       code << "  const CeedInt ncomp_in_"<<i<<" = "<<ncomp<<";\n";
910     }
911 
912     // Load basis data
913     code << "  // EvalMode: "<<CeedEvalModes[emode]<<"\n";
914     switch (emode) {
915     case CEED_EVAL_NONE:
916       break;
917     case CEED_EVAL_INTERP:
918       ierr = CeedBasisGetData(basis, &basis_data); CeedChkBackend(ierr);
919       data->B.in[i] = basis_data->d_interp1d;
920       code << "  __shared__ double s_B_in_"<<i<<"["<<P1d*Q1d<<"];\n";
921       code << "  loadMatrix<P_in_"<<i<<",Q1d>(data, B.in["<<i<<"], s_B_in_"<<i<<");\n";
922       break;
923     case CEED_EVAL_GRAD:
924       ierr = CeedBasisGetData(basis, &basis_data); CeedChkBackend(ierr);
925       data->B.in[i] = basis_data->d_interp1d;
926       code << "  __shared__ double s_B_in_"<<i<<"["<<P1d*Q1d<<"];\n";
927       code << "  loadMatrix<P_in_"<<i<<",Q1d>(data, B.in["<<i<<"], s_B_in_"<<i<<");\n";
928       if (useCollograd) {
929         data->G.in[i] = basis_data->d_collograd1d;
930         code << "  __shared__ double s_G_in_"<<i<<"["<<Q1d*Q1d<<"];\n";
931         code << "  loadMatrix<Q1d,Q1d>(data, G.in["<<i<<"], s_G_in_"<<i<<");\n";
932       } else {
933         data->G.in[i] = basis_data->d_grad1d;
934         code << "  __shared__ double s_G_in_"<<i<<"["<<P1d*Q1d<<"];\n";
935         code << "  loadMatrix<P_in_"<<i<<",Q1d>(data, G.in["<<i<<"], s_G_in_"<<i<<");\n";
936       }
937       break;
938     case CEED_EVAL_WEIGHT:
939       break; // No action
940     case CEED_EVAL_DIV:
941       break; // TODO: Not implemented
942     case CEED_EVAL_CURL:
943       break; // TODO: Not implemented
944     }
945   }
946 
947   code << "\n  // -- Output field constants and basis data --\n";
948   for (CeedInt i = 0; i < numoutputfields; i++) {
949     code << "  // ---- Output field "<<i<<" ----\n";
950     // Get elemsize, emode, ncomp
951     ierr = CeedOperatorFieldGetElemRestriction(opoutputfields[i], &Erestrict);
952     CeedChkBackend(ierr);
953     ierr = CeedElemRestrictionGetElementSize(Erestrict, &elemsize);
954     CeedChkBackend(ierr);
955     ierr = CeedQFunctionFieldGetEvalMode(qfoutputfields[i], &emode);
956     CeedChkBackend(ierr);
957     ierr = CeedElemRestrictionGetNumComponents(Erestrict, &ncomp);
958     CeedChkBackend(ierr);
959 
960     // Set field constants
961     ierr = CeedOperatorFieldGetBasis(opoutputfields[i], &basis); CeedChkBackend(ierr);
962     if (basis != CEED_BASIS_COLLOCATED) {
963       ierr = CeedBasisGetNumNodes1D(basis, &P1d); CeedChkBackend(ierr);
964       code << "  const CeedInt P_out_"<<i<<" = "<<P1d<<";\n";
965     } else {
966       code << "  const CeedInt P_out_"<<i<<" = "<<Q1d<<";\n";
967     }
968     code << "  const CeedInt ncomp_out_"<<i<<" = "<<ncomp<<";\n";
969 
970     // Load basis data
971     code << "  // EvalMode: "<<CeedEvalModes[emode]<<"\n";
972     switch (emode) {
973     case CEED_EVAL_NONE:
974       break; // No action
975     case CEED_EVAL_INTERP:
976       ierr = CeedBasisGetData(basis, &basis_data); CeedChkBackend(ierr);
977       data->B.out[i] = basis_data->d_interp1d;
978       code << "  __shared__ double s_B_out_"<<i<<"["<<P1d*Q1d<<"];\n";
979       code << "  loadMatrix<P_out_"<<i<<",Q1d>(data, B.out["<<i<<"], s_B_out_"<<i<<");\n";
980       break;
981     case CEED_EVAL_GRAD:
982       ierr = CeedBasisGetData(basis, &basis_data); CeedChkBackend(ierr);
983       data->B.out[i] = basis_data->d_interp1d;
984       code << "  __shared__ double s_B_out_"<<i<<"["<<P1d*Q1d<<"];\n";
985       code << "  loadMatrix<P_out_"<<i<<",Q1d>(data, B.out["<<i<<"], s_B_out_"<<i<<");\n";
986       if (useCollograd) {
987         data->G.out[i] = basis_data->d_collograd1d;
988         code << "  __shared__ double s_G_out_"<<i<<"["<<Q1d*Q1d<<"];\n";
989         code << "  loadMatrix<Q1d,Q1d>(data, G.out["<<i<<"], s_G_out_"<<i<<");\n";
990       } else {
991         data->G.out[i] = basis_data->d_grad1d;
992         code << "  __shared__ double s_G_out_"<<i<<"["<<P1d*Q1d<<"];\n";
993         code << "  loadMatrix<P_out_"<<i<<",Q1d>(data, G.out["<<i<<"], s_G_out_"<<i<<");\n";
994       }
995       break;
996     // LCOV_EXCL_START
997     case CEED_EVAL_WEIGHT: {
998       Ceed ceed;
999       ierr = CeedOperatorGetCeed(op, &ceed); CeedChkBackend(ierr);
1000       return CeedError(ceed, CEED_ERROR_BACKEND,
1001                        "CEED_EVAL_WEIGHT cannot be an output evaluation mode");
1002       break; // Should not occur
1003     }
1004     case CEED_EVAL_DIV:
1005       break; // TODO: Not implemented
1006     case CEED_EVAL_CURL:
1007       break; // TODO: Not implemented
1008       // LCOV_EXCL_STOP
1009     }
1010   }
1011   code << "\n  // -- Element loop --\n";
1012   code << "  __syncthreads();\n";
1013   code << "  for (CeedInt elem = blockIdx.x*blockDim.z + threadIdx.z; elem < nelem; elem += gridDim.x*blockDim.z) {\n";
1014   // Input basis apply if needed
1015   // Generate the correct eval mode code for each input
1016   code << "    // -- Input field restrictions and basis actions --\n";
1017   for (CeedInt i = 0; i < numinputfields; i++) {
1018     code << "    // ---- Input field "<<i<<" ----\n";
1019     // Get elemsize, emode, ncomp
1020     ierr = CeedOperatorFieldGetElemRestriction(opinputfields[i], &Erestrict);
1021     CeedChkBackend(ierr);
1022     ierr = CeedElemRestrictionGetElementSize(Erestrict, &elemsize);
1023     CeedChkBackend(ierr);
1024     ierr = CeedQFunctionFieldGetEvalMode(qfinputfields[i], &emode);
1025     CeedChkBackend(ierr);
1026     ierr = CeedElemRestrictionGetNumComponents(Erestrict, &ncomp);
1027     CeedChkBackend(ierr);
1028 
1029     // Restriction
1030     if (emode != CEED_EVAL_WEIGHT &&
1031         !((emode == CEED_EVAL_NONE) && useCollograd)) {
1032       code << "    CeedScalar r_u"<<i<<"[ncomp_in_"<<i<<"*P_in_"<<i<<"];\n";
1033 
1034       bool isStrided;
1035       ierr = CeedElemRestrictionIsStrided(Erestrict, &isStrided); CeedChkBackend(ierr);
1036       if (!isStrided) {
1037         ierr = CeedElemRestrictionGetLVectorSize(Erestrict, &lsize);
1038         CeedChkBackend(ierr);
1039         code << "    const CeedInt lsize_in_"<<i<<" = "<<lsize<<";\n";
1040         CeedInt compstride;
1041         ierr = CeedElemRestrictionGetCompStride(Erestrict, &compstride); CeedChkBackend(ierr);
1042         code << "    // CompStride: "<<compstride<<"\n";
1043         ierr = CeedElemRestrictionGetData(Erestrict, &restr_data); CeedChkBackend(ierr);
1044         data->indices.in[i] = restr_data->d_ind;
1045         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";
1046       } else {
1047         bool backendstrides;
1048         ierr = CeedElemRestrictionHasBackendStrides(Erestrict, &backendstrides);
1049         CeedChkBackend(ierr);
1050         CeedInt nelem;
1051         ierr = CeedElemRestrictionGetNumElements(Erestrict, &nelem);
1052         CeedChkBackend(ierr);
1053         CeedInt strides[3] = {1, elemsize*nelem, elemsize};
1054         if (!backendstrides) {
1055           ierr = CeedElemRestrictionGetStrides(Erestrict, &strides);
1056           CeedChkBackend(ierr);
1057         }
1058         code << "    // Strides: {"<<strides[0]<<", "<<strides[1]<<", "<<strides[2]<<"}\n";
1059         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";
1060       }
1061     }
1062 
1063     // Basis action
1064     code << "    // EvalMode: "<<CeedEvalModes[emode]<<"\n";
1065     switch (emode) {
1066     case CEED_EVAL_NONE:
1067       if (!useCollograd) {
1068         code << "    CeedScalar* r_t"<<i<<" = r_u"<<i<<";\n";
1069       }
1070       break;
1071     case CEED_EVAL_INTERP:
1072       code << "    CeedScalar r_t"<<i<<"[ncomp_in_"<<i<<"*Q1d];\n";
1073       code << "    interp"<<dim<<"d<ncomp_in_"<<i<<",P_in_"<<i<<",Q1d>(data, r_u"<<i<<", s_B_in_"<<i<<", r_t"<<i<<");\n";
1074       break;
1075     case CEED_EVAL_GRAD:
1076       if (useCollograd) {
1077         code << "    CeedScalar r_t"<<i<<"[ncomp_in_"<<i<<"*Q1d];\n";
1078         code << "    interp"<<dim<<"d<ncomp_in_"<<i<<",P_in_"<<i<<",Q1d>(data, r_u"<<i<<", s_B_in_"<<i<<", r_t"<<i<<");\n";
1079       } else {
1080         code << "    CeedScalar r_t"<<i<<"[ncomp_in_"<<i<<"*Dim*Q1d];\n";
1081         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";
1082       }
1083       break;
1084     case CEED_EVAL_WEIGHT:
1085       code << "    CeedScalar r_t"<<i<<"[Q1d];\n";
1086       ierr = CeedOperatorFieldGetBasis(opinputfields[i], &basis); CeedChkBackend(ierr);
1087       ierr = CeedBasisGetData(basis, &basis_data); CeedChkBackend(ierr);
1088       data->W = basis_data->d_qweight1d;
1089       code << "    weight"<<dim<<"d<Q1d>(data, W, r_t"<<i<<");\n";
1090       break; // No action
1091     case CEED_EVAL_DIV:
1092       break; // TODO: Not implemented
1093     case CEED_EVAL_CURL:
1094       break; // TODO: Not implemented
1095     }
1096   }
1097 
1098   // Q function
1099   code << "\n    // -- Output field setup --\n";
1100   for (CeedInt i = 0; i < numoutputfields; i++) {
1101       code << "\n    // ---- Output field "<<i<<" ----\n";
1102     ierr = CeedQFunctionFieldGetEvalMode(qfoutputfields[i], &emode);
1103     CeedChkBackend(ierr);
1104     if (emode==CEED_EVAL_GRAD)
1105     {
1106       if (useCollograd) {
1107         //Accumulator for gradient slices
1108         code << "    CeedScalar r_tt"<<i<<"[ncomp_out_"<<i<<"*Q1d];\n";
1109         code << "    for (CeedInt i = 0; i < ncomp_out_"<<i<<"; ++i) {\n";
1110         code << "      for (CeedInt j = 0; j < Q1d; ++j) {\n";
1111         code << "        r_tt"<<i<<"[j + i*Q1d] = 0.0;\n";
1112         code << "      }\n";
1113         code << "    }\n";
1114       } else {
1115         code << "    CeedScalar r_tt"<<i<<"[ncomp_out_"<<i<<"*Dim*Q1d];\n";
1116       }
1117     }
1118     if (emode==CEED_EVAL_NONE || emode==CEED_EVAL_INTERP)
1119     {
1120       code << "    CeedScalar r_tt"<<i<<"[ncomp_out_"<<i<<"*Q1d];\n";
1121     }
1122   }
1123   // We treat quadrature points per slice in 3d to save registers
1124   if (useCollograd) {
1125     code << "\n    // Note: Collocated Gradient\n";
1126     code << "#pragma unroll\n";
1127     code << "    for (CeedInt q=0; q<Q1d; q++) {\n";
1128     code << "      // -- Input fields --\n";
1129     for (CeedInt i = 0; i < numinputfields; i++) {
1130       code << "      // ---- Input field "<<i<<" ----\n";
1131       // Get elemsize, emode, ncomp
1132       ierr = CeedQFunctionFieldGetEvalMode(qfinputfields[i], &emode);
1133       CeedChkBackend(ierr);
1134       // Basis action
1135       code << "      // EvalMode: "<<CeedEvalModes[emode]<<"\n";
1136       switch (emode) {
1137       case CEED_EVAL_NONE:
1138         code << "      CeedScalar r_q"<<i<<"[ncomp_in_"<<i<<"];\n";
1139 
1140         bool isStrided;
1141         ierr = CeedOperatorFieldGetElemRestriction(opinputfields[i], &Erestrict); CeedChkBackend(ierr);
1142         ierr = CeedElemRestrictionGetElementSize(Erestrict, &elemsize); CeedChkBackend(ierr);
1143         ierr = CeedElemRestrictionIsStrided(Erestrict, &isStrided); CeedChkBackend(ierr);
1144         if (!isStrided) {
1145           ierr = CeedElemRestrictionGetLVectorSize(Erestrict, &lsize);
1146           CeedChkBackend(ierr);
1147           code << "      const CeedInt lsize_in_"<<i<<" = "<<lsize<<";\n";
1148           CeedInt compstride;
1149           ierr = CeedElemRestrictionGetCompStride(Erestrict, &compstride); CeedChkBackend(ierr);
1150           code << "      // CompStride: "<<compstride<<"\n";
1151           ierr = CeedElemRestrictionGetData(Erestrict, &restr_data); CeedChkBackend(ierr);
1152           data->indices.in[i] = restr_data->d_ind;
1153           code << "      readSliceQuadsOffset"<<"3d<ncomp_in_"<<i<<", "<<compstride<<", Q1d>(data, lsize_in_"<<i<<", elem, q, indices.in["<<i<<"], d_u"<<i<<", r_q"<<i<<");\n";
1154         } else {
1155           bool backendstrides;
1156           ierr = CeedElemRestrictionHasBackendStrides(Erestrict, &backendstrides);
1157           CeedChkBackend(ierr);
1158           CeedInt nelem;
1159           ierr = CeedElemRestrictionGetNumElements(Erestrict, &nelem);
1160           CeedChkBackend(ierr);
1161           CeedInt strides[3] = {1, elemsize*nelem, elemsize};
1162           if (!backendstrides) {
1163             ierr = CeedElemRestrictionGetStrides(Erestrict, &strides);
1164             CeedChkBackend(ierr);
1165           }
1166           code << "      // Strides: {"<<strides[0]<<", "<<strides[1]<<", "<<strides[2]<<"}\n";
1167           code << "      readSliceQuadsStrided"<<"3d<ncomp_in_"<<i<<",Q1d"","<<strides[0]<<","<<strides[1]<<","<<strides[2]<<">(data, elem, q, d_u"<<i<<", r_q"<<i<<");\n";
1168         }
1169         break;
1170       case CEED_EVAL_INTERP:
1171         code << "      CeedScalar r_q"<<i<<"[ncomp_in_"<<i<<"];\n";
1172         code << "      for (CeedInt j = 0; j < ncomp_in_"<<i<<" ; ++j) {\n";
1173         code << "        r_q"<<i<<"[j] = r_t"<<i<<"[q + j*Q1d];\n";
1174         code << "      }\n";
1175         break;
1176       case CEED_EVAL_GRAD:
1177         code << "      CeedScalar r_q"<<i<<"[ncomp_in_"<<i<<"*Dim];\n";
1178         code << "      gradCollo3d<ncomp_in_"<<i<<",Q1d>(data, q, r_t"<<i<<", s_G_in_"<<i<<", r_q"<<i<<");\n";
1179         break;
1180       case CEED_EVAL_WEIGHT:
1181         code << "      CeedScalar r_q"<<i<<"[1];\n";
1182         code << "      r_q"<<i<<"[0] = r_t"<<i<<"[q];\n";
1183         break; // No action
1184       case CEED_EVAL_DIV:
1185         break; // TODO: Not implemented
1186       case CEED_EVAL_CURL:
1187         break; // TODO: Not implemented
1188       }
1189     }
1190     code << "\n      // -- Output fields --\n";
1191     for (CeedInt i = 0; i < numoutputfields; i++) {
1192       code << "      // ---- Output field "<<i<<" ----\n";
1193       ierr = CeedQFunctionFieldGetEvalMode(qfoutputfields[i], &emode);
1194       CeedChkBackend(ierr);
1195       // Basis action
1196       switch (emode) {
1197       case CEED_EVAL_NONE:
1198         code << "      CeedScalar r_qq"<<i<<"[ncomp_out_"<<i<<"];\n";
1199         break; // No action
1200       case CEED_EVAL_INTERP:
1201         code << "      CeedScalar r_qq"<<i<<"[ncomp_out_"<<i<<"];\n";
1202         break;
1203       case CEED_EVAL_GRAD:
1204         code << "      CeedScalar r_qq"<<i<<"[ncomp_out_"<<i<<"*Dim];\n";
1205         break;
1206       case CEED_EVAL_WEIGHT:
1207         break; // Should not occur
1208       case CEED_EVAL_DIV:
1209         break; // TODO: Not implemented
1210       case CEED_EVAL_CURL:
1211         break; // TODO: Not implemented
1212       }
1213     }
1214   } else {
1215     code << "\n      // Note: No Collocated Gradient\n";
1216     code << "      // -- Input fields --\n";
1217     for (CeedInt i = 0; i < numinputfields; i++) {
1218       code << "      // ---- Input field "<<i<<" ----\n";
1219       code << "      CeedScalar* r_q"<<i<<" = r_t"<<i<<";\n";
1220     }
1221     code << "      // -- Output fields --\n";
1222     for (CeedInt i = 0; i < numoutputfields; i++) {
1223       code << "      // ---- Output field "<<i<<" ----\n";
1224       code << "      CeedScalar* r_qq"<<i<<" = r_tt"<<i<<";\n";
1225     }
1226   }
1227   code << "\n      // -- QFunction Inputs and outputs --\n";
1228   code << "      CeedScalar* in["<<numinputfields<<"];\n";
1229   for (CeedInt i = 0; i < numinputfields; i++) {
1230     code << "      // ---- Input field "<<i<<" ----\n";
1231     code << "      in["<<i<<"] = r_q"<<i<<";\n";
1232   }
1233   code << "      CeedScalar* out["<<numoutputfields<<"];\n";
1234   for (CeedInt i = 0; i < numoutputfields; i++) {
1235     code << "      // ---- Output field "<<i<<" ----\n";
1236     code << "      out["<<i<<"] = r_qq"<<i<<";\n";
1237   }
1238   code << "\n      // -- Apply QFunction --\n";
1239   code << "      "<<qFunctionName<<"(ctx, ";
1240   if (dim != 3 || useCollograd) {
1241     code << "1";
1242   } else {
1243     code << "Q1d";
1244   }
1245   code << ", in, out);\n";
1246   if (useCollograd) {
1247     code << "\n      // Note: Collocated Gradient\n";
1248     code << "      // -- Output fields --\n";
1249     for (CeedInt i = 0; i < numoutputfields; i++) {
1250       code << "      // ---- Output field "<<i<<" ----\n";
1251       ierr = CeedQFunctionFieldGetEvalMode(qfoutputfields[i], &emode);
1252       CeedChkBackend(ierr);
1253       // Basis action
1254       code << "      // EvalMode: "<<CeedEvalModes[emode]<<"\n";
1255       switch (emode) {
1256       case CEED_EVAL_NONE:
1257         code << "      for (CeedInt j = 0; j < ncomp_out_"<<i<<" ; ++j) {\n";
1258         code << "        r_tt"<<i<<"[q + j*Q1d] = r_qq"<<i<<"[j];\n";
1259         code << "      }\n";
1260         break; // No action
1261       case CEED_EVAL_INTERP:
1262         code << "      for (CeedInt j = 0; j < ncomp_out_"<<i<<" ; ++j) {\n";
1263         code << "        r_tt"<<i<<"[q + j*Q1d] = r_qq"<<i<<"[j];\n";
1264         code << "      }\n";
1265         break;
1266       case CEED_EVAL_GRAD:
1267         code << "      gradColloTranspose3d<ncomp_out_"<<i<<",Q1d>(data, q, r_qq"<<i<<", s_G_out_"<<i<<", r_tt"<<i<<");\n";
1268         break;
1269       case CEED_EVAL_WEIGHT:
1270         break; // Should not occur
1271       case CEED_EVAL_DIV:
1272         break; // TODO: Not implemented
1273       case CEED_EVAL_CURL:
1274         break; // TODO: Not implemented
1275       }
1276     }
1277     code << "    }\n";
1278   }
1279 
1280   // Output basis apply if needed
1281   // Generate the correct eval mode code for each output
1282   code << "\n    // -- Output field basis action and restrictions --\n";
1283   for (CeedInt i = 0; i < numoutputfields; i++) {
1284     code << "    // ---- Output field "<<i<<" ----\n";
1285     // Get elemsize, emode, ncomp
1286     ierr = CeedOperatorFieldGetElemRestriction(opoutputfields[i], &Erestrict);
1287     CeedChkBackend(ierr);
1288     ierr = CeedElemRestrictionGetElementSize(Erestrict, &elemsize);
1289     CeedChkBackend(ierr);
1290     ierr = CeedQFunctionFieldGetEvalMode(qfoutputfields[i], &emode);
1291     CeedChkBackend(ierr);
1292     ierr = CeedElemRestrictionGetNumComponents(Erestrict, &ncomp);
1293     CeedChkBackend(ierr);
1294     // Basis action
1295     code << "    // EvalMode: "<<CeedEvalModes[emode]<<"\n";
1296     switch (emode) {
1297     case CEED_EVAL_NONE:
1298       code << "    CeedScalar* r_v"<<i<<" = r_tt"<<i<<";\n";
1299       break; // No action
1300     case CEED_EVAL_INTERP:
1301       code << "    CeedScalar r_v"<<i<<"[ncomp_out_"<<i<<"*P_out_"<<i<<"];\n";
1302       code << "    interpTranspose"<<dim<<"d<ncomp_out_"<<i<<",P_out_"<<i<<",Q1d>(data, r_tt"<<i<<", s_B_out_"<<i<<", r_v"<<i<<");\n";
1303       break;
1304     case CEED_EVAL_GRAD:
1305       code << "    CeedScalar r_v"<<i<<"[ncomp_out_"<<i<<"*P_out_"<<i<<"];\n";
1306       if (useCollograd) {
1307         code << "    interpTranspose"<<dim<<"d<ncomp_out_"<<i<<",P_out_"<<i<<",Q1d>(data, r_tt"<<i<<", s_B_out_"<<i<<", r_v"<<i<<");\n";
1308       } else {
1309         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";
1310       }
1311       break;
1312     // LCOV_EXCL_START
1313     case CEED_EVAL_WEIGHT: {
1314       Ceed ceed;
1315       ierr = CeedOperatorGetCeed(op, &ceed); CeedChkBackend(ierr);
1316       return CeedError(ceed, CEED_ERROR_BACKEND,
1317                        "CEED_EVAL_WEIGHT cannot be an output evaluation mode");
1318       break; // Should not occur
1319     }
1320     case CEED_EVAL_DIV:
1321       break; // TODO: Not implemented
1322     case CEED_EVAL_CURL:
1323       break; // TODO: Not implemented
1324       // LCOV_EXCL_STOP
1325     }
1326     // Restriction
1327       bool isStrided;
1328       ierr = CeedElemRestrictionIsStrided(Erestrict, &isStrided); CeedChkBackend(ierr);
1329     if (!isStrided) {
1330       ierr = CeedElemRestrictionGetLVectorSize(Erestrict, &lsize);
1331       CeedChkBackend(ierr);
1332       code << "    const CeedInt lsize_out_"<<i<<" = "<<lsize<<";\n";
1333       CeedInt compstride;
1334       ierr = CeedElemRestrictionGetCompStride(Erestrict, &compstride); CeedChkBackend(ierr);
1335       code << "    // CompStride: "<<compstride<<"\n";
1336       ierr = CeedElemRestrictionGetData(Erestrict, &restr_data); CeedChkBackend(ierr);
1337       data->indices.out[i] = restr_data->d_ind;
1338       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";
1339     } else {
1340       bool backendstrides;
1341       ierr = CeedElemRestrictionHasBackendStrides(Erestrict, &backendstrides);
1342       CeedChkBackend(ierr);
1343       CeedInt nelem;
1344       ierr = CeedElemRestrictionGetNumElements(Erestrict, &nelem);
1345       CeedChkBackend(ierr);
1346       CeedInt strides[3] = {1, elemsize*nelem, elemsize};
1347       if (!backendstrides) {
1348         ierr = CeedElemRestrictionGetStrides(Erestrict, &strides);
1349         CeedChkBackend(ierr);
1350       }
1351       code << "    // Strides: {"<<strides[0]<<", "<<strides[1]<<", "<<strides[2]<<"}\n";
1352       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";
1353     }
1354   }
1355 
1356   code << "  }\n";
1357   code << "}\n";
1358   code << "// -----------------------------------------------------------------------------\n\n";
1359 
1360   // View kernel for debugging
1361   CeedDebug(code.str().c_str());
1362 
1363   ierr = CeedCompileHip(ceed, code.str().c_str(), &data->module, 1,
1364                          "T1d", CeedIntMax(Q1d, data->maxP1d));
1365   CeedChkBackend(ierr);
1366   ierr = CeedGetKernelHip(ceed, data->module, oper.c_str(), &data->op);
1367   CeedChkBackend(ierr);
1368 
1369   ierr = CeedOperatorSetSetupDone(op); CeedChkBackend(ierr);
1370   return CEED_ERROR_SUCCESS;
1371 }
1372 //------------------------------------------------------------------------------
1373