Actual source code: ex48.c
petsc-3.13.4 2020-08-01
2: static char help[] = "Solves a tridiagonal linear system with KSP.\n\n";
4: /*
5: Test example that demonstrates how MINRES can produce a dp of zero
6: but still converge.
8: Provided by: Mark Filipiak <mjf@staffmail.ed.ac.uk>
9: */
10: #include <petscksp.h>
12: int main(int argc,char **args)
13: {
14: Vec x, b, u; /* approx solution, RHS, exact solution */
15: Mat A; /* linear system matrix */
16: KSP ksp; /* linear solver context */
17: PC pc; /* preconditioner context */
18: PetscReal norm;
20: PetscInt i,n = 2,col[3],its;
21: PetscMPIInt size;
22: PetscScalar one = 1.0,value[3];
23: PetscBool nonzeroguess = PETSC_FALSE;
25: PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
26: MPI_Comm_size(PETSC_COMM_WORLD,&size);
27: if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
28: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
29: PetscOptionsGetBool(NULL,NULL,"-nonzero_guess",&nonzeroguess,NULL);
32: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
33: Compute the matrix and right-hand-side vector that define
34: the linear system, Ax = b.
35: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
37: /*
38: Create vectors. Note that we form 1 vector from scratch and
39: then duplicate as needed.
40: */
41: VecCreate(PETSC_COMM_WORLD,&x);
42: PetscObjectSetName((PetscObject) x, "Solution");
43: VecSetSizes(x,PETSC_DECIDE,n);
44: VecSetFromOptions(x);
45: VecDuplicate(x,&b);
46: VecDuplicate(x,&u);
48: /*
49: Create matrix. When using MatCreate(), the matrix format can
50: be specified at runtime.
52: Performance tuning note: For problems of substantial size,
53: preallocation of matrix memory is crucial for attaining good
54: performance. See the matrix chapter of the users manual for details.
55: */
56: MatCreate(PETSC_COMM_WORLD,&A);
57: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
58: MatSetFromOptions(A);
59: MatSetUp(A);
61: /*
62: Assemble matrix
63: */
64: value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
65: for (i=1; i<n-1; i++) {
66: col[0] = i-1; col[1] = i; col[2] = i+1;
67: MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
68: }
69: i = n - 1; col[0] = n - 2; col[1] = n - 1;
70: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
71: i = 0; col[0] = 0; col[1] = 1; value[0] = 2.0; value[1] = -1.0;
72: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
73: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
74: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
76: /*
77: Set constant right-hand-side vector.
78: */
79: VecSet(b,one);
80: /*
81: Solution = RHS for the matrix [[2 -1] [-1 2]] and RHS [1 1]
82: */
83: VecSet(u,one);
85: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
86: Create the linear solver and set various options
87: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
88: /*
89: Create linear solver context
90: */
91: KSPCreate(PETSC_COMM_WORLD,&ksp);
93: /*
94: Set operators. Here the matrix that defines the linear system
95: also serves as the preconditioning matrix.
96: */
97: KSPSetOperators(ksp,A,A);
99: /*
100: Set linear solver defaults for this problem (optional).
101: - By extracting the KSP and PC contexts from the KSP context,
102: we can then directly call any KSP and PC routines to set
103: various options.
104: - The following four statements are optional; all of these
105: parameters could alternatively be specified at runtime via
106: KSPSetFromOptions();
107: */
108: KSPGetPC(ksp,&pc);
109: PCSetType(pc,PCJACOBI);
110: KSPSetTolerances(ksp,1.e-5,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
112: /*
113: Set runtime options, e.g.,
114: -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
115: These options will override those specified above as long as
116: KSPSetFromOptions() is called _after_ any other customization
117: routines.
118: */
119: KSPSetFromOptions(ksp);
121: if (nonzeroguess) {
122: PetscScalar p = .5;
123: VecSet(x,p);
124: KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);
125: }
127: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
128: Solve the linear system
129: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
130: /*
131: Solve linear system
132: */
133: KSPSolve(ksp,b,x);
135: /*
136: View solver info; we could instead use the option -ksp_view to
137: print this info to the screen at the conclusion of KSPSolve().
138: */
139: KSPView(ksp,PETSC_VIEWER_STDOUT_WORLD);
141: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
142: Check solution and clean up
143: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
144: /*
145: Check the error
146: */
147: VecAXPY(x,-1.0,u);
148: VecNorm(x,NORM_2,&norm);
149: KSPGetIterationNumber(ksp,&its);
150: PetscPrintf(PETSC_COMM_WORLD,"Norm of error %g, Iterations %D\n",(double)norm,its);
153: /*
154: Free work space. All PETSc objects should be destroyed when they
155: are no longer needed.
156: */
157: VecDestroy(&x); VecDestroy(&u);
158: VecDestroy(&b); MatDestroy(&A);
159: KSPDestroy(&ksp);
161: /*
162: Always call PetscFinalize() before exiting a program. This routine
163: - finalizes the PETSc libraries as well as MPI
164: - provides summary and diagnostic information if certain runtime
165: options are chosen (e.g., -log_view).
166: */
167: PetscFinalize();
168: return ierr;
169: }
172: /*TEST
174: test:
175: args: -ksp_monitor_short -ksp_converged_reason -ksp_type minres -pc_type jacobi -ksp_error_if_not_converged
177: TEST*/