1*3ac85a22SJunchao Zhang static char help[] = "Test matmat products with matdiagonal on gpus \n\n"; 2*3ac85a22SJunchao Zhang 3*3ac85a22SJunchao Zhang // Contributed by: Steven Dargaville 4*3ac85a22SJunchao Zhang 5*3ac85a22SJunchao Zhang #include <petscmat.h> 6*3ac85a22SJunchao Zhang 7*3ac85a22SJunchao Zhang int main(int argc, char **args) 8*3ac85a22SJunchao Zhang { 9*3ac85a22SJunchao Zhang const PetscInt inds[] = {0, 1}; 10*3ac85a22SJunchao Zhang PetscScalar avals[] = {2, 3, 5, 7}; 11*3ac85a22SJunchao Zhang Mat A, B_diag, B_aij_diag, result, result_diag; 12*3ac85a22SJunchao Zhang Vec diag; 13*3ac85a22SJunchao Zhang PetscBool equal = PETSC_FALSE; 14*3ac85a22SJunchao Zhang 15*3ac85a22SJunchao Zhang PetscFunctionBeginUser; 16*3ac85a22SJunchao Zhang PetscCall(PetscInitialize(&argc, &args, NULL, help)); 17*3ac85a22SJunchao Zhang 18*3ac85a22SJunchao Zhang // Create matrix to start 19*3ac85a22SJunchao Zhang PetscCall(MatCreateFromOptions(PETSC_COMM_WORLD, NULL, 1, 2, 2, 2, 2, &A)); 20*3ac85a22SJunchao Zhang PetscCall(MatSetUp(A)); 21*3ac85a22SJunchao Zhang PetscCall(MatSetValues(A, 2, inds, 2, inds, avals, INSERT_VALUES)); 22*3ac85a22SJunchao Zhang PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 23*3ac85a22SJunchao Zhang PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 24*3ac85a22SJunchao Zhang 25*3ac85a22SJunchao Zhang // Create a matdiagonal matrix 26*3ac85a22SJunchao Zhang // Will be the matching vec type as A 27*3ac85a22SJunchao Zhang PetscCall(MatCreateVecs(A, &diag, NULL)); 28*3ac85a22SJunchao Zhang PetscCall(VecSet(diag, 2.0)); 29*3ac85a22SJunchao Zhang PetscCall(MatCreateDiagonal(diag, &B_diag)); 30*3ac85a22SJunchao Zhang 31*3ac85a22SJunchao Zhang // Create the same matrix as the matdiagonal but in aij format 32*3ac85a22SJunchao Zhang PetscCall(MatCreateFromOptions(PETSC_COMM_WORLD, NULL, 1, 2, 2, 2, 2, &B_aij_diag)); 33*3ac85a22SJunchao Zhang PetscCall(MatSetUp(B_aij_diag)); 34*3ac85a22SJunchao Zhang PetscCall(MatDiagonalSet(B_aij_diag, diag, INSERT_VALUES)); 35*3ac85a22SJunchao Zhang PetscCall(MatAssemblyBegin(B_aij_diag, MAT_FINAL_ASSEMBLY)); 36*3ac85a22SJunchao Zhang PetscCall(MatAssemblyEnd(B_aij_diag, MAT_FINAL_ASSEMBLY)); 37*3ac85a22SJunchao Zhang PetscCall(VecDestroy(&diag)); 38*3ac85a22SJunchao Zhang 39*3ac85a22SJunchao Zhang // ~~~~~~~~~~~~~ 40*3ac85a22SJunchao Zhang // Do an initial matmatmult 41*3ac85a22SJunchao Zhang // A * B_aij_diag 42*3ac85a22SJunchao Zhang // and then 43*3ac85a22SJunchao Zhang // A * B_diag but just using MatDiagonalScale 44*3ac85a22SJunchao Zhang // ~~~~~~~~~~~~~ 45*3ac85a22SJunchao Zhang 46*3ac85a22SJunchao Zhang // aij * aij 47*3ac85a22SJunchao Zhang PetscCall(MatMatMult(A, B_aij_diag, MAT_INITIAL_MATRIX, 1.5, &result)); 48*3ac85a22SJunchao Zhang // PetscCall(MatView(result, PETSC_VIEWER_STDOUT_WORLD)); 49*3ac85a22SJunchao Zhang 50*3ac85a22SJunchao Zhang // aij * diagonal 51*3ac85a22SJunchao Zhang PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &result_diag)); 52*3ac85a22SJunchao Zhang PetscCall(MatDiagonalGetDiagonal(B_diag, &diag)); 53*3ac85a22SJunchao Zhang PetscCall(MatDiagonalScale(result_diag, NULL, diag)); 54*3ac85a22SJunchao Zhang PetscCall(MatDiagonalRestoreDiagonal(B_diag, &diag)); 55*3ac85a22SJunchao Zhang // PetscCall(MatView(result_diag, PETSC_VIEWER_STDOUT_WORLD)); 56*3ac85a22SJunchao Zhang 57*3ac85a22SJunchao Zhang PetscCall(MatEqual(result, result_diag, &equal)); 58*3ac85a22SJunchao Zhang PetscCheck(equal, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "MatMatMult and MatDiagonalScale do not give the same result"); 59*3ac85a22SJunchao Zhang 60*3ac85a22SJunchao Zhang // ~~~~~~~~~~~~~ 61*3ac85a22SJunchao Zhang // Now let's modify the diagonal and do it again with "reuse" 62*3ac85a22SJunchao Zhang // ~~~~~~~~~~~~~ 63*3ac85a22SJunchao Zhang PetscCall(MatDiagonalGetDiagonal(B_diag, &diag)); 64*3ac85a22SJunchao Zhang PetscCall(VecSet(diag, 3.0)); 65*3ac85a22SJunchao Zhang PetscCall(MatDiagonalSet(B_aij_diag, diag, INSERT_VALUES)); 66*3ac85a22SJunchao Zhang PetscCall(MatDiagonalRestoreDiagonal(B_diag, &diag)); 67*3ac85a22SJunchao Zhang 68*3ac85a22SJunchao Zhang // aij * aij 69*3ac85a22SJunchao Zhang PetscCall(MatMatMult(A, B_aij_diag, MAT_REUSE_MATRIX, 1.5, &result)); 70*3ac85a22SJunchao Zhang 71*3ac85a22SJunchao Zhang // aij * diagonal 72*3ac85a22SJunchao Zhang PetscCall(MatCopy(A, result_diag, SAME_NONZERO_PATTERN)); 73*3ac85a22SJunchao Zhang PetscCall(MatDiagonalGetDiagonal(B_diag, &diag)); 74*3ac85a22SJunchao Zhang PetscCall(MatDiagonalScale(result_diag, NULL, diag)); 75*3ac85a22SJunchao Zhang PetscCall(MatDiagonalRestoreDiagonal(B_diag, &diag)); 76*3ac85a22SJunchao Zhang 77*3ac85a22SJunchao Zhang PetscCall(MatEqual(result, result_diag, &equal)); 78*3ac85a22SJunchao Zhang PetscCheck(equal, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "MatMatMult and MatDiagonalScale do not give the same result"); 79*3ac85a22SJunchao Zhang 80*3ac85a22SJunchao Zhang PetscCall(MatDestroy(&A)); 81*3ac85a22SJunchao Zhang PetscCall(MatDestroy(&B_diag)); 82*3ac85a22SJunchao Zhang PetscCall(MatDestroy(&B_aij_diag)); 83*3ac85a22SJunchao Zhang PetscCall(MatDestroy(&result)); 84*3ac85a22SJunchao Zhang PetscCall(MatDestroy(&result_diag)); 85*3ac85a22SJunchao Zhang PetscCall(VecDestroy(&diag)); 86*3ac85a22SJunchao Zhang 87*3ac85a22SJunchao Zhang PetscCall(PetscFinalize()); 88*3ac85a22SJunchao Zhang return 0; 89*3ac85a22SJunchao Zhang } 90*3ac85a22SJunchao Zhang 91*3ac85a22SJunchao Zhang /*TEST 92*3ac85a22SJunchao Zhang test: 93*3ac85a22SJunchao Zhang requires: kokkos_kernels 94*3ac85a22SJunchao Zhang args: -mat_type aijkokkos 95*3ac85a22SJunchao Zhang output_file: output/empty.out 96*3ac85a22SJunchao Zhang TEST*/ 97