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