xref: /libCEED/backends/cuda-gen/ceed-cuda-gen-operator-build.cpp (revision 6b92dc4b38fb4a2ad1e306f024e5961058ad12b3)
1 // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors.
2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3 //
4 // SPDX-License-Identifier: BSD-2-Clause
5 //
6 // This file is part of CEED:  http://github.com/ceed
7 
8 #define CEED_DEBUG_COLOR 12
9 
10 #include <ceed.h>
11 #include <ceed/backend.h>
12 #include <ceed/jit-tools.h>
13 #include <cuda_runtime.h>
14 
15 #include <iostream>
16 #include <sstream>
17 #include <string>
18 
19 #include "../cuda-ref/ceed-cuda-ref.h"
20 #include "../cuda-shared/ceed-cuda-shared.h"
21 #include "../cuda/ceed-cuda-common.h"
22 #include "../cuda/ceed-cuda-compile.h"
23 #include "ceed-cuda-gen.h"
24 
25 struct FieldReuse_Cuda {
26   CeedInt      index;
27   bool         is_input;
28   CeedEvalMode eval_mode;
29 };
30 
31 //------------------------------------------------------------------------------
32 // Determine type of operator
33 //------------------------------------------------------------------------------
34 static int CeedOperatorBuildKernelData_Cuda_gen(Ceed ceed, CeedInt num_input_fields, CeedOperatorField *op_input_fields,
35                                                 CeedQFunctionField *qf_input_fields, CeedInt num_output_fields, CeedOperatorField *op_output_fields,
36                                                 CeedQFunctionField *qf_output_fields, CeedInt *max_P_1d, CeedInt *Q_1d, CeedInt *dim, bool *is_tensor,
37                                                 bool *use_3d_slices) {
38   // Find dim, P_1d, Q_1d
39   *max_P_1d  = 0;
40   *Q_1d      = 0;
41   *dim       = 0;
42   *is_tensor = true;
43 
44   for (CeedInt i = 0; i < num_input_fields; i++) {
45     CeedBasis basis;
46 
47     CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
48     if (basis != CEED_BASIS_NONE) {
49       bool    is_field_tensor;
50       CeedInt field_P_1d = 0, field_Q_1d = 0, field_dim = 0;
51 
52       // Collect dim, P_1d, and Q_1d
53       CeedCallBackend(CeedBasisIsTensor(basis, &is_field_tensor));
54       *is_tensor = *is_tensor && is_field_tensor;
55       if (is_field_tensor) CeedCallBackend(CeedBasisGetNumNodes1D(basis, &field_P_1d));
56       else CeedCallBackend(CeedBasisGetNumNodes(basis, &field_P_1d));
57       *max_P_1d = CeedIntMax(*max_P_1d, field_P_1d);
58       CeedCallBackend(CeedBasisGetDimension(basis, &field_dim));
59       CeedCheck(*dim == 0 || field_dim == *dim, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible");
60       *dim = field_dim;
61       if (is_field_tensor) CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &field_Q_1d));
62       else CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis, &field_Q_1d));
63       CeedCheck(*Q_1d == 0 || field_Q_1d == *Q_1d, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible");
64       *Q_1d = field_Q_1d;
65     }
66     CeedCallBackend(CeedBasisDestroy(&basis));
67   }
68   for (CeedInt i = 0; i < num_output_fields; i++) {
69     CeedBasis basis;
70 
71     CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
72     if (basis != CEED_BASIS_NONE) {
73       bool    is_field_tensor;
74       CeedInt field_P_1d = 0, field_Q_1d = 0, field_dim = 0;
75 
76       // Collect dim, P_1d, and Q_1d
77       CeedCallBackend(CeedBasisIsTensor(basis, &is_field_tensor));
78       *is_tensor = *is_tensor && is_field_tensor;
79       if (is_field_tensor) CeedCallBackend(CeedBasisGetNumNodes1D(basis, &field_P_1d));
80       else CeedCallBackend(CeedBasisGetNumNodes(basis, &field_P_1d));
81       *max_P_1d = CeedIntMax(*max_P_1d, field_P_1d);
82       CeedCallBackend(CeedBasisGetDimension(basis, &field_dim));
83       CeedCheck(*dim == 0 || field_dim == *dim, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible");
84       *dim = field_dim;
85       if (is_field_tensor) CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &field_Q_1d));
86       else CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis, &field_Q_1d));
87       CeedCheck(*Q_1d == 0 || field_Q_1d == *Q_1d, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible");
88       *Q_1d = field_Q_1d;
89     }
90     CeedCallBackend(CeedBasisDestroy(&basis));
91   }
92 
93   // Only use 3D collocated gradient parallelization strategy when gradient is computed
94   *use_3d_slices = false;
95   if (*dim == 3) {
96     bool was_grad_found = false;
97 
98     for (CeedInt i = 0; i < num_input_fields; i++) {
99       CeedEvalMode eval_mode;
100 
101       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
102       if (eval_mode == CEED_EVAL_GRAD) {
103         CeedBasis_Cuda_shared *basis_data;
104         CeedBasis              basis;
105 
106         CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
107         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
108         *use_3d_slices = basis_data->d_collo_grad_1d && (was_grad_found ? *use_3d_slices : true);
109         was_grad_found = true;
110         CeedCallBackend(CeedBasisDestroy(&basis));
111       }
112     }
113     for (CeedInt i = 0; i < num_output_fields; i++) {
114       CeedEvalMode eval_mode;
115 
116       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
117       if (eval_mode == CEED_EVAL_GRAD) {
118         CeedBasis_Cuda_shared *basis_data;
119         CeedBasis              basis;
120 
121         CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
122         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
123         *use_3d_slices = basis_data->d_collo_grad_1d && (was_grad_found ? *use_3d_slices : true);
124         was_grad_found = true;
125         CeedCallBackend(CeedBasisDestroy(&basis));
126       }
127     }
128   }
129   return CEED_ERROR_SUCCESS;
130 }
131 
132 //------------------------------------------------------------------------------
133 // Setup fields
134 //------------------------------------------------------------------------------
135 static int CeedOperatorBuildKernelFieldData_Cuda_gen(std::ostringstream &code, CeedOperator_Cuda_gen *data, CeedInt i, CeedOperatorField op_field,
136                                                      CeedQFunctionField qf_field, FieldReuse_Cuda field_reuse, CeedInt Q_1d, bool is_input,
137                                                      bool is_tensor, bool is_at_points, bool use_3d_slices) {
138   const char            *field_name;
139   std::string            var_suffix = (is_input ? "_in_" : "_out_") + std::to_string(i);
140   std::string            P_name = (is_tensor ? "P_1d" : "P") + var_suffix, Q_name = is_tensor ? "Q_1d" : "Q";
141   std::string            option_name = (is_input ? "inputs" : "outputs");
142   CeedEvalMode           eval_mode   = CEED_EVAL_NONE;
143   CeedInt                elem_size = 0, num_comp = 0, P_1d = 0;
144   CeedElemRestriction    elem_rstr;
145   CeedBasis_Cuda_shared *basis_data;
146   CeedBasis              basis;
147 
148   // Field reuse info
149   bool use_previous_field = field_reuse.index != -1;
150 
151   CeedCallBackend(CeedOperatorFieldGetName(op_field, &field_name));
152   code << "  // -- " << (is_input ? "Input" : "Output") << " field " << i << ": " << field_name << "\n";
153 
154   // Get field data
155   CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_field, &elem_rstr));
156   if (elem_rstr != CEED_ELEMRESTRICTION_NONE) {
157     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
158     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
159   }
160   CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
161   CeedCallBackend(CeedOperatorFieldGetBasis(op_field, &basis));
162   if (basis != CEED_BASIS_NONE) {
163     CeedCallBackend(CeedBasisGetData(basis, &basis_data));
164     if (is_tensor) CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
165     else CeedCallBackend(CeedBasisGetNumNodes(basis, &P_1d));
166   }
167   CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_field, &eval_mode));
168 
169   // Set field constants
170   if (eval_mode != CEED_EVAL_WEIGHT) {
171     code << "  const CeedInt " << P_name << " = " << (basis == CEED_BASIS_NONE ? Q_1d : P_1d) << ";\n";
172     code << "  const CeedInt num_comp" << var_suffix << " = " << num_comp << ";\n";
173   }
174 
175   // Load basis data
176   code << "  // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
177   switch (eval_mode) {
178     case CEED_EVAL_NONE:
179       break;
180     case CEED_EVAL_INTERP:
181       if (is_at_points) {
182         // AtPoints
183         if (!basis_data->d_chebyshev_interp_1d) {
184           CeedSize    interp_bytes;
185           CeedScalar *chebyshev_interp_1d;
186 
187           interp_bytes = P_1d * Q_1d * sizeof(CeedScalar);
188           CeedCallBackend(CeedCalloc(P_1d * Q_1d, &chebyshev_interp_1d));
189           CeedCallBackend(CeedBasisGetChebyshevInterp1D(basis, chebyshev_interp_1d));
190           CeedCallCuda(CeedBasisReturnCeed(basis), cudaMalloc((void **)&basis_data->d_chebyshev_interp_1d, interp_bytes));
191           CeedCallCuda(CeedBasisReturnCeed(basis),
192                        cudaMemcpy(basis_data->d_chebyshev_interp_1d, chebyshev_interp_1d, interp_bytes, cudaMemcpyHostToDevice));
193           CeedCallBackend(CeedFree(&chebyshev_interp_1d));
194         }
195         if (is_input) data->B.inputs[i] = basis_data->d_chebyshev_interp_1d;
196         else data->B.outputs[i] = basis_data->d_chebyshev_interp_1d;
197       } else {
198         // Standard quadrature
199         if (is_input) data->B.inputs[i] = basis_data->d_interp_1d;
200         else data->B.outputs[i] = basis_data->d_interp_1d;
201       }
202       if (use_previous_field) {
203         std::string reuse_var = "s_B" + ((field_reuse.is_input ? "_in_" : "_out_") + std::to_string(field_reuse.index));
204 
205         code << "  CeedScalar *s_B" << var_suffix << " = " << reuse_var << ";\n";
206       } else {
207         code << "  __shared__ CeedScalar s_B" << var_suffix << "[" << P_name << "*" << Q_name << "];\n";
208         code << "  LoadMatrix<" << P_name << ", " << Q_name << ">(data, B." << option_name << "[" << i << "], s_B" << var_suffix << ");\n";
209       }
210       break;
211     case CEED_EVAL_GRAD:
212       if (is_at_points) {
213         // AtPoints
214         if (!basis_data->d_chebyshev_interp_1d) {
215           CeedSize    interp_bytes;
216           CeedScalar *chebyshev_interp_1d;
217 
218           interp_bytes = P_1d * Q_1d * sizeof(CeedScalar);
219           CeedCallBackend(CeedCalloc(P_1d * Q_1d, &chebyshev_interp_1d));
220           CeedCallBackend(CeedBasisGetChebyshevInterp1D(basis, chebyshev_interp_1d));
221           CeedCallCuda(CeedBasisReturnCeed(basis), cudaMalloc((void **)&basis_data->d_chebyshev_interp_1d, interp_bytes));
222           CeedCallCuda(CeedBasisReturnCeed(basis),
223                        cudaMemcpy(basis_data->d_chebyshev_interp_1d, chebyshev_interp_1d, interp_bytes, cudaMemcpyHostToDevice));
224           CeedCallBackend(CeedFree(&chebyshev_interp_1d));
225         }
226         if (is_input) data->B.inputs[i] = basis_data->d_chebyshev_interp_1d;
227         else data->B.outputs[i] = basis_data->d_chebyshev_interp_1d;
228       } else {
229         // Standard quadrature
230         if (is_input) data->B.inputs[i] = basis_data->d_interp_1d;
231         else data->B.outputs[i] = basis_data->d_interp_1d;
232       }
233       if (is_tensor) {
234         if (use_previous_field) {
235           std::string reuse_var = "s_B" + ((field_reuse.is_input ? "_in_" : "_out_") + std::to_string(field_reuse.index));
236 
237           code << "  CeedScalar *s_B" << var_suffix << " = " << reuse_var << ";\n";
238         } else {
239           code << "  __shared__ CeedScalar s_B" << var_suffix << "[" << P_name << "*" << Q_name << "];\n";
240           code << "  LoadMatrix<" << P_name << ", " << Q_name << ">(data, B." << option_name << "[" << i << "], s_B" << var_suffix << ");\n";
241         }
242       }
243       if (is_at_points) break;  // No G mat for AtPoints
244       if (use_3d_slices) {
245         if (is_input) data->G.inputs[i] = basis_data->d_collo_grad_1d;
246         else data->G.outputs[i] = basis_data->d_collo_grad_1d;
247         if (use_previous_field && field_reuse.eval_mode == CEED_EVAL_GRAD) {
248           std::string reuse_var = "s_G" + ((field_reuse.is_input ? "_in_" : "_out_") + std::to_string(field_reuse.index));
249 
250           code << "  CeedScalar *s_G" << var_suffix << " = " << reuse_var << ";\n";
251         } else {
252           code << "  __shared__ CeedScalar s_G" << var_suffix << "[" << Q_name << "*" << Q_name << "];\n";
253           code << "  LoadMatrix<" << Q_name << ", " << Q_name << ">(data, G." << option_name << "[" << i << "], s_G" << var_suffix << ");\n";
254         }
255       } else {
256         bool has_collo_grad = basis_data->d_collo_grad_1d;
257 
258         if (is_input) data->G.inputs[i] = has_collo_grad ? basis_data->d_collo_grad_1d : basis_data->d_grad_1d;
259         else data->G.outputs[i] = has_collo_grad ? basis_data->d_collo_grad_1d : basis_data->d_grad_1d;
260         if (has_collo_grad) {
261           if (use_previous_field && field_reuse.eval_mode == CEED_EVAL_GRAD) {
262             std::string reuse_var = "s_G" + ((field_reuse.is_input ? "_in_" : "_out_") + std::to_string(field_reuse.index));
263 
264             code << "  CeedScalar *s_G" << var_suffix << " = " << reuse_var << ";\n";
265           } else {
266             code << "  __shared__ CeedScalar s_G" << var_suffix << "[" << Q_name << "*" << Q_name << "];\n";
267             code << "  LoadMatrix<" << Q_name << ", " << Q_name << ">(data, G." << option_name << "[" << i << "], s_G" << var_suffix << ");\n";
268           }
269         } else {
270           if (use_previous_field && field_reuse.eval_mode == CEED_EVAL_GRAD) {
271             std::string reuse_var = "s_G" + ((field_reuse.is_input ? "_in_" : "_out_") + std::to_string(field_reuse.index));
272 
273             code << "  CeedScalar *s_G" << var_suffix << " = " << reuse_var << ";\n";
274           } else {
275             code << "  __shared__ CeedScalar s_G" << var_suffix << "[" << P_name << "*" << Q_name << (is_tensor ? "" : "*dim") << "];\n";
276             code << "  LoadMatrix<" << P_name << ", " << Q_name << (is_tensor ? "" : "*dim") << ">(data, G." << option_name << "[" << i << "], s_G"
277                  << var_suffix << ");\n";
278           }
279         }
280       }
281       break;
282     case CEED_EVAL_WEIGHT:
283       break;  // No action
284       // LCOV_EXCL_START
285     case CEED_EVAL_DIV:
286     case CEED_EVAL_CURL:
287       break;  // TODO: Not implemented
288               // LCOV_EXCL_STOP
289   }
290   CeedCallBackend(CeedBasisDestroy(&basis));
291   return CEED_ERROR_SUCCESS;
292 }
293 
294 //------------------------------------------------------------------------------
295 // Restriction
296 //------------------------------------------------------------------------------
297 static int CeedOperatorBuildKernelRestriction_Cuda_gen(std::ostringstream &code, CeedOperator_Cuda_gen *data, CeedInt i, CeedInt dim,
298                                                        CeedInt field_input_buffer[], CeedOperatorField op_field, CeedQFunctionField qf_field,
299                                                        CeedInt Q_1d, bool is_input, bool is_tensor, bool is_at_points, bool use_3d_slices) {
300   std::string               var_suffix = (is_input ? "_in_" : "_out_") + std::to_string(i);
301   std::string               P_name     = (is_tensor ? "P_1d" : "P") + var_suffix;
302   CeedEvalMode              eval_mode  = CEED_EVAL_NONE;
303   CeedInt                   elem_size = 0, num_comp = 0, P_1d = 0;
304   CeedSize                  l_size;
305   CeedRestrictionType       rstr_type = CEED_RESTRICTION_STANDARD;
306   CeedElemRestriction_Cuda *rstr_data;
307   CeedElemRestriction       elem_rstr;
308   CeedBasis                 basis;
309 
310   // Get field data
311   CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_field, &elem_rstr));
312   if (elem_rstr != CEED_ELEMRESTRICTION_NONE) {
313     CeedCallBackend(CeedElemRestrictionGetType(elem_rstr, &rstr_type));
314     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
315     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
316     CeedCallBackend(CeedElemRestrictionGetData(elem_rstr, &rstr_data));
317   }
318   CeedCallBackend(CeedOperatorFieldGetBasis(op_field, &basis));
319   if (basis != CEED_BASIS_NONE) {
320     if (is_tensor) CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
321     else CeedCallBackend(CeedBasisGetNumNodes(basis, &P_1d));
322   }
323   CeedCallBackend(CeedBasisDestroy(&basis));
324   CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_field, &eval_mode));
325 
326   // Restriction
327   if (is_input) {
328     // Input
329     if (field_input_buffer[i] != i) {
330       std::string buffer_name = "r_e_in_" + std::to_string(field_input_buffer[i]);
331 
332       // Restriction was already done for previous input
333       code << "    CeedScalar *r_e" << var_suffix << " = " << buffer_name << ";\n";
334     } else if (eval_mode != CEED_EVAL_WEIGHT && !((eval_mode == CEED_EVAL_NONE) && use_3d_slices && is_at_points)) {
335       if (eval_mode == CEED_EVAL_NONE && rstr_type != CEED_RESTRICTION_POINTS) {
336         // No basis action, so r_e_in_* in also r_q_in_* and needs to be allocated
337         code << "    CeedScalar r_e" << var_suffix << "[num_comp" << var_suffix << "*" << P_name << "];\n";
338       } else if (rstr_type != CEED_RESTRICTION_POINTS) {
339         // Otherwise we're using the scratch space
340         code << "    CeedScalar *r_e" << var_suffix << " = r_e_scratch;\n";
341       }
342       switch (rstr_type) {
343         case CEED_RESTRICTION_STANDARD: {
344           CeedInt comp_stride;
345 
346           CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size));
347           code << "    const CeedInt l_size" << var_suffix << " = " << l_size << ";\n";
348           CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
349           code << "    // CompStride: " << comp_stride << "\n";
350           data->indices.inputs[i] = (CeedInt *)rstr_data->d_offsets;
351           code << "    ReadLVecStandard" << (is_tensor ? dim : 1) << "d<num_comp" << var_suffix << ", " << comp_stride << ", " << P_name
352                << ">(data, l_size" << var_suffix << ", elem, indices.inputs[" << i << "], d" << var_suffix << ", r_e" << var_suffix << ");\n";
353           break;
354         }
355         case CEED_RESTRICTION_STRIDED: {
356           bool    has_backend_strides;
357           CeedInt num_elem;
358 
359           CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides));
360           CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem));
361           CeedInt strides[3] = {1, elem_size * num_elem, elem_size};
362 
363           if (!has_backend_strides) {
364             CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides));
365           }
366           code << "    // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n";
367           code << "    ReadLVecStrided" << (is_tensor ? dim : 1) << "d<num_comp" << var_suffix << ", " << P_name << ", " << strides[0] << ", "
368                << strides[1] << ", " << strides[2] << ">(data, elem, d" << var_suffix << ", r_e" << var_suffix << ");\n";
369           break;
370         }
371         case CEED_RESTRICTION_POINTS: {
372           CeedInt comp_stride;
373 
374           CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
375           code << "    const CeedInt comp_stride" << var_suffix << " = " << comp_stride << ";\n";
376           data->indices.inputs[i] = (CeedInt *)rstr_data->d_offsets;
377           break;
378         }
379         // LCOV_EXCL_START
380         case CEED_RESTRICTION_ORIENTED:
381         case CEED_RESTRICTION_CURL_ORIENTED:
382           break;  // TODO: Not implemented
383                   // LCOV_EXCL_STOP
384       }
385     }
386   } else {
387     // Output
388     switch (rstr_type) {
389       case CEED_RESTRICTION_STANDARD: {
390         CeedInt comp_stride;
391 
392         CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size));
393         code << "    const CeedInt l_size" << var_suffix << " = " << l_size << ";\n";
394         CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
395         code << "    // CompStride: " << comp_stride << "\n";
396         data->indices.outputs[i] = (CeedInt *)rstr_data->d_offsets;
397         code << "    WriteLVecStandard" << (is_tensor ? dim : 1) << "d<num_comp" << var_suffix << ", " << comp_stride << ", " << P_name
398              << ">(data, l_size" << var_suffix << ", elem, indices.outputs[" << i << "], r_e" << var_suffix << ", d" << var_suffix << ");\n";
399         break;
400       }
401       case CEED_RESTRICTION_STRIDED: {
402         bool    has_backend_strides;
403         CeedInt num_elem;
404 
405         CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides));
406         CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem));
407         CeedInt strides[3] = {1, elem_size * num_elem, elem_size};
408 
409         if (!has_backend_strides) {
410           CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides));
411         }
412         code << "    // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n";
413         code << "    WriteLVecStrided" << (is_tensor ? dim : 1) << "d<num_comp" << var_suffix << ", " << P_name << ", " << strides[0] << ", "
414              << strides[1] << ", " << strides[2] << ">(data, elem, r_e" << var_suffix << ", d" << var_suffix << ");\n";
415         break;
416       }
417       case CEED_RESTRICTION_POINTS:
418         data->indices.outputs[i] = (CeedInt *)rstr_data->d_offsets;
419         break;
420       // LCOV_EXCL_START
421       case CEED_RESTRICTION_ORIENTED:
422       case CEED_RESTRICTION_CURL_ORIENTED:
423         break;  // TODO: Not implemented
424                 // LCOV_EXCL_STOP
425     }
426   }
427   CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
428   return CEED_ERROR_SUCCESS;
429 }
430 
431 //------------------------------------------------------------------------------
432 // Basis
433 //------------------------------------------------------------------------------
434 static int CeedOperatorBuildKernelBasis_Cuda_gen(std::ostringstream &code, CeedOperator_Cuda_gen *data, CeedInt i, CeedInt dim,
435                                                  CeedOperatorField op_field, CeedQFunctionField qf_field, CeedInt Q_1d, bool is_input, bool is_tensor,
436                                                  bool is_at_points, bool use_3d_slices) {
437   std::string         var_suffix = (is_input ? "_in_" : "_out_") + std::to_string(i);
438   std::string         P_name = (is_tensor ? "P_1d" : "P") + var_suffix, Q_name = is_tensor ? "Q_1d" : "Q";
439   CeedEvalMode        eval_mode = CEED_EVAL_NONE;
440   CeedInt             elem_size = 0, num_comp = 0, P_1d = 0;
441   CeedElemRestriction elem_rstr;
442   CeedBasis           basis;
443 
444   // Get field data
445   CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_field, &elem_rstr));
446   if (elem_rstr != CEED_ELEMRESTRICTION_NONE) {
447     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
448     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
449   }
450   CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
451   CeedCallBackend(CeedOperatorFieldGetBasis(op_field, &basis));
452   if (basis != CEED_BASIS_NONE) {
453     if (is_tensor) CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
454     else CeedCallBackend(CeedBasisGetNumNodes(basis, &P_1d));
455   }
456   CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_field, &eval_mode));
457 
458   // Basis
459   code << "    // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
460   if (is_input) {
461     switch (eval_mode) {
462       case CEED_EVAL_NONE:
463         if (!use_3d_slices && !is_at_points) {
464           code << "    CeedScalar *r_q" << var_suffix << " = r_e" << var_suffix << ";\n";
465         }
466         break;
467       case CEED_EVAL_INTERP:
468         if (is_at_points) {
469           std::string function_name = (dim == 1 ? "Interp" : "InterpTensor") + std::to_string(dim) + "d";
470 
471           code << "    CeedScalar r_c" << var_suffix << "[num_comp" << var_suffix << "*" << (dim >= 3 ? Q_name : "1") << "];\n";
472           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ", OP_T_1D>(data, r_e" << var_suffix
473                << ", s_B" << var_suffix << ", r_c" << var_suffix << ");\n";
474         } else {
475           std::string function_name = is_tensor ? ((dim == 1 ? "Interp" : "InterpTensor") + std::to_string(dim) + "d") : "InterpNonTensor";
476 
477           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << (is_tensor && (dim >= 3) ? Q_name : "1") << "];\n";
478           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ", OP_T_1D>(data, r_e" << var_suffix
479                << ", s_B" << var_suffix << ", r_q" << var_suffix << ");\n";
480         }
481         break;
482       case CEED_EVAL_GRAD:
483         if (is_at_points) {
484           std::string function_name = (dim == 1 ? "Interp" : "InterpTensor") + std::to_string(dim) + "d";
485 
486           code << "    CeedScalar r_c" << var_suffix << "[num_comp" << var_suffix << "*" << (dim >= 3 ? Q_name : "1") << "];\n";
487           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ", OP_T_1D>(data, r_e" << var_suffix
488                << ", s_B" << var_suffix << ", r_c" << var_suffix << ");\n";
489         } else if (use_3d_slices) {
490           std::string function_name = (dim > 1 ? "InterpTensor" : "Interp") + std::to_string(dim) + "d";
491 
492           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << Q_name << "];\n";
493           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ", OP_T_1D>(data, r_e" << var_suffix
494                << ", s_B" << var_suffix << ", r_q" << var_suffix << ");\n";
495         } else if (is_tensor) {
496           bool        is_collocated = dim == 3 && Q_1d >= P_1d;
497           std::string function_name = (dim == 1 ? "Grad" : (is_collocated ? "GradTensorCollocated" : "GradTensor")) + std::to_string(dim) + "d";
498 
499           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*dim*" << (dim >= 3 ? Q_name : "1") << "];\n";
500           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ", OP_T_1D>(data, r_e" << var_suffix
501                << ", s_B" << var_suffix << ", s_G" << var_suffix << ", r_q" << var_suffix << ");\n";
502         } else {
503           std::string function_name = "GradNonTensor";
504 
505           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*dim];\n";
506           code << "    " << function_name << "<num_comp" << var_suffix << ", dim, " << P_name << ", " << Q_name << ", OP_T_1D>(data, r_e"
507                << var_suffix << ", s_G" << var_suffix << ", r_q" << var_suffix << ");\n";
508         }
509         break;
510       case CEED_EVAL_WEIGHT: {
511         if (is_at_points) {
512           code << "    // Nothing to do AtPoints\n";
513         } else {
514           CeedBasis_Cuda_shared *basis_data;
515           std::string            function_name = is_tensor ? ((dim == 1 ? "Weight" : "WeightTensor") + std::to_string(dim) + "d") : "WeightNonTensor";
516 
517           code << "    CeedScalar r_q" << var_suffix << "[" << (is_tensor && (dim >= 3) ? Q_name : "1") << "];\n";
518           CeedCallBackend(CeedBasisGetData(basis, &basis_data));
519           data->W = basis_data->d_q_weight_1d;
520           code << "    " << function_name << "<" << Q_name << ">(data, W, r_q" << var_suffix << ");\n";
521         }
522         break;
523       }
524       // LCOV_EXCL_START
525       case CEED_EVAL_DIV:
526       case CEED_EVAL_CURL:
527         break;  // TODO: Not implemented
528                 // LCOV_EXCL_STOP
529     }
530   } else {
531     switch (eval_mode) {
532       case CEED_EVAL_NONE:
533         code << "    CeedScalar *r_e" << var_suffix << " = r_q" << var_suffix << ";\n";
534         break;  // No action
535       case CEED_EVAL_INTERP:
536         code << "    CeedScalar *r_e" << var_suffix << " = r_e_scratch;\n";
537         if (is_at_points) {
538           std::string function_name = (dim == 1 ? "InterpTranspose" : "InterpTransposeTensor") + std::to_string(dim) + "d";
539 
540           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ", OP_T_1D>(data, r_c" << var_suffix
541                << ", s_B" << var_suffix << ", r_e" << var_suffix << ");\n";
542         } else {
543           std::string function_name =
544               is_tensor ? ((dim == 1 ? "InterpTranspose" : "InterpTransposeTensor") + std::to_string(dim) + "d") : "InterpTransposeNonTensor";
545 
546           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ", OP_T_1D>(data, r_q" << var_suffix
547                << ", s_B" << var_suffix << ", r_e" << var_suffix << ");\n";
548         }
549         break;
550       case CEED_EVAL_GRAD:
551         code << "    CeedScalar *r_e" << var_suffix << " = r_e_scratch;\n";
552         if (is_at_points) {
553           std::string function_name = (dim == 1 ? "InterpTranspose" : "InterpTransposeTensor") + std::to_string(dim) + "d";
554 
555           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ", OP_T_1D>(data, r_c" << var_suffix
556                << ", s_B" << var_suffix << ", r_e" << var_suffix << ");\n";
557         } else if (use_3d_slices) {
558           std::string function_name = (dim == 1 ? "InterpTranspose" : "InterpTransposeTensor") + std::to_string(dim) + "d";
559 
560           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ", OP_T_1D>(data, r_q" << var_suffix
561                << ", s_B" << var_suffix << ", r_e" << var_suffix << ");\n";
562         } else if (is_tensor) {
563           bool        is_collocated = dim == 3 && Q_1d >= P_1d;
564           std::string function_name =
565               (dim == 1 ? "GradTranspose" : (is_collocated ? "GradTransposeTensorCollocated" : "GradTransposeTensor")) + std::to_string(dim) + "d";
566 
567           code << "    " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ", OP_T_1D>(data, r_q" << var_suffix
568                << ", s_B" << var_suffix << ", s_G" << var_suffix << ", r_e" << var_suffix << ");\n";
569         } else {
570           std::string function_name = "GradTransposeNonTensor";
571 
572           code << "    " << function_name << "<num_comp" << var_suffix << ", dim, " << P_name << ", " << Q_name << ", OP_T_1D>(data, r_q"
573                << var_suffix << ", s_G" << var_suffix << ", r_e" << var_suffix << ");\n";
574         }
575         break;
576       // LCOV_EXCL_START
577       case CEED_EVAL_WEIGHT:
578         break;  // Should not occur
579       case CEED_EVAL_DIV:
580       case CEED_EVAL_CURL:
581         break;  // TODO: Not implemented
582                 // LCOV_EXCL_STOP
583     }
584   }
585   CeedCallBackend(CeedBasisDestroy(&basis));
586   return CEED_ERROR_SUCCESS;
587 }
588 
589 //------------------------------------------------------------------------------
590 // QFunction
591 //------------------------------------------------------------------------------
592 static int CeedOperatorBuildKernelQFunction_Cuda_gen(std::ostringstream &code, CeedOperator_Cuda_gen *data, CeedInt dim, CeedInt max_num_points,
593                                                      CeedInt num_input_fields, CeedOperatorField *op_input_fields,
594                                                      CeedQFunctionField *qf_input_fields, CeedInt num_output_fields,
595                                                      CeedOperatorField *op_output_fields, CeedQFunctionField *qf_output_fields,
596                                                      std::string qfunction_name, CeedInt Q_1d, bool is_tensor, bool is_at_points,
597                                                      bool use_3d_slices) {
598   std::string         Q_name    = is_tensor ? "Q_1d" : "Q";
599   CeedEvalMode        eval_mode = CEED_EVAL_NONE;
600   CeedElemRestriction elem_rstr;
601 
602   // Setup output arrays
603   code << "\n    // -- Output field setup\n";
604   for (CeedInt i = 0; i < num_output_fields; i++) {
605     const char *field_name;
606     std::string var_suffix = "_out_" + std::to_string(i);
607 
608     CeedCallBackend(CeedOperatorFieldGetName(op_output_fields[i], &field_name));
609     code << "    // ---- Output field " << i << ": " << field_name << "\n";
610     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
611     switch (eval_mode) {
612       case CEED_EVAL_NONE:
613         if (is_at_points) {
614           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "];\n";
615         } else {
616           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << (is_tensor && (dim >= 3) ? Q_name : "1") << "];\n";
617         }
618         break;
619       case CEED_EVAL_INTERP:
620         if (is_at_points) {
621           // Accumulator for point data
622           code << "    CeedScalar r_c" << var_suffix << "[num_comp" << var_suffix << "*" << (dim >= 3 ? Q_name : "1") << "];\n";
623           code << "    for (CeedInt i = 0; i < num_comp" << var_suffix << "*" << (dim >= 3 ? Q_name : "1") << "; i++) {\n";
624           code << "      r_c" << var_suffix << "[i] = 0.0;\n";
625           code << "    }\n";
626         } else {
627           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << (is_tensor && (dim >= 3) ? Q_name : "1") << "];\n";
628         }
629         break;
630       case CEED_EVAL_GRAD:
631         if (is_at_points) {
632           // Accumulator for point data
633           code << "    CeedScalar r_c" << var_suffix << "[num_comp" << var_suffix << "*" << (dim >= 3 ? Q_name : "1") << "*dim];\n";
634           code << "    for (CeedInt i = 0; i < num_comp" << var_suffix << "*" << (dim >= 3 ? Q_name : "1") << "; i++) {\n";
635           code << "      r_c" << var_suffix << "[i] = 0.0;\n";
636           code << "    }\n";
637         } else if (use_3d_slices) {
638           // Accumulator for gradient slices
639           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << Q_name << "];\n";
640           code << "    for (CeedInt i = 0; i < num_comp" << var_suffix << "*" << Q_name << "; i++) {\n";
641           code << "      r_q" << var_suffix << "[i] = 0.0;\n";
642           code << "    }\n";
643         } else {
644           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*dim*" << (is_tensor && (dim >= 3) ? Q_name : "1") << "];\n";
645         }
646         break;
647       case CEED_EVAL_WEIGHT:
648         break;
649         // LCOV_EXCL_START
650       case CEED_EVAL_DIV:
651       case CEED_EVAL_CURL:
652         break;  // TODO: Not implemented
653                 // LCOV_EXCL_STOP
654     }
655   }
656 
657   if (is_at_points) {
658     // We need to handle batches of points
659     code << "\n    // Note: Using batches of points\n";
660     code << "    const CeedInt point_loop_bound = (blockDim.x * blockDim.y) * ceil(1.0 * max_num_points / (blockDim.x * blockDim.y));\n\n";
661     code << "    #pragma unroll\n";
662     code << "    for (CeedInt i = threadIdx.x + threadIdx.y * blockDim.x; i < point_loop_bound; i += blockDim.x * blockDim.y) {\n";
663     code << "      const CeedInt p = i % max_num_points;\n\n";
664 
665     code << "      // -- Coordinates\n";
666     code << "      CeedScalar r_x[dim];\n";
667     code << "      ReadPoint<dim, coords_comp_stride, max_num_points>(data, elem, p, max_num_points, points.indices, points.coords, r_x);\n\n";
668 
669     code << "      // -- Input fields\n";
670     for (CeedInt i = 0; i < num_input_fields; i++) {
671       const char *field_name;
672       std::string var_suffix = "_in_" + std::to_string(i);
673 
674       CeedCallBackend(CeedOperatorFieldGetName(op_input_fields[i], &field_name));
675       code << "      // ---- Input field " << i << ": " << field_name << "\n";
676       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
677       // Basis action
678       code << "      // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
679       switch (eval_mode) {
680         case CEED_EVAL_NONE:
681           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
682           code << "      ReadPoint<num_comp" << var_suffix << ", comp_stride" << var_suffix
683                << ", max_num_points>(data, elem, p, max_num_points, indices.inputs[" << i << "], d" << var_suffix << ", r_s" << var_suffix << ");\n";
684           break;
685         case CEED_EVAL_INTERP:
686           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
687           code << "      InterpAtPoints" << dim << "d<num_comp" << var_suffix << ", max_num_points, " << Q_name << ">(data, i, r_c" << var_suffix
688                << ", r_x, r_s" << var_suffix << ");\n";
689           break;
690         case CEED_EVAL_GRAD:
691           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim];\n";
692           code << "      GradAtPoints" << dim << "d<num_comp" << var_suffix << ", max_num_points, " << Q_name << ">(data, i, r_c" << var_suffix
693                << ", r_x, r_s" << var_suffix << ");\n";
694           break;
695         case CEED_EVAL_WEIGHT:
696           code << "      CeedScalar r_s" << var_suffix << "[1];\n";
697           code << "      r_s" << var_suffix << "[0] = 1.0;\n";
698           break;
699           // LCOV_EXCL_START
700         case CEED_EVAL_DIV:
701         case CEED_EVAL_CURL:
702           break;  // TODO: Not implemented
703                   // LCOV_EXCL_STOP
704       }
705     }
706     code << "\n      // -- Output fields\n";
707     for (CeedInt i = 0; i < num_output_fields; i++) {
708       const char *field_name;
709       std::string var_suffix = "_out_" + std::to_string(i);
710 
711       CeedCallBackend(CeedOperatorFieldGetName(op_output_fields[i], &field_name));
712       code << "      // ---- Output field " << i << ": " << field_name << "\n";
713       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
714       // Basis action
715       switch (eval_mode) {
716         case CEED_EVAL_NONE:
717           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
718           break;
719         case CEED_EVAL_INTERP:
720           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
721           break;
722         case CEED_EVAL_GRAD:
723           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim];\n";
724           break;
725           // LCOV_EXCL_START
726         case CEED_EVAL_WEIGHT:
727           break;  // Should not occur
728         case CEED_EVAL_DIV:
729         case CEED_EVAL_CURL:
730           break;  // TODO: Not implemented
731                   // LCOV_EXCL_STOP
732       }
733     }
734 
735   } else if (use_3d_slices) {
736     // We treat quadrature points per slice in 3d to save registers
737     code << "\n    // Note: Using planes of 3D elements\n";
738     code << "    #pragma unroll\n";
739     code << "    for (CeedInt q = 0; q < " << Q_name << "; q++) {\n";
740     code << "      // -- Input fields\n";
741     for (CeedInt i = 0; i < num_input_fields; i++) {
742       const char *field_name;
743       std::string var_suffix = "_in_" + std::to_string(i);
744 
745       CeedCallBackend(CeedOperatorFieldGetName(op_input_fields[i], &field_name));
746       code << "      // ---- Input field " << i << ": " << field_name << "\n";
747       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
748       // Basis action
749       code << "      // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
750       switch (eval_mode) {
751         case CEED_EVAL_NONE:
752           bool is_strided;
753 
754           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
755 
756           CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr));
757           CeedCallBackend(CeedElemRestrictionIsStrided(elem_rstr, &is_strided));
758           if (is_strided) {
759             bool    has_backend_strides;
760             CeedInt num_elem, elem_size;
761 
762             CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
763             CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides));
764             CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem));
765             CeedInt strides[3] = {1, elem_size * num_elem, elem_size};
766 
767             if (!has_backend_strides) {
768               CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides));
769             }
770             code << "      // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n";
771             code << "      ReadEVecSliceStrided3d<num_comp" << var_suffix << ", " << Q_name << ", " << strides[0] << ", " << strides[1] << ", "
772                  << strides[2] << ">(data, elem, q, d" << var_suffix << ", r_s" << var_suffix << ");\n";
773           } else {
774             CeedSize                  l_size = 0;
775             CeedInt                   comp_stride;
776             CeedElemRestriction_Cuda *rstr_data;
777 
778             CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size));
779             code << "      const CeedInt l_size" << var_suffix << " = " << l_size << ";\n";
780             CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
781             code << "      // CompStride: " << comp_stride << "\n";
782             CeedCallBackend(CeedElemRestrictionGetData(elem_rstr, &rstr_data));
783             data->indices.inputs[i] = (CeedInt *)rstr_data->d_offsets;
784             code << "      ReadEVecSliceStandard3d<num_comp" << var_suffix << ", " << comp_stride << ", " << Q_name << ">(data, l_size" << var_suffix
785                  << ", elem, q, indices.inputs[" << i << "], d" << var_suffix << ", r_s" << var_suffix << ");\n";
786           }
787           CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
788           break;
789         case CEED_EVAL_INTERP:
790           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
791           code << "      for (CeedInt j = 0; j < num_comp" << var_suffix << "; j++) {\n";
792           code << "        r_s" << var_suffix << "[j] = r_q" << var_suffix << "[q + j*" << Q_name << "];\n";
793           code << "      }\n";
794           break;
795         case CEED_EVAL_GRAD:
796           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim];\n";
797           code << "      GradColloSlice3d<num_comp" << var_suffix << ", " << Q_name << ", OP_T_1D>(data, q, r_q" << var_suffix << ", s_G"
798                << var_suffix << ", r_s" << var_suffix << ");\n";
799           break;
800         case CEED_EVAL_WEIGHT:
801           code << "      CeedScalar r_s" << var_suffix << "[1];\n";
802           code << "      r_s" << var_suffix << "[0] = r_q" << var_suffix << "[q];\n";
803           break;
804           // LCOV_EXCL_START
805         case CEED_EVAL_DIV:
806         case CEED_EVAL_CURL:
807           break;  // TODO: Not implemented
808                   // LCOV_EXCL_STOP
809       }
810     }
811     code << "\n      // -- Output fields\n";
812     for (CeedInt i = 0; i < num_output_fields; i++) {
813       const char *field_name;
814       std::string var_suffix = "_out_" + std::to_string(i);
815 
816       CeedCallBackend(CeedOperatorFieldGetName(op_output_fields[i], &field_name));
817       code << "      // ---- Output field " << i << ": " << field_name << "\n";
818       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
819       // Basis action
820       switch (eval_mode) {
821         case CEED_EVAL_NONE:
822           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
823           break;
824         case CEED_EVAL_INTERP:
825           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
826           break;
827         case CEED_EVAL_GRAD:
828           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim];\n";
829           break;
830           // LCOV_EXCL_START
831         case CEED_EVAL_WEIGHT:
832           break;  // Should not occur
833         case CEED_EVAL_DIV:
834         case CEED_EVAL_CURL:
835           break;  // TODO: Not implemented
836                   // LCOV_EXCL_STOP
837       }
838     }
839   } else {
840     code << "\n    // Note: Using full elements\n";
841     code << "    {\n";
842     code << "      // -- Input fields\n";
843     for (CeedInt i = 0; i < num_input_fields; i++) {
844       const char *field_name;
845 
846       CeedCallBackend(CeedOperatorFieldGetName(op_input_fields[i], &field_name));
847       code << "      // ---- Input field " << i << ": " << field_name << "\n";
848       code << "      CeedScalar *r_s_in_" << i << " = r_q_in_" << i << ";\n";
849     }
850     code << "      // -- Output fields\n";
851     for (CeedInt i = 0; i < num_output_fields; i++) {
852       const char *field_name;
853 
854       CeedCallBackend(CeedOperatorFieldGetName(op_output_fields[i], &field_name));
855       code << "      // ---- Output field " << i << ": " << field_name << "\n";
856       code << "      CeedScalar *r_s_out_" << i << " = r_q_out_" << i << ";\n";
857     }
858   }
859 
860   // Input and output buffers
861   code << "\n      // -- QFunction inputs and outputs\n";
862   code << "      // ---- Inputs\n";
863   code << "      CeedScalar *inputs[" << CeedIntMax(num_input_fields, 1) << "];\n";
864   for (CeedInt i = 0; i < num_input_fields; i++) {
865     const char *field_name;
866 
867     CeedCallBackend(CeedOperatorFieldGetName(op_input_fields[i], &field_name));
868     code << "      // ------ Input field " << i << ": " << field_name << "\n";
869     code << "      inputs[" << i << "] = r_s_in_" << i << ";\n";
870   }
871   code << "      // ---- Outputs\n";
872   code << "      CeedScalar *outputs[" << CeedIntMax(num_output_fields, 1) << "];\n";
873   for (CeedInt i = 0; i < num_output_fields; i++) {
874     const char *field_name;
875 
876     CeedCallBackend(CeedOperatorFieldGetName(op_output_fields[i], &field_name));
877     code << "      // ------ Output field " << i << ": " << field_name << "\n";
878     code << "      outputs[" << i << "] = r_s_out_" << i << ";\n";
879   }
880 
881   // Apply QFunction
882   code << "\n      // -- Apply QFunction\n";
883   code << "      " << qfunction_name << "(ctx, ";
884   if (dim != 3 || is_at_points || use_3d_slices || !is_tensor) {
885     code << "1";
886   } else {
887     code << Q_name;
888   }
889   code << ", inputs, outputs);\n";
890 
891   if (is_at_points) {
892     // Map back to coefficients
893     code << "\n      // -- Output fields\n";
894     for (CeedInt i = 0; i < num_output_fields; i++) {
895       const char *field_name;
896       std::string var_suffix = "_out_" + std::to_string(i);
897       std::string P_name     = "P_1d" + var_suffix;
898 
899       CeedCallBackend(CeedOperatorFieldGetName(op_output_fields[i], &field_name));
900       code << "      // ---- Output field " << i << ": " << field_name << "\n";
901       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
902       // Basis action
903       code << "      // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
904       switch (eval_mode) {
905         case CEED_EVAL_NONE: {
906           CeedInt             comp_stride;
907           CeedElemRestriction elem_rstr;
908 
909           CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &elem_rstr));
910           CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
911           CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
912           code << "      const CeedInt comp_stride" << var_suffix << " = " << comp_stride << ";\n";
913           code << "      WritePoint<num_comp" << var_suffix << ", comp_stride" << var_suffix
914                << ", max_num_points>(data, elem, i, points.num_per_elem[elem], indices.outputs[" << i << "]"
915                << ", r_s" << var_suffix << ", d" << var_suffix << ");\n";
916           break;
917         }
918         case CEED_EVAL_INTERP:
919           code << "      if (i >= points.num_per_elem[elem]) {\n";
920           code << "        for (CeedInt j = 0; j < num_comp" << var_suffix << "; j++) r_s" << var_suffix << "[j] = 0.0;\n";
921           code << "      }\n";
922           code << "      InterpTransposeAtPoints" << dim << "d<num_comp" << var_suffix << ", max_num_points, " << Q_name << ">(data, i, r_s"
923                << var_suffix << ", r_x, r_c" << var_suffix << ");\n";
924           break;
925         case CEED_EVAL_GRAD:
926           code << "      if (i >= points.num_per_elem[elem]) {\n";
927           code << "        for (CeedInt j = 0; j < num_comp" << var_suffix << "*dim; j++) r_s" << var_suffix << "[j] = 0.0;\n";
928           code << "      }\n";
929           code << "      GradTransposeAtPoints" << dim << "d<num_comp" << var_suffix << ", max_num_points, " << Q_name << ">(data, i, r_s"
930                << var_suffix << ", r_x, r_c" << var_suffix << ");\n";
931           break;
932           // LCOV_EXCL_START
933         case CEED_EVAL_WEIGHT:
934           break;  // Should not occur
935         case CEED_EVAL_DIV:
936         case CEED_EVAL_CURL:
937           break;  // TODO: Not implemented
938                   // LCOV_EXCL_STOP
939       }
940     }
941   } else if (use_3d_slices) {
942     // Copy or apply transpose grad, if needed
943     code << "\n      // -- Output fields\n";
944     for (CeedInt i = 0; i < num_output_fields; i++) {
945       const char *field_name;
946       std::string var_suffix = "_out_" + std::to_string(i);
947       std::string P_name     = "P_1d" + var_suffix;
948 
949       CeedCallBackend(CeedOperatorFieldGetName(op_output_fields[i], &field_name));
950       code << "      // ---- Output field " << i << ": " << field_name << "\n";
951       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
952       // Basis action
953       code << "      // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
954       switch (eval_mode) {
955         case CEED_EVAL_NONE:
956           code << "      for (CeedInt j = 0; j < num_comp" << var_suffix << " ; j++) {\n";
957           code << "        r_q" << var_suffix << "[q + j*" << Q_name << "] = r_s" << var_suffix << "[j];\n";
958           code << "      }\n";
959           break;
960         case CEED_EVAL_INTERP:
961           code << "      for (CeedInt j = 0; j < num_comp" << var_suffix << " ; j++) {\n";
962           code << "        r_q" << var_suffix << "[q + j*" << Q_name << "] = r_s" << var_suffix << "[j];\n";
963           code << "      }\n";
964           break;
965         case CEED_EVAL_GRAD:
966           code << "      GradColloSliceTranspose3d<num_comp" << var_suffix << ", " << Q_name << ", OP_T_1D>(data, q, r_s" << var_suffix << ", s_G"
967                << var_suffix << ", r_q" << var_suffix << ");\n";
968           break;
969           // LCOV_EXCL_START
970         case CEED_EVAL_WEIGHT:
971           break;  // Should not occur
972         case CEED_EVAL_DIV:
973         case CEED_EVAL_CURL:
974           break;  // TODO: Not implemented
975                   // LCOV_EXCL_STOP
976       }
977     }
978   }
979   code << "    }\n";
980   return CEED_ERROR_SUCCESS;
981 }
982 
983 //------------------------------------------------------------------------------
984 // Build single operator kernel
985 //------------------------------------------------------------------------------
986 extern "C" int CeedOperatorBuildKernel_Cuda_gen(CeedOperator op, bool *is_good_build) {
987   bool                    is_tensor = true, is_at_points = false, use_3d_slices = false;
988   Ceed                    ceed;
989   CeedInt                 Q_1d, num_input_fields, num_output_fields, dim = 1, max_num_points = 0, coords_comp_stride = 0;
990   CeedQFunctionField     *qf_input_fields, *qf_output_fields;
991   CeedQFunction_Cuda_gen *qf_data;
992   CeedQFunction           qf;
993   CeedOperatorField      *op_input_fields, *op_output_fields;
994   CeedOperator_Cuda_gen  *data;
995   std::ostringstream      code;
996 
997   CeedCallBackend(CeedOperatorGetData(op, &data));
998   {
999     bool is_setup_done;
1000 
1001     CeedCallBackend(CeedOperatorIsSetupDone(op, &is_setup_done));
1002     if (is_setup_done) {
1003       *is_good_build = !data->use_fallback;
1004       return CEED_ERROR_SUCCESS;
1005     }
1006   }
1007 
1008   // Check field compatibility
1009   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
1010   {
1011     bool has_shared_bases = true, is_all_tensor = true, is_all_nontensor = true;
1012 
1013     for (CeedInt i = 0; i < num_input_fields; i++) {
1014       CeedBasis basis;
1015 
1016       CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
1017       if (basis != CEED_BASIS_NONE) {
1018         bool        is_tensor = true;
1019         const char *resource;
1020         char       *resource_root;
1021         Ceed        basis_ceed;
1022 
1023         CeedCallBackend(CeedBasisIsTensor(basis, &is_tensor));
1024         is_all_tensor    = is_all_tensor && is_tensor;
1025         is_all_nontensor = is_all_nontensor && !is_tensor;
1026         CeedCallBackend(CeedBasisGetCeed(basis, &basis_ceed));
1027         CeedCallBackend(CeedGetResource(basis_ceed, &resource));
1028         CeedCallBackend(CeedGetResourceRoot(basis_ceed, resource, ":", &resource_root));
1029         has_shared_bases = has_shared_bases && !strcmp(resource_root, "/gpu/cuda/shared");
1030         CeedCallBackend(CeedFree(&resource_root));
1031         CeedCallBackend(CeedDestroy(&basis_ceed));
1032       }
1033       CeedCallBackend(CeedBasisDestroy(&basis));
1034     }
1035 
1036     for (CeedInt i = 0; i < num_output_fields; i++) {
1037       CeedBasis basis;
1038 
1039       CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
1040       if (basis != CEED_BASIS_NONE) {
1041         bool        is_tensor = true;
1042         const char *resource;
1043         char       *resource_root;
1044         Ceed        basis_ceed;
1045 
1046         CeedCallBackend(CeedBasisIsTensor(basis, &is_tensor));
1047         is_all_tensor    = is_all_tensor && is_tensor;
1048         is_all_nontensor = is_all_nontensor && !is_tensor;
1049 
1050         CeedCallBackend(CeedBasisGetCeed(basis, &basis_ceed));
1051         CeedCallBackend(CeedGetResource(basis_ceed, &resource));
1052         CeedCallBackend(CeedGetResourceRoot(basis_ceed, resource, ":", &resource_root));
1053         has_shared_bases = has_shared_bases && !strcmp(resource_root, "/gpu/cuda/shared");
1054         CeedCallBackend(CeedFree(&resource_root));
1055         CeedCallBackend(CeedDestroy(&basis_ceed));
1056       }
1057       CeedCallBackend(CeedBasisDestroy(&basis));
1058     }
1059     // -- Fallback to ref if not all bases are shared
1060     if (!has_shared_bases || (!is_all_tensor && !is_all_nontensor)) {
1061       *is_good_build = false;
1062       return CEED_ERROR_SUCCESS;
1063     }
1064   }
1065   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
1066   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
1067   CeedCallBackend(CeedQFunctionGetData(qf, &qf_data));
1068   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
1069 
1070   // Get operator data
1071   CeedCallBackend(CeedOperatorIsAtPoints(op, &is_at_points));
1072   CeedCallBackend(CeedOperatorBuildKernelData_Cuda_gen(ceed, num_input_fields, op_input_fields, qf_input_fields, num_output_fields, op_output_fields,
1073                                                        qf_output_fields, &data->max_P_1d, &Q_1d, &dim, &is_tensor, &use_3d_slices));
1074   if (dim == 0) dim = 1;
1075   data->dim = dim;
1076   if (is_at_points) {
1077     CeedElemRestriction_Cuda *rstr_data;
1078     CeedElemRestriction       rstr_points = NULL;
1079 
1080     CeedCallBackend(CeedOperatorAtPointsGetPoints(op, &rstr_points, NULL));
1081     CeedCallBackend(CeedElemRestrictionGetMaxPointsInElement(rstr_points, &max_num_points));
1082     CeedCallBackend(CeedElemRestrictionGetCompStride(rstr_points, &coords_comp_stride));
1083     CeedCallBackend(CeedElemRestrictionGetData(rstr_points, &rstr_data));
1084     data->points.indices = (CeedInt *)rstr_data->d_offsets;
1085     CeedCallBackend(CeedElemRestrictionDestroy(&rstr_points));
1086   }
1087   if (is_at_points) use_3d_slices = false;
1088   if (Q_1d == 0) {
1089     if (is_at_points) Q_1d = max_num_points;
1090     else CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q_1d));
1091   }
1092   data->Q_1d = Q_1d;
1093 
1094   // Check for restriction only identity operator
1095   {
1096     bool is_identity_qf;
1097 
1098     CeedCallBackend(CeedQFunctionIsIdentity(qf, &is_identity_qf));
1099     if (is_identity_qf) {
1100       CeedEvalMode eval_mode_in, eval_mode_out;
1101 
1102       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[0], &eval_mode_in));
1103       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[0], &eval_mode_out));
1104       CeedCheck(eval_mode_in != CEED_EVAL_NONE || eval_mode_out != CEED_EVAL_NONE, ceed, CEED_ERROR_BACKEND,
1105                 "Backend does not implement restriction only identity operators");
1106     }
1107   }
1108 
1109   // Add atomicAdd function for old NVidia architectures
1110   {
1111     Ceed_Cuda            *ceed_data;
1112     struct cudaDeviceProp prop;
1113 
1114     CeedCallBackend(CeedGetData(ceed, &ceed_data));
1115     CeedCallBackend(cudaGetDeviceProperties(&prop, ceed_data->device_id));
1116     if ((prop.major < 6) && (CEED_SCALAR_TYPE != CEED_SCALAR_FP32)) {
1117       code << "// AtomicAdd fallback source\n";
1118       code << "#include <ceed/jit-source/cuda/cuda-atomic-add-fallback.h>\n\n";
1119     }
1120   }
1121 
1122   // Load basis source files
1123   if (is_tensor) {
1124     code << "// Tensor basis source\n";
1125     code << "#include <ceed/jit-source/cuda/cuda-shared-basis-tensor-templates.h>\n\n";
1126   } else {
1127     code << "// Non-tensor basis source\n";
1128     code << "#include <ceed/jit-source/cuda/cuda-shared-basis-nontensor-templates.h>\n\n";
1129   }
1130   if (is_at_points) {
1131     code << "// AtPoints basis source\n";
1132     code << "#include <ceed/jit-source/cuda/cuda-shared-basis-tensor-at-points-templates.h>\n\n";
1133   }
1134   code << "// CodeGen operator source\n";
1135   code << "#include <ceed/jit-source/cuda/cuda-gen-templates.h>\n\n";
1136 
1137   // Get QFunction name
1138   std::string qfunction_name(qf_data->qfunction_name);
1139   std::string operator_name;
1140 
1141   operator_name = "CeedKernelCudaGenOperator_" + qfunction_name;
1142 
1143   // Define CEED_Q_VLA
1144   code << "\n#undef CEED_Q_VLA\n";
1145   if (dim != 3 || is_at_points || use_3d_slices || !is_tensor) {
1146     code << "#define CEED_Q_VLA 1\n\n";
1147   } else {
1148     code << "#define CEED_Q_VLA " << Q_1d << "\n\n";
1149   }
1150 
1151   // Add user QFunction source
1152   {
1153     const char *source_path;
1154 
1155     CeedCallBackend(CeedQFunctionGetSourcePath(qf, &source_path));
1156     CeedCheck(source_path, ceed, CEED_ERROR_UNSUPPORTED, "/gpu/cuda/gen backend requires QFunction source code file");
1157 
1158     code << "// User QFunction source\n";
1159     code << "#include \"" << source_path << "\"\n\n";
1160   }
1161 
1162   // Setup
1163   code << "\n// -----------------------------------------------------------------------------\n";
1164   code << "// Operator Kernel\n";
1165   code << "// \n";
1166   code << "// d_[in,out]_i:   CeedVector device array\n";
1167   code << "// r_[in,out]_e_i: Element vector register\n";
1168   code << "// r_[in,out]_q_i: Quadrature space vector register\n";
1169   code << "// r_[in,out]_c_i: AtPoints Chebyshev coefficients register\n";
1170   code << "// r_[in,out]_s_i: Quadrature space slice vector register\n";
1171   code << "// \n";
1172   code << "// s_B_[in,out]_i: Interpolation matrix, shared memory\n";
1173   code << "// s_G_[in,out]_i: Gradient matrix, shared memory\n";
1174   code << "// -----------------------------------------------------------------------------\n";
1175   code << "extern \"C\" __global__ void " << operator_name
1176        << "(CeedInt num_elem, void* ctx, FieldsInt_Cuda indices, Fields_Cuda fields, Fields_Cuda B, Fields_Cuda G, CeedScalar *W, Points_Cuda "
1177           "points) {\n";
1178 
1179   // Scratch buffers
1180   for (CeedInt i = 0; i < num_input_fields; i++) {
1181     CeedEvalMode eval_mode;
1182 
1183     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
1184     if (eval_mode != CEED_EVAL_WEIGHT) {  // Skip CEED_EVAL_WEIGHT
1185       code << "  const CeedScalar *__restrict__ d_in_" << i << " = fields.inputs[" << i << "];\n";
1186     }
1187   }
1188   for (CeedInt i = 0; i < num_output_fields; i++) {
1189     code << "  CeedScalar *__restrict__ d_out_" << i << " = fields.outputs[" << i << "];\n";
1190   }
1191 
1192   code << "  const CeedInt dim = " << dim << ";\n";
1193   code << "  const CeedInt " << (is_tensor ? "Q_1d" : "Q") << " = " << Q_1d << ";\n";
1194   if (is_at_points) {
1195     code << "  const CeedInt max_num_points = " << max_num_points << ";\n";
1196     code << "  const CeedInt coords_comp_stride = " << coords_comp_stride << ";\n";
1197   }
1198 
1199   // Shared data
1200   code << "  extern __shared__ CeedScalar slice[];\n";
1201   code << "  SharedData_Cuda data;\n";
1202   code << "  data.t_id_x = threadIdx.x;\n";
1203   code << "  data.t_id_y = threadIdx.y;\n";
1204   code << "  data.t_id_z = threadIdx.z;\n";
1205   code << "  data.t_id  = threadIdx.x + threadIdx.y*blockDim.x + threadIdx.z*blockDim.y*blockDim.x;\n";
1206   code << "  data.slice = slice + data.t_id_z*OP_T_1D" << ((!is_tensor || dim == 1) ? "" : "*OP_T_1D") << ";\n";
1207 
1208   // -- Determine input mat reuse
1209   FieldReuse_Cuda input_matrix_reuse[CEED_FIELD_MAX];
1210 
1211   for (CeedInt i = 0; i < num_input_fields; i++) {
1212     input_matrix_reuse[i].index = -1;
1213   }
1214   for (CeedInt i = 0; i < num_input_fields; i++) {
1215     CeedEvalMode eval_mode_i;
1216     CeedBasis    basis_i;
1217 
1218     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode_i));
1219     if (eval_mode_i == CEED_EVAL_WEIGHT) continue;
1220     CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis_i));
1221     for (CeedInt j = 0; (input_matrix_reuse[i].index == -1) && (j < i); j++) {
1222       CeedEvalMode eval_mode_j;
1223       CeedBasis    basis_j;
1224 
1225       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[j], &eval_mode_j));
1226       if (eval_mode_j == CEED_EVAL_WEIGHT) continue;
1227       CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[j], &basis_j));
1228       if (basis_i == basis_j) {
1229         if (is_tensor) {
1230           input_matrix_reuse[i].index     = j;
1231           input_matrix_reuse[i].is_input  = true;
1232           input_matrix_reuse[i].eval_mode = eval_mode_j;
1233         } else {
1234           // For non-tensor can only re-use with the same eval mode
1235           if (eval_mode_i == eval_mode_j) {
1236             input_matrix_reuse[i].index     = j;
1237             input_matrix_reuse[i].is_input  = true;
1238             input_matrix_reuse[i].eval_mode = eval_mode_j;
1239           }
1240         }
1241       }
1242       CeedCallBackend(CeedBasisDestroy(&basis_j));
1243     }
1244     CeedCallBackend(CeedBasisDestroy(&basis_i));
1245   }
1246 
1247   // -- Determine output mat reuse
1248   FieldReuse_Cuda output_matrix_reuse[CEED_FIELD_MAX];
1249 
1250   for (CeedInt i = 0; i < num_output_fields; i++) {
1251     output_matrix_reuse[i].index = -1;
1252   }
1253   for (CeedInt i = 0; i < num_output_fields; i++) {
1254     CeedEvalMode eval_mode_i;
1255     CeedBasis    basis_i;
1256 
1257     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode_i));
1258     CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis_i));
1259     for (CeedInt j = 0; (output_matrix_reuse[i].index == -1) && (j < num_input_fields); j++) {
1260       CeedEvalMode eval_mode_j;
1261       CeedBasis    basis_j;
1262 
1263       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[j], &eval_mode_j));
1264       if (eval_mode_j == CEED_EVAL_WEIGHT) continue;
1265       CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[j], &basis_j));
1266       if (basis_i == basis_j) {
1267         if (is_tensor) {
1268           output_matrix_reuse[i].index     = j;
1269           output_matrix_reuse[i].is_input  = true;
1270           output_matrix_reuse[i].eval_mode = eval_mode_j;
1271         } else {
1272           // For non-tensor can only re-use with the same eval mode
1273           if (eval_mode_i == eval_mode_j) {
1274             output_matrix_reuse[i].index     = j;
1275             output_matrix_reuse[i].is_input  = true;
1276             output_matrix_reuse[i].eval_mode = eval_mode_j;
1277           }
1278         }
1279       }
1280       CeedCallBackend(CeedBasisDestroy(&basis_j));
1281     }
1282     for (CeedInt j = 0; (output_matrix_reuse[i].index == -1) && (j < i); j++) {
1283       CeedEvalMode eval_mode_j;
1284       CeedBasis    basis_j;
1285 
1286       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[j], &eval_mode_j));
1287       if (eval_mode_j == CEED_EVAL_WEIGHT) continue;
1288       CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[j], &basis_j));
1289       if (basis_i == basis_j) {
1290         if (is_tensor) {
1291           output_matrix_reuse[i].index     = j;
1292           output_matrix_reuse[i].is_input  = false;
1293           output_matrix_reuse[i].eval_mode = eval_mode_j;
1294         } else {
1295           // For non-tensor can only re-use with the same eval mode
1296           if (eval_mode_i == eval_mode_j) {
1297             output_matrix_reuse[i].index     = j;
1298             output_matrix_reuse[i].is_input  = false;
1299             output_matrix_reuse[i].eval_mode = eval_mode_j;
1300           }
1301         }
1302       }
1303       CeedCallBackend(CeedBasisDestroy(&basis_j));
1304     }
1305     CeedCallBackend(CeedBasisDestroy(&basis_i));
1306   }
1307 
1308   // Initialize constants, and matrices B and G
1309   code << "\n  // Input field constants and basis data\n";
1310   for (CeedInt i = 0; i < num_input_fields; i++) {
1311     CeedCallBackend(CeedOperatorBuildKernelFieldData_Cuda_gen(code, data, i, op_input_fields[i], qf_input_fields[i], input_matrix_reuse[i], Q_1d,
1312                                                               true, is_tensor, is_at_points, use_3d_slices));
1313   }
1314   code << "\n  // Output field constants and basis data\n";
1315   for (CeedInt i = 0; i < num_output_fields; i++) {
1316     CeedCallBackend(CeedOperatorBuildKernelFieldData_Cuda_gen(code, data, i, op_output_fields[i], qf_output_fields[i], output_matrix_reuse[i], Q_1d,
1317                                                               false, is_tensor, is_at_points, use_3d_slices));
1318   }
1319 
1320   // Loop over all elements
1321   code << "\n  // Element loop\n";
1322   code << "  __syncthreads();\n";
1323   code << "  for (CeedInt elem = blockIdx.x*blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x*blockDim.z) {\n";
1324 
1325   // -- Compute minimum buffer space needed
1326   CeedInt max_rstr_buffer_size = 1;
1327 
1328   for (CeedInt i = 0; i < num_input_fields; i++) {
1329     CeedInt             num_comp, elem_size;
1330     CeedElemRestriction elem_rstr;
1331 
1332     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr));
1333     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
1334     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
1335     max_rstr_buffer_size = CeedIntMax(max_rstr_buffer_size, num_comp * (is_tensor && (dim >= 3) ? elem_size : 1));
1336     CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
1337   }
1338   for (CeedInt i = 0; i < num_output_fields; i++) {
1339     CeedInt             num_comp, elem_size;
1340     CeedElemRestriction elem_rstr;
1341 
1342     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &elem_rstr));
1343     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
1344     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
1345     max_rstr_buffer_size = CeedIntMax(max_rstr_buffer_size, num_comp * (is_tensor && (dim >= 3) ? elem_size : 1));
1346     CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
1347   }
1348   code << "    // Scratch restriction buffer space\n";
1349   code << "    CeedScalar r_e_scratch[" << max_rstr_buffer_size << "];\n";
1350 
1351   // -- Determine best input field processing order
1352   CeedInt field_rstr_in_buffer[CEED_FIELD_MAX], input_field_order[CEED_FIELD_MAX];
1353 
1354   for (CeedInt i = 0; i < num_input_fields; i++) {
1355     field_rstr_in_buffer[i] = -1;
1356     input_field_order[i]    = -1;
1357   }
1358   {
1359     bool    is_ordered[CEED_FIELD_MAX];
1360     CeedInt curr_index = 0;
1361 
1362     for (CeedInt i = 0; i < num_input_fields; i++) is_ordered[i] = false;
1363     for (CeedInt i = 0; i < num_input_fields; i++) {
1364       CeedVector          vec_i;
1365       CeedElemRestriction rstr_i;
1366 
1367       if (is_ordered[i]) continue;
1368       field_rstr_in_buffer[i]       = i;
1369       is_ordered[i]                 = true;
1370       input_field_order[curr_index] = i;
1371       curr_index++;
1372       CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec_i));
1373       if (vec_i == CEED_VECTOR_NONE) continue;  // CEED_EVAL_WEIGHT
1374       CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &rstr_i));
1375       for (CeedInt j = i + 1; j < num_input_fields; j++) {
1376         CeedVector          vec_j;
1377         CeedElemRestriction rstr_j;
1378 
1379         CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[j], &vec_j));
1380         CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[j], &rstr_j));
1381         if (rstr_i == rstr_j && vec_i == vec_j) {
1382           field_rstr_in_buffer[j]       = i;
1383           is_ordered[j]                 = true;
1384           input_field_order[curr_index] = j;
1385           curr_index++;
1386         }
1387         CeedCallBackend(CeedVectorDestroy(&vec_j));
1388         CeedCallBackend(CeedElemRestrictionDestroy(&rstr_j));
1389       }
1390       CeedCallBackend(CeedVectorDestroy(&vec_i));
1391       CeedCallBackend(CeedElemRestrictionDestroy(&rstr_i));
1392     }
1393   }
1394 
1395   // -- Input restriction and basis
1396   code << "\n    // -- Input field restrictions and basis actions\n";
1397   for (CeedInt i = 0; i < num_input_fields; i++) {
1398     const char   *field_name;
1399     const CeedInt f = input_field_order[i];
1400 
1401     CeedCallBackend(CeedOperatorFieldGetName(op_input_fields[f], &field_name));
1402     code << "    // ---- Input field " << f << ": " << field_name << "\n";
1403 
1404     // ---- Restriction
1405     CeedCallBackend(CeedOperatorBuildKernelRestriction_Cuda_gen(code, data, f, dim, field_rstr_in_buffer, op_input_fields[f], qf_input_fields[f],
1406                                                                 Q_1d, true, is_tensor, is_at_points, use_3d_slices));
1407 
1408     // ---- Basis action
1409     CeedCallBackend(CeedOperatorBuildKernelBasis_Cuda_gen(code, data, f, dim, op_input_fields[f], qf_input_fields[f], Q_1d, true, is_tensor,
1410                                                           is_at_points, use_3d_slices));
1411   }
1412 
1413   // -- Q function
1414   CeedCallBackend(CeedOperatorBuildKernelQFunction_Cuda_gen(code, data, dim, max_num_points, num_input_fields, op_input_fields, qf_input_fields,
1415                                                             num_output_fields, op_output_fields, qf_output_fields, qfunction_name, Q_1d, is_tensor,
1416                                                             is_at_points, use_3d_slices));
1417 
1418   // -- Output basis and restriction
1419   code << "\n    // -- Output field basis action and restrictions\n";
1420   for (CeedInt i = 0; i < num_output_fields; i++) {
1421     const char *field_name;
1422 
1423     CeedCallBackend(CeedOperatorFieldGetName(op_output_fields[i], &field_name));
1424     code << "    // ---- Output field " << i << ": " << field_name << "\n";
1425 
1426     // ---- Basis action
1427     CeedCallBackend(CeedOperatorBuildKernelBasis_Cuda_gen(code, data, i, dim, op_output_fields[i], qf_output_fields[i], Q_1d, false, is_tensor,
1428                                                           is_at_points, use_3d_slices));
1429 
1430     // ---- Restriction
1431     CeedCallBackend(CeedOperatorBuildKernelRestriction_Cuda_gen(code, data, i, dim, NULL, op_output_fields[i], qf_output_fields[i], Q_1d, false,
1432                                                                 is_tensor, is_at_points, use_3d_slices));
1433   }
1434 
1435   // Close loop and function
1436   code << "  }\n";
1437   code << "}\n";
1438   code << "// -----------------------------------------------------------------------------\n\n";
1439 
1440   // Compile
1441   {
1442     bool is_compile_good = false;
1443 
1444     CeedCallBackend(CeedTryCompile_Cuda(ceed, code.str().c_str(), &is_compile_good, &data->module, 1, "OP_T_1D", CeedIntMax(Q_1d, data->max_P_1d)));
1445     if (is_compile_good) {
1446       *is_good_build = true;
1447       CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, operator_name.c_str(), &data->op));
1448     } else {
1449       *is_good_build     = false;
1450       data->use_fallback = true;
1451     }
1452   }
1453   CeedCallBackend(CeedOperatorSetSetupDone(op));
1454   CeedCallBackend(CeedDestroy(&ceed));
1455   CeedCallBackend(CeedQFunctionDestroy(&qf));
1456   return CEED_ERROR_SUCCESS;
1457 }
1458 
1459 //------------------------------------------------------------------------------
1460