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