Actual source code: ex58.c
petsc-3.9.4 2018-09-11
2: static char help[] = "Solves a tridiagonal linear system with KSP.\n\n";
4: /*T
5: Concepts: KSP^solving a system of linear equations
6: Processors: 1
7: T*/
9: /*
10: Modified from ex1.c for testing matrix operations when matrix structure is changed.
11: Contributed by Jose E. Roman, Feb. 2012.
12: */
13: #include <petscksp.h>
15: int main(int argc,char **args)
16: {
17: Vec x, b, u; /* approx solution, RHS, exact solution */
18: Mat A,B,C; /* linear system matrix */
19: KSP ksp; /* linear solver context */
20: PC pc; /* preconditioner context */
21: PetscReal norm; /* norm of solution error */
23: PetscInt i,n = 20,col[3],its;
24: PetscMPIInt size;
25: PetscScalar one = 1.0,value[3];
26: PetscBool nonzeroguess = PETSC_FALSE;
28: PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
29: MPI_Comm_size(PETSC_COMM_WORLD,&size);
30: if (size != 1) SETERRQ(PETSC_COMM_WORLD,1,"This is a uniprocessor example only!");
31: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
32: PetscOptionsGetBool(NULL,NULL,"-nonzero_guess",&nonzeroguess,NULL);
35: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
36: Compute the matrix and right-hand-side vector that define
37: the linear system, Ax = b.
38: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
40: /*
41: Create vectors. Note that we form 1 vector from scratch and
42: then duplicate as needed.
43: */
44: VecCreate(PETSC_COMM_WORLD,&x);
45: PetscObjectSetName((PetscObject) x, "Solution");
46: VecSetSizes(x,PETSC_DECIDE,n);
47: VecSetFromOptions(x);
48: VecDuplicate(x,&b);
49: VecDuplicate(x,&u);
51: /*
52: Create matrix. When using MatCreate(), the matrix format can
53: be specified at runtime.
55: Performance tuning note: For problems of substantial size,
56: preallocation of matrix memory is crucial for attaining good
57: performance. See the matrix chapter of the users manual for details.
58: */
59: MatCreate(PETSC_COMM_WORLD,&A);
60: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
61: MatSetFromOptions(A);
62: MatSetUp(A);
64: /*
65: Assemble matrix
66: */
67: value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
68: for (i=1; i<n-1; i++) {
69: col[0] = i-1; col[1] = i; col[2] = i+1;
70: MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
71: }
72: i = n - 1; col[0] = n - 2; col[1] = n - 1;
73: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
74: i = 0; col[0] = 0; col[1] = 1; value[0] = 2.0; value[1] = -1.0;
75: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
76: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
77: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
79: /*
80: jroman: added matrices
81: */
82: MatCreate(PETSC_COMM_WORLD,&B);
83: MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n);
84: MatSetFromOptions(B);
85: MatSetUp(B);
86: for (i=0; i<n; i++) {
87: MatSetValue(B,i,i,value[1],INSERT_VALUES);
88: if (n-i+n/3<n) {
89: MatSetValue(B,n-i+n/3,i,value[0],INSERT_VALUES);
90: MatSetValue(B,i,n-i+n/3,value[0],INSERT_VALUES);
91: }
92: }
93: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
94: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
95: MatDuplicate(A,MAT_COPY_VALUES,&C);
96: MatAXPY(C,2.0,B,DIFFERENT_NONZERO_PATTERN);
98: /*
99: Set exact solution; then compute right-hand-side vector.
100: */
101: VecSet(u,one);
102: MatMult(C,u,b);
104: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
105: Create the linear solver and set various options
106: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
107: /*
108: Create linear solver context
109: */
110: KSPCreate(PETSC_COMM_WORLD,&ksp);
112: /*
113: Set operators. Here the matrix that defines the linear system
114: also serves as the preconditioning matrix.
115: */
116: KSPSetOperators(ksp,C,C);
118: /*
119: Set linear solver defaults for this problem (optional).
120: - By extracting the KSP and PC contexts from the KSP context,
121: we can then directly call any KSP and PC routines to set
122: various options.
123: - The following four statements are optional; all of these
124: parameters could alternatively be specified at runtime via
125: KSPSetFromOptions();
126: */
127: KSPGetPC(ksp,&pc);
128: PCSetType(pc,PCJACOBI);
129: KSPSetTolerances(ksp,1.e-5,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
131: /*
132: Set runtime options, e.g.,
133: -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
134: These options will override those specified above as long as
135: KSPSetFromOptions() is called _after_ any other customization
136: routines.
137: */
138: KSPSetFromOptions(ksp);
140: if (nonzeroguess) {
141: PetscScalar p = .5;
142: VecSet(x,p);
143: KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);
144: }
146: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
147: Solve the linear system
148: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
149: /*
150: Solve linear system
151: */
152: KSPSolve(ksp,b,x);
154: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
155: Check solution and clean up
156: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
157: /*
158: Check the error
159: */
160: VecAXPY(x,-1.0,u);
161: VecNorm(x,NORM_2,&norm);
162: KSPGetIterationNumber(ksp,&its);
163: PetscPrintf(PETSC_COMM_WORLD,"Norm of error %g, Iterations %D\n",(double)norm,its);
165: /*
166: Free work space. All PETSc objects should be destroyed when they
167: are no longer needed.
168: */
169: VecDestroy(&x); VecDestroy(&u);
170: VecDestroy(&b); MatDestroy(&A);
171: MatDestroy(&B);
172: MatDestroy(&C);
173: KSPDestroy(&ksp);
175: /*
176: Always call PetscFinalize() before exiting a program. This routine
177: - finalizes the PETSc libraries as well as MPI
178: - provides summary and diagnostic information if certain runtime
179: options are chosen (e.g., -log_view).
180: */
181: PetscFinalize();
182: return ierr;
183: }
186: /*TEST
188: test:
189: args: -mat_type aij
190: output_file: output/ex58.out
192: test:
193: suffix: baij
194: args: -mat_type baij
195: output_file: output/ex58.out
197: test:
198: suffix: sbaij
199: args: -mat_type sbaij
200: output_file: output/ex58.out
202: TEST*/