xref: /libCEED/rust/libceed-sys/c-src/backends/hip-ref/ceed-hip-ref-operator.c (revision 506b1a0c535297a01950817d4491440b9ddb5287)
13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
30d0321e0SJeremy L Thompson //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
50d0321e0SJeremy L Thompson //
63d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
70d0321e0SJeremy L Thompson 
849aac155SJeremy L Thompson #include <ceed.h>
90d0321e0SJeremy L Thompson #include <ceed/backend.h>
1007b31e0eSJeremy L Thompson #include <ceed/jit-tools.h>
11c85e8640SSebastian Grimberg #include <assert.h>
120d0321e0SJeremy L Thompson #include <stdbool.h>
130d0321e0SJeremy L Thompson #include <string.h>
14c85e8640SSebastian Grimberg #include <hip/hip_runtime.h>
152b730f8bSJeremy L Thompson 
1649aac155SJeremy L Thompson #include "../hip/ceed-hip-common.h"
170d0321e0SJeremy L Thompson #include "../hip/ceed-hip-compile.h"
182b730f8bSJeremy L Thompson #include "ceed-hip-ref.h"
190d0321e0SJeremy L Thompson 
200d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
210d0321e0SJeremy L Thompson // Destroy operator
220d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
230d0321e0SJeremy L Thompson static int CeedOperatorDestroy_Hip(CeedOperator op) {
240d0321e0SJeremy L Thompson   CeedOperator_Hip *impl;
25b7453713SJeremy L Thompson 
262b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &impl));
270d0321e0SJeremy L Thompson 
280d0321e0SJeremy L Thompson   // Apply data
29b7453713SJeremy L Thompson   for (CeedInt i = 0; i < impl->num_inputs + impl->num_outputs; i++) {
30b7453713SJeremy L Thompson     CeedCallBackend(CeedVectorDestroy(&impl->e_vecs[i]));
310d0321e0SJeremy L Thompson   }
32b7453713SJeremy L Thompson   CeedCallBackend(CeedFree(&impl->e_vecs));
330d0321e0SJeremy L Thompson 
34b7453713SJeremy L Thompson   for (CeedInt i = 0; i < impl->num_inputs; i++) {
35b7453713SJeremy L Thompson     CeedCallBackend(CeedVectorDestroy(&impl->q_vecs_in[i]));
360d0321e0SJeremy L Thompson   }
37b7453713SJeremy L Thompson   CeedCallBackend(CeedFree(&impl->q_vecs_in));
380d0321e0SJeremy L Thompson 
39b7453713SJeremy L Thompson   for (CeedInt i = 0; i < impl->num_outputs; i++) {
40b7453713SJeremy L Thompson     CeedCallBackend(CeedVectorDestroy(&impl->q_vecs_out[i]));
410d0321e0SJeremy L Thompson   }
42b7453713SJeremy L Thompson   CeedCallBackend(CeedFree(&impl->q_vecs_out));
430d0321e0SJeremy L Thompson 
44b2165e7aSSebastian Grimberg   // QFunction assembly data
45b7453713SJeremy L Thompson   for (CeedInt i = 0; i < impl->num_active_in; i++) {
46b7453713SJeremy L Thompson     CeedCallBackend(CeedVectorDestroy(&impl->qf_active_in[i]));
470d0321e0SJeremy L Thompson   }
48b7453713SJeremy L Thompson   CeedCallBackend(CeedFree(&impl->qf_active_in));
490d0321e0SJeremy L Thompson 
500d0321e0SJeremy L Thompson   // Diag data
510d0321e0SJeremy L Thompson   if (impl->diag) {
520d0321e0SJeremy L Thompson     Ceed ceed;
53b7453713SJeremy L Thompson 
542b730f8bSJeremy L Thompson     CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
552b730f8bSJeremy L Thompson     CeedCallHip(ceed, hipModuleUnload(impl->diag->module));
56b7453713SJeremy L Thompson     CeedCallBackend(CeedFree(&impl->diag->h_e_mode_in));
57b7453713SJeremy L Thompson     CeedCallBackend(CeedFree(&impl->diag->h_e_mode_out));
58b7453713SJeremy L Thompson     CeedCallHip(ceed, hipFree(impl->diag->d_e_mode_in));
59b7453713SJeremy L Thompson     CeedCallHip(ceed, hipFree(impl->diag->d_e_mode_out));
602b730f8bSJeremy L Thompson     CeedCallHip(ceed, hipFree(impl->diag->d_identity));
61b7453713SJeremy L Thompson     CeedCallHip(ceed, hipFree(impl->diag->d_interp_in));
62b7453713SJeremy L Thompson     CeedCallHip(ceed, hipFree(impl->diag->d_interp_out));
63b7453713SJeremy L Thompson     CeedCallHip(ceed, hipFree(impl->diag->d_grad_in));
64b7453713SJeremy L Thompson     CeedCallHip(ceed, hipFree(impl->diag->d_grad_out));
65b7453713SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionDestroy(&impl->diag->point_block_diag_rstr));
66b7453713SJeremy L Thompson     CeedCallBackend(CeedVectorDestroy(&impl->diag->elem_diag));
67b7453713SJeremy L Thompson     CeedCallBackend(CeedVectorDestroy(&impl->diag->point_block_elem_diag));
680d0321e0SJeremy L Thompson   }
692b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&impl->diag));
700d0321e0SJeremy L Thompson 
71a835093fSnbeams   if (impl->asmb) {
72a835093fSnbeams     Ceed ceed;
73b7453713SJeremy L Thompson 
742b730f8bSJeremy L Thompson     CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
752b730f8bSJeremy L Thompson     CeedCallHip(ceed, hipModuleUnload(impl->asmb->module));
762b730f8bSJeremy L Thompson     CeedCallHip(ceed, hipFree(impl->asmb->d_B_in));
772b730f8bSJeremy L Thompson     CeedCallHip(ceed, hipFree(impl->asmb->d_B_out));
78a835093fSnbeams   }
792b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&impl->asmb));
80a835093fSnbeams 
812b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&impl));
820d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
830d0321e0SJeremy L Thompson }
840d0321e0SJeremy L Thompson 
850d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
860d0321e0SJeremy L Thompson // Setup infields or outfields
870d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
88b7453713SJeremy L Thompson static int CeedOperatorSetupFields_Hip(CeedQFunction qf, CeedOperator op, bool is_input, CeedVector *e_vecs, CeedVector *q_vecs, CeedInt start_e,
89b7453713SJeremy L Thompson                                        CeedInt num_fields, CeedInt Q, CeedInt num_elem) {
900d0321e0SJeremy L Thompson   Ceed                ceed;
91b7453713SJeremy L Thompson   CeedQFunctionField *qf_fields;
92b7453713SJeremy L Thompson   CeedOperatorField  *op_fields;
930d0321e0SJeremy L Thompson 
94b7453713SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
95b7453713SJeremy L Thompson   if (is_input) {
96b7453713SJeremy L Thompson     CeedCallBackend(CeedOperatorGetFields(op, NULL, &op_fields, NULL, NULL));
97b7453713SJeremy L Thompson     CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_fields, NULL, NULL));
980d0321e0SJeremy L Thompson   } else {
99b7453713SJeremy L Thompson     CeedCallBackend(CeedOperatorGetFields(op, NULL, NULL, NULL, &op_fields));
100b7453713SJeremy L Thompson     CeedCallBackend(CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qf_fields));
1010d0321e0SJeremy L Thompson   }
1020d0321e0SJeremy L Thompson 
1030d0321e0SJeremy L Thompson   // Loop over fields
104b7453713SJeremy L Thompson   for (CeedInt i = 0; i < num_fields; i++) {
105b7453713SJeremy L Thompson     bool                is_strided, skip_restriction;
106b7453713SJeremy L Thompson     CeedSize            q_size;
107b7453713SJeremy L Thompson     CeedInt             dim, size;
108b7453713SJeremy L Thompson     CeedEvalMode        e_mode;
109b7453713SJeremy L Thompson     CeedVector          vec;
110b7453713SJeremy L Thompson     CeedElemRestriction elem_rstr;
111b7453713SJeremy L Thompson     CeedBasis           basis;
1120d0321e0SJeremy L Thompson 
113b7453713SJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &e_mode));
114b7453713SJeremy L Thompson     is_strided       = false;
115b7453713SJeremy L Thompson     skip_restriction = false;
116b7453713SJeremy L Thompson     if (e_mode != CEED_EVAL_WEIGHT) {
117b7453713SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_fields[i], &elem_rstr));
1180d0321e0SJeremy L Thompson 
1190d0321e0SJeremy L Thompson       // Check whether this field can skip the element restriction:
120b7453713SJeremy L Thompson       // must be passive input, with e_mode NONE, and have a strided restriction with CEED_STRIDES_BACKEND.
1210d0321e0SJeremy L Thompson 
1220d0321e0SJeremy L Thompson       // First, check whether the field is input or output:
123b7453713SJeremy L Thompson       if (is_input) {
1240d0321e0SJeremy L Thompson         // Check for passive input:
125b7453713SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec));
126b7453713SJeremy L Thompson         if (vec != CEED_VECTOR_ACTIVE) {
127b7453713SJeremy L Thompson           // Check e_mode
128b7453713SJeremy L Thompson           if (e_mode == CEED_EVAL_NONE) {
1290d0321e0SJeremy L Thompson             // Check for strided restriction
130b7453713SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionIsStrided(elem_rstr, &is_strided));
131b7453713SJeremy L Thompson             if (is_strided) {
1320d0321e0SJeremy L Thompson               // Check if vector is already in preferred backend ordering
133b7453713SJeremy L Thompson               CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &skip_restriction));
1340d0321e0SJeremy L Thompson             }
1350d0321e0SJeremy L Thompson           }
1360d0321e0SJeremy L Thompson         }
1370d0321e0SJeremy L Thompson       }
138b7453713SJeremy L Thompson       if (skip_restriction) {
139ea61e9acSJeremy L Thompson         // We do not need an E-Vector, but will use the input field vector's data directly in the operator application.
140b7453713SJeremy L Thompson         e_vecs[i + start_e] = NULL;
1410d0321e0SJeremy L Thompson       } else {
142b7453713SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionCreateVector(elem_rstr, NULL, &e_vecs[i + start_e]));
1430d0321e0SJeremy L Thompson       }
1440d0321e0SJeremy L Thompson     }
1450d0321e0SJeremy L Thompson 
146b7453713SJeremy L Thompson     switch (e_mode) {
1470d0321e0SJeremy L Thompson       case CEED_EVAL_NONE:
148b7453713SJeremy L Thompson         CeedCallBackend(CeedQFunctionFieldGetSize(qf_fields[i], &size));
149b7453713SJeremy L Thompson         q_size = (CeedSize)num_elem * Q * size;
150b7453713SJeremy L Thompson         CeedCallBackend(CeedVectorCreate(ceed, q_size, &q_vecs[i]));
1510d0321e0SJeremy L Thompson         break;
1520d0321e0SJeremy L Thompson       case CEED_EVAL_INTERP:
153b7453713SJeremy L Thompson         CeedCallBackend(CeedQFunctionFieldGetSize(qf_fields[i], &size));
154b7453713SJeremy L Thompson         q_size = (CeedSize)num_elem * Q * size;
155b7453713SJeremy L Thompson         CeedCallBackend(CeedVectorCreate(ceed, q_size, &q_vecs[i]));
1560d0321e0SJeremy L Thompson         break;
1570d0321e0SJeremy L Thompson       case CEED_EVAL_GRAD:
158b7453713SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis));
159b7453713SJeremy L Thompson         CeedCallBackend(CeedQFunctionFieldGetSize(qf_fields[i], &size));
1602b730f8bSJeremy L Thompson         CeedCallBackend(CeedBasisGetDimension(basis, &dim));
161b7453713SJeremy L Thompson         q_size = (CeedSize)num_elem * Q * size;
162b7453713SJeremy L Thompson         CeedCallBackend(CeedVectorCreate(ceed, q_size, &q_vecs[i]));
1630d0321e0SJeremy L Thompson         break;
1640d0321e0SJeremy L Thompson       case CEED_EVAL_WEIGHT:  // Only on input fields
165b7453713SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis));
166b7453713SJeremy L Thompson         q_size = (CeedSize)num_elem * Q;
167b7453713SJeremy L Thompson         CeedCallBackend(CeedVectorCreate(ceed, q_size, &q_vecs[i]));
168b7453713SJeremy L Thompson         CeedCallBackend(CeedBasisApply(basis, num_elem, CEED_NOTRANSPOSE, CEED_EVAL_WEIGHT, CEED_VECTOR_NONE, q_vecs[i]));
1690d0321e0SJeremy L Thompson         break;
1700d0321e0SJeremy L Thompson       case CEED_EVAL_DIV:
1710d0321e0SJeremy L Thompson         break;  // TODO: Not implemented
1720d0321e0SJeremy L Thompson       case CEED_EVAL_CURL:
1730d0321e0SJeremy L Thompson         break;  // TODO: Not implemented
1740d0321e0SJeremy L Thompson     }
1750d0321e0SJeremy L Thompson   }
1760d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1770d0321e0SJeremy L Thompson }
1780d0321e0SJeremy L Thompson 
1790d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
180ea61e9acSJeremy L Thompson // CeedOperator needs to connect all the named fields (be they active or passive) to the named inputs and outputs of its CeedQFunction.
1810d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
1820d0321e0SJeremy L Thompson static int CeedOperatorSetup_Hip(CeedOperator op) {
1830d0321e0SJeremy L Thompson   Ceed                ceed;
184b7453713SJeremy L Thompson   bool                is_setup_done;
185b7453713SJeremy L Thompson   CeedInt             Q, num_elem, num_input_fields, num_output_fields;
186b7453713SJeremy L Thompson   CeedQFunctionField *qf_input_fields, *qf_output_fields;
1870d0321e0SJeremy L Thompson   CeedQFunction       qf;
188b7453713SJeremy L Thompson   CeedOperatorField  *op_input_fields, *op_output_fields;
189b7453713SJeremy L Thompson   CeedOperator_Hip   *impl;
190b7453713SJeremy L Thompson 
191b7453713SJeremy L Thompson   CeedCallBackend(CeedOperatorIsSetupDone(op, &is_setup_done));
192b7453713SJeremy L Thompson   if (is_setup_done) return CEED_ERROR_SUCCESS;
193b7453713SJeremy L Thompson 
194b7453713SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
195b7453713SJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &impl));
1962b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
1972b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q));
198b7453713SJeremy L Thompson   CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem));
199b7453713SJeremy L Thompson   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
200b7453713SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
2010d0321e0SJeremy L Thompson 
2020d0321e0SJeremy L Thompson   // Allocate
203b7453713SJeremy L Thompson   CeedCallBackend(CeedCalloc(num_input_fields + num_output_fields, &impl->e_vecs));
2040d0321e0SJeremy L Thompson 
205b7453713SJeremy L Thompson   CeedCallBackend(CeedCalloc(CEED_FIELD_MAX, &impl->q_vecs_in));
206b7453713SJeremy L Thompson   CeedCallBackend(CeedCalloc(CEED_FIELD_MAX, &impl->q_vecs_out));
2070d0321e0SJeremy L Thompson 
208b7453713SJeremy L Thompson   impl->num_inputs  = num_input_fields;
209b7453713SJeremy L Thompson   impl->num_outputs = num_output_fields;
2100d0321e0SJeremy L Thompson 
211b7453713SJeremy L Thompson   // Set up infield and outfield e_vecs and q_vecs
2120d0321e0SJeremy L Thompson   // Infields
213b7453713SJeremy L Thompson   CeedCallBackend(CeedOperatorSetupFields_Hip(qf, op, true, impl->e_vecs, impl->q_vecs_in, 0, num_input_fields, Q, num_elem));
2140d0321e0SJeremy L Thompson 
2150d0321e0SJeremy L Thompson   // Outfields
216b7453713SJeremy L Thompson   CeedCallBackend(CeedOperatorSetupFields_Hip(qf, op, false, impl->e_vecs, impl->q_vecs_out, num_input_fields, num_output_fields, Q, num_elem));
2170d0321e0SJeremy L Thompson 
2182b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorSetSetupDone(op));
2190d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2200d0321e0SJeremy L Thompson }
2210d0321e0SJeremy L Thompson 
2220d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
2230d0321e0SJeremy L Thompson // Setup Operator Inputs
2240d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
225b7453713SJeremy L Thompson static inline int CeedOperatorSetupInputs_Hip(CeedInt num_input_fields, CeedQFunctionField *qf_input_fields, CeedOperatorField *op_input_fields,
226b7453713SJeremy L Thompson                                               CeedVector in_vec, const bool skip_active, CeedScalar *e_data[2 * CEED_FIELD_MAX],
227b7453713SJeremy L Thompson                                               CeedOperator_Hip *impl, CeedRequest *request) {
228b7453713SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
229b7453713SJeremy L Thompson     CeedEvalMode        e_mode;
2300d0321e0SJeremy L Thompson     CeedVector          vec;
231b7453713SJeremy L Thompson     CeedElemRestriction elem_rstr;
2320d0321e0SJeremy L Thompson 
2330d0321e0SJeremy L Thompson     // Get input vector
234b7453713SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec));
2350d0321e0SJeremy L Thompson     if (vec == CEED_VECTOR_ACTIVE) {
236b7453713SJeremy L Thompson       if (skip_active) continue;
237b7453713SJeremy L Thompson       else vec = in_vec;
2380d0321e0SJeremy L Thompson     }
2390d0321e0SJeremy L Thompson 
240b7453713SJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &e_mode));
241b7453713SJeremy L Thompson     if (e_mode == CEED_EVAL_WEIGHT) {  // Skip
2420d0321e0SJeremy L Thompson     } else {
2430d0321e0SJeremy L Thompson       // Get input vector
244b7453713SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec));
2450d0321e0SJeremy L Thompson       // Get input element restriction
246b7453713SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr));
247b7453713SJeremy L Thompson       if (vec == CEED_VECTOR_ACTIVE) vec = in_vec;
2480d0321e0SJeremy L Thompson       // Restrict, if necessary
249b7453713SJeremy L Thompson       if (!impl->e_vecs[i]) {
2500d0321e0SJeremy L Thompson         // No restriction for this field; read data directly from vec.
251b7453713SJeremy L Thompson         CeedCallBackend(CeedVectorGetArrayRead(vec, CEED_MEM_DEVICE, (const CeedScalar **)&e_data[i]));
2520d0321e0SJeremy L Thompson       } else {
253b7453713SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionApply(elem_rstr, CEED_NOTRANSPOSE, vec, impl->e_vecs[i], request));
2540d0321e0SJeremy L Thompson         // Get evec
255b7453713SJeremy L Thompson         CeedCallBackend(CeedVectorGetArrayRead(impl->e_vecs[i], CEED_MEM_DEVICE, (const CeedScalar **)&e_data[i]));
2560d0321e0SJeremy L Thompson       }
2570d0321e0SJeremy L Thompson     }
2580d0321e0SJeremy L Thompson   }
2590d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2600d0321e0SJeremy L Thompson }
2610d0321e0SJeremy L Thompson 
2620d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
2630d0321e0SJeremy L Thompson // Input Basis Action
2640d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
265b7453713SJeremy L Thompson static inline int CeedOperatorInputBasis_Hip(CeedInt num_elem, CeedQFunctionField *qf_input_fields, CeedOperatorField *op_input_fields,
266b7453713SJeremy L Thompson                                              CeedInt num_input_fields, const bool skip_active, CeedScalar *e_data[2 * CEED_FIELD_MAX],
2672b730f8bSJeremy L Thompson                                              CeedOperator_Hip *impl) {
268b7453713SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
269b7453713SJeremy L Thompson     CeedInt             elem_size, size;
270b7453713SJeremy L Thompson     CeedEvalMode        e_mode;
271b7453713SJeremy L Thompson     CeedElemRestriction elem_rstr;
2720d0321e0SJeremy L Thompson     CeedBasis           basis;
2730d0321e0SJeremy L Thompson 
2740d0321e0SJeremy L Thompson     // Skip active input
275b7453713SJeremy L Thompson     if (skip_active) {
2760d0321e0SJeremy L Thompson       CeedVector vec;
277b7453713SJeremy L Thompson 
278b7453713SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec));
2792b730f8bSJeremy L Thompson       if (vec == CEED_VECTOR_ACTIVE) continue;
2800d0321e0SJeremy L Thompson     }
281b7453713SJeremy L Thompson     // Get elem_size, e_mode, size
282b7453713SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr));
283b7453713SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
284b7453713SJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &e_mode));
285b7453713SJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetSize(qf_input_fields[i], &size));
2860d0321e0SJeremy L Thompson     // Basis action
287b7453713SJeremy L Thompson     switch (e_mode) {
2880d0321e0SJeremy L Thompson       case CEED_EVAL_NONE:
289b7453713SJeremy L Thompson         CeedCallBackend(CeedVectorSetArray(impl->q_vecs_in[i], CEED_MEM_DEVICE, CEED_USE_POINTER, e_data[i]));
2900d0321e0SJeremy L Thompson         break;
2910d0321e0SJeremy L Thompson       case CEED_EVAL_INTERP:
292b7453713SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
293b7453713SJeremy L Thompson         CeedCallBackend(CeedBasisApply(basis, num_elem, CEED_NOTRANSPOSE, CEED_EVAL_INTERP, impl->e_vecs[i], impl->q_vecs_in[i]));
2940d0321e0SJeremy L Thompson         break;
2950d0321e0SJeremy L Thompson       case CEED_EVAL_GRAD:
296b7453713SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
297b7453713SJeremy L Thompson         CeedCallBackend(CeedBasisApply(basis, num_elem, CEED_NOTRANSPOSE, CEED_EVAL_GRAD, impl->e_vecs[i], impl->q_vecs_in[i]));
2980d0321e0SJeremy L Thompson         break;
2990d0321e0SJeremy L Thompson       case CEED_EVAL_WEIGHT:
3000d0321e0SJeremy L Thompson         break;  // No action
3010d0321e0SJeremy L Thompson       case CEED_EVAL_DIV:
3020d0321e0SJeremy L Thompson         break;  // TODO: Not implemented
3030d0321e0SJeremy L Thompson       case CEED_EVAL_CURL:
3040d0321e0SJeremy L Thompson         break;  // TODO: Not implemented
3050d0321e0SJeremy L Thompson     }
3060d0321e0SJeremy L Thompson   }
3070d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3080d0321e0SJeremy L Thompson }
3090d0321e0SJeremy L Thompson 
3100d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
3110d0321e0SJeremy L Thompson // Restore Input Vectors
3120d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
313b7453713SJeremy L Thompson static inline int CeedOperatorRestoreInputs_Hip(CeedInt num_input_fields, CeedQFunctionField *qf_input_fields, CeedOperatorField *op_input_fields,
314b7453713SJeremy L Thompson                                                 const bool skip_active, CeedScalar *e_data[2 * CEED_FIELD_MAX], CeedOperator_Hip *impl) {
315b7453713SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
316b7453713SJeremy L Thompson     CeedEvalMode e_mode;
3170d0321e0SJeremy L Thompson     CeedVector   vec;
3180d0321e0SJeremy L Thompson     // Skip active input
319b7453713SJeremy L Thompson     if (skip_active) {
320b7453713SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec));
3212b730f8bSJeremy L Thompson       if (vec == CEED_VECTOR_ACTIVE) continue;
3220d0321e0SJeremy L Thompson     }
323b7453713SJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &e_mode));
324b7453713SJeremy L Thompson     if (e_mode == CEED_EVAL_WEIGHT) {  // Skip
3250d0321e0SJeremy L Thompson     } else {
326b7453713SJeremy L Thompson       if (!impl->e_vecs[i]) {  // This was a skip_restriction case
327b7453713SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec));
328b7453713SJeremy L Thompson         CeedCallBackend(CeedVectorRestoreArrayRead(vec, (const CeedScalar **)&e_data[i]));
3290d0321e0SJeremy L Thompson       } else {
330b7453713SJeremy L Thompson         CeedCallBackend(CeedVectorRestoreArrayRead(impl->e_vecs[i], (const CeedScalar **)&e_data[i]));
3310d0321e0SJeremy L Thompson       }
3320d0321e0SJeremy L Thompson     }
3330d0321e0SJeremy L Thompson   }
3340d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3350d0321e0SJeremy L Thompson }
3360d0321e0SJeremy L Thompson 
3370d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
3380d0321e0SJeremy L Thompson // Apply and add to output
3390d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
340b7453713SJeremy L Thompson static int CeedOperatorApplyAdd_Hip(CeedOperator op, CeedVector in_vec, CeedVector out_vec, CeedRequest *request) {
341b7453713SJeremy L Thompson   CeedInt             Q, num_elem, elem_size, num_input_fields, num_output_fields, size;
342b7453713SJeremy L Thompson   CeedScalar         *e_data[2 * CEED_FIELD_MAX] = {NULL};
343b7453713SJeremy L Thompson   CeedQFunctionField *qf_input_fields, *qf_output_fields;
3440d0321e0SJeremy L Thompson   CeedQFunction       qf;
345b7453713SJeremy L Thompson   CeedOperatorField  *op_input_fields, *op_output_fields;
346b7453713SJeremy L Thompson   CeedOperator_Hip   *impl;
347b7453713SJeremy L Thompson 
348b7453713SJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &impl));
3492b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
3502b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q));
351b7453713SJeremy L Thompson   CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem));
352b7453713SJeremy L Thompson   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
353b7453713SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
3540d0321e0SJeremy L Thompson 
3550d0321e0SJeremy L Thompson   // Setup
3562b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorSetup_Hip(op));
3570d0321e0SJeremy L Thompson 
3580d0321e0SJeremy L Thompson   // Input Evecs and Restriction
359b7453713SJeremy L Thompson   CeedCallBackend(CeedOperatorSetupInputs_Hip(num_input_fields, qf_input_fields, op_input_fields, in_vec, false, e_data, impl, request));
3600d0321e0SJeremy L Thompson 
3610d0321e0SJeremy L Thompson   // Input basis apply if needed
362b7453713SJeremy L Thompson   CeedCallBackend(CeedOperatorInputBasis_Hip(num_elem, qf_input_fields, op_input_fields, num_input_fields, false, e_data, impl));
3630d0321e0SJeremy L Thompson 
3640d0321e0SJeremy L Thompson   // Output pointers, as necessary
365b7453713SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
366b7453713SJeremy L Thompson     CeedEvalMode e_mode;
367b7453713SJeremy L Thompson 
368b7453713SJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &e_mode));
369b7453713SJeremy L Thompson     if (e_mode == CEED_EVAL_NONE) {
3700d0321e0SJeremy L Thompson       // Set the output Q-Vector to use the E-Vector data directly.
371b7453713SJeremy L Thompson       CeedCallBackend(CeedVectorGetArrayWrite(impl->e_vecs[i + impl->num_inputs], CEED_MEM_DEVICE, &e_data[i + num_input_fields]));
372b7453713SJeremy L Thompson       CeedCallBackend(CeedVectorSetArray(impl->q_vecs_out[i], CEED_MEM_DEVICE, CEED_USE_POINTER, e_data[i + num_input_fields]));
3730d0321e0SJeremy L Thompson     }
3740d0321e0SJeremy L Thompson   }
3750d0321e0SJeremy L Thompson 
3760d0321e0SJeremy L Thompson   // Q function
377b7453713SJeremy L Thompson   CeedCallBackend(CeedQFunctionApply(qf, num_elem * Q, impl->q_vecs_in, impl->q_vecs_out));
3780d0321e0SJeremy L Thompson 
3790d0321e0SJeremy L Thompson   // Output basis apply if needed
380b7453713SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
381b7453713SJeremy L Thompson     CeedEvalMode        e_mode;
382b7453713SJeremy L Thompson     CeedElemRestriction elem_rstr;
383b7453713SJeremy L Thompson     CeedBasis           basis;
384b7453713SJeremy L Thompson 
385b7453713SJeremy L Thompson     // Get elem_size, e_mode, size
386b7453713SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &elem_rstr));
387b7453713SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
388b7453713SJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &e_mode));
389b7453713SJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetSize(qf_output_fields[i], &size));
3900d0321e0SJeremy L Thompson     // Basis action
391b7453713SJeremy L Thompson     switch (e_mode) {
3920d0321e0SJeremy L Thompson       case CEED_EVAL_NONE:
3930d0321e0SJeremy L Thompson         break;
3940d0321e0SJeremy L Thompson       case CEED_EVAL_INTERP:
395b7453713SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
396b7453713SJeremy L Thompson         CeedCallBackend(CeedBasisApply(basis, num_elem, CEED_TRANSPOSE, CEED_EVAL_INTERP, impl->q_vecs_out[i], impl->e_vecs[i + impl->num_inputs]));
3970d0321e0SJeremy L Thompson         break;
3980d0321e0SJeremy L Thompson       case CEED_EVAL_GRAD:
399b7453713SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
400b7453713SJeremy L Thompson         CeedCallBackend(CeedBasisApply(basis, num_elem, CEED_TRANSPOSE, CEED_EVAL_GRAD, impl->q_vecs_out[i], impl->e_vecs[i + impl->num_inputs]));
4010d0321e0SJeremy L Thompson         break;
4020d0321e0SJeremy L Thompson       // LCOV_EXCL_START
4030d0321e0SJeremy L Thompson       case CEED_EVAL_WEIGHT: {
4040d0321e0SJeremy L Thompson         Ceed ceed;
405b7453713SJeremy L Thompson 
4062b730f8bSJeremy L Thompson         CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
4072b730f8bSJeremy L Thompson         return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_WEIGHT cannot be an output evaluation mode");
4080d0321e0SJeremy L Thompson         break;  // Should not occur
4090d0321e0SJeremy L Thompson       }
4100d0321e0SJeremy L Thompson       case CEED_EVAL_DIV:
4110d0321e0SJeremy L Thompson         break;  // TODO: Not implemented
4120d0321e0SJeremy L Thompson       case CEED_EVAL_CURL:
4130d0321e0SJeremy L Thompson         break;  // TODO: Not implemented
4140d0321e0SJeremy L Thompson                 // LCOV_EXCL_STOP
4150d0321e0SJeremy L Thompson     }
4160d0321e0SJeremy L Thompson   }
4170d0321e0SJeremy L Thompson 
4180d0321e0SJeremy L Thompson   // Output restriction
419b7453713SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
420b7453713SJeremy L Thompson     CeedEvalMode        e_mode;
421b7453713SJeremy L Thompson     CeedVector          vec;
422b7453713SJeremy L Thompson     CeedElemRestriction elem_rstr;
423b7453713SJeremy L Thompson 
4240d0321e0SJeremy L Thompson     // Restore evec
425b7453713SJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &e_mode));
426b7453713SJeremy L Thompson     if (e_mode == CEED_EVAL_NONE) {
427b7453713SJeremy L Thompson       CeedCallBackend(CeedVectorRestoreArray(impl->e_vecs[i + impl->num_inputs], &e_data[i + num_input_fields]));
4280d0321e0SJeremy L Thompson     }
4290d0321e0SJeremy L Thompson     // Get output vector
430b7453713SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[i], &vec));
4310d0321e0SJeremy L Thompson     // Restrict
432b7453713SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &elem_rstr));
4330d0321e0SJeremy L Thompson     // Active
434b7453713SJeremy L Thompson     if (vec == CEED_VECTOR_ACTIVE) vec = out_vec;
4350d0321e0SJeremy L Thompson 
436b7453713SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionApply(elem_rstr, CEED_TRANSPOSE, impl->e_vecs[i + impl->num_inputs], vec, request));
4370d0321e0SJeremy L Thompson   }
4380d0321e0SJeremy L Thompson 
4390d0321e0SJeremy L Thompson   // Restore input arrays
440b7453713SJeremy L Thompson   CeedCallBackend(CeedOperatorRestoreInputs_Hip(num_input_fields, qf_input_fields, op_input_fields, false, e_data, impl));
4410d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
4420d0321e0SJeremy L Thompson }
4430d0321e0SJeremy L Thompson 
4440d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
4450d0321e0SJeremy L Thompson // Core code for assembling linear QFunction
4460d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
4472b730f8bSJeremy L Thompson static inline int CeedOperatorLinearAssembleQFunctionCore_Hip(CeedOperator op, bool build_objects, CeedVector *assembled, CeedElemRestriction *rstr,
4480d0321e0SJeremy L Thompson                                                               CeedRequest *request) {
449b7453713SJeremy L Thompson   Ceed                ceed, ceed_parent;
450e984cf9aSJeremy L Thompson   CeedSize            q_size;
451b7453713SJeremy L Thompson   CeedInt             num_active_in, num_active_out, Q, num_elem, num_input_fields, num_output_fields, size;
452b7453713SJeremy L Thompson   CeedScalar         *assembled_array, *e_data[2 * CEED_FIELD_MAX] = {NULL};
453b7453713SJeremy L Thompson   CeedVector         *active_in;
454b7453713SJeremy L Thompson   CeedQFunctionField *qf_input_fields, *qf_output_fields;
455b7453713SJeremy L Thompson   CeedQFunction       qf;
456b7453713SJeremy L Thompson   CeedOperatorField  *op_input_fields, *op_output_fields;
457b7453713SJeremy L Thompson   CeedOperator_Hip   *impl;
458b7453713SJeremy L Thompson 
4592b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
460b7453713SJeremy L Thompson   CeedCallBackend(CeedOperatorGetFallbackParentCeed(op, &ceed_parent));
461e984cf9aSJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &impl));
462e984cf9aSJeremy L Thompson   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
463e984cf9aSJeremy L Thompson   CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q));
464b7453713SJeremy L Thompson   CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem));
465b7453713SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
466b7453713SJeremy L Thompson   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
467b7453713SJeremy L Thompson   active_in      = impl->qf_active_in;
468b7453713SJeremy L Thompson   num_active_in  = impl->num_active_in;
469b7453713SJeremy L Thompson   num_active_out = impl->num_active_out;
4700d0321e0SJeremy L Thompson 
4710d0321e0SJeremy L Thompson   // Setup
4722b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorSetup_Hip(op));
4730d0321e0SJeremy L Thompson 
4740d0321e0SJeremy L Thompson   // Input Evecs and Restriction
475b7453713SJeremy L Thompson   CeedCallBackend(CeedOperatorSetupInputs_Hip(num_input_fields, qf_input_fields, op_input_fields, NULL, true, e_data, impl, request));
4760d0321e0SJeremy L Thompson 
4770d0321e0SJeremy L Thompson   // Count number of active input fields
478b7453713SJeremy L Thompson   if (!num_active_in) {
479b7453713SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
480b7453713SJeremy L Thompson       CeedScalar *q_vec_array;
481b7453713SJeremy L Thompson       CeedVector  vec;
482b7453713SJeremy L Thompson 
4830d0321e0SJeremy L Thompson       // Get input vector
484b7453713SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec));
4850d0321e0SJeremy L Thompson       // Check if active input
4860d0321e0SJeremy L Thompson       if (vec == CEED_VECTOR_ACTIVE) {
487b7453713SJeremy L Thompson         CeedCallBackend(CeedQFunctionFieldGetSize(qf_input_fields[i], &size));
488b7453713SJeremy L Thompson         CeedCallBackend(CeedVectorSetValue(impl->q_vecs_in[i], 0.0));
489b7453713SJeremy L Thompson         CeedCallBackend(CeedVectorGetArray(impl->q_vecs_in[i], CEED_MEM_DEVICE, &q_vec_array));
490b7453713SJeremy L Thompson         CeedCallBackend(CeedRealloc(num_active_in + size, &active_in));
4910d0321e0SJeremy L Thompson         for (CeedInt field = 0; field < size; field++) {
492b7453713SJeremy L Thompson           q_size = (CeedSize)Q * num_elem;
493b7453713SJeremy L Thompson           CeedCallBackend(CeedVectorCreate(ceed, q_size, &active_in[num_active_in + field]));
494b7453713SJeremy L Thompson           CeedCallBackend(
495b7453713SJeremy L Thompson               CeedVectorSetArray(active_in[num_active_in + field], CEED_MEM_DEVICE, CEED_USE_POINTER, &q_vec_array[field * Q * num_elem]));
4960d0321e0SJeremy L Thompson         }
497b7453713SJeremy L Thompson         num_active_in += size;
498b7453713SJeremy L Thompson         CeedCallBackend(CeedVectorRestoreArray(impl->q_vecs_in[i], &q_vec_array));
4990d0321e0SJeremy L Thompson       }
5000d0321e0SJeremy L Thompson     }
501b7453713SJeremy L Thompson     impl->num_active_in = num_active_in;
502b7453713SJeremy L Thompson     impl->qf_active_in  = active_in;
5030d0321e0SJeremy L Thompson   }
5040d0321e0SJeremy L Thompson 
5050d0321e0SJeremy L Thompson   // Count number of active output fields
506b7453713SJeremy L Thompson   if (!num_active_out) {
507b7453713SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
508b7453713SJeremy L Thompson       CeedVector vec;
509b7453713SJeremy L Thompson 
5100d0321e0SJeremy L Thompson       // Get output vector
511b7453713SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[i], &vec));
5120d0321e0SJeremy L Thompson       // Check if active output
5130d0321e0SJeremy L Thompson       if (vec == CEED_VECTOR_ACTIVE) {
514b7453713SJeremy L Thompson         CeedCallBackend(CeedQFunctionFieldGetSize(qf_output_fields[i], &size));
515b7453713SJeremy L Thompson         num_active_out += size;
5160d0321e0SJeremy L Thompson       }
5170d0321e0SJeremy L Thompson     }
518b7453713SJeremy L Thompson     impl->num_active_out = num_active_out;
5190d0321e0SJeremy L Thompson   }
5200d0321e0SJeremy L Thompson 
5210d0321e0SJeremy L Thompson   // Check sizes
522b7453713SJeremy L Thompson   CeedCheck(num_active_in > 0 && num_active_out > 0, ceed, CEED_ERROR_BACKEND, "Cannot assemble QFunction without active inputs and outputs");
5230d0321e0SJeremy L Thompson 
5240d0321e0SJeremy L Thompson   // Build objects if needed
5250d0321e0SJeremy L Thompson   if (build_objects) {
5260d0321e0SJeremy L Thompson     // Create output restriction
527b7453713SJeremy L Thompson     CeedSize l_size     = (CeedSize)num_elem * Q * num_active_in * num_active_out;
528b7453713SJeremy L Thompson     CeedInt  strides[3] = {1, num_elem * Q, Q}; /* *NOPAD* */
529b7453713SJeremy L Thompson 
530b7453713SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionCreateStrided(ceed_parent, num_elem, Q, num_active_in * num_active_out,
531b7453713SJeremy L Thompson                                                      num_active_in * num_active_out * num_elem * Q, strides, rstr));
5320d0321e0SJeremy L Thompson     // Create assembled vector
533b7453713SJeremy L Thompson     CeedCallBackend(CeedVectorCreate(ceed_parent, l_size, assembled));
5340d0321e0SJeremy L Thompson   }
5352b730f8bSJeremy L Thompson   CeedCallBackend(CeedVectorSetValue(*assembled, 0.0));
536b7453713SJeremy L Thompson   CeedCallBackend(CeedVectorGetArray(*assembled, CEED_MEM_DEVICE, &assembled_array));
5370d0321e0SJeremy L Thompson 
5380d0321e0SJeremy L Thompson   // Input basis apply
539b7453713SJeremy L Thompson   CeedCallBackend(CeedOperatorInputBasis_Hip(num_elem, qf_input_fields, op_input_fields, num_input_fields, true, e_data, impl));
5400d0321e0SJeremy L Thompson 
5410d0321e0SJeremy L Thompson   // Assemble QFunction
542b7453713SJeremy L Thompson   for (CeedInt in = 0; in < num_active_in; in++) {
5430d0321e0SJeremy L Thompson     // Set Inputs
544b7453713SJeremy L Thompson     CeedCallBackend(CeedVectorSetValue(active_in[in], 1.0));
545b7453713SJeremy L Thompson     if (num_active_in > 1) {
546b7453713SJeremy L Thompson       CeedCallBackend(CeedVectorSetValue(active_in[(in + num_active_in - 1) % num_active_in], 0.0));
5470d0321e0SJeremy L Thompson     }
5480d0321e0SJeremy L Thompson     // Set Outputs
549b7453713SJeremy L Thompson     for (CeedInt out = 0; out < num_output_fields; out++) {
550b7453713SJeremy L Thompson       CeedVector vec;
551b7453713SJeremy L Thompson 
5520d0321e0SJeremy L Thompson       // Get output vector
553b7453713SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[out], &vec));
5540d0321e0SJeremy L Thompson       // Check if active output
5550d0321e0SJeremy L Thompson       if (vec == CEED_VECTOR_ACTIVE) {
556b7453713SJeremy L Thompson         CeedCallBackend(CeedVectorSetArray(impl->q_vecs_out[out], CEED_MEM_DEVICE, CEED_USE_POINTER, assembled_array));
557b7453713SJeremy L Thompson         CeedCallBackend(CeedQFunctionFieldGetSize(qf_output_fields[out], &size));
558b7453713SJeremy L Thompson         assembled_array += size * Q * num_elem;  // Advance the pointer by the size of the output
5590d0321e0SJeremy L Thompson       }
5600d0321e0SJeremy L Thompson     }
5610d0321e0SJeremy L Thompson     // Apply QFunction
562b7453713SJeremy L Thompson     CeedCallBackend(CeedQFunctionApply(qf, Q * num_elem, impl->q_vecs_in, impl->q_vecs_out));
5630d0321e0SJeremy L Thompson   }
5640d0321e0SJeremy L Thompson 
5650d0321e0SJeremy L Thompson   // Un-set output Qvecs to prevent accidental overwrite of Assembled
566b7453713SJeremy L Thompson   for (CeedInt out = 0; out < num_output_fields; out++) {
567b7453713SJeremy L Thompson     CeedVector vec;
568b7453713SJeremy L Thompson 
5690d0321e0SJeremy L Thompson     // Get output vector
570b7453713SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[out], &vec));
5710d0321e0SJeremy L Thompson     // Check if active output
5720d0321e0SJeremy L Thompson     if (vec == CEED_VECTOR_ACTIVE) {
573b7453713SJeremy L Thompson       CeedCallBackend(CeedVectorTakeArray(impl->q_vecs_out[out], CEED_MEM_DEVICE, NULL));
5740d0321e0SJeremy L Thompson     }
5750d0321e0SJeremy L Thompson   }
5760d0321e0SJeremy L Thompson 
5770d0321e0SJeremy L Thompson   // Restore input arrays
578b7453713SJeremy L Thompson   CeedCallBackend(CeedOperatorRestoreInputs_Hip(num_input_fields, qf_input_fields, op_input_fields, true, e_data, impl));
5790d0321e0SJeremy L Thompson 
5800d0321e0SJeremy L Thompson   // Restore output
581b7453713SJeremy L Thompson   CeedCallBackend(CeedVectorRestoreArray(*assembled, &assembled_array));
5820d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
5830d0321e0SJeremy L Thompson }
5840d0321e0SJeremy L Thompson 
5850d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
5860d0321e0SJeremy L Thompson // Assemble Linear QFunction
5870d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
5882b730f8bSJeremy L Thompson static int CeedOperatorLinearAssembleQFunction_Hip(CeedOperator op, CeedVector *assembled, CeedElemRestriction *rstr, CeedRequest *request) {
5892b730f8bSJeremy L Thompson   return CeedOperatorLinearAssembleQFunctionCore_Hip(op, true, assembled, rstr, request);
5900d0321e0SJeremy L Thompson }
5910d0321e0SJeremy L Thompson 
5920d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
593b2165e7aSSebastian Grimberg // Update Assembled Linear QFunction
5940d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
5952b730f8bSJeremy L Thompson static int CeedOperatorLinearAssembleQFunctionUpdate_Hip(CeedOperator op, CeedVector assembled, CeedElemRestriction rstr, CeedRequest *request) {
5962b730f8bSJeremy L Thompson   return CeedOperatorLinearAssembleQFunctionCore_Hip(op, false, &assembled, &rstr, request);
5970d0321e0SJeremy L Thompson }
5980d0321e0SJeremy L Thompson 
5990d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
6000d0321e0SJeremy L Thompson // Assemble diagonal setup
6010d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
602*506b1a0cSSebastian Grimberg static inline int CeedOperatorAssembleDiagonalSetup_Hip(CeedOperator op, CeedInt use_ceedsize_idx) {
6030d0321e0SJeremy L Thompson   Ceed                ceed;
604b7453713SJeremy L Thompson   char               *diagonal_kernel_path, *diagonal_kernel_source;
605b7453713SJeremy L Thompson   CeedInt             num_input_fields, num_output_fields, num_e_mode_in = 0, num_comp = 0, dim = 1, num_e_mode_out = 0;
606b7453713SJeremy L Thompson   CeedEvalMode       *e_mode_in = NULL, *e_mode_out = NULL;
607b7453713SJeremy L Thompson   CeedElemRestriction rstr_in = NULL, rstr_out = NULL;
608b7453713SJeremy L Thompson   CeedBasis           basis_in = NULL, basis_out = NULL;
609b7453713SJeremy L Thompson   CeedQFunctionField *qf_fields;
6100d0321e0SJeremy L Thompson   CeedQFunction       qf;
611b7453713SJeremy L Thompson   CeedOperatorField  *op_fields;
612b7453713SJeremy L Thompson   CeedOperator_Hip   *impl;
613b7453713SJeremy L Thompson 
614b7453713SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
6152b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
616b7453713SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetNumArgs(qf, &num_input_fields, &num_output_fields));
6170d0321e0SJeremy L Thompson 
6180d0321e0SJeremy L Thompson   // Determine active input basis
619b7453713SJeremy L Thompson   CeedCallBackend(CeedOperatorGetFields(op, NULL, &op_fields, NULL, NULL));
620b7453713SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_fields, NULL, NULL));
621b7453713SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
6220d0321e0SJeremy L Thompson     CeedVector vec;
623b7453713SJeremy L Thompson 
624b7453713SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec));
6250d0321e0SJeremy L Thompson     if (vec == CEED_VECTOR_ACTIVE) {
626b7453713SJeremy L Thompson       CeedEvalMode        e_mode;
6270d0321e0SJeremy L Thompson       CeedElemRestriction rstr;
628b7453713SJeremy L Thompson 
629b7453713SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis_in));
630b7453713SJeremy L Thompson       CeedCallBackend(CeedBasisGetNumComponents(basis_in, &num_comp));
631b7453713SJeremy L Thompson       CeedCallBackend(CeedBasisGetDimension(basis_in, &dim));
632b7453713SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_fields[i], &rstr));
633b7453713SJeremy L Thompson       CeedCheck(!rstr_in || rstr_in == rstr, ceed, CEED_ERROR_BACKEND,
6346574a04fSJeremy L Thompson                 "Backend does not implement multi-field non-composite operator diagonal assembly");
635b7453713SJeremy L Thompson       rstr_in = rstr;
636b7453713SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &e_mode));
637b7453713SJeremy L Thompson       switch (e_mode) {
6380d0321e0SJeremy L Thompson         case CEED_EVAL_NONE:
6390d0321e0SJeremy L Thompson         case CEED_EVAL_INTERP:
640b7453713SJeremy L Thompson           CeedCallBackend(CeedRealloc(num_e_mode_in + 1, &e_mode_in));
641b7453713SJeremy L Thompson           e_mode_in[num_e_mode_in] = e_mode;
642b7453713SJeremy L Thompson           num_e_mode_in += 1;
6430d0321e0SJeremy L Thompson           break;
6440d0321e0SJeremy L Thompson         case CEED_EVAL_GRAD:
645b7453713SJeremy L Thompson           CeedCallBackend(CeedRealloc(num_e_mode_in + dim, &e_mode_in));
646b7453713SJeremy L Thompson           for (CeedInt d = 0; d < dim; d++) e_mode_in[num_e_mode_in + d] = e_mode;
647b7453713SJeremy L Thompson           num_e_mode_in += dim;
6480d0321e0SJeremy L Thompson           break;
6490d0321e0SJeremy L Thompson         case CEED_EVAL_WEIGHT:
6500d0321e0SJeremy L Thompson         case CEED_EVAL_DIV:
6510d0321e0SJeremy L Thompson         case CEED_EVAL_CURL:
6520d0321e0SJeremy L Thompson           break;  // Caught by QF Assembly
6530d0321e0SJeremy L Thompson       }
6540d0321e0SJeremy L Thompson     }
6550d0321e0SJeremy L Thompson   }
6560d0321e0SJeremy L Thompson 
6570d0321e0SJeremy L Thompson   // Determine active output basis
658b7453713SJeremy L Thompson   CeedCallBackend(CeedOperatorGetFields(op, NULL, NULL, NULL, &op_fields));
659b7453713SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qf_fields));
660b7453713SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
6610d0321e0SJeremy L Thompson     CeedVector vec;
662b7453713SJeremy L Thompson 
663b7453713SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec));
6640d0321e0SJeremy L Thompson     if (vec == CEED_VECTOR_ACTIVE) {
665b7453713SJeremy L Thompson       CeedEvalMode        e_mode;
6660d0321e0SJeremy L Thompson       CeedElemRestriction rstr;
667b7453713SJeremy L Thompson 
668b7453713SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis_out));
669b7453713SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_fields[i], &rstr));
670b7453713SJeremy L Thompson       CeedCheck(!rstr_out || rstr_out == rstr, ceed, CEED_ERROR_BACKEND,
6716574a04fSJeremy L Thompson                 "Backend does not implement multi-field non-composite operator diagonal assembly");
672b7453713SJeremy L Thompson       rstr_out = rstr;
673b7453713SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &e_mode));
674b7453713SJeremy L Thompson       switch (e_mode) {
6750d0321e0SJeremy L Thompson         case CEED_EVAL_NONE:
6760d0321e0SJeremy L Thompson         case CEED_EVAL_INTERP:
677b7453713SJeremy L Thompson           CeedCallBackend(CeedRealloc(num_e_mode_out + 1, &e_mode_out));
678b7453713SJeremy L Thompson           e_mode_out[num_e_mode_out] = e_mode;
679b7453713SJeremy L Thompson           num_e_mode_out += 1;
6800d0321e0SJeremy L Thompson           break;
6810d0321e0SJeremy L Thompson         case CEED_EVAL_GRAD:
682b7453713SJeremy L Thompson           CeedCallBackend(CeedRealloc(num_e_mode_out + dim, &e_mode_out));
683b7453713SJeremy L Thompson           for (CeedInt d = 0; d < dim; d++) e_mode_out[num_e_mode_out + d] = e_mode;
684b7453713SJeremy L Thompson           num_e_mode_out += dim;
6850d0321e0SJeremy L Thompson           break;
6860d0321e0SJeremy L Thompson         case CEED_EVAL_WEIGHT:
6870d0321e0SJeremy L Thompson         case CEED_EVAL_DIV:
6880d0321e0SJeremy L Thompson         case CEED_EVAL_CURL:
6890d0321e0SJeremy L Thompson           break;  // Caught by QF Assembly
6900d0321e0SJeremy L Thompson       }
6910d0321e0SJeremy L Thompson     }
6920d0321e0SJeremy L Thompson   }
6930d0321e0SJeremy L Thompson 
6940d0321e0SJeremy L Thompson   // Operator data struct
6952b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &impl));
6962b730f8bSJeremy L Thompson   CeedCallBackend(CeedCalloc(1, &impl->diag));
6970d0321e0SJeremy L Thompson   CeedOperatorDiag_Hip *diag = impl->diag;
698b7453713SJeremy L Thompson 
699b7453713SJeremy L Thompson   diag->basis_in       = basis_in;
700b7453713SJeremy L Thompson   diag->basis_out      = basis_out;
701b7453713SJeremy L Thompson   diag->h_e_mode_in    = e_mode_in;
702b7453713SJeremy L Thompson   diag->h_e_mode_out   = e_mode_out;
703b7453713SJeremy L Thompson   diag->num_e_mode_in  = num_e_mode_in;
704b7453713SJeremy L Thompson   diag->num_e_mode_out = num_e_mode_out;
7050d0321e0SJeremy L Thompson 
7060d0321e0SJeremy L Thompson   // Assemble kernel
7072b730f8bSJeremy L Thompson   CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/hip/hip-ref-operator-assemble-diagonal.h", &diagonal_kernel_path));
70823d4529eSJeremy L Thompson   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Diagonal Assembly Kernel Source -----\n");
7092b730f8bSJeremy L Thompson   CeedCallBackend(CeedLoadSourceToBuffer(ceed, diagonal_kernel_path, &diagonal_kernel_source));
71023d4529eSJeremy L Thompson   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Diagonal Assembly Source Complete! -----\n");
711b7453713SJeremy L Thompson   CeedInt num_modes, num_qpts;
712b7453713SJeremy L Thompson   CeedCallBackend(CeedBasisGetNumNodes(basis_in, &num_modes));
713b7453713SJeremy L Thompson   CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis_in, &num_qpts));
714b7453713SJeremy L Thompson   diag->num_modes = num_modes;
715b7453713SJeremy L Thompson   CeedCallBackend(CeedCompile_Hip(ceed, diagonal_kernel_source, &diag->module, 6, "NUMEMODEIN", num_e_mode_in, "NUMEMODEOUT", num_e_mode_out,
716b7453713SJeremy L Thompson                                   "NNODES", num_modes, "NQPTS", num_qpts, "NCOMP", num_comp, "CEEDSIZE", use_ceedsize_idx));
717eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedGetKernel_Hip(ceed, diag->module, "linearDiagonal", &diag->linearDiagonal));
718eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedGetKernel_Hip(ceed, diag->module, "linearPointBlockDiagonal", &diag->linearPointBlock));
7192b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&diagonal_kernel_path));
7202b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&diagonal_kernel_source));
7210d0321e0SJeremy L Thompson 
7220d0321e0SJeremy L Thompson   // Basis matrices
723b7453713SJeremy L Thompson   const CeedInt     q_bytes      = num_qpts * sizeof(CeedScalar);
724b7453713SJeremy L Thompson   const CeedInt     interp_bytes = q_bytes * num_modes;
725b7453713SJeremy L Thompson   const CeedInt     grad_bytes   = q_bytes * num_modes * dim;
726b7453713SJeremy L Thompson   const CeedInt     e_mode_bytes = sizeof(CeedEvalMode);
727b7453713SJeremy L Thompson   const CeedScalar *interp_in, *interp_out, *grad_in, *grad_out;
7280d0321e0SJeremy L Thompson 
7290d0321e0SJeremy L Thompson   // CEED_EVAL_NONE
7300d0321e0SJeremy L Thompson   CeedScalar *identity     = NULL;
731b7453713SJeremy L Thompson   bool        is_eval_none = false;
732b7453713SJeremy L Thompson 
733b7453713SJeremy L Thompson   for (CeedInt i = 0; i < num_e_mode_in; i++) is_eval_none = is_eval_none || (e_mode_in[i] == CEED_EVAL_NONE);
734b7453713SJeremy L Thompson   for (CeedInt i = 0; i < num_e_mode_out; i++) is_eval_none = is_eval_none || (e_mode_out[i] == CEED_EVAL_NONE);
735b7453713SJeremy L Thompson   if (is_eval_none) {
736b7453713SJeremy L Thompson     CeedCallBackend(CeedCalloc(num_qpts * num_modes, &identity));
737b7453713SJeremy L Thompson     for (CeedInt i = 0; i < (num_modes < num_qpts ? num_modes : num_qpts); i++) identity[i * num_modes + i] = 1.0;
738b7453713SJeremy L Thompson     CeedCallHip(ceed, hipMalloc((void **)&diag->d_identity, interp_bytes));
739b7453713SJeremy L Thompson     CeedCallHip(ceed, hipMemcpy(diag->d_identity, identity, interp_bytes, hipMemcpyHostToDevice));
7400d0321e0SJeremy L Thompson   }
7410d0321e0SJeremy L Thompson 
7420d0321e0SJeremy L Thompson   // CEED_EVAL_INTERP
743b7453713SJeremy L Thompson   CeedCallBackend(CeedBasisGetInterp(basis_in, &interp_in));
744b7453713SJeremy L Thompson   CeedCallHip(ceed, hipMalloc((void **)&diag->d_interp_in, interp_bytes));
745b7453713SJeremy L Thompson   CeedCallHip(ceed, hipMemcpy(diag->d_interp_in, interp_in, interp_bytes, hipMemcpyHostToDevice));
746b7453713SJeremy L Thompson   CeedCallBackend(CeedBasisGetInterp(basis_out, &interp_out));
747b7453713SJeremy L Thompson   CeedCallHip(ceed, hipMalloc((void **)&diag->d_interp_out, interp_bytes));
748b7453713SJeremy L Thompson   CeedCallHip(ceed, hipMemcpy(diag->d_interp_out, interp_out, interp_bytes, hipMemcpyHostToDevice));
7490d0321e0SJeremy L Thompson 
7500d0321e0SJeremy L Thompson   // CEED_EVAL_GRAD
751b7453713SJeremy L Thompson   CeedCallBackend(CeedBasisGetGrad(basis_in, &grad_in));
752b7453713SJeremy L Thompson   CeedCallHip(ceed, hipMalloc((void **)&diag->d_grad_in, grad_bytes));
753b7453713SJeremy L Thompson   CeedCallHip(ceed, hipMemcpy(diag->d_grad_in, grad_in, grad_bytes, hipMemcpyHostToDevice));
754b7453713SJeremy L Thompson   CeedCallBackend(CeedBasisGetGrad(basis_out, &grad_out));
755b7453713SJeremy L Thompson   CeedCallHip(ceed, hipMalloc((void **)&diag->d_grad_out, grad_bytes));
756b7453713SJeremy L Thompson   CeedCallHip(ceed, hipMemcpy(diag->d_grad_out, grad_out, grad_bytes, hipMemcpyHostToDevice));
7570d0321e0SJeremy L Thompson 
758b7453713SJeremy L Thompson   // Arrays of e_modes
759b7453713SJeremy L Thompson   CeedCallHip(ceed, hipMalloc((void **)&diag->d_e_mode_in, num_e_mode_in * e_mode_bytes));
760b7453713SJeremy L Thompson   CeedCallHip(ceed, hipMemcpy(diag->d_e_mode_in, e_mode_in, num_e_mode_in * e_mode_bytes, hipMemcpyHostToDevice));
761b7453713SJeremy L Thompson   CeedCallHip(ceed, hipMalloc((void **)&diag->d_e_mode_out, num_e_mode_out * e_mode_bytes));
762b7453713SJeremy L Thompson   CeedCallHip(ceed, hipMemcpy(diag->d_e_mode_out, e_mode_out, num_e_mode_out * e_mode_bytes, hipMemcpyHostToDevice));
7630d0321e0SJeremy L Thompson 
7640d0321e0SJeremy L Thompson   // Restriction
765b7453713SJeremy L Thompson   diag->diag_rstr = rstr_out;
7660d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
7670d0321e0SJeremy L Thompson }
7680d0321e0SJeremy L Thompson 
7690d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
7700d0321e0SJeremy L Thompson // Assemble diagonal common code
7710d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
772b7453713SJeremy L Thompson static inline int CeedOperatorAssembleDiagonalCore_Hip(CeedOperator op, CeedVector assembled, CeedRequest *request, const bool is_point_block) {
7730d0321e0SJeremy L Thompson   Ceed                ceed;
774b7453713SJeremy L Thompson   CeedSize            assembled_length = 0, assembled_qf_length = 0;
775b7453713SJeremy L Thompson   CeedInt             use_ceedsize_idx = 0, num_elem;
776b7453713SJeremy L Thompson   CeedScalar         *elem_diag_array;
777b7453713SJeremy L Thompson   const CeedScalar   *assembled_qf_array;
778b7453713SJeremy L Thompson   CeedVector          assembled_qf = NULL;
779b7453713SJeremy L Thompson   CeedElemRestriction rstr         = NULL;
7800d0321e0SJeremy L Thompson   CeedOperator_Hip   *impl;
781b7453713SJeremy L Thompson 
782b7453713SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
7832b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &impl));
7840d0321e0SJeremy L Thompson 
7850d0321e0SJeremy L Thompson   // Assemble QFunction
786b7453713SJeremy L Thompson   CeedCallBackend(CeedOperatorLinearAssembleQFunctionBuildOrUpdate(op, &assembled_qf, &rstr, request));
7872b730f8bSJeremy L Thompson   CeedCallBackend(CeedElemRestrictionDestroy(&rstr));
7880d0321e0SJeremy L Thompson 
7899330daecSnbeams   CeedCallBackend(CeedVectorGetLength(assembled, &assembled_length));
790b7453713SJeremy L Thompson   CeedCallBackend(CeedVectorGetLength(assembled_qf, &assembled_qf_length));
791b7453713SJeremy L Thompson   if ((assembled_length > INT_MAX) || (assembled_qf_length > INT_MAX)) use_ceedsize_idx = 1;
7929330daecSnbeams 
7930d0321e0SJeremy L Thompson   // Setup
794*506b1a0cSSebastian Grimberg   if (!impl->diag) CeedCallBackend(CeedOperatorAssembleDiagonalSetup_Hip(op, use_ceedsize_idx));
7950d0321e0SJeremy L Thompson   CeedOperatorDiag_Hip *diag = impl->diag;
796*506b1a0cSSebastian Grimberg 
7970d0321e0SJeremy L Thompson   assert(diag != NULL);
7980d0321e0SJeremy L Thompson 
7990d0321e0SJeremy L Thompson   // Restriction
800b7453713SJeremy L Thompson   if (is_point_block && !diag->point_block_diag_rstr) {
801*506b1a0cSSebastian Grimberg     CeedCallBackend(CeedOperatorCreateActivePointBlockRestriction(diag->diag_rstr, &diag->point_block_diag_rstr));
8020d0321e0SJeremy L Thompson   }
803b7453713SJeremy L Thompson   CeedElemRestriction diag_rstr = is_point_block ? diag->point_block_diag_rstr : diag->diag_rstr;
8040d0321e0SJeremy L Thompson 
8050d0321e0SJeremy L Thompson   // Create diagonal vector
806b7453713SJeremy L Thompson   CeedVector elem_diag = is_point_block ? diag->point_block_elem_diag : diag->elem_diag;
807b7453713SJeremy L Thompson 
808b7453713SJeremy L Thompson   if (!elem_diag) {
809b7453713SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionCreateVector(diag_rstr, NULL, &elem_diag));
810b7453713SJeremy L Thompson     if (is_point_block) diag->point_block_elem_diag = elem_diag;
811b7453713SJeremy L Thompson     else diag->elem_diag = elem_diag;
8120d0321e0SJeremy L Thompson   }
813b7453713SJeremy L Thompson   CeedCallBackend(CeedVectorSetValue(elem_diag, 0.0));
8140d0321e0SJeremy L Thompson 
8150d0321e0SJeremy L Thompson   // Assemble element operator diagonals
816b7453713SJeremy L Thompson   CeedCallBackend(CeedVectorGetArray(elem_diag, CEED_MEM_DEVICE, &elem_diag_array));
817b7453713SJeremy L Thompson   CeedCallBackend(CeedVectorGetArrayRead(assembled_qf, CEED_MEM_DEVICE, &assembled_qf_array));
818b7453713SJeremy L Thompson   CeedCallBackend(CeedElemRestrictionGetNumElements(diag_rstr, &num_elem));
8190d0321e0SJeremy L Thompson 
8200d0321e0SJeremy L Thompson   // Compute the diagonal of B^T D B
821b7453713SJeremy L Thompson   int   elem_per_block = 1;
822b7453713SJeremy L Thompson   int   grid           = num_elem / elem_per_block + ((num_elem / elem_per_block * elem_per_block < num_elem) ? 1 : 0);
823b7453713SJeremy L Thompson   void *args[]         = {(void *)&num_elem, &diag->d_identity,  &diag->d_interp_in,  &diag->d_grad_in,    &diag->d_interp_out,
824b7453713SJeremy L Thompson                           &diag->d_grad_out, &diag->d_e_mode_in, &diag->d_e_mode_out, &assembled_qf_array, &elem_diag_array};
825b7453713SJeremy L Thompson 
826b7453713SJeremy L Thompson   if (is_point_block) {
827b7453713SJeremy L Thompson     CeedCallBackend(CeedRunKernelDim_Hip(ceed, diag->linearPointBlock, grid, diag->num_modes, 1, elem_per_block, args));
8280d0321e0SJeremy L Thompson   } else {
829b7453713SJeremy L Thompson     CeedCallBackend(CeedRunKernelDim_Hip(ceed, diag->linearDiagonal, grid, diag->num_modes, 1, elem_per_block, args));
8300d0321e0SJeremy L Thompson   }
8310d0321e0SJeremy L Thompson 
8320d0321e0SJeremy L Thompson   // Restore arrays
833b7453713SJeremy L Thompson   CeedCallBackend(CeedVectorRestoreArray(elem_diag, &elem_diag_array));
834b7453713SJeremy L Thompson   CeedCallBackend(CeedVectorRestoreArrayRead(assembled_qf, &assembled_qf_array));
8350d0321e0SJeremy L Thompson 
8360d0321e0SJeremy L Thompson   // Assemble local operator diagonal
837b7453713SJeremy L Thompson   CeedCallBackend(CeedElemRestrictionApply(diag_rstr, CEED_TRANSPOSE, elem_diag, assembled, request));
8380d0321e0SJeremy L Thompson 
8390d0321e0SJeremy L Thompson   // Cleanup
840b7453713SJeremy L Thompson   CeedCallBackend(CeedVectorDestroy(&assembled_qf));
8410d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
8420d0321e0SJeremy L Thompson }
8430d0321e0SJeremy L Thompson 
8440d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
8450d0321e0SJeremy L Thompson // Assemble Linear Diagonal
8460d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
8472b730f8bSJeremy L Thompson static int CeedOperatorLinearAssembleAddDiagonal_Hip(CeedOperator op, CeedVector assembled, CeedRequest *request) {
8482b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorAssembleDiagonalCore_Hip(op, assembled, request, false));
8496aa95790SJeremy L Thompson   return CEED_ERROR_SUCCESS;
8500d0321e0SJeremy L Thompson }
8510d0321e0SJeremy L Thompson 
8520d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
8530d0321e0SJeremy L Thompson // Assemble Linear Point Block Diagonal
8540d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
8552b730f8bSJeremy L Thompson static int CeedOperatorLinearAssembleAddPointBlockDiagonal_Hip(CeedOperator op, CeedVector assembled, CeedRequest *request) {
8562b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorAssembleDiagonalCore_Hip(op, assembled, request, true));
8576aa95790SJeremy L Thompson   return CEED_ERROR_SUCCESS;
8580d0321e0SJeremy L Thompson }
8590d0321e0SJeremy L Thompson 
860a835093fSnbeams //------------------------------------------------------------------------------
861a835093fSnbeams // Single operator assembly setup
862a835093fSnbeams //------------------------------------------------------------------------------
8639330daecSnbeams static int CeedSingleOperatorAssembleSetup_Hip(CeedOperator op, CeedInt use_ceedsize_idx) {
864a835093fSnbeams   Ceed    ceed;
865b7453713SJeremy L Thompson   CeedInt num_input_fields, num_output_fields, num_e_mode_in = 0, dim = 1, num_B_in_mats_to_load = 0, size_B_in = 0, num_qpts = 0, elem_size = 0,
866b7453713SJeremy L Thompson                                                num_e_mode_out = 0, num_B_out_mats_to_load = 0, size_B_out = 0, num_elem, num_comp;
867b7453713SJeremy L Thompson   CeedEvalMode       *eval_mode_in = NULL, *eval_mode_out = NULL;
868b7453713SJeremy L Thompson   CeedElemRestriction rstr_in = NULL, rstr_out = NULL;
869b7453713SJeremy L Thompson   CeedBasis           basis_in = NULL, basis_out = NULL;
870b7453713SJeremy L Thompson   CeedQFunctionField *qf_fields;
871b7453713SJeremy L Thompson   CeedQFunction       qf;
872b7453713SJeremy L Thompson   CeedOperatorField  *input_fields, *output_fields;
873a835093fSnbeams   CeedOperator_Hip   *impl;
874b7453713SJeremy L Thompson 
875b7453713SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
8762b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &impl));
877a835093fSnbeams 
878a835093fSnbeams   // Get intput and output fields
8792b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &input_fields, &num_output_fields, &output_fields));
880a835093fSnbeams 
881a835093fSnbeams   // Determine active input basis eval mode
8822b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
8832b730f8bSJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_fields, NULL, NULL));
884a835093fSnbeams   // Note that the kernel will treat each dimension of a gradient action separately;
885b7453713SJeremy L Thompson   // i.e., when an active input has a CEED_EVAL_GRAD mode, num_e_mode_in will increment by dim.
886ea61e9acSJeremy L Thompson   // However, for the purposes of loading the B matrices, it will be treated as one mode, and we will load/copy the entire gradient matrix at once, so
887a835093fSnbeams   // num_B_in_mats_to_load will be incremented by 1.
888a835093fSnbeams   for (CeedInt i = 0; i < num_input_fields; i++) {
889a835093fSnbeams     CeedVector vec;
890b7453713SJeremy L Thompson 
8912b730f8bSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetVector(input_fields[i], &vec));
892a835093fSnbeams     if (vec == CEED_VECTOR_ACTIVE) {
893b7453713SJeremy L Thompson       CeedEvalMode eval_mode;
894b7453713SJeremy L Thompson 
8952b730f8bSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetBasis(input_fields[i], &basis_in));
8962b730f8bSJeremy L Thompson       CeedCallBackend(CeedBasisGetDimension(basis_in, &dim));
897b7453713SJeremy L Thompson       CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis_in, &num_qpts));
8982b730f8bSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetElemRestriction(input_fields[i], &rstr_in));
899b7453713SJeremy L Thompson       CeedCallBackend(CeedElemRestrictionGetElementSize(rstr_in, &elem_size));
9002b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode));
901a835093fSnbeams       if (eval_mode != CEED_EVAL_NONE) {
9022b730f8bSJeremy L Thompson         CeedCallBackend(CeedRealloc(num_B_in_mats_to_load + 1, &eval_mode_in));
903a835093fSnbeams         eval_mode_in[num_B_in_mats_to_load] = eval_mode;
904a835093fSnbeams         num_B_in_mats_to_load += 1;
905a835093fSnbeams         if (eval_mode == CEED_EVAL_GRAD) {
906b7453713SJeremy L Thompson           num_e_mode_in += dim;
907b7453713SJeremy L Thompson           size_B_in += dim * elem_size * num_qpts;
908a835093fSnbeams         } else {
909b7453713SJeremy L Thompson           num_e_mode_in += 1;
910b7453713SJeremy L Thompson           size_B_in += elem_size * num_qpts;
911a835093fSnbeams         }
912a835093fSnbeams       }
913a835093fSnbeams     }
914a835093fSnbeams   }
915a835093fSnbeams 
916a835093fSnbeams   // Determine active output basis; basis_out and rstr_out only used if same as input, TODO
9172b730f8bSJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qf_fields));
918a835093fSnbeams   for (CeedInt i = 0; i < num_output_fields; i++) {
919a835093fSnbeams     CeedVector vec;
920b7453713SJeremy L Thompson 
9212b730f8bSJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetVector(output_fields[i], &vec));
922a835093fSnbeams     if (vec == CEED_VECTOR_ACTIVE) {
923b7453713SJeremy L Thompson       CeedEvalMode eval_mode;
924b7453713SJeremy L Thompson 
9252b730f8bSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetBasis(output_fields[i], &basis_out));
9262b730f8bSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetElemRestriction(output_fields[i], &rstr_out));
9276574a04fSJeremy L Thompson       CeedCheck(!rstr_out || rstr_out == rstr_in, ceed, CEED_ERROR_BACKEND, "Backend does not implement multi-field non-composite operator assembly");
9282b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode));
929a835093fSnbeams       if (eval_mode != CEED_EVAL_NONE) {
9302b730f8bSJeremy L Thompson         CeedCallBackend(CeedRealloc(num_B_out_mats_to_load + 1, &eval_mode_out));
931a835093fSnbeams         eval_mode_out[num_B_out_mats_to_load] = eval_mode;
932a835093fSnbeams         num_B_out_mats_to_load += 1;
933a835093fSnbeams         if (eval_mode == CEED_EVAL_GRAD) {
934b7453713SJeremy L Thompson           num_e_mode_out += dim;
935b7453713SJeremy L Thompson           size_B_out += dim * elem_size * num_qpts;
936a835093fSnbeams         } else {
937b7453713SJeremy L Thompson           num_e_mode_out += 1;
938b7453713SJeremy L Thompson           size_B_out += elem_size * num_qpts;
939a835093fSnbeams         }
940a835093fSnbeams       }
941a835093fSnbeams     }
942a835093fSnbeams   }
943a835093fSnbeams 
944b7453713SJeremy L Thompson   CeedCheck(num_e_mode_in > 0 && num_e_mode_out > 0, ceed, CEED_ERROR_UNSUPPORTED, "Cannot assemble operator without inputs/outputs");
945a835093fSnbeams 
946b7453713SJeremy L Thompson   CeedCallBackend(CeedElemRestrictionGetNumElements(rstr_in, &num_elem));
947b7453713SJeremy L Thompson   CeedCallBackend(CeedElemRestrictionGetNumComponents(rstr_in, &num_comp));
948a835093fSnbeams 
9492b730f8bSJeremy L Thompson   CeedCallBackend(CeedCalloc(1, &impl->asmb));
950a835093fSnbeams   CeedOperatorAssemble_Hip *asmb = impl->asmb;
951b7453713SJeremy L Thompson   asmb->num_elem                 = num_elem;
952a835093fSnbeams 
953a835093fSnbeams   // Compile kernels
954b7453713SJeremy L Thompson   int elem_per_block   = 1;
955b7453713SJeremy L Thompson   asmb->elem_per_block = elem_per_block;
956b7453713SJeremy L Thompson   CeedInt block_size   = elem_size * elem_size * elem_per_block;
95707b31e0eSJeremy L Thompson   char   *assembly_kernel_path, *assembly_kernel_source;
9582b730f8bSJeremy L Thompson   CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/hip/hip-ref-operator-assemble.h", &assembly_kernel_path));
95923d4529eSJeremy L Thompson   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Assembly Kernel Source -----\n");
9602b730f8bSJeremy L Thompson   CeedCallBackend(CeedLoadSourceToBuffer(ceed, assembly_kernel_path, &assembly_kernel_source));
96123d4529eSJeremy L Thompson   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Assembly Source Complete! -----\n");
96207b31e0eSJeremy L Thompson   bool fallback = block_size > 1024;
96307b31e0eSJeremy L Thompson   if (fallback) {  // Use fallback kernel with 1D threadblock
964b7453713SJeremy L Thompson     block_size         = elem_size * elem_per_block;
965b7453713SJeremy L Thompson     asmb->block_size_x = elem_size;
96659ad764aSnbeams     asmb->block_size_y = 1;
96759ad764aSnbeams   } else {  // Use kernel with 2D threadblock
968b7453713SJeremy L Thompson     asmb->block_size_x = elem_size;
969b7453713SJeremy L Thompson     asmb->block_size_y = elem_size;
97007b31e0eSJeremy L Thompson   }
971b7453713SJeremy L Thompson   CeedCallBackend(CeedCompile_Hip(ceed, assembly_kernel_source, &asmb->module, 8, "NELEM", num_elem, "NUMEMODEIN", num_e_mode_in, "NUMEMODEOUT",
972b7453713SJeremy L Thompson                                   num_e_mode_out, "NQPTS", num_qpts, "NNODES", elem_size, "BLOCK_SIZE", block_size, "NCOMP", num_comp, "CEEDSIZE",
9739330daecSnbeams                                   use_ceedsize_idx));
974eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedGetKernel_Hip(ceed, asmb->module, fallback ? "linearAssembleFallback" : "linearAssemble", &asmb->linearAssemble));
9752b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&assembly_kernel_path));
9762b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&assembly_kernel_source));
977a835093fSnbeams 
978a835093fSnbeams   // Build 'full' B matrices (not 1D arrays used for tensor-product matrices)
979a835093fSnbeams   const CeedScalar *interp_in, *grad_in;
9802b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetInterp(basis_in, &interp_in));
9812b730f8bSJeremy L Thompson   CeedCallBackend(CeedBasisGetGrad(basis_in, &grad_in));
982a835093fSnbeams 
983a835093fSnbeams   // Load into B_in, in order that they will be used in eval_mode
984b7453713SJeremy L Thompson   const CeedInt in_bytes  = size_B_in * sizeof(CeedScalar);
985a835093fSnbeams   CeedInt       mat_start = 0;
986b7453713SJeremy L Thompson 
987b7453713SJeremy L Thompson   CeedCallHip(ceed, hipMalloc((void **)&asmb->d_B_in, in_bytes));
988a835093fSnbeams   for (int i = 0; i < num_B_in_mats_to_load; i++) {
989a835093fSnbeams     CeedEvalMode eval_mode = eval_mode_in[i];
990a835093fSnbeams     if (eval_mode == CEED_EVAL_INTERP) {
991b7453713SJeremy L Thompson       CeedCallHip(ceed, hipMemcpy(&asmb->d_B_in[mat_start], interp_in, elem_size * num_qpts * sizeof(CeedScalar), hipMemcpyHostToDevice));
992b7453713SJeremy L Thompson       mat_start += elem_size * num_qpts;
993a835093fSnbeams     } else if (eval_mode == CEED_EVAL_GRAD) {
994b7453713SJeremy L Thompson       CeedCallHip(ceed, hipMemcpy(&asmb->d_B_in[mat_start], grad_in, dim * elem_size * num_qpts * sizeof(CeedScalar), hipMemcpyHostToDevice));
995b7453713SJeremy L Thompson       mat_start += dim * elem_size * num_qpts;
996a835093fSnbeams     }
997a835093fSnbeams   }
998a835093fSnbeams 
999a835093fSnbeams   const CeedScalar *interp_out, *grad_out;
1000b7453713SJeremy L Thompson 
1001b7453713SJeremy L Thompson   // Note that this function currently assumes 1 basis, so this should always be true for now
1002a835093fSnbeams   if (basis_out == basis_in) {
1003a835093fSnbeams     interp_out = interp_in;
1004a835093fSnbeams     grad_out   = grad_in;
1005a835093fSnbeams   } else {
10062b730f8bSJeremy L Thompson     CeedCallBackend(CeedBasisGetInterp(basis_out, &interp_out));
10072b730f8bSJeremy L Thompson     CeedCallBackend(CeedBasisGetGrad(basis_out, &grad_out));
1008a835093fSnbeams   }
1009a835093fSnbeams 
1010a835093fSnbeams   // Load into B_out, in order that they will be used in eval_mode
1011b7453713SJeremy L Thompson   const CeedInt out_bytes = size_B_out * sizeof(CeedScalar);
1012b7453713SJeremy L Thompson 
1013a835093fSnbeams   mat_start = 0;
1014b7453713SJeremy L Thompson   CeedCallHip(ceed, hipMalloc((void **)&asmb->d_B_out, out_bytes));
1015a835093fSnbeams   for (int i = 0; i < num_B_out_mats_to_load; i++) {
1016a835093fSnbeams     CeedEvalMode eval_mode = eval_mode_out[i];
1017a835093fSnbeams     if (eval_mode == CEED_EVAL_INTERP) {
1018b7453713SJeremy L Thompson       CeedCallHip(ceed, hipMemcpy(&asmb->d_B_out[mat_start], interp_out, elem_size * num_qpts * sizeof(CeedScalar), hipMemcpyHostToDevice));
1019b7453713SJeremy L Thompson       mat_start += elem_size * num_qpts;
1020a835093fSnbeams     } else if (eval_mode == CEED_EVAL_GRAD) {
1021b7453713SJeremy L Thompson       CeedCallHip(ceed, hipMemcpy(&asmb->d_B_out[mat_start], grad_out, dim * elem_size * num_qpts * sizeof(CeedScalar), hipMemcpyHostToDevice));
1022b7453713SJeremy L Thompson       mat_start += dim * elem_size * num_qpts;
1023a835093fSnbeams     }
1024a835093fSnbeams   }
1025a835093fSnbeams   return CEED_ERROR_SUCCESS;
1026a835093fSnbeams }
1027a835093fSnbeams 
1028a835093fSnbeams //------------------------------------------------------------------------------
1029cefa2673SJeremy L Thompson // Assemble matrix data for COO matrix of assembled operator.
1030cefa2673SJeremy L Thompson // The sparsity pattern is set by CeedOperatorLinearAssembleSymbolic.
1031cefa2673SJeremy L Thompson //
1032ea61e9acSJeremy L Thompson // Note that this (and other assembly routines) currently assume only one active input restriction/basis per operator (could have multiple basis eval
1033cefa2673SJeremy L Thompson // modes).
1034cefa2673SJeremy L Thompson // TODO: allow multiple active input restrictions/basis objects
1035a835093fSnbeams //------------------------------------------------------------------------------
10362b730f8bSJeremy L Thompson static int CeedSingleOperatorAssemble_Hip(CeedOperator op, CeedInt offset, CeedVector values) {
1037a835093fSnbeams   Ceed                ceed;
1038b7453713SJeremy L Thompson   CeedSize            values_length = 0, assembled_qf_length = 0;
1039b7453713SJeremy L Thompson   CeedInt             use_ceedsize_idx = 0;
1040b7453713SJeremy L Thompson   CeedScalar         *values_array;
1041b7453713SJeremy L Thompson   const CeedScalar   *qf_array;
1042b7453713SJeremy L Thompson   CeedVector          assembled_qf = NULL;
1043b7453713SJeremy L Thompson   CeedElemRestriction rstr_q       = NULL;
1044a835093fSnbeams   CeedOperator_Hip   *impl;
1045b7453713SJeremy L Thompson 
1046b7453713SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
10472b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &impl));
1048a835093fSnbeams 
1049a835093fSnbeams   // Assemble QFunction
10502b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorLinearAssembleQFunctionBuildOrUpdate(op, &assembled_qf, &rstr_q, CEED_REQUEST_IMMEDIATE));
10512b730f8bSJeremy L Thompson   CeedCallBackend(CeedElemRestrictionDestroy(&rstr_q));
105228ec399dSJeremy L Thompson   CeedCallBackend(CeedVectorGetArray(values, CEED_MEM_DEVICE, &values_array));
1053a835093fSnbeams   values_array += offset;
10542b730f8bSJeremy L Thompson   CeedCallBackend(CeedVectorGetArrayRead(assembled_qf, CEED_MEM_DEVICE, &qf_array));
1055a835093fSnbeams 
10569330daecSnbeams   CeedCallBackend(CeedVectorGetLength(values, &values_length));
10579330daecSnbeams   CeedCallBackend(CeedVectorGetLength(assembled_qf, &assembled_qf_length));
10589330daecSnbeams   if ((values_length > INT_MAX) || (assembled_qf_length > INT_MAX)) use_ceedsize_idx = 1;
10599330daecSnbeams   // Setup
10609330daecSnbeams   if (!impl->asmb) {
10619330daecSnbeams     CeedCallBackend(CeedSingleOperatorAssembleSetup_Hip(op, use_ceedsize_idx));
10629330daecSnbeams     assert(impl->asmb != NULL);
10639330daecSnbeams   }
10649330daecSnbeams 
1065a835093fSnbeams   // Compute B^T D B
1066b7453713SJeremy L Thompson   const CeedInt num_elem       = impl->asmb->num_elem;
1067b7453713SJeremy L Thompson   const CeedInt elem_per_block = impl->asmb->elem_per_block;
1068b7453713SJeremy L Thompson   const CeedInt grid           = num_elem / elem_per_block + ((num_elem / elem_per_block * elem_per_block < num_elem) ? 1 : 0);
10692b730f8bSJeremy L Thompson   void         *args[]         = {&impl->asmb->d_B_in, &impl->asmb->d_B_out, &qf_array, &values_array};
1070b7453713SJeremy L Thompson 
10712b730f8bSJeremy L Thompson   CeedCallBackend(
1072b7453713SJeremy L Thompson       CeedRunKernelDim_Hip(ceed, impl->asmb->linearAssemble, grid, impl->asmb->block_size_x, impl->asmb->block_size_y, elem_per_block, args));
1073a835093fSnbeams 
1074a835093fSnbeams   // Restore arrays
10752b730f8bSJeremy L Thompson   CeedCallBackend(CeedVectorRestoreArray(values, &values_array));
10762b730f8bSJeremy L Thompson   CeedCallBackend(CeedVectorRestoreArrayRead(assembled_qf, &qf_array));
1077a835093fSnbeams 
1078a835093fSnbeams   // Cleanup
10792b730f8bSJeremy L Thompson   CeedCallBackend(CeedVectorDestroy(&assembled_qf));
1080a835093fSnbeams   return CEED_ERROR_SUCCESS;
1081a835093fSnbeams }
1082a835093fSnbeams 
1083a835093fSnbeams //------------------------------------------------------------------------------
10840d0321e0SJeremy L Thompson // Create operator
10850d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
10860d0321e0SJeremy L Thompson int CeedOperatorCreate_Hip(CeedOperator op) {
10870d0321e0SJeremy L Thompson   Ceed              ceed;
10880d0321e0SJeremy L Thompson   CeedOperator_Hip *impl;
10890d0321e0SJeremy L Thompson 
1090b7453713SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
10912b730f8bSJeremy L Thompson   CeedCallBackend(CeedCalloc(1, &impl));
10922b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorSetData(op, impl));
10932b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleQFunction", CeedOperatorLinearAssembleQFunction_Hip));
10942b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleQFunctionUpdate", CeedOperatorLinearAssembleQFunctionUpdate_Hip));
10952b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleAddDiagonal", CeedOperatorLinearAssembleAddDiagonal_Hip));
10962b730f8bSJeremy L Thompson   CeedCallBackend(
10972b730f8bSJeremy L Thompson       CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleAddPointBlockDiagonal", CeedOperatorLinearAssembleAddPointBlockDiagonal_Hip));
10982b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleSingle", CeedSingleOperatorAssemble_Hip));
10992b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "ApplyAdd", CeedOperatorApplyAdd_Hip));
11002b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "Destroy", CeedOperatorDestroy_Hip));
11010d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
11020d0321e0SJeremy L Thompson }
11030d0321e0SJeremy L Thompson 
11040d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
1105