xref: /libCEED/rust/libceed-sys/c-src/backends/sycl-ref/ceed-sycl-ref-operator.sycl.cpp (revision 681d0ea73ee05192cf73f31e6e4a886b41175395)
1 // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other
2 // CEED contributors. All Rights Reserved. See the top-level LICENSE and NOTICE
3 // files for details.
4 //
5 // SPDX-License-Identifier: BSD-2-Clause
6 //
7 // This file is part of CEED:  http://github.com/ceed
8 
9 #include <ceed/backend.h>
10 #include <ceed/ceed.h>
11 
12 #include <cassert>
13 #include <string>
14 #include <sycl/sycl.hpp>
15 
16 #include "../sycl/ceed-sycl-compile.hpp"
17 #include "ceed-sycl-ref.hpp"
18 
19 class CeedOperatorSyclLinearDiagonal;
20 class CeedOperatorSyclLinearAssemble;
21 class CeedOperatorSyclLinearAssembleFallback;
22 
23 //------------------------------------------------------------------------------
24 //  Get Basis Emode Pointer
25 //------------------------------------------------------------------------------
26 void CeedOperatorGetBasisPointer_Sycl(const CeedScalar **basis_ptr, CeedEvalMode eval_mode, const CeedScalar *identity, const CeedScalar *interp,
27                                       const CeedScalar *grad) {
28   switch (eval_mode) {
29     case CEED_EVAL_NONE:
30       *basis_ptr = identity;
31       break;
32     case CEED_EVAL_INTERP:
33       *basis_ptr = interp;
34       break;
35     case CEED_EVAL_GRAD:
36       *basis_ptr = grad;
37       break;
38     case CEED_EVAL_WEIGHT:
39     case CEED_EVAL_DIV:
40     case CEED_EVAL_CURL:
41       break;  // Caught by QF Assembly
42   }
43 }
44 
45 //------------------------------------------------------------------------------
46 // Destroy operator
47 //------------------------------------------------------------------------------
48 static int CeedOperatorDestroy_Sycl(CeedOperator op) {
49   Ceed               ceed;
50   Ceed_Sycl         *sycl_data;
51   CeedOperator_Sycl *impl;
52 
53   CeedCallBackend(CeedOperatorGetData(op, &impl));
54   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
55   CeedCallBackend(CeedGetData(ceed, &sycl_data));
56 
57   // Apply data
58   for (CeedInt i = 0; i < impl->num_e_in + impl->num_e_out; i++) {
59     CeedCallBackend(CeedVectorDestroy(&impl->e_vecs[i]));
60   }
61   CeedCallBackend(CeedFree(&impl->e_vecs));
62 
63   for (CeedInt i = 0; i < impl->num_e_in; i++) {
64     CeedCallBackend(CeedVectorDestroy(&impl->q_vecs_in[i]));
65   }
66   CeedCallBackend(CeedFree(&impl->q_vecs_in));
67 
68   for (CeedInt i = 0; i < impl->num_e_out; i++) {
69     CeedCallBackend(CeedVectorDestroy(&impl->q_vecs_out[i]));
70   }
71   CeedCallBackend(CeedFree(&impl->q_vecs_out));
72 
73   // QFunction assembly dataf
74   for (CeedInt i = 0; i < impl->num_active_in; i++) {
75     CeedCallBackend(CeedVectorDestroy(&impl->qf_active_in[i]));
76   }
77   CeedCallBackend(CeedFree(&impl->qf_active_in));
78 
79   // Diag data
80   if (impl->diag) {
81     CeedCallBackend(CeedFree(&impl->diag->h_eval_mode_in));
82     CeedCallBackend(CeedFree(&impl->diag->h_eval_mode_out));
83 
84     CeedCallSycl(ceed, sycl_data->sycl_queue.wait_and_throw());
85     CeedCallSycl(ceed, sycl::free(impl->diag->d_eval_mode_in, sycl_data->sycl_context));
86     CeedCallSycl(ceed, sycl::free(impl->diag->d_eval_mode_out, sycl_data->sycl_context));
87     CeedCallSycl(ceed, sycl::free(impl->diag->d_identity, sycl_data->sycl_context));
88     CeedCallSycl(ceed, sycl::free(impl->diag->d_interp_in, sycl_data->sycl_context));
89     CeedCallSycl(ceed, sycl::free(impl->diag->d_interp_out, sycl_data->sycl_context));
90     CeedCallSycl(ceed, sycl::free(impl->diag->d_grad_in, sycl_data->sycl_context));
91     CeedCallSycl(ceed, sycl::free(impl->diag->d_grad_out, sycl_data->sycl_context));
92     CeedCallBackend(CeedElemRestrictionDestroy(&impl->diag->point_block_diag_rstr));
93 
94     CeedCallBackend(CeedVectorDestroy(&impl->diag->elem_diag));
95     CeedCallBackend(CeedVectorDestroy(&impl->diag->point_block_elem_diag));
96   }
97   CeedCallBackend(CeedFree(&impl->diag));
98 
99   if (impl->asmb) {
100     CeedCallSycl(ceed, sycl_data->sycl_queue.wait_and_throw());
101     CeedCallSycl(ceed, sycl::free(impl->asmb->d_B_in, sycl_data->sycl_context));
102     CeedCallSycl(ceed, sycl::free(impl->asmb->d_B_out, sycl_data->sycl_context));
103   }
104   CeedCallBackend(CeedFree(&impl->asmb));
105 
106   CeedCallBackend(CeedFree(&impl));
107   return CEED_ERROR_SUCCESS;
108 }
109 
110 //------------------------------------------------------------------------------
111 // Setup infields or outfields
112 //------------------------------------------------------------------------------
113 static int CeedOperatorSetupFields_Sycl(CeedQFunction qf, CeedOperator op, bool is_input, CeedVector *e_vecs, CeedVector *q_vecs, CeedInt start_e,
114                                         CeedInt num_fields, CeedInt Q, CeedInt num_elem) {
115   Ceed                ceed;
116   CeedSize            q_size;
117   bool                is_strided, skip_restriction;
118   CeedInt             dim, size;
119   CeedOperatorField  *op_fields;
120   CeedQFunctionField *qf_fields;
121 
122   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
123   if (is_input) {
124     CeedCallBackend(CeedOperatorGetFields(op, NULL, &op_fields, NULL, NULL));
125     CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_fields, NULL, NULL));
126   } else {
127     CeedCallBackend(CeedOperatorGetFields(op, NULL, NULL, NULL, &op_fields));
128     CeedCallBackend(CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qf_fields));
129   }
130 
131   // Loop over fields
132   for (CeedInt i = 0; i < num_fields; i++) {
133     CeedEvalMode        eval_mode;
134     CeedVector          vec;
135     CeedElemRestriction rstr;
136     CeedBasis           basis;
137 
138     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode));
139 
140     is_strided       = false;
141     skip_restriction = false;
142     if (eval_mode != CEED_EVAL_WEIGHT) {
143       CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_fields[i], &rstr));
144 
145       // Check whether this field can skip the element restriction:
146       // must be passive input, with  eval_mode NONE, and have a strided restriction with CEED_STRIDES_BACKEND.
147 
148       // First, check whether the field is input or output:
149       if (is_input) {
150         // Check for passive input:
151         CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec));
152         if (vec != CEED_VECTOR_ACTIVE) {
153           // Check  eval_mode
154           if (eval_mode == CEED_EVAL_NONE) {
155             // Check for  is_strided restriction
156             CeedCallBackend(CeedElemRestrictionIsStrided(rstr, &is_strided));
157             if (is_strided) {
158               // Check if vector is already in preferred backend ordering
159               CeedCallBackend(CeedElemRestrictionHasBackendStrides(rstr, &skip_restriction));
160             }
161           }
162         }
163         CeedCallBackend(CeedVectorDestroy(&vec));
164       }
165       if (skip_restriction) {
166         // We do not need an E-Vector, but will use the input field vector's data directly in the operator application
167         e_vecs[i + start_e] = NULL;
168       } else {
169         CeedCallBackend(CeedElemRestrictionCreateVector(rstr, NULL, &e_vecs[i + start_e]));
170       }
171       CeedCallBackend(CeedElemRestrictionDestroy(&rstr));
172     }
173 
174     switch (eval_mode) {
175       case CEED_EVAL_NONE:
176         CeedCallBackend(CeedQFunctionFieldGetSize(qf_fields[i], &size));
177         q_size = (CeedSize)num_elem * Q * size;
178         CeedCallBackend(CeedVectorCreate(ceed, q_size, &q_vecs[i]));
179         break;
180       case CEED_EVAL_INTERP:
181         CeedCallBackend(CeedQFunctionFieldGetSize(qf_fields[i], &size));
182         q_size = (CeedSize)num_elem * Q * size;
183         CeedCallBackend(CeedVectorCreate(ceed, q_size, &q_vecs[i]));
184         break;
185       case CEED_EVAL_GRAD:
186         CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis));
187         CeedCallBackend(CeedQFunctionFieldGetSize(qf_fields[i], &size));
188         CeedCallBackend(CeedBasisGetDimension(basis, &dim));
189         CeedCallBackend(CeedBasisDestroy(&basis));
190         q_size = (CeedSize)num_elem * Q * size;
191         CeedCallBackend(CeedVectorCreate(ceed, q_size, &q_vecs[i]));
192         break;
193       case CEED_EVAL_WEIGHT:  // Only on input fields
194         CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis));
195         q_size = (CeedSize)num_elem * Q;
196         CeedCallBackend(CeedVectorCreate(ceed, q_size, &q_vecs[i]));
197         CeedCallBackend(CeedBasisApply(basis, num_elem, CEED_NOTRANSPOSE, CEED_EVAL_WEIGHT, CEED_VECTOR_NONE, q_vecs[i]));
198         CeedCallBackend(CeedBasisDestroy(&basis));
199         break;
200       case CEED_EVAL_DIV:
201         break;  // TODO: Not implemented
202       case CEED_EVAL_CURL:
203         break;  // TODO: Not implemented
204     }
205   }
206   return CEED_ERROR_SUCCESS;
207 }
208 
209 //------------------------------------------------------------------------------
210 // CeedOperator needs to connect all the named fields (be they active or
211 // passive) to the named inputs and outputs of its CeedQFunction.
212 //------------------------------------------------------------------------------
213 static int CeedOperatorSetup_Sycl(CeedOperator op) {
214   Ceed                ceed;
215   bool                is_setup_done;
216   CeedInt             Q, num_elem, num_input_fields, num_output_fields;
217   CeedQFunctionField *qf_input_fields, *qf_output_fields;
218   CeedQFunction       qf;
219   CeedOperatorField  *op_input_fields, *op_output_fields;
220   CeedOperator_Sycl  *impl;
221 
222   CeedCallBackend(CeedOperatorIsSetupDone(op, &is_setup_done));
223   if (is_setup_done) return CEED_ERROR_SUCCESS;
224 
225   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
226   CeedCallBackend(CeedOperatorGetData(op, &impl));
227   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
228   CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q));
229   CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem));
230   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
231   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
232 
233   // Allocate
234   CeedCallBackend(CeedCalloc(num_input_fields + num_output_fields, &impl->e_vecs));
235 
236   CeedCallBackend(CeedCalloc(CEED_FIELD_MAX, &impl->q_vecs_in));
237   CeedCallBackend(CeedCalloc(CEED_FIELD_MAX, &impl->q_vecs_out));
238 
239   impl->num_e_in  = num_input_fields;
240   impl->num_e_out = num_output_fields;
241 
242   // Set up infield and outfield  e_vecs and  q_vecs
243   // Infields
244   CeedCallBackend(CeedOperatorSetupFields_Sycl(qf, op, true, impl->e_vecs, impl->q_vecs_in, 0, num_input_fields, Q, num_elem));
245   // Outfields
246   CeedCallBackend(CeedOperatorSetupFields_Sycl(qf, op, false, impl->e_vecs, impl->q_vecs_out, num_input_fields, num_output_fields, Q, num_elem));
247 
248   CeedCallBackend(CeedOperatorSetSetupDone(op));
249   return CEED_ERROR_SUCCESS;
250 }
251 
252 //------------------------------------------------------------------------------
253 // Setup Operator Inputs
254 //------------------------------------------------------------------------------
255 static inline int CeedOperatorSetupInputs_Sycl(CeedInt num_input_fields, CeedQFunctionField *qf_input_fields, CeedOperatorField *op_input_fields,
256                                                CeedVector in_vec, const bool skip_active, CeedScalar *e_data[2 * CEED_FIELD_MAX],
257                                                CeedOperator_Sycl *impl, CeedRequest *request) {
258   for (CeedInt i = 0; i < num_input_fields; i++) {
259     bool         is_active;
260     CeedEvalMode eval_mode;
261     CeedVector   vec;
262 
263     // Get input vector
264     CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec));
265     is_active = vec == CEED_VECTOR_ACTIVE;
266     if (is_active) {
267       if (skip_active) continue;
268       else vec = in_vec;
269     }
270 
271     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
272     if (eval_mode == CEED_EVAL_WEIGHT) {  // Skip
273     } else {
274       // Restrict, if necessary
275       if (!impl->e_vecs[i]) {
276         // No restriction for this field; read data directly from vec.
277         CeedCallBackend(CeedVectorGetArrayRead(vec, CEED_MEM_DEVICE, (const CeedScalar **)&e_data[i]));
278       } else {
279         CeedElemRestriction rstr;
280 
281         CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &rstr));
282         CeedCallBackend(CeedElemRestrictionApply(rstr, CEED_NOTRANSPOSE, vec, impl->e_vecs[i], request));
283         CeedCallBackend(CeedElemRestrictionDestroy(&rstr));
284         CeedCallBackend(CeedVectorGetArrayRead(impl->e_vecs[i], CEED_MEM_DEVICE, (const CeedScalar **)&e_data[i]));
285       }
286     }
287     if (!is_active) CeedCallBackend(CeedVectorDestroy(&vec));
288   }
289   return CEED_ERROR_SUCCESS;
290 }
291 
292 //------------------------------------------------------------------------------
293 // Input Basis Action
294 //------------------------------------------------------------------------------
295 static inline int CeedOperatorInputBasis_Sycl(CeedInt num_elem, CeedQFunctionField *qf_input_fields, CeedOperatorField *op_input_fields,
296                                               CeedInt num_input_fields, const bool skip_active, CeedScalar *e_data[2 * CEED_FIELD_MAX],
297                                               CeedOperator_Sycl *impl) {
298   for (CeedInt i = 0; i < num_input_fields; i++) {
299     CeedEvalMode eval_mode;
300     CeedBasis    basis;
301 
302     // Skip active input
303     if (skip_active) {
304       bool       is_active;
305       CeedVector vec;
306 
307       CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec));
308       is_active = vec == CEED_VECTOR_ACTIVE;
309       CeedCallBackend(CeedVectorDestroy(&vec));
310       if (is_active) continue;
311     }
312     // Basis action
313     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
314     switch (eval_mode) {
315       case CEED_EVAL_NONE:
316         CeedCallBackend(CeedVectorSetArray(impl->q_vecs_in[i], CEED_MEM_DEVICE, CEED_USE_POINTER, e_data[i]));
317         break;
318       case CEED_EVAL_INTERP:
319       case CEED_EVAL_GRAD:
320         CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
321         CeedCallBackend(CeedBasisApply(basis, num_elem, CEED_NOTRANSPOSE, eval_mode, impl->e_vecs[i], impl->q_vecs_in[i]));
322         CeedCallBackend(CeedBasisDestroy(&basis));
323         break;
324       case CEED_EVAL_WEIGHT:
325         break;  // No action
326       case CEED_EVAL_DIV:
327         break;  // TODO: Not implemented
328       case CEED_EVAL_CURL:
329         break;  // TODO: Not implemented
330     }
331   }
332   return CEED_ERROR_SUCCESS;
333 }
334 
335 //------------------------------------------------------------------------------
336 // Restore Input Vectors
337 //------------------------------------------------------------------------------
338 static inline int CeedOperatorRestoreInputs_Sycl(CeedInt num_input_fields, CeedQFunctionField *qf_input_fields, CeedOperatorField *op_input_fields,
339                                                  const bool skip_active, CeedScalar *e_data[2 * CEED_FIELD_MAX], CeedOperator_Sycl *impl) {
340   for (CeedInt i = 0; i < num_input_fields; i++) {
341     bool         is_active;
342     CeedEvalMode eval_mode;
343     CeedVector   vec;
344 
345     CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec));
346     is_active = vec == CEED_VECTOR_ACTIVE;
347     // Skip active input
348     if (skip_active) {
349       if (is_active) continue;
350     }
351     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
352     if (eval_mode == CEED_EVAL_WEIGHT) {  // Skip
353     } else {
354       if (!impl->e_vecs[i]) {  // This was a  skip_restriction case
355         CeedCallBackend(CeedVectorRestoreArrayRead(vec, (const CeedScalar **)&e_data[i]));
356       } else {
357         CeedCallBackend(CeedVectorRestoreArrayRead(impl->e_vecs[i], (const CeedScalar **)&e_data[i]));
358       }
359     }
360     if (!is_active) CeedCallBackend(CeedVectorDestroy(&vec));
361   }
362   return CEED_ERROR_SUCCESS;
363 }
364 
365 //------------------------------------------------------------------------------
366 // Apply and add to output
367 //------------------------------------------------------------------------------
368 static int CeedOperatorApplyAdd_Sycl(CeedOperator op, CeedVector in_vec, CeedVector out_vec, CeedRequest *request) {
369   CeedInt             Q, num_elem, elem_size, num_input_fields, num_output_fields, size;
370   CeedEvalMode        eval_mode;
371   CeedScalar         *e_data[2 * CEED_FIELD_MAX] = {0};
372   CeedQFunctionField *qf_input_fields, *qf_output_fields;
373   CeedQFunction       qf;
374   CeedOperatorField  *op_input_fields, *op_output_fields;
375   CeedOperator_Sycl  *impl;
376 
377   CeedCallBackend(CeedOperatorGetData(op, &impl));
378   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
379   CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q));
380   CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem));
381   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
382   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
383 
384   // Setup
385   CeedCallBackend(CeedOperatorSetup_Sycl(op));
386 
387   // Input Evecs and Restriction
388   CeedCallBackend(CeedOperatorSetupInputs_Sycl(num_input_fields, qf_input_fields, op_input_fields, in_vec, false, e_data, impl, request));
389 
390   // Input basis apply if needed
391   CeedCallBackend(CeedOperatorInputBasis_Sycl(num_elem, qf_input_fields, op_input_fields, num_input_fields, false, e_data, impl));
392 
393   // Output pointers, as necessary
394   for (CeedInt i = 0; i < num_output_fields; i++) {
395     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
396     if (eval_mode == CEED_EVAL_NONE) {
397       // Set the output Q-Vector to use the E-Vector data directly
398       CeedCallBackend(CeedVectorGetArrayWrite(impl->e_vecs[i + impl->num_e_in], CEED_MEM_DEVICE, &e_data[i + num_input_fields]));
399       CeedCallBackend(CeedVectorSetArray(impl->q_vecs_out[i], CEED_MEM_DEVICE, CEED_USE_POINTER, e_data[i + num_input_fields]));
400     }
401   }
402 
403   // Q function
404   CeedCallBackend(CeedQFunctionApply(qf, num_elem * Q, impl->q_vecs_in, impl->q_vecs_out));
405 
406   // Output basis apply if needed
407   for (CeedInt i = 0; i < num_output_fields; i++) {
408     CeedElemRestriction rstr;
409     CeedBasis           basis;
410 
411     // Get elem_size, eval_mode, size
412     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &rstr));
413     CeedCallBackend(CeedElemRestrictionGetElementSize(rstr, &elem_size));
414     CeedCallBackend(CeedElemRestrictionDestroy(&rstr));
415     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
416     CeedCallBackend(CeedQFunctionFieldGetSize(qf_output_fields[i], &size));
417     // Basis action
418     switch (eval_mode) {
419       case CEED_EVAL_NONE:
420         break;
421       case CEED_EVAL_INTERP:
422       case CEED_EVAL_GRAD:
423         CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
424         CeedCallBackend(CeedBasisApply(basis, num_elem, CEED_TRANSPOSE, eval_mode, impl->q_vecs_out[i], impl->e_vecs[i + impl->num_e_in]));
425         CeedCallBackend(CeedBasisDestroy(&basis));
426         break;
427       // LCOV_EXCL_START
428       case CEED_EVAL_WEIGHT:
429         Ceed ceed;
430 
431         CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
432         return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_WEIGHT cannot be an output evaluation mode");
433         break;  // Should not occur
434       case CEED_EVAL_DIV:
435       case CEED_EVAL_CURL: {
436         Ceed ceed;
437 
438         CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
439         return CeedError(ceed, CEED_ERROR_BACKEND, "%s not supported", CeedEvalModes[eval_mode]);
440         break;  // Should not occur
441       }
442         // LCOV_EXCL_STOP
443     }
444   }
445 
446   // Output restriction
447   for (CeedInt i = 0; i < num_output_fields; i++) {
448     bool                is_active;
449     CeedEvalMode        eval_mode;
450     CeedVector          vec;
451     CeedElemRestriction rstr;
452 
453     // Restore evec
454     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
455     if (eval_mode == CEED_EVAL_NONE) {
456       CeedCallBackend(CeedVectorRestoreArray(impl->e_vecs[i + impl->num_e_in], &e_data[i + num_input_fields]));
457     }
458     // Restrict
459     CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[i], &vec));
460     is_active = vec == CEED_VECTOR_ACTIVE;
461     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &rstr));
462     if (is_active) vec = out_vec;
463     CeedCallBackend(CeedElemRestrictionApply(rstr, CEED_TRANSPOSE, impl->e_vecs[i + impl->num_inputs], vec, request));
464     if (!is_active) CeedCallBackend(CeedVectorDestroy(&vec));
465     CeedCallBackend(CeedElemRestrictionDestroy(&rstr));
466   }
467 
468   // Restore input arrays
469   CeedCallBackend(CeedOperatorRestoreInputs_Sycl(num_input_fields, qf_input_fields, op_input_fields, false, e_data, impl));
470   return CEED_ERROR_SUCCESS;
471 }
472 
473 //------------------------------------------------------------------------------
474 // Core code for assembling linear QFunction
475 //------------------------------------------------------------------------------
476 static inline int CeedOperatorLinearAssembleQFunctionCore_Sycl(CeedOperator op, bool build_objects, CeedVector *assembled, CeedElemRestriction *rstr,
477                                                                CeedRequest *request) {
478   Ceed                ceed, ceed_parent;
479   CeedSize            q_size;
480   CeedInt             num_active_in, num_active_out, Q, num_elem, num_input_fields, num_output_fields, size;
481   CeedScalar         *assembled_array, *e_data[2 * CEED_FIELD_MAX] = {NULL};
482   CeedVector         *active_in;
483   CeedQFunctionField *qf_input_fields, *qf_output_fields;
484   CeedQFunction       qf;
485   CeedOperatorField  *op_input_fields, *op_output_fields;
486   CeedOperator_Sycl  *impl;
487 
488   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
489   CeedCallBackend(CeedOperatorGetFallbackParentCeed(op, &ceed_parent));
490   CeedCallBackend(CeedOperatorGetData(op, &impl));
491   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
492   CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q));
493   CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem));
494   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
495   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
496   active_in     = impl->qf_active_in;
497   num_active_in = impl->num_active_in, num_active_out = impl->num_active_out;
498 
499   // Setup
500   CeedCallBackend(CeedOperatorSetup_Sycl(op));
501 
502   // Input Evecs and Restriction
503   CeedCallBackend(CeedOperatorSetupInputs_Sycl(num_input_fields, qf_input_fields, op_input_fields, NULL, true, e_data, impl, request));
504 
505   // Count number of active input fields
506   if (!num_active_in) {
507     for (CeedInt i = 0; i < num_input_fields; i++) {
508       CeedScalar *q_vec_array;
509       CeedVector  vec;
510 
511       // Check if active input
512       CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec));
513       if (vec == CEED_VECTOR_ACTIVE) {
514         CeedCallBackend(CeedQFunctionFieldGetSize(qf_input_fields[i], &size));
515         CeedCallBackend(CeedVectorSetValue(impl->q_vecs_in[i], 0.0));
516         CeedCallBackend(CeedVectorGetArray(impl->q_vecs_in[i], CEED_MEM_DEVICE, &q_vec_array));
517         CeedCallBackend(CeedRealloc(num_active_in + size, &active_in));
518         for (CeedInt field = 0; field < size; field++) {
519           q_size = (CeedSize)Q * num_elem;
520           CeedCallBackend(CeedVectorCreate(ceed, q_size, &active_in[num_active_in + field]));
521           CeedCallBackend(
522               CeedVectorSetArray(active_in[num_active_in + field], CEED_MEM_DEVICE, CEED_USE_POINTER, &q_vec_array[field * Q * num_elem]));
523         }
524         num_active_in += size;
525         CeedCallBackend(CeedVectorRestoreArray(impl->q_vecs_in[i], &q_vec_array));
526       }
527       CeedCallBackend(CeedVectorDestroy(&vec));
528     }
529     impl->num_active_in = num_active_in;
530     impl->qf_active_in  = active_in;
531   }
532 
533   // Count number of active output fields
534   if (!num_active_out) {
535     for (CeedInt i = 0; i < num_output_fields; i++) {
536       CeedVector vec;
537 
538       // Check if active output
539       CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[i], &vec));
540       if (vec == CEED_VECTOR_ACTIVE) {
541         CeedCallBackend(CeedQFunctionFieldGetSize(qf_output_fields[i], &size));
542         num_active_out += size;
543       }
544       CeedCallBackend(CeedVectorDestroy(&vec));
545     }
546     impl->num_active_out = num_active_out;
547   }
548 
549   // Check sizes
550   CeedCheck(num_active_in > 0 && num_active_out > 0, ceed, CEED_ERROR_BACKEND, "Cannot assemble QFunction without active inputs and outputs");
551 
552   // Build objects if needed
553   if (build_objects) {
554     CeedSize l_size     = (CeedSize)num_elem * Q * num_active_in * num_active_out;
555     CeedInt  strides[3] = {1, num_elem * Q, Q}; /* *NOPAD* */
556 
557     // Create output restriction
558     CeedCallBackend(CeedElemRestrictionCreateStrided(ceed_parent, num_elem, Q, num_active_in * num_active_out, l_size, strides, rstr));
559     // Create assembled vector
560     CeedCallBackend(CeedVectorCreate(ceed_parent, l_size, assembled));
561   }
562   CeedCallBackend(CeedVectorSetValue(*assembled, 0.0));
563   CeedCallBackend(CeedVectorGetArray(*assembled, CEED_MEM_DEVICE, &assembled_array));
564 
565   // Input basis apply
566   CeedCallBackend(CeedOperatorInputBasis_Sycl(num_elem, qf_input_fields, op_input_fields, num_input_fields, true, e_data, impl));
567 
568   // Assemble QFunction
569   for (CeedInt in = 0; in < num_active_in; in++) {
570     // Set Inputs
571     CeedCallBackend(CeedVectorSetValue(active_in[in], 1.0));
572     if (num_active_in > 1) {
573       CeedCallBackend(CeedVectorSetValue(active_in[(in + num_active_in - 1) % num_active_in], 0.0));
574     }
575     // Set Outputs
576     for (CeedInt out = 0; out < num_output_fields; out++) {
577       CeedVector vec;
578 
579       // Check if active output
580       CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[out], &vec));
581       if (vec == CEED_VECTOR_ACTIVE) {
582         CeedCallBackend(CeedVectorSetArray(impl->q_vecs_out[out], CEED_MEM_DEVICE, CEED_USE_POINTER, assembled_array));
583         CeedCallBackend(CeedQFunctionFieldGetSize(qf_output_fields[out], &size));
584         assembled_array += size * Q * num_elem;  // Advance the pointer by the size of the output
585       }
586       CeedCallBackend(CeedVectorDestroy(&vec));
587     }
588     // Apply QFunction
589     CeedCallBackend(CeedQFunctionApply(qf, Q * num_elem, impl->q_vecs_in, impl->q_vecs_out));
590   }
591 
592   // Un-set output Qvecs to prevent accidental overwrite of Assembled
593   for (CeedInt out = 0; out < num_output_fields; out++) {
594     CeedVector vec;
595 
596     // Check if active output
597     CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[out], &vec));
598     if (vec == CEED_VECTOR_ACTIVE) {
599       CeedCallBackend(CeedVectorTakeArray(impl->q_vecs_out[out], CEED_MEM_DEVICE, NULL));
600     }
601     CeedCallBackend(CeedVectorDestroy(&vec));
602   }
603 
604   // Restore input arrays
605   CeedCallBackend(CeedOperatorRestoreInputs_Sycl(num_input_fields, qf_input_fields, op_input_fields, true, e_data, impl));
606 
607   // Restore output
608   CeedCallBackend(CeedVectorRestoreArray(*assembled, &assembled_array));
609   return CEED_ERROR_SUCCESS;
610 }
611 
612 //------------------------------------------------------------------------------
613 // Assemble Linear QFunction
614 //------------------------------------------------------------------------------
615 static int CeedOperatorLinearAssembleQFunction_Sycl(CeedOperator op, CeedVector *assembled, CeedElemRestriction *rstr, CeedRequest *request) {
616   return CeedOperatorLinearAssembleQFunctionCore_Sycl(op, true, assembled, rstr, request);
617 }
618 
619 //------------------------------------------------------------------------------
620 // Update Assembled Linear QFunction
621 //------------------------------------------------------------------------------
622 static int CeedOperatorLinearAssembleQFunctionUpdate_Sycl(CeedOperator op, CeedVector assembled, CeedElemRestriction rstr, CeedRequest *request) {
623   return CeedOperatorLinearAssembleQFunctionCore_Sycl(op, false, &assembled, &rstr, request);
624 }
625 
626 //------------------------------------------------------------------------------
627 // Assemble diagonal setup
628 //------------------------------------------------------------------------------
629 static inline int CeedOperatorAssembleDiagonalSetup_Sycl(CeedOperator op) {
630   Ceed                ceed;
631   Ceed_Sycl          *sycl_data;
632   CeedInt             num_input_fields, num_output_fields, num_eval_mode_in = 0, num_comp = 0, dim = 1, num_eval_mode_out = 0;
633   CeedEvalMode       *eval_mode_in = NULL, *eval_mode_out = NULL;
634   CeedBasis           basis_in = NULL, basis_out = NULL;
635   CeedQFunctionField *qf_fields;
636   CeedQFunction       qf;
637   CeedOperatorField  *op_fields;
638   CeedOperator_Sycl  *impl;
639 
640   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
641   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
642   CeedCallBackend(CeedQFunctionGetNumArgs(qf, &num_input_fields, &num_output_fields));
643 
644   // Determine active input basis
645   CeedCallBackend(CeedOperatorGetFields(op, NULL, &op_fields, NULL, NULL));
646   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_fields, NULL, NULL));
647   for (CeedInt i = 0; i < num_input_fields; i++) {
648     CeedVector vec;
649 
650     CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec));
651     if (vec == CEED_VECTOR_ACTIVE) {
652       CeedEvalMode eval_mode;
653       CeedBasis    basis;
654 
655       CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis));
656       CeedCheck(!basis_in || basis_in == basis, ceed, CEED_ERROR_BACKEND,
657                 "Backend does not implement operator diagonal assembly with multiple active bases");
658       if (!basis_in) CeedCallBackend(CeedBasisReferenceCopy(basis, &basis_in));
659       CeedCallBackend(CeedBasisDestroy(&basis));
660       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode));
661       switch (eval_mode) {
662         case CEED_EVAL_NONE:
663         case CEED_EVAL_INTERP:
664           CeedCallBackend(CeedRealloc(num_eval_mode_in + 1, &eval_mode_in));
665           eval_mode_in[num_eval_mode_in] = eval_mode;
666           num_eval_mode_in += 1;
667           break;
668         case CEED_EVAL_GRAD:
669           CeedCallBackend(CeedRealloc(num_eval_mode_in + dim, &eval_mode_in));
670           for (CeedInt d = 0; d < dim; d++) eval_mode_in[num_eval_mode_in + d] = eval_mode;
671           num_eval_mode_in += dim;
672           break;
673         case CEED_EVAL_WEIGHT:
674         case CEED_EVAL_DIV:
675         case CEED_EVAL_CURL:
676           break;  // Caught by QF Assembly
677       }
678     }
679     CeedCallBackend(CeedVectorDestroy(&vec));
680   }
681 
682   // Determine active output basis
683   CeedCallBackend(CeedOperatorGetFields(op, NULL, NULL, NULL, &op_fields));
684   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qf_fields));
685   for (CeedInt i = 0; i < num_output_fields; i++) {
686     CeedVector vec;
687 
688     CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec));
689     if (vec == CEED_VECTOR_ACTIVE) {
690       CeedEvalMode eval_mode;
691       CeedBasis    basis;
692 
693       CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis));
694       CeedCheck(!basis_out || basis_out == basis, ceed, CEED_ERROR_BACKEND,
695                 "Backend does not implement operator diagonal assembly with multiple active bases");
696       if (!basis_out) CeedCallBackend(CeedBasisReferenceCopy(basis, &basis_out));
697       CeedCallBackend(CeedBasisDestroy(&basis));
698       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode));
699       switch (eval_mode) {
700         case CEED_EVAL_NONE:
701         case CEED_EVAL_INTERP:
702           CeedCallBackend(CeedRealloc(num_eval_mode_in + 1, &eval_mode_in));
703           eval_mode_in[num_eval_mode_in] = eval_mode;
704           num_eval_mode_in += 1;
705           break;
706         case CEED_EVAL_GRAD:
707           CeedCallBackend(CeedRealloc(num_eval_mode_in + dim, &eval_mode_in));
708           for (CeedInt d = 0; d < dim; d++) eval_mode_in[num_eval_mode_in + d] = eval_mode;
709           num_eval_mode_in += dim;
710           break;
711         case CEED_EVAL_WEIGHT:
712         case CEED_EVAL_DIV:
713         case CEED_EVAL_CURL:
714           break;  // Caught by QF Assembly
715       }
716     }
717     CeedCallBackend(CeedVectorDestroy(&vec));
718   }
719 
720   // Operator data struct
721   CeedCallBackend(CeedOperatorGetData(op, &impl));
722   CeedCallBackend(CeedGetData(ceed, &sycl_data));
723   CeedCallBackend(CeedCalloc(1, &impl->diag));
724   CeedOperatorDiag_Sycl *diag = impl->diag;
725 
726   diag->basis_in          = basis_in;
727   diag->basis_out         = basis_out;
728   diag->h_eval_mode_in    = eval_mode_in;
729   diag->h_eval_mode_out   = eval_mode_out;
730   diag->num_eval_mode_in  = num_eval_mode_in;
731   diag->num_eval_mode_out = num_eval_mode_out;
732 
733   // Kernel parameters
734   CeedInt num_nodes, num_qpts;
735   CeedCallBackend(CeedBasisGetNumNodes(basis_in, &num_nodes));
736   CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis_in, &num_qpts));
737   diag->num_nodes = num_nodes;
738   diag->num_qpts  = num_qpts;
739   diag->num_comp  = num_comp;
740 
741   // Basis matrices
742   const CeedInt     i_len = num_qpts * num_nodes;
743   const CeedInt     g_len = num_qpts * num_nodes * dim;
744   const CeedScalar *interp_in, *interp_out, *grad_in, *grad_out;
745 
746   // CEED_EVAL_NONE
747   CeedScalar *identity      = NULL;
748   bool        has_eval_none = false;
749   for (CeedInt i = 0; i < num_eval_mode_in; i++) has_eval_none = has_eval_none || (eval_mode_in[i] == CEED_EVAL_NONE);
750   for (CeedInt i = 0; i < num_eval_mode_out; i++) has_eval_none = has_eval_none || (eval_mode_out[i] == CEED_EVAL_NONE);
751 
752   std::vector<sycl::event> e;
753 
754   if (!sycl_data->sycl_queue.is_in_order()) e = {sycl_data->sycl_queue.ext_oneapi_submit_barrier()};
755 
756   std::vector<sycl::event> copy_events;
757 
758   if (has_eval_none) {
759     CeedCallBackend(CeedCalloc(num_qpts * num_nodes, &identity));
760     for (CeedSize i = 0; i < (num_nodes < num_qpts ? num_nodes : num_qpts); i++) identity[i * num_nodes + i] = 1.0;
761     CeedCallSycl(ceed, diag->d_identity = sycl::malloc_device<CeedScalar>(i_len, sycl_data->sycl_device, sycl_data->sycl_context));
762     sycl::event identity_copy = sycl_data->sycl_queue.copy<CeedScalar>(identity, diag->d_identity, i_len, e);
763     copy_events.push_back(identity_copy);
764   }
765 
766   // CEED_EVAL_INTERP
767   CeedCallBackend(CeedBasisGetInterp(basis_in, &interp_in));
768   CeedCallSycl(ceed, diag->d_interp_in = sycl::malloc_device<CeedScalar>(i_len, sycl_data->sycl_device, sycl_data->sycl_context));
769   sycl::event interp_in_copy = sycl_data->sycl_queue.copy<CeedScalar>(interp_in, diag->d_interp_in, i_len, e);
770   copy_events.push_back(interp_in_copy);
771 
772   CeedCallBackend(CeedBasisGetInterp(basis_out, &interp_out));
773   CeedCallSycl(ceed, diag->d_interp_out = sycl::malloc_device<CeedScalar>(i_len, sycl_data->sycl_device, sycl_data->sycl_context));
774   sycl::event interp_out_copy = sycl_data->sycl_queue.copy<CeedScalar>(interp_out, diag->d_interp_out, i_len, e);
775   copy_events.push_back(interp_out_copy);
776 
777   // CEED_EVAL_GRAD
778   CeedCallBackend(CeedBasisGetGrad(basis_in, &grad_in));
779   CeedCallSycl(ceed, diag->d_grad_in = sycl::malloc_device<CeedScalar>(g_len, sycl_data->sycl_device, sycl_data->sycl_context));
780   sycl::event grad_in_copy = sycl_data->sycl_queue.copy<CeedScalar>(grad_in, diag->d_grad_in, g_len, e);
781   copy_events.push_back(grad_in_copy);
782 
783   CeedCallBackend(CeedBasisGetGrad(basis_out, &grad_out));
784   CeedCallSycl(ceed, diag->d_grad_out = sycl::malloc_device<CeedScalar>(g_len, sycl_data->sycl_device, sycl_data->sycl_context));
785   sycl::event grad_out_copy = sycl_data->sycl_queue.copy<CeedScalar>(grad_out, diag->d_grad_out, g_len, e);
786   copy_events.push_back(grad_out_copy);
787 
788   // Arrays of  eval_modes
789   CeedCallSycl(ceed, diag->d_eval_mode_in = sycl::malloc_device<CeedEvalMode>(num_eval_mode_in, sycl_data->sycl_device, sycl_data->sycl_context));
790   sycl::event eval_mode_in_copy = sycl_data->sycl_queue.copy<CeedEvalMode>(eval_mode_in, diag->d_eval_mode_in, num_eval_mode_in, e);
791   copy_events.push_back(eval_mode_in_copy);
792 
793   CeedCallSycl(ceed, diag->d_eval_mode_out = sycl::malloc_device<CeedEvalMode>(num_eval_mode_out, sycl_data->sycl_device, sycl_data->sycl_context));
794   sycl::event eval_mode_out_copy = sycl_data->sycl_queue.copy<CeedEvalMode>(eval_mode_out, diag->d_eval_mode_out, num_eval_mode_out, e);
795   copy_events.push_back(eval_mode_out_copy);
796 
797   // Restriction
798   {
799     CeedElemRestriction rstr_out;
800 
801     CeedCallBackend(CeedOperatorGetActiveElemRestrictions(op, NULL, &rstr_out));
802     diag->diag_rstr = rstr_out;
803     CeedCallBackend(CeedElemRestrictionDestroy(&rstr_out));
804   }
805 
806   // Wait for all copies to complete and handle exceptions
807   CeedCallSycl(ceed, sycl::event::wait_and_throw(copy_events));
808   return CEED_ERROR_SUCCESS;
809 }
810 
811 //------------------------------------------------------------------------------
812 //  Kernel for diagonal assembly
813 //------------------------------------------------------------------------------
814 static int CeedOperatorLinearDiagonal_Sycl(sycl::queue &sycl_queue, const bool is_point_block, const CeedInt num_elem,
815                                            const CeedOperatorDiag_Sycl *diag, const CeedScalar *assembled_qf_array, CeedScalar *elem_diag_array) {
816   const CeedSize      num_nodes         = diag->num_nodes;
817   const CeedSize      num_qpts          = diag->num_qpts;
818   const CeedSize      num_comp          = diag->num_comp;
819   const CeedSize      num_eval_mode_in  = diag->num_eval_mode_in;
820   const CeedSize      num_eval_mode_out = diag->num_eval_mode_out;
821   const CeedScalar   *identity          = diag->d_identity;
822   const CeedScalar   *interp_in         = diag->d_interp_in;
823   const CeedScalar   *grad_in           = diag->d_grad_in;
824   const CeedScalar   *interp_out        = diag->d_interp_out;
825   const CeedScalar   *grad_out          = diag->d_grad_out;
826   const CeedEvalMode *eval_mode_in      = diag->d_eval_mode_in;
827   const CeedEvalMode *eval_mode_out     = diag->d_eval_mode_out;
828 
829   sycl::range<1> kernel_range(num_elem * num_nodes);
830 
831   std::vector<sycl::event> e;
832 
833   if (!sycl_queue.is_in_order()) e = {sycl_queue.ext_oneapi_submit_barrier()};
834 
835   sycl_queue.parallel_for<CeedOperatorSyclLinearDiagonal>(kernel_range, e, [=](sycl::id<1> idx) {
836     const CeedInt tid = idx % num_nodes;
837     const CeedInt e   = idx / num_nodes;
838 
839     // Compute the diagonal of B^T D B
840     // Each element
841     CeedInt d_out = -1;
842     // Each basis eval mode pair
843     for (CeedSize e_out = 0; e_out < num_eval_mode_out; e_out++) {
844       const CeedScalar *bt = NULL;
845 
846       if (eval_mode_out[e_out] == CEED_EVAL_GRAD) ++d_out;
847       CeedOperatorGetBasisPointer_Sycl(&bt, eval_mode_out[e_out], identity, interp_out, &grad_out[d_out * num_qpts * num_nodes]);
848       CeedInt d_in = -1;
849 
850       for (CeedSize e_in = 0; e_in < num_eval_mode_in; e_in++) {
851         const CeedScalar *b = NULL;
852 
853         if (eval_mode_in[e_in] == CEED_EVAL_GRAD) ++d_in;
854         CeedOperatorGetBasisPointer_Sycl(&b, eval_mode_in[e_in], identity, interp_in, &grad_in[d_in * num_qpts * num_nodes]);
855         // Each component
856         for (CeedSize comp_out = 0; comp_out < num_comp; comp_out++) {
857           // Each qpoint/node pair
858           if (is_point_block) {
859             // Point Block Diagonal
860             for (CeedInt comp_in = 0; comp_in < num_comp; comp_in++) {
861               CeedScalar e_value = 0.0;
862 
863               for (CeedSize q = 0; q < num_qpts; q++) {
864                 const CeedScalar qf_value =
865                     assembled_qf_array[((((e_in * num_comp + comp_in) * num_eval_mode_out + e_out) * num_comp + comp_out) * num_elem + e) * num_qpts +
866                                        q];
867 
868                 e_value += bt[q * num_nodes + tid] * qf_value * b[q * num_nodes + tid];
869               }
870               elem_diag_array[((comp_out * num_comp + comp_in) * num_elem + e) * num_nodes + tid] += e_value;
871             }
872           } else {
873             // Diagonal Only
874             CeedScalar e_value = 0.0;
875 
876             for (CeedSize q = 0; q < num_qpts; q++) {
877               const CeedScalar qf_value =
878                   assembled_qf_array[((((e_in * num_comp + comp_out) * num_eval_mode_out + e_out) * num_comp + comp_out) * num_elem + e) * num_qpts +
879                                      q];
880               e_value += bt[q * num_nodes + tid] * qf_value * b[q * num_nodes + tid];
881             }
882             elem_diag_array[(comp_out * num_elem + e) * num_nodes + tid] += e_value;
883           }
884         }
885       }
886     }
887   });
888   return CEED_ERROR_SUCCESS;
889 }
890 
891 //------------------------------------------------------------------------------
892 // Assemble diagonal common code
893 //------------------------------------------------------------------------------
894 static inline int CeedOperatorAssembleDiagonalCore_Sycl(CeedOperator op, CeedVector assembled, CeedRequest *request, const bool is_point_block) {
895   Ceed                ceed;
896   Ceed_Sycl          *sycl_data;
897   CeedInt             num_elem;
898   CeedScalar         *elem_diag_array;
899   const CeedScalar   *assembled_qf_array;
900   CeedVector          assembled_qf = NULL;
901   CeedElemRestriction rstr         = NULL;
902   CeedOperator_Sycl  *impl;
903 
904   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
905   CeedCallBackend(CeedOperatorGetData(op, &impl));
906   CeedCallBackend(CeedGetData(ceed, &sycl_data));
907 
908   // Assemble QFunction
909   CeedCallBackend(CeedOperatorLinearAssembleQFunctionBuildOrUpdate(op, &assembled_qf, &rstr, request));
910   CeedCallBackend(CeedElemRestrictionDestroy(&rstr));
911 
912   // Setup
913   if (!impl->diag) {
914     CeedCallBackend(CeedOperatorAssembleDiagonalSetup_Sycl(op));
915   }
916   CeedOperatorDiag_Sycl *diag = impl->diag;
917 
918   assert(diag != NULL);
919 
920   // Restriction
921   if (is_point_block && !diag->point_block_diag_rstr) {
922     CeedCallBackend(CeedOperatorCreateActivePointBlockRestriction(diag->diag_rstr, &diag->point_block_diag_rstr));
923   }
924   CeedElemRestriction diag_rstr = is_point_block ? diag->point_block_diag_rstr : diag->diag_rstr;
925 
926   // Create diagonal vector
927   CeedVector elem_diag = is_point_block ? diag->point_block_elem_diag : diag->elem_diag;
928 
929   if (!elem_diag) {
930     CeedCallBackend(CeedElemRestrictionCreateVector(diag_rstr, NULL, &elem_diag));
931     if (is_point_block) diag->point_block_elem_diag = elem_diag;
932     else diag->elem_diag = elem_diag;
933   }
934   CeedCallBackend(CeedVectorSetValue(elem_diag, 0.0));
935 
936   // Assemble element operator diagonals
937   CeedCallBackend(CeedVectorGetArray(elem_diag, CEED_MEM_DEVICE, &elem_diag_array));
938   CeedCallBackend(CeedVectorGetArrayRead(assembled_qf, CEED_MEM_DEVICE, &assembled_qf_array));
939   CeedCallBackend(CeedElemRestrictionGetNumElements(diag_rstr, &num_elem));
940 
941   // Compute the diagonal of B^T D B
942   CeedCallBackend(CeedOperatorLinearDiagonal_Sycl(sycl_data->sycl_queue, is_point_block, num_elem, diag, assembled_qf_array, elem_diag_array));
943 
944   // Wait for queue to complete and handle exceptions
945   sycl_data->sycl_queue.wait_and_throw();
946 
947   // Restore arrays
948   CeedCallBackend(CeedVectorRestoreArray(elem_diag, &elem_diag_array));
949   CeedCallBackend(CeedVectorRestoreArrayRead(assembled_qf, &assembled_qf_array));
950 
951   // Assemble local operator diagonal
952   CeedCallBackend(CeedElemRestrictionApply(diag_rstr, CEED_TRANSPOSE, elem_diag, assembled, request));
953 
954   // Cleanup
955   CeedCallBackend(CeedVectorDestroy(&assembled_qf));
956   return CEED_ERROR_SUCCESS;
957 }
958 
959 //------------------------------------------------------------------------------
960 // Assemble Linear Diagonal
961 //------------------------------------------------------------------------------
962 static int CeedOperatorLinearAssembleAddDiagonal_Sycl(CeedOperator op, CeedVector assembled, CeedRequest *request) {
963   CeedCallBackend(CeedOperatorAssembleDiagonalCore_Sycl(op, assembled, request, false));
964   return CEED_ERROR_SUCCESS;
965 }
966 
967 //------------------------------------------------------------------------------
968 // Assemble Linear Point Block Diagonal
969 //------------------------------------------------------------------------------
970 static int CeedOperatorLinearAssembleAddPointBlockDiagonal_Sycl(CeedOperator op, CeedVector assembled, CeedRequest *request) {
971   CeedCallBackend(CeedOperatorAssembleDiagonalCore_Sycl(op, assembled, request, true));
972   return CEED_ERROR_SUCCESS;
973 }
974 
975 //------------------------------------------------------------------------------
976 // Single operator assembly setup
977 //------------------------------------------------------------------------------
978 static int CeedSingleOperatorAssembleSetup_Sycl(CeedOperator op) {
979   Ceed    ceed;
980   CeedInt num_input_fields, num_output_fields, num_eval_mode_in = 0, dim = 1, num_B_in_mats_to_load = 0, size_B_in = 0, num_eval_mode_out = 0,
981                                                num_B_out_mats_to_load = 0, size_B_out = 0, num_qpts = 0, elem_size = 0, num_elem, num_comp,
982                                                mat_start = 0;
983   CeedEvalMode       *eval_mode_in = NULL, *eval_mode_out = NULL;
984   const CeedScalar   *interp_in, *grad_in;
985   CeedElemRestriction rstr_in = NULL, rstr_out = NULL;
986   CeedBasis           basis_in = NULL, basis_out = NULL;
987   CeedQFunctionField *qf_fields;
988   CeedQFunction       qf;
989   CeedOperatorField  *input_fields, *output_fields;
990   CeedOperator_Sycl  *impl;
991 
992   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
993   CeedCallBackend(CeedOperatorGetData(op, &impl));
994 
995   // Get input and output fields
996   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &input_fields, &num_output_fields, &output_fields));
997 
998   // Determine active input basis eval mode
999   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
1000   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_fields, NULL, NULL));
1001   // Note that the kernel will treat each dimension of a gradient action separately;
1002   // i.e., when an active input has a CEED_EVAL_GRAD mode, num_ eval_mode_in will increment by dim.
1003   // 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,
1004   // so num_B_in_mats_to_load will be incremented by 1.
1005   for (CeedInt i = 0; i < num_input_fields; i++) {
1006     CeedVector vec;
1007 
1008     CeedCallBackend(CeedOperatorFieldGetVector(input_fields[i], &vec));
1009     if (vec == CEED_VECTOR_ACTIVE) {
1010       CeedEvalMode        eval_mode;
1011       CeedElemRestriction elem_rstr;
1012       CeedBasis           basis;
1013 
1014       CeedCallBackend(CeedOperatorFieldGetBasis(input_fields[i], &basis));
1015       CeedCheck(!basis_in || basis_in == basis, ceed, CEED_ERROR_BACKEND, "Backend does not implement operator assembly with multiple active bases");
1016       if (!basis_in) CeedCallBackend(CeedBasisReferenceCopy(basis, &basis_in));
1017       CeedCallBackend(CeedBasisDestroy(&basis));
1018       CeedCallBackend(CeedBasisGetDimension(basis_in, &dim));
1019       CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis_in, &num_qpts));
1020       CeedCallBackend(CeedOperatorFieldGetElemRestriction(input_fields[i], &elem_rstr));
1021       if (!rstr_in) CeedCallBackend(CeedElemRestrictionReferenceCopy(elem_rstr, &rstr_in));
1022       CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
1023       CeedCallBackend(CeedElemRestrictionGetElementSize(rstr_in, &elem_size));
1024       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode));
1025       if (eval_mode != CEED_EVAL_NONE) {
1026         CeedCallBackend(CeedRealloc(num_B_in_mats_to_load + 1, &eval_mode_in));
1027         eval_mode_in[num_B_in_mats_to_load] = eval_mode;
1028         num_B_in_mats_to_load += 1;
1029         if (eval_mode == CEED_EVAL_GRAD) {
1030           num_eval_mode_in += dim;
1031           size_B_in += dim * elem_size * num_qpts;
1032         } else {
1033           num_eval_mode_in += 1;
1034           size_B_in += elem_size * num_qpts;
1035         }
1036       }
1037     }
1038     CeedCallBackend(CeedVectorDestroy(&vec));
1039   }
1040 
1041   // Determine active output basis; basis_out and rstr_out only used if same as input, TODO
1042   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qf_fields));
1043   for (CeedInt i = 0; i < num_output_fields; i++) {
1044     CeedVector vec;
1045 
1046     CeedCallBackend(CeedOperatorFieldGetVector(output_fields[i], &vec));
1047     if (vec == CEED_VECTOR_ACTIVE) {
1048       CeedEvalMode        eval_mode;
1049       CeedElemRestriction elem_rstr;
1050       CeedBasis           basis;
1051 
1052       CeedCallBackend(CeedOperatorFieldGetBasis(output_fields[i], &basis));
1053       CeedCheck(!basis_out || basis_out == basis, ceed, CEED_ERROR_BACKEND,
1054                 "Backend does not implement operator assembly with multiple active bases");
1055       if (!basis_out) CeedCallBackend(CeedBasisReferenceCopy(basis, &basis_out));
1056       CeedCallBackend(CeedBasisDestroy(&basis));
1057       CeedCallBackend(CeedOperatorFieldGetElemRestriction(output_fields[i], &elem_rstr));
1058       if (!rstr_out) CeedCallBackend(CeedElemRestrictionReferenceCopy(elem_rstr, &rstr_out));
1059       CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
1060       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode));
1061       if (eval_mode != CEED_EVAL_NONE) {
1062         CeedCallBackend(CeedRealloc(num_B_out_mats_to_load + 1, &eval_mode_out));
1063         eval_mode_out[num_B_out_mats_to_load] = eval_mode;
1064         num_B_out_mats_to_load += 1;
1065         if (eval_mode == CEED_EVAL_GRAD) {
1066           num_eval_mode_out += dim;
1067           size_B_out += dim * elem_size * num_qpts;
1068         } else {
1069           num_eval_mode_out += 1;
1070           size_B_out += elem_size * num_qpts;
1071         }
1072       }
1073     }
1074     CeedCallBackend(CeedVectorDestroy(&vec));
1075   }
1076   CeedCheck(num_eval_mode_in > 0 && num_eval_mode_out > 0, ceed, CEED_ERROR_UNSUPPORTED, "Cannot assemble operator without inputs/outputs");
1077 
1078   CeedCallBackend(CeedElemRestrictionGetNumElements(rstr_in, &num_elem));
1079   CeedCallBackend(CeedElemRestrictionGetNumComponents(rstr_in, &num_comp));
1080 
1081   CeedCallBackend(CeedCalloc(1, &impl->asmb));
1082   CeedOperatorAssemble_Sycl *asmb = impl->asmb;
1083   asmb->num_elem                  = num_elem;
1084 
1085   Ceed_Sycl *sycl_data;
1086   CeedCallBackend(CeedGetData(ceed, &sycl_data));
1087 
1088   // Kernel setup
1089   int elems_per_block     = 1;
1090   asmb->elems_per_block   = elems_per_block;
1091   asmb->block_size_x      = elem_size;
1092   asmb->block_size_y      = elem_size;
1093   asmb->num_eval_mode_in  = num_eval_mode_in;
1094   asmb->num_eval_mode_out = num_eval_mode_out;
1095   asmb->num_qpts          = num_qpts;
1096   asmb->num_nodes         = elem_size;
1097   asmb->block_size        = elem_size * elem_size * elems_per_block;
1098   asmb->num_comp          = num_comp;
1099 
1100   // Build 'full' B matrices (not 1D arrays used for tensor-product matrices
1101   CeedCallBackend(CeedBasisGetInterp(basis_in, &interp_in));
1102   CeedCallBackend(CeedBasisGetGrad(basis_in, &grad_in));
1103 
1104   // Load into B_in, in order that they will be used in eval_mode
1105   CeedCallSycl(ceed, asmb->d_B_in = sycl::malloc_device<CeedScalar>(size_B_in, sycl_data->sycl_device, sycl_data->sycl_context));
1106   for (int i = 0; i < num_B_in_mats_to_load; i++) {
1107     CeedEvalMode eval_mode = eval_mode_in[i];
1108 
1109     if (eval_mode == CEED_EVAL_INTERP) {
1110       std::vector<sycl::event> e;
1111 
1112       if (!sycl_data->sycl_queue.is_in_order()) e = {sycl_data->sycl_queue.ext_oneapi_submit_barrier()};
1113       sycl_data->sycl_queue.copy<CeedScalar>(interp_in, &asmb->d_B_in[mat_start], elem_size * num_qpts, e);
1114       mat_start += elem_size * num_qpts;
1115     } else if (eval_mode == CEED_EVAL_GRAD) {
1116       std::vector<sycl::event> e;
1117 
1118       if (!sycl_data->sycl_queue.is_in_order()) e = {sycl_data->sycl_queue.ext_oneapi_submit_barrier()};
1119       sycl_data->sycl_queue.copy<CeedScalar>(grad_in, &asmb->d_B_in[mat_start], dim * elem_size * num_qpts, e);
1120       mat_start += dim * elem_size * num_qpts;
1121     }
1122   }
1123 
1124   const CeedScalar *interp_out, *grad_out;
1125   // Note that this function currently assumes 1 basis, so this should always be true
1126   // for now
1127   if (basis_out == basis_in) {
1128     interp_out = interp_in;
1129     grad_out   = grad_in;
1130   } else {
1131     CeedCallBackend(CeedBasisGetInterp(basis_out, &interp_out));
1132     CeedCallBackend(CeedBasisGetGrad(basis_out, &grad_out));
1133   }
1134 
1135   // Load into B_out, in order that they will be used in eval_mode
1136   mat_start = 0;
1137   CeedCallSycl(ceed, asmb->d_B_out = sycl::malloc_device<CeedScalar>(size_B_out, sycl_data->sycl_device, sycl_data->sycl_context));
1138   for (int i = 0; i < num_B_out_mats_to_load; i++) {
1139     CeedEvalMode eval_mode = eval_mode_out[i];
1140 
1141     if (eval_mode == CEED_EVAL_INTERP) {
1142       std::vector<sycl::event> e;
1143 
1144       if (!sycl_data->sycl_queue.is_in_order()) e = {sycl_data->sycl_queue.ext_oneapi_submit_barrier()};
1145       sycl_data->sycl_queue.copy<CeedScalar>(interp_out, &asmb->d_B_out[mat_start], elem_size * num_qpts, e);
1146       mat_start += elem_size * num_qpts;
1147     } else if (eval_mode == CEED_EVAL_GRAD) {
1148       std::vector<sycl::event> e;
1149 
1150       if (!sycl_data->sycl_queue.is_in_order()) e = {sycl_data->sycl_queue.ext_oneapi_submit_barrier()};
1151       sycl_data->sycl_queue.copy<CeedScalar>(grad_out, &asmb->d_B_out[mat_start], dim * elem_size * num_qpts, e);
1152       mat_start += dim * elem_size * num_qpts;
1153     }
1154   }
1155   CeedCallBackend(CeedElemRestrictionDestroy(&rstr_in));
1156   CeedCallBackend(CeedElemRestrictionDestroy(&rstr_out));
1157   CeedCallBackend(CeedBasisDestroy(&basis_in));
1158   CeedCallBackend(CeedBasisDestroy(&basis_out));
1159   return CEED_ERROR_SUCCESS;
1160 }
1161 
1162 //------------------------------------------------------------------------------
1163 // Matrix assembly kernel for low-order elements (3D thread block)
1164 //------------------------------------------------------------------------------
1165 static int CeedOperatorLinearAssemble_Sycl(sycl::queue &sycl_queue, const CeedOperator_Sycl *impl, const CeedScalar *qf_array,
1166                                            CeedScalar *values_array) {
1167   // This kernels assumes B_in and B_out have the same number of quadrature points and basis points.
1168   // TODO: expand to more general cases
1169   CeedOperatorAssemble_Sycl *asmb              = impl->asmb;
1170   const CeedInt              num_elem          = asmb->num_elem;
1171   const CeedSize             num_nodes         = asmb->num_nodes;
1172   const CeedSize             num_comp          = asmb->num_comp;
1173   const CeedSize             num_qpts          = asmb->num_qpts;
1174   const CeedSize             num_eval_mode_in  = asmb->num_eval_mode_in;
1175   const CeedSize             num_eval_mode_out = asmb->num_eval_mode_out;
1176 
1177   // Strides for final output ordering, determined by the reference (inference) implementation of the symbolic assembly, slowest --> fastest: element,
1178   // comp_in, comp_out, node_row, node_col
1179   const CeedSize comp_out_stride = num_nodes * num_nodes;
1180   const CeedSize comp_in_stride  = comp_out_stride * num_comp;
1181   const CeedSize e_stride        = comp_in_stride * num_comp;
1182   // Strides for QF array, slowest --> fastest:  eval_mode_in, comp_in,  eval_mode_out, comp_out, elem, qpt
1183   const CeedSize q_e_stride             = num_qpts;
1184   const CeedSize q_comp_out_stride      = num_elem * q_e_stride;
1185   const CeedSize q_eval_mode_out_stride = q_comp_out_stride * num_comp;
1186   const CeedSize q_comp_in_stride       = q_eval_mode_out_stride * num_eval_mode_out;
1187   const CeedSize q_eval_mode_in_stride  = q_comp_in_stride * num_comp;
1188 
1189   CeedScalar *B_in, *B_out;
1190   B_in                       = asmb->d_B_in;
1191   B_out                      = asmb->d_B_out;
1192   const CeedInt block_size_x = asmb->block_size_x;
1193   const CeedInt block_size_y = asmb->block_size_y;
1194 
1195   sycl::range<3> kernel_range(num_elem, block_size_y, block_size_x);
1196 
1197   std::vector<sycl::event> e;
1198 
1199   if (!sycl_queue.is_in_order()) e = {sycl_queue.ext_oneapi_submit_barrier()};
1200   sycl_queue.parallel_for<CeedOperatorSyclLinearAssemble>(kernel_range, e, [=](sycl::id<3> idx) {
1201     const int e = idx.get(0);  // Element index
1202     const int l = idx.get(1);  // The output column index of each B^TDB operation
1203     const int i = idx.get(2);  // The output row index of each B^TDB operation
1204                                // such that we have (Bout^T)_ij D_jk Bin_kl = C_il
1205     for (CeedSize comp_in = 0; comp_in < num_comp; comp_in++) {
1206       for (CeedSize comp_out = 0; comp_out < num_comp; comp_out++) {
1207         CeedScalar result        = 0.0;
1208         CeedSize   qf_index_comp = q_comp_in_stride * comp_in + q_comp_out_stride * comp_out + q_e_stride * e;
1209 
1210         for (CeedSize eval_mode_in = 0; eval_mode_in < num_eval_mode_in; eval_mode_in++) {
1211           CeedSize b_in_index = eval_mode_in * num_qpts * num_nodes;
1212 
1213           for (CeedSize eval_mode_out = 0; eval_mode_out < num_eval_mode_out; eval_mode_out++) {
1214             CeedSize b_out_index = eval_mode_out * num_qpts * num_nodes;
1215             CeedSize qf_index    = qf_index_comp + q_eval_mode_out_stride * eval_mode_out + q_eval_mode_in_stride * eval_mode_in;
1216 
1217             // Perform the B^T D B operation for this 'chunk' of D (the qf_array)
1218             for (CeedSize j = 0; j < num_qpts; j++) {
1219               result += B_out[b_out_index + j * num_nodes + i] * qf_array[qf_index + j] * B_in[b_in_index + j * num_nodes + l];
1220             }
1221           }  // end of  eval_mode_out
1222         }  // end of  eval_mode_in
1223         CeedSize val_index = comp_in_stride * comp_in + comp_out_stride * comp_out + e_stride * e + num_nodes * i + l;
1224 
1225         values_array[val_index] = result;
1226       }  // end of out component
1227     }  // end of in component
1228   });
1229   return CEED_ERROR_SUCCESS;
1230 }
1231 
1232 //------------------------------------------------------------------------------
1233 // Fallback kernel for larger orders (1D thread block)
1234 //------------------------------------------------------------------------------
1235 /*
1236 static int CeedOperatorLinearAssembleFallback_Sycl(sycl::queue &sycl_queue, const CeedOperator_Sycl *impl, const CeedScalar *qf_array,
1237                                                    CeedScalar *values_array) {
1238   // This kernel assumes B_in and B_out have the same number of quadrature points and basis points.
1239   // TODO: expand to more general cases
1240   CeedOperatorAssemble_Sycl *asmb        = impl->asmb;
1241   const CeedInt              num_elem       = asmb->num_elem;
1242   const CeedInt              num_nodes      = asmb->num_nodes;
1243   const CeedInt              num_comp       = asmb->num_comp;
1244   const CeedInt              num_qpts       = asmb->num_qpts;
1245   const CeedInt              num_eval_mode_in  = asmb->num_eval_mode_in;
1246   const CeedInt              num_eval_mode_out = asmb->num_eval_mode_out;
1247 
1248   // Strides for final output ordering, determined by the reference (interface) implementation of the symbolic assembly, slowest --> fastest: elememt,
1249   // comp_in, comp_out, node_row, node_col
1250   const CeedInt comp_out_stride = num_nodes * num_nodes;
1251   const CeedInt comp_in_stride  = comp_out_stride * num_comp;
1252   const CeedInt e_stride        = comp_in_stride * num_comp;
1253   // Strides for QF array, slowest --> fastest:  eval_mode_in, comp_in,  eval_mode_out, comp_out, elem, qpt
1254   const CeedInt q_e_stride         = num_qpts;
1255   const CeedInt q_comp_out_stride  = num_elem * q_e_stride;
1256   const CeedInt q_eval_mode_out_stride = q_comp_out_stride * num_comp;
1257   const CeedInt q_comp_in_stride   = q_eval_mode_out_stride * num_eval_mode_out;
1258   const CeedInt q_eval_mode_in_stride  = q_comp_in_stride * num_comp;
1259 
1260   CeedScalar *B_in, *B_out;
1261   B_in                        = asmb->d_B_in;
1262   B_out                       = asmb->d_B_out;
1263   const CeedInt elems_per_block = asmb->elems_per_block;
1264   const CeedInt block_size_x  = asmb->block_size_x;
1265   const CeedInt block_size_y  = asmb->block_size_y;  // This will be 1 for the fallback kernel
1266 
1267   const CeedInt     grid = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
1268   sycl::range<3>    local_range(block_size_x, block_size_y, elems_per_block);
1269   sycl::range<3>    global_range(grid * block_size_x, block_size_y, elems_per_block);
1270   sycl::nd_range<3> kernel_range(global_range, local_range);
1271 
1272   sycl_queue.parallel_for<CeedOperatorSyclLinearAssembleFallback>(kernel_range, [=](sycl::nd_item<3> work_item) {
1273     const CeedInt blockIdx  = work_item.get_group(0);
1274     const CeedInt gridDimx  = work_item.get_group_range(0);
1275     const CeedInt threadIdx = work_item.get_local_id(0);
1276     const CeedInt threadIdz = work_item.get_local_id(2);
1277     const CeedInt blockDimz = work_item.get_local_range(2);
1278 
1279     const int l = threadIdx;  // The output column index of each B^TDB operation
1280                               // such that we have (Bout^T)_ij D_jk Bin_kl = C_il
1281     for (CeedInt e = blockIdx * blockDimz + threadIdz; e < num_elem; e += gridDimx * blockDimz) {
1282       for (CeedInt comp_in = 0; comp_in < num_comp; comp_in++) {
1283         for (CeedInt comp_out = 0; comp_out < num_comp; comp_out++) {
1284           for (CeedInt i = 0; i < num_nodes; i++) {
1285             CeedScalar result        = 0.0;
1286             CeedInt    qf_index_comp = q_comp_in_stride * comp_in + q_comp_out_stride * comp_out + q_e_stride * e;
1287             for (CeedInt  eval_mode_in = 0;  eval_mode_in < num_eval_mode_in;  eval_mode_in++) {
1288               CeedInt b_in_index =  eval_mode_in * num_qpts * num_nodes;
1289               for (CeedInt  eval_mode_out = 0;  eval_mode_out < num_eval_mode_out;  eval_mode_out++) {
1290                 CeedInt b_out_index =  eval_mode_out * num_qpts * num_nodes;
1291                 CeedInt qf_index    = qf_index_comp + q_eval_mode_out_stride *  eval_mode_out + q_eval_mode_in_stride *  eval_mode_in;
1292                 // Perform the B^T D B operation for this 'chunk' of D (the qf_array)
1293                 for (CeedInt j = 0; j < num_qpts; j++) {
1294                   result += B_out[b_out_index + j * num_nodes + i] * qf_array[qf_index + j] * B_in[b_in_index + j * num_nodes + l];
1295                 }
1296               }  // end of  eval_mode_out
1297             }    // end of  eval_mode_in
1298             CeedInt val_index       = comp_in_stride * comp_in + comp_out_stride * comp_out + e_stride * e + num_nodes * i + l;
1299             values_array[val_index] = result;
1300           }  // end of loop over element node index, i
1301         }    // end of out component
1302       }      // end of in component
1303     }        // end of element loop
1304   });
1305   return CEED_ERROR_SUCCESS;
1306 }*/
1307 
1308 //------------------------------------------------------------------------------
1309 // Assemble matrix data for COO matrix of assembled operator.
1310 // The sparsity pattern is set by CeedOperatorLinearAssembleSymbolic.
1311 //
1312 // Note that this (and other assembly routines) currently assume only one active
1313 // input restriction/basis per operator (could have multiple basis eval modes).
1314 // TODO: allow multiple active input restrictions/basis objects
1315 //------------------------------------------------------------------------------
1316 static int CeedSingleOperatorAssemble_Sycl(CeedOperator op, CeedInt offset, CeedVector values) {
1317   Ceed                ceed;
1318   Ceed_Sycl          *sycl_data;
1319   CeedScalar         *values_array;
1320   const CeedScalar   *qf_array;
1321   CeedVector          assembled_qf = NULL;
1322   CeedElemRestriction rstr_q       = NULL;
1323   CeedOperator_Sycl  *impl;
1324 
1325   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
1326   CeedCallBackend(CeedOperatorGetData(op, &impl));
1327   CeedCallBackend(CeedGetData(ceed, &sycl_data));
1328 
1329   // Setup
1330   if (!impl->asmb) {
1331     CeedCallBackend(CeedSingleOperatorAssembleSetup_Sycl(op));
1332     assert(impl->asmb != NULL);
1333   }
1334 
1335   // Assemble QFunction
1336   CeedCallBackend(CeedOperatorLinearAssembleQFunctionBuildOrUpdate(op, &assembled_qf, &rstr_q, CEED_REQUEST_IMMEDIATE));
1337   CeedCallBackend(CeedElemRestrictionDestroy(&rstr_q));
1338   CeedCallBackend(CeedVectorGetArrayWrite(values, CEED_MEM_DEVICE, &values_array));
1339   values_array += offset;
1340   CeedCallBackend(CeedVectorGetArrayRead(assembled_qf, CEED_MEM_DEVICE, &qf_array));
1341 
1342   // Compute B^T D B
1343   CeedCallBackend(CeedOperatorLinearAssemble_Sycl(sycl_data->sycl_queue, impl, qf_array, values_array));
1344 
1345   // Wait for kernels to be completed
1346   // Kris: Review if this is necessary -- enqueing an async barrier may be sufficient
1347   sycl_data->sycl_queue.wait_and_throw();
1348 
1349   // Restore arrays
1350   CeedCallBackend(CeedVectorRestoreArray(values, &values_array));
1351   CeedCallBackend(CeedVectorRestoreArrayRead(assembled_qf, &qf_array));
1352 
1353   // Cleanup
1354   CeedCallBackend(CeedVectorDestroy(&assembled_qf));
1355   return CEED_ERROR_SUCCESS;
1356 }
1357 
1358 //------------------------------------------------------------------------------
1359 // Create operator
1360 //------------------------------------------------------------------------------
1361 int CeedOperatorCreate_Sycl(CeedOperator op) {
1362   Ceed               ceed;
1363   CeedOperator_Sycl *impl;
1364 
1365   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
1366 
1367   CeedCallBackend(CeedCalloc(1, &impl));
1368   CeedCallBackend(CeedOperatorSetData(op, impl));
1369 
1370   CeedCallBackend(CeedSetBackendFunctionCpp(ceed, "Operator", op, "LinearAssembleQFunction", CeedOperatorLinearAssembleQFunction_Sycl));
1371   CeedCallBackend(CeedSetBackendFunctionCpp(ceed, "Operator", op, "LinearAssembleQFunctionUpdate", CeedOperatorLinearAssembleQFunctionUpdate_Sycl));
1372   CeedCallBackend(CeedSetBackendFunctionCpp(ceed, "Operator", op, "LinearAssembleAddDiagonal", CeedOperatorLinearAssembleAddDiagonal_Sycl));
1373   CeedCallBackend(
1374       CeedSetBackendFunctionCpp(ceed, "Operator", op, "LinearAssembleAddPointBlockDiagonal", CeedOperatorLinearAssembleAddPointBlockDiagonal_Sycl));
1375   CeedCallBackend(CeedSetBackendFunctionCpp(ceed, "Operator", op, "LinearAssembleSingle", CeedSingleOperatorAssemble_Sycl));
1376   CeedCallBackend(CeedSetBackendFunctionCpp(ceed, "Operator", op, "ApplyAdd", CeedOperatorApplyAdd_Sycl));
1377   CeedCallBackend(CeedSetBackendFunctionCpp(ceed, "Operator", op, "Destroy", CeedOperatorDestroy_Sycl));
1378   return CEED_ERROR_SUCCESS;
1379 }
1380 
1381 //------------------------------------------------------------------------------
1382