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