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