1 // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors. 2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3 // 4 // SPDX-License-Identifier: BSD-2-Clause 5 // 6 // This file is part of CEED: http://github.com/ceed 7 8 #include <ceed.h> 9 #include <ceed/backend.h> 10 #include <ceed/jit-tools.h> 11 #include <assert.h> 12 #include <cuda.h> 13 #include <cuda_runtime.h> 14 #include <stdbool.h> 15 #include <string.h> 16 17 #include "../cuda/ceed-cuda-common.h" 18 #include "../cuda/ceed-cuda-compile.h" 19 #include "ceed-cuda-ref.h" 20 21 //------------------------------------------------------------------------------ 22 // Destroy operator 23 //------------------------------------------------------------------------------ 24 static int CeedOperatorDestroy_Cuda(CeedOperator op) { 25 CeedOperator_Cuda *impl; 26 27 CeedCallBackend(CeedOperatorGetData(op, &impl)); 28 29 // Apply data 30 CeedCallBackend(CeedFree(&impl->num_points)); 31 CeedCallBackend(CeedFree(&impl->skip_rstr_in)); 32 CeedCallBackend(CeedFree(&impl->skip_rstr_out)); 33 CeedCallBackend(CeedFree(&impl->apply_add_basis_out)); 34 CeedCallBackend(CeedFree(&impl->input_field_order)); 35 CeedCallBackend(CeedFree(&impl->output_field_order)); 36 CeedCallBackend(CeedFree(&impl->input_states)); 37 38 for (CeedInt i = 0; i < impl->num_inputs; i++) { 39 CeedCallBackend(CeedVectorDestroy(&impl->e_vecs_in[i])); 40 CeedCallBackend(CeedVectorDestroy(&impl->q_vecs_in[i])); 41 } 42 CeedCallBackend(CeedFree(&impl->e_vecs_in)); 43 CeedCallBackend(CeedFree(&impl->q_vecs_in)); 44 45 for (CeedInt i = 0; i < impl->num_outputs; i++) { 46 CeedCallBackend(CeedVectorDestroy(&impl->e_vecs_out[i])); 47 CeedCallBackend(CeedVectorDestroy(&impl->q_vecs_out[i])); 48 } 49 CeedCallBackend(CeedFree(&impl->e_vecs_out)); 50 CeedCallBackend(CeedFree(&impl->q_vecs_out)); 51 CeedCallBackend(CeedVectorDestroy(&impl->point_coords_elem)); 52 53 // QFunction assembly data 54 for (CeedInt i = 0; i < impl->num_active_in; i++) { 55 CeedCallBackend(CeedVectorDestroy(&impl->qf_active_in[i])); 56 } 57 CeedCallBackend(CeedFree(&impl->qf_active_in)); 58 59 // Diag data 60 if (impl->diag) { 61 Ceed ceed; 62 63 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 64 if (impl->diag->module) { 65 CeedCallCuda(ceed, cuModuleUnload(impl->diag->module)); 66 } 67 if (impl->diag->module_point_block) { 68 CeedCallCuda(ceed, cuModuleUnload(impl->diag->module_point_block)); 69 } 70 CeedCallCuda(ceed, cudaFree(impl->diag->d_eval_modes_in)); 71 CeedCallCuda(ceed, cudaFree(impl->diag->d_eval_modes_out)); 72 CeedCallCuda(ceed, cudaFree(impl->diag->d_identity)); 73 CeedCallCuda(ceed, cudaFree(impl->diag->d_interp_in)); 74 CeedCallCuda(ceed, cudaFree(impl->diag->d_interp_out)); 75 CeedCallCuda(ceed, cudaFree(impl->diag->d_grad_in)); 76 CeedCallCuda(ceed, cudaFree(impl->diag->d_grad_out)); 77 CeedCallCuda(ceed, cudaFree(impl->diag->d_div_in)); 78 CeedCallCuda(ceed, cudaFree(impl->diag->d_div_out)); 79 CeedCallCuda(ceed, cudaFree(impl->diag->d_curl_in)); 80 CeedCallCuda(ceed, cudaFree(impl->diag->d_curl_out)); 81 CeedCallBackend(CeedElemRestrictionDestroy(&impl->diag->diag_rstr)); 82 CeedCallBackend(CeedElemRestrictionDestroy(&impl->diag->point_block_diag_rstr)); 83 CeedCallBackend(CeedVectorDestroy(&impl->diag->elem_diag)); 84 CeedCallBackend(CeedVectorDestroy(&impl->diag->point_block_elem_diag)); 85 } 86 CeedCallBackend(CeedFree(&impl->diag)); 87 88 if (impl->asmb) { 89 Ceed ceed; 90 91 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 92 CeedCallCuda(ceed, cuModuleUnload(impl->asmb->module)); 93 CeedCallCuda(ceed, cudaFree(impl->asmb->d_B_in)); 94 CeedCallCuda(ceed, cudaFree(impl->asmb->d_B_out)); 95 } 96 CeedCallBackend(CeedFree(&impl->asmb)); 97 98 CeedCallBackend(CeedFree(&impl)); 99 return CEED_ERROR_SUCCESS; 100 } 101 102 //------------------------------------------------------------------------------ 103 // Setup infields or outfields 104 //------------------------------------------------------------------------------ 105 static int CeedOperatorSetupFields_Cuda(CeedQFunction qf, CeedOperator op, bool is_input, bool is_at_points, bool *skip_rstr, bool *apply_add_basis, 106 CeedVector *e_vecs, CeedVector *q_vecs, CeedInt num_fields, CeedInt Q, CeedInt num_elem) { 107 Ceed ceed; 108 CeedQFunctionField *qf_fields; 109 CeedOperatorField *op_fields; 110 111 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 112 if (is_input) { 113 CeedCallBackend(CeedOperatorGetFields(op, NULL, &op_fields, NULL, NULL)); 114 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_fields, NULL, NULL)); 115 } else { 116 CeedCallBackend(CeedOperatorGetFields(op, NULL, NULL, NULL, &op_fields)); 117 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qf_fields)); 118 } 119 120 // Loop over fields 121 for (CeedInt i = 0; i < num_fields; i++) { 122 bool is_strided = false, skip_e_vec = false; 123 CeedSize q_size; 124 CeedInt size; 125 CeedEvalMode eval_mode; 126 CeedBasis basis; 127 128 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode)); 129 if (eval_mode != CEED_EVAL_WEIGHT) { 130 CeedElemRestriction elem_rstr; 131 132 // Check whether this field can skip the element restriction: 133 // Must be passive input, with eval_mode NONE, and have a strided restriction with CEED_STRIDES_BACKEND. 134 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_fields[i], &elem_rstr)); 135 136 // First, check whether the field is input or output: 137 if (is_input) { 138 CeedVector l_vec; 139 140 // Check for passive input 141 CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &l_vec)); 142 if (l_vec != CEED_VECTOR_ACTIVE && eval_mode == CEED_EVAL_NONE) { 143 // Check for strided restriction 144 CeedCallBackend(CeedElemRestrictionIsStrided(elem_rstr, &is_strided)); 145 if (is_strided) { 146 // Check if vector is already in preferred backend ordering 147 CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &skip_e_vec)); 148 } 149 } 150 } 151 if (skip_e_vec) { 152 // Either an active field or strided local vec in backend ordering 153 e_vecs[i] = NULL; 154 } else { 155 CeedCallBackend(CeedElemRestrictionCreateVector(elem_rstr, NULL, &e_vecs[i])); 156 } 157 } 158 159 switch (eval_mode) { 160 case CEED_EVAL_NONE: 161 CeedCallBackend(CeedQFunctionFieldGetSize(qf_fields[i], &size)); 162 q_size = (CeedSize)num_elem * Q * size; 163 CeedCallBackend(CeedVectorCreate(ceed, q_size, &q_vecs[i])); 164 break; 165 case CEED_EVAL_INTERP: 166 case CEED_EVAL_GRAD: 167 case CEED_EVAL_DIV: 168 case CEED_EVAL_CURL: 169 CeedCallBackend(CeedQFunctionFieldGetSize(qf_fields[i], &size)); 170 q_size = (CeedSize)num_elem * Q * size; 171 CeedCallBackend(CeedVectorCreate(ceed, q_size, &q_vecs[i])); 172 break; 173 case CEED_EVAL_WEIGHT: // Only on input fields 174 CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis)); 175 q_size = (CeedSize)num_elem * Q; 176 CeedCallBackend(CeedVectorCreate(ceed, q_size, &q_vecs[i])); 177 if (is_at_points) { 178 CeedInt num_points[num_elem]; 179 180 for (CeedInt i = 0; i < num_elem; i++) num_points[i] = Q; 181 CeedCallBackend( 182 CeedBasisApplyAtPoints(basis, num_elem, num_points, CEED_NOTRANSPOSE, CEED_EVAL_WEIGHT, CEED_VECTOR_NONE, CEED_VECTOR_NONE, q_vecs[i])); 183 } else { 184 CeedCallBackend(CeedBasisApply(basis, num_elem, CEED_NOTRANSPOSE, CEED_EVAL_WEIGHT, CEED_VECTOR_NONE, q_vecs[i])); 185 } 186 break; 187 } 188 } 189 // Drop duplicate restrictions 190 if (is_input) { 191 for (CeedInt i = 0; i < num_fields; i++) { 192 CeedVector vec_i; 193 CeedElemRestriction rstr_i; 194 195 CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec_i)); 196 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_fields[i], &rstr_i)); 197 for (CeedInt j = i + 1; j < num_fields; j++) { 198 CeedVector vec_j; 199 CeedElemRestriction rstr_j; 200 201 CeedCallBackend(CeedOperatorFieldGetVector(op_fields[j], &vec_j)); 202 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_fields[j], &rstr_j)); 203 if (vec_i == vec_j && rstr_i == rstr_j) { 204 CeedCallBackend(CeedVectorReferenceCopy(e_vecs[i], &e_vecs[j])); 205 skip_rstr[j] = true; 206 } 207 } 208 } 209 } else { 210 for (CeedInt i = num_fields - 1; i >= 0; i--) { 211 CeedVector vec_i; 212 CeedElemRestriction rstr_i; 213 214 CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec_i)); 215 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_fields[i], &rstr_i)); 216 for (CeedInt j = i - 1; j >= 0; j--) { 217 CeedVector vec_j; 218 CeedElemRestriction rstr_j; 219 220 CeedCallBackend(CeedOperatorFieldGetVector(op_fields[j], &vec_j)); 221 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_fields[j], &rstr_j)); 222 if (vec_i == vec_j && rstr_i == rstr_j) { 223 CeedCallBackend(CeedVectorReferenceCopy(e_vecs[i], &e_vecs[j])); 224 skip_rstr[j] = true; 225 apply_add_basis[i] = true; 226 } 227 } 228 } 229 } 230 return CEED_ERROR_SUCCESS; 231 } 232 233 //------------------------------------------------------------------------------ 234 // CeedOperator needs to connect all the named fields (be they active or passive) to the named inputs and outputs of its CeedQFunction. 235 //------------------------------------------------------------------------------ 236 static int CeedOperatorSetup_Cuda(CeedOperator op) { 237 Ceed ceed; 238 bool is_setup_done; 239 CeedInt Q, num_elem, num_input_fields, num_output_fields; 240 CeedQFunctionField *qf_input_fields, *qf_output_fields; 241 CeedQFunction qf; 242 CeedOperatorField *op_input_fields, *op_output_fields; 243 CeedOperator_Cuda *impl; 244 245 CeedCallBackend(CeedOperatorIsSetupDone(op, &is_setup_done)); 246 if (is_setup_done) return CEED_ERROR_SUCCESS; 247 248 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 249 CeedCallBackend(CeedOperatorGetData(op, &impl)); 250 CeedCallBackend(CeedOperatorGetQFunction(op, &qf)); 251 CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q)); 252 CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem)); 253 CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields)); 254 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields)); 255 256 // Allocate 257 CeedCallBackend(CeedCalloc(num_input_fields, &impl->e_vecs_in)); 258 CeedCallBackend(CeedCalloc(num_output_fields, &impl->e_vecs_out)); 259 CeedCallBackend(CeedCalloc(num_input_fields, &impl->skip_rstr_in)); 260 CeedCallBackend(CeedCalloc(num_output_fields, &impl->skip_rstr_out)); 261 CeedCallBackend(CeedCalloc(num_output_fields, &impl->apply_add_basis_out)); 262 CeedCallBackend(CeedCalloc(num_input_fields, &impl->input_field_order)); 263 CeedCallBackend(CeedCalloc(num_output_fields, &impl->output_field_order)); 264 CeedCallBackend(CeedCalloc(num_input_fields, &impl->input_states)); 265 CeedCallBackend(CeedCalloc(num_input_fields, &impl->q_vecs_in)); 266 CeedCallBackend(CeedCalloc(num_output_fields, &impl->q_vecs_out)); 267 impl->num_inputs = num_input_fields; 268 impl->num_outputs = num_output_fields; 269 270 // Set up infield and outfield e-vecs and q-vecs 271 CeedCallBackend( 272 CeedOperatorSetupFields_Cuda(qf, op, true, false, impl->skip_rstr_in, NULL, impl->e_vecs_in, impl->q_vecs_in, num_input_fields, Q, num_elem)); 273 CeedCallBackend(CeedOperatorSetupFields_Cuda(qf, op, false, false, impl->skip_rstr_out, impl->apply_add_basis_out, impl->e_vecs_out, 274 impl->q_vecs_out, num_output_fields, Q, num_elem)); 275 276 // Reorder fields to allow reuse of buffers 277 impl->max_active_e_vec_len = 0; 278 { 279 bool is_ordered[CEED_FIELD_MAX]; 280 CeedInt curr_index = 0; 281 282 for (CeedInt i = 0; i < num_input_fields; i++) is_ordered[i] = false; 283 for (CeedInt i = 0; i < num_input_fields; i++) { 284 CeedSize e_vec_len_i; 285 CeedVector vec_i; 286 CeedElemRestriction rstr_i; 287 288 if (is_ordered[i]) continue; 289 is_ordered[i] = true; 290 impl->input_field_order[curr_index] = i; 291 curr_index++; 292 CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec_i)); 293 if (vec_i == CEED_VECTOR_NONE) continue; // CEED_EVAL_WEIGHT 294 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &rstr_i)); 295 CeedCallBackend(CeedElemRestrictionGetEVectorSize(rstr_i, &e_vec_len_i)); 296 impl->max_active_e_vec_len = e_vec_len_i > impl->max_active_e_vec_len ? e_vec_len_i : impl->max_active_e_vec_len; 297 for (CeedInt j = i + 1; j < num_input_fields; j++) { 298 CeedVector vec_j; 299 CeedElemRestriction rstr_j; 300 301 CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[j], &vec_j)); 302 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[j], &rstr_j)); 303 if (rstr_i == rstr_j && vec_i == vec_j) { 304 is_ordered[j] = true; 305 impl->input_field_order[curr_index] = j; 306 curr_index++; 307 } 308 } 309 } 310 } 311 { 312 bool is_ordered[CEED_FIELD_MAX]; 313 CeedInt curr_index = 0; 314 315 for (CeedInt i = 0; i < num_output_fields; i++) is_ordered[i] = false; 316 for (CeedInt i = 0; i < num_output_fields; i++) { 317 CeedSize e_vec_len_i; 318 CeedVector vec_i; 319 CeedElemRestriction rstr_i; 320 321 if (is_ordered[i]) continue; 322 is_ordered[i] = true; 323 impl->output_field_order[curr_index] = i; 324 curr_index++; 325 CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[i], &vec_i)); 326 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &rstr_i)); 327 CeedCallBackend(CeedElemRestrictionGetEVectorSize(rstr_i, &e_vec_len_i)); 328 impl->max_active_e_vec_len = e_vec_len_i > impl->max_active_e_vec_len ? e_vec_len_i : impl->max_active_e_vec_len; 329 for (CeedInt j = i + 1; j < num_output_fields; j++) { 330 CeedVector vec_j; 331 CeedElemRestriction rstr_j; 332 333 CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[j], &vec_j)); 334 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[j], &rstr_j)); 335 if (rstr_i == rstr_j && vec_i == vec_j) { 336 is_ordered[j] = true; 337 impl->output_field_order[curr_index] = j; 338 curr_index++; 339 } 340 } 341 } 342 } 343 CeedCallBackend(CeedOperatorSetSetupDone(op)); 344 return CEED_ERROR_SUCCESS; 345 } 346 347 //------------------------------------------------------------------------------ 348 // Restrict Operator Inputs 349 //------------------------------------------------------------------------------ 350 static inline int CeedOperatorInputRestrict_Cuda(CeedOperatorField op_input_field, CeedQFunctionField qf_input_field, CeedInt input_field, 351 CeedVector in_vec, const bool skip_active, CeedScalar **e_data, CeedOperator_Cuda *impl, 352 CeedRequest *request) { 353 CeedEvalMode eval_mode; 354 CeedVector l_vec, e_vec = impl->e_vecs_in[input_field]; 355 356 // Get input vector 357 CeedCallBackend(CeedOperatorFieldGetVector(op_input_field, &l_vec)); 358 if (l_vec == CEED_VECTOR_ACTIVE) { 359 if (skip_active) return CEED_ERROR_SUCCESS; 360 else l_vec = in_vec; 361 } 362 363 // Restriction action 364 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_field, &eval_mode)); 365 if (eval_mode == CEED_EVAL_WEIGHT) { // Skip 366 } else { 367 if (!e_vec) { 368 // No restriction for this field; read data directly from vec. 369 CeedCallBackend(CeedVectorGetArrayRead(l_vec, CEED_MEM_DEVICE, (const CeedScalar **)e_data)); 370 } else { 371 // Restrict, if necessary 372 if (!impl->skip_rstr_in[input_field]) { 373 uint64_t state; 374 375 CeedCallBackend(CeedVectorGetState(l_vec, &state)); 376 if (state != impl->input_states[input_field] || l_vec == in_vec) { 377 CeedElemRestriction elem_rstr; 378 379 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_field, &elem_rstr)); 380 CeedCallBackend(CeedElemRestrictionApply(elem_rstr, CEED_NOTRANSPOSE, l_vec, e_vec, request)); 381 } 382 impl->input_states[input_field] = state; 383 } 384 // Get e-vec 385 CeedCallBackend(CeedVectorGetArrayRead(e_vec, CEED_MEM_DEVICE, (const CeedScalar **)e_data)); 386 } 387 } 388 return CEED_ERROR_SUCCESS; 389 } 390 391 //------------------------------------------------------------------------------ 392 // Input Basis Action 393 //------------------------------------------------------------------------------ 394 static inline int CeedOperatorInputBasis_Cuda(CeedOperatorField op_input_field, CeedQFunctionField qf_input_field, CeedInt input_field, 395 CeedInt num_elem, const bool skip_active, CeedScalar *e_data, CeedOperator_Cuda *impl) { 396 CeedEvalMode eval_mode; 397 CeedVector e_vec = impl->e_vecs_in[input_field], q_vec = impl->q_vecs_in[input_field]; 398 399 // Skip active input 400 if (skip_active) { 401 CeedVector l_vec; 402 403 CeedCallBackend(CeedOperatorFieldGetVector(op_input_field, &l_vec)); 404 if (l_vec == CEED_VECTOR_ACTIVE) return CEED_ERROR_SUCCESS; 405 } 406 407 // Basis action 408 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_field, &eval_mode)); 409 switch (eval_mode) { 410 case CEED_EVAL_NONE: 411 CeedCallBackend(CeedVectorSetArray(q_vec, CEED_MEM_DEVICE, CEED_USE_POINTER, e_data)); 412 break; 413 case CEED_EVAL_INTERP: 414 case CEED_EVAL_GRAD: 415 case CEED_EVAL_DIV: 416 case CEED_EVAL_CURL: { 417 CeedBasis basis; 418 419 CeedCallBackend(CeedOperatorFieldGetBasis(op_input_field, &basis)); 420 CeedCallBackend(CeedBasisApply(basis, num_elem, CEED_NOTRANSPOSE, eval_mode, e_vec, q_vec)); 421 break; 422 } 423 case CEED_EVAL_WEIGHT: 424 break; // No action 425 } 426 return CEED_ERROR_SUCCESS; 427 } 428 429 //------------------------------------------------------------------------------ 430 // Restore Input Vectors 431 //------------------------------------------------------------------------------ 432 static inline int CeedOperatorInputRestore_Cuda(CeedOperatorField op_input_field, CeedQFunctionField qf_input_field, CeedInt input_field, 433 const bool skip_active, CeedScalar **e_data, CeedOperator_Cuda *impl) { 434 CeedEvalMode eval_mode; 435 CeedVector l_vec, e_vec = impl->e_vecs_in[input_field]; 436 437 // Skip active input 438 CeedCallBackend(CeedOperatorFieldGetVector(op_input_field, &l_vec)); 439 if (skip_active && l_vec == CEED_VECTOR_ACTIVE) return CEED_ERROR_SUCCESS; 440 441 // Restore e-vec 442 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_field, &eval_mode)); 443 if (eval_mode == CEED_EVAL_WEIGHT) { // Skip 444 } else { 445 if (!e_vec) { // This was a skip_restriction case 446 CeedCallBackend(CeedVectorRestoreArrayRead(l_vec, (const CeedScalar **)e_data)); 447 } else { 448 CeedCallBackend(CeedVectorRestoreArrayRead(e_vec, (const CeedScalar **)e_data)); 449 } 450 } 451 return CEED_ERROR_SUCCESS; 452 } 453 454 //------------------------------------------------------------------------------ 455 // Apply and add to output 456 //------------------------------------------------------------------------------ 457 static int CeedOperatorApplyAdd_Cuda(CeedOperator op, CeedVector in_vec, CeedVector out_vec, CeedRequest *request) { 458 CeedInt Q, num_elem, num_input_fields, num_output_fields; 459 CeedScalar *e_data_in[CEED_FIELD_MAX] = {NULL}, *e_data_out[CEED_FIELD_MAX] = {NULL}; 460 CeedQFunctionField *qf_input_fields, *qf_output_fields; 461 CeedQFunction qf; 462 CeedOperatorField *op_input_fields, *op_output_fields; 463 CeedOperator_Cuda *impl; 464 465 CeedCallBackend(CeedOperatorGetData(op, &impl)); 466 CeedCallBackend(CeedOperatorGetQFunction(op, &qf)); 467 CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q)); 468 CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem)); 469 CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields)); 470 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields)); 471 472 // Setup 473 CeedCallBackend(CeedOperatorSetup_Cuda(op)); 474 475 // Process inputs 476 for (CeedInt i = 0; i < num_input_fields; i++) { 477 CeedInt field = impl->input_field_order[i]; 478 479 CeedCallBackend( 480 CeedOperatorInputRestrict_Cuda(op_input_fields[field], qf_input_fields[field], field, in_vec, false, &e_data_in[field], impl, request)); 481 CeedCallBackend(CeedOperatorInputBasis_Cuda(op_input_fields[field], qf_input_fields[field], field, num_elem, false, e_data_in[field], impl)); 482 } 483 484 // Output pointers, as necessary 485 for (CeedInt i = 0; i < num_output_fields; i++) { 486 CeedEvalMode eval_mode; 487 488 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 489 if (eval_mode == CEED_EVAL_NONE) { 490 // Set the output Q-Vector to use the E-Vector data directly. 491 CeedCallBackend(CeedVectorGetArrayWrite(impl->e_vecs_out[i], CEED_MEM_DEVICE, &e_data_out[i])); 492 CeedCallBackend(CeedVectorSetArray(impl->q_vecs_out[i], CEED_MEM_DEVICE, CEED_USE_POINTER, e_data_out[i])); 493 } 494 } 495 496 // Q function 497 CeedCallBackend(CeedQFunctionApply(qf, num_elem * Q, impl->q_vecs_in, impl->q_vecs_out)); 498 499 // Restore input arrays 500 for (CeedInt i = 0; i < num_input_fields; i++) { 501 CeedCallBackend(CeedOperatorInputRestore_Cuda(op_input_fields[i], qf_input_fields[i], i, false, &e_data_in[i], impl)); 502 } 503 504 // Output basis apply if needed 505 for (CeedInt i = 0; i < num_output_fields; i++) { 506 CeedInt field = impl->output_field_order[i]; 507 CeedEvalMode eval_mode; 508 CeedVector l_vec, e_vec = impl->e_vecs_out[field], q_vec = impl->q_vecs_out[field]; 509 CeedElemRestriction elem_rstr; 510 CeedBasis basis; 511 512 // Output vector 513 CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[field], &l_vec)); 514 if (l_vec == CEED_VECTOR_ACTIVE) l_vec = out_vec; 515 516 // Basis action 517 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[field], &eval_mode)); 518 switch (eval_mode) { 519 case CEED_EVAL_NONE: 520 break; // No action 521 case CEED_EVAL_INTERP: 522 case CEED_EVAL_GRAD: 523 case CEED_EVAL_DIV: 524 case CEED_EVAL_CURL: 525 CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[field], &basis)); 526 if (impl->apply_add_basis_out[field]) { 527 CeedCallBackend(CeedBasisApplyAdd(basis, num_elem, CEED_TRANSPOSE, eval_mode, q_vec, e_vec)); 528 } else { 529 CeedCallBackend(CeedBasisApply(basis, num_elem, CEED_TRANSPOSE, eval_mode, q_vec, e_vec)); 530 } 531 break; 532 // LCOV_EXCL_START 533 case CEED_EVAL_WEIGHT: { 534 return CeedError(CeedOperatorReturnCeed(op), CEED_ERROR_BACKEND, "CEED_EVAL_WEIGHT cannot be an output evaluation mode"); 535 // LCOV_EXCL_STOP 536 } 537 } 538 539 // Restore evec 540 if (eval_mode == CEED_EVAL_NONE) { 541 CeedCallBackend(CeedVectorRestoreArray(e_vec, &e_data_out[field])); 542 } 543 544 // Restrict 545 if (impl->skip_rstr_out[field]) continue; 546 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[field], &elem_rstr)); 547 CeedCallBackend(CeedElemRestrictionApply(elem_rstr, CEED_TRANSPOSE, e_vec, l_vec, request)); 548 } 549 return CEED_ERROR_SUCCESS; 550 } 551 552 //------------------------------------------------------------------------------ 553 // CeedOperator needs to connect all the named fields (be they active or passive) to the named inputs and outputs of its CeedQFunction. 554 //------------------------------------------------------------------------------ 555 static int CeedOperatorSetupAtPoints_Cuda(CeedOperator op) { 556 Ceed ceed; 557 bool is_setup_done; 558 CeedInt max_num_points = -1, num_elem, num_input_fields, num_output_fields; 559 CeedQFunctionField *qf_input_fields, *qf_output_fields; 560 CeedQFunction qf; 561 CeedOperatorField *op_input_fields, *op_output_fields; 562 CeedOperator_Cuda *impl; 563 564 CeedCallBackend(CeedOperatorIsSetupDone(op, &is_setup_done)); 565 if (is_setup_done) return CEED_ERROR_SUCCESS; 566 567 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 568 CeedCallBackend(CeedOperatorGetData(op, &impl)); 569 CeedCallBackend(CeedOperatorGetQFunction(op, &qf)); 570 CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem)); 571 CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields)); 572 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields)); 573 { 574 CeedElemRestriction rstr_points = NULL; 575 576 CeedCallBackend(CeedOperatorAtPointsGetPoints(op, &rstr_points, NULL)); 577 CeedCallBackend(CeedElemRestrictionGetMaxPointsInElement(rstr_points, &max_num_points)); 578 CeedCallBackend(CeedCalloc(num_elem, &impl->num_points)); 579 for (CeedInt e = 0; e < num_elem; e++) { 580 CeedInt num_points_elem; 581 582 CeedCallBackend(CeedElemRestrictionGetNumPointsInElement(rstr_points, e, &num_points_elem)); 583 impl->num_points[e] = num_points_elem; 584 } 585 } 586 impl->max_num_points = max_num_points; 587 588 // Allocate 589 CeedCallBackend(CeedCalloc(num_input_fields, &impl->e_vecs_in)); 590 CeedCallBackend(CeedCalloc(num_output_fields, &impl->e_vecs_out)); 591 CeedCallBackend(CeedCalloc(num_input_fields, &impl->skip_rstr_in)); 592 CeedCallBackend(CeedCalloc(num_output_fields, &impl->skip_rstr_out)); 593 CeedCallBackend(CeedCalloc(num_output_fields, &impl->apply_add_basis_out)); 594 CeedCallBackend(CeedCalloc(num_input_fields, &impl->input_field_order)); 595 CeedCallBackend(CeedCalloc(num_output_fields, &impl->output_field_order)); 596 CeedCallBackend(CeedCalloc(num_input_fields, &impl->input_states)); 597 CeedCallBackend(CeedCalloc(num_input_fields, &impl->q_vecs_in)); 598 CeedCallBackend(CeedCalloc(num_output_fields, &impl->q_vecs_out)); 599 impl->num_inputs = num_input_fields; 600 impl->num_outputs = num_output_fields; 601 602 // Set up infield and outfield e-vecs and q-vecs 603 CeedCallBackend(CeedOperatorSetupFields_Cuda(qf, op, true, true, impl->skip_rstr_in, NULL, impl->e_vecs_in, impl->q_vecs_in, num_input_fields, 604 max_num_points, num_elem)); 605 CeedCallBackend(CeedOperatorSetupFields_Cuda(qf, op, false, true, impl->skip_rstr_out, impl->apply_add_basis_out, impl->e_vecs_out, 606 impl->q_vecs_out, num_output_fields, max_num_points, num_elem)); 607 608 // Reorder fields to allow reuse of buffers 609 { 610 bool is_ordered[CEED_FIELD_MAX]; 611 CeedInt curr_index = 0; 612 613 for (CeedInt i = 0; i < num_input_fields; i++) is_ordered[i] = false; 614 for (CeedInt i = 0; i < num_input_fields; i++) { 615 CeedVector vec_i; 616 CeedElemRestriction rstr_i; 617 618 if (is_ordered[i]) continue; 619 is_ordered[i] = true; 620 impl->input_field_order[curr_index] = i; 621 curr_index++; 622 CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec_i)); 623 if (vec_i == CEED_VECTOR_NONE) continue; // CEED_EVAL_WEIGHT 624 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &rstr_i)); 625 for (CeedInt j = i + 1; j < num_input_fields; j++) { 626 CeedVector vec_j; 627 CeedElemRestriction rstr_j; 628 629 CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[j], &vec_j)); 630 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[j], &rstr_j)); 631 if (rstr_i == rstr_j && vec_i == vec_j) { 632 is_ordered[j] = true; 633 impl->input_field_order[curr_index] = j; 634 curr_index++; 635 } 636 } 637 } 638 } 639 { 640 bool is_ordered[CEED_FIELD_MAX]; 641 CeedInt curr_index = 0; 642 643 for (CeedInt i = 0; i < num_output_fields; i++) is_ordered[i] = false; 644 for (CeedInt i = 0; i < num_output_fields; i++) { 645 CeedVector vec_i; 646 CeedElemRestriction rstr_i; 647 648 if (is_ordered[i]) continue; 649 is_ordered[i] = true; 650 impl->output_field_order[curr_index] = i; 651 curr_index++; 652 CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[i], &vec_i)); 653 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &rstr_i)); 654 for (CeedInt j = i + 1; j < num_output_fields; j++) { 655 CeedVector vec_j; 656 CeedElemRestriction rstr_j; 657 658 CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[j], &vec_j)); 659 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[j], &rstr_j)); 660 if (rstr_i == rstr_j && vec_i == vec_j) { 661 is_ordered[j] = true; 662 impl->output_field_order[curr_index] = j; 663 curr_index++; 664 } 665 } 666 } 667 } 668 CeedCallBackend(CeedOperatorSetSetupDone(op)); 669 return CEED_ERROR_SUCCESS; 670 } 671 672 //------------------------------------------------------------------------------ 673 // Input Basis Action AtPoints 674 //------------------------------------------------------------------------------ 675 static inline int CeedOperatorInputBasisAtPoints_Cuda(CeedOperatorField op_input_field, CeedQFunctionField qf_input_field, CeedInt input_field, 676 CeedInt num_elem, const CeedInt *num_points, const bool skip_active, CeedScalar *e_data, 677 CeedOperator_Cuda *impl) { 678 CeedEvalMode eval_mode; 679 CeedVector e_vec = impl->e_vecs_in[input_field], q_vec = impl->q_vecs_in[input_field]; 680 681 // Skip active input 682 if (skip_active) { 683 CeedVector l_vec; 684 685 CeedCallBackend(CeedOperatorFieldGetVector(op_input_field, &l_vec)); 686 if (l_vec == CEED_VECTOR_ACTIVE) return CEED_ERROR_SUCCESS; 687 } 688 689 // Basis action 690 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_field, &eval_mode)); 691 switch (eval_mode) { 692 case CEED_EVAL_NONE: 693 CeedCallBackend(CeedVectorSetArray(q_vec, CEED_MEM_DEVICE, CEED_USE_POINTER, e_data)); 694 break; 695 case CEED_EVAL_INTERP: 696 case CEED_EVAL_GRAD: 697 case CEED_EVAL_DIV: 698 case CEED_EVAL_CURL: { 699 CeedBasis basis; 700 701 CeedCallBackend(CeedOperatorFieldGetBasis(op_input_field, &basis)); 702 CeedCallBackend(CeedBasisApplyAtPoints(basis, num_elem, num_points, CEED_NOTRANSPOSE, eval_mode, impl->point_coords_elem, e_vec, q_vec)); 703 break; 704 } 705 case CEED_EVAL_WEIGHT: 706 break; // No action 707 } 708 return CEED_ERROR_SUCCESS; 709 } 710 711 //------------------------------------------------------------------------------ 712 // Apply and add to output AtPoints 713 //------------------------------------------------------------------------------ 714 static int CeedOperatorApplyAddAtPoints_Cuda(CeedOperator op, CeedVector in_vec, CeedVector out_vec, CeedRequest *request) { 715 CeedInt max_num_points, *num_points, num_elem, num_input_fields, num_output_fields; 716 CeedScalar *e_data_in[CEED_FIELD_MAX] = {NULL}, *e_data_out[CEED_FIELD_MAX] = {NULL}; 717 CeedQFunctionField *qf_input_fields, *qf_output_fields; 718 CeedQFunction qf; 719 CeedOperatorField *op_input_fields, *op_output_fields; 720 CeedOperator_Cuda *impl; 721 722 CeedCallBackend(CeedOperatorGetData(op, &impl)); 723 CeedCallBackend(CeedOperatorGetQFunction(op, &qf)); 724 CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem)); 725 CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields)); 726 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields)); 727 728 // Setup 729 CeedCallBackend(CeedOperatorSetupAtPoints_Cuda(op)); 730 num_points = impl->num_points; 731 max_num_points = impl->max_num_points; 732 733 // Get point coordinates 734 if (!impl->point_coords_elem) { 735 CeedVector point_coords = NULL; 736 CeedElemRestriction rstr_points = NULL; 737 738 CeedCallBackend(CeedOperatorAtPointsGetPoints(op, &rstr_points, &point_coords)); 739 CeedCallBackend(CeedElemRestrictionCreateVector(rstr_points, NULL, &impl->point_coords_elem)); 740 CeedCallBackend(CeedElemRestrictionApply(rstr_points, CEED_NOTRANSPOSE, point_coords, impl->point_coords_elem, request)); 741 } 742 743 // Process inputs 744 for (CeedInt i = 0; i < num_input_fields; i++) { 745 CeedInt field = impl->input_field_order[i]; 746 747 CeedCallBackend( 748 CeedOperatorInputRestrict_Cuda(op_input_fields[field], qf_input_fields[field], field, in_vec, false, &e_data_in[field], impl, request)); 749 CeedCallBackend(CeedOperatorInputBasisAtPoints_Cuda(op_input_fields[field], qf_input_fields[field], field, num_elem, num_points, false, 750 e_data_in[field], impl)); 751 } 752 753 // Output pointers, as necessary 754 for (CeedInt i = 0; i < num_output_fields; i++) { 755 CeedEvalMode eval_mode; 756 757 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 758 if (eval_mode == CEED_EVAL_NONE) { 759 // Set the output Q-Vector to use the E-Vector data directly. 760 CeedCallBackend(CeedVectorGetArrayWrite(impl->e_vecs_out[i], CEED_MEM_DEVICE, &e_data_out[i])); 761 CeedCallBackend(CeedVectorSetArray(impl->q_vecs_out[i], CEED_MEM_DEVICE, CEED_USE_POINTER, e_data_out[i])); 762 } 763 } 764 765 // Q function 766 CeedCallBackend(CeedQFunctionApply(qf, num_elem * max_num_points, impl->q_vecs_in, impl->q_vecs_out)); 767 768 // Restore input arrays 769 for (CeedInt i = 0; i < num_input_fields; i++) { 770 CeedCallBackend(CeedOperatorInputRestore_Cuda(op_input_fields[i], qf_input_fields[i], i, false, &e_data_in[i], impl)); 771 } 772 773 // Output basis apply if needed 774 for (CeedInt i = 0; i < num_output_fields; i++) { 775 CeedInt field = impl->output_field_order[i]; 776 CeedEvalMode eval_mode; 777 CeedVector l_vec, e_vec = impl->e_vecs_out[field], q_vec = impl->q_vecs_out[field]; 778 CeedElemRestriction elem_rstr; 779 CeedBasis basis; 780 781 // Output vector 782 CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[field], &l_vec)); 783 if (l_vec == CEED_VECTOR_ACTIVE) l_vec = out_vec; 784 785 // Basis action 786 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[field], &eval_mode)); 787 switch (eval_mode) { 788 case CEED_EVAL_NONE: 789 break; // No action 790 case CEED_EVAL_INTERP: 791 case CEED_EVAL_GRAD: 792 case CEED_EVAL_DIV: 793 case CEED_EVAL_CURL: 794 CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[field], &basis)); 795 if (impl->apply_add_basis_out[field]) { 796 CeedCallBackend(CeedBasisApplyAddAtPoints(basis, num_elem, num_points, CEED_TRANSPOSE, eval_mode, impl->point_coords_elem, q_vec, e_vec)); 797 } else { 798 CeedCallBackend(CeedBasisApplyAtPoints(basis, num_elem, num_points, CEED_TRANSPOSE, eval_mode, impl->point_coords_elem, q_vec, e_vec)); 799 } 800 break; 801 // LCOV_EXCL_START 802 case CEED_EVAL_WEIGHT: { 803 return CeedError(CeedOperatorReturnCeed(op), CEED_ERROR_BACKEND, "CEED_EVAL_WEIGHT cannot be an output evaluation mode"); 804 // LCOV_EXCL_STOP 805 } 806 } 807 808 // Restore evec 809 if (eval_mode == CEED_EVAL_NONE) { 810 CeedCallBackend(CeedVectorRestoreArray(e_vec, &e_data_out[field])); 811 } 812 813 // Restrict 814 if (impl->skip_rstr_out[field]) continue; 815 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[field], &elem_rstr)); 816 CeedCallBackend(CeedElemRestrictionApply(elem_rstr, CEED_TRANSPOSE, e_vec, l_vec, request)); 817 } 818 return CEED_ERROR_SUCCESS; 819 } 820 821 //------------------------------------------------------------------------------ 822 // Linear QFunction Assembly Core 823 //------------------------------------------------------------------------------ 824 static inline int CeedOperatorLinearAssembleQFunctionCore_Cuda(CeedOperator op, bool build_objects, CeedVector *assembled, CeedElemRestriction *rstr, 825 CeedRequest *request) { 826 Ceed ceed, ceed_parent; 827 CeedInt num_active_in, num_active_out, Q, num_elem, num_input_fields, num_output_fields, size; 828 CeedScalar *assembled_array, *e_data[2 * CEED_FIELD_MAX] = {NULL}; 829 CeedVector *active_inputs; 830 CeedQFunctionField *qf_input_fields, *qf_output_fields; 831 CeedQFunction qf; 832 CeedOperatorField *op_input_fields, *op_output_fields; 833 CeedOperator_Cuda *impl; 834 835 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 836 CeedCallBackend(CeedOperatorGetFallbackParentCeed(op, &ceed_parent)); 837 CeedCallBackend(CeedOperatorGetData(op, &impl)); 838 CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q)); 839 CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem)); 840 CeedCallBackend(CeedOperatorGetQFunction(op, &qf)); 841 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields)); 842 CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields)); 843 active_inputs = impl->qf_active_in; 844 num_active_in = impl->num_active_in, num_active_out = impl->num_active_out; 845 846 // Setup 847 CeedCallBackend(CeedOperatorSetup_Cuda(op)); 848 849 // Process inputs 850 for (CeedInt i = 0; i < num_input_fields; i++) { 851 CeedCallBackend(CeedOperatorInputRestrict_Cuda(op_input_fields[i], qf_input_fields[i], i, NULL, true, &e_data[i], impl, request)); 852 CeedCallBackend(CeedOperatorInputBasis_Cuda(op_input_fields[i], qf_input_fields[i], i, num_elem, true, e_data[i], impl)); 853 } 854 855 // Count number of active input fields 856 if (!num_active_in) { 857 for (CeedInt i = 0; i < num_input_fields; i++) { 858 CeedScalar *q_vec_array; 859 CeedVector l_vec; 860 861 // Check if active input 862 CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &l_vec)); 863 if (l_vec == CEED_VECTOR_ACTIVE) { 864 CeedCallBackend(CeedQFunctionFieldGetSize(qf_input_fields[i], &size)); 865 CeedCallBackend(CeedVectorSetValue(impl->q_vecs_in[i], 0.0)); 866 CeedCallBackend(CeedVectorGetArray(impl->q_vecs_in[i], CEED_MEM_DEVICE, &q_vec_array)); 867 CeedCallBackend(CeedRealloc(num_active_in + size, &active_inputs)); 868 for (CeedInt field = 0; field < size; field++) { 869 CeedSize q_size = (CeedSize)Q * num_elem; 870 871 CeedCallBackend(CeedVectorCreate(ceed, q_size, &active_inputs[num_active_in + field])); 872 CeedCallBackend( 873 CeedVectorSetArray(active_inputs[num_active_in + field], CEED_MEM_DEVICE, CEED_USE_POINTER, &q_vec_array[field * Q * num_elem])); 874 } 875 num_active_in += size; 876 CeedCallBackend(CeedVectorRestoreArray(impl->q_vecs_in[i], &q_vec_array)); 877 } 878 } 879 impl->num_active_in = num_active_in; 880 impl->qf_active_in = active_inputs; 881 } 882 883 // Count number of active output fields 884 if (!num_active_out) { 885 for (CeedInt i = 0; i < num_output_fields; i++) { 886 CeedVector l_vec; 887 888 // Check if active output 889 CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[i], &l_vec)); 890 if (l_vec == CEED_VECTOR_ACTIVE) { 891 CeedCallBackend(CeedQFunctionFieldGetSize(qf_output_fields[i], &size)); 892 num_active_out += size; 893 } 894 } 895 impl->num_active_out = num_active_out; 896 } 897 898 // Check sizes 899 CeedCheck(num_active_in > 0 && num_active_out > 0, ceed, CEED_ERROR_BACKEND, "Cannot assemble QFunction without active inputs and outputs"); 900 901 // Build objects if needed 902 if (build_objects) { 903 CeedSize l_size = (CeedSize)num_elem * Q * num_active_in * num_active_out; 904 CeedInt strides[3] = {1, num_elem * Q, Q}; /* *NOPAD* */ 905 906 // Create output restriction 907 CeedCallBackend(CeedElemRestrictionCreateStrided(ceed_parent, num_elem, Q, num_active_in * num_active_out, 908 (CeedSize)num_active_in * (CeedSize)num_active_out * (CeedSize)num_elem * (CeedSize)Q, strides, 909 rstr)); 910 // Create assembled vector 911 CeedCallBackend(CeedVectorCreate(ceed_parent, l_size, assembled)); 912 } 913 CeedCallBackend(CeedVectorSetValue(*assembled, 0.0)); 914 CeedCallBackend(CeedVectorGetArray(*assembled, CEED_MEM_DEVICE, &assembled_array)); 915 916 // Assemble QFunction 917 for (CeedInt in = 0; in < num_active_in; in++) { 918 // Set Inputs 919 CeedCallBackend(CeedVectorSetValue(active_inputs[in], 1.0)); 920 if (num_active_in > 1) { 921 CeedCallBackend(CeedVectorSetValue(active_inputs[(in + num_active_in - 1) % num_active_in], 0.0)); 922 } 923 // Set Outputs 924 for (CeedInt out = 0; out < num_output_fields; out++) { 925 CeedVector l_vec; 926 927 // Get output vector 928 CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[out], &l_vec)); 929 // Check if active output 930 if (l_vec == CEED_VECTOR_ACTIVE) { 931 CeedCallBackend(CeedVectorSetArray(impl->q_vecs_out[out], CEED_MEM_DEVICE, CEED_USE_POINTER, assembled_array)); 932 CeedCallBackend(CeedQFunctionFieldGetSize(qf_output_fields[out], &size)); 933 assembled_array += size * Q * num_elem; // Advance the pointer by the size of the output 934 } 935 } 936 // Apply QFunction 937 CeedCallBackend(CeedQFunctionApply(qf, Q * num_elem, impl->q_vecs_in, impl->q_vecs_out)); 938 } 939 940 // Un-set output q-vecs to prevent accidental overwrite of Assembled 941 for (CeedInt out = 0; out < num_output_fields; out++) { 942 CeedVector l_vec; 943 944 CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[out], &l_vec)); 945 if (l_vec == CEED_VECTOR_ACTIVE) { 946 CeedCallBackend(CeedVectorTakeArray(impl->q_vecs_out[out], CEED_MEM_DEVICE, NULL)); 947 } 948 } 949 950 // Restore input arrays 951 for (CeedInt i = 0; i < num_input_fields; i++) { 952 CeedCallBackend(CeedOperatorInputRestore_Cuda(op_input_fields[i], qf_input_fields[i], i, true, &e_data[i], impl)); 953 } 954 955 // Restore output 956 CeedCallBackend(CeedVectorRestoreArray(*assembled, &assembled_array)); 957 return CEED_ERROR_SUCCESS; 958 } 959 960 //------------------------------------------------------------------------------ 961 // Assemble Linear QFunction 962 //------------------------------------------------------------------------------ 963 static int CeedOperatorLinearAssembleQFunction_Cuda(CeedOperator op, CeedVector *assembled, CeedElemRestriction *rstr, CeedRequest *request) { 964 return CeedOperatorLinearAssembleQFunctionCore_Cuda(op, true, assembled, rstr, request); 965 } 966 967 //------------------------------------------------------------------------------ 968 // Update Assembled Linear QFunction 969 //------------------------------------------------------------------------------ 970 static int CeedOperatorLinearAssembleQFunctionUpdate_Cuda(CeedOperator op, CeedVector assembled, CeedElemRestriction rstr, CeedRequest *request) { 971 return CeedOperatorLinearAssembleQFunctionCore_Cuda(op, false, &assembled, &rstr, request); 972 } 973 974 //------------------------------------------------------------------------------ 975 // Assemble Diagonal Setup 976 //------------------------------------------------------------------------------ 977 static inline int CeedOperatorAssembleDiagonalSetup_Cuda(CeedOperator op) { 978 Ceed ceed; 979 CeedInt num_input_fields, num_output_fields, num_eval_modes_in = 0, num_eval_modes_out = 0; 980 CeedInt q_comp, num_nodes, num_qpts; 981 CeedEvalMode *eval_modes_in = NULL, *eval_modes_out = NULL; 982 CeedBasis basis_in = NULL, basis_out = NULL; 983 CeedQFunctionField *qf_fields; 984 CeedQFunction qf; 985 CeedOperatorField *op_fields; 986 CeedOperator_Cuda *impl; 987 988 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 989 CeedCallBackend(CeedOperatorGetQFunction(op, &qf)); 990 CeedCallBackend(CeedQFunctionGetNumArgs(qf, &num_input_fields, &num_output_fields)); 991 992 // Determine active input basis 993 CeedCallBackend(CeedOperatorGetFields(op, NULL, &op_fields, NULL, NULL)); 994 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_fields, NULL, NULL)); 995 for (CeedInt i = 0; i < num_input_fields; i++) { 996 CeedVector vec; 997 998 CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec)); 999 if (vec == CEED_VECTOR_ACTIVE) { 1000 CeedBasis basis; 1001 CeedEvalMode eval_mode; 1002 1003 CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis)); 1004 CeedCheck(!basis_in || basis_in == basis, ceed, CEED_ERROR_BACKEND, 1005 "Backend does not implement operator diagonal assembly with multiple active bases"); 1006 basis_in = basis; 1007 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode)); 1008 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_in, eval_mode, &q_comp)); 1009 if (eval_mode != CEED_EVAL_WEIGHT) { 1010 // q_comp = 1 if CEED_EVAL_NONE, CEED_EVAL_WEIGHT caught by QF assembly 1011 CeedCallBackend(CeedRealloc(num_eval_modes_in + q_comp, &eval_modes_in)); 1012 for (CeedInt d = 0; d < q_comp; d++) eval_modes_in[num_eval_modes_in + d] = eval_mode; 1013 num_eval_modes_in += q_comp; 1014 } 1015 } 1016 } 1017 1018 // Determine active output basis 1019 CeedCallBackend(CeedOperatorGetFields(op, NULL, NULL, NULL, &op_fields)); 1020 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qf_fields)); 1021 for (CeedInt i = 0; i < num_output_fields; i++) { 1022 CeedVector vec; 1023 1024 CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec)); 1025 if (vec == CEED_VECTOR_ACTIVE) { 1026 CeedBasis basis; 1027 CeedEvalMode eval_mode; 1028 1029 CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis)); 1030 CeedCheck(!basis_out || basis_out == basis, ceed, CEED_ERROR_BACKEND, 1031 "Backend does not implement operator diagonal assembly with multiple active bases"); 1032 basis_out = basis; 1033 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode)); 1034 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_out, eval_mode, &q_comp)); 1035 if (eval_mode != CEED_EVAL_WEIGHT) { 1036 // q_comp = 1 if CEED_EVAL_NONE, CEED_EVAL_WEIGHT caught by QF assembly 1037 CeedCallBackend(CeedRealloc(num_eval_modes_out + q_comp, &eval_modes_out)); 1038 for (CeedInt d = 0; d < q_comp; d++) eval_modes_out[num_eval_modes_out + d] = eval_mode; 1039 num_eval_modes_out += q_comp; 1040 } 1041 } 1042 } 1043 1044 // Operator data struct 1045 CeedCallBackend(CeedOperatorGetData(op, &impl)); 1046 CeedCallBackend(CeedCalloc(1, &impl->diag)); 1047 CeedOperatorDiag_Cuda *diag = impl->diag; 1048 1049 // Basis matrices 1050 CeedCallBackend(CeedBasisGetNumNodes(basis_in, &num_nodes)); 1051 if (basis_in == CEED_BASIS_NONE) num_qpts = num_nodes; 1052 else CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis_in, &num_qpts)); 1053 const CeedInt interp_bytes = num_nodes * num_qpts * sizeof(CeedScalar); 1054 const CeedInt eval_modes_bytes = sizeof(CeedEvalMode); 1055 bool has_eval_none = false; 1056 1057 // CEED_EVAL_NONE 1058 for (CeedInt i = 0; i < num_eval_modes_in; i++) has_eval_none = has_eval_none || (eval_modes_in[i] == CEED_EVAL_NONE); 1059 for (CeedInt i = 0; i < num_eval_modes_out; i++) has_eval_none = has_eval_none || (eval_modes_out[i] == CEED_EVAL_NONE); 1060 if (has_eval_none) { 1061 CeedScalar *identity = NULL; 1062 1063 CeedCallBackend(CeedCalloc(num_nodes * num_qpts, &identity)); 1064 for (CeedInt i = 0; i < (num_nodes < num_qpts ? num_nodes : num_qpts); i++) identity[i * num_nodes + i] = 1.0; 1065 CeedCallCuda(ceed, cudaMalloc((void **)&diag->d_identity, interp_bytes)); 1066 CeedCallCuda(ceed, cudaMemcpy(diag->d_identity, identity, interp_bytes, cudaMemcpyHostToDevice)); 1067 CeedCallBackend(CeedFree(&identity)); 1068 } 1069 1070 // CEED_EVAL_INTERP, CEED_EVAL_GRAD, CEED_EVAL_DIV, and CEED_EVAL_CURL 1071 for (CeedInt in = 0; in < 2; in++) { 1072 CeedFESpace fespace; 1073 CeedBasis basis = in ? basis_in : basis_out; 1074 1075 CeedCallBackend(CeedBasisGetFESpace(basis, &fespace)); 1076 switch (fespace) { 1077 case CEED_FE_SPACE_H1: { 1078 CeedInt q_comp_interp, q_comp_grad; 1079 const CeedScalar *interp, *grad; 1080 CeedScalar *d_interp, *d_grad; 1081 1082 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp)); 1083 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_GRAD, &q_comp_grad)); 1084 1085 CeedCallBackend(CeedBasisGetInterp(basis, &interp)); 1086 CeedCallCuda(ceed, cudaMalloc((void **)&d_interp, interp_bytes * q_comp_interp)); 1087 CeedCallCuda(ceed, cudaMemcpy(d_interp, interp, interp_bytes * q_comp_interp, cudaMemcpyHostToDevice)); 1088 CeedCallBackend(CeedBasisGetGrad(basis, &grad)); 1089 CeedCallCuda(ceed, cudaMalloc((void **)&d_grad, interp_bytes * q_comp_grad)); 1090 CeedCallCuda(ceed, cudaMemcpy(d_grad, grad, interp_bytes * q_comp_grad, cudaMemcpyHostToDevice)); 1091 if (in) { 1092 diag->d_interp_in = d_interp; 1093 diag->d_grad_in = d_grad; 1094 } else { 1095 diag->d_interp_out = d_interp; 1096 diag->d_grad_out = d_grad; 1097 } 1098 } break; 1099 case CEED_FE_SPACE_HDIV: { 1100 CeedInt q_comp_interp, q_comp_div; 1101 const CeedScalar *interp, *div; 1102 CeedScalar *d_interp, *d_div; 1103 1104 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp)); 1105 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_DIV, &q_comp_div)); 1106 1107 CeedCallBackend(CeedBasisGetInterp(basis, &interp)); 1108 CeedCallCuda(ceed, cudaMalloc((void **)&d_interp, interp_bytes * q_comp_interp)); 1109 CeedCallCuda(ceed, cudaMemcpy(d_interp, interp, interp_bytes * q_comp_interp, cudaMemcpyHostToDevice)); 1110 CeedCallBackend(CeedBasisGetDiv(basis, &div)); 1111 CeedCallCuda(ceed, cudaMalloc((void **)&d_div, interp_bytes * q_comp_div)); 1112 CeedCallCuda(ceed, cudaMemcpy(d_div, div, interp_bytes * q_comp_div, cudaMemcpyHostToDevice)); 1113 if (in) { 1114 diag->d_interp_in = d_interp; 1115 diag->d_div_in = d_div; 1116 } else { 1117 diag->d_interp_out = d_interp; 1118 diag->d_div_out = d_div; 1119 } 1120 } break; 1121 case CEED_FE_SPACE_HCURL: { 1122 CeedInt q_comp_interp, q_comp_curl; 1123 const CeedScalar *interp, *curl; 1124 CeedScalar *d_interp, *d_curl; 1125 1126 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp)); 1127 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_CURL, &q_comp_curl)); 1128 1129 CeedCallBackend(CeedBasisGetInterp(basis, &interp)); 1130 CeedCallCuda(ceed, cudaMalloc((void **)&d_interp, interp_bytes * q_comp_interp)); 1131 CeedCallCuda(ceed, cudaMemcpy(d_interp, interp, interp_bytes * q_comp_interp, cudaMemcpyHostToDevice)); 1132 CeedCallBackend(CeedBasisGetCurl(basis, &curl)); 1133 CeedCallCuda(ceed, cudaMalloc((void **)&d_curl, interp_bytes * q_comp_curl)); 1134 CeedCallCuda(ceed, cudaMemcpy(d_curl, curl, interp_bytes * q_comp_curl, cudaMemcpyHostToDevice)); 1135 if (in) { 1136 diag->d_interp_in = d_interp; 1137 diag->d_curl_in = d_curl; 1138 } else { 1139 diag->d_interp_out = d_interp; 1140 diag->d_curl_out = d_curl; 1141 } 1142 } break; 1143 } 1144 } 1145 1146 // Arrays of eval_modes 1147 CeedCallCuda(ceed, cudaMalloc((void **)&diag->d_eval_modes_in, num_eval_modes_in * eval_modes_bytes)); 1148 CeedCallCuda(ceed, cudaMemcpy(diag->d_eval_modes_in, eval_modes_in, num_eval_modes_in * eval_modes_bytes, cudaMemcpyHostToDevice)); 1149 CeedCallCuda(ceed, cudaMalloc((void **)&diag->d_eval_modes_out, num_eval_modes_out * eval_modes_bytes)); 1150 CeedCallCuda(ceed, cudaMemcpy(diag->d_eval_modes_out, eval_modes_out, num_eval_modes_out * eval_modes_bytes, cudaMemcpyHostToDevice)); 1151 CeedCallBackend(CeedFree(&eval_modes_in)); 1152 CeedCallBackend(CeedFree(&eval_modes_out)); 1153 return CEED_ERROR_SUCCESS; 1154 } 1155 1156 //------------------------------------------------------------------------------ 1157 // Assemble Diagonal Setup (Compilation) 1158 //------------------------------------------------------------------------------ 1159 static inline int CeedOperatorAssembleDiagonalSetupCompile_Cuda(CeedOperator op, CeedInt use_ceedsize_idx, const bool is_point_block) { 1160 Ceed ceed; 1161 char *diagonal_kernel_source; 1162 const char *diagonal_kernel_path; 1163 CeedInt num_input_fields, num_output_fields, num_eval_modes_in = 0, num_eval_modes_out = 0; 1164 CeedInt num_comp, q_comp, num_nodes, num_qpts; 1165 CeedBasis basis_in = NULL, basis_out = NULL; 1166 CeedQFunctionField *qf_fields; 1167 CeedQFunction qf; 1168 CeedOperatorField *op_fields; 1169 CeedOperator_Cuda *impl; 1170 1171 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 1172 CeedCallBackend(CeedOperatorGetQFunction(op, &qf)); 1173 CeedCallBackend(CeedQFunctionGetNumArgs(qf, &num_input_fields, &num_output_fields)); 1174 1175 // Determine active input basis 1176 CeedCallBackend(CeedOperatorGetFields(op, NULL, &op_fields, NULL, NULL)); 1177 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_fields, NULL, NULL)); 1178 for (CeedInt i = 0; i < num_input_fields; i++) { 1179 CeedVector vec; 1180 1181 CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec)); 1182 if (vec == CEED_VECTOR_ACTIVE) { 1183 CeedEvalMode eval_mode; 1184 1185 CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis_in)); 1186 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode)); 1187 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_in, eval_mode, &q_comp)); 1188 if (eval_mode != CEED_EVAL_WEIGHT) { 1189 num_eval_modes_in += q_comp; 1190 } 1191 } 1192 } 1193 1194 // Determine active output basis 1195 CeedCallBackend(CeedOperatorGetFields(op, NULL, NULL, NULL, &op_fields)); 1196 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qf_fields)); 1197 for (CeedInt i = 0; i < num_output_fields; i++) { 1198 CeedVector vec; 1199 1200 CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec)); 1201 if (vec == CEED_VECTOR_ACTIVE) { 1202 CeedEvalMode eval_mode; 1203 1204 CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis_out)); 1205 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode)); 1206 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_out, eval_mode, &q_comp)); 1207 if (eval_mode != CEED_EVAL_WEIGHT) { 1208 num_eval_modes_out += q_comp; 1209 } 1210 } 1211 } 1212 1213 // Operator data struct 1214 CeedCallBackend(CeedOperatorGetData(op, &impl)); 1215 CeedOperatorDiag_Cuda *diag = impl->diag; 1216 1217 // Assemble kernel 1218 CUmodule *module = is_point_block ? &diag->module_point_block : &diag->module; 1219 CeedInt elems_per_block = 1; 1220 CeedCallBackend(CeedBasisGetNumNodes(basis_in, &num_nodes)); 1221 CeedCallBackend(CeedBasisGetNumComponents(basis_in, &num_comp)); 1222 if (basis_in == CEED_BASIS_NONE) num_qpts = num_nodes; 1223 else CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis_in, &num_qpts)); 1224 CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-operator-assemble-diagonal.h", &diagonal_kernel_path)); 1225 CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Diagonal Assembly Kernel Source -----\n"); 1226 CeedCallBackend(CeedLoadSourceToBuffer(ceed, diagonal_kernel_path, &diagonal_kernel_source)); 1227 CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Diagonal Assembly Source Complete! -----\n"); 1228 CeedCallCuda(ceed, CeedCompile_Cuda(ceed, diagonal_kernel_source, module, 8, "NUM_EVAL_MODES_IN", num_eval_modes_in, "NUM_EVAL_MODES_OUT", 1229 num_eval_modes_out, "NUM_COMP", num_comp, "NUM_NODES", num_nodes, "NUM_QPTS", num_qpts, "USE_CEEDSIZE", 1230 use_ceedsize_idx, "USE_POINT_BLOCK", is_point_block ? 1 : 0, "BLOCK_SIZE", num_nodes * elems_per_block)); 1231 CeedCallCuda(ceed, CeedGetKernel_Cuda(ceed, *module, "LinearDiagonal", is_point_block ? &diag->LinearPointBlock : &diag->LinearDiagonal)); 1232 CeedCallBackend(CeedFree(&diagonal_kernel_path)); 1233 CeedCallBackend(CeedFree(&diagonal_kernel_source)); 1234 return CEED_ERROR_SUCCESS; 1235 } 1236 1237 //------------------------------------------------------------------------------ 1238 // Assemble Diagonal Core 1239 //------------------------------------------------------------------------------ 1240 static inline int CeedOperatorAssembleDiagonalCore_Cuda(CeedOperator op, CeedVector assembled, CeedRequest *request, const bool is_point_block) { 1241 Ceed ceed; 1242 CeedInt num_elem, num_nodes; 1243 CeedScalar *elem_diag_array; 1244 const CeedScalar *assembled_qf_array; 1245 CeedVector assembled_qf = NULL, elem_diag; 1246 CeedElemRestriction assembled_rstr = NULL, rstr_in, rstr_out, diag_rstr; 1247 CeedOperator_Cuda *impl; 1248 1249 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 1250 CeedCallBackend(CeedOperatorGetData(op, &impl)); 1251 1252 // Assemble QFunction 1253 CeedCallBackend(CeedOperatorLinearAssembleQFunctionBuildOrUpdate(op, &assembled_qf, &assembled_rstr, request)); 1254 CeedCallBackend(CeedElemRestrictionDestroy(&assembled_rstr)); 1255 CeedCallBackend(CeedVectorGetArrayRead(assembled_qf, CEED_MEM_DEVICE, &assembled_qf_array)); 1256 1257 // Setup 1258 if (!impl->diag) CeedCallBackend(CeedOperatorAssembleDiagonalSetup_Cuda(op)); 1259 CeedOperatorDiag_Cuda *diag = impl->diag; 1260 1261 assert(diag != NULL); 1262 1263 // Assemble kernel if needed 1264 if ((!is_point_block && !diag->LinearDiagonal) || (is_point_block && !diag->LinearPointBlock)) { 1265 CeedSize assembled_length, assembled_qf_length; 1266 CeedInt use_ceedsize_idx = 0; 1267 CeedCallBackend(CeedVectorGetLength(assembled, &assembled_length)); 1268 CeedCallBackend(CeedVectorGetLength(assembled_qf, &assembled_qf_length)); 1269 if ((assembled_length > INT_MAX) || (assembled_qf_length > INT_MAX)) use_ceedsize_idx = 1; 1270 1271 CeedCallBackend(CeedOperatorAssembleDiagonalSetupCompile_Cuda(op, use_ceedsize_idx, is_point_block)); 1272 } 1273 1274 // Restriction and diagonal vector 1275 CeedCallBackend(CeedOperatorGetActiveElemRestrictions(op, &rstr_in, &rstr_out)); 1276 CeedCheck(rstr_in == rstr_out, ceed, CEED_ERROR_BACKEND, 1277 "Cannot assemble operator diagonal with different input and output active element restrictions"); 1278 if (!is_point_block && !diag->diag_rstr) { 1279 CeedCallBackend(CeedElemRestrictionCreateUnsignedCopy(rstr_out, &diag->diag_rstr)); 1280 CeedCallBackend(CeedElemRestrictionCreateVector(diag->diag_rstr, NULL, &diag->elem_diag)); 1281 } else if (is_point_block && !diag->point_block_diag_rstr) { 1282 CeedCallBackend(CeedOperatorCreateActivePointBlockRestriction(rstr_out, &diag->point_block_diag_rstr)); 1283 CeedCallBackend(CeedElemRestrictionCreateVector(diag->point_block_diag_rstr, NULL, &diag->point_block_elem_diag)); 1284 } 1285 diag_rstr = is_point_block ? diag->point_block_diag_rstr : diag->diag_rstr; 1286 elem_diag = is_point_block ? diag->point_block_elem_diag : diag->elem_diag; 1287 CeedCallBackend(CeedVectorSetValue(elem_diag, 0.0)); 1288 1289 // Only assemble diagonal if the basis has nodes, otherwise inputs are null pointers 1290 CeedCallBackend(CeedElemRestrictionGetElementSize(diag_rstr, &num_nodes)); 1291 if (num_nodes > 0) { 1292 // Assemble element operator diagonals 1293 CeedCallBackend(CeedElemRestrictionGetNumElements(diag_rstr, &num_elem)); 1294 CeedCallBackend(CeedVectorGetArray(elem_diag, CEED_MEM_DEVICE, &elem_diag_array)); 1295 1296 // Compute the diagonal of B^T D B 1297 CeedInt elems_per_block = 1; 1298 CeedInt grid = CeedDivUpInt(num_elem, elems_per_block); 1299 void *args[] = {(void *)&num_elem, &diag->d_identity, &diag->d_interp_in, &diag->d_grad_in, &diag->d_div_in, 1300 &diag->d_curl_in, &diag->d_interp_out, &diag->d_grad_out, &diag->d_div_out, &diag->d_curl_out, 1301 &diag->d_eval_modes_in, &diag->d_eval_modes_out, &assembled_qf_array, &elem_diag_array}; 1302 1303 if (is_point_block) { 1304 CeedCallBackend(CeedRunKernelDim_Cuda(ceed, diag->LinearPointBlock, grid, num_nodes, 1, elems_per_block, args)); 1305 } else { 1306 CeedCallBackend(CeedRunKernelDim_Cuda(ceed, diag->LinearDiagonal, grid, num_nodes, 1, elems_per_block, args)); 1307 } 1308 1309 // Restore arrays 1310 CeedCallBackend(CeedVectorRestoreArray(elem_diag, &elem_diag_array)); 1311 CeedCallBackend(CeedVectorRestoreArrayRead(assembled_qf, &assembled_qf_array)); 1312 } 1313 1314 // Assemble local operator diagonal 1315 CeedCallBackend(CeedElemRestrictionApply(diag_rstr, CEED_TRANSPOSE, elem_diag, assembled, request)); 1316 1317 // Cleanup 1318 CeedCallBackend(CeedVectorDestroy(&assembled_qf)); 1319 return CEED_ERROR_SUCCESS; 1320 } 1321 1322 //------------------------------------------------------------------------------ 1323 // Assemble Linear Diagonal 1324 //------------------------------------------------------------------------------ 1325 static int CeedOperatorLinearAssembleAddDiagonal_Cuda(CeedOperator op, CeedVector assembled, CeedRequest *request) { 1326 CeedCallBackend(CeedOperatorAssembleDiagonalCore_Cuda(op, assembled, request, false)); 1327 return CEED_ERROR_SUCCESS; 1328 } 1329 1330 //------------------------------------------------------------------------------ 1331 // Assemble Linear Point Block Diagonal 1332 //------------------------------------------------------------------------------ 1333 static int CeedOperatorLinearAssembleAddPointBlockDiagonal_Cuda(CeedOperator op, CeedVector assembled, CeedRequest *request) { 1334 CeedCallBackend(CeedOperatorAssembleDiagonalCore_Cuda(op, assembled, request, true)); 1335 return CEED_ERROR_SUCCESS; 1336 } 1337 1338 //------------------------------------------------------------------------------ 1339 // Single Operator Assembly Setup 1340 //------------------------------------------------------------------------------ 1341 static int CeedSingleOperatorAssembleSetup_Cuda(CeedOperator op, CeedInt use_ceedsize_idx) { 1342 Ceed ceed; 1343 Ceed_Cuda *cuda_data; 1344 char *assembly_kernel_source; 1345 const char *assembly_kernel_path; 1346 CeedInt num_input_fields, num_output_fields, num_eval_modes_in = 0, num_eval_modes_out = 0; 1347 CeedInt elem_size_in, num_qpts_in = 0, num_comp_in, elem_size_out, num_qpts_out, num_comp_out, q_comp; 1348 CeedEvalMode *eval_modes_in = NULL, *eval_modes_out = NULL; 1349 CeedElemRestriction rstr_in = NULL, rstr_out = NULL; 1350 CeedBasis basis_in = NULL, basis_out = NULL; 1351 CeedQFunctionField *qf_fields; 1352 CeedQFunction qf; 1353 CeedOperatorField *input_fields, *output_fields; 1354 CeedOperator_Cuda *impl; 1355 1356 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 1357 CeedCallBackend(CeedOperatorGetData(op, &impl)); 1358 1359 // Get intput and output fields 1360 CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &input_fields, &num_output_fields, &output_fields)); 1361 1362 // Determine active input basis eval mode 1363 CeedCallBackend(CeedOperatorGetQFunction(op, &qf)); 1364 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_fields, NULL, NULL)); 1365 for (CeedInt i = 0; i < num_input_fields; i++) { 1366 CeedVector vec; 1367 1368 CeedCallBackend(CeedOperatorFieldGetVector(input_fields[i], &vec)); 1369 if (vec == CEED_VECTOR_ACTIVE) { 1370 CeedBasis basis; 1371 CeedEvalMode eval_mode; 1372 1373 CeedCallBackend(CeedOperatorFieldGetBasis(input_fields[i], &basis)); 1374 CeedCheck(!basis_in || basis_in == basis, ceed, CEED_ERROR_BACKEND, "Backend does not implement operator assembly with multiple active bases"); 1375 basis_in = basis; 1376 CeedCallBackend(CeedOperatorFieldGetElemRestriction(input_fields[i], &rstr_in)); 1377 CeedCallBackend(CeedElemRestrictionGetElementSize(rstr_in, &elem_size_in)); 1378 if (basis_in == CEED_BASIS_NONE) num_qpts_in = elem_size_in; 1379 else CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis_in, &num_qpts_in)); 1380 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode)); 1381 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_in, eval_mode, &q_comp)); 1382 if (eval_mode != CEED_EVAL_WEIGHT) { 1383 // q_comp = 1 if CEED_EVAL_NONE, CEED_EVAL_WEIGHT caught by QF Assembly 1384 CeedCallBackend(CeedRealloc(num_eval_modes_in + q_comp, &eval_modes_in)); 1385 for (CeedInt d = 0; d < q_comp; d++) { 1386 eval_modes_in[num_eval_modes_in + d] = eval_mode; 1387 } 1388 num_eval_modes_in += q_comp; 1389 } 1390 } 1391 } 1392 1393 // Determine active output basis; basis_out and rstr_out only used if same as input, TODO 1394 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qf_fields)); 1395 for (CeedInt i = 0; i < num_output_fields; i++) { 1396 CeedVector vec; 1397 1398 CeedCallBackend(CeedOperatorFieldGetVector(output_fields[i], &vec)); 1399 if (vec == CEED_VECTOR_ACTIVE) { 1400 CeedBasis basis; 1401 CeedEvalMode eval_mode; 1402 1403 CeedCallBackend(CeedOperatorFieldGetBasis(output_fields[i], &basis)); 1404 CeedCheck(!basis_out || basis_out == basis, ceed, CEED_ERROR_BACKEND, 1405 "Backend does not implement operator assembly with multiple active bases"); 1406 basis_out = basis; 1407 CeedCallBackend(CeedOperatorFieldGetElemRestriction(output_fields[i], &rstr_out)); 1408 CeedCallBackend(CeedElemRestrictionGetElementSize(rstr_out, &elem_size_out)); 1409 if (basis_out == CEED_BASIS_NONE) num_qpts_out = elem_size_out; 1410 else CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis_out, &num_qpts_out)); 1411 CeedCheck(num_qpts_in == num_qpts_out, ceed, CEED_ERROR_UNSUPPORTED, 1412 "Active input and output bases must have the same number of quadrature points"); 1413 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode)); 1414 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_out, eval_mode, &q_comp)); 1415 if (eval_mode != CEED_EVAL_WEIGHT) { 1416 // q_comp = 1 if CEED_EVAL_NONE, CEED_EVAL_WEIGHT caught by QF Assembly 1417 CeedCallBackend(CeedRealloc(num_eval_modes_out + q_comp, &eval_modes_out)); 1418 for (CeedInt d = 0; d < q_comp; d++) { 1419 eval_modes_out[num_eval_modes_out + d] = eval_mode; 1420 } 1421 num_eval_modes_out += q_comp; 1422 } 1423 } 1424 } 1425 CeedCheck(num_eval_modes_in > 0 && num_eval_modes_out > 0, ceed, CEED_ERROR_UNSUPPORTED, "Cannot assemble operator without inputs/outputs"); 1426 1427 CeedCallBackend(CeedCalloc(1, &impl->asmb)); 1428 CeedOperatorAssemble_Cuda *asmb = impl->asmb; 1429 asmb->elems_per_block = 1; 1430 asmb->block_size_x = elem_size_in; 1431 asmb->block_size_y = elem_size_out; 1432 1433 CeedCallBackend(CeedGetData(ceed, &cuda_data)); 1434 bool fallback = asmb->block_size_x * asmb->block_size_y * asmb->elems_per_block > cuda_data->device_prop.maxThreadsPerBlock; 1435 1436 if (fallback) { 1437 // Use fallback kernel with 1D threadblock 1438 asmb->block_size_y = 1; 1439 } 1440 1441 // Compile kernels 1442 CeedCallBackend(CeedElemRestrictionGetNumComponents(rstr_in, &num_comp_in)); 1443 CeedCallBackend(CeedElemRestrictionGetNumComponents(rstr_out, &num_comp_out)); 1444 CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-operator-assemble.h", &assembly_kernel_path)); 1445 CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Assembly Kernel Source -----\n"); 1446 CeedCallBackend(CeedLoadSourceToBuffer(ceed, assembly_kernel_path, &assembly_kernel_source)); 1447 CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Assembly Source Complete! -----\n"); 1448 CeedCallBackend(CeedCompile_Cuda(ceed, assembly_kernel_source, &asmb->module, 10, "NUM_EVAL_MODES_IN", num_eval_modes_in, "NUM_EVAL_MODES_OUT", 1449 num_eval_modes_out, "NUM_COMP_IN", num_comp_in, "NUM_COMP_OUT", num_comp_out, "NUM_NODES_IN", elem_size_in, 1450 "NUM_NODES_OUT", elem_size_out, "NUM_QPTS", num_qpts_in, "BLOCK_SIZE", 1451 asmb->block_size_x * asmb->block_size_y * asmb->elems_per_block, "BLOCK_SIZE_Y", asmb->block_size_y, 1452 "USE_CEEDSIZE", use_ceedsize_idx)); 1453 CeedCallBackend(CeedGetKernel_Cuda(ceed, asmb->module, "LinearAssemble", &asmb->LinearAssemble)); 1454 CeedCallBackend(CeedFree(&assembly_kernel_path)); 1455 CeedCallBackend(CeedFree(&assembly_kernel_source)); 1456 1457 // Load into B_in, in order that they will be used in eval_modes_in 1458 { 1459 const CeedInt in_bytes = elem_size_in * num_qpts_in * num_eval_modes_in * sizeof(CeedScalar); 1460 CeedInt d_in = 0; 1461 CeedEvalMode eval_modes_in_prev = CEED_EVAL_NONE; 1462 bool has_eval_none = false; 1463 CeedScalar *identity = NULL; 1464 1465 for (CeedInt i = 0; i < num_eval_modes_in; i++) { 1466 has_eval_none = has_eval_none || (eval_modes_in[i] == CEED_EVAL_NONE); 1467 } 1468 if (has_eval_none) { 1469 CeedCallBackend(CeedCalloc(elem_size_in * num_qpts_in, &identity)); 1470 for (CeedInt i = 0; i < (elem_size_in < num_qpts_in ? elem_size_in : num_qpts_in); i++) identity[i * elem_size_in + i] = 1.0; 1471 } 1472 1473 CeedCallCuda(ceed, cudaMalloc((void **)&asmb->d_B_in, in_bytes)); 1474 for (CeedInt i = 0; i < num_eval_modes_in; i++) { 1475 const CeedScalar *h_B_in; 1476 1477 CeedCallBackend(CeedOperatorGetBasisPointer(basis_in, eval_modes_in[i], identity, &h_B_in)); 1478 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_in, eval_modes_in[i], &q_comp)); 1479 if (q_comp > 1) { 1480 if (i == 0 || eval_modes_in[i] != eval_modes_in_prev) d_in = 0; 1481 else h_B_in = &h_B_in[(++d_in) * elem_size_in * num_qpts_in]; 1482 } 1483 eval_modes_in_prev = eval_modes_in[i]; 1484 1485 CeedCallCuda(ceed, cudaMemcpy(&asmb->d_B_in[i * elem_size_in * num_qpts_in], h_B_in, elem_size_in * num_qpts_in * sizeof(CeedScalar), 1486 cudaMemcpyHostToDevice)); 1487 } 1488 1489 if (identity) { 1490 CeedCallBackend(CeedFree(&identity)); 1491 } 1492 } 1493 1494 // Load into B_out, in order that they will be used in eval_modes_out 1495 { 1496 const CeedInt out_bytes = elem_size_out * num_qpts_out * num_eval_modes_out * sizeof(CeedScalar); 1497 CeedInt d_out = 0; 1498 CeedEvalMode eval_modes_out_prev = CEED_EVAL_NONE; 1499 bool has_eval_none = false; 1500 CeedScalar *identity = NULL; 1501 1502 for (CeedInt i = 0; i < num_eval_modes_out; i++) { 1503 has_eval_none = has_eval_none || (eval_modes_out[i] == CEED_EVAL_NONE); 1504 } 1505 if (has_eval_none) { 1506 CeedCallBackend(CeedCalloc(elem_size_out * num_qpts_out, &identity)); 1507 for (CeedInt i = 0; i < (elem_size_out < num_qpts_out ? elem_size_out : num_qpts_out); i++) identity[i * elem_size_out + i] = 1.0; 1508 } 1509 1510 CeedCallCuda(ceed, cudaMalloc((void **)&asmb->d_B_out, out_bytes)); 1511 for (CeedInt i = 0; i < num_eval_modes_out; i++) { 1512 const CeedScalar *h_B_out; 1513 1514 CeedCallBackend(CeedOperatorGetBasisPointer(basis_out, eval_modes_out[i], identity, &h_B_out)); 1515 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_out, eval_modes_out[i], &q_comp)); 1516 if (q_comp > 1) { 1517 if (i == 0 || eval_modes_out[i] != eval_modes_out_prev) d_out = 0; 1518 else h_B_out = &h_B_out[(++d_out) * elem_size_out * num_qpts_out]; 1519 } 1520 eval_modes_out_prev = eval_modes_out[i]; 1521 1522 CeedCallCuda(ceed, cudaMemcpy(&asmb->d_B_out[i * elem_size_out * num_qpts_out], h_B_out, elem_size_out * num_qpts_out * sizeof(CeedScalar), 1523 cudaMemcpyHostToDevice)); 1524 } 1525 1526 if (identity) { 1527 CeedCallBackend(CeedFree(&identity)); 1528 } 1529 } 1530 return CEED_ERROR_SUCCESS; 1531 } 1532 1533 //------------------------------------------------------------------------------ 1534 // Assemble matrix data for COO matrix of assembled operator. 1535 // The sparsity pattern is set by CeedOperatorLinearAssembleSymbolic. 1536 // 1537 // Note that this (and other assembly routines) currently assume only one active input restriction/basis per operator 1538 // (could have multiple basis eval modes). 1539 // TODO: allow multiple active input restrictions/basis objects 1540 //------------------------------------------------------------------------------ 1541 static int CeedSingleOperatorAssemble_Cuda(CeedOperator op, CeedInt offset, CeedVector values) { 1542 Ceed ceed; 1543 CeedSize values_length = 0, assembled_qf_length = 0; 1544 CeedInt use_ceedsize_idx = 0, num_elem_in, num_elem_out, elem_size_in, elem_size_out; 1545 CeedScalar *values_array; 1546 const CeedScalar *assembled_qf_array; 1547 CeedVector assembled_qf = NULL; 1548 CeedElemRestriction assembled_rstr = NULL, rstr_in, rstr_out; 1549 CeedRestrictionType rstr_type_in, rstr_type_out; 1550 const bool *orients_in = NULL, *orients_out = NULL; 1551 const CeedInt8 *curl_orients_in = NULL, *curl_orients_out = NULL; 1552 CeedOperator_Cuda *impl; 1553 1554 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 1555 CeedCallBackend(CeedOperatorGetData(op, &impl)); 1556 1557 // Assemble QFunction 1558 CeedCallBackend(CeedOperatorLinearAssembleQFunctionBuildOrUpdate(op, &assembled_qf, &assembled_rstr, CEED_REQUEST_IMMEDIATE)); 1559 CeedCallBackend(CeedElemRestrictionDestroy(&assembled_rstr)); 1560 CeedCallBackend(CeedVectorGetArrayRead(assembled_qf, CEED_MEM_DEVICE, &assembled_qf_array)); 1561 1562 CeedCallBackend(CeedVectorGetLength(values, &values_length)); 1563 CeedCallBackend(CeedVectorGetLength(assembled_qf, &assembled_qf_length)); 1564 if ((values_length > INT_MAX) || (assembled_qf_length > INT_MAX)) use_ceedsize_idx = 1; 1565 1566 // Setup 1567 if (!impl->asmb) CeedCallBackend(CeedSingleOperatorAssembleSetup_Cuda(op, use_ceedsize_idx)); 1568 CeedOperatorAssemble_Cuda *asmb = impl->asmb; 1569 1570 assert(asmb != NULL); 1571 1572 // Assemble element operator 1573 CeedCallBackend(CeedVectorGetArray(values, CEED_MEM_DEVICE, &values_array)); 1574 values_array += offset; 1575 1576 CeedCallBackend(CeedOperatorGetActiveElemRestrictions(op, &rstr_in, &rstr_out)); 1577 CeedCallBackend(CeedElemRestrictionGetNumElements(rstr_in, &num_elem_in)); 1578 CeedCallBackend(CeedElemRestrictionGetElementSize(rstr_in, &elem_size_in)); 1579 1580 CeedCallBackend(CeedElemRestrictionGetType(rstr_in, &rstr_type_in)); 1581 if (rstr_type_in == CEED_RESTRICTION_ORIENTED) { 1582 CeedCallBackend(CeedElemRestrictionGetOrientations(rstr_in, CEED_MEM_DEVICE, &orients_in)); 1583 } else if (rstr_type_in == CEED_RESTRICTION_CURL_ORIENTED) { 1584 CeedCallBackend(CeedElemRestrictionGetCurlOrientations(rstr_in, CEED_MEM_DEVICE, &curl_orients_in)); 1585 } 1586 1587 if (rstr_in != rstr_out) { 1588 CeedCallBackend(CeedElemRestrictionGetNumElements(rstr_out, &num_elem_out)); 1589 CeedCheck(num_elem_in == num_elem_out, ceed, CEED_ERROR_UNSUPPORTED, 1590 "Active input and output operator restrictions must have the same number of elements"); 1591 CeedCallBackend(CeedElemRestrictionGetElementSize(rstr_out, &elem_size_out)); 1592 1593 CeedCallBackend(CeedElemRestrictionGetType(rstr_out, &rstr_type_out)); 1594 if (rstr_type_out == CEED_RESTRICTION_ORIENTED) { 1595 CeedCallBackend(CeedElemRestrictionGetOrientations(rstr_out, CEED_MEM_DEVICE, &orients_out)); 1596 } else if (rstr_type_out == CEED_RESTRICTION_CURL_ORIENTED) { 1597 CeedCallBackend(CeedElemRestrictionGetCurlOrientations(rstr_out, CEED_MEM_DEVICE, &curl_orients_out)); 1598 } 1599 } else { 1600 elem_size_out = elem_size_in; 1601 orients_out = orients_in; 1602 curl_orients_out = curl_orients_in; 1603 } 1604 1605 // Compute B^T D B 1606 CeedInt shared_mem = 1607 ((curl_orients_in || curl_orients_out ? elem_size_in * elem_size_out : 0) + (curl_orients_in ? elem_size_in * asmb->block_size_y : 0)) * 1608 sizeof(CeedScalar); 1609 CeedInt grid = CeedDivUpInt(num_elem_in, asmb->elems_per_block); 1610 void *args[] = {(void *)&num_elem_in, &asmb->d_B_in, &asmb->d_B_out, &orients_in, &curl_orients_in, 1611 &orients_out, &curl_orients_out, &assembled_qf_array, &values_array}; 1612 1613 CeedCallBackend( 1614 CeedRunKernelDimShared_Cuda(ceed, asmb->LinearAssemble, grid, asmb->block_size_x, asmb->block_size_y, asmb->elems_per_block, shared_mem, args)); 1615 1616 // Restore arrays 1617 CeedCallBackend(CeedVectorRestoreArray(values, &values_array)); 1618 CeedCallBackend(CeedVectorRestoreArrayRead(assembled_qf, &assembled_qf_array)); 1619 1620 // Cleanup 1621 CeedCallBackend(CeedVectorDestroy(&assembled_qf)); 1622 if (rstr_type_in == CEED_RESTRICTION_ORIENTED) { 1623 CeedCallBackend(CeedElemRestrictionRestoreOrientations(rstr_in, &orients_in)); 1624 } else if (rstr_type_in == CEED_RESTRICTION_CURL_ORIENTED) { 1625 CeedCallBackend(CeedElemRestrictionRestoreCurlOrientations(rstr_in, &curl_orients_in)); 1626 } 1627 if (rstr_in != rstr_out) { 1628 if (rstr_type_out == CEED_RESTRICTION_ORIENTED) { 1629 CeedCallBackend(CeedElemRestrictionRestoreOrientations(rstr_out, &orients_out)); 1630 } else if (rstr_type_out == CEED_RESTRICTION_CURL_ORIENTED) { 1631 CeedCallBackend(CeedElemRestrictionRestoreCurlOrientations(rstr_out, &curl_orients_out)); 1632 } 1633 } 1634 return CEED_ERROR_SUCCESS; 1635 } 1636 1637 //------------------------------------------------------------------------------ 1638 // Assemble Linear QFunction AtPoints 1639 //------------------------------------------------------------------------------ 1640 static int CeedOperatorLinearAssembleQFunctionAtPoints_Cuda(CeedOperator op, CeedVector *assembled, CeedElemRestriction *rstr, CeedRequest *request) { 1641 return CeedError(CeedOperatorReturnCeed(op), CEED_ERROR_BACKEND, "Backend does not implement CeedOperatorLinearAssembleQFunction"); 1642 } 1643 1644 //------------------------------------------------------------------------------ 1645 // Assemble Linear Diagonal AtPoints 1646 //------------------------------------------------------------------------------ 1647 static int CeedOperatorLinearAssembleAddDiagonalAtPoints_Cuda(CeedOperator op, CeedVector assembled, CeedRequest *request) { 1648 CeedInt max_num_points, *num_points, num_elem, num_input_fields, num_output_fields; 1649 CeedScalar *e_data_in[CEED_FIELD_MAX] = {NULL}, *e_data_out[CEED_FIELD_MAX] = {NULL}; 1650 CeedQFunctionField *qf_input_fields, *qf_output_fields; 1651 CeedQFunction qf; 1652 CeedOperatorField *op_input_fields, *op_output_fields; 1653 CeedOperator_Cuda *impl; 1654 1655 CeedCallBackend(CeedOperatorGetData(op, &impl)); 1656 CeedCallBackend(CeedOperatorGetQFunction(op, &qf)); 1657 CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem)); 1658 CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields)); 1659 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields)); 1660 1661 // Setup 1662 CeedCallBackend(CeedOperatorSetupAtPoints_Cuda(op)); 1663 num_points = impl->num_points; 1664 max_num_points = impl->max_num_points; 1665 1666 // Create separate output e-vecs 1667 if (impl->has_shared_e_vecs) { 1668 for (CeedInt i = 0; i < impl->num_outputs; i++) { 1669 CeedCallBackend(CeedVectorDestroy(&impl->q_vecs_out[i])); 1670 CeedCallBackend(CeedVectorDestroy(&impl->e_vecs_out[i])); 1671 } 1672 CeedCallBackend(CeedOperatorSetupFields_Cuda(qf, op, false, true, impl->skip_rstr_out, impl->apply_add_basis_out, impl->e_vecs_out, 1673 impl->q_vecs_out, num_output_fields, max_num_points, num_elem)); 1674 } 1675 impl->has_shared_e_vecs = false; 1676 1677 // Get point coordinates 1678 if (!impl->point_coords_elem) { 1679 CeedVector point_coords = NULL; 1680 CeedElemRestriction rstr_points = NULL; 1681 1682 CeedCallBackend(CeedOperatorAtPointsGetPoints(op, &rstr_points, &point_coords)); 1683 CeedCallBackend(CeedElemRestrictionCreateVector(rstr_points, NULL, &impl->point_coords_elem)); 1684 CeedCallBackend(CeedElemRestrictionApply(rstr_points, CEED_NOTRANSPOSE, point_coords, impl->point_coords_elem, request)); 1685 } 1686 1687 // Process inputs 1688 for (CeedInt i = 0; i < num_input_fields; i++) { 1689 CeedCallBackend(CeedOperatorInputRestrict_Cuda(op_input_fields[i], qf_input_fields[i], i, NULL, true, &e_data_in[i], impl, request)); 1690 CeedCallBackend(CeedOperatorInputBasisAtPoints_Cuda(op_input_fields[i], qf_input_fields[i], i, num_elem, num_points, true, e_data_in[i], impl)); 1691 } 1692 1693 // Clear active input Qvecs 1694 for (CeedInt i = 0; i < num_input_fields; i++) { 1695 CeedVector vec; 1696 1697 CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec)); 1698 if (vec != CEED_VECTOR_ACTIVE) continue; 1699 CeedCallBackend(CeedVectorSetValue(impl->q_vecs_in[i], 0.0)); 1700 } 1701 1702 // Output pointers, as necessary 1703 for (CeedInt i = 0; i < num_output_fields; i++) { 1704 CeedEvalMode eval_mode; 1705 1706 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 1707 if (eval_mode == CEED_EVAL_NONE) { 1708 // Set the output Q-Vector to use the E-Vector data directly. 1709 CeedCallBackend(CeedVectorGetArrayWrite(impl->e_vecs_out[i], CEED_MEM_DEVICE, &e_data_out[i])); 1710 CeedCallBackend(CeedVectorSetArray(impl->q_vecs_out[i], CEED_MEM_DEVICE, CEED_USE_POINTER, e_data_out[i])); 1711 } 1712 } 1713 1714 // Loop over active fields 1715 for (CeedInt i = 0; i < num_input_fields; i++) { 1716 bool is_active_at_points = true; 1717 CeedInt elem_size = 1, num_comp_active = 1, e_vec_size = 0; 1718 CeedRestrictionType rstr_type; 1719 CeedVector l_vec; 1720 CeedElemRestriction elem_rstr; 1721 1722 CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &l_vec)); 1723 // -- Skip non-active input 1724 if (l_vec != CEED_VECTOR_ACTIVE) continue; 1725 1726 // -- Get active restriction type 1727 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr)); 1728 CeedCallBackend(CeedElemRestrictionGetType(elem_rstr, &rstr_type)); 1729 is_active_at_points = rstr_type == CEED_RESTRICTION_POINTS; 1730 if (!is_active_at_points) CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size)); 1731 else elem_size = max_num_points; 1732 CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp_active)); 1733 1734 e_vec_size = elem_size * num_comp_active; 1735 for (CeedInt s = 0; s < e_vec_size; s++) { 1736 bool is_active_input = false; 1737 CeedEvalMode eval_mode; 1738 CeedVector vec; 1739 CeedBasis basis; 1740 1741 CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec)); 1742 // Skip non-active input 1743 is_active_input = vec == CEED_VECTOR_ACTIVE; 1744 if (!is_active_input) continue; 1745 1746 // Update unit vector 1747 if (s == 0) CeedCallBackend(CeedVectorSetValue(impl->e_vecs_in[i], 0.0)); 1748 else CeedCallBackend(CeedVectorSetValueStrided(impl->e_vecs_in[i], s - 1, e_vec_size, 0.0)); 1749 CeedCallBackend(CeedVectorSetValueStrided(impl->e_vecs_in[i], s, e_vec_size, 1.0)); 1750 1751 // Basis action 1752 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode)); 1753 switch (eval_mode) { 1754 case CEED_EVAL_NONE: 1755 CeedCallBackend(CeedVectorSetArray(impl->q_vecs_in[i], CEED_MEM_DEVICE, CEED_USE_POINTER, e_data_in[i])); 1756 break; 1757 case CEED_EVAL_INTERP: 1758 case CEED_EVAL_GRAD: 1759 case CEED_EVAL_DIV: 1760 case CEED_EVAL_CURL: 1761 CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis)); 1762 CeedCallBackend(CeedBasisApplyAtPoints(basis, num_elem, num_points, CEED_NOTRANSPOSE, eval_mode, impl->point_coords_elem, 1763 impl->e_vecs_in[i], impl->q_vecs_in[i])); 1764 break; 1765 case CEED_EVAL_WEIGHT: 1766 break; // No action 1767 } 1768 1769 // Q function 1770 CeedCallBackend(CeedQFunctionApply(qf, num_elem * max_num_points, impl->q_vecs_in, impl->q_vecs_out)); 1771 1772 // Output basis apply if needed 1773 for (CeedInt j = 0; j < num_output_fields; j++) { 1774 bool is_active_output = false; 1775 CeedInt elem_size = 0; 1776 CeedRestrictionType rstr_type; 1777 CeedEvalMode eval_mode; 1778 CeedVector l_vec; 1779 CeedElemRestriction elem_rstr; 1780 CeedBasis basis; 1781 1782 CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[j], &l_vec)); 1783 // ---- Skip non-active output 1784 is_active_output = l_vec == CEED_VECTOR_ACTIVE; 1785 if (!is_active_output) continue; 1786 1787 // ---- Check if elem size matches 1788 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[j], &elem_rstr)); 1789 CeedCallBackend(CeedElemRestrictionGetType(elem_rstr, &rstr_type)); 1790 if (is_active_at_points && rstr_type != CEED_RESTRICTION_POINTS) continue; 1791 if (rstr_type == CEED_RESTRICTION_POINTS) { 1792 CeedCallBackend(CeedElemRestrictionGetMaxPointsInElement(elem_rstr, &elem_size)); 1793 } else { 1794 CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size)); 1795 } 1796 { 1797 CeedInt num_comp = 0; 1798 1799 CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp)); 1800 if (e_vec_size != num_comp * elem_size) continue; 1801 } 1802 1803 // Basis action 1804 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[j], &eval_mode)); 1805 switch (eval_mode) { 1806 case CEED_EVAL_NONE: 1807 CeedCallBackend(CeedVectorRestoreArray(impl->e_vecs_out[j], &e_data_out[j])); 1808 break; 1809 case CEED_EVAL_INTERP: 1810 case CEED_EVAL_GRAD: 1811 case CEED_EVAL_DIV: 1812 case CEED_EVAL_CURL: 1813 CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[j], &basis)); 1814 CeedCallBackend(CeedBasisApplyAtPoints(basis, num_elem, num_points, CEED_TRANSPOSE, eval_mode, impl->point_coords_elem, 1815 impl->q_vecs_out[j], impl->e_vecs_out[j])); 1816 break; 1817 // LCOV_EXCL_START 1818 case CEED_EVAL_WEIGHT: { 1819 return CeedError(CeedOperatorReturnCeed(op), CEED_ERROR_BACKEND, "CEED_EVAL_WEIGHT cannot be an output evaluation mode"); 1820 // LCOV_EXCL_STOP 1821 } 1822 } 1823 1824 // Mask output e-vec 1825 CeedCallBackend(CeedVectorPointwiseMult(impl->e_vecs_out[j], impl->e_vecs_in[i], impl->e_vecs_out[j])); 1826 1827 // Restrict 1828 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[j], &elem_rstr)); 1829 CeedCallBackend(CeedElemRestrictionApply(elem_rstr, CEED_TRANSPOSE, impl->e_vecs_out[j], assembled, request)); 1830 1831 // Reset q_vec for 1832 if (eval_mode == CEED_EVAL_NONE) { 1833 CeedCallBackend(CeedVectorGetArrayWrite(impl->e_vecs_out[j], CEED_MEM_DEVICE, &e_data_out[j])); 1834 CeedCallBackend(CeedVectorSetArray(impl->q_vecs_out[j], CEED_MEM_DEVICE, CEED_USE_POINTER, e_data_out[j])); 1835 } 1836 } 1837 1838 // Reset vec 1839 if (s == e_vec_size - 1 && i != num_input_fields - 1) CeedCallBackend(CeedVectorSetValue(impl->q_vecs_in[i], 0.0)); 1840 } 1841 } 1842 1843 // Restore CEED_EVAL_NONE 1844 for (CeedInt i = 0; i < num_output_fields; i++) { 1845 CeedEvalMode eval_mode; 1846 1847 // Get eval_mode 1848 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 1849 1850 // Restore evec 1851 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 1852 if (eval_mode == CEED_EVAL_NONE) { 1853 CeedCallBackend(CeedVectorRestoreArray(impl->e_vecs_in[i], &e_data_in[i])); 1854 } 1855 } 1856 1857 // Restore input arrays 1858 for (CeedInt i = 0; i < num_input_fields; i++) { 1859 CeedCallBackend(CeedOperatorInputRestore_Cuda(op_input_fields[i], qf_input_fields[i], i, true, &e_data_in[i], impl)); 1860 } 1861 return CEED_ERROR_SUCCESS; 1862 } 1863 1864 //------------------------------------------------------------------------------ 1865 // Create operator 1866 //------------------------------------------------------------------------------ 1867 int CeedOperatorCreate_Cuda(CeedOperator op) { 1868 Ceed ceed; 1869 CeedOperator_Cuda *impl; 1870 1871 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 1872 CeedCallBackend(CeedCalloc(1, &impl)); 1873 CeedCallBackend(CeedOperatorSetData(op, impl)); 1874 1875 CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleQFunction", CeedOperatorLinearAssembleQFunction_Cuda)); 1876 CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleQFunctionUpdate", CeedOperatorLinearAssembleQFunctionUpdate_Cuda)); 1877 CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleAddDiagonal", CeedOperatorLinearAssembleAddDiagonal_Cuda)); 1878 CeedCallBackend( 1879 CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleAddPointBlockDiagonal", CeedOperatorLinearAssembleAddPointBlockDiagonal_Cuda)); 1880 CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleSingle", CeedSingleOperatorAssemble_Cuda)); 1881 CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "ApplyAdd", CeedOperatorApplyAdd_Cuda)); 1882 CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "Destroy", CeedOperatorDestroy_Cuda)); 1883 return CEED_ERROR_SUCCESS; 1884 } 1885 1886 //------------------------------------------------------------------------------ 1887 // Create operator AtPoints 1888 //------------------------------------------------------------------------------ 1889 int CeedOperatorCreateAtPoints_Cuda(CeedOperator op) { 1890 Ceed ceed; 1891 CeedOperator_Cuda *impl; 1892 1893 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 1894 CeedCallBackend(CeedCalloc(1, &impl)); 1895 CeedCallBackend(CeedOperatorSetData(op, impl)); 1896 1897 CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleQFunction", CeedOperatorLinearAssembleQFunctionAtPoints_Cuda)); 1898 CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleAddDiagonal", CeedOperatorLinearAssembleAddDiagonalAtPoints_Cuda)); 1899 CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "ApplyAdd", CeedOperatorApplyAddAtPoints_Cuda)); 1900 CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "Destroy", CeedOperatorDestroy_Cuda)); 1901 return CEED_ERROR_SUCCESS; 1902 } 1903 1904 //------------------------------------------------------------------------------ 1905