1 // Copyright (c) 2017-2022, 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 CeedCallBackend(CeedOperatorGetData(op, &impl)); 27 28 // Apply data 29 for (CeedInt i = 0; i < impl->numein + impl->numeout; i++) { 30 CeedCallBackend(CeedVectorDestroy(&impl->evecs[i])); 31 } 32 CeedCallBackend(CeedFree(&impl->evecs)); 33 34 for (CeedInt i = 0; i < impl->numein; i++) { 35 CeedCallBackend(CeedVectorDestroy(&impl->qvecsin[i])); 36 } 37 CeedCallBackend(CeedFree(&impl->qvecsin)); 38 39 for (CeedInt i = 0; i < impl->numeout; i++) { 40 CeedCallBackend(CeedVectorDestroy(&impl->qvecsout[i])); 41 } 42 CeedCallBackend(CeedFree(&impl->qvecsout)); 43 44 // QFunction assembly data 45 for (CeedInt i = 0; i < impl->qfnumactivein; i++) { 46 CeedCallBackend(CeedVectorDestroy(&impl->qfactivein[i])); 47 } 48 CeedCallBackend(CeedFree(&impl->qfactivein)); 49 50 // Diag data 51 if (impl->diag) { 52 Ceed ceed; 53 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 54 CeedCallCuda(ceed, cuModuleUnload(impl->diag->module)); 55 CeedCallBackend(CeedFree(&impl->diag->h_emodein)); 56 CeedCallBackend(CeedFree(&impl->diag->h_emodeout)); 57 CeedCallCuda(ceed, cudaFree(impl->diag->d_emodein)); 58 CeedCallCuda(ceed, cudaFree(impl->diag->d_emodeout)); 59 CeedCallCuda(ceed, cudaFree(impl->diag->d_identity)); 60 CeedCallCuda(ceed, cudaFree(impl->diag->d_interpin)); 61 CeedCallCuda(ceed, cudaFree(impl->diag->d_interpout)); 62 CeedCallCuda(ceed, cudaFree(impl->diag->d_gradin)); 63 CeedCallCuda(ceed, cudaFree(impl->diag->d_gradout)); 64 CeedCallBackend(CeedElemRestrictionDestroy(&impl->diag->pbdiagrstr)); 65 CeedCallBackend(CeedVectorDestroy(&impl->diag->elemdiag)); 66 CeedCallBackend(CeedVectorDestroy(&impl->diag->pbelemdiag)); 67 } 68 CeedCallBackend(CeedFree(&impl->diag)); 69 70 if (impl->asmb) { 71 Ceed ceed; 72 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 73 CeedCallCuda(ceed, cuModuleUnload(impl->asmb->module)); 74 CeedCallCuda(ceed, cudaFree(impl->asmb->d_B_in)); 75 CeedCallCuda(ceed, cudaFree(impl->asmb->d_B_out)); 76 } 77 CeedCallBackend(CeedFree(&impl->asmb)); 78 79 CeedCallBackend(CeedFree(&impl)); 80 return CEED_ERROR_SUCCESS; 81 } 82 83 //------------------------------------------------------------------------------ 84 // Setup infields or outfields 85 //------------------------------------------------------------------------------ 86 static int CeedOperatorSetupFields_Cuda(CeedQFunction qf, CeedOperator op, bool isinput, CeedVector *evecs, CeedVector *qvecs, CeedInt starte, 87 CeedInt numfields, CeedInt Q, CeedInt numelements) { 88 CeedInt dim, size; 89 CeedSize q_size; 90 Ceed ceed; 91 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 92 CeedBasis basis; 93 CeedElemRestriction Erestrict; 94 CeedOperatorField *opfields; 95 CeedQFunctionField *qffields; 96 CeedVector fieldvec; 97 bool strided; 98 bool skiprestrict; 99 100 if (isinput) { 101 CeedCallBackend(CeedOperatorGetFields(op, NULL, &opfields, NULL, NULL)); 102 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qffields, NULL, NULL)); 103 } else { 104 CeedCallBackend(CeedOperatorGetFields(op, NULL, NULL, NULL, &opfields)); 105 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qffields)); 106 } 107 108 // Loop over fields 109 for (CeedInt i = 0; i < numfields; i++) { 110 CeedEvalMode emode; 111 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qffields[i], &emode)); 112 113 strided = false; 114 skiprestrict = false; 115 if (emode != CEED_EVAL_WEIGHT) { 116 CeedCallBackend(CeedOperatorFieldGetElemRestriction(opfields[i], &Erestrict)); 117 118 // Check whether this field can skip the element restriction: 119 // must be passive input, with emode NONE, and have a strided restriction with CEED_STRIDES_BACKEND. 120 121 // First, check whether the field is input or output: 122 if (isinput) { 123 // Check for passive input: 124 CeedCallBackend(CeedOperatorFieldGetVector(opfields[i], &fieldvec)); 125 if (fieldvec != CEED_VECTOR_ACTIVE) { 126 // Check emode 127 if (emode == CEED_EVAL_NONE) { 128 // Check for strided restriction 129 CeedCallBackend(CeedElemRestrictionIsStrided(Erestrict, &strided)); 130 if (strided) { 131 // Check if vector is already in preferred backend ordering 132 CeedCallBackend(CeedElemRestrictionHasBackendStrides(Erestrict, &skiprestrict)); 133 } 134 } 135 } 136 } 137 if (skiprestrict) { 138 // We do not need an E-Vector, but will use the input field vector's data directly in the operator application. 139 evecs[i + starte] = NULL; 140 } else { 141 CeedCallBackend(CeedElemRestrictionCreateVector(Erestrict, NULL, &evecs[i + starte])); 142 } 143 } 144 145 switch (emode) { 146 case CEED_EVAL_NONE: 147 CeedCallBackend(CeedQFunctionFieldGetSize(qffields[i], &size)); 148 q_size = (CeedSize)numelements * Q * size; 149 CeedCallBackend(CeedVectorCreate(ceed, q_size, &qvecs[i])); 150 break; 151 case CEED_EVAL_INTERP: 152 CeedCallBackend(CeedQFunctionFieldGetSize(qffields[i], &size)); 153 q_size = (CeedSize)numelements * Q * size; 154 CeedCallBackend(CeedVectorCreate(ceed, q_size, &qvecs[i])); 155 break; 156 case CEED_EVAL_GRAD: 157 CeedCallBackend(CeedOperatorFieldGetBasis(opfields[i], &basis)); 158 CeedCallBackend(CeedQFunctionFieldGetSize(qffields[i], &size)); 159 CeedCallBackend(CeedBasisGetDimension(basis, &dim)); 160 q_size = (CeedSize)numelements * Q * size; 161 CeedCallBackend(CeedVectorCreate(ceed, q_size, &qvecs[i])); 162 break; 163 case CEED_EVAL_WEIGHT: // Only on input fields 164 CeedCallBackend(CeedOperatorFieldGetBasis(opfields[i], &basis)); 165 q_size = (CeedSize)numelements * Q; 166 CeedCallBackend(CeedVectorCreate(ceed, q_size, &qvecs[i])); 167 CeedCallBackend(CeedBasisApply(basis, numelements, CEED_NOTRANSPOSE, CEED_EVAL_WEIGHT, CEED_VECTOR_NONE, qvecs[i])); 168 break; 169 case CEED_EVAL_DIV: 170 break; // TODO: Not implemented 171 case CEED_EVAL_CURL: 172 break; // TODO: Not implemented 173 } 174 } 175 return CEED_ERROR_SUCCESS; 176 } 177 178 //------------------------------------------------------------------------------ 179 // CeedOperator needs to connect all the named fields (be they active or passive) to the named inputs and outputs of its CeedQFunction. 180 //------------------------------------------------------------------------------ 181 static int CeedOperatorSetup_Cuda(CeedOperator op) { 182 bool setupdone; 183 CeedCallBackend(CeedOperatorIsSetupDone(op, &setupdone)); 184 if (setupdone) return CEED_ERROR_SUCCESS; 185 Ceed ceed; 186 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 187 CeedOperator_Cuda *impl; 188 CeedCallBackend(CeedOperatorGetData(op, &impl)); 189 CeedQFunction qf; 190 CeedCallBackend(CeedOperatorGetQFunction(op, &qf)); 191 CeedInt Q, numelements, numinputfields, numoutputfields; 192 CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q)); 193 CeedCallBackend(CeedOperatorGetNumElements(op, &numelements)); 194 CeedOperatorField *opinputfields, *opoutputfields; 195 CeedCallBackend(CeedOperatorGetFields(op, &numinputfields, &opinputfields, &numoutputfields, &opoutputfields)); 196 CeedQFunctionField *qfinputfields, *qfoutputfields; 197 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qfinputfields, NULL, &qfoutputfields)); 198 199 // Allocate 200 CeedCallBackend(CeedCalloc(numinputfields + numoutputfields, &impl->evecs)); 201 202 CeedCallBackend(CeedCalloc(CEED_FIELD_MAX, &impl->qvecsin)); 203 CeedCallBackend(CeedCalloc(CEED_FIELD_MAX, &impl->qvecsout)); 204 205 impl->numein = numinputfields; 206 impl->numeout = numoutputfields; 207 208 // Set up infield and outfield evecs and qvecs 209 // Infields 210 CeedCallBackend(CeedOperatorSetupFields_Cuda(qf, op, true, impl->evecs, impl->qvecsin, 0, numinputfields, Q, numelements)); 211 212 // Outfields 213 CeedCallBackend(CeedOperatorSetupFields_Cuda(qf, op, false, impl->evecs, impl->qvecsout, numinputfields, numoutputfields, Q, numelements)); 214 215 CeedCallBackend(CeedOperatorSetSetupDone(op)); 216 return CEED_ERROR_SUCCESS; 217 } 218 219 //------------------------------------------------------------------------------ 220 // Setup Operator Inputs 221 //------------------------------------------------------------------------------ 222 static inline int CeedOperatorSetupInputs_Cuda(CeedInt numinputfields, CeedQFunctionField *qfinputfields, CeedOperatorField *opinputfields, 223 CeedVector invec, const bool skipactive, CeedScalar *edata[2 * CEED_FIELD_MAX], 224 CeedOperator_Cuda *impl, CeedRequest *request) { 225 CeedEvalMode emode; 226 CeedVector vec; 227 CeedElemRestriction Erestrict; 228 229 for (CeedInt i = 0; i < numinputfields; i++) { 230 // Get input vector 231 CeedCallBackend(CeedOperatorFieldGetVector(opinputfields[i], &vec)); 232 if (vec == CEED_VECTOR_ACTIVE) { 233 if (skipactive) continue; 234 else vec = invec; 235 } 236 237 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qfinputfields[i], &emode)); 238 if (emode == CEED_EVAL_WEIGHT) { // Skip 239 } else { 240 // Get input vector 241 CeedCallBackend(CeedOperatorFieldGetVector(opinputfields[i], &vec)); 242 // Get input element restriction 243 CeedCallBackend(CeedOperatorFieldGetElemRestriction(opinputfields[i], &Erestrict)); 244 if (vec == CEED_VECTOR_ACTIVE) vec = invec; 245 // Restrict, if necessary 246 if (!impl->evecs[i]) { 247 // No restriction for this field; read data directly from vec. 248 CeedCallBackend(CeedVectorGetArrayRead(vec, CEED_MEM_DEVICE, (const CeedScalar **)&edata[i])); 249 } else { 250 CeedCallBackend(CeedElemRestrictionApply(Erestrict, CEED_NOTRANSPOSE, vec, impl->evecs[i], request)); 251 // Get evec 252 CeedCallBackend(CeedVectorGetArrayRead(impl->evecs[i], CEED_MEM_DEVICE, (const CeedScalar **)&edata[i])); 253 } 254 } 255 } 256 return CEED_ERROR_SUCCESS; 257 } 258 259 //------------------------------------------------------------------------------ 260 // Input Basis Action 261 //------------------------------------------------------------------------------ 262 static inline int CeedOperatorInputBasis_Cuda(CeedInt numelements, CeedQFunctionField *qfinputfields, CeedOperatorField *opinputfields, 263 CeedInt numinputfields, const bool skipactive, CeedScalar *edata[2 * CEED_FIELD_MAX], 264 CeedOperator_Cuda *impl) { 265 CeedInt elemsize, size; 266 CeedElemRestriction Erestrict; 267 CeedEvalMode emode; 268 CeedBasis basis; 269 270 for (CeedInt i = 0; i < numinputfields; i++) { 271 // Skip active input 272 if (skipactive) { 273 CeedVector vec; 274 CeedCallBackend(CeedOperatorFieldGetVector(opinputfields[i], &vec)); 275 if (vec == CEED_VECTOR_ACTIVE) continue; 276 } 277 // Get elemsize, emode, size 278 CeedCallBackend(CeedOperatorFieldGetElemRestriction(opinputfields[i], &Erestrict)); 279 CeedCallBackend(CeedElemRestrictionGetElementSize(Erestrict, &elemsize)); 280 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qfinputfields[i], &emode)); 281 CeedCallBackend(CeedQFunctionFieldGetSize(qfinputfields[i], &size)); 282 // Basis action 283 switch (emode) { 284 case CEED_EVAL_NONE: 285 CeedCallBackend(CeedVectorSetArray(impl->qvecsin[i], CEED_MEM_DEVICE, CEED_USE_POINTER, edata[i])); 286 break; 287 case CEED_EVAL_INTERP: 288 CeedCallBackend(CeedOperatorFieldGetBasis(opinputfields[i], &basis)); 289 CeedCallBackend(CeedBasisApply(basis, numelements, CEED_NOTRANSPOSE, CEED_EVAL_INTERP, impl->evecs[i], impl->qvecsin[i])); 290 break; 291 case CEED_EVAL_GRAD: 292 CeedCallBackend(CeedOperatorFieldGetBasis(opinputfields[i], &basis)); 293 CeedCallBackend(CeedBasisApply(basis, numelements, CEED_NOTRANSPOSE, CEED_EVAL_GRAD, impl->evecs[i], impl->qvecsin[i])); 294 break; 295 case CEED_EVAL_WEIGHT: 296 break; // No action 297 case CEED_EVAL_DIV: 298 break; // TODO: Not implemented 299 case CEED_EVAL_CURL: 300 break; // TODO: Not implemented 301 } 302 } 303 return CEED_ERROR_SUCCESS; 304 } 305 306 //------------------------------------------------------------------------------ 307 // Restore Input Vectors 308 //------------------------------------------------------------------------------ 309 static inline int CeedOperatorRestoreInputs_Cuda(CeedInt numinputfields, CeedQFunctionField *qfinputfields, CeedOperatorField *opinputfields, 310 const bool skipactive, CeedScalar *edata[2 * CEED_FIELD_MAX], CeedOperator_Cuda *impl) { 311 CeedEvalMode emode; 312 CeedVector vec; 313 314 for (CeedInt i = 0; i < numinputfields; i++) { 315 // Skip active input 316 if (skipactive) { 317 CeedCallBackend(CeedOperatorFieldGetVector(opinputfields[i], &vec)); 318 if (vec == CEED_VECTOR_ACTIVE) continue; 319 } 320 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qfinputfields[i], &emode)); 321 if (emode == CEED_EVAL_WEIGHT) { // Skip 322 } else { 323 if (!impl->evecs[i]) { // This was a skiprestrict case 324 CeedCallBackend(CeedOperatorFieldGetVector(opinputfields[i], &vec)); 325 CeedCallBackend(CeedVectorRestoreArrayRead(vec, (const CeedScalar **)&edata[i])); 326 } else { 327 CeedCallBackend(CeedVectorRestoreArrayRead(impl->evecs[i], (const CeedScalar **)&edata[i])); 328 } 329 } 330 } 331 return CEED_ERROR_SUCCESS; 332 } 333 334 //------------------------------------------------------------------------------ 335 // Apply and add to output 336 //------------------------------------------------------------------------------ 337 static int CeedOperatorApplyAdd_Cuda(CeedOperator op, CeedVector invec, CeedVector outvec, CeedRequest *request) { 338 CeedOperator_Cuda *impl; 339 CeedCallBackend(CeedOperatorGetData(op, &impl)); 340 CeedQFunction qf; 341 CeedCallBackend(CeedOperatorGetQFunction(op, &qf)); 342 CeedInt Q, numelements, elemsize, numinputfields, numoutputfields, size; 343 CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q)); 344 CeedCallBackend(CeedOperatorGetNumElements(op, &numelements)); 345 CeedOperatorField *opinputfields, *opoutputfields; 346 CeedCallBackend(CeedOperatorGetFields(op, &numinputfields, &opinputfields, &numoutputfields, &opoutputfields)); 347 CeedQFunctionField *qfinputfields, *qfoutputfields; 348 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qfinputfields, NULL, &qfoutputfields)); 349 CeedEvalMode emode; 350 CeedVector vec; 351 CeedBasis basis; 352 CeedElemRestriction Erestrict; 353 CeedScalar *edata[2 * CEED_FIELD_MAX] = {NULL}; 354 355 // Setup 356 CeedCallBackend(CeedOperatorSetup_Cuda(op)); 357 358 // Input Evecs and Restriction 359 CeedCallBackend(CeedOperatorSetupInputs_Cuda(numinputfields, qfinputfields, opinputfields, invec, false, edata, impl, request)); 360 361 // Input basis apply if needed 362 CeedCallBackend(CeedOperatorInputBasis_Cuda(numelements, qfinputfields, opinputfields, numinputfields, false, edata, impl)); 363 364 // Output pointers, as necessary 365 for (CeedInt i = 0; i < numoutputfields; i++) { 366 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qfoutputfields[i], &emode)); 367 if (emode == CEED_EVAL_NONE) { 368 // Set the output Q-Vector to use the E-Vector data directly. 369 CeedCallBackend(CeedVectorGetArrayWrite(impl->evecs[i + impl->numein], CEED_MEM_DEVICE, &edata[i + numinputfields])); 370 CeedCallBackend(CeedVectorSetArray(impl->qvecsout[i], CEED_MEM_DEVICE, CEED_USE_POINTER, edata[i + numinputfields])); 371 } 372 } 373 374 // Q function 375 CeedCallBackend(CeedQFunctionApply(qf, numelements * Q, impl->qvecsin, impl->qvecsout)); 376 377 // Output basis apply if needed 378 for (CeedInt i = 0; i < numoutputfields; i++) { 379 // Get elemsize, emode, size 380 CeedCallBackend(CeedOperatorFieldGetElemRestriction(opoutputfields[i], &Erestrict)); 381 CeedCallBackend(CeedElemRestrictionGetElementSize(Erestrict, &elemsize)); 382 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qfoutputfields[i], &emode)); 383 CeedCallBackend(CeedQFunctionFieldGetSize(qfoutputfields[i], &size)); 384 // Basis action 385 switch (emode) { 386 case CEED_EVAL_NONE: 387 break; 388 case CEED_EVAL_INTERP: 389 CeedCallBackend(CeedOperatorFieldGetBasis(opoutputfields[i], &basis)); 390 CeedCallBackend(CeedBasisApply(basis, numelements, CEED_TRANSPOSE, CEED_EVAL_INTERP, impl->qvecsout[i], impl->evecs[i + impl->numein])); 391 break; 392 case CEED_EVAL_GRAD: 393 CeedCallBackend(CeedOperatorFieldGetBasis(opoutputfields[i], &basis)); 394 CeedCallBackend(CeedBasisApply(basis, numelements, CEED_TRANSPOSE, CEED_EVAL_GRAD, impl->qvecsout[i], impl->evecs[i + impl->numein])); 395 break; 396 // LCOV_EXCL_START 397 case CEED_EVAL_WEIGHT: { 398 Ceed ceed; 399 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 400 return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_WEIGHT cannot be an output evaluation mode"); 401 break; // Should not occur 402 } 403 case CEED_EVAL_DIV: 404 break; // TODO: Not implemented 405 case CEED_EVAL_CURL: 406 break; // TODO: Not implemented 407 // LCOV_EXCL_STOP 408 } 409 } 410 411 // Output restriction 412 for (CeedInt i = 0; i < numoutputfields; i++) { 413 // Restore evec 414 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qfoutputfields[i], &emode)); 415 if (emode == CEED_EVAL_NONE) { 416 CeedCallBackend(CeedVectorRestoreArray(impl->evecs[i + impl->numein], &edata[i + numinputfields])); 417 } 418 // Get output vector 419 CeedCallBackend(CeedOperatorFieldGetVector(opoutputfields[i], &vec)); 420 // Restrict 421 CeedCallBackend(CeedOperatorFieldGetElemRestriction(opoutputfields[i], &Erestrict)); 422 // Active 423 if (vec == CEED_VECTOR_ACTIVE) vec = outvec; 424 425 CeedCallBackend(CeedElemRestrictionApply(Erestrict, CEED_TRANSPOSE, impl->evecs[i + impl->numein], vec, request)); 426 } 427 428 // Restore input arrays 429 CeedCallBackend(CeedOperatorRestoreInputs_Cuda(numinputfields, qfinputfields, opinputfields, false, edata, impl)); 430 return CEED_ERROR_SUCCESS; 431 } 432 433 //------------------------------------------------------------------------------ 434 // Core code for assembling linear QFunction 435 //------------------------------------------------------------------------------ 436 static inline int CeedOperatorLinearAssembleQFunctionCore_Cuda(CeedOperator op, bool build_objects, CeedVector *assembled, CeedElemRestriction *rstr, 437 CeedRequest *request) { 438 Ceed ceed, ceedparent; 439 CeedOperator_Cuda *impl; 440 CeedQFunction qf; 441 CeedQFunctionField *qfinputfields, *qfoutputfields; 442 CeedOperatorField *opinputfields, *opoutputfields; 443 CeedVector vec, *activein; 444 CeedInt numactivein, numactiveout, Q, numelements, numinputfields, numoutputfields, size; 445 CeedSize q_size; 446 CeedScalar *a, *tmp, *edata[2 * CEED_FIELD_MAX] = {NULL}; 447 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 448 CeedCallBackend(CeedOperatorGetFallbackParentCeed(op, &ceedparent)); 449 CeedCallBackend(CeedOperatorGetData(op, &impl)); 450 activein = impl->qfactivein; 451 numactivein = impl->qfnumactivein, numactiveout = impl->qfnumactiveout; 452 CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q)); 453 CeedCallBackend(CeedOperatorGetNumElements(op, &numelements)); 454 CeedCallBackend(CeedOperatorGetQFunction(op, &qf)); 455 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qfinputfields, NULL, &qfoutputfields)); 456 CeedCallBackend(CeedOperatorGetFields(op, &numinputfields, &opinputfields, &numoutputfields, &opoutputfields)); 457 458 // Setup 459 CeedCallBackend(CeedOperatorSetup_Cuda(op)); 460 461 // Check for identity 462 bool identityqf; 463 CeedCallBackend(CeedQFunctionIsIdentity(qf, &identityqf)); 464 CeedCheck(!identityqf, ceed, CEED_ERROR_BACKEND, "Assembling identity QFunctions not supported"); 465 466 // Input Evecs and Restriction 467 CeedCallBackend(CeedOperatorSetupInputs_Cuda(numinputfields, qfinputfields, opinputfields, NULL, true, edata, impl, request)); 468 469 // Count number of active input fields 470 if (!numactivein) { 471 for (CeedInt i = 0; i < numinputfields; i++) { 472 // Get input vector 473 CeedCallBackend(CeedOperatorFieldGetVector(opinputfields[i], &vec)); 474 // Check if active input 475 if (vec == CEED_VECTOR_ACTIVE) { 476 CeedCallBackend(CeedQFunctionFieldGetSize(qfinputfields[i], &size)); 477 CeedCallBackend(CeedVectorSetValue(impl->qvecsin[i], 0.0)); 478 CeedCallBackend(CeedVectorGetArray(impl->qvecsin[i], CEED_MEM_DEVICE, &tmp)); 479 CeedCallBackend(CeedRealloc(numactivein + size, &activein)); 480 for (CeedInt field = 0; field < size; field++) { 481 q_size = (CeedSize)Q * numelements; 482 CeedCallBackend(CeedVectorCreate(ceed, q_size, &activein[numactivein + field])); 483 CeedCallBackend(CeedVectorSetArray(activein[numactivein + field], CEED_MEM_DEVICE, CEED_USE_POINTER, &tmp[field * Q * numelements])); 484 } 485 numactivein += size; 486 CeedCallBackend(CeedVectorRestoreArray(impl->qvecsin[i], &tmp)); 487 } 488 } 489 impl->qfnumactivein = numactivein; 490 impl->qfactivein = activein; 491 } 492 493 // Count number of active output fields 494 if (!numactiveout) { 495 for (CeedInt i = 0; i < numoutputfields; i++) { 496 // Get output vector 497 CeedCallBackend(CeedOperatorFieldGetVector(opoutputfields[i], &vec)); 498 // Check if active output 499 if (vec == CEED_VECTOR_ACTIVE) { 500 CeedCallBackend(CeedQFunctionFieldGetSize(qfoutputfields[i], &size)); 501 numactiveout += size; 502 } 503 } 504 impl->qfnumactiveout = numactiveout; 505 } 506 507 // Check sizes 508 CeedCheck(numactivein > 0 && numactiveout > 0, ceed, CEED_ERROR_BACKEND, "Cannot assemble QFunction without active inputs and outputs"); 509 510 // Build objects if needed 511 if (build_objects) { 512 // Create output restriction 513 CeedInt strides[3] = {1, numelements * Q, Q}; /* *NOPAD* */ 514 CeedCallBackend(CeedElemRestrictionCreateStrided(ceedparent, numelements, Q, numactivein * numactiveout, 515 numactivein * numactiveout * numelements * Q, strides, rstr)); 516 // Create assembled vector 517 CeedSize l_size = (CeedSize)numelements * Q * numactivein * numactiveout; 518 CeedCallBackend(CeedVectorCreate(ceedparent, l_size, assembled)); 519 } 520 CeedCallBackend(CeedVectorSetValue(*assembled, 0.0)); 521 CeedCallBackend(CeedVectorGetArray(*assembled, CEED_MEM_DEVICE, &a)); 522 523 // Input basis apply 524 CeedCallBackend(CeedOperatorInputBasis_Cuda(numelements, qfinputfields, opinputfields, numinputfields, true, edata, impl)); 525 526 // Assemble QFunction 527 for (CeedInt in = 0; in < numactivein; in++) { 528 // Set Inputs 529 CeedCallBackend(CeedVectorSetValue(activein[in], 1.0)); 530 if (numactivein > 1) { 531 CeedCallBackend(CeedVectorSetValue(activein[(in + numactivein - 1) % numactivein], 0.0)); 532 } 533 // Set Outputs 534 for (CeedInt out = 0; out < numoutputfields; out++) { 535 // Get output vector 536 CeedCallBackend(CeedOperatorFieldGetVector(opoutputfields[out], &vec)); 537 // Check if active output 538 if (vec == CEED_VECTOR_ACTIVE) { 539 CeedCallBackend(CeedVectorSetArray(impl->qvecsout[out], CEED_MEM_DEVICE, CEED_USE_POINTER, a)); 540 CeedCallBackend(CeedQFunctionFieldGetSize(qfoutputfields[out], &size)); 541 a += size * Q * numelements; // Advance the pointer by the size of the output 542 } 543 } 544 // Apply QFunction 545 CeedCallBackend(CeedQFunctionApply(qf, Q * numelements, impl->qvecsin, impl->qvecsout)); 546 } 547 548 // Un-set output Qvecs to prevent accidental overwrite of Assembled 549 for (CeedInt out = 0; out < numoutputfields; out++) { 550 // Get output vector 551 CeedCallBackend(CeedOperatorFieldGetVector(opoutputfields[out], &vec)); 552 // Check if active output 553 if (vec == CEED_VECTOR_ACTIVE) { 554 CeedCallBackend(CeedVectorTakeArray(impl->qvecsout[out], CEED_MEM_DEVICE, NULL)); 555 } 556 } 557 558 // Restore input arrays 559 CeedCallBackend(CeedOperatorRestoreInputs_Cuda(numinputfields, qfinputfields, opinputfields, true, edata, impl)); 560 561 // Restore output 562 CeedCallBackend(CeedVectorRestoreArray(*assembled, &a)); 563 564 return CEED_ERROR_SUCCESS; 565 } 566 567 //------------------------------------------------------------------------------ 568 // Assemble Linear QFunction 569 //------------------------------------------------------------------------------ 570 static int CeedOperatorLinearAssembleQFunction_Cuda(CeedOperator op, CeedVector *assembled, CeedElemRestriction *rstr, CeedRequest *request) { 571 return CeedOperatorLinearAssembleQFunctionCore_Cuda(op, true, assembled, rstr, request); 572 } 573 574 //------------------------------------------------------------------------------ 575 // Update Assembled Linear QFunction 576 //------------------------------------------------------------------------------ 577 static int CeedOperatorLinearAssembleQFunctionUpdate_Cuda(CeedOperator op, CeedVector assembled, CeedElemRestriction rstr, CeedRequest *request) { 578 return CeedOperatorLinearAssembleQFunctionCore_Cuda(op, false, &assembled, &rstr, request); 579 } 580 581 //------------------------------------------------------------------------------ 582 // Create point block restriction 583 //------------------------------------------------------------------------------ 584 static int CreatePBRestriction(CeedElemRestriction rstr, CeedElemRestriction *pbRstr) { 585 Ceed ceed; 586 CeedCallBackend(CeedElemRestrictionGetCeed(rstr, &ceed)); 587 const CeedInt *offsets; 588 CeedCallBackend(CeedElemRestrictionGetOffsets(rstr, CEED_MEM_HOST, &offsets)); 589 590 // Expand offsets 591 CeedInt nelem, ncomp, elemsize, compstride, *pbOffsets; 592 CeedSize l_size; 593 CeedCallBackend(CeedElemRestrictionGetNumElements(rstr, &nelem)); 594 CeedCallBackend(CeedElemRestrictionGetNumComponents(rstr, &ncomp)); 595 CeedCallBackend(CeedElemRestrictionGetElementSize(rstr, &elemsize)); 596 CeedCallBackend(CeedElemRestrictionGetCompStride(rstr, &compstride)); 597 CeedCallBackend(CeedElemRestrictionGetLVectorSize(rstr, &l_size)); 598 CeedInt shift = ncomp; 599 if (compstride != 1) shift *= ncomp; 600 CeedCallBackend(CeedCalloc(nelem * elemsize, &pbOffsets)); 601 for (CeedInt i = 0; i < nelem * elemsize; i++) { 602 pbOffsets[i] = offsets[i] * shift; 603 } 604 605 // Create new restriction 606 CeedCallBackend( 607 CeedElemRestrictionCreate(ceed, nelem, elemsize, ncomp * ncomp, 1, l_size * ncomp, CEED_MEM_HOST, CEED_OWN_POINTER, pbOffsets, pbRstr)); 608 609 // Cleanup 610 CeedCallBackend(CeedElemRestrictionRestoreOffsets(rstr, &offsets)); 611 612 return CEED_ERROR_SUCCESS; 613 } 614 615 //------------------------------------------------------------------------------ 616 // Assemble diagonal setup 617 //------------------------------------------------------------------------------ 618 static inline int CeedOperatorAssembleDiagonalSetup_Cuda(CeedOperator op, const bool pointBlock, CeedInt use_ceedsize_idx) { 619 Ceed ceed; 620 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 621 CeedQFunction qf; 622 CeedCallBackend(CeedOperatorGetQFunction(op, &qf)); 623 CeedInt numinputfields, numoutputfields; 624 CeedCallBackend(CeedQFunctionGetNumArgs(qf, &numinputfields, &numoutputfields)); 625 626 // Determine active input basis 627 CeedOperatorField *opfields; 628 CeedQFunctionField *qffields; 629 CeedCallBackend(CeedOperatorGetFields(op, NULL, &opfields, NULL, NULL)); 630 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qffields, NULL, NULL)); 631 CeedInt numemodein = 0, ncomp = 0, dim = 1; 632 CeedEvalMode *emodein = NULL; 633 CeedBasis basisin = NULL; 634 CeedElemRestriction rstrin = NULL; 635 for (CeedInt i = 0; i < numinputfields; i++) { 636 CeedVector vec; 637 CeedCallBackend(CeedOperatorFieldGetVector(opfields[i], &vec)); 638 if (vec == CEED_VECTOR_ACTIVE) { 639 CeedElemRestriction rstr; 640 CeedCallBackend(CeedOperatorFieldGetBasis(opfields[i], &basisin)); 641 CeedCallBackend(CeedBasisGetNumComponents(basisin, &ncomp)); 642 CeedCallBackend(CeedBasisGetDimension(basisin, &dim)); 643 CeedCallBackend(CeedOperatorFieldGetElemRestriction(opfields[i], &rstr)); 644 CeedCheck(!rstrin || rstrin == rstr, ceed, CEED_ERROR_BACKEND, 645 "Backend does not implement multi-field non-composite operator diagonal assembly"); 646 rstrin = rstr; 647 CeedEvalMode emode; 648 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qffields[i], &emode)); 649 switch (emode) { 650 case CEED_EVAL_NONE: 651 case CEED_EVAL_INTERP: 652 CeedCallBackend(CeedRealloc(numemodein + 1, &emodein)); 653 emodein[numemodein] = emode; 654 numemodein += 1; 655 break; 656 case CEED_EVAL_GRAD: 657 CeedCallBackend(CeedRealloc(numemodein + dim, &emodein)); 658 for (CeedInt d = 0; d < dim; d++) emodein[numemodein + d] = emode; 659 numemodein += dim; 660 break; 661 case CEED_EVAL_WEIGHT: 662 case CEED_EVAL_DIV: 663 case CEED_EVAL_CURL: 664 break; // Caught by QF Assembly 665 } 666 } 667 } 668 669 // Determine active output basis 670 CeedCallBackend(CeedOperatorGetFields(op, NULL, NULL, NULL, &opfields)); 671 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qffields)); 672 CeedInt numemodeout = 0; 673 CeedEvalMode *emodeout = NULL; 674 CeedBasis basisout = NULL; 675 CeedElemRestriction rstrout = NULL; 676 for (CeedInt i = 0; i < numoutputfields; i++) { 677 CeedVector vec; 678 CeedCallBackend(CeedOperatorFieldGetVector(opfields[i], &vec)); 679 if (vec == CEED_VECTOR_ACTIVE) { 680 CeedElemRestriction rstr; 681 CeedCallBackend(CeedOperatorFieldGetBasis(opfields[i], &basisout)); 682 CeedCallBackend(CeedOperatorFieldGetElemRestriction(opfields[i], &rstr)); 683 CeedCheck(!rstrout || rstrout == rstr, ceed, CEED_ERROR_BACKEND, 684 "Backend does not implement multi-field non-composite operator diagonal assembly"); 685 rstrout = rstr; 686 CeedEvalMode emode; 687 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qffields[i], &emode)); 688 switch (emode) { 689 case CEED_EVAL_NONE: 690 case CEED_EVAL_INTERP: 691 CeedCallBackend(CeedRealloc(numemodeout + 1, &emodeout)); 692 emodeout[numemodeout] = emode; 693 numemodeout += 1; 694 break; 695 case CEED_EVAL_GRAD: 696 CeedCallBackend(CeedRealloc(numemodeout + dim, &emodeout)); 697 for (CeedInt d = 0; d < dim; d++) emodeout[numemodeout + d] = emode; 698 numemodeout += dim; 699 break; 700 case CEED_EVAL_WEIGHT: 701 case CEED_EVAL_DIV: 702 case CEED_EVAL_CURL: 703 break; // Caught by QF Assembly 704 } 705 } 706 } 707 708 // Operator data struct 709 CeedOperator_Cuda *impl; 710 CeedCallBackend(CeedOperatorGetData(op, &impl)); 711 CeedCallBackend(CeedCalloc(1, &impl->diag)); 712 CeedOperatorDiag_Cuda *diag = impl->diag; 713 diag->basisin = basisin; 714 diag->basisout = basisout; 715 diag->h_emodein = emodein; 716 diag->h_emodeout = emodeout; 717 diag->numemodein = numemodein; 718 diag->numemodeout = numemodeout; 719 720 // Assemble kernel 721 char *diagonal_kernel_path, *diagonal_kernel_source; 722 CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-operator-assemble-diagonal.h", &diagonal_kernel_path)); 723 CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Diagonal Assembly Kernel Source -----\n"); 724 CeedCallBackend(CeedLoadSourceToBuffer(ceed, diagonal_kernel_path, &diagonal_kernel_source)); 725 CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Diagonal Assembly Source Complete! -----\n"); 726 CeedInt nnodes, nqpts; 727 CeedCallBackend(CeedBasisGetNumNodes(basisin, &nnodes)); 728 CeedCallBackend(CeedBasisGetNumQuadraturePoints(basisin, &nqpts)); 729 diag->nnodes = nnodes; 730 CeedCallCuda(ceed, CeedCompile_Cuda(ceed, diagonal_kernel_source, &diag->module, 6, "NUMEMODEIN", numemodein, "NUMEMODEOUT", numemodeout, "NNODES", 731 nnodes, "NQPTS", nqpts, "NCOMP", ncomp, "CEEDSIZE", use_ceedsize_idx)); 732 CeedCallCuda(ceed, CeedGetKernel_Cuda(ceed, diag->module, "linearDiagonal", &diag->linearDiagonal)); 733 CeedCallCuda(ceed, CeedGetKernel_Cuda(ceed, diag->module, "linearPointBlockDiagonal", &diag->linearPointBlock)); 734 CeedCallBackend(CeedFree(&diagonal_kernel_path)); 735 CeedCallBackend(CeedFree(&diagonal_kernel_source)); 736 737 // Basis matrices 738 const CeedInt qBytes = nqpts * sizeof(CeedScalar); 739 const CeedInt iBytes = qBytes * nnodes; 740 const CeedInt gBytes = qBytes * nnodes * dim; 741 const CeedInt eBytes = sizeof(CeedEvalMode); 742 const CeedScalar *interpin, *interpout, *gradin, *gradout; 743 744 // CEED_EVAL_NONE 745 CeedScalar *identity = NULL; 746 bool evalNone = false; 747 for (CeedInt i = 0; i < numemodein; i++) evalNone = evalNone || (emodein[i] == CEED_EVAL_NONE); 748 for (CeedInt i = 0; i < numemodeout; i++) evalNone = evalNone || (emodeout[i] == CEED_EVAL_NONE); 749 if (evalNone) { 750 CeedCallBackend(CeedCalloc(nqpts * nnodes, &identity)); 751 for (CeedInt i = 0; i < (nnodes < nqpts ? nnodes : nqpts); i++) identity[i * nnodes + i] = 1.0; 752 CeedCallCuda(ceed, cudaMalloc((void **)&diag->d_identity, iBytes)); 753 CeedCallCuda(ceed, cudaMemcpy(diag->d_identity, identity, iBytes, cudaMemcpyHostToDevice)); 754 } 755 756 // CEED_EVAL_INTERP 757 CeedCallBackend(CeedBasisGetInterp(basisin, &interpin)); 758 CeedCallCuda(ceed, cudaMalloc((void **)&diag->d_interpin, iBytes)); 759 CeedCallCuda(ceed, cudaMemcpy(diag->d_interpin, interpin, iBytes, cudaMemcpyHostToDevice)); 760 CeedCallBackend(CeedBasisGetInterp(basisout, &interpout)); 761 CeedCallCuda(ceed, cudaMalloc((void **)&diag->d_interpout, iBytes)); 762 CeedCallCuda(ceed, cudaMemcpy(diag->d_interpout, interpout, iBytes, cudaMemcpyHostToDevice)); 763 764 // CEED_EVAL_GRAD 765 CeedCallBackend(CeedBasisGetGrad(basisin, &gradin)); 766 CeedCallCuda(ceed, cudaMalloc((void **)&diag->d_gradin, gBytes)); 767 CeedCallCuda(ceed, cudaMemcpy(diag->d_gradin, gradin, gBytes, cudaMemcpyHostToDevice)); 768 CeedCallBackend(CeedBasisGetGrad(basisout, &gradout)); 769 CeedCallCuda(ceed, cudaMalloc((void **)&diag->d_gradout, gBytes)); 770 CeedCallCuda(ceed, cudaMemcpy(diag->d_gradout, gradout, gBytes, cudaMemcpyHostToDevice)); 771 772 // Arrays of emodes 773 CeedCallCuda(ceed, cudaMalloc((void **)&diag->d_emodein, numemodein * eBytes)); 774 CeedCallCuda(ceed, cudaMemcpy(diag->d_emodein, emodein, numemodein * eBytes, cudaMemcpyHostToDevice)); 775 CeedCallCuda(ceed, cudaMalloc((void **)&diag->d_emodeout, numemodeout * eBytes)); 776 CeedCallCuda(ceed, cudaMemcpy(diag->d_emodeout, emodeout, numemodeout * eBytes, cudaMemcpyHostToDevice)); 777 778 // Restriction 779 diag->diagrstr = rstrout; 780 781 return CEED_ERROR_SUCCESS; 782 } 783 784 //------------------------------------------------------------------------------ 785 // Assemble diagonal common code 786 //------------------------------------------------------------------------------ 787 static inline int CeedOperatorAssembleDiagonalCore_Cuda(CeedOperator op, CeedVector assembled, CeedRequest *request, const bool pointBlock) { 788 Ceed ceed; 789 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 790 CeedOperator_Cuda *impl; 791 CeedCallBackend(CeedOperatorGetData(op, &impl)); 792 793 // Assemble QFunction 794 CeedVector assembledqf = NULL; 795 CeedElemRestriction rstr = NULL; 796 CeedCallBackend(CeedOperatorLinearAssembleQFunctionBuildOrUpdate(op, &assembledqf, &rstr, request)); 797 CeedCallBackend(CeedElemRestrictionDestroy(&rstr)); 798 799 CeedSize assembled_length = 0, assembledqf_length = 0; 800 CeedCallBackend(CeedVectorGetLength(assembled, &assembled_length)); 801 CeedCallBackend(CeedVectorGetLength(assembledqf, &assembledqf_length)); 802 CeedInt use_ceedsize_idx = 0; 803 if ((assembled_length > INT_MAX) || (assembledqf_length > INT_MAX)) use_ceedsize_idx = 1; 804 805 // Setup 806 if (!impl->diag) CeedCallBackend(CeedOperatorAssembleDiagonalSetup_Cuda(op, pointBlock, use_ceedsize_idx)); 807 CeedOperatorDiag_Cuda *diag = impl->diag; 808 assert(diag != NULL); 809 810 // Restriction 811 if (pointBlock && !diag->pbdiagrstr) { 812 CeedElemRestriction pbdiagrstr; 813 CeedCallBackend(CreatePBRestriction(diag->diagrstr, &pbdiagrstr)); 814 diag->pbdiagrstr = pbdiagrstr; 815 } 816 CeedElemRestriction diagrstr = pointBlock ? diag->pbdiagrstr : diag->diagrstr; 817 818 // Create diagonal vector 819 CeedVector elemdiag = pointBlock ? diag->pbelemdiag : diag->elemdiag; 820 if (!elemdiag) { 821 CeedCallBackend(CeedElemRestrictionCreateVector(diagrstr, NULL, &elemdiag)); 822 if (pointBlock) diag->pbelemdiag = elemdiag; 823 else diag->elemdiag = elemdiag; 824 } 825 CeedCallBackend(CeedVectorSetValue(elemdiag, 0.0)); 826 827 // Assemble element operator diagonals 828 CeedScalar *elemdiagarray; 829 const CeedScalar *assembledqfarray; 830 CeedCallBackend(CeedVectorGetArray(elemdiag, CEED_MEM_DEVICE, &elemdiagarray)); 831 CeedCallBackend(CeedVectorGetArrayRead(assembledqf, CEED_MEM_DEVICE, &assembledqfarray)); 832 CeedInt nelem; 833 CeedCallBackend(CeedElemRestrictionGetNumElements(diagrstr, &nelem)); 834 835 // Compute the diagonal of B^T D B 836 int elemsPerBlock = 1; 837 int grid = nelem / elemsPerBlock + ((nelem / elemsPerBlock * elemsPerBlock < nelem) ? 1 : 0); 838 void *args[] = {(void *)&nelem, &diag->d_identity, &diag->d_interpin, &diag->d_gradin, &diag->d_interpout, 839 &diag->d_gradout, &diag->d_emodein, &diag->d_emodeout, &assembledqfarray, &elemdiagarray}; 840 if (pointBlock) { 841 CeedCallBackend(CeedRunKernelDim_Cuda(ceed, diag->linearPointBlock, grid, diag->nnodes, 1, elemsPerBlock, args)); 842 } else { 843 CeedCallBackend(CeedRunKernelDim_Cuda(ceed, diag->linearDiagonal, grid, diag->nnodes, 1, elemsPerBlock, args)); 844 } 845 846 // Restore arrays 847 CeedCallBackend(CeedVectorRestoreArray(elemdiag, &elemdiagarray)); 848 CeedCallBackend(CeedVectorRestoreArrayRead(assembledqf, &assembledqfarray)); 849 850 // Assemble local operator diagonal 851 CeedCallBackend(CeedElemRestrictionApply(diagrstr, CEED_TRANSPOSE, elemdiag, assembled, request)); 852 853 // Cleanup 854 CeedCallBackend(CeedVectorDestroy(&assembledqf)); 855 856 return CEED_ERROR_SUCCESS; 857 } 858 859 //------------------------------------------------------------------------------ 860 // Assemble Linear Diagonal 861 //------------------------------------------------------------------------------ 862 static int CeedOperatorLinearAssembleAddDiagonal_Cuda(CeedOperator op, CeedVector assembled, CeedRequest *request) { 863 CeedCallBackend(CeedOperatorAssembleDiagonalCore_Cuda(op, assembled, request, false)); 864 return CEED_ERROR_SUCCESS; 865 } 866 867 //------------------------------------------------------------------------------ 868 // Assemble Linear Point Block Diagonal 869 //------------------------------------------------------------------------------ 870 static int CeedOperatorLinearAssembleAddPointBlockDiagonal_Cuda(CeedOperator op, CeedVector assembled, CeedRequest *request) { 871 CeedCallBackend(CeedOperatorAssembleDiagonalCore_Cuda(op, assembled, request, true)); 872 return CEED_ERROR_SUCCESS; 873 } 874 875 //------------------------------------------------------------------------------ 876 // Single operator assembly setup 877 //------------------------------------------------------------------------------ 878 static int CeedSingleOperatorAssembleSetup_Cuda(CeedOperator op, CeedInt use_ceedsize_idx) { 879 Ceed ceed; 880 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 881 CeedOperator_Cuda *impl; 882 CeedCallBackend(CeedOperatorGetData(op, &impl)); 883 884 // Get intput and output fields 885 CeedInt num_input_fields, num_output_fields; 886 CeedOperatorField *input_fields; 887 CeedOperatorField *output_fields; 888 CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &input_fields, &num_output_fields, &output_fields)); 889 890 // Determine active input basis eval mode 891 CeedQFunction qf; 892 CeedCallBackend(CeedOperatorGetQFunction(op, &qf)); 893 CeedQFunctionField *qf_fields; 894 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_fields, NULL, NULL)); 895 // Note that the kernel will treat each dimension of a gradient action separately; 896 // i.e., when an active input has a CEED_EVAL_GRAD mode, num_emode_in will increment by dim. 897 // However, for the purposes of loading the B matrices, it will be treated as one mode, and we will load/copy the entire gradient matrix at once, so 898 // num_B_in_mats_to_load will be incremented by 1. 899 CeedInt num_emode_in = 0, dim = 1, num_B_in_mats_to_load = 0, size_B_in = 0; 900 CeedEvalMode *eval_mode_in = NULL; // will be of size num_B_in_mats_load 901 CeedBasis basis_in = NULL; 902 CeedInt nqpts = 0, esize = 0; 903 CeedElemRestriction rstr_in = NULL; 904 for (CeedInt i = 0; i < num_input_fields; i++) { 905 CeedVector vec; 906 CeedCallBackend(CeedOperatorFieldGetVector(input_fields[i], &vec)); 907 if (vec == CEED_VECTOR_ACTIVE) { 908 CeedCallBackend(CeedOperatorFieldGetBasis(input_fields[i], &basis_in)); 909 CeedCallBackend(CeedBasisGetDimension(basis_in, &dim)); 910 CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis_in, &nqpts)); 911 CeedCallBackend(CeedOperatorFieldGetElemRestriction(input_fields[i], &rstr_in)); 912 CeedCallBackend(CeedElemRestrictionGetElementSize(rstr_in, &esize)); 913 CeedEvalMode eval_mode; 914 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode)); 915 if (eval_mode != CEED_EVAL_NONE) { 916 CeedCallBackend(CeedRealloc(num_B_in_mats_to_load + 1, &eval_mode_in)); 917 eval_mode_in[num_B_in_mats_to_load] = eval_mode; 918 num_B_in_mats_to_load += 1; 919 if (eval_mode == CEED_EVAL_GRAD) { 920 num_emode_in += dim; 921 size_B_in += dim * esize * nqpts; 922 } else { 923 num_emode_in += 1; 924 size_B_in += esize * nqpts; 925 } 926 } 927 } 928 } 929 930 // Determine active output basis; basis_out and rstr_out only used if same as input, TODO 931 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qf_fields)); 932 CeedInt num_emode_out = 0, num_B_out_mats_to_load = 0, size_B_out = 0; 933 CeedEvalMode *eval_mode_out = NULL; 934 CeedBasis basis_out = NULL; 935 CeedElemRestriction rstr_out = NULL; 936 for (CeedInt i = 0; i < num_output_fields; i++) { 937 CeedVector vec; 938 CeedCallBackend(CeedOperatorFieldGetVector(output_fields[i], &vec)); 939 if (vec == CEED_VECTOR_ACTIVE) { 940 CeedCallBackend(CeedOperatorFieldGetBasis(output_fields[i], &basis_out)); 941 CeedCallBackend(CeedOperatorFieldGetElemRestriction(output_fields[i], &rstr_out)); 942 CeedCheck(!rstr_out || rstr_out == rstr_in, ceed, CEED_ERROR_BACKEND, "Backend does not implement multi-field non-composite operator assembly"); 943 CeedEvalMode eval_mode; 944 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode)); 945 if (eval_mode != CEED_EVAL_NONE) { 946 CeedCallBackend(CeedRealloc(num_B_out_mats_to_load + 1, &eval_mode_out)); 947 eval_mode_out[num_B_out_mats_to_load] = eval_mode; 948 num_B_out_mats_to_load += 1; 949 if (eval_mode == CEED_EVAL_GRAD) { 950 num_emode_out += dim; 951 size_B_out += dim * esize * nqpts; 952 } else { 953 num_emode_out += 1; 954 size_B_out += esize * nqpts; 955 } 956 } 957 } 958 } 959 CeedCheck(num_emode_in > 0 && num_emode_out > 0, ceed, CEED_ERROR_UNSUPPORTED, "Cannot assemble operator without inputs/outputs"); 960 961 CeedInt nelem, ncomp; 962 CeedCallBackend(CeedElemRestrictionGetNumElements(rstr_in, &nelem)); 963 CeedCallBackend(CeedElemRestrictionGetNumComponents(rstr_in, &ncomp)); 964 965 CeedCallBackend(CeedCalloc(1, &impl->asmb)); 966 CeedOperatorAssemble_Cuda *asmb = impl->asmb; 967 asmb->nelem = nelem; 968 969 // Compile kernels 970 int elemsPerBlock = 1; 971 asmb->elemsPerBlock = elemsPerBlock; 972 CeedInt block_size = esize * esize * elemsPerBlock; 973 Ceed_Cuda *cuda_data; 974 CeedCallBackend(CeedGetData(ceed, &cuda_data)); 975 char *assembly_kernel_path, *assembly_kernel_source; 976 CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-operator-assemble.h", &assembly_kernel_path)); 977 CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Assembly Kernel Source -----\n"); 978 CeedCallBackend(CeedLoadSourceToBuffer(ceed, assembly_kernel_path, &assembly_kernel_source)); 979 CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Assembly Source Complete! -----\n"); 980 bool fallback = block_size > cuda_data->device_prop.maxThreadsPerBlock; 981 if (fallback) { 982 // Use fallback kernel with 1D threadblock 983 block_size = esize * elemsPerBlock; 984 asmb->block_size_x = esize; 985 asmb->block_size_y = 1; 986 } else { // Use kernel with 2D threadblock 987 asmb->block_size_x = esize; 988 asmb->block_size_y = esize; 989 } 990 CeedCallBackend(CeedCompile_Cuda(ceed, assembly_kernel_source, &asmb->module, 8, "NELEM", nelem, "NUMEMODEIN", num_emode_in, "NUMEMODEOUT", 991 num_emode_out, "NQPTS", nqpts, "NNODES", esize, "BLOCK_SIZE", block_size, "NCOMP", ncomp, "CEEDSIZE", 992 use_ceedsize_idx)); 993 CeedCallBackend(CeedGetKernel_Cuda(ceed, asmb->module, fallback ? "linearAssembleFallback" : "linearAssemble", &asmb->linearAssemble)); 994 CeedCallBackend(CeedFree(&assembly_kernel_path)); 995 CeedCallBackend(CeedFree(&assembly_kernel_source)); 996 997 // Build 'full' B matrices (not 1D arrays used for tensor-product matrices) 998 const CeedScalar *interp_in, *grad_in; 999 CeedCallBackend(CeedBasisGetInterp(basis_in, &interp_in)); 1000 CeedCallBackend(CeedBasisGetGrad(basis_in, &grad_in)); 1001 1002 // Load into B_in, in order that they will be used in eval_mode 1003 const CeedInt inBytes = size_B_in * sizeof(CeedScalar); 1004 CeedInt mat_start = 0; 1005 CeedCallCuda(ceed, cudaMalloc((void **)&asmb->d_B_in, inBytes)); 1006 for (int i = 0; i < num_B_in_mats_to_load; i++) { 1007 CeedEvalMode eval_mode = eval_mode_in[i]; 1008 if (eval_mode == CEED_EVAL_INTERP) { 1009 CeedCallCuda(ceed, cudaMemcpy(&asmb->d_B_in[mat_start], interp_in, esize * nqpts * sizeof(CeedScalar), cudaMemcpyHostToDevice)); 1010 mat_start += esize * nqpts; 1011 } else if (eval_mode == CEED_EVAL_GRAD) { 1012 CeedCallCuda(ceed, cudaMemcpy(&asmb->d_B_in[mat_start], grad_in, dim * esize * nqpts * sizeof(CeedScalar), cudaMemcpyHostToDevice)); 1013 mat_start += dim * esize * nqpts; 1014 } 1015 } 1016 1017 const CeedScalar *interp_out, *grad_out; 1018 // Note that this function currently assumes 1 basis, so this should always be true 1019 // for now 1020 if (basis_out == basis_in) { 1021 interp_out = interp_in; 1022 grad_out = grad_in; 1023 } else { 1024 CeedCallBackend(CeedBasisGetInterp(basis_out, &interp_out)); 1025 CeedCallBackend(CeedBasisGetGrad(basis_out, &grad_out)); 1026 } 1027 1028 // Load into B_out, in order that they will be used in eval_mode 1029 const CeedInt outBytes = size_B_out * sizeof(CeedScalar); 1030 mat_start = 0; 1031 CeedCallCuda(ceed, cudaMalloc((void **)&asmb->d_B_out, outBytes)); 1032 for (int i = 0; i < num_B_out_mats_to_load; i++) { 1033 CeedEvalMode eval_mode = eval_mode_out[i]; 1034 if (eval_mode == CEED_EVAL_INTERP) { 1035 CeedCallCuda(ceed, cudaMemcpy(&asmb->d_B_out[mat_start], interp_out, esize * nqpts * sizeof(CeedScalar), cudaMemcpyHostToDevice)); 1036 mat_start += esize * nqpts; 1037 } else if (eval_mode == CEED_EVAL_GRAD) { 1038 CeedCallCuda(ceed, cudaMemcpy(&asmb->d_B_out[mat_start], grad_out, dim * esize * nqpts * sizeof(CeedScalar), cudaMemcpyHostToDevice)); 1039 mat_start += dim * esize * nqpts; 1040 } 1041 } 1042 return CEED_ERROR_SUCCESS; 1043 } 1044 1045 //------------------------------------------------------------------------------ 1046 // Assemble matrix data for COO matrix of assembled operator. 1047 // The sparsity pattern is set by CeedOperatorLinearAssembleSymbolic. 1048 // 1049 // Note that this (and other assembly routines) currently assume only one active input restriction/basis per operator (could have multiple basis eval 1050 // modes). 1051 // TODO: allow multiple active input restrictions/basis objects 1052 //------------------------------------------------------------------------------ 1053 static int CeedSingleOperatorAssemble_Cuda(CeedOperator op, CeedInt offset, CeedVector values) { 1054 Ceed ceed; 1055 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 1056 CeedOperator_Cuda *impl; 1057 CeedCallBackend(CeedOperatorGetData(op, &impl)); 1058 1059 // Assemble QFunction 1060 CeedVector assembled_qf = NULL; 1061 CeedElemRestriction rstr_q = NULL; 1062 CeedCallBackend(CeedOperatorLinearAssembleQFunctionBuildOrUpdate(op, &assembled_qf, &rstr_q, CEED_REQUEST_IMMEDIATE)); 1063 CeedCallBackend(CeedElemRestrictionDestroy(&rstr_q)); 1064 CeedScalar *values_array; 1065 CeedCallBackend(CeedVectorGetArray(values, CEED_MEM_DEVICE, &values_array)); 1066 values_array += offset; 1067 const CeedScalar *qf_array; 1068 CeedCallBackend(CeedVectorGetArrayRead(assembled_qf, CEED_MEM_DEVICE, &qf_array)); 1069 1070 CeedSize values_length = 0, assembled_qf_length = 0; 1071 CeedCallBackend(CeedVectorGetLength(values, &values_length)); 1072 CeedCallBackend(CeedVectorGetLength(assembled_qf, &assembled_qf_length)); 1073 CeedInt use_ceedsize_idx = 0; 1074 if ((values_length > INT_MAX) || (assembled_qf_length > INT_MAX)) use_ceedsize_idx = 1; 1075 // Setup 1076 if (!impl->asmb) { 1077 CeedCallBackend(CeedSingleOperatorAssembleSetup_Cuda(op, use_ceedsize_idx)); 1078 assert(impl->asmb != NULL); 1079 } 1080 1081 // Compute B^T D B 1082 const CeedInt nelem = impl->asmb->nelem; 1083 const CeedInt elemsPerBlock = impl->asmb->elemsPerBlock; 1084 const CeedInt grid = nelem / elemsPerBlock + ((nelem / elemsPerBlock * elemsPerBlock < nelem) ? 1 : 0); 1085 void *args[] = {&impl->asmb->d_B_in, &impl->asmb->d_B_out, &qf_array, &values_array}; 1086 CeedCallBackend( 1087 CeedRunKernelDim_Cuda(ceed, impl->asmb->linearAssemble, grid, impl->asmb->block_size_x, impl->asmb->block_size_y, elemsPerBlock, args)); 1088 1089 // Restore arrays 1090 CeedCallBackend(CeedVectorRestoreArray(values, &values_array)); 1091 CeedCallBackend(CeedVectorRestoreArrayRead(assembled_qf, &qf_array)); 1092 1093 // Cleanup 1094 CeedCallBackend(CeedVectorDestroy(&assembled_qf)); 1095 1096 return CEED_ERROR_SUCCESS; 1097 } 1098 1099 //------------------------------------------------------------------------------ 1100 // Create operator 1101 //------------------------------------------------------------------------------ 1102 int CeedOperatorCreate_Cuda(CeedOperator op) { 1103 Ceed ceed; 1104 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 1105 CeedOperator_Cuda *impl; 1106 1107 CeedCallBackend(CeedCalloc(1, &impl)); 1108 CeedCallBackend(CeedOperatorSetData(op, impl)); 1109 1110 CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleQFunction", CeedOperatorLinearAssembleQFunction_Cuda)); 1111 CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleQFunctionUpdate", CeedOperatorLinearAssembleQFunctionUpdate_Cuda)); 1112 CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleAddDiagonal", CeedOperatorLinearAssembleAddDiagonal_Cuda)); 1113 CeedCallBackend( 1114 CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleAddPointBlockDiagonal", CeedOperatorLinearAssembleAddPointBlockDiagonal_Cuda)); 1115 CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleSingle", CeedSingleOperatorAssemble_Cuda)); 1116 CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "ApplyAdd", CeedOperatorApplyAdd_Cuda)); 1117 CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "Destroy", CeedOperatorDestroy_Cuda)); 1118 return CEED_ERROR_SUCCESS; 1119 } 1120 1121 //------------------------------------------------------------------------------ 1122