xref: /libCEED/backends/cuda-ref/ceed-cuda-ref-operator.c (revision c85e864000754611328ed43fe87e29f99787d593)
1 // Copyright (c) 2017-2022, 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   CeedCallBackend(CeedOperatorGetData(op, &impl));
27 
28   // Apply data
29   for (CeedInt i = 0; i < impl->numein + impl->numeout; i++) {
30     CeedCallBackend(CeedVectorDestroy(&impl->evecs[i]));
31   }
32   CeedCallBackend(CeedFree(&impl->evecs));
33 
34   for (CeedInt i = 0; i < impl->numein; i++) {
35     CeedCallBackend(CeedVectorDestroy(&impl->qvecsin[i]));
36   }
37   CeedCallBackend(CeedFree(&impl->qvecsin));
38 
39   for (CeedInt i = 0; i < impl->numeout; i++) {
40     CeedCallBackend(CeedVectorDestroy(&impl->qvecsout[i]));
41   }
42   CeedCallBackend(CeedFree(&impl->qvecsout));
43 
44   // QFunction assembly data
45   for (CeedInt i = 0; i < impl->qfnumactivein; i++) {
46     CeedCallBackend(CeedVectorDestroy(&impl->qfactivein[i]));
47   }
48   CeedCallBackend(CeedFree(&impl->qfactivein));
49 
50   // Diag data
51   if (impl->diag) {
52     Ceed ceed;
53     CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
54     CeedCallCuda(ceed, cuModuleUnload(impl->diag->module));
55     CeedCallBackend(CeedFree(&impl->diag->h_emodein));
56     CeedCallBackend(CeedFree(&impl->diag->h_emodeout));
57     CeedCallCuda(ceed, cudaFree(impl->diag->d_emodein));
58     CeedCallCuda(ceed, cudaFree(impl->diag->d_emodeout));
59     CeedCallCuda(ceed, cudaFree(impl->diag->d_identity));
60     CeedCallCuda(ceed, cudaFree(impl->diag->d_interpin));
61     CeedCallCuda(ceed, cudaFree(impl->diag->d_interpout));
62     CeedCallCuda(ceed, cudaFree(impl->diag->d_gradin));
63     CeedCallCuda(ceed, cudaFree(impl->diag->d_gradout));
64     CeedCallBackend(CeedElemRestrictionDestroy(&impl->diag->pbdiagrstr));
65     CeedCallBackend(CeedVectorDestroy(&impl->diag->elemdiag));
66     CeedCallBackend(CeedVectorDestroy(&impl->diag->pbelemdiag));
67   }
68   CeedCallBackend(CeedFree(&impl->diag));
69 
70   if (impl->asmb) {
71     Ceed ceed;
72     CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
73     CeedCallCuda(ceed, cuModuleUnload(impl->asmb->module));
74     CeedCallCuda(ceed, cudaFree(impl->asmb->d_B_in));
75     CeedCallCuda(ceed, cudaFree(impl->asmb->d_B_out));
76   }
77   CeedCallBackend(CeedFree(&impl->asmb));
78 
79   CeedCallBackend(CeedFree(&impl));
80   return CEED_ERROR_SUCCESS;
81 }
82 
83 //------------------------------------------------------------------------------
84 // Setup infields or outfields
85 //------------------------------------------------------------------------------
86 static int CeedOperatorSetupFields_Cuda(CeedQFunction qf, CeedOperator op, bool isinput, CeedVector *evecs, CeedVector *qvecs, CeedInt starte,
87                                         CeedInt numfields, CeedInt Q, CeedInt numelements) {
88   CeedInt  dim, size;
89   CeedSize q_size;
90   Ceed     ceed;
91   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
92   CeedBasis           basis;
93   CeedElemRestriction Erestrict;
94   CeedOperatorField  *opfields;
95   CeedQFunctionField *qffields;
96   CeedVector          fieldvec;
97   bool                strided;
98   bool                skiprestrict;
99 
100   if (isinput) {
101     CeedCallBackend(CeedOperatorGetFields(op, NULL, &opfields, NULL, NULL));
102     CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qffields, NULL, NULL));
103   } else {
104     CeedCallBackend(CeedOperatorGetFields(op, NULL, NULL, NULL, &opfields));
105     CeedCallBackend(CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qffields));
106   }
107 
108   // Loop over fields
109   for (CeedInt i = 0; i < numfields; i++) {
110     CeedEvalMode emode;
111     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qffields[i], &emode));
112 
113     strided      = false;
114     skiprestrict = false;
115     if (emode != CEED_EVAL_WEIGHT) {
116       CeedCallBackend(CeedOperatorFieldGetElemRestriction(opfields[i], &Erestrict));
117 
118       // Check whether this field can skip the element restriction:
119       // must be passive input, with emode NONE, and have a strided restriction with CEED_STRIDES_BACKEND.
120 
121       // First, check whether the field is input or output:
122       if (isinput) {
123         // Check for passive input:
124         CeedCallBackend(CeedOperatorFieldGetVector(opfields[i], &fieldvec));
125         if (fieldvec != CEED_VECTOR_ACTIVE) {
126           // Check emode
127           if (emode == CEED_EVAL_NONE) {
128             // Check for strided restriction
129             CeedCallBackend(CeedElemRestrictionIsStrided(Erestrict, &strided));
130             if (strided) {
131               // Check if vector is already in preferred backend ordering
132               CeedCallBackend(CeedElemRestrictionHasBackendStrides(Erestrict, &skiprestrict));
133             }
134           }
135         }
136       }
137       if (skiprestrict) {
138         // We do not need an E-Vector, but will use the input field vector's data directly in the operator application.
139         evecs[i + starte] = NULL;
140       } else {
141         CeedCallBackend(CeedElemRestrictionCreateVector(Erestrict, NULL, &evecs[i + starte]));
142       }
143     }
144 
145     switch (emode) {
146       case CEED_EVAL_NONE:
147         CeedCallBackend(CeedQFunctionFieldGetSize(qffields[i], &size));
148         q_size = (CeedSize)numelements * Q * size;
149         CeedCallBackend(CeedVectorCreate(ceed, q_size, &qvecs[i]));
150         break;
151       case CEED_EVAL_INTERP:
152         CeedCallBackend(CeedQFunctionFieldGetSize(qffields[i], &size));
153         q_size = (CeedSize)numelements * Q * size;
154         CeedCallBackend(CeedVectorCreate(ceed, q_size, &qvecs[i]));
155         break;
156       case CEED_EVAL_GRAD:
157         CeedCallBackend(CeedOperatorFieldGetBasis(opfields[i], &basis));
158         CeedCallBackend(CeedQFunctionFieldGetSize(qffields[i], &size));
159         CeedCallBackend(CeedBasisGetDimension(basis, &dim));
160         q_size = (CeedSize)numelements * Q * size;
161         CeedCallBackend(CeedVectorCreate(ceed, q_size, &qvecs[i]));
162         break;
163       case CEED_EVAL_WEIGHT:  // Only on input fields
164         CeedCallBackend(CeedOperatorFieldGetBasis(opfields[i], &basis));
165         q_size = (CeedSize)numelements * Q;
166         CeedCallBackend(CeedVectorCreate(ceed, q_size, &qvecs[i]));
167         CeedCallBackend(CeedBasisApply(basis, numelements, CEED_NOTRANSPOSE, CEED_EVAL_WEIGHT, NULL, qvecs[i]));
168         break;
169       case CEED_EVAL_DIV:
170         break;  // TODO: Not implemented
171       case CEED_EVAL_CURL:
172         break;  // TODO: Not implemented
173     }
174   }
175   return CEED_ERROR_SUCCESS;
176 }
177 
178 //------------------------------------------------------------------------------
179 // CeedOperator needs to connect all the named fields (be they active or passive) to the named inputs and outputs of its CeedQFunction.
180 //------------------------------------------------------------------------------
181 static int CeedOperatorSetup_Cuda(CeedOperator op) {
182   bool setupdone;
183   CeedCallBackend(CeedOperatorIsSetupDone(op, &setupdone));
184   if (setupdone) return CEED_ERROR_SUCCESS;
185   Ceed ceed;
186   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
187   CeedOperator_Cuda *impl;
188   CeedCallBackend(CeedOperatorGetData(op, &impl));
189   CeedQFunction qf;
190   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
191   CeedInt Q, numelements, numinputfields, numoutputfields;
192   CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q));
193   CeedCallBackend(CeedOperatorGetNumElements(op, &numelements));
194   CeedOperatorField *opinputfields, *opoutputfields;
195   CeedCallBackend(CeedOperatorGetFields(op, &numinputfields, &opinputfields, &numoutputfields, &opoutputfields));
196   CeedQFunctionField *qfinputfields, *qfoutputfields;
197   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qfinputfields, NULL, &qfoutputfields));
198 
199   // Allocate
200   CeedCallBackend(CeedCalloc(numinputfields + numoutputfields, &impl->evecs));
201 
202   CeedCallBackend(CeedCalloc(CEED_FIELD_MAX, &impl->qvecsin));
203   CeedCallBackend(CeedCalloc(CEED_FIELD_MAX, &impl->qvecsout));
204 
205   impl->numein  = numinputfields;
206   impl->numeout = numoutputfields;
207 
208   // Set up infield and outfield evecs and qvecs
209   // Infields
210   CeedCallBackend(CeedOperatorSetupFields_Cuda(qf, op, true, impl->evecs, impl->qvecsin, 0, numinputfields, Q, numelements));
211 
212   // Outfields
213   CeedCallBackend(CeedOperatorSetupFields_Cuda(qf, op, false, impl->evecs, impl->qvecsout, numinputfields, numoutputfields, Q, numelements));
214 
215   CeedCallBackend(CeedOperatorSetSetupDone(op));
216   return CEED_ERROR_SUCCESS;
217 }
218 
219 //------------------------------------------------------------------------------
220 // Setup Operator Inputs
221 //------------------------------------------------------------------------------
222 static inline int CeedOperatorSetupInputs_Cuda(CeedInt numinputfields, CeedQFunctionField *qfinputfields, CeedOperatorField *opinputfields,
223                                                CeedVector invec, const bool skipactive, CeedScalar *edata[2 * CEED_FIELD_MAX],
224                                                CeedOperator_Cuda *impl, CeedRequest *request) {
225   CeedEvalMode        emode;
226   CeedVector          vec;
227   CeedElemRestriction Erestrict;
228 
229   for (CeedInt i = 0; i < numinputfields; i++) {
230     // Get input vector
231     CeedCallBackend(CeedOperatorFieldGetVector(opinputfields[i], &vec));
232     if (vec == CEED_VECTOR_ACTIVE) {
233       if (skipactive) continue;
234       else vec = invec;
235     }
236 
237     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qfinputfields[i], &emode));
238     if (emode == CEED_EVAL_WEIGHT) {  // Skip
239     } else {
240       // Get input vector
241       CeedCallBackend(CeedOperatorFieldGetVector(opinputfields[i], &vec));
242       // Get input element restriction
243       CeedCallBackend(CeedOperatorFieldGetElemRestriction(opinputfields[i], &Erestrict));
244       if (vec == CEED_VECTOR_ACTIVE) vec = invec;
245       // Restrict, if necessary
246       if (!impl->evecs[i]) {
247         // No restriction for this field; read data directly from vec.
248         CeedCallBackend(CeedVectorGetArrayRead(vec, CEED_MEM_DEVICE, (const CeedScalar **)&edata[i]));
249       } else {
250         CeedCallBackend(CeedElemRestrictionApply(Erestrict, CEED_NOTRANSPOSE, vec, impl->evecs[i], request));
251         // Get evec
252         CeedCallBackend(CeedVectorGetArrayRead(impl->evecs[i], CEED_MEM_DEVICE, (const CeedScalar **)&edata[i]));
253       }
254     }
255   }
256   return CEED_ERROR_SUCCESS;
257 }
258 
259 //------------------------------------------------------------------------------
260 // Input Basis Action
261 //------------------------------------------------------------------------------
262 static inline int CeedOperatorInputBasis_Cuda(CeedInt numelements, CeedQFunctionField *qfinputfields, CeedOperatorField *opinputfields,
263                                               CeedInt numinputfields, const bool skipactive, CeedScalar *edata[2 * CEED_FIELD_MAX],
264                                               CeedOperator_Cuda *impl) {
265   CeedInt             elemsize, size;
266   CeedElemRestriction Erestrict;
267   CeedEvalMode        emode;
268   CeedBasis           basis;
269 
270   for (CeedInt i = 0; i < numinputfields; i++) {
271     // Skip active input
272     if (skipactive) {
273       CeedVector vec;
274       CeedCallBackend(CeedOperatorFieldGetVector(opinputfields[i], &vec));
275       if (vec == CEED_VECTOR_ACTIVE) continue;
276     }
277     // Get elemsize, emode, size
278     CeedCallBackend(CeedOperatorFieldGetElemRestriction(opinputfields[i], &Erestrict));
279     CeedCallBackend(CeedElemRestrictionGetElementSize(Erestrict, &elemsize));
280     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qfinputfields[i], &emode));
281     CeedCallBackend(CeedQFunctionFieldGetSize(qfinputfields[i], &size));
282     // Basis action
283     switch (emode) {
284       case CEED_EVAL_NONE:
285         CeedCallBackend(CeedVectorSetArray(impl->qvecsin[i], CEED_MEM_DEVICE, CEED_USE_POINTER, edata[i]));
286         break;
287       case CEED_EVAL_INTERP:
288         CeedCallBackend(CeedOperatorFieldGetBasis(opinputfields[i], &basis));
289         CeedCallBackend(CeedBasisApply(basis, numelements, CEED_NOTRANSPOSE, CEED_EVAL_INTERP, impl->evecs[i], impl->qvecsin[i]));
290         break;
291       case CEED_EVAL_GRAD:
292         CeedCallBackend(CeedOperatorFieldGetBasis(opinputfields[i], &basis));
293         CeedCallBackend(CeedBasisApply(basis, numelements, CEED_NOTRANSPOSE, CEED_EVAL_GRAD, impl->evecs[i], impl->qvecsin[i]));
294         break;
295       case CEED_EVAL_WEIGHT:
296         break;  // No action
297       case CEED_EVAL_DIV:
298         break;  // TODO: Not implemented
299       case CEED_EVAL_CURL:
300         break;  // TODO: Not implemented
301     }
302   }
303   return CEED_ERROR_SUCCESS;
304 }
305 
306 //------------------------------------------------------------------------------
307 // Restore Input Vectors
308 //------------------------------------------------------------------------------
309 static inline int CeedOperatorRestoreInputs_Cuda(CeedInt numinputfields, CeedQFunctionField *qfinputfields, CeedOperatorField *opinputfields,
310                                                  const bool skipactive, CeedScalar *edata[2 * CEED_FIELD_MAX], CeedOperator_Cuda *impl) {
311   CeedEvalMode emode;
312   CeedVector   vec;
313 
314   for (CeedInt i = 0; i < numinputfields; i++) {
315     // Skip active input
316     if (skipactive) {
317       CeedCallBackend(CeedOperatorFieldGetVector(opinputfields[i], &vec));
318       if (vec == CEED_VECTOR_ACTIVE) continue;
319     }
320     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qfinputfields[i], &emode));
321     if (emode == CEED_EVAL_WEIGHT) {  // Skip
322     } else {
323       if (!impl->evecs[i]) {  // This was a skiprestrict case
324         CeedCallBackend(CeedOperatorFieldGetVector(opinputfields[i], &vec));
325         CeedCallBackend(CeedVectorRestoreArrayRead(vec, (const CeedScalar **)&edata[i]));
326       } else {
327         CeedCallBackend(CeedVectorRestoreArrayRead(impl->evecs[i], (const CeedScalar **)&edata[i]));
328       }
329     }
330   }
331   return CEED_ERROR_SUCCESS;
332 }
333 
334 //------------------------------------------------------------------------------
335 // Apply and add to output
336 //------------------------------------------------------------------------------
337 static int CeedOperatorApplyAdd_Cuda(CeedOperator op, CeedVector invec, CeedVector outvec, CeedRequest *request) {
338   CeedOperator_Cuda *impl;
339   CeedCallBackend(CeedOperatorGetData(op, &impl));
340   CeedQFunction qf;
341   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
342   CeedInt Q, numelements, elemsize, numinputfields, numoutputfields, size;
343   CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q));
344   CeedCallBackend(CeedOperatorGetNumElements(op, &numelements));
345   CeedOperatorField *opinputfields, *opoutputfields;
346   CeedCallBackend(CeedOperatorGetFields(op, &numinputfields, &opinputfields, &numoutputfields, &opoutputfields));
347   CeedQFunctionField *qfinputfields, *qfoutputfields;
348   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qfinputfields, NULL, &qfoutputfields));
349   CeedEvalMode        emode;
350   CeedVector          vec;
351   CeedBasis           basis;
352   CeedElemRestriction Erestrict;
353   CeedScalar         *edata[2 * CEED_FIELD_MAX] = {0};
354 
355   // Setup
356   CeedCallBackend(CeedOperatorSetup_Cuda(op));
357 
358   // Input Evecs and Restriction
359   CeedCallBackend(CeedOperatorSetupInputs_Cuda(numinputfields, qfinputfields, opinputfields, invec, false, edata, impl, request));
360 
361   // Input basis apply if needed
362   CeedCallBackend(CeedOperatorInputBasis_Cuda(numelements, qfinputfields, opinputfields, numinputfields, false, edata, impl));
363 
364   // Output pointers, as necessary
365   for (CeedInt i = 0; i < numoutputfields; i++) {
366     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qfoutputfields[i], &emode));
367     if (emode == CEED_EVAL_NONE) {
368       // Set the output Q-Vector to use the E-Vector data directly.
369       CeedCallBackend(CeedVectorGetArrayWrite(impl->evecs[i + impl->numein], CEED_MEM_DEVICE, &edata[i + numinputfields]));
370       CeedCallBackend(CeedVectorSetArray(impl->qvecsout[i], CEED_MEM_DEVICE, CEED_USE_POINTER, edata[i + numinputfields]));
371     }
372   }
373 
374   // Q function
375   CeedCallBackend(CeedQFunctionApply(qf, numelements * Q, impl->qvecsin, impl->qvecsout));
376 
377   // Output basis apply if needed
378   for (CeedInt i = 0; i < numoutputfields; i++) {
379     // Get elemsize, emode, size
380     CeedCallBackend(CeedOperatorFieldGetElemRestriction(opoutputfields[i], &Erestrict));
381     CeedCallBackend(CeedElemRestrictionGetElementSize(Erestrict, &elemsize));
382     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qfoutputfields[i], &emode));
383     CeedCallBackend(CeedQFunctionFieldGetSize(qfoutputfields[i], &size));
384     // Basis action
385     switch (emode) {
386       case CEED_EVAL_NONE:
387         break;
388       case CEED_EVAL_INTERP:
389         CeedCallBackend(CeedOperatorFieldGetBasis(opoutputfields[i], &basis));
390         CeedCallBackend(CeedBasisApply(basis, numelements, CEED_TRANSPOSE, CEED_EVAL_INTERP, impl->qvecsout[i], impl->evecs[i + impl->numein]));
391         break;
392       case CEED_EVAL_GRAD:
393         CeedCallBackend(CeedOperatorFieldGetBasis(opoutputfields[i], &basis));
394         CeedCallBackend(CeedBasisApply(basis, numelements, CEED_TRANSPOSE, CEED_EVAL_GRAD, impl->qvecsout[i], impl->evecs[i + impl->numein]));
395         break;
396       // LCOV_EXCL_START
397       case CEED_EVAL_WEIGHT: {
398         Ceed ceed;
399         CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
400         return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_WEIGHT cannot be an output evaluation mode");
401         break;  // Should not occur
402       }
403       case CEED_EVAL_DIV:
404         break;  // TODO: Not implemented
405       case CEED_EVAL_CURL:
406         break;  // TODO: Not implemented
407                 // LCOV_EXCL_STOP
408     }
409   }
410 
411   // Output restriction
412   for (CeedInt i = 0; i < numoutputfields; i++) {
413     // Restore evec
414     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qfoutputfields[i], &emode));
415     if (emode == CEED_EVAL_NONE) {
416       CeedCallBackend(CeedVectorRestoreArray(impl->evecs[i + impl->numein], &edata[i + numinputfields]));
417     }
418     // Get output vector
419     CeedCallBackend(CeedOperatorFieldGetVector(opoutputfields[i], &vec));
420     // Restrict
421     CeedCallBackend(CeedOperatorFieldGetElemRestriction(opoutputfields[i], &Erestrict));
422     // Active
423     if (vec == CEED_VECTOR_ACTIVE) vec = outvec;
424 
425     CeedCallBackend(CeedElemRestrictionApply(Erestrict, CEED_TRANSPOSE, impl->evecs[i + impl->numein], vec, request));
426   }
427 
428   // Restore input arrays
429   CeedCallBackend(CeedOperatorRestoreInputs_Cuda(numinputfields, qfinputfields, opinputfields, false, edata, impl));
430   return CEED_ERROR_SUCCESS;
431 }
432 
433 //------------------------------------------------------------------------------
434 // Core code for assembling linear QFunction
435 //------------------------------------------------------------------------------
436 static inline int CeedOperatorLinearAssembleQFunctionCore_Cuda(CeedOperator op, bool build_objects, CeedVector *assembled, CeedElemRestriction *rstr,
437                                                                CeedRequest *request) {
438   CeedOperator_Cuda *impl;
439   CeedCallBackend(CeedOperatorGetData(op, &impl));
440   CeedQFunction qf;
441   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
442   CeedInt  Q, numelements, numinputfields, numoutputfields, size;
443   CeedSize q_size;
444   CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q));
445   CeedCallBackend(CeedOperatorGetNumElements(op, &numelements));
446   CeedOperatorField *opinputfields, *opoutputfields;
447   CeedCallBackend(CeedOperatorGetFields(op, &numinputfields, &opinputfields, &numoutputfields, &opoutputfields));
448   CeedQFunctionField *qfinputfields, *qfoutputfields;
449   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qfinputfields, NULL, &qfoutputfields));
450   CeedVector  vec;
451   CeedInt     numactivein = impl->qfnumactivein, numactiveout = impl->qfnumactiveout;
452   CeedVector *activein = impl->qfactivein;
453   CeedScalar *a, *tmp;
454   Ceed        ceed, ceedparent;
455   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
456   CeedCallBackend(CeedGetOperatorFallbackParentCeed(ceed, &ceedparent));
457   ceedparent = ceedparent ? ceedparent : ceed;
458   CeedScalar *edata[2 * CEED_FIELD_MAX];
459 
460   // Setup
461   CeedCallBackend(CeedOperatorSetup_Cuda(op));
462 
463   // Check for identity
464   bool identityqf;
465   CeedCallBackend(CeedQFunctionIsIdentity(qf, &identityqf));
466   if (identityqf) {
467     // LCOV_EXCL_START
468     return CeedError(ceed, CEED_ERROR_BACKEND, "Assembling identity QFunctions not supported");
469     // LCOV_EXCL_STOP
470   }
471 
472   // Input Evecs and Restriction
473   CeedCallBackend(CeedOperatorSetupInputs_Cuda(numinputfields, qfinputfields, opinputfields, NULL, true, edata, impl, request));
474 
475   // Count number of active input fields
476   if (!numactivein) {
477     for (CeedInt i = 0; i < numinputfields; i++) {
478       // Get input vector
479       CeedCallBackend(CeedOperatorFieldGetVector(opinputfields[i], &vec));
480       // Check if active input
481       if (vec == CEED_VECTOR_ACTIVE) {
482         CeedCallBackend(CeedQFunctionFieldGetSize(qfinputfields[i], &size));
483         CeedCallBackend(CeedVectorSetValue(impl->qvecsin[i], 0.0));
484         CeedCallBackend(CeedVectorGetArray(impl->qvecsin[i], CEED_MEM_DEVICE, &tmp));
485         CeedCallBackend(CeedRealloc(numactivein + size, &activein));
486         for (CeedInt field = 0; field < size; field++) {
487           q_size = (CeedSize)Q * numelements;
488           CeedCallBackend(CeedVectorCreate(ceed, q_size, &activein[numactivein + field]));
489           CeedCallBackend(CeedVectorSetArray(activein[numactivein + field], CEED_MEM_DEVICE, CEED_USE_POINTER, &tmp[field * Q * numelements]));
490         }
491         numactivein += size;
492         CeedCallBackend(CeedVectorRestoreArray(impl->qvecsin[i], &tmp));
493       }
494     }
495     impl->qfnumactivein = numactivein;
496     impl->qfactivein    = activein;
497   }
498 
499   // Count number of active output fields
500   if (!numactiveout) {
501     for (CeedInt i = 0; i < numoutputfields; i++) {
502       // Get output vector
503       CeedCallBackend(CeedOperatorFieldGetVector(opoutputfields[i], &vec));
504       // Check if active output
505       if (vec == CEED_VECTOR_ACTIVE) {
506         CeedCallBackend(CeedQFunctionFieldGetSize(qfoutputfields[i], &size));
507         numactiveout += size;
508       }
509     }
510     impl->qfnumactiveout = numactiveout;
511   }
512 
513   // Check sizes
514   if (!numactivein || !numactiveout) {
515     // LCOV_EXCL_START
516     return CeedError(ceed, CEED_ERROR_BACKEND, "Cannot assemble QFunction without active inputs and outputs");
517     // LCOV_EXCL_STOP
518   }
519 
520   // Build objects if needed
521   if (build_objects) {
522     // Create output restriction
523     CeedInt strides[3] = {1, numelements * Q, Q}; /* *NOPAD* */
524     CeedCallBackend(CeedElemRestrictionCreateStrided(ceedparent, numelements, Q, numactivein * numactiveout,
525                                                      numactivein * numactiveout * numelements * Q, strides, rstr));
526     // Create assembled vector
527     CeedSize l_size = (CeedSize)numelements * Q * numactivein * numactiveout;
528     CeedCallBackend(CeedVectorCreate(ceedparent, l_size, assembled));
529   }
530   CeedCallBackend(CeedVectorSetValue(*assembled, 0.0));
531   CeedCallBackend(CeedVectorGetArray(*assembled, CEED_MEM_DEVICE, &a));
532 
533   // Input basis apply
534   CeedCallBackend(CeedOperatorInputBasis_Cuda(numelements, qfinputfields, opinputfields, numinputfields, true, edata, impl));
535 
536   // Assemble QFunction
537   for (CeedInt in = 0; in < numactivein; in++) {
538     // Set Inputs
539     CeedCallBackend(CeedVectorSetValue(activein[in], 1.0));
540     if (numactivein > 1) {
541       CeedCallBackend(CeedVectorSetValue(activein[(in + numactivein - 1) % numactivein], 0.0));
542     }
543     // Set Outputs
544     for (CeedInt out = 0; out < numoutputfields; out++) {
545       // Get output vector
546       CeedCallBackend(CeedOperatorFieldGetVector(opoutputfields[out], &vec));
547       // Check if active output
548       if (vec == CEED_VECTOR_ACTIVE) {
549         CeedCallBackend(CeedVectorSetArray(impl->qvecsout[out], CEED_MEM_DEVICE, CEED_USE_POINTER, a));
550         CeedCallBackend(CeedQFunctionFieldGetSize(qfoutputfields[out], &size));
551         a += size * Q * numelements;  // Advance the pointer by the size of the output
552       }
553     }
554     // Apply QFunction
555     CeedCallBackend(CeedQFunctionApply(qf, Q * numelements, impl->qvecsin, impl->qvecsout));
556   }
557 
558   // Un-set output Qvecs to prevent accidental overwrite of Assembled
559   for (CeedInt out = 0; out < numoutputfields; out++) {
560     // Get output vector
561     CeedCallBackend(CeedOperatorFieldGetVector(opoutputfields[out], &vec));
562     // Check if active output
563     if (vec == CEED_VECTOR_ACTIVE) {
564       CeedCallBackend(CeedVectorTakeArray(impl->qvecsout[out], CEED_MEM_DEVICE, NULL));
565     }
566   }
567 
568   // Restore input arrays
569   CeedCallBackend(CeedOperatorRestoreInputs_Cuda(numinputfields, qfinputfields, opinputfields, true, edata, impl));
570 
571   // Restore output
572   CeedCallBackend(CeedVectorRestoreArray(*assembled, &a));
573 
574   return CEED_ERROR_SUCCESS;
575 }
576 
577 //------------------------------------------------------------------------------
578 // Assemble Linear QFunction
579 //------------------------------------------------------------------------------
580 static int CeedOperatorLinearAssembleQFunction_Cuda(CeedOperator op, CeedVector *assembled, CeedElemRestriction *rstr, CeedRequest *request) {
581   return CeedOperatorLinearAssembleQFunctionCore_Cuda(op, true, assembled, rstr, request);
582 }
583 
584 //------------------------------------------------------------------------------
585 // Update Assembled Linear QFunction
586 //------------------------------------------------------------------------------
587 static int CeedOperatorLinearAssembleQFunctionUpdate_Cuda(CeedOperator op, CeedVector assembled, CeedElemRestriction rstr, CeedRequest *request) {
588   return CeedOperatorLinearAssembleQFunctionCore_Cuda(op, false, &assembled, &rstr, request);
589 }
590 
591 //------------------------------------------------------------------------------
592 // Create point block restriction
593 //------------------------------------------------------------------------------
594 static int CreatePBRestriction(CeedElemRestriction rstr, CeedElemRestriction *pbRstr) {
595   Ceed ceed;
596   CeedCallBackend(CeedElemRestrictionGetCeed(rstr, &ceed));
597   const CeedInt *offsets;
598   CeedCallBackend(CeedElemRestrictionGetOffsets(rstr, CEED_MEM_HOST, &offsets));
599 
600   // Expand offsets
601   CeedInt  nelem, ncomp, elemsize, compstride, *pbOffsets;
602   CeedSize l_size;
603   CeedCallBackend(CeedElemRestrictionGetNumElements(rstr, &nelem));
604   CeedCallBackend(CeedElemRestrictionGetNumComponents(rstr, &ncomp));
605   CeedCallBackend(CeedElemRestrictionGetElementSize(rstr, &elemsize));
606   CeedCallBackend(CeedElemRestrictionGetCompStride(rstr, &compstride));
607   CeedCallBackend(CeedElemRestrictionGetLVectorSize(rstr, &l_size));
608   CeedInt shift = ncomp;
609   if (compstride != 1) shift *= ncomp;
610   CeedCallBackend(CeedCalloc(nelem * elemsize, &pbOffsets));
611   for (CeedInt i = 0; i < nelem * elemsize; i++) {
612     pbOffsets[i] = offsets[i] * shift;
613   }
614 
615   // Create new restriction
616   CeedCallBackend(
617       CeedElemRestrictionCreate(ceed, nelem, elemsize, ncomp * ncomp, 1, l_size * ncomp, CEED_MEM_HOST, CEED_OWN_POINTER, pbOffsets, pbRstr));
618 
619   // Cleanup
620   CeedCallBackend(CeedElemRestrictionRestoreOffsets(rstr, &offsets));
621 
622   return CEED_ERROR_SUCCESS;
623 }
624 
625 //------------------------------------------------------------------------------
626 // Assemble diagonal setup
627 //------------------------------------------------------------------------------
628 static inline int CeedOperatorAssembleDiagonalSetup_Cuda(CeedOperator op, const bool pointBlock) {
629   Ceed ceed;
630   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
631   CeedQFunction qf;
632   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
633   CeedInt numinputfields, numoutputfields;
634   CeedCallBackend(CeedQFunctionGetNumArgs(qf, &numinputfields, &numoutputfields));
635 
636   // Determine active input basis
637   CeedOperatorField  *opfields;
638   CeedQFunctionField *qffields;
639   CeedCallBackend(CeedOperatorGetFields(op, NULL, &opfields, NULL, NULL));
640   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qffields, NULL, NULL));
641   CeedInt             numemodein = 0, ncomp = 0, dim = 1;
642   CeedEvalMode       *emodein = NULL;
643   CeedBasis           basisin = NULL;
644   CeedElemRestriction rstrin  = NULL;
645   for (CeedInt i = 0; i < numinputfields; i++) {
646     CeedVector vec;
647     CeedCallBackend(CeedOperatorFieldGetVector(opfields[i], &vec));
648     if (vec == CEED_VECTOR_ACTIVE) {
649       CeedElemRestriction rstr;
650       CeedCallBackend(CeedOperatorFieldGetBasis(opfields[i], &basisin));
651       CeedCallBackend(CeedBasisGetNumComponents(basisin, &ncomp));
652       CeedCallBackend(CeedBasisGetDimension(basisin, &dim));
653       CeedCallBackend(CeedOperatorFieldGetElemRestriction(opfields[i], &rstr));
654       if (rstrin && rstrin != rstr) {
655         // LCOV_EXCL_START
656         return CeedError(ceed, CEED_ERROR_BACKEND, "Backend does not implement multi-field non-composite operator diagonal assembly");
657         // LCOV_EXCL_STOP
658       }
659       rstrin = rstr;
660       CeedEvalMode emode;
661       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qffields[i], &emode));
662       switch (emode) {
663         case CEED_EVAL_NONE:
664         case CEED_EVAL_INTERP:
665           CeedCallBackend(CeedRealloc(numemodein + 1, &emodein));
666           emodein[numemodein] = emode;
667           numemodein += 1;
668           break;
669         case CEED_EVAL_GRAD:
670           CeedCallBackend(CeedRealloc(numemodein + dim, &emodein));
671           for (CeedInt d = 0; d < dim; d++) emodein[numemodein + d] = emode;
672           numemodein += dim;
673           break;
674         case CEED_EVAL_WEIGHT:
675         case CEED_EVAL_DIV:
676         case CEED_EVAL_CURL:
677           break;  // Caught by QF Assembly
678       }
679     }
680   }
681 
682   // Determine active output basis
683   CeedCallBackend(CeedOperatorGetFields(op, NULL, NULL, NULL, &opfields));
684   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qffields));
685   CeedInt             numemodeout = 0;
686   CeedEvalMode       *emodeout    = NULL;
687   CeedBasis           basisout    = NULL;
688   CeedElemRestriction rstrout     = NULL;
689   for (CeedInt i = 0; i < numoutputfields; i++) {
690     CeedVector vec;
691     CeedCallBackend(CeedOperatorFieldGetVector(opfields[i], &vec));
692     if (vec == CEED_VECTOR_ACTIVE) {
693       CeedElemRestriction rstr;
694       CeedCallBackend(CeedOperatorFieldGetBasis(opfields[i], &basisout));
695       CeedCallBackend(CeedOperatorFieldGetElemRestriction(opfields[i], &rstr));
696       if (rstrout && rstrout != rstr) {
697         // LCOV_EXCL_START
698         return CeedError(ceed, CEED_ERROR_BACKEND, "Backend does not implement multi-field non-composite operator diagonal assembly");
699         // LCOV_EXCL_STOP
700       }
701       rstrout = rstr;
702       CeedEvalMode emode;
703       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qffields[i], &emode));
704       switch (emode) {
705         case CEED_EVAL_NONE:
706         case CEED_EVAL_INTERP:
707           CeedCallBackend(CeedRealloc(numemodeout + 1, &emodeout));
708           emodeout[numemodeout] = emode;
709           numemodeout += 1;
710           break;
711         case CEED_EVAL_GRAD:
712           CeedCallBackend(CeedRealloc(numemodeout + dim, &emodeout));
713           for (CeedInt d = 0; d < dim; d++) emodeout[numemodeout + d] = emode;
714           numemodeout += dim;
715           break;
716         case CEED_EVAL_WEIGHT:
717         case CEED_EVAL_DIV:
718         case CEED_EVAL_CURL:
719           break;  // Caught by QF Assembly
720       }
721     }
722   }
723 
724   // Operator data struct
725   CeedOperator_Cuda *impl;
726   CeedCallBackend(CeedOperatorGetData(op, &impl));
727   CeedCallBackend(CeedCalloc(1, &impl->diag));
728   CeedOperatorDiag_Cuda *diag = impl->diag;
729   diag->basisin               = basisin;
730   diag->basisout              = basisout;
731   diag->h_emodein             = emodein;
732   diag->h_emodeout            = emodeout;
733   diag->numemodein            = numemodein;
734   diag->numemodeout           = numemodeout;
735 
736   // Assemble kernel
737   char *diagonal_kernel_path, *diagonal_kernel_source;
738   CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-operator-assemble-diagonal.h", &diagonal_kernel_path));
739   CeedDebug256(ceed, 2, "----- Loading Diagonal Assembly Kernel Source -----\n");
740   CeedCallBackend(CeedLoadSourceToBuffer(ceed, diagonal_kernel_path, &diagonal_kernel_source));
741   CeedDebug256(ceed, 2, "----- Loading Diagonal Assembly Source Complete! -----\n");
742   CeedInt nnodes, nqpts;
743   CeedCallBackend(CeedBasisGetNumNodes(basisin, &nnodes));
744   CeedCallBackend(CeedBasisGetNumQuadraturePoints(basisin, &nqpts));
745   diag->nnodes = nnodes;
746   CeedCallCuda(ceed, CeedCompileCuda(ceed, diagonal_kernel_source, &diag->module, 5, "NUMEMODEIN", numemodein, "NUMEMODEOUT", numemodeout, "NNODES",
747                                      nnodes, "NQPTS", nqpts, "NCOMP", ncomp));
748   CeedCallCuda(ceed, CeedGetKernelCuda(ceed, diag->module, "linearDiagonal", &diag->linearDiagonal));
749   CeedCallCuda(ceed, CeedGetKernelCuda(ceed, diag->module, "linearPointBlockDiagonal", &diag->linearPointBlock));
750   CeedCallBackend(CeedFree(&diagonal_kernel_path));
751   CeedCallBackend(CeedFree(&diagonal_kernel_source));
752 
753   // Basis matrices
754   const CeedInt     qBytes = nqpts * sizeof(CeedScalar);
755   const CeedInt     iBytes = qBytes * nnodes;
756   const CeedInt     gBytes = qBytes * nnodes * dim;
757   const CeedInt     eBytes = sizeof(CeedEvalMode);
758   const CeedScalar *interpin, *interpout, *gradin, *gradout;
759 
760   // CEED_EVAL_NONE
761   CeedScalar *identity = NULL;
762   bool        evalNone = false;
763   for (CeedInt i = 0; i < numemodein; i++) evalNone = evalNone || (emodein[i] == CEED_EVAL_NONE);
764   for (CeedInt i = 0; i < numemodeout; i++) evalNone = evalNone || (emodeout[i] == CEED_EVAL_NONE);
765   if (evalNone) {
766     CeedCallBackend(CeedCalloc(nqpts * nnodes, &identity));
767     for (CeedInt i = 0; i < (nnodes < nqpts ? nnodes : nqpts); i++) identity[i * nnodes + i] = 1.0;
768     CeedCallCuda(ceed, cudaMalloc((void **)&diag->d_identity, iBytes));
769     CeedCallCuda(ceed, cudaMemcpy(diag->d_identity, identity, iBytes, cudaMemcpyHostToDevice));
770   }
771 
772   // CEED_EVAL_INTERP
773   CeedCallBackend(CeedBasisGetInterp(basisin, &interpin));
774   CeedCallCuda(ceed, cudaMalloc((void **)&diag->d_interpin, iBytes));
775   CeedCallCuda(ceed, cudaMemcpy(diag->d_interpin, interpin, iBytes, cudaMemcpyHostToDevice));
776   CeedCallBackend(CeedBasisGetInterp(basisout, &interpout));
777   CeedCallCuda(ceed, cudaMalloc((void **)&diag->d_interpout, iBytes));
778   CeedCallCuda(ceed, cudaMemcpy(diag->d_interpout, interpout, iBytes, cudaMemcpyHostToDevice));
779 
780   // CEED_EVAL_GRAD
781   CeedCallBackend(CeedBasisGetGrad(basisin, &gradin));
782   CeedCallCuda(ceed, cudaMalloc((void **)&diag->d_gradin, gBytes));
783   CeedCallCuda(ceed, cudaMemcpy(diag->d_gradin, gradin, gBytes, cudaMemcpyHostToDevice));
784   CeedCallBackend(CeedBasisGetGrad(basisout, &gradout));
785   CeedCallCuda(ceed, cudaMalloc((void **)&diag->d_gradout, gBytes));
786   CeedCallCuda(ceed, cudaMemcpy(diag->d_gradout, gradout, gBytes, cudaMemcpyHostToDevice));
787 
788   // Arrays of emodes
789   CeedCallCuda(ceed, cudaMalloc((void **)&diag->d_emodein, numemodein * eBytes));
790   CeedCallCuda(ceed, cudaMemcpy(diag->d_emodein, emodein, numemodein * eBytes, cudaMemcpyHostToDevice));
791   CeedCallCuda(ceed, cudaMalloc((void **)&diag->d_emodeout, numemodeout * eBytes));
792   CeedCallCuda(ceed, cudaMemcpy(diag->d_emodeout, emodeout, numemodeout * eBytes, cudaMemcpyHostToDevice));
793 
794   // Restriction
795   diag->diagrstr = rstrout;
796 
797   return CEED_ERROR_SUCCESS;
798 }
799 
800 //------------------------------------------------------------------------------
801 // Assemble diagonal common code
802 //------------------------------------------------------------------------------
803 static inline int CeedOperatorAssembleDiagonalCore_Cuda(CeedOperator op, CeedVector assembled, CeedRequest *request, const bool pointBlock) {
804   Ceed ceed;
805   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
806   CeedOperator_Cuda *impl;
807   CeedCallBackend(CeedOperatorGetData(op, &impl));
808 
809   // Assemble QFunction
810   CeedVector          assembledqf;
811   CeedElemRestriction rstr;
812   CeedCallBackend(CeedOperatorLinearAssembleQFunctionBuildOrUpdate(op, &assembledqf, &rstr, request));
813   CeedCallBackend(CeedElemRestrictionDestroy(&rstr));
814 
815   // Setup
816   if (!impl->diag) {
817     CeedCallBackend(CeedOperatorAssembleDiagonalSetup_Cuda(op, pointBlock));
818   }
819   CeedOperatorDiag_Cuda *diag = impl->diag;
820   assert(diag != NULL);
821 
822   // Restriction
823   if (pointBlock && !diag->pbdiagrstr) {
824     CeedElemRestriction pbdiagrstr;
825     CeedCallBackend(CreatePBRestriction(diag->diagrstr, &pbdiagrstr));
826     diag->pbdiagrstr = pbdiagrstr;
827   }
828   CeedElemRestriction diagrstr = pointBlock ? diag->pbdiagrstr : diag->diagrstr;
829 
830   // Create diagonal vector
831   CeedVector elemdiag = pointBlock ? diag->pbelemdiag : diag->elemdiag;
832   if (!elemdiag) {
833     CeedCallBackend(CeedElemRestrictionCreateVector(diagrstr, NULL, &elemdiag));
834     if (pointBlock) diag->pbelemdiag = elemdiag;
835     else diag->elemdiag = elemdiag;
836   }
837   CeedCallBackend(CeedVectorSetValue(elemdiag, 0.0));
838 
839   // Assemble element operator diagonals
840   CeedScalar       *elemdiagarray;
841   const CeedScalar *assembledqfarray;
842   CeedCallBackend(CeedVectorGetArray(elemdiag, CEED_MEM_DEVICE, &elemdiagarray));
843   CeedCallBackend(CeedVectorGetArrayRead(assembledqf, CEED_MEM_DEVICE, &assembledqfarray));
844   CeedInt nelem;
845   CeedCallBackend(CeedElemRestrictionGetNumElements(diagrstr, &nelem));
846 
847   // Compute the diagonal of B^T D B
848   int   elemsPerBlock = 1;
849   int   grid          = nelem / elemsPerBlock + ((nelem / elemsPerBlock * elemsPerBlock < nelem) ? 1 : 0);
850   void *args[]        = {(void *)&nelem,   &diag->d_identity, &diag->d_interpin, &diag->d_gradin,   &diag->d_interpout,
851                          &diag->d_gradout, &diag->d_emodein,  &diag->d_emodeout, &assembledqfarray, &elemdiagarray};
852   if (pointBlock) {
853     CeedCallBackend(CeedRunKernelDimCuda(ceed, diag->linearPointBlock, grid, diag->nnodes, 1, elemsPerBlock, args));
854   } else {
855     CeedCallBackend(CeedRunKernelDimCuda(ceed, diag->linearDiagonal, grid, diag->nnodes, 1, elemsPerBlock, args));
856   }
857 
858   // Restore arrays
859   CeedCallBackend(CeedVectorRestoreArray(elemdiag, &elemdiagarray));
860   CeedCallBackend(CeedVectorRestoreArrayRead(assembledqf, &assembledqfarray));
861 
862   // Assemble local operator diagonal
863   CeedCallBackend(CeedElemRestrictionApply(diagrstr, CEED_TRANSPOSE, elemdiag, assembled, request));
864 
865   // Cleanup
866   CeedCallBackend(CeedVectorDestroy(&assembledqf));
867 
868   return CEED_ERROR_SUCCESS;
869 }
870 
871 //------------------------------------------------------------------------------
872 // Assemble Linear Diagonal
873 //------------------------------------------------------------------------------
874 static int CeedOperatorLinearAssembleAddDiagonal_Cuda(CeedOperator op, CeedVector assembled, CeedRequest *request) {
875   CeedCallBackend(CeedOperatorAssembleDiagonalCore_Cuda(op, assembled, request, false));
876   return CEED_ERROR_SUCCESS;
877 }
878 
879 //------------------------------------------------------------------------------
880 // Assemble Linear Point Block Diagonal
881 //------------------------------------------------------------------------------
882 static int CeedOperatorLinearAssembleAddPointBlockDiagonal_Cuda(CeedOperator op, CeedVector assembled, CeedRequest *request) {
883   CeedCallBackend(CeedOperatorAssembleDiagonalCore_Cuda(op, assembled, request, true));
884   return CEED_ERROR_SUCCESS;
885 }
886 
887 //------------------------------------------------------------------------------
888 // Single operator assembly setup
889 //------------------------------------------------------------------------------
890 static int CeedSingleOperatorAssembleSetup_Cuda(CeedOperator op) {
891   Ceed ceed;
892   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
893   CeedOperator_Cuda *impl;
894   CeedCallBackend(CeedOperatorGetData(op, &impl));
895 
896   // Get intput and output fields
897   CeedInt            num_input_fields, num_output_fields;
898   CeedOperatorField *input_fields;
899   CeedOperatorField *output_fields;
900   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &input_fields, &num_output_fields, &output_fields));
901 
902   // Determine active input basis eval mode
903   CeedQFunction qf;
904   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
905   CeedQFunctionField *qf_fields;
906   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_fields, NULL, NULL));
907   // Note that the kernel will treat each dimension of a gradient action separately;
908   // i.e., when an active input has a CEED_EVAL_GRAD mode, num_emode_in will increment by dim.
909   // However, for the purposes of loading the B matrices, it will be treated as one mode, and we will load/copy the entire gradient matrix at once, so
910   // num_B_in_mats_to_load will be incremented by 1.
911   CeedInt             num_emode_in = 0, dim = 1, num_B_in_mats_to_load = 0, size_B_in = 0;
912   CeedEvalMode       *eval_mode_in = NULL;  // will be of size num_B_in_mats_load
913   CeedBasis           basis_in     = NULL;
914   CeedInt             nqpts = 0, esize = 0;
915   CeedElemRestriction rstr_in = NULL;
916   for (CeedInt i = 0; i < num_input_fields; i++) {
917     CeedVector vec;
918     CeedCallBackend(CeedOperatorFieldGetVector(input_fields[i], &vec));
919     if (vec == CEED_VECTOR_ACTIVE) {
920       CeedCallBackend(CeedOperatorFieldGetBasis(input_fields[i], &basis_in));
921       CeedCallBackend(CeedBasisGetDimension(basis_in, &dim));
922       CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis_in, &nqpts));
923       CeedCallBackend(CeedOperatorFieldGetElemRestriction(input_fields[i], &rstr_in));
924       CeedCallBackend(CeedElemRestrictionGetElementSize(rstr_in, &esize));
925       CeedEvalMode eval_mode;
926       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode));
927       if (eval_mode != CEED_EVAL_NONE) {
928         CeedCallBackend(CeedRealloc(num_B_in_mats_to_load + 1, &eval_mode_in));
929         eval_mode_in[num_B_in_mats_to_load] = eval_mode;
930         num_B_in_mats_to_load += 1;
931         if (eval_mode == CEED_EVAL_GRAD) {
932           num_emode_in += dim;
933           size_B_in += dim * esize * nqpts;
934         } else {
935           num_emode_in += 1;
936           size_B_in += esize * nqpts;
937         }
938       }
939     }
940   }
941 
942   // Determine active output basis; basis_out and rstr_out only used if same as input, TODO
943   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qf_fields));
944   CeedInt             num_emode_out = 0, num_B_out_mats_to_load = 0, size_B_out = 0;
945   CeedEvalMode       *eval_mode_out = NULL;
946   CeedBasis           basis_out     = NULL;
947   CeedElemRestriction rstr_out      = NULL;
948   for (CeedInt i = 0; i < num_output_fields; i++) {
949     CeedVector vec;
950     CeedCallBackend(CeedOperatorFieldGetVector(output_fields[i], &vec));
951     if (vec == CEED_VECTOR_ACTIVE) {
952       CeedCallBackend(CeedOperatorFieldGetBasis(output_fields[i], &basis_out));
953       CeedCallBackend(CeedOperatorFieldGetElemRestriction(output_fields[i], &rstr_out));
954       if (rstr_out && rstr_out != rstr_in) {
955         // LCOV_EXCL_START
956         return CeedError(ceed, CEED_ERROR_BACKEND, "Backend does not implement multi-field non-composite operator assembly");
957         // LCOV_EXCL_STOP
958       }
959       CeedEvalMode eval_mode;
960       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode));
961       if (eval_mode != CEED_EVAL_NONE) {
962         CeedCallBackend(CeedRealloc(num_B_out_mats_to_load + 1, &eval_mode_out));
963         eval_mode_out[num_B_out_mats_to_load] = eval_mode;
964         num_B_out_mats_to_load += 1;
965         if (eval_mode == CEED_EVAL_GRAD) {
966           num_emode_out += dim;
967           size_B_out += dim * esize * nqpts;
968         } else {
969           num_emode_out += 1;
970           size_B_out += esize * nqpts;
971         }
972       }
973     }
974   }
975 
976   if (num_emode_in == 0 || num_emode_out == 0) {
977     // LCOV_EXCL_START
978     return CeedError(ceed, CEED_ERROR_UNSUPPORTED, "Cannot assemble operator without inputs/outputs");
979     // LCOV_EXCL_STOP
980   }
981 
982   CeedInt nelem, ncomp;
983   CeedCallBackend(CeedElemRestrictionGetNumElements(rstr_in, &nelem));
984   CeedCallBackend(CeedElemRestrictionGetNumComponents(rstr_in, &ncomp));
985 
986   CeedCallBackend(CeedCalloc(1, &impl->asmb));
987   CeedOperatorAssemble_Cuda *asmb = impl->asmb;
988   asmb->nelem                     = nelem;
989 
990   // Compile kernels
991   int elemsPerBlock     = 1;
992   asmb->elemsPerBlock   = elemsPerBlock;
993   CeedInt    block_size = esize * esize * elemsPerBlock;
994   Ceed_Cuda *cuda_data;
995   CeedCallBackend(CeedGetData(ceed, &cuda_data));
996   char *assembly_kernel_path, *assembly_kernel_source;
997   CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-operator-assemble.h", &assembly_kernel_path));
998   CeedDebug256(ceed, 2, "----- Loading Assembly Kernel Source -----\n");
999   CeedCallBackend(CeedLoadSourceToBuffer(ceed, assembly_kernel_path, &assembly_kernel_source));
1000   CeedDebug256(ceed, 2, "----- Loading Assembly Source Complete! -----\n");
1001   bool fallback = block_size > cuda_data->device_prop.maxThreadsPerBlock;
1002   if (fallback) {
1003     // Use fallback kernel with 1D threadblock
1004     block_size         = esize * elemsPerBlock;
1005     asmb->block_size_x = esize;
1006     asmb->block_size_y = 1;
1007   } else {  // Use kernel with 2D threadblock
1008     asmb->block_size_x = esize;
1009     asmb->block_size_y = esize;
1010   }
1011   CeedCallCuda(ceed, CeedCompileCuda(ceed, assembly_kernel_source, &asmb->module, 7, "NELEM", nelem, "NUMEMODEIN", num_emode_in, "NUMEMODEOUT",
1012                                      num_emode_out, "NQPTS", nqpts, "NNODES", esize, "BLOCK_SIZE", block_size, "NCOMP", ncomp));
1013   CeedCallCuda(ceed, CeedGetKernelCuda(ceed, asmb->module, fallback ? "linearAssembleFallback" : "linearAssemble", &asmb->linearAssemble));
1014   CeedCallBackend(CeedFree(&assembly_kernel_path));
1015   CeedCallBackend(CeedFree(&assembly_kernel_source));
1016 
1017   // Build 'full' B matrices (not 1D arrays used for tensor-product matrices)
1018   const CeedScalar *interp_in, *grad_in;
1019   CeedCallBackend(CeedBasisGetInterp(basis_in, &interp_in));
1020   CeedCallBackend(CeedBasisGetGrad(basis_in, &grad_in));
1021 
1022   // Load into B_in, in order that they will be used in eval_mode
1023   const CeedInt inBytes   = size_B_in * sizeof(CeedScalar);
1024   CeedInt       mat_start = 0;
1025   CeedCallCuda(ceed, cudaMalloc((void **)&asmb->d_B_in, inBytes));
1026   for (int i = 0; i < num_B_in_mats_to_load; i++) {
1027     CeedEvalMode eval_mode = eval_mode_in[i];
1028     if (eval_mode == CEED_EVAL_INTERP) {
1029       CeedCallCuda(ceed, cudaMemcpy(&asmb->d_B_in[mat_start], interp_in, esize * nqpts * sizeof(CeedScalar), cudaMemcpyHostToDevice));
1030       mat_start += esize * nqpts;
1031     } else if (eval_mode == CEED_EVAL_GRAD) {
1032       CeedCallCuda(ceed, cudaMemcpy(&asmb->d_B_in[mat_start], grad_in, dim * esize * nqpts * sizeof(CeedScalar), cudaMemcpyHostToDevice));
1033       mat_start += dim * esize * nqpts;
1034     }
1035   }
1036 
1037   const CeedScalar *interp_out, *grad_out;
1038   // Note that this function currently assumes 1 basis, so this should always be true
1039   // for now
1040   if (basis_out == basis_in) {
1041     interp_out = interp_in;
1042     grad_out   = grad_in;
1043   } else {
1044     CeedCallBackend(CeedBasisGetInterp(basis_out, &interp_out));
1045     CeedCallBackend(CeedBasisGetGrad(basis_out, &grad_out));
1046   }
1047 
1048   // Load into B_out, in order that they will be used in eval_mode
1049   const CeedInt outBytes = size_B_out * sizeof(CeedScalar);
1050   mat_start              = 0;
1051   CeedCallCuda(ceed, cudaMalloc((void **)&asmb->d_B_out, outBytes));
1052   for (int i = 0; i < num_B_out_mats_to_load; i++) {
1053     CeedEvalMode eval_mode = eval_mode_out[i];
1054     if (eval_mode == CEED_EVAL_INTERP) {
1055       CeedCallCuda(ceed, cudaMemcpy(&asmb->d_B_out[mat_start], interp_out, esize * nqpts * sizeof(CeedScalar), cudaMemcpyHostToDevice));
1056       mat_start += esize * nqpts;
1057     } else if (eval_mode == CEED_EVAL_GRAD) {
1058       CeedCallCuda(ceed, cudaMemcpy(&asmb->d_B_out[mat_start], grad_out, dim * esize * nqpts * sizeof(CeedScalar), cudaMemcpyHostToDevice));
1059       mat_start += dim * esize * nqpts;
1060     }
1061   }
1062   return CEED_ERROR_SUCCESS;
1063 }
1064 
1065 //------------------------------------------------------------------------------
1066 // Assemble matrix data for COO matrix of assembled operator.
1067 // The sparsity pattern is set by CeedOperatorLinearAssembleSymbolic.
1068 //
1069 // Note that this (and other assembly routines) currently assume only one active input restriction/basis per operator (could have multiple basis eval
1070 // modes).
1071 // TODO: allow multiple active input restrictions/basis objects
1072 //------------------------------------------------------------------------------
1073 static int CeedSingleOperatorAssemble_Cuda(CeedOperator op, CeedInt offset, CeedVector values) {
1074   Ceed ceed;
1075   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
1076   CeedOperator_Cuda *impl;
1077   CeedCallBackend(CeedOperatorGetData(op, &impl));
1078 
1079   // Setup
1080   if (!impl->asmb) {
1081     CeedCallBackend(CeedSingleOperatorAssembleSetup_Cuda(op));
1082     assert(impl->asmb != NULL);
1083   }
1084 
1085   // Assemble QFunction
1086   CeedVector          assembled_qf;
1087   CeedElemRestriction rstr_q;
1088   CeedCallBackend(CeedOperatorLinearAssembleQFunctionBuildOrUpdate(op, &assembled_qf, &rstr_q, CEED_REQUEST_IMMEDIATE));
1089   CeedCallBackend(CeedElemRestrictionDestroy(&rstr_q));
1090   CeedScalar *values_array;
1091   CeedCallBackend(CeedVectorGetArray(values, CEED_MEM_DEVICE, &values_array));
1092   values_array += offset;
1093   const CeedScalar *qf_array;
1094   CeedCallBackend(CeedVectorGetArrayRead(assembled_qf, CEED_MEM_DEVICE, &qf_array));
1095 
1096   // Compute B^T D B
1097   const CeedInt nelem         = impl->asmb->nelem;
1098   const CeedInt elemsPerBlock = impl->asmb->elemsPerBlock;
1099   const CeedInt grid          = nelem / elemsPerBlock + ((nelem / elemsPerBlock * elemsPerBlock < nelem) ? 1 : 0);
1100   void         *args[]        = {&impl->asmb->d_B_in, &impl->asmb->d_B_out, &qf_array, &values_array};
1101   CeedCallBackend(
1102       CeedRunKernelDimCuda(ceed, impl->asmb->linearAssemble, grid, impl->asmb->block_size_x, impl->asmb->block_size_y, elemsPerBlock, args));
1103 
1104   // Restore arrays
1105   CeedCallBackend(CeedVectorRestoreArray(values, &values_array));
1106   CeedCallBackend(CeedVectorRestoreArrayRead(assembled_qf, &qf_array));
1107 
1108   // Cleanup
1109   CeedCallBackend(CeedVectorDestroy(&assembled_qf));
1110 
1111   return CEED_ERROR_SUCCESS;
1112 }
1113 
1114 //------------------------------------------------------------------------------
1115 // Create operator
1116 //------------------------------------------------------------------------------
1117 int CeedOperatorCreate_Cuda(CeedOperator op) {
1118   Ceed ceed;
1119   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
1120   CeedOperator_Cuda *impl;
1121 
1122   CeedCallBackend(CeedCalloc(1, &impl));
1123   CeedCallBackend(CeedOperatorSetData(op, impl));
1124 
1125   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleQFunction", CeedOperatorLinearAssembleQFunction_Cuda));
1126   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleQFunctionUpdate", CeedOperatorLinearAssembleQFunctionUpdate_Cuda));
1127   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleAddDiagonal", CeedOperatorLinearAssembleAddDiagonal_Cuda));
1128   CeedCallBackend(
1129       CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleAddPointBlockDiagonal", CeedOperatorLinearAssembleAddPointBlockDiagonal_Cuda));
1130   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleSingle", CeedSingleOperatorAssemble_Cuda));
1131   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "ApplyAdd", CeedOperatorApplyAdd_Cuda));
1132   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "Destroy", CeedOperatorDestroy_Cuda));
1133   return CEED_ERROR_SUCCESS;
1134 }
1135 
1136 //------------------------------------------------------------------------------
1137