xref: /petsc/src/mat/tests/ex214.c (revision c6a7a37075f8bf8d34d92c4910d42445b7a3482d)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests MatMatSolve() and MatMatTransposeSolve() for computing inv(A) with MUMPS.\n\
3c4762a1bSJed Brown Example: mpiexec -n <np> ./ex214 -displ \n\n";
4c4762a1bSJed Brown 
5c4762a1bSJed Brown #include <petscmat.h>
6c4762a1bSJed Brown 
7d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
8d71ae5a4SJacob Faibussowitsch {
9c4762a1bSJed Brown   PetscMPIInt size, rank;
10c4762a1bSJed Brown #if defined(PETSC_HAVE_MUMPS)
11c4762a1bSJed Brown   Mat         A, RHS, C, F, X, AX, spRHST;
12c4762a1bSJed Brown   PetscInt    m, n, nrhs, M, N, i, Istart, Iend, Ii, j, J, test;
13c4762a1bSJed Brown   PetscScalar v;
14c4762a1bSJed Brown   PetscReal   norm, tol = PETSC_SQRT_MACHINE_EPSILON;
15c4762a1bSJed Brown   PetscRandom rand;
16c4762a1bSJed Brown   PetscBool   displ = PETSC_FALSE;
17c4762a1bSJed Brown   char        solver[256];
18c4762a1bSJed Brown #endif
19c4762a1bSJed Brown 
20327415f7SBarry Smith   PetscFunctionBeginUser;
219566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
229566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
239566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
24c4762a1bSJed Brown 
25c4762a1bSJed Brown #if !defined(PETSC_HAVE_MUMPS)
269566063dSJacob Faibussowitsch   if (rank == 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, "This example requires MUMPS, exit...\n"));
279566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
28b122ec5aSJacob Faibussowitsch   return 0;
29c4762a1bSJed Brown #else
30c4762a1bSJed Brown 
319566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-displ", &displ, NULL));
32c4762a1bSJed Brown 
33c4762a1bSJed Brown   /* Create matrix A */
349371c9d4SSatish Balay   m = 4;
359371c9d4SSatish Balay   n = 4;
369566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-m", &m, NULL));
379566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL));
38c4762a1bSJed Brown 
399566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
409566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, m * n, m * n));
419566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
429566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(A, 5, NULL, 5, NULL));
439566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(A, 5, NULL));
44c4762a1bSJed Brown 
459566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A, &Istart, &Iend));
46c4762a1bSJed Brown   for (Ii = Istart; Ii < Iend; Ii++) {
479371c9d4SSatish Balay     v = -1.0;
489371c9d4SSatish Balay     i = Ii / n;
499371c9d4SSatish Balay     j = Ii - i * n;
509371c9d4SSatish Balay     if (i > 0) {
519371c9d4SSatish Balay       J = Ii - n;
529371c9d4SSatish Balay       PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, ADD_VALUES));
539371c9d4SSatish Balay     }
549371c9d4SSatish Balay     if (i < m - 1) {
559371c9d4SSatish Balay       J = Ii + n;
569371c9d4SSatish Balay       PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, ADD_VALUES));
579371c9d4SSatish Balay     }
589371c9d4SSatish Balay     if (j > 0) {
599371c9d4SSatish Balay       J = Ii - 1;
609371c9d4SSatish Balay       PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, ADD_VALUES));
619371c9d4SSatish Balay     }
629371c9d4SSatish Balay     if (j < n - 1) {
639371c9d4SSatish Balay       J = Ii + 1;
649371c9d4SSatish Balay       PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, ADD_VALUES));
659371c9d4SSatish Balay     }
669371c9d4SSatish Balay     v = 4.0;
679371c9d4SSatish Balay     PetscCall(MatSetValues(A, 1, &Ii, 1, &Ii, &v, ADD_VALUES));
68c4762a1bSJed Brown   }
699566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
709566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
71c4762a1bSJed Brown 
729566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &m, &n));
739566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, &M, &N));
7408401ef6SPierre Jolivet   PetscCheck(m == n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "This example is not intended for rectangular matrices (%" PetscInt_FMT ", %" PetscInt_FMT ")", m, n);
75c4762a1bSJed Brown 
76a5b23f4aSJose E. Roman   /* Create dense matrix C and X; C holds true solution with identical columns */
77c4762a1bSJed Brown   nrhs = N;
789566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-nrhs", &nrhs, NULL));
799566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &C));
809566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(C, m, PETSC_DECIDE, PETSC_DECIDE, nrhs));
819566063dSJacob Faibussowitsch   PetscCall(MatSetType(C, MATDENSE));
829566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(C));
839566063dSJacob Faibussowitsch   PetscCall(MatSetUp(C));
84c4762a1bSJed Brown 
859566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rand));
869566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rand));
879566063dSJacob Faibussowitsch   PetscCall(MatSetRandom(C, rand));
889566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(C, MAT_DO_NOT_COPY_VALUES, &X));
89c4762a1bSJed Brown 
90*c6a7a370SJeremy L Thompson   PetscCall(PetscStrncpy(solver, MATSOLVERMUMPS, sizeof(solver)));
919566063dSJacob Faibussowitsch   if (rank == 0 && displ) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Solving with %s: nrhs %" PetscInt_FMT ", size mat %" PetscInt_FMT " x %" PetscInt_FMT "\n", solver, nrhs, M, N));
92c4762a1bSJed Brown 
93c4762a1bSJed Brown   for (test = 0; test < 2; test++) {
94c4762a1bSJed Brown     if (test == 0) {
95c4762a1bSJed Brown       /* Test LU Factorization */
969566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "test LU factorization\n"));
979566063dSJacob Faibussowitsch       PetscCall(MatGetFactor(A, solver, MAT_FACTOR_LU, &F));
989566063dSJacob Faibussowitsch       PetscCall(MatLUFactorSymbolic(F, A, NULL, NULL, NULL));
999566063dSJacob Faibussowitsch       PetscCall(MatLUFactorNumeric(F, A, NULL));
100c4762a1bSJed Brown     } else {
101c4762a1bSJed Brown       /* Test Cholesky Factorization */
102c4762a1bSJed Brown       PetscBool flg;
1039566063dSJacob Faibussowitsch       PetscCall(MatIsSymmetric(A, 0.0, &flg));
10428b400f6SJacob Faibussowitsch       PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "A must be symmetric for Cholesky factorization");
105c4762a1bSJed Brown 
1069566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "test Cholesky factorization\n"));
1079566063dSJacob Faibussowitsch       PetscCall(MatGetFactor(A, solver, MAT_FACTOR_CHOLESKY, &F));
1089566063dSJacob Faibussowitsch       PetscCall(MatCholeskyFactorSymbolic(F, A, NULL, NULL));
1099566063dSJacob Faibussowitsch       PetscCall(MatCholeskyFactorNumeric(F, A, NULL));
110c4762a1bSJed Brown     }
111c4762a1bSJed Brown 
112c4762a1bSJed Brown     /* (1) Test MatMatSolve(): dense RHS = A*C, C: true solutions */
113c4762a1bSJed Brown     /* ---------------------------------------------------------- */
1149566063dSJacob Faibussowitsch     PetscCall(MatMatMult(A, C, MAT_INITIAL_MATRIX, 2.0, &RHS));
1159566063dSJacob Faibussowitsch     PetscCall(MatMatSolve(F, RHS, X));
116c4762a1bSJed Brown     /* Check the error */
1179566063dSJacob Faibussowitsch     PetscCall(MatAXPY(X, -1.0, C, SAME_NONZERO_PATTERN));
1189566063dSJacob Faibussowitsch     PetscCall(MatNorm(X, NORM_FROBENIUS, &norm));
11948a46eb9SPierre Jolivet     if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "(1) MatMatSolve: Norm of error %g\n", norm));
120c4762a1bSJed Brown 
121c4762a1bSJed Brown     /* Test X=RHS */
1229566063dSJacob Faibussowitsch     PetscCall(MatMatSolve(F, RHS, RHS));
123c4762a1bSJed Brown     /* Check the error */
1249566063dSJacob Faibussowitsch     PetscCall(MatAXPY(RHS, -1.0, C, SAME_NONZERO_PATTERN));
1259566063dSJacob Faibussowitsch     PetscCall(MatNorm(RHS, NORM_FROBENIUS, &norm));
12648a46eb9SPierre Jolivet     if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "(1.1) MatMatSolve: Norm of error %g\n", norm));
127c4762a1bSJed Brown 
128c4762a1bSJed Brown     /* (2) Test MatMatSolve() for inv(A) with dense RHS:
129c4762a1bSJed Brown      RHS = [e[0],...,e[nrhs-1]], dense X holds first nrhs columns of inv(A) */
130c4762a1bSJed Brown     /* -------------------------------------------------------------------- */
1319566063dSJacob Faibussowitsch     PetscCall(MatZeroEntries(RHS));
132c4762a1bSJed Brown     for (i = 0; i < nrhs; i++) {
133c4762a1bSJed Brown       v = 1.0;
1349566063dSJacob Faibussowitsch       PetscCall(MatSetValues(RHS, 1, &i, 1, &i, &v, INSERT_VALUES));
135c4762a1bSJed Brown     }
1369566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(RHS, MAT_FINAL_ASSEMBLY));
1379566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(RHS, MAT_FINAL_ASSEMBLY));
138c4762a1bSJed Brown 
1399566063dSJacob Faibussowitsch     PetscCall(MatMatSolve(F, RHS, X));
140c4762a1bSJed Brown     if (displ) {
1419566063dSJacob Faibussowitsch       if (rank == 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, " \n(2) first %" PetscInt_FMT " columns of inv(A) with dense RHS:\n", nrhs));
1429566063dSJacob Faibussowitsch       PetscCall(MatView(X, PETSC_VIEWER_STDOUT_WORLD));
143c4762a1bSJed Brown     }
144c4762a1bSJed Brown 
145c4762a1bSJed Brown     /* Check the residual */
1469566063dSJacob Faibussowitsch     PetscCall(MatMatMult(A, X, MAT_INITIAL_MATRIX, 2.0, &AX));
1479566063dSJacob Faibussowitsch     PetscCall(MatAXPY(AX, -1.0, RHS, SAME_NONZERO_PATTERN));
1489566063dSJacob Faibussowitsch     PetscCall(MatNorm(AX, NORM_INFINITY, &norm));
14948a46eb9SPierre Jolivet     if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "(2) MatMatSolve: Norm of residual %g\n", norm));
1509566063dSJacob Faibussowitsch     PetscCall(MatZeroEntries(X));
151c4762a1bSJed Brown 
152c4762a1bSJed Brown     /* (3) Test MatMatTransposeSolve() for inv(A) with sparse RHS stored in the host:
153c4762a1bSJed Brown      spRHST = [e[0],...,e[nrhs-1]]^T, dense X holds first nrhs columns of inv(A) */
154c4762a1bSJed Brown     /* --------------------------------------------------------------------------- */
155c4762a1bSJed Brown     /* Create spRHST: PETSc does not support compressed column format which is required by MUMPS for sparse RHS matrix,
156c4762a1bSJed Brown      thus user must create spRHST=spRHS^T and call MatMatTransposeSolve() */
1579566063dSJacob Faibussowitsch     PetscCall(MatCreate(PETSC_COMM_WORLD, &spRHST));
158dd400576SPatrick Sanan     if (rank == 0) {
159c4762a1bSJed Brown       /* MUMPS requirs RHS be centralized on the host! */
1609566063dSJacob Faibussowitsch       PetscCall(MatSetSizes(spRHST, nrhs, M, PETSC_DECIDE, PETSC_DECIDE));
161c4762a1bSJed Brown     } else {
1629566063dSJacob Faibussowitsch       PetscCall(MatSetSizes(spRHST, 0, 0, PETSC_DECIDE, PETSC_DECIDE));
163c4762a1bSJed Brown     }
1649566063dSJacob Faibussowitsch     PetscCall(MatSetType(spRHST, MATAIJ));
1659566063dSJacob Faibussowitsch     PetscCall(MatSetFromOptions(spRHST));
1669566063dSJacob Faibussowitsch     PetscCall(MatSetUp(spRHST));
167dd400576SPatrick Sanan     if (rank == 0) {
168c4762a1bSJed Brown       v = 1.0;
16948a46eb9SPierre Jolivet       for (i = 0; i < nrhs; i++) PetscCall(MatSetValues(spRHST, 1, &i, 1, &i, &v, INSERT_VALUES));
170c4762a1bSJed Brown     }
1719566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(spRHST, MAT_FINAL_ASSEMBLY));
1729566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(spRHST, MAT_FINAL_ASSEMBLY));
173c4762a1bSJed Brown 
1749566063dSJacob Faibussowitsch     PetscCall(MatMatTransposeSolve(F, spRHST, X));
175c4762a1bSJed Brown 
176c4762a1bSJed Brown     if (displ) {
1779566063dSJacob Faibussowitsch       if (rank == 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, " \n(3) first %" PetscInt_FMT " columns of inv(A) with sparse RHS:\n", nrhs));
1789566063dSJacob Faibussowitsch       PetscCall(MatView(X, PETSC_VIEWER_STDOUT_WORLD));
179c4762a1bSJed Brown     }
180c4762a1bSJed Brown 
181c4762a1bSJed Brown     /* Check the residual: R = A*X - RHS */
1829566063dSJacob Faibussowitsch     PetscCall(MatMatMult(A, X, MAT_REUSE_MATRIX, 2.0, &AX));
183c4762a1bSJed Brown 
1849566063dSJacob Faibussowitsch     PetscCall(MatAXPY(AX, -1.0, RHS, SAME_NONZERO_PATTERN));
1859566063dSJacob Faibussowitsch     PetscCall(MatNorm(AX, NORM_INFINITY, &norm));
18648a46eb9SPierre Jolivet     if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "(3) MatMatSolve: Norm of residual %g\n", norm));
187c4762a1bSJed Brown 
188c4762a1bSJed Brown     /* (4) Test MatMatSolve() for inv(A) with selected entries:
189c4762a1bSJed Brown      input: spRHS gives selected indices; output: spRHS holds selected entries of inv(A) */
190c4762a1bSJed Brown     /* --------------------------------------------------------------------------------- */
191c4762a1bSJed Brown     if (nrhs == N) { /* mumps requires nrhs = n */
192c4762a1bSJed Brown       /* Create spRHS on proc[0] */
193c4762a1bSJed Brown       Mat spRHS = NULL;
194c4762a1bSJed Brown 
195c4762a1bSJed Brown       /* Create spRHS = spRHST^T. Two matrices share internal matrix data structure */
1969566063dSJacob Faibussowitsch       PetscCall(MatCreateTranspose(spRHST, &spRHS));
1979566063dSJacob Faibussowitsch       PetscCall(MatMumpsGetInverse(F, spRHS));
1989566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&spRHS));
199c4762a1bSJed Brown 
2009566063dSJacob Faibussowitsch       PetscCall(MatMumpsGetInverseTranspose(F, spRHST));
201c4762a1bSJed Brown       if (displ) {
2029566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nSelected entries of inv(A^T):\n"));
2039566063dSJacob Faibussowitsch         PetscCall(MatView(spRHST, PETSC_VIEWER_STDOUT_WORLD));
204c4762a1bSJed Brown       }
2059566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&spRHS));
206c4762a1bSJed Brown     }
2079566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&AX));
2089566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&F));
2099566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&RHS));
2109566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&spRHST));
211c4762a1bSJed Brown   }
212c4762a1bSJed Brown 
213c4762a1bSJed Brown   /* Free data structures */
2149566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
2159566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
2169566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&X));
2179566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rand));
2189566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
219b122ec5aSJacob Faibussowitsch   return 0;
220c4762a1bSJed Brown #endif
221c4762a1bSJed Brown }
222c4762a1bSJed Brown 
223c4762a1bSJed Brown /*TEST
224c4762a1bSJed Brown 
225c4762a1bSJed Brown    test:
226c4762a1bSJed Brown      requires: mumps double !complex
227c4762a1bSJed Brown 
228c4762a1bSJed Brown    test:
229c4762a1bSJed Brown      suffix: 2
230c4762a1bSJed Brown      requires: mumps double !complex
231c4762a1bSJed Brown      nsize: 2
232c4762a1bSJed Brown 
233c4762a1bSJed Brown TEST*/
234