xref: /libCEED/backends/cuda-ref/ceed-cuda-ref-operator.c (revision 021a32ddaeb3f5ae90defc18a1194608e9df3b1b)
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 #include <ceed.h>
9 #include <ceed/backend.h>
10 #include <ceed/jit-tools.h>
11 #include <assert.h>
12 #include <cuda.h>
13 #include <cuda_runtime.h>
14 #include <stdbool.h>
15 #include <string.h>
16 
17 #include "../cuda/ceed-cuda-common.h"
18 #include "../cuda/ceed-cuda-compile.h"
19 #include "ceed-cuda-ref.h"
20 
21 //------------------------------------------------------------------------------
22 // Destroy operator
23 //------------------------------------------------------------------------------
24 static int CeedOperatorDestroy_Cuda(CeedOperator op) {
25   CeedOperator_Cuda *impl;
26 
27   CeedCallBackend(CeedOperatorGetData(op, &impl));
28 
29   // Apply data
30   CeedCallBackend(CeedFree(&impl->num_points));
31   CeedCallBackend(CeedFree(&impl->skip_rstr_in));
32   CeedCallBackend(CeedFree(&impl->skip_rstr_out));
33   CeedCallBackend(CeedFree(&impl->apply_add_basis_out));
34   CeedCallBackend(CeedFree(&impl->input_field_order));
35   CeedCallBackend(CeedFree(&impl->output_field_order));
36   CeedCallBackend(CeedFree(&impl->input_states));
37 
38   for (CeedInt i = 0; i < impl->num_inputs; i++) {
39     CeedCallBackend(CeedVectorDestroy(&impl->e_vecs_in[i]));
40     CeedCallBackend(CeedVectorDestroy(&impl->q_vecs_in[i]));
41   }
42   CeedCallBackend(CeedFree(&impl->e_vecs_in));
43   CeedCallBackend(CeedFree(&impl->q_vecs_in));
44 
45   for (CeedInt i = 0; i < impl->num_outputs; i++) {
46     CeedCallBackend(CeedVectorDestroy(&impl->e_vecs_out[i]));
47     CeedCallBackend(CeedVectorDestroy(&impl->q_vecs_out[i]));
48   }
49   CeedCallBackend(CeedFree(&impl->e_vecs_out));
50   CeedCallBackend(CeedFree(&impl->q_vecs_out));
51   CeedCallBackend(CeedVectorDestroy(&impl->point_coords_elem));
52 
53   // QFunction assembly data
54   for (CeedInt i = 0; i < impl->num_active_in; i++) {
55     CeedCallBackend(CeedVectorDestroy(&impl->qf_active_in[i]));
56   }
57   CeedCallBackend(CeedFree(&impl->qf_active_in));
58 
59   // Diag data
60   if (impl->diag) {
61     Ceed ceed;
62 
63     CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
64     if (impl->diag->module) {
65       CeedCallCuda(ceed, cuModuleUnload(impl->diag->module));
66     }
67     if (impl->diag->module_point_block) {
68       CeedCallCuda(ceed, cuModuleUnload(impl->diag->module_point_block));
69     }
70     CeedCallCuda(ceed, cudaFree(impl->diag->d_eval_modes_in));
71     CeedCallCuda(ceed, cudaFree(impl->diag->d_eval_modes_out));
72     CeedCallCuda(ceed, cudaFree(impl->diag->d_identity));
73     CeedCallCuda(ceed, cudaFree(impl->diag->d_interp_in));
74     CeedCallCuda(ceed, cudaFree(impl->diag->d_interp_out));
75     CeedCallCuda(ceed, cudaFree(impl->diag->d_grad_in));
76     CeedCallCuda(ceed, cudaFree(impl->diag->d_grad_out));
77     CeedCallCuda(ceed, cudaFree(impl->diag->d_div_in));
78     CeedCallCuda(ceed, cudaFree(impl->diag->d_div_out));
79     CeedCallCuda(ceed, cudaFree(impl->diag->d_curl_in));
80     CeedCallCuda(ceed, cudaFree(impl->diag->d_curl_out));
81     CeedCallBackend(CeedDestroy(&ceed));
82     CeedCallBackend(CeedVectorDestroy(&impl->diag->elem_diag));
83     CeedCallBackend(CeedVectorDestroy(&impl->diag->point_block_elem_diag));
84     CeedCallBackend(CeedElemRestrictionDestroy(&impl->diag->diag_rstr));
85     CeedCallBackend(CeedElemRestrictionDestroy(&impl->diag->point_block_diag_rstr));
86   }
87   CeedCallBackend(CeedFree(&impl->diag));
88 
89   if (impl->asmb) {
90     Ceed ceed;
91 
92     CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
93     CeedCallCuda(ceed, cuModuleUnload(impl->asmb->module));
94     CeedCallCuda(ceed, cudaFree(impl->asmb->d_B_in));
95     CeedCallCuda(ceed, cudaFree(impl->asmb->d_B_out));
96     CeedCallBackend(CeedDestroy(&ceed));
97   }
98   CeedCallBackend(CeedFree(&impl->asmb));
99 
100   CeedCallBackend(CeedFree(&impl));
101   return CEED_ERROR_SUCCESS;
102 }
103 
104 //------------------------------------------------------------------------------
105 // Setup infields or outfields
106 //------------------------------------------------------------------------------
107 static int CeedOperatorSetupFields_Cuda(CeedQFunction qf, CeedOperator op, bool is_input, bool is_at_points, bool *skip_rstr, bool *apply_add_basis,
108                                         CeedVector *e_vecs, CeedVector *q_vecs, CeedInt num_fields, CeedInt Q, CeedInt num_elem) {
109   Ceed                ceed;
110   CeedQFunctionField *qf_fields;
111   CeedOperatorField  *op_fields;
112 
113   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
114   if (is_input) {
115     CeedCallBackend(CeedOperatorGetFields(op, NULL, &op_fields, NULL, NULL));
116     CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_fields, NULL, NULL));
117   } else {
118     CeedCallBackend(CeedOperatorGetFields(op, NULL, NULL, NULL, &op_fields));
119     CeedCallBackend(CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qf_fields));
120   }
121 
122   // Loop over fields
123   for (CeedInt i = 0; i < num_fields; i++) {
124     bool                is_active = false, is_strided = false, skip_e_vec = false;
125     CeedSize            q_size;
126     CeedInt             size;
127     CeedEvalMode        eval_mode;
128     CeedVector          l_vec;
129     CeedElemRestriction elem_rstr;
130 
131     // Check whether this field can skip the element restriction:
132     // Input CEED_VECTOR_ACTIVE
133     // Output CEED_VECTOR_ACTIVE without CEED_EVAL_NONE
134     // Input CEED_VECTOR_NONE with CEED_EVAL_WEIGHT
135     // Input passive vector with CEED_EVAL_NONE and strided restriction with CEED_STRIDES_BACKEND
136     CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &l_vec));
137     is_active = l_vec == CEED_VECTOR_ACTIVE;
138     CeedCallBackend(CeedVectorDestroy(&l_vec));
139     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_fields[i], &elem_rstr));
140     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode));
141     skip_e_vec = (is_input && is_active) || (is_active && eval_mode != CEED_EVAL_NONE) || (eval_mode == CEED_EVAL_WEIGHT);
142     if (!skip_e_vec && is_input && !is_active && eval_mode == CEED_EVAL_NONE) {
143       CeedCallBackend(CeedElemRestrictionIsStrided(elem_rstr, &is_strided));
144       if (is_strided) CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &skip_e_vec));
145     }
146     if (skip_e_vec) {
147       e_vecs[i] = NULL;
148     } else {
149       CeedCallBackend(CeedElemRestrictionCreateVector(elem_rstr, NULL, &e_vecs[i]));
150     }
151     CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
152 
153     switch (eval_mode) {
154       case CEED_EVAL_NONE:
155       case CEED_EVAL_INTERP:
156       case CEED_EVAL_GRAD:
157       case CEED_EVAL_DIV:
158       case CEED_EVAL_CURL:
159         CeedCallBackend(CeedQFunctionFieldGetSize(qf_fields[i], &size));
160         q_size = (CeedSize)num_elem * (CeedSize)Q * (CeedSize)size;
161         CeedCallBackend(CeedVectorCreate(ceed, q_size, &q_vecs[i]));
162         break;
163       case CEED_EVAL_WEIGHT: {
164         CeedBasis basis;
165 
166         CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis));
167         q_size = (CeedSize)num_elem * (CeedSize)Q;
168         CeedCallBackend(CeedVectorCreate(ceed, q_size, &q_vecs[i]));
169         if (is_at_points) {
170           CeedInt num_points[num_elem];
171 
172           for (CeedInt i = 0; i < num_elem; i++) num_points[i] = Q;
173           CeedCallBackend(
174               CeedBasisApplyAtPoints(basis, num_elem, num_points, CEED_NOTRANSPOSE, CEED_EVAL_WEIGHT, CEED_VECTOR_NONE, CEED_VECTOR_NONE, q_vecs[i]));
175         } else {
176           CeedCallBackend(CeedBasisApply(basis, num_elem, CEED_NOTRANSPOSE, CEED_EVAL_WEIGHT, CEED_VECTOR_NONE, q_vecs[i]));
177         }
178         CeedCallBackend(CeedBasisDestroy(&basis));
179         break;
180       }
181     }
182   }
183   // Drop duplicate restrictions
184   if (is_input) {
185     for (CeedInt i = 0; i < num_fields; i++) {
186       CeedVector          vec_i;
187       CeedElemRestriction rstr_i;
188 
189       CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec_i));
190       CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_fields[i], &rstr_i));
191       for (CeedInt j = i + 1; j < num_fields; j++) {
192         CeedVector          vec_j;
193         CeedElemRestriction rstr_j;
194 
195         CeedCallBackend(CeedOperatorFieldGetVector(op_fields[j], &vec_j));
196         CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_fields[j], &rstr_j));
197         if (vec_i == vec_j && rstr_i == rstr_j) {
198           if (e_vecs[i]) CeedCallBackend(CeedVectorReferenceCopy(e_vecs[i], &e_vecs[j]));
199           skip_rstr[j] = true;
200         }
201         CeedCallBackend(CeedVectorDestroy(&vec_j));
202         CeedCallBackend(CeedElemRestrictionDestroy(&rstr_j));
203       }
204       CeedCallBackend(CeedVectorDestroy(&vec_i));
205       CeedCallBackend(CeedElemRestrictionDestroy(&rstr_i));
206     }
207   } else {
208     for (CeedInt i = num_fields - 1; i >= 0; i--) {
209       CeedVector          vec_i;
210       CeedElemRestriction rstr_i;
211 
212       CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec_i));
213       CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_fields[i], &rstr_i));
214       for (CeedInt j = i - 1; j >= 0; j--) {
215         CeedVector          vec_j;
216         CeedElemRestriction rstr_j;
217 
218         CeedCallBackend(CeedOperatorFieldGetVector(op_fields[j], &vec_j));
219         CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_fields[j], &rstr_j));
220         if (vec_i == vec_j && rstr_i == rstr_j) {
221           if (e_vecs[i]) CeedCallBackend(CeedVectorReferenceCopy(e_vecs[i], &e_vecs[j]));
222           skip_rstr[j]       = true;
223           apply_add_basis[i] = true;
224         }
225         CeedCallBackend(CeedVectorDestroy(&vec_j));
226         CeedCallBackend(CeedElemRestrictionDestroy(&rstr_j));
227       }
228       CeedCallBackend(CeedVectorDestroy(&vec_i));
229       CeedCallBackend(CeedElemRestrictionDestroy(&rstr_i));
230     }
231   }
232   CeedCallBackend(CeedDestroy(&ceed));
233   return CEED_ERROR_SUCCESS;
234 }
235 
236 //------------------------------------------------------------------------------
237 // CeedOperator needs to connect all the named fields (be they active or passive) to the named inputs and outputs of its CeedQFunction.
238 //------------------------------------------------------------------------------
239 static int CeedOperatorSetup_Cuda(CeedOperator op) {
240   bool                is_setup_done;
241   CeedInt             Q, num_elem, num_input_fields, num_output_fields;
242   CeedQFunctionField *qf_input_fields, *qf_output_fields;
243   CeedQFunction       qf;
244   CeedOperatorField  *op_input_fields, *op_output_fields;
245   CeedOperator_Cuda  *impl;
246 
247   CeedCallBackend(CeedOperatorIsSetupDone(op, &is_setup_done));
248   if (is_setup_done) return CEED_ERROR_SUCCESS;
249 
250   CeedCallBackend(CeedOperatorGetData(op, &impl));
251   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
252   CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q));
253   CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem));
254   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
255   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
256 
257   // Allocate
258   CeedCallBackend(CeedCalloc(num_input_fields, &impl->e_vecs_in));
259   CeedCallBackend(CeedCalloc(num_output_fields, &impl->e_vecs_out));
260   CeedCallBackend(CeedCalloc(num_input_fields, &impl->skip_rstr_in));
261   CeedCallBackend(CeedCalloc(num_output_fields, &impl->skip_rstr_out));
262   CeedCallBackend(CeedCalloc(num_output_fields, &impl->apply_add_basis_out));
263   CeedCallBackend(CeedCalloc(num_input_fields, &impl->input_field_order));
264   CeedCallBackend(CeedCalloc(num_output_fields, &impl->output_field_order));
265   CeedCallBackend(CeedCalloc(num_input_fields, &impl->input_states));
266   CeedCallBackend(CeedCalloc(num_input_fields, &impl->q_vecs_in));
267   CeedCallBackend(CeedCalloc(num_output_fields, &impl->q_vecs_out));
268   impl->num_inputs  = num_input_fields;
269   impl->num_outputs = num_output_fields;
270 
271   // Set up infield and outfield e-vecs and q-vecs
272   CeedCallBackend(
273       CeedOperatorSetupFields_Cuda(qf, op, true, false, impl->skip_rstr_in, NULL, impl->e_vecs_in, impl->q_vecs_in, num_input_fields, Q, num_elem));
274   CeedCallBackend(CeedOperatorSetupFields_Cuda(qf, op, false, false, impl->skip_rstr_out, impl->apply_add_basis_out, impl->e_vecs_out,
275                                                impl->q_vecs_out, num_output_fields, Q, num_elem));
276 
277   // Reorder fields to allow reuse of buffers
278   impl->max_active_e_vec_len = 0;
279   {
280     bool    is_ordered[CEED_FIELD_MAX];
281     CeedInt curr_index = 0;
282 
283     for (CeedInt i = 0; i < num_input_fields; i++) is_ordered[i] = false;
284     for (CeedInt i = 0; i < num_input_fields; i++) {
285       CeedSize            e_vec_len_i;
286       CeedVector          vec_i;
287       CeedElemRestriction rstr_i;
288 
289       if (is_ordered[i]) continue;
290       is_ordered[i]                       = true;
291       impl->input_field_order[curr_index] = i;
292       curr_index++;
293       CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec_i));
294       if (vec_i == CEED_VECTOR_NONE) {
295         // CEED_EVAL_WEIGHT
296         CeedCallBackend(CeedVectorDestroy(&vec_i));
297         continue;
298       };
299       CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &rstr_i));
300       CeedCallBackend(CeedElemRestrictionGetEVectorSize(rstr_i, &e_vec_len_i));
301       impl->max_active_e_vec_len = e_vec_len_i > impl->max_active_e_vec_len ? e_vec_len_i : impl->max_active_e_vec_len;
302       for (CeedInt j = i + 1; j < num_input_fields; j++) {
303         CeedVector          vec_j;
304         CeedElemRestriction rstr_j;
305 
306         CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[j], &vec_j));
307         CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[j], &rstr_j));
308         if (rstr_i == rstr_j && vec_i == vec_j) {
309           is_ordered[j]                       = true;
310           impl->input_field_order[curr_index] = j;
311           curr_index++;
312         }
313         CeedCallBackend(CeedVectorDestroy(&vec_j));
314         CeedCallBackend(CeedElemRestrictionDestroy(&rstr_j));
315       }
316       CeedCallBackend(CeedVectorDestroy(&vec_i));
317       CeedCallBackend(CeedElemRestrictionDestroy(&rstr_i));
318     }
319   }
320   {
321     bool    is_ordered[CEED_FIELD_MAX];
322     CeedInt curr_index = 0;
323 
324     for (CeedInt i = 0; i < num_output_fields; i++) is_ordered[i] = false;
325     for (CeedInt i = 0; i < num_output_fields; i++) {
326       CeedSize            e_vec_len_i;
327       CeedVector          vec_i;
328       CeedElemRestriction rstr_i;
329 
330       if (is_ordered[i]) continue;
331       is_ordered[i]                        = true;
332       impl->output_field_order[curr_index] = i;
333       curr_index++;
334       CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[i], &vec_i));
335       CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &rstr_i));
336       CeedCallBackend(CeedElemRestrictionGetEVectorSize(rstr_i, &e_vec_len_i));
337       impl->max_active_e_vec_len = e_vec_len_i > impl->max_active_e_vec_len ? e_vec_len_i : impl->max_active_e_vec_len;
338       for (CeedInt j = i + 1; j < num_output_fields; j++) {
339         CeedVector          vec_j;
340         CeedElemRestriction rstr_j;
341 
342         CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[j], &vec_j));
343         CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[j], &rstr_j));
344         if (rstr_i == rstr_j && vec_i == vec_j) {
345           is_ordered[j]                        = true;
346           impl->output_field_order[curr_index] = j;
347           curr_index++;
348         }
349         CeedCallBackend(CeedVectorDestroy(&vec_j));
350         CeedCallBackend(CeedElemRestrictionDestroy(&rstr_j));
351       }
352       CeedCallBackend(CeedVectorDestroy(&vec_i));
353       CeedCallBackend(CeedElemRestrictionDestroy(&rstr_i));
354     }
355   }
356 
357   CeedCallBackend(CeedClearWorkVectors(CeedOperatorReturnCeed(op), impl->max_active_e_vec_len));
358   CeedCallBackend(CeedOperatorSetSetupDone(op));
359   CeedCallBackend(CeedQFunctionDestroy(&qf));
360   return CEED_ERROR_SUCCESS;
361 }
362 
363 //------------------------------------------------------------------------------
364 // Restrict Operator Inputs
365 //------------------------------------------------------------------------------
366 static inline int CeedOperatorInputRestrict_Cuda(CeedOperatorField op_input_field, CeedQFunctionField qf_input_field, CeedInt input_field,
367                                                  CeedVector in_vec, CeedVector active_e_vec, const bool skip_active, CeedOperator_Cuda *impl,
368                                                  CeedRequest *request) {
369   bool       is_active = false;
370   CeedVector l_vec, e_vec = impl->e_vecs_in[input_field];
371 
372   // Get input vector
373   CeedCallBackend(CeedOperatorFieldGetVector(op_input_field, &l_vec));
374   is_active = l_vec == CEED_VECTOR_ACTIVE;
375   if (is_active && skip_active) return CEED_ERROR_SUCCESS;
376   if (is_active) {
377     l_vec = in_vec;
378     if (!e_vec) e_vec = active_e_vec;
379   }
380 
381   // Restriction action
382   if (e_vec) {
383     // Restrict, if necessary
384     if (!impl->skip_rstr_in[input_field]) {
385       uint64_t state;
386 
387       CeedCallBackend(CeedVectorGetState(l_vec, &state));
388       if (is_active || state != impl->input_states[input_field]) {
389         CeedElemRestriction elem_rstr;
390 
391         CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_field, &elem_rstr));
392         CeedCallBackend(CeedElemRestrictionApply(elem_rstr, CEED_NOTRANSPOSE, l_vec, e_vec, request));
393         CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
394       }
395       impl->input_states[input_field] = state;
396     }
397   }
398   if (!is_active) CeedCallBackend(CeedVectorDestroy(&l_vec));
399   return CEED_ERROR_SUCCESS;
400 }
401 
402 //------------------------------------------------------------------------------
403 // Input Basis Action
404 //------------------------------------------------------------------------------
405 static inline int CeedOperatorInputBasis_Cuda(CeedOperatorField op_input_field, CeedQFunctionField qf_input_field, CeedInt input_field,
406                                               CeedVector in_vec, CeedVector active_e_vec, CeedInt num_elem, const bool skip_active,
407                                               CeedOperator_Cuda *impl) {
408   bool         is_active = false;
409   CeedEvalMode eval_mode;
410   CeedVector   l_vec, e_vec = impl->e_vecs_in[input_field], q_vec = impl->q_vecs_in[input_field];
411 
412   // Skip active input
413   CeedCallBackend(CeedOperatorFieldGetVector(op_input_field, &l_vec));
414   is_active = l_vec == CEED_VECTOR_ACTIVE;
415   if (is_active && skip_active) return CEED_ERROR_SUCCESS;
416   if (is_active) {
417     l_vec = in_vec;
418     if (!e_vec) e_vec = active_e_vec;
419   }
420 
421   // Basis action
422   CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_field, &eval_mode));
423   switch (eval_mode) {
424     case CEED_EVAL_NONE: {
425       const CeedScalar *e_vec_array;
426 
427       if (e_vec) {
428         CeedCallBackend(CeedVectorGetArrayRead(e_vec, CEED_MEM_DEVICE, &e_vec_array));
429       } else {
430         CeedCallBackend(CeedVectorGetArrayRead(l_vec, CEED_MEM_DEVICE, &e_vec_array));
431       }
432       CeedCallBackend(CeedVectorSetArray(q_vec, CEED_MEM_DEVICE, CEED_USE_POINTER, (CeedScalar *)e_vec_array));
433       break;
434     }
435     case CEED_EVAL_INTERP:
436     case CEED_EVAL_GRAD:
437     case CEED_EVAL_DIV:
438     case CEED_EVAL_CURL: {
439       CeedBasis basis;
440 
441       CeedCallBackend(CeedOperatorFieldGetBasis(op_input_field, &basis));
442       CeedCallBackend(CeedBasisApply(basis, num_elem, CEED_NOTRANSPOSE, eval_mode, e_vec, q_vec));
443       CeedCallBackend(CeedBasisDestroy(&basis));
444       break;
445     }
446     case CEED_EVAL_WEIGHT:
447       break;  // No action
448   }
449   if (!is_active) CeedCallBackend(CeedVectorDestroy(&l_vec));
450   return CEED_ERROR_SUCCESS;
451 }
452 
453 //------------------------------------------------------------------------------
454 // Restore Input Vectors
455 //------------------------------------------------------------------------------
456 static inline int CeedOperatorInputRestore_Cuda(CeedOperatorField op_input_field, CeedQFunctionField qf_input_field, CeedInt input_field,
457                                                 CeedVector in_vec, CeedVector active_e_vec, const bool skip_active, CeedOperator_Cuda *impl) {
458   bool         is_active = false;
459   CeedEvalMode eval_mode;
460   CeedVector   l_vec, e_vec = impl->e_vecs_in[input_field];
461 
462   // Skip active input
463   CeedCallBackend(CeedOperatorFieldGetVector(op_input_field, &l_vec));
464   is_active = l_vec == CEED_VECTOR_ACTIVE;
465   if (is_active && skip_active) return CEED_ERROR_SUCCESS;
466   if (is_active) {
467     l_vec = in_vec;
468     if (!e_vec) e_vec = active_e_vec;
469   }
470 
471   // Restore e-vec
472   CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_field, &eval_mode));
473   if (eval_mode == CEED_EVAL_NONE) {
474     const CeedScalar *e_vec_array;
475 
476     CeedCallBackend(CeedVectorTakeArray(impl->q_vecs_in[input_field], CEED_MEM_DEVICE, (CeedScalar **)&e_vec_array));
477     if (e_vec) {
478       CeedCallBackend(CeedVectorRestoreArrayRead(e_vec, &e_vec_array));
479     } else {
480       CeedCallBackend(CeedVectorRestoreArrayRead(l_vec, &e_vec_array));
481     }
482   }
483   if (!is_active) CeedCallBackend(CeedVectorDestroy(&l_vec));
484   return CEED_ERROR_SUCCESS;
485 }
486 
487 //------------------------------------------------------------------------------
488 // Apply and add to output
489 //------------------------------------------------------------------------------
490 static int CeedOperatorApplyAdd_Cuda(CeedOperator op, CeedVector in_vec, CeedVector out_vec, CeedRequest *request) {
491   CeedInt             Q, num_elem, num_input_fields, num_output_fields;
492   Ceed                ceed;
493   CeedVector          active_e_vec;
494   CeedQFunctionField *qf_input_fields, *qf_output_fields;
495   CeedQFunction       qf;
496   CeedOperatorField  *op_input_fields, *op_output_fields;
497   CeedOperator_Cuda  *impl;
498 
499   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
500   CeedCallBackend(CeedOperatorGetData(op, &impl));
501   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
502   CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q));
503   CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem));
504   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
505   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
506 
507   // Setup
508   CeedCallBackend(CeedOperatorSetup_Cuda(op));
509 
510   // Work vector
511   CeedCallBackend(CeedGetWorkVector(ceed, impl->max_active_e_vec_len, &active_e_vec));
512 
513   // Process inputs
514   for (CeedInt i = 0; i < num_input_fields; i++) {
515     CeedInt field = impl->input_field_order[i];
516 
517     CeedCallBackend(
518         CeedOperatorInputRestrict_Cuda(op_input_fields[field], qf_input_fields[field], field, in_vec, active_e_vec, false, impl, request));
519     CeedCallBackend(CeedOperatorInputBasis_Cuda(op_input_fields[field], qf_input_fields[field], field, in_vec, active_e_vec, num_elem, false, impl));
520   }
521 
522   // Output pointers, as necessary
523   for (CeedInt i = 0; i < num_output_fields; i++) {
524     CeedEvalMode eval_mode;
525 
526     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
527     if (eval_mode == CEED_EVAL_NONE) {
528       CeedScalar *e_vec_array;
529 
530       CeedCallBackend(CeedVectorGetArrayWrite(impl->e_vecs_out[i], CEED_MEM_DEVICE, &e_vec_array));
531       CeedCallBackend(CeedVectorSetArray(impl->q_vecs_out[i], CEED_MEM_DEVICE, CEED_USE_POINTER, e_vec_array));
532     }
533   }
534 
535   // Q function
536   CeedCallBackend(CeedQFunctionApply(qf, num_elem * Q, impl->q_vecs_in, impl->q_vecs_out));
537 
538   // Restore input arrays
539   for (CeedInt i = 0; i < num_input_fields; i++) {
540     CeedCallBackend(CeedOperatorInputRestore_Cuda(op_input_fields[i], qf_input_fields[i], i, in_vec, active_e_vec, false, impl));
541   }
542 
543   // Output basis and restriction
544   for (CeedInt i = 0; i < num_output_fields; i++) {
545     bool         is_active = false;
546     CeedInt      field     = impl->output_field_order[i];
547     CeedEvalMode eval_mode;
548     CeedVector   l_vec, e_vec = impl->e_vecs_out[field], q_vec = impl->q_vecs_out[field];
549 
550     // Output vector
551     CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[field], &l_vec));
552     is_active = l_vec == CEED_VECTOR_ACTIVE;
553     if (is_active) {
554       l_vec = out_vec;
555       if (!e_vec) e_vec = active_e_vec;
556     }
557 
558     // Basis action
559     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[field], &eval_mode));
560     switch (eval_mode) {
561       case CEED_EVAL_NONE:
562         break;  // No action
563       case CEED_EVAL_INTERP:
564       case CEED_EVAL_GRAD:
565       case CEED_EVAL_DIV:
566       case CEED_EVAL_CURL: {
567         CeedBasis basis;
568 
569         CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[field], &basis));
570         if (impl->apply_add_basis_out[field]) {
571           CeedCallBackend(CeedBasisApplyAdd(basis, num_elem, CEED_TRANSPOSE, eval_mode, q_vec, e_vec));
572         } else {
573           CeedCallBackend(CeedBasisApply(basis, num_elem, CEED_TRANSPOSE, eval_mode, q_vec, e_vec));
574         }
575         CeedCallBackend(CeedBasisDestroy(&basis));
576         break;
577       }
578       // LCOV_EXCL_START
579       case CEED_EVAL_WEIGHT: {
580         return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_WEIGHT cannot be an output evaluation mode");
581         // LCOV_EXCL_STOP
582       }
583     }
584 
585     // Restore evec
586     if (eval_mode == CEED_EVAL_NONE) {
587       CeedScalar *e_vec_array;
588 
589       CeedCallBackend(CeedVectorTakeArray(impl->q_vecs_out[i], CEED_MEM_DEVICE, &e_vec_array));
590       CeedCallBackend(CeedVectorRestoreArray(e_vec, &e_vec_array));
591     }
592 
593     // Restrict
594     if (!impl->skip_rstr_out[field]) {
595       CeedElemRestriction elem_rstr;
596 
597       CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[field], &elem_rstr));
598       CeedCallBackend(CeedElemRestrictionApply(elem_rstr, CEED_TRANSPOSE, e_vec, l_vec, request));
599       CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
600     }
601     if (!is_active) CeedCallBackend(CeedVectorDestroy(&l_vec));
602   }
603 
604   // Return work vector
605   CeedCallBackend(CeedRestoreWorkVector(ceed, &active_e_vec));
606   CeedCallBackend(CeedDestroy(&ceed));
607   CeedCallBackend(CeedQFunctionDestroy(&qf));
608   return CEED_ERROR_SUCCESS;
609 }
610 
611 //------------------------------------------------------------------------------
612 // CeedOperator needs to connect all the named fields (be they active or passive) to the named inputs and outputs of its CeedQFunction.
613 //------------------------------------------------------------------------------
614 static int CeedOperatorSetupAtPoints_Cuda(CeedOperator op) {
615   bool                is_setup_done;
616   CeedInt             max_num_points = -1, num_elem, num_input_fields, num_output_fields;
617   CeedQFunctionField *qf_input_fields, *qf_output_fields;
618   CeedQFunction       qf;
619   CeedOperatorField  *op_input_fields, *op_output_fields;
620   CeedOperator_Cuda  *impl;
621 
622   CeedCallBackend(CeedOperatorIsSetupDone(op, &is_setup_done));
623   if (is_setup_done) return CEED_ERROR_SUCCESS;
624 
625   CeedCallBackend(CeedOperatorGetData(op, &impl));
626   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
627   CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem));
628   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
629   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
630   {
631     CeedElemRestriction rstr_points = NULL;
632 
633     CeedCallBackend(CeedOperatorAtPointsGetPoints(op, &rstr_points, NULL));
634     CeedCallBackend(CeedElemRestrictionGetMaxPointsInElement(rstr_points, &max_num_points));
635     CeedCallBackend(CeedCalloc(num_elem, &impl->num_points));
636     for (CeedInt e = 0; e < num_elem; e++) {
637       CeedInt num_points_elem;
638 
639       CeedCallBackend(CeedElemRestrictionGetNumPointsInElement(rstr_points, e, &num_points_elem));
640       impl->num_points[e] = num_points_elem;
641     }
642     CeedCallBackend(CeedElemRestrictionDestroy(&rstr_points));
643   }
644   impl->max_num_points = max_num_points;
645 
646   // Allocate
647   CeedCallBackend(CeedCalloc(num_input_fields, &impl->e_vecs_in));
648   CeedCallBackend(CeedCalloc(num_output_fields, &impl->e_vecs_out));
649   CeedCallBackend(CeedCalloc(num_input_fields, &impl->skip_rstr_in));
650   CeedCallBackend(CeedCalloc(num_output_fields, &impl->skip_rstr_out));
651   CeedCallBackend(CeedCalloc(num_output_fields, &impl->apply_add_basis_out));
652   CeedCallBackend(CeedCalloc(num_input_fields, &impl->input_field_order));
653   CeedCallBackend(CeedCalloc(num_output_fields, &impl->output_field_order));
654   CeedCallBackend(CeedCalloc(num_input_fields, &impl->input_states));
655   CeedCallBackend(CeedCalloc(num_input_fields, &impl->q_vecs_in));
656   CeedCallBackend(CeedCalloc(num_output_fields, &impl->q_vecs_out));
657   impl->num_inputs  = num_input_fields;
658   impl->num_outputs = num_output_fields;
659 
660   // Set up infield and outfield e-vecs and q-vecs
661   CeedCallBackend(CeedOperatorSetupFields_Cuda(qf, op, true, true, impl->skip_rstr_in, NULL, impl->e_vecs_in, impl->q_vecs_in, num_input_fields,
662                                                max_num_points, num_elem));
663   CeedCallBackend(CeedOperatorSetupFields_Cuda(qf, op, false, true, impl->skip_rstr_out, impl->apply_add_basis_out, impl->e_vecs_out,
664                                                impl->q_vecs_out, num_output_fields, max_num_points, num_elem));
665 
666   // Reorder fields to allow reuse of buffers
667   impl->max_active_e_vec_len = 0;
668   {
669     bool    is_ordered[CEED_FIELD_MAX];
670     CeedInt curr_index = 0;
671 
672     for (CeedInt i = 0; i < num_input_fields; i++) is_ordered[i] = false;
673     for (CeedInt i = 0; i < num_input_fields; i++) {
674       CeedSize            e_vec_len_i;
675       CeedVector          vec_i;
676       CeedElemRestriction rstr_i;
677 
678       if (is_ordered[i]) continue;
679       is_ordered[i]                       = true;
680       impl->input_field_order[curr_index] = i;
681       curr_index++;
682       CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec_i));
683       if (vec_i == CEED_VECTOR_NONE) {
684         // CEED_EVAL_WEIGHT
685         CeedCallBackend(CeedVectorDestroy(&vec_i));
686         continue;
687       };
688       CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &rstr_i));
689       CeedCallBackend(CeedElemRestrictionGetEVectorSize(rstr_i, &e_vec_len_i));
690       impl->max_active_e_vec_len = e_vec_len_i > impl->max_active_e_vec_len ? e_vec_len_i : impl->max_active_e_vec_len;
691       for (CeedInt j = i + 1; j < num_input_fields; j++) {
692         CeedVector          vec_j;
693         CeedElemRestriction rstr_j;
694 
695         CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[j], &vec_j));
696         CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[j], &rstr_j));
697         if (rstr_i == rstr_j && vec_i == vec_j) {
698           is_ordered[j]                       = true;
699           impl->input_field_order[curr_index] = j;
700           curr_index++;
701         }
702         CeedCallBackend(CeedVectorDestroy(&vec_j));
703         CeedCallBackend(CeedElemRestrictionDestroy(&rstr_j));
704       }
705       CeedCallBackend(CeedVectorDestroy(&vec_i));
706       CeedCallBackend(CeedElemRestrictionDestroy(&rstr_i));
707     }
708   }
709   {
710     bool    is_ordered[CEED_FIELD_MAX];
711     CeedInt curr_index = 0;
712 
713     for (CeedInt i = 0; i < num_output_fields; i++) is_ordered[i] = false;
714     for (CeedInt i = 0; i < num_output_fields; i++) {
715       CeedSize            e_vec_len_i;
716       CeedVector          vec_i;
717       CeedElemRestriction rstr_i;
718 
719       if (is_ordered[i]) continue;
720       is_ordered[i]                        = true;
721       impl->output_field_order[curr_index] = i;
722       curr_index++;
723       CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[i], &vec_i));
724       CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &rstr_i));
725       CeedCallBackend(CeedElemRestrictionGetEVectorSize(rstr_i, &e_vec_len_i));
726       impl->max_active_e_vec_len = e_vec_len_i > impl->max_active_e_vec_len ? e_vec_len_i : impl->max_active_e_vec_len;
727       for (CeedInt j = i + 1; j < num_output_fields; j++) {
728         CeedVector          vec_j;
729         CeedElemRestriction rstr_j;
730 
731         CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[j], &vec_j));
732         CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[j], &rstr_j));
733         if (rstr_i == rstr_j && vec_i == vec_j) {
734           is_ordered[j]                        = true;
735           impl->output_field_order[curr_index] = j;
736           curr_index++;
737         }
738         CeedCallBackend(CeedVectorDestroy(&vec_j));
739         CeedCallBackend(CeedElemRestrictionDestroy(&rstr_j));
740       }
741       CeedCallBackend(CeedVectorDestroy(&vec_i));
742       CeedCallBackend(CeedElemRestrictionDestroy(&rstr_i));
743     }
744   }
745 
746   CeedCallBackend(CeedClearWorkVectors(CeedOperatorReturnCeed(op), impl->max_active_e_vec_len));
747 
748   CeedCallBackend(CeedOperatorSetSetupDone(op));
749   CeedCallBackend(CeedQFunctionDestroy(&qf));
750   return CEED_ERROR_SUCCESS;
751 }
752 
753 //------------------------------------------------------------------------------
754 // Input Basis Action AtPoints
755 //------------------------------------------------------------------------------
756 static inline int CeedOperatorInputBasisAtPoints_Cuda(CeedOperatorField op_input_field, CeedQFunctionField qf_input_field, CeedInt input_field,
757                                                       CeedVector in_vec, CeedVector active_e_vec, CeedInt num_elem, const CeedInt *num_points,
758                                                       const bool skip_active, const bool skip_passive, CeedOperator_Cuda *impl) {
759   bool         is_active = false;
760   CeedEvalMode eval_mode;
761   CeedVector   l_vec, e_vec = impl->e_vecs_in[input_field], q_vec = impl->q_vecs_in[input_field];
762 
763   // Skip active input
764   CeedCallBackend(CeedOperatorFieldGetVector(op_input_field, &l_vec));
765   is_active = l_vec == CEED_VECTOR_ACTIVE;
766   if (skip_active && is_active) return CEED_ERROR_SUCCESS;
767   if (skip_passive && !is_active) {
768     CeedCallBackend(CeedVectorDestroy(&l_vec));
769     return CEED_ERROR_SUCCESS;
770   }
771   if (is_active) {
772     l_vec = in_vec;
773     if (!e_vec) e_vec = active_e_vec;
774   }
775 
776   // Basis action
777   CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_field, &eval_mode));
778   switch (eval_mode) {
779     case CEED_EVAL_NONE: {
780       const CeedScalar *e_vec_array;
781 
782       if (e_vec) {
783         CeedCallBackend(CeedVectorGetArrayRead(e_vec, CEED_MEM_DEVICE, &e_vec_array));
784       } else {
785         CeedCallBackend(CeedVectorGetArrayRead(l_vec, CEED_MEM_DEVICE, &e_vec_array));
786       }
787       CeedCallBackend(CeedVectorSetArray(q_vec, CEED_MEM_DEVICE, CEED_USE_POINTER, (CeedScalar *)e_vec_array));
788       break;
789     }
790     case CEED_EVAL_INTERP:
791     case CEED_EVAL_GRAD:
792     case CEED_EVAL_DIV:
793     case CEED_EVAL_CURL: {
794       CeedBasis basis;
795 
796       CeedCallBackend(CeedOperatorFieldGetBasis(op_input_field, &basis));
797       CeedCallBackend(CeedBasisApplyAtPoints(basis, num_elem, num_points, CEED_NOTRANSPOSE, eval_mode, impl->point_coords_elem, e_vec, q_vec));
798       CeedCallBackend(CeedBasisDestroy(&basis));
799       break;
800     }
801     case CEED_EVAL_WEIGHT:
802       break;  // No action
803   }
804   if (!is_active) CeedCallBackend(CeedVectorDestroy(&l_vec));
805   return CEED_ERROR_SUCCESS;
806 }
807 
808 //------------------------------------------------------------------------------
809 // Apply and add to output AtPoints
810 //------------------------------------------------------------------------------
811 static int CeedOperatorApplyAddAtPoints_Cuda(CeedOperator op, CeedVector in_vec, CeedVector out_vec, CeedRequest *request) {
812   CeedInt             max_num_points, *num_points, num_elem, num_input_fields, num_output_fields;
813   Ceed                ceed;
814   CeedVector          active_e_vec;
815   CeedQFunctionField *qf_input_fields, *qf_output_fields;
816   CeedQFunction       qf;
817   CeedOperatorField  *op_input_fields, *op_output_fields;
818   CeedOperator_Cuda  *impl;
819 
820   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
821   CeedCallBackend(CeedOperatorGetData(op, &impl));
822   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
823   CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem));
824   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
825   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
826 
827   // Setup
828   CeedCallBackend(CeedOperatorSetupAtPoints_Cuda(op));
829   num_points     = impl->num_points;
830   max_num_points = impl->max_num_points;
831 
832   // Work vector
833   CeedCallBackend(CeedGetWorkVector(ceed, impl->max_active_e_vec_len, &active_e_vec));
834 
835   // Get point coordinates
836   if (!impl->point_coords_elem) {
837     CeedVector          point_coords = NULL;
838     CeedElemRestriction rstr_points  = NULL;
839 
840     CeedCallBackend(CeedOperatorAtPointsGetPoints(op, &rstr_points, &point_coords));
841     CeedCallBackend(CeedElemRestrictionCreateVector(rstr_points, NULL, &impl->point_coords_elem));
842     CeedCallBackend(CeedElemRestrictionApply(rstr_points, CEED_NOTRANSPOSE, point_coords, impl->point_coords_elem, request));
843     CeedCallBackend(CeedVectorDestroy(&point_coords));
844     CeedCallBackend(CeedElemRestrictionDestroy(&rstr_points));
845   }
846 
847   // Process inputs
848   for (CeedInt i = 0; i < num_input_fields; i++) {
849     CeedInt field = impl->input_field_order[i];
850 
851     CeedCallBackend(
852         CeedOperatorInputRestrict_Cuda(op_input_fields[field], qf_input_fields[field], field, in_vec, active_e_vec, false, impl, request));
853     CeedCallBackend(CeedOperatorInputBasisAtPoints_Cuda(op_input_fields[field], qf_input_fields[field], field, in_vec, active_e_vec, num_elem,
854                                                         num_points, false, false, impl));
855   }
856 
857   // Output pointers, as necessary
858   for (CeedInt i = 0; i < num_output_fields; i++) {
859     CeedEvalMode eval_mode;
860 
861     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
862     if (eval_mode == CEED_EVAL_NONE) {
863       CeedScalar *e_vec_array;
864 
865       CeedCallBackend(CeedVectorGetArrayWrite(impl->e_vecs_out[i], CEED_MEM_DEVICE, &e_vec_array));
866       CeedCallBackend(CeedVectorSetArray(impl->q_vecs_out[i], CEED_MEM_DEVICE, CEED_USE_POINTER, e_vec_array));
867     }
868   }
869 
870   // Q function
871   CeedCallBackend(CeedQFunctionApply(qf, num_elem * max_num_points, impl->q_vecs_in, impl->q_vecs_out));
872 
873   // Restore input arrays
874   for (CeedInt i = 0; i < num_input_fields; i++) {
875     CeedCallBackend(CeedOperatorInputRestore_Cuda(op_input_fields[i], qf_input_fields[i], i, in_vec, active_e_vec, false, impl));
876   }
877 
878   // Output basis and restriction
879   for (CeedInt i = 0; i < num_output_fields; i++) {
880     bool         is_active = false;
881     CeedInt      field     = impl->output_field_order[i];
882     CeedEvalMode eval_mode;
883     CeedVector   l_vec, e_vec = impl->e_vecs_out[field], q_vec = impl->q_vecs_out[field];
884 
885     // Output vector
886     CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[field], &l_vec));
887     is_active = l_vec == CEED_VECTOR_ACTIVE;
888     if (is_active) {
889       l_vec = out_vec;
890       if (!e_vec) e_vec = active_e_vec;
891     }
892 
893     // Basis action
894     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[field], &eval_mode));
895     switch (eval_mode) {
896       case CEED_EVAL_NONE:
897         break;  // No action
898       case CEED_EVAL_INTERP:
899       case CEED_EVAL_GRAD:
900       case CEED_EVAL_DIV:
901       case CEED_EVAL_CURL: {
902         CeedBasis basis;
903 
904         CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[field], &basis));
905         if (impl->apply_add_basis_out[field]) {
906           CeedCallBackend(CeedBasisApplyAddAtPoints(basis, num_elem, num_points, CEED_TRANSPOSE, eval_mode, impl->point_coords_elem, q_vec, e_vec));
907         } else {
908           CeedCallBackend(CeedBasisApplyAtPoints(basis, num_elem, num_points, CEED_TRANSPOSE, eval_mode, impl->point_coords_elem, q_vec, e_vec));
909         }
910         CeedCallBackend(CeedBasisDestroy(&basis));
911         break;
912       }
913       // LCOV_EXCL_START
914       case CEED_EVAL_WEIGHT: {
915         return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_WEIGHT cannot be an output evaluation mode");
916         // LCOV_EXCL_STOP
917       }
918     }
919 
920     // Restore evec
921     if (eval_mode == CEED_EVAL_NONE) {
922       CeedScalar *e_vec_array;
923 
924       CeedCallBackend(CeedVectorTakeArray(impl->q_vecs_out[i], CEED_MEM_DEVICE, &e_vec_array));
925       CeedCallBackend(CeedVectorRestoreArray(e_vec, &e_vec_array));
926     }
927 
928     // Restrict
929     if (!impl->skip_rstr_out[field]) {
930       CeedElemRestriction elem_rstr;
931 
932       CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[field], &elem_rstr));
933       CeedCallBackend(CeedElemRestrictionApply(elem_rstr, CEED_TRANSPOSE, e_vec, l_vec, request));
934       CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
935     }
936     if (!is_active) CeedCallBackend(CeedVectorDestroy(&l_vec));
937   }
938 
939   // Restore work vector
940   CeedCallBackend(CeedRestoreWorkVector(ceed, &active_e_vec));
941   CeedCallBackend(CeedDestroy(&ceed));
942   CeedCallBackend(CeedQFunctionDestroy(&qf));
943   return CEED_ERROR_SUCCESS;
944 }
945 
946 //------------------------------------------------------------------------------
947 // Linear QFunction Assembly Core
948 //------------------------------------------------------------------------------
949 static inline int CeedOperatorLinearAssembleQFunctionCore_Cuda(CeedOperator op, bool build_objects, CeedVector *assembled, CeedElemRestriction *rstr,
950                                                                CeedRequest *request) {
951   Ceed                ceed, ceed_parent;
952   CeedInt             num_active_in, num_active_out, Q, num_elem, num_input_fields, num_output_fields, size;
953   CeedScalar         *assembled_array;
954   CeedVector         *active_inputs;
955   CeedQFunctionField *qf_input_fields, *qf_output_fields;
956   CeedQFunction       qf;
957   CeedOperatorField  *op_input_fields, *op_output_fields;
958   CeedOperator_Cuda  *impl;
959 
960   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
961   CeedCallBackend(CeedOperatorGetFallbackParentCeed(op, &ceed_parent));
962   CeedCallBackend(CeedOperatorGetData(op, &impl));
963   CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q));
964   CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem));
965   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
966   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
967   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
968   active_inputs = impl->qf_active_in;
969   num_active_in = impl->num_active_in, num_active_out = impl->num_active_out;
970 
971   // Setup
972   CeedCallBackend(CeedOperatorSetup_Cuda(op));
973 
974   // Process inputs
975   for (CeedInt i = 0; i < num_input_fields; i++) {
976     CeedCallBackend(CeedOperatorInputRestrict_Cuda(op_input_fields[i], qf_input_fields[i], i, NULL, NULL, true, impl, request));
977     CeedCallBackend(CeedOperatorInputBasis_Cuda(op_input_fields[i], qf_input_fields[i], i, NULL, NULL, num_elem, true, impl));
978   }
979 
980   // Count number of active input fields
981   if (!num_active_in) {
982     for (CeedInt i = 0; i < num_input_fields; i++) {
983       CeedScalar *q_vec_array;
984       CeedVector  l_vec;
985 
986       // Check if active input
987       CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &l_vec));
988       if (l_vec == CEED_VECTOR_ACTIVE) {
989         CeedCallBackend(CeedQFunctionFieldGetSize(qf_input_fields[i], &size));
990         CeedCallBackend(CeedVectorSetValue(impl->q_vecs_in[i], 0.0));
991         CeedCallBackend(CeedVectorGetArray(impl->q_vecs_in[i], CEED_MEM_DEVICE, &q_vec_array));
992         CeedCallBackend(CeedRealloc(num_active_in + size, &active_inputs));
993         for (CeedInt field = 0; field < size; field++) {
994           CeedSize q_size = (CeedSize)Q * num_elem;
995 
996           CeedCallBackend(CeedVectorCreate(ceed, q_size, &active_inputs[num_active_in + field]));
997           CeedCallBackend(
998               CeedVectorSetArray(active_inputs[num_active_in + field], CEED_MEM_DEVICE, CEED_USE_POINTER, &q_vec_array[field * Q * num_elem]));
999         }
1000         num_active_in += size;
1001         CeedCallBackend(CeedVectorRestoreArray(impl->q_vecs_in[i], &q_vec_array));
1002       }
1003       CeedCallBackend(CeedVectorDestroy(&l_vec));
1004     }
1005     impl->num_active_in = num_active_in;
1006     impl->qf_active_in  = active_inputs;
1007   }
1008 
1009   // Count number of active output fields
1010   if (!num_active_out) {
1011     for (CeedInt i = 0; i < num_output_fields; i++) {
1012       CeedVector l_vec;
1013 
1014       // Check if active output
1015       CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[i], &l_vec));
1016       if (l_vec == CEED_VECTOR_ACTIVE) {
1017         CeedCallBackend(CeedQFunctionFieldGetSize(qf_output_fields[i], &size));
1018         num_active_out += size;
1019       }
1020       CeedCallBackend(CeedVectorDestroy(&l_vec));
1021     }
1022     impl->num_active_out = num_active_out;
1023   }
1024 
1025   // Check sizes
1026   CeedCheck(num_active_in > 0 && num_active_out > 0, ceed, CEED_ERROR_BACKEND, "Cannot assemble QFunction without active inputs and outputs");
1027 
1028   // Build objects if needed
1029   if (build_objects) {
1030     CeedSize l_size     = (CeedSize)num_elem * Q * num_active_in * num_active_out;
1031     CeedInt  strides[3] = {1, num_elem * Q, Q}; /* *NOPAD* */
1032 
1033     // Create output restriction
1034     CeedCallBackend(CeedElemRestrictionCreateStrided(ceed_parent, num_elem, Q, num_active_in * num_active_out,
1035                                                      (CeedSize)num_active_in * (CeedSize)num_active_out * (CeedSize)num_elem * (CeedSize)Q, strides,
1036                                                      rstr));
1037     // Create assembled vector
1038     CeedCallBackend(CeedVectorCreate(ceed_parent, l_size, assembled));
1039   }
1040   CeedCallBackend(CeedVectorSetValue(*assembled, 0.0));
1041   CeedCallBackend(CeedVectorGetArray(*assembled, CEED_MEM_DEVICE, &assembled_array));
1042 
1043   // Assemble QFunction
1044   for (CeedInt in = 0; in < num_active_in; in++) {
1045     // Set Inputs
1046     CeedCallBackend(CeedVectorSetValue(active_inputs[in], 1.0));
1047     if (num_active_in > 1) {
1048       CeedCallBackend(CeedVectorSetValue(active_inputs[(in + num_active_in - 1) % num_active_in], 0.0));
1049     }
1050     // Set Outputs
1051     for (CeedInt out = 0; out < num_output_fields; out++) {
1052       CeedVector l_vec;
1053 
1054       // Check if active output
1055       CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[out], &l_vec));
1056       if (l_vec == CEED_VECTOR_ACTIVE) {
1057         CeedCallBackend(CeedVectorSetArray(impl->q_vecs_out[out], CEED_MEM_DEVICE, CEED_USE_POINTER, assembled_array));
1058         CeedCallBackend(CeedQFunctionFieldGetSize(qf_output_fields[out], &size));
1059         assembled_array += size * Q * num_elem;  // Advance the pointer by the size of the output
1060       }
1061       CeedCallBackend(CeedVectorDestroy(&l_vec));
1062     }
1063     // Apply QFunction
1064     CeedCallBackend(CeedQFunctionApply(qf, Q * num_elem, impl->q_vecs_in, impl->q_vecs_out));
1065   }
1066 
1067   // Un-set output q-vecs to prevent accidental overwrite of Assembled
1068   for (CeedInt out = 0; out < num_output_fields; out++) {
1069     CeedVector l_vec;
1070 
1071     CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[out], &l_vec));
1072     if (l_vec == CEED_VECTOR_ACTIVE) {
1073       CeedCallBackend(CeedVectorTakeArray(impl->q_vecs_out[out], CEED_MEM_DEVICE, NULL));
1074     }
1075     CeedCallBackend(CeedVectorDestroy(&l_vec));
1076   }
1077 
1078   // Restore input arrays
1079   for (CeedInt i = 0; i < num_input_fields; i++) {
1080     CeedCallBackend(CeedOperatorInputRestore_Cuda(op_input_fields[i], qf_input_fields[i], i, NULL, NULL, true, impl));
1081   }
1082 
1083   // Restore output
1084   CeedCallBackend(CeedVectorRestoreArray(*assembled, &assembled_array));
1085   CeedCallBackend(CeedDestroy(&ceed));
1086   CeedCallBackend(CeedDestroy(&ceed_parent));
1087   CeedCallBackend(CeedQFunctionDestroy(&qf));
1088   return CEED_ERROR_SUCCESS;
1089 }
1090 
1091 //------------------------------------------------------------------------------
1092 // Assemble Linear QFunction
1093 //------------------------------------------------------------------------------
1094 static int CeedOperatorLinearAssembleQFunction_Cuda(CeedOperator op, CeedVector *assembled, CeedElemRestriction *rstr, CeedRequest *request) {
1095   return CeedOperatorLinearAssembleQFunctionCore_Cuda(op, true, assembled, rstr, request);
1096 }
1097 
1098 //------------------------------------------------------------------------------
1099 // Update Assembled Linear QFunction
1100 //------------------------------------------------------------------------------
1101 static int CeedOperatorLinearAssembleQFunctionUpdate_Cuda(CeedOperator op, CeedVector assembled, CeedElemRestriction rstr, CeedRequest *request) {
1102   return CeedOperatorLinearAssembleQFunctionCore_Cuda(op, false, &assembled, &rstr, request);
1103 }
1104 
1105 //------------------------------------------------------------------------------
1106 // Assemble Diagonal Setup
1107 //------------------------------------------------------------------------------
1108 static inline int CeedOperatorAssembleDiagonalSetup_Cuda(CeedOperator op) {
1109   Ceed                ceed;
1110   CeedInt             num_input_fields, num_output_fields, num_eval_modes_in = 0, num_eval_modes_out = 0;
1111   CeedInt             q_comp, num_nodes, num_qpts;
1112   CeedEvalMode       *eval_modes_in = NULL, *eval_modes_out = NULL;
1113   CeedBasis           basis_in = NULL, basis_out = NULL;
1114   CeedQFunctionField *qf_fields;
1115   CeedQFunction       qf;
1116   CeedOperatorField  *op_fields;
1117   CeedOperator_Cuda  *impl;
1118 
1119   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
1120   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
1121   CeedCallBackend(CeedQFunctionGetNumArgs(qf, &num_input_fields, &num_output_fields));
1122 
1123   // Determine active input basis
1124   CeedCallBackend(CeedOperatorGetFields(op, NULL, &op_fields, NULL, NULL));
1125   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_fields, NULL, NULL));
1126   for (CeedInt i = 0; i < num_input_fields; i++) {
1127     CeedVector vec;
1128 
1129     CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec));
1130     if (vec == CEED_VECTOR_ACTIVE) {
1131       CeedEvalMode eval_mode;
1132       CeedBasis    basis;
1133 
1134       CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis));
1135       CeedCheck(!basis_in || basis_in == basis, ceed, CEED_ERROR_BACKEND,
1136                 "Backend does not implement operator diagonal assembly with multiple active bases");
1137       if (!basis_in) CeedCallBackend(CeedBasisReferenceCopy(basis, &basis_in));
1138       CeedCallBackend(CeedBasisDestroy(&basis));
1139       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode));
1140       CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_in, eval_mode, &q_comp));
1141       if (eval_mode != CEED_EVAL_WEIGHT) {
1142         // q_comp = 1 if CEED_EVAL_NONE, CEED_EVAL_WEIGHT caught by QF assembly
1143         CeedCallBackend(CeedRealloc(num_eval_modes_in + q_comp, &eval_modes_in));
1144         for (CeedInt d = 0; d < q_comp; d++) eval_modes_in[num_eval_modes_in + d] = eval_mode;
1145         num_eval_modes_in += q_comp;
1146       }
1147     }
1148     CeedCallBackend(CeedVectorDestroy(&vec));
1149   }
1150 
1151   // Determine active output basis
1152   CeedCallBackend(CeedOperatorGetFields(op, NULL, NULL, NULL, &op_fields));
1153   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qf_fields));
1154   for (CeedInt i = 0; i < num_output_fields; i++) {
1155     CeedVector vec;
1156 
1157     CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec));
1158     if (vec == CEED_VECTOR_ACTIVE) {
1159       CeedBasis    basis;
1160       CeedEvalMode eval_mode;
1161 
1162       CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis));
1163       CeedCheck(!basis_out || basis_out == basis, ceed, CEED_ERROR_BACKEND,
1164                 "Backend does not implement operator diagonal assembly with multiple active bases");
1165       if (!basis_out) CeedCallBackend(CeedBasisReferenceCopy(basis, &basis_out));
1166       CeedCallBackend(CeedBasisDestroy(&basis));
1167       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode));
1168       CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_out, eval_mode, &q_comp));
1169       if (eval_mode != CEED_EVAL_WEIGHT) {
1170         // q_comp = 1 if CEED_EVAL_NONE, CEED_EVAL_WEIGHT caught by QF assembly
1171         CeedCallBackend(CeedRealloc(num_eval_modes_out + q_comp, &eval_modes_out));
1172         for (CeedInt d = 0; d < q_comp; d++) eval_modes_out[num_eval_modes_out + d] = eval_mode;
1173         num_eval_modes_out += q_comp;
1174       }
1175     }
1176     CeedCallBackend(CeedVectorDestroy(&vec));
1177   }
1178 
1179   // Operator data struct
1180   CeedCallBackend(CeedOperatorGetData(op, &impl));
1181   CeedCallBackend(CeedCalloc(1, &impl->diag));
1182   CeedOperatorDiag_Cuda *diag = impl->diag;
1183 
1184   // Basis matrices
1185   CeedCallBackend(CeedBasisGetNumNodes(basis_in, &num_nodes));
1186   if (basis_in == CEED_BASIS_NONE) num_qpts = num_nodes;
1187   else CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis_in, &num_qpts));
1188   const CeedInt interp_bytes     = num_nodes * num_qpts * sizeof(CeedScalar);
1189   const CeedInt eval_modes_bytes = sizeof(CeedEvalMode);
1190   bool          has_eval_none    = false;
1191 
1192   // CEED_EVAL_NONE
1193   for (CeedInt i = 0; i < num_eval_modes_in; i++) has_eval_none = has_eval_none || (eval_modes_in[i] == CEED_EVAL_NONE);
1194   for (CeedInt i = 0; i < num_eval_modes_out; i++) has_eval_none = has_eval_none || (eval_modes_out[i] == CEED_EVAL_NONE);
1195   if (has_eval_none) {
1196     CeedScalar *identity = NULL;
1197 
1198     CeedCallBackend(CeedCalloc(num_nodes * num_qpts, &identity));
1199     for (CeedInt i = 0; i < (num_nodes < num_qpts ? num_nodes : num_qpts); i++) identity[i * num_nodes + i] = 1.0;
1200     CeedCallCuda(ceed, cudaMalloc((void **)&diag->d_identity, interp_bytes));
1201     CeedCallCuda(ceed, cudaMemcpy(diag->d_identity, identity, interp_bytes, cudaMemcpyHostToDevice));
1202     CeedCallBackend(CeedFree(&identity));
1203   }
1204 
1205   // CEED_EVAL_INTERP, CEED_EVAL_GRAD, CEED_EVAL_DIV, and CEED_EVAL_CURL
1206   for (CeedInt in = 0; in < 2; in++) {
1207     CeedFESpace fespace;
1208     CeedBasis   basis = in ? basis_in : basis_out;
1209 
1210     CeedCallBackend(CeedBasisGetFESpace(basis, &fespace));
1211     switch (fespace) {
1212       case CEED_FE_SPACE_H1: {
1213         CeedInt           q_comp_interp, q_comp_grad;
1214         const CeedScalar *interp, *grad;
1215         CeedScalar       *d_interp, *d_grad;
1216 
1217         CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp));
1218         CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_GRAD, &q_comp_grad));
1219 
1220         CeedCallBackend(CeedBasisGetInterp(basis, &interp));
1221         CeedCallCuda(ceed, cudaMalloc((void **)&d_interp, interp_bytes * q_comp_interp));
1222         CeedCallCuda(ceed, cudaMemcpy(d_interp, interp, interp_bytes * q_comp_interp, cudaMemcpyHostToDevice));
1223         CeedCallBackend(CeedBasisGetGrad(basis, &grad));
1224         CeedCallCuda(ceed, cudaMalloc((void **)&d_grad, interp_bytes * q_comp_grad));
1225         CeedCallCuda(ceed, cudaMemcpy(d_grad, grad, interp_bytes * q_comp_grad, cudaMemcpyHostToDevice));
1226         if (in) {
1227           diag->d_interp_in = d_interp;
1228           diag->d_grad_in   = d_grad;
1229         } else {
1230           diag->d_interp_out = d_interp;
1231           diag->d_grad_out   = d_grad;
1232         }
1233       } break;
1234       case CEED_FE_SPACE_HDIV: {
1235         CeedInt           q_comp_interp, q_comp_div;
1236         const CeedScalar *interp, *div;
1237         CeedScalar       *d_interp, *d_div;
1238 
1239         CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp));
1240         CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_DIV, &q_comp_div));
1241 
1242         CeedCallBackend(CeedBasisGetInterp(basis, &interp));
1243         CeedCallCuda(ceed, cudaMalloc((void **)&d_interp, interp_bytes * q_comp_interp));
1244         CeedCallCuda(ceed, cudaMemcpy(d_interp, interp, interp_bytes * q_comp_interp, cudaMemcpyHostToDevice));
1245         CeedCallBackend(CeedBasisGetDiv(basis, &div));
1246         CeedCallCuda(ceed, cudaMalloc((void **)&d_div, interp_bytes * q_comp_div));
1247         CeedCallCuda(ceed, cudaMemcpy(d_div, div, interp_bytes * q_comp_div, cudaMemcpyHostToDevice));
1248         if (in) {
1249           diag->d_interp_in = d_interp;
1250           diag->d_div_in    = d_div;
1251         } else {
1252           diag->d_interp_out = d_interp;
1253           diag->d_div_out    = d_div;
1254         }
1255       } break;
1256       case CEED_FE_SPACE_HCURL: {
1257         CeedInt           q_comp_interp, q_comp_curl;
1258         const CeedScalar *interp, *curl;
1259         CeedScalar       *d_interp, *d_curl;
1260 
1261         CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp));
1262         CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_CURL, &q_comp_curl));
1263 
1264         CeedCallBackend(CeedBasisGetInterp(basis, &interp));
1265         CeedCallCuda(ceed, cudaMalloc((void **)&d_interp, interp_bytes * q_comp_interp));
1266         CeedCallCuda(ceed, cudaMemcpy(d_interp, interp, interp_bytes * q_comp_interp, cudaMemcpyHostToDevice));
1267         CeedCallBackend(CeedBasisGetCurl(basis, &curl));
1268         CeedCallCuda(ceed, cudaMalloc((void **)&d_curl, interp_bytes * q_comp_curl));
1269         CeedCallCuda(ceed, cudaMemcpy(d_curl, curl, interp_bytes * q_comp_curl, cudaMemcpyHostToDevice));
1270         if (in) {
1271           diag->d_interp_in = d_interp;
1272           diag->d_curl_in   = d_curl;
1273         } else {
1274           diag->d_interp_out = d_interp;
1275           diag->d_curl_out   = d_curl;
1276         }
1277       } break;
1278     }
1279   }
1280 
1281   // Arrays of eval_modes
1282   CeedCallCuda(ceed, cudaMalloc((void **)&diag->d_eval_modes_in, num_eval_modes_in * eval_modes_bytes));
1283   CeedCallCuda(ceed, cudaMemcpy(diag->d_eval_modes_in, eval_modes_in, num_eval_modes_in * eval_modes_bytes, cudaMemcpyHostToDevice));
1284   CeedCallCuda(ceed, cudaMalloc((void **)&diag->d_eval_modes_out, num_eval_modes_out * eval_modes_bytes));
1285   CeedCallCuda(ceed, cudaMemcpy(diag->d_eval_modes_out, eval_modes_out, num_eval_modes_out * eval_modes_bytes, cudaMemcpyHostToDevice));
1286   CeedCallBackend(CeedFree(&eval_modes_in));
1287   CeedCallBackend(CeedFree(&eval_modes_out));
1288   CeedCallBackend(CeedDestroy(&ceed));
1289   CeedCallBackend(CeedBasisDestroy(&basis_in));
1290   CeedCallBackend(CeedBasisDestroy(&basis_out));
1291   CeedCallBackend(CeedQFunctionDestroy(&qf));
1292   return CEED_ERROR_SUCCESS;
1293 }
1294 
1295 //------------------------------------------------------------------------------
1296 // Assemble Diagonal Setup (Compilation)
1297 //------------------------------------------------------------------------------
1298 static inline int CeedOperatorAssembleDiagonalSetupCompile_Cuda(CeedOperator op, CeedInt use_ceedsize_idx, const bool is_point_block) {
1299   Ceed                ceed;
1300   CeedInt             num_input_fields, num_output_fields, num_eval_modes_in = 0, num_eval_modes_out = 0;
1301   CeedInt             num_comp, q_comp, num_nodes, num_qpts;
1302   CeedBasis           basis_in = NULL, basis_out = NULL;
1303   CeedQFunctionField *qf_fields;
1304   CeedQFunction       qf;
1305   CeedOperatorField  *op_fields;
1306   CeedOperator_Cuda  *impl;
1307 
1308   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
1309   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
1310   CeedCallBackend(CeedQFunctionGetNumArgs(qf, &num_input_fields, &num_output_fields));
1311 
1312   // Determine active input basis
1313   CeedCallBackend(CeedOperatorGetFields(op, NULL, &op_fields, NULL, NULL));
1314   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_fields, NULL, NULL));
1315   for (CeedInt i = 0; i < num_input_fields; i++) {
1316     CeedVector vec;
1317 
1318     CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec));
1319     if (vec == CEED_VECTOR_ACTIVE) {
1320       CeedEvalMode eval_mode;
1321       CeedBasis    basis;
1322 
1323       CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis));
1324       if (!basis_in) CeedCallBackend(CeedBasisReferenceCopy(basis, &basis_in));
1325       CeedCallBackend(CeedBasisDestroy(&basis));
1326       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode));
1327       CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_in, eval_mode, &q_comp));
1328       if (eval_mode != CEED_EVAL_WEIGHT) {
1329         num_eval_modes_in += q_comp;
1330       }
1331     }
1332     CeedCallBackend(CeedVectorDestroy(&vec));
1333   }
1334 
1335   // Determine active output basis
1336   CeedCallBackend(CeedOperatorGetFields(op, NULL, NULL, NULL, &op_fields));
1337   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qf_fields));
1338   for (CeedInt i = 0; i < num_output_fields; i++) {
1339     CeedVector vec;
1340 
1341     CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec));
1342     if (vec == CEED_VECTOR_ACTIVE) {
1343       CeedEvalMode eval_mode;
1344       CeedBasis    basis;
1345 
1346       CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis));
1347       if (!basis_out) CeedCallBackend(CeedBasisReferenceCopy(basis, &basis_out));
1348       CeedCallBackend(CeedBasisDestroy(&basis));
1349       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode));
1350       CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_out, eval_mode, &q_comp));
1351       if (eval_mode != CEED_EVAL_WEIGHT) {
1352         num_eval_modes_out += q_comp;
1353       }
1354     }
1355     CeedCallBackend(CeedVectorDestroy(&vec));
1356   }
1357 
1358   // Operator data struct
1359   CeedCallBackend(CeedOperatorGetData(op, &impl));
1360   CeedOperatorDiag_Cuda *diag = impl->diag;
1361 
1362   // Assemble kernel
1363   const char diagonal_kernel_source[] = "// Diagonal assembly source\n#include <ceed/jit-source/cuda/cuda-ref-operator-assemble-diagonal.h>\n";
1364   CUmodule  *module                   = is_point_block ? &diag->module_point_block : &diag->module;
1365   CeedInt    elems_per_block          = 1;
1366 
1367   CeedCallBackend(CeedBasisGetNumNodes(basis_in, &num_nodes));
1368   CeedCallBackend(CeedBasisGetNumComponents(basis_in, &num_comp));
1369   if (basis_in == CEED_BASIS_NONE) num_qpts = num_nodes;
1370   else CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis_in, &num_qpts));
1371   CeedCallCuda(ceed, CeedCompile_Cuda(ceed, diagonal_kernel_source, module, 8, "NUM_EVAL_MODES_IN", num_eval_modes_in, "NUM_EVAL_MODES_OUT",
1372                                       num_eval_modes_out, "NUM_COMP", num_comp, "NUM_NODES", num_nodes, "NUM_QPTS", num_qpts, "USE_CEEDSIZE",
1373                                       use_ceedsize_idx, "USE_POINT_BLOCK", is_point_block ? 1 : 0, "BLOCK_SIZE", num_nodes * elems_per_block));
1374   CeedCallCuda(ceed, CeedGetKernel_Cuda(ceed, *module, "LinearDiagonal", is_point_block ? &diag->LinearPointBlock : &diag->LinearDiagonal));
1375   CeedCallBackend(CeedDestroy(&ceed));
1376   CeedCallBackend(CeedBasisDestroy(&basis_in));
1377   CeedCallBackend(CeedBasisDestroy(&basis_out));
1378   CeedCallBackend(CeedQFunctionDestroy(&qf));
1379   return CEED_ERROR_SUCCESS;
1380 }
1381 
1382 //------------------------------------------------------------------------------
1383 // Assemble Diagonal Core
1384 //------------------------------------------------------------------------------
1385 static inline int CeedOperatorAssembleDiagonalCore_Cuda(CeedOperator op, CeedVector assembled, CeedRequest *request, const bool is_point_block) {
1386   Ceed                ceed;
1387   CeedInt             num_elem, num_nodes;
1388   CeedScalar         *elem_diag_array;
1389   const CeedScalar   *assembled_qf_array;
1390   CeedVector          assembled_qf   = NULL, elem_diag;
1391   CeedElemRestriction assembled_rstr = NULL, rstr_in, rstr_out, diag_rstr;
1392   CeedOperator_Cuda  *impl;
1393 
1394   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
1395   CeedCallBackend(CeedOperatorGetData(op, &impl));
1396 
1397   // Assemble QFunction
1398   CeedCallBackend(CeedOperatorLinearAssembleQFunctionBuildOrUpdate(op, &assembled_qf, &assembled_rstr, request));
1399   CeedCallBackend(CeedElemRestrictionDestroy(&assembled_rstr));
1400   CeedCallBackend(CeedVectorGetArrayRead(assembled_qf, CEED_MEM_DEVICE, &assembled_qf_array));
1401 
1402   // Setup
1403   if (!impl->diag) CeedCallBackend(CeedOperatorAssembleDiagonalSetup_Cuda(op));
1404   CeedOperatorDiag_Cuda *diag = impl->diag;
1405 
1406   assert(diag != NULL);
1407 
1408   // Assemble kernel if needed
1409   if ((!is_point_block && !diag->LinearDiagonal) || (is_point_block && !diag->LinearPointBlock)) {
1410     CeedSize assembled_length, assembled_qf_length;
1411     CeedInt  use_ceedsize_idx = 0;
1412     CeedCallBackend(CeedVectorGetLength(assembled, &assembled_length));
1413     CeedCallBackend(CeedVectorGetLength(assembled_qf, &assembled_qf_length));
1414     if ((assembled_length > INT_MAX) || (assembled_qf_length > INT_MAX)) use_ceedsize_idx = 1;
1415 
1416     CeedCallBackend(CeedOperatorAssembleDiagonalSetupCompile_Cuda(op, use_ceedsize_idx, is_point_block));
1417   }
1418 
1419   // Restriction and diagonal vector
1420   CeedCallBackend(CeedOperatorGetActiveElemRestrictions(op, &rstr_in, &rstr_out));
1421   CeedCheck(rstr_in == rstr_out, ceed, CEED_ERROR_BACKEND,
1422             "Cannot assemble operator diagonal with different input and output active element restrictions");
1423   if (!is_point_block && !diag->diag_rstr) {
1424     CeedCallBackend(CeedElemRestrictionCreateUnsignedCopy(rstr_out, &diag->diag_rstr));
1425     CeedCallBackend(CeedElemRestrictionCreateVector(diag->diag_rstr, NULL, &diag->elem_diag));
1426   } else if (is_point_block && !diag->point_block_diag_rstr) {
1427     CeedCallBackend(CeedOperatorCreateActivePointBlockRestriction(rstr_out, &diag->point_block_diag_rstr));
1428     CeedCallBackend(CeedElemRestrictionCreateVector(diag->point_block_diag_rstr, NULL, &diag->point_block_elem_diag));
1429   }
1430   CeedCallBackend(CeedElemRestrictionDestroy(&rstr_in));
1431   CeedCallBackend(CeedElemRestrictionDestroy(&rstr_out));
1432   diag_rstr = is_point_block ? diag->point_block_diag_rstr : diag->diag_rstr;
1433   elem_diag = is_point_block ? diag->point_block_elem_diag : diag->elem_diag;
1434   CeedCallBackend(CeedVectorSetValue(elem_diag, 0.0));
1435 
1436   // Only assemble diagonal if the basis has nodes, otherwise inputs are null pointers
1437   CeedCallBackend(CeedElemRestrictionGetElementSize(diag_rstr, &num_nodes));
1438   if (num_nodes > 0) {
1439     // Assemble element operator diagonals
1440     CeedCallBackend(CeedElemRestrictionGetNumElements(diag_rstr, &num_elem));
1441     CeedCallBackend(CeedVectorGetArray(elem_diag, CEED_MEM_DEVICE, &elem_diag_array));
1442 
1443     // Compute the diagonal of B^T D B
1444     CeedInt elems_per_block = 1;
1445     CeedInt grid            = CeedDivUpInt(num_elem, elems_per_block);
1446     void   *args[]          = {(void *)&num_elem,      &diag->d_identity,       &diag->d_interp_in,  &diag->d_grad_in, &diag->d_div_in,
1447                                &diag->d_curl_in,       &diag->d_interp_out,     &diag->d_grad_out,   &diag->d_div_out, &diag->d_curl_out,
1448                                &diag->d_eval_modes_in, &diag->d_eval_modes_out, &assembled_qf_array, &elem_diag_array};
1449 
1450     if (is_point_block) {
1451       CeedCallBackend(CeedRunKernelDim_Cuda(ceed, diag->LinearPointBlock, grid, num_nodes, 1, elems_per_block, args));
1452     } else {
1453       CeedCallBackend(CeedRunKernelDim_Cuda(ceed, diag->LinearDiagonal, grid, num_nodes, 1, elems_per_block, args));
1454     }
1455 
1456     // Restore arrays
1457     CeedCallBackend(CeedVectorRestoreArray(elem_diag, &elem_diag_array));
1458     CeedCallBackend(CeedVectorRestoreArrayRead(assembled_qf, &assembled_qf_array));
1459   }
1460 
1461   // Assemble local operator diagonal
1462   CeedCallBackend(CeedElemRestrictionApply(diag_rstr, CEED_TRANSPOSE, elem_diag, assembled, request));
1463 
1464   // Cleanup
1465   CeedCallBackend(CeedDestroy(&ceed));
1466   CeedCallBackend(CeedVectorDestroy(&assembled_qf));
1467   return CEED_ERROR_SUCCESS;
1468 }
1469 
1470 //------------------------------------------------------------------------------
1471 // Assemble Linear Diagonal
1472 //------------------------------------------------------------------------------
1473 static int CeedOperatorLinearAssembleAddDiagonal_Cuda(CeedOperator op, CeedVector assembled, CeedRequest *request) {
1474   CeedCallBackend(CeedOperatorAssembleDiagonalCore_Cuda(op, assembled, request, false));
1475   return CEED_ERROR_SUCCESS;
1476 }
1477 
1478 //------------------------------------------------------------------------------
1479 // Assemble Linear Point Block Diagonal
1480 //------------------------------------------------------------------------------
1481 static int CeedOperatorLinearAssembleAddPointBlockDiagonal_Cuda(CeedOperator op, CeedVector assembled, CeedRequest *request) {
1482   CeedCallBackend(CeedOperatorAssembleDiagonalCore_Cuda(op, assembled, request, true));
1483   return CEED_ERROR_SUCCESS;
1484 }
1485 
1486 //------------------------------------------------------------------------------
1487 // Single Operator Assembly Setup
1488 //------------------------------------------------------------------------------
1489 static int CeedSingleOperatorAssembleSetup_Cuda(CeedOperator op, CeedInt use_ceedsize_idx) {
1490   Ceed                ceed;
1491   Ceed_Cuda          *cuda_data;
1492   CeedInt             num_input_fields, num_output_fields, num_eval_modes_in = 0, num_eval_modes_out = 0;
1493   CeedInt             elem_size_in, num_qpts_in = 0, num_comp_in, elem_size_out, num_qpts_out, num_comp_out, q_comp;
1494   CeedEvalMode       *eval_modes_in = NULL, *eval_modes_out = NULL;
1495   CeedElemRestriction rstr_in = NULL, rstr_out = NULL;
1496   CeedBasis           basis_in = NULL, basis_out = NULL;
1497   CeedQFunctionField *qf_fields;
1498   CeedQFunction       qf;
1499   CeedOperatorField  *input_fields, *output_fields;
1500   CeedOperator_Cuda  *impl;
1501 
1502   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
1503   CeedCallBackend(CeedOperatorGetData(op, &impl));
1504 
1505   // Get intput and output fields
1506   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &input_fields, &num_output_fields, &output_fields));
1507 
1508   // Determine active input basis eval mode
1509   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
1510   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_fields, NULL, NULL));
1511   for (CeedInt i = 0; i < num_input_fields; i++) {
1512     CeedVector vec;
1513 
1514     CeedCallBackend(CeedOperatorFieldGetVector(input_fields[i], &vec));
1515     if (vec == CEED_VECTOR_ACTIVE) {
1516       CeedEvalMode        eval_mode;
1517       CeedElemRestriction elem_rstr;
1518       CeedBasis           basis;
1519 
1520       CeedCallBackend(CeedOperatorFieldGetBasis(input_fields[i], &basis));
1521       CeedCheck(!basis_in || basis_in == basis, ceed, CEED_ERROR_BACKEND, "Backend does not implement operator assembly with multiple active bases");
1522       if (!basis_in) CeedCallBackend(CeedBasisReferenceCopy(basis, &basis_in));
1523       CeedCallBackend(CeedBasisDestroy(&basis));
1524       CeedCallBackend(CeedOperatorFieldGetElemRestriction(input_fields[i], &elem_rstr));
1525       if (!rstr_in) CeedCallBackend(CeedElemRestrictionReferenceCopy(elem_rstr, &rstr_in));
1526       CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
1527       CeedCallBackend(CeedElemRestrictionGetElementSize(rstr_in, &elem_size_in));
1528       if (basis_in == CEED_BASIS_NONE) num_qpts_in = elem_size_in;
1529       else CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis_in, &num_qpts_in));
1530       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode));
1531       CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_in, eval_mode, &q_comp));
1532       if (eval_mode != CEED_EVAL_WEIGHT) {
1533         // q_comp = 1 if CEED_EVAL_NONE, CEED_EVAL_WEIGHT caught by QF Assembly
1534         CeedCallBackend(CeedRealloc(num_eval_modes_in + q_comp, &eval_modes_in));
1535         for (CeedInt d = 0; d < q_comp; d++) {
1536           eval_modes_in[num_eval_modes_in + d] = eval_mode;
1537         }
1538         num_eval_modes_in += q_comp;
1539       }
1540     }
1541     CeedCallBackend(CeedVectorDestroy(&vec));
1542   }
1543 
1544   // Determine active output basis; basis_out and rstr_out only used if same as input, TODO
1545   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qf_fields));
1546   for (CeedInt i = 0; i < num_output_fields; i++) {
1547     CeedVector vec;
1548 
1549     CeedCallBackend(CeedOperatorFieldGetVector(output_fields[i], &vec));
1550     if (vec == CEED_VECTOR_ACTIVE) {
1551       CeedEvalMode        eval_mode;
1552       CeedElemRestriction elem_rstr;
1553       CeedBasis           basis;
1554 
1555       CeedCallBackend(CeedOperatorFieldGetBasis(output_fields[i], &basis));
1556       CeedCheck(!basis_out || basis_out == basis, ceed, CEED_ERROR_BACKEND,
1557                 "Backend does not implement operator assembly with multiple active bases");
1558       if (!basis_out) CeedCallBackend(CeedBasisReferenceCopy(basis, &basis_out));
1559       CeedCallBackend(CeedBasisDestroy(&basis));
1560       CeedCallBackend(CeedOperatorFieldGetElemRestriction(output_fields[i], &elem_rstr));
1561       if (!rstr_out) CeedCallBackend(CeedElemRestrictionReferenceCopy(elem_rstr, &rstr_out));
1562       CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
1563       CeedCallBackend(CeedElemRestrictionGetElementSize(rstr_out, &elem_size_out));
1564       if (basis_out == CEED_BASIS_NONE) num_qpts_out = elem_size_out;
1565       else CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis_out, &num_qpts_out));
1566       CeedCheck(num_qpts_in == num_qpts_out, ceed, CEED_ERROR_UNSUPPORTED,
1567                 "Active input and output bases must have the same number of quadrature points");
1568       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode));
1569       CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_out, eval_mode, &q_comp));
1570       if (eval_mode != CEED_EVAL_WEIGHT) {
1571         // q_comp = 1 if CEED_EVAL_NONE, CEED_EVAL_WEIGHT caught by QF Assembly
1572         CeedCallBackend(CeedRealloc(num_eval_modes_out + q_comp, &eval_modes_out));
1573         for (CeedInt d = 0; d < q_comp; d++) {
1574           eval_modes_out[num_eval_modes_out + d] = eval_mode;
1575         }
1576         num_eval_modes_out += q_comp;
1577       }
1578     }
1579     CeedCallBackend(CeedVectorDestroy(&vec));
1580   }
1581   CeedCheck(num_eval_modes_in > 0 && num_eval_modes_out > 0, ceed, CEED_ERROR_UNSUPPORTED, "Cannot assemble operator without inputs/outputs");
1582 
1583   CeedCallBackend(CeedCalloc(1, &impl->asmb));
1584   CeedOperatorAssemble_Cuda *asmb = impl->asmb;
1585   asmb->elems_per_block           = 1;
1586   asmb->block_size_x              = elem_size_in;
1587   asmb->block_size_y              = elem_size_out;
1588 
1589   CeedCallBackend(CeedGetData(ceed, &cuda_data));
1590   bool fallback = asmb->block_size_x * asmb->block_size_y * asmb->elems_per_block > cuda_data->device_prop.maxThreadsPerBlock;
1591 
1592   if (fallback) {
1593     // Use fallback kernel with 1D threadblock
1594     asmb->block_size_y = 1;
1595   }
1596 
1597   // Compile kernels
1598   const char assembly_kernel_source[] = "// Full assembly source\n#include <ceed/jit-source/cuda/cuda-ref-operator-assemble.h>\n";
1599 
1600   CeedCallBackend(CeedElemRestrictionGetNumComponents(rstr_in, &num_comp_in));
1601   CeedCallBackend(CeedElemRestrictionGetNumComponents(rstr_out, &num_comp_out));
1602   CeedCallBackend(CeedCompile_Cuda(ceed, assembly_kernel_source, &asmb->module, 10, "NUM_EVAL_MODES_IN", num_eval_modes_in, "NUM_EVAL_MODES_OUT",
1603                                    num_eval_modes_out, "NUM_COMP_IN", num_comp_in, "NUM_COMP_OUT", num_comp_out, "NUM_NODES_IN", elem_size_in,
1604                                    "NUM_NODES_OUT", elem_size_out, "NUM_QPTS", num_qpts_in, "BLOCK_SIZE",
1605                                    asmb->block_size_x * asmb->block_size_y * asmb->elems_per_block, "BLOCK_SIZE_Y", asmb->block_size_y,
1606                                    "USE_CEEDSIZE", use_ceedsize_idx));
1607   CeedCallBackend(CeedGetKernel_Cuda(ceed, asmb->module, "LinearAssemble", &asmb->LinearAssemble));
1608 
1609   // Load into B_in, in order that they will be used in eval_modes_in
1610   {
1611     const CeedInt in_bytes           = elem_size_in * num_qpts_in * num_eval_modes_in * sizeof(CeedScalar);
1612     CeedInt       d_in               = 0;
1613     CeedEvalMode  eval_modes_in_prev = CEED_EVAL_NONE;
1614     bool          has_eval_none      = false;
1615     CeedScalar   *identity           = NULL;
1616 
1617     for (CeedInt i = 0; i < num_eval_modes_in; i++) {
1618       has_eval_none = has_eval_none || (eval_modes_in[i] == CEED_EVAL_NONE);
1619     }
1620     if (has_eval_none) {
1621       CeedCallBackend(CeedCalloc(elem_size_in * num_qpts_in, &identity));
1622       for (CeedInt i = 0; i < (elem_size_in < num_qpts_in ? elem_size_in : num_qpts_in); i++) identity[i * elem_size_in + i] = 1.0;
1623     }
1624 
1625     CeedCallCuda(ceed, cudaMalloc((void **)&asmb->d_B_in, in_bytes));
1626     for (CeedInt i = 0; i < num_eval_modes_in; i++) {
1627       const CeedScalar *h_B_in;
1628 
1629       CeedCallBackend(CeedOperatorGetBasisPointer(basis_in, eval_modes_in[i], identity, &h_B_in));
1630       CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_in, eval_modes_in[i], &q_comp));
1631       if (q_comp > 1) {
1632         if (i == 0 || eval_modes_in[i] != eval_modes_in_prev) d_in = 0;
1633         else h_B_in = &h_B_in[(++d_in) * elem_size_in * num_qpts_in];
1634       }
1635       eval_modes_in_prev = eval_modes_in[i];
1636 
1637       CeedCallCuda(ceed, cudaMemcpy(&asmb->d_B_in[i * elem_size_in * num_qpts_in], h_B_in, elem_size_in * num_qpts_in * sizeof(CeedScalar),
1638                                     cudaMemcpyHostToDevice));
1639     }
1640     CeedCallBackend(CeedFree(&identity));
1641   }
1642   CeedCallBackend(CeedFree(&eval_modes_in));
1643 
1644   // Load into B_out, in order that they will be used in eval_modes_out
1645   {
1646     const CeedInt out_bytes           = elem_size_out * num_qpts_out * num_eval_modes_out * sizeof(CeedScalar);
1647     CeedInt       d_out               = 0;
1648     CeedEvalMode  eval_modes_out_prev = CEED_EVAL_NONE;
1649     bool          has_eval_none       = false;
1650     CeedScalar   *identity            = NULL;
1651 
1652     for (CeedInt i = 0; i < num_eval_modes_out; i++) {
1653       has_eval_none = has_eval_none || (eval_modes_out[i] == CEED_EVAL_NONE);
1654     }
1655     if (has_eval_none) {
1656       CeedCallBackend(CeedCalloc(elem_size_out * num_qpts_out, &identity));
1657       for (CeedInt i = 0; i < (elem_size_out < num_qpts_out ? elem_size_out : num_qpts_out); i++) identity[i * elem_size_out + i] = 1.0;
1658     }
1659 
1660     CeedCallCuda(ceed, cudaMalloc((void **)&asmb->d_B_out, out_bytes));
1661     for (CeedInt i = 0; i < num_eval_modes_out; i++) {
1662       const CeedScalar *h_B_out;
1663 
1664       CeedCallBackend(CeedOperatorGetBasisPointer(basis_out, eval_modes_out[i], identity, &h_B_out));
1665       CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_out, eval_modes_out[i], &q_comp));
1666       if (q_comp > 1) {
1667         if (i == 0 || eval_modes_out[i] != eval_modes_out_prev) d_out = 0;
1668         else h_B_out = &h_B_out[(++d_out) * elem_size_out * num_qpts_out];
1669       }
1670       eval_modes_out_prev = eval_modes_out[i];
1671 
1672       CeedCallCuda(ceed, cudaMemcpy(&asmb->d_B_out[i * elem_size_out * num_qpts_out], h_B_out, elem_size_out * num_qpts_out * sizeof(CeedScalar),
1673                                     cudaMemcpyHostToDevice));
1674     }
1675     CeedCallBackend(CeedFree(&identity));
1676   }
1677   CeedCallBackend(CeedFree(&eval_modes_out));
1678   CeedCallBackend(CeedDestroy(&ceed));
1679   CeedCallBackend(CeedElemRestrictionDestroy(&rstr_in));
1680   CeedCallBackend(CeedElemRestrictionDestroy(&rstr_out));
1681   CeedCallBackend(CeedBasisDestroy(&basis_in));
1682   CeedCallBackend(CeedBasisDestroy(&basis_out));
1683   CeedCallBackend(CeedQFunctionDestroy(&qf));
1684   return CEED_ERROR_SUCCESS;
1685 }
1686 
1687 //------------------------------------------------------------------------------
1688 // Assemble matrix data for COO matrix of assembled operator.
1689 // The sparsity pattern is set by CeedOperatorLinearAssembleSymbolic.
1690 //
1691 // Note that this (and other assembly routines) currently assume only one active input restriction/basis per operator
1692 // (could have multiple basis eval modes).
1693 // TODO: allow multiple active input restrictions/basis objects
1694 //------------------------------------------------------------------------------
1695 static int CeedSingleOperatorAssemble_Cuda(CeedOperator op, CeedInt offset, CeedVector values) {
1696   Ceed                ceed;
1697   CeedSize            values_length = 0, assembled_qf_length = 0;
1698   CeedInt             use_ceedsize_idx = 0, num_elem_in, num_elem_out, elem_size_in, elem_size_out;
1699   CeedScalar         *values_array;
1700   const CeedScalar   *assembled_qf_array;
1701   CeedVector          assembled_qf   = NULL;
1702   CeedElemRestriction assembled_rstr = NULL, rstr_in, rstr_out;
1703   CeedRestrictionType rstr_type_in, rstr_type_out;
1704   const bool         *orients_in = NULL, *orients_out = NULL;
1705   const CeedInt8     *curl_orients_in = NULL, *curl_orients_out = NULL;
1706   CeedOperator_Cuda  *impl;
1707 
1708   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
1709   CeedCallBackend(CeedOperatorGetData(op, &impl));
1710 
1711   // Assemble QFunction
1712   CeedCallBackend(CeedOperatorLinearAssembleQFunctionBuildOrUpdate(op, &assembled_qf, &assembled_rstr, CEED_REQUEST_IMMEDIATE));
1713   CeedCallBackend(CeedElemRestrictionDestroy(&assembled_rstr));
1714   CeedCallBackend(CeedVectorGetArrayRead(assembled_qf, CEED_MEM_DEVICE, &assembled_qf_array));
1715 
1716   CeedCallBackend(CeedVectorGetLength(values, &values_length));
1717   CeedCallBackend(CeedVectorGetLength(assembled_qf, &assembled_qf_length));
1718   if ((values_length > INT_MAX) || (assembled_qf_length > INT_MAX)) use_ceedsize_idx = 1;
1719 
1720   // Setup
1721   if (!impl->asmb) CeedCallBackend(CeedSingleOperatorAssembleSetup_Cuda(op, use_ceedsize_idx));
1722   CeedOperatorAssemble_Cuda *asmb = impl->asmb;
1723 
1724   assert(asmb != NULL);
1725 
1726   // Assemble element operator
1727   CeedCallBackend(CeedVectorGetArray(values, CEED_MEM_DEVICE, &values_array));
1728   values_array += offset;
1729 
1730   CeedCallBackend(CeedOperatorGetActiveElemRestrictions(op, &rstr_in, &rstr_out));
1731   CeedCallBackend(CeedElemRestrictionGetNumElements(rstr_in, &num_elem_in));
1732   CeedCallBackend(CeedElemRestrictionGetElementSize(rstr_in, &elem_size_in));
1733 
1734   CeedCallBackend(CeedElemRestrictionGetType(rstr_in, &rstr_type_in));
1735   if (rstr_type_in == CEED_RESTRICTION_ORIENTED) {
1736     CeedCallBackend(CeedElemRestrictionGetOrientations(rstr_in, CEED_MEM_DEVICE, &orients_in));
1737   } else if (rstr_type_in == CEED_RESTRICTION_CURL_ORIENTED) {
1738     CeedCallBackend(CeedElemRestrictionGetCurlOrientations(rstr_in, CEED_MEM_DEVICE, &curl_orients_in));
1739   }
1740 
1741   if (rstr_in != rstr_out) {
1742     CeedCallBackend(CeedElemRestrictionGetNumElements(rstr_out, &num_elem_out));
1743     CeedCheck(num_elem_in == num_elem_out, ceed, CEED_ERROR_UNSUPPORTED,
1744               "Active input and output operator restrictions must have the same number of elements");
1745     CeedCallBackend(CeedElemRestrictionGetElementSize(rstr_out, &elem_size_out));
1746 
1747     CeedCallBackend(CeedElemRestrictionGetType(rstr_out, &rstr_type_out));
1748     if (rstr_type_out == CEED_RESTRICTION_ORIENTED) {
1749       CeedCallBackend(CeedElemRestrictionGetOrientations(rstr_out, CEED_MEM_DEVICE, &orients_out));
1750     } else if (rstr_type_out == CEED_RESTRICTION_CURL_ORIENTED) {
1751       CeedCallBackend(CeedElemRestrictionGetCurlOrientations(rstr_out, CEED_MEM_DEVICE, &curl_orients_out));
1752     }
1753   } else {
1754     elem_size_out    = elem_size_in;
1755     orients_out      = orients_in;
1756     curl_orients_out = curl_orients_in;
1757   }
1758 
1759   // Compute B^T D B
1760   CeedInt shared_mem =
1761       ((curl_orients_in || curl_orients_out ? elem_size_in * elem_size_out : 0) + (curl_orients_in ? elem_size_in * asmb->block_size_y : 0)) *
1762       sizeof(CeedScalar);
1763   CeedInt grid   = CeedDivUpInt(num_elem_in, asmb->elems_per_block);
1764   void   *args[] = {(void *)&num_elem_in, &asmb->d_B_in,     &asmb->d_B_out,      &orients_in,  &curl_orients_in,
1765                     &orients_out,         &curl_orients_out, &assembled_qf_array, &values_array};
1766 
1767   CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, asmb->LinearAssemble, NULL, grid, asmb->block_size_x, asmb->block_size_y, asmb->elems_per_block,
1768                                               shared_mem, args));
1769 
1770   // Restore arrays
1771   CeedCallBackend(CeedVectorRestoreArray(values, &values_array));
1772   CeedCallBackend(CeedVectorRestoreArrayRead(assembled_qf, &assembled_qf_array));
1773 
1774   // Cleanup
1775   CeedCallBackend(CeedVectorDestroy(&assembled_qf));
1776   if (rstr_type_in == CEED_RESTRICTION_ORIENTED) {
1777     CeedCallBackend(CeedElemRestrictionRestoreOrientations(rstr_in, &orients_in));
1778   } else if (rstr_type_in == CEED_RESTRICTION_CURL_ORIENTED) {
1779     CeedCallBackend(CeedElemRestrictionRestoreCurlOrientations(rstr_in, &curl_orients_in));
1780   }
1781   if (rstr_in != rstr_out) {
1782     if (rstr_type_out == CEED_RESTRICTION_ORIENTED) {
1783       CeedCallBackend(CeedElemRestrictionRestoreOrientations(rstr_out, &orients_out));
1784     } else if (rstr_type_out == CEED_RESTRICTION_CURL_ORIENTED) {
1785       CeedCallBackend(CeedElemRestrictionRestoreCurlOrientations(rstr_out, &curl_orients_out));
1786     }
1787   }
1788   CeedCallBackend(CeedDestroy(&ceed));
1789   CeedCallBackend(CeedElemRestrictionDestroy(&rstr_in));
1790   CeedCallBackend(CeedElemRestrictionDestroy(&rstr_out));
1791   return CEED_ERROR_SUCCESS;
1792 }
1793 
1794 //------------------------------------------------------------------------------
1795 // Assemble Linear QFunction AtPoints
1796 //------------------------------------------------------------------------------
1797 static int CeedOperatorLinearAssembleQFunctionAtPoints_Cuda(CeedOperator op, CeedVector *assembled, CeedElemRestriction *rstr, CeedRequest *request) {
1798   return CeedError(CeedOperatorReturnCeed(op), CEED_ERROR_BACKEND, "Backend does not implement CeedOperatorLinearAssembleQFunction");
1799 }
1800 
1801 //------------------------------------------------------------------------------
1802 // Assemble Linear Diagonal AtPoints
1803 //------------------------------------------------------------------------------
1804 static int CeedOperatorLinearAssembleAddDiagonalAtPoints_Cuda(CeedOperator op, CeedVector assembled, CeedRequest *request) {
1805   CeedInt             max_num_points, *num_points, num_elem, num_input_fields, num_output_fields;
1806   Ceed                ceed;
1807   CeedVector          active_e_vec_in, active_e_vec_out;
1808   CeedQFunctionField *qf_input_fields, *qf_output_fields;
1809   CeedQFunction       qf;
1810   CeedOperatorField  *op_input_fields, *op_output_fields;
1811   CeedOperator_Cuda  *impl;
1812 
1813   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
1814   CeedCallBackend(CeedOperatorGetData(op, &impl));
1815   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
1816   CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem));
1817   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
1818   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
1819 
1820   // Setup
1821   CeedCallBackend(CeedOperatorSetupAtPoints_Cuda(op));
1822   num_points     = impl->num_points;
1823   max_num_points = impl->max_num_points;
1824 
1825   // Work vector
1826   CeedCallBackend(CeedGetWorkVector(ceed, impl->max_active_e_vec_len, &active_e_vec_in));
1827   CeedCallBackend(CeedGetWorkVector(ceed, impl->max_active_e_vec_len, &active_e_vec_out));
1828   {
1829     CeedSize length_in, length_out;
1830 
1831     CeedCallBackend(CeedVectorGetLength(active_e_vec_in, &length_in));
1832     CeedCallBackend(CeedVectorGetLength(active_e_vec_out, &length_out));
1833     // Need input e_vec to be longer
1834     if (length_in < length_out) {
1835       CeedVector temp = active_e_vec_in;
1836 
1837       active_e_vec_in  = active_e_vec_out;
1838       active_e_vec_out = temp;
1839     }
1840   }
1841 
1842   // Get point coordinates
1843   if (!impl->point_coords_elem) {
1844     CeedVector          point_coords = NULL;
1845     CeedElemRestriction rstr_points  = NULL;
1846 
1847     CeedCallBackend(CeedOperatorAtPointsGetPoints(op, &rstr_points, &point_coords));
1848     CeedCallBackend(CeedElemRestrictionCreateVector(rstr_points, NULL, &impl->point_coords_elem));
1849     CeedCallBackend(CeedElemRestrictionApply(rstr_points, CEED_NOTRANSPOSE, point_coords, impl->point_coords_elem, request));
1850     CeedCallBackend(CeedVectorDestroy(&point_coords));
1851     CeedCallBackend(CeedElemRestrictionDestroy(&rstr_points));
1852   }
1853 
1854   // Process inputs
1855   for (CeedInt i = 0; i < num_input_fields; i++) {
1856     CeedCallBackend(CeedOperatorInputRestrict_Cuda(op_input_fields[i], qf_input_fields[i], i, NULL, NULL, true, impl, request));
1857     CeedCallBackend(
1858         CeedOperatorInputBasisAtPoints_Cuda(op_input_fields[i], qf_input_fields[i], i, NULL, NULL, num_elem, num_points, true, false, impl));
1859   }
1860 
1861   // Output pointers, as necessary
1862   for (CeedInt i = 0; i < num_output_fields; i++) {
1863     CeedEvalMode eval_mode;
1864 
1865     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
1866     if (eval_mode == CEED_EVAL_NONE) {
1867       CeedScalar *e_vec_array;
1868 
1869       CeedCallBackend(CeedVectorGetArrayWrite(impl->e_vecs_out[i], CEED_MEM_DEVICE, &e_vec_array));
1870       CeedCallBackend(CeedVectorSetArray(impl->q_vecs_out[i], CEED_MEM_DEVICE, CEED_USE_POINTER, e_vec_array));
1871     }
1872   }
1873 
1874   // Loop over active fields
1875   for (CeedInt i = 0; i < num_input_fields; i++) {
1876     bool                is_active = false, is_active_at_points = true;
1877     CeedInt             elem_size = 1, num_comp_active = 1, e_vec_size = 0, field_in = impl->input_field_order[i];
1878     CeedRestrictionType rstr_type;
1879     CeedVector          l_vec;
1880     CeedElemRestriction elem_rstr;
1881 
1882     // -- Skip non-active input
1883     CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[field_in], &l_vec));
1884     is_active = l_vec == CEED_VECTOR_ACTIVE;
1885     CeedCallBackend(CeedVectorDestroy(&l_vec));
1886     if (!is_active || impl->skip_rstr_in[field_in]) continue;
1887 
1888     // -- Get active restriction type
1889     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[field_in], &elem_rstr));
1890     CeedCallBackend(CeedElemRestrictionGetType(elem_rstr, &rstr_type));
1891     is_active_at_points = rstr_type == CEED_RESTRICTION_POINTS;
1892     if (!is_active_at_points) CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
1893     else elem_size = max_num_points;
1894     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp_active));
1895     CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
1896 
1897     e_vec_size = elem_size * num_comp_active;
1898     CeedCallBackend(CeedVectorSetValue(active_e_vec_in, 0.0));
1899     for (CeedInt s = 0; s < e_vec_size; s++) {
1900       CeedVector q_vec = impl->q_vecs_in[field_in];
1901 
1902       // Update unit vector
1903       {
1904         // Note: E-vec strides are node * (1) + comp * (elem_size * num_elem) + elem * (elem_size)
1905         CeedInt  node = (s - 1) % elem_size, comp = (s - 1) / elem_size;
1906         CeedSize start = node * 1 + comp * (elem_size * num_elem);
1907         CeedSize stop  = (comp + 1) * (elem_size * num_elem);
1908 
1909         if (s != 0) CeedCallBackend(CeedVectorSetValueStrided(active_e_vec_in, start, stop, elem_size, 0.0));
1910 
1911         node = s % elem_size, comp = s / elem_size;
1912         start = node * 1 + comp * (elem_size * num_elem);
1913         stop  = (comp + 1) * (elem_size * num_elem);
1914         CeedCallBackend(CeedVectorSetValueStrided(active_e_vec_in, start, stop, elem_size, 1.0));
1915       }
1916 
1917       // Basis action
1918       for (CeedInt j = 0; j < num_input_fields; j++) {
1919         CeedInt field = impl->input_field_order[j];
1920 
1921         CeedCallBackend(CeedOperatorInputBasisAtPoints_Cuda(op_input_fields[field], qf_input_fields[field], field, NULL, active_e_vec_in, num_elem,
1922                                                             num_points, false, true, impl));
1923       }
1924 
1925       // Q function
1926       CeedCallBackend(CeedQFunctionApply(qf, num_elem * max_num_points, impl->q_vecs_in, impl->q_vecs_out));
1927 
1928       // Output basis apply if needed
1929       for (CeedInt j = 0; j < num_output_fields; j++) {
1930         bool                is_active = false;
1931         CeedInt             elem_size = 0;
1932         CeedInt             field_out = impl->output_field_order[j];
1933         CeedRestrictionType rstr_type;
1934         CeedEvalMode        eval_mode;
1935         CeedVector          l_vec, e_vec = impl->e_vecs_out[field_out], q_vec = impl->q_vecs_out[field_out];
1936         CeedElemRestriction elem_rstr;
1937 
1938         // ---- Skip non-active output
1939         CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[field_out], &l_vec));
1940         is_active = l_vec == CEED_VECTOR_ACTIVE;
1941         CeedCallBackend(CeedVectorDestroy(&l_vec));
1942         if (!is_active) continue;
1943         if (!e_vec) e_vec = active_e_vec_out;
1944 
1945         // ---- Check if elem size matches
1946         CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[field_out], &elem_rstr));
1947         CeedCallBackend(CeedElemRestrictionGetType(elem_rstr, &rstr_type));
1948         if (is_active_at_points && rstr_type != CEED_RESTRICTION_POINTS) continue;
1949         if (rstr_type == CEED_RESTRICTION_POINTS) {
1950           CeedCallBackend(CeedElemRestrictionGetMaxPointsInElement(elem_rstr, &elem_size));
1951         } else {
1952           CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
1953         }
1954         {
1955           CeedInt num_comp = 0;
1956 
1957           CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
1958           if (e_vec_size != num_comp * elem_size) continue;
1959         }
1960 
1961         // Basis action
1962         CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[field_out], &eval_mode));
1963         switch (eval_mode) {
1964           case CEED_EVAL_NONE: {
1965             CeedScalar *e_vec_array;
1966 
1967             CeedCallBackend(CeedVectorTakeArray(q_vec, CEED_MEM_DEVICE, &e_vec_array));
1968             CeedCallBackend(CeedVectorRestoreArray(e_vec, &e_vec_array));
1969             break;
1970           }
1971           case CEED_EVAL_INTERP:
1972           case CEED_EVAL_GRAD:
1973           case CEED_EVAL_DIV:
1974           case CEED_EVAL_CURL: {
1975             CeedBasis basis;
1976 
1977             CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[field_out], &basis));
1978             if (impl->apply_add_basis_out[field_out]) {
1979               CeedCallBackend(
1980                   CeedBasisApplyAddAtPoints(basis, num_elem, num_points, CEED_TRANSPOSE, eval_mode, impl->point_coords_elem, q_vec, e_vec));
1981             } else {
1982               CeedCallBackend(CeedBasisApplyAtPoints(basis, num_elem, num_points, CEED_TRANSPOSE, eval_mode, impl->point_coords_elem, q_vec, e_vec));
1983             }
1984             CeedCallBackend(CeedBasisDestroy(&basis));
1985             break;
1986           }
1987           // LCOV_EXCL_START
1988           case CEED_EVAL_WEIGHT: {
1989             return CeedError(CeedOperatorReturnCeed(op), CEED_ERROR_BACKEND, "CEED_EVAL_WEIGHT cannot be an output evaluation mode");
1990             // LCOV_EXCL_STOP
1991           }
1992         }
1993 
1994         // Mask output e-vec
1995         if (impl->skip_rstr_out[field_out]) {
1996           CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
1997           continue;
1998         }
1999         CeedCallBackend(CeedVectorPointwiseMult(e_vec, active_e_vec_in, e_vec));
2000 
2001         // Restrict
2002         CeedCallBackend(CeedElemRestrictionApply(elem_rstr, CEED_TRANSPOSE, e_vec, assembled, request));
2003         CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
2004 
2005         // Reset q_vec for
2006         if (eval_mode == CEED_EVAL_NONE) {
2007           CeedScalar *e_vec_array;
2008 
2009           CeedCallBackend(CeedVectorGetArrayWrite(e_vec, CEED_MEM_DEVICE, &e_vec_array));
2010           CeedCallBackend(CeedVectorSetArray(q_vec, CEED_MEM_DEVICE, CEED_USE_POINTER, e_vec_array));
2011         }
2012       }
2013 
2014       // Reset vec
2015       if (s == e_vec_size - 1 && i != num_input_fields - 1) CeedCallBackend(CeedVectorSetValue(q_vec, 0.0));
2016     }
2017   }
2018 
2019   // Restore CEED_EVAL_NONE
2020   for (CeedInt i = 0; i < num_output_fields; i++) {
2021     CeedEvalMode eval_mode;
2022 
2023     // Get eval_mode
2024     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
2025 
2026     // Restore evec
2027     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
2028     if (eval_mode == CEED_EVAL_NONE) {
2029       CeedScalar *e_vec_array;
2030 
2031       CeedCallBackend(CeedVectorTakeArray(impl->q_vecs_in[i], CEED_MEM_DEVICE, &e_vec_array));
2032       CeedCallBackend(CeedVectorRestoreArray(impl->e_vecs_in[i], &e_vec_array));
2033     }
2034   }
2035 
2036   // Restore input arrays
2037   for (CeedInt i = 0; i < num_input_fields; i++) {
2038     CeedCallBackend(CeedOperatorInputRestore_Cuda(op_input_fields[i], qf_input_fields[i], i, NULL, NULL, true, impl));
2039   }
2040 
2041   // Restore work vector
2042   CeedCallBackend(CeedRestoreWorkVector(ceed, &active_e_vec_in));
2043   CeedCallBackend(CeedRestoreWorkVector(ceed, &active_e_vec_out));
2044   CeedCallBackend(CeedDestroy(&ceed));
2045   CeedCallBackend(CeedQFunctionDestroy(&qf));
2046   return CEED_ERROR_SUCCESS;
2047 }
2048 
2049 //------------------------------------------------------------------------------
2050 // Create operator
2051 //------------------------------------------------------------------------------
2052 int CeedOperatorCreate_Cuda(CeedOperator op) {
2053   Ceed               ceed;
2054   CeedOperator_Cuda *impl;
2055 
2056   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
2057   CeedCallBackend(CeedCalloc(1, &impl));
2058   CeedCallBackend(CeedOperatorSetData(op, impl));
2059 
2060   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleQFunction", CeedOperatorLinearAssembleQFunction_Cuda));
2061   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleQFunctionUpdate", CeedOperatorLinearAssembleQFunctionUpdate_Cuda));
2062   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleAddDiagonal", CeedOperatorLinearAssembleAddDiagonal_Cuda));
2063   CeedCallBackend(
2064       CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleAddPointBlockDiagonal", CeedOperatorLinearAssembleAddPointBlockDiagonal_Cuda));
2065   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleSingle", CeedSingleOperatorAssemble_Cuda));
2066   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "ApplyAdd", CeedOperatorApplyAdd_Cuda));
2067   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "Destroy", CeedOperatorDestroy_Cuda));
2068   CeedCallBackend(CeedDestroy(&ceed));
2069   return CEED_ERROR_SUCCESS;
2070 }
2071 
2072 //------------------------------------------------------------------------------
2073 // Create operator AtPoints
2074 //------------------------------------------------------------------------------
2075 int CeedOperatorCreateAtPoints_Cuda(CeedOperator op) {
2076   Ceed               ceed;
2077   CeedOperator_Cuda *impl;
2078 
2079   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
2080   CeedCallBackend(CeedCalloc(1, &impl));
2081   CeedCallBackend(CeedOperatorSetData(op, impl));
2082 
2083   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleQFunction", CeedOperatorLinearAssembleQFunctionAtPoints_Cuda));
2084   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleAddDiagonal", CeedOperatorLinearAssembleAddDiagonalAtPoints_Cuda));
2085   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "ApplyAdd", CeedOperatorApplyAddAtPoints_Cuda));
2086   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "Destroy", CeedOperatorDestroy_Cuda));
2087   CeedCallBackend(CeedDestroy(&ceed));
2088   return CEED_ERROR_SUCCESS;
2089 }
2090 
2091 //------------------------------------------------------------------------------
2092