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