xref: /petsc/src/ksp/pc/impls/asm/asm.c (revision 963223946dee96861353a954add71b45b080f26e)
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 
13910cf402Sprj- #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);}
35ffc4695bSBarry Smith     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRMPI(ierr);
3692bb6962SLisandro Dalcin     if (osm->same_local_solves) {
374d219a6aSLisandro Dalcin       if (osm->ksp) {
388d76b567SPatrick Sanan         ierr = PetscViewerASCIIPrintf(viewer,"  Local solver is the same for all blocks, as in the following KSP and PC objects on rank 0:\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;
88ffc4695bSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size);CHKERRMPI(ierr);
89ffc4695bSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank);CHKERRMPI(ierr);
9092bb6962SLisandro Dalcin   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
91589a23caSBarry Smith   ierr = PetscOptionsGetString(NULL,prefix,"-pc_asm_print_subdomains",fname,sizeof(fname),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 
212ffc4695bSBarry Smith       ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRMPI(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);
2821dd8081eSeaulisa 
2834b9ad928SBarry Smith     scall = MAT_INITIAL_MATRIX;
2844b9ad928SBarry Smith   } else {
2854b9ad928SBarry Smith     /*
2864b9ad928SBarry Smith        Destroy the blocks from the previous iteration
2874b9ad928SBarry Smith     */
288e09e08ccSBarry Smith     if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
2894b9ad928SBarry Smith       ierr = MatDestroyMatrices(osm->n_local_true,&osm->pmat);CHKERRQ(ierr);
2904b9ad928SBarry Smith       scall = MAT_INITIAL_MATRIX;
2914b9ad928SBarry Smith     }
2924b9ad928SBarry Smith   }
2934b9ad928SBarry Smith 
294b58cb649SBarry Smith   /* Destroy previous submatrices of a different type than pc->pmat since MAT_REUSE_MATRIX won't work in that case */
295b58cb649SBarry Smith   if ((scall == MAT_REUSE_MATRIX) && osm->sub_mat_type) {
296b58cb649SBarry Smith     if (osm->n_local_true > 0) {
297b58cb649SBarry Smith       ierr = MatDestroySubMatrices(osm->n_local_true,&osm->pmat);CHKERRQ(ierr);
298b58cb649SBarry Smith     }
299b58cb649SBarry Smith     scall = MAT_INITIAL_MATRIX;
300b58cb649SBarry Smith   }
301b58cb649SBarry Smith 
30292bb6962SLisandro Dalcin   /*
30392bb6962SLisandro Dalcin      Extract out the submatrices
30492bb6962SLisandro Dalcin   */
3057dae84e0SHong Zhang   ierr = MatCreateSubMatrices(pc->pmat,osm->n_local_true,osm->is,osm->is,scall,&osm->pmat);CHKERRQ(ierr);
30692bb6962SLisandro Dalcin   if (scall == MAT_INITIAL_MATRIX) {
30792bb6962SLisandro Dalcin     ierr = PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix);CHKERRQ(ierr);
30892bb6962SLisandro Dalcin     for (i=0; i<osm->n_local_true; i++) {
3093bb1ff40SBarry Smith       ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)osm->pmat[i]);CHKERRQ(ierr);
31078904715SLisandro Dalcin       ierr = PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i],pprefix);CHKERRQ(ierr);
31192bb6962SLisandro Dalcin     }
31292bb6962SLisandro Dalcin   }
31380ec0b7dSPatrick Sanan 
31480ec0b7dSPatrick Sanan   /* Convert the types of the submatrices (if needbe) */
31580ec0b7dSPatrick Sanan   if (osm->sub_mat_type) {
31680ec0b7dSPatrick Sanan     for (i=0; i<osm->n_local_true; i++) {
31780ec0b7dSPatrick Sanan       ierr = MatConvert(osm->pmat[i],osm->sub_mat_type,MAT_INPLACE_MATRIX,&(osm->pmat[i]));CHKERRQ(ierr);
31880ec0b7dSPatrick Sanan     }
31980ec0b7dSPatrick Sanan   }
32080ec0b7dSPatrick Sanan 
32180ec0b7dSPatrick Sanan   if (!pc->setupcalled) {
32256ea09ceSStefano Zampini     VecType vtype;
32356ea09ceSStefano Zampini 
32480ec0b7dSPatrick Sanan     /* Create the local work vectors (from the local matrices) and scatter contexts */
3250a545947SLisandro Dalcin     ierr = MatCreateVecs(pc->pmat,&vec,NULL);CHKERRQ(ierr);
3265b3c0d42Seaulisa 
3271dd8081eSeaulisa     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()");
3281dd8081eSeaulisa     if (osm->is_local && osm->type == PC_ASM_RESTRICT && osm->loctype == PC_COMPOSITE_ADDITIVE) {
3291dd8081eSeaulisa       ierr = PetscMalloc1(osm->n_local_true,&osm->lprolongation);CHKERRQ(ierr);
3301dd8081eSeaulisa     }
3311dd8081eSeaulisa     ierr = PetscMalloc1(osm->n_local_true,&osm->lrestriction);CHKERRQ(ierr);
3321dd8081eSeaulisa     ierr = PetscMalloc1(osm->n_local_true,&osm->x);CHKERRQ(ierr);
3331dd8081eSeaulisa     ierr = PetscMalloc1(osm->n_local_true,&osm->y);CHKERRQ(ierr);
334347574c9Seaulisa 
335347574c9Seaulisa     ierr = ISGetLocalSize(osm->lis,&m);CHKERRQ(ierr);
336347574c9Seaulisa     ierr = ISCreateStride(PETSC_COMM_SELF,m,0,1,&isl);CHKERRQ(ierr);
33756ea09ceSStefano Zampini     ierr = MatGetVecType(osm->pmat[0],&vtype);CHKERRQ(ierr);
33856ea09ceSStefano Zampini     ierr = VecCreate(PETSC_COMM_SELF,&osm->lx);CHKERRQ(ierr);
33956ea09ceSStefano Zampini     ierr = VecSetSizes(osm->lx,m,m);CHKERRQ(ierr);
34056ea09ceSStefano Zampini     ierr = VecSetType(osm->lx,vtype);CHKERRQ(ierr);
34156ea09ceSStefano Zampini     ierr = VecDuplicate(osm->lx, &osm->ly);CHKERRQ(ierr);
3429448b7f1SJunchao Zhang     ierr = VecScatterCreate(vec,osm->lis,osm->lx,isl,&osm->restriction);CHKERRQ(ierr);
343347574c9Seaulisa     ierr = ISDestroy(&isl);CHKERRQ(ierr);
344347574c9Seaulisa 
34580ec0b7dSPatrick Sanan     for (i=0; i<osm->n_local_true; ++i) {
3465b3c0d42Seaulisa       ISLocalToGlobalMapping ltog;
3475b3c0d42Seaulisa       IS                     isll;
3485b3c0d42Seaulisa       const PetscInt         *idx_is;
3495b3c0d42Seaulisa       PetscInt               *idx_lis,nout;
3505b3c0d42Seaulisa 
35155817e79SBarry Smith       ierr = ISGetLocalSize(osm->is[i],&m);CHKERRQ(ierr);
35255817e79SBarry Smith       ierr = MatCreateVecs(osm->pmat[i],&osm->x[i],NULL);CHKERRQ(ierr);
35355817e79SBarry Smith       ierr = VecDuplicate(osm->x[i],&osm->y[i]);CHKERRQ(ierr);
35455817e79SBarry Smith 
355b0de9f23SBarry Smith       /* generate a scatter from ly to y[i] picking all the overlapping is[i] entries */
3565b3c0d42Seaulisa       ierr = ISLocalToGlobalMappingCreateIS(osm->lis,&ltog);CHKERRQ(ierr);
3575b3c0d42Seaulisa       ierr = ISGetLocalSize(osm->is[i],&m);CHKERRQ(ierr);
3585b3c0d42Seaulisa       ierr = ISGetIndices(osm->is[i], &idx_is);CHKERRQ(ierr);
3595b3c0d42Seaulisa       ierr = PetscMalloc1(m,&idx_lis);CHKERRQ(ierr);
3605b3c0d42Seaulisa       ierr = ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,m,idx_is,&nout,idx_lis);CHKERRQ(ierr);
3615b3c0d42Seaulisa       if (nout != m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"is not a subset of lis");
3625b3c0d42Seaulisa       ierr = ISRestoreIndices(osm->is[i], &idx_is);CHKERRQ(ierr);
3635b3c0d42Seaulisa       ierr = ISCreateGeneral(PETSC_COMM_SELF,m,idx_lis,PETSC_OWN_POINTER,&isll);CHKERRQ(ierr);
364d8b3b5e3Seaulisa       ierr = ISLocalToGlobalMappingDestroy(&ltog);CHKERRQ(ierr);
3655b3c0d42Seaulisa       ierr = ISCreateStride(PETSC_COMM_SELF,m,0,1,&isl);CHKERRQ(ierr);
3669448b7f1SJunchao Zhang       ierr = VecScatterCreate(osm->ly,isll,osm->y[i],isl,&osm->lrestriction[i]);CHKERRQ(ierr);
3675b3c0d42Seaulisa       ierr = ISDestroy(&isll);CHKERRQ(ierr);
3685b3c0d42Seaulisa       ierr = ISDestroy(&isl);CHKERRQ(ierr);
369910cf402Sprj-       if (osm->lprolongation) { /* generate a scatter from y[i] to ly picking only the the non-overlapping is_local[i] entries */
370d8b3b5e3Seaulisa         ISLocalToGlobalMapping ltog;
371d8b3b5e3Seaulisa         IS                     isll,isll_local;
372d8b3b5e3Seaulisa         const PetscInt         *idx_local;
373d8b3b5e3Seaulisa         PetscInt               *idx1, *idx2, nout;
374d8b3b5e3Seaulisa 
375d8b3b5e3Seaulisa         ierr = ISGetLocalSize(osm->is_local[i],&m_local);CHKERRQ(ierr);
376d8b3b5e3Seaulisa         ierr = ISGetIndices(osm->is_local[i], &idx_local);CHKERRQ(ierr);
377d8b3b5e3Seaulisa 
378d8b3b5e3Seaulisa         ierr = ISLocalToGlobalMappingCreateIS(osm->is[i],&ltog);CHKERRQ(ierr);
379d8b3b5e3Seaulisa         ierr = PetscMalloc1(m_local,&idx1);CHKERRQ(ierr);
380d8b3b5e3Seaulisa         ierr = ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,m_local,idx_local,&nout,idx1);CHKERRQ(ierr);
381d8b3b5e3Seaulisa         ierr = ISLocalToGlobalMappingDestroy(&ltog);CHKERRQ(ierr);
382d8b3b5e3Seaulisa         if (nout != m_local) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"is_local not a subset of is");
383d8b3b5e3Seaulisa         ierr = ISCreateGeneral(PETSC_COMM_SELF,m_local,idx1,PETSC_OWN_POINTER,&isll);CHKERRQ(ierr);
384d8b3b5e3Seaulisa 
385d8b3b5e3Seaulisa         ierr = ISLocalToGlobalMappingCreateIS(osm->lis,&ltog);CHKERRQ(ierr);
386d8b3b5e3Seaulisa         ierr = PetscMalloc1(m_local,&idx2);CHKERRQ(ierr);
387d8b3b5e3Seaulisa         ierr = ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,m_local,idx_local,&nout,idx2);CHKERRQ(ierr);
388d8b3b5e3Seaulisa         ierr = ISLocalToGlobalMappingDestroy(&ltog);CHKERRQ(ierr);
389d8b3b5e3Seaulisa         if (nout != m_local) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"is_local not a subset of lis");
390d8b3b5e3Seaulisa         ierr = ISCreateGeneral(PETSC_COMM_SELF,m_local,idx2,PETSC_OWN_POINTER,&isll_local);CHKERRQ(ierr);
391d8b3b5e3Seaulisa 
392d8b3b5e3Seaulisa         ierr = ISRestoreIndices(osm->is_local[i], &idx_local);CHKERRQ(ierr);
3939448b7f1SJunchao Zhang         ierr = VecScatterCreate(osm->y[i],isll,osm->ly,isll_local,&osm->lprolongation[i]);CHKERRQ(ierr);
394d8b3b5e3Seaulisa 
395d8b3b5e3Seaulisa         ierr = ISDestroy(&isll);CHKERRQ(ierr);
396d8b3b5e3Seaulisa         ierr = ISDestroy(&isll_local);CHKERRQ(ierr);
397d8b3b5e3Seaulisa       }
39880ec0b7dSPatrick Sanan     }
39980ec0b7dSPatrick Sanan     ierr = VecDestroy(&vec);CHKERRQ(ierr);
40080ec0b7dSPatrick Sanan   }
40180ec0b7dSPatrick Sanan 
402fb745f2cSMatthew G. Knepley   if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) {
403235cc4ceSMatthew G. Knepley     IS      *cis;
404235cc4ceSMatthew G. Knepley     PetscInt c;
405235cc4ceSMatthew G. Knepley 
406235cc4ceSMatthew G. Knepley     ierr = PetscMalloc1(osm->n_local_true, &cis);CHKERRQ(ierr);
407235cc4ceSMatthew G. Knepley     for (c = 0; c < osm->n_local_true; ++c) cis[c] = osm->lis;
4087dae84e0SHong Zhang     ierr = MatCreateSubMatrices(pc->pmat, osm->n_local_true, osm->is, cis, scall, &osm->lmats);CHKERRQ(ierr);
409235cc4ceSMatthew G. Knepley     ierr = PetscFree(cis);CHKERRQ(ierr);
410fb745f2cSMatthew G. Knepley   }
4114b9ad928SBarry Smith 
4124b9ad928SBarry Smith   /* Return control to the user so that the submatrices can be modified (e.g., to apply
4134b9ad928SBarry Smith      different boundary conditions for the submatrices than for the global problem) */
41492bb6962SLisandro Dalcin   ierr = PCModifySubMatrices(pc,osm->n_local_true,osm->is,osm->is,osm->pmat,pc->modifysubmatricesP);CHKERRQ(ierr);
4154b9ad928SBarry Smith 
41692bb6962SLisandro Dalcin   /*
41792bb6962SLisandro Dalcin      Loop over subdomains putting them into local ksp
41892bb6962SLisandro Dalcin   */
41983f9b43bSPierre Jolivet   ierr = KSPGetOptionsPrefix(osm->ksp[0],&prefix);CHKERRQ(ierr);
42092bb6962SLisandro Dalcin   for (i=0; i<osm->n_local_true; i++) {
42123ee1639SBarry Smith     ierr = KSPSetOperators(osm->ksp[i],osm->pmat[i],osm->pmat[i]);CHKERRQ(ierr);
42283f9b43bSPierre Jolivet     ierr = MatSetOptionsPrefix(osm->pmat[i],prefix);CHKERRQ(ierr);
42392bb6962SLisandro Dalcin     if (!pc->setupcalled) {
424bf108f30SBarry Smith       ierr = KSPSetFromOptions(osm->ksp[i]);CHKERRQ(ierr);
4254b9ad928SBarry Smith     }
42692bb6962SLisandro Dalcin   }
4274b9ad928SBarry Smith   PetscFunctionReturn(0);
4284b9ad928SBarry Smith }
4294b9ad928SBarry Smith 
4306849ba73SBarry Smith static PetscErrorCode PCSetUpOnBlocks_ASM(PC pc)
4314b9ad928SBarry Smith {
4324b9ad928SBarry Smith   PC_ASM             *osm = (PC_ASM*)pc->data;
4336849ba73SBarry Smith   PetscErrorCode     ierr;
43413f74950SBarry Smith   PetscInt           i;
435539c167fSBarry Smith   KSPConvergedReason reason;
4364b9ad928SBarry Smith 
4374b9ad928SBarry Smith   PetscFunctionBegin;
4384b9ad928SBarry Smith   for (i=0; i<osm->n_local_true; i++) {
4394b9ad928SBarry Smith     ierr = KSPSetUp(osm->ksp[i]);CHKERRQ(ierr);
440539c167fSBarry Smith     ierr = KSPGetConvergedReason(osm->ksp[i],&reason);CHKERRQ(ierr);
441c0decd05SBarry Smith     if (reason == KSP_DIVERGED_PC_FAILED) {
442261222b2SHong Zhang       pc->failedreason = PC_SUBPC_ERROR;
443e0eafd54SHong Zhang     }
4444b9ad928SBarry Smith   }
4454b9ad928SBarry Smith   PetscFunctionReturn(0);
4464b9ad928SBarry Smith }
4474b9ad928SBarry Smith 
4486849ba73SBarry Smith static PetscErrorCode PCApply_ASM(PC pc,Vec x,Vec y)
4494b9ad928SBarry Smith {
4504b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
4516849ba73SBarry Smith   PetscErrorCode ierr;
4521dd8081eSeaulisa   PetscInt       i,n_local_true = osm->n_local_true;
4534b9ad928SBarry Smith   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;
4544b9ad928SBarry Smith 
4554b9ad928SBarry Smith   PetscFunctionBegin;
4564b9ad928SBarry Smith   /*
45748e38a8aSPierre Jolivet      support for limiting the restriction or interpolation to only local
4584b9ad928SBarry Smith      subdomain values (leaving the other values 0).
4594b9ad928SBarry Smith   */
4604b9ad928SBarry Smith   if (!(osm->type & PC_ASM_RESTRICT)) {
4614b9ad928SBarry Smith     forward = SCATTER_FORWARD_LOCAL;
4624b9ad928SBarry Smith     /* have to zero the work RHS since scatter may leave some slots empty */
4631dd8081eSeaulisa     ierr = VecSet(osm->lx, 0.0);CHKERRQ(ierr);
4644b9ad928SBarry Smith   }
465347574c9Seaulisa   if (!(osm->type & PC_ASM_INTERPOLATE)) {
466347574c9Seaulisa     reverse = SCATTER_REVERSE_LOCAL;
467347574c9Seaulisa   }
4684b9ad928SBarry Smith 
469347574c9Seaulisa   if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE || osm->loctype == PC_COMPOSITE_ADDITIVE) {
470b0de9f23SBarry Smith     /* zero the global and the local solutions */
471910cf402Sprj-     ierr = VecSet(y, 0.0);CHKERRQ(ierr);
4725b3c0d42Seaulisa     ierr = VecSet(osm->ly, 0.0);CHKERRQ(ierr);
473347574c9Seaulisa 
47448e38a8aSPierre Jolivet     /* copy the global RHS to local RHS including the ghost nodes */
4751dd8081eSeaulisa     ierr = VecScatterBegin(osm->restriction, x, osm->lx, INSERT_VALUES, forward);CHKERRQ(ierr);
4761dd8081eSeaulisa     ierr = VecScatterEnd(osm->restriction, x, osm->lx, INSERT_VALUES, forward);CHKERRQ(ierr);
477347574c9Seaulisa 
47848e38a8aSPierre Jolivet     /* restrict local RHS to the overlapping 0-block RHS */
4791dd8081eSeaulisa     ierr = VecScatterBegin(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward);CHKERRQ(ierr);
4801dd8081eSeaulisa     ierr = VecScatterEnd(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward);CHKERRQ(ierr);
481d8b3b5e3Seaulisa 
48212cd4985SMatthew G. Knepley     /* do the local solves */
48312cd4985SMatthew G. Knepley     for (i = 0; i < n_local_true; ++i) {
484347574c9Seaulisa 
485b0de9f23SBarry Smith       /* solve the overlapping i-block */
486fa2bb9feSLisandro Dalcin       ierr = PetscLogEventBegin(PC_ApplyOnBlocks, osm->ksp[i], osm->x[i], osm->y[i],0);CHKERRQ(ierr);
48712cd4985SMatthew G. Knepley       ierr = KSPSolve(osm->ksp[i], osm->x[i], osm->y[i]);CHKERRQ(ierr);
488c0decd05SBarry Smith       ierr = KSPCheckSolve(osm->ksp[i], pc, osm->y[i]);CHKERRQ(ierr);
489fa2bb9feSLisandro Dalcin       ierr = PetscLogEventEnd(PC_ApplyOnBlocks, osm->ksp[i], osm->x[i], osm->y[i], 0);CHKERRQ(ierr);
490d8b3b5e3Seaulisa 
491910cf402Sprj-       if (osm->lprolongation) { /* interpolate the non-overlapping i-block solution to the local solution (only for restrictive additive) */
4921dd8081eSeaulisa         ierr = VecScatterBegin(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward);CHKERRQ(ierr);
4931dd8081eSeaulisa         ierr = VecScatterEnd(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward);CHKERRQ(ierr);
49448e38a8aSPierre Jolivet       } else { /* interpolate the overlapping i-block solution to the local solution */
4951dd8081eSeaulisa         ierr = VecScatterBegin(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr);
4961dd8081eSeaulisa         ierr = VecScatterEnd(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr);
497d8b3b5e3Seaulisa       }
498347574c9Seaulisa 
499347574c9Seaulisa       if (i < n_local_true-1) {
50048e38a8aSPierre Jolivet         /* restrict local RHS to the overlapping (i+1)-block RHS */
5011dd8081eSeaulisa         ierr = VecScatterBegin(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward);CHKERRQ(ierr);
5021dd8081eSeaulisa         ierr = VecScatterEnd(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward);CHKERRQ(ierr);
503347574c9Seaulisa 
504347574c9Seaulisa         if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) {
5058966356dSPierre Jolivet           /* update the overlapping (i+1)-block RHS using the current local solution */
506347574c9Seaulisa           ierr = MatMult(osm->lmats[i+1], osm->ly, osm->y[i+1]);CHKERRQ(ierr);
507347574c9Seaulisa           ierr = VecAXPBY(osm->x[i+1],-1.,1., osm->y[i+1]);CHKERRQ(ierr);
5087c3d802fSMatthew G. Knepley         }
50912cd4985SMatthew G. Knepley       }
51012cd4985SMatthew G. Knepley     }
51148e38a8aSPierre Jolivet     /* add the local solution to the global solution including the ghost nodes */
5121dd8081eSeaulisa     ierr = VecScatterBegin(osm->restriction, osm->ly, y, ADD_VALUES, reverse);CHKERRQ(ierr);
5131dd8081eSeaulisa     ierr = VecScatterEnd(osm->restriction, osm->ly, y, ADD_VALUES, reverse);CHKERRQ(ierr);
51489c010cfSBarry Smith   } else SETERRQ1(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Invalid local composition type: %s", PCCompositeTypes[osm->loctype]);
5154b9ad928SBarry Smith   PetscFunctionReturn(0);
5164b9ad928SBarry Smith }
5174b9ad928SBarry Smith 
51848e38a8aSPierre Jolivet static PetscErrorCode PCMatApply_ASM(PC pc,Mat X,Mat Y)
51948e38a8aSPierre Jolivet {
52048e38a8aSPierre Jolivet   PC_ASM         *osm = (PC_ASM*)pc->data;
52148e38a8aSPierre Jolivet   Mat            Z,W;
52248e38a8aSPierre Jolivet   Vec            x;
52348e38a8aSPierre Jolivet   PetscInt       i,m,N;
52448e38a8aSPierre Jolivet   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;
52548e38a8aSPierre Jolivet   PetscErrorCode ierr;
52648e38a8aSPierre Jolivet 
52748e38a8aSPierre Jolivet   PetscFunctionBegin;
52848e38a8aSPierre Jolivet   if (osm->n_local_true > 1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Not yet implemented");
52948e38a8aSPierre Jolivet   /*
53048e38a8aSPierre Jolivet      support for limiting the restriction or interpolation to only local
53148e38a8aSPierre Jolivet      subdomain values (leaving the other values 0).
53248e38a8aSPierre Jolivet   */
53348e38a8aSPierre Jolivet   if (!(osm->type & PC_ASM_RESTRICT)) {
53448e38a8aSPierre Jolivet     forward = SCATTER_FORWARD_LOCAL;
53548e38a8aSPierre Jolivet     /* have to zero the work RHS since scatter may leave some slots empty */
53648e38a8aSPierre Jolivet     ierr = VecSet(osm->lx, 0.0);CHKERRQ(ierr);
53748e38a8aSPierre Jolivet   }
53848e38a8aSPierre Jolivet   if (!(osm->type & PC_ASM_INTERPOLATE)) {
53948e38a8aSPierre Jolivet     reverse = SCATTER_REVERSE_LOCAL;
54048e38a8aSPierre Jolivet   }
54148e38a8aSPierre Jolivet   ierr = VecGetLocalSize(osm->x[0], &m);CHKERRQ(ierr);
54248e38a8aSPierre Jolivet   ierr = MatGetSize(X, NULL, &N);CHKERRQ(ierr);
54348e38a8aSPierre Jolivet   ierr = MatCreateSeqDense(PETSC_COMM_SELF, m, N, NULL, &Z);CHKERRQ(ierr);
54448e38a8aSPierre Jolivet   if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE || osm->loctype == PC_COMPOSITE_ADDITIVE) {
54548e38a8aSPierre Jolivet     /* zero the global and the local solutions */
54648e38a8aSPierre Jolivet     ierr = MatZeroEntries(Y);CHKERRQ(ierr);
54748e38a8aSPierre Jolivet     ierr = VecSet(osm->ly, 0.0);CHKERRQ(ierr);
54848e38a8aSPierre Jolivet 
54948e38a8aSPierre Jolivet     for (i = 0; i < N; ++i) {
55048e38a8aSPierre Jolivet       ierr = MatDenseGetColumnVecRead(X, i, &x);
55148e38a8aSPierre Jolivet       /* copy the global RHS to local RHS including the ghost nodes */
55248e38a8aSPierre Jolivet       ierr = VecScatterBegin(osm->restriction, x, osm->lx, INSERT_VALUES, forward);CHKERRQ(ierr);
55348e38a8aSPierre Jolivet       ierr = VecScatterEnd(osm->restriction, x, osm->lx, INSERT_VALUES, forward);CHKERRQ(ierr);
55448e38a8aSPierre Jolivet       ierr = MatDenseRestoreColumnVecRead(X, i, &x);
55548e38a8aSPierre Jolivet 
55648e38a8aSPierre Jolivet       ierr = MatDenseGetColumnVecWrite(Z, i, &x);
55748e38a8aSPierre Jolivet       /* restrict local RHS to the overlapping 0-block RHS */
55848e38a8aSPierre Jolivet       ierr = VecScatterBegin(osm->lrestriction[0], osm->lx, x, INSERT_VALUES, forward);CHKERRQ(ierr);
55948e38a8aSPierre Jolivet       ierr = VecScatterEnd(osm->lrestriction[0], osm->lx, x, INSERT_VALUES, forward);CHKERRQ(ierr);
56048e38a8aSPierre Jolivet       ierr = MatDenseRestoreColumnVecWrite(Z, i, &x);
56148e38a8aSPierre Jolivet     }
56248e38a8aSPierre Jolivet     ierr = MatCreateSeqDense(PETSC_COMM_SELF, m, N, NULL, &W);CHKERRQ(ierr);
56348e38a8aSPierre Jolivet     /* solve the overlapping 0-block */
56448e38a8aSPierre Jolivet     ierr = PetscLogEventBegin(PC_ApplyOnBlocks, osm->ksp[0], Z, W, 0);CHKERRQ(ierr);
56548e38a8aSPierre Jolivet     ierr = KSPMatSolve(osm->ksp[0], Z, W);CHKERRQ(ierr);
56648e38a8aSPierre Jolivet     ierr = KSPCheckSolve(osm->ksp[0], pc, NULL);CHKERRQ(ierr);
56748e38a8aSPierre Jolivet     ierr = PetscLogEventEnd(PC_ApplyOnBlocks, osm->ksp[0], Z, W,0);CHKERRQ(ierr);
56848e38a8aSPierre Jolivet     ierr = MatDestroy(&Z);CHKERRQ(ierr);
56948e38a8aSPierre Jolivet 
57048e38a8aSPierre Jolivet     for (i = 0; i < N; ++i) {
57148e38a8aSPierre Jolivet       ierr = VecSet(osm->ly, 0.0);CHKERRQ(ierr);
57248e38a8aSPierre Jolivet       ierr = MatDenseGetColumnVecRead(W, i, &x);
57348e38a8aSPierre Jolivet       if (osm->lprolongation) { /* interpolate the non-overlapping 0-block solution to the local solution (only for restrictive additive) */
57448e38a8aSPierre Jolivet         ierr = VecScatterBegin(osm->lprolongation[0], x, osm->ly, ADD_VALUES, forward);CHKERRQ(ierr);
57548e38a8aSPierre Jolivet         ierr = VecScatterEnd(osm->lprolongation[0], x, osm->ly, ADD_VALUES, forward);CHKERRQ(ierr);
57648e38a8aSPierre Jolivet       } else { /* interpolate the overlapping 0-block solution to the local solution */
57748e38a8aSPierre Jolivet         ierr = VecScatterBegin(osm->lrestriction[0], x, osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr);
57848e38a8aSPierre Jolivet         ierr = VecScatterEnd(osm->lrestriction[0], x, osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr);
57948e38a8aSPierre Jolivet       }
58048e38a8aSPierre Jolivet       ierr = MatDenseRestoreColumnVecRead(W, i, &x);
58148e38a8aSPierre Jolivet 
58248e38a8aSPierre Jolivet       ierr = MatDenseGetColumnVecWrite(Y, i, &x);
58348e38a8aSPierre Jolivet       /* add the local solution to the global solution including the ghost nodes */
58448e38a8aSPierre Jolivet       ierr = VecScatterBegin(osm->restriction, osm->ly, x, ADD_VALUES, reverse);CHKERRQ(ierr);
58548e38a8aSPierre Jolivet       ierr = VecScatterEnd(osm->restriction, osm->ly, x, ADD_VALUES, reverse);CHKERRQ(ierr);
58648e38a8aSPierre Jolivet       ierr = MatDenseRestoreColumnVecWrite(Y, i, &x);
58748e38a8aSPierre Jolivet     }
58848e38a8aSPierre Jolivet     ierr = MatDestroy(&W);CHKERRQ(ierr);
58989c010cfSBarry Smith   } else SETERRQ1(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Invalid local composition type: %s", PCCompositeTypes[osm->loctype]);
59048e38a8aSPierre Jolivet   PetscFunctionReturn(0);
59148e38a8aSPierre Jolivet }
59248e38a8aSPierre Jolivet 
5936849ba73SBarry Smith static PetscErrorCode PCApplyTranspose_ASM(PC pc,Vec x,Vec y)
5944b9ad928SBarry Smith {
5954b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
5966849ba73SBarry Smith   PetscErrorCode ierr;
5971dd8081eSeaulisa   PetscInt       i,n_local_true = osm->n_local_true;
5984b9ad928SBarry Smith   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;
5994b9ad928SBarry Smith 
6004b9ad928SBarry Smith   PetscFunctionBegin;
6014b9ad928SBarry Smith   /*
6024b9ad928SBarry Smith      Support for limiting the restriction or interpolation to only local
6034b9ad928SBarry Smith      subdomain values (leaving the other values 0).
6044b9ad928SBarry Smith 
6054b9ad928SBarry Smith      Note: these are reversed from the PCApply_ASM() because we are applying the
6064b9ad928SBarry Smith      transpose of the three terms
6074b9ad928SBarry Smith   */
608d8b3b5e3Seaulisa 
6094b9ad928SBarry Smith   if (!(osm->type & PC_ASM_INTERPOLATE)) {
6104b9ad928SBarry Smith     forward = SCATTER_FORWARD_LOCAL;
6114b9ad928SBarry Smith     /* have to zero the work RHS since scatter may leave some slots empty */
6121dd8081eSeaulisa     ierr = VecSet(osm->lx, 0.0);CHKERRQ(ierr);
6134b9ad928SBarry Smith   }
6142fa5cd67SKarl Rupp   if (!(osm->type & PC_ASM_RESTRICT)) reverse = SCATTER_REVERSE_LOCAL;
6154b9ad928SBarry Smith 
616b0de9f23SBarry Smith   /* zero the global and the local solutions */
617910cf402Sprj-   ierr = VecSet(y, 0.0);CHKERRQ(ierr);
618d8b3b5e3Seaulisa   ierr = VecSet(osm->ly, 0.0);CHKERRQ(ierr);
619d8b3b5e3Seaulisa 
620b0de9f23SBarry Smith   /* Copy the global RHS to local RHS including the ghost nodes */
6211dd8081eSeaulisa   ierr = VecScatterBegin(osm->restriction, x, osm->lx, INSERT_VALUES, forward);CHKERRQ(ierr);
6221dd8081eSeaulisa   ierr = VecScatterEnd(osm->restriction, x, osm->lx, INSERT_VALUES, forward);CHKERRQ(ierr);
623d8b3b5e3Seaulisa 
624b0de9f23SBarry Smith   /* Restrict local RHS to the overlapping 0-block RHS */
6251dd8081eSeaulisa   ierr = VecScatterBegin(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward);CHKERRQ(ierr);
6261dd8081eSeaulisa   ierr = VecScatterEnd(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward);CHKERRQ(ierr);
627d8b3b5e3Seaulisa 
6284b9ad928SBarry Smith   /* do the local solves */
629d8b3b5e3Seaulisa   for (i = 0; i < n_local_true; ++i) {
630d8b3b5e3Seaulisa 
631b0de9f23SBarry Smith     /* solve the overlapping i-block */
632fa2bb9feSLisandro Dalcin     ierr = PetscLogEventBegin(PC_ApplyOnBlocks,osm->ksp[i],osm->x[i],osm->y[i],0);CHKERRQ(ierr);
633e1bcd54cSBarry Smith     ierr = KSPSolveTranspose(osm->ksp[i], osm->x[i], osm->y[i]);CHKERRQ(ierr);
634c0decd05SBarry Smith     ierr = KSPCheckSolve(osm->ksp[i],pc,osm->y[i]);CHKERRQ(ierr);
635fa2bb9feSLisandro Dalcin     ierr = PetscLogEventEnd(PC_ApplyOnBlocks,osm->ksp[i],osm->x[i],osm->y[i],0);CHKERRQ(ierr);
636d8b3b5e3Seaulisa 
637910cf402Sprj-     if (osm->lprolongation) { /* interpolate the non-overlapping i-block solution to the local solution */
6381dd8081eSeaulisa       ierr = VecScatterBegin(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward);CHKERRQ(ierr);
6391dd8081eSeaulisa       ierr = VecScatterEnd(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward);CHKERRQ(ierr);
640910cf402Sprj-     } else { /* interpolate the overlapping i-block solution to the local solution */
6411dd8081eSeaulisa       ierr = VecScatterBegin(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr);
6421dd8081eSeaulisa       ierr = VecScatterEnd(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr);
6434b9ad928SBarry Smith     }
644d8b3b5e3Seaulisa 
645d8b3b5e3Seaulisa     if (i < n_local_true-1) {
646b0de9f23SBarry Smith       /* Restrict local RHS to the overlapping (i+1)-block RHS */
6471dd8081eSeaulisa       ierr = VecScatterBegin(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward);CHKERRQ(ierr);
6481dd8081eSeaulisa       ierr = VecScatterEnd(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward);CHKERRQ(ierr);
6494b9ad928SBarry Smith     }
6504b9ad928SBarry Smith   }
651b0de9f23SBarry Smith   /* Add the local solution to the global solution including the ghost nodes */
6521dd8081eSeaulisa   ierr = VecScatterBegin(osm->restriction, osm->ly, y, ADD_VALUES, reverse);CHKERRQ(ierr);
6531dd8081eSeaulisa   ierr = VecScatterEnd(osm->restriction, osm->ly, y, ADD_VALUES, reverse);CHKERRQ(ierr);
654d8b3b5e3Seaulisa 
6554b9ad928SBarry Smith   PetscFunctionReturn(0);
656d8b3b5e3Seaulisa 
6574b9ad928SBarry Smith }
6584b9ad928SBarry Smith 
659e91c6855SBarry Smith static PetscErrorCode PCReset_ASM(PC pc)
6604b9ad928SBarry Smith {
6614b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
6626849ba73SBarry Smith   PetscErrorCode ierr;
66313f74950SBarry Smith   PetscInt       i;
6644b9ad928SBarry Smith 
6654b9ad928SBarry Smith   PetscFunctionBegin;
66692bb6962SLisandro Dalcin   if (osm->ksp) {
66792bb6962SLisandro Dalcin     for (i=0; i<osm->n_local_true; i++) {
668e91c6855SBarry Smith       ierr = KSPReset(osm->ksp[i]);CHKERRQ(ierr);
66992bb6962SLisandro Dalcin     }
67092bb6962SLisandro Dalcin   }
671e09e08ccSBarry Smith   if (osm->pmat) {
67292bb6962SLisandro Dalcin     if (osm->n_local_true > 0) {
67330a70a9aSHong Zhang       ierr = MatDestroySubMatrices(osm->n_local_true,&osm->pmat);CHKERRQ(ierr);
67492bb6962SLisandro Dalcin     }
67592bb6962SLisandro Dalcin   }
6761dd8081eSeaulisa   if (osm->lrestriction) {
6771dd8081eSeaulisa     ierr = VecScatterDestroy(&osm->restriction);CHKERRQ(ierr);
6781dd8081eSeaulisa     for (i=0; i<osm->n_local_true; i++) {
6791dd8081eSeaulisa       ierr = VecScatterDestroy(&osm->lrestriction[i]);CHKERRQ(ierr);
6801dd8081eSeaulisa       if (osm->lprolongation) {ierr = VecScatterDestroy(&osm->lprolongation[i]);CHKERRQ(ierr);}
681fcfd50ebSBarry Smith       ierr = VecDestroy(&osm->x[i]);CHKERRQ(ierr);
682fcfd50ebSBarry Smith       ierr = VecDestroy(&osm->y[i]);CHKERRQ(ierr);
6834b9ad928SBarry Smith     }
6841dd8081eSeaulisa     ierr = PetscFree(osm->lrestriction);CHKERRQ(ierr);
6851dd8081eSeaulisa     if (osm->lprolongation) {ierr = PetscFree(osm->lprolongation);CHKERRQ(ierr);}
68605b42c5fSBarry Smith     ierr = PetscFree(osm->x);CHKERRQ(ierr);
68778904715SLisandro Dalcin     ierr = PetscFree(osm->y);CHKERRQ(ierr);
6881dd8081eSeaulisa 
68992bb6962SLisandro Dalcin   }
6902b691e39SMatthew Knepley   ierr = PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);CHKERRQ(ierr);
691fb745f2cSMatthew G. Knepley   ierr = ISDestroy(&osm->lis);CHKERRQ(ierr);
692fb745f2cSMatthew G. Knepley   ierr = VecDestroy(&osm->lx);CHKERRQ(ierr);
693fb745f2cSMatthew G. Knepley   ierr = VecDestroy(&osm->ly);CHKERRQ(ierr);
694347574c9Seaulisa   if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) {
695347574c9Seaulisa     ierr = MatDestroyMatrices(osm->n_local_true, &osm->lmats);CHKERRQ(ierr);
696fb745f2cSMatthew G. Knepley   }
6972fa5cd67SKarl Rupp 
69880ec0b7dSPatrick Sanan   ierr = PetscFree(osm->sub_mat_type);CHKERRQ(ierr);
69980ec0b7dSPatrick Sanan 
7000a545947SLisandro Dalcin   osm->is       = NULL;
7010a545947SLisandro Dalcin   osm->is_local = NULL;
702e91c6855SBarry Smith   PetscFunctionReturn(0);
703e91c6855SBarry Smith }
704e91c6855SBarry Smith 
705e91c6855SBarry Smith static PetscErrorCode PCDestroy_ASM(PC pc)
706e91c6855SBarry Smith {
707e91c6855SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
708e91c6855SBarry Smith   PetscErrorCode ierr;
709e91c6855SBarry Smith   PetscInt       i;
710e91c6855SBarry Smith 
711e91c6855SBarry Smith   PetscFunctionBegin;
712e91c6855SBarry Smith   ierr = PCReset_ASM(pc);CHKERRQ(ierr);
713e91c6855SBarry Smith   if (osm->ksp) {
714e91c6855SBarry Smith     for (i=0; i<osm->n_local_true; i++) {
715fcfd50ebSBarry Smith       ierr = KSPDestroy(&osm->ksp[i]);CHKERRQ(ierr);
716e91c6855SBarry Smith     }
717e91c6855SBarry Smith     ierr = PetscFree(osm->ksp);CHKERRQ(ierr);
718e91c6855SBarry Smith   }
719e91c6855SBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
720*96322394SPierre Jolivet 
721*96322394SPierre Jolivet   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetLocalSubdomains_C",NULL);CHKERRQ(ierr);
722*96322394SPierre Jolivet   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetTotalSubdomains_C",NULL);CHKERRQ(ierr);
723*96322394SPierre Jolivet   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetOverlap_C",NULL);CHKERRQ(ierr);
724*96322394SPierre Jolivet   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetType_C",NULL);CHKERRQ(ierr);
725*96322394SPierre Jolivet   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetType_C",NULL);CHKERRQ(ierr);
726*96322394SPierre Jolivet   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetLocalType_C",NULL);CHKERRQ(ierr);
727*96322394SPierre Jolivet   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetLocalType_C",NULL);CHKERRQ(ierr);
728*96322394SPierre Jolivet   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetSortIndices_C",NULL);CHKERRQ(ierr);
729*96322394SPierre Jolivet   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetSubKSP_C",NULL);CHKERRQ(ierr);
730*96322394SPierre Jolivet   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetSubMatType_C",NULL);CHKERRQ(ierr);
731*96322394SPierre Jolivet   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetSubMatType_C",NULL);CHKERRQ(ierr);
7324b9ad928SBarry Smith   PetscFunctionReturn(0);
7334b9ad928SBarry Smith }
7344b9ad928SBarry Smith 
7354416b707SBarry Smith static PetscErrorCode PCSetFromOptions_ASM(PetscOptionItems *PetscOptionsObject,PC pc)
7364b9ad928SBarry Smith {
7374b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
7386849ba73SBarry Smith   PetscErrorCode ierr;
7399dcbbd2bSBarry Smith   PetscInt       blocks,ovl;
74057501b6eSBarry Smith   PetscBool      flg;
74192bb6962SLisandro Dalcin   PCASMType      asmtype;
74212cd4985SMatthew G. Knepley   PCCompositeType loctype;
74380ec0b7dSPatrick Sanan   char           sub_mat_type[256];
7444b9ad928SBarry Smith 
7454b9ad928SBarry Smith   PetscFunctionBegin;
746e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Additive Schwarz options");CHKERRQ(ierr);
747d709ea83SDmitry Karpeev   ierr = PetscOptionsBool("-pc_asm_dm_subdomains","Use DMCreateDomainDecomposition() to define subdomains","PCASMSetDMSubdomains",osm->dm_subdomains,&osm->dm_subdomains,&flg);CHKERRQ(ierr);
74890d69ab7SBarry Smith   ierr = PetscOptionsInt("-pc_asm_blocks","Number of subdomains","PCASMSetTotalSubdomains",osm->n,&blocks,&flg);CHKERRQ(ierr);
74965db9045SDmitry Karpeev   if (flg) {
750f77a5242SKarl Rupp     ierr = PCASMSetTotalSubdomains(pc,blocks,NULL,NULL);CHKERRQ(ierr);
751d709ea83SDmitry Karpeev     osm->dm_subdomains = PETSC_FALSE;
75265db9045SDmitry Karpeev   }
753342c94f9SMatthew G. Knepley   ierr = PetscOptionsInt("-pc_asm_local_blocks","Number of local subdomains","PCASMSetLocalSubdomains",osm->n_local_true,&blocks,&flg);CHKERRQ(ierr);
754342c94f9SMatthew G. Knepley   if (flg) {
755342c94f9SMatthew G. Knepley     ierr = PCASMSetLocalSubdomains(pc,blocks,NULL,NULL);CHKERRQ(ierr);
756342c94f9SMatthew G. Knepley     osm->dm_subdomains = PETSC_FALSE;
757342c94f9SMatthew G. Knepley   }
75890d69ab7SBarry Smith   ierr = PetscOptionsInt("-pc_asm_overlap","Number of grid points overlap","PCASMSetOverlap",osm->overlap,&ovl,&flg);CHKERRQ(ierr);
75965db9045SDmitry Karpeev   if (flg) {
76065db9045SDmitry Karpeev     ierr = PCASMSetOverlap(pc,ovl);CHKERRQ(ierr);
761d709ea83SDmitry Karpeev     osm->dm_subdomains = PETSC_FALSE;
76265db9045SDmitry Karpeev   }
76390d69ab7SBarry Smith   flg  = PETSC_FALSE;
76490d69ab7SBarry Smith   ierr = PetscOptionsEnum("-pc_asm_type","Type of restriction/extension","PCASMSetType",PCASMTypes,(PetscEnum)osm->type,(PetscEnum*)&asmtype,&flg);CHKERRQ(ierr);
76592bb6962SLisandro Dalcin   if (flg) {ierr = PCASMSetType(pc,asmtype);CHKERRQ(ierr); }
76612cd4985SMatthew G. Knepley   flg  = PETSC_FALSE;
76712cd4985SMatthew G. Knepley   ierr = PetscOptionsEnum("-pc_asm_local_type","Type of local solver composition","PCASMSetLocalType",PCCompositeTypes,(PetscEnum)osm->loctype,(PetscEnum*)&loctype,&flg);CHKERRQ(ierr);
76812cd4985SMatthew G. Knepley   if (flg) {ierr = PCASMSetLocalType(pc,loctype);CHKERRQ(ierr); }
76980ec0b7dSPatrick Sanan   ierr = PetscOptionsFList("-pc_asm_sub_mat_type","Subsolve Matrix Type","PCASMSetSubMatType",MatList,NULL,sub_mat_type,256,&flg);CHKERRQ(ierr);
77080ec0b7dSPatrick Sanan   if (flg) {
771459726d8SSatish Balay     ierr = PCASMSetSubMatType(pc,sub_mat_type);CHKERRQ(ierr);
77280ec0b7dSPatrick Sanan   }
7734b9ad928SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
7744b9ad928SBarry Smith   PetscFunctionReturn(0);
7754b9ad928SBarry Smith }
7764b9ad928SBarry Smith 
7774b9ad928SBarry Smith /*------------------------------------------------------------------------------------*/
7784b9ad928SBarry Smith 
7791e6b0712SBarry Smith static PetscErrorCode  PCASMSetLocalSubdomains_ASM(PC pc,PetscInt n,IS is[],IS is_local[])
7804b9ad928SBarry Smith {
7814b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
78292bb6962SLisandro Dalcin   PetscErrorCode ierr;
78392bb6962SLisandro Dalcin   PetscInt       i;
7844b9ad928SBarry Smith 
7854b9ad928SBarry Smith   PetscFunctionBegin;
786e32f2f54SBarry Smith   if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Each process must have 1 or more blocks, n = %D",n);
787ce94432eSBarry 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().");
788e7e72b3dSBarry Smith 
7894b9ad928SBarry Smith   if (!pc->setupcalled) {
79092bb6962SLisandro Dalcin     if (is) {
79192bb6962SLisandro Dalcin       for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)is[i]);CHKERRQ(ierr);}
79292bb6962SLisandro Dalcin     }
793832fc9a5SMatthew Knepley     if (is_local) {
794832fc9a5SMatthew Knepley       for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)is_local[i]);CHKERRQ(ierr);}
795832fc9a5SMatthew Knepley     }
7962b691e39SMatthew Knepley     ierr = PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);CHKERRQ(ierr);
7972fa5cd67SKarl Rupp 
7984b9ad928SBarry Smith     osm->n_local_true = n;
7990a545947SLisandro Dalcin     osm->is           = NULL;
8000a545947SLisandro Dalcin     osm->is_local     = NULL;
80192bb6962SLisandro Dalcin     if (is) {
802785e854fSJed Brown       ierr = PetscMalloc1(n,&osm->is);CHKERRQ(ierr);
8032fa5cd67SKarl Rupp       for (i=0; i<n; i++) osm->is[i] = is[i];
8043d03bb48SJed Brown       /* Flag indicating that the user has set overlapping subdomains so PCASM should not increase their size. */
8053d03bb48SJed Brown       osm->overlap = -1;
80692bb6962SLisandro Dalcin     }
8072b691e39SMatthew Knepley     if (is_local) {
808785e854fSJed Brown       ierr = PetscMalloc1(n,&osm->is_local);CHKERRQ(ierr);
8092fa5cd67SKarl Rupp       for (i=0; i<n; i++) osm->is_local[i] = is_local[i];
810a35b7b57SDmitry Karpeev       if (!is) {
811785e854fSJed Brown         ierr = PetscMalloc1(osm->n_local_true,&osm->is);CHKERRQ(ierr);
812a35b7b57SDmitry Karpeev         for (i=0; i<osm->n_local_true; i++) {
813a35b7b57SDmitry Karpeev           if (osm->overlap > 0) { /* With positive overlap, osm->is[i] will be modified */
814a35b7b57SDmitry Karpeev             ierr = ISDuplicate(osm->is_local[i],&osm->is[i]);CHKERRQ(ierr);
815a35b7b57SDmitry Karpeev             ierr = ISCopy(osm->is_local[i],osm->is[i]);CHKERRQ(ierr);
816a35b7b57SDmitry Karpeev           } else {
817a35b7b57SDmitry Karpeev             ierr = PetscObjectReference((PetscObject)osm->is_local[i]);CHKERRQ(ierr);
818a35b7b57SDmitry Karpeev             osm->is[i] = osm->is_local[i];
819a35b7b57SDmitry Karpeev           }
820a35b7b57SDmitry Karpeev         }
821a35b7b57SDmitry Karpeev       }
8222b691e39SMatthew Knepley     }
8234b9ad928SBarry Smith   }
8244b9ad928SBarry Smith   PetscFunctionReturn(0);
8254b9ad928SBarry Smith }
8264b9ad928SBarry Smith 
8271e6b0712SBarry Smith static PetscErrorCode  PCASMSetTotalSubdomains_ASM(PC pc,PetscInt N,IS *is,IS *is_local)
8284b9ad928SBarry Smith {
8294b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
8306849ba73SBarry Smith   PetscErrorCode ierr;
83113f74950SBarry Smith   PetscMPIInt    rank,size;
83278904715SLisandro Dalcin   PetscInt       n;
8334b9ad928SBarry Smith 
8344b9ad928SBarry Smith   PetscFunctionBegin;
835ce94432eSBarry Smith   if (N < 1) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Number of total blocks must be > 0, N = %D",N);
836ce94432eSBarry 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.");
8374b9ad928SBarry Smith 
8384b9ad928SBarry Smith   /*
839880770e9SJed Brown      Split the subdomains equally among all processors
8404b9ad928SBarry Smith   */
841ffc4695bSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRMPI(ierr);
842ffc4695bSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRMPI(ierr);
8434b9ad928SBarry Smith   n    = N/size + ((N % size) > rank);
844e32f2f54SBarry 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);
845e32f2f54SBarry Smith   if (pc->setupcalled && n != osm->n_local_true) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"PCASMSetTotalSubdomains() should be called before PCSetUp().");
8464b9ad928SBarry Smith   if (!pc->setupcalled) {
8472b691e39SMatthew Knepley     ierr = PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);CHKERRQ(ierr);
8482fa5cd67SKarl Rupp 
8494b9ad928SBarry Smith     osm->n_local_true = n;
8500a545947SLisandro Dalcin     osm->is           = NULL;
8510a545947SLisandro Dalcin     osm->is_local     = NULL;
8524b9ad928SBarry Smith   }
8534b9ad928SBarry Smith   PetscFunctionReturn(0);
8544b9ad928SBarry Smith }
8554b9ad928SBarry Smith 
8561e6b0712SBarry Smith static PetscErrorCode  PCASMSetOverlap_ASM(PC pc,PetscInt ovl)
8574b9ad928SBarry Smith {
85892bb6962SLisandro Dalcin   PC_ASM *osm = (PC_ASM*)pc->data;
8594b9ad928SBarry Smith 
8604b9ad928SBarry Smith   PetscFunctionBegin;
861ce94432eSBarry Smith   if (ovl < 0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested");
862ce94432eSBarry Smith   if (pc->setupcalled && ovl != osm->overlap) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCASMSetOverlap() should be called before PCSetUp().");
8632fa5cd67SKarl Rupp   if (!pc->setupcalled) osm->overlap = ovl;
8644b9ad928SBarry Smith   PetscFunctionReturn(0);
8654b9ad928SBarry Smith }
8664b9ad928SBarry Smith 
8671e6b0712SBarry Smith static PetscErrorCode  PCASMSetType_ASM(PC pc,PCASMType type)
8684b9ad928SBarry Smith {
86992bb6962SLisandro Dalcin   PC_ASM *osm = (PC_ASM*)pc->data;
8704b9ad928SBarry Smith 
8714b9ad928SBarry Smith   PetscFunctionBegin;
8724b9ad928SBarry Smith   osm->type     = type;
873bf108f30SBarry Smith   osm->type_set = PETSC_TRUE;
8744b9ad928SBarry Smith   PetscFunctionReturn(0);
8754b9ad928SBarry Smith }
8764b9ad928SBarry Smith 
877c60c7ad4SBarry Smith static PetscErrorCode  PCASMGetType_ASM(PC pc,PCASMType *type)
878c60c7ad4SBarry Smith {
879c60c7ad4SBarry Smith   PC_ASM *osm = (PC_ASM*)pc->data;
880c60c7ad4SBarry Smith 
881c60c7ad4SBarry Smith   PetscFunctionBegin;
882c60c7ad4SBarry Smith   *type = osm->type;
883c60c7ad4SBarry Smith   PetscFunctionReturn(0);
884c60c7ad4SBarry Smith }
885c60c7ad4SBarry Smith 
88612cd4985SMatthew G. Knepley static PetscErrorCode  PCASMSetLocalType_ASM(PC pc, PCCompositeType type)
88712cd4985SMatthew G. Knepley {
88812cd4985SMatthew G. Knepley   PC_ASM *osm = (PC_ASM *) pc->data;
88912cd4985SMatthew G. Knepley 
89012cd4985SMatthew G. Knepley   PetscFunctionBegin;
891d0ecd4dfSBarry 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");
89212cd4985SMatthew G. Knepley   osm->loctype = type;
89312cd4985SMatthew G. Knepley   PetscFunctionReturn(0);
89412cd4985SMatthew G. Knepley }
89512cd4985SMatthew G. Knepley 
89612cd4985SMatthew G. Knepley static PetscErrorCode  PCASMGetLocalType_ASM(PC pc, PCCompositeType *type)
89712cd4985SMatthew G. Knepley {
89812cd4985SMatthew G. Knepley   PC_ASM *osm = (PC_ASM *) pc->data;
89912cd4985SMatthew G. Knepley 
90012cd4985SMatthew G. Knepley   PetscFunctionBegin;
90112cd4985SMatthew G. Knepley   *type = osm->loctype;
90212cd4985SMatthew G. Knepley   PetscFunctionReturn(0);
90312cd4985SMatthew G. Knepley }
90412cd4985SMatthew G. Knepley 
9051e6b0712SBarry Smith static PetscErrorCode  PCASMSetSortIndices_ASM(PC pc,PetscBool  doSort)
9066ed231c7SMatthew Knepley {
9076ed231c7SMatthew Knepley   PC_ASM *osm = (PC_ASM*)pc->data;
9086ed231c7SMatthew Knepley 
9096ed231c7SMatthew Knepley   PetscFunctionBegin;
9106ed231c7SMatthew Knepley   osm->sort_indices = doSort;
9116ed231c7SMatthew Knepley   PetscFunctionReturn(0);
9126ed231c7SMatthew Knepley }
9136ed231c7SMatthew Knepley 
9141e6b0712SBarry Smith static PetscErrorCode  PCASMGetSubKSP_ASM(PC pc,PetscInt *n_local,PetscInt *first_local,KSP **ksp)
9154b9ad928SBarry Smith {
91692bb6962SLisandro Dalcin   PC_ASM         *osm = (PC_ASM*)pc->data;
917dfbe8321SBarry Smith   PetscErrorCode ierr;
9184b9ad928SBarry Smith 
9194b9ad928SBarry Smith   PetscFunctionBegin;
92034a84908Sprj-   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");
9214b9ad928SBarry Smith 
9222fa5cd67SKarl Rupp   if (n_local) *n_local = osm->n_local_true;
92392bb6962SLisandro Dalcin   if (first_local) {
924ce94432eSBarry Smith     ierr          = MPI_Scan(&osm->n_local_true,first_local,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
92592bb6962SLisandro Dalcin     *first_local -= osm->n_local_true;
92692bb6962SLisandro Dalcin   }
92792bb6962SLisandro Dalcin   if (ksp) {
92892bb6962SLisandro Dalcin     /* Assume that local solves are now different; not necessarily
92992bb6962SLisandro Dalcin        true though!  This flag is used only for PCView_ASM() */
93092bb6962SLisandro Dalcin     *ksp                   = osm->ksp;
93192bb6962SLisandro Dalcin     osm->same_local_solves = PETSC_FALSE;
93292bb6962SLisandro Dalcin   }
9334b9ad928SBarry Smith   PetscFunctionReturn(0);
9344b9ad928SBarry Smith }
9354b9ad928SBarry Smith 
93680ec0b7dSPatrick Sanan static PetscErrorCode  PCASMGetSubMatType_ASM(PC pc,MatType *sub_mat_type)
93780ec0b7dSPatrick Sanan {
93880ec0b7dSPatrick Sanan   PC_ASM         *osm = (PC_ASM*)pc->data;
93980ec0b7dSPatrick Sanan 
94080ec0b7dSPatrick Sanan   PetscFunctionBegin;
94180ec0b7dSPatrick Sanan   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
94280ec0b7dSPatrick Sanan   PetscValidPointer(sub_mat_type,2);
94380ec0b7dSPatrick Sanan   *sub_mat_type = osm->sub_mat_type;
94480ec0b7dSPatrick Sanan   PetscFunctionReturn(0);
94580ec0b7dSPatrick Sanan }
94680ec0b7dSPatrick Sanan 
94780ec0b7dSPatrick Sanan static PetscErrorCode PCASMSetSubMatType_ASM(PC pc,MatType sub_mat_type)
94880ec0b7dSPatrick Sanan {
94980ec0b7dSPatrick Sanan   PetscErrorCode    ierr;
95080ec0b7dSPatrick Sanan   PC_ASM            *osm = (PC_ASM*)pc->data;
95180ec0b7dSPatrick Sanan 
95280ec0b7dSPatrick Sanan   PetscFunctionBegin;
95380ec0b7dSPatrick Sanan   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
95480ec0b7dSPatrick Sanan   ierr = PetscFree(osm->sub_mat_type);CHKERRQ(ierr);
95580ec0b7dSPatrick Sanan   ierr = PetscStrallocpy(sub_mat_type,(char**)&osm->sub_mat_type);CHKERRQ(ierr);
95680ec0b7dSPatrick Sanan   PetscFunctionReturn(0);
95780ec0b7dSPatrick Sanan }
95880ec0b7dSPatrick Sanan 
9594b9ad928SBarry Smith /*@C
9601093a601SBarry Smith     PCASMSetLocalSubdomains - Sets the local subdomains (for this processor only) for the additive Schwarz preconditioner.
9614b9ad928SBarry Smith 
962d083f849SBarry Smith     Collective on pc
9634b9ad928SBarry Smith 
9644b9ad928SBarry Smith     Input Parameters:
9654b9ad928SBarry Smith +   pc - the preconditioner context
9664b9ad928SBarry Smith .   n - the number of subdomains for this processor (default value = 1)
9678c03b21aSDmitry Karpeev .   is - the index set that defines the subdomains for this processor
9680298fd71SBarry Smith          (or NULL for PETSc to determine subdomains)
969f1ee410cSBarry 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
970f1ee410cSBarry Smith          (or NULL to not provide these)
9714b9ad928SBarry Smith 
972342c94f9SMatthew G. Knepley     Options Database Key:
973342c94f9SMatthew G. Knepley     To set the total number of subdomain blocks rather than specify the
974342c94f9SMatthew G. Knepley     index sets, use the option
975342c94f9SMatthew G. Knepley .    -pc_asm_local_blocks <blks> - Sets local blocks
976342c94f9SMatthew G. Knepley 
9774b9ad928SBarry Smith     Notes:
9781093a601SBarry Smith     The IS numbering is in the parallel, global numbering of the vector for both is and is_local
9794b9ad928SBarry Smith 
9804b9ad928SBarry Smith     By default the ASM preconditioner uses 1 block per processor.
9814b9ad928SBarry Smith 
9824b9ad928SBarry Smith     Use PCASMSetTotalSubdomains() to set the subdomains for all processors.
9834b9ad928SBarry Smith 
984f1ee410cSBarry Smith     If is_local is provided and PCASMType is PC_ASM_RESTRICT then the solution only over the is_local region is interpolated
985f1ee410cSBarry Smith     back to form the global solution (this is the standard restricted additive Schwarz method)
986f1ee410cSBarry Smith 
987f1ee410cSBarry 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
988f1ee410cSBarry Smith     no code to handle that case.
989f1ee410cSBarry Smith 
9904b9ad928SBarry Smith     Level: advanced
9914b9ad928SBarry Smith 
9924b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(),
993f1ee410cSBarry Smith           PCASMCreateSubdomains2D(), PCASMGetLocalSubdomains(), PCASMType, PCASMSetType()
9944b9ad928SBarry Smith @*/
9957087cfbeSBarry Smith PetscErrorCode  PCASMSetLocalSubdomains(PC pc,PetscInt n,IS is[],IS is_local[])
9964b9ad928SBarry Smith {
9977bb14e67SBarry Smith   PetscErrorCode ierr;
9984b9ad928SBarry Smith 
9994b9ad928SBarry Smith   PetscFunctionBegin;
10000700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
10017bb14e67SBarry Smith   ierr = PetscTryMethod(pc,"PCASMSetLocalSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,n,is,is_local));CHKERRQ(ierr);
10024b9ad928SBarry Smith   PetscFunctionReturn(0);
10034b9ad928SBarry Smith }
10044b9ad928SBarry Smith 
10054b9ad928SBarry Smith /*@C
1006feb221f8SDmitry Karpeev     PCASMSetTotalSubdomains - Sets the subdomains for all processors for the
10074b9ad928SBarry Smith     additive Schwarz preconditioner.  Either all or no processors in the
10084b9ad928SBarry Smith     PC communicator must call this routine, with the same index sets.
10094b9ad928SBarry Smith 
1010d083f849SBarry Smith     Collective on pc
10114b9ad928SBarry Smith 
10124b9ad928SBarry Smith     Input Parameters:
10134b9ad928SBarry Smith +   pc - the preconditioner context
1014feb221f8SDmitry Karpeev .   N  - the number of subdomains for all processors
1015feb221f8SDmitry Karpeev .   is - the index sets that define the subdomains for all processors
1016dfaa0c5fSBarry Smith          (or NULL to ask PETSc to determine the subdomains)
10172b691e39SMatthew Knepley -   is_local - the index sets that define the local part of the subdomains for this processor
1018f1ee410cSBarry Smith          (or NULL to not provide this information)
10194b9ad928SBarry Smith 
10204b9ad928SBarry Smith     Options Database Key:
10214b9ad928SBarry Smith     To set the total number of subdomain blocks rather than specify the
10224b9ad928SBarry Smith     index sets, use the option
10234b9ad928SBarry Smith .    -pc_asm_blocks <blks> - Sets total blocks
10244b9ad928SBarry Smith 
10254b9ad928SBarry Smith     Notes:
1026f1ee410cSBarry Smith     Currently you cannot use this to set the actual subdomains with the argument is or is_local.
10274b9ad928SBarry Smith 
10284b9ad928SBarry Smith     By default the ASM preconditioner uses 1 block per processor.
10294b9ad928SBarry Smith 
10304b9ad928SBarry Smith     These index sets cannot be destroyed until after completion of the
10314b9ad928SBarry Smith     linear solves for which the ASM preconditioner is being used.
10324b9ad928SBarry Smith 
10334b9ad928SBarry Smith     Use PCASMSetLocalSubdomains() to set local subdomains.
10344b9ad928SBarry Smith 
10351093a601SBarry Smith     The IS numbering is in the parallel, global numbering of the vector for both is and is_local
10361093a601SBarry Smith 
10374b9ad928SBarry Smith     Level: advanced
10384b9ad928SBarry Smith 
10394b9ad928SBarry Smith .seealso: PCASMSetLocalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(),
10404b9ad928SBarry Smith           PCASMCreateSubdomains2D()
10414b9ad928SBarry Smith @*/
10427087cfbeSBarry Smith PetscErrorCode  PCASMSetTotalSubdomains(PC pc,PetscInt N,IS is[],IS is_local[])
10434b9ad928SBarry Smith {
10447bb14e67SBarry Smith   PetscErrorCode ierr;
10454b9ad928SBarry Smith 
10464b9ad928SBarry Smith   PetscFunctionBegin;
10470700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
10487bb14e67SBarry Smith   ierr = PetscTryMethod(pc,"PCASMSetTotalSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,N,is,is_local));CHKERRQ(ierr);
10494b9ad928SBarry Smith   PetscFunctionReturn(0);
10504b9ad928SBarry Smith }
10514b9ad928SBarry Smith 
10524b9ad928SBarry Smith /*@
10534b9ad928SBarry Smith     PCASMSetOverlap - Sets the overlap between a pair of subdomains for the
10544b9ad928SBarry Smith     additive Schwarz preconditioner.  Either all or no processors in the
1055f1ee410cSBarry Smith     PC communicator must call this routine.
10564b9ad928SBarry Smith 
1057d083f849SBarry Smith     Logically Collective on pc
10584b9ad928SBarry Smith 
10594b9ad928SBarry Smith     Input Parameters:
10604b9ad928SBarry Smith +   pc  - the preconditioner context
10614b9ad928SBarry Smith -   ovl - the amount of overlap between subdomains (ovl >= 0, default value = 1)
10624b9ad928SBarry Smith 
10634b9ad928SBarry Smith     Options Database Key:
10644b9ad928SBarry Smith .   -pc_asm_overlap <ovl> - Sets overlap
10654b9ad928SBarry Smith 
10664b9ad928SBarry Smith     Notes:
10674b9ad928SBarry Smith     By default the ASM preconditioner uses 1 block per processor.  To use
10684b9ad928SBarry Smith     multiple blocks per perocessor, see PCASMSetTotalSubdomains() and
10694b9ad928SBarry Smith     PCASMSetLocalSubdomains() (and the option -pc_asm_blocks <blks>).
10704b9ad928SBarry Smith 
10714b9ad928SBarry Smith     The overlap defaults to 1, so if one desires that no additional
10724b9ad928SBarry Smith     overlap be computed beyond what may have been set with a call to
10734b9ad928SBarry Smith     PCASMSetTotalSubdomains() or PCASMSetLocalSubdomains(), then ovl
10744b9ad928SBarry Smith     must be set to be 0.  In particular, if one does not explicitly set
10754b9ad928SBarry Smith     the subdomains an application code, then all overlap would be computed
10764b9ad928SBarry Smith     internally by PETSc, and using an overlap of 0 would result in an ASM
10774b9ad928SBarry Smith     variant that is equivalent to the block Jacobi preconditioner.
10784b9ad928SBarry Smith 
1079f1ee410cSBarry Smith     The default algorithm used by PETSc to increase overlap is fast, but not scalable,
1080f1ee410cSBarry Smith     use the option -mat_increase_overlap_scalable when the problem and number of processes is large.
1081f1ee410cSBarry Smith 
10824b9ad928SBarry Smith     Note that one can define initial index sets with any overlap via
1083f1ee410cSBarry Smith     PCASMSetLocalSubdomains(); the routine
10844b9ad928SBarry Smith     PCASMSetOverlap() merely allows PETSc to extend that overlap further
10854b9ad928SBarry Smith     if desired.
10864b9ad928SBarry Smith 
10874b9ad928SBarry Smith     Level: intermediate
10884b9ad928SBarry Smith 
10894b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetLocalSubdomains(), PCASMGetSubKSP(),
1090f1ee410cSBarry Smith           PCASMCreateSubdomains2D(), PCASMGetLocalSubdomains(), MatIncreaseOverlap()
10914b9ad928SBarry Smith @*/
10927087cfbeSBarry Smith PetscErrorCode  PCASMSetOverlap(PC pc,PetscInt ovl)
10934b9ad928SBarry Smith {
10947bb14e67SBarry Smith   PetscErrorCode ierr;
10954b9ad928SBarry Smith 
10964b9ad928SBarry Smith   PetscFunctionBegin;
10970700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1098c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,ovl,2);
10997bb14e67SBarry Smith   ierr = PetscTryMethod(pc,"PCASMSetOverlap_C",(PC,PetscInt),(pc,ovl));CHKERRQ(ierr);
11004b9ad928SBarry Smith   PetscFunctionReturn(0);
11014b9ad928SBarry Smith }
11024b9ad928SBarry Smith 
11034b9ad928SBarry Smith /*@
11044b9ad928SBarry Smith     PCASMSetType - Sets the type of restriction and interpolation used
11054b9ad928SBarry Smith     for local problems in the additive Schwarz method.
11064b9ad928SBarry Smith 
1107d083f849SBarry Smith     Logically Collective on pc
11084b9ad928SBarry Smith 
11094b9ad928SBarry Smith     Input Parameters:
11104b9ad928SBarry Smith +   pc  - the preconditioner context
11114b9ad928SBarry Smith -   type - variant of ASM, one of
11124b9ad928SBarry Smith .vb
11134b9ad928SBarry Smith       PC_ASM_BASIC       - full interpolation and restriction
11144b9ad928SBarry Smith       PC_ASM_RESTRICT    - full restriction, local processor interpolation
11154b9ad928SBarry Smith       PC_ASM_INTERPOLATE - full interpolation, local processor restriction
11164b9ad928SBarry Smith       PC_ASM_NONE        - local processor restriction and interpolation
11174b9ad928SBarry Smith .ve
11184b9ad928SBarry Smith 
11194b9ad928SBarry Smith     Options Database Key:
11204b9ad928SBarry Smith .   -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type
11214b9ad928SBarry Smith 
112295452b02SPatrick Sanan     Notes:
112395452b02SPatrick Sanan     if the is_local arguments are passed to PCASMSetLocalSubdomains() then they are used when PC_ASM_RESTRICT has been selected
1124f1ee410cSBarry Smith     to limit the local processor interpolation
1125f1ee410cSBarry Smith 
11264b9ad928SBarry Smith     Level: intermediate
11274b9ad928SBarry Smith 
11284b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(),
1129f1ee410cSBarry Smith           PCASMCreateSubdomains2D(), PCASMType, PCASMSetLocalType(), PCASMGetLocalType()
11304b9ad928SBarry Smith @*/
11317087cfbeSBarry Smith PetscErrorCode  PCASMSetType(PC pc,PCASMType type)
11324b9ad928SBarry Smith {
11337bb14e67SBarry Smith   PetscErrorCode ierr;
11344b9ad928SBarry Smith 
11354b9ad928SBarry Smith   PetscFunctionBegin;
11360700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1137c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,type,2);
11387bb14e67SBarry Smith   ierr = PetscTryMethod(pc,"PCASMSetType_C",(PC,PCASMType),(pc,type));CHKERRQ(ierr);
11394b9ad928SBarry Smith   PetscFunctionReturn(0);
11404b9ad928SBarry Smith }
11414b9ad928SBarry Smith 
1142c60c7ad4SBarry Smith /*@
1143c60c7ad4SBarry Smith     PCASMGetType - Gets the type of restriction and interpolation used
1144c60c7ad4SBarry Smith     for local problems in the additive Schwarz method.
1145c60c7ad4SBarry Smith 
1146d083f849SBarry Smith     Logically Collective on pc
1147c60c7ad4SBarry Smith 
1148c60c7ad4SBarry Smith     Input Parameter:
1149c60c7ad4SBarry Smith .   pc  - the preconditioner context
1150c60c7ad4SBarry Smith 
1151c60c7ad4SBarry Smith     Output Parameter:
1152c60c7ad4SBarry Smith .   type - variant of ASM, one of
1153c60c7ad4SBarry Smith 
1154c60c7ad4SBarry Smith .vb
1155c60c7ad4SBarry Smith       PC_ASM_BASIC       - full interpolation and restriction
1156c60c7ad4SBarry Smith       PC_ASM_RESTRICT    - full restriction, local processor interpolation
1157c60c7ad4SBarry Smith       PC_ASM_INTERPOLATE - full interpolation, local processor restriction
1158c60c7ad4SBarry Smith       PC_ASM_NONE        - local processor restriction and interpolation
1159c60c7ad4SBarry Smith .ve
1160c60c7ad4SBarry Smith 
1161c60c7ad4SBarry Smith     Options Database Key:
1162c60c7ad4SBarry Smith .   -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type
1163c60c7ad4SBarry Smith 
1164c60c7ad4SBarry Smith     Level: intermediate
1165c60c7ad4SBarry Smith 
1166c60c7ad4SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(),
1167f1ee410cSBarry Smith           PCASMCreateSubdomains2D(), PCASMType, PCASMSetType(), PCASMSetLocalType(), PCASMGetLocalType()
1168c60c7ad4SBarry Smith @*/
1169c60c7ad4SBarry Smith PetscErrorCode  PCASMGetType(PC pc,PCASMType *type)
1170c60c7ad4SBarry Smith {
1171c60c7ad4SBarry Smith   PetscErrorCode ierr;
1172c60c7ad4SBarry Smith 
1173c60c7ad4SBarry Smith   PetscFunctionBegin;
1174c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1175c60c7ad4SBarry Smith   ierr = PetscUseMethod(pc,"PCASMGetType_C",(PC,PCASMType*),(pc,type));CHKERRQ(ierr);
1176c60c7ad4SBarry Smith   PetscFunctionReturn(0);
1177c60c7ad4SBarry Smith }
1178c60c7ad4SBarry Smith 
117912cd4985SMatthew G. Knepley /*@
118012cd4985SMatthew G. Knepley   PCASMSetLocalType - Sets the type of composition used for local problems in the additive Schwarz method.
118112cd4985SMatthew G. Knepley 
1182d083f849SBarry Smith   Logically Collective on pc
118312cd4985SMatthew G. Knepley 
118412cd4985SMatthew G. Knepley   Input Parameters:
118512cd4985SMatthew G. Knepley + pc  - the preconditioner context
118612cd4985SMatthew G. Knepley - type - type of composition, one of
118712cd4985SMatthew G. Knepley .vb
118812cd4985SMatthew G. Knepley   PC_COMPOSITE_ADDITIVE       - local additive combination
118912cd4985SMatthew G. Knepley   PC_COMPOSITE_MULTIPLICATIVE - local multiplicative combination
119012cd4985SMatthew G. Knepley .ve
119112cd4985SMatthew G. Knepley 
119212cd4985SMatthew G. Knepley   Options Database Key:
119312cd4985SMatthew G. Knepley . -pc_asm_local_type [additive,multiplicative] - Sets local solver composition type
119412cd4985SMatthew G. Knepley 
119512cd4985SMatthew G. Knepley   Level: intermediate
119612cd4985SMatthew G. Knepley 
1197f1ee410cSBarry Smith .seealso: PCASMSetType(), PCASMGetType(), PCASMGetLocalType(), PCASM, PCASMType, PCASMSetType(), PCASMGetType(), PCCompositeType
119812cd4985SMatthew G. Knepley @*/
119912cd4985SMatthew G. Knepley PetscErrorCode PCASMSetLocalType(PC pc, PCCompositeType type)
120012cd4985SMatthew G. Knepley {
120112cd4985SMatthew G. Knepley   PetscErrorCode ierr;
120212cd4985SMatthew G. Knepley 
120312cd4985SMatthew G. Knepley   PetscFunctionBegin;
120412cd4985SMatthew G. Knepley   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
120512cd4985SMatthew G. Knepley   PetscValidLogicalCollectiveEnum(pc, type, 2);
120612cd4985SMatthew G. Knepley   ierr = PetscTryMethod(pc, "PCASMSetLocalType_C", (PC, PCCompositeType), (pc, type));CHKERRQ(ierr);
120712cd4985SMatthew G. Knepley   PetscFunctionReturn(0);
120812cd4985SMatthew G. Knepley }
120912cd4985SMatthew G. Knepley 
121012cd4985SMatthew G. Knepley /*@
121112cd4985SMatthew G. Knepley   PCASMGetLocalType - Gets the type of composition used for local problems in the additive Schwarz method.
121212cd4985SMatthew G. Knepley 
1213d083f849SBarry Smith   Logically Collective on pc
121412cd4985SMatthew G. Knepley 
121512cd4985SMatthew G. Knepley   Input Parameter:
121612cd4985SMatthew G. Knepley . pc  - the preconditioner context
121712cd4985SMatthew G. Knepley 
121812cd4985SMatthew G. Knepley   Output Parameter:
121912cd4985SMatthew G. Knepley . type - type of composition, one of
122012cd4985SMatthew G. Knepley .vb
122112cd4985SMatthew G. Knepley   PC_COMPOSITE_ADDITIVE       - local additive combination
122212cd4985SMatthew G. Knepley   PC_COMPOSITE_MULTIPLICATIVE - local multiplicative combination
122312cd4985SMatthew G. Knepley .ve
122412cd4985SMatthew G. Knepley 
122512cd4985SMatthew G. Knepley   Options Database Key:
122612cd4985SMatthew G. Knepley . -pc_asm_local_type [additive,multiplicative] - Sets local solver composition type
122712cd4985SMatthew G. Knepley 
122812cd4985SMatthew G. Knepley   Level: intermediate
122912cd4985SMatthew G. Knepley 
1230f1ee410cSBarry Smith .seealso: PCASMSetType(), PCASMGetType(), PCASMSetLocalType(), PCASMCreate(), PCASMType, PCASMSetType(), PCASMGetType(), PCCompositeType
123112cd4985SMatthew G. Knepley @*/
123212cd4985SMatthew G. Knepley PetscErrorCode PCASMGetLocalType(PC pc, PCCompositeType *type)
123312cd4985SMatthew G. Knepley {
123412cd4985SMatthew G. Knepley   PetscErrorCode ierr;
123512cd4985SMatthew G. Knepley 
123612cd4985SMatthew G. Knepley   PetscFunctionBegin;
123712cd4985SMatthew G. Knepley   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
123812cd4985SMatthew G. Knepley   PetscValidPointer(type, 2);
123912cd4985SMatthew G. Knepley   ierr = PetscUseMethod(pc, "PCASMGetLocalType_C", (PC, PCCompositeType *), (pc, type));CHKERRQ(ierr);
124012cd4985SMatthew G. Knepley   PetscFunctionReturn(0);
124112cd4985SMatthew G. Knepley }
124212cd4985SMatthew G. Knepley 
12436ed231c7SMatthew Knepley /*@
12446ed231c7SMatthew Knepley     PCASMSetSortIndices - Determines whether subdomain indices are sorted.
12456ed231c7SMatthew Knepley 
1246d083f849SBarry Smith     Logically Collective on pc
12476ed231c7SMatthew Knepley 
12486ed231c7SMatthew Knepley     Input Parameters:
12496ed231c7SMatthew Knepley +   pc  - the preconditioner context
12506ed231c7SMatthew Knepley -   doSort - sort the subdomain indices
12516ed231c7SMatthew Knepley 
12526ed231c7SMatthew Knepley     Level: intermediate
12536ed231c7SMatthew Knepley 
12546ed231c7SMatthew Knepley .seealso: PCASMSetLocalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(),
12556ed231c7SMatthew Knepley           PCASMCreateSubdomains2D()
12566ed231c7SMatthew Knepley @*/
12577087cfbeSBarry Smith PetscErrorCode  PCASMSetSortIndices(PC pc,PetscBool doSort)
12586ed231c7SMatthew Knepley {
12597bb14e67SBarry Smith   PetscErrorCode ierr;
12606ed231c7SMatthew Knepley 
12616ed231c7SMatthew Knepley   PetscFunctionBegin;
12620700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1263acfcf0e5SJed Brown   PetscValidLogicalCollectiveBool(pc,doSort,2);
12647bb14e67SBarry Smith   ierr = PetscTryMethod(pc,"PCASMSetSortIndices_C",(PC,PetscBool),(pc,doSort));CHKERRQ(ierr);
12656ed231c7SMatthew Knepley   PetscFunctionReturn(0);
12666ed231c7SMatthew Knepley }
12676ed231c7SMatthew Knepley 
12684b9ad928SBarry Smith /*@C
12694b9ad928SBarry Smith    PCASMGetSubKSP - Gets the local KSP contexts for all blocks on
12704b9ad928SBarry Smith    this processor.
12714b9ad928SBarry Smith 
1272d083f849SBarry Smith    Collective on pc iff first_local is requested
12734b9ad928SBarry Smith 
12744b9ad928SBarry Smith    Input Parameter:
12754b9ad928SBarry Smith .  pc - the preconditioner context
12764b9ad928SBarry Smith 
12774b9ad928SBarry Smith    Output Parameters:
12780298fd71SBarry Smith +  n_local - the number of blocks on this processor or NULL
12790298fd71SBarry Smith .  first_local - the global number of the first block on this processor or NULL,
12800298fd71SBarry Smith                  all processors must request or all must pass NULL
12814b9ad928SBarry Smith -  ksp - the array of KSP contexts
12824b9ad928SBarry Smith 
12834b9ad928SBarry Smith    Note:
1284d29017ddSJed Brown    After PCASMGetSubKSP() the array of KSPes is not to be freed.
12854b9ad928SBarry Smith 
12864b9ad928SBarry Smith    You must call KSPSetUp() before calling PCASMGetSubKSP().
12874b9ad928SBarry Smith 
1288d29017ddSJed Brown    Fortran note:
12892bf68e3eSBarry 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.
1290d29017ddSJed Brown 
12914b9ad928SBarry Smith    Level: advanced
12924b9ad928SBarry Smith 
12934b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap(),
12944b9ad928SBarry Smith           PCASMCreateSubdomains2D(),
12954b9ad928SBarry Smith @*/
12967087cfbeSBarry Smith PetscErrorCode  PCASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[])
12974b9ad928SBarry Smith {
12987bb14e67SBarry Smith   PetscErrorCode ierr;
12994b9ad928SBarry Smith 
13004b9ad928SBarry Smith   PetscFunctionBegin;
13010700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
13027bb14e67SBarry Smith   ierr = PetscUseMethod(pc,"PCASMGetSubKSP_C",(PC,PetscInt*,PetscInt*,KSP **),(pc,n_local,first_local,ksp));CHKERRQ(ierr);
13034b9ad928SBarry Smith   PetscFunctionReturn(0);
13044b9ad928SBarry Smith }
13054b9ad928SBarry Smith 
13064b9ad928SBarry Smith /* -------------------------------------------------------------------------------------*/
13074b9ad928SBarry Smith /*MC
13083b09bd56SBarry Smith    PCASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with
13094b9ad928SBarry Smith            its own KSP object.
13104b9ad928SBarry Smith 
13114b9ad928SBarry Smith    Options Database Keys:
131249517cdeSBarry Smith +  -pc_asm_blocks <blks> - Sets total blocks
13134b9ad928SBarry Smith .  -pc_asm_overlap <ovl> - Sets overlap
1314f1ee410cSBarry Smith .  -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type, default is restrict
1315f1ee410cSBarry Smith -  -pc_asm_local_type [additive, multiplicative] - Sets ASM type, default is additive
13164b9ad928SBarry Smith 
13173b09bd56SBarry Smith      IMPORTANT: If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you
13183b09bd56SBarry Smith       will get a different convergence rate due to the default option of -pc_asm_type restrict. Use
13193b09bd56SBarry Smith       -pc_asm_type basic to use the standard ASM.
13203b09bd56SBarry Smith 
132195452b02SPatrick Sanan    Notes:
132295452b02SPatrick Sanan     Each processor can have one or more blocks, but a block cannot be shared by more
1323f1ee410cSBarry Smith      than one processor. Use PCGASM for subdomains shared by multiple processes. Defaults to one block per processor.
13244b9ad928SBarry Smith 
13253b09bd56SBarry Smith      To set options on the solvers for each block append -sub_ to all the KSP, and PC
1326d7ee0231SBarry Smith         options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly
13274b9ad928SBarry Smith 
1328a8c7a070SBarry Smith      To set the options on the solvers separate for each block call PCASMGetSubKSP()
13294b9ad928SBarry Smith          and set the options directly on the resulting KSP object (you can access its PC
13304b9ad928SBarry Smith          with KSPGetPC())
13314b9ad928SBarry Smith 
13324b9ad928SBarry Smith    Level: beginner
13334b9ad928SBarry Smith 
1334c582cd25SBarry Smith     References:
133596a0c994SBarry Smith +   1. - M Dryja, OB Widlund, An additive variant of the Schwarz alternating method for the case of many subregions
133696a0c994SBarry Smith      Courant Institute, New York University Technical report
13376d33885cSprj- -   2. - Barry Smith, Petter Bjorstad, and William Gropp, Domain Decompositions: Parallel Multilevel Methods for Elliptic Partial Differential Equations,
133896a0c994SBarry Smith     Cambridge University Press.
1339c582cd25SBarry Smith 
13404b9ad928SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
1341f1ee410cSBarry Smith            PCBJACOBI, PCASMGetSubKSP(), PCASMSetLocalSubdomains(), PCASMType, PCASMGetType(), PCASMSetLocalType(), PCASMGetLocalType()
134234a84908Sprj-            PCASMSetTotalSubdomains(), PCSetModifySubMatrices(), PCASMSetOverlap(), PCASMSetType(), PCCompositeType
1343e09e08ccSBarry Smith 
13444b9ad928SBarry Smith M*/
13454b9ad928SBarry Smith 
13468cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_ASM(PC pc)
13474b9ad928SBarry Smith {
1348dfbe8321SBarry Smith   PetscErrorCode ierr;
13494b9ad928SBarry Smith   PC_ASM         *osm;
13504b9ad928SBarry Smith 
13514b9ad928SBarry Smith   PetscFunctionBegin;
1352b00a9115SJed Brown   ierr = PetscNewLog(pc,&osm);CHKERRQ(ierr);
13532fa5cd67SKarl Rupp 
13544b9ad928SBarry Smith   osm->n                 = PETSC_DECIDE;
13554b9ad928SBarry Smith   osm->n_local           = 0;
13562b837212SDmitry Karpeev   osm->n_local_true      = PETSC_DECIDE;
13574b9ad928SBarry Smith   osm->overlap           = 1;
13580a545947SLisandro Dalcin   osm->ksp               = NULL;
13590a545947SLisandro Dalcin   osm->restriction       = NULL;
13600a545947SLisandro Dalcin   osm->lprolongation     = NULL;
13610a545947SLisandro Dalcin   osm->lrestriction      = NULL;
13620a545947SLisandro Dalcin   osm->x                 = NULL;
13630a545947SLisandro Dalcin   osm->y                 = NULL;
13640a545947SLisandro Dalcin   osm->is                = NULL;
13650a545947SLisandro Dalcin   osm->is_local          = NULL;
13660a545947SLisandro Dalcin   osm->mat               = NULL;
13670a545947SLisandro Dalcin   osm->pmat              = NULL;
13684b9ad928SBarry Smith   osm->type              = PC_ASM_RESTRICT;
136912cd4985SMatthew G. Knepley   osm->loctype           = PC_COMPOSITE_ADDITIVE;
13704b9ad928SBarry Smith   osm->same_local_solves = PETSC_TRUE;
13716ed231c7SMatthew Knepley   osm->sort_indices      = PETSC_TRUE;
1372d709ea83SDmitry Karpeev   osm->dm_subdomains     = PETSC_FALSE;
137380ec0b7dSPatrick Sanan   osm->sub_mat_type      = NULL;
13744b9ad928SBarry Smith 
137592bb6962SLisandro Dalcin   pc->data                 = (void*)osm;
13764b9ad928SBarry Smith   pc->ops->apply           = PCApply_ASM;
137748e38a8aSPierre Jolivet   pc->ops->matapply        = PCMatApply_ASM;
13784b9ad928SBarry Smith   pc->ops->applytranspose  = PCApplyTranspose_ASM;
13794b9ad928SBarry Smith   pc->ops->setup           = PCSetUp_ASM;
1380e91c6855SBarry Smith   pc->ops->reset           = PCReset_ASM;
13814b9ad928SBarry Smith   pc->ops->destroy         = PCDestroy_ASM;
13824b9ad928SBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_ASM;
13834b9ad928SBarry Smith   pc->ops->setuponblocks   = PCSetUpOnBlocks_ASM;
13844b9ad928SBarry Smith   pc->ops->view            = PCView_ASM;
13850a545947SLisandro Dalcin   pc->ops->applyrichardson = NULL;
13864b9ad928SBarry Smith 
1387bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetLocalSubdomains_C",PCASMSetLocalSubdomains_ASM);CHKERRQ(ierr);
1388bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetTotalSubdomains_C",PCASMSetTotalSubdomains_ASM);CHKERRQ(ierr);
1389bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetOverlap_C",PCASMSetOverlap_ASM);CHKERRQ(ierr);
1390bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetType_C",PCASMSetType_ASM);CHKERRQ(ierr);
1391c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetType_C",PCASMGetType_ASM);CHKERRQ(ierr);
139212cd4985SMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetLocalType_C",PCASMSetLocalType_ASM);CHKERRQ(ierr);
139312cd4985SMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetLocalType_C",PCASMGetLocalType_ASM);CHKERRQ(ierr);
1394bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetSortIndices_C",PCASMSetSortIndices_ASM);CHKERRQ(ierr);
1395bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetSubKSP_C",PCASMGetSubKSP_ASM);CHKERRQ(ierr);
139680ec0b7dSPatrick Sanan   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetSubMatType_C",PCASMGetSubMatType_ASM);CHKERRQ(ierr);
139780ec0b7dSPatrick Sanan   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetSubMatType_C",PCASMSetSubMatType_ASM);CHKERRQ(ierr);
13984b9ad928SBarry Smith   PetscFunctionReturn(0);
13994b9ad928SBarry Smith }
14004b9ad928SBarry Smith 
140192bb6962SLisandro Dalcin /*@C
140292bb6962SLisandro Dalcin    PCASMCreateSubdomains - Creates the index sets for the overlapping Schwarz
140392bb6962SLisandro Dalcin    preconditioner for a any problem on a general grid.
140492bb6962SLisandro Dalcin 
140592bb6962SLisandro Dalcin    Collective
140692bb6962SLisandro Dalcin 
140792bb6962SLisandro Dalcin    Input Parameters:
140892bb6962SLisandro Dalcin +  A - The global matrix operator
140992bb6962SLisandro Dalcin -  n - the number of local blocks
141092bb6962SLisandro Dalcin 
141192bb6962SLisandro Dalcin    Output Parameters:
141292bb6962SLisandro Dalcin .  outis - the array of index sets defining the subdomains
141392bb6962SLisandro Dalcin 
141492bb6962SLisandro Dalcin    Level: advanced
141592bb6962SLisandro Dalcin 
14167d6bfa3bSBarry Smith    Note: this generates nonoverlapping subdomains; the PCASM will generate the overlap
14177d6bfa3bSBarry Smith     from these if you use PCASMSetLocalSubdomains()
14187d6bfa3bSBarry Smith 
14197d6bfa3bSBarry Smith     In the Fortran version you must provide the array outis[] already allocated of length n.
14207d6bfa3bSBarry Smith 
142192bb6962SLisandro Dalcin .seealso: PCASMSetLocalSubdomains(), PCASMDestroySubdomains()
142292bb6962SLisandro Dalcin @*/
14237087cfbeSBarry Smith PetscErrorCode  PCASMCreateSubdomains(Mat A, PetscInt n, IS* outis[])
142492bb6962SLisandro Dalcin {
142592bb6962SLisandro Dalcin   MatPartitioning mpart;
142692bb6962SLisandro Dalcin   const char      *prefix;
142792bb6962SLisandro Dalcin   PetscInt        i,j,rstart,rend,bs;
1428976e8c5aSLisandro Dalcin   PetscBool       hasop, isbaij = PETSC_FALSE,foundpart = PETSC_FALSE;
14290298fd71SBarry Smith   Mat             Ad     = NULL, adj;
143092bb6962SLisandro Dalcin   IS              ispart,isnumb,*is;
143192bb6962SLisandro Dalcin   PetscErrorCode  ierr;
143292bb6962SLisandro Dalcin 
143392bb6962SLisandro Dalcin   PetscFunctionBegin;
14340700a824SBarry Smith   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
143592bb6962SLisandro Dalcin   PetscValidPointer(outis,3);
1436c1235816SBarry Smith   if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of local blocks must be > 0, n = %D",n);
143792bb6962SLisandro Dalcin 
143892bb6962SLisandro Dalcin   /* Get prefix, row distribution, and block size */
143992bb6962SLisandro Dalcin   ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr);
144092bb6962SLisandro Dalcin   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
144192bb6962SLisandro Dalcin   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
144265e19b50SBarry 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);
144365e19b50SBarry Smith 
144492bb6962SLisandro Dalcin   /* Get diagonal block from matrix if possible */
1445976e8c5aSLisandro Dalcin   ierr = MatHasOperation(A,MATOP_GET_DIAGONAL_BLOCK,&hasop);CHKERRQ(ierr);
1446976e8c5aSLisandro Dalcin   if (hasop) {
144711bd1e4dSLisandro Dalcin     ierr = MatGetDiagonalBlock(A,&Ad);CHKERRQ(ierr);
144892bb6962SLisandro Dalcin   }
144992bb6962SLisandro Dalcin   if (Ad) {
1450b9e7e5c1SBarry Smith     ierr = PetscObjectBaseTypeCompare((PetscObject)Ad,MATSEQBAIJ,&isbaij);CHKERRQ(ierr);
1451b9e7e5c1SBarry Smith     if (!isbaij) {ierr = PetscObjectBaseTypeCompare((PetscObject)Ad,MATSEQSBAIJ,&isbaij);CHKERRQ(ierr);}
145292bb6962SLisandro Dalcin   }
145392bb6962SLisandro Dalcin   if (Ad && n > 1) {
1454ace3abfcSBarry Smith     PetscBool match,done;
145592bb6962SLisandro Dalcin     /* Try to setup a good matrix partitioning if available */
145692bb6962SLisandro Dalcin     ierr = MatPartitioningCreate(PETSC_COMM_SELF,&mpart);CHKERRQ(ierr);
145792bb6962SLisandro Dalcin     ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr);
145892bb6962SLisandro Dalcin     ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr);
1459251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGCURRENT,&match);CHKERRQ(ierr);
146092bb6962SLisandro Dalcin     if (!match) {
1461251f4c67SDmitry Karpeev       ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGSQUARE,&match);CHKERRQ(ierr);
146292bb6962SLisandro Dalcin     }
146392bb6962SLisandro Dalcin     if (!match) { /* assume a "good" partitioner is available */
14641a83f524SJed Brown       PetscInt       na;
14651a83f524SJed Brown       const PetscInt *ia,*ja;
146692bb6962SLisandro Dalcin       ierr = MatGetRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr);
146792bb6962SLisandro Dalcin       if (done) {
146892bb6962SLisandro Dalcin         /* Build adjacency matrix by hand. Unfortunately a call to
146992bb6962SLisandro Dalcin            MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will
147092bb6962SLisandro Dalcin            remove the block-aij structure and we cannot expect
147192bb6962SLisandro Dalcin            MatPartitioning to split vertices as we need */
14720a545947SLisandro Dalcin         PetscInt       i,j,len,nnz,cnt,*iia=NULL,*jja=NULL;
14731a83f524SJed Brown         const PetscInt *row;
147492bb6962SLisandro Dalcin         nnz = 0;
147592bb6962SLisandro Dalcin         for (i=0; i<na; i++) { /* count number of nonzeros */
147692bb6962SLisandro Dalcin           len = ia[i+1] - ia[i];
147792bb6962SLisandro Dalcin           row = ja + ia[i];
147892bb6962SLisandro Dalcin           for (j=0; j<len; j++) {
147992bb6962SLisandro Dalcin             if (row[j] == i) { /* don't count diagonal */
148092bb6962SLisandro Dalcin               len--; break;
148192bb6962SLisandro Dalcin             }
148292bb6962SLisandro Dalcin           }
148392bb6962SLisandro Dalcin           nnz += len;
148492bb6962SLisandro Dalcin         }
1485854ce69bSBarry Smith         ierr   = PetscMalloc1(na+1,&iia);CHKERRQ(ierr);
1486854ce69bSBarry Smith         ierr   = PetscMalloc1(nnz,&jja);CHKERRQ(ierr);
148792bb6962SLisandro Dalcin         nnz    = 0;
148892bb6962SLisandro Dalcin         iia[0] = 0;
148992bb6962SLisandro Dalcin         for (i=0; i<na; i++) { /* fill adjacency */
149092bb6962SLisandro Dalcin           cnt = 0;
149192bb6962SLisandro Dalcin           len = ia[i+1] - ia[i];
149292bb6962SLisandro Dalcin           row = ja + ia[i];
149392bb6962SLisandro Dalcin           for (j=0; j<len; j++) {
149492bb6962SLisandro Dalcin             if (row[j] != i) { /* if not diagonal */
149592bb6962SLisandro Dalcin               jja[nnz+cnt++] = row[j];
149692bb6962SLisandro Dalcin             }
149792bb6962SLisandro Dalcin           }
149892bb6962SLisandro Dalcin           nnz     += cnt;
149992bb6962SLisandro Dalcin           iia[i+1] = nnz;
150092bb6962SLisandro Dalcin         }
150192bb6962SLisandro Dalcin         /* Partitioning of the adjacency matrix */
15020298fd71SBarry Smith         ierr      = MatCreateMPIAdj(PETSC_COMM_SELF,na,na,iia,jja,NULL,&adj);CHKERRQ(ierr);
150392bb6962SLisandro Dalcin         ierr      = MatPartitioningSetAdjacency(mpart,adj);CHKERRQ(ierr);
150492bb6962SLisandro Dalcin         ierr      = MatPartitioningSetNParts(mpart,n);CHKERRQ(ierr);
150592bb6962SLisandro Dalcin         ierr      = MatPartitioningApply(mpart,&ispart);CHKERRQ(ierr);
150692bb6962SLisandro Dalcin         ierr      = ISPartitioningToNumbering(ispart,&isnumb);CHKERRQ(ierr);
1507fcfd50ebSBarry Smith         ierr      = MatDestroy(&adj);CHKERRQ(ierr);
150892bb6962SLisandro Dalcin         foundpart = PETSC_TRUE;
150992bb6962SLisandro Dalcin       }
151092bb6962SLisandro Dalcin       ierr = MatRestoreRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr);
151192bb6962SLisandro Dalcin     }
1512fcfd50ebSBarry Smith     ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr);
151392bb6962SLisandro Dalcin   }
151492bb6962SLisandro Dalcin 
1515785e854fSJed Brown   ierr   = PetscMalloc1(n,&is);CHKERRQ(ierr);
151692bb6962SLisandro Dalcin   *outis = is;
151792bb6962SLisandro Dalcin 
151892bb6962SLisandro Dalcin   if (!foundpart) {
151992bb6962SLisandro Dalcin 
152092bb6962SLisandro Dalcin     /* Partitioning by contiguous chunks of rows */
152192bb6962SLisandro Dalcin 
152292bb6962SLisandro Dalcin     PetscInt mbs   = (rend-rstart)/bs;
152392bb6962SLisandro Dalcin     PetscInt start = rstart;
152492bb6962SLisandro Dalcin     for (i=0; i<n; i++) {
152592bb6962SLisandro Dalcin       PetscInt count = (mbs/n + ((mbs % n) > i)) * bs;
152692bb6962SLisandro Dalcin       ierr   = ISCreateStride(PETSC_COMM_SELF,count,start,1,&is[i]);CHKERRQ(ierr);
152792bb6962SLisandro Dalcin       start += count;
152892bb6962SLisandro Dalcin     }
152992bb6962SLisandro Dalcin 
153092bb6962SLisandro Dalcin   } else {
153192bb6962SLisandro Dalcin 
153292bb6962SLisandro Dalcin     /* Partitioning by adjacency of diagonal block  */
153392bb6962SLisandro Dalcin 
153492bb6962SLisandro Dalcin     const PetscInt *numbering;
153592bb6962SLisandro Dalcin     PetscInt       *count,nidx,*indices,*newidx,start=0;
153692bb6962SLisandro Dalcin     /* Get node count in each partition */
1537785e854fSJed Brown     ierr = PetscMalloc1(n,&count);CHKERRQ(ierr);
153892bb6962SLisandro Dalcin     ierr = ISPartitioningCount(ispart,n,count);CHKERRQ(ierr);
153992bb6962SLisandro Dalcin     if (isbaij && bs > 1) { /* adjust for the block-aij case */
154092bb6962SLisandro Dalcin       for (i=0; i<n; i++) count[i] *= bs;
154192bb6962SLisandro Dalcin     }
154292bb6962SLisandro Dalcin     /* Build indices from node numbering */
154392bb6962SLisandro Dalcin     ierr = ISGetLocalSize(isnumb,&nidx);CHKERRQ(ierr);
1544785e854fSJed Brown     ierr = PetscMalloc1(nidx,&indices);CHKERRQ(ierr);
154592bb6962SLisandro Dalcin     for (i=0; i<nidx; i++) indices[i] = i; /* needs to be initialized */
154692bb6962SLisandro Dalcin     ierr = ISGetIndices(isnumb,&numbering);CHKERRQ(ierr);
154792bb6962SLisandro Dalcin     ierr = PetscSortIntWithPermutation(nidx,numbering,indices);CHKERRQ(ierr);
154892bb6962SLisandro Dalcin     ierr = ISRestoreIndices(isnumb,&numbering);CHKERRQ(ierr);
154992bb6962SLisandro Dalcin     if (isbaij && bs > 1) { /* adjust for the block-aij case */
1550785e854fSJed Brown       ierr = PetscMalloc1(nidx*bs,&newidx);CHKERRQ(ierr);
15512fa5cd67SKarl Rupp       for (i=0; i<nidx; i++) {
15522fa5cd67SKarl Rupp         for (j=0; j<bs; j++) newidx[i*bs+j] = indices[i]*bs + j;
15532fa5cd67SKarl Rupp       }
155492bb6962SLisandro Dalcin       ierr    = PetscFree(indices);CHKERRQ(ierr);
155592bb6962SLisandro Dalcin       nidx   *= bs;
155692bb6962SLisandro Dalcin       indices = newidx;
155792bb6962SLisandro Dalcin     }
155892bb6962SLisandro Dalcin     /* Shift to get global indices */
155992bb6962SLisandro Dalcin     for (i=0; i<nidx; i++) indices[i] += rstart;
156092bb6962SLisandro Dalcin 
156192bb6962SLisandro Dalcin     /* Build the index sets for each block */
156292bb6962SLisandro Dalcin     for (i=0; i<n; i++) {
156370b3c8c7SBarry Smith       ierr   = ISCreateGeneral(PETSC_COMM_SELF,count[i],&indices[start],PETSC_COPY_VALUES,&is[i]);CHKERRQ(ierr);
156492bb6962SLisandro Dalcin       ierr   = ISSort(is[i]);CHKERRQ(ierr);
156592bb6962SLisandro Dalcin       start += count[i];
156692bb6962SLisandro Dalcin     }
156792bb6962SLisandro Dalcin 
15683bf036e2SBarry Smith     ierr = PetscFree(count);CHKERRQ(ierr);
15693bf036e2SBarry Smith     ierr = PetscFree(indices);CHKERRQ(ierr);
1570fcfd50ebSBarry Smith     ierr = ISDestroy(&isnumb);CHKERRQ(ierr);
1571fcfd50ebSBarry Smith     ierr = ISDestroy(&ispart);CHKERRQ(ierr);
157292bb6962SLisandro Dalcin 
157392bb6962SLisandro Dalcin   }
157492bb6962SLisandro Dalcin   PetscFunctionReturn(0);
157592bb6962SLisandro Dalcin }
157692bb6962SLisandro Dalcin 
157792bb6962SLisandro Dalcin /*@C
157892bb6962SLisandro Dalcin    PCASMDestroySubdomains - Destroys the index sets created with
157992bb6962SLisandro Dalcin    PCASMCreateSubdomains(). Should be called after setting subdomains
158092bb6962SLisandro Dalcin    with PCASMSetLocalSubdomains().
158192bb6962SLisandro Dalcin 
158292bb6962SLisandro Dalcin    Collective
158392bb6962SLisandro Dalcin 
158492bb6962SLisandro Dalcin    Input Parameters:
158592bb6962SLisandro Dalcin +  n - the number of index sets
15862b691e39SMatthew Knepley .  is - the array of index sets
15870298fd71SBarry Smith -  is_local - the array of local index sets, can be NULL
158892bb6962SLisandro Dalcin 
158992bb6962SLisandro Dalcin    Level: advanced
159092bb6962SLisandro Dalcin 
159192bb6962SLisandro Dalcin .seealso: PCASMCreateSubdomains(), PCASMSetLocalSubdomains()
159292bb6962SLisandro Dalcin @*/
15937087cfbeSBarry Smith PetscErrorCode  PCASMDestroySubdomains(PetscInt n, IS is[], IS is_local[])
159492bb6962SLisandro Dalcin {
159592bb6962SLisandro Dalcin   PetscInt       i;
159692bb6962SLisandro Dalcin   PetscErrorCode ierr;
15975fd66863SKarl Rupp 
159892bb6962SLisandro Dalcin   PetscFunctionBegin;
1599a35b7b57SDmitry Karpeev   if (n <= 0) PetscFunctionReturn(0);
1600a35b7b57SDmitry Karpeev   if (is) {
160192bb6962SLisandro Dalcin     PetscValidPointer(is,2);
1602fcfd50ebSBarry Smith     for (i=0; i<n; i++) { ierr = ISDestroy(&is[i]);CHKERRQ(ierr); }
160392bb6962SLisandro Dalcin     ierr = PetscFree(is);CHKERRQ(ierr);
1604a35b7b57SDmitry Karpeev   }
16052b691e39SMatthew Knepley   if (is_local) {
16062b691e39SMatthew Knepley     PetscValidPointer(is_local,3);
1607fcfd50ebSBarry Smith     for (i=0; i<n; i++) { ierr = ISDestroy(&is_local[i]);CHKERRQ(ierr); }
16082b691e39SMatthew Knepley     ierr = PetscFree(is_local);CHKERRQ(ierr);
16092b691e39SMatthew Knepley   }
161092bb6962SLisandro Dalcin   PetscFunctionReturn(0);
161192bb6962SLisandro Dalcin }
161292bb6962SLisandro Dalcin 
16134b9ad928SBarry Smith /*@
16144b9ad928SBarry Smith    PCASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz
16154b9ad928SBarry Smith    preconditioner for a two-dimensional problem on a regular grid.
16164b9ad928SBarry Smith 
16174b9ad928SBarry Smith    Not Collective
16184b9ad928SBarry Smith 
16194b9ad928SBarry Smith    Input Parameters:
16204b9ad928SBarry Smith +  m, n - the number of mesh points in the x and y directions
16214b9ad928SBarry Smith .  M, N - the number of subdomains in the x and y directions
16224b9ad928SBarry Smith .  dof - degrees of freedom per node
16234b9ad928SBarry Smith -  overlap - overlap in mesh lines
16244b9ad928SBarry Smith 
16254b9ad928SBarry Smith    Output Parameters:
16264b9ad928SBarry Smith +  Nsub - the number of subdomains created
16273d03bb48SJed Brown .  is - array of index sets defining overlapping (if overlap > 0) subdomains
16283d03bb48SJed Brown -  is_local - array of index sets defining non-overlapping subdomains
16294b9ad928SBarry Smith 
16304b9ad928SBarry Smith    Note:
16314b9ad928SBarry Smith    Presently PCAMSCreateSubdomains2d() is valid only for sequential
16324b9ad928SBarry Smith    preconditioners.  More general related routines are
16334b9ad928SBarry Smith    PCASMSetTotalSubdomains() and PCASMSetLocalSubdomains().
16344b9ad928SBarry Smith 
16354b9ad928SBarry Smith    Level: advanced
16364b9ad928SBarry Smith 
16374b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetLocalSubdomains(), PCASMGetSubKSP(),
16384b9ad928SBarry Smith           PCASMSetOverlap()
16394b9ad928SBarry Smith @*/
16407087cfbeSBarry Smith PetscErrorCode  PCASMCreateSubdomains2D(PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt dof,PetscInt overlap,PetscInt *Nsub,IS **is,IS **is_local)
16414b9ad928SBarry Smith {
16423d03bb48SJed Brown   PetscInt       i,j,height,width,ystart,xstart,yleft,yright,xleft,xright,loc_outer;
16436849ba73SBarry Smith   PetscErrorCode ierr;
164413f74950SBarry Smith   PetscInt       nidx,*idx,loc,ii,jj,count;
16454b9ad928SBarry Smith 
16464b9ad928SBarry Smith   PetscFunctionBegin;
1647e32f2f54SBarry Smith   if (dof != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP," ");
16484b9ad928SBarry Smith 
16494b9ad928SBarry Smith   *Nsub     = N*M;
1650854ce69bSBarry Smith   ierr      = PetscMalloc1(*Nsub,is);CHKERRQ(ierr);
1651854ce69bSBarry Smith   ierr      = PetscMalloc1(*Nsub,is_local);CHKERRQ(ierr);
16524b9ad928SBarry Smith   ystart    = 0;
16533d03bb48SJed Brown   loc_outer = 0;
16544b9ad928SBarry Smith   for (i=0; i<N; i++) {
16554b9ad928SBarry Smith     height = n/N + ((n % N) > i); /* height of subdomain */
1656e32f2f54SBarry Smith     if (height < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many N subdomains for mesh dimension n");
16574b9ad928SBarry Smith     yleft  = ystart - overlap; if (yleft < 0) yleft = 0;
16584b9ad928SBarry Smith     yright = ystart + height + overlap; if (yright > n) yright = n;
16594b9ad928SBarry Smith     xstart = 0;
16604b9ad928SBarry Smith     for (j=0; j<M; j++) {
16614b9ad928SBarry Smith       width = m/M + ((m % M) > j); /* width of subdomain */
1662e32f2f54SBarry Smith       if (width < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many M subdomains for mesh dimension m");
16634b9ad928SBarry Smith       xleft  = xstart - overlap; if (xleft < 0) xleft = 0;
16644b9ad928SBarry Smith       xright = xstart + width + overlap; if (xright > m) xright = m;
16654b9ad928SBarry Smith       nidx   = (xright - xleft)*(yright - yleft);
1666785e854fSJed Brown       ierr   = PetscMalloc1(nidx,&idx);CHKERRQ(ierr);
16674b9ad928SBarry Smith       loc    = 0;
16684b9ad928SBarry Smith       for (ii=yleft; ii<yright; ii++) {
16694b9ad928SBarry Smith         count = m*ii + xleft;
16702fa5cd67SKarl Rupp         for (jj=xleft; jj<xright; jj++) idx[loc++] = count++;
16714b9ad928SBarry Smith       }
167270b3c8c7SBarry Smith       ierr = ISCreateGeneral(PETSC_COMM_SELF,nidx,idx,PETSC_COPY_VALUES,(*is)+loc_outer);CHKERRQ(ierr);
16733d03bb48SJed Brown       if (overlap == 0) {
16743d03bb48SJed Brown         ierr = PetscObjectReference((PetscObject)(*is)[loc_outer]);CHKERRQ(ierr);
16752fa5cd67SKarl Rupp 
16763d03bb48SJed Brown         (*is_local)[loc_outer] = (*is)[loc_outer];
16773d03bb48SJed Brown       } else {
16783d03bb48SJed Brown         for (loc=0,ii=ystart; ii<ystart+height; ii++) {
16793d03bb48SJed Brown           for (jj=xstart; jj<xstart+width; jj++) {
16803d03bb48SJed Brown             idx[loc++] = m*ii + jj;
16813d03bb48SJed Brown           }
16823d03bb48SJed Brown         }
168370b3c8c7SBarry Smith         ierr = ISCreateGeneral(PETSC_COMM_SELF,loc,idx,PETSC_COPY_VALUES,*is_local+loc_outer);CHKERRQ(ierr);
16843d03bb48SJed Brown       }
16854b9ad928SBarry Smith       ierr    = PetscFree(idx);CHKERRQ(ierr);
16864b9ad928SBarry Smith       xstart += width;
16873d03bb48SJed Brown       loc_outer++;
16884b9ad928SBarry Smith     }
16894b9ad928SBarry Smith     ystart += height;
16904b9ad928SBarry Smith   }
16914b9ad928SBarry Smith   for (i=0; i<*Nsub; i++) { ierr = ISSort((*is)[i]);CHKERRQ(ierr); }
16924b9ad928SBarry Smith   PetscFunctionReturn(0);
16934b9ad928SBarry Smith }
16944b9ad928SBarry Smith 
16954b9ad928SBarry Smith /*@C
16964b9ad928SBarry Smith     PCASMGetLocalSubdomains - Gets the local subdomains (for this processor
16974b9ad928SBarry Smith     only) for the additive Schwarz preconditioner.
16984b9ad928SBarry Smith 
1699ad4df100SBarry Smith     Not Collective
17004b9ad928SBarry Smith 
17014b9ad928SBarry Smith     Input Parameter:
17024b9ad928SBarry Smith .   pc - the preconditioner context
17034b9ad928SBarry Smith 
17044b9ad928SBarry Smith     Output Parameters:
17054b9ad928SBarry Smith +   n - the number of subdomains for this processor (default value = 1)
17062b691e39SMatthew Knepley .   is - the index sets that define the subdomains for this processor
17070298fd71SBarry Smith -   is_local - the index sets that define the local part of the subdomains for this processor (can be NULL)
17084b9ad928SBarry Smith 
17094b9ad928SBarry Smith 
17104b9ad928SBarry Smith     Notes:
17114b9ad928SBarry Smith     The IS numbering is in the parallel, global numbering of the vector.
17124b9ad928SBarry Smith 
17134b9ad928SBarry Smith     Level: advanced
17144b9ad928SBarry Smith 
17154b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(),
17164b9ad928SBarry Smith           PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubmatrices()
17174b9ad928SBarry Smith @*/
17187087cfbeSBarry Smith PetscErrorCode  PCASMGetLocalSubdomains(PC pc,PetscInt *n,IS *is[],IS *is_local[])
17194b9ad928SBarry Smith {
17202a808120SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
172192bb6962SLisandro Dalcin   PetscErrorCode ierr;
1722ace3abfcSBarry Smith   PetscBool      match;
17234b9ad928SBarry Smith 
17244b9ad928SBarry Smith   PetscFunctionBegin;
17250700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
17264482741eSBarry Smith   PetscValidIntPointer(n,2);
172792bb6962SLisandro Dalcin   if (is) PetscValidPointer(is,3);
1728251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr);
17292a808120SBarry Smith   if (!match) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"PC is not a PCASM");
17304b9ad928SBarry Smith   if (n) *n = osm->n_local_true;
17314b9ad928SBarry Smith   if (is) *is = osm->is;
17322b691e39SMatthew Knepley   if (is_local) *is_local = osm->is_local;
17334b9ad928SBarry Smith   PetscFunctionReturn(0);
17344b9ad928SBarry Smith }
17354b9ad928SBarry Smith 
17364b9ad928SBarry Smith /*@C
17374b9ad928SBarry Smith     PCASMGetLocalSubmatrices - Gets the local submatrices (for this processor
17384b9ad928SBarry Smith     only) for the additive Schwarz preconditioner.
17394b9ad928SBarry Smith 
1740ad4df100SBarry Smith     Not Collective
17414b9ad928SBarry Smith 
17424b9ad928SBarry Smith     Input Parameter:
17434b9ad928SBarry Smith .   pc - the preconditioner context
17444b9ad928SBarry Smith 
17454b9ad928SBarry Smith     Output Parameters:
17464b9ad928SBarry Smith +   n - the number of matrices for this processor (default value = 1)
17474b9ad928SBarry Smith -   mat - the matrices
17484b9ad928SBarry Smith 
17494b9ad928SBarry Smith     Level: advanced
17504b9ad928SBarry Smith 
175195452b02SPatrick Sanan     Notes:
175234a84908Sprj-     Call after PCSetUp() (or KSPSetUp()) but before PCApply() and before PCSetUpOnBlocks())
1753cf739d55SBarry Smith 
175434a84908Sprj-            Usually one would use PCSetModifySubMatrices() to change the submatrices in building the preconditioner.
1755cf739d55SBarry Smith 
17564b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(),
175734a84908Sprj-           PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains(), PCSetModifySubMatrices()
17584b9ad928SBarry Smith @*/
17597087cfbeSBarry Smith PetscErrorCode  PCASMGetLocalSubmatrices(PC pc,PetscInt *n,Mat *mat[])
17604b9ad928SBarry Smith {
17614b9ad928SBarry Smith   PC_ASM         *osm;
176292bb6962SLisandro Dalcin   PetscErrorCode ierr;
1763ace3abfcSBarry Smith   PetscBool      match;
17644b9ad928SBarry Smith 
17654b9ad928SBarry Smith   PetscFunctionBegin;
17660700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
176792bb6962SLisandro Dalcin   PetscValidIntPointer(n,2);
176892bb6962SLisandro Dalcin   if (mat) PetscValidPointer(mat,3);
176934a84908Sprj-   if (!pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call after KSPSetUp() or PCSetUp().");
1770251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr);
177192bb6962SLisandro Dalcin   if (!match) {
177292bb6962SLisandro Dalcin     if (n) *n = 0;
17730298fd71SBarry Smith     if (mat) *mat = NULL;
177492bb6962SLisandro Dalcin   } else {
17754b9ad928SBarry Smith     osm = (PC_ASM*)pc->data;
17764b9ad928SBarry Smith     if (n) *n = osm->n_local_true;
17774b9ad928SBarry Smith     if (mat) *mat = osm->pmat;
177892bb6962SLisandro Dalcin   }
17794b9ad928SBarry Smith   PetscFunctionReturn(0);
17804b9ad928SBarry Smith }
1781d709ea83SDmitry Karpeev 
1782d709ea83SDmitry Karpeev /*@
1783d709ea83SDmitry Karpeev     PCASMSetDMSubdomains - Indicates whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible.
1784f1ee410cSBarry Smith 
1785d709ea83SDmitry Karpeev     Logically Collective
1786d709ea83SDmitry Karpeev 
1787d709ea83SDmitry Karpeev     Input Parameter:
1788d709ea83SDmitry Karpeev +   pc  - the preconditioner
1789d709ea83SDmitry Karpeev -   flg - boolean indicating whether to use subdomains defined by the DM
1790d709ea83SDmitry Karpeev 
1791d709ea83SDmitry Karpeev     Options Database Key:
1792d709ea83SDmitry Karpeev .   -pc_asm_dm_subdomains
1793d709ea83SDmitry Karpeev 
1794d709ea83SDmitry Karpeev     Level: intermediate
1795d709ea83SDmitry Karpeev 
1796d709ea83SDmitry Karpeev     Notes:
1797d709ea83SDmitry Karpeev     PCASMSetTotalSubdomains() and PCASMSetOverlap() take precedence over PCASMSetDMSubdomains(),
1798d709ea83SDmitry Karpeev     so setting either of the first two effectively turns the latter off.
1799d709ea83SDmitry Karpeev 
1800d709ea83SDmitry Karpeev .seealso: PCASMGetDMSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap()
1801d709ea83SDmitry Karpeev           PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains()
1802d709ea83SDmitry Karpeev @*/
1803d709ea83SDmitry Karpeev PetscErrorCode  PCASMSetDMSubdomains(PC pc,PetscBool flg)
1804d709ea83SDmitry Karpeev {
1805d709ea83SDmitry Karpeev   PC_ASM         *osm = (PC_ASM*)pc->data;
1806d709ea83SDmitry Karpeev   PetscErrorCode ierr;
1807d709ea83SDmitry Karpeev   PetscBool      match;
1808d709ea83SDmitry Karpeev 
1809d709ea83SDmitry Karpeev   PetscFunctionBegin;
1810d709ea83SDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1811d709ea83SDmitry Karpeev   PetscValidLogicalCollectiveBool(pc,flg,2);
1812d709ea83SDmitry Karpeev   if (pc->setupcalled) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Not for a setup PC.");
1813d709ea83SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr);
1814d709ea83SDmitry Karpeev   if (match) {
1815d709ea83SDmitry Karpeev     osm->dm_subdomains = flg;
1816d709ea83SDmitry Karpeev   }
1817d709ea83SDmitry Karpeev   PetscFunctionReturn(0);
1818d709ea83SDmitry Karpeev }
1819d709ea83SDmitry Karpeev 
1820d709ea83SDmitry Karpeev /*@
1821d709ea83SDmitry Karpeev     PCASMGetDMSubdomains - Returns flag indicating whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible.
1822d709ea83SDmitry Karpeev     Not Collective
1823d709ea83SDmitry Karpeev 
1824d709ea83SDmitry Karpeev     Input Parameter:
1825d709ea83SDmitry Karpeev .   pc  - the preconditioner
1826d709ea83SDmitry Karpeev 
1827d709ea83SDmitry Karpeev     Output Parameter:
1828d709ea83SDmitry Karpeev .   flg - boolean indicating whether to use subdomains defined by the DM
1829d709ea83SDmitry Karpeev 
1830d709ea83SDmitry Karpeev     Level: intermediate
1831d709ea83SDmitry Karpeev 
1832d709ea83SDmitry Karpeev .seealso: PCASMSetDMSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap()
1833d709ea83SDmitry Karpeev           PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains()
1834d709ea83SDmitry Karpeev @*/
1835d709ea83SDmitry Karpeev PetscErrorCode  PCASMGetDMSubdomains(PC pc,PetscBool* flg)
1836d709ea83SDmitry Karpeev {
1837d709ea83SDmitry Karpeev   PC_ASM         *osm = (PC_ASM*)pc->data;
1838d709ea83SDmitry Karpeev   PetscErrorCode ierr;
1839d709ea83SDmitry Karpeev   PetscBool      match;
1840d709ea83SDmitry Karpeev 
1841d709ea83SDmitry Karpeev   PetscFunctionBegin;
1842d709ea83SDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1843534a8f05SLisandro Dalcin   PetscValidBoolPointer(flg,2);
1844d709ea83SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr);
184556ea09ceSStefano Zampini   if (match) *flg = osm->dm_subdomains;
184656ea09ceSStefano Zampini   else *flg = PETSC_FALSE;
1847d709ea83SDmitry Karpeev   PetscFunctionReturn(0);
1848d709ea83SDmitry Karpeev }
184980ec0b7dSPatrick Sanan 
185080ec0b7dSPatrick Sanan /*@
185180ec0b7dSPatrick Sanan      PCASMGetSubMatType - Gets the matrix type used for ASM subsolves, as a string.
185280ec0b7dSPatrick Sanan 
185380ec0b7dSPatrick Sanan    Not Collective
185480ec0b7dSPatrick Sanan 
185580ec0b7dSPatrick Sanan    Input Parameter:
185680ec0b7dSPatrick Sanan .  pc - the PC
185780ec0b7dSPatrick Sanan 
185880ec0b7dSPatrick Sanan    Output Parameter:
1859f1ee410cSBarry Smith .  -pc_asm_sub_mat_type - name of matrix type
186080ec0b7dSPatrick Sanan 
186180ec0b7dSPatrick Sanan    Level: advanced
186280ec0b7dSPatrick Sanan 
186380ec0b7dSPatrick Sanan .seealso: PCASMSetSubMatType(), PCASM, PCSetType(), VecSetType(), MatType, Mat
186480ec0b7dSPatrick Sanan @*/
186556ea09ceSStefano Zampini PetscErrorCode  PCASMGetSubMatType(PC pc,MatType *sub_mat_type)
186656ea09ceSStefano Zampini {
186780ec0b7dSPatrick Sanan   PetscErrorCode ierr;
186880ec0b7dSPatrick Sanan 
186956ea09ceSStefano Zampini   PetscFunctionBegin;
187056ea09ceSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
187180ec0b7dSPatrick Sanan   ierr = PetscTryMethod(pc,"PCASMGetSubMatType_C",(PC,MatType*),(pc,sub_mat_type));CHKERRQ(ierr);
187280ec0b7dSPatrick Sanan   PetscFunctionReturn(0);
187380ec0b7dSPatrick Sanan }
187480ec0b7dSPatrick Sanan 
187580ec0b7dSPatrick Sanan /*@
187680ec0b7dSPatrick Sanan      PCASMSetSubMatType - Set the type of matrix used for ASM subsolves
187780ec0b7dSPatrick Sanan 
187880ec0b7dSPatrick Sanan    Collective on Mat
187980ec0b7dSPatrick Sanan 
188080ec0b7dSPatrick Sanan    Input Parameters:
188180ec0b7dSPatrick Sanan +  pc             - the PC object
188280ec0b7dSPatrick Sanan -  sub_mat_type   - matrix type
188380ec0b7dSPatrick Sanan 
188480ec0b7dSPatrick Sanan    Options Database Key:
188580ec0b7dSPatrick 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.
188680ec0b7dSPatrick Sanan 
188780ec0b7dSPatrick Sanan    Notes:
188880ec0b7dSPatrick Sanan    See "${PETSC_DIR}/include/petscmat.h" for available types
188980ec0b7dSPatrick Sanan 
189080ec0b7dSPatrick Sanan   Level: advanced
189180ec0b7dSPatrick Sanan 
189280ec0b7dSPatrick Sanan .seealso: PCASMGetSubMatType(), PCASM, PCSetType(), VecSetType(), MatType, Mat
189380ec0b7dSPatrick Sanan @*/
189480ec0b7dSPatrick Sanan PetscErrorCode PCASMSetSubMatType(PC pc,MatType sub_mat_type)
189580ec0b7dSPatrick Sanan {
189680ec0b7dSPatrick Sanan   PetscErrorCode ierr;
189780ec0b7dSPatrick Sanan 
189856ea09ceSStefano Zampini   PetscFunctionBegin;
189956ea09ceSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
190080ec0b7dSPatrick Sanan   ierr = PetscTryMethod(pc,"PCASMSetSubMatType_C",(PC,MatType),(pc,sub_mat_type));CHKERRQ(ierr);
190180ec0b7dSPatrick Sanan   PetscFunctionReturn(0);
190280ec0b7dSPatrick Sanan }
1903