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