1bbfc976bSStefano Zampini static char help[] = "Tests MatCreateDenseCUDA(), MatDenseCUDAPlaceArray(), MatDenseCUDAReplaceArray(), MatDenseCUDAResetArray()\n"; 2bbfc976bSStefano Zampini 3bbfc976bSStefano Zampini #include <petscmat.h> 4f5ae8b2eSJunchao Zhang #include <petscpkg_version.h> 5bbfc976bSStefano Zampini 6d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMult_S(Mat S, Vec x, Vec y) 7d71ae5a4SJacob Faibussowitsch { 8bbfc976bSStefano Zampini Mat A; 9bbfc976bSStefano Zampini 10d9b9dddeSStefano Zampini PetscFunctionBeginUser; 119566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(S, &A)); 129566063dSJacob Faibussowitsch PetscCall(MatMult(A, x, y)); 133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14bbfc976bSStefano Zampini } 15bbfc976bSStefano Zampini 1654da937aSStefano Zampini static PetscBool test_cusparse_transgen = PETSC_FALSE; 1754da937aSStefano Zampini 18d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMultTranspose_S(Mat S, Vec x, Vec y) 19d71ae5a4SJacob Faibussowitsch { 20d9b9dddeSStefano Zampini Mat A; 21d9b9dddeSStefano Zampini 22d9b9dddeSStefano Zampini PetscFunctionBeginUser; 239566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(S, &A)); 249566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A, x, y)); 2554da937aSStefano Zampini 2654da937aSStefano Zampini /* alternate transgen true and false to test code logic */ 279566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_FORM_EXPLICIT_TRANSPOSE, test_cusparse_transgen)); 2854da937aSStefano Zampini test_cusparse_transgen = (PetscBool)!test_cusparse_transgen; 293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 30d9b9dddeSStefano Zampini } 31d9b9dddeSStefano Zampini 32d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 33d71ae5a4SJacob Faibussowitsch { 34bbfc976bSStefano Zampini Mat A, B, C, S; 35bbfc976bSStefano Zampini Vec t, v; 36bbfc976bSStefano Zampini PetscScalar *vv, *aa; 37f5ae8b2eSJunchao Zhang // We met a mysterious cudaErrorMisalignedAddress error on some systems with cuda-12.0,1 but not 38f5ae8b2eSJunchao Zhang // with prior cuda-11.2,3,7,8 versions. Making nloc an even number somehow 'fixes' the problem. 39f5ae8b2eSJunchao Zhang // See more at https://gitlab.com/petsc/petsc/-/merge_requests/6225 40f5ae8b2eSJunchao Zhang #if PETSC_PKG_CUDA_VERSION_GE(12, 0, 0) 41f5ae8b2eSJunchao Zhang PetscInt n = 32; 42f5ae8b2eSJunchao Zhang #else 43f5ae8b2eSJunchao Zhang PetscInt n = 30; 44f5ae8b2eSJunchao Zhang #endif 45f5ae8b2eSJunchao Zhang PetscInt k = 6, l = 0, i, Istart, Iend, nloc, bs, test = 1; 46bbfc976bSStefano Zampini PetscBool flg, reset, use_shell = PETSC_FALSE; 47bbfc976bSStefano Zampini VecType vtype; 48bbfc976bSStefano Zampini 49327415f7SBarry Smith PetscFunctionBeginUser; 509566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 519566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL)); 529566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-k", &k, NULL)); 539566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-l", &l, NULL)); 549566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-test", &test, NULL)); 559566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-use_shell", &use_shell, NULL)); 5608401ef6SPierre Jolivet PetscCheck(k >= 0, PETSC_COMM_WORLD, PETSC_ERR_USER, "k %" PetscInt_FMT " must be positive", k); 5708401ef6SPierre Jolivet PetscCheck(l >= 0, PETSC_COMM_WORLD, PETSC_ERR_USER, "l %" PetscInt_FMT " must be positive", l); 5808401ef6SPierre Jolivet PetscCheck(l <= k, PETSC_COMM_WORLD, PETSC_ERR_USER, "l %" PetscInt_FMT " must be smaller or equal than k %" PetscInt_FMT, l, k); 59bbfc976bSStefano Zampini 60bbfc976bSStefano Zampini /* sparse matrix */ 619566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 629566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, n, n)); 639566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATAIJCUSPARSE)); 649566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(A, "A_")); 659566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 669566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 67bbfc976bSStefano Zampini 68d9b9dddeSStefano Zampini /* test special case for SeqAIJCUSPARSE to generate explicit transpose (not default) */ 699566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_FORM_EXPLICIT_TRANSPOSE, PETSC_TRUE)); 70d9b9dddeSStefano Zampini 719566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &Istart, &Iend)); 72bbfc976bSStefano Zampini for (i = Istart; i < Iend; i++) { 739566063dSJacob Faibussowitsch if (i > 0) PetscCall(MatSetValue(A, i, i - 1, -1.0, INSERT_VALUES)); 749566063dSJacob Faibussowitsch if (i < n - 1) PetscCall(MatSetValue(A, i, i + 1, -1.0, INSERT_VALUES)); 759566063dSJacob Faibussowitsch PetscCall(MatSetValue(A, i, i, 2.0, INSERT_VALUES)); 76bbfc976bSStefano Zampini } 779566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 789566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 79bbfc976bSStefano Zampini 80bbfc976bSStefano Zampini /* template vector */ 819566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, NULL, &t)); 829566063dSJacob Faibussowitsch PetscCall(VecGetType(t, &vtype)); 83bbfc976bSStefano Zampini 84bbfc976bSStefano Zampini /* long vector, contains the stacked columns of an nxk dense matrix */ 859566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(t, &nloc)); 869566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(t, &bs)); 879566063dSJacob Faibussowitsch PetscCall(VecCreate(PetscObjectComm((PetscObject)t), &v)); 889566063dSJacob Faibussowitsch PetscCall(VecSetType(v, vtype)); 899566063dSJacob Faibussowitsch PetscCall(VecSetSizes(v, k * nloc, k * n)); 909566063dSJacob Faibussowitsch PetscCall(VecSetBlockSize(v, bs)); 919566063dSJacob Faibussowitsch PetscCall(VecSetRandom(v, NULL)); 92bbfc976bSStefano Zampini 93bbfc976bSStefano Zampini /* dense matrix that contains the columns of v */ 949566063dSJacob Faibussowitsch PetscCall(VecCUDAGetArray(v, &vv)); 95bbfc976bSStefano Zampini 96bbfc976bSStefano Zampini /* test few cases for MatDenseCUDA handling pointers */ 97bbfc976bSStefano Zampini switch (test) { 98bbfc976bSStefano Zampini case 1: 999566063dSJacob Faibussowitsch PetscCall(MatCreateDenseCUDA(PetscObjectComm((PetscObject)v), nloc, PETSC_DECIDE, n, k - l, vv, &B)); /* pass a pointer to avoid allocation of storage */ 1009566063dSJacob Faibussowitsch PetscCall(MatDenseCUDAReplaceArray(B, NULL)); /* replace with a null pointer, the value after BVRestoreMat */ 1019566063dSJacob Faibussowitsch PetscCall(MatDenseCUDAPlaceArray(B, vv + l * nloc)); /* set the actual pointer */ 102bbfc976bSStefano Zampini reset = PETSC_TRUE; 103bbfc976bSStefano Zampini break; 104bbfc976bSStefano Zampini case 2: 1059566063dSJacob Faibussowitsch PetscCall(MatCreateDenseCUDA(PetscObjectComm((PetscObject)v), nloc, PETSC_DECIDE, n, k - l, NULL, &B)); 1069566063dSJacob Faibussowitsch PetscCall(MatDenseCUDAPlaceArray(B, vv + l * nloc)); /* set the actual pointer */ 107bbfc976bSStefano Zampini reset = PETSC_TRUE; 108bbfc976bSStefano Zampini break; 109bbfc976bSStefano Zampini default: 1109566063dSJacob Faibussowitsch PetscCall(MatCreateDenseCUDA(PetscObjectComm((PetscObject)v), nloc, PETSC_DECIDE, n, k - l, vv + l * nloc, &B)); 111bbfc976bSStefano Zampini reset = PETSC_FALSE; 112bbfc976bSStefano Zampini break; 113bbfc976bSStefano Zampini } 114bbfc976bSStefano Zampini 115bbfc976bSStefano Zampini /* Test MatMatMult */ 116bbfc976bSStefano Zampini if (use_shell) { 117*f332b1cbSPierre Jolivet /* we could have called the general converter below, but we explicitly set the operations 11854da937aSStefano Zampini ourselves to test MatProductSymbolic_X_Dense, MatProductNumeric_X_Dense code */ 1199566063dSJacob Faibussowitsch /* PetscCall(MatConvert(A,MATSHELL,MAT_INITIAL_MATRIX,&S)); */ 1209566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PetscObjectComm((PetscObject)v), nloc, nloc, n, n, A, &S)); 1219566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(S, MATOP_MULT, (void (*)(void))MatMult_S)); 1229566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(S, MATOP_MULT_TRANSPOSE, (void (*)(void))MatMultTranspose_S)); 1239566063dSJacob Faibussowitsch PetscCall(MatShellSetVecType(S, vtype)); 124bbfc976bSStefano Zampini } else { 1259566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)A)); 126bbfc976bSStefano Zampini S = A; 127bbfc976bSStefano Zampini } 128bbfc976bSStefano Zampini 1299566063dSJacob Faibussowitsch PetscCall(MatCreateDenseCUDA(PetscObjectComm((PetscObject)v), nloc, PETSC_DECIDE, n, k - l, NULL, &C)); 130d9b9dddeSStefano Zampini 131d9b9dddeSStefano Zampini /* test MatMatMult */ 1329566063dSJacob Faibussowitsch PetscCall(MatProductCreateWithMat(S, B, NULL, C)); 1339566063dSJacob Faibussowitsch PetscCall(MatProductSetType(C, MATPRODUCT_AB)); 1349566063dSJacob Faibussowitsch PetscCall(MatProductSetFromOptions(C)); 1359566063dSJacob Faibussowitsch PetscCall(MatProductSymbolic(C)); 1369566063dSJacob Faibussowitsch PetscCall(MatProductNumeric(C)); 1379566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual(S, B, C, 10, &flg)); 1389566063dSJacob Faibussowitsch if (!flg) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error MatMatMult\n")); 139bbfc976bSStefano Zampini 140d9b9dddeSStefano Zampini /* test MatTransposeMatMult */ 1419566063dSJacob Faibussowitsch PetscCall(MatProductCreateWithMat(S, B, NULL, C)); 1429566063dSJacob Faibussowitsch PetscCall(MatProductSetType(C, MATPRODUCT_AtB)); 1439566063dSJacob Faibussowitsch PetscCall(MatProductSetFromOptions(C)); 1449566063dSJacob Faibussowitsch PetscCall(MatProductSymbolic(C)); 1459566063dSJacob Faibussowitsch PetscCall(MatProductNumeric(C)); 1469566063dSJacob Faibussowitsch PetscCall(MatTransposeMatMultEqual(S, B, C, 10, &flg)); 1479566063dSJacob Faibussowitsch if (!flg) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error MatTransposeMatMult\n")); 148d9b9dddeSStefano Zampini 1499566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 1509566063dSJacob Faibussowitsch PetscCall(MatDestroy(&S)); 151bbfc976bSStefano Zampini 152bbfc976bSStefano Zampini /* finished using B */ 153cd3f9d89SJunchao Zhang PetscCall(MatDenseGetArrayAndMemType(B, &aa, NULL)); 15408401ef6SPierre Jolivet PetscCheck(vv == aa - l * nloc, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Wrong array"); 155cd3f9d89SJunchao Zhang PetscCall(MatDenseRestoreArrayAndMemType(B, &aa)); 1561baa6e33SBarry Smith if (reset) PetscCall(MatDenseCUDAResetArray(B)); 1579566063dSJacob Faibussowitsch PetscCall(VecCUDARestoreArray(v, &vv)); 158bbfc976bSStefano Zampini 159bbfc976bSStefano Zampini if (test == 1) { 160cd3f9d89SJunchao Zhang PetscCall(MatDenseGetArrayAndMemType(B, &aa, NULL)); 16128b400f6SJacob Faibussowitsch PetscCheck(!aa, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Expected a null pointer"); 162cd3f9d89SJunchao Zhang PetscCall(MatDenseRestoreArrayAndMemType(B, &aa)); 163bbfc976bSStefano Zampini } 164bbfc976bSStefano Zampini 165bbfc976bSStefano Zampini /* free work space */ 1669566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 1679566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1689566063dSJacob Faibussowitsch PetscCall(VecDestroy(&t)); 1699566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v)); 1709566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 171b122ec5aSJacob Faibussowitsch return 0; 172bbfc976bSStefano Zampini } 173bbfc976bSStefano Zampini 174bbfc976bSStefano Zampini /*TEST 175bbfc976bSStefano Zampini 176bbfc976bSStefano Zampini build: 177bbfc976bSStefano Zampini requires: cuda 178bbfc976bSStefano Zampini 179bbfc976bSStefano Zampini test: 180bbfc976bSStefano Zampini requires: cuda 181bbfc976bSStefano Zampini suffix: 1 182bbfc976bSStefano Zampini nsize: {{1 2}} 1830a7ee96eSStefano Zampini args: -A_mat_type {{aij aijcusparse}} -test {{0 1 2}} -k 6 -l {{0 5}} -use_shell {{0 1}} 184bbfc976bSStefano Zampini 185bbfc976bSStefano Zampini TEST*/ 186