xref: /petsc/src/ksp/pc/impls/gasm/gasm.c (revision b6b0974b5449c9bbf556d702c62f36941bd51fa0)
1b1a0a8a3SJed Brown 
2b1a0a8a3SJed Brown /*
3f746d493SDmitry Karpeev   This file defines an "generalized" additive Schwarz preconditioner for any Mat implementation.
444fe18e5SDmitry Karpeev   In this version each processor may have any number of subdomains and a subdomain may live on multiple
5f746d493SDmitry Karpeev   processors.
6b1a0a8a3SJed Brown 
7af538404SDmitry Karpeev        N    - total number of subdomains on all processors
8f746d493SDmitry Karpeev        n    - actual number of subdomains on this processor
9f746d493SDmitry Karpeev        nmax - maximum number of subdomains per processor
10b1a0a8a3SJed Brown */
11b45d2f2cSJed Brown #include <petsc-private/pcimpl.h>     /*I "petscpc.h" I*/
12b1a0a8a3SJed Brown 
13b1a0a8a3SJed Brown typedef struct {
14f746d493SDmitry Karpeev   PetscInt   N,n,nmax;
15b1a0a8a3SJed Brown   PetscInt   overlap;             /* overlap requested by user */
16b1a0a8a3SJed Brown   KSP        *ksp;                /* linear solvers for each block */
17f746d493SDmitry Karpeev   Vec        gx,gy;               /* Merged work vectors */
18f746d493SDmitry Karpeev   Vec        *x,*y;               /* Split work vectors; storage aliases pieces of storage of the above merged vectors. */
19f746d493SDmitry Karpeev   IS         gis, gis_local;      /* merged ISs */
20f746d493SDmitry Karpeev   VecScatter grestriction;        /* merged restriction */
21f746d493SDmitry Karpeev   VecScatter gprolongation;       /* merged prolongation */
22b1a0a8a3SJed Brown   IS         *is;                 /* index set that defines each overlapping subdomain */
23f746d493SDmitry Karpeev   IS         *is_local;           /* index set that defines each local subdomain (same as subdomain with the overlap removed); may be NULL */
24f746d493SDmitry Karpeev   Mat        *pmat;               /* subdomain block matrices */
25f746d493SDmitry Karpeev   PCGASMType  type;               /* use reduced interpolation, restriction or both */
26b1a0a8a3SJed Brown   PetscBool  type_set;            /* if user set this value (so won't change it for symmetric problems) */
27b1a0a8a3SJed Brown   PetscBool  same_local_solves;   /* flag indicating whether all local solvers are same */
28b1a0a8a3SJed Brown   PetscBool  sort_indices;        /* flag to sort subdomain indices */
29f746d493SDmitry Karpeev } PC_GASM;
30b1a0a8a3SJed Brown 
31b1a0a8a3SJed Brown #undef __FUNCT__
32af538404SDmitry Karpeev #define __FUNCT__ "PCGASMPrintSubdomains"
33af538404SDmitry Karpeev static PetscErrorCode  PCGASMPrintSubdomains(PC pc)
34af538404SDmitry Karpeev {
35af538404SDmitry Karpeev   PC_GASM         *osm  = (PC_GASM*)pc->data;
36af538404SDmitry Karpeev   const char     *prefix;
37af538404SDmitry Karpeev   char           fname[PETSC_MAX_PATH_LEN+1];
38af538404SDmitry Karpeev   PetscViewer    viewer, sviewer;
39af538404SDmitry Karpeev   PetscInt       i,j,nidx;
40af538404SDmitry Karpeev   const PetscInt *idx;
41af538404SDmitry Karpeev   char           *cidx;
42af538404SDmitry Karpeev   PetscBool      found;
43eec7e3faSDmitry Karpeev   PetscMPIInt    size, rank;
44af538404SDmitry Karpeev   PetscErrorCode ierr;
45af538404SDmitry Karpeev 
46af538404SDmitry Karpeev   PetscFunctionBegin;
47eec7e3faSDmitry Karpeev   ierr = MPI_Comm_size(((PetscObject)pc)->comm, &size);CHKERRQ(ierr);
48eec7e3faSDmitry Karpeev   ierr = MPI_Comm_rank(((PetscObject)pc)->comm, &rank);CHKERRQ(ierr);
49af538404SDmitry Karpeev   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
50af538404SDmitry Karpeev   ierr = PetscOptionsGetString(prefix,"-pc_gasm_print_subdomains",fname,PETSC_MAX_PATH_LEN,&found);CHKERRQ(ierr);
51af538404SDmitry Karpeev   if (!found) { ierr = PetscStrcpy(fname,"stdout");CHKERRQ(ierr); };
52af538404SDmitry Karpeev   for (i=0;i<osm->n;++i) {
53af538404SDmitry Karpeev     ierr = PetscViewerASCIIOpen(((PetscObject)((osm->is)[i]))->comm,fname,&viewer);CHKERRQ(ierr);
54af538404SDmitry Karpeev     ierr = ISGetLocalSize(osm->is[i], &nidx);CHKERRQ(ierr);
55eec7e3faSDmitry Karpeev     /*
56eec7e3faSDmitry Karpeev      No more than 15 characters per index plus a space.
57eec7e3faSDmitry Karpeev      PetscViewerStringSPrintf requires a string of size at least 2, so use (nidx+1) instead of nidx,
58eec7e3faSDmitry Karpeev      in case nidx == 0. That will take care of the space for the trailing '\0' as well.
59eec7e3faSDmitry Karpeev      For nidx == 0, the whole string 16 '\0'.
60eec7e3faSDmitry Karpeev      */
61eec7e3faSDmitry Karpeev     ierr = PetscMalloc(sizeof(char)*(16*(nidx+1)+1), &cidx);CHKERRQ(ierr);
62af538404SDmitry Karpeev     ierr = ISGetIndices(osm->is[i], &idx);CHKERRQ(ierr);
63eec7e3faSDmitry Karpeev     ierr = PetscViewerStringOpen(((PetscObject)(osm->is[i]))->comm, cidx, 16*(nidx+1), &sviewer);CHKERRQ(ierr);
64af538404SDmitry Karpeev     for(j = 0; j < nidx; ++j) {
65af538404SDmitry Karpeev       ierr = PetscViewerStringSPrintf(sviewer, "%D ", idx[j]);CHKERRQ(ierr);
66af538404SDmitry Karpeev     }
676bf464f9SBarry Smith     ierr = PetscViewerDestroy(&sviewer);CHKERRQ(ierr);
68af538404SDmitry Karpeev     ierr = ISRestoreIndices(osm->is[i],&idx);CHKERRQ(ierr);
69eec7e3faSDmitry Karpeev 
707b23a99aSBarry Smith     ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
71af538404SDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer, "Subdomain with overlap\n");CHKERRQ(ierr);
72af538404SDmitry Karpeev     ierr = PetscViewerASCIISynchronizedPrintf(viewer, "%s", cidx);CHKERRQ(ierr);
73af538404SDmitry Karpeev     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
747b23a99aSBarry Smith     ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
75af538404SDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr);
76af538404SDmitry Karpeev     ierr = PetscFree(cidx);CHKERRQ(ierr);
776bf464f9SBarry Smith     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
78af538404SDmitry Karpeev     if (osm->is_local) {
79af538404SDmitry Karpeev       ierr = PetscViewerASCIIOpen(((PetscObject)((osm->is)[i]))->comm,fname,&viewer);CHKERRQ(ierr);
80af538404SDmitry Karpeev       ierr = ISGetLocalSize(osm->is_local[i], &nidx);CHKERRQ(ierr);
81eec7e3faSDmitry Karpeev       /*
82eec7e3faSDmitry Karpeev        No more than 15 characters per index plus a space.
83eec7e3faSDmitry Karpeev        PetscViewerStringSPrintf requires a string of size at least 2, so use (nidx+1) instead of nidx,
84eec7e3faSDmitry Karpeev        in case nidx == 0. That will take care of the space for the trailing '\0' as well.
85eec7e3faSDmitry Karpeev        For nidx == 0, the whole string 16 '\0'.
86eec7e3faSDmitry Karpeev        */
87eec7e3faSDmitry Karpeev       ierr = PetscMalloc(sizeof(char)*(16*(nidx+1)+1), &cidx);CHKERRQ(ierr);
88af538404SDmitry Karpeev       ierr = ISGetIndices(osm->is_local[i], &idx);CHKERRQ(ierr);
89eec7e3faSDmitry Karpeev       ierr = PetscViewerStringOpen(((PetscObject)(osm->is_local[i]))->comm, cidx, 16*(nidx+1), &sviewer);CHKERRQ(ierr);
90af538404SDmitry Karpeev       for(j = 0; j < nidx; ++j) {
91af538404SDmitry Karpeev         ierr = PetscViewerStringSPrintf(sviewer, "%D ", idx[j]);CHKERRQ(ierr);
92af538404SDmitry Karpeev       }
936bf464f9SBarry Smith       ierr = PetscViewerDestroy(&sviewer);CHKERRQ(ierr);
94af538404SDmitry Karpeev       ierr = ISRestoreIndices(osm->is_local[i],&idx);CHKERRQ(ierr);
95eec7e3faSDmitry Karpeev 
967b23a99aSBarry Smith       ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
97af538404SDmitry Karpeev       ierr = PetscViewerASCIIPrintf(viewer, "Subdomain without overlap\n");CHKERRQ(ierr);
98af538404SDmitry Karpeev       ierr = PetscViewerASCIISynchronizedPrintf(viewer, "%s", cidx);CHKERRQ(ierr);
99af538404SDmitry Karpeev       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1007b23a99aSBarry Smith       ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
101af538404SDmitry Karpeev       ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr);
102af538404SDmitry Karpeev       ierr = PetscFree(cidx);CHKERRQ(ierr);
1036bf464f9SBarry Smith       ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
104af538404SDmitry Karpeev     }
105af538404SDmitry Karpeev   }
106af538404SDmitry Karpeev   PetscFunctionReturn(0);
107af538404SDmitry Karpeev }
108af538404SDmitry Karpeev 
109af538404SDmitry Karpeev 
110af538404SDmitry Karpeev #undef __FUNCT__
111f746d493SDmitry Karpeev #define __FUNCT__ "PCView_GASM"
112f746d493SDmitry Karpeev static PetscErrorCode PCView_GASM(PC pc,PetscViewer viewer)
113b1a0a8a3SJed Brown {
114f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
115af538404SDmitry Karpeev   const char     *prefix;
116b1a0a8a3SJed Brown   PetscErrorCode ierr;
117af538404SDmitry Karpeev   PetscMPIInt    rank, size;
118b1a0a8a3SJed Brown   PetscInt       i,bsz;
119e4b3f054SSatish Balay   PetscBool      iascii,isstring, print_subdomains=PETSC_FALSE;
120b1a0a8a3SJed Brown   PetscViewer    sviewer;
121b1a0a8a3SJed Brown 
122b1a0a8a3SJed Brown 
123b1a0a8a3SJed Brown   PetscFunctionBegin;
124af538404SDmitry Karpeev   ierr = MPI_Comm_size(((PetscObject)pc)->comm, &size);CHKERRQ(ierr);
125af538404SDmitry Karpeev   ierr = MPI_Comm_rank(((PetscObject)pc)->comm, &rank);CHKERRQ(ierr);
126b1a0a8a3SJed Brown   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
127b1a0a8a3SJed Brown   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
128af538404SDmitry Karpeev   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
129af538404SDmitry Karpeev   ierr = PetscOptionsGetBool(prefix,"-pc_gasm_print_subdomains",&print_subdomains,PETSC_NULL);CHKERRQ(ierr);
130b1a0a8a3SJed Brown   if (iascii) {
131af538404SDmitry Karpeev     char overlaps[256] = "user-defined overlap",subdomains[256] = "total subdomains set";
132b1a0a8a3SJed Brown     if (osm->overlap >= 0) {ierr = PetscSNPrintf(overlaps,sizeof overlaps,"amount of overlap = %D",osm->overlap);CHKERRQ(ierr);}
133af538404SDmitry Karpeev     if (osm->nmax > 0)     {ierr = PetscSNPrintf(subdomains,sizeof subdomains,"max number of local subdomains = %D",osm->nmax);CHKERRQ(ierr);}
1347b23a99aSBarry Smith     ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
135eec7e3faSDmitry Karpeev     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"  [%d:%d] number of locally-supported subdomains = %D\n",(int)rank,(int)size,osm->n);CHKERRQ(ierr);
136af538404SDmitry Karpeev     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1377b23a99aSBarry Smith     ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
138af538404SDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"  Generalized additive Schwarz: %s, %s\n",subdomains,overlaps);CHKERRQ(ierr);
139f746d493SDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"  Generalized additive Schwarz: restriction/interpolation type - %s\n",PCGASMTypes[osm->type]);CHKERRQ(ierr);
140b1a0a8a3SJed Brown     ierr = MPI_Comm_rank(((PetscObject)pc)->comm,&rank);CHKERRQ(ierr);
141b1a0a8a3SJed Brown     if (osm->same_local_solves) {
142b1a0a8a3SJed Brown       if (osm->ksp) {
143af538404SDmitry Karpeev         ierr = PetscViewerASCIIPrintf(viewer,"  Local solve is same for all subdomains, in the following KSP and PC objects:\n");CHKERRQ(ierr);
144b1a0a8a3SJed Brown         ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
145b1a0a8a3SJed Brown         if (!rank) {
146b1a0a8a3SJed Brown           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
147b1a0a8a3SJed Brown           ierr = KSPView(osm->ksp[0],sviewer);CHKERRQ(ierr);
148b1a0a8a3SJed Brown           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
149b1a0a8a3SJed Brown         }
150b1a0a8a3SJed Brown         ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
151b1a0a8a3SJed Brown       }
152b1a0a8a3SJed Brown     } else {
153b1a0a8a3SJed Brown       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
154af538404SDmitry Karpeev       ierr = PetscViewerASCIIPrintf(viewer,"  Local solve info for each subdomain is in the following KSP and PC objects:\n");CHKERRQ(ierr);
155b1a0a8a3SJed Brown       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
156b1a0a8a3SJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr);
157f746d493SDmitry Karpeev       for (i=0; i<osm->nmax; i++) {
158b1a0a8a3SJed Brown         ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
159f746d493SDmitry Karpeev         if (i < osm->n) {
160b1a0a8a3SJed Brown           ierr = ISGetLocalSize(osm->is[i],&bsz);CHKERRQ(ierr);
1617b23a99aSBarry Smith           ierr = PetscViewerASCIISynchronizedAllow(sviewer,PETSC_TRUE);CHKERRQ(ierr);
162af538404SDmitry Karpeev           ierr = PetscViewerASCIISynchronizedPrintf(sviewer,"[%d:%d] local subdomain number %D, size = %D\n",(int)rank,(int)size,i,bsz);CHKERRQ(ierr);
163b1a0a8a3SJed Brown           ierr = KSPView(osm->ksp[i],sviewer);CHKERRQ(ierr);
164b1a0a8a3SJed Brown           ierr = PetscViewerASCIISynchronizedPrintf(sviewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr);
1657b23a99aSBarry Smith           ierr = PetscViewerFlush(sviewer);CHKERRQ(ierr);
1667b23a99aSBarry Smith           ierr = PetscViewerASCIISynchronizedAllow(sviewer,PETSC_FALSE);CHKERRQ(ierr);
167b1a0a8a3SJed Brown         }
168b1a0a8a3SJed Brown         ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
169b1a0a8a3SJed Brown       }
170b1a0a8a3SJed Brown       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
171b1a0a8a3SJed Brown       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
172b1a0a8a3SJed Brown     }
173b1a0a8a3SJed Brown   } else if (isstring) {
174af538404SDmitry Karpeev     ierr = PetscViewerStringSPrintf(viewer," subdomains=%D, overlap=%D, type=%s",osm->n,osm->overlap,PCGASMTypes[osm->type]);CHKERRQ(ierr);
175b1a0a8a3SJed Brown     ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
176b1a0a8a3SJed Brown     if (osm->ksp) {ierr = KSPView(osm->ksp[0],sviewer);CHKERRQ(ierr);}
177b1a0a8a3SJed Brown     ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
178b1a0a8a3SJed Brown   } else {
179f746d493SDmitry Karpeev     SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCGASM",((PetscObject)viewer)->type_name);
180b1a0a8a3SJed Brown   }
181b1a0a8a3SJed Brown   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
182af538404SDmitry Karpeev   if(print_subdomains) {
183af538404SDmitry Karpeev     ierr = PCGASMPrintSubdomains(pc);CHKERRQ(ierr);
184af538404SDmitry Karpeev   }
185b1a0a8a3SJed Brown   PetscFunctionReturn(0);
186b1a0a8a3SJed Brown }
187b1a0a8a3SJed Brown 
188af538404SDmitry Karpeev 
189af538404SDmitry Karpeev 
190af538404SDmitry Karpeev 
191af538404SDmitry Karpeev 
192b1a0a8a3SJed Brown #undef __FUNCT__
193f746d493SDmitry Karpeev #define __FUNCT__ "PCSetUp_GASM"
194f746d493SDmitry Karpeev static PetscErrorCode PCSetUp_GASM(PC pc)
195b1a0a8a3SJed Brown {
196f746d493SDmitry Karpeev   PC_GASM         *osm  = (PC_GASM*)pc->data;
197b1a0a8a3SJed Brown   PetscErrorCode ierr;
198b1a0a8a3SJed Brown   PetscBool      symset,flg;
19988389755SJed Brown   PetscInt       i;
200c10bc1a1SDmitry Karpeev   PetscMPIInt    rank, size;
201b1a0a8a3SJed Brown   MatReuse       scall = MAT_REUSE_MATRIX;
202b1a0a8a3SJed Brown   KSP            ksp;
203b1a0a8a3SJed Brown   PC             subpc;
204b1a0a8a3SJed Brown   const char     *prefix,*pprefix;
205f746d493SDmitry Karpeev   PetscInt       dn;       /* Number of indices in a single subdomain assigned to this processor. */
206f746d493SDmitry Karpeev   const PetscInt *didx;    /* Indices from a single subdomain assigned to this processor. */
207f746d493SDmitry Karpeev   PetscInt       ddn;      /* Number of indices in all subdomains assigned to this processor. */
208f746d493SDmitry Karpeev   PetscInt       *ddidx;   /* Indices of all subdomains assigned to this processor. */
209f746d493SDmitry Karpeev   IS             gid;      /* Identity IS of the size of all subdomains assigned to this processor. */
210f746d493SDmitry Karpeev   Vec            x,y;
211f746d493SDmitry Karpeev   PetscScalar    *gxarray, *gyarray;
2126d98aee9SDmitry Karpeev   PetscInt       gfirst, glast;
213fdc48646SDmitry Karpeev   DM             *domain_dm = PETSC_NULL;
214b1a0a8a3SJed Brown 
215b1a0a8a3SJed Brown   PetscFunctionBegin;
216c10bc1a1SDmitry Karpeev   ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr);
217c10bc1a1SDmitry Karpeev   ierr = MPI_Comm_rank(((PetscObject)pc)->comm,&rank);CHKERRQ(ierr);
218b1a0a8a3SJed Brown   if (!pc->setupcalled) {
219b1a0a8a3SJed Brown 
220b1a0a8a3SJed Brown     if (!osm->type_set) {
221b1a0a8a3SJed Brown       ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr);
222f746d493SDmitry Karpeev       if (symset && flg) { osm->type = PC_GASM_BASIC; }
223b1a0a8a3SJed Brown     }
224b1a0a8a3SJed Brown 
22544fe18e5SDmitry Karpeev     if (osm->N == PETSC_DECIDE && osm->n < 1) {
226fdc48646SDmitry Karpeev       /* no subdomains given, try pc->dm */
227fdc48646SDmitry Karpeev       if(pc->dm) {
228fdc48646SDmitry Karpeev         char      ddm_name[1024];
229fdc48646SDmitry Karpeev         DM        ddm;
230fdc48646SDmitry Karpeev         PetscBool flg;
231fdc48646SDmitry Karpeev         PetscInt     num_domains, d;
232fdc48646SDmitry Karpeev         char         **domain_names;
233fdc48646SDmitry Karpeev         IS           *domain_is;
234fdc48646SDmitry Karpeev         /* Allow the user to request a decomposition DM by name */
235fdc48646SDmitry Karpeev         ierr = PetscStrncpy(ddm_name, "", 1024); CHKERRQ(ierr);
236fdc48646SDmitry Karpeev         ierr = PetscOptionsString("-pc_asm_decomposition_dm", "Name of the DM defining the composition", "PCSetDM", ddm_name, ddm_name,1024,&flg); CHKERRQ(ierr);
237fdc48646SDmitry Karpeev         if(flg) {
238fdc48646SDmitry Karpeev           ierr = DMCreateDecompositionDM(pc->dm, ddm_name, &ddm); CHKERRQ(ierr);
239fdc48646SDmitry Karpeev           if(!ddm) {
240fdc48646SDmitry Karpeev             SETERRQ1(((PetscObject)pc)->comm, PETSC_ERR_ARG_WRONGSTATE, "Uknown DM decomposition name %s", ddm_name);
241fdc48646SDmitry Karpeev           }
242fdc48646SDmitry Karpeev           ierr = PetscInfo(pc,"Using decomposition DM defined using options database\n");CHKERRQ(ierr);
243fdc48646SDmitry Karpeev           ierr = PCSetDM(pc,ddm); CHKERRQ(ierr);
244fdc48646SDmitry Karpeev         }
245fdc48646SDmitry Karpeev         ierr = DMCreateDecomposition(pc->dm, &num_domains, &domain_names, &domain_is, &domain_dm);    CHKERRQ(ierr);
246fdc48646SDmitry Karpeev         ierr = PCGASMSetLocalSubdomains(pc, num_domains, domain_is, PETSC_NULL);CHKERRQ(ierr);
247fdc48646SDmitry Karpeev         for(d = 0; d < num_domains; ++d) {
248fdc48646SDmitry Karpeev           ierr = PetscFree(domain_names[d]); CHKERRQ(ierr);
249fdc48646SDmitry Karpeev           ierr = ISDestroy(&domain_is[d]);   CHKERRQ(ierr);
250fdc48646SDmitry Karpeev         }
251fdc48646SDmitry Karpeev         ierr = PetscFree(domain_names);CHKERRQ(ierr);
252fdc48646SDmitry Karpeev         ierr = PetscFree(domain_is);CHKERRQ(ierr);
253fdc48646SDmitry Karpeev       }
254fdc48646SDmitry Karpeev       else { /* no dm, use one per processor */
25544fe18e5SDmitry Karpeev         osm->nmax = osm->n = 1;
256b1a0a8a3SJed Brown         ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr);
257f746d493SDmitry Karpeev         osm->N = size;
258fdc48646SDmitry Karpeev       }
259f746d493SDmitry Karpeev     } else if (osm->N == PETSC_DECIDE) {
260b1a0a8a3SJed Brown       PetscInt inwork[2], outwork[2];
261f746d493SDmitry Karpeev       /* determine global number of subdomains and the max number of local subdomains */
262f746d493SDmitry Karpeev       inwork[0] = inwork[1] = osm->n;
263b1a0a8a3SJed Brown       ierr = MPI_Allreduce(inwork,outwork,1,MPIU_2INT,PetscMaxSum_Op,((PetscObject)pc)->comm);CHKERRQ(ierr);
264f746d493SDmitry Karpeev       osm->nmax = outwork[0];
265f746d493SDmitry Karpeev       osm->N    = outwork[1];
266b1a0a8a3SJed Brown     }
267b1a0a8a3SJed Brown     if (!osm->is){ /* create the index sets */
268f746d493SDmitry Karpeev       ierr = PCGASMCreateSubdomains(pc->pmat,osm->n,&osm->is);CHKERRQ(ierr);
269b1a0a8a3SJed Brown     }
27093c1ec32SDmitry Karpeev     if (!osm->is_local) {
27193c1ec32SDmitry Karpeev       /*
27293c1ec32SDmitry Karpeev 	 This indicates that osm->is should define a nonoverlapping decomposition
27393c1ec32SDmitry Karpeev 	 (there is no way to really guarantee that if subdomains are set by the user through PCGASMSetLocalSubdomains,
27493c1ec32SDmitry Karpeev 	  but the assumption is that either the user does the right thing, or subdomains in ossm->is have been created
27593c1ec32SDmitry Karpeev 	  via PCGASMCreateSubdomains, which guarantees a nonoverlapping decomposition).
27693c1ec32SDmitry Karpeev 	 Therefore, osm->is will be used to define osm->is_local.
27793c1ec32SDmitry Karpeev 	 If a nonzero overlap has been requested by the user, then osm->is will be expanded and will overlap,
27893c1ec32SDmitry Karpeev 	 so osm->is_local should obtain a copy of osm->is while they are still (presumably) nonoverlapping.
27993c1ec32SDmitry Karpeev 	 Otherwise (no overlap has been requested), osm->is_local are simply aliases for osm->is.
28093c1ec32SDmitry Karpeev       */
281f746d493SDmitry Karpeev       ierr = PetscMalloc(osm->n*sizeof(IS),&osm->is_local);CHKERRQ(ierr);
282f746d493SDmitry Karpeev       for (i=0; i<osm->n; i++) {
283b1a0a8a3SJed Brown         if (osm->overlap > 0) { /* With positive overlap, osm->is[i] will be modified */
284b1a0a8a3SJed Brown           ierr = ISDuplicate(osm->is[i],&osm->is_local[i]);CHKERRQ(ierr);
285b1a0a8a3SJed Brown           ierr = ISCopy(osm->is[i],osm->is_local[i]);CHKERRQ(ierr);
286b1a0a8a3SJed Brown         } else {
287b1a0a8a3SJed Brown           ierr = PetscObjectReference((PetscObject)osm->is[i]);CHKERRQ(ierr);
288b1a0a8a3SJed Brown           osm->is_local[i] = osm->is[i];
289b1a0a8a3SJed Brown         }
290b1a0a8a3SJed Brown       }
291b1a0a8a3SJed Brown     }
29293c1ec32SDmitry Karpeev     /* Beyond this point osm->is_local is not null. */
29393c1ec32SDmitry Karpeev     if (osm->overlap > 0) {
29493c1ec32SDmitry Karpeev       /* Extend the "overlapping" regions by a number of steps */
29593c1ec32SDmitry Karpeev       ierr = MatIncreaseOverlap(pc->pmat,osm->n,osm->is,osm->overlap);CHKERRQ(ierr);
29693c1ec32SDmitry Karpeev     }
297b1a0a8a3SJed Brown     ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
298b1a0a8a3SJed Brown     flg  = PETSC_FALSE;
299f746d493SDmitry Karpeev     ierr = PetscOptionsGetBool(prefix,"-pc_gasm_print_subdomains",&flg,PETSC_NULL);CHKERRQ(ierr);
300f746d493SDmitry Karpeev     if (flg) { ierr = PCGASMPrintSubdomains(pc);CHKERRQ(ierr); }
301b1a0a8a3SJed Brown 
302b1a0a8a3SJed Brown     if (osm->sort_indices) {
303f746d493SDmitry Karpeev       for (i=0; i<osm->n; i++) {
304b1a0a8a3SJed Brown         ierr = ISSort(osm->is[i]);CHKERRQ(ierr);
305b1a0a8a3SJed Brown 	ierr = ISSort(osm->is_local[i]);CHKERRQ(ierr);
306b1a0a8a3SJed Brown       }
307b1a0a8a3SJed Brown     }
308f746d493SDmitry Karpeev     /* Merge the ISs, create merged vectors and scatter contexts. */
30993c1ec32SDmitry Karpeev     /* Restriction ISs. */
310f746d493SDmitry Karpeev     ddn = 0;
311f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
312f746d493SDmitry Karpeev       ierr = ISGetLocalSize(osm->is[i],&dn);CHKERRQ(ierr);
313f746d493SDmitry Karpeev       ddn += dn;
314f746d493SDmitry Karpeev     }
315f746d493SDmitry Karpeev     ierr = PetscMalloc(ddn*sizeof(PetscInt), &ddidx);CHKERRQ(ierr);
316f746d493SDmitry Karpeev     ddn = 0;
317f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
318f746d493SDmitry Karpeev       ierr = ISGetLocalSize(osm->is[i],&dn);CHKERRQ(ierr);
319f746d493SDmitry Karpeev       ierr = ISGetIndices(osm->is[i],&didx);CHKERRQ(ierr);
320f746d493SDmitry Karpeev       ierr = PetscMemcpy(ddidx+ddn, didx, sizeof(PetscInt)*dn);CHKERRQ(ierr);
321f746d493SDmitry Karpeev       ierr = ISRestoreIndices(osm->is[i], &didx);CHKERRQ(ierr);
32293c1ec32SDmitry Karpeev       ddn += dn;
323f746d493SDmitry Karpeev     }
324f746d493SDmitry Karpeev     ierr = ISCreateGeneral(((PetscObject)(pc))->comm, ddn, ddidx, PETSC_OWN_POINTER, &osm->gis);CHKERRQ(ierr);
325f746d493SDmitry Karpeev     ierr = MatGetVecs(pc->pmat,&x,&y);CHKERRQ(ierr);
326f746d493SDmitry Karpeev     ierr = VecCreateMPI(((PetscObject)pc)->comm, ddn, PETSC_DECIDE, &osm->gx);CHKERRQ(ierr);
327f746d493SDmitry Karpeev     ierr = VecDuplicate(osm->gx,&osm->gy);CHKERRQ(ierr);
3286d98aee9SDmitry Karpeev     ierr = VecGetOwnershipRange(osm->gx, &gfirst, &glast);CHKERRQ(ierr);
3296d98aee9SDmitry Karpeev     ierr = ISCreateStride(((PetscObject)pc)->comm,ddn,gfirst,1, &gid);CHKERRQ(ierr);
330f746d493SDmitry Karpeev     ierr = VecScatterCreate(x,osm->gis,osm->gx,gid, &(osm->grestriction));CHKERRQ(ierr);
331fcfd50ebSBarry Smith     ierr = ISDestroy(&gid);CHKERRQ(ierr);
33293c1ec32SDmitry Karpeev     /* Prolongation ISs */
33381a5c4aaSDmitry Karpeev     { PetscInt       dn_local;       /* Number of indices in the local part of a single domain assigned to this processor. */
334f746d493SDmitry Karpeev       const PetscInt *didx_local;    /* Global indices from the local part of a single domain assigned to this processor. */
335f746d493SDmitry Karpeev       PetscInt       ddn_local;      /* Number of indices in the local part of the disjoint union all domains assigned to this processor. */
336f746d493SDmitry Karpeev       PetscInt       *ddidx_local;   /* Global indices of the local part of the disjoint union of all domains assigned to this processor. */
337f746d493SDmitry Karpeev       /**/
33893c1ec32SDmitry Karpeev       ISLocalToGlobalMapping ltog;          /* Map from global to local indices on the disjoint union of subdomains: "local" ind's run from 0 to ddn-1. */
33993c1ec32SDmitry Karpeev       PetscInt              *ddidx_llocal;  /* Mapped local indices of the disjoint union of local parts of subdomains. */
340f746d493SDmitry 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. */
341f746d493SDmitry Karpeev       IS                     gis_llocal;    /* IS with ddidx_llocal indices. */
342f746d493SDmitry Karpeev       PetscInt               j;
343f746d493SDmitry Karpeev       ddn_local = 0;
344f746d493SDmitry Karpeev       for (i=0; i<osm->n; i++) {
345f746d493SDmitry Karpeev 	ierr = ISGetLocalSize(osm->is_local[i],&dn_local);CHKERRQ(ierr);
346f746d493SDmitry Karpeev 	ddn_local += dn_local;
347f746d493SDmitry Karpeev       }
348f746d493SDmitry Karpeev       ierr = PetscMalloc(ddn_local*sizeof(PetscInt), &ddidx_local);CHKERRQ(ierr);
349f746d493SDmitry Karpeev       ddn_local = 0;
350f746d493SDmitry Karpeev       for (i=0; i<osm->n; i++) {
351f746d493SDmitry Karpeev 	ierr = ISGetLocalSize(osm->is_local[i],&dn_local);CHKERRQ(ierr);
352f746d493SDmitry Karpeev 	ierr = ISGetIndices(osm->is_local[i],&didx_local);CHKERRQ(ierr);
353f746d493SDmitry Karpeev 	ierr = PetscMemcpy(ddidx_local+ddn_local, didx_local, sizeof(PetscInt)*dn_local);CHKERRQ(ierr);
354f746d493SDmitry Karpeev 	ierr = ISRestoreIndices(osm->is_local[i], &didx_local);CHKERRQ(ierr);
355f746d493SDmitry Karpeev 	ddn_local += dn_local;
356f746d493SDmitry Karpeev       }
357ab3e923fSDmitry Karpeev       ierr = PetscMalloc(sizeof(PetscInt)*ddn_local, &ddidx_llocal);CHKERRQ(ierr);
358f746d493SDmitry Karpeev       ierr = ISCreateGeneral(((PetscObject)pc)->comm, ddn_local, ddidx_local, PETSC_OWN_POINTER, &(osm->gis_local));CHKERRQ(ierr);
359f746d493SDmitry Karpeev       ierr = ISLocalToGlobalMappingCreateIS(osm->gis,&ltog);CHKERRQ(ierr);
36081a5c4aaSDmitry Karpeev       ierr = ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,ddn_local,ddidx_local,&ddn_llocal,ddidx_llocal);CHKERRQ(ierr);
3616bf464f9SBarry Smith       ierr = ISLocalToGlobalMappingDestroy(&ltog);CHKERRQ(ierr);
36293c1ec32SDmitry 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);
363f746d493SDmitry Karpeev       /* Now convert these localized indices into the global indices into the merged output vector. */
3646d98aee9SDmitry Karpeev       ierr = VecGetOwnershipRange(osm->gy, &gfirst, &glast);CHKERRQ(ierr);
365f746d493SDmitry Karpeev       for(j=0; j < ddn_llocal; ++j) {
3666d98aee9SDmitry Karpeev 	ddidx_llocal[j] += gfirst;
367f746d493SDmitry Karpeev       }
368f746d493SDmitry Karpeev       ierr = ISCreateGeneral(PETSC_COMM_SELF,ddn_llocal,ddidx_llocal,PETSC_OWN_POINTER,&gis_llocal);CHKERRQ(ierr);
369f746d493SDmitry Karpeev       ierr = VecScatterCreate(y,osm->gis_local,osm->gy,gis_llocal,&osm->gprolongation);CHKERRQ(ierr);
370fcfd50ebSBarry Smith       ierr = ISDestroy(&gis_llocal);CHKERRQ(ierr);
371b1a0a8a3SJed Brown     }
372f746d493SDmitry Karpeev     /* Create the subdomain work vectors. */
373f746d493SDmitry Karpeev     ierr = PetscMalloc(osm->n*sizeof(Vec),&osm->x);CHKERRQ(ierr);
374f746d493SDmitry Karpeev     ierr = PetscMalloc(osm->n*sizeof(Vec),&osm->y);CHKERRQ(ierr);
3756d98aee9SDmitry Karpeev     ierr = VecGetOwnershipRange(osm->gx, &gfirst, &glast);CHKERRQ(ierr);
376f746d493SDmitry Karpeev     ierr = VecGetArray(osm->gx, &gxarray);CHKERRQ(ierr);
377f746d493SDmitry Karpeev     ierr = VecGetArray(osm->gy, &gyarray);CHKERRQ(ierr);
378f746d493SDmitry Karpeev     for (i=0, ddn=0; i<osm->n; ++i, ddn += dn) {
37986b96d36SDmitry Karpeev       PetscInt dN;
380f746d493SDmitry Karpeev       ierr = ISGetLocalSize(osm->is[i],&dn);CHKERRQ(ierr);
38186b96d36SDmitry Karpeev       ierr = ISGetSize(osm->is[i],&dN);CHKERRQ(ierr);
38286b96d36SDmitry Karpeev       ierr = VecCreateMPIWithArray(((PetscObject)(osm->is[i]))->comm,dn,dN,gxarray+ddn,&osm->x[i]);CHKERRQ(ierr);
38386b96d36SDmitry Karpeev       ierr = VecCreateMPIWithArray(((PetscObject)(osm->is[i]))->comm,dn,dN,gyarray+ddn,&osm->y[i]);CHKERRQ(ierr);
384b1a0a8a3SJed Brown     }
385f746d493SDmitry Karpeev     ierr = VecRestoreArray(osm->gx, &gxarray);CHKERRQ(ierr);
386f746d493SDmitry Karpeev     ierr = VecRestoreArray(osm->gy, &gyarray);CHKERRQ(ierr);
3876bf464f9SBarry Smith     ierr = VecDestroy(&x);CHKERRQ(ierr);
3886bf464f9SBarry Smith     ierr = VecDestroy(&y);CHKERRQ(ierr);
389b1a0a8a3SJed Brown     /* Create the local solvers */
390f746d493SDmitry Karpeev     ierr = PetscMalloc(osm->n*sizeof(KSP *),&osm->ksp);CHKERRQ(ierr);
39144fe18e5SDmitry Karpeev     for (i=0; i<osm->n; i++) { /* KSPs are local */
39286b96d36SDmitry Karpeev       ierr = KSPCreate(((PetscObject)(osm->is[i]))->comm,&ksp);CHKERRQ(ierr);
393b1a0a8a3SJed Brown       ierr = PetscLogObjectParent(pc,ksp);CHKERRQ(ierr);
394b1a0a8a3SJed Brown       ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);CHKERRQ(ierr);
395b1a0a8a3SJed Brown       ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr);
396b1a0a8a3SJed Brown       ierr = KSPGetPC(ksp,&subpc);CHKERRQ(ierr);
397b1a0a8a3SJed Brown       ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
398b1a0a8a3SJed Brown       ierr = KSPSetOptionsPrefix(ksp,prefix);CHKERRQ(ierr);
399b1a0a8a3SJed Brown       ierr = KSPAppendOptionsPrefix(ksp,"sub_");CHKERRQ(ierr);
400b1a0a8a3SJed Brown       osm->ksp[i] = ksp;
401b1a0a8a3SJed Brown     }
402b1a0a8a3SJed Brown     scall = MAT_INITIAL_MATRIX;
403b1a0a8a3SJed Brown 
404b1a0a8a3SJed Brown   } else {
405b1a0a8a3SJed Brown     /*
406b1a0a8a3SJed Brown        Destroy the blocks from the previous iteration
407b1a0a8a3SJed Brown     */
408b1a0a8a3SJed Brown     if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
409f746d493SDmitry Karpeev       ierr = MatDestroyMatrices(osm->n,&osm->pmat);CHKERRQ(ierr);
410b1a0a8a3SJed Brown       scall = MAT_INITIAL_MATRIX;
411b1a0a8a3SJed Brown     }
412b1a0a8a3SJed Brown   }
413b1a0a8a3SJed Brown 
414b1a0a8a3SJed Brown   /*
415f746d493SDmitry Karpeev      Extract out the submatrices.
416b1a0a8a3SJed Brown   */
41781a5c4aaSDmitry Karpeev   if(size > 1) {
4186d42056aSDmitry Karpeev     ierr = MatGetSubMatricesParallel(pc->pmat,osm->n,osm->is, osm->is,scall,&osm->pmat);CHKERRQ(ierr);
41981a5c4aaSDmitry Karpeev   }
42081a5c4aaSDmitry Karpeev   else {
42181a5c4aaSDmitry Karpeev     ierr = MatGetSubMatrices(pc->pmat,osm->n,osm->is, osm->is,scall,&osm->pmat);CHKERRQ(ierr);
42281a5c4aaSDmitry Karpeev   }
423b1a0a8a3SJed Brown   if (scall == MAT_INITIAL_MATRIX) {
424b1a0a8a3SJed Brown     ierr = PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix);CHKERRQ(ierr);
425f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
426b1a0a8a3SJed Brown       ierr = PetscLogObjectParent(pc,osm->pmat[i]);CHKERRQ(ierr);
427b1a0a8a3SJed Brown       ierr = PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i],pprefix);CHKERRQ(ierr);
428b1a0a8a3SJed Brown     }
429b1a0a8a3SJed Brown   }
430b1a0a8a3SJed Brown 
431b1a0a8a3SJed Brown   /* Return control to the user so that the submatrices can be modified (e.g., to apply
432b1a0a8a3SJed Brown      different boundary conditions for the submatrices than for the global problem) */
433c10bc1a1SDmitry Karpeev   ierr = PCModifySubMatrices(pc,osm->n,osm->is,osm->is,osm->pmat,pc->modifysubmatricesP);CHKERRQ(ierr);
434b1a0a8a3SJed Brown 
435b1a0a8a3SJed Brown   /*
436f746d493SDmitry Karpeev      Loop over submatrices putting them into local ksp
437b1a0a8a3SJed Brown   */
438f746d493SDmitry Karpeev   for (i=0; i<osm->n; i++) {
439b1a0a8a3SJed Brown     ierr = KSPSetOperators(osm->ksp[i],osm->pmat[i],osm->pmat[i],pc->flag);CHKERRQ(ierr);
440b1a0a8a3SJed Brown     if (!pc->setupcalled) {
441b1a0a8a3SJed Brown       ierr = KSPSetFromOptions(osm->ksp[i]);CHKERRQ(ierr);
442b1a0a8a3SJed Brown     }
443b1a0a8a3SJed Brown   }
444b1a0a8a3SJed Brown 
445b1a0a8a3SJed Brown   PetscFunctionReturn(0);
446b1a0a8a3SJed Brown }
447b1a0a8a3SJed Brown 
448b1a0a8a3SJed Brown #undef __FUNCT__
449f746d493SDmitry Karpeev #define __FUNCT__ "PCSetUpOnBlocks_GASM"
450f746d493SDmitry Karpeev static PetscErrorCode PCSetUpOnBlocks_GASM(PC pc)
451b1a0a8a3SJed Brown {
452f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
453b1a0a8a3SJed Brown   PetscErrorCode ierr;
454b1a0a8a3SJed Brown   PetscInt       i;
455b1a0a8a3SJed Brown 
456b1a0a8a3SJed Brown   PetscFunctionBegin;
457f746d493SDmitry Karpeev   for (i=0; i<osm->n; i++) {
458b1a0a8a3SJed Brown     ierr = KSPSetUp(osm->ksp[i]);CHKERRQ(ierr);
459b1a0a8a3SJed Brown   }
460b1a0a8a3SJed Brown   PetscFunctionReturn(0);
461b1a0a8a3SJed Brown }
462b1a0a8a3SJed Brown 
463b1a0a8a3SJed Brown #undef __FUNCT__
464f746d493SDmitry Karpeev #define __FUNCT__ "PCApply_GASM"
465f746d493SDmitry Karpeev static PetscErrorCode PCApply_GASM(PC pc,Vec x,Vec y)
466b1a0a8a3SJed Brown {
467f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
468b1a0a8a3SJed Brown   PetscErrorCode ierr;
469f746d493SDmitry Karpeev   PetscInt       i;
470b1a0a8a3SJed Brown   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;
471b1a0a8a3SJed Brown 
472b1a0a8a3SJed Brown   PetscFunctionBegin;
473b1a0a8a3SJed Brown   /*
474b1a0a8a3SJed Brown      Support for limiting the restriction or interpolation to only local
475b1a0a8a3SJed Brown      subdomain values (leaving the other values 0).
476b1a0a8a3SJed Brown   */
477f746d493SDmitry Karpeev   if (!(osm->type & PC_GASM_RESTRICT)) {
478b1a0a8a3SJed Brown     forward = SCATTER_FORWARD_LOCAL;
479b1a0a8a3SJed Brown     /* have to zero the work RHS since scatter may leave some slots empty */
480f746d493SDmitry Karpeev     ierr = VecZeroEntries(osm->gx);CHKERRQ(ierr);
481b1a0a8a3SJed Brown   }
482f746d493SDmitry Karpeev   if (!(osm->type & PC_GASM_INTERPOLATE)) {
483b1a0a8a3SJed Brown     reverse = SCATTER_REVERSE_LOCAL;
484b1a0a8a3SJed Brown   }
485b1a0a8a3SJed Brown 
486f746d493SDmitry Karpeev   ierr = VecScatterBegin(osm->grestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
487b1a0a8a3SJed Brown   ierr = VecZeroEntries(y);CHKERRQ(ierr);
488f746d493SDmitry Karpeev   ierr = VecScatterEnd(osm->grestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
48986b96d36SDmitry Karpeev   /* do the subdomain solves */
49086b96d36SDmitry Karpeev   for (i=0; i<osm->n; ++i) {
491b1a0a8a3SJed Brown     ierr = KSPSolve(osm->ksp[i],osm->x[i],osm->y[i]);CHKERRQ(ierr);
492b1a0a8a3SJed Brown   }
493f746d493SDmitry Karpeev   ierr = VecScatterBegin(osm->gprolongation,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
494f746d493SDmitry Karpeev   ierr = VecScatterEnd(osm->gprolongation,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
495b1a0a8a3SJed Brown   PetscFunctionReturn(0);
496b1a0a8a3SJed Brown }
497b1a0a8a3SJed Brown 
498b1a0a8a3SJed Brown #undef __FUNCT__
499f746d493SDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_GASM"
500f746d493SDmitry Karpeev static PetscErrorCode PCApplyTranspose_GASM(PC pc,Vec x,Vec y)
501b1a0a8a3SJed Brown {
502f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
503b1a0a8a3SJed Brown   PetscErrorCode ierr;
504f746d493SDmitry Karpeev   PetscInt       i;
505b1a0a8a3SJed Brown   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;
506b1a0a8a3SJed Brown 
507b1a0a8a3SJed Brown   PetscFunctionBegin;
508b1a0a8a3SJed Brown   /*
509b1a0a8a3SJed Brown      Support for limiting the restriction or interpolation to only local
510b1a0a8a3SJed Brown      subdomain values (leaving the other values 0).
511b1a0a8a3SJed Brown 
512f746d493SDmitry Karpeev      Note: these are reversed from the PCApply_GASM() because we are applying the
513b1a0a8a3SJed Brown      transpose of the three terms
514b1a0a8a3SJed Brown   */
515f746d493SDmitry Karpeev   if (!(osm->type & PC_GASM_INTERPOLATE)) {
516b1a0a8a3SJed Brown     forward = SCATTER_FORWARD_LOCAL;
517b1a0a8a3SJed Brown     /* have to zero the work RHS since scatter may leave some slots empty */
518f746d493SDmitry Karpeev     ierr = VecZeroEntries(osm->gx);CHKERRQ(ierr);
519b1a0a8a3SJed Brown   }
520f746d493SDmitry Karpeev   if (!(osm->type & PC_GASM_RESTRICT)) {
521b1a0a8a3SJed Brown     reverse = SCATTER_REVERSE_LOCAL;
522b1a0a8a3SJed Brown   }
523b1a0a8a3SJed Brown 
524ab3e923fSDmitry Karpeev   ierr = VecScatterBegin(osm->grestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
525b1a0a8a3SJed Brown   ierr = VecZeroEntries(y);CHKERRQ(ierr);
526ab3e923fSDmitry Karpeev   ierr = VecScatterEnd(osm->grestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
527b1a0a8a3SJed Brown   /* do the local solves */
528ab3e923fSDmitry 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. */
529b1a0a8a3SJed Brown     ierr = KSPSolveTranspose(osm->ksp[i],osm->x[i],osm->y[i]);CHKERRQ(ierr);
530b1a0a8a3SJed Brown   }
531ab3e923fSDmitry Karpeev   ierr = VecScatterBegin(osm->gprolongation,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
532ab3e923fSDmitry Karpeev   ierr = VecScatterEnd(osm->gprolongation,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
533b1a0a8a3SJed Brown   PetscFunctionReturn(0);
534b1a0a8a3SJed Brown }
535b1a0a8a3SJed Brown 
536b1a0a8a3SJed Brown #undef __FUNCT__
537a06653b4SBarry Smith #define __FUNCT__ "PCReset_GASM"
538a06653b4SBarry Smith static PetscErrorCode PCReset_GASM(PC pc)
539b1a0a8a3SJed Brown {
540f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
541b1a0a8a3SJed Brown   PetscErrorCode ierr;
542b1a0a8a3SJed Brown   PetscInt       i;
543b1a0a8a3SJed Brown 
544b1a0a8a3SJed Brown   PetscFunctionBegin;
545b1a0a8a3SJed Brown   if (osm->ksp) {
546f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
547a06653b4SBarry Smith       ierr = KSPReset(osm->ksp[i]);CHKERRQ(ierr);
548b1a0a8a3SJed Brown     }
549b1a0a8a3SJed Brown   }
550b1a0a8a3SJed Brown   if (osm->pmat) {
551f746d493SDmitry Karpeev     if (osm->n > 0) {
552ab3e923fSDmitry Karpeev       ierr = MatDestroyMatrices(osm->n,&osm->pmat);CHKERRQ(ierr);
553b1a0a8a3SJed Brown     }
554b1a0a8a3SJed Brown   }
555d34fcf5fSBarry Smith   if (osm->x) {
556f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
557fcfd50ebSBarry Smith       ierr = VecDestroy(&osm->x[i]);CHKERRQ(ierr);
558fcfd50ebSBarry Smith       ierr = VecDestroy(&osm->y[i]);CHKERRQ(ierr);
559b1a0a8a3SJed Brown     }
560d34fcf5fSBarry Smith   }
5616bf464f9SBarry Smith   ierr = VecDestroy(&osm->gx);CHKERRQ(ierr);
5626bf464f9SBarry Smith   ierr = VecDestroy(&osm->gy);CHKERRQ(ierr);
563ab3e923fSDmitry Karpeev 
5646bf464f9SBarry Smith   ierr = VecScatterDestroy(&osm->grestriction);CHKERRQ(ierr);
5656bf464f9SBarry Smith   ierr = VecScatterDestroy(&osm->gprolongation);CHKERRQ(ierr);
566*b6b0974bSDmitry Karpeev   if (osm->is) {
567*b6b0974bSDmitry Karpeev     ierr = PCGASMDestroySubdomains(osm->n,osm->is,osm->is_local);CHKERRQ(ierr);
568*b6b0974bSDmitry Karpeev     osm->is = 0;
569*b6b0974bSDmitry Karpeev     osm->is_local = 0;
570*b6b0974bSDmitry Karpeev   }
571fcfd50ebSBarry Smith   ierr = ISDestroy(&osm->gis);CHKERRQ(ierr);
572fcfd50ebSBarry Smith   ierr = ISDestroy(&osm->gis_local);CHKERRQ(ierr);
573a06653b4SBarry Smith   PetscFunctionReturn(0);
574a06653b4SBarry Smith }
575a06653b4SBarry Smith 
576a06653b4SBarry Smith #undef __FUNCT__
577a06653b4SBarry Smith #define __FUNCT__ "PCDestroy_GASM"
578a06653b4SBarry Smith static PetscErrorCode PCDestroy_GASM(PC pc)
579a06653b4SBarry Smith {
580a06653b4SBarry Smith   PC_GASM         *osm = (PC_GASM*)pc->data;
581a06653b4SBarry Smith   PetscErrorCode ierr;
582a06653b4SBarry Smith   PetscInt       i;
583a06653b4SBarry Smith 
584a06653b4SBarry Smith   PetscFunctionBegin;
585a06653b4SBarry Smith   ierr = PCReset_GASM(pc);CHKERRQ(ierr);
586a06653b4SBarry Smith   if (osm->ksp) {
587a06653b4SBarry Smith     for (i=0; i<osm->n; i++) {
5886bf464f9SBarry Smith       ierr = KSPDestroy(&osm->ksp[i]);CHKERRQ(ierr);
589a06653b4SBarry Smith     }
590a06653b4SBarry Smith     ierr = PetscFree(osm->ksp);CHKERRQ(ierr);
591a06653b4SBarry Smith   }
592a06653b4SBarry Smith   ierr = PetscFree(osm->x);CHKERRQ(ierr);
593a06653b4SBarry Smith   ierr = PetscFree(osm->y);CHKERRQ(ierr);
594c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
595b1a0a8a3SJed Brown   PetscFunctionReturn(0);
596b1a0a8a3SJed Brown }
597b1a0a8a3SJed Brown 
598b1a0a8a3SJed Brown #undef __FUNCT__
599f746d493SDmitry Karpeev #define __FUNCT__ "PCSetFromOptions_GASM"
600ab3e923fSDmitry Karpeev static PetscErrorCode PCSetFromOptions_GASM(PC pc) {
601f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
602b1a0a8a3SJed Brown   PetscErrorCode ierr;
603b1a0a8a3SJed Brown   PetscInt       blocks,ovl;
604b1a0a8a3SJed Brown   PetscBool      symset,flg;
605f746d493SDmitry Karpeev   PCGASMType      gasmtype;
606b1a0a8a3SJed Brown 
607b1a0a8a3SJed Brown   PetscFunctionBegin;
608b1a0a8a3SJed Brown   /* set the type to symmetric if matrix is symmetric */
609b1a0a8a3SJed Brown   if (!osm->type_set && pc->pmat) {
610b1a0a8a3SJed Brown     ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr);
611f746d493SDmitry Karpeev     if (symset && flg) { osm->type = PC_GASM_BASIC; }
612b1a0a8a3SJed Brown   }
61344fe18e5SDmitry Karpeev   ierr = PetscOptionsHead("Generalized additive Schwarz options");CHKERRQ(ierr);
614f746d493SDmitry Karpeev     ierr = PetscOptionsInt("-pc_gasm_blocks","Number of subdomains","PCGASMSetTotalSubdomains",osm->n,&blocks,&flg);CHKERRQ(ierr);
615ab3e923fSDmitry Karpeev     if (flg) {ierr = PCGASMSetTotalSubdomains(pc,blocks);CHKERRQ(ierr); }
616f746d493SDmitry Karpeev     ierr = PetscOptionsInt("-pc_gasm_overlap","Number of grid points overlap","PCGASMSetOverlap",osm->overlap,&ovl,&flg);CHKERRQ(ierr);
617f746d493SDmitry Karpeev     if (flg) {ierr = PCGASMSetOverlap(pc,ovl);CHKERRQ(ierr); }
618b1a0a8a3SJed Brown     flg  = PETSC_FALSE;
619f746d493SDmitry Karpeev     ierr = PetscOptionsEnum("-pc_gasm_type","Type of restriction/extension","PCGASMSetType",PCGASMTypes,(PetscEnum)osm->type,(PetscEnum*)&gasmtype,&flg);CHKERRQ(ierr);
620f746d493SDmitry Karpeev     if (flg) {ierr = PCGASMSetType(pc,gasmtype);CHKERRQ(ierr); }
621b1a0a8a3SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
622b1a0a8a3SJed Brown   PetscFunctionReturn(0);
623b1a0a8a3SJed Brown }
624b1a0a8a3SJed Brown 
625b1a0a8a3SJed Brown /*------------------------------------------------------------------------------------*/
626b1a0a8a3SJed Brown 
627b1a0a8a3SJed Brown EXTERN_C_BEGIN
628b1a0a8a3SJed Brown #undef __FUNCT__
629f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetLocalSubdomains_GASM"
6307087cfbeSBarry Smith PetscErrorCode  PCGASMSetLocalSubdomains_GASM(PC pc,PetscInt n,IS is[],IS is_local[])
631b1a0a8a3SJed Brown {
632f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
633b1a0a8a3SJed Brown   PetscErrorCode ierr;
634b1a0a8a3SJed Brown   PetscInt       i;
635b1a0a8a3SJed Brown 
636b1a0a8a3SJed Brown   PetscFunctionBegin;
637b1a0a8a3SJed Brown   if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Each process must have 1 or more blocks, n = %D",n);
638f746d493SDmitry Karpeev   if (pc->setupcalled && (n != osm->n || is)) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetLocalSubdomains() should be called before calling PCSetUp().");
639b1a0a8a3SJed Brown 
640b1a0a8a3SJed Brown   if (!pc->setupcalled) {
64193c1ec32SDmitry Karpeev     osm->n            = n;
64293c1ec32SDmitry Karpeev     osm->is           = 0;
64393c1ec32SDmitry Karpeev     osm->is_local     = 0;
644b1a0a8a3SJed Brown     if (is) {
645b1a0a8a3SJed Brown       for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)is[i]);CHKERRQ(ierr);}
646b1a0a8a3SJed Brown     }
647b1a0a8a3SJed Brown     if (is_local) {
648b1a0a8a3SJed Brown       for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)is_local[i]);CHKERRQ(ierr);}
649b1a0a8a3SJed Brown     }
650b1a0a8a3SJed Brown     if (osm->is) {
651ab3e923fSDmitry Karpeev       ierr = PCGASMDestroySubdomains(osm->n,osm->is,osm->is_local);CHKERRQ(ierr);
652b1a0a8a3SJed Brown     }
653b1a0a8a3SJed Brown     if (is) {
654b1a0a8a3SJed Brown       ierr = PetscMalloc(n*sizeof(IS),&osm->is);CHKERRQ(ierr);
655b1a0a8a3SJed Brown       for (i=0; i<n; i++) { osm->is[i] = is[i]; }
656f746d493SDmitry Karpeev       /* Flag indicating that the user has set overlapping subdomains so PCGASM should not increase their size. */
657b1a0a8a3SJed Brown       osm->overlap = -1;
658b1a0a8a3SJed Brown     }
659b1a0a8a3SJed Brown     if (is_local) {
660b1a0a8a3SJed Brown       ierr = PetscMalloc(n*sizeof(IS),&osm->is_local);CHKERRQ(ierr);
661b1a0a8a3SJed Brown       for (i=0; i<n; i++) { osm->is_local[i] = is_local[i]; }
662b1a0a8a3SJed Brown     }
663b1a0a8a3SJed Brown   }
664b1a0a8a3SJed Brown   PetscFunctionReturn(0);
665b1a0a8a3SJed Brown }
666b1a0a8a3SJed Brown EXTERN_C_END
667b1a0a8a3SJed Brown 
668b1a0a8a3SJed Brown EXTERN_C_BEGIN
669b1a0a8a3SJed Brown #undef __FUNCT__
670f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetTotalSubdomains_GASM"
6717087cfbeSBarry Smith PetscErrorCode  PCGASMSetTotalSubdomains_GASM(PC pc,PetscInt N) {
672f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
673b1a0a8a3SJed Brown   PetscErrorCode ierr;
674b1a0a8a3SJed Brown   PetscMPIInt    rank,size;
675b1a0a8a3SJed Brown   PetscInt       n;
676b1a0a8a3SJed Brown 
677b1a0a8a3SJed Brown   PetscFunctionBegin;
678b1a0a8a3SJed Brown   if (N < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of total blocks must be > 0, N = %D",N);
679b1a0a8a3SJed Brown 
680b1a0a8a3SJed Brown   /*
681b1a0a8a3SJed Brown      Split the subdomains equally among all processors
682b1a0a8a3SJed Brown   */
683b1a0a8a3SJed Brown   ierr = MPI_Comm_rank(((PetscObject)pc)->comm,&rank);CHKERRQ(ierr);
684b1a0a8a3SJed Brown   ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr);
685b1a0a8a3SJed Brown   n = N/size + ((N % size) > rank);
686b1a0a8a3SJed 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);
687f746d493SDmitry Karpeev   if (pc->setupcalled && n != osm->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetTotalSubdomains() should be called before PCSetUp().");
688b1a0a8a3SJed Brown   if (!pc->setupcalled) {
689b1a0a8a3SJed Brown     if (osm->is) {
690ab3e923fSDmitry Karpeev       ierr = PCGASMDestroySubdomains(osm->n,osm->is,osm->is_local);CHKERRQ(ierr);
691b1a0a8a3SJed Brown     }
69293c1ec32SDmitry Karpeev     osm->N            = N;
693f746d493SDmitry Karpeev     osm->n            = n;
694b1a0a8a3SJed Brown     osm->is           = 0;
695b1a0a8a3SJed Brown     osm->is_local     = 0;
696b1a0a8a3SJed Brown   }
697b1a0a8a3SJed Brown   PetscFunctionReturn(0);
698ab3e923fSDmitry Karpeev }/* PCGASMSetTotalSubdomains_GASM() */
699b1a0a8a3SJed Brown EXTERN_C_END
700b1a0a8a3SJed Brown 
701b1a0a8a3SJed Brown EXTERN_C_BEGIN
702b1a0a8a3SJed Brown #undef __FUNCT__
703f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetOverlap_GASM"
7047087cfbeSBarry Smith PetscErrorCode  PCGASMSetOverlap_GASM(PC pc,PetscInt ovl)
705b1a0a8a3SJed Brown {
706f746d493SDmitry Karpeev   PC_GASM *osm = (PC_GASM*)pc->data;
707b1a0a8a3SJed Brown 
708b1a0a8a3SJed Brown   PetscFunctionBegin;
709b1a0a8a3SJed Brown   if (ovl < 0) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested");
710f746d493SDmitry Karpeev   if (pc->setupcalled && ovl != osm->overlap) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetOverlap() should be called before PCSetUp().");
711b1a0a8a3SJed Brown   if (!pc->setupcalled) {
712b1a0a8a3SJed Brown     osm->overlap = ovl;
713b1a0a8a3SJed Brown   }
714b1a0a8a3SJed Brown   PetscFunctionReturn(0);
715b1a0a8a3SJed Brown }
716b1a0a8a3SJed Brown EXTERN_C_END
717b1a0a8a3SJed Brown 
718b1a0a8a3SJed Brown EXTERN_C_BEGIN
719b1a0a8a3SJed Brown #undef __FUNCT__
720f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetType_GASM"
7217087cfbeSBarry Smith PetscErrorCode  PCGASMSetType_GASM(PC pc,PCGASMType type)
722b1a0a8a3SJed Brown {
723f746d493SDmitry Karpeev   PC_GASM *osm = (PC_GASM*)pc->data;
724b1a0a8a3SJed Brown 
725b1a0a8a3SJed Brown   PetscFunctionBegin;
726b1a0a8a3SJed Brown   osm->type     = type;
727b1a0a8a3SJed Brown   osm->type_set = PETSC_TRUE;
728b1a0a8a3SJed Brown   PetscFunctionReturn(0);
729b1a0a8a3SJed Brown }
730b1a0a8a3SJed Brown EXTERN_C_END
731b1a0a8a3SJed Brown 
732b1a0a8a3SJed Brown EXTERN_C_BEGIN
733b1a0a8a3SJed Brown #undef __FUNCT__
734f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetSortIndices_GASM"
7357087cfbeSBarry Smith PetscErrorCode  PCGASMSetSortIndices_GASM(PC pc,PetscBool  doSort)
736b1a0a8a3SJed Brown {
737f746d493SDmitry Karpeev   PC_GASM *osm = (PC_GASM*)pc->data;
738b1a0a8a3SJed Brown 
739b1a0a8a3SJed Brown   PetscFunctionBegin;
740b1a0a8a3SJed Brown   osm->sort_indices = doSort;
741b1a0a8a3SJed Brown   PetscFunctionReturn(0);
742b1a0a8a3SJed Brown }
743b1a0a8a3SJed Brown EXTERN_C_END
744b1a0a8a3SJed Brown 
745b1a0a8a3SJed Brown EXTERN_C_BEGIN
746b1a0a8a3SJed Brown #undef __FUNCT__
747f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMGetSubKSP_GASM"
74844fe18e5SDmitry Karpeev /*
74944fe18e5SDmitry Karpeev    FIX: This routine might need to be modified once multiple ranks per subdomain are allowed.
75044fe18e5SDmitry Karpeev         In particular, it would upset the global subdomain number calculation.
75144fe18e5SDmitry Karpeev */
7527087cfbeSBarry Smith PetscErrorCode  PCGASMGetSubKSP_GASM(PC pc,PetscInt *n,PetscInt *first,KSP **ksp)
753b1a0a8a3SJed Brown {
754f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
755b1a0a8a3SJed Brown   PetscErrorCode ierr;
756b1a0a8a3SJed Brown 
757b1a0a8a3SJed Brown   PetscFunctionBegin;
758ab3e923fSDmitry 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");
759b1a0a8a3SJed Brown 
760ab3e923fSDmitry Karpeev   if (n) {
761ab3e923fSDmitry Karpeev     *n = osm->n;
762b1a0a8a3SJed Brown   }
763ab3e923fSDmitry Karpeev   if (first) {
764ab3e923fSDmitry Karpeev     ierr = MPI_Scan(&osm->n,first,1,MPIU_INT,MPI_SUM,((PetscObject)pc)->comm);CHKERRQ(ierr);
765ab3e923fSDmitry Karpeev     *first -= osm->n;
766b1a0a8a3SJed Brown   }
767b1a0a8a3SJed Brown   if (ksp) {
768b1a0a8a3SJed Brown     /* Assume that local solves are now different; not necessarily
769f746d493SDmitry Karpeev        true though!  This flag is used only for PCView_GASM() */
770b1a0a8a3SJed Brown     *ksp                   = osm->ksp;
771b1a0a8a3SJed Brown     osm->same_local_solves = PETSC_FALSE;
772b1a0a8a3SJed Brown   }
773b1a0a8a3SJed Brown   PetscFunctionReturn(0);
774ab3e923fSDmitry Karpeev }/* PCGASMGetSubKSP_GASM() */
775b1a0a8a3SJed Brown EXTERN_C_END
776b1a0a8a3SJed Brown 
777b1a0a8a3SJed Brown 
778b1a0a8a3SJed Brown #undef __FUNCT__
779f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetLocalSubdomains"
780b1a0a8a3SJed Brown /*@C
781f746d493SDmitry Karpeev     PCGASMSetLocalSubdomains - Sets the local subdomains (for this processor
782b1a0a8a3SJed Brown     only) for the additive Schwarz preconditioner.
783b1a0a8a3SJed Brown 
784b1a0a8a3SJed Brown     Collective on PC
785b1a0a8a3SJed Brown 
786b1a0a8a3SJed Brown     Input Parameters:
787b1a0a8a3SJed Brown +   pc - the preconditioner context
788b1a0a8a3SJed Brown .   n - the number of subdomains for this processor (default value = 1)
789b1a0a8a3SJed Brown .   is - the index set that defines the subdomains for this processor
790b1a0a8a3SJed Brown          (or PETSC_NULL for PETSc to determine subdomains)
791b1a0a8a3SJed Brown -   is_local - the index sets that define the local part of the subdomains for this processor
792b1a0a8a3SJed Brown          (or PETSC_NULL to use the default of 1 subdomain per process)
793b1a0a8a3SJed Brown 
794b1a0a8a3SJed Brown     Notes:
795b1a0a8a3SJed Brown     The IS numbering is in the parallel, global numbering of the vector.
796b1a0a8a3SJed Brown 
797f746d493SDmitry Karpeev     By default the GASM preconditioner uses 1 block per processor.
798b1a0a8a3SJed Brown 
799f746d493SDmitry Karpeev     Use PCGASMSetTotalSubdomains() to set the subdomains for all processors.
800b1a0a8a3SJed Brown 
801b1a0a8a3SJed Brown     Level: advanced
802b1a0a8a3SJed Brown 
803f746d493SDmitry Karpeev .keywords: PC, GASM, set, local, subdomains, additive Schwarz
804b1a0a8a3SJed Brown 
805f746d493SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(),
806f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D(), PCGASMGetLocalSubdomains()
807b1a0a8a3SJed Brown @*/
8087087cfbeSBarry Smith PetscErrorCode  PCGASMSetLocalSubdomains(PC pc,PetscInt n,IS is[],IS is_local[])
809b1a0a8a3SJed Brown {
810b1a0a8a3SJed Brown   PetscErrorCode ierr;
811b1a0a8a3SJed Brown 
812b1a0a8a3SJed Brown   PetscFunctionBegin;
813b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
814f746d493SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetLocalSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,n,is,is_local));CHKERRQ(ierr);
815b1a0a8a3SJed Brown   PetscFunctionReturn(0);
816b1a0a8a3SJed Brown }
817b1a0a8a3SJed Brown 
818b1a0a8a3SJed Brown #undef __FUNCT__
819f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetTotalSubdomains"
820b1a0a8a3SJed Brown /*@C
821f746d493SDmitry Karpeev     PCGASMSetTotalSubdomains - Sets the subdomains for all processor for the
822b1a0a8a3SJed Brown     additive Schwarz preconditioner.  Either all or no processors in the
823b1a0a8a3SJed Brown     PC communicator must call this routine, with the same index sets.
824b1a0a8a3SJed Brown 
825b1a0a8a3SJed Brown     Collective on PC
826b1a0a8a3SJed Brown 
827b1a0a8a3SJed Brown     Input Parameters:
828b1a0a8a3SJed Brown +   pc - the preconditioner context
829b1a0a8a3SJed Brown .   n - the number of subdomains for all processors
830b1a0a8a3SJed Brown .   is - the index sets that define the subdomains for all processor
831b1a0a8a3SJed Brown          (or PETSC_NULL for PETSc to determine subdomains)
832b1a0a8a3SJed Brown -   is_local - the index sets that define the local part of the subdomains for this processor
833b1a0a8a3SJed Brown          (or PETSC_NULL to use the default of 1 subdomain per process)
834b1a0a8a3SJed Brown 
835b1a0a8a3SJed Brown     Options Database Key:
836b1a0a8a3SJed Brown     To set the total number of subdomain blocks rather than specify the
837b1a0a8a3SJed Brown     index sets, use the option
838f746d493SDmitry Karpeev .    -pc_gasm_blocks <blks> - Sets total blocks
839b1a0a8a3SJed Brown 
840b1a0a8a3SJed Brown     Notes:
841b1a0a8a3SJed Brown     Currently you cannot use this to set the actual subdomains with the argument is.
842b1a0a8a3SJed Brown 
843f746d493SDmitry Karpeev     By default the GASM preconditioner uses 1 block per processor.
844b1a0a8a3SJed Brown 
845b1a0a8a3SJed Brown     These index sets cannot be destroyed until after completion of the
846f746d493SDmitry Karpeev     linear solves for which the GASM preconditioner is being used.
847b1a0a8a3SJed Brown 
848f746d493SDmitry Karpeev     Use PCGASMSetLocalSubdomains() to set local subdomains.
849b1a0a8a3SJed Brown 
850b1a0a8a3SJed Brown     Level: advanced
851b1a0a8a3SJed Brown 
852f746d493SDmitry Karpeev .keywords: PC, GASM, set, total, global, subdomains, additive Schwarz
853b1a0a8a3SJed Brown 
854f746d493SDmitry Karpeev .seealso: PCGASMSetLocalSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(),
855f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D()
856b1a0a8a3SJed Brown @*/
8577087cfbeSBarry Smith PetscErrorCode  PCGASMSetTotalSubdomains(PC pc,PetscInt N)
858b1a0a8a3SJed Brown {
859b1a0a8a3SJed Brown   PetscErrorCode ierr;
860b1a0a8a3SJed Brown 
861b1a0a8a3SJed Brown   PetscFunctionBegin;
862b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
863f746d493SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetTotalSubdomains_C",(PC,PetscInt),(pc,N));CHKERRQ(ierr);
864b1a0a8a3SJed Brown   PetscFunctionReturn(0);
865b1a0a8a3SJed Brown }
866b1a0a8a3SJed Brown 
867b1a0a8a3SJed Brown #undef __FUNCT__
868f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetOverlap"
869b1a0a8a3SJed Brown /*@
870f746d493SDmitry Karpeev     PCGASMSetOverlap - Sets the overlap between a pair of subdomains for the
871b1a0a8a3SJed Brown     additive Schwarz preconditioner.  Either all or no processors in the
872b1a0a8a3SJed Brown     PC communicator must call this routine.
873b1a0a8a3SJed Brown 
874b1a0a8a3SJed Brown     Logically Collective on PC
875b1a0a8a3SJed Brown 
876b1a0a8a3SJed Brown     Input Parameters:
877b1a0a8a3SJed Brown +   pc  - the preconditioner context
878b1a0a8a3SJed Brown -   ovl - the amount of overlap between subdomains (ovl >= 0, default value = 1)
879b1a0a8a3SJed Brown 
880b1a0a8a3SJed Brown     Options Database Key:
881f746d493SDmitry Karpeev .   -pc_gasm_overlap <ovl> - Sets overlap
882b1a0a8a3SJed Brown 
883b1a0a8a3SJed Brown     Notes:
884f746d493SDmitry Karpeev     By default the GASM preconditioner uses 1 block per processor.  To use
885f746d493SDmitry Karpeev     multiple blocks per perocessor, see PCGASMSetTotalSubdomains() and
886f746d493SDmitry Karpeev     PCGASMSetLocalSubdomains() (and the option -pc_gasm_blocks <blks>).
887b1a0a8a3SJed Brown 
888b1a0a8a3SJed Brown     The overlap defaults to 1, so if one desires that no additional
889b1a0a8a3SJed Brown     overlap be computed beyond what may have been set with a call to
890f746d493SDmitry Karpeev     PCGASMSetTotalSubdomains() or PCGASMSetLocalSubdomains(), then ovl
891b1a0a8a3SJed Brown     must be set to be 0.  In particular, if one does not explicitly set
892b1a0a8a3SJed Brown     the subdomains an application code, then all overlap would be computed
893f746d493SDmitry Karpeev     internally by PETSc, and using an overlap of 0 would result in an GASM
894b1a0a8a3SJed Brown     variant that is equivalent to the block Jacobi preconditioner.
895b1a0a8a3SJed Brown 
896b1a0a8a3SJed Brown     Note that one can define initial index sets with any overlap via
897f746d493SDmitry Karpeev     PCGASMSetTotalSubdomains() or PCGASMSetLocalSubdomains(); the routine
898f746d493SDmitry Karpeev     PCGASMSetOverlap() merely allows PETSc to extend that overlap further
899b1a0a8a3SJed Brown     if desired.
900b1a0a8a3SJed Brown 
901b1a0a8a3SJed Brown     Level: intermediate
902b1a0a8a3SJed Brown 
903f746d493SDmitry Karpeev .keywords: PC, GASM, set, overlap
904b1a0a8a3SJed Brown 
905f746d493SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetLocalSubdomains(), PCGASMGetSubKSP(),
906f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D(), PCGASMGetLocalSubdomains()
907b1a0a8a3SJed Brown @*/
9087087cfbeSBarry Smith PetscErrorCode  PCGASMSetOverlap(PC pc,PetscInt ovl)
909b1a0a8a3SJed Brown {
910b1a0a8a3SJed Brown   PetscErrorCode ierr;
911b1a0a8a3SJed Brown 
912b1a0a8a3SJed Brown   PetscFunctionBegin;
913b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
914b1a0a8a3SJed Brown   PetscValidLogicalCollectiveInt(pc,ovl,2);
915f746d493SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetOverlap_C",(PC,PetscInt),(pc,ovl));CHKERRQ(ierr);
916b1a0a8a3SJed Brown   PetscFunctionReturn(0);
917b1a0a8a3SJed Brown }
918b1a0a8a3SJed Brown 
919b1a0a8a3SJed Brown #undef __FUNCT__
920f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetType"
921b1a0a8a3SJed Brown /*@
922f746d493SDmitry Karpeev     PCGASMSetType - Sets the type of restriction and interpolation used
923b1a0a8a3SJed Brown     for local problems in the additive Schwarz method.
924b1a0a8a3SJed Brown 
925b1a0a8a3SJed Brown     Logically Collective on PC
926b1a0a8a3SJed Brown 
927b1a0a8a3SJed Brown     Input Parameters:
928b1a0a8a3SJed Brown +   pc  - the preconditioner context
929f746d493SDmitry Karpeev -   type - variant of GASM, one of
930b1a0a8a3SJed Brown .vb
931f746d493SDmitry Karpeev       PC_GASM_BASIC       - full interpolation and restriction
932f746d493SDmitry Karpeev       PC_GASM_RESTRICT    - full restriction, local processor interpolation
933f746d493SDmitry Karpeev       PC_GASM_INTERPOLATE - full interpolation, local processor restriction
934f746d493SDmitry Karpeev       PC_GASM_NONE        - local processor restriction and interpolation
935b1a0a8a3SJed Brown .ve
936b1a0a8a3SJed Brown 
937b1a0a8a3SJed Brown     Options Database Key:
938f746d493SDmitry Karpeev .   -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type
939b1a0a8a3SJed Brown 
940b1a0a8a3SJed Brown     Level: intermediate
941b1a0a8a3SJed Brown 
942f746d493SDmitry Karpeev .keywords: PC, GASM, set, type
943b1a0a8a3SJed Brown 
944f746d493SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetTotalSubdomains(), PCGASMGetSubKSP(),
945f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D()
946b1a0a8a3SJed Brown @*/
9477087cfbeSBarry Smith PetscErrorCode  PCGASMSetType(PC pc,PCGASMType type)
948b1a0a8a3SJed Brown {
949b1a0a8a3SJed Brown   PetscErrorCode ierr;
950b1a0a8a3SJed Brown 
951b1a0a8a3SJed Brown   PetscFunctionBegin;
952b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
953b1a0a8a3SJed Brown   PetscValidLogicalCollectiveEnum(pc,type,2);
954f746d493SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetType_C",(PC,PCGASMType),(pc,type));CHKERRQ(ierr);
955b1a0a8a3SJed Brown   PetscFunctionReturn(0);
956b1a0a8a3SJed Brown }
957b1a0a8a3SJed Brown 
958b1a0a8a3SJed Brown #undef __FUNCT__
959f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetSortIndices"
960b1a0a8a3SJed Brown /*@
961f746d493SDmitry Karpeev     PCGASMSetSortIndices - Determines whether subdomain indices are sorted.
962b1a0a8a3SJed Brown 
963b1a0a8a3SJed Brown     Logically Collective on PC
964b1a0a8a3SJed Brown 
965b1a0a8a3SJed Brown     Input Parameters:
966b1a0a8a3SJed Brown +   pc  - the preconditioner context
967b1a0a8a3SJed Brown -   doSort - sort the subdomain indices
968b1a0a8a3SJed Brown 
969b1a0a8a3SJed Brown     Level: intermediate
970b1a0a8a3SJed Brown 
971f746d493SDmitry Karpeev .keywords: PC, GASM, set, type
972b1a0a8a3SJed Brown 
973f746d493SDmitry Karpeev .seealso: PCGASMSetLocalSubdomains(), PCGASMSetTotalSubdomains(), PCGASMGetSubKSP(),
974f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D()
975b1a0a8a3SJed Brown @*/
9767087cfbeSBarry Smith PetscErrorCode  PCGASMSetSortIndices(PC pc,PetscBool  doSort)
977b1a0a8a3SJed Brown {
978b1a0a8a3SJed Brown   PetscErrorCode ierr;
979b1a0a8a3SJed Brown 
980b1a0a8a3SJed Brown   PetscFunctionBegin;
981b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
982b1a0a8a3SJed Brown   PetscValidLogicalCollectiveBool(pc,doSort,2);
983f746d493SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetSortIndices_C",(PC,PetscBool),(pc,doSort));CHKERRQ(ierr);
984b1a0a8a3SJed Brown   PetscFunctionReturn(0);
985b1a0a8a3SJed Brown }
986b1a0a8a3SJed Brown 
987b1a0a8a3SJed Brown #undef __FUNCT__
988f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMGetSubKSP"
989b1a0a8a3SJed Brown /*@C
990f746d493SDmitry Karpeev    PCGASMGetSubKSP - Gets the local KSP contexts for all blocks on
991b1a0a8a3SJed Brown    this processor.
992b1a0a8a3SJed Brown 
993b1a0a8a3SJed Brown    Collective on PC iff first_local is requested
994b1a0a8a3SJed Brown 
995b1a0a8a3SJed Brown    Input Parameter:
996b1a0a8a3SJed Brown .  pc - the preconditioner context
997b1a0a8a3SJed Brown 
998b1a0a8a3SJed Brown    Output Parameters:
999b1a0a8a3SJed Brown +  n_local - the number of blocks on this processor or PETSC_NULL
1000b1a0a8a3SJed Brown .  first_local - the global number of the first block on this processor or PETSC_NULL,
1001b1a0a8a3SJed Brown                  all processors must request or all must pass PETSC_NULL
1002b1a0a8a3SJed Brown -  ksp - the array of KSP contexts
1003b1a0a8a3SJed Brown 
1004b1a0a8a3SJed Brown    Note:
1005f746d493SDmitry Karpeev    After PCGASMGetSubKSP() the array of KSPes is not to be freed
1006b1a0a8a3SJed Brown 
1007b1a0a8a3SJed Brown    Currently for some matrix implementations only 1 block per processor
1008b1a0a8a3SJed Brown    is supported.
1009b1a0a8a3SJed Brown 
1010f746d493SDmitry Karpeev    You must call KSPSetUp() before calling PCGASMGetSubKSP().
1011b1a0a8a3SJed Brown 
1012b1a0a8a3SJed Brown    Level: advanced
1013b1a0a8a3SJed Brown 
1014f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, get, sub, KSP, context
1015b1a0a8a3SJed Brown 
1016f746d493SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetTotalSubdomains(), PCGASMSetOverlap(),
1017f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D(),
1018b1a0a8a3SJed Brown @*/
10197087cfbeSBarry Smith PetscErrorCode  PCGASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[])
1020b1a0a8a3SJed Brown {
1021b1a0a8a3SJed Brown   PetscErrorCode ierr;
1022b1a0a8a3SJed Brown 
1023b1a0a8a3SJed Brown   PetscFunctionBegin;
1024b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1025f746d493SDmitry Karpeev   ierr = PetscUseMethod(pc,"PCGASMGetSubKSP_C",(PC,PetscInt*,PetscInt*,KSP **),(pc,n_local,first_local,ksp));CHKERRQ(ierr);
1026b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1027b1a0a8a3SJed Brown }
1028b1a0a8a3SJed Brown 
1029b1a0a8a3SJed Brown /* -------------------------------------------------------------------------------------*/
1030b1a0a8a3SJed Brown /*MC
1031f746d493SDmitry Karpeev    PCGASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with
1032b1a0a8a3SJed Brown            its own KSP object.
1033b1a0a8a3SJed Brown 
1034b1a0a8a3SJed Brown    Options Database Keys:
1035f746d493SDmitry Karpeev +  -pc_gasm_truelocal - Activates PCGGASMSetUseTrueLocal()
1036f746d493SDmitry Karpeev .  -pc_gasm_blocks <blks> - Sets total blocks
1037f746d493SDmitry Karpeev .  -pc_gasm_overlap <ovl> - Sets overlap
1038f746d493SDmitry Karpeev -  -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type
1039b1a0a8a3SJed Brown 
1040b1a0a8a3SJed Brown      IMPORTANT: If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you
1041f746d493SDmitry Karpeev       will get a different convergence rate due to the default option of -pc_gasm_type restrict. Use
1042f746d493SDmitry Karpeev       -pc_gasm_type basic to use the standard GASM.
1043b1a0a8a3SJed Brown 
1044b1a0a8a3SJed Brown    Notes: Each processor can have one or more blocks, but a block cannot be shared by more
1045b1a0a8a3SJed Brown      than one processor. Defaults to one block per processor.
1046b1a0a8a3SJed Brown 
1047b1a0a8a3SJed Brown      To set options on the solvers for each block append -sub_ to all the KSP, and PC
1048b1a0a8a3SJed Brown         options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly
1049b1a0a8a3SJed Brown 
1050f746d493SDmitry Karpeev      To set the options on the solvers separate for each block call PCGASMGetSubKSP()
1051b1a0a8a3SJed Brown          and set the options directly on the resulting KSP object (you can access its PC
1052b1a0a8a3SJed Brown          with KSPGetPC())
1053b1a0a8a3SJed Brown 
1054b1a0a8a3SJed Brown 
1055b1a0a8a3SJed Brown    Level: beginner
1056b1a0a8a3SJed Brown 
1057b1a0a8a3SJed Brown    Concepts: additive Schwarz method
1058b1a0a8a3SJed Brown 
1059b1a0a8a3SJed Brown     References:
1060b1a0a8a3SJed Brown     An additive variant of the Schwarz alternating method for the case of many subregions
1061b1a0a8a3SJed Brown     M Dryja, OB Widlund - Courant Institute, New York University Technical report
1062b1a0a8a3SJed Brown 
1063b1a0a8a3SJed Brown     Domain Decompositions: Parallel Multilevel Methods for Elliptic Partial Differential Equations,
1064b1a0a8a3SJed Brown     Barry Smith, Petter Bjorstad, and William Gropp, Cambridge University Press, ISBN 0-521-49589-X.
1065b1a0a8a3SJed Brown 
1066b1a0a8a3SJed Brown .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
1067f746d493SDmitry Karpeev            PCBJACOBI, PCGASMSetUseTrueLocal(), PCGASMGetSubKSP(), PCGASMSetLocalSubdomains(),
1068f746d493SDmitry Karpeev            PCGASMSetTotalSubdomains(), PCSetModifySubmatrices(), PCGASMSetOverlap(), PCGASMSetType()
1069b1a0a8a3SJed Brown 
1070b1a0a8a3SJed Brown M*/
1071b1a0a8a3SJed Brown 
1072b1a0a8a3SJed Brown EXTERN_C_BEGIN
1073b1a0a8a3SJed Brown #undef __FUNCT__
1074f746d493SDmitry Karpeev #define __FUNCT__ "PCCreate_GASM"
10757087cfbeSBarry Smith PetscErrorCode  PCCreate_GASM(PC pc)
1076b1a0a8a3SJed Brown {
1077b1a0a8a3SJed Brown   PetscErrorCode ierr;
1078f746d493SDmitry Karpeev   PC_GASM         *osm;
1079b1a0a8a3SJed Brown 
1080b1a0a8a3SJed Brown   PetscFunctionBegin;
1081f746d493SDmitry Karpeev   ierr = PetscNewLog(pc,PC_GASM,&osm);CHKERRQ(ierr);
1082ab3e923fSDmitry Karpeev   osm->N                 = PETSC_DECIDE;
1083ab3e923fSDmitry Karpeev   osm->n                 = 0;
1084ab3e923fSDmitry Karpeev   osm->nmax              = 0;
1085b1a0a8a3SJed Brown   osm->overlap           = 1;
1086b1a0a8a3SJed Brown   osm->ksp               = 0;
1087ab3e923fSDmitry Karpeev   osm->grestriction      = 0;
1088ab3e923fSDmitry Karpeev   osm->gprolongation     = 0;
1089ab3e923fSDmitry Karpeev   osm->gx                = 0;
1090ab3e923fSDmitry Karpeev   osm->gy                = 0;
1091b1a0a8a3SJed Brown   osm->x                 = 0;
1092b1a0a8a3SJed Brown   osm->y                 = 0;
1093b1a0a8a3SJed Brown   osm->is                = 0;
1094b1a0a8a3SJed Brown   osm->is_local          = 0;
1095b1a0a8a3SJed Brown   osm->pmat              = 0;
1096f746d493SDmitry Karpeev   osm->type              = PC_GASM_RESTRICT;
1097b1a0a8a3SJed Brown   osm->same_local_solves = PETSC_TRUE;
1098b1a0a8a3SJed Brown   osm->sort_indices      = PETSC_TRUE;
1099b1a0a8a3SJed Brown 
1100b1a0a8a3SJed Brown   pc->data                   = (void*)osm;
1101f746d493SDmitry Karpeev   pc->ops->apply             = PCApply_GASM;
1102f746d493SDmitry Karpeev   pc->ops->applytranspose    = PCApplyTranspose_GASM;
1103f746d493SDmitry Karpeev   pc->ops->setup             = PCSetUp_GASM;
1104a06653b4SBarry Smith   pc->ops->reset             = PCReset_GASM;
1105f746d493SDmitry Karpeev   pc->ops->destroy           = PCDestroy_GASM;
1106f746d493SDmitry Karpeev   pc->ops->setfromoptions    = PCSetFromOptions_GASM;
1107f746d493SDmitry Karpeev   pc->ops->setuponblocks     = PCSetUpOnBlocks_GASM;
1108f746d493SDmitry Karpeev   pc->ops->view              = PCView_GASM;
1109b1a0a8a3SJed Brown   pc->ops->applyrichardson   = 0;
1110b1a0a8a3SJed Brown 
1111f746d493SDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetLocalSubdomains_C","PCGASMSetLocalSubdomains_GASM",
1112f746d493SDmitry Karpeev                     PCGASMSetLocalSubdomains_GASM);CHKERRQ(ierr);
1113f746d493SDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetTotalSubdomains_C","PCGASMSetTotalSubdomains_GASM",
1114f746d493SDmitry Karpeev                     PCGASMSetTotalSubdomains_GASM);CHKERRQ(ierr);
1115f746d493SDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetOverlap_C","PCGASMSetOverlap_GASM",
1116f746d493SDmitry Karpeev                     PCGASMSetOverlap_GASM);CHKERRQ(ierr);
1117f746d493SDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetType_C","PCGASMSetType_GASM",
1118f746d493SDmitry Karpeev                     PCGASMSetType_GASM);CHKERRQ(ierr);
1119f746d493SDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetSortIndices_C","PCGASMSetSortIndices_GASM",
1120f746d493SDmitry Karpeev                     PCGASMSetSortIndices_GASM);CHKERRQ(ierr);
1121f746d493SDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMGetSubKSP_C","PCGASMGetSubKSP_GASM",
1122f746d493SDmitry Karpeev                     PCGASMGetSubKSP_GASM);CHKERRQ(ierr);
1123b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1124b1a0a8a3SJed Brown }
1125b1a0a8a3SJed Brown EXTERN_C_END
1126b1a0a8a3SJed Brown 
1127b1a0a8a3SJed Brown 
1128b1a0a8a3SJed Brown #undef __FUNCT__
1129f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMCreateSubdomains"
1130b1a0a8a3SJed Brown /*@C
1131f746d493SDmitry Karpeev    PCGASMCreateSubdomains - Creates the index sets for the overlapping Schwarz
1132b1a0a8a3SJed Brown    preconditioner for a any problem on a general grid.
1133b1a0a8a3SJed Brown 
1134b1a0a8a3SJed Brown    Collective
1135b1a0a8a3SJed Brown 
1136b1a0a8a3SJed Brown    Input Parameters:
1137b1a0a8a3SJed Brown +  A - The global matrix operator
1138b1a0a8a3SJed Brown -  n - the number of local blocks
1139b1a0a8a3SJed Brown 
1140b1a0a8a3SJed Brown    Output Parameters:
1141b1a0a8a3SJed Brown .  outis - the array of index sets defining the subdomains
1142b1a0a8a3SJed Brown 
1143b1a0a8a3SJed Brown    Level: advanced
1144b1a0a8a3SJed Brown 
1145f746d493SDmitry Karpeev    Note: this generates nonoverlapping subdomains; the PCGASM will generate the overlap
1146f746d493SDmitry Karpeev     from these if you use PCGASMSetLocalSubdomains()
1147b1a0a8a3SJed Brown 
1148b1a0a8a3SJed Brown     In the Fortran version you must provide the array outis[] already allocated of length n.
1149b1a0a8a3SJed Brown 
1150f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, unstructured grid
1151b1a0a8a3SJed Brown 
1152f746d493SDmitry Karpeev .seealso: PCGASMSetLocalSubdomains(), PCGASMDestroySubdomains()
1153b1a0a8a3SJed Brown @*/
11547087cfbeSBarry Smith PetscErrorCode  PCGASMCreateSubdomains(Mat A, PetscInt n, IS* outis[])
1155b1a0a8a3SJed Brown {
1156b1a0a8a3SJed Brown   MatPartitioning           mpart;
1157b1a0a8a3SJed Brown   const char                *prefix;
115811bd1e4dSLisandro Dalcin   PetscErrorCode            (*f)(Mat,MatReuse,Mat*);
1159b1a0a8a3SJed Brown   PetscMPIInt               size;
1160b1a0a8a3SJed Brown   PetscInt                  i,j,rstart,rend,bs;
116111bd1e4dSLisandro Dalcin   PetscBool                 isbaij = PETSC_FALSE,foundpart = PETSC_FALSE;
1162b1a0a8a3SJed Brown   Mat                       Ad = PETSC_NULL, adj;
1163b1a0a8a3SJed Brown   IS                        ispart,isnumb,*is;
1164b1a0a8a3SJed Brown   PetscErrorCode            ierr;
1165b1a0a8a3SJed Brown 
1166b1a0a8a3SJed Brown   PetscFunctionBegin;
1167b1a0a8a3SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1168b1a0a8a3SJed Brown   PetscValidPointer(outis,3);
1169b1a0a8a3SJed Brown   if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of local blocks must be > 0, n = %D",n);
1170b1a0a8a3SJed Brown 
1171b1a0a8a3SJed Brown   /* Get prefix, row distribution, and block size */
1172b1a0a8a3SJed Brown   ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr);
1173b1a0a8a3SJed Brown   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
1174b1a0a8a3SJed Brown   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
1175b1a0a8a3SJed 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);
1176b1a0a8a3SJed Brown 
1177b1a0a8a3SJed Brown   /* Get diagonal block from matrix if possible */
1178b1a0a8a3SJed Brown   ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr);
1179b1a0a8a3SJed Brown   ierr = PetscObjectQueryFunction((PetscObject)A,"MatGetDiagonalBlock_C",(void (**)(void))&f);CHKERRQ(ierr);
1180b1a0a8a3SJed Brown   if (f) {
118111bd1e4dSLisandro Dalcin     ierr = MatGetDiagonalBlock(A,&Ad);CHKERRQ(ierr);
1182b1a0a8a3SJed Brown   } else if (size == 1) {
118311bd1e4dSLisandro Dalcin     Ad = A;
1184b1a0a8a3SJed Brown   }
1185b1a0a8a3SJed Brown   if (Ad) {
1186b1a0a8a3SJed Brown     ierr = PetscTypeCompare((PetscObject)Ad,MATSEQBAIJ,&isbaij);CHKERRQ(ierr);
1187b1a0a8a3SJed Brown     if (!isbaij) {ierr = PetscTypeCompare((PetscObject)Ad,MATSEQSBAIJ,&isbaij);CHKERRQ(ierr);}
1188b1a0a8a3SJed Brown   }
1189b1a0a8a3SJed Brown   if (Ad && n > 1) {
1190b1a0a8a3SJed Brown     PetscBool  match,done;
1191b1a0a8a3SJed Brown     /* Try to setup a good matrix partitioning if available */
1192b1a0a8a3SJed Brown     ierr = MatPartitioningCreate(PETSC_COMM_SELF,&mpart);CHKERRQ(ierr);
1193b1a0a8a3SJed Brown     ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr);
1194b1a0a8a3SJed Brown     ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr);
1195b1a0a8a3SJed Brown     ierr = PetscTypeCompare((PetscObject)mpart,MATPARTITIONINGCURRENT,&match);CHKERRQ(ierr);
1196b1a0a8a3SJed Brown     if (!match) {
1197b1a0a8a3SJed Brown       ierr = PetscTypeCompare((PetscObject)mpart,MATPARTITIONINGSQUARE,&match);CHKERRQ(ierr);
1198b1a0a8a3SJed Brown     }
1199b1a0a8a3SJed Brown     if (!match) { /* assume a "good" partitioner is available */
1200b1a0a8a3SJed Brown       PetscInt na,*ia,*ja;
1201b1a0a8a3SJed Brown       ierr = MatGetRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr);
1202b1a0a8a3SJed Brown       if (done) {
1203b1a0a8a3SJed Brown         /* Build adjacency matrix by hand. Unfortunately a call to
1204b1a0a8a3SJed Brown            MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will
1205b1a0a8a3SJed Brown            remove the block-aij structure and we cannot expect
1206b1a0a8a3SJed Brown            MatPartitioning to split vertices as we need */
1207b1a0a8a3SJed Brown         PetscInt i,j,*row,len,nnz,cnt,*iia=0,*jja=0;
1208b1a0a8a3SJed Brown         nnz = 0;
1209b1a0a8a3SJed Brown         for (i=0; i<na; i++) { /* count number of nonzeros */
1210b1a0a8a3SJed Brown           len = ia[i+1] - ia[i];
1211b1a0a8a3SJed Brown           row = ja + ia[i];
1212b1a0a8a3SJed Brown           for (j=0; j<len; j++) {
1213b1a0a8a3SJed Brown             if (row[j] == i) { /* don't count diagonal */
1214b1a0a8a3SJed Brown               len--; break;
1215b1a0a8a3SJed Brown             }
1216b1a0a8a3SJed Brown           }
1217b1a0a8a3SJed Brown           nnz += len;
1218b1a0a8a3SJed Brown         }
1219b1a0a8a3SJed Brown         ierr = PetscMalloc((na+1)*sizeof(PetscInt),&iia);CHKERRQ(ierr);
1220b1a0a8a3SJed Brown         ierr = PetscMalloc((nnz)*sizeof(PetscInt),&jja);CHKERRQ(ierr);
1221b1a0a8a3SJed Brown         nnz    = 0;
1222b1a0a8a3SJed Brown         iia[0] = 0;
1223b1a0a8a3SJed Brown         for (i=0; i<na; i++) { /* fill adjacency */
1224b1a0a8a3SJed Brown           cnt = 0;
1225b1a0a8a3SJed Brown           len = ia[i+1] - ia[i];
1226b1a0a8a3SJed Brown           row = ja + ia[i];
1227b1a0a8a3SJed Brown           for (j=0; j<len; j++) {
1228b1a0a8a3SJed Brown             if (row[j] != i) { /* if not diagonal */
1229b1a0a8a3SJed Brown               jja[nnz+cnt++] = row[j];
1230b1a0a8a3SJed Brown             }
1231b1a0a8a3SJed Brown           }
1232b1a0a8a3SJed Brown           nnz += cnt;
1233b1a0a8a3SJed Brown           iia[i+1] = nnz;
1234b1a0a8a3SJed Brown         }
1235b1a0a8a3SJed Brown         /* Partitioning of the adjacency matrix */
1236b1a0a8a3SJed Brown         ierr = MatCreateMPIAdj(PETSC_COMM_SELF,na,na,iia,jja,PETSC_NULL,&adj);CHKERRQ(ierr);
1237b1a0a8a3SJed Brown         ierr = MatPartitioningSetAdjacency(mpart,adj);CHKERRQ(ierr);
1238b1a0a8a3SJed Brown         ierr = MatPartitioningSetNParts(mpart,n);CHKERRQ(ierr);
1239b1a0a8a3SJed Brown         ierr = MatPartitioningApply(mpart,&ispart);CHKERRQ(ierr);
1240b1a0a8a3SJed Brown         ierr = ISPartitioningToNumbering(ispart,&isnumb);CHKERRQ(ierr);
12416bf464f9SBarry Smith         ierr = MatDestroy(&adj);CHKERRQ(ierr);
1242b1a0a8a3SJed Brown         foundpart = PETSC_TRUE;
1243b1a0a8a3SJed Brown       }
1244b1a0a8a3SJed Brown       ierr = MatRestoreRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr);
1245b1a0a8a3SJed Brown     }
12466bf464f9SBarry Smith     ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr);
1247b1a0a8a3SJed Brown   }
1248b1a0a8a3SJed Brown 
1249b1a0a8a3SJed Brown   ierr = PetscMalloc(n*sizeof(IS),&is);CHKERRQ(ierr);
1250b1a0a8a3SJed Brown   *outis = is;
1251b1a0a8a3SJed Brown 
1252b1a0a8a3SJed Brown   if (!foundpart) {
1253b1a0a8a3SJed Brown 
1254b1a0a8a3SJed Brown     /* Partitioning by contiguous chunks of rows */
1255b1a0a8a3SJed Brown 
1256b1a0a8a3SJed Brown     PetscInt mbs   = (rend-rstart)/bs;
1257b1a0a8a3SJed Brown     PetscInt start = rstart;
1258b1a0a8a3SJed Brown     for (i=0; i<n; i++) {
1259b1a0a8a3SJed Brown       PetscInt count = (mbs/n + ((mbs % n) > i)) * bs;
1260b1a0a8a3SJed Brown       ierr   = ISCreateStride(PETSC_COMM_SELF,count,start,1,&is[i]);CHKERRQ(ierr);
1261b1a0a8a3SJed Brown       start += count;
1262b1a0a8a3SJed Brown     }
1263b1a0a8a3SJed Brown 
1264b1a0a8a3SJed Brown   } else {
1265b1a0a8a3SJed Brown 
1266b1a0a8a3SJed Brown     /* Partitioning by adjacency of diagonal block  */
1267b1a0a8a3SJed Brown 
1268b1a0a8a3SJed Brown     const PetscInt *numbering;
1269b1a0a8a3SJed Brown     PetscInt       *count,nidx,*indices,*newidx,start=0;
1270b1a0a8a3SJed Brown     /* Get node count in each partition */
1271b1a0a8a3SJed Brown     ierr = PetscMalloc(n*sizeof(PetscInt),&count);CHKERRQ(ierr);
1272b1a0a8a3SJed Brown     ierr = ISPartitioningCount(ispart,n,count);CHKERRQ(ierr);
1273b1a0a8a3SJed Brown     if (isbaij && bs > 1) { /* adjust for the block-aij case */
1274b1a0a8a3SJed Brown       for (i=0; i<n; i++) count[i] *= bs;
1275b1a0a8a3SJed Brown     }
1276b1a0a8a3SJed Brown     /* Build indices from node numbering */
1277b1a0a8a3SJed Brown     ierr = ISGetLocalSize(isnumb,&nidx);CHKERRQ(ierr);
1278b1a0a8a3SJed Brown     ierr = PetscMalloc(nidx*sizeof(PetscInt),&indices);CHKERRQ(ierr);
1279b1a0a8a3SJed Brown     for (i=0; i<nidx; i++) indices[i] = i; /* needs to be initialized */
1280b1a0a8a3SJed Brown     ierr = ISGetIndices(isnumb,&numbering);CHKERRQ(ierr);
1281b1a0a8a3SJed Brown     ierr = PetscSortIntWithPermutation(nidx,numbering,indices);CHKERRQ(ierr);
1282b1a0a8a3SJed Brown     ierr = ISRestoreIndices(isnumb,&numbering);CHKERRQ(ierr);
1283b1a0a8a3SJed Brown     if (isbaij && bs > 1) { /* adjust for the block-aij case */
1284b1a0a8a3SJed Brown       ierr = PetscMalloc(nidx*bs*sizeof(PetscInt),&newidx);CHKERRQ(ierr);
1285b1a0a8a3SJed Brown       for (i=0; i<nidx; i++)
1286b1a0a8a3SJed Brown         for (j=0; j<bs; j++)
1287b1a0a8a3SJed Brown           newidx[i*bs+j] = indices[i]*bs + j;
1288b1a0a8a3SJed Brown       ierr = PetscFree(indices);CHKERRQ(ierr);
1289b1a0a8a3SJed Brown       nidx   *= bs;
1290b1a0a8a3SJed Brown       indices = newidx;
1291b1a0a8a3SJed Brown     }
1292b1a0a8a3SJed Brown     /* Shift to get global indices */
1293b1a0a8a3SJed Brown     for (i=0; i<nidx; i++) indices[i] += rstart;
1294b1a0a8a3SJed Brown 
1295b1a0a8a3SJed Brown     /* Build the index sets for each block */
1296b1a0a8a3SJed Brown     for (i=0; i<n; i++) {
1297b1a0a8a3SJed Brown       ierr   = ISCreateGeneral(PETSC_COMM_SELF,count[i],&indices[start],PETSC_COPY_VALUES,&is[i]);CHKERRQ(ierr);
1298b1a0a8a3SJed Brown       ierr   = ISSort(is[i]);CHKERRQ(ierr);
1299b1a0a8a3SJed Brown       start += count[i];
1300b1a0a8a3SJed Brown     }
1301b1a0a8a3SJed Brown 
1302b1a0a8a3SJed Brown     ierr = PetscFree(count);
1303b1a0a8a3SJed Brown     ierr = PetscFree(indices);
1304fcfd50ebSBarry Smith     ierr = ISDestroy(&isnumb);CHKERRQ(ierr);
1305fcfd50ebSBarry Smith     ierr = ISDestroy(&ispart);CHKERRQ(ierr);
1306b1a0a8a3SJed Brown   }
1307b1a0a8a3SJed Brown 
1308b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1309b1a0a8a3SJed Brown }
1310b1a0a8a3SJed Brown 
1311b1a0a8a3SJed Brown #undef __FUNCT__
1312f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMDestroySubdomains"
1313b1a0a8a3SJed Brown /*@C
1314f746d493SDmitry Karpeev    PCGASMDestroySubdomains - Destroys the index sets created with
1315f746d493SDmitry Karpeev    PCGASMCreateSubdomains(). Should be called after setting subdomains
1316f746d493SDmitry Karpeev    with PCGASMSetLocalSubdomains().
1317b1a0a8a3SJed Brown 
1318b1a0a8a3SJed Brown    Collective
1319b1a0a8a3SJed Brown 
1320b1a0a8a3SJed Brown    Input Parameters:
1321b1a0a8a3SJed Brown +  n - the number of index sets
1322b1a0a8a3SJed Brown .  is - the array of index sets
1323b1a0a8a3SJed Brown -  is_local - the array of local index sets, can be PETSC_NULL
1324b1a0a8a3SJed Brown 
1325b1a0a8a3SJed Brown    Level: advanced
1326b1a0a8a3SJed Brown 
1327f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, unstructured grid
1328b1a0a8a3SJed Brown 
1329f746d493SDmitry Karpeev .seealso: PCGASMCreateSubdomains(), PCGASMSetLocalSubdomains()
1330b1a0a8a3SJed Brown @*/
13317087cfbeSBarry Smith PetscErrorCode  PCGASMDestroySubdomains(PetscInt n, IS is[], IS is_local[])
1332b1a0a8a3SJed Brown {
1333b1a0a8a3SJed Brown   PetscInt       i;
1334b1a0a8a3SJed Brown   PetscErrorCode ierr;
1335b1a0a8a3SJed Brown   PetscFunctionBegin;
1336b1a0a8a3SJed Brown   if (n <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"n must be > 0: n = %D",n);
1337b1a0a8a3SJed Brown   PetscValidPointer(is,2);
1338fcfd50ebSBarry Smith   for (i=0; i<n; i++) { ierr = ISDestroy(&is[i]);CHKERRQ(ierr); }
1339b1a0a8a3SJed Brown   ierr = PetscFree(is);CHKERRQ(ierr);
1340b1a0a8a3SJed Brown   if (is_local) {
1341b1a0a8a3SJed Brown     PetscValidPointer(is_local,3);
1342fcfd50ebSBarry Smith     for (i=0; i<n; i++) { ierr = ISDestroy(&is_local[i]);CHKERRQ(ierr); }
1343b1a0a8a3SJed Brown     ierr = PetscFree(is_local);CHKERRQ(ierr);
1344b1a0a8a3SJed Brown   }
1345b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1346b1a0a8a3SJed Brown }
1347b1a0a8a3SJed Brown 
1348af538404SDmitry Karpeev 
1349af538404SDmitry Karpeev #define PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,xleft_loc,ylow_loc,xright_loc,yhigh_loc,n) \
1350af538404SDmitry Karpeev {                                                                                                       \
1351af538404SDmitry Karpeev  PetscInt first_row = first/M, last_row = last/M+1;                                                     \
1352af538404SDmitry Karpeev   /*                                                                                                    \
1353af538404SDmitry Karpeev    Compute ylow_loc and yhigh_loc so that (ylow_loc,xleft) and (yhigh_loc,xright) are the corners       \
1354af538404SDmitry Karpeev    of the bounding box of the intersection of the subdomain with the local ownership range (local       \
1355af538404SDmitry Karpeev    subdomain).                                                                                          \
1356af538404SDmitry Karpeev    Also compute xleft_loc and xright_loc as the lower and upper bounds on the first and last rows       \
1357af538404SDmitry Karpeev    of the intersection.                                                                                 \
1358af538404SDmitry Karpeev   */                                                                                                    \
1359af538404SDmitry Karpeev   /* ylow_loc is the grid row containing the first element of the local sumbdomain */                   \
1360eec7e3faSDmitry Karpeev   *ylow_loc = PetscMax(first_row,ylow);                                                                    \
1361af538404SDmitry Karpeev   /* xleft_loc is the offset of first element of the local subdomain within its grid row (might actually be outside the local subdomain) */ \
1362eec7e3faSDmitry Karpeev   *xleft_loc = *ylow_loc==first_row?PetscMax(first%M,xleft):xleft;                                                                            \
1363af538404SDmitry Karpeev   /* yhigh_loc is the grid row above the last local subdomain element */                                                                    \
1364eec7e3faSDmitry Karpeev   *yhigh_loc = PetscMin(last_row,yhigh);                                                                                                     \
1365af538404SDmitry Karpeev   /* xright is the offset of the end of the  local subdomain within its grid row (might actually be outside the local subdomain) */         \
1366eec7e3faSDmitry Karpeev   *xright_loc = *yhigh_loc==last_row?PetscMin(xright,last%M):xright;                                                                          \
1367af538404SDmitry Karpeev   /* Now compute the size of the local subdomain n. */ \
1368c3518bceSDmitry Karpeev   *n = 0;                                               \
1369eec7e3faSDmitry Karpeev   if(*ylow_loc < *yhigh_loc) {                           \
1370af538404SDmitry Karpeev     PetscInt width = xright-xleft;                     \
1371eec7e3faSDmitry Karpeev     *n += width*(*yhigh_loc-*ylow_loc-1);                 \
1372eec7e3faSDmitry Karpeev     *n += PetscMin(PetscMax(*xright_loc-xleft,0),width); \
1373eec7e3faSDmitry Karpeev     *n -= PetscMin(PetscMax(*xleft_loc-xleft,0), width); \
1374af538404SDmitry Karpeev   }\
1375af538404SDmitry Karpeev }
1376af538404SDmitry Karpeev 
1377c3518bceSDmitry Karpeev 
1378eec7e3faSDmitry Karpeev 
1379eec7e3faSDmitry Karpeev #undef __FUNCT__
1380f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMCreateSubdomains2D"
1381b1a0a8a3SJed Brown /*@
1382f746d493SDmitry Karpeev    PCGASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz
1383b1a0a8a3SJed Brown    preconditioner for a two-dimensional problem on a regular grid.
1384b1a0a8a3SJed Brown 
1385af538404SDmitry Karpeev    Collective
1386b1a0a8a3SJed Brown 
1387b1a0a8a3SJed Brown    Input Parameters:
1388af538404SDmitry Karpeev +  M, N - the global number of mesh points in the x and y directions
1389af538404SDmitry Karpeev .  Mdomains, Ndomains - the global number of subdomains in the x and y directions
1390b1a0a8a3SJed Brown .  dof - degrees of freedom per node
1391b1a0a8a3SJed Brown -  overlap - overlap in mesh lines
1392b1a0a8a3SJed Brown 
1393b1a0a8a3SJed Brown    Output Parameters:
1394af538404SDmitry Karpeev +  Nsub - the number of local subdomains created
1395b1a0a8a3SJed Brown .  is - array of index sets defining overlapping (if overlap > 0) subdomains
1396b1a0a8a3SJed Brown -  is_local - array of index sets defining non-overlapping subdomains
1397b1a0a8a3SJed Brown 
1398b1a0a8a3SJed Brown    Note:
1399b1a0a8a3SJed Brown    Presently PCAMSCreateSubdomains2d() is valid only for sequential
1400b1a0a8a3SJed Brown    preconditioners.  More general related routines are
1401f746d493SDmitry Karpeev    PCGASMSetTotalSubdomains() and PCGASMSetLocalSubdomains().
1402b1a0a8a3SJed Brown 
1403b1a0a8a3SJed Brown    Level: advanced
1404b1a0a8a3SJed Brown 
1405f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, 2D, regular grid
1406b1a0a8a3SJed Brown 
1407f746d493SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetLocalSubdomains(), PCGASMGetSubKSP(),
1408f746d493SDmitry Karpeev           PCGASMSetOverlap()
1409b1a0a8a3SJed Brown @*/
14107087cfbeSBarry Smith PetscErrorCode  PCGASMCreateSubdomains2D(PC pc, PetscInt M,PetscInt N,PetscInt Mdomains,PetscInt Ndomains,PetscInt dof,PetscInt overlap, PetscInt *nsub,IS **is,IS **is_local)
1411b1a0a8a3SJed Brown {
1412b1a0a8a3SJed Brown   PetscErrorCode ierr;
14132bbb417fSDmitry Karpeev   PetscMPIInt    size, rank;
14142bbb417fSDmitry Karpeev   PetscInt       i, j;
14152bbb417fSDmitry Karpeev   PetscInt       maxheight, maxwidth;
14162bbb417fSDmitry Karpeev   PetscInt       xstart, xleft, xright, xleft_loc, xright_loc;
14172bbb417fSDmitry Karpeev   PetscInt       ystart, ylow,  yhigh,  ylow_loc,  yhigh_loc;
1418eec7e3faSDmitry Karpeev   PetscInt       x[2][2], y[2][2], n[2];
14192bbb417fSDmitry Karpeev   PetscInt       first, last;
14202bbb417fSDmitry Karpeev   PetscInt       nidx, *idx;
14212bbb417fSDmitry Karpeev   PetscInt       ii,jj,s,q,d;
1422af538404SDmitry Karpeev   PetscInt       k,kk;
14232bbb417fSDmitry Karpeev   PetscMPIInt    color;
14242bbb417fSDmitry Karpeev   MPI_Comm       comm, subcomm;
14255d38c402SBarry Smith   IS             **iis = 0;
1426d34fcf5fSBarry Smith 
1427b1a0a8a3SJed Brown   PetscFunctionBegin;
14282bbb417fSDmitry Karpeev   ierr = PetscObjectGetComm((PetscObject)pc, &comm);CHKERRQ(ierr);
14292bbb417fSDmitry Karpeev   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
14302bbb417fSDmitry Karpeev   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
14312bbb417fSDmitry Karpeev   ierr = MatGetOwnershipRange(pc->pmat, &first, &last);CHKERRQ(ierr);
1432d34fcf5fSBarry Smith   if (first%dof || last%dof) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Matrix row partitioning unsuitable for domain decomposition: local row range (%D,%D) "
14332bbb417fSDmitry Karpeev 	     "does not respect the number of degrees of freedom per grid point %D", first, last, dof);
1434d34fcf5fSBarry Smith 
1435af538404SDmitry Karpeev   /* Determine the number of domains with nonzero intersections with the local ownership range. */
14362bbb417fSDmitry Karpeev   s = 0;
1437af538404SDmitry Karpeev   ystart = 0;
1438af538404SDmitry Karpeev   for (j=0; j<Ndomains; ++j) {
1439af538404SDmitry Karpeev     maxheight = N/Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */
1440af538404SDmitry Karpeev     if (maxheight < 2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many %D subdomains in the vertical directon for mesh height %D", Ndomains, N);
1441eec7e3faSDmitry Karpeev     /* Vertical domain limits with an overlap. */
1442eec7e3faSDmitry Karpeev     ylow = PetscMax(ystart - overlap,0);
1443eec7e3faSDmitry Karpeev     yhigh = PetscMin(ystart + maxheight + overlap,N);
1444b1a0a8a3SJed Brown     xstart = 0;
1445af538404SDmitry Karpeev     for (i=0; i<Mdomains; ++i) {
1446af538404SDmitry Karpeev       maxwidth = M/Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */
1447af538404SDmitry Karpeev       if (maxwidth < 2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many %D subdomains in the horizontal direction for mesh width %D", Mdomains, M);
1448eec7e3faSDmitry Karpeev       /* Horizontal domain limits with an overlap. */
1449eec7e3faSDmitry Karpeev       xleft   = PetscMax(xstart - overlap,0);
1450eec7e3faSDmitry Karpeev       xright  = PetscMin(xstart + maxwidth + overlap,M);
1451af538404SDmitry Karpeev       /*
1452af538404SDmitry Karpeev 	 Determine whether this subdomain intersects this processor's ownership range of pc->pmat.
1453af538404SDmitry Karpeev       */
1454c3518bceSDmitry Karpeev       PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,(&xleft_loc),(&ylow_loc),(&xright_loc),(&yhigh_loc),(&nidx));
1455af538404SDmitry Karpeev       if(nidx) {
1456af538404SDmitry Karpeev         ++s;
1457af538404SDmitry Karpeev       }
1458af538404SDmitry Karpeev       xstart += maxwidth;
1459af538404SDmitry Karpeev     }/* for(i = 0; i < Mdomains; ++i) */
1460af538404SDmitry Karpeev     ystart += maxheight;
1461af538404SDmitry Karpeev   }/* for(j = 0; j < Ndomains; ++j) */
1462af538404SDmitry Karpeev   /* Now we can allocate the necessary number of ISs. */
1463af538404SDmitry Karpeev   *nsub = s;
1464af538404SDmitry Karpeev   ierr = PetscMalloc((*nsub)*sizeof(IS*),is);CHKERRQ(ierr);
1465af538404SDmitry Karpeev   ierr = PetscMalloc((*nsub)*sizeof(IS*),is_local);CHKERRQ(ierr);
1466af538404SDmitry Karpeev   s = 0;
1467af538404SDmitry Karpeev   ystart = 0;
1468af538404SDmitry Karpeev   for (j=0; j<Ndomains; ++j) {
1469af538404SDmitry Karpeev     maxheight = N/Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */
1470af538404SDmitry Karpeev     if (maxheight < 2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many %D subdomains in the vertical directon for mesh height %D", Ndomains, N);
1471af538404SDmitry Karpeev     /* Vertical domain limits with an overlap. */
1472eec7e3faSDmitry Karpeev     y[0][0] = PetscMax(ystart - overlap,0);
1473eec7e3faSDmitry Karpeev     y[0][1] = PetscMin(ystart + maxheight + overlap,N);
1474eec7e3faSDmitry Karpeev     /* Vertical domain limits without an overlap. */
1475eec7e3faSDmitry Karpeev     y[1][0] = ystart;
1476eec7e3faSDmitry Karpeev     y[1][1] = PetscMin(ystart + maxheight,N);
1477eec7e3faSDmitry Karpeev     xstart = 0;
1478af538404SDmitry Karpeev     for (i=0; i<Mdomains; ++i) {
1479af538404SDmitry Karpeev       maxwidth = M/Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */
1480af538404SDmitry Karpeev       if (maxwidth < 2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many %D subdomains in the horizontal direction for mesh width %D", Mdomains, M);
1481af538404SDmitry Karpeev       /* Horizontal domain limits with an overlap. */
1482eec7e3faSDmitry Karpeev       x[0][0]  = PetscMax(xstart - overlap,0);
1483eec7e3faSDmitry Karpeev       x[0][1]  = PetscMin(xstart + maxwidth + overlap,M);
1484eec7e3faSDmitry Karpeev       /* Horizontal domain limits without an overlap. */
1485eec7e3faSDmitry Karpeev       x[1][0] = xstart;
1486eec7e3faSDmitry Karpeev       x[1][1] = PetscMin(xstart+maxwidth,M);
14872bbb417fSDmitry Karpeev       /*
14882bbb417fSDmitry Karpeev 	 Determine whether this domain intersects this processor's ownership range of pc->pmat.
14892bbb417fSDmitry Karpeev 	 Do this twice: first for the domains with overlaps, and once without.
14902bbb417fSDmitry Karpeev 	 During the first pass create the subcommunicators, and use them on the second pass as well.
14912bbb417fSDmitry Karpeev       */
14922bbb417fSDmitry Karpeev       for(q = 0; q < 2; ++q) {
14932bbb417fSDmitry Karpeev 	/*
14942bbb417fSDmitry Karpeev 	  domain limits, (xleft, xright) and (ylow, yheigh) are adjusted
14952bbb417fSDmitry Karpeev 	  according to whether the domain with an overlap or without is considered.
14962bbb417fSDmitry Karpeev 	*/
14972bbb417fSDmitry Karpeev 	xleft = x[q][0]; xright = x[q][1];
14982bbb417fSDmitry Karpeev 	ylow  = y[q][0]; yhigh  = y[q][1];
1499c3518bceSDmitry Karpeev         PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,(&xleft_loc),(&ylow_loc),(&xright_loc),(&yhigh_loc),(&nidx));
1500af538404SDmitry Karpeev 	nidx *= dof;
1501eec7e3faSDmitry Karpeev         n[q] = nidx;
15022bbb417fSDmitry Karpeev         /*
1503eec7e3faSDmitry Karpeev          Based on the counted number of indices in the local domain *with an overlap*,
1504af538404SDmitry Karpeev          construct a subcommunicator of all the processors supporting this domain.
1505eec7e3faSDmitry Karpeev          Observe that a domain with an overlap might have nontrivial local support,
1506eec7e3faSDmitry Karpeev          while the domain without an overlap might not.  Hence, the decision to participate
1507eec7e3faSDmitry Karpeev          in the subcommunicator must be based on the domain with an overlap.
15082bbb417fSDmitry Karpeev          */
15092bbb417fSDmitry Karpeev 	if (q == 0) {
15102bbb417fSDmitry Karpeev 	  if(nidx) {
15112bbb417fSDmitry Karpeev 	    color = 1;
1512d34fcf5fSBarry Smith 	  } else {
15132bbb417fSDmitry Karpeev 	    color = MPI_UNDEFINED;
15142bbb417fSDmitry Karpeev 	  }
15152bbb417fSDmitry Karpeev 	  ierr = MPI_Comm_split(comm, color, rank, &subcomm);CHKERRQ(ierr);
15162bbb417fSDmitry Karpeev 	}
1517af538404SDmitry Karpeev         /*
1518eec7e3faSDmitry Karpeev          Proceed only if the number of local indices *with an overlap* is nonzero.
1519af538404SDmitry Karpeev          */
1520eec7e3faSDmitry Karpeev         if (n[0]) {
1521af538404SDmitry Karpeev           if(q == 0) {
1522eec7e3faSDmitry Karpeev             iis = is;
1523af538404SDmitry Karpeev           }
1524af538404SDmitry Karpeev           if (q == 1) {
1525af538404SDmitry Karpeev             /*
1526eec7e3faSDmitry Karpeev              The IS for the no-overlap subdomain shares a communicator with the overlapping domain.
1527af538404SDmitry Karpeev              Moreover, if the overlap is zero, the two ISs are identical.
1528af538404SDmitry Karpeev              */
1529b1a0a8a3SJed Brown             if (overlap == 0) {
1530eec7e3faSDmitry Karpeev               (*is_local)[s] = (*is)[s];
15312bbb417fSDmitry Karpeev               ierr = PetscObjectReference((PetscObject)(*is)[s]);CHKERRQ(ierr);
15322bbb417fSDmitry Karpeev               continue;
1533d34fcf5fSBarry Smith             } else {
1534eec7e3faSDmitry Karpeev               iis = is_local;
1535eec7e3faSDmitry Karpeev               subcomm = ((PetscObject)(*is)[s])->comm;
15362bbb417fSDmitry Karpeev             }
1537af538404SDmitry Karpeev           }/* if(q == 1) */
1538eec7e3faSDmitry Karpeev           idx = PETSC_NULL;
15392bbb417fSDmitry Karpeev 	  ierr = PetscMalloc(nidx*sizeof(PetscInt),&idx);CHKERRQ(ierr);
1540eec7e3faSDmitry Karpeev           if(nidx) {
15412bbb417fSDmitry Karpeev             k    = 0;
15422bbb417fSDmitry Karpeev             for (jj=ylow_loc; jj<yhigh_loc; ++jj) {
1543af538404SDmitry Karpeev               PetscInt x0 = (jj==ylow_loc)?xleft_loc:xleft;
1544af538404SDmitry Karpeev               PetscInt x1 = (jj==yhigh_loc-1)?xright_loc:xright;
1545af538404SDmitry Karpeev               kk = dof*(M*jj + x0);
15462bbb417fSDmitry Karpeev               for (ii=x0; ii<x1; ++ii) {
15472bbb417fSDmitry Karpeev                 for(d = 0; d < dof; ++d) {
15482bbb417fSDmitry Karpeev                   idx[k++] = kk++;
1549b1a0a8a3SJed Brown                 }
1550b1a0a8a3SJed Brown               }
1551b1a0a8a3SJed Brown             }
1552eec7e3faSDmitry Karpeev           }
15532bbb417fSDmitry Karpeev 	  ierr = ISCreateGeneral(subcomm,nidx,idx,PETSC_OWN_POINTER,(*iis)+s);CHKERRQ(ierr);
1554eec7e3faSDmitry Karpeev 	}/* if(n[0]) */
15552bbb417fSDmitry Karpeev       }/* for(q = 0; q < 2; ++q) */
1556eec7e3faSDmitry Karpeev       if(n[0]) {
15572bbb417fSDmitry Karpeev         ++s;
1558b1a0a8a3SJed Brown       }
1559af538404SDmitry Karpeev       xstart += maxwidth;
1560af538404SDmitry Karpeev     }/* for(i = 0; i < Mdomains; ++i) */
15612bbb417fSDmitry Karpeev     ystart += maxheight;
1562af538404SDmitry Karpeev   }/* for(j = 0; j < Ndomains; ++j) */
1563b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1564b1a0a8a3SJed Brown }
1565b1a0a8a3SJed Brown 
1566b1a0a8a3SJed Brown #undef __FUNCT__
1567f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMGetLocalSubdomains"
1568b1a0a8a3SJed Brown /*@C
1569f746d493SDmitry Karpeev     PCGASMGetLocalSubdomains - Gets the local subdomains (for this processor
1570b1a0a8a3SJed Brown     only) for the additive Schwarz preconditioner.
1571b1a0a8a3SJed Brown 
1572b1a0a8a3SJed Brown     Not Collective
1573b1a0a8a3SJed Brown 
1574b1a0a8a3SJed Brown     Input Parameter:
1575b1a0a8a3SJed Brown .   pc - the preconditioner context
1576b1a0a8a3SJed Brown 
1577b1a0a8a3SJed Brown     Output Parameters:
1578b1a0a8a3SJed Brown +   n - the number of subdomains for this processor (default value = 1)
1579b1a0a8a3SJed Brown .   is - the index sets that define the subdomains for this processor
1580b1a0a8a3SJed Brown -   is_local - the index sets that define the local part of the subdomains for this processor (can be PETSC_NULL)
1581b1a0a8a3SJed Brown 
1582b1a0a8a3SJed Brown 
1583b1a0a8a3SJed Brown     Notes:
1584b1a0a8a3SJed Brown     The IS numbering is in the parallel, global numbering of the vector.
1585b1a0a8a3SJed Brown 
1586b1a0a8a3SJed Brown     Level: advanced
1587b1a0a8a3SJed Brown 
1588f746d493SDmitry Karpeev .keywords: PC, GASM, set, local, subdomains, additive Schwarz
1589b1a0a8a3SJed Brown 
1590f746d493SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(),
1591f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D(), PCGASMSetLocalSubdomains(), PCGASMGetLocalSubmatrices()
1592b1a0a8a3SJed Brown @*/
15937087cfbeSBarry Smith PetscErrorCode  PCGASMGetLocalSubdomains(PC pc,PetscInt *n,IS *is[],IS *is_local[])
1594b1a0a8a3SJed Brown {
1595f746d493SDmitry Karpeev   PC_GASM         *osm;
1596b1a0a8a3SJed Brown   PetscErrorCode ierr;
1597b1a0a8a3SJed Brown   PetscBool      match;
1598b1a0a8a3SJed Brown 
1599b1a0a8a3SJed Brown   PetscFunctionBegin;
1600b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1601b1a0a8a3SJed Brown   PetscValidIntPointer(n,2);
1602b1a0a8a3SJed Brown   if (is) PetscValidPointer(is,3);
1603f746d493SDmitry Karpeev   ierr = PetscTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr);
1604b1a0a8a3SJed Brown   if (!match) {
1605b1a0a8a3SJed Brown     if (n)  *n  = 0;
1606b1a0a8a3SJed Brown     if (is) *is = PETSC_NULL;
1607b1a0a8a3SJed Brown   } else {
1608f746d493SDmitry Karpeev     osm = (PC_GASM*)pc->data;
1609ab3e923fSDmitry Karpeev     if (n)  *n  = osm->n;
1610b1a0a8a3SJed Brown     if (is) *is = osm->is;
1611b1a0a8a3SJed Brown     if (is_local) *is_local = osm->is_local;
1612b1a0a8a3SJed Brown   }
1613b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1614b1a0a8a3SJed Brown }
1615b1a0a8a3SJed Brown 
1616b1a0a8a3SJed Brown #undef __FUNCT__
1617f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMGetLocalSubmatrices"
1618b1a0a8a3SJed Brown /*@C
1619f746d493SDmitry Karpeev     PCGASMGetLocalSubmatrices - Gets the local submatrices (for this processor
1620b1a0a8a3SJed Brown     only) for the additive Schwarz preconditioner.
1621b1a0a8a3SJed Brown 
1622b1a0a8a3SJed Brown     Not Collective
1623b1a0a8a3SJed Brown 
1624b1a0a8a3SJed Brown     Input Parameter:
1625b1a0a8a3SJed Brown .   pc - the preconditioner context
1626b1a0a8a3SJed Brown 
1627b1a0a8a3SJed Brown     Output Parameters:
1628b1a0a8a3SJed Brown +   n - the number of matrices for this processor (default value = 1)
1629b1a0a8a3SJed Brown -   mat - the matrices
1630b1a0a8a3SJed Brown 
1631b1a0a8a3SJed Brown 
1632b1a0a8a3SJed Brown     Level: advanced
1633b1a0a8a3SJed Brown 
1634f746d493SDmitry Karpeev .keywords: PC, GASM, set, local, subdomains, additive Schwarz, block Jacobi
1635b1a0a8a3SJed Brown 
1636f746d493SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(),
1637f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D(), PCGASMSetLocalSubdomains(), PCGASMGetLocalSubdomains()
1638b1a0a8a3SJed Brown @*/
16397087cfbeSBarry Smith PetscErrorCode  PCGASMGetLocalSubmatrices(PC pc,PetscInt *n,Mat *mat[])
1640b1a0a8a3SJed Brown {
1641f746d493SDmitry Karpeev   PC_GASM         *osm;
1642b1a0a8a3SJed Brown   PetscErrorCode ierr;
1643b1a0a8a3SJed Brown   PetscBool      match;
1644b1a0a8a3SJed Brown 
1645b1a0a8a3SJed Brown   PetscFunctionBegin;
1646b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1647b1a0a8a3SJed Brown   PetscValidIntPointer(n,2);
1648b1a0a8a3SJed Brown   if (mat) PetscValidPointer(mat,3);
1649b1a0a8a3SJed Brown   if (!pc->setupcalled) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must call after KSPSetUP() or PCSetUp().");
1650f746d493SDmitry Karpeev   ierr = PetscTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr);
1651b1a0a8a3SJed Brown   if (!match) {
1652b1a0a8a3SJed Brown     if (n)   *n   = 0;
1653b1a0a8a3SJed Brown     if (mat) *mat = PETSC_NULL;
1654b1a0a8a3SJed Brown   } else {
1655f746d493SDmitry Karpeev     osm = (PC_GASM*)pc->data;
1656ab3e923fSDmitry Karpeev     if (n)   *n   = osm->n;
1657b1a0a8a3SJed Brown     if (mat) *mat = osm->pmat;
1658b1a0a8a3SJed Brown   }
1659b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1660b1a0a8a3SJed Brown }
1661