1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests various routines in MatSeqBAIJ format.\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, B, C, D, Fact; 9c4762a1bSJed Brown Vec xx, s1, s2, yy; 10c4762a1bSJed Brown PetscInt m = 45, rows[2], cols[2], bs = 1, i, row, col, *idx, M; 11c4762a1bSJed Brown PetscScalar rval, vals1[4], vals2[4]; 12c4762a1bSJed Brown PetscRandom rdm; 13c4762a1bSJed Brown IS is1, is2; 14c4762a1bSJed Brown PetscReal s1norm, s2norm, rnorm, tol = 1.e-4; 15c4762a1bSJed Brown PetscBool flg; 16c4762a1bSJed Brown MatFactorInfo info; 17c4762a1bSJed Brown 18327415f7SBarry Smith PetscFunctionBeginUser; 199566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 20c4762a1bSJed Brown /* Test MatSetValues() and MatGetValues() */ 219566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-mat_block_size", &bs, NULL)); 229566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-mat_size", &m, NULL)); 23c4762a1bSJed Brown M = m * bs; 249566063dSJacob Faibussowitsch PetscCall(MatCreateSeqBAIJ(PETSC_COMM_SELF, bs, M, M, 1, NULL, &A)); 259566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE)); 269566063dSJacob Faibussowitsch PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, M, M, 15, NULL, &B)); 279566063dSJacob Faibussowitsch PetscCall(MatSetOption(B, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE)); 289566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rdm)); 299566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rdm)); 309566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, M, &xx)); 319566063dSJacob Faibussowitsch PetscCall(VecDuplicate(xx, &s1)); 329566063dSJacob Faibussowitsch PetscCall(VecDuplicate(xx, &s2)); 339566063dSJacob Faibussowitsch PetscCall(VecDuplicate(xx, &yy)); 34c4762a1bSJed Brown 35c4762a1bSJed Brown /* For each row add at least 15 elements */ 36c4762a1bSJed Brown for (row = 0; row < M; row++) { 37c4762a1bSJed Brown for (i = 0; i < 25 * bs; i++) { 389566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rdm, &rval)); 39c4762a1bSJed Brown col = PetscMin(M - 1, (PetscInt)(PetscRealPart(rval) * M)); 409566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 1, &col, &rval, INSERT_VALUES)); 419566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 1, &row, 1, &col, &rval, INSERT_VALUES)); 42c4762a1bSJed Brown } 43c4762a1bSJed Brown } 44c4762a1bSJed Brown 45c4762a1bSJed Brown /* Now set blocks of values */ 46c4762a1bSJed Brown for (i = 0; i < 20 * bs; i++) { 479566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rdm, &rval)); 48c4762a1bSJed Brown cols[0] = PetscMin(M - 1, (PetscInt)(PetscRealPart(rval) * M)); 49c4762a1bSJed Brown vals1[0] = rval; 509566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rdm, &rval)); 51c4762a1bSJed Brown cols[1] = PetscMin(M - 1, (PetscInt)(PetscRealPart(rval) * M)); 52c4762a1bSJed Brown vals1[1] = rval; 539566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rdm, &rval)); 54c4762a1bSJed Brown rows[0] = PetscMin(M - 1, (PetscInt)(PetscRealPart(rval) * M)); 55c4762a1bSJed Brown vals1[2] = rval; 569566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rdm, &rval)); 57c4762a1bSJed Brown rows[1] = PetscMin(M - 1, (PetscInt)(PetscRealPart(rval) * M)); 58c4762a1bSJed Brown vals1[3] = rval; 599566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 2, rows, 2, cols, vals1, INSERT_VALUES)); 609566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 2, rows, 2, cols, vals1, INSERT_VALUES)); 61c4762a1bSJed Brown } 62c4762a1bSJed Brown 639566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 649566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 659566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 669566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 67c4762a1bSJed Brown 68c4762a1bSJed Brown /* Test MatNorm() */ 699566063dSJacob Faibussowitsch PetscCall(MatNorm(A, NORM_FROBENIUS, &s1norm)); 709566063dSJacob Faibussowitsch PetscCall(MatNorm(B, NORM_FROBENIUS, &s2norm)); 71c4762a1bSJed Brown rnorm = PetscAbsReal(s2norm - s1norm) / s2norm; 7248a46eb9SPierre Jolivet if (rnorm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error: MatNorm_FROBENIUS()- NormA=%16.14e NormB=%16.14e bs = %" PetscInt_FMT "\n", (double)s1norm, (double)s2norm, bs)); 739566063dSJacob Faibussowitsch PetscCall(MatNorm(A, NORM_INFINITY, &s1norm)); 749566063dSJacob Faibussowitsch PetscCall(MatNorm(B, NORM_INFINITY, &s2norm)); 75c4762a1bSJed Brown rnorm = PetscAbsReal(s2norm - s1norm) / s2norm; 7648a46eb9SPierre Jolivet if (rnorm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error: MatNorm_INFINITY()- NormA=%16.14e NormB=%16.14e bs = %" PetscInt_FMT "\n", (double)s1norm, (double)s2norm, bs)); 779566063dSJacob Faibussowitsch PetscCall(MatNorm(A, NORM_1, &s1norm)); 789566063dSJacob Faibussowitsch PetscCall(MatNorm(B, NORM_1, &s2norm)); 79c4762a1bSJed Brown rnorm = PetscAbsReal(s2norm - s1norm) / s2norm; 8048a46eb9SPierre Jolivet if (rnorm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error: MatNorm_NORM_1()- NormA=%16.14e NormB=%16.14e bs = %" PetscInt_FMT "\n", (double)s1norm, (double)s2norm, bs)); 81c4762a1bSJed Brown 82c4762a1bSJed Brown /* MatShift() */ 83c4762a1bSJed Brown rval = 10 * s1norm; 849566063dSJacob Faibussowitsch PetscCall(MatShift(A, rval)); 859566063dSJacob Faibussowitsch PetscCall(MatShift(B, rval)); 86c4762a1bSJed Brown 87c4762a1bSJed Brown /* Test MatTranspose() */ 889566063dSJacob Faibussowitsch PetscCall(MatTranspose(A, MAT_INPLACE_MATRIX, &A)); 899566063dSJacob Faibussowitsch PetscCall(MatTranspose(B, MAT_INPLACE_MATRIX, &B)); 90c4762a1bSJed Brown 91c4762a1bSJed Brown /* Now do MatGetValues() */ 92c4762a1bSJed Brown for (i = 0; i < 30; i++) { 939566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rdm, &rval)); 94c4762a1bSJed Brown cols[0] = PetscMin(M - 1, (PetscInt)(PetscRealPart(rval) * M)); 959566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rdm, &rval)); 96c4762a1bSJed Brown cols[1] = PetscMin(M - 1, (PetscInt)(PetscRealPart(rval) * M)); 979566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rdm, &rval)); 98c4762a1bSJed Brown rows[0] = PetscMin(M - 1, (PetscInt)(PetscRealPart(rval) * M)); 999566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rdm, &rval)); 100c4762a1bSJed Brown rows[1] = PetscMin(M - 1, (PetscInt)(PetscRealPart(rval) * M)); 1019566063dSJacob Faibussowitsch PetscCall(MatGetValues(A, 2, rows, 2, cols, vals1)); 1029566063dSJacob Faibussowitsch PetscCall(MatGetValues(B, 2, rows, 2, cols, vals2)); 1039566063dSJacob Faibussowitsch PetscCall(PetscArraycmp(vals1, vals2, 4, &flg)); 10448a46eb9SPierre Jolivet if (!flg) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error: MatGetValues bs = %" PetscInt_FMT "\n", bs)); 105c4762a1bSJed Brown } 106c4762a1bSJed Brown 107c4762a1bSJed Brown /* Test MatMult(), MatMultAdd() */ 108c4762a1bSJed Brown for (i = 0; i < 40; i++) { 1099566063dSJacob Faibussowitsch PetscCall(VecSetRandom(xx, rdm)); 1109566063dSJacob Faibussowitsch PetscCall(VecSet(s2, 0.0)); 1119566063dSJacob Faibussowitsch PetscCall(MatMult(A, xx, s1)); 1129566063dSJacob Faibussowitsch PetscCall(MatMultAdd(A, xx, s2, s2)); 1139566063dSJacob Faibussowitsch PetscCall(VecNorm(s1, NORM_2, &s1norm)); 1149566063dSJacob Faibussowitsch PetscCall(VecNorm(s2, NORM_2, &s2norm)); 115c4762a1bSJed Brown rnorm = s2norm - s1norm; 11648a46eb9SPierre Jolivet if (rnorm < -tol || rnorm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "MatMult not equal to MatMultAdd Norm1=%e Norm2=%e bs = %" PetscInt_FMT "\n", (double)s1norm, (double)s2norm, bs)); 117c4762a1bSJed Brown } 118c4762a1bSJed Brown 119c4762a1bSJed Brown /* Test MatMult() */ 1209566063dSJacob Faibussowitsch PetscCall(MatMultEqual(A, B, 10, &flg)); 12148a46eb9SPierre Jolivet if (!flg) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error: MatMult()\n")); 122c4762a1bSJed Brown 123c4762a1bSJed Brown /* Test MatMultAdd() */ 1249566063dSJacob Faibussowitsch PetscCall(MatMultAddEqual(A, B, 10, &flg)); 12548a46eb9SPierre Jolivet if (!flg) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error: MatMultAdd()\n")); 126c4762a1bSJed Brown 127c4762a1bSJed Brown /* Test MatMultTranspose() */ 1289566063dSJacob Faibussowitsch PetscCall(MatMultTransposeEqual(A, B, 10, &flg)); 12948a46eb9SPierre Jolivet if (!flg) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error: MatMultTranspose()\n")); 130c4762a1bSJed Brown 131c4762a1bSJed Brown /* Test MatMultTransposeAdd() */ 1329566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAddEqual(A, B, 10, &flg)); 13348a46eb9SPierre Jolivet if (!flg) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error: MatMultTransposeAdd()\n")); 134c4762a1bSJed Brown 135c4762a1bSJed Brown /* Test MatMatMult() */ 1369566063dSJacob Faibussowitsch PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, M, 40, NULL, &C)); 1379566063dSJacob Faibussowitsch PetscCall(MatSetRandom(C, rdm)); 1389566063dSJacob Faibussowitsch PetscCall(MatMatMult(A, C, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &D)); 1399566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual(A, C, D, 40, &flg)); 14048a46eb9SPierre Jolivet if (!flg) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error: MatMatMult()\n")); 1419566063dSJacob Faibussowitsch PetscCall(MatDestroy(&D)); 1429566063dSJacob Faibussowitsch PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, M, 40, NULL, &D)); 1439566063dSJacob Faibussowitsch PetscCall(MatMatMult(A, C, MAT_REUSE_MATRIX, PETSC_DEFAULT, &D)); 1449566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual(A, C, D, 40, &flg)); 14548a46eb9SPierre Jolivet if (!flg) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error: MatMatMult()\n")); 146c4762a1bSJed Brown 147c4762a1bSJed Brown /* Do LUFactor() on both the matrices */ 1489566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(M, &idx)); 149c4762a1bSJed Brown for (i = 0; i < M; i++) idx[i] = i; 1509566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, M, idx, PETSC_COPY_VALUES, &is1)); 1519566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, M, idx, PETSC_COPY_VALUES, &is2)); 1529566063dSJacob Faibussowitsch PetscCall(PetscFree(idx)); 1539566063dSJacob Faibussowitsch PetscCall(ISSetPermutation(is1)); 1549566063dSJacob Faibussowitsch PetscCall(ISSetPermutation(is2)); 155c4762a1bSJed Brown 1569566063dSJacob Faibussowitsch PetscCall(MatFactorInfoInitialize(&info)); 157c4762a1bSJed Brown 158c4762a1bSJed Brown info.fill = 2.0; 159c4762a1bSJed Brown info.dtcol = 0.0; 160c4762a1bSJed Brown info.zeropivot = 1.e-14; 161c4762a1bSJed Brown info.pivotinblocks = 1.0; 162c4762a1bSJed Brown 163c4762a1bSJed Brown if (bs < 4) { 1649566063dSJacob Faibussowitsch PetscCall(MatGetFactor(A, "petsc", MAT_FACTOR_LU, &Fact)); 1659566063dSJacob Faibussowitsch PetscCall(MatLUFactorSymbolic(Fact, A, is1, is2, &info)); 1669566063dSJacob Faibussowitsch PetscCall(MatLUFactorNumeric(Fact, A, &info)); 1679566063dSJacob Faibussowitsch PetscCall(VecSetRandom(yy, rdm)); 1689566063dSJacob Faibussowitsch PetscCall(MatForwardSolve(Fact, yy, xx)); 1699566063dSJacob Faibussowitsch PetscCall(MatBackwardSolve(Fact, xx, s1)); 1709566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Fact)); 1719566063dSJacob Faibussowitsch PetscCall(VecScale(s1, -1.0)); 1729566063dSJacob Faibussowitsch PetscCall(MatMultAdd(A, s1, yy, yy)); 1739566063dSJacob Faibussowitsch PetscCall(VecNorm(yy, NORM_2, &rnorm)); 17448a46eb9SPierre Jolivet if (rnorm < -tol || rnorm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error:MatForwardSolve/MatBackwardSolve - Norm1=%16.14e bs = %" PetscInt_FMT "\n", (double)rnorm, bs)); 175c4762a1bSJed Brown } 176c4762a1bSJed Brown 1779566063dSJacob Faibussowitsch PetscCall(MatLUFactor(B, is1, is2, &info)); 1789566063dSJacob Faibussowitsch PetscCall(MatLUFactor(A, is1, is2, &info)); 179c4762a1bSJed Brown 180c4762a1bSJed Brown /* Test MatSolveAdd() */ 181c4762a1bSJed Brown for (i = 0; i < 10; i++) { 1829566063dSJacob Faibussowitsch PetscCall(VecSetRandom(xx, rdm)); 1839566063dSJacob Faibussowitsch PetscCall(VecSetRandom(yy, rdm)); 1849566063dSJacob Faibussowitsch PetscCall(MatSolveAdd(B, xx, yy, s2)); 1859566063dSJacob Faibussowitsch PetscCall(MatSolveAdd(A, xx, yy, s1)); 1869566063dSJacob Faibussowitsch PetscCall(VecNorm(s1, NORM_2, &s1norm)); 1879566063dSJacob Faibussowitsch PetscCall(VecNorm(s2, NORM_2, &s2norm)); 188c4762a1bSJed Brown rnorm = s2norm - s1norm; 18948a46eb9SPierre Jolivet if (rnorm < -tol || rnorm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error:MatSolveAdd - Norm1=%16.14e Norm2=%16.14e bs = %" PetscInt_FMT "\n", (double)s1norm, (double)s2norm, bs)); 190c4762a1bSJed Brown } 191c4762a1bSJed Brown 192c4762a1bSJed Brown /* Test MatSolveAdd() when x = A'b +x */ 193c4762a1bSJed Brown for (i = 0; i < 10; i++) { 1949566063dSJacob Faibussowitsch PetscCall(VecSetRandom(xx, rdm)); 1959566063dSJacob Faibussowitsch PetscCall(VecSetRandom(s1, rdm)); 1969566063dSJacob Faibussowitsch PetscCall(VecCopy(s2, s1)); 1979566063dSJacob Faibussowitsch PetscCall(MatSolveAdd(B, xx, s2, s2)); 1989566063dSJacob Faibussowitsch PetscCall(MatSolveAdd(A, xx, s1, s1)); 1999566063dSJacob Faibussowitsch PetscCall(VecNorm(s1, NORM_2, &s1norm)); 2009566063dSJacob Faibussowitsch PetscCall(VecNorm(s2, NORM_2, &s2norm)); 201c4762a1bSJed Brown rnorm = s2norm - s1norm; 20248a46eb9SPierre Jolivet if (rnorm < -tol || rnorm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error:MatSolveAdd(same) - Norm1=%16.14e Norm2=%16.14e bs = %" PetscInt_FMT "\n", (double)s1norm, (double)s2norm, bs)); 203c4762a1bSJed Brown } 204c4762a1bSJed Brown 205c4762a1bSJed Brown /* Test MatSolve() */ 206c4762a1bSJed Brown for (i = 0; i < 10; i++) { 2079566063dSJacob Faibussowitsch PetscCall(VecSetRandom(xx, rdm)); 2089566063dSJacob Faibussowitsch PetscCall(MatSolve(B, xx, s2)); 2099566063dSJacob Faibussowitsch PetscCall(MatSolve(A, xx, s1)); 2109566063dSJacob Faibussowitsch PetscCall(VecNorm(s1, NORM_2, &s1norm)); 2119566063dSJacob Faibussowitsch PetscCall(VecNorm(s2, NORM_2, &s2norm)); 212c4762a1bSJed Brown rnorm = s2norm - s1norm; 21348a46eb9SPierre Jolivet if (rnorm < -tol || rnorm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error:MatSolve - Norm1=%16.14e Norm2=%16.14e bs = %" PetscInt_FMT "\n", (double)s1norm, (double)s2norm, bs)); 214c4762a1bSJed Brown } 215c4762a1bSJed Brown 216c4762a1bSJed Brown /* Test MatSolveTranspose() */ 217c4762a1bSJed Brown if (bs < 8) { 218c4762a1bSJed Brown for (i = 0; i < 10; i++) { 2199566063dSJacob Faibussowitsch PetscCall(VecSetRandom(xx, rdm)); 2209566063dSJacob Faibussowitsch PetscCall(MatSolveTranspose(B, xx, s2)); 2219566063dSJacob Faibussowitsch PetscCall(MatSolveTranspose(A, xx, s1)); 2229566063dSJacob Faibussowitsch PetscCall(VecNorm(s1, NORM_2, &s1norm)); 2239566063dSJacob Faibussowitsch PetscCall(VecNorm(s2, NORM_2, &s2norm)); 224c4762a1bSJed Brown rnorm = s2norm - s1norm; 22548a46eb9SPierre Jolivet if (rnorm < -tol || rnorm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error:MatSolveTranspose - Norm1=%16.14e Norm2=%16.14e bs = %" PetscInt_FMT "\n", (double)s1norm, (double)s2norm, bs)); 226c4762a1bSJed Brown } 227c4762a1bSJed Brown } 228c4762a1bSJed Brown 2299566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 2309566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 2319566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 2329566063dSJacob Faibussowitsch PetscCall(MatDestroy(&D)); 2339566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xx)); 2349566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s1)); 2359566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s2)); 2369566063dSJacob Faibussowitsch PetscCall(VecDestroy(&yy)); 2379566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is1)); 2389566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is2)); 2399566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rdm)); 2409566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 241b122ec5aSJacob Faibussowitsch return 0; 242c4762a1bSJed Brown } 243c4762a1bSJed Brown 244c4762a1bSJed Brown /*TEST 245c4762a1bSJed Brown 246c4762a1bSJed Brown test: 247c4762a1bSJed Brown args: -mat_block_size {{1 2 3 4 5 6 7 8}} 248c4762a1bSJed Brown 249c4762a1bSJed Brown TEST*/ 250