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 6*d71ae5a4SJacob Faibussowitsch int main(int argc, char **args) 7*d71ae5a4SJacob Faibussowitsch { 8c4762a1bSJed Brown Vec x, y, b, s1, s2; 9c4762a1bSJed Brown Mat A; /* linear system matrix */ 10c4762a1bSJed Brown Mat sA; /* symmetric part of the matrices */ 11c4762a1bSJed Brown PetscInt n, mbs = 16, bs = 1, nz = 3, prob = 2, i, j, col[3], row, Ii, J, n1; 12c4762a1bSJed Brown const PetscInt *ip_ptr; 13c4762a1bSJed Brown PetscScalar neg_one = -1.0, value[3], alpha = 0.1; 14c4762a1bSJed Brown PetscMPIInt size; 15c4762a1bSJed Brown IS ip, isrow, iscol; 16c4762a1bSJed Brown PetscRandom rdm; 17c4762a1bSJed Brown PetscBool reorder = PETSC_FALSE; 18c4762a1bSJed Brown MatInfo minfo1, minfo2; 19c4762a1bSJed Brown PetscReal norm1, norm2, tol = 10 * PETSC_SMALL; 20c4762a1bSJed Brown 21327415f7SBarry Smith PetscFunctionBeginUser; 229566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 239566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 24be096a46SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!"); 259566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-bs", &bs, NULL)); 269566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-mbs", &mbs, NULL)); 27c4762a1bSJed Brown 28c4762a1bSJed Brown n = mbs * bs; 299566063dSJacob Faibussowitsch PetscCall(MatCreateSeqBAIJ(PETSC_COMM_WORLD, bs, n, n, nz, NULL, &A)); 309566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE)); 319566063dSJacob Faibussowitsch PetscCall(MatCreateSeqSBAIJ(PETSC_COMM_WORLD, bs, n, n, nz, NULL, &sA)); 329566063dSJacob Faibussowitsch PetscCall(MatSetOption(sA, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE)); 33c4762a1bSJed Brown 34c4762a1bSJed Brown /* Test MatGetOwnershipRange() */ 359566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &Ii, &J)); 369566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(sA, &i, &j)); 3748a46eb9SPierre Jolivet if (i - Ii || j - J) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error: MatGetOwnershipRange() in MatSBAIJ format\n")); 38c4762a1bSJed Brown 39c4762a1bSJed Brown /* Assemble matrix */ 40c4762a1bSJed Brown if (bs == 1) { 419566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-test_problem", &prob, NULL)); 42c4762a1bSJed Brown if (prob == 1) { /* tridiagonal matrix */ 439371c9d4SSatish Balay value[0] = -1.0; 449371c9d4SSatish Balay value[1] = 2.0; 459371c9d4SSatish Balay value[2] = -1.0; 46c4762a1bSJed Brown for (i = 1; i < n - 1; i++) { 479371c9d4SSatish Balay col[0] = i - 1; 489371c9d4SSatish Balay col[1] = i; 499371c9d4SSatish Balay col[2] = i + 1; 509566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 3, col, value, INSERT_VALUES)); 519566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA, 1, &i, 3, col, value, INSERT_VALUES)); 52c4762a1bSJed Brown } 539371c9d4SSatish Balay i = n - 1; 549371c9d4SSatish Balay col[0] = 0; 559371c9d4SSatish Balay col[1] = n - 2; 569371c9d4SSatish Balay col[2] = n - 1; 57c4762a1bSJed Brown 589371c9d4SSatish Balay value[0] = 0.1; 599371c9d4SSatish Balay value[1] = -1; 609371c9d4SSatish Balay value[2] = 2; 619566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 3, col, value, INSERT_VALUES)); 629566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA, 1, &i, 3, col, value, INSERT_VALUES)); 63c4762a1bSJed Brown 649371c9d4SSatish Balay i = 0; 659371c9d4SSatish Balay col[0] = 0; 669371c9d4SSatish Balay col[1] = 1; 679371c9d4SSatish Balay col[2] = n - 1; 68c4762a1bSJed Brown 699371c9d4SSatish Balay value[0] = 2.0; 709371c9d4SSatish Balay value[1] = -1.0; 719371c9d4SSatish Balay value[2] = 0.1; 729566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 3, col, value, INSERT_VALUES)); 739566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA, 1, &i, 3, col, value, INSERT_VALUES)); 74c4762a1bSJed Brown } else if (prob == 2) { /* matrix for the five point stencil */ 75c4762a1bSJed Brown n1 = (PetscInt)(PetscSqrtReal((PetscReal)n) + 0.001); 76bc30f867SBarry Smith PetscCheck(n1 * n1 == n, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "sqrt(n) must be a positive integer!"); 77c4762a1bSJed Brown for (i = 0; i < n1; i++) { 78c4762a1bSJed Brown for (j = 0; j < n1; j++) { 79c4762a1bSJed Brown Ii = j + n1 * i; 80c4762a1bSJed Brown if (i > 0) { 81c4762a1bSJed Brown J = Ii - n1; 829566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &neg_one, INSERT_VALUES)); 839566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA, 1, &Ii, 1, &J, &neg_one, INSERT_VALUES)); 84c4762a1bSJed Brown } 85c4762a1bSJed Brown if (i < n1 - 1) { 86c4762a1bSJed Brown J = Ii + n1; 879566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &neg_one, INSERT_VALUES)); 889566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA, 1, &Ii, 1, &J, &neg_one, INSERT_VALUES)); 89c4762a1bSJed Brown } 90c4762a1bSJed Brown if (j > 0) { 91c4762a1bSJed Brown J = Ii - 1; 929566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &neg_one, INSERT_VALUES)); 939566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA, 1, &Ii, 1, &J, &neg_one, INSERT_VALUES)); 94c4762a1bSJed Brown } 95c4762a1bSJed Brown if (j < n1 - 1) { 96c4762a1bSJed Brown J = Ii + 1; 979566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &neg_one, INSERT_VALUES)); 989566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA, 1, &Ii, 1, &J, &neg_one, INSERT_VALUES)); 99c4762a1bSJed Brown } 100c4762a1bSJed Brown } 101c4762a1bSJed Brown } 102c4762a1bSJed Brown } 103c4762a1bSJed Brown } else { /* bs > 1 */ 104c4762a1bSJed Brown #if defined(DIAGB) 105c4762a1bSJed Brown for (block = 0; block < n / bs; block++) { 106c4762a1bSJed Brown /* diagonal blocks */ 1079371c9d4SSatish Balay value[0] = -1.0; 1089371c9d4SSatish Balay value[1] = 4.0; 1099371c9d4SSatish Balay value[2] = -1.0; 110c4762a1bSJed Brown for (i = 1 + block * bs; i < bs - 1 + block * bs; i++) { 1119371c9d4SSatish Balay col[0] = i - 1; 1129371c9d4SSatish Balay col[1] = i; 1139371c9d4SSatish Balay col[2] = i + 1; 1149566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 3, col, value, INSERT_VALUES)); 1159566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA, 1, &i, 3, col, value, INSERT_VALUES)); 116c4762a1bSJed Brown } 1179371c9d4SSatish Balay i = bs - 1 + block * bs; 1189371c9d4SSatish Balay col[0] = bs - 2 + block * bs; 1199371c9d4SSatish Balay col[1] = bs - 1 + block * bs; 120c4762a1bSJed Brown 1219371c9d4SSatish Balay value[0] = -1.0; 1229371c9d4SSatish Balay value[1] = 4.0; 1239566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 2, col, value, INSERT_VALUES)); 1249566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA, 1, &i, 2, col, value, INSERT_VALUES)); 125c4762a1bSJed Brown 1269371c9d4SSatish Balay i = 0 + block * bs; 1279371c9d4SSatish Balay col[0] = 0 + block * bs; 1289371c9d4SSatish Balay col[1] = 1 + block * bs; 129c4762a1bSJed Brown 1309371c9d4SSatish Balay value[0] = 4.0; 1319371c9d4SSatish Balay value[1] = -1.0; 1329566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 2, col, value, INSERT_VALUES)); 1339566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA, 1, &i, 2, col, value, INSERT_VALUES)); 134c4762a1bSJed Brown } 135c4762a1bSJed Brown #endif 136c4762a1bSJed Brown /* off-diagonal blocks */ 137c4762a1bSJed Brown value[0] = -1.0; 138c4762a1bSJed Brown for (i = 0; i < (n / bs - 1) * bs; i++) { 139c4762a1bSJed Brown col[0] = i + bs; 1409566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 1, col, value, INSERT_VALUES)); 1419566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA, 1, &i, 1, col, value, INSERT_VALUES)); 1429371c9d4SSatish Balay col[0] = i; 1439371c9d4SSatish Balay row = i + bs; 1449566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 1, col, value, INSERT_VALUES)); 1459566063dSJacob Faibussowitsch PetscCall(MatSetValues(sA, 1, &row, 1, col, value, INSERT_VALUES)); 146c4762a1bSJed Brown } 147c4762a1bSJed Brown } 1489566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 1499566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 150c4762a1bSJed Brown 1519566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(sA, MAT_FINAL_ASSEMBLY)); 1529566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(sA, MAT_FINAL_ASSEMBLY)); 153c4762a1bSJed Brown 154c4762a1bSJed Brown /* Test MatNorm() */ 1559566063dSJacob Faibussowitsch PetscCall(MatNorm(A, NORM_FROBENIUS, &norm1)); 1569566063dSJacob Faibussowitsch PetscCall(MatNorm(sA, NORM_FROBENIUS, &norm2)); 157c4762a1bSJed Brown norm1 -= norm2; 15848a46eb9SPierre Jolivet if (norm1 < -tol || norm1 > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error: MatNorm(), fnorm1-fnorm2=%16.14e\n", (double)norm1)); 1599566063dSJacob Faibussowitsch PetscCall(MatNorm(A, NORM_INFINITY, &norm1)); 1609566063dSJacob Faibussowitsch PetscCall(MatNorm(sA, NORM_INFINITY, &norm2)); 161c4762a1bSJed Brown norm1 -= norm2; 16248a46eb9SPierre Jolivet if (norm1 < -tol || norm1 > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error: MatNorm(), inf_norm1-inf_norm2=%16.14e\n", (double)norm1)); 163c4762a1bSJed Brown 164c4762a1bSJed Brown /* Test MatGetInfo(), MatGetSize(), MatGetBlockSize() */ 1659566063dSJacob Faibussowitsch PetscCall(MatGetInfo(A, MAT_LOCAL, &minfo1)); 1669566063dSJacob Faibussowitsch PetscCall(MatGetInfo(sA, MAT_LOCAL, &minfo2)); 167c4762a1bSJed Brown i = (int)(minfo1.nz_used - minfo2.nz_used); 168c4762a1bSJed Brown j = (int)(minfo1.nz_allocated - minfo2.nz_allocated); 16948a46eb9SPierre Jolivet if (i < 0 || j < 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error: MatGetInfo()\n")); 170c4762a1bSJed Brown 1719566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &Ii, &J)); 1729566063dSJacob Faibussowitsch PetscCall(MatGetSize(sA, &i, &j)); 17348a46eb9SPierre Jolivet if (i - Ii || j - J) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error: MatGetSize()\n")); 174c4762a1bSJed Brown 1759566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(A, &Ii)); 1769566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(sA, &i)); 17748a46eb9SPierre Jolivet if (i - Ii) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error: MatGetBlockSize()\n")); 178c4762a1bSJed Brown 179c4762a1bSJed Brown /* Test MatDiagonalScale(), MatGetDiagonal(), MatScale() */ 1809566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rdm)); 1819566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rdm)); 1829566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, n, &x)); 1839566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &s1)); 1849566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &s2)); 1859566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &y)); 1869566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &b)); 187c4762a1bSJed Brown 1889566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x, rdm)); 189c4762a1bSJed Brown 1909566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(A, x, x)); 1919566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(sA, x, x)); 192c4762a1bSJed Brown 1939566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(A, s1)); 1949566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(sA, s2)); 1959566063dSJacob Faibussowitsch PetscCall(VecNorm(s1, NORM_1, &norm1)); 1969566063dSJacob Faibussowitsch PetscCall(VecNorm(s2, NORM_1, &norm2)); 197c4762a1bSJed Brown norm1 -= norm2; 19848a46eb9SPierre Jolivet if (norm1 < -tol || norm1 > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error:MatGetDiagonal() \n")); 199c4762a1bSJed Brown 2009566063dSJacob Faibussowitsch PetscCall(MatScale(A, alpha)); 2019566063dSJacob Faibussowitsch PetscCall(MatScale(sA, alpha)); 202c4762a1bSJed Brown 203c4762a1bSJed Brown /* Test MatMult(), MatMultAdd() */ 204c4762a1bSJed Brown for (i = 0; i < 40; i++) { 2059566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x, rdm)); 2069566063dSJacob Faibussowitsch PetscCall(MatMult(A, x, s1)); 2079566063dSJacob Faibussowitsch PetscCall(MatMult(sA, x, s2)); 2089566063dSJacob Faibussowitsch PetscCall(VecNorm(s1, NORM_1, &norm1)); 2099566063dSJacob Faibussowitsch PetscCall(VecNorm(s2, NORM_1, &norm2)); 210c4762a1bSJed Brown norm1 -= norm2; 21148a46eb9SPierre Jolivet if (norm1 < -tol || norm1 > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error: MatMult(), MatDiagonalScale() or MatScale()\n")); 212c4762a1bSJed Brown } 213c4762a1bSJed Brown 214c4762a1bSJed Brown for (i = 0; i < 40; i++) { 2159566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x, rdm)); 2169566063dSJacob Faibussowitsch PetscCall(VecSetRandom(y, rdm)); 2179566063dSJacob Faibussowitsch PetscCall(MatMultAdd(A, x, y, s1)); 2189566063dSJacob Faibussowitsch PetscCall(MatMultAdd(sA, x, y, s2)); 2199566063dSJacob Faibussowitsch PetscCall(VecNorm(s1, NORM_1, &norm1)); 2209566063dSJacob Faibussowitsch PetscCall(VecNorm(s2, NORM_1, &norm2)); 221c4762a1bSJed Brown norm1 -= norm2; 22248a46eb9SPierre Jolivet if (norm1 < -tol || norm1 > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error:MatMultAdd(), MatDiagonalScale() or MatScale() \n")); 223c4762a1bSJed Brown } 224c4762a1bSJed Brown 225c4762a1bSJed Brown /* Test MatReordering() */ 2269566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(A, MATORDERINGNATURAL, &isrow, &iscol)); 227c4762a1bSJed Brown ip = isrow; 228c4762a1bSJed Brown 229c4762a1bSJed Brown if (reorder) { 230c4762a1bSJed Brown IS nip; 231c4762a1bSJed Brown PetscInt *nip_ptr; 2329566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(mbs, &nip_ptr)); 2339566063dSJacob Faibussowitsch PetscCall(ISGetIndices(ip, &ip_ptr)); 2349566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(nip_ptr, ip_ptr, mbs)); 2359371c9d4SSatish Balay i = nip_ptr[1]; 2369371c9d4SSatish Balay nip_ptr[1] = nip_ptr[mbs - 2]; 2379371c9d4SSatish Balay nip_ptr[mbs - 2] = i; 2389371c9d4SSatish Balay i = nip_ptr[0]; 2399371c9d4SSatish Balay nip_ptr[0] = nip_ptr[mbs - 1]; 2409371c9d4SSatish Balay nip_ptr[mbs - 1] = i; 2419566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(ip, &ip_ptr)); 2429566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, mbs, nip_ptr, PETSC_COPY_VALUES, &nip)); 2439566063dSJacob Faibussowitsch PetscCall(PetscFree(nip_ptr)); 244c4762a1bSJed Brown 2459566063dSJacob Faibussowitsch PetscCall(MatReorderingSeqSBAIJ(sA, ip)); 2469566063dSJacob Faibussowitsch PetscCall(ISDestroy(&nip)); 247c4762a1bSJed Brown } 248c4762a1bSJed Brown 2499566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iscol)); 2509566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isrow)); 2519566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 2529566063dSJacob Faibussowitsch PetscCall(MatDestroy(&sA)); 2539566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 2549566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 2559566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s1)); 2569566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s2)); 2579566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b)); 2589566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rdm)); 259c4762a1bSJed Brown 2609566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 261b122ec5aSJacob Faibussowitsch return 0; 262c4762a1bSJed Brown } 263c4762a1bSJed Brown 264c4762a1bSJed Brown /*TEST 265c4762a1bSJed Brown 266c4762a1bSJed Brown test: 267c4762a1bSJed Brown args: -bs {{1 2 3 4 5 6 7 8}} 268c4762a1bSJed Brown 269c4762a1bSJed Brown TEST*/ 270