Actual source code: ex1.c
2: /* Program usage: mpirun 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: petsc.h - base PETSc routines petscvec.h - vectors
15: petscsys.h - system routines 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];
37: PetscInitialize(&argc,&args,(char *)0,help);
38: MPI_Comm_size(PETSC_COMM_WORLD,&size);
39: if (size != 1) SETERRQ(1,"This is a uniprocessor example only!");
40: PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
42: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
43: Compute the matrix and right-hand-side vector that define
44: the linear system, Ax = b.
45: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
47: /*
48: Create vectors. Note that we form 1 vector from scratch and
49: then duplicate as needed.
50: */
51: VecCreate(PETSC_COMM_WORLD,&x);
52: PetscObjectSetName((PetscObject) x, "Solution");
53: VecSetSizes(x,PETSC_DECIDE,n);
54: VecSetFromOptions(x);
55: VecDuplicate(x,&b);
56: VecDuplicate(x,&u);
58: /*
59: Create matrix. When using MatCreate(), the matrix format can
60: be specified at runtime.
62: Performance tuning note: For problems of substantial size,
63: preallocation of matrix memory is crucial for attaining good
64: performance. Since preallocation is not possible via the generic
65: matrix creation routine MatCreate(), we recommend for practical
66: problems instead to use the creation routine for a particular matrix
67: format, e.g.,
68: MatCreateSeqAIJ() - sequential AIJ (compressed sparse row)
69: MatCreateSeqBAIJ() - block AIJ
70: See the matrix chapter of the users manual for details.
71: */
72: MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n,&A);
73: MatSetFromOptions(A);
75: /*
76: Assemble matrix
77: */
78: value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
79: for (i=1; i<n-1; i++) {
80: col[0] = i-1; col[1] = i; col[2] = i+1;
81: MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
82: }
83: i = n - 1; col[0] = n - 2; col[1] = n - 1;
84: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
85: i = 0; col[0] = 0; col[1] = 1; value[0] = 2.0; value[1] = -1.0;
86: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
87: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
88: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
90: /*
91: Set exact solution; then compute right-hand-side vector.
92: */
93: VecSet(&one,u);
94: MatMult(A,u,b);
96: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
97: Create the linear solver and set various options
98: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
99: /*
100: Create linear solver context
101: */
102: KSPCreate(PETSC_COMM_WORLD,&ksp);
104: /*
105: Set operators. Here the matrix that defines the linear system
106: also serves as the preconditioning matrix.
107: */
108: KSPSetOperators(ksp,A,A,DIFFERENT_NONZERO_PATTERN);
110: /*
111: Set linear solver defaults for this problem (optional).
112: - By extracting the KSP and PC contexts from the KSP context,
113: we can then directly call any KSP and PC routines to set
114: various options.
115: - The following four statements are optional; all of these
116: parameters could alternatively be specified at runtime via
117: KSPSetFromOptions();
118: */
119: KSPGetPC(ksp,&pc);
120: PCSetType(pc,PCJACOBI);
121: KSPSetTolerances(ksp,1.e-7,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
123: /*
124: Set runtime options, e.g.,
125: -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
126: These options will override those specified above as long as
127: KSPSetFromOptions() is called _after_ any other customization
128: routines.
129: */
130: KSPSetFromOptions(ksp);
131:
132: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
133: Solve the linear system
134: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
135: /*
136: Solve linear system
137: */
138: KSPSolve(ksp,b,x);
140: /*
141: View solver info; we could instead use the option -ksp_view to
142: print this info to the screen at the conclusion of KSPSolve().
143: */
144: KSPView(ksp,PETSC_VIEWER_STDOUT_WORLD);
146: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
147: Check solution and clean up
148: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
149: /*
150: Check the error
151: */
152: VecAXPY(&neg_one,u,x);
153: VecNorm(x,NORM_2,&norm);
154: KSPGetIterationNumber(ksp,&its);
155: PetscPrintf(PETSC_COMM_WORLD,"Norm of error %A, Iterations %D\n",
156: norm,its);
157: /*
158: Free work space. All PETSc objects should be destroyed when they
159: are no longer needed.
160: */
161: VecDestroy(x); VecDestroy(u);
162: VecDestroy(b); MatDestroy(A);
163: KSPDestroy(ksp);
165: /*
166: Always call PetscFinalize() before exiting a program. This routine
167: - finalizes the PETSc libraries as well as MPI
168: - provides summary and diagnostic information if certain runtime
169: options are chosen (e.g., -log_summary).
170: */
171: PetscFinalize();
172: return 0;
173: }