1 2 static char help[] = "Tests MatIncreaseOverlap(), MatCreateSubMatrices() for parallel AIJ and BAIJ formats.\n"; 3 4 #include <petscmat.h> 5 6 int main(int argc, char **args) 7 { 8 Mat A, B, *submatA, *submatB; 9 PetscInt bs = 1, m = 11, ov = 1, i, j, k, *rows, *cols, nd = 5, *idx, rstart, rend, sz, mm, nn, M, N, Mbs; 10 PetscMPIInt size, rank; 11 PetscScalar *vals, rval; 12 IS *is1, *is2; 13 PetscRandom rdm; 14 Vec xx, s1, s2; 15 PetscReal s1norm, s2norm, rnorm, tol = 100 * PETSC_SMALL; 16 PetscBool flg, test_nd0 = PETSC_FALSE, emptynd; 17 18 PetscFunctionBeginUser; 19 PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 20 PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 21 PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 22 23 PetscCall(PetscOptionsGetInt(NULL, NULL, "-mat_block_size", &bs, NULL)); 24 PetscCall(PetscOptionsGetInt(NULL, NULL, "-mat_size", &m, NULL)); 25 PetscCall(PetscOptionsGetInt(NULL, NULL, "-ov", &ov, NULL)); 26 PetscCall(PetscOptionsGetInt(NULL, NULL, "-nd", &nd, NULL)); 27 PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_nd0", &test_nd0, NULL)); 28 29 /* Create a AIJ matrix A */ 30 PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 31 PetscCall(MatSetSizes(A, m * bs, m * bs, PETSC_DECIDE, PETSC_DECIDE)); 32 PetscCall(MatSetType(A, MATAIJ)); 33 PetscCall(MatSeqAIJSetPreallocation(A, PETSC_DEFAULT, NULL)); 34 PetscCall(MatMPIAIJSetPreallocation(A, PETSC_DEFAULT, NULL, PETSC_DEFAULT, NULL)); 35 PetscCall(MatSetFromOptions(A)); 36 PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE)); 37 38 /* Create a BAIJ matrix B */ 39 PetscCall(MatCreate(PETSC_COMM_WORLD, &B)); 40 PetscCall(MatSetSizes(B, m * bs, m * bs, PETSC_DECIDE, PETSC_DECIDE)); 41 PetscCall(MatSetType(B, MATBAIJ)); 42 PetscCall(MatSeqBAIJSetPreallocation(B, bs, PETSC_DEFAULT, NULL)); 43 PetscCall(MatMPIBAIJSetPreallocation(B, bs, PETSC_DEFAULT, NULL, PETSC_DEFAULT, NULL)); 44 PetscCall(MatSetFromOptions(B)); 45 PetscCall(MatSetOption(B, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE)); 46 47 PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rdm)); 48 PetscCall(PetscRandomSetFromOptions(rdm)); 49 50 PetscCall(MatGetOwnershipRange(A, &rstart, &rend)); 51 PetscCall(MatGetSize(A, &M, &N)); 52 Mbs = M / bs; 53 54 PetscCall(PetscMalloc1(bs, &rows)); 55 PetscCall(PetscMalloc1(bs, &cols)); 56 PetscCall(PetscMalloc1(bs * bs, &vals)); 57 PetscCall(PetscMalloc1(M, &idx)); 58 59 /* Now set blocks of values */ 60 for (i = 0; i < 40 * bs; i++) { 61 PetscCall(PetscRandomGetValue(rdm, &rval)); 62 cols[0] = bs * (int)(PetscRealPart(rval) * Mbs); 63 PetscCall(PetscRandomGetValue(rdm, &rval)); 64 rows[0] = rstart + bs * (int)(PetscRealPart(rval) * m); 65 for (j = 1; j < bs; j++) { 66 rows[j] = rows[j - 1] + 1; 67 cols[j] = cols[j - 1] + 1; 68 } 69 70 for (j = 0; j < bs * bs; j++) { 71 PetscCall(PetscRandomGetValue(rdm, &rval)); 72 vals[j] = rval; 73 } 74 PetscCall(MatSetValues(A, bs, rows, bs, cols, vals, ADD_VALUES)); 75 PetscCall(MatSetValues(B, bs, rows, bs, cols, vals, ADD_VALUES)); 76 } 77 PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 78 PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 79 PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 80 PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 81 82 /* Test MatIncreaseOverlap() */ 83 PetscCall(PetscMalloc1(nd, &is1)); 84 PetscCall(PetscMalloc1(nd, &is2)); 85 86 emptynd = PETSC_FALSE; 87 if (rank == 0 && test_nd0) emptynd = PETSC_TRUE; /* test case */ 88 89 for (i = 0; i < nd; i++) { 90 PetscCall(PetscRandomGetValue(rdm, &rval)); 91 sz = (int)(PetscRealPart(rval) * m); 92 for (j = 0; j < sz; j++) { 93 PetscCall(PetscRandomGetValue(rdm, &rval)); 94 idx[j * bs] = bs * (int)(PetscRealPart(rval) * Mbs); 95 for (k = 1; k < bs; k++) idx[j * bs + k] = idx[j * bs] + k; 96 } 97 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, emptynd ? 0 : sz * bs, idx, PETSC_COPY_VALUES, is1 + i)); 98 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, emptynd ? 0 : sz * bs, idx, PETSC_COPY_VALUES, is2 + i)); 99 } 100 PetscCall(MatIncreaseOverlap(A, nd, is1, ov)); 101 PetscCall(MatIncreaseOverlap(B, nd, is2, ov)); 102 103 for (i = 0; i < nd; ++i) { 104 PetscCall(ISEqual(is1[i], is2[i], &flg)); 105 106 if (!flg) PetscCall(PetscPrintf(PETSC_COMM_SELF, "i=%" PetscInt_FMT ", flg=%d :bs=%" PetscInt_FMT " m=%" PetscInt_FMT " ov=%" PetscInt_FMT " nd=%" PetscInt_FMT " np=%d\n", i, flg, bs, m, ov, nd, size)); 107 } 108 109 for (i = 0; i < nd; ++i) { 110 PetscCall(ISSort(is1[i])); 111 PetscCall(ISSort(is2[i])); 112 } 113 114 PetscCall(MatCreateSubMatrices(B, nd, is2, is2, MAT_INITIAL_MATRIX, &submatB)); 115 PetscCall(MatCreateSubMatrices(A, nd, is1, is1, MAT_INITIAL_MATRIX, &submatA)); 116 117 /* Test MatMult() */ 118 for (i = 0; i < nd; i++) { 119 PetscCall(MatGetSize(submatA[i], &mm, &nn)); 120 PetscCall(VecCreateSeq(PETSC_COMM_SELF, mm, &xx)); 121 PetscCall(VecDuplicate(xx, &s1)); 122 PetscCall(VecDuplicate(xx, &s2)); 123 for (j = 0; j < 3; j++) { 124 PetscCall(VecSetRandom(xx, rdm)); 125 PetscCall(MatMult(submatA[i], xx, s1)); 126 PetscCall(MatMult(submatB[i], xx, s2)); 127 PetscCall(VecNorm(s1, NORM_2, &s1norm)); 128 PetscCall(VecNorm(s2, NORM_2, &s2norm)); 129 rnorm = s2norm - s1norm; 130 if (rnorm < -tol || rnorm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%d]Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n", rank, (double)s1norm, (double)s2norm)); 131 } 132 PetscCall(VecDestroy(&xx)); 133 PetscCall(VecDestroy(&s1)); 134 PetscCall(VecDestroy(&s2)); 135 } 136 137 /* Now test MatCreateSubmatrices with MAT_REUSE_MATRIX option */ 138 PetscCall(MatCreateSubMatrices(A, nd, is1, is1, MAT_REUSE_MATRIX, &submatA)); 139 PetscCall(MatCreateSubMatrices(B, nd, is2, is2, MAT_REUSE_MATRIX, &submatB)); 140 141 /* Test MatMult() */ 142 for (i = 0; i < nd; i++) { 143 PetscCall(MatGetSize(submatA[i], &mm, &nn)); 144 PetscCall(VecCreateSeq(PETSC_COMM_SELF, mm, &xx)); 145 PetscCall(VecDuplicate(xx, &s1)); 146 PetscCall(VecDuplicate(xx, &s2)); 147 for (j = 0; j < 3; j++) { 148 PetscCall(VecSetRandom(xx, rdm)); 149 PetscCall(MatMult(submatA[i], xx, s1)); 150 PetscCall(MatMult(submatB[i], xx, s2)); 151 PetscCall(VecNorm(s1, NORM_2, &s1norm)); 152 PetscCall(VecNorm(s2, NORM_2, &s2norm)); 153 rnorm = s2norm - s1norm; 154 if (rnorm < -tol || rnorm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%d]Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n", rank, (double)s1norm, (double)s2norm)); 155 } 156 PetscCall(VecDestroy(&xx)); 157 PetscCall(VecDestroy(&s1)); 158 PetscCall(VecDestroy(&s2)); 159 } 160 161 /* Free allocated memory */ 162 for (i = 0; i < nd; ++i) { 163 PetscCall(ISDestroy(&is1[i])); 164 PetscCall(ISDestroy(&is2[i])); 165 } 166 PetscCall(MatDestroySubMatrices(nd, &submatA)); 167 PetscCall(MatDestroySubMatrices(nd, &submatB)); 168 169 PetscCall(PetscFree(is1)); 170 PetscCall(PetscFree(is2)); 171 PetscCall(PetscFree(idx)); 172 PetscCall(PetscFree(rows)); 173 PetscCall(PetscFree(cols)); 174 PetscCall(PetscFree(vals)); 175 PetscCall(MatDestroy(&A)); 176 PetscCall(MatDestroy(&B)); 177 PetscCall(PetscRandomDestroy(&rdm)); 178 PetscCall(PetscFinalize()); 179 return 0; 180 } 181 182 /*TEST 183 184 test: 185 nsize: {{1 3}} 186 args: -mat_block_size {{1 3 4 6 8}} -ov {{1 3}} -mat_size {{11 13}} -nd 7 187 output_file: output/ex54.out 188 189 test: 190 suffix: 2 191 args: -nd 2 -test_nd0 192 output_file: output/ex54.out 193 194 test: 195 suffix: 3 196 nsize: 3 197 args: -nd 2 -test_nd0 198 output_file: output/ex54.out 199 200 TEST*/ 201