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