1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests MatIncreaseOverlap(), MatCreateSubMatrices() for parallel AIJ and BAIJ formats.\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown 69371c9d4SSatish Balay int main(int argc, char **args) { 7c4762a1bSJed Brown Mat A, B, *submatA, *submatB; 8c4762a1bSJed Brown PetscInt bs = 1, m = 11, ov = 1, i, j, k, *rows, *cols, nd = 5, *idx, rstart, rend, sz, mm, nn, M, N, Mbs; 9c4762a1bSJed Brown PetscMPIInt size, rank; 10c4762a1bSJed Brown PetscScalar *vals, rval; 11c4762a1bSJed Brown IS *is1, *is2; 12c4762a1bSJed Brown PetscRandom rdm; 13c4762a1bSJed Brown Vec xx, s1, s2; 14c4762a1bSJed Brown PetscReal s1norm, s2norm, rnorm, tol = 100 * PETSC_SMALL; 1582b5ce2aSStefano Zampini PetscBool flg, test_nd0 = PETSC_FALSE, emptynd; 16c4762a1bSJed Brown 17327415f7SBarry Smith PetscFunctionBeginUser; 189566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 199566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 209566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 21c4762a1bSJed Brown 229566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-mat_block_size", &bs, NULL)); 239566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-mat_size", &m, NULL)); 249566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-ov", &ov, NULL)); 259566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-nd", &nd, NULL)); 269566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_nd0", &test_nd0, NULL)); 27c4762a1bSJed Brown 28c4762a1bSJed Brown /* Create a AIJ matrix A */ 299566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 309566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, m * bs, m * bs, PETSC_DECIDE, PETSC_DECIDE)); 319566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATAIJ)); 329566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(A, PETSC_DEFAULT, NULL)); 339566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(A, PETSC_DEFAULT, NULL, PETSC_DEFAULT, NULL)); 349566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 359566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE)); 36c4762a1bSJed Brown 37c4762a1bSJed Brown /* Create a BAIJ matrix B */ 389566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &B)); 399566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B, m * bs, m * bs, PETSC_DECIDE, PETSC_DECIDE)); 409566063dSJacob Faibussowitsch PetscCall(MatSetType(B, MATBAIJ)); 419566063dSJacob Faibussowitsch PetscCall(MatSeqBAIJSetPreallocation(B, bs, PETSC_DEFAULT, NULL)); 429566063dSJacob Faibussowitsch PetscCall(MatMPIBAIJSetPreallocation(B, bs, PETSC_DEFAULT, NULL, PETSC_DEFAULT, NULL)); 439566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(B)); 449566063dSJacob Faibussowitsch PetscCall(MatSetOption(B, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE)); 45c4762a1bSJed Brown 469566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rdm)); 479566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rdm)); 48c4762a1bSJed Brown 499566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &rstart, &rend)); 509566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &M, &N)); 51c4762a1bSJed Brown Mbs = M / bs; 52c4762a1bSJed Brown 539566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bs, &rows)); 549566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bs, &cols)); 559566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bs * bs, &vals)); 569566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(M, &idx)); 57c4762a1bSJed Brown 58c4762a1bSJed Brown /* Now set blocks of values */ 59c4762a1bSJed Brown for (i = 0; i < 40 * bs; i++) { 609566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rdm, &rval)); 61c4762a1bSJed Brown cols[0] = bs * (int)(PetscRealPart(rval) * Mbs); 629566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rdm, &rval)); 63c4762a1bSJed Brown rows[0] = rstart + bs * (int)(PetscRealPart(rval) * m); 64c4762a1bSJed Brown for (j = 1; j < bs; j++) { 65c4762a1bSJed Brown rows[j] = rows[j - 1] + 1; 66c4762a1bSJed Brown cols[j] = cols[j - 1] + 1; 67c4762a1bSJed Brown } 68c4762a1bSJed Brown 69c4762a1bSJed Brown for (j = 0; j < bs * bs; j++) { 709566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rdm, &rval)); 71c4762a1bSJed Brown vals[j] = rval; 72c4762a1bSJed Brown } 739566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, bs, rows, bs, cols, vals, ADD_VALUES)); 749566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, bs, rows, bs, cols, vals, ADD_VALUES)); 75c4762a1bSJed Brown } 769566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 779566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 789566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 799566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 80c4762a1bSJed Brown 81c4762a1bSJed Brown /* Test MatIncreaseOverlap() */ 829566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nd, &is1)); 839566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nd, &is2)); 84c4762a1bSJed Brown 8582b5ce2aSStefano Zampini emptynd = PETSC_FALSE; 8682b5ce2aSStefano Zampini if (rank == 0 && test_nd0) emptynd = PETSC_TRUE; /* test case */ 87c4762a1bSJed Brown 88c4762a1bSJed Brown for (i = 0; i < nd; i++) { 899566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rdm, &rval)); 90c4762a1bSJed Brown sz = (int)(PetscRealPart(rval) * m); 91c4762a1bSJed Brown for (j = 0; j < sz; j++) { 929566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rdm, &rval)); 93c4762a1bSJed Brown idx[j * bs] = bs * (int)(PetscRealPart(rval) * Mbs); 94c4762a1bSJed Brown for (k = 1; k < bs; k++) idx[j * bs + k] = idx[j * bs] + k; 95c4762a1bSJed Brown } 9682b5ce2aSStefano Zampini PetscCall(ISCreateGeneral(PETSC_COMM_SELF, emptynd ? 0 : sz * bs, idx, PETSC_COPY_VALUES, is1 + i)); 9782b5ce2aSStefano Zampini PetscCall(ISCreateGeneral(PETSC_COMM_SELF, emptynd ? 0 : sz * bs, idx, PETSC_COPY_VALUES, is2 + i)); 98c4762a1bSJed Brown } 999566063dSJacob Faibussowitsch PetscCall(MatIncreaseOverlap(A, nd, is1, ov)); 1009566063dSJacob Faibussowitsch PetscCall(MatIncreaseOverlap(B, nd, is2, ov)); 101c4762a1bSJed Brown 102c4762a1bSJed Brown for (i = 0; i < nd; ++i) { 1039566063dSJacob Faibussowitsch PetscCall(ISEqual(is1[i], is2[i], &flg)); 104c4762a1bSJed Brown 105*48a46eb9SPierre Jolivet 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)); 106c4762a1bSJed Brown } 107c4762a1bSJed Brown 108c4762a1bSJed Brown for (i = 0; i < nd; ++i) { 1099566063dSJacob Faibussowitsch PetscCall(ISSort(is1[i])); 1109566063dSJacob Faibussowitsch PetscCall(ISSort(is2[i])); 111c4762a1bSJed Brown } 112c4762a1bSJed Brown 1139566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(B, nd, is2, is2, MAT_INITIAL_MATRIX, &submatB)); 1149566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(A, nd, is1, is1, MAT_INITIAL_MATRIX, &submatA)); 115c4762a1bSJed Brown 116c4762a1bSJed Brown /* Test MatMult() */ 117c4762a1bSJed Brown for (i = 0; i < nd; i++) { 1189566063dSJacob Faibussowitsch PetscCall(MatGetSize(submatA[i], &mm, &nn)); 1199566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, mm, &xx)); 1209566063dSJacob Faibussowitsch PetscCall(VecDuplicate(xx, &s1)); 1219566063dSJacob Faibussowitsch PetscCall(VecDuplicate(xx, &s2)); 122c4762a1bSJed Brown for (j = 0; j < 3; j++) { 1239566063dSJacob Faibussowitsch PetscCall(VecSetRandom(xx, rdm)); 1249566063dSJacob Faibussowitsch PetscCall(MatMult(submatA[i], xx, s1)); 1259566063dSJacob Faibussowitsch PetscCall(MatMult(submatB[i], xx, s2)); 1269566063dSJacob Faibussowitsch PetscCall(VecNorm(s1, NORM_2, &s1norm)); 1279566063dSJacob Faibussowitsch PetscCall(VecNorm(s2, NORM_2, &s2norm)); 128c4762a1bSJed Brown rnorm = s2norm - s1norm; 129*48a46eb9SPierre Jolivet 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)); 130c4762a1bSJed Brown } 1319566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xx)); 1329566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s1)); 1339566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s2)); 134c4762a1bSJed Brown } 135c4762a1bSJed Brown 136c4762a1bSJed Brown /* Now test MatCreateSubmatrices with MAT_REUSE_MATRIX option */ 1379566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(A, nd, is1, is1, MAT_REUSE_MATRIX, &submatA)); 1389566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(B, nd, is2, is2, MAT_REUSE_MATRIX, &submatB)); 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, rdm)); 1489566063dSJacob Faibussowitsch PetscCall(MatMult(submatA[i], xx, s1)); 1499566063dSJacob Faibussowitsch PetscCall(MatMult(submatB[i], xx, s2)); 1509566063dSJacob Faibussowitsch PetscCall(VecNorm(s1, NORM_2, &s1norm)); 1519566063dSJacob Faibussowitsch PetscCall(VecNorm(s2, NORM_2, &s2norm)); 152c4762a1bSJed Brown rnorm = s2norm - s1norm; 153*48a46eb9SPierre Jolivet 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)); 154c4762a1bSJed Brown } 1559566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xx)); 1569566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s1)); 1579566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s2)); 158c4762a1bSJed Brown } 159c4762a1bSJed Brown 160c4762a1bSJed Brown /* Free allocated memory */ 161c4762a1bSJed Brown for (i = 0; i < nd; ++i) { 1629566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is1[i])); 1639566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is2[i])); 164c4762a1bSJed Brown } 1659566063dSJacob Faibussowitsch PetscCall(MatDestroySubMatrices(nd, &submatA)); 1669566063dSJacob Faibussowitsch PetscCall(MatDestroySubMatrices(nd, &submatB)); 167c4762a1bSJed Brown 1689566063dSJacob Faibussowitsch PetscCall(PetscFree(is1)); 1699566063dSJacob Faibussowitsch PetscCall(PetscFree(is2)); 1709566063dSJacob Faibussowitsch PetscCall(PetscFree(idx)); 1719566063dSJacob Faibussowitsch PetscCall(PetscFree(rows)); 1729566063dSJacob Faibussowitsch PetscCall(PetscFree(cols)); 1739566063dSJacob Faibussowitsch PetscCall(PetscFree(vals)); 1749566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1759566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 1769566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rdm)); 1779566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 178b122ec5aSJacob Faibussowitsch return 0; 179c4762a1bSJed Brown } 180c4762a1bSJed Brown 181c4762a1bSJed Brown /*TEST 182c4762a1bSJed Brown 183c4762a1bSJed Brown test: 184c4762a1bSJed Brown nsize: {{1 3}} 18582b5ce2aSStefano Zampini args: -mat_block_size {{1 3 4 6 8}} -ov {{1 3}} -mat_size {{11 13}} -nd 7 186c4762a1bSJed Brown output_file: output/ex54.out 187c4762a1bSJed Brown 188c4762a1bSJed Brown test: 189c4762a1bSJed Brown suffix: 2 190c4762a1bSJed Brown args: -nd 2 -test_nd0 191c4762a1bSJed Brown output_file: output/ex54.out 192c4762a1bSJed Brown 193c4762a1bSJed Brown test: 194c4762a1bSJed Brown suffix: 3 195c4762a1bSJed Brown nsize: 3 196c4762a1bSJed Brown args: -nd 2 -test_nd0 197c4762a1bSJed Brown output_file: output/ex54.out 198c4762a1bSJed Brown 199c4762a1bSJed Brown TEST*/ 200