Actual source code: ex1.c
2: /* Program usage: mpiexec ex1 [-help] [all PETSc options] */
4: static char help[] = "Solves a tridiagonal linear system with KSP.\n\n";
6: /*T
7: Concepts: KSP^solving a system of linear equations
8: Processors: 1
9: T*/
11: /*
12: Include "petscksp.h" so that we can use KSP solvers. Note that this file
13: automatically includes:
14: petscsys.h - base PETSc routines petscvec.h - vectors
15: petscmat.h - matrices
16: petscis.h - index sets petscksp.h - Krylov subspace methods
17: petscviewer.h - viewers petscpc.h - preconditioners
19: Note: The corresponding parallel example is ex23.c
20: */
21: #include <petscksp.h>
25: int main(int argc,char **args)
26: {
27: Vec x, b, u; /* approx solution, RHS, exact solution */
28: Mat A; /* linear system matrix */
29: KSP ksp; /* linear solver context */
30: PC pc; /* preconditioner context */
31: PetscReal norm; /* norm of solution error */
33: PetscInt i,n = 10,col[3],its;
34: PetscMPIInt size;
35: PetscScalar neg_one = -1.0,one = 1.0,value[3];
36: PetscBool nonzeroguess = PETSC_FALSE;
38: PetscInitialize(&argc,&args,(char *)0,help);
39: MPI_Comm_size(PETSC_COMM_WORLD,&size);
40: if (size != 1) SETERRQ(PETSC_COMM_WORLD,1,"This is a uniprocessor example only!");
41: PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
42: PetscOptionsGetBool(PETSC_NULL,"-nonzero_guess",&nonzeroguess,PETSC_NULL);
45: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
46: Compute the matrix and right-hand-side vector that define
47: the linear system, Ax = b.
48: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
50: /*
51: Create vectors. Note that we form 1 vector from scratch and
52: then duplicate as needed.
53: */
54: VecCreate(PETSC_COMM_WORLD,&x);
55: PetscObjectSetName((PetscObject) x, "Solution");
56: VecSetSizes(x,PETSC_DECIDE,n);
57: VecSetFromOptions(x);
58: VecDuplicate(x,&b);
59: VecDuplicate(x,&u);
61: /*
62: Create matrix. When using MatCreate(), the matrix format can
63: be specified at runtime.
65: Performance tuning note: For problems of substantial size,
66: preallocation of matrix memory is crucial for attaining good
67: performance. See the matrix chapter of the users manual for details.
68: */
69: MatCreate(PETSC_COMM_WORLD,&A);
70: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
71: MatSetFromOptions(A);
73: /*
74: Assemble matrix
75: */
76: value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
77: for (i=1; i<n-1; i++) {
78: col[0] = i-1; col[1] = i; col[2] = i+1;
79: MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
80: }
81: i = n - 1; col[0] = n - 2; col[1] = n - 1;
82: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
83: i = 0; col[0] = 0; col[1] = 1; value[0] = 2.0; value[1] = -1.0;
84: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
85: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
86: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
88: /*
89: Set exact solution; then compute right-hand-side vector.
90: */
91: VecSet(u,one);
92: MatMult(A,u,b);
94: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
95: Create the linear solver and set various options
96: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
97: /*
98: Create linear solver context
99: */
100: KSPCreate(PETSC_COMM_WORLD,&ksp);
102: /*
103: Set operators. Here the matrix that defines the linear system
104: also serves as the preconditioning matrix.
105: */
106: KSPSetOperators(ksp,A,A,DIFFERENT_NONZERO_PATTERN);
108: /*
109: Set linear solver defaults for this problem (optional).
110: - By extracting the KSP and PC contexts from the KSP context,
111: we can then directly call any KSP and PC routines to set
112: various options.
113: - The following four statements are optional; all of these
114: parameters could alternatively be specified at runtime via
115: KSPSetFromOptions();
116: */
117: KSPGetPC(ksp,&pc);
118: PCSetType(pc,PCJACOBI);
119: KSPSetTolerances(ksp,1.e-5,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
121: /*
122: Set runtime options, e.g.,
123: -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
124: These options will override those specified above as long as
125: KSPSetFromOptions() is called _after_ any other customization
126: routines.
127: */
128: KSPSetFromOptions(ksp);
130: if (nonzeroguess) {
131: PetscScalar p = .5;
132: VecSet(x,p);
133: KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);
134: }
135:
136: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
137: Solve the linear system
138: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
139: /*
140: Solve linear system
141: */
142: KSPSolve(ksp,b,x);
144: /*
145: View solver info; we could instead use the option -ksp_view to
146: print this info to the screen at the conclusion of KSPSolve().
147: */
148: KSPView(ksp,PETSC_VIEWER_STDOUT_WORLD);
150: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
151: Check solution and clean up
152: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
153: /*
154: Check the error
155: */
156: VecAXPY(x,neg_one,u);
157: VecNorm(x,NORM_2,&norm);
158: KSPGetIterationNumber(ksp,&its);
159: PetscPrintf(PETSC_COMM_WORLD,"Norm of error %A, Iterations %D\n",
160: norm,its);
162: /*
163: Free work space. All PETSc objects should be destroyed when they
164: are no longer needed.
165: */
166: VecDestroy(&x); VecDestroy(&u);
167: VecDestroy(&b); MatDestroy(&A);
168: KSPDestroy(&ksp);
170: /*
171: Always call PetscFinalize() before exiting a program. This routine
172: - finalizes the PETSc libraries as well as MPI
173: - provides summary and diagnostic information if certain runtime
174: options are chosen (e.g., -log_summary).
175: */
176: PetscFinalize();
177: return 0;
178: }