xref: /petsc/src/ksp/pc/impls/gasm/gasm.c (revision c10bc1a15c3326c73cfa28ba5063292eba8d33d9)
1b1a0a8a3SJed Brown #define PETSCKSP_DLL
2b1a0a8a3SJed Brown 
3b1a0a8a3SJed Brown /*
4f746d493SDmitry Karpeev   This file defines an "generalized" additive Schwarz preconditioner for any Mat implementation.
544fe18e5SDmitry Karpeev   In this version each processor may have any number of subdomains and a subdomain may live on multiple
6f746d493SDmitry Karpeev   processors.
7b1a0a8a3SJed Brown 
8f746d493SDmitry Karpeev        N    - total number of true subdomains on all processors
9f746d493SDmitry Karpeev        n    - actual number of subdomains on this processor
10f746d493SDmitry Karpeev        nmax - maximum number of subdomains per processor
11b1a0a8a3SJed Brown */
12b1a0a8a3SJed Brown #include "private/pcimpl.h"     /*I "petscpc.h" I*/
13b1a0a8a3SJed Brown 
14b1a0a8a3SJed Brown typedef struct {
15f746d493SDmitry Karpeev   PetscInt   N,n,nmax;
16b1a0a8a3SJed Brown   PetscInt   overlap;             /* overlap requested by user */
17b1a0a8a3SJed Brown   KSP        *ksp;                /* linear solvers for each block */
18f746d493SDmitry Karpeev   Vec        gx,gy;               /* Merged work vectors */
19f746d493SDmitry Karpeev   Vec        *x,*y;               /* Split work vectors; storage aliases pieces of storage of the above merged vectors. */
20f746d493SDmitry Karpeev   IS         gis, gis_local;      /* merged ISs */
21f746d493SDmitry Karpeev   VecScatter grestriction;        /* merged restriction */
22f746d493SDmitry Karpeev   VecScatter gprolongation;       /* merged prolongation */
23b1a0a8a3SJed Brown   IS         *is;                 /* index set that defines each overlapping subdomain */
24f746d493SDmitry Karpeev   IS         *is_local;           /* index set that defines each local subdomain (same as subdomain with the overlap removed); may be NULL */
25f746d493SDmitry Karpeev   Mat        *pmat;               /* subdomain block matrices */
26f746d493SDmitry Karpeev   PCGASMType  type;               /* use reduced interpolation, restriction or both */
27b1a0a8a3SJed Brown   PetscBool  type_set;            /* if user set this value (so won't change it for symmetric problems) */
28b1a0a8a3SJed Brown   PetscBool  same_local_solves;   /* flag indicating whether all local solvers are same */
29b1a0a8a3SJed Brown   PetscBool  sort_indices;        /* flag to sort subdomain indices */
30f746d493SDmitry Karpeev } PC_GASM;
31b1a0a8a3SJed Brown 
32b1a0a8a3SJed Brown #undef __FUNCT__
33f746d493SDmitry Karpeev #define __FUNCT__ "PCView_GASM"
34f746d493SDmitry Karpeev static PetscErrorCode PCView_GASM(PC pc,PetscViewer viewer)
35b1a0a8a3SJed Brown {
36f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
37b1a0a8a3SJed Brown   PetscErrorCode ierr;
38b1a0a8a3SJed Brown   PetscMPIInt    rank;
39b1a0a8a3SJed Brown   PetscInt       i,bsz;
40b1a0a8a3SJed Brown   PetscBool      iascii,isstring;
41b1a0a8a3SJed Brown   PetscViewer    sviewer;
42b1a0a8a3SJed Brown 
43b1a0a8a3SJed Brown 
44b1a0a8a3SJed Brown   PetscFunctionBegin;
45b1a0a8a3SJed Brown   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
46b1a0a8a3SJed Brown   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
47b1a0a8a3SJed Brown   if (iascii) {
48b1a0a8a3SJed Brown     char overlaps[256] = "user-defined overlap",blocks[256] = "total subdomain blocks not yet set";
49b1a0a8a3SJed Brown     if (osm->overlap >= 0) {ierr = PetscSNPrintf(overlaps,sizeof overlaps,"amount of overlap = %D",osm->overlap);CHKERRQ(ierr);}
5044fe18e5SDmitry Karpeev     if (osm->nmax > 0)     {ierr = PetscSNPrintf(blocks,sizeof blocks,"max number of  subdomain blocks = %D",osm->nmax);CHKERRQ(ierr);}
51f746d493SDmitry Karpeev     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"  [%d] number of local blocks = %D\n",(int)rank,osm->n);CHKERRQ(ierr);
52f746d493SDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"  Generalized additive Schwarz: %s, %s\n",blocks,overlaps);CHKERRQ(ierr);
53f746d493SDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"  Generalized additive Schwarz: restriction/interpolation type - %s\n",PCGASMTypes[osm->type]);CHKERRQ(ierr);
54b1a0a8a3SJed Brown     ierr = MPI_Comm_rank(((PetscObject)pc)->comm,&rank);CHKERRQ(ierr);
55b1a0a8a3SJed Brown     if (osm->same_local_solves) {
56b1a0a8a3SJed Brown       if (osm->ksp) {
57b1a0a8a3SJed Brown         ierr = PetscViewerASCIIPrintf(viewer,"  Local solve is same for all blocks, in the following KSP and PC objects:\n");CHKERRQ(ierr);
58b1a0a8a3SJed Brown         ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
59b1a0a8a3SJed Brown         if (!rank) {
60b1a0a8a3SJed Brown           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
61b1a0a8a3SJed Brown           ierr = KSPView(osm->ksp[0],sviewer);CHKERRQ(ierr);
62b1a0a8a3SJed Brown           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
63b1a0a8a3SJed Brown         }
64b1a0a8a3SJed Brown         ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
65b1a0a8a3SJed Brown       }
66b1a0a8a3SJed Brown     } else {
67b1a0a8a3SJed Brown       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
68b1a0a8a3SJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  Local solve info for each block is in the following KSP and PC objects:\n");CHKERRQ(ierr);
69b1a0a8a3SJed Brown       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
70b1a0a8a3SJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr);
71f746d493SDmitry Karpeev       for (i=0; i<osm->nmax; i++) {
72b1a0a8a3SJed Brown         ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
73f746d493SDmitry Karpeev         if (i < osm->n) {
74b1a0a8a3SJed Brown           ierr = ISGetLocalSize(osm->is[i],&bsz);CHKERRQ(ierr);
75b1a0a8a3SJed Brown           ierr = PetscViewerASCIISynchronizedPrintf(sviewer,"[%d] local block number %D, size = %D\n",(int)rank,i,bsz);CHKERRQ(ierr);
76b1a0a8a3SJed Brown           ierr = KSPView(osm->ksp[i],sviewer);CHKERRQ(ierr);
77b1a0a8a3SJed Brown           ierr = PetscViewerASCIISynchronizedPrintf(sviewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr);
78b1a0a8a3SJed Brown         }
79b1a0a8a3SJed Brown         ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
80b1a0a8a3SJed Brown       }
81b1a0a8a3SJed Brown       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
82b1a0a8a3SJed Brown       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
83b1a0a8a3SJed Brown     }
84b1a0a8a3SJed Brown   } else if (isstring) {
85f746d493SDmitry Karpeev     ierr = PetscViewerStringSPrintf(viewer," blocks=%D, overlap=%D, type=%s",osm->n,osm->overlap,PCGASMTypes[osm->type]);CHKERRQ(ierr);
86b1a0a8a3SJed Brown     ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
87b1a0a8a3SJed Brown     if (osm->ksp) {ierr = KSPView(osm->ksp[0],sviewer);CHKERRQ(ierr);}
88b1a0a8a3SJed Brown     ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
89b1a0a8a3SJed Brown   } else {
90f746d493SDmitry Karpeev     SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCGASM",((PetscObject)viewer)->type_name);
91b1a0a8a3SJed Brown   }
92b1a0a8a3SJed Brown   PetscFunctionReturn(0);
93b1a0a8a3SJed Brown }
94b1a0a8a3SJed Brown 
95b1a0a8a3SJed Brown #undef __FUNCT__
96f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMPrintSubdomains"
97f746d493SDmitry Karpeev static PetscErrorCode PCGASMPrintSubdomains(PC pc)
98b1a0a8a3SJed Brown {
99f746d493SDmitry Karpeev   PC_GASM         *osm  = (PC_GASM*)pc->data;
100b1a0a8a3SJed Brown   const char     *prefix;
101b1a0a8a3SJed Brown   char           fname[PETSC_MAX_PATH_LEN+1];
102b1a0a8a3SJed Brown   PetscViewer    viewer;
103b1a0a8a3SJed Brown   PetscInt       i,j,nidx;
104b1a0a8a3SJed Brown   const PetscInt *idx;
105b1a0a8a3SJed Brown   PetscErrorCode ierr;
106b1a0a8a3SJed Brown 
107b1a0a8a3SJed Brown   PetscFunctionBegin;
108b1a0a8a3SJed Brown   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
109f746d493SDmitry Karpeev   ierr = PetscOptionsGetString(prefix,"-pc_gasm_print_subdomains",fname,PETSC_MAX_PATH_LEN,PETSC_NULL);CHKERRQ(ierr);
110b1a0a8a3SJed Brown   if (fname[0] == 0) { ierr = PetscStrcpy(fname,"stdout");CHKERRQ(ierr); };
111b1a0a8a3SJed Brown   ierr = PetscViewerASCIIOpen(((PetscObject)pc)->comm,fname,&viewer);CHKERRQ(ierr);
112f746d493SDmitry Karpeev   for (i=0;i<osm->n;++i) {
113b1a0a8a3SJed Brown     ierr = ISGetLocalSize(osm->is[i],&nidx);CHKERRQ(ierr);
114b1a0a8a3SJed Brown     ierr = ISGetIndices(osm->is[i],&idx);CHKERRQ(ierr);
115b1a0a8a3SJed Brown     for (j=0; j<nidx; j++) {
116b1a0a8a3SJed Brown       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"%D ",idx[j]);CHKERRQ(ierr);
117b1a0a8a3SJed Brown     }
118b1a0a8a3SJed Brown     ierr = ISRestoreIndices(osm->is[i],&idx);CHKERRQ(ierr);
119b1a0a8a3SJed Brown     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr);
120b1a0a8a3SJed Brown     if (osm->is_local) {
121b1a0a8a3SJed Brown       ierr = ISGetLocalSize(osm->is_local[i],&nidx);CHKERRQ(ierr);
122b1a0a8a3SJed Brown       ierr = ISGetIndices(osm->is_local[i],&idx);CHKERRQ(ierr);
123b1a0a8a3SJed Brown       for (j=0; j<nidx; j++) {
124b1a0a8a3SJed Brown         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"%D ",idx[j]);CHKERRQ(ierr);
125b1a0a8a3SJed Brown       }
126b1a0a8a3SJed Brown       ierr = ISRestoreIndices(osm->is_local[i],&idx);CHKERRQ(ierr);
127b1a0a8a3SJed Brown       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr);
128b1a0a8a3SJed Brown     }
129b1a0a8a3SJed Brown   }
130b1a0a8a3SJed Brown   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
131b1a0a8a3SJed Brown   ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr);
132b1a0a8a3SJed Brown   PetscFunctionReturn(0);
133b1a0a8a3SJed Brown }
134b1a0a8a3SJed Brown 
135b1a0a8a3SJed Brown #undef __FUNCT__
136f746d493SDmitry Karpeev #define __FUNCT__ "PCSetUp_GASM"
137f746d493SDmitry Karpeev static PetscErrorCode PCSetUp_GASM(PC pc)
138b1a0a8a3SJed Brown {
139f746d493SDmitry Karpeev   PC_GASM         *osm  = (PC_GASM*)pc->data;
140b1a0a8a3SJed Brown   PetscErrorCode ierr;
141b1a0a8a3SJed Brown   PetscBool      symset,flg;
142f746d493SDmitry Karpeev   PetscInt       i,firstRow,lastRow;
143*c10bc1a1SDmitry Karpeev   PetscMPIInt    rank, size;
144b1a0a8a3SJed Brown   MatReuse       scall = MAT_REUSE_MATRIX;
145b1a0a8a3SJed Brown   KSP            ksp;
146b1a0a8a3SJed Brown   PC             subpc;
147b1a0a8a3SJed Brown   const char     *prefix,*pprefix;
148f746d493SDmitry Karpeev   PetscInt       dn;       /* Number of indices in a single subdomain assigned to this processor. */
149f746d493SDmitry Karpeev   const PetscInt *didx;    /* Indices from a single subdomain assigned to this processor. */
150f746d493SDmitry Karpeev   PetscInt       ddn;      /* Number of indices in all subdomains assigned to this processor. */
151f746d493SDmitry Karpeev   PetscInt       *ddidx;   /* Indices of all subdomains assigned to this processor. */
152f746d493SDmitry Karpeev   IS             gid;      /* Identity IS of the size of all subdomains assigned to this processor. */
153f746d493SDmitry Karpeev   Vec            x,y;
154f746d493SDmitry Karpeev   PetscScalar    *gxarray, *gyarray;
155b1a0a8a3SJed Brown 
156b1a0a8a3SJed Brown   PetscFunctionBegin;
157*c10bc1a1SDmitry Karpeev   ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr);
158*c10bc1a1SDmitry Karpeev   ierr = MPI_Comm_rank(((PetscObject)pc)->comm,&rank);CHKERRQ(ierr);
159b1a0a8a3SJed Brown   if (!pc->setupcalled) {
160b1a0a8a3SJed Brown 
161b1a0a8a3SJed Brown     if (!osm->type_set) {
162b1a0a8a3SJed Brown       ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr);
163f746d493SDmitry Karpeev       if (symset && flg) { osm->type = PC_GASM_BASIC; }
164b1a0a8a3SJed Brown     }
165b1a0a8a3SJed Brown 
16644fe18e5SDmitry Karpeev     if (osm->N == PETSC_DECIDE && osm->n < 1) {
167b1a0a8a3SJed Brown       /* no subdomains given, use one per processor */
16844fe18e5SDmitry Karpeev       osm->nmax = osm->n = 1;
169b1a0a8a3SJed Brown       ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr);
170f746d493SDmitry Karpeev       osm->N = size;
171f746d493SDmitry Karpeev     } else if (osm->N == PETSC_DECIDE) {
172b1a0a8a3SJed Brown       PetscInt inwork[2], outwork[2];
173f746d493SDmitry Karpeev       /* determine global number of subdomains and the max number of local subdomains */
174f746d493SDmitry Karpeev       inwork[0] = inwork[1] = osm->n;
175b1a0a8a3SJed Brown       ierr = MPI_Allreduce(inwork,outwork,1,MPIU_2INT,PetscMaxSum_Op,((PetscObject)pc)->comm);CHKERRQ(ierr);
176f746d493SDmitry Karpeev       osm->nmax = outwork[0];
177f746d493SDmitry Karpeev       osm->N    = outwork[1];
178b1a0a8a3SJed Brown     }
179b1a0a8a3SJed Brown     if (!osm->is){ /* create the index sets */
180f746d493SDmitry Karpeev       ierr = PCGASMCreateSubdomains(pc->pmat,osm->n,&osm->is);CHKERRQ(ierr);
181b1a0a8a3SJed Brown     }
18293c1ec32SDmitry Karpeev     if (!osm->is_local) {
18393c1ec32SDmitry Karpeev       /*
18493c1ec32SDmitry Karpeev 	 This indicates that osm->is should define a nonoverlapping decomposition
18593c1ec32SDmitry Karpeev 	 (there is no way to really guarantee that if subdomains are set by the user through PCGASMSetLocalSubdomains,
18693c1ec32SDmitry Karpeev 	  but the assumption is that either the user does the right thing, or subdomains in ossm->is have been created
18793c1ec32SDmitry Karpeev 	  via PCGASMCreateSubdomains, which guarantees a nonoverlapping decomposition).
18893c1ec32SDmitry Karpeev 	 Therefore, osm->is will be used to define osm->is_local.
18993c1ec32SDmitry Karpeev 	 If a nonzero overlap has been requested by the user, then osm->is will be expanded and will overlap,
19093c1ec32SDmitry Karpeev 	 so osm->is_local should obtain a copy of osm->is while they are still (presumably) nonoverlapping.
19193c1ec32SDmitry Karpeev 	 Otherwise (no overlap has been requested), osm->is_local are simply aliases for osm->is.
19293c1ec32SDmitry Karpeev       */
193f746d493SDmitry Karpeev       ierr = PetscMalloc(osm->n*sizeof(IS),&osm->is_local);CHKERRQ(ierr);
194f746d493SDmitry Karpeev       for (i=0; i<osm->n; i++) {
195b1a0a8a3SJed Brown         if (osm->overlap > 0) { /* With positive overlap, osm->is[i] will be modified */
196b1a0a8a3SJed Brown           ierr = ISDuplicate(osm->is[i],&osm->is_local[i]);CHKERRQ(ierr);
197b1a0a8a3SJed Brown           ierr = ISCopy(osm->is[i],osm->is_local[i]);CHKERRQ(ierr);
198b1a0a8a3SJed Brown         } else {
199b1a0a8a3SJed Brown           ierr = PetscObjectReference((PetscObject)osm->is[i]);CHKERRQ(ierr);
200b1a0a8a3SJed Brown           osm->is_local[i] = osm->is[i];
201b1a0a8a3SJed Brown         }
202b1a0a8a3SJed Brown       }
203b1a0a8a3SJed Brown     }
20493c1ec32SDmitry Karpeev     /* Beyond this point osm->is_local is not null. */
20593c1ec32SDmitry Karpeev     if (osm->overlap > 0) {
20693c1ec32SDmitry Karpeev       /* Extend the "overlapping" regions by a number of steps */
20793c1ec32SDmitry Karpeev       ierr = MatIncreaseOverlap(pc->pmat,osm->n,osm->is,osm->overlap);CHKERRQ(ierr);
20893c1ec32SDmitry Karpeev     }
209b1a0a8a3SJed Brown     ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
210b1a0a8a3SJed Brown     flg  = PETSC_FALSE;
211f746d493SDmitry Karpeev     ierr = PetscOptionsGetBool(prefix,"-pc_gasm_print_subdomains",&flg,PETSC_NULL);CHKERRQ(ierr);
212f746d493SDmitry Karpeev     if (flg) { ierr = PCGASMPrintSubdomains(pc);CHKERRQ(ierr); }
213b1a0a8a3SJed Brown 
214b1a0a8a3SJed Brown     if (osm->sort_indices) {
215f746d493SDmitry Karpeev       for (i=0; i<osm->n; i++) {
216b1a0a8a3SJed Brown         ierr = ISSort(osm->is[i]);CHKERRQ(ierr);
217b1a0a8a3SJed Brown 	ierr = ISSort(osm->is_local[i]);CHKERRQ(ierr);
218b1a0a8a3SJed Brown       }
219b1a0a8a3SJed Brown     }
220f746d493SDmitry Karpeev     /* Merge the ISs, create merged vectors and scatter contexts. */
22193c1ec32SDmitry Karpeev     /* Restriction ISs. */
222f746d493SDmitry Karpeev     ddn = 0;
223f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
224f746d493SDmitry Karpeev       ierr = ISGetLocalSize(osm->is[i],&dn);          CHKERRQ(ierr);
225f746d493SDmitry Karpeev       ddn += dn;
226f746d493SDmitry Karpeev     }
227f746d493SDmitry Karpeev     ierr = PetscMalloc(ddn*sizeof(PetscInt), &ddidx); CHKERRQ(ierr);
228f746d493SDmitry Karpeev     ddn = 0;
229f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
230f746d493SDmitry Karpeev       ierr = ISGetLocalSize(osm->is[i],&dn); CHKERRQ(ierr);
231f746d493SDmitry Karpeev       ierr = ISGetIndices(osm->is[i],&didx); CHKERRQ(ierr);
232f746d493SDmitry Karpeev       ierr = PetscMemcpy(ddidx+ddn, didx, sizeof(PetscInt)*dn); CHKERRQ(ierr);
233f746d493SDmitry Karpeev       ierr = ISRestoreIndices(osm->is[i], &didx); CHKERRQ(ierr);
23493c1ec32SDmitry Karpeev       ddn += dn;
235f746d493SDmitry Karpeev     }
236f746d493SDmitry Karpeev     ierr = ISCreateGeneral(((PetscObject)(pc))->comm, ddn, ddidx, PETSC_OWN_POINTER, &osm->gis); CHKERRQ(ierr);
237f746d493SDmitry Karpeev     ierr = MatGetVecs(pc->pmat,&x,&y);CHKERRQ(ierr);
238f746d493SDmitry Karpeev     ierr = VecCreateMPI(((PetscObject)pc)->comm, ddn, PETSC_DECIDE, &osm->gx); CHKERRQ(ierr);
239f746d493SDmitry Karpeev     ierr = VecDuplicate(osm->gx,&osm->gy);                                     CHKERRQ(ierr);
240f746d493SDmitry Karpeev     ierr = ISCreateStride(((PetscObject)pc)->comm,ddn,0,1, &gid);              CHKERRQ(ierr);
241f746d493SDmitry Karpeev     ierr = VecScatterCreate(x,osm->gis,osm->gx,gid, &(osm->grestriction));   CHKERRQ(ierr);
242f746d493SDmitry Karpeev     ierr = ISDestroy(gid);                                                     CHKERRQ(ierr);
24393c1ec32SDmitry Karpeev     /* Prolongation ISs */
24493c1ec32SDmitry Karpeev     { PetscInt       dn_local;       /* Number of indices in the local part of single domain assigned to this processor. */
245f746d493SDmitry Karpeev       const PetscInt *didx_local;    /* Global indices from the local part of a single domain assigned to this processor. */
246f746d493SDmitry Karpeev       PetscInt       ddn_local;      /* Number of indices in the local part of the disjoint union all domains assigned to this processor. */
247f746d493SDmitry Karpeev       PetscInt       *ddidx_local;   /* Global indices of the local part of the disjoint union of all domains assigned to this processor. */
248f746d493SDmitry Karpeev       /**/
24993c1ec32SDmitry Karpeev       ISLocalToGlobalMapping ltog;          /* Map from global to local indices on the disjoint union of subdomains: "local" ind's run from 0 to ddn-1. */
25093c1ec32SDmitry Karpeev       PetscInt              *ddidx_llocal;  /* Mapped local indices of the disjoint union of local parts of subdomains. */
251f746d493SDmitry Karpeev       PetscInt               ddn_llocal;    /* Number of indices in ddidx_llocal; must equal ddn_local, or else gis_local is not a sub-IS of gis. */
252f746d493SDmitry Karpeev       IS                     gis_llocal;    /* IS with ddidx_llocal indices. */
253f746d493SDmitry Karpeev       PetscInt               j;
254f746d493SDmitry Karpeev       ddn_local = 0;
255f746d493SDmitry Karpeev       for (i=0; i<osm->n; i++) {
256f746d493SDmitry Karpeev 	ierr = ISGetLocalSize(osm->is_local[i],&dn_local); CHKERRQ(ierr);
257f746d493SDmitry Karpeev 	ddn_local += dn_local;
258f746d493SDmitry Karpeev       }
259f746d493SDmitry Karpeev       ierr = PetscMalloc(ddn_local*sizeof(PetscInt), &ddidx_local); CHKERRQ(ierr);
260f746d493SDmitry Karpeev       ddn_local = 0;
261f746d493SDmitry Karpeev       for (i=0; i<osm->n; i++) {
262f746d493SDmitry Karpeev 	ierr = ISGetLocalSize(osm->is_local[i],&dn_local); CHKERRQ(ierr);
263f746d493SDmitry Karpeev 	ierr = ISGetIndices(osm->is_local[i],&didx_local); CHKERRQ(ierr);
264f746d493SDmitry Karpeev 	ierr = PetscMemcpy(ddidx_local+ddn_local, didx_local, sizeof(PetscInt)*dn_local); CHKERRQ(ierr);
265f746d493SDmitry Karpeev 	ierr = ISRestoreIndices(osm->is_local[i], &didx_local); CHKERRQ(ierr);
266f746d493SDmitry Karpeev 	ddn_local += dn_local;
267f746d493SDmitry Karpeev       }
268ab3e923fSDmitry Karpeev       ierr = PetscMalloc(sizeof(PetscInt)*ddn_local, &ddidx_llocal); CHKERRQ(ierr);
269f746d493SDmitry Karpeev       ierr = ISCreateGeneral(((PetscObject)pc)->comm, ddn_local, ddidx_local, PETSC_OWN_POINTER, &(osm->gis_local)); CHKERRQ(ierr);
270f746d493SDmitry Karpeev       ierr = ISLocalToGlobalMappingCreateIS(osm->gis,&ltog);CHKERRQ(ierr);
271f746d493SDmitry Karpeev       ierr = ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,dn_local,ddidx_local,&ddn_llocal,ddidx_llocal);CHKERRQ(ierr);
272b1a0a8a3SJed Brown       ierr = ISLocalToGlobalMappingDestroy(ltog);CHKERRQ(ierr);
27393c1ec32SDmitry Karpeev       if (ddn_llocal != ddn_local) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"gis_local contains %D indices outside of gis", ddn_llocal - ddn_local);
274f746d493SDmitry Karpeev       /* Now convert these localized indices into the global indices into the merged output vector. */
275f746d493SDmitry Karpeev       ierr = VecGetOwnershipRange(osm->gy, &firstRow, &lastRow); CHKERRQ(ierr);
276f746d493SDmitry Karpeev       for(j=0; j < ddn_llocal; ++j) {
277f746d493SDmitry Karpeev 	ddidx_llocal[j] += firstRow;
278f746d493SDmitry Karpeev       }
279f746d493SDmitry Karpeev       ierr = ISCreateGeneral(PETSC_COMM_SELF,ddn_llocal,ddidx_llocal,PETSC_OWN_POINTER,&gis_llocal); CHKERRQ(ierr);
280f746d493SDmitry Karpeev       ierr = VecScatterCreate(y,osm->gis_local,osm->gy,gis_llocal,&osm->gprolongation);              CHKERRQ(ierr);
281f746d493SDmitry Karpeev       ierr = ISDestroy(gis_llocal);CHKERRQ(ierr);
282b1a0a8a3SJed Brown     }
283f746d493SDmitry Karpeev     /* Create the subdomain work vectors. */
284f746d493SDmitry Karpeev     ierr = PetscMalloc(osm->n*sizeof(Vec),&osm->x);CHKERRQ(ierr);
285f746d493SDmitry Karpeev     ierr = PetscMalloc(osm->n*sizeof(Vec),&osm->y);CHKERRQ(ierr);
286f746d493SDmitry Karpeev     ierr = VecGetOwnershipRange(osm->gx, &firstRow, &lastRow);CHKERRQ(ierr);
287f746d493SDmitry Karpeev     ierr = VecGetArray(osm->gx, &gxarray);                                     CHKERRQ(ierr);
288f746d493SDmitry Karpeev     ierr = VecGetArray(osm->gy, &gyarray);                                     CHKERRQ(ierr);
289f746d493SDmitry Karpeev     for (i=0, ddn=0; i<osm->n; ++i, ddn += dn) {
290f746d493SDmitry Karpeev       ierr = ISGetLocalSize(osm->is[i],&dn);CHKERRQ(ierr);
291f746d493SDmitry Karpeev       ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,dn,gxarray+ddn,&osm->x[i]);CHKERRQ(ierr);
292f746d493SDmitry Karpeev       ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,dn,gyarray+ddn,&osm->y[i]);CHKERRQ(ierr);
293b1a0a8a3SJed Brown     }
294f746d493SDmitry Karpeev     ierr = VecRestoreArray(osm->gx, &gxarray); CHKERRQ(ierr);
295f746d493SDmitry Karpeev     ierr = VecRestoreArray(osm->gy, &gyarray); CHKERRQ(ierr);
296f746d493SDmitry Karpeev     ierr = VecDestroy(x); CHKERRQ(ierr);
297f746d493SDmitry Karpeev     ierr = VecDestroy(y); CHKERRQ(ierr);
298b1a0a8a3SJed Brown     /* Create the local solvers */
299f746d493SDmitry Karpeev     ierr = PetscMalloc(osm->n*sizeof(KSP *),&osm->ksp);CHKERRQ(ierr);
30044fe18e5SDmitry Karpeev     for (i=0; i<osm->n; i++) { /* KSPs are local */
301b1a0a8a3SJed Brown       ierr = KSPCreate(PETSC_COMM_SELF,&ksp);CHKERRQ(ierr);
302b1a0a8a3SJed Brown       ierr = PetscLogObjectParent(pc,ksp);CHKERRQ(ierr);
303b1a0a8a3SJed Brown       ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);CHKERRQ(ierr);
304b1a0a8a3SJed Brown       ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr);
305b1a0a8a3SJed Brown       ierr = KSPGetPC(ksp,&subpc);CHKERRQ(ierr);
306b1a0a8a3SJed Brown       ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
307b1a0a8a3SJed Brown       ierr = KSPSetOptionsPrefix(ksp,prefix);CHKERRQ(ierr);
308b1a0a8a3SJed Brown       ierr = KSPAppendOptionsPrefix(ksp,"sub_");CHKERRQ(ierr);
309b1a0a8a3SJed Brown       osm->ksp[i] = ksp;
310b1a0a8a3SJed Brown     }
311b1a0a8a3SJed Brown     scall = MAT_INITIAL_MATRIX;
312b1a0a8a3SJed Brown 
313b1a0a8a3SJed Brown   } else {
314b1a0a8a3SJed Brown     /*
315b1a0a8a3SJed Brown        Destroy the blocks from the previous iteration
316b1a0a8a3SJed Brown     */
317b1a0a8a3SJed Brown     if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
318f746d493SDmitry Karpeev       ierr = MatDestroyMatrices(osm->n,&osm->pmat);CHKERRQ(ierr);
319b1a0a8a3SJed Brown       scall = MAT_INITIAL_MATRIX;
320b1a0a8a3SJed Brown     }
321b1a0a8a3SJed Brown   }
322b1a0a8a3SJed Brown 
323b1a0a8a3SJed Brown   /*
324f746d493SDmitry Karpeev      Extract out the submatrices.
325b1a0a8a3SJed Brown   */
326*c10bc1a1SDmitry Karpeev   ierr = MatGetSubMatrices(pc->pmat,osm->n,osm->is, osm->is,scall,&osm->pmat);CHKERRQ(ierr);
327b1a0a8a3SJed Brown   if (scall == MAT_INITIAL_MATRIX) {
328b1a0a8a3SJed Brown     ierr = PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix);CHKERRQ(ierr);
329f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
330b1a0a8a3SJed Brown       ierr = PetscLogObjectParent(pc,osm->pmat[i]);CHKERRQ(ierr);
331b1a0a8a3SJed Brown       ierr = PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i],pprefix);CHKERRQ(ierr);
332b1a0a8a3SJed Brown     }
333b1a0a8a3SJed Brown   }
334b1a0a8a3SJed Brown 
335b1a0a8a3SJed Brown   /* Return control to the user so that the submatrices can be modified (e.g., to apply
336b1a0a8a3SJed Brown      different boundary conditions for the submatrices than for the global problem) */
337*c10bc1a1SDmitry Karpeev   ierr = PCModifySubMatrices(pc,osm->n,osm->is,osm->is,osm->pmat,pc->modifysubmatricesP);CHKERRQ(ierr);
338b1a0a8a3SJed Brown 
339b1a0a8a3SJed Brown   /*
340f746d493SDmitry Karpeev      Loop over submatrices putting them into local ksp
341b1a0a8a3SJed Brown   */
342f746d493SDmitry Karpeev   for (i=0; i<osm->n; i++) {
343b1a0a8a3SJed Brown     ierr = KSPSetOperators(osm->ksp[i],osm->pmat[i],osm->pmat[i],pc->flag);CHKERRQ(ierr);
344b1a0a8a3SJed Brown     if (!pc->setupcalled) {
345b1a0a8a3SJed Brown       ierr = KSPSetFromOptions(osm->ksp[i]);CHKERRQ(ierr);
346b1a0a8a3SJed Brown     }
347b1a0a8a3SJed Brown   }
348b1a0a8a3SJed Brown 
349b1a0a8a3SJed Brown   PetscFunctionReturn(0);
350b1a0a8a3SJed Brown }
351b1a0a8a3SJed Brown 
352b1a0a8a3SJed Brown #undef __FUNCT__
353f746d493SDmitry Karpeev #define __FUNCT__ "PCSetUpOnBlocks_GASM"
354f746d493SDmitry Karpeev static PetscErrorCode PCSetUpOnBlocks_GASM(PC pc)
355b1a0a8a3SJed Brown {
356f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
357b1a0a8a3SJed Brown   PetscErrorCode ierr;
358b1a0a8a3SJed Brown   PetscInt       i;
359b1a0a8a3SJed Brown 
360b1a0a8a3SJed Brown   PetscFunctionBegin;
361f746d493SDmitry Karpeev   for (i=0; i<osm->n; i++) {
362b1a0a8a3SJed Brown     ierr = KSPSetUp(osm->ksp[i]);CHKERRQ(ierr);
363b1a0a8a3SJed Brown   }
364b1a0a8a3SJed Brown   PetscFunctionReturn(0);
365b1a0a8a3SJed Brown }
366b1a0a8a3SJed Brown 
367b1a0a8a3SJed Brown #undef __FUNCT__
368f746d493SDmitry Karpeev #define __FUNCT__ "PCApply_GASM"
369f746d493SDmitry Karpeev static PetscErrorCode PCApply_GASM(PC pc,Vec x,Vec y)
370b1a0a8a3SJed Brown {
371f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
372b1a0a8a3SJed Brown   PetscErrorCode ierr;
373f746d493SDmitry Karpeev   PetscInt       i;
374b1a0a8a3SJed Brown   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;
375b1a0a8a3SJed Brown 
376b1a0a8a3SJed Brown   PetscFunctionBegin;
377b1a0a8a3SJed Brown   /*
378b1a0a8a3SJed Brown      Support for limiting the restriction or interpolation to only local
379b1a0a8a3SJed Brown      subdomain values (leaving the other values 0).
380b1a0a8a3SJed Brown   */
381f746d493SDmitry Karpeev   if (!(osm->type & PC_GASM_RESTRICT)) {
382b1a0a8a3SJed Brown     forward = SCATTER_FORWARD_LOCAL;
383b1a0a8a3SJed Brown     /* have to zero the work RHS since scatter may leave some slots empty */
384f746d493SDmitry Karpeev     ierr = VecZeroEntries(osm->gx); CHKERRQ(ierr);
385b1a0a8a3SJed Brown   }
386f746d493SDmitry Karpeev   if (!(osm->type & PC_GASM_INTERPOLATE)) {
387b1a0a8a3SJed Brown     reverse = SCATTER_REVERSE_LOCAL;
388b1a0a8a3SJed Brown   }
389b1a0a8a3SJed Brown 
390f746d493SDmitry Karpeev   ierr = VecScatterBegin(osm->grestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
391b1a0a8a3SJed Brown   ierr = VecZeroEntries(y);CHKERRQ(ierr);
392f746d493SDmitry Karpeev   ierr = VecScatterEnd(osm->grestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
393b1a0a8a3SJed Brown   /* do the local solves */
394f746d493SDmitry Karpeev   for (i=0; i<osm->n; ++i) { /* Note that the solves are local, so we can go to osm->n, rather than osm->nmax. */
395b1a0a8a3SJed Brown     ierr = KSPSolve(osm->ksp[i],osm->x[i],osm->y[i]);CHKERRQ(ierr);
396b1a0a8a3SJed Brown   }
397f746d493SDmitry Karpeev   ierr = VecScatterBegin(osm->gprolongation,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
398f746d493SDmitry Karpeev   ierr = VecScatterEnd(osm->gprolongation,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
399b1a0a8a3SJed Brown   PetscFunctionReturn(0);
400b1a0a8a3SJed Brown }
401b1a0a8a3SJed Brown 
402b1a0a8a3SJed Brown #undef __FUNCT__
403f746d493SDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_GASM"
404f746d493SDmitry Karpeev static PetscErrorCode PCApplyTranspose_GASM(PC pc,Vec x,Vec y)
405b1a0a8a3SJed Brown {
406f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
407b1a0a8a3SJed Brown   PetscErrorCode ierr;
408f746d493SDmitry Karpeev   PetscInt       i;
409b1a0a8a3SJed Brown   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;
410b1a0a8a3SJed Brown 
411b1a0a8a3SJed Brown   PetscFunctionBegin;
412b1a0a8a3SJed Brown   /*
413b1a0a8a3SJed Brown      Support for limiting the restriction or interpolation to only local
414b1a0a8a3SJed Brown      subdomain values (leaving the other values 0).
415b1a0a8a3SJed Brown 
416f746d493SDmitry Karpeev      Note: these are reversed from the PCApply_GASM() because we are applying the
417b1a0a8a3SJed Brown      transpose of the three terms
418b1a0a8a3SJed Brown   */
419f746d493SDmitry Karpeev   if (!(osm->type & PC_GASM_INTERPOLATE)) {
420b1a0a8a3SJed Brown     forward = SCATTER_FORWARD_LOCAL;
421b1a0a8a3SJed Brown     /* have to zero the work RHS since scatter may leave some slots empty */
422f746d493SDmitry Karpeev     ierr = VecZeroEntries(osm->gx);CHKERRQ(ierr);
423b1a0a8a3SJed Brown   }
424f746d493SDmitry Karpeev   if (!(osm->type & PC_GASM_RESTRICT)) {
425b1a0a8a3SJed Brown     reverse = SCATTER_REVERSE_LOCAL;
426b1a0a8a3SJed Brown   }
427b1a0a8a3SJed Brown 
428ab3e923fSDmitry Karpeev   ierr = VecScatterBegin(osm->grestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
429b1a0a8a3SJed Brown   ierr = VecZeroEntries(y);CHKERRQ(ierr);
430ab3e923fSDmitry Karpeev   ierr = VecScatterEnd(osm->grestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
431b1a0a8a3SJed Brown   /* do the local solves */
432ab3e923fSDmitry Karpeev   for (i=0; i<osm->n; ++i) { /* Note that the solves are local, so we can go to osm->n, rather than osm->nmax. */
433b1a0a8a3SJed Brown     ierr = KSPSolveTranspose(osm->ksp[i],osm->x[i],osm->y[i]);CHKERRQ(ierr);
434b1a0a8a3SJed Brown   }
435ab3e923fSDmitry Karpeev   ierr = VecScatterBegin(osm->gprolongation,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
436ab3e923fSDmitry Karpeev   ierr = VecScatterEnd(osm->gprolongation,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
437b1a0a8a3SJed Brown   PetscFunctionReturn(0);
438b1a0a8a3SJed Brown }
439b1a0a8a3SJed Brown 
440b1a0a8a3SJed Brown #undef __FUNCT__
441f746d493SDmitry Karpeev #define __FUNCT__ "PCDestroy_GASM"
442f746d493SDmitry Karpeev static PetscErrorCode PCDestroy_GASM(PC pc)
443b1a0a8a3SJed Brown {
444f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
445b1a0a8a3SJed Brown   PetscErrorCode ierr;
446b1a0a8a3SJed Brown   PetscInt       i;
447b1a0a8a3SJed Brown 
448b1a0a8a3SJed Brown   PetscFunctionBegin;
449b1a0a8a3SJed Brown   if (osm->ksp) {
450f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
451b1a0a8a3SJed Brown       ierr = KSPDestroy(osm->ksp[i]);CHKERRQ(ierr);
452b1a0a8a3SJed Brown     }
453b1a0a8a3SJed Brown     ierr = PetscFree(osm->ksp);CHKERRQ(ierr);
454b1a0a8a3SJed Brown   }
455b1a0a8a3SJed Brown   if (osm->pmat) {
456f746d493SDmitry Karpeev     if (osm->n > 0) {
457ab3e923fSDmitry Karpeev       ierr = MatDestroyMatrices(osm->n,&osm->pmat);CHKERRQ(ierr);
458b1a0a8a3SJed Brown     }
459b1a0a8a3SJed Brown   }
460f746d493SDmitry Karpeev   for (i=0; i<osm->n; i++) {
461b1a0a8a3SJed Brown     ierr = VecDestroy(osm->x[i]);CHKERRQ(ierr);
462b1a0a8a3SJed Brown     ierr = VecDestroy(osm->y[i]);CHKERRQ(ierr);
463b1a0a8a3SJed Brown   }
464b1a0a8a3SJed Brown   ierr = PetscFree(osm->x);CHKERRQ(ierr);
465b1a0a8a3SJed Brown   ierr = PetscFree(osm->y);CHKERRQ(ierr);
466f746d493SDmitry Karpeev   ierr = VecDestroy(osm->gx); CHKERRQ(ierr);
467ab3e923fSDmitry Karpeev   ierr = VecDestroy(osm->gy); CHKERRQ(ierr);
468ab3e923fSDmitry Karpeev 
469ab3e923fSDmitry Karpeev   ierr = VecScatterDestroy(osm->grestriction); CHKERRQ(ierr);
470ab3e923fSDmitry Karpeev   ierr = VecScatterDestroy(osm->gprolongation); CHKERRQ(ierr);
471f746d493SDmitry Karpeev 
472b1a0a8a3SJed Brown   if (osm->is) {
473ab3e923fSDmitry Karpeev     ierr = PCGASMDestroySubdomains(osm->n,osm->is,osm->is_local);CHKERRQ(ierr);
474b1a0a8a3SJed Brown   }
475b1a0a8a3SJed Brown   ierr = PetscFree(osm);CHKERRQ(ierr);
476b1a0a8a3SJed Brown   PetscFunctionReturn(0);
477b1a0a8a3SJed Brown }
478b1a0a8a3SJed Brown 
479b1a0a8a3SJed Brown #undef __FUNCT__
480f746d493SDmitry Karpeev #define __FUNCT__ "PCSetFromOptions_GASM"
481ab3e923fSDmitry Karpeev static PetscErrorCode PCSetFromOptions_GASM(PC pc) {
482f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
483b1a0a8a3SJed Brown   PetscErrorCode ierr;
484b1a0a8a3SJed Brown   PetscInt       blocks,ovl;
485b1a0a8a3SJed Brown   PetscBool      symset,flg;
486f746d493SDmitry Karpeev   PCGASMType      gasmtype;
487b1a0a8a3SJed Brown 
488b1a0a8a3SJed Brown   PetscFunctionBegin;
489b1a0a8a3SJed Brown   /* set the type to symmetric if matrix is symmetric */
490b1a0a8a3SJed Brown   if (!osm->type_set && pc->pmat) {
491b1a0a8a3SJed Brown     ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr);
492f746d493SDmitry Karpeev     if (symset && flg) { osm->type = PC_GASM_BASIC; }
493b1a0a8a3SJed Brown   }
49444fe18e5SDmitry Karpeev   ierr = PetscOptionsHead("Generalized additive Schwarz options");CHKERRQ(ierr);
495f746d493SDmitry Karpeev     ierr = PetscOptionsInt("-pc_gasm_blocks","Number of subdomains","PCGASMSetTotalSubdomains",osm->n,&blocks,&flg);CHKERRQ(ierr);
496ab3e923fSDmitry Karpeev     if (flg) {ierr = PCGASMSetTotalSubdomains(pc,blocks);CHKERRQ(ierr); }
497f746d493SDmitry Karpeev     ierr = PetscOptionsInt("-pc_gasm_overlap","Number of grid points overlap","PCGASMSetOverlap",osm->overlap,&ovl,&flg);CHKERRQ(ierr);
498f746d493SDmitry Karpeev     if (flg) {ierr = PCGASMSetOverlap(pc,ovl);CHKERRQ(ierr); }
499b1a0a8a3SJed Brown     flg  = PETSC_FALSE;
500f746d493SDmitry Karpeev     ierr = PetscOptionsEnum("-pc_gasm_type","Type of restriction/extension","PCGASMSetType",PCGASMTypes,(PetscEnum)osm->type,(PetscEnum*)&gasmtype,&flg);CHKERRQ(ierr);
501f746d493SDmitry Karpeev     if (flg) {ierr = PCGASMSetType(pc,gasmtype);CHKERRQ(ierr); }
502b1a0a8a3SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
503b1a0a8a3SJed Brown   PetscFunctionReturn(0);
504b1a0a8a3SJed Brown }
505b1a0a8a3SJed Brown 
506b1a0a8a3SJed Brown /*------------------------------------------------------------------------------------*/
507b1a0a8a3SJed Brown 
508b1a0a8a3SJed Brown EXTERN_C_BEGIN
509b1a0a8a3SJed Brown #undef __FUNCT__
510f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetLocalSubdomains_GASM"
511f746d493SDmitry Karpeev PetscErrorCode PETSCKSP_DLLEXPORT PCGASMSetLocalSubdomains_GASM(PC pc,PetscInt n,IS is[],IS is_local[])
512b1a0a8a3SJed Brown {
513f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
514b1a0a8a3SJed Brown   PetscErrorCode ierr;
515b1a0a8a3SJed Brown   PetscInt       i;
516b1a0a8a3SJed Brown 
517b1a0a8a3SJed Brown   PetscFunctionBegin;
518b1a0a8a3SJed Brown   if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Each process must have 1 or more blocks, n = %D",n);
519f746d493SDmitry Karpeev   if (pc->setupcalled && (n != osm->n || is)) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetLocalSubdomains() should be called before calling PCSetUp().");
520b1a0a8a3SJed Brown 
521b1a0a8a3SJed Brown   if (!pc->setupcalled) {
52293c1ec32SDmitry Karpeev     osm->n            = n;
52393c1ec32SDmitry Karpeev     osm->is           = 0;
52493c1ec32SDmitry Karpeev     osm->is_local     = 0;
525b1a0a8a3SJed Brown     if (is) {
526b1a0a8a3SJed Brown       for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)is[i]);CHKERRQ(ierr);}
527b1a0a8a3SJed Brown     }
528b1a0a8a3SJed Brown     if (is_local) {
529b1a0a8a3SJed Brown       for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)is_local[i]);CHKERRQ(ierr);}
530b1a0a8a3SJed Brown     }
531b1a0a8a3SJed Brown     if (osm->is) {
532ab3e923fSDmitry Karpeev       ierr = PCGASMDestroySubdomains(osm->n,osm->is,osm->is_local);CHKERRQ(ierr);
533b1a0a8a3SJed Brown     }
534b1a0a8a3SJed Brown     if (is) {
535b1a0a8a3SJed Brown       ierr = PetscMalloc(n*sizeof(IS),&osm->is);CHKERRQ(ierr);
536b1a0a8a3SJed Brown       for (i=0; i<n; i++) { osm->is[i] = is[i]; }
537f746d493SDmitry Karpeev       /* Flag indicating that the user has set overlapping subdomains so PCGASM should not increase their size. */
538b1a0a8a3SJed Brown       osm->overlap = -1;
539b1a0a8a3SJed Brown     }
540b1a0a8a3SJed Brown     if (is_local) {
541b1a0a8a3SJed Brown       ierr = PetscMalloc(n*sizeof(IS),&osm->is_local);CHKERRQ(ierr);
542b1a0a8a3SJed Brown       for (i=0; i<n; i++) { osm->is_local[i] = is_local[i]; }
543b1a0a8a3SJed Brown     }
544b1a0a8a3SJed Brown   }
545b1a0a8a3SJed Brown   PetscFunctionReturn(0);
546b1a0a8a3SJed Brown }
547b1a0a8a3SJed Brown EXTERN_C_END
548b1a0a8a3SJed Brown 
549b1a0a8a3SJed Brown EXTERN_C_BEGIN
550b1a0a8a3SJed Brown #undef __FUNCT__
551f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetTotalSubdomains_GASM"
552ab3e923fSDmitry Karpeev PetscErrorCode PETSCKSP_DLLEXPORT PCGASMSetTotalSubdomains_GASM(PC pc,PetscInt N) {
553f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
554b1a0a8a3SJed Brown   PetscErrorCode ierr;
555b1a0a8a3SJed Brown   PetscMPIInt    rank,size;
556b1a0a8a3SJed Brown   PetscInt       n;
557b1a0a8a3SJed Brown 
558b1a0a8a3SJed Brown   PetscFunctionBegin;
559b1a0a8a3SJed Brown   if (N < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of total blocks must be > 0, N = %D",N);
560b1a0a8a3SJed Brown 
561b1a0a8a3SJed Brown   /*
562b1a0a8a3SJed Brown      Split the subdomains equally among all processors
563b1a0a8a3SJed Brown   */
564b1a0a8a3SJed Brown   ierr = MPI_Comm_rank(((PetscObject)pc)->comm,&rank);CHKERRQ(ierr);
565b1a0a8a3SJed Brown   ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr);
566b1a0a8a3SJed Brown   n = N/size + ((N % size) > rank);
567b1a0a8a3SJed Brown   if (!n) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Process %d must have at least one block: total processors %d total blocks %D",(int)rank,(int)size,N);
568f746d493SDmitry Karpeev   if (pc->setupcalled && n != osm->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetTotalSubdomains() should be called before PCSetUp().");
569b1a0a8a3SJed Brown   if (!pc->setupcalled) {
570b1a0a8a3SJed Brown     if (osm->is) {
571ab3e923fSDmitry Karpeev       ierr = PCGASMDestroySubdomains(osm->n,osm->is,osm->is_local);CHKERRQ(ierr);
572b1a0a8a3SJed Brown     }
57393c1ec32SDmitry Karpeev     osm->N            = N;
574f746d493SDmitry Karpeev     osm->n            = n;
575b1a0a8a3SJed Brown     osm->is           = 0;
576b1a0a8a3SJed Brown     osm->is_local     = 0;
577b1a0a8a3SJed Brown   }
578b1a0a8a3SJed Brown   PetscFunctionReturn(0);
579ab3e923fSDmitry Karpeev }/* PCGASMSetTotalSubdomains_GASM() */
580b1a0a8a3SJed Brown EXTERN_C_END
581b1a0a8a3SJed Brown 
582b1a0a8a3SJed Brown EXTERN_C_BEGIN
583b1a0a8a3SJed Brown #undef __FUNCT__
584f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetOverlap_GASM"
585f746d493SDmitry Karpeev PetscErrorCode PETSCKSP_DLLEXPORT PCGASMSetOverlap_GASM(PC pc,PetscInt ovl)
586b1a0a8a3SJed Brown {
587f746d493SDmitry Karpeev   PC_GASM *osm = (PC_GASM*)pc->data;
588b1a0a8a3SJed Brown 
589b1a0a8a3SJed Brown   PetscFunctionBegin;
590b1a0a8a3SJed Brown   if (ovl < 0) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested");
591f746d493SDmitry Karpeev   if (pc->setupcalled && ovl != osm->overlap) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetOverlap() should be called before PCSetUp().");
592b1a0a8a3SJed Brown   if (!pc->setupcalled) {
593b1a0a8a3SJed Brown     osm->overlap = ovl;
594b1a0a8a3SJed Brown   }
595b1a0a8a3SJed Brown   PetscFunctionReturn(0);
596b1a0a8a3SJed Brown }
597b1a0a8a3SJed Brown EXTERN_C_END
598b1a0a8a3SJed Brown 
599b1a0a8a3SJed Brown EXTERN_C_BEGIN
600b1a0a8a3SJed Brown #undef __FUNCT__
601f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetType_GASM"
602f746d493SDmitry Karpeev PetscErrorCode PETSCKSP_DLLEXPORT PCGASMSetType_GASM(PC pc,PCGASMType type)
603b1a0a8a3SJed Brown {
604f746d493SDmitry Karpeev   PC_GASM *osm = (PC_GASM*)pc->data;
605b1a0a8a3SJed Brown 
606b1a0a8a3SJed Brown   PetscFunctionBegin;
607b1a0a8a3SJed Brown   osm->type     = type;
608b1a0a8a3SJed Brown   osm->type_set = PETSC_TRUE;
609b1a0a8a3SJed Brown   PetscFunctionReturn(0);
610b1a0a8a3SJed Brown }
611b1a0a8a3SJed Brown EXTERN_C_END
612b1a0a8a3SJed Brown 
613b1a0a8a3SJed Brown EXTERN_C_BEGIN
614b1a0a8a3SJed Brown #undef __FUNCT__
615f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetSortIndices_GASM"
616f746d493SDmitry Karpeev PetscErrorCode PETSCKSP_DLLEXPORT PCGASMSetSortIndices_GASM(PC pc,PetscBool  doSort)
617b1a0a8a3SJed Brown {
618f746d493SDmitry Karpeev   PC_GASM *osm = (PC_GASM*)pc->data;
619b1a0a8a3SJed Brown 
620b1a0a8a3SJed Brown   PetscFunctionBegin;
621b1a0a8a3SJed Brown   osm->sort_indices = doSort;
622b1a0a8a3SJed Brown   PetscFunctionReturn(0);
623b1a0a8a3SJed Brown }
624b1a0a8a3SJed Brown EXTERN_C_END
625b1a0a8a3SJed Brown 
626b1a0a8a3SJed Brown EXTERN_C_BEGIN
627b1a0a8a3SJed Brown #undef __FUNCT__
628f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMGetSubKSP_GASM"
62944fe18e5SDmitry Karpeev /*
63044fe18e5SDmitry Karpeev    FIX: This routine might need to be modified once multiple ranks per subdomain are allowed.
63144fe18e5SDmitry Karpeev         In particular, it would upset the global subdomain number calculation.
63244fe18e5SDmitry Karpeev */
633ab3e923fSDmitry Karpeev PetscErrorCode PETSCKSP_DLLEXPORT PCGASMGetSubKSP_GASM(PC pc,PetscInt *n,PetscInt *first,KSP **ksp)
634b1a0a8a3SJed Brown {
635f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
636b1a0a8a3SJed Brown   PetscErrorCode ierr;
637b1a0a8a3SJed Brown 
638b1a0a8a3SJed Brown   PetscFunctionBegin;
639ab3e923fSDmitry Karpeev   if (osm->n < 1) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ORDER,"Need to call PCSetUP() on PC (or KSPSetUp() on the outer KSP object) before calling here");
640b1a0a8a3SJed Brown 
641ab3e923fSDmitry Karpeev   if (n) {
642ab3e923fSDmitry Karpeev     *n = osm->n;
643b1a0a8a3SJed Brown   }
644ab3e923fSDmitry Karpeev   if (first) {
645ab3e923fSDmitry Karpeev     ierr = MPI_Scan(&osm->n,first,1,MPIU_INT,MPI_SUM,((PetscObject)pc)->comm);CHKERRQ(ierr);
646ab3e923fSDmitry Karpeev     *first -= osm->n;
647b1a0a8a3SJed Brown   }
648b1a0a8a3SJed Brown   if (ksp) {
649b1a0a8a3SJed Brown     /* Assume that local solves are now different; not necessarily
650f746d493SDmitry Karpeev        true though!  This flag is used only for PCView_GASM() */
651b1a0a8a3SJed Brown     *ksp                   = osm->ksp;
652b1a0a8a3SJed Brown     osm->same_local_solves = PETSC_FALSE;
653b1a0a8a3SJed Brown   }
654b1a0a8a3SJed Brown   PetscFunctionReturn(0);
655ab3e923fSDmitry Karpeev }/* PCGASMGetSubKSP_GASM() */
656b1a0a8a3SJed Brown EXTERN_C_END
657b1a0a8a3SJed Brown 
658b1a0a8a3SJed Brown 
659b1a0a8a3SJed Brown #undef __FUNCT__
660f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetLocalSubdomains"
661b1a0a8a3SJed Brown /*@C
662f746d493SDmitry Karpeev     PCGASMSetLocalSubdomains - Sets the local subdomains (for this processor
663b1a0a8a3SJed Brown     only) for the additive Schwarz preconditioner.
664b1a0a8a3SJed Brown 
665b1a0a8a3SJed Brown     Collective on PC
666b1a0a8a3SJed Brown 
667b1a0a8a3SJed Brown     Input Parameters:
668b1a0a8a3SJed Brown +   pc - the preconditioner context
669b1a0a8a3SJed Brown .   n - the number of subdomains for this processor (default value = 1)
670b1a0a8a3SJed Brown .   is - the index set that defines the subdomains for this processor
671b1a0a8a3SJed Brown          (or PETSC_NULL for PETSc to determine subdomains)
672b1a0a8a3SJed Brown -   is_local - the index sets that define the local part of the subdomains for this processor
673b1a0a8a3SJed Brown          (or PETSC_NULL to use the default of 1 subdomain per process)
674b1a0a8a3SJed Brown 
675b1a0a8a3SJed Brown     Notes:
676b1a0a8a3SJed Brown     The IS numbering is in the parallel, global numbering of the vector.
677b1a0a8a3SJed Brown 
678f746d493SDmitry Karpeev     By default the GASM preconditioner uses 1 block per processor.
679b1a0a8a3SJed Brown 
680f746d493SDmitry Karpeev     Use PCGASMSetTotalSubdomains() to set the subdomains for all processors.
681b1a0a8a3SJed Brown 
682b1a0a8a3SJed Brown     Level: advanced
683b1a0a8a3SJed Brown 
684f746d493SDmitry Karpeev .keywords: PC, GASM, set, local, subdomains, additive Schwarz
685b1a0a8a3SJed Brown 
686f746d493SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(),
687f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D(), PCGASMGetLocalSubdomains()
688b1a0a8a3SJed Brown @*/
689f746d493SDmitry Karpeev PetscErrorCode PETSCKSP_DLLEXPORT PCGASMSetLocalSubdomains(PC pc,PetscInt n,IS is[],IS is_local[])
690b1a0a8a3SJed Brown {
691b1a0a8a3SJed Brown   PetscErrorCode ierr;
692b1a0a8a3SJed Brown 
693b1a0a8a3SJed Brown   PetscFunctionBegin;
694b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
695f746d493SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetLocalSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,n,is,is_local));CHKERRQ(ierr);
696b1a0a8a3SJed Brown   PetscFunctionReturn(0);
697b1a0a8a3SJed Brown }
698b1a0a8a3SJed Brown 
699b1a0a8a3SJed Brown #undef __FUNCT__
700f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetTotalSubdomains"
701b1a0a8a3SJed Brown /*@C
702f746d493SDmitry Karpeev     PCGASMSetTotalSubdomains - Sets the subdomains for all processor for the
703b1a0a8a3SJed Brown     additive Schwarz preconditioner.  Either all or no processors in the
704b1a0a8a3SJed Brown     PC communicator must call this routine, with the same index sets.
705b1a0a8a3SJed Brown 
706b1a0a8a3SJed Brown     Collective on PC
707b1a0a8a3SJed Brown 
708b1a0a8a3SJed Brown     Input Parameters:
709b1a0a8a3SJed Brown +   pc - the preconditioner context
710b1a0a8a3SJed Brown .   n - the number of subdomains for all processors
711b1a0a8a3SJed Brown .   is - the index sets that define the subdomains for all processor
712b1a0a8a3SJed Brown          (or PETSC_NULL for PETSc to determine subdomains)
713b1a0a8a3SJed Brown -   is_local - the index sets that define the local part of the subdomains for this processor
714b1a0a8a3SJed Brown          (or PETSC_NULL to use the default of 1 subdomain per process)
715b1a0a8a3SJed Brown 
716b1a0a8a3SJed Brown     Options Database Key:
717b1a0a8a3SJed Brown     To set the total number of subdomain blocks rather than specify the
718b1a0a8a3SJed Brown     index sets, use the option
719f746d493SDmitry Karpeev .    -pc_gasm_blocks <blks> - Sets total blocks
720b1a0a8a3SJed Brown 
721b1a0a8a3SJed Brown     Notes:
722b1a0a8a3SJed Brown     Currently you cannot use this to set the actual subdomains with the argument is.
723b1a0a8a3SJed Brown 
724f746d493SDmitry Karpeev     By default the GASM preconditioner uses 1 block per processor.
725b1a0a8a3SJed Brown 
726b1a0a8a3SJed Brown     These index sets cannot be destroyed until after completion of the
727f746d493SDmitry Karpeev     linear solves for which the GASM preconditioner is being used.
728b1a0a8a3SJed Brown 
729f746d493SDmitry Karpeev     Use PCGASMSetLocalSubdomains() to set local subdomains.
730b1a0a8a3SJed Brown 
731b1a0a8a3SJed Brown     Level: advanced
732b1a0a8a3SJed Brown 
733f746d493SDmitry Karpeev .keywords: PC, GASM, set, total, global, subdomains, additive Schwarz
734b1a0a8a3SJed Brown 
735f746d493SDmitry Karpeev .seealso: PCGASMSetLocalSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(),
736f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D()
737b1a0a8a3SJed Brown @*/
738f746d493SDmitry Karpeev PetscErrorCode PETSCKSP_DLLEXPORT PCGASMSetTotalSubdomains(PC pc,PetscInt N)
739b1a0a8a3SJed Brown {
740b1a0a8a3SJed Brown   PetscErrorCode ierr;
741b1a0a8a3SJed Brown 
742b1a0a8a3SJed Brown   PetscFunctionBegin;
743b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
744f746d493SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetTotalSubdomains_C",(PC,PetscInt),(pc,N));CHKERRQ(ierr);
745b1a0a8a3SJed Brown   PetscFunctionReturn(0);
746b1a0a8a3SJed Brown }
747b1a0a8a3SJed Brown 
748b1a0a8a3SJed Brown #undef __FUNCT__
749f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetOverlap"
750b1a0a8a3SJed Brown /*@
751f746d493SDmitry Karpeev     PCGASMSetOverlap - Sets the overlap between a pair of subdomains for the
752b1a0a8a3SJed Brown     additive Schwarz preconditioner.  Either all or no processors in the
753b1a0a8a3SJed Brown     PC communicator must call this routine.
754b1a0a8a3SJed Brown 
755b1a0a8a3SJed Brown     Logically Collective on PC
756b1a0a8a3SJed Brown 
757b1a0a8a3SJed Brown     Input Parameters:
758b1a0a8a3SJed Brown +   pc  - the preconditioner context
759b1a0a8a3SJed Brown -   ovl - the amount of overlap between subdomains (ovl >= 0, default value = 1)
760b1a0a8a3SJed Brown 
761b1a0a8a3SJed Brown     Options Database Key:
762f746d493SDmitry Karpeev .   -pc_gasm_overlap <ovl> - Sets overlap
763b1a0a8a3SJed Brown 
764b1a0a8a3SJed Brown     Notes:
765f746d493SDmitry Karpeev     By default the GASM preconditioner uses 1 block per processor.  To use
766f746d493SDmitry Karpeev     multiple blocks per perocessor, see PCGASMSetTotalSubdomains() and
767f746d493SDmitry Karpeev     PCGASMSetLocalSubdomains() (and the option -pc_gasm_blocks <blks>).
768b1a0a8a3SJed Brown 
769b1a0a8a3SJed Brown     The overlap defaults to 1, so if one desires that no additional
770b1a0a8a3SJed Brown     overlap be computed beyond what may have been set with a call to
771f746d493SDmitry Karpeev     PCGASMSetTotalSubdomains() or PCGASMSetLocalSubdomains(), then ovl
772b1a0a8a3SJed Brown     must be set to be 0.  In particular, if one does not explicitly set
773b1a0a8a3SJed Brown     the subdomains an application code, then all overlap would be computed
774f746d493SDmitry Karpeev     internally by PETSc, and using an overlap of 0 would result in an GASM
775b1a0a8a3SJed Brown     variant that is equivalent to the block Jacobi preconditioner.
776b1a0a8a3SJed Brown 
777b1a0a8a3SJed Brown     Note that one can define initial index sets with any overlap via
778f746d493SDmitry Karpeev     PCGASMSetTotalSubdomains() or PCGASMSetLocalSubdomains(); the routine
779f746d493SDmitry Karpeev     PCGASMSetOverlap() merely allows PETSc to extend that overlap further
780b1a0a8a3SJed Brown     if desired.
781b1a0a8a3SJed Brown 
782b1a0a8a3SJed Brown     Level: intermediate
783b1a0a8a3SJed Brown 
784f746d493SDmitry Karpeev .keywords: PC, GASM, set, overlap
785b1a0a8a3SJed Brown 
786f746d493SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetLocalSubdomains(), PCGASMGetSubKSP(),
787f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D(), PCGASMGetLocalSubdomains()
788b1a0a8a3SJed Brown @*/
789f746d493SDmitry Karpeev PetscErrorCode PETSCKSP_DLLEXPORT PCGASMSetOverlap(PC pc,PetscInt ovl)
790b1a0a8a3SJed Brown {
791b1a0a8a3SJed Brown   PetscErrorCode ierr;
792b1a0a8a3SJed Brown 
793b1a0a8a3SJed Brown   PetscFunctionBegin;
794b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
795b1a0a8a3SJed Brown   PetscValidLogicalCollectiveInt(pc,ovl,2);
796f746d493SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetOverlap_C",(PC,PetscInt),(pc,ovl));CHKERRQ(ierr);
797b1a0a8a3SJed Brown   PetscFunctionReturn(0);
798b1a0a8a3SJed Brown }
799b1a0a8a3SJed Brown 
800b1a0a8a3SJed Brown #undef __FUNCT__
801f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetType"
802b1a0a8a3SJed Brown /*@
803f746d493SDmitry Karpeev     PCGASMSetType - Sets the type of restriction and interpolation used
804b1a0a8a3SJed Brown     for local problems in the additive Schwarz method.
805b1a0a8a3SJed Brown 
806b1a0a8a3SJed Brown     Logically Collective on PC
807b1a0a8a3SJed Brown 
808b1a0a8a3SJed Brown     Input Parameters:
809b1a0a8a3SJed Brown +   pc  - the preconditioner context
810f746d493SDmitry Karpeev -   type - variant of GASM, one of
811b1a0a8a3SJed Brown .vb
812f746d493SDmitry Karpeev       PC_GASM_BASIC       - full interpolation and restriction
813f746d493SDmitry Karpeev       PC_GASM_RESTRICT    - full restriction, local processor interpolation
814f746d493SDmitry Karpeev       PC_GASM_INTERPOLATE - full interpolation, local processor restriction
815f746d493SDmitry Karpeev       PC_GASM_NONE        - local processor restriction and interpolation
816b1a0a8a3SJed Brown .ve
817b1a0a8a3SJed Brown 
818b1a0a8a3SJed Brown     Options Database Key:
819f746d493SDmitry Karpeev .   -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type
820b1a0a8a3SJed Brown 
821b1a0a8a3SJed Brown     Level: intermediate
822b1a0a8a3SJed Brown 
823f746d493SDmitry Karpeev .keywords: PC, GASM, set, type
824b1a0a8a3SJed Brown 
825f746d493SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetTotalSubdomains(), PCGASMGetSubKSP(),
826f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D()
827b1a0a8a3SJed Brown @*/
828f746d493SDmitry Karpeev PetscErrorCode PETSCKSP_DLLEXPORT PCGASMSetType(PC pc,PCGASMType type)
829b1a0a8a3SJed Brown {
830b1a0a8a3SJed Brown   PetscErrorCode ierr;
831b1a0a8a3SJed Brown 
832b1a0a8a3SJed Brown   PetscFunctionBegin;
833b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
834b1a0a8a3SJed Brown   PetscValidLogicalCollectiveEnum(pc,type,2);
835f746d493SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetType_C",(PC,PCGASMType),(pc,type));CHKERRQ(ierr);
836b1a0a8a3SJed Brown   PetscFunctionReturn(0);
837b1a0a8a3SJed Brown }
838b1a0a8a3SJed Brown 
839b1a0a8a3SJed Brown #undef __FUNCT__
840f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetSortIndices"
841b1a0a8a3SJed Brown /*@
842f746d493SDmitry Karpeev     PCGASMSetSortIndices - Determines whether subdomain indices are sorted.
843b1a0a8a3SJed Brown 
844b1a0a8a3SJed Brown     Logically Collective on PC
845b1a0a8a3SJed Brown 
846b1a0a8a3SJed Brown     Input Parameters:
847b1a0a8a3SJed Brown +   pc  - the preconditioner context
848b1a0a8a3SJed Brown -   doSort - sort the subdomain indices
849b1a0a8a3SJed Brown 
850b1a0a8a3SJed Brown     Level: intermediate
851b1a0a8a3SJed Brown 
852f746d493SDmitry Karpeev .keywords: PC, GASM, set, type
853b1a0a8a3SJed Brown 
854f746d493SDmitry Karpeev .seealso: PCGASMSetLocalSubdomains(), PCGASMSetTotalSubdomains(), PCGASMGetSubKSP(),
855f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D()
856b1a0a8a3SJed Brown @*/
857f746d493SDmitry Karpeev PetscErrorCode PETSCKSP_DLLEXPORT PCGASMSetSortIndices(PC pc,PetscBool  doSort)
858b1a0a8a3SJed Brown {
859b1a0a8a3SJed Brown   PetscErrorCode ierr;
860b1a0a8a3SJed Brown 
861b1a0a8a3SJed Brown   PetscFunctionBegin;
862b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
863b1a0a8a3SJed Brown   PetscValidLogicalCollectiveBool(pc,doSort,2);
864f746d493SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetSortIndices_C",(PC,PetscBool),(pc,doSort));CHKERRQ(ierr);
865b1a0a8a3SJed Brown   PetscFunctionReturn(0);
866b1a0a8a3SJed Brown }
867b1a0a8a3SJed Brown 
868b1a0a8a3SJed Brown #undef __FUNCT__
869f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMGetSubKSP"
870b1a0a8a3SJed Brown /*@C
871f746d493SDmitry Karpeev    PCGASMGetSubKSP - Gets the local KSP contexts for all blocks on
872b1a0a8a3SJed Brown    this processor.
873b1a0a8a3SJed Brown 
874b1a0a8a3SJed Brown    Collective on PC iff first_local is requested
875b1a0a8a3SJed Brown 
876b1a0a8a3SJed Brown    Input Parameter:
877b1a0a8a3SJed Brown .  pc - the preconditioner context
878b1a0a8a3SJed Brown 
879b1a0a8a3SJed Brown    Output Parameters:
880b1a0a8a3SJed Brown +  n_local - the number of blocks on this processor or PETSC_NULL
881b1a0a8a3SJed Brown .  first_local - the global number of the first block on this processor or PETSC_NULL,
882b1a0a8a3SJed Brown                  all processors must request or all must pass PETSC_NULL
883b1a0a8a3SJed Brown -  ksp - the array of KSP contexts
884b1a0a8a3SJed Brown 
885b1a0a8a3SJed Brown    Note:
886f746d493SDmitry Karpeev    After PCGASMGetSubKSP() the array of KSPes is not to be freed
887b1a0a8a3SJed Brown 
888b1a0a8a3SJed Brown    Currently for some matrix implementations only 1 block per processor
889b1a0a8a3SJed Brown    is supported.
890b1a0a8a3SJed Brown 
891f746d493SDmitry Karpeev    You must call KSPSetUp() before calling PCGASMGetSubKSP().
892b1a0a8a3SJed Brown 
893b1a0a8a3SJed Brown    Level: advanced
894b1a0a8a3SJed Brown 
895f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, get, sub, KSP, context
896b1a0a8a3SJed Brown 
897f746d493SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetTotalSubdomains(), PCGASMSetOverlap(),
898f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D(),
899b1a0a8a3SJed Brown @*/
900f746d493SDmitry Karpeev PetscErrorCode PETSCKSP_DLLEXPORT PCGASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[])
901b1a0a8a3SJed Brown {
902b1a0a8a3SJed Brown   PetscErrorCode ierr;
903b1a0a8a3SJed Brown 
904b1a0a8a3SJed Brown   PetscFunctionBegin;
905b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
906f746d493SDmitry Karpeev   ierr = PetscUseMethod(pc,"PCGASMGetSubKSP_C",(PC,PetscInt*,PetscInt*,KSP **),(pc,n_local,first_local,ksp));CHKERRQ(ierr);
907b1a0a8a3SJed Brown   PetscFunctionReturn(0);
908b1a0a8a3SJed Brown }
909b1a0a8a3SJed Brown 
910b1a0a8a3SJed Brown /* -------------------------------------------------------------------------------------*/
911b1a0a8a3SJed Brown /*MC
912f746d493SDmitry Karpeev    PCGASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with
913b1a0a8a3SJed Brown            its own KSP object.
914b1a0a8a3SJed Brown 
915b1a0a8a3SJed Brown    Options Database Keys:
916f746d493SDmitry Karpeev +  -pc_gasm_truelocal - Activates PCGGASMSetUseTrueLocal()
917f746d493SDmitry Karpeev .  -pc_gasm_blocks <blks> - Sets total blocks
918f746d493SDmitry Karpeev .  -pc_gasm_overlap <ovl> - Sets overlap
919f746d493SDmitry Karpeev -  -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type
920b1a0a8a3SJed Brown 
921b1a0a8a3SJed Brown      IMPORTANT: If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you
922f746d493SDmitry Karpeev       will get a different convergence rate due to the default option of -pc_gasm_type restrict. Use
923f746d493SDmitry Karpeev       -pc_gasm_type basic to use the standard GASM.
924b1a0a8a3SJed Brown 
925b1a0a8a3SJed Brown    Notes: Each processor can have one or more blocks, but a block cannot be shared by more
926b1a0a8a3SJed Brown      than one processor. Defaults to one block per processor.
927b1a0a8a3SJed Brown 
928b1a0a8a3SJed Brown      To set options on the solvers for each block append -sub_ to all the KSP, and PC
929b1a0a8a3SJed Brown         options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly
930b1a0a8a3SJed Brown 
931f746d493SDmitry Karpeev      To set the options on the solvers separate for each block call PCGASMGetSubKSP()
932b1a0a8a3SJed Brown          and set the options directly on the resulting KSP object (you can access its PC
933b1a0a8a3SJed Brown          with KSPGetPC())
934b1a0a8a3SJed Brown 
935b1a0a8a3SJed Brown 
936b1a0a8a3SJed Brown    Level: beginner
937b1a0a8a3SJed Brown 
938b1a0a8a3SJed Brown    Concepts: additive Schwarz method
939b1a0a8a3SJed Brown 
940b1a0a8a3SJed Brown     References:
941b1a0a8a3SJed Brown     An additive variant of the Schwarz alternating method for the case of many subregions
942b1a0a8a3SJed Brown     M Dryja, OB Widlund - Courant Institute, New York University Technical report
943b1a0a8a3SJed Brown 
944b1a0a8a3SJed Brown     Domain Decompositions: Parallel Multilevel Methods for Elliptic Partial Differential Equations,
945b1a0a8a3SJed Brown     Barry Smith, Petter Bjorstad, and William Gropp, Cambridge University Press, ISBN 0-521-49589-X.
946b1a0a8a3SJed Brown 
947b1a0a8a3SJed Brown .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
948f746d493SDmitry Karpeev            PCBJACOBI, PCGASMSetUseTrueLocal(), PCGASMGetSubKSP(), PCGASMSetLocalSubdomains(),
949f746d493SDmitry Karpeev            PCGASMSetTotalSubdomains(), PCSetModifySubmatrices(), PCGASMSetOverlap(), PCGASMSetType()
950b1a0a8a3SJed Brown 
951b1a0a8a3SJed Brown M*/
952b1a0a8a3SJed Brown 
953b1a0a8a3SJed Brown EXTERN_C_BEGIN
954b1a0a8a3SJed Brown #undef __FUNCT__
955f746d493SDmitry Karpeev #define __FUNCT__ "PCCreate_GASM"
956f746d493SDmitry Karpeev PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_GASM(PC pc)
957b1a0a8a3SJed Brown {
958b1a0a8a3SJed Brown   PetscErrorCode ierr;
959f746d493SDmitry Karpeev   PC_GASM         *osm;
960b1a0a8a3SJed Brown 
961b1a0a8a3SJed Brown   PetscFunctionBegin;
962f746d493SDmitry Karpeev   ierr = PetscNewLog(pc,PC_GASM,&osm);CHKERRQ(ierr);
963ab3e923fSDmitry Karpeev   osm->N                 = PETSC_DECIDE;
964ab3e923fSDmitry Karpeev   osm->n                 = 0;
965ab3e923fSDmitry Karpeev   osm->nmax              = 0;
966b1a0a8a3SJed Brown   osm->overlap           = 1;
967b1a0a8a3SJed Brown   osm->ksp               = 0;
968ab3e923fSDmitry Karpeev   osm->grestriction      = 0;
969ab3e923fSDmitry Karpeev   osm->gprolongation     = 0;
970ab3e923fSDmitry Karpeev   osm->gx                = 0;
971ab3e923fSDmitry Karpeev   osm->gy                = 0;
972b1a0a8a3SJed Brown   osm->x                 = 0;
973b1a0a8a3SJed Brown   osm->y                 = 0;
974b1a0a8a3SJed Brown   osm->is                = 0;
975b1a0a8a3SJed Brown   osm->is_local          = 0;
976b1a0a8a3SJed Brown   osm->pmat              = 0;
977f746d493SDmitry Karpeev   osm->type              = PC_GASM_RESTRICT;
978b1a0a8a3SJed Brown   osm->same_local_solves = PETSC_TRUE;
979b1a0a8a3SJed Brown   osm->sort_indices      = PETSC_TRUE;
980b1a0a8a3SJed Brown 
981b1a0a8a3SJed Brown   pc->data                   = (void*)osm;
982f746d493SDmitry Karpeev   pc->ops->apply             = PCApply_GASM;
983f746d493SDmitry Karpeev   pc->ops->applytranspose    = PCApplyTranspose_GASM;
984f746d493SDmitry Karpeev   pc->ops->setup             = PCSetUp_GASM;
985f746d493SDmitry Karpeev   pc->ops->destroy           = PCDestroy_GASM;
986f746d493SDmitry Karpeev   pc->ops->setfromoptions    = PCSetFromOptions_GASM;
987f746d493SDmitry Karpeev   pc->ops->setuponblocks     = PCSetUpOnBlocks_GASM;
988f746d493SDmitry Karpeev   pc->ops->view              = PCView_GASM;
989b1a0a8a3SJed Brown   pc->ops->applyrichardson   = 0;
990b1a0a8a3SJed Brown 
991f746d493SDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetLocalSubdomains_C","PCGASMSetLocalSubdomains_GASM",
992f746d493SDmitry Karpeev                     PCGASMSetLocalSubdomains_GASM);CHKERRQ(ierr);
993f746d493SDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetTotalSubdomains_C","PCGASMSetTotalSubdomains_GASM",
994f746d493SDmitry Karpeev                     PCGASMSetTotalSubdomains_GASM);CHKERRQ(ierr);
995f746d493SDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetOverlap_C","PCGASMSetOverlap_GASM",
996f746d493SDmitry Karpeev                     PCGASMSetOverlap_GASM);CHKERRQ(ierr);
997f746d493SDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetType_C","PCGASMSetType_GASM",
998f746d493SDmitry Karpeev                     PCGASMSetType_GASM);CHKERRQ(ierr);
999f746d493SDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetSortIndices_C","PCGASMSetSortIndices_GASM",
1000f746d493SDmitry Karpeev                     PCGASMSetSortIndices_GASM);CHKERRQ(ierr);
1001f746d493SDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMGetSubKSP_C","PCGASMGetSubKSP_GASM",
1002f746d493SDmitry Karpeev                     PCGASMGetSubKSP_GASM);CHKERRQ(ierr);
1003b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1004b1a0a8a3SJed Brown }
1005b1a0a8a3SJed Brown EXTERN_C_END
1006b1a0a8a3SJed Brown 
1007b1a0a8a3SJed Brown 
1008b1a0a8a3SJed Brown #undef __FUNCT__
1009f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMCreateSubdomains"
1010b1a0a8a3SJed Brown /*@C
1011f746d493SDmitry Karpeev    PCGASMCreateSubdomains - Creates the index sets for the overlapping Schwarz
1012b1a0a8a3SJed Brown    preconditioner for a any problem on a general grid.
1013b1a0a8a3SJed Brown 
1014b1a0a8a3SJed Brown    Collective
1015b1a0a8a3SJed Brown 
1016b1a0a8a3SJed Brown    Input Parameters:
1017b1a0a8a3SJed Brown +  A - The global matrix operator
1018b1a0a8a3SJed Brown -  n - the number of local blocks
1019b1a0a8a3SJed Brown 
1020b1a0a8a3SJed Brown    Output Parameters:
1021b1a0a8a3SJed Brown .  outis - the array of index sets defining the subdomains
1022b1a0a8a3SJed Brown 
1023b1a0a8a3SJed Brown    Level: advanced
1024b1a0a8a3SJed Brown 
1025f746d493SDmitry Karpeev    Note: this generates nonoverlapping subdomains; the PCGASM will generate the overlap
1026f746d493SDmitry Karpeev     from these if you use PCGASMSetLocalSubdomains()
1027b1a0a8a3SJed Brown 
1028b1a0a8a3SJed Brown     In the Fortran version you must provide the array outis[] already allocated of length n.
1029b1a0a8a3SJed Brown 
1030f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, unstructured grid
1031b1a0a8a3SJed Brown 
1032f746d493SDmitry Karpeev .seealso: PCGASMSetLocalSubdomains(), PCGASMDestroySubdomains()
1033b1a0a8a3SJed Brown @*/
1034f746d493SDmitry Karpeev PetscErrorCode PETSCKSP_DLLEXPORT PCGASMCreateSubdomains(Mat A, PetscInt n, IS* outis[])
1035b1a0a8a3SJed Brown {
1036b1a0a8a3SJed Brown   MatPartitioning           mpart;
1037b1a0a8a3SJed Brown   const char                *prefix;
1038b1a0a8a3SJed Brown   PetscErrorCode            (*f)(Mat,PetscBool *,MatReuse,Mat*);
1039b1a0a8a3SJed Brown   PetscMPIInt               size;
1040b1a0a8a3SJed Brown   PetscInt                  i,j,rstart,rend,bs;
1041b1a0a8a3SJed Brown   PetscBool                 iscopy = PETSC_FALSE,isbaij = PETSC_FALSE,foundpart = PETSC_FALSE;
1042b1a0a8a3SJed Brown   Mat                       Ad = PETSC_NULL, adj;
1043b1a0a8a3SJed Brown   IS                        ispart,isnumb,*is;
1044b1a0a8a3SJed Brown   PetscErrorCode            ierr;
1045b1a0a8a3SJed Brown 
1046b1a0a8a3SJed Brown   PetscFunctionBegin;
1047b1a0a8a3SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1048b1a0a8a3SJed Brown   PetscValidPointer(outis,3);
1049b1a0a8a3SJed Brown   if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of local blocks must be > 0, n = %D",n);
1050b1a0a8a3SJed Brown 
1051b1a0a8a3SJed Brown   /* Get prefix, row distribution, and block size */
1052b1a0a8a3SJed Brown   ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr);
1053b1a0a8a3SJed Brown   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
1054b1a0a8a3SJed Brown   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
1055b1a0a8a3SJed Brown   if (rstart/bs*bs != rstart || rend/bs*bs != rend) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"bad row distribution [%D,%D) for matrix block size %D",rstart,rend,bs);
1056b1a0a8a3SJed Brown 
1057b1a0a8a3SJed Brown   /* Get diagonal block from matrix if possible */
1058b1a0a8a3SJed Brown   ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr);
1059b1a0a8a3SJed Brown   ierr = PetscObjectQueryFunction((PetscObject)A,"MatGetDiagonalBlock_C",(void (**)(void))&f);CHKERRQ(ierr);
1060b1a0a8a3SJed Brown   if (f) {
1061b1a0a8a3SJed Brown     ierr = (*f)(A,&iscopy,MAT_INITIAL_MATRIX,&Ad);CHKERRQ(ierr);
1062b1a0a8a3SJed Brown   } else if (size == 1) {
1063b1a0a8a3SJed Brown     iscopy = PETSC_FALSE; Ad = A;
1064b1a0a8a3SJed Brown   } else {
1065b1a0a8a3SJed Brown     iscopy = PETSC_FALSE; Ad = PETSC_NULL;
1066b1a0a8a3SJed Brown   }
1067b1a0a8a3SJed Brown   if (Ad) {
1068b1a0a8a3SJed Brown     ierr = PetscTypeCompare((PetscObject)Ad,MATSEQBAIJ,&isbaij);CHKERRQ(ierr);
1069b1a0a8a3SJed Brown     if (!isbaij) {ierr = PetscTypeCompare((PetscObject)Ad,MATSEQSBAIJ,&isbaij);CHKERRQ(ierr);}
1070b1a0a8a3SJed Brown   }
1071b1a0a8a3SJed Brown   if (Ad && n > 1) {
1072b1a0a8a3SJed Brown     PetscBool  match,done;
1073b1a0a8a3SJed Brown     /* Try to setup a good matrix partitioning if available */
1074b1a0a8a3SJed Brown     ierr = MatPartitioningCreate(PETSC_COMM_SELF,&mpart);CHKERRQ(ierr);
1075b1a0a8a3SJed Brown     ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr);
1076b1a0a8a3SJed Brown     ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr);
1077b1a0a8a3SJed Brown     ierr = PetscTypeCompare((PetscObject)mpart,MATPARTITIONINGCURRENT,&match);CHKERRQ(ierr);
1078b1a0a8a3SJed Brown     if (!match) {
1079b1a0a8a3SJed Brown       ierr = PetscTypeCompare((PetscObject)mpart,MATPARTITIONINGSQUARE,&match);CHKERRQ(ierr);
1080b1a0a8a3SJed Brown     }
1081b1a0a8a3SJed Brown     if (!match) { /* assume a "good" partitioner is available */
1082b1a0a8a3SJed Brown       PetscInt na,*ia,*ja;
1083b1a0a8a3SJed Brown       ierr = MatGetRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr);
1084b1a0a8a3SJed Brown       if (done) {
1085b1a0a8a3SJed Brown 	/* Build adjacency matrix by hand. Unfortunately a call to
1086b1a0a8a3SJed Brown 	   MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will
1087b1a0a8a3SJed Brown 	   remove the block-aij structure and we cannot expect
1088b1a0a8a3SJed Brown 	   MatPartitioning to split vertices as we need */
1089b1a0a8a3SJed Brown 	PetscInt i,j,*row,len,nnz,cnt,*iia=0,*jja=0;
1090b1a0a8a3SJed Brown 	nnz = 0;
1091b1a0a8a3SJed Brown 	for (i=0; i<na; i++) { /* count number of nonzeros */
1092b1a0a8a3SJed Brown 	  len = ia[i+1] - ia[i];
1093b1a0a8a3SJed Brown 	  row = ja + ia[i];
1094b1a0a8a3SJed Brown 	  for (j=0; j<len; j++) {
1095b1a0a8a3SJed Brown 	    if (row[j] == i) { /* don't count diagonal */
1096b1a0a8a3SJed Brown 	      len--; break;
1097b1a0a8a3SJed Brown 	    }
1098b1a0a8a3SJed Brown 	  }
1099b1a0a8a3SJed Brown 	  nnz += len;
1100b1a0a8a3SJed Brown 	}
1101b1a0a8a3SJed Brown 	ierr = PetscMalloc((na+1)*sizeof(PetscInt),&iia);CHKERRQ(ierr);
1102b1a0a8a3SJed Brown 	ierr = PetscMalloc((nnz)*sizeof(PetscInt),&jja);CHKERRQ(ierr);
1103b1a0a8a3SJed Brown 	nnz    = 0;
1104b1a0a8a3SJed Brown 	iia[0] = 0;
1105b1a0a8a3SJed Brown 	for (i=0; i<na; i++) { /* fill adjacency */
1106b1a0a8a3SJed Brown 	  cnt = 0;
1107b1a0a8a3SJed Brown 	  len = ia[i+1] - ia[i];
1108b1a0a8a3SJed Brown 	  row = ja + ia[i];
1109b1a0a8a3SJed Brown 	  for (j=0; j<len; j++) {
1110b1a0a8a3SJed Brown 	    if (row[j] != i) { /* if not diagonal */
1111b1a0a8a3SJed Brown 	      jja[nnz+cnt++] = row[j];
1112b1a0a8a3SJed Brown 	    }
1113b1a0a8a3SJed Brown 	  }
1114b1a0a8a3SJed Brown 	  nnz += cnt;
1115b1a0a8a3SJed Brown 	  iia[i+1] = nnz;
1116b1a0a8a3SJed Brown 	}
1117b1a0a8a3SJed Brown 	/* Partitioning of the adjacency matrix */
1118b1a0a8a3SJed Brown 	ierr = MatCreateMPIAdj(PETSC_COMM_SELF,na,na,iia,jja,PETSC_NULL,&adj);CHKERRQ(ierr);
1119b1a0a8a3SJed Brown 	ierr = MatPartitioningSetAdjacency(mpart,adj);CHKERRQ(ierr);
1120b1a0a8a3SJed Brown 	ierr = MatPartitioningSetNParts(mpart,n);CHKERRQ(ierr);
1121b1a0a8a3SJed Brown 	ierr = MatPartitioningApply(mpart,&ispart);CHKERRQ(ierr);
1122b1a0a8a3SJed Brown 	ierr = ISPartitioningToNumbering(ispart,&isnumb);CHKERRQ(ierr);
1123b1a0a8a3SJed Brown 	ierr = MatDestroy(adj);CHKERRQ(ierr);
1124b1a0a8a3SJed Brown 	foundpart = PETSC_TRUE;
1125b1a0a8a3SJed Brown       }
1126b1a0a8a3SJed Brown       ierr = MatRestoreRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr);
1127b1a0a8a3SJed Brown     }
1128b1a0a8a3SJed Brown     ierr = MatPartitioningDestroy(mpart);CHKERRQ(ierr);
1129b1a0a8a3SJed Brown   }
1130b1a0a8a3SJed Brown   if (iscopy) {ierr = MatDestroy(Ad);CHKERRQ(ierr);}
1131b1a0a8a3SJed Brown 
1132b1a0a8a3SJed Brown   ierr = PetscMalloc(n*sizeof(IS),&is);CHKERRQ(ierr);
1133b1a0a8a3SJed Brown   *outis = is;
1134b1a0a8a3SJed Brown 
1135b1a0a8a3SJed Brown   if (!foundpart) {
1136b1a0a8a3SJed Brown 
1137b1a0a8a3SJed Brown     /* Partitioning by contiguous chunks of rows */
1138b1a0a8a3SJed Brown 
1139b1a0a8a3SJed Brown     PetscInt mbs   = (rend-rstart)/bs;
1140b1a0a8a3SJed Brown     PetscInt start = rstart;
1141b1a0a8a3SJed Brown     for (i=0; i<n; i++) {
1142b1a0a8a3SJed Brown       PetscInt count = (mbs/n + ((mbs % n) > i)) * bs;
1143b1a0a8a3SJed Brown       ierr   = ISCreateStride(PETSC_COMM_SELF,count,start,1,&is[i]);CHKERRQ(ierr);
1144b1a0a8a3SJed Brown       start += count;
1145b1a0a8a3SJed Brown     }
1146b1a0a8a3SJed Brown 
1147b1a0a8a3SJed Brown   } else {
1148b1a0a8a3SJed Brown 
1149b1a0a8a3SJed Brown     /* Partitioning by adjacency of diagonal block  */
1150b1a0a8a3SJed Brown 
1151b1a0a8a3SJed Brown     const PetscInt *numbering;
1152b1a0a8a3SJed Brown     PetscInt       *count,nidx,*indices,*newidx,start=0;
1153b1a0a8a3SJed Brown     /* Get node count in each partition */
1154b1a0a8a3SJed Brown     ierr = PetscMalloc(n*sizeof(PetscInt),&count);CHKERRQ(ierr);
1155b1a0a8a3SJed Brown     ierr = ISPartitioningCount(ispart,n,count);CHKERRQ(ierr);
1156b1a0a8a3SJed Brown     if (isbaij && bs > 1) { /* adjust for the block-aij case */
1157b1a0a8a3SJed Brown       for (i=0; i<n; i++) count[i] *= bs;
1158b1a0a8a3SJed Brown     }
1159b1a0a8a3SJed Brown     /* Build indices from node numbering */
1160b1a0a8a3SJed Brown     ierr = ISGetLocalSize(isnumb,&nidx);CHKERRQ(ierr);
1161b1a0a8a3SJed Brown     ierr = PetscMalloc(nidx*sizeof(PetscInt),&indices);CHKERRQ(ierr);
1162b1a0a8a3SJed Brown     for (i=0; i<nidx; i++) indices[i] = i; /* needs to be initialized */
1163b1a0a8a3SJed Brown     ierr = ISGetIndices(isnumb,&numbering);CHKERRQ(ierr);
1164b1a0a8a3SJed Brown     ierr = PetscSortIntWithPermutation(nidx,numbering,indices);CHKERRQ(ierr);
1165b1a0a8a3SJed Brown     ierr = ISRestoreIndices(isnumb,&numbering);CHKERRQ(ierr);
1166b1a0a8a3SJed Brown     if (isbaij && bs > 1) { /* adjust for the block-aij case */
1167b1a0a8a3SJed Brown       ierr = PetscMalloc(nidx*bs*sizeof(PetscInt),&newidx);CHKERRQ(ierr);
1168b1a0a8a3SJed Brown       for (i=0; i<nidx; i++)
1169b1a0a8a3SJed Brown 	for (j=0; j<bs; j++)
1170b1a0a8a3SJed Brown 	  newidx[i*bs+j] = indices[i]*bs + j;
1171b1a0a8a3SJed Brown       ierr = PetscFree(indices);CHKERRQ(ierr);
1172b1a0a8a3SJed Brown       nidx   *= bs;
1173b1a0a8a3SJed Brown       indices = newidx;
1174b1a0a8a3SJed Brown     }
1175b1a0a8a3SJed Brown     /* Shift to get global indices */
1176b1a0a8a3SJed Brown     for (i=0; i<nidx; i++) indices[i] += rstart;
1177b1a0a8a3SJed Brown 
1178b1a0a8a3SJed Brown     /* Build the index sets for each block */
1179b1a0a8a3SJed Brown     for (i=0; i<n; i++) {
1180b1a0a8a3SJed Brown       ierr   = ISCreateGeneral(PETSC_COMM_SELF,count[i],&indices[start],PETSC_COPY_VALUES,&is[i]);CHKERRQ(ierr);
1181b1a0a8a3SJed Brown       ierr   = ISSort(is[i]);CHKERRQ(ierr);
1182b1a0a8a3SJed Brown       start += count[i];
1183b1a0a8a3SJed Brown     }
1184b1a0a8a3SJed Brown 
1185b1a0a8a3SJed Brown     ierr = PetscFree(count);
1186b1a0a8a3SJed Brown     ierr = PetscFree(indices);
1187b1a0a8a3SJed Brown     ierr = ISDestroy(isnumb);CHKERRQ(ierr);
1188b1a0a8a3SJed Brown     ierr = ISDestroy(ispart);CHKERRQ(ierr);
1189b1a0a8a3SJed Brown 
1190b1a0a8a3SJed Brown   }
1191b1a0a8a3SJed Brown 
1192b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1193b1a0a8a3SJed Brown }
1194b1a0a8a3SJed Brown 
1195b1a0a8a3SJed Brown #undef __FUNCT__
1196f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMDestroySubdomains"
1197b1a0a8a3SJed Brown /*@C
1198f746d493SDmitry Karpeev    PCGASMDestroySubdomains - Destroys the index sets created with
1199f746d493SDmitry Karpeev    PCGASMCreateSubdomains(). Should be called after setting subdomains
1200f746d493SDmitry Karpeev    with PCGASMSetLocalSubdomains().
1201b1a0a8a3SJed Brown 
1202b1a0a8a3SJed Brown    Collective
1203b1a0a8a3SJed Brown 
1204b1a0a8a3SJed Brown    Input Parameters:
1205b1a0a8a3SJed Brown +  n - the number of index sets
1206b1a0a8a3SJed Brown .  is - the array of index sets
1207b1a0a8a3SJed Brown -  is_local - the array of local index sets, can be PETSC_NULL
1208b1a0a8a3SJed Brown 
1209b1a0a8a3SJed Brown    Level: advanced
1210b1a0a8a3SJed Brown 
1211f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, unstructured grid
1212b1a0a8a3SJed Brown 
1213f746d493SDmitry Karpeev .seealso: PCGASMCreateSubdomains(), PCGASMSetLocalSubdomains()
1214b1a0a8a3SJed Brown @*/
1215f746d493SDmitry Karpeev PetscErrorCode PETSCKSP_DLLEXPORT PCGASMDestroySubdomains(PetscInt n, IS is[], IS is_local[])
1216b1a0a8a3SJed Brown {
1217b1a0a8a3SJed Brown   PetscInt       i;
1218b1a0a8a3SJed Brown   PetscErrorCode ierr;
1219b1a0a8a3SJed Brown   PetscFunctionBegin;
1220b1a0a8a3SJed Brown   if (n <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"n must be > 0: n = %D",n);
1221b1a0a8a3SJed Brown   PetscValidPointer(is,2);
1222b1a0a8a3SJed Brown   for (i=0; i<n; i++) { ierr = ISDestroy(is[i]);CHKERRQ(ierr); }
1223b1a0a8a3SJed Brown   ierr = PetscFree(is);CHKERRQ(ierr);
1224b1a0a8a3SJed Brown   if (is_local) {
1225b1a0a8a3SJed Brown     PetscValidPointer(is_local,3);
1226b1a0a8a3SJed Brown     for (i=0; i<n; i++) { ierr = ISDestroy(is_local[i]);CHKERRQ(ierr); }
1227b1a0a8a3SJed Brown     ierr = PetscFree(is_local);CHKERRQ(ierr);
1228b1a0a8a3SJed Brown   }
1229b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1230b1a0a8a3SJed Brown }
1231b1a0a8a3SJed Brown 
1232b1a0a8a3SJed Brown #undef __FUNCT__
1233f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMCreateSubdomains2D"
1234b1a0a8a3SJed Brown /*@
1235f746d493SDmitry Karpeev    PCGASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz
1236b1a0a8a3SJed Brown    preconditioner for a two-dimensional problem on a regular grid.
1237b1a0a8a3SJed Brown 
1238b1a0a8a3SJed Brown    Not Collective
1239b1a0a8a3SJed Brown 
1240b1a0a8a3SJed Brown    Input Parameters:
1241b1a0a8a3SJed Brown +  m, n - the number of mesh points in the x and y directions
1242b1a0a8a3SJed Brown .  M, N - the number of subdomains in the x and y directions
1243b1a0a8a3SJed Brown .  dof - degrees of freedom per node
1244b1a0a8a3SJed Brown -  overlap - overlap in mesh lines
1245b1a0a8a3SJed Brown 
1246b1a0a8a3SJed Brown    Output Parameters:
1247b1a0a8a3SJed Brown +  Nsub - the number of subdomains created
1248b1a0a8a3SJed Brown .  is - array of index sets defining overlapping (if overlap > 0) subdomains
1249b1a0a8a3SJed Brown -  is_local - array of index sets defining non-overlapping subdomains
1250b1a0a8a3SJed Brown 
1251b1a0a8a3SJed Brown    Note:
1252b1a0a8a3SJed Brown    Presently PCAMSCreateSubdomains2d() is valid only for sequential
1253b1a0a8a3SJed Brown    preconditioners.  More general related routines are
1254f746d493SDmitry Karpeev    PCGASMSetTotalSubdomains() and PCGASMSetLocalSubdomains().
1255b1a0a8a3SJed Brown 
1256b1a0a8a3SJed Brown    Level: advanced
1257b1a0a8a3SJed Brown 
1258f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, 2D, regular grid
1259b1a0a8a3SJed Brown 
1260f746d493SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetLocalSubdomains(), PCGASMGetSubKSP(),
1261f746d493SDmitry Karpeev           PCGASMSetOverlap()
1262b1a0a8a3SJed Brown @*/
1263f746d493SDmitry Karpeev PetscErrorCode PETSCKSP_DLLEXPORT PCGASMCreateSubdomains2D(PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt dof,PetscInt overlap,PetscInt *Nsub,IS **is,IS **is_local)
1264b1a0a8a3SJed Brown {
1265b1a0a8a3SJed Brown   PetscInt       i,j,height,width,ystart,xstart,yleft,yright,xleft,xright,loc_outer;
1266b1a0a8a3SJed Brown   PetscErrorCode ierr;
1267b1a0a8a3SJed Brown   PetscInt       nidx,*idx,loc,ii,jj,count;
1268b1a0a8a3SJed Brown 
1269b1a0a8a3SJed Brown   PetscFunctionBegin;
1270b1a0a8a3SJed Brown   if (dof != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP," ");
1271b1a0a8a3SJed Brown 
1272b1a0a8a3SJed Brown   *Nsub = N*M;
1273b1a0a8a3SJed Brown   ierr = PetscMalloc((*Nsub)*sizeof(IS*),is);CHKERRQ(ierr);
1274b1a0a8a3SJed Brown   ierr = PetscMalloc((*Nsub)*sizeof(IS*),is_local);CHKERRQ(ierr);
1275b1a0a8a3SJed Brown   ystart = 0;
1276b1a0a8a3SJed Brown   loc_outer = 0;
1277b1a0a8a3SJed Brown   for (i=0; i<N; i++) {
1278b1a0a8a3SJed Brown     height = n/N + ((n % N) > i); /* height of subdomain */
1279b1a0a8a3SJed Brown     if (height < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many N subdomains for mesh dimension n");
1280b1a0a8a3SJed Brown     yleft  = ystart - overlap; if (yleft < 0) yleft = 0;
1281b1a0a8a3SJed Brown     yright = ystart + height + overlap; if (yright > n) yright = n;
1282b1a0a8a3SJed Brown     xstart = 0;
1283b1a0a8a3SJed Brown     for (j=0; j<M; j++) {
1284b1a0a8a3SJed Brown       width = m/M + ((m % M) > j); /* width of subdomain */
1285b1a0a8a3SJed Brown       if (width < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many M subdomains for mesh dimension m");
1286b1a0a8a3SJed Brown       xleft  = xstart - overlap; if (xleft < 0) xleft = 0;
1287b1a0a8a3SJed Brown       xright = xstart + width + overlap; if (xright > m) xright = m;
1288b1a0a8a3SJed Brown       nidx   = (xright - xleft)*(yright - yleft);
1289b1a0a8a3SJed Brown       ierr = PetscMalloc(nidx*sizeof(PetscInt),&idx);CHKERRQ(ierr);
1290b1a0a8a3SJed Brown       loc    = 0;
1291b1a0a8a3SJed Brown       for (ii=yleft; ii<yright; ii++) {
1292b1a0a8a3SJed Brown         count = m*ii + xleft;
1293b1a0a8a3SJed Brown         for (jj=xleft; jj<xright; jj++) {
1294b1a0a8a3SJed Brown           idx[loc++] = count++;
1295b1a0a8a3SJed Brown         }
1296b1a0a8a3SJed Brown       }
1297b1a0a8a3SJed Brown       ierr = ISCreateGeneral(PETSC_COMM_SELF,nidx,idx,PETSC_COPY_VALUES,(*is)+loc_outer);CHKERRQ(ierr);
1298b1a0a8a3SJed Brown       if (overlap == 0) {
1299b1a0a8a3SJed Brown         ierr = PetscObjectReference((PetscObject)(*is)[loc_outer]);CHKERRQ(ierr);
1300b1a0a8a3SJed Brown         (*is_local)[loc_outer] = (*is)[loc_outer];
1301b1a0a8a3SJed Brown       } else {
1302b1a0a8a3SJed Brown         for (loc=0,ii=ystart; ii<ystart+height; ii++) {
1303b1a0a8a3SJed Brown           for (jj=xstart; jj<xstart+width; jj++) {
1304b1a0a8a3SJed Brown             idx[loc++] = m*ii + jj;
1305b1a0a8a3SJed Brown           }
1306b1a0a8a3SJed Brown         }
1307b1a0a8a3SJed Brown         ierr = ISCreateGeneral(PETSC_COMM_SELF,loc,idx,PETSC_COPY_VALUES,*is_local+loc_outer);CHKERRQ(ierr);
1308b1a0a8a3SJed Brown       }
1309b1a0a8a3SJed Brown       ierr = PetscFree(idx);CHKERRQ(ierr);
1310b1a0a8a3SJed Brown       xstart += width;
1311b1a0a8a3SJed Brown       loc_outer++;
1312b1a0a8a3SJed Brown     }
1313b1a0a8a3SJed Brown     ystart += height;
1314b1a0a8a3SJed Brown   }
1315b1a0a8a3SJed Brown   for (i=0; i<*Nsub; i++) { ierr = ISSort((*is)[i]);CHKERRQ(ierr); }
1316b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1317b1a0a8a3SJed Brown }
1318b1a0a8a3SJed Brown 
1319b1a0a8a3SJed Brown #undef __FUNCT__
1320f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMGetLocalSubdomains"
1321b1a0a8a3SJed Brown /*@C
1322f746d493SDmitry Karpeev     PCGASMGetLocalSubdomains - Gets the local subdomains (for this processor
1323b1a0a8a3SJed Brown     only) for the additive Schwarz preconditioner.
1324b1a0a8a3SJed Brown 
1325b1a0a8a3SJed Brown     Not Collective
1326b1a0a8a3SJed Brown 
1327b1a0a8a3SJed Brown     Input Parameter:
1328b1a0a8a3SJed Brown .   pc - the preconditioner context
1329b1a0a8a3SJed Brown 
1330b1a0a8a3SJed Brown     Output Parameters:
1331b1a0a8a3SJed Brown +   n - the number of subdomains for this processor (default value = 1)
1332b1a0a8a3SJed Brown .   is - the index sets that define the subdomains for this processor
1333b1a0a8a3SJed Brown -   is_local - the index sets that define the local part of the subdomains for this processor (can be PETSC_NULL)
1334b1a0a8a3SJed Brown 
1335b1a0a8a3SJed Brown 
1336b1a0a8a3SJed Brown     Notes:
1337b1a0a8a3SJed Brown     The IS numbering is in the parallel, global numbering of the vector.
1338b1a0a8a3SJed Brown 
1339b1a0a8a3SJed Brown     Level: advanced
1340b1a0a8a3SJed Brown 
1341f746d493SDmitry Karpeev .keywords: PC, GASM, set, local, subdomains, additive Schwarz
1342b1a0a8a3SJed Brown 
1343f746d493SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(),
1344f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D(), PCGASMSetLocalSubdomains(), PCGASMGetLocalSubmatrices()
1345b1a0a8a3SJed Brown @*/
1346f746d493SDmitry Karpeev PetscErrorCode PETSCKSP_DLLEXPORT PCGASMGetLocalSubdomains(PC pc,PetscInt *n,IS *is[],IS *is_local[])
1347b1a0a8a3SJed Brown {
1348f746d493SDmitry Karpeev   PC_GASM         *osm;
1349b1a0a8a3SJed Brown   PetscErrorCode ierr;
1350b1a0a8a3SJed Brown   PetscBool      match;
1351b1a0a8a3SJed Brown 
1352b1a0a8a3SJed Brown   PetscFunctionBegin;
1353b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1354b1a0a8a3SJed Brown   PetscValidIntPointer(n,2);
1355b1a0a8a3SJed Brown   if (is) PetscValidPointer(is,3);
1356f746d493SDmitry Karpeev   ierr = PetscTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr);
1357b1a0a8a3SJed Brown   if (!match) {
1358b1a0a8a3SJed Brown     if (n)  *n  = 0;
1359b1a0a8a3SJed Brown     if (is) *is = PETSC_NULL;
1360b1a0a8a3SJed Brown   } else {
1361f746d493SDmitry Karpeev     osm = (PC_GASM*)pc->data;
1362ab3e923fSDmitry Karpeev     if (n)  *n  = osm->n;
1363b1a0a8a3SJed Brown     if (is) *is = osm->is;
1364b1a0a8a3SJed Brown     if (is_local) *is_local = osm->is_local;
1365b1a0a8a3SJed Brown   }
1366b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1367b1a0a8a3SJed Brown }
1368b1a0a8a3SJed Brown 
1369b1a0a8a3SJed Brown #undef __FUNCT__
1370f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMGetLocalSubmatrices"
1371b1a0a8a3SJed Brown /*@C
1372f746d493SDmitry Karpeev     PCGASMGetLocalSubmatrices - Gets the local submatrices (for this processor
1373b1a0a8a3SJed Brown     only) for the additive Schwarz preconditioner.
1374b1a0a8a3SJed Brown 
1375b1a0a8a3SJed Brown     Not Collective
1376b1a0a8a3SJed Brown 
1377b1a0a8a3SJed Brown     Input Parameter:
1378b1a0a8a3SJed Brown .   pc - the preconditioner context
1379b1a0a8a3SJed Brown 
1380b1a0a8a3SJed Brown     Output Parameters:
1381b1a0a8a3SJed Brown +   n - the number of matrices for this processor (default value = 1)
1382b1a0a8a3SJed Brown -   mat - the matrices
1383b1a0a8a3SJed Brown 
1384b1a0a8a3SJed Brown 
1385b1a0a8a3SJed Brown     Level: advanced
1386b1a0a8a3SJed Brown 
1387f746d493SDmitry Karpeev .keywords: PC, GASM, set, local, subdomains, additive Schwarz, block Jacobi
1388b1a0a8a3SJed Brown 
1389f746d493SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(),
1390f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D(), PCGASMSetLocalSubdomains(), PCGASMGetLocalSubdomains()
1391b1a0a8a3SJed Brown @*/
1392f746d493SDmitry Karpeev PetscErrorCode PETSCKSP_DLLEXPORT PCGASMGetLocalSubmatrices(PC pc,PetscInt *n,Mat *mat[])
1393b1a0a8a3SJed Brown {
1394f746d493SDmitry Karpeev   PC_GASM         *osm;
1395b1a0a8a3SJed Brown   PetscErrorCode ierr;
1396b1a0a8a3SJed Brown   PetscBool      match;
1397b1a0a8a3SJed Brown 
1398b1a0a8a3SJed Brown   PetscFunctionBegin;
1399b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1400b1a0a8a3SJed Brown   PetscValidIntPointer(n,2);
1401b1a0a8a3SJed Brown   if (mat) PetscValidPointer(mat,3);
1402b1a0a8a3SJed Brown   if (!pc->setupcalled) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must call after KSPSetUP() or PCSetUp().");
1403f746d493SDmitry Karpeev   ierr = PetscTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr);
1404b1a0a8a3SJed Brown   if (!match) {
1405b1a0a8a3SJed Brown     if (n)   *n   = 0;
1406b1a0a8a3SJed Brown     if (mat) *mat = PETSC_NULL;
1407b1a0a8a3SJed Brown   } else {
1408f746d493SDmitry Karpeev     osm = (PC_GASM*)pc->data;
1409ab3e923fSDmitry Karpeev     if (n)   *n   = osm->n;
1410b1a0a8a3SJed Brown     if (mat) *mat = osm->pmat;
1411b1a0a8a3SJed Brown   }
1412b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1413b1a0a8a3SJed Brown }
1414