Actual source code: asm.c
1: /*$Id: asm.c,v 1.129 2001/04/10 22:35:59 bsmith Exp $*/
2: /*
3: This file defines an additive Schwarz preconditioner for any Mat implementation.
5: Note that each processor may have any number of subdomains. But in order to
6: deal easily with the VecScatter(), we treat each processor as if it has the
7: same number of subdomains.
9: n - total number of true subdomains on all processors
10: n_local_true - actual number of subdomains on this processor
11: n_local = maximum over all processors of n_local_true
12: */
13: #include src/sles/pc/pcimpl.h
14: #include petscsles.h
16: typedef struct {
17: int n,n_local,n_local_true;
18: PetscTruth is_flg; /* flg set to 1 if the IS created in pcsetup */
19: int overlap; /* overlap requested by user */
20: SLES *sles; /* linear solvers for each block */
21: VecScatter *scat; /* mapping to subregion */
22: Vec *x,*y;
23: IS *is; /* index set that defines each subdomain */
24: Mat *mat,*pmat; /* mat is not currently used */
25: PCASMType type; /* use reduced interpolation, restriction or both */
26: PetscTruth same_local_solves; /* flag indicating whether all local solvers are same */
27: PetscTruth inplace; /* indicates that the sub-matrices are deleted after
28: PCSetUpOnBlocks() is done. Similar to inplace
29: factorization in the case of LU and ILU */
30: } PC_ASM;
32: static int PCView_ASM(PC pc,PetscViewer viewer)
33: {
34: PC_ASM *jac = (PC_ASM*)pc->data;
35: int rank,ierr,i;
36: char *cstring = 0;
37: PetscTruth isascii,isstring;
38: PetscViewer sviewer;
42: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
43: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);
44: if (isascii) {
45: PetscViewerASCIIPrintf(viewer," Additive Schwarz: total subdomain blocks = %d, amount of overlap = %dn",jac->n,jac->overlap);
46: if (jac->type == PC_ASM_NONE) cstring = "limited restriction and interpolation (PC_ASM_NONE)";
47: else if (jac->type == PC_ASM_RESTRICT) cstring = "full restriction (PC_ASM_RESTRICT)";
48: else if (jac->type == PC_ASM_INTERPOLATE) cstring = "full interpolation (PC_ASM_INTERPOLATE)";
49: else if (jac->type == PC_ASM_BASIC) cstring = "full restriction and interpolation (PC_ASM_BASIC)";
50: else cstring = "Unknown ASM type";
51: PetscViewerASCIIPrintf(viewer," Additive Schwarz: type - %sn",cstring);
52: MPI_Comm_rank(pc->comm,&rank);
53: if (jac->same_local_solves) {
54: PetscViewerASCIIPrintf(viewer," Local solve is same for all blocks, in the following KSP and PC objects:n");
55: PetscViewerGetSingleton(viewer,&sviewer);
56: if (!rank && jac->sles) {
57: PetscViewerASCIIPushTab(viewer);
58: SLESView(jac->sles[0],sviewer);
59: PetscViewerASCIIPopTab(viewer);
60: }
61: PetscViewerRestoreSingleton(viewer,&sviewer);
62: } else {
63: PetscViewerASCIIPrintf(viewer," Local solve info for each block is in the following KSP and PC objects:n");
64: PetscViewerASCIISynchronizedPrintf(viewer,"Proc %d: number of local blocks = %dn",rank,jac->n_local);
65: PetscViewerASCIIPushTab(viewer);
66: for (i=0; i<jac->n_local; i++) {
67: PetscViewerASCIISynchronizedPrintf(viewer,"Proc %d: local block number %dn",rank,i);
68: PetscViewerGetSingleton(viewer,&sviewer);
69: SLESView(jac->sles[i],sviewer);
70: PetscViewerRestoreSingleton(viewer,&sviewer);
71: if (i != jac->n_local-1) {
72: PetscViewerASCIISynchronizedPrintf(viewer,"- - - - - - - - - - - - - - - - - -n");
73: }
74: }
75: PetscViewerASCIIPopTab(viewer);
76: PetscViewerFlush(viewer);
77: }
78: } else if (isstring) {
79: PetscViewerStringSPrintf(viewer," blks=%d, overlap=%d, type=%d",jac->n,jac->overlap,jac->type);
80: PetscViewerGetSingleton(viewer,&sviewer);
81: if (jac->sles) {SLESView(jac->sles[0],sviewer);}
82: PetscViewerGetSingleton(viewer,&sviewer);
83: } else {
84: SETERRQ1(1,"Viewer type %s not supported for PCASM",((PetscObject)viewer)->type_name);
85: }
86: return(0);
87: }
89: static int PCSetUp_ASM(PC pc)
90: {
91: PC_ASM *osm = (PC_ASM*)pc->data;
92: int i,ierr,m,n_local = osm->n_local,n_local_true = osm->n_local_true;
93: int start,start_val,end_val,size,sz,bs;
94: MatReuse scall = MAT_REUSE_MATRIX;
95: IS isl;
96: SLES sles;
97: KSP subksp;
98: PC subpc;
99: char *prefix,*pprefix;
102: if (!pc->setupcalled) {
103: if (osm->n == PETSC_DECIDE && osm->n_local_true == PETSC_DECIDE) {
104: /* no subdomains given, use one per processor */
105: osm->n_local_true = osm->n_local = 1;
106: MPI_Comm_size(pc->comm,&size);
107: osm->n = size;
108: } else if (osm->n == PETSC_DECIDE) { /* determine global number of subdomains */
109: int inwork[2],outwork[2];
110: inwork[0] = inwork[1] = osm->n_local_true;
111: MPI_Allreduce(inwork,outwork,2,MPI_INT,PetscMaxSum_Op,pc->comm);
112: osm->n_local = outwork[0];
113: osm->n = outwork[1];
114: }
115: n_local = osm->n_local;
116: n_local_true = osm->n_local_true;
117: if (!osm->is){ /* build the index sets */
118: ierr = PetscMalloc((n_local_true+1)*sizeof(IS **),&osm->is);
119: ierr = MatGetOwnershipRange(pc->pmat,&start_val,&end_val);
120: ierr = MatGetBlockSize(pc->pmat,&bs);
121: sz = end_val - start_val;
122: start = start_val;
123: if (end_val/bs*bs != end_val || start_val/bs*bs != start_val) {
124: SETERRQ(PETSC_ERR_ARG_WRONG,"Bad distribution for matrix block size");
125: }
126: for (i=0; i<n_local_true; i++){
127: size = ((sz/bs)/n_local_true + (((sz/bs) % n_local_true) > i))*bs;
128: ierr = ISCreateStride(PETSC_COMM_SELF,size,start,1,&isl);
129: start += size;
130: osm->is[i] = isl;
131: }
132: osm->is_flg = PETSC_TRUE;
133: }
135: ierr = PetscMalloc((n_local_true+1)*sizeof(SLES **),&osm->sles);
136: ierr = PetscMalloc(n_local*sizeof(VecScatter **),&osm->scat);
137: ierr = PetscMalloc(2*n_local*sizeof(Vec **),&osm->x);
138: osm->y = osm->x + n_local;
140: /* Extend the "overlapping" regions by a number of steps */
141: MatIncreaseOverlap(pc->pmat,n_local_true,osm->is,osm->overlap);
142: for (i=0; i<n_local_true; i++) {
143: ISSort(osm->is[i]);
144: }
146: /* create the local work vectors and scatter contexts */
147: for (i=0; i<n_local_true; i++) {
148: ISGetLocalSize(osm->is[i],&m);
149: VecCreateSeq(PETSC_COMM_SELF,m,&osm->x[i]);
150: VecDuplicate(osm->x[i],&osm->y[i]);
151: ISCreateStride(PETSC_COMM_SELF,m,0,1,&isl);
152: VecScatterCreate(pc->vec,osm->is[i],osm->x[i],isl,&osm->scat[i]);
153: ISDestroy(isl);
154: }
155: for (i=n_local_true; i<n_local; i++) {
156: VecCreateSeq(PETSC_COMM_SELF,0,&osm->x[i]);
157: VecDuplicate(osm->x[i],&osm->y[i]);
158: ISCreateStride(PETSC_COMM_SELF,0,0,1,&isl);
159: VecScatterCreate(pc->vec,isl,osm->x[i],isl,&osm->scat[i]);
160: ISDestroy(isl);
161: }
163: /*
164: Create the local solvers.
165: */
166: for (i=0; i<n_local_true; i++) {
167: SLESCreate(PETSC_COMM_SELF,&sles);
168: PetscLogObjectParent(pc,sles);
169: SLESGetKSP(sles,&subksp);
170: KSPSetType(subksp,KSPPREONLY);
171: SLESGetPC(sles,&subpc);
172: PCSetType(subpc,PCILU);
173: PCGetOptionsPrefix(pc,&prefix);
174: SLESSetOptionsPrefix(sles,prefix);
175: SLESAppendOptionsPrefix(sles,"sub_");
176: SLESSetFromOptions(sles);
177: osm->sles[i] = sles;
178: }
179: scall = MAT_INITIAL_MATRIX;
180: } else {
181: /*
182: Destroy the blocks from the previous iteration
183: */
184: if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
185: MatDestroyMatrices(osm->n_local_true,&osm->pmat);
186: scall = MAT_INITIAL_MATRIX;
187: }
188: }
190: /* extract out the submatrices */
191: MatGetSubMatrices(pc->pmat,osm->n_local_true,osm->is,osm->is,scall,&osm->pmat);
193: /* Return control to the user so that the submatrices can be modified (e.g., to apply
194: different boundary conditions for the submatrices than for the global problem) */
195: PCModifySubMatrices(pc,osm->n_local,osm->is,osm->is,osm->pmat,pc->modifysubmatricesP);
197: /* loop over subdomains putting them into local sles */
198: PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix);
199: for (i=0; i<n_local_true; i++) {
200: PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i],pprefix);
201: PetscLogObjectParent(pc,osm->pmat[i]);
202: SLESSetOperators(osm->sles[i],osm->pmat[i],osm->pmat[i],pc->flag);
203: }
204: return(0);
205: }
207: static int PCSetUpOnBlocks_ASM(PC pc)
208: {
209: PC_ASM *osm = (PC_ASM*)pc->data;
210: int i,ierr;
213: for (i=0; i<osm->n_local_true; i++) {
214: SLESSetUp(osm->sles[i],osm->x[i],osm->y[i]);
215: }
216: /*
217: If inplace flag is set, then destroy the matrix after the setup
218: on blocks is done.
219: */
220: if (osm->inplace && osm->n_local_true > 0) {
221: MatDestroyMatrices(osm->n_local_true,&osm->pmat);
222: }
223: return(0);
224: }
226: static int PCApply_ASM(PC pc,Vec x,Vec y)
227: {
228: PC_ASM *osm = (PC_ASM*)pc->data;
229: int i,n_local = osm->n_local,n_local_true = osm->n_local_true,ierr,its;
230: Scalar zero = 0.0;
231: ScatterMode forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;
234: /*
235: Support for limiting the restriction or interpolation to only local
236: subdomain values (leaving the other values 0).
237: */
238: if (!(osm->type & PC_ASM_RESTRICT)) {
239: forward = SCATTER_FORWARD_LOCAL;
240: /* have to zero the work RHS since scatter may leave some slots empty */
241: for (i=0; i<n_local; i++) {
242: VecSet(&zero,osm->x[i]);
243: }
244: }
245: if (!(osm->type & PC_ASM_INTERPOLATE)) {
246: reverse = SCATTER_REVERSE_LOCAL;
247: }
249: for (i=0; i<n_local; i++) {
250: VecScatterBegin(x,osm->x[i],INSERT_VALUES,forward,osm->scat[i]);
251: }
252: VecSet(&zero,y);
253: /* do the local solves */
254: for (i=0; i<n_local_true; i++) {
255: VecScatterEnd(x,osm->x[i],INSERT_VALUES,forward,osm->scat[i]);
256: SLESSolve(osm->sles[i],osm->x[i],osm->y[i],&its);
257: VecScatterBegin(osm->y[i],y,ADD_VALUES,reverse,osm->scat[i]);
258: }
259: /* handle the rest of the scatters that do not have local solves */
260: for (i=n_local_true; i<n_local; i++) {
261: VecScatterEnd(x,osm->x[i],INSERT_VALUES,forward,osm->scat[i]);
262: VecScatterBegin(osm->y[i],y,ADD_VALUES,reverse,osm->scat[i]);
263: }
264: for (i=0; i<n_local; i++) {
265: VecScatterEnd(osm->y[i],y,ADD_VALUES,reverse,osm->scat[i]);
266: }
267: return(0);
268: }
270: static int PCApplyTranspose_ASM(PC pc,Vec x,Vec y)
271: {
272: PC_ASM *osm = (PC_ASM*)pc->data;
273: int i,n_local = osm->n_local,n_local_true = osm->n_local_true,ierr,its;
274: Scalar zero = 0.0;
275: ScatterMode forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;
278: /*
279: Support for limiting the restriction or interpolation to only local
280: subdomain values (leaving the other values 0).
282: Note: these are reversed from the PCApply_ASM() because we are applying the
283: transpose of the three terms
284: */
285: if (!(osm->type & PC_ASM_INTERPOLATE)) {
286: forward = SCATTER_FORWARD_LOCAL;
287: /* have to zero the work RHS since scatter may leave some slots empty */
288: for (i=0; i<n_local; i++) {
289: VecSet(&zero,osm->x[i]);
290: }
291: }
292: if (!(osm->type & PC_ASM_RESTRICT)) {
293: reverse = SCATTER_REVERSE_LOCAL;
294: }
296: for (i=0; i<n_local; i++) {
297: VecScatterBegin(x,osm->x[i],INSERT_VALUES,forward,osm->scat[i]);
298: }
299: VecSet(&zero,y);
300: /* do the local solves */
301: for (i=0; i<n_local_true; i++) {
302: VecScatterEnd(x,osm->x[i],INSERT_VALUES,forward,osm->scat[i]);
303: SLESSolveTranspose(osm->sles[i],osm->x[i],osm->y[i],&its);
304: VecScatterBegin(osm->y[i],y,ADD_VALUES,reverse,osm->scat[i]);
305: }
306: /* handle the rest of the scatters that do not have local solves */
307: for (i=n_local_true; i<n_local; i++) {
308: VecScatterEnd(x,osm->x[i],INSERT_VALUES,forward,osm->scat[i]);
309: VecScatterBegin(osm->y[i],y,ADD_VALUES,reverse,osm->scat[i]);
310: }
311: for (i=0; i<n_local; i++) {
312: VecScatterEnd(osm->y[i],y,ADD_VALUES,reverse,osm->scat[i]);
313: }
314: return(0);
315: }
317: static int PCDestroy_ASM(PC pc)
318: {
319: PC_ASM *osm = (PC_ASM*)pc->data;
320: int i,ierr;
323: for (i=0; i<osm->n_local; i++) {
324: VecScatterDestroy(osm->scat[i]);
325: VecDestroy(osm->x[i]);
326: VecDestroy(osm->y[i]);
327: }
328: if (osm->n_local_true > 0 && !osm->inplace && osm->pmat) {
329: MatDestroyMatrices(osm->n_local_true,&osm->pmat);
330: }
331: for (i=0; i<osm->n_local_true; i++) {
332: SLESDestroy(osm->sles[i]);
333: }
334: if (osm->is_flg) {
335: for (i=0; i<osm->n_local_true; i++) {ISDestroy(osm->is[i]);}
336: PetscFree(osm->is);
337: }
338: if (osm->sles) {PetscFree(osm->sles);}
339: if (osm->scat) {PetscFree(osm->scat);}
340: if (osm->x) {PetscFree(osm->x);}
341: PetscFree(osm);
342: return(0);
343: }
345: static int PCSetFromOptions_ASM(PC pc)
346: {
347: PC_ASM *osm = (PC_ASM*)pc->data;
348: int blocks,ovl,ierr;
349: PetscTruth flg;
350: char buff[16],*type[] = {"basic","restrict","interpolate","none"};
353: PetscOptionsHead("Additive Schwarz options");
354: PetscOptionsInt("-pc_asm_blocks","Number of subdomains","PCASMSetTotalSubdomains",osm->n,&blocks,&flg);
355: if (flg) {PCASMSetTotalSubdomains(pc,blocks,PETSC_NULL); }
356: PetscOptionsInt("-pc_asm_overlap","Number of grid points overlap","PCASMSetOverlap",osm->overlap,&ovl,&flg);
357: if (flg) {PCASMSetOverlap(pc,ovl); }
358: PetscOptionsName("-pc_asm_in_place","Perform matrix factorization inplace","PCASMSetUseInPlace",&flg);
359: if (flg) {PCASMSetUseInPlace(pc); }
360: PetscOptionsEList("-pc_asm_type","Type of restriction/extension","PCASMSetType",type,4,"restrict",buff,15,&flg);
361: if (flg) {
362: PCASMType atype;
363: PetscTruth isbasic,isrestrict,isinterpolate,isnone;
365: PetscStrcmp(buff,type[0],&isbasic);
366: PetscStrcmp(buff,type[1],&isrestrict);
367: PetscStrcmp(buff,type[2],&isinterpolate);
368: PetscStrcmp(buff,type[3],&isnone);
370: if (isbasic) atype = PC_ASM_BASIC;
371: else if (isrestrict) atype = PC_ASM_RESTRICT;
372: else if (isinterpolate) atype = PC_ASM_INTERPOLATE;
373: else if (isnone) atype = PC_ASM_NONE;
374: else SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Unknown type");
375: PCASMSetType(pc,atype);
376: }
377: PetscOptionsTail();
378: return(0);
379: }
381: /*------------------------------------------------------------------------------------*/
383: EXTERN_C_BEGIN
384: int PCASMSetLocalSubdomains_ASM(PC pc,int n,IS *is)
385: {
386: PC_ASM *osm;
390: if (pc->setupcalled) {
391: SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"PCASMSetLocalSubdomains() should be called before calling PCSetup().");
392: }
394: if (n < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Each process must have 0 or more blocks");
395: osm = (PC_ASM*)pc->data;
396: osm->n_local_true = n;
397: osm->is = is;
398: return(0);
399: }
400: EXTERN_C_END
402: EXTERN_C_BEGIN
403: int PCASMSetTotalSubdomains_ASM(PC pc,int N,IS *is)
404: {
405: PC_ASM *osm;
406: int rank,size,ierr;
409: if (pc->setupcalled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"PCASMSetTotalSubdomains() should be called before PCSetup().");
411: if (is) SETERRQ(PETSC_ERR_SUP,"Use PCASMSetLocalSubdomains to set specific index setsn
412: they cannot be set globally yet.");
414: osm = (PC_ASM*)pc->data;
415: /*
416: Split the subdomains equally amoung all processors
417: */
418: MPI_Comm_rank(pc->comm,&rank);
419: MPI_Comm_size(pc->comm,&size);
420: osm->n_local_true = N/size + ((N % size) > rank);
421: if (osm->n_local_true < 0) {
422: SETERRQ(PETSC_ERR_SUP,"Each process must have 0 or more blocks");
423: }
424: osm->is = 0;
425: return(0);
426: }
427: EXTERN_C_END
429: EXTERN_C_BEGIN
430: int PCASMSetOverlap_ASM(PC pc,int ovl)
431: {
432: PC_ASM *osm;
435: if (ovl < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested");
437: osm = (PC_ASM*)pc->data;
438: osm->overlap = ovl;
439: return(0);
440: }
441: EXTERN_C_END
443: EXTERN_C_BEGIN
444: int PCASMSetType_ASM(PC pc,PCASMType type)
445: {
446: PC_ASM *osm;
449: osm = (PC_ASM*)pc->data;
450: osm->type = type;
451: return(0);
452: }
453: EXTERN_C_END
455: EXTERN_C_BEGIN
456: int PCASMGetSubSLES_ASM(PC pc,int *n_local,int *first_local,SLES **sles)
457: {
458: PC_ASM *jac = (PC_ASM*)pc->data;
459: int ierr;
462: if (jac->n_local_true < 0) {
463: SETERRQ(1,"Need to call PCSetUP() on PC (or SLESSetUp() on the outer SLES object) before calling here");
464: }
466: if (n_local) *n_local = jac->n_local_true;
467: if (first_local) {
468: MPI_Scan(&jac->n_local_true,first_local,1,MPI_INT,MPI_SUM,pc->comm);
469: *first_local -= jac->n_local_true;
470: }
471: *sles = jac->sles;
472: jac->same_local_solves = PETSC_FALSE; /* Assume that local solves are now different;
473: not necessarily true though! This flag is
474: used only for PCView_ASM() */
475: return(0);
476: }
477: EXTERN_C_END
479: EXTERN_C_BEGIN
480: int PCASMSetUseInPlace_ASM(PC pc)
481: {
482: PC_ASM *dir;
485: dir = (PC_ASM*)pc->data;
486: dir->inplace = PETSC_TRUE;
487: return(0);
488: }
489: EXTERN_C_END
491: /*----------------------------------------------------------------------------*/
492: /*@
493: PCASMSetUseInPlace - Tells the system to destroy the matrix after setup is done.
495: Collective on PC
497: Input Parameters:
498: . pc - the preconditioner context
500: Options Database Key:
501: . -pc_asm_in_place - Activates in-place factorization
503: Note:
504: PCASMSetUseInplace() can only be used with the KSP method KSPPREONLY, and
505: when the original matrix is not required during the Solve process.
506: This destroys the matrix, early thus, saving on memory usage.
508: Level: intermediate
510: .keywords: PC, set, factorization, direct, inplace, in-place, ASM
512: .seealso: PCILUSetUseInPlace(), PCLUSetUseInPlace ()
513: @*/
514: int PCASMSetUseInPlace(PC pc)
515: {
516: int ierr,(*f)(PC);
520: PetscObjectQueryFunction((PetscObject)pc,"PCASMSetUseInPlace_C",(void (**)())&f);
521: if (f) {
522: (*f)(pc);
523: }
524: return(0);
525: }
526: /*----------------------------------------------------------------------------*/
528: /*@C
529: PCASMSetLocalSubdomains - Sets the local subdomains (for this processor
530: only) for the additive Schwarz preconditioner.
532: Collective on PC
534: Input Parameters:
535: + pc - the preconditioner context
536: . n - the number of subdomains for this processor (default value = 1)
537: - is - the index sets that define the subdomains for this processor
538: (or PETSC_NULL for PETSc to determine subdomains)
540: Notes:
541: The IS numbering is in the parallel, global numbering of the vector.
543: By default the ASM preconditioner uses 1 block per processor.
545: These index sets cannot be destroyed until after completion of the
546: linear solves for which the ASM preconditioner is being used.
548: Use PCASMSetTotalSubdomains() to set the subdomains for all processors.
550: Level: advanced
552: .keywords: PC, ASM, set, local, subdomains, additive Schwarz
554: .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubSLES(),
555: PCASMCreateSubdomains2D(), PCASMGetLocalSubdomains()
556: @*/
557: int PCASMSetLocalSubdomains(PC pc,int n,IS *is)
558: {
559: int ierr,(*f)(PC,int,IS *);
563: PetscObjectQueryFunction((PetscObject)pc,"PCASMSetLocalSubdomains_C",(void (**)())&f);
564: if (f) {
565: (*f)(pc,n,is);
566: }
567: return(0);
568: }
570: /*@C
571: PCASMSetTotalSubdomains - Sets the subdomains for all processor for the
572: additive Schwarz preconditioner. Either all or no processors in the
573: PC communicator must call this routine, with the same index sets.
575: Collective on PC
577: Input Parameters:
578: + pc - the preconditioner context
579: . n - the number of subdomains for all processors
580: - is - the index sets that define the subdomains for all processor
581: (or PETSC_NULL for PETSc to determine subdomains)
583: Options Database Key:
584: To set the total number of subdomain blocks rather than specify the
585: index sets, use the option
586: . -pc_asm_blocks <blks> - Sets total blocks
588: Notes:
589: Currently you cannot use this to set the actual subdomains with the argument is.
591: By default the ASM preconditioner uses 1 block per processor.
593: These index sets cannot be destroyed until after completion of the
594: linear solves for which the ASM preconditioner is being used.
596: Use PCASMSetLocalSubdomains() to set local subdomains.
598: Level: advanced
600: .keywords: PC, ASM, set, total, global, subdomains, additive Schwarz
602: .seealso: PCASMSetLocalSubdomains(), PCASMSetOverlap(), PCASMGetSubSLES(),
603: PCASMCreateSubdomains2D()
604: @*/
605: int PCASMSetTotalSubdomains(PC pc,int N,IS *is)
606: {
607: int ierr,(*f)(PC,int,IS *);
611: PetscObjectQueryFunction((PetscObject)pc,"PCASMSetTotalSubdomains_C",(void (**)())&f);
612: if (f) {
613: (*f)(pc,N,is);
614: }
615: return(0);
616: }
618: /*@
619: PCASMSetOverlap - Sets the overlap between a pair of subdomains for the
620: additive Schwarz preconditioner. Either all or no processors in the
621: PC communicator must call this routine.
623: Collective on PC
625: Input Parameters:
626: + pc - the preconditioner context
627: - ovl - the amount of overlap between subdomains (ovl >= 0, default value = 1)
629: Options Database Key:
630: . -pc_asm_overlap <ovl> - Sets overlap
632: Notes:
633: By default the ASM preconditioner uses 1 block per processor. To use
634: multiple blocks per perocessor, see PCASMSetTotalSubdomains() and
635: PCASMSetLocalSubdomains() (and the option -pc_asm_blocks <blks>).
637: The overlap defaults to 1, so if one desires that no additional
638: overlap be computed beyond what may have been set with a call to
639: PCASMSetTotalSubdomains() or PCASMSetLocalSubdomains(), then ovl
640: must be set to be 0. In particular, if one does not explicitly set
641: the subdomains an application code, then all overlap would be computed
642: internally by PETSc, and using an overlap of 0 would result in an ASM
643: variant that is equivalent to the block Jacobi preconditioner.
645: Note that one can define initial index sets with any overlap via
646: PCASMSetTotalSubdomains() or PCASMSetLocalSubdomains(); the routine
647: PCASMSetOverlap() merely allows PETSc to extend that overlap further
648: if desired.
650: Level: intermediate
652: .keywords: PC, ASM, set, overlap
654: .seealso: PCASMSetTotalSubdomains(), PCASMSetLocalSubdomains(), PCASMGetSubSLES(),
655: PCASMCreateSubdomains2D(), PCASMGetLocalSubdomains()
656: @*/
657: int PCASMSetOverlap(PC pc,int ovl)
658: {
659: int ierr,(*f)(PC,int);
663: PetscObjectQueryFunction((PetscObject)pc,"PCASMSetOverlap_C",(void (**)())&f);
664: if (f) {
665: (*f)(pc,ovl);
666: }
667: return(0);
668: }
670: /*@
671: PCASMSetType - Sets the type of restriction and interpolation used
672: for local problems in the additive Schwarz method.
674: Collective on PC
676: Input Parameters:
677: + pc - the preconditioner context
678: - type - variant of ASM, one of
679: .vb
680: PC_ASM_BASIC - full interpolation and restriction
681: PC_ASM_RESTRICT - full restriction, local processor interpolation
682: PC_ASM_INTERPOLATE - full interpolation, local processor restriction
683: PC_ASM_NONE - local processor restriction and interpolation
684: .ve
686: Options Database Key:
687: $ -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type
689: Level: intermediate
691: .keywords: PC, ASM, set, type
693: .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubSLES(),
694: PCASMCreateSubdomains2D()
695: @*/
696: int PCASMSetType(PC pc,PCASMType type)
697: {
698: int ierr,(*f)(PC,PCASMType);
702: PetscObjectQueryFunction((PetscObject)pc,"PCASMSetType_C",(void (**)())&f);
703: if (f) {
704: (*f)(pc,type);
705: }
706: return(0);
707: }
709: /*@C
710: PCASMGetSubSLES - Gets the local SLES contexts for all blocks on
711: this processor.
712:
713: Collective on PC iff first_local is requested
715: Input Parameter:
716: . pc - the preconditioner context
718: Output Parameters:
719: + n_local - the number of blocks on this processor or PETSC_NULL
720: . first_local - the global number of the first block on this processor or PETSC_NULL,
721: all processors must request or all must pass PETSC_NULL
722: - sles - the array of SLES contexts
724: Note:
725: After PCASMGetSubSLES() the array of SLESes is not to be freed
727: Currently for some matrix implementations only 1 block per processor
728: is supported.
729:
730: You must call SLESSetUp() before calling PCASMGetSubSLES().
732: Level: advanced
734: .keywords: PC, ASM, additive Schwarz, get, sub, SLES, context
736: .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap(),
737: PCASMCreateSubdomains2D(),
738: @*/
739: int PCASMGetSubSLES(PC pc,int *n_local,int *first_local,SLES **sles)
740: {
741: int ierr,(*f)(PC,int*,int*,SLES **);
745: PetscObjectQueryFunction((PetscObject)pc,"PCASMGetSubSLES_C",(void (**)())&f);
746: if (f) {
747: (*f)(pc,n_local,first_local,sles);
748: } else {
749: SETERRQ(1,"Cannot get subsles for this type of PC");
750: }
752: return(0);
753: }
755: /* -------------------------------------------------------------------------------------*/
756: EXTERN_C_BEGIN
757: int PCCreate_ASM(PC pc)
758: {
759: int ierr;
760: PC_ASM *osm;
763: PetscNew(PC_ASM,&osm);
764: PetscLogObjectMemory(pc,sizeof(PC_ASM));
765: PetscMemzero(osm,sizeof(PC_ASM));
766: osm->n = PETSC_DECIDE;
767: osm->n_local = 0;
768: osm->n_local_true = PETSC_DECIDE;
769: osm->overlap = 1;
770: osm->is_flg = PETSC_FALSE;
771: osm->sles = 0;
772: osm->scat = 0;
773: osm->is = 0;
774: osm->mat = 0;
775: osm->pmat = 0;
776: osm->type = PC_ASM_RESTRICT;
777: osm->same_local_solves = PETSC_TRUE;
778: osm->inplace = PETSC_FALSE;
779: pc->data = (void*)osm;
781: pc->ops->apply = PCApply_ASM;
782: pc->ops->applytranspose = PCApplyTranspose_ASM;
783: pc->ops->setup = PCSetUp_ASM;
784: pc->ops->destroy = PCDestroy_ASM;
785: pc->ops->setfromoptions = PCSetFromOptions_ASM;
786: pc->ops->setuponblocks = PCSetUpOnBlocks_ASM;
787: pc->ops->view = PCView_ASM;
788: pc->ops->applyrichardson = 0;
790: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCASMSetLocalSubdomains_C","PCASMSetLocalSubdomains_ASM",
791: PCASMSetLocalSubdomains_ASM);
792: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCASMSetTotalSubdomains_C","PCASMSetTotalSubdomains_ASM",
793: PCASMSetTotalSubdomains_ASM);
794: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCASMSetOverlap_C","PCASMSetOverlap_ASM",
795: PCASMSetOverlap_ASM);
796: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCASMSetType_C","PCASMSetType_ASM",
797: PCASMSetType_ASM);
798: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCASMGetSubSLES_C","PCASMGetSubSLES_ASM",
799: PCASMGetSubSLES_ASM);
800: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCASMSetUseInPlace_C","PCASMSetUseInPlace_ASM",
801: PCASMSetUseInPlace_ASM);
802: return(0);
803: }
804: EXTERN_C_END
807: /*@
808: PCASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz
809: preconditioner for a two-dimensional problem on a regular grid.
811: Not Collective
813: Input Parameters:
814: + m, n - the number of mesh points in the x and y directions
815: . M, N - the number of subdomains in the x and y directions
816: . dof - degrees of freedom per node
817: - overlap - overlap in mesh lines
819: Output Parameters:
820: + Nsub - the number of subdomains created
821: - is - the array of index sets defining the subdomains
823: Note:
824: Presently PCAMSCreateSubdomains2d() is valid only for sequential
825: preconditioners. More general related routines are
826: PCASMSetTotalSubdomains() and PCASMSetLocalSubdomains().
828: Level: advanced
830: .keywords: PC, ASM, additive Schwarz, create, subdomains, 2D, regular grid
832: .seealso: PCASMSetTotalSubdomains(), PCASMSetLocalSubdomains(), PCASMGetSubSLES(),
833: PCASMSetOverlap()
834: @*/
835: int PCASMCreateSubdomains2D(int m,int n,int M,int N,int dof,int overlap,int *Nsub,IS **is)
836: {
837: int i,j,height,width,ystart,xstart,yleft,yright,xleft,xright,loc_outter;
838: int nidx,*idx,loc,ii,jj,ierr,count;
841: if (dof != 1) SETERRQ(PETSC_ERR_SUP,"");
843: *Nsub = N*M;
844: PetscMalloc((*Nsub)*sizeof(IS **),is);
845: ystart = 0;
846: loc_outter = 0;
847: for (i=0; i<N; i++) {
848: height = n/N + ((n % N) > i); /* height of subdomain */
849: if (height < 2) SETERRQ(1,"Too many N subdomains for mesh dimension n");
850: yleft = ystart - overlap; if (yleft < 0) yleft = 0;
851: yright = ystart + height + overlap; if (yright > n) yright = n;
852: xstart = 0;
853: for (j=0; j<M; j++) {
854: width = m/M + ((m % M) > j); /* width of subdomain */
855: if (width < 2) SETERRQ(1,"Too many M subdomains for mesh dimension m");
856: xleft = xstart - overlap; if (xleft < 0) xleft = 0;
857: xright = xstart + width + overlap; if (xright > m) xright = m;
858: nidx = (xright - xleft)*(yright - yleft);
859: PetscMalloc(nidx*sizeof(int),&idx);
860: loc = 0;
861: for (ii=yleft; ii<yright; ii++) {
862: count = m*ii + xleft;
863: for (jj=xleft; jj<xright; jj++) {
864: idx[loc++] = count++;
865: }
866: }
867: ISCreateGeneral(PETSC_COMM_SELF,nidx,idx,(*is)+loc_outter++);
868: PetscFree(idx);
869: xstart += width;
870: }
871: ystart += height;
872: }
873: for (i=0; i<*Nsub; i++) { ISSort((*is)[i]); }
874: return(0);
875: }
877: /*@C
878: PCASMGetLocalSubdomains - Gets the local subdomains (for this processor
879: only) for the additive Schwarz preconditioner.
881: Collective on PC
883: Input Parameter:
884: . pc - the preconditioner context
886: Output Parameters:
887: + n - the number of subdomains for this processor (default value = 1)
888: - is - the index sets that define the subdomains for this processor
889:
891: Notes:
892: The IS numbering is in the parallel, global numbering of the vector.
894: Level: advanced
896: .keywords: PC, ASM, set, local, subdomains, additive Schwarz
898: .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubSLES(),
899: PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains()
900: @*/
901: int PCASMGetLocalSubdomains(PC pc,int *n,IS **is)
902: {
903: PC_ASM *osm;
907: if (!pc->setupcalled) {
908: SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call after SLESSetUP() or PCSetUp().");
909: }
911: osm = (PC_ASM*)pc->data;
912: if (n) *n = osm->n_local_true;
913: if (is) *is = osm->is;
914: return(0);
915: }