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