xref: /libCEED/backends/sycl-ref/ceed-sycl-ref-operator.sycl.cpp (revision cde30410219da4e6c581f011cb04f3bf412e880c)
1 // Copyright (c) 2017-2022, 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, NULL, 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,
557                                                      num_active_in * num_active_out * num_elem * Q, strides, rstr));
558     // Create assembled vector
559     CeedCallBackend(CeedVectorCreate(ceed_parent, l_size, assembled));
560   }
561   CeedCallBackend(CeedVectorSetValue(*assembled, 0.0));
562   CeedCallBackend(CeedVectorGetArray(*assembled, CEED_MEM_DEVICE, &assembled_array));
563 
564   // Input basis apply
565   CeedCallBackend(CeedOperatorInputBasis_Sycl(num_elem, qf_input_fields, op_input_fields, num_input_fields, true, e_data, impl));
566 
567   // Assemble QFunction
568   for (CeedInt in = 0; in < num_active_in; in++) {
569     // Set Inputs
570     CeedCallBackend(CeedVectorSetValue(active_in[in], 1.0));
571     if (num_active_in > 1) {
572       CeedCallBackend(CeedVectorSetValue(active_in[(in + num_active_in - 1) % num_active_in], 0.0));
573     }
574     // Set Outputs
575     for (CeedInt out = 0; out < num_output_fields; out++) {
576       CeedVector vec;
577 
578       // Get output vector
579       CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[out], &vec));
580       // Check if active output
581       if (vec == CEED_VECTOR_ACTIVE) {
582         CeedCallBackend(CeedVectorSetArray(impl->q_vecs_out[out], CEED_MEM_DEVICE, CEED_USE_POINTER, assembled_array));
583         CeedCallBackend(CeedQFunctionFieldGetSize(qf_output_fields[out], &size));
584         assembled_array += size * Q * num_elem;  // Advance the pointer by the size of the output
585       }
586     }
587     // Apply QFunction
588     CeedCallBackend(CeedQFunctionApply(qf, Q * num_elem, impl->q_vecs_in, impl->q_vecs_out));
589   }
590 
591   // Un-set output Qvecs to prevent accidental overwrite of Assembled
592   for (CeedInt out = 0; out < num_output_fields; out++) {
593     CeedVector vec;
594 
595     // Get output vector
596     CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[out], &vec));
597     // Check if active output
598     if (vec == CEED_VECTOR_ACTIVE) {
599       CeedCallBackend(CeedVectorTakeArray(impl->q_vecs_out[out], CEED_MEM_DEVICE, NULL));
600     }
601   }
602 
603   // Restore input arrays
604   CeedCallBackend(CeedOperatorRestoreInputs_Sycl(num_input_fields, qf_input_fields, op_input_fields, true, e_data, impl));
605 
606   // Restore output
607   CeedCallBackend(CeedVectorRestoreArray(*assembled, &assembled_array));
608   return CEED_ERROR_SUCCESS;
609 }
610 
611 //------------------------------------------------------------------------------
612 // Assemble Linear QFunction
613 //------------------------------------------------------------------------------
614 static int CeedOperatorLinearAssembleQFunction_Sycl(CeedOperator op, CeedVector *assembled, CeedElemRestriction *rstr, CeedRequest *request) {
615   return CeedOperatorLinearAssembleQFunctionCore_Sycl(op, true, assembled, rstr, request);
616 }
617 
618 //------------------------------------------------------------------------------
619 // Update Assembled Linear QFunction
620 //------------------------------------------------------------------------------
621 static int CeedOperatorLinearAssembleQFunctionUpdate_Sycl(CeedOperator op, CeedVector assembled, CeedElemRestriction rstr, CeedRequest *request) {
622   return CeedOperatorLinearAssembleQFunctionCore_Sycl(op, false, &assembled, &rstr, request);
623 }
624 
625 //------------------------------------------------------------------------------
626 // Assemble diagonal setup
627 //------------------------------------------------------------------------------
628 static inline int CeedOperatorAssembleDiagonalSetup_Sycl(CeedOperator op) {
629   Ceed                ceed;
630   Ceed_Sycl          *sycl_data;
631   CeedInt             num_input_fields, num_output_fields, num_e_mode_in = 0, num_comp = 0, dim = 1, num_e_mode_out = 0;
632   CeedEvalMode       *e_mode_in = NULL, *e_mode_out = NULL;
633   CeedBasis           basis_in = NULL, basis_out = NULL;
634   CeedElemRestriction rstr_in = NULL, rstr_out = NULL;
635   CeedQFunctionField *qf_fields;
636   CeedQFunction       qf;
637   CeedOperatorField  *op_fields;
638   CeedOperator_Sycl  *impl;
639 
640   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
641   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
642   CeedCallBackend(CeedQFunctionGetNumArgs(qf, &num_input_fields, &num_output_fields));
643 
644   // Determine active input basis
645   CeedCallBackend(CeedOperatorGetFields(op, NULL, &op_fields, NULL, NULL));
646   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_fields, NULL, NULL));
647   for (CeedInt i = 0; i < num_input_fields; i++) {
648     CeedVector vec;
649 
650     CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec));
651     if (vec == CEED_VECTOR_ACTIVE) {
652       CeedEvalMode        e_mode;
653       CeedElemRestriction rstr;
654 
655       CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis_in));
656       CeedCallBackend(CeedBasisGetNumComponents(basis_in, &num_comp));
657       CeedCallBackend(CeedBasisGetDimension(basis_in, &dim));
658       CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_fields[i], &rstr));
659       CeedCheck(!rstr_in || rstr_in == rstr, ceed, CEED_ERROR_BACKEND,
660                 "Backend does not implement multi-field non-composite operator diagonal assembly");
661       rstr_in = rstr;
662       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &e_mode));
663       switch (e_mode) {
664         case CEED_EVAL_NONE:
665         case CEED_EVAL_INTERP:
666           CeedCallBackend(CeedRealloc(num_e_mode_in + 1, &e_mode_in));
667           e_mode_in[num_e_mode_in] = e_mode;
668           num_e_mode_in += 1;
669           break;
670         case CEED_EVAL_GRAD:
671           CeedCallBackend(CeedRealloc(num_e_mode_in + dim, &e_mode_in));
672           for (CeedInt d = 0; d < dim; d++) e_mode_in[num_e_mode_in + d] = e_mode;
673           num_e_mode_in += dim;
674           break;
675         case CEED_EVAL_WEIGHT:
676         case CEED_EVAL_DIV:
677         case CEED_EVAL_CURL:
678           break;  // Caught by QF Assembly
679       }
680     }
681   }
682 
683   // Determine active output basis
684   CeedCallBackend(CeedOperatorGetFields(op, NULL, NULL, NULL, &op_fields));
685   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qf_fields));
686   for (CeedInt i = 0; i < num_output_fields; i++) {
687     CeedVector vec;
688 
689     CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec));
690     if (vec == CEED_VECTOR_ACTIVE) {
691       CeedEvalMode        e_mode;
692       CeedElemRestriction rstr;
693 
694       CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis_out));
695       CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_fields[i], &rstr));
696       CeedCheck(!rstr_out || rstr_out == rstr, ceed, CEED_ERROR_BACKEND,
697                 "Backend does not implement multi-field non-composite operator diagonal assembly");
698       rstr_out = rstr;
699       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &e_mode));
700       switch (e_mode) {
701         case CEED_EVAL_NONE:
702         case CEED_EVAL_INTERP:
703           CeedCallBackend(CeedRealloc(num_e_mode_out + 1, &e_mode_out));
704           e_mode_out[num_e_mode_out] = e_mode;
705           num_e_mode_out += 1;
706           break;
707         case CEED_EVAL_GRAD:
708           CeedCallBackend(CeedRealloc(num_e_mode_out + dim, &e_mode_out));
709           for (CeedInt d = 0; d < dim; d++) e_mode_out[num_e_mode_out + d] = e_mode;
710           num_e_mode_out += dim;
711           break;
712         case CEED_EVAL_WEIGHT:
713         case CEED_EVAL_DIV:
714         case CEED_EVAL_CURL:
715           break;  // Caught by QF Assembly
716       }
717     }
718   }
719 
720   // Operator data struct
721   CeedCallBackend(CeedOperatorGetData(op, &impl));
722   CeedCallBackend(CeedGetData(ceed, &sycl_data));
723   CeedCallBackend(CeedCalloc(1, &impl->diag));
724   CeedOperatorDiag_Sycl *diag = impl->diag;
725 
726   diag->basis_in       = basis_in;
727   diag->basis_out      = basis_out;
728   diag->h_e_mode_in    = e_mode_in;
729   diag->h_e_mode_out   = e_mode_out;
730   diag->num_e_mode_in  = num_e_mode_in;
731   diag->num_e_mode_out = num_e_mode_out;
732 
733   // Kernel parameters
734   CeedInt num_nodes, num_qpts;
735   CeedCallBackend(CeedBasisGetNumNodes(basis_in, &num_nodes));
736   CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis_in, &num_qpts));
737   diag->num_nodes = num_nodes;
738   diag->num_qpts  = num_qpts;
739   diag->num_comp  = num_comp;
740 
741   // Basis matrices
742   const CeedInt     i_len = num_qpts * num_nodes;
743   const CeedInt     g_len = num_qpts * num_nodes * dim;
744   const CeedScalar *interp_in, *interp_out, *grad_in, *grad_out;
745 
746   // CEED_EVAL_NONE
747   CeedScalar *identity      = NULL;
748   bool        has_eval_none = false;
749   for (CeedInt i = 0; i < num_e_mode_in; i++) has_eval_none = has_eval_none || (e_mode_in[i] == CEED_EVAL_NONE);
750   for (CeedInt i = 0; i < num_e_mode_out; i++) has_eval_none = has_eval_none || (e_mode_out[i] == CEED_EVAL_NONE);
751 
752   // Order queue
753   sycl::event e = sycl_data->sycl_queue.ext_oneapi_submit_barrier();
754 
755   std::vector<sycl::event> copy_events;
756   if (has_eval_none) {
757     CeedCallBackend(CeedCalloc(num_qpts * num_nodes, &identity));
758     for (CeedSize i = 0; i < (num_nodes < num_qpts ? num_nodes : num_qpts); i++) identity[i * num_nodes + i] = 1.0;
759     CeedCallSycl(ceed, diag->d_identity = sycl::malloc_device<CeedScalar>(i_len, sycl_data->sycl_device, sycl_data->sycl_context));
760     sycl::event identity_copy = sycl_data->sycl_queue.copy<CeedScalar>(identity, diag->d_identity, i_len, {e});
761     copy_events.push_back(identity_copy);
762   }
763 
764   // CEED_EVAL_INTERP
765   CeedCallBackend(CeedBasisGetInterp(basis_in, &interp_in));
766   CeedCallSycl(ceed, diag->d_interp_in = sycl::malloc_device<CeedScalar>(i_len, sycl_data->sycl_device, sycl_data->sycl_context));
767   sycl::event interp_in_copy = sycl_data->sycl_queue.copy<CeedScalar>(interp_in, diag->d_interp_in, i_len, {e});
768   copy_events.push_back(interp_in_copy);
769 
770   CeedCallBackend(CeedBasisGetInterp(basis_out, &interp_out));
771   CeedCallSycl(ceed, diag->d_interp_out = sycl::malloc_device<CeedScalar>(i_len, sycl_data->sycl_device, sycl_data->sycl_context));
772   sycl::event interp_out_copy = sycl_data->sycl_queue.copy<CeedScalar>(interp_out, diag->d_interp_out, i_len, {e});
773   copy_events.push_back(interp_out_copy);
774 
775   // CEED_EVAL_GRAD
776   CeedCallBackend(CeedBasisGetGrad(basis_in, &grad_in));
777   CeedCallSycl(ceed, diag->d_grad_in = sycl::malloc_device<CeedScalar>(g_len, sycl_data->sycl_device, sycl_data->sycl_context));
778   sycl::event grad_in_copy = sycl_data->sycl_queue.copy<CeedScalar>(grad_in, diag->d_grad_in, g_len, {e});
779   copy_events.push_back(grad_in_copy);
780 
781   CeedCallBackend(CeedBasisGetGrad(basis_out, &grad_out));
782   CeedCallSycl(ceed, diag->d_grad_out = sycl::malloc_device<CeedScalar>(g_len, sycl_data->sycl_device, sycl_data->sycl_context));
783   sycl::event grad_out_copy = sycl_data->sycl_queue.copy<CeedScalar>(grad_out, diag->d_grad_out, g_len, {e});
784   copy_events.push_back(grad_out_copy);
785 
786   // Arrays of  e_modes
787   CeedCallSycl(ceed, diag->d_e_mode_in = sycl::malloc_device<CeedEvalMode>(num_e_mode_in, sycl_data->sycl_device, sycl_data->sycl_context));
788   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});
789   copy_events.push_back(e_mode_in_copy);
790 
791   CeedCallSycl(ceed, diag->d_e_mode_out = sycl::malloc_device<CeedEvalMode>(num_e_mode_out, sycl_data->sycl_device, sycl_data->sycl_context));
792   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});
793   copy_events.push_back(e_mode_out_copy);
794 
795   // Restriction
796   diag->diag_rstr = rstr_out;
797 
798   // Wait for all copies to complete and handle exceptions
799   CeedCallSycl(ceed, sycl::event::wait_and_throw(copy_events));
800   return CEED_ERROR_SUCCESS;
801 }
802 
803 //------------------------------------------------------------------------------
804 //  Kernel for diagonal assembly
805 //------------------------------------------------------------------------------
806 static int CeedOperatorLinearDiagonal_Sycl(sycl::queue &sycl_queue, const bool is_point_block, const CeedInt num_elem,
807                                            const CeedOperatorDiag_Sycl *diag, const CeedScalar *assembled_qf_array, CeedScalar *elem_diag_array) {
808   const CeedSize      num_nodes      = diag->num_nodes;
809   const CeedSize      num_qpts       = diag->num_qpts;
810   const CeedSize      num_comp       = diag->num_comp;
811   const CeedSize      num_e_mode_in  = diag->num_e_mode_in;
812   const CeedSize      num_e_mode_out = diag->num_e_mode_out;
813   const CeedScalar   *identity       = diag->d_identity;
814   const CeedScalar   *interp_in      = diag->d_interp_in;
815   const CeedScalar   *grad_in        = diag->d_grad_in;
816   const CeedScalar   *interp_out     = diag->d_interp_out;
817   const CeedScalar   *grad_out       = diag->d_grad_out;
818   const CeedEvalMode *e_mode_in      = diag->d_e_mode_in;
819   const CeedEvalMode *e_mode_out     = diag->d_e_mode_out;
820 
821   sycl::range<1> kernel_range(num_elem * num_nodes);
822 
823   // Order queue
824   sycl::event e = sycl_queue.ext_oneapi_submit_barrier();
825   sycl_queue.parallel_for<CeedOperatorSyclLinearDiagonal>(kernel_range, {e}, [=](sycl::id<1> idx) {
826     const CeedInt tid = idx % num_nodes;
827     const CeedInt e   = idx / num_nodes;
828 
829     // Compute the diagonal of B^T D B
830     // Each element
831     CeedInt d_out = -1;
832     // Each basis eval mode pair
833     for (CeedSize e_out = 0; e_out < num_e_mode_out; e_out++) {
834       const CeedScalar *bt = NULL;
835 
836       if (e_mode_out[e_out] == CEED_EVAL_GRAD) ++d_out;
837       CeedOperatorGetBasisPointer_Sycl(&bt, e_mode_out[e_out], identity, interp_out, &grad_out[d_out * num_qpts * num_nodes]);
838       CeedInt d_in = -1;
839 
840       for (CeedSize e_in = 0; e_in < num_e_mode_in; e_in++) {
841         const CeedScalar *b = NULL;
842 
843         if (e_mode_in[e_in] == CEED_EVAL_GRAD) ++d_in;
844         CeedOperatorGetBasisPointer_Sycl(&b, e_mode_in[e_in], identity, interp_in, &grad_in[d_in * num_qpts * num_nodes]);
845         // Each component
846         for (CeedSize comp_out = 0; comp_out < num_comp; comp_out++) {
847           // Each qpoint/node pair
848           if (is_point_block) {
849             // Point Block Diagonal
850             for (CeedInt comp_in = 0; comp_in < num_comp; comp_in++) {
851               CeedScalar e_value = 0.0;
852 
853               for (CeedSize q = 0; q < num_qpts; q++) {
854                 const CeedScalar qf_value =
855                     assembled_qf_array[((((e_in * num_comp + comp_in) * num_e_mode_out + e_out) * num_comp + comp_out) * num_elem + e) * num_qpts +
856                                        q];
857 
858                 e_value += bt[q * num_nodes + tid] * qf_value * b[q * num_nodes + tid];
859               }
860               elem_diag_array[((comp_out * num_comp + comp_in) * num_elem + e) * num_nodes + tid] += e_value;
861             }
862           } else {
863             // Diagonal Only
864             CeedScalar e_value = 0.0;
865 
866             for (CeedSize q = 0; q < num_qpts; q++) {
867               const CeedScalar qf_value =
868                   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];
869               e_value += bt[q * num_nodes + tid] * qf_value * b[q * num_nodes + tid];
870             }
871             elem_diag_array[(comp_out * num_elem + e) * num_nodes + tid] += e_value;
872           }
873         }
874       }
875     }
876   });
877   return CEED_ERROR_SUCCESS;
878 }
879 
880 //------------------------------------------------------------------------------
881 // Assemble diagonal common code
882 //------------------------------------------------------------------------------
883 static inline int CeedOperatorAssembleDiagonalCore_Sycl(CeedOperator op, CeedVector assembled, CeedRequest *request, const bool is_point_block) {
884   Ceed                ceed;
885   Ceed_Sycl          *sycl_data;
886   CeedInt             num_elem;
887   CeedScalar         *elem_diag_array;
888   const CeedScalar   *assembled_qf_array;
889   CeedVector          assembled_qf = NULL;
890   CeedElemRestriction rstr         = NULL;
891   CeedOperator_Sycl  *impl;
892 
893   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
894   CeedCallBackend(CeedOperatorGetData(op, &impl));
895   CeedCallBackend(CeedGetData(ceed, &sycl_data));
896 
897   // Assemble QFunction
898   CeedCallBackend(CeedOperatorLinearAssembleQFunctionBuildOrUpdate(op, &assembled_qf, &rstr, request));
899   CeedCallBackend(CeedElemRestrictionDestroy(&rstr));
900 
901   // Setup
902   if (!impl->diag) {
903     CeedCallBackend(CeedOperatorAssembleDiagonalSetup_Sycl(op));
904   }
905   CeedOperatorDiag_Sycl *diag = impl->diag;
906 
907   assert(diag != NULL);
908 
909   // Restriction
910   if (is_point_block && !diag->point_block_diag_rstr) {
911     CeedCallBackend(CeedOperatorCreateActivePointBlockRestriction(diag->diag_rstr, &diag->point_block_diag_rstr));
912   }
913   CeedElemRestriction diag_rstr = is_point_block ? diag->point_block_diag_rstr : diag->diag_rstr;
914 
915   // Create diagonal vector
916   CeedVector elem_diag = is_point_block ? diag->point_block_elem_diag : diag->elem_diag;
917 
918   if (!elem_diag) {
919     CeedCallBackend(CeedElemRestrictionCreateVector(diag_rstr, NULL, &elem_diag));
920     if (is_point_block) diag->point_block_elem_diag = elem_diag;
921     else diag->elem_diag = elem_diag;
922   }
923   CeedCallBackend(CeedVectorSetValue(elem_diag, 0.0));
924 
925   // Assemble element operator diagonals
926   CeedCallBackend(CeedVectorGetArray(elem_diag, CEED_MEM_DEVICE, &elem_diag_array));
927   CeedCallBackend(CeedVectorGetArrayRead(assembled_qf, CEED_MEM_DEVICE, &assembled_qf_array));
928   CeedCallBackend(CeedElemRestrictionGetNumElements(diag_rstr, &num_elem));
929 
930   // Compute the diagonal of B^T D B
931   CeedCallBackend(CeedOperatorLinearDiagonal_Sycl(sycl_data->sycl_queue, is_point_block, num_elem, diag, assembled_qf_array, elem_diag_array));
932 
933   // Wait for queue to complete and handle exceptions
934   sycl_data->sycl_queue.wait_and_throw();
935 
936   // Restore arrays
937   CeedCallBackend(CeedVectorRestoreArray(elem_diag, &elem_diag_array));
938   CeedCallBackend(CeedVectorRestoreArrayRead(assembled_qf, &assembled_qf_array));
939 
940   // Assemble local operator diagonal
941   CeedCallBackend(CeedElemRestrictionApply(diag_rstr, CEED_TRANSPOSE, elem_diag, assembled, request));
942 
943   // Cleanup
944   CeedCallBackend(CeedVectorDestroy(&assembled_qf));
945   return CEED_ERROR_SUCCESS;
946 }
947 
948 //------------------------------------------------------------------------------
949 // Assemble Linear Diagonal
950 //------------------------------------------------------------------------------
951 static int CeedOperatorLinearAssembleAddDiagonal_Sycl(CeedOperator op, CeedVector assembled, CeedRequest *request) {
952   CeedCallBackend(CeedOperatorAssembleDiagonalCore_Sycl(op, assembled, request, false));
953   return CEED_ERROR_SUCCESS;
954 }
955 
956 //------------------------------------------------------------------------------
957 // Assemble Linear Point Block Diagonal
958 //------------------------------------------------------------------------------
959 static int CeedOperatorLinearAssembleAddPointBlockDiagonal_Sycl(CeedOperator op, CeedVector assembled, CeedRequest *request) {
960   CeedCallBackend(CeedOperatorAssembleDiagonalCore_Sycl(op, assembled, request, true));
961   return CEED_ERROR_SUCCESS;
962 }
963 
964 //------------------------------------------------------------------------------
965 // Single operator assembly setup
966 //------------------------------------------------------------------------------
967 static int CeedSingleOperatorAssembleSetup_Sycl(CeedOperator op) {
968   Ceed    ceed;
969   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,
970                                                num_B_out_mats_to_load = 0, size_B_out = 0, num_qpts = 0, elem_size = 0, num_elem, num_comp,
971                                                mat_start = 0;
972   CeedEvalMode       *eval_mode_in = NULL, *eval_mode_out = NULL;
973   const CeedScalar   *interp_in, *grad_in;
974   CeedElemRestriction rstr_in = NULL, rstr_out = NULL;
975   CeedBasis           basis_in = NULL, basis_out = NULL;
976   CeedQFunctionField *qf_fields;
977   CeedQFunction       qf;
978   CeedOperatorField  *input_fields, *output_fields;
979   CeedOperator_Sycl  *impl;
980 
981   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
982   CeedCallBackend(CeedOperatorGetData(op, &impl));
983 
984   // Get input and output fields
985   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &input_fields, &num_output_fields, &output_fields));
986 
987   // Determine active input basis eval mode
988   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
989   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_fields, NULL, NULL));
990   // Note that the kernel will treat each dimension of a gradient action separately;
991   // i.e., when an active input has a CEED_EVAL_GRAD mode, num_ e_mode_in will increment by dim.
992   // 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,
993   // so num_B_in_mats_to_load will be incremented by 1.
994   for (CeedInt i = 0; i < num_input_fields; i++) {
995     CeedEvalMode eval_mode;
996     CeedVector   vec;
997 
998     CeedCallBackend(CeedOperatorFieldGetVector(input_fields[i], &vec));
999     if (vec == CEED_VECTOR_ACTIVE) {
1000       CeedCallBackend(CeedOperatorFieldGetBasis(input_fields[i], &basis_in));
1001       CeedCallBackend(CeedBasisGetDimension(basis_in, &dim));
1002       CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis_in, &num_qpts));
1003       CeedCallBackend(CeedOperatorFieldGetElemRestriction(input_fields[i], &rstr_in));
1004       CeedCallBackend(CeedElemRestrictionGetElementSize(rstr_in, &elem_size));
1005       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode));
1006       if (eval_mode != CEED_EVAL_NONE) {
1007         CeedCallBackend(CeedRealloc(num_B_in_mats_to_load + 1, &eval_mode_in));
1008         eval_mode_in[num_B_in_mats_to_load] = eval_mode;
1009         num_B_in_mats_to_load += 1;
1010         if (eval_mode == CEED_EVAL_GRAD) {
1011           num_e_mode_in += dim;
1012           size_B_in += dim * elem_size * num_qpts;
1013         } else {
1014           num_e_mode_in += 1;
1015           size_B_in += elem_size * num_qpts;
1016         }
1017       }
1018     }
1019   }
1020 
1021   // Determine active output basis; basis_out and rstr_out only used if same as input, TODO
1022   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qf_fields));
1023   for (CeedInt i = 0; i < num_output_fields; i++) {
1024     CeedEvalMode eval_mode;
1025     CeedVector   vec;
1026 
1027     CeedCallBackend(CeedOperatorFieldGetVector(output_fields[i], &vec));
1028     if (vec == CEED_VECTOR_ACTIVE) {
1029       CeedCallBackend(CeedOperatorFieldGetBasis(output_fields[i], &basis_out));
1030       CeedCallBackend(CeedOperatorFieldGetElemRestriction(output_fields[i], &rstr_out));
1031       CeedCheck(!rstr_out || rstr_out == rstr_in, ceed, CEED_ERROR_BACKEND, "Backend does not implement multi-field non-composite operator assembly");
1032       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode));
1033       if (eval_mode != CEED_EVAL_NONE) {
1034         CeedCallBackend(CeedRealloc(num_B_out_mats_to_load + 1, &eval_mode_out));
1035         eval_mode_out[num_B_out_mats_to_load] = eval_mode;
1036         num_B_out_mats_to_load += 1;
1037         if (eval_mode == CEED_EVAL_GRAD) {
1038           num_e_mode_out += dim;
1039           size_B_out += dim * elem_size * num_qpts;
1040         } else {
1041           num_e_mode_out += 1;
1042           size_B_out += elem_size * num_qpts;
1043         }
1044       }
1045     }
1046   }
1047   CeedCheck(num_e_mode_in > 0 && num_e_mode_out > 0, ceed, CEED_ERROR_UNSUPPORTED, "Cannot assemble operator without inputs/outputs");
1048 
1049   CeedCallBackend(CeedElemRestrictionGetNumElements(rstr_in, &num_elem));
1050   CeedCallBackend(CeedElemRestrictionGetNumComponents(rstr_in, &num_comp));
1051 
1052   CeedCallBackend(CeedCalloc(1, &impl->asmb));
1053   CeedOperatorAssemble_Sycl *asmb = impl->asmb;
1054   asmb->num_elem                  = num_elem;
1055 
1056   Ceed_Sycl *sycl_data;
1057   CeedCallBackend(CeedGetData(ceed, &sycl_data));
1058 
1059   // Kernel setup
1060   int elems_per_block   = 1;
1061   asmb->elems_per_block = elems_per_block;
1062   asmb->block_size_x    = elem_size;
1063   asmb->block_size_y    = elem_size;
1064   asmb->num_e_mode_in   = num_e_mode_in;
1065   asmb->num_e_mode_out  = num_e_mode_out;
1066   asmb->num_qpts        = num_qpts;
1067   asmb->num_nodes       = elem_size;
1068   asmb->block_size      = elem_size * elem_size * elems_per_block;
1069   asmb->num_comp        = num_comp;
1070 
1071   // Build 'full' B matrices (not 1D arrays used for tensor-product matrices
1072   CeedCallBackend(CeedBasisGetInterp(basis_in, &interp_in));
1073   CeedCallBackend(CeedBasisGetGrad(basis_in, &grad_in));
1074 
1075   // Load into B_in, in order that they will be used in eval_mode
1076   CeedCallSycl(ceed, asmb->d_B_in = sycl::malloc_device<CeedScalar>(size_B_in, sycl_data->sycl_device, sycl_data->sycl_context));
1077   for (int i = 0; i < num_B_in_mats_to_load; i++) {
1078     CeedEvalMode eval_mode = eval_mode_in[i];
1079 
1080     if (eval_mode == CEED_EVAL_INTERP) {
1081       // Order queue
1082       sycl::event e = sycl_data->sycl_queue.ext_oneapi_submit_barrier();
1083       sycl_data->sycl_queue.copy<CeedScalar>(interp_in, &asmb->d_B_in[mat_start], elem_size * num_qpts, {e});
1084       mat_start += elem_size * num_qpts;
1085     } else if (eval_mode == CEED_EVAL_GRAD) {
1086       // Order queue
1087       sycl::event e = sycl_data->sycl_queue.ext_oneapi_submit_barrier();
1088       sycl_data->sycl_queue.copy<CeedScalar>(grad_in, &asmb->d_B_in[mat_start], dim * elem_size * num_qpts, {e});
1089       mat_start += dim * elem_size * num_qpts;
1090     }
1091   }
1092 
1093   const CeedScalar *interp_out, *grad_out;
1094   // Note that this function currently assumes 1 basis, so this should always be true
1095   // for now
1096   if (basis_out == basis_in) {
1097     interp_out = interp_in;
1098     grad_out   = grad_in;
1099   } else {
1100     CeedCallBackend(CeedBasisGetInterp(basis_out, &interp_out));
1101     CeedCallBackend(CeedBasisGetGrad(basis_out, &grad_out));
1102   }
1103 
1104   // Load into B_out, in order that they will be used in eval_mode
1105   mat_start = 0;
1106   CeedCallSycl(ceed, asmb->d_B_out = sycl::malloc_device<CeedScalar>(size_B_out, sycl_data->sycl_device, sycl_data->sycl_context));
1107   for (int i = 0; i < num_B_out_mats_to_load; i++) {
1108     CeedEvalMode eval_mode = eval_mode_out[i];
1109 
1110     if (eval_mode == CEED_EVAL_INTERP) {
1111       // Order queue
1112       sycl::event e = sycl_data->sycl_queue.ext_oneapi_submit_barrier();
1113       sycl_data->sycl_queue.copy<CeedScalar>(interp_out, &asmb->d_B_out[mat_start], elem_size * num_qpts, {e});
1114       mat_start += elem_size * num_qpts;
1115     } else if (eval_mode == CEED_EVAL_GRAD) {
1116       // Order queue
1117       sycl::event e = sycl_data->sycl_queue.ext_oneapi_submit_barrier();
1118       sycl_data->sycl_queue.copy<CeedScalar>(grad_out, &asmb->d_B_out[mat_start], dim * elem_size * num_qpts, {e});
1119       mat_start += dim * elem_size * num_qpts;
1120     }
1121   }
1122   return CEED_ERROR_SUCCESS;
1123 }
1124 
1125 //------------------------------------------------------------------------------
1126 // Matrix assembly kernel for low-order elements (3D thread block)
1127 //------------------------------------------------------------------------------
1128 static int CeedOperatorLinearAssemble_Sycl(sycl::queue &sycl_queue, const CeedOperator_Sycl *impl, const CeedScalar *qf_array,
1129                                            CeedScalar *values_array) {
1130   // This kernels assumes B_in and B_out have the same number of quadrature points and basis points.
1131   // TODO: expand to more general cases
1132   CeedOperatorAssemble_Sycl *asmb           = impl->asmb;
1133   const CeedInt              num_elem       = asmb->num_elem;
1134   const CeedSize             num_nodes      = asmb->num_nodes;
1135   const CeedSize             num_comp       = asmb->num_comp;
1136   const CeedSize             num_qpts       = asmb->num_qpts;
1137   const CeedSize             num_e_mode_in  = asmb->num_e_mode_in;
1138   const CeedSize             num_e_mode_out = asmb->num_e_mode_out;
1139 
1140   // Strides for final output ordering, determined by the reference (inference) implementation of the symbolic assembly, slowest --> fastest: element,
1141   // comp_in, comp_out, node_row, node_col
1142   const CeedSize comp_out_stride = num_nodes * num_nodes;
1143   const CeedSize comp_in_stride  = comp_out_stride * num_comp;
1144   const CeedSize e_stride        = comp_in_stride * num_comp;
1145   // Strides for QF array, slowest --> fastest:  e_mode_in, comp_in,  e_mode_out, comp_out, elem, qpt
1146   const CeedSize q_e_stride          = num_qpts;
1147   const CeedSize q_comp_out_stride   = num_elem * q_e_stride;
1148   const CeedSize q_e_mode_out_stride = q_comp_out_stride * num_comp;
1149   const CeedSize q_comp_in_stride    = q_e_mode_out_stride * num_e_mode_out;
1150   const CeedSize q_e_mode_in_stride  = q_comp_in_stride * num_comp;
1151 
1152   CeedScalar *B_in, *B_out;
1153   B_in                       = asmb->d_B_in;
1154   B_out                      = asmb->d_B_out;
1155   const CeedInt block_size_x = asmb->block_size_x;
1156   const CeedInt block_size_y = asmb->block_size_y;
1157 
1158   sycl::range<3> kernel_range(num_elem, block_size_y, block_size_x);
1159 
1160   // Order queue
1161   sycl::event e = sycl_queue.ext_oneapi_submit_barrier();
1162   sycl_queue.parallel_for<CeedOperatorSyclLinearAssemble>(kernel_range, {e}, [=](sycl::id<3> idx) {
1163     const int e = idx.get(0);  // Element index
1164     const int l = idx.get(1);  // The output column index of each B^TDB operation
1165     const int i = idx.get(2);  // The output row index of each B^TDB operation
1166                                // such that we have (Bout^T)_ij D_jk Bin_kl = C_il
1167     for (CeedSize comp_in = 0; comp_in < num_comp; comp_in++) {
1168       for (CeedSize comp_out = 0; comp_out < num_comp; comp_out++) {
1169         CeedScalar result        = 0.0;
1170         CeedSize   qf_index_comp = q_comp_in_stride * comp_in + q_comp_out_stride * comp_out + q_e_stride * e;
1171 
1172         for (CeedSize e_mode_in = 0; e_mode_in < num_e_mode_in; e_mode_in++) {
1173           CeedSize b_in_index = e_mode_in * num_qpts * num_nodes;
1174 
1175           for (CeedSize e_mode_out = 0; e_mode_out < num_e_mode_out; e_mode_out++) {
1176             CeedSize b_out_index = e_mode_out * num_qpts * num_nodes;
1177             CeedSize qf_index    = qf_index_comp + q_e_mode_out_stride * e_mode_out + q_e_mode_in_stride * e_mode_in;
1178 
1179             // Perform the B^T D B operation for this 'chunk' of D (the qf_array)
1180             for (CeedSize j = 0; j < num_qpts; j++) {
1181               result += B_out[b_out_index + j * num_nodes + i] * qf_array[qf_index + j] * B_in[b_in_index + j * num_nodes + l];
1182             }
1183           }  // end of  e_mode_out
1184         }    // end of  e_mode_in
1185         CeedSize val_index = comp_in_stride * comp_in + comp_out_stride * comp_out + e_stride * e + num_nodes * i + l;
1186 
1187         values_array[val_index] = result;
1188       }  // end of out component
1189     }    // end of in component
1190   });
1191   return CEED_ERROR_SUCCESS;
1192 }
1193 
1194 //------------------------------------------------------------------------------
1195 // Fallback kernel for larger orders (1D thread block)
1196 //------------------------------------------------------------------------------
1197 /*
1198 static int CeedOperatorLinearAssembleFallback_Sycl(sycl::queue &sycl_queue, const CeedOperator_Sycl *impl, const CeedScalar *qf_array,
1199                                                    CeedScalar *values_array) {
1200   // This kernel assumes B_in and B_out have the same number of quadrature points and basis points.
1201   // TODO: expand to more general cases
1202   CeedOperatorAssemble_Sycl *asmb        = impl->asmb;
1203   const CeedInt              num_elem       = asmb->num_elem;
1204   const CeedInt              num_nodes      = asmb->num_nodes;
1205   const CeedInt              num_comp       = asmb->num_comp;
1206   const CeedInt              num_qpts       = asmb->num_qpts;
1207   const CeedInt              num_e_mode_in  = asmb->num_e_mode_in;
1208   const CeedInt              num_e_mode_out = asmb->num_e_mode_out;
1209 
1210   // Strides for final output ordering, determined by the reference (interface) implementation of the symbolic assembly, slowest --> fastest: elememt,
1211   // comp_in, comp_out, node_row, node_col
1212   const CeedInt comp_out_stride = num_nodes * num_nodes;
1213   const CeedInt comp_in_stride  = comp_out_stride * num_comp;
1214   const CeedInt e_stride        = comp_in_stride * num_comp;
1215   // Strides for QF array, slowest --> fastest:  e_mode_in, comp_in,  e_mode_out, comp_out, elem, qpt
1216   const CeedInt q_e_stride         = num_qpts;
1217   const CeedInt q_comp_out_stride  = num_elem * q_e_stride;
1218   const CeedInt q_e_mode_out_stride = q_comp_out_stride * num_comp;
1219   const CeedInt q_comp_in_stride   = q_e_mode_out_stride * num_e_mode_out;
1220   const CeedInt q_e_mode_in_stride  = q_comp_in_stride * num_comp;
1221 
1222   CeedScalar *B_in, *B_out;
1223   B_in                        = asmb->d_B_in;
1224   B_out                       = asmb->d_B_out;
1225   const CeedInt elems_per_block = asmb->elems_per_block;
1226   const CeedInt block_size_x  = asmb->block_size_x;
1227   const CeedInt block_size_y  = asmb->block_size_y;  // This will be 1 for the fallback kernel
1228 
1229   const CeedInt     grid = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
1230   sycl::range<3>    local_range(block_size_x, block_size_y, elems_per_block);
1231   sycl::range<3>    global_range(grid * block_size_x, block_size_y, elems_per_block);
1232   sycl::nd_range<3> kernel_range(global_range, local_range);
1233 
1234   sycl_queue.parallel_for<CeedOperatorSyclLinearAssembleFallback>(kernel_range, [=](sycl::nd_item<3> work_item) {
1235     const CeedInt blockIdx  = work_item.get_group(0);
1236     const CeedInt gridDimx  = work_item.get_group_range(0);
1237     const CeedInt threadIdx = work_item.get_local_id(0);
1238     const CeedInt threadIdz = work_item.get_local_id(2);
1239     const CeedInt blockDimz = work_item.get_local_range(2);
1240 
1241     const int l = threadIdx;  // The output column index of each B^TDB operation
1242                               // such that we have (Bout^T)_ij D_jk Bin_kl = C_il
1243     for (CeedInt e = blockIdx * blockDimz + threadIdz; e < num_elem; e += gridDimx * blockDimz) {
1244       for (CeedInt comp_in = 0; comp_in < num_comp; comp_in++) {
1245         for (CeedInt comp_out = 0; comp_out < num_comp; comp_out++) {
1246           for (CeedInt i = 0; i < num_nodes; i++) {
1247             CeedScalar result        = 0.0;
1248             CeedInt    qf_index_comp = q_comp_in_stride * comp_in + q_comp_out_stride * comp_out + q_e_stride * e;
1249             for (CeedInt  e_mode_in = 0;  e_mode_in < num_e_mode_in;  e_mode_in++) {
1250               CeedInt b_in_index =  e_mode_in * num_qpts * num_nodes;
1251               for (CeedInt  e_mode_out = 0;  e_mode_out < num_e_mode_out;  e_mode_out++) {
1252                 CeedInt b_out_index =  e_mode_out * num_qpts * num_nodes;
1253                 CeedInt qf_index    = qf_index_comp + q_e_mode_out_stride *  e_mode_out + q_e_mode_in_stride *  e_mode_in;
1254                 // Perform the B^T D B operation for this 'chunk' of D (the qf_array)
1255                 for (CeedInt j = 0; j < num_qpts; j++) {
1256                   result += B_out[b_out_index + j * num_nodes + i] * qf_array[qf_index + j] * B_in[b_in_index + j * num_nodes + l];
1257                 }
1258               }  // end of  e_mode_out
1259             }    // end of  e_mode_in
1260             CeedInt val_index       = comp_in_stride * comp_in + comp_out_stride * comp_out + e_stride * e + num_nodes * i + l;
1261             values_array[val_index] = result;
1262           }  // end of loop over element node index, i
1263         }    // end of out component
1264       }      // end of in component
1265     }        // end of element loop
1266   });
1267   return CEED_ERROR_SUCCESS;
1268 }*/
1269 
1270 //------------------------------------------------------------------------------
1271 // Assemble matrix data for COO matrix of assembled operator.
1272 // The sparsity pattern is set by CeedOperatorLinearAssembleSymbolic.
1273 //
1274 // Note that this (and other assembly routines) currently assume only one active
1275 // input restriction/basis per operator (could have multiple basis eval modes).
1276 // TODO: allow multiple active input restrictions/basis objects
1277 //------------------------------------------------------------------------------
1278 static int CeedSingleOperatorAssemble_Sycl(CeedOperator op, CeedInt offset, CeedVector values) {
1279   Ceed                ceed;
1280   Ceed_Sycl          *sycl_data;
1281   CeedScalar         *values_array;
1282   const CeedScalar   *qf_array;
1283   CeedVector          assembled_qf = NULL;
1284   CeedElemRestriction rstr_q       = NULL;
1285   CeedOperator_Sycl  *impl;
1286 
1287   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
1288   CeedCallBackend(CeedOperatorGetData(op, &impl));
1289   CeedCallBackend(CeedGetData(ceed, &sycl_data));
1290 
1291   // Setup
1292   if (!impl->asmb) {
1293     CeedCallBackend(CeedSingleOperatorAssembleSetup_Sycl(op));
1294     assert(impl->asmb != NULL);
1295   }
1296 
1297   // Assemble QFunction
1298   CeedCallBackend(CeedOperatorLinearAssembleQFunctionBuildOrUpdate(op, &assembled_qf, &rstr_q, CEED_REQUEST_IMMEDIATE));
1299   CeedCallBackend(CeedElemRestrictionDestroy(&rstr_q));
1300   CeedCallBackend(CeedVectorGetArrayWrite(values, CEED_MEM_DEVICE, &values_array));
1301   values_array += offset;
1302   CeedCallBackend(CeedVectorGetArrayRead(assembled_qf, CEED_MEM_DEVICE, &qf_array));
1303 
1304   // Compute B^T D B
1305   CeedCallBackend(CeedOperatorLinearAssemble_Sycl(sycl_data->sycl_queue, impl, qf_array, values_array));
1306 
1307   // Wait for kernels to be completed
1308   // Kris: Review if this is necessary -- enqueing an async barrier may be sufficient
1309   sycl_data->sycl_queue.wait_and_throw();
1310 
1311   // Restore arrays
1312   CeedCallBackend(CeedVectorRestoreArray(values, &values_array));
1313   CeedCallBackend(CeedVectorRestoreArrayRead(assembled_qf, &qf_array));
1314 
1315   // Cleanup
1316   CeedCallBackend(CeedVectorDestroy(&assembled_qf));
1317   return CEED_ERROR_SUCCESS;
1318 }
1319 
1320 //------------------------------------------------------------------------------
1321 // Create operator
1322 //------------------------------------------------------------------------------
1323 int CeedOperatorCreate_Sycl(CeedOperator op) {
1324   Ceed               ceed;
1325   CeedOperator_Sycl *impl;
1326 
1327   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
1328 
1329   CeedCallBackend(CeedCalloc(1, &impl));
1330   CeedCallBackend(CeedOperatorSetData(op, impl));
1331 
1332   CeedCallBackend(CeedSetBackendFunctionCpp(ceed, "Operator", op, "LinearAssembleQFunction", CeedOperatorLinearAssembleQFunction_Sycl));
1333   CeedCallBackend(CeedSetBackendFunctionCpp(ceed, "Operator", op, "LinearAssembleQFunctionUpdate", CeedOperatorLinearAssembleQFunctionUpdate_Sycl));
1334   CeedCallBackend(CeedSetBackendFunctionCpp(ceed, "Operator", op, "LinearAssembleAddDiagonal", CeedOperatorLinearAssembleAddDiagonal_Sycl));
1335   CeedCallBackend(
1336       CeedSetBackendFunctionCpp(ceed, "Operator", op, "LinearAssembleAddPointBlockDiagonal", CeedOperatorLinearAssembleAddPointBlockDiagonal_Sycl));
1337   CeedCallBackend(CeedSetBackendFunctionCpp(ceed, "Operator", op, "LinearAssembleSingle", CeedSingleOperatorAssemble_Sycl));
1338   CeedCallBackend(CeedSetBackendFunctionCpp(ceed, "Operator", op, "ApplyAdd", CeedOperatorApplyAdd_Sycl));
1339   CeedCallBackend(CeedSetBackendFunctionCpp(ceed, "Operator", op, "Destroy", CeedOperatorDestroy_Sycl));
1340   return CEED_ERROR_SUCCESS;
1341 }
1342 
1343 //------------------------------------------------------------------------------
1344