1 2 static char help[] = "Tests MatIncreaseOverlap(), MatCreateSubMatrices() for MatBAIJ format.\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 = 43, ov = 1, i, j, k, *rows, *cols, M, nd = 5, *idx, mm, nn, lsize; 10 PetscScalar *vals, rval; 11 IS *is1, *is2; 12 PetscRandom rdm; 13 Vec xx, s1, s2; 14 PetscReal s1norm, s2norm, rnorm, tol = PETSC_SQRT_MACHINE_EPSILON; 15 PetscBool flg; 16 17 PetscFunctionBeginUser; 18 PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 19 PetscCall(PetscOptionsGetInt(NULL, NULL, "-mat_block_size", &bs, NULL)); 20 PetscCall(PetscOptionsGetInt(NULL, NULL, "-mat_size", &m, NULL)); 21 PetscCall(PetscOptionsGetInt(NULL, NULL, "-ov", &ov, NULL)); 22 PetscCall(PetscOptionsGetInt(NULL, NULL, "-nd", &nd, NULL)); 23 M = m * bs; 24 25 PetscCall(MatCreateSeqBAIJ(PETSC_COMM_SELF, bs, M, M, 1, NULL, &A)); 26 PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE)); 27 PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, M, M, 15, NULL, &B)); 28 PetscCall(MatSetOption(B, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE)); 29 PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rdm)); 30 PetscCall(PetscRandomSetFromOptions(rdm)); 31 32 PetscCall(PetscMalloc1(bs, &rows)); 33 PetscCall(PetscMalloc1(bs, &cols)); 34 PetscCall(PetscMalloc1(bs * bs, &vals)); 35 PetscCall(PetscMalloc1(M, &idx)); 36 37 /* Now set blocks of values */ 38 for (i = 0; i < 20 * bs; i++) { 39 PetscCall(PetscRandomGetValue(rdm, &rval)); 40 cols[0] = bs * (int)(PetscRealPart(rval) * m); 41 PetscCall(PetscRandomGetValue(rdm, &rval)); 42 rows[0] = bs * (int)(PetscRealPart(rval) * m); 43 for (j = 1; j < bs; j++) { 44 rows[j] = rows[j - 1] + 1; 45 cols[j] = cols[j - 1] + 1; 46 } 47 48 for (j = 0; j < bs * bs; j++) { 49 PetscCall(PetscRandomGetValue(rdm, &rval)); 50 vals[j] = rval; 51 } 52 PetscCall(MatSetValues(A, bs, rows, bs, cols, vals, ADD_VALUES)); 53 PetscCall(MatSetValues(B, bs, rows, bs, cols, vals, ADD_VALUES)); 54 } 55 56 PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 57 PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 58 PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 59 PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 60 61 /* Test MatIncreaseOverlap() */ 62 PetscCall(PetscMalloc1(nd, &is1)); 63 PetscCall(PetscMalloc1(nd, &is2)); 64 65 for (i = 0; i < nd; i++) { 66 PetscCall(PetscRandomGetValue(rdm, &rval)); 67 lsize = (int)(PetscRealPart(rval) * m); 68 for (j = 0; j < lsize; j++) { 69 PetscCall(PetscRandomGetValue(rdm, &rval)); 70 idx[j * bs] = bs * (int)(PetscRealPart(rval) * m); 71 for (k = 1; k < bs; k++) idx[j * bs + k] = idx[j * bs] + k; 72 } 73 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, lsize * bs, idx, PETSC_COPY_VALUES, is1 + i)); 74 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, lsize * bs, idx, PETSC_COPY_VALUES, is2 + i)); 75 } 76 PetscCall(MatIncreaseOverlap(A, nd, is1, ov)); 77 PetscCall(MatIncreaseOverlap(B, nd, is2, ov)); 78 79 for (i = 0; i < nd; ++i) { 80 PetscCall(ISEqual(is1[i], is2[i], &flg)); 81 PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_PLIB, "i=%" PetscInt_FMT ", flg =%d", i, (int)flg); 82 } 83 84 for (i = 0; i < nd; ++i) { 85 PetscCall(ISSort(is1[i])); 86 PetscCall(ISSort(is2[i])); 87 } 88 89 PetscCall(MatCreateSubMatrices(A, nd, is1, is1, MAT_INITIAL_MATRIX, &submatA)); 90 PetscCall(MatCreateSubMatrices(B, nd, is2, is2, MAT_INITIAL_MATRIX, &submatB)); 91 92 /* Test MatMult() */ 93 for (i = 0; i < nd; i++) { 94 PetscCall(MatGetSize(submatA[i], &mm, &nn)); 95 PetscCall(VecCreateSeq(PETSC_COMM_SELF, mm, &xx)); 96 PetscCall(VecDuplicate(xx, &s1)); 97 PetscCall(VecDuplicate(xx, &s2)); 98 for (j = 0; j < 3; j++) { 99 PetscCall(VecSetRandom(xx, rdm)); 100 PetscCall(MatMult(submatA[i], xx, s1)); 101 PetscCall(MatMult(submatB[i], xx, s2)); 102 PetscCall(VecNorm(s1, NORM_2, &s1norm)); 103 PetscCall(VecNorm(s2, NORM_2, &s2norm)); 104 rnorm = s2norm - s1norm; 105 if (rnorm < -tol || rnorm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n", (double)s1norm, (double)s2norm)); 106 } 107 PetscCall(VecDestroy(&xx)); 108 PetscCall(VecDestroy(&s1)); 109 PetscCall(VecDestroy(&s2)); 110 } 111 /* Now test MatCreateSubmatrices with MAT_REUSE_MATRIX option */ 112 PetscCall(MatCreateSubMatrices(A, nd, is1, is1, MAT_REUSE_MATRIX, &submatA)); 113 PetscCall(MatCreateSubMatrices(B, nd, is2, is2, MAT_REUSE_MATRIX, &submatB)); 114 115 /* Test MatMult() */ 116 for (i = 0; i < nd; i++) { 117 PetscCall(MatGetSize(submatA[i], &mm, &nn)); 118 PetscCall(VecCreateSeq(PETSC_COMM_SELF, mm, &xx)); 119 PetscCall(VecDuplicate(xx, &s1)); 120 PetscCall(VecDuplicate(xx, &s2)); 121 for (j = 0; j < 3; j++) { 122 PetscCall(VecSetRandom(xx, rdm)); 123 PetscCall(MatMult(submatA[i], xx, s1)); 124 PetscCall(MatMult(submatB[i], xx, s2)); 125 PetscCall(VecNorm(s1, NORM_2, &s1norm)); 126 PetscCall(VecNorm(s2, NORM_2, &s2norm)); 127 rnorm = s2norm - s1norm; 128 if (rnorm < -tol || rnorm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n", (double)s1norm, (double)s2norm)); 129 } 130 PetscCall(VecDestroy(&xx)); 131 PetscCall(VecDestroy(&s1)); 132 PetscCall(VecDestroy(&s2)); 133 } 134 135 /* Free allocated memory */ 136 for (i = 0; i < nd; ++i) { 137 PetscCall(ISDestroy(&is1[i])); 138 PetscCall(ISDestroy(&is2[i])); 139 } 140 PetscCall(MatDestroySubMatrices(nd, &submatA)); 141 PetscCall(MatDestroySubMatrices(nd, &submatB)); 142 PetscCall(PetscFree(is1)); 143 PetscCall(PetscFree(is2)); 144 PetscCall(PetscFree(idx)); 145 PetscCall(PetscFree(rows)); 146 PetscCall(PetscFree(cols)); 147 PetscCall(PetscFree(vals)); 148 PetscCall(MatDestroy(&A)); 149 PetscCall(MatDestroy(&B)); 150 PetscCall(PetscRandomDestroy(&rdm)); 151 PetscCall(PetscFinalize()); 152 return 0; 153 } 154 155 /*TEST 156 157 test: 158 args: -mat_block_size {{1 2 5 7 8}} -ov {{1 3}} -mat_size {{11 13}} -nd {{7}} 159 160 TEST*/ 161