Actual source code: ex53.c
petsc-3.4.5 2014-06-29
2: static char help[] = "Solves a tridiagonal linear system with KSP. \n\
3: Modified from ex1.c to illustrate reuse of preconditioner \n\
4: Written as requested by [petsc-maint #63875] \n\n";
6: #include <petscksp.h>
10: int main(int argc,char **args)
11: {
12: Vec x,x2,b,u; /* approx solution, RHS, exact solution */
13: Mat A; /* linear system matrix */
14: KSP ksp; /* linear solver context */
15: PC pc; /* preconditioner context */
16: PetscReal norm,tol=1.e-14; /* norm of solution error */
18: PetscInt i,n = 10,col[3],its;
19: PetscMPIInt rank;
20: PetscScalar neg_one = -1.0,one = 1.0,value[3];
22: PetscInitialize(&argc,&args,(char*)0,help);
23: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
24: PetscOptionsGetInt(NULL,"-n",&n,NULL);
26: /* Create vectors.*/
27: VecCreate(PETSC_COMM_WORLD,&x);
28: PetscObjectSetName((PetscObject) x, "Solution");
29: VecSetSizes(x,PETSC_DECIDE,n);
30: VecSetFromOptions(x);
31: VecDuplicate(x,&b);
32: VecDuplicate(x,&u);
33: VecDuplicate(x,&x2);
35: /* Create matrix. Only proc[0] sets values - not efficient for parallel processing!
36: See ex23.c for efficient parallel assembly matrix */
37: MatCreate(PETSC_COMM_WORLD,&A);
38: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
39: MatSetFromOptions(A);
40: MatSetUp(A);
42: if (!rank) {
43: value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
44: for (i=1; i<n-1; i++) {
45: col[0] = i-1; col[1] = i; col[2] = i+1;
46: MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
47: }
48: i = n - 1; col[0] = n - 2; col[1] = n - 1;
49: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
50: i = 0; col[0] = 0; col[1] = 1; value[0] = 2.0; value[1] = -1.0;
51: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
53: i = 0; col[0] = n-1; value[0] = 0.5; /* make A non-symmetric */
54: MatSetValues(A,1,&i,1,col,value,INSERT_VALUES);
55: }
56: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
57: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
59: /* Set exact solution */
60: VecSet(u,one);
62: /* Create linear solver context */
63: KSPCreate(PETSC_COMM_WORLD,&ksp);
64: KSPSetOperators(ksp,A,A,SAME_PRECONDITIONER);
65: KSPGetPC(ksp,&pc);
66: PCSetType(pc,PCLU);
67: #if defined(PETSC_HAVE_MUMPS)
68: PCFactorSetMatSolverPackage(pc,MATSOLVERMUMPS);
69: #endif
70: KSPSetFromOptions(ksp);
72: /* 1. Solve linear system A x = b */
73: MatMult(A,u,b);
74: KSPSolve(ksp,b,x);
76: /* Check the error */
77: VecAXPY(x,neg_one,u);
78: VecNorm(x,NORM_2,&norm);
79: KSPGetIterationNumber(ksp,&its);
80: if (norm > tol) {
81: PetscPrintf(PETSC_COMM_WORLD,"1. Norm of error for Ax=b: %G, Iterations %D\n",
82: norm,its);
83: }
85: /* 2. Solve linear system A^T x = b*/
86: MatMultTranspose(A,u,b);
87: KSPSolveTranspose(ksp,b,x2);
89: /* Check the error */
90: VecAXPY(x2,neg_one,u);
91: VecNorm(x2,NORM_2,&norm);
92: KSPGetIterationNumber(ksp,&its);
93: if (norm > tol) {
94: PetscPrintf(PETSC_COMM_WORLD,"2. Norm of error for A^T x=b: %G, Iterations %D\n",
95: norm,its);
96: }
98: /* 3. Change A and solve A x = b with an iterative solver using A=LU as a preconditioner*/
99: if (!rank) {
100: i = 0; col[0] = n-1; value[0] = 1.e-2;
101: MatSetValues(A,1,&i,1,col,value,ADD_VALUES);
102: }
103: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
104: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
106: MatMult(A,u,b);
107: KSPSolve(ksp,b,x);
109: /* Check the error */
110: VecAXPY(x,neg_one,u);
111: VecNorm(x,NORM_2,&norm);
112: KSPGetIterationNumber(ksp,&its);
113: if (norm > tol) {
114: PetscPrintf(PETSC_COMM_WORLD,"3. Norm of error for (A+Delta) x=b: %G, Iterations %D\n",norm,its);
115: }
117: /* Free work space. */
118: VecDestroy(&x);
119: VecDestroy(&u);
120: VecDestroy(&x2);
121: VecDestroy(&b);
122: MatDestroy(&A);
123: KSPDestroy(&ksp);
125: PetscFinalize();
126: return 0;
127: }