Actual source code: ex1.c
1: /*$Id: ex1.c,v 1.90 2001/08/07 21:30:54 bsmith Exp $*/
3: /* Program usage: mpirun ex1 [-help] [all PETSc options] */
5: static char help[] = "Solves a tridiagonal linear system with SLES.nn";
7: /*T
8: Concepts: SLES^solving a system of linear equations
9: Processors: 1
10: T*/
12: /*
13: Include "petscsles.h" so that we can use SLES solvers. Note that this file
14: automatically includes:
15: petsc.h - base PETSc routines petscvec.h - vectors
16: petscsys.h - system routines petscmat.h - matrices
17: petscis.h - index sets petscksp.h - Krylov subspace methods
18: petscviewer.h - viewers petscpc.h - preconditioners
20: Note: The corresponding parallel example is ex23.c
21: */
22: #include petscsles.h
24: int main(int argc,char **args)
25: {
26: Vec x, b, u; /* approx solution, RHS, exact solution */
27: Mat A; /* linear system matrix */
28: SLES sles; /* linear solver context */
29: PC pc; /* preconditioner context */
30: KSP ksp; /* Krylov subspace method context */
31: PetscReal norm; /* norm of solution error */
32: int ierr,i,n = 10,col[3],its,size;
33: PetscScalar neg_one = -1.0,one = 1.0,value[3];
35: PetscInitialize(&argc,&args,(char *)0,help);
36: MPI_Comm_size(PETSC_COMM_WORLD,&size);
37: if (size != 1) SETERRQ(1,"This is a uniprocessor example only!");
38: PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
40: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
41: Compute the matrix and right-hand-side vector that define
42: the linear system, Ax = b.
43: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
45: /*
46: Create vectors. Note that we form 1 vector from scratch and
47: then duplicate as needed.
48: */
49: VecCreate(PETSC_COMM_WORLD,&x);
50: VecSetSizes(x,PETSC_DECIDE,n);
51: VecSetFromOptions(x);
52: VecDuplicate(x,&b);
53: VecDuplicate(x,&u);
55: /*
56: Create matrix. When using MatCreate(), the matrix format can
57: be specified at runtime.
59: Performance tuning note: For problems of substantial size,
60: preallocation of matrix memory is crucial for attaining good
61: performance. Since preallocation is not possible via the generic
62: matrix creation routine MatCreate(), we recommend for practical
63: problems instead to use the creation routine for a particular matrix
64: format, e.g.,
65: MatCreateSeqAIJ() - sequential AIJ (compressed sparse row)
66: MatCreateSeqBAIJ() - block AIJ
67: See the matrix chapter of the users manual for details.
68: */
69: MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n,&A);
70: MatSetFromOptions(A);
72: /*
73: Assemble matrix
74: */
75: value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
76: for (i=1; i<n-1; i++) {
77: col[0] = i-1; col[1] = i; col[2] = i+1;
78: MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
79: }
80: i = n - 1; col[0] = n - 2; col[1] = n - 1;
81: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
82: i = 0; col[0] = 0; col[1] = 1; value[0] = 2.0; value[1] = -1.0;
83: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
84: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
85: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
87: /*
88: Set exact solution; then compute right-hand-side vector.
89: */
90: VecSet(&one,u);
91: MatMult(A,u,b);
93: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
94: Create the linear solver and set various options
95: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
96: /*
97: Create linear solver context
98: */
99: SLESCreate(PETSC_COMM_WORLD,&sles);
101: /*
102: Set operators. Here the matrix that defines the linear system
103: also serves as the preconditioning matrix.
104: */
105: SLESSetOperators(sles,A,A,DIFFERENT_NONZERO_PATTERN);
107: /*
108: Set linear solver defaults for this problem (optional).
109: - By extracting the KSP and PC contexts from the SLES context,
110: we can then directly call any KSP and PC routines to set
111: various options.
112: - The following four statements are optional; all of these
113: parameters could alternatively be specified at runtime via
114: SLESSetFromOptions();
115: */
116: SLESGetKSP(sles,&ksp);
117: SLESGetPC(sles,&pc);
118: PCSetType(pc,PCJACOBI);
119: KSPSetTolerances(ksp,1.e-7,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: SLESSetFromOptions() is called _after_ any other customization
126: routines.
127: */
128: SLESSetFromOptions(sles);
129:
130: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
131: Solve the linear system
132: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
133: /*
134: Solve linear system
135: */
136: SLESSolve(sles,b,x,&its);
138: /*
139: View solver info; we could instead use the option -sles_view to
140: print this info to the screen at the conclusion of SLESSolve().
141: */
142: SLESView(sles,PETSC_VIEWER_STDOUT_WORLD);
144: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
145: Check solution and clean up
146: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
147: /*
148: Check the error
149: */
150: VecAXPY(&neg_one,u,x);
151: ierr = VecNorm(x,NORM_2,&norm);
152: PetscPrintf(PETSC_COMM_WORLD,"Norm of error %A, Iterations %dn",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: SLESDestroy(sles);
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_summary).
166: */
167: PetscFinalize();
168: return 0;
169: }