xref: /libCEED/rust/libceed-sys/c-src/backends/cuda-ref/ceed-cuda-ref-operator.c (revision b216396c6ed5d64bcae95e98b720fbdd56ea5be1)
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>
100d0321e0SJeremy L Thompson #include <assert.h>
110d0321e0SJeremy L Thompson #include <cuda.h>
120d0321e0SJeremy L Thompson #include <cuda_runtime.h>
130d0321e0SJeremy L Thompson #include <stdbool.h>
140d0321e0SJeremy L Thompson #include <string.h>
150d0321e0SJeremy L Thompson #include "ceed-cuda-ref.h"
160d0321e0SJeremy L Thompson #include "../cuda/ceed-cuda-compile.h"
170d0321e0SJeremy L Thompson 
180d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
190d0321e0SJeremy L Thompson // Destroy operator
200d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
210d0321e0SJeremy L Thompson static int CeedOperatorDestroy_Cuda(CeedOperator op) {
220d0321e0SJeremy L Thompson   int ierr;
230d0321e0SJeremy L Thompson   CeedOperator_Cuda *impl;
240d0321e0SJeremy L Thompson   ierr = CeedOperatorGetData(op, &impl); CeedChkBackend(ierr);
250d0321e0SJeremy L Thompson 
260d0321e0SJeremy L Thompson   // Apply data
270d0321e0SJeremy L Thompson   for (CeedInt i = 0; i < impl->numein + impl->numeout; i++) {
280d0321e0SJeremy L Thompson     ierr = CeedVectorDestroy(&impl->evecs[i]); CeedChkBackend(ierr);
290d0321e0SJeremy L Thompson   }
300d0321e0SJeremy L Thompson   ierr = CeedFree(&impl->evecs); CeedChkBackend(ierr);
310d0321e0SJeremy L Thompson 
320d0321e0SJeremy L Thompson   for (CeedInt i = 0; i < impl->numein; i++) {
330d0321e0SJeremy L Thompson     ierr = CeedVectorDestroy(&impl->qvecsin[i]); CeedChkBackend(ierr);
340d0321e0SJeremy L Thompson   }
350d0321e0SJeremy L Thompson   ierr = CeedFree(&impl->qvecsin); CeedChkBackend(ierr);
360d0321e0SJeremy L Thompson 
370d0321e0SJeremy L Thompson   for (CeedInt i = 0; i < impl->numeout; i++) {
380d0321e0SJeremy L Thompson     ierr = CeedVectorDestroy(&impl->qvecsout[i]); CeedChkBackend(ierr);
390d0321e0SJeremy L Thompson   }
400d0321e0SJeremy L Thompson   ierr = CeedFree(&impl->qvecsout); CeedChkBackend(ierr);
410d0321e0SJeremy L Thompson 
420d0321e0SJeremy L Thompson   // QFunction assembly data
430d0321e0SJeremy L Thompson   for (CeedInt i=0; i<impl->qfnumactivein; i++) {
440d0321e0SJeremy L Thompson     ierr = CeedVectorDestroy(&impl->qfactivein[i]); CeedChkBackend(ierr);
450d0321e0SJeremy L Thompson   }
460d0321e0SJeremy L Thompson   ierr = CeedFree(&impl->qfactivein); CeedChkBackend(ierr);
470d0321e0SJeremy L Thompson 
480d0321e0SJeremy L Thompson   // Diag data
490d0321e0SJeremy L Thompson   if (impl->diag) {
500d0321e0SJeremy L Thompson     Ceed ceed;
510d0321e0SJeremy L Thompson     ierr = CeedOperatorGetCeed(op, &ceed); CeedChkBackend(ierr);
520d0321e0SJeremy L Thompson     CeedChk_Cu(ceed, cuModuleUnload(impl->diag->module));
530d0321e0SJeremy L Thompson     ierr = CeedFree(&impl->diag->h_emodein); CeedChkBackend(ierr);
540d0321e0SJeremy L Thompson     ierr = CeedFree(&impl->diag->h_emodeout); CeedChkBackend(ierr);
550d0321e0SJeremy L Thompson     ierr = cudaFree(impl->diag->d_emodein); CeedChk_Cu(ceed, ierr);
560d0321e0SJeremy L Thompson     ierr = cudaFree(impl->diag->d_emodeout); CeedChk_Cu(ceed, ierr);
570d0321e0SJeremy L Thompson     ierr = cudaFree(impl->diag->d_identity); CeedChk_Cu(ceed, ierr);
580d0321e0SJeremy L Thompson     ierr = cudaFree(impl->diag->d_interpin); CeedChk_Cu(ceed, ierr);
590d0321e0SJeremy L Thompson     ierr = cudaFree(impl->diag->d_interpout); CeedChk_Cu(ceed, ierr);
600d0321e0SJeremy L Thompson     ierr = cudaFree(impl->diag->d_gradin); CeedChk_Cu(ceed, ierr);
610d0321e0SJeremy L Thompson     ierr = cudaFree(impl->diag->d_gradout); CeedChk_Cu(ceed, ierr);
620d0321e0SJeremy L Thompson     ierr = CeedElemRestrictionDestroy(&impl->diag->pbdiagrstr);
630d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
640d0321e0SJeremy L Thompson     ierr = CeedVectorDestroy(&impl->diag->elemdiag); CeedChkBackend(ierr);
650d0321e0SJeremy L Thompson     ierr = CeedVectorDestroy(&impl->diag->pbelemdiag); CeedChkBackend(ierr);
660d0321e0SJeremy L Thompson   }
670d0321e0SJeremy L Thompson   ierr = CeedFree(&impl->diag); CeedChkBackend(ierr);
680d0321e0SJeremy L Thompson 
69cc132f9aSnbeams   if (impl->asmb) {
70cc132f9aSnbeams     Ceed ceed;
71cc132f9aSnbeams     ierr = CeedOperatorGetCeed(op, &ceed); CeedChkBackend(ierr);
72cc132f9aSnbeams     CeedChk_Cu(ceed, cuModuleUnload(impl->asmb->module));
73cc132f9aSnbeams     ierr = cudaFree(impl->asmb->d_B_in); CeedChk_Cu(ceed, ierr);
74cc132f9aSnbeams     ierr = cudaFree(impl->asmb->d_B_out); CeedChk_Cu(ceed, ierr);
75cc132f9aSnbeams   }
76cc132f9aSnbeams   ierr = CeedFree(&impl->asmb); CeedChkBackend(ierr);
77cc132f9aSnbeams 
780d0321e0SJeremy L Thompson   ierr = CeedFree(&impl); CeedChkBackend(ierr);
790d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
800d0321e0SJeremy L Thompson }
810d0321e0SJeremy L Thompson 
820d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
830d0321e0SJeremy L Thompson // Setup infields or outfields
840d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
850d0321e0SJeremy L Thompson static int CeedOperatorSetupFields_Cuda(CeedQFunction qf, CeedOperator op,
860d0321e0SJeremy L Thompson                                         bool isinput, CeedVector *evecs,
870d0321e0SJeremy L Thompson                                         CeedVector *qvecs, CeedInt starte,
880d0321e0SJeremy L Thompson                                         CeedInt numfields, CeedInt Q,
890d0321e0SJeremy L Thompson                                         CeedInt numelements) {
900d0321e0SJeremy L Thompson   CeedInt dim, ierr, size;
91d2643443SJeremy L Thompson   CeedSize q_size;
920d0321e0SJeremy L Thompson   Ceed ceed;
930d0321e0SJeremy L Thompson   ierr = CeedOperatorGetCeed(op, &ceed); CeedChkBackend(ierr);
940d0321e0SJeremy L Thompson   CeedBasis basis;
950d0321e0SJeremy L Thompson   CeedElemRestriction Erestrict;
960d0321e0SJeremy L Thompson   CeedOperatorField *opfields;
970d0321e0SJeremy L Thompson   CeedQFunctionField *qffields;
980d0321e0SJeremy L Thompson   CeedVector fieldvec;
990d0321e0SJeremy L Thompson   bool strided;
1000d0321e0SJeremy L Thompson   bool skiprestrict;
1010d0321e0SJeremy L Thompson 
1020d0321e0SJeremy L Thompson   if (isinput) {
1030d0321e0SJeremy L Thompson     ierr = CeedOperatorGetFields(op, NULL, &opfields, NULL, NULL);
1040d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
1050d0321e0SJeremy L Thompson     ierr = CeedQFunctionGetFields(qf, NULL, &qffields, NULL, NULL);
1060d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
1070d0321e0SJeremy L Thompson   } else {
1080d0321e0SJeremy L Thompson     ierr = CeedOperatorGetFields(op, NULL, NULL, NULL, &opfields);
1090d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
1100d0321e0SJeremy L Thompson     ierr = CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qffields);
1110d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
1120d0321e0SJeremy L Thompson   }
1130d0321e0SJeremy L Thompson 
1140d0321e0SJeremy L Thompson   // Loop over fields
1150d0321e0SJeremy L Thompson   for (CeedInt i = 0; i < numfields; i++) {
1160d0321e0SJeremy L Thompson     CeedEvalMode emode;
1170d0321e0SJeremy L Thompson     ierr = CeedQFunctionFieldGetEvalMode(qffields[i], &emode); CeedChkBackend(ierr);
1180d0321e0SJeremy L Thompson 
1190d0321e0SJeremy L Thompson     strided = false;
1200d0321e0SJeremy L Thompson     skiprestrict = false;
1210d0321e0SJeremy L Thompson     if (emode != CEED_EVAL_WEIGHT) {
1220d0321e0SJeremy L Thompson       ierr = CeedOperatorFieldGetElemRestriction(opfields[i], &Erestrict);
1230d0321e0SJeremy L Thompson       CeedChkBackend(ierr);
1240d0321e0SJeremy L Thompson 
1250d0321e0SJeremy L Thompson       // Check whether this field can skip the element restriction:
1260d0321e0SJeremy L Thompson       // must be passive input, with emode NONE, and have a strided restriction with
1270d0321e0SJeremy L Thompson       // CEED_STRIDES_BACKEND.
1280d0321e0SJeremy L Thompson 
1290d0321e0SJeremy L Thompson       // First, check whether the field is input or output:
1300d0321e0SJeremy L Thompson       if (isinput) {
1310d0321e0SJeremy L Thompson         // Check for passive input:
1320d0321e0SJeremy L Thompson         ierr = CeedOperatorFieldGetVector(opfields[i], &fieldvec); CeedChkBackend(ierr);
1330d0321e0SJeremy L Thompson         if (fieldvec != CEED_VECTOR_ACTIVE) {
1340d0321e0SJeremy L Thompson           // Check emode
1350d0321e0SJeremy L Thompson           if (emode == CEED_EVAL_NONE) {
1360d0321e0SJeremy L Thompson             // Check for strided restriction
1370d0321e0SJeremy L Thompson             ierr = CeedElemRestrictionIsStrided(Erestrict, &strided);
1380d0321e0SJeremy L Thompson             CeedChkBackend(ierr);
1390d0321e0SJeremy L Thompson             if (strided) {
1400d0321e0SJeremy L Thompson               // Check if vector is already in preferred backend ordering
1410d0321e0SJeremy L Thompson               ierr = CeedElemRestrictionHasBackendStrides(Erestrict,
1420d0321e0SJeremy L Thompson                      &skiprestrict); CeedChkBackend(ierr);
1430d0321e0SJeremy L Thompson             }
1440d0321e0SJeremy L Thompson           }
1450d0321e0SJeremy L Thompson         }
1460d0321e0SJeremy L Thompson       }
1470d0321e0SJeremy L Thompson       if (skiprestrict) {
1480d0321e0SJeremy L Thompson         // We do not need an E-Vector, but will use the input field vector's data
1490d0321e0SJeremy L Thompson         // directly in the operator application.
1500d0321e0SJeremy L Thompson         evecs[i + starte] = NULL;
1510d0321e0SJeremy L Thompson       } else {
1520d0321e0SJeremy L Thompson         ierr = CeedElemRestrictionCreateVector(Erestrict, NULL,
1530d0321e0SJeremy L Thompson                                                &evecs[i + starte]);
1540d0321e0SJeremy L Thompson         CeedChkBackend(ierr);
1550d0321e0SJeremy L Thompson       }
1560d0321e0SJeremy L Thompson     }
1570d0321e0SJeremy L Thompson 
1580d0321e0SJeremy L Thompson     switch (emode) {
1590d0321e0SJeremy L Thompson     case CEED_EVAL_NONE:
1600d0321e0SJeremy L Thompson       ierr = CeedQFunctionFieldGetSize(qffields[i], &size); CeedChkBackend(ierr);
161d2643443SJeremy L Thompson       q_size = (CeedSize)numelements * Q * size;
162d2643443SJeremy L Thompson       ierr = CeedVectorCreate(ceed, q_size, &qvecs[i]); CeedChkBackend(ierr);
1630d0321e0SJeremy L Thompson       break;
1640d0321e0SJeremy L Thompson     case CEED_EVAL_INTERP:
1650d0321e0SJeremy L Thompson       ierr = CeedQFunctionFieldGetSize(qffields[i], &size); CeedChkBackend(ierr);
166d2643443SJeremy L Thompson       q_size = (CeedSize)numelements * Q * size;
167d2643443SJeremy L Thompson       ierr = CeedVectorCreate(ceed, q_size, &qvecs[i]); CeedChkBackend(ierr);
1680d0321e0SJeremy L Thompson       break;
1690d0321e0SJeremy L Thompson     case CEED_EVAL_GRAD:
1700d0321e0SJeremy L Thompson       ierr = CeedOperatorFieldGetBasis(opfields[i], &basis); CeedChkBackend(ierr);
1710d0321e0SJeremy L Thompson       ierr = CeedQFunctionFieldGetSize(qffields[i], &size); CeedChkBackend(ierr);
1720d0321e0SJeremy L Thompson       ierr = CeedBasisGetDimension(basis, &dim); CeedChkBackend(ierr);
173d2643443SJeremy L Thompson       q_size = (CeedSize)numelements * Q * size;
174d2643443SJeremy L Thompson       ierr = CeedVectorCreate(ceed, q_size, &qvecs[i]); CeedChkBackend(ierr);
1750d0321e0SJeremy L Thompson       break;
1760d0321e0SJeremy L Thompson     case CEED_EVAL_WEIGHT: // Only on input fields
1770d0321e0SJeremy L Thompson       ierr = CeedOperatorFieldGetBasis(opfields[i], &basis); CeedChkBackend(ierr);
178d2643443SJeremy L Thompson       q_size = (CeedSize)numelements * Q;
179d2643443SJeremy L Thompson       ierr = CeedVectorCreate(ceed, q_size, &qvecs[i]); CeedChkBackend(ierr);
1800d0321e0SJeremy L Thompson       ierr = CeedBasisApply(basis, numelements, CEED_NOTRANSPOSE,
1810d0321e0SJeremy L Thompson                             CEED_EVAL_WEIGHT, NULL, qvecs[i]); CeedChkBackend(ierr);
1820d0321e0SJeremy L Thompson       break;
1830d0321e0SJeremy L Thompson     case CEED_EVAL_DIV:
1840d0321e0SJeremy L Thompson       break; // TODO: Not implemented
1850d0321e0SJeremy L Thompson     case CEED_EVAL_CURL:
1860d0321e0SJeremy L Thompson       break; // TODO: Not implemented
1870d0321e0SJeremy L Thompson     }
1880d0321e0SJeremy L Thompson   }
1890d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1900d0321e0SJeremy L Thompson }
1910d0321e0SJeremy L Thompson 
1920d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
1930d0321e0SJeremy L Thompson // CeedOperator needs to connect all the named fields (be they active or passive)
1940d0321e0SJeremy L Thompson //   to the named inputs and outputs of its CeedQFunction.
1950d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
1960d0321e0SJeremy L Thompson static int CeedOperatorSetup_Cuda(CeedOperator op) {
1970d0321e0SJeremy L Thompson   int ierr;
1980d0321e0SJeremy L Thompson   bool setupdone;
1990d0321e0SJeremy L Thompson   ierr = CeedOperatorIsSetupDone(op, &setupdone); CeedChkBackend(ierr);
2000d0321e0SJeremy L Thompson   if (setupdone)
2010d0321e0SJeremy L Thompson     return CEED_ERROR_SUCCESS;
2020d0321e0SJeremy L Thompson   Ceed ceed;
2030d0321e0SJeremy L Thompson   ierr = CeedOperatorGetCeed(op, &ceed); CeedChkBackend(ierr);
2040d0321e0SJeremy L Thompson   CeedOperator_Cuda *impl;
2050d0321e0SJeremy L Thompson   ierr = CeedOperatorGetData(op, &impl); CeedChkBackend(ierr);
2060d0321e0SJeremy L Thompson   CeedQFunction qf;
2070d0321e0SJeremy L Thompson   ierr = CeedOperatorGetQFunction(op, &qf); CeedChkBackend(ierr);
2080d0321e0SJeremy L Thompson   CeedInt Q, numelements, numinputfields, numoutputfields;
2090d0321e0SJeremy L Thompson   ierr = CeedOperatorGetNumQuadraturePoints(op, &Q); CeedChkBackend(ierr);
2100d0321e0SJeremy L Thompson   ierr = CeedOperatorGetNumElements(op, &numelements); CeedChkBackend(ierr);
2110d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
2120d0321e0SJeremy L Thompson   CeedOperatorField *opinputfields, *opoutputfields;
2130d0321e0SJeremy L Thompson   ierr = CeedOperatorGetFields(op, &numinputfields, &opinputfields,
2140d0321e0SJeremy L Thompson                                &numoutputfields, &opoutputfields);
2150d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
2160d0321e0SJeremy L Thompson   CeedQFunctionField *qfinputfields, *qfoutputfields;
2170d0321e0SJeremy L Thompson   ierr = CeedQFunctionGetFields(qf, NULL, &qfinputfields, NULL, &qfoutputfields);
2180d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
2190d0321e0SJeremy L Thompson 
2200d0321e0SJeremy L Thompson   // Allocate
2210d0321e0SJeremy L Thompson   ierr = CeedCalloc(numinputfields + numoutputfields, &impl->evecs);
2220d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
2230d0321e0SJeremy L Thompson 
2240d0321e0SJeremy L Thompson   ierr = CeedCalloc(CEED_FIELD_MAX, &impl->qvecsin); CeedChkBackend(ierr);
2250d0321e0SJeremy L Thompson   ierr = CeedCalloc(CEED_FIELD_MAX, &impl->qvecsout); CeedChkBackend(ierr);
2260d0321e0SJeremy L Thompson 
2270d0321e0SJeremy L Thompson   impl->numein = numinputfields; impl->numeout = numoutputfields;
2280d0321e0SJeremy L Thompson 
2290d0321e0SJeremy L Thompson   // Set up infield and outfield evecs and qvecs
2300d0321e0SJeremy L Thompson   // Infields
2310d0321e0SJeremy L Thompson   ierr = CeedOperatorSetupFields_Cuda(qf, op, true,
2320d0321e0SJeremy L Thompson                                       impl->evecs, impl->qvecsin, 0,
2330d0321e0SJeremy L Thompson                                       numinputfields, Q, numelements);
2340d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
2350d0321e0SJeremy L Thompson 
2360d0321e0SJeremy L Thompson   // Outfields
2370d0321e0SJeremy L Thompson   ierr = CeedOperatorSetupFields_Cuda(qf, op, false,
2380d0321e0SJeremy L Thompson                                       impl->evecs, impl->qvecsout,
2390d0321e0SJeremy L Thompson                                       numinputfields, numoutputfields, Q,
2400d0321e0SJeremy L Thompson                                       numelements); CeedChkBackend(ierr);
2410d0321e0SJeremy L Thompson 
2420d0321e0SJeremy L Thompson   ierr = CeedOperatorSetSetupDone(op); CeedChkBackend(ierr);
2430d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2440d0321e0SJeremy L Thompson }
2450d0321e0SJeremy L Thompson 
2460d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
2470d0321e0SJeremy L Thompson // Setup Operator Inputs
2480d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
2490d0321e0SJeremy L Thompson static inline int CeedOperatorSetupInputs_Cuda(CeedInt numinputfields,
2500d0321e0SJeremy L Thompson     CeedQFunctionField *qfinputfields, CeedOperatorField *opinputfields,
2510d0321e0SJeremy L Thompson     CeedVector invec, const bool skipactive, CeedScalar *edata[2*CEED_FIELD_MAX],
2520d0321e0SJeremy L Thompson     CeedOperator_Cuda *impl, CeedRequest *request) {
2530d0321e0SJeremy L Thompson   CeedInt ierr;
2540d0321e0SJeremy L Thompson   CeedEvalMode emode;
2550d0321e0SJeremy L Thompson   CeedVector vec;
2560d0321e0SJeremy L Thompson   CeedElemRestriction Erestrict;
2570d0321e0SJeremy L Thompson 
2580d0321e0SJeremy L Thompson   for (CeedInt i = 0; i < numinputfields; i++) {
2590d0321e0SJeremy L Thompson     // Get input vector
2600d0321e0SJeremy L Thompson     ierr = CeedOperatorFieldGetVector(opinputfields[i], &vec); CeedChkBackend(ierr);
2610d0321e0SJeremy L Thompson     if (vec == CEED_VECTOR_ACTIVE) {
2620d0321e0SJeremy L Thompson       if (skipactive)
2630d0321e0SJeremy L Thompson         continue;
2640d0321e0SJeremy L Thompson       else
2650d0321e0SJeremy L Thompson         vec = invec;
2660d0321e0SJeremy L Thompson     }
2670d0321e0SJeremy L Thompson 
2680d0321e0SJeremy L Thompson     ierr = CeedQFunctionFieldGetEvalMode(qfinputfields[i], &emode);
2690d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
2700d0321e0SJeremy L Thompson     if (emode == CEED_EVAL_WEIGHT) { // Skip
2710d0321e0SJeremy L Thompson     } else {
2720d0321e0SJeremy L Thompson       // Get input vector
2730d0321e0SJeremy L Thompson       ierr = CeedOperatorFieldGetVector(opinputfields[i], &vec); CeedChkBackend(ierr);
2740d0321e0SJeremy L Thompson       // Get input element restriction
2750d0321e0SJeremy L Thompson       ierr = CeedOperatorFieldGetElemRestriction(opinputfields[i], &Erestrict);
2760d0321e0SJeremy L Thompson       CeedChkBackend(ierr);
2770d0321e0SJeremy L Thompson       if (vec == CEED_VECTOR_ACTIVE)
2780d0321e0SJeremy L Thompson         vec = invec;
2790d0321e0SJeremy L Thompson       // Restrict, if necessary
2800d0321e0SJeremy L Thompson       if (!impl->evecs[i]) {
2810d0321e0SJeremy L Thompson         // No restriction for this field; read data directly from vec.
2820d0321e0SJeremy L Thompson         ierr = CeedVectorGetArrayRead(vec, CEED_MEM_DEVICE,
2830d0321e0SJeremy L Thompson                                       (const CeedScalar **) &edata[i]);
2840d0321e0SJeremy L Thompson         CeedChkBackend(ierr);
2850d0321e0SJeremy L Thompson       } else {
2860d0321e0SJeremy L Thompson         ierr = CeedElemRestrictionApply(Erestrict, CEED_NOTRANSPOSE, vec,
2870d0321e0SJeremy L Thompson                                         impl->evecs[i], request); CeedChkBackend(ierr);
2880d0321e0SJeremy L Thompson         // Get evec
2890d0321e0SJeremy L Thompson         ierr = CeedVectorGetArrayRead(impl->evecs[i], CEED_MEM_DEVICE,
2900d0321e0SJeremy L Thompson                                       (const CeedScalar **) &edata[i]);
2910d0321e0SJeremy L Thompson         CeedChkBackend(ierr);
2920d0321e0SJeremy L Thompson       }
2930d0321e0SJeremy L Thompson     }
2940d0321e0SJeremy L Thompson   }
2950d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2960d0321e0SJeremy L Thompson }
2970d0321e0SJeremy L Thompson 
2980d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
2990d0321e0SJeremy L Thompson // Input Basis Action
3000d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
3010d0321e0SJeremy L Thompson static inline int CeedOperatorInputBasis_Cuda(CeedInt numelements,
3020d0321e0SJeremy L Thompson     CeedQFunctionField *qfinputfields, CeedOperatorField *opinputfields,
3030d0321e0SJeremy L Thompson     CeedInt numinputfields, const bool skipactive,
3040d0321e0SJeremy L Thompson     CeedScalar *edata[2*CEED_FIELD_MAX],CeedOperator_Cuda *impl) {
3050d0321e0SJeremy L Thompson   CeedInt ierr;
3060d0321e0SJeremy L Thompson   CeedInt elemsize, size;
3070d0321e0SJeremy L Thompson   CeedElemRestriction Erestrict;
3080d0321e0SJeremy L Thompson   CeedEvalMode emode;
3090d0321e0SJeremy L Thompson   CeedBasis basis;
3100d0321e0SJeremy L Thompson 
3110d0321e0SJeremy L Thompson   for (CeedInt i=0; i<numinputfields; i++) {
3120d0321e0SJeremy L Thompson     // Skip active input
3130d0321e0SJeremy L Thompson     if (skipactive) {
3140d0321e0SJeremy L Thompson       CeedVector vec;
3150d0321e0SJeremy L Thompson       ierr = CeedOperatorFieldGetVector(opinputfields[i], &vec); CeedChkBackend(ierr);
3160d0321e0SJeremy L Thompson       if (vec == CEED_VECTOR_ACTIVE)
3170d0321e0SJeremy L Thompson         continue;
3180d0321e0SJeremy L Thompson     }
3190d0321e0SJeremy L Thompson     // Get elemsize, emode, size
3200d0321e0SJeremy L Thompson     ierr = CeedOperatorFieldGetElemRestriction(opinputfields[i], &Erestrict);
3210d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
3220d0321e0SJeremy L Thompson     ierr = CeedElemRestrictionGetElementSize(Erestrict, &elemsize);
3230d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
3240d0321e0SJeremy L Thompson     ierr = CeedQFunctionFieldGetEvalMode(qfinputfields[i], &emode);
3250d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
3260d0321e0SJeremy L Thompson     ierr = CeedQFunctionFieldGetSize(qfinputfields[i], &size); CeedChkBackend(ierr);
3270d0321e0SJeremy L Thompson     // Basis action
3280d0321e0SJeremy L Thompson     switch (emode) {
3290d0321e0SJeremy L Thompson     case CEED_EVAL_NONE:
3300d0321e0SJeremy L Thompson       ierr = CeedVectorSetArray(impl->qvecsin[i], CEED_MEM_DEVICE,
3310d0321e0SJeremy L Thompson                                 CEED_USE_POINTER, edata[i]); CeedChkBackend(ierr);
3320d0321e0SJeremy L Thompson       break;
3330d0321e0SJeremy L Thompson     case CEED_EVAL_INTERP:
3340d0321e0SJeremy L Thompson       ierr = CeedOperatorFieldGetBasis(opinputfields[i], &basis);
3350d0321e0SJeremy L Thompson       CeedChkBackend(ierr);
3360d0321e0SJeremy L Thompson       ierr = CeedBasisApply(basis, numelements, CEED_NOTRANSPOSE,
3370d0321e0SJeremy L Thompson                             CEED_EVAL_INTERP, impl->evecs[i],
3380d0321e0SJeremy L Thompson                             impl->qvecsin[i]); CeedChkBackend(ierr);
3390d0321e0SJeremy L Thompson       break;
3400d0321e0SJeremy L Thompson     case CEED_EVAL_GRAD:
3410d0321e0SJeremy L Thompson       ierr = CeedOperatorFieldGetBasis(opinputfields[i], &basis);
3420d0321e0SJeremy L Thompson       CeedChkBackend(ierr);
3430d0321e0SJeremy L Thompson       ierr = CeedBasisApply(basis, numelements, CEED_NOTRANSPOSE,
3440d0321e0SJeremy L Thompson                             CEED_EVAL_GRAD, impl->evecs[i],
3450d0321e0SJeremy L Thompson                             impl->qvecsin[i]); CeedChkBackend(ierr);
3460d0321e0SJeremy L Thompson       break;
3470d0321e0SJeremy L Thompson     case CEED_EVAL_WEIGHT:
3480d0321e0SJeremy L Thompson       break; // No action
3490d0321e0SJeremy L Thompson     case CEED_EVAL_DIV:
3500d0321e0SJeremy L Thompson       break; // TODO: Not implemented
3510d0321e0SJeremy L Thompson     case CEED_EVAL_CURL:
3520d0321e0SJeremy L Thompson       break; // TODO: Not implemented
3530d0321e0SJeremy L Thompson     }
3540d0321e0SJeremy L Thompson   }
3550d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3560d0321e0SJeremy L Thompson }
3570d0321e0SJeremy L Thompson 
3580d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
3590d0321e0SJeremy L Thompson // Restore Input Vectors
3600d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
3610d0321e0SJeremy L Thompson static inline int CeedOperatorRestoreInputs_Cuda(CeedInt numinputfields,
3620d0321e0SJeremy L Thompson     CeedQFunctionField *qfinputfields, CeedOperatorField *opinputfields,
3630d0321e0SJeremy L Thompson     const bool skipactive, CeedScalar *edata[2*CEED_FIELD_MAX],
3640d0321e0SJeremy L Thompson     CeedOperator_Cuda *impl) {
3650d0321e0SJeremy L Thompson   CeedInt ierr;
3660d0321e0SJeremy L Thompson   CeedEvalMode emode;
3670d0321e0SJeremy L Thompson   CeedVector vec;
3680d0321e0SJeremy L Thompson 
3690d0321e0SJeremy L Thompson   for (CeedInt i = 0; i < numinputfields; i++) {
3700d0321e0SJeremy L Thompson     // Skip active input
3710d0321e0SJeremy L Thompson     if (skipactive) {
3720d0321e0SJeremy L Thompson       ierr = CeedOperatorFieldGetVector(opinputfields[i], &vec); CeedChkBackend(ierr);
3730d0321e0SJeremy L Thompson       if (vec == CEED_VECTOR_ACTIVE)
3740d0321e0SJeremy L Thompson         continue;
3750d0321e0SJeremy L Thompson     }
3760d0321e0SJeremy L Thompson     ierr = CeedQFunctionFieldGetEvalMode(qfinputfields[i], &emode);
3770d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
3780d0321e0SJeremy L Thompson     if (emode == CEED_EVAL_WEIGHT) { // Skip
3790d0321e0SJeremy L Thompson     } else {
3800d0321e0SJeremy L Thompson       if (!impl->evecs[i]) {  // This was a skiprestrict case
3810d0321e0SJeremy L Thompson         ierr = CeedOperatorFieldGetVector(opinputfields[i], &vec); CeedChkBackend(ierr);
3820d0321e0SJeremy L Thompson         ierr = CeedVectorRestoreArrayRead(vec,
3830d0321e0SJeremy L Thompson                                           (const CeedScalar **)&edata[i]);
3840d0321e0SJeremy L Thompson         CeedChkBackend(ierr);
3850d0321e0SJeremy L Thompson       } else {
3860d0321e0SJeremy L Thompson         ierr = CeedVectorRestoreArrayRead(impl->evecs[i],
3870d0321e0SJeremy L Thompson                                           (const CeedScalar **) &edata[i]);
3880d0321e0SJeremy L Thompson         CeedChkBackend(ierr);
3890d0321e0SJeremy L Thompson       }
3900d0321e0SJeremy L Thompson     }
3910d0321e0SJeremy L Thompson   }
3920d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3930d0321e0SJeremy L Thompson }
3940d0321e0SJeremy L Thompson 
3950d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
3960d0321e0SJeremy L Thompson // Apply and add to output
3970d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
3980d0321e0SJeremy L Thompson static int CeedOperatorApplyAdd_Cuda(CeedOperator op, CeedVector invec,
3990d0321e0SJeremy L Thompson                                      CeedVector outvec, CeedRequest *request) {
4000d0321e0SJeremy L Thompson   int ierr;
4010d0321e0SJeremy L Thompson   CeedOperator_Cuda *impl;
4020d0321e0SJeremy L Thompson   ierr = CeedOperatorGetData(op, &impl); CeedChkBackend(ierr);
4030d0321e0SJeremy L Thompson   CeedQFunction qf;
4040d0321e0SJeremy L Thompson   ierr = CeedOperatorGetQFunction(op, &qf); CeedChkBackend(ierr);
4050d0321e0SJeremy L Thompson   CeedInt Q, numelements, elemsize, numinputfields, numoutputfields, size;
4060d0321e0SJeremy L Thompson   ierr = CeedOperatorGetNumQuadraturePoints(op, &Q); CeedChkBackend(ierr);
4070d0321e0SJeremy L Thompson   ierr = CeedOperatorGetNumElements(op, &numelements); CeedChkBackend(ierr);
4080d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
4090d0321e0SJeremy L Thompson   CeedOperatorField *opinputfields, *opoutputfields;
4100d0321e0SJeremy L Thompson   ierr = CeedOperatorGetFields(op, &numinputfields, &opinputfields,
4110d0321e0SJeremy L Thompson                                &numoutputfields, &opoutputfields);
4120d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
4130d0321e0SJeremy L Thompson   CeedQFunctionField *qfinputfields, *qfoutputfields;
4140d0321e0SJeremy L Thompson   ierr = CeedQFunctionGetFields(qf, NULL, &qfinputfields, NULL, &qfoutputfields);
4150d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
4160d0321e0SJeremy L Thompson   CeedEvalMode emode;
4170d0321e0SJeremy L Thompson   CeedVector vec;
4180d0321e0SJeremy L Thompson   CeedBasis basis;
4190d0321e0SJeremy L Thompson   CeedElemRestriction Erestrict;
4200d0321e0SJeremy L Thompson   CeedScalar *edata[2*CEED_FIELD_MAX] = {0};
4210d0321e0SJeremy L Thompson 
4220d0321e0SJeremy L Thompson   // Setup
4230d0321e0SJeremy L Thompson   ierr = CeedOperatorSetup_Cuda(op); CeedChkBackend(ierr);
4240d0321e0SJeremy L Thompson 
4250d0321e0SJeremy L Thompson   // Input Evecs and Restriction
4260d0321e0SJeremy L Thompson   ierr = CeedOperatorSetupInputs_Cuda(numinputfields, qfinputfields,
4270d0321e0SJeremy L Thompson                                       opinputfields, invec, false, edata,
4280d0321e0SJeremy L Thompson                                       impl, request); CeedChkBackend(ierr);
4290d0321e0SJeremy L Thompson 
4300d0321e0SJeremy L Thompson   // Input basis apply if needed
4310d0321e0SJeremy L Thompson   ierr = CeedOperatorInputBasis_Cuda(numelements, qfinputfields, opinputfields,
4320d0321e0SJeremy L Thompson                                      numinputfields, false, edata, impl);
4330d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
4340d0321e0SJeremy L Thompson 
4350d0321e0SJeremy L Thompson   // Output pointers, as necessary
4360d0321e0SJeremy L Thompson   for (CeedInt i = 0; i < numoutputfields; i++) {
4370d0321e0SJeremy L Thompson     ierr = CeedQFunctionFieldGetEvalMode(qfoutputfields[i], &emode);
4380d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
4390d0321e0SJeremy L Thompson     if (emode == CEED_EVAL_NONE) {
4400d0321e0SJeremy L Thompson       // Set the output Q-Vector to use the E-Vector data directly.
4410d0321e0SJeremy L Thompson       ierr = CeedVectorGetArrayWrite(impl->evecs[i + impl->numein], CEED_MEM_DEVICE,
4420d0321e0SJeremy L Thompson                                      &edata[i + numinputfields]); CeedChkBackend(ierr);
4430d0321e0SJeremy L Thompson       ierr = CeedVectorSetArray(impl->qvecsout[i], CEED_MEM_DEVICE,
4440d0321e0SJeremy L Thompson                                 CEED_USE_POINTER, edata[i + numinputfields]);
4450d0321e0SJeremy L Thompson       CeedChkBackend(ierr);
4460d0321e0SJeremy L Thompson     }
4470d0321e0SJeremy L Thompson   }
4480d0321e0SJeremy L Thompson 
4490d0321e0SJeremy L Thompson   // Q function
4500d0321e0SJeremy L Thompson   ierr = CeedQFunctionApply(qf, numelements * Q, impl->qvecsin, impl->qvecsout);
4510d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
4520d0321e0SJeremy L Thompson 
4530d0321e0SJeremy L Thompson   // Output basis apply if needed
4540d0321e0SJeremy L Thompson   for (CeedInt i = 0; i < numoutputfields; i++) {
4550d0321e0SJeremy L Thompson     // Get elemsize, emode, size
4560d0321e0SJeremy L Thompson     ierr = CeedOperatorFieldGetElemRestriction(opoutputfields[i], &Erestrict);
4570d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
4580d0321e0SJeremy L Thompson     ierr = CeedElemRestrictionGetElementSize(Erestrict, &elemsize);
4590d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
4600d0321e0SJeremy L Thompson     ierr = CeedQFunctionFieldGetEvalMode(qfoutputfields[i], &emode);
4610d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
4620d0321e0SJeremy L Thompson     ierr = CeedQFunctionFieldGetSize(qfoutputfields[i], &size);
4630d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
4640d0321e0SJeremy L Thompson     // Basis action
4650d0321e0SJeremy L Thompson     switch (emode) {
4660d0321e0SJeremy L Thompson     case CEED_EVAL_NONE:
4670d0321e0SJeremy L Thompson       break;
4680d0321e0SJeremy L Thompson     case CEED_EVAL_INTERP:
4690d0321e0SJeremy L Thompson       ierr = CeedOperatorFieldGetBasis(opoutputfields[i], &basis);
4700d0321e0SJeremy L Thompson       CeedChkBackend(ierr);
4710d0321e0SJeremy L Thompson       ierr = CeedBasisApply(basis, numelements, CEED_TRANSPOSE,
4720d0321e0SJeremy L Thompson                             CEED_EVAL_INTERP, impl->qvecsout[i],
4730d0321e0SJeremy L Thompson                             impl->evecs[i + impl->numein]); CeedChkBackend(ierr);
4740d0321e0SJeremy L Thompson       break;
4750d0321e0SJeremy L Thompson     case CEED_EVAL_GRAD:
4760d0321e0SJeremy L Thompson       ierr = CeedOperatorFieldGetBasis(opoutputfields[i], &basis);
4770d0321e0SJeremy L Thompson       CeedChkBackend(ierr);
4780d0321e0SJeremy L Thompson       ierr = CeedBasisApply(basis, numelements, CEED_TRANSPOSE,
4790d0321e0SJeremy L Thompson                             CEED_EVAL_GRAD, impl->qvecsout[i],
4800d0321e0SJeremy L Thompson                             impl->evecs[i + impl->numein]); CeedChkBackend(ierr);
4810d0321e0SJeremy L Thompson       break;
4820d0321e0SJeremy L Thompson     // LCOV_EXCL_START
4830d0321e0SJeremy L Thompson     case CEED_EVAL_WEIGHT: {
4840d0321e0SJeremy L Thompson       Ceed ceed;
4850d0321e0SJeremy L Thompson       ierr = CeedOperatorGetCeed(op, &ceed); CeedChkBackend(ierr);
4860d0321e0SJeremy L Thompson       return CeedError(ceed, CEED_ERROR_BACKEND,
4870d0321e0SJeremy L Thompson                        "CEED_EVAL_WEIGHT cannot be an output evaluation mode");
4880d0321e0SJeremy L Thompson       break; // Should not occur
4890d0321e0SJeremy L Thompson     }
4900d0321e0SJeremy L Thompson     case CEED_EVAL_DIV:
4910d0321e0SJeremy L Thompson       break; // TODO: Not implemented
4920d0321e0SJeremy L Thompson     case CEED_EVAL_CURL:
4930d0321e0SJeremy L Thompson       break; // TODO: Not implemented
4940d0321e0SJeremy L Thompson       // LCOV_EXCL_STOP
4950d0321e0SJeremy L Thompson     }
4960d0321e0SJeremy L Thompson   }
4970d0321e0SJeremy L Thompson 
4980d0321e0SJeremy L Thompson   // Output restriction
4990d0321e0SJeremy L Thompson   for (CeedInt i = 0; i < numoutputfields; i++) {
5000d0321e0SJeremy L Thompson     // Restore evec
5010d0321e0SJeremy L Thompson     ierr = CeedQFunctionFieldGetEvalMode(qfoutputfields[i], &emode);
5020d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
5030d0321e0SJeremy L Thompson     if (emode == CEED_EVAL_NONE) {
5040d0321e0SJeremy L Thompson       ierr = CeedVectorRestoreArray(impl->evecs[i+impl->numein],
5050d0321e0SJeremy L Thompson                                     &edata[i + numinputfields]);
5060d0321e0SJeremy L Thompson       CeedChkBackend(ierr);
5070d0321e0SJeremy L Thompson     }
5080d0321e0SJeremy L Thompson     // Get output vector
5090d0321e0SJeremy L Thompson     ierr = CeedOperatorFieldGetVector(opoutputfields[i], &vec);
5100d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
5110d0321e0SJeremy L Thompson     // Restrict
5120d0321e0SJeremy L Thompson     ierr = CeedOperatorFieldGetElemRestriction(opoutputfields[i], &Erestrict);
5130d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
5140d0321e0SJeremy L Thompson     // Active
5150d0321e0SJeremy L Thompson     if (vec == CEED_VECTOR_ACTIVE)
5160d0321e0SJeremy L Thompson       vec = outvec;
5170d0321e0SJeremy L Thompson 
5180d0321e0SJeremy L Thompson     ierr = CeedElemRestrictionApply(Erestrict, CEED_TRANSPOSE,
5190d0321e0SJeremy L Thompson                                     impl->evecs[i + impl->numein], vec,
5200d0321e0SJeremy L Thompson                                     request); CeedChkBackend(ierr);
5210d0321e0SJeremy L Thompson   }
5220d0321e0SJeremy L Thompson 
5230d0321e0SJeremy L Thompson   // Restore input arrays
5240d0321e0SJeremy L Thompson   ierr = CeedOperatorRestoreInputs_Cuda(numinputfields, qfinputfields,
5250d0321e0SJeremy L Thompson                                         opinputfields, false, edata, impl);
5260d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
5270d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
5280d0321e0SJeremy L Thompson }
5290d0321e0SJeremy L Thompson 
5300d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
5310d0321e0SJeremy L Thompson // Core code for assembling linear QFunction
5320d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
5330d0321e0SJeremy L Thompson static inline int CeedOperatorLinearAssembleQFunctionCore_Cuda(CeedOperator op,
5340d0321e0SJeremy L Thompson     bool build_objects, CeedVector *assembled, CeedElemRestriction *rstr,
5350d0321e0SJeremy L Thompson     CeedRequest *request) {
5360d0321e0SJeremy L Thompson   int ierr;
5370d0321e0SJeremy L Thompson   CeedOperator_Cuda *impl;
5380d0321e0SJeremy L Thompson   ierr = CeedOperatorGetData(op, &impl); CeedChkBackend(ierr);
5390d0321e0SJeremy L Thompson   CeedQFunction qf;
5400d0321e0SJeremy L Thompson   ierr = CeedOperatorGetQFunction(op, &qf); CeedChkBackend(ierr);
5410d0321e0SJeremy L Thompson   CeedInt Q, numelements, numinputfields, numoutputfields, size;
542d2643443SJeremy L Thompson   CeedSize q_size;
5430d0321e0SJeremy L Thompson   ierr = CeedOperatorGetNumQuadraturePoints(op, &Q); CeedChkBackend(ierr);
5440d0321e0SJeremy L Thompson   ierr = CeedOperatorGetNumElements(op, &numelements); CeedChkBackend(ierr);
5450d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
5460d0321e0SJeremy L Thompson   CeedOperatorField *opinputfields, *opoutputfields;
5470d0321e0SJeremy L Thompson   ierr = CeedOperatorGetFields(op, &numinputfields, &opinputfields,
5480d0321e0SJeremy L Thompson                                &numoutputfields, &opoutputfields);
5490d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
5500d0321e0SJeremy L Thompson   CeedQFunctionField *qfinputfields, *qfoutputfields;
5510d0321e0SJeremy L Thompson   ierr = CeedQFunctionGetFields(qf, NULL, &qfinputfields, NULL, &qfoutputfields);
5520d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
5530d0321e0SJeremy L Thompson   CeedVector vec;
5540d0321e0SJeremy L Thompson   CeedInt numactivein = impl->qfnumactivein, numactiveout = impl->qfnumactiveout;
5550d0321e0SJeremy L Thompson   CeedVector *activein = impl->qfactivein;
5560d0321e0SJeremy L Thompson   CeedScalar *a, *tmp;
5570d0321e0SJeremy L Thompson   Ceed ceed, ceedparent;
5580d0321e0SJeremy L Thompson   ierr = CeedOperatorGetCeed(op, &ceed); CeedChkBackend(ierr);
5590d0321e0SJeremy L Thompson   ierr = CeedGetOperatorFallbackParentCeed(ceed, &ceedparent);
5600d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
5610d0321e0SJeremy L Thompson   ceedparent = ceedparent ? ceedparent : ceed;
5620d0321e0SJeremy L Thompson   CeedScalar *edata[2*CEED_FIELD_MAX];
5630d0321e0SJeremy L Thompson 
5640d0321e0SJeremy L Thompson   // Setup
5650d0321e0SJeremy L Thompson   ierr = CeedOperatorSetup_Cuda(op); CeedChkBackend(ierr);
5660d0321e0SJeremy L Thompson 
5670d0321e0SJeremy L Thompson   // Check for identity
5680d0321e0SJeremy L Thompson   bool identityqf;
5690d0321e0SJeremy L Thompson   ierr = CeedQFunctionIsIdentity(qf, &identityqf); CeedChkBackend(ierr);
5700d0321e0SJeremy L Thompson   if (identityqf)
5710d0321e0SJeremy L Thompson     // LCOV_EXCL_START
5720d0321e0SJeremy L Thompson     return CeedError(ceed, CEED_ERROR_BACKEND,
5730d0321e0SJeremy L Thompson                      "Assembling identity QFunctions not supported");
5740d0321e0SJeremy L Thompson   // LCOV_EXCL_STOP
5750d0321e0SJeremy L Thompson 
5760d0321e0SJeremy L Thompson   // Input Evecs and Restriction
5770d0321e0SJeremy L Thompson   ierr = CeedOperatorSetupInputs_Cuda(numinputfields, qfinputfields,
5780d0321e0SJeremy L Thompson                                       opinputfields, NULL, true, edata,
5790d0321e0SJeremy L Thompson                                       impl, request); CeedChkBackend(ierr);
5800d0321e0SJeremy L Thompson 
5810d0321e0SJeremy L Thompson   // Count number of active input fields
5820d0321e0SJeremy L Thompson   if (!numactivein) {
5830d0321e0SJeremy L Thompson     for (CeedInt i=0; i<numinputfields; i++) {
5840d0321e0SJeremy L Thompson       // Get input vector
5850d0321e0SJeremy L Thompson       ierr = CeedOperatorFieldGetVector(opinputfields[i], &vec); CeedChkBackend(ierr);
5860d0321e0SJeremy L Thompson       // Check if active input
5870d0321e0SJeremy L Thompson       if (vec == CEED_VECTOR_ACTIVE) {
5880d0321e0SJeremy L Thompson         ierr = CeedQFunctionFieldGetSize(qfinputfields[i], &size); CeedChkBackend(ierr);
5890d0321e0SJeremy L Thompson         ierr = CeedVectorSetValue(impl->qvecsin[i], 0.0); CeedChkBackend(ierr);
5900d0321e0SJeremy L Thompson         ierr = CeedVectorGetArray(impl->qvecsin[i], CEED_MEM_DEVICE, &tmp);
5910d0321e0SJeremy L Thompson         CeedChkBackend(ierr);
5920d0321e0SJeremy L Thompson         ierr = CeedRealloc(numactivein + size, &activein); CeedChkBackend(ierr);
5930d0321e0SJeremy L Thompson         for (CeedInt field = 0; field < size; field++) {
594d2643443SJeremy L Thompson           q_size = (CeedSize)Q*numelements;
595d2643443SJeremy L Thompson           ierr = CeedVectorCreate(ceed, q_size, &activein[numactivein+field]);
596d2643443SJeremy L Thompson           CeedChkBackend(ierr);
5970d0321e0SJeremy L Thompson           ierr = CeedVectorSetArray(activein[numactivein+field], CEED_MEM_DEVICE,
5980d0321e0SJeremy L Thompson                                     CEED_USE_POINTER, &tmp[field*Q*numelements]);
5990d0321e0SJeremy L Thompson           CeedChkBackend(ierr);
6000d0321e0SJeremy L Thompson         }
6010d0321e0SJeremy L Thompson         numactivein += size;
6020d0321e0SJeremy L Thompson         ierr = CeedVectorRestoreArray(impl->qvecsin[i], &tmp); CeedChkBackend(ierr);
6030d0321e0SJeremy L Thompson       }
6040d0321e0SJeremy L Thompson     }
6050d0321e0SJeremy L Thompson     impl->qfnumactivein = numactivein;
6060d0321e0SJeremy L Thompson     impl->qfactivein = activein;
6070d0321e0SJeremy L Thompson   }
6080d0321e0SJeremy L Thompson 
6090d0321e0SJeremy L Thompson   // Count number of active output fields
6100d0321e0SJeremy L Thompson   if (!numactiveout) {
6110d0321e0SJeremy L Thompson     for (CeedInt i=0; i<numoutputfields; i++) {
6120d0321e0SJeremy L Thompson       // Get output vector
6130d0321e0SJeremy L Thompson       ierr = CeedOperatorFieldGetVector(opoutputfields[i], &vec);
6140d0321e0SJeremy L Thompson       CeedChkBackend(ierr);
6150d0321e0SJeremy L Thompson       // Check if active output
6160d0321e0SJeremy L Thompson       if (vec == CEED_VECTOR_ACTIVE) {
6170d0321e0SJeremy L Thompson         ierr = CeedQFunctionFieldGetSize(qfoutputfields[i], &size);
6180d0321e0SJeremy L Thompson         CeedChkBackend(ierr);
6190d0321e0SJeremy L Thompson         numactiveout += size;
6200d0321e0SJeremy L Thompson       }
6210d0321e0SJeremy L Thompson     }
6220d0321e0SJeremy L Thompson     impl->qfnumactiveout = numactiveout;
6230d0321e0SJeremy L Thompson   }
6240d0321e0SJeremy L Thompson 
6250d0321e0SJeremy L Thompson   // Check sizes
6260d0321e0SJeremy L Thompson   if (!numactivein || !numactiveout)
6270d0321e0SJeremy L Thompson     // LCOV_EXCL_START
6280d0321e0SJeremy L Thompson     return CeedError(ceed, CEED_ERROR_BACKEND,
6290d0321e0SJeremy L Thompson                      "Cannot assemble QFunction without active inputs "
6300d0321e0SJeremy L Thompson                      "and outputs");
6310d0321e0SJeremy L Thompson   // LCOV_EXCL_STOP
6320d0321e0SJeremy L Thompson 
6330d0321e0SJeremy L Thompson   // Build objects if needed
6340d0321e0SJeremy L Thompson   if (build_objects) {
6350d0321e0SJeremy L Thompson     // Create output restriction
6360d0321e0SJeremy L Thompson     CeedInt strides[3] = {1, numelements*Q, Q}; /* *NOPAD* */
6370d0321e0SJeremy L Thompson     ierr = CeedElemRestrictionCreateStrided(ceedparent, numelements, Q,
6380d0321e0SJeremy L Thompson                                             numactivein*numactiveout,
6390d0321e0SJeremy L Thompson                                             numactivein*numactiveout*numelements*Q,
6400d0321e0SJeremy L Thompson                                             strides, rstr); CeedChkBackend(ierr);
6410d0321e0SJeremy L Thompson     // Create assembled vector
642d2643443SJeremy L Thompson     CeedSize l_size = (CeedSize)numelements*Q*numactivein*numactiveout;
643d2643443SJeremy L Thompson     ierr = CeedVectorCreate(ceedparent, l_size, assembled); CeedChkBackend(ierr);
6440d0321e0SJeremy L Thompson   }
6450d0321e0SJeremy L Thompson   ierr = CeedVectorSetValue(*assembled, 0.0); CeedChkBackend(ierr);
6460d0321e0SJeremy L Thompson   ierr = CeedVectorGetArray(*assembled, CEED_MEM_DEVICE, &a);
6470d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
6480d0321e0SJeremy L Thompson 
6490d0321e0SJeremy L Thompson   // Input basis apply
6500d0321e0SJeremy L Thompson   ierr = CeedOperatorInputBasis_Cuda(numelements, qfinputfields, opinputfields,
6510d0321e0SJeremy L Thompson                                      numinputfields, true, edata, impl);
6520d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
6530d0321e0SJeremy L Thompson 
6540d0321e0SJeremy L Thompson   // Assemble QFunction
6550d0321e0SJeremy L Thompson   for (CeedInt in=0; in<numactivein; in++) {
6560d0321e0SJeremy L Thompson     // Set Inputs
6570d0321e0SJeremy L Thompson     ierr = CeedVectorSetValue(activein[in], 1.0); CeedChkBackend(ierr);
6580d0321e0SJeremy L Thompson     if (numactivein > 1) {
6590d0321e0SJeremy L Thompson       ierr = CeedVectorSetValue(activein[(in+numactivein-1)%numactivein],
6600d0321e0SJeremy L Thompson                                 0.0); CeedChkBackend(ierr);
6610d0321e0SJeremy L Thompson     }
6620d0321e0SJeremy L Thompson     // Set Outputs
6630d0321e0SJeremy L Thompson     for (CeedInt out=0; out<numoutputfields; out++) {
6640d0321e0SJeremy L Thompson       // Get output vector
6650d0321e0SJeremy L Thompson       ierr = CeedOperatorFieldGetVector(opoutputfields[out], &vec);
6660d0321e0SJeremy L Thompson       CeedChkBackend(ierr);
6670d0321e0SJeremy L Thompson       // Check if active output
6680d0321e0SJeremy L Thompson       if (vec == CEED_VECTOR_ACTIVE) {
6690d0321e0SJeremy L Thompson         CeedVectorSetArray(impl->qvecsout[out], CEED_MEM_DEVICE,
6700d0321e0SJeremy L Thompson                            CEED_USE_POINTER, a); CeedChkBackend(ierr);
6710d0321e0SJeremy L Thompson         ierr = CeedQFunctionFieldGetSize(qfoutputfields[out], &size);
6720d0321e0SJeremy L Thompson         CeedChkBackend(ierr);
6730d0321e0SJeremy L Thompson         a += size*Q*numelements; // Advance the pointer by the size of the output
6740d0321e0SJeremy L Thompson       }
6750d0321e0SJeremy L Thompson     }
6760d0321e0SJeremy L Thompson     // Apply QFunction
6770d0321e0SJeremy L Thompson     ierr = CeedQFunctionApply(qf, Q*numelements, impl->qvecsin, impl->qvecsout);
6780d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
6790d0321e0SJeremy L Thompson   }
6800d0321e0SJeremy L Thompson 
6810d0321e0SJeremy L Thompson   // Un-set output Qvecs to prevent accidental overwrite of Assembled
6820d0321e0SJeremy L Thompson   for (CeedInt out=0; out<numoutputfields; out++) {
6830d0321e0SJeremy L Thompson     // Get output vector
6840d0321e0SJeremy L Thompson     ierr = CeedOperatorFieldGetVector(opoutputfields[out], &vec);
6850d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
6860d0321e0SJeremy L Thompson     // Check if active output
6870d0321e0SJeremy L Thompson     if (vec == CEED_VECTOR_ACTIVE) {
6880d0321e0SJeremy L Thompson       ierr = CeedVectorTakeArray(impl->qvecsout[out], CEED_MEM_DEVICE, NULL);
6890d0321e0SJeremy L Thompson       CeedChkBackend(ierr);
6900d0321e0SJeremy L Thompson     }
6910d0321e0SJeremy L Thompson   }
6920d0321e0SJeremy L Thompson 
6930d0321e0SJeremy L Thompson   // Restore input arrays
6940d0321e0SJeremy L Thompson   ierr = CeedOperatorRestoreInputs_Cuda(numinputfields, qfinputfields,
6950d0321e0SJeremy L Thompson                                         opinputfields, true, edata, impl);
6960d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
6970d0321e0SJeremy L Thompson 
6980d0321e0SJeremy L Thompson   // Restore output
6990d0321e0SJeremy L Thompson   ierr = CeedVectorRestoreArray(*assembled, &a); CeedChkBackend(ierr);
7000d0321e0SJeremy L Thompson 
7010d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
7020d0321e0SJeremy L Thompson }
7030d0321e0SJeremy L Thompson 
7040d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
7050d0321e0SJeremy L Thompson // Assemble Linear QFunction
7060d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
7070d0321e0SJeremy L Thompson static int CeedOperatorLinearAssembleQFunction_Cuda(CeedOperator op,
7080d0321e0SJeremy L Thompson     CeedVector *assembled, CeedElemRestriction *rstr, CeedRequest *request) {
7090d0321e0SJeremy L Thompson   return CeedOperatorLinearAssembleQFunctionCore_Cuda(op, true, assembled, rstr,
7100d0321e0SJeremy L Thompson          request);
7110d0321e0SJeremy L Thompson }
7120d0321e0SJeremy L Thompson 
7130d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
7140d0321e0SJeremy L Thompson // Update Assembled Linear QFunction
7150d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
7160d0321e0SJeremy L Thompson static int CeedOperatorLinearAssembleQFunctionUpdate_Cuda(CeedOperator op,
7170d0321e0SJeremy L Thompson     CeedVector assembled, CeedElemRestriction rstr, CeedRequest *request) {
7180d0321e0SJeremy L Thompson   return CeedOperatorLinearAssembleQFunctionCore_Cuda(op, false, &assembled,
7190d0321e0SJeremy L Thompson          &rstr, request);
7200d0321e0SJeremy L Thompson }
7210d0321e0SJeremy L Thompson 
7220d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
7230d0321e0SJeremy L Thompson // Diagonal assembly kernels
7240d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
7250d0321e0SJeremy L Thompson // *INDENT-OFF*
7260d0321e0SJeremy L Thompson static const char *diagonalkernels = QUOTE(
7270d0321e0SJeremy L Thompson 
7280d0321e0SJeremy L Thompson typedef enum {
7290d0321e0SJeremy L Thompson   /// Perform no evaluation (either because there is no data or it is already at
7300d0321e0SJeremy L Thompson   /// quadrature points)
7310d0321e0SJeremy L Thompson   CEED_EVAL_NONE   = 0,
7320d0321e0SJeremy L Thompson   /// Interpolate from nodes to quadrature points
7330d0321e0SJeremy L Thompson   CEED_EVAL_INTERP = 1,
7340d0321e0SJeremy L Thompson   /// Evaluate gradients at quadrature points from input in a nodal basis
7350d0321e0SJeremy L Thompson   CEED_EVAL_GRAD   = 2,
7360d0321e0SJeremy L Thompson   /// Evaluate divergence at quadrature points from input in a nodal basis
7370d0321e0SJeremy L Thompson   CEED_EVAL_DIV    = 4,
7380d0321e0SJeremy L Thompson   /// Evaluate curl at quadrature points from input in a nodal basis
7390d0321e0SJeremy L Thompson   CEED_EVAL_CURL   = 8,
7400d0321e0SJeremy L Thompson   /// Using no input, evaluate quadrature weights on the reference element
7410d0321e0SJeremy L Thompson   CEED_EVAL_WEIGHT = 16,
7420d0321e0SJeremy L Thompson } CeedEvalMode;
7430d0321e0SJeremy L Thompson 
7440d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
7450d0321e0SJeremy L Thompson // Get Basis Emode Pointer
7460d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
7470d0321e0SJeremy L Thompson extern "C" __device__ void CeedOperatorGetBasisPointer_Cuda(const CeedScalar **basisptr,
7480d0321e0SJeremy L Thompson     CeedEvalMode emode, const CeedScalar *identity, const CeedScalar *interp,
7490d0321e0SJeremy L Thompson     const CeedScalar *grad) {
7500d0321e0SJeremy L Thompson   switch (emode) {
7510d0321e0SJeremy L Thompson   case CEED_EVAL_NONE:
7520d0321e0SJeremy L Thompson     *basisptr = identity;
7530d0321e0SJeremy L Thompson     break;
7540d0321e0SJeremy L Thompson   case CEED_EVAL_INTERP:
7550d0321e0SJeremy L Thompson     *basisptr = interp;
7560d0321e0SJeremy L Thompson     break;
7570d0321e0SJeremy L Thompson   case CEED_EVAL_GRAD:
7580d0321e0SJeremy L Thompson     *basisptr = grad;
7590d0321e0SJeremy L Thompson     break;
7600d0321e0SJeremy L Thompson   case CEED_EVAL_WEIGHT:
7610d0321e0SJeremy L Thompson   case CEED_EVAL_DIV:
7620d0321e0SJeremy L Thompson   case CEED_EVAL_CURL:
7630d0321e0SJeremy L Thompson     break; // Caught by QF Assembly
7640d0321e0SJeremy L Thompson   }
7650d0321e0SJeremy L Thompson }
7660d0321e0SJeremy L Thompson 
7670d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
7680d0321e0SJeremy L Thompson // Core code for diagonal assembly
7690d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
7700d0321e0SJeremy L Thompson __device__ void diagonalCore(const CeedInt nelem,
7710d0321e0SJeremy L Thompson     const CeedScalar maxnorm, const bool pointBlock,
7720d0321e0SJeremy L Thompson     const CeedScalar *identity,
7730d0321e0SJeremy L Thompson     const CeedScalar *interpin, const CeedScalar *gradin,
7740d0321e0SJeremy L Thompson     const CeedScalar *interpout, const CeedScalar *gradout,
7750d0321e0SJeremy L Thompson     const CeedEvalMode *emodein, const CeedEvalMode *emodeout,
7760d0321e0SJeremy L Thompson     const CeedScalar *__restrict__ assembledqfarray,
7770d0321e0SJeremy L Thompson     CeedScalar *__restrict__ elemdiagarray) {
7780d0321e0SJeremy L Thompson   const int tid = threadIdx.x; // running with P threads, tid is evec node
7790d0321e0SJeremy L Thompson   const CeedScalar qfvaluebound = maxnorm*1e-12;
7800d0321e0SJeremy L Thompson 
7810d0321e0SJeremy L Thompson   // Compute the diagonal of B^T D B
7820d0321e0SJeremy L Thompson   // Each element
7830d0321e0SJeremy L Thompson   for (CeedInt e = blockIdx.x*blockDim.z + threadIdx.z; e < nelem;
7840d0321e0SJeremy L Thompson        e += gridDim.x*blockDim.z) {
7850d0321e0SJeremy L Thompson     CeedInt dout = -1;
7860d0321e0SJeremy L Thompson     // Each basis eval mode pair
7870d0321e0SJeremy L Thompson     for (CeedInt eout = 0; eout < NUMEMODEOUT; eout++) {
7880d0321e0SJeremy L Thompson       const CeedScalar *bt = NULL;
7890d0321e0SJeremy L Thompson       if (emodeout[eout] == CEED_EVAL_GRAD)
7900d0321e0SJeremy L Thompson         dout += 1;
7910d0321e0SJeremy L Thompson       CeedOperatorGetBasisPointer_Cuda(&bt, emodeout[eout], identity, interpout,
7920d0321e0SJeremy L Thompson                                       &gradout[dout*NQPTS*NNODES]);
7930d0321e0SJeremy L Thompson       CeedInt din = -1;
7940d0321e0SJeremy L Thompson       for (CeedInt ein = 0; ein < NUMEMODEIN; ein++) {
7950d0321e0SJeremy L Thompson         const CeedScalar *b = NULL;
7960d0321e0SJeremy L Thompson         if (emodein[ein] == CEED_EVAL_GRAD)
7970d0321e0SJeremy L Thompson           din += 1;
7980d0321e0SJeremy L Thompson         CeedOperatorGetBasisPointer_Cuda(&b, emodein[ein], identity, interpin,
7990d0321e0SJeremy L Thompson                                         &gradin[din*NQPTS*NNODES]);
8000d0321e0SJeremy L Thompson         // Each component
8010d0321e0SJeremy L Thompson         for (CeedInt compOut = 0; compOut < NCOMP; compOut++) {
8020d0321e0SJeremy L Thompson           // Each qpoint/node pair
8030d0321e0SJeremy L Thompson           if (pointBlock) {
8040d0321e0SJeremy L Thompson             // Point Block Diagonal
8050d0321e0SJeremy L Thompson             for (CeedInt compIn = 0; compIn < NCOMP; compIn++) {
8060d0321e0SJeremy L Thompson               CeedScalar evalue = 0.;
8070d0321e0SJeremy L Thompson               for (CeedInt q = 0; q < NQPTS; q++) {
8080d0321e0SJeremy L Thompson                 const CeedScalar qfvalue =
8090d0321e0SJeremy L Thompson                   assembledqfarray[((((ein*NCOMP+compIn)*NUMEMODEOUT+eout)*
8100d0321e0SJeremy L Thompson                                      NCOMP+compOut)*nelem+e)*NQPTS+q];
8110d0321e0SJeremy L Thompson                 if (abs(qfvalue) > qfvaluebound)
8120d0321e0SJeremy L Thompson                   evalue += bt[q*NNODES+tid] * qfvalue * b[q*NNODES+tid];
8130d0321e0SJeremy L Thompson               }
8140d0321e0SJeremy L Thompson               elemdiagarray[((compOut*NCOMP+compIn)*nelem+e)*NNODES+tid] += evalue;
8150d0321e0SJeremy L Thompson             }
8160d0321e0SJeremy L Thompson           } else {
8170d0321e0SJeremy L Thompson             // Diagonal Only
8180d0321e0SJeremy L Thompson             CeedScalar evalue = 0.;
8190d0321e0SJeremy L Thompson             for (CeedInt q = 0; q < NQPTS; q++) {
8200d0321e0SJeremy L Thompson               const CeedScalar qfvalue =
8210d0321e0SJeremy L Thompson                 assembledqfarray[((((ein*NCOMP+compOut)*NUMEMODEOUT+eout)*
8220d0321e0SJeremy L Thompson                                    NCOMP+compOut)*nelem+e)*NQPTS+q];
8230d0321e0SJeremy L Thompson               if (abs(qfvalue) > qfvaluebound)
8240d0321e0SJeremy L Thompson                 evalue += bt[q*NNODES+tid] * qfvalue * b[q*NNODES+tid];
8250d0321e0SJeremy L Thompson             }
8260d0321e0SJeremy L Thompson             elemdiagarray[(compOut*nelem+e)*NNODES+tid] += evalue;
8270d0321e0SJeremy L Thompson           }
8280d0321e0SJeremy L Thompson         }
8290d0321e0SJeremy L Thompson       }
8300d0321e0SJeremy L Thompson     }
8310d0321e0SJeremy L Thompson   }
8320d0321e0SJeremy L Thompson }
8330d0321e0SJeremy L Thompson 
8340d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
8350d0321e0SJeremy L Thompson // Linear diagonal
8360d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
8370d0321e0SJeremy L Thompson extern "C" __global__ void linearDiagonal(const CeedInt nelem,
8380d0321e0SJeremy L Thompson     const CeedScalar maxnorm, const CeedScalar *identity,
8390d0321e0SJeremy L Thompson     const CeedScalar *interpin, const CeedScalar *gradin,
8400d0321e0SJeremy L Thompson     const CeedScalar *interpout, const CeedScalar *gradout,
8410d0321e0SJeremy L Thompson     const CeedEvalMode *emodein, const CeedEvalMode *emodeout,
8420d0321e0SJeremy L Thompson     const CeedScalar *__restrict__ assembledqfarray,
8430d0321e0SJeremy L Thompson     CeedScalar *__restrict__ elemdiagarray) {
8440d0321e0SJeremy L Thompson   diagonalCore(nelem, maxnorm, false, identity, interpin, gradin, interpout,
8450d0321e0SJeremy L Thompson                gradout, emodein, emodeout, assembledqfarray, elemdiagarray);
8460d0321e0SJeremy L Thompson }
8470d0321e0SJeremy L Thompson 
8480d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
8490d0321e0SJeremy L Thompson // Linear point block diagonal
8500d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
8510d0321e0SJeremy L Thompson extern "C" __global__ void linearPointBlockDiagonal(const CeedInt nelem,
8520d0321e0SJeremy L Thompson     const CeedScalar maxnorm, const CeedScalar *identity,
8530d0321e0SJeremy L Thompson     const CeedScalar *interpin, const CeedScalar *gradin,
8540d0321e0SJeremy L Thompson     const CeedScalar *interpout, const CeedScalar *gradout,
8550d0321e0SJeremy L Thompson     const CeedEvalMode *emodein, const CeedEvalMode *emodeout,
8560d0321e0SJeremy L Thompson     const CeedScalar *__restrict__ assembledqfarray,
8570d0321e0SJeremy L Thompson     CeedScalar *__restrict__ elemdiagarray) {
8580d0321e0SJeremy L Thompson   diagonalCore(nelem, maxnorm, true, identity, interpin, gradin, interpout,
8590d0321e0SJeremy L Thompson                gradout, emodein, emodeout, assembledqfarray, elemdiagarray);
8600d0321e0SJeremy L Thompson }
8610d0321e0SJeremy L Thompson 
8620d0321e0SJeremy L Thompson );
8630d0321e0SJeremy L Thompson // *INDENT-ON*
8640d0321e0SJeremy L Thompson 
8650d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
8660d0321e0SJeremy L Thompson // Create point block restriction
8670d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
8680d0321e0SJeremy L Thompson static int CreatePBRestriction(CeedElemRestriction rstr,
8690d0321e0SJeremy L Thompson                                CeedElemRestriction *pbRstr) {
8700d0321e0SJeremy L Thompson   int ierr;
8710d0321e0SJeremy L Thompson   Ceed ceed;
8720d0321e0SJeremy L Thompson   ierr = CeedElemRestrictionGetCeed(rstr, &ceed); CeedChkBackend(ierr);
8730d0321e0SJeremy L Thompson   const CeedInt *offsets;
8740d0321e0SJeremy L Thompson   ierr = CeedElemRestrictionGetOffsets(rstr, CEED_MEM_HOST, &offsets);
8750d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
8760d0321e0SJeremy L Thompson 
8770d0321e0SJeremy L Thompson   // Expand offsets
8780d0321e0SJeremy L Thompson   CeedInt nelem, ncomp, elemsize, compstride, max = 1, *pbOffsets;
8790d0321e0SJeremy L Thompson   ierr = CeedElemRestrictionGetNumElements(rstr, &nelem); CeedChkBackend(ierr);
8800d0321e0SJeremy L Thompson   ierr = CeedElemRestrictionGetNumComponents(rstr, &ncomp); CeedChkBackend(ierr);
8810d0321e0SJeremy L Thompson   ierr = CeedElemRestrictionGetElementSize(rstr, &elemsize); CeedChkBackend(ierr);
8820d0321e0SJeremy L Thompson   ierr = CeedElemRestrictionGetCompStride(rstr, &compstride);
8830d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
8840d0321e0SJeremy L Thompson   CeedInt shift = ncomp;
8850d0321e0SJeremy L Thompson   if (compstride != 1)
8860d0321e0SJeremy L Thompson     shift *= ncomp;
8870d0321e0SJeremy L Thompson   ierr = CeedCalloc(nelem*elemsize, &pbOffsets); CeedChkBackend(ierr);
8880d0321e0SJeremy L Thompson   for (CeedInt i = 0; i < nelem*elemsize; i++) {
8890d0321e0SJeremy L Thompson     pbOffsets[i] = offsets[i]*shift;
8900d0321e0SJeremy L Thompson     if (pbOffsets[i] > max)
8910d0321e0SJeremy L Thompson       max = pbOffsets[i];
8920d0321e0SJeremy L Thompson   }
8930d0321e0SJeremy L Thompson 
8940d0321e0SJeremy L Thompson   // Create new restriction
8950d0321e0SJeremy L Thompson   ierr = CeedElemRestrictionCreate(ceed, nelem, elemsize, ncomp*ncomp, 1,
8960d0321e0SJeremy L Thompson                                    max + ncomp*ncomp, CEED_MEM_HOST,
8970d0321e0SJeremy L Thompson                                    CEED_OWN_POINTER, pbOffsets, pbRstr);
8980d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
8990d0321e0SJeremy L Thompson 
9000d0321e0SJeremy L Thompson   // Cleanup
9010d0321e0SJeremy L Thompson   ierr = CeedElemRestrictionRestoreOffsets(rstr, &offsets); CeedChkBackend(ierr);
9020d0321e0SJeremy L Thompson 
9030d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
9040d0321e0SJeremy L Thompson }
9050d0321e0SJeremy L Thompson 
9060d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
9070d0321e0SJeremy L Thompson // Assemble diagonal setup
9080d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
9090d0321e0SJeremy L Thompson static inline int CeedOperatorAssembleDiagonalSetup_Cuda(CeedOperator op,
9100d0321e0SJeremy L Thompson     const bool pointBlock) {
9110d0321e0SJeremy L Thompson   int ierr;
9120d0321e0SJeremy L Thompson   Ceed ceed;
9130d0321e0SJeremy L Thompson   ierr = CeedOperatorGetCeed(op, &ceed); CeedChkBackend(ierr);
9140d0321e0SJeremy L Thompson   CeedQFunction qf;
9150d0321e0SJeremy L Thompson   ierr = CeedOperatorGetQFunction(op, &qf); CeedChkBackend(ierr);
9160d0321e0SJeremy L Thompson   CeedInt numinputfields, numoutputfields;
9170d0321e0SJeremy L Thompson   ierr = CeedQFunctionGetNumArgs(qf, &numinputfields, &numoutputfields);
9180d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
9190d0321e0SJeremy L Thompson 
9200d0321e0SJeremy L Thompson   // Determine active input basis
9210d0321e0SJeremy L Thompson   CeedOperatorField *opfields;
9220d0321e0SJeremy L Thompson   CeedQFunctionField *qffields;
9230d0321e0SJeremy L Thompson   ierr = CeedOperatorGetFields(op, NULL, &opfields, NULL, NULL);
9240d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
9250d0321e0SJeremy L Thompson   ierr = CeedQFunctionGetFields(qf, NULL, &qffields, NULL, NULL);
9260d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
9270d0321e0SJeremy L Thompson   CeedInt numemodein = 0, ncomp = 0, dim = 1;
9280d0321e0SJeremy L Thompson   CeedEvalMode *emodein = NULL;
9290d0321e0SJeremy L Thompson   CeedBasis basisin = NULL;
9300d0321e0SJeremy L Thompson   CeedElemRestriction rstrin = NULL;
9310d0321e0SJeremy L Thompson   for (CeedInt i = 0; i < numinputfields; i++) {
9320d0321e0SJeremy L Thompson     CeedVector vec;
9330d0321e0SJeremy L Thompson     ierr = CeedOperatorFieldGetVector(opfields[i], &vec); CeedChkBackend(ierr);
9340d0321e0SJeremy L Thompson     if (vec == CEED_VECTOR_ACTIVE) {
9350d0321e0SJeremy L Thompson       CeedElemRestriction rstr;
9360d0321e0SJeremy L Thompson       ierr = CeedOperatorFieldGetBasis(opfields[i], &basisin); CeedChkBackend(ierr);
9370d0321e0SJeremy L Thompson       ierr = CeedBasisGetNumComponents(basisin, &ncomp); CeedChkBackend(ierr);
9380d0321e0SJeremy L Thompson       ierr = CeedBasisGetDimension(basisin, &dim); CeedChkBackend(ierr);
9390d0321e0SJeremy L Thompson       ierr = CeedOperatorFieldGetElemRestriction(opfields[i], &rstr);
9400d0321e0SJeremy L Thompson       CeedChkBackend(ierr);
9410d0321e0SJeremy L Thompson       if (rstrin && rstrin != rstr)
9420d0321e0SJeremy L Thompson         // LCOV_EXCL_START
9430d0321e0SJeremy L Thompson         return CeedError(ceed, CEED_ERROR_BACKEND,
9440d0321e0SJeremy L Thompson                          "Multi-field non-composite operator diagonal assembly not supported");
9450d0321e0SJeremy L Thompson       // LCOV_EXCL_STOP
9460d0321e0SJeremy L Thompson       rstrin = rstr;
9470d0321e0SJeremy L Thompson       CeedEvalMode emode;
9480d0321e0SJeremy L Thompson       ierr = CeedQFunctionFieldGetEvalMode(qffields[i], &emode);
9490d0321e0SJeremy L Thompson       CeedChkBackend(ierr);
9500d0321e0SJeremy L Thompson       switch (emode) {
9510d0321e0SJeremy L Thompson       case CEED_EVAL_NONE:
9520d0321e0SJeremy L Thompson       case CEED_EVAL_INTERP:
9530d0321e0SJeremy L Thompson         ierr = CeedRealloc(numemodein + 1, &emodein); CeedChkBackend(ierr);
9540d0321e0SJeremy L Thompson         emodein[numemodein] = emode;
9550d0321e0SJeremy L Thompson         numemodein += 1;
9560d0321e0SJeremy L Thompson         break;
9570d0321e0SJeremy L Thompson       case CEED_EVAL_GRAD:
9580d0321e0SJeremy L Thompson         ierr = CeedRealloc(numemodein + dim, &emodein); CeedChkBackend(ierr);
9590d0321e0SJeremy L Thompson         for (CeedInt d = 0; d < dim; d++)
9600d0321e0SJeremy L Thompson           emodein[numemodein+d] = emode;
9610d0321e0SJeremy L Thompson         numemodein += dim;
9620d0321e0SJeremy L Thompson         break;
9630d0321e0SJeremy L Thompson       case CEED_EVAL_WEIGHT:
9640d0321e0SJeremy L Thompson       case CEED_EVAL_DIV:
9650d0321e0SJeremy L Thompson       case CEED_EVAL_CURL:
9660d0321e0SJeremy L Thompson         break; // Caught by QF Assembly
9670d0321e0SJeremy L Thompson       }
9680d0321e0SJeremy L Thompson     }
9690d0321e0SJeremy L Thompson   }
9700d0321e0SJeremy L Thompson 
9710d0321e0SJeremy L Thompson   // Determine active output basis
9720d0321e0SJeremy L Thompson   ierr = CeedOperatorGetFields(op, NULL, NULL, NULL, &opfields);
9730d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
9740d0321e0SJeremy L Thompson   ierr = CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qffields);
9750d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
9760d0321e0SJeremy L Thompson   CeedInt numemodeout = 0;
9770d0321e0SJeremy L Thompson   CeedEvalMode *emodeout = NULL;
9780d0321e0SJeremy L Thompson   CeedBasis basisout = NULL;
9790d0321e0SJeremy L Thompson   CeedElemRestriction rstrout = NULL;
9800d0321e0SJeremy L Thompson   for (CeedInt i = 0; i < numoutputfields; i++) {
9810d0321e0SJeremy L Thompson     CeedVector vec;
9820d0321e0SJeremy L Thompson     ierr = CeedOperatorFieldGetVector(opfields[i], &vec); CeedChkBackend(ierr);
9830d0321e0SJeremy L Thompson     if (vec == CEED_VECTOR_ACTIVE) {
9840d0321e0SJeremy L Thompson       CeedElemRestriction rstr;
9850d0321e0SJeremy L Thompson       ierr = CeedOperatorFieldGetBasis(opfields[i], &basisout); CeedChkBackend(ierr);
9860d0321e0SJeremy L Thompson       ierr = CeedOperatorFieldGetElemRestriction(opfields[i], &rstr);
9870d0321e0SJeremy L Thompson       CeedChkBackend(ierr);
9880d0321e0SJeremy L Thompson       if (rstrout && rstrout != rstr)
9890d0321e0SJeremy L Thompson         // LCOV_EXCL_START
9900d0321e0SJeremy L Thompson         return CeedError(ceed, CEED_ERROR_BACKEND,
9910d0321e0SJeremy L Thompson                          "Multi-field non-composite operator diagonal assembly not supported");
9920d0321e0SJeremy L Thompson       // LCOV_EXCL_STOP
9930d0321e0SJeremy L Thompson       rstrout = rstr;
9940d0321e0SJeremy L Thompson       CeedEvalMode emode;
9950d0321e0SJeremy L Thompson       ierr = CeedQFunctionFieldGetEvalMode(qffields[i], &emode); CeedChkBackend(ierr);
9960d0321e0SJeremy L Thompson       switch (emode) {
9970d0321e0SJeremy L Thompson       case CEED_EVAL_NONE:
9980d0321e0SJeremy L Thompson       case CEED_EVAL_INTERP:
9990d0321e0SJeremy L Thompson         ierr = CeedRealloc(numemodeout + 1, &emodeout); CeedChkBackend(ierr);
10000d0321e0SJeremy L Thompson         emodeout[numemodeout] = emode;
10010d0321e0SJeremy L Thompson         numemodeout += 1;
10020d0321e0SJeremy L Thompson         break;
10030d0321e0SJeremy L Thompson       case CEED_EVAL_GRAD:
10040d0321e0SJeremy L Thompson         ierr = CeedRealloc(numemodeout + dim, &emodeout); CeedChkBackend(ierr);
10050d0321e0SJeremy L Thompson         for (CeedInt d = 0; d < dim; d++)
10060d0321e0SJeremy L Thompson           emodeout[numemodeout+d] = emode;
10070d0321e0SJeremy L Thompson         numemodeout += dim;
10080d0321e0SJeremy L Thompson         break;
10090d0321e0SJeremy L Thompson       case CEED_EVAL_WEIGHT:
10100d0321e0SJeremy L Thompson       case CEED_EVAL_DIV:
10110d0321e0SJeremy L Thompson       case CEED_EVAL_CURL:
10120d0321e0SJeremy L Thompson         break; // Caught by QF Assembly
10130d0321e0SJeremy L Thompson       }
10140d0321e0SJeremy L Thompson     }
10150d0321e0SJeremy L Thompson   }
10160d0321e0SJeremy L Thompson 
10170d0321e0SJeremy L Thompson   // Operator data struct
10180d0321e0SJeremy L Thompson   CeedOperator_Cuda *impl;
10190d0321e0SJeremy L Thompson   ierr = CeedOperatorGetData(op, &impl); CeedChkBackend(ierr);
10200d0321e0SJeremy L Thompson   ierr = CeedCalloc(1, &impl->diag); CeedChkBackend(ierr);
10210d0321e0SJeremy L Thompson   CeedOperatorDiag_Cuda *diag = impl->diag;
10220d0321e0SJeremy L Thompson   diag->basisin = basisin;
10230d0321e0SJeremy L Thompson   diag->basisout = basisout;
10240d0321e0SJeremy L Thompson   diag->h_emodein = emodein;
10250d0321e0SJeremy L Thompson   diag->h_emodeout = emodeout;
10260d0321e0SJeremy L Thompson   diag->numemodein = numemodein;
10270d0321e0SJeremy L Thompson   diag->numemodeout = numemodeout;
10280d0321e0SJeremy L Thompson 
10290d0321e0SJeremy L Thompson   // Assemble kernel
10300d0321e0SJeremy L Thompson   CeedInt nnodes, nqpts;
10310d0321e0SJeremy L Thompson   ierr = CeedBasisGetNumNodes(basisin, &nnodes); CeedChkBackend(ierr);
10320d0321e0SJeremy L Thompson   ierr = CeedBasisGetNumQuadraturePoints(basisin, &nqpts); CeedChkBackend(ierr);
10330d0321e0SJeremy L Thompson   diag->nnodes = nnodes;
10340d0321e0SJeremy L Thompson   ierr = CeedCompileCuda(ceed, diagonalkernels, &diag->module, 5,
10350d0321e0SJeremy L Thompson                          "NUMEMODEIN", numemodein,
10360d0321e0SJeremy L Thompson                          "NUMEMODEOUT", numemodeout,
10370d0321e0SJeremy L Thompson                          "NNODES", nnodes,
10380d0321e0SJeremy L Thompson                          "NQPTS", nqpts,
10390d0321e0SJeremy L Thompson                          "NCOMP", ncomp
10400d0321e0SJeremy L Thompson                         ); CeedChk_Cu(ceed, ierr);
10410d0321e0SJeremy L Thompson   ierr = CeedGetKernelCuda(ceed, diag->module, "linearDiagonal",
10420d0321e0SJeremy L Thompson                            &diag->linearDiagonal); CeedChk_Cu(ceed, ierr);
10430d0321e0SJeremy L Thompson   ierr = CeedGetKernelCuda(ceed, diag->module, "linearPointBlockDiagonal",
10440d0321e0SJeremy L Thompson                            &diag->linearPointBlock);
10450d0321e0SJeremy L Thompson   CeedChk_Cu(ceed, ierr);
10460d0321e0SJeremy L Thompson 
10470d0321e0SJeremy L Thompson   // Basis matrices
10480d0321e0SJeremy L Thompson   const CeedInt qBytes = nqpts * sizeof(CeedScalar);
10490d0321e0SJeremy L Thompson   const CeedInt iBytes = qBytes * nnodes;
10500d0321e0SJeremy L Thompson   const CeedInt gBytes = qBytes * nnodes * dim;
10510d0321e0SJeremy L Thompson   const CeedInt eBytes = sizeof(CeedEvalMode);
10520d0321e0SJeremy L Thompson   const CeedScalar *interpin, *interpout, *gradin, *gradout;
10530d0321e0SJeremy L Thompson 
10540d0321e0SJeremy L Thompson   // CEED_EVAL_NONE
10550d0321e0SJeremy L Thompson   CeedScalar *identity = NULL;
10560d0321e0SJeremy L Thompson   bool evalNone = false;
10570d0321e0SJeremy L Thompson   for (CeedInt i=0; i<numemodein; i++)
10580d0321e0SJeremy L Thompson     evalNone = evalNone || (emodein[i] == CEED_EVAL_NONE);
10590d0321e0SJeremy L Thompson   for (CeedInt i=0; i<numemodeout; i++)
10600d0321e0SJeremy L Thompson     evalNone = evalNone || (emodeout[i] == CEED_EVAL_NONE);
10610d0321e0SJeremy L Thompson   if (evalNone) {
10620d0321e0SJeremy L Thompson     ierr = CeedCalloc(nqpts*nnodes, &identity); CeedChkBackend(ierr);
10630d0321e0SJeremy L Thompson     for (CeedInt i=0; i<(nnodes<nqpts?nnodes:nqpts); i++)
10640d0321e0SJeremy L Thompson       identity[i*nnodes+i] = 1.0;
10650d0321e0SJeremy L Thompson     ierr = cudaMalloc((void **)&diag->d_identity, iBytes); CeedChk_Cu(ceed, ierr);
10660d0321e0SJeremy L Thompson     ierr = cudaMemcpy(diag->d_identity, identity, iBytes,
10670d0321e0SJeremy L Thompson                       cudaMemcpyHostToDevice); CeedChk_Cu(ceed, ierr);
10680d0321e0SJeremy L Thompson   }
10690d0321e0SJeremy L Thompson 
10700d0321e0SJeremy L Thompson   // CEED_EVAL_INTERP
10710d0321e0SJeremy L Thompson   ierr = CeedBasisGetInterp(basisin, &interpin); CeedChkBackend(ierr);
10720d0321e0SJeremy L Thompson   ierr = cudaMalloc((void **)&diag->d_interpin, iBytes); CeedChk_Cu(ceed, ierr);
10730d0321e0SJeremy L Thompson   ierr = cudaMemcpy(diag->d_interpin, interpin, iBytes,
10740d0321e0SJeremy L Thompson                     cudaMemcpyHostToDevice); CeedChk_Cu(ceed, ierr);
10750d0321e0SJeremy L Thompson   ierr = CeedBasisGetInterp(basisout, &interpout); CeedChkBackend(ierr);
10760d0321e0SJeremy L Thompson   ierr = cudaMalloc((void **)&diag->d_interpout, iBytes); CeedChk_Cu(ceed, ierr);
10770d0321e0SJeremy L Thompson   ierr = cudaMemcpy(diag->d_interpout, interpout, iBytes,
10780d0321e0SJeremy L Thompson                     cudaMemcpyHostToDevice); CeedChk_Cu(ceed, ierr);
10790d0321e0SJeremy L Thompson 
10800d0321e0SJeremy L Thompson   // CEED_EVAL_GRAD
10810d0321e0SJeremy L Thompson   ierr = CeedBasisGetGrad(basisin, &gradin); CeedChkBackend(ierr);
10820d0321e0SJeremy L Thompson   ierr = cudaMalloc((void **)&diag->d_gradin, gBytes); CeedChk_Cu(ceed, ierr);
10830d0321e0SJeremy L Thompson   ierr = cudaMemcpy(diag->d_gradin, gradin, gBytes,
10840d0321e0SJeremy L Thompson                     cudaMemcpyHostToDevice); CeedChk_Cu(ceed, ierr);
10850d0321e0SJeremy L Thompson   ierr = CeedBasisGetGrad(basisout, &gradout); CeedChkBackend(ierr);
10860d0321e0SJeremy L Thompson   ierr = cudaMalloc((void **)&diag->d_gradout, gBytes); CeedChk_Cu(ceed, ierr);
10870d0321e0SJeremy L Thompson   ierr = cudaMemcpy(diag->d_gradout, gradout, gBytes,
10880d0321e0SJeremy L Thompson                     cudaMemcpyHostToDevice); CeedChk_Cu(ceed, ierr);
10890d0321e0SJeremy L Thompson 
10900d0321e0SJeremy L Thompson   // Arrays of emodes
10910d0321e0SJeremy L Thompson   ierr = cudaMalloc((void **)&diag->d_emodein, numemodein * eBytes);
10920d0321e0SJeremy L Thompson   CeedChk_Cu(ceed, ierr);
10930d0321e0SJeremy L Thompson   ierr = cudaMemcpy(diag->d_emodein, emodein, numemodein * eBytes,
10940d0321e0SJeremy L Thompson                     cudaMemcpyHostToDevice); CeedChk_Cu(ceed, ierr);
10950d0321e0SJeremy L Thompson   ierr = cudaMalloc((void **)&diag->d_emodeout, numemodeout * eBytes);
10960d0321e0SJeremy L Thompson   CeedChk_Cu(ceed, ierr);
10970d0321e0SJeremy L Thompson   ierr = cudaMemcpy(diag->d_emodeout, emodeout, numemodeout * eBytes,
10980d0321e0SJeremy L Thompson                     cudaMemcpyHostToDevice); CeedChk_Cu(ceed, ierr);
10990d0321e0SJeremy L Thompson 
11000d0321e0SJeremy L Thompson   // Restriction
11010d0321e0SJeremy L Thompson   diag->diagrstr = rstrout;
11020d0321e0SJeremy L Thompson 
11030d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
11040d0321e0SJeremy L Thompson }
11050d0321e0SJeremy L Thompson 
11060d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
11070d0321e0SJeremy L Thompson // Assemble diagonal common code
11080d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
11090d0321e0SJeremy L Thompson static inline int CeedOperatorAssembleDiagonalCore_Cuda(CeedOperator op,
11100d0321e0SJeremy L Thompson     CeedVector assembled, CeedRequest *request, const bool pointBlock) {
11110d0321e0SJeremy L Thompson   int ierr;
11120d0321e0SJeremy L Thompson   Ceed ceed;
11130d0321e0SJeremy L Thompson   ierr = CeedOperatorGetCeed(op, &ceed); CeedChkBackend(ierr);
11140d0321e0SJeremy L Thompson   CeedOperator_Cuda *impl;
11150d0321e0SJeremy L Thompson   ierr = CeedOperatorGetData(op, &impl); CeedChkBackend(ierr);
11160d0321e0SJeremy L Thompson 
11170d0321e0SJeremy L Thompson   // Assemble QFunction
11180d0321e0SJeremy L Thompson   CeedVector assembledqf;
11190d0321e0SJeremy L Thompson   CeedElemRestriction rstr;
11200d0321e0SJeremy L Thompson   ierr = CeedOperatorLinearAssembleQFunctionBuildOrUpdate(op, &assembledqf,
11210d0321e0SJeremy L Thompson          &rstr, request); CeedChkBackend(ierr);
11220d0321e0SJeremy L Thompson   ierr = CeedElemRestrictionDestroy(&rstr); CeedChkBackend(ierr);
11230d0321e0SJeremy L Thompson   CeedScalar maxnorm = 0;
11240d0321e0SJeremy L Thompson   ierr = CeedVectorNorm(assembledqf, CEED_NORM_MAX, &maxnorm);
11250d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
11260d0321e0SJeremy L Thompson 
11270d0321e0SJeremy L Thompson   // Setup
11280d0321e0SJeremy L Thompson   if (!impl->diag) {
11290d0321e0SJeremy L Thompson     ierr = CeedOperatorAssembleDiagonalSetup_Cuda(op, pointBlock);
11300d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
11310d0321e0SJeremy L Thompson   }
11320d0321e0SJeremy L Thompson   CeedOperatorDiag_Cuda *diag = impl->diag;
11330d0321e0SJeremy L Thompson   assert(diag != NULL);
11340d0321e0SJeremy L Thompson 
11350d0321e0SJeremy L Thompson   // Restriction
11360d0321e0SJeremy L Thompson   if (pointBlock && !diag->pbdiagrstr) {
11370d0321e0SJeremy L Thompson     CeedElemRestriction pbdiagrstr;
11380d0321e0SJeremy L Thompson     ierr = CreatePBRestriction(diag->diagrstr, &pbdiagrstr); CeedChkBackend(ierr);
11390d0321e0SJeremy L Thompson     diag->pbdiagrstr = pbdiagrstr;
11400d0321e0SJeremy L Thompson   }
11410d0321e0SJeremy L Thompson   CeedElemRestriction diagrstr = pointBlock ? diag->pbdiagrstr : diag->diagrstr;
11420d0321e0SJeremy L Thompson 
11430d0321e0SJeremy L Thompson   // Create diagonal vector
11440d0321e0SJeremy L Thompson   CeedVector elemdiag = pointBlock ? diag->pbelemdiag : diag->elemdiag;
11450d0321e0SJeremy L Thompson   if (!elemdiag) {
11460d0321e0SJeremy L Thompson     ierr = CeedElemRestrictionCreateVector(diagrstr, NULL, &elemdiag);
11470d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
11480d0321e0SJeremy L Thompson     if (pointBlock)
11490d0321e0SJeremy L Thompson       diag->pbelemdiag = elemdiag;
11500d0321e0SJeremy L Thompson     else
11510d0321e0SJeremy L Thompson       diag->elemdiag = elemdiag;
11520d0321e0SJeremy L Thompson   }
11530d0321e0SJeremy L Thompson   ierr = CeedVectorSetValue(elemdiag, 0.0); CeedChkBackend(ierr);
11540d0321e0SJeremy L Thompson 
11550d0321e0SJeremy L Thompson   // Assemble element operator diagonals
11560d0321e0SJeremy L Thompson   CeedScalar *elemdiagarray;
11570d0321e0SJeremy L Thompson   const CeedScalar *assembledqfarray;
11580d0321e0SJeremy L Thompson   ierr = CeedVectorGetArray(elemdiag, CEED_MEM_DEVICE, &elemdiagarray);
11590d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
11600d0321e0SJeremy L Thompson   ierr = CeedVectorGetArrayRead(assembledqf, CEED_MEM_DEVICE, &assembledqfarray);
11610d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
11620d0321e0SJeremy L Thompson   CeedInt nelem;
11630d0321e0SJeremy L Thompson   ierr = CeedElemRestrictionGetNumElements(diagrstr, &nelem);
11640d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
11650d0321e0SJeremy L Thompson 
11660d0321e0SJeremy L Thompson   // Compute the diagonal of B^T D B
11670d0321e0SJeremy L Thompson   int elemsPerBlock = 1;
11680d0321e0SJeremy L Thompson   int grid = nelem/elemsPerBlock+((nelem/elemsPerBlock*elemsPerBlock<nelem)?1:0);
11690d0321e0SJeremy L Thompson   void *args[] = {(void *) &nelem, (void *) &maxnorm, &diag->d_identity,
11700d0321e0SJeremy L Thompson                   &diag->d_interpin, &diag->d_gradin, &diag->d_interpout,
11710d0321e0SJeremy L Thompson                   &diag->d_gradout, &diag->d_emodein, &diag->d_emodeout,
11720d0321e0SJeremy L Thompson                   &assembledqfarray, &elemdiagarray
11730d0321e0SJeremy L Thompson                  };
11740d0321e0SJeremy L Thompson   if (pointBlock) {
11750d0321e0SJeremy L Thompson     ierr = CeedRunKernelDimCuda(ceed, diag->linearPointBlock, grid,
11760d0321e0SJeremy L Thompson                                 diag->nnodes, 1, elemsPerBlock, args);
11770d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
11780d0321e0SJeremy L Thompson   } else {
11790d0321e0SJeremy L Thompson     ierr = CeedRunKernelDimCuda(ceed, diag->linearDiagonal, grid,
11800d0321e0SJeremy L Thompson                                 diag->nnodes, 1, elemsPerBlock, args);
11810d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
11820d0321e0SJeremy L Thompson   }
11830d0321e0SJeremy L Thompson 
11840d0321e0SJeremy L Thompson   // Restore arrays
11850d0321e0SJeremy L Thompson   ierr = CeedVectorRestoreArray(elemdiag, &elemdiagarray); CeedChkBackend(ierr);
11860d0321e0SJeremy L Thompson   ierr = CeedVectorRestoreArrayRead(assembledqf, &assembledqfarray);
11870d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
11880d0321e0SJeremy L Thompson 
11890d0321e0SJeremy L Thompson   // Assemble local operator diagonal
11900d0321e0SJeremy L Thompson   ierr = CeedElemRestrictionApply(diagrstr, CEED_TRANSPOSE, elemdiag,
11910d0321e0SJeremy L Thompson                                   assembled, request); CeedChkBackend(ierr);
11920d0321e0SJeremy L Thompson 
11930d0321e0SJeremy L Thompson   // Cleanup
11940d0321e0SJeremy L Thompson   ierr = CeedVectorDestroy(&assembledqf); CeedChkBackend(ierr);
11950d0321e0SJeremy L Thompson 
11960d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
11970d0321e0SJeremy L Thompson }
11980d0321e0SJeremy L Thompson 
11990d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
12000d0321e0SJeremy L Thompson // Assemble composite diagonal common code
12010d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
12020d0321e0SJeremy L Thompson static inline int CeedOperatorLinearAssembleAddDiagonalCompositeCore_Cuda(
12030d0321e0SJeremy L Thompson   CeedOperator op, CeedVector assembled, CeedRequest *request,
12040d0321e0SJeremy L Thompson   const bool pointBlock) {
12050d0321e0SJeremy L Thompson   int ierr;
12060d0321e0SJeremy L Thompson   CeedInt numSub;
12070d0321e0SJeremy L Thompson   CeedOperator *subOperators;
12080d0321e0SJeremy L Thompson   ierr = CeedOperatorGetNumSub(op, &numSub); CeedChkBackend(ierr);
12090d0321e0SJeremy L Thompson   ierr = CeedOperatorGetSubList(op, &subOperators); CeedChkBackend(ierr);
12100d0321e0SJeremy L Thompson   for (CeedInt i = 0; i < numSub; i++) {
12110d0321e0SJeremy L Thompson     ierr = CeedOperatorAssembleDiagonalCore_Cuda(subOperators[i], assembled,
12120d0321e0SJeremy L Thompson            request, pointBlock); CeedChkBackend(ierr);
12130d0321e0SJeremy L Thompson   }
12140d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
12150d0321e0SJeremy L Thompson }
12160d0321e0SJeremy L Thompson 
12170d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
12180d0321e0SJeremy L Thompson // Assemble Linear Diagonal
12190d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
12200d0321e0SJeremy L Thompson static int CeedOperatorLinearAssembleAddDiagonal_Cuda(CeedOperator op,
12210d0321e0SJeremy L Thompson     CeedVector assembled, CeedRequest *request) {
12220d0321e0SJeremy L Thompson   int ierr;
12230d0321e0SJeremy L Thompson   bool isComposite;
12240d0321e0SJeremy L Thompson   ierr = CeedOperatorIsComposite(op, &isComposite); CeedChkBackend(ierr);
12250d0321e0SJeremy L Thompson   if (isComposite) {
12260d0321e0SJeremy L Thompson     return CeedOperatorLinearAssembleAddDiagonalCompositeCore_Cuda(op, assembled,
12270d0321e0SJeremy L Thompson            request, false);
12280d0321e0SJeremy L Thompson   } else {
12290d0321e0SJeremy L Thompson     return CeedOperatorAssembleDiagonalCore_Cuda(op, assembled, request, false);
12300d0321e0SJeremy L Thompson   }
12310d0321e0SJeremy L Thompson }
12320d0321e0SJeremy L Thompson 
12330d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
12340d0321e0SJeremy L Thompson // Assemble Linear Point Block Diagonal
12350d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
12360d0321e0SJeremy L Thompson static int CeedOperatorLinearAssembleAddPointBlockDiagonal_Cuda(CeedOperator op,
12370d0321e0SJeremy L Thompson     CeedVector assembled, CeedRequest *request) {
12380d0321e0SJeremy L Thompson   int ierr;
12390d0321e0SJeremy L Thompson   bool isComposite;
12400d0321e0SJeremy L Thompson   ierr = CeedOperatorIsComposite(op, &isComposite); CeedChkBackend(ierr);
12410d0321e0SJeremy L Thompson   if (isComposite) {
12420d0321e0SJeremy L Thompson     return CeedOperatorLinearAssembleAddDiagonalCompositeCore_Cuda(op, assembled,
12430d0321e0SJeremy L Thompson            request, true);
12440d0321e0SJeremy L Thompson   } else {
12450d0321e0SJeremy L Thompson     return CeedOperatorAssembleDiagonalCore_Cuda(op, assembled, request, true);
12460d0321e0SJeremy L Thompson   }
12470d0321e0SJeremy L Thompson }
12480d0321e0SJeremy L Thompson 
12490d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
125059ad764aSnbeams // Matrix assembly kernel for low-order elements (2D thread block)
1251cc132f9aSnbeams //------------------------------------------------------------------------------
1252cc132f9aSnbeams // *INDENT-OFF*
125359ad764aSnbeams static const char *assemblykernel = QUOTE(
1254cc132f9aSnbeams extern "C" __launch_bounds__(BLOCK_SIZE)
1255cc132f9aSnbeams            __global__ void linearAssemble(const CeedScalar *B_in, const CeedScalar *B_out,
1256cc132f9aSnbeams                    const CeedScalar *__restrict__ qf_array,
1257cc132f9aSnbeams                    CeedScalar *__restrict__ values_array) {
1258cc132f9aSnbeams 
1259cc132f9aSnbeams   // This kernel assumes B_in and B_out have the same number of quadrature points and
1260cc132f9aSnbeams   // basis points.
1261cc132f9aSnbeams   // TODO: expand to more general cases
1262cc132f9aSnbeams   const int i = threadIdx.x; // The output row index of each B^TDB operation
1263cc132f9aSnbeams   const int l = threadIdx.y; // The output column index of each B^TDB operation
1264cc132f9aSnbeams 			     // such that we have (Bout^T)_ij D_jk Bin_kl = C_il
1265cc132f9aSnbeams 
1266cc132f9aSnbeams   // Strides for final output ordering, determined by the reference (interface) implementation of
1267cc132f9aSnbeams   // the symbolic assembly, slowest --> fastest: element, comp_in, comp_out, node_row, node_col
1268cc132f9aSnbeams   const CeedInt comp_out_stride = NNODES * NNODES;
1269cc132f9aSnbeams   const CeedInt comp_in_stride = comp_out_stride * NCOMP;
1270cc132f9aSnbeams   const CeedInt e_stride = comp_in_stride * NCOMP;
1271cc132f9aSnbeams   // Strides for QF array, slowest --> fastest:  emode_in, comp_in, emode_out, comp_out, elem, qpt
1272cc132f9aSnbeams   const CeedInt qe_stride = NQPTS;
1273cc132f9aSnbeams   const CeedInt qcomp_out_stride = NELEM * qe_stride;
1274cc132f9aSnbeams   const CeedInt qemode_out_stride = qcomp_out_stride * NCOMP;
1275cc132f9aSnbeams   const CeedInt qcomp_in_stride = qemode_out_stride * NUMEMODEOUT;
1276cc132f9aSnbeams   const CeedInt qemode_in_stride = qcomp_in_stride * NCOMP;
1277cc132f9aSnbeams 
1278cc132f9aSnbeams   // Loop over each element (if necessary)
1279cc132f9aSnbeams   for (CeedInt e = blockIdx.x*blockDim.z + threadIdx.z; e < NELEM;
1280cc132f9aSnbeams          e += gridDim.x*blockDim.z) {
1281cc132f9aSnbeams     for (CeedInt comp_in = 0; comp_in < NCOMP; comp_in++) {
1282cc132f9aSnbeams       for (CeedInt comp_out = 0; comp_out < NCOMP; comp_out++) {
1283cc132f9aSnbeams         CeedScalar result = 0.0;
1284cc132f9aSnbeams         CeedInt qf_index_comp = qcomp_in_stride * comp_in + qcomp_out_stride * comp_out + qe_stride * e;
1285cc132f9aSnbeams         for (CeedInt emode_in = 0; emode_in < NUMEMODEIN; emode_in++) {
1286cc132f9aSnbeams           CeedInt b_in_index = emode_in * NQPTS * NNODES;
1287cc132f9aSnbeams       	  for (CeedInt emode_out = 0; emode_out < NUMEMODEOUT; emode_out++) {
1288cc132f9aSnbeams              CeedInt b_out_index = emode_out * NQPTS * NNODES;
1289cc132f9aSnbeams              CeedInt qf_index = qf_index_comp + qemode_out_stride * emode_out + qemode_in_stride * emode_in;
1290cc132f9aSnbeams  	     // Perform the B^T D B operation for this 'chunk' of D (the qf_array)
1291cc132f9aSnbeams             for (CeedInt j = 0; j < NQPTS; j++) {
1292cc132f9aSnbeams      	      result += B_out[b_out_index + j * NNODES  + i] * qf_array[qf_index + j] * B_in[b_in_index + j * NNODES + l];
1293cc132f9aSnbeams 	    }
1294cc132f9aSnbeams 
1295cc132f9aSnbeams           }// end of emode_out
1296cc132f9aSnbeams         } // end of emode_in
1297cc132f9aSnbeams         CeedInt val_index = comp_in_stride * comp_in + comp_out_stride * comp_out + e_stride * e + NNODES * i + l;
1298cc132f9aSnbeams    	values_array[val_index] = result;
1299cc132f9aSnbeams       } // end of out component
1300cc132f9aSnbeams     } // end of in component
1301cc132f9aSnbeams   } // end of element loop
1302cc132f9aSnbeams }
1303cc132f9aSnbeams );
130459ad764aSnbeams 
130559ad764aSnbeams //------------------------------------------------------------------------------
130659ad764aSnbeams // Fallback kernel for larger orders (1D thread block)
130759ad764aSnbeams //------------------------------------------------------------------------------
130859ad764aSnbeams static const char *assemblykernelbigelem = QUOTE(
130959ad764aSnbeams extern "C" __launch_bounds__(BLOCK_SIZE)
131059ad764aSnbeams            __global__ void linearAssemble(const CeedScalar *B_in, const CeedScalar *B_out,
131159ad764aSnbeams                    const CeedScalar *__restrict__ qf_array,
131259ad764aSnbeams                    CeedScalar *__restrict__ values_array) {
131359ad764aSnbeams 
131459ad764aSnbeams   // This kernel assumes B_in and B_out have the same number of quadrature points and
131559ad764aSnbeams   // basis points.
131659ad764aSnbeams   // TODO: expand to more general cases
131759ad764aSnbeams   const int l = threadIdx.x; // The output column index of each B^TDB operation
131859ad764aSnbeams 			     // such that we have (Bout^T)_ij D_jk Bin_kl = C_il
131959ad764aSnbeams 
132059ad764aSnbeams   // Strides for final output ordering, determined by the reference (interface) implementation of
132159ad764aSnbeams   // the symbolic assembly, slowest --> fastest: element, comp_in, comp_out, node_row, node_col
132259ad764aSnbeams   const CeedInt comp_out_stride = NNODES * NNODES;
132359ad764aSnbeams   const CeedInt comp_in_stride = comp_out_stride * NCOMP;
132459ad764aSnbeams   const CeedInt e_stride = comp_in_stride * NCOMP;
132559ad764aSnbeams   // Strides for QF array, slowest --> fastest:  emode_in, comp_in, emode_out, comp_out, elem, qpt
132659ad764aSnbeams   const CeedInt qe_stride = NQPTS;
132759ad764aSnbeams   const CeedInt qcomp_out_stride = NELEM * qe_stride;
132859ad764aSnbeams   const CeedInt qemode_out_stride = qcomp_out_stride * NCOMP;
132959ad764aSnbeams   const CeedInt qcomp_in_stride = qemode_out_stride * NUMEMODEOUT;
133059ad764aSnbeams   const CeedInt qemode_in_stride = qcomp_in_stride * NCOMP;
133159ad764aSnbeams 
133259ad764aSnbeams     // Loop over each element (if necessary)
133359ad764aSnbeams   for (CeedInt e = blockIdx.x*blockDim.z + threadIdx.z; e < NELEM;
133459ad764aSnbeams          e += gridDim.x*blockDim.z) {
133559ad764aSnbeams     for (CeedInt comp_in = 0; comp_in < NCOMP; comp_in++) {
133659ad764aSnbeams       for (CeedInt comp_out = 0; comp_out < NCOMP; comp_out++) {
133759ad764aSnbeams         for (CeedInt i = 0; i < NNODES; i++) {
133859ad764aSnbeams           CeedScalar result = 0.0;
133959ad764aSnbeams           CeedInt qf_index_comp = qcomp_in_stride * comp_in + qcomp_out_stride * comp_out + qe_stride * e;
134059ad764aSnbeams           for (CeedInt emode_in = 0; emode_in < NUMEMODEIN; emode_in++) {
134159ad764aSnbeams             CeedInt b_in_index = emode_in * NQPTS * NNODES;
134259ad764aSnbeams         	  for (CeedInt emode_out = 0; emode_out < NUMEMODEOUT; emode_out++) {
134359ad764aSnbeams                CeedInt b_out_index = emode_out * NQPTS * NNODES;
134459ad764aSnbeams                CeedInt qf_index = qf_index_comp + qemode_out_stride * emode_out + qemode_in_stride * emode_in;
134559ad764aSnbeams    	     // Perform the B^T D B operation for this 'chunk' of D (the qf_array)
134659ad764aSnbeams               for (CeedInt j = 0; j < NQPTS; j++) {
134759ad764aSnbeams        	      result += B_out[b_out_index + j * NNODES  + i] * qf_array[qf_index + j] * B_in[b_in_index + j * NNODES + l];
134859ad764aSnbeams   	    }
134959ad764aSnbeams 
135059ad764aSnbeams             }// end of emode_out
135159ad764aSnbeams           } // end of emode_in
135259ad764aSnbeams           CeedInt val_index = comp_in_stride * comp_in + comp_out_stride * comp_out + e_stride * e + NNODES * i + l;
135359ad764aSnbeams      	  values_array[val_index] = result;
135459ad764aSnbeams         } // end of loop over element node index, i
135559ad764aSnbeams       } // end of out component
135659ad764aSnbeams     } // end of in component
135759ad764aSnbeams   } // end of element loop
135859ad764aSnbeams }
135959ad764aSnbeams );
1360cc132f9aSnbeams // *INDENT-ON*
1361cc132f9aSnbeams 
1362cc132f9aSnbeams //------------------------------------------------------------------------------
1363cc132f9aSnbeams // Single operator assembly setup
1364cc132f9aSnbeams //------------------------------------------------------------------------------
1365cc132f9aSnbeams static int CeedSingleOperatorAssembleSetup_Cuda(CeedOperator op) {
1366cc132f9aSnbeams   int ierr;
1367cc132f9aSnbeams   Ceed ceed;
1368cc132f9aSnbeams   ierr = CeedOperatorGetCeed(op, &ceed); CeedChkBackend(ierr);
1369cc132f9aSnbeams   CeedOperator_Cuda *impl;
1370cc132f9aSnbeams   ierr = CeedOperatorGetData(op, &impl); CeedChkBackend(ierr);
1371cc132f9aSnbeams 
1372cc132f9aSnbeams   // Get intput and output fields
1373cc132f9aSnbeams   CeedInt num_input_fields, num_output_fields;
1374cc132f9aSnbeams   CeedOperatorField *input_fields;
1375cc132f9aSnbeams   CeedOperatorField *output_fields;
1376cc132f9aSnbeams   ierr = CeedOperatorGetFields(op, &num_input_fields, &input_fields,
1377*b216396cSJeremy L Thompson                                &num_output_fields, &output_fields); CeedChkBackend(ierr);
1378cc132f9aSnbeams 
1379cc132f9aSnbeams   // Determine active input basis eval mode
1380cc132f9aSnbeams   CeedQFunction qf;
1381*b216396cSJeremy L Thompson   ierr = CeedOperatorGetQFunction(op, &qf); CeedChkBackend(ierr);
1382cc132f9aSnbeams   CeedQFunctionField *qf_fields;
1383*b216396cSJeremy L Thompson   ierr = CeedQFunctionGetFields(qf, NULL, &qf_fields, NULL, NULL);
1384*b216396cSJeremy L Thompson   CeedChkBackend(ierr);
1385cc132f9aSnbeams   // Note that the kernel will treat each dimension of a gradient action separately;
1386cc132f9aSnbeams   // i.e., when an active input has a CEED_EVAL_GRAD mode, num_emode_in will increment
1387cc132f9aSnbeams   // by dim.  However, for the purposes of loading the B matrices, it will be treated
1388cc132f9aSnbeams   // as one mode, and we will load/copy the entire gradient matrix at once, so
1389cc132f9aSnbeams   // num_B_in_mats_to_load will be incremented by 1.
1390cc132f9aSnbeams   CeedInt num_emode_in = 0, dim = 1, num_B_in_mats_to_load = 0, size_B_in = 0;
1391cc132f9aSnbeams   CeedEvalMode *eval_mode_in = NULL; //will be of size num_B_in_mats_load
1392cc132f9aSnbeams   CeedBasis basis_in = NULL;
1393*b216396cSJeremy L Thompson   CeedInt nqpts = 0, esize = 0;
1394cc132f9aSnbeams   CeedElemRestriction rstr_in = NULL;
1395cc132f9aSnbeams   for (CeedInt i=0; i<num_input_fields; i++) {
1396cc132f9aSnbeams     CeedVector vec;
1397*b216396cSJeremy L Thompson     ierr = CeedOperatorFieldGetVector(input_fields[i], &vec); CeedChkBackend(ierr);
1398cc132f9aSnbeams     if (vec == CEED_VECTOR_ACTIVE) {
1399cc132f9aSnbeams       ierr = CeedOperatorFieldGetBasis(input_fields[i], &basis_in);
1400*b216396cSJeremy L Thompson       CeedChkBackend(ierr);
1401*b216396cSJeremy L Thompson       ierr = CeedBasisGetDimension(basis_in, &dim); CeedChkBackend(ierr);
1402*b216396cSJeremy L Thompson       ierr = CeedBasisGetNumQuadraturePoints(basis_in, &nqpts); CeedChkBackend(ierr);
1403cc132f9aSnbeams       ierr = CeedOperatorFieldGetElemRestriction(input_fields[i], &rstr_in);
1404*b216396cSJeremy L Thompson       CeedChkBackend(ierr);
1405*b216396cSJeremy L Thompson       ierr = CeedElemRestrictionGetElementSize(rstr_in, &esize); CeedChkBackend(ierr);
1406cc132f9aSnbeams       CeedEvalMode eval_mode;
1407cc132f9aSnbeams       ierr = CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode);
1408*b216396cSJeremy L Thompson       CeedChkBackend(ierr);
1409cc132f9aSnbeams       if (eval_mode != CEED_EVAL_NONE) {
1410*b216396cSJeremy L Thompson         ierr = CeedRealloc(num_B_in_mats_to_load + 1, &eval_mode_in);
1411*b216396cSJeremy L Thompson         CeedChkBackend(ierr);
1412cc132f9aSnbeams         eval_mode_in[num_B_in_mats_to_load] = eval_mode;
1413cc132f9aSnbeams         num_B_in_mats_to_load += 1;
1414cc132f9aSnbeams         if (eval_mode == CEED_EVAL_GRAD) {
1415cc132f9aSnbeams           num_emode_in += dim;
1416cc132f9aSnbeams           size_B_in += dim * esize * nqpts;
1417cc132f9aSnbeams         } else {
1418cc132f9aSnbeams           num_emode_in +=1;
1419cc132f9aSnbeams           size_B_in += esize * nqpts;
1420cc132f9aSnbeams         }
1421cc132f9aSnbeams       }
1422cc132f9aSnbeams     }
1423cc132f9aSnbeams   }
1424cc132f9aSnbeams 
1425cc132f9aSnbeams   // Determine active output basis; basis_out and rstr_out only used if same as input, TODO
1426*b216396cSJeremy L Thompson   ierr = CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qf_fields);
1427*b216396cSJeremy L Thompson   CeedChkBackend(ierr);
1428cc132f9aSnbeams   CeedInt num_emode_out = 0, num_B_out_mats_to_load = 0, size_B_out = 0;
1429cc132f9aSnbeams   CeedEvalMode *eval_mode_out = NULL;
1430cc132f9aSnbeams   CeedBasis basis_out = NULL;
1431cc132f9aSnbeams   CeedElemRestriction rstr_out = NULL;
1432cc132f9aSnbeams   for (CeedInt i=0; i<num_output_fields; i++) {
1433cc132f9aSnbeams     CeedVector vec;
1434*b216396cSJeremy L Thompson     ierr = CeedOperatorFieldGetVector(output_fields[i], &vec); CeedChkBackend(ierr);
1435cc132f9aSnbeams     if (vec == CEED_VECTOR_ACTIVE) {
1436cc132f9aSnbeams       ierr = CeedOperatorFieldGetBasis(output_fields[i], &basis_out);
1437*b216396cSJeremy L Thompson       CeedChkBackend(ierr);
1438cc132f9aSnbeams       ierr = CeedOperatorFieldGetElemRestriction(output_fields[i], &rstr_out);
1439*b216396cSJeremy L Thompson       CeedChkBackend(ierr);
1440cc132f9aSnbeams       if (rstr_out && rstr_out != rstr_in)
1441cc132f9aSnbeams         // LCOV_EXCL_START
1442cc132f9aSnbeams         return CeedError(ceed, CEED_ERROR_BACKEND,
1443cc132f9aSnbeams                          "Multi-field non-composite operator assembly not supported");
1444cc132f9aSnbeams       // LCOV_EXCL_STOP
1445cc132f9aSnbeams       CeedEvalMode eval_mode;
1446cc132f9aSnbeams       ierr = CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode);
1447*b216396cSJeremy L Thompson       CeedChkBackend(ierr);
1448cc132f9aSnbeams       if (eval_mode != CEED_EVAL_NONE) {
1449*b216396cSJeremy L Thompson         ierr = CeedRealloc(num_B_out_mats_to_load + 1, &eval_mode_out);
1450*b216396cSJeremy L Thompson         CeedChkBackend(ierr);
1451cc132f9aSnbeams         eval_mode_out[num_B_out_mats_to_load] = eval_mode;
1452cc132f9aSnbeams         num_B_out_mats_to_load += 1;
1453cc132f9aSnbeams         if (eval_mode == CEED_EVAL_GRAD) {
1454cc132f9aSnbeams           num_emode_out += dim;
1455cc132f9aSnbeams           size_B_out += dim * esize * nqpts;
1456cc132f9aSnbeams         } else {
1457cc132f9aSnbeams           num_emode_out +=1;
1458cc132f9aSnbeams           size_B_out += esize * nqpts;
1459cc132f9aSnbeams         }
1460cc132f9aSnbeams       }
1461cc132f9aSnbeams     }
1462cc132f9aSnbeams   }
1463cc132f9aSnbeams 
1464cc132f9aSnbeams   if (num_emode_in == 0 || num_emode_out == 0)
1465cc132f9aSnbeams     // LCOV_EXCL_START
1466cc132f9aSnbeams     return CeedError(ceed, CEED_ERROR_UNSUPPORTED,
1467cc132f9aSnbeams                      "Cannot assemble operator without inputs/outputs");
1468cc132f9aSnbeams   // LCOV_EXCL_STOP
1469cc132f9aSnbeams 
1470cc132f9aSnbeams   CeedInt nelem, ncomp;
1471*b216396cSJeremy L Thompson   ierr = CeedElemRestrictionGetNumElements(rstr_in, &nelem); CeedChkBackend(ierr);
1472*b216396cSJeremy L Thompson   ierr = CeedElemRestrictionGetNumComponents(rstr_in, &ncomp);
1473*b216396cSJeremy L Thompson   CeedChkBackend(ierr);
1474cc132f9aSnbeams 
1475cc132f9aSnbeams   ierr = CeedCalloc(1, &impl->asmb); CeedChkBackend(ierr);
1476cc132f9aSnbeams   CeedOperatorAssemble_Cuda *asmb = impl->asmb;
1477cc132f9aSnbeams   asmb->nelem = nelem;
1478cc132f9aSnbeams 
1479cc132f9aSnbeams   // Compile kernels
1480cc132f9aSnbeams   int elemsPerBlock = 1;
1481cc132f9aSnbeams   asmb->elemsPerBlock = elemsPerBlock;
148259ad764aSnbeams   CeedInt block_size = esize * esize * elemsPerBlock;
1483cc132f9aSnbeams   Ceed_Cuda *cuda_data;
1484cc132f9aSnbeams   ierr = CeedGetData(ceed, &cuda_data); CeedChkBackend(ierr);
148559ad764aSnbeams   if (block_size > cuda_data->device_prop.maxThreadsPerBlock) {
148659ad764aSnbeams     // Use fallback kernel with 1D threadblock
148759ad764aSnbeams     block_size = esize * elemsPerBlock;
148859ad764aSnbeams     asmb->block_size_x = esize;
148959ad764aSnbeams     asmb->block_size_y = 1;
149059ad764aSnbeams     ierr = CeedCompileCuda(ceed, assemblykernelbigelem, &asmb->module, 7,
1491cc132f9aSnbeams                            "NELEM", nelem,
1492cc132f9aSnbeams                            "NUMEMODEIN", num_emode_in,
1493cc132f9aSnbeams                            "NUMEMODEOUT", num_emode_out,
1494cc132f9aSnbeams                            "NQPTS", nqpts,
1495cc132f9aSnbeams                            "NNODES", esize,
149659ad764aSnbeams                            "BLOCK_SIZE", block_size,
1497cc132f9aSnbeams                            "NCOMP", ncomp
1498cc132f9aSnbeams                           ); CeedChk_Cu(ceed, ierr);
149959ad764aSnbeams   } else {  // Use kernel with 2D threadblock
150059ad764aSnbeams     asmb->block_size_x = esize;
150159ad764aSnbeams     asmb->block_size_y = esize;
150259ad764aSnbeams     ierr = CeedCompileCuda(ceed, assemblykernel, &asmb->module, 7,
150359ad764aSnbeams                            "NELEM", nelem,
150459ad764aSnbeams                            "NUMEMODEIN", num_emode_in,
150559ad764aSnbeams                            "NUMEMODEOUT", num_emode_out,
150659ad764aSnbeams                            "NQPTS", nqpts,
150759ad764aSnbeams                            "NNODES", esize,
150859ad764aSnbeams                            "BLOCK_SIZE", block_size,
150959ad764aSnbeams                            "NCOMP", ncomp
151059ad764aSnbeams                           ); CeedChk_Cu(ceed, ierr);
151159ad764aSnbeams   }
1512cc132f9aSnbeams   ierr = CeedGetKernelCuda(ceed, asmb->module, "linearAssemble",
1513cc132f9aSnbeams                            &asmb->linearAssemble); CeedChk_Cu(ceed, ierr);
1514cc132f9aSnbeams 
1515cc132f9aSnbeams   // Build 'full' B matrices (not 1D arrays used for tensor-product matrices)
1516cc132f9aSnbeams   const CeedScalar *interp_in, *grad_in;
1517cc132f9aSnbeams   ierr = CeedBasisGetInterp(basis_in, &interp_in); CeedChkBackend(ierr);
1518cc132f9aSnbeams   ierr = CeedBasisGetGrad(basis_in, &grad_in); CeedChkBackend(ierr);
1519cc132f9aSnbeams 
1520cc132f9aSnbeams   // Load into B_in, in order that they will be used in eval_mode
1521cc132f9aSnbeams   const CeedInt inBytes = size_B_in * sizeof(CeedScalar);
1522cc132f9aSnbeams   CeedInt mat_start = 0;
1523cc132f9aSnbeams   ierr = cudaMalloc((void **) &asmb->d_B_in, inBytes); CeedChk_Cu(ceed, ierr);
1524cc132f9aSnbeams   for (int i = 0; i < num_B_in_mats_to_load; i++) {
1525cc132f9aSnbeams     CeedEvalMode eval_mode = eval_mode_in[i];
1526cc132f9aSnbeams     if (eval_mode == CEED_EVAL_INTERP) {
1527cc132f9aSnbeams       ierr = cudaMemcpy(&asmb->d_B_in[mat_start], interp_in,
1528cc132f9aSnbeams                         esize * nqpts * sizeof(CeedScalar),
1529cc132f9aSnbeams                         cudaMemcpyHostToDevice); CeedChk_Cu(ceed, ierr);
1530cc132f9aSnbeams       mat_start += esize * nqpts;
1531cc132f9aSnbeams     } else if (eval_mode == CEED_EVAL_GRAD) {
1532cc132f9aSnbeams       ierr = cudaMemcpy(asmb->d_B_in, grad_in,
1533cc132f9aSnbeams                         dim * esize * nqpts * sizeof(CeedScalar),
1534cc132f9aSnbeams                         cudaMemcpyHostToDevice); CeedChk_Cu(ceed, ierr);
1535cc132f9aSnbeams       mat_start += dim * esize * nqpts;
1536cc132f9aSnbeams     }
1537cc132f9aSnbeams   }
1538cc132f9aSnbeams 
1539cc132f9aSnbeams   const CeedScalar *interp_out, *grad_out;
1540cc132f9aSnbeams   // Note that this function currently assumes 1 basis, so this should always be true
1541cc132f9aSnbeams   // for now
1542cc132f9aSnbeams   if (basis_out == basis_in) {
1543cc132f9aSnbeams     interp_out = interp_in;
1544cc132f9aSnbeams     grad_out = grad_in;
1545cc132f9aSnbeams   } else {
1546cc132f9aSnbeams     ierr = CeedBasisGetInterp(basis_out, &interp_out); CeedChkBackend(ierr);
1547cc132f9aSnbeams     ierr = CeedBasisGetGrad(basis_out, &grad_out); CeedChkBackend(ierr);
1548cc132f9aSnbeams   }
1549cc132f9aSnbeams 
1550cc132f9aSnbeams   // Load into B_out, in order that they will be used in eval_mode
1551cc132f9aSnbeams   const CeedInt outBytes = size_B_out * sizeof(CeedScalar);
1552cc132f9aSnbeams   mat_start = 0;
1553cc132f9aSnbeams   ierr = cudaMalloc((void **) &asmb->d_B_out, outBytes); CeedChk_Cu(ceed, ierr);
1554cc132f9aSnbeams   for (int i = 0; i < num_B_out_mats_to_load; i++) {
1555cc132f9aSnbeams     CeedEvalMode eval_mode = eval_mode_out[i];
1556cc132f9aSnbeams     if (eval_mode == CEED_EVAL_INTERP) {
1557cc132f9aSnbeams       ierr = cudaMemcpy(&asmb->d_B_out[mat_start], interp_out,
1558cc132f9aSnbeams                         esize * nqpts * sizeof(CeedScalar),
1559cc132f9aSnbeams                         cudaMemcpyHostToDevice); CeedChk_Cu(ceed, ierr);
1560cc132f9aSnbeams       mat_start += esize * nqpts;
1561cc132f9aSnbeams     } else if (eval_mode == CEED_EVAL_GRAD) {
1562cc132f9aSnbeams       ierr = cudaMemcpy(&asmb->d_B_out[mat_start], grad_out,
1563cc132f9aSnbeams                         dim * esize * nqpts * sizeof(CeedScalar),
1564cc132f9aSnbeams                         cudaMemcpyHostToDevice); CeedChk_Cu(ceed, ierr);
1565cc132f9aSnbeams       mat_start += dim * esize * nqpts;
1566cc132f9aSnbeams     }
1567cc132f9aSnbeams   }
1568cc132f9aSnbeams   return CEED_ERROR_SUCCESS;
1569cc132f9aSnbeams }
1570cc132f9aSnbeams 
1571cc132f9aSnbeams //------------------------------------------------------------------------------
1572cc132f9aSnbeams // Single operator assembly
1573cc132f9aSnbeams //------------------------------------------------------------------------------
1574cc132f9aSnbeams static int CeedSingleOperatorAssemble_Cuda(CeedOperator op, CeedInt offset,
1575cc132f9aSnbeams     CeedVector values) {
1576cc132f9aSnbeams 
1577cc132f9aSnbeams   int ierr;
1578cc132f9aSnbeams   Ceed ceed;
1579cc132f9aSnbeams   ierr = CeedOperatorGetCeed(op, &ceed); CeedChkBackend(ierr);
1580cc132f9aSnbeams   CeedOperator_Cuda *impl;
1581cc132f9aSnbeams   ierr = CeedOperatorGetData(op, &impl); CeedChkBackend(ierr);
1582cc132f9aSnbeams 
1583cc132f9aSnbeams   // Setup
1584cc132f9aSnbeams   if (!impl->asmb) {
1585cc132f9aSnbeams     ierr = CeedSingleOperatorAssembleSetup_Cuda(op);
1586cc132f9aSnbeams     CeedChkBackend(ierr);
1587*b216396cSJeremy L Thompson     assert(impl->asmb != NULL);
1588cc132f9aSnbeams   }
1589cc132f9aSnbeams 
1590cc132f9aSnbeams   // Assemble QFunction
1591cc132f9aSnbeams   CeedVector assembled_qf;
1592cc132f9aSnbeams   CeedElemRestriction rstr_q;
1593cc132f9aSnbeams   ierr = CeedOperatorLinearAssembleQFunctionBuildOrUpdate(
1594*b216396cSJeremy L Thompson            op, &assembled_qf, &rstr_q, CEED_REQUEST_IMMEDIATE); CeedChkBackend(ierr);
1595cc132f9aSnbeams   ierr = CeedElemRestrictionDestroy(&rstr_q); CeedChkBackend(ierr);
1596cc132f9aSnbeams   CeedScalar *values_array;
1597cc132f9aSnbeams   ierr = CeedVectorGetArrayWrite(values, CEED_MEM_DEVICE, &values_array);
1598cc132f9aSnbeams   CeedChkBackend(ierr);
1599cc132f9aSnbeams   values_array += offset;
1600cc132f9aSnbeams   const CeedScalar *qf_array;
1601cc132f9aSnbeams   ierr = CeedVectorGetArrayRead(assembled_qf, CEED_MEM_DEVICE, &qf_array);
1602cc132f9aSnbeams   CeedChkBackend(ierr);
1603cc132f9aSnbeams 
1604cc132f9aSnbeams   // Compute B^T D B
1605cc132f9aSnbeams   const CeedInt nelem = impl->asmb->nelem;
1606cc132f9aSnbeams   const CeedInt elemsPerBlock = impl->asmb->elemsPerBlock;
1607cc132f9aSnbeams   const CeedInt grid = nelem/elemsPerBlock+((
1608cc132f9aSnbeams                          nelem/elemsPerBlock*elemsPerBlock<nelem)?1:0);
1609cc132f9aSnbeams   void *args[] = {&impl->asmb->d_B_in, &impl->asmb->d_B_out,
1610cc132f9aSnbeams                   &qf_array, &values_array
1611cc132f9aSnbeams                  };
1612cc132f9aSnbeams   ierr = CeedRunKernelDimCuda(ceed, impl->asmb->linearAssemble, grid,
161359ad764aSnbeams                               impl->asmb->block_size_x, impl->asmb->block_size_y,
161459ad764aSnbeams                               elemsPerBlock, args);
1615cc132f9aSnbeams   CeedChkBackend(ierr);
1616cc132f9aSnbeams 
1617cc132f9aSnbeams 
1618cc132f9aSnbeams   // Restore arrays
1619cc132f9aSnbeams   ierr = CeedVectorRestoreArray(values, &values_array); CeedChkBackend(ierr);
1620cc132f9aSnbeams   ierr = CeedVectorRestoreArrayRead(assembled_qf, &qf_array);
1621cc132f9aSnbeams   CeedChkBackend(ierr);
1622cc132f9aSnbeams 
1623cc132f9aSnbeams   // Cleanup
1624cc132f9aSnbeams   ierr = CeedVectorDestroy(&assembled_qf); CeedChkBackend(ierr);
1625cc132f9aSnbeams 
1626cc132f9aSnbeams   return CEED_ERROR_SUCCESS;
1627cc132f9aSnbeams }
1628cc132f9aSnbeams 
1629cc132f9aSnbeams //------------------------------------------------------------------------------
1630cc132f9aSnbeams // Assemble matrix data for COO matrix of assembled operator.
1631cc132f9aSnbeams // The sparsity pattern is set by CeedOperatorLinearAssembleSymbolic.
1632cc132f9aSnbeams //
1633cc132f9aSnbeams // Note that this (and other assembly routines) currently assume only one
1634cc132f9aSnbeams // active input restriction/basis per operator (could have multiple basis eval
1635cc132f9aSnbeams // modes).
1636cc132f9aSnbeams // TODO: allow multiple active input restrictions/basis objects
1637cc132f9aSnbeams //------------------------------------------------------------------------------
1638cc132f9aSnbeams int CeedOperatorLinearAssemble_Cuda(CeedOperator op, CeedVector values) {
1639cc132f9aSnbeams 
1640cc132f9aSnbeams   // As done in the default implementation, loop through suboperators
1641cc132f9aSnbeams   // for composite operators, or call single operator assembly otherwise
1642cc132f9aSnbeams   bool is_composite;
1643cc132f9aSnbeams   CeedInt ierr;
1644*b216396cSJeremy L Thompson   ierr = CeedOperatorIsComposite(op, &is_composite); CeedChkBackend(ierr);
1645cc132f9aSnbeams 
1646cc132f9aSnbeams   CeedElemRestriction rstr;
1647cc132f9aSnbeams   CeedInt num_elem, elem_size, num_comp;
1648cc132f9aSnbeams 
1649cc132f9aSnbeams   CeedInt offset = 0;
1650cc132f9aSnbeams   if (is_composite) {
1651cc132f9aSnbeams     CeedInt num_suboperators;
1652*b216396cSJeremy L Thompson     ierr = CeedOperatorGetNumSub(op, &num_suboperators); CeedChkBackend(ierr);
1653cc132f9aSnbeams     CeedOperator *sub_operators;
1654*b216396cSJeremy L Thompson     ierr = CeedOperatorGetSubList(op, &sub_operators); CeedChkBackend(ierr);
1655cc132f9aSnbeams     for (int k = 0; k < num_suboperators; ++k) {
1656cc132f9aSnbeams       ierr = CeedSingleOperatorAssemble_Cuda(sub_operators[k], offset, values);
1657*b216396cSJeremy L Thompson       CeedChkBackend(ierr);
1658cc132f9aSnbeams       ierr = CeedOperatorGetActiveElemRestriction(sub_operators[k], &rstr);
1659*b216396cSJeremy L Thompson       CeedChkBackend(ierr);
1660*b216396cSJeremy L Thompson       ierr = CeedElemRestrictionGetNumElements(rstr, &num_elem); CeedChkBackend(ierr);
1661*b216396cSJeremy L Thompson       ierr = CeedElemRestrictionGetElementSize(rstr, &elem_size);
1662*b216396cSJeremy L Thompson       CeedChkBackend(ierr);
1663*b216396cSJeremy L Thompson       ierr = CeedElemRestrictionGetNumComponents(rstr, &num_comp);
1664*b216396cSJeremy L Thompson       CeedChkBackend(ierr);
1665cc132f9aSnbeams       offset += elem_size*num_comp * elem_size*num_comp * num_elem;
1666cc132f9aSnbeams     }
1667cc132f9aSnbeams   } else {
1668*b216396cSJeremy L Thompson     ierr = CeedSingleOperatorAssemble_Cuda(op, offset, values);
1669*b216396cSJeremy L Thompson     CeedChkBackend(ierr);
1670cc132f9aSnbeams   }
1671cc132f9aSnbeams 
1672cc132f9aSnbeams   return CEED_ERROR_SUCCESS;
1673cc132f9aSnbeams }
1674cc132f9aSnbeams //------------------------------------------------------------------------------
16750d0321e0SJeremy L Thompson // Create operator
16760d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
16770d0321e0SJeremy L Thompson int CeedOperatorCreate_Cuda(CeedOperator op) {
16780d0321e0SJeremy L Thompson   int ierr;
16790d0321e0SJeremy L Thompson   Ceed ceed;
16800d0321e0SJeremy L Thompson   ierr = CeedOperatorGetCeed(op, &ceed); CeedChkBackend(ierr);
16810d0321e0SJeremy L Thompson   CeedOperator_Cuda *impl;
16820d0321e0SJeremy L Thompson 
16830d0321e0SJeremy L Thompson   ierr = CeedCalloc(1, &impl); CeedChkBackend(ierr);
16840d0321e0SJeremy L Thompson   ierr = CeedOperatorSetData(op, impl); CeedChkBackend(ierr);
16850d0321e0SJeremy L Thompson 
16860d0321e0SJeremy L Thompson   ierr = CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleQFunction",
16870d0321e0SJeremy L Thompson                                 CeedOperatorLinearAssembleQFunction_Cuda);
16880d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
16890d0321e0SJeremy L Thompson   ierr = CeedSetBackendFunction(ceed, "Operator", op,
16900d0321e0SJeremy L Thompson                                 "LinearAssembleQFunctionUpdate",
16910d0321e0SJeremy L Thompson                                 CeedOperatorLinearAssembleQFunctionUpdate_Cuda);
16920d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
16930d0321e0SJeremy L Thompson   ierr = CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleAddDiagonal",
16940d0321e0SJeremy L Thompson                                 CeedOperatorLinearAssembleAddDiagonal_Cuda);
16950d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
16960d0321e0SJeremy L Thompson   ierr = CeedSetBackendFunction(ceed, "Operator", op,
16970d0321e0SJeremy L Thompson                                 "LinearAssembleAddPointBlockDiagonal",
16980d0321e0SJeremy L Thompson                                 CeedOperatorLinearAssembleAddPointBlockDiagonal_Cuda);
16990d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
1700cc132f9aSnbeams   ierr = CeedSetBackendFunction(ceed, "Operator", op,
1701cc132f9aSnbeams                                 "LinearAssemble", CeedOperatorLinearAssemble_Cuda);
1702cc132f9aSnbeams   CeedChkBackend(ierr);
17030d0321e0SJeremy L Thompson   ierr = CeedSetBackendFunction(ceed, "Operator", op, "ApplyAdd",
17040d0321e0SJeremy L Thompson                                 CeedOperatorApplyAdd_Cuda); CeedChkBackend(ierr);
17050d0321e0SJeremy L Thompson   ierr = CeedSetBackendFunction(ceed, "Operator", op, "Destroy",
17060d0321e0SJeremy L Thompson                                 CeedOperatorDestroy_Cuda); CeedChkBackend(ierr);
17070d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
17080d0321e0SJeremy L Thompson }
17090d0321e0SJeremy L Thompson 
17100d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
17110d0321e0SJeremy L Thompson // Composite Operator Create
17120d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
17130d0321e0SJeremy L Thompson int CeedCompositeOperatorCreate_Cuda(CeedOperator op) {
17140d0321e0SJeremy L Thompson   int ierr;
17150d0321e0SJeremy L Thompson   Ceed ceed;
17160d0321e0SJeremy L Thompson   ierr = CeedOperatorGetCeed(op, &ceed); CeedChkBackend(ierr);
17170d0321e0SJeremy L Thompson 
17180d0321e0SJeremy L Thompson   ierr = CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleAddDiagonal",
17190d0321e0SJeremy L Thompson                                 CeedOperatorLinearAssembleAddDiagonal_Cuda);
17200d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
17210d0321e0SJeremy L Thompson   ierr = CeedSetBackendFunction(ceed, "Operator", op,
17220d0321e0SJeremy L Thompson                                 "LinearAssembleAddPointBlockDiagonal",
17230d0321e0SJeremy L Thompson                                 CeedOperatorLinearAssembleAddPointBlockDiagonal_Cuda);
17240d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
1725cc132f9aSnbeams   ierr = CeedSetBackendFunction(ceed, "Operator", op,
1726cc132f9aSnbeams                                 "LinearAssemble", CeedOperatorLinearAssemble_Cuda);
1727cc132f9aSnbeams   CeedChkBackend(ierr);
17280d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
17290d0321e0SJeremy L Thompson }
17300d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
1731