xref: /libCEED/rust/libceed-sys/c-src/backends/cuda-ref/ceed-cuda-ref-operator.c (revision 8b0f73481df5f9f5ba2e5a396f0652e2438a7896)
15aed82e4SJeremy L Thompson // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors.
23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
30d0321e0SJeremy L Thompson //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
50d0321e0SJeremy L Thompson //
63d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
70d0321e0SJeremy L Thompson 
849aac155SJeremy L Thompson #include <ceed.h>
92b730f8bSJeremy L Thompson #include <ceed/backend.h>
102b730f8bSJeremy L Thompson #include <ceed/jit-tools.h>
11c85e8640SSebastian Grimberg #include <assert.h>
120d0321e0SJeremy L Thompson #include <cuda.h>
130d0321e0SJeremy L Thompson #include <cuda_runtime.h>
140d0321e0SJeremy L Thompson #include <stdbool.h>
150d0321e0SJeremy L Thompson #include <string.h>
162b730f8bSJeremy L Thompson 
1749aac155SJeremy L Thompson #include "../cuda/ceed-cuda-common.h"
180d0321e0SJeremy L Thompson #include "../cuda/ceed-cuda-compile.h"
192b730f8bSJeremy L Thompson #include "ceed-cuda-ref.h"
200d0321e0SJeremy L Thompson 
210d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
220d0321e0SJeremy L Thompson // Destroy operator
230d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
240d0321e0SJeremy L Thompson static int CeedOperatorDestroy_Cuda(CeedOperator op) {
250d0321e0SJeremy L Thompson   CeedOperator_Cuda *impl;
26ca735530SJeremy L Thompson 
272b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &impl));
280d0321e0SJeremy L Thompson 
290d0321e0SJeremy L Thompson   // Apply data
30111870feSJeremy L Thompson   CeedCallBackend(CeedFree(&impl->num_points));
313aab95c0SJeremy L Thompson   CeedCallBackend(CeedFree(&impl->skip_rstr_in));
32f8a0df59SJeremy L Thompson   CeedCallBackend(CeedFree(&impl->skip_rstr_out));
33f8a0df59SJeremy L Thompson   CeedCallBackend(CeedFree(&impl->apply_add_basis_out));
34034f99fdSJeremy L Thompson   CeedCallBackend(CeedFree(&impl->input_field_order));
35034f99fdSJeremy L Thompson   CeedCallBackend(CeedFree(&impl->output_field_order));
36c1222711SJeremy L Thompson   CeedCallBackend(CeedFree(&impl->input_states));
370d0321e0SJeremy L Thompson 
38ca735530SJeremy L Thompson   for (CeedInt i = 0; i < impl->num_inputs; i++) {
39034f99fdSJeremy L Thompson     CeedCallBackend(CeedVectorDestroy(&impl->e_vecs_in[i]));
40ca735530SJeremy L Thompson     CeedCallBackend(CeedVectorDestroy(&impl->q_vecs_in[i]));
410d0321e0SJeremy L Thompson   }
42034f99fdSJeremy L Thompson   CeedCallBackend(CeedFree(&impl->e_vecs_in));
43ca735530SJeremy L Thompson   CeedCallBackend(CeedFree(&impl->q_vecs_in));
440d0321e0SJeremy L Thompson 
45ca735530SJeremy L Thompson   for (CeedInt i = 0; i < impl->num_outputs; i++) {
46034f99fdSJeremy L Thompson     CeedCallBackend(CeedVectorDestroy(&impl->e_vecs_out[i]));
47ca735530SJeremy L Thompson     CeedCallBackend(CeedVectorDestroy(&impl->q_vecs_out[i]));
480d0321e0SJeremy L Thompson   }
49034f99fdSJeremy L Thompson   CeedCallBackend(CeedFree(&impl->e_vecs_out));
50ca735530SJeremy L Thompson   CeedCallBackend(CeedFree(&impl->q_vecs_out));
51756ca9e9SJeremy L Thompson   CeedCallBackend(CeedVectorDestroy(&impl->point_coords_elem));
520d0321e0SJeremy L Thompson 
530d0321e0SJeremy L Thompson   // QFunction assembly data
54ca735530SJeremy L Thompson   for (CeedInt i = 0; i < impl->num_active_in; i++) {
55ca735530SJeremy L Thompson     CeedCallBackend(CeedVectorDestroy(&impl->qf_active_in[i]));
560d0321e0SJeremy L Thompson   }
57ca735530SJeremy L Thompson   CeedCallBackend(CeedFree(&impl->qf_active_in));
580d0321e0SJeremy L Thompson 
590d0321e0SJeremy L Thompson   // Diag data
600d0321e0SJeremy L Thompson   if (impl->diag) {
610d0321e0SJeremy L Thompson     Ceed ceed;
62ca735530SJeremy L Thompson 
632b730f8bSJeremy L Thompson     CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
64cbfe683aSSebastian Grimberg     if (impl->diag->module) {
652b730f8bSJeremy L Thompson       CeedCallCuda(ceed, cuModuleUnload(impl->diag->module));
66cbfe683aSSebastian Grimberg     }
67cbfe683aSSebastian Grimberg     if (impl->diag->module_point_block) {
68cbfe683aSSebastian Grimberg       CeedCallCuda(ceed, cuModuleUnload(impl->diag->module_point_block));
69cbfe683aSSebastian Grimberg     }
70004e4986SSebastian Grimberg     CeedCallCuda(ceed, cudaFree(impl->diag->d_eval_modes_in));
71004e4986SSebastian Grimberg     CeedCallCuda(ceed, cudaFree(impl->diag->d_eval_modes_out));
722b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaFree(impl->diag->d_identity));
73ca735530SJeremy L Thompson     CeedCallCuda(ceed, cudaFree(impl->diag->d_interp_in));
74ca735530SJeremy L Thompson     CeedCallCuda(ceed, cudaFree(impl->diag->d_interp_out));
75ca735530SJeremy L Thompson     CeedCallCuda(ceed, cudaFree(impl->diag->d_grad_in));
76ca735530SJeremy L Thompson     CeedCallCuda(ceed, cudaFree(impl->diag->d_grad_out));
77004e4986SSebastian Grimberg     CeedCallCuda(ceed, cudaFree(impl->diag->d_div_in));
78004e4986SSebastian Grimberg     CeedCallCuda(ceed, cudaFree(impl->diag->d_div_out));
79004e4986SSebastian Grimberg     CeedCallCuda(ceed, cudaFree(impl->diag->d_curl_in));
80004e4986SSebastian Grimberg     CeedCallCuda(ceed, cudaFree(impl->diag->d_curl_out));
81004e4986SSebastian Grimberg     CeedCallBackend(CeedElemRestrictionDestroy(&impl->diag->diag_rstr));
82506b1a0cSSebastian Grimberg     CeedCallBackend(CeedElemRestrictionDestroy(&impl->diag->point_block_diag_rstr));
83ca735530SJeremy L Thompson     CeedCallBackend(CeedVectorDestroy(&impl->diag->elem_diag));
84ca735530SJeremy L Thompson     CeedCallBackend(CeedVectorDestroy(&impl->diag->point_block_elem_diag));
850d0321e0SJeremy L Thompson   }
862b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&impl->diag));
870d0321e0SJeremy L Thompson 
88cc132f9aSnbeams   if (impl->asmb) {
89cc132f9aSnbeams     Ceed ceed;
90ca735530SJeremy L Thompson 
912b730f8bSJeremy L Thompson     CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
922b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cuModuleUnload(impl->asmb->module));
932b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaFree(impl->asmb->d_B_in));
942b730f8bSJeremy L Thompson     CeedCallCuda(ceed, cudaFree(impl->asmb->d_B_out));
95cc132f9aSnbeams   }
962b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&impl->asmb));
97cc132f9aSnbeams 
982b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&impl));
990d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1000d0321e0SJeremy L Thompson }
1010d0321e0SJeremy L Thompson 
1020d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
1030d0321e0SJeremy L Thompson // Setup infields or outfields
1040d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
105f8a0df59SJeremy L Thompson static int CeedOperatorSetupFields_Cuda(CeedQFunction qf, CeedOperator op, bool is_input, bool is_at_points, bool *skip_rstr, bool *apply_add_basis,
106034f99fdSJeremy L Thompson                                         CeedVector *e_vecs, CeedVector *q_vecs, CeedInt num_fields, CeedInt Q, CeedInt num_elem) {
1070d0321e0SJeremy L Thompson   Ceed                ceed;
108ca735530SJeremy L Thompson   CeedQFunctionField *qf_fields;
109ca735530SJeremy L Thompson   CeedOperatorField  *op_fields;
1100d0321e0SJeremy L Thompson 
111ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
112ca735530SJeremy L Thompson   if (is_input) {
113ca735530SJeremy L Thompson     CeedCallBackend(CeedOperatorGetFields(op, NULL, &op_fields, NULL, NULL));
114ca735530SJeremy L Thompson     CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_fields, NULL, NULL));
1150d0321e0SJeremy L Thompson   } else {
116ca735530SJeremy L Thompson     CeedCallBackend(CeedOperatorGetFields(op, NULL, NULL, NULL, &op_fields));
117ca735530SJeremy L Thompson     CeedCallBackend(CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qf_fields));
1180d0321e0SJeremy L Thompson   }
1190d0321e0SJeremy L Thompson 
1200d0321e0SJeremy L Thompson   // Loop over fields
121ca735530SJeremy L Thompson   for (CeedInt i = 0; i < num_fields; i++) {
1228bbba8cdSJeremy L Thompson     bool                is_active = false, is_strided = false, skip_e_vec = false;
123004e4986SSebastian Grimberg     CeedSize            q_size;
124004e4986SSebastian Grimberg     CeedInt             size;
125004e4986SSebastian Grimberg     CeedEvalMode        eval_mode;
1268bbba8cdSJeremy L Thompson     CeedVector          l_vec;
127edb2538eSJeremy L Thompson     CeedElemRestriction elem_rstr;
128ca735530SJeremy L Thompson 
1290d0321e0SJeremy L Thompson     // Check whether this field can skip the element restriction:
1308bbba8cdSJeremy L Thompson     // Input CEED_VECTOR_ACTIVE
1318bbba8cdSJeremy L Thompson     // Output CEED_VECTOR_ACTIVE without CEED_EVAL_NONE
1328bbba8cdSJeremy L Thompson     // Input CEED_VECTOR_NONE with CEED_EVAL_WEIGHT
1338bbba8cdSJeremy L Thompson     // Input passive vectorr with CEED_EVAL_NONE and strided restriction with CEED_STRIDES_BACKEND
134034f99fdSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &l_vec));
1358bbba8cdSJeremy L Thompson     is_active = l_vec == CEED_VECTOR_ACTIVE;
1368bbba8cdSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_fields[i], &elem_rstr));
1378bbba8cdSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode));
1388bbba8cdSJeremy L Thompson     skip_e_vec = (is_input && is_active) || (is_active && eval_mode != CEED_EVAL_NONE) || (eval_mode == CEED_EVAL_WEIGHT);
1398bbba8cdSJeremy L Thompson     if (!skip_e_vec && is_input && !is_active && eval_mode == CEED_EVAL_NONE) {
140edb2538eSJeremy L Thompson       CeedCallBackend(CeedElemRestrictionIsStrided(elem_rstr, &is_strided));
1418bbba8cdSJeremy L Thompson       if (is_strided) CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &skip_e_vec));
1420d0321e0SJeremy L Thompson     }
143034f99fdSJeremy L Thompson     if (skip_e_vec) {
144034f99fdSJeremy L Thompson       e_vecs[i] = NULL;
1450d0321e0SJeremy L Thompson     } else {
146034f99fdSJeremy L Thompson       CeedCallBackend(CeedElemRestrictionCreateVector(elem_rstr, NULL, &e_vecs[i]));
1470d0321e0SJeremy L Thompson     }
1480d0321e0SJeremy L Thompson 
149004e4986SSebastian Grimberg     switch (eval_mode) {
1500d0321e0SJeremy L Thompson       case CEED_EVAL_NONE:
1510d0321e0SJeremy L Thompson       case CEED_EVAL_INTERP:
1520d0321e0SJeremy L Thompson       case CEED_EVAL_GRAD:
153004e4986SSebastian Grimberg       case CEED_EVAL_DIV:
154004e4986SSebastian Grimberg       case CEED_EVAL_CURL:
155ca735530SJeremy L Thompson         CeedCallBackend(CeedQFunctionFieldGetSize(qf_fields[i], &size));
1568bbba8cdSJeremy L Thompson         q_size = (CeedSize)num_elem * (CeedSize)Q * (CeedSize)size;
157ca735530SJeremy L Thompson         CeedCallBackend(CeedVectorCreate(ceed, q_size, &q_vecs[i]));
1580d0321e0SJeremy L Thompson         break;
1598bbba8cdSJeremy L Thompson       case CEED_EVAL_WEIGHT: {
1608bbba8cdSJeremy L Thompson         CeedBasis basis;
1618bbba8cdSJeremy L Thompson 
162ca735530SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis));
1638bbba8cdSJeremy L Thompson         q_size = (CeedSize)num_elem * (CeedSize)Q;
164ca735530SJeremy L Thompson         CeedCallBackend(CeedVectorCreate(ceed, q_size, &q_vecs[i]));
165756ca9e9SJeremy L Thompson         if (is_at_points) {
166756ca9e9SJeremy L Thompson           CeedInt num_points[num_elem];
167756ca9e9SJeremy L Thompson 
168756ca9e9SJeremy L Thompson           for (CeedInt i = 0; i < num_elem; i++) num_points[i] = Q;
169756ca9e9SJeremy L Thompson           CeedCallBackend(
170756ca9e9SJeremy L Thompson               CeedBasisApplyAtPoints(basis, num_elem, num_points, CEED_NOTRANSPOSE, CEED_EVAL_WEIGHT, CEED_VECTOR_NONE, CEED_VECTOR_NONE, q_vecs[i]));
171756ca9e9SJeremy L Thompson         } else {
172ca735530SJeremy L Thompson           CeedCallBackend(CeedBasisApply(basis, num_elem, CEED_NOTRANSPOSE, CEED_EVAL_WEIGHT, CEED_VECTOR_NONE, q_vecs[i]));
173756ca9e9SJeremy L Thompson         }
1740d0321e0SJeremy L Thompson         break;
1750d0321e0SJeremy L Thompson       }
1760d0321e0SJeremy L Thompson     }
1778bbba8cdSJeremy L Thompson   }
178f8a0df59SJeremy L Thompson   // Drop duplicate restrictions
1793aab95c0SJeremy L Thompson   if (is_input) {
1803aab95c0SJeremy L Thompson     for (CeedInt i = 0; i < num_fields; i++) {
1813aab95c0SJeremy L Thompson       CeedVector          vec_i;
1823aab95c0SJeremy L Thompson       CeedElemRestriction rstr_i;
1833aab95c0SJeremy L Thompson 
1843aab95c0SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec_i));
1853aab95c0SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_fields[i], &rstr_i));
1863aab95c0SJeremy L Thompson       for (CeedInt j = i + 1; j < num_fields; j++) {
1873aab95c0SJeremy L Thompson         CeedVector          vec_j;
1883aab95c0SJeremy L Thompson         CeedElemRestriction rstr_j;
1893aab95c0SJeremy L Thompson 
1903aab95c0SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetVector(op_fields[j], &vec_j));
1913aab95c0SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_fields[j], &rstr_j));
1923aab95c0SJeremy L Thompson         if (vec_i == vec_j && rstr_i == rstr_j) {
1938bbba8cdSJeremy L Thompson           if (e_vecs[i]) CeedCallBackend(CeedVectorReferenceCopy(e_vecs[i], &e_vecs[j]));
1943aab95c0SJeremy L Thompson           skip_rstr[j] = true;
1953aab95c0SJeremy L Thompson         }
1963aab95c0SJeremy L Thompson       }
1973aab95c0SJeremy L Thompson     }
198f8a0df59SJeremy L Thompson   } else {
199f8a0df59SJeremy L Thompson     for (CeedInt i = num_fields - 1; i >= 0; i--) {
200f8a0df59SJeremy L Thompson       CeedVector          vec_i;
201f8a0df59SJeremy L Thompson       CeedElemRestriction rstr_i;
202f8a0df59SJeremy L Thompson 
203f8a0df59SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec_i));
204f8a0df59SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_fields[i], &rstr_i));
205f8a0df59SJeremy L Thompson       for (CeedInt j = i - 1; j >= 0; j--) {
206f8a0df59SJeremy L Thompson         CeedVector          vec_j;
207f8a0df59SJeremy L Thompson         CeedElemRestriction rstr_j;
208f8a0df59SJeremy L Thompson 
209f8a0df59SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetVector(op_fields[j], &vec_j));
210f8a0df59SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_fields[j], &rstr_j));
211f8a0df59SJeremy L Thompson         if (vec_i == vec_j && rstr_i == rstr_j) {
2128bbba8cdSJeremy L Thompson           if (e_vecs[i]) CeedCallBackend(CeedVectorReferenceCopy(e_vecs[i], &e_vecs[j]));
213f8a0df59SJeremy L Thompson           skip_rstr[j]       = true;
214f8a0df59SJeremy L Thompson           apply_add_basis[i] = true;
215f8a0df59SJeremy L Thompson         }
216f8a0df59SJeremy L Thompson       }
217f8a0df59SJeremy L Thompson     }
2183aab95c0SJeremy L Thompson   }
2190d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2200d0321e0SJeremy L Thompson }
2210d0321e0SJeremy L Thompson 
2220d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
223ea61e9acSJeremy L Thompson // CeedOperator needs to connect all the named fields (be they active or passive) to the named inputs and outputs of its CeedQFunction.
2240d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
2250d0321e0SJeremy L Thompson static int CeedOperatorSetup_Cuda(CeedOperator op) {
2260d0321e0SJeremy L Thompson   Ceed                ceed;
227ca735530SJeremy L Thompson   bool                is_setup_done;
228ca735530SJeremy L Thompson   CeedInt             Q, num_elem, num_input_fields, num_output_fields;
229ca735530SJeremy L Thompson   CeedQFunctionField *qf_input_fields, *qf_output_fields;
2300d0321e0SJeremy L Thompson   CeedQFunction       qf;
231ca735530SJeremy L Thompson   CeedOperatorField  *op_input_fields, *op_output_fields;
232ca735530SJeremy L Thompson   CeedOperator_Cuda  *impl;
233ca735530SJeremy L Thompson 
234ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorIsSetupDone(op, &is_setup_done));
235ca735530SJeremy L Thompson   if (is_setup_done) return CEED_ERROR_SUCCESS;
236ca735530SJeremy L Thompson 
237ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
238ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &impl));
2392b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
2402b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q));
241ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem));
242ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
243ca735530SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
2440d0321e0SJeremy L Thompson 
2450d0321e0SJeremy L Thompson   // Allocate
246034f99fdSJeremy L Thompson   CeedCallBackend(CeedCalloc(num_input_fields, &impl->e_vecs_in));
247034f99fdSJeremy L Thompson   CeedCallBackend(CeedCalloc(num_output_fields, &impl->e_vecs_out));
248034f99fdSJeremy L Thompson   CeedCallBackend(CeedCalloc(num_input_fields, &impl->skip_rstr_in));
249034f99fdSJeremy L Thompson   CeedCallBackend(CeedCalloc(num_output_fields, &impl->skip_rstr_out));
250034f99fdSJeremy L Thompson   CeedCallBackend(CeedCalloc(num_output_fields, &impl->apply_add_basis_out));
251034f99fdSJeremy L Thompson   CeedCallBackend(CeedCalloc(num_input_fields, &impl->input_field_order));
252034f99fdSJeremy L Thompson   CeedCallBackend(CeedCalloc(num_output_fields, &impl->output_field_order));
253034f99fdSJeremy L Thompson   CeedCallBackend(CeedCalloc(num_input_fields, &impl->input_states));
254034f99fdSJeremy L Thompson   CeedCallBackend(CeedCalloc(num_input_fields, &impl->q_vecs_in));
255034f99fdSJeremy L Thompson   CeedCallBackend(CeedCalloc(num_output_fields, &impl->q_vecs_out));
256ca735530SJeremy L Thompson   impl->num_inputs  = num_input_fields;
257ca735530SJeremy L Thompson   impl->num_outputs = num_output_fields;
2580d0321e0SJeremy L Thompson 
25941655a23SJeremy L Thompson   // Set up infield and outfield e-vecs and q-vecs
2603aab95c0SJeremy L Thompson   CeedCallBackend(
261034f99fdSJeremy L Thompson       CeedOperatorSetupFields_Cuda(qf, op, true, false, impl->skip_rstr_in, NULL, impl->e_vecs_in, impl->q_vecs_in, num_input_fields, Q, num_elem));
262034f99fdSJeremy L Thompson   CeedCallBackend(CeedOperatorSetupFields_Cuda(qf, op, false, false, impl->skip_rstr_out, impl->apply_add_basis_out, impl->e_vecs_out,
263034f99fdSJeremy L Thompson                                                impl->q_vecs_out, num_output_fields, Q, num_elem));
2640d0321e0SJeremy L Thompson 
265034f99fdSJeremy L Thompson   // Reorder fields to allow reuse of buffers
266034f99fdSJeremy L Thompson   impl->max_active_e_vec_len = 0;
26741655a23SJeremy L Thompson   {
268034f99fdSJeremy L Thompson     bool    is_ordered[CEED_FIELD_MAX];
269034f99fdSJeremy L Thompson     CeedInt curr_index = 0;
27041655a23SJeremy L Thompson 
271034f99fdSJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) is_ordered[i] = false;
27241655a23SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
273034f99fdSJeremy L Thompson       CeedSize            e_vec_len_i;
27441655a23SJeremy L Thompson       CeedVector          vec_i;
27541655a23SJeremy L Thompson       CeedElemRestriction rstr_i;
27641655a23SJeremy L Thompson 
277034f99fdSJeremy L Thompson       if (is_ordered[i]) continue;
278034f99fdSJeremy L Thompson       is_ordered[i]                       = true;
279034f99fdSJeremy L Thompson       impl->input_field_order[curr_index] = i;
280034f99fdSJeremy L Thompson       curr_index++;
28141655a23SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec_i));
282034f99fdSJeremy L Thompson       if (vec_i == CEED_VECTOR_NONE) continue;  // CEED_EVAL_WEIGHT
28341655a23SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &rstr_i));
284034f99fdSJeremy L Thompson       CeedCallBackend(CeedElemRestrictionGetEVectorSize(rstr_i, &e_vec_len_i));
285034f99fdSJeremy L Thompson       impl->max_active_e_vec_len = e_vec_len_i > impl->max_active_e_vec_len ? e_vec_len_i : impl->max_active_e_vec_len;
286034f99fdSJeremy L Thompson       for (CeedInt j = i + 1; j < num_input_fields; j++) {
287034f99fdSJeremy L Thompson         CeedVector          vec_j;
288034f99fdSJeremy L Thompson         CeedElemRestriction rstr_j;
289034f99fdSJeremy L Thompson 
290034f99fdSJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[j], &vec_j));
291034f99fdSJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[j], &rstr_j));
292034f99fdSJeremy L Thompson         if (rstr_i == rstr_j && vec_i == vec_j) {
293034f99fdSJeremy L Thompson           is_ordered[j]                       = true;
294034f99fdSJeremy L Thompson           impl->input_field_order[curr_index] = j;
295034f99fdSJeremy L Thompson           curr_index++;
29641655a23SJeremy L Thompson         }
297034f99fdSJeremy L Thompson       }
298034f99fdSJeremy L Thompson     }
299034f99fdSJeremy L Thompson   }
300034f99fdSJeremy L Thompson   {
301034f99fdSJeremy L Thompson     bool    is_ordered[CEED_FIELD_MAX];
302034f99fdSJeremy L Thompson     CeedInt curr_index = 0;
303034f99fdSJeremy L Thompson 
304034f99fdSJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) is_ordered[i] = false;
305034f99fdSJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
306034f99fdSJeremy L Thompson       CeedSize            e_vec_len_i;
307034f99fdSJeremy L Thompson       CeedVector          vec_i;
308034f99fdSJeremy L Thompson       CeedElemRestriction rstr_i;
309034f99fdSJeremy L Thompson 
310034f99fdSJeremy L Thompson       if (is_ordered[i]) continue;
311034f99fdSJeremy L Thompson       is_ordered[i]                        = true;
312034f99fdSJeremy L Thompson       impl->output_field_order[curr_index] = i;
313034f99fdSJeremy L Thompson       curr_index++;
314034f99fdSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[i], &vec_i));
315034f99fdSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &rstr_i));
316034f99fdSJeremy L Thompson       CeedCallBackend(CeedElemRestrictionGetEVectorSize(rstr_i, &e_vec_len_i));
317034f99fdSJeremy L Thompson       impl->max_active_e_vec_len = e_vec_len_i > impl->max_active_e_vec_len ? e_vec_len_i : impl->max_active_e_vec_len;
318034f99fdSJeremy L Thompson       for (CeedInt j = i + 1; j < num_output_fields; j++) {
31941655a23SJeremy L Thompson         CeedVector          vec_j;
32041655a23SJeremy L Thompson         CeedElemRestriction rstr_j;
32141655a23SJeremy L Thompson 
32241655a23SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[j], &vec_j));
32341655a23SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[j], &rstr_j));
324034f99fdSJeremy L Thompson         if (rstr_i == rstr_j && vec_i == vec_j) {
325034f99fdSJeremy L Thompson           is_ordered[j]                        = true;
326034f99fdSJeremy L Thompson           impl->output_field_order[curr_index] = j;
327034f99fdSJeremy L Thompson           curr_index++;
32841655a23SJeremy L Thompson         }
32941655a23SJeremy L Thompson       }
33041655a23SJeremy L Thompson     }
33141655a23SJeremy L Thompson   }
3322b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorSetSetupDone(op));
3330d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3340d0321e0SJeremy L Thompson }
3350d0321e0SJeremy L Thompson 
3360d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
33743e13feeSJeremy L Thompson // Restrict Operator Inputs
3380d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
33943e13feeSJeremy L Thompson static inline int CeedOperatorInputRestrict_Cuda(CeedOperatorField op_input_field, CeedQFunctionField qf_input_field, CeedInt input_field,
3408bbba8cdSJeremy L Thompson                                                  CeedVector in_vec, CeedVector active_e_vec, const bool skip_active, CeedOperator_Cuda *impl,
34143e13feeSJeremy L Thompson                                                  CeedRequest *request) {
3428bbba8cdSJeremy L Thompson   bool       is_active = false;
343034f99fdSJeremy L Thompson   CeedVector l_vec, e_vec = impl->e_vecs_in[input_field];
3440d0321e0SJeremy L Thompson 
3450d0321e0SJeremy L Thompson   // Get input vector
346034f99fdSJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetVector(op_input_field, &l_vec));
3478bbba8cdSJeremy L Thompson   is_active = l_vec == CEED_VECTOR_ACTIVE;
3488bbba8cdSJeremy L Thompson   if (is_active && skip_active) return CEED_ERROR_SUCCESS;
3498bbba8cdSJeremy L Thompson   if (is_active) {
3508bbba8cdSJeremy L Thompson     l_vec = in_vec;
3518bbba8cdSJeremy L Thompson     if (!e_vec) e_vec = active_e_vec;
3520d0321e0SJeremy L Thompson   }
3530d0321e0SJeremy L Thompson 
35443e13feeSJeremy L Thompson   // Restriction action
3558bbba8cdSJeremy L Thompson   if (e_vec) {
356034f99fdSJeremy L Thompson     // Restrict, if necessary
357034f99fdSJeremy L Thompson     if (!impl->skip_rstr_in[input_field]) {
358c1222711SJeremy L Thompson       uint64_t state;
359c1222711SJeremy L Thompson 
360034f99fdSJeremy L Thompson       CeedCallBackend(CeedVectorGetState(l_vec, &state));
3618bbba8cdSJeremy L Thompson       if (is_active || state != impl->input_states[input_field]) {
362034f99fdSJeremy L Thompson         CeedElemRestriction elem_rstr;
363034f99fdSJeremy L Thompson 
364034f99fdSJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_field, &elem_rstr));
365034f99fdSJeremy L Thompson         CeedCallBackend(CeedElemRestrictionApply(elem_rstr, CEED_NOTRANSPOSE, l_vec, e_vec, request));
366c1222711SJeremy L Thompson       }
36743e13feeSJeremy L Thompson       impl->input_states[input_field] = state;
368034f99fdSJeremy L Thompson     }
3690d0321e0SJeremy L Thompson   }
3700d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3710d0321e0SJeremy L Thompson }
3720d0321e0SJeremy L Thompson 
3730d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
3740d0321e0SJeremy L Thompson // Input Basis Action
3750d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
37643e13feeSJeremy L Thompson static inline int CeedOperatorInputBasis_Cuda(CeedOperatorField op_input_field, CeedQFunctionField qf_input_field, CeedInt input_field,
3778bbba8cdSJeremy L Thompson                                               CeedVector in_vec, CeedVector active_e_vec, CeedInt num_elem, const bool skip_active,
3788bbba8cdSJeremy L Thompson                                               CeedOperator_Cuda *impl) {
3798bbba8cdSJeremy L Thompson   bool         is_active = false;
380004e4986SSebastian Grimberg   CeedEvalMode eval_mode;
3818bbba8cdSJeremy L Thompson   CeedVector   l_vec, e_vec = impl->e_vecs_in[input_field], q_vec = impl->q_vecs_in[input_field];
3820d0321e0SJeremy L Thompson 
3830d0321e0SJeremy L Thompson   // Skip active input
384034f99fdSJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetVector(op_input_field, &l_vec));
3858bbba8cdSJeremy L Thompson   is_active = l_vec == CEED_VECTOR_ACTIVE;
3868bbba8cdSJeremy L Thompson   if (is_active && skip_active) return CEED_ERROR_SUCCESS;
3878bbba8cdSJeremy L Thompson   if (is_active) {
3888bbba8cdSJeremy L Thompson     l_vec = in_vec;
3898bbba8cdSJeremy L Thompson     if (!e_vec) e_vec = active_e_vec;
3900d0321e0SJeremy L Thompson   }
39143e13feeSJeremy L Thompson 
3920d0321e0SJeremy L Thompson   // Basis action
39343e13feeSJeremy L Thompson   CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_field, &eval_mode));
394004e4986SSebastian Grimberg   switch (eval_mode) {
3958bbba8cdSJeremy L Thompson     case CEED_EVAL_NONE: {
3968bbba8cdSJeremy L Thompson       const CeedScalar *e_vec_array;
3978bbba8cdSJeremy L Thompson 
3988bbba8cdSJeremy L Thompson       if (e_vec) {
3998bbba8cdSJeremy L Thompson         CeedCallBackend(CeedVectorGetArrayRead(e_vec, CEED_MEM_DEVICE, &e_vec_array));
4008bbba8cdSJeremy L Thompson       } else {
4018bbba8cdSJeremy L Thompson         CeedCallBackend(CeedVectorGetArrayRead(l_vec, CEED_MEM_DEVICE, &e_vec_array));
4028bbba8cdSJeremy L Thompson       }
4038bbba8cdSJeremy L Thompson       CeedCallBackend(CeedVectorSetArray(q_vec, CEED_MEM_DEVICE, CEED_USE_POINTER, (CeedScalar *)e_vec_array));
4040d0321e0SJeremy L Thompson       break;
4058bbba8cdSJeremy L Thompson     }
4060d0321e0SJeremy L Thompson     case CEED_EVAL_INTERP:
4070d0321e0SJeremy L Thompson     case CEED_EVAL_GRAD:
408004e4986SSebastian Grimberg     case CEED_EVAL_DIV:
40943e13feeSJeremy L Thompson     case CEED_EVAL_CURL: {
41043e13feeSJeremy L Thompson       CeedBasis basis;
41143e13feeSJeremy L Thompson 
41243e13feeSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetBasis(op_input_field, &basis));
413034f99fdSJeremy L Thompson       CeedCallBackend(CeedBasisApply(basis, num_elem, CEED_NOTRANSPOSE, eval_mode, e_vec, q_vec));
4140d0321e0SJeremy L Thompson       break;
41543e13feeSJeremy L Thompson     }
4160d0321e0SJeremy L Thompson     case CEED_EVAL_WEIGHT:
4170d0321e0SJeremy L Thompson       break;  // No action
4180d0321e0SJeremy L Thompson   }
4190d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
4200d0321e0SJeremy L Thompson }
4210d0321e0SJeremy L Thompson 
4220d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
4230d0321e0SJeremy L Thompson // Restore Input Vectors
4240d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
42543e13feeSJeremy L Thompson static inline int CeedOperatorInputRestore_Cuda(CeedOperatorField op_input_field, CeedQFunctionField qf_input_field, CeedInt input_field,
4268bbba8cdSJeremy L Thompson                                                 CeedVector in_vec, CeedVector active_e_vec, const bool skip_active, CeedOperator_Cuda *impl) {
4278bbba8cdSJeremy L Thompson   bool         is_active = false;
428004e4986SSebastian Grimberg   CeedEvalMode eval_mode;
429034f99fdSJeremy L Thompson   CeedVector   l_vec, e_vec = impl->e_vecs_in[input_field];
4300d0321e0SJeremy L Thompson 
4310d0321e0SJeremy L Thompson   // Skip active input
432034f99fdSJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetVector(op_input_field, &l_vec));
4338bbba8cdSJeremy L Thompson   is_active = l_vec == CEED_VECTOR_ACTIVE;
4348bbba8cdSJeremy L Thompson   if (is_active && skip_active) return CEED_ERROR_SUCCESS;
4358bbba8cdSJeremy L Thompson   if (is_active) {
4368bbba8cdSJeremy L Thompson     l_vec = in_vec;
4378bbba8cdSJeremy L Thompson     if (!e_vec) e_vec = active_e_vec;
4388bbba8cdSJeremy L Thompson   }
43943e13feeSJeremy L Thompson 
44043e13feeSJeremy L Thompson   // Restore e-vec
44143e13feeSJeremy L Thompson   CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_field, &eval_mode));
4428bbba8cdSJeremy L Thompson   if (eval_mode == CEED_EVAL_NONE) {
4438bbba8cdSJeremy L Thompson     const CeedScalar *e_vec_array;
4448bbba8cdSJeremy L Thompson 
4458bbba8cdSJeremy L Thompson     CeedCallBackend(CeedVectorTakeArray(impl->q_vecs_in[input_field], CEED_MEM_DEVICE, (CeedScalar **)&e_vec_array));
4468bbba8cdSJeremy L Thompson     if (e_vec) {
4478bbba8cdSJeremy L Thompson       CeedCallBackend(CeedVectorRestoreArrayRead(e_vec, &e_vec_array));
4480d0321e0SJeremy L Thompson     } else {
4498bbba8cdSJeremy L Thompson       CeedCallBackend(CeedVectorRestoreArrayRead(l_vec, &e_vec_array));
4500d0321e0SJeremy L Thompson     }
4510d0321e0SJeremy L Thompson   }
4520d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
4530d0321e0SJeremy L Thompson }
4540d0321e0SJeremy L Thompson 
4550d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
4560d0321e0SJeremy L Thompson // Apply and add to output
4570d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
458ca735530SJeremy L Thompson static int CeedOperatorApplyAdd_Cuda(CeedOperator op, CeedVector in_vec, CeedVector out_vec, CeedRequest *request) {
459034f99fdSJeremy L Thompson   CeedInt             Q, num_elem, num_input_fields, num_output_fields;
4608bbba8cdSJeremy L Thompson   Ceed                ceed;
4618bbba8cdSJeremy L Thompson   CeedVector          active_e_vec;
462ca735530SJeremy L Thompson   CeedQFunctionField *qf_input_fields, *qf_output_fields;
4630d0321e0SJeremy L Thompson   CeedQFunction       qf;
464004e4986SSebastian Grimberg   CeedOperatorField  *op_input_fields, *op_output_fields;
465004e4986SSebastian Grimberg   CeedOperator_Cuda  *impl;
466ca735530SJeremy L Thompson 
4678bbba8cdSJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
468ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &impl));
4692b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
4702b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q));
471ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem));
472ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
473ca735530SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
4740d0321e0SJeremy L Thompson 
4750d0321e0SJeremy L Thompson   // Setup
4762b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorSetup_Cuda(op));
4770d0321e0SJeremy L Thompson 
4788bbba8cdSJeremy L Thompson   // Work vector
4798bbba8cdSJeremy L Thompson   CeedCallBackend(CeedGetWorkVector(ceed, impl->max_active_e_vec_len, &active_e_vec));
4808bbba8cdSJeremy L Thompson 
48143e13feeSJeremy L Thompson   // Process inputs
48243e13feeSJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
483034f99fdSJeremy L Thompson     CeedInt field = impl->input_field_order[i];
484034f99fdSJeremy L Thompson 
485034f99fdSJeremy L Thompson     CeedCallBackend(
4868bbba8cdSJeremy L Thompson         CeedOperatorInputRestrict_Cuda(op_input_fields[field], qf_input_fields[field], field, in_vec, active_e_vec, false, impl, request));
4878bbba8cdSJeremy L Thompson     CeedCallBackend(CeedOperatorInputBasis_Cuda(op_input_fields[field], qf_input_fields[field], field, in_vec, active_e_vec, num_elem, false, impl));
48843e13feeSJeremy L Thompson   }
4890d0321e0SJeremy L Thompson 
4900d0321e0SJeremy L Thompson   // Output pointers, as necessary
491ca735530SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
492004e4986SSebastian Grimberg     CeedEvalMode eval_mode;
493004e4986SSebastian Grimberg 
494004e4986SSebastian Grimberg     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
495004e4986SSebastian Grimberg     if (eval_mode == CEED_EVAL_NONE) {
4968bbba8cdSJeremy L Thompson       CeedScalar *e_vec_array;
4978bbba8cdSJeremy L Thompson 
4988bbba8cdSJeremy L Thompson       CeedCallBackend(CeedVectorGetArrayWrite(impl->e_vecs_out[i], CEED_MEM_DEVICE, &e_vec_array));
4998bbba8cdSJeremy L Thompson       CeedCallBackend(CeedVectorSetArray(impl->q_vecs_out[i], CEED_MEM_DEVICE, CEED_USE_POINTER, e_vec_array));
5000d0321e0SJeremy L Thompson     }
5010d0321e0SJeremy L Thompson   }
5020d0321e0SJeremy L Thompson 
5030d0321e0SJeremy L Thompson   // Q function
504ca735530SJeremy L Thompson   CeedCallBackend(CeedQFunctionApply(qf, num_elem * Q, impl->q_vecs_in, impl->q_vecs_out));
5050d0321e0SJeremy L Thompson 
50641655a23SJeremy L Thompson   // Restore input arrays
50743e13feeSJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
5088bbba8cdSJeremy L Thompson     CeedCallBackend(CeedOperatorInputRestore_Cuda(op_input_fields[i], qf_input_fields[i], i, in_vec, active_e_vec, false, impl));
50943e13feeSJeremy L Thompson   }
51041655a23SJeremy L Thompson 
5118bbba8cdSJeremy L Thompson   // Output basis and restriction
512ca735530SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
5138bbba8cdSJeremy L Thompson     bool                is_active = false;
514034f99fdSJeremy L Thompson     CeedInt             field     = impl->output_field_order[i];
515004e4986SSebastian Grimberg     CeedEvalMode        eval_mode;
516034f99fdSJeremy L Thompson     CeedVector          l_vec, e_vec = impl->e_vecs_out[field], q_vec = impl->q_vecs_out[field];
517edb2538eSJeremy L Thompson     CeedElemRestriction elem_rstr;
518ca735530SJeremy L Thompson     CeedBasis           basis;
519ca735530SJeremy L Thompson 
520034f99fdSJeremy L Thompson     // Output vector
521034f99fdSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[field], &l_vec));
5228bbba8cdSJeremy L Thompson     is_active = l_vec == CEED_VECTOR_ACTIVE;
5238bbba8cdSJeremy L Thompson     if (is_active) {
5248bbba8cdSJeremy L Thompson       l_vec = out_vec;
5258bbba8cdSJeremy L Thompson       if (!e_vec) e_vec = active_e_vec;
5268bbba8cdSJeremy L Thompson     }
527034f99fdSJeremy L Thompson 
5280d0321e0SJeremy L Thompson     // Basis action
529034f99fdSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[field], &eval_mode));
530004e4986SSebastian Grimberg     switch (eval_mode) {
5310d0321e0SJeremy L Thompson       case CEED_EVAL_NONE:
532004e4986SSebastian Grimberg         break;  // No action
5330d0321e0SJeremy L Thompson       case CEED_EVAL_INTERP:
5340d0321e0SJeremy L Thompson       case CEED_EVAL_GRAD:
535004e4986SSebastian Grimberg       case CEED_EVAL_DIV:
536004e4986SSebastian Grimberg       case CEED_EVAL_CURL:
537034f99fdSJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[field], &basis));
538034f99fdSJeremy L Thompson         if (impl->apply_add_basis_out[field]) {
539034f99fdSJeremy L Thompson           CeedCallBackend(CeedBasisApplyAdd(basis, num_elem, CEED_TRANSPOSE, eval_mode, q_vec, e_vec));
540f8a0df59SJeremy L Thompson         } else {
541034f99fdSJeremy L Thompson           CeedCallBackend(CeedBasisApply(basis, num_elem, CEED_TRANSPOSE, eval_mode, q_vec, e_vec));
542f8a0df59SJeremy L Thompson         }
5430d0321e0SJeremy L Thompson         break;
5440d0321e0SJeremy L Thompson       // LCOV_EXCL_START
5450d0321e0SJeremy L Thompson       case CEED_EVAL_WEIGHT: {
5468bbba8cdSJeremy L Thompson         return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_WEIGHT cannot be an output evaluation mode");
5470d0321e0SJeremy L Thompson         // LCOV_EXCL_STOP
5480d0321e0SJeremy L Thompson       }
5490d0321e0SJeremy L Thompson     }
550ca735530SJeremy L Thompson 
5510d0321e0SJeremy L Thompson     // Restore evec
552004e4986SSebastian Grimberg     if (eval_mode == CEED_EVAL_NONE) {
5538bbba8cdSJeremy L Thompson       CeedScalar *e_vec_array;
5548bbba8cdSJeremy L Thompson 
5558bbba8cdSJeremy L Thompson       CeedCallBackend(CeedVectorTakeArray(impl->q_vecs_out[i], CEED_MEM_DEVICE, &e_vec_array));
5568bbba8cdSJeremy L Thompson       CeedCallBackend(CeedVectorRestoreArray(e_vec, &e_vec_array));
5570d0321e0SJeremy L Thompson     }
5580d0321e0SJeremy L Thompson 
559034f99fdSJeremy L Thompson     // Restrict
560034f99fdSJeremy L Thompson     if (impl->skip_rstr_out[field]) continue;
561034f99fdSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[field], &elem_rstr));
562034f99fdSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionApply(elem_rstr, CEED_TRANSPOSE, e_vec, l_vec, request));
5630d0321e0SJeremy L Thompson   }
5648bbba8cdSJeremy L Thompson 
5658bbba8cdSJeremy L Thompson   // Return work vector
5668bbba8cdSJeremy L Thompson   CeedCallBackend(CeedRestoreWorkVector(ceed, &active_e_vec));
5670d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
5680d0321e0SJeremy L Thompson }
5690d0321e0SJeremy L Thompson 
5700d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
571756ca9e9SJeremy L Thompson // CeedOperator needs to connect all the named fields (be they active or passive) to the named inputs and outputs of its CeedQFunction.
572756ca9e9SJeremy L Thompson //------------------------------------------------------------------------------
573756ca9e9SJeremy L Thompson static int CeedOperatorSetupAtPoints_Cuda(CeedOperator op) {
574756ca9e9SJeremy L Thompson   Ceed                ceed;
575756ca9e9SJeremy L Thompson   bool                is_setup_done;
576756ca9e9SJeremy L Thompson   CeedInt             max_num_points = -1, num_elem, num_input_fields, num_output_fields;
577756ca9e9SJeremy L Thompson   CeedQFunctionField *qf_input_fields, *qf_output_fields;
578756ca9e9SJeremy L Thompson   CeedQFunction       qf;
579756ca9e9SJeremy L Thompson   CeedOperatorField  *op_input_fields, *op_output_fields;
580756ca9e9SJeremy L Thompson   CeedOperator_Cuda  *impl;
581756ca9e9SJeremy L Thompson 
582756ca9e9SJeremy L Thompson   CeedCallBackend(CeedOperatorIsSetupDone(op, &is_setup_done));
583756ca9e9SJeremy L Thompson   if (is_setup_done) return CEED_ERROR_SUCCESS;
584756ca9e9SJeremy L Thompson 
585756ca9e9SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
586756ca9e9SJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &impl));
587756ca9e9SJeremy L Thompson   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
588756ca9e9SJeremy L Thompson   CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem));
589756ca9e9SJeremy L Thompson   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
590756ca9e9SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
591756ca9e9SJeremy L Thompson   {
592111870feSJeremy L Thompson     CeedElemRestriction rstr_points = NULL;
593756ca9e9SJeremy L Thompson 
594111870feSJeremy L Thompson     CeedCallBackend(CeedOperatorAtPointsGetPoints(op, &rstr_points, NULL));
595111870feSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetMaxPointsInElement(rstr_points, &max_num_points));
596111870feSJeremy L Thompson     CeedCallBackend(CeedCalloc(num_elem, &impl->num_points));
597111870feSJeremy L Thompson     for (CeedInt e = 0; e < num_elem; e++) {
598111870feSJeremy L Thompson       CeedInt num_points_elem;
599111870feSJeremy L Thompson 
600111870feSJeremy L Thompson       CeedCallBackend(CeedElemRestrictionGetNumPointsInElement(rstr_points, e, &num_points_elem));
601111870feSJeremy L Thompson       impl->num_points[e] = num_points_elem;
602111870feSJeremy L Thompson     }
603*8b0f7348SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionDestroy(&rstr_points));
604756ca9e9SJeremy L Thompson   }
605756ca9e9SJeremy L Thompson   impl->max_num_points = max_num_points;
606756ca9e9SJeremy L Thompson 
607756ca9e9SJeremy L Thompson   // Allocate
608034f99fdSJeremy L Thompson   CeedCallBackend(CeedCalloc(num_input_fields, &impl->e_vecs_in));
609034f99fdSJeremy L Thompson   CeedCallBackend(CeedCalloc(num_output_fields, &impl->e_vecs_out));
610034f99fdSJeremy L Thompson   CeedCallBackend(CeedCalloc(num_input_fields, &impl->skip_rstr_in));
611034f99fdSJeremy L Thompson   CeedCallBackend(CeedCalloc(num_output_fields, &impl->skip_rstr_out));
612034f99fdSJeremy L Thompson   CeedCallBackend(CeedCalloc(num_output_fields, &impl->apply_add_basis_out));
613034f99fdSJeremy L Thompson   CeedCallBackend(CeedCalloc(num_input_fields, &impl->input_field_order));
614034f99fdSJeremy L Thompson   CeedCallBackend(CeedCalloc(num_output_fields, &impl->output_field_order));
615034f99fdSJeremy L Thompson   CeedCallBackend(CeedCalloc(num_input_fields, &impl->input_states));
616034f99fdSJeremy L Thompson   CeedCallBackend(CeedCalloc(num_input_fields, &impl->q_vecs_in));
617034f99fdSJeremy L Thompson   CeedCallBackend(CeedCalloc(num_output_fields, &impl->q_vecs_out));
618756ca9e9SJeremy L Thompson   impl->num_inputs  = num_input_fields;
619756ca9e9SJeremy L Thompson   impl->num_outputs = num_output_fields;
620756ca9e9SJeremy L Thompson 
6218a213570SJeremy L Thompson   // Set up infield and outfield e-vecs and q-vecs
622034f99fdSJeremy L Thompson   CeedCallBackend(CeedOperatorSetupFields_Cuda(qf, op, true, true, impl->skip_rstr_in, NULL, impl->e_vecs_in, impl->q_vecs_in, num_input_fields,
6233aab95c0SJeremy L Thompson                                                max_num_points, num_elem));
624034f99fdSJeremy L Thompson   CeedCallBackend(CeedOperatorSetupFields_Cuda(qf, op, false, true, impl->skip_rstr_out, impl->apply_add_basis_out, impl->e_vecs_out,
625034f99fdSJeremy L Thompson                                                impl->q_vecs_out, num_output_fields, max_num_points, num_elem));
626756ca9e9SJeremy L Thompson 
627034f99fdSJeremy L Thompson   // Reorder fields to allow reuse of buffers
6288bbba8cdSJeremy L Thompson   impl->max_active_e_vec_len = 0;
6298a213570SJeremy L Thompson   {
630034f99fdSJeremy L Thompson     bool    is_ordered[CEED_FIELD_MAX];
631034f99fdSJeremy L Thompson     CeedInt curr_index = 0;
6328a213570SJeremy L Thompson 
633034f99fdSJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) is_ordered[i] = false;
6348a213570SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
6358bbba8cdSJeremy L Thompson       CeedSize            e_vec_len_i;
6368a213570SJeremy L Thompson       CeedVector          vec_i;
6378a213570SJeremy L Thompson       CeedElemRestriction rstr_i;
6388a213570SJeremy L Thompson 
639034f99fdSJeremy L Thompson       if (is_ordered[i]) continue;
640034f99fdSJeremy L Thompson       is_ordered[i]                       = true;
641034f99fdSJeremy L Thompson       impl->input_field_order[curr_index] = i;
642034f99fdSJeremy L Thompson       curr_index++;
6438a213570SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec_i));
644034f99fdSJeremy L Thompson       if (vec_i == CEED_VECTOR_NONE) continue;  // CEED_EVAL_WEIGHT
6458a213570SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &rstr_i));
6468bbba8cdSJeremy L Thompson       CeedCallBackend(CeedElemRestrictionGetEVectorSize(rstr_i, &e_vec_len_i));
6478bbba8cdSJeremy L Thompson       impl->max_active_e_vec_len = e_vec_len_i > impl->max_active_e_vec_len ? e_vec_len_i : impl->max_active_e_vec_len;
648034f99fdSJeremy L Thompson       for (CeedInt j = i + 1; j < num_input_fields; j++) {
649034f99fdSJeremy L Thompson         CeedVector          vec_j;
650034f99fdSJeremy L Thompson         CeedElemRestriction rstr_j;
651034f99fdSJeremy L Thompson 
652034f99fdSJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[j], &vec_j));
653034f99fdSJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[j], &rstr_j));
654034f99fdSJeremy L Thompson         if (rstr_i == rstr_j && vec_i == vec_j) {
655034f99fdSJeremy L Thompson           is_ordered[j]                       = true;
656034f99fdSJeremy L Thompson           impl->input_field_order[curr_index] = j;
657034f99fdSJeremy L Thompson           curr_index++;
6588a213570SJeremy L Thompson         }
659034f99fdSJeremy L Thompson       }
660034f99fdSJeremy L Thompson     }
661034f99fdSJeremy L Thompson   }
662034f99fdSJeremy L Thompson   {
663034f99fdSJeremy L Thompson     bool    is_ordered[CEED_FIELD_MAX];
664034f99fdSJeremy L Thompson     CeedInt curr_index = 0;
665034f99fdSJeremy L Thompson 
666034f99fdSJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) is_ordered[i] = false;
667034f99fdSJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
6688bbba8cdSJeremy L Thompson       CeedSize            e_vec_len_i;
669034f99fdSJeremy L Thompson       CeedVector          vec_i;
670034f99fdSJeremy L Thompson       CeedElemRestriction rstr_i;
671034f99fdSJeremy L Thompson 
672034f99fdSJeremy L Thompson       if (is_ordered[i]) continue;
673034f99fdSJeremy L Thompson       is_ordered[i]                        = true;
674034f99fdSJeremy L Thompson       impl->output_field_order[curr_index] = i;
675034f99fdSJeremy L Thompson       curr_index++;
676034f99fdSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[i], &vec_i));
677034f99fdSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &rstr_i));
6788bbba8cdSJeremy L Thompson       CeedCallBackend(CeedElemRestrictionGetEVectorSize(rstr_i, &e_vec_len_i));
6798bbba8cdSJeremy L Thompson       impl->max_active_e_vec_len = e_vec_len_i > impl->max_active_e_vec_len ? e_vec_len_i : impl->max_active_e_vec_len;
680034f99fdSJeremy L Thompson       for (CeedInt j = i + 1; j < num_output_fields; j++) {
6818a213570SJeremy L Thompson         CeedVector          vec_j;
6828a213570SJeremy L Thompson         CeedElemRestriction rstr_j;
6838a213570SJeremy L Thompson 
6848a213570SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[j], &vec_j));
6858a213570SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[j], &rstr_j));
686034f99fdSJeremy L Thompson         if (rstr_i == rstr_j && vec_i == vec_j) {
687034f99fdSJeremy L Thompson           is_ordered[j]                        = true;
688034f99fdSJeremy L Thompson           impl->output_field_order[curr_index] = j;
689034f99fdSJeremy L Thompson           curr_index++;
6908a213570SJeremy L Thompson         }
6918a213570SJeremy L Thompson       }
6928a213570SJeremy L Thompson     }
6938a213570SJeremy L Thompson   }
6948bbba8cdSJeremy L Thompson 
695756ca9e9SJeremy L Thompson   CeedCallBackend(CeedOperatorSetSetupDone(op));
696756ca9e9SJeremy L Thompson   return CEED_ERROR_SUCCESS;
697756ca9e9SJeremy L Thompson }
698756ca9e9SJeremy L Thompson 
699756ca9e9SJeremy L Thompson //------------------------------------------------------------------------------
70067d9480aSJeremy L Thompson // Input Basis Action AtPoints
701756ca9e9SJeremy L Thompson //------------------------------------------------------------------------------
70243e13feeSJeremy L Thompson static inline int CeedOperatorInputBasisAtPoints_Cuda(CeedOperatorField op_input_field, CeedQFunctionField qf_input_field, CeedInt input_field,
7038bbba8cdSJeremy L Thompson                                                       CeedVector in_vec, CeedVector active_e_vec, CeedInt num_elem, const CeedInt *num_points,
7048bbba8cdSJeremy L Thompson                                                       const bool skip_active, CeedOperator_Cuda *impl) {
7058bbba8cdSJeremy L Thompson   bool         is_active = false;
706756ca9e9SJeremy L Thompson   CeedEvalMode eval_mode;
7078bbba8cdSJeremy L Thompson   CeedVector   l_vec, e_vec = impl->e_vecs_in[input_field], q_vec = impl->q_vecs_in[input_field];
708756ca9e9SJeremy L Thompson 
709756ca9e9SJeremy L Thompson   // Skip active input
710034f99fdSJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetVector(op_input_field, &l_vec));
7118bbba8cdSJeremy L Thompson   is_active = l_vec == CEED_VECTOR_ACTIVE;
7128bbba8cdSJeremy L Thompson   if (is_active && skip_active) return CEED_ERROR_SUCCESS;
7138bbba8cdSJeremy L Thompson   if (is_active) {
7148bbba8cdSJeremy L Thompson     l_vec = in_vec;
7158bbba8cdSJeremy L Thompson     if (!e_vec) e_vec = active_e_vec;
716756ca9e9SJeremy L Thompson   }
71743e13feeSJeremy L Thompson 
718756ca9e9SJeremy L Thompson   // Basis action
71943e13feeSJeremy L Thompson   CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_field, &eval_mode));
720756ca9e9SJeremy L Thompson   switch (eval_mode) {
7218bbba8cdSJeremy L Thompson     case CEED_EVAL_NONE: {
7228bbba8cdSJeremy L Thompson       const CeedScalar *e_vec_array;
7238bbba8cdSJeremy L Thompson 
7248bbba8cdSJeremy L Thompson       if (e_vec) {
7258bbba8cdSJeremy L Thompson         CeedCallBackend(CeedVectorGetArrayRead(e_vec, CEED_MEM_DEVICE, &e_vec_array));
7268bbba8cdSJeremy L Thompson       } else {
7278bbba8cdSJeremy L Thompson         CeedCallBackend(CeedVectorGetArrayRead(l_vec, CEED_MEM_DEVICE, &e_vec_array));
7288bbba8cdSJeremy L Thompson       }
7298bbba8cdSJeremy L Thompson       CeedCallBackend(CeedVectorSetArray(q_vec, CEED_MEM_DEVICE, CEED_USE_POINTER, (CeedScalar *)e_vec_array));
730756ca9e9SJeremy L Thompson       break;
7318bbba8cdSJeremy L Thompson     }
732756ca9e9SJeremy L Thompson     case CEED_EVAL_INTERP:
733756ca9e9SJeremy L Thompson     case CEED_EVAL_GRAD:
734756ca9e9SJeremy L Thompson     case CEED_EVAL_DIV:
73543e13feeSJeremy L Thompson     case CEED_EVAL_CURL: {
73643e13feeSJeremy L Thompson       CeedBasis basis;
73743e13feeSJeremy L Thompson 
73843e13feeSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetBasis(op_input_field, &basis));
739034f99fdSJeremy L Thompson       CeedCallBackend(CeedBasisApplyAtPoints(basis, num_elem, num_points, CEED_NOTRANSPOSE, eval_mode, impl->point_coords_elem, e_vec, q_vec));
740756ca9e9SJeremy L Thompson       break;
74143e13feeSJeremy L Thompson     }
742756ca9e9SJeremy L Thompson     case CEED_EVAL_WEIGHT:
743756ca9e9SJeremy L Thompson       break;  // No action
744756ca9e9SJeremy L Thompson   }
745756ca9e9SJeremy L Thompson   return CEED_ERROR_SUCCESS;
746756ca9e9SJeremy L Thompson }
747756ca9e9SJeremy L Thompson 
748756ca9e9SJeremy L Thompson //------------------------------------------------------------------------------
74967d9480aSJeremy L Thompson // Apply and add to output AtPoints
750756ca9e9SJeremy L Thompson //------------------------------------------------------------------------------
751756ca9e9SJeremy L Thompson static int CeedOperatorApplyAddAtPoints_Cuda(CeedOperator op, CeedVector in_vec, CeedVector out_vec, CeedRequest *request) {
752034f99fdSJeremy L Thompson   CeedInt             max_num_points, *num_points, num_elem, num_input_fields, num_output_fields;
7538bbba8cdSJeremy L Thompson   Ceed                ceed;
7548bbba8cdSJeremy L Thompson   CeedVector          active_e_vec;
755756ca9e9SJeremy L Thompson   CeedQFunctionField *qf_input_fields, *qf_output_fields;
756756ca9e9SJeremy L Thompson   CeedQFunction       qf;
757756ca9e9SJeremy L Thompson   CeedOperatorField  *op_input_fields, *op_output_fields;
758756ca9e9SJeremy L Thompson   CeedOperator_Cuda  *impl;
759756ca9e9SJeremy L Thompson 
7608bbba8cdSJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
761756ca9e9SJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &impl));
762756ca9e9SJeremy L Thompson   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
763756ca9e9SJeremy L Thompson   CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem));
764756ca9e9SJeremy L Thompson   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
765756ca9e9SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
766756ca9e9SJeremy L Thompson 
767756ca9e9SJeremy L Thompson   // Setup
768756ca9e9SJeremy L Thompson   CeedCallBackend(CeedOperatorSetupAtPoints_Cuda(op));
769111870feSJeremy L Thompson   num_points     = impl->num_points;
770756ca9e9SJeremy L Thompson   max_num_points = impl->max_num_points;
771756ca9e9SJeremy L Thompson 
7728bbba8cdSJeremy L Thompson   // Work vector
7738bbba8cdSJeremy L Thompson   CeedCallBackend(CeedGetWorkVector(ceed, impl->max_active_e_vec_len, &active_e_vec));
7748bbba8cdSJeremy L Thompson 
775756ca9e9SJeremy L Thompson   // Get point coordinates
776756ca9e9SJeremy L Thompson   if (!impl->point_coords_elem) {
777756ca9e9SJeremy L Thompson     CeedVector          point_coords = NULL;
778756ca9e9SJeremy L Thompson     CeedElemRestriction rstr_points  = NULL;
779756ca9e9SJeremy L Thompson 
780756ca9e9SJeremy L Thompson     CeedCallBackend(CeedOperatorAtPointsGetPoints(op, &rstr_points, &point_coords));
781756ca9e9SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionCreateVector(rstr_points, NULL, &impl->point_coords_elem));
782756ca9e9SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionApply(rstr_points, CEED_NOTRANSPOSE, point_coords, impl->point_coords_elem, request));
783*8b0f7348SJeremy L Thompson     CeedCallBackend(CeedVectorDestroy(&point_coords));
784*8b0f7348SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionDestroy(&rstr_points));
785756ca9e9SJeremy L Thompson   }
786756ca9e9SJeremy L Thompson 
78743e13feeSJeremy L Thompson   // Process inputs
78843e13feeSJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
789034f99fdSJeremy L Thompson     CeedInt field = impl->input_field_order[i];
790034f99fdSJeremy L Thompson 
791034f99fdSJeremy L Thompson     CeedCallBackend(
7928bbba8cdSJeremy L Thompson         CeedOperatorInputRestrict_Cuda(op_input_fields[field], qf_input_fields[field], field, in_vec, active_e_vec, false, impl, request));
7938bbba8cdSJeremy L Thompson     CeedCallBackend(CeedOperatorInputBasisAtPoints_Cuda(op_input_fields[field], qf_input_fields[field], field, in_vec, active_e_vec, num_elem,
7948bbba8cdSJeremy L Thompson                                                         num_points, false, impl));
79543e13feeSJeremy L Thompson   }
796756ca9e9SJeremy L Thompson 
797756ca9e9SJeremy L Thompson   // Output pointers, as necessary
798756ca9e9SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
799756ca9e9SJeremy L Thompson     CeedEvalMode eval_mode;
800756ca9e9SJeremy L Thompson 
801756ca9e9SJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
802756ca9e9SJeremy L Thompson     if (eval_mode == CEED_EVAL_NONE) {
8038bbba8cdSJeremy L Thompson       CeedScalar *e_vec_array;
8048bbba8cdSJeremy L Thompson 
8058bbba8cdSJeremy L Thompson       CeedCallBackend(CeedVectorGetArrayWrite(impl->e_vecs_out[i], CEED_MEM_DEVICE, &e_vec_array));
8068bbba8cdSJeremy L Thompson       CeedCallBackend(CeedVectorSetArray(impl->q_vecs_out[i], CEED_MEM_DEVICE, CEED_USE_POINTER, e_vec_array));
807756ca9e9SJeremy L Thompson     }
808756ca9e9SJeremy L Thompson   }
809756ca9e9SJeremy L Thompson 
810756ca9e9SJeremy L Thompson   // Q function
811756ca9e9SJeremy L Thompson   CeedCallBackend(CeedQFunctionApply(qf, num_elem * max_num_points, impl->q_vecs_in, impl->q_vecs_out));
812756ca9e9SJeremy L Thompson 
8138a213570SJeremy L Thompson   // Restore input arrays
81443e13feeSJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
8158bbba8cdSJeremy L Thompson     CeedCallBackend(CeedOperatorInputRestore_Cuda(op_input_fields[i], qf_input_fields[i], i, in_vec, active_e_vec, false, impl));
81643e13feeSJeremy L Thompson   }
8178a213570SJeremy L Thompson 
8188bbba8cdSJeremy L Thompson   // Output basis and restriction
819756ca9e9SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
8208bbba8cdSJeremy L Thompson     bool                is_active = false;
821034f99fdSJeremy L Thompson     CeedInt             field     = impl->output_field_order[i];
822756ca9e9SJeremy L Thompson     CeedEvalMode        eval_mode;
823034f99fdSJeremy L Thompson     CeedVector          l_vec, e_vec = impl->e_vecs_out[field], q_vec = impl->q_vecs_out[field];
824756ca9e9SJeremy L Thompson     CeedElemRestriction elem_rstr;
825756ca9e9SJeremy L Thompson     CeedBasis           basis;
826756ca9e9SJeremy L Thompson 
827034f99fdSJeremy L Thompson     // Output vector
828034f99fdSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[field], &l_vec));
8298bbba8cdSJeremy L Thompson     is_active = l_vec == CEED_VECTOR_ACTIVE;
8308bbba8cdSJeremy L Thompson     if (is_active) {
8318bbba8cdSJeremy L Thompson       l_vec = out_vec;
8328bbba8cdSJeremy L Thompson       if (!e_vec) e_vec = active_e_vec;
8338bbba8cdSJeremy L Thompson     }
834034f99fdSJeremy L Thompson 
835756ca9e9SJeremy L Thompson     // Basis action
836034f99fdSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[field], &eval_mode));
837756ca9e9SJeremy L Thompson     switch (eval_mode) {
838756ca9e9SJeremy L Thompson       case CEED_EVAL_NONE:
839756ca9e9SJeremy L Thompson         break;  // No action
840756ca9e9SJeremy L Thompson       case CEED_EVAL_INTERP:
841756ca9e9SJeremy L Thompson       case CEED_EVAL_GRAD:
842756ca9e9SJeremy L Thompson       case CEED_EVAL_DIV:
843756ca9e9SJeremy L Thompson       case CEED_EVAL_CURL:
844034f99fdSJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[field], &basis));
845034f99fdSJeremy L Thompson         if (impl->apply_add_basis_out[field]) {
846034f99fdSJeremy L Thompson           CeedCallBackend(CeedBasisApplyAddAtPoints(basis, num_elem, num_points, CEED_TRANSPOSE, eval_mode, impl->point_coords_elem, q_vec, e_vec));
847f8a0df59SJeremy L Thompson         } else {
848034f99fdSJeremy L Thompson           CeedCallBackend(CeedBasisApplyAtPoints(basis, num_elem, num_points, CEED_TRANSPOSE, eval_mode, impl->point_coords_elem, q_vec, e_vec));
849f8a0df59SJeremy L Thompson         }
850756ca9e9SJeremy L Thompson         break;
851756ca9e9SJeremy L Thompson       // LCOV_EXCL_START
852756ca9e9SJeremy L Thompson       case CEED_EVAL_WEIGHT: {
8538bbba8cdSJeremy L Thompson         return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_WEIGHT cannot be an output evaluation mode");
854756ca9e9SJeremy L Thompson         // LCOV_EXCL_STOP
855756ca9e9SJeremy L Thompson       }
856756ca9e9SJeremy L Thompson     }
857756ca9e9SJeremy L Thompson 
858756ca9e9SJeremy L Thompson     // Restore evec
859756ca9e9SJeremy L Thompson     if (eval_mode == CEED_EVAL_NONE) {
8608bbba8cdSJeremy L Thompson       CeedScalar *e_vec_array;
8618bbba8cdSJeremy L Thompson 
8628bbba8cdSJeremy L Thompson       CeedCallBackend(CeedVectorTakeArray(impl->q_vecs_out[i], CEED_MEM_DEVICE, &e_vec_array));
8638bbba8cdSJeremy L Thompson       CeedCallBackend(CeedVectorRestoreArray(e_vec, &e_vec_array));
864756ca9e9SJeremy L Thompson     }
865756ca9e9SJeremy L Thompson 
866034f99fdSJeremy L Thompson     // Restrict
867034f99fdSJeremy L Thompson     if (impl->skip_rstr_out[field]) continue;
868034f99fdSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[field], &elem_rstr));
869034f99fdSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionApply(elem_rstr, CEED_TRANSPOSE, e_vec, l_vec, request));
870756ca9e9SJeremy L Thompson   }
8718bbba8cdSJeremy L Thompson 
8728bbba8cdSJeremy L Thompson   // Restore work vector
8738bbba8cdSJeremy L Thompson   CeedCallBackend(CeedRestoreWorkVector(ceed, &active_e_vec));
874756ca9e9SJeremy L Thompson   return CEED_ERROR_SUCCESS;
875756ca9e9SJeremy L Thompson }
876756ca9e9SJeremy L Thompson 
877756ca9e9SJeremy L Thompson //------------------------------------------------------------------------------
878004e4986SSebastian Grimberg // Linear QFunction Assembly Core
8790d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
8802b730f8bSJeremy L Thompson static inline int CeedOperatorLinearAssembleQFunctionCore_Cuda(CeedOperator op, bool build_objects, CeedVector *assembled, CeedElemRestriction *rstr,
8810d0321e0SJeremy L Thompson                                                                CeedRequest *request) {
882ca735530SJeremy L Thompson   Ceed                ceed, ceed_parent;
883ca735530SJeremy L Thompson   CeedInt             num_active_in, num_active_out, Q, num_elem, num_input_fields, num_output_fields, size;
8848bbba8cdSJeremy L Thompson   CeedScalar         *assembled_array;
885ca735530SJeremy L Thompson   CeedVector         *active_inputs;
886ca735530SJeremy L Thompson   CeedQFunctionField *qf_input_fields, *qf_output_fields;
887ca735530SJeremy L Thompson   CeedQFunction       qf;
888ca735530SJeremy L Thompson   CeedOperatorField  *op_input_fields, *op_output_fields;
889ca735530SJeremy L Thompson   CeedOperator_Cuda  *impl;
890ca735530SJeremy L Thompson 
8912b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
892ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetFallbackParentCeed(op, &ceed_parent));
893e984cf9aSJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &impl));
894e984cf9aSJeremy L Thompson   CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q));
895ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem));
896e984cf9aSJeremy L Thompson   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
897ca735530SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
898ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
899ca735530SJeremy L Thompson   active_inputs = impl->qf_active_in;
900ca735530SJeremy L Thompson   num_active_in = impl->num_active_in, num_active_out = impl->num_active_out;
9010d0321e0SJeremy L Thompson 
9020d0321e0SJeremy L Thompson   // Setup
9032b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorSetup_Cuda(op));
9040d0321e0SJeremy L Thompson 
905034f99fdSJeremy L Thompson   // Process inputs
90643e13feeSJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
9078bbba8cdSJeremy L Thompson     CeedCallBackend(CeedOperatorInputRestrict_Cuda(op_input_fields[i], qf_input_fields[i], i, NULL, NULL, true, impl, request));
9088bbba8cdSJeremy L Thompson     CeedCallBackend(CeedOperatorInputBasis_Cuda(op_input_fields[i], qf_input_fields[i], i, NULL, NULL, num_elem, true, impl));
90943e13feeSJeremy L Thompson   }
9100d0321e0SJeremy L Thompson 
9110d0321e0SJeremy L Thompson   // Count number of active input fields
912ca735530SJeremy L Thompson   if (!num_active_in) {
913ca735530SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
914ca735530SJeremy L Thompson       CeedScalar *q_vec_array;
915034f99fdSJeremy L Thompson       CeedVector  l_vec;
916ca735530SJeremy L Thompson 
9170d0321e0SJeremy L Thompson       // Check if active input
918034f99fdSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &l_vec));
919034f99fdSJeremy L Thompson       if (l_vec == CEED_VECTOR_ACTIVE) {
920ca735530SJeremy L Thompson         CeedCallBackend(CeedQFunctionFieldGetSize(qf_input_fields[i], &size));
921ca735530SJeremy L Thompson         CeedCallBackend(CeedVectorSetValue(impl->q_vecs_in[i], 0.0));
922ca735530SJeremy L Thompson         CeedCallBackend(CeedVectorGetArray(impl->q_vecs_in[i], CEED_MEM_DEVICE, &q_vec_array));
923ca735530SJeremy L Thompson         CeedCallBackend(CeedRealloc(num_active_in + size, &active_inputs));
9240d0321e0SJeremy L Thompson         for (CeedInt field = 0; field < size; field++) {
925004e4986SSebastian Grimberg           CeedSize q_size = (CeedSize)Q * num_elem;
926004e4986SSebastian Grimberg 
927ca735530SJeremy L Thompson           CeedCallBackend(CeedVectorCreate(ceed, q_size, &active_inputs[num_active_in + field]));
928ca735530SJeremy L Thompson           CeedCallBackend(
929ca735530SJeremy L Thompson               CeedVectorSetArray(active_inputs[num_active_in + field], CEED_MEM_DEVICE, CEED_USE_POINTER, &q_vec_array[field * Q * num_elem]));
9300d0321e0SJeremy L Thompson         }
931ca735530SJeremy L Thompson         num_active_in += size;
932ca735530SJeremy L Thompson         CeedCallBackend(CeedVectorRestoreArray(impl->q_vecs_in[i], &q_vec_array));
9330d0321e0SJeremy L Thompson       }
9340d0321e0SJeremy L Thompson     }
935ca735530SJeremy L Thompson     impl->num_active_in = num_active_in;
936ca735530SJeremy L Thompson     impl->qf_active_in  = active_inputs;
9370d0321e0SJeremy L Thompson   }
9380d0321e0SJeremy L Thompson 
9390d0321e0SJeremy L Thompson   // Count number of active output fields
940ca735530SJeremy L Thompson   if (!num_active_out) {
941ca735530SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
942034f99fdSJeremy L Thompson       CeedVector l_vec;
943ca735530SJeremy L Thompson 
9440d0321e0SJeremy L Thompson       // Check if active output
945034f99fdSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[i], &l_vec));
946034f99fdSJeremy L Thompson       if (l_vec == CEED_VECTOR_ACTIVE) {
947ca735530SJeremy L Thompson         CeedCallBackend(CeedQFunctionFieldGetSize(qf_output_fields[i], &size));
948ca735530SJeremy L Thompson         num_active_out += size;
9490d0321e0SJeremy L Thompson       }
9500d0321e0SJeremy L Thompson     }
951ca735530SJeremy L Thompson     impl->num_active_out = num_active_out;
9520d0321e0SJeremy L Thompson   }
9530d0321e0SJeremy L Thompson 
9540d0321e0SJeremy L Thompson   // Check sizes
955ca735530SJeremy L Thompson   CeedCheck(num_active_in > 0 && num_active_out > 0, ceed, CEED_ERROR_BACKEND, "Cannot assemble QFunction without active inputs and outputs");
9560d0321e0SJeremy L Thompson 
9570d0321e0SJeremy L Thompson   // Build objects if needed
9580d0321e0SJeremy L Thompson   if (build_objects) {
959004e4986SSebastian Grimberg     CeedSize l_size     = (CeedSize)num_elem * Q * num_active_in * num_active_out;
960ca735530SJeremy L Thompson     CeedInt  strides[3] = {1, num_elem * Q, Q}; /* *NOPAD* */
961004e4986SSebastian Grimberg 
962004e4986SSebastian Grimberg     // Create output restriction
963ca735530SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionCreateStrided(ceed_parent, num_elem, Q, num_active_in * num_active_out,
9640a5597ceSJeremy L Thompson                                                      (CeedSize)num_active_in * (CeedSize)num_active_out * (CeedSize)num_elem * (CeedSize)Q, strides,
9650a5597ceSJeremy L Thompson                                                      rstr));
9660d0321e0SJeremy L Thompson     // Create assembled vector
967ca735530SJeremy L Thompson     CeedCallBackend(CeedVectorCreate(ceed_parent, l_size, assembled));
9680d0321e0SJeremy L Thompson   }
9692b730f8bSJeremy L Thompson   CeedCallBackend(CeedVectorSetValue(*assembled, 0.0));
970ca735530SJeremy L Thompson   CeedCallBackend(CeedVectorGetArray(*assembled, CEED_MEM_DEVICE, &assembled_array));
9710d0321e0SJeremy L Thompson 
9720d0321e0SJeremy L Thompson   // Assemble QFunction
973ca735530SJeremy L Thompson   for (CeedInt in = 0; in < num_active_in; in++) {
9740d0321e0SJeremy L Thompson     // Set Inputs
975ca735530SJeremy L Thompson     CeedCallBackend(CeedVectorSetValue(active_inputs[in], 1.0));
976ca735530SJeremy L Thompson     if (num_active_in > 1) {
977ca735530SJeremy L Thompson       CeedCallBackend(CeedVectorSetValue(active_inputs[(in + num_active_in - 1) % num_active_in], 0.0));
9780d0321e0SJeremy L Thompson     }
9790d0321e0SJeremy L Thompson     // Set Outputs
980ca735530SJeremy L Thompson     for (CeedInt out = 0; out < num_output_fields; out++) {
981034f99fdSJeremy L Thompson       CeedVector l_vec;
982ca735530SJeremy L Thompson 
9830d0321e0SJeremy L Thompson       // Get output vector
984034f99fdSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[out], &l_vec));
9850d0321e0SJeremy L Thompson       // Check if active output
986034f99fdSJeremy L Thompson       if (l_vec == CEED_VECTOR_ACTIVE) {
987ca735530SJeremy L Thompson         CeedCallBackend(CeedVectorSetArray(impl->q_vecs_out[out], CEED_MEM_DEVICE, CEED_USE_POINTER, assembled_array));
988ca735530SJeremy L Thompson         CeedCallBackend(CeedQFunctionFieldGetSize(qf_output_fields[out], &size));
989ca735530SJeremy L Thompson         assembled_array += size * Q * num_elem;  // Advance the pointer by the size of the output
9900d0321e0SJeremy L Thompson       }
9910d0321e0SJeremy L Thompson     }
9920d0321e0SJeremy L Thompson     // Apply QFunction
993ca735530SJeremy L Thompson     CeedCallBackend(CeedQFunctionApply(qf, Q * num_elem, impl->q_vecs_in, impl->q_vecs_out));
9940d0321e0SJeremy L Thompson   }
9950d0321e0SJeremy L Thompson 
9968a213570SJeremy L Thompson   // Un-set output q-vecs to prevent accidental overwrite of Assembled
997ca735530SJeremy L Thompson   for (CeedInt out = 0; out < num_output_fields; out++) {
998034f99fdSJeremy L Thompson     CeedVector l_vec;
999ca735530SJeremy L Thompson 
1000034f99fdSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[out], &l_vec));
1001034f99fdSJeremy L Thompson     if (l_vec == CEED_VECTOR_ACTIVE) {
1002ca735530SJeremy L Thompson       CeedCallBackend(CeedVectorTakeArray(impl->q_vecs_out[out], CEED_MEM_DEVICE, NULL));
10030d0321e0SJeremy L Thompson     }
10040d0321e0SJeremy L Thompson   }
10050d0321e0SJeremy L Thompson 
10060d0321e0SJeremy L Thompson   // Restore input arrays
100743e13feeSJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
10088bbba8cdSJeremy L Thompson     CeedCallBackend(CeedOperatorInputRestore_Cuda(op_input_fields[i], qf_input_fields[i], i, NULL, NULL, true, impl));
100943e13feeSJeremy L Thompson   }
10100d0321e0SJeremy L Thompson 
10110d0321e0SJeremy L Thompson   // Restore output
1012ca735530SJeremy L Thompson   CeedCallBackend(CeedVectorRestoreArray(*assembled, &assembled_array));
10130d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
10140d0321e0SJeremy L Thompson }
10150d0321e0SJeremy L Thompson 
10160d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
10170d0321e0SJeremy L Thompson // Assemble Linear QFunction
10180d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
10192b730f8bSJeremy L Thompson static int CeedOperatorLinearAssembleQFunction_Cuda(CeedOperator op, CeedVector *assembled, CeedElemRestriction *rstr, CeedRequest *request) {
10202b730f8bSJeremy L Thompson   return CeedOperatorLinearAssembleQFunctionCore_Cuda(op, true, assembled, rstr, request);
10210d0321e0SJeremy L Thompson }
10220d0321e0SJeremy L Thompson 
10230d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
10240d0321e0SJeremy L Thompson // Update Assembled Linear QFunction
10250d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
10262b730f8bSJeremy L Thompson static int CeedOperatorLinearAssembleQFunctionUpdate_Cuda(CeedOperator op, CeedVector assembled, CeedElemRestriction rstr, CeedRequest *request) {
10272b730f8bSJeremy L Thompson   return CeedOperatorLinearAssembleQFunctionCore_Cuda(op, false, &assembled, &rstr, request);
10280d0321e0SJeremy L Thompson }
10290d0321e0SJeremy L Thompson 
10300d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
1031004e4986SSebastian Grimberg // Assemble Diagonal Setup
10320d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
1033cbfe683aSSebastian Grimberg static inline int CeedOperatorAssembleDiagonalSetup_Cuda(CeedOperator op) {
10340d0321e0SJeremy L Thompson   Ceed                ceed;
1035004e4986SSebastian Grimberg   CeedInt             num_input_fields, num_output_fields, num_eval_modes_in = 0, num_eval_modes_out = 0;
1036cbfe683aSSebastian Grimberg   CeedInt             q_comp, num_nodes, num_qpts;
1037004e4986SSebastian Grimberg   CeedEvalMode       *eval_modes_in = NULL, *eval_modes_out = NULL;
1038ca735530SJeremy L Thompson   CeedBasis           basis_in = NULL, basis_out = NULL;
1039ca735530SJeremy L Thompson   CeedQFunctionField *qf_fields;
10400d0321e0SJeremy L Thompson   CeedQFunction       qf;
1041ca735530SJeremy L Thompson   CeedOperatorField  *op_fields;
1042ca735530SJeremy L Thompson   CeedOperator_Cuda  *impl;
1043ca735530SJeremy L Thompson 
1044ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
10452b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
1046ca735530SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetNumArgs(qf, &num_input_fields, &num_output_fields));
10470d0321e0SJeremy L Thompson 
10480d0321e0SJeremy L Thompson   // Determine active input basis
1049ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetFields(op, NULL, &op_fields, NULL, NULL));
1050ca735530SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_fields, NULL, NULL));
1051ca735530SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
10520d0321e0SJeremy L Thompson     CeedVector vec;
1053ca735530SJeremy L Thompson 
1054ca735530SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec));
10550d0321e0SJeremy L Thompson     if (vec == CEED_VECTOR_ACTIVE) {
1056004e4986SSebastian Grimberg       CeedBasis    basis;
1057004e4986SSebastian Grimberg       CeedEvalMode eval_mode;
1058ca735530SJeremy L Thompson 
1059004e4986SSebastian Grimberg       CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis));
1060004e4986SSebastian Grimberg       CeedCheck(!basis_in || basis_in == basis, ceed, CEED_ERROR_BACKEND,
1061004e4986SSebastian Grimberg                 "Backend does not implement operator diagonal assembly with multiple active bases");
1062004e4986SSebastian Grimberg       basis_in = basis;
1063004e4986SSebastian Grimberg       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode));
1064004e4986SSebastian Grimberg       CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_in, eval_mode, &q_comp));
1065004e4986SSebastian Grimberg       if (eval_mode != CEED_EVAL_WEIGHT) {
1066004e4986SSebastian Grimberg         // q_comp = 1 if CEED_EVAL_NONE, CEED_EVAL_WEIGHT caught by QF assembly
1067004e4986SSebastian Grimberg         CeedCallBackend(CeedRealloc(num_eval_modes_in + q_comp, &eval_modes_in));
1068004e4986SSebastian Grimberg         for (CeedInt d = 0; d < q_comp; d++) eval_modes_in[num_eval_modes_in + d] = eval_mode;
1069004e4986SSebastian Grimberg         num_eval_modes_in += q_comp;
10700d0321e0SJeremy L Thompson       }
10710d0321e0SJeremy L Thompson     }
10720d0321e0SJeremy L Thompson   }
10730d0321e0SJeremy L Thompson 
10740d0321e0SJeremy L Thompson   // Determine active output basis
1075ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetFields(op, NULL, NULL, NULL, &op_fields));
1076ca735530SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qf_fields));
1077ca735530SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
10780d0321e0SJeremy L Thompson     CeedVector vec;
1079ca735530SJeremy L Thompson 
1080ca735530SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec));
10810d0321e0SJeremy L Thompson     if (vec == CEED_VECTOR_ACTIVE) {
1082004e4986SSebastian Grimberg       CeedBasis    basis;
1083004e4986SSebastian Grimberg       CeedEvalMode eval_mode;
1084ca735530SJeremy L Thompson 
1085004e4986SSebastian Grimberg       CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis));
1086004e4986SSebastian Grimberg       CeedCheck(!basis_out || basis_out == basis, ceed, CEED_ERROR_BACKEND,
1087004e4986SSebastian Grimberg                 "Backend does not implement operator diagonal assembly with multiple active bases");
1088004e4986SSebastian Grimberg       basis_out = basis;
1089004e4986SSebastian Grimberg       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode));
1090004e4986SSebastian Grimberg       CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_out, eval_mode, &q_comp));
1091004e4986SSebastian Grimberg       if (eval_mode != CEED_EVAL_WEIGHT) {
1092004e4986SSebastian Grimberg         // q_comp = 1 if CEED_EVAL_NONE, CEED_EVAL_WEIGHT caught by QF assembly
1093004e4986SSebastian Grimberg         CeedCallBackend(CeedRealloc(num_eval_modes_out + q_comp, &eval_modes_out));
1094004e4986SSebastian Grimberg         for (CeedInt d = 0; d < q_comp; d++) eval_modes_out[num_eval_modes_out + d] = eval_mode;
1095004e4986SSebastian Grimberg         num_eval_modes_out += q_comp;
10960d0321e0SJeremy L Thompson       }
10970d0321e0SJeremy L Thompson     }
10980d0321e0SJeremy L Thompson   }
10990d0321e0SJeremy L Thompson 
11000d0321e0SJeremy L Thompson   // Operator data struct
11012b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &impl));
11022b730f8bSJeremy L Thompson   CeedCallBackend(CeedCalloc(1, &impl->diag));
11030d0321e0SJeremy L Thompson   CeedOperatorDiag_Cuda *diag = impl->diag;
1104ca735530SJeremy L Thompson 
1105cbfe683aSSebastian Grimberg   // Basis matrices
1106004e4986SSebastian Grimberg   CeedCallBackend(CeedBasisGetNumNodes(basis_in, &num_nodes));
1107004e4986SSebastian Grimberg   if (basis_in == CEED_BASIS_NONE) num_qpts = num_nodes;
1108004e4986SSebastian Grimberg   else CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis_in, &num_qpts));
1109004e4986SSebastian Grimberg   const CeedInt interp_bytes     = num_nodes * num_qpts * sizeof(CeedScalar);
1110004e4986SSebastian Grimberg   const CeedInt eval_modes_bytes = sizeof(CeedEvalMode);
1111004e4986SSebastian Grimberg   bool          has_eval_none    = false;
11120d0321e0SJeremy L Thompson 
11130d0321e0SJeremy L Thompson   // CEED_EVAL_NONE
1114004e4986SSebastian Grimberg   for (CeedInt i = 0; i < num_eval_modes_in; i++) has_eval_none = has_eval_none || (eval_modes_in[i] == CEED_EVAL_NONE);
1115004e4986SSebastian Grimberg   for (CeedInt i = 0; i < num_eval_modes_out; i++) has_eval_none = has_eval_none || (eval_modes_out[i] == CEED_EVAL_NONE);
1116004e4986SSebastian Grimberg   if (has_eval_none) {
11170d0321e0SJeremy L Thompson     CeedScalar *identity = NULL;
1118ca735530SJeremy L Thompson 
1119004e4986SSebastian Grimberg     CeedCallBackend(CeedCalloc(num_nodes * num_qpts, &identity));
1120ca735530SJeremy L Thompson     for (CeedInt i = 0; i < (num_nodes < num_qpts ? num_nodes : num_qpts); i++) identity[i * num_nodes + i] = 1.0;
1121ca735530SJeremy L Thompson     CeedCallCuda(ceed, cudaMalloc((void **)&diag->d_identity, interp_bytes));
1122ca735530SJeremy L Thompson     CeedCallCuda(ceed, cudaMemcpy(diag->d_identity, identity, interp_bytes, cudaMemcpyHostToDevice));
1123004e4986SSebastian Grimberg     CeedCallBackend(CeedFree(&identity));
11240d0321e0SJeremy L Thompson   }
11250d0321e0SJeremy L Thompson 
1126004e4986SSebastian Grimberg   // CEED_EVAL_INTERP, CEED_EVAL_GRAD, CEED_EVAL_DIV, and CEED_EVAL_CURL
1127004e4986SSebastian Grimberg   for (CeedInt in = 0; in < 2; in++) {
1128004e4986SSebastian Grimberg     CeedFESpace fespace;
1129004e4986SSebastian Grimberg     CeedBasis   basis = in ? basis_in : basis_out;
11300d0321e0SJeremy L Thompson 
1131004e4986SSebastian Grimberg     CeedCallBackend(CeedBasisGetFESpace(basis, &fespace));
1132004e4986SSebastian Grimberg     switch (fespace) {
1133004e4986SSebastian Grimberg       case CEED_FE_SPACE_H1: {
1134004e4986SSebastian Grimberg         CeedInt           q_comp_interp, q_comp_grad;
1135004e4986SSebastian Grimberg         const CeedScalar *interp, *grad;
1136004e4986SSebastian Grimberg         CeedScalar       *d_interp, *d_grad;
11370d0321e0SJeremy L Thompson 
1138004e4986SSebastian Grimberg         CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp));
1139004e4986SSebastian Grimberg         CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_GRAD, &q_comp_grad));
11400d0321e0SJeremy L Thompson 
1141004e4986SSebastian Grimberg         CeedCallBackend(CeedBasisGetInterp(basis, &interp));
1142004e4986SSebastian Grimberg         CeedCallCuda(ceed, cudaMalloc((void **)&d_interp, interp_bytes * q_comp_interp));
1143004e4986SSebastian Grimberg         CeedCallCuda(ceed, cudaMemcpy(d_interp, interp, interp_bytes * q_comp_interp, cudaMemcpyHostToDevice));
1144004e4986SSebastian Grimberg         CeedCallBackend(CeedBasisGetGrad(basis, &grad));
1145004e4986SSebastian Grimberg         CeedCallCuda(ceed, cudaMalloc((void **)&d_grad, interp_bytes * q_comp_grad));
1146004e4986SSebastian Grimberg         CeedCallCuda(ceed, cudaMemcpy(d_grad, grad, interp_bytes * q_comp_grad, cudaMemcpyHostToDevice));
1147004e4986SSebastian Grimberg         if (in) {
1148004e4986SSebastian Grimberg           diag->d_interp_in = d_interp;
1149004e4986SSebastian Grimberg           diag->d_grad_in   = d_grad;
1150004e4986SSebastian Grimberg         } else {
1151004e4986SSebastian Grimberg           diag->d_interp_out = d_interp;
1152004e4986SSebastian Grimberg           diag->d_grad_out   = d_grad;
1153004e4986SSebastian Grimberg         }
1154004e4986SSebastian Grimberg       } break;
1155004e4986SSebastian Grimberg       case CEED_FE_SPACE_HDIV: {
1156004e4986SSebastian Grimberg         CeedInt           q_comp_interp, q_comp_div;
1157004e4986SSebastian Grimberg         const CeedScalar *interp, *div;
1158004e4986SSebastian Grimberg         CeedScalar       *d_interp, *d_div;
1159004e4986SSebastian Grimberg 
1160004e4986SSebastian Grimberg         CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp));
1161004e4986SSebastian Grimberg         CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_DIV, &q_comp_div));
1162004e4986SSebastian Grimberg 
1163004e4986SSebastian Grimberg         CeedCallBackend(CeedBasisGetInterp(basis, &interp));
1164004e4986SSebastian Grimberg         CeedCallCuda(ceed, cudaMalloc((void **)&d_interp, interp_bytes * q_comp_interp));
1165004e4986SSebastian Grimberg         CeedCallCuda(ceed, cudaMemcpy(d_interp, interp, interp_bytes * q_comp_interp, cudaMemcpyHostToDevice));
1166004e4986SSebastian Grimberg         CeedCallBackend(CeedBasisGetDiv(basis, &div));
1167004e4986SSebastian Grimberg         CeedCallCuda(ceed, cudaMalloc((void **)&d_div, interp_bytes * q_comp_div));
1168004e4986SSebastian Grimberg         CeedCallCuda(ceed, cudaMemcpy(d_div, div, interp_bytes * q_comp_div, cudaMemcpyHostToDevice));
1169004e4986SSebastian Grimberg         if (in) {
1170004e4986SSebastian Grimberg           diag->d_interp_in = d_interp;
1171004e4986SSebastian Grimberg           diag->d_div_in    = d_div;
1172004e4986SSebastian Grimberg         } else {
1173004e4986SSebastian Grimberg           diag->d_interp_out = d_interp;
1174004e4986SSebastian Grimberg           diag->d_div_out    = d_div;
1175004e4986SSebastian Grimberg         }
1176004e4986SSebastian Grimberg       } break;
1177004e4986SSebastian Grimberg       case CEED_FE_SPACE_HCURL: {
1178004e4986SSebastian Grimberg         CeedInt           q_comp_interp, q_comp_curl;
1179004e4986SSebastian Grimberg         const CeedScalar *interp, *curl;
1180004e4986SSebastian Grimberg         CeedScalar       *d_interp, *d_curl;
1181004e4986SSebastian Grimberg 
1182004e4986SSebastian Grimberg         CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp));
1183004e4986SSebastian Grimberg         CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_CURL, &q_comp_curl));
1184004e4986SSebastian Grimberg 
1185004e4986SSebastian Grimberg         CeedCallBackend(CeedBasisGetInterp(basis, &interp));
1186004e4986SSebastian Grimberg         CeedCallCuda(ceed, cudaMalloc((void **)&d_interp, interp_bytes * q_comp_interp));
1187004e4986SSebastian Grimberg         CeedCallCuda(ceed, cudaMemcpy(d_interp, interp, interp_bytes * q_comp_interp, cudaMemcpyHostToDevice));
1188004e4986SSebastian Grimberg         CeedCallBackend(CeedBasisGetCurl(basis, &curl));
1189004e4986SSebastian Grimberg         CeedCallCuda(ceed, cudaMalloc((void **)&d_curl, interp_bytes * q_comp_curl));
1190004e4986SSebastian Grimberg         CeedCallCuda(ceed, cudaMemcpy(d_curl, curl, interp_bytes * q_comp_curl, cudaMemcpyHostToDevice));
1191004e4986SSebastian Grimberg         if (in) {
1192004e4986SSebastian Grimberg           diag->d_interp_in = d_interp;
1193004e4986SSebastian Grimberg           diag->d_curl_in   = d_curl;
1194004e4986SSebastian Grimberg         } else {
1195004e4986SSebastian Grimberg           diag->d_interp_out = d_interp;
1196004e4986SSebastian Grimberg           diag->d_curl_out   = d_curl;
1197004e4986SSebastian Grimberg         }
1198004e4986SSebastian Grimberg       } break;
1199004e4986SSebastian Grimberg     }
1200004e4986SSebastian Grimberg   }
1201004e4986SSebastian Grimberg 
1202004e4986SSebastian Grimberg   // Arrays of eval_modes
1203004e4986SSebastian Grimberg   CeedCallCuda(ceed, cudaMalloc((void **)&diag->d_eval_modes_in, num_eval_modes_in * eval_modes_bytes));
1204004e4986SSebastian Grimberg   CeedCallCuda(ceed, cudaMemcpy(diag->d_eval_modes_in, eval_modes_in, num_eval_modes_in * eval_modes_bytes, cudaMemcpyHostToDevice));
1205004e4986SSebastian Grimberg   CeedCallCuda(ceed, cudaMalloc((void **)&diag->d_eval_modes_out, num_eval_modes_out * eval_modes_bytes));
1206004e4986SSebastian Grimberg   CeedCallCuda(ceed, cudaMemcpy(diag->d_eval_modes_out, eval_modes_out, num_eval_modes_out * eval_modes_bytes, cudaMemcpyHostToDevice));
1207004e4986SSebastian Grimberg   CeedCallBackend(CeedFree(&eval_modes_in));
1208004e4986SSebastian Grimberg   CeedCallBackend(CeedFree(&eval_modes_out));
12090d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
12100d0321e0SJeremy L Thompson }
12110d0321e0SJeremy L Thompson 
12120d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
1213cbfe683aSSebastian Grimberg // Assemble Diagonal Setup (Compilation)
1214cbfe683aSSebastian Grimberg //------------------------------------------------------------------------------
1215cbfe683aSSebastian Grimberg static inline int CeedOperatorAssembleDiagonalSetupCompile_Cuda(CeedOperator op, CeedInt use_ceedsize_idx, const bool is_point_block) {
1216cbfe683aSSebastian Grimberg   Ceed                ceed;
121722070f95SJeremy L Thompson   char               *diagonal_kernel_source;
121822070f95SJeremy L Thompson   const char         *diagonal_kernel_path;
1219cbfe683aSSebastian Grimberg   CeedInt             num_input_fields, num_output_fields, num_eval_modes_in = 0, num_eval_modes_out = 0;
1220cbfe683aSSebastian Grimberg   CeedInt             num_comp, q_comp, num_nodes, num_qpts;
1221cbfe683aSSebastian Grimberg   CeedBasis           basis_in = NULL, basis_out = NULL;
1222cbfe683aSSebastian Grimberg   CeedQFunctionField *qf_fields;
1223cbfe683aSSebastian Grimberg   CeedQFunction       qf;
1224cbfe683aSSebastian Grimberg   CeedOperatorField  *op_fields;
1225cbfe683aSSebastian Grimberg   CeedOperator_Cuda  *impl;
1226cbfe683aSSebastian Grimberg 
1227cbfe683aSSebastian Grimberg   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
1228cbfe683aSSebastian Grimberg   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
1229cbfe683aSSebastian Grimberg   CeedCallBackend(CeedQFunctionGetNumArgs(qf, &num_input_fields, &num_output_fields));
1230cbfe683aSSebastian Grimberg 
1231cbfe683aSSebastian Grimberg   // Determine active input basis
1232cbfe683aSSebastian Grimberg   CeedCallBackend(CeedOperatorGetFields(op, NULL, &op_fields, NULL, NULL));
1233cbfe683aSSebastian Grimberg   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_fields, NULL, NULL));
1234cbfe683aSSebastian Grimberg   for (CeedInt i = 0; i < num_input_fields; i++) {
1235cbfe683aSSebastian Grimberg     CeedVector vec;
1236cbfe683aSSebastian Grimberg 
1237cbfe683aSSebastian Grimberg     CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec));
1238cbfe683aSSebastian Grimberg     if (vec == CEED_VECTOR_ACTIVE) {
1239cbfe683aSSebastian Grimberg       CeedEvalMode eval_mode;
1240cbfe683aSSebastian Grimberg 
1241cbfe683aSSebastian Grimberg       CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis_in));
1242cbfe683aSSebastian Grimberg       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode));
1243cbfe683aSSebastian Grimberg       CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_in, eval_mode, &q_comp));
1244cbfe683aSSebastian Grimberg       if (eval_mode != CEED_EVAL_WEIGHT) {
1245cbfe683aSSebastian Grimberg         num_eval_modes_in += q_comp;
1246cbfe683aSSebastian Grimberg       }
1247cbfe683aSSebastian Grimberg     }
1248cbfe683aSSebastian Grimberg   }
1249cbfe683aSSebastian Grimberg 
1250cbfe683aSSebastian Grimberg   // Determine active output basis
1251cbfe683aSSebastian Grimberg   CeedCallBackend(CeedOperatorGetFields(op, NULL, NULL, NULL, &op_fields));
1252cbfe683aSSebastian Grimberg   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qf_fields));
1253cbfe683aSSebastian Grimberg   for (CeedInt i = 0; i < num_output_fields; i++) {
1254cbfe683aSSebastian Grimberg     CeedVector vec;
1255cbfe683aSSebastian Grimberg 
1256cbfe683aSSebastian Grimberg     CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec));
1257cbfe683aSSebastian Grimberg     if (vec == CEED_VECTOR_ACTIVE) {
1258cbfe683aSSebastian Grimberg       CeedEvalMode eval_mode;
1259cbfe683aSSebastian Grimberg 
1260cbfe683aSSebastian Grimberg       CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis_out));
1261cbfe683aSSebastian Grimberg       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode));
1262cbfe683aSSebastian Grimberg       CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_out, eval_mode, &q_comp));
1263cbfe683aSSebastian Grimberg       if (eval_mode != CEED_EVAL_WEIGHT) {
1264cbfe683aSSebastian Grimberg         num_eval_modes_out += q_comp;
1265cbfe683aSSebastian Grimberg       }
1266cbfe683aSSebastian Grimberg     }
1267cbfe683aSSebastian Grimberg   }
1268cbfe683aSSebastian Grimberg 
1269cbfe683aSSebastian Grimberg   // Operator data struct
1270cbfe683aSSebastian Grimberg   CeedCallBackend(CeedOperatorGetData(op, &impl));
1271cbfe683aSSebastian Grimberg   CeedOperatorDiag_Cuda *diag = impl->diag;
1272cbfe683aSSebastian Grimberg 
1273cbfe683aSSebastian Grimberg   // Assemble kernel
1274cbfe683aSSebastian Grimberg   CUmodule *module          = is_point_block ? &diag->module_point_block : &diag->module;
1275cbfe683aSSebastian Grimberg   CeedInt   elems_per_block = 1;
1276cbfe683aSSebastian Grimberg   CeedCallBackend(CeedBasisGetNumNodes(basis_in, &num_nodes));
1277cbfe683aSSebastian Grimberg   CeedCallBackend(CeedBasisGetNumComponents(basis_in, &num_comp));
1278cbfe683aSSebastian Grimberg   if (basis_in == CEED_BASIS_NONE) num_qpts = num_nodes;
1279cbfe683aSSebastian Grimberg   else CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis_in, &num_qpts));
1280cbfe683aSSebastian Grimberg   CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-operator-assemble-diagonal.h", &diagonal_kernel_path));
1281cbfe683aSSebastian Grimberg   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Diagonal Assembly Kernel Source -----\n");
1282cbfe683aSSebastian Grimberg   CeedCallBackend(CeedLoadSourceToBuffer(ceed, diagonal_kernel_path, &diagonal_kernel_source));
1283cbfe683aSSebastian Grimberg   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Diagonal Assembly Source Complete! -----\n");
1284cbfe683aSSebastian Grimberg   CeedCallCuda(ceed, CeedCompile_Cuda(ceed, diagonal_kernel_source, module, 8, "NUM_EVAL_MODES_IN", num_eval_modes_in, "NUM_EVAL_MODES_OUT",
1285cbfe683aSSebastian Grimberg                                       num_eval_modes_out, "NUM_COMP", num_comp, "NUM_NODES", num_nodes, "NUM_QPTS", num_qpts, "USE_CEEDSIZE",
1286cbfe683aSSebastian Grimberg                                       use_ceedsize_idx, "USE_POINT_BLOCK", is_point_block ? 1 : 0, "BLOCK_SIZE", num_nodes * elems_per_block));
1287cbfe683aSSebastian Grimberg   CeedCallCuda(ceed, CeedGetKernel_Cuda(ceed, *module, "LinearDiagonal", is_point_block ? &diag->LinearPointBlock : &diag->LinearDiagonal));
1288cbfe683aSSebastian Grimberg   CeedCallBackend(CeedFree(&diagonal_kernel_path));
1289cbfe683aSSebastian Grimberg   CeedCallBackend(CeedFree(&diagonal_kernel_source));
1290cbfe683aSSebastian Grimberg   return CEED_ERROR_SUCCESS;
1291cbfe683aSSebastian Grimberg }
1292cbfe683aSSebastian Grimberg 
1293cbfe683aSSebastian Grimberg //------------------------------------------------------------------------------
1294004e4986SSebastian Grimberg // Assemble Diagonal Core
12950d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
1296ca735530SJeremy L Thompson static inline int CeedOperatorAssembleDiagonalCore_Cuda(CeedOperator op, CeedVector assembled, CeedRequest *request, const bool is_point_block) {
12970d0321e0SJeremy L Thompson   Ceed                ceed;
1298cbfe683aSSebastian Grimberg   CeedInt             num_elem, num_nodes;
1299ca735530SJeremy L Thompson   CeedScalar         *elem_diag_array;
1300ca735530SJeremy L Thompson   const CeedScalar   *assembled_qf_array;
1301004e4986SSebastian Grimberg   CeedVector          assembled_qf   = NULL, elem_diag;
1302004e4986SSebastian Grimberg   CeedElemRestriction assembled_rstr = NULL, rstr_in, rstr_out, diag_rstr;
13030d0321e0SJeremy L Thompson   CeedOperator_Cuda  *impl;
1304ca735530SJeremy L Thompson 
1305ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
13062b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &impl));
13070d0321e0SJeremy L Thompson 
13080d0321e0SJeremy L Thompson   // Assemble QFunction
1309004e4986SSebastian Grimberg   CeedCallBackend(CeedOperatorLinearAssembleQFunctionBuildOrUpdate(op, &assembled_qf, &assembled_rstr, request));
1310004e4986SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionDestroy(&assembled_rstr));
1311004e4986SSebastian Grimberg   CeedCallBackend(CeedVectorGetArrayRead(assembled_qf, CEED_MEM_DEVICE, &assembled_qf_array));
13120d0321e0SJeremy L Thompson 
1313cbfe683aSSebastian Grimberg   // Setup
1314cbfe683aSSebastian Grimberg   if (!impl->diag) CeedCallBackend(CeedOperatorAssembleDiagonalSetup_Cuda(op));
1315cbfe683aSSebastian Grimberg   CeedOperatorDiag_Cuda *diag = impl->diag;
1316cbfe683aSSebastian Grimberg 
1317cbfe683aSSebastian Grimberg   assert(diag != NULL);
1318cbfe683aSSebastian Grimberg 
1319cbfe683aSSebastian Grimberg   // Assemble kernel if needed
1320cbfe683aSSebastian Grimberg   if ((!is_point_block && !diag->LinearDiagonal) || (is_point_block && !diag->LinearPointBlock)) {
1321cbfe683aSSebastian Grimberg     CeedSize assembled_length, assembled_qf_length;
1322cbfe683aSSebastian Grimberg     CeedInt  use_ceedsize_idx = 0;
1323f7c1b517Snbeams     CeedCallBackend(CeedVectorGetLength(assembled, &assembled_length));
1324ca735530SJeremy L Thompson     CeedCallBackend(CeedVectorGetLength(assembled_qf, &assembled_qf_length));
1325ca735530SJeremy L Thompson     if ((assembled_length > INT_MAX) || (assembled_qf_length > INT_MAX)) use_ceedsize_idx = 1;
1326f7c1b517Snbeams 
1327cbfe683aSSebastian Grimberg     CeedCallBackend(CeedOperatorAssembleDiagonalSetupCompile_Cuda(op, use_ceedsize_idx, is_point_block));
1328cbfe683aSSebastian Grimberg   }
13290d0321e0SJeremy L Thompson 
1330004e4986SSebastian Grimberg   // Restriction and diagonal vector
1331004e4986SSebastian Grimberg   CeedCallBackend(CeedOperatorGetActiveElemRestrictions(op, &rstr_in, &rstr_out));
1332004e4986SSebastian Grimberg   CeedCheck(rstr_in == rstr_out, ceed, CEED_ERROR_BACKEND,
1333004e4986SSebastian Grimberg             "Cannot assemble operator diagonal with different input and output active element restrictions");
1334004e4986SSebastian Grimberg   if (!is_point_block && !diag->diag_rstr) {
1335004e4986SSebastian Grimberg     CeedCallBackend(CeedElemRestrictionCreateUnsignedCopy(rstr_out, &diag->diag_rstr));
1336004e4986SSebastian Grimberg     CeedCallBackend(CeedElemRestrictionCreateVector(diag->diag_rstr, NULL, &diag->elem_diag));
1337004e4986SSebastian Grimberg   } else if (is_point_block && !diag->point_block_diag_rstr) {
1338004e4986SSebastian Grimberg     CeedCallBackend(CeedOperatorCreateActivePointBlockRestriction(rstr_out, &diag->point_block_diag_rstr));
1339004e4986SSebastian Grimberg     CeedCallBackend(CeedElemRestrictionCreateVector(diag->point_block_diag_rstr, NULL, &diag->point_block_elem_diag));
13400d0321e0SJeremy L Thompson   }
1341004e4986SSebastian Grimberg   diag_rstr = is_point_block ? diag->point_block_diag_rstr : diag->diag_rstr;
1342004e4986SSebastian Grimberg   elem_diag = is_point_block ? diag->point_block_elem_diag : diag->elem_diag;
1343ca735530SJeremy L Thompson   CeedCallBackend(CeedVectorSetValue(elem_diag, 0.0));
13440d0321e0SJeremy L Thompson 
134591db28b6SZach Atkins   // Only assemble diagonal if the basis has nodes, otherwise inputs are null pointers
1346004e4986SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionGetElementSize(diag_rstr, &num_nodes));
1347004e4986SSebastian Grimberg   if (num_nodes > 0) {
13480d0321e0SJeremy L Thompson     // Assemble element operator diagonals
1349ca735530SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumElements(diag_rstr, &num_elem));
1350004e4986SSebastian Grimberg     CeedCallBackend(CeedVectorGetArray(elem_diag, CEED_MEM_DEVICE, &elem_diag_array));
13510d0321e0SJeremy L Thompson 
13520d0321e0SJeremy L Thompson     // Compute the diagonal of B^T D B
1353004e4986SSebastian Grimberg     CeedInt elems_per_block = 1;
1354004e4986SSebastian Grimberg     CeedInt grid            = CeedDivUpInt(num_elem, elems_per_block);
1355004e4986SSebastian Grimberg     void   *args[]          = {(void *)&num_elem,      &diag->d_identity,       &diag->d_interp_in,  &diag->d_grad_in, &diag->d_div_in,
1356004e4986SSebastian Grimberg                                &diag->d_curl_in,       &diag->d_interp_out,     &diag->d_grad_out,   &diag->d_div_out, &diag->d_curl_out,
1357004e4986SSebastian Grimberg                                &diag->d_eval_modes_in, &diag->d_eval_modes_out, &assembled_qf_array, &elem_diag_array};
1358004e4986SSebastian Grimberg 
1359ca735530SJeremy L Thompson     if (is_point_block) {
1360004e4986SSebastian Grimberg       CeedCallBackend(CeedRunKernelDim_Cuda(ceed, diag->LinearPointBlock, grid, num_nodes, 1, elems_per_block, args));
13610d0321e0SJeremy L Thompson     } else {
1362004e4986SSebastian Grimberg       CeedCallBackend(CeedRunKernelDim_Cuda(ceed, diag->LinearDiagonal, grid, num_nodes, 1, elems_per_block, args));
13630d0321e0SJeremy L Thompson     }
13640d0321e0SJeremy L Thompson 
13650d0321e0SJeremy L Thompson     // Restore arrays
1366ca735530SJeremy L Thompson     CeedCallBackend(CeedVectorRestoreArray(elem_diag, &elem_diag_array));
1367ca735530SJeremy L Thompson     CeedCallBackend(CeedVectorRestoreArrayRead(assembled_qf, &assembled_qf_array));
136891db28b6SZach Atkins   }
13690d0321e0SJeremy L Thompson 
13700d0321e0SJeremy L Thompson   // Assemble local operator diagonal
1371ca735530SJeremy L Thompson   CeedCallBackend(CeedElemRestrictionApply(diag_rstr, CEED_TRANSPOSE, elem_diag, assembled, request));
13720d0321e0SJeremy L Thompson 
13730d0321e0SJeremy L Thompson   // Cleanup
1374ca735530SJeremy L Thompson   CeedCallBackend(CeedVectorDestroy(&assembled_qf));
13750d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
13760d0321e0SJeremy L Thompson }
13770d0321e0SJeremy L Thompson 
13780d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
13790d0321e0SJeremy L Thompson // Assemble Linear Diagonal
13800d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
13812b730f8bSJeremy L Thompson static int CeedOperatorLinearAssembleAddDiagonal_Cuda(CeedOperator op, CeedVector assembled, CeedRequest *request) {
13822b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorAssembleDiagonalCore_Cuda(op, assembled, request, false));
13836aa95790SJeremy L Thompson   return CEED_ERROR_SUCCESS;
13840d0321e0SJeremy L Thompson }
13850d0321e0SJeremy L Thompson 
13860d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
13870d0321e0SJeremy L Thompson // Assemble Linear Point Block Diagonal
13880d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
13892b730f8bSJeremy L Thompson static int CeedOperatorLinearAssembleAddPointBlockDiagonal_Cuda(CeedOperator op, CeedVector assembled, CeedRequest *request) {
13902b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorAssembleDiagonalCore_Cuda(op, assembled, request, true));
13916aa95790SJeremy L Thompson   return CEED_ERROR_SUCCESS;
13920d0321e0SJeremy L Thompson }
13930d0321e0SJeremy L Thompson 
13940d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
1395004e4986SSebastian Grimberg // Single Operator Assembly Setup
1396cc132f9aSnbeams //------------------------------------------------------------------------------
1397f7c1b517Snbeams static int CeedSingleOperatorAssembleSetup_Cuda(CeedOperator op, CeedInt use_ceedsize_idx) {
1398cc132f9aSnbeams   Ceed                ceed;
1399004e4986SSebastian Grimberg   Ceed_Cuda          *cuda_data;
140022070f95SJeremy L Thompson   char               *assembly_kernel_source;
140122070f95SJeremy L Thompson   const char         *assembly_kernel_path;
1402004e4986SSebastian Grimberg   CeedInt             num_input_fields, num_output_fields, num_eval_modes_in = 0, num_eval_modes_out = 0;
14033b38d1dfSJeremy L Thompson   CeedInt             elem_size_in, num_qpts_in = 0, num_comp_in, elem_size_out, num_qpts_out, num_comp_out, q_comp;
1404004e4986SSebastian Grimberg   CeedEvalMode       *eval_modes_in = NULL, *eval_modes_out = NULL;
1405ca735530SJeremy L Thompson   CeedElemRestriction rstr_in = NULL, rstr_out = NULL;
1406ca735530SJeremy L Thompson   CeedBasis           basis_in = NULL, basis_out = NULL;
1407ca735530SJeremy L Thompson   CeedQFunctionField *qf_fields;
1408ca735530SJeremy L Thompson   CeedQFunction       qf;
1409ca735530SJeremy L Thompson   CeedOperatorField  *input_fields, *output_fields;
1410cc132f9aSnbeams   CeedOperator_Cuda  *impl;
1411ca735530SJeremy L Thompson 
1412ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
14132b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &impl));
1414cc132f9aSnbeams 
1415cc132f9aSnbeams   // Get intput and output fields
14162b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &input_fields, &num_output_fields, &output_fields));
1417cc132f9aSnbeams 
1418cc132f9aSnbeams   // Determine active input basis eval mode
14192b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
14202b730f8bSJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_fields, NULL, NULL));
1421cc132f9aSnbeams   for (CeedInt i = 0; i < num_input_fields; i++) {
1422cc132f9aSnbeams     CeedVector vec;
1423ca735530SJeremy L Thompson 
14242b730f8bSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetVector(input_fields[i], &vec));
1425cc132f9aSnbeams     if (vec == CEED_VECTOR_ACTIVE) {
1426004e4986SSebastian Grimberg       CeedBasis    basis;
1427ca735530SJeremy L Thompson       CeedEvalMode eval_mode;
1428ca735530SJeremy L Thompson 
1429004e4986SSebastian Grimberg       CeedCallBackend(CeedOperatorFieldGetBasis(input_fields[i], &basis));
1430004e4986SSebastian Grimberg       CeedCheck(!basis_in || basis_in == basis, ceed, CEED_ERROR_BACKEND, "Backend does not implement operator assembly with multiple active bases");
1431004e4986SSebastian Grimberg       basis_in = basis;
14322b730f8bSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetElemRestriction(input_fields[i], &rstr_in));
1433004e4986SSebastian Grimberg       CeedCallBackend(CeedElemRestrictionGetElementSize(rstr_in, &elem_size_in));
1434004e4986SSebastian Grimberg       if (basis_in == CEED_BASIS_NONE) num_qpts_in = elem_size_in;
1435004e4986SSebastian Grimberg       else CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis_in, &num_qpts_in));
14362b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode));
1437004e4986SSebastian Grimberg       CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_in, eval_mode, &q_comp));
1438004e4986SSebastian Grimberg       if (eval_mode != CEED_EVAL_WEIGHT) {
1439004e4986SSebastian Grimberg         // q_comp = 1 if CEED_EVAL_NONE, CEED_EVAL_WEIGHT caught by QF Assembly
1440004e4986SSebastian Grimberg         CeedCallBackend(CeedRealloc(num_eval_modes_in + q_comp, &eval_modes_in));
1441004e4986SSebastian Grimberg         for (CeedInt d = 0; d < q_comp; d++) {
1442004e4986SSebastian Grimberg           eval_modes_in[num_eval_modes_in + d] = eval_mode;
1443cc132f9aSnbeams         }
1444004e4986SSebastian Grimberg         num_eval_modes_in += q_comp;
1445cc132f9aSnbeams       }
1446cc132f9aSnbeams     }
1447cc132f9aSnbeams   }
1448cc132f9aSnbeams 
1449cc132f9aSnbeams   // Determine active output basis; basis_out and rstr_out only used if same as input, TODO
14502b730f8bSJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qf_fields));
1451cc132f9aSnbeams   for (CeedInt i = 0; i < num_output_fields; i++) {
1452cc132f9aSnbeams     CeedVector vec;
1453ca735530SJeremy L Thompson 
14542b730f8bSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetVector(output_fields[i], &vec));
1455cc132f9aSnbeams     if (vec == CEED_VECTOR_ACTIVE) {
1456004e4986SSebastian Grimberg       CeedBasis    basis;
1457ca735530SJeremy L Thompson       CeedEvalMode eval_mode;
1458ca735530SJeremy L Thompson 
1459004e4986SSebastian Grimberg       CeedCallBackend(CeedOperatorFieldGetBasis(output_fields[i], &basis));
1460004e4986SSebastian Grimberg       CeedCheck(!basis_out || basis_out == basis, ceed, CEED_ERROR_BACKEND,
1461004e4986SSebastian Grimberg                 "Backend does not implement operator assembly with multiple active bases");
1462004e4986SSebastian Grimberg       basis_out = basis;
14632b730f8bSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetElemRestriction(output_fields[i], &rstr_out));
1464004e4986SSebastian Grimberg       CeedCallBackend(CeedElemRestrictionGetElementSize(rstr_out, &elem_size_out));
1465004e4986SSebastian Grimberg       if (basis_out == CEED_BASIS_NONE) num_qpts_out = elem_size_out;
1466004e4986SSebastian Grimberg       else CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis_out, &num_qpts_out));
1467004e4986SSebastian Grimberg       CeedCheck(num_qpts_in == num_qpts_out, ceed, CEED_ERROR_UNSUPPORTED,
1468004e4986SSebastian Grimberg                 "Active input and output bases must have the same number of quadrature points");
14692b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode));
1470004e4986SSebastian Grimberg       CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_out, eval_mode, &q_comp));
1471004e4986SSebastian Grimberg       if (eval_mode != CEED_EVAL_WEIGHT) {
1472004e4986SSebastian Grimberg         // q_comp = 1 if CEED_EVAL_NONE, CEED_EVAL_WEIGHT caught by QF Assembly
1473004e4986SSebastian Grimberg         CeedCallBackend(CeedRealloc(num_eval_modes_out + q_comp, &eval_modes_out));
1474004e4986SSebastian Grimberg         for (CeedInt d = 0; d < q_comp; d++) {
1475004e4986SSebastian Grimberg           eval_modes_out[num_eval_modes_out + d] = eval_mode;
1476004e4986SSebastian Grimberg         }
1477004e4986SSebastian Grimberg         num_eval_modes_out += q_comp;
1478cc132f9aSnbeams       }
1479cc132f9aSnbeams     }
1480cc132f9aSnbeams   }
1481004e4986SSebastian Grimberg   CeedCheck(num_eval_modes_in > 0 && num_eval_modes_out > 0, ceed, CEED_ERROR_UNSUPPORTED, "Cannot assemble operator without inputs/outputs");
1482cc132f9aSnbeams 
14832b730f8bSJeremy L Thompson   CeedCallBackend(CeedCalloc(1, &impl->asmb));
1484cc132f9aSnbeams   CeedOperatorAssemble_Cuda *asmb = impl->asmb;
1485004e4986SSebastian Grimberg   asmb->elems_per_block           = 1;
1486004e4986SSebastian Grimberg   asmb->block_size_x              = elem_size_in;
1487004e4986SSebastian Grimberg   asmb->block_size_y              = elem_size_out;
1488ca735530SJeremy L Thompson 
14892b730f8bSJeremy L Thompson   CeedCallBackend(CeedGetData(ceed, &cuda_data));
1490004e4986SSebastian Grimberg   bool fallback = asmb->block_size_x * asmb->block_size_y * asmb->elems_per_block > cuda_data->device_prop.maxThreadsPerBlock;
1491004e4986SSebastian Grimberg 
1492004e4986SSebastian Grimberg   if (fallback) {
1493004e4986SSebastian Grimberg     // Use fallback kernel with 1D threadblock
1494004e4986SSebastian Grimberg     asmb->block_size_y = 1;
1495004e4986SSebastian Grimberg   }
1496004e4986SSebastian Grimberg 
1497004e4986SSebastian Grimberg   // Compile kernels
1498004e4986SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionGetNumComponents(rstr_in, &num_comp_in));
1499004e4986SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionGetNumComponents(rstr_out, &num_comp_out));
15002b730f8bSJeremy L Thompson   CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-operator-assemble.h", &assembly_kernel_path));
150123d4529eSJeremy L Thompson   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Assembly Kernel Source -----\n");
15022b730f8bSJeremy L Thompson   CeedCallBackend(CeedLoadSourceToBuffer(ceed, assembly_kernel_path, &assembly_kernel_source));
150323d4529eSJeremy L Thompson   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Assembly Source Complete! -----\n");
1504004e4986SSebastian Grimberg   CeedCallBackend(CeedCompile_Cuda(ceed, assembly_kernel_source, &asmb->module, 10, "NUM_EVAL_MODES_IN", num_eval_modes_in, "NUM_EVAL_MODES_OUT",
1505004e4986SSebastian Grimberg                                    num_eval_modes_out, "NUM_COMP_IN", num_comp_in, "NUM_COMP_OUT", num_comp_out, "NUM_NODES_IN", elem_size_in,
1506004e4986SSebastian Grimberg                                    "NUM_NODES_OUT", elem_size_out, "NUM_QPTS", num_qpts_in, "BLOCK_SIZE",
1507cbfe683aSSebastian Grimberg                                    asmb->block_size_x * asmb->block_size_y * asmb->elems_per_block, "BLOCK_SIZE_Y", asmb->block_size_y,
1508cbfe683aSSebastian Grimberg                                    "USE_CEEDSIZE", use_ceedsize_idx));
1509004e4986SSebastian Grimberg   CeedCallBackend(CeedGetKernel_Cuda(ceed, asmb->module, "LinearAssemble", &asmb->LinearAssemble));
15102b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&assembly_kernel_path));
15112b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&assembly_kernel_source));
1512cc132f9aSnbeams 
1513004e4986SSebastian Grimberg   // Load into B_in, in order that they will be used in eval_modes_in
1514004e4986SSebastian Grimberg   {
1515004e4986SSebastian Grimberg     const CeedInt in_bytes           = elem_size_in * num_qpts_in * num_eval_modes_in * sizeof(CeedScalar);
1516004e4986SSebastian Grimberg     CeedInt       d_in               = 0;
1517004e4986SSebastian Grimberg     CeedEvalMode  eval_modes_in_prev = CEED_EVAL_NONE;
1518004e4986SSebastian Grimberg     bool          has_eval_none      = false;
1519004e4986SSebastian Grimberg     CeedScalar   *identity           = NULL;
1520ca735530SJeremy L Thompson 
1521004e4986SSebastian Grimberg     for (CeedInt i = 0; i < num_eval_modes_in; i++) {
1522004e4986SSebastian Grimberg       has_eval_none = has_eval_none || (eval_modes_in[i] == CEED_EVAL_NONE);
1523004e4986SSebastian Grimberg     }
1524004e4986SSebastian Grimberg     if (has_eval_none) {
1525004e4986SSebastian Grimberg       CeedCallBackend(CeedCalloc(elem_size_in * num_qpts_in, &identity));
1526004e4986SSebastian 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;
1527004e4986SSebastian Grimberg     }
1528cc132f9aSnbeams 
1529004e4986SSebastian Grimberg     CeedCallCuda(ceed, cudaMalloc((void **)&asmb->d_B_in, in_bytes));
1530004e4986SSebastian Grimberg     for (CeedInt i = 0; i < num_eval_modes_in; i++) {
1531004e4986SSebastian Grimberg       const CeedScalar *h_B_in;
1532ca735530SJeremy L Thompson 
1533004e4986SSebastian Grimberg       CeedCallBackend(CeedOperatorGetBasisPointer(basis_in, eval_modes_in[i], identity, &h_B_in));
1534004e4986SSebastian Grimberg       CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_in, eval_modes_in[i], &q_comp));
1535004e4986SSebastian Grimberg       if (q_comp > 1) {
1536004e4986SSebastian Grimberg         if (i == 0 || eval_modes_in[i] != eval_modes_in_prev) d_in = 0;
1537004e4986SSebastian Grimberg         else h_B_in = &h_B_in[(++d_in) * elem_size_in * num_qpts_in];
1538004e4986SSebastian Grimberg       }
1539004e4986SSebastian Grimberg       eval_modes_in_prev = eval_modes_in[i];
1540ca735530SJeremy L Thompson 
1541004e4986SSebastian 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),
1542004e4986SSebastian Grimberg                                     cudaMemcpyHostToDevice));
1543004e4986SSebastian Grimberg     }
1544004e4986SSebastian Grimberg 
1545004e4986SSebastian Grimberg     if (identity) {
1546004e4986SSebastian Grimberg       CeedCallBackend(CeedFree(&identity));
1547cc132f9aSnbeams     }
1548cc132f9aSnbeams   }
1549*8b0f7348SJeremy L Thompson   CeedCallBackend(CeedFree(&eval_modes_in));
1550cc132f9aSnbeams 
1551004e4986SSebastian Grimberg   // Load into B_out, in order that they will be used in eval_modes_out
1552004e4986SSebastian Grimberg   {
1553004e4986SSebastian Grimberg     const CeedInt out_bytes           = elem_size_out * num_qpts_out * num_eval_modes_out * sizeof(CeedScalar);
1554004e4986SSebastian Grimberg     CeedInt       d_out               = 0;
1555004e4986SSebastian Grimberg     CeedEvalMode  eval_modes_out_prev = CEED_EVAL_NONE;
1556004e4986SSebastian Grimberg     bool          has_eval_none       = false;
1557004e4986SSebastian Grimberg     CeedScalar   *identity            = NULL;
1558ca735530SJeremy L Thompson 
1559004e4986SSebastian Grimberg     for (CeedInt i = 0; i < num_eval_modes_out; i++) {
1560004e4986SSebastian Grimberg       has_eval_none = has_eval_none || (eval_modes_out[i] == CEED_EVAL_NONE);
1561004e4986SSebastian Grimberg     }
1562004e4986SSebastian Grimberg     if (has_eval_none) {
1563004e4986SSebastian Grimberg       CeedCallBackend(CeedCalloc(elem_size_out * num_qpts_out, &identity));
1564004e4986SSebastian 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;
1565cc132f9aSnbeams     }
1566cc132f9aSnbeams 
1567004e4986SSebastian Grimberg     CeedCallCuda(ceed, cudaMalloc((void **)&asmb->d_B_out, out_bytes));
1568004e4986SSebastian Grimberg     for (CeedInt i = 0; i < num_eval_modes_out; i++) {
1569004e4986SSebastian Grimberg       const CeedScalar *h_B_out;
1570ca735530SJeremy L Thompson 
1571004e4986SSebastian Grimberg       CeedCallBackend(CeedOperatorGetBasisPointer(basis_out, eval_modes_out[i], identity, &h_B_out));
1572004e4986SSebastian Grimberg       CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_out, eval_modes_out[i], &q_comp));
1573004e4986SSebastian Grimberg       if (q_comp > 1) {
1574004e4986SSebastian Grimberg         if (i == 0 || eval_modes_out[i] != eval_modes_out_prev) d_out = 0;
1575004e4986SSebastian Grimberg         else h_B_out = &h_B_out[(++d_out) * elem_size_out * num_qpts_out];
1576004e4986SSebastian Grimberg       }
1577004e4986SSebastian Grimberg       eval_modes_out_prev = eval_modes_out[i];
1578ca735530SJeremy L Thompson 
1579004e4986SSebastian 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),
1580004e4986SSebastian Grimberg                                     cudaMemcpyHostToDevice));
1581004e4986SSebastian Grimberg     }
1582004e4986SSebastian Grimberg 
1583004e4986SSebastian Grimberg     if (identity) {
1584004e4986SSebastian Grimberg       CeedCallBackend(CeedFree(&identity));
1585cc132f9aSnbeams     }
1586cc132f9aSnbeams   }
1587*8b0f7348SJeremy L Thompson   CeedCallBackend(CeedFree(&eval_modes_out));
1588cc132f9aSnbeams   return CEED_ERROR_SUCCESS;
1589cc132f9aSnbeams }
1590cc132f9aSnbeams 
1591cc132f9aSnbeams //------------------------------------------------------------------------------
1592cefa2673SJeremy L Thompson // Assemble matrix data for COO matrix of assembled operator.
1593cefa2673SJeremy L Thompson // The sparsity pattern is set by CeedOperatorLinearAssembleSymbolic.
1594cefa2673SJeremy L Thompson //
1595034f99fdSJeremy L Thompson // Note that this (and other assembly routines) currently assume only one active input restriction/basis per operator
1596034f99fdSJeremy L Thompson // (could have multiple basis eval modes).
1597cefa2673SJeremy L Thompson // TODO: allow multiple active input restrictions/basis objects
1598cc132f9aSnbeams //------------------------------------------------------------------------------
15992b730f8bSJeremy L Thompson static int CeedSingleOperatorAssemble_Cuda(CeedOperator op, CeedInt offset, CeedVector values) {
1600cc132f9aSnbeams   Ceed                ceed;
1601ca735530SJeremy L Thompson   CeedSize            values_length = 0, assembled_qf_length = 0;
1602004e4986SSebastian Grimberg   CeedInt             use_ceedsize_idx = 0, num_elem_in, num_elem_out, elem_size_in, elem_size_out;
1603ca735530SJeremy L Thompson   CeedScalar         *values_array;
1604004e4986SSebastian Grimberg   const CeedScalar   *assembled_qf_array;
1605ca735530SJeremy L Thompson   CeedVector          assembled_qf   = NULL;
1606004e4986SSebastian Grimberg   CeedElemRestriction assembled_rstr = NULL, rstr_in, rstr_out;
1607004e4986SSebastian Grimberg   CeedRestrictionType rstr_type_in, rstr_type_out;
1608004e4986SSebastian Grimberg   const bool         *orients_in = NULL, *orients_out = NULL;
1609004e4986SSebastian Grimberg   const CeedInt8     *curl_orients_in = NULL, *curl_orients_out = NULL;
1610cc132f9aSnbeams   CeedOperator_Cuda  *impl;
1611ca735530SJeremy L Thompson 
1612ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
16132b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &impl));
1614cc132f9aSnbeams 
1615cc132f9aSnbeams   // Assemble QFunction
1616004e4986SSebastian Grimberg   CeedCallBackend(CeedOperatorLinearAssembleQFunctionBuildOrUpdate(op, &assembled_qf, &assembled_rstr, CEED_REQUEST_IMMEDIATE));
1617004e4986SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionDestroy(&assembled_rstr));
1618004e4986SSebastian Grimberg   CeedCallBackend(CeedVectorGetArrayRead(assembled_qf, CEED_MEM_DEVICE, &assembled_qf_array));
1619cc132f9aSnbeams 
1620f7c1b517Snbeams   CeedCallBackend(CeedVectorGetLength(values, &values_length));
1621f7c1b517Snbeams   CeedCallBackend(CeedVectorGetLength(assembled_qf, &assembled_qf_length));
1622f7c1b517Snbeams   if ((values_length > INT_MAX) || (assembled_qf_length > INT_MAX)) use_ceedsize_idx = 1;
1623004e4986SSebastian Grimberg 
1624f7c1b517Snbeams   // Setup
1625004e4986SSebastian Grimberg   if (!impl->asmb) CeedCallBackend(CeedSingleOperatorAssembleSetup_Cuda(op, use_ceedsize_idx));
1626004e4986SSebastian Grimberg   CeedOperatorAssemble_Cuda *asmb = impl->asmb;
1627004e4986SSebastian Grimberg 
1628004e4986SSebastian Grimberg   assert(asmb != NULL);
1629004e4986SSebastian Grimberg 
1630004e4986SSebastian Grimberg   // Assemble element operator
1631004e4986SSebastian Grimberg   CeedCallBackend(CeedVectorGetArray(values, CEED_MEM_DEVICE, &values_array));
1632004e4986SSebastian Grimberg   values_array += offset;
1633004e4986SSebastian Grimberg 
1634004e4986SSebastian Grimberg   CeedCallBackend(CeedOperatorGetActiveElemRestrictions(op, &rstr_in, &rstr_out));
1635004e4986SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionGetNumElements(rstr_in, &num_elem_in));
1636004e4986SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionGetElementSize(rstr_in, &elem_size_in));
1637004e4986SSebastian Grimberg 
1638004e4986SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionGetType(rstr_in, &rstr_type_in));
1639004e4986SSebastian Grimberg   if (rstr_type_in == CEED_RESTRICTION_ORIENTED) {
1640004e4986SSebastian Grimberg     CeedCallBackend(CeedElemRestrictionGetOrientations(rstr_in, CEED_MEM_DEVICE, &orients_in));
1641004e4986SSebastian Grimberg   } else if (rstr_type_in == CEED_RESTRICTION_CURL_ORIENTED) {
1642004e4986SSebastian Grimberg     CeedCallBackend(CeedElemRestrictionGetCurlOrientations(rstr_in, CEED_MEM_DEVICE, &curl_orients_in));
1643004e4986SSebastian Grimberg   }
1644004e4986SSebastian Grimberg 
1645004e4986SSebastian Grimberg   if (rstr_in != rstr_out) {
1646004e4986SSebastian Grimberg     CeedCallBackend(CeedElemRestrictionGetNumElements(rstr_out, &num_elem_out));
1647004e4986SSebastian Grimberg     CeedCheck(num_elem_in == num_elem_out, ceed, CEED_ERROR_UNSUPPORTED,
1648004e4986SSebastian Grimberg               "Active input and output operator restrictions must have the same number of elements");
1649004e4986SSebastian Grimberg     CeedCallBackend(CeedElemRestrictionGetElementSize(rstr_out, &elem_size_out));
1650004e4986SSebastian Grimberg 
1651004e4986SSebastian Grimberg     CeedCallBackend(CeedElemRestrictionGetType(rstr_out, &rstr_type_out));
1652004e4986SSebastian Grimberg     if (rstr_type_out == CEED_RESTRICTION_ORIENTED) {
1653004e4986SSebastian Grimberg       CeedCallBackend(CeedElemRestrictionGetOrientations(rstr_out, CEED_MEM_DEVICE, &orients_out));
1654004e4986SSebastian Grimberg     } else if (rstr_type_out == CEED_RESTRICTION_CURL_ORIENTED) {
1655004e4986SSebastian Grimberg       CeedCallBackend(CeedElemRestrictionGetCurlOrientations(rstr_out, CEED_MEM_DEVICE, &curl_orients_out));
1656004e4986SSebastian Grimberg     }
1657004e4986SSebastian Grimberg   } else {
1658004e4986SSebastian Grimberg     elem_size_out    = elem_size_in;
1659004e4986SSebastian Grimberg     orients_out      = orients_in;
1660004e4986SSebastian Grimberg     curl_orients_out = curl_orients_in;
1661f7c1b517Snbeams   }
1662f7c1b517Snbeams 
1663cc132f9aSnbeams   // Compute B^T D B
1664004e4986SSebastian Grimberg   CeedInt shared_mem =
1665004e4986SSebastian 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)) *
1666004e4986SSebastian Grimberg       sizeof(CeedScalar);
1667004e4986SSebastian Grimberg   CeedInt grid   = CeedDivUpInt(num_elem_in, asmb->elems_per_block);
1668004e4986SSebastian Grimberg   void   *args[] = {(void *)&num_elem_in, &asmb->d_B_in,     &asmb->d_B_out,      &orients_in,  &curl_orients_in,
1669004e4986SSebastian Grimberg                     &orients_out,         &curl_orients_out, &assembled_qf_array, &values_array};
1670ca735530SJeremy L Thompson 
16712b730f8bSJeremy L Thompson   CeedCallBackend(
1672004e4986SSebastian Grimberg       CeedRunKernelDimShared_Cuda(ceed, asmb->LinearAssemble, grid, asmb->block_size_x, asmb->block_size_y, asmb->elems_per_block, shared_mem, args));
1673cc132f9aSnbeams 
1674cc132f9aSnbeams   // Restore arrays
16752b730f8bSJeremy L Thompson   CeedCallBackend(CeedVectorRestoreArray(values, &values_array));
1676004e4986SSebastian Grimberg   CeedCallBackend(CeedVectorRestoreArrayRead(assembled_qf, &assembled_qf_array));
1677cc132f9aSnbeams 
1678cc132f9aSnbeams   // Cleanup
16792b730f8bSJeremy L Thompson   CeedCallBackend(CeedVectorDestroy(&assembled_qf));
1680004e4986SSebastian Grimberg   if (rstr_type_in == CEED_RESTRICTION_ORIENTED) {
1681004e4986SSebastian Grimberg     CeedCallBackend(CeedElemRestrictionRestoreOrientations(rstr_in, &orients_in));
1682004e4986SSebastian Grimberg   } else if (rstr_type_in == CEED_RESTRICTION_CURL_ORIENTED) {
1683004e4986SSebastian Grimberg     CeedCallBackend(CeedElemRestrictionRestoreCurlOrientations(rstr_in, &curl_orients_in));
1684004e4986SSebastian Grimberg   }
1685004e4986SSebastian Grimberg   if (rstr_in != rstr_out) {
1686004e4986SSebastian Grimberg     if (rstr_type_out == CEED_RESTRICTION_ORIENTED) {
1687004e4986SSebastian Grimberg       CeedCallBackend(CeedElemRestrictionRestoreOrientations(rstr_out, &orients_out));
1688004e4986SSebastian Grimberg     } else if (rstr_type_out == CEED_RESTRICTION_CURL_ORIENTED) {
1689004e4986SSebastian Grimberg       CeedCallBackend(CeedElemRestrictionRestoreCurlOrientations(rstr_out, &curl_orients_out));
1690004e4986SSebastian Grimberg     }
1691004e4986SSebastian Grimberg   }
1692cc132f9aSnbeams   return CEED_ERROR_SUCCESS;
1693cc132f9aSnbeams }
1694cc132f9aSnbeams 
1695cc132f9aSnbeams //------------------------------------------------------------------------------
1696756ca9e9SJeremy L Thompson // Assemble Linear QFunction AtPoints
1697756ca9e9SJeremy L Thompson //------------------------------------------------------------------------------
1698756ca9e9SJeremy L Thompson static int CeedOperatorLinearAssembleQFunctionAtPoints_Cuda(CeedOperator op, CeedVector *assembled, CeedElemRestriction *rstr, CeedRequest *request) {
1699756ca9e9SJeremy L Thompson   return CeedError(CeedOperatorReturnCeed(op), CEED_ERROR_BACKEND, "Backend does not implement CeedOperatorLinearAssembleQFunction");
1700756ca9e9SJeremy L Thompson }
1701756ca9e9SJeremy L Thompson 
1702756ca9e9SJeremy L Thompson //------------------------------------------------------------------------------
1703756ca9e9SJeremy L Thompson // Assemble Linear Diagonal AtPoints
1704756ca9e9SJeremy L Thompson //------------------------------------------------------------------------------
1705756ca9e9SJeremy L Thompson static int CeedOperatorLinearAssembleAddDiagonalAtPoints_Cuda(CeedOperator op, CeedVector assembled, CeedRequest *request) {
1706111870feSJeremy L Thompson   CeedInt             max_num_points, *num_points, num_elem, num_input_fields, num_output_fields;
17078bbba8cdSJeremy L Thompson   Ceed                ceed;
17088bbba8cdSJeremy L Thompson   CeedVector          active_e_vec_in, active_e_vec_out;
1709349fb27dSJeremy L Thompson   CeedQFunctionField *qf_input_fields, *qf_output_fields;
1710349fb27dSJeremy L Thompson   CeedQFunction       qf;
1711349fb27dSJeremy L Thompson   CeedOperatorField  *op_input_fields, *op_output_fields;
1712349fb27dSJeremy L Thompson   CeedOperator_Cuda  *impl;
1713349fb27dSJeremy L Thompson 
17148bbba8cdSJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
1715349fb27dSJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &impl));
1716349fb27dSJeremy L Thompson   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
1717349fb27dSJeremy L Thompson   CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem));
1718349fb27dSJeremy L Thompson   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
1719349fb27dSJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
1720349fb27dSJeremy L Thompson 
1721349fb27dSJeremy L Thompson   // Setup
1722349fb27dSJeremy L Thompson   CeedCallBackend(CeedOperatorSetupAtPoints_Cuda(op));
1723111870feSJeremy L Thompson   num_points     = impl->num_points;
1724349fb27dSJeremy L Thompson   max_num_points = impl->max_num_points;
1725349fb27dSJeremy L Thompson 
17268bbba8cdSJeremy L Thompson   // Work vector
17278bbba8cdSJeremy L Thompson   CeedCallBackend(CeedGetWorkVector(ceed, impl->max_active_e_vec_len, &active_e_vec_in));
17288bbba8cdSJeremy L Thompson   CeedCallBackend(CeedGetWorkVector(ceed, impl->max_active_e_vec_len, &active_e_vec_out));
172954404f0bSJeremy L Thompson   {
173054404f0bSJeremy L Thompson     CeedSize length_in, length_out;
173154404f0bSJeremy L Thompson 
173254404f0bSJeremy L Thompson     CeedCallBackend(CeedVectorGetLength(active_e_vec_in, &length_in));
173354404f0bSJeremy L Thompson     CeedCallBackend(CeedVectorGetLength(active_e_vec_out, &length_out));
173454404f0bSJeremy L Thompson     // Need input e_vec to be longer
173554404f0bSJeremy L Thompson     if (length_in < length_out) {
173654404f0bSJeremy L Thompson       CeedVector temp = active_e_vec_in;
173754404f0bSJeremy L Thompson 
173854404f0bSJeremy L Thompson       active_e_vec_in  = active_e_vec_out;
173954404f0bSJeremy L Thompson       active_e_vec_out = temp;
174054404f0bSJeremy L Thompson     }
174154404f0bSJeremy L Thompson   }
17428a213570SJeremy L Thompson 
1743349fb27dSJeremy L Thompson   // Get point coordinates
1744349fb27dSJeremy L Thompson   if (!impl->point_coords_elem) {
1745349fb27dSJeremy L Thompson     CeedVector          point_coords = NULL;
1746349fb27dSJeremy L Thompson     CeedElemRestriction rstr_points  = NULL;
1747349fb27dSJeremy L Thompson 
1748349fb27dSJeremy L Thompson     CeedCallBackend(CeedOperatorAtPointsGetPoints(op, &rstr_points, &point_coords));
1749349fb27dSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionCreateVector(rstr_points, NULL, &impl->point_coords_elem));
1750349fb27dSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionApply(rstr_points, CEED_NOTRANSPOSE, point_coords, impl->point_coords_elem, request));
1751*8b0f7348SJeremy L Thompson     CeedCallBackend(CeedVectorDestroy(&point_coords));
1752*8b0f7348SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionDestroy(&rstr_points));
1753349fb27dSJeremy L Thompson   }
1754349fb27dSJeremy L Thompson 
1755034f99fdSJeremy L Thompson   // Process inputs
175643e13feeSJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
17578bbba8cdSJeremy L Thompson     CeedCallBackend(CeedOperatorInputRestrict_Cuda(op_input_fields[i], qf_input_fields[i], i, NULL, NULL, true, impl, request));
17588bbba8cdSJeremy L Thompson     CeedCallBackend(CeedOperatorInputBasisAtPoints_Cuda(op_input_fields[i], qf_input_fields[i], i, NULL, NULL, num_elem, num_points, true, impl));
175943e13feeSJeremy L Thompson   }
176043e13feeSJeremy L Thompson 
1761382e9c83SJeremy L Thompson   // Clear active input Qvecs
1762382e9c83SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
1763382e9c83SJeremy L Thompson     CeedVector vec;
1764382e9c83SJeremy L Thompson 
1765382e9c83SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec));
1766382e9c83SJeremy L Thompson     if (vec != CEED_VECTOR_ACTIVE) continue;
1767382e9c83SJeremy L Thompson     CeedCallBackend(CeedVectorSetValue(impl->q_vecs_in[i], 0.0));
1768382e9c83SJeremy L Thompson   }
1769382e9c83SJeremy L Thompson 
1770349fb27dSJeremy L Thompson   // Output pointers, as necessary
1771349fb27dSJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
1772349fb27dSJeremy L Thompson     CeedEvalMode eval_mode;
1773349fb27dSJeremy L Thompson 
1774349fb27dSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
1775349fb27dSJeremy L Thompson     if (eval_mode == CEED_EVAL_NONE) {
17768bbba8cdSJeremy L Thompson       CeedScalar *e_vec_array;
17778bbba8cdSJeremy L Thompson 
17788bbba8cdSJeremy L Thompson       CeedCallBackend(CeedVectorGetArrayWrite(impl->e_vecs_out[i], CEED_MEM_DEVICE, &e_vec_array));
17798bbba8cdSJeremy L Thompson       CeedCallBackend(CeedVectorSetArray(impl->q_vecs_out[i], CEED_MEM_DEVICE, CEED_USE_POINTER, e_vec_array));
1780349fb27dSJeremy L Thompson     }
1781349fb27dSJeremy L Thompson   }
1782349fb27dSJeremy L Thompson 
1783349fb27dSJeremy L Thompson   // Loop over active fields
1784349fb27dSJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
1785382e9c83SJeremy L Thompson     bool                is_active_at_points = true;
1786382e9c83SJeremy L Thompson     CeedInt             elem_size = 1, num_comp_active = 1, e_vec_size = 0;
1787382e9c83SJeremy L Thompson     CeedRestrictionType rstr_type;
1788034f99fdSJeremy L Thompson     CeedVector          l_vec;
1789382e9c83SJeremy L Thompson     CeedElemRestriction elem_rstr;
1790382e9c83SJeremy L Thompson 
1791034f99fdSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &l_vec));
1792382e9c83SJeremy L Thompson     // -- Skip non-active input
1793034f99fdSJeremy L Thompson     if (l_vec != CEED_VECTOR_ACTIVE) continue;
1794382e9c83SJeremy L Thompson 
1795382e9c83SJeremy L Thompson     // -- Get active restriction type
1796382e9c83SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr));
1797382e9c83SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetType(elem_rstr, &rstr_type));
1798382e9c83SJeremy L Thompson     is_active_at_points = rstr_type == CEED_RESTRICTION_POINTS;
1799382e9c83SJeremy L Thompson     if (!is_active_at_points) CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
1800382e9c83SJeremy L Thompson     else elem_size = max_num_points;
1801382e9c83SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp_active));
1802382e9c83SJeremy L Thompson 
1803382e9c83SJeremy L Thompson     e_vec_size = elem_size * num_comp_active;
1804382e9c83SJeremy L Thompson     for (CeedInt s = 0; s < e_vec_size; s++) {
1805349fb27dSJeremy L Thompson       bool         is_active_input = false;
1806349fb27dSJeremy L Thompson       CeedEvalMode eval_mode;
18078bbba8cdSJeremy L Thompson       CeedVector   l_vec, q_vec = impl->q_vecs_in[i];
1808349fb27dSJeremy L Thompson       CeedBasis    basis;
1809349fb27dSJeremy L Thompson 
18108bbba8cdSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &l_vec));
1811349fb27dSJeremy L Thompson       // Skip non-active input
18128bbba8cdSJeremy L Thompson       is_active_input = l_vec == CEED_VECTOR_ACTIVE;
1813349fb27dSJeremy L Thompson       if (!is_active_input) continue;
1814349fb27dSJeremy L Thompson 
1815349fb27dSJeremy L Thompson       // Update unit vector
18168bbba8cdSJeremy L Thompson       if (s == 0) CeedCallBackend(CeedVectorSetValue(active_e_vec_in, 0.0));
18178bbba8cdSJeremy L Thompson       else CeedCallBackend(CeedVectorSetValueStrided(active_e_vec_in, s - 1, e_vec_size, 0.0));
18188bbba8cdSJeremy L Thompson       CeedCallBackend(CeedVectorSetValueStrided(active_e_vec_in, s, e_vec_size, 1.0));
1819349fb27dSJeremy L Thompson 
1820349fb27dSJeremy L Thompson       // Basis action
1821349fb27dSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
1822349fb27dSJeremy L Thompson       switch (eval_mode) {
18238bbba8cdSJeremy L Thompson         case CEED_EVAL_NONE: {
18248bbba8cdSJeremy L Thompson           const CeedScalar *e_vec_array;
18258bbba8cdSJeremy L Thompson 
18268bbba8cdSJeremy L Thompson           CeedCallBackend(CeedVectorGetArrayRead(active_e_vec_in, CEED_MEM_DEVICE, &e_vec_array));
182754404f0bSJeremy L Thompson           CeedCallBackend(CeedVectorSetArray(q_vec, CEED_MEM_DEVICE, CEED_USE_POINTER, (CeedScalar *)e_vec_array));
1828349fb27dSJeremy L Thompson           break;
18298bbba8cdSJeremy L Thompson         }
1830349fb27dSJeremy L Thompson         case CEED_EVAL_INTERP:
1831349fb27dSJeremy L Thompson         case CEED_EVAL_GRAD:
1832349fb27dSJeremy L Thompson         case CEED_EVAL_DIV:
1833349fb27dSJeremy L Thompson         case CEED_EVAL_CURL:
1834349fb27dSJeremy L Thompson           CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
18358bbba8cdSJeremy L Thompson           CeedCallBackend(
18368bbba8cdSJeremy L Thompson               CeedBasisApplyAtPoints(basis, num_elem, num_points, CEED_NOTRANSPOSE, eval_mode, impl->point_coords_elem, active_e_vec_in, q_vec));
1837349fb27dSJeremy L Thompson           break;
1838349fb27dSJeremy L Thompson         case CEED_EVAL_WEIGHT:
1839349fb27dSJeremy L Thompson           break;  // No action
1840349fb27dSJeremy L Thompson       }
1841349fb27dSJeremy L Thompson 
1842349fb27dSJeremy L Thompson       // Q function
1843349fb27dSJeremy L Thompson       CeedCallBackend(CeedQFunctionApply(qf, num_elem * max_num_points, impl->q_vecs_in, impl->q_vecs_out));
1844349fb27dSJeremy L Thompson 
1845349fb27dSJeremy L Thompson       // Output basis apply if needed
1846382e9c83SJeremy L Thompson       for (CeedInt j = 0; j < num_output_fields; j++) {
1847349fb27dSJeremy L Thompson         bool                is_active_output = false;
1848382e9c83SJeremy L Thompson         CeedInt             elem_size        = 0;
1849382e9c83SJeremy L Thompson         CeedRestrictionType rstr_type;
1850349fb27dSJeremy L Thompson         CeedEvalMode        eval_mode;
18518bbba8cdSJeremy L Thompson         CeedVector          l_vec, e_vec = impl->e_vecs_out[j], q_vec = impl->q_vecs_out[j];
1852349fb27dSJeremy L Thompson         CeedElemRestriction elem_rstr;
1853349fb27dSJeremy L Thompson         CeedBasis           basis;
1854349fb27dSJeremy L Thompson 
1855034f99fdSJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[j], &l_vec));
1856382e9c83SJeremy L Thompson         // ---- Skip non-active output
1857034f99fdSJeremy L Thompson         is_active_output = l_vec == CEED_VECTOR_ACTIVE;
1858349fb27dSJeremy L Thompson         if (!is_active_output) continue;
18598bbba8cdSJeremy L Thompson         if (!e_vec) e_vec = active_e_vec_out;
1860349fb27dSJeremy L Thompson 
1861382e9c83SJeremy L Thompson         // ---- Check if elem size matches
1862382e9c83SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[j], &elem_rstr));
1863382e9c83SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetType(elem_rstr, &rstr_type));
1864382e9c83SJeremy L Thompson         if (is_active_at_points && rstr_type != CEED_RESTRICTION_POINTS) continue;
1865382e9c83SJeremy L Thompson         if (rstr_type == CEED_RESTRICTION_POINTS) {
1866382e9c83SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetMaxPointsInElement(elem_rstr, &elem_size));
1867382e9c83SJeremy L Thompson         } else {
1868382e9c83SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
1869382e9c83SJeremy L Thompson         }
1870382e9c83SJeremy L Thompson         {
1871382e9c83SJeremy L Thompson           CeedInt num_comp = 0;
1872382e9c83SJeremy L Thompson 
1873382e9c83SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
1874382e9c83SJeremy L Thompson           if (e_vec_size != num_comp * elem_size) continue;
1875382e9c83SJeremy L Thompson         }
1876382e9c83SJeremy L Thompson 
1877349fb27dSJeremy L Thompson         // Basis action
1878382e9c83SJeremy L Thompson         CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[j], &eval_mode));
1879349fb27dSJeremy L Thompson         switch (eval_mode) {
18808bbba8cdSJeremy L Thompson           case CEED_EVAL_NONE: {
18818bbba8cdSJeremy L Thompson             CeedScalar *e_vec_array;
18828bbba8cdSJeremy L Thompson 
18838bbba8cdSJeremy L Thompson             CeedCallBackend(CeedVectorTakeArray(q_vec, CEED_MEM_DEVICE, &e_vec_array));
18848bbba8cdSJeremy L Thompson             CeedCallBackend(CeedVectorRestoreArray(e_vec, &e_vec_array));
1885349fb27dSJeremy L Thompson             break;
18868bbba8cdSJeremy L Thompson           }
1887349fb27dSJeremy L Thompson           case CEED_EVAL_INTERP:
1888349fb27dSJeremy L Thompson           case CEED_EVAL_GRAD:
1889349fb27dSJeremy L Thompson           case CEED_EVAL_DIV:
1890349fb27dSJeremy L Thompson           case CEED_EVAL_CURL:
1891382e9c83SJeremy L Thompson             CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[j], &basis));
18928bbba8cdSJeremy L Thompson             CeedCallBackend(CeedBasisApplyAtPoints(basis, num_elem, num_points, CEED_TRANSPOSE, eval_mode, impl->point_coords_elem, q_vec, e_vec));
1893349fb27dSJeremy L Thompson             break;
1894349fb27dSJeremy L Thompson           // LCOV_EXCL_START
1895349fb27dSJeremy L Thompson           case CEED_EVAL_WEIGHT: {
1896349fb27dSJeremy L Thompson             return CeedError(CeedOperatorReturnCeed(op), CEED_ERROR_BACKEND, "CEED_EVAL_WEIGHT cannot be an output evaluation mode");
1897349fb27dSJeremy L Thompson             // LCOV_EXCL_STOP
1898349fb27dSJeremy L Thompson           }
1899349fb27dSJeremy L Thompson         }
1900349fb27dSJeremy L Thompson 
1901349fb27dSJeremy L Thompson         // Mask output e-vec
19028bbba8cdSJeremy L Thompson         CeedCallBackend(CeedVectorPointwiseMult(e_vec, active_e_vec_in, e_vec));
1903349fb27dSJeremy L Thompson 
1904349fb27dSJeremy L Thompson         // Restrict
1905382e9c83SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[j], &elem_rstr));
19068bbba8cdSJeremy L Thompson         CeedCallBackend(CeedElemRestrictionApply(elem_rstr, CEED_TRANSPOSE, e_vec, assembled, request));
1907349fb27dSJeremy L Thompson 
1908349fb27dSJeremy L Thompson         // Reset q_vec for
1909349fb27dSJeremy L Thompson         if (eval_mode == CEED_EVAL_NONE) {
19108bbba8cdSJeremy L Thompson           CeedScalar *e_vec_array;
19118bbba8cdSJeremy L Thompson 
19128bbba8cdSJeremy L Thompson           CeedCallBackend(CeedVectorGetArrayWrite(e_vec, CEED_MEM_DEVICE, &e_vec_array));
19138bbba8cdSJeremy L Thompson           CeedCallBackend(CeedVectorSetArray(q_vec, CEED_MEM_DEVICE, CEED_USE_POINTER, e_vec_array));
1914349fb27dSJeremy L Thompson         }
1915349fb27dSJeremy L Thompson       }
1916382e9c83SJeremy L Thompson 
1917382e9c83SJeremy L Thompson       // Reset vec
19188bbba8cdSJeremy L Thompson       if (s == e_vec_size - 1 && i != num_input_fields - 1) CeedCallBackend(CeedVectorSetValue(q_vec, 0.0));
191986e10729SJeremy L Thompson     }
1920349fb27dSJeremy L Thompson   }
1921349fb27dSJeremy L Thompson 
1922349fb27dSJeremy L Thompson   // Restore CEED_EVAL_NONE
1923349fb27dSJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
1924349fb27dSJeremy L Thompson     CeedEvalMode eval_mode;
1925349fb27dSJeremy L Thompson 
1926349fb27dSJeremy L Thompson     // Get eval_mode
1927349fb27dSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
1928349fb27dSJeremy L Thompson 
1929349fb27dSJeremy L Thompson     // Restore evec
1930349fb27dSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
1931349fb27dSJeremy L Thompson     if (eval_mode == CEED_EVAL_NONE) {
19328bbba8cdSJeremy L Thompson       CeedScalar *e_vec_array;
19338bbba8cdSJeremy L Thompson 
19348bbba8cdSJeremy L Thompson       CeedCallBackend(CeedVectorTakeArray(impl->q_vecs_in[i], CEED_MEM_DEVICE, &e_vec_array));
19358bbba8cdSJeremy L Thompson       CeedCallBackend(CeedVectorRestoreArray(impl->e_vecs_in[i], &e_vec_array));
1936349fb27dSJeremy L Thompson     }
1937349fb27dSJeremy L Thompson   }
1938349fb27dSJeremy L Thompson 
1939349fb27dSJeremy L Thompson   // Restore input arrays
194043e13feeSJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
19418bbba8cdSJeremy L Thompson     CeedCallBackend(CeedOperatorInputRestore_Cuda(op_input_fields[i], qf_input_fields[i], i, NULL, NULL, true, impl));
194243e13feeSJeremy L Thompson   }
1943*8b0f7348SJeremy L Thompson 
1944*8b0f7348SJeremy L Thompson   // Restore work vector
1945*8b0f7348SJeremy L Thompson   CeedCallBackend(CeedRestoreWorkVector(ceed, &active_e_vec_in));
1946*8b0f7348SJeremy L Thompson   CeedCallBackend(CeedRestoreWorkVector(ceed, &active_e_vec_out));
1947349fb27dSJeremy L Thompson   return CEED_ERROR_SUCCESS;
1948756ca9e9SJeremy L Thompson }
1949756ca9e9SJeremy L Thompson 
1950756ca9e9SJeremy L Thompson //------------------------------------------------------------------------------
19510d0321e0SJeremy L Thompson // Create operator
19520d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
19530d0321e0SJeremy L Thompson int CeedOperatorCreate_Cuda(CeedOperator op) {
19540d0321e0SJeremy L Thompson   Ceed               ceed;
19550d0321e0SJeremy L Thompson   CeedOperator_Cuda *impl;
19560d0321e0SJeremy L Thompson 
1957ca735530SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
19582b730f8bSJeremy L Thompson   CeedCallBackend(CeedCalloc(1, &impl));
19592b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorSetData(op, impl));
19600d0321e0SJeremy L Thompson 
19612b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleQFunction", CeedOperatorLinearAssembleQFunction_Cuda));
19622b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleQFunctionUpdate", CeedOperatorLinearAssembleQFunctionUpdate_Cuda));
19632b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleAddDiagonal", CeedOperatorLinearAssembleAddDiagonal_Cuda));
19642b730f8bSJeremy L Thompson   CeedCallBackend(
19652b730f8bSJeremy L Thompson       CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleAddPointBlockDiagonal", CeedOperatorLinearAssembleAddPointBlockDiagonal_Cuda));
19662b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleSingle", CeedSingleOperatorAssemble_Cuda));
19672b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "ApplyAdd", CeedOperatorApplyAdd_Cuda));
19682b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "Destroy", CeedOperatorDestroy_Cuda));
19690d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
19700d0321e0SJeremy L Thompson }
19710d0321e0SJeremy L Thompson 
19720d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
197367d9480aSJeremy L Thompson // Create operator AtPoints
1974756ca9e9SJeremy L Thompson //------------------------------------------------------------------------------
1975756ca9e9SJeremy L Thompson int CeedOperatorCreateAtPoints_Cuda(CeedOperator op) {
1976756ca9e9SJeremy L Thompson   Ceed               ceed;
1977756ca9e9SJeremy L Thompson   CeedOperator_Cuda *impl;
1978756ca9e9SJeremy L Thompson 
1979756ca9e9SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
1980756ca9e9SJeremy L Thompson   CeedCallBackend(CeedCalloc(1, &impl));
1981756ca9e9SJeremy L Thompson   CeedCallBackend(CeedOperatorSetData(op, impl));
1982756ca9e9SJeremy L Thompson 
1983756ca9e9SJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleQFunction", CeedOperatorLinearAssembleQFunctionAtPoints_Cuda));
1984756ca9e9SJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleAddDiagonal", CeedOperatorLinearAssembleAddDiagonalAtPoints_Cuda));
1985756ca9e9SJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "ApplyAdd", CeedOperatorApplyAddAtPoints_Cuda));
1986756ca9e9SJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "Destroy", CeedOperatorDestroy_Cuda));
1987756ca9e9SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1988756ca9e9SJeremy L Thompson }
1989756ca9e9SJeremy L Thompson 
1990756ca9e9SJeremy L Thompson //------------------------------------------------------------------------------
1991