xref: /petsc/src/ksp/pc/impls/asm/asm.c (revision 28b400f66ebc7ae0049166a2294dfcd3df27e64b)
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 
1334eca32bSPierre Jolivet #include <petsc/private/pcasmimpl.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;
1813f74950SBarry Smith   PetscMPIInt       rank;
194d219a6aSLisandro Dalcin   PetscInt          i,bsz;
20ace3abfcSBarry Smith   PetscBool         iascii,isstring;
214b9ad928SBarry Smith   PetscViewer       sviewer;
22ed155784SPierre Jolivet   PetscViewerFormat format;
23ed155784SPierre Jolivet   const char        *prefix;
244b9ad928SBarry Smith 
254b9ad928SBarry Smith   PetscFunctionBegin;
265f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
275f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring));
2832077d6dSBarry Smith   if (iascii) {
293d03bb48SJed Brown     char overlaps[256] = "user-defined overlap",blocks[256] = "total subdomain blocks not yet set";
305f80ce2aSJacob Faibussowitsch     if (osm->overlap >= 0) CHKERRQ(PetscSNPrintf(overlaps,sizeof(overlaps),"amount of overlap = %D",osm->overlap));
315f80ce2aSJacob Faibussowitsch     if (osm->n > 0) CHKERRQ(PetscSNPrintf(blocks,sizeof(blocks),"total subdomain blocks = %D",osm->n));
325f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPrintf(viewer,"  %s, %s\n",blocks,overlaps));
335f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPrintf(viewer,"  restriction/interpolation type - %s\n",PCASMTypes[osm->type]));
345f80ce2aSJacob Faibussowitsch     if (osm->dm_subdomains) CHKERRQ(PetscViewerASCIIPrintf(viewer,"  Additive Schwarz: using DM to define subdomains\n"));
355f80ce2aSJacob Faibussowitsch     if (osm->loctype != PC_COMPOSITE_ADDITIVE) CHKERRQ(PetscViewerASCIIPrintf(viewer,"  Additive Schwarz: local solve composition type - %s\n",PCCompositeTypes[osm->loctype]));
365f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank));
375f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerGetFormat(viewer,&format));
38ed155784SPierre Jolivet     if (format != PETSC_VIEWER_ASCII_INFO_DETAIL) {
394d219a6aSLisandro Dalcin       if (osm->ksp) {
405f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscViewerASCIIPrintf(viewer,"  Local solver information for first block is in the following KSP and PC objects on rank 0:\n"));
415f80ce2aSJacob Faibussowitsch         CHKERRQ(PCGetOptionsPrefix(pc,&prefix));
425f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscViewerASCIIPrintf(viewer,"  Use -%sksp_view ::ascii_info_detail to display information for all blocks\n",prefix?prefix:""));
435f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer));
44dd400576SPatrick Sanan         if (rank == 0) {
455f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscViewerASCIIPushTab(viewer));
465f80ce2aSJacob Faibussowitsch           CHKERRQ(KSPView(osm->ksp[0],sviewer));
475f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscViewerASCIIPopTab(viewer));
484b9ad928SBarry Smith         }
495f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer));
504d219a6aSLisandro Dalcin       }
514b9ad928SBarry Smith     } else {
525f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerASCIIPushSynchronized(viewer));
535f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerASCIISynchronizedPrintf(viewer,"  [%d] number of local blocks = %D\n",(int)rank,osm->n_local_true));
545f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerFlush(viewer));
555f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerASCIIPrintf(viewer,"  Local solver information for each block is in the following KSP and PC objects:\n"));
565f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerASCIIPushTab(viewer));
575f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerASCIIPrintf(viewer,"- - - - - - - - - - - - - - - - - -\n"));
585f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer));
59dd2fa690SBarry Smith       for (i=0; i<osm->n_local_true; i++) {
605f80ce2aSJacob Faibussowitsch         CHKERRQ(ISGetLocalSize(osm->is[i],&bsz));
615f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscViewerASCIISynchronizedPrintf(sviewer,"[%d] local block number %D, size = %D\n",(int)rank,i,bsz));
625f80ce2aSJacob Faibussowitsch         CHKERRQ(KSPView(osm->ksp[i],sviewer));
635f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscViewerASCIISynchronizedPrintf(sviewer,"- - - - - - - - - - - - - - - - - -\n"));
644b9ad928SBarry Smith       }
655f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer));
665f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerASCIIPopTab(viewer));
675f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerFlush(viewer));
685f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerASCIIPopSynchronized(viewer));
694b9ad928SBarry Smith     }
704b9ad928SBarry Smith   } else if (isstring) {
715f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerStringSPrintf(viewer," blocks=%D, overlap=%D, type=%s",osm->n,osm->overlap,PCASMTypes[osm->type]));
725f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer));
735f80ce2aSJacob Faibussowitsch     if (osm->ksp) CHKERRQ(KSPView(osm->ksp[0],sviewer));
745f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer));
754b9ad928SBarry Smith   }
764b9ad928SBarry Smith   PetscFunctionReturn(0);
774b9ad928SBarry Smith }
784b9ad928SBarry Smith 
7992bb6962SLisandro Dalcin static PetscErrorCode PCASMPrintSubdomains(PC pc)
8092bb6962SLisandro Dalcin {
8192bb6962SLisandro Dalcin   PC_ASM         *osm = (PC_ASM*)pc->data;
8292bb6962SLisandro Dalcin   const char     *prefix;
8392bb6962SLisandro Dalcin   char           fname[PETSC_MAX_PATH_LEN+1];
84643535aeSDmitry Karpeev   PetscViewer    viewer, sviewer;
85643535aeSDmitry Karpeev   char           *s;
8692bb6962SLisandro Dalcin   PetscInt       i,j,nidx;
8792bb6962SLisandro Dalcin   const PetscInt *idx;
88643535aeSDmitry Karpeev   PetscMPIInt    rank, size;
89846783a0SBarry Smith 
9092bb6962SLisandro Dalcin   PetscFunctionBegin;
915f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size));
925f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank));
935f80ce2aSJacob Faibussowitsch   CHKERRQ(PCGetOptionsPrefix(pc,&prefix));
945f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetString(NULL,prefix,"-pc_asm_print_subdomains",fname,sizeof(fname),NULL));
955f80ce2aSJacob Faibussowitsch   if (fname[0] == 0) CHKERRQ(PetscStrcpy(fname,"stdout"));
965f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIOpen(PetscObjectComm((PetscObject)pc),fname,&viewer));
97643535aeSDmitry Karpeev   for (i=0; i<osm->n_local; i++) {
98643535aeSDmitry Karpeev     if (i < osm->n_local_true) {
995f80ce2aSJacob Faibussowitsch       CHKERRQ(ISGetLocalSize(osm->is[i],&nidx));
1005f80ce2aSJacob Faibussowitsch       CHKERRQ(ISGetIndices(osm->is[i],&idx));
101643535aeSDmitry Karpeev       /* Print to a string viewer; no more than 15 characters per index plus 512 char for the header.*/
10236a9e3b9SBarry Smith #define len  16*(nidx+1)+512
1035f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscMalloc1(len,&s));
1045f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerStringOpen(PETSC_COMM_SELF, s, len, &sviewer));
10536a9e3b9SBarry Smith #undef len
1065f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerStringSPrintf(sviewer, "[%D:%D] Subdomain %D with overlap:\n", rank, size, i));
10792bb6962SLisandro Dalcin       for (j=0; j<nidx; j++) {
1085f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscViewerStringSPrintf(sviewer,"%D ",idx[j]));
10992bb6962SLisandro Dalcin       }
1105f80ce2aSJacob Faibussowitsch       CHKERRQ(ISRestoreIndices(osm->is[i],&idx));
1115f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerStringSPrintf(sviewer,"\n"));
1125f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerDestroy(&sviewer));
1135f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerASCIIPushSynchronized(viewer));
1145f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerASCIISynchronizedPrintf(viewer, s));
1155f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerFlush(viewer));
1165f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerASCIIPopSynchronized(viewer));
1175f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFree(s));
1182b691e39SMatthew Knepley       if (osm->is_local) {
119643535aeSDmitry Karpeev         /* Print to a string viewer; no more than 15 characters per index plus 512 char for the header.*/
12036a9e3b9SBarry Smith #define len  16*(nidx+1)+512
1215f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscMalloc1(len, &s));
1225f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscViewerStringOpen(PETSC_COMM_SELF, s, len, &sviewer));
12336a9e3b9SBarry Smith #undef len
1245f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscViewerStringSPrintf(sviewer, "[%D:%D] Subdomain %D without overlap:\n", rank, size, i));
1255f80ce2aSJacob Faibussowitsch         CHKERRQ(ISGetLocalSize(osm->is_local[i],&nidx));
1265f80ce2aSJacob Faibussowitsch         CHKERRQ(ISGetIndices(osm->is_local[i],&idx));
1272b691e39SMatthew Knepley         for (j=0; j<nidx; j++) {
1285f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscViewerStringSPrintf(sviewer,"%D ",idx[j]));
1292b691e39SMatthew Knepley         }
1305f80ce2aSJacob Faibussowitsch         CHKERRQ(ISRestoreIndices(osm->is_local[i],&idx));
1315f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscViewerStringSPrintf(sviewer,"\n"));
1325f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscViewerDestroy(&sviewer));
1335f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscViewerASCIIPushSynchronized(viewer));
1345f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscViewerASCIISynchronizedPrintf(viewer, s));
1355f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscViewerFlush(viewer));
1365f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscViewerASCIIPopSynchronized(viewer));
1375f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscFree(s));
138643535aeSDmitry Karpeev       }
1392fa5cd67SKarl Rupp     } else {
140643535aeSDmitry Karpeev       /* Participate in collective viewer calls. */
1415f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerASCIIPushSynchronized(viewer));
1425f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerFlush(viewer));
1435f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerASCIIPopSynchronized(viewer));
144643535aeSDmitry Karpeev       /* Assume either all ranks have is_local or none do. */
145643535aeSDmitry Karpeev       if (osm->is_local) {
1465f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscViewerASCIIPushSynchronized(viewer));
1475f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscViewerFlush(viewer));
1485f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscViewerASCIIPopSynchronized(viewer));
149643535aeSDmitry Karpeev       }
150fdc48646SDmitry Karpeev     }
15192bb6962SLisandro Dalcin   }
1525f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerFlush(viewer));
1535f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&viewer));
15492bb6962SLisandro Dalcin   PetscFunctionReturn(0);
15592bb6962SLisandro Dalcin }
15692bb6962SLisandro Dalcin 
1576849ba73SBarry Smith static PetscErrorCode PCSetUp_ASM(PC pc)
1584b9ad928SBarry Smith {
1594b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
16057501b6eSBarry Smith   PetscBool      flg;
16187e86712SBarry Smith   PetscInt       i,m,m_local;
1624b9ad928SBarry Smith   MatReuse       scall = MAT_REUSE_MATRIX;
1634b9ad928SBarry Smith   IS             isl;
1644b9ad928SBarry Smith   KSP            ksp;
1654b9ad928SBarry Smith   PC             subpc;
1662dcb1b2aSMatthew Knepley   const char     *prefix,*pprefix;
16723ce1328SBarry Smith   Vec            vec;
1680298fd71SBarry Smith   DM             *domain_dm = NULL;
1694b9ad928SBarry Smith 
1704b9ad928SBarry Smith   PetscFunctionBegin;
1714b9ad928SBarry Smith   if (!pc->setupcalled) {
172265a2120SBarry Smith     PetscInt m;
17392bb6962SLisandro Dalcin 
1742b837212SDmitry Karpeev     /* Note: if subdomains have been set either via PCASMSetTotalSubdomains() or via PCASMSetLocalSubdomains(), osm->n_local_true will not be PETSC_DECIDE */
1752b837212SDmitry Karpeev     if (osm->n_local_true == PETSC_DECIDE) {
17669ca1f37SDmitry Karpeev       /* no subdomains given */
17765db9045SDmitry Karpeev       /* try pc->dm first, if allowed */
178d709ea83SDmitry Karpeev       if (osm->dm_subdomains && pc->dm) {
179feb221f8SDmitry Karpeev         PetscInt  num_domains, d;
180feb221f8SDmitry Karpeev         char      **domain_names;
1818d4ac253SDmitry Karpeev         IS        *inner_domain_is, *outer_domain_is;
1825f80ce2aSJacob Faibussowitsch         CHKERRQ(DMCreateDomainDecomposition(pc->dm, &num_domains, &domain_names, &inner_domain_is, &outer_domain_is, &domain_dm));
183704f0395SPatrick Sanan         osm->overlap = -1; /* We do not want to increase the overlap of the IS.
184704f0395SPatrick Sanan                               A future improvement of this code might allow one to use
185704f0395SPatrick Sanan                               DM-defined subdomains and also increase the overlap,
186704f0395SPatrick Sanan                               but that is not currently supported */
18769ca1f37SDmitry Karpeev         if (num_domains) {
1885f80ce2aSJacob Faibussowitsch           CHKERRQ(PCASMSetLocalSubdomains(pc, num_domains, outer_domain_is, inner_domain_is));
18969ca1f37SDmitry Karpeev         }
190feb221f8SDmitry Karpeev         for (d = 0; d < num_domains; ++d) {
1915f80ce2aSJacob Faibussowitsch           if (domain_names)    CHKERRQ(PetscFree(domain_names[d]));
1925f80ce2aSJacob Faibussowitsch           if (inner_domain_is) CHKERRQ(ISDestroy(&inner_domain_is[d]));
1935f80ce2aSJacob Faibussowitsch           if (outer_domain_is) CHKERRQ(ISDestroy(&outer_domain_is[d]));
194feb221f8SDmitry Karpeev         }
1955f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscFree(domain_names));
1965f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscFree(inner_domain_is));
1975f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscFree(outer_domain_is));
198feb221f8SDmitry Karpeev       }
1992b837212SDmitry Karpeev       if (osm->n_local_true == PETSC_DECIDE) {
20069ca1f37SDmitry Karpeev         /* still no subdomains; use one subdomain per processor */
2012b837212SDmitry Karpeev         osm->n_local_true = 1;
202feb221f8SDmitry Karpeev       }
2032b837212SDmitry Karpeev     }
2042b837212SDmitry Karpeev     { /* determine the global and max number of subdomains */
2056ac3741eSJed Brown       struct {PetscInt max,sum;} inwork,outwork;
206c10200c1SHong Zhang       PetscMPIInt size;
207c10200c1SHong Zhang 
2086ac3741eSJed Brown       inwork.max   = osm->n_local_true;
2096ac3741eSJed Brown       inwork.sum   = osm->n_local_true;
2105f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPIU_Allreduce(&inwork,&outwork,1,MPIU_2INT,MPIU_MAXSUM_OP,PetscObjectComm((PetscObject)pc)));
2116ac3741eSJed Brown       osm->n_local = outwork.max;
2126ac3741eSJed Brown       osm->n       = outwork.sum;
213c10200c1SHong Zhang 
2145f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size));
215c10200c1SHong Zhang       if (outwork.max == 1 && outwork.sum == size) {
2167dae84e0SHong Zhang         /* osm->n_local_true = 1 on all processes, set this option may enable use of optimized MatCreateSubMatrices() implementation */
2175f80ce2aSJacob Faibussowitsch         CHKERRQ(MatSetOption(pc->pmat,MAT_SUBMAT_SINGLEIS,PETSC_TRUE));
218c10200c1SHong Zhang       }
2194b9ad928SBarry Smith     }
22078904715SLisandro Dalcin     if (!osm->is) { /* create the index sets */
2215f80ce2aSJacob Faibussowitsch       CHKERRQ(PCASMCreateSubdomains(pc->pmat,osm->n_local_true,&osm->is));
2224b9ad928SBarry Smith     }
223f5234e35SJed Brown     if (osm->n_local_true > 1 && !osm->is_local) {
2245f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscMalloc1(osm->n_local_true,&osm->is_local));
225f5234e35SJed Brown       for (i=0; i<osm->n_local_true; i++) {
226f5234e35SJed Brown         if (osm->overlap > 0) { /* With positive overlap, osm->is[i] will be modified */
2275f80ce2aSJacob Faibussowitsch           CHKERRQ(ISDuplicate(osm->is[i],&osm->is_local[i]));
2285f80ce2aSJacob Faibussowitsch           CHKERRQ(ISCopy(osm->is[i],osm->is_local[i]));
229f5234e35SJed Brown         } else {
2305f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscObjectReference((PetscObject)osm->is[i]));
231f5234e35SJed Brown           osm->is_local[i] = osm->is[i];
232f5234e35SJed Brown         }
233f5234e35SJed Brown       }
234f5234e35SJed Brown     }
2355f80ce2aSJacob Faibussowitsch     CHKERRQ(PCGetOptionsPrefix(pc,&prefix));
23690d69ab7SBarry Smith     flg  = PETSC_FALSE;
2375f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsHasName(NULL,prefix,"-pc_asm_print_subdomains",&flg));
2385f80ce2aSJacob Faibussowitsch     if (flg) CHKERRQ(PCASMPrintSubdomains(pc));
2394b9ad928SBarry Smith 
2403d03bb48SJed Brown     if (osm->overlap > 0) {
2414b9ad928SBarry Smith       /* Extend the "overlapping" regions by a number of steps */
2425f80ce2aSJacob Faibussowitsch       CHKERRQ(MatIncreaseOverlap(pc->pmat,osm->n_local_true,osm->is,osm->overlap));
2433d03bb48SJed Brown     }
2446ed231c7SMatthew Knepley     if (osm->sort_indices) {
24592bb6962SLisandro Dalcin       for (i=0; i<osm->n_local_true; i++) {
2465f80ce2aSJacob Faibussowitsch         CHKERRQ(ISSort(osm->is[i]));
2472b691e39SMatthew Knepley         if (osm->is_local) {
2485f80ce2aSJacob Faibussowitsch           CHKERRQ(ISSort(osm->is_local[i]));
2492b691e39SMatthew Knepley         }
2504b9ad928SBarry Smith       }
2516ed231c7SMatthew Knepley     }
2524b9ad928SBarry Smith 
253f6991133SBarry Smith     if (!osm->ksp) {
25478904715SLisandro Dalcin       /* Create the local solvers */
2555f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscMalloc1(osm->n_local_true,&osm->ksp));
256feb221f8SDmitry Karpeev       if (domain_dm) {
2575f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscInfo(pc,"Setting up ASM subproblems using the embedded DM\n"));
258feb221f8SDmitry Karpeev       }
25992bb6962SLisandro Dalcin       for (i=0; i<osm->n_local_true; i++) {
2605f80ce2aSJacob Faibussowitsch         CHKERRQ(KSPCreate(PETSC_COMM_SELF,&ksp));
2615f80ce2aSJacob Faibussowitsch         CHKERRQ(KSPSetErrorIfNotConverged(ksp,pc->erroriffailure));
2625f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscLogObjectParent((PetscObject)pc,(PetscObject)ksp));
2635f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1));
2645f80ce2aSJacob Faibussowitsch         CHKERRQ(KSPSetType(ksp,KSPPREONLY));
2655f80ce2aSJacob Faibussowitsch         CHKERRQ(KSPGetPC(ksp,&subpc));
2665f80ce2aSJacob Faibussowitsch         CHKERRQ(PCGetOptionsPrefix(pc,&prefix));
2675f80ce2aSJacob Faibussowitsch         CHKERRQ(KSPSetOptionsPrefix(ksp,prefix));
2685f80ce2aSJacob Faibussowitsch         CHKERRQ(KSPAppendOptionsPrefix(ksp,"sub_"));
269feb221f8SDmitry Karpeev         if (domain_dm) {
2705f80ce2aSJacob Faibussowitsch           CHKERRQ(KSPSetDM(ksp, domain_dm[i]));
2715f80ce2aSJacob Faibussowitsch           CHKERRQ(KSPSetDMActive(ksp, PETSC_FALSE));
2725f80ce2aSJacob Faibussowitsch           CHKERRQ(DMDestroy(&domain_dm[i]));
273feb221f8SDmitry Karpeev         }
2744b9ad928SBarry Smith         osm->ksp[i] = ksp;
2754b9ad928SBarry Smith       }
276feb221f8SDmitry Karpeev       if (domain_dm) {
2775f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscFree(domain_dm));
278feb221f8SDmitry Karpeev       }
279f6991133SBarry Smith     }
2801dd8081eSeaulisa 
2815f80ce2aSJacob Faibussowitsch     CHKERRQ(ISConcatenate(PETSC_COMM_SELF, osm->n_local_true, osm->is, &osm->lis));
2825f80ce2aSJacob Faibussowitsch     CHKERRQ(ISSortRemoveDups(osm->lis));
2835f80ce2aSJacob Faibussowitsch     CHKERRQ(ISGetLocalSize(osm->lis, &m));
2841dd8081eSeaulisa 
2854b9ad928SBarry Smith     scall = MAT_INITIAL_MATRIX;
2864b9ad928SBarry Smith   } else {
2874b9ad928SBarry Smith     /*
2884b9ad928SBarry Smith        Destroy the blocks from the previous iteration
2894b9ad928SBarry Smith     */
290e09e08ccSBarry Smith     if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
2915f80ce2aSJacob Faibussowitsch       CHKERRQ(MatDestroyMatrices(osm->n_local_true,&osm->pmat));
2924b9ad928SBarry Smith       scall = MAT_INITIAL_MATRIX;
2934b9ad928SBarry Smith     }
2944b9ad928SBarry Smith   }
2954b9ad928SBarry Smith 
296b58cb649SBarry Smith   /* Destroy previous submatrices of a different type than pc->pmat since MAT_REUSE_MATRIX won't work in that case */
297b58cb649SBarry Smith   if ((scall == MAT_REUSE_MATRIX) && osm->sub_mat_type) {
298b58cb649SBarry Smith     if (osm->n_local_true > 0) {
2995f80ce2aSJacob Faibussowitsch       CHKERRQ(MatDestroySubMatrices(osm->n_local_true,&osm->pmat));
300b58cb649SBarry Smith     }
301b58cb649SBarry Smith     scall = MAT_INITIAL_MATRIX;
302b58cb649SBarry Smith   }
303b58cb649SBarry Smith 
30492bb6962SLisandro Dalcin   /*
30592bb6962SLisandro Dalcin      Extract out the submatrices
30692bb6962SLisandro Dalcin   */
3075f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSubMatrices(pc->pmat,osm->n_local_true,osm->is,osm->is,scall,&osm->pmat));
30892bb6962SLisandro Dalcin   if (scall == MAT_INITIAL_MATRIX) {
3095f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix));
31092bb6962SLisandro Dalcin     for (i=0; i<osm->n_local_true; i++) {
3115f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscLogObjectParent((PetscObject)pc,(PetscObject)osm->pmat[i]));
3125f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i],pprefix));
31392bb6962SLisandro Dalcin     }
31492bb6962SLisandro Dalcin   }
31580ec0b7dSPatrick Sanan 
31680ec0b7dSPatrick Sanan   /* Convert the types of the submatrices (if needbe) */
31780ec0b7dSPatrick Sanan   if (osm->sub_mat_type) {
31880ec0b7dSPatrick Sanan     for (i=0; i<osm->n_local_true; i++) {
3195f80ce2aSJacob Faibussowitsch       CHKERRQ(MatConvert(osm->pmat[i],osm->sub_mat_type,MAT_INPLACE_MATRIX,&(osm->pmat[i])));
32080ec0b7dSPatrick Sanan     }
32180ec0b7dSPatrick Sanan   }
32280ec0b7dSPatrick Sanan 
32380ec0b7dSPatrick Sanan   if (!pc->setupcalled) {
32456ea09ceSStefano Zampini     VecType vtype;
32556ea09ceSStefano Zampini 
32680ec0b7dSPatrick Sanan     /* Create the local work vectors (from the local matrices) and scatter contexts */
3275f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateVecs(pc->pmat,&vec,NULL));
3285b3c0d42Seaulisa 
3292c71b3e2SJacob Faibussowitsch     PetscCheckFalse(osm->is_local && (osm->type == PC_ASM_INTERPOLATE || osm->type == PC_ASM_NONE),PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot use interpolate or none PCASMType if is_local was provided to PCASMSetLocalSubdomains()");
3301dd8081eSeaulisa     if (osm->is_local && osm->type == PC_ASM_RESTRICT && osm->loctype == PC_COMPOSITE_ADDITIVE) {
3315f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscMalloc1(osm->n_local_true,&osm->lprolongation));
3321dd8081eSeaulisa     }
3335f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(osm->n_local_true,&osm->lrestriction));
3345f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(osm->n_local_true,&osm->x));
3355f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(osm->n_local_true,&osm->y));
336347574c9Seaulisa 
3375f80ce2aSJacob Faibussowitsch     CHKERRQ(ISGetLocalSize(osm->lis,&m));
3385f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateStride(PETSC_COMM_SELF,m,0,1,&isl));
3395f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetVecType(osm->pmat[0],&vtype));
3405f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCreate(PETSC_COMM_SELF,&osm->lx));
3415f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetSizes(osm->lx,m,m));
3425f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetType(osm->lx,vtype));
3435f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDuplicate(osm->lx, &osm->ly));
3445f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterCreate(vec,osm->lis,osm->lx,isl,&osm->restriction));
3455f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&isl));
346347574c9Seaulisa 
34780ec0b7dSPatrick Sanan     for (i=0; i<osm->n_local_true; ++i) {
3485b3c0d42Seaulisa       ISLocalToGlobalMapping ltog;
3495b3c0d42Seaulisa       IS                     isll;
3505b3c0d42Seaulisa       const PetscInt         *idx_is;
3515b3c0d42Seaulisa       PetscInt               *idx_lis,nout;
3525b3c0d42Seaulisa 
3535f80ce2aSJacob Faibussowitsch       CHKERRQ(ISGetLocalSize(osm->is[i],&m));
3545f80ce2aSJacob Faibussowitsch       CHKERRQ(MatCreateVecs(osm->pmat[i],&osm->x[i],NULL));
3555f80ce2aSJacob Faibussowitsch       CHKERRQ(VecDuplicate(osm->x[i],&osm->y[i]));
35655817e79SBarry Smith 
357b0de9f23SBarry Smith       /* generate a scatter from ly to y[i] picking all the overlapping is[i] entries */
3585f80ce2aSJacob Faibussowitsch       CHKERRQ(ISLocalToGlobalMappingCreateIS(osm->lis,&ltog));
3595f80ce2aSJacob Faibussowitsch       CHKERRQ(ISGetLocalSize(osm->is[i],&m));
3605f80ce2aSJacob Faibussowitsch       CHKERRQ(ISGetIndices(osm->is[i], &idx_is));
3615f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscMalloc1(m,&idx_lis));
3625f80ce2aSJacob Faibussowitsch       CHKERRQ(ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,m,idx_is,&nout,idx_lis));
3632c71b3e2SJacob Faibussowitsch       PetscCheckFalse(nout != m,PETSC_COMM_SELF,PETSC_ERR_PLIB,"is not a subset of lis");
3645f80ce2aSJacob Faibussowitsch       CHKERRQ(ISRestoreIndices(osm->is[i], &idx_is));
3655f80ce2aSJacob Faibussowitsch       CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,m,idx_lis,PETSC_OWN_POINTER,&isll));
3665f80ce2aSJacob Faibussowitsch       CHKERRQ(ISLocalToGlobalMappingDestroy(&ltog));
3675f80ce2aSJacob Faibussowitsch       CHKERRQ(ISCreateStride(PETSC_COMM_SELF,m,0,1,&isl));
3685f80ce2aSJacob Faibussowitsch       CHKERRQ(VecScatterCreate(osm->ly,isll,osm->y[i],isl,&osm->lrestriction[i]));
3695f80ce2aSJacob Faibussowitsch       CHKERRQ(ISDestroy(&isll));
3705f80ce2aSJacob Faibussowitsch       CHKERRQ(ISDestroy(&isl));
371910cf402Sprj-       if (osm->lprolongation) { /* generate a scatter from y[i] to ly picking only the the non-overlapping is_local[i] entries */
372d8b3b5e3Seaulisa         ISLocalToGlobalMapping ltog;
373d8b3b5e3Seaulisa         IS                     isll,isll_local;
374d8b3b5e3Seaulisa         const PetscInt         *idx_local;
375d8b3b5e3Seaulisa         PetscInt               *idx1, *idx2, nout;
376d8b3b5e3Seaulisa 
3775f80ce2aSJacob Faibussowitsch         CHKERRQ(ISGetLocalSize(osm->is_local[i],&m_local));
3785f80ce2aSJacob Faibussowitsch         CHKERRQ(ISGetIndices(osm->is_local[i], &idx_local));
379d8b3b5e3Seaulisa 
3805f80ce2aSJacob Faibussowitsch         CHKERRQ(ISLocalToGlobalMappingCreateIS(osm->is[i],&ltog));
3815f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscMalloc1(m_local,&idx1));
3825f80ce2aSJacob Faibussowitsch         CHKERRQ(ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,m_local,idx_local,&nout,idx1));
3835f80ce2aSJacob Faibussowitsch         CHKERRQ(ISLocalToGlobalMappingDestroy(&ltog));
3842c71b3e2SJacob Faibussowitsch         PetscCheckFalse(nout != m_local,PETSC_COMM_SELF,PETSC_ERR_PLIB,"is_local not a subset of is");
3855f80ce2aSJacob Faibussowitsch         CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,m_local,idx1,PETSC_OWN_POINTER,&isll));
386d8b3b5e3Seaulisa 
3875f80ce2aSJacob Faibussowitsch         CHKERRQ(ISLocalToGlobalMappingCreateIS(osm->lis,&ltog));
3885f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscMalloc1(m_local,&idx2));
3895f80ce2aSJacob Faibussowitsch         CHKERRQ(ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,m_local,idx_local,&nout,idx2));
3905f80ce2aSJacob Faibussowitsch         CHKERRQ(ISLocalToGlobalMappingDestroy(&ltog));
3912c71b3e2SJacob Faibussowitsch         PetscCheckFalse(nout != m_local,PETSC_COMM_SELF,PETSC_ERR_PLIB,"is_local not a subset of lis");
3925f80ce2aSJacob Faibussowitsch         CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,m_local,idx2,PETSC_OWN_POINTER,&isll_local));
393d8b3b5e3Seaulisa 
3945f80ce2aSJacob Faibussowitsch         CHKERRQ(ISRestoreIndices(osm->is_local[i], &idx_local));
3955f80ce2aSJacob Faibussowitsch         CHKERRQ(VecScatterCreate(osm->y[i],isll,osm->ly,isll_local,&osm->lprolongation[i]));
396d8b3b5e3Seaulisa 
3975f80ce2aSJacob Faibussowitsch         CHKERRQ(ISDestroy(&isll));
3985f80ce2aSJacob Faibussowitsch         CHKERRQ(ISDestroy(&isll_local));
399d8b3b5e3Seaulisa       }
40080ec0b7dSPatrick Sanan     }
4015f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&vec));
40280ec0b7dSPatrick Sanan   }
40380ec0b7dSPatrick Sanan 
404fb745f2cSMatthew G. Knepley   if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) {
405235cc4ceSMatthew G. Knepley     IS      *cis;
406235cc4ceSMatthew G. Knepley     PetscInt c;
407235cc4ceSMatthew G. Knepley 
4085f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(osm->n_local_true, &cis));
409235cc4ceSMatthew G. Knepley     for (c = 0; c < osm->n_local_true; ++c) cis[c] = osm->lis;
4105f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateSubMatrices(pc->pmat, osm->n_local_true, osm->is, cis, scall, &osm->lmats));
4115f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(cis));
412fb745f2cSMatthew G. Knepley   }
4134b9ad928SBarry Smith 
4144b9ad928SBarry Smith   /* Return control to the user so that the submatrices can be modified (e.g., to apply
4154b9ad928SBarry Smith      different boundary conditions for the submatrices than for the global problem) */
4165f80ce2aSJacob Faibussowitsch   CHKERRQ(PCModifySubMatrices(pc,osm->n_local_true,osm->is,osm->is,osm->pmat,pc->modifysubmatricesP));
4174b9ad928SBarry Smith 
41892bb6962SLisandro Dalcin   /*
41992bb6962SLisandro Dalcin      Loop over subdomains putting them into local ksp
42092bb6962SLisandro Dalcin   */
4215f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPGetOptionsPrefix(osm->ksp[0],&prefix));
42292bb6962SLisandro Dalcin   for (i=0; i<osm->n_local_true; i++) {
4235f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPSetOperators(osm->ksp[i],osm->pmat[i],osm->pmat[i]));
4245f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetOptionsPrefix(osm->pmat[i],prefix));
42592bb6962SLisandro Dalcin     if (!pc->setupcalled) {
4265f80ce2aSJacob Faibussowitsch       CHKERRQ(KSPSetFromOptions(osm->ksp[i]));
4274b9ad928SBarry Smith     }
42892bb6962SLisandro Dalcin   }
4294b9ad928SBarry Smith   PetscFunctionReturn(0);
4304b9ad928SBarry Smith }
4314b9ad928SBarry Smith 
4326849ba73SBarry Smith static PetscErrorCode PCSetUpOnBlocks_ASM(PC pc)
4334b9ad928SBarry Smith {
4344b9ad928SBarry Smith   PC_ASM             *osm = (PC_ASM*)pc->data;
43513f74950SBarry Smith   PetscInt           i;
436539c167fSBarry Smith   KSPConvergedReason reason;
4374b9ad928SBarry Smith 
4384b9ad928SBarry Smith   PetscFunctionBegin;
4394b9ad928SBarry Smith   for (i=0; i<osm->n_local_true; i++) {
4405f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPSetUp(osm->ksp[i]));
4415f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPGetConvergedReason(osm->ksp[i],&reason));
442c0decd05SBarry Smith     if (reason == KSP_DIVERGED_PC_FAILED) {
443261222b2SHong Zhang       pc->failedreason = PC_SUBPC_ERROR;
444e0eafd54SHong Zhang     }
4454b9ad928SBarry Smith   }
4464b9ad928SBarry Smith   PetscFunctionReturn(0);
4474b9ad928SBarry Smith }
4484b9ad928SBarry Smith 
4496849ba73SBarry Smith static PetscErrorCode PCApply_ASM(PC pc,Vec x,Vec y)
4504b9ad928SBarry Smith {
4514b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
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 */
4635f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSet(osm->lx, 0.0));
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 */
4715f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSet(y, 0.0));
4725f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSet(osm->ly, 0.0));
473347574c9Seaulisa 
47448e38a8aSPierre Jolivet     /* copy the global RHS to local RHS including the ghost nodes */
4755f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterBegin(osm->restriction, x, osm->lx, INSERT_VALUES, forward));
4765f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterEnd(osm->restriction, x, osm->lx, INSERT_VALUES, forward));
477347574c9Seaulisa 
47848e38a8aSPierre Jolivet     /* restrict local RHS to the overlapping 0-block RHS */
4795f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterBegin(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward));
4805f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterEnd(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward));
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 */
4865f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscLogEventBegin(PC_ApplyOnBlocks, osm->ksp[i], osm->x[i], osm->y[i],0));
4875f80ce2aSJacob Faibussowitsch       CHKERRQ(KSPSolve(osm->ksp[i], osm->x[i], osm->y[i]));
4885f80ce2aSJacob Faibussowitsch       CHKERRQ(KSPCheckSolve(osm->ksp[i], pc, osm->y[i]));
4895f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscLogEventEnd(PC_ApplyOnBlocks, osm->ksp[i], osm->x[i], osm->y[i], 0));
490d8b3b5e3Seaulisa 
491910cf402Sprj-       if (osm->lprolongation) { /* interpolate the non-overlapping i-block solution to the local solution (only for restrictive additive) */
4925f80ce2aSJacob Faibussowitsch         CHKERRQ(VecScatterBegin(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward));
4935f80ce2aSJacob Faibussowitsch         CHKERRQ(VecScatterEnd(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward));
49448e38a8aSPierre Jolivet       } else { /* interpolate the overlapping i-block solution to the local solution */
4955f80ce2aSJacob Faibussowitsch         CHKERRQ(VecScatterBegin(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse));
4965f80ce2aSJacob Faibussowitsch         CHKERRQ(VecScatterEnd(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse));
497d8b3b5e3Seaulisa       }
498347574c9Seaulisa 
499347574c9Seaulisa       if (i < n_local_true-1) {
50048e38a8aSPierre Jolivet         /* restrict local RHS to the overlapping (i+1)-block RHS */
5015f80ce2aSJacob Faibussowitsch         CHKERRQ(VecScatterBegin(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward));
5025f80ce2aSJacob Faibussowitsch         CHKERRQ(VecScatterEnd(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward));
503347574c9Seaulisa 
504347574c9Seaulisa         if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) {
5058966356dSPierre Jolivet           /* update the overlapping (i+1)-block RHS using the current local solution */
5065f80ce2aSJacob Faibussowitsch           CHKERRQ(MatMult(osm->lmats[i+1], osm->ly, osm->y[i+1]));
5075f80ce2aSJacob Faibussowitsch           CHKERRQ(VecAXPBY(osm->x[i+1],-1.,1., osm->y[i+1]));
5087c3d802fSMatthew G. Knepley         }
50912cd4985SMatthew G. Knepley       }
51012cd4985SMatthew G. Knepley     }
51148e38a8aSPierre Jolivet     /* add the local solution to the global solution including the ghost nodes */
5125f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterBegin(osm->restriction, osm->ly, y, ADD_VALUES, reverse));
5135f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterEnd(osm->restriction, osm->ly, y, ADD_VALUES, reverse));
51498921bdaSJacob Faibussowitsch   } else SETERRQ(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 
52648e38a8aSPierre Jolivet   PetscFunctionBegin;
5272c71b3e2SJacob Faibussowitsch   PetscCheckFalse(osm->n_local_true > 1,PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Not yet implemented");
52848e38a8aSPierre Jolivet   /*
52948e38a8aSPierre Jolivet      support for limiting the restriction or interpolation to only local
53048e38a8aSPierre Jolivet      subdomain values (leaving the other values 0).
53148e38a8aSPierre Jolivet   */
53248e38a8aSPierre Jolivet   if (!(osm->type & PC_ASM_RESTRICT)) {
53348e38a8aSPierre Jolivet     forward = SCATTER_FORWARD_LOCAL;
53448e38a8aSPierre Jolivet     /* have to zero the work RHS since scatter may leave some slots empty */
5355f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSet(osm->lx, 0.0));
53648e38a8aSPierre Jolivet   }
53748e38a8aSPierre Jolivet   if (!(osm->type & PC_ASM_INTERPOLATE)) {
53848e38a8aSPierre Jolivet     reverse = SCATTER_REVERSE_LOCAL;
53948e38a8aSPierre Jolivet   }
5405f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetLocalSize(osm->x[0], &m));
5415f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetSize(X, NULL, &N));
5425f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSeqDense(PETSC_COMM_SELF, m, N, NULL, &Z));
54348e38a8aSPierre Jolivet   if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE || osm->loctype == PC_COMPOSITE_ADDITIVE) {
54448e38a8aSPierre Jolivet     /* zero the global and the local solutions */
5455f80ce2aSJacob Faibussowitsch     CHKERRQ(MatZeroEntries(Y));
5465f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSet(osm->ly, 0.0));
54748e38a8aSPierre Jolivet 
54848e38a8aSPierre Jolivet     for (i = 0; i < N; ++i) {
5495f80ce2aSJacob Faibussowitsch       CHKERRQ(MatDenseGetColumnVecRead(X, i, &x));
55048e38a8aSPierre Jolivet       /* copy the global RHS to local RHS including the ghost nodes */
5515f80ce2aSJacob Faibussowitsch       CHKERRQ(VecScatterBegin(osm->restriction, x, osm->lx, INSERT_VALUES, forward));
5525f80ce2aSJacob Faibussowitsch       CHKERRQ(VecScatterEnd(osm->restriction, x, osm->lx, INSERT_VALUES, forward));
5535f80ce2aSJacob Faibussowitsch       CHKERRQ(MatDenseRestoreColumnVecRead(X, i, &x));
55448e38a8aSPierre Jolivet 
5555f80ce2aSJacob Faibussowitsch       CHKERRQ(MatDenseGetColumnVecWrite(Z, i, &x));
55648e38a8aSPierre Jolivet       /* restrict local RHS to the overlapping 0-block RHS */
5575f80ce2aSJacob Faibussowitsch       CHKERRQ(VecScatterBegin(osm->lrestriction[0], osm->lx, x, INSERT_VALUES, forward));
5585f80ce2aSJacob Faibussowitsch       CHKERRQ(VecScatterEnd(osm->lrestriction[0], osm->lx, x, INSERT_VALUES, forward));
5595f80ce2aSJacob Faibussowitsch       CHKERRQ(MatDenseRestoreColumnVecWrite(Z, i, &x));
56048e38a8aSPierre Jolivet     }
5615f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateSeqDense(PETSC_COMM_SELF, m, N, NULL, &W));
56248e38a8aSPierre Jolivet     /* solve the overlapping 0-block */
5635f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscLogEventBegin(PC_ApplyOnBlocks, osm->ksp[0], Z, W, 0));
5645f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPMatSolve(osm->ksp[0], Z, W));
5655f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPCheckSolve(osm->ksp[0], pc, NULL));
5665f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscLogEventEnd(PC_ApplyOnBlocks, osm->ksp[0], Z, W,0));
5675f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&Z));
56848e38a8aSPierre Jolivet 
56948e38a8aSPierre Jolivet     for (i = 0; i < N; ++i) {
5705f80ce2aSJacob Faibussowitsch       CHKERRQ(VecSet(osm->ly, 0.0));
5715f80ce2aSJacob Faibussowitsch       CHKERRQ(MatDenseGetColumnVecRead(W, i, &x));
57248e38a8aSPierre Jolivet       if (osm->lprolongation) { /* interpolate the non-overlapping 0-block solution to the local solution (only for restrictive additive) */
5735f80ce2aSJacob Faibussowitsch         CHKERRQ(VecScatterBegin(osm->lprolongation[0], x, osm->ly, ADD_VALUES, forward));
5745f80ce2aSJacob Faibussowitsch         CHKERRQ(VecScatterEnd(osm->lprolongation[0], x, osm->ly, ADD_VALUES, forward));
57548e38a8aSPierre Jolivet       } else { /* interpolate the overlapping 0-block solution to the local solution */
5765f80ce2aSJacob Faibussowitsch         CHKERRQ(VecScatterBegin(osm->lrestriction[0], x, osm->ly, ADD_VALUES, reverse));
5775f80ce2aSJacob Faibussowitsch         CHKERRQ(VecScatterEnd(osm->lrestriction[0], x, osm->ly, ADD_VALUES, reverse));
57848e38a8aSPierre Jolivet       }
5795f80ce2aSJacob Faibussowitsch       CHKERRQ(MatDenseRestoreColumnVecRead(W, i, &x));
58048e38a8aSPierre Jolivet 
5815f80ce2aSJacob Faibussowitsch       CHKERRQ(MatDenseGetColumnVecWrite(Y, i, &x));
58248e38a8aSPierre Jolivet       /* add the local solution to the global solution including the ghost nodes */
5835f80ce2aSJacob Faibussowitsch       CHKERRQ(VecScatterBegin(osm->restriction, osm->ly, x, ADD_VALUES, reverse));
5845f80ce2aSJacob Faibussowitsch       CHKERRQ(VecScatterEnd(osm->restriction, osm->ly, x, ADD_VALUES, reverse));
5855f80ce2aSJacob Faibussowitsch       CHKERRQ(MatDenseRestoreColumnVecWrite(Y, i, &x));
58648e38a8aSPierre Jolivet     }
5875f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&W));
58898921bdaSJacob Faibussowitsch   } else SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Invalid local composition type: %s", PCCompositeTypes[osm->loctype]);
58948e38a8aSPierre Jolivet   PetscFunctionReturn(0);
59048e38a8aSPierre Jolivet }
59148e38a8aSPierre Jolivet 
5926849ba73SBarry Smith static PetscErrorCode PCApplyTranspose_ASM(PC pc,Vec x,Vec y)
5934b9ad928SBarry Smith {
5944b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
5951dd8081eSeaulisa   PetscInt       i,n_local_true = osm->n_local_true;
5964b9ad928SBarry Smith   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;
5974b9ad928SBarry Smith 
5984b9ad928SBarry Smith   PetscFunctionBegin;
5994b9ad928SBarry Smith   /*
6004b9ad928SBarry Smith      Support for limiting the restriction or interpolation to only local
6014b9ad928SBarry Smith      subdomain values (leaving the other values 0).
6024b9ad928SBarry Smith 
6034b9ad928SBarry Smith      Note: these are reversed from the PCApply_ASM() because we are applying the
6044b9ad928SBarry Smith      transpose of the three terms
6054b9ad928SBarry Smith   */
606d8b3b5e3Seaulisa 
6074b9ad928SBarry Smith   if (!(osm->type & PC_ASM_INTERPOLATE)) {
6084b9ad928SBarry Smith     forward = SCATTER_FORWARD_LOCAL;
6094b9ad928SBarry Smith     /* have to zero the work RHS since scatter may leave some slots empty */
6105f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSet(osm->lx, 0.0));
6114b9ad928SBarry Smith   }
6122fa5cd67SKarl Rupp   if (!(osm->type & PC_ASM_RESTRICT)) reverse = SCATTER_REVERSE_LOCAL;
6134b9ad928SBarry Smith 
614b0de9f23SBarry Smith   /* zero the global and the local solutions */
6155f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(y, 0.0));
6165f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(osm->ly, 0.0));
617d8b3b5e3Seaulisa 
618b0de9f23SBarry Smith   /* Copy the global RHS to local RHS including the ghost nodes */
6195f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(osm->restriction, x, osm->lx, INSERT_VALUES, forward));
6205f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(osm->restriction, x, osm->lx, INSERT_VALUES, forward));
621d8b3b5e3Seaulisa 
622b0de9f23SBarry Smith   /* Restrict local RHS to the overlapping 0-block RHS */
6235f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward));
6245f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward));
625d8b3b5e3Seaulisa 
6264b9ad928SBarry Smith   /* do the local solves */
627d8b3b5e3Seaulisa   for (i = 0; i < n_local_true; ++i) {
628d8b3b5e3Seaulisa 
629b0de9f23SBarry Smith     /* solve the overlapping i-block */
6305f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscLogEventBegin(PC_ApplyOnBlocks,osm->ksp[i],osm->x[i],osm->y[i],0));
6315f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPSolveTranspose(osm->ksp[i], osm->x[i], osm->y[i]));
6325f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPCheckSolve(osm->ksp[i],pc,osm->y[i]));
6335f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscLogEventEnd(PC_ApplyOnBlocks,osm->ksp[i],osm->x[i],osm->y[i],0));
634d8b3b5e3Seaulisa 
635910cf402Sprj-     if (osm->lprolongation) { /* interpolate the non-overlapping i-block solution to the local solution */
6365f80ce2aSJacob Faibussowitsch       CHKERRQ(VecScatterBegin(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward));
6375f80ce2aSJacob Faibussowitsch       CHKERRQ(VecScatterEnd(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward));
638910cf402Sprj-     } else { /* interpolate the overlapping i-block solution to the local solution */
6395f80ce2aSJacob Faibussowitsch       CHKERRQ(VecScatterBegin(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse));
6405f80ce2aSJacob Faibussowitsch       CHKERRQ(VecScatterEnd(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse));
6414b9ad928SBarry Smith     }
642d8b3b5e3Seaulisa 
643d8b3b5e3Seaulisa     if (i < n_local_true-1) {
644b0de9f23SBarry Smith       /* Restrict local RHS to the overlapping (i+1)-block RHS */
6455f80ce2aSJacob Faibussowitsch       CHKERRQ(VecScatterBegin(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward));
6465f80ce2aSJacob Faibussowitsch       CHKERRQ(VecScatterEnd(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward));
6474b9ad928SBarry Smith     }
6484b9ad928SBarry Smith   }
649b0de9f23SBarry Smith   /* Add the local solution to the global solution including the ghost nodes */
6505f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(osm->restriction, osm->ly, y, ADD_VALUES, reverse));
6515f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(osm->restriction, osm->ly, y, ADD_VALUES, reverse));
652d8b3b5e3Seaulisa 
6534b9ad928SBarry Smith   PetscFunctionReturn(0);
654d8b3b5e3Seaulisa 
6554b9ad928SBarry Smith }
6564b9ad928SBarry Smith 
657e91c6855SBarry Smith static PetscErrorCode PCReset_ASM(PC pc)
6584b9ad928SBarry Smith {
6594b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
66013f74950SBarry Smith   PetscInt       i;
6614b9ad928SBarry Smith 
6624b9ad928SBarry Smith   PetscFunctionBegin;
66392bb6962SLisandro Dalcin   if (osm->ksp) {
66492bb6962SLisandro Dalcin     for (i=0; i<osm->n_local_true; i++) {
6655f80ce2aSJacob Faibussowitsch       CHKERRQ(KSPReset(osm->ksp[i]));
66692bb6962SLisandro Dalcin     }
66792bb6962SLisandro Dalcin   }
668e09e08ccSBarry Smith   if (osm->pmat) {
66992bb6962SLisandro Dalcin     if (osm->n_local_true > 0) {
6705f80ce2aSJacob Faibussowitsch       CHKERRQ(MatDestroySubMatrices(osm->n_local_true,&osm->pmat));
67192bb6962SLisandro Dalcin     }
67292bb6962SLisandro Dalcin   }
6731dd8081eSeaulisa   if (osm->lrestriction) {
6745f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterDestroy(&osm->restriction));
6751dd8081eSeaulisa     for (i=0; i<osm->n_local_true; i++) {
6765f80ce2aSJacob Faibussowitsch       CHKERRQ(VecScatterDestroy(&osm->lrestriction[i]));
6775f80ce2aSJacob Faibussowitsch       if (osm->lprolongation) CHKERRQ(VecScatterDestroy(&osm->lprolongation[i]));
6785f80ce2aSJacob Faibussowitsch       CHKERRQ(VecDestroy(&osm->x[i]));
6795f80ce2aSJacob Faibussowitsch       CHKERRQ(VecDestroy(&osm->y[i]));
6804b9ad928SBarry Smith     }
6815f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(osm->lrestriction));
6825f80ce2aSJacob Faibussowitsch     if (osm->lprolongation) CHKERRQ(PetscFree(osm->lprolongation));
6835f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(osm->x));
6845f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(osm->y));
6851dd8081eSeaulisa 
68692bb6962SLisandro Dalcin   }
6875f80ce2aSJacob Faibussowitsch   CHKERRQ(PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local));
6885f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&osm->lis));
6895f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&osm->lx));
6905f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&osm->ly));
691347574c9Seaulisa   if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) {
6925f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroyMatrices(osm->n_local_true, &osm->lmats));
693fb745f2cSMatthew G. Knepley   }
6942fa5cd67SKarl Rupp 
6955f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(osm->sub_mat_type));
69680ec0b7dSPatrick Sanan 
6970a545947SLisandro Dalcin   osm->is       = NULL;
6980a545947SLisandro Dalcin   osm->is_local = NULL;
699e91c6855SBarry Smith   PetscFunctionReturn(0);
700e91c6855SBarry Smith }
701e91c6855SBarry Smith 
702e91c6855SBarry Smith static PetscErrorCode PCDestroy_ASM(PC pc)
703e91c6855SBarry Smith {
704e91c6855SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
705e91c6855SBarry Smith   PetscInt       i;
706e91c6855SBarry Smith 
707e91c6855SBarry Smith   PetscFunctionBegin;
7085f80ce2aSJacob Faibussowitsch   CHKERRQ(PCReset_ASM(pc));
709e91c6855SBarry Smith   if (osm->ksp) {
710e91c6855SBarry Smith     for (i=0; i<osm->n_local_true; i++) {
7115f80ce2aSJacob Faibussowitsch       CHKERRQ(KSPDestroy(&osm->ksp[i]));
712e91c6855SBarry Smith     }
7135f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(osm->ksp));
714e91c6855SBarry Smith   }
7155f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(pc->data));
71696322394SPierre Jolivet 
7175f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectComposeFunction((PetscObject)pc,"PCASMSetLocalSubdomains_C",NULL));
7185f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectComposeFunction((PetscObject)pc,"PCASMSetTotalSubdomains_C",NULL));
7195f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectComposeFunction((PetscObject)pc,"PCASMSetOverlap_C",NULL));
7205f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectComposeFunction((PetscObject)pc,"PCASMSetType_C",NULL));
7215f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectComposeFunction((PetscObject)pc,"PCASMGetType_C",NULL));
7225f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectComposeFunction((PetscObject)pc,"PCASMSetLocalType_C",NULL));
7235f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectComposeFunction((PetscObject)pc,"PCASMGetLocalType_C",NULL));
7245f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectComposeFunction((PetscObject)pc,"PCASMSetSortIndices_C",NULL));
7255f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectComposeFunction((PetscObject)pc,"PCASMGetSubKSP_C",NULL));
7265f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectComposeFunction((PetscObject)pc,"PCASMGetSubMatType_C",NULL));
7275f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectComposeFunction((PetscObject)pc,"PCASMSetSubMatType_C",NULL));
7284b9ad928SBarry Smith   PetscFunctionReturn(0);
7294b9ad928SBarry Smith }
7304b9ad928SBarry Smith 
7314416b707SBarry Smith static PetscErrorCode PCSetFromOptions_ASM(PetscOptionItems *PetscOptionsObject,PC pc)
7324b9ad928SBarry Smith {
7334b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
7349dcbbd2bSBarry Smith   PetscInt       blocks,ovl;
73557501b6eSBarry Smith   PetscBool      flg;
73692bb6962SLisandro Dalcin   PCASMType      asmtype;
73712cd4985SMatthew G. Knepley   PCCompositeType loctype;
73880ec0b7dSPatrick Sanan   char           sub_mat_type[256];
7394b9ad928SBarry Smith 
7404b9ad928SBarry Smith   PetscFunctionBegin;
7415f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHead(PetscOptionsObject,"Additive Schwarz options"));
7425f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBool("-pc_asm_dm_subdomains","Use DMCreateDomainDecomposition() to define subdomains","PCASMSetDMSubdomains",osm->dm_subdomains,&osm->dm_subdomains,&flg));
7435f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsInt("-pc_asm_blocks","Number of subdomains","PCASMSetTotalSubdomains",osm->n,&blocks,&flg));
74465db9045SDmitry Karpeev   if (flg) {
7455f80ce2aSJacob Faibussowitsch     CHKERRQ(PCASMSetTotalSubdomains(pc,blocks,NULL,NULL));
746d709ea83SDmitry Karpeev     osm->dm_subdomains = PETSC_FALSE;
74765db9045SDmitry Karpeev   }
7485f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsInt("-pc_asm_local_blocks","Number of local subdomains","PCASMSetLocalSubdomains",osm->n_local_true,&blocks,&flg));
749342c94f9SMatthew G. Knepley   if (flg) {
7505f80ce2aSJacob Faibussowitsch     CHKERRQ(PCASMSetLocalSubdomains(pc,blocks,NULL,NULL));
751342c94f9SMatthew G. Knepley     osm->dm_subdomains = PETSC_FALSE;
752342c94f9SMatthew G. Knepley   }
7535f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsInt("-pc_asm_overlap","Number of grid points overlap","PCASMSetOverlap",osm->overlap,&ovl,&flg));
75465db9045SDmitry Karpeev   if (flg) {
7555f80ce2aSJacob Faibussowitsch     CHKERRQ(PCASMSetOverlap(pc,ovl));
756d709ea83SDmitry Karpeev     osm->dm_subdomains = PETSC_FALSE;
75765db9045SDmitry Karpeev   }
75890d69ab7SBarry Smith   flg  = PETSC_FALSE;
7595f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsEnum("-pc_asm_type","Type of restriction/extension","PCASMSetType",PCASMTypes,(PetscEnum)osm->type,(PetscEnum*)&asmtype,&flg));
7605f80ce2aSJacob Faibussowitsch   if (flg) CHKERRQ(PCASMSetType(pc,asmtype));
76112cd4985SMatthew G. Knepley   flg  = PETSC_FALSE;
7625f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsEnum("-pc_asm_local_type","Type of local solver composition","PCASMSetLocalType",PCCompositeTypes,(PetscEnum)osm->loctype,(PetscEnum*)&loctype,&flg));
7635f80ce2aSJacob Faibussowitsch   if (flg) CHKERRQ(PCASMSetLocalType(pc,loctype));
7645f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsFList("-pc_asm_sub_mat_type","Subsolve Matrix Type","PCASMSetSubMatType",MatList,NULL,sub_mat_type,256,&flg));
76580ec0b7dSPatrick Sanan   if (flg) {
7665f80ce2aSJacob Faibussowitsch     CHKERRQ(PCASMSetSubMatType(pc,sub_mat_type));
76780ec0b7dSPatrick Sanan   }
7685f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsTail());
7694b9ad928SBarry Smith   PetscFunctionReturn(0);
7704b9ad928SBarry Smith }
7714b9ad928SBarry Smith 
7724b9ad928SBarry Smith /*------------------------------------------------------------------------------------*/
7734b9ad928SBarry Smith 
7741e6b0712SBarry Smith static PetscErrorCode  PCASMSetLocalSubdomains_ASM(PC pc,PetscInt n,IS is[],IS is_local[])
7754b9ad928SBarry Smith {
7764b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
77792bb6962SLisandro Dalcin   PetscInt       i;
7784b9ad928SBarry Smith 
7794b9ad928SBarry Smith   PetscFunctionBegin;
7802c71b3e2SJacob Faibussowitsch   PetscCheckFalse(n < 1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Each process must have 1 or more blocks, n = %D",n);
7812c71b3e2SJacob Faibussowitsch   PetscCheckFalse(pc->setupcalled && (n != osm->n_local_true || is),PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCASMSetLocalSubdomains() should be called before calling PCSetUp().");
782e7e72b3dSBarry Smith 
7834b9ad928SBarry Smith   if (!pc->setupcalled) {
78492bb6962SLisandro Dalcin     if (is) {
7855f80ce2aSJacob Faibussowitsch       for (i=0; i<n; i++) CHKERRQ(PetscObjectReference((PetscObject)is[i]));
78692bb6962SLisandro Dalcin     }
787832fc9a5SMatthew Knepley     if (is_local) {
7885f80ce2aSJacob Faibussowitsch       for (i=0; i<n; i++) CHKERRQ(PetscObjectReference((PetscObject)is_local[i]));
789832fc9a5SMatthew Knepley     }
7905f80ce2aSJacob Faibussowitsch     CHKERRQ(PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local));
7912fa5cd67SKarl Rupp 
7924b9ad928SBarry Smith     osm->n_local_true = n;
7930a545947SLisandro Dalcin     osm->is           = NULL;
7940a545947SLisandro Dalcin     osm->is_local     = NULL;
79592bb6962SLisandro Dalcin     if (is) {
7965f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscMalloc1(n,&osm->is));
7972fa5cd67SKarl Rupp       for (i=0; i<n; i++) osm->is[i] = is[i];
7983d03bb48SJed Brown       /* Flag indicating that the user has set overlapping subdomains so PCASM should not increase their size. */
7993d03bb48SJed Brown       osm->overlap = -1;
80092bb6962SLisandro Dalcin     }
8012b691e39SMatthew Knepley     if (is_local) {
8025f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscMalloc1(n,&osm->is_local));
8032fa5cd67SKarl Rupp       for (i=0; i<n; i++) osm->is_local[i] = is_local[i];
804a35b7b57SDmitry Karpeev       if (!is) {
8055f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscMalloc1(osm->n_local_true,&osm->is));
806a35b7b57SDmitry Karpeev         for (i=0; i<osm->n_local_true; i++) {
807a35b7b57SDmitry Karpeev           if (osm->overlap > 0) { /* With positive overlap, osm->is[i] will be modified */
8085f80ce2aSJacob Faibussowitsch             CHKERRQ(ISDuplicate(osm->is_local[i],&osm->is[i]));
8095f80ce2aSJacob Faibussowitsch             CHKERRQ(ISCopy(osm->is_local[i],osm->is[i]));
810a35b7b57SDmitry Karpeev           } else {
8115f80ce2aSJacob Faibussowitsch             CHKERRQ(PetscObjectReference((PetscObject)osm->is_local[i]));
812a35b7b57SDmitry Karpeev             osm->is[i] = osm->is_local[i];
813a35b7b57SDmitry Karpeev           }
814a35b7b57SDmitry Karpeev         }
815a35b7b57SDmitry Karpeev       }
8162b691e39SMatthew Knepley     }
8174b9ad928SBarry Smith   }
8184b9ad928SBarry Smith   PetscFunctionReturn(0);
8194b9ad928SBarry Smith }
8204b9ad928SBarry Smith 
8211e6b0712SBarry Smith static PetscErrorCode  PCASMSetTotalSubdomains_ASM(PC pc,PetscInt N,IS *is,IS *is_local)
8224b9ad928SBarry Smith {
8234b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
82413f74950SBarry Smith   PetscMPIInt    rank,size;
82578904715SLisandro Dalcin   PetscInt       n;
8264b9ad928SBarry Smith 
8274b9ad928SBarry Smith   PetscFunctionBegin;
8282c71b3e2SJacob Faibussowitsch   PetscCheckFalse(N < 1,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Number of total blocks must be > 0, N = %D",N);
8292c71b3e2SJacob Faibussowitsch   PetscCheckFalse(is || is_local,PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Use PCASMSetLocalSubdomains() to set specific index sets\n\they cannot be set globally yet.");
8304b9ad928SBarry Smith 
8314b9ad928SBarry Smith   /*
832880770e9SJed Brown      Split the subdomains equally among all processors
8334b9ad928SBarry Smith   */
8345f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank));
8355f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size));
8364b9ad928SBarry Smith   n    = N/size + ((N % size) > rank);
837*28b400f6SJacob Faibussowitsch   PetscCheck(n,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);
8382c71b3e2SJacob Faibussowitsch   PetscCheckFalse(pc->setupcalled && n != osm->n_local_true,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"PCASMSetTotalSubdomains() should be called before PCSetUp().");
8394b9ad928SBarry Smith   if (!pc->setupcalled) {
8405f80ce2aSJacob Faibussowitsch     CHKERRQ(PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local));
8412fa5cd67SKarl Rupp 
8424b9ad928SBarry Smith     osm->n_local_true = n;
8430a545947SLisandro Dalcin     osm->is           = NULL;
8440a545947SLisandro Dalcin     osm->is_local     = NULL;
8454b9ad928SBarry Smith   }
8464b9ad928SBarry Smith   PetscFunctionReturn(0);
8474b9ad928SBarry Smith }
8484b9ad928SBarry Smith 
8491e6b0712SBarry Smith static PetscErrorCode  PCASMSetOverlap_ASM(PC pc,PetscInt ovl)
8504b9ad928SBarry Smith {
85192bb6962SLisandro Dalcin   PC_ASM *osm = (PC_ASM*)pc->data;
8524b9ad928SBarry Smith 
8534b9ad928SBarry Smith   PetscFunctionBegin;
8542c71b3e2SJacob Faibussowitsch   PetscCheckFalse(ovl < 0,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested");
8552c71b3e2SJacob Faibussowitsch   PetscCheckFalse(pc->setupcalled && ovl != osm->overlap,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCASMSetOverlap() should be called before PCSetUp().");
8562fa5cd67SKarl Rupp   if (!pc->setupcalled) osm->overlap = ovl;
8574b9ad928SBarry Smith   PetscFunctionReturn(0);
8584b9ad928SBarry Smith }
8594b9ad928SBarry Smith 
8601e6b0712SBarry Smith static PetscErrorCode  PCASMSetType_ASM(PC pc,PCASMType type)
8614b9ad928SBarry Smith {
86292bb6962SLisandro Dalcin   PC_ASM *osm = (PC_ASM*)pc->data;
8634b9ad928SBarry Smith 
8644b9ad928SBarry Smith   PetscFunctionBegin;
8654b9ad928SBarry Smith   osm->type     = type;
866bf108f30SBarry Smith   osm->type_set = PETSC_TRUE;
8674b9ad928SBarry Smith   PetscFunctionReturn(0);
8684b9ad928SBarry Smith }
8694b9ad928SBarry Smith 
870c60c7ad4SBarry Smith static PetscErrorCode  PCASMGetType_ASM(PC pc,PCASMType *type)
871c60c7ad4SBarry Smith {
872c60c7ad4SBarry Smith   PC_ASM *osm = (PC_ASM*)pc->data;
873c60c7ad4SBarry Smith 
874c60c7ad4SBarry Smith   PetscFunctionBegin;
875c60c7ad4SBarry Smith   *type = osm->type;
876c60c7ad4SBarry Smith   PetscFunctionReturn(0);
877c60c7ad4SBarry Smith }
878c60c7ad4SBarry Smith 
87912cd4985SMatthew G. Knepley static PetscErrorCode  PCASMSetLocalType_ASM(PC pc, PCCompositeType type)
88012cd4985SMatthew G. Knepley {
88112cd4985SMatthew G. Knepley   PC_ASM *osm = (PC_ASM *) pc->data;
88212cd4985SMatthew G. Knepley 
88312cd4985SMatthew G. Knepley   PetscFunctionBegin;
8842c71b3e2SJacob Faibussowitsch   PetscCheckFalse(type != PC_COMPOSITE_ADDITIVE && type != PC_COMPOSITE_MULTIPLICATIVE,PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Only supports additive or multiplicative as the local type");
88512cd4985SMatthew G. Knepley   osm->loctype = type;
88612cd4985SMatthew G. Knepley   PetscFunctionReturn(0);
88712cd4985SMatthew G. Knepley }
88812cd4985SMatthew G. Knepley 
88912cd4985SMatthew G. Knepley static PetscErrorCode  PCASMGetLocalType_ASM(PC pc, PCCompositeType *type)
89012cd4985SMatthew G. Knepley {
89112cd4985SMatthew G. Knepley   PC_ASM *osm = (PC_ASM *) pc->data;
89212cd4985SMatthew G. Knepley 
89312cd4985SMatthew G. Knepley   PetscFunctionBegin;
89412cd4985SMatthew G. Knepley   *type = osm->loctype;
89512cd4985SMatthew G. Knepley   PetscFunctionReturn(0);
89612cd4985SMatthew G. Knepley }
89712cd4985SMatthew G. Knepley 
8981e6b0712SBarry Smith static PetscErrorCode  PCASMSetSortIndices_ASM(PC pc,PetscBool  doSort)
8996ed231c7SMatthew Knepley {
9006ed231c7SMatthew Knepley   PC_ASM *osm = (PC_ASM*)pc->data;
9016ed231c7SMatthew Knepley 
9026ed231c7SMatthew Knepley   PetscFunctionBegin;
9036ed231c7SMatthew Knepley   osm->sort_indices = doSort;
9046ed231c7SMatthew Knepley   PetscFunctionReturn(0);
9056ed231c7SMatthew Knepley }
9066ed231c7SMatthew Knepley 
9071e6b0712SBarry Smith static PetscErrorCode  PCASMGetSubKSP_ASM(PC pc,PetscInt *n_local,PetscInt *first_local,KSP **ksp)
9084b9ad928SBarry Smith {
90992bb6962SLisandro Dalcin   PC_ASM         *osm = (PC_ASM*)pc->data;
9104b9ad928SBarry Smith 
9114b9ad928SBarry Smith   PetscFunctionBegin;
9122c71b3e2SJacob Faibussowitsch   PetscCheckFalse(osm->n_local_true < 1,PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Need to call PCSetUp() on PC (or KSPSetUp() on the outer KSP object) before calling here");
9134b9ad928SBarry Smith 
9142fa5cd67SKarl Rupp   if (n_local) *n_local = osm->n_local_true;
91592bb6962SLisandro Dalcin   if (first_local) {
9165f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Scan(&osm->n_local_true,first_local,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)pc)));
91792bb6962SLisandro Dalcin     *first_local -= osm->n_local_true;
91892bb6962SLisandro Dalcin   }
919ed155784SPierre Jolivet   if (ksp) *ksp   = osm->ksp;
9204b9ad928SBarry Smith   PetscFunctionReturn(0);
9214b9ad928SBarry Smith }
9224b9ad928SBarry Smith 
92380ec0b7dSPatrick Sanan static PetscErrorCode  PCASMGetSubMatType_ASM(PC pc,MatType *sub_mat_type)
92480ec0b7dSPatrick Sanan {
92580ec0b7dSPatrick Sanan   PC_ASM         *osm = (PC_ASM*)pc->data;
92680ec0b7dSPatrick Sanan 
92780ec0b7dSPatrick Sanan   PetscFunctionBegin;
92880ec0b7dSPatrick Sanan   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
92980ec0b7dSPatrick Sanan   PetscValidPointer(sub_mat_type,2);
93080ec0b7dSPatrick Sanan   *sub_mat_type = osm->sub_mat_type;
93180ec0b7dSPatrick Sanan   PetscFunctionReturn(0);
93280ec0b7dSPatrick Sanan }
93380ec0b7dSPatrick Sanan 
93480ec0b7dSPatrick Sanan static PetscErrorCode PCASMSetSubMatType_ASM(PC pc,MatType sub_mat_type)
93580ec0b7dSPatrick Sanan {
93680ec0b7dSPatrick Sanan   PC_ASM            *osm = (PC_ASM*)pc->data;
93780ec0b7dSPatrick Sanan 
93880ec0b7dSPatrick Sanan   PetscFunctionBegin;
93980ec0b7dSPatrick Sanan   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
9405f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(osm->sub_mat_type));
9415f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrallocpy(sub_mat_type,(char**)&osm->sub_mat_type));
94280ec0b7dSPatrick Sanan   PetscFunctionReturn(0);
94380ec0b7dSPatrick Sanan }
94480ec0b7dSPatrick Sanan 
9454b9ad928SBarry Smith /*@C
9461093a601SBarry Smith     PCASMSetLocalSubdomains - Sets the local subdomains (for this processor only) for the additive Schwarz preconditioner.
9474b9ad928SBarry Smith 
948d083f849SBarry Smith     Collective on pc
9494b9ad928SBarry Smith 
9504b9ad928SBarry Smith     Input Parameters:
9514b9ad928SBarry Smith +   pc - the preconditioner context
9524b9ad928SBarry Smith .   n - the number of subdomains for this processor (default value = 1)
9538c03b21aSDmitry Karpeev .   is - the index set that defines the subdomains for this processor
9540298fd71SBarry Smith          (or NULL for PETSc to determine subdomains)
955f1ee410cSBarry 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
956f1ee410cSBarry Smith          (or NULL to not provide these)
9574b9ad928SBarry Smith 
958342c94f9SMatthew G. Knepley     Options Database Key:
959342c94f9SMatthew G. Knepley     To set the total number of subdomain blocks rather than specify the
960342c94f9SMatthew G. Knepley     index sets, use the option
961342c94f9SMatthew G. Knepley .    -pc_asm_local_blocks <blks> - Sets local blocks
962342c94f9SMatthew G. Knepley 
9634b9ad928SBarry Smith     Notes:
9641093a601SBarry Smith     The IS numbering is in the parallel, global numbering of the vector for both is and is_local
9654b9ad928SBarry Smith 
9664b9ad928SBarry Smith     By default the ASM preconditioner uses 1 block per processor.
9674b9ad928SBarry Smith 
9684b9ad928SBarry Smith     Use PCASMSetTotalSubdomains() to set the subdomains for all processors.
9694b9ad928SBarry Smith 
970f1ee410cSBarry Smith     If is_local is provided and PCASMType is PC_ASM_RESTRICT then the solution only over the is_local region is interpolated
971f1ee410cSBarry Smith     back to form the global solution (this is the standard restricted additive Schwarz method)
972f1ee410cSBarry Smith 
973f1ee410cSBarry 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
974f1ee410cSBarry Smith     no code to handle that case.
975f1ee410cSBarry Smith 
9764b9ad928SBarry Smith     Level: advanced
9774b9ad928SBarry Smith 
9784b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(),
979f1ee410cSBarry Smith           PCASMCreateSubdomains2D(), PCASMGetLocalSubdomains(), PCASMType, PCASMSetType()
9804b9ad928SBarry Smith @*/
9817087cfbeSBarry Smith PetscErrorCode  PCASMSetLocalSubdomains(PC pc,PetscInt n,IS is[],IS is_local[])
9824b9ad928SBarry Smith {
9834b9ad928SBarry Smith   PetscFunctionBegin;
9840700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
9855f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscTryMethod(pc,"PCASMSetLocalSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,n,is,is_local)));
9864b9ad928SBarry Smith   PetscFunctionReturn(0);
9874b9ad928SBarry Smith }
9884b9ad928SBarry Smith 
9894b9ad928SBarry Smith /*@C
990feb221f8SDmitry Karpeev     PCASMSetTotalSubdomains - Sets the subdomains for all processors for the
9914b9ad928SBarry Smith     additive Schwarz preconditioner.  Either all or no processors in the
9924b9ad928SBarry Smith     PC communicator must call this routine, with the same index sets.
9934b9ad928SBarry Smith 
994d083f849SBarry Smith     Collective on pc
9954b9ad928SBarry Smith 
9964b9ad928SBarry Smith     Input Parameters:
9974b9ad928SBarry Smith +   pc - the preconditioner context
998feb221f8SDmitry Karpeev .   N  - the number of subdomains for all processors
999feb221f8SDmitry Karpeev .   is - the index sets that define the subdomains for all processors
1000dfaa0c5fSBarry Smith          (or NULL to ask PETSc to determine the subdomains)
10012b691e39SMatthew Knepley -   is_local - the index sets that define the local part of the subdomains for this processor
1002f1ee410cSBarry Smith          (or NULL to not provide this information)
10034b9ad928SBarry Smith 
10044b9ad928SBarry Smith     Options Database Key:
10054b9ad928SBarry Smith     To set the total number of subdomain blocks rather than specify the
10064b9ad928SBarry Smith     index sets, use the option
10074b9ad928SBarry Smith .    -pc_asm_blocks <blks> - Sets total blocks
10084b9ad928SBarry Smith 
10094b9ad928SBarry Smith     Notes:
1010f1ee410cSBarry Smith     Currently you cannot use this to set the actual subdomains with the argument is or is_local.
10114b9ad928SBarry Smith 
10124b9ad928SBarry Smith     By default the ASM preconditioner uses 1 block per processor.
10134b9ad928SBarry Smith 
10144b9ad928SBarry Smith     These index sets cannot be destroyed until after completion of the
10154b9ad928SBarry Smith     linear solves for which the ASM preconditioner is being used.
10164b9ad928SBarry Smith 
10174b9ad928SBarry Smith     Use PCASMSetLocalSubdomains() to set local subdomains.
10184b9ad928SBarry Smith 
10191093a601SBarry Smith     The IS numbering is in the parallel, global numbering of the vector for both is and is_local
10201093a601SBarry Smith 
10214b9ad928SBarry Smith     Level: advanced
10224b9ad928SBarry Smith 
10234b9ad928SBarry Smith .seealso: PCASMSetLocalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(),
10244b9ad928SBarry Smith           PCASMCreateSubdomains2D()
10254b9ad928SBarry Smith @*/
10267087cfbeSBarry Smith PetscErrorCode  PCASMSetTotalSubdomains(PC pc,PetscInt N,IS is[],IS is_local[])
10274b9ad928SBarry Smith {
10284b9ad928SBarry Smith   PetscFunctionBegin;
10290700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
10305f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscTryMethod(pc,"PCASMSetTotalSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,N,is,is_local)));
10314b9ad928SBarry Smith   PetscFunctionReturn(0);
10324b9ad928SBarry Smith }
10334b9ad928SBarry Smith 
10344b9ad928SBarry Smith /*@
10354b9ad928SBarry Smith     PCASMSetOverlap - Sets the overlap between a pair of subdomains for the
10364b9ad928SBarry Smith     additive Schwarz preconditioner.  Either all or no processors in the
1037f1ee410cSBarry Smith     PC communicator must call this routine.
10384b9ad928SBarry Smith 
1039d083f849SBarry Smith     Logically Collective on pc
10404b9ad928SBarry Smith 
10414b9ad928SBarry Smith     Input Parameters:
10424b9ad928SBarry Smith +   pc  - the preconditioner context
10434b9ad928SBarry Smith -   ovl - the amount of overlap between subdomains (ovl >= 0, default value = 1)
10444b9ad928SBarry Smith 
10454b9ad928SBarry Smith     Options Database Key:
10464b9ad928SBarry Smith .   -pc_asm_overlap <ovl> - Sets overlap
10474b9ad928SBarry Smith 
10484b9ad928SBarry Smith     Notes:
10494b9ad928SBarry Smith     By default the ASM preconditioner uses 1 block per processor.  To use
10504b9ad928SBarry Smith     multiple blocks per perocessor, see PCASMSetTotalSubdomains() and
10514b9ad928SBarry Smith     PCASMSetLocalSubdomains() (and the option -pc_asm_blocks <blks>).
10524b9ad928SBarry Smith 
10534b9ad928SBarry Smith     The overlap defaults to 1, so if one desires that no additional
10544b9ad928SBarry Smith     overlap be computed beyond what may have been set with a call to
10554b9ad928SBarry Smith     PCASMSetTotalSubdomains() or PCASMSetLocalSubdomains(), then ovl
10564b9ad928SBarry Smith     must be set to be 0.  In particular, if one does not explicitly set
10574b9ad928SBarry Smith     the subdomains an application code, then all overlap would be computed
10584b9ad928SBarry Smith     internally by PETSc, and using an overlap of 0 would result in an ASM
10594b9ad928SBarry Smith     variant that is equivalent to the block Jacobi preconditioner.
10604b9ad928SBarry Smith 
1061f1ee410cSBarry Smith     The default algorithm used by PETSc to increase overlap is fast, but not scalable,
1062f1ee410cSBarry Smith     use the option -mat_increase_overlap_scalable when the problem and number of processes is large.
1063f1ee410cSBarry Smith 
10644b9ad928SBarry Smith     Note that one can define initial index sets with any overlap via
1065f1ee410cSBarry Smith     PCASMSetLocalSubdomains(); the routine
10664b9ad928SBarry Smith     PCASMSetOverlap() merely allows PETSc to extend that overlap further
10674b9ad928SBarry Smith     if desired.
10684b9ad928SBarry Smith 
10694b9ad928SBarry Smith     Level: intermediate
10704b9ad928SBarry Smith 
10714b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetLocalSubdomains(), PCASMGetSubKSP(),
1072f1ee410cSBarry Smith           PCASMCreateSubdomains2D(), PCASMGetLocalSubdomains(), MatIncreaseOverlap()
10734b9ad928SBarry Smith @*/
10747087cfbeSBarry Smith PetscErrorCode  PCASMSetOverlap(PC pc,PetscInt ovl)
10754b9ad928SBarry Smith {
10764b9ad928SBarry Smith   PetscFunctionBegin;
10770700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1078c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,ovl,2);
10795f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscTryMethod(pc,"PCASMSetOverlap_C",(PC,PetscInt),(pc,ovl)));
10804b9ad928SBarry Smith   PetscFunctionReturn(0);
10814b9ad928SBarry Smith }
10824b9ad928SBarry Smith 
10834b9ad928SBarry Smith /*@
10844b9ad928SBarry Smith     PCASMSetType - Sets the type of restriction and interpolation used
10854b9ad928SBarry Smith     for local problems in the additive Schwarz method.
10864b9ad928SBarry Smith 
1087d083f849SBarry Smith     Logically Collective on pc
10884b9ad928SBarry Smith 
10894b9ad928SBarry Smith     Input Parameters:
10904b9ad928SBarry Smith +   pc  - the preconditioner context
10914b9ad928SBarry Smith -   type - variant of ASM, one of
10924b9ad928SBarry Smith .vb
10934b9ad928SBarry Smith       PC_ASM_BASIC       - full interpolation and restriction
10944b9ad928SBarry Smith       PC_ASM_RESTRICT    - full restriction, local processor interpolation
10954b9ad928SBarry Smith       PC_ASM_INTERPOLATE - full interpolation, local processor restriction
10964b9ad928SBarry Smith       PC_ASM_NONE        - local processor restriction and interpolation
10974b9ad928SBarry Smith .ve
10984b9ad928SBarry Smith 
10994b9ad928SBarry Smith     Options Database Key:
11004b9ad928SBarry Smith .   -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type
11014b9ad928SBarry Smith 
110295452b02SPatrick Sanan     Notes:
110395452b02SPatrick Sanan     if the is_local arguments are passed to PCASMSetLocalSubdomains() then they are used when PC_ASM_RESTRICT has been selected
1104f1ee410cSBarry Smith     to limit the local processor interpolation
1105f1ee410cSBarry Smith 
11064b9ad928SBarry Smith     Level: intermediate
11074b9ad928SBarry Smith 
11084b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(),
1109f1ee410cSBarry Smith           PCASMCreateSubdomains2D(), PCASMType, PCASMSetLocalType(), PCASMGetLocalType()
11104b9ad928SBarry Smith @*/
11117087cfbeSBarry Smith PetscErrorCode  PCASMSetType(PC pc,PCASMType type)
11124b9ad928SBarry Smith {
11134b9ad928SBarry Smith   PetscFunctionBegin;
11140700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1115c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,type,2);
11165f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscTryMethod(pc,"PCASMSetType_C",(PC,PCASMType),(pc,type)));
11174b9ad928SBarry Smith   PetscFunctionReturn(0);
11184b9ad928SBarry Smith }
11194b9ad928SBarry Smith 
1120c60c7ad4SBarry Smith /*@
1121c60c7ad4SBarry Smith     PCASMGetType - Gets the type of restriction and interpolation used
1122c60c7ad4SBarry Smith     for local problems in the additive Schwarz method.
1123c60c7ad4SBarry Smith 
1124d083f849SBarry Smith     Logically Collective on pc
1125c60c7ad4SBarry Smith 
1126c60c7ad4SBarry Smith     Input Parameter:
1127c60c7ad4SBarry Smith .   pc  - the preconditioner context
1128c60c7ad4SBarry Smith 
1129c60c7ad4SBarry Smith     Output Parameter:
1130c60c7ad4SBarry Smith .   type - variant of ASM, one of
1131c60c7ad4SBarry Smith 
1132c60c7ad4SBarry Smith .vb
1133c60c7ad4SBarry Smith       PC_ASM_BASIC       - full interpolation and restriction
1134c60c7ad4SBarry Smith       PC_ASM_RESTRICT    - full restriction, local processor interpolation
1135c60c7ad4SBarry Smith       PC_ASM_INTERPOLATE - full interpolation, local processor restriction
1136c60c7ad4SBarry Smith       PC_ASM_NONE        - local processor restriction and interpolation
1137c60c7ad4SBarry Smith .ve
1138c60c7ad4SBarry Smith 
1139c60c7ad4SBarry Smith     Options Database Key:
1140c60c7ad4SBarry Smith .   -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type
1141c60c7ad4SBarry Smith 
1142c60c7ad4SBarry Smith     Level: intermediate
1143c60c7ad4SBarry Smith 
1144c60c7ad4SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(),
1145f1ee410cSBarry Smith           PCASMCreateSubdomains2D(), PCASMType, PCASMSetType(), PCASMSetLocalType(), PCASMGetLocalType()
1146c60c7ad4SBarry Smith @*/
1147c60c7ad4SBarry Smith PetscErrorCode  PCASMGetType(PC pc,PCASMType *type)
1148c60c7ad4SBarry Smith {
1149c60c7ad4SBarry Smith   PetscFunctionBegin;
1150c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
11515f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscUseMethod(pc,"PCASMGetType_C",(PC,PCASMType*),(pc,type)));
1152c60c7ad4SBarry Smith   PetscFunctionReturn(0);
1153c60c7ad4SBarry Smith }
1154c60c7ad4SBarry Smith 
115512cd4985SMatthew G. Knepley /*@
115612cd4985SMatthew G. Knepley   PCASMSetLocalType - Sets the type of composition used for local problems in the additive Schwarz method.
115712cd4985SMatthew G. Knepley 
1158d083f849SBarry Smith   Logically Collective on pc
115912cd4985SMatthew G. Knepley 
116012cd4985SMatthew G. Knepley   Input Parameters:
116112cd4985SMatthew G. Knepley + pc  - the preconditioner context
116212cd4985SMatthew G. Knepley - type - type of composition, one of
116312cd4985SMatthew G. Knepley .vb
116412cd4985SMatthew G. Knepley   PC_COMPOSITE_ADDITIVE       - local additive combination
116512cd4985SMatthew G. Knepley   PC_COMPOSITE_MULTIPLICATIVE - local multiplicative combination
116612cd4985SMatthew G. Knepley .ve
116712cd4985SMatthew G. Knepley 
116812cd4985SMatthew G. Knepley   Options Database Key:
116912cd4985SMatthew G. Knepley . -pc_asm_local_type [additive,multiplicative] - Sets local solver composition type
117012cd4985SMatthew G. Knepley 
117112cd4985SMatthew G. Knepley   Level: intermediate
117212cd4985SMatthew G. Knepley 
1173f1ee410cSBarry Smith .seealso: PCASMSetType(), PCASMGetType(), PCASMGetLocalType(), PCASM, PCASMType, PCASMSetType(), PCASMGetType(), PCCompositeType
117412cd4985SMatthew G. Knepley @*/
117512cd4985SMatthew G. Knepley PetscErrorCode PCASMSetLocalType(PC pc, PCCompositeType type)
117612cd4985SMatthew G. Knepley {
117712cd4985SMatthew G. Knepley   PetscFunctionBegin;
117812cd4985SMatthew G. Knepley   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
117912cd4985SMatthew G. Knepley   PetscValidLogicalCollectiveEnum(pc, type, 2);
11805f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscTryMethod(pc, "PCASMSetLocalType_C", (PC, PCCompositeType), (pc, type)));
118112cd4985SMatthew G. Knepley   PetscFunctionReturn(0);
118212cd4985SMatthew G. Knepley }
118312cd4985SMatthew G. Knepley 
118412cd4985SMatthew G. Knepley /*@
118512cd4985SMatthew G. Knepley   PCASMGetLocalType - Gets the type of composition used for local problems in the additive Schwarz method.
118612cd4985SMatthew G. Knepley 
1187d083f849SBarry Smith   Logically Collective on pc
118812cd4985SMatthew G. Knepley 
118912cd4985SMatthew G. Knepley   Input Parameter:
119012cd4985SMatthew G. Knepley . pc  - the preconditioner context
119112cd4985SMatthew G. Knepley 
119212cd4985SMatthew G. Knepley   Output Parameter:
119312cd4985SMatthew G. Knepley . type - type of composition, one of
119412cd4985SMatthew G. Knepley .vb
119512cd4985SMatthew G. Knepley   PC_COMPOSITE_ADDITIVE       - local additive combination
119612cd4985SMatthew G. Knepley   PC_COMPOSITE_MULTIPLICATIVE - local multiplicative combination
119712cd4985SMatthew G. Knepley .ve
119812cd4985SMatthew G. Knepley 
119912cd4985SMatthew G. Knepley   Options Database Key:
120012cd4985SMatthew G. Knepley . -pc_asm_local_type [additive,multiplicative] - Sets local solver composition type
120112cd4985SMatthew G. Knepley 
120212cd4985SMatthew G. Knepley   Level: intermediate
120312cd4985SMatthew G. Knepley 
1204f1ee410cSBarry Smith .seealso: PCASMSetType(), PCASMGetType(), PCASMSetLocalType(), PCASMCreate(), PCASMType, PCASMSetType(), PCASMGetType(), PCCompositeType
120512cd4985SMatthew G. Knepley @*/
120612cd4985SMatthew G. Knepley PetscErrorCode PCASMGetLocalType(PC pc, PCCompositeType *type)
120712cd4985SMatthew G. Knepley {
120812cd4985SMatthew G. Knepley   PetscFunctionBegin;
120912cd4985SMatthew G. Knepley   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
121012cd4985SMatthew G. Knepley   PetscValidPointer(type, 2);
12115f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscUseMethod(pc, "PCASMGetLocalType_C", (PC, PCCompositeType *), (pc, type)));
121212cd4985SMatthew G. Knepley   PetscFunctionReturn(0);
121312cd4985SMatthew G. Knepley }
121412cd4985SMatthew G. Knepley 
12156ed231c7SMatthew Knepley /*@
12166ed231c7SMatthew Knepley     PCASMSetSortIndices - Determines whether subdomain indices are sorted.
12176ed231c7SMatthew Knepley 
1218d083f849SBarry Smith     Logically Collective on pc
12196ed231c7SMatthew Knepley 
12206ed231c7SMatthew Knepley     Input Parameters:
12216ed231c7SMatthew Knepley +   pc  - the preconditioner context
12226ed231c7SMatthew Knepley -   doSort - sort the subdomain indices
12236ed231c7SMatthew Knepley 
12246ed231c7SMatthew Knepley     Level: intermediate
12256ed231c7SMatthew Knepley 
12266ed231c7SMatthew Knepley .seealso: PCASMSetLocalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(),
12276ed231c7SMatthew Knepley           PCASMCreateSubdomains2D()
12286ed231c7SMatthew Knepley @*/
12297087cfbeSBarry Smith PetscErrorCode  PCASMSetSortIndices(PC pc,PetscBool doSort)
12306ed231c7SMatthew Knepley {
12316ed231c7SMatthew Knepley   PetscFunctionBegin;
12320700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1233acfcf0e5SJed Brown   PetscValidLogicalCollectiveBool(pc,doSort,2);
12345f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscTryMethod(pc,"PCASMSetSortIndices_C",(PC,PetscBool),(pc,doSort)));
12356ed231c7SMatthew Knepley   PetscFunctionReturn(0);
12366ed231c7SMatthew Knepley }
12376ed231c7SMatthew Knepley 
12384b9ad928SBarry Smith /*@C
12394b9ad928SBarry Smith    PCASMGetSubKSP - Gets the local KSP contexts for all blocks on
12404b9ad928SBarry Smith    this processor.
12414b9ad928SBarry Smith 
1242d083f849SBarry Smith    Collective on pc iff first_local is requested
12434b9ad928SBarry Smith 
12444b9ad928SBarry Smith    Input Parameter:
12454b9ad928SBarry Smith .  pc - the preconditioner context
12464b9ad928SBarry Smith 
12474b9ad928SBarry Smith    Output Parameters:
12480298fd71SBarry Smith +  n_local - the number of blocks on this processor or NULL
12490298fd71SBarry Smith .  first_local - the global number of the first block on this processor or NULL,
12500298fd71SBarry Smith                  all processors must request or all must pass NULL
12514b9ad928SBarry Smith -  ksp - the array of KSP contexts
12524b9ad928SBarry Smith 
12534b9ad928SBarry Smith    Note:
1254d29017ddSJed Brown    After PCASMGetSubKSP() the array of KSPes is not to be freed.
12554b9ad928SBarry Smith 
12564b9ad928SBarry Smith    You must call KSPSetUp() before calling PCASMGetSubKSP().
12574b9ad928SBarry Smith 
1258d29017ddSJed Brown    Fortran note:
12592bf68e3eSBarry 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.
1260d29017ddSJed Brown 
12614b9ad928SBarry Smith    Level: advanced
12624b9ad928SBarry Smith 
12634b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap(),
12644b9ad928SBarry Smith           PCASMCreateSubdomains2D(),
12654b9ad928SBarry Smith @*/
12667087cfbeSBarry Smith PetscErrorCode  PCASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[])
12674b9ad928SBarry Smith {
12684b9ad928SBarry Smith   PetscFunctionBegin;
12690700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
12705f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscUseMethod(pc,"PCASMGetSubKSP_C",(PC,PetscInt*,PetscInt*,KSP **),(pc,n_local,first_local,ksp)));
12714b9ad928SBarry Smith   PetscFunctionReturn(0);
12724b9ad928SBarry Smith }
12734b9ad928SBarry Smith 
12744b9ad928SBarry Smith /* -------------------------------------------------------------------------------------*/
12754b9ad928SBarry Smith /*MC
12763b09bd56SBarry Smith    PCASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with
12774b9ad928SBarry Smith            its own KSP object.
12784b9ad928SBarry Smith 
12794b9ad928SBarry Smith    Options Database Keys:
128049517cdeSBarry Smith +  -pc_asm_blocks <blks> - Sets total blocks
12814b9ad928SBarry Smith .  -pc_asm_overlap <ovl> - Sets overlap
1282f1ee410cSBarry Smith .  -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type, default is restrict
1283f1ee410cSBarry Smith -  -pc_asm_local_type [additive, multiplicative] - Sets ASM type, default is additive
12844b9ad928SBarry Smith 
12853b09bd56SBarry Smith      IMPORTANT: If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you
12863b09bd56SBarry Smith       will get a different convergence rate due to the default option of -pc_asm_type restrict. Use
12873b09bd56SBarry Smith       -pc_asm_type basic to use the standard ASM.
12883b09bd56SBarry Smith 
128995452b02SPatrick Sanan    Notes:
129095452b02SPatrick Sanan     Each processor can have one or more blocks, but a block cannot be shared by more
1291f1ee410cSBarry Smith      than one processor. Use PCGASM for subdomains shared by multiple processes. Defaults to one block per processor.
12924b9ad928SBarry Smith 
12933b09bd56SBarry Smith      To set options on the solvers for each block append -sub_ to all the KSP, and PC
1294d7ee0231SBarry Smith         options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly
12954b9ad928SBarry Smith 
1296a8c7a070SBarry Smith      To set the options on the solvers separate for each block call PCASMGetSubKSP()
12974b9ad928SBarry Smith          and set the options directly on the resulting KSP object (you can access its PC
12984b9ad928SBarry Smith          with KSPGetPC())
12994b9ad928SBarry Smith 
13004b9ad928SBarry Smith    Level: beginner
13014b9ad928SBarry Smith 
1302c582cd25SBarry Smith     References:
1303606c0280SSatish Balay +   * - M Dryja, OB Widlund, An additive variant of the Schwarz alternating method for the case of many subregions
130496a0c994SBarry Smith      Courant Institute, New York University Technical report
1305606c0280SSatish Balay -   * - Barry Smith, Petter Bjorstad, and William Gropp, Domain Decompositions: Parallel Multilevel Methods for Elliptic Partial Differential Equations,
130696a0c994SBarry Smith     Cambridge University Press.
1307c582cd25SBarry Smith 
13084b9ad928SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
1309f1ee410cSBarry Smith            PCBJACOBI, PCASMGetSubKSP(), PCASMSetLocalSubdomains(), PCASMType, PCASMGetType(), PCASMSetLocalType(), PCASMGetLocalType()
131034a84908Sprj-            PCASMSetTotalSubdomains(), PCSetModifySubMatrices(), PCASMSetOverlap(), PCASMSetType(), PCCompositeType
1311e09e08ccSBarry Smith 
13124b9ad928SBarry Smith M*/
13134b9ad928SBarry Smith 
13148cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_ASM(PC pc)
13154b9ad928SBarry Smith {
13164b9ad928SBarry Smith   PC_ASM         *osm;
13174b9ad928SBarry Smith 
13184b9ad928SBarry Smith   PetscFunctionBegin;
13195f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscNewLog(pc,&osm));
13202fa5cd67SKarl Rupp 
13214b9ad928SBarry Smith   osm->n                 = PETSC_DECIDE;
13224b9ad928SBarry Smith   osm->n_local           = 0;
13232b837212SDmitry Karpeev   osm->n_local_true      = PETSC_DECIDE;
13244b9ad928SBarry Smith   osm->overlap           = 1;
13250a545947SLisandro Dalcin   osm->ksp               = NULL;
13260a545947SLisandro Dalcin   osm->restriction       = NULL;
13270a545947SLisandro Dalcin   osm->lprolongation     = NULL;
13280a545947SLisandro Dalcin   osm->lrestriction      = NULL;
13290a545947SLisandro Dalcin   osm->x                 = NULL;
13300a545947SLisandro Dalcin   osm->y                 = NULL;
13310a545947SLisandro Dalcin   osm->is                = NULL;
13320a545947SLisandro Dalcin   osm->is_local          = NULL;
13330a545947SLisandro Dalcin   osm->mat               = NULL;
13340a545947SLisandro Dalcin   osm->pmat              = NULL;
13354b9ad928SBarry Smith   osm->type              = PC_ASM_RESTRICT;
133612cd4985SMatthew G. Knepley   osm->loctype           = PC_COMPOSITE_ADDITIVE;
13376ed231c7SMatthew Knepley   osm->sort_indices      = PETSC_TRUE;
1338d709ea83SDmitry Karpeev   osm->dm_subdomains     = PETSC_FALSE;
133980ec0b7dSPatrick Sanan   osm->sub_mat_type      = NULL;
13404b9ad928SBarry Smith 
134192bb6962SLisandro Dalcin   pc->data                 = (void*)osm;
13424b9ad928SBarry Smith   pc->ops->apply           = PCApply_ASM;
134348e38a8aSPierre Jolivet   pc->ops->matapply        = PCMatApply_ASM;
13444b9ad928SBarry Smith   pc->ops->applytranspose  = PCApplyTranspose_ASM;
13454b9ad928SBarry Smith   pc->ops->setup           = PCSetUp_ASM;
1346e91c6855SBarry Smith   pc->ops->reset           = PCReset_ASM;
13474b9ad928SBarry Smith   pc->ops->destroy         = PCDestroy_ASM;
13484b9ad928SBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_ASM;
13494b9ad928SBarry Smith   pc->ops->setuponblocks   = PCSetUpOnBlocks_ASM;
13504b9ad928SBarry Smith   pc->ops->view            = PCView_ASM;
13510a545947SLisandro Dalcin   pc->ops->applyrichardson = NULL;
13524b9ad928SBarry Smith 
13535f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectComposeFunction((PetscObject)pc,"PCASMSetLocalSubdomains_C",PCASMSetLocalSubdomains_ASM));
13545f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectComposeFunction((PetscObject)pc,"PCASMSetTotalSubdomains_C",PCASMSetTotalSubdomains_ASM));
13555f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectComposeFunction((PetscObject)pc,"PCASMSetOverlap_C",PCASMSetOverlap_ASM));
13565f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectComposeFunction((PetscObject)pc,"PCASMSetType_C",PCASMSetType_ASM));
13575f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectComposeFunction((PetscObject)pc,"PCASMGetType_C",PCASMGetType_ASM));
13585f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectComposeFunction((PetscObject)pc,"PCASMSetLocalType_C",PCASMSetLocalType_ASM));
13595f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectComposeFunction((PetscObject)pc,"PCASMGetLocalType_C",PCASMGetLocalType_ASM));
13605f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectComposeFunction((PetscObject)pc,"PCASMSetSortIndices_C",PCASMSetSortIndices_ASM));
13615f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectComposeFunction((PetscObject)pc,"PCASMGetSubKSP_C",PCASMGetSubKSP_ASM));
13625f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectComposeFunction((PetscObject)pc,"PCASMGetSubMatType_C",PCASMGetSubMatType_ASM));
13635f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectComposeFunction((PetscObject)pc,"PCASMSetSubMatType_C",PCASMSetSubMatType_ASM));
13644b9ad928SBarry Smith   PetscFunctionReturn(0);
13654b9ad928SBarry Smith }
13664b9ad928SBarry Smith 
136792bb6962SLisandro Dalcin /*@C
136892bb6962SLisandro Dalcin    PCASMCreateSubdomains - Creates the index sets for the overlapping Schwarz
13693b8d980fSPierre Jolivet    preconditioner for any problem on a general grid.
137092bb6962SLisandro Dalcin 
137192bb6962SLisandro Dalcin    Collective
137292bb6962SLisandro Dalcin 
137392bb6962SLisandro Dalcin    Input Parameters:
137492bb6962SLisandro Dalcin +  A - The global matrix operator
137592bb6962SLisandro Dalcin -  n - the number of local blocks
137692bb6962SLisandro Dalcin 
137792bb6962SLisandro Dalcin    Output Parameters:
137892bb6962SLisandro Dalcin .  outis - the array of index sets defining the subdomains
137992bb6962SLisandro Dalcin 
138092bb6962SLisandro Dalcin    Level: advanced
138192bb6962SLisandro Dalcin 
13827d6bfa3bSBarry Smith    Note: this generates nonoverlapping subdomains; the PCASM will generate the overlap
13837d6bfa3bSBarry Smith     from these if you use PCASMSetLocalSubdomains()
13847d6bfa3bSBarry Smith 
13857d6bfa3bSBarry Smith     In the Fortran version you must provide the array outis[] already allocated of length n.
13867d6bfa3bSBarry Smith 
138792bb6962SLisandro Dalcin .seealso: PCASMSetLocalSubdomains(), PCASMDestroySubdomains()
138892bb6962SLisandro Dalcin @*/
13897087cfbeSBarry Smith PetscErrorCode  PCASMCreateSubdomains(Mat A, PetscInt n, IS* outis[])
139092bb6962SLisandro Dalcin {
139192bb6962SLisandro Dalcin   MatPartitioning mpart;
139292bb6962SLisandro Dalcin   const char      *prefix;
139392bb6962SLisandro Dalcin   PetscInt        i,j,rstart,rend,bs;
1394976e8c5aSLisandro Dalcin   PetscBool       hasop, isbaij = PETSC_FALSE,foundpart = PETSC_FALSE;
13950298fd71SBarry Smith   Mat             Ad     = NULL, adj;
139692bb6962SLisandro Dalcin   IS              ispart,isnumb,*is;
139792bb6962SLisandro Dalcin 
139892bb6962SLisandro Dalcin   PetscFunctionBegin;
13990700a824SBarry Smith   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
140092bb6962SLisandro Dalcin   PetscValidPointer(outis,3);
14012c71b3e2SJacob Faibussowitsch   PetscCheckFalse(n < 1,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of local blocks must be > 0, n = %D",n);
140292bb6962SLisandro Dalcin 
140392bb6962SLisandro Dalcin   /* Get prefix, row distribution, and block size */
14045f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOptionsPrefix(A,&prefix));
14055f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOwnershipRange(A,&rstart,&rend));
14065f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetBlockSize(A,&bs));
14072c71b3e2SJacob Faibussowitsch   PetscCheckFalse(rstart/bs*bs != rstart || rend/bs*bs != rend,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"bad row distribution [%D,%D) for matrix block size %D",rstart,rend,bs);
140865e19b50SBarry Smith 
140992bb6962SLisandro Dalcin   /* Get diagonal block from matrix if possible */
14105f80ce2aSJacob Faibussowitsch   CHKERRQ(MatHasOperation(A,MATOP_GET_DIAGONAL_BLOCK,&hasop));
1411976e8c5aSLisandro Dalcin   if (hasop) {
14125f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetDiagonalBlock(A,&Ad));
141392bb6962SLisandro Dalcin   }
141492bb6962SLisandro Dalcin   if (Ad) {
14155f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectBaseTypeCompare((PetscObject)Ad,MATSEQBAIJ,&isbaij));
14165f80ce2aSJacob Faibussowitsch     if (!isbaij) CHKERRQ(PetscObjectBaseTypeCompare((PetscObject)Ad,MATSEQSBAIJ,&isbaij));
141792bb6962SLisandro Dalcin   }
141892bb6962SLisandro Dalcin   if (Ad && n > 1) {
1419ace3abfcSBarry Smith     PetscBool match,done;
142092bb6962SLisandro Dalcin     /* Try to setup a good matrix partitioning if available */
14215f80ce2aSJacob Faibussowitsch     CHKERRQ(MatPartitioningCreate(PETSC_COMM_SELF,&mpart));
14225f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix));
14235f80ce2aSJacob Faibussowitsch     CHKERRQ(MatPartitioningSetFromOptions(mpart));
14245f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGCURRENT,&match));
142592bb6962SLisandro Dalcin     if (!match) {
14265f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGSQUARE,&match));
142792bb6962SLisandro Dalcin     }
142892bb6962SLisandro Dalcin     if (!match) { /* assume a "good" partitioner is available */
14291a83f524SJed Brown       PetscInt       na;
14301a83f524SJed Brown       const PetscInt *ia,*ja;
14315f80ce2aSJacob Faibussowitsch       CHKERRQ(MatGetRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done));
143292bb6962SLisandro Dalcin       if (done) {
143392bb6962SLisandro Dalcin         /* Build adjacency matrix by hand. Unfortunately a call to
143492bb6962SLisandro Dalcin            MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will
143592bb6962SLisandro Dalcin            remove the block-aij structure and we cannot expect
143692bb6962SLisandro Dalcin            MatPartitioning to split vertices as we need */
14370a545947SLisandro Dalcin         PetscInt       i,j,len,nnz,cnt,*iia=NULL,*jja=NULL;
14381a83f524SJed Brown         const PetscInt *row;
143992bb6962SLisandro Dalcin         nnz = 0;
144092bb6962SLisandro Dalcin         for (i=0; i<na; i++) { /* count number of nonzeros */
144192bb6962SLisandro Dalcin           len = ia[i+1] - ia[i];
144292bb6962SLisandro Dalcin           row = ja + ia[i];
144392bb6962SLisandro Dalcin           for (j=0; j<len; j++) {
144492bb6962SLisandro Dalcin             if (row[j] == i) { /* don't count diagonal */
144592bb6962SLisandro Dalcin               len--; break;
144692bb6962SLisandro Dalcin             }
144792bb6962SLisandro Dalcin           }
144892bb6962SLisandro Dalcin           nnz += len;
144992bb6962SLisandro Dalcin         }
14505f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscMalloc1(na+1,&iia));
14515f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscMalloc1(nnz,&jja));
145292bb6962SLisandro Dalcin         nnz    = 0;
145392bb6962SLisandro Dalcin         iia[0] = 0;
145492bb6962SLisandro Dalcin         for (i=0; i<na; i++) { /* fill adjacency */
145592bb6962SLisandro Dalcin           cnt = 0;
145692bb6962SLisandro Dalcin           len = ia[i+1] - ia[i];
145792bb6962SLisandro Dalcin           row = ja + ia[i];
145892bb6962SLisandro Dalcin           for (j=0; j<len; j++) {
145992bb6962SLisandro Dalcin             if (row[j] != i) { /* if not diagonal */
146092bb6962SLisandro Dalcin               jja[nnz+cnt++] = row[j];
146192bb6962SLisandro Dalcin             }
146292bb6962SLisandro Dalcin           }
146392bb6962SLisandro Dalcin           nnz     += cnt;
146492bb6962SLisandro Dalcin           iia[i+1] = nnz;
146592bb6962SLisandro Dalcin         }
146692bb6962SLisandro Dalcin         /* Partitioning of the adjacency matrix */
14675f80ce2aSJacob Faibussowitsch         CHKERRQ(MatCreateMPIAdj(PETSC_COMM_SELF,na,na,iia,jja,NULL,&adj));
14685f80ce2aSJacob Faibussowitsch         CHKERRQ(MatPartitioningSetAdjacency(mpart,adj));
14695f80ce2aSJacob Faibussowitsch         CHKERRQ(MatPartitioningSetNParts(mpart,n));
14705f80ce2aSJacob Faibussowitsch         CHKERRQ(MatPartitioningApply(mpart,&ispart));
14715f80ce2aSJacob Faibussowitsch         CHKERRQ(ISPartitioningToNumbering(ispart,&isnumb));
14725f80ce2aSJacob Faibussowitsch         CHKERRQ(MatDestroy(&adj));
147392bb6962SLisandro Dalcin         foundpart = PETSC_TRUE;
147492bb6962SLisandro Dalcin       }
14755f80ce2aSJacob Faibussowitsch       CHKERRQ(MatRestoreRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done));
147692bb6962SLisandro Dalcin     }
14775f80ce2aSJacob Faibussowitsch     CHKERRQ(MatPartitioningDestroy(&mpart));
147892bb6962SLisandro Dalcin   }
147992bb6962SLisandro Dalcin 
14805f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(n,&is));
148192bb6962SLisandro Dalcin   *outis = is;
148292bb6962SLisandro Dalcin 
148392bb6962SLisandro Dalcin   if (!foundpart) {
148492bb6962SLisandro Dalcin 
148592bb6962SLisandro Dalcin     /* Partitioning by contiguous chunks of rows */
148692bb6962SLisandro Dalcin 
148792bb6962SLisandro Dalcin     PetscInt mbs   = (rend-rstart)/bs;
148892bb6962SLisandro Dalcin     PetscInt start = rstart;
148992bb6962SLisandro Dalcin     for (i=0; i<n; i++) {
149092bb6962SLisandro Dalcin       PetscInt count = (mbs/n + ((mbs % n) > i)) * bs;
14915f80ce2aSJacob Faibussowitsch       CHKERRQ(ISCreateStride(PETSC_COMM_SELF,count,start,1,&is[i]));
149292bb6962SLisandro Dalcin       start += count;
149392bb6962SLisandro Dalcin     }
149492bb6962SLisandro Dalcin 
149592bb6962SLisandro Dalcin   } else {
149692bb6962SLisandro Dalcin 
149792bb6962SLisandro Dalcin     /* Partitioning by adjacency of diagonal block  */
149892bb6962SLisandro Dalcin 
149992bb6962SLisandro Dalcin     const PetscInt *numbering;
150092bb6962SLisandro Dalcin     PetscInt       *count,nidx,*indices,*newidx,start=0;
150192bb6962SLisandro Dalcin     /* Get node count in each partition */
15025f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(n,&count));
15035f80ce2aSJacob Faibussowitsch     CHKERRQ(ISPartitioningCount(ispart,n,count));
150492bb6962SLisandro Dalcin     if (isbaij && bs > 1) { /* adjust for the block-aij case */
150592bb6962SLisandro Dalcin       for (i=0; i<n; i++) count[i] *= bs;
150692bb6962SLisandro Dalcin     }
150792bb6962SLisandro Dalcin     /* Build indices from node numbering */
15085f80ce2aSJacob Faibussowitsch     CHKERRQ(ISGetLocalSize(isnumb,&nidx));
15095f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(nidx,&indices));
151092bb6962SLisandro Dalcin     for (i=0; i<nidx; i++) indices[i] = i; /* needs to be initialized */
15115f80ce2aSJacob Faibussowitsch     CHKERRQ(ISGetIndices(isnumb,&numbering));
15125f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSortIntWithPermutation(nidx,numbering,indices));
15135f80ce2aSJacob Faibussowitsch     CHKERRQ(ISRestoreIndices(isnumb,&numbering));
151492bb6962SLisandro Dalcin     if (isbaij && bs > 1) { /* adjust for the block-aij case */
15155f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscMalloc1(nidx*bs,&newidx));
15162fa5cd67SKarl Rupp       for (i=0; i<nidx; i++) {
15172fa5cd67SKarl Rupp         for (j=0; j<bs; j++) newidx[i*bs+j] = indices[i]*bs + j;
15182fa5cd67SKarl Rupp       }
15195f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFree(indices));
152092bb6962SLisandro Dalcin       nidx   *= bs;
152192bb6962SLisandro Dalcin       indices = newidx;
152292bb6962SLisandro Dalcin     }
152392bb6962SLisandro Dalcin     /* Shift to get global indices */
152492bb6962SLisandro Dalcin     for (i=0; i<nidx; i++) indices[i] += rstart;
152592bb6962SLisandro Dalcin 
152692bb6962SLisandro Dalcin     /* Build the index sets for each block */
152792bb6962SLisandro Dalcin     for (i=0; i<n; i++) {
15285f80ce2aSJacob Faibussowitsch       CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,count[i],&indices[start],PETSC_COPY_VALUES,&is[i]));
15295f80ce2aSJacob Faibussowitsch       CHKERRQ(ISSort(is[i]));
153092bb6962SLisandro Dalcin       start += count[i];
153192bb6962SLisandro Dalcin     }
153292bb6962SLisandro Dalcin 
15335f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(count));
15345f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(indices));
15355f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&isnumb));
15365f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&ispart));
153792bb6962SLisandro Dalcin 
153892bb6962SLisandro Dalcin   }
153992bb6962SLisandro Dalcin   PetscFunctionReturn(0);
154092bb6962SLisandro Dalcin }
154192bb6962SLisandro Dalcin 
154292bb6962SLisandro Dalcin /*@C
154392bb6962SLisandro Dalcin    PCASMDestroySubdomains - Destroys the index sets created with
154492bb6962SLisandro Dalcin    PCASMCreateSubdomains(). Should be called after setting subdomains
154592bb6962SLisandro Dalcin    with PCASMSetLocalSubdomains().
154692bb6962SLisandro Dalcin 
154792bb6962SLisandro Dalcin    Collective
154892bb6962SLisandro Dalcin 
154992bb6962SLisandro Dalcin    Input Parameters:
155092bb6962SLisandro Dalcin +  n - the number of index sets
15512b691e39SMatthew Knepley .  is - the array of index sets
15520298fd71SBarry Smith -  is_local - the array of local index sets, can be NULL
155392bb6962SLisandro Dalcin 
155492bb6962SLisandro Dalcin    Level: advanced
155592bb6962SLisandro Dalcin 
155692bb6962SLisandro Dalcin .seealso: PCASMCreateSubdomains(), PCASMSetLocalSubdomains()
155792bb6962SLisandro Dalcin @*/
15587087cfbeSBarry Smith PetscErrorCode  PCASMDestroySubdomains(PetscInt n, IS is[], IS is_local[])
155992bb6962SLisandro Dalcin {
156092bb6962SLisandro Dalcin   PetscInt       i;
15615fd66863SKarl Rupp 
156292bb6962SLisandro Dalcin   PetscFunctionBegin;
1563a35b7b57SDmitry Karpeev   if (n <= 0) PetscFunctionReturn(0);
1564a35b7b57SDmitry Karpeev   if (is) {
156592bb6962SLisandro Dalcin     PetscValidPointer(is,2);
15665f80ce2aSJacob Faibussowitsch     for (i=0; i<n; i++) CHKERRQ(ISDestroy(&is[i]));
15675f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(is));
1568a35b7b57SDmitry Karpeev   }
15692b691e39SMatthew Knepley   if (is_local) {
15702b691e39SMatthew Knepley     PetscValidPointer(is_local,3);
15715f80ce2aSJacob Faibussowitsch     for (i=0; i<n; i++) CHKERRQ(ISDestroy(&is_local[i]));
15725f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(is_local));
15732b691e39SMatthew Knepley   }
157492bb6962SLisandro Dalcin   PetscFunctionReturn(0);
157592bb6962SLisandro Dalcin }
157692bb6962SLisandro Dalcin 
15774b9ad928SBarry Smith /*@
15784b9ad928SBarry Smith    PCASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz
15794b9ad928SBarry Smith    preconditioner for a two-dimensional problem on a regular grid.
15804b9ad928SBarry Smith 
15814b9ad928SBarry Smith    Not Collective
15824b9ad928SBarry Smith 
15834b9ad928SBarry Smith    Input Parameters:
15846b867d5aSJose E. Roman +  m   - the number of mesh points in the x direction
15856b867d5aSJose E. Roman .  n   - the number of mesh points in the y direction
15866b867d5aSJose E. Roman .  M   - the number of subdomains in the x direction
15876b867d5aSJose E. Roman .  N   - the number of subdomains in the y direction
15884b9ad928SBarry Smith .  dof - degrees of freedom per node
15894b9ad928SBarry Smith -  overlap - overlap in mesh lines
15904b9ad928SBarry Smith 
15914b9ad928SBarry Smith    Output Parameters:
15924b9ad928SBarry Smith +  Nsub - the number of subdomains created
15933d03bb48SJed Brown .  is - array of index sets defining overlapping (if overlap > 0) subdomains
15943d03bb48SJed Brown -  is_local - array of index sets defining non-overlapping subdomains
15954b9ad928SBarry Smith 
15964b9ad928SBarry Smith    Note:
15974b9ad928SBarry Smith    Presently PCAMSCreateSubdomains2d() is valid only for sequential
15984b9ad928SBarry Smith    preconditioners.  More general related routines are
15994b9ad928SBarry Smith    PCASMSetTotalSubdomains() and PCASMSetLocalSubdomains().
16004b9ad928SBarry Smith 
16014b9ad928SBarry Smith    Level: advanced
16024b9ad928SBarry Smith 
16034b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetLocalSubdomains(), PCASMGetSubKSP(),
16044b9ad928SBarry Smith           PCASMSetOverlap()
16054b9ad928SBarry Smith @*/
16067087cfbeSBarry Smith PetscErrorCode  PCASMCreateSubdomains2D(PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt dof,PetscInt overlap,PetscInt *Nsub,IS **is,IS **is_local)
16074b9ad928SBarry Smith {
16083d03bb48SJed Brown   PetscInt       i,j,height,width,ystart,xstart,yleft,yright,xleft,xright,loc_outer;
160913f74950SBarry Smith   PetscInt       nidx,*idx,loc,ii,jj,count;
16104b9ad928SBarry Smith 
16114b9ad928SBarry Smith   PetscFunctionBegin;
16122c71b3e2SJacob Faibussowitsch   PetscCheckFalse(dof != 1,PETSC_COMM_SELF,PETSC_ERR_SUP," ");
16134b9ad928SBarry Smith 
16144b9ad928SBarry Smith   *Nsub     = N*M;
16155f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(*Nsub,is));
16165f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(*Nsub,is_local));
16174b9ad928SBarry Smith   ystart    = 0;
16183d03bb48SJed Brown   loc_outer = 0;
16194b9ad928SBarry Smith   for (i=0; i<N; i++) {
16204b9ad928SBarry Smith     height = n/N + ((n % N) > i); /* height of subdomain */
16212c71b3e2SJacob Faibussowitsch     PetscCheckFalse(height < 2,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many N subdomains for mesh dimension n");
16224b9ad928SBarry Smith     yleft  = ystart - overlap; if (yleft < 0) yleft = 0;
16234b9ad928SBarry Smith     yright = ystart + height + overlap; if (yright > n) yright = n;
16244b9ad928SBarry Smith     xstart = 0;
16254b9ad928SBarry Smith     for (j=0; j<M; j++) {
16264b9ad928SBarry Smith       width = m/M + ((m % M) > j); /* width of subdomain */
16272c71b3e2SJacob Faibussowitsch       PetscCheckFalse(width < 2,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many M subdomains for mesh dimension m");
16284b9ad928SBarry Smith       xleft  = xstart - overlap; if (xleft < 0) xleft = 0;
16294b9ad928SBarry Smith       xright = xstart + width + overlap; if (xright > m) xright = m;
16304b9ad928SBarry Smith       nidx   = (xright - xleft)*(yright - yleft);
16315f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscMalloc1(nidx,&idx));
16324b9ad928SBarry Smith       loc    = 0;
16334b9ad928SBarry Smith       for (ii=yleft; ii<yright; ii++) {
16344b9ad928SBarry Smith         count = m*ii + xleft;
16352fa5cd67SKarl Rupp         for (jj=xleft; jj<xright; jj++) idx[loc++] = count++;
16364b9ad928SBarry Smith       }
16375f80ce2aSJacob Faibussowitsch       CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,nidx,idx,PETSC_COPY_VALUES,(*is)+loc_outer));
16383d03bb48SJed Brown       if (overlap == 0) {
16395f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscObjectReference((PetscObject)(*is)[loc_outer]));
16402fa5cd67SKarl Rupp 
16413d03bb48SJed Brown         (*is_local)[loc_outer] = (*is)[loc_outer];
16423d03bb48SJed Brown       } else {
16433d03bb48SJed Brown         for (loc=0,ii=ystart; ii<ystart+height; ii++) {
16443d03bb48SJed Brown           for (jj=xstart; jj<xstart+width; jj++) {
16453d03bb48SJed Brown             idx[loc++] = m*ii + jj;
16463d03bb48SJed Brown           }
16473d03bb48SJed Brown         }
16485f80ce2aSJacob Faibussowitsch         CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,loc,idx,PETSC_COPY_VALUES,*is_local+loc_outer));
16493d03bb48SJed Brown       }
16505f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFree(idx));
16514b9ad928SBarry Smith       xstart += width;
16523d03bb48SJed Brown       loc_outer++;
16534b9ad928SBarry Smith     }
16544b9ad928SBarry Smith     ystart += height;
16554b9ad928SBarry Smith   }
16565f80ce2aSJacob Faibussowitsch   for (i=0; i<*Nsub; i++) CHKERRQ(ISSort((*is)[i]));
16574b9ad928SBarry Smith   PetscFunctionReturn(0);
16584b9ad928SBarry Smith }
16594b9ad928SBarry Smith 
16604b9ad928SBarry Smith /*@C
16614b9ad928SBarry Smith     PCASMGetLocalSubdomains - Gets the local subdomains (for this processor
16624b9ad928SBarry Smith     only) for the additive Schwarz preconditioner.
16634b9ad928SBarry Smith 
1664ad4df100SBarry Smith     Not Collective
16654b9ad928SBarry Smith 
16664b9ad928SBarry Smith     Input Parameter:
16674b9ad928SBarry Smith .   pc - the preconditioner context
16684b9ad928SBarry Smith 
16694b9ad928SBarry Smith     Output Parameters:
1670f17a6784SPierre Jolivet +   n - if requested, the number of subdomains for this processor (default value = 1)
1671f17a6784SPierre Jolivet .   is - if requested, the index sets that define the subdomains for this processor
1672f17a6784SPierre Jolivet -   is_local - if requested, the index sets that define the local part of the subdomains for this processor (can be NULL)
16734b9ad928SBarry Smith 
16744b9ad928SBarry Smith     Notes:
16754b9ad928SBarry Smith     The IS numbering is in the parallel, global numbering of the vector.
16764b9ad928SBarry Smith 
16774b9ad928SBarry Smith     Level: advanced
16784b9ad928SBarry Smith 
16794b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(),
16804b9ad928SBarry Smith           PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubmatrices()
16814b9ad928SBarry Smith @*/
16827087cfbeSBarry Smith PetscErrorCode  PCASMGetLocalSubdomains(PC pc,PetscInt *n,IS *is[],IS *is_local[])
16834b9ad928SBarry Smith {
16842a808120SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
1685ace3abfcSBarry Smith   PetscBool      match;
16864b9ad928SBarry Smith 
16874b9ad928SBarry Smith   PetscFunctionBegin;
16880700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1689f17a6784SPierre Jolivet   if (n) PetscValidIntPointer(n,2);
169092bb6962SLisandro Dalcin   if (is) PetscValidPointer(is,3);
1691f17a6784SPierre Jolivet   if (is_local) PetscValidPointer(is_local,4);
16925f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectTypeCompare((PetscObject)pc,PCASM,&match));
1693*28b400f6SJacob Faibussowitsch   PetscCheck(match,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"PC is not a PCASM");
16944b9ad928SBarry Smith   if (n) *n = osm->n_local_true;
16954b9ad928SBarry Smith   if (is) *is = osm->is;
16962b691e39SMatthew Knepley   if (is_local) *is_local = osm->is_local;
16974b9ad928SBarry Smith   PetscFunctionReturn(0);
16984b9ad928SBarry Smith }
16994b9ad928SBarry Smith 
17004b9ad928SBarry Smith /*@C
17014b9ad928SBarry Smith     PCASMGetLocalSubmatrices - Gets the local submatrices (for this processor
17024b9ad928SBarry Smith     only) for the additive Schwarz preconditioner.
17034b9ad928SBarry Smith 
1704ad4df100SBarry Smith     Not Collective
17054b9ad928SBarry Smith 
17064b9ad928SBarry Smith     Input Parameter:
17074b9ad928SBarry Smith .   pc - the preconditioner context
17084b9ad928SBarry Smith 
17094b9ad928SBarry Smith     Output Parameters:
1710f17a6784SPierre Jolivet +   n - if requested, the number of matrices for this processor (default value = 1)
1711f17a6784SPierre Jolivet -   mat - if requested, the matrices
17124b9ad928SBarry Smith 
17134b9ad928SBarry Smith     Level: advanced
17144b9ad928SBarry Smith 
171595452b02SPatrick Sanan     Notes:
171634a84908Sprj-     Call after PCSetUp() (or KSPSetUp()) but before PCApply() and before PCSetUpOnBlocks())
1717cf739d55SBarry Smith 
171834a84908Sprj-            Usually one would use PCSetModifySubMatrices() to change the submatrices in building the preconditioner.
1719cf739d55SBarry Smith 
17204b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(),
172134a84908Sprj-           PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains(), PCSetModifySubMatrices()
17224b9ad928SBarry Smith @*/
17237087cfbeSBarry Smith PetscErrorCode  PCASMGetLocalSubmatrices(PC pc,PetscInt *n,Mat *mat[])
17244b9ad928SBarry Smith {
17254b9ad928SBarry Smith   PC_ASM         *osm;
1726ace3abfcSBarry Smith   PetscBool      match;
17274b9ad928SBarry Smith 
17284b9ad928SBarry Smith   PetscFunctionBegin;
17290700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1730f17a6784SPierre Jolivet   if (n) PetscValidIntPointer(n,2);
173192bb6962SLisandro Dalcin   if (mat) PetscValidPointer(mat,3);
1732*28b400f6SJacob Faibussowitsch   PetscCheck(pc->setupcalled,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call after KSPSetUp() or PCSetUp().");
17335f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectTypeCompare((PetscObject)pc,PCASM,&match));
173492bb6962SLisandro Dalcin   if (!match) {
173592bb6962SLisandro Dalcin     if (n) *n = 0;
17360298fd71SBarry Smith     if (mat) *mat = NULL;
173792bb6962SLisandro Dalcin   } else {
17384b9ad928SBarry Smith     osm = (PC_ASM*)pc->data;
17394b9ad928SBarry Smith     if (n) *n = osm->n_local_true;
17404b9ad928SBarry Smith     if (mat) *mat = osm->pmat;
174192bb6962SLisandro Dalcin   }
17424b9ad928SBarry Smith   PetscFunctionReturn(0);
17434b9ad928SBarry Smith }
1744d709ea83SDmitry Karpeev 
1745d709ea83SDmitry Karpeev /*@
1746d709ea83SDmitry Karpeev     PCASMSetDMSubdomains - Indicates whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible.
1747f1ee410cSBarry Smith 
1748d709ea83SDmitry Karpeev     Logically Collective
1749d709ea83SDmitry Karpeev 
1750d8d19677SJose E. Roman     Input Parameters:
1751d709ea83SDmitry Karpeev +   pc  - the preconditioner
1752d709ea83SDmitry Karpeev -   flg - boolean indicating whether to use subdomains defined by the DM
1753d709ea83SDmitry Karpeev 
1754d709ea83SDmitry Karpeev     Options Database Key:
1755147403d9SBarry Smith .   -pc_asm_dm_subdomains <bool> - use subdomains defined by the DM
1756d709ea83SDmitry Karpeev 
1757d709ea83SDmitry Karpeev     Level: intermediate
1758d709ea83SDmitry Karpeev 
1759d709ea83SDmitry Karpeev     Notes:
1760d709ea83SDmitry Karpeev     PCASMSetTotalSubdomains() and PCASMSetOverlap() take precedence over PCASMSetDMSubdomains(),
1761d709ea83SDmitry Karpeev     so setting either of the first two effectively turns the latter off.
1762d709ea83SDmitry Karpeev 
1763d709ea83SDmitry Karpeev .seealso: PCASMGetDMSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap()
1764d709ea83SDmitry Karpeev           PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains()
1765d709ea83SDmitry Karpeev @*/
1766d709ea83SDmitry Karpeev PetscErrorCode  PCASMSetDMSubdomains(PC pc,PetscBool flg)
1767d709ea83SDmitry Karpeev {
1768d709ea83SDmitry Karpeev   PC_ASM         *osm = (PC_ASM*)pc->data;
1769d709ea83SDmitry Karpeev   PetscBool      match;
1770d709ea83SDmitry Karpeev 
1771d709ea83SDmitry Karpeev   PetscFunctionBegin;
1772d709ea83SDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1773d709ea83SDmitry Karpeev   PetscValidLogicalCollectiveBool(pc,flg,2);
1774*28b400f6SJacob Faibussowitsch   PetscCheck(!pc->setupcalled,((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Not for a setup PC.");
17755f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectTypeCompare((PetscObject)pc,PCASM,&match));
1776d709ea83SDmitry Karpeev   if (match) {
1777d709ea83SDmitry Karpeev     osm->dm_subdomains = flg;
1778d709ea83SDmitry Karpeev   }
1779d709ea83SDmitry Karpeev   PetscFunctionReturn(0);
1780d709ea83SDmitry Karpeev }
1781d709ea83SDmitry Karpeev 
1782d709ea83SDmitry Karpeev /*@
1783d709ea83SDmitry Karpeev     PCASMGetDMSubdomains - Returns flag indicating whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible.
1784d709ea83SDmitry Karpeev     Not Collective
1785d709ea83SDmitry Karpeev 
1786d709ea83SDmitry Karpeev     Input Parameter:
1787d709ea83SDmitry Karpeev .   pc  - the preconditioner
1788d709ea83SDmitry Karpeev 
1789d709ea83SDmitry Karpeev     Output Parameter:
1790d709ea83SDmitry Karpeev .   flg - boolean indicating whether to use subdomains defined by the DM
1791d709ea83SDmitry Karpeev 
1792d709ea83SDmitry Karpeev     Level: intermediate
1793d709ea83SDmitry Karpeev 
1794d709ea83SDmitry Karpeev .seealso: PCASMSetDMSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap()
1795d709ea83SDmitry Karpeev           PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains()
1796d709ea83SDmitry Karpeev @*/
1797d709ea83SDmitry Karpeev PetscErrorCode  PCASMGetDMSubdomains(PC pc,PetscBool* flg)
1798d709ea83SDmitry Karpeev {
1799d709ea83SDmitry Karpeev   PC_ASM         *osm = (PC_ASM*)pc->data;
1800d709ea83SDmitry Karpeev   PetscBool      match;
1801d709ea83SDmitry Karpeev 
1802d709ea83SDmitry Karpeev   PetscFunctionBegin;
1803d709ea83SDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1804534a8f05SLisandro Dalcin   PetscValidBoolPointer(flg,2);
18055f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectTypeCompare((PetscObject)pc,PCASM,&match));
180656ea09ceSStefano Zampini   if (match) *flg = osm->dm_subdomains;
180756ea09ceSStefano Zampini   else *flg = PETSC_FALSE;
1808d709ea83SDmitry Karpeev   PetscFunctionReturn(0);
1809d709ea83SDmitry Karpeev }
181080ec0b7dSPatrick Sanan 
181180ec0b7dSPatrick Sanan /*@
181280ec0b7dSPatrick Sanan      PCASMGetSubMatType - Gets the matrix type used for ASM subsolves, as a string.
181380ec0b7dSPatrick Sanan 
181480ec0b7dSPatrick Sanan    Not Collective
181580ec0b7dSPatrick Sanan 
181680ec0b7dSPatrick Sanan    Input Parameter:
181780ec0b7dSPatrick Sanan .  pc - the PC
181880ec0b7dSPatrick Sanan 
181980ec0b7dSPatrick Sanan    Output Parameter:
1820f1ee410cSBarry Smith .  -pc_asm_sub_mat_type - name of matrix type
182180ec0b7dSPatrick Sanan 
182280ec0b7dSPatrick Sanan    Level: advanced
182380ec0b7dSPatrick Sanan 
182480ec0b7dSPatrick Sanan .seealso: PCASMSetSubMatType(), PCASM, PCSetType(), VecSetType(), MatType, Mat
182580ec0b7dSPatrick Sanan @*/
182656ea09ceSStefano Zampini PetscErrorCode  PCASMGetSubMatType(PC pc,MatType *sub_mat_type)
182756ea09ceSStefano Zampini {
182856ea09ceSStefano Zampini   PetscFunctionBegin;
182956ea09ceSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
18305f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscTryMethod(pc,"PCASMGetSubMatType_C",(PC,MatType*),(pc,sub_mat_type)));
183180ec0b7dSPatrick Sanan   PetscFunctionReturn(0);
183280ec0b7dSPatrick Sanan }
183380ec0b7dSPatrick Sanan 
183480ec0b7dSPatrick Sanan /*@
183580ec0b7dSPatrick Sanan      PCASMSetSubMatType - Set the type of matrix used for ASM subsolves
183680ec0b7dSPatrick Sanan 
183780ec0b7dSPatrick Sanan    Collective on Mat
183880ec0b7dSPatrick Sanan 
183980ec0b7dSPatrick Sanan    Input Parameters:
184080ec0b7dSPatrick Sanan +  pc             - the PC object
184180ec0b7dSPatrick Sanan -  sub_mat_type   - matrix type
184280ec0b7dSPatrick Sanan 
184380ec0b7dSPatrick Sanan    Options Database Key:
184480ec0b7dSPatrick 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.
184580ec0b7dSPatrick Sanan 
184680ec0b7dSPatrick Sanan    Notes:
184780ec0b7dSPatrick Sanan    See "${PETSC_DIR}/include/petscmat.h" for available types
184880ec0b7dSPatrick Sanan 
184980ec0b7dSPatrick Sanan   Level: advanced
185080ec0b7dSPatrick Sanan 
185180ec0b7dSPatrick Sanan .seealso: PCASMGetSubMatType(), PCASM, PCSetType(), VecSetType(), MatType, Mat
185280ec0b7dSPatrick Sanan @*/
185380ec0b7dSPatrick Sanan PetscErrorCode PCASMSetSubMatType(PC pc,MatType sub_mat_type)
185480ec0b7dSPatrick Sanan {
185556ea09ceSStefano Zampini   PetscFunctionBegin;
185656ea09ceSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
18575f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscTryMethod(pc,"PCASMSetSubMatType_C",(PC,MatType),(pc,sub_mat_type)));
185880ec0b7dSPatrick Sanan   PetscFunctionReturn(0);
185980ec0b7dSPatrick Sanan }
1860