1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests MatIncreaseOverlap(), MatCreateSubMatrices() for sequential MatSBAIJ format. Derived from ex51.c\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown 6*d71ae5a4SJacob Faibussowitsch int main(int argc, char **args) 7*d71ae5a4SJacob Faibussowitsch { 8c4762a1bSJed Brown Mat A, Atrans, sA, *submatA, *submatsA; 9c4762a1bSJed Brown PetscInt bs = 1, m = 43, ov = 1, i, j, k, *rows, *cols, M, nd = 5, *idx, mm, nn; 10c4762a1bSJed Brown PetscMPIInt size; 11c4762a1bSJed Brown PetscScalar *vals, rval, one = 1.0; 12c4762a1bSJed Brown IS *is1, *is2; 13c4762a1bSJed Brown PetscRandom rand; 14c4762a1bSJed Brown Vec xx, s1, s2; 15c4762a1bSJed Brown PetscReal s1norm, s2norm, rnorm, tol = 10 * PETSC_SMALL; 16c4762a1bSJed Brown PetscBool flg; 17c4762a1bSJed Brown 18327415f7SBarry Smith PetscFunctionBeginUser; 199566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 209566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-mat_block_size", &bs, NULL)); 219566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-mat_size", &m, NULL)); 229566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-ov", &ov, NULL)); 239566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-nd", &nd, NULL)); 24c4762a1bSJed Brown 25c4762a1bSJed Brown /* create a SeqBAIJ matrix A */ 26c4762a1bSJed Brown M = m * bs; 279566063dSJacob Faibussowitsch PetscCall(MatCreateSeqBAIJ(PETSC_COMM_SELF, bs, M, M, 1, NULL, &A)); 289566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE)); 299566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rand)); 309566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rand)); 31c4762a1bSJed Brown 329566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bs, &rows)); 339566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bs, &cols)); 349566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bs * bs, &vals)); 359566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(M, &idx)); 36c4762a1bSJed Brown 37c4762a1bSJed Brown /* Now set blocks of random values */ 38c4762a1bSJed Brown /* first, set diagonal blocks as zero */ 39c4762a1bSJed Brown for (j = 0; j < bs * bs; j++) vals[j] = 0.0; 40c4762a1bSJed Brown for (i = 0; i < m; i++) { 419371c9d4SSatish Balay cols[0] = i * bs; 429371c9d4SSatish Balay rows[0] = i * bs; 43c4762a1bSJed Brown for (j = 1; j < bs; j++) { 44c4762a1bSJed Brown rows[j] = rows[j - 1] + 1; 45c4762a1bSJed Brown cols[j] = cols[j - 1] + 1; 46c4762a1bSJed Brown } 479566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, bs, rows, bs, cols, vals, ADD_VALUES)); 48c4762a1bSJed Brown } 49c4762a1bSJed Brown /* second, add random blocks */ 50c4762a1bSJed Brown for (i = 0; i < 20 * bs; i++) { 519566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rand, &rval)); 52c4762a1bSJed Brown cols[0] = bs * (int)(PetscRealPart(rval) * m); 539566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rand, &rval)); 54c4762a1bSJed Brown rows[0] = bs * (int)(PetscRealPart(rval) * m); 55c4762a1bSJed Brown for (j = 1; j < bs; j++) { 56c4762a1bSJed Brown rows[j] = rows[j - 1] + 1; 57c4762a1bSJed Brown cols[j] = cols[j - 1] + 1; 58c4762a1bSJed Brown } 59c4762a1bSJed Brown 60c4762a1bSJed Brown for (j = 0; j < bs * bs; j++) { 619566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rand, &rval)); 62c4762a1bSJed Brown vals[j] = rval; 63c4762a1bSJed Brown } 649566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, bs, rows, bs, cols, vals, ADD_VALUES)); 65c4762a1bSJed Brown } 66c4762a1bSJed Brown 679566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 689566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 69c4762a1bSJed Brown 70c4762a1bSJed Brown /* make A a symmetric matrix: A <- A^T + A */ 719566063dSJacob Faibussowitsch PetscCall(MatTranspose(A, MAT_INITIAL_MATRIX, &Atrans)); 729566063dSJacob Faibussowitsch PetscCall(MatAXPY(A, one, Atrans, DIFFERENT_NONZERO_PATTERN)); 739566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Atrans)); 749566063dSJacob Faibussowitsch PetscCall(MatTranspose(A, MAT_INITIAL_MATRIX, &Atrans)); 759566063dSJacob Faibussowitsch PetscCall(MatEqual(A, Atrans, &flg)); 7628b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_PLIB, "A+A^T is non-symmetric"); 779566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Atrans)); 78c4762a1bSJed Brown 79c4762a1bSJed Brown /* create a SeqSBAIJ matrix sA (= A) */ 809566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_SYMMETRIC, PETSC_TRUE)); 819566063dSJacob Faibussowitsch PetscCall(MatConvert(A, MATSEQSBAIJ, MAT_INITIAL_MATRIX, &sA)); 82c4762a1bSJed Brown 83c4762a1bSJed Brown /* Test sA==A through MatMult() */ 84c4762a1bSJed Brown for (i = 0; i < nd; i++) { 859566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &mm, &nn)); 869566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, mm, &xx)); 879566063dSJacob Faibussowitsch PetscCall(VecDuplicate(xx, &s1)); 889566063dSJacob Faibussowitsch PetscCall(VecDuplicate(xx, &s2)); 89c4762a1bSJed Brown for (j = 0; j < 3; j++) { 909566063dSJacob Faibussowitsch PetscCall(VecSetRandom(xx, rand)); 919566063dSJacob Faibussowitsch PetscCall(MatMult(A, xx, s1)); 929566063dSJacob Faibussowitsch PetscCall(MatMult(sA, xx, s2)); 939566063dSJacob Faibussowitsch PetscCall(VecNorm(s1, NORM_2, &s1norm)); 949566063dSJacob Faibussowitsch PetscCall(VecNorm(s2, NORM_2, &s2norm)); 95c4762a1bSJed Brown rnorm = s2norm - s1norm; 9648a46eb9SPierre Jolivet if (rnorm < -tol || rnorm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n", (double)s1norm, (double)s2norm)); 97c4762a1bSJed Brown } 989566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xx)); 999566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s1)); 1009566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s2)); 101c4762a1bSJed Brown } 102c4762a1bSJed Brown 103c4762a1bSJed Brown /* Test MatIncreaseOverlap() */ 1049566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nd, &is1)); 1059566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nd, &is2)); 106c4762a1bSJed Brown 107c4762a1bSJed Brown for (i = 0; i < nd; i++) { 1089566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rand, &rval)); 109c4762a1bSJed Brown size = (int)(PetscRealPart(rval) * m); 110c4762a1bSJed Brown for (j = 0; j < size; j++) { 1119566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rand, &rval)); 112c4762a1bSJed Brown idx[j * bs] = bs * (int)(PetscRealPart(rval) * m); 113c4762a1bSJed Brown for (k = 1; k < bs; k++) idx[j * bs + k] = idx[j * bs] + k; 114c4762a1bSJed Brown } 1159566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, size * bs, idx, PETSC_COPY_VALUES, is1 + i)); 1169566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, size * bs, idx, PETSC_COPY_VALUES, is2 + i)); 117c4762a1bSJed Brown } 118c4762a1bSJed Brown /* for debugging */ 119c4762a1bSJed Brown /* 1209566063dSJacob Faibussowitsch PetscCall(MatView(A,PETSC_VIEWER_STDOUT_SELF)); 1219566063dSJacob Faibussowitsch PetscCall(MatView(sA,PETSC_VIEWER_STDOUT_SELF)); 122c4762a1bSJed Brown */ 123c4762a1bSJed Brown 1249566063dSJacob Faibussowitsch PetscCall(MatIncreaseOverlap(A, nd, is1, ov)); 1259566063dSJacob Faibussowitsch PetscCall(MatIncreaseOverlap(sA, nd, is2, ov)); 126c4762a1bSJed Brown 127c4762a1bSJed Brown for (i = 0; i < nd; ++i) { 1289566063dSJacob Faibussowitsch PetscCall(ISSort(is1[i])); 1299566063dSJacob Faibussowitsch PetscCall(ISSort(is2[i])); 130c4762a1bSJed Brown } 131c4762a1bSJed Brown 132c4762a1bSJed Brown for (i = 0; i < nd; ++i) { 1339566063dSJacob Faibussowitsch PetscCall(ISEqual(is1[i], is2[i], &flg)); 13428b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_PLIB, "i=%" PetscInt_FMT ", is1 != is2", i); 135c4762a1bSJed Brown } 136c4762a1bSJed Brown 1379566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(A, nd, is1, is1, MAT_INITIAL_MATRIX, &submatA)); 1389566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(sA, nd, is2, is2, MAT_INITIAL_MATRIX, &submatsA)); 139c4762a1bSJed Brown 140c4762a1bSJed Brown /* Test MatMult() */ 141c4762a1bSJed Brown for (i = 0; i < nd; i++) { 1429566063dSJacob Faibussowitsch PetscCall(MatGetSize(submatA[i], &mm, &nn)); 1439566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, mm, &xx)); 1449566063dSJacob Faibussowitsch PetscCall(VecDuplicate(xx, &s1)); 1459566063dSJacob Faibussowitsch PetscCall(VecDuplicate(xx, &s2)); 146c4762a1bSJed Brown for (j = 0; j < 3; j++) { 1479566063dSJacob Faibussowitsch PetscCall(VecSetRandom(xx, rand)); 1489566063dSJacob Faibussowitsch PetscCall(MatMult(submatA[i], xx, s1)); 1499566063dSJacob Faibussowitsch PetscCall(MatMult(submatsA[i], xx, s2)); 1509566063dSJacob Faibussowitsch PetscCall(VecNorm(s1, NORM_2, &s1norm)); 1519566063dSJacob Faibussowitsch PetscCall(VecNorm(s2, NORM_2, &s2norm)); 152c4762a1bSJed Brown rnorm = s2norm - s1norm; 15348a46eb9SPierre Jolivet if (rnorm < -tol || rnorm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n", (double)s1norm, (double)s2norm)); 154c4762a1bSJed Brown } 1559566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xx)); 1569566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s1)); 1579566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s2)); 158c4762a1bSJed Brown } 159c4762a1bSJed Brown 160c4762a1bSJed Brown /* Now test MatCreateSubmatrices with MAT_REUSE_MATRIX option */ 1619566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(A, nd, is1, is1, MAT_REUSE_MATRIX, &submatA)); 1629566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(sA, nd, is2, is2, MAT_REUSE_MATRIX, &submatsA)); 163c4762a1bSJed Brown 164c4762a1bSJed Brown /* Test MatMult() */ 165c4762a1bSJed Brown for (i = 0; i < nd; i++) { 1669566063dSJacob Faibussowitsch PetscCall(MatGetSize(submatA[i], &mm, &nn)); 1679566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, mm, &xx)); 1689566063dSJacob Faibussowitsch PetscCall(VecDuplicate(xx, &s1)); 1699566063dSJacob Faibussowitsch PetscCall(VecDuplicate(xx, &s2)); 170c4762a1bSJed Brown for (j = 0; j < 3; j++) { 1719566063dSJacob Faibussowitsch PetscCall(VecSetRandom(xx, rand)); 1729566063dSJacob Faibussowitsch PetscCall(MatMult(submatA[i], xx, s1)); 1739566063dSJacob Faibussowitsch PetscCall(MatMult(submatsA[i], xx, s2)); 1749566063dSJacob Faibussowitsch PetscCall(VecNorm(s1, NORM_2, &s1norm)); 1759566063dSJacob Faibussowitsch PetscCall(VecNorm(s2, NORM_2, &s2norm)); 176c4762a1bSJed Brown rnorm = s2norm - s1norm; 17748a46eb9SPierre Jolivet if (rnorm < -tol || rnorm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n", (double)s1norm, (double)s2norm)); 178c4762a1bSJed Brown } 1799566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xx)); 1809566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s1)); 1819566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s2)); 182c4762a1bSJed Brown } 183c4762a1bSJed Brown 184c4762a1bSJed Brown /* Free allocated memory */ 185c4762a1bSJed Brown for (i = 0; i < nd; ++i) { 1869566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is1[i])); 1879566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is2[i])); 188c4762a1bSJed Brown } 1899566063dSJacob Faibussowitsch PetscCall(MatDestroySubMatrices(nd, &submatA)); 1909566063dSJacob Faibussowitsch PetscCall(MatDestroySubMatrices(nd, &submatsA)); 191c4762a1bSJed Brown 1929566063dSJacob Faibussowitsch PetscCall(PetscFree(is1)); 1939566063dSJacob Faibussowitsch PetscCall(PetscFree(is2)); 1949566063dSJacob Faibussowitsch PetscCall(PetscFree(idx)); 1959566063dSJacob Faibussowitsch PetscCall(PetscFree(rows)); 1969566063dSJacob Faibussowitsch PetscCall(PetscFree(cols)); 1979566063dSJacob Faibussowitsch PetscCall(PetscFree(vals)); 1989566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1999566063dSJacob Faibussowitsch PetscCall(MatDestroy(&sA)); 2009566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rand)); 2019566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 202b122ec5aSJacob Faibussowitsch return 0; 203c4762a1bSJed Brown } 204c4762a1bSJed Brown 205c4762a1bSJed Brown /*TEST 206c4762a1bSJed Brown 207c4762a1bSJed Brown test: 208c4762a1bSJed Brown args: -ov 2 209c4762a1bSJed Brown 210c4762a1bSJed Brown TEST*/ 211