18a9c020eSBarry Smith static char help[] = "Tests MatInvertVariableBlockEnvelope()\n\n"; 28a9c020eSBarry Smith 38a9c020eSBarry Smith #include <petscmat.h> 48a9c020eSBarry Smith extern PetscErrorCode MatIsDiagonal(Mat); 58a9c020eSBarry Smith extern PetscErrorCode BuildMatrix(const PetscInt *, PetscInt, const PetscInt *, Mat *); 68a9c020eSBarry Smith 7d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 8d71ae5a4SJacob Faibussowitsch { 98a9c020eSBarry Smith Mat A, C, D, F; 108a9c020eSBarry Smith PetscInt i, j, rows[2], *parts, cnt, N = 21, nblocks, *blocksizes; 118a9c020eSBarry Smith PetscScalar values[2][2]; 128a9c020eSBarry Smith PetscReal rand; 138a9c020eSBarry Smith PetscRandom rctx; 148a9c020eSBarry Smith PetscMPIInt size; 158a9c020eSBarry Smith 16327415f7SBarry Smith PetscFunctionBeginUser; 17c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 188a9c020eSBarry Smith PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_DENSE)); 198a9c020eSBarry Smith 208a9c020eSBarry Smith PetscCall(MatCreate(PETSC_COMM_WORLD, &C)); 218a9c020eSBarry Smith PetscCall(MatSetSizes(C, PETSC_DECIDE, PETSC_DECIDE, 6, 18)); 228a9c020eSBarry Smith PetscCall(MatSetFromOptions(C)); 238a9c020eSBarry Smith PetscCall(MatSetUp(C)); 249371c9d4SSatish Balay values[0][0] = 2; 259371c9d4SSatish Balay values[0][1] = 1; 269371c9d4SSatish Balay values[1][0] = 1; 279371c9d4SSatish Balay values[1][1] = 2; 288a9c020eSBarry Smith for (i = 0; i < 3; i++) { 299371c9d4SSatish Balay rows[0] = 2 * i; 309371c9d4SSatish Balay rows[1] = 2 * i + 1; 318a9c020eSBarry Smith PetscCall(MatSetValues(C, 2, rows, 2, rows, (PetscScalar *)values, INSERT_VALUES)); 328a9c020eSBarry Smith } 338a9c020eSBarry Smith PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY)); 348a9c020eSBarry Smith PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY)); 358a9c020eSBarry Smith PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD)); 368a9c020eSBarry Smith 378a9c020eSBarry Smith PetscCall(MatMatTransposeMult(C, C, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &A)); 388a9c020eSBarry Smith PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD)); 398a9c020eSBarry Smith 408a9c020eSBarry Smith PetscCall(MatInvertVariableBlockEnvelope(A, MAT_INITIAL_MATRIX, &D)); 418a9c020eSBarry Smith PetscCall(MatView(D, PETSC_VIEWER_STDOUT_WORLD)); 428a9c020eSBarry Smith 438a9c020eSBarry Smith PetscCall(MatMatMult(A, D, MAT_INITIAL_MATRIX, 1.0, &F)); 448a9c020eSBarry Smith PetscCall(MatView(F, PETSC_VIEWER_STDOUT_WORLD)); 458a9c020eSBarry Smith PetscCall(MatIsDiagonal(F)); 468a9c020eSBarry Smith 478a9c020eSBarry Smith PetscCall(MatDestroy(&A)); 488a9c020eSBarry Smith PetscCall(MatDestroy(&D)); 498a9c020eSBarry Smith PetscCall(MatDestroy(&C)); 508a9c020eSBarry Smith PetscCall(MatDestroy(&F)); 518a9c020eSBarry Smith 528a9c020eSBarry Smith PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rctx)); 538a9c020eSBarry Smith PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 548a9c020eSBarry Smith PetscCall(PetscMalloc1(size, &parts)); 558a9c020eSBarry Smith 568a9c020eSBarry Smith for (j = 0; j < 128; j++) { 578a9c020eSBarry Smith cnt = 0; 588a9c020eSBarry Smith for (i = 0; i < size - 1; i++) { 598a9c020eSBarry Smith PetscCall(PetscRandomGetValueReal(rctx, &rand)); 60*300f1712SStefano Zampini parts[i] = (PetscInt)(N * rand); 618a9c020eSBarry Smith parts[i] = PetscMin(parts[i], N - cnt); 628a9c020eSBarry Smith cnt += parts[i]; 638a9c020eSBarry Smith } 648a9c020eSBarry Smith parts[size - 1] = N - cnt; 658a9c020eSBarry Smith 668a9c020eSBarry Smith PetscCall(PetscRandomGetValueReal(rctx, &rand)); 678a9c020eSBarry Smith nblocks = rand * 10; 688a9c020eSBarry Smith nblocks = PetscMax(nblocks, 2); 698a9c020eSBarry Smith cnt = 0; 708a9c020eSBarry Smith PetscCall(PetscMalloc1(nblocks, &blocksizes)); 718a9c020eSBarry Smith for (i = 0; i < nblocks - 1; i++) { 728a9c020eSBarry Smith PetscCall(PetscRandomGetValueReal(rctx, &rand)); 73*300f1712SStefano Zampini blocksizes[i] = PetscMax(1, (PetscInt)(N * rand)); 748a9c020eSBarry Smith blocksizes[i] = PetscMin(blocksizes[i], N - cnt); 758a9c020eSBarry Smith cnt += blocksizes[i]; 768a9c020eSBarry Smith if (cnt == N) { 778a9c020eSBarry Smith nblocks = i + 1; 788a9c020eSBarry Smith break; 798a9c020eSBarry Smith } 808a9c020eSBarry Smith } 81ad540459SPierre Jolivet if (cnt < N) blocksizes[nblocks - 1] = N - cnt; 828a9c020eSBarry Smith 838a9c020eSBarry Smith PetscCall(BuildMatrix(parts, nblocks, blocksizes, &A)); 848a9c020eSBarry Smith PetscCall(PetscFree(blocksizes)); 858a9c020eSBarry Smith 868a9c020eSBarry Smith PetscCall(MatInvertVariableBlockEnvelope(A, MAT_INITIAL_MATRIX, &D)); 878a9c020eSBarry Smith 888a9c020eSBarry Smith PetscCall(MatMatMult(A, D, MAT_INITIAL_MATRIX, 1.0, &F)); 898a9c020eSBarry Smith PetscCall(MatIsDiagonal(F)); 908a9c020eSBarry Smith 918a9c020eSBarry Smith PetscCall(MatDestroy(&A)); 928a9c020eSBarry Smith PetscCall(MatDestroy(&D)); 938a9c020eSBarry Smith PetscCall(MatDestroy(&F)); 948a9c020eSBarry Smith } 958a9c020eSBarry Smith PetscCall(PetscFree(parts)); 968a9c020eSBarry Smith PetscCall(PetscRandomDestroy(&rctx)); 978a9c020eSBarry Smith 988a9c020eSBarry Smith PetscCall(PetscFinalize()); 998a9c020eSBarry Smith return 0; 1008a9c020eSBarry Smith } 1018a9c020eSBarry Smith 102d71ae5a4SJacob Faibussowitsch PetscErrorCode MatIsDiagonal(Mat A) 103d71ae5a4SJacob Faibussowitsch { 1048a9c020eSBarry Smith PetscInt ncols, i, j, rstart, rend; 1058a9c020eSBarry Smith const PetscInt *cols; 1068a9c020eSBarry Smith const PetscScalar *vals; 1078a9c020eSBarry Smith PetscBool founddiag; 1088a9c020eSBarry Smith 1098a9c020eSBarry Smith PetscFunctionBeginUser; 1108a9c020eSBarry Smith PetscCall(MatGetOwnershipRange(A, &rstart, &rend)); 1118a9c020eSBarry Smith for (i = rstart; i < rend; i++) { 1128a9c020eSBarry Smith founddiag = PETSC_FALSE; 1138a9c020eSBarry Smith PetscCall(MatGetRow(A, i, &ncols, &cols, &vals)); 1148a9c020eSBarry Smith for (j = 0; j < ncols; j++) { 1158a9c020eSBarry Smith if (cols[j] == i) { 1168a9c020eSBarry Smith PetscCheck(PetscAbsScalar(vals[j] - 1) < PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Row %" PetscInt_FMT " does not have 1 on the diagonal, it has %g", i, (double)PetscAbsScalar(vals[j])); 1178a9c020eSBarry Smith founddiag = PETSC_TRUE; 1188a9c020eSBarry Smith } else { 119e978a55eSPierre Jolivet PetscCheck(PetscAbsScalar(vals[j]) < PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Row %" PetscInt_FMT " has off-diagonal value %g at %" PetscInt_FMT, i, (double)PetscAbsScalar(vals[j]), cols[j]); 1208a9c020eSBarry Smith } 1218a9c020eSBarry Smith } 1228a9c020eSBarry Smith PetscCall(MatRestoreRow(A, i, &ncols, &cols, &vals)); 123aaa8cc7dSPierre Jolivet PetscCheck(founddiag, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Row %" PetscInt_FMT " does not have diagonal entry", i); 1248a9c020eSBarry Smith } 1253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1268a9c020eSBarry Smith } 1278a9c020eSBarry Smith 1288a9c020eSBarry Smith /* 1298a9c020eSBarry Smith All processes receive all the block information 1308a9c020eSBarry Smith */ 131d71ae5a4SJacob Faibussowitsch PetscErrorCode BuildMatrix(const PetscInt *parts, PetscInt nblocks, const PetscInt *blocksizes, Mat *A) 132d71ae5a4SJacob Faibussowitsch { 1338a9c020eSBarry Smith PetscInt i, cnt = 0; 1348a9c020eSBarry Smith PetscMPIInt rank; 1358a9c020eSBarry Smith 1368a9c020eSBarry Smith PetscFunctionBeginUser; 1378a9c020eSBarry Smith PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 1388a9c020eSBarry Smith PetscCall(MatCreateAIJ(PETSC_COMM_WORLD, parts[rank], parts[rank], PETSC_DETERMINE, PETSC_DETERMINE, 0, NULL, 0, NULL, A)); 1398a9c020eSBarry Smith PetscCall(MatSetOption(*A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE)); 140c5853193SPierre Jolivet if (rank == 0) { 1418a9c020eSBarry Smith for (i = 0; i < nblocks; i++) { 1428a9c020eSBarry Smith PetscCall(MatSetValue(*A, cnt, cnt + blocksizes[i] - 1, 1.0, INSERT_VALUES)); 1438a9c020eSBarry Smith PetscCall(MatSetValue(*A, cnt + blocksizes[i] - 1, cnt, 1.0, INSERT_VALUES)); 1448a9c020eSBarry Smith cnt += blocksizes[i]; 1458a9c020eSBarry Smith } 1468a9c020eSBarry Smith } 1478a9c020eSBarry Smith PetscCall(MatAssemblyBegin(*A, MAT_FINAL_ASSEMBLY)); 1488a9c020eSBarry Smith PetscCall(MatAssemblyEnd(*A, MAT_FINAL_ASSEMBLY)); 1498a9c020eSBarry Smith PetscCall(MatShift(*A, 10)); 1503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1518a9c020eSBarry Smith } 1528a9c020eSBarry Smith 1538a9c020eSBarry Smith /*TEST 1548a9c020eSBarry Smith 1558a9c020eSBarry Smith test: 1568a9c020eSBarry Smith 1578a9c020eSBarry Smith test: 1588a9c020eSBarry Smith suffix: 2 1598a9c020eSBarry Smith nsize: 2 1608a9c020eSBarry Smith 1618a9c020eSBarry Smith test: 1628a9c020eSBarry Smith suffix: 5 1638a9c020eSBarry Smith nsize: 5 1648a9c020eSBarry Smith 1658a9c020eSBarry Smith TEST*/ 166