xref: /libCEED/backends/cuda-ref/ceed-cuda-ref-operator.c (revision 6beac0efa77054b88d7ece2abf46f259be72ade1)
13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
30d0321e0SJeremy L Thompson //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
50d0321e0SJeremy L Thompson //
63d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
70d0321e0SJeremy L Thompson 
80d0321e0SJeremy L Thompson #include <ceed/ceed.h>
90d0321e0SJeremy L Thompson #include <ceed/backend.h>
1007b31e0eSJeremy L Thompson #include <ceed/jit-tools.h>
110d0321e0SJeremy L Thompson #include <assert.h>
120d0321e0SJeremy L Thompson #include <cuda.h>
130d0321e0SJeremy L Thompson #include <cuda_runtime.h>
140d0321e0SJeremy L Thompson #include <stdbool.h>
150d0321e0SJeremy L Thompson #include <string.h>
160d0321e0SJeremy L Thompson #include "ceed-cuda-ref.h"
170d0321e0SJeremy L Thompson #include "../cuda/ceed-cuda-compile.h"
180d0321e0SJeremy L Thompson 
190d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
200d0321e0SJeremy L Thompson // Destroy operator
210d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
220d0321e0SJeremy L Thompson static int CeedOperatorDestroy_Cuda(CeedOperator op) {
230d0321e0SJeremy L Thompson   int ierr;
240d0321e0SJeremy L Thompson   CeedOperator_Cuda *impl;
250d0321e0SJeremy L Thompson   ierr = CeedOperatorGetData(op, &impl); CeedChkBackend(ierr);
260d0321e0SJeremy L Thompson 
270d0321e0SJeremy L Thompson   // Apply data
280d0321e0SJeremy L Thompson   for (CeedInt i = 0; i < impl->numein + impl->numeout; i++) {
290d0321e0SJeremy L Thompson     ierr = CeedVectorDestroy(&impl->evecs[i]); CeedChkBackend(ierr);
300d0321e0SJeremy L Thompson   }
310d0321e0SJeremy L Thompson   ierr = CeedFree(&impl->evecs); CeedChkBackend(ierr);
320d0321e0SJeremy L Thompson 
330d0321e0SJeremy L Thompson   for (CeedInt i = 0; i < impl->numein; i++) {
340d0321e0SJeremy L Thompson     ierr = CeedVectorDestroy(&impl->qvecsin[i]); CeedChkBackend(ierr);
350d0321e0SJeremy L Thompson   }
360d0321e0SJeremy L Thompson   ierr = CeedFree(&impl->qvecsin); CeedChkBackend(ierr);
370d0321e0SJeremy L Thompson 
380d0321e0SJeremy L Thompson   for (CeedInt i = 0; i < impl->numeout; i++) {
390d0321e0SJeremy L Thompson     ierr = CeedVectorDestroy(&impl->qvecsout[i]); CeedChkBackend(ierr);
400d0321e0SJeremy L Thompson   }
410d0321e0SJeremy L Thompson   ierr = CeedFree(&impl->qvecsout); CeedChkBackend(ierr);
420d0321e0SJeremy L Thompson 
430d0321e0SJeremy L Thompson   // QFunction assembly data
440d0321e0SJeremy L Thompson   for (CeedInt i=0; i<impl->qfnumactivein; i++) {
450d0321e0SJeremy L Thompson     ierr = CeedVectorDestroy(&impl->qfactivein[i]); CeedChkBackend(ierr);
460d0321e0SJeremy L Thompson   }
470d0321e0SJeremy L Thompson   ierr = CeedFree(&impl->qfactivein); CeedChkBackend(ierr);
480d0321e0SJeremy L Thompson 
490d0321e0SJeremy L Thompson   // Diag data
500d0321e0SJeremy L Thompson   if (impl->diag) {
510d0321e0SJeremy L Thompson     Ceed ceed;
520d0321e0SJeremy L Thompson     ierr = CeedOperatorGetCeed(op, &ceed); CeedChkBackend(ierr);
530d0321e0SJeremy L Thompson     CeedChk_Cu(ceed, cuModuleUnload(impl->diag->module));
540d0321e0SJeremy L Thompson     ierr = CeedFree(&impl->diag->h_emodein); CeedChkBackend(ierr);
550d0321e0SJeremy L Thompson     ierr = CeedFree(&impl->diag->h_emodeout); CeedChkBackend(ierr);
560d0321e0SJeremy L Thompson     ierr = cudaFree(impl->diag->d_emodein); CeedChk_Cu(ceed, ierr);
570d0321e0SJeremy L Thompson     ierr = cudaFree(impl->diag->d_emodeout); CeedChk_Cu(ceed, ierr);
580d0321e0SJeremy L Thompson     ierr = cudaFree(impl->diag->d_identity); CeedChk_Cu(ceed, ierr);
590d0321e0SJeremy L Thompson     ierr = cudaFree(impl->diag->d_interpin); CeedChk_Cu(ceed, ierr);
600d0321e0SJeremy L Thompson     ierr = cudaFree(impl->diag->d_interpout); CeedChk_Cu(ceed, ierr);
610d0321e0SJeremy L Thompson     ierr = cudaFree(impl->diag->d_gradin); CeedChk_Cu(ceed, ierr);
620d0321e0SJeremy L Thompson     ierr = cudaFree(impl->diag->d_gradout); CeedChk_Cu(ceed, ierr);
630d0321e0SJeremy L Thompson     ierr = CeedElemRestrictionDestroy(&impl->diag->pbdiagrstr);
640d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
650d0321e0SJeremy L Thompson     ierr = CeedVectorDestroy(&impl->diag->elemdiag); CeedChkBackend(ierr);
660d0321e0SJeremy L Thompson     ierr = CeedVectorDestroy(&impl->diag->pbelemdiag); CeedChkBackend(ierr);
670d0321e0SJeremy L Thompson   }
680d0321e0SJeremy L Thompson   ierr = CeedFree(&impl->diag); CeedChkBackend(ierr);
690d0321e0SJeremy L Thompson 
70cc132f9aSnbeams   if (impl->asmb) {
71cc132f9aSnbeams     Ceed ceed;
72cc132f9aSnbeams     ierr = CeedOperatorGetCeed(op, &ceed); CeedChkBackend(ierr);
73cc132f9aSnbeams     CeedChk_Cu(ceed, cuModuleUnload(impl->asmb->module));
74cc132f9aSnbeams     ierr = cudaFree(impl->asmb->d_B_in); CeedChk_Cu(ceed, ierr);
75cc132f9aSnbeams     ierr = cudaFree(impl->asmb->d_B_out); CeedChk_Cu(ceed, ierr);
76cc132f9aSnbeams   }
77cc132f9aSnbeams   ierr = CeedFree(&impl->asmb); CeedChkBackend(ierr);
78cc132f9aSnbeams 
790d0321e0SJeremy L Thompson   ierr = CeedFree(&impl); CeedChkBackend(ierr);
800d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
810d0321e0SJeremy L Thompson }
820d0321e0SJeremy L Thompson 
830d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
840d0321e0SJeremy L Thompson // Setup infields or outfields
850d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
860d0321e0SJeremy L Thompson static int CeedOperatorSetupFields_Cuda(CeedQFunction qf, CeedOperator op,
870d0321e0SJeremy L Thompson                                         bool isinput, CeedVector *evecs,
880d0321e0SJeremy L Thompson                                         CeedVector *qvecs, CeedInt starte,
890d0321e0SJeremy L Thompson                                         CeedInt numfields, CeedInt Q,
900d0321e0SJeremy L Thompson                                         CeedInt numelements) {
910d0321e0SJeremy L Thompson   CeedInt dim, ierr, size;
92d2643443SJeremy L Thompson   CeedSize q_size;
930d0321e0SJeremy L Thompson   Ceed ceed;
940d0321e0SJeremy L Thompson   ierr = CeedOperatorGetCeed(op, &ceed); CeedChkBackend(ierr);
950d0321e0SJeremy L Thompson   CeedBasis basis;
960d0321e0SJeremy L Thompson   CeedElemRestriction Erestrict;
970d0321e0SJeremy L Thompson   CeedOperatorField *opfields;
980d0321e0SJeremy L Thompson   CeedQFunctionField *qffields;
990d0321e0SJeremy L Thompson   CeedVector fieldvec;
1000d0321e0SJeremy L Thompson   bool strided;
1010d0321e0SJeremy L Thompson   bool skiprestrict;
1020d0321e0SJeremy L Thompson 
1030d0321e0SJeremy L Thompson   if (isinput) {
1040d0321e0SJeremy L Thompson     ierr = CeedOperatorGetFields(op, NULL, &opfields, NULL, NULL);
1050d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
1060d0321e0SJeremy L Thompson     ierr = CeedQFunctionGetFields(qf, NULL, &qffields, NULL, NULL);
1070d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
1080d0321e0SJeremy L Thompson   } else {
1090d0321e0SJeremy L Thompson     ierr = CeedOperatorGetFields(op, NULL, NULL, NULL, &opfields);
1100d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
1110d0321e0SJeremy L Thompson     ierr = CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qffields);
1120d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
1130d0321e0SJeremy L Thompson   }
1140d0321e0SJeremy L Thompson 
1150d0321e0SJeremy L Thompson   // Loop over fields
1160d0321e0SJeremy L Thompson   for (CeedInt i = 0; i < numfields; i++) {
1170d0321e0SJeremy L Thompson     CeedEvalMode emode;
1180d0321e0SJeremy L Thompson     ierr = CeedQFunctionFieldGetEvalMode(qffields[i], &emode); CeedChkBackend(ierr);
1190d0321e0SJeremy L Thompson 
1200d0321e0SJeremy L Thompson     strided = false;
1210d0321e0SJeremy L Thompson     skiprestrict = false;
1220d0321e0SJeremy L Thompson     if (emode != CEED_EVAL_WEIGHT) {
1230d0321e0SJeremy L Thompson       ierr = CeedOperatorFieldGetElemRestriction(opfields[i], &Erestrict);
1240d0321e0SJeremy L Thompson       CeedChkBackend(ierr);
1250d0321e0SJeremy L Thompson 
1260d0321e0SJeremy L Thompson       // Check whether this field can skip the element restriction:
1270d0321e0SJeremy L Thompson       // must be passive input, with emode NONE, and have a strided restriction with
1280d0321e0SJeremy L Thompson       // CEED_STRIDES_BACKEND.
1290d0321e0SJeremy L Thompson 
1300d0321e0SJeremy L Thompson       // First, check whether the field is input or output:
1310d0321e0SJeremy L Thompson       if (isinput) {
1320d0321e0SJeremy L Thompson         // Check for passive input:
1330d0321e0SJeremy L Thompson         ierr = CeedOperatorFieldGetVector(opfields[i], &fieldvec); CeedChkBackend(ierr);
1340d0321e0SJeremy L Thompson         if (fieldvec != CEED_VECTOR_ACTIVE) {
1350d0321e0SJeremy L Thompson           // Check emode
1360d0321e0SJeremy L Thompson           if (emode == CEED_EVAL_NONE) {
1370d0321e0SJeremy L Thompson             // Check for strided restriction
1380d0321e0SJeremy L Thompson             ierr = CeedElemRestrictionIsStrided(Erestrict, &strided);
1390d0321e0SJeremy L Thompson             CeedChkBackend(ierr);
1400d0321e0SJeremy L Thompson             if (strided) {
1410d0321e0SJeremy L Thompson               // Check if vector is already in preferred backend ordering
1420d0321e0SJeremy L Thompson               ierr = CeedElemRestrictionHasBackendStrides(Erestrict,
1430d0321e0SJeremy L Thompson                      &skiprestrict); CeedChkBackend(ierr);
1440d0321e0SJeremy L Thompson             }
1450d0321e0SJeremy L Thompson           }
1460d0321e0SJeremy L Thompson         }
1470d0321e0SJeremy L Thompson       }
1480d0321e0SJeremy L Thompson       if (skiprestrict) {
1490d0321e0SJeremy L Thompson         // We do not need an E-Vector, but will use the input field vector's data
1500d0321e0SJeremy L Thompson         // directly in the operator application.
1510d0321e0SJeremy L Thompson         evecs[i + starte] = NULL;
1520d0321e0SJeremy L Thompson       } else {
1530d0321e0SJeremy L Thompson         ierr = CeedElemRestrictionCreateVector(Erestrict, NULL,
1540d0321e0SJeremy L Thompson                                                &evecs[i + starte]);
1550d0321e0SJeremy L Thompson         CeedChkBackend(ierr);
1560d0321e0SJeremy L Thompson       }
1570d0321e0SJeremy L Thompson     }
1580d0321e0SJeremy L Thompson 
1590d0321e0SJeremy L Thompson     switch (emode) {
1600d0321e0SJeremy L Thompson     case CEED_EVAL_NONE:
1610d0321e0SJeremy L Thompson       ierr = CeedQFunctionFieldGetSize(qffields[i], &size); CeedChkBackend(ierr);
162d2643443SJeremy L Thompson       q_size = (CeedSize)numelements * Q * size;
163d2643443SJeremy L Thompson       ierr = CeedVectorCreate(ceed, q_size, &qvecs[i]); CeedChkBackend(ierr);
1640d0321e0SJeremy L Thompson       break;
1650d0321e0SJeremy L Thompson     case CEED_EVAL_INTERP:
1660d0321e0SJeremy L Thompson       ierr = CeedQFunctionFieldGetSize(qffields[i], &size); CeedChkBackend(ierr);
167d2643443SJeremy L Thompson       q_size = (CeedSize)numelements * Q * size;
168d2643443SJeremy L Thompson       ierr = CeedVectorCreate(ceed, q_size, &qvecs[i]); CeedChkBackend(ierr);
1690d0321e0SJeremy L Thompson       break;
1700d0321e0SJeremy L Thompson     case CEED_EVAL_GRAD:
1710d0321e0SJeremy L Thompson       ierr = CeedOperatorFieldGetBasis(opfields[i], &basis); CeedChkBackend(ierr);
1720d0321e0SJeremy L Thompson       ierr = CeedQFunctionFieldGetSize(qffields[i], &size); CeedChkBackend(ierr);
1730d0321e0SJeremy L Thompson       ierr = CeedBasisGetDimension(basis, &dim); CeedChkBackend(ierr);
174d2643443SJeremy L Thompson       q_size = (CeedSize)numelements * Q * size;
175d2643443SJeremy L Thompson       ierr = CeedVectorCreate(ceed, q_size, &qvecs[i]); CeedChkBackend(ierr);
1760d0321e0SJeremy L Thompson       break;
1770d0321e0SJeremy L Thompson     case CEED_EVAL_WEIGHT: // Only on input fields
1780d0321e0SJeremy L Thompson       ierr = CeedOperatorFieldGetBasis(opfields[i], &basis); CeedChkBackend(ierr);
179d2643443SJeremy L Thompson       q_size = (CeedSize)numelements * Q;
180d2643443SJeremy L Thompson       ierr = CeedVectorCreate(ceed, q_size, &qvecs[i]); CeedChkBackend(ierr);
1810d0321e0SJeremy L Thompson       ierr = CeedBasisApply(basis, numelements, CEED_NOTRANSPOSE,
1820d0321e0SJeremy L Thompson                             CEED_EVAL_WEIGHT, NULL, qvecs[i]); CeedChkBackend(ierr);
1830d0321e0SJeremy L Thompson       break;
1840d0321e0SJeremy L Thompson     case CEED_EVAL_DIV:
1850d0321e0SJeremy L Thompson       break; // TODO: Not implemented
1860d0321e0SJeremy L Thompson     case CEED_EVAL_CURL:
1870d0321e0SJeremy L Thompson       break; // TODO: Not implemented
1880d0321e0SJeremy L Thompson     }
1890d0321e0SJeremy L Thompson   }
1900d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1910d0321e0SJeremy L Thompson }
1920d0321e0SJeremy L Thompson 
1930d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
1940d0321e0SJeremy L Thompson // CeedOperator needs to connect all the named fields (be they active or passive)
1950d0321e0SJeremy L Thompson //   to the named inputs and outputs of its CeedQFunction.
1960d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
1970d0321e0SJeremy L Thompson static int CeedOperatorSetup_Cuda(CeedOperator op) {
1980d0321e0SJeremy L Thompson   int ierr;
1990d0321e0SJeremy L Thompson   bool setupdone;
2000d0321e0SJeremy L Thompson   ierr = CeedOperatorIsSetupDone(op, &setupdone); CeedChkBackend(ierr);
2010d0321e0SJeremy L Thompson   if (setupdone)
2020d0321e0SJeremy L Thompson     return CEED_ERROR_SUCCESS;
2030d0321e0SJeremy L Thompson   Ceed ceed;
2040d0321e0SJeremy L Thompson   ierr = CeedOperatorGetCeed(op, &ceed); CeedChkBackend(ierr);
2050d0321e0SJeremy L Thompson   CeedOperator_Cuda *impl;
2060d0321e0SJeremy L Thompson   ierr = CeedOperatorGetData(op, &impl); CeedChkBackend(ierr);
2070d0321e0SJeremy L Thompson   CeedQFunction qf;
2080d0321e0SJeremy L Thompson   ierr = CeedOperatorGetQFunction(op, &qf); CeedChkBackend(ierr);
2090d0321e0SJeremy L Thompson   CeedInt Q, numelements, numinputfields, numoutputfields;
2100d0321e0SJeremy L Thompson   ierr = CeedOperatorGetNumQuadraturePoints(op, &Q); CeedChkBackend(ierr);
2110d0321e0SJeremy L Thompson   ierr = CeedOperatorGetNumElements(op, &numelements); CeedChkBackend(ierr);
2120d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
2130d0321e0SJeremy L Thompson   CeedOperatorField *opinputfields, *opoutputfields;
2140d0321e0SJeremy L Thompson   ierr = CeedOperatorGetFields(op, &numinputfields, &opinputfields,
2150d0321e0SJeremy L Thompson                                &numoutputfields, &opoutputfields);
2160d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
2170d0321e0SJeremy L Thompson   CeedQFunctionField *qfinputfields, *qfoutputfields;
2180d0321e0SJeremy L Thompson   ierr = CeedQFunctionGetFields(qf, NULL, &qfinputfields, NULL, &qfoutputfields);
2190d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
2200d0321e0SJeremy L Thompson 
2210d0321e0SJeremy L Thompson   // Allocate
2220d0321e0SJeremy L Thompson   ierr = CeedCalloc(numinputfields + numoutputfields, &impl->evecs);
2230d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
2240d0321e0SJeremy L Thompson 
2250d0321e0SJeremy L Thompson   ierr = CeedCalloc(CEED_FIELD_MAX, &impl->qvecsin); CeedChkBackend(ierr);
2260d0321e0SJeremy L Thompson   ierr = CeedCalloc(CEED_FIELD_MAX, &impl->qvecsout); CeedChkBackend(ierr);
2270d0321e0SJeremy L Thompson 
2280d0321e0SJeremy L Thompson   impl->numein = numinputfields; impl->numeout = numoutputfields;
2290d0321e0SJeremy L Thompson 
2300d0321e0SJeremy L Thompson   // Set up infield and outfield evecs and qvecs
2310d0321e0SJeremy L Thompson   // Infields
2320d0321e0SJeremy L Thompson   ierr = CeedOperatorSetupFields_Cuda(qf, op, true,
2330d0321e0SJeremy L Thompson                                       impl->evecs, impl->qvecsin, 0,
2340d0321e0SJeremy L Thompson                                       numinputfields, Q, numelements);
2350d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
2360d0321e0SJeremy L Thompson 
2370d0321e0SJeremy L Thompson   // Outfields
2380d0321e0SJeremy L Thompson   ierr = CeedOperatorSetupFields_Cuda(qf, op, false,
2390d0321e0SJeremy L Thompson                                       impl->evecs, impl->qvecsout,
2400d0321e0SJeremy L Thompson                                       numinputfields, numoutputfields, Q,
2410d0321e0SJeremy L Thompson                                       numelements); CeedChkBackend(ierr);
2420d0321e0SJeremy L Thompson 
2430d0321e0SJeremy L Thompson   ierr = CeedOperatorSetSetupDone(op); CeedChkBackend(ierr);
2440d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2450d0321e0SJeremy L Thompson }
2460d0321e0SJeremy L Thompson 
2470d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
2480d0321e0SJeremy L Thompson // Setup Operator Inputs
2490d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
2500d0321e0SJeremy L Thompson static inline int CeedOperatorSetupInputs_Cuda(CeedInt numinputfields,
2510d0321e0SJeremy L Thompson     CeedQFunctionField *qfinputfields, CeedOperatorField *opinputfields,
2520d0321e0SJeremy L Thompson     CeedVector invec, const bool skipactive, CeedScalar *edata[2*CEED_FIELD_MAX],
2530d0321e0SJeremy L Thompson     CeedOperator_Cuda *impl, CeedRequest *request) {
2540d0321e0SJeremy L Thompson   CeedInt ierr;
2550d0321e0SJeremy L Thompson   CeedEvalMode emode;
2560d0321e0SJeremy L Thompson   CeedVector vec;
2570d0321e0SJeremy L Thompson   CeedElemRestriction Erestrict;
2580d0321e0SJeremy L Thompson 
2590d0321e0SJeremy L Thompson   for (CeedInt i = 0; i < numinputfields; i++) {
2600d0321e0SJeremy L Thompson     // Get input vector
2610d0321e0SJeremy L Thompson     ierr = CeedOperatorFieldGetVector(opinputfields[i], &vec); CeedChkBackend(ierr);
2620d0321e0SJeremy L Thompson     if (vec == CEED_VECTOR_ACTIVE) {
2630d0321e0SJeremy L Thompson       if (skipactive)
2640d0321e0SJeremy L Thompson         continue;
2650d0321e0SJeremy L Thompson       else
2660d0321e0SJeremy L Thompson         vec = invec;
2670d0321e0SJeremy L Thompson     }
2680d0321e0SJeremy L Thompson 
2690d0321e0SJeremy L Thompson     ierr = CeedQFunctionFieldGetEvalMode(qfinputfields[i], &emode);
2700d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
2710d0321e0SJeremy L Thompson     if (emode == CEED_EVAL_WEIGHT) { // Skip
2720d0321e0SJeremy L Thompson     } else {
2730d0321e0SJeremy L Thompson       // Get input vector
2740d0321e0SJeremy L Thompson       ierr = CeedOperatorFieldGetVector(opinputfields[i], &vec); CeedChkBackend(ierr);
2750d0321e0SJeremy L Thompson       // Get input element restriction
2760d0321e0SJeremy L Thompson       ierr = CeedOperatorFieldGetElemRestriction(opinputfields[i], &Erestrict);
2770d0321e0SJeremy L Thompson       CeedChkBackend(ierr);
2780d0321e0SJeremy L Thompson       if (vec == CEED_VECTOR_ACTIVE)
2790d0321e0SJeremy L Thompson         vec = invec;
2800d0321e0SJeremy L Thompson       // Restrict, if necessary
2810d0321e0SJeremy L Thompson       if (!impl->evecs[i]) {
2820d0321e0SJeremy L Thompson         // No restriction for this field; read data directly from vec.
2830d0321e0SJeremy L Thompson         ierr = CeedVectorGetArrayRead(vec, CEED_MEM_DEVICE,
2840d0321e0SJeremy L Thompson                                       (const CeedScalar **) &edata[i]);
2850d0321e0SJeremy L Thompson         CeedChkBackend(ierr);
2860d0321e0SJeremy L Thompson       } else {
2870d0321e0SJeremy L Thompson         ierr = CeedElemRestrictionApply(Erestrict, CEED_NOTRANSPOSE, vec,
2880d0321e0SJeremy L Thompson                                         impl->evecs[i], request); CeedChkBackend(ierr);
2890d0321e0SJeremy L Thompson         // Get evec
2900d0321e0SJeremy L Thompson         ierr = CeedVectorGetArrayRead(impl->evecs[i], CEED_MEM_DEVICE,
2910d0321e0SJeremy L Thompson                                       (const CeedScalar **) &edata[i]);
2920d0321e0SJeremy L Thompson         CeedChkBackend(ierr);
2930d0321e0SJeremy L Thompson       }
2940d0321e0SJeremy L Thompson     }
2950d0321e0SJeremy L Thompson   }
2960d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2970d0321e0SJeremy L Thompson }
2980d0321e0SJeremy L Thompson 
2990d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
3000d0321e0SJeremy L Thompson // Input Basis Action
3010d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
3020d0321e0SJeremy L Thompson static inline int CeedOperatorInputBasis_Cuda(CeedInt numelements,
3030d0321e0SJeremy L Thompson     CeedQFunctionField *qfinputfields, CeedOperatorField *opinputfields,
3040d0321e0SJeremy L Thompson     CeedInt numinputfields, const bool skipactive,
3050d0321e0SJeremy L Thompson     CeedScalar *edata[2*CEED_FIELD_MAX],CeedOperator_Cuda *impl) {
3060d0321e0SJeremy L Thompson   CeedInt ierr;
3070d0321e0SJeremy L Thompson   CeedInt elemsize, size;
3080d0321e0SJeremy L Thompson   CeedElemRestriction Erestrict;
3090d0321e0SJeremy L Thompson   CeedEvalMode emode;
3100d0321e0SJeremy L Thompson   CeedBasis basis;
3110d0321e0SJeremy L Thompson 
3120d0321e0SJeremy L Thompson   for (CeedInt i=0; i<numinputfields; i++) {
3130d0321e0SJeremy L Thompson     // Skip active input
3140d0321e0SJeremy L Thompson     if (skipactive) {
3150d0321e0SJeremy L Thompson       CeedVector vec;
3160d0321e0SJeremy L Thompson       ierr = CeedOperatorFieldGetVector(opinputfields[i], &vec); CeedChkBackend(ierr);
3170d0321e0SJeremy L Thompson       if (vec == CEED_VECTOR_ACTIVE)
3180d0321e0SJeremy L Thompson         continue;
3190d0321e0SJeremy L Thompson     }
3200d0321e0SJeremy L Thompson     // Get elemsize, emode, size
3210d0321e0SJeremy L Thompson     ierr = CeedOperatorFieldGetElemRestriction(opinputfields[i], &Erestrict);
3220d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
3230d0321e0SJeremy L Thompson     ierr = CeedElemRestrictionGetElementSize(Erestrict, &elemsize);
3240d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
3250d0321e0SJeremy L Thompson     ierr = CeedQFunctionFieldGetEvalMode(qfinputfields[i], &emode);
3260d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
3270d0321e0SJeremy L Thompson     ierr = CeedQFunctionFieldGetSize(qfinputfields[i], &size); CeedChkBackend(ierr);
3280d0321e0SJeremy L Thompson     // Basis action
3290d0321e0SJeremy L Thompson     switch (emode) {
3300d0321e0SJeremy L Thompson     case CEED_EVAL_NONE:
3310d0321e0SJeremy L Thompson       ierr = CeedVectorSetArray(impl->qvecsin[i], CEED_MEM_DEVICE,
3320d0321e0SJeremy L Thompson                                 CEED_USE_POINTER, edata[i]); CeedChkBackend(ierr);
3330d0321e0SJeremy L Thompson       break;
3340d0321e0SJeremy L Thompson     case CEED_EVAL_INTERP:
3350d0321e0SJeremy L Thompson       ierr = CeedOperatorFieldGetBasis(opinputfields[i], &basis);
3360d0321e0SJeremy L Thompson       CeedChkBackend(ierr);
3370d0321e0SJeremy L Thompson       ierr = CeedBasisApply(basis, numelements, CEED_NOTRANSPOSE,
3380d0321e0SJeremy L Thompson                             CEED_EVAL_INTERP, impl->evecs[i],
3390d0321e0SJeremy L Thompson                             impl->qvecsin[i]); CeedChkBackend(ierr);
3400d0321e0SJeremy L Thompson       break;
3410d0321e0SJeremy L Thompson     case CEED_EVAL_GRAD:
3420d0321e0SJeremy L Thompson       ierr = CeedOperatorFieldGetBasis(opinputfields[i], &basis);
3430d0321e0SJeremy L Thompson       CeedChkBackend(ierr);
3440d0321e0SJeremy L Thompson       ierr = CeedBasisApply(basis, numelements, CEED_NOTRANSPOSE,
3450d0321e0SJeremy L Thompson                             CEED_EVAL_GRAD, impl->evecs[i],
3460d0321e0SJeremy L Thompson                             impl->qvecsin[i]); CeedChkBackend(ierr);
3470d0321e0SJeremy L Thompson       break;
3480d0321e0SJeremy L Thompson     case CEED_EVAL_WEIGHT:
3490d0321e0SJeremy L Thompson       break; // No action
3500d0321e0SJeremy L Thompson     case CEED_EVAL_DIV:
3510d0321e0SJeremy L Thompson       break; // TODO: Not implemented
3520d0321e0SJeremy L Thompson     case CEED_EVAL_CURL:
3530d0321e0SJeremy L Thompson       break; // TODO: Not implemented
3540d0321e0SJeremy L Thompson     }
3550d0321e0SJeremy L Thompson   }
3560d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3570d0321e0SJeremy L Thompson }
3580d0321e0SJeremy L Thompson 
3590d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
3600d0321e0SJeremy L Thompson // Restore Input Vectors
3610d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
3620d0321e0SJeremy L Thompson static inline int CeedOperatorRestoreInputs_Cuda(CeedInt numinputfields,
3630d0321e0SJeremy L Thompson     CeedQFunctionField *qfinputfields, CeedOperatorField *opinputfields,
3640d0321e0SJeremy L Thompson     const bool skipactive, CeedScalar *edata[2*CEED_FIELD_MAX],
3650d0321e0SJeremy L Thompson     CeedOperator_Cuda *impl) {
3660d0321e0SJeremy L Thompson   CeedInt ierr;
3670d0321e0SJeremy L Thompson   CeedEvalMode emode;
3680d0321e0SJeremy L Thompson   CeedVector vec;
3690d0321e0SJeremy L Thompson 
3700d0321e0SJeremy L Thompson   for (CeedInt i = 0; i < numinputfields; i++) {
3710d0321e0SJeremy L Thompson     // Skip active input
3720d0321e0SJeremy L Thompson     if (skipactive) {
3730d0321e0SJeremy L Thompson       ierr = CeedOperatorFieldGetVector(opinputfields[i], &vec); CeedChkBackend(ierr);
3740d0321e0SJeremy L Thompson       if (vec == CEED_VECTOR_ACTIVE)
3750d0321e0SJeremy L Thompson         continue;
3760d0321e0SJeremy L Thompson     }
3770d0321e0SJeremy L Thompson     ierr = CeedQFunctionFieldGetEvalMode(qfinputfields[i], &emode);
3780d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
3790d0321e0SJeremy L Thompson     if (emode == CEED_EVAL_WEIGHT) { // Skip
3800d0321e0SJeremy L Thompson     } else {
3810d0321e0SJeremy L Thompson       if (!impl->evecs[i]) {  // This was a skiprestrict case
3820d0321e0SJeremy L Thompson         ierr = CeedOperatorFieldGetVector(opinputfields[i], &vec); CeedChkBackend(ierr);
3830d0321e0SJeremy L Thompson         ierr = CeedVectorRestoreArrayRead(vec,
3840d0321e0SJeremy L Thompson                                           (const CeedScalar **)&edata[i]);
3850d0321e0SJeremy L Thompson         CeedChkBackend(ierr);
3860d0321e0SJeremy L Thompson       } else {
3870d0321e0SJeremy L Thompson         ierr = CeedVectorRestoreArrayRead(impl->evecs[i],
3880d0321e0SJeremy L Thompson                                           (const CeedScalar **) &edata[i]);
3890d0321e0SJeremy L Thompson         CeedChkBackend(ierr);
3900d0321e0SJeremy L Thompson       }
3910d0321e0SJeremy L Thompson     }
3920d0321e0SJeremy L Thompson   }
3930d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3940d0321e0SJeremy L Thompson }
3950d0321e0SJeremy L Thompson 
3960d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
3970d0321e0SJeremy L Thompson // Apply and add to output
3980d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
3990d0321e0SJeremy L Thompson static int CeedOperatorApplyAdd_Cuda(CeedOperator op, CeedVector invec,
4000d0321e0SJeremy L Thompson                                      CeedVector outvec, CeedRequest *request) {
4010d0321e0SJeremy L Thompson   int ierr;
4020d0321e0SJeremy L Thompson   CeedOperator_Cuda *impl;
4030d0321e0SJeremy L Thompson   ierr = CeedOperatorGetData(op, &impl); CeedChkBackend(ierr);
4040d0321e0SJeremy L Thompson   CeedQFunction qf;
4050d0321e0SJeremy L Thompson   ierr = CeedOperatorGetQFunction(op, &qf); CeedChkBackend(ierr);
4060d0321e0SJeremy L Thompson   CeedInt Q, numelements, elemsize, numinputfields, numoutputfields, size;
4070d0321e0SJeremy L Thompson   ierr = CeedOperatorGetNumQuadraturePoints(op, &Q); CeedChkBackend(ierr);
4080d0321e0SJeremy L Thompson   ierr = CeedOperatorGetNumElements(op, &numelements); CeedChkBackend(ierr);
4090d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
4100d0321e0SJeremy L Thompson   CeedOperatorField *opinputfields, *opoutputfields;
4110d0321e0SJeremy L Thompson   ierr = CeedOperatorGetFields(op, &numinputfields, &opinputfields,
4120d0321e0SJeremy L Thompson                                &numoutputfields, &opoutputfields);
4130d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
4140d0321e0SJeremy L Thompson   CeedQFunctionField *qfinputfields, *qfoutputfields;
4150d0321e0SJeremy L Thompson   ierr = CeedQFunctionGetFields(qf, NULL, &qfinputfields, NULL, &qfoutputfields);
4160d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
4170d0321e0SJeremy L Thompson   CeedEvalMode emode;
4180d0321e0SJeremy L Thompson   CeedVector vec;
4190d0321e0SJeremy L Thompson   CeedBasis basis;
4200d0321e0SJeremy L Thompson   CeedElemRestriction Erestrict;
4210d0321e0SJeremy L Thompson   CeedScalar *edata[2*CEED_FIELD_MAX] = {0};
4220d0321e0SJeremy L Thompson 
4230d0321e0SJeremy L Thompson   // Setup
4240d0321e0SJeremy L Thompson   ierr = CeedOperatorSetup_Cuda(op); CeedChkBackend(ierr);
4250d0321e0SJeremy L Thompson 
4260d0321e0SJeremy L Thompson   // Input Evecs and Restriction
4270d0321e0SJeremy L Thompson   ierr = CeedOperatorSetupInputs_Cuda(numinputfields, qfinputfields,
4280d0321e0SJeremy L Thompson                                       opinputfields, invec, false, edata,
4290d0321e0SJeremy L Thompson                                       impl, request); CeedChkBackend(ierr);
4300d0321e0SJeremy L Thompson 
4310d0321e0SJeremy L Thompson   // Input basis apply if needed
4320d0321e0SJeremy L Thompson   ierr = CeedOperatorInputBasis_Cuda(numelements, qfinputfields, opinputfields,
4330d0321e0SJeremy L Thompson                                      numinputfields, false, edata, impl);
4340d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
4350d0321e0SJeremy L Thompson 
4360d0321e0SJeremy L Thompson   // Output pointers, as necessary
4370d0321e0SJeremy L Thompson   for (CeedInt i = 0; i < numoutputfields; i++) {
4380d0321e0SJeremy L Thompson     ierr = CeedQFunctionFieldGetEvalMode(qfoutputfields[i], &emode);
4390d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
4400d0321e0SJeremy L Thompson     if (emode == CEED_EVAL_NONE) {
4410d0321e0SJeremy L Thompson       // Set the output Q-Vector to use the E-Vector data directly.
4420d0321e0SJeremy L Thompson       ierr = CeedVectorGetArrayWrite(impl->evecs[i + impl->numein], CEED_MEM_DEVICE,
4430d0321e0SJeremy L Thompson                                      &edata[i + numinputfields]); CeedChkBackend(ierr);
4440d0321e0SJeremy L Thompson       ierr = CeedVectorSetArray(impl->qvecsout[i], CEED_MEM_DEVICE,
4450d0321e0SJeremy L Thompson                                 CEED_USE_POINTER, edata[i + numinputfields]);
4460d0321e0SJeremy L Thompson       CeedChkBackend(ierr);
4470d0321e0SJeremy L Thompson     }
4480d0321e0SJeremy L Thompson   }
4490d0321e0SJeremy L Thompson 
4500d0321e0SJeremy L Thompson   // Q function
4510d0321e0SJeremy L Thompson   ierr = CeedQFunctionApply(qf, numelements * Q, impl->qvecsin, impl->qvecsout);
4520d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
4530d0321e0SJeremy L Thompson 
4540d0321e0SJeremy L Thompson   // Output basis apply if needed
4550d0321e0SJeremy L Thompson   for (CeedInt i = 0; i < numoutputfields; i++) {
4560d0321e0SJeremy L Thompson     // Get elemsize, emode, size
4570d0321e0SJeremy L Thompson     ierr = CeedOperatorFieldGetElemRestriction(opoutputfields[i], &Erestrict);
4580d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
4590d0321e0SJeremy L Thompson     ierr = CeedElemRestrictionGetElementSize(Erestrict, &elemsize);
4600d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
4610d0321e0SJeremy L Thompson     ierr = CeedQFunctionFieldGetEvalMode(qfoutputfields[i], &emode);
4620d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
4630d0321e0SJeremy L Thompson     ierr = CeedQFunctionFieldGetSize(qfoutputfields[i], &size);
4640d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
4650d0321e0SJeremy L Thompson     // Basis action
4660d0321e0SJeremy L Thompson     switch (emode) {
4670d0321e0SJeremy L Thompson     case CEED_EVAL_NONE:
4680d0321e0SJeremy L Thompson       break;
4690d0321e0SJeremy L Thompson     case CEED_EVAL_INTERP:
4700d0321e0SJeremy L Thompson       ierr = CeedOperatorFieldGetBasis(opoutputfields[i], &basis);
4710d0321e0SJeremy L Thompson       CeedChkBackend(ierr);
4720d0321e0SJeremy L Thompson       ierr = CeedBasisApply(basis, numelements, CEED_TRANSPOSE,
4730d0321e0SJeremy L Thompson                             CEED_EVAL_INTERP, impl->qvecsout[i],
4740d0321e0SJeremy L Thompson                             impl->evecs[i + impl->numein]); CeedChkBackend(ierr);
4750d0321e0SJeremy L Thompson       break;
4760d0321e0SJeremy L Thompson     case CEED_EVAL_GRAD:
4770d0321e0SJeremy L Thompson       ierr = CeedOperatorFieldGetBasis(opoutputfields[i], &basis);
4780d0321e0SJeremy L Thompson       CeedChkBackend(ierr);
4790d0321e0SJeremy L Thompson       ierr = CeedBasisApply(basis, numelements, CEED_TRANSPOSE,
4800d0321e0SJeremy L Thompson                             CEED_EVAL_GRAD, impl->qvecsout[i],
4810d0321e0SJeremy L Thompson                             impl->evecs[i + impl->numein]); CeedChkBackend(ierr);
4820d0321e0SJeremy L Thompson       break;
4830d0321e0SJeremy L Thompson     // LCOV_EXCL_START
4840d0321e0SJeremy L Thompson     case CEED_EVAL_WEIGHT: {
4850d0321e0SJeremy L Thompson       Ceed ceed;
4860d0321e0SJeremy L Thompson       ierr = CeedOperatorGetCeed(op, &ceed); CeedChkBackend(ierr);
4870d0321e0SJeremy L Thompson       return CeedError(ceed, CEED_ERROR_BACKEND,
4880d0321e0SJeremy L Thompson                        "CEED_EVAL_WEIGHT cannot be an output evaluation mode");
4890d0321e0SJeremy L Thompson       break; // Should not occur
4900d0321e0SJeremy L Thompson     }
4910d0321e0SJeremy L Thompson     case CEED_EVAL_DIV:
4920d0321e0SJeremy L Thompson       break; // TODO: Not implemented
4930d0321e0SJeremy L Thompson     case CEED_EVAL_CURL:
4940d0321e0SJeremy L Thompson       break; // TODO: Not implemented
4950d0321e0SJeremy L Thompson       // LCOV_EXCL_STOP
4960d0321e0SJeremy L Thompson     }
4970d0321e0SJeremy L Thompson   }
4980d0321e0SJeremy L Thompson 
4990d0321e0SJeremy L Thompson   // Output restriction
5000d0321e0SJeremy L Thompson   for (CeedInt i = 0; i < numoutputfields; i++) {
5010d0321e0SJeremy L Thompson     // Restore evec
5020d0321e0SJeremy L Thompson     ierr = CeedQFunctionFieldGetEvalMode(qfoutputfields[i], &emode);
5030d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
5040d0321e0SJeremy L Thompson     if (emode == CEED_EVAL_NONE) {
5050d0321e0SJeremy L Thompson       ierr = CeedVectorRestoreArray(impl->evecs[i+impl->numein],
5060d0321e0SJeremy L Thompson                                     &edata[i + numinputfields]);
5070d0321e0SJeremy L Thompson       CeedChkBackend(ierr);
5080d0321e0SJeremy L Thompson     }
5090d0321e0SJeremy L Thompson     // Get output vector
5100d0321e0SJeremy L Thompson     ierr = CeedOperatorFieldGetVector(opoutputfields[i], &vec);
5110d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
5120d0321e0SJeremy L Thompson     // Restrict
5130d0321e0SJeremy L Thompson     ierr = CeedOperatorFieldGetElemRestriction(opoutputfields[i], &Erestrict);
5140d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
5150d0321e0SJeremy L Thompson     // Active
5160d0321e0SJeremy L Thompson     if (vec == CEED_VECTOR_ACTIVE)
5170d0321e0SJeremy L Thompson       vec = outvec;
5180d0321e0SJeremy L Thompson 
5190d0321e0SJeremy L Thompson     ierr = CeedElemRestrictionApply(Erestrict, CEED_TRANSPOSE,
5200d0321e0SJeremy L Thompson                                     impl->evecs[i + impl->numein], vec,
5210d0321e0SJeremy L Thompson                                     request); CeedChkBackend(ierr);
5220d0321e0SJeremy L Thompson   }
5230d0321e0SJeremy L Thompson 
5240d0321e0SJeremy L Thompson   // Restore input arrays
5250d0321e0SJeremy L Thompson   ierr = CeedOperatorRestoreInputs_Cuda(numinputfields, qfinputfields,
5260d0321e0SJeremy L Thompson                                         opinputfields, false, edata, impl);
5270d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
5280d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
5290d0321e0SJeremy L Thompson }
5300d0321e0SJeremy L Thompson 
5310d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
5320d0321e0SJeremy L Thompson // Core code for assembling linear QFunction
5330d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
5340d0321e0SJeremy L Thompson static inline int CeedOperatorLinearAssembleQFunctionCore_Cuda(CeedOperator op,
5350d0321e0SJeremy L Thompson     bool build_objects, CeedVector *assembled, CeedElemRestriction *rstr,
5360d0321e0SJeremy L Thompson     CeedRequest *request) {
5370d0321e0SJeremy L Thompson   int ierr;
5380d0321e0SJeremy L Thompson   CeedOperator_Cuda *impl;
5390d0321e0SJeremy L Thompson   ierr = CeedOperatorGetData(op, &impl); CeedChkBackend(ierr);
5400d0321e0SJeremy L Thompson   CeedQFunction qf;
5410d0321e0SJeremy L Thompson   ierr = CeedOperatorGetQFunction(op, &qf); CeedChkBackend(ierr);
5420d0321e0SJeremy L Thompson   CeedInt Q, numelements, numinputfields, numoutputfields, size;
543d2643443SJeremy L Thompson   CeedSize q_size;
5440d0321e0SJeremy L Thompson   ierr = CeedOperatorGetNumQuadraturePoints(op, &Q); CeedChkBackend(ierr);
5450d0321e0SJeremy L Thompson   ierr = CeedOperatorGetNumElements(op, &numelements); CeedChkBackend(ierr);
5460d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
5470d0321e0SJeremy L Thompson   CeedOperatorField *opinputfields, *opoutputfields;
5480d0321e0SJeremy L Thompson   ierr = CeedOperatorGetFields(op, &numinputfields, &opinputfields,
5490d0321e0SJeremy L Thompson                                &numoutputfields, &opoutputfields);
5500d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
5510d0321e0SJeremy L Thompson   CeedQFunctionField *qfinputfields, *qfoutputfields;
5520d0321e0SJeremy L Thompson   ierr = CeedQFunctionGetFields(qf, NULL, &qfinputfields, NULL, &qfoutputfields);
5530d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
5540d0321e0SJeremy L Thompson   CeedVector vec;
5550d0321e0SJeremy L Thompson   CeedInt numactivein = impl->qfnumactivein, numactiveout = impl->qfnumactiveout;
5560d0321e0SJeremy L Thompson   CeedVector *activein = impl->qfactivein;
5570d0321e0SJeremy L Thompson   CeedScalar *a, *tmp;
5580d0321e0SJeremy L Thompson   Ceed ceed, ceedparent;
5590d0321e0SJeremy L Thompson   ierr = CeedOperatorGetCeed(op, &ceed); CeedChkBackend(ierr);
5600d0321e0SJeremy L Thompson   ierr = CeedGetOperatorFallbackParentCeed(ceed, &ceedparent);
5610d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
5620d0321e0SJeremy L Thompson   ceedparent = ceedparent ? ceedparent : ceed;
5630d0321e0SJeremy L Thompson   CeedScalar *edata[2*CEED_FIELD_MAX];
5640d0321e0SJeremy L Thompson 
5650d0321e0SJeremy L Thompson   // Setup
5660d0321e0SJeremy L Thompson   ierr = CeedOperatorSetup_Cuda(op); CeedChkBackend(ierr);
5670d0321e0SJeremy L Thompson 
5680d0321e0SJeremy L Thompson   // Check for identity
5690d0321e0SJeremy L Thompson   bool identityqf;
5700d0321e0SJeremy L Thompson   ierr = CeedQFunctionIsIdentity(qf, &identityqf); CeedChkBackend(ierr);
5710d0321e0SJeremy L Thompson   if (identityqf)
5720d0321e0SJeremy L Thompson     // LCOV_EXCL_START
5730d0321e0SJeremy L Thompson     return CeedError(ceed, CEED_ERROR_BACKEND,
5740d0321e0SJeremy L Thompson                      "Assembling identity QFunctions not supported");
5750d0321e0SJeremy L Thompson   // LCOV_EXCL_STOP
5760d0321e0SJeremy L Thompson 
5770d0321e0SJeremy L Thompson   // Input Evecs and Restriction
5780d0321e0SJeremy L Thompson   ierr = CeedOperatorSetupInputs_Cuda(numinputfields, qfinputfields,
5790d0321e0SJeremy L Thompson                                       opinputfields, NULL, true, edata,
5800d0321e0SJeremy L Thompson                                       impl, request); CeedChkBackend(ierr);
5810d0321e0SJeremy L Thompson 
5820d0321e0SJeremy L Thompson   // Count number of active input fields
5830d0321e0SJeremy L Thompson   if (!numactivein) {
5840d0321e0SJeremy L Thompson     for (CeedInt i=0; i<numinputfields; i++) {
5850d0321e0SJeremy L Thompson       // Get input vector
5860d0321e0SJeremy L Thompson       ierr = CeedOperatorFieldGetVector(opinputfields[i], &vec); CeedChkBackend(ierr);
5870d0321e0SJeremy L Thompson       // Check if active input
5880d0321e0SJeremy L Thompson       if (vec == CEED_VECTOR_ACTIVE) {
5890d0321e0SJeremy L Thompson         ierr = CeedQFunctionFieldGetSize(qfinputfields[i], &size); CeedChkBackend(ierr);
5900d0321e0SJeremy L Thompson         ierr = CeedVectorSetValue(impl->qvecsin[i], 0.0); CeedChkBackend(ierr);
5910d0321e0SJeremy L Thompson         ierr = CeedVectorGetArray(impl->qvecsin[i], CEED_MEM_DEVICE, &tmp);
5920d0321e0SJeremy L Thompson         CeedChkBackend(ierr);
5930d0321e0SJeremy L Thompson         ierr = CeedRealloc(numactivein + size, &activein); CeedChkBackend(ierr);
5940d0321e0SJeremy L Thompson         for (CeedInt field = 0; field < size; field++) {
595d2643443SJeremy L Thompson           q_size = (CeedSize)Q*numelements;
596d2643443SJeremy L Thompson           ierr = CeedVectorCreate(ceed, q_size, &activein[numactivein+field]);
597d2643443SJeremy L Thompson           CeedChkBackend(ierr);
5980d0321e0SJeremy L Thompson           ierr = CeedVectorSetArray(activein[numactivein+field], CEED_MEM_DEVICE,
5990d0321e0SJeremy L Thompson                                     CEED_USE_POINTER, &tmp[field*Q*numelements]);
6000d0321e0SJeremy L Thompson           CeedChkBackend(ierr);
6010d0321e0SJeremy L Thompson         }
6020d0321e0SJeremy L Thompson         numactivein += size;
6030d0321e0SJeremy L Thompson         ierr = CeedVectorRestoreArray(impl->qvecsin[i], &tmp); CeedChkBackend(ierr);
6040d0321e0SJeremy L Thompson       }
6050d0321e0SJeremy L Thompson     }
6060d0321e0SJeremy L Thompson     impl->qfnumactivein = numactivein;
6070d0321e0SJeremy L Thompson     impl->qfactivein = activein;
6080d0321e0SJeremy L Thompson   }
6090d0321e0SJeremy L Thompson 
6100d0321e0SJeremy L Thompson   // Count number of active output fields
6110d0321e0SJeremy L Thompson   if (!numactiveout) {
6120d0321e0SJeremy L Thompson     for (CeedInt i=0; i<numoutputfields; i++) {
6130d0321e0SJeremy L Thompson       // Get output vector
6140d0321e0SJeremy L Thompson       ierr = CeedOperatorFieldGetVector(opoutputfields[i], &vec);
6150d0321e0SJeremy L Thompson       CeedChkBackend(ierr);
6160d0321e0SJeremy L Thompson       // Check if active output
6170d0321e0SJeremy L Thompson       if (vec == CEED_VECTOR_ACTIVE) {
6180d0321e0SJeremy L Thompson         ierr = CeedQFunctionFieldGetSize(qfoutputfields[i], &size);
6190d0321e0SJeremy L Thompson         CeedChkBackend(ierr);
6200d0321e0SJeremy L Thompson         numactiveout += size;
6210d0321e0SJeremy L Thompson       }
6220d0321e0SJeremy L Thompson     }
6230d0321e0SJeremy L Thompson     impl->qfnumactiveout = numactiveout;
6240d0321e0SJeremy L Thompson   }
6250d0321e0SJeremy L Thompson 
6260d0321e0SJeremy L Thompson   // Check sizes
6270d0321e0SJeremy L Thompson   if (!numactivein || !numactiveout)
6280d0321e0SJeremy L Thompson     // LCOV_EXCL_START
6290d0321e0SJeremy L Thompson     return CeedError(ceed, CEED_ERROR_BACKEND,
6300d0321e0SJeremy L Thompson                      "Cannot assemble QFunction without active inputs "
6310d0321e0SJeremy L Thompson                      "and outputs");
6320d0321e0SJeremy L Thompson   // LCOV_EXCL_STOP
6330d0321e0SJeremy L Thompson 
6340d0321e0SJeremy L Thompson   // Build objects if needed
6350d0321e0SJeremy L Thompson   if (build_objects) {
6360d0321e0SJeremy L Thompson     // Create output restriction
6370d0321e0SJeremy L Thompson     CeedInt strides[3] = {1, numelements*Q, Q}; /* *NOPAD* */
6380d0321e0SJeremy L Thompson     ierr = CeedElemRestrictionCreateStrided(ceedparent, numelements, Q,
6390d0321e0SJeremy L Thompson                                             numactivein*numactiveout,
6400d0321e0SJeremy L Thompson                                             numactivein*numactiveout*numelements*Q,
6410d0321e0SJeremy L Thompson                                             strides, rstr); CeedChkBackend(ierr);
6420d0321e0SJeremy L Thompson     // Create assembled vector
643d2643443SJeremy L Thompson     CeedSize l_size = (CeedSize)numelements*Q*numactivein*numactiveout;
644d2643443SJeremy L Thompson     ierr = CeedVectorCreate(ceedparent, l_size, assembled); CeedChkBackend(ierr);
6450d0321e0SJeremy L Thompson   }
6460d0321e0SJeremy L Thompson   ierr = CeedVectorSetValue(*assembled, 0.0); CeedChkBackend(ierr);
6470d0321e0SJeremy L Thompson   ierr = CeedVectorGetArray(*assembled, CEED_MEM_DEVICE, &a);
6480d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
6490d0321e0SJeremy L Thompson 
6500d0321e0SJeremy L Thompson   // Input basis apply
6510d0321e0SJeremy L Thompson   ierr = CeedOperatorInputBasis_Cuda(numelements, qfinputfields, opinputfields,
6520d0321e0SJeremy L Thompson                                      numinputfields, true, edata, impl);
6530d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
6540d0321e0SJeremy L Thompson 
6550d0321e0SJeremy L Thompson   // Assemble QFunction
6560d0321e0SJeremy L Thompson   for (CeedInt in=0; in<numactivein; in++) {
6570d0321e0SJeremy L Thompson     // Set Inputs
6580d0321e0SJeremy L Thompson     ierr = CeedVectorSetValue(activein[in], 1.0); CeedChkBackend(ierr);
6590d0321e0SJeremy L Thompson     if (numactivein > 1) {
6600d0321e0SJeremy L Thompson       ierr = CeedVectorSetValue(activein[(in+numactivein-1)%numactivein],
6610d0321e0SJeremy L Thompson                                 0.0); CeedChkBackend(ierr);
6620d0321e0SJeremy L Thompson     }
6630d0321e0SJeremy L Thompson     // Set Outputs
6640d0321e0SJeremy L Thompson     for (CeedInt out=0; out<numoutputfields; out++) {
6650d0321e0SJeremy L Thompson       // Get output vector
6660d0321e0SJeremy L Thompson       ierr = CeedOperatorFieldGetVector(opoutputfields[out], &vec);
6670d0321e0SJeremy L Thompson       CeedChkBackend(ierr);
6680d0321e0SJeremy L Thompson       // Check if active output
6690d0321e0SJeremy L Thompson       if (vec == CEED_VECTOR_ACTIVE) {
6700d0321e0SJeremy L Thompson         CeedVectorSetArray(impl->qvecsout[out], CEED_MEM_DEVICE,
6710d0321e0SJeremy L Thompson                            CEED_USE_POINTER, a); CeedChkBackend(ierr);
6720d0321e0SJeremy L Thompson         ierr = CeedQFunctionFieldGetSize(qfoutputfields[out], &size);
6730d0321e0SJeremy L Thompson         CeedChkBackend(ierr);
6740d0321e0SJeremy L Thompson         a += size*Q*numelements; // Advance the pointer by the size of the output
6750d0321e0SJeremy L Thompson       }
6760d0321e0SJeremy L Thompson     }
6770d0321e0SJeremy L Thompson     // Apply QFunction
6780d0321e0SJeremy L Thompson     ierr = CeedQFunctionApply(qf, Q*numelements, impl->qvecsin, impl->qvecsout);
6790d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
6800d0321e0SJeremy L Thompson   }
6810d0321e0SJeremy L Thompson 
6820d0321e0SJeremy L Thompson   // Un-set output Qvecs to prevent accidental overwrite of Assembled
6830d0321e0SJeremy L Thompson   for (CeedInt out=0; out<numoutputfields; out++) {
6840d0321e0SJeremy L Thompson     // Get output vector
6850d0321e0SJeremy L Thompson     ierr = CeedOperatorFieldGetVector(opoutputfields[out], &vec);
6860d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
6870d0321e0SJeremy L Thompson     // Check if active output
6880d0321e0SJeremy L Thompson     if (vec == CEED_VECTOR_ACTIVE) {
6890d0321e0SJeremy L Thompson       ierr = CeedVectorTakeArray(impl->qvecsout[out], CEED_MEM_DEVICE, NULL);
6900d0321e0SJeremy L Thompson       CeedChkBackend(ierr);
6910d0321e0SJeremy L Thompson     }
6920d0321e0SJeremy L Thompson   }
6930d0321e0SJeremy L Thompson 
6940d0321e0SJeremy L Thompson   // Restore input arrays
6950d0321e0SJeremy L Thompson   ierr = CeedOperatorRestoreInputs_Cuda(numinputfields, qfinputfields,
6960d0321e0SJeremy L Thompson                                         opinputfields, true, edata, impl);
6970d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
6980d0321e0SJeremy L Thompson 
6990d0321e0SJeremy L Thompson   // Restore output
7000d0321e0SJeremy L Thompson   ierr = CeedVectorRestoreArray(*assembled, &a); CeedChkBackend(ierr);
7010d0321e0SJeremy L Thompson 
7020d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
7030d0321e0SJeremy L Thompson }
7040d0321e0SJeremy L Thompson 
7050d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
7060d0321e0SJeremy L Thompson // Assemble Linear QFunction
7070d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
7080d0321e0SJeremy L Thompson static int CeedOperatorLinearAssembleQFunction_Cuda(CeedOperator op,
7090d0321e0SJeremy L Thompson     CeedVector *assembled, CeedElemRestriction *rstr, CeedRequest *request) {
7100d0321e0SJeremy L Thompson   return CeedOperatorLinearAssembleQFunctionCore_Cuda(op, true, assembled, rstr,
7110d0321e0SJeremy L Thompson          request);
7120d0321e0SJeremy L Thompson }
7130d0321e0SJeremy L Thompson 
7140d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
7150d0321e0SJeremy L Thompson // Update Assembled Linear QFunction
7160d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
7170d0321e0SJeremy L Thompson static int CeedOperatorLinearAssembleQFunctionUpdate_Cuda(CeedOperator op,
7180d0321e0SJeremy L Thompson     CeedVector assembled, CeedElemRestriction rstr, CeedRequest *request) {
7190d0321e0SJeremy L Thompson   return CeedOperatorLinearAssembleQFunctionCore_Cuda(op, false, &assembled,
7200d0321e0SJeremy L Thompson          &rstr, request);
7210d0321e0SJeremy L Thompson }
7220d0321e0SJeremy L Thompson 
7230d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
7240d0321e0SJeremy L Thompson // Create point block restriction
7250d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
7260d0321e0SJeremy L Thompson static int CreatePBRestriction(CeedElemRestriction rstr,
7270d0321e0SJeremy L Thompson                                CeedElemRestriction *pbRstr) {
7280d0321e0SJeremy L Thompson   int ierr;
7290d0321e0SJeremy L Thompson   Ceed ceed;
7300d0321e0SJeremy L Thompson   ierr = CeedElemRestrictionGetCeed(rstr, &ceed); CeedChkBackend(ierr);
7310d0321e0SJeremy L Thompson   const CeedInt *offsets;
7320d0321e0SJeremy L Thompson   ierr = CeedElemRestrictionGetOffsets(rstr, CEED_MEM_HOST, &offsets);
7330d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
7340d0321e0SJeremy L Thompson 
7350d0321e0SJeremy L Thompson   // Expand offsets
7360d0321e0SJeremy L Thompson   CeedInt nelem, ncomp, elemsize, compstride, max = 1, *pbOffsets;
7370d0321e0SJeremy L Thompson   ierr = CeedElemRestrictionGetNumElements(rstr, &nelem); CeedChkBackend(ierr);
7380d0321e0SJeremy L Thompson   ierr = CeedElemRestrictionGetNumComponents(rstr, &ncomp); CeedChkBackend(ierr);
7390d0321e0SJeremy L Thompson   ierr = CeedElemRestrictionGetElementSize(rstr, &elemsize); CeedChkBackend(ierr);
7400d0321e0SJeremy L Thompson   ierr = CeedElemRestrictionGetCompStride(rstr, &compstride);
7410d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
7420d0321e0SJeremy L Thompson   CeedInt shift = ncomp;
7430d0321e0SJeremy L Thompson   if (compstride != 1)
7440d0321e0SJeremy L Thompson     shift *= ncomp;
7450d0321e0SJeremy L Thompson   ierr = CeedCalloc(nelem*elemsize, &pbOffsets); CeedChkBackend(ierr);
7460d0321e0SJeremy L Thompson   for (CeedInt i = 0; i < nelem*elemsize; i++) {
7470d0321e0SJeremy L Thompson     pbOffsets[i] = offsets[i]*shift;
7480d0321e0SJeremy L Thompson     if (pbOffsets[i] > max)
7490d0321e0SJeremy L Thompson       max = pbOffsets[i];
7500d0321e0SJeremy L Thompson   }
7510d0321e0SJeremy L Thompson 
7520d0321e0SJeremy L Thompson   // Create new restriction
7530d0321e0SJeremy L Thompson   ierr = CeedElemRestrictionCreate(ceed, nelem, elemsize, ncomp*ncomp, 1,
7540d0321e0SJeremy L Thompson                                    max + ncomp*ncomp, CEED_MEM_HOST,
7550d0321e0SJeremy L Thompson                                    CEED_OWN_POINTER, pbOffsets, pbRstr);
7560d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
7570d0321e0SJeremy L Thompson 
7580d0321e0SJeremy L Thompson   // Cleanup
7590d0321e0SJeremy L Thompson   ierr = CeedElemRestrictionRestoreOffsets(rstr, &offsets); CeedChkBackend(ierr);
7600d0321e0SJeremy L Thompson 
7610d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
7620d0321e0SJeremy L Thompson }
7630d0321e0SJeremy L Thompson 
7640d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
7650d0321e0SJeremy L Thompson // Assemble diagonal setup
7660d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
7670d0321e0SJeremy L Thompson static inline int CeedOperatorAssembleDiagonalSetup_Cuda(CeedOperator op,
7680d0321e0SJeremy L Thompson     const bool pointBlock) {
7690d0321e0SJeremy L Thompson   int ierr;
7700d0321e0SJeremy L Thompson   Ceed ceed;
7710d0321e0SJeremy L Thompson   ierr = CeedOperatorGetCeed(op, &ceed); CeedChkBackend(ierr);
7720d0321e0SJeremy L Thompson   CeedQFunction qf;
7730d0321e0SJeremy L Thompson   ierr = CeedOperatorGetQFunction(op, &qf); CeedChkBackend(ierr);
7740d0321e0SJeremy L Thompson   CeedInt numinputfields, numoutputfields;
7750d0321e0SJeremy L Thompson   ierr = CeedQFunctionGetNumArgs(qf, &numinputfields, &numoutputfields);
7760d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
7770d0321e0SJeremy L Thompson 
7780d0321e0SJeremy L Thompson   // Determine active input basis
7790d0321e0SJeremy L Thompson   CeedOperatorField *opfields;
7800d0321e0SJeremy L Thompson   CeedQFunctionField *qffields;
7810d0321e0SJeremy L Thompson   ierr = CeedOperatorGetFields(op, NULL, &opfields, NULL, NULL);
7820d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
7830d0321e0SJeremy L Thompson   ierr = CeedQFunctionGetFields(qf, NULL, &qffields, NULL, NULL);
7840d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
7850d0321e0SJeremy L Thompson   CeedInt numemodein = 0, ncomp = 0, dim = 1;
7860d0321e0SJeremy L Thompson   CeedEvalMode *emodein = NULL;
7870d0321e0SJeremy L Thompson   CeedBasis basisin = NULL;
7880d0321e0SJeremy L Thompson   CeedElemRestriction rstrin = NULL;
7890d0321e0SJeremy L Thompson   for (CeedInt i = 0; i < numinputfields; i++) {
7900d0321e0SJeremy L Thompson     CeedVector vec;
7910d0321e0SJeremy L Thompson     ierr = CeedOperatorFieldGetVector(opfields[i], &vec); CeedChkBackend(ierr);
7920d0321e0SJeremy L Thompson     if (vec == CEED_VECTOR_ACTIVE) {
7930d0321e0SJeremy L Thompson       CeedElemRestriction rstr;
7940d0321e0SJeremy L Thompson       ierr = CeedOperatorFieldGetBasis(opfields[i], &basisin); CeedChkBackend(ierr);
7950d0321e0SJeremy L Thompson       ierr = CeedBasisGetNumComponents(basisin, &ncomp); CeedChkBackend(ierr);
7960d0321e0SJeremy L Thompson       ierr = CeedBasisGetDimension(basisin, &dim); CeedChkBackend(ierr);
7970d0321e0SJeremy L Thompson       ierr = CeedOperatorFieldGetElemRestriction(opfields[i], &rstr);
7980d0321e0SJeremy L Thompson       CeedChkBackend(ierr);
7990d0321e0SJeremy L Thompson       if (rstrin && rstrin != rstr)
8000d0321e0SJeremy L Thompson         // LCOV_EXCL_START
8010d0321e0SJeremy L Thompson         return CeedError(ceed, CEED_ERROR_BACKEND,
8020d0321e0SJeremy L Thompson                          "Multi-field non-composite operator diagonal assembly not supported");
8030d0321e0SJeremy L Thompson       // LCOV_EXCL_STOP
8040d0321e0SJeremy L Thompson       rstrin = rstr;
8050d0321e0SJeremy L Thompson       CeedEvalMode emode;
8060d0321e0SJeremy L Thompson       ierr = CeedQFunctionFieldGetEvalMode(qffields[i], &emode);
8070d0321e0SJeremy L Thompson       CeedChkBackend(ierr);
8080d0321e0SJeremy L Thompson       switch (emode) {
8090d0321e0SJeremy L Thompson       case CEED_EVAL_NONE:
8100d0321e0SJeremy L Thompson       case CEED_EVAL_INTERP:
8110d0321e0SJeremy L Thompson         ierr = CeedRealloc(numemodein + 1, &emodein); CeedChkBackend(ierr);
8120d0321e0SJeremy L Thompson         emodein[numemodein] = emode;
8130d0321e0SJeremy L Thompson         numemodein += 1;
8140d0321e0SJeremy L Thompson         break;
8150d0321e0SJeremy L Thompson       case CEED_EVAL_GRAD:
8160d0321e0SJeremy L Thompson         ierr = CeedRealloc(numemodein + dim, &emodein); CeedChkBackend(ierr);
8170d0321e0SJeremy L Thompson         for (CeedInt d = 0; d < dim; d++)
8180d0321e0SJeremy L Thompson           emodein[numemodein+d] = emode;
8190d0321e0SJeremy L Thompson         numemodein += dim;
8200d0321e0SJeremy L Thompson         break;
8210d0321e0SJeremy L Thompson       case CEED_EVAL_WEIGHT:
8220d0321e0SJeremy L Thompson       case CEED_EVAL_DIV:
8230d0321e0SJeremy L Thompson       case CEED_EVAL_CURL:
8240d0321e0SJeremy L Thompson         break; // Caught by QF Assembly
8250d0321e0SJeremy L Thompson       }
8260d0321e0SJeremy L Thompson     }
8270d0321e0SJeremy L Thompson   }
8280d0321e0SJeremy L Thompson 
8290d0321e0SJeremy L Thompson   // Determine active output basis
8300d0321e0SJeremy L Thompson   ierr = CeedOperatorGetFields(op, NULL, NULL, NULL, &opfields);
8310d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
8320d0321e0SJeremy L Thompson   ierr = CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qffields);
8330d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
8340d0321e0SJeremy L Thompson   CeedInt numemodeout = 0;
8350d0321e0SJeremy L Thompson   CeedEvalMode *emodeout = NULL;
8360d0321e0SJeremy L Thompson   CeedBasis basisout = NULL;
8370d0321e0SJeremy L Thompson   CeedElemRestriction rstrout = NULL;
8380d0321e0SJeremy L Thompson   for (CeedInt i = 0; i < numoutputfields; i++) {
8390d0321e0SJeremy L Thompson     CeedVector vec;
8400d0321e0SJeremy L Thompson     ierr = CeedOperatorFieldGetVector(opfields[i], &vec); CeedChkBackend(ierr);
8410d0321e0SJeremy L Thompson     if (vec == CEED_VECTOR_ACTIVE) {
8420d0321e0SJeremy L Thompson       CeedElemRestriction rstr;
8430d0321e0SJeremy L Thompson       ierr = CeedOperatorFieldGetBasis(opfields[i], &basisout); CeedChkBackend(ierr);
8440d0321e0SJeremy L Thompson       ierr = CeedOperatorFieldGetElemRestriction(opfields[i], &rstr);
8450d0321e0SJeremy L Thompson       CeedChkBackend(ierr);
8460d0321e0SJeremy L Thompson       if (rstrout && rstrout != rstr)
8470d0321e0SJeremy L Thompson         // LCOV_EXCL_START
8480d0321e0SJeremy L Thompson         return CeedError(ceed, CEED_ERROR_BACKEND,
8490d0321e0SJeremy L Thompson                          "Multi-field non-composite operator diagonal assembly not supported");
8500d0321e0SJeremy L Thompson       // LCOV_EXCL_STOP
8510d0321e0SJeremy L Thompson       rstrout = rstr;
8520d0321e0SJeremy L Thompson       CeedEvalMode emode;
8530d0321e0SJeremy L Thompson       ierr = CeedQFunctionFieldGetEvalMode(qffields[i], &emode); CeedChkBackend(ierr);
8540d0321e0SJeremy L Thompson       switch (emode) {
8550d0321e0SJeremy L Thompson       case CEED_EVAL_NONE:
8560d0321e0SJeremy L Thompson       case CEED_EVAL_INTERP:
8570d0321e0SJeremy L Thompson         ierr = CeedRealloc(numemodeout + 1, &emodeout); CeedChkBackend(ierr);
8580d0321e0SJeremy L Thompson         emodeout[numemodeout] = emode;
8590d0321e0SJeremy L Thompson         numemodeout += 1;
8600d0321e0SJeremy L Thompson         break;
8610d0321e0SJeremy L Thompson       case CEED_EVAL_GRAD:
8620d0321e0SJeremy L Thompson         ierr = CeedRealloc(numemodeout + dim, &emodeout); CeedChkBackend(ierr);
8630d0321e0SJeremy L Thompson         for (CeedInt d = 0; d < dim; d++)
8640d0321e0SJeremy L Thompson           emodeout[numemodeout+d] = emode;
8650d0321e0SJeremy L Thompson         numemodeout += dim;
8660d0321e0SJeremy L Thompson         break;
8670d0321e0SJeremy L Thompson       case CEED_EVAL_WEIGHT:
8680d0321e0SJeremy L Thompson       case CEED_EVAL_DIV:
8690d0321e0SJeremy L Thompson       case CEED_EVAL_CURL:
8700d0321e0SJeremy L Thompson         break; // Caught by QF Assembly
8710d0321e0SJeremy L Thompson       }
8720d0321e0SJeremy L Thompson     }
8730d0321e0SJeremy L Thompson   }
8740d0321e0SJeremy L Thompson 
8750d0321e0SJeremy L Thompson   // Operator data struct
8760d0321e0SJeremy L Thompson   CeedOperator_Cuda *impl;
8770d0321e0SJeremy L Thompson   ierr = CeedOperatorGetData(op, &impl); CeedChkBackend(ierr);
8780d0321e0SJeremy L Thompson   ierr = CeedCalloc(1, &impl->diag); CeedChkBackend(ierr);
8790d0321e0SJeremy L Thompson   CeedOperatorDiag_Cuda *diag = impl->diag;
8800d0321e0SJeremy L Thompson   diag->basisin = basisin;
8810d0321e0SJeremy L Thompson   diag->basisout = basisout;
8820d0321e0SJeremy L Thompson   diag->h_emodein = emodein;
8830d0321e0SJeremy L Thompson   diag->h_emodeout = emodeout;
8840d0321e0SJeremy L Thompson   diag->numemodein = numemodein;
8850d0321e0SJeremy L Thompson   diag->numemodeout = numemodeout;
8860d0321e0SJeremy L Thompson 
8870d0321e0SJeremy L Thompson   // Assemble kernel
88807b31e0eSJeremy L Thompson   char *diagonal_kernel_path, *diagonal_kernel_source;
88907b31e0eSJeremy L Thompson   ierr = CeedGetJitAbsolutePath(ceed,
89007b31e0eSJeremy L Thompson                                 "ceed/jit-source/cuda/cuda-ref-operator-assemble-diagonal.h",
89107b31e0eSJeremy L Thompson                                 &diagonal_kernel_path); CeedChkBackend(ierr);
89207b31e0eSJeremy L Thompson   CeedDebug256(ceed, 2, "----- Loading Diagonal Assembly Kernel Source -----\n");
89307b31e0eSJeremy L Thompson   ierr = CeedLoadSourceToBuffer(ceed, diagonal_kernel_path,
89407b31e0eSJeremy L Thompson                                 &diagonal_kernel_source);
89507b31e0eSJeremy L Thompson   CeedChkBackend(ierr);
89607b31e0eSJeremy L Thompson   CeedDebug256(ceed, 2,
89707b31e0eSJeremy L Thompson                "----- Loading Diagonal Assembly Source Complete! -----\n");
8980d0321e0SJeremy L Thompson   CeedInt nnodes, nqpts;
8990d0321e0SJeremy L Thompson   ierr = CeedBasisGetNumNodes(basisin, &nnodes); CeedChkBackend(ierr);
9000d0321e0SJeremy L Thompson   ierr = CeedBasisGetNumQuadraturePoints(basisin, &nqpts); CeedChkBackend(ierr);
9010d0321e0SJeremy L Thompson   diag->nnodes = nnodes;
90207b31e0eSJeremy L Thompson   ierr = CeedCompileCuda(ceed, diagonal_kernel_source, &diag->module, 5,
9030d0321e0SJeremy L Thompson                          "NUMEMODEIN", numemodein,
9040d0321e0SJeremy L Thompson                          "NUMEMODEOUT", numemodeout,
9050d0321e0SJeremy L Thompson                          "NNODES", nnodes,
9060d0321e0SJeremy L Thompson                          "NQPTS", nqpts,
9070d0321e0SJeremy L Thompson                          "NCOMP", ncomp
9080d0321e0SJeremy L Thompson                         ); CeedChk_Cu(ceed, ierr);
9090d0321e0SJeremy L Thompson   ierr = CeedGetKernelCuda(ceed, diag->module, "linearDiagonal",
9100d0321e0SJeremy L Thompson                            &diag->linearDiagonal); CeedChk_Cu(ceed, ierr);
9110d0321e0SJeremy L Thompson   ierr = CeedGetKernelCuda(ceed, diag->module, "linearPointBlockDiagonal",
9120d0321e0SJeremy L Thompson                            &diag->linearPointBlock);
9130d0321e0SJeremy L Thompson   CeedChk_Cu(ceed, ierr);
91407b31e0eSJeremy L Thompson   ierr = CeedFree(&diagonal_kernel_path); CeedChkBackend(ierr);
91507b31e0eSJeremy L Thompson   ierr = CeedFree(&diagonal_kernel_source); CeedChkBackend(ierr);
9160d0321e0SJeremy L Thompson 
9170d0321e0SJeremy L Thompson   // Basis matrices
9180d0321e0SJeremy L Thompson   const CeedInt qBytes = nqpts * sizeof(CeedScalar);
9190d0321e0SJeremy L Thompson   const CeedInt iBytes = qBytes * nnodes;
9200d0321e0SJeremy L Thompson   const CeedInt gBytes = qBytes * nnodes * dim;
9210d0321e0SJeremy L Thompson   const CeedInt eBytes = sizeof(CeedEvalMode);
9220d0321e0SJeremy L Thompson   const CeedScalar *interpin, *interpout, *gradin, *gradout;
9230d0321e0SJeremy L Thompson 
9240d0321e0SJeremy L Thompson   // CEED_EVAL_NONE
9250d0321e0SJeremy L Thompson   CeedScalar *identity = NULL;
9260d0321e0SJeremy L Thompson   bool evalNone = false;
9270d0321e0SJeremy L Thompson   for (CeedInt i=0; i<numemodein; i++)
9280d0321e0SJeremy L Thompson     evalNone = evalNone || (emodein[i] == CEED_EVAL_NONE);
9290d0321e0SJeremy L Thompson   for (CeedInt i=0; i<numemodeout; i++)
9300d0321e0SJeremy L Thompson     evalNone = evalNone || (emodeout[i] == CEED_EVAL_NONE);
9310d0321e0SJeremy L Thompson   if (evalNone) {
9320d0321e0SJeremy L Thompson     ierr = CeedCalloc(nqpts*nnodes, &identity); CeedChkBackend(ierr);
9330d0321e0SJeremy L Thompson     for (CeedInt i=0; i<(nnodes<nqpts?nnodes:nqpts); i++)
9340d0321e0SJeremy L Thompson       identity[i*nnodes+i] = 1.0;
9350d0321e0SJeremy L Thompson     ierr = cudaMalloc((void **)&diag->d_identity, iBytes); CeedChk_Cu(ceed, ierr);
9360d0321e0SJeremy L Thompson     ierr = cudaMemcpy(diag->d_identity, identity, iBytes,
9370d0321e0SJeremy L Thompson                       cudaMemcpyHostToDevice); CeedChk_Cu(ceed, ierr);
9380d0321e0SJeremy L Thompson   }
9390d0321e0SJeremy L Thompson 
9400d0321e0SJeremy L Thompson   // CEED_EVAL_INTERP
9410d0321e0SJeremy L Thompson   ierr = CeedBasisGetInterp(basisin, &interpin); CeedChkBackend(ierr);
9420d0321e0SJeremy L Thompson   ierr = cudaMalloc((void **)&diag->d_interpin, iBytes); CeedChk_Cu(ceed, ierr);
9430d0321e0SJeremy L Thompson   ierr = cudaMemcpy(diag->d_interpin, interpin, iBytes,
9440d0321e0SJeremy L Thompson                     cudaMemcpyHostToDevice); CeedChk_Cu(ceed, ierr);
9450d0321e0SJeremy L Thompson   ierr = CeedBasisGetInterp(basisout, &interpout); CeedChkBackend(ierr);
9460d0321e0SJeremy L Thompson   ierr = cudaMalloc((void **)&diag->d_interpout, iBytes); CeedChk_Cu(ceed, ierr);
9470d0321e0SJeremy L Thompson   ierr = cudaMemcpy(diag->d_interpout, interpout, iBytes,
9480d0321e0SJeremy L Thompson                     cudaMemcpyHostToDevice); CeedChk_Cu(ceed, ierr);
9490d0321e0SJeremy L Thompson 
9500d0321e0SJeremy L Thompson   // CEED_EVAL_GRAD
9510d0321e0SJeremy L Thompson   ierr = CeedBasisGetGrad(basisin, &gradin); CeedChkBackend(ierr);
9520d0321e0SJeremy L Thompson   ierr = cudaMalloc((void **)&diag->d_gradin, gBytes); CeedChk_Cu(ceed, ierr);
9530d0321e0SJeremy L Thompson   ierr = cudaMemcpy(diag->d_gradin, gradin, gBytes,
9540d0321e0SJeremy L Thompson                     cudaMemcpyHostToDevice); CeedChk_Cu(ceed, ierr);
9550d0321e0SJeremy L Thompson   ierr = CeedBasisGetGrad(basisout, &gradout); CeedChkBackend(ierr);
9560d0321e0SJeremy L Thompson   ierr = cudaMalloc((void **)&diag->d_gradout, gBytes); CeedChk_Cu(ceed, ierr);
9570d0321e0SJeremy L Thompson   ierr = cudaMemcpy(diag->d_gradout, gradout, gBytes,
9580d0321e0SJeremy L Thompson                     cudaMemcpyHostToDevice); CeedChk_Cu(ceed, ierr);
9590d0321e0SJeremy L Thompson 
9600d0321e0SJeremy L Thompson   // Arrays of emodes
9610d0321e0SJeremy L Thompson   ierr = cudaMalloc((void **)&diag->d_emodein, numemodein * eBytes);
9620d0321e0SJeremy L Thompson   CeedChk_Cu(ceed, ierr);
9630d0321e0SJeremy L Thompson   ierr = cudaMemcpy(diag->d_emodein, emodein, numemodein * eBytes,
9640d0321e0SJeremy L Thompson                     cudaMemcpyHostToDevice); CeedChk_Cu(ceed, ierr);
9650d0321e0SJeremy L Thompson   ierr = cudaMalloc((void **)&diag->d_emodeout, numemodeout * eBytes);
9660d0321e0SJeremy L Thompson   CeedChk_Cu(ceed, ierr);
9670d0321e0SJeremy L Thompson   ierr = cudaMemcpy(diag->d_emodeout, emodeout, numemodeout * eBytes,
9680d0321e0SJeremy L Thompson                     cudaMemcpyHostToDevice); CeedChk_Cu(ceed, ierr);
9690d0321e0SJeremy L Thompson 
9700d0321e0SJeremy L Thompson   // Restriction
9710d0321e0SJeremy L Thompson   diag->diagrstr = rstrout;
9720d0321e0SJeremy L Thompson 
9730d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
9740d0321e0SJeremy L Thompson }
9750d0321e0SJeremy L Thompson 
9760d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
9770d0321e0SJeremy L Thompson // Assemble diagonal common code
9780d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
9790d0321e0SJeremy L Thompson static inline int CeedOperatorAssembleDiagonalCore_Cuda(CeedOperator op,
9800d0321e0SJeremy L Thompson     CeedVector assembled, CeedRequest *request, const bool pointBlock) {
9810d0321e0SJeremy L Thompson   int ierr;
9820d0321e0SJeremy L Thompson   Ceed ceed;
9830d0321e0SJeremy L Thompson   ierr = CeedOperatorGetCeed(op, &ceed); CeedChkBackend(ierr);
9840d0321e0SJeremy L Thompson   CeedOperator_Cuda *impl;
9850d0321e0SJeremy L Thompson   ierr = CeedOperatorGetData(op, &impl); CeedChkBackend(ierr);
9860d0321e0SJeremy L Thompson 
9870d0321e0SJeremy L Thompson   // Assemble QFunction
9880d0321e0SJeremy L Thompson   CeedVector assembledqf;
9890d0321e0SJeremy L Thompson   CeedElemRestriction rstr;
9900d0321e0SJeremy L Thompson   ierr = CeedOperatorLinearAssembleQFunctionBuildOrUpdate(op, &assembledqf,
9910d0321e0SJeremy L Thompson          &rstr, request); CeedChkBackend(ierr);
9920d0321e0SJeremy L Thompson   ierr = CeedElemRestrictionDestroy(&rstr); CeedChkBackend(ierr);
9930d0321e0SJeremy L Thompson   CeedScalar maxnorm = 0;
9940d0321e0SJeremy L Thompson   ierr = CeedVectorNorm(assembledqf, CEED_NORM_MAX, &maxnorm);
9950d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
9960d0321e0SJeremy L Thompson 
9970d0321e0SJeremy L Thompson   // Setup
9980d0321e0SJeremy L Thompson   if (!impl->diag) {
9990d0321e0SJeremy L Thompson     ierr = CeedOperatorAssembleDiagonalSetup_Cuda(op, pointBlock);
10000d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
10010d0321e0SJeremy L Thompson   }
10020d0321e0SJeremy L Thompson   CeedOperatorDiag_Cuda *diag = impl->diag;
10030d0321e0SJeremy L Thompson   assert(diag != NULL);
10040d0321e0SJeremy L Thompson 
10050d0321e0SJeremy L Thompson   // Restriction
10060d0321e0SJeremy L Thompson   if (pointBlock && !diag->pbdiagrstr) {
10070d0321e0SJeremy L Thompson     CeedElemRestriction pbdiagrstr;
10080d0321e0SJeremy L Thompson     ierr = CreatePBRestriction(diag->diagrstr, &pbdiagrstr); CeedChkBackend(ierr);
10090d0321e0SJeremy L Thompson     diag->pbdiagrstr = pbdiagrstr;
10100d0321e0SJeremy L Thompson   }
10110d0321e0SJeremy L Thompson   CeedElemRestriction diagrstr = pointBlock ? diag->pbdiagrstr : diag->diagrstr;
10120d0321e0SJeremy L Thompson 
10130d0321e0SJeremy L Thompson   // Create diagonal vector
10140d0321e0SJeremy L Thompson   CeedVector elemdiag = pointBlock ? diag->pbelemdiag : diag->elemdiag;
10150d0321e0SJeremy L Thompson   if (!elemdiag) {
10160d0321e0SJeremy L Thompson     ierr = CeedElemRestrictionCreateVector(diagrstr, NULL, &elemdiag);
10170d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
10180d0321e0SJeremy L Thompson     if (pointBlock)
10190d0321e0SJeremy L Thompson       diag->pbelemdiag = elemdiag;
10200d0321e0SJeremy L Thompson     else
10210d0321e0SJeremy L Thompson       diag->elemdiag = elemdiag;
10220d0321e0SJeremy L Thompson   }
10230d0321e0SJeremy L Thompson   ierr = CeedVectorSetValue(elemdiag, 0.0); CeedChkBackend(ierr);
10240d0321e0SJeremy L Thompson 
10250d0321e0SJeremy L Thompson   // Assemble element operator diagonals
10260d0321e0SJeremy L Thompson   CeedScalar *elemdiagarray;
10270d0321e0SJeremy L Thompson   const CeedScalar *assembledqfarray;
10280d0321e0SJeremy L Thompson   ierr = CeedVectorGetArray(elemdiag, CEED_MEM_DEVICE, &elemdiagarray);
10290d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
10300d0321e0SJeremy L Thompson   ierr = CeedVectorGetArrayRead(assembledqf, CEED_MEM_DEVICE, &assembledqfarray);
10310d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
10320d0321e0SJeremy L Thompson   CeedInt nelem;
10330d0321e0SJeremy L Thompson   ierr = CeedElemRestrictionGetNumElements(diagrstr, &nelem);
10340d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
10350d0321e0SJeremy L Thompson 
10360d0321e0SJeremy L Thompson   // Compute the diagonal of B^T D B
10370d0321e0SJeremy L Thompson   int elemsPerBlock = 1;
10380d0321e0SJeremy L Thompson   int grid = nelem/elemsPerBlock+((nelem/elemsPerBlock*elemsPerBlock<nelem)?1:0);
10390d0321e0SJeremy L Thompson   void *args[] = {(void *) &nelem, (void *) &maxnorm, &diag->d_identity,
10400d0321e0SJeremy L Thompson                   &diag->d_interpin, &diag->d_gradin, &diag->d_interpout,
10410d0321e0SJeremy L Thompson                   &diag->d_gradout, &diag->d_emodein, &diag->d_emodeout,
10420d0321e0SJeremy L Thompson                   &assembledqfarray, &elemdiagarray
10430d0321e0SJeremy L Thompson                  };
10440d0321e0SJeremy L Thompson   if (pointBlock) {
10450d0321e0SJeremy L Thompson     ierr = CeedRunKernelDimCuda(ceed, diag->linearPointBlock, grid,
10460d0321e0SJeremy L Thompson                                 diag->nnodes, 1, elemsPerBlock, args);
10470d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
10480d0321e0SJeremy L Thompson   } else {
10490d0321e0SJeremy L Thompson     ierr = CeedRunKernelDimCuda(ceed, diag->linearDiagonal, grid,
10500d0321e0SJeremy L Thompson                                 diag->nnodes, 1, elemsPerBlock, args);
10510d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
10520d0321e0SJeremy L Thompson   }
10530d0321e0SJeremy L Thompson 
10540d0321e0SJeremy L Thompson   // Restore arrays
10550d0321e0SJeremy L Thompson   ierr = CeedVectorRestoreArray(elemdiag, &elemdiagarray); CeedChkBackend(ierr);
10560d0321e0SJeremy L Thompson   ierr = CeedVectorRestoreArrayRead(assembledqf, &assembledqfarray);
10570d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
10580d0321e0SJeremy L Thompson 
10590d0321e0SJeremy L Thompson   // Assemble local operator diagonal
10600d0321e0SJeremy L Thompson   ierr = CeedElemRestrictionApply(diagrstr, CEED_TRANSPOSE, elemdiag,
10610d0321e0SJeremy L Thompson                                   assembled, request); CeedChkBackend(ierr);
10620d0321e0SJeremy L Thompson 
10630d0321e0SJeremy L Thompson   // Cleanup
10640d0321e0SJeremy L Thompson   ierr = CeedVectorDestroy(&assembledqf); CeedChkBackend(ierr);
10650d0321e0SJeremy L Thompson 
10660d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
10670d0321e0SJeremy L Thompson }
10680d0321e0SJeremy L Thompson 
10690d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
10700d0321e0SJeremy L Thompson // Assemble composite diagonal common code
10710d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
10720d0321e0SJeremy L Thompson static inline int CeedOperatorLinearAssembleAddDiagonalCompositeCore_Cuda(
10730d0321e0SJeremy L Thompson   CeedOperator op, CeedVector assembled, CeedRequest *request,
10740d0321e0SJeremy L Thompson   const bool pointBlock) {
10750d0321e0SJeremy L Thompson   int ierr;
10760d0321e0SJeremy L Thompson   CeedInt numSub;
10770d0321e0SJeremy L Thompson   CeedOperator *subOperators;
10780d0321e0SJeremy L Thompson   ierr = CeedOperatorGetNumSub(op, &numSub); CeedChkBackend(ierr);
10790d0321e0SJeremy L Thompson   ierr = CeedOperatorGetSubList(op, &subOperators); CeedChkBackend(ierr);
10800d0321e0SJeremy L Thompson   for (CeedInt i = 0; i < numSub; i++) {
10810d0321e0SJeremy L Thompson     ierr = CeedOperatorAssembleDiagonalCore_Cuda(subOperators[i], assembled,
10820d0321e0SJeremy L Thompson            request, pointBlock); CeedChkBackend(ierr);
10830d0321e0SJeremy L Thompson   }
10840d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
10850d0321e0SJeremy L Thompson }
10860d0321e0SJeremy L Thompson 
10870d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
10880d0321e0SJeremy L Thompson // Assemble Linear Diagonal
10890d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
10900d0321e0SJeremy L Thompson static int CeedOperatorLinearAssembleAddDiagonal_Cuda(CeedOperator op,
10910d0321e0SJeremy L Thompson     CeedVector assembled, CeedRequest *request) {
10920d0321e0SJeremy L Thompson   int ierr;
10930d0321e0SJeremy L Thompson   bool isComposite;
10940d0321e0SJeremy L Thompson   ierr = CeedOperatorIsComposite(op, &isComposite); CeedChkBackend(ierr);
10950d0321e0SJeremy L Thompson   if (isComposite) {
10960d0321e0SJeremy L Thompson     return CeedOperatorLinearAssembleAddDiagonalCompositeCore_Cuda(op, assembled,
10970d0321e0SJeremy L Thompson            request, false);
10980d0321e0SJeremy L Thompson   } else {
10990d0321e0SJeremy L Thompson     return CeedOperatorAssembleDiagonalCore_Cuda(op, assembled, request, false);
11000d0321e0SJeremy L Thompson   }
11010d0321e0SJeremy L Thompson }
11020d0321e0SJeremy L Thompson 
11030d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
11040d0321e0SJeremy L Thompson // Assemble Linear Point Block Diagonal
11050d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
11060d0321e0SJeremy L Thompson static int CeedOperatorLinearAssembleAddPointBlockDiagonal_Cuda(CeedOperator op,
11070d0321e0SJeremy L Thompson     CeedVector assembled, CeedRequest *request) {
11080d0321e0SJeremy L Thompson   int ierr;
11090d0321e0SJeremy L Thompson   bool isComposite;
11100d0321e0SJeremy L Thompson   ierr = CeedOperatorIsComposite(op, &isComposite); CeedChkBackend(ierr);
11110d0321e0SJeremy L Thompson   if (isComposite) {
11120d0321e0SJeremy L Thompson     return CeedOperatorLinearAssembleAddDiagonalCompositeCore_Cuda(op, assembled,
11130d0321e0SJeremy L Thompson            request, true);
11140d0321e0SJeremy L Thompson   } else {
11150d0321e0SJeremy L Thompson     return CeedOperatorAssembleDiagonalCore_Cuda(op, assembled, request, true);
11160d0321e0SJeremy L Thompson   }
11170d0321e0SJeremy L Thompson }
11180d0321e0SJeremy L Thompson 
11190d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
1120cc132f9aSnbeams // Single operator assembly setup
1121cc132f9aSnbeams //------------------------------------------------------------------------------
1122cc132f9aSnbeams static int CeedSingleOperatorAssembleSetup_Cuda(CeedOperator op) {
1123cc132f9aSnbeams   int ierr;
1124cc132f9aSnbeams   Ceed ceed;
1125cc132f9aSnbeams   ierr = CeedOperatorGetCeed(op, &ceed); CeedChkBackend(ierr);
1126cc132f9aSnbeams   CeedOperator_Cuda *impl;
1127cc132f9aSnbeams   ierr = CeedOperatorGetData(op, &impl); CeedChkBackend(ierr);
1128cc132f9aSnbeams 
1129cc132f9aSnbeams   // Get intput and output fields
1130cc132f9aSnbeams   CeedInt num_input_fields, num_output_fields;
1131cc132f9aSnbeams   CeedOperatorField *input_fields;
1132cc132f9aSnbeams   CeedOperatorField *output_fields;
1133cc132f9aSnbeams   ierr = CeedOperatorGetFields(op, &num_input_fields, &input_fields,
1134b216396cSJeremy L Thompson                                &num_output_fields, &output_fields); CeedChkBackend(ierr);
1135cc132f9aSnbeams 
1136cc132f9aSnbeams   // Determine active input basis eval mode
1137cc132f9aSnbeams   CeedQFunction qf;
1138b216396cSJeremy L Thompson   ierr = CeedOperatorGetQFunction(op, &qf); CeedChkBackend(ierr);
1139cc132f9aSnbeams   CeedQFunctionField *qf_fields;
1140b216396cSJeremy L Thompson   ierr = CeedQFunctionGetFields(qf, NULL, &qf_fields, NULL, NULL);
1141b216396cSJeremy L Thompson   CeedChkBackend(ierr);
1142cc132f9aSnbeams   // Note that the kernel will treat each dimension of a gradient action separately;
1143cc132f9aSnbeams   // i.e., when an active input has a CEED_EVAL_GRAD mode, num_emode_in will increment
1144cc132f9aSnbeams   // by dim.  However, for the purposes of loading the B matrices, it will be treated
1145cc132f9aSnbeams   // as one mode, and we will load/copy the entire gradient matrix at once, so
1146cc132f9aSnbeams   // num_B_in_mats_to_load will be incremented by 1.
1147cc132f9aSnbeams   CeedInt num_emode_in = 0, dim = 1, num_B_in_mats_to_load = 0, size_B_in = 0;
1148cc132f9aSnbeams   CeedEvalMode *eval_mode_in = NULL; //will be of size num_B_in_mats_load
1149cc132f9aSnbeams   CeedBasis basis_in = NULL;
1150b216396cSJeremy L Thompson   CeedInt nqpts = 0, esize = 0;
1151cc132f9aSnbeams   CeedElemRestriction rstr_in = NULL;
1152cc132f9aSnbeams   for (CeedInt i=0; i<num_input_fields; i++) {
1153cc132f9aSnbeams     CeedVector vec;
1154b216396cSJeremy L Thompson     ierr = CeedOperatorFieldGetVector(input_fields[i], &vec); CeedChkBackend(ierr);
1155cc132f9aSnbeams     if (vec == CEED_VECTOR_ACTIVE) {
1156cc132f9aSnbeams       ierr = CeedOperatorFieldGetBasis(input_fields[i], &basis_in);
1157b216396cSJeremy L Thompson       CeedChkBackend(ierr);
1158b216396cSJeremy L Thompson       ierr = CeedBasisGetDimension(basis_in, &dim); CeedChkBackend(ierr);
1159b216396cSJeremy L Thompson       ierr = CeedBasisGetNumQuadraturePoints(basis_in, &nqpts); CeedChkBackend(ierr);
1160cc132f9aSnbeams       ierr = CeedOperatorFieldGetElemRestriction(input_fields[i], &rstr_in);
1161b216396cSJeremy L Thompson       CeedChkBackend(ierr);
1162b216396cSJeremy L Thompson       ierr = CeedElemRestrictionGetElementSize(rstr_in, &esize); CeedChkBackend(ierr);
1163cc132f9aSnbeams       CeedEvalMode eval_mode;
1164cc132f9aSnbeams       ierr = CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode);
1165b216396cSJeremy L Thompson       CeedChkBackend(ierr);
1166cc132f9aSnbeams       if (eval_mode != CEED_EVAL_NONE) {
1167b216396cSJeremy L Thompson         ierr = CeedRealloc(num_B_in_mats_to_load + 1, &eval_mode_in);
1168b216396cSJeremy L Thompson         CeedChkBackend(ierr);
1169cc132f9aSnbeams         eval_mode_in[num_B_in_mats_to_load] = eval_mode;
1170cc132f9aSnbeams         num_B_in_mats_to_load += 1;
1171cc132f9aSnbeams         if (eval_mode == CEED_EVAL_GRAD) {
1172cc132f9aSnbeams           num_emode_in += dim;
1173cc132f9aSnbeams           size_B_in += dim * esize * nqpts;
1174cc132f9aSnbeams         } else {
1175cc132f9aSnbeams           num_emode_in +=1;
1176cc132f9aSnbeams           size_B_in += esize * nqpts;
1177cc132f9aSnbeams         }
1178cc132f9aSnbeams       }
1179cc132f9aSnbeams     }
1180cc132f9aSnbeams   }
1181cc132f9aSnbeams 
1182cc132f9aSnbeams   // Determine active output basis; basis_out and rstr_out only used if same as input, TODO
1183b216396cSJeremy L Thompson   ierr = CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qf_fields);
1184b216396cSJeremy L Thompson   CeedChkBackend(ierr);
1185cc132f9aSnbeams   CeedInt num_emode_out = 0, num_B_out_mats_to_load = 0, size_B_out = 0;
1186cc132f9aSnbeams   CeedEvalMode *eval_mode_out = NULL;
1187cc132f9aSnbeams   CeedBasis basis_out = NULL;
1188cc132f9aSnbeams   CeedElemRestriction rstr_out = NULL;
1189cc132f9aSnbeams   for (CeedInt i=0; i<num_output_fields; i++) {
1190cc132f9aSnbeams     CeedVector vec;
1191b216396cSJeremy L Thompson     ierr = CeedOperatorFieldGetVector(output_fields[i], &vec); CeedChkBackend(ierr);
1192cc132f9aSnbeams     if (vec == CEED_VECTOR_ACTIVE) {
1193cc132f9aSnbeams       ierr = CeedOperatorFieldGetBasis(output_fields[i], &basis_out);
1194b216396cSJeremy L Thompson       CeedChkBackend(ierr);
1195cc132f9aSnbeams       ierr = CeedOperatorFieldGetElemRestriction(output_fields[i], &rstr_out);
1196b216396cSJeremy L Thompson       CeedChkBackend(ierr);
1197cc132f9aSnbeams       if (rstr_out && rstr_out != rstr_in)
1198cc132f9aSnbeams         // LCOV_EXCL_START
1199cc132f9aSnbeams         return CeedError(ceed, CEED_ERROR_BACKEND,
1200cc132f9aSnbeams                          "Multi-field non-composite operator assembly not supported");
1201cc132f9aSnbeams       // LCOV_EXCL_STOP
1202cc132f9aSnbeams       CeedEvalMode eval_mode;
1203cc132f9aSnbeams       ierr = CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode);
1204b216396cSJeremy L Thompson       CeedChkBackend(ierr);
1205cc132f9aSnbeams       if (eval_mode != CEED_EVAL_NONE) {
1206b216396cSJeremy L Thompson         ierr = CeedRealloc(num_B_out_mats_to_load + 1, &eval_mode_out);
1207b216396cSJeremy L Thompson         CeedChkBackend(ierr);
1208cc132f9aSnbeams         eval_mode_out[num_B_out_mats_to_load] = eval_mode;
1209cc132f9aSnbeams         num_B_out_mats_to_load += 1;
1210cc132f9aSnbeams         if (eval_mode == CEED_EVAL_GRAD) {
1211cc132f9aSnbeams           num_emode_out += dim;
1212cc132f9aSnbeams           size_B_out += dim * esize * nqpts;
1213cc132f9aSnbeams         } else {
1214cc132f9aSnbeams           num_emode_out +=1;
1215cc132f9aSnbeams           size_B_out += esize * nqpts;
1216cc132f9aSnbeams         }
1217cc132f9aSnbeams       }
1218cc132f9aSnbeams     }
1219cc132f9aSnbeams   }
1220cc132f9aSnbeams 
1221cc132f9aSnbeams   if (num_emode_in == 0 || num_emode_out == 0)
1222cc132f9aSnbeams     // LCOV_EXCL_START
1223cc132f9aSnbeams     return CeedError(ceed, CEED_ERROR_UNSUPPORTED,
1224cc132f9aSnbeams                      "Cannot assemble operator without inputs/outputs");
1225cc132f9aSnbeams   // LCOV_EXCL_STOP
1226cc132f9aSnbeams 
1227cc132f9aSnbeams   CeedInt nelem, ncomp;
1228b216396cSJeremy L Thompson   ierr = CeedElemRestrictionGetNumElements(rstr_in, &nelem); CeedChkBackend(ierr);
1229b216396cSJeremy L Thompson   ierr = CeedElemRestrictionGetNumComponents(rstr_in, &ncomp);
1230b216396cSJeremy L Thompson   CeedChkBackend(ierr);
1231cc132f9aSnbeams 
1232cc132f9aSnbeams   ierr = CeedCalloc(1, &impl->asmb); CeedChkBackend(ierr);
1233cc132f9aSnbeams   CeedOperatorAssemble_Cuda *asmb = impl->asmb;
1234cc132f9aSnbeams   asmb->nelem = nelem;
1235cc132f9aSnbeams 
1236cc132f9aSnbeams   // Compile kernels
1237cc132f9aSnbeams   int elemsPerBlock = 1;
1238cc132f9aSnbeams   asmb->elemsPerBlock = elemsPerBlock;
123959ad764aSnbeams   CeedInt block_size = esize * esize * elemsPerBlock;
1240cc132f9aSnbeams   Ceed_Cuda *cuda_data;
1241cc132f9aSnbeams   ierr = CeedGetData(ceed, &cuda_data); CeedChkBackend(ierr);
124207b31e0eSJeremy L Thompson   char *assembly_kernel_path, *assembly_kernel_source;
124307b31e0eSJeremy L Thompson   ierr = CeedGetJitAbsolutePath(ceed,
124407b31e0eSJeremy L Thompson                                 "ceed/jit-source/cuda/cuda-ref-operator-assemble.h",
124507b31e0eSJeremy L Thompson                                 &assembly_kernel_path); CeedChkBackend(ierr);
124607b31e0eSJeremy L Thompson   CeedDebug256(ceed, 2, "----- Loading Assembly Kernel Source -----\n");
124707b31e0eSJeremy L Thompson   ierr = CeedLoadSourceToBuffer(ceed, assembly_kernel_path,
124807b31e0eSJeremy L Thompson                                 &assembly_kernel_source);
124907b31e0eSJeremy L Thompson   CeedChkBackend(ierr);
125007b31e0eSJeremy L Thompson   CeedDebug256(ceed, 2, "----- Loading Assembly Source Complete! -----\n");
125107b31e0eSJeremy L Thompson   bool fallback = block_size > cuda_data->device_prop.maxThreadsPerBlock;
125207b31e0eSJeremy L Thompson   if (fallback) {
125359ad764aSnbeams     // Use fallback kernel with 1D threadblock
125459ad764aSnbeams     block_size = esize * elemsPerBlock;
125559ad764aSnbeams     asmb->block_size_x = esize;
125659ad764aSnbeams     asmb->block_size_y = 1;
125759ad764aSnbeams   } else {  // Use kernel with 2D threadblock
125859ad764aSnbeams     asmb->block_size_x = esize;
125959ad764aSnbeams     asmb->block_size_y = esize;
126007b31e0eSJeremy L Thompson   }
126107b31e0eSJeremy L Thompson   ierr = CeedCompileCuda(ceed, assembly_kernel_source, &asmb->module, 7,
126259ad764aSnbeams                          "NELEM", nelem,
126359ad764aSnbeams                          "NUMEMODEIN", num_emode_in,
126459ad764aSnbeams                          "NUMEMODEOUT", num_emode_out,
126559ad764aSnbeams                          "NQPTS", nqpts,
126659ad764aSnbeams                          "NNODES", esize,
126759ad764aSnbeams                          "BLOCK_SIZE", block_size,
126859ad764aSnbeams                          "NCOMP", ncomp
126959ad764aSnbeams                         ); CeedChk_Cu(ceed, ierr);
127007b31e0eSJeremy L Thompson   ierr = CeedGetKernelCuda(ceed, asmb->module,
127107b31e0eSJeremy L Thompson                            fallback ? "linearAssembleFallback" : "linearAssemble",
1272cc132f9aSnbeams                            &asmb->linearAssemble); CeedChk_Cu(ceed, ierr);
127307b31e0eSJeremy L Thompson   ierr = CeedFree(&assembly_kernel_path); CeedChkBackend(ierr);
127407b31e0eSJeremy L Thompson   ierr = CeedFree(&assembly_kernel_source); CeedChkBackend(ierr);
1275cc132f9aSnbeams 
1276cc132f9aSnbeams   // Build 'full' B matrices (not 1D arrays used for tensor-product matrices)
1277cc132f9aSnbeams   const CeedScalar *interp_in, *grad_in;
1278cc132f9aSnbeams   ierr = CeedBasisGetInterp(basis_in, &interp_in); CeedChkBackend(ierr);
1279cc132f9aSnbeams   ierr = CeedBasisGetGrad(basis_in, &grad_in); CeedChkBackend(ierr);
1280cc132f9aSnbeams 
1281cc132f9aSnbeams   // Load into B_in, in order that they will be used in eval_mode
1282cc132f9aSnbeams   const CeedInt inBytes = size_B_in * sizeof(CeedScalar);
1283cc132f9aSnbeams   CeedInt mat_start = 0;
1284cc132f9aSnbeams   ierr = cudaMalloc((void **) &asmb->d_B_in, inBytes); CeedChk_Cu(ceed, ierr);
1285cc132f9aSnbeams   for (int i = 0; i < num_B_in_mats_to_load; i++) {
1286cc132f9aSnbeams     CeedEvalMode eval_mode = eval_mode_in[i];
1287cc132f9aSnbeams     if (eval_mode == CEED_EVAL_INTERP) {
1288cc132f9aSnbeams       ierr = cudaMemcpy(&asmb->d_B_in[mat_start], interp_in,
1289cc132f9aSnbeams                         esize * nqpts * sizeof(CeedScalar),
1290cc132f9aSnbeams                         cudaMemcpyHostToDevice); CeedChk_Cu(ceed, ierr);
1291cc132f9aSnbeams       mat_start += esize * nqpts;
1292cc132f9aSnbeams     } else if (eval_mode == CEED_EVAL_GRAD) {
1293*6beac0efSnbeams       ierr = cudaMemcpy(&asmb->d_B_in[mat_start], grad_in,
1294cc132f9aSnbeams                         dim * esize * nqpts * sizeof(CeedScalar),
1295cc132f9aSnbeams                         cudaMemcpyHostToDevice); CeedChk_Cu(ceed, ierr);
1296cc132f9aSnbeams       mat_start += dim * esize * nqpts;
1297cc132f9aSnbeams     }
1298cc132f9aSnbeams   }
1299cc132f9aSnbeams 
1300cc132f9aSnbeams   const CeedScalar *interp_out, *grad_out;
1301cc132f9aSnbeams   // Note that this function currently assumes 1 basis, so this should always be true
1302cc132f9aSnbeams   // for now
1303cc132f9aSnbeams   if (basis_out == basis_in) {
1304cc132f9aSnbeams     interp_out = interp_in;
1305cc132f9aSnbeams     grad_out = grad_in;
1306cc132f9aSnbeams   } else {
1307cc132f9aSnbeams     ierr = CeedBasisGetInterp(basis_out, &interp_out); CeedChkBackend(ierr);
1308cc132f9aSnbeams     ierr = CeedBasisGetGrad(basis_out, &grad_out); CeedChkBackend(ierr);
1309cc132f9aSnbeams   }
1310cc132f9aSnbeams 
1311cc132f9aSnbeams   // Load into B_out, in order that they will be used in eval_mode
1312cc132f9aSnbeams   const CeedInt outBytes = size_B_out * sizeof(CeedScalar);
1313cc132f9aSnbeams   mat_start = 0;
1314cc132f9aSnbeams   ierr = cudaMalloc((void **) &asmb->d_B_out, outBytes); CeedChk_Cu(ceed, ierr);
1315cc132f9aSnbeams   for (int i = 0; i < num_B_out_mats_to_load; i++) {
1316cc132f9aSnbeams     CeedEvalMode eval_mode = eval_mode_out[i];
1317cc132f9aSnbeams     if (eval_mode == CEED_EVAL_INTERP) {
1318cc132f9aSnbeams       ierr = cudaMemcpy(&asmb->d_B_out[mat_start], interp_out,
1319cc132f9aSnbeams                         esize * nqpts * sizeof(CeedScalar),
1320cc132f9aSnbeams                         cudaMemcpyHostToDevice); CeedChk_Cu(ceed, ierr);
1321cc132f9aSnbeams       mat_start += esize * nqpts;
1322cc132f9aSnbeams     } else if (eval_mode == CEED_EVAL_GRAD) {
1323cc132f9aSnbeams       ierr = cudaMemcpy(&asmb->d_B_out[mat_start], grad_out,
1324cc132f9aSnbeams                         dim * esize * nqpts * sizeof(CeedScalar),
1325cc132f9aSnbeams                         cudaMemcpyHostToDevice); CeedChk_Cu(ceed, ierr);
1326cc132f9aSnbeams       mat_start += dim * esize * nqpts;
1327cc132f9aSnbeams     }
1328cc132f9aSnbeams   }
1329cc132f9aSnbeams   return CEED_ERROR_SUCCESS;
1330cc132f9aSnbeams }
1331cc132f9aSnbeams 
1332cc132f9aSnbeams //------------------------------------------------------------------------------
1333cc132f9aSnbeams // Single operator assembly
1334cc132f9aSnbeams //------------------------------------------------------------------------------
1335cc132f9aSnbeams static int CeedSingleOperatorAssemble_Cuda(CeedOperator op, CeedInt offset,
1336cc132f9aSnbeams     CeedVector values) {
1337cc132f9aSnbeams 
1338cc132f9aSnbeams   int ierr;
1339cc132f9aSnbeams   Ceed ceed;
1340cc132f9aSnbeams   ierr = CeedOperatorGetCeed(op, &ceed); CeedChkBackend(ierr);
1341cc132f9aSnbeams   CeedOperator_Cuda *impl;
1342cc132f9aSnbeams   ierr = CeedOperatorGetData(op, &impl); CeedChkBackend(ierr);
1343cc132f9aSnbeams 
1344cc132f9aSnbeams   // Setup
1345cc132f9aSnbeams   if (!impl->asmb) {
1346cc132f9aSnbeams     ierr = CeedSingleOperatorAssembleSetup_Cuda(op);
1347cc132f9aSnbeams     CeedChkBackend(ierr);
1348b216396cSJeremy L Thompson     assert(impl->asmb != NULL);
1349cc132f9aSnbeams   }
1350cc132f9aSnbeams 
1351cc132f9aSnbeams   // Assemble QFunction
1352cc132f9aSnbeams   CeedVector assembled_qf;
1353cc132f9aSnbeams   CeedElemRestriction rstr_q;
1354cc132f9aSnbeams   ierr = CeedOperatorLinearAssembleQFunctionBuildOrUpdate(
1355b216396cSJeremy L Thompson            op, &assembled_qf, &rstr_q, CEED_REQUEST_IMMEDIATE); CeedChkBackend(ierr);
1356cc132f9aSnbeams   ierr = CeedElemRestrictionDestroy(&rstr_q); CeedChkBackend(ierr);
1357cc132f9aSnbeams   CeedScalar *values_array;
1358cc132f9aSnbeams   ierr = CeedVectorGetArrayWrite(values, CEED_MEM_DEVICE, &values_array);
1359cc132f9aSnbeams   CeedChkBackend(ierr);
1360cc132f9aSnbeams   values_array += offset;
1361cc132f9aSnbeams   const CeedScalar *qf_array;
1362cc132f9aSnbeams   ierr = CeedVectorGetArrayRead(assembled_qf, CEED_MEM_DEVICE, &qf_array);
1363cc132f9aSnbeams   CeedChkBackend(ierr);
1364cc132f9aSnbeams 
1365cc132f9aSnbeams   // Compute B^T D B
1366cc132f9aSnbeams   const CeedInt nelem = impl->asmb->nelem;
1367cc132f9aSnbeams   const CeedInt elemsPerBlock = impl->asmb->elemsPerBlock;
1368cc132f9aSnbeams   const CeedInt grid = nelem/elemsPerBlock+((
1369cc132f9aSnbeams                          nelem/elemsPerBlock*elemsPerBlock<nelem)?1:0);
1370cc132f9aSnbeams   void *args[] = {&impl->asmb->d_B_in, &impl->asmb->d_B_out,
1371cc132f9aSnbeams                   &qf_array, &values_array
1372cc132f9aSnbeams                  };
1373cc132f9aSnbeams   ierr = CeedRunKernelDimCuda(ceed, impl->asmb->linearAssemble, grid,
137459ad764aSnbeams                               impl->asmb->block_size_x, impl->asmb->block_size_y,
137559ad764aSnbeams                               elemsPerBlock, args);
1376cc132f9aSnbeams   CeedChkBackend(ierr);
1377cc132f9aSnbeams 
1378cc132f9aSnbeams 
1379cc132f9aSnbeams   // Restore arrays
1380cc132f9aSnbeams   ierr = CeedVectorRestoreArray(values, &values_array); CeedChkBackend(ierr);
1381cc132f9aSnbeams   ierr = CeedVectorRestoreArrayRead(assembled_qf, &qf_array);
1382cc132f9aSnbeams   CeedChkBackend(ierr);
1383cc132f9aSnbeams 
1384cc132f9aSnbeams   // Cleanup
1385cc132f9aSnbeams   ierr = CeedVectorDestroy(&assembled_qf); CeedChkBackend(ierr);
1386cc132f9aSnbeams 
1387cc132f9aSnbeams   return CEED_ERROR_SUCCESS;
1388cc132f9aSnbeams }
1389cc132f9aSnbeams 
1390cc132f9aSnbeams //------------------------------------------------------------------------------
1391cc132f9aSnbeams // Assemble matrix data for COO matrix of assembled operator.
1392cc132f9aSnbeams // The sparsity pattern is set by CeedOperatorLinearAssembleSymbolic.
1393cc132f9aSnbeams //
1394cc132f9aSnbeams // Note that this (and other assembly routines) currently assume only one
1395cc132f9aSnbeams // active input restriction/basis per operator (could have multiple basis eval
1396cc132f9aSnbeams // modes).
1397cc132f9aSnbeams // TODO: allow multiple active input restrictions/basis objects
1398cc132f9aSnbeams //------------------------------------------------------------------------------
1399cc132f9aSnbeams int CeedOperatorLinearAssemble_Cuda(CeedOperator op, CeedVector values) {
1400cc132f9aSnbeams 
1401cc132f9aSnbeams   // As done in the default implementation, loop through suboperators
1402cc132f9aSnbeams   // for composite operators, or call single operator assembly otherwise
1403cc132f9aSnbeams   bool is_composite;
1404cc132f9aSnbeams   CeedInt ierr;
1405b216396cSJeremy L Thompson   ierr = CeedOperatorIsComposite(op, &is_composite); CeedChkBackend(ierr);
1406cc132f9aSnbeams 
1407cc132f9aSnbeams   CeedElemRestriction rstr;
1408cc132f9aSnbeams   CeedInt num_elem, elem_size, num_comp;
1409cc132f9aSnbeams 
1410cc132f9aSnbeams   CeedInt offset = 0;
1411cc132f9aSnbeams   if (is_composite) {
1412cc132f9aSnbeams     CeedInt num_suboperators;
1413b216396cSJeremy L Thompson     ierr = CeedOperatorGetNumSub(op, &num_suboperators); CeedChkBackend(ierr);
1414cc132f9aSnbeams     CeedOperator *sub_operators;
1415b216396cSJeremy L Thompson     ierr = CeedOperatorGetSubList(op, &sub_operators); CeedChkBackend(ierr);
1416cc132f9aSnbeams     for (int k = 0; k < num_suboperators; ++k) {
1417cc132f9aSnbeams       ierr = CeedSingleOperatorAssemble_Cuda(sub_operators[k], offset, values);
1418b216396cSJeremy L Thompson       CeedChkBackend(ierr);
1419cc132f9aSnbeams       ierr = CeedOperatorGetActiveElemRestriction(sub_operators[k], &rstr);
1420b216396cSJeremy L Thompson       CeedChkBackend(ierr);
1421b216396cSJeremy L Thompson       ierr = CeedElemRestrictionGetNumElements(rstr, &num_elem); CeedChkBackend(ierr);
1422b216396cSJeremy L Thompson       ierr = CeedElemRestrictionGetElementSize(rstr, &elem_size);
1423b216396cSJeremy L Thompson       CeedChkBackend(ierr);
1424b216396cSJeremy L Thompson       ierr = CeedElemRestrictionGetNumComponents(rstr, &num_comp);
1425b216396cSJeremy L Thompson       CeedChkBackend(ierr);
1426cc132f9aSnbeams       offset += elem_size*num_comp * elem_size*num_comp * num_elem;
1427cc132f9aSnbeams     }
1428cc132f9aSnbeams   } else {
1429b216396cSJeremy L Thompson     ierr = CeedSingleOperatorAssemble_Cuda(op, offset, values);
1430b216396cSJeremy L Thompson     CeedChkBackend(ierr);
1431cc132f9aSnbeams   }
1432cc132f9aSnbeams 
1433cc132f9aSnbeams   return CEED_ERROR_SUCCESS;
1434cc132f9aSnbeams }
1435cc132f9aSnbeams //------------------------------------------------------------------------------
14360d0321e0SJeremy L Thompson // Create operator
14370d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
14380d0321e0SJeremy L Thompson int CeedOperatorCreate_Cuda(CeedOperator op) {
14390d0321e0SJeremy L Thompson   int ierr;
14400d0321e0SJeremy L Thompson   Ceed ceed;
14410d0321e0SJeremy L Thompson   ierr = CeedOperatorGetCeed(op, &ceed); CeedChkBackend(ierr);
14420d0321e0SJeremy L Thompson   CeedOperator_Cuda *impl;
14430d0321e0SJeremy L Thompson 
14440d0321e0SJeremy L Thompson   ierr = CeedCalloc(1, &impl); CeedChkBackend(ierr);
14450d0321e0SJeremy L Thompson   ierr = CeedOperatorSetData(op, impl); CeedChkBackend(ierr);
14460d0321e0SJeremy L Thompson 
14470d0321e0SJeremy L Thompson   ierr = CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleQFunction",
14480d0321e0SJeremy L Thompson                                 CeedOperatorLinearAssembleQFunction_Cuda);
14490d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
14500d0321e0SJeremy L Thompson   ierr = CeedSetBackendFunction(ceed, "Operator", op,
14510d0321e0SJeremy L Thompson                                 "LinearAssembleQFunctionUpdate",
14520d0321e0SJeremy L Thompson                                 CeedOperatorLinearAssembleQFunctionUpdate_Cuda);
14530d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
14540d0321e0SJeremy L Thompson   ierr = CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleAddDiagonal",
14550d0321e0SJeremy L Thompson                                 CeedOperatorLinearAssembleAddDiagonal_Cuda);
14560d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
14570d0321e0SJeremy L Thompson   ierr = CeedSetBackendFunction(ceed, "Operator", op,
14580d0321e0SJeremy L Thompson                                 "LinearAssembleAddPointBlockDiagonal",
14590d0321e0SJeremy L Thompson                                 CeedOperatorLinearAssembleAddPointBlockDiagonal_Cuda);
14600d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
1461cc132f9aSnbeams   ierr = CeedSetBackendFunction(ceed, "Operator", op,
1462cc132f9aSnbeams                                 "LinearAssemble", CeedOperatorLinearAssemble_Cuda);
1463cc132f9aSnbeams   CeedChkBackend(ierr);
14640d0321e0SJeremy L Thompson   ierr = CeedSetBackendFunction(ceed, "Operator", op, "ApplyAdd",
14650d0321e0SJeremy L Thompson                                 CeedOperatorApplyAdd_Cuda); CeedChkBackend(ierr);
14660d0321e0SJeremy L Thompson   ierr = CeedSetBackendFunction(ceed, "Operator", op, "Destroy",
14670d0321e0SJeremy L Thompson                                 CeedOperatorDestroy_Cuda); CeedChkBackend(ierr);
14680d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
14690d0321e0SJeremy L Thompson }
14700d0321e0SJeremy L Thompson 
14710d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
14720d0321e0SJeremy L Thompson // Composite Operator Create
14730d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
14740d0321e0SJeremy L Thompson int CeedCompositeOperatorCreate_Cuda(CeedOperator op) {
14750d0321e0SJeremy L Thompson   int ierr;
14760d0321e0SJeremy L Thompson   Ceed ceed;
14770d0321e0SJeremy L Thompson   ierr = CeedOperatorGetCeed(op, &ceed); CeedChkBackend(ierr);
14780d0321e0SJeremy L Thompson 
14790d0321e0SJeremy L Thompson   ierr = CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleAddDiagonal",
14800d0321e0SJeremy L Thompson                                 CeedOperatorLinearAssembleAddDiagonal_Cuda);
14810d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
14820d0321e0SJeremy L Thompson   ierr = CeedSetBackendFunction(ceed, "Operator", op,
14830d0321e0SJeremy L Thompson                                 "LinearAssembleAddPointBlockDiagonal",
14840d0321e0SJeremy L Thompson                                 CeedOperatorLinearAssembleAddPointBlockDiagonal_Cuda);
14850d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
1486cc132f9aSnbeams   ierr = CeedSetBackendFunction(ceed, "Operator", op,
1487cc132f9aSnbeams                                 "LinearAssemble", CeedOperatorLinearAssemble_Cuda);
1488cc132f9aSnbeams   CeedChkBackend(ierr);
14890d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
14900d0321e0SJeremy L Thompson }
14910d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
1492