Actual source code: ex2.c

petsc-3.7.2 2016-06-05
Report Typos and Errors
  2: static char help[] = "Tests repeated solving linear system on 2 by 2 matrix provided by MUMPS developer, Dec 17, 2012.\n\n";
  3: /*
  4: We have investigated the problem further, and we have
  5: been able to reproduce it and obtain an erroneous
  6: solution with an even smaller, 2x2, matrix:
  7:     [1 2]
  8:     [2 3]
  9: and a right-hand side vector with all ones (1,1)
 10: The correct solution is the vector (-1,1), in both solves.

 12: mpiexec -n 2 ./ex2 -ksp_type preonly -pc_type lu  -pc_factor_mat_solver_package mumps  -mat_mumps_icntl_7 6 -mat_mumps_cntl_1 0.99

 14: With this combination of options, I get off-diagonal pivots during the
 15: factorization, which is the cause of the problem (different isol_loc
 16: returned in the second solve, whereas, as I understand it, Petsc expects
 17: isol_loc not to change between successive solves).
 18: */

 20: #include <petscksp.h>

 24: int main(int argc,char **args)
 25: {
 26:   Mat            C;
 28:   PetscInt       N = 2,rowidx,colidx;
 29:   Vec            u,b,r;
 30:   KSP            ksp;
 31:   PetscReal      norm;
 32:   PetscMPIInt    rank,size;
 33:   PetscScalar    v;

 35:   PetscInitialize(&argc,&args,(char*)0,help);
 36:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 37:   MPI_Comm_size(PETSC_COMM_WORLD,&size);

 39:   /* create stiffness matrix C = [1 2; 2 3] */
 40:   MatCreate(PETSC_COMM_WORLD,&C);
 41:   MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,N,N);
 42:   MatSetFromOptions(C);
 43:   MatSetUp(C);
 44:   if (!rank) {
 45:     rowidx = 0; colidx = 0; v = 1.0;
 46:     MatSetValues(C,1,&rowidx,1,&colidx,&v,INSERT_VALUES);
 47:     rowidx = 0; colidx = 1; v = 2.0;
 48:     MatSetValues(C,1,&rowidx,1,&colidx,&v,INSERT_VALUES);

 50:     rowidx = 1; colidx = 0; v = 2.0;
 51:     MatSetValues(C,1,&rowidx,1,&colidx,&v,INSERT_VALUES);
 52:     rowidx = 1; colidx = 1; v = 3.0;
 53:     MatSetValues(C,1,&rowidx,1,&colidx,&v,INSERT_VALUES);
 54:   }
 55:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
 56:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

 58:   /* create right hand side and solution */
 59:   VecCreate(PETSC_COMM_WORLD,&u);
 60:   VecSetSizes(u,PETSC_DECIDE,N);
 61:   VecSetFromOptions(u);
 62:   VecDuplicate(u,&b);
 63:   VecDuplicate(u,&r);
 64:   VecSet(u,0.0);
 65:   VecSet(b,1.0);

 67:   /* solve linear system C*u = b */
 68:   KSPCreate(PETSC_COMM_WORLD,&ksp);
 69:   KSPSetOperators(ksp,C,C);
 70:   KSPSetFromOptions(ksp);
 71:   KSPSolve(ksp,b,u);

 73:   /* check residual r = C*u - b */
 74:   MatMult(C,u,r);
 75:   VecAXPY(r,-1.0,b);
 76:   VecNorm(r,NORM_2,&norm);
 77:   PetscPrintf(PETSC_COMM_WORLD,"|| C*u - b|| = %g\n",(double)norm);

 79:   /* solve C^T*u = b twice */
 80:   KSPSolveTranspose(ksp,b,u);
 81:   /* check residual r = C^T*u - b */
 82:   MatMultTranspose(C,u,r);
 83:   VecAXPY(r,-1.0,b);
 84:   VecNorm(r,NORM_2,&norm);
 85:   PetscPrintf(PETSC_COMM_WORLD,"|| C^T*u - b|| =  %g\n",(double)norm);

 87:   KSPSolveTranspose(ksp,b,u);
 88:   MatMultTranspose(C,u,r);
 89:   VecAXPY(r,-1.0,b);
 90:   VecNorm(r,NORM_2,&norm);
 91:   PetscPrintf(PETSC_COMM_WORLD,"|| C^T*u - b|| =  %g\n",(double)norm);

 93:   /* solve C*u = b again */
 94:   KSPSolve(ksp,b,u);
 95:   MatMult(C,u,r);
 96:   VecAXPY(r,-1.0,b);
 97:   VecNorm(r,NORM_2,&norm);
 98:   PetscPrintf(PETSC_COMM_WORLD,"|| C*u - b|| = %g\n",(double)norm);

100:   KSPDestroy(&ksp);
101:   VecDestroy(&u);
102:   VecDestroy(&r);
103:   VecDestroy(&b);
104:   MatDestroy(&C);
105:   PetscFinalize();
106:   return 0;
107: }