Actual source code: ex53.c
petsc-3.7.5 2017-01-01
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=100.*PETSC_MACHINE_EPSILON; /* norm of solution error */
18: PetscInt i,n = 10,col[3],its;
19: PetscMPIInt rank,size;
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: MPI_Comm_size(PETSC_COMM_WORLD,&size);
25: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
27: /* Create vectors.*/
28: VecCreate(PETSC_COMM_WORLD,&x);
29: PetscObjectSetName((PetscObject) x, "Solution");
30: VecSetSizes(x,PETSC_DECIDE,n);
31: VecSetFromOptions(x);
32: VecDuplicate(x,&b);
33: VecDuplicate(x,&u);
34: VecDuplicate(x,&x2);
36: /* Create matrix. Only proc[0] sets values - not efficient for parallel processing!
37: See ex23.c for efficient parallel assembly matrix */
38: MatCreate(PETSC_COMM_WORLD,&A);
39: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
40: MatSetFromOptions(A);
41: MatSetUp(A);
43: if (!rank) {
44: value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
45: for (i=1; i<n-1; i++) {
46: col[0] = i-1; col[1] = i; col[2] = i+1;
47: MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
48: }
49: i = n - 1; col[0] = n - 2; col[1] = n - 1;
50: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
51: i = 0; col[0] = 0; col[1] = 1; value[0] = 2.0; value[1] = -1.0;
52: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
54: i = 0; col[0] = n-1; value[0] = 0.5; /* make A non-symmetric */
55: MatSetValues(A,1,&i,1,col,value,INSERT_VALUES);
56: }
57: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
58: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
60: /* Set exact solution */
61: VecSet(u,one);
63: /* Create linear solver context */
64: KSPCreate(PETSC_COMM_WORLD,&ksp);
65: KSPSetOperators(ksp,A,A);
66: KSPGetPC(ksp,&pc);
67: PCSetType(pc,PCLU);
68: #if defined(PETSC_HAVE_MUMPS)
69: if (size > 1) {
70: PCFactorSetMatSolverPackage(pc,MATSOLVERMUMPS);
71: }
72: #endif
73: KSPSetFromOptions(ksp);
75: /* 1. Solve linear system A x = b */
76: MatMult(A,u,b);
77: KSPSolve(ksp,b,x);
79: /* Check the error */
80: VecAXPY(x,neg_one,u);
81: VecNorm(x,NORM_2,&norm);
82: KSPGetIterationNumber(ksp,&its);
83: if (norm > tol) {
84: PetscPrintf(PETSC_COMM_WORLD,"1. Norm of error for Ax=b: %g, Iterations %D\n",(double)norm,its);
85: }
87: /* 2. Solve linear system A^T x = b*/
88: MatMultTranspose(A,u,b);
89: KSPSolveTranspose(ksp,b,x2);
91: /* Check the error */
92: VecAXPY(x2,neg_one,u);
93: VecNorm(x2,NORM_2,&norm);
94: KSPGetIterationNumber(ksp,&its);
95: if (norm > tol) {
96: PetscPrintf(PETSC_COMM_WORLD,"2. Norm of error for A^T x=b: %g, Iterations %D\n",(double)norm,its);
97: }
99: /* 3. Change A and solve A x = b with an iterative solver using A=LU as a preconditioner*/
100: if (!rank) {
101: i = 0; col[0] = n-1; value[0] = 1.e-2;
102: MatSetValues(A,1,&i,1,col,value,ADD_VALUES);
103: }
104: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
105: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
107: MatMult(A,u,b);
108: KSPSolve(ksp,b,x);
110: /* Check the error */
111: VecAXPY(x,neg_one,u);
112: VecNorm(x,NORM_2,&norm);
113: KSPGetIterationNumber(ksp,&its);
114: if (norm > tol) {
115: PetscPrintf(PETSC_COMM_WORLD,"3. Norm of error for (A+Delta) x=b: %g, Iterations %D\n",(double)norm,its);
116: }
118: /* Free work space. */
119: VecDestroy(&x);
120: VecDestroy(&u);
121: VecDestroy(&x2);
122: VecDestroy(&b);
123: MatDestroy(&A);
124: KSPDestroy(&ksp);
126: PetscFinalize();
127: return 0;
128: }