xref: /libCEED/backends/cuda-ref/ceed-cuda-ref-operator.c (revision 6e536b992ff6bc401c55631f1bc4464446496b52)
13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
30d0321e0SJeremy L Thompson //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
50d0321e0SJeremy L Thompson //
63d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
70d0321e0SJeremy L Thompson 
849aac155SJeremy L Thompson #include <ceed.h>
92b730f8bSJeremy L Thompson #include <ceed/backend.h>
102b730f8bSJeremy L Thompson #include <ceed/jit-tools.h>
11c85e8640SSebastian Grimberg #include <assert.h>
120d0321e0SJeremy L Thompson #include <cuda.h>
130d0321e0SJeremy L Thompson #include <cuda_runtime.h>
140d0321e0SJeremy L Thompson #include <stdbool.h>
150d0321e0SJeremy L Thompson #include <string.h>
162b730f8bSJeremy L Thompson 
1749aac155SJeremy L Thompson #include "../cuda/ceed-cuda-common.h"
180d0321e0SJeremy L Thompson #include "../cuda/ceed-cuda-compile.h"
192b730f8bSJeremy L Thompson #include "ceed-cuda-ref.h"
200d0321e0SJeremy L Thompson 
210d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
220d0321e0SJeremy L Thompson // Destroy operator
230d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
240d0321e0SJeremy L Thompson static int CeedOperatorDestroy_Cuda(CeedOperator op) {
250d0321e0SJeremy L Thompson   CeedOperator_Cuda *impl;
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));
340d0321e0SJeremy L Thompson 
35ca735530SJeremy L Thompson   for (CeedInt i = 0; i < impl->num_inputs; i++) {
36ca735530SJeremy L Thompson     CeedCallBackend(CeedVectorDestroy(&impl->q_vecs_in[i]));
370d0321e0SJeremy L Thompson   }
38ca735530SJeremy L Thompson   CeedCallBackend(CeedFree(&impl->q_vecs_in));
390d0321e0SJeremy L Thompson 
40ca735530SJeremy L Thompson   for (CeedInt i = 0; i < impl->num_outputs; i++) {
41ca735530SJeremy L Thompson     CeedCallBackend(CeedVectorDestroy(&impl->q_vecs_out[i]));
420d0321e0SJeremy L Thompson   }
43ca735530SJeremy L Thompson   CeedCallBackend(CeedFree(&impl->q_vecs_out));
440d0321e0SJeremy L Thompson 
450d0321e0SJeremy L Thompson   // QFunction assembly data
46ca735530SJeremy L Thompson   for (CeedInt i = 0; i < impl->num_active_in; i++) {
47ca735530SJeremy L Thompson     CeedCallBackend(CeedVectorDestroy(&impl->qf_active_in[i]));
480d0321e0SJeremy L Thompson   }
49ca735530SJeremy L Thompson   CeedCallBackend(CeedFree(&impl->qf_active_in));
500d0321e0SJeremy L Thompson 
510d0321e0SJeremy L Thompson   // Diag data
520d0321e0SJeremy L Thompson   if (impl->diag) {
530d0321e0SJeremy L Thompson     Ceed ceed;
54ca735530SJeremy L Thompson 
552b730f8bSJeremy L Thompson     CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
56cbfe683aSSebastian Grimberg     if (impl->diag->module) {
572b730f8bSJeremy L Thompson       CeedCallCuda(ceed, cuModuleUnload(impl->diag->module));
58cbfe683aSSebastian Grimberg     }
59cbfe683aSSebastian Grimberg     if (impl->diag->module_point_block) {
60cbfe683aSSebastian Grimberg       CeedCallCuda(ceed, cuModuleUnload(impl->diag->module_point_block));
61cbfe683aSSebastian Grimberg     }
62004e4986SSebastian Grimberg     CeedCallCuda(ceed, cudaFree(impl->diag->d_eval_modes_in));
63004e4986SSebastian Grimberg     CeedCallCuda(ceed, cudaFree(impl->diag->d_eval_modes_out));
642b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaFree(impl->diag->d_identity));
65ca735530SJeremy L Thompson     CeedCallCuda(ceed, cudaFree(impl->diag->d_interp_in));
66ca735530SJeremy L Thompson     CeedCallCuda(ceed, cudaFree(impl->diag->d_interp_out));
67ca735530SJeremy L Thompson     CeedCallCuda(ceed, cudaFree(impl->diag->d_grad_in));
68ca735530SJeremy L Thompson     CeedCallCuda(ceed, cudaFree(impl->diag->d_grad_out));
69004e4986SSebastian Grimberg     CeedCallCuda(ceed, cudaFree(impl->diag->d_div_in));
70004e4986SSebastian Grimberg     CeedCallCuda(ceed, cudaFree(impl->diag->d_div_out));
71004e4986SSebastian Grimberg     CeedCallCuda(ceed, cudaFree(impl->diag->d_curl_in));
72004e4986SSebastian Grimberg     CeedCallCuda(ceed, cudaFree(impl->diag->d_curl_out));
73004e4986SSebastian Grimberg     CeedCallBackend(CeedElemRestrictionDestroy(&impl->diag->diag_rstr));
74506b1a0cSSebastian Grimberg     CeedCallBackend(CeedElemRestrictionDestroy(&impl->diag->point_block_diag_rstr));
75ca735530SJeremy L Thompson     CeedCallBackend(CeedVectorDestroy(&impl->diag->elem_diag));
76ca735530SJeremy L Thompson     CeedCallBackend(CeedVectorDestroy(&impl->diag->point_block_elem_diag));
770d0321e0SJeremy L Thompson   }
782b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&impl->diag));
790d0321e0SJeremy L Thompson 
80cc132f9aSnbeams   if (impl->asmb) {
81cc132f9aSnbeams     Ceed ceed;
82ca735530SJeremy L Thompson 
832b730f8bSJeremy L Thompson     CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
842b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cuModuleUnload(impl->asmb->module));
852b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaFree(impl->asmb->d_B_in));
862b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaFree(impl->asmb->d_B_out));
87cc132f9aSnbeams   }
882b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&impl->asmb));
89cc132f9aSnbeams 
902b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&impl));
910d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
920d0321e0SJeremy L Thompson }
930d0321e0SJeremy L Thompson 
940d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
950d0321e0SJeremy L Thompson // Setup infields or outfields
960d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
97004e4986SSebastian Grimberg static int CeedOperatorSetupFields_Cuda(CeedQFunction qf, CeedOperator op, bool is_input, CeedVector *e_vecs, CeedVector *q_vecs, CeedInt start_e,
98ca735530SJeremy L Thompson                                         CeedInt num_fields, CeedInt Q, CeedInt num_elem) {
990d0321e0SJeremy L Thompson   Ceed                ceed;
100ca735530SJeremy L Thompson   CeedQFunctionField *qf_fields;
101ca735530SJeremy L Thompson   CeedOperatorField  *op_fields;
1020d0321e0SJeremy L Thompson 
103ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
104ca735530SJeremy L Thompson   if (is_input) {
105ca735530SJeremy L Thompson     CeedCallBackend(CeedOperatorGetFields(op, NULL, &op_fields, NULL, NULL));
106ca735530SJeremy L Thompson     CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_fields, NULL, NULL));
1070d0321e0SJeremy L Thompson   } else {
108ca735530SJeremy L Thompson     CeedCallBackend(CeedOperatorGetFields(op, NULL, NULL, NULL, &op_fields));
109ca735530SJeremy L Thompson     CeedCallBackend(CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qf_fields));
1100d0321e0SJeremy L Thompson   }
1110d0321e0SJeremy L Thompson 
1120d0321e0SJeremy L Thompson   // Loop over fields
113ca735530SJeremy L Thompson   for (CeedInt i = 0; i < num_fields; i++) {
114004e4986SSebastian Grimberg     bool         is_strided = false, skip_restriction = false;
115004e4986SSebastian Grimberg     CeedSize     q_size;
116004e4986SSebastian Grimberg     CeedInt      size;
117004e4986SSebastian Grimberg     CeedEvalMode eval_mode;
118ca735530SJeremy L Thompson     CeedBasis    basis;
1190d0321e0SJeremy L Thompson 
120004e4986SSebastian Grimberg     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode));
121004e4986SSebastian Grimberg     if (eval_mode != CEED_EVAL_WEIGHT) {
122edb2538eSJeremy L Thompson       CeedElemRestriction elem_rstr;
123ca735530SJeremy L Thompson 
1240d0321e0SJeremy L Thompson       // Check whether this field can skip the element restriction:
125004e4986SSebastian Grimberg       // Must be passive input, with eval_mode NONE, and have a strided restriction with CEED_STRIDES_BACKEND.
126004e4986SSebastian Grimberg       CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_fields[i], &elem_rstr));
1270d0321e0SJeremy L Thompson 
1280d0321e0SJeremy L Thompson       // First, check whether the field is input or output:
129ca735530SJeremy L Thompson       if (is_input) {
130ca735530SJeremy L Thompson         CeedVector vec;
131ca735530SJeremy L Thompson 
132004e4986SSebastian Grimberg         // Check for passive input
133ca735530SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec));
134ca735530SJeremy L Thompson         if (vec != CEED_VECTOR_ACTIVE) {
135004e4986SSebastian Grimberg           // Check eval_mode
136004e4986SSebastian Grimberg           if (eval_mode == CEED_EVAL_NONE) {
1370d0321e0SJeremy L Thompson             // Check for strided restriction
138edb2538eSJeremy L Thompson             CeedCallBackend(CeedElemRestrictionIsStrided(elem_rstr, &is_strided));
139ca735530SJeremy L Thompson             if (is_strided) {
1400d0321e0SJeremy L Thompson               // Check if vector is already in preferred backend ordering
141edb2538eSJeremy L Thompson               CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &skip_restriction));
1420d0321e0SJeremy L Thompson             }
1430d0321e0SJeremy L Thompson           }
1440d0321e0SJeremy L Thompson         }
1450d0321e0SJeremy L Thompson       }
146ca735530SJeremy L Thompson       if (skip_restriction) {
147ea61e9acSJeremy L Thompson         // We do not need an E-Vector, but will use the input field vector's data directly in the operator application.
148004e4986SSebastian Grimberg         e_vecs[i + start_e] = NULL;
1490d0321e0SJeremy L Thompson       } else {
150004e4986SSebastian Grimberg         CeedCallBackend(CeedElemRestrictionCreateVector(elem_rstr, NULL, &e_vecs[i + start_e]));
1510d0321e0SJeremy L Thompson       }
1520d0321e0SJeremy L Thompson     }
1530d0321e0SJeremy L Thompson 
154004e4986SSebastian Grimberg     switch (eval_mode) {
1550d0321e0SJeremy L Thompson       case CEED_EVAL_NONE:
156ca735530SJeremy L Thompson         CeedCallBackend(CeedQFunctionFieldGetSize(qf_fields[i], &size));
157ca735530SJeremy L Thompson         q_size = (CeedSize)num_elem * Q * size;
158ca735530SJeremy L Thompson         CeedCallBackend(CeedVectorCreate(ceed, q_size, &q_vecs[i]));
1590d0321e0SJeremy L Thompson         break;
1600d0321e0SJeremy L Thompson       case CEED_EVAL_INTERP:
1610d0321e0SJeremy L Thompson       case CEED_EVAL_GRAD:
162004e4986SSebastian Grimberg       case CEED_EVAL_DIV:
163004e4986SSebastian Grimberg       case CEED_EVAL_CURL:
164ca735530SJeremy L Thompson         CeedCallBackend(CeedQFunctionFieldGetSize(qf_fields[i], &size));
165ca735530SJeremy L Thompson         q_size = (CeedSize)num_elem * Q * size;
166ca735530SJeremy L Thompson         CeedCallBackend(CeedVectorCreate(ceed, q_size, &q_vecs[i]));
1670d0321e0SJeremy L Thompson         break;
1680d0321e0SJeremy L Thompson       case CEED_EVAL_WEIGHT:  // Only on input fields
169ca735530SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis));
170ca735530SJeremy L Thompson         q_size = (CeedSize)num_elem * Q;
171ca735530SJeremy L Thompson         CeedCallBackend(CeedVectorCreate(ceed, q_size, &q_vecs[i]));
172ca735530SJeremy L Thompson         CeedCallBackend(CeedBasisApply(basis, num_elem, CEED_NOTRANSPOSE, CEED_EVAL_WEIGHT, CEED_VECTOR_NONE, q_vecs[i]));
1730d0321e0SJeremy L Thompson         break;
1740d0321e0SJeremy L Thompson     }
1750d0321e0SJeremy L Thompson   }
1760d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1770d0321e0SJeremy L Thompson }
1780d0321e0SJeremy L Thompson 
1790d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
180ea61e9acSJeremy L Thompson // CeedOperator needs to connect all the named fields (be they active or passive) to the named inputs and outputs of its CeedQFunction.
1810d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
1820d0321e0SJeremy L Thompson static int CeedOperatorSetup_Cuda(CeedOperator op) {
1830d0321e0SJeremy L Thompson   Ceed                ceed;
184ca735530SJeremy L Thompson   bool                is_setup_done;
185ca735530SJeremy L Thompson   CeedInt             Q, num_elem, num_input_fields, num_output_fields;
186ca735530SJeremy L Thompson   CeedQFunctionField *qf_input_fields, *qf_output_fields;
1870d0321e0SJeremy L Thompson   CeedQFunction       qf;
188ca735530SJeremy L Thompson   CeedOperatorField  *op_input_fields, *op_output_fields;
189ca735530SJeremy L Thompson   CeedOperator_Cuda  *impl;
190ca735530SJeremy L Thompson 
191ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorIsSetupDone(op, &is_setup_done));
192ca735530SJeremy L Thompson   if (is_setup_done) return CEED_ERROR_SUCCESS;
193ca735530SJeremy L Thompson 
194ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
195ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &impl));
1962b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
1972b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q));
198ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem));
199ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
200ca735530SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
2010d0321e0SJeremy L Thompson 
2020d0321e0SJeremy L Thompson   // Allocate
203ca735530SJeremy L Thompson   CeedCallBackend(CeedCalloc(num_input_fields + num_output_fields, &impl->e_vecs));
204ca735530SJeremy L Thompson   CeedCallBackend(CeedCalloc(CEED_FIELD_MAX, &impl->q_vecs_in));
205ca735530SJeremy L Thompson   CeedCallBackend(CeedCalloc(CEED_FIELD_MAX, &impl->q_vecs_out));
206ca735530SJeremy L Thompson   impl->num_inputs  = num_input_fields;
207ca735530SJeremy L Thompson   impl->num_outputs = num_output_fields;
2080d0321e0SJeremy L Thompson 
209ca735530SJeremy L Thompson   // Set up infield and outfield e_vecs and q_vecs
2100d0321e0SJeremy L Thompson   // Infields
211ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorSetupFields_Cuda(qf, op, true, impl->e_vecs, impl->q_vecs_in, 0, num_input_fields, Q, num_elem));
2120d0321e0SJeremy L Thompson   // Outfields
213ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorSetupFields_Cuda(qf, op, false, impl->e_vecs, impl->q_vecs_out, num_input_fields, num_output_fields, Q, num_elem));
2140d0321e0SJeremy L Thompson 
2152b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorSetSetupDone(op));
2160d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2170d0321e0SJeremy L Thompson }
2180d0321e0SJeremy L Thompson 
2190d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
2200d0321e0SJeremy L Thompson // Setup Operator Inputs
2210d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
222ca735530SJeremy L Thompson static inline int CeedOperatorSetupInputs_Cuda(CeedInt num_input_fields, CeedQFunctionField *qf_input_fields, CeedOperatorField *op_input_fields,
223004e4986SSebastian Grimberg                                                CeedVector in_vec, const bool skip_active, CeedScalar *e_data[2 * CEED_FIELD_MAX],
2240d0321e0SJeremy L Thompson                                                CeedOperator_Cuda *impl, CeedRequest *request) {
225ca735530SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
226004e4986SSebastian Grimberg     CeedEvalMode        eval_mode;
2270d0321e0SJeremy L Thompson     CeedVector          vec;
228edb2538eSJeremy L Thompson     CeedElemRestriction elem_rstr;
2290d0321e0SJeremy L Thompson 
2300d0321e0SJeremy L Thompson     // Get input vector
231ca735530SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec));
2320d0321e0SJeremy L Thompson     if (vec == CEED_VECTOR_ACTIVE) {
233004e4986SSebastian Grimberg       if (skip_active) continue;
234ca735530SJeremy L Thompson       else vec = in_vec;
2350d0321e0SJeremy L Thompson     }
2360d0321e0SJeremy L Thompson 
237004e4986SSebastian Grimberg     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
238004e4986SSebastian Grimberg     if (eval_mode == CEED_EVAL_WEIGHT) {  // Skip
2390d0321e0SJeremy L Thompson     } else {
240004e4986SSebastian Grimberg       // Get input vector
241004e4986SSebastian Grimberg       CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec));
2420d0321e0SJeremy L Thompson       // Get input element restriction
243edb2538eSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr));
244ca735530SJeremy L Thompson       if (vec == CEED_VECTOR_ACTIVE) vec = in_vec;
2450d0321e0SJeremy L Thompson       // Restrict, if necessary
246ca735530SJeremy L Thompson       if (!impl->e_vecs[i]) {
2470d0321e0SJeremy L Thompson         // No restriction for this field; read data directly from vec.
248ca735530SJeremy L Thompson         CeedCallBackend(CeedVectorGetArrayRead(vec, CEED_MEM_DEVICE, (const CeedScalar **)&e_data[i]));
2490d0321e0SJeremy L Thompson       } else {
250edb2538eSJeremy L Thompson         CeedCallBackend(CeedElemRestrictionApply(elem_rstr, CEED_NOTRANSPOSE, vec, impl->e_vecs[i], request));
2510d0321e0SJeremy L Thompson         // Get evec
252ca735530SJeremy L Thompson         CeedCallBackend(CeedVectorGetArrayRead(impl->e_vecs[i], CEED_MEM_DEVICE, (const CeedScalar **)&e_data[i]));
2530d0321e0SJeremy L Thompson       }
2540d0321e0SJeremy L Thompson     }
2550d0321e0SJeremy L Thompson   }
2560d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2570d0321e0SJeremy L Thompson }
2580d0321e0SJeremy L Thompson 
2590d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
2600d0321e0SJeremy L Thompson // Input Basis Action
2610d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
262ca735530SJeremy L Thompson static inline int CeedOperatorInputBasis_Cuda(CeedInt num_elem, CeedQFunctionField *qf_input_fields, CeedOperatorField *op_input_fields,
263004e4986SSebastian Grimberg                                               CeedInt num_input_fields, const bool skip_active, CeedScalar *e_data[2 * CEED_FIELD_MAX],
2642b730f8bSJeremy L Thompson                                               CeedOperator_Cuda *impl) {
265ca735530SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
266ca735530SJeremy L Thompson     CeedInt             elem_size, size;
267004e4986SSebastian Grimberg     CeedEvalMode        eval_mode;
268edb2538eSJeremy L Thompson     CeedElemRestriction elem_rstr;
2690d0321e0SJeremy L Thompson     CeedBasis           basis;
2700d0321e0SJeremy L Thompson 
2710d0321e0SJeremy L Thompson     // Skip active input
272004e4986SSebastian Grimberg     if (skip_active) {
2730d0321e0SJeremy L Thompson       CeedVector vec;
274ca735530SJeremy L Thompson 
275ca735530SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec));
2762b730f8bSJeremy L Thompson       if (vec == CEED_VECTOR_ACTIVE) continue;
2770d0321e0SJeremy L Thompson     }
278004e4986SSebastian Grimberg     // Get elem_size, eval_mode, size
279edb2538eSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr));
280edb2538eSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
281004e4986SSebastian Grimberg     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
282ca735530SJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetSize(qf_input_fields[i], &size));
2830d0321e0SJeremy L Thompson     // Basis action
284004e4986SSebastian Grimberg     switch (eval_mode) {
2850d0321e0SJeremy L Thompson       case CEED_EVAL_NONE:
286ca735530SJeremy L Thompson         CeedCallBackend(CeedVectorSetArray(impl->q_vecs_in[i], CEED_MEM_DEVICE, CEED_USE_POINTER, e_data[i]));
2870d0321e0SJeremy L Thompson         break;
2880d0321e0SJeremy L Thompson       case CEED_EVAL_INTERP:
2890d0321e0SJeremy L Thompson       case CEED_EVAL_GRAD:
290004e4986SSebastian Grimberg       case CEED_EVAL_DIV:
291004e4986SSebastian Grimberg       case CEED_EVAL_CURL:
292ca735530SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
293004e4986SSebastian Grimberg         CeedCallBackend(CeedBasisApply(basis, num_elem, CEED_NOTRANSPOSE, eval_mode, impl->e_vecs[i], impl->q_vecs_in[i]));
2940d0321e0SJeremy L Thompson         break;
2950d0321e0SJeremy L Thompson       case CEED_EVAL_WEIGHT:
2960d0321e0SJeremy L Thompson         break;  // No action
2970d0321e0SJeremy L Thompson     }
2980d0321e0SJeremy L Thompson   }
2990d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3000d0321e0SJeremy L Thompson }
3010d0321e0SJeremy L Thompson 
3020d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
3030d0321e0SJeremy L Thompson // Restore Input Vectors
3040d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
305ca735530SJeremy L Thompson static inline int CeedOperatorRestoreInputs_Cuda(CeedInt num_input_fields, CeedQFunctionField *qf_input_fields, CeedOperatorField *op_input_fields,
306004e4986SSebastian Grimberg                                                  const bool skip_active, CeedScalar *e_data[2 * CEED_FIELD_MAX], CeedOperator_Cuda *impl) {
307ca735530SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
308004e4986SSebastian Grimberg     CeedEvalMode eval_mode;
3090d0321e0SJeremy L Thompson     CeedVector   vec;
3100d0321e0SJeremy L Thompson 
3110d0321e0SJeremy L Thompson     // Skip active input
312004e4986SSebastian Grimberg     if (skip_active) {
313ca735530SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec));
3142b730f8bSJeremy L Thompson       if (vec == CEED_VECTOR_ACTIVE) continue;
3150d0321e0SJeremy L Thompson     }
316004e4986SSebastian Grimberg     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
317004e4986SSebastian Grimberg     if (eval_mode == CEED_EVAL_WEIGHT) {  // Skip
3180d0321e0SJeremy L Thompson     } else {
319ca735530SJeremy L Thompson       if (!impl->e_vecs[i]) {  // This was a skip_restriction case
320ca735530SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec));
321ca735530SJeremy L Thompson         CeedCallBackend(CeedVectorRestoreArrayRead(vec, (const CeedScalar **)&e_data[i]));
3220d0321e0SJeremy L Thompson       } else {
323ca735530SJeremy L Thompson         CeedCallBackend(CeedVectorRestoreArrayRead(impl->e_vecs[i], (const CeedScalar **)&e_data[i]));
3240d0321e0SJeremy L Thompson       }
3250d0321e0SJeremy L Thompson     }
3260d0321e0SJeremy L Thompson   }
3270d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3280d0321e0SJeremy L Thompson }
3290d0321e0SJeremy L Thompson 
3300d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
3310d0321e0SJeremy L Thompson // Apply and add to output
3320d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
333ca735530SJeremy L Thompson static int CeedOperatorApplyAdd_Cuda(CeedOperator op, CeedVector in_vec, CeedVector out_vec, CeedRequest *request) {
334ca735530SJeremy L Thompson   CeedInt             Q, num_elem, elem_size, num_input_fields, num_output_fields, size;
335ca735530SJeremy L Thompson   CeedScalar         *e_data[2 * CEED_FIELD_MAX] = {NULL};
336ca735530SJeremy L Thompson   CeedQFunctionField *qf_input_fields, *qf_output_fields;
3370d0321e0SJeremy L Thompson   CeedQFunction       qf;
338004e4986SSebastian Grimberg   CeedOperatorField  *op_input_fields, *op_output_fields;
339004e4986SSebastian Grimberg   CeedOperator_Cuda  *impl;
340ca735530SJeremy L Thompson 
341ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &impl));
3422b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
3432b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q));
344ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem));
345ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
346ca735530SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
3470d0321e0SJeremy L Thompson 
3480d0321e0SJeremy L Thompson   // Setup
3492b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorSetup_Cuda(op));
3500d0321e0SJeremy L Thompson 
351004e4986SSebastian Grimberg   // Input Evecs and Restriction
352ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorSetupInputs_Cuda(num_input_fields, qf_input_fields, op_input_fields, in_vec, false, e_data, impl, request));
3530d0321e0SJeremy L Thompson 
3540d0321e0SJeremy L Thompson   // Input basis apply if needed
355ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorInputBasis_Cuda(num_elem, qf_input_fields, op_input_fields, num_input_fields, false, e_data, impl));
3560d0321e0SJeremy L Thompson 
3570d0321e0SJeremy L Thompson   // Output pointers, as necessary
358ca735530SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
359004e4986SSebastian Grimberg     CeedEvalMode eval_mode;
360004e4986SSebastian Grimberg 
361004e4986SSebastian Grimberg     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
362004e4986SSebastian Grimberg     if (eval_mode == CEED_EVAL_NONE) {
3630d0321e0SJeremy L Thompson       // Set the output Q-Vector to use the E-Vector data directly.
364ca735530SJeremy L Thompson       CeedCallBackend(CeedVectorGetArrayWrite(impl->e_vecs[i + impl->num_inputs], CEED_MEM_DEVICE, &e_data[i + num_input_fields]));
365ca735530SJeremy L Thompson       CeedCallBackend(CeedVectorSetArray(impl->q_vecs_out[i], CEED_MEM_DEVICE, CEED_USE_POINTER, e_data[i + num_input_fields]));
3660d0321e0SJeremy L Thompson     }
3670d0321e0SJeremy L Thompson   }
3680d0321e0SJeremy L Thompson 
3690d0321e0SJeremy L Thompson   // Q function
370ca735530SJeremy L Thompson   CeedCallBackend(CeedQFunctionApply(qf, num_elem * Q, impl->q_vecs_in, impl->q_vecs_out));
3710d0321e0SJeremy L Thompson 
3720d0321e0SJeremy L Thompson   // Output basis apply if needed
373ca735530SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
374004e4986SSebastian Grimberg     CeedEvalMode        eval_mode;
375edb2538eSJeremy L Thompson     CeedElemRestriction elem_rstr;
376ca735530SJeremy L Thompson     CeedBasis           basis;
377ca735530SJeremy L Thompson 
378004e4986SSebastian Grimberg     // Get elem_size, eval_mode, size
379edb2538eSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &elem_rstr));
380edb2538eSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
381004e4986SSebastian Grimberg     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
382ca735530SJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetSize(qf_output_fields[i], &size));
3830d0321e0SJeremy L Thompson     // Basis action
384004e4986SSebastian Grimberg     switch (eval_mode) {
3850d0321e0SJeremy L Thompson       case CEED_EVAL_NONE:
386004e4986SSebastian Grimberg         break;  // No action
3870d0321e0SJeremy L Thompson       case CEED_EVAL_INTERP:
3880d0321e0SJeremy L Thompson       case CEED_EVAL_GRAD:
389004e4986SSebastian Grimberg       case CEED_EVAL_DIV:
390004e4986SSebastian Grimberg       case CEED_EVAL_CURL:
391ca735530SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
392004e4986SSebastian Grimberg         CeedCallBackend(CeedBasisApply(basis, num_elem, CEED_TRANSPOSE, eval_mode, impl->q_vecs_out[i], impl->e_vecs[i + impl->num_inputs]));
3930d0321e0SJeremy L Thompson         break;
3940d0321e0SJeremy L Thompson       // LCOV_EXCL_START
3950d0321e0SJeremy L Thompson       case CEED_EVAL_WEIGHT: {
396*6e536b99SJeremy L Thompson         return CeedError(CeedOperatorReturnCeed(op), CEED_ERROR_BACKEND, "CEED_EVAL_WEIGHT cannot be an output evaluation mode");
3970d0321e0SJeremy L Thompson         // LCOV_EXCL_STOP
3980d0321e0SJeremy L Thompson       }
3990d0321e0SJeremy L Thompson     }
400004e4986SSebastian Grimberg   }
4010d0321e0SJeremy L Thompson 
4020d0321e0SJeremy L Thompson   // Output restriction
403ca735530SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
404004e4986SSebastian Grimberg     CeedEvalMode        eval_mode;
405ca735530SJeremy L Thompson     CeedVector          vec;
406edb2538eSJeremy L Thompson     CeedElemRestriction elem_rstr;
407ca735530SJeremy L Thompson 
4080d0321e0SJeremy L Thompson     // Restore evec
409004e4986SSebastian Grimberg     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
410004e4986SSebastian Grimberg     if (eval_mode == CEED_EVAL_NONE) {
411ca735530SJeremy L Thompson       CeedCallBackend(CeedVectorRestoreArray(impl->e_vecs[i + impl->num_inputs], &e_data[i + num_input_fields]));
4120d0321e0SJeremy L Thompson     }
4130d0321e0SJeremy L Thompson     // Get output vector
414ca735530SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[i], &vec));
4150d0321e0SJeremy L Thompson     // Restrict
416edb2538eSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &elem_rstr));
4170d0321e0SJeremy L Thompson     // Active
418ca735530SJeremy L Thompson     if (vec == CEED_VECTOR_ACTIVE) vec = out_vec;
4190d0321e0SJeremy L Thompson 
420edb2538eSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionApply(elem_rstr, CEED_TRANSPOSE, impl->e_vecs[i + impl->num_inputs], vec, request));
4210d0321e0SJeremy L Thompson   }
4220d0321e0SJeremy L Thompson 
4230d0321e0SJeremy L Thompson   // Restore input arrays
424ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorRestoreInputs_Cuda(num_input_fields, qf_input_fields, op_input_fields, false, e_data, impl));
4250d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
4260d0321e0SJeremy L Thompson }
4270d0321e0SJeremy L Thompson 
4280d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
429004e4986SSebastian Grimberg // Linear QFunction Assembly Core
4300d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
4312b730f8bSJeremy L Thompson static inline int CeedOperatorLinearAssembleQFunctionCore_Cuda(CeedOperator op, bool build_objects, CeedVector *assembled, CeedElemRestriction *rstr,
4320d0321e0SJeremy L Thompson                                                                CeedRequest *request) {
433ca735530SJeremy L Thompson   Ceed                ceed, ceed_parent;
434ca735530SJeremy L Thompson   CeedInt             num_active_in, num_active_out, Q, num_elem, num_input_fields, num_output_fields, size;
435ca735530SJeremy L Thompson   CeedScalar         *assembled_array, *e_data[2 * CEED_FIELD_MAX] = {NULL};
436ca735530SJeremy L Thompson   CeedVector         *active_inputs;
437ca735530SJeremy L Thompson   CeedQFunctionField *qf_input_fields, *qf_output_fields;
438ca735530SJeremy L Thompson   CeedQFunction       qf;
439ca735530SJeremy L Thompson   CeedOperatorField  *op_input_fields, *op_output_fields;
440ca735530SJeremy L Thompson   CeedOperator_Cuda  *impl;
441ca735530SJeremy L Thompson 
4422b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
443ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetFallbackParentCeed(op, &ceed_parent));
444e984cf9aSJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &impl));
445e984cf9aSJeremy L Thompson   CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q));
446ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem));
447e984cf9aSJeremy L Thompson   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
448ca735530SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
449ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
450ca735530SJeremy L Thompson   active_inputs = impl->qf_active_in;
451ca735530SJeremy L Thompson   num_active_in = impl->num_active_in, num_active_out = impl->num_active_out;
4520d0321e0SJeremy L Thompson 
4530d0321e0SJeremy L Thompson   // Setup
4542b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorSetup_Cuda(op));
4550d0321e0SJeremy L Thompson 
456004e4986SSebastian Grimberg   // Input Evecs and Restriction
457ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorSetupInputs_Cuda(num_input_fields, qf_input_fields, op_input_fields, NULL, true, e_data, impl, request));
4580d0321e0SJeremy L Thompson 
4590d0321e0SJeremy L Thompson   // Count number of active input fields
460ca735530SJeremy L Thompson   if (!num_active_in) {
461ca735530SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
462ca735530SJeremy L Thompson       CeedScalar *q_vec_array;
463ca735530SJeremy L Thompson       CeedVector  vec;
464ca735530SJeremy L Thompson 
4650d0321e0SJeremy L Thompson       // Get input vector
466ca735530SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec));
4670d0321e0SJeremy L Thompson       // Check if active input
4680d0321e0SJeremy L Thompson       if (vec == CEED_VECTOR_ACTIVE) {
469ca735530SJeremy L Thompson         CeedCallBackend(CeedQFunctionFieldGetSize(qf_input_fields[i], &size));
470ca735530SJeremy L Thompson         CeedCallBackend(CeedVectorSetValue(impl->q_vecs_in[i], 0.0));
471ca735530SJeremy L Thompson         CeedCallBackend(CeedVectorGetArray(impl->q_vecs_in[i], CEED_MEM_DEVICE, &q_vec_array));
472ca735530SJeremy L Thompson         CeedCallBackend(CeedRealloc(num_active_in + size, &active_inputs));
4730d0321e0SJeremy L Thompson         for (CeedInt field = 0; field < size; field++) {
474004e4986SSebastian Grimberg           CeedSize q_size = (CeedSize)Q * num_elem;
475004e4986SSebastian Grimberg 
476ca735530SJeremy L Thompson           CeedCallBackend(CeedVectorCreate(ceed, q_size, &active_inputs[num_active_in + field]));
477ca735530SJeremy L Thompson           CeedCallBackend(
478ca735530SJeremy L Thompson               CeedVectorSetArray(active_inputs[num_active_in + field], CEED_MEM_DEVICE, CEED_USE_POINTER, &q_vec_array[field * Q * num_elem]));
4790d0321e0SJeremy L Thompson         }
480ca735530SJeremy L Thompson         num_active_in += size;
481ca735530SJeremy L Thompson         CeedCallBackend(CeedVectorRestoreArray(impl->q_vecs_in[i], &q_vec_array));
4820d0321e0SJeremy L Thompson       }
4830d0321e0SJeremy L Thompson     }
484ca735530SJeremy L Thompson     impl->num_active_in = num_active_in;
485ca735530SJeremy L Thompson     impl->qf_active_in  = active_inputs;
4860d0321e0SJeremy L Thompson   }
4870d0321e0SJeremy L Thompson 
4880d0321e0SJeremy L Thompson   // Count number of active output fields
489ca735530SJeremy L Thompson   if (!num_active_out) {
490ca735530SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
491ca735530SJeremy L Thompson       CeedVector vec;
492ca735530SJeremy L Thompson 
4930d0321e0SJeremy L Thompson       // Get output vector
494ca735530SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[i], &vec));
4950d0321e0SJeremy L Thompson       // Check if active output
4960d0321e0SJeremy L Thompson       if (vec == CEED_VECTOR_ACTIVE) {
497ca735530SJeremy L Thompson         CeedCallBackend(CeedQFunctionFieldGetSize(qf_output_fields[i], &size));
498ca735530SJeremy L Thompson         num_active_out += size;
4990d0321e0SJeremy L Thompson       }
5000d0321e0SJeremy L Thompson     }
501ca735530SJeremy L Thompson     impl->num_active_out = num_active_out;
5020d0321e0SJeremy L Thompson   }
5030d0321e0SJeremy L Thompson 
5040d0321e0SJeremy L Thompson   // Check sizes
505ca735530SJeremy L Thompson   CeedCheck(num_active_in > 0 && num_active_out > 0, ceed, CEED_ERROR_BACKEND, "Cannot assemble QFunction without active inputs and outputs");
5060d0321e0SJeremy L Thompson 
5070d0321e0SJeremy L Thompson   // Build objects if needed
5080d0321e0SJeremy L Thompson   if (build_objects) {
509004e4986SSebastian Grimberg     CeedSize l_size     = (CeedSize)num_elem * Q * num_active_in * num_active_out;
510ca735530SJeremy L Thompson     CeedInt  strides[3] = {1, num_elem * Q, Q}; /* *NOPAD* */
511004e4986SSebastian Grimberg 
512004e4986SSebastian Grimberg     // Create output restriction
513ca735530SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionCreateStrided(ceed_parent, num_elem, Q, num_active_in * num_active_out,
514ca735530SJeremy L Thompson                                                      num_active_in * num_active_out * num_elem * Q, strides, rstr));
5150d0321e0SJeremy L Thompson     // Create assembled vector
516ca735530SJeremy L Thompson     CeedCallBackend(CeedVectorCreate(ceed_parent, l_size, assembled));
5170d0321e0SJeremy L Thompson   }
5182b730f8bSJeremy L Thompson   CeedCallBackend(CeedVectorSetValue(*assembled, 0.0));
519ca735530SJeremy L Thompson   CeedCallBackend(CeedVectorGetArray(*assembled, CEED_MEM_DEVICE, &assembled_array));
5200d0321e0SJeremy L Thompson 
5210d0321e0SJeremy L Thompson   // Input basis apply
522ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorInputBasis_Cuda(num_elem, qf_input_fields, op_input_fields, num_input_fields, true, e_data, impl));
5230d0321e0SJeremy L Thompson 
5240d0321e0SJeremy L Thompson   // Assemble QFunction
525ca735530SJeremy L Thompson   for (CeedInt in = 0; in < num_active_in; in++) {
5260d0321e0SJeremy L Thompson     // Set Inputs
527ca735530SJeremy L Thompson     CeedCallBackend(CeedVectorSetValue(active_inputs[in], 1.0));
528ca735530SJeremy L Thompson     if (num_active_in > 1) {
529ca735530SJeremy L Thompson       CeedCallBackend(CeedVectorSetValue(active_inputs[(in + num_active_in - 1) % num_active_in], 0.0));
5300d0321e0SJeremy L Thompson     }
5310d0321e0SJeremy L Thompson     // Set Outputs
532ca735530SJeremy L Thompson     for (CeedInt out = 0; out < num_output_fields; out++) {
533ca735530SJeremy L Thompson       CeedVector vec;
534ca735530SJeremy L Thompson 
5350d0321e0SJeremy L Thompson       // Get output vector
536ca735530SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[out], &vec));
5370d0321e0SJeremy L Thompson       // Check if active output
5380d0321e0SJeremy L Thompson       if (vec == CEED_VECTOR_ACTIVE) {
539ca735530SJeremy L Thompson         CeedCallBackend(CeedVectorSetArray(impl->q_vecs_out[out], CEED_MEM_DEVICE, CEED_USE_POINTER, assembled_array));
540ca735530SJeremy L Thompson         CeedCallBackend(CeedQFunctionFieldGetSize(qf_output_fields[out], &size));
541ca735530SJeremy L Thompson         assembled_array += size * Q * num_elem;  // Advance the pointer by the size of the output
5420d0321e0SJeremy L Thompson       }
5430d0321e0SJeremy L Thompson     }
5440d0321e0SJeremy L Thompson     // Apply QFunction
545ca735530SJeremy L Thompson     CeedCallBackend(CeedQFunctionApply(qf, Q * num_elem, impl->q_vecs_in, impl->q_vecs_out));
5460d0321e0SJeremy L Thompson   }
5470d0321e0SJeremy L Thompson 
548ca735530SJeremy L Thompson   // Un-set output q_vecs to prevent accidental overwrite of Assembled
549ca735530SJeremy L Thompson   for (CeedInt out = 0; out < num_output_fields; out++) {
550ca735530SJeremy L Thompson     CeedVector vec;
551ca735530SJeremy L Thompson 
5520d0321e0SJeremy L Thompson     // Get output vector
553ca735530SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[out], &vec));
5540d0321e0SJeremy L Thompson     // Check if active output
5550d0321e0SJeremy L Thompson     if (vec == CEED_VECTOR_ACTIVE) {
556ca735530SJeremy L Thompson       CeedCallBackend(CeedVectorTakeArray(impl->q_vecs_out[out], CEED_MEM_DEVICE, NULL));
5570d0321e0SJeremy L Thompson     }
5580d0321e0SJeremy L Thompson   }
5590d0321e0SJeremy L Thompson 
5600d0321e0SJeremy L Thompson   // Restore input arrays
561ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorRestoreInputs_Cuda(num_input_fields, qf_input_fields, op_input_fields, true, e_data, impl));
5620d0321e0SJeremy L Thompson 
5630d0321e0SJeremy L Thompson   // Restore output
564ca735530SJeremy L Thompson   CeedCallBackend(CeedVectorRestoreArray(*assembled, &assembled_array));
5650d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
5660d0321e0SJeremy L Thompson }
5670d0321e0SJeremy L Thompson 
5680d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
5690d0321e0SJeremy L Thompson // Assemble Linear QFunction
5700d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
5712b730f8bSJeremy L Thompson static int CeedOperatorLinearAssembleQFunction_Cuda(CeedOperator op, CeedVector *assembled, CeedElemRestriction *rstr, CeedRequest *request) {
5722b730f8bSJeremy L Thompson   return CeedOperatorLinearAssembleQFunctionCore_Cuda(op, true, assembled, rstr, request);
5730d0321e0SJeremy L Thompson }
5740d0321e0SJeremy L Thompson 
5750d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
5760d0321e0SJeremy L Thompson // Update Assembled Linear QFunction
5770d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
5782b730f8bSJeremy L Thompson static int CeedOperatorLinearAssembleQFunctionUpdate_Cuda(CeedOperator op, CeedVector assembled, CeedElemRestriction rstr, CeedRequest *request) {
5792b730f8bSJeremy L Thompson   return CeedOperatorLinearAssembleQFunctionCore_Cuda(op, false, &assembled, &rstr, request);
5800d0321e0SJeremy L Thompson }
5810d0321e0SJeremy L Thompson 
5820d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
583004e4986SSebastian Grimberg // Assemble Diagonal Setup
5840d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
585cbfe683aSSebastian Grimberg static inline int CeedOperatorAssembleDiagonalSetup_Cuda(CeedOperator op) {
5860d0321e0SJeremy L Thompson   Ceed                ceed;
587004e4986SSebastian Grimberg   CeedInt             num_input_fields, num_output_fields, num_eval_modes_in = 0, num_eval_modes_out = 0;
588cbfe683aSSebastian Grimberg   CeedInt             q_comp, num_nodes, num_qpts;
589004e4986SSebastian Grimberg   CeedEvalMode       *eval_modes_in = NULL, *eval_modes_out = NULL;
590ca735530SJeremy L Thompson   CeedBasis           basis_in = NULL, basis_out = NULL;
591ca735530SJeremy L Thompson   CeedQFunctionField *qf_fields;
5920d0321e0SJeremy L Thompson   CeedQFunction       qf;
593ca735530SJeremy L Thompson   CeedOperatorField  *op_fields;
594ca735530SJeremy L Thompson   CeedOperator_Cuda  *impl;
595ca735530SJeremy L Thompson 
596ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
5972b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
598ca735530SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetNumArgs(qf, &num_input_fields, &num_output_fields));
5990d0321e0SJeremy L Thompson 
6000d0321e0SJeremy L Thompson   // Determine active input basis
601ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetFields(op, NULL, &op_fields, NULL, NULL));
602ca735530SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_fields, NULL, NULL));
603ca735530SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
6040d0321e0SJeremy L Thompson     CeedVector vec;
605ca735530SJeremy L Thompson 
606ca735530SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec));
6070d0321e0SJeremy L Thompson     if (vec == CEED_VECTOR_ACTIVE) {
608004e4986SSebastian Grimberg       CeedBasis    basis;
609004e4986SSebastian Grimberg       CeedEvalMode eval_mode;
610ca735530SJeremy L Thompson 
611004e4986SSebastian Grimberg       CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis));
612004e4986SSebastian Grimberg       CeedCheck(!basis_in || basis_in == basis, ceed, CEED_ERROR_BACKEND,
613004e4986SSebastian Grimberg                 "Backend does not implement operator diagonal assembly with multiple active bases");
614004e4986SSebastian Grimberg       basis_in = basis;
615004e4986SSebastian Grimberg       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode));
616004e4986SSebastian Grimberg       CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_in, eval_mode, &q_comp));
617004e4986SSebastian Grimberg       if (eval_mode != CEED_EVAL_WEIGHT) {
618004e4986SSebastian Grimberg         // q_comp = 1 if CEED_EVAL_NONE, CEED_EVAL_WEIGHT caught by QF assembly
619004e4986SSebastian Grimberg         CeedCallBackend(CeedRealloc(num_eval_modes_in + q_comp, &eval_modes_in));
620004e4986SSebastian Grimberg         for (CeedInt d = 0; d < q_comp; d++) eval_modes_in[num_eval_modes_in + d] = eval_mode;
621004e4986SSebastian Grimberg         num_eval_modes_in += q_comp;
6220d0321e0SJeremy L Thompson       }
6230d0321e0SJeremy L Thompson     }
6240d0321e0SJeremy L Thompson   }
6250d0321e0SJeremy L Thompson 
6260d0321e0SJeremy L Thompson   // Determine active output basis
627ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetFields(op, NULL, NULL, NULL, &op_fields));
628ca735530SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qf_fields));
629ca735530SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
6300d0321e0SJeremy L Thompson     CeedVector vec;
631ca735530SJeremy L Thompson 
632ca735530SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec));
6330d0321e0SJeremy L Thompson     if (vec == CEED_VECTOR_ACTIVE) {
634004e4986SSebastian Grimberg       CeedBasis    basis;
635004e4986SSebastian Grimberg       CeedEvalMode eval_mode;
636ca735530SJeremy L Thompson 
637004e4986SSebastian Grimberg       CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis));
638004e4986SSebastian Grimberg       CeedCheck(!basis_out || basis_out == basis, ceed, CEED_ERROR_BACKEND,
639004e4986SSebastian Grimberg                 "Backend does not implement operator diagonal assembly with multiple active bases");
640004e4986SSebastian Grimberg       basis_out = basis;
641004e4986SSebastian Grimberg       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode));
642004e4986SSebastian Grimberg       CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_out, eval_mode, &q_comp));
643004e4986SSebastian Grimberg       if (eval_mode != CEED_EVAL_WEIGHT) {
644004e4986SSebastian Grimberg         // q_comp = 1 if CEED_EVAL_NONE, CEED_EVAL_WEIGHT caught by QF assembly
645004e4986SSebastian Grimberg         CeedCallBackend(CeedRealloc(num_eval_modes_out + q_comp, &eval_modes_out));
646004e4986SSebastian Grimberg         for (CeedInt d = 0; d < q_comp; d++) eval_modes_out[num_eval_modes_out + d] = eval_mode;
647004e4986SSebastian Grimberg         num_eval_modes_out += q_comp;
6480d0321e0SJeremy L Thompson       }
6490d0321e0SJeremy L Thompson     }
6500d0321e0SJeremy L Thompson   }
6510d0321e0SJeremy L Thompson 
6520d0321e0SJeremy L Thompson   // Operator data struct
6532b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &impl));
6542b730f8bSJeremy L Thompson   CeedCallBackend(CeedCalloc(1, &impl->diag));
6550d0321e0SJeremy L Thompson   CeedOperatorDiag_Cuda *diag = impl->diag;
656ca735530SJeremy L Thompson 
657cbfe683aSSebastian Grimberg   // Basis matrices
658004e4986SSebastian Grimberg   CeedCallBackend(CeedBasisGetNumNodes(basis_in, &num_nodes));
659004e4986SSebastian Grimberg   if (basis_in == CEED_BASIS_NONE) num_qpts = num_nodes;
660004e4986SSebastian Grimberg   else CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis_in, &num_qpts));
661004e4986SSebastian Grimberg   const CeedInt interp_bytes     = num_nodes * num_qpts * sizeof(CeedScalar);
662004e4986SSebastian Grimberg   const CeedInt eval_modes_bytes = sizeof(CeedEvalMode);
663004e4986SSebastian Grimberg   bool          has_eval_none    = false;
6640d0321e0SJeremy L Thompson 
6650d0321e0SJeremy L Thompson   // CEED_EVAL_NONE
666004e4986SSebastian Grimberg   for (CeedInt i = 0; i < num_eval_modes_in; i++) has_eval_none = has_eval_none || (eval_modes_in[i] == CEED_EVAL_NONE);
667004e4986SSebastian Grimberg   for (CeedInt i = 0; i < num_eval_modes_out; i++) has_eval_none = has_eval_none || (eval_modes_out[i] == CEED_EVAL_NONE);
668004e4986SSebastian Grimberg   if (has_eval_none) {
6690d0321e0SJeremy L Thompson     CeedScalar *identity = NULL;
670ca735530SJeremy L Thompson 
671004e4986SSebastian Grimberg     CeedCallBackend(CeedCalloc(num_nodes * num_qpts, &identity));
672ca735530SJeremy L Thompson     for (CeedInt i = 0; i < (num_nodes < num_qpts ? num_nodes : num_qpts); i++) identity[i * num_nodes + i] = 1.0;
673ca735530SJeremy L Thompson     CeedCallCuda(ceed, cudaMalloc((void **)&diag->d_identity, interp_bytes));
674ca735530SJeremy L Thompson     CeedCallCuda(ceed, cudaMemcpy(diag->d_identity, identity, interp_bytes, cudaMemcpyHostToDevice));
675004e4986SSebastian Grimberg     CeedCallBackend(CeedFree(&identity));
6760d0321e0SJeremy L Thompson   }
6770d0321e0SJeremy L Thompson 
678004e4986SSebastian Grimberg   // CEED_EVAL_INTERP, CEED_EVAL_GRAD, CEED_EVAL_DIV, and CEED_EVAL_CURL
679004e4986SSebastian Grimberg   for (CeedInt in = 0; in < 2; in++) {
680004e4986SSebastian Grimberg     CeedFESpace fespace;
681004e4986SSebastian Grimberg     CeedBasis   basis = in ? basis_in : basis_out;
6820d0321e0SJeremy L Thompson 
683004e4986SSebastian Grimberg     CeedCallBackend(CeedBasisGetFESpace(basis, &fespace));
684004e4986SSebastian Grimberg     switch (fespace) {
685004e4986SSebastian Grimberg       case CEED_FE_SPACE_H1: {
686004e4986SSebastian Grimberg         CeedInt           q_comp_interp, q_comp_grad;
687004e4986SSebastian Grimberg         const CeedScalar *interp, *grad;
688004e4986SSebastian Grimberg         CeedScalar       *d_interp, *d_grad;
6890d0321e0SJeremy L Thompson 
690004e4986SSebastian Grimberg         CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp));
691004e4986SSebastian Grimberg         CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_GRAD, &q_comp_grad));
6920d0321e0SJeremy L Thompson 
693004e4986SSebastian Grimberg         CeedCallBackend(CeedBasisGetInterp(basis, &interp));
694004e4986SSebastian Grimberg         CeedCallCuda(ceed, cudaMalloc((void **)&d_interp, interp_bytes * q_comp_interp));
695004e4986SSebastian Grimberg         CeedCallCuda(ceed, cudaMemcpy(d_interp, interp, interp_bytes * q_comp_interp, cudaMemcpyHostToDevice));
696004e4986SSebastian Grimberg         CeedCallBackend(CeedBasisGetGrad(basis, &grad));
697004e4986SSebastian Grimberg         CeedCallCuda(ceed, cudaMalloc((void **)&d_grad, interp_bytes * q_comp_grad));
698004e4986SSebastian Grimberg         CeedCallCuda(ceed, cudaMemcpy(d_grad, grad, interp_bytes * q_comp_grad, cudaMemcpyHostToDevice));
699004e4986SSebastian Grimberg         if (in) {
700004e4986SSebastian Grimberg           diag->d_interp_in = d_interp;
701004e4986SSebastian Grimberg           diag->d_grad_in   = d_grad;
702004e4986SSebastian Grimberg         } else {
703004e4986SSebastian Grimberg           diag->d_interp_out = d_interp;
704004e4986SSebastian Grimberg           diag->d_grad_out   = d_grad;
705004e4986SSebastian Grimberg         }
706004e4986SSebastian Grimberg       } break;
707004e4986SSebastian Grimberg       case CEED_FE_SPACE_HDIV: {
708004e4986SSebastian Grimberg         CeedInt           q_comp_interp, q_comp_div;
709004e4986SSebastian Grimberg         const CeedScalar *interp, *div;
710004e4986SSebastian Grimberg         CeedScalar       *d_interp, *d_div;
711004e4986SSebastian Grimberg 
712004e4986SSebastian Grimberg         CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp));
713004e4986SSebastian Grimberg         CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_DIV, &q_comp_div));
714004e4986SSebastian Grimberg 
715004e4986SSebastian Grimberg         CeedCallBackend(CeedBasisGetInterp(basis, &interp));
716004e4986SSebastian Grimberg         CeedCallCuda(ceed, cudaMalloc((void **)&d_interp, interp_bytes * q_comp_interp));
717004e4986SSebastian Grimberg         CeedCallCuda(ceed, cudaMemcpy(d_interp, interp, interp_bytes * q_comp_interp, cudaMemcpyHostToDevice));
718004e4986SSebastian Grimberg         CeedCallBackend(CeedBasisGetDiv(basis, &div));
719004e4986SSebastian Grimberg         CeedCallCuda(ceed, cudaMalloc((void **)&d_div, interp_bytes * q_comp_div));
720004e4986SSebastian Grimberg         CeedCallCuda(ceed, cudaMemcpy(d_div, div, interp_bytes * q_comp_div, cudaMemcpyHostToDevice));
721004e4986SSebastian Grimberg         if (in) {
722004e4986SSebastian Grimberg           diag->d_interp_in = d_interp;
723004e4986SSebastian Grimberg           diag->d_div_in    = d_div;
724004e4986SSebastian Grimberg         } else {
725004e4986SSebastian Grimberg           diag->d_interp_out = d_interp;
726004e4986SSebastian Grimberg           diag->d_div_out    = d_div;
727004e4986SSebastian Grimberg         }
728004e4986SSebastian Grimberg       } break;
729004e4986SSebastian Grimberg       case CEED_FE_SPACE_HCURL: {
730004e4986SSebastian Grimberg         CeedInt           q_comp_interp, q_comp_curl;
731004e4986SSebastian Grimberg         const CeedScalar *interp, *curl;
732004e4986SSebastian Grimberg         CeedScalar       *d_interp, *d_curl;
733004e4986SSebastian Grimberg 
734004e4986SSebastian Grimberg         CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp));
735004e4986SSebastian Grimberg         CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_CURL, &q_comp_curl));
736004e4986SSebastian Grimberg 
737004e4986SSebastian Grimberg         CeedCallBackend(CeedBasisGetInterp(basis, &interp));
738004e4986SSebastian Grimberg         CeedCallCuda(ceed, cudaMalloc((void **)&d_interp, interp_bytes * q_comp_interp));
739004e4986SSebastian Grimberg         CeedCallCuda(ceed, cudaMemcpy(d_interp, interp, interp_bytes * q_comp_interp, cudaMemcpyHostToDevice));
740004e4986SSebastian Grimberg         CeedCallBackend(CeedBasisGetCurl(basis, &curl));
741004e4986SSebastian Grimberg         CeedCallCuda(ceed, cudaMalloc((void **)&d_curl, interp_bytes * q_comp_curl));
742004e4986SSebastian Grimberg         CeedCallCuda(ceed, cudaMemcpy(d_curl, curl, interp_bytes * q_comp_curl, cudaMemcpyHostToDevice));
743004e4986SSebastian Grimberg         if (in) {
744004e4986SSebastian Grimberg           diag->d_interp_in = d_interp;
745004e4986SSebastian Grimberg           diag->d_curl_in   = d_curl;
746004e4986SSebastian Grimberg         } else {
747004e4986SSebastian Grimberg           diag->d_interp_out = d_interp;
748004e4986SSebastian Grimberg           diag->d_curl_out   = d_curl;
749004e4986SSebastian Grimberg         }
750004e4986SSebastian Grimberg       } break;
751004e4986SSebastian Grimberg     }
752004e4986SSebastian Grimberg   }
753004e4986SSebastian Grimberg 
754004e4986SSebastian Grimberg   // Arrays of eval_modes
755004e4986SSebastian Grimberg   CeedCallCuda(ceed, cudaMalloc((void **)&diag->d_eval_modes_in, num_eval_modes_in * eval_modes_bytes));
756004e4986SSebastian Grimberg   CeedCallCuda(ceed, cudaMemcpy(diag->d_eval_modes_in, eval_modes_in, num_eval_modes_in * eval_modes_bytes, cudaMemcpyHostToDevice));
757004e4986SSebastian Grimberg   CeedCallCuda(ceed, cudaMalloc((void **)&diag->d_eval_modes_out, num_eval_modes_out * eval_modes_bytes));
758004e4986SSebastian Grimberg   CeedCallCuda(ceed, cudaMemcpy(diag->d_eval_modes_out, eval_modes_out, num_eval_modes_out * eval_modes_bytes, cudaMemcpyHostToDevice));
759004e4986SSebastian Grimberg   CeedCallBackend(CeedFree(&eval_modes_in));
760004e4986SSebastian Grimberg   CeedCallBackend(CeedFree(&eval_modes_out));
7610d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
7620d0321e0SJeremy L Thompson }
7630d0321e0SJeremy L Thompson 
7640d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
765cbfe683aSSebastian Grimberg // Assemble Diagonal Setup (Compilation)
766cbfe683aSSebastian Grimberg //------------------------------------------------------------------------------
767cbfe683aSSebastian Grimberg static inline int CeedOperatorAssembleDiagonalSetupCompile_Cuda(CeedOperator op, CeedInt use_ceedsize_idx, const bool is_point_block) {
768cbfe683aSSebastian Grimberg   Ceed                ceed;
76922070f95SJeremy L Thompson   char               *diagonal_kernel_source;
77022070f95SJeremy L Thompson   const char         *diagonal_kernel_path;
771cbfe683aSSebastian Grimberg   CeedInt             num_input_fields, num_output_fields, num_eval_modes_in = 0, num_eval_modes_out = 0;
772cbfe683aSSebastian Grimberg   CeedInt             num_comp, q_comp, num_nodes, num_qpts;
773cbfe683aSSebastian Grimberg   CeedBasis           basis_in = NULL, basis_out = NULL;
774cbfe683aSSebastian Grimberg   CeedQFunctionField *qf_fields;
775cbfe683aSSebastian Grimberg   CeedQFunction       qf;
776cbfe683aSSebastian Grimberg   CeedOperatorField  *op_fields;
777cbfe683aSSebastian Grimberg   CeedOperator_Cuda  *impl;
778cbfe683aSSebastian Grimberg 
779cbfe683aSSebastian Grimberg   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
780cbfe683aSSebastian Grimberg   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
781cbfe683aSSebastian Grimberg   CeedCallBackend(CeedQFunctionGetNumArgs(qf, &num_input_fields, &num_output_fields));
782cbfe683aSSebastian Grimberg 
783cbfe683aSSebastian Grimberg   // Determine active input basis
784cbfe683aSSebastian Grimberg   CeedCallBackend(CeedOperatorGetFields(op, NULL, &op_fields, NULL, NULL));
785cbfe683aSSebastian Grimberg   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_fields, NULL, NULL));
786cbfe683aSSebastian Grimberg   for (CeedInt i = 0; i < num_input_fields; i++) {
787cbfe683aSSebastian Grimberg     CeedVector vec;
788cbfe683aSSebastian Grimberg 
789cbfe683aSSebastian Grimberg     CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec));
790cbfe683aSSebastian Grimberg     if (vec == CEED_VECTOR_ACTIVE) {
791cbfe683aSSebastian Grimberg       CeedEvalMode eval_mode;
792cbfe683aSSebastian Grimberg 
793cbfe683aSSebastian Grimberg       CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis_in));
794cbfe683aSSebastian Grimberg       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode));
795cbfe683aSSebastian Grimberg       CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_in, eval_mode, &q_comp));
796cbfe683aSSebastian Grimberg       if (eval_mode != CEED_EVAL_WEIGHT) {
797cbfe683aSSebastian Grimberg         num_eval_modes_in += q_comp;
798cbfe683aSSebastian Grimberg       }
799cbfe683aSSebastian Grimberg     }
800cbfe683aSSebastian Grimberg   }
801cbfe683aSSebastian Grimberg 
802cbfe683aSSebastian Grimberg   // Determine active output basis
803cbfe683aSSebastian Grimberg   CeedCallBackend(CeedOperatorGetFields(op, NULL, NULL, NULL, &op_fields));
804cbfe683aSSebastian Grimberg   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qf_fields));
805cbfe683aSSebastian Grimberg   for (CeedInt i = 0; i < num_output_fields; i++) {
806cbfe683aSSebastian Grimberg     CeedVector vec;
807cbfe683aSSebastian Grimberg 
808cbfe683aSSebastian Grimberg     CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec));
809cbfe683aSSebastian Grimberg     if (vec == CEED_VECTOR_ACTIVE) {
810cbfe683aSSebastian Grimberg       CeedEvalMode eval_mode;
811cbfe683aSSebastian Grimberg 
812cbfe683aSSebastian Grimberg       CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis_out));
813cbfe683aSSebastian Grimberg       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode));
814cbfe683aSSebastian Grimberg       CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_out, eval_mode, &q_comp));
815cbfe683aSSebastian Grimberg       if (eval_mode != CEED_EVAL_WEIGHT) {
816cbfe683aSSebastian Grimberg         num_eval_modes_out += q_comp;
817cbfe683aSSebastian Grimberg       }
818cbfe683aSSebastian Grimberg     }
819cbfe683aSSebastian Grimberg   }
820cbfe683aSSebastian Grimberg 
821cbfe683aSSebastian Grimberg   // Operator data struct
822cbfe683aSSebastian Grimberg   CeedCallBackend(CeedOperatorGetData(op, &impl));
823cbfe683aSSebastian Grimberg   CeedOperatorDiag_Cuda *diag = impl->diag;
824cbfe683aSSebastian Grimberg 
825cbfe683aSSebastian Grimberg   // Assemble kernel
826cbfe683aSSebastian Grimberg   CUmodule *module          = is_point_block ? &diag->module_point_block : &diag->module;
827cbfe683aSSebastian Grimberg   CeedInt   elems_per_block = 1;
828cbfe683aSSebastian Grimberg   CeedCallBackend(CeedBasisGetNumNodes(basis_in, &num_nodes));
829cbfe683aSSebastian Grimberg   CeedCallBackend(CeedBasisGetNumComponents(basis_in, &num_comp));
830cbfe683aSSebastian Grimberg   if (basis_in == CEED_BASIS_NONE) num_qpts = num_nodes;
831cbfe683aSSebastian Grimberg   else CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis_in, &num_qpts));
832cbfe683aSSebastian Grimberg   CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-operator-assemble-diagonal.h", &diagonal_kernel_path));
833cbfe683aSSebastian Grimberg   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Diagonal Assembly Kernel Source -----\n");
834cbfe683aSSebastian Grimberg   CeedCallBackend(CeedLoadSourceToBuffer(ceed, diagonal_kernel_path, &diagonal_kernel_source));
835cbfe683aSSebastian Grimberg   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Diagonal Assembly Source Complete! -----\n");
836cbfe683aSSebastian Grimberg   CeedCallCuda(ceed, CeedCompile_Cuda(ceed, diagonal_kernel_source, module, 8, "NUM_EVAL_MODES_IN", num_eval_modes_in, "NUM_EVAL_MODES_OUT",
837cbfe683aSSebastian Grimberg                                       num_eval_modes_out, "NUM_COMP", num_comp, "NUM_NODES", num_nodes, "NUM_QPTS", num_qpts, "USE_CEEDSIZE",
838cbfe683aSSebastian Grimberg                                       use_ceedsize_idx, "USE_POINT_BLOCK", is_point_block ? 1 : 0, "BLOCK_SIZE", num_nodes * elems_per_block));
839cbfe683aSSebastian Grimberg   CeedCallCuda(ceed, CeedGetKernel_Cuda(ceed, *module, "LinearDiagonal", is_point_block ? &diag->LinearPointBlock : &diag->LinearDiagonal));
840cbfe683aSSebastian Grimberg   CeedCallBackend(CeedFree(&diagonal_kernel_path));
841cbfe683aSSebastian Grimberg   CeedCallBackend(CeedFree(&diagonal_kernel_source));
842cbfe683aSSebastian Grimberg   return CEED_ERROR_SUCCESS;
843cbfe683aSSebastian Grimberg }
844cbfe683aSSebastian Grimberg 
845cbfe683aSSebastian Grimberg //------------------------------------------------------------------------------
846004e4986SSebastian Grimberg // Assemble Diagonal Core
8470d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
848ca735530SJeremy L Thompson static inline int CeedOperatorAssembleDiagonalCore_Cuda(CeedOperator op, CeedVector assembled, CeedRequest *request, const bool is_point_block) {
8490d0321e0SJeremy L Thompson   Ceed                ceed;
850cbfe683aSSebastian Grimberg   CeedInt             num_elem, num_nodes;
851ca735530SJeremy L Thompson   CeedScalar         *elem_diag_array;
852ca735530SJeremy L Thompson   const CeedScalar   *assembled_qf_array;
853004e4986SSebastian Grimberg   CeedVector          assembled_qf   = NULL, elem_diag;
854004e4986SSebastian Grimberg   CeedElemRestriction assembled_rstr = NULL, rstr_in, rstr_out, diag_rstr;
8550d0321e0SJeremy L Thompson   CeedOperator_Cuda  *impl;
856ca735530SJeremy L Thompson 
857ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
8582b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &impl));
8590d0321e0SJeremy L Thompson 
8600d0321e0SJeremy L Thompson   // Assemble QFunction
861004e4986SSebastian Grimberg   CeedCallBackend(CeedOperatorLinearAssembleQFunctionBuildOrUpdate(op, &assembled_qf, &assembled_rstr, request));
862004e4986SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionDestroy(&assembled_rstr));
863004e4986SSebastian Grimberg   CeedCallBackend(CeedVectorGetArrayRead(assembled_qf, CEED_MEM_DEVICE, &assembled_qf_array));
8640d0321e0SJeremy L Thompson 
865cbfe683aSSebastian Grimberg   // Setup
866cbfe683aSSebastian Grimberg   if (!impl->diag) CeedCallBackend(CeedOperatorAssembleDiagonalSetup_Cuda(op));
867cbfe683aSSebastian Grimberg   CeedOperatorDiag_Cuda *diag = impl->diag;
868cbfe683aSSebastian Grimberg 
869cbfe683aSSebastian Grimberg   assert(diag != NULL);
870cbfe683aSSebastian Grimberg 
871cbfe683aSSebastian Grimberg   // Assemble kernel if needed
872cbfe683aSSebastian Grimberg   if ((!is_point_block && !diag->LinearDiagonal) || (is_point_block && !diag->LinearPointBlock)) {
873cbfe683aSSebastian Grimberg     CeedSize assembled_length, assembled_qf_length;
874cbfe683aSSebastian Grimberg     CeedInt  use_ceedsize_idx = 0;
875f7c1b517Snbeams     CeedCallBackend(CeedVectorGetLength(assembled, &assembled_length));
876ca735530SJeremy L Thompson     CeedCallBackend(CeedVectorGetLength(assembled_qf, &assembled_qf_length));
877ca735530SJeremy L Thompson     if ((assembled_length > INT_MAX) || (assembled_qf_length > INT_MAX)) use_ceedsize_idx = 1;
878f7c1b517Snbeams 
879cbfe683aSSebastian Grimberg     CeedCallBackend(CeedOperatorAssembleDiagonalSetupCompile_Cuda(op, use_ceedsize_idx, is_point_block));
880cbfe683aSSebastian Grimberg   }
8810d0321e0SJeremy L Thompson 
882004e4986SSebastian Grimberg   // Restriction and diagonal vector
883004e4986SSebastian Grimberg   CeedCallBackend(CeedOperatorGetActiveElemRestrictions(op, &rstr_in, &rstr_out));
884004e4986SSebastian Grimberg   CeedCheck(rstr_in == rstr_out, ceed, CEED_ERROR_BACKEND,
885004e4986SSebastian Grimberg             "Cannot assemble operator diagonal with different input and output active element restrictions");
886004e4986SSebastian Grimberg   if (!is_point_block && !diag->diag_rstr) {
887004e4986SSebastian Grimberg     CeedCallBackend(CeedElemRestrictionCreateUnsignedCopy(rstr_out, &diag->diag_rstr));
888004e4986SSebastian Grimberg     CeedCallBackend(CeedElemRestrictionCreateVector(diag->diag_rstr, NULL, &diag->elem_diag));
889004e4986SSebastian Grimberg   } else if (is_point_block && !diag->point_block_diag_rstr) {
890004e4986SSebastian Grimberg     CeedCallBackend(CeedOperatorCreateActivePointBlockRestriction(rstr_out, &diag->point_block_diag_rstr));
891004e4986SSebastian Grimberg     CeedCallBackend(CeedElemRestrictionCreateVector(diag->point_block_diag_rstr, NULL, &diag->point_block_elem_diag));
8920d0321e0SJeremy L Thompson   }
893004e4986SSebastian Grimberg   diag_rstr = is_point_block ? diag->point_block_diag_rstr : diag->diag_rstr;
894004e4986SSebastian Grimberg   elem_diag = is_point_block ? diag->point_block_elem_diag : diag->elem_diag;
895ca735530SJeremy L Thompson   CeedCallBackend(CeedVectorSetValue(elem_diag, 0.0));
8960d0321e0SJeremy L Thompson 
89791db28b6SZach Atkins   // Only assemble diagonal if the basis has nodes, otherwise inputs are null pointers
898004e4986SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionGetElementSize(diag_rstr, &num_nodes));
899004e4986SSebastian Grimberg   if (num_nodes > 0) {
9000d0321e0SJeremy L Thompson     // Assemble element operator diagonals
901ca735530SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumElements(diag_rstr, &num_elem));
902004e4986SSebastian Grimberg     CeedCallBackend(CeedVectorGetArray(elem_diag, CEED_MEM_DEVICE, &elem_diag_array));
9030d0321e0SJeremy L Thompson 
9040d0321e0SJeremy L Thompson     // Compute the diagonal of B^T D B
905004e4986SSebastian Grimberg     CeedInt elems_per_block = 1;
906004e4986SSebastian Grimberg     CeedInt grid            = CeedDivUpInt(num_elem, elems_per_block);
907004e4986SSebastian Grimberg     void   *args[]          = {(void *)&num_elem,      &diag->d_identity,       &diag->d_interp_in,  &diag->d_grad_in, &diag->d_div_in,
908004e4986SSebastian Grimberg                                &diag->d_curl_in,       &diag->d_interp_out,     &diag->d_grad_out,   &diag->d_div_out, &diag->d_curl_out,
909004e4986SSebastian Grimberg                                &diag->d_eval_modes_in, &diag->d_eval_modes_out, &assembled_qf_array, &elem_diag_array};
910004e4986SSebastian Grimberg 
911ca735530SJeremy L Thompson     if (is_point_block) {
912004e4986SSebastian Grimberg       CeedCallBackend(CeedRunKernelDim_Cuda(ceed, diag->LinearPointBlock, grid, num_nodes, 1, elems_per_block, args));
9130d0321e0SJeremy L Thompson     } else {
914004e4986SSebastian Grimberg       CeedCallBackend(CeedRunKernelDim_Cuda(ceed, diag->LinearDiagonal, grid, num_nodes, 1, elems_per_block, args));
9150d0321e0SJeremy L Thompson     }
9160d0321e0SJeremy L Thompson 
9170d0321e0SJeremy L Thompson     // Restore arrays
918ca735530SJeremy L Thompson     CeedCallBackend(CeedVectorRestoreArray(elem_diag, &elem_diag_array));
919ca735530SJeremy L Thompson     CeedCallBackend(CeedVectorRestoreArrayRead(assembled_qf, &assembled_qf_array));
92091db28b6SZach Atkins   }
9210d0321e0SJeremy L Thompson 
9220d0321e0SJeremy L Thompson   // Assemble local operator diagonal
923ca735530SJeremy L Thompson   CeedCallBackend(CeedElemRestrictionApply(diag_rstr, CEED_TRANSPOSE, elem_diag, assembled, request));
9240d0321e0SJeremy L Thompson 
9250d0321e0SJeremy L Thompson   // Cleanup
926ca735530SJeremy L Thompson   CeedCallBackend(CeedVectorDestroy(&assembled_qf));
9270d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
9280d0321e0SJeremy L Thompson }
9290d0321e0SJeremy L Thompson 
9300d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
9310d0321e0SJeremy L Thompson // Assemble Linear Diagonal
9320d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
9332b730f8bSJeremy L Thompson static int CeedOperatorLinearAssembleAddDiagonal_Cuda(CeedOperator op, CeedVector assembled, CeedRequest *request) {
9342b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorAssembleDiagonalCore_Cuda(op, assembled, request, false));
9356aa95790SJeremy L Thompson   return CEED_ERROR_SUCCESS;
9360d0321e0SJeremy L Thompson }
9370d0321e0SJeremy L Thompson 
9380d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
9390d0321e0SJeremy L Thompson // Assemble Linear Point Block Diagonal
9400d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
9412b730f8bSJeremy L Thompson static int CeedOperatorLinearAssembleAddPointBlockDiagonal_Cuda(CeedOperator op, CeedVector assembled, CeedRequest *request) {
9422b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorAssembleDiagonalCore_Cuda(op, assembled, request, true));
9436aa95790SJeremy L Thompson   return CEED_ERROR_SUCCESS;
9440d0321e0SJeremy L Thompson }
9450d0321e0SJeremy L Thompson 
9460d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
947004e4986SSebastian Grimberg // Single Operator Assembly Setup
948cc132f9aSnbeams //------------------------------------------------------------------------------
949f7c1b517Snbeams static int CeedSingleOperatorAssembleSetup_Cuda(CeedOperator op, CeedInt use_ceedsize_idx) {
950cc132f9aSnbeams   Ceed                ceed;
951004e4986SSebastian Grimberg   Ceed_Cuda          *cuda_data;
95222070f95SJeremy L Thompson   char               *assembly_kernel_source;
95322070f95SJeremy L Thompson   const char         *assembly_kernel_path;
954004e4986SSebastian Grimberg   CeedInt             num_input_fields, num_output_fields, num_eval_modes_in = 0, num_eval_modes_out = 0;
9553b38d1dfSJeremy L Thompson   CeedInt             elem_size_in, num_qpts_in = 0, num_comp_in, elem_size_out, num_qpts_out, num_comp_out, q_comp;
956004e4986SSebastian Grimberg   CeedEvalMode       *eval_modes_in = NULL, *eval_modes_out = NULL;
957ca735530SJeremy L Thompson   CeedElemRestriction rstr_in = NULL, rstr_out = NULL;
958ca735530SJeremy L Thompson   CeedBasis           basis_in = NULL, basis_out = NULL;
959ca735530SJeremy L Thompson   CeedQFunctionField *qf_fields;
960ca735530SJeremy L Thompson   CeedQFunction       qf;
961ca735530SJeremy L Thompson   CeedOperatorField  *input_fields, *output_fields;
962cc132f9aSnbeams   CeedOperator_Cuda  *impl;
963ca735530SJeremy L Thompson 
964ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
9652b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &impl));
966cc132f9aSnbeams 
967cc132f9aSnbeams   // Get intput and output fields
9682b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &input_fields, &num_output_fields, &output_fields));
969cc132f9aSnbeams 
970cc132f9aSnbeams   // Determine active input basis eval mode
9712b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
9722b730f8bSJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_fields, NULL, NULL));
973cc132f9aSnbeams   for (CeedInt i = 0; i < num_input_fields; i++) {
974cc132f9aSnbeams     CeedVector vec;
975ca735530SJeremy L Thompson 
9762b730f8bSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetVector(input_fields[i], &vec));
977cc132f9aSnbeams     if (vec == CEED_VECTOR_ACTIVE) {
978004e4986SSebastian Grimberg       CeedBasis    basis;
979ca735530SJeremy L Thompson       CeedEvalMode eval_mode;
980ca735530SJeremy L Thompson 
981004e4986SSebastian Grimberg       CeedCallBackend(CeedOperatorFieldGetBasis(input_fields[i], &basis));
982004e4986SSebastian Grimberg       CeedCheck(!basis_in || basis_in == basis, ceed, CEED_ERROR_BACKEND, "Backend does not implement operator assembly with multiple active bases");
983004e4986SSebastian Grimberg       basis_in = basis;
9842b730f8bSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetElemRestriction(input_fields[i], &rstr_in));
985004e4986SSebastian Grimberg       CeedCallBackend(CeedElemRestrictionGetElementSize(rstr_in, &elem_size_in));
986004e4986SSebastian Grimberg       if (basis_in == CEED_BASIS_NONE) num_qpts_in = elem_size_in;
987004e4986SSebastian Grimberg       else CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis_in, &num_qpts_in));
9882b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode));
989004e4986SSebastian Grimberg       CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_in, eval_mode, &q_comp));
990004e4986SSebastian Grimberg       if (eval_mode != CEED_EVAL_WEIGHT) {
991004e4986SSebastian Grimberg         // q_comp = 1 if CEED_EVAL_NONE, CEED_EVAL_WEIGHT caught by QF Assembly
992004e4986SSebastian Grimberg         CeedCallBackend(CeedRealloc(num_eval_modes_in + q_comp, &eval_modes_in));
993004e4986SSebastian Grimberg         for (CeedInt d = 0; d < q_comp; d++) {
994004e4986SSebastian Grimberg           eval_modes_in[num_eval_modes_in + d] = eval_mode;
995cc132f9aSnbeams         }
996004e4986SSebastian Grimberg         num_eval_modes_in += q_comp;
997cc132f9aSnbeams       }
998cc132f9aSnbeams     }
999cc132f9aSnbeams   }
1000cc132f9aSnbeams 
1001cc132f9aSnbeams   // Determine active output basis; basis_out and rstr_out only used if same as input, TODO
10022b730f8bSJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qf_fields));
1003cc132f9aSnbeams   for (CeedInt i = 0; i < num_output_fields; i++) {
1004cc132f9aSnbeams     CeedVector vec;
1005ca735530SJeremy L Thompson 
10062b730f8bSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetVector(output_fields[i], &vec));
1007cc132f9aSnbeams     if (vec == CEED_VECTOR_ACTIVE) {
1008004e4986SSebastian Grimberg       CeedBasis    basis;
1009ca735530SJeremy L Thompson       CeedEvalMode eval_mode;
1010ca735530SJeremy L Thompson 
1011004e4986SSebastian Grimberg       CeedCallBackend(CeedOperatorFieldGetBasis(output_fields[i], &basis));
1012004e4986SSebastian Grimberg       CeedCheck(!basis_out || basis_out == basis, ceed, CEED_ERROR_BACKEND,
1013004e4986SSebastian Grimberg                 "Backend does not implement operator assembly with multiple active bases");
1014004e4986SSebastian Grimberg       basis_out = basis;
10152b730f8bSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetElemRestriction(output_fields[i], &rstr_out));
1016004e4986SSebastian Grimberg       CeedCallBackend(CeedElemRestrictionGetElementSize(rstr_out, &elem_size_out));
1017004e4986SSebastian Grimberg       if (basis_out == CEED_BASIS_NONE) num_qpts_out = elem_size_out;
1018004e4986SSebastian Grimberg       else CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis_out, &num_qpts_out));
1019004e4986SSebastian Grimberg       CeedCheck(num_qpts_in == num_qpts_out, ceed, CEED_ERROR_UNSUPPORTED,
1020004e4986SSebastian Grimberg                 "Active input and output bases must have the same number of quadrature points");
10212b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode));
1022004e4986SSebastian Grimberg       CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_out, eval_mode, &q_comp));
1023004e4986SSebastian Grimberg       if (eval_mode != CEED_EVAL_WEIGHT) {
1024004e4986SSebastian Grimberg         // q_comp = 1 if CEED_EVAL_NONE, CEED_EVAL_WEIGHT caught by QF Assembly
1025004e4986SSebastian Grimberg         CeedCallBackend(CeedRealloc(num_eval_modes_out + q_comp, &eval_modes_out));
1026004e4986SSebastian Grimberg         for (CeedInt d = 0; d < q_comp; d++) {
1027004e4986SSebastian Grimberg           eval_modes_out[num_eval_modes_out + d] = eval_mode;
1028004e4986SSebastian Grimberg         }
1029004e4986SSebastian Grimberg         num_eval_modes_out += q_comp;
1030cc132f9aSnbeams       }
1031cc132f9aSnbeams     }
1032cc132f9aSnbeams   }
1033004e4986SSebastian Grimberg   CeedCheck(num_eval_modes_in > 0 && num_eval_modes_out > 0, ceed, CEED_ERROR_UNSUPPORTED, "Cannot assemble operator without inputs/outputs");
1034cc132f9aSnbeams 
10352b730f8bSJeremy L Thompson   CeedCallBackend(CeedCalloc(1, &impl->asmb));
1036cc132f9aSnbeams   CeedOperatorAssemble_Cuda *asmb = impl->asmb;
1037004e4986SSebastian Grimberg   asmb->elems_per_block           = 1;
1038004e4986SSebastian Grimberg   asmb->block_size_x              = elem_size_in;
1039004e4986SSebastian Grimberg   asmb->block_size_y              = elem_size_out;
1040ca735530SJeremy L Thompson 
10412b730f8bSJeremy L Thompson   CeedCallBackend(CeedGetData(ceed, &cuda_data));
1042004e4986SSebastian Grimberg   bool fallback = asmb->block_size_x * asmb->block_size_y * asmb->elems_per_block > cuda_data->device_prop.maxThreadsPerBlock;
1043004e4986SSebastian Grimberg 
1044004e4986SSebastian Grimberg   if (fallback) {
1045004e4986SSebastian Grimberg     // Use fallback kernel with 1D threadblock
1046004e4986SSebastian Grimberg     asmb->block_size_y = 1;
1047004e4986SSebastian Grimberg   }
1048004e4986SSebastian Grimberg 
1049004e4986SSebastian Grimberg   // Compile kernels
1050004e4986SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionGetNumComponents(rstr_in, &num_comp_in));
1051004e4986SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionGetNumComponents(rstr_out, &num_comp_out));
10522b730f8bSJeremy L Thompson   CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-operator-assemble.h", &assembly_kernel_path));
105323d4529eSJeremy L Thompson   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Assembly Kernel Source -----\n");
10542b730f8bSJeremy L Thompson   CeedCallBackend(CeedLoadSourceToBuffer(ceed, assembly_kernel_path, &assembly_kernel_source));
105523d4529eSJeremy L Thompson   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Assembly Source Complete! -----\n");
1056004e4986SSebastian Grimberg   CeedCallBackend(CeedCompile_Cuda(ceed, assembly_kernel_source, &asmb->module, 10, "NUM_EVAL_MODES_IN", num_eval_modes_in, "NUM_EVAL_MODES_OUT",
1057004e4986SSebastian Grimberg                                    num_eval_modes_out, "NUM_COMP_IN", num_comp_in, "NUM_COMP_OUT", num_comp_out, "NUM_NODES_IN", elem_size_in,
1058004e4986SSebastian Grimberg                                    "NUM_NODES_OUT", elem_size_out, "NUM_QPTS", num_qpts_in, "BLOCK_SIZE",
1059cbfe683aSSebastian Grimberg                                    asmb->block_size_x * asmb->block_size_y * asmb->elems_per_block, "BLOCK_SIZE_Y", asmb->block_size_y,
1060cbfe683aSSebastian Grimberg                                    "USE_CEEDSIZE", use_ceedsize_idx));
1061004e4986SSebastian Grimberg   CeedCallBackend(CeedGetKernel_Cuda(ceed, asmb->module, "LinearAssemble", &asmb->LinearAssemble));
10622b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&assembly_kernel_path));
10632b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&assembly_kernel_source));
1064cc132f9aSnbeams 
1065004e4986SSebastian Grimberg   // Load into B_in, in order that they will be used in eval_modes_in
1066004e4986SSebastian Grimberg   {
1067004e4986SSebastian Grimberg     const CeedInt in_bytes           = elem_size_in * num_qpts_in * num_eval_modes_in * sizeof(CeedScalar);
1068004e4986SSebastian Grimberg     CeedInt       d_in               = 0;
1069004e4986SSebastian Grimberg     CeedEvalMode  eval_modes_in_prev = CEED_EVAL_NONE;
1070004e4986SSebastian Grimberg     bool          has_eval_none      = false;
1071004e4986SSebastian Grimberg     CeedScalar   *identity           = NULL;
1072ca735530SJeremy L Thompson 
1073004e4986SSebastian Grimberg     for (CeedInt i = 0; i < num_eval_modes_in; i++) {
1074004e4986SSebastian Grimberg       has_eval_none = has_eval_none || (eval_modes_in[i] == CEED_EVAL_NONE);
1075004e4986SSebastian Grimberg     }
1076004e4986SSebastian Grimberg     if (has_eval_none) {
1077004e4986SSebastian Grimberg       CeedCallBackend(CeedCalloc(elem_size_in * num_qpts_in, &identity));
1078004e4986SSebastian 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;
1079004e4986SSebastian Grimberg     }
1080cc132f9aSnbeams 
1081004e4986SSebastian Grimberg     CeedCallCuda(ceed, cudaMalloc((void **)&asmb->d_B_in, in_bytes));
1082004e4986SSebastian Grimberg     for (CeedInt i = 0; i < num_eval_modes_in; i++) {
1083004e4986SSebastian Grimberg       const CeedScalar *h_B_in;
1084ca735530SJeremy L Thompson 
1085004e4986SSebastian Grimberg       CeedCallBackend(CeedOperatorGetBasisPointer(basis_in, eval_modes_in[i], identity, &h_B_in));
1086004e4986SSebastian Grimberg       CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_in, eval_modes_in[i], &q_comp));
1087004e4986SSebastian Grimberg       if (q_comp > 1) {
1088004e4986SSebastian Grimberg         if (i == 0 || eval_modes_in[i] != eval_modes_in_prev) d_in = 0;
1089004e4986SSebastian Grimberg         else h_B_in = &h_B_in[(++d_in) * elem_size_in * num_qpts_in];
1090004e4986SSebastian Grimberg       }
1091004e4986SSebastian Grimberg       eval_modes_in_prev = eval_modes_in[i];
1092ca735530SJeremy L Thompson 
1093004e4986SSebastian 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),
1094004e4986SSebastian Grimberg                                     cudaMemcpyHostToDevice));
1095004e4986SSebastian Grimberg     }
1096004e4986SSebastian Grimberg 
1097004e4986SSebastian Grimberg     if (identity) {
1098004e4986SSebastian Grimberg       CeedCallBackend(CeedFree(&identity));
1099cc132f9aSnbeams     }
1100cc132f9aSnbeams   }
1101cc132f9aSnbeams 
1102004e4986SSebastian Grimberg   // Load into B_out, in order that they will be used in eval_modes_out
1103004e4986SSebastian Grimberg   {
1104004e4986SSebastian Grimberg     const CeedInt out_bytes           = elem_size_out * num_qpts_out * num_eval_modes_out * sizeof(CeedScalar);
1105004e4986SSebastian Grimberg     CeedInt       d_out               = 0;
1106004e4986SSebastian Grimberg     CeedEvalMode  eval_modes_out_prev = CEED_EVAL_NONE;
1107004e4986SSebastian Grimberg     bool          has_eval_none       = false;
1108004e4986SSebastian Grimberg     CeedScalar   *identity            = NULL;
1109ca735530SJeremy L Thompson 
1110004e4986SSebastian Grimberg     for (CeedInt i = 0; i < num_eval_modes_out; i++) {
1111004e4986SSebastian Grimberg       has_eval_none = has_eval_none || (eval_modes_out[i] == CEED_EVAL_NONE);
1112004e4986SSebastian Grimberg     }
1113004e4986SSebastian Grimberg     if (has_eval_none) {
1114004e4986SSebastian Grimberg       CeedCallBackend(CeedCalloc(elem_size_out * num_qpts_out, &identity));
1115004e4986SSebastian 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;
1116cc132f9aSnbeams     }
1117cc132f9aSnbeams 
1118004e4986SSebastian Grimberg     CeedCallCuda(ceed, cudaMalloc((void **)&asmb->d_B_out, out_bytes));
1119004e4986SSebastian Grimberg     for (CeedInt i = 0; i < num_eval_modes_out; i++) {
1120004e4986SSebastian Grimberg       const CeedScalar *h_B_out;
1121ca735530SJeremy L Thompson 
1122004e4986SSebastian Grimberg       CeedCallBackend(CeedOperatorGetBasisPointer(basis_out, eval_modes_out[i], identity, &h_B_out));
1123004e4986SSebastian Grimberg       CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_out, eval_modes_out[i], &q_comp));
1124004e4986SSebastian Grimberg       if (q_comp > 1) {
1125004e4986SSebastian Grimberg         if (i == 0 || eval_modes_out[i] != eval_modes_out_prev) d_out = 0;
1126004e4986SSebastian Grimberg         else h_B_out = &h_B_out[(++d_out) * elem_size_out * num_qpts_out];
1127004e4986SSebastian Grimberg       }
1128004e4986SSebastian Grimberg       eval_modes_out_prev = eval_modes_out[i];
1129ca735530SJeremy L Thompson 
1130004e4986SSebastian 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),
1131004e4986SSebastian Grimberg                                     cudaMemcpyHostToDevice));
1132004e4986SSebastian Grimberg     }
1133004e4986SSebastian Grimberg 
1134004e4986SSebastian Grimberg     if (identity) {
1135004e4986SSebastian Grimberg       CeedCallBackend(CeedFree(&identity));
1136cc132f9aSnbeams     }
1137cc132f9aSnbeams   }
1138cc132f9aSnbeams   return CEED_ERROR_SUCCESS;
1139cc132f9aSnbeams }
1140cc132f9aSnbeams 
1141cc132f9aSnbeams //------------------------------------------------------------------------------
1142cefa2673SJeremy L Thompson // Assemble matrix data for COO matrix of assembled operator.
1143cefa2673SJeremy L Thompson // The sparsity pattern is set by CeedOperatorLinearAssembleSymbolic.
1144cefa2673SJeremy L Thompson //
1145ea61e9acSJeremy L Thompson // Note that this (and other assembly routines) currently assume only one active input restriction/basis per operator (could have multiple basis eval
1146cefa2673SJeremy L Thompson // modes).
1147cefa2673SJeremy L Thompson // TODO: allow multiple active input restrictions/basis objects
1148cc132f9aSnbeams //------------------------------------------------------------------------------
11492b730f8bSJeremy L Thompson static int CeedSingleOperatorAssemble_Cuda(CeedOperator op, CeedInt offset, CeedVector values) {
1150cc132f9aSnbeams   Ceed                ceed;
1151ca735530SJeremy L Thompson   CeedSize            values_length = 0, assembled_qf_length = 0;
1152004e4986SSebastian Grimberg   CeedInt             use_ceedsize_idx = 0, num_elem_in, num_elem_out, elem_size_in, elem_size_out;
1153ca735530SJeremy L Thompson   CeedScalar         *values_array;
1154004e4986SSebastian Grimberg   const CeedScalar   *assembled_qf_array;
1155ca735530SJeremy L Thompson   CeedVector          assembled_qf   = NULL;
1156004e4986SSebastian Grimberg   CeedElemRestriction assembled_rstr = NULL, rstr_in, rstr_out;
1157004e4986SSebastian Grimberg   CeedRestrictionType rstr_type_in, rstr_type_out;
1158004e4986SSebastian Grimberg   const bool         *orients_in = NULL, *orients_out = NULL;
1159004e4986SSebastian Grimberg   const CeedInt8     *curl_orients_in = NULL, *curl_orients_out = NULL;
1160cc132f9aSnbeams   CeedOperator_Cuda  *impl;
1161ca735530SJeremy L Thompson 
1162ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
11632b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &impl));
1164cc132f9aSnbeams 
1165cc132f9aSnbeams   // Assemble QFunction
1166004e4986SSebastian Grimberg   CeedCallBackend(CeedOperatorLinearAssembleQFunctionBuildOrUpdate(op, &assembled_qf, &assembled_rstr, CEED_REQUEST_IMMEDIATE));
1167004e4986SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionDestroy(&assembled_rstr));
1168004e4986SSebastian Grimberg   CeedCallBackend(CeedVectorGetArrayRead(assembled_qf, CEED_MEM_DEVICE, &assembled_qf_array));
1169cc132f9aSnbeams 
1170f7c1b517Snbeams   CeedCallBackend(CeedVectorGetLength(values, &values_length));
1171f7c1b517Snbeams   CeedCallBackend(CeedVectorGetLength(assembled_qf, &assembled_qf_length));
1172f7c1b517Snbeams   if ((values_length > INT_MAX) || (assembled_qf_length > INT_MAX)) use_ceedsize_idx = 1;
1173004e4986SSebastian Grimberg 
1174f7c1b517Snbeams   // Setup
1175004e4986SSebastian Grimberg   if (!impl->asmb) CeedCallBackend(CeedSingleOperatorAssembleSetup_Cuda(op, use_ceedsize_idx));
1176004e4986SSebastian Grimberg   CeedOperatorAssemble_Cuda *asmb = impl->asmb;
1177004e4986SSebastian Grimberg 
1178004e4986SSebastian Grimberg   assert(asmb != NULL);
1179004e4986SSebastian Grimberg 
1180004e4986SSebastian Grimberg   // Assemble element operator
1181004e4986SSebastian Grimberg   CeedCallBackend(CeedVectorGetArray(values, CEED_MEM_DEVICE, &values_array));
1182004e4986SSebastian Grimberg   values_array += offset;
1183004e4986SSebastian Grimberg 
1184004e4986SSebastian Grimberg   CeedCallBackend(CeedOperatorGetActiveElemRestrictions(op, &rstr_in, &rstr_out));
1185004e4986SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionGetNumElements(rstr_in, &num_elem_in));
1186004e4986SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionGetElementSize(rstr_in, &elem_size_in));
1187004e4986SSebastian Grimberg 
1188004e4986SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionGetType(rstr_in, &rstr_type_in));
1189004e4986SSebastian Grimberg   if (rstr_type_in == CEED_RESTRICTION_ORIENTED) {
1190004e4986SSebastian Grimberg     CeedCallBackend(CeedElemRestrictionGetOrientations(rstr_in, CEED_MEM_DEVICE, &orients_in));
1191004e4986SSebastian Grimberg   } else if (rstr_type_in == CEED_RESTRICTION_CURL_ORIENTED) {
1192004e4986SSebastian Grimberg     CeedCallBackend(CeedElemRestrictionGetCurlOrientations(rstr_in, CEED_MEM_DEVICE, &curl_orients_in));
1193004e4986SSebastian Grimberg   }
1194004e4986SSebastian Grimberg 
1195004e4986SSebastian Grimberg   if (rstr_in != rstr_out) {
1196004e4986SSebastian Grimberg     CeedCallBackend(CeedElemRestrictionGetNumElements(rstr_out, &num_elem_out));
1197004e4986SSebastian Grimberg     CeedCheck(num_elem_in == num_elem_out, ceed, CEED_ERROR_UNSUPPORTED,
1198004e4986SSebastian Grimberg               "Active input and output operator restrictions must have the same number of elements");
1199004e4986SSebastian Grimberg     CeedCallBackend(CeedElemRestrictionGetElementSize(rstr_out, &elem_size_out));
1200004e4986SSebastian Grimberg 
1201004e4986SSebastian Grimberg     CeedCallBackend(CeedElemRestrictionGetType(rstr_out, &rstr_type_out));
1202004e4986SSebastian Grimberg     if (rstr_type_out == CEED_RESTRICTION_ORIENTED) {
1203004e4986SSebastian Grimberg       CeedCallBackend(CeedElemRestrictionGetOrientations(rstr_out, CEED_MEM_DEVICE, &orients_out));
1204004e4986SSebastian Grimberg     } else if (rstr_type_out == CEED_RESTRICTION_CURL_ORIENTED) {
1205004e4986SSebastian Grimberg       CeedCallBackend(CeedElemRestrictionGetCurlOrientations(rstr_out, CEED_MEM_DEVICE, &curl_orients_out));
1206004e4986SSebastian Grimberg     }
1207004e4986SSebastian Grimberg   } else {
1208004e4986SSebastian Grimberg     elem_size_out    = elem_size_in;
1209004e4986SSebastian Grimberg     orients_out      = orients_in;
1210004e4986SSebastian Grimberg     curl_orients_out = curl_orients_in;
1211f7c1b517Snbeams   }
1212f7c1b517Snbeams 
1213cc132f9aSnbeams   // Compute B^T D B
1214004e4986SSebastian Grimberg   CeedInt shared_mem =
1215004e4986SSebastian 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)) *
1216004e4986SSebastian Grimberg       sizeof(CeedScalar);
1217004e4986SSebastian Grimberg   CeedInt grid   = CeedDivUpInt(num_elem_in, asmb->elems_per_block);
1218004e4986SSebastian Grimberg   void   *args[] = {(void *)&num_elem_in, &asmb->d_B_in,     &asmb->d_B_out,      &orients_in,  &curl_orients_in,
1219004e4986SSebastian Grimberg                     &orients_out,         &curl_orients_out, &assembled_qf_array, &values_array};
1220ca735530SJeremy L Thompson 
12212b730f8bSJeremy L Thompson   CeedCallBackend(
1222004e4986SSebastian Grimberg       CeedRunKernelDimShared_Cuda(ceed, asmb->LinearAssemble, grid, asmb->block_size_x, asmb->block_size_y, asmb->elems_per_block, shared_mem, args));
1223cc132f9aSnbeams 
1224cc132f9aSnbeams   // Restore arrays
12252b730f8bSJeremy L Thompson   CeedCallBackend(CeedVectorRestoreArray(values, &values_array));
1226004e4986SSebastian Grimberg   CeedCallBackend(CeedVectorRestoreArrayRead(assembled_qf, &assembled_qf_array));
1227cc132f9aSnbeams 
1228cc132f9aSnbeams   // Cleanup
12292b730f8bSJeremy L Thompson   CeedCallBackend(CeedVectorDestroy(&assembled_qf));
1230004e4986SSebastian Grimberg   if (rstr_type_in == CEED_RESTRICTION_ORIENTED) {
1231004e4986SSebastian Grimberg     CeedCallBackend(CeedElemRestrictionRestoreOrientations(rstr_in, &orients_in));
1232004e4986SSebastian Grimberg   } else if (rstr_type_in == CEED_RESTRICTION_CURL_ORIENTED) {
1233004e4986SSebastian Grimberg     CeedCallBackend(CeedElemRestrictionRestoreCurlOrientations(rstr_in, &curl_orients_in));
1234004e4986SSebastian Grimberg   }
1235004e4986SSebastian Grimberg   if (rstr_in != rstr_out) {
1236004e4986SSebastian Grimberg     if (rstr_type_out == CEED_RESTRICTION_ORIENTED) {
1237004e4986SSebastian Grimberg       CeedCallBackend(CeedElemRestrictionRestoreOrientations(rstr_out, &orients_out));
1238004e4986SSebastian Grimberg     } else if (rstr_type_out == CEED_RESTRICTION_CURL_ORIENTED) {
1239004e4986SSebastian Grimberg       CeedCallBackend(CeedElemRestrictionRestoreCurlOrientations(rstr_out, &curl_orients_out));
1240004e4986SSebastian Grimberg     }
1241004e4986SSebastian Grimberg   }
1242cc132f9aSnbeams   return CEED_ERROR_SUCCESS;
1243cc132f9aSnbeams }
1244cc132f9aSnbeams 
1245cc132f9aSnbeams //------------------------------------------------------------------------------
12460d0321e0SJeremy L Thompson // Create operator
12470d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
12480d0321e0SJeremy L Thompson int CeedOperatorCreate_Cuda(CeedOperator op) {
12490d0321e0SJeremy L Thompson   Ceed               ceed;
12500d0321e0SJeremy L Thompson   CeedOperator_Cuda *impl;
12510d0321e0SJeremy L Thompson 
1252ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
12532b730f8bSJeremy L Thompson   CeedCallBackend(CeedCalloc(1, &impl));
12542b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorSetData(op, impl));
12550d0321e0SJeremy L Thompson 
12562b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleQFunction", CeedOperatorLinearAssembleQFunction_Cuda));
12572b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleQFunctionUpdate", CeedOperatorLinearAssembleQFunctionUpdate_Cuda));
12582b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleAddDiagonal", CeedOperatorLinearAssembleAddDiagonal_Cuda));
12592b730f8bSJeremy L Thompson   CeedCallBackend(
12602b730f8bSJeremy L Thompson       CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleAddPointBlockDiagonal", CeedOperatorLinearAssembleAddPointBlockDiagonal_Cuda));
12612b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleSingle", CeedSingleOperatorAssemble_Cuda));
12622b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "ApplyAdd", CeedOperatorApplyAdd_Cuda));
12632b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "Destroy", CeedOperatorDestroy_Cuda));
12640d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
12650d0321e0SJeremy L Thompson }
12660d0321e0SJeremy L Thompson 
12670d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
1268