Actual source code: asm.c
1: /*$Id: asm.c,v 1.131 2001/08/07 03:03:38 balay 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: PetscScalar 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: PetscScalar 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: if (osm->sles) {
332: for (i=0; i<osm->n_local_true; i++) {
333: SLESDestroy(osm->sles[i]);
334: }
335: }
336: if (osm->is_flg) {
337: for (i=0; i<osm->n_local_true; i++) {ISDestroy(osm->is[i]);}
338: PetscFree(osm->is);
339: }
340: if (osm->sles) {PetscFree(osm->sles);}
341: if (osm->scat) {PetscFree(osm->scat);}
342: if (osm->x) {PetscFree(osm->x);}
343: PetscFree(osm);
344: return(0);
345: }
347: static int PCSetFromOptions_ASM(PC pc)
348: {
349: PC_ASM *osm = (PC_ASM*)pc->data;
350: int blocks,ovl,ierr;
351: PetscTruth flg;
352: char buff[16],*type[] = {"basic","restrict","interpolate","none"};
355: PetscOptionsHead("Additive Schwarz options");
356: PetscOptionsInt("-pc_asm_blocks","Number of subdomains","PCASMSetTotalSubdomains",osm->n,&blocks,&flg);
357: if (flg) {PCASMSetTotalSubdomains(pc,blocks,PETSC_NULL); }
358: PetscOptionsInt("-pc_asm_overlap","Number of grid points overlap","PCASMSetOverlap",osm->overlap,&ovl,&flg);
359: if (flg) {PCASMSetOverlap(pc,ovl); }
360: PetscOptionsName("-pc_asm_in_place","Perform matrix factorization inplace","PCASMSetUseInPlace",&flg);
361: if (flg) {PCASMSetUseInPlace(pc); }
362: PetscOptionsEList("-pc_asm_type","Type of restriction/extension","PCASMSetType",type,4,"restrict",buff,15,&flg);
363: if (flg) {
364: PCASMType atype;
365: PetscTruth isbasic,isrestrict,isinterpolate,isnone;
367: PetscStrcmp(buff,type[0],&isbasic);
368: PetscStrcmp(buff,type[1],&isrestrict);
369: PetscStrcmp(buff,type[2],&isinterpolate);
370: PetscStrcmp(buff,type[3],&isnone);
372: if (isbasic) atype = PC_ASM_BASIC;
373: else if (isrestrict) atype = PC_ASM_RESTRICT;
374: else if (isinterpolate) atype = PC_ASM_INTERPOLATE;
375: else if (isnone) atype = PC_ASM_NONE;
376: else SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Unknown type");
377: PCASMSetType(pc,atype);
378: }
379: PetscOptionsTail();
380: return(0);
381: }
383: /*------------------------------------------------------------------------------------*/
385: EXTERN_C_BEGIN
386: int PCASMSetLocalSubdomains_ASM(PC pc,int n,IS *is)
387: {
388: PC_ASM *osm;
392: if (pc->setupcalled) {
393: SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"PCASMSetLocalSubdomains() should be called before calling PCSetup().");
394: }
396: if (n < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Each process must have 0 or more blocks");
397: osm = (PC_ASM*)pc->data;
398: osm->n_local_true = n;
399: osm->is = is;
400: return(0);
401: }
402: EXTERN_C_END
404: EXTERN_C_BEGIN
405: int PCASMSetTotalSubdomains_ASM(PC pc,int N,IS *is)
406: {
407: PC_ASM *osm;
408: int rank,size,ierr;
411: if (pc->setupcalled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"PCASMSetTotalSubdomains() should be called before PCSetup().");
413: if (is) SETERRQ(PETSC_ERR_SUP,"Use PCASMSetLocalSubdomains to set specific index setsn
414: they cannot be set globally yet.");
416: osm = (PC_ASM*)pc->data;
417: /*
418: Split the subdomains equally amoung all processors
419: */
420: MPI_Comm_rank(pc->comm,&rank);
421: MPI_Comm_size(pc->comm,&size);
422: osm->n_local_true = N/size + ((N % size) > rank);
423: if (osm->n_local_true < 0) {
424: SETERRQ(PETSC_ERR_SUP,"Each process must have 0 or more blocks");
425: }
426: osm->is = 0;
427: return(0);
428: }
429: EXTERN_C_END
431: EXTERN_C_BEGIN
432: int PCASMSetOverlap_ASM(PC pc,int ovl)
433: {
434: PC_ASM *osm;
437: if (ovl < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested");
439: osm = (PC_ASM*)pc->data;
440: osm->overlap = ovl;
441: return(0);
442: }
443: EXTERN_C_END
445: EXTERN_C_BEGIN
446: int PCASMSetType_ASM(PC pc,PCASMType type)
447: {
448: PC_ASM *osm;
451: osm = (PC_ASM*)pc->data;
452: osm->type = type;
453: return(0);
454: }
455: EXTERN_C_END
457: EXTERN_C_BEGIN
458: int PCASMGetSubSLES_ASM(PC pc,int *n_local,int *first_local,SLES **sles)
459: {
460: PC_ASM *jac = (PC_ASM*)pc->data;
461: int ierr;
464: if (jac->n_local_true < 0) {
465: SETERRQ(1,"Need to call PCSetUP() on PC (or SLESSetUp() on the outer SLES object) before calling here");
466: }
468: if (n_local) *n_local = jac->n_local_true;
469: if (first_local) {
470: MPI_Scan(&jac->n_local_true,first_local,1,MPI_INT,MPI_SUM,pc->comm);
471: *first_local -= jac->n_local_true;
472: }
473: *sles = jac->sles;
474: jac->same_local_solves = PETSC_FALSE; /* Assume that local solves are now different;
475: not necessarily true though! This flag is
476: used only for PCView_ASM() */
477: return(0);
478: }
479: EXTERN_C_END
481: EXTERN_C_BEGIN
482: int PCASMSetUseInPlace_ASM(PC pc)
483: {
484: PC_ASM *dir;
487: dir = (PC_ASM*)pc->data;
488: dir->inplace = PETSC_TRUE;
489: return(0);
490: }
491: EXTERN_C_END
493: /*----------------------------------------------------------------------------*/
494: /*@
495: PCASMSetUseInPlace - Tells the system to destroy the matrix after setup is done.
497: Collective on PC
499: Input Parameters:
500: . pc - the preconditioner context
502: Options Database Key:
503: . -pc_asm_in_place - Activates in-place factorization
505: Note:
506: PCASMSetUseInplace() can only be used with the KSP method KSPPREONLY, and
507: when the original matrix is not required during the Solve process.
508: This destroys the matrix, early thus, saving on memory usage.
510: Level: intermediate
512: .keywords: PC, set, factorization, direct, inplace, in-place, ASM
514: .seealso: PCILUSetUseInPlace(), PCLUSetUseInPlace ()
515: @*/
516: int PCASMSetUseInPlace(PC pc)
517: {
518: int ierr,(*f)(PC);
522: PetscObjectQueryFunction((PetscObject)pc,"PCASMSetUseInPlace_C",(void (**)(void))&f);
523: if (f) {
524: (*f)(pc);
525: }
526: return(0);
527: }
528: /*----------------------------------------------------------------------------*/
530: /*@C
531: PCASMSetLocalSubdomains - Sets the local subdomains (for this processor
532: only) for the additive Schwarz preconditioner.
534: Collective on PC
536: Input Parameters:
537: + pc - the preconditioner context
538: . n - the number of subdomains for this processor (default value = 1)
539: - is - the index sets that define the subdomains for this processor
540: (or PETSC_NULL for PETSc to determine subdomains)
542: Notes:
543: The IS numbering is in the parallel, global numbering of the vector.
545: By default the ASM preconditioner uses 1 block per processor.
547: These index sets cannot be destroyed until after completion of the
548: linear solves for which the ASM preconditioner is being used.
550: Use PCASMSetTotalSubdomains() to set the subdomains for all processors.
552: Level: advanced
554: .keywords: PC, ASM, set, local, subdomains, additive Schwarz
556: .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubSLES(),
557: PCASMCreateSubdomains2D(), PCASMGetLocalSubdomains()
558: @*/
559: int PCASMSetLocalSubdomains(PC pc,int n,IS *is)
560: {
561: int ierr,(*f)(PC,int,IS *);
565: PetscObjectQueryFunction((PetscObject)pc,"PCASMSetLocalSubdomains_C",(void (**)(void))&f);
566: if (f) {
567: (*f)(pc,n,is);
568: }
569: return(0);
570: }
572: /*@C
573: PCASMSetTotalSubdomains - Sets the subdomains for all processor for the
574: additive Schwarz preconditioner. Either all or no processors in the
575: PC communicator must call this routine, with the same index sets.
577: Collective on PC
579: Input Parameters:
580: + pc - the preconditioner context
581: . n - the number of subdomains for all processors
582: - is - the index sets that define the subdomains for all processor
583: (or PETSC_NULL for PETSc to determine subdomains)
585: Options Database Key:
586: To set the total number of subdomain blocks rather than specify the
587: index sets, use the option
588: . -pc_asm_blocks <blks> - Sets total blocks
590: Notes:
591: Currently you cannot use this to set the actual subdomains with the argument is.
593: By default the ASM preconditioner uses 1 block per processor.
595: These index sets cannot be destroyed until after completion of the
596: linear solves for which the ASM preconditioner is being used.
598: Use PCASMSetLocalSubdomains() to set local subdomains.
600: Level: advanced
602: .keywords: PC, ASM, set, total, global, subdomains, additive Schwarz
604: .seealso: PCASMSetLocalSubdomains(), PCASMSetOverlap(), PCASMGetSubSLES(),
605: PCASMCreateSubdomains2D()
606: @*/
607: int PCASMSetTotalSubdomains(PC pc,int N,IS *is)
608: {
609: int ierr,(*f)(PC,int,IS *);
613: PetscObjectQueryFunction((PetscObject)pc,"PCASMSetTotalSubdomains_C",(void (**)(void))&f);
614: if (f) {
615: (*f)(pc,N,is);
616: }
617: return(0);
618: }
620: /*@
621: PCASMSetOverlap - Sets the overlap between a pair of subdomains for the
622: additive Schwarz preconditioner. Either all or no processors in the
623: PC communicator must call this routine.
625: Collective on PC
627: Input Parameters:
628: + pc - the preconditioner context
629: - ovl - the amount of overlap between subdomains (ovl >= 0, default value = 1)
631: Options Database Key:
632: . -pc_asm_overlap <ovl> - Sets overlap
634: Notes:
635: By default the ASM preconditioner uses 1 block per processor. To use
636: multiple blocks per perocessor, see PCASMSetTotalSubdomains() and
637: PCASMSetLocalSubdomains() (and the option -pc_asm_blocks <blks>).
639: The overlap defaults to 1, so if one desires that no additional
640: overlap be computed beyond what may have been set with a call to
641: PCASMSetTotalSubdomains() or PCASMSetLocalSubdomains(), then ovl
642: must be set to be 0. In particular, if one does not explicitly set
643: the subdomains an application code, then all overlap would be computed
644: internally by PETSc, and using an overlap of 0 would result in an ASM
645: variant that is equivalent to the block Jacobi preconditioner.
647: Note that one can define initial index sets with any overlap via
648: PCASMSetTotalSubdomains() or PCASMSetLocalSubdomains(); the routine
649: PCASMSetOverlap() merely allows PETSc to extend that overlap further
650: if desired.
652: Level: intermediate
654: .keywords: PC, ASM, set, overlap
656: .seealso: PCASMSetTotalSubdomains(), PCASMSetLocalSubdomains(), PCASMGetSubSLES(),
657: PCASMCreateSubdomains2D(), PCASMGetLocalSubdomains()
658: @*/
659: int PCASMSetOverlap(PC pc,int ovl)
660: {
661: int ierr,(*f)(PC,int);
665: PetscObjectQueryFunction((PetscObject)pc,"PCASMSetOverlap_C",(void (**)(void))&f);
666: if (f) {
667: (*f)(pc,ovl);
668: }
669: return(0);
670: }
672: /*@
673: PCASMSetType - Sets the type of restriction and interpolation used
674: for local problems in the additive Schwarz method.
676: Collective on PC
678: Input Parameters:
679: + pc - the preconditioner context
680: - type - variant of ASM, one of
681: .vb
682: PC_ASM_BASIC - full interpolation and restriction
683: PC_ASM_RESTRICT - full restriction, local processor interpolation
684: PC_ASM_INTERPOLATE - full interpolation, local processor restriction
685: PC_ASM_NONE - local processor restriction and interpolation
686: .ve
688: Options Database Key:
689: $ -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type
691: Level: intermediate
693: .keywords: PC, ASM, set, type
695: .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubSLES(),
696: PCASMCreateSubdomains2D()
697: @*/
698: int PCASMSetType(PC pc,PCASMType type)
699: {
700: int ierr,(*f)(PC,PCASMType);
704: PetscObjectQueryFunction((PetscObject)pc,"PCASMSetType_C",(void (**)(void))&f);
705: if (f) {
706: (*f)(pc,type);
707: }
708: return(0);
709: }
711: /*@C
712: PCASMGetSubSLES - Gets the local SLES contexts for all blocks on
713: this processor.
714:
715: Collective on PC iff first_local is requested
717: Input Parameter:
718: . pc - the preconditioner context
720: Output Parameters:
721: + n_local - the number of blocks on this processor or PETSC_NULL
722: . first_local - the global number of the first block on this processor or PETSC_NULL,
723: all processors must request or all must pass PETSC_NULL
724: - sles - the array of SLES contexts
726: Note:
727: After PCASMGetSubSLES() the array of SLESes is not to be freed
729: Currently for some matrix implementations only 1 block per processor
730: is supported.
731:
732: You must call SLESSetUp() before calling PCASMGetSubSLES().
734: Level: advanced
736: .keywords: PC, ASM, additive Schwarz, get, sub, SLES, context
738: .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap(),
739: PCASMCreateSubdomains2D(),
740: @*/
741: int PCASMGetSubSLES(PC pc,int *n_local,int *first_local,SLES **sles)
742: {
743: int ierr,(*f)(PC,int*,int*,SLES **);
747: PetscObjectQueryFunction((PetscObject)pc,"PCASMGetSubSLES_C",(void (**)(void))&f);
748: if (f) {
749: (*f)(pc,n_local,first_local,sles);
750: } else {
751: SETERRQ(1,"Cannot get subsles for this type of PC");
752: }
754: return(0);
755: }
757: /* -------------------------------------------------------------------------------------*/
758: EXTERN_C_BEGIN
759: int PCCreate_ASM(PC pc)
760: {
761: int ierr;
762: PC_ASM *osm;
765: PetscNew(PC_ASM,&osm);
766: PetscLogObjectMemory(pc,sizeof(PC_ASM));
767: PetscMemzero(osm,sizeof(PC_ASM));
768: osm->n = PETSC_DECIDE;
769: osm->n_local = 0;
770: osm->n_local_true = PETSC_DECIDE;
771: osm->overlap = 1;
772: osm->is_flg = PETSC_FALSE;
773: osm->sles = 0;
774: osm->scat = 0;
775: osm->is = 0;
776: osm->mat = 0;
777: osm->pmat = 0;
778: osm->type = PC_ASM_RESTRICT;
779: osm->same_local_solves = PETSC_TRUE;
780: osm->inplace = PETSC_FALSE;
781: pc->data = (void*)osm;
783: pc->ops->apply = PCApply_ASM;
784: pc->ops->applytranspose = PCApplyTranspose_ASM;
785: pc->ops->setup = PCSetUp_ASM;
786: pc->ops->destroy = PCDestroy_ASM;
787: pc->ops->setfromoptions = PCSetFromOptions_ASM;
788: pc->ops->setuponblocks = PCSetUpOnBlocks_ASM;
789: pc->ops->view = PCView_ASM;
790: pc->ops->applyrichardson = 0;
792: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCASMSetLocalSubdomains_C","PCASMSetLocalSubdomains_ASM",
793: PCASMSetLocalSubdomains_ASM);
794: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCASMSetTotalSubdomains_C","PCASMSetTotalSubdomains_ASM",
795: PCASMSetTotalSubdomains_ASM);
796: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCASMSetOverlap_C","PCASMSetOverlap_ASM",
797: PCASMSetOverlap_ASM);
798: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCASMSetType_C","PCASMSetType_ASM",
799: PCASMSetType_ASM);
800: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCASMGetSubSLES_C","PCASMGetSubSLES_ASM",
801: PCASMGetSubSLES_ASM);
802: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCASMSetUseInPlace_C","PCASMSetUseInPlace_ASM",
803: PCASMSetUseInPlace_ASM);
804: return(0);
805: }
806: EXTERN_C_END
809: /*@
810: PCASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz
811: preconditioner for a two-dimensional problem on a regular grid.
813: Not Collective
815: Input Parameters:
816: + m, n - the number of mesh points in the x and y directions
817: . M, N - the number of subdomains in the x and y directions
818: . dof - degrees of freedom per node
819: - overlap - overlap in mesh lines
821: Output Parameters:
822: + Nsub - the number of subdomains created
823: - is - the array of index sets defining the subdomains
825: Note:
826: Presently PCAMSCreateSubdomains2d() is valid only for sequential
827: preconditioners. More general related routines are
828: PCASMSetTotalSubdomains() and PCASMSetLocalSubdomains().
830: Level: advanced
832: .keywords: PC, ASM, additive Schwarz, create, subdomains, 2D, regular grid
834: .seealso: PCASMSetTotalSubdomains(), PCASMSetLocalSubdomains(), PCASMGetSubSLES(),
835: PCASMSetOverlap()
836: @*/
837: int PCASMCreateSubdomains2D(int m,int n,int M,int N,int dof,int overlap,int *Nsub,IS **is)
838: {
839: int i,j,height,width,ystart,xstart,yleft,yright,xleft,xright,loc_outter;
840: int nidx,*idx,loc,ii,jj,ierr,count;
843: if (dof != 1) SETERRQ(PETSC_ERR_SUP," ");
845: *Nsub = N*M;
846: PetscMalloc((*Nsub)*sizeof(IS **),is);
847: ystart = 0;
848: loc_outter = 0;
849: for (i=0; i<N; i++) {
850: height = n/N + ((n % N) > i); /* height of subdomain */
851: if (height < 2) SETERRQ(1,"Too many N subdomains for mesh dimension n");
852: yleft = ystart - overlap; if (yleft < 0) yleft = 0;
853: yright = ystart + height + overlap; if (yright > n) yright = n;
854: xstart = 0;
855: for (j=0; j<M; j++) {
856: width = m/M + ((m % M) > j); /* width of subdomain */
857: if (width < 2) SETERRQ(1,"Too many M subdomains for mesh dimension m");
858: xleft = xstart - overlap; if (xleft < 0) xleft = 0;
859: xright = xstart + width + overlap; if (xright > m) xright = m;
860: nidx = (xright - xleft)*(yright - yleft);
861: PetscMalloc(nidx*sizeof(int),&idx);
862: loc = 0;
863: for (ii=yleft; ii<yright; ii++) {
864: count = m*ii + xleft;
865: for (jj=xleft; jj<xright; jj++) {
866: idx[loc++] = count++;
867: }
868: }
869: ISCreateGeneral(PETSC_COMM_SELF,nidx,idx,(*is)+loc_outter++);
870: PetscFree(idx);
871: xstart += width;
872: }
873: ystart += height;
874: }
875: for (i=0; i<*Nsub; i++) { ISSort((*is)[i]); }
876: return(0);
877: }
879: /*@C
880: PCASMGetLocalSubdomains - Gets the local subdomains (for this processor
881: only) for the additive Schwarz preconditioner.
883: Collective on PC
885: Input Parameter:
886: . pc - the preconditioner context
888: Output Parameters:
889: + n - the number of subdomains for this processor (default value = 1)
890: - is - the index sets that define the subdomains for this processor
891:
893: Notes:
894: The IS numbering is in the parallel, global numbering of the vector.
896: Level: advanced
898: .keywords: PC, ASM, set, local, subdomains, additive Schwarz
900: .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubSLES(),
901: PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubmatrices()
902: @*/
903: int PCASMGetLocalSubdomains(PC pc,int *n,IS **is)
904: {
905: PC_ASM *osm;
909: if (!pc->setupcalled) {
910: SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call after SLESSetUP() or PCSetUp().");
911: }
913: osm = (PC_ASM*)pc->data;
914: if (n) *n = osm->n_local_true;
915: if (is) *is = osm->is;
916: return(0);
917: }
919: /*@C
920: PCASMGetLocalSubmatrices - Gets the local submatrices (for this processor
921: only) for the additive Schwarz preconditioner.
923: Collective on PC
925: Input Parameter:
926: . pc - the preconditioner context
928: Output Parameters:
929: + n - the number of matrices for this processor (default value = 1)
930: - mat - the matrices
931:
933: Level: advanced
935: .keywords: PC, ASM, set, local, subdomains, additive Schwarz, block Jacobi
937: .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubSLES(),
938: PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains()
939: @*/
940: int PCASMGetLocalSubmatrices(PC pc,int *n,Mat **mat)
941: {
942: PC_ASM *osm;
946: if (!pc->setupcalled) {
947: SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call after SLESSetUP() or PCSetUp().");
948: }
950: osm = (PC_ASM*)pc->data;
951: if (n) *n = osm->n_local_true;
952: if (mat) *mat = osm->pmat;
953: return(0);
954: }