10d0321e0SJeremy L Thompson // Copyright (c) 2017-2018, Lawrence Livermore National Security, LLC. 20d0321e0SJeremy L Thompson // Produced at the Lawrence Livermore National Laboratory. LLNL-CODE-734707. 30d0321e0SJeremy L Thompson // All Rights reserved. See files LICENSE and NOTICE for details. 40d0321e0SJeremy L Thompson // 50d0321e0SJeremy L Thompson // This file is part of CEED, a collection of benchmarks, miniapps, software 60d0321e0SJeremy L Thompson // libraries and APIs for efficient high-order finite element and spectral 70d0321e0SJeremy L Thompson // element discretizations for exascale applications. For more information and 80d0321e0SJeremy L Thompson // source code availability see http://github.com/ceed. 90d0321e0SJeremy L Thompson // 100d0321e0SJeremy L Thompson // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC, 110d0321e0SJeremy L Thompson // a collaborative effort of two U.S. Department of Energy organizations (Office 120d0321e0SJeremy L Thompson // of Science and the National Nuclear Security Administration) responsible for 130d0321e0SJeremy L Thompson // the planning and preparation of a capable exascale ecosystem, including 140d0321e0SJeremy L Thompson // software, applications, hardware, advanced system engineering and early 150d0321e0SJeremy L Thompson // testbed platforms, in support of the nation's exascale computing imperative. 160d0321e0SJeremy L Thompson 170d0321e0SJeremy L Thompson #include <ceed/ceed.h> 180d0321e0SJeremy L Thompson #include <ceed/backend.h> 190d0321e0SJeremy L Thompson #include <hip/hip_runtime.h> 200d0321e0SJeremy L Thompson #include <assert.h> 210d0321e0SJeremy L Thompson #include <stdbool.h> 220d0321e0SJeremy L Thompson #include <string.h> 230d0321e0SJeremy L Thompson #include "ceed-hip-ref.h" 240d0321e0SJeremy L Thompson #include "../hip/ceed-hip-compile.h" 250d0321e0SJeremy L Thompson 260d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 270d0321e0SJeremy L Thompson // Destroy operator 280d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 290d0321e0SJeremy L Thompson static int CeedOperatorDestroy_Hip(CeedOperator op) { 300d0321e0SJeremy L Thompson int ierr; 310d0321e0SJeremy L Thompson CeedOperator_Hip *impl; 320d0321e0SJeremy L Thompson ierr = CeedOperatorGetData(op, &impl); CeedChkBackend(ierr); 330d0321e0SJeremy L Thompson 340d0321e0SJeremy L Thompson // Apply data 350d0321e0SJeremy L Thompson for (CeedInt i = 0; i < impl->numein + impl->numeout; i++) { 360d0321e0SJeremy L Thompson ierr = CeedVectorDestroy(&impl->evecs[i]); CeedChkBackend(ierr); 370d0321e0SJeremy L Thompson } 380d0321e0SJeremy L Thompson ierr = CeedFree(&impl->evecs); CeedChkBackend(ierr); 390d0321e0SJeremy L Thompson 400d0321e0SJeremy L Thompson for (CeedInt i = 0; i < impl->numein; i++) { 410d0321e0SJeremy L Thompson ierr = CeedVectorDestroy(&impl->qvecsin[i]); CeedChkBackend(ierr); 420d0321e0SJeremy L Thompson } 430d0321e0SJeremy L Thompson ierr = CeedFree(&impl->qvecsin); CeedChkBackend(ierr); 440d0321e0SJeremy L Thompson 450d0321e0SJeremy L Thompson for (CeedInt i = 0; i < impl->numeout; i++) { 460d0321e0SJeremy L Thompson ierr = CeedVectorDestroy(&impl->qvecsout[i]); CeedChkBackend(ierr); 470d0321e0SJeremy L Thompson } 480d0321e0SJeremy L Thompson ierr = CeedFree(&impl->qvecsout); CeedChkBackend(ierr); 490d0321e0SJeremy L Thompson 500d0321e0SJeremy L Thompson // QFunction diagonal assembly data 510d0321e0SJeremy L Thompson for (CeedInt i=0; i<impl->qfnumactivein; i++) { 520d0321e0SJeremy L Thompson ierr = CeedVectorDestroy(&impl->qfactivein[i]); CeedChkBackend(ierr); 530d0321e0SJeremy L Thompson } 540d0321e0SJeremy L Thompson ierr = CeedFree(&impl->qfactivein); CeedChkBackend(ierr); 550d0321e0SJeremy L Thompson 560d0321e0SJeremy L Thompson // Diag data 570d0321e0SJeremy L Thompson if (impl->diag) { 580d0321e0SJeremy L Thompson Ceed ceed; 590d0321e0SJeremy L Thompson ierr = CeedOperatorGetCeed(op, &ceed); CeedChkBackend(ierr); 600d0321e0SJeremy L Thompson CeedChk_Hip(ceed, hipModuleUnload(impl->diag->module)); 610d0321e0SJeremy L Thompson ierr = CeedFree(&impl->diag->h_emodein); CeedChkBackend(ierr); 620d0321e0SJeremy L Thompson ierr = CeedFree(&impl->diag->h_emodeout); CeedChkBackend(ierr); 630d0321e0SJeremy L Thompson ierr = hipFree(impl->diag->d_emodein); CeedChk_Hip(ceed, ierr); 640d0321e0SJeremy L Thompson ierr = hipFree(impl->diag->d_emodeout); CeedChk_Hip(ceed, ierr); 650d0321e0SJeremy L Thompson ierr = hipFree(impl->diag->d_identity); CeedChk_Hip(ceed, ierr); 660d0321e0SJeremy L Thompson ierr = hipFree(impl->diag->d_interpin); CeedChk_Hip(ceed, ierr); 670d0321e0SJeremy L Thompson ierr = hipFree(impl->diag->d_interpout); CeedChk_Hip(ceed, ierr); 680d0321e0SJeremy L Thompson ierr = hipFree(impl->diag->d_gradin); CeedChk_Hip(ceed, ierr); 690d0321e0SJeremy L Thompson ierr = hipFree(impl->diag->d_gradout); CeedChk_Hip(ceed, ierr); 700d0321e0SJeremy L Thompson ierr = CeedElemRestrictionDestroy(&impl->diag->pbdiagrstr); 710d0321e0SJeremy L Thompson CeedChkBackend(ierr); 720d0321e0SJeremy L Thompson ierr = CeedVectorDestroy(&impl->diag->elemdiag); CeedChkBackend(ierr); 730d0321e0SJeremy L Thompson ierr = CeedVectorDestroy(&impl->diag->pbelemdiag); CeedChkBackend(ierr); 740d0321e0SJeremy L Thompson } 750d0321e0SJeremy L Thompson ierr = CeedFree(&impl->diag); CeedChkBackend(ierr); 760d0321e0SJeremy L Thompson 77*a835093fSnbeams if (impl->asmb) { 78*a835093fSnbeams Ceed ceed; 79*a835093fSnbeams ierr = CeedOperatorGetCeed(op, &ceed); CeedChkBackend(ierr); 80*a835093fSnbeams CeedChk_Hip(ceed, hipModuleUnload(impl->asmb->module)); 81*a835093fSnbeams ierr = hipFree(impl->asmb->d_B_in); CeedChk_Hip(ceed, ierr); 82*a835093fSnbeams ierr = hipFree(impl->asmb->d_B_out); CeedChk_Hip(ceed, ierr); 83*a835093fSnbeams } 84*a835093fSnbeams ierr = CeedFree(&impl->asmb); CeedChkBackend(ierr); 85*a835093fSnbeams 860d0321e0SJeremy L Thompson ierr = CeedFree(&impl); CeedChkBackend(ierr); 870d0321e0SJeremy L Thompson return CEED_ERROR_SUCCESS; 880d0321e0SJeremy L Thompson } 890d0321e0SJeremy L Thompson 900d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 910d0321e0SJeremy L Thompson // Setup infields or outfields 920d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 930d0321e0SJeremy L Thompson static int CeedOperatorSetupFields_Hip(CeedQFunction qf, CeedOperator op, 940d0321e0SJeremy L Thompson bool isinput, CeedVector *evecs, 950d0321e0SJeremy L Thompson CeedVector *qvecs, CeedInt starte, 960d0321e0SJeremy L Thompson CeedInt numfields, CeedInt Q, 970d0321e0SJeremy L Thompson CeedInt numelements) { 980d0321e0SJeremy L Thompson CeedInt dim, ierr, size; 990d0321e0SJeremy L Thompson Ceed ceed; 1000d0321e0SJeremy L Thompson ierr = CeedOperatorGetCeed(op, &ceed); CeedChkBackend(ierr); 1010d0321e0SJeremy L Thompson CeedBasis basis; 1020d0321e0SJeremy L Thompson CeedElemRestriction Erestrict; 1030d0321e0SJeremy L Thompson CeedOperatorField *opfields; 1040d0321e0SJeremy L Thompson CeedQFunctionField *qffields; 1050d0321e0SJeremy L Thompson CeedVector fieldvec; 1060d0321e0SJeremy L Thompson bool strided; 1070d0321e0SJeremy L Thompson bool skiprestrict; 1080d0321e0SJeremy L Thompson 1090d0321e0SJeremy L Thompson if (isinput) { 1100d0321e0SJeremy L Thompson ierr = CeedOperatorGetFields(op, NULL, &opfields, NULL, NULL); 1110d0321e0SJeremy L Thompson CeedChkBackend(ierr); 1120d0321e0SJeremy L Thompson ierr = CeedQFunctionGetFields(qf, NULL, &qffields, NULL, NULL); 1130d0321e0SJeremy L Thompson CeedChkBackend(ierr); 1140d0321e0SJeremy L Thompson } else { 1150d0321e0SJeremy L Thompson ierr = CeedOperatorGetFields(op, NULL, NULL, NULL, &opfields); 1160d0321e0SJeremy L Thompson CeedChkBackend(ierr); 1170d0321e0SJeremy L Thompson ierr = CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qffields); 1180d0321e0SJeremy L Thompson CeedChkBackend(ierr); 1190d0321e0SJeremy L Thompson } 1200d0321e0SJeremy L Thompson 1210d0321e0SJeremy L Thompson // Loop over fields 1220d0321e0SJeremy L Thompson for (CeedInt i = 0; i < numfields; i++) { 1230d0321e0SJeremy L Thompson CeedEvalMode emode; 1240d0321e0SJeremy L Thompson ierr = CeedQFunctionFieldGetEvalMode(qffields[i], &emode); CeedChkBackend(ierr); 1250d0321e0SJeremy L Thompson 1260d0321e0SJeremy L Thompson strided = false; 1270d0321e0SJeremy L Thompson skiprestrict = false; 1280d0321e0SJeremy L Thompson if (emode != CEED_EVAL_WEIGHT) { 1290d0321e0SJeremy L Thompson ierr = CeedOperatorFieldGetElemRestriction(opfields[i], &Erestrict); 1300d0321e0SJeremy L Thompson CeedChkBackend(ierr); 1310d0321e0SJeremy L Thompson 1320d0321e0SJeremy L Thompson // Check whether this field can skip the element restriction: 1330d0321e0SJeremy L Thompson // must be passive input, with emode NONE, and have a strided restriction with 1340d0321e0SJeremy L Thompson // CEED_STRIDES_BACKEND. 1350d0321e0SJeremy L Thompson 1360d0321e0SJeremy L Thompson // First, check whether the field is input or output: 1370d0321e0SJeremy L Thompson if (isinput) { 1380d0321e0SJeremy L Thompson // Check for passive input: 1390d0321e0SJeremy L Thompson ierr = CeedOperatorFieldGetVector(opfields[i], &fieldvec); CeedChkBackend(ierr); 1400d0321e0SJeremy L Thompson if (fieldvec != CEED_VECTOR_ACTIVE) { 1410d0321e0SJeremy L Thompson // Check emode 1420d0321e0SJeremy L Thompson if (emode == CEED_EVAL_NONE) { 1430d0321e0SJeremy L Thompson // Check for strided restriction 1440d0321e0SJeremy L Thompson ierr = CeedElemRestrictionIsStrided(Erestrict, &strided); 1450d0321e0SJeremy L Thompson CeedChkBackend(ierr); 1460d0321e0SJeremy L Thompson if (strided) { 1470d0321e0SJeremy L Thompson // Check if vector is already in preferred backend ordering 1480d0321e0SJeremy L Thompson ierr = CeedElemRestrictionHasBackendStrides(Erestrict, 1490d0321e0SJeremy L Thompson &skiprestrict); CeedChkBackend(ierr); 1500d0321e0SJeremy L Thompson } 1510d0321e0SJeremy L Thompson } 1520d0321e0SJeremy L Thompson } 1530d0321e0SJeremy L Thompson } 1540d0321e0SJeremy L Thompson if (skiprestrict) { 1550d0321e0SJeremy L Thompson // We do not need an E-Vector, but will use the input field vector's data 1560d0321e0SJeremy L Thompson // directly in the operator application. 1570d0321e0SJeremy L Thompson evecs[i + starte] = NULL; 1580d0321e0SJeremy L Thompson } else { 1590d0321e0SJeremy L Thompson ierr = CeedElemRestrictionCreateVector(Erestrict, NULL, 1600d0321e0SJeremy L Thompson &evecs[i + starte]); 1610d0321e0SJeremy L Thompson CeedChkBackend(ierr); 1620d0321e0SJeremy L Thompson } 1630d0321e0SJeremy L Thompson } 1640d0321e0SJeremy L Thompson 1650d0321e0SJeremy L Thompson switch (emode) { 1660d0321e0SJeremy L Thompson case CEED_EVAL_NONE: 1670d0321e0SJeremy L Thompson ierr = CeedQFunctionFieldGetSize(qffields[i], &size); CeedChkBackend(ierr); 1680d0321e0SJeremy L Thompson ierr = CeedVectorCreate(ceed, numelements * Q * size, &qvecs[i]); 1690d0321e0SJeremy L Thompson CeedChkBackend(ierr); 1700d0321e0SJeremy L Thompson break; 1710d0321e0SJeremy L Thompson case CEED_EVAL_INTERP: 1720d0321e0SJeremy L Thompson ierr = CeedQFunctionFieldGetSize(qffields[i], &size); CeedChkBackend(ierr); 1730d0321e0SJeremy L Thompson ierr = CeedVectorCreate(ceed, numelements * Q * size, &qvecs[i]); 1740d0321e0SJeremy L Thompson CeedChkBackend(ierr); 1750d0321e0SJeremy L Thompson break; 1760d0321e0SJeremy L Thompson case CEED_EVAL_GRAD: 1770d0321e0SJeremy L Thompson ierr = CeedOperatorFieldGetBasis(opfields[i], &basis); CeedChkBackend(ierr); 1780d0321e0SJeremy L Thompson ierr = CeedQFunctionFieldGetSize(qffields[i], &size); CeedChkBackend(ierr); 1790d0321e0SJeremy L Thompson ierr = CeedBasisGetDimension(basis, &dim); CeedChkBackend(ierr); 1800d0321e0SJeremy L Thompson ierr = CeedVectorCreate(ceed, numelements * Q * size, &qvecs[i]); 1810d0321e0SJeremy L Thompson CeedChkBackend(ierr); 1820d0321e0SJeremy L Thompson break; 1830d0321e0SJeremy L Thompson case CEED_EVAL_WEIGHT: // Only on input fields 1840d0321e0SJeremy L Thompson ierr = CeedOperatorFieldGetBasis(opfields[i], &basis); CeedChkBackend(ierr); 1850d0321e0SJeremy L Thompson ierr = CeedVectorCreate(ceed, numelements * Q, &qvecs[i]); CeedChkBackend(ierr); 1860d0321e0SJeremy L Thompson ierr = CeedBasisApply(basis, numelements, CEED_NOTRANSPOSE, 1870d0321e0SJeremy L Thompson CEED_EVAL_WEIGHT, NULL, qvecs[i]); CeedChkBackend(ierr); 1880d0321e0SJeremy L Thompson break; 1890d0321e0SJeremy L Thompson case CEED_EVAL_DIV: 1900d0321e0SJeremy L Thompson break; // TODO: Not implemented 1910d0321e0SJeremy L Thompson case CEED_EVAL_CURL: 1920d0321e0SJeremy L Thompson break; // TODO: Not implemented 1930d0321e0SJeremy L Thompson } 1940d0321e0SJeremy L Thompson } 1950d0321e0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1960d0321e0SJeremy L Thompson } 1970d0321e0SJeremy L Thompson 1980d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 1990d0321e0SJeremy L Thompson // CeedOperator needs to connect all the named fields (be they active or passive) 2000d0321e0SJeremy L Thompson // to the named inputs and outputs of its CeedQFunction. 2010d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 2020d0321e0SJeremy L Thompson static int CeedOperatorSetup_Hip(CeedOperator op) { 2030d0321e0SJeremy L Thompson int ierr; 2040d0321e0SJeremy L Thompson bool setupdone; 2050d0321e0SJeremy L Thompson ierr = CeedOperatorIsSetupDone(op, &setupdone); CeedChkBackend(ierr); 2060d0321e0SJeremy L Thompson if (setupdone) 2070d0321e0SJeremy L Thompson return CEED_ERROR_SUCCESS; 2080d0321e0SJeremy L Thompson Ceed ceed; 2090d0321e0SJeremy L Thompson ierr = CeedOperatorGetCeed(op, &ceed); CeedChkBackend(ierr); 2100d0321e0SJeremy L Thompson CeedOperator_Hip *impl; 2110d0321e0SJeremy L Thompson ierr = CeedOperatorGetData(op, &impl); CeedChkBackend(ierr); 2120d0321e0SJeremy L Thompson CeedQFunction qf; 2130d0321e0SJeremy L Thompson ierr = CeedOperatorGetQFunction(op, &qf); CeedChkBackend(ierr); 2140d0321e0SJeremy L Thompson CeedInt Q, numelements, numinputfields, numoutputfields; 2150d0321e0SJeremy L Thompson ierr = CeedOperatorGetNumQuadraturePoints(op, &Q); CeedChkBackend(ierr); 2160d0321e0SJeremy L Thompson ierr = CeedOperatorGetNumElements(op, &numelements); CeedChkBackend(ierr); 2170d0321e0SJeremy L Thompson CeedOperatorField *opinputfields, *opoutputfields; 2180d0321e0SJeremy L Thompson ierr = CeedOperatorGetFields(op, &numinputfields, &opinputfields, 2190d0321e0SJeremy L Thompson &numoutputfields, &opoutputfields); 2200d0321e0SJeremy L Thompson CeedChkBackend(ierr); 2210d0321e0SJeremy L Thompson CeedQFunctionField *qfinputfields, *qfoutputfields; 2220d0321e0SJeremy L Thompson ierr = CeedQFunctionGetFields(qf, NULL, &qfinputfields, NULL, &qfoutputfields); 2230d0321e0SJeremy L Thompson CeedChkBackend(ierr); 2240d0321e0SJeremy L Thompson 2250d0321e0SJeremy L Thompson // Allocate 2260d0321e0SJeremy L Thompson ierr = CeedCalloc(numinputfields + numoutputfields, &impl->evecs); 2270d0321e0SJeremy L Thompson CeedChkBackend(ierr); 2280d0321e0SJeremy L Thompson 2290d0321e0SJeremy L Thompson ierr = CeedCalloc(CEED_FIELD_MAX, &impl->qvecsin); CeedChkBackend(ierr); 2300d0321e0SJeremy L Thompson ierr = CeedCalloc(CEED_FIELD_MAX, &impl->qvecsout); CeedChkBackend(ierr); 2310d0321e0SJeremy L Thompson 2320d0321e0SJeremy L Thompson impl->numein = numinputfields; impl->numeout = numoutputfields; 2330d0321e0SJeremy L Thompson 2340d0321e0SJeremy L Thompson // Set up infield and outfield evecs and qvecs 2350d0321e0SJeremy L Thompson // Infields 2360d0321e0SJeremy L Thompson ierr = CeedOperatorSetupFields_Hip(qf, op, true, 2370d0321e0SJeremy L Thompson impl->evecs, impl->qvecsin, 0, 2380d0321e0SJeremy L Thompson numinputfields, Q, numelements); 2390d0321e0SJeremy L Thompson CeedChkBackend(ierr); 2400d0321e0SJeremy L Thompson 2410d0321e0SJeremy L Thompson // Outfields 2420d0321e0SJeremy L Thompson ierr = CeedOperatorSetupFields_Hip(qf, op, false, 2430d0321e0SJeremy L Thompson impl->evecs, impl->qvecsout, 2440d0321e0SJeremy L Thompson numinputfields, numoutputfields, Q, 2450d0321e0SJeremy L Thompson numelements); CeedChkBackend(ierr); 2460d0321e0SJeremy L Thompson 2470d0321e0SJeremy L Thompson ierr = CeedOperatorSetSetupDone(op); CeedChkBackend(ierr); 2480d0321e0SJeremy L Thompson return CEED_ERROR_SUCCESS; 2490d0321e0SJeremy L Thompson } 2500d0321e0SJeremy L Thompson 2510d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 2520d0321e0SJeremy L Thompson // Setup Operator Inputs 2530d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 2540d0321e0SJeremy L Thompson static inline int CeedOperatorSetupInputs_Hip(CeedInt numinputfields, 2550d0321e0SJeremy L Thompson CeedQFunctionField *qfinputfields, CeedOperatorField *opinputfields, 2560d0321e0SJeremy L Thompson CeedVector invec, const bool skipactive, CeedScalar *edata[2*CEED_FIELD_MAX], 2570d0321e0SJeremy L Thompson CeedOperator_Hip *impl, CeedRequest *request) { 2580d0321e0SJeremy L Thompson CeedInt ierr; 2590d0321e0SJeremy L Thompson CeedEvalMode emode; 2600d0321e0SJeremy L Thompson CeedVector vec; 2610d0321e0SJeremy L Thompson CeedElemRestriction Erestrict; 2620d0321e0SJeremy L Thompson 2630d0321e0SJeremy L Thompson for (CeedInt i = 0; i < numinputfields; i++) { 2640d0321e0SJeremy L Thompson // Get input vector 2650d0321e0SJeremy L Thompson ierr = CeedOperatorFieldGetVector(opinputfields[i], &vec); CeedChkBackend(ierr); 2660d0321e0SJeremy L Thompson if (vec == CEED_VECTOR_ACTIVE) { 2670d0321e0SJeremy L Thompson if (skipactive) 2680d0321e0SJeremy L Thompson continue; 2690d0321e0SJeremy L Thompson else 2700d0321e0SJeremy L Thompson vec = invec; 2710d0321e0SJeremy L Thompson } 2720d0321e0SJeremy L Thompson 2730d0321e0SJeremy L Thompson ierr = CeedQFunctionFieldGetEvalMode(qfinputfields[i], &emode); 2740d0321e0SJeremy L Thompson CeedChkBackend(ierr); 2750d0321e0SJeremy L Thompson if (emode == CEED_EVAL_WEIGHT) { // Skip 2760d0321e0SJeremy L Thompson } else { 2770d0321e0SJeremy L Thompson // Get input vector 2780d0321e0SJeremy L Thompson ierr = CeedOperatorFieldGetVector(opinputfields[i], &vec); CeedChkBackend(ierr); 2790d0321e0SJeremy L Thompson // Get input element restriction 2800d0321e0SJeremy L Thompson ierr = CeedOperatorFieldGetElemRestriction(opinputfields[i], &Erestrict); 2810d0321e0SJeremy L Thompson CeedChkBackend(ierr); 2820d0321e0SJeremy L Thompson if (vec == CEED_VECTOR_ACTIVE) 2830d0321e0SJeremy L Thompson vec = invec; 2840d0321e0SJeremy L Thompson // Restrict, if necessary 2850d0321e0SJeremy L Thompson if (!impl->evecs[i]) { 2860d0321e0SJeremy L Thompson // No restriction for this field; read data directly from vec. 2870d0321e0SJeremy L Thompson ierr = CeedVectorGetArrayRead(vec, CEED_MEM_DEVICE, 2880d0321e0SJeremy L Thompson (const CeedScalar **) &edata[i]); 2890d0321e0SJeremy L Thompson CeedChkBackend(ierr); 2900d0321e0SJeremy L Thompson } else { 2910d0321e0SJeremy L Thompson ierr = CeedElemRestrictionApply(Erestrict, CEED_NOTRANSPOSE, vec, 2920d0321e0SJeremy L Thompson impl->evecs[i], request); CeedChkBackend(ierr); 2930d0321e0SJeremy L Thompson // Get evec 2940d0321e0SJeremy L Thompson ierr = CeedVectorGetArrayRead(impl->evecs[i], CEED_MEM_DEVICE, 2950d0321e0SJeremy L Thompson (const CeedScalar **) &edata[i]); 2960d0321e0SJeremy L Thompson CeedChkBackend(ierr); 2970d0321e0SJeremy L Thompson } 2980d0321e0SJeremy L Thompson } 2990d0321e0SJeremy L Thompson } 3000d0321e0SJeremy L Thompson return CEED_ERROR_SUCCESS; 3010d0321e0SJeremy L Thompson } 3020d0321e0SJeremy L Thompson 3030d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 3040d0321e0SJeremy L Thompson // Input Basis Action 3050d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 3060d0321e0SJeremy L Thompson static inline int CeedOperatorInputBasis_Hip(CeedInt numelements, 3070d0321e0SJeremy L Thompson CeedQFunctionField *qfinputfields, CeedOperatorField *opinputfields, 3080d0321e0SJeremy L Thompson CeedInt numinputfields, const bool skipactive, 3090d0321e0SJeremy L Thompson CeedScalar *edata[2*CEED_FIELD_MAX], CeedOperator_Hip *impl) { 3100d0321e0SJeremy L Thompson CeedInt ierr; 3110d0321e0SJeremy L Thompson CeedInt elemsize, size; 3120d0321e0SJeremy L Thompson CeedElemRestriction Erestrict; 3130d0321e0SJeremy L Thompson CeedEvalMode emode; 3140d0321e0SJeremy L Thompson CeedBasis basis; 3150d0321e0SJeremy L Thompson 3160d0321e0SJeremy L Thompson for (CeedInt i=0; i<numinputfields; i++) { 3170d0321e0SJeremy L Thompson // Skip active input 3180d0321e0SJeremy L Thompson if (skipactive) { 3190d0321e0SJeremy L Thompson CeedVector vec; 3200d0321e0SJeremy L Thompson ierr = CeedOperatorFieldGetVector(opinputfields[i], &vec); CeedChkBackend(ierr); 3210d0321e0SJeremy L Thompson if (vec == CEED_VECTOR_ACTIVE) 3220d0321e0SJeremy L Thompson continue; 3230d0321e0SJeremy L Thompson } 3240d0321e0SJeremy L Thompson // Get elemsize, emode, size 3250d0321e0SJeremy L Thompson ierr = CeedOperatorFieldGetElemRestriction(opinputfields[i], &Erestrict); 3260d0321e0SJeremy L Thompson CeedChkBackend(ierr); 3270d0321e0SJeremy L Thompson ierr = CeedElemRestrictionGetElementSize(Erestrict, &elemsize); 3280d0321e0SJeremy L Thompson CeedChkBackend(ierr); 3290d0321e0SJeremy L Thompson ierr = CeedQFunctionFieldGetEvalMode(qfinputfields[i], &emode); 3300d0321e0SJeremy L Thompson CeedChkBackend(ierr); 3310d0321e0SJeremy L Thompson ierr = CeedQFunctionFieldGetSize(qfinputfields[i], &size); CeedChkBackend(ierr); 3320d0321e0SJeremy L Thompson // Basis action 3330d0321e0SJeremy L Thompson switch (emode) { 3340d0321e0SJeremy L Thompson case CEED_EVAL_NONE: 3350d0321e0SJeremy L Thompson ierr = CeedVectorSetArray(impl->qvecsin[i], CEED_MEM_DEVICE, 3360d0321e0SJeremy L Thompson CEED_USE_POINTER, edata[i]); CeedChkBackend(ierr); 3370d0321e0SJeremy L Thompson break; 3380d0321e0SJeremy L Thompson case CEED_EVAL_INTERP: 3390d0321e0SJeremy L Thompson ierr = CeedOperatorFieldGetBasis(opinputfields[i], &basis); 3400d0321e0SJeremy L Thompson CeedChkBackend(ierr); 3410d0321e0SJeremy L Thompson ierr = CeedBasisApply(basis, numelements, CEED_NOTRANSPOSE, 3420d0321e0SJeremy L Thompson CEED_EVAL_INTERP, impl->evecs[i], 3430d0321e0SJeremy L Thompson impl->qvecsin[i]); CeedChkBackend(ierr); 3440d0321e0SJeremy L Thompson break; 3450d0321e0SJeremy L Thompson case CEED_EVAL_GRAD: 3460d0321e0SJeremy L Thompson ierr = CeedOperatorFieldGetBasis(opinputfields[i], &basis); 3470d0321e0SJeremy L Thompson CeedChkBackend(ierr); 3480d0321e0SJeremy L Thompson ierr = CeedBasisApply(basis, numelements, CEED_NOTRANSPOSE, 3490d0321e0SJeremy L Thompson CEED_EVAL_GRAD, impl->evecs[i], 3500d0321e0SJeremy L Thompson impl->qvecsin[i]); CeedChkBackend(ierr); 3510d0321e0SJeremy L Thompson break; 3520d0321e0SJeremy L Thompson case CEED_EVAL_WEIGHT: 3530d0321e0SJeremy L Thompson break; // No action 3540d0321e0SJeremy L Thompson case CEED_EVAL_DIV: 3550d0321e0SJeremy L Thompson break; // TODO: Not implemented 3560d0321e0SJeremy L Thompson case CEED_EVAL_CURL: 3570d0321e0SJeremy L Thompson break; // TODO: Not implemented 3580d0321e0SJeremy L Thompson } 3590d0321e0SJeremy L Thompson } 3600d0321e0SJeremy L Thompson return CEED_ERROR_SUCCESS; 3610d0321e0SJeremy L Thompson } 3620d0321e0SJeremy L Thompson 3630d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 3640d0321e0SJeremy L Thompson // Restore Input Vectors 3650d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 3660d0321e0SJeremy L Thompson static inline int CeedOperatorRestoreInputs_Hip(CeedInt numinputfields, 3670d0321e0SJeremy L Thompson CeedQFunctionField *qfinputfields, CeedOperatorField *opinputfields, 3680d0321e0SJeremy L Thompson const bool skipactive, CeedScalar *edata[2*CEED_FIELD_MAX], 3690d0321e0SJeremy L Thompson CeedOperator_Hip *impl) { 3700d0321e0SJeremy L Thompson CeedInt ierr; 3710d0321e0SJeremy L Thompson CeedEvalMode emode; 3720d0321e0SJeremy L Thompson CeedVector vec; 3730d0321e0SJeremy L Thompson 3740d0321e0SJeremy L Thompson for (CeedInt i = 0; i < numinputfields; i++) { 3750d0321e0SJeremy L Thompson // Skip active input 3760d0321e0SJeremy L Thompson if (skipactive) { 3770d0321e0SJeremy L Thompson ierr = CeedOperatorFieldGetVector(opinputfields[i], &vec); CeedChkBackend(ierr); 3780d0321e0SJeremy L Thompson if (vec == CEED_VECTOR_ACTIVE) 3790d0321e0SJeremy L Thompson continue; 3800d0321e0SJeremy L Thompson } 3810d0321e0SJeremy L Thompson ierr = CeedQFunctionFieldGetEvalMode(qfinputfields[i], &emode); 3820d0321e0SJeremy L Thompson CeedChkBackend(ierr); 3830d0321e0SJeremy L Thompson if (emode == CEED_EVAL_WEIGHT) { // Skip 3840d0321e0SJeremy L Thompson } else { 3850d0321e0SJeremy L Thompson if (!impl->evecs[i]) { // This was a skiprestrict case 3860d0321e0SJeremy L Thompson ierr = CeedOperatorFieldGetVector(opinputfields[i], &vec); CeedChkBackend(ierr); 3870d0321e0SJeremy L Thompson ierr = CeedVectorRestoreArrayRead(vec, 3880d0321e0SJeremy L Thompson (const CeedScalar **)&edata[i]); 3890d0321e0SJeremy L Thompson CeedChkBackend(ierr); 3900d0321e0SJeremy L Thompson } else { 3910d0321e0SJeremy L Thompson ierr = CeedVectorRestoreArrayRead(impl->evecs[i], 3920d0321e0SJeremy L Thompson (const CeedScalar **) &edata[i]); 3930d0321e0SJeremy L Thompson CeedChkBackend(ierr); 3940d0321e0SJeremy L Thompson } 3950d0321e0SJeremy L Thompson } 3960d0321e0SJeremy L Thompson } 3970d0321e0SJeremy L Thompson return CEED_ERROR_SUCCESS; 3980d0321e0SJeremy L Thompson } 3990d0321e0SJeremy L Thompson 4000d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 4010d0321e0SJeremy L Thompson // Apply and add to output 4020d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 4030d0321e0SJeremy L Thompson static int CeedOperatorApplyAdd_Hip(CeedOperator op, CeedVector invec, 4040d0321e0SJeremy L Thompson CeedVector outvec, CeedRequest *request) { 4050d0321e0SJeremy L Thompson int ierr; 4060d0321e0SJeremy L Thompson CeedOperator_Hip *impl; 4070d0321e0SJeremy L Thompson ierr = CeedOperatorGetData(op, &impl); CeedChkBackend(ierr); 4080d0321e0SJeremy L Thompson CeedQFunction qf; 4090d0321e0SJeremy L Thompson ierr = CeedOperatorGetQFunction(op, &qf); CeedChkBackend(ierr); 4100d0321e0SJeremy L Thompson CeedInt Q, numelements, elemsize, numinputfields, numoutputfields, size; 4110d0321e0SJeremy L Thompson ierr = CeedOperatorGetNumQuadraturePoints(op, &Q); CeedChkBackend(ierr); 4120d0321e0SJeremy L Thompson ierr = CeedOperatorGetNumElements(op, &numelements); CeedChkBackend(ierr); 4130d0321e0SJeremy L Thompson CeedOperatorField *opinputfields, *opoutputfields; 4140d0321e0SJeremy L Thompson ierr = CeedOperatorGetFields(op, &numinputfields, &opinputfields, 4150d0321e0SJeremy L Thompson &numoutputfields, &opoutputfields); 4160d0321e0SJeremy L Thompson CeedChkBackend(ierr); 4170d0321e0SJeremy L Thompson CeedQFunctionField *qfinputfields, *qfoutputfields; 4180d0321e0SJeremy L Thompson ierr = CeedQFunctionGetFields(qf, NULL, &qfinputfields, NULL, &qfoutputfields); 4190d0321e0SJeremy L Thompson CeedChkBackend(ierr); 4200d0321e0SJeremy L Thompson CeedEvalMode emode; 4210d0321e0SJeremy L Thompson CeedVector vec; 4220d0321e0SJeremy L Thompson CeedBasis basis; 4230d0321e0SJeremy L Thompson CeedElemRestriction Erestrict; 4240d0321e0SJeremy L Thompson CeedScalar *edata[2*CEED_FIELD_MAX]; 4250d0321e0SJeremy L Thompson 4260d0321e0SJeremy L Thompson // Setup 4270d0321e0SJeremy L Thompson ierr = CeedOperatorSetup_Hip(op); CeedChkBackend(ierr); 4280d0321e0SJeremy L Thompson 4290d0321e0SJeremy L Thompson // Input Evecs and Restriction 4300d0321e0SJeremy L Thompson ierr = CeedOperatorSetupInputs_Hip(numinputfields, qfinputfields, 4310d0321e0SJeremy L Thompson opinputfields, invec, false, edata, 4320d0321e0SJeremy L Thompson impl, request); CeedChkBackend(ierr); 4330d0321e0SJeremy L Thompson 4340d0321e0SJeremy L Thompson // Input basis apply if needed 4350d0321e0SJeremy L Thompson ierr = CeedOperatorInputBasis_Hip(numelements, qfinputfields, opinputfields, 4360d0321e0SJeremy L Thompson numinputfields, false, edata, impl); 4370d0321e0SJeremy L Thompson CeedChkBackend(ierr); 4380d0321e0SJeremy L Thompson 4390d0321e0SJeremy L Thompson // Output pointers, as necessary 4400d0321e0SJeremy L Thompson for (CeedInt i = 0; i < numoutputfields; i++) { 4410d0321e0SJeremy L Thompson ierr = CeedQFunctionFieldGetEvalMode(qfoutputfields[i], &emode); 4420d0321e0SJeremy L Thompson CeedChkBackend(ierr); 4430d0321e0SJeremy L Thompson if (emode == CEED_EVAL_NONE) { 4440d0321e0SJeremy L Thompson // Set the output Q-Vector to use the E-Vector data directly. 4450d0321e0SJeremy L Thompson ierr = CeedVectorGetArrayWrite(impl->evecs[i + impl->numein], CEED_MEM_DEVICE, 4460d0321e0SJeremy L Thompson &edata[i + numinputfields]); CeedChkBackend(ierr); 4470d0321e0SJeremy L Thompson ierr = CeedVectorSetArray(impl->qvecsout[i], CEED_MEM_DEVICE, 4480d0321e0SJeremy L Thompson CEED_USE_POINTER, edata[i + numinputfields]); 4490d0321e0SJeremy L Thompson CeedChkBackend(ierr); 4500d0321e0SJeremy L Thompson } 4510d0321e0SJeremy L Thompson } 4520d0321e0SJeremy L Thompson 4530d0321e0SJeremy L Thompson // Q function 4540d0321e0SJeremy L Thompson ierr = CeedQFunctionApply(qf, numelements * Q, impl->qvecsin, impl->qvecsout); 4550d0321e0SJeremy L Thompson CeedChkBackend(ierr); 4560d0321e0SJeremy L Thompson 4570d0321e0SJeremy L Thompson // Output basis apply if needed 4580d0321e0SJeremy L Thompson for (CeedInt i = 0; i < numoutputfields; i++) { 4590d0321e0SJeremy L Thompson // Get elemsize, emode, size 4600d0321e0SJeremy L Thompson ierr = CeedOperatorFieldGetElemRestriction(opoutputfields[i], &Erestrict); 4610d0321e0SJeremy L Thompson CeedChkBackend(ierr); 4620d0321e0SJeremy L Thompson ierr = CeedElemRestrictionGetElementSize(Erestrict, &elemsize); 4630d0321e0SJeremy L Thompson CeedChkBackend(ierr); 4640d0321e0SJeremy L Thompson ierr = CeedQFunctionFieldGetEvalMode(qfoutputfields[i], &emode); 4650d0321e0SJeremy L Thompson CeedChkBackend(ierr); 4660d0321e0SJeremy L Thompson ierr = CeedQFunctionFieldGetSize(qfoutputfields[i], &size); 4670d0321e0SJeremy L Thompson CeedChkBackend(ierr); 4680d0321e0SJeremy L Thompson // Basis action 4690d0321e0SJeremy L Thompson switch (emode) { 4700d0321e0SJeremy L Thompson case CEED_EVAL_NONE: 4710d0321e0SJeremy L Thompson break; 4720d0321e0SJeremy L Thompson case CEED_EVAL_INTERP: 4730d0321e0SJeremy L Thompson ierr = CeedOperatorFieldGetBasis(opoutputfields[i], &basis); 4740d0321e0SJeremy L Thompson CeedChkBackend(ierr); 4750d0321e0SJeremy L Thompson ierr = CeedBasisApply(basis, numelements, CEED_TRANSPOSE, 4760d0321e0SJeremy L Thompson CEED_EVAL_INTERP, impl->qvecsout[i], 4770d0321e0SJeremy L Thompson impl->evecs[i + impl->numein]); CeedChkBackend(ierr); 4780d0321e0SJeremy L Thompson break; 4790d0321e0SJeremy L Thompson case CEED_EVAL_GRAD: 4800d0321e0SJeremy L Thompson ierr = CeedOperatorFieldGetBasis(opoutputfields[i], &basis); 4810d0321e0SJeremy L Thompson CeedChkBackend(ierr); 4820d0321e0SJeremy L Thompson ierr = CeedBasisApply(basis, numelements, CEED_TRANSPOSE, 4830d0321e0SJeremy L Thompson CEED_EVAL_GRAD, impl->qvecsout[i], 4840d0321e0SJeremy L Thompson impl->evecs[i + impl->numein]); CeedChkBackend(ierr); 4850d0321e0SJeremy L Thompson break; 4860d0321e0SJeremy L Thompson // LCOV_EXCL_START 4870d0321e0SJeremy L Thompson case CEED_EVAL_WEIGHT: { 4880d0321e0SJeremy L Thompson Ceed ceed; 4890d0321e0SJeremy L Thompson ierr = CeedOperatorGetCeed(op, &ceed); CeedChkBackend(ierr); 4900d0321e0SJeremy L Thompson return CeedError(ceed, CEED_ERROR_BACKEND, 4910d0321e0SJeremy L Thompson "CEED_EVAL_WEIGHT cannot be an output evaluation mode"); 4920d0321e0SJeremy L Thompson break; // Should not occur 4930d0321e0SJeremy L Thompson } 4940d0321e0SJeremy L Thompson case CEED_EVAL_DIV: 4950d0321e0SJeremy L Thompson break; // TODO: Not implemented 4960d0321e0SJeremy L Thompson case CEED_EVAL_CURL: 4970d0321e0SJeremy L Thompson break; // TODO: Not implemented 4980d0321e0SJeremy L Thompson // LCOV_EXCL_STOP 4990d0321e0SJeremy L Thompson } 5000d0321e0SJeremy L Thompson } 5010d0321e0SJeremy L Thompson 5020d0321e0SJeremy L Thompson // Output restriction 5030d0321e0SJeremy L Thompson for (CeedInt i = 0; i < numoutputfields; i++) { 5040d0321e0SJeremy L Thompson // Restore evec 5050d0321e0SJeremy L Thompson ierr = CeedQFunctionFieldGetEvalMode(qfoutputfields[i], &emode); 5060d0321e0SJeremy L Thompson CeedChkBackend(ierr); 5070d0321e0SJeremy L Thompson if (emode == CEED_EVAL_NONE) { 5080d0321e0SJeremy L Thompson ierr = CeedVectorRestoreArray(impl->evecs[i+impl->numein], 5090d0321e0SJeremy L Thompson &edata[i + numinputfields]); 5100d0321e0SJeremy L Thompson CeedChkBackend(ierr); 5110d0321e0SJeremy L Thompson } 5120d0321e0SJeremy L Thompson // Get output vector 5130d0321e0SJeremy L Thompson ierr = CeedOperatorFieldGetVector(opoutputfields[i], &vec); 5140d0321e0SJeremy L Thompson CeedChkBackend(ierr); 5150d0321e0SJeremy L Thompson // Restrict 5160d0321e0SJeremy L Thompson ierr = CeedOperatorFieldGetElemRestriction(opoutputfields[i], &Erestrict); 5170d0321e0SJeremy L Thompson CeedChkBackend(ierr); 5180d0321e0SJeremy L Thompson // Active 5190d0321e0SJeremy L Thompson if (vec == CEED_VECTOR_ACTIVE) 5200d0321e0SJeremy L Thompson vec = outvec; 5210d0321e0SJeremy L Thompson 5220d0321e0SJeremy L Thompson ierr = CeedElemRestrictionApply(Erestrict, CEED_TRANSPOSE, 5230d0321e0SJeremy L Thompson impl->evecs[i + impl->numein], vec, 5240d0321e0SJeremy L Thompson request); CeedChkBackend(ierr); 5250d0321e0SJeremy L Thompson } 5260d0321e0SJeremy L Thompson 5270d0321e0SJeremy L Thompson // Restore input arrays 5280d0321e0SJeremy L Thompson ierr = CeedOperatorRestoreInputs_Hip(numinputfields, qfinputfields, 5290d0321e0SJeremy L Thompson opinputfields, false, edata, impl); 5300d0321e0SJeremy L Thompson CeedChkBackend(ierr); 5310d0321e0SJeremy L Thompson return CEED_ERROR_SUCCESS; 5320d0321e0SJeremy L Thompson } 5330d0321e0SJeremy L Thompson 5340d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 5350d0321e0SJeremy L Thompson // Core code for assembling linear QFunction 5360d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 5370d0321e0SJeremy L Thompson static inline int CeedOperatorLinearAssembleQFunctionCore_Hip(CeedOperator op, 5380d0321e0SJeremy L Thompson bool build_objects, CeedVector *assembled, CeedElemRestriction *rstr, 5390d0321e0SJeremy L Thompson CeedRequest *request) { 5400d0321e0SJeremy L Thompson int ierr; 5410d0321e0SJeremy L Thompson CeedOperator_Hip *impl; 5420d0321e0SJeremy L Thompson ierr = CeedOperatorGetData(op, &impl); CeedChkBackend(ierr); 5430d0321e0SJeremy L Thompson CeedQFunction qf; 5440d0321e0SJeremy L Thompson ierr = CeedOperatorGetQFunction(op, &qf); CeedChkBackend(ierr); 5450d0321e0SJeremy L Thompson CeedInt Q, numelements, numinputfields, numoutputfields, size; 5460d0321e0SJeremy L Thompson ierr = CeedOperatorGetNumQuadraturePoints(op, &Q); CeedChkBackend(ierr); 5470d0321e0SJeremy L Thompson ierr = CeedOperatorGetNumElements(op, &numelements); CeedChkBackend(ierr); 5480d0321e0SJeremy L Thompson CeedOperatorField *opinputfields, *opoutputfields; 5490d0321e0SJeremy L Thompson ierr = CeedOperatorGetFields(op, &numinputfields, &opinputfields, 5500d0321e0SJeremy L Thompson &numoutputfields, &opoutputfields); 5510d0321e0SJeremy L Thompson CeedChkBackend(ierr); 5520d0321e0SJeremy L Thompson CeedQFunctionField *qfinputfields, *qfoutputfields; 5530d0321e0SJeremy L Thompson ierr = CeedQFunctionGetFields(qf, NULL, &qfinputfields, NULL, &qfoutputfields); 5540d0321e0SJeremy L Thompson CeedChkBackend(ierr); 5550d0321e0SJeremy L Thompson CeedVector vec; 5560d0321e0SJeremy L Thompson CeedInt numactivein = impl->qfnumactivein, numactiveout = impl->qfnumactiveout; 5570d0321e0SJeremy L Thompson CeedVector *activein = impl->qfactivein; 5580d0321e0SJeremy L Thompson CeedScalar *a, *tmp; 5590d0321e0SJeremy L Thompson Ceed ceed, ceedparent; 5600d0321e0SJeremy L Thompson ierr = CeedOperatorGetCeed(op, &ceed); CeedChkBackend(ierr); 5610d0321e0SJeremy L Thompson ierr = CeedGetOperatorFallbackParentCeed(ceed, &ceedparent); 5620d0321e0SJeremy L Thompson CeedChkBackend(ierr); 5630d0321e0SJeremy L Thompson ceedparent = ceedparent ? ceedparent : ceed; 5640d0321e0SJeremy L Thompson CeedScalar *edata[2*CEED_FIELD_MAX]; 5650d0321e0SJeremy L Thompson 5660d0321e0SJeremy L Thompson // Setup 5670d0321e0SJeremy L Thompson ierr = CeedOperatorSetup_Hip(op); CeedChkBackend(ierr); 5680d0321e0SJeremy L Thompson 5690d0321e0SJeremy L Thompson // Check for identity 5700d0321e0SJeremy L Thompson bool identityqf; 5710d0321e0SJeremy L Thompson ierr = CeedQFunctionIsIdentity(qf, &identityqf); CeedChkBackend(ierr); 5720d0321e0SJeremy L Thompson if (identityqf) 5730d0321e0SJeremy L Thompson // LCOV_EXCL_START 5740d0321e0SJeremy L Thompson return CeedError(ceed, CEED_ERROR_BACKEND, 5750d0321e0SJeremy L Thompson "Assembling identity QFunctions not supported"); 5760d0321e0SJeremy L Thompson // LCOV_EXCL_STOP 5770d0321e0SJeremy L Thompson 5780d0321e0SJeremy L Thompson // Input Evecs and Restriction 5790d0321e0SJeremy L Thompson ierr = CeedOperatorSetupInputs_Hip(numinputfields, qfinputfields, 5800d0321e0SJeremy L Thompson opinputfields, NULL, true, edata, 5810d0321e0SJeremy L Thompson impl, request); CeedChkBackend(ierr); 5820d0321e0SJeremy L Thompson 5830d0321e0SJeremy L Thompson // Count number of active input fields 5840d0321e0SJeremy L Thompson if (!numactivein) { 5850d0321e0SJeremy L Thompson for (CeedInt i=0; i<numinputfields; i++) { 5860d0321e0SJeremy L Thompson // Get input vector 5870d0321e0SJeremy L Thompson ierr = CeedOperatorFieldGetVector(opinputfields[i], &vec); CeedChkBackend(ierr); 5880d0321e0SJeremy L Thompson // Check if active input 5890d0321e0SJeremy L Thompson if (vec == CEED_VECTOR_ACTIVE) { 5900d0321e0SJeremy L Thompson ierr = CeedQFunctionFieldGetSize(qfinputfields[i], &size); CeedChkBackend(ierr); 5910d0321e0SJeremy L Thompson ierr = CeedVectorSetValue(impl->qvecsin[i], 0.0); CeedChkBackend(ierr); 5920d0321e0SJeremy L Thompson ierr = CeedVectorGetArray(impl->qvecsin[i], CEED_MEM_DEVICE, &tmp); 5930d0321e0SJeremy L Thompson CeedChkBackend(ierr); 5940d0321e0SJeremy L Thompson ierr = CeedRealloc(numactivein + size, &activein); CeedChkBackend(ierr); 5950d0321e0SJeremy L Thompson for (CeedInt field = 0; field < size; field++) { 5960d0321e0SJeremy L Thompson ierr = CeedVectorCreate(ceed, Q*numelements, 5970d0321e0SJeremy L Thompson &activein[numactivein+field]); CeedChkBackend(ierr); 5980d0321e0SJeremy L Thompson ierr = CeedVectorSetArray(activein[numactivein+field], CEED_MEM_DEVICE, 5990d0321e0SJeremy L Thompson CEED_USE_POINTER, &tmp[field*Q*numelements]); 6000d0321e0SJeremy L Thompson CeedChkBackend(ierr); 6010d0321e0SJeremy L Thompson } 6020d0321e0SJeremy L Thompson numactivein += size; 6030d0321e0SJeremy L Thompson ierr = CeedVectorRestoreArray(impl->qvecsin[i], &tmp); CeedChkBackend(ierr); 6040d0321e0SJeremy L Thompson } 6050d0321e0SJeremy L Thompson } 6060d0321e0SJeremy L Thompson impl->qfnumactivein = numactivein; 6070d0321e0SJeremy L Thompson impl->qfactivein = activein; 6080d0321e0SJeremy L Thompson } 6090d0321e0SJeremy L Thompson 6100d0321e0SJeremy L Thompson // Count number of active output fields 6110d0321e0SJeremy L Thompson if (!numactiveout) { 6120d0321e0SJeremy L Thompson for (CeedInt i=0; i<numoutputfields; i++) { 6130d0321e0SJeremy L Thompson // Get output vector 6140d0321e0SJeremy L Thompson ierr = CeedOperatorFieldGetVector(opoutputfields[i], &vec); 6150d0321e0SJeremy L Thompson CeedChkBackend(ierr); 6160d0321e0SJeremy L Thompson // Check if active output 6170d0321e0SJeremy L Thompson if (vec == CEED_VECTOR_ACTIVE) { 6180d0321e0SJeremy L Thompson ierr = CeedQFunctionFieldGetSize(qfoutputfields[i], &size); 6190d0321e0SJeremy L Thompson CeedChkBackend(ierr); 6200d0321e0SJeremy L Thompson numactiveout += size; 6210d0321e0SJeremy L Thompson } 6220d0321e0SJeremy L Thompson } 6230d0321e0SJeremy L Thompson impl->qfnumactiveout = numactiveout; 6240d0321e0SJeremy L Thompson } 6250d0321e0SJeremy L Thompson 6260d0321e0SJeremy L Thompson // Check sizes 6270d0321e0SJeremy L Thompson if (!numactivein || !numactiveout) 6280d0321e0SJeremy L Thompson // LCOV_EXCL_START 6290d0321e0SJeremy L Thompson return CeedError(ceed, CEED_ERROR_BACKEND, 6300d0321e0SJeremy L Thompson "Cannot assemble QFunction without active inputs " 6310d0321e0SJeremy L Thompson "and outputs"); 6320d0321e0SJeremy L Thompson // LCOV_EXCL_STOP 6330d0321e0SJeremy L Thompson 6340d0321e0SJeremy L Thompson // Build objects if needed 6350d0321e0SJeremy L Thompson if (build_objects) { 6360d0321e0SJeremy L Thompson // Create output restriction 6370d0321e0SJeremy L Thompson CeedInt strides[3] = {1, numelements*Q, Q}; /* *NOPAD* */ 6380d0321e0SJeremy L Thompson ierr = CeedElemRestrictionCreateStrided(ceedparent, numelements, Q, 6390d0321e0SJeremy L Thompson numactivein*numactiveout, 6400d0321e0SJeremy L Thompson numactivein*numactiveout*numelements*Q, 6410d0321e0SJeremy L Thompson strides, rstr); CeedChkBackend(ierr); 6420d0321e0SJeremy L Thompson // Create assembled vector 6430d0321e0SJeremy L Thompson ierr = CeedVectorCreate(ceedparent, numelements*Q*numactivein*numactiveout, 6440d0321e0SJeremy L Thompson assembled); CeedChkBackend(ierr); 6450d0321e0SJeremy L Thompson } 6460d0321e0SJeremy L Thompson ierr = CeedVectorSetValue(*assembled, 0.0); CeedChkBackend(ierr); 6470d0321e0SJeremy L Thompson ierr = CeedVectorGetArray(*assembled, CEED_MEM_DEVICE, &a); 6480d0321e0SJeremy L Thompson CeedChkBackend(ierr); 6490d0321e0SJeremy L Thompson 6500d0321e0SJeremy L Thompson // Input basis apply 6510d0321e0SJeremy L Thompson ierr = CeedOperatorInputBasis_Hip(numelements, qfinputfields, opinputfields, 6520d0321e0SJeremy L Thompson numinputfields, true, edata, impl); 6530d0321e0SJeremy L Thompson CeedChkBackend(ierr); 6540d0321e0SJeremy L Thompson 6550d0321e0SJeremy L Thompson // Assemble QFunction 6560d0321e0SJeremy L Thompson for (CeedInt in=0; in<numactivein; in++) { 6570d0321e0SJeremy L Thompson // Set Inputs 6580d0321e0SJeremy L Thompson ierr = CeedVectorSetValue(activein[in], 1.0); CeedChkBackend(ierr); 6590d0321e0SJeremy L Thompson if (numactivein > 1) { 6600d0321e0SJeremy L Thompson ierr = CeedVectorSetValue(activein[(in+numactivein-1)%numactivein], 6610d0321e0SJeremy L Thompson 0.0); CeedChkBackend(ierr); 6620d0321e0SJeremy L Thompson } 6630d0321e0SJeremy L Thompson // Set Outputs 6640d0321e0SJeremy L Thompson for (CeedInt out=0; out<numoutputfields; out++) { 6650d0321e0SJeremy L Thompson // Get output vector 6660d0321e0SJeremy L Thompson ierr = CeedOperatorFieldGetVector(opoutputfields[out], &vec); 6670d0321e0SJeremy L Thompson CeedChkBackend(ierr); 6680d0321e0SJeremy L Thompson // Check if active output 6690d0321e0SJeremy L Thompson if (vec == CEED_VECTOR_ACTIVE) { 6700d0321e0SJeremy L Thompson CeedVectorSetArray(impl->qvecsout[out], CEED_MEM_DEVICE, 6710d0321e0SJeremy L Thompson CEED_USE_POINTER, a); CeedChkBackend(ierr); 6720d0321e0SJeremy L Thompson ierr = CeedQFunctionFieldGetSize(qfoutputfields[out], &size); 6730d0321e0SJeremy L Thompson CeedChkBackend(ierr); 6740d0321e0SJeremy L Thompson a += size*Q*numelements; // Advance the pointer by the size of the output 6750d0321e0SJeremy L Thompson } 6760d0321e0SJeremy L Thompson } 6770d0321e0SJeremy L Thompson // Apply QFunction 6780d0321e0SJeremy L Thompson ierr = CeedQFunctionApply(qf, Q*numelements, impl->qvecsin, impl->qvecsout); 6790d0321e0SJeremy L Thompson CeedChkBackend(ierr); 6800d0321e0SJeremy L Thompson } 6810d0321e0SJeremy L Thompson 6820d0321e0SJeremy L Thompson // Un-set output Qvecs to prevent accidental overwrite of Assembled 6830d0321e0SJeremy L Thompson for (CeedInt out=0; out<numoutputfields; out++) { 6840d0321e0SJeremy L Thompson // Get output vector 6850d0321e0SJeremy L Thompson ierr = CeedOperatorFieldGetVector(opoutputfields[out], &vec); 6860d0321e0SJeremy L Thompson CeedChkBackend(ierr); 6870d0321e0SJeremy L Thompson // Check if active output 6880d0321e0SJeremy L Thompson if (vec == CEED_VECTOR_ACTIVE) { 6890d0321e0SJeremy L Thompson ierr = CeedVectorTakeArray(impl->qvecsout[out], CEED_MEM_DEVICE, NULL); 6900d0321e0SJeremy L Thompson CeedChkBackend(ierr); 6910d0321e0SJeremy L Thompson } 6920d0321e0SJeremy L Thompson } 6930d0321e0SJeremy L Thompson 6940d0321e0SJeremy L Thompson // Restore input arrays 6950d0321e0SJeremy L Thompson ierr = CeedOperatorRestoreInputs_Hip(numinputfields, qfinputfields, 6960d0321e0SJeremy L Thompson opinputfields, true, edata, impl); 6970d0321e0SJeremy L Thompson CeedChkBackend(ierr); 6980d0321e0SJeremy L Thompson 6990d0321e0SJeremy L Thompson // Restore output 7000d0321e0SJeremy L Thompson ierr = CeedVectorRestoreArray(*assembled, &a); CeedChkBackend(ierr); 7010d0321e0SJeremy L Thompson 7020d0321e0SJeremy L Thompson return CEED_ERROR_SUCCESS; 7030d0321e0SJeremy L Thompson } 7040d0321e0SJeremy L Thompson 7050d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 7060d0321e0SJeremy L Thompson // Assemble Linear QFunction 7070d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 7080d0321e0SJeremy L Thompson static int CeedOperatorLinearAssembleQFunction_Hip(CeedOperator op, 7090d0321e0SJeremy L Thompson CeedVector *assembled, CeedElemRestriction *rstr, CeedRequest *request) { 7100d0321e0SJeremy L Thompson return CeedOperatorLinearAssembleQFunctionCore_Hip(op, true, assembled, rstr, 7110d0321e0SJeremy L Thompson request); 7120d0321e0SJeremy L Thompson } 7130d0321e0SJeremy L Thompson 7140d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 7150d0321e0SJeremy L Thompson // Assemble Linear QFunction 7160d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 7170d0321e0SJeremy L Thompson static int CeedOperatorLinearAssembleQFunctionUpdate_Hip(CeedOperator op, 7180d0321e0SJeremy L Thompson CeedVector assembled, CeedElemRestriction rstr, CeedRequest *request) { 7190d0321e0SJeremy L Thompson return CeedOperatorLinearAssembleQFunctionCore_Hip(op, false, &assembled, &rstr, 7200d0321e0SJeremy L Thompson request); 7210d0321e0SJeremy L Thompson } 7220d0321e0SJeremy L Thompson 7230d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 7240d0321e0SJeremy L Thompson // Diagonal assembly kernels 7250d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 7260d0321e0SJeremy L Thompson // *INDENT-OFF* 7270d0321e0SJeremy L Thompson static const char *diagonalkernels = QUOTE( 7280d0321e0SJeremy L Thompson 7290d0321e0SJeremy L Thompson typedef enum { 7300d0321e0SJeremy L Thompson /// Perform no evaluation (either because there is no data or it is already at 7310d0321e0SJeremy L Thompson /// quadrature points) 7320d0321e0SJeremy L Thompson CEED_EVAL_NONE = 0, 7330d0321e0SJeremy L Thompson /// Interpolate from nodes to quadrature points 7340d0321e0SJeremy L Thompson CEED_EVAL_INTERP = 1, 7350d0321e0SJeremy L Thompson /// Evaluate gradients at quadrature points from input in a nodal basis 7360d0321e0SJeremy L Thompson CEED_EVAL_GRAD = 2, 7370d0321e0SJeremy L Thompson /// Evaluate divergence at quadrature points from input in a nodal basis 7380d0321e0SJeremy L Thompson CEED_EVAL_DIV = 4, 7390d0321e0SJeremy L Thompson /// Evaluate curl at quadrature points from input in a nodal basis 7400d0321e0SJeremy L Thompson CEED_EVAL_CURL = 8, 7410d0321e0SJeremy L Thompson /// Using no input, evaluate quadrature weights on the reference element 7420d0321e0SJeremy L Thompson CEED_EVAL_WEIGHT = 16, 7430d0321e0SJeremy L Thompson } CeedEvalMode; 7440d0321e0SJeremy L Thompson 7450d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 7460d0321e0SJeremy L Thompson // Get Basis Emode Pointer 7470d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 7480d0321e0SJeremy L Thompson extern "C" __device__ void CeedOperatorGetBasisPointer_Hip(const CeedScalar **basisptr, 7490d0321e0SJeremy L Thompson CeedEvalMode emode, const CeedScalar *identity, const CeedScalar *interp, 7500d0321e0SJeremy L Thompson const CeedScalar *grad) { 7510d0321e0SJeremy L Thompson switch (emode) { 7520d0321e0SJeremy L Thompson case CEED_EVAL_NONE: 7530d0321e0SJeremy L Thompson *basisptr = identity; 7540d0321e0SJeremy L Thompson break; 7550d0321e0SJeremy L Thompson case CEED_EVAL_INTERP: 7560d0321e0SJeremy L Thompson *basisptr = interp; 7570d0321e0SJeremy L Thompson break; 7580d0321e0SJeremy L Thompson case CEED_EVAL_GRAD: 7590d0321e0SJeremy L Thompson *basisptr = grad; 7600d0321e0SJeremy L Thompson break; 7610d0321e0SJeremy L Thompson case CEED_EVAL_WEIGHT: 7620d0321e0SJeremy L Thompson case CEED_EVAL_DIV: 7630d0321e0SJeremy L Thompson case CEED_EVAL_CURL: 7640d0321e0SJeremy L Thompson break; // Caught by QF Assembly 7650d0321e0SJeremy L Thompson } 7660d0321e0SJeremy L Thompson } 7670d0321e0SJeremy L Thompson 7680d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 7690d0321e0SJeremy L Thompson // Core code for diagonal assembly 7700d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 7710d0321e0SJeremy L Thompson __device__ void diagonalCore(const CeedInt nelem, 7720d0321e0SJeremy L Thompson const CeedScalar maxnorm, const bool pointBlock, 7730d0321e0SJeremy L Thompson const CeedScalar *identity, 7740d0321e0SJeremy L Thompson const CeedScalar *interpin, const CeedScalar *gradin, 7750d0321e0SJeremy L Thompson const CeedScalar *interpout, const CeedScalar *gradout, 7760d0321e0SJeremy L Thompson const CeedEvalMode *emodein, const CeedEvalMode *emodeout, 7770d0321e0SJeremy L Thompson const CeedScalar *__restrict__ assembledqfarray, 7780d0321e0SJeremy L Thompson CeedScalar *__restrict__ elemdiagarray) { 7790d0321e0SJeremy L Thompson const int tid = threadIdx.x; // running with P threads, tid is evec node 7800d0321e0SJeremy L Thompson const CeedScalar qfvaluebound = maxnorm*1e-12; 7810d0321e0SJeremy L Thompson 7820d0321e0SJeremy L Thompson // Compute the diagonal of B^T D B 7830d0321e0SJeremy L Thompson // Each element 7840d0321e0SJeremy L Thompson for (CeedInt e = blockIdx.x*blockDim.z + threadIdx.z; e < nelem; 7850d0321e0SJeremy L Thompson e += gridDim.x*blockDim.z) { 7860d0321e0SJeremy L Thompson CeedInt dout = -1; 7870d0321e0SJeremy L Thompson // Each basis eval mode pair 7880d0321e0SJeremy L Thompson for (CeedInt eout = 0; eout < NUMEMODEOUT; eout++) { 7890d0321e0SJeremy L Thompson const CeedScalar *bt = NULL; 7900d0321e0SJeremy L Thompson if (emodeout[eout] == CEED_EVAL_GRAD) 7910d0321e0SJeremy L Thompson dout += 1; 7920d0321e0SJeremy L Thompson CeedOperatorGetBasisPointer_Hip(&bt, emodeout[eout], identity, interpout, 7930d0321e0SJeremy L Thompson &gradout[dout*NQPTS*NNODES]); 7940d0321e0SJeremy L Thompson CeedInt din = -1; 7950d0321e0SJeremy L Thompson for (CeedInt ein = 0; ein < NUMEMODEIN; ein++) { 7960d0321e0SJeremy L Thompson const CeedScalar *b = NULL; 7970d0321e0SJeremy L Thompson if (emodein[ein] == CEED_EVAL_GRAD) 7980d0321e0SJeremy L Thompson din += 1; 7990d0321e0SJeremy L Thompson CeedOperatorGetBasisPointer_Hip(&b, emodein[ein], identity, interpin, 8000d0321e0SJeremy L Thompson &gradin[din*NQPTS*NNODES]); 8010d0321e0SJeremy L Thompson // Each component 8020d0321e0SJeremy L Thompson for (CeedInt compOut = 0; compOut < NCOMP; compOut++) { 8030d0321e0SJeremy L Thompson // Each qpoint/node pair 8040d0321e0SJeremy L Thompson if (pointBlock) { 8050d0321e0SJeremy L Thompson // Point Block Diagonal 8060d0321e0SJeremy L Thompson for (CeedInt compIn = 0; compIn < NCOMP; compIn++) { 8070d0321e0SJeremy L Thompson CeedScalar evalue = 0.; 8080d0321e0SJeremy L Thompson for (CeedInt q = 0; q < NQPTS; q++) { 8090d0321e0SJeremy L Thompson const CeedScalar qfvalue = 8100d0321e0SJeremy L Thompson assembledqfarray[((((ein*NCOMP+compIn)*NUMEMODEOUT+eout)* 8110d0321e0SJeremy L Thompson NCOMP+compOut)*nelem+e)*NQPTS+q]; 8120d0321e0SJeremy L Thompson if (abs(qfvalue) > qfvaluebound) 8130d0321e0SJeremy L Thompson evalue += bt[q*NNODES+tid] * qfvalue * b[q*NNODES+tid]; 8140d0321e0SJeremy L Thompson } 8150d0321e0SJeremy L Thompson elemdiagarray[((compOut*NCOMP+compIn)*nelem+e)*NNODES+tid] += evalue; 8160d0321e0SJeremy L Thompson } 8170d0321e0SJeremy L Thompson } else { 8180d0321e0SJeremy L Thompson // Diagonal Only 8190d0321e0SJeremy L Thompson CeedScalar evalue = 0.; 8200d0321e0SJeremy L Thompson for (CeedInt q = 0; q < NQPTS; q++) { 8210d0321e0SJeremy L Thompson const CeedScalar qfvalue = 8220d0321e0SJeremy L Thompson assembledqfarray[((((ein*NCOMP+compOut)*NUMEMODEOUT+eout)* 8230d0321e0SJeremy L Thompson NCOMP+compOut)*nelem+e)*NQPTS+q]; 8240d0321e0SJeremy L Thompson if (abs(qfvalue) > qfvaluebound) 8250d0321e0SJeremy L Thompson evalue += bt[q*NNODES+tid] * qfvalue * b[q*NNODES+tid]; 8260d0321e0SJeremy L Thompson } 8270d0321e0SJeremy L Thompson elemdiagarray[(compOut*nelem+e)*NNODES+tid] += evalue; 8280d0321e0SJeremy L Thompson } 8290d0321e0SJeremy L Thompson } 8300d0321e0SJeremy L Thompson } 8310d0321e0SJeremy L Thompson } 8320d0321e0SJeremy L Thompson } 8330d0321e0SJeremy L Thompson } 8340d0321e0SJeremy L Thompson 8350d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 8360d0321e0SJeremy L Thompson // Linear diagonal 8370d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 8380d0321e0SJeremy L Thompson extern "C" __global__ void linearDiagonal(const CeedInt nelem, 8390d0321e0SJeremy L Thompson const CeedScalar maxnorm, const CeedScalar *identity, 8400d0321e0SJeremy L Thompson const CeedScalar *interpin, const CeedScalar *gradin, 8410d0321e0SJeremy L Thompson const CeedScalar *interpout, const CeedScalar *gradout, 8420d0321e0SJeremy L Thompson const CeedEvalMode *emodein, const CeedEvalMode *emodeout, 8430d0321e0SJeremy L Thompson const CeedScalar *__restrict__ assembledqfarray, 8440d0321e0SJeremy L Thompson CeedScalar *__restrict__ elemdiagarray) { 8450d0321e0SJeremy L Thompson diagonalCore(nelem, maxnorm, false, identity, interpin, gradin, interpout, 8460d0321e0SJeremy L Thompson gradout, emodein, emodeout, assembledqfarray, elemdiagarray); 8470d0321e0SJeremy L Thompson } 8480d0321e0SJeremy L Thompson 8490d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 8500d0321e0SJeremy L Thompson // Linear point block diagonal 8510d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 8520d0321e0SJeremy L Thompson extern "C" __global__ void linearPointBlockDiagonal(const CeedInt nelem, 8530d0321e0SJeremy L Thompson const CeedScalar maxnorm, const CeedScalar *identity, 8540d0321e0SJeremy L Thompson const CeedScalar *interpin, const CeedScalar *gradin, 8550d0321e0SJeremy L Thompson const CeedScalar *interpout, const CeedScalar *gradout, 8560d0321e0SJeremy L Thompson const CeedEvalMode *emodein, const CeedEvalMode *emodeout, 8570d0321e0SJeremy L Thompson const CeedScalar *__restrict__ assembledqfarray, 8580d0321e0SJeremy L Thompson CeedScalar *__restrict__ elemdiagarray) { 8590d0321e0SJeremy L Thompson diagonalCore(nelem, maxnorm, true, identity, interpin, gradin, interpout, 8600d0321e0SJeremy L Thompson gradout, emodein, emodeout, assembledqfarray, elemdiagarray); 8610d0321e0SJeremy L Thompson } 8620d0321e0SJeremy L Thompson 8630d0321e0SJeremy L Thompson ); 8640d0321e0SJeremy L Thompson // *INDENT-ON* 8650d0321e0SJeremy L Thompson 8660d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 8670d0321e0SJeremy L Thompson // Create point block restriction 8680d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 8690d0321e0SJeremy L Thompson static int CreatePBRestriction(CeedElemRestriction rstr, 8700d0321e0SJeremy L Thompson CeedElemRestriction *pbRstr) { 8710d0321e0SJeremy L Thompson int ierr; 8720d0321e0SJeremy L Thompson Ceed ceed; 8730d0321e0SJeremy L Thompson ierr = CeedElemRestrictionGetCeed(rstr, &ceed); CeedChkBackend(ierr); 8740d0321e0SJeremy L Thompson const CeedInt *offsets; 8750d0321e0SJeremy L Thompson ierr = CeedElemRestrictionGetOffsets(rstr, CEED_MEM_HOST, &offsets); 8760d0321e0SJeremy L Thompson CeedChkBackend(ierr); 8770d0321e0SJeremy L Thompson 8780d0321e0SJeremy L Thompson // Expand offsets 8790d0321e0SJeremy L Thompson CeedInt nelem, ncomp, elemsize, compstride, max = 1, *pbOffsets; 8800d0321e0SJeremy L Thompson ierr = CeedElemRestrictionGetNumElements(rstr, &nelem); CeedChkBackend(ierr); 8810d0321e0SJeremy L Thompson ierr = CeedElemRestrictionGetNumComponents(rstr, &ncomp); CeedChkBackend(ierr); 8820d0321e0SJeremy L Thompson ierr = CeedElemRestrictionGetElementSize(rstr, &elemsize); CeedChkBackend(ierr); 8830d0321e0SJeremy L Thompson ierr = CeedElemRestrictionGetCompStride(rstr, &compstride); 8840d0321e0SJeremy L Thompson CeedChkBackend(ierr); 8850d0321e0SJeremy L Thompson CeedInt shift = ncomp; 8860d0321e0SJeremy L Thompson if (compstride != 1) 8870d0321e0SJeremy L Thompson shift *= ncomp; 8880d0321e0SJeremy L Thompson ierr = CeedCalloc(nelem*elemsize, &pbOffsets); CeedChkBackend(ierr); 8890d0321e0SJeremy L Thompson for (CeedInt i = 0; i < nelem*elemsize; i++) { 8900d0321e0SJeremy L Thompson pbOffsets[i] = offsets[i]*shift; 8910d0321e0SJeremy L Thompson if (pbOffsets[i] > max) 8920d0321e0SJeremy L Thompson max = pbOffsets[i]; 8930d0321e0SJeremy L Thompson } 8940d0321e0SJeremy L Thompson 8950d0321e0SJeremy L Thompson // Create new restriction 8960d0321e0SJeremy L Thompson ierr = CeedElemRestrictionCreate(ceed, nelem, elemsize, ncomp*ncomp, 1, 8970d0321e0SJeremy L Thompson max + ncomp*ncomp, CEED_MEM_HOST, 8980d0321e0SJeremy L Thompson CEED_OWN_POINTER, pbOffsets, pbRstr); 8990d0321e0SJeremy L Thompson CeedChkBackend(ierr); 9000d0321e0SJeremy L Thompson 9010d0321e0SJeremy L Thompson // Cleanup 9020d0321e0SJeremy L Thompson ierr = CeedElemRestrictionRestoreOffsets(rstr, &offsets); CeedChkBackend(ierr); 9030d0321e0SJeremy L Thompson 9040d0321e0SJeremy L Thompson return CEED_ERROR_SUCCESS; 9050d0321e0SJeremy L Thompson } 9060d0321e0SJeremy L Thompson 9070d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 9080d0321e0SJeremy L Thompson // Assemble diagonal setup 9090d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 9100d0321e0SJeremy L Thompson static inline int CeedOperatorAssembleDiagonalSetup_Hip(CeedOperator op, 9110d0321e0SJeremy L Thompson const bool pointBlock) { 9120d0321e0SJeremy L Thompson int ierr; 9130d0321e0SJeremy L Thompson Ceed ceed; 9140d0321e0SJeremy L Thompson ierr = CeedOperatorGetCeed(op, &ceed); CeedChkBackend(ierr); 9150d0321e0SJeremy L Thompson CeedQFunction qf; 9160d0321e0SJeremy L Thompson ierr = CeedOperatorGetQFunction(op, &qf); CeedChkBackend(ierr); 9170d0321e0SJeremy L Thompson CeedInt numinputfields, numoutputfields; 9180d0321e0SJeremy L Thompson ierr = CeedQFunctionGetNumArgs(qf, &numinputfields, &numoutputfields); 9190d0321e0SJeremy L Thompson CeedChkBackend(ierr); 9200d0321e0SJeremy L Thompson 9210d0321e0SJeremy L Thompson // Determine active input basis 9220d0321e0SJeremy L Thompson CeedOperatorField *opfields; 9230d0321e0SJeremy L Thompson CeedQFunctionField *qffields; 9240d0321e0SJeremy L Thompson ierr = CeedOperatorGetFields(op, NULL, &opfields, NULL, NULL); 9250d0321e0SJeremy L Thompson CeedChkBackend(ierr); 9260d0321e0SJeremy L Thompson ierr = CeedQFunctionGetFields(qf, NULL, &qffields, NULL, NULL); 9270d0321e0SJeremy L Thompson CeedChkBackend(ierr); 9280d0321e0SJeremy L Thompson CeedInt numemodein = 0, ncomp = 0, dim = 1; 9290d0321e0SJeremy L Thompson CeedEvalMode *emodein = NULL; 9300d0321e0SJeremy L Thompson CeedBasis basisin = NULL; 9310d0321e0SJeremy L Thompson CeedElemRestriction rstrin = NULL; 9320d0321e0SJeremy L Thompson for (CeedInt i = 0; i < numinputfields; i++) { 9330d0321e0SJeremy L Thompson CeedVector vec; 9340d0321e0SJeremy L Thompson ierr = CeedOperatorFieldGetVector(opfields[i], &vec); CeedChkBackend(ierr); 9350d0321e0SJeremy L Thompson if (vec == CEED_VECTOR_ACTIVE) { 9360d0321e0SJeremy L Thompson CeedElemRestriction rstr; 9370d0321e0SJeremy L Thompson ierr = CeedOperatorFieldGetBasis(opfields[i], &basisin); CeedChkBackend(ierr); 9380d0321e0SJeremy L Thompson ierr = CeedBasisGetNumComponents(basisin, &ncomp); CeedChkBackend(ierr); 9390d0321e0SJeremy L Thompson ierr = CeedBasisGetDimension(basisin, &dim); CeedChkBackend(ierr); 9400d0321e0SJeremy L Thompson ierr = CeedOperatorFieldGetElemRestriction(opfields[i], &rstr); 9410d0321e0SJeremy L Thompson CeedChkBackend(ierr); 9420d0321e0SJeremy L Thompson if (rstrin && rstrin != rstr) 9430d0321e0SJeremy L Thompson // LCOV_EXCL_START 9440d0321e0SJeremy L Thompson return CeedError(ceed, CEED_ERROR_BACKEND, 9450d0321e0SJeremy L Thompson "Multi-field non-composite operator diagonal assembly not supported"); 9460d0321e0SJeremy L Thompson // LCOV_EXCL_STOP 9470d0321e0SJeremy L Thompson rstrin = rstr; 9480d0321e0SJeremy L Thompson CeedEvalMode emode; 9490d0321e0SJeremy L Thompson ierr = CeedQFunctionFieldGetEvalMode(qffields[i], &emode); 9500d0321e0SJeremy L Thompson CeedChkBackend(ierr); 9510d0321e0SJeremy L Thompson switch (emode) { 9520d0321e0SJeremy L Thompson case CEED_EVAL_NONE: 9530d0321e0SJeremy L Thompson case CEED_EVAL_INTERP: 9540d0321e0SJeremy L Thompson ierr = CeedRealloc(numemodein + 1, &emodein); CeedChkBackend(ierr); 9550d0321e0SJeremy L Thompson emodein[numemodein] = emode; 9560d0321e0SJeremy L Thompson numemodein += 1; 9570d0321e0SJeremy L Thompson break; 9580d0321e0SJeremy L Thompson case CEED_EVAL_GRAD: 9590d0321e0SJeremy L Thompson ierr = CeedRealloc(numemodein + dim, &emodein); CeedChkBackend(ierr); 9600d0321e0SJeremy L Thompson for (CeedInt d = 0; d < dim; d++) 9610d0321e0SJeremy L Thompson emodein[numemodein+d] = emode; 9620d0321e0SJeremy L Thompson numemodein += dim; 9630d0321e0SJeremy L Thompson break; 9640d0321e0SJeremy L Thompson case CEED_EVAL_WEIGHT: 9650d0321e0SJeremy L Thompson case CEED_EVAL_DIV: 9660d0321e0SJeremy L Thompson case CEED_EVAL_CURL: 9670d0321e0SJeremy L Thompson break; // Caught by QF Assembly 9680d0321e0SJeremy L Thompson } 9690d0321e0SJeremy L Thompson } 9700d0321e0SJeremy L Thompson } 9710d0321e0SJeremy L Thompson 9720d0321e0SJeremy L Thompson // Determine active output basis 9730d0321e0SJeremy L Thompson ierr = CeedOperatorGetFields(op, NULL, NULL, NULL, &opfields); 9740d0321e0SJeremy L Thompson CeedChkBackend(ierr); 9750d0321e0SJeremy L Thompson ierr = CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qffields); 9760d0321e0SJeremy L Thompson CeedChkBackend(ierr); 9770d0321e0SJeremy L Thompson CeedInt numemodeout = 0; 9780d0321e0SJeremy L Thompson CeedEvalMode *emodeout = NULL; 9790d0321e0SJeremy L Thompson CeedBasis basisout = NULL; 9800d0321e0SJeremy L Thompson CeedElemRestriction rstrout = NULL; 9810d0321e0SJeremy L Thompson for (CeedInt i = 0; i < numoutputfields; i++) { 9820d0321e0SJeremy L Thompson CeedVector vec; 9830d0321e0SJeremy L Thompson ierr = CeedOperatorFieldGetVector(opfields[i], &vec); CeedChkBackend(ierr); 9840d0321e0SJeremy L Thompson if (vec == CEED_VECTOR_ACTIVE) { 9850d0321e0SJeremy L Thompson CeedElemRestriction rstr; 9860d0321e0SJeremy L Thompson ierr = CeedOperatorFieldGetBasis(opfields[i], &basisout); CeedChkBackend(ierr); 9870d0321e0SJeremy L Thompson ierr = CeedOperatorFieldGetElemRestriction(opfields[i], &rstr); 9880d0321e0SJeremy L Thompson CeedChkBackend(ierr); 9890d0321e0SJeremy L Thompson if (rstrout && rstrout != rstr) 9900d0321e0SJeremy L Thompson // LCOV_EXCL_START 9910d0321e0SJeremy L Thompson return CeedError(ceed, CEED_ERROR_BACKEND, 9920d0321e0SJeremy L Thompson "Multi-field non-composite operator diagonal assembly not supported"); 9930d0321e0SJeremy L Thompson // LCOV_EXCL_STOP 9940d0321e0SJeremy L Thompson rstrout = rstr; 9950d0321e0SJeremy L Thompson CeedEvalMode emode; 9960d0321e0SJeremy L Thompson ierr = CeedQFunctionFieldGetEvalMode(qffields[i], &emode); CeedChkBackend(ierr); 9970d0321e0SJeremy L Thompson switch (emode) { 9980d0321e0SJeremy L Thompson case CEED_EVAL_NONE: 9990d0321e0SJeremy L Thompson case CEED_EVAL_INTERP: 10000d0321e0SJeremy L Thompson ierr = CeedRealloc(numemodeout + 1, &emodeout); CeedChkBackend(ierr); 10010d0321e0SJeremy L Thompson emodeout[numemodeout] = emode; 10020d0321e0SJeremy L Thompson numemodeout += 1; 10030d0321e0SJeremy L Thompson break; 10040d0321e0SJeremy L Thompson case CEED_EVAL_GRAD: 10050d0321e0SJeremy L Thompson ierr = CeedRealloc(numemodeout + dim, &emodeout); CeedChkBackend(ierr); 10060d0321e0SJeremy L Thompson for (CeedInt d = 0; d < dim; d++) 10070d0321e0SJeremy L Thompson emodeout[numemodeout+d] = emode; 10080d0321e0SJeremy L Thompson numemodeout += dim; 10090d0321e0SJeremy L Thompson break; 10100d0321e0SJeremy L Thompson case CEED_EVAL_WEIGHT: 10110d0321e0SJeremy L Thompson case CEED_EVAL_DIV: 10120d0321e0SJeremy L Thompson case CEED_EVAL_CURL: 10130d0321e0SJeremy L Thompson break; // Caught by QF Assembly 10140d0321e0SJeremy L Thompson } 10150d0321e0SJeremy L Thompson } 10160d0321e0SJeremy L Thompson } 10170d0321e0SJeremy L Thompson 10180d0321e0SJeremy L Thompson // Operator data struct 10190d0321e0SJeremy L Thompson CeedOperator_Hip *impl; 10200d0321e0SJeremy L Thompson ierr = CeedOperatorGetData(op, &impl); CeedChkBackend(ierr); 10210d0321e0SJeremy L Thompson ierr = CeedCalloc(1, &impl->diag); CeedChkBackend(ierr); 10220d0321e0SJeremy L Thompson CeedOperatorDiag_Hip *diag = impl->diag; 10230d0321e0SJeremy L Thompson diag->basisin = basisin; 10240d0321e0SJeremy L Thompson diag->basisout = basisout; 10250d0321e0SJeremy L Thompson diag->h_emodein = emodein; 10260d0321e0SJeremy L Thompson diag->h_emodeout = emodeout; 10270d0321e0SJeremy L Thompson diag->numemodein = numemodein; 10280d0321e0SJeremy L Thompson diag->numemodeout = numemodeout; 10290d0321e0SJeremy L Thompson 10300d0321e0SJeremy L Thompson // Assemble kernel 10310d0321e0SJeremy L Thompson CeedInt nnodes, nqpts; 10320d0321e0SJeremy L Thompson ierr = CeedBasisGetNumNodes(basisin, &nnodes); CeedChkBackend(ierr); 10330d0321e0SJeremy L Thompson ierr = CeedBasisGetNumQuadraturePoints(basisin, &nqpts); CeedChkBackend(ierr); 10340d0321e0SJeremy L Thompson diag->nnodes = nnodes; 10350d0321e0SJeremy L Thompson ierr = CeedCompileHip(ceed, diagonalkernels, &diag->module, 5, 10360d0321e0SJeremy L Thompson "NUMEMODEIN", numemodein, 10370d0321e0SJeremy L Thompson "NUMEMODEOUT", numemodeout, 10380d0321e0SJeremy L Thompson "NNODES", nnodes, 10390d0321e0SJeremy L Thompson "NQPTS", nqpts, 10400d0321e0SJeremy L Thompson "NCOMP", ncomp 10410d0321e0SJeremy L Thompson ); CeedChk_Hip(ceed, ierr); 10420d0321e0SJeremy L Thompson ierr = CeedGetKernelHip(ceed, diag->module, "linearDiagonal", 10430d0321e0SJeremy L Thompson &diag->linearDiagonal); CeedChk_Hip(ceed, ierr); 10440d0321e0SJeremy L Thompson ierr = CeedGetKernelHip(ceed, diag->module, "linearPointBlockDiagonal", 10450d0321e0SJeremy L Thompson &diag->linearPointBlock); 10460d0321e0SJeremy L Thompson CeedChk_Hip(ceed, ierr); 10470d0321e0SJeremy L Thompson 10480d0321e0SJeremy L Thompson // Basis matrices 10490d0321e0SJeremy L Thompson const CeedInt qBytes = nqpts * sizeof(CeedScalar); 10500d0321e0SJeremy L Thompson const CeedInt iBytes = qBytes * nnodes; 10510d0321e0SJeremy L Thompson const CeedInt gBytes = qBytes * nnodes * dim; 10520d0321e0SJeremy L Thompson const CeedInt eBytes = sizeof(CeedEvalMode); 10530d0321e0SJeremy L Thompson const CeedScalar *interpin, *interpout, *gradin, *gradout; 10540d0321e0SJeremy L Thompson 10550d0321e0SJeremy L Thompson // CEED_EVAL_NONE 10560d0321e0SJeremy L Thompson CeedScalar *identity = NULL; 10570d0321e0SJeremy L Thompson bool evalNone = false; 10580d0321e0SJeremy L Thompson for (CeedInt i=0; i<numemodein; i++) 10590d0321e0SJeremy L Thompson evalNone = evalNone || (emodein[i] == CEED_EVAL_NONE); 10600d0321e0SJeremy L Thompson for (CeedInt i=0; i<numemodeout; i++) 10610d0321e0SJeremy L Thompson evalNone = evalNone || (emodeout[i] == CEED_EVAL_NONE); 10620d0321e0SJeremy L Thompson if (evalNone) { 10630d0321e0SJeremy L Thompson ierr = CeedCalloc(nqpts*nnodes, &identity); CeedChkBackend(ierr); 10640d0321e0SJeremy L Thompson for (CeedInt i=0; i<(nnodes<nqpts?nnodes:nqpts); i++) 10650d0321e0SJeremy L Thompson identity[i*nnodes+i] = 1.0; 10660d0321e0SJeremy L Thompson ierr = hipMalloc((void **)&diag->d_identity, iBytes); CeedChk_Hip(ceed, ierr); 10670d0321e0SJeremy L Thompson ierr = hipMemcpy(diag->d_identity, identity, iBytes, 10680d0321e0SJeremy L Thompson hipMemcpyHostToDevice); CeedChk_Hip(ceed, ierr); 10690d0321e0SJeremy L Thompson } 10700d0321e0SJeremy L Thompson 10710d0321e0SJeremy L Thompson // CEED_EVAL_INTERP 10720d0321e0SJeremy L Thompson ierr = CeedBasisGetInterp(basisin, &interpin); CeedChkBackend(ierr); 10730d0321e0SJeremy L Thompson ierr = hipMalloc((void **)&diag->d_interpin, iBytes); CeedChk_Hip(ceed, ierr); 10740d0321e0SJeremy L Thompson ierr = hipMemcpy(diag->d_interpin, interpin, iBytes, 10750d0321e0SJeremy L Thompson hipMemcpyHostToDevice); CeedChk_Hip(ceed, ierr); 10760d0321e0SJeremy L Thompson ierr = CeedBasisGetInterp(basisout, &interpout); CeedChkBackend(ierr); 10770d0321e0SJeremy L Thompson ierr = hipMalloc((void **)&diag->d_interpout, iBytes); CeedChk_Hip(ceed, ierr); 10780d0321e0SJeremy L Thompson ierr = hipMemcpy(diag->d_interpout, interpout, iBytes, 10790d0321e0SJeremy L Thompson hipMemcpyHostToDevice); CeedChk_Hip(ceed, ierr); 10800d0321e0SJeremy L Thompson 10810d0321e0SJeremy L Thompson // CEED_EVAL_GRAD 10820d0321e0SJeremy L Thompson ierr = CeedBasisGetGrad(basisin, &gradin); CeedChkBackend(ierr); 10830d0321e0SJeremy L Thompson ierr = hipMalloc((void **)&diag->d_gradin, gBytes); CeedChk_Hip(ceed, ierr); 10840d0321e0SJeremy L Thompson ierr = hipMemcpy(diag->d_gradin, gradin, gBytes, 10850d0321e0SJeremy L Thompson hipMemcpyHostToDevice); CeedChk_Hip(ceed, ierr); 10860d0321e0SJeremy L Thompson ierr = CeedBasisGetGrad(basisout, &gradout); CeedChkBackend(ierr); 10870d0321e0SJeremy L Thompson ierr = hipMalloc((void **)&diag->d_gradout, gBytes); CeedChk_Hip(ceed, ierr); 10880d0321e0SJeremy L Thompson ierr = hipMemcpy(diag->d_gradout, gradout, gBytes, 10890d0321e0SJeremy L Thompson hipMemcpyHostToDevice); CeedChk_Hip(ceed, ierr); 10900d0321e0SJeremy L Thompson 10910d0321e0SJeremy L Thompson // Arrays of emodes 10920d0321e0SJeremy L Thompson ierr = hipMalloc((void **)&diag->d_emodein, numemodein * eBytes); 10930d0321e0SJeremy L Thompson CeedChk_Hip(ceed, ierr); 10940d0321e0SJeremy L Thompson ierr = hipMemcpy(diag->d_emodein, emodein, numemodein * eBytes, 10950d0321e0SJeremy L Thompson hipMemcpyHostToDevice); CeedChk_Hip(ceed, ierr); 10960d0321e0SJeremy L Thompson ierr = hipMalloc((void **)&diag->d_emodeout, numemodeout * eBytes); 10970d0321e0SJeremy L Thompson CeedChk_Hip(ceed, ierr); 10980d0321e0SJeremy L Thompson ierr = hipMemcpy(diag->d_emodeout, emodeout, numemodeout * eBytes, 10990d0321e0SJeremy L Thompson hipMemcpyHostToDevice); CeedChk_Hip(ceed, ierr); 11000d0321e0SJeremy L Thompson 11010d0321e0SJeremy L Thompson // Restriction 11020d0321e0SJeremy L Thompson diag->diagrstr = rstrout; 11030d0321e0SJeremy L Thompson 11040d0321e0SJeremy L Thompson return CEED_ERROR_SUCCESS; 11050d0321e0SJeremy L Thompson } 11060d0321e0SJeremy L Thompson 11070d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 11080d0321e0SJeremy L Thompson // Assemble diagonal common code 11090d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 11100d0321e0SJeremy L Thompson static inline int CeedOperatorAssembleDiagonalCore_Hip(CeedOperator op, 11110d0321e0SJeremy L Thompson CeedVector assembled, CeedRequest *request, const bool pointBlock) { 11120d0321e0SJeremy L Thompson int ierr; 11130d0321e0SJeremy L Thompson Ceed ceed; 11140d0321e0SJeremy L Thompson ierr = CeedOperatorGetCeed(op, &ceed); CeedChkBackend(ierr); 11150d0321e0SJeremy L Thompson CeedOperator_Hip *impl; 11160d0321e0SJeremy L Thompson ierr = CeedOperatorGetData(op, &impl); CeedChkBackend(ierr); 11170d0321e0SJeremy L Thompson 11180d0321e0SJeremy L Thompson // Assemble QFunction 11190d0321e0SJeremy L Thompson CeedVector assembledqf; 11200d0321e0SJeremy L Thompson CeedElemRestriction rstr; 11210d0321e0SJeremy L Thompson ierr = CeedOperatorLinearAssembleQFunctionBuildOrUpdate(op, &assembledqf, 11220d0321e0SJeremy L Thompson &rstr, request); CeedChkBackend(ierr); 11230d0321e0SJeremy L Thompson ierr = CeedElemRestrictionDestroy(&rstr); CeedChkBackend(ierr); 11240d0321e0SJeremy L Thompson CeedScalar maxnorm = 0; 11250d0321e0SJeremy L Thompson ierr = CeedVectorNorm(assembledqf, CEED_NORM_MAX, &maxnorm); 11260d0321e0SJeremy L Thompson CeedChkBackend(ierr); 11270d0321e0SJeremy L Thompson 11280d0321e0SJeremy L Thompson // Setup 11290d0321e0SJeremy L Thompson if (!impl->diag) { 11300d0321e0SJeremy L Thompson ierr = CeedOperatorAssembleDiagonalSetup_Hip(op, pointBlock); 11310d0321e0SJeremy L Thompson CeedChkBackend(ierr); 11320d0321e0SJeremy L Thompson } 11330d0321e0SJeremy L Thompson CeedOperatorDiag_Hip *diag = impl->diag; 11340d0321e0SJeremy L Thompson assert(diag != NULL); 11350d0321e0SJeremy L Thompson 11360d0321e0SJeremy L Thompson // Restriction 11370d0321e0SJeremy L Thompson if (pointBlock && !diag->pbdiagrstr) { 11380d0321e0SJeremy L Thompson CeedElemRestriction pbdiagrstr; 11390d0321e0SJeremy L Thompson ierr = CreatePBRestriction(diag->diagrstr, &pbdiagrstr); CeedChkBackend(ierr); 11400d0321e0SJeremy L Thompson diag->pbdiagrstr = pbdiagrstr; 11410d0321e0SJeremy L Thompson } 11420d0321e0SJeremy L Thompson CeedElemRestriction diagrstr = pointBlock ? diag->pbdiagrstr : diag->diagrstr; 11430d0321e0SJeremy L Thompson 11440d0321e0SJeremy L Thompson // Create diagonal vector 11450d0321e0SJeremy L Thompson CeedVector elemdiag = pointBlock ? diag->pbelemdiag : diag->elemdiag; 11460d0321e0SJeremy L Thompson if (!elemdiag) { 11470d0321e0SJeremy L Thompson // Element diagonal vector 11480d0321e0SJeremy L Thompson ierr = CeedElemRestrictionCreateVector(diagrstr, NULL, &elemdiag); 11490d0321e0SJeremy L Thompson CeedChkBackend(ierr); 11500d0321e0SJeremy L Thompson if (pointBlock) 11510d0321e0SJeremy L Thompson diag->pbelemdiag = elemdiag; 11520d0321e0SJeremy L Thompson else 11530d0321e0SJeremy L Thompson diag->elemdiag = elemdiag; 11540d0321e0SJeremy L Thompson } 11550d0321e0SJeremy L Thompson ierr = CeedVectorSetValue(elemdiag, 0.0); CeedChkBackend(ierr); 11560d0321e0SJeremy L Thompson 11570d0321e0SJeremy L Thompson // Assemble element operator diagonals 11580d0321e0SJeremy L Thompson CeedScalar *elemdiagarray; 11590d0321e0SJeremy L Thompson const CeedScalar *assembledqfarray; 11600d0321e0SJeremy L Thompson ierr = CeedVectorGetArray(elemdiag, CEED_MEM_DEVICE, &elemdiagarray); 11610d0321e0SJeremy L Thompson CeedChkBackend(ierr); 11620d0321e0SJeremy L Thompson ierr = CeedVectorGetArrayRead(assembledqf, CEED_MEM_DEVICE, &assembledqfarray); 11630d0321e0SJeremy L Thompson CeedChkBackend(ierr); 11640d0321e0SJeremy L Thompson CeedInt nelem; 11650d0321e0SJeremy L Thompson ierr = CeedElemRestrictionGetNumElements(diagrstr, &nelem); 11660d0321e0SJeremy L Thompson CeedChkBackend(ierr); 11670d0321e0SJeremy L Thompson 11680d0321e0SJeremy L Thompson // Compute the diagonal of B^T D B 11690d0321e0SJeremy L Thompson int elemsPerBlock = 1; 11700d0321e0SJeremy L Thompson int grid = nelem/elemsPerBlock+((nelem/elemsPerBlock*elemsPerBlock<nelem)?1:0); 11710d0321e0SJeremy L Thompson void *args[] = {(void *) &nelem, (void *) &maxnorm, &diag->d_identity, 11720d0321e0SJeremy L Thompson &diag->d_interpin, &diag->d_gradin, &diag->d_interpout, 11730d0321e0SJeremy L Thompson &diag->d_gradout, &diag->d_emodein, &diag->d_emodeout, 11740d0321e0SJeremy L Thompson &assembledqfarray, &elemdiagarray 11750d0321e0SJeremy L Thompson }; 11760d0321e0SJeremy L Thompson if (pointBlock) { 11770d0321e0SJeremy L Thompson ierr = CeedRunKernelDimHip(ceed, diag->linearPointBlock, grid, 11780d0321e0SJeremy L Thompson diag->nnodes, 1, elemsPerBlock, args); 11790d0321e0SJeremy L Thompson CeedChkBackend(ierr); 11800d0321e0SJeremy L Thompson } else { 11810d0321e0SJeremy L Thompson ierr = CeedRunKernelDimHip(ceed, diag->linearDiagonal, grid, 11820d0321e0SJeremy L Thompson diag->nnodes, 1, elemsPerBlock, args); 11830d0321e0SJeremy L Thompson CeedChkBackend(ierr); 11840d0321e0SJeremy L Thompson } 11850d0321e0SJeremy L Thompson 11860d0321e0SJeremy L Thompson // Restore arrays 11870d0321e0SJeremy L Thompson ierr = CeedVectorRestoreArray(elemdiag, &elemdiagarray); CeedChkBackend(ierr); 11880d0321e0SJeremy L Thompson ierr = CeedVectorRestoreArrayRead(assembledqf, &assembledqfarray); 11890d0321e0SJeremy L Thompson CeedChkBackend(ierr); 11900d0321e0SJeremy L Thompson 11910d0321e0SJeremy L Thompson // Assemble local operator diagonal 11920d0321e0SJeremy L Thompson ierr = CeedElemRestrictionApply(diagrstr, CEED_TRANSPOSE, elemdiag, 11930d0321e0SJeremy L Thompson assembled, request); CeedChkBackend(ierr); 11940d0321e0SJeremy L Thompson 11950d0321e0SJeremy L Thompson // Cleanup 11960d0321e0SJeremy L Thompson ierr = CeedVectorDestroy(&assembledqf); CeedChkBackend(ierr); 11970d0321e0SJeremy L Thompson 11980d0321e0SJeremy L Thompson return CEED_ERROR_SUCCESS; 11990d0321e0SJeremy L Thompson } 12000d0321e0SJeremy L Thompson 12010d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 12020d0321e0SJeremy L Thompson // Assemble composite diagonal common code 12030d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 12040d0321e0SJeremy L Thompson static inline int CeedOperatorLinearAssembleAddDiagonalCompositeCore_Hip( 12050d0321e0SJeremy L Thompson CeedOperator op, CeedVector assembled, CeedRequest *request, 12060d0321e0SJeremy L Thompson const bool pointBlock) { 12070d0321e0SJeremy L Thompson int ierr; 12080d0321e0SJeremy L Thompson CeedInt numSub; 12090d0321e0SJeremy L Thompson CeedOperator *subOperators; 12100d0321e0SJeremy L Thompson ierr = CeedOperatorGetNumSub(op, &numSub); CeedChkBackend(ierr); 12110d0321e0SJeremy L Thompson ierr = CeedOperatorGetSubList(op, &subOperators); CeedChkBackend(ierr); 12120d0321e0SJeremy L Thompson for (CeedInt i = 0; i < numSub; i++) { 12130d0321e0SJeremy L Thompson ierr = CeedOperatorAssembleDiagonalCore_Hip(subOperators[i], assembled, 12140d0321e0SJeremy L Thompson request, pointBlock); CeedChkBackend(ierr); 12150d0321e0SJeremy L Thompson } 12160d0321e0SJeremy L Thompson return CEED_ERROR_SUCCESS; 12170d0321e0SJeremy L Thompson } 12180d0321e0SJeremy L Thompson 12190d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 12200d0321e0SJeremy L Thompson // Assemble Linear Diagonal 12210d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 12220d0321e0SJeremy L Thompson static int CeedOperatorLinearAssembleAddDiagonal_Hip(CeedOperator op, 12230d0321e0SJeremy L Thompson CeedVector assembled, CeedRequest *request) { 12240d0321e0SJeremy L Thompson int ierr; 12250d0321e0SJeremy L Thompson bool isComposite; 12260d0321e0SJeremy L Thompson ierr = CeedOperatorIsComposite(op, &isComposite); CeedChkBackend(ierr); 12270d0321e0SJeremy L Thompson if (isComposite) { 12280d0321e0SJeremy L Thompson return CeedOperatorLinearAssembleAddDiagonalCompositeCore_Hip(op, assembled, 12290d0321e0SJeremy L Thompson request, false); 12300d0321e0SJeremy L Thompson } else { 12310d0321e0SJeremy L Thompson return CeedOperatorAssembleDiagonalCore_Hip(op, assembled, request, false); 12320d0321e0SJeremy L Thompson } 12330d0321e0SJeremy L Thompson } 12340d0321e0SJeremy L Thompson 12350d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 12360d0321e0SJeremy L Thompson // Assemble Linear Point Block Diagonal 12370d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 12380d0321e0SJeremy L Thompson static int CeedOperatorLinearAssembleAddPointBlockDiagonal_Hip(CeedOperator op, 12390d0321e0SJeremy L Thompson CeedVector assembled, CeedRequest *request) { 12400d0321e0SJeremy L Thompson int ierr; 12410d0321e0SJeremy L Thompson bool isComposite; 12420d0321e0SJeremy L Thompson ierr = CeedOperatorIsComposite(op, &isComposite); CeedChkBackend(ierr); 12430d0321e0SJeremy L Thompson if (isComposite) { 12440d0321e0SJeremy L Thompson return CeedOperatorLinearAssembleAddDiagonalCompositeCore_Hip(op, assembled, 12450d0321e0SJeremy L Thompson request, true); 12460d0321e0SJeremy L Thompson } else { 12470d0321e0SJeremy L Thompson return CeedOperatorAssembleDiagonalCore_Hip(op, assembled, request, true); 12480d0321e0SJeremy L Thompson } 12490d0321e0SJeremy L Thompson } 12500d0321e0SJeremy L Thompson 1251*a835093fSnbeams //------------------------------------------------------------------------------ 1252*a835093fSnbeams // Matrix assembly kernel 1253*a835093fSnbeams //------------------------------------------------------------------------------ 1254*a835093fSnbeams // *INDENT-OFF* 1255*a835093fSnbeams static const char *assemblykernels = QUOTE( 1256*a835093fSnbeams extern "C" __launch_bounds__(BLOCK_SIZE) 1257*a835093fSnbeams __global__ void linearAssemble(const CeedScalar *B_in, const CeedScalar *B_out, 1258*a835093fSnbeams const CeedScalar *__restrict__ qf_array, 1259*a835093fSnbeams CeedScalar *__restrict__ values_array) { 1260*a835093fSnbeams 1261*a835093fSnbeams // This kernel assumes B_in and B_out have the same number of quadrature points and 1262*a835093fSnbeams // basis points. 1263*a835093fSnbeams // TODO: expand to more general cases 1264*a835093fSnbeams const int i = threadIdx.x; // The output row index of each B^TDB operation 1265*a835093fSnbeams const int l = threadIdx.y; // The output column index of each B^TDB operation 1266*a835093fSnbeams // such that we have (Bout^T)_ij D_jk Bin_kl = C_il 1267*a835093fSnbeams 1268*a835093fSnbeams // Strides for final output ordering, determined by the reference (interface) implementation of 1269*a835093fSnbeams // the symbolic assembly, slowest --> fastest: element, comp_in, comp_out, node_row, node_col 1270*a835093fSnbeams const CeedInt comp_out_stride = NNODES * NNODES; 1271*a835093fSnbeams const CeedInt comp_in_stride = comp_out_stride * NCOMP; 1272*a835093fSnbeams const CeedInt e_stride = comp_in_stride * NCOMP; 1273*a835093fSnbeams // Strides for QF array, slowest --> fastest: emode_in, comp_in, emode_out, comp_out, elem, qpt 1274*a835093fSnbeams const CeedInt qe_stride = NQPTS; 1275*a835093fSnbeams const CeedInt qcomp_out_stride = NELEM * qe_stride; 1276*a835093fSnbeams const CeedInt qemode_out_stride = qcomp_out_stride * NCOMP; 1277*a835093fSnbeams const CeedInt qcomp_in_stride = qemode_out_stride * NUMEMODEOUT; 1278*a835093fSnbeams const CeedInt qemode_in_stride = qcomp_in_stride * NCOMP; 1279*a835093fSnbeams 1280*a835093fSnbeams // Loop over each element (if necessary) 1281*a835093fSnbeams for (CeedInt e = blockIdx.x*blockDim.z + threadIdx.z; e < NELEM; 1282*a835093fSnbeams e += gridDim.x*blockDim.z) { 1283*a835093fSnbeams for (CeedInt comp_in = 0; comp_in < NCOMP; comp_in++) { 1284*a835093fSnbeams for (CeedInt comp_out = 0; comp_out < NCOMP; comp_out++) { 1285*a835093fSnbeams CeedScalar result = 0.0; 1286*a835093fSnbeams CeedInt qf_index_comp = qcomp_in_stride * comp_in + qcomp_out_stride * comp_out + qe_stride * e; 1287*a835093fSnbeams for (CeedInt emode_in = 0; emode_in < NUMEMODEIN; emode_in++) { 1288*a835093fSnbeams CeedInt b_in_index = emode_in * NQPTS * NNODES; 1289*a835093fSnbeams for (CeedInt emode_out = 0; emode_out < NUMEMODEOUT; emode_out++) { 1290*a835093fSnbeams CeedInt b_out_index = emode_out * NQPTS * NNODES; 1291*a835093fSnbeams CeedInt qf_index = qf_index_comp + qemode_out_stride * emode_out + qemode_in_stride * emode_in; 1292*a835093fSnbeams // Perform the B^T D B operation for this 'chunk' of D (the qf_array) 1293*a835093fSnbeams for (CeedInt j = 0; j < NQPTS; j++) { 1294*a835093fSnbeams result += B_out[b_out_index + j * NNODES + i] * qf_array[qf_index + j] * B_in[b_in_index + j * NNODES + l]; 1295*a835093fSnbeams } 1296*a835093fSnbeams 1297*a835093fSnbeams }// end of emode_out 1298*a835093fSnbeams } // end of emode_in 1299*a835093fSnbeams CeedInt val_index = comp_in_stride * comp_in + comp_out_stride * comp_out + e_stride * e + NNODES * i + l; 1300*a835093fSnbeams values_array[val_index] = result; 1301*a835093fSnbeams } // end of out component 1302*a835093fSnbeams } // end of in component 1303*a835093fSnbeams } // end of element loop 1304*a835093fSnbeams } 1305*a835093fSnbeams ); 1306*a835093fSnbeams // *INDENT-ON* 1307*a835093fSnbeams 1308*a835093fSnbeams //------------------------------------------------------------------------------ 1309*a835093fSnbeams // Single operator assembly setup 1310*a835093fSnbeams //------------------------------------------------------------------------------ 1311*a835093fSnbeams static int CeedSingleOperatorAssembleSetup_Hip(CeedOperator op) { 1312*a835093fSnbeams int ierr; 1313*a835093fSnbeams Ceed ceed; 1314*a835093fSnbeams ierr = CeedOperatorGetCeed(op, &ceed); CeedChkBackend(ierr); 1315*a835093fSnbeams CeedOperator_Hip *impl; 1316*a835093fSnbeams ierr = CeedOperatorGetData(op, &impl); CeedChkBackend(ierr); 1317*a835093fSnbeams 1318*a835093fSnbeams // Get intput and output fields 1319*a835093fSnbeams CeedInt num_input_fields, num_output_fields; 1320*a835093fSnbeams CeedOperatorField *input_fields; 1321*a835093fSnbeams CeedOperatorField *output_fields; 1322*a835093fSnbeams ierr = CeedOperatorGetFields(op, &num_input_fields, &input_fields, 1323*a835093fSnbeams &num_output_fields, &output_fields); CeedChk(ierr); 1324*a835093fSnbeams 1325*a835093fSnbeams // Determine active input basis eval mode 1326*a835093fSnbeams CeedQFunction qf; 1327*a835093fSnbeams ierr = CeedOperatorGetQFunction(op, &qf); CeedChk(ierr); 1328*a835093fSnbeams CeedQFunctionField *qf_fields; 1329*a835093fSnbeams ierr = CeedQFunctionGetFields(qf, NULL, &qf_fields, NULL, NULL); CeedChk(ierr); 1330*a835093fSnbeams // Note that the kernel will treat each dimension of a gradient action separately; 1331*a835093fSnbeams // i.e., when an active input has a CEED_EVAL_GRAD mode, num_emode_in will increment 1332*a835093fSnbeams // by dim. However, for the purposes of loading the B matrices, it will be treated 1333*a835093fSnbeams // as one mode, and we will load/copy the entire gradient matrix at once, so 1334*a835093fSnbeams // num_B_in_mats_to_load will be incremented by 1. 1335*a835093fSnbeams CeedInt num_emode_in = 0, dim = 1, num_B_in_mats_to_load = 0, size_B_in = 0; 1336*a835093fSnbeams CeedEvalMode *eval_mode_in = NULL; //will be of size num_B_in_mats_load 1337*a835093fSnbeams CeedBasis basis_in = NULL; 1338*a835093fSnbeams CeedInt nqpts, esize; 1339*a835093fSnbeams CeedElemRestriction rstr_in = NULL; 1340*a835093fSnbeams for (CeedInt i=0; i<num_input_fields; i++) { 1341*a835093fSnbeams CeedVector vec; 1342*a835093fSnbeams ierr = CeedOperatorFieldGetVector(input_fields[i], &vec); CeedChk(ierr); 1343*a835093fSnbeams if (vec == CEED_VECTOR_ACTIVE) { 1344*a835093fSnbeams ierr = CeedOperatorFieldGetBasis(input_fields[i], &basis_in); 1345*a835093fSnbeams CeedChk(ierr); 1346*a835093fSnbeams ierr = CeedBasisGetDimension(basis_in, &dim); CeedChk(ierr); 1347*a835093fSnbeams ierr = CeedBasisGetNumQuadraturePoints(basis_in, &nqpts); CeedChk(ierr); 1348*a835093fSnbeams ierr = CeedOperatorFieldGetElemRestriction(input_fields[i], &rstr_in); 1349*a835093fSnbeams ierr = CeedElemRestrictionGetElementSize(rstr_in, &esize); CeedChk(ierr); 1350*a835093fSnbeams CeedChk(ierr); 1351*a835093fSnbeams CeedEvalMode eval_mode; 1352*a835093fSnbeams ierr = CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode); 1353*a835093fSnbeams CeedChk(ierr); 1354*a835093fSnbeams if (eval_mode != CEED_EVAL_NONE) { 1355*a835093fSnbeams ierr = CeedRealloc(num_B_in_mats_to_load + 1, &eval_mode_in); CeedChk(ierr); 1356*a835093fSnbeams eval_mode_in[num_B_in_mats_to_load] = eval_mode; 1357*a835093fSnbeams num_B_in_mats_to_load += 1; 1358*a835093fSnbeams if (eval_mode == CEED_EVAL_GRAD) { 1359*a835093fSnbeams num_emode_in += dim; 1360*a835093fSnbeams size_B_in += dim * esize * nqpts; 1361*a835093fSnbeams } else { 1362*a835093fSnbeams num_emode_in +=1; 1363*a835093fSnbeams size_B_in += esize * nqpts; 1364*a835093fSnbeams } 1365*a835093fSnbeams } 1366*a835093fSnbeams } 1367*a835093fSnbeams } 1368*a835093fSnbeams 1369*a835093fSnbeams // Determine active output basis; basis_out and rstr_out only used if same as input, TODO 1370*a835093fSnbeams ierr = CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qf_fields); CeedChk(ierr); 1371*a835093fSnbeams CeedInt num_emode_out = 0, num_B_out_mats_to_load = 0, size_B_out = 0; 1372*a835093fSnbeams CeedEvalMode *eval_mode_out = NULL; 1373*a835093fSnbeams CeedBasis basis_out = NULL; 1374*a835093fSnbeams CeedElemRestriction rstr_out = NULL; 1375*a835093fSnbeams for (CeedInt i=0; i<num_output_fields; i++) { 1376*a835093fSnbeams CeedVector vec; 1377*a835093fSnbeams ierr = CeedOperatorFieldGetVector(output_fields[i], &vec); CeedChk(ierr); 1378*a835093fSnbeams if (vec == CEED_VECTOR_ACTIVE) { 1379*a835093fSnbeams ierr = CeedOperatorFieldGetBasis(output_fields[i], &basis_out); 1380*a835093fSnbeams CeedChk(ierr); 1381*a835093fSnbeams ierr = CeedOperatorFieldGetElemRestriction(output_fields[i], &rstr_out); 1382*a835093fSnbeams CeedChk(ierr); 1383*a835093fSnbeams if (rstr_out && rstr_out != rstr_in) 1384*a835093fSnbeams // LCOV_EXCL_START 1385*a835093fSnbeams return CeedError(ceed, CEED_ERROR_BACKEND, 1386*a835093fSnbeams "Multi-field non-composite operator assembly not supported"); 1387*a835093fSnbeams // LCOV_EXCL_STOP 1388*a835093fSnbeams CeedEvalMode eval_mode; 1389*a835093fSnbeams ierr = CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode); 1390*a835093fSnbeams CeedChk(ierr); 1391*a835093fSnbeams if (eval_mode != CEED_EVAL_NONE) { 1392*a835093fSnbeams ierr = CeedRealloc(num_B_out_mats_to_load + 1, &eval_mode_out); CeedChk(ierr); 1393*a835093fSnbeams eval_mode_out[num_B_out_mats_to_load] = eval_mode; 1394*a835093fSnbeams num_B_out_mats_to_load += 1; 1395*a835093fSnbeams if (eval_mode == CEED_EVAL_GRAD) { 1396*a835093fSnbeams num_emode_out += dim; 1397*a835093fSnbeams size_B_out += dim * esize * nqpts; 1398*a835093fSnbeams } else { 1399*a835093fSnbeams num_emode_out +=1; 1400*a835093fSnbeams size_B_out += esize * nqpts; 1401*a835093fSnbeams } 1402*a835093fSnbeams } 1403*a835093fSnbeams } 1404*a835093fSnbeams } 1405*a835093fSnbeams 1406*a835093fSnbeams if (num_emode_in == 0 || num_emode_out == 0) 1407*a835093fSnbeams // LCOV_EXCL_START 1408*a835093fSnbeams return CeedError(ceed, CEED_ERROR_UNSUPPORTED, 1409*a835093fSnbeams "Cannot assemble operator without inputs/outputs"); 1410*a835093fSnbeams // LCOV_EXCL_STOP 1411*a835093fSnbeams 1412*a835093fSnbeams CeedInt nelem, ncomp; 1413*a835093fSnbeams ierr = CeedElemRestrictionGetNumElements(rstr_in, &nelem); CeedChk(ierr); 1414*a835093fSnbeams ierr = CeedElemRestrictionGetNumComponents(rstr_in, &ncomp); CeedChk(ierr); 1415*a835093fSnbeams 1416*a835093fSnbeams ierr = CeedCalloc(1, &impl->asmb); CeedChkBackend(ierr); 1417*a835093fSnbeams CeedOperatorAssemble_Hip *asmb = impl->asmb; 1418*a835093fSnbeams asmb->nelem = nelem; 1419*a835093fSnbeams asmb->nnodes = esize; 1420*a835093fSnbeams 1421*a835093fSnbeams // Compile kernels 1422*a835093fSnbeams int elemsPerBlock = 1; 1423*a835093fSnbeams asmb->elemsPerBlock = elemsPerBlock; 1424*a835093fSnbeams // TODO: 1D threadblock version for larger elements 1425*a835093fSnbeams if (esize * esize * elemsPerBlock > 1024) 1426*a835093fSnbeams // LCOV_EXCL_START 1427*a835093fSnbeams return CeedError(ceed, CEED_ERROR_UNSUPPORTED, 1428*a835093fSnbeams "GPU assembly not currently available for elements of this size"); 1429*a835093fSnbeams // LCOV_EXCL_STOP 1430*a835093fSnbeams 1431*a835093fSnbeams ierr = CeedCompileHip(ceed, assemblykernels, &asmb->module, 7, 1432*a835093fSnbeams "NELEM", nelem, 1433*a835093fSnbeams "NUMEMODEIN", num_emode_in, 1434*a835093fSnbeams "NUMEMODEOUT", num_emode_out, 1435*a835093fSnbeams "NQPTS", nqpts, 1436*a835093fSnbeams "NNODES", esize, 1437*a835093fSnbeams "BLOCK_SIZE", esize * esize * elemsPerBlock, 1438*a835093fSnbeams "NCOMP", ncomp 1439*a835093fSnbeams ); CeedChk_Hip(ceed, ierr); 1440*a835093fSnbeams ierr = CeedGetKernelHip(ceed, asmb->module, "linearAssemble", 1441*a835093fSnbeams &asmb->linearAssemble); CeedChk_Hip(ceed, ierr); 1442*a835093fSnbeams 1443*a835093fSnbeams // Build 'full' B matrices (not 1D arrays used for tensor-product matrices) 1444*a835093fSnbeams const CeedScalar *interp_in, *grad_in; 1445*a835093fSnbeams ierr = CeedBasisGetInterp(basis_in, &interp_in); CeedChkBackend(ierr); 1446*a835093fSnbeams ierr = CeedBasisGetGrad(basis_in, &grad_in); CeedChkBackend(ierr); 1447*a835093fSnbeams 1448*a835093fSnbeams // Load into B_in, in order that they will be used in eval_mode 1449*a835093fSnbeams const CeedInt inBytes = size_B_in * sizeof(CeedScalar); 1450*a835093fSnbeams CeedInt mat_start = 0; 1451*a835093fSnbeams ierr = hipMalloc((void **) &asmb->d_B_in, inBytes); CeedChk_Hip(ceed, ierr); 1452*a835093fSnbeams for (int i = 0; i < num_B_in_mats_to_load; i++) { 1453*a835093fSnbeams CeedEvalMode eval_mode = eval_mode_in[i]; 1454*a835093fSnbeams if (eval_mode == CEED_EVAL_INTERP) { 1455*a835093fSnbeams ierr = hipMemcpy(&asmb->d_B_in[mat_start], interp_in, 1456*a835093fSnbeams esize * nqpts * sizeof(CeedScalar), 1457*a835093fSnbeams hipMemcpyHostToDevice); CeedChk_Hip(ceed, ierr); 1458*a835093fSnbeams mat_start += esize * nqpts; 1459*a835093fSnbeams } else if (eval_mode == CEED_EVAL_GRAD) { 1460*a835093fSnbeams ierr = hipMemcpy(asmb->d_B_in, grad_in, 1461*a835093fSnbeams dim * esize * nqpts * sizeof(CeedScalar), 1462*a835093fSnbeams hipMemcpyHostToDevice); CeedChk_Hip(ceed, ierr); 1463*a835093fSnbeams mat_start += dim * esize * nqpts; 1464*a835093fSnbeams } 1465*a835093fSnbeams } 1466*a835093fSnbeams 1467*a835093fSnbeams const CeedScalar *interp_out, *grad_out; 1468*a835093fSnbeams // Note that this function currently assumes 1 basis, so this should always be true 1469*a835093fSnbeams // for now 1470*a835093fSnbeams if (basis_out == basis_in) { 1471*a835093fSnbeams interp_out = interp_in; 1472*a835093fSnbeams grad_out = grad_in; 1473*a835093fSnbeams } else { 1474*a835093fSnbeams ierr = CeedBasisGetInterp(basis_out, &interp_out); CeedChkBackend(ierr); 1475*a835093fSnbeams ierr = CeedBasisGetGrad(basis_out, &grad_out); CeedChkBackend(ierr); 1476*a835093fSnbeams } 1477*a835093fSnbeams 1478*a835093fSnbeams // Load into B_out, in order that they will be used in eval_mode 1479*a835093fSnbeams const CeedInt outBytes = size_B_out * sizeof(CeedScalar); 1480*a835093fSnbeams mat_start = 0; 1481*a835093fSnbeams ierr = hipMalloc((void **) &asmb->d_B_out, outBytes); CeedChk_Hip(ceed, ierr); 1482*a835093fSnbeams for (int i = 0; i < num_B_out_mats_to_load; i++) { 1483*a835093fSnbeams CeedEvalMode eval_mode = eval_mode_out[i]; 1484*a835093fSnbeams if (eval_mode == CEED_EVAL_INTERP) { 1485*a835093fSnbeams ierr = hipMemcpy(&asmb->d_B_out[mat_start], interp_out, 1486*a835093fSnbeams esize * nqpts * sizeof(CeedScalar), 1487*a835093fSnbeams hipMemcpyHostToDevice); CeedChk_Hip(ceed, ierr); 1488*a835093fSnbeams mat_start += esize * nqpts; 1489*a835093fSnbeams } else if (eval_mode == CEED_EVAL_GRAD) { 1490*a835093fSnbeams ierr = hipMemcpy(&asmb->d_B_out[mat_start], grad_out, 1491*a835093fSnbeams dim * esize * nqpts * sizeof(CeedScalar), 1492*a835093fSnbeams hipMemcpyHostToDevice); CeedChk_Hip(ceed, ierr); 1493*a835093fSnbeams mat_start += dim * esize * nqpts; 1494*a835093fSnbeams } 1495*a835093fSnbeams } 1496*a835093fSnbeams return CEED_ERROR_SUCCESS; 1497*a835093fSnbeams } 1498*a835093fSnbeams 1499*a835093fSnbeams //------------------------------------------------------------------------------ 1500*a835093fSnbeams // Single operator assembly 1501*a835093fSnbeams //------------------------------------------------------------------------------ 1502*a835093fSnbeams static int CeedSingleOperatorAssemble_Hip(CeedOperator op, CeedInt offset, 1503*a835093fSnbeams CeedVector values) { 1504*a835093fSnbeams 1505*a835093fSnbeams int ierr; 1506*a835093fSnbeams Ceed ceed; 1507*a835093fSnbeams ierr = CeedOperatorGetCeed(op, &ceed); CeedChkBackend(ierr); 1508*a835093fSnbeams CeedOperator_Hip *impl; 1509*a835093fSnbeams ierr = CeedOperatorGetData(op, &impl); CeedChkBackend(ierr); 1510*a835093fSnbeams 1511*a835093fSnbeams // Setup 1512*a835093fSnbeams if (!impl->asmb) { 1513*a835093fSnbeams ierr = CeedSingleOperatorAssembleSetup_Hip(op); 1514*a835093fSnbeams CeedChkBackend(ierr); 1515*a835093fSnbeams } 1516*a835093fSnbeams 1517*a835093fSnbeams // Assemble QFunction 1518*a835093fSnbeams CeedVector assembled_qf; 1519*a835093fSnbeams CeedElemRestriction rstr_q; 1520*a835093fSnbeams ierr = CeedOperatorLinearAssembleQFunctionBuildOrUpdate( 1521*a835093fSnbeams op, &assembled_qf, &rstr_q, CEED_REQUEST_IMMEDIATE); CeedChk(ierr); 1522*a835093fSnbeams ierr = CeedElemRestrictionDestroy(&rstr_q); CeedChkBackend(ierr); 1523*a835093fSnbeams CeedScalar *values_array; 1524*a835093fSnbeams ierr = CeedVectorGetArrayWrite(values, CEED_MEM_DEVICE, &values_array); 1525*a835093fSnbeams CeedChkBackend(ierr); 1526*a835093fSnbeams values_array += offset; 1527*a835093fSnbeams const CeedScalar *qf_array; 1528*a835093fSnbeams ierr = CeedVectorGetArrayRead(assembled_qf, CEED_MEM_DEVICE, &qf_array); 1529*a835093fSnbeams CeedChkBackend(ierr); 1530*a835093fSnbeams 1531*a835093fSnbeams // Compute B^T D B 1532*a835093fSnbeams const CeedInt nelem = impl->asmb->nelem; 1533*a835093fSnbeams const CeedInt elemsPerBlock = impl->asmb->elemsPerBlock; 1534*a835093fSnbeams const CeedInt nnodes = impl->asmb->nnodes; 1535*a835093fSnbeams const CeedInt grid = nelem/elemsPerBlock+(( 1536*a835093fSnbeams nelem/elemsPerBlock*elemsPerBlock<nelem)?1:0); 1537*a835093fSnbeams void *args[] = {&impl->asmb->d_B_in, &impl->asmb->d_B_out, 1538*a835093fSnbeams &qf_array, &values_array 1539*a835093fSnbeams }; 1540*a835093fSnbeams ierr = CeedRunKernelDimHip(ceed, impl->asmb->linearAssemble, grid, 1541*a835093fSnbeams nnodes, nnodes, elemsPerBlock, args); 1542*a835093fSnbeams CeedChkBackend(ierr); 1543*a835093fSnbeams 1544*a835093fSnbeams 1545*a835093fSnbeams // Restore arrays 1546*a835093fSnbeams ierr = CeedVectorRestoreArray(values, &values_array); CeedChkBackend(ierr); 1547*a835093fSnbeams ierr = CeedVectorRestoreArrayRead(assembled_qf, &qf_array); 1548*a835093fSnbeams CeedChkBackend(ierr); 1549*a835093fSnbeams 1550*a835093fSnbeams // Cleanup 1551*a835093fSnbeams ierr = CeedVectorDestroy(&assembled_qf); CeedChkBackend(ierr); 1552*a835093fSnbeams 1553*a835093fSnbeams return CEED_ERROR_SUCCESS; 1554*a835093fSnbeams } 1555*a835093fSnbeams 1556*a835093fSnbeams //------------------------------------------------------------------------------ 1557*a835093fSnbeams // Assemble matrix data for COO matrix of assembled operator. 1558*a835093fSnbeams // The sparsity pattern is set by CeedOperatorLinearAssembleSymbolic. 1559*a835093fSnbeams // 1560*a835093fSnbeams // Note that this (and other assembly routines) currently assume only one 1561*a835093fSnbeams // active input restriction/basis per operator (could have multiple basis eval 1562*a835093fSnbeams // modes). 1563*a835093fSnbeams // TODO: allow multiple active input restrictions/basis objects 1564*a835093fSnbeams //------------------------------------------------------------------------------ 1565*a835093fSnbeams int CeedOperatorLinearAssemble_Hip(CeedOperator op, CeedVector values) { 1566*a835093fSnbeams 1567*a835093fSnbeams // As done in the default implementation, loop through suboperators 1568*a835093fSnbeams // for composite operators, or call single operator assembly otherwise 1569*a835093fSnbeams bool is_composite; 1570*a835093fSnbeams CeedInt ierr; 1571*a835093fSnbeams ierr = CeedOperatorIsComposite(op, &is_composite); CeedChk(ierr); 1572*a835093fSnbeams 1573*a835093fSnbeams CeedElemRestriction rstr; 1574*a835093fSnbeams CeedInt num_elem, elem_size, num_comp; 1575*a835093fSnbeams 1576*a835093fSnbeams CeedInt offset = 0; 1577*a835093fSnbeams if (is_composite) { 1578*a835093fSnbeams CeedInt num_suboperators; 1579*a835093fSnbeams ierr = CeedOperatorGetNumSub(op, &num_suboperators); CeedChk(ierr); 1580*a835093fSnbeams CeedOperator *sub_operators; 1581*a835093fSnbeams ierr = CeedOperatorGetSubList(op, &sub_operators); CeedChk(ierr); 1582*a835093fSnbeams for (int k = 0; k < num_suboperators; ++k) { 1583*a835093fSnbeams ierr = CeedSingleOperatorAssemble_Hip(sub_operators[k], offset, values); 1584*a835093fSnbeams CeedChk(ierr); 1585*a835093fSnbeams ierr = CeedOperatorGetActiveElemRestriction(sub_operators[k], &rstr); 1586*a835093fSnbeams CeedChk(ierr); 1587*a835093fSnbeams ierr = CeedElemRestrictionGetNumElements(rstr, &num_elem); CeedChk(ierr); 1588*a835093fSnbeams ierr = CeedElemRestrictionGetElementSize(rstr, &elem_size); CeedChk(ierr); 1589*a835093fSnbeams ierr = CeedElemRestrictionGetNumComponents(rstr, &num_comp); CeedChk(ierr); 1590*a835093fSnbeams offset += elem_size*num_comp * elem_size*num_comp * num_elem; 1591*a835093fSnbeams } 1592*a835093fSnbeams } else { 1593*a835093fSnbeams ierr = CeedSingleOperatorAssemble_Hip(op, offset, values); CeedChk(ierr); 1594*a835093fSnbeams } 1595*a835093fSnbeams 1596*a835093fSnbeams return CEED_ERROR_SUCCESS; 1597*a835093fSnbeams } 15980d0321e0SJeremy L Thompson 15990d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 16000d0321e0SJeremy L Thompson // Create operator 16010d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 16020d0321e0SJeremy L Thompson int CeedOperatorCreate_Hip(CeedOperator op) { 16030d0321e0SJeremy L Thompson int ierr; 16040d0321e0SJeremy L Thompson Ceed ceed; 16050d0321e0SJeremy L Thompson ierr = CeedOperatorGetCeed(op, &ceed); CeedChkBackend(ierr); 16060d0321e0SJeremy L Thompson CeedOperator_Hip *impl; 16070d0321e0SJeremy L Thompson 16080d0321e0SJeremy L Thompson ierr = CeedCalloc(1, &impl); CeedChkBackend(ierr); 16090d0321e0SJeremy L Thompson ierr = CeedOperatorSetData(op, impl); CeedChkBackend(ierr); 16100d0321e0SJeremy L Thompson 16110d0321e0SJeremy L Thompson ierr = CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleQFunction", 16120d0321e0SJeremy L Thompson CeedOperatorLinearAssembleQFunction_Hip); 16130d0321e0SJeremy L Thompson CeedChkBackend(ierr); 16140d0321e0SJeremy L Thompson ierr = CeedSetBackendFunction(ceed, "Operator", op, 16150d0321e0SJeremy L Thompson "LinearAssembleQFunctionUpdate", 16160d0321e0SJeremy L Thompson CeedOperatorLinearAssembleQFunctionUpdate_Hip); 16170d0321e0SJeremy L Thompson CeedChkBackend(ierr); 16180d0321e0SJeremy L Thompson ierr = CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleAddDiagonal", 16190d0321e0SJeremy L Thompson CeedOperatorLinearAssembleAddDiagonal_Hip); 16200d0321e0SJeremy L Thompson CeedChkBackend(ierr); 16210d0321e0SJeremy L Thompson ierr = CeedSetBackendFunction(ceed, "Operator", op, 16220d0321e0SJeremy L Thompson "LinearAssembleAddPointBlockDiagonal", 16230d0321e0SJeremy L Thompson CeedOperatorLinearAssembleAddPointBlockDiagonal_Hip); 16240d0321e0SJeremy L Thompson CeedChkBackend(ierr); 1625*a835093fSnbeams ierr = CeedSetBackendFunction(ceed, "Operator", op, 1626*a835093fSnbeams "LinearAssemble", CeedOperatorLinearAssemble_Hip); 1627*a835093fSnbeams CeedChkBackend(ierr); 16280d0321e0SJeremy L Thompson ierr = CeedSetBackendFunction(ceed, "Operator", op, "ApplyAdd", 16290d0321e0SJeremy L Thompson CeedOperatorApplyAdd_Hip); CeedChkBackend(ierr); 16300d0321e0SJeremy L Thompson ierr = CeedSetBackendFunction(ceed, "Operator", op, "Destroy", 16310d0321e0SJeremy L Thompson CeedOperatorDestroy_Hip); CeedChkBackend(ierr); 16320d0321e0SJeremy L Thompson return CEED_ERROR_SUCCESS; 16330d0321e0SJeremy L Thompson } 16340d0321e0SJeremy L Thompson 16350d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 16360d0321e0SJeremy L Thompson // Composite Operator Create 16370d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 16380d0321e0SJeremy L Thompson int CeedCompositeOperatorCreate_Hip(CeedOperator op) { 16390d0321e0SJeremy L Thompson int ierr; 16400d0321e0SJeremy L Thompson Ceed ceed; 16410d0321e0SJeremy L Thompson ierr = CeedOperatorGetCeed(op, &ceed); CeedChkBackend(ierr); 16420d0321e0SJeremy L Thompson 16430d0321e0SJeremy L Thompson ierr = CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleAddDiagonal", 16440d0321e0SJeremy L Thompson CeedOperatorLinearAssembleAddDiagonal_Hip); 16450d0321e0SJeremy L Thompson CeedChkBackend(ierr); 16460d0321e0SJeremy L Thompson ierr = CeedSetBackendFunction(ceed, "Operator", op, 16470d0321e0SJeremy L Thompson "LinearAssembleAddPointBlockDiagonal", 16480d0321e0SJeremy L Thompson CeedOperatorLinearAssembleAddPointBlockDiagonal_Hip); 16490d0321e0SJeremy L Thompson CeedChkBackend(ierr); 1650*a835093fSnbeams ierr = CeedSetBackendFunction(ceed, "Operator", op, 1651*a835093fSnbeams "LinearAssemble", CeedOperatorLinearAssemble_Hip); 1652*a835093fSnbeams CeedChkBackend(ierr); 16530d0321e0SJeremy L Thompson return CEED_ERROR_SUCCESS; 16540d0321e0SJeremy L Thompson } 16550d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 1656