Actual source code: ex53.c
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>
8: int main(int argc,char **args)
9: {
10: Vec x,x2,b,u; /* approx solution, RHS, exact solution */
11: Mat A; /* linear system matrix */
12: KSP ksp; /* linear solver context */
13: PC pc; /* preconditioner context */
14: PetscReal norm,tol=100.*PETSC_MACHINE_EPSILON; /* norm of solution error */
15: PetscInt i,n = 10,col[3],its;
16: PetscMPIInt rank,size;
17: PetscScalar one = 1.0,value[3];
19: PetscInitialize(&argc,&args,(char*)0,help);
20: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
21: MPI_Comm_size(PETSC_COMM_WORLD,&size);
22: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
24: /* Create vectors.*/
25: VecCreate(PETSC_COMM_WORLD,&x);
26: PetscObjectSetName((PetscObject) x, "Solution");
27: VecSetSizes(x,PETSC_DECIDE,n);
28: VecSetFromOptions(x);
29: VecDuplicate(x,&b);
30: VecDuplicate(x,&u);
31: VecDuplicate(x,&x2);
33: /* Create matrix. Only proc[0] sets values - not efficient for parallel processing!
34: See ex23.c for efficient parallel assembly matrix */
35: MatCreate(PETSC_COMM_WORLD,&A);
36: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
37: MatSetFromOptions(A);
38: MatSetUp(A);
40: if (rank == 0) {
41: value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
42: for (i=1; i<n-1; i++) {
43: col[0] = i-1; col[1] = i; col[2] = i+1;
44: MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
45: }
46: i = n - 1; col[0] = n - 2; col[1] = n - 1;
47: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
48: i = 0; col[0] = 0; col[1] = 1; value[0] = 2.0; value[1] = -1.0;
49: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
51: i = 0; col[0] = n-1; value[0] = 0.5; /* make A non-symmetric */
52: MatSetValues(A,1,&i,1,col,value,INSERT_VALUES);
53: }
54: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
55: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
57: /* Set exact solution */
58: VecSet(u,one);
60: /* Create linear solver context */
61: KSPCreate(PETSC_COMM_WORLD,&ksp);
62: KSPSetOperators(ksp,A,A);
63: KSPGetPC(ksp,&pc);
64: PCSetType(pc,PCLU);
65: #if defined(PETSC_HAVE_MUMPS)
66: if (size > 1) {
67: PCFactorSetMatSolverType(pc,MATSOLVERMUMPS);
68: }
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,-1.0,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",(double)norm,its);
82: }
84: /* 2. Solve linear system A^T x = b*/
85: MatMultTranspose(A,u,b);
86: KSPSolveTranspose(ksp,b,x2);
88: /* Check the error */
89: VecAXPY(x2,-1.0,u);
90: VecNorm(x2,NORM_2,&norm);
91: KSPGetIterationNumber(ksp,&its);
92: if (norm > tol) {
93: PetscPrintf(PETSC_COMM_WORLD,"2. Norm of error for A^T x=b: %g, Iterations %D\n",(double)norm,its);
94: }
96: /* 3. Change A and solve A x = b with an iterative solver using A=LU as a preconditioner*/
97: if (rank == 0) {
98: i = 0; col[0] = n-1; value[0] = 1.e-2;
99: MatSetValues(A,1,&i,1,col,value,ADD_VALUES);
100: }
101: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
102: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
104: MatMult(A,u,b);
105: KSPSolve(ksp,b,x);
107: /* Check the error */
108: VecAXPY(x,-1.0,u);
109: VecNorm(x,NORM_2,&norm);
110: KSPGetIterationNumber(ksp,&its);
111: if (norm > tol) {
112: PetscPrintf(PETSC_COMM_WORLD,"3. Norm of error for (A+Delta) x=b: %g, Iterations %D\n",(double)norm,its);
113: }
115: /* Free work space. */
116: VecDestroy(&x);
117: VecDestroy(&u);
118: VecDestroy(&x2);
119: VecDestroy(&b);
120: MatDestroy(&A);
121: KSPDestroy(&ksp);
123: PetscFinalize();
124: return 0;
125: }
127: /*TEST
129: test:
130: requires: mumps
132: test:
133: suffix: 2
134: nsize: 2
135: requires: mumps
136: output_file: output/ex53.out
138: TEST*/