xref: /petsc/src/mat/tests/ex178.c (revision aaa8cc7d2a5c3913edcbb923e20f154fe9c4aa65)
18a9c020eSBarry Smith 
28a9c020eSBarry Smith static char help[] = "Tests MatInvertVariableBlockEnvelope()\n\n";
38a9c020eSBarry Smith 
48a9c020eSBarry Smith #include <petscmat.h>
58a9c020eSBarry Smith extern PetscErrorCode MatIsDiagonal(Mat);
68a9c020eSBarry Smith extern PetscErrorCode BuildMatrix(const PetscInt *, PetscInt, const PetscInt *, Mat *);
78a9c020eSBarry Smith 
8d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
9d71ae5a4SJacob Faibussowitsch {
108a9c020eSBarry Smith   Mat         A, C, D, F;
118a9c020eSBarry Smith   PetscInt    i, j, rows[2], *parts, cnt, N = 21, nblocks, *blocksizes;
128a9c020eSBarry Smith   PetscScalar values[2][2];
138a9c020eSBarry Smith   PetscReal   rand;
148a9c020eSBarry Smith   PetscRandom rctx;
158a9c020eSBarry Smith   PetscMPIInt size;
168a9c020eSBarry Smith 
17327415f7SBarry Smith   PetscFunctionBeginUser;
188a9c020eSBarry Smith   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
198a9c020eSBarry Smith   PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_DENSE));
208a9c020eSBarry Smith 
218a9c020eSBarry Smith   PetscCall(MatCreate(PETSC_COMM_WORLD, &C));
228a9c020eSBarry Smith   PetscCall(MatSetSizes(C, PETSC_DECIDE, PETSC_DECIDE, 6, 18));
238a9c020eSBarry Smith   PetscCall(MatSetFromOptions(C));
248a9c020eSBarry Smith   PetscCall(MatSetUp(C));
259371c9d4SSatish Balay   values[0][0] = 2;
269371c9d4SSatish Balay   values[0][1] = 1;
279371c9d4SSatish Balay   values[1][0] = 1;
289371c9d4SSatish Balay   values[1][1] = 2;
298a9c020eSBarry Smith   for (i = 0; i < 3; i++) {
309371c9d4SSatish Balay     rows[0] = 2 * i;
319371c9d4SSatish Balay     rows[1] = 2 * i + 1;
328a9c020eSBarry Smith     PetscCall(MatSetValues(C, 2, rows, 2, rows, (PetscScalar *)values, INSERT_VALUES));
338a9c020eSBarry Smith   }
348a9c020eSBarry Smith   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
358a9c020eSBarry Smith   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
368a9c020eSBarry Smith   PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD));
378a9c020eSBarry Smith 
388a9c020eSBarry Smith   PetscCall(MatMatTransposeMult(C, C, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &A));
398a9c020eSBarry Smith   PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD));
408a9c020eSBarry Smith 
418a9c020eSBarry Smith   PetscCall(MatInvertVariableBlockEnvelope(A, MAT_INITIAL_MATRIX, &D));
428a9c020eSBarry Smith   PetscCall(MatView(D, PETSC_VIEWER_STDOUT_WORLD));
438a9c020eSBarry Smith 
448a9c020eSBarry Smith   PetscCall(MatMatMult(A, D, MAT_INITIAL_MATRIX, 1.0, &F));
458a9c020eSBarry Smith   PetscCall(MatView(F, PETSC_VIEWER_STDOUT_WORLD));
468a9c020eSBarry Smith   PetscCall(MatIsDiagonal(F));
478a9c020eSBarry Smith 
488a9c020eSBarry Smith   PetscCall(MatDestroy(&A));
498a9c020eSBarry Smith   PetscCall(MatDestroy(&D));
508a9c020eSBarry Smith   PetscCall(MatDestroy(&C));
518a9c020eSBarry Smith   PetscCall(MatDestroy(&F));
528a9c020eSBarry Smith 
538a9c020eSBarry Smith   PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rctx));
548a9c020eSBarry Smith   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
558a9c020eSBarry Smith   PetscCall(PetscMalloc1(size, &parts));
568a9c020eSBarry Smith 
578a9c020eSBarry Smith   for (j = 0; j < 128; j++) {
588a9c020eSBarry Smith     cnt = 0;
598a9c020eSBarry Smith     for (i = 0; i < size - 1; i++) {
608a9c020eSBarry Smith       PetscCall(PetscRandomGetValueReal(rctx, &rand));
618a9c020eSBarry Smith       parts[i] = (PetscInt)N * rand;
628a9c020eSBarry Smith       parts[i] = PetscMin(parts[i], N - cnt);
638a9c020eSBarry Smith       cnt += parts[i];
648a9c020eSBarry Smith     }
658a9c020eSBarry Smith     parts[size - 1] = N - cnt;
668a9c020eSBarry Smith 
678a9c020eSBarry Smith     PetscCall(PetscRandomGetValueReal(rctx, &rand));
688a9c020eSBarry Smith     nblocks = rand * 10;
698a9c020eSBarry Smith     nblocks = PetscMax(nblocks, 2);
708a9c020eSBarry Smith     cnt     = 0;
718a9c020eSBarry Smith     PetscCall(PetscMalloc1(nblocks, &blocksizes));
728a9c020eSBarry Smith     for (i = 0; i < nblocks - 1; i++) {
738a9c020eSBarry Smith       PetscCall(PetscRandomGetValueReal(rctx, &rand));
748a9c020eSBarry Smith       blocksizes[i] = PetscMax(1, (PetscInt)N * rand);
758a9c020eSBarry Smith       blocksizes[i] = PetscMin(blocksizes[i], N - cnt);
768a9c020eSBarry Smith       cnt += blocksizes[i];
778a9c020eSBarry Smith       if (cnt == N) {
788a9c020eSBarry Smith         nblocks = i + 1;
798a9c020eSBarry Smith         break;
808a9c020eSBarry Smith       }
818a9c020eSBarry Smith     }
82ad540459SPierre Jolivet     if (cnt < N) blocksizes[nblocks - 1] = N - cnt;
838a9c020eSBarry Smith 
848a9c020eSBarry Smith     PetscCall(BuildMatrix(parts, nblocks, blocksizes, &A));
858a9c020eSBarry Smith     PetscCall(PetscFree(blocksizes));
868a9c020eSBarry Smith 
878a9c020eSBarry Smith     PetscCall(MatInvertVariableBlockEnvelope(A, MAT_INITIAL_MATRIX, &D));
888a9c020eSBarry Smith 
898a9c020eSBarry Smith     PetscCall(MatMatMult(A, D, MAT_INITIAL_MATRIX, 1.0, &F));
908a9c020eSBarry Smith     PetscCall(MatIsDiagonal(F));
918a9c020eSBarry Smith 
928a9c020eSBarry Smith     PetscCall(MatDestroy(&A));
938a9c020eSBarry Smith     PetscCall(MatDestroy(&D));
948a9c020eSBarry Smith     PetscCall(MatDestroy(&F));
958a9c020eSBarry Smith   }
968a9c020eSBarry Smith   PetscCall(PetscFree(parts));
978a9c020eSBarry Smith   PetscCall(PetscRandomDestroy(&rctx));
988a9c020eSBarry Smith 
998a9c020eSBarry Smith   PetscCall(PetscFinalize());
1008a9c020eSBarry Smith   return 0;
1018a9c020eSBarry Smith }
1028a9c020eSBarry Smith 
103d71ae5a4SJacob Faibussowitsch PetscErrorCode MatIsDiagonal(Mat A)
104d71ae5a4SJacob Faibussowitsch {
1058a9c020eSBarry Smith   PetscInt           ncols, i, j, rstart, rend;
1068a9c020eSBarry Smith   const PetscInt    *cols;
1078a9c020eSBarry Smith   const PetscScalar *vals;
1088a9c020eSBarry Smith   PetscBool          founddiag;
1098a9c020eSBarry Smith 
1108a9c020eSBarry Smith   PetscFunctionBeginUser;
1118a9c020eSBarry Smith   PetscCall(MatGetOwnershipRange(A, &rstart, &rend));
1128a9c020eSBarry Smith   for (i = rstart; i < rend; i++) {
1138a9c020eSBarry Smith     founddiag = PETSC_FALSE;
1148a9c020eSBarry Smith     PetscCall(MatGetRow(A, i, &ncols, &cols, &vals));
1158a9c020eSBarry Smith     for (j = 0; j < ncols; j++) {
1168a9c020eSBarry Smith       if (cols[j] == i) {
1178a9c020eSBarry 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]));
1188a9c020eSBarry Smith         founddiag = PETSC_TRUE;
1198a9c020eSBarry Smith       } else {
1208a9c020eSBarry Smith         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]);
1218a9c020eSBarry Smith       }
1228a9c020eSBarry Smith     }
1238a9c020eSBarry Smith     PetscCall(MatRestoreRow(A, i, &ncols, &cols, &vals));
124*aaa8cc7dSPierre Jolivet     PetscCheck(founddiag, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Row %" PetscInt_FMT " does not have diagonal entry", i);
1258a9c020eSBarry Smith   }
1263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1278a9c020eSBarry Smith }
1288a9c020eSBarry Smith 
1298a9c020eSBarry Smith /*
1308a9c020eSBarry Smith     All processes receive all the block information
1318a9c020eSBarry Smith */
132d71ae5a4SJacob Faibussowitsch PetscErrorCode BuildMatrix(const PetscInt *parts, PetscInt nblocks, const PetscInt *blocksizes, Mat *A)
133d71ae5a4SJacob Faibussowitsch {
1348a9c020eSBarry Smith   PetscInt    i, cnt = 0;
1358a9c020eSBarry Smith   PetscMPIInt rank;
1368a9c020eSBarry Smith 
1378a9c020eSBarry Smith   PetscFunctionBeginUser;
1388a9c020eSBarry Smith   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
1398a9c020eSBarry Smith   PetscCall(MatCreateAIJ(PETSC_COMM_WORLD, parts[rank], parts[rank], PETSC_DETERMINE, PETSC_DETERMINE, 0, NULL, 0, NULL, A));
1408a9c020eSBarry Smith   PetscCall(MatSetOption(*A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE));
141c5853193SPierre Jolivet   if (rank == 0) {
1428a9c020eSBarry Smith     for (i = 0; i < nblocks; i++) {
1438a9c020eSBarry Smith       PetscCall(MatSetValue(*A, cnt, cnt + blocksizes[i] - 1, 1.0, INSERT_VALUES));
1448a9c020eSBarry Smith       PetscCall(MatSetValue(*A, cnt + blocksizes[i] - 1, cnt, 1.0, INSERT_VALUES));
1458a9c020eSBarry Smith       cnt += blocksizes[i];
1468a9c020eSBarry Smith     }
1478a9c020eSBarry Smith   }
1488a9c020eSBarry Smith   PetscCall(MatAssemblyBegin(*A, MAT_FINAL_ASSEMBLY));
1498a9c020eSBarry Smith   PetscCall(MatAssemblyEnd(*A, MAT_FINAL_ASSEMBLY));
1508a9c020eSBarry Smith   PetscCall(MatShift(*A, 10));
1513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1528a9c020eSBarry Smith }
1538a9c020eSBarry Smith 
1548a9c020eSBarry Smith /*TEST
1558a9c020eSBarry Smith 
1568a9c020eSBarry Smith    test:
1578a9c020eSBarry Smith 
1588a9c020eSBarry Smith    test:
1598a9c020eSBarry Smith      suffix: 2
1608a9c020eSBarry Smith      nsize: 2
1618a9c020eSBarry Smith 
1628a9c020eSBarry Smith    test:
1638a9c020eSBarry Smith      suffix: 5
1648a9c020eSBarry Smith      nsize: 5
1658a9c020eSBarry Smith 
1668a9c020eSBarry Smith TEST*/
167