Actual source code: ex8.c
petsc-3.9.4 2018-09-11
2: static char help[] = "Illustrates use of the preconditioner ASM.\n\
3: The Additive Schwarz Method for solving a linear system in parallel with KSP. The\n\
4: code indicates the procedure for setting user-defined subdomains. Input\n\
5: parameters include:\n\
6: -user_set_subdomain_solvers: User explicitly sets subdomain solvers\n\
7: -user_set_subdomains: Activate user-defined subdomains\n\n";
9: /*
10: Note: This example focuses on setting the subdomains for the ASM
11: preconditioner for a problem on a 2D rectangular grid. See ex1.c
12: and ex2.c for more detailed comments on the basic usage of KSP
13: (including working with matrices and vectors).
15: The ASM preconditioner is fully parallel, but currently the routine
16: PCASMCreateSubdomains2D(), which is used in this example to demonstrate
17: user-defined subdomains (activated via -user_set_subdomains), is
18: uniprocessor only.
20: This matrix in this linear system arises from the discretized Laplacian,
21: and thus is not very interesting in terms of experimenting with variants
22: of the ASM preconditioner.
23: */
25: /*T
26: Concepts: KSP^Additive Schwarz Method (ASM) with user-defined subdomains
27: Processors: n
28: requires: x TODO: Need to determine if deprecated
29: T*/
33: /*
34: Include "petscksp.h" so that we can use KSP solvers. Note that this file
35: automatically includes:
36: petscsys.h - base PETSc routines petscvec.h - vectors
37: petscmat.h - matrices
38: petscis.h - index sets petscksp.h - Krylov subspace methods
39: petscviewer.h - viewers petscpc.h - preconditioners
40: */
41: #include <petscksp.h>
43: int main(int argc,char **args)
44: {
45: Vec x,b,u; /* approx solution, RHS, exact solution */
46: Mat A; /* linear system matrix */
47: KSP ksp; /* linear solver context */
48: PC pc; /* PC context */
49: IS *is,*is_local; /* array of index sets that define the subdomains */
50: PetscInt overlap = 1; /* width of subdomain overlap */
51: PetscInt Nsub; /* number of subdomains */
52: PetscInt m = 15,n = 17; /* mesh dimensions in x- and y- directions */
53: PetscInt M = 2,N = 1; /* number of subdomains in x- and y- directions */
54: PetscInt i,j,Ii,J,Istart,Iend;
56: PetscMPIInt size;
57: PetscBool flg;
58: PetscBool user_subdomains = PETSC_FALSE;
59: PetscScalar v, one = 1.0;
60: PetscReal e;
62: PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
63: MPI_Comm_size(PETSC_COMM_WORLD,&size);
64: PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
65: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
66: PetscOptionsGetInt(NULL,NULL,"-Mdomains",&M,NULL);
67: PetscOptionsGetInt(NULL,NULL,"-Ndomains",&N,NULL);
68: PetscOptionsGetInt(NULL,NULL,"-overlap",&overlap,NULL);
69: PetscOptionsGetBool(NULL,NULL,"-user_set_subdomains",&user_subdomains,NULL);
71: /* -------------------------------------------------------------------
72: Compute the matrix and right-hand-side vector that define
73: the linear system, Ax = b.
74: ------------------------------------------------------------------- */
76: /*
77: Assemble the matrix for the five point stencil, YET AGAIN
78: */
79: MatCreate(PETSC_COMM_WORLD,&A);
80: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n);
81: MatSetFromOptions(A);
82: MatSetUp(A);
83: MatGetOwnershipRange(A,&Istart,&Iend);
84: for (Ii=Istart; Ii<Iend; Ii++) {
85: v = -1.0; i = Ii/n; j = Ii - i*n;
86: if (i>0) {J = Ii - n; MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);}
87: if (i<m-1) {J = Ii + n; MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);}
88: if (j>0) {J = Ii - 1; MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);}
89: if (j<n-1) {J = Ii + 1; MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);}
90: v = 4.0; MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES);
91: }
92: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
93: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
95: /*
96: Create and set vectors
97: */
98: VecCreate(PETSC_COMM_WORLD,&b);
99: VecSetSizes(b,PETSC_DECIDE,m*n);
100: VecSetFromOptions(b);
101: VecDuplicate(b,&u);
102: VecDuplicate(b,&x);
103: VecSet(u,one);
104: MatMult(A,u,b);
106: /*
107: Create linear solver context
108: */
109: KSPCreate(PETSC_COMM_WORLD,&ksp);
111: /*
112: Set operators. Here the matrix that defines the linear system
113: also serves as the preconditioning matrix.
114: */
115: KSPSetOperators(ksp,A,A);
117: /*
118: Set the default preconditioner for this program to be ASM
119: */
120: KSPGetPC(ksp,&pc);
121: PCSetType(pc,PCASM);
123: /* -------------------------------------------------------------------
124: Define the problem decomposition
125: ------------------------------------------------------------------- */
127: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
128: Basic method, should be sufficient for the needs of many users.
129: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
131: Set the overlap, using the default PETSc decomposition via
132: PCASMSetOverlap(pc,overlap);
133: Could instead use the option -pc_asm_overlap <ovl>
135: Set the total number of blocks via -pc_asm_blocks <blks>
136: Note: The ASM default is to use 1 block per processor. To
137: experiment on a single processor with various overlaps, you
138: must specify use of multiple blocks!
139: */
141: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
142: More advanced method, setting user-defined subdomains
143: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
145: Firstly, create index sets that define the subdomains. The utility
146: routine PCASMCreateSubdomains2D() is a simple example (that currently
147: supports 1 processor only!). More generally, the user should write
148: a custom routine for a particular problem geometry.
150: Then call either PCASMSetLocalSubdomains() or PCASMSetTotalSubdomains()
151: to set the subdomains for the ASM preconditioner.
152: */
154: if (!user_subdomains) { /* basic version */
155: PCASMSetOverlap(pc,overlap);
156: } else { /* advanced version */
157: if (size != 1) SETERRQ(PETSC_COMM_WORLD,1,"PCASMCreateSubdomains2D() is currently a uniprocessor routine only!");
158: PCASMCreateSubdomains2D(m,n,M,N,1,overlap,&Nsub,&is,&is_local);
159: PCASMSetLocalSubdomains(pc,Nsub,is,is_local);
160: flg = PETSC_FALSE;
161: PetscOptionsGetBool(NULL,NULL,"-subdomain_view",&flg,NULL);
162: if (flg) {
163: PetscPrintf(PETSC_COMM_SELF,"Nmesh points: %D x %D; subdomain partition: %D x %D; overlap: %D; Nsub: %D\n",m,n,M,N,overlap,Nsub);
164: PetscPrintf(PETSC_COMM_SELF,"IS:\n");
165: for (i=0; i<Nsub; i++) {
166: PetscPrintf(PETSC_COMM_SELF," IS[%D]\n",i);
167: ISView(is[i],PETSC_VIEWER_STDOUT_SELF);
168: }
169: PetscPrintf(PETSC_COMM_SELF,"IS_local:\n");
170: for (i=0; i<Nsub; i++) {
171: PetscPrintf(PETSC_COMM_SELF," IS_local[%D]\n",i);
172: ISView(is_local[i],PETSC_VIEWER_STDOUT_SELF);
173: }
174: }
175: }
177: /* -------------------------------------------------------------------
178: Set the linear solvers for the subblocks
179: ------------------------------------------------------------------- */
181: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
182: Basic method, should be sufficient for the needs of most users.
183: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
185: By default, the ASM preconditioner uses the same solver on each
186: block of the problem. To set the same solver options on all blocks,
187: use the prefix -sub before the usual PC and KSP options, e.g.,
188: -sub_pc_type <pc> -sub_ksp_type <ksp> -sub_ksp_rtol 1.e-4
190: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
191: Advanced method, setting different solvers for various blocks.
192: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
194: Note that each block's KSP context is completely independent of
195: the others, and the full range of uniprocessor KSP options is
196: available for each block.
198: - Use PCASMGetSubKSP() to extract the array of KSP contexts for
199: the local blocks.
200: - See ex7.c for a simple example of setting different linear solvers
201: for the individual blocks for the block Jacobi method (which is
202: equivalent to the ASM method with zero overlap).
203: */
205: flg = PETSC_FALSE;
206: PetscOptionsGetBool(NULL,NULL,"-user_set_subdomain_solvers",&flg,NULL);
207: if (flg) {
208: KSP *subksp; /* array of KSP contexts for local subblocks */
209: PetscInt nlocal,first; /* number of local subblocks, first local subblock */
210: PC subpc; /* PC context for subblock */
211: PetscBool isasm;
213: PetscPrintf(PETSC_COMM_WORLD,"User explicitly sets subdomain solvers.\n");
215: /*
216: Set runtime options
217: */
218: KSPSetFromOptions(ksp);
220: /*
221: Flag an error if PCTYPE is changed from the runtime options
222: */
223: PetscObjectTypeCompare((PetscObject)pc,PCASM,&isasm);
224: if (!isasm) SETERRQ(PETSC_COMM_WORLD,1,"Cannot Change the PCTYPE when manually changing the subdomain solver settings");
226: /*
227: Call KSPSetUp() to set the block Jacobi data structures (including
228: creation of an internal KSP context for each block).
230: Note: KSPSetUp() MUST be called before PCASMGetSubKSP().
231: */
232: KSPSetUp(ksp);
234: /*
235: Extract the array of KSP contexts for the local blocks
236: */
237: PCASMGetSubKSP(pc,&nlocal,&first,&subksp);
239: /*
240: Loop over the local blocks, setting various KSP options
241: for each block.
242: */
243: for (i=0; i<nlocal; i++) {
244: KSPGetPC(subksp[i],&subpc);
245: PCSetType(subpc,PCILU);
246: KSPSetType(subksp[i],KSPGMRES);
247: KSPSetTolerances(subksp[i],1.e-7,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
248: }
249: } else {
250: /*
251: Set runtime options
252: */
253: KSPSetFromOptions(ksp);
254: }
256: /* -------------------------------------------------------------------
257: Solve the linear system
258: ------------------------------------------------------------------- */
260: KSPSolve(ksp,b,x);
262: /* -------------------------------------------------------------------
263: Compare result to the exact solution
264: ------------------------------------------------------------------- */
265: VecAXPY(x,-1.0,u);
266: VecNorm(x,NORM_INFINITY, &e);
268: flg = PETSC_FALSE;
269: PetscOptionsGetBool(NULL,NULL,"-print_error",&flg,NULL);
270: if (flg) {
271: PetscPrintf(PETSC_COMM_WORLD, "Infinity norm of the error: %g\n",(double) e);
272: }
274: /*
275: Free work space. All PETSc objects should be destroyed when they
276: are no longer needed.
277: */
279: if (user_subdomains) {
280: for (i=0; i<Nsub; i++) {
281: ISDestroy(&is[i]);
282: ISDestroy(&is_local[i]);
283: }
284: PetscFree(is);
285: PetscFree(is_local);
286: }
287: KSPDestroy(&ksp);
288: VecDestroy(&u);
289: VecDestroy(&x);
290: VecDestroy(&b);
291: MatDestroy(&A);
292: PetscFinalize();
293: return ierr;
294: }
297: /*TEST
299: test:
300: suffix: 1
301: args: -print_error
303: TEST*/