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