xref: /petsc/src/ksp/pc/impls/asm/asm.c (revision 910cf402deff9031709e244fa7c1e84ef2fdc3da)
14b9ad928SBarry Smith /*
24b9ad928SBarry Smith   This file defines an additive Schwarz preconditioner for any Mat implementation.
34b9ad928SBarry Smith 
44b9ad928SBarry Smith   Note that each processor may have any number of subdomains. But in order to
54b9ad928SBarry Smith   deal easily with the VecScatter(), we treat each processor as if it has the
64b9ad928SBarry Smith   same number of subdomains.
74b9ad928SBarry Smith 
84b9ad928SBarry Smith        n - total number of true subdomains on all processors
94b9ad928SBarry Smith        n_local_true - actual number of subdomains on this processor
104b9ad928SBarry Smith        n_local = maximum over all processors of n_local_true
114b9ad928SBarry Smith */
124b9ad928SBarry Smith 
13*910cf402Sprj- #include <../src/ksp/pc/impls/asm/asm.h> /*I "petscpc.h" I*/
144b9ad928SBarry Smith 
156849ba73SBarry Smith static PetscErrorCode PCView_ASM(PC pc,PetscViewer viewer)
164b9ad928SBarry Smith {
1792bb6962SLisandro Dalcin   PC_ASM         *osm = (PC_ASM*)pc->data;
186849ba73SBarry Smith   PetscErrorCode ierr;
1913f74950SBarry Smith   PetscMPIInt    rank;
204d219a6aSLisandro Dalcin   PetscInt       i,bsz;
21ace3abfcSBarry Smith   PetscBool      iascii,isstring;
224b9ad928SBarry Smith   PetscViewer    sviewer;
234b9ad928SBarry Smith 
244b9ad928SBarry Smith   PetscFunctionBegin;
25251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
26251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
2732077d6dSBarry Smith   if (iascii) {
283d03bb48SJed Brown     char overlaps[256] = "user-defined overlap",blocks[256] = "total subdomain blocks not yet set";
298caf3d72SBarry Smith     if (osm->overlap >= 0) {ierr = PetscSNPrintf(overlaps,sizeof(overlaps),"amount of overlap = %D",osm->overlap);CHKERRQ(ierr);}
308caf3d72SBarry Smith     if (osm->n > 0) {ierr = PetscSNPrintf(blocks,sizeof(blocks),"total subdomain blocks = %D",osm->n);CHKERRQ(ierr);}
31efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  %s, %s\n",blocks,overlaps);CHKERRQ(ierr);
32efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  restriction/interpolation type - %s\n",PCASMTypes[osm->type]);CHKERRQ(ierr);
33e64f0791SPatrick Sanan     if (osm->dm_subdomains) {ierr = PetscViewerASCIIPrintf(viewer,"  Additive Schwarz: using DM to define subdomains\n");CHKERRQ(ierr);}
3412cd4985SMatthew G. Knepley     if (osm->loctype != PC_COMPOSITE_ADDITIVE) {ierr = PetscViewerASCIIPrintf(viewer,"  Additive Schwarz: local solve composition type - %s\n",PCCompositeTypes[osm->loctype]);CHKERRQ(ierr);}
35ce94432eSBarry Smith     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr);
3692bb6962SLisandro Dalcin     if (osm->same_local_solves) {
374d219a6aSLisandro Dalcin       if (osm->ksp) {
384b9ad928SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"  Local solve is same for all blocks, in the following KSP and PC objects:\n");CHKERRQ(ierr);
39c5e4d11fSDmitry Karpeev         ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
404d219a6aSLisandro Dalcin         if (!rank) {
414b9ad928SBarry Smith           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
4292bb6962SLisandro Dalcin           ierr = KSPView(osm->ksp[0],sviewer);CHKERRQ(ierr);
434b9ad928SBarry Smith           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
444b9ad928SBarry Smith         }
45c5e4d11fSDmitry Karpeev         ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
464d219a6aSLisandro Dalcin       }
474b9ad928SBarry Smith     } else {
48c5e4d11fSDmitry Karpeev       ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
494d219a6aSLisandro Dalcin       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"  [%d] number of local blocks = %D\n",(int)rank,osm->n_local_true);CHKERRQ(ierr);
504d219a6aSLisandro Dalcin       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
514b9ad928SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  Local solve info for each block is in the following KSP and PC objects:\n");CHKERRQ(ierr);
524b9ad928SBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
534d219a6aSLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr);
54c5e4d11fSDmitry Karpeev       ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
55dd2fa690SBarry Smith       for (i=0; i<osm->n_local_true; i++) {
564d219a6aSLisandro Dalcin         ierr = ISGetLocalSize(osm->is[i],&bsz);CHKERRQ(ierr);
574d219a6aSLisandro Dalcin         ierr = PetscViewerASCIISynchronizedPrintf(sviewer,"[%d] local block number %D, size = %D\n",(int)rank,i,bsz);CHKERRQ(ierr);
5892bb6962SLisandro Dalcin         ierr = KSPView(osm->ksp[i],sviewer);CHKERRQ(ierr);
594d219a6aSLisandro Dalcin         ierr = PetscViewerASCIISynchronizedPrintf(sviewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr);
604b9ad928SBarry Smith       }
61c5e4d11fSDmitry Karpeev       ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
624b9ad928SBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
634b9ad928SBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
64c5e4d11fSDmitry Karpeev       ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
654b9ad928SBarry Smith     }
664b9ad928SBarry Smith   } else if (isstring) {
674d219a6aSLisandro Dalcin     ierr = PetscViewerStringSPrintf(viewer," blocks=%D, overlap=%D, type=%s",osm->n,osm->overlap,PCASMTypes[osm->type]);CHKERRQ(ierr);
68c5e4d11fSDmitry Karpeev     ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
6992bb6962SLisandro Dalcin     if (osm->ksp) {ierr = KSPView(osm->ksp[0],sviewer);CHKERRQ(ierr);}
70c5e4d11fSDmitry Karpeev     ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
714b9ad928SBarry Smith   }
724b9ad928SBarry Smith   PetscFunctionReturn(0);
734b9ad928SBarry Smith }
744b9ad928SBarry Smith 
7592bb6962SLisandro Dalcin static PetscErrorCode PCASMPrintSubdomains(PC pc)
7692bb6962SLisandro Dalcin {
7792bb6962SLisandro Dalcin   PC_ASM         *osm = (PC_ASM*)pc->data;
7892bb6962SLisandro Dalcin   const char     *prefix;
7992bb6962SLisandro Dalcin   char           fname[PETSC_MAX_PATH_LEN+1];
80643535aeSDmitry Karpeev   PetscViewer    viewer, sviewer;
81643535aeSDmitry Karpeev   char           *s;
8292bb6962SLisandro Dalcin   PetscInt       i,j,nidx;
8392bb6962SLisandro Dalcin   const PetscInt *idx;
84643535aeSDmitry Karpeev   PetscMPIInt    rank, size;
8592bb6962SLisandro Dalcin   PetscErrorCode ierr;
86846783a0SBarry Smith 
8792bb6962SLisandro Dalcin   PetscFunctionBegin;
88ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size);CHKERRQ(ierr);
89ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank);CHKERRQ(ierr);
9092bb6962SLisandro Dalcin   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
91c5929fdfSBarry Smith   ierr = PetscOptionsGetString(NULL,prefix,"-pc_asm_print_subdomains",fname,PETSC_MAX_PATH_LEN,NULL);CHKERRQ(ierr);
9292bb6962SLisandro Dalcin   if (fname[0] == 0) { ierr = PetscStrcpy(fname,"stdout");CHKERRQ(ierr); };
93ce94432eSBarry Smith   ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)pc),fname,&viewer);CHKERRQ(ierr);
94643535aeSDmitry Karpeev   for (i=0; i<osm->n_local; i++) {
95643535aeSDmitry Karpeev     if (i < osm->n_local_true) {
9692bb6962SLisandro Dalcin       ierr = ISGetLocalSize(osm->is[i],&nidx);CHKERRQ(ierr);
9792bb6962SLisandro Dalcin       ierr = ISGetIndices(osm->is[i],&idx);CHKERRQ(ierr);
98643535aeSDmitry Karpeev       /* Print to a string viewer; no more than 15 characters per index plus 512 char for the header.*/
9936a9e3b9SBarry Smith #define len  16*(nidx+1)+512
10036a9e3b9SBarry Smith       ierr = PetscMalloc1(len,&s);CHKERRQ(ierr);
10136a9e3b9SBarry Smith       ierr = PetscViewerStringOpen(PETSC_COMM_SELF, s, len, &sviewer);CHKERRQ(ierr);
10236a9e3b9SBarry Smith #undef len
103643535aeSDmitry Karpeev       ierr = PetscViewerStringSPrintf(sviewer, "[%D:%D] Subdomain %D with overlap:\n", rank, size, i);CHKERRQ(ierr);
10492bb6962SLisandro Dalcin       for (j=0; j<nidx; j++) {
105643535aeSDmitry Karpeev         ierr = PetscViewerStringSPrintf(sviewer,"%D ",idx[j]);CHKERRQ(ierr);
10692bb6962SLisandro Dalcin       }
10792bb6962SLisandro Dalcin       ierr = ISRestoreIndices(osm->is[i],&idx);CHKERRQ(ierr);
108643535aeSDmitry Karpeev       ierr = PetscViewerStringSPrintf(sviewer,"\n");CHKERRQ(ierr);
109643535aeSDmitry Karpeev       ierr = PetscViewerDestroy(&sviewer);CHKERRQ(ierr);
110c5e4d11fSDmitry Karpeev       ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
111643535aeSDmitry Karpeev       ierr = PetscViewerASCIISynchronizedPrintf(viewer, s);CHKERRQ(ierr);
112643535aeSDmitry Karpeev       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
113c5e4d11fSDmitry Karpeev       ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
114643535aeSDmitry Karpeev       ierr = PetscFree(s);CHKERRQ(ierr);
1152b691e39SMatthew Knepley       if (osm->is_local) {
116643535aeSDmitry Karpeev         /* Print to a string viewer; no more than 15 characters per index plus 512 char for the header.*/
11736a9e3b9SBarry Smith #define len  16*(nidx+1)+512
11836a9e3b9SBarry Smith         ierr = PetscMalloc1(len, &s);CHKERRQ(ierr);
11936a9e3b9SBarry Smith         ierr = PetscViewerStringOpen(PETSC_COMM_SELF, s, len, &sviewer);CHKERRQ(ierr);
12036a9e3b9SBarry Smith #undef len
12109d011a0SDmitry Karpeev         ierr = PetscViewerStringSPrintf(sviewer, "[%D:%D] Subdomain %D without overlap:\n", rank, size, i);CHKERRQ(ierr);
1222b691e39SMatthew Knepley         ierr = ISGetLocalSize(osm->is_local[i],&nidx);CHKERRQ(ierr);
1232b691e39SMatthew Knepley         ierr = ISGetIndices(osm->is_local[i],&idx);CHKERRQ(ierr);
1242b691e39SMatthew Knepley         for (j=0; j<nidx; j++) {
12509d011a0SDmitry Karpeev           ierr = PetscViewerStringSPrintf(sviewer,"%D ",idx[j]);CHKERRQ(ierr);
1262b691e39SMatthew Knepley         }
1272b691e39SMatthew Knepley         ierr = ISRestoreIndices(osm->is_local[i],&idx);CHKERRQ(ierr);
12809d011a0SDmitry Karpeev         ierr = PetscViewerStringSPrintf(sviewer,"\n");CHKERRQ(ierr);
129643535aeSDmitry Karpeev         ierr = PetscViewerDestroy(&sviewer);CHKERRQ(ierr);
130c5e4d11fSDmitry Karpeev         ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
131643535aeSDmitry Karpeev         ierr = PetscViewerASCIISynchronizedPrintf(viewer, s);CHKERRQ(ierr);
132643535aeSDmitry Karpeev         ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
133c5e4d11fSDmitry Karpeev         ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
134643535aeSDmitry Karpeev         ierr = PetscFree(s);CHKERRQ(ierr);
135643535aeSDmitry Karpeev       }
1362fa5cd67SKarl Rupp     } else {
137643535aeSDmitry Karpeev       /* Participate in collective viewer calls. */
138c5e4d11fSDmitry Karpeev       ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
139643535aeSDmitry Karpeev       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
140c5e4d11fSDmitry Karpeev       ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
141643535aeSDmitry Karpeev       /* Assume either all ranks have is_local or none do. */
142643535aeSDmitry Karpeev       if (osm->is_local) {
143c5e4d11fSDmitry Karpeev         ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
144643535aeSDmitry Karpeev         ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
145c5e4d11fSDmitry Karpeev         ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
146643535aeSDmitry Karpeev       }
147fdc48646SDmitry Karpeev     }
14892bb6962SLisandro Dalcin   }
14992bb6962SLisandro Dalcin   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
150fcfd50ebSBarry Smith   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
15192bb6962SLisandro Dalcin   PetscFunctionReturn(0);
15292bb6962SLisandro Dalcin }
15392bb6962SLisandro Dalcin 
1546849ba73SBarry Smith static PetscErrorCode PCSetUp_ASM(PC pc)
1554b9ad928SBarry Smith {
1564b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
1576849ba73SBarry Smith   PetscErrorCode ierr;
15857501b6eSBarry Smith   PetscBool      flg;
15987e86712SBarry Smith   PetscInt       i,m,m_local;
1604b9ad928SBarry Smith   MatReuse       scall = MAT_REUSE_MATRIX;
1614b9ad928SBarry Smith   IS             isl;
1624b9ad928SBarry Smith   KSP            ksp;
1634b9ad928SBarry Smith   PC             subpc;
1642dcb1b2aSMatthew Knepley   const char     *prefix,*pprefix;
16523ce1328SBarry Smith   Vec            vec;
1660298fd71SBarry Smith   DM             *domain_dm = NULL;
1674b9ad928SBarry Smith 
1684b9ad928SBarry Smith   PetscFunctionBegin;
1694b9ad928SBarry Smith   if (!pc->setupcalled) {
170265a2120SBarry Smith     PetscInt m;
17192bb6962SLisandro Dalcin 
1722b837212SDmitry Karpeev     /* Note: if subdomains have been set either via PCASMSetTotalSubdomains() or via PCASMSetLocalSubdomains(), osm->n_local_true will not be PETSC_DECIDE */
1732b837212SDmitry Karpeev     if (osm->n_local_true == PETSC_DECIDE) {
17469ca1f37SDmitry Karpeev       /* no subdomains given */
17565db9045SDmitry Karpeev       /* try pc->dm first, if allowed */
176d709ea83SDmitry Karpeev       if (osm->dm_subdomains && pc->dm) {
177feb221f8SDmitry Karpeev         PetscInt  num_domains, d;
178feb221f8SDmitry Karpeev         char      **domain_names;
1798d4ac253SDmitry Karpeev         IS        *inner_domain_is, *outer_domain_is;
1808d4ac253SDmitry Karpeev         ierr = DMCreateDomainDecomposition(pc->dm, &num_domains, &domain_names, &inner_domain_is, &outer_domain_is, &domain_dm);CHKERRQ(ierr);
181704f0395SPatrick Sanan         osm->overlap = -1; /* We do not want to increase the overlap of the IS.
182704f0395SPatrick Sanan                               A future improvement of this code might allow one to use
183704f0395SPatrick Sanan                               DM-defined subdomains and also increase the overlap,
184704f0395SPatrick Sanan                               but that is not currently supported */
18569ca1f37SDmitry Karpeev         if (num_domains) {
1868d4ac253SDmitry Karpeev           ierr = PCASMSetLocalSubdomains(pc, num_domains, outer_domain_is, inner_domain_is);CHKERRQ(ierr);
18769ca1f37SDmitry Karpeev         }
188feb221f8SDmitry Karpeev         for (d = 0; d < num_domains; ++d) {
189a35b7b57SDmitry Karpeev           if (domain_names)    {ierr = PetscFree(domain_names[d]);CHKERRQ(ierr);}
190a35b7b57SDmitry Karpeev           if (inner_domain_is) {ierr = ISDestroy(&inner_domain_is[d]);CHKERRQ(ierr);}
191a35b7b57SDmitry Karpeev           if (outer_domain_is) {ierr = ISDestroy(&outer_domain_is[d]);CHKERRQ(ierr);}
192feb221f8SDmitry Karpeev         }
193feb221f8SDmitry Karpeev         ierr = PetscFree(domain_names);CHKERRQ(ierr);
1948d4ac253SDmitry Karpeev         ierr = PetscFree(inner_domain_is);CHKERRQ(ierr);
1958d4ac253SDmitry Karpeev         ierr = PetscFree(outer_domain_is);CHKERRQ(ierr);
196feb221f8SDmitry Karpeev       }
1972b837212SDmitry Karpeev       if (osm->n_local_true == PETSC_DECIDE) {
19869ca1f37SDmitry Karpeev         /* still no subdomains; use one subdomain per processor */
1992b837212SDmitry Karpeev         osm->n_local_true = 1;
200feb221f8SDmitry Karpeev       }
2012b837212SDmitry Karpeev     }
2022b837212SDmitry Karpeev     { /* determine the global and max number of subdomains */
2036ac3741eSJed Brown       struct {PetscInt max,sum;} inwork,outwork;
204c10200c1SHong Zhang       PetscMPIInt size;
205c10200c1SHong Zhang 
2066ac3741eSJed Brown       inwork.max   = osm->n_local_true;
2076ac3741eSJed Brown       inwork.sum   = osm->n_local_true;
208367daffbSBarry Smith       ierr         = MPIU_Allreduce(&inwork,&outwork,1,MPIU_2INT,MPIU_MAXSUM_OP,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
2096ac3741eSJed Brown       osm->n_local = outwork.max;
2106ac3741eSJed Brown       osm->n       = outwork.sum;
211c10200c1SHong Zhang 
212c10200c1SHong Zhang       ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr);
213c10200c1SHong Zhang       if (outwork.max == 1 && outwork.sum == size) {
2147dae84e0SHong Zhang         /* osm->n_local_true = 1 on all processes, set this option may enable use of optimized MatCreateSubMatrices() implementation */
215c10200c1SHong Zhang         ierr = MatSetOption(pc->pmat,MAT_SUBMAT_SINGLEIS,PETSC_TRUE);CHKERRQ(ierr);
216c10200c1SHong Zhang       }
2174b9ad928SBarry Smith     }
21878904715SLisandro Dalcin     if (!osm->is) { /* create the index sets */
21978904715SLisandro Dalcin       ierr = PCASMCreateSubdomains(pc->pmat,osm->n_local_true,&osm->is);CHKERRQ(ierr);
2204b9ad928SBarry Smith     }
221f5234e35SJed Brown     if (osm->n_local_true > 1 && !osm->is_local) {
222785e854fSJed Brown       ierr = PetscMalloc1(osm->n_local_true,&osm->is_local);CHKERRQ(ierr);
223f5234e35SJed Brown       for (i=0; i<osm->n_local_true; i++) {
224f5234e35SJed Brown         if (osm->overlap > 0) { /* With positive overlap, osm->is[i] will be modified */
225f5234e35SJed Brown           ierr = ISDuplicate(osm->is[i],&osm->is_local[i]);CHKERRQ(ierr);
226f5234e35SJed Brown           ierr = ISCopy(osm->is[i],osm->is_local[i]);CHKERRQ(ierr);
227f5234e35SJed Brown         } else {
228f5234e35SJed Brown           ierr             = PetscObjectReference((PetscObject)osm->is[i]);CHKERRQ(ierr);
229f5234e35SJed Brown           osm->is_local[i] = osm->is[i];
230f5234e35SJed Brown         }
231f5234e35SJed Brown       }
232f5234e35SJed Brown     }
23392bb6962SLisandro Dalcin     ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
23490d69ab7SBarry Smith     flg  = PETSC_FALSE;
235c5929fdfSBarry Smith     ierr = PetscOptionsGetBool(NULL,prefix,"-pc_asm_print_subdomains",&flg,NULL);CHKERRQ(ierr);
23692bb6962SLisandro Dalcin     if (flg) { ierr = PCASMPrintSubdomains(pc);CHKERRQ(ierr); }
2374b9ad928SBarry Smith 
2383d03bb48SJed Brown     if (osm->overlap > 0) {
2394b9ad928SBarry Smith       /* Extend the "overlapping" regions by a number of steps */
24092bb6962SLisandro Dalcin       ierr = MatIncreaseOverlap(pc->pmat,osm->n_local_true,osm->is,osm->overlap);CHKERRQ(ierr);
2413d03bb48SJed Brown     }
2426ed231c7SMatthew Knepley     if (osm->sort_indices) {
24392bb6962SLisandro Dalcin       for (i=0; i<osm->n_local_true; i++) {
2444b9ad928SBarry Smith         ierr = ISSort(osm->is[i]);CHKERRQ(ierr);
2452b691e39SMatthew Knepley         if (osm->is_local) {
2462b691e39SMatthew Knepley           ierr = ISSort(osm->is_local[i]);CHKERRQ(ierr);
2472b691e39SMatthew Knepley         }
2484b9ad928SBarry Smith       }
2496ed231c7SMatthew Knepley     }
2504b9ad928SBarry Smith 
251f6991133SBarry Smith     if (!osm->ksp) {
25278904715SLisandro Dalcin       /* Create the local solvers */
253785e854fSJed Brown       ierr = PetscMalloc1(osm->n_local_true,&osm->ksp);CHKERRQ(ierr);
254feb221f8SDmitry Karpeev       if (domain_dm) {
255feb221f8SDmitry Karpeev         ierr = PetscInfo(pc,"Setting up ASM subproblems using the embedded DM\n");CHKERRQ(ierr);
256feb221f8SDmitry Karpeev       }
25792bb6962SLisandro Dalcin       for (i=0; i<osm->n_local_true; i++) {
2584b9ad928SBarry Smith         ierr = KSPCreate(PETSC_COMM_SELF,&ksp);CHKERRQ(ierr);
259422a814eSBarry Smith         ierr = KSPSetErrorIfNotConverged(ksp,pc->erroriffailure);CHKERRQ(ierr);
2603bb1ff40SBarry Smith         ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ksp);CHKERRQ(ierr);
26192bb6962SLisandro Dalcin         ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);CHKERRQ(ierr);
2624b9ad928SBarry Smith         ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr);
2634b9ad928SBarry Smith         ierr = KSPGetPC(ksp,&subpc);CHKERRQ(ierr);
2644b9ad928SBarry Smith         ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
2654b9ad928SBarry Smith         ierr = KSPSetOptionsPrefix(ksp,prefix);CHKERRQ(ierr);
2664b9ad928SBarry Smith         ierr = KSPAppendOptionsPrefix(ksp,"sub_");CHKERRQ(ierr);
267feb221f8SDmitry Karpeev         if (domain_dm) {
268feb221f8SDmitry Karpeev           ierr = KSPSetDM(ksp, domain_dm[i]);CHKERRQ(ierr);
269feb221f8SDmitry Karpeev           ierr = KSPSetDMActive(ksp, PETSC_FALSE);CHKERRQ(ierr);
270feb221f8SDmitry Karpeev           ierr = DMDestroy(&domain_dm[i]);CHKERRQ(ierr);
271feb221f8SDmitry Karpeev         }
2724b9ad928SBarry Smith         osm->ksp[i] = ksp;
2734b9ad928SBarry Smith       }
274feb221f8SDmitry Karpeev       if (domain_dm) {
275feb221f8SDmitry Karpeev         ierr = PetscFree(domain_dm);CHKERRQ(ierr);
276feb221f8SDmitry Karpeev       }
277f6991133SBarry Smith     }
2781dd8081eSeaulisa 
279fb745f2cSMatthew G. Knepley     ierr = ISConcatenate(PETSC_COMM_SELF, osm->n_local_true, osm->is, &osm->lis);CHKERRQ(ierr);
280fb745f2cSMatthew G. Knepley     ierr = ISSortRemoveDups(osm->lis);CHKERRQ(ierr);
281fb745f2cSMatthew G. Knepley     ierr = ISGetLocalSize(osm->lis, &m);CHKERRQ(ierr);
282fb745f2cSMatthew G. Knepley     ierr = VecCreateSeq(PETSC_COMM_SELF, m, &osm->lx);CHKERRQ(ierr);
283fb745f2cSMatthew G. Knepley     ierr = VecDuplicate(osm->lx, &osm->ly);CHKERRQ(ierr);
2841dd8081eSeaulisa 
2854b9ad928SBarry Smith     scall = MAT_INITIAL_MATRIX;
2864b9ad928SBarry Smith   } else {
2874b9ad928SBarry Smith     /*
2884b9ad928SBarry Smith        Destroy the blocks from the previous iteration
2894b9ad928SBarry Smith     */
290e09e08ccSBarry Smith     if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
2914b9ad928SBarry Smith       ierr = MatDestroyMatrices(osm->n_local_true,&osm->pmat);CHKERRQ(ierr);
2924b9ad928SBarry Smith       scall = MAT_INITIAL_MATRIX;
2934b9ad928SBarry Smith     }
2944b9ad928SBarry Smith   }
2954b9ad928SBarry Smith 
29692bb6962SLisandro Dalcin   /*
29792bb6962SLisandro Dalcin      Extract out the submatrices
29892bb6962SLisandro Dalcin   */
2997dae84e0SHong Zhang   ierr = MatCreateSubMatrices(pc->pmat,osm->n_local_true,osm->is,osm->is,scall,&osm->pmat);CHKERRQ(ierr);
30092bb6962SLisandro Dalcin   if (scall == MAT_INITIAL_MATRIX) {
30192bb6962SLisandro Dalcin     ierr = PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix);CHKERRQ(ierr);
30292bb6962SLisandro Dalcin     for (i=0; i<osm->n_local_true; i++) {
3033bb1ff40SBarry Smith       ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)osm->pmat[i]);CHKERRQ(ierr);
30478904715SLisandro Dalcin       ierr = PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i],pprefix);CHKERRQ(ierr);
30592bb6962SLisandro Dalcin     }
30692bb6962SLisandro Dalcin   }
30780ec0b7dSPatrick Sanan 
30880ec0b7dSPatrick Sanan   /* Convert the types of the submatrices (if needbe) */
30980ec0b7dSPatrick Sanan   if (osm->sub_mat_type) {
31080ec0b7dSPatrick Sanan     for (i=0; i<osm->n_local_true; i++) {
31180ec0b7dSPatrick Sanan       ierr = MatConvert(osm->pmat[i],osm->sub_mat_type,MAT_INPLACE_MATRIX,&(osm->pmat[i]));CHKERRQ(ierr);
31280ec0b7dSPatrick Sanan     }
31380ec0b7dSPatrick Sanan   }
31480ec0b7dSPatrick Sanan 
31580ec0b7dSPatrick Sanan   if (!pc->setupcalled) {
31680ec0b7dSPatrick Sanan     /* Create the local work vectors (from the local matrices) and scatter contexts */
31780ec0b7dSPatrick Sanan     ierr = MatCreateVecs(pc->pmat,&vec,0);CHKERRQ(ierr);
3185b3c0d42Seaulisa 
3191dd8081eSeaulisa     if (osm->is_local && (osm->type == PC_ASM_INTERPOLATE || osm->type == PC_ASM_NONE )) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot use interpolate or none PCASMType if is_local was provided to PCASMSetLocalSubdomains()");
3201dd8081eSeaulisa     if (osm->is_local && osm->type == PC_ASM_RESTRICT && osm->loctype == PC_COMPOSITE_ADDITIVE) {
3211dd8081eSeaulisa       ierr = PetscMalloc1(osm->n_local_true,&osm->lprolongation);CHKERRQ(ierr);
3221dd8081eSeaulisa     }
3231dd8081eSeaulisa     ierr = PetscMalloc1(osm->n_local_true,&osm->lrestriction);CHKERRQ(ierr);
3241dd8081eSeaulisa     ierr = PetscMalloc1(osm->n_local_true,&osm->x);CHKERRQ(ierr);
3251dd8081eSeaulisa     ierr = PetscMalloc1(osm->n_local_true,&osm->y);CHKERRQ(ierr);
326347574c9Seaulisa 
327347574c9Seaulisa     ierr = ISGetLocalSize(osm->lis,&m);CHKERRQ(ierr);
328347574c9Seaulisa     ierr = ISCreateStride(PETSC_COMM_SELF,m,0,1,&isl);CHKERRQ(ierr);
3299448b7f1SJunchao Zhang     ierr = VecScatterCreate(vec,osm->lis,osm->lx,isl,&osm->restriction);CHKERRQ(ierr);
330347574c9Seaulisa     ierr = ISDestroy(&isl);CHKERRQ(ierr);
331347574c9Seaulisa 
332347574c9Seaulisa 
33380ec0b7dSPatrick Sanan     for (i=0; i<osm->n_local_true; ++i) {
3345b3c0d42Seaulisa       ISLocalToGlobalMapping ltog;
3355b3c0d42Seaulisa       IS                     isll;
3365b3c0d42Seaulisa       const PetscInt         *idx_is;
3375b3c0d42Seaulisa       PetscInt               *idx_lis,nout;
3385b3c0d42Seaulisa 
33955817e79SBarry Smith       ierr = ISGetLocalSize(osm->is[i],&m);CHKERRQ(ierr);
34055817e79SBarry Smith       ierr = MatCreateVecs(osm->pmat[i],&osm->x[i],NULL);CHKERRQ(ierr);
34155817e79SBarry Smith       ierr = VecDuplicate(osm->x[i],&osm->y[i]);CHKERRQ(ierr);
34255817e79SBarry Smith 
343b0de9f23SBarry Smith       /* generate a scatter from ly to y[i] picking all the overlapping is[i] entries */
3445b3c0d42Seaulisa       ierr = ISLocalToGlobalMappingCreateIS(osm->lis,&ltog);CHKERRQ(ierr);
3455b3c0d42Seaulisa       ierr = ISGetLocalSize(osm->is[i],&m);CHKERRQ(ierr);
3465b3c0d42Seaulisa       ierr = ISGetIndices(osm->is[i], &idx_is);CHKERRQ(ierr);
3475b3c0d42Seaulisa       ierr = PetscMalloc1(m,&idx_lis);CHKERRQ(ierr);
3485b3c0d42Seaulisa       ierr = ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,m,idx_is,&nout,idx_lis);CHKERRQ(ierr);
3495b3c0d42Seaulisa       if (nout != m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"is not a subset of lis");
3505b3c0d42Seaulisa       ierr = ISRestoreIndices(osm->is[i], &idx_is);CHKERRQ(ierr);
3515b3c0d42Seaulisa       ierr = ISCreateGeneral(PETSC_COMM_SELF,m,idx_lis,PETSC_OWN_POINTER,&isll);CHKERRQ(ierr);
352d8b3b5e3Seaulisa       ierr = ISLocalToGlobalMappingDestroy(&ltog);CHKERRQ(ierr);
3535b3c0d42Seaulisa       ierr = ISCreateStride(PETSC_COMM_SELF,m,0,1,&isl);CHKERRQ(ierr);
3549448b7f1SJunchao Zhang       ierr = VecScatterCreate(osm->ly,isll,osm->y[i],isl,&osm->lrestriction[i]);CHKERRQ(ierr);
3555b3c0d42Seaulisa       ierr = ISDestroy(&isll);CHKERRQ(ierr);
3565b3c0d42Seaulisa       ierr = ISDestroy(&isl);CHKERRQ(ierr);
357*910cf402Sprj-       if (osm->lprolongation) { /* generate a scatter from y[i] to ly picking only the the non-overlapping is_local[i] entries */
358d8b3b5e3Seaulisa         ISLocalToGlobalMapping ltog;
359d8b3b5e3Seaulisa         IS                     isll,isll_local;
360d8b3b5e3Seaulisa         const PetscInt         *idx_local;
361d8b3b5e3Seaulisa         PetscInt               *idx1, *idx2, nout;
362d8b3b5e3Seaulisa 
363d8b3b5e3Seaulisa         ierr = ISGetLocalSize(osm->is_local[i],&m_local);CHKERRQ(ierr);
364d8b3b5e3Seaulisa         ierr = ISGetIndices(osm->is_local[i], &idx_local);CHKERRQ(ierr);
365d8b3b5e3Seaulisa 
366d8b3b5e3Seaulisa         ierr = ISLocalToGlobalMappingCreateIS(osm->is[i],&ltog);CHKERRQ(ierr);
367d8b3b5e3Seaulisa         ierr = PetscMalloc1(m_local,&idx1);CHKERRQ(ierr);
368d8b3b5e3Seaulisa         ierr = ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,m_local,idx_local,&nout,idx1);CHKERRQ(ierr);
369d8b3b5e3Seaulisa         ierr = ISLocalToGlobalMappingDestroy(&ltog);CHKERRQ(ierr);
370d8b3b5e3Seaulisa         if (nout != m_local) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"is_local not a subset of is");
371d8b3b5e3Seaulisa         ierr = ISCreateGeneral(PETSC_COMM_SELF,m_local,idx1,PETSC_OWN_POINTER,&isll);CHKERRQ(ierr);
372d8b3b5e3Seaulisa 
373d8b3b5e3Seaulisa         ierr = ISLocalToGlobalMappingCreateIS(osm->lis,&ltog);CHKERRQ(ierr);
374d8b3b5e3Seaulisa         ierr = PetscMalloc1(m_local,&idx2);CHKERRQ(ierr);
375d8b3b5e3Seaulisa         ierr = ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,m_local,idx_local,&nout,idx2);CHKERRQ(ierr);
376d8b3b5e3Seaulisa         ierr = ISLocalToGlobalMappingDestroy(&ltog);CHKERRQ(ierr);
377d8b3b5e3Seaulisa         if (nout != m_local) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"is_local not a subset of lis");
378d8b3b5e3Seaulisa         ierr = ISCreateGeneral(PETSC_COMM_SELF,m_local,idx2,PETSC_OWN_POINTER,&isll_local);CHKERRQ(ierr);
379d8b3b5e3Seaulisa 
380d8b3b5e3Seaulisa         ierr = ISRestoreIndices(osm->is_local[i], &idx_local);CHKERRQ(ierr);
3819448b7f1SJunchao Zhang         ierr = VecScatterCreate(osm->y[i],isll,osm->ly,isll_local,&osm->lprolongation[i]);CHKERRQ(ierr);
382d8b3b5e3Seaulisa 
383d8b3b5e3Seaulisa         ierr = ISDestroy(&isll);CHKERRQ(ierr);
384d8b3b5e3Seaulisa         ierr = ISDestroy(&isll_local);CHKERRQ(ierr);
385d8b3b5e3Seaulisa       }
38680ec0b7dSPatrick Sanan     }
38780ec0b7dSPatrick Sanan     ierr = VecDestroy(&vec);CHKERRQ(ierr);
38880ec0b7dSPatrick Sanan   }
38980ec0b7dSPatrick Sanan 
390fb745f2cSMatthew G. Knepley   if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) {
391235cc4ceSMatthew G. Knepley     IS      *cis;
392235cc4ceSMatthew G. Knepley     PetscInt c;
393235cc4ceSMatthew G. Knepley 
394235cc4ceSMatthew G. Knepley     ierr = PetscMalloc1(osm->n_local_true, &cis);CHKERRQ(ierr);
395235cc4ceSMatthew G. Knepley     for (c = 0; c < osm->n_local_true; ++c) cis[c] = osm->lis;
3967dae84e0SHong Zhang     ierr = MatCreateSubMatrices(pc->pmat, osm->n_local_true, osm->is, cis, scall, &osm->lmats);CHKERRQ(ierr);
397235cc4ceSMatthew G. Knepley     ierr = PetscFree(cis);CHKERRQ(ierr);
398fb745f2cSMatthew G. Knepley   }
3994b9ad928SBarry Smith 
4004b9ad928SBarry Smith   /* Return control to the user so that the submatrices can be modified (e.g., to apply
4014b9ad928SBarry Smith      different boundary conditions for the submatrices than for the global problem) */
40292bb6962SLisandro Dalcin   ierr = PCModifySubMatrices(pc,osm->n_local_true,osm->is,osm->is,osm->pmat,pc->modifysubmatricesP);CHKERRQ(ierr);
4034b9ad928SBarry Smith 
40492bb6962SLisandro Dalcin   /*
40592bb6962SLisandro Dalcin      Loop over subdomains putting them into local ksp
40692bb6962SLisandro Dalcin   */
40792bb6962SLisandro Dalcin   for (i=0; i<osm->n_local_true; i++) {
40823ee1639SBarry Smith     ierr = KSPSetOperators(osm->ksp[i],osm->pmat[i],osm->pmat[i]);CHKERRQ(ierr);
40992bb6962SLisandro Dalcin     if (!pc->setupcalled) {
410bf108f30SBarry Smith       ierr = KSPSetFromOptions(osm->ksp[i]);CHKERRQ(ierr);
4114b9ad928SBarry Smith     }
41292bb6962SLisandro Dalcin   }
4134b9ad928SBarry Smith   PetscFunctionReturn(0);
4144b9ad928SBarry Smith }
4154b9ad928SBarry Smith 
4166849ba73SBarry Smith static PetscErrorCode PCSetUpOnBlocks_ASM(PC pc)
4174b9ad928SBarry Smith {
4184b9ad928SBarry Smith   PC_ASM             *osm = (PC_ASM*)pc->data;
4196849ba73SBarry Smith   PetscErrorCode     ierr;
42013f74950SBarry Smith   PetscInt           i;
421539c167fSBarry Smith   KSPConvergedReason reason;
4224b9ad928SBarry Smith 
4234b9ad928SBarry Smith   PetscFunctionBegin;
4244b9ad928SBarry Smith   for (i=0; i<osm->n_local_true; i++) {
4254b9ad928SBarry Smith     ierr = KSPSetUp(osm->ksp[i]);CHKERRQ(ierr);
426539c167fSBarry Smith     ierr = KSPGetConvergedReason(osm->ksp[i],&reason);CHKERRQ(ierr);
427c0decd05SBarry Smith     if (reason == KSP_DIVERGED_PC_FAILED) {
428261222b2SHong Zhang       pc->failedreason = PC_SUBPC_ERROR;
429e0eafd54SHong Zhang     }
4304b9ad928SBarry Smith   }
4314b9ad928SBarry Smith   PetscFunctionReturn(0);
4324b9ad928SBarry Smith }
4334b9ad928SBarry Smith 
4346849ba73SBarry Smith static PetscErrorCode PCApply_ASM(PC pc,Vec x,Vec y)
4354b9ad928SBarry Smith {
4364b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
4376849ba73SBarry Smith   PetscErrorCode ierr;
4381dd8081eSeaulisa   PetscInt       i,n_local_true = osm->n_local_true;
4394b9ad928SBarry Smith   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;
4404b9ad928SBarry Smith 
4414b9ad928SBarry Smith   PetscFunctionBegin;
4424b9ad928SBarry Smith   /*
4434b9ad928SBarry Smith      Support for limiting the restriction or interpolation to only local
4444b9ad928SBarry Smith      subdomain values (leaving the other values 0).
4454b9ad928SBarry Smith   */
4464b9ad928SBarry Smith   if (!(osm->type & PC_ASM_RESTRICT)) {
4474b9ad928SBarry Smith     forward = SCATTER_FORWARD_LOCAL;
4484b9ad928SBarry Smith     /* have to zero the work RHS since scatter may leave some slots empty */
4491dd8081eSeaulisa     ierr = VecSet(osm->lx, 0.0);CHKERRQ(ierr);
4504b9ad928SBarry Smith   }
451347574c9Seaulisa   if (!(osm->type & PC_ASM_INTERPOLATE)) {
452347574c9Seaulisa     reverse = SCATTER_REVERSE_LOCAL;
453347574c9Seaulisa   }
4544b9ad928SBarry Smith 
455347574c9Seaulisa   if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE || osm->loctype == PC_COMPOSITE_ADDITIVE) {
456b0de9f23SBarry Smith     /* zero the global and the local solutions */
457*910cf402Sprj-     ierr = VecSet(y, 0.0);CHKERRQ(ierr);
4585b3c0d42Seaulisa     ierr = VecSet(osm->ly, 0.0);CHKERRQ(ierr);
459347574c9Seaulisa 
460b0de9f23SBarry Smith     /* Copy the global RHS to local RHS including the ghost nodes */
4611dd8081eSeaulisa     ierr = VecScatterBegin(osm->restriction, x, osm->lx, INSERT_VALUES, forward);CHKERRQ(ierr);
4621dd8081eSeaulisa     ierr = VecScatterEnd(osm->restriction, x, osm->lx, INSERT_VALUES, forward);CHKERRQ(ierr);
463347574c9Seaulisa 
464b0de9f23SBarry Smith     /* Restrict local RHS to the overlapping 0-block RHS */
4651dd8081eSeaulisa     ierr = VecScatterBegin(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward);CHKERRQ(ierr);
4661dd8081eSeaulisa     ierr = VecScatterEnd(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward);CHKERRQ(ierr);
467d8b3b5e3Seaulisa 
46812cd4985SMatthew G. Knepley     /* do the local solves */
46912cd4985SMatthew G. Knepley     for (i = 0; i < n_local_true; ++i) {
470347574c9Seaulisa 
471b0de9f23SBarry Smith       /* solve the overlapping i-block */
472fa2bb9feSLisandro Dalcin       ierr = PetscLogEventBegin(PC_ApplyOnBlocks,osm->ksp[i],osm->x[i],osm->y[i],0);CHKERRQ(ierr);
47312cd4985SMatthew G. Knepley       ierr = KSPSolve(osm->ksp[i], osm->x[i], osm->y[i]);CHKERRQ(ierr);
474c0decd05SBarry Smith       ierr = KSPCheckSolve(osm->ksp[i],pc,osm->y[i]);CHKERRQ(ierr);
475fa2bb9feSLisandro Dalcin       ierr = PetscLogEventEnd(PC_ApplyOnBlocks,osm->ksp[i],osm->x[i],osm->y[i],0);CHKERRQ(ierr);
476d8b3b5e3Seaulisa 
477*910cf402Sprj-       if (osm->lprolongation) { /* interpolate the non-overlapping i-block solution to the local solution (only for restrictive additive) */
4781dd8081eSeaulisa         ierr = VecScatterBegin(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward);CHKERRQ(ierr);
4791dd8081eSeaulisa         ierr = VecScatterEnd(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward);CHKERRQ(ierr);
480d8b3b5e3Seaulisa       }
481*910cf402Sprj-       else { /* interpolate the overlapping i-block solution to the local solution */
4821dd8081eSeaulisa         ierr = VecScatterBegin(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr);
4831dd8081eSeaulisa         ierr = VecScatterEnd(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr);
484d8b3b5e3Seaulisa       }
485347574c9Seaulisa 
486347574c9Seaulisa       if (i < n_local_true-1) {
487b0de9f23SBarry Smith         /* Restrict local RHS to the overlapping (i+1)-block RHS */
4881dd8081eSeaulisa         ierr = VecScatterBegin(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward);CHKERRQ(ierr);
4891dd8081eSeaulisa         ierr = VecScatterEnd(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward);CHKERRQ(ierr);
490347574c9Seaulisa 
491347574c9Seaulisa         if ( osm->loctype == PC_COMPOSITE_MULTIPLICATIVE){
4928966356dSPierre Jolivet           /* update the overlapping (i+1)-block RHS using the current local solution */
493347574c9Seaulisa           ierr = MatMult(osm->lmats[i+1], osm->ly, osm->y[i+1]);CHKERRQ(ierr);
494347574c9Seaulisa           ierr = VecAXPBY(osm->x[i+1],-1.,1., osm->y[i+1]); CHKERRQ(ierr);
4957c3d802fSMatthew G. Knepley         }
49612cd4985SMatthew G. Knepley       }
49712cd4985SMatthew G. Knepley     }
498b0de9f23SBarry Smith     /* Add the local solution to the global solution including the ghost nodes */
4991dd8081eSeaulisa     ierr = VecScatterBegin(osm->restriction, osm->ly, y, ADD_VALUES, reverse);CHKERRQ(ierr);
5001dd8081eSeaulisa     ierr = VecScatterEnd(osm->restriction, osm->ly, y, ADD_VALUES, reverse);CHKERRQ(ierr);
501347574c9Seaulisa   } else {
502347574c9Seaulisa     SETERRQ1(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Invalid local composition type: %s", PCCompositeTypes[osm->loctype]);
50312cd4985SMatthew G. Knepley   }
5044b9ad928SBarry Smith   PetscFunctionReturn(0);
5054b9ad928SBarry Smith }
5064b9ad928SBarry Smith 
5076849ba73SBarry Smith static PetscErrorCode PCApplyTranspose_ASM(PC pc,Vec x,Vec y)
5084b9ad928SBarry Smith {
5094b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
5106849ba73SBarry Smith   PetscErrorCode ierr;
5111dd8081eSeaulisa   PetscInt       i,n_local_true = osm->n_local_true;
5124b9ad928SBarry Smith   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;
5134b9ad928SBarry Smith 
5144b9ad928SBarry Smith   PetscFunctionBegin;
5154b9ad928SBarry Smith   /*
5164b9ad928SBarry Smith      Support for limiting the restriction or interpolation to only local
5174b9ad928SBarry Smith      subdomain values (leaving the other values 0).
5184b9ad928SBarry Smith 
5194b9ad928SBarry Smith      Note: these are reversed from the PCApply_ASM() because we are applying the
5204b9ad928SBarry Smith      transpose of the three terms
5214b9ad928SBarry Smith   */
522d8b3b5e3Seaulisa 
5234b9ad928SBarry Smith   if (!(osm->type & PC_ASM_INTERPOLATE)) {
5244b9ad928SBarry Smith     forward = SCATTER_FORWARD_LOCAL;
5254b9ad928SBarry Smith     /* have to zero the work RHS since scatter may leave some slots empty */
5261dd8081eSeaulisa     ierr = VecSet(osm->lx, 0.0);CHKERRQ(ierr);
5274b9ad928SBarry Smith   }
5282fa5cd67SKarl Rupp   if (!(osm->type & PC_ASM_RESTRICT)) reverse = SCATTER_REVERSE_LOCAL;
5294b9ad928SBarry Smith 
530b0de9f23SBarry Smith   /* zero the global and the local solutions */
531*910cf402Sprj-   ierr = VecSet(y, 0.0);CHKERRQ(ierr);
532d8b3b5e3Seaulisa   ierr = VecSet(osm->ly, 0.0);CHKERRQ(ierr);
533d8b3b5e3Seaulisa 
534b0de9f23SBarry Smith   /* Copy the global RHS to local RHS including the ghost nodes */
5351dd8081eSeaulisa   ierr = VecScatterBegin(osm->restriction, x, osm->lx, INSERT_VALUES, forward);CHKERRQ(ierr);
5361dd8081eSeaulisa   ierr = VecScatterEnd(osm->restriction, x, osm->lx, INSERT_VALUES, forward);CHKERRQ(ierr);
537d8b3b5e3Seaulisa 
538b0de9f23SBarry Smith   /* Restrict local RHS to the overlapping 0-block RHS */
5391dd8081eSeaulisa   ierr = VecScatterBegin(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward);CHKERRQ(ierr);
5401dd8081eSeaulisa   ierr = VecScatterEnd(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward);CHKERRQ(ierr);
541d8b3b5e3Seaulisa 
5424b9ad928SBarry Smith   /* do the local solves */
543d8b3b5e3Seaulisa   for (i = 0; i < n_local_true; ++i) {
544d8b3b5e3Seaulisa 
545b0de9f23SBarry Smith     /* solve the overlapping i-block */
546fa2bb9feSLisandro Dalcin     ierr = PetscLogEventBegin(PC_ApplyOnBlocks,osm->ksp[i],osm->x[i],osm->y[i],0);CHKERRQ(ierr);
547e1bcd54cSBarry Smith     ierr = KSPSolveTranspose(osm->ksp[i], osm->x[i], osm->y[i]);CHKERRQ(ierr);
548c0decd05SBarry Smith     ierr = KSPCheckSolve(osm->ksp[i],pc,osm->y[i]);CHKERRQ(ierr);
549fa2bb9feSLisandro Dalcin     ierr = PetscLogEventEnd(PC_ApplyOnBlocks,osm->ksp[i],osm->x[i],osm->y[i],0);CHKERRQ(ierr);
550d8b3b5e3Seaulisa 
551*910cf402Sprj-     if (osm->lprolongation) { /* interpolate the non-overlapping i-block solution to the local solution */
5521dd8081eSeaulisa       ierr = VecScatterBegin(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward);CHKERRQ(ierr);
5531dd8081eSeaulisa       ierr = VecScatterEnd(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward);CHKERRQ(ierr);
554*910cf402Sprj-     } else { /* interpolate the overlapping i-block solution to the local solution */
5551dd8081eSeaulisa       ierr = VecScatterBegin(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr);
5561dd8081eSeaulisa       ierr = VecScatterEnd(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr);
5574b9ad928SBarry Smith     }
558d8b3b5e3Seaulisa 
559d8b3b5e3Seaulisa     if (i < n_local_true-1) {
560b0de9f23SBarry Smith       /* Restrict local RHS to the overlapping (i+1)-block RHS */
5611dd8081eSeaulisa       ierr = VecScatterBegin(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward);CHKERRQ(ierr);
5621dd8081eSeaulisa       ierr = VecScatterEnd(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward);CHKERRQ(ierr);
5634b9ad928SBarry Smith     }
5644b9ad928SBarry Smith   }
565b0de9f23SBarry Smith   /* Add the local solution to the global solution including the ghost nodes */
5661dd8081eSeaulisa   ierr = VecScatterBegin(osm->restriction, osm->ly, y, ADD_VALUES, reverse);CHKERRQ(ierr);
5671dd8081eSeaulisa   ierr = VecScatterEnd(osm->restriction, osm->ly, y, ADD_VALUES, reverse);CHKERRQ(ierr);
568d8b3b5e3Seaulisa 
5694b9ad928SBarry Smith   PetscFunctionReturn(0);
570d8b3b5e3Seaulisa 
5714b9ad928SBarry Smith }
5724b9ad928SBarry Smith 
573e91c6855SBarry Smith static PetscErrorCode PCReset_ASM(PC pc)
5744b9ad928SBarry Smith {
5754b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
5766849ba73SBarry Smith   PetscErrorCode ierr;
57713f74950SBarry Smith   PetscInt       i;
5784b9ad928SBarry Smith 
5794b9ad928SBarry Smith   PetscFunctionBegin;
58092bb6962SLisandro Dalcin   if (osm->ksp) {
58192bb6962SLisandro Dalcin     for (i=0; i<osm->n_local_true; i++) {
582e91c6855SBarry Smith       ierr = KSPReset(osm->ksp[i]);CHKERRQ(ierr);
58392bb6962SLisandro Dalcin     }
58492bb6962SLisandro Dalcin   }
585e09e08ccSBarry Smith   if (osm->pmat) {
58692bb6962SLisandro Dalcin     if (osm->n_local_true > 0) {
58730a70a9aSHong Zhang       ierr = MatDestroySubMatrices(osm->n_local_true,&osm->pmat);CHKERRQ(ierr);
58892bb6962SLisandro Dalcin     }
58992bb6962SLisandro Dalcin   }
5901dd8081eSeaulisa   if (osm->lrestriction) {
5911dd8081eSeaulisa     ierr = VecScatterDestroy(&osm->restriction);CHKERRQ(ierr);
5921dd8081eSeaulisa     for (i=0; i<osm->n_local_true; i++) {
5931dd8081eSeaulisa       ierr = VecScatterDestroy(&osm->lrestriction[i]);CHKERRQ(ierr);
5941dd8081eSeaulisa       if (osm->lprolongation) {ierr = VecScatterDestroy(&osm->lprolongation[i]);CHKERRQ(ierr);}
595fcfd50ebSBarry Smith       ierr = VecDestroy(&osm->x[i]);CHKERRQ(ierr);
596fcfd50ebSBarry Smith       ierr = VecDestroy(&osm->y[i]);CHKERRQ(ierr);
5974b9ad928SBarry Smith     }
5981dd8081eSeaulisa     ierr = PetscFree(osm->lrestriction);CHKERRQ(ierr);
5991dd8081eSeaulisa     if (osm->lprolongation) {ierr = PetscFree(osm->lprolongation);CHKERRQ(ierr);}
60005b42c5fSBarry Smith     ierr = PetscFree(osm->x);CHKERRQ(ierr);
60178904715SLisandro Dalcin     ierr = PetscFree(osm->y);CHKERRQ(ierr);
6021dd8081eSeaulisa 
60392bb6962SLisandro Dalcin   }
6042b691e39SMatthew Knepley   ierr = PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);CHKERRQ(ierr);
605fb745f2cSMatthew G. Knepley   ierr = ISDestroy(&osm->lis);CHKERRQ(ierr);
606fb745f2cSMatthew G. Knepley   ierr = VecDestroy(&osm->lx);CHKERRQ(ierr);
607fb745f2cSMatthew G. Knepley   ierr = VecDestroy(&osm->ly);CHKERRQ(ierr);
608347574c9Seaulisa   if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) {
609347574c9Seaulisa     ierr = MatDestroyMatrices(osm->n_local_true, &osm->lmats);CHKERRQ(ierr);
610fb745f2cSMatthew G. Knepley   }
6112fa5cd67SKarl Rupp 
61280ec0b7dSPatrick Sanan   ierr = PetscFree(osm->sub_mat_type);CHKERRQ(ierr);
61380ec0b7dSPatrick Sanan 
614e91c6855SBarry Smith   osm->is       = 0;
615e91c6855SBarry Smith   osm->is_local = 0;
616e91c6855SBarry Smith   PetscFunctionReturn(0);
617e91c6855SBarry Smith }
618e91c6855SBarry Smith 
619e91c6855SBarry Smith static PetscErrorCode PCDestroy_ASM(PC pc)
620e91c6855SBarry Smith {
621e91c6855SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
622e91c6855SBarry Smith   PetscErrorCode ierr;
623e91c6855SBarry Smith   PetscInt       i;
624e91c6855SBarry Smith 
625e91c6855SBarry Smith   PetscFunctionBegin;
626e91c6855SBarry Smith   ierr = PCReset_ASM(pc);CHKERRQ(ierr);
627e91c6855SBarry Smith   if (osm->ksp) {
628e91c6855SBarry Smith     for (i=0; i<osm->n_local_true; i++) {
629fcfd50ebSBarry Smith       ierr = KSPDestroy(&osm->ksp[i]);CHKERRQ(ierr);
630e91c6855SBarry Smith     }
631e91c6855SBarry Smith     ierr = PetscFree(osm->ksp);CHKERRQ(ierr);
632e91c6855SBarry Smith   }
633e91c6855SBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
6344b9ad928SBarry Smith   PetscFunctionReturn(0);
6354b9ad928SBarry Smith }
6364b9ad928SBarry Smith 
6374416b707SBarry Smith static PetscErrorCode PCSetFromOptions_ASM(PetscOptionItems *PetscOptionsObject,PC pc)
6384b9ad928SBarry Smith {
6394b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
6406849ba73SBarry Smith   PetscErrorCode ierr;
6419dcbbd2bSBarry Smith   PetscInt       blocks,ovl;
64257501b6eSBarry Smith   PetscBool      flg;
64392bb6962SLisandro Dalcin   PCASMType      asmtype;
64412cd4985SMatthew G. Knepley   PCCompositeType loctype;
64580ec0b7dSPatrick Sanan   char           sub_mat_type[256];
6464b9ad928SBarry Smith 
6474b9ad928SBarry Smith   PetscFunctionBegin;
648e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Additive Schwarz options");CHKERRQ(ierr);
649d709ea83SDmitry Karpeev   ierr = PetscOptionsBool("-pc_asm_dm_subdomains","Use DMCreateDomainDecomposition() to define subdomains","PCASMSetDMSubdomains",osm->dm_subdomains,&osm->dm_subdomains,&flg);CHKERRQ(ierr);
65090d69ab7SBarry Smith   ierr = PetscOptionsInt("-pc_asm_blocks","Number of subdomains","PCASMSetTotalSubdomains",osm->n,&blocks,&flg);CHKERRQ(ierr);
65165db9045SDmitry Karpeev   if (flg) {
652f77a5242SKarl Rupp     ierr = PCASMSetTotalSubdomains(pc,blocks,NULL,NULL);CHKERRQ(ierr);
653d709ea83SDmitry Karpeev     osm->dm_subdomains = PETSC_FALSE;
65465db9045SDmitry Karpeev   }
655342c94f9SMatthew G. Knepley   ierr = PetscOptionsInt("-pc_asm_local_blocks","Number of local subdomains","PCASMSetLocalSubdomains",osm->n_local_true,&blocks,&flg);CHKERRQ(ierr);
656342c94f9SMatthew G. Knepley   if (flg) {
657342c94f9SMatthew G. Knepley     ierr = PCASMSetLocalSubdomains(pc,blocks,NULL,NULL);CHKERRQ(ierr);
658342c94f9SMatthew G. Knepley     osm->dm_subdomains = PETSC_FALSE;
659342c94f9SMatthew G. Knepley   }
66090d69ab7SBarry Smith   ierr = PetscOptionsInt("-pc_asm_overlap","Number of grid points overlap","PCASMSetOverlap",osm->overlap,&ovl,&flg);CHKERRQ(ierr);
66165db9045SDmitry Karpeev   if (flg) {
66265db9045SDmitry Karpeev     ierr = PCASMSetOverlap(pc,ovl);CHKERRQ(ierr);
663d709ea83SDmitry Karpeev     osm->dm_subdomains = PETSC_FALSE;
66465db9045SDmitry Karpeev   }
66590d69ab7SBarry Smith   flg  = PETSC_FALSE;
66690d69ab7SBarry Smith   ierr = PetscOptionsEnum("-pc_asm_type","Type of restriction/extension","PCASMSetType",PCASMTypes,(PetscEnum)osm->type,(PetscEnum*)&asmtype,&flg);CHKERRQ(ierr);
66792bb6962SLisandro Dalcin   if (flg) {ierr = PCASMSetType(pc,asmtype);CHKERRQ(ierr); }
66812cd4985SMatthew G. Knepley   flg  = PETSC_FALSE;
66912cd4985SMatthew G. Knepley   ierr = PetscOptionsEnum("-pc_asm_local_type","Type of local solver composition","PCASMSetLocalType",PCCompositeTypes,(PetscEnum)osm->loctype,(PetscEnum*)&loctype,&flg);CHKERRQ(ierr);
67012cd4985SMatthew G. Knepley   if (flg) {ierr = PCASMSetLocalType(pc,loctype);CHKERRQ(ierr); }
67180ec0b7dSPatrick Sanan   ierr = PetscOptionsFList("-pc_asm_sub_mat_type","Subsolve Matrix Type","PCASMSetSubMatType",MatList,NULL,sub_mat_type,256,&flg);CHKERRQ(ierr);
67280ec0b7dSPatrick Sanan   if (flg) {
673459726d8SSatish Balay     ierr = PCASMSetSubMatType(pc,sub_mat_type);CHKERRQ(ierr);
67480ec0b7dSPatrick Sanan   }
6754b9ad928SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
6764b9ad928SBarry Smith   PetscFunctionReturn(0);
6774b9ad928SBarry Smith }
6784b9ad928SBarry Smith 
6794b9ad928SBarry Smith /*------------------------------------------------------------------------------------*/
6804b9ad928SBarry Smith 
6811e6b0712SBarry Smith static PetscErrorCode  PCASMSetLocalSubdomains_ASM(PC pc,PetscInt n,IS is[],IS is_local[])
6824b9ad928SBarry Smith {
6834b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
68492bb6962SLisandro Dalcin   PetscErrorCode ierr;
68592bb6962SLisandro Dalcin   PetscInt       i;
6864b9ad928SBarry Smith 
6874b9ad928SBarry Smith   PetscFunctionBegin;
688e32f2f54SBarry Smith   if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Each process must have 1 or more blocks, n = %D",n);
689ce94432eSBarry Smith   if (pc->setupcalled && (n != osm->n_local_true || is)) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCASMSetLocalSubdomains() should be called before calling PCSetUp().");
690e7e72b3dSBarry Smith 
6914b9ad928SBarry Smith   if (!pc->setupcalled) {
69292bb6962SLisandro Dalcin     if (is) {
69392bb6962SLisandro Dalcin       for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)is[i]);CHKERRQ(ierr);}
69492bb6962SLisandro Dalcin     }
695832fc9a5SMatthew Knepley     if (is_local) {
696832fc9a5SMatthew Knepley       for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)is_local[i]);CHKERRQ(ierr);}
697832fc9a5SMatthew Knepley     }
6982b691e39SMatthew Knepley     ierr = PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);CHKERRQ(ierr);
6992fa5cd67SKarl Rupp 
7004b9ad928SBarry Smith     osm->n_local_true = n;
70192bb6962SLisandro Dalcin     osm->is           = 0;
7022b691e39SMatthew Knepley     osm->is_local     = 0;
70392bb6962SLisandro Dalcin     if (is) {
704785e854fSJed Brown       ierr = PetscMalloc1(n,&osm->is);CHKERRQ(ierr);
7052fa5cd67SKarl Rupp       for (i=0; i<n; i++) osm->is[i] = is[i];
7063d03bb48SJed Brown       /* Flag indicating that the user has set overlapping subdomains so PCASM should not increase their size. */
7073d03bb48SJed Brown       osm->overlap = -1;
70892bb6962SLisandro Dalcin     }
7092b691e39SMatthew Knepley     if (is_local) {
710785e854fSJed Brown       ierr = PetscMalloc1(n,&osm->is_local);CHKERRQ(ierr);
7112fa5cd67SKarl Rupp       for (i=0; i<n; i++) osm->is_local[i] = is_local[i];
712a35b7b57SDmitry Karpeev       if (!is) {
713785e854fSJed Brown         ierr = PetscMalloc1(osm->n_local_true,&osm->is);CHKERRQ(ierr);
714a35b7b57SDmitry Karpeev         for (i=0; i<osm->n_local_true; i++) {
715a35b7b57SDmitry Karpeev           if (osm->overlap > 0) { /* With positive overlap, osm->is[i] will be modified */
716a35b7b57SDmitry Karpeev             ierr = ISDuplicate(osm->is_local[i],&osm->is[i]);CHKERRQ(ierr);
717a35b7b57SDmitry Karpeev             ierr = ISCopy(osm->is_local[i],osm->is[i]);CHKERRQ(ierr);
718a35b7b57SDmitry Karpeev           } else {
719a35b7b57SDmitry Karpeev             ierr = PetscObjectReference((PetscObject)osm->is_local[i]);CHKERRQ(ierr);
720a35b7b57SDmitry Karpeev             osm->is[i] = osm->is_local[i];
721a35b7b57SDmitry Karpeev           }
722a35b7b57SDmitry Karpeev         }
723a35b7b57SDmitry Karpeev       }
7242b691e39SMatthew Knepley     }
7254b9ad928SBarry Smith   }
7264b9ad928SBarry Smith   PetscFunctionReturn(0);
7274b9ad928SBarry Smith }
7284b9ad928SBarry Smith 
7291e6b0712SBarry Smith static PetscErrorCode  PCASMSetTotalSubdomains_ASM(PC pc,PetscInt N,IS *is,IS *is_local)
7304b9ad928SBarry Smith {
7314b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
7326849ba73SBarry Smith   PetscErrorCode ierr;
73313f74950SBarry Smith   PetscMPIInt    rank,size;
73478904715SLisandro Dalcin   PetscInt       n;
7354b9ad928SBarry Smith 
7364b9ad928SBarry Smith   PetscFunctionBegin;
737ce94432eSBarry Smith   if (N < 1) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Number of total blocks must be > 0, N = %D",N);
738ce94432eSBarry Smith   if (is || is_local) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Use PCASMSetLocalSubdomains() to set specific index sets\n\they cannot be set globally yet.");
7394b9ad928SBarry Smith 
7404b9ad928SBarry Smith   /*
741880770e9SJed Brown      Split the subdomains equally among all processors
7424b9ad928SBarry Smith   */
743ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr);
744ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr);
7454b9ad928SBarry Smith   n    = N/size + ((N % size) > rank);
746e32f2f54SBarry Smith   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);
747e32f2f54SBarry Smith   if (pc->setupcalled && n != osm->n_local_true) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"PCASMSetTotalSubdomains() should be called before PCSetUp().");
7484b9ad928SBarry Smith   if (!pc->setupcalled) {
7492b691e39SMatthew Knepley     ierr = PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);CHKERRQ(ierr);
7502fa5cd67SKarl Rupp 
7514b9ad928SBarry Smith     osm->n_local_true = n;
7524b9ad928SBarry Smith     osm->is           = 0;
7532b691e39SMatthew Knepley     osm->is_local     = 0;
7544b9ad928SBarry Smith   }
7554b9ad928SBarry Smith   PetscFunctionReturn(0);
7564b9ad928SBarry Smith }
7574b9ad928SBarry Smith 
7581e6b0712SBarry Smith static PetscErrorCode  PCASMSetOverlap_ASM(PC pc,PetscInt ovl)
7594b9ad928SBarry Smith {
76092bb6962SLisandro Dalcin   PC_ASM *osm = (PC_ASM*)pc->data;
7614b9ad928SBarry Smith 
7624b9ad928SBarry Smith   PetscFunctionBegin;
763ce94432eSBarry Smith   if (ovl < 0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested");
764ce94432eSBarry Smith   if (pc->setupcalled && ovl != osm->overlap) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCASMSetOverlap() should be called before PCSetUp().");
7652fa5cd67SKarl Rupp   if (!pc->setupcalled) osm->overlap = ovl;
7664b9ad928SBarry Smith   PetscFunctionReturn(0);
7674b9ad928SBarry Smith }
7684b9ad928SBarry Smith 
7691e6b0712SBarry Smith static PetscErrorCode  PCASMSetType_ASM(PC pc,PCASMType type)
7704b9ad928SBarry Smith {
77192bb6962SLisandro Dalcin   PC_ASM *osm = (PC_ASM*)pc->data;
7724b9ad928SBarry Smith 
7734b9ad928SBarry Smith   PetscFunctionBegin;
7744b9ad928SBarry Smith   osm->type     = type;
775bf108f30SBarry Smith   osm->type_set = PETSC_TRUE;
7764b9ad928SBarry Smith   PetscFunctionReturn(0);
7774b9ad928SBarry Smith }
7784b9ad928SBarry Smith 
779c60c7ad4SBarry Smith static PetscErrorCode  PCASMGetType_ASM(PC pc,PCASMType *type)
780c60c7ad4SBarry Smith {
781c60c7ad4SBarry Smith   PC_ASM *osm = (PC_ASM*)pc->data;
782c60c7ad4SBarry Smith 
783c60c7ad4SBarry Smith   PetscFunctionBegin;
784c60c7ad4SBarry Smith   *type = osm->type;
785c60c7ad4SBarry Smith   PetscFunctionReturn(0);
786c60c7ad4SBarry Smith }
787c60c7ad4SBarry Smith 
78812cd4985SMatthew G. Knepley static PetscErrorCode  PCASMSetLocalType_ASM(PC pc, PCCompositeType type)
78912cd4985SMatthew G. Knepley {
79012cd4985SMatthew G. Knepley   PC_ASM *osm = (PC_ASM *) pc->data;
79112cd4985SMatthew G. Knepley 
79212cd4985SMatthew G. Knepley   PetscFunctionBegin;
793d0ecd4dfSBarry Smith   if (type != PC_COMPOSITE_ADDITIVE && type != PC_COMPOSITE_MULTIPLICATIVE) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Only supports additive or multiplicative as the local type");
79412cd4985SMatthew G. Knepley   osm->loctype = type;
79512cd4985SMatthew G. Knepley   PetscFunctionReturn(0);
79612cd4985SMatthew G. Knepley }
79712cd4985SMatthew G. Knepley 
79812cd4985SMatthew G. Knepley static PetscErrorCode  PCASMGetLocalType_ASM(PC pc, PCCompositeType *type)
79912cd4985SMatthew G. Knepley {
80012cd4985SMatthew G. Knepley   PC_ASM *osm = (PC_ASM *) pc->data;
80112cd4985SMatthew G. Knepley 
80212cd4985SMatthew G. Knepley   PetscFunctionBegin;
80312cd4985SMatthew G. Knepley   *type = osm->loctype;
80412cd4985SMatthew G. Knepley   PetscFunctionReturn(0);
80512cd4985SMatthew G. Knepley }
80612cd4985SMatthew G. Knepley 
8071e6b0712SBarry Smith static PetscErrorCode  PCASMSetSortIndices_ASM(PC pc,PetscBool  doSort)
8086ed231c7SMatthew Knepley {
8096ed231c7SMatthew Knepley   PC_ASM *osm = (PC_ASM*)pc->data;
8106ed231c7SMatthew Knepley 
8116ed231c7SMatthew Knepley   PetscFunctionBegin;
8126ed231c7SMatthew Knepley   osm->sort_indices = doSort;
8136ed231c7SMatthew Knepley   PetscFunctionReturn(0);
8146ed231c7SMatthew Knepley }
8156ed231c7SMatthew Knepley 
8161e6b0712SBarry Smith static PetscErrorCode  PCASMGetSubKSP_ASM(PC pc,PetscInt *n_local,PetscInt *first_local,KSP **ksp)
8174b9ad928SBarry Smith {
81892bb6962SLisandro Dalcin   PC_ASM         *osm = (PC_ASM*)pc->data;
819dfbe8321SBarry Smith   PetscErrorCode ierr;
8204b9ad928SBarry Smith 
8214b9ad928SBarry Smith   PetscFunctionBegin;
82234a84908Sprj-   if (osm->n_local_true < 1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Need to call PCSetUp() on PC (or KSPSetUp() on the outer KSP object) before calling here");
8234b9ad928SBarry Smith 
8242fa5cd67SKarl Rupp   if (n_local) *n_local = osm->n_local_true;
82592bb6962SLisandro Dalcin   if (first_local) {
826ce94432eSBarry Smith     ierr          = MPI_Scan(&osm->n_local_true,first_local,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
82792bb6962SLisandro Dalcin     *first_local -= osm->n_local_true;
82892bb6962SLisandro Dalcin   }
82992bb6962SLisandro Dalcin   if (ksp) {
83092bb6962SLisandro Dalcin     /* Assume that local solves are now different; not necessarily
83192bb6962SLisandro Dalcin        true though!  This flag is used only for PCView_ASM() */
83292bb6962SLisandro Dalcin     *ksp                   = osm->ksp;
83392bb6962SLisandro Dalcin     osm->same_local_solves = PETSC_FALSE;
83492bb6962SLisandro Dalcin   }
8354b9ad928SBarry Smith   PetscFunctionReturn(0);
8364b9ad928SBarry Smith }
8374b9ad928SBarry Smith 
83880ec0b7dSPatrick Sanan static PetscErrorCode  PCASMGetSubMatType_ASM(PC pc,MatType *sub_mat_type)
83980ec0b7dSPatrick Sanan {
84080ec0b7dSPatrick Sanan   PC_ASM         *osm = (PC_ASM*)pc->data;
84180ec0b7dSPatrick Sanan 
84280ec0b7dSPatrick Sanan   PetscFunctionBegin;
84380ec0b7dSPatrick Sanan   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
84480ec0b7dSPatrick Sanan   PetscValidPointer(sub_mat_type,2);
84580ec0b7dSPatrick Sanan   *sub_mat_type = osm->sub_mat_type;
84680ec0b7dSPatrick Sanan   PetscFunctionReturn(0);
84780ec0b7dSPatrick Sanan }
84880ec0b7dSPatrick Sanan 
84980ec0b7dSPatrick Sanan static PetscErrorCode PCASMSetSubMatType_ASM(PC pc,MatType sub_mat_type)
85080ec0b7dSPatrick Sanan {
85180ec0b7dSPatrick Sanan   PetscErrorCode    ierr;
85280ec0b7dSPatrick Sanan   PC_ASM            *osm = (PC_ASM*)pc->data;
85380ec0b7dSPatrick Sanan 
85480ec0b7dSPatrick Sanan   PetscFunctionBegin;
85580ec0b7dSPatrick Sanan   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
85680ec0b7dSPatrick Sanan   ierr = PetscFree(osm->sub_mat_type);CHKERRQ(ierr);
85780ec0b7dSPatrick Sanan   ierr = PetscStrallocpy(sub_mat_type,(char**)&osm->sub_mat_type);CHKERRQ(ierr);
85880ec0b7dSPatrick Sanan   PetscFunctionReturn(0);
85980ec0b7dSPatrick Sanan }
86080ec0b7dSPatrick Sanan 
8614b9ad928SBarry Smith /*@C
8621093a601SBarry Smith     PCASMSetLocalSubdomains - Sets the local subdomains (for this processor only) for the additive Schwarz preconditioner.
8634b9ad928SBarry Smith 
864d083f849SBarry Smith     Collective on pc
8654b9ad928SBarry Smith 
8664b9ad928SBarry Smith     Input Parameters:
8674b9ad928SBarry Smith +   pc - the preconditioner context
8684b9ad928SBarry Smith .   n - the number of subdomains for this processor (default value = 1)
8698c03b21aSDmitry Karpeev .   is - the index set that defines the subdomains for this processor
8700298fd71SBarry Smith          (or NULL for PETSc to determine subdomains)
871f1ee410cSBarry Smith -   is_local - the index sets that define the local part of the subdomains for this processor, not used unless PCASMType is PC_ASM_RESTRICT
872f1ee410cSBarry Smith          (or NULL to not provide these)
8734b9ad928SBarry Smith 
874342c94f9SMatthew G. Knepley     Options Database Key:
875342c94f9SMatthew G. Knepley     To set the total number of subdomain blocks rather than specify the
876342c94f9SMatthew G. Knepley     index sets, use the option
877342c94f9SMatthew G. Knepley .    -pc_asm_local_blocks <blks> - Sets local blocks
878342c94f9SMatthew G. Knepley 
8794b9ad928SBarry Smith     Notes:
8801093a601SBarry Smith     The IS numbering is in the parallel, global numbering of the vector for both is and is_local
8814b9ad928SBarry Smith 
8824b9ad928SBarry Smith     By default the ASM preconditioner uses 1 block per processor.
8834b9ad928SBarry Smith 
8844b9ad928SBarry Smith     Use PCASMSetTotalSubdomains() to set the subdomains for all processors.
8854b9ad928SBarry Smith 
886f1ee410cSBarry Smith     If is_local is provided and PCASMType is PC_ASM_RESTRICT then the solution only over the is_local region is interpolated
887f1ee410cSBarry Smith     back to form the global solution (this is the standard restricted additive Schwarz method)
888f1ee410cSBarry Smith 
889f1ee410cSBarry Smith     If the is_local is provided and PCASMType is PC_ASM_INTERPOLATE or PC_ASM_NONE then an error is generated since there is
890f1ee410cSBarry Smith     no code to handle that case.
891f1ee410cSBarry Smith 
8924b9ad928SBarry Smith     Level: advanced
8934b9ad928SBarry Smith 
8944b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(),
895f1ee410cSBarry Smith           PCASMCreateSubdomains2D(), PCASMGetLocalSubdomains(), PCASMType, PCASMSetType()
8964b9ad928SBarry Smith @*/
8977087cfbeSBarry Smith PetscErrorCode  PCASMSetLocalSubdomains(PC pc,PetscInt n,IS is[],IS is_local[])
8984b9ad928SBarry Smith {
8997bb14e67SBarry Smith   PetscErrorCode ierr;
9004b9ad928SBarry Smith 
9014b9ad928SBarry Smith   PetscFunctionBegin;
9020700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
9037bb14e67SBarry Smith   ierr = PetscTryMethod(pc,"PCASMSetLocalSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,n,is,is_local));CHKERRQ(ierr);
9044b9ad928SBarry Smith   PetscFunctionReturn(0);
9054b9ad928SBarry Smith }
9064b9ad928SBarry Smith 
9074b9ad928SBarry Smith /*@C
908feb221f8SDmitry Karpeev     PCASMSetTotalSubdomains - Sets the subdomains for all processors for the
9094b9ad928SBarry Smith     additive Schwarz preconditioner.  Either all or no processors in the
9104b9ad928SBarry Smith     PC communicator must call this routine, with the same index sets.
9114b9ad928SBarry Smith 
912d083f849SBarry Smith     Collective on pc
9134b9ad928SBarry Smith 
9144b9ad928SBarry Smith     Input Parameters:
9154b9ad928SBarry Smith +   pc - the preconditioner context
916feb221f8SDmitry Karpeev .   N  - the number of subdomains for all processors
917feb221f8SDmitry Karpeev .   is - the index sets that define the subdomains for all processors
918dfaa0c5fSBarry Smith          (or NULL to ask PETSc to determine the subdomains)
9192b691e39SMatthew Knepley -   is_local - the index sets that define the local part of the subdomains for this processor
920f1ee410cSBarry Smith          (or NULL to not provide this information)
9214b9ad928SBarry Smith 
9224b9ad928SBarry Smith     Options Database Key:
9234b9ad928SBarry Smith     To set the total number of subdomain blocks rather than specify the
9244b9ad928SBarry Smith     index sets, use the option
9254b9ad928SBarry Smith .    -pc_asm_blocks <blks> - Sets total blocks
9264b9ad928SBarry Smith 
9274b9ad928SBarry Smith     Notes:
928f1ee410cSBarry Smith     Currently you cannot use this to set the actual subdomains with the argument is or is_local.
9294b9ad928SBarry Smith 
9304b9ad928SBarry Smith     By default the ASM preconditioner uses 1 block per processor.
9314b9ad928SBarry Smith 
9324b9ad928SBarry Smith     These index sets cannot be destroyed until after completion of the
9334b9ad928SBarry Smith     linear solves for which the ASM preconditioner is being used.
9344b9ad928SBarry Smith 
9354b9ad928SBarry Smith     Use PCASMSetLocalSubdomains() to set local subdomains.
9364b9ad928SBarry Smith 
9371093a601SBarry Smith     The IS numbering is in the parallel, global numbering of the vector for both is and is_local
9381093a601SBarry Smith 
9394b9ad928SBarry Smith     Level: advanced
9404b9ad928SBarry Smith 
9414b9ad928SBarry Smith .seealso: PCASMSetLocalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(),
9424b9ad928SBarry Smith           PCASMCreateSubdomains2D()
9434b9ad928SBarry Smith @*/
9447087cfbeSBarry Smith PetscErrorCode  PCASMSetTotalSubdomains(PC pc,PetscInt N,IS is[],IS is_local[])
9454b9ad928SBarry Smith {
9467bb14e67SBarry Smith   PetscErrorCode ierr;
9474b9ad928SBarry Smith 
9484b9ad928SBarry Smith   PetscFunctionBegin;
9490700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
9507bb14e67SBarry Smith   ierr = PetscTryMethod(pc,"PCASMSetTotalSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,N,is,is_local));CHKERRQ(ierr);
9514b9ad928SBarry Smith   PetscFunctionReturn(0);
9524b9ad928SBarry Smith }
9534b9ad928SBarry Smith 
9544b9ad928SBarry Smith /*@
9554b9ad928SBarry Smith     PCASMSetOverlap - Sets the overlap between a pair of subdomains for the
9564b9ad928SBarry Smith     additive Schwarz preconditioner.  Either all or no processors in the
957f1ee410cSBarry Smith     PC communicator must call this routine.
9584b9ad928SBarry Smith 
959d083f849SBarry Smith     Logically Collective on pc
9604b9ad928SBarry Smith 
9614b9ad928SBarry Smith     Input Parameters:
9624b9ad928SBarry Smith +   pc  - the preconditioner context
9634b9ad928SBarry Smith -   ovl - the amount of overlap between subdomains (ovl >= 0, default value = 1)
9644b9ad928SBarry Smith 
9654b9ad928SBarry Smith     Options Database Key:
9664b9ad928SBarry Smith .   -pc_asm_overlap <ovl> - Sets overlap
9674b9ad928SBarry Smith 
9684b9ad928SBarry Smith     Notes:
9694b9ad928SBarry Smith     By default the ASM preconditioner uses 1 block per processor.  To use
9704b9ad928SBarry Smith     multiple blocks per perocessor, see PCASMSetTotalSubdomains() and
9714b9ad928SBarry Smith     PCASMSetLocalSubdomains() (and the option -pc_asm_blocks <blks>).
9724b9ad928SBarry Smith 
9734b9ad928SBarry Smith     The overlap defaults to 1, so if one desires that no additional
9744b9ad928SBarry Smith     overlap be computed beyond what may have been set with a call to
9754b9ad928SBarry Smith     PCASMSetTotalSubdomains() or PCASMSetLocalSubdomains(), then ovl
9764b9ad928SBarry Smith     must be set to be 0.  In particular, if one does not explicitly set
9774b9ad928SBarry Smith     the subdomains an application code, then all overlap would be computed
9784b9ad928SBarry Smith     internally by PETSc, and using an overlap of 0 would result in an ASM
9794b9ad928SBarry Smith     variant that is equivalent to the block Jacobi preconditioner.
9804b9ad928SBarry Smith 
981f1ee410cSBarry Smith     The default algorithm used by PETSc to increase overlap is fast, but not scalable,
982f1ee410cSBarry Smith     use the option -mat_increase_overlap_scalable when the problem and number of processes is large.
983f1ee410cSBarry Smith 
9844b9ad928SBarry Smith     Note that one can define initial index sets with any overlap via
985f1ee410cSBarry Smith     PCASMSetLocalSubdomains(); the routine
9864b9ad928SBarry Smith     PCASMSetOverlap() merely allows PETSc to extend that overlap further
9874b9ad928SBarry Smith     if desired.
9884b9ad928SBarry Smith 
9894b9ad928SBarry Smith     Level: intermediate
9904b9ad928SBarry Smith 
9914b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetLocalSubdomains(), PCASMGetSubKSP(),
992f1ee410cSBarry Smith           PCASMCreateSubdomains2D(), PCASMGetLocalSubdomains(), MatIncreaseOverlap()
9934b9ad928SBarry Smith @*/
9947087cfbeSBarry Smith PetscErrorCode  PCASMSetOverlap(PC pc,PetscInt ovl)
9954b9ad928SBarry Smith {
9967bb14e67SBarry Smith   PetscErrorCode ierr;
9974b9ad928SBarry Smith 
9984b9ad928SBarry Smith   PetscFunctionBegin;
9990700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1000c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,ovl,2);
10017bb14e67SBarry Smith   ierr = PetscTryMethod(pc,"PCASMSetOverlap_C",(PC,PetscInt),(pc,ovl));CHKERRQ(ierr);
10024b9ad928SBarry Smith   PetscFunctionReturn(0);
10034b9ad928SBarry Smith }
10044b9ad928SBarry Smith 
10054b9ad928SBarry Smith /*@
10064b9ad928SBarry Smith     PCASMSetType - Sets the type of restriction and interpolation used
10074b9ad928SBarry Smith     for local problems in the additive Schwarz method.
10084b9ad928SBarry Smith 
1009d083f849SBarry Smith     Logically Collective on pc
10104b9ad928SBarry Smith 
10114b9ad928SBarry Smith     Input Parameters:
10124b9ad928SBarry Smith +   pc  - the preconditioner context
10134b9ad928SBarry Smith -   type - variant of ASM, one of
10144b9ad928SBarry Smith .vb
10154b9ad928SBarry Smith       PC_ASM_BASIC       - full interpolation and restriction
10164b9ad928SBarry Smith       PC_ASM_RESTRICT    - full restriction, local processor interpolation
10174b9ad928SBarry Smith       PC_ASM_INTERPOLATE - full interpolation, local processor restriction
10184b9ad928SBarry Smith       PC_ASM_NONE        - local processor restriction and interpolation
10194b9ad928SBarry Smith .ve
10204b9ad928SBarry Smith 
10214b9ad928SBarry Smith     Options Database Key:
10224b9ad928SBarry Smith .   -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type
10234b9ad928SBarry Smith 
102495452b02SPatrick Sanan     Notes:
102595452b02SPatrick Sanan     if the is_local arguments are passed to PCASMSetLocalSubdomains() then they are used when PC_ASM_RESTRICT has been selected
1026f1ee410cSBarry Smith     to limit the local processor interpolation
1027f1ee410cSBarry Smith 
10284b9ad928SBarry Smith     Level: intermediate
10294b9ad928SBarry Smith 
10304b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(),
1031f1ee410cSBarry Smith           PCASMCreateSubdomains2D(), PCASMType, PCASMSetLocalType(), PCASMGetLocalType()
10324b9ad928SBarry Smith @*/
10337087cfbeSBarry Smith PetscErrorCode  PCASMSetType(PC pc,PCASMType type)
10344b9ad928SBarry Smith {
10357bb14e67SBarry Smith   PetscErrorCode ierr;
10364b9ad928SBarry Smith 
10374b9ad928SBarry Smith   PetscFunctionBegin;
10380700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1039c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,type,2);
10407bb14e67SBarry Smith   ierr = PetscTryMethod(pc,"PCASMSetType_C",(PC,PCASMType),(pc,type));CHKERRQ(ierr);
10414b9ad928SBarry Smith   PetscFunctionReturn(0);
10424b9ad928SBarry Smith }
10434b9ad928SBarry Smith 
1044c60c7ad4SBarry Smith /*@
1045c60c7ad4SBarry Smith     PCASMGetType - Gets the type of restriction and interpolation used
1046c60c7ad4SBarry Smith     for local problems in the additive Schwarz method.
1047c60c7ad4SBarry Smith 
1048d083f849SBarry Smith     Logically Collective on pc
1049c60c7ad4SBarry Smith 
1050c60c7ad4SBarry Smith     Input Parameter:
1051c60c7ad4SBarry Smith .   pc  - the preconditioner context
1052c60c7ad4SBarry Smith 
1053c60c7ad4SBarry Smith     Output Parameter:
1054c60c7ad4SBarry Smith .   type - variant of ASM, one of
1055c60c7ad4SBarry Smith 
1056c60c7ad4SBarry Smith .vb
1057c60c7ad4SBarry Smith       PC_ASM_BASIC       - full interpolation and restriction
1058c60c7ad4SBarry Smith       PC_ASM_RESTRICT    - full restriction, local processor interpolation
1059c60c7ad4SBarry Smith       PC_ASM_INTERPOLATE - full interpolation, local processor restriction
1060c60c7ad4SBarry Smith       PC_ASM_NONE        - local processor restriction and interpolation
1061c60c7ad4SBarry Smith .ve
1062c60c7ad4SBarry Smith 
1063c60c7ad4SBarry Smith     Options Database Key:
1064c60c7ad4SBarry Smith .   -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type
1065c60c7ad4SBarry Smith 
1066c60c7ad4SBarry Smith     Level: intermediate
1067c60c7ad4SBarry Smith 
1068c60c7ad4SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(),
1069f1ee410cSBarry Smith           PCASMCreateSubdomains2D(), PCASMType, PCASMSetType(), PCASMSetLocalType(), PCASMGetLocalType()
1070c60c7ad4SBarry Smith @*/
1071c60c7ad4SBarry Smith PetscErrorCode  PCASMGetType(PC pc,PCASMType *type)
1072c60c7ad4SBarry Smith {
1073c60c7ad4SBarry Smith   PetscErrorCode ierr;
1074c60c7ad4SBarry Smith 
1075c60c7ad4SBarry Smith   PetscFunctionBegin;
1076c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1077c60c7ad4SBarry Smith   ierr = PetscUseMethod(pc,"PCASMGetType_C",(PC,PCASMType*),(pc,type));CHKERRQ(ierr);
1078c60c7ad4SBarry Smith   PetscFunctionReturn(0);
1079c60c7ad4SBarry Smith }
1080c60c7ad4SBarry Smith 
108112cd4985SMatthew G. Knepley /*@
108212cd4985SMatthew G. Knepley   PCASMSetLocalType - Sets the type of composition used for local problems in the additive Schwarz method.
108312cd4985SMatthew G. Knepley 
1084d083f849SBarry Smith   Logically Collective on pc
108512cd4985SMatthew G. Knepley 
108612cd4985SMatthew G. Knepley   Input Parameters:
108712cd4985SMatthew G. Knepley + pc  - the preconditioner context
108812cd4985SMatthew G. Knepley - type - type of composition, one of
108912cd4985SMatthew G. Knepley .vb
109012cd4985SMatthew G. Knepley   PC_COMPOSITE_ADDITIVE       - local additive combination
109112cd4985SMatthew G. Knepley   PC_COMPOSITE_MULTIPLICATIVE - local multiplicative combination
109212cd4985SMatthew G. Knepley .ve
109312cd4985SMatthew G. Knepley 
109412cd4985SMatthew G. Knepley   Options Database Key:
109512cd4985SMatthew G. Knepley . -pc_asm_local_type [additive,multiplicative] - Sets local solver composition type
109612cd4985SMatthew G. Knepley 
109712cd4985SMatthew G. Knepley   Level: intermediate
109812cd4985SMatthew G. Knepley 
1099f1ee410cSBarry Smith .seealso: PCASMSetType(), PCASMGetType(), PCASMGetLocalType(), PCASM, PCASMType, PCASMSetType(), PCASMGetType(), PCCompositeType
110012cd4985SMatthew G. Knepley @*/
110112cd4985SMatthew G. Knepley PetscErrorCode PCASMSetLocalType(PC pc, PCCompositeType type)
110212cd4985SMatthew G. Knepley {
110312cd4985SMatthew G. Knepley   PetscErrorCode ierr;
110412cd4985SMatthew G. Knepley 
110512cd4985SMatthew G. Knepley   PetscFunctionBegin;
110612cd4985SMatthew G. Knepley   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
110712cd4985SMatthew G. Knepley   PetscValidLogicalCollectiveEnum(pc, type, 2);
110812cd4985SMatthew G. Knepley   ierr = PetscTryMethod(pc, "PCASMSetLocalType_C", (PC, PCCompositeType), (pc, type));CHKERRQ(ierr);
110912cd4985SMatthew G. Knepley   PetscFunctionReturn(0);
111012cd4985SMatthew G. Knepley }
111112cd4985SMatthew G. Knepley 
111212cd4985SMatthew G. Knepley /*@
111312cd4985SMatthew G. Knepley   PCASMGetLocalType - Gets the type of composition used for local problems in the additive Schwarz method.
111412cd4985SMatthew G. Knepley 
1115d083f849SBarry Smith   Logically Collective on pc
111612cd4985SMatthew G. Knepley 
111712cd4985SMatthew G. Knepley   Input Parameter:
111812cd4985SMatthew G. Knepley . pc  - the preconditioner context
111912cd4985SMatthew G. Knepley 
112012cd4985SMatthew G. Knepley   Output Parameter:
112112cd4985SMatthew G. Knepley . type - type of composition, one of
112212cd4985SMatthew G. Knepley .vb
112312cd4985SMatthew G. Knepley   PC_COMPOSITE_ADDITIVE       - local additive combination
112412cd4985SMatthew G. Knepley   PC_COMPOSITE_MULTIPLICATIVE - local multiplicative combination
112512cd4985SMatthew G. Knepley .ve
112612cd4985SMatthew G. Knepley 
112712cd4985SMatthew G. Knepley   Options Database Key:
112812cd4985SMatthew G. Knepley . -pc_asm_local_type [additive,multiplicative] - Sets local solver composition type
112912cd4985SMatthew G. Knepley 
113012cd4985SMatthew G. Knepley   Level: intermediate
113112cd4985SMatthew G. Knepley 
1132f1ee410cSBarry Smith .seealso: PCASMSetType(), PCASMGetType(), PCASMSetLocalType(), PCASMCreate(), PCASMType, PCASMSetType(), PCASMGetType(), PCCompositeType
113312cd4985SMatthew G. Knepley @*/
113412cd4985SMatthew G. Knepley PetscErrorCode PCASMGetLocalType(PC pc, PCCompositeType *type)
113512cd4985SMatthew G. Knepley {
113612cd4985SMatthew G. Knepley   PetscErrorCode ierr;
113712cd4985SMatthew G. Knepley 
113812cd4985SMatthew G. Knepley   PetscFunctionBegin;
113912cd4985SMatthew G. Knepley   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
114012cd4985SMatthew G. Knepley   PetscValidPointer(type, 2);
114112cd4985SMatthew G. Knepley   ierr = PetscUseMethod(pc, "PCASMGetLocalType_C", (PC, PCCompositeType *), (pc, type));CHKERRQ(ierr);
114212cd4985SMatthew G. Knepley   PetscFunctionReturn(0);
114312cd4985SMatthew G. Knepley }
114412cd4985SMatthew G. Knepley 
11456ed231c7SMatthew Knepley /*@
11466ed231c7SMatthew Knepley     PCASMSetSortIndices - Determines whether subdomain indices are sorted.
11476ed231c7SMatthew Knepley 
1148d083f849SBarry Smith     Logically Collective on pc
11496ed231c7SMatthew Knepley 
11506ed231c7SMatthew Knepley     Input Parameters:
11516ed231c7SMatthew Knepley +   pc  - the preconditioner context
11526ed231c7SMatthew Knepley -   doSort - sort the subdomain indices
11536ed231c7SMatthew Knepley 
11546ed231c7SMatthew Knepley     Level: intermediate
11556ed231c7SMatthew Knepley 
11566ed231c7SMatthew Knepley .seealso: PCASMSetLocalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(),
11576ed231c7SMatthew Knepley           PCASMCreateSubdomains2D()
11586ed231c7SMatthew Knepley @*/
11597087cfbeSBarry Smith PetscErrorCode  PCASMSetSortIndices(PC pc,PetscBool doSort)
11606ed231c7SMatthew Knepley {
11617bb14e67SBarry Smith   PetscErrorCode ierr;
11626ed231c7SMatthew Knepley 
11636ed231c7SMatthew Knepley   PetscFunctionBegin;
11640700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1165acfcf0e5SJed Brown   PetscValidLogicalCollectiveBool(pc,doSort,2);
11667bb14e67SBarry Smith   ierr = PetscTryMethod(pc,"PCASMSetSortIndices_C",(PC,PetscBool),(pc,doSort));CHKERRQ(ierr);
11676ed231c7SMatthew Knepley   PetscFunctionReturn(0);
11686ed231c7SMatthew Knepley }
11696ed231c7SMatthew Knepley 
11704b9ad928SBarry Smith /*@C
11714b9ad928SBarry Smith    PCASMGetSubKSP - Gets the local KSP contexts for all blocks on
11724b9ad928SBarry Smith    this processor.
11734b9ad928SBarry Smith 
1174d083f849SBarry Smith    Collective on pc iff first_local is requested
11754b9ad928SBarry Smith 
11764b9ad928SBarry Smith    Input Parameter:
11774b9ad928SBarry Smith .  pc - the preconditioner context
11784b9ad928SBarry Smith 
11794b9ad928SBarry Smith    Output Parameters:
11800298fd71SBarry Smith +  n_local - the number of blocks on this processor or NULL
11810298fd71SBarry Smith .  first_local - the global number of the first block on this processor or NULL,
11820298fd71SBarry Smith                  all processors must request or all must pass NULL
11834b9ad928SBarry Smith -  ksp - the array of KSP contexts
11844b9ad928SBarry Smith 
11854b9ad928SBarry Smith    Note:
1186d29017ddSJed Brown    After PCASMGetSubKSP() the array of KSPes is not to be freed.
11874b9ad928SBarry Smith 
11884b9ad928SBarry Smith    You must call KSPSetUp() before calling PCASMGetSubKSP().
11894b9ad928SBarry Smith 
1190d29017ddSJed Brown    Fortran note:
11912bf68e3eSBarry Smith    The output argument 'ksp' must be an array of sufficient length or PETSC_NULL_KSP. The latter can be used to learn the necessary length.
1192d29017ddSJed Brown 
11934b9ad928SBarry Smith    Level: advanced
11944b9ad928SBarry Smith 
11954b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap(),
11964b9ad928SBarry Smith           PCASMCreateSubdomains2D(),
11974b9ad928SBarry Smith @*/
11987087cfbeSBarry Smith PetscErrorCode  PCASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[])
11994b9ad928SBarry Smith {
12007bb14e67SBarry Smith   PetscErrorCode ierr;
12014b9ad928SBarry Smith 
12024b9ad928SBarry Smith   PetscFunctionBegin;
12030700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
12047bb14e67SBarry Smith   ierr = PetscUseMethod(pc,"PCASMGetSubKSP_C",(PC,PetscInt*,PetscInt*,KSP **),(pc,n_local,first_local,ksp));CHKERRQ(ierr);
12054b9ad928SBarry Smith   PetscFunctionReturn(0);
12064b9ad928SBarry Smith }
12074b9ad928SBarry Smith 
12084b9ad928SBarry Smith /* -------------------------------------------------------------------------------------*/
12094b9ad928SBarry Smith /*MC
12103b09bd56SBarry Smith    PCASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with
12114b9ad928SBarry Smith            its own KSP object.
12124b9ad928SBarry Smith 
12134b9ad928SBarry Smith    Options Database Keys:
121449517cdeSBarry Smith +  -pc_asm_blocks <blks> - Sets total blocks
12154b9ad928SBarry Smith .  -pc_asm_overlap <ovl> - Sets overlap
1216f1ee410cSBarry Smith .  -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type, default is restrict
1217f1ee410cSBarry Smith -  -pc_asm_local_type [additive, multiplicative] - Sets ASM type, default is additive
12184b9ad928SBarry Smith 
12193b09bd56SBarry Smith      IMPORTANT: If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you
12203b09bd56SBarry Smith       will get a different convergence rate due to the default option of -pc_asm_type restrict. Use
12213b09bd56SBarry Smith       -pc_asm_type basic to use the standard ASM.
12223b09bd56SBarry Smith 
122395452b02SPatrick Sanan    Notes:
122495452b02SPatrick Sanan     Each processor can have one or more blocks, but a block cannot be shared by more
1225f1ee410cSBarry Smith      than one processor. Use PCGASM for subdomains shared by multiple processes. Defaults to one block per processor.
12264b9ad928SBarry Smith 
12273b09bd56SBarry Smith      To set options on the solvers for each block append -sub_ to all the KSP, and PC
1228d7ee0231SBarry Smith         options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly
12294b9ad928SBarry Smith 
1230a8c7a070SBarry Smith      To set the options on the solvers separate for each block call PCASMGetSubKSP()
12314b9ad928SBarry Smith          and set the options directly on the resulting KSP object (you can access its PC
12324b9ad928SBarry Smith          with KSPGetPC())
12334b9ad928SBarry Smith 
12344b9ad928SBarry Smith    Level: beginner
12354b9ad928SBarry Smith 
1236c582cd25SBarry Smith     References:
123796a0c994SBarry Smith +   1. - M Dryja, OB Widlund, An additive variant of the Schwarz alternating method for the case of many subregions
123896a0c994SBarry Smith      Courant Institute, New York University Technical report
12396d33885cSprj- -   2. - Barry Smith, Petter Bjorstad, and William Gropp, Domain Decompositions: Parallel Multilevel Methods for Elliptic Partial Differential Equations,
124096a0c994SBarry Smith     Cambridge University Press.
1241c582cd25SBarry Smith 
12424b9ad928SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
1243f1ee410cSBarry Smith            PCBJACOBI, PCASMGetSubKSP(), PCASMSetLocalSubdomains(), PCASMType, PCASMGetType(), PCASMSetLocalType(), PCASMGetLocalType()
124434a84908Sprj-            PCASMSetTotalSubdomains(), PCSetModifySubMatrices(), PCASMSetOverlap(), PCASMSetType(), PCCompositeType
1245e09e08ccSBarry Smith 
12464b9ad928SBarry Smith M*/
12474b9ad928SBarry Smith 
12488cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_ASM(PC pc)
12494b9ad928SBarry Smith {
1250dfbe8321SBarry Smith   PetscErrorCode ierr;
12514b9ad928SBarry Smith   PC_ASM         *osm;
12524b9ad928SBarry Smith 
12534b9ad928SBarry Smith   PetscFunctionBegin;
1254b00a9115SJed Brown   ierr = PetscNewLog(pc,&osm);CHKERRQ(ierr);
12552fa5cd67SKarl Rupp 
12564b9ad928SBarry Smith   osm->n                 = PETSC_DECIDE;
12574b9ad928SBarry Smith   osm->n_local           = 0;
12582b837212SDmitry Karpeev   osm->n_local_true      = PETSC_DECIDE;
12594b9ad928SBarry Smith   osm->overlap           = 1;
12604b9ad928SBarry Smith   osm->ksp               = 0;
12612b691e39SMatthew Knepley   osm->restriction       = 0;
12625b3c0d42Seaulisa   osm->lprolongation     = 0;
12631dd8081eSeaulisa   osm->lrestriction      = 0;
126492bb6962SLisandro Dalcin   osm->x                 = 0;
126592bb6962SLisandro Dalcin   osm->y                 = 0;
12664b9ad928SBarry Smith   osm->is                = 0;
12672b691e39SMatthew Knepley   osm->is_local          = 0;
12684b9ad928SBarry Smith   osm->mat               = 0;
12694b9ad928SBarry Smith   osm->pmat              = 0;
12704b9ad928SBarry Smith   osm->type              = PC_ASM_RESTRICT;
127112cd4985SMatthew G. Knepley   osm->loctype           = PC_COMPOSITE_ADDITIVE;
12724b9ad928SBarry Smith   osm->same_local_solves = PETSC_TRUE;
12736ed231c7SMatthew Knepley   osm->sort_indices      = PETSC_TRUE;
1274d709ea83SDmitry Karpeev   osm->dm_subdomains     = PETSC_FALSE;
127580ec0b7dSPatrick Sanan   osm->sub_mat_type      = NULL;
12764b9ad928SBarry Smith 
127792bb6962SLisandro Dalcin   pc->data                 = (void*)osm;
12784b9ad928SBarry Smith   pc->ops->apply           = PCApply_ASM;
12794b9ad928SBarry Smith   pc->ops->applytranspose  = PCApplyTranspose_ASM;
12804b9ad928SBarry Smith   pc->ops->setup           = PCSetUp_ASM;
1281e91c6855SBarry Smith   pc->ops->reset           = PCReset_ASM;
12824b9ad928SBarry Smith   pc->ops->destroy         = PCDestroy_ASM;
12834b9ad928SBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_ASM;
12844b9ad928SBarry Smith   pc->ops->setuponblocks   = PCSetUpOnBlocks_ASM;
12854b9ad928SBarry Smith   pc->ops->view            = PCView_ASM;
12864b9ad928SBarry Smith   pc->ops->applyrichardson = 0;
12874b9ad928SBarry Smith 
1288bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetLocalSubdomains_C",PCASMSetLocalSubdomains_ASM);CHKERRQ(ierr);
1289bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetTotalSubdomains_C",PCASMSetTotalSubdomains_ASM);CHKERRQ(ierr);
1290bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetOverlap_C",PCASMSetOverlap_ASM);CHKERRQ(ierr);
1291bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetType_C",PCASMSetType_ASM);CHKERRQ(ierr);
1292c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetType_C",PCASMGetType_ASM);CHKERRQ(ierr);
129312cd4985SMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetLocalType_C",PCASMSetLocalType_ASM);CHKERRQ(ierr);
129412cd4985SMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetLocalType_C",PCASMGetLocalType_ASM);CHKERRQ(ierr);
1295bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetSortIndices_C",PCASMSetSortIndices_ASM);CHKERRQ(ierr);
1296bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetSubKSP_C",PCASMGetSubKSP_ASM);CHKERRQ(ierr);
129780ec0b7dSPatrick Sanan   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetSubMatType_C",PCASMGetSubMatType_ASM);CHKERRQ(ierr);
129880ec0b7dSPatrick Sanan   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetSubMatType_C",PCASMSetSubMatType_ASM);CHKERRQ(ierr);
12994b9ad928SBarry Smith   PetscFunctionReturn(0);
13004b9ad928SBarry Smith }
13014b9ad928SBarry Smith 
130292bb6962SLisandro Dalcin /*@C
130392bb6962SLisandro Dalcin    PCASMCreateSubdomains - Creates the index sets for the overlapping Schwarz
130492bb6962SLisandro Dalcin    preconditioner for a any problem on a general grid.
130592bb6962SLisandro Dalcin 
130692bb6962SLisandro Dalcin    Collective
130792bb6962SLisandro Dalcin 
130892bb6962SLisandro Dalcin    Input Parameters:
130992bb6962SLisandro Dalcin +  A - The global matrix operator
131092bb6962SLisandro Dalcin -  n - the number of local blocks
131192bb6962SLisandro Dalcin 
131292bb6962SLisandro Dalcin    Output Parameters:
131392bb6962SLisandro Dalcin .  outis - the array of index sets defining the subdomains
131492bb6962SLisandro Dalcin 
131592bb6962SLisandro Dalcin    Level: advanced
131692bb6962SLisandro Dalcin 
13177d6bfa3bSBarry Smith    Note: this generates nonoverlapping subdomains; the PCASM will generate the overlap
13187d6bfa3bSBarry Smith     from these if you use PCASMSetLocalSubdomains()
13197d6bfa3bSBarry Smith 
13207d6bfa3bSBarry Smith     In the Fortran version you must provide the array outis[] already allocated of length n.
13217d6bfa3bSBarry Smith 
132292bb6962SLisandro Dalcin .seealso: PCASMSetLocalSubdomains(), PCASMDestroySubdomains()
132392bb6962SLisandro Dalcin @*/
13247087cfbeSBarry Smith PetscErrorCode  PCASMCreateSubdomains(Mat A, PetscInt n, IS* outis[])
132592bb6962SLisandro Dalcin {
132692bb6962SLisandro Dalcin   MatPartitioning mpart;
132792bb6962SLisandro Dalcin   const char      *prefix;
132892bb6962SLisandro Dalcin   PetscInt        i,j,rstart,rend,bs;
1329976e8c5aSLisandro Dalcin   PetscBool       hasop, isbaij = PETSC_FALSE,foundpart = PETSC_FALSE;
13300298fd71SBarry Smith   Mat             Ad     = NULL, adj;
133192bb6962SLisandro Dalcin   IS              ispart,isnumb,*is;
133292bb6962SLisandro Dalcin   PetscErrorCode  ierr;
133392bb6962SLisandro Dalcin 
133492bb6962SLisandro Dalcin   PetscFunctionBegin;
13350700a824SBarry Smith   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
133692bb6962SLisandro Dalcin   PetscValidPointer(outis,3);
1337c1235816SBarry Smith   if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of local blocks must be > 0, n = %D",n);
133892bb6962SLisandro Dalcin 
133992bb6962SLisandro Dalcin   /* Get prefix, row distribution, and block size */
134092bb6962SLisandro Dalcin   ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr);
134192bb6962SLisandro Dalcin   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
134292bb6962SLisandro Dalcin   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
134365e19b50SBarry Smith   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);
134465e19b50SBarry Smith 
134592bb6962SLisandro Dalcin   /* Get diagonal block from matrix if possible */
1346976e8c5aSLisandro Dalcin   ierr = MatHasOperation(A,MATOP_GET_DIAGONAL_BLOCK,&hasop);CHKERRQ(ierr);
1347976e8c5aSLisandro Dalcin   if (hasop) {
134811bd1e4dSLisandro Dalcin     ierr = MatGetDiagonalBlock(A,&Ad);CHKERRQ(ierr);
134992bb6962SLisandro Dalcin   }
135092bb6962SLisandro Dalcin   if (Ad) {
1351b9e7e5c1SBarry Smith     ierr = PetscObjectBaseTypeCompare((PetscObject)Ad,MATSEQBAIJ,&isbaij);CHKERRQ(ierr);
1352b9e7e5c1SBarry Smith     if (!isbaij) {ierr = PetscObjectBaseTypeCompare((PetscObject)Ad,MATSEQSBAIJ,&isbaij);CHKERRQ(ierr);}
135392bb6962SLisandro Dalcin   }
135492bb6962SLisandro Dalcin   if (Ad && n > 1) {
1355ace3abfcSBarry Smith     PetscBool match,done;
135692bb6962SLisandro Dalcin     /* Try to setup a good matrix partitioning if available */
135792bb6962SLisandro Dalcin     ierr = MatPartitioningCreate(PETSC_COMM_SELF,&mpart);CHKERRQ(ierr);
135892bb6962SLisandro Dalcin     ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr);
135992bb6962SLisandro Dalcin     ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr);
1360251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGCURRENT,&match);CHKERRQ(ierr);
136192bb6962SLisandro Dalcin     if (!match) {
1362251f4c67SDmitry Karpeev       ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGSQUARE,&match);CHKERRQ(ierr);
136392bb6962SLisandro Dalcin     }
136492bb6962SLisandro Dalcin     if (!match) { /* assume a "good" partitioner is available */
13651a83f524SJed Brown       PetscInt       na;
13661a83f524SJed Brown       const PetscInt *ia,*ja;
136792bb6962SLisandro Dalcin       ierr = MatGetRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr);
136892bb6962SLisandro Dalcin       if (done) {
136992bb6962SLisandro Dalcin         /* Build adjacency matrix by hand. Unfortunately a call to
137092bb6962SLisandro Dalcin            MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will
137192bb6962SLisandro Dalcin            remove the block-aij structure and we cannot expect
137292bb6962SLisandro Dalcin            MatPartitioning to split vertices as we need */
13731a83f524SJed Brown         PetscInt       i,j,len,nnz,cnt,*iia=0,*jja=0;
13741a83f524SJed Brown         const PetscInt *row;
137592bb6962SLisandro Dalcin         nnz = 0;
137692bb6962SLisandro Dalcin         for (i=0; i<na; i++) { /* count number of nonzeros */
137792bb6962SLisandro Dalcin           len = ia[i+1] - ia[i];
137892bb6962SLisandro Dalcin           row = ja + ia[i];
137992bb6962SLisandro Dalcin           for (j=0; j<len; j++) {
138092bb6962SLisandro Dalcin             if (row[j] == i) { /* don't count diagonal */
138192bb6962SLisandro Dalcin               len--; break;
138292bb6962SLisandro Dalcin             }
138392bb6962SLisandro Dalcin           }
138492bb6962SLisandro Dalcin           nnz += len;
138592bb6962SLisandro Dalcin         }
1386854ce69bSBarry Smith         ierr   = PetscMalloc1(na+1,&iia);CHKERRQ(ierr);
1387854ce69bSBarry Smith         ierr   = PetscMalloc1(nnz,&jja);CHKERRQ(ierr);
138892bb6962SLisandro Dalcin         nnz    = 0;
138992bb6962SLisandro Dalcin         iia[0] = 0;
139092bb6962SLisandro Dalcin         for (i=0; i<na; i++) { /* fill adjacency */
139192bb6962SLisandro Dalcin           cnt = 0;
139292bb6962SLisandro Dalcin           len = ia[i+1] - ia[i];
139392bb6962SLisandro Dalcin           row = ja + ia[i];
139492bb6962SLisandro Dalcin           for (j=0; j<len; j++) {
139592bb6962SLisandro Dalcin             if (row[j] != i) { /* if not diagonal */
139692bb6962SLisandro Dalcin               jja[nnz+cnt++] = row[j];
139792bb6962SLisandro Dalcin             }
139892bb6962SLisandro Dalcin           }
139992bb6962SLisandro Dalcin           nnz     += cnt;
140092bb6962SLisandro Dalcin           iia[i+1] = nnz;
140192bb6962SLisandro Dalcin         }
140292bb6962SLisandro Dalcin         /* Partitioning of the adjacency matrix */
14030298fd71SBarry Smith         ierr      = MatCreateMPIAdj(PETSC_COMM_SELF,na,na,iia,jja,NULL,&adj);CHKERRQ(ierr);
140492bb6962SLisandro Dalcin         ierr      = MatPartitioningSetAdjacency(mpart,adj);CHKERRQ(ierr);
140592bb6962SLisandro Dalcin         ierr      = MatPartitioningSetNParts(mpart,n);CHKERRQ(ierr);
140692bb6962SLisandro Dalcin         ierr      = MatPartitioningApply(mpart,&ispart);CHKERRQ(ierr);
140792bb6962SLisandro Dalcin         ierr      = ISPartitioningToNumbering(ispart,&isnumb);CHKERRQ(ierr);
1408fcfd50ebSBarry Smith         ierr      = MatDestroy(&adj);CHKERRQ(ierr);
140992bb6962SLisandro Dalcin         foundpart = PETSC_TRUE;
141092bb6962SLisandro Dalcin       }
141192bb6962SLisandro Dalcin       ierr = MatRestoreRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr);
141292bb6962SLisandro Dalcin     }
1413fcfd50ebSBarry Smith     ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr);
141492bb6962SLisandro Dalcin   }
141592bb6962SLisandro Dalcin 
1416785e854fSJed Brown   ierr   = PetscMalloc1(n,&is);CHKERRQ(ierr);
141792bb6962SLisandro Dalcin   *outis = is;
141892bb6962SLisandro Dalcin 
141992bb6962SLisandro Dalcin   if (!foundpart) {
142092bb6962SLisandro Dalcin 
142192bb6962SLisandro Dalcin     /* Partitioning by contiguous chunks of rows */
142292bb6962SLisandro Dalcin 
142392bb6962SLisandro Dalcin     PetscInt mbs   = (rend-rstart)/bs;
142492bb6962SLisandro Dalcin     PetscInt start = rstart;
142592bb6962SLisandro Dalcin     for (i=0; i<n; i++) {
142692bb6962SLisandro Dalcin       PetscInt count = (mbs/n + ((mbs % n) > i)) * bs;
142792bb6962SLisandro Dalcin       ierr   = ISCreateStride(PETSC_COMM_SELF,count,start,1,&is[i]);CHKERRQ(ierr);
142892bb6962SLisandro Dalcin       start += count;
142992bb6962SLisandro Dalcin     }
143092bb6962SLisandro Dalcin 
143192bb6962SLisandro Dalcin   } else {
143292bb6962SLisandro Dalcin 
143392bb6962SLisandro Dalcin     /* Partitioning by adjacency of diagonal block  */
143492bb6962SLisandro Dalcin 
143592bb6962SLisandro Dalcin     const PetscInt *numbering;
143692bb6962SLisandro Dalcin     PetscInt       *count,nidx,*indices,*newidx,start=0;
143792bb6962SLisandro Dalcin     /* Get node count in each partition */
1438785e854fSJed Brown     ierr = PetscMalloc1(n,&count);CHKERRQ(ierr);
143992bb6962SLisandro Dalcin     ierr = ISPartitioningCount(ispart,n,count);CHKERRQ(ierr);
144092bb6962SLisandro Dalcin     if (isbaij && bs > 1) { /* adjust for the block-aij case */
144192bb6962SLisandro Dalcin       for (i=0; i<n; i++) count[i] *= bs;
144292bb6962SLisandro Dalcin     }
144392bb6962SLisandro Dalcin     /* Build indices from node numbering */
144492bb6962SLisandro Dalcin     ierr = ISGetLocalSize(isnumb,&nidx);CHKERRQ(ierr);
1445785e854fSJed Brown     ierr = PetscMalloc1(nidx,&indices);CHKERRQ(ierr);
144692bb6962SLisandro Dalcin     for (i=0; i<nidx; i++) indices[i] = i; /* needs to be initialized */
144792bb6962SLisandro Dalcin     ierr = ISGetIndices(isnumb,&numbering);CHKERRQ(ierr);
144892bb6962SLisandro Dalcin     ierr = PetscSortIntWithPermutation(nidx,numbering,indices);CHKERRQ(ierr);
144992bb6962SLisandro Dalcin     ierr = ISRestoreIndices(isnumb,&numbering);CHKERRQ(ierr);
145092bb6962SLisandro Dalcin     if (isbaij && bs > 1) { /* adjust for the block-aij case */
1451785e854fSJed Brown       ierr = PetscMalloc1(nidx*bs,&newidx);CHKERRQ(ierr);
14522fa5cd67SKarl Rupp       for (i=0; i<nidx; i++) {
14532fa5cd67SKarl Rupp         for (j=0; j<bs; j++) newidx[i*bs+j] = indices[i]*bs + j;
14542fa5cd67SKarl Rupp       }
145592bb6962SLisandro Dalcin       ierr    = PetscFree(indices);CHKERRQ(ierr);
145692bb6962SLisandro Dalcin       nidx   *= bs;
145792bb6962SLisandro Dalcin       indices = newidx;
145892bb6962SLisandro Dalcin     }
145992bb6962SLisandro Dalcin     /* Shift to get global indices */
146092bb6962SLisandro Dalcin     for (i=0; i<nidx; i++) indices[i] += rstart;
146192bb6962SLisandro Dalcin 
146292bb6962SLisandro Dalcin     /* Build the index sets for each block */
146392bb6962SLisandro Dalcin     for (i=0; i<n; i++) {
146470b3c8c7SBarry Smith       ierr   = ISCreateGeneral(PETSC_COMM_SELF,count[i],&indices[start],PETSC_COPY_VALUES,&is[i]);CHKERRQ(ierr);
146592bb6962SLisandro Dalcin       ierr   = ISSort(is[i]);CHKERRQ(ierr);
146692bb6962SLisandro Dalcin       start += count[i];
146792bb6962SLisandro Dalcin     }
146892bb6962SLisandro Dalcin 
14693bf036e2SBarry Smith     ierr = PetscFree(count);CHKERRQ(ierr);
14703bf036e2SBarry Smith     ierr = PetscFree(indices);CHKERRQ(ierr);
1471fcfd50ebSBarry Smith     ierr = ISDestroy(&isnumb);CHKERRQ(ierr);
1472fcfd50ebSBarry Smith     ierr = ISDestroy(&ispart);CHKERRQ(ierr);
147392bb6962SLisandro Dalcin 
147492bb6962SLisandro Dalcin   }
147592bb6962SLisandro Dalcin   PetscFunctionReturn(0);
147692bb6962SLisandro Dalcin }
147792bb6962SLisandro Dalcin 
147892bb6962SLisandro Dalcin /*@C
147992bb6962SLisandro Dalcin    PCASMDestroySubdomains - Destroys the index sets created with
148092bb6962SLisandro Dalcin    PCASMCreateSubdomains(). Should be called after setting subdomains
148192bb6962SLisandro Dalcin    with PCASMSetLocalSubdomains().
148292bb6962SLisandro Dalcin 
148392bb6962SLisandro Dalcin    Collective
148492bb6962SLisandro Dalcin 
148592bb6962SLisandro Dalcin    Input Parameters:
148692bb6962SLisandro Dalcin +  n - the number of index sets
14872b691e39SMatthew Knepley .  is - the array of index sets
14880298fd71SBarry Smith -  is_local - the array of local index sets, can be NULL
148992bb6962SLisandro Dalcin 
149092bb6962SLisandro Dalcin    Level: advanced
149192bb6962SLisandro Dalcin 
149292bb6962SLisandro Dalcin .seealso: PCASMCreateSubdomains(), PCASMSetLocalSubdomains()
149392bb6962SLisandro Dalcin @*/
14947087cfbeSBarry Smith PetscErrorCode  PCASMDestroySubdomains(PetscInt n, IS is[], IS is_local[])
149592bb6962SLisandro Dalcin {
149692bb6962SLisandro Dalcin   PetscInt       i;
149792bb6962SLisandro Dalcin   PetscErrorCode ierr;
14985fd66863SKarl Rupp 
149992bb6962SLisandro Dalcin   PetscFunctionBegin;
1500a35b7b57SDmitry Karpeev   if (n <= 0) PetscFunctionReturn(0);
1501a35b7b57SDmitry Karpeev   if (is) {
150292bb6962SLisandro Dalcin     PetscValidPointer(is,2);
1503fcfd50ebSBarry Smith     for (i=0; i<n; i++) { ierr = ISDestroy(&is[i]);CHKERRQ(ierr); }
150492bb6962SLisandro Dalcin     ierr = PetscFree(is);CHKERRQ(ierr);
1505a35b7b57SDmitry Karpeev   }
15062b691e39SMatthew Knepley   if (is_local) {
15072b691e39SMatthew Knepley     PetscValidPointer(is_local,3);
1508fcfd50ebSBarry Smith     for (i=0; i<n; i++) { ierr = ISDestroy(&is_local[i]);CHKERRQ(ierr); }
15092b691e39SMatthew Knepley     ierr = PetscFree(is_local);CHKERRQ(ierr);
15102b691e39SMatthew Knepley   }
151192bb6962SLisandro Dalcin   PetscFunctionReturn(0);
151292bb6962SLisandro Dalcin }
151392bb6962SLisandro Dalcin 
15144b9ad928SBarry Smith /*@
15154b9ad928SBarry Smith    PCASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz
15164b9ad928SBarry Smith    preconditioner for a two-dimensional problem on a regular grid.
15174b9ad928SBarry Smith 
15184b9ad928SBarry Smith    Not Collective
15194b9ad928SBarry Smith 
15204b9ad928SBarry Smith    Input Parameters:
15214b9ad928SBarry Smith +  m, n - the number of mesh points in the x and y directions
15224b9ad928SBarry Smith .  M, N - the number of subdomains in the x and y directions
15234b9ad928SBarry Smith .  dof - degrees of freedom per node
15244b9ad928SBarry Smith -  overlap - overlap in mesh lines
15254b9ad928SBarry Smith 
15264b9ad928SBarry Smith    Output Parameters:
15274b9ad928SBarry Smith +  Nsub - the number of subdomains created
15283d03bb48SJed Brown .  is - array of index sets defining overlapping (if overlap > 0) subdomains
15293d03bb48SJed Brown -  is_local - array of index sets defining non-overlapping subdomains
15304b9ad928SBarry Smith 
15314b9ad928SBarry Smith    Note:
15324b9ad928SBarry Smith    Presently PCAMSCreateSubdomains2d() is valid only for sequential
15334b9ad928SBarry Smith    preconditioners.  More general related routines are
15344b9ad928SBarry Smith    PCASMSetTotalSubdomains() and PCASMSetLocalSubdomains().
15354b9ad928SBarry Smith 
15364b9ad928SBarry Smith    Level: advanced
15374b9ad928SBarry Smith 
15384b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetLocalSubdomains(), PCASMGetSubKSP(),
15394b9ad928SBarry Smith           PCASMSetOverlap()
15404b9ad928SBarry Smith @*/
15417087cfbeSBarry Smith PetscErrorCode  PCASMCreateSubdomains2D(PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt dof,PetscInt overlap,PetscInt *Nsub,IS **is,IS **is_local)
15424b9ad928SBarry Smith {
15433d03bb48SJed Brown   PetscInt       i,j,height,width,ystart,xstart,yleft,yright,xleft,xright,loc_outer;
15446849ba73SBarry Smith   PetscErrorCode ierr;
154513f74950SBarry Smith   PetscInt       nidx,*idx,loc,ii,jj,count;
15464b9ad928SBarry Smith 
15474b9ad928SBarry Smith   PetscFunctionBegin;
1548e32f2f54SBarry Smith   if (dof != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP," ");
15494b9ad928SBarry Smith 
15504b9ad928SBarry Smith   *Nsub     = N*M;
1551854ce69bSBarry Smith   ierr      = PetscMalloc1(*Nsub,is);CHKERRQ(ierr);
1552854ce69bSBarry Smith   ierr      = PetscMalloc1(*Nsub,is_local);CHKERRQ(ierr);
15534b9ad928SBarry Smith   ystart    = 0;
15543d03bb48SJed Brown   loc_outer = 0;
15554b9ad928SBarry Smith   for (i=0; i<N; i++) {
15564b9ad928SBarry Smith     height = n/N + ((n % N) > i); /* height of subdomain */
1557e32f2f54SBarry Smith     if (height < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many N subdomains for mesh dimension n");
15584b9ad928SBarry Smith     yleft  = ystart - overlap; if (yleft < 0) yleft = 0;
15594b9ad928SBarry Smith     yright = ystart + height + overlap; if (yright > n) yright = n;
15604b9ad928SBarry Smith     xstart = 0;
15614b9ad928SBarry Smith     for (j=0; j<M; j++) {
15624b9ad928SBarry Smith       width = m/M + ((m % M) > j); /* width of subdomain */
1563e32f2f54SBarry Smith       if (width < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many M subdomains for mesh dimension m");
15644b9ad928SBarry Smith       xleft  = xstart - overlap; if (xleft < 0) xleft = 0;
15654b9ad928SBarry Smith       xright = xstart + width + overlap; if (xright > m) xright = m;
15664b9ad928SBarry Smith       nidx   = (xright - xleft)*(yright - yleft);
1567785e854fSJed Brown       ierr   = PetscMalloc1(nidx,&idx);CHKERRQ(ierr);
15684b9ad928SBarry Smith       loc    = 0;
15694b9ad928SBarry Smith       for (ii=yleft; ii<yright; ii++) {
15704b9ad928SBarry Smith         count = m*ii + xleft;
15712fa5cd67SKarl Rupp         for (jj=xleft; jj<xright; jj++) idx[loc++] = count++;
15724b9ad928SBarry Smith       }
157370b3c8c7SBarry Smith       ierr = ISCreateGeneral(PETSC_COMM_SELF,nidx,idx,PETSC_COPY_VALUES,(*is)+loc_outer);CHKERRQ(ierr);
15743d03bb48SJed Brown       if (overlap == 0) {
15753d03bb48SJed Brown         ierr = PetscObjectReference((PetscObject)(*is)[loc_outer]);CHKERRQ(ierr);
15762fa5cd67SKarl Rupp 
15773d03bb48SJed Brown         (*is_local)[loc_outer] = (*is)[loc_outer];
15783d03bb48SJed Brown       } else {
15793d03bb48SJed Brown         for (loc=0,ii=ystart; ii<ystart+height; ii++) {
15803d03bb48SJed Brown           for (jj=xstart; jj<xstart+width; jj++) {
15813d03bb48SJed Brown             idx[loc++] = m*ii + jj;
15823d03bb48SJed Brown           }
15833d03bb48SJed Brown         }
158470b3c8c7SBarry Smith         ierr = ISCreateGeneral(PETSC_COMM_SELF,loc,idx,PETSC_COPY_VALUES,*is_local+loc_outer);CHKERRQ(ierr);
15853d03bb48SJed Brown       }
15864b9ad928SBarry Smith       ierr    = PetscFree(idx);CHKERRQ(ierr);
15874b9ad928SBarry Smith       xstart += width;
15883d03bb48SJed Brown       loc_outer++;
15894b9ad928SBarry Smith     }
15904b9ad928SBarry Smith     ystart += height;
15914b9ad928SBarry Smith   }
15924b9ad928SBarry Smith   for (i=0; i<*Nsub; i++) { ierr = ISSort((*is)[i]);CHKERRQ(ierr); }
15934b9ad928SBarry Smith   PetscFunctionReturn(0);
15944b9ad928SBarry Smith }
15954b9ad928SBarry Smith 
15964b9ad928SBarry Smith /*@C
15974b9ad928SBarry Smith     PCASMGetLocalSubdomains - Gets the local subdomains (for this processor
15984b9ad928SBarry Smith     only) for the additive Schwarz preconditioner.
15994b9ad928SBarry Smith 
1600ad4df100SBarry Smith     Not Collective
16014b9ad928SBarry Smith 
16024b9ad928SBarry Smith     Input Parameter:
16034b9ad928SBarry Smith .   pc - the preconditioner context
16044b9ad928SBarry Smith 
16054b9ad928SBarry Smith     Output Parameters:
16064b9ad928SBarry Smith +   n - the number of subdomains for this processor (default value = 1)
16072b691e39SMatthew Knepley .   is - the index sets that define the subdomains for this processor
16080298fd71SBarry Smith -   is_local - the index sets that define the local part of the subdomains for this processor (can be NULL)
16094b9ad928SBarry Smith 
16104b9ad928SBarry Smith 
16114b9ad928SBarry Smith     Notes:
16124b9ad928SBarry Smith     The IS numbering is in the parallel, global numbering of the vector.
16134b9ad928SBarry Smith 
16144b9ad928SBarry Smith     Level: advanced
16154b9ad928SBarry Smith 
16164b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(),
16174b9ad928SBarry Smith           PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubmatrices()
16184b9ad928SBarry Smith @*/
16197087cfbeSBarry Smith PetscErrorCode  PCASMGetLocalSubdomains(PC pc,PetscInt *n,IS *is[],IS *is_local[])
16204b9ad928SBarry Smith {
16212a808120SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
162292bb6962SLisandro Dalcin   PetscErrorCode ierr;
1623ace3abfcSBarry Smith   PetscBool      match;
16244b9ad928SBarry Smith 
16254b9ad928SBarry Smith   PetscFunctionBegin;
16260700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
16274482741eSBarry Smith   PetscValidIntPointer(n,2);
162892bb6962SLisandro Dalcin   if (is) PetscValidPointer(is,3);
1629251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr);
16302a808120SBarry Smith   if (!match) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"PC is not a PCASM");
16314b9ad928SBarry Smith   if (n) *n = osm->n_local_true;
16324b9ad928SBarry Smith   if (is) *is = osm->is;
16332b691e39SMatthew Knepley   if (is_local) *is_local = osm->is_local;
16344b9ad928SBarry Smith   PetscFunctionReturn(0);
16354b9ad928SBarry Smith }
16364b9ad928SBarry Smith 
16374b9ad928SBarry Smith /*@C
16384b9ad928SBarry Smith     PCASMGetLocalSubmatrices - Gets the local submatrices (for this processor
16394b9ad928SBarry Smith     only) for the additive Schwarz preconditioner.
16404b9ad928SBarry Smith 
1641ad4df100SBarry Smith     Not Collective
16424b9ad928SBarry Smith 
16434b9ad928SBarry Smith     Input Parameter:
16444b9ad928SBarry Smith .   pc - the preconditioner context
16454b9ad928SBarry Smith 
16464b9ad928SBarry Smith     Output Parameters:
16474b9ad928SBarry Smith +   n - the number of matrices for this processor (default value = 1)
16484b9ad928SBarry Smith -   mat - the matrices
16494b9ad928SBarry Smith 
16504b9ad928SBarry Smith     Level: advanced
16514b9ad928SBarry Smith 
165295452b02SPatrick Sanan     Notes:
165334a84908Sprj-     Call after PCSetUp() (or KSPSetUp()) but before PCApply() and before PCSetUpOnBlocks())
1654cf739d55SBarry Smith 
165534a84908Sprj-            Usually one would use PCSetModifySubMatrices() to change the submatrices in building the preconditioner.
1656cf739d55SBarry Smith 
16574b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(),
165834a84908Sprj-           PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains(), PCSetModifySubMatrices()
16594b9ad928SBarry Smith @*/
16607087cfbeSBarry Smith PetscErrorCode  PCASMGetLocalSubmatrices(PC pc,PetscInt *n,Mat *mat[])
16614b9ad928SBarry Smith {
16624b9ad928SBarry Smith   PC_ASM         *osm;
166392bb6962SLisandro Dalcin   PetscErrorCode ierr;
1664ace3abfcSBarry Smith   PetscBool      match;
16654b9ad928SBarry Smith 
16664b9ad928SBarry Smith   PetscFunctionBegin;
16670700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
166892bb6962SLisandro Dalcin   PetscValidIntPointer(n,2);
166992bb6962SLisandro Dalcin   if (mat) PetscValidPointer(mat,3);
167034a84908Sprj-   if (!pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call after KSPSetUp() or PCSetUp().");
1671251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr);
167292bb6962SLisandro Dalcin   if (!match) {
167392bb6962SLisandro Dalcin     if (n) *n = 0;
16740298fd71SBarry Smith     if (mat) *mat = NULL;
167592bb6962SLisandro Dalcin   } else {
16764b9ad928SBarry Smith     osm = (PC_ASM*)pc->data;
16774b9ad928SBarry Smith     if (n) *n = osm->n_local_true;
16784b9ad928SBarry Smith     if (mat) *mat = osm->pmat;
167992bb6962SLisandro Dalcin   }
16804b9ad928SBarry Smith   PetscFunctionReturn(0);
16814b9ad928SBarry Smith }
1682d709ea83SDmitry Karpeev 
1683d709ea83SDmitry Karpeev /*@
1684d709ea83SDmitry Karpeev     PCASMSetDMSubdomains - Indicates whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible.
1685f1ee410cSBarry Smith 
1686d709ea83SDmitry Karpeev     Logically Collective
1687d709ea83SDmitry Karpeev 
1688d709ea83SDmitry Karpeev     Input Parameter:
1689d709ea83SDmitry Karpeev +   pc  - the preconditioner
1690d709ea83SDmitry Karpeev -   flg - boolean indicating whether to use subdomains defined by the DM
1691d709ea83SDmitry Karpeev 
1692d709ea83SDmitry Karpeev     Options Database Key:
1693d709ea83SDmitry Karpeev .   -pc_asm_dm_subdomains
1694d709ea83SDmitry Karpeev 
1695d709ea83SDmitry Karpeev     Level: intermediate
1696d709ea83SDmitry Karpeev 
1697d709ea83SDmitry Karpeev     Notes:
1698d709ea83SDmitry Karpeev     PCASMSetTotalSubdomains() and PCASMSetOverlap() take precedence over PCASMSetDMSubdomains(),
1699d709ea83SDmitry Karpeev     so setting either of the first two effectively turns the latter off.
1700d709ea83SDmitry Karpeev 
1701d709ea83SDmitry Karpeev .seealso: PCASMGetDMSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap()
1702d709ea83SDmitry Karpeev           PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains()
1703d709ea83SDmitry Karpeev @*/
1704d709ea83SDmitry Karpeev PetscErrorCode  PCASMSetDMSubdomains(PC pc,PetscBool flg)
1705d709ea83SDmitry Karpeev {
1706d709ea83SDmitry Karpeev   PC_ASM         *osm = (PC_ASM*)pc->data;
1707d709ea83SDmitry Karpeev   PetscErrorCode ierr;
1708d709ea83SDmitry Karpeev   PetscBool      match;
1709d709ea83SDmitry Karpeev 
1710d709ea83SDmitry Karpeev   PetscFunctionBegin;
1711d709ea83SDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1712d709ea83SDmitry Karpeev   PetscValidLogicalCollectiveBool(pc,flg,2);
1713d709ea83SDmitry Karpeev   if (pc->setupcalled) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Not for a setup PC.");
1714d709ea83SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr);
1715d709ea83SDmitry Karpeev   if (match) {
1716d709ea83SDmitry Karpeev     osm->dm_subdomains = flg;
1717d709ea83SDmitry Karpeev   }
1718d709ea83SDmitry Karpeev   PetscFunctionReturn(0);
1719d709ea83SDmitry Karpeev }
1720d709ea83SDmitry Karpeev 
1721d709ea83SDmitry Karpeev /*@
1722d709ea83SDmitry Karpeev     PCASMGetDMSubdomains - Returns flag indicating whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible.
1723d709ea83SDmitry Karpeev     Not Collective
1724d709ea83SDmitry Karpeev 
1725d709ea83SDmitry Karpeev     Input Parameter:
1726d709ea83SDmitry Karpeev .   pc  - the preconditioner
1727d709ea83SDmitry Karpeev 
1728d709ea83SDmitry Karpeev     Output Parameter:
1729d709ea83SDmitry Karpeev .   flg - boolean indicating whether to use subdomains defined by the DM
1730d709ea83SDmitry Karpeev 
1731d709ea83SDmitry Karpeev     Level: intermediate
1732d709ea83SDmitry Karpeev 
1733d709ea83SDmitry Karpeev .seealso: PCASMSetDMSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap()
1734d709ea83SDmitry Karpeev           PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains()
1735d709ea83SDmitry Karpeev @*/
1736d709ea83SDmitry Karpeev PetscErrorCode  PCASMGetDMSubdomains(PC pc,PetscBool* flg)
1737d709ea83SDmitry Karpeev {
1738d709ea83SDmitry Karpeev   PC_ASM         *osm = (PC_ASM*)pc->data;
1739d709ea83SDmitry Karpeev   PetscErrorCode ierr;
1740d709ea83SDmitry Karpeev   PetscBool      match;
1741d709ea83SDmitry Karpeev 
1742d709ea83SDmitry Karpeev   PetscFunctionBegin;
1743d709ea83SDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1744534a8f05SLisandro Dalcin   PetscValidBoolPointer(flg,2);
1745d709ea83SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr);
1746d709ea83SDmitry Karpeev   if (match) {
1747d709ea83SDmitry Karpeev     if (flg) *flg = osm->dm_subdomains;
1748d709ea83SDmitry Karpeev   }
1749d709ea83SDmitry Karpeev   PetscFunctionReturn(0);
1750d709ea83SDmitry Karpeev }
175180ec0b7dSPatrick Sanan 
175280ec0b7dSPatrick Sanan /*@
175380ec0b7dSPatrick Sanan      PCASMGetSubMatType - Gets the matrix type used for ASM subsolves, as a string.
175480ec0b7dSPatrick Sanan 
175580ec0b7dSPatrick Sanan    Not Collective
175680ec0b7dSPatrick Sanan 
175780ec0b7dSPatrick Sanan    Input Parameter:
175880ec0b7dSPatrick Sanan .  pc - the PC
175980ec0b7dSPatrick Sanan 
176080ec0b7dSPatrick Sanan    Output Parameter:
1761f1ee410cSBarry Smith .  -pc_asm_sub_mat_type - name of matrix type
176280ec0b7dSPatrick Sanan 
176380ec0b7dSPatrick Sanan    Level: advanced
176480ec0b7dSPatrick Sanan 
176580ec0b7dSPatrick Sanan .seealso: PCASMSetSubMatType(), PCASM, PCSetType(), VecSetType(), MatType, Mat
176680ec0b7dSPatrick Sanan @*/
176780ec0b7dSPatrick Sanan PetscErrorCode  PCASMGetSubMatType(PC pc,MatType *sub_mat_type){
176880ec0b7dSPatrick Sanan   PetscErrorCode ierr;
176980ec0b7dSPatrick Sanan 
177080ec0b7dSPatrick Sanan   ierr = PetscTryMethod(pc,"PCASMGetSubMatType_C",(PC,MatType*),(pc,sub_mat_type));CHKERRQ(ierr);
177180ec0b7dSPatrick Sanan   PetscFunctionReturn(0);
177280ec0b7dSPatrick Sanan }
177380ec0b7dSPatrick Sanan 
177480ec0b7dSPatrick Sanan /*@
177580ec0b7dSPatrick Sanan      PCASMSetSubMatType - Set the type of matrix used for ASM subsolves
177680ec0b7dSPatrick Sanan 
177780ec0b7dSPatrick Sanan    Collective on Mat
177880ec0b7dSPatrick Sanan 
177980ec0b7dSPatrick Sanan    Input Parameters:
178080ec0b7dSPatrick Sanan +  pc             - the PC object
178180ec0b7dSPatrick Sanan -  sub_mat_type   - matrix type
178280ec0b7dSPatrick Sanan 
178380ec0b7dSPatrick Sanan    Options Database Key:
178480ec0b7dSPatrick Sanan .  -pc_asm_sub_mat_type  <sub_mat_type> - Sets the matrix type used for subsolves, for example, seqaijviennacl. If you specify a base name like aijviennacl, the corresponding sequential type is assumed.
178580ec0b7dSPatrick Sanan 
178680ec0b7dSPatrick Sanan    Notes:
178780ec0b7dSPatrick Sanan    See "${PETSC_DIR}/include/petscmat.h" for available types
178880ec0b7dSPatrick Sanan 
178980ec0b7dSPatrick Sanan   Level: advanced
179080ec0b7dSPatrick Sanan 
179180ec0b7dSPatrick Sanan .seealso: PCASMGetSubMatType(), PCASM, PCSetType(), VecSetType(), MatType, Mat
179280ec0b7dSPatrick Sanan @*/
179380ec0b7dSPatrick Sanan PetscErrorCode PCASMSetSubMatType(PC pc,MatType sub_mat_type)
179480ec0b7dSPatrick Sanan {
179580ec0b7dSPatrick Sanan   PetscErrorCode ierr;
179680ec0b7dSPatrick Sanan 
179780ec0b7dSPatrick Sanan   ierr = PetscTryMethod(pc,"PCASMSetSubMatType_C",(PC,MatType),(pc,sub_mat_type));CHKERRQ(ierr);
179880ec0b7dSPatrick Sanan   PetscFunctionReturn(0);
179980ec0b7dSPatrick Sanan }
1800