xref: /libCEED/backends/cuda-gen/ceed-cuda-gen-operator-build.cpp (revision 1a8516d00062e8132c3db0515cc9f5fa064f6664)
1 // Copyright (c) 2017-2025, 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/gen-tools.h>
13 #include <ceed/jit-tools.h>
14 #include <cuda_runtime.h>
15 
16 #include <iostream>
17 #include <sstream>
18 #include <string>
19 
20 #include "../cuda-ref/ceed-cuda-ref.h"
21 #include "../cuda-shared/ceed-cuda-shared.h"
22 #include "../cuda/ceed-cuda-common.h"
23 #include "../cuda/ceed-cuda-compile.h"
24 #include "ceed-cuda-gen.h"
25 
26 struct FieldReuse_Cuda {
27   CeedInt      index;
28   bool         is_input;
29   CeedEvalMode eval_mode;
30 };
31 
32 //------------------------------------------------------------------------------
33 // Determine type of operator
34 //------------------------------------------------------------------------------
35 static int CeedOperatorBuildKernelData_Cuda_gen(Ceed ceed, CeedInt num_input_fields, CeedOperatorField *op_input_fields,
36                                                 CeedQFunctionField *qf_input_fields, CeedInt num_output_fields, CeedOperatorField *op_output_fields,
37                                                 CeedQFunctionField *qf_output_fields, CeedInt *max_P, CeedInt *max_P_1d, CeedInt *Q, CeedInt *Q_1d,
38                                                 CeedInt *max_dim, bool *is_all_tensor, bool *use_3d_slices) {
39   // Check if all are tensor
40   *is_all_tensor = true;
41   for (CeedInt i = 0; i < num_input_fields; i++) {
42     CeedBasis basis;
43 
44     CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
45     if (basis != CEED_BASIS_NONE) {
46       bool is_field_tensor;
47 
48       CeedCallBackend(CeedBasisIsTensor(basis, &is_field_tensor));
49       *is_all_tensor = *is_all_tensor && is_field_tensor;
50     }
51     CeedCallBackend(CeedBasisDestroy(&basis));
52   }
53   for (CeedInt i = 0; i < num_output_fields; i++) {
54     CeedBasis basis;
55 
56     CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
57     if (basis != CEED_BASIS_NONE) {
58       bool is_field_tensor;
59 
60       CeedCallBackend(CeedBasisIsTensor(basis, &is_field_tensor));
61       *is_all_tensor = *is_all_tensor && is_field_tensor;
62     }
63     CeedCallBackend(CeedBasisDestroy(&basis));
64   }
65 
66   // Find max_P, max_P_1d, Q, and Q_1d
67   bool is_all_3d = true;
68 
69   *max_P    = 0;
70   *max_P_1d = 0;
71   *Q        = 0;
72   *Q_1d     = 0;
73   for (CeedInt i = 0; i < num_input_fields; i++) {
74     CeedBasis basis;
75 
76     CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
77     if (basis != CEED_BASIS_NONE) {
78       bool    is_field_tensor;
79       CeedInt field_dim = 0, field_P = 0, field_P_1d = 0, field_Q = 0, field_Q_1d = 0;
80 
81       // Check if 3D
82       CeedCallBackend(CeedBasisGetDimension(basis, &field_dim));
83       is_all_3d = is_all_3d && (field_dim == 3);
84       *max_dim  = CeedIntMax(*max_dim, field_dim);
85 
86       // Collect P, P_1d, Q, and Q_1d
87       CeedCallBackend(CeedBasisGetNumNodes(basis, &field_P));
88       *max_P = CeedIntMax(*max_P, field_P);
89       CeedCallBackend(CeedBasisIsTensor(basis, &is_field_tensor));
90       if (is_field_tensor) {
91         CeedCallBackend(CeedBasisGetNumNodes1D(basis, &field_P_1d));
92         *max_P_1d = CeedIntMax(*max_P_1d, field_P_1d);
93       }
94       CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis, &field_Q));
95       CeedCheck(*Q == 0 || field_Q == *Q, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible");
96       *Q = field_Q;
97       if (is_field_tensor) {
98         CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &field_Q_1d));
99         CeedCheck(*Q_1d == 0 || field_Q_1d == *Q_1d, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible");
100         *Q_1d = field_Q_1d;
101       }
102     }
103     CeedCallBackend(CeedBasisDestroy(&basis));
104   }
105   for (CeedInt i = 0; i < num_output_fields; i++) {
106     CeedBasis basis;
107 
108     CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
109     if (basis != CEED_BASIS_NONE) {
110       bool    is_field_tensor;
111       CeedInt field_dim = 0, field_P = 0, field_P_1d = 0, field_Q = 0, field_Q_1d = 0;
112 
113       // Check if 3D
114       CeedCallBackend(CeedBasisGetDimension(basis, &field_dim));
115       is_all_3d = is_all_3d && (field_dim == 3);
116       *max_dim  = CeedIntMax(*max_dim, field_dim);
117 
118       // Collect P, P_1d, Q, and Q_1d
119       CeedCallBackend(CeedBasisGetNumNodes(basis, &field_P));
120       *max_P = CeedIntMax(*max_P, field_P);
121       CeedCallBackend(CeedBasisIsTensor(basis, &is_field_tensor));
122       if (is_field_tensor) {
123         CeedCallBackend(CeedBasisGetNumNodes1D(basis, &field_P_1d));
124         *max_P_1d = CeedIntMax(*max_P_1d, field_P_1d);
125       }
126       CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis, &field_Q));
127       CeedCheck(*Q == 0 || field_Q == *Q, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible");
128       *Q = field_Q;
129       if (is_field_tensor) {
130         CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &field_Q_1d));
131         CeedCheck(*Q_1d == 0 || field_Q_1d == *Q_1d, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible");
132         *Q_1d = field_Q_1d;
133       }
134     }
135     CeedCallBackend(CeedBasisDestroy(&basis));
136   }
137 
138   // Only use 3D collocated gradient parallelization strategy when gradient is computed
139   *use_3d_slices = false;
140   if (is_all_3d && *is_all_tensor) {
141     bool was_grad_found = false;
142 
143     for (CeedInt i = 0; i < num_input_fields; i++) {
144       CeedEvalMode eval_mode;
145 
146       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
147       if (eval_mode == CEED_EVAL_GRAD) {
148         CeedBasis_Cuda_shared *basis_data;
149         CeedBasis              basis;
150 
151         CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
152         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
153         *use_3d_slices = basis_data->d_collo_grad_1d && (was_grad_found ? *use_3d_slices : true);
154         was_grad_found = true;
155         CeedCallBackend(CeedBasisDestroy(&basis));
156       }
157     }
158     for (CeedInt i = 0; i < num_output_fields; i++) {
159       CeedEvalMode eval_mode;
160 
161       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
162       if (eval_mode == CEED_EVAL_GRAD) {
163         CeedBasis_Cuda_shared *basis_data;
164         CeedBasis              basis;
165 
166         CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
167         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
168         *use_3d_slices = basis_data->d_collo_grad_1d && (was_grad_found ? *use_3d_slices : true);
169         was_grad_found = true;
170         CeedCallBackend(CeedBasisDestroy(&basis));
171       }
172     }
173   }
174   return CEED_ERROR_SUCCESS;
175 }
176 
177 //------------------------------------------------------------------------------
178 // Setup fields
179 //------------------------------------------------------------------------------
180 static int CeedOperatorBuildKernelFieldData_Cuda_gen(std::ostringstream &code, CeedOperator_Cuda_gen *data, Tab &tab, CeedInt i,
181                                                      CeedOperatorField op_field, CeedQFunctionField qf_field, FieldReuse_Cuda field_reuse,
182                                                      CeedInt max_dim, CeedInt Q, CeedInt Q_1d, bool is_input, bool is_all_tensor, bool is_at_points,
183                                                      bool use_3d_slices, bool skip_active_load) {
184   bool      is_tensor = true, is_active = true;
185   CeedBasis basis;
186 
187   CeedCallBackend(CeedOperatorFieldGetBasis(op_field, &basis));
188   if (basis != CEED_BASIS_NONE) CeedCallBackend(CeedBasisIsTensor(basis, &is_tensor));
189   {
190     CeedVector vec;
191 
192     CeedCallBackend(CeedOperatorFieldGetVector(op_field, &vec));
193     is_active = vec == CEED_VECTOR_ACTIVE;
194     CeedCallBackend(CeedVectorDestroy(&vec));
195   }
196 
197   const char            *field_name;
198   std::string            var_suffix = (is_input ? "_in_" : "_out_") + std::to_string(i);
199   std::string            P_name = (is_tensor ? "P_1d" : "P") + var_suffix, Q_name = is_tensor ? "Q_1d" : "Q";
200   std::string            option_name = (is_input ? "inputs" : "outputs");
201   CeedEvalMode           eval_mode   = CEED_EVAL_NONE;
202   CeedInt                elem_size = 0, num_comp = 0, dim = max_dim, P_1d = 0;
203   CeedElemRestriction    elem_rstr;
204   CeedBasis_Cuda_shared *basis_data;
205 
206   // Field reuse info
207   bool use_previous_field = field_reuse.index != -1;
208 
209   CeedCallBackend(CeedOperatorFieldGetName(op_field, &field_name));
210   code << tab << "// -- " << (is_input ? "Input" : "Output") << " field " << i << ": " << field_name << "\n";
211 
212   // Get field data
213   CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_field, &elem_rstr));
214   if (elem_rstr != CEED_ELEMRESTRICTION_NONE) {
215     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
216     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
217   }
218   CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
219   if (basis != CEED_BASIS_NONE) {
220     CeedCallBackend(CeedBasisGetData(basis, &basis_data));
221     CeedCallBackend(CeedBasisGetDimension(basis, &dim));
222     if (is_tensor) CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
223     else CeedCallBackend(CeedBasisGetNumNodes(basis, &P_1d));
224   }
225   CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_field, &eval_mode));
226 
227   // Set field constants
228   code << tab << "const CeedInt dim" << var_suffix << " = " << dim << ";\n";
229   if (is_tensor && !is_all_tensor) {
230     CeedInt P = 0;
231 
232     CeedCallBackend(CeedBasisGetNumNodes(basis, &P));
233     code << tab << "const CeedInt P" << var_suffix << " = " << (basis == CEED_BASIS_NONE ? Q : P) << ";\n";
234   }
235   code << tab << "const CeedInt " << P_name << " = " << (basis == CEED_BASIS_NONE ? Q_1d : P_1d) << ";\n";
236   if (eval_mode != CEED_EVAL_WEIGHT) {
237     code << tab << "const CeedInt num_comp" << var_suffix << " = " << num_comp << ";\n";
238   }
239 
240   // Load basis data
241   code << tab << "// EvalMode: " << CeedEvalModes[eval_mode] << "\n";
242   switch (eval_mode) {
243     case CEED_EVAL_NONE:
244       break;
245     case CEED_EVAL_INTERP:
246       if (is_at_points) {
247         // AtPoints
248         if (!basis_data->d_chebyshev_interp_1d) {
249           CeedSize    interp_bytes;
250           CeedScalar *chebyshev_interp_1d;
251 
252           interp_bytes = P_1d * Q_1d * sizeof(CeedScalar);
253           CeedCallBackend(CeedCalloc(P_1d * Q_1d, &chebyshev_interp_1d));
254           CeedCallBackend(CeedBasisGetChebyshevInterp1D(basis, chebyshev_interp_1d));
255           CeedCallCuda(CeedBasisReturnCeed(basis), cudaMalloc((void **)&basis_data->d_chebyshev_interp_1d, interp_bytes));
256           CeedCallCuda(CeedBasisReturnCeed(basis),
257                        cudaMemcpy(basis_data->d_chebyshev_interp_1d, chebyshev_interp_1d, interp_bytes, cudaMemcpyHostToDevice));
258           CeedCallBackend(CeedFree(&chebyshev_interp_1d));
259         }
260         if (is_input) data->B.inputs[i] = basis_data->d_chebyshev_interp_1d;
261         else data->B.outputs[i] = basis_data->d_chebyshev_interp_1d;
262       } else {
263         // Standard quadrature
264         if (is_input) data->B.inputs[i] = basis_data->d_interp_1d;
265         else data->B.outputs[i] = basis_data->d_interp_1d;
266       }
267       if (use_previous_field && !skip_active_load) {
268         std::string reuse_var = "s_B" + ((field_reuse.is_input ? "_in_" : "_out_") + std::to_string(field_reuse.index));
269 
270         code << tab << "CeedScalar *s_B" << var_suffix << " = " << reuse_var << ";\n";
271       } else {
272         bool is_collocated = false;
273 
274         CeedCallBackend(CeedBasisIsCollocated(basis, &is_collocated));
275         if ((is_active && skip_active_load) || (is_collocated && !is_at_points)) {
276           code << tab << "CeedScalar *s_B" << var_suffix << " = NULL;\n";
277         } else {
278           code << tab << "__shared__ CeedScalar s_B" << var_suffix << "[" << P_name << "*" << Q_name << "];\n";
279           code << tab << "LoadMatrix<" << P_name << ", " << Q_name << ">(data, B." << option_name << "[" << i << "], s_B" << var_suffix << ");\n";
280         }
281       }
282       break;
283     case CEED_EVAL_GRAD:
284       if (is_at_points) {
285         // AtPoints
286         if (!basis_data->d_chebyshev_interp_1d) {
287           CeedSize    interp_bytes;
288           CeedScalar *chebyshev_interp_1d;
289 
290           interp_bytes = P_1d * Q_1d * sizeof(CeedScalar);
291           CeedCallBackend(CeedCalloc(P_1d * Q_1d, &chebyshev_interp_1d));
292           CeedCallBackend(CeedBasisGetChebyshevInterp1D(basis, chebyshev_interp_1d));
293           CeedCallCuda(CeedBasisReturnCeed(basis), cudaMalloc((void **)&basis_data->d_chebyshev_interp_1d, interp_bytes));
294           CeedCallCuda(CeedBasisReturnCeed(basis),
295                        cudaMemcpy(basis_data->d_chebyshev_interp_1d, chebyshev_interp_1d, interp_bytes, cudaMemcpyHostToDevice));
296           CeedCallBackend(CeedFree(&chebyshev_interp_1d));
297         }
298         if (is_input) data->B.inputs[i] = basis_data->d_chebyshev_interp_1d;
299         else data->B.outputs[i] = basis_data->d_chebyshev_interp_1d;
300       } else {
301         // Standard quadrature
302         if (is_input) data->B.inputs[i] = basis_data->d_interp_1d;
303         else data->B.outputs[i] = basis_data->d_interp_1d;
304       }
305       if (is_tensor) {
306         if (use_previous_field && !skip_active_load) {
307           std::string reuse_var = "s_B" + ((field_reuse.is_input ? "_in_" : "_out_") + std::to_string(field_reuse.index));
308 
309           code << tab << "CeedScalar *s_B" << var_suffix << " = " << reuse_var << ";\n";
310         } else {
311           bool is_collocated = false;
312 
313           CeedCallBackend(CeedBasisIsCollocated(basis, &is_collocated));
314           if ((is_active && skip_active_load) || (is_collocated && !is_at_points)) {
315             code << tab << "CeedScalar *s_B" << var_suffix << " = NULL;\n";
316           } else {
317             code << tab << "__shared__ CeedScalar s_B" << var_suffix << "[" << P_name << "*" << Q_name << "];\n";
318             code << tab << "LoadMatrix<" << P_name << ", " << Q_name << ">(data, B." << option_name << "[" << i << "], s_B" << var_suffix << ");\n";
319           }
320         }
321       }
322       if (is_at_points) break;  // No G mat for AtPoints
323       if (use_3d_slices) {
324         if (is_input) data->G.inputs[i] = basis_data->d_collo_grad_1d;
325         else data->G.outputs[i] = basis_data->d_collo_grad_1d;
326         if (use_previous_field && field_reuse.eval_mode == CEED_EVAL_GRAD && !skip_active_load) {
327           std::string reuse_var = "s_G" + ((field_reuse.is_input ? "_in_" : "_out_") + std::to_string(field_reuse.index));
328 
329           code << tab << "CeedScalar *s_G" << var_suffix << " = " << reuse_var << ";\n";
330         } else if (is_active && skip_active_load) {
331           code << tab << "CeedScalar *s_G" << var_suffix << " = NULL;\n";
332         } else {
333           code << tab << "__shared__ CeedScalar s_G" << var_suffix << "[" << Q_name << "*" << Q_name << "];\n";
334           code << tab << "LoadMatrix<" << Q_name << ", " << Q_name << ">(data, G." << option_name << "[" << i << "], s_G" << var_suffix << ");\n";
335         }
336       } else {
337         bool has_collo_grad = basis_data->d_collo_grad_1d;
338 
339         if (is_input) data->G.inputs[i] = has_collo_grad ? basis_data->d_collo_grad_1d : basis_data->d_grad_1d;
340         else data->G.outputs[i] = has_collo_grad ? basis_data->d_collo_grad_1d : basis_data->d_grad_1d;
341         if (has_collo_grad) {
342           if (use_previous_field && field_reuse.eval_mode == CEED_EVAL_GRAD && !skip_active_load) {
343             std::string reuse_var = "s_G" + ((field_reuse.is_input ? "_in_" : "_out_") + std::to_string(field_reuse.index));
344 
345             code << tab << "CeedScalar *s_G" << var_suffix << " = " << reuse_var << ";\n";
346           } else if (is_active && skip_active_load) {
347             code << tab << "CeedScalar *s_G" << var_suffix << " = NULL;\n";
348           } else {
349             code << tab << "__shared__ CeedScalar s_G" << var_suffix << "[" << Q_name << "*" << Q_name << "];\n";
350             code << tab << "LoadMatrix<" << Q_name << ", " << Q_name << ">(data, G." << option_name << "[" << i << "], s_G" << var_suffix << ");\n";
351           }
352         } else {
353           if (use_previous_field && field_reuse.eval_mode == CEED_EVAL_GRAD && !skip_active_load) {
354             std::string reuse_var = "s_G" + ((field_reuse.is_input ? "_in_" : "_out_") + std::to_string(field_reuse.index));
355 
356             code << tab << "CeedScalar *s_G" << var_suffix << " = " << reuse_var << ";\n";
357           } else if (is_active && skip_active_load) {
358             code << tab << "CeedScalar *s_G" << var_suffix << " = NULL;\n";
359           } else {
360             code << tab << "__shared__ CeedScalar s_G" << var_suffix << "[" << P_name << "*" << Q_name << (is_tensor ? "" : "*dim")
361                  << (is_tensor ? "" : var_suffix) << "];\n";
362             code << tab << "LoadMatrix<" << P_name << ", " << Q_name << (is_tensor ? "" : "*dim") << (is_tensor ? "" : var_suffix) << ">(data, G."
363                  << option_name << "[" << i << "], s_G" << var_suffix << ");\n";
364           }
365         }
366       }
367       break;
368     case CEED_EVAL_WEIGHT:
369       break;  // No action
370       // LCOV_EXCL_START
371     case CEED_EVAL_DIV:
372     case CEED_EVAL_CURL:
373       break;  // TODO: Not implemented
374               // LCOV_EXCL_STOP
375   }
376   CeedCallBackend(CeedBasisDestroy(&basis));
377   return CEED_ERROR_SUCCESS;
378 }
379 
380 //------------------------------------------------------------------------------
381 // Restriction
382 //------------------------------------------------------------------------------
383 static int CeedOperatorBuildKernelRestriction_Cuda_gen(std::ostringstream &code, CeedOperator_Cuda_gen *data, Tab &tab, CeedInt i,
384                                                        CeedInt field_input_buffer[], CeedOperatorField op_field, CeedQFunctionField qf_field,
385                                                        CeedInt max_dim, CeedInt Q_1d, bool is_input, bool is_all_tensor, bool is_at_points,
386                                                        bool use_3d_slices) {
387   std::string               var_suffix = (is_input ? "_in_" : "_out_") + std::to_string(i);
388   std::string               P_name     = (is_all_tensor ? "P_1d" : "P") + var_suffix;
389   CeedEvalMode              eval_mode  = CEED_EVAL_NONE;
390   CeedInt                   elem_size = 0, num_comp = 0;
391   CeedSize                  l_size;
392   CeedRestrictionType       rstr_type = CEED_RESTRICTION_STANDARD;
393   CeedElemRestriction_Cuda *rstr_data;
394   CeedElemRestriction       elem_rstr;
395 
396   // Get field data
397   CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_field, &elem_rstr));
398   if (elem_rstr != CEED_ELEMRESTRICTION_NONE) {
399     CeedCallBackend(CeedElemRestrictionGetType(elem_rstr, &rstr_type));
400     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
401     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
402     CeedCallBackend(CeedElemRestrictionGetData(elem_rstr, &rstr_data));
403   }
404   CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_field, &eval_mode));
405 
406   // Restriction
407   if (is_input) {
408     // Input
409     if (field_input_buffer[i] != i) {
410       std::string buffer_name = "r_e_in_" + std::to_string(field_input_buffer[i]);
411 
412       // Restriction was already done for previous input
413       code << tab << "CeedScalar *r_e" << var_suffix << " = " << buffer_name << ";\n";
414     } else if (eval_mode != CEED_EVAL_WEIGHT && !((eval_mode == CEED_EVAL_NONE) && use_3d_slices && is_at_points)) {
415       if (eval_mode == CEED_EVAL_NONE && rstr_type != CEED_RESTRICTION_POINTS) {
416         // No basis action, so r_e_in_* in also r_q_in_* and needs to be allocated
417         code << tab << "CeedScalar r_e" << var_suffix << "[num_comp" << var_suffix << "*" << P_name << "];\n";
418       } else if (rstr_type != CEED_RESTRICTION_POINTS) {
419         // Otherwise we're using the scratch space
420         code << tab << "CeedScalar *r_e" << var_suffix << " = r_e_scratch;\n";
421       }
422       switch (rstr_type) {
423         case CEED_RESTRICTION_STANDARD: {
424           CeedInt comp_stride;
425 
426           CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size));
427           code << tab << "const CeedInt l_size" << var_suffix << " = " << l_size << ";\n";
428           CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
429           code << tab << "const CeedInt comp_stride" << var_suffix << " = " << comp_stride << ";\n";
430           data->indices.inputs[i] = (CeedInt *)rstr_data->d_offsets;
431           code << tab << "ReadLVecStandard" << (is_all_tensor ? max_dim : 1) << "d<num_comp" << var_suffix << ", comp_stride" << var_suffix << ", "
432                << P_name << ">(data, l_size" << var_suffix << ", elem, indices.inputs[" << i << "], d" << var_suffix << ", r_e" << var_suffix
433                << ");\n";
434           break;
435         }
436         case CEED_RESTRICTION_STRIDED: {
437           bool    has_backend_strides;
438           CeedInt num_elem;
439 
440           CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides));
441           CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem));
442           CeedInt strides[3] = {1, elem_size * num_elem, elem_size};
443 
444           if (!has_backend_strides) {
445             CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides));
446           }
447           code << tab << "const CeedInt strides" << var_suffix << "_0 = " << strides[0] << ", strides" << var_suffix << "_1 = " << strides[1]
448                << ", strides" << var_suffix << "_2 = " << strides[2] << ";\n";
449           code << tab << "ReadLVecStrided" << (is_all_tensor ? max_dim : 1) << "d<num_comp" << var_suffix << ", " << P_name << ", strides"
450                << var_suffix << "_0, strides" << var_suffix << "_1, strides" << var_suffix << "_2>(data, elem, d" << var_suffix << ", r_e"
451                << var_suffix << ");\n";
452           break;
453         }
454         case CEED_RESTRICTION_POINTS: {
455           CeedInt comp_stride;
456 
457           CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
458           code << tab << "const CeedInt comp_stride" << var_suffix << " = " << comp_stride << ";\n";
459           data->indices.inputs[i] = (CeedInt *)rstr_data->d_offsets;
460           break;
461         }
462         // LCOV_EXCL_START
463         case CEED_RESTRICTION_ORIENTED:
464         case CEED_RESTRICTION_CURL_ORIENTED:
465           break;  // TODO: Not implemented
466                   // LCOV_EXCL_STOP
467       }
468     }
469   } else {
470     // Output
471     switch (rstr_type) {
472       case CEED_RESTRICTION_STANDARD: {
473         CeedInt comp_stride;
474 
475         CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size));
476         code << tab << "const CeedInt l_size" << var_suffix << " = " << l_size << ";\n";
477         CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
478         code << tab << "const CeedInt comp_stride" << var_suffix << " = " << comp_stride << ";\n";
479         data->indices.outputs[i] = (CeedInt *)rstr_data->d_offsets;
480         code << tab << "WriteLVecStandard" << (is_all_tensor ? max_dim : 1) << "d<num_comp" << var_suffix << ", comp_stride" << var_suffix << ", "
481              << P_name << ">(data, l_size" << var_suffix << ", elem, indices.outputs[" << i << "], r_e" << var_suffix << ", d" << var_suffix
482              << ");\n";
483         break;
484       }
485       case CEED_RESTRICTION_STRIDED: {
486         bool    has_backend_strides;
487         CeedInt num_elem;
488 
489         CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides));
490         CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem));
491         CeedInt strides[3] = {1, elem_size * num_elem, elem_size};
492 
493         if (!has_backend_strides) {
494           CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides));
495         }
496         code << tab << "const CeedInt strides" << var_suffix << "_0 = " << strides[0] << ", strides" << var_suffix << "_1 = " << strides[1]
497              << ", strides" << var_suffix << "_2 = " << strides[2] << ";\n";
498         code << tab << "WriteLVecStrided" << (is_all_tensor ? max_dim : 1) << "d<num_comp" << var_suffix << ", " << P_name << ", strides"
499              << var_suffix << "_0, strides" << var_suffix << "_1, strides" << var_suffix << "_2>(data, elem, r_e" << var_suffix << ", d" << var_suffix
500              << ");\n";
501         break;
502       }
503       case CEED_RESTRICTION_POINTS:
504         data->indices.outputs[i] = (CeedInt *)rstr_data->d_offsets;
505         break;
506       // LCOV_EXCL_START
507       case CEED_RESTRICTION_ORIENTED:
508       case CEED_RESTRICTION_CURL_ORIENTED:
509         break;  // TODO: Not implemented
510                 // LCOV_EXCL_STOP
511     }
512   }
513   CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
514   return CEED_ERROR_SUCCESS;
515 }
516 
517 //------------------------------------------------------------------------------
518 // Basis
519 //------------------------------------------------------------------------------
520 static int CeedOperatorBuildKernelBasis_Cuda_gen(std::ostringstream &code, CeedOperator_Cuda_gen *data, Tab &tab, CeedInt i,
521                                                  CeedOperatorField op_field, CeedQFunctionField qf_field, CeedInt max_dim, CeedInt Q_1d,
522                                                  bool is_input, bool is_all_tensor, bool is_at_points, bool use_3d_slices) {
523   bool      is_tensor = true, is_collocated = true;
524   CeedBasis basis;
525   CeedCallBackend(CeedOperatorFieldGetBasis(op_field, &basis));
526   CeedCallBackend(CeedBasisIsTensor(basis, &is_tensor));
527   CeedCallBackend(CeedBasisIsCollocated(basis, &is_collocated));
528 
529   std::string         var_suffix = (is_input ? "_in_" : "_out_") + std::to_string(i);
530   std::string         P_name = (is_tensor ? "P_1d" : "P") + var_suffix, Q_name = is_tensor ? "Q_1d" : "Q";
531   CeedEvalMode        eval_mode = CEED_EVAL_NONE;
532   CeedInt             dim = max_dim, elem_size = 0, num_comp = 0, P_1d = 0;
533   CeedElemRestriction elem_rstr;
534 
535   // Get field data
536   CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_field, &elem_rstr));
537   if (elem_rstr != CEED_ELEMRESTRICTION_NONE) {
538     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
539     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
540   }
541   CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
542   if (basis != CEED_BASIS_NONE) {
543     CeedCallBackend(CeedBasisGetDimension(basis, &dim));
544     if (is_tensor) CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
545     else CeedCallBackend(CeedBasisGetNumNodes(basis, &P_1d));
546   }
547   CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_field, &eval_mode));
548 
549   // Basis
550   code << tab << "// EvalMode: " << CeedEvalModes[eval_mode] << "\n";
551   if (is_input) {
552     switch (eval_mode) {
553       case CEED_EVAL_NONE:
554         if (!use_3d_slices && !is_at_points) {
555           code << tab << "CeedScalar *r_q" << var_suffix << " = r_e" << var_suffix << ";\n";
556         }
557         break;
558       case CEED_EVAL_INTERP:
559         if (is_at_points) {
560           std::string function_name = (dim == 1 ? "Interp" : "InterpTensor") + std::to_string(dim) + "d";
561 
562           code << tab << "CeedScalar r_c" << var_suffix << "[num_comp" << var_suffix << "*" << (dim >= 3 ? Q_name : "1") << "];\n";
563           code << tab << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ", OP_T_1D>(data, r_e" << var_suffix
564                << ", s_B" << var_suffix << ", r_c" << var_suffix << ");\n";
565         } else {
566           std::string function_name = is_tensor ? ((dim == 1 ? "Interp" : "InterpTensor") + std::string(is_collocated ? "CollocatedNodes" : "") +
567                                                    std::to_string(dim) + "d" + (is_all_tensor ? "" : "Flattened"))
568                                                 : "InterpNonTensor";
569           std::string op_t_1d_name  = (is_all_tensor || !is_tensor) ? "OP_T_1D" : (P_1d > Q_1d ? P_name : Q_name);
570 
571           code << tab << "CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << (is_all_tensor && (dim >= 3) ? Q_name : "1") << "];\n";
572           code << tab << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ", " << op_t_1d_name << ">(data, r_e"
573                << var_suffix << ", s_B" << var_suffix << ", r_q" << var_suffix << ");\n";
574         }
575         break;
576       case CEED_EVAL_GRAD:
577         if (is_at_points) {
578           std::string function_name = (dim == 1 ? "Interp" : "InterpTensor") + std::to_string(dim) + "d";
579 
580           code << tab << "CeedScalar r_c" << var_suffix << "[num_comp" << var_suffix << "*" << (dim >= 3 ? Q_name : "1") << "];\n";
581           code << tab << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ", OP_T_1D>(data, r_e" << var_suffix
582                << ", s_B" << var_suffix << ", r_c" << var_suffix << ");\n";
583         } else if (use_3d_slices) {
584           std::string function_name =
585               (dim > 1 ? "InterpTensor" : "Interp") + std::string(is_collocated ? "CollocatedNodes" : "") + std::to_string(dim) + "d";
586 
587           code << tab << "CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << Q_name << "];\n";
588           code << tab << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ", OP_T_1D>(data, r_e" << var_suffix
589                << ", s_B" << var_suffix << ", r_q" << var_suffix << ");\n";
590         } else if (is_tensor) {
591           bool        is_collocated_grad = dim == 3 && Q_1d >= P_1d;
592           std::string function_name =
593               (dim == 1 ? "Grad" : ("GradTensor" + std::string(is_collocated ? "CollocatedNodes" : (is_collocated_grad ? "Collocated" : "")))) +
594               std::to_string(dim) + "d" + (is_all_tensor ? "" : "Flattened");
595           std::string op_t_1d_name = is_all_tensor ? "OP_T_1D" : (P_1d > Q_1d ? P_name : Q_name);
596 
597           code << tab << "CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*dim" << var_suffix << "*"
598                << (is_all_tensor && dim >= 3 ? Q_name : "1") << "];\n";
599           code << tab << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ", " << op_t_1d_name << ">(data, r_e"
600                << var_suffix << ", s_B" << var_suffix << ", s_G" << var_suffix << ", r_q" << var_suffix << ");\n";
601         } else {
602           std::string function_name = "GradNonTensor";
603 
604           code << tab << "CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*dim" << var_suffix << "];\n";
605           code << tab << function_name << "<num_comp" << var_suffix << ", dim" << var_suffix << ", " << P_name << ", " << Q_name
606                << ", OP_T_1D>(data, r_e" << var_suffix << ", s_G" << var_suffix << ", r_q" << var_suffix << ");\n";
607         }
608         break;
609       case CEED_EVAL_WEIGHT: {
610         if (is_at_points) {
611           code << tab << "// Nothing to do AtPoints\n";
612         } else {
613           CeedBasis_Cuda_shared *basis_data;
614           std::string            function_name = is_tensor
615                                                      ? ((dim == 1 ? "Weight" : "WeightTensor") + std::to_string(dim) + "d" + (is_all_tensor ? "" : "Flattened"))
616                                                      : "WeightNonTensor";
617 
618           code << tab << "CeedScalar r_q" << var_suffix << "[" << (is_all_tensor && (dim >= 3) ? Q_name : "1") << "];\n";
619           CeedCallBackend(CeedBasisGetData(basis, &basis_data));
620           data->W = basis_data->d_q_weight_1d;
621           code << tab << function_name << "<" << P_name << ", " << Q_name << ">(data, W, r_q" << var_suffix << ");\n";
622         }
623         break;
624       }
625       // LCOV_EXCL_START
626       case CEED_EVAL_DIV:
627       case CEED_EVAL_CURL:
628         break;  // TODO: Not implemented
629                 // LCOV_EXCL_STOP
630     }
631   } else {
632     switch (eval_mode) {
633       case CEED_EVAL_NONE:
634         code << tab << "CeedScalar *r_e" << var_suffix << " = r_q" << var_suffix << ";\n";
635         break;  // No action
636       case CEED_EVAL_INTERP:
637         code << tab << "CeedScalar *r_e" << var_suffix << " = r_e_scratch;\n";
638         if (is_at_points) {
639           std::string function_name = (dim == 1 ? "InterpTranspose" : "InterpTransposeTensor") + std::to_string(dim) + "d";
640 
641           code << tab << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ", OP_T_1D>(data, r_c" << var_suffix
642                << ", s_B" << var_suffix << ", r_e" << var_suffix << ");\n";
643         } else {
644           std::string function_name =
645               is_tensor ? ((dim == 1 ? "InterpTranspose" : "InterpTransposeTensor") + std::string(is_collocated ? "CollocatedNodes" : "") +
646                            std::to_string(dim) + "d" + (is_all_tensor ? "" : "Flattened"))
647                         : "InterpTransposeNonTensor";
648           std::string op_t_1d_name = (is_all_tensor || !is_tensor) ? "OP_T_1D" : (P_1d > Q_1d ? P_name : Q_name);
649 
650           code << tab << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ", " << op_t_1d_name << ">(data, r_q"
651                << var_suffix << ", s_B" << var_suffix << ", r_e" << var_suffix << ");\n";
652         }
653         break;
654       case CEED_EVAL_GRAD:
655         code << tab << "CeedScalar *r_e" << var_suffix << " = r_e_scratch;\n";
656         if (is_at_points) {
657           std::string function_name = (dim == 1 ? "InterpTranspose" : "InterpTransposeTensor") + std::to_string(dim) + "d";
658 
659           code << tab << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ", OP_T_1D>(data, r_c" << var_suffix
660                << ", s_B" << var_suffix << ", r_e" << var_suffix << ");\n";
661         } else if (use_3d_slices) {
662           std::string function_name = (dim == 1 ? "InterpTranspose" : "InterpTransposeTensor") + std::string(is_collocated ? "CollocatedNodes" : "") +
663                                       std::to_string(dim) + "d";
664 
665           code << tab << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ", OP_T_1D>(data, r_q" << var_suffix
666                << ", s_B" << var_suffix << ", r_e" << var_suffix << ");\n";
667         } else if (is_tensor) {
668           bool        is_collocated_grad = dim == 3 && Q_1d >= P_1d;
669           std::string function_name =
670               (dim == 1 ? "GradTranspose"
671                         : ("GradTransposeTensor" + std::string(is_collocated ? "CollocatedNodes" : (is_collocated_grad ? "Collocated" : "")))) +
672               std::to_string(dim) + "d" + (is_all_tensor ? "" : "Flattened");
673           std::string op_t_1d_name = is_all_tensor ? "OP_T_1D" : (P_1d > Q_1d ? P_name : Q_name);
674 
675           code << tab << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ", " << op_t_1d_name << ">(data, r_q"
676                << var_suffix << ", s_B" << var_suffix << ", s_G" << var_suffix << ", r_e" << var_suffix << ");\n";
677         } else {
678           std::string function_name = "GradTransposeNonTensor";
679 
680           code << tab << function_name << "<num_comp" << var_suffix << ", dim" << var_suffix << ", " << P_name << ", " << Q_name
681                << ", OP_T_1D>(data, r_q" << var_suffix << ", s_G" << var_suffix << ", r_e" << var_suffix << ");\n";
682         }
683         break;
684       // LCOV_EXCL_START
685       case CEED_EVAL_WEIGHT:
686         break;  // Should not occur
687       case CEED_EVAL_DIV:
688       case CEED_EVAL_CURL:
689         break;  // TODO: Not implemented
690                 // LCOV_EXCL_STOP
691     }
692   }
693   CeedCallBackend(CeedBasisDestroy(&basis));
694   return CEED_ERROR_SUCCESS;
695 }
696 
697 //------------------------------------------------------------------------------
698 // QFunction
699 //------------------------------------------------------------------------------
700 static int CeedOperatorBuildKernelQFunction_Cuda_gen(std::ostringstream &code, CeedOperator_Cuda_gen *data, Tab &tab, CeedInt max_dim,
701                                                      CeedInt max_num_points, CeedInt num_input_fields, CeedOperatorField *op_input_fields,
702                                                      CeedQFunctionField *qf_input_fields, CeedInt num_output_fields,
703                                                      CeedOperatorField *op_output_fields, CeedQFunctionField *qf_output_fields,
704                                                      std::string qfunction_name, CeedInt Q_1d, bool is_all_tensor, bool is_at_points,
705                                                      bool use_3d_slices) {
706   std::string         Q_name    = is_all_tensor ? "Q_1d" : "Q";
707   CeedEvalMode        eval_mode = CEED_EVAL_NONE;
708   CeedElemRestriction elem_rstr;
709 
710   // Setup output arrays
711   code << "\n";
712   code << tab << "// -- Output field setup\n";
713   for (CeedInt i = 0; i < num_output_fields; i++) {
714     const char *field_name;
715     std::string var_suffix = "_out_" + std::to_string(i);
716 
717     CeedCallBackend(CeedOperatorFieldGetName(op_output_fields[i], &field_name));
718     code << tab << "// ---- Output field " << i << ": " << field_name << "\n";
719     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
720     switch (eval_mode) {
721       case CEED_EVAL_NONE:
722         if (is_at_points) {
723           code << tab << "CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "];\n";
724         } else {
725           code << tab << "CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << (is_all_tensor && (max_dim >= 3) ? Q_name : "1")
726                << "];\n";
727         }
728         break;
729       case CEED_EVAL_INTERP:
730         if (is_at_points) {
731           // Accumulator for point data
732           code << tab << "CeedScalar r_c" << var_suffix << "[num_comp" << var_suffix << "*" << (max_dim >= 3 ? Q_name : "1") << "];\n";
733           code << tab << "for (CeedInt i = 0; i < num_comp" << var_suffix << "*" << (max_dim >= 3 ? Q_name : "1") << "; i++) r_c" << var_suffix
734                << "[i] = 0.0;\n";
735         } else {
736           code << tab << "CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << (is_all_tensor && (max_dim >= 3) ? Q_name : "1")
737                << "];\n";
738         }
739         break;
740       case CEED_EVAL_GRAD:
741         if (is_at_points) {
742           // Accumulator for point data
743           code << tab << "CeedScalar r_c" << var_suffix << "[num_comp" << var_suffix << "*" << (max_dim >= 3 ? Q_name : "1") << "];\n";
744           code << tab << "for (CeedInt i = 0; i < num_comp" << var_suffix << "*" << (max_dim >= 3 ? Q_name : "1") << "; i++) r_c" << var_suffix
745                << "[i] = 0.0;\n";
746         } else if (use_3d_slices) {
747           // Accumulator for gradient slices
748           code << tab << "CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << Q_name << "];\n";
749           code << tab << "for (CeedInt i = 0; i < num_comp" << var_suffix << "*" << Q_name << "; i++) r_q" << var_suffix << "[i] = 0.0;\n";
750         } else {
751           code << tab << "CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*dim" << var_suffix << "*"
752                << (is_all_tensor && (max_dim >= 3) ? Q_name : "1") << "];\n";
753         }
754         break;
755       case CEED_EVAL_WEIGHT:
756         break;
757         // LCOV_EXCL_START
758       case CEED_EVAL_DIV:
759       case CEED_EVAL_CURL:
760         break;  // TODO: Not implemented
761                 // LCOV_EXCL_STOP
762     }
763   }
764 
765   if (is_at_points) {
766     // We need to handle batches of points
767     code << "\n";
768     code << tab << "// Note: Using batches of points\n";
769     code << tab << "const CeedInt point_loop_bound = (blockDim.x*blockDim.y) * ceil((1.0*max_num_points) / (blockDim.x*blockDim.y));\n\n";
770     code << tab << "#pragma unroll\n";
771     code << tab << "for (CeedInt i = threadIdx.x + threadIdx.y*blockDim.x; i < point_loop_bound; i += blockDim.x*blockDim.y) {\n";
772     tab.push();
773     code << tab << "const CeedInt p = i % max_num_points;\n\n";
774 
775     code << tab << "// -- Coordinates\n";
776     code << tab << "CeedScalar r_x[max_dim];\n";
777     code << tab << "ReadPoint<max_dim, coords_comp_stride, max_num_points>(data, elem, p, max_num_points, points.indices, points.coords, r_x);\n\n";
778 
779     code << tab << "// -- Input fields\n";
780     for (CeedInt i = 0; i < num_input_fields; i++) {
781       const char *field_name;
782       std::string var_suffix = "_in_" + std::to_string(i);
783       std::string P_name     = "P_1d" + var_suffix;
784 
785       CeedCallBackend(CeedOperatorFieldGetName(op_input_fields[i], &field_name));
786       code << tab << "// ---- Input field " << i << ": " << field_name << "\n";
787       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
788       // Basis action
789       code << tab << "// EvalMode: " << CeedEvalModes[eval_mode] << "\n";
790       switch (eval_mode) {
791         case CEED_EVAL_NONE:
792           code << tab << "CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
793           code << tab << "ReadPoint<num_comp" << var_suffix << ", comp_stride" << var_suffix
794                << ", max_num_points>(data, elem, p, max_num_points, indices.inputs[" << i << "], d" << var_suffix << ", r_s" << var_suffix << ");\n";
795           break;
796         case CEED_EVAL_INTERP:
797           code << tab << "CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
798           code << tab << "InterpAtPoints" << max_dim << "d<num_comp" << var_suffix << ", max_num_points, " << P_name << ", " << Q_name
799                << ">(data, i, r_c" << var_suffix << ", r_x, r_s" << var_suffix << ");\n";
800           break;
801         case CEED_EVAL_GRAD:
802           code << tab << "CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim" << var_suffix << "];\n";
803           code << tab << "GradAtPoints" << max_dim << "d<num_comp" << var_suffix << ", max_num_points, " << P_name << ", " << Q_name
804                << ">(data, i, r_c" << var_suffix << ", r_x, r_s" << var_suffix << ");\n";
805           break;
806         case CEED_EVAL_WEIGHT:
807           code << tab << "CeedScalar r_s" << var_suffix << "[1];\n";
808           code << tab << "r_s" << var_suffix << "[0] = 1.0;\n";
809           break;
810           // LCOV_EXCL_START
811         case CEED_EVAL_DIV:
812         case CEED_EVAL_CURL:
813           break;  // TODO: Not implemented
814                   // LCOV_EXCL_STOP
815       }
816     }
817     code << "\n";
818     code << tab << "// -- Output fields\n";
819     for (CeedInt i = 0; i < num_output_fields; i++) {
820       const char *field_name;
821       std::string var_suffix = "_out_" + std::to_string(i);
822 
823       CeedCallBackend(CeedOperatorFieldGetName(op_output_fields[i], &field_name));
824       code << tab << "// ---- Output field " << i << ": " << field_name << "\n";
825       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
826       // Basis action
827       switch (eval_mode) {
828         case CEED_EVAL_NONE:
829           code << tab << "CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
830           break;
831         case CEED_EVAL_INTERP:
832           code << tab << "CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
833           break;
834         case CEED_EVAL_GRAD:
835           code << tab << "CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim" << var_suffix << "];\n";
836           break;
837           // LCOV_EXCL_START
838         case CEED_EVAL_WEIGHT:
839           break;  // Should not occur
840         case CEED_EVAL_DIV:
841         case CEED_EVAL_CURL:
842           break;  // TODO: Not implemented
843                   // LCOV_EXCL_STOP
844       }
845     }
846 
847   } else if (use_3d_slices) {
848     // We treat quadrature points per slice in 3d to save registers
849     code << "\n";
850     code << tab << "// Note: Using planes of 3D elements\n";
851     code << tab << "#pragma unroll\n";
852     code << tab << "for (CeedInt q = 0; q < " << Q_name << "; q++) {\n";
853     tab.push();
854     code << tab << "// -- Input fields\n";
855     for (CeedInt i = 0; i < num_input_fields; i++) {
856       const char *field_name;
857       std::string var_suffix = "_in_" + std::to_string(i);
858 
859       CeedCallBackend(CeedOperatorFieldGetName(op_input_fields[i], &field_name));
860       code << tab << "// ---- Input field " << i << ": " << field_name << "\n";
861       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
862       // Basis action
863       code << tab << "// EvalMode: " << CeedEvalModes[eval_mode] << "\n";
864       switch (eval_mode) {
865         case CEED_EVAL_NONE:
866           bool is_strided;
867 
868           code << tab << "CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
869 
870           CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr));
871           CeedCallBackend(CeedElemRestrictionIsStrided(elem_rstr, &is_strided));
872           if (is_strided) {
873             bool    has_backend_strides;
874             CeedInt num_elem, elem_size;
875 
876             CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
877             CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides));
878             CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem));
879             CeedInt strides[3] = {1, elem_size * num_elem, elem_size};
880 
881             if (!has_backend_strides) {
882               CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides));
883             }
884             code << tab << "const CeedInt strides" << var_suffix << "_0 = " << strides[0] << ", strides" << var_suffix << "_1 = " << strides[1]
885                  << ", strides" << var_suffix << "_2 = " << strides[2] << ";\n";
886             code << tab << "ReadEVecSliceStrided3d<num_comp" << var_suffix << ", " << Q_name << ", strides" << var_suffix << "_0, strides"
887                  << var_suffix << "_1, strides" << var_suffix << "_2>(data, elem, q, d" << var_suffix << ", r_s" << var_suffix << ");\n";
888           } else {
889             CeedSize                  l_size = 0;
890             CeedInt                   comp_stride;
891             CeedElemRestriction_Cuda *rstr_data;
892 
893             CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size));
894             code << tab << "const CeedInt l_size" << var_suffix << " = " << l_size << ";\n";
895             CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
896             code << tab << "const CeedInt comp_stride" << var_suffix << " = " << comp_stride << ";\n";
897             CeedCallBackend(CeedElemRestrictionGetData(elem_rstr, &rstr_data));
898             data->indices.inputs[i] = (CeedInt *)rstr_data->d_offsets;
899             code << tab << "ReadEVecSliceStandard3d<num_comp" << var_suffix << ", comp_stride" << var_suffix << ", " << Q_name << ">(data, l_size"
900                  << var_suffix << ", elem, q, indices.inputs[" << i << "], d" << var_suffix << ", r_s" << var_suffix << ");\n";
901           }
902           CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
903           break;
904         case CEED_EVAL_INTERP:
905           code << tab << "CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
906           code << tab << "for (CeedInt j = 0; j < num_comp" << var_suffix << "; j++) {\n";
907           tab.push();
908           code << tab << "r_s" << var_suffix << "[j] = r_q" << var_suffix << "[q + j*" << Q_name << "];\n";
909           tab.pop();
910           code << tab << "}\n";
911           break;
912         case CEED_EVAL_GRAD:
913           code << tab << "CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim" << var_suffix << "];\n";
914           code << tab << "GradColloSlice3d<num_comp" << var_suffix << ", " << Q_name << ", OP_T_1D>(data, q, r_q" << var_suffix << ", s_G"
915                << var_suffix << ", r_s" << var_suffix << ");\n";
916           break;
917         case CEED_EVAL_WEIGHT:
918           code << tab << "CeedScalar r_s" << var_suffix << "[1];\n";
919           code << tab << "r_s" << var_suffix << "[0] = r_q" << var_suffix << "[q];\n";
920           break;
921           // LCOV_EXCL_START
922         case CEED_EVAL_DIV:
923         case CEED_EVAL_CURL:
924           break;  // TODO: Not implemented
925                   // LCOV_EXCL_STOP
926       }
927     }
928     code << "\n";
929     code << tab << "// -- Output fields\n";
930     for (CeedInt i = 0; i < num_output_fields; i++) {
931       const char *field_name;
932       std::string var_suffix = "_out_" + std::to_string(i);
933 
934       CeedCallBackend(CeedOperatorFieldGetName(op_output_fields[i], &field_name));
935       code << tab << "// ---- Output field " << i << ": " << field_name << "\n";
936       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
937       // Basis action
938       switch (eval_mode) {
939         case CEED_EVAL_NONE:
940           code << tab << "CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
941           break;
942         case CEED_EVAL_INTERP:
943           code << tab << "CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
944           break;
945         case CEED_EVAL_GRAD:
946           code << tab << "CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim" << var_suffix << "];\n";
947           break;
948           // LCOV_EXCL_START
949         case CEED_EVAL_WEIGHT:
950           break;  // Should not occur
951         case CEED_EVAL_DIV:
952         case CEED_EVAL_CURL:
953           break;  // TODO: Not implemented
954                   // LCOV_EXCL_STOP
955       }
956     }
957   } else {
958     code << "\n";
959     code << tab << "// Note: Using full elements\n";
960     code << tab << "{\n";
961     tab.push();
962     code << tab << "// -- Input fields\n";
963     for (CeedInt i = 0; i < num_input_fields; i++) {
964       const char *field_name;
965 
966       CeedCallBackend(CeedOperatorFieldGetName(op_input_fields[i], &field_name));
967       code << tab << "// ---- Input field " << i << ": " << field_name << "\n";
968       code << tab << "CeedScalar *r_s_in_" << i << " = r_q_in_" << i << ";\n";
969     }
970     code << tab << "// -- Output fields\n";
971     for (CeedInt i = 0; i < num_output_fields; i++) {
972       const char *field_name;
973 
974       CeedCallBackend(CeedOperatorFieldGetName(op_output_fields[i], &field_name));
975       code << tab << "// ---- Output field " << i << ": " << field_name << "\n";
976       code << tab << "CeedScalar *r_s_out_" << i << " = r_q_out_" << i << ";\n";
977     }
978   }
979 
980   // Input and output buffers
981   code << "\n";
982   code << tab << "// -- QFunction inputs and outputs\n";
983   code << tab << "// ---- Inputs\n";
984   code << tab << "CeedScalar *inputs[" << CeedIntMax(num_input_fields, 1) << "];\n";
985   for (CeedInt i = 0; i < num_input_fields; i++) {
986     const char *field_name;
987 
988     CeedCallBackend(CeedOperatorFieldGetName(op_input_fields[i], &field_name));
989     code << tab << "// ------ Input field " << i << ": " << field_name << "\n";
990     code << tab << "inputs[" << i << "] = r_s_in_" << i << ";\n";
991   }
992   code << tab << "// ---- Outputs\n";
993   code << tab << "CeedScalar *outputs[" << CeedIntMax(num_output_fields, 1) << "];\n";
994   for (CeedInt i = 0; i < num_output_fields; i++) {
995     const char *field_name;
996 
997     CeedCallBackend(CeedOperatorFieldGetName(op_output_fields[i], &field_name));
998     code << tab << "// ------ Output field " << i << ": " << field_name << "\n";
999     code << tab << "outputs[" << i << "] = r_s_out_" << i << ";\n";
1000   }
1001 
1002   // Apply QFunction
1003   code << "\n";
1004   code << tab << "// -- Apply QFunction\n";
1005   code << tab << "" << qfunction_name << "(ctx, ";
1006   if (max_dim != 3 || is_at_points || use_3d_slices || !is_all_tensor) {
1007     code << "1";
1008   } else {
1009     code << Q_name;
1010   }
1011   code << ", inputs, outputs);\n";
1012 
1013   if (is_at_points) {
1014     // Map back to coefficients
1015     code << "\n";
1016     code << tab << "// -- Output fields\n";
1017     for (CeedInt i = 0; i < num_output_fields; i++) {
1018       const char *field_name;
1019       std::string var_suffix = "_out_" + std::to_string(i);
1020       std::string P_name     = "P_1d" + var_suffix;
1021 
1022       CeedCallBackend(CeedOperatorFieldGetName(op_output_fields[i], &field_name));
1023       code << tab << "// ---- Output field " << i << ": " << field_name << "\n";
1024       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
1025       // Basis action
1026       code << tab << "// EvalMode: " << CeedEvalModes[eval_mode] << "\n";
1027       switch (eval_mode) {
1028         case CEED_EVAL_NONE: {
1029           CeedInt             comp_stride;
1030           CeedElemRestriction elem_rstr;
1031 
1032           CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &elem_rstr));
1033           CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
1034           CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
1035           code << tab << "const CeedInt comp_stride" << var_suffix << " = " << comp_stride << ";\n";
1036           code << tab << "WritePoint<num_comp" << var_suffix << ", comp_stride" << var_suffix
1037                << ", max_num_points>(data, elem, i, points.num_per_elem[elem], indices.outputs[" << i << "]"
1038                << ", r_s" << var_suffix << ", d" << var_suffix << ");\n";
1039           break;
1040         }
1041         case CEED_EVAL_INTERP:
1042           code << tab << "if (i >= points.num_per_elem[elem]) {\n";
1043           tab.push();
1044           code << tab << "for (CeedInt j = 0; j < num_comp" << var_suffix << "; j++) r_s" << var_suffix << "[j] = 0.0;\n";
1045           tab.pop();
1046           code << tab << "}\n";
1047           code << tab << "InterpTransposeAtPoints" << max_dim << "d<num_comp" << var_suffix << ", max_num_points, " << P_name << ", " << Q_name
1048                << ">(data, i, r_s" << var_suffix << ", r_x, r_c" << var_suffix << ");\n";
1049           break;
1050         case CEED_EVAL_GRAD:
1051           code << tab << "if (i >= points.num_per_elem[elem]) {\n";
1052           tab.push();
1053           code << tab << "for (CeedInt j = 0; j < num_comp" << var_suffix << "*dim" << var_suffix << "; j++) r_s" << var_suffix << "[j] = 0.0;\n";
1054           tab.pop();
1055           code << tab << "}\n";
1056           code << tab << "GradTransposeAtPoints" << max_dim << "d<num_comp" << var_suffix << ", max_num_points, " << P_name << ", " << Q_name
1057                << ">(data, i, r_s" << var_suffix << ", r_x, r_c" << var_suffix << ");\n";
1058           break;
1059           // LCOV_EXCL_START
1060         case CEED_EVAL_WEIGHT:
1061           break;  // Should not occur
1062         case CEED_EVAL_DIV:
1063         case CEED_EVAL_CURL:
1064           break;  // TODO: Not implemented
1065                   // LCOV_EXCL_STOP
1066       }
1067     }
1068   } else if (use_3d_slices) {
1069     // Copy or apply transpose grad, if needed
1070     code << "\n";
1071     code << tab << "// -- Output fields\n";
1072     for (CeedInt i = 0; i < num_output_fields; i++) {
1073       const char *field_name;
1074       std::string var_suffix = "_out_" + std::to_string(i);
1075       std::string P_name     = "P_1d" + var_suffix;
1076 
1077       CeedCallBackend(CeedOperatorFieldGetName(op_output_fields[i], &field_name));
1078       code << tab << "// ---- Output field " << i << ": " << field_name << "\n";
1079       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
1080       // Basis action
1081       code << tab << "// EvalMode: " << CeedEvalModes[eval_mode] << "\n";
1082       switch (eval_mode) {
1083         case CEED_EVAL_NONE:
1084           code << tab << "for (CeedInt j = 0; j < num_comp" << var_suffix << " ; j++) {\n";
1085           tab.push();
1086           code << tab << "r_q" << var_suffix << "[q + j*" << Q_name << "] = r_s" << var_suffix << "[j];\n";
1087           tab.pop();
1088           code << tab << "}\n";
1089           break;
1090         case CEED_EVAL_INTERP:
1091           code << tab << "for (CeedInt j = 0; j < num_comp" << var_suffix << " ; j++) {\n";
1092           tab.push();
1093           code << tab << "r_q" << var_suffix << "[q + j*" << Q_name << "] = r_s" << var_suffix << "[j];\n";
1094           tab.pop();
1095           code << tab << "}\n";
1096           break;
1097         case CEED_EVAL_GRAD:
1098           code << tab << "GradColloSliceTranspose3d<num_comp" << var_suffix << ", " << Q_name << ", OP_T_1D>(data, q, r_s" << var_suffix << ", s_G"
1099                << var_suffix << ", r_q" << var_suffix << ");\n";
1100           break;
1101           // LCOV_EXCL_START
1102         case CEED_EVAL_WEIGHT:
1103           break;  // Should not occur
1104         case CEED_EVAL_DIV:
1105         case CEED_EVAL_CURL:
1106           break;  // TODO: Not implemented
1107                   // LCOV_EXCL_STOP
1108       }
1109     }
1110   }
1111   tab.pop();
1112   code << tab << "}\n";
1113   return CEED_ERROR_SUCCESS;
1114 }
1115 
1116 //------------------------------------------------------------------------------
1117 // Build single operator kernel
1118 //------------------------------------------------------------------------------
1119 extern "C" int CeedOperatorBuildKernel_Cuda_gen(CeedOperator op, bool *is_good_build) {
1120   bool                    is_all_tensor = true, is_all_nontensor = true, is_at_points = false, use_3d_slices = false;
1121   Ceed                    ceed;
1122   CeedInt                 Q = 0, Q_1d = 0, num_input_fields, num_output_fields, max_dim = 1, max_num_points = 0, coords_comp_stride = 0;
1123   CeedQFunctionField     *qf_input_fields, *qf_output_fields;
1124   CeedQFunction_Cuda_gen *qf_data;
1125   CeedQFunction           qf;
1126   CeedOperatorField      *op_input_fields, *op_output_fields;
1127   CeedOperator_Cuda_gen  *data;
1128   std::ostringstream      code;
1129   Tab                     tab;
1130 
1131   CeedCallBackend(CeedOperatorGetData(op, &data));
1132   {
1133     bool is_setup_done;
1134 
1135     CeedCallBackend(CeedOperatorIsSetupDone(op, &is_setup_done));
1136     if (is_setup_done) {
1137       *is_good_build = !data->use_fallback;
1138       return CEED_ERROR_SUCCESS;
1139     }
1140   }
1141 
1142   // Check field compatibility
1143   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
1144   {
1145     bool has_shared_bases = true;
1146 
1147     for (CeedInt i = 0; i < num_input_fields; i++) {
1148       CeedBasis basis;
1149 
1150       CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
1151       if (basis != CEED_BASIS_NONE) {
1152         bool        is_tensor = true;
1153         const char *resource;
1154         char       *resource_root;
1155         Ceed        basis_ceed;
1156 
1157         CeedCallBackend(CeedBasisIsTensor(basis, &is_tensor));
1158         is_all_tensor    = is_all_tensor && is_tensor;
1159         is_all_nontensor = is_all_nontensor && !is_tensor;
1160         CeedCallBackend(CeedBasisGetCeed(basis, &basis_ceed));
1161         CeedCallBackend(CeedGetResource(basis_ceed, &resource));
1162         CeedCallBackend(CeedGetResourceRoot(basis_ceed, resource, ":", &resource_root));
1163         has_shared_bases = has_shared_bases && !strcmp(resource_root, "/gpu/cuda/shared");
1164         CeedCallBackend(CeedFree(&resource_root));
1165         CeedCallBackend(CeedDestroy(&basis_ceed));
1166       }
1167       CeedCallBackend(CeedBasisDestroy(&basis));
1168     }
1169 
1170     for (CeedInt i = 0; i < num_output_fields; i++) {
1171       CeedBasis basis;
1172 
1173       CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
1174       if (basis != CEED_BASIS_NONE) {
1175         bool        is_tensor = true;
1176         const char *resource;
1177         char       *resource_root;
1178         Ceed        basis_ceed;
1179 
1180         CeedCallBackend(CeedBasisIsTensor(basis, &is_tensor));
1181         is_all_tensor    = is_all_tensor && is_tensor;
1182         is_all_nontensor = is_all_nontensor && !is_tensor;
1183 
1184         CeedCallBackend(CeedBasisGetCeed(basis, &basis_ceed));
1185         CeedCallBackend(CeedGetResource(basis_ceed, &resource));
1186         CeedCallBackend(CeedGetResourceRoot(basis_ceed, resource, ":", &resource_root));
1187         has_shared_bases = has_shared_bases && !strcmp(resource_root, "/gpu/cuda/shared");
1188         CeedCallBackend(CeedFree(&resource_root));
1189         CeedCallBackend(CeedDestroy(&basis_ceed));
1190       }
1191       CeedCallBackend(CeedBasisDestroy(&basis));
1192     }
1193     // -- Fallback to ref if not all bases are shared
1194     if (!has_shared_bases) {
1195       *is_good_build = false;
1196       return CEED_ERROR_SUCCESS;
1197     }
1198   }
1199   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
1200   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
1201   CeedCallBackend(CeedQFunctionGetData(qf, &qf_data));
1202   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
1203 
1204   // Get operator data
1205   CeedCallBackend(CeedOperatorIsAtPoints(op, &is_at_points));
1206   {
1207     CeedInt max_P = 0, max_P_1d = 0;
1208 
1209     CeedCallBackend(CeedOperatorBuildKernelData_Cuda_gen(ceed, num_input_fields, op_input_fields, qf_input_fields, num_output_fields,
1210                                                          op_output_fields, qf_output_fields, &max_P, &max_P_1d, &Q, &Q_1d, &max_dim, &is_all_tensor,
1211                                                          &use_3d_slices));
1212     data->max_P_1d = is_all_tensor ? max_P_1d : max_P;
1213   }
1214   if (max_dim == 0) max_dim = 1;
1215   data->dim = max_dim;
1216   if (is_at_points) {
1217     CeedElemRestriction_Cuda *rstr_data;
1218     CeedElemRestriction       rstr_points = NULL;
1219 
1220     CeedCallBackend(CeedOperatorAtPointsGetPoints(op, &rstr_points, NULL));
1221     CeedCallBackend(CeedElemRestrictionGetMaxPointsInElement(rstr_points, &max_num_points));
1222     CeedCallBackend(CeedElemRestrictionGetCompStride(rstr_points, &coords_comp_stride));
1223     CeedCallBackend(CeedElemRestrictionGetData(rstr_points, &rstr_data));
1224     data->points.indices = (CeedInt *)rstr_data->d_offsets;
1225     CeedCallBackend(CeedElemRestrictionDestroy(&rstr_points));
1226   }
1227   if (is_at_points) use_3d_slices = false;
1228   if (Q_1d == 0) {
1229     if (is_at_points) Q_1d = max_num_points;
1230     else CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q_1d));
1231   }
1232   if (Q == 0) Q = Q_1d;
1233   data->Q    = Q;
1234   data->Q_1d = Q_1d;
1235 
1236   // Check for restriction only identity operator
1237   {
1238     bool is_identity_qf;
1239 
1240     CeedCallBackend(CeedQFunctionIsIdentity(qf, &is_identity_qf));
1241     if (is_identity_qf) {
1242       CeedEvalMode eval_mode_in, eval_mode_out;
1243 
1244       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[0], &eval_mode_in));
1245       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[0], &eval_mode_out));
1246       CeedCheck(eval_mode_in != CEED_EVAL_NONE || eval_mode_out != CEED_EVAL_NONE, ceed, CEED_ERROR_BACKEND,
1247                 "Backend does not implement restriction only identity operators");
1248     }
1249   }
1250 
1251   // Add atomicAdd function for old NVidia architectures
1252   {
1253     Ceed_Cuda            *ceed_data;
1254     struct cudaDeviceProp prop;
1255 
1256     CeedCallBackend(CeedGetData(ceed, &ceed_data));
1257     CeedCallBackend(cudaGetDeviceProperties(&prop, ceed_data->device_id));
1258     if ((prop.major < 6) && (CEED_SCALAR_TYPE != CEED_SCALAR_FP32)) {
1259       code << tab << "// AtomicAdd fallback source\n";
1260       code << tab << "#include <ceed/jit-source/cuda/cuda-atomic-add-fallback.h>\n\n";
1261     }
1262   }
1263 
1264   // Load basis source files
1265   if (!is_all_nontensor) {
1266     code << tab << "// Tensor basis source\n";
1267     code << tab << "#include <ceed/jit-source/cuda/cuda-shared-basis-tensor-templates.h>\n\n";
1268   }
1269   if (!is_all_tensor) {
1270     code << tab << "// Non-tensor basis source\n";
1271     code << tab << "#include <ceed/jit-source/cuda/cuda-shared-basis-nontensor-templates.h>\n\n";
1272   }
1273   if (!is_all_tensor && !is_all_nontensor) {
1274     code << "// Tensor basis source\n";
1275     code << "#include <ceed/jit-source/cuda/cuda-shared-basis-tensor-flattened-templates.h>\n\n";
1276   }
1277   if (is_at_points) {
1278     code << "// AtPoints basis source\n";
1279     code << "#include <ceed/jit-source/cuda/cuda-shared-basis-tensor-at-points-templates.h>\n\n";
1280   }
1281   code << "// CodeGen operator source\n";
1282   code << "#include <ceed/jit-source/cuda/cuda-gen-templates.h>\n\n";
1283 
1284   // Get QFunction name
1285   std::string qfunction_name(qf_data->qfunction_name);
1286   std::string operator_name;
1287 
1288   operator_name = "CeedKernelCudaGenOperator_" + qfunction_name;
1289 
1290   // Define CEED_Q_VLA
1291   code << "\n" << tab << "#undef CEED_Q_VLA\n";
1292   if (max_dim != 3 || is_at_points || use_3d_slices || !is_all_tensor) {
1293     code << tab << "#define CEED_Q_VLA 1\n\n";
1294   } else {
1295     code << tab << "#define CEED_Q_VLA " << Q_1d << "\n\n";
1296   }
1297 
1298   // Add user QFunction source
1299   {
1300     const char *source_path;
1301 
1302     CeedCallBackend(CeedQFunctionGetSourcePath(qf, &source_path));
1303     CeedCheck(source_path, ceed, CEED_ERROR_UNSUPPORTED, "/gpu/cuda/gen backend requires QFunction source code file");
1304 
1305     code << tab << "// User QFunction source\n";
1306     code << tab << "#include \"" << source_path << "\"\n\n";
1307   }
1308 
1309   // Setup
1310   code << "\n" << tab << "// -----------------------------------------------------------------------------\n";
1311   code << tab << "// Operator Kernel\n";
1312   code << tab << "// \n";
1313   code << tab << "// d_[in,out]_i:   CeedVector device array\n";
1314   code << tab << "// r_[in,out]_e_i: Element vector register\n";
1315   code << tab << "// r_[in,out]_q_i: Quadrature space vector register\n";
1316   code << tab << "// r_[in,out]_c_i: AtPoints Chebyshev coefficients register\n";
1317   code << tab << "// r_[in,out]_s_i: Quadrature space slice vector register\n";
1318   code << tab << "// \n";
1319   code << tab << "// s_B_[in,out]_i: Interpolation matrix, shared memory\n";
1320   code << tab << "// s_G_[in,out]_i: Gradient matrix, shared memory\n";
1321   code << tab << "// -----------------------------------------------------------------------------\n";
1322   code << tab << "extern \"C\" __global__ void " << operator_name
1323        << "(CeedInt num_elem, void* ctx, FieldsInt_Cuda indices, Fields_Cuda fields, Fields_Cuda B, Fields_Cuda G, CeedScalar *W, Points_Cuda "
1324           "points) {\n";
1325   tab.push();
1326 
1327   // Scratch buffers
1328   for (CeedInt i = 0; i < num_input_fields; i++) {
1329     CeedEvalMode eval_mode;
1330 
1331     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
1332     if (eval_mode != CEED_EVAL_WEIGHT) {  // Skip CEED_EVAL_WEIGHT
1333       code << tab << "const CeedScalar *__restrict__ d_in_" << i << " = fields.inputs[" << i << "];\n";
1334     }
1335   }
1336   for (CeedInt i = 0; i < num_output_fields; i++) {
1337     code << tab << "CeedScalar *__restrict__ d_out_" << i << " = fields.outputs[" << i << "];\n";
1338   }
1339 
1340   code << tab << "const CeedInt max_dim = " << max_dim << ";\n";
1341   if (!is_all_tensor) {
1342     code << tab << "const CeedInt Q = " << Q << ";\n";
1343   }
1344   if (!is_all_nontensor) {
1345     code << tab << "const CeedInt Q_1d = " << Q_1d << ";\n";
1346   }
1347   if (is_at_points) {
1348     code << tab << "const CeedInt max_num_points = " << max_num_points << ";\n";
1349     code << tab << "const CeedInt coords_comp_stride = " << coords_comp_stride << ";\n";
1350   }
1351 
1352   // Shared data
1353   code << tab << "extern __shared__ CeedScalar slice[];\n";
1354   code << tab << "SharedData_Cuda data;\n";
1355   code << tab << "data.t_id_x = threadIdx.x;\n";
1356   code << tab << "data.t_id_y = threadIdx.y;\n";
1357   code << tab << "data.t_id_z = threadIdx.z;\n";
1358   code << tab << "data.t_id   = threadIdx.x + threadIdx.y*blockDim.x + threadIdx.z*blockDim.y*blockDim.x;\n";
1359   code << tab << "data.slice  = slice + data.t_id_z*OP_T_1D" << ((!is_all_tensor || max_dim == 1) ? "" : "*OP_T_1D") << ";\n";
1360 
1361   // -- Determine input mat reuse
1362   FieldReuse_Cuda input_matrix_reuse[CEED_FIELD_MAX];
1363 
1364   for (CeedInt i = 0; i < num_input_fields; i++) {
1365     input_matrix_reuse[i].index = -1;
1366   }
1367   for (CeedInt i = 0; i < num_input_fields; i++) {
1368     bool         is_tensor = true;
1369     CeedEvalMode eval_mode_i;
1370     CeedBasis    basis_i;
1371 
1372     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode_i));
1373     if (eval_mode_i == CEED_EVAL_WEIGHT) continue;
1374     CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis_i));
1375     CeedCallBackend(CeedBasisIsTensor(basis_i, &is_tensor));
1376     for (CeedInt j = 0; (input_matrix_reuse[i].index == -1) && (j < i); j++) {
1377       CeedEvalMode eval_mode_j;
1378       CeedBasis    basis_j;
1379 
1380       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[j], &eval_mode_j));
1381       if (eval_mode_j == CEED_EVAL_WEIGHT) continue;
1382       CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[j], &basis_j));
1383       if (basis_i == basis_j) {
1384         if (is_tensor) {
1385           input_matrix_reuse[i].index     = j;
1386           input_matrix_reuse[i].is_input  = true;
1387           input_matrix_reuse[i].eval_mode = eval_mode_j;
1388         } else {
1389           // For non-tensor can only re-use with the same eval mode
1390           if (eval_mode_i == eval_mode_j) {
1391             input_matrix_reuse[i].index     = j;
1392             input_matrix_reuse[i].is_input  = true;
1393             input_matrix_reuse[i].eval_mode = eval_mode_j;
1394           }
1395         }
1396       }
1397       CeedCallBackend(CeedBasisDestroy(&basis_j));
1398     }
1399     CeedCallBackend(CeedBasisDestroy(&basis_i));
1400   }
1401 
1402   // -- Determine output mat reuse
1403   FieldReuse_Cuda output_matrix_reuse[CEED_FIELD_MAX];
1404 
1405   for (CeedInt i = 0; i < num_output_fields; i++) {
1406     output_matrix_reuse[i].index = -1;
1407   }
1408   for (CeedInt i = 0; i < num_output_fields; i++) {
1409     bool         is_tensor = true;
1410     CeedEvalMode eval_mode_i;
1411     CeedBasis    basis_i;
1412 
1413     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode_i));
1414     CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis_i));
1415     CeedCallBackend(CeedBasisIsTensor(basis_i, &is_tensor));
1416     for (CeedInt j = 0; (output_matrix_reuse[i].index == -1) && (j < num_input_fields); j++) {
1417       CeedEvalMode eval_mode_j;
1418       CeedBasis    basis_j;
1419 
1420       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[j], &eval_mode_j));
1421       if (eval_mode_j == CEED_EVAL_WEIGHT) continue;
1422       CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[j], &basis_j));
1423       if (basis_i == basis_j) {
1424         if (is_tensor) {
1425           output_matrix_reuse[i].index     = j;
1426           output_matrix_reuse[i].is_input  = true;
1427           output_matrix_reuse[i].eval_mode = eval_mode_j;
1428         } else {
1429           // For non-tensor can only re-use with the same eval mode
1430           if (eval_mode_i == eval_mode_j) {
1431             output_matrix_reuse[i].index     = j;
1432             output_matrix_reuse[i].is_input  = true;
1433             output_matrix_reuse[i].eval_mode = eval_mode_j;
1434           }
1435         }
1436       }
1437       CeedCallBackend(CeedBasisDestroy(&basis_j));
1438     }
1439     for (CeedInt j = 0; (output_matrix_reuse[i].index == -1) && (j < i); j++) {
1440       CeedEvalMode eval_mode_j;
1441       CeedBasis    basis_j;
1442 
1443       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[j], &eval_mode_j));
1444       if (eval_mode_j == CEED_EVAL_WEIGHT) continue;
1445       CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[j], &basis_j));
1446       if (basis_i == basis_j) {
1447         if (is_tensor) {
1448           output_matrix_reuse[i].index     = j;
1449           output_matrix_reuse[i].is_input  = false;
1450           output_matrix_reuse[i].eval_mode = eval_mode_j;
1451         } else {
1452           // For non-tensor can only re-use with the same eval mode
1453           if (eval_mode_i == eval_mode_j) {
1454             output_matrix_reuse[i].index     = j;
1455             output_matrix_reuse[i].is_input  = false;
1456             output_matrix_reuse[i].eval_mode = eval_mode_j;
1457           }
1458         }
1459       }
1460       CeedCallBackend(CeedBasisDestroy(&basis_j));
1461     }
1462     CeedCallBackend(CeedBasisDestroy(&basis_i));
1463   }
1464 
1465   // Initialize constants, and matrices B and G
1466   code << "\n" << tab << "// Input field constants and basis data\n";
1467   for (CeedInt i = 0; i < num_input_fields; i++) {
1468     CeedCallBackend(CeedOperatorBuildKernelFieldData_Cuda_gen(code, data, tab, i, op_input_fields[i], qf_input_fields[i], input_matrix_reuse[i],
1469                                                               max_dim, Q, Q_1d, true, is_all_tensor, is_at_points, use_3d_slices, false));
1470   }
1471   code << "\n" << tab << "// Output field constants and basis data\n";
1472   for (CeedInt i = 0; i < num_output_fields; i++) {
1473     CeedCallBackend(CeedOperatorBuildKernelFieldData_Cuda_gen(code, data, tab, i, op_output_fields[i], qf_output_fields[i], output_matrix_reuse[i],
1474                                                               max_dim, Q, Q_1d, false, is_all_tensor, is_at_points, use_3d_slices, false));
1475   }
1476 
1477   // Loop over all elements
1478   code << "\n" << tab << "// Element loop\n";
1479   code << tab << "__syncthreads();\n";
1480   code << tab << "for (CeedInt elem = blockIdx.x*blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x*blockDim.z) {\n";
1481   tab.push();
1482 
1483   // -- Compute minimum buffer space needed
1484   CeedInt max_rstr_buffer_size = 1;
1485 
1486   for (CeedInt i = 0; i < num_input_fields; i++) {
1487     CeedEvalMode eval_mode;
1488 
1489     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
1490     if (eval_mode != CEED_EVAL_NONE && eval_mode != CEED_EVAL_WEIGHT) {
1491       CeedInt             num_comp;
1492       CeedElemRestriction elem_rstr;
1493 
1494       CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr));
1495       CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
1496       max_rstr_buffer_size = CeedIntMax(max_rstr_buffer_size, num_comp * (is_all_tensor && (max_dim >= 3) ? Q_1d : 1));
1497       CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
1498     }
1499   }
1500   for (CeedInt i = 0; i < num_output_fields; i++) {
1501     CeedEvalMode eval_mode;
1502 
1503     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
1504     if (eval_mode != CEED_EVAL_NONE) {
1505       CeedInt             num_comp;
1506       CeedElemRestriction elem_rstr;
1507 
1508       CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &elem_rstr));
1509       CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
1510       max_rstr_buffer_size = CeedIntMax(max_rstr_buffer_size, num_comp * (is_all_tensor && (max_dim >= 3) ? Q_1d : 1));
1511       CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
1512     }
1513   }
1514   code << tab << "// Scratch restriction buffer space\n";
1515   code << tab << "CeedScalar r_e_scratch[" << max_rstr_buffer_size << "];\n";
1516 
1517   // -- Determine best input field processing order
1518   CeedInt field_rstr_in_buffer[CEED_FIELD_MAX], input_field_order[CEED_FIELD_MAX];
1519 
1520   for (CeedInt i = 0; i < num_input_fields; i++) {
1521     field_rstr_in_buffer[i] = -1;
1522     input_field_order[i]    = -1;
1523   }
1524   {
1525     bool    is_ordered[CEED_FIELD_MAX];
1526     CeedInt curr_index = 0;
1527 
1528     for (CeedInt i = 0; i < num_input_fields; i++) is_ordered[i] = false;
1529     for (CeedInt i = 0; i < num_input_fields; i++) {
1530       CeedVector          vec_i;
1531       CeedElemRestriction rstr_i;
1532 
1533       if (is_ordered[i]) continue;
1534       field_rstr_in_buffer[i]       = i;
1535       is_ordered[i]                 = true;
1536       input_field_order[curr_index] = i;
1537       curr_index++;
1538       CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec_i));
1539       if (vec_i == CEED_VECTOR_NONE) continue;  // CEED_EVAL_WEIGHT
1540       CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &rstr_i));
1541       for (CeedInt j = i + 1; j < num_input_fields; j++) {
1542         CeedVector          vec_j;
1543         CeedElemRestriction rstr_j;
1544 
1545         CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[j], &vec_j));
1546         CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[j], &rstr_j));
1547         if (rstr_i == rstr_j && vec_i == vec_j) {
1548           field_rstr_in_buffer[j]       = i;
1549           is_ordered[j]                 = true;
1550           input_field_order[curr_index] = j;
1551           curr_index++;
1552         }
1553         CeedCallBackend(CeedVectorDestroy(&vec_j));
1554         CeedCallBackend(CeedElemRestrictionDestroy(&rstr_j));
1555       }
1556       CeedCallBackend(CeedVectorDestroy(&vec_i));
1557       CeedCallBackend(CeedElemRestrictionDestroy(&rstr_i));
1558     }
1559   }
1560 
1561   // -- Input restriction and basis
1562   code << "\n" << tab << "// -- Input field restrictions and basis actions\n";
1563   for (CeedInt i = 0; i < num_input_fields; i++) {
1564     const char   *field_name;
1565     const CeedInt f = input_field_order[i];
1566 
1567     CeedCallBackend(CeedOperatorFieldGetName(op_input_fields[f], &field_name));
1568     code << tab << "// ---- Input field " << f << ": " << field_name << "\n";
1569 
1570     // ---- Restriction
1571     CeedCallBackend(CeedOperatorBuildKernelRestriction_Cuda_gen(code, data, tab, f, field_rstr_in_buffer, op_input_fields[f], qf_input_fields[f],
1572                                                                 max_dim, Q_1d, true, is_all_tensor, is_at_points, use_3d_slices));
1573 
1574     // ---- Basis action
1575     CeedCallBackend(CeedOperatorBuildKernelBasis_Cuda_gen(code, data, tab, f, op_input_fields[f], qf_input_fields[f], max_dim, Q_1d, true,
1576                                                           is_all_tensor, is_at_points, use_3d_slices));
1577   }
1578 
1579   // -- Q function
1580   CeedCallBackend(CeedOperatorBuildKernelQFunction_Cuda_gen(code, data, tab, max_dim, max_num_points, num_input_fields, op_input_fields,
1581                                                             qf_input_fields, num_output_fields, op_output_fields, qf_output_fields, qfunction_name,
1582                                                             Q_1d, is_all_tensor, is_at_points, use_3d_slices));
1583 
1584   // -- Output basis and restriction
1585   code << "\n" << tab << "// -- Output field basis action and restrictions\n";
1586   for (CeedInt i = 0; i < num_output_fields; i++) {
1587     const char *field_name;
1588 
1589     CeedCallBackend(CeedOperatorFieldGetName(op_output_fields[i], &field_name));
1590     code << tab << "// ---- Output field " << i << ": " << field_name << "\n";
1591 
1592     // ---- Basis action
1593     CeedCallBackend(CeedOperatorBuildKernelBasis_Cuda_gen(code, data, tab, i, op_output_fields[i], qf_output_fields[i], max_dim, Q_1d, false,
1594                                                           is_all_tensor, is_at_points, use_3d_slices));
1595 
1596     // ---- Restriction
1597     CeedCallBackend(CeedOperatorBuildKernelRestriction_Cuda_gen(code, data, tab, i, NULL, op_output_fields[i], qf_output_fields[i], max_dim, Q_1d,
1598                                                                 false, is_all_tensor, is_at_points, use_3d_slices));
1599   }
1600 
1601   // Close loop and function
1602   tab.pop();
1603   code << tab << "}\n";
1604   tab.pop();
1605   code << tab << "}\n";
1606   code << tab << "// -----------------------------------------------------------------------------\n\n";
1607 
1608   // Compile
1609   {
1610     bool          is_compile_good = false;
1611     const CeedInt T_1d            = CeedIntMax(is_all_tensor ? Q_1d : Q, data->max_P_1d);
1612 
1613     data->thread_1d = T_1d;
1614     CeedCallBackend(CeedTryCompile_Cuda(ceed, code.str().c_str(), &is_compile_good, &data->module, 1, "OP_T_1D", T_1d));
1615     if (is_compile_good) {
1616       *is_good_build = true;
1617       CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, operator_name.c_str(), &data->op));
1618     } else {
1619       *is_good_build     = false;
1620       data->use_fallback = true;
1621     }
1622   }
1623   CeedCallBackend(CeedOperatorSetSetupDone(op));
1624   CeedCallBackend(CeedDestroy(&ceed));
1625   CeedCallBackend(CeedQFunctionDestroy(&qf));
1626   return CEED_ERROR_SUCCESS;
1627 }
1628 
1629 //------------------------------------------------------------------------------
1630 // Build AtPoints assembly operator kernel
1631 //------------------------------------------------------------------------------
1632 static int CeedOperatorBuildKernelAssemblyAtPoints_Cuda_gen(CeedOperator op, bool is_full, bool *is_good_build) {
1633   bool                    is_all_tensor = true, is_at_points = false, use_3d_slices = false;
1634   Ceed                    ceed;
1635   CeedInt                 Q, Q_1d, num_input_fields, num_output_fields, max_dim = 1, max_num_points = 0, coords_comp_stride = 0;
1636   CeedQFunctionField     *qf_input_fields, *qf_output_fields;
1637   CeedQFunction_Cuda_gen *qf_data;
1638   CeedQFunction           qf;
1639   CeedOperatorField      *op_input_fields, *op_output_fields;
1640   CeedOperator_Cuda_gen  *data;
1641   std::ostringstream      code;
1642   Tab                     tab;
1643 
1644   // Check compatibility
1645   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
1646   CeedCallBackend(CeedOperatorIsAtPoints(op, &is_at_points));
1647   CeedCheck(is_at_points, ceed, CEED_ERROR_BACKEND, "Only AtPoints operator assembly supported");
1648 
1649   // Retrieve operator data
1650   CeedCallBackend(CeedOperatorGetData(op, &data));
1651   Q       = data->Q;
1652   Q_1d    = data->Q_1d;
1653   max_dim = data->dim;
1654   {
1655     CeedElemRestriction rstr_points = NULL;
1656 
1657     CeedCallBackend(CeedOperatorAtPointsGetPoints(op, &rstr_points, NULL));
1658     CeedCallBackend(CeedElemRestrictionGetMaxPointsInElement(rstr_points, &max_num_points));
1659     CeedCallBackend(CeedElemRestrictionGetCompStride(rstr_points, &coords_comp_stride));
1660     CeedCallBackend(CeedElemRestrictionDestroy(&rstr_points));
1661   }
1662   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
1663   CeedCallBackend(CeedQFunctionGetData(qf, &qf_data));
1664   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
1665   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
1666 
1667   // Add atomicAdd function for old NVidia architectures
1668   {
1669     Ceed_Cuda            *ceed_data;
1670     struct cudaDeviceProp prop;
1671 
1672     CeedCallBackend(CeedGetData(ceed, &ceed_data));
1673     CeedCallBackend(cudaGetDeviceProperties(&prop, ceed_data->device_id));
1674     if ((prop.major < 6) && (CEED_SCALAR_TYPE != CEED_SCALAR_FP32)) {
1675       code << tab << "// AtomicAdd fallback source\n";
1676       code << tab << "#include <ceed/jit-source/cuda/cuda-atomic-add-fallback.h>\n\n";
1677     }
1678   }
1679 
1680   // Load basis source files
1681   code << tab << "// Tensor basis source\n";
1682   code << tab << "#include <ceed/jit-source/cuda/cuda-shared-basis-tensor-templates.h>\n\n";
1683   code << tab << "// AtPoints basis source\n";
1684   code << tab << "#include <ceed/jit-source/cuda/cuda-shared-basis-tensor-at-points-templates.h>\n\n";
1685   code << tab << "// CodeGen operator source\n";
1686   code << tab << "#include <ceed/jit-source/cuda/cuda-gen-templates.h>\n\n";
1687 
1688   // Get QFunction name
1689   std::string qfunction_name(qf_data->qfunction_name);
1690   std::string operator_name;
1691 
1692   if (is_full) {
1693     operator_name = "CeedKernelCudaGenOperatorFullAssembly_" + qfunction_name;
1694   } else {
1695     operator_name = "CeedKernelCudaGenOperatorDiagonalAssembly_" + qfunction_name;
1696   }
1697 
1698   // Define CEED_Q_VLA
1699   code << "\n" << tab << "#undef CEED_Q_VLA\n";
1700   code << tab << "#define CEED_Q_VLA 1\n\n";
1701 
1702   // Add user QFunction source
1703   {
1704     const char *source_path;
1705 
1706     CeedCallBackend(CeedQFunctionGetSourcePath(qf, &source_path));
1707     CeedCheck(source_path, ceed, CEED_ERROR_UNSUPPORTED, "/gpu/cuda/gen backend requires QFunction source code file");
1708 
1709     code << tab << "// User QFunction source\n";
1710     code << tab << "#include \"" << source_path << "\"\n\n";
1711   }
1712 
1713   // Setup
1714   code << "\n" << tab << "// -----------------------------------------------------------------------------\n";
1715   code << tab << "// Operator Assembly Kernel\n";
1716   code << tab << "// \n";
1717   code << tab << "// d_[in,out]_i:   CeedVector device array\n";
1718   code << tab << "// r_[in,out]_e_i: Element vector register\n";
1719   code << tab << "// r_[in,out]_q_i: Quadrature space vector register\n";
1720   code << tab << "// r_[in,out]_c_i: AtPoints Chebyshev coefficients register\n";
1721   code << tab << "// r_[in,out]_s_i: Quadrature space slice vector register\n";
1722   code << tab << "// \n";
1723   code << tab << "// s_B_[in,out]_i: Interpolation matrix, shared memory\n";
1724   code << tab << "// s_G_[in,out]_i: Gradient matrix, shared memory\n";
1725   code << tab << "// -----------------------------------------------------------------------------\n";
1726   code << tab << "extern \"C\" __global__ void " << operator_name
1727        << "(CeedInt num_elem, void* ctx, FieldsInt_Cuda indices, Fields_Cuda fields, Fields_Cuda B, Fields_Cuda G, CeedScalar *W, Points_Cuda "
1728           "points, CeedScalar *__restrict__ values_array) {\n";
1729   tab.push();
1730 
1731   // Scratch buffers
1732   for (CeedInt i = 0; i < num_input_fields; i++) {
1733     CeedEvalMode eval_mode;
1734 
1735     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
1736     if (eval_mode != CEED_EVAL_WEIGHT) {  // Skip CEED_EVAL_WEIGHT
1737       code << tab << "const CeedScalar *__restrict__ d_in_" << i << " = fields.inputs[" << i << "];\n";
1738     }
1739   }
1740   for (CeedInt i = 0; i < num_output_fields; i++) {
1741     code << tab << "CeedScalar *__restrict__ d_out_" << i << " = fields.outputs[" << i << "];\n";
1742   }
1743 
1744   code << tab << "const CeedInt max_dim = " << max_dim << ";\n";
1745   code << tab << "const CeedInt Q_1d = " << Q_1d << ";\n";
1746   code << tab << "const CeedInt max_num_points = " << max_num_points << ";\n";
1747   code << tab << "const CeedInt coords_comp_stride = " << coords_comp_stride << ";\n";
1748 
1749   // Shared data
1750   code << tab << "extern __shared__ CeedScalar slice[];\n";
1751   code << tab << "SharedData_Cuda data;\n";
1752   code << tab << "data.t_id_x = threadIdx.x;\n";
1753   code << tab << "data.t_id_y = threadIdx.y;\n";
1754   code << tab << "data.t_id_z = threadIdx.z;\n";
1755   code << tab << "data.t_id   = threadIdx.x + threadIdx.y*blockDim.x + threadIdx.z*blockDim.y*blockDim.x;\n";
1756   code << tab << "data.slice  = slice + data.t_id_z*OP_T_1D" << ((!is_all_tensor || max_dim == 1) ? "" : "*OP_T_1D") << ";\n";
1757 
1758   // -- Determine input mat reuse
1759   FieldReuse_Cuda input_matrix_reuse[CEED_FIELD_MAX];
1760 
1761   for (CeedInt i = 0; i < num_input_fields; i++) {
1762     input_matrix_reuse[i].index = -1;
1763   }
1764   for (CeedInt i = 0; i < num_input_fields; i++) {
1765     CeedEvalMode eval_mode_i;
1766     CeedBasis    basis_i;
1767 
1768     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode_i));
1769     if (eval_mode_i == CEED_EVAL_WEIGHT) continue;
1770     CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis_i));
1771     for (CeedInt j = 0; (input_matrix_reuse[i].index == -1) && (j < i); j++) {
1772       CeedEvalMode eval_mode_j;
1773       CeedBasis    basis_j;
1774 
1775       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[j], &eval_mode_j));
1776       if (eval_mode_j == CEED_EVAL_WEIGHT) continue;
1777       CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[j], &basis_j));
1778       if (basis_i == basis_j) {
1779         input_matrix_reuse[i].index     = j;
1780         input_matrix_reuse[i].is_input  = true;
1781         input_matrix_reuse[i].eval_mode = eval_mode_j;
1782       }
1783       CeedCallBackend(CeedBasisDestroy(&basis_j));
1784     }
1785     CeedCallBackend(CeedBasisDestroy(&basis_i));
1786   }
1787 
1788   // -- Determine output mat reuse
1789   FieldReuse_Cuda output_matrix_reuse[CEED_FIELD_MAX];
1790 
1791   for (CeedInt i = 0; i < num_output_fields; i++) {
1792     output_matrix_reuse[i].index = -1;
1793   }
1794   for (CeedInt i = 0; i < num_output_fields; i++) {
1795     CeedEvalMode eval_mode_i;
1796     CeedBasis    basis_i;
1797 
1798     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode_i));
1799     CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis_i));
1800     for (CeedInt j = 0; (output_matrix_reuse[i].index == -1) && (j < num_input_fields); j++) {
1801       CeedEvalMode eval_mode_j;
1802       CeedBasis    basis_j;
1803 
1804       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[j], &eval_mode_j));
1805       if (eval_mode_j == CEED_EVAL_WEIGHT) continue;
1806       CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[j], &basis_j));
1807       if (basis_i == basis_j) {
1808         output_matrix_reuse[i].index     = j;
1809         output_matrix_reuse[i].is_input  = true;
1810         output_matrix_reuse[i].eval_mode = eval_mode_j;
1811       }
1812       CeedCallBackend(CeedBasisDestroy(&basis_j));
1813     }
1814     for (CeedInt j = 0; (output_matrix_reuse[i].index == -1) && (j < i); j++) {
1815       CeedEvalMode eval_mode_j;
1816       CeedBasis    basis_j;
1817 
1818       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[j], &eval_mode_j));
1819       if (eval_mode_j == CEED_EVAL_WEIGHT) continue;
1820       CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[j], &basis_j));
1821       if (basis_i == basis_j) {
1822         output_matrix_reuse[i].index     = j;
1823         output_matrix_reuse[i].is_input  = false;
1824         output_matrix_reuse[i].eval_mode = eval_mode_j;
1825       }
1826       CeedCallBackend(CeedBasisDestroy(&basis_j));
1827     }
1828     CeedCallBackend(CeedBasisDestroy(&basis_i));
1829   }
1830 
1831   // Initialize constants, and matrices B and G
1832   code << "\n" << tab << "// Input field constants and basis data\n";
1833   for (CeedInt i = 0; i < num_input_fields; i++) {
1834     CeedCallBackend(CeedOperatorBuildKernelFieldData_Cuda_gen(code, data, tab, i, op_input_fields[i], qf_input_fields[i], input_matrix_reuse[i],
1835                                                               max_dim, Q, Q_1d, true, is_all_tensor, is_at_points, use_3d_slices, false));
1836   }
1837   code << "\n" << tab << "// Output field constants and basis data\n";
1838   for (CeedInt i = 0; i < num_output_fields; i++) {
1839     CeedCallBackend(CeedOperatorBuildKernelFieldData_Cuda_gen(code, data, tab, i, op_output_fields[i], qf_output_fields[i], output_matrix_reuse[i],
1840                                                               max_dim, Q, Q_1d, false, is_all_tensor, is_at_points, use_3d_slices, false));
1841   }
1842 
1843   // Loop over all elements
1844   code << "\n" << tab << "// Element loop\n";
1845   code << tab << "__syncthreads();\n";
1846   code << tab << "for (CeedInt elem = blockIdx.x*blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x*blockDim.z) {\n";
1847   tab.push();
1848 
1849   // -- Compute minimum buffer space needed
1850   CeedInt max_rstr_buffer_size = 1;
1851 
1852   for (CeedInt i = 0; i < num_input_fields; i++) {
1853     CeedEvalMode eval_mode;
1854 
1855     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
1856     if (eval_mode != CEED_EVAL_NONE && eval_mode != CEED_EVAL_WEIGHT) {
1857       CeedInt             num_comp;
1858       CeedElemRestriction elem_rstr;
1859 
1860       CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr));
1861       CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
1862       max_rstr_buffer_size = CeedIntMax(max_rstr_buffer_size, num_comp * (is_all_tensor && (max_dim >= 3) ? Q_1d : 1));
1863       CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
1864     }
1865   }
1866   for (CeedInt i = 0; i < num_output_fields; i++) {
1867     CeedEvalMode eval_mode;
1868 
1869     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
1870     if (eval_mode != CEED_EVAL_NONE) {
1871       CeedInt             num_comp;
1872       CeedElemRestriction elem_rstr;
1873 
1874       CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &elem_rstr));
1875       CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
1876       max_rstr_buffer_size = CeedIntMax(max_rstr_buffer_size, num_comp * (is_all_tensor && (max_dim >= 3) ? Q_1d : 1));
1877       CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
1878     }
1879   }
1880   code << tab << "// Scratch restriction buffer space\n";
1881   code << tab << "CeedScalar r_e_scratch[" << max_rstr_buffer_size << "];\n";
1882 
1883   // -- Determine best input field processing order
1884   CeedInt field_rstr_in_buffer[CEED_FIELD_MAX], input_field_order[CEED_FIELD_MAX];
1885 
1886   for (CeedInt i = 0; i < num_input_fields; i++) {
1887     field_rstr_in_buffer[i] = -1;
1888     input_field_order[i]    = -1;
1889   }
1890   {
1891     bool    is_ordered[CEED_FIELD_MAX];
1892     CeedInt curr_index = 0;
1893 
1894     for (CeedInt i = 0; i < num_input_fields; i++) is_ordered[i] = false;
1895     for (CeedInt i = 0; i < num_input_fields; i++) {
1896       CeedVector          vec_i;
1897       CeedElemRestriction rstr_i;
1898 
1899       if (is_ordered[i]) continue;
1900       field_rstr_in_buffer[i]       = i;
1901       is_ordered[i]                 = true;
1902       input_field_order[curr_index] = i;
1903       curr_index++;
1904       CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec_i));
1905       if (vec_i == CEED_VECTOR_NONE) continue;  // CEED_EVAL_WEIGHT
1906       CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &rstr_i));
1907       for (CeedInt j = i + 1; j < num_input_fields; j++) {
1908         CeedVector          vec_j;
1909         CeedElemRestriction rstr_j;
1910 
1911         CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[j], &vec_j));
1912         CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[j], &rstr_j));
1913         if (rstr_i == rstr_j && vec_i == vec_j) {
1914           field_rstr_in_buffer[j]       = i;
1915           is_ordered[j]                 = true;
1916           input_field_order[curr_index] = j;
1917           curr_index++;
1918         }
1919         CeedCallBackend(CeedVectorDestroy(&vec_j));
1920         CeedCallBackend(CeedElemRestrictionDestroy(&rstr_j));
1921       }
1922       CeedCallBackend(CeedVectorDestroy(&vec_i));
1923       CeedCallBackend(CeedElemRestrictionDestroy(&rstr_i));
1924     }
1925   }
1926 
1927   // -- Input restriction and basis
1928   code << "\n" << tab << "// -- Input field restrictions and basis actions\n";
1929   CeedInt active_field_index = -1;
1930 
1931   for (CeedInt i = 0; i < num_input_fields; i++) {
1932     bool          is_active = false;
1933     const char   *field_name;
1934     const CeedInt f = input_field_order[i];
1935 
1936     {
1937       CeedVector vec;
1938 
1939       CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[f], &vec));
1940       is_active = vec == CEED_VECTOR_ACTIVE;
1941       CeedCallBackend(CeedVectorDestroy(&vec));
1942     }
1943 
1944     CeedCallBackend(CeedOperatorFieldGetName(op_input_fields[f], &field_name));
1945     code << tab << "// ---- Input field " << f << ": " << field_name << "\n";
1946 
1947     if (is_active) {
1948       std::string var_suffix = "_in_" + std::to_string(f);
1949 
1950       code << tab << "// Active field - no restriction or basis action here\n";
1951       if (active_field_index == -1) {
1952         active_field_index = f;
1953         code << tab << "CeedScalar r_e" << var_suffix << "[num_comp" << var_suffix << "*" << (max_dim >= 3 ? "P_1d" + var_suffix : "1")
1954              << "] = {0.0};\n";
1955       } else {
1956         code << tab << "CeedScalar *r_e" << var_suffix << " = r_e_in_" << active_field_index << ";\n";
1957       }
1958     } else {
1959       // ---- Restriction
1960       CeedCallBackend(CeedOperatorBuildKernelRestriction_Cuda_gen(code, data, tab, f, field_rstr_in_buffer, op_input_fields[f], qf_input_fields[f],
1961                                                                   max_dim, Q_1d, true, is_all_tensor, is_at_points, use_3d_slices));
1962 
1963       // ---- Basis action
1964       CeedCallBackend(CeedOperatorBuildKernelBasis_Cuda_gen(code, data, tab, f, op_input_fields[f], qf_input_fields[f], max_dim, Q_1d, true,
1965                                                             is_all_tensor, is_at_points, use_3d_slices));
1966     }
1967   }
1968 
1969   // -- Loop over active field
1970   std::string active_var_suffix = "_in_" + std::to_string(active_field_index);
1971 
1972   code << "\n" << tab << "// Loop over nodes in active field\n";
1973   code << tab << "for (CeedInt n = 0; n < num_comp" << active_var_suffix << "*P_1d" << active_var_suffix
1974        << (max_dim > 1 ? "*P_1d" + active_var_suffix : "") << (max_dim > 2 ? "*P_1d" + active_var_suffix : "") << "; n++) {\n";
1975   tab.push();
1976 
1977   // -- Set current active node and component to 1
1978   code << tab << "// Set current active node and component to 1.0\n";
1979   code << tab << "SetEVecStandard" << max_dim << "d_Single<num_comp" << active_var_suffix << ", P_1d" << active_var_suffix << ">(data, n, 1.0, r_e"
1980        << active_var_suffix << ");\n\n";
1981 
1982   for (CeedInt i = 0; i < num_input_fields; i++) {
1983     bool          is_active = false;
1984     const char   *field_name;
1985     const CeedInt f = input_field_order[i];
1986 
1987     {
1988       CeedVector vec;
1989 
1990       CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[f], &vec));
1991       is_active = vec == CEED_VECTOR_ACTIVE;
1992       CeedCallBackend(CeedVectorDestroy(&vec));
1993     }
1994     if (!is_active) continue;
1995 
1996     CeedCallBackend(CeedOperatorFieldGetName(op_input_fields[f], &field_name));
1997     code << tab << "// ---- Input field " << f << ": " << field_name << "\n";
1998 
1999     // ---- Basis action
2000     CeedCallBackend(CeedOperatorBuildKernelBasis_Cuda_gen(code, data, tab, f, op_input_fields[f], qf_input_fields[f], max_dim, Q_1d, true,
2001                                                           is_all_tensor, is_at_points, use_3d_slices));
2002   }
2003 
2004   // -- Q function
2005   CeedCallBackend(CeedOperatorBuildKernelQFunction_Cuda_gen(code, data, tab, max_dim, max_num_points, num_input_fields, op_input_fields,
2006                                                             qf_input_fields, num_output_fields, op_output_fields, qf_output_fields, qfunction_name,
2007                                                             Q_1d, is_all_tensor, is_at_points, use_3d_slices));
2008 
2009   // -- Output basis and restriction
2010   code << "\n" << tab << "// -- Output field basis action and restrictions\n";
2011   for (CeedInt i = 0; i < num_output_fields; i++) {
2012     bool        is_active = false;
2013     const char *field_name;
2014 
2015     {
2016       CeedVector vec;
2017 
2018       CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[i], &vec));
2019       is_active = vec == CEED_VECTOR_ACTIVE;
2020       CeedCallBackend(CeedVectorDestroy(&vec));
2021     }
2022     if (!is_active) continue;
2023 
2024     CeedCallBackend(CeedOperatorFieldGetName(op_output_fields[i], &field_name));
2025     code << tab << "// ---- Output field " << i << ": " << field_name << "\n";
2026 
2027     // ---- Basis action
2028     CeedCallBackend(CeedOperatorBuildKernelBasis_Cuda_gen(code, data, tab, i, op_output_fields[i], qf_output_fields[i], max_dim, Q_1d, false,
2029                                                           is_all_tensor, is_at_points, use_3d_slices));
2030 
2031     // ---- Restriction
2032     if (is_full) {
2033       std::string         var_suffix = "_out_" + std::to_string(i);
2034       CeedInt             comp_stride;
2035       CeedSize            l_size;
2036       CeedElemRestriction elem_rstr;
2037 
2038       CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &elem_rstr));
2039       CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size));
2040       code << tab << "const CeedInt l_size" << var_suffix << " = " << l_size << ";\n";
2041       CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
2042       code << tab << "const CeedInt comp_stride" << var_suffix << " = " << comp_stride << ";\n";
2043       code << tab << "WriteLVecStandard" << max_dim << "d_Assembly<num_comp" << var_suffix << ", comp_stride" << var_suffix << ", P_1d" + var_suffix
2044            << ">(data, l_size" << var_suffix << ", elem, n, r_e" << var_suffix << ", values_array);\n";
2045       CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
2046     } else {
2047       std::string         var_suffix = "_out_" + std::to_string(i);
2048       CeedInt             comp_stride;
2049       CeedSize            l_size;
2050       CeedElemRestriction elem_rstr;
2051 
2052       CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &elem_rstr));
2053       CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size));
2054       code << tab << "const CeedInt l_size" << var_suffix << " = " << l_size << ";\n";
2055       CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
2056       code << tab << "const CeedInt comp_stride" << var_suffix << " = " << comp_stride << ";\n";
2057       code << tab << "WriteLVecStandard" << max_dim << "d_Single<num_comp" << var_suffix << ", comp_stride" << var_suffix << ", P_1d" + var_suffix
2058            << ">(data, l_size" << var_suffix << ", elem, n, indices.outputs[" << i << "], r_e" << var_suffix << ", values_array);\n";
2059       CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
2060     }
2061   }
2062 
2063   // -- Reset current active node and component
2064   code << "\n" << tab << "// Reset current active node and component to 0.0\n";
2065   code << tab << "SetEVecStandard" << max_dim << "d_Single<num_comp" << active_var_suffix << ", P_1d" << active_var_suffix << ">(data, n, 0.0, r_e"
2066        << active_var_suffix << ");\n";
2067 
2068   // -- End of loop over active field
2069   tab.pop();
2070   code << tab << "}\n";
2071 
2072   // Close loop and function
2073   tab.pop();
2074   code << tab << "}\n";
2075   tab.pop();
2076   code << tab << "}\n";
2077   code << tab << "// -----------------------------------------------------------------------------\n\n";
2078 
2079   // Compile
2080   {
2081     bool          is_compile_good = false;
2082     const CeedInt T_1d            = CeedIntMax(is_all_tensor ? Q_1d : Q, data->max_P_1d);
2083 
2084     data->thread_1d = T_1d;
2085     CeedCallBackend(CeedTryCompile_Cuda(ceed, code.str().c_str(), &is_compile_good,
2086                                         is_full ? &data->module_assemble_full : &data->module_assemble_diagonal, 1, "OP_T_1D", T_1d));
2087     if (is_compile_good) {
2088       *is_good_build = true;
2089       CeedCallBackend(CeedGetKernel_Cuda(ceed, is_full ? data->module_assemble_full : data->module_assemble_diagonal, operator_name.c_str(),
2090                                          is_full ? &data->assemble_full : &data->assemble_diagonal));
2091     } else {
2092       *is_good_build              = false;
2093       data->use_assembly_fallback = true;
2094     }
2095   }
2096   CeedCallBackend(CeedDestroy(&ceed));
2097   CeedCallBackend(CeedQFunctionDestroy(&qf));
2098   return CEED_ERROR_SUCCESS;
2099 }
2100 
2101 extern "C" int CeedOperatorBuildKernelDiagonalAssemblyAtPoints_Cuda_gen(CeedOperator op, bool *is_good_build) {
2102   return CeedOperatorBuildKernelAssemblyAtPoints_Cuda_gen(op, false, is_good_build);
2103 }
2104 
2105 extern "C" int CeedOperatorBuildKernelFullAssemblyAtPoints_Cuda_gen(CeedOperator op, bool *is_good_build) {
2106   return CeedOperatorBuildKernelAssemblyAtPoints_Cuda_gen(op, true, is_good_build);
2107 }
2108 
2109 //------------------------------------------------------------------------------
2110 // Build QFunction assembly operator kernel
2111 //------------------------------------------------------------------------------
2112 extern "C" int CeedOperatorBuildKernelLinearAssembleQFunction_Cuda_gen(CeedOperator op, bool *is_good_build) {
2113   bool                    is_all_tensor = true, is_all_nontensor = true, is_at_points = false, use_3d_slices = false;
2114   Ceed                    ceed;
2115   CeedInt                 Q, Q_1d, num_input_fields, num_output_fields, max_dim = 1, max_num_points = 0;
2116   CeedQFunctionField     *qf_input_fields, *qf_output_fields;
2117   CeedQFunction_Cuda_gen *qf_data;
2118   CeedQFunction           qf;
2119   CeedOperatorField      *op_input_fields, *op_output_fields;
2120   CeedOperator_Cuda_gen  *data;
2121   std::ostringstream      code;
2122   Tab                     tab;
2123 
2124   // Check compatibility
2125   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
2126   CeedCallBackend(CeedOperatorIsAtPoints(op, &is_at_points));
2127   CeedCheck(!is_at_points, ceed, CEED_ERROR_BACKEND, "AtPoints QFunction assembly is not supported");
2128 
2129   // Check field compatibility
2130   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
2131   {
2132     bool has_shared_bases = true;
2133 
2134     for (CeedInt i = 0; i < num_input_fields; i++) {
2135       CeedBasis basis;
2136 
2137       CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
2138       if (basis != CEED_BASIS_NONE) {
2139         bool        is_tensor = true;
2140         const char *resource;
2141         char       *resource_root;
2142         Ceed        basis_ceed;
2143 
2144         CeedCallBackend(CeedBasisIsTensor(basis, &is_tensor));
2145         is_all_tensor    = is_all_tensor && is_tensor;
2146         is_all_nontensor = is_all_nontensor && !is_tensor;
2147         CeedCallBackend(CeedBasisGetCeed(basis, &basis_ceed));
2148         CeedCallBackend(CeedGetResource(basis_ceed, &resource));
2149         CeedCallBackend(CeedGetResourceRoot(basis_ceed, resource, ":", &resource_root));
2150         has_shared_bases = has_shared_bases && !strcmp(resource_root, "/gpu/cuda/shared");
2151         CeedCallBackend(CeedFree(&resource_root));
2152         CeedCallBackend(CeedDestroy(&basis_ceed));
2153       }
2154       CeedCallBackend(CeedBasisDestroy(&basis));
2155     }
2156 
2157     for (CeedInt i = 0; i < num_output_fields; i++) {
2158       CeedBasis basis;
2159 
2160       CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
2161       if (basis != CEED_BASIS_NONE) {
2162         bool        is_tensor = true;
2163         const char *resource;
2164         char       *resource_root;
2165         Ceed        basis_ceed;
2166 
2167         CeedCallBackend(CeedBasisIsTensor(basis, &is_tensor));
2168         is_all_tensor    = is_all_tensor && is_tensor;
2169         is_all_nontensor = is_all_nontensor && !is_tensor;
2170 
2171         CeedCallBackend(CeedBasisGetCeed(basis, &basis_ceed));
2172         CeedCallBackend(CeedGetResource(basis_ceed, &resource));
2173         CeedCallBackend(CeedGetResourceRoot(basis_ceed, resource, ":", &resource_root));
2174         has_shared_bases = has_shared_bases && !strcmp(resource_root, "/gpu/cuda/shared");
2175         CeedCallBackend(CeedFree(&resource_root));
2176         CeedCallBackend(CeedDestroy(&basis_ceed));
2177       }
2178       CeedCallBackend(CeedBasisDestroy(&basis));
2179     }
2180   }
2181 
2182   // Retrieve operator data
2183   CeedCallBackend(CeedOperatorGetData(op, &data));
2184   Q       = data->Q;
2185   Q_1d    = data->Q_1d;
2186   max_dim = data->dim;
2187   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
2188   CeedCallBackend(CeedQFunctionGetData(qf, &qf_data));
2189   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
2190 
2191   // Add atomicAdd function for old NVidia architectures
2192   {
2193     Ceed_Cuda            *ceed_data;
2194     struct cudaDeviceProp prop;
2195 
2196     CeedCallBackend(CeedGetData(ceed, &ceed_data));
2197     CeedCallBackend(cudaGetDeviceProperties(&prop, ceed_data->device_id));
2198     if ((prop.major < 6) && (CEED_SCALAR_TYPE != CEED_SCALAR_FP32)) {
2199       code << tab << "// AtomicAdd fallback source\n";
2200       code << tab << "#include <ceed/jit-source/cuda/cuda-atomic-add-fallback.h>\n\n";
2201     }
2202   }
2203 
2204   // Load basis source files
2205   if (!is_all_nontensor) {
2206     code << tab << "// Tensor basis source\n";
2207     code << tab << "#include <ceed/jit-source/cuda/cuda-shared-basis-tensor-templates.h>\n\n";
2208   }
2209   if (!is_all_tensor) {
2210     code << tab << "// Non-tensor basis source\n";
2211     code << tab << "#include <ceed/jit-source/cuda/cuda-shared-basis-nontensor-templates.h>\n\n";
2212   }
2213   if (!is_all_tensor && !is_all_nontensor) {
2214     code << "// Tensor basis source\n";
2215     code << "#include <ceed/jit-source/cuda/cuda-shared-basis-tensor-flattened-templates.h>\n\n";
2216   }
2217   code << "// CodeGen operator source\n";
2218   code << "#include <ceed/jit-source/cuda/cuda-gen-templates.h>\n\n";
2219 
2220   // Get QFunction name
2221   std::string qfunction_name(qf_data->qfunction_name);
2222   std::string operator_name;
2223 
2224   operator_name = "CeedKernelCudaGenQFunctionAssembly_" + qfunction_name;
2225 
2226   // Define CEED_Q_VLA
2227   code << "\n" << tab << "#undef CEED_Q_VLA\n";
2228   if (max_dim != 3 || is_at_points || use_3d_slices || !is_all_tensor) {
2229     code << tab << "#define CEED_Q_VLA 1\n\n";
2230   } else {
2231     code << tab << "#define CEED_Q_VLA " << Q_1d << "\n\n";
2232   }
2233 
2234   // Add user QFunction source
2235   {
2236     const char *source_path;
2237 
2238     CeedCallBackend(CeedQFunctionGetSourcePath(qf, &source_path));
2239     CeedCheck(source_path, ceed, CEED_ERROR_UNSUPPORTED, "/gpu/cuda/gen backend requires QFunction source code file");
2240 
2241     code << tab << "// User QFunction source\n";
2242     code << tab << "#include \"" << source_path << "\"\n\n";
2243   }
2244 
2245   // Setup
2246   code << "\n" << tab << "// -----------------------------------------------------------------------------\n";
2247   code << tab << "// Operator Assembly Kernel\n";
2248   code << tab << "// \n";
2249   code << tab << "// d_[in,out]_i:   CeedVector device array\n";
2250   code << tab << "// r_[in,out]_e_i: Element vector register\n";
2251   code << tab << "// r_[in,out]_q_i: Quadrature space vector register\n";
2252   code << tab << "// r_[in,out]_c_i: AtPoints Chebyshev coefficients register\n";
2253   code << tab << "// r_[in,out]_s_i: Quadrature space slice vector register\n";
2254   code << tab << "// \n";
2255   code << tab << "// s_B_[in,out]_i: Interpolation matrix, shared memory\n";
2256   code << tab << "// s_G_[in,out]_i: Gradient matrix, shared memory\n";
2257   code << tab << "// -----------------------------------------------------------------------------\n";
2258   code << tab << "extern \"C\" __global__ void " << operator_name
2259        << "(CeedInt num_elem, void* ctx, FieldsInt_Cuda indices, Fields_Cuda fields, Fields_Cuda B, Fields_Cuda G, CeedScalar *W, Points_Cuda "
2260           "points, CeedScalar *__restrict__ values_array) {\n";
2261   tab.push();
2262 
2263   // Scratch buffers
2264   for (CeedInt i = 0; i < num_input_fields; i++) {
2265     CeedEvalMode eval_mode;
2266 
2267     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
2268     if (eval_mode != CEED_EVAL_WEIGHT) {  // Skip CEED_EVAL_WEIGHT
2269       code << tab << "const CeedScalar *__restrict__ d_in_" << i << " = fields.inputs[" << i << "];\n";
2270     }
2271   }
2272   for (CeedInt i = 0; i < num_output_fields; i++) {
2273     code << tab << "CeedScalar *__restrict__ d_out_" << i << " = fields.outputs[" << i << "];\n";
2274   }
2275 
2276   code << tab << "const CeedInt max_dim = " << max_dim << ";\n";
2277   if (!is_all_tensor) {
2278     code << tab << "const CeedInt Q = " << Q << ";\n";
2279   }
2280   if (!is_all_nontensor) {
2281     code << tab << "const CeedInt Q_1d = " << Q_1d << ";\n";
2282   }
2283 
2284   // Shared data
2285   code << tab << "extern __shared__ CeedScalar slice[];\n";
2286   code << tab << "SharedData_Cuda data;\n";
2287   code << tab << "data.t_id_x = threadIdx.x;\n";
2288   code << tab << "data.t_id_y = threadIdx.y;\n";
2289   code << tab << "data.t_id_z = threadIdx.z;\n";
2290   code << tab << "data.t_id   = threadIdx.x + threadIdx.y*blockDim.x + threadIdx.z*blockDim.y*blockDim.x;\n";
2291   code << tab << "data.slice  = slice + data.t_id_z*OP_T_1D" << ((!is_all_tensor || max_dim == 1) ? "" : "*OP_T_1D") << ";\n";
2292 
2293   // -- Determine input mat reuse
2294   FieldReuse_Cuda input_matrix_reuse[CEED_FIELD_MAX];
2295 
2296   for (CeedInt i = 0; i < num_input_fields; i++) {
2297     input_matrix_reuse[i].index = -1;
2298   }
2299   for (CeedInt i = 0; i < num_input_fields; i++) {
2300     bool         is_tensor = true;
2301     CeedEvalMode eval_mode_i;
2302     CeedBasis    basis_i;
2303 
2304     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode_i));
2305     if (eval_mode_i == CEED_EVAL_WEIGHT) continue;
2306     CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis_i));
2307     CeedCallBackend(CeedBasisIsTensor(basis_i, &is_tensor));
2308     for (CeedInt j = 0; (input_matrix_reuse[i].index == -1) && (j < i); j++) {
2309       CeedEvalMode eval_mode_j;
2310       CeedBasis    basis_j;
2311 
2312       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[j], &eval_mode_j));
2313       if (eval_mode_j == CEED_EVAL_WEIGHT) continue;
2314       CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[j], &basis_j));
2315       if (basis_i == basis_j) {
2316         if (is_tensor) {
2317           input_matrix_reuse[i].index     = j;
2318           input_matrix_reuse[i].is_input  = true;
2319           input_matrix_reuse[i].eval_mode = eval_mode_j;
2320         } else {
2321           // For non-tensor can only re-use with the same eval mode
2322           if (eval_mode_i == eval_mode_j) {
2323             input_matrix_reuse[i].index     = j;
2324             input_matrix_reuse[i].is_input  = true;
2325             input_matrix_reuse[i].eval_mode = eval_mode_j;
2326           }
2327         }
2328       }
2329       CeedCallBackend(CeedBasisDestroy(&basis_j));
2330     }
2331     CeedCallBackend(CeedBasisDestroy(&basis_i));
2332   }
2333 
2334   // -- Determine output mat reuse
2335   FieldReuse_Cuda output_matrix_reuse[CEED_FIELD_MAX];
2336 
2337   for (CeedInt i = 0; i < num_output_fields; i++) {
2338     output_matrix_reuse[i].index = -1;
2339   }
2340   for (CeedInt i = 0; i < num_output_fields; i++) {
2341     bool         is_tensor = true;
2342     CeedEvalMode eval_mode_i;
2343     CeedBasis    basis_i;
2344 
2345     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode_i));
2346     CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis_i));
2347     CeedCallBackend(CeedBasisIsTensor(basis_i, &is_tensor));
2348     for (CeedInt j = 0; (output_matrix_reuse[i].index == -1) && (j < num_input_fields); j++) {
2349       CeedEvalMode eval_mode_j;
2350       CeedBasis    basis_j;
2351 
2352       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[j], &eval_mode_j));
2353       if (eval_mode_j == CEED_EVAL_WEIGHT) continue;
2354       CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[j], &basis_j));
2355       if (basis_i == basis_j) {
2356         if (is_tensor) {
2357           output_matrix_reuse[i].index     = j;
2358           output_matrix_reuse[i].is_input  = true;
2359           output_matrix_reuse[i].eval_mode = eval_mode_j;
2360         } else {
2361           // For non-tensor can only re-use with the same eval mode
2362           if (eval_mode_i == eval_mode_j) {
2363             output_matrix_reuse[i].index     = j;
2364             output_matrix_reuse[i].is_input  = true;
2365             output_matrix_reuse[i].eval_mode = eval_mode_j;
2366           }
2367         }
2368       }
2369       CeedCallBackend(CeedBasisDestroy(&basis_j));
2370     }
2371     for (CeedInt j = 0; (output_matrix_reuse[i].index == -1) && (j < i); j++) {
2372       CeedEvalMode eval_mode_j;
2373       CeedBasis    basis_j;
2374 
2375       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[j], &eval_mode_j));
2376       if (eval_mode_j == CEED_EVAL_WEIGHT) continue;
2377       CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[j], &basis_j));
2378       if (basis_i == basis_j) {
2379         if (is_tensor) {
2380           output_matrix_reuse[i].index     = j;
2381           output_matrix_reuse[i].is_input  = false;
2382           output_matrix_reuse[i].eval_mode = eval_mode_j;
2383         } else {
2384           // For non-tensor can only re-use with the same eval mode
2385           if (eval_mode_i == eval_mode_j) {
2386             output_matrix_reuse[i].index     = j;
2387             output_matrix_reuse[i].is_input  = false;
2388             output_matrix_reuse[i].eval_mode = eval_mode_j;
2389           }
2390         }
2391       }
2392       CeedCallBackend(CeedBasisDestroy(&basis_j));
2393     }
2394     CeedCallBackend(CeedBasisDestroy(&basis_i));
2395   }
2396 
2397   // Initialize constants, and matrices B and G
2398   code << "\n" << tab << "// Input field constants and basis data\n";
2399   for (CeedInt i = 0; i < num_input_fields; i++) {
2400     CeedCallBackend(CeedOperatorBuildKernelFieldData_Cuda_gen(code, data, tab, i, op_input_fields[i], qf_input_fields[i], input_matrix_reuse[i],
2401                                                               max_dim, Q, Q_1d, true, is_all_tensor, is_at_points, use_3d_slices, true));
2402   }
2403   code << "\n" << tab << "// Output field constants and basis data\n";
2404   for (CeedInt i = 0; i < num_output_fields; i++) {
2405     CeedCallBackend(CeedOperatorBuildKernelFieldData_Cuda_gen(code, data, tab, i, op_output_fields[i], qf_output_fields[i], output_matrix_reuse[i],
2406                                                               max_dim, Q, Q_1d, false, is_all_tensor, is_at_points, use_3d_slices, true));
2407   }
2408 
2409   // Loop over all elements
2410   code << "\n" << tab << "// Element loop\n";
2411   code << tab << "__syncthreads();\n";
2412   code << tab << "for (CeedInt elem = blockIdx.x*blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x*blockDim.z) {\n";
2413   tab.push();
2414 
2415   // -- Compute minimum buffer space needed
2416   CeedInt max_rstr_buffer_size = 1;
2417 
2418   for (CeedInt i = 0; i < num_input_fields; i++) {
2419     CeedEvalMode eval_mode;
2420 
2421     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
2422     if (eval_mode != CEED_EVAL_NONE && eval_mode != CEED_EVAL_WEIGHT) {
2423       CeedInt             num_comp;
2424       CeedElemRestriction elem_rstr;
2425 
2426       CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr));
2427       CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
2428       max_rstr_buffer_size = CeedIntMax(max_rstr_buffer_size, num_comp * (is_all_tensor && (max_dim >= 3) ? Q_1d : 1));
2429       CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
2430     }
2431   }
2432   for (CeedInt i = 0; i < num_output_fields; i++) {
2433     CeedEvalMode eval_mode;
2434 
2435     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
2436     if (eval_mode != CEED_EVAL_NONE) {
2437       CeedInt             num_comp;
2438       CeedElemRestriction elem_rstr;
2439 
2440       CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &elem_rstr));
2441       CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
2442       max_rstr_buffer_size = CeedIntMax(max_rstr_buffer_size, num_comp * (is_all_tensor && (max_dim >= 3) ? Q_1d : 1));
2443       CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
2444     }
2445   }
2446   code << tab << "// Scratch restriction buffer space\n";
2447   code << tab << "CeedScalar r_e_scratch[" << max_rstr_buffer_size << "];\n";
2448 
2449   // -- Determine best input field processing order
2450   CeedInt field_rstr_in_buffer[CEED_FIELD_MAX], input_field_order[CEED_FIELD_MAX];
2451 
2452   for (CeedInt i = 0; i < num_input_fields; i++) {
2453     field_rstr_in_buffer[i] = -1;
2454     input_field_order[i]    = -1;
2455   }
2456   {
2457     bool    is_ordered[CEED_FIELD_MAX];
2458     CeedInt curr_index = 0;
2459 
2460     for (CeedInt i = 0; i < num_input_fields; i++) is_ordered[i] = false;
2461     for (CeedInt i = 0; i < num_input_fields; i++) {
2462       CeedVector          vec_i;
2463       CeedElemRestriction rstr_i;
2464 
2465       if (is_ordered[i]) continue;
2466       field_rstr_in_buffer[i]       = i;
2467       is_ordered[i]                 = true;
2468       input_field_order[curr_index] = i;
2469       curr_index++;
2470       CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec_i));
2471       if (vec_i == CEED_VECTOR_NONE) continue;  // CEED_EVAL_WEIGHT
2472       CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &rstr_i));
2473       for (CeedInt j = i + 1; j < num_input_fields; j++) {
2474         CeedVector          vec_j;
2475         CeedElemRestriction rstr_j;
2476 
2477         CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[j], &vec_j));
2478         CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[j], &rstr_j));
2479         if (rstr_i == rstr_j && vec_i == vec_j) {
2480           field_rstr_in_buffer[j]       = i;
2481           is_ordered[j]                 = true;
2482           input_field_order[curr_index] = j;
2483           curr_index++;
2484         }
2485         CeedCallBackend(CeedVectorDestroy(&vec_j));
2486         CeedCallBackend(CeedElemRestrictionDestroy(&rstr_j));
2487       }
2488       CeedCallBackend(CeedVectorDestroy(&vec_i));
2489       CeedCallBackend(CeedElemRestrictionDestroy(&rstr_i));
2490     }
2491   }
2492 
2493   // -- Input restriction and basis
2494   code << "\n" << tab << "// -- Input field restrictions and basis actions\n";
2495   CeedInt num_active_in = 0, num_active_out = 0, qf_assembly_size_out = 0;
2496   CeedInt active_fields_in[CEED_FIELD_MAX], active_fields_out[CEED_FIELD_MAX];
2497 
2498   for (CeedInt i = 0; i < num_input_fields; i++) {
2499     bool          is_active = false;
2500     const char   *field_name;
2501     const CeedInt f = input_field_order[i];
2502 
2503     {
2504       CeedVector vec;
2505 
2506       CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[f], &vec));
2507       is_active = vec == CEED_VECTOR_ACTIVE;
2508       CeedCallBackend(CeedVectorDestroy(&vec));
2509     }
2510 
2511     CeedCallBackend(CeedOperatorFieldGetName(op_input_fields[f], &field_name));
2512     code << tab << "// ---- Input field " << f << ": " << field_name << "\n";
2513 
2514     if (is_active) {
2515       CeedEvalMode eval_mode;
2516       CeedInt      field_size;
2517 
2518       active_fields_in[num_active_in] = f;
2519       num_active_in++;
2520       CeedCallBackend(CeedQFunctionFieldGetSize(qf_input_fields[f], &field_size));
2521       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[f], &eval_mode));
2522       if (eval_mode == CEED_EVAL_GRAD) {
2523         code << tab << "CeedScalar r_q_in_" << f << "[num_comp_in_" << f << "*" << "dim_in_" << f << "*"
2524              << (is_all_tensor && (max_dim >= 3) ? "Q_1d" : "1") << "] = {0.};\n";
2525       } else {
2526         code << tab << "CeedScalar r_q_in_" << f << "[num_comp_in_" << f << "*" << (is_all_tensor && (max_dim >= 3) ? "Q_1d" : "1") << "] = {0.};\n";
2527       }
2528       code << tab << "const CeedInt field_size_in_" << f << " = " << field_size << ";\n";
2529     } else {
2530       // ---- Restriction
2531       CeedCallBackend(CeedOperatorBuildKernelRestriction_Cuda_gen(code, data, tab, f, field_rstr_in_buffer, op_input_fields[f], qf_input_fields[f],
2532                                                                   max_dim, Q_1d, true, is_all_tensor, is_at_points, use_3d_slices));
2533 
2534       // ---- Basis action
2535       CeedCallBackend(CeedOperatorBuildKernelBasis_Cuda_gen(code, data, tab, f, op_input_fields[f], qf_input_fields[f], max_dim, Q_1d, true,
2536                                                             is_all_tensor, is_at_points, use_3d_slices));
2537     }
2538   }
2539   code << tab << "const CeedInt field_sizes_in[" << num_active_in << "] = {";
2540   for (CeedInt i = 0; i < num_active_in; i++) {
2541     code << "field_size_in_" << active_fields_in[i] << (i < num_active_in - 1 ? ", " : "");
2542   }
2543   code << "};\n";
2544   code << tab << "CeedScalar * r_q_in[" << num_active_in << "] = {";
2545   for (CeedInt i = 0; i < num_active_in; i++) {
2546     code << "r_q_in_" << active_fields_in[i] << (i < num_active_in - 1 ? ", " : "");
2547   }
2548   code << "};\n";
2549 
2550   for (CeedInt i = 0; i < num_output_fields; i++) {
2551     bool is_active = false;
2552 
2553     {
2554       CeedVector vec;
2555 
2556       CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[i], &vec));
2557       is_active = vec == CEED_VECTOR_ACTIVE;
2558       CeedCallBackend(CeedVectorDestroy(&vec));
2559     }
2560     if (is_active) {
2561       const char *field_name;
2562       CeedInt     field_size;
2563 
2564       active_fields_out[num_active_out] = i;
2565       num_active_out++;
2566       CeedCallBackend(CeedQFunctionFieldGetSize(qf_output_fields[i], &field_size));
2567       qf_assembly_size_out += field_size;
2568       CeedCallBackend(CeedOperatorFieldGetName(op_output_fields[i], &field_name));
2569       code << tab << "// ---- Output field " << i << ": " << field_name << "\n";
2570       code << tab << "const CeedInt field_size_out_" << i << " = " << field_size << ";\n";
2571     }
2572   }
2573   code << tab << "const CeedInt field_sizes_out[" << num_active_out << "] = {";
2574   for (CeedInt i = 0; i < num_active_out; i++) {
2575     code << "field_size_out_" << active_fields_out[i] << (i < num_active_out - 1 ? ", " : "");
2576   }
2577   code << "};\n";
2578   code << tab << "const CeedInt total_size_out = " << qf_assembly_size_out << ";\n";
2579 
2580   // -- Loop over active field
2581   code << "\n" << tab << "CeedInt input_offset = 0;\n";
2582   code << tab << "// Loop over active QFunction input fields\n";
2583   code << tab << "const CeedInt num_active_in = " << num_active_in << ";\n";
2584   code << tab << "for (CeedInt a = 0; a < num_active_in; a++) {\n";
2585   tab.push();
2586 
2587   // -- Loop over size of active field
2588   code << "\n" << tab << "// Loop over current active input field size\n";
2589   code << tab << "const CeedInt field_size_in = field_sizes_in[a];\n";
2590   code << tab << "for (CeedInt s = 0; s < field_size_in; s++) {\n";
2591   tab.push();
2592 
2593   // -- Set current active point and component to 1
2594   code << tab << "// Set current active point and component to 1.0\n";
2595   if (is_all_tensor && (max_dim >= 3)) {
2596     code << tab << "for (CeedInt i = 0; i < Q_1d; i++) r_q_in[a][i + s * Q_1d] = 1.0;\n";
2597   } else {
2598     code << tab << "r_q_in[a][s] = 1.0;\n";
2599   }
2600 
2601   // -- Q function
2602   CeedCallBackend(CeedOperatorBuildKernelQFunction_Cuda_gen(code, data, tab, max_dim, max_num_points, num_input_fields, op_input_fields,
2603                                                             qf_input_fields, num_output_fields, op_output_fields, qf_output_fields, qfunction_name,
2604                                                             Q_1d, is_all_tensor, is_at_points, use_3d_slices));
2605 
2606   // -- Output basis and restriction
2607   code << "\n" << tab << "// -- Output field basis action and restrictions\n";
2608   CeedScalar offset = 0;
2609 
2610   for (CeedInt i = 0; i < num_output_fields; i++) {
2611     bool        is_active = false;
2612     const char *field_name;
2613 
2614     {
2615       CeedVector vec;
2616 
2617       CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[i], &vec));
2618       is_active = vec == CEED_VECTOR_ACTIVE;
2619       CeedCallBackend(CeedVectorDestroy(&vec));
2620     }
2621     if (!is_active) continue;
2622 
2623     CeedCallBackend(CeedOperatorFieldGetName(op_output_fields[i], &field_name));
2624     code << tab << "// ---- Output field " << i << ": " << field_name << "\n";
2625 
2626     // ---- Restriction
2627     CeedInt field_size;
2628 
2629     code << tab << "WriteLVecStandard" << (is_all_tensor ? max_dim : 1) << "d_QFAssembly<total_size_out, field_size_out_" << i << ", "
2630          << (is_all_tensor ? "Q_1d" : "Q") << ">(data, num_elem, elem, input_offset + s, " << offset << ", r_q_out_" << i << ", values_array);\n";
2631     CeedCallBackend(CeedQFunctionFieldGetSize(qf_output_fields[i], &field_size));
2632     offset += field_size;
2633   }
2634 
2635   // -- Reset current active node and component
2636   code << "\n" << tab << "// Reset current active node and component to 0.0\n";
2637   if (is_all_tensor && (max_dim >= 3)) {
2638     code << tab << "for (CeedInt i = 0; i < Q_1d; i++) r_q_in[a][i + s * Q_1d] = 0.0;\n";
2639   } else {
2640     code << tab << "r_q_in[a][s] = 0.0;\n";
2641   }
2642 
2643   // -- End of loop over size of active field
2644   tab.pop();
2645   code << tab << "}\n";
2646   code << tab << "input_offset += field_size_in;\n";
2647 
2648   // -- End of loop over active field
2649   tab.pop();
2650   code << tab << "}\n";
2651 
2652   // Close loop and function
2653   tab.pop();
2654   code << tab << "}\n";
2655   tab.pop();
2656   code << tab << "}\n";
2657   code << tab << "// -----------------------------------------------------------------------------\n\n";
2658 
2659   // Compile
2660   {
2661     bool          is_compile_good = false;
2662     const CeedInt T_1d            = CeedIntMax(is_all_tensor ? Q_1d : Q, data->max_P_1d);
2663 
2664     data->thread_1d = T_1d;
2665     CeedCallBackend(CeedTryCompile_Cuda(ceed, code.str().c_str(), &is_compile_good, &data->module_assemble_qfunction, 1, "OP_T_1D", T_1d));
2666     if (is_compile_good) {
2667       *is_good_build = true;
2668       CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module_assemble_qfunction, operator_name.c_str(), &data->assemble_qfunction));
2669     } else {
2670       *is_good_build              = false;
2671       data->use_assembly_fallback = true;
2672     }
2673   }
2674   CeedCallBackend(CeedDestroy(&ceed));
2675   CeedCallBackend(CeedQFunctionDestroy(&qf));
2676   return CEED_ERROR_SUCCESS;
2677 }
2678 
2679 //------------------------------------------------------------------------------
2680