xref: /libCEED/backends/cuda-ref/ceed-cuda-ref-operator.c (revision ca7355308a14805110a7d98dabf1f658d5ec16d5)
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 
849aac155SJeremy L Thompson #include <ceed.h>
92b730f8bSJeremy L Thompson #include <ceed/backend.h>
102b730f8bSJeremy L Thompson #include <ceed/jit-tools.h>
11c85e8640SSebastian Grimberg #include <assert.h>
120d0321e0SJeremy L Thompson #include <cuda.h>
130d0321e0SJeremy L Thompson #include <cuda_runtime.h>
140d0321e0SJeremy L Thompson #include <stdbool.h>
150d0321e0SJeremy L Thompson #include <string.h>
162b730f8bSJeremy L Thompson 
1749aac155SJeremy L Thompson #include "../cuda/ceed-cuda-common.h"
180d0321e0SJeremy L Thompson #include "../cuda/ceed-cuda-compile.h"
192b730f8bSJeremy L Thompson #include "ceed-cuda-ref.h"
200d0321e0SJeremy L Thompson 
210d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
220d0321e0SJeremy L Thompson // Destroy operator
230d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
240d0321e0SJeremy L Thompson static int CeedOperatorDestroy_Cuda(CeedOperator op) {
250d0321e0SJeremy L Thompson   CeedOperator_Cuda *impl;
26*ca735530SJeremy L Thompson 
272b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &impl));
280d0321e0SJeremy L Thompson 
290d0321e0SJeremy L Thompson   // Apply data
30*ca735530SJeremy L Thompson   for (CeedInt i = 0; i < impl->num_inputs + impl->num_outputs; i++) {
31*ca735530SJeremy L Thompson     CeedCallBackend(CeedVectorDestroy(&impl->e_vecs[i]));
320d0321e0SJeremy L Thompson   }
33*ca735530SJeremy L Thompson   CeedCallBackend(CeedFree(&impl->e_vecs));
340d0321e0SJeremy L Thompson 
35*ca735530SJeremy L Thompson   for (CeedInt i = 0; i < impl->num_inputs; i++) {
36*ca735530SJeremy L Thompson     CeedCallBackend(CeedVectorDestroy(&impl->q_vecs_in[i]));
370d0321e0SJeremy L Thompson   }
38*ca735530SJeremy L Thompson   CeedCallBackend(CeedFree(&impl->q_vecs_in));
390d0321e0SJeremy L Thompson 
40*ca735530SJeremy L Thompson   for (CeedInt i = 0; i < impl->num_outputs; i++) {
41*ca735530SJeremy L Thompson     CeedCallBackend(CeedVectorDestroy(&impl->q_vecs_out[i]));
420d0321e0SJeremy L Thompson   }
43*ca735530SJeremy L Thompson   CeedCallBackend(CeedFree(&impl->q_vecs_out));
440d0321e0SJeremy L Thompson 
450d0321e0SJeremy L Thompson   // QFunction assembly data
46*ca735530SJeremy L Thompson   for (CeedInt i = 0; i < impl->num_active_in; i++) {
47*ca735530SJeremy L Thompson     CeedCallBackend(CeedVectorDestroy(&impl->qf_active_in[i]));
480d0321e0SJeremy L Thompson   }
49*ca735530SJeremy L Thompson   CeedCallBackend(CeedFree(&impl->qf_active_in));
500d0321e0SJeremy L Thompson 
510d0321e0SJeremy L Thompson   // Diag data
520d0321e0SJeremy L Thompson   if (impl->diag) {
530d0321e0SJeremy L Thompson     Ceed ceed;
54*ca735530SJeremy L Thompson 
552b730f8bSJeremy L Thompson     CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
562b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cuModuleUnload(impl->diag->module));
57*ca735530SJeremy L Thompson     CeedCallBackend(CeedFree(&impl->diag->h_e_mode_in));
58*ca735530SJeremy L Thompson     CeedCallBackend(CeedFree(&impl->diag->h_e_mode_out));
59*ca735530SJeremy L Thompson     CeedCallCuda(ceed, cudaFree(impl->diag->d_e_mode_in));
60*ca735530SJeremy L Thompson     CeedCallCuda(ceed, cudaFree(impl->diag->d_e_mode_out));
612b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaFree(impl->diag->d_identity));
62*ca735530SJeremy L Thompson     CeedCallCuda(ceed, cudaFree(impl->diag->d_interp_in));
63*ca735530SJeremy L Thompson     CeedCallCuda(ceed, cudaFree(impl->diag->d_interp_out));
64*ca735530SJeremy L Thompson     CeedCallCuda(ceed, cudaFree(impl->diag->d_grad_in));
65*ca735530SJeremy L Thompson     CeedCallCuda(ceed, cudaFree(impl->diag->d_grad_out));
66*ca735530SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionDestroy(&impl->diag->point_block_rstr));
67*ca735530SJeremy L Thompson     CeedCallBackend(CeedVectorDestroy(&impl->diag->elem_diag));
68*ca735530SJeremy L Thompson     CeedCallBackend(CeedVectorDestroy(&impl->diag->point_block_elem_diag));
690d0321e0SJeremy L Thompson   }
702b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&impl->diag));
710d0321e0SJeremy L Thompson 
72cc132f9aSnbeams   if (impl->asmb) {
73cc132f9aSnbeams     Ceed ceed;
74*ca735530SJeremy L Thompson 
752b730f8bSJeremy L Thompson     CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
762b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cuModuleUnload(impl->asmb->module));
772b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaFree(impl->asmb->d_B_in));
782b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaFree(impl->asmb->d_B_out));
79cc132f9aSnbeams   }
802b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&impl->asmb));
81cc132f9aSnbeams 
822b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&impl));
830d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
840d0321e0SJeremy L Thompson }
850d0321e0SJeremy L Thompson 
860d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
870d0321e0SJeremy L Thompson // Setup infields or outfields
880d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
89*ca735530SJeremy L Thompson static int CeedOperatorSetupFields_Cuda(CeedQFunction qf, CeedOperator op, bool is_input, CeedVector *e_vecs, CeedVector *q_vecs, CeedInt e_start,
90*ca735530SJeremy L Thompson                                         CeedInt num_fields, CeedInt Q, CeedInt num_elem) {
910d0321e0SJeremy L Thompson   Ceed                ceed;
92*ca735530SJeremy L Thompson   bool                is_strided, skip_restriction;
93*ca735530SJeremy L Thompson   CeedSize            q_size;
94*ca735530SJeremy L Thompson   CeedInt             dim, size;
95*ca735530SJeremy L Thompson   CeedQFunctionField *qf_fields;
96*ca735530SJeremy L Thompson   CeedOperatorField  *op_fields;
970d0321e0SJeremy L Thompson 
98*ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
99*ca735530SJeremy L Thompson 
100*ca735530SJeremy L Thompson   if (is_input) {
101*ca735530SJeremy L Thompson     CeedCallBackend(CeedOperatorGetFields(op, NULL, &op_fields, NULL, NULL));
102*ca735530SJeremy L Thompson     CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_fields, NULL, NULL));
1030d0321e0SJeremy L Thompson   } else {
104*ca735530SJeremy L Thompson     CeedCallBackend(CeedOperatorGetFields(op, NULL, NULL, NULL, &op_fields));
105*ca735530SJeremy L Thompson     CeedCallBackend(CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qf_fields));
1060d0321e0SJeremy L Thompson   }
1070d0321e0SJeremy L Thompson 
1080d0321e0SJeremy L Thompson   // Loop over fields
109*ca735530SJeremy L Thompson   for (CeedInt i = 0; i < num_fields; i++) {
110*ca735530SJeremy L Thompson     CeedEvalMode e_mode;
111*ca735530SJeremy L Thompson     CeedBasis    basis;
1120d0321e0SJeremy L Thompson 
113*ca735530SJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &e_mode));
114*ca735530SJeremy L Thompson 
115*ca735530SJeremy L Thompson     is_strided       = false;
116*ca735530SJeremy L Thompson     skip_restriction = false;
117*ca735530SJeremy L Thompson     if (e_mode != CEED_EVAL_WEIGHT) {
118*ca735530SJeremy L Thompson       CeedElemRestriction elem_restr;
119*ca735530SJeremy L Thompson 
120*ca735530SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_fields[i], &elem_restr));
1210d0321e0SJeremy L Thompson 
1220d0321e0SJeremy L Thompson       // Check whether this field can skip the element restriction:
123*ca735530SJeremy L Thompson       // must be passive input, with e_mode NONE, and have a strided restriction with CEED_STRIDES_BACKEND.
1240d0321e0SJeremy L Thompson 
1250d0321e0SJeremy L Thompson       // First, check whether the field is input or output:
126*ca735530SJeremy L Thompson       if (is_input) {
127*ca735530SJeremy L Thompson         CeedVector vec;
128*ca735530SJeremy L Thompson 
1290d0321e0SJeremy L Thompson         // Check for passive input:
130*ca735530SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec));
131*ca735530SJeremy L Thompson         if (vec != CEED_VECTOR_ACTIVE) {
132*ca735530SJeremy L Thompson           // Check e_mode
133*ca735530SJeremy L Thompson           if (e_mode == CEED_EVAL_NONE) {
1340d0321e0SJeremy L Thompson             // Check for strided restriction
135*ca735530SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionIsStrided(elem_restr, &is_strided));
136*ca735530SJeremy L Thompson             if (is_strided) {
1370d0321e0SJeremy L Thompson               // Check if vector is already in preferred backend ordering
138*ca735530SJeremy L Thompson               CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_restr, &skip_restriction));
1390d0321e0SJeremy L Thompson             }
1400d0321e0SJeremy L Thompson           }
1410d0321e0SJeremy L Thompson         }
1420d0321e0SJeremy L Thompson       }
143*ca735530SJeremy L Thompson       if (skip_restriction) {
144ea61e9acSJeremy L Thompson         // We do not need an E-Vector, but will use the input field vector's data directly in the operator application.
145*ca735530SJeremy L Thompson         e_vecs[i + e_start] = NULL;
1460d0321e0SJeremy L Thompson       } else {
147*ca735530SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionCreateVector(elem_restr, NULL, &e_vecs[i + e_start]));
1480d0321e0SJeremy L Thompson       }
1490d0321e0SJeremy L Thompson     }
1500d0321e0SJeremy L Thompson 
151*ca735530SJeremy L Thompson     switch (e_mode) {
1520d0321e0SJeremy L Thompson       case CEED_EVAL_NONE:
153*ca735530SJeremy L Thompson         CeedCallBackend(CeedQFunctionFieldGetSize(qf_fields[i], &size));
154*ca735530SJeremy L Thompson         q_size = (CeedSize)num_elem * Q * size;
155*ca735530SJeremy L Thompson         CeedCallBackend(CeedVectorCreate(ceed, q_size, &q_vecs[i]));
1560d0321e0SJeremy L Thompson         break;
1570d0321e0SJeremy L Thompson       case CEED_EVAL_INTERP:
158*ca735530SJeremy L Thompson         CeedCallBackend(CeedQFunctionFieldGetSize(qf_fields[i], &size));
159*ca735530SJeremy L Thompson         q_size = (CeedSize)num_elem * Q * size;
160*ca735530SJeremy L Thompson         CeedCallBackend(CeedVectorCreate(ceed, q_size, &q_vecs[i]));
1610d0321e0SJeremy L Thompson         break;
1620d0321e0SJeremy L Thompson       case CEED_EVAL_GRAD:
163*ca735530SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis));
164*ca735530SJeremy L Thompson         CeedCallBackend(CeedQFunctionFieldGetSize(qf_fields[i], &size));
1652b730f8bSJeremy L Thompson         CeedCallBackend(CeedBasisGetDimension(basis, &dim));
166*ca735530SJeremy L Thompson         q_size = (CeedSize)num_elem * Q * size;
167*ca735530SJeremy L Thompson         CeedCallBackend(CeedVectorCreate(ceed, q_size, &q_vecs[i]));
1680d0321e0SJeremy L Thompson         break;
1690d0321e0SJeremy L Thompson       case CEED_EVAL_WEIGHT:  // Only on input fields
170*ca735530SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis));
171*ca735530SJeremy L Thompson         q_size = (CeedSize)num_elem * Q;
172*ca735530SJeremy L Thompson         CeedCallBackend(CeedVectorCreate(ceed, q_size, &q_vecs[i]));
173*ca735530SJeremy L Thompson         CeedCallBackend(CeedBasisApply(basis, num_elem, CEED_NOTRANSPOSE, CEED_EVAL_WEIGHT, CEED_VECTOR_NONE, q_vecs[i]));
1740d0321e0SJeremy L Thompson         break;
1750d0321e0SJeremy L Thompson       case CEED_EVAL_DIV:
1760d0321e0SJeremy L Thompson         break;  // TODO: Not implemented
1770d0321e0SJeremy L Thompson       case CEED_EVAL_CURL:
1780d0321e0SJeremy L Thompson         break;  // TODO: Not implemented
1790d0321e0SJeremy L Thompson     }
1800d0321e0SJeremy L Thompson   }
1810d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1820d0321e0SJeremy L Thompson }
1830d0321e0SJeremy L Thompson 
1840d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
185ea61e9acSJeremy L Thompson // CeedOperator needs to connect all the named fields (be they active or passive) to the named inputs and outputs of its CeedQFunction.
1860d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
1870d0321e0SJeremy L Thompson static int CeedOperatorSetup_Cuda(CeedOperator op) {
1880d0321e0SJeremy L Thompson   Ceed                ceed;
189*ca735530SJeremy L Thompson   bool                is_setup_done;
190*ca735530SJeremy L Thompson   CeedInt             Q, num_elem, num_input_fields, num_output_fields;
191*ca735530SJeremy L Thompson   CeedQFunctionField *qf_input_fields, *qf_output_fields;
1920d0321e0SJeremy L Thompson   CeedQFunction       qf;
193*ca735530SJeremy L Thompson   CeedOperatorField  *op_input_fields, *op_output_fields;
194*ca735530SJeremy L Thompson   CeedOperator_Cuda  *impl;
195*ca735530SJeremy L Thompson 
196*ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorIsSetupDone(op, &is_setup_done));
197*ca735530SJeremy L Thompson   if (is_setup_done) return CEED_ERROR_SUCCESS;
198*ca735530SJeremy L Thompson 
199*ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
200*ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &impl));
2012b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
2022b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q));
203*ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem));
204*ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
205*ca735530SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
2060d0321e0SJeremy L Thompson 
2070d0321e0SJeremy L Thompson   // Allocate
208*ca735530SJeremy L Thompson   CeedCallBackend(CeedCalloc(num_input_fields + num_output_fields, &impl->e_vecs));
2090d0321e0SJeremy L Thompson 
210*ca735530SJeremy L Thompson   CeedCallBackend(CeedCalloc(CEED_FIELD_MAX, &impl->q_vecs_in));
211*ca735530SJeremy L Thompson   CeedCallBackend(CeedCalloc(CEED_FIELD_MAX, &impl->q_vecs_out));
2120d0321e0SJeremy L Thompson 
213*ca735530SJeremy L Thompson   impl->num_inputs  = num_input_fields;
214*ca735530SJeremy L Thompson   impl->num_outputs = num_output_fields;
2150d0321e0SJeremy L Thompson 
216*ca735530SJeremy L Thompson   // Set up infield and outfield e_vecs and q_vecs
2170d0321e0SJeremy L Thompson   // Infields
218*ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorSetupFields_Cuda(qf, op, true, impl->e_vecs, impl->q_vecs_in, 0, num_input_fields, Q, num_elem));
2190d0321e0SJeremy L Thompson   // Outfields
220*ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorSetupFields_Cuda(qf, op, false, impl->e_vecs, impl->q_vecs_out, num_input_fields, num_output_fields, Q, num_elem));
2210d0321e0SJeremy L Thompson 
2222b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorSetSetupDone(op));
2230d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2240d0321e0SJeremy L Thompson }
2250d0321e0SJeremy L Thompson 
2260d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
2270d0321e0SJeremy L Thompson // Setup Operator Inputs
2280d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
229*ca735530SJeremy L Thompson static inline int CeedOperatorSetupInputs_Cuda(CeedInt num_input_fields, CeedQFunctionField *qf_input_fields, CeedOperatorField *op_input_fields,
230*ca735530SJeremy L Thompson                                                CeedVector in_vec, const bool skip_active_in, CeedScalar *e_data[2 * CEED_FIELD_MAX],
2310d0321e0SJeremy L Thompson                                                CeedOperator_Cuda *impl, CeedRequest *request) {
232*ca735530SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
233*ca735530SJeremy L Thompson     CeedEvalMode        e_mode;
2340d0321e0SJeremy L Thompson     CeedVector          vec;
235*ca735530SJeremy L Thompson     CeedElemRestriction elem_restr;
2360d0321e0SJeremy L Thompson 
2370d0321e0SJeremy L Thompson     // Get input vector
238*ca735530SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec));
2390d0321e0SJeremy L Thompson     if (vec == CEED_VECTOR_ACTIVE) {
240*ca735530SJeremy L Thompson       if (skip_active_in) continue;
241*ca735530SJeremy L Thompson       else vec = in_vec;
2420d0321e0SJeremy L Thompson     }
2430d0321e0SJeremy L Thompson 
244*ca735530SJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &e_mode));
245*ca735530SJeremy L Thompson     if (e_mode == CEED_EVAL_WEIGHT) {  // Skip
2460d0321e0SJeremy L Thompson     } else {
2470d0321e0SJeremy L Thompson       // Get input element restriction
248*ca735530SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_restr));
249*ca735530SJeremy L Thompson       if (vec == CEED_VECTOR_ACTIVE) vec = in_vec;
2500d0321e0SJeremy L Thompson       // Restrict, if necessary
251*ca735530SJeremy L Thompson       if (!impl->e_vecs[i]) {
2520d0321e0SJeremy L Thompson         // No restriction for this field; read data directly from vec.
253*ca735530SJeremy L Thompson         CeedCallBackend(CeedVectorGetArrayRead(vec, CEED_MEM_DEVICE, (const CeedScalar **)&e_data[i]));
2540d0321e0SJeremy L Thompson       } else {
255*ca735530SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionApply(elem_restr, CEED_NOTRANSPOSE, vec, impl->e_vecs[i], request));
2560d0321e0SJeremy L Thompson         // Get evec
257*ca735530SJeremy L Thompson         CeedCallBackend(CeedVectorGetArrayRead(impl->e_vecs[i], CEED_MEM_DEVICE, (const CeedScalar **)&e_data[i]));
2580d0321e0SJeremy L Thompson       }
2590d0321e0SJeremy L Thompson     }
2600d0321e0SJeremy L Thompson   }
2610d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2620d0321e0SJeremy L Thompson }
2630d0321e0SJeremy L Thompson 
2640d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
2650d0321e0SJeremy L Thompson // Input Basis Action
2660d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
267*ca735530SJeremy L Thompson static inline int CeedOperatorInputBasis_Cuda(CeedInt num_elem, CeedQFunctionField *qf_input_fields, CeedOperatorField *op_input_fields,
268*ca735530SJeremy L Thompson                                               CeedInt num_input_fields, const bool skip_active_in, CeedScalar *e_data[2 * CEED_FIELD_MAX],
2692b730f8bSJeremy L Thompson                                               CeedOperator_Cuda *impl) {
270*ca735530SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
271*ca735530SJeremy L Thompson     CeedInt             elem_size, size;
272*ca735530SJeremy L Thompson     CeedEvalMode        e_mode;
273*ca735530SJeremy L Thompson     CeedElemRestriction elem_restr;
2740d0321e0SJeremy L Thompson     CeedBasis           basis;
2750d0321e0SJeremy L Thompson 
2760d0321e0SJeremy L Thompson     // Skip active input
277*ca735530SJeremy L Thompson     if (skip_active_in) {
2780d0321e0SJeremy L Thompson       CeedVector vec;
279*ca735530SJeremy L Thompson 
280*ca735530SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec));
2812b730f8bSJeremy L Thompson       if (vec == CEED_VECTOR_ACTIVE) continue;
2820d0321e0SJeremy L Thompson     }
283*ca735530SJeremy L Thompson     // Get elem_size, e_mode, size
284*ca735530SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_restr));
285*ca735530SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_restr, &elem_size));
286*ca735530SJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &e_mode));
287*ca735530SJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetSize(qf_input_fields[i], &size));
2880d0321e0SJeremy L Thompson     // Basis action
289*ca735530SJeremy L Thompson     switch (e_mode) {
2900d0321e0SJeremy L Thompson       case CEED_EVAL_NONE:
291*ca735530SJeremy L Thompson         CeedCallBackend(CeedVectorSetArray(impl->q_vecs_in[i], CEED_MEM_DEVICE, CEED_USE_POINTER, e_data[i]));
2920d0321e0SJeremy L Thompson         break;
2930d0321e0SJeremy L Thompson       case CEED_EVAL_INTERP:
294*ca735530SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
295*ca735530SJeremy L Thompson         CeedCallBackend(CeedBasisApply(basis, num_elem, CEED_NOTRANSPOSE, CEED_EVAL_INTERP, impl->e_vecs[i], impl->q_vecs_in[i]));
2960d0321e0SJeremy L Thompson         break;
2970d0321e0SJeremy L Thompson       case CEED_EVAL_GRAD:
298*ca735530SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
299*ca735530SJeremy L Thompson         CeedCallBackend(CeedBasisApply(basis, num_elem, CEED_NOTRANSPOSE, CEED_EVAL_GRAD, impl->e_vecs[i], impl->q_vecs_in[i]));
3000d0321e0SJeremy L Thompson         break;
3010d0321e0SJeremy L Thompson       case CEED_EVAL_WEIGHT:
3020d0321e0SJeremy L Thompson         break;  // No action
3030d0321e0SJeremy L Thompson       case CEED_EVAL_DIV:
3040d0321e0SJeremy L Thompson         break;  // TODO: Not implemented
3050d0321e0SJeremy L Thompson       case CEED_EVAL_CURL:
3060d0321e0SJeremy L Thompson         break;  // TODO: Not implemented
3070d0321e0SJeremy L Thompson     }
3080d0321e0SJeremy L Thompson   }
3090d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3100d0321e0SJeremy L Thompson }
3110d0321e0SJeremy L Thompson 
3120d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
3130d0321e0SJeremy L Thompson // Restore Input Vectors
3140d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
315*ca735530SJeremy L Thompson static inline int CeedOperatorRestoreInputs_Cuda(CeedInt num_input_fields, CeedQFunctionField *qf_input_fields, CeedOperatorField *op_input_fields,
316*ca735530SJeremy L Thompson                                                  const bool skip_active_in, CeedScalar *e_data[2 * CEED_FIELD_MAX], CeedOperator_Cuda *impl) {
317*ca735530SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
318*ca735530SJeremy L Thompson     CeedEvalMode e_mode;
3190d0321e0SJeremy L Thompson     CeedVector   vec;
3200d0321e0SJeremy L Thompson 
3210d0321e0SJeremy L Thompson     // Skip active input
322*ca735530SJeremy L Thompson     if (skip_active_in) {
323*ca735530SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec));
3242b730f8bSJeremy L Thompson       if (vec == CEED_VECTOR_ACTIVE) continue;
3250d0321e0SJeremy L Thompson     }
326*ca735530SJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &e_mode));
327*ca735530SJeremy L Thompson     if (e_mode == CEED_EVAL_WEIGHT) {  // Skip
3280d0321e0SJeremy L Thompson     } else {
329*ca735530SJeremy L Thompson       if (!impl->e_vecs[i]) {  // This was a skip_restriction case
330*ca735530SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec));
331*ca735530SJeremy L Thompson         CeedCallBackend(CeedVectorRestoreArrayRead(vec, (const CeedScalar **)&e_data[i]));
3320d0321e0SJeremy L Thompson       } else {
333*ca735530SJeremy L Thompson         CeedCallBackend(CeedVectorRestoreArrayRead(impl->e_vecs[i], (const CeedScalar **)&e_data[i]));
3340d0321e0SJeremy L Thompson       }
3350d0321e0SJeremy L Thompson     }
3360d0321e0SJeremy L Thompson   }
3370d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3380d0321e0SJeremy L Thompson }
3390d0321e0SJeremy L Thompson 
3400d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
3410d0321e0SJeremy L Thompson // Apply and add to output
3420d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
343*ca735530SJeremy L Thompson static int CeedOperatorApplyAdd_Cuda(CeedOperator op, CeedVector in_vec, CeedVector out_vec, CeedRequest *request) {
3440d0321e0SJeremy L Thompson   CeedOperator_Cuda  *impl;
345*ca735530SJeremy L Thompson   CeedInt             Q, num_elem, elem_size, num_input_fields, num_output_fields, size;
346*ca735530SJeremy L Thompson   CeedEvalMode        e_mode;
347*ca735530SJeremy L Thompson   CeedScalar         *e_data[2 * CEED_FIELD_MAX] = {NULL};
348*ca735530SJeremy L Thompson   CeedOperatorField  *op_input_fields, *op_output_fields;
349*ca735530SJeremy L Thompson   CeedQFunctionField *qf_input_fields, *qf_output_fields;
3500d0321e0SJeremy L Thompson   CeedQFunction       qf;
351*ca735530SJeremy L Thompson 
352*ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &impl));
3532b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
3542b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q));
355*ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem));
356*ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
357*ca735530SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
3580d0321e0SJeremy L Thompson 
3590d0321e0SJeremy L Thompson   // Setup
3602b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorSetup_Cuda(op));
3610d0321e0SJeremy L Thompson 
362*ca735530SJeremy L Thompson   // Input e_vecs and Restriction
363*ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorSetupInputs_Cuda(num_input_fields, qf_input_fields, op_input_fields, in_vec, false, e_data, impl, request));
3640d0321e0SJeremy L Thompson 
3650d0321e0SJeremy L Thompson   // Input basis apply if needed
366*ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorInputBasis_Cuda(num_elem, qf_input_fields, op_input_fields, num_input_fields, false, e_data, impl));
3670d0321e0SJeremy L Thompson 
3680d0321e0SJeremy L Thompson   // Output pointers, as necessary
369*ca735530SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
370*ca735530SJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &e_mode));
371*ca735530SJeremy L Thompson     if (e_mode == CEED_EVAL_NONE) {
3720d0321e0SJeremy L Thompson       // Set the output Q-Vector to use the E-Vector data directly.
373*ca735530SJeremy L Thompson       CeedCallBackend(CeedVectorGetArrayWrite(impl->e_vecs[i + impl->num_inputs], CEED_MEM_DEVICE, &e_data[i + num_input_fields]));
374*ca735530SJeremy L Thompson       CeedCallBackend(CeedVectorSetArray(impl->q_vecs_out[i], CEED_MEM_DEVICE, CEED_USE_POINTER, e_data[i + num_input_fields]));
3750d0321e0SJeremy L Thompson     }
3760d0321e0SJeremy L Thompson   }
3770d0321e0SJeremy L Thompson 
3780d0321e0SJeremy L Thompson   // Q function
379*ca735530SJeremy L Thompson   CeedCallBackend(CeedQFunctionApply(qf, num_elem * Q, impl->q_vecs_in, impl->q_vecs_out));
3800d0321e0SJeremy L Thompson 
3810d0321e0SJeremy L Thompson   // Output basis apply if needed
382*ca735530SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
383*ca735530SJeremy L Thompson     CeedElemRestriction elem_restr;
384*ca735530SJeremy L Thompson     CeedBasis           basis;
385*ca735530SJeremy L Thompson 
386*ca735530SJeremy L Thompson     // Get elem_size, e_mode, size
387*ca735530SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &elem_restr));
388*ca735530SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_restr, &elem_size));
389*ca735530SJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &e_mode));
390*ca735530SJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetSize(qf_output_fields[i], &size));
3910d0321e0SJeremy L Thompson     // Basis action
392*ca735530SJeremy L Thompson     switch (e_mode) {
3930d0321e0SJeremy L Thompson       case CEED_EVAL_NONE:
3940d0321e0SJeremy L Thompson         break;
3950d0321e0SJeremy L Thompson       case CEED_EVAL_INTERP:
396*ca735530SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
397*ca735530SJeremy L Thompson         CeedCallBackend(CeedBasisApply(basis, num_elem, CEED_TRANSPOSE, CEED_EVAL_INTERP, impl->q_vecs_out[i], impl->e_vecs[i + impl->num_inputs]));
3980d0321e0SJeremy L Thompson         break;
3990d0321e0SJeremy L Thompson       case CEED_EVAL_GRAD:
400*ca735530SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
401*ca735530SJeremy L Thompson         CeedCallBackend(CeedBasisApply(basis, num_elem, CEED_TRANSPOSE, CEED_EVAL_GRAD, impl->q_vecs_out[i], impl->e_vecs[i + impl->num_inputs]));
4020d0321e0SJeremy L Thompson         break;
4030d0321e0SJeremy L Thompson       // LCOV_EXCL_START
4040d0321e0SJeremy L Thompson       case CEED_EVAL_WEIGHT: {
4050d0321e0SJeremy L Thompson         Ceed ceed;
4062b730f8bSJeremy L Thompson         CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
4072b730f8bSJeremy L Thompson         return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_WEIGHT cannot be an output evaluation mode");
4080d0321e0SJeremy L Thompson         break;  // Should not occur
4090d0321e0SJeremy L Thompson       }
4100d0321e0SJeremy L Thompson       case CEED_EVAL_DIV:
4110d0321e0SJeremy L Thompson         break;  // TODO: Not implemented
4120d0321e0SJeremy L Thompson       case CEED_EVAL_CURL:
4130d0321e0SJeremy L Thompson         break;  // TODO: Not implemented
4140d0321e0SJeremy L Thompson                 // LCOV_EXCL_STOP
4150d0321e0SJeremy L Thompson     }
4160d0321e0SJeremy L Thompson   }
4170d0321e0SJeremy L Thompson 
4180d0321e0SJeremy L Thompson   // Output restriction
419*ca735530SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
420*ca735530SJeremy L Thompson     CeedVector          vec;
421*ca735530SJeremy L Thompson     CeedElemRestriction elem_restr;
422*ca735530SJeremy L Thompson 
4230d0321e0SJeremy L Thompson     // Restore evec
424*ca735530SJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &e_mode));
425*ca735530SJeremy L Thompson     if (e_mode == CEED_EVAL_NONE) {
426*ca735530SJeremy L Thompson       CeedCallBackend(CeedVectorRestoreArray(impl->e_vecs[i + impl->num_inputs], &e_data[i + num_input_fields]));
4270d0321e0SJeremy L Thompson     }
4280d0321e0SJeremy L Thompson     // Get output vector
429*ca735530SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[i], &vec));
4300d0321e0SJeremy L Thompson     // Restrict
431*ca735530SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &elem_restr));
4320d0321e0SJeremy L Thompson     // Active
433*ca735530SJeremy L Thompson     if (vec == CEED_VECTOR_ACTIVE) vec = out_vec;
4340d0321e0SJeremy L Thompson 
435*ca735530SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionApply(elem_restr, CEED_TRANSPOSE, impl->e_vecs[i + impl->num_inputs], vec, request));
4360d0321e0SJeremy L Thompson   }
4370d0321e0SJeremy L Thompson 
4380d0321e0SJeremy L Thompson   // Restore input arrays
439*ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorRestoreInputs_Cuda(num_input_fields, qf_input_fields, op_input_fields, false, e_data, impl));
4400d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
4410d0321e0SJeremy L Thompson }
4420d0321e0SJeremy L Thompson 
4430d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
4440d0321e0SJeremy L Thompson // Core code for assembling linear QFunction
4450d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
4462b730f8bSJeremy L Thompson static inline int CeedOperatorLinearAssembleQFunctionCore_Cuda(CeedOperator op, bool build_objects, CeedVector *assembled, CeedElemRestriction *rstr,
4470d0321e0SJeremy L Thompson                                                                CeedRequest *request) {
448*ca735530SJeremy L Thompson   Ceed                ceed, ceed_parent;
449*ca735530SJeremy L Thompson   bool                is_identity_qf;
450*ca735530SJeremy L Thompson   CeedInt             num_active_in, num_active_out, Q, num_elem, num_input_fields, num_output_fields, size;
451e984cf9aSJeremy L Thompson   CeedSize            q_size;
452*ca735530SJeremy L Thompson   CeedScalar         *assembled_array, *e_data[2 * CEED_FIELD_MAX] = {NULL};
453*ca735530SJeremy L Thompson   CeedVector         *active_inputs;
454*ca735530SJeremy L Thompson   CeedQFunctionField *qf_input_fields, *qf_output_fields;
455*ca735530SJeremy L Thompson   CeedQFunction       qf;
456*ca735530SJeremy L Thompson   CeedOperatorField  *op_input_fields, *op_output_fields;
457*ca735530SJeremy L Thompson   CeedOperator_Cuda  *impl;
458*ca735530SJeremy L Thompson 
4592b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
460*ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetFallbackParentCeed(op, &ceed_parent));
461e984cf9aSJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &impl));
462e984cf9aSJeremy L Thompson   CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q));
463*ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem));
464e984cf9aSJeremy L Thompson   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
465*ca735530SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
466*ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
467*ca735530SJeremy L Thompson   active_inputs = impl->qf_active_in;
468*ca735530SJeremy L Thompson   num_active_in = impl->num_active_in, num_active_out = impl->num_active_out;
4690d0321e0SJeremy L Thompson 
4700d0321e0SJeremy L Thompson   // Setup
4712b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorSetup_Cuda(op));
4720d0321e0SJeremy L Thompson 
4730d0321e0SJeremy L Thompson   // Check for identity
474*ca735530SJeremy L Thompson   CeedCallBackend(CeedQFunctionIsIdentity(qf, &is_identity_qf));
475*ca735530SJeremy L Thompson   CeedCheck(!is_identity_qf, ceed, CEED_ERROR_BACKEND, "Assembling identity QFunctions not supported");
4760d0321e0SJeremy L Thompson 
477*ca735530SJeremy L Thompson   // Input e_vecs and Restriction
478*ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorSetupInputs_Cuda(num_input_fields, qf_input_fields, op_input_fields, NULL, true, e_data, impl, request));
4790d0321e0SJeremy L Thompson 
4800d0321e0SJeremy L Thompson   // Count number of active input fields
481*ca735530SJeremy L Thompson   if (!num_active_in) {
482*ca735530SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
483*ca735530SJeremy L Thompson       CeedScalar *q_vec_array;
484*ca735530SJeremy L Thompson       CeedVector  vec;
485*ca735530SJeremy L Thompson 
4860d0321e0SJeremy L Thompson       // Get input vector
487*ca735530SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec));
4880d0321e0SJeremy L Thompson       // Check if active input
4890d0321e0SJeremy L Thompson       if (vec == CEED_VECTOR_ACTIVE) {
490*ca735530SJeremy L Thompson         CeedCallBackend(CeedQFunctionFieldGetSize(qf_input_fields[i], &size));
491*ca735530SJeremy L Thompson         CeedCallBackend(CeedVectorSetValue(impl->q_vecs_in[i], 0.0));
492*ca735530SJeremy L Thompson         CeedCallBackend(CeedVectorGetArray(impl->q_vecs_in[i], CEED_MEM_DEVICE, &q_vec_array));
493*ca735530SJeremy L Thompson         CeedCallBackend(CeedRealloc(num_active_in + size, &active_inputs));
4940d0321e0SJeremy L Thompson         for (CeedInt field = 0; field < size; field++) {
495*ca735530SJeremy L Thompson           q_size = (CeedSize)Q * num_elem;
496*ca735530SJeremy L Thompson           CeedCallBackend(CeedVectorCreate(ceed, q_size, &active_inputs[num_active_in + field]));
497*ca735530SJeremy L Thompson           CeedCallBackend(
498*ca735530SJeremy L Thompson               CeedVectorSetArray(active_inputs[num_active_in + field], CEED_MEM_DEVICE, CEED_USE_POINTER, &q_vec_array[field * Q * num_elem]));
4990d0321e0SJeremy L Thompson         }
500*ca735530SJeremy L Thompson         num_active_in += size;
501*ca735530SJeremy L Thompson         CeedCallBackend(CeedVectorRestoreArray(impl->q_vecs_in[i], &q_vec_array));
5020d0321e0SJeremy L Thompson       }
5030d0321e0SJeremy L Thompson     }
504*ca735530SJeremy L Thompson     impl->num_active_in = num_active_in;
505*ca735530SJeremy L Thompson     impl->qf_active_in  = active_inputs;
5060d0321e0SJeremy L Thompson   }
5070d0321e0SJeremy L Thompson 
5080d0321e0SJeremy L Thompson   // Count number of active output fields
509*ca735530SJeremy L Thompson   if (!num_active_out) {
510*ca735530SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
511*ca735530SJeremy L Thompson       CeedVector vec;
512*ca735530SJeremy L Thompson 
5130d0321e0SJeremy L Thompson       // Get output vector
514*ca735530SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[i], &vec));
5150d0321e0SJeremy L Thompson       // Check if active output
5160d0321e0SJeremy L Thompson       if (vec == CEED_VECTOR_ACTIVE) {
517*ca735530SJeremy L Thompson         CeedCallBackend(CeedQFunctionFieldGetSize(qf_output_fields[i], &size));
518*ca735530SJeremy L Thompson         num_active_out += size;
5190d0321e0SJeremy L Thompson       }
5200d0321e0SJeremy L Thompson     }
521*ca735530SJeremy L Thompson     impl->num_active_out = num_active_out;
5220d0321e0SJeremy L Thompson   }
5230d0321e0SJeremy L Thompson 
5240d0321e0SJeremy L Thompson   // Check sizes
525*ca735530SJeremy L Thompson   CeedCheck(num_active_in > 0 && num_active_out > 0, ceed, CEED_ERROR_BACKEND, "Cannot assemble QFunction without active inputs and outputs");
5260d0321e0SJeremy L Thompson 
5270d0321e0SJeremy L Thompson   // Build objects if needed
5280d0321e0SJeremy L Thompson   if (build_objects) {
5290d0321e0SJeremy L Thompson     // Create output restriction
530*ca735530SJeremy L Thompson     CeedInt strides[3] = {1, num_elem * Q, Q}; /* *NOPAD* */
531*ca735530SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionCreateStrided(ceed_parent, num_elem, Q, num_active_in * num_active_out,
532*ca735530SJeremy L Thompson                                                      num_active_in * num_active_out * num_elem * Q, strides, rstr));
5330d0321e0SJeremy L Thompson     // Create assembled vector
534*ca735530SJeremy L Thompson     CeedSize l_size = (CeedSize)num_elem * Q * num_active_in * num_active_out;
535*ca735530SJeremy L Thompson     CeedCallBackend(CeedVectorCreate(ceed_parent, l_size, assembled));
5360d0321e0SJeremy L Thompson   }
5372b730f8bSJeremy L Thompson   CeedCallBackend(CeedVectorSetValue(*assembled, 0.0));
538*ca735530SJeremy L Thompson   CeedCallBackend(CeedVectorGetArray(*assembled, CEED_MEM_DEVICE, &assembled_array));
5390d0321e0SJeremy L Thompson 
5400d0321e0SJeremy L Thompson   // Input basis apply
541*ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorInputBasis_Cuda(num_elem, qf_input_fields, op_input_fields, num_input_fields, true, e_data, impl));
5420d0321e0SJeremy L Thompson 
5430d0321e0SJeremy L Thompson   // Assemble QFunction
544*ca735530SJeremy L Thompson   for (CeedInt in = 0; in < num_active_in; in++) {
5450d0321e0SJeremy L Thompson     // Set Inputs
546*ca735530SJeremy L Thompson     CeedCallBackend(CeedVectorSetValue(active_inputs[in], 1.0));
547*ca735530SJeremy L Thompson     if (num_active_in > 1) {
548*ca735530SJeremy L Thompson       CeedCallBackend(CeedVectorSetValue(active_inputs[(in + num_active_in - 1) % num_active_in], 0.0));
5490d0321e0SJeremy L Thompson     }
5500d0321e0SJeremy L Thompson     // Set Outputs
551*ca735530SJeremy L Thompson     for (CeedInt out = 0; out < num_output_fields; out++) {
552*ca735530SJeremy L Thompson       CeedVector vec;
553*ca735530SJeremy L Thompson 
5540d0321e0SJeremy L Thompson       // Get output vector
555*ca735530SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[out], &vec));
5560d0321e0SJeremy L Thompson       // Check if active output
5570d0321e0SJeremy L Thompson       if (vec == CEED_VECTOR_ACTIVE) {
558*ca735530SJeremy L Thompson         CeedCallBackend(CeedVectorSetArray(impl->q_vecs_out[out], CEED_MEM_DEVICE, CEED_USE_POINTER, assembled_array));
559*ca735530SJeremy L Thompson         CeedCallBackend(CeedQFunctionFieldGetSize(qf_output_fields[out], &size));
560*ca735530SJeremy L Thompson         assembled_array += size * Q * num_elem;  // Advance the pointer by the size of the output
5610d0321e0SJeremy L Thompson       }
5620d0321e0SJeremy L Thompson     }
5630d0321e0SJeremy L Thompson     // Apply QFunction
564*ca735530SJeremy L Thompson     CeedCallBackend(CeedQFunctionApply(qf, Q * num_elem, impl->q_vecs_in, impl->q_vecs_out));
5650d0321e0SJeremy L Thompson   }
5660d0321e0SJeremy L Thompson 
567*ca735530SJeremy L Thompson   // Un-set output q_vecs to prevent accidental overwrite of Assembled
568*ca735530SJeremy L Thompson   for (CeedInt out = 0; out < num_output_fields; out++) {
569*ca735530SJeremy L Thompson     CeedVector vec;
570*ca735530SJeremy L Thompson 
5710d0321e0SJeremy L Thompson     // Get output vector
572*ca735530SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[out], &vec));
5730d0321e0SJeremy L Thompson     // Check if active output
5740d0321e0SJeremy L Thompson     if (vec == CEED_VECTOR_ACTIVE) {
575*ca735530SJeremy L Thompson       CeedCallBackend(CeedVectorTakeArray(impl->q_vecs_out[out], CEED_MEM_DEVICE, NULL));
5760d0321e0SJeremy L Thompson     }
5770d0321e0SJeremy L Thompson   }
5780d0321e0SJeremy L Thompson 
5790d0321e0SJeremy L Thompson   // Restore input arrays
580*ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorRestoreInputs_Cuda(num_input_fields, qf_input_fields, op_input_fields, true, e_data, impl));
5810d0321e0SJeremy L Thompson 
5820d0321e0SJeremy L Thompson   // Restore output
583*ca735530SJeremy L Thompson   CeedCallBackend(CeedVectorRestoreArray(*assembled, &assembled_array));
5840d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
5850d0321e0SJeremy L Thompson }
5860d0321e0SJeremy L Thompson 
5870d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
5880d0321e0SJeremy L Thompson // Assemble Linear QFunction
5890d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
5902b730f8bSJeremy L Thompson static int CeedOperatorLinearAssembleQFunction_Cuda(CeedOperator op, CeedVector *assembled, CeedElemRestriction *rstr, CeedRequest *request) {
5912b730f8bSJeremy L Thompson   return CeedOperatorLinearAssembleQFunctionCore_Cuda(op, true, assembled, rstr, request);
5920d0321e0SJeremy L Thompson }
5930d0321e0SJeremy L Thompson 
5940d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
5950d0321e0SJeremy L Thompson // Update Assembled Linear QFunction
5960d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
5972b730f8bSJeremy L Thompson static int CeedOperatorLinearAssembleQFunctionUpdate_Cuda(CeedOperator op, CeedVector assembled, CeedElemRestriction rstr, CeedRequest *request) {
5982b730f8bSJeremy L Thompson   return CeedOperatorLinearAssembleQFunctionCore_Cuda(op, false, &assembled, &rstr, request);
5990d0321e0SJeremy L Thompson }
6000d0321e0SJeremy L Thompson 
6010d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
6020d0321e0SJeremy L Thompson // Create point block restriction
6030d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
604*ca735530SJeremy L Thompson static int CreatePointBlockRestriction(CeedElemRestriction rstr, CeedElemRestriction *point_block_rstr) {
6050d0321e0SJeremy L Thompson   Ceed           ceed;
606*ca735530SJeremy L Thompson   CeedSize       l_size;
607*ca735530SJeremy L Thompson   CeedInt        num_elem, num_comp, elem_size, comp_stride, *point_block_offsets;
6080d0321e0SJeremy L Thompson   const CeedInt *offsets;
609*ca735530SJeremy L Thompson 
610*ca735530SJeremy L Thompson   CeedCallBackend(CeedElemRestrictionGetCeed(rstr, &ceed));
6112b730f8bSJeremy L Thompson   CeedCallBackend(CeedElemRestrictionGetOffsets(rstr, CEED_MEM_HOST, &offsets));
6120d0321e0SJeremy L Thompson 
6130d0321e0SJeremy L Thompson   // Expand offsets
614*ca735530SJeremy L Thompson   CeedCallBackend(CeedElemRestrictionGetNumElements(rstr, &num_elem));
615*ca735530SJeremy L Thompson   CeedCallBackend(CeedElemRestrictionGetNumComponents(rstr, &num_comp));
616*ca735530SJeremy L Thompson   CeedCallBackend(CeedElemRestrictionGetElementSize(rstr, &elem_size));
617*ca735530SJeremy L Thompson   CeedCallBackend(CeedElemRestrictionGetCompStride(rstr, &comp_stride));
6182b730f8bSJeremy L Thompson   CeedCallBackend(CeedElemRestrictionGetLVectorSize(rstr, &l_size));
619*ca735530SJeremy L Thompson   CeedInt shift = num_comp;
620*ca735530SJeremy L Thompson 
621*ca735530SJeremy L Thompson   if (comp_stride != 1) shift *= num_comp;
622*ca735530SJeremy L Thompson   CeedCallBackend(CeedCalloc(num_elem * elem_size, &point_block_offsets));
623*ca735530SJeremy L Thompson   for (CeedInt i = 0; i < num_elem * elem_size; i++) {
624*ca735530SJeremy L Thompson     point_block_offsets[i] = offsets[i] * shift;
6250d0321e0SJeremy L Thompson   }
6260d0321e0SJeremy L Thompson 
6270d0321e0SJeremy L Thompson   // Create new restriction
628*ca735530SJeremy L Thompson   CeedCallBackend(CeedElemRestrictionCreate(ceed, num_elem, elem_size, num_comp * num_comp, 1, l_size * num_comp, CEED_MEM_HOST, CEED_OWN_POINTER,
629*ca735530SJeremy L Thompson                                             point_block_offsets, point_block_rstr));
6300d0321e0SJeremy L Thompson 
6310d0321e0SJeremy L Thompson   // Cleanup
6322b730f8bSJeremy L Thompson   CeedCallBackend(CeedElemRestrictionRestoreOffsets(rstr, &offsets));
6330d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
6340d0321e0SJeremy L Thompson }
6350d0321e0SJeremy L Thompson 
6360d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
6370d0321e0SJeremy L Thompson // Assemble diagonal setup
6380d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
639*ca735530SJeremy L Thompson static inline int CeedOperatorAssembleDiagonalSetup_Cuda(CeedOperator op, const bool is_point_block, CeedInt use_ceedsize_idx) {
6400d0321e0SJeremy L Thompson   Ceed                ceed;
641*ca735530SJeremy L Thompson   char               *diagonal_kernel_path, *diagonal_kernel_source;
642*ca735530SJeremy L Thompson   CeedInt             num_input_fields, num_output_fields, num_e_mode_in = 0, num_comp = 0, dim = 1, num_e_mode_out = 0, num_nodes, num_qpts;
643*ca735530SJeremy L Thompson   CeedEvalMode       *e_mode_in = NULL, *e_mode_out = NULL;
644*ca735530SJeremy L Thompson   CeedElemRestriction rstr_in = NULL, rstr_out = NULL;
645*ca735530SJeremy L Thompson   CeedBasis           basis_in = NULL, basis_out = NULL;
646*ca735530SJeremy L Thompson   CeedQFunctionField *qf_fields;
6470d0321e0SJeremy L Thompson   CeedQFunction       qf;
648*ca735530SJeremy L Thompson   CeedOperatorField  *op_fields;
649*ca735530SJeremy L Thompson   CeedOperator_Cuda  *impl;
650*ca735530SJeremy L Thompson 
651*ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
6522b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
653*ca735530SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetNumArgs(qf, &num_input_fields, &num_output_fields));
6540d0321e0SJeremy L Thompson 
6550d0321e0SJeremy L Thompson   // Determine active input basis
656*ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetFields(op, NULL, &op_fields, NULL, NULL));
657*ca735530SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_fields, NULL, NULL));
658*ca735530SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
6590d0321e0SJeremy L Thompson     CeedVector vec;
660*ca735530SJeremy L Thompson 
661*ca735530SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec));
6620d0321e0SJeremy L Thompson     if (vec == CEED_VECTOR_ACTIVE) {
663*ca735530SJeremy L Thompson       CeedEvalMode        e_mode;
6640d0321e0SJeremy L Thompson       CeedElemRestriction rstr;
665*ca735530SJeremy L Thompson 
666*ca735530SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis_in));
667*ca735530SJeremy L Thompson       CeedCallBackend(CeedBasisGetNumComponents(basis_in, &num_comp));
668*ca735530SJeremy L Thompson       CeedCallBackend(CeedBasisGetDimension(basis_in, &dim));
669*ca735530SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_fields[i], &rstr));
670*ca735530SJeremy L Thompson       CeedCheck(!rstr_in || rstr_in == rstr, ceed, CEED_ERROR_BACKEND,
6716574a04fSJeremy L Thompson                 "Backend does not implement multi-field non-composite operator diagonal assembly");
672*ca735530SJeremy L Thompson       rstr_in = rstr;
673*ca735530SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &e_mode));
674*ca735530SJeremy L Thompson       switch (e_mode) {
6750d0321e0SJeremy L Thompson         case CEED_EVAL_NONE:
6760d0321e0SJeremy L Thompson         case CEED_EVAL_INTERP:
677*ca735530SJeremy L Thompson           CeedCallBackend(CeedRealloc(num_e_mode_in + 1, &e_mode_in));
678*ca735530SJeremy L Thompson           e_mode_in[num_e_mode_in] = e_mode;
679*ca735530SJeremy L Thompson           num_e_mode_in += 1;
6800d0321e0SJeremy L Thompson           break;
6810d0321e0SJeremy L Thompson         case CEED_EVAL_GRAD:
682*ca735530SJeremy L Thompson           CeedCallBackend(CeedRealloc(num_e_mode_in + dim, &e_mode_in));
683*ca735530SJeremy L Thompson           for (CeedInt d = 0; d < dim; d++) e_mode_in[num_e_mode_in + d] = e_mode;
684*ca735530SJeremy L Thompson           num_e_mode_in += dim;
6850d0321e0SJeremy L Thompson           break;
6860d0321e0SJeremy L Thompson         case CEED_EVAL_WEIGHT:
6870d0321e0SJeremy L Thompson         case CEED_EVAL_DIV:
6880d0321e0SJeremy L Thompson         case CEED_EVAL_CURL:
6890d0321e0SJeremy L Thompson           break;  // Caught by QF Assembly
6900d0321e0SJeremy L Thompson       }
6910d0321e0SJeremy L Thompson     }
6920d0321e0SJeremy L Thompson   }
6930d0321e0SJeremy L Thompson 
6940d0321e0SJeremy L Thompson   // Determine active output basis
695*ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetFields(op, NULL, NULL, NULL, &op_fields));
696*ca735530SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qf_fields));
697*ca735530SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
6980d0321e0SJeremy L Thompson     CeedVector vec;
699*ca735530SJeremy L Thompson 
700*ca735530SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec));
7010d0321e0SJeremy L Thompson     if (vec == CEED_VECTOR_ACTIVE) {
702*ca735530SJeremy L Thompson       CeedEvalMode        e_mode;
7030d0321e0SJeremy L Thompson       CeedElemRestriction rstr;
704*ca735530SJeremy L Thompson 
705*ca735530SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis_out));
706*ca735530SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_fields[i], &rstr));
707*ca735530SJeremy L Thompson       CeedCheck(!rstr_out || rstr_out == rstr, ceed, CEED_ERROR_BACKEND,
7086574a04fSJeremy L Thompson                 "Backend does not implement multi-field non-composite operator diagonal assembly");
709*ca735530SJeremy L Thompson       rstr_out = rstr;
710*ca735530SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &e_mode));
711*ca735530SJeremy L Thompson       switch (e_mode) {
7120d0321e0SJeremy L Thompson         case CEED_EVAL_NONE:
7130d0321e0SJeremy L Thompson         case CEED_EVAL_INTERP:
714*ca735530SJeremy L Thompson           CeedCallBackend(CeedRealloc(num_e_mode_out + 1, &e_mode_out));
715*ca735530SJeremy L Thompson           e_mode_out[num_e_mode_out] = e_mode;
716*ca735530SJeremy L Thompson           num_e_mode_out += 1;
7170d0321e0SJeremy L Thompson           break;
7180d0321e0SJeremy L Thompson         case CEED_EVAL_GRAD:
719*ca735530SJeremy L Thompson           CeedCallBackend(CeedRealloc(num_e_mode_out + dim, &e_mode_out));
720*ca735530SJeremy L Thompson           for (CeedInt d = 0; d < dim; d++) e_mode_out[num_e_mode_out + d] = e_mode;
721*ca735530SJeremy L Thompson           num_e_mode_out += dim;
7220d0321e0SJeremy L Thompson           break;
7230d0321e0SJeremy L Thompson         case CEED_EVAL_WEIGHT:
7240d0321e0SJeremy L Thompson         case CEED_EVAL_DIV:
7250d0321e0SJeremy L Thompson         case CEED_EVAL_CURL:
7260d0321e0SJeremy L Thompson           break;  // Caught by QF Assembly
7270d0321e0SJeremy L Thompson       }
7280d0321e0SJeremy L Thompson     }
7290d0321e0SJeremy L Thompson   }
7300d0321e0SJeremy L Thompson 
7310d0321e0SJeremy L Thompson   // Operator data struct
7322b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &impl));
7332b730f8bSJeremy L Thompson   CeedCallBackend(CeedCalloc(1, &impl->diag));
7340d0321e0SJeremy L Thompson   CeedOperatorDiag_Cuda *diag = impl->diag;
735*ca735530SJeremy L Thompson 
736*ca735530SJeremy L Thompson   diag->basis_in       = basis_in;
737*ca735530SJeremy L Thompson   diag->basis_out      = basis_out;
738*ca735530SJeremy L Thompson   diag->h_e_mode_in    = e_mode_in;
739*ca735530SJeremy L Thompson   diag->h_e_mode_out   = e_mode_out;
740*ca735530SJeremy L Thompson   diag->num_e_mode_in  = num_e_mode_in;
741*ca735530SJeremy L Thompson   diag->num_e_mode_out = num_e_mode_out;
7420d0321e0SJeremy L Thompson 
7430d0321e0SJeremy L Thompson   // Assemble kernel
7442b730f8bSJeremy L Thompson   CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-operator-assemble-diagonal.h", &diagonal_kernel_path));
74523d4529eSJeremy L Thompson   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Diagonal Assembly Kernel Source -----\n");
7462b730f8bSJeremy L Thompson   CeedCallBackend(CeedLoadSourceToBuffer(ceed, diagonal_kernel_path, &diagonal_kernel_source));
74723d4529eSJeremy L Thompson   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Diagonal Assembly Source Complete! -----\n");
748*ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetNumNodes(basis_in, &num_nodes));
749*ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis_in, &num_qpts));
750*ca735530SJeremy L Thompson   diag->num_nodes = num_nodes;
751*ca735530SJeremy L Thompson   CeedCallCuda(ceed,
752*ca735530SJeremy L Thompson                CeedCompile_Cuda(ceed, diagonal_kernel_source, &diag->module, 6, "NUM_E_MODE_IN", num_e_mode_in, "NUM_E_MODE_OUT", num_e_mode_out,
753*ca735530SJeremy L Thompson                                 "NUM_NODES", num_nodes, "NUM_QPTS", num_qpts, "NUM_COMP", num_comp, "USE_CEEDSIZE", use_ceedsize_idx));
754eb7e6cafSJeremy L Thompson   CeedCallCuda(ceed, CeedGetKernel_Cuda(ceed, diag->module, "linearDiagonal", &diag->linearDiagonal));
755eb7e6cafSJeremy L Thompson   CeedCallCuda(ceed, CeedGetKernel_Cuda(ceed, diag->module, "linearPointBlockDiagonal", &diag->linearPointBlock));
7562b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&diagonal_kernel_path));
7572b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&diagonal_kernel_source));
7580d0321e0SJeremy L Thompson 
7590d0321e0SJeremy L Thompson   // Basis matrices
760*ca735530SJeremy L Thompson   const CeedInt     q_bytes      = num_qpts * sizeof(CeedScalar);
761*ca735530SJeremy L Thompson   const CeedInt     interp_bytes = q_bytes * num_nodes;
762*ca735530SJeremy L Thompson   const CeedInt     grad_bytes   = q_bytes * num_nodes * dim;
763*ca735530SJeremy L Thompson   const CeedInt     e_mode_bytes = sizeof(CeedEvalMode);
764*ca735530SJeremy L Thompson   const CeedScalar *interp_in, *interp_out, *grad_in, *grad_out;
7650d0321e0SJeremy L Thompson 
7660d0321e0SJeremy L Thompson   // CEED_EVAL_NONE
7670d0321e0SJeremy L Thompson   CeedScalar *identity     = NULL;
768*ca735530SJeremy L Thompson   bool        is_eval_none = false;
769*ca735530SJeremy L Thompson 
770*ca735530SJeremy L Thompson   for (CeedInt i = 0; i < num_e_mode_in; i++) is_eval_none = is_eval_none || (e_mode_in[i] == CEED_EVAL_NONE);
771*ca735530SJeremy L Thompson   for (CeedInt i = 0; i < num_e_mode_out; i++) is_eval_none = is_eval_none || (e_mode_out[i] == CEED_EVAL_NONE);
772*ca735530SJeremy L Thompson   if (is_eval_none) {
773*ca735530SJeremy L Thompson     CeedCallBackend(CeedCalloc(num_qpts * num_nodes, &identity));
774*ca735530SJeremy L Thompson     for (CeedInt i = 0; i < (num_nodes < num_qpts ? num_nodes : num_qpts); i++) identity[i * num_nodes + i] = 1.0;
775*ca735530SJeremy L Thompson     CeedCallCuda(ceed, cudaMalloc((void **)&diag->d_identity, interp_bytes));
776*ca735530SJeremy L Thompson     CeedCallCuda(ceed, cudaMemcpy(diag->d_identity, identity, interp_bytes, cudaMemcpyHostToDevice));
7770d0321e0SJeremy L Thompson   }
7780d0321e0SJeremy L Thompson 
7790d0321e0SJeremy L Thompson   // CEED_EVAL_INTERP
780*ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetInterp(basis_in, &interp_in));
781*ca735530SJeremy L Thompson   CeedCallCuda(ceed, cudaMalloc((void **)&diag->d_interp_in, interp_bytes));
782*ca735530SJeremy L Thompson   CeedCallCuda(ceed, cudaMemcpy(diag->d_interp_in, interp_in, interp_bytes, cudaMemcpyHostToDevice));
783*ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetInterp(basis_out, &interp_out));
784*ca735530SJeremy L Thompson   CeedCallCuda(ceed, cudaMalloc((void **)&diag->d_interp_out, interp_bytes));
785*ca735530SJeremy L Thompson   CeedCallCuda(ceed, cudaMemcpy(diag->d_interp_out, interp_out, interp_bytes, cudaMemcpyHostToDevice));
7860d0321e0SJeremy L Thompson 
7870d0321e0SJeremy L Thompson   // CEED_EVAL_GRAD
788*ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetGrad(basis_in, &grad_in));
789*ca735530SJeremy L Thompson   CeedCallCuda(ceed, cudaMalloc((void **)&diag->d_grad_in, grad_bytes));
790*ca735530SJeremy L Thompson   CeedCallCuda(ceed, cudaMemcpy(diag->d_grad_in, grad_in, grad_bytes, cudaMemcpyHostToDevice));
791*ca735530SJeremy L Thompson   CeedCallBackend(CeedBasisGetGrad(basis_out, &grad_out));
792*ca735530SJeremy L Thompson   CeedCallCuda(ceed, cudaMalloc((void **)&diag->d_grad_out, grad_bytes));
793*ca735530SJeremy L Thompson   CeedCallCuda(ceed, cudaMemcpy(diag->d_grad_out, grad_out, grad_bytes, cudaMemcpyHostToDevice));
7940d0321e0SJeremy L Thompson 
795*ca735530SJeremy L Thompson   // Arrays of e_modes
796*ca735530SJeremy L Thompson   CeedCallCuda(ceed, cudaMalloc((void **)&diag->d_e_mode_in, num_e_mode_in * e_mode_bytes));
797*ca735530SJeremy L Thompson   CeedCallCuda(ceed, cudaMemcpy(diag->d_e_mode_in, e_mode_in, num_e_mode_in * e_mode_bytes, cudaMemcpyHostToDevice));
798*ca735530SJeremy L Thompson   CeedCallCuda(ceed, cudaMalloc((void **)&diag->d_e_mode_out, num_e_mode_out * e_mode_bytes));
799*ca735530SJeremy L Thompson   CeedCallCuda(ceed, cudaMemcpy(diag->d_e_mode_out, e_mode_out, num_e_mode_out * e_mode_bytes, cudaMemcpyHostToDevice));
8000d0321e0SJeremy L Thompson 
8010d0321e0SJeremy L Thompson   // Restriction
802*ca735530SJeremy L Thompson   diag->diag_rstr = rstr_out;
8030d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
8040d0321e0SJeremy L Thompson }
8050d0321e0SJeremy L Thompson 
8060d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
8070d0321e0SJeremy L Thompson // Assemble diagonal common code
8080d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
809*ca735530SJeremy L Thompson static inline int CeedOperatorAssembleDiagonalCore_Cuda(CeedOperator op, CeedVector assembled, CeedRequest *request, const bool is_point_block) {
8100d0321e0SJeremy L Thompson   Ceed                ceed;
811*ca735530SJeremy L Thompson   CeedSize            assembled_length = 0, assembled_qf_length = 0;
812*ca735530SJeremy L Thompson   CeedInt             use_ceedsize_idx = 0, num_elem;
813*ca735530SJeremy L Thompson   CeedScalar         *elem_diag_array;
814*ca735530SJeremy L Thompson   const CeedScalar   *assembled_qf_array;
815*ca735530SJeremy L Thompson   CeedVector          assembled_qf = NULL;
816*ca735530SJeremy L Thompson   CeedElemRestriction rstr         = NULL;
8170d0321e0SJeremy L Thompson   CeedOperator_Cuda  *impl;
818*ca735530SJeremy L Thompson 
819*ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
8202b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &impl));
8210d0321e0SJeremy L Thompson 
8220d0321e0SJeremy L Thompson   // Assemble QFunction
823*ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorLinearAssembleQFunctionBuildOrUpdate(op, &assembled_qf, &rstr, request));
8242b730f8bSJeremy L Thompson   CeedCallBackend(CeedElemRestrictionDestroy(&rstr));
8250d0321e0SJeremy L Thompson 
826f7c1b517Snbeams   CeedCallBackend(CeedVectorGetLength(assembled, &assembled_length));
827*ca735530SJeremy L Thompson   CeedCallBackend(CeedVectorGetLength(assembled_qf, &assembled_qf_length));
828*ca735530SJeremy L Thompson   if ((assembled_length > INT_MAX) || (assembled_qf_length > INT_MAX)) use_ceedsize_idx = 1;
829f7c1b517Snbeams 
8300d0321e0SJeremy L Thompson   // Setup
831*ca735530SJeremy L Thompson   if (!impl->diag) CeedCallBackend(CeedOperatorAssembleDiagonalSetup_Cuda(op, is_point_block, use_ceedsize_idx));
8320d0321e0SJeremy L Thompson   CeedOperatorDiag_Cuda *diag = impl->diag;
833*ca735530SJeremy L Thompson 
8340d0321e0SJeremy L Thompson   assert(diag != NULL);
8350d0321e0SJeremy L Thompson 
8360d0321e0SJeremy L Thompson   // Restriction
837*ca735530SJeremy L Thompson   if (is_point_block && !diag->point_block_rstr) {
838*ca735530SJeremy L Thompson     CeedElemRestriction point_block_rstr;
839*ca735530SJeremy L Thompson 
840*ca735530SJeremy L Thompson     CeedCallBackend(CreatePointBlockRestriction(diag->diag_rstr, &point_block_rstr));
841*ca735530SJeremy L Thompson     diag->point_block_rstr = point_block_rstr;
8420d0321e0SJeremy L Thompson   }
843*ca735530SJeremy L Thompson   CeedElemRestriction diag_rstr = is_point_block ? diag->point_block_rstr : diag->diag_rstr;
8440d0321e0SJeremy L Thompson 
8450d0321e0SJeremy L Thompson   // Create diagonal vector
846*ca735530SJeremy L Thompson   CeedVector elem_diag = is_point_block ? diag->point_block_elem_diag : diag->elem_diag;
847*ca735530SJeremy L Thompson 
848*ca735530SJeremy L Thompson   if (!elem_diag) {
849*ca735530SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionCreateVector(diag_rstr, NULL, &elem_diag));
850*ca735530SJeremy L Thompson     if (is_point_block) diag->point_block_elem_diag = elem_diag;
851*ca735530SJeremy L Thompson     else diag->elem_diag = elem_diag;
8520d0321e0SJeremy L Thompson   }
853*ca735530SJeremy L Thompson   CeedCallBackend(CeedVectorSetValue(elem_diag, 0.0));
8540d0321e0SJeremy L Thompson 
8550d0321e0SJeremy L Thompson   // Assemble element operator diagonals
856*ca735530SJeremy L Thompson   CeedCallBackend(CeedVectorGetArray(elem_diag, CEED_MEM_DEVICE, &elem_diag_array));
857*ca735530SJeremy L Thompson   CeedCallBackend(CeedVectorGetArrayRead(assembled_qf, CEED_MEM_DEVICE, &assembled_qf_array));
858*ca735530SJeremy L Thompson   CeedCallBackend(CeedElemRestrictionGetNumElements(diag_rstr, &num_elem));
8590d0321e0SJeremy L Thompson 
8600d0321e0SJeremy L Thompson   // Compute the diagonal of B^T D B
861*ca735530SJeremy L Thompson   int   elem_per_block = 1;
862*ca735530SJeremy L Thompson   int   grid           = num_elem / elem_per_block + ((num_elem / elem_per_block * elem_per_block < num_elem) ? 1 : 0);
863*ca735530SJeremy L Thompson   void *args[]         = {(void *)&num_elem, &diag->d_identity,  &diag->d_interp_in,  &diag->d_grad_in,    &diag->d_interp_out,
864*ca735530SJeremy L Thompson                           &diag->d_grad_out, &diag->d_e_mode_in, &diag->d_e_mode_out, &assembled_qf_array, &elem_diag_array};
865*ca735530SJeremy L Thompson   if (is_point_block) {
866*ca735530SJeremy L Thompson     CeedCallBackend(CeedRunKernelDim_Cuda(ceed, diag->linearPointBlock, grid, diag->num_nodes, 1, elem_per_block, args));
8670d0321e0SJeremy L Thompson   } else {
868*ca735530SJeremy L Thompson     CeedCallBackend(CeedRunKernelDim_Cuda(ceed, diag->linearDiagonal, grid, diag->num_nodes, 1, elem_per_block, args));
8690d0321e0SJeremy L Thompson   }
8700d0321e0SJeremy L Thompson 
8710d0321e0SJeremy L Thompson   // Restore arrays
872*ca735530SJeremy L Thompson   CeedCallBackend(CeedVectorRestoreArray(elem_diag, &elem_diag_array));
873*ca735530SJeremy L Thompson   CeedCallBackend(CeedVectorRestoreArrayRead(assembled_qf, &assembled_qf_array));
8740d0321e0SJeremy L Thompson 
8750d0321e0SJeremy L Thompson   // Assemble local operator diagonal
876*ca735530SJeremy L Thompson   CeedCallBackend(CeedElemRestrictionApply(diag_rstr, CEED_TRANSPOSE, elem_diag, assembled, request));
8770d0321e0SJeremy L Thompson 
8780d0321e0SJeremy L Thompson   // Cleanup
879*ca735530SJeremy L Thompson   CeedCallBackend(CeedVectorDestroy(&assembled_qf));
8800d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
8810d0321e0SJeremy L Thompson }
8820d0321e0SJeremy L Thompson 
8830d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
8840d0321e0SJeremy L Thompson // Assemble Linear Diagonal
8850d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
8862b730f8bSJeremy L Thompson static int CeedOperatorLinearAssembleAddDiagonal_Cuda(CeedOperator op, CeedVector assembled, CeedRequest *request) {
8872b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorAssembleDiagonalCore_Cuda(op, assembled, request, false));
8886aa95790SJeremy L Thompson   return CEED_ERROR_SUCCESS;
8890d0321e0SJeremy L Thompson }
8900d0321e0SJeremy L Thompson 
8910d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
8920d0321e0SJeremy L Thompson // Assemble Linear Point Block Diagonal
8930d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
8942b730f8bSJeremy L Thompson static int CeedOperatorLinearAssembleAddPointBlockDiagonal_Cuda(CeedOperator op, CeedVector assembled, CeedRequest *request) {
8952b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorAssembleDiagonalCore_Cuda(op, assembled, request, true));
8966aa95790SJeremy L Thompson   return CEED_ERROR_SUCCESS;
8970d0321e0SJeremy L Thompson }
8980d0321e0SJeremy L Thompson 
8990d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
900cc132f9aSnbeams // Single operator assembly setup
901cc132f9aSnbeams //------------------------------------------------------------------------------
902f7c1b517Snbeams static int CeedSingleOperatorAssembleSetup_Cuda(CeedOperator op, CeedInt use_ceedsize_idx) {
903cc132f9aSnbeams   Ceed    ceed;
904*ca735530SJeremy L Thompson   char   *assembly_kernel_path, *assembly_kernel_source;
905*ca735530SJeremy L Thompson   CeedInt num_input_fields, num_output_fields, num_e_mode_in = 0, dim = 1, num_B_in_mats_to_load = 0, size_B_in = 0, num_qpts = 0, elem_size = 0,
906*ca735530SJeremy L Thompson                                                num_e_mode_out = 0, num_B_out_mats_to_load = 0, size_B_out = 0, num_elem, num_comp;
907*ca735530SJeremy L Thompson   CeedEvalMode       *eval_mode_in = NULL, *eval_mode_out = NULL;
908*ca735530SJeremy L Thompson   CeedElemRestriction rstr_in = NULL, rstr_out = NULL;
909*ca735530SJeremy L Thompson   CeedBasis           basis_in = NULL, basis_out = NULL;
910*ca735530SJeremy L Thompson   CeedQFunctionField *qf_fields;
911*ca735530SJeremy L Thompson   CeedQFunction       qf;
912*ca735530SJeremy L Thompson   CeedOperatorField  *input_fields, *output_fields;
913cc132f9aSnbeams   CeedOperator_Cuda  *impl;
914*ca735530SJeremy L Thompson 
915*ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
9162b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &impl));
917cc132f9aSnbeams 
918cc132f9aSnbeams   // Get intput and output fields
9192b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &input_fields, &num_output_fields, &output_fields));
920cc132f9aSnbeams 
921cc132f9aSnbeams   // Determine active input basis eval mode
9222b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
9232b730f8bSJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_fields, NULL, NULL));
924cc132f9aSnbeams   // Note that the kernel will treat each dimension of a gradient action separately;
925*ca735530SJeremy L Thompson   // i.e., when an active input has a CEED_EVAL_GRAD mode, num_e_mode_in will increment by dim.
926ea61e9acSJeremy L Thompson   // However, for the purposes of loading the B matrices, it will be treated as one mode, and we will load/copy the entire gradient matrix at once, so
927cc132f9aSnbeams   // num_B_in_mats_to_load will be incremented by 1.
928cc132f9aSnbeams   for (CeedInt i = 0; i < num_input_fields; i++) {
929cc132f9aSnbeams     CeedVector vec;
930*ca735530SJeremy L Thompson 
9312b730f8bSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetVector(input_fields[i], &vec));
932cc132f9aSnbeams     if (vec == CEED_VECTOR_ACTIVE) {
933*ca735530SJeremy L Thompson       CeedEvalMode eval_mode;
934*ca735530SJeremy L Thompson 
9352b730f8bSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetBasis(input_fields[i], &basis_in));
9362b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetDimension(basis_in, &dim));
937*ca735530SJeremy L Thompson       CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis_in, &num_qpts));
9382b730f8bSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetElemRestriction(input_fields[i], &rstr_in));
939*ca735530SJeremy L Thompson       CeedCallBackend(CeedElemRestrictionGetElementSize(rstr_in, &elem_size));
9402b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode));
941cc132f9aSnbeams       if (eval_mode != CEED_EVAL_NONE) {
9422b730f8bSJeremy L Thompson         CeedCallBackend(CeedRealloc(num_B_in_mats_to_load + 1, &eval_mode_in));
943cc132f9aSnbeams         eval_mode_in[num_B_in_mats_to_load] = eval_mode;
944cc132f9aSnbeams         num_B_in_mats_to_load += 1;
945cc132f9aSnbeams         if (eval_mode == CEED_EVAL_GRAD) {
946*ca735530SJeremy L Thompson           num_e_mode_in += dim;
947*ca735530SJeremy L Thompson           size_B_in += dim * elem_size * num_qpts;
948cc132f9aSnbeams         } else {
949*ca735530SJeremy L Thompson           num_e_mode_in += 1;
950*ca735530SJeremy L Thompson           size_B_in += elem_size * num_qpts;
951cc132f9aSnbeams         }
952cc132f9aSnbeams       }
953cc132f9aSnbeams     }
954cc132f9aSnbeams   }
955cc132f9aSnbeams 
956cc132f9aSnbeams   // Determine active output basis; basis_out and rstr_out only used if same as input, TODO
9572b730f8bSJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qf_fields));
958cc132f9aSnbeams   for (CeedInt i = 0; i < num_output_fields; i++) {
959cc132f9aSnbeams     CeedVector vec;
960*ca735530SJeremy L Thompson 
9612b730f8bSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetVector(output_fields[i], &vec));
962cc132f9aSnbeams     if (vec == CEED_VECTOR_ACTIVE) {
963*ca735530SJeremy L Thompson       CeedEvalMode eval_mode;
964*ca735530SJeremy L Thompson 
9652b730f8bSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetBasis(output_fields[i], &basis_out));
9662b730f8bSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetElemRestriction(output_fields[i], &rstr_out));
9676574a04fSJeremy L Thompson       CeedCheck(!rstr_out || rstr_out == rstr_in, ceed, CEED_ERROR_BACKEND, "Backend does not implement multi-field non-composite operator assembly");
9682b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode));
969cc132f9aSnbeams       if (eval_mode != CEED_EVAL_NONE) {
9702b730f8bSJeremy L Thompson         CeedCallBackend(CeedRealloc(num_B_out_mats_to_load + 1, &eval_mode_out));
971cc132f9aSnbeams         eval_mode_out[num_B_out_mats_to_load] = eval_mode;
972cc132f9aSnbeams         num_B_out_mats_to_load += 1;
973cc132f9aSnbeams         if (eval_mode == CEED_EVAL_GRAD) {
974*ca735530SJeremy L Thompson           num_e_mode_out += dim;
975*ca735530SJeremy L Thompson           size_B_out += dim * elem_size * num_qpts;
976cc132f9aSnbeams         } else {
977*ca735530SJeremy L Thompson           num_e_mode_out += 1;
978*ca735530SJeremy L Thompson           size_B_out += elem_size * num_qpts;
979cc132f9aSnbeams         }
980cc132f9aSnbeams       }
981cc132f9aSnbeams     }
982cc132f9aSnbeams   }
983*ca735530SJeremy L Thompson   CeedCheck(num_e_mode_in > 0 && num_e_mode_out > 0, ceed, CEED_ERROR_UNSUPPORTED, "Cannot assemble operator without inputs/outputs");
984cc132f9aSnbeams 
985*ca735530SJeremy L Thompson   CeedCallBackend(CeedElemRestrictionGetNumElements(rstr_in, &num_elem));
986*ca735530SJeremy L Thompson   CeedCallBackend(CeedElemRestrictionGetNumComponents(rstr_in, &num_comp));
987cc132f9aSnbeams 
9882b730f8bSJeremy L Thompson   CeedCallBackend(CeedCalloc(1, &impl->asmb));
989cc132f9aSnbeams   CeedOperatorAssemble_Cuda *asmb = impl->asmb;
990*ca735530SJeremy L Thompson   asmb->num_elem                  = num_elem;
991cc132f9aSnbeams 
992cc132f9aSnbeams   // Compile kernels
993*ca735530SJeremy L Thompson   int elem_per_block    = 1;
994*ca735530SJeremy L Thompson   asmb->elem_per_block  = elem_per_block;
995*ca735530SJeremy L Thompson   CeedInt    block_size = elem_size * elem_size * elem_per_block;
996cc132f9aSnbeams   Ceed_Cuda *cuda_data;
997*ca735530SJeremy L Thompson 
9982b730f8bSJeremy L Thompson   CeedCallBackend(CeedGetData(ceed, &cuda_data));
9992b730f8bSJeremy L Thompson   CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-operator-assemble.h", &assembly_kernel_path));
100023d4529eSJeremy L Thompson   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Assembly Kernel Source -----\n");
10012b730f8bSJeremy L Thompson   CeedCallBackend(CeedLoadSourceToBuffer(ceed, assembly_kernel_path, &assembly_kernel_source));
100223d4529eSJeremy L Thompson   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Assembly Source Complete! -----\n");
100307b31e0eSJeremy L Thompson   bool fallback = block_size > cuda_data->device_prop.maxThreadsPerBlock;
1004*ca735530SJeremy L Thompson 
100507b31e0eSJeremy L Thompson   if (fallback) {
100659ad764aSnbeams     // Use fallback kernel with 1D threadblock
1007*ca735530SJeremy L Thompson     block_size         = elem_size * elem_per_block;
1008*ca735530SJeremy L Thompson     asmb->block_size_x = elem_size;
100959ad764aSnbeams     asmb->block_size_y = 1;
101059ad764aSnbeams   } else {  // Use kernel with 2D threadblock
1011*ca735530SJeremy L Thompson     asmb->block_size_x = elem_size;
1012*ca735530SJeremy L Thompson     asmb->block_size_y = elem_size;
101307b31e0eSJeremy L Thompson   }
1014*ca735530SJeremy L Thompson   CeedCallBackend(CeedCompile_Cuda(ceed, assembly_kernel_source, &asmb->module, 8, "NUM_ELEM", num_elem, "NUM_E_MODE_IN", num_e_mode_in,
1015*ca735530SJeremy L Thompson                                    "NUM_E_MODE_OUT", num_e_mode_out, "NUM_QPTS", num_qpts, "NUM_NODES", elem_size, "BLOCK_SIZE", block_size,
1016*ca735530SJeremy L Thompson                                    "NUM_COMP", num_comp, "USE_CEEDSIZE", use_ceedsize_idx));
1017b2165e7aSSebastian Grimberg   CeedCallBackend(CeedGetKernel_Cuda(ceed, asmb->module, fallback ? "linearAssembleFallback" : "linearAssemble", &asmb->linearAssemble));
10182b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&assembly_kernel_path));
10192b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&assembly_kernel_source));
1020cc132f9aSnbeams 
1021cc132f9aSnbeams   // Build 'full' B matrices (not 1D arrays used for tensor-product matrices)
1022cc132f9aSnbeams   const CeedScalar *interp_in, *grad_in;
1023*ca735530SJeremy L Thompson 
10242b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetInterp(basis_in, &interp_in));
10252b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetGrad(basis_in, &grad_in));
1026cc132f9aSnbeams 
1027cc132f9aSnbeams   // Load into B_in, in order that they will be used in eval_mode
1028cc132f9aSnbeams   const CeedInt inBytes   = size_B_in * sizeof(CeedScalar);
1029cc132f9aSnbeams   CeedInt       mat_start = 0;
1030*ca735530SJeremy L Thompson 
10312b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMalloc((void **)&asmb->d_B_in, inBytes));
1032cc132f9aSnbeams   for (int i = 0; i < num_B_in_mats_to_load; i++) {
1033cc132f9aSnbeams     CeedEvalMode eval_mode = eval_mode_in[i];
1034*ca735530SJeremy L Thompson 
1035cc132f9aSnbeams     if (eval_mode == CEED_EVAL_INTERP) {
1036*ca735530SJeremy L Thompson       CeedCallCuda(ceed, cudaMemcpy(&asmb->d_B_in[mat_start], interp_in, elem_size * num_qpts * sizeof(CeedScalar), cudaMemcpyHostToDevice));
1037*ca735530SJeremy L Thompson       mat_start += elem_size * num_qpts;
1038cc132f9aSnbeams     } else if (eval_mode == CEED_EVAL_GRAD) {
1039*ca735530SJeremy L Thompson       CeedCallCuda(ceed, cudaMemcpy(&asmb->d_B_in[mat_start], grad_in, dim * elem_size * num_qpts * sizeof(CeedScalar), cudaMemcpyHostToDevice));
1040*ca735530SJeremy L Thompson       mat_start += dim * elem_size * num_qpts;
1041cc132f9aSnbeams     }
1042cc132f9aSnbeams   }
1043cc132f9aSnbeams 
1044cc132f9aSnbeams   const CeedScalar *interp_out, *grad_out;
1045*ca735530SJeremy L Thompson 
1046*ca735530SJeremy L Thompson   // Note that this function currently assumes 1 basis, so this should always be true for now
1047cc132f9aSnbeams   if (basis_out == basis_in) {
1048cc132f9aSnbeams     interp_out = interp_in;
1049cc132f9aSnbeams     grad_out   = grad_in;
1050cc132f9aSnbeams   } else {
10512b730f8bSJeremy L Thompson     CeedCallBackend(CeedBasisGetInterp(basis_out, &interp_out));
10522b730f8bSJeremy L Thompson     CeedCallBackend(CeedBasisGetGrad(basis_out, &grad_out));
1053cc132f9aSnbeams   }
1054cc132f9aSnbeams 
1055cc132f9aSnbeams   // Load into B_out, in order that they will be used in eval_mode
1056cc132f9aSnbeams   const CeedInt outBytes = size_B_out * sizeof(CeedScalar);
1057cc132f9aSnbeams   mat_start              = 0;
1058*ca735530SJeremy L Thompson 
10592b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cudaMalloc((void **)&asmb->d_B_out, outBytes));
1060cc132f9aSnbeams   for (int i = 0; i < num_B_out_mats_to_load; i++) {
1061cc132f9aSnbeams     CeedEvalMode eval_mode = eval_mode_out[i];
1062*ca735530SJeremy L Thompson 
1063cc132f9aSnbeams     if (eval_mode == CEED_EVAL_INTERP) {
1064*ca735530SJeremy L Thompson       CeedCallCuda(ceed, cudaMemcpy(&asmb->d_B_out[mat_start], interp_out, elem_size * num_qpts * sizeof(CeedScalar), cudaMemcpyHostToDevice));
1065*ca735530SJeremy L Thompson       mat_start += elem_size * num_qpts;
1066cc132f9aSnbeams     } else if (eval_mode == CEED_EVAL_GRAD) {
1067*ca735530SJeremy L Thompson       CeedCallCuda(ceed, cudaMemcpy(&asmb->d_B_out[mat_start], grad_out, dim * elem_size * num_qpts * sizeof(CeedScalar), cudaMemcpyHostToDevice));
1068*ca735530SJeremy L Thompson       mat_start += dim * elem_size * num_qpts;
1069cc132f9aSnbeams     }
1070cc132f9aSnbeams   }
1071cc132f9aSnbeams   return CEED_ERROR_SUCCESS;
1072cc132f9aSnbeams }
1073cc132f9aSnbeams 
1074cc132f9aSnbeams //------------------------------------------------------------------------------
1075cefa2673SJeremy L Thompson // Assemble matrix data for COO matrix of assembled operator.
1076cefa2673SJeremy L Thompson // The sparsity pattern is set by CeedOperatorLinearAssembleSymbolic.
1077cefa2673SJeremy L Thompson //
1078ea61e9acSJeremy L Thompson // Note that this (and other assembly routines) currently assume only one active input restriction/basis per operator (could have multiple basis eval
1079cefa2673SJeremy L Thompson // modes).
1080cefa2673SJeremy L Thompson // TODO: allow multiple active input restrictions/basis objects
1081cc132f9aSnbeams //------------------------------------------------------------------------------
10822b730f8bSJeremy L Thompson static int CeedSingleOperatorAssemble_Cuda(CeedOperator op, CeedInt offset, CeedVector values) {
1083cc132f9aSnbeams   Ceed                ceed;
1084*ca735530SJeremy L Thompson   CeedSize            values_length = 0, assembled_qf_length = 0;
1085*ca735530SJeremy L Thompson   CeedInt             use_ceedsize_idx = 0;
1086*ca735530SJeremy L Thompson   CeedScalar         *values_array;
1087*ca735530SJeremy L Thompson   const CeedScalar   *qf_array;
1088*ca735530SJeremy L Thompson   CeedVector          assembled_qf = NULL;
1089*ca735530SJeremy L Thompson   CeedElemRestriction rstr_q       = NULL;
1090cc132f9aSnbeams   CeedOperator_Cuda  *impl;
1091*ca735530SJeremy L Thompson 
1092*ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
10932b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &impl));
1094cc132f9aSnbeams 
1095cc132f9aSnbeams   // Assemble QFunction
10962b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorLinearAssembleQFunctionBuildOrUpdate(op, &assembled_qf, &rstr_q, CEED_REQUEST_IMMEDIATE));
10972b730f8bSJeremy L Thompson   CeedCallBackend(CeedElemRestrictionDestroy(&rstr_q));
109828ec399dSJeremy L Thompson   CeedCallBackend(CeedVectorGetArray(values, CEED_MEM_DEVICE, &values_array));
1099cc132f9aSnbeams   values_array += offset;
11002b730f8bSJeremy L Thompson   CeedCallBackend(CeedVectorGetArrayRead(assembled_qf, CEED_MEM_DEVICE, &qf_array));
1101cc132f9aSnbeams 
1102f7c1b517Snbeams   CeedCallBackend(CeedVectorGetLength(values, &values_length));
1103f7c1b517Snbeams   CeedCallBackend(CeedVectorGetLength(assembled_qf, &assembled_qf_length));
1104f7c1b517Snbeams   if ((values_length > INT_MAX) || (assembled_qf_length > INT_MAX)) use_ceedsize_idx = 1;
1105f7c1b517Snbeams   // Setup
1106f7c1b517Snbeams   if (!impl->asmb) {
1107f7c1b517Snbeams     CeedCallBackend(CeedSingleOperatorAssembleSetup_Cuda(op, use_ceedsize_idx));
1108f7c1b517Snbeams     assert(impl->asmb != NULL);
1109f7c1b517Snbeams   }
1110f7c1b517Snbeams 
1111cc132f9aSnbeams   // Compute B^T D B
1112*ca735530SJeremy L Thompson   const CeedInt num_elem       = impl->asmb->num_elem;
1113*ca735530SJeremy L Thompson   const CeedInt elem_per_block = impl->asmb->elem_per_block;
1114*ca735530SJeremy L Thompson   const CeedInt grid           = num_elem / elem_per_block + ((num_elem / elem_per_block * elem_per_block < num_elem) ? 1 : 0);
11152b730f8bSJeremy L Thompson   void         *args[]         = {&impl->asmb->d_B_in, &impl->asmb->d_B_out, &qf_array, &values_array};
1116*ca735530SJeremy L Thompson 
11172b730f8bSJeremy L Thompson   CeedCallBackend(
1118*ca735530SJeremy L Thompson       CeedRunKernelDim_Cuda(ceed, impl->asmb->linearAssemble, grid, impl->asmb->block_size_x, impl->asmb->block_size_y, elem_per_block, args));
1119cc132f9aSnbeams 
1120cc132f9aSnbeams   // Restore arrays
11212b730f8bSJeremy L Thompson   CeedCallBackend(CeedVectorRestoreArray(values, &values_array));
11222b730f8bSJeremy L Thompson   CeedCallBackend(CeedVectorRestoreArrayRead(assembled_qf, &qf_array));
1123cc132f9aSnbeams 
1124cc132f9aSnbeams   // Cleanup
11252b730f8bSJeremy L Thompson   CeedCallBackend(CeedVectorDestroy(&assembled_qf));
1126cc132f9aSnbeams   return CEED_ERROR_SUCCESS;
1127cc132f9aSnbeams }
1128cc132f9aSnbeams 
1129cc132f9aSnbeams //------------------------------------------------------------------------------
11300d0321e0SJeremy L Thompson // Create operator
11310d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
11320d0321e0SJeremy L Thompson int CeedOperatorCreate_Cuda(CeedOperator op) {
11330d0321e0SJeremy L Thompson   Ceed               ceed;
11340d0321e0SJeremy L Thompson   CeedOperator_Cuda *impl;
11350d0321e0SJeremy L Thompson 
1136*ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
11372b730f8bSJeremy L Thompson   CeedCallBackend(CeedCalloc(1, &impl));
11382b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorSetData(op, impl));
11390d0321e0SJeremy L Thompson 
11402b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleQFunction", CeedOperatorLinearAssembleQFunction_Cuda));
11412b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleQFunctionUpdate", CeedOperatorLinearAssembleQFunctionUpdate_Cuda));
11422b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleAddDiagonal", CeedOperatorLinearAssembleAddDiagonal_Cuda));
11432b730f8bSJeremy L Thompson   CeedCallBackend(
11442b730f8bSJeremy L Thompson       CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleAddPointBlockDiagonal", CeedOperatorLinearAssembleAddPointBlockDiagonal_Cuda));
11452b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleSingle", CeedSingleOperatorAssemble_Cuda));
11462b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "ApplyAdd", CeedOperatorApplyAdd_Cuda));
11472b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "Destroy", CeedOperatorDestroy_Cuda));
11480d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
11490d0321e0SJeremy L Thompson }
11500d0321e0SJeremy L Thompson 
11510d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
1152