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