xref: /libCEED/rust/libceed-sys/c-src/backends/cuda-ref/ceed-cuda-ref-operator.c (revision c1222711bbaa154065e44a2caf05c765092ea957)
15aed82e4SJeremy L Thompson // Copyright (c) 2017-2024, 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;
26ca735530SJeremy L Thompson 
272b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &impl));
280d0321e0SJeremy L Thompson 
290d0321e0SJeremy L Thompson   // Apply data
30ca735530SJeremy L Thompson   for (CeedInt i = 0; i < impl->num_inputs + impl->num_outputs; i++) {
31ca735530SJeremy L Thompson     CeedCallBackend(CeedVectorDestroy(&impl->e_vecs[i]));
320d0321e0SJeremy L Thompson   }
33ca735530SJeremy L Thompson   CeedCallBackend(CeedFree(&impl->e_vecs));
34*c1222711SJeremy L Thompson   CeedCallBackend(CeedFree(&impl->input_states));
350d0321e0SJeremy L Thompson 
36ca735530SJeremy L Thompson   for (CeedInt i = 0; i < impl->num_inputs; i++) {
37ca735530SJeremy L Thompson     CeedCallBackend(CeedVectorDestroy(&impl->q_vecs_in[i]));
380d0321e0SJeremy L Thompson   }
39ca735530SJeremy L Thompson   CeedCallBackend(CeedFree(&impl->q_vecs_in));
400d0321e0SJeremy L Thompson 
41ca735530SJeremy L Thompson   for (CeedInt i = 0; i < impl->num_outputs; i++) {
42ca735530SJeremy L Thompson     CeedCallBackend(CeedVectorDestroy(&impl->q_vecs_out[i]));
430d0321e0SJeremy L Thompson   }
44ca735530SJeremy L Thompson   CeedCallBackend(CeedFree(&impl->q_vecs_out));
450d0321e0SJeremy L Thompson 
460d0321e0SJeremy L Thompson   // QFunction assembly data
47ca735530SJeremy L Thompson   for (CeedInt i = 0; i < impl->num_active_in; i++) {
48ca735530SJeremy L Thompson     CeedCallBackend(CeedVectorDestroy(&impl->qf_active_in[i]));
490d0321e0SJeremy L Thompson   }
50ca735530SJeremy L Thompson   CeedCallBackend(CeedFree(&impl->qf_active_in));
510d0321e0SJeremy L Thompson 
520d0321e0SJeremy L Thompson   // Diag data
530d0321e0SJeremy L Thompson   if (impl->diag) {
540d0321e0SJeremy L Thompson     Ceed ceed;
55ca735530SJeremy L Thompson 
562b730f8bSJeremy L Thompson     CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
57cbfe683aSSebastian Grimberg     if (impl->diag->module) {
582b730f8bSJeremy L Thompson       CeedCallCuda(ceed, cuModuleUnload(impl->diag->module));
59cbfe683aSSebastian Grimberg     }
60cbfe683aSSebastian Grimberg     if (impl->diag->module_point_block) {
61cbfe683aSSebastian Grimberg       CeedCallCuda(ceed, cuModuleUnload(impl->diag->module_point_block));
62cbfe683aSSebastian Grimberg     }
63004e4986SSebastian Grimberg     CeedCallCuda(ceed, cudaFree(impl->diag->d_eval_modes_in));
64004e4986SSebastian Grimberg     CeedCallCuda(ceed, cudaFree(impl->diag->d_eval_modes_out));
652b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaFree(impl->diag->d_identity));
66ca735530SJeremy L Thompson     CeedCallCuda(ceed, cudaFree(impl->diag->d_interp_in));
67ca735530SJeremy L Thompson     CeedCallCuda(ceed, cudaFree(impl->diag->d_interp_out));
68ca735530SJeremy L Thompson     CeedCallCuda(ceed, cudaFree(impl->diag->d_grad_in));
69ca735530SJeremy L Thompson     CeedCallCuda(ceed, cudaFree(impl->diag->d_grad_out));
70004e4986SSebastian Grimberg     CeedCallCuda(ceed, cudaFree(impl->diag->d_div_in));
71004e4986SSebastian Grimberg     CeedCallCuda(ceed, cudaFree(impl->diag->d_div_out));
72004e4986SSebastian Grimberg     CeedCallCuda(ceed, cudaFree(impl->diag->d_curl_in));
73004e4986SSebastian Grimberg     CeedCallCuda(ceed, cudaFree(impl->diag->d_curl_out));
74004e4986SSebastian Grimberg     CeedCallBackend(CeedElemRestrictionDestroy(&impl->diag->diag_rstr));
75506b1a0cSSebastian Grimberg     CeedCallBackend(CeedElemRestrictionDestroy(&impl->diag->point_block_diag_rstr));
76ca735530SJeremy L Thompson     CeedCallBackend(CeedVectorDestroy(&impl->diag->elem_diag));
77ca735530SJeremy L Thompson     CeedCallBackend(CeedVectorDestroy(&impl->diag->point_block_elem_diag));
780d0321e0SJeremy L Thompson   }
792b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&impl->diag));
800d0321e0SJeremy L Thompson 
81cc132f9aSnbeams   if (impl->asmb) {
82cc132f9aSnbeams     Ceed ceed;
83ca735530SJeremy L Thompson 
842b730f8bSJeremy L Thompson     CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
852b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cuModuleUnload(impl->asmb->module));
862b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaFree(impl->asmb->d_B_in));
872b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaFree(impl->asmb->d_B_out));
88cc132f9aSnbeams   }
892b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&impl->asmb));
90cc132f9aSnbeams 
912b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&impl));
920d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
930d0321e0SJeremy L Thompson }
940d0321e0SJeremy L Thompson 
950d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
960d0321e0SJeremy L Thompson // Setup infields or outfields
970d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
98004e4986SSebastian Grimberg static int CeedOperatorSetupFields_Cuda(CeedQFunction qf, CeedOperator op, bool is_input, CeedVector *e_vecs, CeedVector *q_vecs, CeedInt start_e,
99ca735530SJeremy L Thompson                                         CeedInt num_fields, CeedInt Q, CeedInt num_elem) {
1000d0321e0SJeremy L Thompson   Ceed                ceed;
101ca735530SJeremy L Thompson   CeedQFunctionField *qf_fields;
102ca735530SJeremy L Thompson   CeedOperatorField  *op_fields;
1030d0321e0SJeremy L Thompson 
104ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
105ca735530SJeremy L Thompson   if (is_input) {
106ca735530SJeremy L Thompson     CeedCallBackend(CeedOperatorGetFields(op, NULL, &op_fields, NULL, NULL));
107ca735530SJeremy L Thompson     CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_fields, NULL, NULL));
1080d0321e0SJeremy L Thompson   } else {
109ca735530SJeremy L Thompson     CeedCallBackend(CeedOperatorGetFields(op, NULL, NULL, NULL, &op_fields));
110ca735530SJeremy L Thompson     CeedCallBackend(CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qf_fields));
1110d0321e0SJeremy L Thompson   }
1120d0321e0SJeremy L Thompson 
1130d0321e0SJeremy L Thompson   // Loop over fields
114ca735530SJeremy L Thompson   for (CeedInt i = 0; i < num_fields; i++) {
115004e4986SSebastian Grimberg     bool         is_strided = false, skip_restriction = false;
116004e4986SSebastian Grimberg     CeedSize     q_size;
117004e4986SSebastian Grimberg     CeedInt      size;
118004e4986SSebastian Grimberg     CeedEvalMode eval_mode;
119ca735530SJeremy L Thompson     CeedBasis    basis;
1200d0321e0SJeremy L Thompson 
121004e4986SSebastian Grimberg     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode));
122004e4986SSebastian Grimberg     if (eval_mode != CEED_EVAL_WEIGHT) {
123edb2538eSJeremy L Thompson       CeedElemRestriction elem_rstr;
124ca735530SJeremy L Thompson 
1250d0321e0SJeremy L Thompson       // Check whether this field can skip the element restriction:
126004e4986SSebastian Grimberg       // Must be passive input, with eval_mode NONE, and have a strided restriction with CEED_STRIDES_BACKEND.
127004e4986SSebastian Grimberg       CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_fields[i], &elem_rstr));
1280d0321e0SJeremy L Thompson 
1290d0321e0SJeremy L Thompson       // First, check whether the field is input or output:
130ca735530SJeremy L Thompson       if (is_input) {
131ca735530SJeremy L Thompson         CeedVector vec;
132ca735530SJeremy L Thompson 
133004e4986SSebastian Grimberg         // Check for passive input
134ca735530SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec));
135ca735530SJeremy L Thompson         if (vec != CEED_VECTOR_ACTIVE) {
136004e4986SSebastian Grimberg           // Check eval_mode
137004e4986SSebastian Grimberg           if (eval_mode == CEED_EVAL_NONE) {
1380d0321e0SJeremy L Thompson             // Check for strided restriction
139edb2538eSJeremy L Thompson             CeedCallBackend(CeedElemRestrictionIsStrided(elem_rstr, &is_strided));
140ca735530SJeremy L Thompson             if (is_strided) {
1410d0321e0SJeremy L Thompson               // Check if vector is already in preferred backend ordering
142edb2538eSJeremy L Thompson               CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &skip_restriction));
1430d0321e0SJeremy L Thompson             }
1440d0321e0SJeremy L Thompson           }
1450d0321e0SJeremy L Thompson         }
1460d0321e0SJeremy L Thompson       }
147ca735530SJeremy L Thompson       if (skip_restriction) {
148ea61e9acSJeremy L Thompson         // We do not need an E-Vector, but will use the input field vector's data directly in the operator application.
149004e4986SSebastian Grimberg         e_vecs[i + start_e] = NULL;
1500d0321e0SJeremy L Thompson       } else {
151004e4986SSebastian Grimberg         CeedCallBackend(CeedElemRestrictionCreateVector(elem_rstr, NULL, &e_vecs[i + start_e]));
1520d0321e0SJeremy L Thompson       }
1530d0321e0SJeremy L Thompson     }
1540d0321e0SJeremy L Thompson 
155004e4986SSebastian Grimberg     switch (eval_mode) {
1560d0321e0SJeremy L Thompson       case CEED_EVAL_NONE:
157ca735530SJeremy L Thompson         CeedCallBackend(CeedQFunctionFieldGetSize(qf_fields[i], &size));
158ca735530SJeremy L Thompson         q_size = (CeedSize)num_elem * Q * size;
159ca735530SJeremy L Thompson         CeedCallBackend(CeedVectorCreate(ceed, q_size, &q_vecs[i]));
1600d0321e0SJeremy L Thompson         break;
1610d0321e0SJeremy L Thompson       case CEED_EVAL_INTERP:
1620d0321e0SJeremy L Thompson       case CEED_EVAL_GRAD:
163004e4986SSebastian Grimberg       case CEED_EVAL_DIV:
164004e4986SSebastian Grimberg       case CEED_EVAL_CURL:
165ca735530SJeremy L Thompson         CeedCallBackend(CeedQFunctionFieldGetSize(qf_fields[i], &size));
166ca735530SJeremy L Thompson         q_size = (CeedSize)num_elem * Q * size;
167ca735530SJeremy 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
170ca735530SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis));
171ca735530SJeremy L Thompson         q_size = (CeedSize)num_elem * Q;
172ca735530SJeremy L Thompson         CeedCallBackend(CeedVectorCreate(ceed, q_size, &q_vecs[i]));
173ca735530SJeremy 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     }
1760d0321e0SJeremy L Thompson   }
1770d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1780d0321e0SJeremy L Thompson }
1790d0321e0SJeremy L Thompson 
1800d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
181ea61e9acSJeremy L Thompson // CeedOperator needs to connect all the named fields (be they active or passive) to the named inputs and outputs of its CeedQFunction.
1820d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
1830d0321e0SJeremy L Thompson static int CeedOperatorSetup_Cuda(CeedOperator op) {
1840d0321e0SJeremy L Thompson   Ceed                ceed;
185ca735530SJeremy L Thompson   bool                is_setup_done;
186ca735530SJeremy L Thompson   CeedInt             Q, num_elem, num_input_fields, num_output_fields;
187ca735530SJeremy L Thompson   CeedQFunctionField *qf_input_fields, *qf_output_fields;
1880d0321e0SJeremy L Thompson   CeedQFunction       qf;
189ca735530SJeremy L Thompson   CeedOperatorField  *op_input_fields, *op_output_fields;
190ca735530SJeremy L Thompson   CeedOperator_Cuda  *impl;
191ca735530SJeremy L Thompson 
192ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorIsSetupDone(op, &is_setup_done));
193ca735530SJeremy L Thompson   if (is_setup_done) return CEED_ERROR_SUCCESS;
194ca735530SJeremy L Thompson 
195ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
196ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &impl));
1972b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
1982b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q));
199ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem));
200ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
201ca735530SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
2020d0321e0SJeremy L Thompson 
2030d0321e0SJeremy L Thompson   // Allocate
204ca735530SJeremy L Thompson   CeedCallBackend(CeedCalloc(num_input_fields + num_output_fields, &impl->e_vecs));
205*c1222711SJeremy L Thompson   CeedCallBackend(CeedCalloc(CEED_FIELD_MAX, &impl->input_states));
206ca735530SJeremy L Thompson   CeedCallBackend(CeedCalloc(CEED_FIELD_MAX, &impl->q_vecs_in));
207ca735530SJeremy L Thompson   CeedCallBackend(CeedCalloc(CEED_FIELD_MAX, &impl->q_vecs_out));
208ca735530SJeremy L Thompson   impl->num_inputs  = num_input_fields;
209ca735530SJeremy L Thompson   impl->num_outputs = num_output_fields;
2100d0321e0SJeremy L Thompson 
211ca735530SJeremy L Thompson   // Set up infield and outfield e_vecs and q_vecs
2120d0321e0SJeremy L Thompson   // Infields
213ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorSetupFields_Cuda(qf, op, true, impl->e_vecs, impl->q_vecs_in, 0, num_input_fields, Q, num_elem));
2140d0321e0SJeremy L Thompson   // Outfields
215ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorSetupFields_Cuda(qf, op, false, impl->e_vecs, impl->q_vecs_out, num_input_fields, num_output_fields, Q, num_elem));
2160d0321e0SJeremy L Thompson 
2172b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorSetSetupDone(op));
2180d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2190d0321e0SJeremy L Thompson }
2200d0321e0SJeremy L Thompson 
2210d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
2220d0321e0SJeremy L Thompson // Setup Operator Inputs
2230d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
224ca735530SJeremy L Thompson static inline int CeedOperatorSetupInputs_Cuda(CeedInt num_input_fields, CeedQFunctionField *qf_input_fields, CeedOperatorField *op_input_fields,
225004e4986SSebastian Grimberg                                                CeedVector in_vec, const bool skip_active, CeedScalar *e_data[2 * CEED_FIELD_MAX],
2260d0321e0SJeremy L Thompson                                                CeedOperator_Cuda *impl, CeedRequest *request) {
227ca735530SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
228004e4986SSebastian Grimberg     CeedEvalMode        eval_mode;
2290d0321e0SJeremy L Thompson     CeedVector          vec;
230edb2538eSJeremy L Thompson     CeedElemRestriction elem_rstr;
2310d0321e0SJeremy L Thompson 
2320d0321e0SJeremy L Thompson     // Get input vector
233ca735530SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec));
2340d0321e0SJeremy L Thompson     if (vec == CEED_VECTOR_ACTIVE) {
235004e4986SSebastian Grimberg       if (skip_active) continue;
236ca735530SJeremy L Thompson       else vec = in_vec;
2370d0321e0SJeremy L Thompson     }
2380d0321e0SJeremy L Thompson 
239004e4986SSebastian Grimberg     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
240004e4986SSebastian Grimberg     if (eval_mode == CEED_EVAL_WEIGHT) {  // Skip
2410d0321e0SJeremy L Thompson     } else {
242004e4986SSebastian Grimberg       // Get input vector
243004e4986SSebastian Grimberg       CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec));
2440d0321e0SJeremy L Thompson       // Get input element restriction
245edb2538eSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr));
246ca735530SJeremy L Thompson       if (vec == CEED_VECTOR_ACTIVE) vec = in_vec;
2470d0321e0SJeremy L Thompson       // Restrict, if necessary
248ca735530SJeremy L Thompson       if (!impl->e_vecs[i]) {
2490d0321e0SJeremy L Thompson         // No restriction for this field; read data directly from vec.
250ca735530SJeremy L Thompson         CeedCallBackend(CeedVectorGetArrayRead(vec, CEED_MEM_DEVICE, (const CeedScalar **)&e_data[i]));
2510d0321e0SJeremy L Thompson       } else {
252*c1222711SJeremy L Thompson         uint64_t state;
253*c1222711SJeremy L Thompson 
254*c1222711SJeremy L Thompson         CeedCallBackend(CeedVectorGetState(vec, &state));
255*c1222711SJeremy L Thompson         if (state != impl->input_states[i]) {
256edb2538eSJeremy L Thompson           CeedCallBackend(CeedElemRestrictionApply(elem_rstr, CEED_NOTRANSPOSE, vec, impl->e_vecs[i], request));
257*c1222711SJeremy L Thompson           impl->input_states[i] = state;
258*c1222711SJeremy L Thompson         }
2590d0321e0SJeremy L Thompson         // Get evec
260ca735530SJeremy L Thompson         CeedCallBackend(CeedVectorGetArrayRead(impl->e_vecs[i], CEED_MEM_DEVICE, (const CeedScalar **)&e_data[i]));
2610d0321e0SJeremy L Thompson       }
2620d0321e0SJeremy L Thompson     }
2630d0321e0SJeremy L Thompson   }
2640d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2650d0321e0SJeremy L Thompson }
2660d0321e0SJeremy L Thompson 
2670d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
2680d0321e0SJeremy L Thompson // Input Basis Action
2690d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
270ca735530SJeremy L Thompson static inline int CeedOperatorInputBasis_Cuda(CeedInt num_elem, CeedQFunctionField *qf_input_fields, CeedOperatorField *op_input_fields,
271004e4986SSebastian Grimberg                                               CeedInt num_input_fields, const bool skip_active, CeedScalar *e_data[2 * CEED_FIELD_MAX],
2722b730f8bSJeremy L Thompson                                               CeedOperator_Cuda *impl) {
273ca735530SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
274ca735530SJeremy L Thompson     CeedInt             elem_size, size;
275004e4986SSebastian Grimberg     CeedEvalMode        eval_mode;
276edb2538eSJeremy L Thompson     CeedElemRestriction elem_rstr;
2770d0321e0SJeremy L Thompson     CeedBasis           basis;
2780d0321e0SJeremy L Thompson 
2790d0321e0SJeremy L Thompson     // Skip active input
280004e4986SSebastian Grimberg     if (skip_active) {
2810d0321e0SJeremy L Thompson       CeedVector vec;
282ca735530SJeremy L Thompson 
283ca735530SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec));
2842b730f8bSJeremy L Thompson       if (vec == CEED_VECTOR_ACTIVE) continue;
2850d0321e0SJeremy L Thompson     }
286004e4986SSebastian Grimberg     // Get elem_size, eval_mode, size
287edb2538eSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr));
288edb2538eSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
289004e4986SSebastian Grimberg     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
290ca735530SJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetSize(qf_input_fields[i], &size));
2910d0321e0SJeremy L Thompson     // Basis action
292004e4986SSebastian Grimberg     switch (eval_mode) {
2930d0321e0SJeremy L Thompson       case CEED_EVAL_NONE:
294ca735530SJeremy L Thompson         CeedCallBackend(CeedVectorSetArray(impl->q_vecs_in[i], CEED_MEM_DEVICE, CEED_USE_POINTER, e_data[i]));
2950d0321e0SJeremy L Thompson         break;
2960d0321e0SJeremy L Thompson       case CEED_EVAL_INTERP:
2970d0321e0SJeremy L Thompson       case CEED_EVAL_GRAD:
298004e4986SSebastian Grimberg       case CEED_EVAL_DIV:
299004e4986SSebastian Grimberg       case CEED_EVAL_CURL:
300ca735530SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
301004e4986SSebastian Grimberg         CeedCallBackend(CeedBasisApply(basis, num_elem, CEED_NOTRANSPOSE, eval_mode, impl->e_vecs[i], impl->q_vecs_in[i]));
3020d0321e0SJeremy L Thompson         break;
3030d0321e0SJeremy L Thompson       case CEED_EVAL_WEIGHT:
3040d0321e0SJeremy L Thompson         break;  // No action
3050d0321e0SJeremy L Thompson     }
3060d0321e0SJeremy L Thompson   }
3070d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3080d0321e0SJeremy L Thompson }
3090d0321e0SJeremy L Thompson 
3100d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
3110d0321e0SJeremy L Thompson // Restore Input Vectors
3120d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
313ca735530SJeremy L Thompson static inline int CeedOperatorRestoreInputs_Cuda(CeedInt num_input_fields, CeedQFunctionField *qf_input_fields, CeedOperatorField *op_input_fields,
314004e4986SSebastian Grimberg                                                  const bool skip_active, CeedScalar *e_data[2 * CEED_FIELD_MAX], CeedOperator_Cuda *impl) {
315ca735530SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
316004e4986SSebastian Grimberg     CeedEvalMode eval_mode;
3170d0321e0SJeremy L Thompson     CeedVector   vec;
3180d0321e0SJeremy L Thompson 
3190d0321e0SJeremy L Thompson     // Skip active input
320004e4986SSebastian Grimberg     if (skip_active) {
321ca735530SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec));
3222b730f8bSJeremy L Thompson       if (vec == CEED_VECTOR_ACTIVE) continue;
3230d0321e0SJeremy L Thompson     }
324004e4986SSebastian Grimberg     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
325004e4986SSebastian Grimberg     if (eval_mode == CEED_EVAL_WEIGHT) {  // Skip
3260d0321e0SJeremy L Thompson     } else {
327ca735530SJeremy L Thompson       if (!impl->e_vecs[i]) {  // This was a skip_restriction case
328ca735530SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec));
329ca735530SJeremy L Thompson         CeedCallBackend(CeedVectorRestoreArrayRead(vec, (const CeedScalar **)&e_data[i]));
3300d0321e0SJeremy L Thompson       } else {
331ca735530SJeremy L Thompson         CeedCallBackend(CeedVectorRestoreArrayRead(impl->e_vecs[i], (const CeedScalar **)&e_data[i]));
3320d0321e0SJeremy L Thompson       }
3330d0321e0SJeremy L Thompson     }
3340d0321e0SJeremy L Thompson   }
3350d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3360d0321e0SJeremy L Thompson }
3370d0321e0SJeremy L Thompson 
3380d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
3390d0321e0SJeremy L Thompson // Apply and add to output
3400d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
341ca735530SJeremy L Thompson static int CeedOperatorApplyAdd_Cuda(CeedOperator op, CeedVector in_vec, CeedVector out_vec, CeedRequest *request) {
342ca735530SJeremy L Thompson   CeedInt             Q, num_elem, elem_size, num_input_fields, num_output_fields, size;
343ca735530SJeremy L Thompson   CeedScalar         *e_data[2 * CEED_FIELD_MAX] = {NULL};
344ca735530SJeremy L Thompson   CeedQFunctionField *qf_input_fields, *qf_output_fields;
3450d0321e0SJeremy L Thompson   CeedQFunction       qf;
346004e4986SSebastian Grimberg   CeedOperatorField  *op_input_fields, *op_output_fields;
347004e4986SSebastian Grimberg   CeedOperator_Cuda  *impl;
348ca735530SJeremy L Thompson 
349ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &impl));
3502b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
3512b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q));
352ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem));
353ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
354ca735530SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
3550d0321e0SJeremy L Thompson 
3560d0321e0SJeremy L Thompson   // Setup
3572b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorSetup_Cuda(op));
3580d0321e0SJeremy L Thompson 
359004e4986SSebastian Grimberg   // Input Evecs and Restriction
360ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorSetupInputs_Cuda(num_input_fields, qf_input_fields, op_input_fields, in_vec, false, e_data, impl, request));
3610d0321e0SJeremy L Thompson 
3620d0321e0SJeremy L Thompson   // Input basis apply if needed
363ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorInputBasis_Cuda(num_elem, qf_input_fields, op_input_fields, num_input_fields, false, e_data, impl));
3640d0321e0SJeremy L Thompson 
3650d0321e0SJeremy L Thompson   // Output pointers, as necessary
366ca735530SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
367004e4986SSebastian Grimberg     CeedEvalMode eval_mode;
368004e4986SSebastian Grimberg 
369004e4986SSebastian Grimberg     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
370004e4986SSebastian Grimberg     if (eval_mode == CEED_EVAL_NONE) {
3710d0321e0SJeremy L Thompson       // Set the output Q-Vector to use the E-Vector data directly.
372ca735530SJeremy L Thompson       CeedCallBackend(CeedVectorGetArrayWrite(impl->e_vecs[i + impl->num_inputs], CEED_MEM_DEVICE, &e_data[i + num_input_fields]));
373ca735530SJeremy L Thompson       CeedCallBackend(CeedVectorSetArray(impl->q_vecs_out[i], CEED_MEM_DEVICE, CEED_USE_POINTER, e_data[i + num_input_fields]));
3740d0321e0SJeremy L Thompson     }
3750d0321e0SJeremy L Thompson   }
3760d0321e0SJeremy L Thompson 
3770d0321e0SJeremy L Thompson   // Q function
378ca735530SJeremy L Thompson   CeedCallBackend(CeedQFunctionApply(qf, num_elem * Q, impl->q_vecs_in, impl->q_vecs_out));
3790d0321e0SJeremy L Thompson 
3800d0321e0SJeremy L Thompson   // Output basis apply if needed
381ca735530SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
382004e4986SSebastian Grimberg     CeedEvalMode        eval_mode;
383edb2538eSJeremy L Thompson     CeedElemRestriction elem_rstr;
384ca735530SJeremy L Thompson     CeedBasis           basis;
385ca735530SJeremy L Thompson 
386004e4986SSebastian Grimberg     // Get elem_size, eval_mode, size
387edb2538eSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &elem_rstr));
388edb2538eSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
389004e4986SSebastian Grimberg     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
390ca735530SJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetSize(qf_output_fields[i], &size));
3910d0321e0SJeremy L Thompson     // Basis action
392004e4986SSebastian Grimberg     switch (eval_mode) {
3930d0321e0SJeremy L Thompson       case CEED_EVAL_NONE:
394004e4986SSebastian Grimberg         break;  // No action
3950d0321e0SJeremy L Thompson       case CEED_EVAL_INTERP:
3960d0321e0SJeremy L Thompson       case CEED_EVAL_GRAD:
397004e4986SSebastian Grimberg       case CEED_EVAL_DIV:
398004e4986SSebastian Grimberg       case CEED_EVAL_CURL:
399ca735530SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
400004e4986SSebastian Grimberg         CeedCallBackend(CeedBasisApply(basis, num_elem, CEED_TRANSPOSE, eval_mode, impl->q_vecs_out[i], impl->e_vecs[i + impl->num_inputs]));
4010d0321e0SJeremy L Thompson         break;
4020d0321e0SJeremy L Thompson       // LCOV_EXCL_START
4030d0321e0SJeremy L Thompson       case CEED_EVAL_WEIGHT: {
4046e536b99SJeremy L Thompson         return CeedError(CeedOperatorReturnCeed(op), CEED_ERROR_BACKEND, "CEED_EVAL_WEIGHT cannot be an output evaluation mode");
4050d0321e0SJeremy L Thompson         // LCOV_EXCL_STOP
4060d0321e0SJeremy L Thompson       }
4070d0321e0SJeremy L Thompson     }
408004e4986SSebastian Grimberg   }
4090d0321e0SJeremy L Thompson 
4100d0321e0SJeremy L Thompson   // Output restriction
411ca735530SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
412004e4986SSebastian Grimberg     CeedEvalMode        eval_mode;
413ca735530SJeremy L Thompson     CeedVector          vec;
414edb2538eSJeremy L Thompson     CeedElemRestriction elem_rstr;
415ca735530SJeremy L Thompson 
4160d0321e0SJeremy L Thompson     // Restore evec
417004e4986SSebastian Grimberg     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
418004e4986SSebastian Grimberg     if (eval_mode == CEED_EVAL_NONE) {
419ca735530SJeremy L Thompson       CeedCallBackend(CeedVectorRestoreArray(impl->e_vecs[i + impl->num_inputs], &e_data[i + num_input_fields]));
4200d0321e0SJeremy L Thompson     }
4210d0321e0SJeremy L Thompson     // Get output vector
422ca735530SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[i], &vec));
4230d0321e0SJeremy L Thompson     // Restrict
424edb2538eSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &elem_rstr));
4250d0321e0SJeremy L Thompson     // Active
426ca735530SJeremy L Thompson     if (vec == CEED_VECTOR_ACTIVE) vec = out_vec;
4270d0321e0SJeremy L Thompson 
428edb2538eSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionApply(elem_rstr, CEED_TRANSPOSE, impl->e_vecs[i + impl->num_inputs], vec, request));
4290d0321e0SJeremy L Thompson   }
4300d0321e0SJeremy L Thompson 
4310d0321e0SJeremy L Thompson   // Restore input arrays
432ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorRestoreInputs_Cuda(num_input_fields, qf_input_fields, op_input_fields, false, e_data, impl));
4330d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
4340d0321e0SJeremy L Thompson }
4350d0321e0SJeremy L Thompson 
4360d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
437004e4986SSebastian Grimberg // Linear QFunction Assembly Core
4380d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
4392b730f8bSJeremy L Thompson static inline int CeedOperatorLinearAssembleQFunctionCore_Cuda(CeedOperator op, bool build_objects, CeedVector *assembled, CeedElemRestriction *rstr,
4400d0321e0SJeremy L Thompson                                                                CeedRequest *request) {
441ca735530SJeremy L Thompson   Ceed                ceed, ceed_parent;
442ca735530SJeremy L Thompson   CeedInt             num_active_in, num_active_out, Q, num_elem, num_input_fields, num_output_fields, size;
443ca735530SJeremy L Thompson   CeedScalar         *assembled_array, *e_data[2 * CEED_FIELD_MAX] = {NULL};
444ca735530SJeremy L Thompson   CeedVector         *active_inputs;
445ca735530SJeremy L Thompson   CeedQFunctionField *qf_input_fields, *qf_output_fields;
446ca735530SJeremy L Thompson   CeedQFunction       qf;
447ca735530SJeremy L Thompson   CeedOperatorField  *op_input_fields, *op_output_fields;
448ca735530SJeremy L Thompson   CeedOperator_Cuda  *impl;
449ca735530SJeremy L Thompson 
4502b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
451ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetFallbackParentCeed(op, &ceed_parent));
452e984cf9aSJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &impl));
453e984cf9aSJeremy L Thompson   CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q));
454ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem));
455e984cf9aSJeremy L Thompson   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
456ca735530SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
457ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
458ca735530SJeremy L Thompson   active_inputs = impl->qf_active_in;
459ca735530SJeremy L Thompson   num_active_in = impl->num_active_in, num_active_out = impl->num_active_out;
4600d0321e0SJeremy L Thompson 
4610d0321e0SJeremy L Thompson   // Setup
4622b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorSetup_Cuda(op));
4630d0321e0SJeremy L Thompson 
464004e4986SSebastian Grimberg   // Input Evecs and Restriction
465ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorSetupInputs_Cuda(num_input_fields, qf_input_fields, op_input_fields, NULL, true, e_data, impl, request));
4660d0321e0SJeremy L Thompson 
4670d0321e0SJeremy L Thompson   // Count number of active input fields
468ca735530SJeremy L Thompson   if (!num_active_in) {
469ca735530SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
470ca735530SJeremy L Thompson       CeedScalar *q_vec_array;
471ca735530SJeremy L Thompson       CeedVector  vec;
472ca735530SJeremy L Thompson 
4730d0321e0SJeremy L Thompson       // Get input vector
474ca735530SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec));
4750d0321e0SJeremy L Thompson       // Check if active input
4760d0321e0SJeremy L Thompson       if (vec == CEED_VECTOR_ACTIVE) {
477ca735530SJeremy L Thompson         CeedCallBackend(CeedQFunctionFieldGetSize(qf_input_fields[i], &size));
478ca735530SJeremy L Thompson         CeedCallBackend(CeedVectorSetValue(impl->q_vecs_in[i], 0.0));
479ca735530SJeremy L Thompson         CeedCallBackend(CeedVectorGetArray(impl->q_vecs_in[i], CEED_MEM_DEVICE, &q_vec_array));
480ca735530SJeremy L Thompson         CeedCallBackend(CeedRealloc(num_active_in + size, &active_inputs));
4810d0321e0SJeremy L Thompson         for (CeedInt field = 0; field < size; field++) {
482004e4986SSebastian Grimberg           CeedSize q_size = (CeedSize)Q * num_elem;
483004e4986SSebastian Grimberg 
484ca735530SJeremy L Thompson           CeedCallBackend(CeedVectorCreate(ceed, q_size, &active_inputs[num_active_in + field]));
485ca735530SJeremy L Thompson           CeedCallBackend(
486ca735530SJeremy L Thompson               CeedVectorSetArray(active_inputs[num_active_in + field], CEED_MEM_DEVICE, CEED_USE_POINTER, &q_vec_array[field * Q * num_elem]));
4870d0321e0SJeremy L Thompson         }
488ca735530SJeremy L Thompson         num_active_in += size;
489ca735530SJeremy L Thompson         CeedCallBackend(CeedVectorRestoreArray(impl->q_vecs_in[i], &q_vec_array));
4900d0321e0SJeremy L Thompson       }
4910d0321e0SJeremy L Thompson     }
492ca735530SJeremy L Thompson     impl->num_active_in = num_active_in;
493ca735530SJeremy L Thompson     impl->qf_active_in  = active_inputs;
4940d0321e0SJeremy L Thompson   }
4950d0321e0SJeremy L Thompson 
4960d0321e0SJeremy L Thompson   // Count number of active output fields
497ca735530SJeremy L Thompson   if (!num_active_out) {
498ca735530SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
499ca735530SJeremy L Thompson       CeedVector vec;
500ca735530SJeremy L Thompson 
5010d0321e0SJeremy L Thompson       // Get output vector
502ca735530SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[i], &vec));
5030d0321e0SJeremy L Thompson       // Check if active output
5040d0321e0SJeremy L Thompson       if (vec == CEED_VECTOR_ACTIVE) {
505ca735530SJeremy L Thompson         CeedCallBackend(CeedQFunctionFieldGetSize(qf_output_fields[i], &size));
506ca735530SJeremy L Thompson         num_active_out += size;
5070d0321e0SJeremy L Thompson       }
5080d0321e0SJeremy L Thompson     }
509ca735530SJeremy L Thompson     impl->num_active_out = num_active_out;
5100d0321e0SJeremy L Thompson   }
5110d0321e0SJeremy L Thompson 
5120d0321e0SJeremy L Thompson   // Check sizes
513ca735530SJeremy L Thompson   CeedCheck(num_active_in > 0 && num_active_out > 0, ceed, CEED_ERROR_BACKEND, "Cannot assemble QFunction without active inputs and outputs");
5140d0321e0SJeremy L Thompson 
5150d0321e0SJeremy L Thompson   // Build objects if needed
5160d0321e0SJeremy L Thompson   if (build_objects) {
517004e4986SSebastian Grimberg     CeedSize l_size     = (CeedSize)num_elem * Q * num_active_in * num_active_out;
518ca735530SJeremy L Thompson     CeedInt  strides[3] = {1, num_elem * Q, Q}; /* *NOPAD* */
519004e4986SSebastian Grimberg 
520004e4986SSebastian Grimberg     // Create output restriction
521ca735530SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionCreateStrided(ceed_parent, num_elem, Q, num_active_in * num_active_out,
522ca735530SJeremy L Thompson                                                      num_active_in * num_active_out * num_elem * Q, strides, rstr));
5230d0321e0SJeremy L Thompson     // Create assembled vector
524ca735530SJeremy L Thompson     CeedCallBackend(CeedVectorCreate(ceed_parent, l_size, assembled));
5250d0321e0SJeremy L Thompson   }
5262b730f8bSJeremy L Thompson   CeedCallBackend(CeedVectorSetValue(*assembled, 0.0));
527ca735530SJeremy L Thompson   CeedCallBackend(CeedVectorGetArray(*assembled, CEED_MEM_DEVICE, &assembled_array));
5280d0321e0SJeremy L Thompson 
5290d0321e0SJeremy L Thompson   // Input basis apply
530ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorInputBasis_Cuda(num_elem, qf_input_fields, op_input_fields, num_input_fields, true, e_data, impl));
5310d0321e0SJeremy L Thompson 
5320d0321e0SJeremy L Thompson   // Assemble QFunction
533ca735530SJeremy L Thompson   for (CeedInt in = 0; in < num_active_in; in++) {
5340d0321e0SJeremy L Thompson     // Set Inputs
535ca735530SJeremy L Thompson     CeedCallBackend(CeedVectorSetValue(active_inputs[in], 1.0));
536ca735530SJeremy L Thompson     if (num_active_in > 1) {
537ca735530SJeremy L Thompson       CeedCallBackend(CeedVectorSetValue(active_inputs[(in + num_active_in - 1) % num_active_in], 0.0));
5380d0321e0SJeremy L Thompson     }
5390d0321e0SJeremy L Thompson     // Set Outputs
540ca735530SJeremy L Thompson     for (CeedInt out = 0; out < num_output_fields; out++) {
541ca735530SJeremy L Thompson       CeedVector vec;
542ca735530SJeremy L Thompson 
5430d0321e0SJeremy L Thompson       // Get output vector
544ca735530SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[out], &vec));
5450d0321e0SJeremy L Thompson       // Check if active output
5460d0321e0SJeremy L Thompson       if (vec == CEED_VECTOR_ACTIVE) {
547ca735530SJeremy L Thompson         CeedCallBackend(CeedVectorSetArray(impl->q_vecs_out[out], CEED_MEM_DEVICE, CEED_USE_POINTER, assembled_array));
548ca735530SJeremy L Thompson         CeedCallBackend(CeedQFunctionFieldGetSize(qf_output_fields[out], &size));
549ca735530SJeremy L Thompson         assembled_array += size * Q * num_elem;  // Advance the pointer by the size of the output
5500d0321e0SJeremy L Thompson       }
5510d0321e0SJeremy L Thompson     }
5520d0321e0SJeremy L Thompson     // Apply QFunction
553ca735530SJeremy L Thompson     CeedCallBackend(CeedQFunctionApply(qf, Q * num_elem, impl->q_vecs_in, impl->q_vecs_out));
5540d0321e0SJeremy L Thompson   }
5550d0321e0SJeremy L Thompson 
556ca735530SJeremy L Thompson   // Un-set output q_vecs to prevent accidental overwrite of Assembled
557ca735530SJeremy L Thompson   for (CeedInt out = 0; out < num_output_fields; out++) {
558ca735530SJeremy L Thompson     CeedVector vec;
559ca735530SJeremy L Thompson 
5600d0321e0SJeremy L Thompson     // Get output vector
561ca735530SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[out], &vec));
5620d0321e0SJeremy L Thompson     // Check if active output
5630d0321e0SJeremy L Thompson     if (vec == CEED_VECTOR_ACTIVE) {
564ca735530SJeremy L Thompson       CeedCallBackend(CeedVectorTakeArray(impl->q_vecs_out[out], CEED_MEM_DEVICE, NULL));
5650d0321e0SJeremy L Thompson     }
5660d0321e0SJeremy L Thompson   }
5670d0321e0SJeremy L Thompson 
5680d0321e0SJeremy L Thompson   // Restore input arrays
569ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorRestoreInputs_Cuda(num_input_fields, qf_input_fields, op_input_fields, true, e_data, impl));
5700d0321e0SJeremy L Thompson 
5710d0321e0SJeremy L Thompson   // Restore output
572ca735530SJeremy L Thompson   CeedCallBackend(CeedVectorRestoreArray(*assembled, &assembled_array));
5730d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
5740d0321e0SJeremy L Thompson }
5750d0321e0SJeremy L Thompson 
5760d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
5770d0321e0SJeremy L Thompson // Assemble Linear QFunction
5780d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
5792b730f8bSJeremy L Thompson static int CeedOperatorLinearAssembleQFunction_Cuda(CeedOperator op, CeedVector *assembled, CeedElemRestriction *rstr, CeedRequest *request) {
5802b730f8bSJeremy L Thompson   return CeedOperatorLinearAssembleQFunctionCore_Cuda(op, true, assembled, rstr, request);
5810d0321e0SJeremy L Thompson }
5820d0321e0SJeremy L Thompson 
5830d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
5840d0321e0SJeremy L Thompson // Update Assembled Linear QFunction
5850d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
5862b730f8bSJeremy L Thompson static int CeedOperatorLinearAssembleQFunctionUpdate_Cuda(CeedOperator op, CeedVector assembled, CeedElemRestriction rstr, CeedRequest *request) {
5872b730f8bSJeremy L Thompson   return CeedOperatorLinearAssembleQFunctionCore_Cuda(op, false, &assembled, &rstr, request);
5880d0321e0SJeremy L Thompson }
5890d0321e0SJeremy L Thompson 
5900d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
591004e4986SSebastian Grimberg // Assemble Diagonal Setup
5920d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
593cbfe683aSSebastian Grimberg static inline int CeedOperatorAssembleDiagonalSetup_Cuda(CeedOperator op) {
5940d0321e0SJeremy L Thompson   Ceed                ceed;
595004e4986SSebastian Grimberg   CeedInt             num_input_fields, num_output_fields, num_eval_modes_in = 0, num_eval_modes_out = 0;
596cbfe683aSSebastian Grimberg   CeedInt             q_comp, num_nodes, num_qpts;
597004e4986SSebastian Grimberg   CeedEvalMode       *eval_modes_in = NULL, *eval_modes_out = NULL;
598ca735530SJeremy L Thompson   CeedBasis           basis_in = NULL, basis_out = NULL;
599ca735530SJeremy L Thompson   CeedQFunctionField *qf_fields;
6000d0321e0SJeremy L Thompson   CeedQFunction       qf;
601ca735530SJeremy L Thompson   CeedOperatorField  *op_fields;
602ca735530SJeremy L Thompson   CeedOperator_Cuda  *impl;
603ca735530SJeremy L Thompson 
604ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
6052b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
606ca735530SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetNumArgs(qf, &num_input_fields, &num_output_fields));
6070d0321e0SJeremy L Thompson 
6080d0321e0SJeremy L Thompson   // Determine active input basis
609ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetFields(op, NULL, &op_fields, NULL, NULL));
610ca735530SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_fields, NULL, NULL));
611ca735530SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
6120d0321e0SJeremy L Thompson     CeedVector vec;
613ca735530SJeremy L Thompson 
614ca735530SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec));
6150d0321e0SJeremy L Thompson     if (vec == CEED_VECTOR_ACTIVE) {
616004e4986SSebastian Grimberg       CeedBasis    basis;
617004e4986SSebastian Grimberg       CeedEvalMode eval_mode;
618ca735530SJeremy L Thompson 
619004e4986SSebastian Grimberg       CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis));
620004e4986SSebastian Grimberg       CeedCheck(!basis_in || basis_in == basis, ceed, CEED_ERROR_BACKEND,
621004e4986SSebastian Grimberg                 "Backend does not implement operator diagonal assembly with multiple active bases");
622004e4986SSebastian Grimberg       basis_in = basis;
623004e4986SSebastian Grimberg       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode));
624004e4986SSebastian Grimberg       CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_in, eval_mode, &q_comp));
625004e4986SSebastian Grimberg       if (eval_mode != CEED_EVAL_WEIGHT) {
626004e4986SSebastian Grimberg         // q_comp = 1 if CEED_EVAL_NONE, CEED_EVAL_WEIGHT caught by QF assembly
627004e4986SSebastian Grimberg         CeedCallBackend(CeedRealloc(num_eval_modes_in + q_comp, &eval_modes_in));
628004e4986SSebastian Grimberg         for (CeedInt d = 0; d < q_comp; d++) eval_modes_in[num_eval_modes_in + d] = eval_mode;
629004e4986SSebastian Grimberg         num_eval_modes_in += q_comp;
6300d0321e0SJeremy L Thompson       }
6310d0321e0SJeremy L Thompson     }
6320d0321e0SJeremy L Thompson   }
6330d0321e0SJeremy L Thompson 
6340d0321e0SJeremy L Thompson   // Determine active output basis
635ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetFields(op, NULL, NULL, NULL, &op_fields));
636ca735530SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qf_fields));
637ca735530SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
6380d0321e0SJeremy L Thompson     CeedVector vec;
639ca735530SJeremy L Thompson 
640ca735530SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec));
6410d0321e0SJeremy L Thompson     if (vec == CEED_VECTOR_ACTIVE) {
642004e4986SSebastian Grimberg       CeedBasis    basis;
643004e4986SSebastian Grimberg       CeedEvalMode eval_mode;
644ca735530SJeremy L Thompson 
645004e4986SSebastian Grimberg       CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis));
646004e4986SSebastian Grimberg       CeedCheck(!basis_out || basis_out == basis, ceed, CEED_ERROR_BACKEND,
647004e4986SSebastian Grimberg                 "Backend does not implement operator diagonal assembly with multiple active bases");
648004e4986SSebastian Grimberg       basis_out = basis;
649004e4986SSebastian Grimberg       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode));
650004e4986SSebastian Grimberg       CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_out, eval_mode, &q_comp));
651004e4986SSebastian Grimberg       if (eval_mode != CEED_EVAL_WEIGHT) {
652004e4986SSebastian Grimberg         // q_comp = 1 if CEED_EVAL_NONE, CEED_EVAL_WEIGHT caught by QF assembly
653004e4986SSebastian Grimberg         CeedCallBackend(CeedRealloc(num_eval_modes_out + q_comp, &eval_modes_out));
654004e4986SSebastian Grimberg         for (CeedInt d = 0; d < q_comp; d++) eval_modes_out[num_eval_modes_out + d] = eval_mode;
655004e4986SSebastian Grimberg         num_eval_modes_out += q_comp;
6560d0321e0SJeremy L Thompson       }
6570d0321e0SJeremy L Thompson     }
6580d0321e0SJeremy L Thompson   }
6590d0321e0SJeremy L Thompson 
6600d0321e0SJeremy L Thompson   // Operator data struct
6612b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &impl));
6622b730f8bSJeremy L Thompson   CeedCallBackend(CeedCalloc(1, &impl->diag));
6630d0321e0SJeremy L Thompson   CeedOperatorDiag_Cuda *diag = impl->diag;
664ca735530SJeremy L Thompson 
665cbfe683aSSebastian Grimberg   // Basis matrices
666004e4986SSebastian Grimberg   CeedCallBackend(CeedBasisGetNumNodes(basis_in, &num_nodes));
667004e4986SSebastian Grimberg   if (basis_in == CEED_BASIS_NONE) num_qpts = num_nodes;
668004e4986SSebastian Grimberg   else CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis_in, &num_qpts));
669004e4986SSebastian Grimberg   const CeedInt interp_bytes     = num_nodes * num_qpts * sizeof(CeedScalar);
670004e4986SSebastian Grimberg   const CeedInt eval_modes_bytes = sizeof(CeedEvalMode);
671004e4986SSebastian Grimberg   bool          has_eval_none    = false;
6720d0321e0SJeremy L Thompson 
6730d0321e0SJeremy L Thompson   // CEED_EVAL_NONE
674004e4986SSebastian Grimberg   for (CeedInt i = 0; i < num_eval_modes_in; i++) has_eval_none = has_eval_none || (eval_modes_in[i] == CEED_EVAL_NONE);
675004e4986SSebastian Grimberg   for (CeedInt i = 0; i < num_eval_modes_out; i++) has_eval_none = has_eval_none || (eval_modes_out[i] == CEED_EVAL_NONE);
676004e4986SSebastian Grimberg   if (has_eval_none) {
6770d0321e0SJeremy L Thompson     CeedScalar *identity = NULL;
678ca735530SJeremy L Thompson 
679004e4986SSebastian Grimberg     CeedCallBackend(CeedCalloc(num_nodes * num_qpts, &identity));
680ca735530SJeremy L Thompson     for (CeedInt i = 0; i < (num_nodes < num_qpts ? num_nodes : num_qpts); i++) identity[i * num_nodes + i] = 1.0;
681ca735530SJeremy L Thompson     CeedCallCuda(ceed, cudaMalloc((void **)&diag->d_identity, interp_bytes));
682ca735530SJeremy L Thompson     CeedCallCuda(ceed, cudaMemcpy(diag->d_identity, identity, interp_bytes, cudaMemcpyHostToDevice));
683004e4986SSebastian Grimberg     CeedCallBackend(CeedFree(&identity));
6840d0321e0SJeremy L Thompson   }
6850d0321e0SJeremy L Thompson 
686004e4986SSebastian Grimberg   // CEED_EVAL_INTERP, CEED_EVAL_GRAD, CEED_EVAL_DIV, and CEED_EVAL_CURL
687004e4986SSebastian Grimberg   for (CeedInt in = 0; in < 2; in++) {
688004e4986SSebastian Grimberg     CeedFESpace fespace;
689004e4986SSebastian Grimberg     CeedBasis   basis = in ? basis_in : basis_out;
6900d0321e0SJeremy L Thompson 
691004e4986SSebastian Grimberg     CeedCallBackend(CeedBasisGetFESpace(basis, &fespace));
692004e4986SSebastian Grimberg     switch (fespace) {
693004e4986SSebastian Grimberg       case CEED_FE_SPACE_H1: {
694004e4986SSebastian Grimberg         CeedInt           q_comp_interp, q_comp_grad;
695004e4986SSebastian Grimberg         const CeedScalar *interp, *grad;
696004e4986SSebastian Grimberg         CeedScalar       *d_interp, *d_grad;
6970d0321e0SJeremy L Thompson 
698004e4986SSebastian Grimberg         CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp));
699004e4986SSebastian Grimberg         CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_GRAD, &q_comp_grad));
7000d0321e0SJeremy L Thompson 
701004e4986SSebastian Grimberg         CeedCallBackend(CeedBasisGetInterp(basis, &interp));
702004e4986SSebastian Grimberg         CeedCallCuda(ceed, cudaMalloc((void **)&d_interp, interp_bytes * q_comp_interp));
703004e4986SSebastian Grimberg         CeedCallCuda(ceed, cudaMemcpy(d_interp, interp, interp_bytes * q_comp_interp, cudaMemcpyHostToDevice));
704004e4986SSebastian Grimberg         CeedCallBackend(CeedBasisGetGrad(basis, &grad));
705004e4986SSebastian Grimberg         CeedCallCuda(ceed, cudaMalloc((void **)&d_grad, interp_bytes * q_comp_grad));
706004e4986SSebastian Grimberg         CeedCallCuda(ceed, cudaMemcpy(d_grad, grad, interp_bytes * q_comp_grad, cudaMemcpyHostToDevice));
707004e4986SSebastian Grimberg         if (in) {
708004e4986SSebastian Grimberg           diag->d_interp_in = d_interp;
709004e4986SSebastian Grimberg           diag->d_grad_in   = d_grad;
710004e4986SSebastian Grimberg         } else {
711004e4986SSebastian Grimberg           diag->d_interp_out = d_interp;
712004e4986SSebastian Grimberg           diag->d_grad_out   = d_grad;
713004e4986SSebastian Grimberg         }
714004e4986SSebastian Grimberg       } break;
715004e4986SSebastian Grimberg       case CEED_FE_SPACE_HDIV: {
716004e4986SSebastian Grimberg         CeedInt           q_comp_interp, q_comp_div;
717004e4986SSebastian Grimberg         const CeedScalar *interp, *div;
718004e4986SSebastian Grimberg         CeedScalar       *d_interp, *d_div;
719004e4986SSebastian Grimberg 
720004e4986SSebastian Grimberg         CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp));
721004e4986SSebastian Grimberg         CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_DIV, &q_comp_div));
722004e4986SSebastian Grimberg 
723004e4986SSebastian Grimberg         CeedCallBackend(CeedBasisGetInterp(basis, &interp));
724004e4986SSebastian Grimberg         CeedCallCuda(ceed, cudaMalloc((void **)&d_interp, interp_bytes * q_comp_interp));
725004e4986SSebastian Grimberg         CeedCallCuda(ceed, cudaMemcpy(d_interp, interp, interp_bytes * q_comp_interp, cudaMemcpyHostToDevice));
726004e4986SSebastian Grimberg         CeedCallBackend(CeedBasisGetDiv(basis, &div));
727004e4986SSebastian Grimberg         CeedCallCuda(ceed, cudaMalloc((void **)&d_div, interp_bytes * q_comp_div));
728004e4986SSebastian Grimberg         CeedCallCuda(ceed, cudaMemcpy(d_div, div, interp_bytes * q_comp_div, cudaMemcpyHostToDevice));
729004e4986SSebastian Grimberg         if (in) {
730004e4986SSebastian Grimberg           diag->d_interp_in = d_interp;
731004e4986SSebastian Grimberg           diag->d_div_in    = d_div;
732004e4986SSebastian Grimberg         } else {
733004e4986SSebastian Grimberg           diag->d_interp_out = d_interp;
734004e4986SSebastian Grimberg           diag->d_div_out    = d_div;
735004e4986SSebastian Grimberg         }
736004e4986SSebastian Grimberg       } break;
737004e4986SSebastian Grimberg       case CEED_FE_SPACE_HCURL: {
738004e4986SSebastian Grimberg         CeedInt           q_comp_interp, q_comp_curl;
739004e4986SSebastian Grimberg         const CeedScalar *interp, *curl;
740004e4986SSebastian Grimberg         CeedScalar       *d_interp, *d_curl;
741004e4986SSebastian Grimberg 
742004e4986SSebastian Grimberg         CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp));
743004e4986SSebastian Grimberg         CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_CURL, &q_comp_curl));
744004e4986SSebastian Grimberg 
745004e4986SSebastian Grimberg         CeedCallBackend(CeedBasisGetInterp(basis, &interp));
746004e4986SSebastian Grimberg         CeedCallCuda(ceed, cudaMalloc((void **)&d_interp, interp_bytes * q_comp_interp));
747004e4986SSebastian Grimberg         CeedCallCuda(ceed, cudaMemcpy(d_interp, interp, interp_bytes * q_comp_interp, cudaMemcpyHostToDevice));
748004e4986SSebastian Grimberg         CeedCallBackend(CeedBasisGetCurl(basis, &curl));
749004e4986SSebastian Grimberg         CeedCallCuda(ceed, cudaMalloc((void **)&d_curl, interp_bytes * q_comp_curl));
750004e4986SSebastian Grimberg         CeedCallCuda(ceed, cudaMemcpy(d_curl, curl, interp_bytes * q_comp_curl, cudaMemcpyHostToDevice));
751004e4986SSebastian Grimberg         if (in) {
752004e4986SSebastian Grimberg           diag->d_interp_in = d_interp;
753004e4986SSebastian Grimberg           diag->d_curl_in   = d_curl;
754004e4986SSebastian Grimberg         } else {
755004e4986SSebastian Grimberg           diag->d_interp_out = d_interp;
756004e4986SSebastian Grimberg           diag->d_curl_out   = d_curl;
757004e4986SSebastian Grimberg         }
758004e4986SSebastian Grimberg       } break;
759004e4986SSebastian Grimberg     }
760004e4986SSebastian Grimberg   }
761004e4986SSebastian Grimberg 
762004e4986SSebastian Grimberg   // Arrays of eval_modes
763004e4986SSebastian Grimberg   CeedCallCuda(ceed, cudaMalloc((void **)&diag->d_eval_modes_in, num_eval_modes_in * eval_modes_bytes));
764004e4986SSebastian Grimberg   CeedCallCuda(ceed, cudaMemcpy(diag->d_eval_modes_in, eval_modes_in, num_eval_modes_in * eval_modes_bytes, cudaMemcpyHostToDevice));
765004e4986SSebastian Grimberg   CeedCallCuda(ceed, cudaMalloc((void **)&diag->d_eval_modes_out, num_eval_modes_out * eval_modes_bytes));
766004e4986SSebastian Grimberg   CeedCallCuda(ceed, cudaMemcpy(diag->d_eval_modes_out, eval_modes_out, num_eval_modes_out * eval_modes_bytes, cudaMemcpyHostToDevice));
767004e4986SSebastian Grimberg   CeedCallBackend(CeedFree(&eval_modes_in));
768004e4986SSebastian Grimberg   CeedCallBackend(CeedFree(&eval_modes_out));
7690d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
7700d0321e0SJeremy L Thompson }
7710d0321e0SJeremy L Thompson 
7720d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
773cbfe683aSSebastian Grimberg // Assemble Diagonal Setup (Compilation)
774cbfe683aSSebastian Grimberg //------------------------------------------------------------------------------
775cbfe683aSSebastian Grimberg static inline int CeedOperatorAssembleDiagonalSetupCompile_Cuda(CeedOperator op, CeedInt use_ceedsize_idx, const bool is_point_block) {
776cbfe683aSSebastian Grimberg   Ceed                ceed;
77722070f95SJeremy L Thompson   char               *diagonal_kernel_source;
77822070f95SJeremy L Thompson   const char         *diagonal_kernel_path;
779cbfe683aSSebastian Grimberg   CeedInt             num_input_fields, num_output_fields, num_eval_modes_in = 0, num_eval_modes_out = 0;
780cbfe683aSSebastian Grimberg   CeedInt             num_comp, q_comp, num_nodes, num_qpts;
781cbfe683aSSebastian Grimberg   CeedBasis           basis_in = NULL, basis_out = NULL;
782cbfe683aSSebastian Grimberg   CeedQFunctionField *qf_fields;
783cbfe683aSSebastian Grimberg   CeedQFunction       qf;
784cbfe683aSSebastian Grimberg   CeedOperatorField  *op_fields;
785cbfe683aSSebastian Grimberg   CeedOperator_Cuda  *impl;
786cbfe683aSSebastian Grimberg 
787cbfe683aSSebastian Grimberg   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
788cbfe683aSSebastian Grimberg   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
789cbfe683aSSebastian Grimberg   CeedCallBackend(CeedQFunctionGetNumArgs(qf, &num_input_fields, &num_output_fields));
790cbfe683aSSebastian Grimberg 
791cbfe683aSSebastian Grimberg   // Determine active input basis
792cbfe683aSSebastian Grimberg   CeedCallBackend(CeedOperatorGetFields(op, NULL, &op_fields, NULL, NULL));
793cbfe683aSSebastian Grimberg   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_fields, NULL, NULL));
794cbfe683aSSebastian Grimberg   for (CeedInt i = 0; i < num_input_fields; i++) {
795cbfe683aSSebastian Grimberg     CeedVector vec;
796cbfe683aSSebastian Grimberg 
797cbfe683aSSebastian Grimberg     CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec));
798cbfe683aSSebastian Grimberg     if (vec == CEED_VECTOR_ACTIVE) {
799cbfe683aSSebastian Grimberg       CeedEvalMode eval_mode;
800cbfe683aSSebastian Grimberg 
801cbfe683aSSebastian Grimberg       CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis_in));
802cbfe683aSSebastian Grimberg       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode));
803cbfe683aSSebastian Grimberg       CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_in, eval_mode, &q_comp));
804cbfe683aSSebastian Grimberg       if (eval_mode != CEED_EVAL_WEIGHT) {
805cbfe683aSSebastian Grimberg         num_eval_modes_in += q_comp;
806cbfe683aSSebastian Grimberg       }
807cbfe683aSSebastian Grimberg     }
808cbfe683aSSebastian Grimberg   }
809cbfe683aSSebastian Grimberg 
810cbfe683aSSebastian Grimberg   // Determine active output basis
811cbfe683aSSebastian Grimberg   CeedCallBackend(CeedOperatorGetFields(op, NULL, NULL, NULL, &op_fields));
812cbfe683aSSebastian Grimberg   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qf_fields));
813cbfe683aSSebastian Grimberg   for (CeedInt i = 0; i < num_output_fields; i++) {
814cbfe683aSSebastian Grimberg     CeedVector vec;
815cbfe683aSSebastian Grimberg 
816cbfe683aSSebastian Grimberg     CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec));
817cbfe683aSSebastian Grimberg     if (vec == CEED_VECTOR_ACTIVE) {
818cbfe683aSSebastian Grimberg       CeedEvalMode eval_mode;
819cbfe683aSSebastian Grimberg 
820cbfe683aSSebastian Grimberg       CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis_out));
821cbfe683aSSebastian Grimberg       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode));
822cbfe683aSSebastian Grimberg       CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_out, eval_mode, &q_comp));
823cbfe683aSSebastian Grimberg       if (eval_mode != CEED_EVAL_WEIGHT) {
824cbfe683aSSebastian Grimberg         num_eval_modes_out += q_comp;
825cbfe683aSSebastian Grimberg       }
826cbfe683aSSebastian Grimberg     }
827cbfe683aSSebastian Grimberg   }
828cbfe683aSSebastian Grimberg 
829cbfe683aSSebastian Grimberg   // Operator data struct
830cbfe683aSSebastian Grimberg   CeedCallBackend(CeedOperatorGetData(op, &impl));
831cbfe683aSSebastian Grimberg   CeedOperatorDiag_Cuda *diag = impl->diag;
832cbfe683aSSebastian Grimberg 
833cbfe683aSSebastian Grimberg   // Assemble kernel
834cbfe683aSSebastian Grimberg   CUmodule *module          = is_point_block ? &diag->module_point_block : &diag->module;
835cbfe683aSSebastian Grimberg   CeedInt   elems_per_block = 1;
836cbfe683aSSebastian Grimberg   CeedCallBackend(CeedBasisGetNumNodes(basis_in, &num_nodes));
837cbfe683aSSebastian Grimberg   CeedCallBackend(CeedBasisGetNumComponents(basis_in, &num_comp));
838cbfe683aSSebastian Grimberg   if (basis_in == CEED_BASIS_NONE) num_qpts = num_nodes;
839cbfe683aSSebastian Grimberg   else CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis_in, &num_qpts));
840cbfe683aSSebastian Grimberg   CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-operator-assemble-diagonal.h", &diagonal_kernel_path));
841cbfe683aSSebastian Grimberg   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Diagonal Assembly Kernel Source -----\n");
842cbfe683aSSebastian Grimberg   CeedCallBackend(CeedLoadSourceToBuffer(ceed, diagonal_kernel_path, &diagonal_kernel_source));
843cbfe683aSSebastian Grimberg   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Diagonal Assembly Source Complete! -----\n");
844cbfe683aSSebastian Grimberg   CeedCallCuda(ceed, CeedCompile_Cuda(ceed, diagonal_kernel_source, module, 8, "NUM_EVAL_MODES_IN", num_eval_modes_in, "NUM_EVAL_MODES_OUT",
845cbfe683aSSebastian Grimberg                                       num_eval_modes_out, "NUM_COMP", num_comp, "NUM_NODES", num_nodes, "NUM_QPTS", num_qpts, "USE_CEEDSIZE",
846cbfe683aSSebastian Grimberg                                       use_ceedsize_idx, "USE_POINT_BLOCK", is_point_block ? 1 : 0, "BLOCK_SIZE", num_nodes * elems_per_block));
847cbfe683aSSebastian Grimberg   CeedCallCuda(ceed, CeedGetKernel_Cuda(ceed, *module, "LinearDiagonal", is_point_block ? &diag->LinearPointBlock : &diag->LinearDiagonal));
848cbfe683aSSebastian Grimberg   CeedCallBackend(CeedFree(&diagonal_kernel_path));
849cbfe683aSSebastian Grimberg   CeedCallBackend(CeedFree(&diagonal_kernel_source));
850cbfe683aSSebastian Grimberg   return CEED_ERROR_SUCCESS;
851cbfe683aSSebastian Grimberg }
852cbfe683aSSebastian Grimberg 
853cbfe683aSSebastian Grimberg //------------------------------------------------------------------------------
854004e4986SSebastian Grimberg // Assemble Diagonal Core
8550d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
856ca735530SJeremy L Thompson static inline int CeedOperatorAssembleDiagonalCore_Cuda(CeedOperator op, CeedVector assembled, CeedRequest *request, const bool is_point_block) {
8570d0321e0SJeremy L Thompson   Ceed                ceed;
858cbfe683aSSebastian Grimberg   CeedInt             num_elem, num_nodes;
859ca735530SJeremy L Thompson   CeedScalar         *elem_diag_array;
860ca735530SJeremy L Thompson   const CeedScalar   *assembled_qf_array;
861004e4986SSebastian Grimberg   CeedVector          assembled_qf   = NULL, elem_diag;
862004e4986SSebastian Grimberg   CeedElemRestriction assembled_rstr = NULL, rstr_in, rstr_out, diag_rstr;
8630d0321e0SJeremy L Thompson   CeedOperator_Cuda  *impl;
864ca735530SJeremy L Thompson 
865ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
8662b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &impl));
8670d0321e0SJeremy L Thompson 
8680d0321e0SJeremy L Thompson   // Assemble QFunction
869004e4986SSebastian Grimberg   CeedCallBackend(CeedOperatorLinearAssembleQFunctionBuildOrUpdate(op, &assembled_qf, &assembled_rstr, request));
870004e4986SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionDestroy(&assembled_rstr));
871004e4986SSebastian Grimberg   CeedCallBackend(CeedVectorGetArrayRead(assembled_qf, CEED_MEM_DEVICE, &assembled_qf_array));
8720d0321e0SJeremy L Thompson 
873cbfe683aSSebastian Grimberg   // Setup
874cbfe683aSSebastian Grimberg   if (!impl->diag) CeedCallBackend(CeedOperatorAssembleDiagonalSetup_Cuda(op));
875cbfe683aSSebastian Grimberg   CeedOperatorDiag_Cuda *diag = impl->diag;
876cbfe683aSSebastian Grimberg 
877cbfe683aSSebastian Grimberg   assert(diag != NULL);
878cbfe683aSSebastian Grimberg 
879cbfe683aSSebastian Grimberg   // Assemble kernel if needed
880cbfe683aSSebastian Grimberg   if ((!is_point_block && !diag->LinearDiagonal) || (is_point_block && !diag->LinearPointBlock)) {
881cbfe683aSSebastian Grimberg     CeedSize assembled_length, assembled_qf_length;
882cbfe683aSSebastian Grimberg     CeedInt  use_ceedsize_idx = 0;
883f7c1b517Snbeams     CeedCallBackend(CeedVectorGetLength(assembled, &assembled_length));
884ca735530SJeremy L Thompson     CeedCallBackend(CeedVectorGetLength(assembled_qf, &assembled_qf_length));
885ca735530SJeremy L Thompson     if ((assembled_length > INT_MAX) || (assembled_qf_length > INT_MAX)) use_ceedsize_idx = 1;
886f7c1b517Snbeams 
887cbfe683aSSebastian Grimberg     CeedCallBackend(CeedOperatorAssembleDiagonalSetupCompile_Cuda(op, use_ceedsize_idx, is_point_block));
888cbfe683aSSebastian Grimberg   }
8890d0321e0SJeremy L Thompson 
890004e4986SSebastian Grimberg   // Restriction and diagonal vector
891004e4986SSebastian Grimberg   CeedCallBackend(CeedOperatorGetActiveElemRestrictions(op, &rstr_in, &rstr_out));
892004e4986SSebastian Grimberg   CeedCheck(rstr_in == rstr_out, ceed, CEED_ERROR_BACKEND,
893004e4986SSebastian Grimberg             "Cannot assemble operator diagonal with different input and output active element restrictions");
894004e4986SSebastian Grimberg   if (!is_point_block && !diag->diag_rstr) {
895004e4986SSebastian Grimberg     CeedCallBackend(CeedElemRestrictionCreateUnsignedCopy(rstr_out, &diag->diag_rstr));
896004e4986SSebastian Grimberg     CeedCallBackend(CeedElemRestrictionCreateVector(diag->diag_rstr, NULL, &diag->elem_diag));
897004e4986SSebastian Grimberg   } else if (is_point_block && !diag->point_block_diag_rstr) {
898004e4986SSebastian Grimberg     CeedCallBackend(CeedOperatorCreateActivePointBlockRestriction(rstr_out, &diag->point_block_diag_rstr));
899004e4986SSebastian Grimberg     CeedCallBackend(CeedElemRestrictionCreateVector(diag->point_block_diag_rstr, NULL, &diag->point_block_elem_diag));
9000d0321e0SJeremy L Thompson   }
901004e4986SSebastian Grimberg   diag_rstr = is_point_block ? diag->point_block_diag_rstr : diag->diag_rstr;
902004e4986SSebastian Grimberg   elem_diag = is_point_block ? diag->point_block_elem_diag : diag->elem_diag;
903ca735530SJeremy L Thompson   CeedCallBackend(CeedVectorSetValue(elem_diag, 0.0));
9040d0321e0SJeremy L Thompson 
90591db28b6SZach Atkins   // Only assemble diagonal if the basis has nodes, otherwise inputs are null pointers
906004e4986SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionGetElementSize(diag_rstr, &num_nodes));
907004e4986SSebastian Grimberg   if (num_nodes > 0) {
9080d0321e0SJeremy L Thompson     // Assemble element operator diagonals
909ca735530SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumElements(diag_rstr, &num_elem));
910004e4986SSebastian Grimberg     CeedCallBackend(CeedVectorGetArray(elem_diag, CEED_MEM_DEVICE, &elem_diag_array));
9110d0321e0SJeremy L Thompson 
9120d0321e0SJeremy L Thompson     // Compute the diagonal of B^T D B
913004e4986SSebastian Grimberg     CeedInt elems_per_block = 1;
914004e4986SSebastian Grimberg     CeedInt grid            = CeedDivUpInt(num_elem, elems_per_block);
915004e4986SSebastian Grimberg     void   *args[]          = {(void *)&num_elem,      &diag->d_identity,       &diag->d_interp_in,  &diag->d_grad_in, &diag->d_div_in,
916004e4986SSebastian Grimberg                                &diag->d_curl_in,       &diag->d_interp_out,     &diag->d_grad_out,   &diag->d_div_out, &diag->d_curl_out,
917004e4986SSebastian Grimberg                                &diag->d_eval_modes_in, &diag->d_eval_modes_out, &assembled_qf_array, &elem_diag_array};
918004e4986SSebastian Grimberg 
919ca735530SJeremy L Thompson     if (is_point_block) {
920004e4986SSebastian Grimberg       CeedCallBackend(CeedRunKernelDim_Cuda(ceed, diag->LinearPointBlock, grid, num_nodes, 1, elems_per_block, args));
9210d0321e0SJeremy L Thompson     } else {
922004e4986SSebastian Grimberg       CeedCallBackend(CeedRunKernelDim_Cuda(ceed, diag->LinearDiagonal, grid, num_nodes, 1, elems_per_block, args));
9230d0321e0SJeremy L Thompson     }
9240d0321e0SJeremy L Thompson 
9250d0321e0SJeremy L Thompson     // Restore arrays
926ca735530SJeremy L Thompson     CeedCallBackend(CeedVectorRestoreArray(elem_diag, &elem_diag_array));
927ca735530SJeremy L Thompson     CeedCallBackend(CeedVectorRestoreArrayRead(assembled_qf, &assembled_qf_array));
92891db28b6SZach Atkins   }
9290d0321e0SJeremy L Thompson 
9300d0321e0SJeremy L Thompson   // Assemble local operator diagonal
931ca735530SJeremy L Thompson   CeedCallBackend(CeedElemRestrictionApply(diag_rstr, CEED_TRANSPOSE, elem_diag, assembled, request));
9320d0321e0SJeremy L Thompson 
9330d0321e0SJeremy L Thompson   // Cleanup
934ca735530SJeremy L Thompson   CeedCallBackend(CeedVectorDestroy(&assembled_qf));
9350d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
9360d0321e0SJeremy L Thompson }
9370d0321e0SJeremy L Thompson 
9380d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
9390d0321e0SJeremy L Thompson // Assemble Linear Diagonal
9400d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
9412b730f8bSJeremy L Thompson static int CeedOperatorLinearAssembleAddDiagonal_Cuda(CeedOperator op, CeedVector assembled, CeedRequest *request) {
9422b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorAssembleDiagonalCore_Cuda(op, assembled, request, false));
9436aa95790SJeremy L Thompson   return CEED_ERROR_SUCCESS;
9440d0321e0SJeremy L Thompson }
9450d0321e0SJeremy L Thompson 
9460d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
9470d0321e0SJeremy L Thompson // Assemble Linear Point Block Diagonal
9480d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
9492b730f8bSJeremy L Thompson static int CeedOperatorLinearAssembleAddPointBlockDiagonal_Cuda(CeedOperator op, CeedVector assembled, CeedRequest *request) {
9502b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorAssembleDiagonalCore_Cuda(op, assembled, request, true));
9516aa95790SJeremy L Thompson   return CEED_ERROR_SUCCESS;
9520d0321e0SJeremy L Thompson }
9530d0321e0SJeremy L Thompson 
9540d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
955004e4986SSebastian Grimberg // Single Operator Assembly Setup
956cc132f9aSnbeams //------------------------------------------------------------------------------
957f7c1b517Snbeams static int CeedSingleOperatorAssembleSetup_Cuda(CeedOperator op, CeedInt use_ceedsize_idx) {
958cc132f9aSnbeams   Ceed                ceed;
959004e4986SSebastian Grimberg   Ceed_Cuda          *cuda_data;
96022070f95SJeremy L Thompson   char               *assembly_kernel_source;
96122070f95SJeremy L Thompson   const char         *assembly_kernel_path;
962004e4986SSebastian Grimberg   CeedInt             num_input_fields, num_output_fields, num_eval_modes_in = 0, num_eval_modes_out = 0;
9633b38d1dfSJeremy L Thompson   CeedInt             elem_size_in, num_qpts_in = 0, num_comp_in, elem_size_out, num_qpts_out, num_comp_out, q_comp;
964004e4986SSebastian Grimberg   CeedEvalMode       *eval_modes_in = NULL, *eval_modes_out = NULL;
965ca735530SJeremy L Thompson   CeedElemRestriction rstr_in = NULL, rstr_out = NULL;
966ca735530SJeremy L Thompson   CeedBasis           basis_in = NULL, basis_out = NULL;
967ca735530SJeremy L Thompson   CeedQFunctionField *qf_fields;
968ca735530SJeremy L Thompson   CeedQFunction       qf;
969ca735530SJeremy L Thompson   CeedOperatorField  *input_fields, *output_fields;
970cc132f9aSnbeams   CeedOperator_Cuda  *impl;
971ca735530SJeremy L Thompson 
972ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
9732b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &impl));
974cc132f9aSnbeams 
975cc132f9aSnbeams   // Get intput and output fields
9762b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &input_fields, &num_output_fields, &output_fields));
977cc132f9aSnbeams 
978cc132f9aSnbeams   // Determine active input basis eval mode
9792b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
9802b730f8bSJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_fields, NULL, NULL));
981cc132f9aSnbeams   for (CeedInt i = 0; i < num_input_fields; i++) {
982cc132f9aSnbeams     CeedVector vec;
983ca735530SJeremy L Thompson 
9842b730f8bSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetVector(input_fields[i], &vec));
985cc132f9aSnbeams     if (vec == CEED_VECTOR_ACTIVE) {
986004e4986SSebastian Grimberg       CeedBasis    basis;
987ca735530SJeremy L Thompson       CeedEvalMode eval_mode;
988ca735530SJeremy L Thompson 
989004e4986SSebastian Grimberg       CeedCallBackend(CeedOperatorFieldGetBasis(input_fields[i], &basis));
990004e4986SSebastian Grimberg       CeedCheck(!basis_in || basis_in == basis, ceed, CEED_ERROR_BACKEND, "Backend does not implement operator assembly with multiple active bases");
991004e4986SSebastian Grimberg       basis_in = basis;
9922b730f8bSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetElemRestriction(input_fields[i], &rstr_in));
993004e4986SSebastian Grimberg       CeedCallBackend(CeedElemRestrictionGetElementSize(rstr_in, &elem_size_in));
994004e4986SSebastian Grimberg       if (basis_in == CEED_BASIS_NONE) num_qpts_in = elem_size_in;
995004e4986SSebastian Grimberg       else CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis_in, &num_qpts_in));
9962b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode));
997004e4986SSebastian Grimberg       CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_in, eval_mode, &q_comp));
998004e4986SSebastian Grimberg       if (eval_mode != CEED_EVAL_WEIGHT) {
999004e4986SSebastian Grimberg         // q_comp = 1 if CEED_EVAL_NONE, CEED_EVAL_WEIGHT caught by QF Assembly
1000004e4986SSebastian Grimberg         CeedCallBackend(CeedRealloc(num_eval_modes_in + q_comp, &eval_modes_in));
1001004e4986SSebastian Grimberg         for (CeedInt d = 0; d < q_comp; d++) {
1002004e4986SSebastian Grimberg           eval_modes_in[num_eval_modes_in + d] = eval_mode;
1003cc132f9aSnbeams         }
1004004e4986SSebastian Grimberg         num_eval_modes_in += q_comp;
1005cc132f9aSnbeams       }
1006cc132f9aSnbeams     }
1007cc132f9aSnbeams   }
1008cc132f9aSnbeams 
1009cc132f9aSnbeams   // Determine active output basis; basis_out and rstr_out only used if same as input, TODO
10102b730f8bSJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qf_fields));
1011cc132f9aSnbeams   for (CeedInt i = 0; i < num_output_fields; i++) {
1012cc132f9aSnbeams     CeedVector vec;
1013ca735530SJeremy L Thompson 
10142b730f8bSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetVector(output_fields[i], &vec));
1015cc132f9aSnbeams     if (vec == CEED_VECTOR_ACTIVE) {
1016004e4986SSebastian Grimberg       CeedBasis    basis;
1017ca735530SJeremy L Thompson       CeedEvalMode eval_mode;
1018ca735530SJeremy L Thompson 
1019004e4986SSebastian Grimberg       CeedCallBackend(CeedOperatorFieldGetBasis(output_fields[i], &basis));
1020004e4986SSebastian Grimberg       CeedCheck(!basis_out || basis_out == basis, ceed, CEED_ERROR_BACKEND,
1021004e4986SSebastian Grimberg                 "Backend does not implement operator assembly with multiple active bases");
1022004e4986SSebastian Grimberg       basis_out = basis;
10232b730f8bSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetElemRestriction(output_fields[i], &rstr_out));
1024004e4986SSebastian Grimberg       CeedCallBackend(CeedElemRestrictionGetElementSize(rstr_out, &elem_size_out));
1025004e4986SSebastian Grimberg       if (basis_out == CEED_BASIS_NONE) num_qpts_out = elem_size_out;
1026004e4986SSebastian Grimberg       else CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis_out, &num_qpts_out));
1027004e4986SSebastian Grimberg       CeedCheck(num_qpts_in == num_qpts_out, ceed, CEED_ERROR_UNSUPPORTED,
1028004e4986SSebastian Grimberg                 "Active input and output bases must have the same number of quadrature points");
10292b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode));
1030004e4986SSebastian Grimberg       CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_out, eval_mode, &q_comp));
1031004e4986SSebastian Grimberg       if (eval_mode != CEED_EVAL_WEIGHT) {
1032004e4986SSebastian Grimberg         // q_comp = 1 if CEED_EVAL_NONE, CEED_EVAL_WEIGHT caught by QF Assembly
1033004e4986SSebastian Grimberg         CeedCallBackend(CeedRealloc(num_eval_modes_out + q_comp, &eval_modes_out));
1034004e4986SSebastian Grimberg         for (CeedInt d = 0; d < q_comp; d++) {
1035004e4986SSebastian Grimberg           eval_modes_out[num_eval_modes_out + d] = eval_mode;
1036004e4986SSebastian Grimberg         }
1037004e4986SSebastian Grimberg         num_eval_modes_out += q_comp;
1038cc132f9aSnbeams       }
1039cc132f9aSnbeams     }
1040cc132f9aSnbeams   }
1041004e4986SSebastian Grimberg   CeedCheck(num_eval_modes_in > 0 && num_eval_modes_out > 0, ceed, CEED_ERROR_UNSUPPORTED, "Cannot assemble operator without inputs/outputs");
1042cc132f9aSnbeams 
10432b730f8bSJeremy L Thompson   CeedCallBackend(CeedCalloc(1, &impl->asmb));
1044cc132f9aSnbeams   CeedOperatorAssemble_Cuda *asmb = impl->asmb;
1045004e4986SSebastian Grimberg   asmb->elems_per_block           = 1;
1046004e4986SSebastian Grimberg   asmb->block_size_x              = elem_size_in;
1047004e4986SSebastian Grimberg   asmb->block_size_y              = elem_size_out;
1048ca735530SJeremy L Thompson 
10492b730f8bSJeremy L Thompson   CeedCallBackend(CeedGetData(ceed, &cuda_data));
1050004e4986SSebastian Grimberg   bool fallback = asmb->block_size_x * asmb->block_size_y * asmb->elems_per_block > cuda_data->device_prop.maxThreadsPerBlock;
1051004e4986SSebastian Grimberg 
1052004e4986SSebastian Grimberg   if (fallback) {
1053004e4986SSebastian Grimberg     // Use fallback kernel with 1D threadblock
1054004e4986SSebastian Grimberg     asmb->block_size_y = 1;
1055004e4986SSebastian Grimberg   }
1056004e4986SSebastian Grimberg 
1057004e4986SSebastian Grimberg   // Compile kernels
1058004e4986SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionGetNumComponents(rstr_in, &num_comp_in));
1059004e4986SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionGetNumComponents(rstr_out, &num_comp_out));
10602b730f8bSJeremy L Thompson   CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-operator-assemble.h", &assembly_kernel_path));
106123d4529eSJeremy L Thompson   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Assembly Kernel Source -----\n");
10622b730f8bSJeremy L Thompson   CeedCallBackend(CeedLoadSourceToBuffer(ceed, assembly_kernel_path, &assembly_kernel_source));
106323d4529eSJeremy L Thompson   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Assembly Source Complete! -----\n");
1064004e4986SSebastian Grimberg   CeedCallBackend(CeedCompile_Cuda(ceed, assembly_kernel_source, &asmb->module, 10, "NUM_EVAL_MODES_IN", num_eval_modes_in, "NUM_EVAL_MODES_OUT",
1065004e4986SSebastian Grimberg                                    num_eval_modes_out, "NUM_COMP_IN", num_comp_in, "NUM_COMP_OUT", num_comp_out, "NUM_NODES_IN", elem_size_in,
1066004e4986SSebastian Grimberg                                    "NUM_NODES_OUT", elem_size_out, "NUM_QPTS", num_qpts_in, "BLOCK_SIZE",
1067cbfe683aSSebastian Grimberg                                    asmb->block_size_x * asmb->block_size_y * asmb->elems_per_block, "BLOCK_SIZE_Y", asmb->block_size_y,
1068cbfe683aSSebastian Grimberg                                    "USE_CEEDSIZE", use_ceedsize_idx));
1069004e4986SSebastian Grimberg   CeedCallBackend(CeedGetKernel_Cuda(ceed, asmb->module, "LinearAssemble", &asmb->LinearAssemble));
10702b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&assembly_kernel_path));
10712b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&assembly_kernel_source));
1072cc132f9aSnbeams 
1073004e4986SSebastian Grimberg   // Load into B_in, in order that they will be used in eval_modes_in
1074004e4986SSebastian Grimberg   {
1075004e4986SSebastian Grimberg     const CeedInt in_bytes           = elem_size_in * num_qpts_in * num_eval_modes_in * sizeof(CeedScalar);
1076004e4986SSebastian Grimberg     CeedInt       d_in               = 0;
1077004e4986SSebastian Grimberg     CeedEvalMode  eval_modes_in_prev = CEED_EVAL_NONE;
1078004e4986SSebastian Grimberg     bool          has_eval_none      = false;
1079004e4986SSebastian Grimberg     CeedScalar   *identity           = NULL;
1080ca735530SJeremy L Thompson 
1081004e4986SSebastian Grimberg     for (CeedInt i = 0; i < num_eval_modes_in; i++) {
1082004e4986SSebastian Grimberg       has_eval_none = has_eval_none || (eval_modes_in[i] == CEED_EVAL_NONE);
1083004e4986SSebastian Grimberg     }
1084004e4986SSebastian Grimberg     if (has_eval_none) {
1085004e4986SSebastian Grimberg       CeedCallBackend(CeedCalloc(elem_size_in * num_qpts_in, &identity));
1086004e4986SSebastian Grimberg       for (CeedInt i = 0; i < (elem_size_in < num_qpts_in ? elem_size_in : num_qpts_in); i++) identity[i * elem_size_in + i] = 1.0;
1087004e4986SSebastian Grimberg     }
1088cc132f9aSnbeams 
1089004e4986SSebastian Grimberg     CeedCallCuda(ceed, cudaMalloc((void **)&asmb->d_B_in, in_bytes));
1090004e4986SSebastian Grimberg     for (CeedInt i = 0; i < num_eval_modes_in; i++) {
1091004e4986SSebastian Grimberg       const CeedScalar *h_B_in;
1092ca735530SJeremy L Thompson 
1093004e4986SSebastian Grimberg       CeedCallBackend(CeedOperatorGetBasisPointer(basis_in, eval_modes_in[i], identity, &h_B_in));
1094004e4986SSebastian Grimberg       CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_in, eval_modes_in[i], &q_comp));
1095004e4986SSebastian Grimberg       if (q_comp > 1) {
1096004e4986SSebastian Grimberg         if (i == 0 || eval_modes_in[i] != eval_modes_in_prev) d_in = 0;
1097004e4986SSebastian Grimberg         else h_B_in = &h_B_in[(++d_in) * elem_size_in * num_qpts_in];
1098004e4986SSebastian Grimberg       }
1099004e4986SSebastian Grimberg       eval_modes_in_prev = eval_modes_in[i];
1100ca735530SJeremy L Thompson 
1101004e4986SSebastian Grimberg       CeedCallCuda(ceed, cudaMemcpy(&asmb->d_B_in[i * elem_size_in * num_qpts_in], h_B_in, elem_size_in * num_qpts_in * sizeof(CeedScalar),
1102004e4986SSebastian Grimberg                                     cudaMemcpyHostToDevice));
1103004e4986SSebastian Grimberg     }
1104004e4986SSebastian Grimberg 
1105004e4986SSebastian Grimberg     if (identity) {
1106004e4986SSebastian Grimberg       CeedCallBackend(CeedFree(&identity));
1107cc132f9aSnbeams     }
1108cc132f9aSnbeams   }
1109cc132f9aSnbeams 
1110004e4986SSebastian Grimberg   // Load into B_out, in order that they will be used in eval_modes_out
1111004e4986SSebastian Grimberg   {
1112004e4986SSebastian Grimberg     const CeedInt out_bytes           = elem_size_out * num_qpts_out * num_eval_modes_out * sizeof(CeedScalar);
1113004e4986SSebastian Grimberg     CeedInt       d_out               = 0;
1114004e4986SSebastian Grimberg     CeedEvalMode  eval_modes_out_prev = CEED_EVAL_NONE;
1115004e4986SSebastian Grimberg     bool          has_eval_none       = false;
1116004e4986SSebastian Grimberg     CeedScalar   *identity            = NULL;
1117ca735530SJeremy L Thompson 
1118004e4986SSebastian Grimberg     for (CeedInt i = 0; i < num_eval_modes_out; i++) {
1119004e4986SSebastian Grimberg       has_eval_none = has_eval_none || (eval_modes_out[i] == CEED_EVAL_NONE);
1120004e4986SSebastian Grimberg     }
1121004e4986SSebastian Grimberg     if (has_eval_none) {
1122004e4986SSebastian Grimberg       CeedCallBackend(CeedCalloc(elem_size_out * num_qpts_out, &identity));
1123004e4986SSebastian Grimberg       for (CeedInt i = 0; i < (elem_size_out < num_qpts_out ? elem_size_out : num_qpts_out); i++) identity[i * elem_size_out + i] = 1.0;
1124cc132f9aSnbeams     }
1125cc132f9aSnbeams 
1126004e4986SSebastian Grimberg     CeedCallCuda(ceed, cudaMalloc((void **)&asmb->d_B_out, out_bytes));
1127004e4986SSebastian Grimberg     for (CeedInt i = 0; i < num_eval_modes_out; i++) {
1128004e4986SSebastian Grimberg       const CeedScalar *h_B_out;
1129ca735530SJeremy L Thompson 
1130004e4986SSebastian Grimberg       CeedCallBackend(CeedOperatorGetBasisPointer(basis_out, eval_modes_out[i], identity, &h_B_out));
1131004e4986SSebastian Grimberg       CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_out, eval_modes_out[i], &q_comp));
1132004e4986SSebastian Grimberg       if (q_comp > 1) {
1133004e4986SSebastian Grimberg         if (i == 0 || eval_modes_out[i] != eval_modes_out_prev) d_out = 0;
1134004e4986SSebastian Grimberg         else h_B_out = &h_B_out[(++d_out) * elem_size_out * num_qpts_out];
1135004e4986SSebastian Grimberg       }
1136004e4986SSebastian Grimberg       eval_modes_out_prev = eval_modes_out[i];
1137ca735530SJeremy L Thompson 
1138004e4986SSebastian Grimberg       CeedCallCuda(ceed, cudaMemcpy(&asmb->d_B_out[i * elem_size_out * num_qpts_out], h_B_out, elem_size_out * num_qpts_out * sizeof(CeedScalar),
1139004e4986SSebastian Grimberg                                     cudaMemcpyHostToDevice));
1140004e4986SSebastian Grimberg     }
1141004e4986SSebastian Grimberg 
1142004e4986SSebastian Grimberg     if (identity) {
1143004e4986SSebastian Grimberg       CeedCallBackend(CeedFree(&identity));
1144cc132f9aSnbeams     }
1145cc132f9aSnbeams   }
1146cc132f9aSnbeams   return CEED_ERROR_SUCCESS;
1147cc132f9aSnbeams }
1148cc132f9aSnbeams 
1149cc132f9aSnbeams //------------------------------------------------------------------------------
1150cefa2673SJeremy L Thompson // Assemble matrix data for COO matrix of assembled operator.
1151cefa2673SJeremy L Thompson // The sparsity pattern is set by CeedOperatorLinearAssembleSymbolic.
1152cefa2673SJeremy L Thompson //
1153ea61e9acSJeremy L Thompson // Note that this (and other assembly routines) currently assume only one active input restriction/basis per operator (could have multiple basis eval
1154cefa2673SJeremy L Thompson // modes).
1155cefa2673SJeremy L Thompson // TODO: allow multiple active input restrictions/basis objects
1156cc132f9aSnbeams //------------------------------------------------------------------------------
11572b730f8bSJeremy L Thompson static int CeedSingleOperatorAssemble_Cuda(CeedOperator op, CeedInt offset, CeedVector values) {
1158cc132f9aSnbeams   Ceed                ceed;
1159ca735530SJeremy L Thompson   CeedSize            values_length = 0, assembled_qf_length = 0;
1160004e4986SSebastian Grimberg   CeedInt             use_ceedsize_idx = 0, num_elem_in, num_elem_out, elem_size_in, elem_size_out;
1161ca735530SJeremy L Thompson   CeedScalar         *values_array;
1162004e4986SSebastian Grimberg   const CeedScalar   *assembled_qf_array;
1163ca735530SJeremy L Thompson   CeedVector          assembled_qf   = NULL;
1164004e4986SSebastian Grimberg   CeedElemRestriction assembled_rstr = NULL, rstr_in, rstr_out;
1165004e4986SSebastian Grimberg   CeedRestrictionType rstr_type_in, rstr_type_out;
1166004e4986SSebastian Grimberg   const bool         *orients_in = NULL, *orients_out = NULL;
1167004e4986SSebastian Grimberg   const CeedInt8     *curl_orients_in = NULL, *curl_orients_out = NULL;
1168cc132f9aSnbeams   CeedOperator_Cuda  *impl;
1169ca735530SJeremy L Thompson 
1170ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
11712b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &impl));
1172cc132f9aSnbeams 
1173cc132f9aSnbeams   // Assemble QFunction
1174004e4986SSebastian Grimberg   CeedCallBackend(CeedOperatorLinearAssembleQFunctionBuildOrUpdate(op, &assembled_qf, &assembled_rstr, CEED_REQUEST_IMMEDIATE));
1175004e4986SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionDestroy(&assembled_rstr));
1176004e4986SSebastian Grimberg   CeedCallBackend(CeedVectorGetArrayRead(assembled_qf, CEED_MEM_DEVICE, &assembled_qf_array));
1177cc132f9aSnbeams 
1178f7c1b517Snbeams   CeedCallBackend(CeedVectorGetLength(values, &values_length));
1179f7c1b517Snbeams   CeedCallBackend(CeedVectorGetLength(assembled_qf, &assembled_qf_length));
1180f7c1b517Snbeams   if ((values_length > INT_MAX) || (assembled_qf_length > INT_MAX)) use_ceedsize_idx = 1;
1181004e4986SSebastian Grimberg 
1182f7c1b517Snbeams   // Setup
1183004e4986SSebastian Grimberg   if (!impl->asmb) CeedCallBackend(CeedSingleOperatorAssembleSetup_Cuda(op, use_ceedsize_idx));
1184004e4986SSebastian Grimberg   CeedOperatorAssemble_Cuda *asmb = impl->asmb;
1185004e4986SSebastian Grimberg 
1186004e4986SSebastian Grimberg   assert(asmb != NULL);
1187004e4986SSebastian Grimberg 
1188004e4986SSebastian Grimberg   // Assemble element operator
1189004e4986SSebastian Grimberg   CeedCallBackend(CeedVectorGetArray(values, CEED_MEM_DEVICE, &values_array));
1190004e4986SSebastian Grimberg   values_array += offset;
1191004e4986SSebastian Grimberg 
1192004e4986SSebastian Grimberg   CeedCallBackend(CeedOperatorGetActiveElemRestrictions(op, &rstr_in, &rstr_out));
1193004e4986SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionGetNumElements(rstr_in, &num_elem_in));
1194004e4986SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionGetElementSize(rstr_in, &elem_size_in));
1195004e4986SSebastian Grimberg 
1196004e4986SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionGetType(rstr_in, &rstr_type_in));
1197004e4986SSebastian Grimberg   if (rstr_type_in == CEED_RESTRICTION_ORIENTED) {
1198004e4986SSebastian Grimberg     CeedCallBackend(CeedElemRestrictionGetOrientations(rstr_in, CEED_MEM_DEVICE, &orients_in));
1199004e4986SSebastian Grimberg   } else if (rstr_type_in == CEED_RESTRICTION_CURL_ORIENTED) {
1200004e4986SSebastian Grimberg     CeedCallBackend(CeedElemRestrictionGetCurlOrientations(rstr_in, CEED_MEM_DEVICE, &curl_orients_in));
1201004e4986SSebastian Grimberg   }
1202004e4986SSebastian Grimberg 
1203004e4986SSebastian Grimberg   if (rstr_in != rstr_out) {
1204004e4986SSebastian Grimberg     CeedCallBackend(CeedElemRestrictionGetNumElements(rstr_out, &num_elem_out));
1205004e4986SSebastian Grimberg     CeedCheck(num_elem_in == num_elem_out, ceed, CEED_ERROR_UNSUPPORTED,
1206004e4986SSebastian Grimberg               "Active input and output operator restrictions must have the same number of elements");
1207004e4986SSebastian Grimberg     CeedCallBackend(CeedElemRestrictionGetElementSize(rstr_out, &elem_size_out));
1208004e4986SSebastian Grimberg 
1209004e4986SSebastian Grimberg     CeedCallBackend(CeedElemRestrictionGetType(rstr_out, &rstr_type_out));
1210004e4986SSebastian Grimberg     if (rstr_type_out == CEED_RESTRICTION_ORIENTED) {
1211004e4986SSebastian Grimberg       CeedCallBackend(CeedElemRestrictionGetOrientations(rstr_out, CEED_MEM_DEVICE, &orients_out));
1212004e4986SSebastian Grimberg     } else if (rstr_type_out == CEED_RESTRICTION_CURL_ORIENTED) {
1213004e4986SSebastian Grimberg       CeedCallBackend(CeedElemRestrictionGetCurlOrientations(rstr_out, CEED_MEM_DEVICE, &curl_orients_out));
1214004e4986SSebastian Grimberg     }
1215004e4986SSebastian Grimberg   } else {
1216004e4986SSebastian Grimberg     elem_size_out    = elem_size_in;
1217004e4986SSebastian Grimberg     orients_out      = orients_in;
1218004e4986SSebastian Grimberg     curl_orients_out = curl_orients_in;
1219f7c1b517Snbeams   }
1220f7c1b517Snbeams 
1221cc132f9aSnbeams   // Compute B^T D B
1222004e4986SSebastian Grimberg   CeedInt shared_mem =
1223004e4986SSebastian Grimberg       ((curl_orients_in || curl_orients_out ? elem_size_in * elem_size_out : 0) + (curl_orients_in ? elem_size_in * asmb->block_size_y : 0)) *
1224004e4986SSebastian Grimberg       sizeof(CeedScalar);
1225004e4986SSebastian Grimberg   CeedInt grid   = CeedDivUpInt(num_elem_in, asmb->elems_per_block);
1226004e4986SSebastian Grimberg   void   *args[] = {(void *)&num_elem_in, &asmb->d_B_in,     &asmb->d_B_out,      &orients_in,  &curl_orients_in,
1227004e4986SSebastian Grimberg                     &orients_out,         &curl_orients_out, &assembled_qf_array, &values_array};
1228ca735530SJeremy L Thompson 
12292b730f8bSJeremy L Thompson   CeedCallBackend(
1230004e4986SSebastian Grimberg       CeedRunKernelDimShared_Cuda(ceed, asmb->LinearAssemble, grid, asmb->block_size_x, asmb->block_size_y, asmb->elems_per_block, shared_mem, args));
1231cc132f9aSnbeams 
1232cc132f9aSnbeams   // Restore arrays
12332b730f8bSJeremy L Thompson   CeedCallBackend(CeedVectorRestoreArray(values, &values_array));
1234004e4986SSebastian Grimberg   CeedCallBackend(CeedVectorRestoreArrayRead(assembled_qf, &assembled_qf_array));
1235cc132f9aSnbeams 
1236cc132f9aSnbeams   // Cleanup
12372b730f8bSJeremy L Thompson   CeedCallBackend(CeedVectorDestroy(&assembled_qf));
1238004e4986SSebastian Grimberg   if (rstr_type_in == CEED_RESTRICTION_ORIENTED) {
1239004e4986SSebastian Grimberg     CeedCallBackend(CeedElemRestrictionRestoreOrientations(rstr_in, &orients_in));
1240004e4986SSebastian Grimberg   } else if (rstr_type_in == CEED_RESTRICTION_CURL_ORIENTED) {
1241004e4986SSebastian Grimberg     CeedCallBackend(CeedElemRestrictionRestoreCurlOrientations(rstr_in, &curl_orients_in));
1242004e4986SSebastian Grimberg   }
1243004e4986SSebastian Grimberg   if (rstr_in != rstr_out) {
1244004e4986SSebastian Grimberg     if (rstr_type_out == CEED_RESTRICTION_ORIENTED) {
1245004e4986SSebastian Grimberg       CeedCallBackend(CeedElemRestrictionRestoreOrientations(rstr_out, &orients_out));
1246004e4986SSebastian Grimberg     } else if (rstr_type_out == CEED_RESTRICTION_CURL_ORIENTED) {
1247004e4986SSebastian Grimberg       CeedCallBackend(CeedElemRestrictionRestoreCurlOrientations(rstr_out, &curl_orients_out));
1248004e4986SSebastian Grimberg     }
1249004e4986SSebastian Grimberg   }
1250cc132f9aSnbeams   return CEED_ERROR_SUCCESS;
1251cc132f9aSnbeams }
1252cc132f9aSnbeams 
1253cc132f9aSnbeams //------------------------------------------------------------------------------
12540d0321e0SJeremy L Thompson // Create operator
12550d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
12560d0321e0SJeremy L Thompson int CeedOperatorCreate_Cuda(CeedOperator op) {
12570d0321e0SJeremy L Thompson   Ceed               ceed;
12580d0321e0SJeremy L Thompson   CeedOperator_Cuda *impl;
12590d0321e0SJeremy L Thompson 
1260ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
12612b730f8bSJeremy L Thompson   CeedCallBackend(CeedCalloc(1, &impl));
12622b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorSetData(op, impl));
12630d0321e0SJeremy L Thompson 
12642b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleQFunction", CeedOperatorLinearAssembleQFunction_Cuda));
12652b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleQFunctionUpdate", CeedOperatorLinearAssembleQFunctionUpdate_Cuda));
12662b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleAddDiagonal", CeedOperatorLinearAssembleAddDiagonal_Cuda));
12672b730f8bSJeremy L Thompson   CeedCallBackend(
12682b730f8bSJeremy L Thompson       CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleAddPointBlockDiagonal", CeedOperatorLinearAssembleAddPointBlockDiagonal_Cuda));
12692b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleSingle", CeedSingleOperatorAssemble_Cuda));
12702b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "ApplyAdd", CeedOperatorApplyAdd_Cuda));
12712b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "Destroy", CeedOperatorDestroy_Cuda));
12720d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
12730d0321e0SJeremy L Thompson }
12740d0321e0SJeremy L Thompson 
12750d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
1276