xref: /petsc/src/mat/tests/ex214.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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 
7c4762a1bSJed Brown int main(int argc,char **args)
8c4762a1bSJed Brown {
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 
20*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help));
215f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
225f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
23c4762a1bSJed Brown 
24c4762a1bSJed Brown #if !defined(PETSC_HAVE_MUMPS)
255f80ce2aSJacob Faibussowitsch   if (rank == 0) CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"This example requires MUMPS, exit...\n"));
26*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
27*b122ec5aSJacob Faibussowitsch   return 0;
28c4762a1bSJed Brown #else
29c4762a1bSJed Brown 
305f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-displ",&displ,NULL));
31c4762a1bSJed Brown 
32c4762a1bSJed Brown   /* Create matrix A */
33c4762a1bSJed Brown   m = 4; n = 4;
345f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
355f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
36c4762a1bSJed Brown 
375f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
385f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n));
395f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(A));
405f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMPIAIJSetPreallocation(A,5,NULL,5,NULL));
415f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSeqAIJSetPreallocation(A,5,NULL));
42c4762a1bSJed Brown 
435f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOwnershipRange(A,&Istart,&Iend));
44c4762a1bSJed Brown   for (Ii=Istart; Ii<Iend; Ii++) {
45c4762a1bSJed Brown     v = -1.0; i = Ii/n; j = Ii - i*n;
465f80ce2aSJacob Faibussowitsch     if (i>0)   {J = Ii - n; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES));}
475f80ce2aSJacob Faibussowitsch     if (i<m-1) {J = Ii + n; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES));}
485f80ce2aSJacob Faibussowitsch     if (j>0)   {J = Ii - 1; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES));}
495f80ce2aSJacob Faibussowitsch     if (j<n-1) {J = Ii + 1; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES));}
505f80ce2aSJacob Faibussowitsch     v = 4.0; CHKERRQ(MatSetValues(A,1,&Ii,1,&Ii,&v,ADD_VALUES));
51c4762a1bSJed Brown   }
525f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
535f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
54c4762a1bSJed Brown 
555f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetLocalSize(A,&m,&n));
565f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetSize(A,&M,&N));
572c71b3e2SJacob Faibussowitsch   PetscCheckFalse(m != n,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ, "This example is not intended for rectangular matrices (%" PetscInt_FMT ", %" PetscInt_FMT ")", m, n);
58c4762a1bSJed Brown 
59a5b23f4aSJose E. Roman   /* Create dense matrix C and X; C holds true solution with identical columns */
60c4762a1bSJed Brown   nrhs = N;
615f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-nrhs",&nrhs,NULL));
625f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&C));
635f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(C,m,PETSC_DECIDE,PETSC_DECIDE,nrhs));
645f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(C,MATDENSE));
655f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(C));
665f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(C));
67c4762a1bSJed Brown 
685f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&rand));
695f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomSetFromOptions(rand));
705f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetRandom(C,rand));
715f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDuplicate(C,MAT_DO_NOT_COPY_VALUES,&X));
72c4762a1bSJed Brown 
735f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrcpy(solver,MATSOLVERMUMPS));
745f80ce2aSJacob Faibussowitsch   if (rank == 0 && displ) CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Solving with %s: nrhs %" PetscInt_FMT ", size mat %" PetscInt_FMT " x %" PetscInt_FMT "\n",solver,nrhs,M,N));
75c4762a1bSJed Brown 
76c4762a1bSJed Brown   for (test=0; test<2; test++) {
77c4762a1bSJed Brown     if (test == 0) {
78c4762a1bSJed Brown       /* Test LU Factorization */
795f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"test LU factorization\n"));
805f80ce2aSJacob Faibussowitsch       CHKERRQ(MatGetFactor(A,solver,MAT_FACTOR_LU,&F));
815f80ce2aSJacob Faibussowitsch       CHKERRQ(MatLUFactorSymbolic(F,A,NULL,NULL,NULL));
825f80ce2aSJacob Faibussowitsch       CHKERRQ(MatLUFactorNumeric(F,A,NULL));
83c4762a1bSJed Brown     } else {
84c4762a1bSJed Brown       /* Test Cholesky Factorization */
85c4762a1bSJed Brown       PetscBool flg;
865f80ce2aSJacob Faibussowitsch       CHKERRQ(MatIsSymmetric(A,0.0,&flg));
8728b400f6SJacob Faibussowitsch       PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"A must be symmetric for Cholesky factorization");
88c4762a1bSJed Brown 
895f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"test Cholesky factorization\n"));
905f80ce2aSJacob Faibussowitsch       CHKERRQ(MatGetFactor(A,solver,MAT_FACTOR_CHOLESKY,&F));
915f80ce2aSJacob Faibussowitsch       CHKERRQ(MatCholeskyFactorSymbolic(F,A,NULL,NULL));
925f80ce2aSJacob Faibussowitsch       CHKERRQ(MatCholeskyFactorNumeric(F,A,NULL));
93c4762a1bSJed Brown     }
94c4762a1bSJed Brown 
95c4762a1bSJed Brown     /* (1) Test MatMatSolve(): dense RHS = A*C, C: true solutions */
96c4762a1bSJed Brown     /* ---------------------------------------------------------- */
975f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMatMult(A,C,MAT_INITIAL_MATRIX,2.0,&RHS));
985f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMatSolve(F,RHS,X));
99c4762a1bSJed Brown     /* Check the error */
1005f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAXPY(X,-1.0,C,SAME_NONZERO_PATTERN));
1015f80ce2aSJacob Faibussowitsch     CHKERRQ(MatNorm(X,NORM_FROBENIUS,&norm));
102c4762a1bSJed Brown     if (norm > tol) {
1035f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"(1) MatMatSolve: Norm of error %g\n",norm));
104c4762a1bSJed Brown     }
105c4762a1bSJed Brown 
106c4762a1bSJed Brown     /* Test X=RHS */
1075f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMatSolve(F,RHS,RHS));
108c4762a1bSJed Brown     /* Check the error */
1095f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAXPY(RHS,-1.0,C,SAME_NONZERO_PATTERN));
1105f80ce2aSJacob Faibussowitsch     CHKERRQ(MatNorm(RHS,NORM_FROBENIUS,&norm));
111c4762a1bSJed Brown     if (norm > tol) {
1125f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"(1.1) MatMatSolve: Norm of error %g\n",norm));
113c4762a1bSJed Brown     }
114c4762a1bSJed Brown 
115c4762a1bSJed Brown     /* (2) Test MatMatSolve() for inv(A) with dense RHS:
116c4762a1bSJed Brown      RHS = [e[0],...,e[nrhs-1]], dense X holds first nrhs columns of inv(A) */
117c4762a1bSJed Brown     /* -------------------------------------------------------------------- */
1185f80ce2aSJacob Faibussowitsch     CHKERRQ(MatZeroEntries(RHS));
119c4762a1bSJed Brown     for (i=0; i<nrhs; i++) {
120c4762a1bSJed Brown       v = 1.0;
1215f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValues(RHS,1,&i,1,&i,&v,INSERT_VALUES));
122c4762a1bSJed Brown     }
1235f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyBegin(RHS,MAT_FINAL_ASSEMBLY));
1245f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyEnd(RHS,MAT_FINAL_ASSEMBLY));
125c4762a1bSJed Brown 
1265f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMatSolve(F,RHS,X));
127c4762a1bSJed Brown     if (displ) {
1285f80ce2aSJacob Faibussowitsch       if (rank == 0) CHKERRQ(PetscPrintf(PETSC_COMM_SELF," \n(2) first %" PetscInt_FMT " columns of inv(A) with dense RHS:\n",nrhs));
1295f80ce2aSJacob Faibussowitsch       CHKERRQ(MatView(X,PETSC_VIEWER_STDOUT_WORLD));
130c4762a1bSJed Brown     }
131c4762a1bSJed Brown 
132c4762a1bSJed Brown     /* Check the residual */
1335f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMatMult(A,X,MAT_INITIAL_MATRIX,2.0,&AX));
1345f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAXPY(AX,-1.0,RHS,SAME_NONZERO_PATTERN));
1355f80ce2aSJacob Faibussowitsch     CHKERRQ(MatNorm(AX,NORM_INFINITY,&norm));
136c4762a1bSJed Brown     if (norm > tol) {
1375f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"(2) MatMatSolve: Norm of residual %g\n",norm));
138c4762a1bSJed Brown     }
1395f80ce2aSJacob Faibussowitsch     CHKERRQ(MatZeroEntries(X));
140c4762a1bSJed Brown 
141c4762a1bSJed Brown     /* (3) Test MatMatTransposeSolve() for inv(A) with sparse RHS stored in the host:
142c4762a1bSJed Brown      spRHST = [e[0],...,e[nrhs-1]]^T, dense X holds first nrhs columns of inv(A) */
143c4762a1bSJed Brown     /* --------------------------------------------------------------------------- */
144c4762a1bSJed Brown     /* Create spRHST: PETSc does not support compressed column format which is required by MUMPS for sparse RHS matrix,
145c4762a1bSJed Brown      thus user must create spRHST=spRHS^T and call MatMatTransposeSolve() */
1465f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreate(PETSC_COMM_WORLD,&spRHST));
147dd400576SPatrick Sanan     if (rank == 0) {
148c4762a1bSJed Brown       /* MUMPS requirs RHS be centralized on the host! */
1495f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetSizes(spRHST,nrhs,M,PETSC_DECIDE,PETSC_DECIDE));
150c4762a1bSJed Brown     } else {
1515f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetSizes(spRHST,0,0,PETSC_DECIDE,PETSC_DECIDE));
152c4762a1bSJed Brown     }
1535f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetType(spRHST,MATAIJ));
1545f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetFromOptions(spRHST));
1555f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetUp(spRHST));
156dd400576SPatrick Sanan     if (rank == 0) {
157c4762a1bSJed Brown       v = 1.0;
158c4762a1bSJed Brown       for (i=0; i<nrhs; i++) {
1595f80ce2aSJacob Faibussowitsch         CHKERRQ(MatSetValues(spRHST,1,&i,1,&i,&v,INSERT_VALUES));
160c4762a1bSJed Brown       }
161c4762a1bSJed Brown     }
1625f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyBegin(spRHST,MAT_FINAL_ASSEMBLY));
1635f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyEnd(spRHST,MAT_FINAL_ASSEMBLY));
164c4762a1bSJed Brown 
1655f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMatTransposeSolve(F,spRHST,X));
166c4762a1bSJed Brown 
167c4762a1bSJed Brown     if (displ) {
1685f80ce2aSJacob Faibussowitsch       if (rank == 0) CHKERRQ(PetscPrintf(PETSC_COMM_SELF," \n(3) first %" PetscInt_FMT " columns of inv(A) with sparse RHS:\n",nrhs));
1695f80ce2aSJacob Faibussowitsch       CHKERRQ(MatView(X,PETSC_VIEWER_STDOUT_WORLD));
170c4762a1bSJed Brown     }
171c4762a1bSJed Brown 
172c4762a1bSJed Brown     /* Check the residual: R = A*X - RHS */
1735f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMatMult(A,X,MAT_REUSE_MATRIX,2.0,&AX));
174c4762a1bSJed Brown 
1755f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAXPY(AX,-1.0,RHS,SAME_NONZERO_PATTERN));
1765f80ce2aSJacob Faibussowitsch     CHKERRQ(MatNorm(AX,NORM_INFINITY,&norm));
177c4762a1bSJed Brown     if (norm > tol) {
1785f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"(3) MatMatSolve: Norm of residual %g\n",norm));
179c4762a1bSJed Brown     }
180c4762a1bSJed Brown 
181c4762a1bSJed Brown     /* (4) Test MatMatSolve() for inv(A) with selected entries:
182c4762a1bSJed Brown      input: spRHS gives selected indices; output: spRHS holds selected entries of inv(A) */
183c4762a1bSJed Brown     /* --------------------------------------------------------------------------------- */
184c4762a1bSJed Brown     if (nrhs == N) { /* mumps requires nrhs = n */
185c4762a1bSJed Brown       /* Create spRHS on proc[0] */
186c4762a1bSJed Brown       Mat spRHS = NULL;
187c4762a1bSJed Brown 
188c4762a1bSJed Brown       /* Create spRHS = spRHST^T. Two matrices share internal matrix data structure */
1895f80ce2aSJacob Faibussowitsch       CHKERRQ(MatCreateTranspose(spRHST,&spRHS));
1905f80ce2aSJacob Faibussowitsch       CHKERRQ(MatMumpsGetInverse(F,spRHS));
1915f80ce2aSJacob Faibussowitsch       CHKERRQ(MatDestroy(&spRHS));
192c4762a1bSJed Brown 
1935f80ce2aSJacob Faibussowitsch       CHKERRQ(MatMumpsGetInverseTranspose(F,spRHST));
194c4762a1bSJed Brown       if (displ) {
1955f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\nSelected entries of inv(A^T):\n"));
1965f80ce2aSJacob Faibussowitsch         CHKERRQ(MatView(spRHST,PETSC_VIEWER_STDOUT_WORLD));
197c4762a1bSJed Brown       }
1985f80ce2aSJacob Faibussowitsch       CHKERRQ(MatDestroy(&spRHS));
199c4762a1bSJed Brown     }
2005f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&AX));
2015f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&F));
2025f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&RHS));
2035f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&spRHST));
204c4762a1bSJed Brown   }
205c4762a1bSJed Brown 
206c4762a1bSJed Brown   /* Free data structures */
2075f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
2085f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C));
2095f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&X));
2105f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomDestroy(&rand));
211*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
212*b122ec5aSJacob Faibussowitsch   return 0;
213c4762a1bSJed Brown #endif
214c4762a1bSJed Brown }
215c4762a1bSJed Brown 
216c4762a1bSJed Brown /*TEST
217c4762a1bSJed Brown 
218c4762a1bSJed Brown    test:
219c4762a1bSJed Brown      requires: mumps double !complex
220c4762a1bSJed Brown 
221c4762a1bSJed Brown    test:
222c4762a1bSJed Brown      suffix: 2
223c4762a1bSJed Brown      requires: mumps double !complex
224c4762a1bSJed Brown      nsize: 2
225c4762a1bSJed Brown 
226c4762a1bSJed Brown TEST*/
227