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 { 6 Mat A, B, C; 7 PetscInt k; 8 const PetscInt M = 18, N = 18; 9 PetscBool equal; 10 PetscScalar *vals; 11 PetscBool flg = PETSC_FALSE, freecoo = PETSC_FALSE; 12 PetscInt ncoos = 1; 13 14 // clang-format off 15 /* Construct 18 x 18 matrices, which are big enough to have complex communication patterns but still small enough for debugging */ 16 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}; 17 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}; 18 19 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}; 20 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}; 21 22 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}; 23 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}; 24 // clang-format on 25 26 typedef struct { 27 PetscInt *i, *j, n; 28 } coo_data; 29 30 coo_data coos[3] = { 31 {i0, j0, sizeof(i0) / sizeof(PetscInt)}, 32 {i1, j1, sizeof(i1) / sizeof(PetscInt)}, 33 {i2, j2, sizeof(i2) / sizeof(PetscInt)} 34 }; 35 coo_data mycoo; 36 37 PetscFunctionBeginUser; 38 PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 39 PetscCall(PetscOptionsGetBool(NULL, NULL, "-ignore_remote", &flg, NULL)); 40 PetscCall(PetscOptionsGetInt(NULL, NULL, "-ncoos", &ncoos, NULL)); 41 42 mycoo.n = 0; 43 if (ncoos > 1) { 44 PetscLayout map; 45 46 freecoo = PETSC_TRUE; 47 PetscCall(PetscLayoutCreate(PETSC_COMM_WORLD, &map)); 48 PetscCall(PetscLayoutSetSize(map, ncoos)); 49 PetscCall(PetscLayoutSetUp(map)); 50 PetscCall(PetscLayoutGetLocalSize(map, &ncoos)); 51 for (PetscInt i = 0; i < ncoos; i++) mycoo.n += coos[i % 3].n; 52 PetscCall(PetscMalloc2(mycoo.n, &mycoo.i, mycoo.n, &mycoo.j)); 53 mycoo.n = 0; 54 for (PetscInt i = 0; i < ncoos; i++) { 55 PetscCall(PetscArraycpy(mycoo.i + mycoo.n, coos[i % 3].i, coos[i % 3].n)); 56 PetscCall(PetscArraycpy(mycoo.j + mycoo.n, coos[i % 3].j, coos[i % 3].n)); 57 mycoo.n += coos[i % 3].n; 58 } 59 PetscCall(PetscLayoutDestroy(&map)); 60 } else if (ncoos == 1 && PetscGlobalRank < 3) mycoo = coos[PetscGlobalRank]; 61 62 PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 63 PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, M, N)); 64 PetscCall(MatSetType(A, MATAIJ)); 65 PetscCall(MatSeqAIJSetPreallocation(A, 2, NULL)); 66 PetscCall(MatMPIAIJSetPreallocation(A, 2, NULL, 2, NULL)); 67 PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE)); 68 PetscCall(MatSetOption(A, MAT_IGNORE_OFF_PROC_ENTRIES, flg)); 69 70 PetscCall(PetscMalloc1(mycoo.n, &vals)); 71 for (k = 0; k < mycoo.n; k++) { 72 vals[k] = mycoo.j[k]; 73 PetscCall(MatSetValue(A, mycoo.i[k], mycoo.j[k], vals[k], ADD_VALUES)); 74 } 75 PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 76 PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 77 PetscCall(MatViewFromOptions(A, NULL, "-a_view")); 78 79 PetscCall(MatCreate(PETSC_COMM_WORLD, &B)); 80 PetscCall(MatSetSizes(B, PETSC_DECIDE, PETSC_DECIDE, M, N)); 81 PetscCall(MatSetFromOptions(B)); 82 PetscCall(MatSetOption(B, MAT_IGNORE_OFF_PROC_ENTRIES, flg)); 83 PetscCall(MatSetPreallocationCOO(B, mycoo.n, mycoo.i, mycoo.j)); 84 85 /* Test with ADD_VALUES on a zeroed matrix */ 86 PetscCall(MatSetValuesCOO(B, vals, ADD_VALUES)); 87 PetscCall(MatMultEqual(A, B, 10, &equal)); 88 if (!equal) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "MatSetValuesCOO() failed\n")); 89 PetscCall(MatViewFromOptions(B, NULL, "-b_view")); 90 91 /* Test with MatDuplicate on a zeroed matrix */ 92 PetscCall(MatDuplicate(B, MAT_DO_NOT_COPY_VALUES, &C)); 93 PetscCall(MatDestroy(&B)); 94 PetscCall(MatSetValuesCOO(C, vals, ADD_VALUES)); 95 PetscCall(MatMultEqual(A, C, 10, &equal)); 96 if (!equal) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "MatSetValuesCOO() on duplicated matrix failed\n")); 97 PetscCall(MatViewFromOptions(C, NULL, "-c_view")); 98 99 PetscCall(PetscFree(vals)); 100 if (freecoo) PetscCall(PetscFree2(mycoo.i, mycoo.j)); 101 PetscCall(MatDestroy(&A)); 102 PetscCall(MatDestroy(&C)); 103 104 PetscCall(PetscFinalize()); 105 return 0; 106 } 107 108 /*TEST 109 110 testset: 111 output_file: output/ex254_1.out 112 nsize: {{1 2 3}} 113 args: -ignore_remote {{0 1}} 114 filter: grep -v type | grep -v "Mat Object" 115 116 test: 117 suffix: kokkos 118 requires: kokkos_kernels 119 args: -mat_type aijkokkos 120 121 test: 122 suffix: cuda 123 requires: cuda 124 args: -mat_type aijcusparse 125 126 test: 127 suffix: hip 128 requires: hip 129 args: -mat_type aijhipsparse 130 131 test: 132 suffix: aij 133 args: -mat_type aij 134 135 test: 136 suffix: hypre 137 requires: hypre 138 args: -mat_type hypre 139 140 testset: 141 output_file: output/ex254_2.out 142 nsize: 1 143 args: -ncoos 3 144 filter: grep -v type | grep -v "Mat Object" 145 146 test: 147 suffix: 2_kokkos 148 requires: kokkos_kernels 149 args: -mat_type aijkokkos 150 151 test: 152 suffix: 2_cuda 153 requires: cuda 154 args: -mat_type aijcusparse 155 156 test: 157 suffix: 2_hip 158 requires: hip 159 args: -mat_type aijhipsparse 160 161 test: 162 suffix: 2_aij 163 args: -mat_type aij 164 165 test: 166 suffix: 2_hypre 167 requires: hypre 168 args: -mat_type hypre 169 170 TEST*/ 171