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 e_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], &e_mode)); 395 if (e_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, e_mode, size 411 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &rstr)); 412 CeedCallBackend(CeedElemRestrictionGetElementSize(rstr, &elem_size)); 413 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &e_mode)); 414 CeedCallBackend(CeedQFunctionFieldGetSize(qf_output_fields[i], &size)); 415 // Basis action 416 switch (e_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 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 431 return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_WEIGHT cannot be an output evaluation mode"); 432 break; // Should not occur 433 case CEED_EVAL_DIV: 434 break; // TODO: Not implemented 435 case CEED_EVAL_CURL: 436 break; // TODO: Not implemented 437 // LCOV_EXCL_STOP 438 } 439 } 440 441 // Output restriction 442 for (CeedInt i = 0; i < num_output_fields; i++) { 443 CeedVector vec; 444 CeedElemRestriction rstr; 445 446 // Restore evec 447 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &e_mode)); 448 if (e_mode == CEED_EVAL_NONE) { 449 CeedCallBackend(CeedVectorRestoreArray(impl->e_vecs[i + impl->num_e_in], &e_data[i + num_input_fields])); 450 } 451 // Get output vector 452 CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[i], &vec)); 453 // Restrict 454 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &rstr)); 455 // Active 456 if (vec == CEED_VECTOR_ACTIVE) vec = out_vec; 457 458 CeedCallBackend(CeedElemRestrictionApply(rstr, CEED_TRANSPOSE, impl->e_vecs[i + impl->num_e_in], vec, request)); 459 } 460 461 // Restore input arrays 462 CeedCallBackend(CeedOperatorRestoreInputs_Sycl(num_input_fields, qf_input_fields, op_input_fields, false, e_data, impl)); 463 return CEED_ERROR_SUCCESS; 464 } 465 466 //------------------------------------------------------------------------------ 467 // Core code for assembling linear QFunction 468 //------------------------------------------------------------------------------ 469 static inline int CeedOperatorLinearAssembleQFunctionCore_Sycl(CeedOperator op, bool build_objects, CeedVector *assembled, CeedElemRestriction *rstr, 470 CeedRequest *request) { 471 Ceed ceed, ceed_parent; 472 CeedSize q_size; 473 CeedInt num_active_in, num_active_out, Q, num_elem, num_input_fields, num_output_fields, size; 474 CeedScalar *assembled_array, *e_data[2 * CEED_FIELD_MAX] = {NULL}; 475 CeedVector *active_in; 476 CeedQFunctionField *qf_input_fields, *qf_output_fields; 477 CeedQFunction qf; 478 CeedOperatorField *op_input_fields, *op_output_fields; 479 CeedOperator_Sycl *impl; 480 481 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 482 CeedCallBackend(CeedOperatorGetFallbackParentCeed(op, &ceed_parent)); 483 CeedCallBackend(CeedOperatorGetData(op, &impl)); 484 CeedCallBackend(CeedOperatorGetQFunction(op, &qf)); 485 CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q)); 486 CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem)); 487 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields)); 488 CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields)); 489 active_in = impl->qf_active_in; 490 num_active_in = impl->num_active_in, num_active_out = impl->num_active_out; 491 492 // Setup 493 CeedCallBackend(CeedOperatorSetup_Sycl(op)); 494 495 // Input Evecs and Restriction 496 CeedCallBackend(CeedOperatorSetupInputs_Sycl(num_input_fields, qf_input_fields, op_input_fields, NULL, true, e_data, impl, request)); 497 498 // Count number of active input fields 499 if (!num_active_in) { 500 for (CeedInt i = 0; i < num_input_fields; i++) { 501 CeedScalar *q_vec_array; 502 CeedVector vec; 503 504 // Get input vector 505 CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec)); 506 // Check if active input 507 if (vec == CEED_VECTOR_ACTIVE) { 508 CeedCallBackend(CeedQFunctionFieldGetSize(qf_input_fields[i], &size)); 509 CeedCallBackend(CeedVectorSetValue(impl->q_vecs_in[i], 0.0)); 510 CeedCallBackend(CeedVectorGetArray(impl->q_vecs_in[i], CEED_MEM_DEVICE, &q_vec_array)); 511 CeedCallBackend(CeedRealloc(num_active_in + size, &active_in)); 512 for (CeedInt field = 0; field < size; field++) { 513 q_size = (CeedSize)Q * num_elem; 514 CeedCallBackend(CeedVectorCreate(ceed, q_size, &active_in[num_active_in + field])); 515 CeedCallBackend( 516 CeedVectorSetArray(active_in[num_active_in + field], CEED_MEM_DEVICE, CEED_USE_POINTER, &q_vec_array[field * Q * num_elem])); 517 } 518 num_active_in += size; 519 CeedCallBackend(CeedVectorRestoreArray(impl->q_vecs_in[i], &q_vec_array)); 520 } 521 } 522 impl->num_active_in = num_active_in; 523 impl->qf_active_in = active_in; 524 } 525 526 // Count number of active output fields 527 if (!num_active_out) { 528 for (CeedInt i = 0; i < num_output_fields; i++) { 529 CeedVector vec; 530 531 // Get output vector 532 CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[i], &vec)); 533 // Check if active output 534 if (vec == CEED_VECTOR_ACTIVE) { 535 CeedCallBackend(CeedQFunctionFieldGetSize(qf_output_fields[i], &size)); 536 num_active_out += size; 537 } 538 } 539 impl->num_active_out = num_active_out; 540 } 541 542 // Check sizes 543 if (!num_active_in || !num_active_out) { 544 // LCOV_EXCL_START 545 return CeedError(ceed, CEED_ERROR_BACKEND, "Cannot assemble QFunction without active inputs and outputs"); 546 // LCOV_EXCL_STOP 547 } 548 549 // Build objects if needed 550 if (build_objects) { 551 CeedSize l_size = (CeedSize)num_elem * Q * num_active_in * num_active_out; 552 CeedInt strides[3] = {1, num_elem * Q, Q}; /* *NOPAD* */ 553 554 // Create output restriction 555 CeedCallBackend(CeedElemRestrictionCreateStrided(ceed_parent, num_elem, Q, num_active_in * num_active_out, 556 num_active_in * num_active_out * num_elem * Q, 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 if (rstr_in && rstr_in != rstr) { 659 // LCOV_EXCL_START 660 return CeedError(ceed, CEED_ERROR_BACKEND, "Backend does not implement multi-field non-composite operator diagonal assembly"); 661 // LCOV_EXCL_STOP 662 } 663 rstr_in = rstr; 664 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &e_mode)); 665 switch (e_mode) { 666 case CEED_EVAL_NONE: 667 case CEED_EVAL_INTERP: 668 CeedCallBackend(CeedRealloc(num_e_mode_in + 1, &e_mode_in)); 669 e_mode_in[num_e_mode_in] = e_mode; 670 num_e_mode_in += 1; 671 break; 672 case CEED_EVAL_GRAD: 673 CeedCallBackend(CeedRealloc(num_e_mode_in + dim, &e_mode_in)); 674 for (CeedInt d = 0; d < dim; d++) e_mode_in[num_e_mode_in + d] = e_mode; 675 num_e_mode_in += dim; 676 break; 677 case CEED_EVAL_WEIGHT: 678 case CEED_EVAL_DIV: 679 case CEED_EVAL_CURL: 680 break; // Caught by QF Assembly 681 } 682 } 683 } 684 685 // Determine active output basis 686 CeedCallBackend(CeedOperatorGetFields(op, NULL, NULL, NULL, &op_fields)); 687 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qf_fields)); 688 for (CeedInt i = 0; i < num_output_fields; i++) { 689 CeedVector vec; 690 691 CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec)); 692 if (vec == CEED_VECTOR_ACTIVE) { 693 CeedEvalMode e_mode; 694 CeedElemRestriction rstr; 695 696 CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis_out)); 697 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_fields[i], &rstr)); 698 if (rstr_out && rstr_out != rstr) { 699 // LCOV_EXCL_START 700 return CeedError(ceed, CEED_ERROR_BACKEND, "Backend does not implement multi-field non-composite operator diagonal assembly"); 701 // LCOV_EXCL_STOP 702 } 703 rstr_out = rstr; 704 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &e_mode)); 705 switch (e_mode) { 706 case CEED_EVAL_NONE: 707 case CEED_EVAL_INTERP: 708 CeedCallBackend(CeedRealloc(num_e_mode_out + 1, &e_mode_out)); 709 e_mode_out[num_e_mode_out] = e_mode; 710 num_e_mode_out += 1; 711 break; 712 case CEED_EVAL_GRAD: 713 CeedCallBackend(CeedRealloc(num_e_mode_out + dim, &e_mode_out)); 714 for (CeedInt d = 0; d < dim; d++) e_mode_out[num_e_mode_out + d] = e_mode; 715 num_e_mode_out += dim; 716 break; 717 case CEED_EVAL_WEIGHT: 718 case CEED_EVAL_DIV: 719 case CEED_EVAL_CURL: 720 break; // Caught by QF Assembly 721 } 722 } 723 } 724 725 // Operator data struct 726 CeedCallBackend(CeedOperatorGetData(op, &impl)); 727 CeedCallBackend(CeedGetData(ceed, &sycl_data)); 728 CeedCallBackend(CeedCalloc(1, &impl->diag)); 729 CeedOperatorDiag_Sycl *diag = impl->diag; 730 731 diag->basis_in = basis_in; 732 diag->basis_out = basis_out; 733 diag->h_e_mode_in = e_mode_in; 734 diag->h_e_mode_out = e_mode_out; 735 diag->num_e_mode_in = num_e_mode_in; 736 diag->num_e_mode_out = num_e_mode_out; 737 738 // Kernel parameters 739 CeedInt num_nodes, num_qpts; 740 CeedCallBackend(CeedBasisGetNumNodes(basis_in, &num_nodes)); 741 CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis_in, &num_qpts)); 742 diag->num_nodes = num_nodes; 743 diag->num_qpts = num_qpts; 744 diag->num_comp = num_comp; 745 746 // Basis matrices 747 const CeedInt i_len = num_qpts * num_nodes; 748 const CeedInt g_len = num_qpts * num_nodes * dim; 749 const CeedScalar *interp_in, *interp_out, *grad_in, *grad_out; 750 751 // CEED_EVAL_NONE 752 CeedScalar *identity = NULL; 753 bool has_eval_none = false; 754 for (CeedInt i = 0; i < num_e_mode_in; i++) has_eval_none = has_eval_none || (e_mode_in[i] == CEED_EVAL_NONE); 755 for (CeedInt i = 0; i < num_e_mode_out; i++) has_eval_none = has_eval_none || (e_mode_out[i] == CEED_EVAL_NONE); 756 757 // Order queue 758 sycl::event e = sycl_data->sycl_queue.ext_oneapi_submit_barrier(); 759 760 std::vector<sycl::event> copy_events; 761 if (has_eval_none) { 762 CeedCallBackend(CeedCalloc(num_qpts * num_nodes, &identity)); 763 for (CeedSize i = 0; i < (num_nodes < num_qpts ? num_nodes : num_qpts); i++) identity[i * num_nodes + i] = 1.0; 764 CeedCallSycl(ceed, diag->d_identity = sycl::malloc_device<CeedScalar>(i_len, sycl_data->sycl_device, sycl_data->sycl_context)); 765 sycl::event identity_copy = sycl_data->sycl_queue.copy<CeedScalar>(identity, diag->d_identity, i_len, {e}); 766 copy_events.push_back(identity_copy); 767 } 768 769 // CEED_EVAL_INTERP 770 CeedCallBackend(CeedBasisGetInterp(basis_in, &interp_in)); 771 CeedCallSycl(ceed, diag->d_interp_in = sycl::malloc_device<CeedScalar>(i_len, sycl_data->sycl_device, sycl_data->sycl_context)); 772 sycl::event interp_in_copy = sycl_data->sycl_queue.copy<CeedScalar>(interp_in, diag->d_interp_in, i_len, {e}); 773 copy_events.push_back(interp_in_copy); 774 775 CeedCallBackend(CeedBasisGetInterp(basis_out, &interp_out)); 776 CeedCallSycl(ceed, diag->d_interp_out = sycl::malloc_device<CeedScalar>(i_len, sycl_data->sycl_device, sycl_data->sycl_context)); 777 sycl::event interp_out_copy = sycl_data->sycl_queue.copy<CeedScalar>(interp_out, diag->d_interp_out, i_len, {e}); 778 copy_events.push_back(interp_out_copy); 779 780 // CEED_EVAL_GRAD 781 CeedCallBackend(CeedBasisGetGrad(basis_in, &grad_in)); 782 CeedCallSycl(ceed, diag->d_grad_in = sycl::malloc_device<CeedScalar>(g_len, sycl_data->sycl_device, sycl_data->sycl_context)); 783 sycl::event grad_in_copy = sycl_data->sycl_queue.copy<CeedScalar>(grad_in, diag->d_grad_in, g_len, {e}); 784 copy_events.push_back(grad_in_copy); 785 786 CeedCallBackend(CeedBasisGetGrad(basis_out, &grad_out)); 787 CeedCallSycl(ceed, diag->d_grad_out = sycl::malloc_device<CeedScalar>(g_len, sycl_data->sycl_device, sycl_data->sycl_context)); 788 sycl::event grad_out_copy = sycl_data->sycl_queue.copy<CeedScalar>(grad_out, diag->d_grad_out, g_len, {e}); 789 copy_events.push_back(grad_out_copy); 790 791 // Arrays of e_modes 792 CeedCallSycl(ceed, diag->d_e_mode_in = sycl::malloc_device<CeedEvalMode>(num_e_mode_in, sycl_data->sycl_device, sycl_data->sycl_context)); 793 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}); 794 copy_events.push_back(e_mode_in_copy); 795 796 CeedCallSycl(ceed, diag->d_e_mode_out = sycl::malloc_device<CeedEvalMode>(num_e_mode_out, sycl_data->sycl_device, sycl_data->sycl_context)); 797 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}); 798 copy_events.push_back(e_mode_out_copy); 799 800 // Restriction 801 diag->diag_rstr = rstr_out; 802 803 // Wait for all copies to complete and handle exceptions 804 CeedCallSycl(ceed, sycl::event::wait_and_throw(copy_events)); 805 return CEED_ERROR_SUCCESS; 806 } 807 808 //------------------------------------------------------------------------------ 809 // Kernel for diagonal assembly 810 //------------------------------------------------------------------------------ 811 static int CeedOperatorLinearDiagonal_Sycl(sycl::queue &sycl_queue, const bool is_point_block, const CeedInt num_elem, 812 const CeedOperatorDiag_Sycl *diag, const CeedScalar *assembled_qf_array, CeedScalar *elem_diag_array) { 813 const CeedSize num_nodes = diag->num_nodes; 814 const CeedSize num_qpts = diag->num_qpts; 815 const CeedSize num_comp = diag->num_comp; 816 const CeedSize num_e_mode_in = diag->num_e_mode_in; 817 const CeedSize num_e_mode_out = diag->num_e_mode_out; 818 const CeedScalar *identity = diag->d_identity; 819 const CeedScalar *interp_in = diag->d_interp_in; 820 const CeedScalar *grad_in = diag->d_grad_in; 821 const CeedScalar *interp_out = diag->d_interp_out; 822 const CeedScalar *grad_out = diag->d_grad_out; 823 const CeedEvalMode *e_mode_in = diag->d_e_mode_in; 824 const CeedEvalMode *e_mode_out = diag->d_e_mode_out; 825 826 sycl::range<1> kernel_range(num_elem * num_nodes); 827 828 // Order queue 829 sycl::event e = sycl_queue.ext_oneapi_submit_barrier(); 830 sycl_queue.parallel_for<CeedOperatorSyclLinearDiagonal>(kernel_range, {e}, [=](sycl::id<1> idx) { 831 const CeedInt tid = idx % num_nodes; 832 const CeedInt e = idx / num_nodes; 833 834 // Compute the diagonal of B^T D B 835 // Each element 836 CeedInt d_out = -1; 837 // Each basis eval mode pair 838 for (CeedSize e_out = 0; e_out < num_e_mode_out; e_out++) { 839 const CeedScalar *bt = NULL; 840 841 if (e_mode_out[e_out] == CEED_EVAL_GRAD) ++d_out; 842 CeedOperatorGetBasisPointer_Sycl(&bt, e_mode_out[e_out], identity, interp_out, &grad_out[d_out * num_qpts * num_nodes]); 843 CeedInt d_in = -1; 844 845 for (CeedSize e_in = 0; e_in < num_e_mode_in; e_in++) { 846 const CeedScalar *b = NULL; 847 848 if (e_mode_in[e_in] == CEED_EVAL_GRAD) ++d_in; 849 CeedOperatorGetBasisPointer_Sycl(&b, e_mode_in[e_in], identity, interp_in, &grad_in[d_in * num_qpts * num_nodes]); 850 // Each component 851 for (CeedSize comp_out = 0; comp_out < num_comp; comp_out++) { 852 // Each qpoint/node pair 853 if (is_point_block) { 854 // Point Block Diagonal 855 for (CeedInt comp_in = 0; comp_in < num_comp; comp_in++) { 856 CeedScalar e_value = 0.0; 857 858 for (CeedSize q = 0; q < num_qpts; q++) { 859 const CeedScalar qf_value = 860 assembled_qf_array[((((e_in * num_comp + comp_in) * num_e_mode_out + e_out) * num_comp + comp_out) * num_elem + e) * num_qpts + 861 q]; 862 863 e_value += bt[q * num_nodes + tid] * qf_value * b[q * num_nodes + tid]; 864 } 865 elem_diag_array[((comp_out * num_comp + comp_in) * num_elem + e) * num_nodes + tid] += e_value; 866 } 867 } else { 868 // Diagonal Only 869 CeedScalar e_value = 0.0; 870 871 for (CeedSize q = 0; q < num_qpts; q++) { 872 const CeedScalar qf_value = 873 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]; 874 e_value += bt[q * num_nodes + tid] * qf_value * b[q * num_nodes + tid]; 875 } 876 elem_diag_array[(comp_out * num_elem + e) * num_nodes + tid] += e_value; 877 } 878 } 879 } 880 } 881 }); 882 return CEED_ERROR_SUCCESS; 883 } 884 885 //------------------------------------------------------------------------------ 886 // Assemble diagonal common code 887 //------------------------------------------------------------------------------ 888 static inline int CeedOperatorAssembleDiagonalCore_Sycl(CeedOperator op, CeedVector assembled, CeedRequest *request, const bool is_point_block) { 889 Ceed ceed; 890 Ceed_Sycl *sycl_data; 891 CeedInt num_elem; 892 CeedScalar *elem_diag_array; 893 const CeedScalar *assembled_qf_array; 894 CeedVector assembled_qf = NULL; 895 CeedElemRestriction rstr = NULL; 896 CeedOperator_Sycl *impl; 897 898 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 899 CeedCallBackend(CeedOperatorGetData(op, &impl)); 900 CeedCallBackend(CeedGetData(ceed, &sycl_data)); 901 902 // Assemble QFunction 903 CeedCallBackend(CeedOperatorLinearAssembleQFunctionBuildOrUpdate(op, &assembled_qf, &rstr, request)); 904 CeedCallBackend(CeedElemRestrictionDestroy(&rstr)); 905 906 // Setup 907 if (!impl->diag) { 908 CeedCallBackend(CeedOperatorAssembleDiagonalSetup_Sycl(op)); 909 } 910 CeedOperatorDiag_Sycl *diag = impl->diag; 911 912 assert(diag != NULL); 913 914 // Restriction 915 if (is_point_block && !diag->point_block_diag_rstr) { 916 CeedCallBackend(CeedOperatorCreateActivePointBlockRestriction(diag->diag_rstr, &diag->point_block_diag_rstr)); 917 } 918 CeedElemRestriction diag_rstr = is_point_block ? diag->point_block_diag_rstr : diag->diag_rstr; 919 920 // Create diagonal vector 921 CeedVector elem_diag = is_point_block ? diag->point_block_elem_diag : diag->elem_diag; 922 923 if (!elem_diag) { 924 CeedCallBackend(CeedElemRestrictionCreateVector(diag_rstr, NULL, &elem_diag)); 925 if (is_point_block) diag->point_block_elem_diag = elem_diag; 926 else diag->elem_diag = elem_diag; 927 } 928 CeedCallBackend(CeedVectorSetValue(elem_diag, 0.0)); 929 930 // Assemble element operator diagonals 931 CeedCallBackend(CeedVectorGetArray(elem_diag, CEED_MEM_DEVICE, &elem_diag_array)); 932 CeedCallBackend(CeedVectorGetArrayRead(assembled_qf, CEED_MEM_DEVICE, &assembled_qf_array)); 933 CeedCallBackend(CeedElemRestrictionGetNumElements(diag_rstr, &num_elem)); 934 935 // Compute the diagonal of B^T D B 936 CeedCallBackend(CeedOperatorLinearDiagonal_Sycl(sycl_data->sycl_queue, is_point_block, num_elem, diag, assembled_qf_array, elem_diag_array)); 937 938 // Wait for queue to complete and handle exceptions 939 sycl_data->sycl_queue.wait_and_throw(); 940 941 // Restore arrays 942 CeedCallBackend(CeedVectorRestoreArray(elem_diag, &elem_diag_array)); 943 CeedCallBackend(CeedVectorRestoreArrayRead(assembled_qf, &assembled_qf_array)); 944 945 // Assemble local operator diagonal 946 CeedCallBackend(CeedElemRestrictionApply(diag_rstr, CEED_TRANSPOSE, elem_diag, assembled, request)); 947 948 // Cleanup 949 CeedCallBackend(CeedVectorDestroy(&assembled_qf)); 950 return CEED_ERROR_SUCCESS; 951 } 952 953 //------------------------------------------------------------------------------ 954 // Assemble Linear Diagonal 955 //------------------------------------------------------------------------------ 956 static int CeedOperatorLinearAssembleAddDiagonal_Sycl(CeedOperator op, CeedVector assembled, CeedRequest *request) { 957 CeedCallBackend(CeedOperatorAssembleDiagonalCore_Sycl(op, assembled, request, false)); 958 return CEED_ERROR_SUCCESS; 959 } 960 961 //------------------------------------------------------------------------------ 962 // Assemble Linear Point Block Diagonal 963 //------------------------------------------------------------------------------ 964 static int CeedOperatorLinearAssembleAddPointBlockDiagonal_Sycl(CeedOperator op, CeedVector assembled, CeedRequest *request) { 965 CeedCallBackend(CeedOperatorAssembleDiagonalCore_Sycl(op, assembled, request, true)); 966 return CEED_ERROR_SUCCESS; 967 } 968 969 //------------------------------------------------------------------------------ 970 // Single operator assembly setup 971 //------------------------------------------------------------------------------ 972 static int CeedSingleOperatorAssembleSetup_Sycl(CeedOperator op) { 973 Ceed ceed; 974 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, 975 num_B_out_mats_to_load = 0, size_B_out = 0, num_qpts = 0, elem_size = 0, num_elem, num_comp, 976 mat_start = 0; 977 CeedEvalMode *eval_mode_in = NULL, *eval_mode_out = NULL; 978 const CeedScalar *interp_in, *grad_in; 979 CeedElemRestriction rstr_in = NULL, rstr_out = NULL; 980 CeedBasis basis_in = NULL, basis_out = NULL; 981 CeedQFunctionField *qf_fields; 982 CeedQFunction qf; 983 CeedOperatorField *input_fields, *output_fields; 984 CeedOperator_Sycl *impl; 985 986 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 987 CeedCallBackend(CeedOperatorGetData(op, &impl)); 988 989 // Get input and output fields 990 CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &input_fields, &num_output_fields, &output_fields)); 991 992 // Determine active input basis eval mode 993 CeedCallBackend(CeedOperatorGetQFunction(op, &qf)); 994 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_fields, NULL, NULL)); 995 // Note that the kernel will treat each dimension of a gradient action separately; 996 // i.e., when an active input has a CEED_EVAL_GRAD mode, num_ e_mode_in will increment by dim. 997 // 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, 998 // so num_B_in_mats_to_load will be incremented by 1. 999 for (CeedInt i = 0; i < num_input_fields; i++) { 1000 CeedEvalMode eval_mode; 1001 CeedVector vec; 1002 1003 CeedCallBackend(CeedOperatorFieldGetVector(input_fields[i], &vec)); 1004 if (vec == CEED_VECTOR_ACTIVE) { 1005 CeedCallBackend(CeedOperatorFieldGetBasis(input_fields[i], &basis_in)); 1006 CeedCallBackend(CeedBasisGetDimension(basis_in, &dim)); 1007 CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis_in, &num_qpts)); 1008 CeedCallBackend(CeedOperatorFieldGetElemRestriction(input_fields[i], &rstr_in)); 1009 CeedCallBackend(CeedElemRestrictionGetElementSize(rstr_in, &elem_size)); 1010 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode)); 1011 if (eval_mode != CEED_EVAL_NONE) { 1012 CeedCallBackend(CeedRealloc(num_B_in_mats_to_load + 1, &eval_mode_in)); 1013 eval_mode_in[num_B_in_mats_to_load] = eval_mode; 1014 num_B_in_mats_to_load += 1; 1015 if (eval_mode == CEED_EVAL_GRAD) { 1016 num_e_mode_in += dim; 1017 size_B_in += dim * elem_size * num_qpts; 1018 } else { 1019 num_e_mode_in += 1; 1020 size_B_in += elem_size * num_qpts; 1021 } 1022 } 1023 } 1024 } 1025 1026 // Determine active output basis; basis_out and rstr_out only used if same as input, TODO 1027 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qf_fields)); 1028 for (CeedInt i = 0; i < num_output_fields; i++) { 1029 CeedEvalMode eval_mode; 1030 CeedVector vec; 1031 1032 CeedCallBackend(CeedOperatorFieldGetVector(output_fields[i], &vec)); 1033 if (vec == CEED_VECTOR_ACTIVE) { 1034 CeedCallBackend(CeedOperatorFieldGetBasis(output_fields[i], &basis_out)); 1035 CeedCallBackend(CeedOperatorFieldGetElemRestriction(output_fields[i], &rstr_out)); 1036 if (rstr_out && rstr_out != rstr_in) { 1037 // LCOV_EXCL_START 1038 return CeedError(ceed, CEED_ERROR_BACKEND, "Backend does not implement multi-field non-composite operator assembly"); 1039 // LCOV_EXCL_STOP 1040 } 1041 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode)); 1042 if (eval_mode != CEED_EVAL_NONE) { 1043 CeedCallBackend(CeedRealloc(num_B_out_mats_to_load + 1, &eval_mode_out)); 1044 eval_mode_out[num_B_out_mats_to_load] = eval_mode; 1045 num_B_out_mats_to_load += 1; 1046 if (eval_mode == CEED_EVAL_GRAD) { 1047 num_e_mode_out += dim; 1048 size_B_out += dim * elem_size * num_qpts; 1049 } else { 1050 num_e_mode_out += 1; 1051 size_B_out += elem_size * num_qpts; 1052 } 1053 } 1054 } 1055 } 1056 1057 if (num_e_mode_in == 0 || num_e_mode_out == 0) { 1058 // LCOV_EXCL_START 1059 return CeedError(ceed, CEED_ERROR_UNSUPPORTED, "Cannot assemble operator without inputs/outputs"); 1060 // LCOV_EXCL_STOP 1061 } 1062 1063 CeedCallBackend(CeedElemRestrictionGetNumElements(rstr_in, &num_elem)); 1064 CeedCallBackend(CeedElemRestrictionGetNumComponents(rstr_in, &num_comp)); 1065 1066 CeedCallBackend(CeedCalloc(1, &impl->asmb)); 1067 CeedOperatorAssemble_Sycl *asmb = impl->asmb; 1068 asmb->num_elem = num_elem; 1069 1070 Ceed_Sycl *sycl_data; 1071 CeedCallBackend(CeedGetData(ceed, &sycl_data)); 1072 1073 // Kernel setup 1074 int elem_per_block = 1; 1075 asmb->elem_per_block = elem_per_block; 1076 CeedInt block_size = elem_size * elem_size * elem_per_block; 1077 1078 /* CeedInt maxThreadsPerBlock = sycl_data->sycl_device.get_info<sycl::info::device::max_work_group_size>(); 1079 bool fallback = block_size > maxThreadsPerBlock; 1080 asmb->fallback = fallback; 1081 if (fallback) { 1082 // Use fallback kernel with 1D threadblock 1083 block_size = elem_size * elem_per_block; 1084 asmb->block_size_x = elem_size; 1085 asmb->block_size_y = 1; 1086 } else { // Use kernel with 2D threadblock 1087 asmb->block_size_x = elem_size; 1088 asmb->block_size_y = elem_size; 1089 }*/ 1090 asmb->block_size_x = elem_size; 1091 asmb->block_size_y = elem_size; 1092 asmb->num_e_mode_in = num_e_mode_in; 1093 asmb->num_e_mode_out = num_e_mode_out; 1094 asmb->num_qpts = num_qpts; 1095 asmb->num_nodes = elem_size; 1096 asmb->block_size = block_size; 1097 asmb->num_comp = num_comp; 1098 1099 // Build 'full' B matrices (not 1D arrays used for tensor-product matrices 1100 CeedCallBackend(CeedBasisGetInterp(basis_in, &interp_in)); 1101 CeedCallBackend(CeedBasisGetGrad(basis_in, &grad_in)); 1102 1103 // Load into B_in, in order that they will be used in eval_mode 1104 CeedCallSycl(ceed, asmb->d_B_in = sycl::malloc_device<CeedScalar>(size_B_in, sycl_data->sycl_device, sycl_data->sycl_context)); 1105 for (int i = 0; i < num_B_in_mats_to_load; i++) { 1106 CeedEvalMode eval_mode = eval_mode_in[i]; 1107 1108 if (eval_mode == CEED_EVAL_INTERP) { 1109 // Order queue 1110 sycl::event e = sycl_data->sycl_queue.ext_oneapi_submit_barrier(); 1111 sycl_data->sycl_queue.copy<CeedScalar>(interp_in, &asmb->d_B_in[mat_start], elem_size * num_qpts, {e}); 1112 mat_start += elem_size * num_qpts; 1113 } else if (eval_mode == CEED_EVAL_GRAD) { 1114 // Order queue 1115 sycl::event e = sycl_data->sycl_queue.ext_oneapi_submit_barrier(); 1116 sycl_data->sycl_queue.copy<CeedScalar>(grad_in, &asmb->d_B_in[mat_start], dim * elem_size * num_qpts, {e}); 1117 mat_start += dim * elem_size * num_qpts; 1118 } 1119 } 1120 1121 const CeedScalar *interp_out, *grad_out; 1122 // Note that this function currently assumes 1 basis, so this should always be true 1123 // for now 1124 if (basis_out == basis_in) { 1125 interp_out = interp_in; 1126 grad_out = grad_in; 1127 } else { 1128 CeedCallBackend(CeedBasisGetInterp(basis_out, &interp_out)); 1129 CeedCallBackend(CeedBasisGetGrad(basis_out, &grad_out)); 1130 } 1131 1132 // Load into B_out, in order that they will be used in eval_mode 1133 mat_start = 0; 1134 CeedCallSycl(ceed, asmb->d_B_out = sycl::malloc_device<CeedScalar>(size_B_out, sycl_data->sycl_device, sycl_data->sycl_context)); 1135 for (int i = 0; i < num_B_out_mats_to_load; i++) { 1136 CeedEvalMode eval_mode = eval_mode_out[i]; 1137 1138 if (eval_mode == CEED_EVAL_INTERP) { 1139 // Order queue 1140 sycl::event e = sycl_data->sycl_queue.ext_oneapi_submit_barrier(); 1141 sycl_data->sycl_queue.copy<CeedScalar>(interp_out, &asmb->d_B_out[mat_start], elem_size * num_qpts, {e}); 1142 mat_start += elem_size * num_qpts; 1143 } else if (eval_mode == CEED_EVAL_GRAD) { 1144 // Order queue 1145 sycl::event e = sycl_data->sycl_queue.ext_oneapi_submit_barrier(); 1146 sycl_data->sycl_queue.copy<CeedScalar>(grad_out, &asmb->d_B_out[mat_start], dim * elem_size * num_qpts, {e}); 1147 mat_start += dim * elem_size * num_qpts; 1148 } 1149 } 1150 return CEED_ERROR_SUCCESS; 1151 } 1152 1153 //------------------------------------------------------------------------------ 1154 // Matrix assembly kernel for low-order elements (3D thread block) 1155 //------------------------------------------------------------------------------ 1156 static int CeedOperatorLinearAssemble_Sycl(sycl::queue &sycl_queue, const CeedOperator_Sycl *impl, const CeedScalar *qf_array, 1157 CeedScalar *values_array) { 1158 // This kernels assumes B_in and B_out have the same number of quadrature points and basis points. 1159 // TODO: expand to more general cases 1160 CeedOperatorAssemble_Sycl *asmb = impl->asmb; 1161 const CeedInt num_elem = asmb->num_elem; 1162 const CeedSize num_nodes = asmb->num_nodes; 1163 const CeedSize num_comp = asmb->num_comp; 1164 const CeedSize num_qpts = asmb->num_qpts; 1165 const CeedSize num_e_mode_in = asmb->num_e_mode_in; 1166 const CeedSize num_e_mode_out = asmb->num_e_mode_out; 1167 1168 // Strides for final output ordering, determined by the reference (inference) implementation of the symbolic assembly, slowest --> fastest: element, 1169 // comp_in, comp_out, node_row, node_col 1170 const CeedSize comp_out_stride = num_nodes * num_nodes; 1171 const CeedSize comp_in_stride = comp_out_stride * num_comp; 1172 const CeedSize e_stride = comp_in_stride * num_comp; 1173 // Strides for QF array, slowest --> fastest: e_mode_in, comp_in, e_mode_out, comp_out, elem, qpt 1174 const CeedSize q_e_stride = num_qpts; 1175 const CeedSize q_comp_out_stride = num_elem * q_e_stride; 1176 const CeedSize q_e_mode_out_stride = q_comp_out_stride * num_comp; 1177 const CeedSize q_comp_in_stride = q_e_mode_out_stride * num_e_mode_out; 1178 const CeedSize q_e_mode_in_stride = q_comp_in_stride * num_comp; 1179 1180 CeedScalar *B_in, *B_out; 1181 B_in = asmb->d_B_in; 1182 B_out = asmb->d_B_out; 1183 const CeedInt block_size_x = asmb->block_size_x; 1184 const CeedInt block_size_y = asmb->block_size_y; 1185 1186 sycl::range<3> kernel_range(num_elem, block_size_y, block_size_x); 1187 1188 // Order queue 1189 sycl::event e = sycl_queue.ext_oneapi_submit_barrier(); 1190 sycl_queue.parallel_for<CeedOperatorSyclLinearAssemble>(kernel_range, {e}, [=](sycl::id<3> idx) { 1191 const int e = idx.get(0); // Element index 1192 const int l = idx.get(1); // The output column index of each B^TDB operation 1193 const int i = idx.get(2); // The output row index of each B^TDB operation 1194 // such that we have (Bout^T)_ij D_jk Bin_kl = C_il 1195 for (CeedSize comp_in = 0; comp_in < num_comp; comp_in++) { 1196 for (CeedSize comp_out = 0; comp_out < num_comp; comp_out++) { 1197 CeedScalar result = 0.0; 1198 CeedSize qf_index_comp = q_comp_in_stride * comp_in + q_comp_out_stride * comp_out + q_e_stride * e; 1199 1200 for (CeedSize e_mode_in = 0; e_mode_in < num_e_mode_in; e_mode_in++) { 1201 CeedSize b_in_index = e_mode_in * num_qpts * num_nodes; 1202 1203 for (CeedSize e_mode_out = 0; e_mode_out < num_e_mode_out; e_mode_out++) { 1204 CeedSize b_out_index = e_mode_out * num_qpts * num_nodes; 1205 CeedSize qf_index = qf_index_comp + q_e_mode_out_stride * e_mode_out + q_e_mode_in_stride * e_mode_in; 1206 1207 // Perform the B^T D B operation for this 'chunk' of D (the qf_array) 1208 for (CeedSize j = 0; j < num_qpts; j++) { 1209 result += B_out[b_out_index + j * num_nodes + i] * qf_array[qf_index + j] * B_in[b_in_index + j * num_nodes + l]; 1210 } 1211 } // end of e_mode_out 1212 } // end of e_mode_in 1213 CeedSize val_index = comp_in_stride * comp_in + comp_out_stride * comp_out + e_stride * e + num_nodes * i + l; 1214 1215 values_array[val_index] = result; 1216 } // end of out component 1217 } // end of in component 1218 }); 1219 return CEED_ERROR_SUCCESS; 1220 } 1221 1222 //------------------------------------------------------------------------------ 1223 // Fallback kernel for larger orders (1D thread block) 1224 //------------------------------------------------------------------------------ 1225 /* 1226 static int CeedOperatorLinearAssembleFallback_Sycl(sycl::queue &sycl_queue, const CeedOperator_Sycl *impl, const CeedScalar *qf_array, 1227 CeedScalar *values_array) { 1228 // This kernel assumes B_in and B_out have the same number of quadrature points and basis points. 1229 // TODO: expand to more general cases 1230 CeedOperatorAssemble_Sycl *asmb = impl->asmb; 1231 const CeedInt num_elem = asmb->num_elem; 1232 const CeedInt num_nodes = asmb->num_nodes; 1233 const CeedInt num_comp = asmb->num_comp; 1234 const CeedInt num_qpts = asmb->num_qpts; 1235 const CeedInt num_e_mode_in = asmb->num_e_mode_in; 1236 const CeedInt num_e_mode_out = asmb->num_e_mode_out; 1237 1238 // Strides for final output ordering, determined by the reference (interface) implementation of the symbolic assembly, slowest --> fastest: elememt, 1239 // comp_in, comp_out, node_row, node_col 1240 const CeedInt comp_out_stride = num_nodes * num_nodes; 1241 const CeedInt comp_in_stride = comp_out_stride * num_comp; 1242 const CeedInt e_stride = comp_in_stride * num_comp; 1243 // Strides for QF array, slowest --> fastest: e_mode_in, comp_in, e_mode_out, comp_out, elem, qpt 1244 const CeedInt q_e_stride = num_qpts; 1245 const CeedInt q_comp_out_stride = num_elem * q_e_stride; 1246 const CeedInt q_e_mode_out_stride = q_comp_out_stride * num_comp; 1247 const CeedInt q_comp_in_stride = q_e_mode_out_stride * num_e_mode_out; 1248 const CeedInt q_e_mode_in_stride = q_comp_in_stride * num_comp; 1249 1250 CeedScalar *B_in, *B_out; 1251 B_in = asmb->d_B_in; 1252 B_out = asmb->d_B_out; 1253 const CeedInt elem_per_block = asmb->elem_per_block; 1254 const CeedInt block_size_x = asmb->block_size_x; 1255 const CeedInt block_size_y = asmb->block_size_y; // This will be 1 for the fallback kernel 1256 1257 const CeedInt grid = num_elem / elem_per_block + ((num_elem / elem_per_block * elem_per_block < num_elem) ? 1 : 0); 1258 sycl::range<3> local_range(block_size_x, block_size_y, elem_per_block); 1259 sycl::range<3> global_range(grid * block_size_x, block_size_y, elem_per_block); 1260 sycl::nd_range<3> kernel_range(global_range, local_range); 1261 1262 sycl_queue.parallel_for<CeedOperatorSyclLinearAssembleFallback>(kernel_range, [=](sycl::nd_item<3> work_item) { 1263 const CeedInt blockIdx = work_item.get_group(0); 1264 const CeedInt gridDimx = work_item.get_group_range(0); 1265 const CeedInt threadIdx = work_item.get_local_id(0); 1266 const CeedInt threadIdz = work_item.get_local_id(2); 1267 const CeedInt blockDimz = work_item.get_local_range(2); 1268 1269 const int l = threadIdx; // The output column index of each B^TDB operation 1270 // such that we have (Bout^T)_ij D_jk Bin_kl = C_il 1271 for (CeedInt e = blockIdx * blockDimz + threadIdz; e < num_elem; e += gridDimx * blockDimz) { 1272 for (CeedInt comp_in = 0; comp_in < num_comp; comp_in++) { 1273 for (CeedInt comp_out = 0; comp_out < num_comp; comp_out++) { 1274 for (CeedInt i = 0; i < num_nodes; i++) { 1275 CeedScalar result = 0.0; 1276 CeedInt qf_index_comp = q_comp_in_stride * comp_in + q_comp_out_stride * comp_out + q_e_stride * e; 1277 for (CeedInt e_mode_in = 0; e_mode_in < num_e_mode_in; e_mode_in++) { 1278 CeedInt b_in_index = e_mode_in * num_qpts * num_nodes; 1279 for (CeedInt e_mode_out = 0; e_mode_out < num_e_mode_out; e_mode_out++) { 1280 CeedInt b_out_index = e_mode_out * num_qpts * num_nodes; 1281 CeedInt qf_index = qf_index_comp + q_e_mode_out_stride * e_mode_out + q_e_mode_in_stride * e_mode_in; 1282 // Perform the B^T D B operation for this 'chunk' of D (the qf_array) 1283 for (CeedInt j = 0; j < num_qpts; j++) { 1284 result += B_out[b_out_index + j * num_nodes + i] * qf_array[qf_index + j] * B_in[b_in_index + j * num_nodes + l]; 1285 } 1286 } // end of e_mode_out 1287 } // end of e_mode_in 1288 CeedInt val_index = comp_in_stride * comp_in + comp_out_stride * comp_out + e_stride * e + num_nodes * i + l; 1289 values_array[val_index] = result; 1290 } // end of loop over element node index, i 1291 } // end of out component 1292 } // end of in component 1293 } // end of element loop 1294 }); 1295 return CEED_ERROR_SUCCESS; 1296 }*/ 1297 1298 //------------------------------------------------------------------------------ 1299 // Assemble matrix data for COO matrix of assembled operator. 1300 // The sparsity pattern is set by CeedOperatorLinearAssembleSymbolic. 1301 // 1302 // Note that this (and other assembly routines) currently assume only one active 1303 // input restriction/basis per operator (could have multiple basis eval modes). 1304 // TODO: allow multiple active input restrictions/basis objects 1305 //------------------------------------------------------------------------------ 1306 static int CeedSingleOperatorAssemble_Sycl(CeedOperator op, CeedInt offset, CeedVector values) { 1307 Ceed ceed; 1308 Ceed_Sycl *sycl_data; 1309 CeedScalar *values_array; 1310 const CeedScalar *qf_array; 1311 CeedVector assembled_qf = NULL; 1312 CeedElemRestriction rstr_q = NULL; 1313 CeedOperator_Sycl *impl; 1314 1315 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 1316 CeedCallBackend(CeedOperatorGetData(op, &impl)); 1317 CeedCallBackend(CeedGetData(ceed, &sycl_data)); 1318 1319 // Setup 1320 if (!impl->asmb) { 1321 CeedCallBackend(CeedSingleOperatorAssembleSetup_Sycl(op)); 1322 assert(impl->asmb != NULL); 1323 } 1324 1325 // Assemble QFunction 1326 CeedCallBackend(CeedOperatorLinearAssembleQFunctionBuildOrUpdate(op, &assembled_qf, &rstr_q, CEED_REQUEST_IMMEDIATE)); 1327 CeedCallBackend(CeedElemRestrictionDestroy(&rstr_q)); 1328 CeedCallBackend(CeedVectorGetArrayWrite(values, CEED_MEM_DEVICE, &values_array)); 1329 values_array += offset; 1330 CeedCallBackend(CeedVectorGetArrayRead(assembled_qf, CEED_MEM_DEVICE, &qf_array)); 1331 1332 // Compute B^T D B 1333 CeedCallBackend(CeedOperatorLinearAssemble_Sycl(sycl_data->sycl_queue, impl, qf_array, values_array)); 1334 1335 // Wait for kernels to be completed 1336 // Kris: Review if this is necessary -- enqueing an async barrier may be sufficient 1337 sycl_data->sycl_queue.wait_and_throw(); 1338 1339 // Restore arrays 1340 CeedCallBackend(CeedVectorRestoreArray(values, &values_array)); 1341 CeedCallBackend(CeedVectorRestoreArrayRead(assembled_qf, &qf_array)); 1342 1343 // Cleanup 1344 CeedCallBackend(CeedVectorDestroy(&assembled_qf)); 1345 return CEED_ERROR_SUCCESS; 1346 } 1347 1348 //------------------------------------------------------------------------------ 1349 // Create operator 1350 //------------------------------------------------------------------------------ 1351 int CeedOperatorCreate_Sycl(CeedOperator op) { 1352 Ceed ceed; 1353 CeedOperator_Sycl *impl; 1354 1355 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 1356 1357 CeedCallBackend(CeedCalloc(1, &impl)); 1358 CeedCallBackend(CeedOperatorSetData(op, impl)); 1359 1360 CeedCallBackend(CeedSetBackendFunctionCpp(ceed, "Operator", op, "LinearAssembleQFunction", CeedOperatorLinearAssembleQFunction_Sycl)); 1361 CeedCallBackend(CeedSetBackendFunctionCpp(ceed, "Operator", op, "LinearAssembleQFunctionUpdate", CeedOperatorLinearAssembleQFunctionUpdate_Sycl)); 1362 CeedCallBackend(CeedSetBackendFunctionCpp(ceed, "Operator", op, "LinearAssembleAddDiagonal", CeedOperatorLinearAssembleAddDiagonal_Sycl)); 1363 CeedCallBackend( 1364 CeedSetBackendFunctionCpp(ceed, "Operator", op, "LinearAssembleAddPointBlockDiagonal", CeedOperatorLinearAssembleAddPointBlockDiagonal_Sycl)); 1365 CeedCallBackend(CeedSetBackendFunctionCpp(ceed, "Operator", op, "LinearAssembleSingle", CeedSingleOperatorAssemble_Sycl)); 1366 CeedCallBackend(CeedSetBackendFunctionCpp(ceed, "Operator", op, "ApplyAdd", CeedOperatorApplyAdd_Sycl)); 1367 CeedCallBackend(CeedSetBackendFunctionCpp(ceed, "Operator", op, "Destroy", CeedOperatorDestroy_Sycl)); 1368 return CEED_ERROR_SUCCESS; 1369 } 1370 1371 //------------------------------------------------------------------------------ 1372