xref: /libCEED/backends/cuda-ref/ceed-cuda-ref-operator.c (revision 3d8e882215d238700cdceb37404f76ca7fa24eaa)
1*3d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2*3d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
30d0321e0SJeremy L Thompson //
4*3d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
50d0321e0SJeremy L Thompson //
6*3d8e8822SJeremy 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 
690d0321e0SJeremy L Thompson   ierr = CeedFree(&impl); CeedChkBackend(ierr);
700d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
710d0321e0SJeremy L Thompson }
720d0321e0SJeremy L Thompson 
730d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
740d0321e0SJeremy L Thompson // Setup infields or outfields
750d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
760d0321e0SJeremy L Thompson static int CeedOperatorSetupFields_Cuda(CeedQFunction qf, CeedOperator op,
770d0321e0SJeremy L Thompson                                         bool isinput, CeedVector *evecs,
780d0321e0SJeremy L Thompson                                         CeedVector *qvecs, CeedInt starte,
790d0321e0SJeremy L Thompson                                         CeedInt numfields, CeedInt Q,
800d0321e0SJeremy L Thompson                                         CeedInt numelements) {
810d0321e0SJeremy L Thompson   CeedInt dim, ierr, size;
820d0321e0SJeremy L Thompson   Ceed ceed;
830d0321e0SJeremy L Thompson   ierr = CeedOperatorGetCeed(op, &ceed); CeedChkBackend(ierr);
840d0321e0SJeremy L Thompson   CeedBasis basis;
850d0321e0SJeremy L Thompson   CeedElemRestriction Erestrict;
860d0321e0SJeremy L Thompson   CeedOperatorField *opfields;
870d0321e0SJeremy L Thompson   CeedQFunctionField *qffields;
880d0321e0SJeremy L Thompson   CeedVector fieldvec;
890d0321e0SJeremy L Thompson   bool strided;
900d0321e0SJeremy L Thompson   bool skiprestrict;
910d0321e0SJeremy L Thompson 
920d0321e0SJeremy L Thompson   if (isinput) {
930d0321e0SJeremy L Thompson     ierr = CeedOperatorGetFields(op, NULL, &opfields, NULL, NULL);
940d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
950d0321e0SJeremy L Thompson     ierr = CeedQFunctionGetFields(qf, NULL, &qffields, NULL, NULL);
960d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
970d0321e0SJeremy L Thompson   } else {
980d0321e0SJeremy L Thompson     ierr = CeedOperatorGetFields(op, NULL, NULL, NULL, &opfields);
990d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
1000d0321e0SJeremy L Thompson     ierr = CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qffields);
1010d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
1020d0321e0SJeremy L Thompson   }
1030d0321e0SJeremy L Thompson 
1040d0321e0SJeremy L Thompson   // Loop over fields
1050d0321e0SJeremy L Thompson   for (CeedInt i = 0; i < numfields; i++) {
1060d0321e0SJeremy L Thompson     CeedEvalMode emode;
1070d0321e0SJeremy L Thompson     ierr = CeedQFunctionFieldGetEvalMode(qffields[i], &emode); CeedChkBackend(ierr);
1080d0321e0SJeremy L Thompson 
1090d0321e0SJeremy L Thompson     strided = false;
1100d0321e0SJeremy L Thompson     skiprestrict = false;
1110d0321e0SJeremy L Thompson     if (emode != CEED_EVAL_WEIGHT) {
1120d0321e0SJeremy L Thompson       ierr = CeedOperatorFieldGetElemRestriction(opfields[i], &Erestrict);
1130d0321e0SJeremy L Thompson       CeedChkBackend(ierr);
1140d0321e0SJeremy L Thompson 
1150d0321e0SJeremy L Thompson       // Check whether this field can skip the element restriction:
1160d0321e0SJeremy L Thompson       // must be passive input, with emode NONE, and have a strided restriction with
1170d0321e0SJeremy L Thompson       // CEED_STRIDES_BACKEND.
1180d0321e0SJeremy L Thompson 
1190d0321e0SJeremy L Thompson       // First, check whether the field is input or output:
1200d0321e0SJeremy L Thompson       if (isinput) {
1210d0321e0SJeremy L Thompson         // Check for passive input:
1220d0321e0SJeremy L Thompson         ierr = CeedOperatorFieldGetVector(opfields[i], &fieldvec); CeedChkBackend(ierr);
1230d0321e0SJeremy L Thompson         if (fieldvec != CEED_VECTOR_ACTIVE) {
1240d0321e0SJeremy L Thompson           // Check emode
1250d0321e0SJeremy L Thompson           if (emode == CEED_EVAL_NONE) {
1260d0321e0SJeremy L Thompson             // Check for strided restriction
1270d0321e0SJeremy L Thompson             ierr = CeedElemRestrictionIsStrided(Erestrict, &strided);
1280d0321e0SJeremy L Thompson             CeedChkBackend(ierr);
1290d0321e0SJeremy L Thompson             if (strided) {
1300d0321e0SJeremy L Thompson               // Check if vector is already in preferred backend ordering
1310d0321e0SJeremy L Thompson               ierr = CeedElemRestrictionHasBackendStrides(Erestrict,
1320d0321e0SJeremy L Thompson                      &skiprestrict); CeedChkBackend(ierr);
1330d0321e0SJeremy L Thompson             }
1340d0321e0SJeremy L Thompson           }
1350d0321e0SJeremy L Thompson         }
1360d0321e0SJeremy L Thompson       }
1370d0321e0SJeremy L Thompson       if (skiprestrict) {
1380d0321e0SJeremy L Thompson         // We do not need an E-Vector, but will use the input field vector's data
1390d0321e0SJeremy L Thompson         // directly in the operator application.
1400d0321e0SJeremy L Thompson         evecs[i + starte] = NULL;
1410d0321e0SJeremy L Thompson       } else {
1420d0321e0SJeremy L Thompson         ierr = CeedElemRestrictionCreateVector(Erestrict, NULL,
1430d0321e0SJeremy L Thompson                                                &evecs[i + starte]);
1440d0321e0SJeremy L Thompson         CeedChkBackend(ierr);
1450d0321e0SJeremy L Thompson       }
1460d0321e0SJeremy L Thompson     }
1470d0321e0SJeremy L Thompson 
1480d0321e0SJeremy L Thompson     switch (emode) {
1490d0321e0SJeremy L Thompson     case CEED_EVAL_NONE:
1500d0321e0SJeremy L Thompson       ierr = CeedQFunctionFieldGetSize(qffields[i], &size); CeedChkBackend(ierr);
1510d0321e0SJeremy L Thompson       ierr = CeedVectorCreate(ceed, numelements * Q * size, &qvecs[i]);
1520d0321e0SJeremy L Thompson       CeedChkBackend(ierr);
1530d0321e0SJeremy L Thompson       break;
1540d0321e0SJeremy L Thompson     case CEED_EVAL_INTERP:
1550d0321e0SJeremy L Thompson       ierr = CeedQFunctionFieldGetSize(qffields[i], &size); CeedChkBackend(ierr);
1560d0321e0SJeremy L Thompson       ierr = CeedVectorCreate(ceed, numelements * Q * size, &qvecs[i]);
1570d0321e0SJeremy L Thompson       CeedChkBackend(ierr);
1580d0321e0SJeremy L Thompson       break;
1590d0321e0SJeremy L Thompson     case CEED_EVAL_GRAD:
1600d0321e0SJeremy L Thompson       ierr = CeedOperatorFieldGetBasis(opfields[i], &basis); CeedChkBackend(ierr);
1610d0321e0SJeremy L Thompson       ierr = CeedQFunctionFieldGetSize(qffields[i], &size); CeedChkBackend(ierr);
1620d0321e0SJeremy L Thompson       ierr = CeedBasisGetDimension(basis, &dim); CeedChkBackend(ierr);
1630d0321e0SJeremy L Thompson       ierr = CeedVectorCreate(ceed, numelements * Q * size, &qvecs[i]);
1640d0321e0SJeremy L Thompson       CeedChkBackend(ierr);
1650d0321e0SJeremy L Thompson       break;
1660d0321e0SJeremy L Thompson     case CEED_EVAL_WEIGHT: // Only on input fields
1670d0321e0SJeremy L Thompson       ierr = CeedOperatorFieldGetBasis(opfields[i], &basis); CeedChkBackend(ierr);
1680d0321e0SJeremy L Thompson       ierr = CeedVectorCreate(ceed, numelements * Q, &qvecs[i]); CeedChkBackend(ierr);
1690d0321e0SJeremy L Thompson       ierr = CeedBasisApply(basis, numelements, CEED_NOTRANSPOSE,
1700d0321e0SJeremy L Thompson                             CEED_EVAL_WEIGHT, NULL, qvecs[i]); CeedChkBackend(ierr);
1710d0321e0SJeremy L Thompson       break;
1720d0321e0SJeremy L Thompson     case CEED_EVAL_DIV:
1730d0321e0SJeremy L Thompson       break; // TODO: Not implemented
1740d0321e0SJeremy L Thompson     case CEED_EVAL_CURL:
1750d0321e0SJeremy L Thompson       break; // TODO: Not implemented
1760d0321e0SJeremy L Thompson     }
1770d0321e0SJeremy L Thompson   }
1780d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1790d0321e0SJeremy L Thompson }
1800d0321e0SJeremy L Thompson 
1810d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
1820d0321e0SJeremy L Thompson // CeedOperator needs to connect all the named fields (be they active or passive)
1830d0321e0SJeremy L Thompson //   to the named inputs and outputs of its CeedQFunction.
1840d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
1850d0321e0SJeremy L Thompson static int CeedOperatorSetup_Cuda(CeedOperator op) {
1860d0321e0SJeremy L Thompson   int ierr;
1870d0321e0SJeremy L Thompson   bool setupdone;
1880d0321e0SJeremy L Thompson   ierr = CeedOperatorIsSetupDone(op, &setupdone); CeedChkBackend(ierr);
1890d0321e0SJeremy L Thompson   if (setupdone)
1900d0321e0SJeremy L Thompson     return CEED_ERROR_SUCCESS;
1910d0321e0SJeremy L Thompson   Ceed ceed;
1920d0321e0SJeremy L Thompson   ierr = CeedOperatorGetCeed(op, &ceed); CeedChkBackend(ierr);
1930d0321e0SJeremy L Thompson   CeedOperator_Cuda *impl;
1940d0321e0SJeremy L Thompson   ierr = CeedOperatorGetData(op, &impl); CeedChkBackend(ierr);
1950d0321e0SJeremy L Thompson   CeedQFunction qf;
1960d0321e0SJeremy L Thompson   ierr = CeedOperatorGetQFunction(op, &qf); CeedChkBackend(ierr);
1970d0321e0SJeremy L Thompson   CeedInt Q, numelements, numinputfields, numoutputfields;
1980d0321e0SJeremy L Thompson   ierr = CeedOperatorGetNumQuadraturePoints(op, &Q); CeedChkBackend(ierr);
1990d0321e0SJeremy L Thompson   ierr = CeedOperatorGetNumElements(op, &numelements); CeedChkBackend(ierr);
2000d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
2010d0321e0SJeremy L Thompson   CeedOperatorField *opinputfields, *opoutputfields;
2020d0321e0SJeremy L Thompson   ierr = CeedOperatorGetFields(op, &numinputfields, &opinputfields,
2030d0321e0SJeremy L Thompson                                &numoutputfields, &opoutputfields);
2040d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
2050d0321e0SJeremy L Thompson   CeedQFunctionField *qfinputfields, *qfoutputfields;
2060d0321e0SJeremy L Thompson   ierr = CeedQFunctionGetFields(qf, NULL, &qfinputfields, NULL, &qfoutputfields);
2070d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
2080d0321e0SJeremy L Thompson 
2090d0321e0SJeremy L Thompson   // Allocate
2100d0321e0SJeremy L Thompson   ierr = CeedCalloc(numinputfields + numoutputfields, &impl->evecs);
2110d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
2120d0321e0SJeremy L Thompson 
2130d0321e0SJeremy L Thompson   ierr = CeedCalloc(CEED_FIELD_MAX, &impl->qvecsin); CeedChkBackend(ierr);
2140d0321e0SJeremy L Thompson   ierr = CeedCalloc(CEED_FIELD_MAX, &impl->qvecsout); CeedChkBackend(ierr);
2150d0321e0SJeremy L Thompson 
2160d0321e0SJeremy L Thompson   impl->numein = numinputfields; impl->numeout = numoutputfields;
2170d0321e0SJeremy L Thompson 
2180d0321e0SJeremy L Thompson   // Set up infield and outfield evecs and qvecs
2190d0321e0SJeremy L Thompson   // Infields
2200d0321e0SJeremy L Thompson   ierr = CeedOperatorSetupFields_Cuda(qf, op, true,
2210d0321e0SJeremy L Thompson                                       impl->evecs, impl->qvecsin, 0,
2220d0321e0SJeremy L Thompson                                       numinputfields, Q, numelements);
2230d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
2240d0321e0SJeremy L Thompson 
2250d0321e0SJeremy L Thompson   // Outfields
2260d0321e0SJeremy L Thompson   ierr = CeedOperatorSetupFields_Cuda(qf, op, false,
2270d0321e0SJeremy L Thompson                                       impl->evecs, impl->qvecsout,
2280d0321e0SJeremy L Thompson                                       numinputfields, numoutputfields, Q,
2290d0321e0SJeremy L Thompson                                       numelements); CeedChkBackend(ierr);
2300d0321e0SJeremy L Thompson 
2310d0321e0SJeremy L Thompson   ierr = CeedOperatorSetSetupDone(op); CeedChkBackend(ierr);
2320d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2330d0321e0SJeremy L Thompson }
2340d0321e0SJeremy L Thompson 
2350d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
2360d0321e0SJeremy L Thompson // Setup Operator Inputs
2370d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
2380d0321e0SJeremy L Thompson static inline int CeedOperatorSetupInputs_Cuda(CeedInt numinputfields,
2390d0321e0SJeremy L Thompson     CeedQFunctionField *qfinputfields, CeedOperatorField *opinputfields,
2400d0321e0SJeremy L Thompson     CeedVector invec, const bool skipactive, CeedScalar *edata[2*CEED_FIELD_MAX],
2410d0321e0SJeremy L Thompson     CeedOperator_Cuda *impl, CeedRequest *request) {
2420d0321e0SJeremy L Thompson   CeedInt ierr;
2430d0321e0SJeremy L Thompson   CeedEvalMode emode;
2440d0321e0SJeremy L Thompson   CeedVector vec;
2450d0321e0SJeremy L Thompson   CeedElemRestriction Erestrict;
2460d0321e0SJeremy L Thompson 
2470d0321e0SJeremy L Thompson   for (CeedInt i = 0; i < numinputfields; i++) {
2480d0321e0SJeremy L Thompson     // Get input vector
2490d0321e0SJeremy L Thompson     ierr = CeedOperatorFieldGetVector(opinputfields[i], &vec); CeedChkBackend(ierr);
2500d0321e0SJeremy L Thompson     if (vec == CEED_VECTOR_ACTIVE) {
2510d0321e0SJeremy L Thompson       if (skipactive)
2520d0321e0SJeremy L Thompson         continue;
2530d0321e0SJeremy L Thompson       else
2540d0321e0SJeremy L Thompson         vec = invec;
2550d0321e0SJeremy L Thompson     }
2560d0321e0SJeremy L Thompson 
2570d0321e0SJeremy L Thompson     ierr = CeedQFunctionFieldGetEvalMode(qfinputfields[i], &emode);
2580d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
2590d0321e0SJeremy L Thompson     if (emode == CEED_EVAL_WEIGHT) { // Skip
2600d0321e0SJeremy L Thompson     } else {
2610d0321e0SJeremy L Thompson       // Get input vector
2620d0321e0SJeremy L Thompson       ierr = CeedOperatorFieldGetVector(opinputfields[i], &vec); CeedChkBackend(ierr);
2630d0321e0SJeremy L Thompson       // Get input element restriction
2640d0321e0SJeremy L Thompson       ierr = CeedOperatorFieldGetElemRestriction(opinputfields[i], &Erestrict);
2650d0321e0SJeremy L Thompson       CeedChkBackend(ierr);
2660d0321e0SJeremy L Thompson       if (vec == CEED_VECTOR_ACTIVE)
2670d0321e0SJeremy L Thompson         vec = invec;
2680d0321e0SJeremy L Thompson       // Restrict, if necessary
2690d0321e0SJeremy L Thompson       if (!impl->evecs[i]) {
2700d0321e0SJeremy L Thompson         // No restriction for this field; read data directly from vec.
2710d0321e0SJeremy L Thompson         ierr = CeedVectorGetArrayRead(vec, CEED_MEM_DEVICE,
2720d0321e0SJeremy L Thompson                                       (const CeedScalar **) &edata[i]);
2730d0321e0SJeremy L Thompson         CeedChkBackend(ierr);
2740d0321e0SJeremy L Thompson       } else {
2750d0321e0SJeremy L Thompson         ierr = CeedElemRestrictionApply(Erestrict, CEED_NOTRANSPOSE, vec,
2760d0321e0SJeremy L Thompson                                         impl->evecs[i], request); CeedChkBackend(ierr);
2770d0321e0SJeremy L Thompson         // Get evec
2780d0321e0SJeremy L Thompson         ierr = CeedVectorGetArrayRead(impl->evecs[i], CEED_MEM_DEVICE,
2790d0321e0SJeremy L Thompson                                       (const CeedScalar **) &edata[i]);
2800d0321e0SJeremy L Thompson         CeedChkBackend(ierr);
2810d0321e0SJeremy L Thompson       }
2820d0321e0SJeremy L Thompson     }
2830d0321e0SJeremy L Thompson   }
2840d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2850d0321e0SJeremy L Thompson }
2860d0321e0SJeremy L Thompson 
2870d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
2880d0321e0SJeremy L Thompson // Input Basis Action
2890d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
2900d0321e0SJeremy L Thompson static inline int CeedOperatorInputBasis_Cuda(CeedInt numelements,
2910d0321e0SJeremy L Thompson     CeedQFunctionField *qfinputfields, CeedOperatorField *opinputfields,
2920d0321e0SJeremy L Thompson     CeedInt numinputfields, const bool skipactive,
2930d0321e0SJeremy L Thompson     CeedScalar *edata[2*CEED_FIELD_MAX],CeedOperator_Cuda *impl) {
2940d0321e0SJeremy L Thompson   CeedInt ierr;
2950d0321e0SJeremy L Thompson   CeedInt elemsize, size;
2960d0321e0SJeremy L Thompson   CeedElemRestriction Erestrict;
2970d0321e0SJeremy L Thompson   CeedEvalMode emode;
2980d0321e0SJeremy L Thompson   CeedBasis basis;
2990d0321e0SJeremy L Thompson 
3000d0321e0SJeremy L Thompson   for (CeedInt i=0; i<numinputfields; i++) {
3010d0321e0SJeremy L Thompson     // Skip active input
3020d0321e0SJeremy L Thompson     if (skipactive) {
3030d0321e0SJeremy L Thompson       CeedVector vec;
3040d0321e0SJeremy L Thompson       ierr = CeedOperatorFieldGetVector(opinputfields[i], &vec); CeedChkBackend(ierr);
3050d0321e0SJeremy L Thompson       if (vec == CEED_VECTOR_ACTIVE)
3060d0321e0SJeremy L Thompson         continue;
3070d0321e0SJeremy L Thompson     }
3080d0321e0SJeremy L Thompson     // Get elemsize, emode, size
3090d0321e0SJeremy L Thompson     ierr = CeedOperatorFieldGetElemRestriction(opinputfields[i], &Erestrict);
3100d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
3110d0321e0SJeremy L Thompson     ierr = CeedElemRestrictionGetElementSize(Erestrict, &elemsize);
3120d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
3130d0321e0SJeremy L Thompson     ierr = CeedQFunctionFieldGetEvalMode(qfinputfields[i], &emode);
3140d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
3150d0321e0SJeremy L Thompson     ierr = CeedQFunctionFieldGetSize(qfinputfields[i], &size); CeedChkBackend(ierr);
3160d0321e0SJeremy L Thompson     // Basis action
3170d0321e0SJeremy L Thompson     switch (emode) {
3180d0321e0SJeremy L Thompson     case CEED_EVAL_NONE:
3190d0321e0SJeremy L Thompson       ierr = CeedVectorSetArray(impl->qvecsin[i], CEED_MEM_DEVICE,
3200d0321e0SJeremy L Thompson                                 CEED_USE_POINTER, edata[i]); CeedChkBackend(ierr);
3210d0321e0SJeremy L Thompson       break;
3220d0321e0SJeremy L Thompson     case CEED_EVAL_INTERP:
3230d0321e0SJeremy L Thompson       ierr = CeedOperatorFieldGetBasis(opinputfields[i], &basis);
3240d0321e0SJeremy L Thompson       CeedChkBackend(ierr);
3250d0321e0SJeremy L Thompson       ierr = CeedBasisApply(basis, numelements, CEED_NOTRANSPOSE,
3260d0321e0SJeremy L Thompson                             CEED_EVAL_INTERP, impl->evecs[i],
3270d0321e0SJeremy L Thompson                             impl->qvecsin[i]); CeedChkBackend(ierr);
3280d0321e0SJeremy L Thompson       break;
3290d0321e0SJeremy L Thompson     case CEED_EVAL_GRAD:
3300d0321e0SJeremy L Thompson       ierr = CeedOperatorFieldGetBasis(opinputfields[i], &basis);
3310d0321e0SJeremy L Thompson       CeedChkBackend(ierr);
3320d0321e0SJeremy L Thompson       ierr = CeedBasisApply(basis, numelements, CEED_NOTRANSPOSE,
3330d0321e0SJeremy L Thompson                             CEED_EVAL_GRAD, impl->evecs[i],
3340d0321e0SJeremy L Thompson                             impl->qvecsin[i]); CeedChkBackend(ierr);
3350d0321e0SJeremy L Thompson       break;
3360d0321e0SJeremy L Thompson     case CEED_EVAL_WEIGHT:
3370d0321e0SJeremy L Thompson       break; // No action
3380d0321e0SJeremy L Thompson     case CEED_EVAL_DIV:
3390d0321e0SJeremy L Thompson       break; // TODO: Not implemented
3400d0321e0SJeremy L Thompson     case CEED_EVAL_CURL:
3410d0321e0SJeremy L Thompson       break; // TODO: Not implemented
3420d0321e0SJeremy L Thompson     }
3430d0321e0SJeremy L Thompson   }
3440d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3450d0321e0SJeremy L Thompson }
3460d0321e0SJeremy L Thompson 
3470d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
3480d0321e0SJeremy L Thompson // Restore Input Vectors
3490d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
3500d0321e0SJeremy L Thompson static inline int CeedOperatorRestoreInputs_Cuda(CeedInt numinputfields,
3510d0321e0SJeremy L Thompson     CeedQFunctionField *qfinputfields, CeedOperatorField *opinputfields,
3520d0321e0SJeremy L Thompson     const bool skipactive, CeedScalar *edata[2*CEED_FIELD_MAX],
3530d0321e0SJeremy L Thompson     CeedOperator_Cuda *impl) {
3540d0321e0SJeremy L Thompson   CeedInt ierr;
3550d0321e0SJeremy L Thompson   CeedEvalMode emode;
3560d0321e0SJeremy L Thompson   CeedVector vec;
3570d0321e0SJeremy L Thompson 
3580d0321e0SJeremy L Thompson   for (CeedInt i = 0; i < numinputfields; i++) {
3590d0321e0SJeremy L Thompson     // Skip active input
3600d0321e0SJeremy L Thompson     if (skipactive) {
3610d0321e0SJeremy L Thompson       ierr = CeedOperatorFieldGetVector(opinputfields[i], &vec); CeedChkBackend(ierr);
3620d0321e0SJeremy L Thompson       if (vec == CEED_VECTOR_ACTIVE)
3630d0321e0SJeremy L Thompson         continue;
3640d0321e0SJeremy L Thompson     }
3650d0321e0SJeremy L Thompson     ierr = CeedQFunctionFieldGetEvalMode(qfinputfields[i], &emode);
3660d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
3670d0321e0SJeremy L Thompson     if (emode == CEED_EVAL_WEIGHT) { // Skip
3680d0321e0SJeremy L Thompson     } else {
3690d0321e0SJeremy L Thompson       if (!impl->evecs[i]) {  // This was a skiprestrict case
3700d0321e0SJeremy L Thompson         ierr = CeedOperatorFieldGetVector(opinputfields[i], &vec); CeedChkBackend(ierr);
3710d0321e0SJeremy L Thompson         ierr = CeedVectorRestoreArrayRead(vec,
3720d0321e0SJeremy L Thompson                                           (const CeedScalar **)&edata[i]);
3730d0321e0SJeremy L Thompson         CeedChkBackend(ierr);
3740d0321e0SJeremy L Thompson       } else {
3750d0321e0SJeremy L Thompson         ierr = CeedVectorRestoreArrayRead(impl->evecs[i],
3760d0321e0SJeremy L Thompson                                           (const CeedScalar **) &edata[i]);
3770d0321e0SJeremy L Thompson         CeedChkBackend(ierr);
3780d0321e0SJeremy L Thompson       }
3790d0321e0SJeremy L Thompson     }
3800d0321e0SJeremy L Thompson   }
3810d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3820d0321e0SJeremy L Thompson }
3830d0321e0SJeremy L Thompson 
3840d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
3850d0321e0SJeremy L Thompson // Apply and add to output
3860d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
3870d0321e0SJeremy L Thompson static int CeedOperatorApplyAdd_Cuda(CeedOperator op, CeedVector invec,
3880d0321e0SJeremy L Thompson                                      CeedVector outvec, CeedRequest *request) {
3890d0321e0SJeremy L Thompson   int ierr;
3900d0321e0SJeremy L Thompson   CeedOperator_Cuda *impl;
3910d0321e0SJeremy L Thompson   ierr = CeedOperatorGetData(op, &impl); CeedChkBackend(ierr);
3920d0321e0SJeremy L Thompson   CeedQFunction qf;
3930d0321e0SJeremy L Thompson   ierr = CeedOperatorGetQFunction(op, &qf); CeedChkBackend(ierr);
3940d0321e0SJeremy L Thompson   CeedInt Q, numelements, elemsize, numinputfields, numoutputfields, size;
3950d0321e0SJeremy L Thompson   ierr = CeedOperatorGetNumQuadraturePoints(op, &Q); CeedChkBackend(ierr);
3960d0321e0SJeremy L Thompson   ierr = CeedOperatorGetNumElements(op, &numelements); CeedChkBackend(ierr);
3970d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
3980d0321e0SJeremy L Thompson   CeedOperatorField *opinputfields, *opoutputfields;
3990d0321e0SJeremy L Thompson   ierr = CeedOperatorGetFields(op, &numinputfields, &opinputfields,
4000d0321e0SJeremy L Thompson                                &numoutputfields, &opoutputfields);
4010d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
4020d0321e0SJeremy L Thompson   CeedQFunctionField *qfinputfields, *qfoutputfields;
4030d0321e0SJeremy L Thompson   ierr = CeedQFunctionGetFields(qf, NULL, &qfinputfields, NULL, &qfoutputfields);
4040d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
4050d0321e0SJeremy L Thompson   CeedEvalMode emode;
4060d0321e0SJeremy L Thompson   CeedVector vec;
4070d0321e0SJeremy L Thompson   CeedBasis basis;
4080d0321e0SJeremy L Thompson   CeedElemRestriction Erestrict;
4090d0321e0SJeremy L Thompson   CeedScalar *edata[2*CEED_FIELD_MAX] = {0};
4100d0321e0SJeremy L Thompson 
4110d0321e0SJeremy L Thompson   // Setup
4120d0321e0SJeremy L Thompson   ierr = CeedOperatorSetup_Cuda(op); CeedChkBackend(ierr);
4130d0321e0SJeremy L Thompson 
4140d0321e0SJeremy L Thompson   // Input Evecs and Restriction
4150d0321e0SJeremy L Thompson   ierr = CeedOperatorSetupInputs_Cuda(numinputfields, qfinputfields,
4160d0321e0SJeremy L Thompson                                       opinputfields, invec, false, edata,
4170d0321e0SJeremy L Thompson                                       impl, request); CeedChkBackend(ierr);
4180d0321e0SJeremy L Thompson 
4190d0321e0SJeremy L Thompson   // Input basis apply if needed
4200d0321e0SJeremy L Thompson   ierr = CeedOperatorInputBasis_Cuda(numelements, qfinputfields, opinputfields,
4210d0321e0SJeremy L Thompson                                      numinputfields, false, edata, impl);
4220d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
4230d0321e0SJeremy L Thompson 
4240d0321e0SJeremy L Thompson   // Output pointers, as necessary
4250d0321e0SJeremy L Thompson   for (CeedInt i = 0; i < numoutputfields; i++) {
4260d0321e0SJeremy L Thompson     ierr = CeedQFunctionFieldGetEvalMode(qfoutputfields[i], &emode);
4270d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
4280d0321e0SJeremy L Thompson     if (emode == CEED_EVAL_NONE) {
4290d0321e0SJeremy L Thompson       // Set the output Q-Vector to use the E-Vector data directly.
4300d0321e0SJeremy L Thompson       ierr = CeedVectorGetArrayWrite(impl->evecs[i + impl->numein], CEED_MEM_DEVICE,
4310d0321e0SJeremy L Thompson                                      &edata[i + numinputfields]); CeedChkBackend(ierr);
4320d0321e0SJeremy L Thompson       ierr = CeedVectorSetArray(impl->qvecsout[i], CEED_MEM_DEVICE,
4330d0321e0SJeremy L Thompson                                 CEED_USE_POINTER, edata[i + numinputfields]);
4340d0321e0SJeremy L Thompson       CeedChkBackend(ierr);
4350d0321e0SJeremy L Thompson     }
4360d0321e0SJeremy L Thompson   }
4370d0321e0SJeremy L Thompson 
4380d0321e0SJeremy L Thompson   // Q function
4390d0321e0SJeremy L Thompson   ierr = CeedQFunctionApply(qf, numelements * Q, impl->qvecsin, impl->qvecsout);
4400d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
4410d0321e0SJeremy L Thompson 
4420d0321e0SJeremy L Thompson   // Output basis apply if needed
4430d0321e0SJeremy L Thompson   for (CeedInt i = 0; i < numoutputfields; i++) {
4440d0321e0SJeremy L Thompson     // Get elemsize, emode, size
4450d0321e0SJeremy L Thompson     ierr = CeedOperatorFieldGetElemRestriction(opoutputfields[i], &Erestrict);
4460d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
4470d0321e0SJeremy L Thompson     ierr = CeedElemRestrictionGetElementSize(Erestrict, &elemsize);
4480d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
4490d0321e0SJeremy L Thompson     ierr = CeedQFunctionFieldGetEvalMode(qfoutputfields[i], &emode);
4500d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
4510d0321e0SJeremy L Thompson     ierr = CeedQFunctionFieldGetSize(qfoutputfields[i], &size);
4520d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
4530d0321e0SJeremy L Thompson     // Basis action
4540d0321e0SJeremy L Thompson     switch (emode) {
4550d0321e0SJeremy L Thompson     case CEED_EVAL_NONE:
4560d0321e0SJeremy L Thompson       break;
4570d0321e0SJeremy L Thompson     case CEED_EVAL_INTERP:
4580d0321e0SJeremy L Thompson       ierr = CeedOperatorFieldGetBasis(opoutputfields[i], &basis);
4590d0321e0SJeremy L Thompson       CeedChkBackend(ierr);
4600d0321e0SJeremy L Thompson       ierr = CeedBasisApply(basis, numelements, CEED_TRANSPOSE,
4610d0321e0SJeremy L Thompson                             CEED_EVAL_INTERP, impl->qvecsout[i],
4620d0321e0SJeremy L Thompson                             impl->evecs[i + impl->numein]); CeedChkBackend(ierr);
4630d0321e0SJeremy L Thompson       break;
4640d0321e0SJeremy L Thompson     case CEED_EVAL_GRAD:
4650d0321e0SJeremy L Thompson       ierr = CeedOperatorFieldGetBasis(opoutputfields[i], &basis);
4660d0321e0SJeremy L Thompson       CeedChkBackend(ierr);
4670d0321e0SJeremy L Thompson       ierr = CeedBasisApply(basis, numelements, CEED_TRANSPOSE,
4680d0321e0SJeremy L Thompson                             CEED_EVAL_GRAD, impl->qvecsout[i],
4690d0321e0SJeremy L Thompson                             impl->evecs[i + impl->numein]); CeedChkBackend(ierr);
4700d0321e0SJeremy L Thompson       break;
4710d0321e0SJeremy L Thompson     // LCOV_EXCL_START
4720d0321e0SJeremy L Thompson     case CEED_EVAL_WEIGHT: {
4730d0321e0SJeremy L Thompson       Ceed ceed;
4740d0321e0SJeremy L Thompson       ierr = CeedOperatorGetCeed(op, &ceed); CeedChkBackend(ierr);
4750d0321e0SJeremy L Thompson       return CeedError(ceed, CEED_ERROR_BACKEND,
4760d0321e0SJeremy L Thompson                        "CEED_EVAL_WEIGHT cannot be an output evaluation mode");
4770d0321e0SJeremy L Thompson       break; // Should not occur
4780d0321e0SJeremy L Thompson     }
4790d0321e0SJeremy L Thompson     case CEED_EVAL_DIV:
4800d0321e0SJeremy L Thompson       break; // TODO: Not implemented
4810d0321e0SJeremy L Thompson     case CEED_EVAL_CURL:
4820d0321e0SJeremy L Thompson       break; // TODO: Not implemented
4830d0321e0SJeremy L Thompson       // LCOV_EXCL_STOP
4840d0321e0SJeremy L Thompson     }
4850d0321e0SJeremy L Thompson   }
4860d0321e0SJeremy L Thompson 
4870d0321e0SJeremy L Thompson   // Output restriction
4880d0321e0SJeremy L Thompson   for (CeedInt i = 0; i < numoutputfields; i++) {
4890d0321e0SJeremy L Thompson     // Restore evec
4900d0321e0SJeremy L Thompson     ierr = CeedQFunctionFieldGetEvalMode(qfoutputfields[i], &emode);
4910d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
4920d0321e0SJeremy L Thompson     if (emode == CEED_EVAL_NONE) {
4930d0321e0SJeremy L Thompson       ierr = CeedVectorRestoreArray(impl->evecs[i+impl->numein],
4940d0321e0SJeremy L Thompson                                     &edata[i + numinputfields]);
4950d0321e0SJeremy L Thompson       CeedChkBackend(ierr);
4960d0321e0SJeremy L Thompson     }
4970d0321e0SJeremy L Thompson     // Get output vector
4980d0321e0SJeremy L Thompson     ierr = CeedOperatorFieldGetVector(opoutputfields[i], &vec);
4990d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
5000d0321e0SJeremy L Thompson     // Restrict
5010d0321e0SJeremy L Thompson     ierr = CeedOperatorFieldGetElemRestriction(opoutputfields[i], &Erestrict);
5020d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
5030d0321e0SJeremy L Thompson     // Active
5040d0321e0SJeremy L Thompson     if (vec == CEED_VECTOR_ACTIVE)
5050d0321e0SJeremy L Thompson       vec = outvec;
5060d0321e0SJeremy L Thompson 
5070d0321e0SJeremy L Thompson     ierr = CeedElemRestrictionApply(Erestrict, CEED_TRANSPOSE,
5080d0321e0SJeremy L Thompson                                     impl->evecs[i + impl->numein], vec,
5090d0321e0SJeremy L Thompson                                     request); CeedChkBackend(ierr);
5100d0321e0SJeremy L Thompson   }
5110d0321e0SJeremy L Thompson 
5120d0321e0SJeremy L Thompson   // Restore input arrays
5130d0321e0SJeremy L Thompson   ierr = CeedOperatorRestoreInputs_Cuda(numinputfields, qfinputfields,
5140d0321e0SJeremy L Thompson                                         opinputfields, false, edata, impl);
5150d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
5160d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
5170d0321e0SJeremy L Thompson }
5180d0321e0SJeremy L Thompson 
5190d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
5200d0321e0SJeremy L Thompson // Core code for assembling linear QFunction
5210d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
5220d0321e0SJeremy L Thompson static inline int CeedOperatorLinearAssembleQFunctionCore_Cuda(CeedOperator op,
5230d0321e0SJeremy L Thompson     bool build_objects, CeedVector *assembled, CeedElemRestriction *rstr,
5240d0321e0SJeremy L Thompson     CeedRequest *request) {
5250d0321e0SJeremy L Thompson   int ierr;
5260d0321e0SJeremy L Thompson   CeedOperator_Cuda *impl;
5270d0321e0SJeremy L Thompson   ierr = CeedOperatorGetData(op, &impl); CeedChkBackend(ierr);
5280d0321e0SJeremy L Thompson   CeedQFunction qf;
5290d0321e0SJeremy L Thompson   ierr = CeedOperatorGetQFunction(op, &qf); CeedChkBackend(ierr);
5300d0321e0SJeremy L Thompson   CeedInt Q, numelements, numinputfields, numoutputfields, size;
5310d0321e0SJeremy L Thompson   ierr = CeedOperatorGetNumQuadraturePoints(op, &Q); CeedChkBackend(ierr);
5320d0321e0SJeremy L Thompson   ierr = CeedOperatorGetNumElements(op, &numelements); CeedChkBackend(ierr);
5330d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
5340d0321e0SJeremy L Thompson   CeedOperatorField *opinputfields, *opoutputfields;
5350d0321e0SJeremy L Thompson   ierr = CeedOperatorGetFields(op, &numinputfields, &opinputfields,
5360d0321e0SJeremy L Thompson                                &numoutputfields, &opoutputfields);
5370d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
5380d0321e0SJeremy L Thompson   CeedQFunctionField *qfinputfields, *qfoutputfields;
5390d0321e0SJeremy L Thompson   ierr = CeedQFunctionGetFields(qf, NULL, &qfinputfields, NULL, &qfoutputfields);
5400d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
5410d0321e0SJeremy L Thompson   CeedVector vec;
5420d0321e0SJeremy L Thompson   CeedInt numactivein = impl->qfnumactivein, numactiveout = impl->qfnumactiveout;
5430d0321e0SJeremy L Thompson   CeedVector *activein = impl->qfactivein;
5440d0321e0SJeremy L Thompson   CeedScalar *a, *tmp;
5450d0321e0SJeremy L Thompson   Ceed ceed, ceedparent;
5460d0321e0SJeremy L Thompson   ierr = CeedOperatorGetCeed(op, &ceed); CeedChkBackend(ierr);
5470d0321e0SJeremy L Thompson   ierr = CeedGetOperatorFallbackParentCeed(ceed, &ceedparent);
5480d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
5490d0321e0SJeremy L Thompson   ceedparent = ceedparent ? ceedparent : ceed;
5500d0321e0SJeremy L Thompson   CeedScalar *edata[2*CEED_FIELD_MAX];
5510d0321e0SJeremy L Thompson 
5520d0321e0SJeremy L Thompson   // Setup
5530d0321e0SJeremy L Thompson   ierr = CeedOperatorSetup_Cuda(op); CeedChkBackend(ierr);
5540d0321e0SJeremy L Thompson 
5550d0321e0SJeremy L Thompson   // Check for identity
5560d0321e0SJeremy L Thompson   bool identityqf;
5570d0321e0SJeremy L Thompson   ierr = CeedQFunctionIsIdentity(qf, &identityqf); CeedChkBackend(ierr);
5580d0321e0SJeremy L Thompson   if (identityqf)
5590d0321e0SJeremy L Thompson     // LCOV_EXCL_START
5600d0321e0SJeremy L Thompson     return CeedError(ceed, CEED_ERROR_BACKEND,
5610d0321e0SJeremy L Thompson                      "Assembling identity QFunctions not supported");
5620d0321e0SJeremy L Thompson   // LCOV_EXCL_STOP
5630d0321e0SJeremy L Thompson 
5640d0321e0SJeremy L Thompson   // Input Evecs and Restriction
5650d0321e0SJeremy L Thompson   ierr = CeedOperatorSetupInputs_Cuda(numinputfields, qfinputfields,
5660d0321e0SJeremy L Thompson                                       opinputfields, NULL, true, edata,
5670d0321e0SJeremy L Thompson                                       impl, request); CeedChkBackend(ierr);
5680d0321e0SJeremy L Thompson 
5690d0321e0SJeremy L Thompson   // Count number of active input fields
5700d0321e0SJeremy L Thompson   if (!numactivein) {
5710d0321e0SJeremy L Thompson     for (CeedInt i=0; i<numinputfields; i++) {
5720d0321e0SJeremy L Thompson       // Get input vector
5730d0321e0SJeremy L Thompson       ierr = CeedOperatorFieldGetVector(opinputfields[i], &vec); CeedChkBackend(ierr);
5740d0321e0SJeremy L Thompson       // Check if active input
5750d0321e0SJeremy L Thompson       if (vec == CEED_VECTOR_ACTIVE) {
5760d0321e0SJeremy L Thompson         ierr = CeedQFunctionFieldGetSize(qfinputfields[i], &size); CeedChkBackend(ierr);
5770d0321e0SJeremy L Thompson         ierr = CeedVectorSetValue(impl->qvecsin[i], 0.0); CeedChkBackend(ierr);
5780d0321e0SJeremy L Thompson         ierr = CeedVectorGetArray(impl->qvecsin[i], CEED_MEM_DEVICE, &tmp);
5790d0321e0SJeremy L Thompson         CeedChkBackend(ierr);
5800d0321e0SJeremy L Thompson         ierr = CeedRealloc(numactivein + size, &activein); CeedChkBackend(ierr);
5810d0321e0SJeremy L Thompson         for (CeedInt field = 0; field < size; field++) {
5820d0321e0SJeremy L Thompson           ierr = CeedVectorCreate(ceed, Q*numelements,
5830d0321e0SJeremy L Thompson                                   &activein[numactivein+field]); CeedChkBackend(ierr);
5840d0321e0SJeremy L Thompson           ierr = CeedVectorSetArray(activein[numactivein+field], CEED_MEM_DEVICE,
5850d0321e0SJeremy L Thompson                                     CEED_USE_POINTER, &tmp[field*Q*numelements]);
5860d0321e0SJeremy L Thompson           CeedChkBackend(ierr);
5870d0321e0SJeremy L Thompson         }
5880d0321e0SJeremy L Thompson         numactivein += size;
5890d0321e0SJeremy L Thompson         ierr = CeedVectorRestoreArray(impl->qvecsin[i], &tmp); CeedChkBackend(ierr);
5900d0321e0SJeremy L Thompson       }
5910d0321e0SJeremy L Thompson     }
5920d0321e0SJeremy L Thompson     impl->qfnumactivein = numactivein;
5930d0321e0SJeremy L Thompson     impl->qfactivein = activein;
5940d0321e0SJeremy L Thompson   }
5950d0321e0SJeremy L Thompson 
5960d0321e0SJeremy L Thompson   // Count number of active output fields
5970d0321e0SJeremy L Thompson   if (!numactiveout) {
5980d0321e0SJeremy L Thompson     for (CeedInt i=0; i<numoutputfields; i++) {
5990d0321e0SJeremy L Thompson       // Get output vector
6000d0321e0SJeremy L Thompson       ierr = CeedOperatorFieldGetVector(opoutputfields[i], &vec);
6010d0321e0SJeremy L Thompson       CeedChkBackend(ierr);
6020d0321e0SJeremy L Thompson       // Check if active output
6030d0321e0SJeremy L Thompson       if (vec == CEED_VECTOR_ACTIVE) {
6040d0321e0SJeremy L Thompson         ierr = CeedQFunctionFieldGetSize(qfoutputfields[i], &size);
6050d0321e0SJeremy L Thompson         CeedChkBackend(ierr);
6060d0321e0SJeremy L Thompson         numactiveout += size;
6070d0321e0SJeremy L Thompson       }
6080d0321e0SJeremy L Thompson     }
6090d0321e0SJeremy L Thompson     impl->qfnumactiveout = numactiveout;
6100d0321e0SJeremy L Thompson   }
6110d0321e0SJeremy L Thompson 
6120d0321e0SJeremy L Thompson   // Check sizes
6130d0321e0SJeremy L Thompson   if (!numactivein || !numactiveout)
6140d0321e0SJeremy L Thompson     // LCOV_EXCL_START
6150d0321e0SJeremy L Thompson     return CeedError(ceed, CEED_ERROR_BACKEND,
6160d0321e0SJeremy L Thompson                      "Cannot assemble QFunction without active inputs "
6170d0321e0SJeremy L Thompson                      "and outputs");
6180d0321e0SJeremy L Thompson   // LCOV_EXCL_STOP
6190d0321e0SJeremy L Thompson 
6200d0321e0SJeremy L Thompson   // Build objects if needed
6210d0321e0SJeremy L Thompson   if (build_objects) {
6220d0321e0SJeremy L Thompson     // Create output restriction
6230d0321e0SJeremy L Thompson     CeedInt strides[3] = {1, numelements*Q, Q}; /* *NOPAD* */
6240d0321e0SJeremy L Thompson     ierr = CeedElemRestrictionCreateStrided(ceedparent, numelements, Q,
6250d0321e0SJeremy L Thompson                                             numactivein*numactiveout,
6260d0321e0SJeremy L Thompson                                             numactivein*numactiveout*numelements*Q,
6270d0321e0SJeremy L Thompson                                             strides, rstr); CeedChkBackend(ierr);
6280d0321e0SJeremy L Thompson     // Create assembled vector
6290d0321e0SJeremy L Thompson     ierr = CeedVectorCreate(ceedparent, numelements*Q*numactivein*numactiveout,
6300d0321e0SJeremy L Thompson                             assembled); CeedChkBackend(ierr);
6310d0321e0SJeremy L Thompson   }
6320d0321e0SJeremy L Thompson   ierr = CeedVectorSetValue(*assembled, 0.0); CeedChkBackend(ierr);
6330d0321e0SJeremy L Thompson   ierr = CeedVectorGetArray(*assembled, CEED_MEM_DEVICE, &a);
6340d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
6350d0321e0SJeremy L Thompson 
6360d0321e0SJeremy L Thompson   // Input basis apply
6370d0321e0SJeremy L Thompson   ierr = CeedOperatorInputBasis_Cuda(numelements, qfinputfields, opinputfields,
6380d0321e0SJeremy L Thompson                                      numinputfields, true, edata, impl);
6390d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
6400d0321e0SJeremy L Thompson 
6410d0321e0SJeremy L Thompson   // Assemble QFunction
6420d0321e0SJeremy L Thompson   for (CeedInt in=0; in<numactivein; in++) {
6430d0321e0SJeremy L Thompson     // Set Inputs
6440d0321e0SJeremy L Thompson     ierr = CeedVectorSetValue(activein[in], 1.0); CeedChkBackend(ierr);
6450d0321e0SJeremy L Thompson     if (numactivein > 1) {
6460d0321e0SJeremy L Thompson       ierr = CeedVectorSetValue(activein[(in+numactivein-1)%numactivein],
6470d0321e0SJeremy L Thompson                                 0.0); CeedChkBackend(ierr);
6480d0321e0SJeremy L Thompson     }
6490d0321e0SJeremy L Thompson     // Set Outputs
6500d0321e0SJeremy L Thompson     for (CeedInt out=0; out<numoutputfields; out++) {
6510d0321e0SJeremy L Thompson       // Get output vector
6520d0321e0SJeremy L Thompson       ierr = CeedOperatorFieldGetVector(opoutputfields[out], &vec);
6530d0321e0SJeremy L Thompson       CeedChkBackend(ierr);
6540d0321e0SJeremy L Thompson       // Check if active output
6550d0321e0SJeremy L Thompson       if (vec == CEED_VECTOR_ACTIVE) {
6560d0321e0SJeremy L Thompson         CeedVectorSetArray(impl->qvecsout[out], CEED_MEM_DEVICE,
6570d0321e0SJeremy L Thompson                            CEED_USE_POINTER, a); CeedChkBackend(ierr);
6580d0321e0SJeremy L Thompson         ierr = CeedQFunctionFieldGetSize(qfoutputfields[out], &size);
6590d0321e0SJeremy L Thompson         CeedChkBackend(ierr);
6600d0321e0SJeremy L Thompson         a += size*Q*numelements; // Advance the pointer by the size of the output
6610d0321e0SJeremy L Thompson       }
6620d0321e0SJeremy L Thompson     }
6630d0321e0SJeremy L Thompson     // Apply QFunction
6640d0321e0SJeremy L Thompson     ierr = CeedQFunctionApply(qf, Q*numelements, impl->qvecsin, impl->qvecsout);
6650d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
6660d0321e0SJeremy L Thompson   }
6670d0321e0SJeremy L Thompson 
6680d0321e0SJeremy L Thompson   // Un-set output Qvecs to prevent accidental overwrite of Assembled
6690d0321e0SJeremy L Thompson   for (CeedInt out=0; out<numoutputfields; out++) {
6700d0321e0SJeremy L Thompson     // Get output vector
6710d0321e0SJeremy L Thompson     ierr = CeedOperatorFieldGetVector(opoutputfields[out], &vec);
6720d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
6730d0321e0SJeremy L Thompson     // Check if active output
6740d0321e0SJeremy L Thompson     if (vec == CEED_VECTOR_ACTIVE) {
6750d0321e0SJeremy L Thompson       ierr = CeedVectorTakeArray(impl->qvecsout[out], CEED_MEM_DEVICE, NULL);
6760d0321e0SJeremy L Thompson       CeedChkBackend(ierr);
6770d0321e0SJeremy L Thompson     }
6780d0321e0SJeremy L Thompson   }
6790d0321e0SJeremy L Thompson 
6800d0321e0SJeremy L Thompson   // Restore input arrays
6810d0321e0SJeremy L Thompson   ierr = CeedOperatorRestoreInputs_Cuda(numinputfields, qfinputfields,
6820d0321e0SJeremy L Thompson                                         opinputfields, true, edata, impl);
6830d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
6840d0321e0SJeremy L Thompson 
6850d0321e0SJeremy L Thompson   // Restore output
6860d0321e0SJeremy L Thompson   ierr = CeedVectorRestoreArray(*assembled, &a); CeedChkBackend(ierr);
6870d0321e0SJeremy L Thompson 
6880d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
6890d0321e0SJeremy L Thompson }
6900d0321e0SJeremy L Thompson 
6910d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
6920d0321e0SJeremy L Thompson // Assemble Linear QFunction
6930d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
6940d0321e0SJeremy L Thompson static int CeedOperatorLinearAssembleQFunction_Cuda(CeedOperator op,
6950d0321e0SJeremy L Thompson     CeedVector *assembled, CeedElemRestriction *rstr, CeedRequest *request) {
6960d0321e0SJeremy L Thompson   return CeedOperatorLinearAssembleQFunctionCore_Cuda(op, true, assembled, rstr,
6970d0321e0SJeremy L Thompson          request);
6980d0321e0SJeremy L Thompson }
6990d0321e0SJeremy L Thompson 
7000d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
7010d0321e0SJeremy L Thompson // Update Assembled Linear QFunction
7020d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
7030d0321e0SJeremy L Thompson static int CeedOperatorLinearAssembleQFunctionUpdate_Cuda(CeedOperator op,
7040d0321e0SJeremy L Thompson     CeedVector assembled, CeedElemRestriction rstr, CeedRequest *request) {
7050d0321e0SJeremy L Thompson   return CeedOperatorLinearAssembleQFunctionCore_Cuda(op, false, &assembled,
7060d0321e0SJeremy L Thompson          &rstr, request);
7070d0321e0SJeremy L Thompson }
7080d0321e0SJeremy L Thompson 
7090d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
7100d0321e0SJeremy L Thompson // Diagonal assembly kernels
7110d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
7120d0321e0SJeremy L Thompson // *INDENT-OFF*
7130d0321e0SJeremy L Thompson static const char *diagonalkernels = QUOTE(
7140d0321e0SJeremy L Thompson 
7150d0321e0SJeremy L Thompson typedef enum {
7160d0321e0SJeremy L Thompson   /// Perform no evaluation (either because there is no data or it is already at
7170d0321e0SJeremy L Thompson   /// quadrature points)
7180d0321e0SJeremy L Thompson   CEED_EVAL_NONE   = 0,
7190d0321e0SJeremy L Thompson   /// Interpolate from nodes to quadrature points
7200d0321e0SJeremy L Thompson   CEED_EVAL_INTERP = 1,
7210d0321e0SJeremy L Thompson   /// Evaluate gradients at quadrature points from input in a nodal basis
7220d0321e0SJeremy L Thompson   CEED_EVAL_GRAD   = 2,
7230d0321e0SJeremy L Thompson   /// Evaluate divergence at quadrature points from input in a nodal basis
7240d0321e0SJeremy L Thompson   CEED_EVAL_DIV    = 4,
7250d0321e0SJeremy L Thompson   /// Evaluate curl at quadrature points from input in a nodal basis
7260d0321e0SJeremy L Thompson   CEED_EVAL_CURL   = 8,
7270d0321e0SJeremy L Thompson   /// Using no input, evaluate quadrature weights on the reference element
7280d0321e0SJeremy L Thompson   CEED_EVAL_WEIGHT = 16,
7290d0321e0SJeremy L Thompson } CeedEvalMode;
7300d0321e0SJeremy L Thompson 
7310d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
7320d0321e0SJeremy L Thompson // Get Basis Emode Pointer
7330d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
7340d0321e0SJeremy L Thompson extern "C" __device__ void CeedOperatorGetBasisPointer_Cuda(const CeedScalar **basisptr,
7350d0321e0SJeremy L Thompson     CeedEvalMode emode, const CeedScalar *identity, const CeedScalar *interp,
7360d0321e0SJeremy L Thompson     const CeedScalar *grad) {
7370d0321e0SJeremy L Thompson   switch (emode) {
7380d0321e0SJeremy L Thompson   case CEED_EVAL_NONE:
7390d0321e0SJeremy L Thompson     *basisptr = identity;
7400d0321e0SJeremy L Thompson     break;
7410d0321e0SJeremy L Thompson   case CEED_EVAL_INTERP:
7420d0321e0SJeremy L Thompson     *basisptr = interp;
7430d0321e0SJeremy L Thompson     break;
7440d0321e0SJeremy L Thompson   case CEED_EVAL_GRAD:
7450d0321e0SJeremy L Thompson     *basisptr = grad;
7460d0321e0SJeremy L Thompson     break;
7470d0321e0SJeremy L Thompson   case CEED_EVAL_WEIGHT:
7480d0321e0SJeremy L Thompson   case CEED_EVAL_DIV:
7490d0321e0SJeremy L Thompson   case CEED_EVAL_CURL:
7500d0321e0SJeremy L Thompson     break; // Caught by QF Assembly
7510d0321e0SJeremy L Thompson   }
7520d0321e0SJeremy L Thompson }
7530d0321e0SJeremy L Thompson 
7540d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
7550d0321e0SJeremy L Thompson // Core code for diagonal assembly
7560d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
7570d0321e0SJeremy L Thompson __device__ void diagonalCore(const CeedInt nelem,
7580d0321e0SJeremy L Thompson     const CeedScalar maxnorm, const bool pointBlock,
7590d0321e0SJeremy L Thompson     const CeedScalar *identity,
7600d0321e0SJeremy L Thompson     const CeedScalar *interpin, const CeedScalar *gradin,
7610d0321e0SJeremy L Thompson     const CeedScalar *interpout, const CeedScalar *gradout,
7620d0321e0SJeremy L Thompson     const CeedEvalMode *emodein, const CeedEvalMode *emodeout,
7630d0321e0SJeremy L Thompson     const CeedScalar *__restrict__ assembledqfarray,
7640d0321e0SJeremy L Thompson     CeedScalar *__restrict__ elemdiagarray) {
7650d0321e0SJeremy L Thompson   const int tid = threadIdx.x; // running with P threads, tid is evec node
7660d0321e0SJeremy L Thompson   const CeedScalar qfvaluebound = maxnorm*1e-12;
7670d0321e0SJeremy L Thompson 
7680d0321e0SJeremy L Thompson   // Compute the diagonal of B^T D B
7690d0321e0SJeremy L Thompson   // Each element
7700d0321e0SJeremy L Thompson   for (CeedInt e = blockIdx.x*blockDim.z + threadIdx.z; e < nelem;
7710d0321e0SJeremy L Thompson        e += gridDim.x*blockDim.z) {
7720d0321e0SJeremy L Thompson     CeedInt dout = -1;
7730d0321e0SJeremy L Thompson     // Each basis eval mode pair
7740d0321e0SJeremy L Thompson     for (CeedInt eout = 0; eout < NUMEMODEOUT; eout++) {
7750d0321e0SJeremy L Thompson       const CeedScalar *bt = NULL;
7760d0321e0SJeremy L Thompson       if (emodeout[eout] == CEED_EVAL_GRAD)
7770d0321e0SJeremy L Thompson         dout += 1;
7780d0321e0SJeremy L Thompson       CeedOperatorGetBasisPointer_Cuda(&bt, emodeout[eout], identity, interpout,
7790d0321e0SJeremy L Thompson                                       &gradout[dout*NQPTS*NNODES]);
7800d0321e0SJeremy L Thompson       CeedInt din = -1;
7810d0321e0SJeremy L Thompson       for (CeedInt ein = 0; ein < NUMEMODEIN; ein++) {
7820d0321e0SJeremy L Thompson         const CeedScalar *b = NULL;
7830d0321e0SJeremy L Thompson         if (emodein[ein] == CEED_EVAL_GRAD)
7840d0321e0SJeremy L Thompson           din += 1;
7850d0321e0SJeremy L Thompson         CeedOperatorGetBasisPointer_Cuda(&b, emodein[ein], identity, interpin,
7860d0321e0SJeremy L Thompson                                         &gradin[din*NQPTS*NNODES]);
7870d0321e0SJeremy L Thompson         // Each component
7880d0321e0SJeremy L Thompson         for (CeedInt compOut = 0; compOut < NCOMP; compOut++) {
7890d0321e0SJeremy L Thompson           // Each qpoint/node pair
7900d0321e0SJeremy L Thompson           if (pointBlock) {
7910d0321e0SJeremy L Thompson             // Point Block Diagonal
7920d0321e0SJeremy L Thompson             for (CeedInt compIn = 0; compIn < NCOMP; compIn++) {
7930d0321e0SJeremy L Thompson               CeedScalar evalue = 0.;
7940d0321e0SJeremy L Thompson               for (CeedInt q = 0; q < NQPTS; q++) {
7950d0321e0SJeremy L Thompson                 const CeedScalar qfvalue =
7960d0321e0SJeremy L Thompson                   assembledqfarray[((((ein*NCOMP+compIn)*NUMEMODEOUT+eout)*
7970d0321e0SJeremy L Thompson                                      NCOMP+compOut)*nelem+e)*NQPTS+q];
7980d0321e0SJeremy L Thompson                 if (abs(qfvalue) > qfvaluebound)
7990d0321e0SJeremy L Thompson                   evalue += bt[q*NNODES+tid] * qfvalue * b[q*NNODES+tid];
8000d0321e0SJeremy L Thompson               }
8010d0321e0SJeremy L Thompson               elemdiagarray[((compOut*NCOMP+compIn)*nelem+e)*NNODES+tid] += evalue;
8020d0321e0SJeremy L Thompson             }
8030d0321e0SJeremy L Thompson           } else {
8040d0321e0SJeremy L Thompson             // Diagonal Only
8050d0321e0SJeremy L Thompson             CeedScalar evalue = 0.;
8060d0321e0SJeremy L Thompson             for (CeedInt q = 0; q < NQPTS; q++) {
8070d0321e0SJeremy L Thompson               const CeedScalar qfvalue =
8080d0321e0SJeremy L Thompson                 assembledqfarray[((((ein*NCOMP+compOut)*NUMEMODEOUT+eout)*
8090d0321e0SJeremy L Thompson                                    NCOMP+compOut)*nelem+e)*NQPTS+q];
8100d0321e0SJeremy L Thompson               if (abs(qfvalue) > qfvaluebound)
8110d0321e0SJeremy L Thompson                 evalue += bt[q*NNODES+tid] * qfvalue * b[q*NNODES+tid];
8120d0321e0SJeremy L Thompson             }
8130d0321e0SJeremy L Thompson             elemdiagarray[(compOut*nelem+e)*NNODES+tid] += evalue;
8140d0321e0SJeremy L Thompson           }
8150d0321e0SJeremy L Thompson         }
8160d0321e0SJeremy L Thompson       }
8170d0321e0SJeremy L Thompson     }
8180d0321e0SJeremy L Thompson   }
8190d0321e0SJeremy L Thompson }
8200d0321e0SJeremy L Thompson 
8210d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
8220d0321e0SJeremy L Thompson // Linear diagonal
8230d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
8240d0321e0SJeremy L Thompson extern "C" __global__ void linearDiagonal(const CeedInt nelem,
8250d0321e0SJeremy L Thompson     const CeedScalar maxnorm, const CeedScalar *identity,
8260d0321e0SJeremy L Thompson     const CeedScalar *interpin, const CeedScalar *gradin,
8270d0321e0SJeremy L Thompson     const CeedScalar *interpout, const CeedScalar *gradout,
8280d0321e0SJeremy L Thompson     const CeedEvalMode *emodein, const CeedEvalMode *emodeout,
8290d0321e0SJeremy L Thompson     const CeedScalar *__restrict__ assembledqfarray,
8300d0321e0SJeremy L Thompson     CeedScalar *__restrict__ elemdiagarray) {
8310d0321e0SJeremy L Thompson   diagonalCore(nelem, maxnorm, false, identity, interpin, gradin, interpout,
8320d0321e0SJeremy L Thompson                gradout, emodein, emodeout, assembledqfarray, elemdiagarray);
8330d0321e0SJeremy L Thompson }
8340d0321e0SJeremy L Thompson 
8350d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
8360d0321e0SJeremy L Thompson // Linear point block diagonal
8370d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
8380d0321e0SJeremy L Thompson extern "C" __global__ void linearPointBlockDiagonal(const CeedInt nelem,
8390d0321e0SJeremy L Thompson     const CeedScalar maxnorm, const CeedScalar *identity,
8400d0321e0SJeremy L Thompson     const CeedScalar *interpin, const CeedScalar *gradin,
8410d0321e0SJeremy L Thompson     const CeedScalar *interpout, const CeedScalar *gradout,
8420d0321e0SJeremy L Thompson     const CeedEvalMode *emodein, const CeedEvalMode *emodeout,
8430d0321e0SJeremy L Thompson     const CeedScalar *__restrict__ assembledqfarray,
8440d0321e0SJeremy L Thompson     CeedScalar *__restrict__ elemdiagarray) {
8450d0321e0SJeremy L Thompson   diagonalCore(nelem, maxnorm, true, identity, interpin, gradin, interpout,
8460d0321e0SJeremy L Thompson                gradout, emodein, emodeout, assembledqfarray, elemdiagarray);
8470d0321e0SJeremy L Thompson }
8480d0321e0SJeremy L Thompson 
8490d0321e0SJeremy L Thompson );
8500d0321e0SJeremy L Thompson // *INDENT-ON*
8510d0321e0SJeremy L Thompson 
8520d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
8530d0321e0SJeremy L Thompson // Create point block restriction
8540d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
8550d0321e0SJeremy L Thompson static int CreatePBRestriction(CeedElemRestriction rstr,
8560d0321e0SJeremy L Thompson                                CeedElemRestriction *pbRstr) {
8570d0321e0SJeremy L Thompson   int ierr;
8580d0321e0SJeremy L Thompson   Ceed ceed;
8590d0321e0SJeremy L Thompson   ierr = CeedElemRestrictionGetCeed(rstr, &ceed); CeedChkBackend(ierr);
8600d0321e0SJeremy L Thompson   const CeedInt *offsets;
8610d0321e0SJeremy L Thompson   ierr = CeedElemRestrictionGetOffsets(rstr, CEED_MEM_HOST, &offsets);
8620d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
8630d0321e0SJeremy L Thompson 
8640d0321e0SJeremy L Thompson   // Expand offsets
8650d0321e0SJeremy L Thompson   CeedInt nelem, ncomp, elemsize, compstride, max = 1, *pbOffsets;
8660d0321e0SJeremy L Thompson   ierr = CeedElemRestrictionGetNumElements(rstr, &nelem); CeedChkBackend(ierr);
8670d0321e0SJeremy L Thompson   ierr = CeedElemRestrictionGetNumComponents(rstr, &ncomp); CeedChkBackend(ierr);
8680d0321e0SJeremy L Thompson   ierr = CeedElemRestrictionGetElementSize(rstr, &elemsize); CeedChkBackend(ierr);
8690d0321e0SJeremy L Thompson   ierr = CeedElemRestrictionGetCompStride(rstr, &compstride);
8700d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
8710d0321e0SJeremy L Thompson   CeedInt shift = ncomp;
8720d0321e0SJeremy L Thompson   if (compstride != 1)
8730d0321e0SJeremy L Thompson     shift *= ncomp;
8740d0321e0SJeremy L Thompson   ierr = CeedCalloc(nelem*elemsize, &pbOffsets); CeedChkBackend(ierr);
8750d0321e0SJeremy L Thompson   for (CeedInt i = 0; i < nelem*elemsize; i++) {
8760d0321e0SJeremy L Thompson     pbOffsets[i] = offsets[i]*shift;
8770d0321e0SJeremy L Thompson     if (pbOffsets[i] > max)
8780d0321e0SJeremy L Thompson       max = pbOffsets[i];
8790d0321e0SJeremy L Thompson   }
8800d0321e0SJeremy L Thompson 
8810d0321e0SJeremy L Thompson   // Create new restriction
8820d0321e0SJeremy L Thompson   ierr = CeedElemRestrictionCreate(ceed, nelem, elemsize, ncomp*ncomp, 1,
8830d0321e0SJeremy L Thompson                                    max + ncomp*ncomp, CEED_MEM_HOST,
8840d0321e0SJeremy L Thompson                                    CEED_OWN_POINTER, pbOffsets, pbRstr);
8850d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
8860d0321e0SJeremy L Thompson 
8870d0321e0SJeremy L Thompson   // Cleanup
8880d0321e0SJeremy L Thompson   ierr = CeedElemRestrictionRestoreOffsets(rstr, &offsets); CeedChkBackend(ierr);
8890d0321e0SJeremy L Thompson 
8900d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
8910d0321e0SJeremy L Thompson }
8920d0321e0SJeremy L Thompson 
8930d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
8940d0321e0SJeremy L Thompson // Assemble diagonal setup
8950d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
8960d0321e0SJeremy L Thompson static inline int CeedOperatorAssembleDiagonalSetup_Cuda(CeedOperator op,
8970d0321e0SJeremy L Thompson     const bool pointBlock) {
8980d0321e0SJeremy L Thompson   int ierr;
8990d0321e0SJeremy L Thompson   Ceed ceed;
9000d0321e0SJeremy L Thompson   ierr = CeedOperatorGetCeed(op, &ceed); CeedChkBackend(ierr);
9010d0321e0SJeremy L Thompson   CeedQFunction qf;
9020d0321e0SJeremy L Thompson   ierr = CeedOperatorGetQFunction(op, &qf); CeedChkBackend(ierr);
9030d0321e0SJeremy L Thompson   CeedInt numinputfields, numoutputfields;
9040d0321e0SJeremy L Thompson   ierr = CeedQFunctionGetNumArgs(qf, &numinputfields, &numoutputfields);
9050d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
9060d0321e0SJeremy L Thompson 
9070d0321e0SJeremy L Thompson   // Determine active input basis
9080d0321e0SJeremy L Thompson   CeedOperatorField *opfields;
9090d0321e0SJeremy L Thompson   CeedQFunctionField *qffields;
9100d0321e0SJeremy L Thompson   ierr = CeedOperatorGetFields(op, NULL, &opfields, NULL, NULL);
9110d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
9120d0321e0SJeremy L Thompson   ierr = CeedQFunctionGetFields(qf, NULL, &qffields, NULL, NULL);
9130d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
9140d0321e0SJeremy L Thompson   CeedInt numemodein = 0, ncomp = 0, dim = 1;
9150d0321e0SJeremy L Thompson   CeedEvalMode *emodein = NULL;
9160d0321e0SJeremy L Thompson   CeedBasis basisin = NULL;
9170d0321e0SJeremy L Thompson   CeedElemRestriction rstrin = NULL;
9180d0321e0SJeremy L Thompson   for (CeedInt i = 0; i < numinputfields; i++) {
9190d0321e0SJeremy L Thompson     CeedVector vec;
9200d0321e0SJeremy L Thompson     ierr = CeedOperatorFieldGetVector(opfields[i], &vec); CeedChkBackend(ierr);
9210d0321e0SJeremy L Thompson     if (vec == CEED_VECTOR_ACTIVE) {
9220d0321e0SJeremy L Thompson       CeedElemRestriction rstr;
9230d0321e0SJeremy L Thompson       ierr = CeedOperatorFieldGetBasis(opfields[i], &basisin); CeedChkBackend(ierr);
9240d0321e0SJeremy L Thompson       ierr = CeedBasisGetNumComponents(basisin, &ncomp); CeedChkBackend(ierr);
9250d0321e0SJeremy L Thompson       ierr = CeedBasisGetDimension(basisin, &dim); CeedChkBackend(ierr);
9260d0321e0SJeremy L Thompson       ierr = CeedOperatorFieldGetElemRestriction(opfields[i], &rstr);
9270d0321e0SJeremy L Thompson       CeedChkBackend(ierr);
9280d0321e0SJeremy L Thompson       if (rstrin && rstrin != rstr)
9290d0321e0SJeremy L Thompson         // LCOV_EXCL_START
9300d0321e0SJeremy L Thompson         return CeedError(ceed, CEED_ERROR_BACKEND,
9310d0321e0SJeremy L Thompson                          "Multi-field non-composite operator diagonal assembly not supported");
9320d0321e0SJeremy L Thompson       // LCOV_EXCL_STOP
9330d0321e0SJeremy L Thompson       rstrin = rstr;
9340d0321e0SJeremy L Thompson       CeedEvalMode emode;
9350d0321e0SJeremy L Thompson       ierr = CeedQFunctionFieldGetEvalMode(qffields[i], &emode);
9360d0321e0SJeremy L Thompson       CeedChkBackend(ierr);
9370d0321e0SJeremy L Thompson       switch (emode) {
9380d0321e0SJeremy L Thompson       case CEED_EVAL_NONE:
9390d0321e0SJeremy L Thompson       case CEED_EVAL_INTERP:
9400d0321e0SJeremy L Thompson         ierr = CeedRealloc(numemodein + 1, &emodein); CeedChkBackend(ierr);
9410d0321e0SJeremy L Thompson         emodein[numemodein] = emode;
9420d0321e0SJeremy L Thompson         numemodein += 1;
9430d0321e0SJeremy L Thompson         break;
9440d0321e0SJeremy L Thompson       case CEED_EVAL_GRAD:
9450d0321e0SJeremy L Thompson         ierr = CeedRealloc(numemodein + dim, &emodein); CeedChkBackend(ierr);
9460d0321e0SJeremy L Thompson         for (CeedInt d = 0; d < dim; d++)
9470d0321e0SJeremy L Thompson           emodein[numemodein+d] = emode;
9480d0321e0SJeremy L Thompson         numemodein += dim;
9490d0321e0SJeremy L Thompson         break;
9500d0321e0SJeremy L Thompson       case CEED_EVAL_WEIGHT:
9510d0321e0SJeremy L Thompson       case CEED_EVAL_DIV:
9520d0321e0SJeremy L Thompson       case CEED_EVAL_CURL:
9530d0321e0SJeremy L Thompson         break; // Caught by QF Assembly
9540d0321e0SJeremy L Thompson       }
9550d0321e0SJeremy L Thompson     }
9560d0321e0SJeremy L Thompson   }
9570d0321e0SJeremy L Thompson 
9580d0321e0SJeremy L Thompson   // Determine active output basis
9590d0321e0SJeremy L Thompson   ierr = CeedOperatorGetFields(op, NULL, NULL, NULL, &opfields);
9600d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
9610d0321e0SJeremy L Thompson   ierr = CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qffields);
9620d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
9630d0321e0SJeremy L Thompson   CeedInt numemodeout = 0;
9640d0321e0SJeremy L Thompson   CeedEvalMode *emodeout = NULL;
9650d0321e0SJeremy L Thompson   CeedBasis basisout = NULL;
9660d0321e0SJeremy L Thompson   CeedElemRestriction rstrout = NULL;
9670d0321e0SJeremy L Thompson   for (CeedInt i = 0; i < numoutputfields; i++) {
9680d0321e0SJeremy L Thompson     CeedVector vec;
9690d0321e0SJeremy L Thompson     ierr = CeedOperatorFieldGetVector(opfields[i], &vec); CeedChkBackend(ierr);
9700d0321e0SJeremy L Thompson     if (vec == CEED_VECTOR_ACTIVE) {
9710d0321e0SJeremy L Thompson       CeedElemRestriction rstr;
9720d0321e0SJeremy L Thompson       ierr = CeedOperatorFieldGetBasis(opfields[i], &basisout); CeedChkBackend(ierr);
9730d0321e0SJeremy L Thompson       ierr = CeedOperatorFieldGetElemRestriction(opfields[i], &rstr);
9740d0321e0SJeremy L Thompson       CeedChkBackend(ierr);
9750d0321e0SJeremy L Thompson       if (rstrout && rstrout != rstr)
9760d0321e0SJeremy L Thompson         // LCOV_EXCL_START
9770d0321e0SJeremy L Thompson         return CeedError(ceed, CEED_ERROR_BACKEND,
9780d0321e0SJeremy L Thompson                          "Multi-field non-composite operator diagonal assembly not supported");
9790d0321e0SJeremy L Thompson       // LCOV_EXCL_STOP
9800d0321e0SJeremy L Thompson       rstrout = rstr;
9810d0321e0SJeremy L Thompson       CeedEvalMode emode;
9820d0321e0SJeremy L Thompson       ierr = CeedQFunctionFieldGetEvalMode(qffields[i], &emode); CeedChkBackend(ierr);
9830d0321e0SJeremy L Thompson       switch (emode) {
9840d0321e0SJeremy L Thompson       case CEED_EVAL_NONE:
9850d0321e0SJeremy L Thompson       case CEED_EVAL_INTERP:
9860d0321e0SJeremy L Thompson         ierr = CeedRealloc(numemodeout + 1, &emodeout); CeedChkBackend(ierr);
9870d0321e0SJeremy L Thompson         emodeout[numemodeout] = emode;
9880d0321e0SJeremy L Thompson         numemodeout += 1;
9890d0321e0SJeremy L Thompson         break;
9900d0321e0SJeremy L Thompson       case CEED_EVAL_GRAD:
9910d0321e0SJeremy L Thompson         ierr = CeedRealloc(numemodeout + dim, &emodeout); CeedChkBackend(ierr);
9920d0321e0SJeremy L Thompson         for (CeedInt d = 0; d < dim; d++)
9930d0321e0SJeremy L Thompson           emodeout[numemodeout+d] = emode;
9940d0321e0SJeremy L Thompson         numemodeout += dim;
9950d0321e0SJeremy L Thompson         break;
9960d0321e0SJeremy L Thompson       case CEED_EVAL_WEIGHT:
9970d0321e0SJeremy L Thompson       case CEED_EVAL_DIV:
9980d0321e0SJeremy L Thompson       case CEED_EVAL_CURL:
9990d0321e0SJeremy L Thompson         break; // Caught by QF Assembly
10000d0321e0SJeremy L Thompson       }
10010d0321e0SJeremy L Thompson     }
10020d0321e0SJeremy L Thompson   }
10030d0321e0SJeremy L Thompson 
10040d0321e0SJeremy L Thompson   // Operator data struct
10050d0321e0SJeremy L Thompson   CeedOperator_Cuda *impl;
10060d0321e0SJeremy L Thompson   ierr = CeedOperatorGetData(op, &impl); CeedChkBackend(ierr);
10070d0321e0SJeremy L Thompson   ierr = CeedCalloc(1, &impl->diag); CeedChkBackend(ierr);
10080d0321e0SJeremy L Thompson   CeedOperatorDiag_Cuda *diag = impl->diag;
10090d0321e0SJeremy L Thompson   diag->basisin = basisin;
10100d0321e0SJeremy L Thompson   diag->basisout = basisout;
10110d0321e0SJeremy L Thompson   diag->h_emodein = emodein;
10120d0321e0SJeremy L Thompson   diag->h_emodeout = emodeout;
10130d0321e0SJeremy L Thompson   diag->numemodein = numemodein;
10140d0321e0SJeremy L Thompson   diag->numemodeout = numemodeout;
10150d0321e0SJeremy L Thompson 
10160d0321e0SJeremy L Thompson   // Assemble kernel
10170d0321e0SJeremy L Thompson   CeedInt nnodes, nqpts;
10180d0321e0SJeremy L Thompson   ierr = CeedBasisGetNumNodes(basisin, &nnodes); CeedChkBackend(ierr);
10190d0321e0SJeremy L Thompson   ierr = CeedBasisGetNumQuadraturePoints(basisin, &nqpts); CeedChkBackend(ierr);
10200d0321e0SJeremy L Thompson   diag->nnodes = nnodes;
10210d0321e0SJeremy L Thompson   ierr = CeedCompileCuda(ceed, diagonalkernels, &diag->module, 5,
10220d0321e0SJeremy L Thompson                          "NUMEMODEIN", numemodein,
10230d0321e0SJeremy L Thompson                          "NUMEMODEOUT", numemodeout,
10240d0321e0SJeremy L Thompson                          "NNODES", nnodes,
10250d0321e0SJeremy L Thompson                          "NQPTS", nqpts,
10260d0321e0SJeremy L Thompson                          "NCOMP", ncomp
10270d0321e0SJeremy L Thompson                         ); CeedChk_Cu(ceed, ierr);
10280d0321e0SJeremy L Thompson   ierr = CeedGetKernelCuda(ceed, diag->module, "linearDiagonal",
10290d0321e0SJeremy L Thompson                            &diag->linearDiagonal); CeedChk_Cu(ceed, ierr);
10300d0321e0SJeremy L Thompson   ierr = CeedGetKernelCuda(ceed, diag->module, "linearPointBlockDiagonal",
10310d0321e0SJeremy L Thompson                            &diag->linearPointBlock);
10320d0321e0SJeremy L Thompson   CeedChk_Cu(ceed, ierr);
10330d0321e0SJeremy L Thompson 
10340d0321e0SJeremy L Thompson   // Basis matrices
10350d0321e0SJeremy L Thompson   const CeedInt qBytes = nqpts * sizeof(CeedScalar);
10360d0321e0SJeremy L Thompson   const CeedInt iBytes = qBytes * nnodes;
10370d0321e0SJeremy L Thompson   const CeedInt gBytes = qBytes * nnodes * dim;
10380d0321e0SJeremy L Thompson   const CeedInt eBytes = sizeof(CeedEvalMode);
10390d0321e0SJeremy L Thompson   const CeedScalar *interpin, *interpout, *gradin, *gradout;
10400d0321e0SJeremy L Thompson 
10410d0321e0SJeremy L Thompson   // CEED_EVAL_NONE
10420d0321e0SJeremy L Thompson   CeedScalar *identity = NULL;
10430d0321e0SJeremy L Thompson   bool evalNone = false;
10440d0321e0SJeremy L Thompson   for (CeedInt i=0; i<numemodein; i++)
10450d0321e0SJeremy L Thompson     evalNone = evalNone || (emodein[i] == CEED_EVAL_NONE);
10460d0321e0SJeremy L Thompson   for (CeedInt i=0; i<numemodeout; i++)
10470d0321e0SJeremy L Thompson     evalNone = evalNone || (emodeout[i] == CEED_EVAL_NONE);
10480d0321e0SJeremy L Thompson   if (evalNone) {
10490d0321e0SJeremy L Thompson     ierr = CeedCalloc(nqpts*nnodes, &identity); CeedChkBackend(ierr);
10500d0321e0SJeremy L Thompson     for (CeedInt i=0; i<(nnodes<nqpts?nnodes:nqpts); i++)
10510d0321e0SJeremy L Thompson       identity[i*nnodes+i] = 1.0;
10520d0321e0SJeremy L Thompson     ierr = cudaMalloc((void **)&diag->d_identity, iBytes); CeedChk_Cu(ceed, ierr);
10530d0321e0SJeremy L Thompson     ierr = cudaMemcpy(diag->d_identity, identity, iBytes,
10540d0321e0SJeremy L Thompson                       cudaMemcpyHostToDevice); CeedChk_Cu(ceed, ierr);
10550d0321e0SJeremy L Thompson   }
10560d0321e0SJeremy L Thompson 
10570d0321e0SJeremy L Thompson   // CEED_EVAL_INTERP
10580d0321e0SJeremy L Thompson   ierr = CeedBasisGetInterp(basisin, &interpin); CeedChkBackend(ierr);
10590d0321e0SJeremy L Thompson   ierr = cudaMalloc((void **)&diag->d_interpin, iBytes); CeedChk_Cu(ceed, ierr);
10600d0321e0SJeremy L Thompson   ierr = cudaMemcpy(diag->d_interpin, interpin, iBytes,
10610d0321e0SJeremy L Thompson                     cudaMemcpyHostToDevice); CeedChk_Cu(ceed, ierr);
10620d0321e0SJeremy L Thompson   ierr = CeedBasisGetInterp(basisout, &interpout); CeedChkBackend(ierr);
10630d0321e0SJeremy L Thompson   ierr = cudaMalloc((void **)&diag->d_interpout, iBytes); CeedChk_Cu(ceed, ierr);
10640d0321e0SJeremy L Thompson   ierr = cudaMemcpy(diag->d_interpout, interpout, iBytes,
10650d0321e0SJeremy L Thompson                     cudaMemcpyHostToDevice); CeedChk_Cu(ceed, ierr);
10660d0321e0SJeremy L Thompson 
10670d0321e0SJeremy L Thompson   // CEED_EVAL_GRAD
10680d0321e0SJeremy L Thompson   ierr = CeedBasisGetGrad(basisin, &gradin); CeedChkBackend(ierr);
10690d0321e0SJeremy L Thompson   ierr = cudaMalloc((void **)&diag->d_gradin, gBytes); CeedChk_Cu(ceed, ierr);
10700d0321e0SJeremy L Thompson   ierr = cudaMemcpy(diag->d_gradin, gradin, gBytes,
10710d0321e0SJeremy L Thompson                     cudaMemcpyHostToDevice); CeedChk_Cu(ceed, ierr);
10720d0321e0SJeremy L Thompson   ierr = CeedBasisGetGrad(basisout, &gradout); CeedChkBackend(ierr);
10730d0321e0SJeremy L Thompson   ierr = cudaMalloc((void **)&diag->d_gradout, gBytes); CeedChk_Cu(ceed, ierr);
10740d0321e0SJeremy L Thompson   ierr = cudaMemcpy(diag->d_gradout, gradout, gBytes,
10750d0321e0SJeremy L Thompson                     cudaMemcpyHostToDevice); CeedChk_Cu(ceed, ierr);
10760d0321e0SJeremy L Thompson 
10770d0321e0SJeremy L Thompson   // Arrays of emodes
10780d0321e0SJeremy L Thompson   ierr = cudaMalloc((void **)&diag->d_emodein, numemodein * eBytes);
10790d0321e0SJeremy L Thompson   CeedChk_Cu(ceed, ierr);
10800d0321e0SJeremy L Thompson   ierr = cudaMemcpy(diag->d_emodein, emodein, numemodein * eBytes,
10810d0321e0SJeremy L Thompson                     cudaMemcpyHostToDevice); CeedChk_Cu(ceed, ierr);
10820d0321e0SJeremy L Thompson   ierr = cudaMalloc((void **)&diag->d_emodeout, numemodeout * eBytes);
10830d0321e0SJeremy L Thompson   CeedChk_Cu(ceed, ierr);
10840d0321e0SJeremy L Thompson   ierr = cudaMemcpy(diag->d_emodeout, emodeout, numemodeout * eBytes,
10850d0321e0SJeremy L Thompson                     cudaMemcpyHostToDevice); CeedChk_Cu(ceed, ierr);
10860d0321e0SJeremy L Thompson 
10870d0321e0SJeremy L Thompson   // Restriction
10880d0321e0SJeremy L Thompson   diag->diagrstr = rstrout;
10890d0321e0SJeremy L Thompson 
10900d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
10910d0321e0SJeremy L Thompson }
10920d0321e0SJeremy L Thompson 
10930d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
10940d0321e0SJeremy L Thompson // Assemble diagonal common code
10950d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
10960d0321e0SJeremy L Thompson static inline int CeedOperatorAssembleDiagonalCore_Cuda(CeedOperator op,
10970d0321e0SJeremy L Thompson     CeedVector assembled, CeedRequest *request, const bool pointBlock) {
10980d0321e0SJeremy L Thompson   int ierr;
10990d0321e0SJeremy L Thompson   Ceed ceed;
11000d0321e0SJeremy L Thompson   ierr = CeedOperatorGetCeed(op, &ceed); CeedChkBackend(ierr);
11010d0321e0SJeremy L Thompson   CeedOperator_Cuda *impl;
11020d0321e0SJeremy L Thompson   ierr = CeedOperatorGetData(op, &impl); CeedChkBackend(ierr);
11030d0321e0SJeremy L Thompson 
11040d0321e0SJeremy L Thompson   // Assemble QFunction
11050d0321e0SJeremy L Thompson   CeedVector assembledqf;
11060d0321e0SJeremy L Thompson   CeedElemRestriction rstr;
11070d0321e0SJeremy L Thompson   ierr = CeedOperatorLinearAssembleQFunctionBuildOrUpdate(op, &assembledqf,
11080d0321e0SJeremy L Thompson          &rstr, request); CeedChkBackend(ierr);
11090d0321e0SJeremy L Thompson   ierr = CeedElemRestrictionDestroy(&rstr); CeedChkBackend(ierr);
11100d0321e0SJeremy L Thompson   CeedScalar maxnorm = 0;
11110d0321e0SJeremy L Thompson   ierr = CeedVectorNorm(assembledqf, CEED_NORM_MAX, &maxnorm);
11120d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
11130d0321e0SJeremy L Thompson 
11140d0321e0SJeremy L Thompson   // Setup
11150d0321e0SJeremy L Thompson   if (!impl->diag) {
11160d0321e0SJeremy L Thompson     ierr = CeedOperatorAssembleDiagonalSetup_Cuda(op, pointBlock);
11170d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
11180d0321e0SJeremy L Thompson   }
11190d0321e0SJeremy L Thompson   CeedOperatorDiag_Cuda *diag = impl->diag;
11200d0321e0SJeremy L Thompson   assert(diag != NULL);
11210d0321e0SJeremy L Thompson 
11220d0321e0SJeremy L Thompson   // Restriction
11230d0321e0SJeremy L Thompson   if (pointBlock && !diag->pbdiagrstr) {
11240d0321e0SJeremy L Thompson     CeedElemRestriction pbdiagrstr;
11250d0321e0SJeremy L Thompson     ierr = CreatePBRestriction(diag->diagrstr, &pbdiagrstr); CeedChkBackend(ierr);
11260d0321e0SJeremy L Thompson     diag->pbdiagrstr = pbdiagrstr;
11270d0321e0SJeremy L Thompson   }
11280d0321e0SJeremy L Thompson   CeedElemRestriction diagrstr = pointBlock ? diag->pbdiagrstr : diag->diagrstr;
11290d0321e0SJeremy L Thompson 
11300d0321e0SJeremy L Thompson   // Create diagonal vector
11310d0321e0SJeremy L Thompson   CeedVector elemdiag = pointBlock ? diag->pbelemdiag : diag->elemdiag;
11320d0321e0SJeremy L Thompson   if (!elemdiag) {
11330d0321e0SJeremy L Thompson     ierr = CeedElemRestrictionCreateVector(diagrstr, NULL, &elemdiag);
11340d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
11350d0321e0SJeremy L Thompson     if (pointBlock)
11360d0321e0SJeremy L Thompson       diag->pbelemdiag = elemdiag;
11370d0321e0SJeremy L Thompson     else
11380d0321e0SJeremy L Thompson       diag->elemdiag = elemdiag;
11390d0321e0SJeremy L Thompson   }
11400d0321e0SJeremy L Thompson   ierr = CeedVectorSetValue(elemdiag, 0.0); CeedChkBackend(ierr);
11410d0321e0SJeremy L Thompson 
11420d0321e0SJeremy L Thompson   // Assemble element operator diagonals
11430d0321e0SJeremy L Thompson   CeedScalar *elemdiagarray;
11440d0321e0SJeremy L Thompson   const CeedScalar *assembledqfarray;
11450d0321e0SJeremy L Thompson   ierr = CeedVectorGetArray(elemdiag, CEED_MEM_DEVICE, &elemdiagarray);
11460d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
11470d0321e0SJeremy L Thompson   ierr = CeedVectorGetArrayRead(assembledqf, CEED_MEM_DEVICE, &assembledqfarray);
11480d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
11490d0321e0SJeremy L Thompson   CeedInt nelem;
11500d0321e0SJeremy L Thompson   ierr = CeedElemRestrictionGetNumElements(diagrstr, &nelem);
11510d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
11520d0321e0SJeremy L Thompson 
11530d0321e0SJeremy L Thompson   // Compute the diagonal of B^T D B
11540d0321e0SJeremy L Thompson   int elemsPerBlock = 1;
11550d0321e0SJeremy L Thompson   int grid = nelem/elemsPerBlock+((nelem/elemsPerBlock*elemsPerBlock<nelem)?1:0);
11560d0321e0SJeremy L Thompson   void *args[] = {(void *) &nelem, (void *) &maxnorm, &diag->d_identity,
11570d0321e0SJeremy L Thompson                   &diag->d_interpin, &diag->d_gradin, &diag->d_interpout,
11580d0321e0SJeremy L Thompson                   &diag->d_gradout, &diag->d_emodein, &diag->d_emodeout,
11590d0321e0SJeremy L Thompson                   &assembledqfarray, &elemdiagarray
11600d0321e0SJeremy L Thompson                  };
11610d0321e0SJeremy L Thompson   if (pointBlock) {
11620d0321e0SJeremy L Thompson     ierr = CeedRunKernelDimCuda(ceed, diag->linearPointBlock, grid,
11630d0321e0SJeremy L Thompson                                 diag->nnodes, 1, elemsPerBlock, args);
11640d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
11650d0321e0SJeremy L Thompson   } else {
11660d0321e0SJeremy L Thompson     ierr = CeedRunKernelDimCuda(ceed, diag->linearDiagonal, grid,
11670d0321e0SJeremy L Thompson                                 diag->nnodes, 1, elemsPerBlock, args);
11680d0321e0SJeremy L Thompson     CeedChkBackend(ierr);
11690d0321e0SJeremy L Thompson   }
11700d0321e0SJeremy L Thompson 
11710d0321e0SJeremy L Thompson   // Restore arrays
11720d0321e0SJeremy L Thompson   ierr = CeedVectorRestoreArray(elemdiag, &elemdiagarray); CeedChkBackend(ierr);
11730d0321e0SJeremy L Thompson   ierr = CeedVectorRestoreArrayRead(assembledqf, &assembledqfarray);
11740d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
11750d0321e0SJeremy L Thompson 
11760d0321e0SJeremy L Thompson   // Assemble local operator diagonal
11770d0321e0SJeremy L Thompson   ierr = CeedElemRestrictionApply(diagrstr, CEED_TRANSPOSE, elemdiag,
11780d0321e0SJeremy L Thompson                                   assembled, request); CeedChkBackend(ierr);
11790d0321e0SJeremy L Thompson 
11800d0321e0SJeremy L Thompson   // Cleanup
11810d0321e0SJeremy L Thompson   ierr = CeedVectorDestroy(&assembledqf); CeedChkBackend(ierr);
11820d0321e0SJeremy L Thompson 
11830d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
11840d0321e0SJeremy L Thompson }
11850d0321e0SJeremy L Thompson 
11860d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
11870d0321e0SJeremy L Thompson // Assemble composite diagonal common code
11880d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
11890d0321e0SJeremy L Thompson static inline int CeedOperatorLinearAssembleAddDiagonalCompositeCore_Cuda(
11900d0321e0SJeremy L Thompson   CeedOperator op, CeedVector assembled, CeedRequest *request,
11910d0321e0SJeremy L Thompson   const bool pointBlock) {
11920d0321e0SJeremy L Thompson   int ierr;
11930d0321e0SJeremy L Thompson   CeedInt numSub;
11940d0321e0SJeremy L Thompson   CeedOperator *subOperators;
11950d0321e0SJeremy L Thompson   ierr = CeedOperatorGetNumSub(op, &numSub); CeedChkBackend(ierr);
11960d0321e0SJeremy L Thompson   ierr = CeedOperatorGetSubList(op, &subOperators); CeedChkBackend(ierr);
11970d0321e0SJeremy L Thompson   for (CeedInt i = 0; i < numSub; i++) {
11980d0321e0SJeremy L Thompson     ierr = CeedOperatorAssembleDiagonalCore_Cuda(subOperators[i], assembled,
11990d0321e0SJeremy L Thompson            request, pointBlock); CeedChkBackend(ierr);
12000d0321e0SJeremy L Thompson   }
12010d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
12020d0321e0SJeremy L Thompson }
12030d0321e0SJeremy L Thompson 
12040d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
12050d0321e0SJeremy L Thompson // Assemble Linear Diagonal
12060d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
12070d0321e0SJeremy L Thompson static int CeedOperatorLinearAssembleAddDiagonal_Cuda(CeedOperator op,
12080d0321e0SJeremy L Thompson     CeedVector assembled, CeedRequest *request) {
12090d0321e0SJeremy L Thompson   int ierr;
12100d0321e0SJeremy L Thompson   bool isComposite;
12110d0321e0SJeremy L Thompson   ierr = CeedOperatorIsComposite(op, &isComposite); CeedChkBackend(ierr);
12120d0321e0SJeremy L Thompson   if (isComposite) {
12130d0321e0SJeremy L Thompson     return CeedOperatorLinearAssembleAddDiagonalCompositeCore_Cuda(op, assembled,
12140d0321e0SJeremy L Thompson            request, false);
12150d0321e0SJeremy L Thompson   } else {
12160d0321e0SJeremy L Thompson     return CeedOperatorAssembleDiagonalCore_Cuda(op, assembled, request, false);
12170d0321e0SJeremy L Thompson   }
12180d0321e0SJeremy L Thompson }
12190d0321e0SJeremy L Thompson 
12200d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
12210d0321e0SJeremy L Thompson // Assemble Linear Point Block Diagonal
12220d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
12230d0321e0SJeremy L Thompson static int CeedOperatorLinearAssembleAddPointBlockDiagonal_Cuda(CeedOperator op,
12240d0321e0SJeremy L Thompson     CeedVector assembled, CeedRequest *request) {
12250d0321e0SJeremy L Thompson   int ierr;
12260d0321e0SJeremy L Thompson   bool isComposite;
12270d0321e0SJeremy L Thompson   ierr = CeedOperatorIsComposite(op, &isComposite); CeedChkBackend(ierr);
12280d0321e0SJeremy L Thompson   if (isComposite) {
12290d0321e0SJeremy L Thompson     return CeedOperatorLinearAssembleAddDiagonalCompositeCore_Cuda(op, assembled,
12300d0321e0SJeremy L Thompson            request, true);
12310d0321e0SJeremy L Thompson   } else {
12320d0321e0SJeremy L Thompson     return CeedOperatorAssembleDiagonalCore_Cuda(op, assembled, request, true);
12330d0321e0SJeremy L Thompson   }
12340d0321e0SJeremy L Thompson }
12350d0321e0SJeremy L Thompson 
12360d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
12370d0321e0SJeremy L Thompson // Create operator
12380d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
12390d0321e0SJeremy L Thompson int CeedOperatorCreate_Cuda(CeedOperator op) {
12400d0321e0SJeremy L Thompson   int ierr;
12410d0321e0SJeremy L Thompson   Ceed ceed;
12420d0321e0SJeremy L Thompson   ierr = CeedOperatorGetCeed(op, &ceed); CeedChkBackend(ierr);
12430d0321e0SJeremy L Thompson   CeedOperator_Cuda *impl;
12440d0321e0SJeremy L Thompson 
12450d0321e0SJeremy L Thompson   ierr = CeedCalloc(1, &impl); CeedChkBackend(ierr);
12460d0321e0SJeremy L Thompson   ierr = CeedOperatorSetData(op, impl); CeedChkBackend(ierr);
12470d0321e0SJeremy L Thompson 
12480d0321e0SJeremy L Thompson   ierr = CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleQFunction",
12490d0321e0SJeremy L Thompson                                 CeedOperatorLinearAssembleQFunction_Cuda);
12500d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
12510d0321e0SJeremy L Thompson   ierr = CeedSetBackendFunction(ceed, "Operator", op,
12520d0321e0SJeremy L Thompson                                 "LinearAssembleQFunctionUpdate",
12530d0321e0SJeremy L Thompson                                 CeedOperatorLinearAssembleQFunctionUpdate_Cuda);
12540d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
12550d0321e0SJeremy L Thompson   ierr = CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleAddDiagonal",
12560d0321e0SJeremy L Thompson                                 CeedOperatorLinearAssembleAddDiagonal_Cuda);
12570d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
12580d0321e0SJeremy L Thompson   ierr = CeedSetBackendFunction(ceed, "Operator", op,
12590d0321e0SJeremy L Thompson                                 "LinearAssembleAddPointBlockDiagonal",
12600d0321e0SJeremy L Thompson                                 CeedOperatorLinearAssembleAddPointBlockDiagonal_Cuda);
12610d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
12620d0321e0SJeremy L Thompson   ierr = CeedSetBackendFunction(ceed, "Operator", op, "ApplyAdd",
12630d0321e0SJeremy L Thompson                                 CeedOperatorApplyAdd_Cuda); CeedChkBackend(ierr);
12640d0321e0SJeremy L Thompson   ierr = CeedSetBackendFunction(ceed, "Operator", op, "Destroy",
12650d0321e0SJeremy L Thompson                                 CeedOperatorDestroy_Cuda); CeedChkBackend(ierr);
12660d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
12670d0321e0SJeremy L Thompson }
12680d0321e0SJeremy L Thompson 
12690d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
12700d0321e0SJeremy L Thompson // Composite Operator Create
12710d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
12720d0321e0SJeremy L Thompson int CeedCompositeOperatorCreate_Cuda(CeedOperator op) {
12730d0321e0SJeremy L Thompson   int ierr;
12740d0321e0SJeremy L Thompson   Ceed ceed;
12750d0321e0SJeremy L Thompson   ierr = CeedOperatorGetCeed(op, &ceed); CeedChkBackend(ierr);
12760d0321e0SJeremy L Thompson 
12770d0321e0SJeremy L Thompson   ierr = CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleAddDiagonal",
12780d0321e0SJeremy L Thompson                                 CeedOperatorLinearAssembleAddDiagonal_Cuda);
12790d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
12800d0321e0SJeremy L Thompson   ierr = CeedSetBackendFunction(ceed, "Operator", op,
12810d0321e0SJeremy L Thompson                                 "LinearAssembleAddPointBlockDiagonal",
12820d0321e0SJeremy L Thompson                                 CeedOperatorLinearAssembleAddPointBlockDiagonal_Cuda);
12830d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
12840d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
12850d0321e0SJeremy L Thompson }
12860d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
1287