1 static char help[] = "Test MatSetValuesCOO for MPIAIJ and its subclasses \n\n"; 2 3 #include <petscmat.h> 4 int main(int argc, char **args) { 5 Mat A, B; 6 PetscInt k; 7 const PetscInt M = 18, N = 18; 8 PetscMPIInt rank, size; 9 PetscBool equal; 10 PetscScalar *vals; 11 PetscBool flg = PETSC_FALSE; 12 13 /* Construct 18 x 18 matrices, which are big enough to have complex communication patterns but still small enough for debugging */ 14 PetscInt i0[] = {7, 7, 8, 8, 9, 16, 17, 9, 10, 1, 1, -2, 2, 3, 3, 14, 4, 5, 10, 13, 9, 9, 10, 1, 0, 0, 5, 5, 6, 6, 13, 13, 14, -14, 4, 4, 5, 11, 11, 12, 15, 15, 16}; 15 PetscInt j0[] = {1, 6, 2, 4, 10, 15, 13, 16, 11, 2, 7, 3, 8, 4, 9, 13, 5, 2, 15, 14, 10, 16, 11, 2, 0, 1, 5, -11, 0, 7, 15, 17, 11, 13, 4, 8, 2, 12, 17, 13, 3, 16, 9}; 16 17 PetscInt i1[] = {8, 5, 15, 16, 6, 13, 4, 17, 8, 9, 9, 10, -6, 12, 7, 3, -4, 1, 1, 2, 5, 5, 6, 14, 17, 8, 9, 9, 10, 4, 5, 10, 11, 1, 2}; 18 PetscInt j1[] = {2, 3, 16, 9, 5, 17, 1, 13, 4, 10, 16, 11, -5, 12, 1, 7, -1, 2, 7, 3, 6, 11, 0, 11, 13, 4, 10, 16, 11, 8, -2, 15, 12, 7, 3}; 19 20 PetscInt i2[] = {3, 4, 1, 10, 0, 1, 1, 2, 1, 1, 2, 2, 3, 3, 4, 4, 1, 2, 5, 5, 6, 4, 17, 0, 1, 1, 8, 5, 5, 6, 4, 7, 8, 5}; 21 PetscInt j2[] = {7, 1, 2, 11, 5, 2, 7, 3, 2, 7, 3, 8, 4, 9, 3, 5, 7, 3, 6, 11, 0, 1, 13, 5, 2, 7, 4, 6, 11, 0, 1, 3, 4, 2}; 22 23 struct { 24 PetscInt *i, *j, n; 25 } coo[3] = { 26 {i0, j0, sizeof(i0) / sizeof(PetscInt)}, 27 {i1, j1, sizeof(i1) / sizeof(PetscInt)}, 28 {i2, j2, sizeof(i2) / sizeof(PetscInt)} 29 }; 30 31 PetscFunctionBeginUser; 32 PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 33 PetscCall(PetscOptionsGetBool(NULL, NULL, "-ignore_remote", &flg, NULL)); 34 PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 35 PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 36 37 PetscCheck(size <= 3, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This test requires at most 3 processes"); 38 39 PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 40 PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, M, N)); 41 PetscCall(MatSetType(A, MATAIJ)); 42 PetscCall(MatSeqAIJSetPreallocation(A, 2, NULL)); 43 PetscCall(MatMPIAIJSetPreallocation(A, 2, NULL, 2, NULL)); 44 PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE)); 45 PetscCall(MatSetOption(A, MAT_IGNORE_OFF_PROC_ENTRIES, flg)); 46 47 PetscCall(PetscMalloc1(coo[rank].n, &vals)); 48 for (k = 0; k < coo[rank].n; k++) { 49 vals[k] = coo[rank].j[k]; 50 PetscCall(MatSetValue(A, coo[rank].i[k], coo[rank].j[k], vals[k], ADD_VALUES)); 51 } 52 PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 53 PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 54 55 PetscCall(MatCreate(PETSC_COMM_WORLD, &B)); 56 PetscCall(MatSetSizes(B, PETSC_DECIDE, PETSC_DECIDE, M, N)); 57 PetscCall(MatSetFromOptions(B)); 58 PetscCall(MatSetOption(B, MAT_IGNORE_OFF_PROC_ENTRIES, flg)); 59 PetscCall(MatSetPreallocationCOO(B, coo[rank].n, coo[rank].i, coo[rank].j)); 60 61 PetscCall(MatSetValuesCOO(B, vals, ADD_VALUES)); 62 PetscCall(MatEqual(A, B, &equal)); 63 64 if (!equal) { 65 PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD)); 66 PetscCall(MatView(B, PETSC_VIEWER_STDOUT_WORLD)); 67 SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_PLIB, "MatSetValuesCOO() failed"); 68 } 69 70 PetscCall(PetscFree(vals)); 71 PetscCall(MatDestroy(&A)); 72 PetscCall(MatDestroy(&B)); 73 74 PetscCall(PetscFinalize()); 75 return 0; 76 } 77 78 /*TEST 79 80 testset: 81 output_file: output/ex254_1.out 82 nsize: {{1 2 3}} 83 args: -ignore_remote {{0 1}} 84 85 test: 86 suffix: kokkos 87 requires: kokkos_kernels 88 args: -mat_type aijkokkos 89 90 test: 91 suffix: cuda 92 requires: cuda 93 args: -mat_type aijcusparse 94 95 test: 96 suffix: aij 97 args: -mat_type aij 98 99 testset: 100 # MATHYPRE does not support MAT_IGNORE_OFF_PROC_ENTRIES yet 101 output_file: output/ex254_1.out 102 nsize: {{1 2 3}} 103 suffix: hypre 104 requires: hypre kokkos_kernels 105 args: -ignore_remote 0 -mat_type hypre 106 107 TEST*/ 108