1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests the various sequential routines in MatSBAIJ format. Same as ex74.c except diagonal entries of the matrices are zeros.\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown 69371c9d4SSatish Balay int main(int argc, char **args) { 7c4762a1bSJed Brown Vec x, y, b, s1, s2; 8c4762a1bSJed Brown Mat A; /* linear system matrix */ 9c4762a1bSJed Brown Mat sA; /* symmetric part of the matrices */ 10c4762a1bSJed Brown PetscInt n, mbs = 16, bs = 1, nz = 3, prob = 2, i, j, col[3], row, Ii, J, n1; 11c4762a1bSJed Brown const PetscInt *ip_ptr; 12c4762a1bSJed Brown PetscScalar neg_one = -1.0, value[3], alpha = 0.1; 13c4762a1bSJed Brown PetscMPIInt size; 14c4762a1bSJed Brown IS ip, isrow, iscol; 15c4762a1bSJed Brown PetscRandom rdm; 16c4762a1bSJed Brown PetscBool reorder = PETSC_FALSE; 17c4762a1bSJed Brown MatInfo minfo1, minfo2; 18c4762a1bSJed Brown PetscReal norm1, norm2, tol = 10 * PETSC_SMALL; 19c4762a1bSJed Brown 20327415f7SBarry Smith PetscFunctionBeginUser; 219566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 229566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 23be096a46SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!"); 249566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-bs", &bs, NULL)); 259566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-mbs", &mbs, NULL)); 26c4762a1bSJed Brown 27c4762a1bSJed Brown n = mbs * bs; 289566063dSJacob Faibussowitsch PetscCall(MatCreateSeqBAIJ(PETSC_COMM_WORLD, bs, n, n, nz, NULL, &A)); 299566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE)); 309566063dSJacob Faibussowitsch PetscCall(MatCreateSeqSBAIJ(PETSC_COMM_WORLD, bs, n, n, nz, NULL, &sA)); 319566063dSJacob Faibussowitsch PetscCall(MatSetOption(sA, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE)); 32c4762a1bSJed Brown 33c4762a1bSJed Brown /* Test MatGetOwnershipRange() */ 349566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &Ii, &J)); 359566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(sA, &i, &j)); 36*48a46eb9SPierre Jolivet if (i - Ii || j - J) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error: MatGetOwnershipRange() in MatSBAIJ format\n")); 37c4762a1bSJed Brown 38c4762a1bSJed Brown /* Assemble matrix */ 39c4762a1bSJed Brown if (bs == 1) { 409566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-test_problem", &prob, NULL)); 41c4762a1bSJed Brown if (prob == 1) { /* tridiagonal matrix */ 429371c9d4SSatish Balay value[0] = -1.0; 439371c9d4SSatish Balay value[1] = 2.0; 449371c9d4SSatish Balay value[2] = -1.0; 45c4762a1bSJed Brown for (i = 1; i < n - 1; i++) { 469371c9d4SSatish Balay col[0] = i - 1; 479371c9d4SSatish Balay col[1] = i; 489371c9d4SSatish Balay col[2] = i + 1; 499566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 3, col, value, INSERT_VALUES)); 509566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA, 1, &i, 3, col, value, INSERT_VALUES)); 51c4762a1bSJed Brown } 529371c9d4SSatish Balay i = n - 1; 539371c9d4SSatish Balay col[0] = 0; 549371c9d4SSatish Balay col[1] = n - 2; 559371c9d4SSatish Balay col[2] = n - 1; 56c4762a1bSJed Brown 579371c9d4SSatish Balay value[0] = 0.1; 589371c9d4SSatish Balay value[1] = -1; 599371c9d4SSatish Balay value[2] = 2; 609566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 3, col, value, INSERT_VALUES)); 619566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA, 1, &i, 3, col, value, INSERT_VALUES)); 62c4762a1bSJed Brown 639371c9d4SSatish Balay i = 0; 649371c9d4SSatish Balay col[0] = 0; 659371c9d4SSatish Balay col[1] = 1; 669371c9d4SSatish Balay col[2] = n - 1; 67c4762a1bSJed Brown 689371c9d4SSatish Balay value[0] = 2.0; 699371c9d4SSatish Balay value[1] = -1.0; 709371c9d4SSatish Balay value[2] = 0.1; 719566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 3, col, value, INSERT_VALUES)); 729566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA, 1, &i, 3, col, value, INSERT_VALUES)); 73c4762a1bSJed Brown } else if (prob == 2) { /* matrix for the five point stencil */ 74c4762a1bSJed Brown n1 = (PetscInt)(PetscSqrtReal((PetscReal)n) + 0.001); 75bc30f867SBarry Smith PetscCheck(n1 * n1 == n, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "sqrt(n) must be a positive integer!"); 76c4762a1bSJed Brown for (i = 0; i < n1; i++) { 77c4762a1bSJed Brown for (j = 0; j < n1; j++) { 78c4762a1bSJed Brown Ii = j + n1 * i; 79c4762a1bSJed Brown if (i > 0) { 80c4762a1bSJed Brown J = Ii - n1; 819566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &neg_one, INSERT_VALUES)); 829566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA, 1, &Ii, 1, &J, &neg_one, INSERT_VALUES)); 83c4762a1bSJed Brown } 84c4762a1bSJed Brown if (i < n1 - 1) { 85c4762a1bSJed Brown J = Ii + n1; 869566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &neg_one, INSERT_VALUES)); 879566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA, 1, &Ii, 1, &J, &neg_one, INSERT_VALUES)); 88c4762a1bSJed Brown } 89c4762a1bSJed Brown if (j > 0) { 90c4762a1bSJed Brown J = Ii - 1; 919566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &neg_one, INSERT_VALUES)); 929566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA, 1, &Ii, 1, &J, &neg_one, INSERT_VALUES)); 93c4762a1bSJed Brown } 94c4762a1bSJed Brown if (j < n1 - 1) { 95c4762a1bSJed Brown J = Ii + 1; 969566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &neg_one, INSERT_VALUES)); 979566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA, 1, &Ii, 1, &J, &neg_one, INSERT_VALUES)); 98c4762a1bSJed Brown } 99c4762a1bSJed Brown } 100c4762a1bSJed Brown } 101c4762a1bSJed Brown } 102c4762a1bSJed Brown } else { /* bs > 1 */ 103c4762a1bSJed Brown #if defined(DIAGB) 104c4762a1bSJed Brown for (block = 0; block < n / bs; block++) { 105c4762a1bSJed Brown /* diagonal blocks */ 1069371c9d4SSatish Balay value[0] = -1.0; 1079371c9d4SSatish Balay value[1] = 4.0; 1089371c9d4SSatish Balay value[2] = -1.0; 109c4762a1bSJed Brown for (i = 1 + block * bs; i < bs - 1 + block * bs; i++) { 1109371c9d4SSatish Balay col[0] = i - 1; 1119371c9d4SSatish Balay col[1] = i; 1129371c9d4SSatish Balay col[2] = i + 1; 1139566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 3, col, value, INSERT_VALUES)); 1149566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA, 1, &i, 3, col, value, INSERT_VALUES)); 115c4762a1bSJed Brown } 1169371c9d4SSatish Balay i = bs - 1 + block * bs; 1179371c9d4SSatish Balay col[0] = bs - 2 + block * bs; 1189371c9d4SSatish Balay col[1] = bs - 1 + block * bs; 119c4762a1bSJed Brown 1209371c9d4SSatish Balay value[0] = -1.0; 1219371c9d4SSatish Balay value[1] = 4.0; 1229566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 2, col, value, INSERT_VALUES)); 1239566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA, 1, &i, 2, col, value, INSERT_VALUES)); 124c4762a1bSJed Brown 1259371c9d4SSatish Balay i = 0 + block * bs; 1269371c9d4SSatish Balay col[0] = 0 + block * bs; 1279371c9d4SSatish Balay col[1] = 1 + block * bs; 128c4762a1bSJed Brown 1299371c9d4SSatish Balay value[0] = 4.0; 1309371c9d4SSatish Balay value[1] = -1.0; 1319566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 2, col, value, INSERT_VALUES)); 1329566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA, 1, &i, 2, col, value, INSERT_VALUES)); 133c4762a1bSJed Brown } 134c4762a1bSJed Brown #endif 135c4762a1bSJed Brown /* off-diagonal blocks */ 136c4762a1bSJed Brown value[0] = -1.0; 137c4762a1bSJed Brown for (i = 0; i < (n / bs - 1) * bs; i++) { 138c4762a1bSJed Brown col[0] = i + bs; 1399566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 1, col, value, INSERT_VALUES)); 1409566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA, 1, &i, 1, col, value, INSERT_VALUES)); 1419371c9d4SSatish Balay col[0] = i; 1429371c9d4SSatish Balay row = i + bs; 1439566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 1, col, value, INSERT_VALUES)); 1449566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA, 1, &row, 1, col, value, INSERT_VALUES)); 145c4762a1bSJed Brown } 146c4762a1bSJed Brown } 1479566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 1489566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 149c4762a1bSJed Brown 1509566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(sA, MAT_FINAL_ASSEMBLY)); 1519566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(sA, MAT_FINAL_ASSEMBLY)); 152c4762a1bSJed Brown 153c4762a1bSJed Brown /* Test MatNorm() */ 1549566063dSJacob Faibussowitsch PetscCall(MatNorm(A, NORM_FROBENIUS, &norm1)); 1559566063dSJacob Faibussowitsch PetscCall(MatNorm(sA, NORM_FROBENIUS, &norm2)); 156c4762a1bSJed Brown norm1 -= norm2; 157*48a46eb9SPierre Jolivet if (norm1 < -tol || norm1 > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error: MatNorm(), fnorm1-fnorm2=%16.14e\n", (double)norm1)); 1589566063dSJacob Faibussowitsch PetscCall(MatNorm(A, NORM_INFINITY, &norm1)); 1599566063dSJacob Faibussowitsch PetscCall(MatNorm(sA, NORM_INFINITY, &norm2)); 160c4762a1bSJed Brown norm1 -= norm2; 161*48a46eb9SPierre Jolivet if (norm1 < -tol || norm1 > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error: MatNorm(), inf_norm1-inf_norm2=%16.14e\n", (double)norm1)); 162c4762a1bSJed Brown 163c4762a1bSJed Brown /* Test MatGetInfo(), MatGetSize(), MatGetBlockSize() */ 1649566063dSJacob Faibussowitsch PetscCall(MatGetInfo(A, MAT_LOCAL, &minfo1)); 1659566063dSJacob Faibussowitsch PetscCall(MatGetInfo(sA, MAT_LOCAL, &minfo2)); 166c4762a1bSJed Brown i = (int)(minfo1.nz_used - minfo2.nz_used); 167c4762a1bSJed Brown j = (int)(minfo1.nz_allocated - minfo2.nz_allocated); 168*48a46eb9SPierre Jolivet if (i < 0 || j < 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error: MatGetInfo()\n")); 169c4762a1bSJed Brown 1709566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &Ii, &J)); 1719566063dSJacob Faibussowitsch PetscCall(MatGetSize(sA, &i, &j)); 172*48a46eb9SPierre Jolivet if (i - Ii || j - J) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error: MatGetSize()\n")); 173c4762a1bSJed Brown 1749566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(A, &Ii)); 1759566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(sA, &i)); 176*48a46eb9SPierre Jolivet if (i - Ii) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error: MatGetBlockSize()\n")); 177c4762a1bSJed Brown 178c4762a1bSJed Brown /* Test MatDiagonalScale(), MatGetDiagonal(), MatScale() */ 1799566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rdm)); 1809566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rdm)); 1819566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, n, &x)); 1829566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &s1)); 1839566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &s2)); 1849566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &y)); 1859566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &b)); 186c4762a1bSJed Brown 1879566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x, rdm)); 188c4762a1bSJed Brown 1899566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(A, x, x)); 1909566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(sA, x, x)); 191c4762a1bSJed Brown 1929566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(A, s1)); 1939566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(sA, s2)); 1949566063dSJacob Faibussowitsch PetscCall(VecNorm(s1, NORM_1, &norm1)); 1959566063dSJacob Faibussowitsch PetscCall(VecNorm(s2, NORM_1, &norm2)); 196c4762a1bSJed Brown norm1 -= norm2; 197*48a46eb9SPierre Jolivet if (norm1 < -tol || norm1 > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error:MatGetDiagonal() \n")); 198c4762a1bSJed Brown 1999566063dSJacob Faibussowitsch PetscCall(MatScale(A, alpha)); 2009566063dSJacob Faibussowitsch PetscCall(MatScale(sA, alpha)); 201c4762a1bSJed Brown 202c4762a1bSJed Brown /* Test MatMult(), MatMultAdd() */ 203c4762a1bSJed Brown for (i = 0; i < 40; i++) { 2049566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x, rdm)); 2059566063dSJacob Faibussowitsch PetscCall(MatMult(A, x, s1)); 2069566063dSJacob Faibussowitsch PetscCall(MatMult(sA, x, s2)); 2079566063dSJacob Faibussowitsch PetscCall(VecNorm(s1, NORM_1, &norm1)); 2089566063dSJacob Faibussowitsch PetscCall(VecNorm(s2, NORM_1, &norm2)); 209c4762a1bSJed Brown norm1 -= norm2; 210*48a46eb9SPierre Jolivet if (norm1 < -tol || norm1 > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error: MatMult(), MatDiagonalScale() or MatScale()\n")); 211c4762a1bSJed Brown } 212c4762a1bSJed Brown 213c4762a1bSJed Brown for (i = 0; i < 40; i++) { 2149566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x, rdm)); 2159566063dSJacob Faibussowitsch PetscCall(VecSetRandom(y, rdm)); 2169566063dSJacob Faibussowitsch PetscCall(MatMultAdd(A, x, y, s1)); 2179566063dSJacob Faibussowitsch PetscCall(MatMultAdd(sA, x, y, s2)); 2189566063dSJacob Faibussowitsch PetscCall(VecNorm(s1, NORM_1, &norm1)); 2199566063dSJacob Faibussowitsch PetscCall(VecNorm(s2, NORM_1, &norm2)); 220c4762a1bSJed Brown norm1 -= norm2; 221*48a46eb9SPierre Jolivet if (norm1 < -tol || norm1 > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error:MatMultAdd(), MatDiagonalScale() or MatScale() \n")); 222c4762a1bSJed Brown } 223c4762a1bSJed Brown 224c4762a1bSJed Brown /* Test MatReordering() */ 2259566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(A, MATORDERINGNATURAL, &isrow, &iscol)); 226c4762a1bSJed Brown ip = isrow; 227c4762a1bSJed Brown 228c4762a1bSJed Brown if (reorder) { 229c4762a1bSJed Brown IS nip; 230c4762a1bSJed Brown PetscInt *nip_ptr; 2319566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(mbs, &nip_ptr)); 2329566063dSJacob Faibussowitsch PetscCall(ISGetIndices(ip, &ip_ptr)); 2339566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(nip_ptr, ip_ptr, mbs)); 2349371c9d4SSatish Balay i = nip_ptr[1]; 2359371c9d4SSatish Balay nip_ptr[1] = nip_ptr[mbs - 2]; 2369371c9d4SSatish Balay nip_ptr[mbs - 2] = i; 2379371c9d4SSatish Balay i = nip_ptr[0]; 2389371c9d4SSatish Balay nip_ptr[0] = nip_ptr[mbs - 1]; 2399371c9d4SSatish Balay nip_ptr[mbs - 1] = i; 2409566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(ip, &ip_ptr)); 2419566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, mbs, nip_ptr, PETSC_COPY_VALUES, &nip)); 2429566063dSJacob Faibussowitsch PetscCall(PetscFree(nip_ptr)); 243c4762a1bSJed Brown 2449566063dSJacob Faibussowitsch PetscCall(MatReorderingSeqSBAIJ(sA, ip)); 2459566063dSJacob Faibussowitsch PetscCall(ISDestroy(&nip)); 246c4762a1bSJed Brown } 247c4762a1bSJed Brown 2489566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iscol)); 2499566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isrow)); 2509566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 2519566063dSJacob Faibussowitsch PetscCall(MatDestroy(&sA)); 2529566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 2539566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 2549566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s1)); 2559566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s2)); 2569566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b)); 2579566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rdm)); 258c4762a1bSJed Brown 2599566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 260b122ec5aSJacob Faibussowitsch return 0; 261c4762a1bSJed Brown } 262c4762a1bSJed Brown 263c4762a1bSJed Brown /*TEST 264c4762a1bSJed Brown 265c4762a1bSJed Brown test: 266c4762a1bSJed Brown args: -bs {{1 2 3 4 5 6 7 8}} 267c4762a1bSJed Brown 268c4762a1bSJed Brown TEST*/ 269