Actual source code: ex7.c
1: /*$Id: ex7.c,v 1.58 2001/08/07 21:30:54 bsmith Exp $*/
3: static char help[] = "Block Jacobi preconditioner for solving a linear system in parallel with SLES.n
4: The code indicates then
5: procedures for setting the particular block sizes and for using differentn
6: linear solvers on the individual blocks.nn";
8: /*
9: Note: This example focuses on ways to customize the block Jacobi
10: preconditioner. See ex1.c and ex2.c for more detailed comments on
11: the basic usage of SLES (including working with matrices and vectors).
13: Recall: The block Jacobi method is equivalent to the ASM preconditioner
14: with zero overlap.
15: */
17: /*T
18: Concepts: SLES^customizing the block Jacobi preconditioner
19: Processors: n
20: T*/
22: /*
23: Include "petscsles.h" so that we can use SLES solvers. Note that this file
24: automatically includes:
25: petsc.h - base PETSc routines petscvec.h - vectors
26: petscsys.h - system routines petscmat.h - matrices
27: petscis.h - index sets petscksp.h - Krylov subspace methods
28: petscviewer.h - viewers petscpc.h - preconditioners
29: */
30: #include petscsles.h
32: #undef __FUNCT__
34: int main(int argc,char **args)
35: {
36: Vec x,b,u; /* approx solution, RHS, exact solution */
37: Mat A; /* linear system matrix */
38: SLES sles; /* SLES context */
39: SLES *subsles; /* array of local SLES contexts on this processor */
40: PC pc; /* PC context */
41: PC subpc; /* PC context for subdomain */
42: KSP subksp; /* KSP context for subdomain */
43: PetscReal norm; /* norm of solution error */
44: int i,j,I,J,ierr,*blks,m = 8,n;
45: int rank,size,its,nlocal,first,Istart,Iend;
46: PetscScalar v,one = 1.0,none = -1.0;
47: PetscTruth isbjacobi,flg;
49: PetscInitialize(&argc,&args,(char *)0,help);
50: PetscOptionsGetInt(PETSC_NULL,"-m",&m,PETSC_NULL);
51: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
52: MPI_Comm_size(PETSC_COMM_WORLD,&size);
53: n = m+2;
55: /* -------------------------------------------------------------------
56: Compute the matrix and right-hand-side vector that define
57: the linear system, Ax = b.
58: ------------------------------------------------------------------- */
60: /*
61: Create and assemble parallel matrix
62: */
63: MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n,&A);
64: MatSetFromOptions(A);
65: MatGetOwnershipRange(A,&Istart,&Iend);
66: for (I=Istart; I<Iend; I++) {
67: v = -1.0; i = I/n; j = I - i*n;
68: if (i>0) {J = I - n; MatSetValues(A,1,&I,1,&J,&v,ADD_VALUES);}
69: if (i<m-1) {J = I + n; MatSetValues(A,1,&I,1,&J,&v,ADD_VALUES);}
70: if (j>0) {J = I - 1; MatSetValues(A,1,&I,1,&J,&v,ADD_VALUES);}
71: if (j<n-1) {J = I + 1; MatSetValues(A,1,&I,1,&J,&v,ADD_VALUES);}
72: v = 4.0; MatSetValues(A,1,&I,1,&I,&v,ADD_VALUES);
73: }
74: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
75: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
77: /*
78: Create parallel vectors
79: */
80: VecCreate(PETSC_COMM_WORLD,&u);
81: VecSetSizes(u,PETSC_DECIDE,m*n);
82: VecSetFromOptions(u);
83: VecDuplicate(u,&b);
84: VecDuplicate(b,&x);
86: /*
87: Set exact solution; then compute right-hand-side vector.
88: */
89: VecSet(&one,u);
90: MatMult(A,u,b);
92: /*
93: Create linear solver context
94: */
95: SLESCreate(PETSC_COMM_WORLD,&sles);
97: /*
98: Set operators. Here the matrix that defines the linear system
99: also serves as the preconditioning matrix.
100: */
101: SLESSetOperators(sles,A,A,DIFFERENT_NONZERO_PATTERN);
103: /*
104: Set default preconditioner for this program to be block Jacobi.
105: This choice can be overridden at runtime with the option
106: -pc_type <type>
107: */
108: SLESGetPC(sles,&pc);
109: PCSetType(pc,PCBJACOBI);
112: /* -------------------------------------------------------------------
113: Define the problem decomposition
114: ------------------------------------------------------------------- */
116: /*
117: Call PCBJacobiSetTotalBlocks() to set individually the size of
118: each block in the preconditioner. This could also be done with
119: the runtime option
120: -pc_bjacobi_blocks <blocks>
121: Also, see the command PCBJacobiSetLocalBlocks() to set the
122: local blocks.
124: Note: The default decomposition is 1 block per processor.
125: */
126: PetscMalloc(m*sizeof(int),&blks);
127: for (i=0; i<m; i++) blks[i] = n;
128: PCBJacobiSetTotalBlocks(pc,m,blks);
129: PetscFree(blks);
132: /* -------------------------------------------------------------------
133: Set the linear solvers for the subblocks
134: ------------------------------------------------------------------- */
136: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
137: Basic method, should be sufficient for the needs of most users.
138: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
140: By default, the block Jacobi method uses the same solver on each
141: block of the problem. To set the same solver options on all blocks,
142: use the prefix -sub before the usual PC and KSP options, e.g.,
143: -sub_pc_type <pc> -sub_ksp_type <ksp> -sub_ksp_rtol 1.e-4
144: */
146: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
147: Advanced method, setting different solvers for various blocks.
148: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
150: Note that each block's SLES context is completely independent of
151: the others, and the full range of uniprocessor SLES options is
152: available for each block. The following section of code is intended
153: to be a simple illustration of setting different linear solvers for
154: the individual blocks. These choices are obviously not recommended
155: for solving this particular problem.
156: */
157: PetscTypeCompare((PetscObject)pc,PCBJACOBI,&isbjacobi);
158: if (isbjacobi) {
159: /*
160: Call SLESSetUp() to set the block Jacobi data structures (including
161: creation of an internal SLES context for each block).
163: Note: SLESSetUp() MUST be called before PCBJacobiGetSubSLES().
164: */
165: SLESSetUp(sles,x,b);
167: /*
168: Extract the array of SLES contexts for the local blocks
169: */
170: PCBJacobiGetSubSLES(pc,&nlocal,&first,&subsles);
172: /*
173: Loop over the local blocks, setting various SLES options
174: for each block.
175: */
176: for (i=0; i<nlocal; i++) {
177: SLESGetPC(subsles[i],&subpc);
178: SLESGetKSP(subsles[i],&subksp);
179: if (!rank) {
180: if (i%2) {
181: PCSetType(subpc,PCILU);
182: } else {
183: PCSetType(subpc,PCNONE);
184: KSPSetType(subksp,KSPBCGS);
185: KSPSetTolerances(subksp,1.e-6,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
186: }
187: } else {
188: PCSetType(subpc,PCJACOBI);
189: KSPSetType(subksp,KSPGMRES);
190: KSPSetTolerances(subksp,1.e-7,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
191: }
192: }
193: }
195: /* -------------------------------------------------------------------
196: Solve the linear system
197: ------------------------------------------------------------------- */
199: /*
200: Set runtime options
201: */
202: SLESSetFromOptions(sles);
204: /*
205: Solve the linear system
206: */
207: SLESSolve(sles,b,x,&its);
209: /*
210: View info about the solver
211: */
212: PetscOptionsHasName(PETSC_NULL,"-noslesview",&flg);
213: if (!flg) {
214: SLESView(sles,PETSC_VIEWER_STDOUT_WORLD);
215: }
217: /* -------------------------------------------------------------------
218: Check solution and clean up
219: ------------------------------------------------------------------- */
221: /*
222: Check the error
223: */
224: VecAXPY(&none,u,x);
225: VecNorm(x,NORM_2,&norm);
226: PetscPrintf(PETSC_COMM_WORLD,"Norm of error %A iterations %dn",norm,its);
228: /*
229: Free work space. All PETSc objects should be destroyed when they
230: are no longer needed.
231: */
232: SLESDestroy(sles);
233: VecDestroy(u); VecDestroy(x);
234: VecDestroy(b); MatDestroy(A);
235: PetscFinalize();
236: return 0;
237: }