1c4762a1bSJed Brown static char help[] = "Tests MatMatSolve() and MatMatTransposeSolve() for computing inv(A) with MUMPS.\n\
2c4762a1bSJed Brown Example: mpiexec -n <np> ./ex214 -displ \n\n";
3c4762a1bSJed Brown
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown
main(int argc,char ** args)6d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
7d71ae5a4SJacob Faibussowitsch {
8c4762a1bSJed Brown PetscMPIInt size, rank;
9c4762a1bSJed Brown #if defined(PETSC_HAVE_MUMPS)
10c4762a1bSJed Brown Mat A, RHS, C, F, X, AX, spRHST;
11c4762a1bSJed Brown PetscInt m, n, nrhs, M, N, i, Istart, Iend, Ii, j, J, test;
12c4762a1bSJed Brown PetscScalar v;
13c4762a1bSJed Brown PetscReal norm, tol = PETSC_SQRT_MACHINE_EPSILON;
14c4762a1bSJed Brown PetscRandom rand;
15c4762a1bSJed Brown PetscBool displ = PETSC_FALSE;
16c4762a1bSJed Brown char solver[256];
17c4762a1bSJed Brown #endif
18c4762a1bSJed Brown
19327415f7SBarry Smith PetscFunctionBeginUser;
20*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, help));
219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
229566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
23c4762a1bSJed Brown
24c4762a1bSJed Brown #if !defined(PETSC_HAVE_MUMPS)
259566063dSJacob Faibussowitsch if (rank == 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, "This example requires MUMPS, exit...\n"));
269566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
27b122ec5aSJacob Faibussowitsch return 0;
28c4762a1bSJed Brown #else
29c4762a1bSJed Brown
309566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-displ", &displ, NULL));
31c4762a1bSJed Brown
32c4762a1bSJed Brown /* Create matrix A */
339371c9d4SSatish Balay m = 4;
349371c9d4SSatish Balay n = 4;
359566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-m", &m, NULL));
369566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL));
37c4762a1bSJed Brown
389566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
399566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, m * n, m * n));
409566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A));
419566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(A, 5, NULL, 5, NULL));
429566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(A, 5, NULL));
43c4762a1bSJed Brown
449566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &Istart, &Iend));
45c4762a1bSJed Brown for (Ii = Istart; Ii < Iend; Ii++) {
469371c9d4SSatish Balay v = -1.0;
479371c9d4SSatish Balay i = Ii / n;
489371c9d4SSatish Balay j = Ii - i * n;
499371c9d4SSatish Balay if (i > 0) {
509371c9d4SSatish Balay J = Ii - n;
519371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, ADD_VALUES));
529371c9d4SSatish Balay }
539371c9d4SSatish Balay if (i < m - 1) {
549371c9d4SSatish Balay J = Ii + n;
559371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, ADD_VALUES));
569371c9d4SSatish Balay }
579371c9d4SSatish Balay if (j > 0) {
589371c9d4SSatish Balay J = Ii - 1;
599371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, ADD_VALUES));
609371c9d4SSatish Balay }
619371c9d4SSatish Balay if (j < n - 1) {
629371c9d4SSatish Balay J = Ii + 1;
639371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, ADD_VALUES));
649371c9d4SSatish Balay }
659371c9d4SSatish Balay v = 4.0;
669371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &Ii, 1, &Ii, &v, ADD_VALUES));
67c4762a1bSJed Brown }
689566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
699566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
70c4762a1bSJed Brown
719566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &m, &n));
729566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &M, &N));
7308401ef6SPierre 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);
74c4762a1bSJed Brown
75a5b23f4aSJose E. Roman /* Create dense matrix C and X; C holds true solution with identical columns */
76c4762a1bSJed Brown nrhs = N;
779566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-nrhs", &nrhs, NULL));
789566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &C));
799566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C, m, PETSC_DECIDE, PETSC_DECIDE, nrhs));
809566063dSJacob Faibussowitsch PetscCall(MatSetType(C, MATDENSE));
819566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(C));
829566063dSJacob Faibussowitsch PetscCall(MatSetUp(C));
83c4762a1bSJed Brown
849566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rand));
859566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rand));
869566063dSJacob Faibussowitsch PetscCall(MatSetRandom(C, rand));
879566063dSJacob Faibussowitsch PetscCall(MatDuplicate(C, MAT_DO_NOT_COPY_VALUES, &X));
88c4762a1bSJed Brown
89c6a7a370SJeremy L Thompson PetscCall(PetscStrncpy(solver, MATSOLVERMUMPS, sizeof(solver)));
909566063dSJacob 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));
91c4762a1bSJed Brown
92c4762a1bSJed Brown for (test = 0; test < 2; test++) {
93c4762a1bSJed Brown if (test == 0) {
94c4762a1bSJed Brown /* Test LU Factorization */
959566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "test LU factorization\n"));
969566063dSJacob Faibussowitsch PetscCall(MatGetFactor(A, solver, MAT_FACTOR_LU, &F));
979566063dSJacob Faibussowitsch PetscCall(MatLUFactorSymbolic(F, A, NULL, NULL, NULL));
989566063dSJacob Faibussowitsch PetscCall(MatLUFactorNumeric(F, A, NULL));
99c4762a1bSJed Brown } else {
100c4762a1bSJed Brown /* Test Cholesky Factorization */
101c4762a1bSJed Brown PetscBool flg;
1029566063dSJacob Faibussowitsch PetscCall(MatIsSymmetric(A, 0.0, &flg));
10328b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "A must be symmetric for Cholesky factorization");
104c4762a1bSJed Brown
1059566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "test Cholesky factorization\n"));
1069566063dSJacob Faibussowitsch PetscCall(MatGetFactor(A, solver, MAT_FACTOR_CHOLESKY, &F));
1079566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorSymbolic(F, A, NULL, NULL));
1089566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorNumeric(F, A, NULL));
109c4762a1bSJed Brown }
110c4762a1bSJed Brown
111c4762a1bSJed Brown /* (1) Test MatMatSolve(): dense RHS = A*C, C: true solutions */
112c4762a1bSJed Brown /* ---------------------------------------------------------- */
1139566063dSJacob Faibussowitsch PetscCall(MatMatMult(A, C, MAT_INITIAL_MATRIX, 2.0, &RHS));
1149566063dSJacob Faibussowitsch PetscCall(MatMatSolve(F, RHS, X));
115c4762a1bSJed Brown /* Check the error */
1169566063dSJacob Faibussowitsch PetscCall(MatAXPY(X, -1.0, C, SAME_NONZERO_PATTERN));
1179566063dSJacob Faibussowitsch PetscCall(MatNorm(X, NORM_FROBENIUS, &norm));
11848a46eb9SPierre Jolivet if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "(1) MatMatSolve: Norm of error %g\n", norm));
119c4762a1bSJed Brown
120c4762a1bSJed Brown /* Test X=RHS */
1219566063dSJacob Faibussowitsch PetscCall(MatMatSolve(F, RHS, RHS));
122c4762a1bSJed Brown /* Check the error */
1239566063dSJacob Faibussowitsch PetscCall(MatAXPY(RHS, -1.0, C, SAME_NONZERO_PATTERN));
1249566063dSJacob Faibussowitsch PetscCall(MatNorm(RHS, NORM_FROBENIUS, &norm));
12548a46eb9SPierre Jolivet if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "(1.1) MatMatSolve: Norm of error %g\n", norm));
126c4762a1bSJed Brown
127c4762a1bSJed Brown /* (2) Test MatMatSolve() for inv(A) with dense RHS:
128c4762a1bSJed Brown RHS = [e[0],...,e[nrhs-1]], dense X holds first nrhs columns of inv(A) */
129c4762a1bSJed Brown /* -------------------------------------------------------------------- */
1309566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(RHS));
131c4762a1bSJed Brown for (i = 0; i < nrhs; i++) {
132c4762a1bSJed Brown v = 1.0;
1339566063dSJacob Faibussowitsch PetscCall(MatSetValues(RHS, 1, &i, 1, &i, &v, INSERT_VALUES));
134c4762a1bSJed Brown }
1359566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(RHS, MAT_FINAL_ASSEMBLY));
1369566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(RHS, MAT_FINAL_ASSEMBLY));
137c4762a1bSJed Brown
1389566063dSJacob Faibussowitsch PetscCall(MatMatSolve(F, RHS, X));
139c4762a1bSJed Brown if (displ) {
1409566063dSJacob Faibussowitsch if (rank == 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, " \n(2) first %" PetscInt_FMT " columns of inv(A) with dense RHS:\n", nrhs));
1419566063dSJacob Faibussowitsch PetscCall(MatView(X, PETSC_VIEWER_STDOUT_WORLD));
142c4762a1bSJed Brown }
143c4762a1bSJed Brown
144c4762a1bSJed Brown /* Check the residual */
1459566063dSJacob Faibussowitsch PetscCall(MatMatMult(A, X, MAT_INITIAL_MATRIX, 2.0, &AX));
1469566063dSJacob Faibussowitsch PetscCall(MatAXPY(AX, -1.0, RHS, SAME_NONZERO_PATTERN));
1479566063dSJacob Faibussowitsch PetscCall(MatNorm(AX, NORM_INFINITY, &norm));
14848a46eb9SPierre Jolivet if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "(2) MatMatSolve: Norm of residual %g\n", norm));
1499566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(X));
150c4762a1bSJed Brown
151c4762a1bSJed Brown /* (3) Test MatMatTransposeSolve() for inv(A) with sparse RHS stored in the host:
152c4762a1bSJed Brown spRHST = [e[0],...,e[nrhs-1]]^T, dense X holds first nrhs columns of inv(A) */
153c4762a1bSJed Brown /* --------------------------------------------------------------------------- */
154c4762a1bSJed Brown /* Create spRHST: PETSc does not support compressed column format which is required by MUMPS for sparse RHS matrix,
155c4762a1bSJed Brown thus user must create spRHST=spRHS^T and call MatMatTransposeSolve() */
1569566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &spRHST));
157dd400576SPatrick Sanan if (rank == 0) {
158145b44c9SPierre Jolivet /* MUMPS requires RHS be centralized on the host! */
1599566063dSJacob Faibussowitsch PetscCall(MatSetSizes(spRHST, nrhs, M, PETSC_DECIDE, PETSC_DECIDE));
160c4762a1bSJed Brown } else {
1619566063dSJacob Faibussowitsch PetscCall(MatSetSizes(spRHST, 0, 0, PETSC_DECIDE, PETSC_DECIDE));
162c4762a1bSJed Brown }
1639566063dSJacob Faibussowitsch PetscCall(MatSetType(spRHST, MATAIJ));
1649566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(spRHST));
1659566063dSJacob Faibussowitsch PetscCall(MatSetUp(spRHST));
166dd400576SPatrick Sanan if (rank == 0) {
167c4762a1bSJed Brown v = 1.0;
16848a46eb9SPierre Jolivet for (i = 0; i < nrhs; i++) PetscCall(MatSetValues(spRHST, 1, &i, 1, &i, &v, INSERT_VALUES));
169c4762a1bSJed Brown }
1709566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(spRHST, MAT_FINAL_ASSEMBLY));
1719566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(spRHST, MAT_FINAL_ASSEMBLY));
172c4762a1bSJed Brown
1739566063dSJacob Faibussowitsch PetscCall(MatMatTransposeSolve(F, spRHST, X));
174c4762a1bSJed Brown
175c4762a1bSJed Brown if (displ) {
1769566063dSJacob Faibussowitsch if (rank == 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, " \n(3) first %" PetscInt_FMT " columns of inv(A) with sparse RHS:\n", nrhs));
1779566063dSJacob Faibussowitsch PetscCall(MatView(X, PETSC_VIEWER_STDOUT_WORLD));
178c4762a1bSJed Brown }
179c4762a1bSJed Brown
180c4762a1bSJed Brown /* Check the residual: R = A*X - RHS */
1819566063dSJacob Faibussowitsch PetscCall(MatMatMult(A, X, MAT_REUSE_MATRIX, 2.0, &AX));
182c4762a1bSJed Brown
1839566063dSJacob Faibussowitsch PetscCall(MatAXPY(AX, -1.0, RHS, SAME_NONZERO_PATTERN));
1849566063dSJacob Faibussowitsch PetscCall(MatNorm(AX, NORM_INFINITY, &norm));
18548a46eb9SPierre Jolivet if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "(3) MatMatSolve: Norm of residual %g\n", norm));
186c4762a1bSJed Brown
187c4762a1bSJed Brown /* (4) Test MatMatSolve() for inv(A) with selected entries:
188c4762a1bSJed Brown input: spRHS gives selected indices; output: spRHS holds selected entries of inv(A) */
189c4762a1bSJed Brown /* --------------------------------------------------------------------------------- */
190c4762a1bSJed Brown if (nrhs == N) { /* mumps requires nrhs = n */
191c4762a1bSJed Brown /* Create spRHS on proc[0] */
192c4762a1bSJed Brown Mat spRHS = NULL;
193c4762a1bSJed Brown
194c4762a1bSJed Brown /* Create spRHS = spRHST^T. Two matrices share internal matrix data structure */
1959566063dSJacob Faibussowitsch PetscCall(MatCreateTranspose(spRHST, &spRHS));
1969566063dSJacob Faibussowitsch PetscCall(MatMumpsGetInverse(F, spRHS));
1979566063dSJacob Faibussowitsch PetscCall(MatDestroy(&spRHS));
198c4762a1bSJed Brown
1999566063dSJacob Faibussowitsch PetscCall(MatMumpsGetInverseTranspose(F, spRHST));
200c4762a1bSJed Brown if (displ) {
2019566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nSelected entries of inv(A^T):\n"));
2029566063dSJacob Faibussowitsch PetscCall(MatView(spRHST, PETSC_VIEWER_STDOUT_WORLD));
203c4762a1bSJed Brown }
2049566063dSJacob Faibussowitsch PetscCall(MatDestroy(&spRHS));
205c4762a1bSJed Brown }
2069566063dSJacob Faibussowitsch PetscCall(MatDestroy(&AX));
2079566063dSJacob Faibussowitsch PetscCall(MatDestroy(&F));
2089566063dSJacob Faibussowitsch PetscCall(MatDestroy(&RHS));
2099566063dSJacob Faibussowitsch PetscCall(MatDestroy(&spRHST));
210c4762a1bSJed Brown }
211c4762a1bSJed Brown
212c4762a1bSJed Brown /* Free data structures */
2139566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A));
2149566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C));
2159566063dSJacob Faibussowitsch PetscCall(MatDestroy(&X));
2169566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rand));
2179566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
218b122ec5aSJacob Faibussowitsch return 0;
219c4762a1bSJed Brown #endif
220c4762a1bSJed Brown }
221c4762a1bSJed Brown
222c4762a1bSJed Brown /*TEST
223c4762a1bSJed Brown
224c4762a1bSJed Brown test:
225c4762a1bSJed Brown requires: mumps double !complex
226c4762a1bSJed Brown
227c4762a1bSJed Brown test:
228c4762a1bSJed Brown suffix: 2
229c4762a1bSJed Brown requires: mumps double !complex
230c4762a1bSJed Brown nsize: 2
231c4762a1bSJed Brown
232c4762a1bSJed Brown TEST*/
233